21#ifdef DEAL_II_WITH_SUNDIALS
29# ifdef DEAL_II_WITH_TRILINOS
33# ifdef DEAL_II_WITH_PETSC
48 <<
"One of the SUNDIALS linear solver internal"
49 <<
" functions returned a negative error code: " << arg1
50 <<
". Please consult SUNDIALS manual.");
52# define AssertSundialsSolver(code) \
53 Assert(code >= 0, ExcSundialsSolverError(code))
60 template <
typename VectorType>
78# if DEAL_II_SUNDIALS_VERSION_GTE(6, 0, 0)
97 template <
typename VectorType>
109 return SUNLINEARSOLVER_ITERATIVE;
114 template <
typename VectorType>
122 auto content = access_content<VectorType>(LS);
124 auto *src_b = unwrap_nvector_const<VectorType>(b);
125 auto *dst_x = unwrap_nvector<VectorType>(x);
137 content->preconditioner_solve,
143 return content->lsolve(op, preconditioner, *dst_x, *src_b, tol);
148 template <
typename VectorType>
152 auto content = access_content<VectorType>(LS);
153 if (content->preconditioner_setup)
154 return content->preconditioner_setup(content->P_data);
160 template <
typename VectorType>
170 template <
typename VectorType>
174 auto content = access_content<VectorType>(LS);
175 content->A_data = A_data;
176 content->a_times_fn = ATimes;
182 template <
typename VectorType>
189 auto content = access_content<VectorType>(LS);
190 content->P_data = P_data;
191 content->preconditioner_setup = p_setup;
192 content->preconditioner_solve = p_solve;
199 template <
typename VectorType>
204 SUNContext linsol_ctx
209# if DEAL_II_SUNDIALS_VERSION_GTE(6, 0, 0)
221 arkode_linsol_set_preconditioner<VectorType>;
224# if DEAL_II_SUNDIALS_VERSION_GTE(6, 0, 0)
234 template <
typename VectorType>
237 SUNLinSolFreeEmpty(sun_linear_solver);
242 template <
typename VectorType>
245 return sun_linear_solver;
251# if DEAL_II_SUNDIALS_VERSION_GTE(6, 0, 0)
252 template <
typename VectorType>
255 SUNContext linsol_ctx)
257 , a_times_fn(a_times_fn)
258 , linsol_ctx(linsol_ctx)
265 template <
typename VectorType>
269 , a_times_fn(a_times_fn)
278 template <
typename VectorType>
281 const VectorType &src)
const
283 auto sun_dst = internal::make_nvector_view(dst
289 auto sun_src = internal::make_nvector_view(src
295 int status = a_times_fn(A_data, sun_src, sun_dst);
302# if DEAL_II_SUNDIALS_VERSION_GTE(6, 0, 0)
303 template <
typename VectorType>
307 SUNContext linsol_ctx,
310 , p_solve_fn(p_solve_fn)
311 , linsol_ctx(linsol_ctx)
315 template <
typename VectorType>
321 , p_solve_fn(p_solve_fn)
328 template <
typename VectorType>
331 const VectorType &src)
const
340 auto sun_dst = internal::make_nvector_view(dst
346 auto sun_src = internal::make_nvector_view(src
355 p_solve_fn(P_data, sun_src, sun_dst, tol, 0 );
380# ifdef DEAL_II_WITH_MPI
381# ifdef DEAL_II_WITH_TRILINOS
393# ifdef DEAL_II_WITH_PETSC
394# ifndef PETSC_USE_COMPLEX
std::unique_ptr< LinearSolverContent< VectorType > > content
SUNLinearSolver sun_linear_solver
LinearSolverWrapper(LinearSolveFunction< VectorType > lsolve, SUNContext linsol_ctx)
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_SUNDIALS_VERSION_GTE(major, minor, patch)
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
static ::ExceptionBase & ExcInternalError()
#define DeclException1(Exception1, type1, outsequence)
static ::ExceptionBase & ExcSundialsSolverError(int arg1)
std::function< int(SundialsOperator< VectorType > &op, SundialsPreconditioner< VectorType > &prec, VectorType &x, const VectorType &b, double tol)> LinearSolveFunction
void vmult(VectorType &dst, const VectorType &src) const
SundialsOperator(void *A_data, ATimesFn a_times_fn, SUNContext linsol_ctx)
SundialsPreconditioner(void *P_data, PSolveFn p_solve_fn, SUNContext linsol_ctx, double tol)
void vmult(VectorType &dst, const VectorType &src) const
LinearSolveFunction< VectorType > lsolve
PSolveFn preconditioner_solve
PSetupFn preconditioner_setup
#define AssertSundialsSolver(code)