Skip to content

Commit

Permalink
Issue #3682 Cleanup, lint fixes
Browse files Browse the repository at this point in the history
  • Loading branch information
jaeandersson committed May 8, 2024
1 parent f645433 commit 8939d61
Show file tree
Hide file tree
Showing 4 changed files with 17 additions and 14 deletions.
3 changes: 2 additions & 1 deletion casadi/core/integrator.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -908,7 +908,8 @@ int Integrator::bdae_sp_forward(SpForwardMem* m, const bvec_t* x, const bvec_t*
}

int Integrator::bquad_sp_forward(SpForwardMem* m, const bvec_t* x, const bvec_t* z,
const bvec_t* p, const bvec_t* u, const bvec_t* adj_ode, const bvec_t* adj_alg, const bvec_t* adj_quad,
const bvec_t* p, const bvec_t* u,
const bvec_t* adj_ode, const bvec_t* adj_alg, const bvec_t* adj_quad,
bvec_t* adj_p, bvec_t* adj_u) const {
// Evaluate nondifferentiated
m->arg[BDYN_T] = nullptr; // t
Expand Down
3 changes: 2 additions & 1 deletion casadi/interfaces/sundials/cvodes_interface.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -545,7 +545,8 @@ int CvodesInterface::psolveB(double t, N_Vector x, N_Vector xB, N_Vector xdotB,
if (s.second_order_correction_) {
// The outputs will double as seeds for calc_daeB
casadi_clear(v + s.nrx1_ * s.nadj_, s.nrx_ - s.nrx1_ * s.nadj_);
if (s.calc_daeB(m, t, NV_DATA_S(x), nullptr, v, nullptr, nullptr, m->tmp2, nullptr)) return 1;
if (s.calc_daeB(m, t, NV_DATA_S(x), nullptr, v, nullptr, nullptr, m->tmp2, nullptr))
return 1;
// Subtract m->tmp2 from m->tmp1, scaled with gammaB
casadi_axpy(s.nrx_ - s.nrx1_ * s.nadj_, -m->gammaB, m->tmp2 + s.nrx1_ * s.nadj_,
m->tmp1 + s.nrx1_ * s.nadj_);
Expand Down
21 changes: 11 additions & 10 deletions casadi/interfaces/sundials/sundials_interface.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -595,7 +595,8 @@ int SundialsInterface::calc_daeF(SundialsMemory* m, double t, const double* x, c
}

int SundialsInterface::calc_daeB(SundialsMemory* m, double t, const double* x, const double* z,
const double* rx, const double* rz, const double* adj_quad, double* adj_x, double* adj_z) const {
const double* adj_ode, const double* adj_alg, const double* adj_quad,
double* adj_x, double* adj_z) const {
// Evaluate nondifferentiated
m->arg[BDYN_T] = &t; // t
m->arg[BDYN_X] = x; // x
Expand All @@ -606,8 +607,8 @@ int SundialsInterface::calc_daeB(SundialsMemory* m, double t, const double* x, c
m->arg[BDYN_OUT_ALG] = nullptr; // out_alg
m->arg[BDYN_OUT_QUAD] = nullptr; // out_quad
m->arg[BDYN_OUT_ZERO] = nullptr; // out_zero
m->arg[BDYN_ADJ_ODE] = rx; // adj_ode
m->arg[BDYN_ADJ_ALG] = rz; // adj_alg
m->arg[BDYN_ADJ_ODE] = adj_ode; // adj_ode
m->arg[BDYN_ADJ_ALG] = adj_alg; // adj_alg
m->arg[BDYN_ADJ_QUAD] = adj_quad; // adj_quad
m->arg[BDYN_ADJ_ZERO] = nullptr; // adj_zero
m->res[BDAE_ADJ_X] = adj_x; // adj_x
Expand All @@ -627,9 +628,9 @@ int SundialsInterface::calc_daeB(SundialsMemory* m, double t, const double* x, c
m->arg[BDYN_NUM_IN + BDAE_NUM_OUT + BDYN_OUT_QUAD] = nullptr; // fwd:out_quad
m->arg[BDYN_NUM_IN + BDAE_NUM_OUT + BDYN_OUT_ZERO] = nullptr; // fwd:out_zero
m->arg[BDYN_NUM_IN + BDAE_NUM_OUT + BDYN_ADJ_ODE] =
rx ? rx + nrx1_ * nadj_ : 0; // fwd:adj_ode
adj_ode ? adj_ode + nrx1_ * nadj_ : 0; // fwd:adj_ode
m->arg[BDYN_NUM_IN + BDAE_NUM_OUT + BDYN_ADJ_ALG] =
rz ? rz + nrz1_ * nadj_ : 0; // fwd:adj_alg
adj_alg ? adj_alg + nrz1_ * nadj_ : 0; // fwd:adj_alg
m->arg[BDYN_NUM_IN + BDAE_NUM_OUT + BDYN_ADJ_QUAD] =
adj_quad ? adj_quad + nrp1_ * nadj_ : 0; // fwd:adj_quad
m->arg[BDYN_NUM_IN + BDAE_NUM_OUT + BDYN_ADJ_ZERO] = nullptr; // fwd:adj_zero
Expand Down Expand Up @@ -664,7 +665,7 @@ int SundialsInterface::calc_quadF(SundialsMemory* m, double t, const double* x,
}

int SundialsInterface::calc_quadB(SundialsMemory* m, double t, const double* x, const double* z,
const double* rx, const double* rz, double* adj_p, double* adj_u) const {
const double* adj_ode, const double* adj_alg, double* adj_p, double* adj_u) const {
// Evaluate nondifferentiated
m->arg[BDYN_T] = &t; // t
m->arg[BDYN_X] = x; // x
Expand All @@ -675,8 +676,8 @@ int SundialsInterface::calc_quadB(SundialsMemory* m, double t, const double* x,
m->arg[BDYN_OUT_ALG] = nullptr; // out_alg
m->arg[BDYN_OUT_QUAD] = nullptr; // out_quad
m->arg[BDYN_OUT_ZERO] = nullptr; // out_zero
m->arg[BDYN_ADJ_ODE] = rx; // adj_ode
m->arg[BDYN_ADJ_ALG] = rz; // adj_alg
m->arg[BDYN_ADJ_ODE] = adj_ode; // adj_ode
m->arg[BDYN_ADJ_ALG] = adj_alg; // adj_alg
m->arg[BDYN_ADJ_QUAD] = m->adj_q; // adj_quad
m->arg[BDYN_ADJ_ZERO] = nullptr; // adj_zero
m->res[BQUAD_ADJ_P] = adj_p; // adj_p
Expand All @@ -696,9 +697,9 @@ int SundialsInterface::calc_quadB(SundialsMemory* m, double t, const double* x,
m->arg[BDYN_NUM_IN + BQUAD_NUM_OUT + BDYN_OUT_QUAD] = nullptr; // fwd:out_quad
m->arg[BDYN_NUM_IN + BQUAD_NUM_OUT + BDYN_OUT_ZERO] = nullptr; // fwd:out_zero
m->arg[BDYN_NUM_IN + BQUAD_NUM_OUT + BDYN_ADJ_ODE] =
rx ? rx + nrx1_ * nadj_ : 0; // fwd:adj_ode
adj_ode ? adj_ode + nrx1_ * nadj_ : 0; // fwd:adj_ode
m->arg[BDYN_NUM_IN + BQUAD_NUM_OUT + BDYN_ADJ_ALG] =
rz ? rz + nrz1_ * nadj_ : 0; // fwd:adj_alg
adj_alg ? adj_alg + nrz1_ * nadj_ : 0; // fwd:adj_alg
m->arg[BDYN_NUM_IN + BQUAD_NUM_OUT + BDYN_ADJ_QUAD] = m->adj_q + nrp1_ * nadj_; // fwd:adj_quad
m->arg[BDYN_NUM_IN + BQUAD_NUM_OUT + BDYN_ADJ_ZERO] = nullptr; // fwd:adj_zero
m->res[BQUAD_ADJ_P] = adj_p + nrq1_ * nadj_; // fwd:adj_p
Expand Down
4 changes: 2 additions & 2 deletions casadi/interfaces/sundials/sundials_interface.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -125,7 +125,7 @@ namespace casadi {

// DAE right-hand-side, forward problem
int calc_daeB(SundialsMemory* m, double t, const double* x, const double* z,
const double* rx, const double* rz, const double* adj_quad,
const double* adj_ode, const double* adj_alg, const double* adj_quad,
double* adj_x, double* adj_z) const;

// Quadrature right-hand-side, forward problem
Expand All @@ -134,7 +134,7 @@ namespace casadi {

// Quadrature right-hand-side, backward problem
int calc_quadB(SundialsMemory* m, double t, const double* x, const double* z,
const double* rx, const double* rz, double* adj_p, double* adj_u) const;
const double* adj_ode, const double* adj_alg, double* adj_p, double* adj_u) const;

// Jacobian of DAE-times-vector function, forward problem
int calc_jtimesF(SundialsMemory* m, double t, const double* x, const double* z,
Expand Down

0 comments on commit 8939d61

Please sign in to comment.