--- /dev/null
+Deprecated: ARKode does no longer need a `reinit_vector` function.
+<br>
+(Sebastian Proell, 2021/02/01)
* $z_j$ where $j < i$). Additional information on the specific predictor
* algorithms implemented in ARKode is provided in ARKode documentation.
*
- * The user has to provide the implementation of the following
- *`std::function`s:
- * - reinit_vector()
- * and either one or both of
+ * The user has to provide the implementation of at least one (or both) of the
+ * following `std::function`s:
* - implicit_function()
* - explicit_function()
*
*
* SUNDIALS::ARKode<VectorType> ode;
*
- * ode.reinit_vector = [] (VectorType&v)
- * {
- * v.reinit(2);
- * };
- *
- * double kappa = 1.0;
+ * const double kappa = 1.0;
*
* ode.explicit_function = [kappa] (double,
* const VectorType &y,
get_arkode_memory() const;
/**
- * A function object that users need to supply and that is intended to
- * reinit the given vector.
+ * A function object that was used to `reinit` the given vector. Setting
+ * this field does no longer have any effect and all auxiliary vectors are
+ * reinit-ed automatically based on the user-supplied vector in solve_ode().
+ *
+ * @deprecated This function is no longer used and can be safely removed in
+ * user code.
*/
+ DEAL_II_DEPRECATED
std::function<void(VectorType &)> reinit_vector;
/**
/**
* Constructor.
*
- * @param solver The ARKode solver that uses this operator
* @param A_data Data required by @p a_times_fn
* @param a_times_fn A function pointer to the function that computes A*v
*/
- SundialsOperator(ARKode<VectorType> &solver,
- void * A_data,
- ATimesFn a_times_fn);
+ SundialsOperator(void *A_data, ATimesFn a_times_fn);
private:
- /**
- * Reference to the ARKode object that uses this SundialsOperator.
- */
- ARKode<VectorType> &solver;
/**
* Data necessary to evaluate a_times_fn.
*/
/**
* Constructor.
*
- * @param solver The ARKode solver that uses this operator
* @param P_data Data required by @p p_solve_fn
* @param p_solve_fn A function pointer to the function that computes A*v
* @param tol Tolerance, that an iterative solver should use to judge
* convergence
*/
- SundialsPreconditioner(ARKode<VectorType> &solver,
- void * P_data,
- PSolveFn p_solve_fn,
- double tol);
+ SundialsPreconditioner(void *P_data, PSolveFn p_solve_fn, double tol);
private:
- /**
- * Reference to the ARKode object that uses this SundialsPreconditioner.
- */
- ARKode<VectorType> &solver;
-
/**
* Data necessary to calls p_solve_fn
*/
LinearSolveFunction<VectorType> lsolve;
- //! solver reference to forward all solver related calls to user-specified
- //! functions
- ARKode<VectorType> *solver;
- void * P_data;
- void * A_data;
+ void *P_data;
+ void *A_data;
};
N_Vector b,
realtype tol)
{
- auto content = access_content<VectorType>(LS);
- ARKode<VectorType> &solver = *(content->solver);
+ auto content = access_content<VectorType>(LS);
auto *src_b = internal::unwrap_nvector_const<VectorType>(b);
auto *dst_x = internal::unwrap_nvector<VectorType>(x);
- SundialsOperator<VectorType> op(solver,
- content->A_data,
- content->a_times_fn);
+ SundialsOperator<VectorType> op(content->A_data, content->a_times_fn);
SundialsPreconditioner<VectorType> preconditioner(
- solver, content->P_data, content->preconditioner_solve, tol);
+ content->P_data, content->preconditioner_solve, tol);
return content->lsolve(op, preconditioner, *dst_x, *src_b, tol);
}
class SundialsLinearSolverWrapper
{
public:
- SundialsLinearSolverWrapper(ARKode<VectorType> & solver,
- LinearSolveFunction<VectorType> lsolve)
+ SundialsLinearSolverWrapper(LinearSolveFunction<VectorType> lsolve)
{
sun_linear_solver = SUNLinSolNewEmpty();
sun_linear_solver->ops->gettype = arkode_linsol_get_type;
sun_linear_solver->ops->setpreconditioner =
arkode_linsol_set_preconditioner<VectorType>;
- content.solver = &solver;
content.lsolve = lsolve;
sun_linear_solver->content = &content;
}
{
linear_solver =
std::make_unique<SundialsLinearSolverWrapper<VectorType>>(
- *this, solve_linearized_system);
+ solve_linearized_system);
sun_linear_solver = linear_solver->get_wrapped_solver();
}
else
{
mass_solver =
std::make_unique<SundialsLinearSolverWrapper<VectorType>>(
- *this, solve_mass);
+ solve_mass);
sun_mass_linear_solver = mass_solver->get_wrapped_solver();
}
else
# if DEAL_II_SUNDIALS_VERSION_GTE(4, 0, 0)
template <typename VectorType>
- SundialsOperator<VectorType>::SundialsOperator(ARKode<VectorType> &solver,
- void * A_data,
- ATimesFn a_times_fn)
- : solver(solver)
- , A_data(A_data)
+ SundialsOperator<VectorType>::SundialsOperator(void * A_data,
+ ATimesFn a_times_fn)
+ : A_data(A_data)
, a_times_fn(a_times_fn)
{
template <typename VectorType>
SundialsPreconditioner<VectorType>::SundialsPreconditioner(
- ARKode<VectorType> &solver,
- void * P_data,
- PSolveFn p_solve_fn,
- double tol)
- : solver(solver)
- , P_data(P_data)
+ void * P_data,
+ PSolveFn p_solve_fn,
+ double tol)
+ : P_data(P_data)
, p_solve_fn(p_solve_fn)
, tol(tol)
{}
SUNDIALS::ARKode<VectorType> ode(data);
- ode.reinit_vector = [&](VectorType &v) { v.reinit(2); };
-
double kappa = 1.0;
ode.explicit_function =
SUNDIALS::ARKode<VectorType> ode(data);
- ode.reinit_vector = [&](VectorType &v) { v.reinit(2); };
-
double kappa = 1.0;
ode.implicit_function =
SUNDIALS::ARKode<VectorType> ode(data);
- ode.reinit_vector = [&](VectorType &v) {
- // Three independent variables
- v.reinit(3);
- };
-
// Parameters
double u0 = 3.9, v0 = 1.1, w0 = 2.8, a = 1.2, b = 2.5, eps = 1e-5;
SUNDIALS::ARKode<VectorType> ode(data);
- ode.reinit_vector = [&](VectorType &v) {
- // Three independent variables
- v.reinit(3);
- };
-
// Parameters
double u0 = 3.9, v0 = 1.1, w0 = 2.8, a = 1.2, b = 2.5, eps = 1e-5;
// Explicit jacobian.
SUNDIALS::ARKode<VectorType> ode(data);
- ode.reinit_vector = [&](VectorType &v) {
- // Three independent variables
- v.reinit(3);
- };
-
// Parameters
double u0 = 3.9, v0 = 1.1, w0 = 2.8, a = 1.2, b = 2.5, eps = 1e-5;
// Explicit jacobian.
SUNDIALS::ARKode<VectorType> ode(data);
- ode.reinit_vector = [&](VectorType &v) {
- // Three independent variables
- v.reinit(3);
- };
-
// Parameters
double u0 = 3.9, v0 = 1.1, w0 = 2.8, a = 1.2, b = 2.5, eps = 1e-5;
// Explicit jacobian.
SUNDIALS::ARKode<VectorType> ode(data);
- ode.reinit_vector = [&](VectorType &v) {
- // Three independent variables
- v.reinit(3);
- };
-
// Parameters
double u0 = 3.9, v0 = 1.1, w0 = 2.8, a = 1.2, b = 2.5, eps = 1e-5;
// Explicit jacobian.
SUNDIALS::ARKode<VectorType> ode(data);
- ode.reinit_vector = [&](VectorType &v) { v.reinit(2); };
-
// Explicit jacobian = stiffness matrix
FullMatrix<double> K(2, 2);
K(0, 0) = K(1, 1) = 0.5;