Reference documentation for deal.II version 9.4.1
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
sunlinsol_wrapper.cc
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2021 - 2022 by the deal.II authors
4//
5// This file is part of the deal.II library.
6//
7// The deal.II library is free software; you can use it, redistribute
8// it, and/or modify it under the terms of the GNU Lesser General
9// Public License as published by the Free Software Foundation; either
10// version 2.1 of the License, or (at your option) any later version.
11// The full text of the license can be found in the file LICENSE.md at
12// the top level directory of deal.II.
13//
14// ---------------------------------------------------------------------
15
16
17#include <deal.II/base/config.h>
18
20
21#ifdef DEAL_II_WITH_SUNDIALS
22
24
28# include <deal.II/lac/vector.h>
29# ifdef DEAL_II_WITH_TRILINOS
32# endif
33# ifdef DEAL_II_WITH_PETSC
36# endif
37
39
41
42
43
44namespace SUNDIALS
45{
47 int,
48 << "One of the SUNDIALS linear solver internal"
49 << " functions returned a negative error code: " << arg1
50 << ". Please consult SUNDIALS manual.");
51
52# define AssertSundialsSolver(code) \
53 Assert(code >= 0, ExcSundialsSolverError(code))
54
55 namespace internal
56 {
60 template <typename VectorType>
62 {
64 : a_times_fn(nullptr)
65 , preconditioner_setup(nullptr)
66 , preconditioner_solve(nullptr)
68 , linsol_ctx(nullptr)
69# endif
70 , P_data(nullptr)
71 , A_data(nullptr)
72 {}
73
77
78# if DEAL_II_SUNDIALS_VERSION_GTE(6, 0, 0)
79 SUNContext linsol_ctx;
80# endif
81
83
84 void *P_data;
85 void *A_data;
86 };
87 } // namespace internal
88
89 namespace
90 {
91 using namespace SUNDIALS::internal;
92
97 template <typename VectorType>
99 access_content(SUNLinearSolver ls)
100 {
101 Assert(ls->content != nullptr, ExcInternalError());
102 return static_cast<LinearSolverContent<VectorType> *>(ls->content);
103 }
104
105
106
107 SUNLinearSolver_Type arkode_linsol_get_type(SUNLinearSolver)
108 {
109 return SUNLINEARSOLVER_ITERATIVE;
110 }
111
112
113
114 template <typename VectorType>
115 int
116 arkode_linsol_solve(SUNLinearSolver LS,
117 SUNMatrix /*ignored*/,
118 N_Vector x,
119 N_Vector b,
120 realtype tol)
121 {
122 auto content = access_content<VectorType>(LS);
123
124 auto *src_b = unwrap_nvector_const<VectorType>(b);
125 auto *dst_x = unwrap_nvector<VectorType>(x);
126
127 SundialsOperator<VectorType> op(content->A_data,
128 content->a_times_fn
130 ,
131 content->linsol_ctx
132# endif
133 );
134
136 content->P_data,
137 content->preconditioner_solve,
139 content->linsol_ctx,
140# endif
141 tol);
142
143 return content->lsolve(op, preconditioner, *dst_x, *src_b, tol);
144 }
145
146
147
148 template <typename VectorType>
149 int
150 arkode_linsol_setup(SUNLinearSolver LS, SUNMatrix /*ignored*/)
151 {
152 auto content = access_content<VectorType>(LS);
153 if (content->preconditioner_setup)
154 return content->preconditioner_setup(content->P_data);
155 return 0;
156 }
157
158
159
160 template <typename VectorType>
161 int arkode_linsol_initialize(SUNLinearSolver)
162 {
163 // this method is currently only provided because SUNDIALS 4.0.0 requires
164 // it - no user-set action is implemented so far
165 return 0;
166 }
167
168
169
170 template <typename VectorType>
171 int
172 arkode_linsol_set_a_times(SUNLinearSolver LS, void *A_data, ATimesFn ATimes)
173 {
174 auto content = access_content<VectorType>(LS);
175 content->A_data = A_data;
176 content->a_times_fn = ATimes;
177 return 0;
178 }
179
180
181
182 template <typename VectorType>
183 int
184 arkode_linsol_set_preconditioner(SUNLinearSolver LS,
185 void * P_data,
186 PSetupFn p_setup,
187 PSolveFn p_solve)
188 {
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;
193 return 0;
194 }
195 } // namespace
196
197
198
199 template <typename VectorType>
203 ,
204 SUNContext linsol_ctx
205# endif
206 )
207 : content(std::make_unique<LinearSolverContent<VectorType>>())
208 {
209# if DEAL_II_SUNDIALS_VERSION_GTE(6, 0, 0)
210 sun_linear_solver = SUNLinSolNewEmpty(linsol_ctx);
211# else
212 sun_linear_solver = SUNLinSolNewEmpty();
213# endif
214
215 sun_linear_solver->ops->gettype = arkode_linsol_get_type;
216 sun_linear_solver->ops->solve = arkode_linsol_solve<VectorType>;
217 sun_linear_solver->ops->setup = arkode_linsol_setup<VectorType>;
218 sun_linear_solver->ops->initialize = arkode_linsol_initialize<VectorType>;
219 sun_linear_solver->ops->setatimes = arkode_linsol_set_a_times<VectorType>;
220 sun_linear_solver->ops->setpreconditioner =
221 arkode_linsol_set_preconditioner<VectorType>;
222
223 content->lsolve = lsolve;
224# if DEAL_II_SUNDIALS_VERSION_GTE(6, 0, 0)
225 content->linsol_ctx = linsol_ctx;
226# endif
227 sun_linear_solver->content = content.get();
228 }
229
231
232 namespace internal
233 {
234 template <typename VectorType>
236 {
237 SUNLinSolFreeEmpty(sun_linear_solver);
238 }
239
240
241
242 template <typename VectorType>
244 {
245 return sun_linear_solver;
246 }
247 } // namespace internal
248
249
250
251# if DEAL_II_SUNDIALS_VERSION_GTE(6, 0, 0)
252 template <typename VectorType>
254 ATimesFn a_times_fn,
255 SUNContext linsol_ctx)
256 : A_data(A_data)
257 , a_times_fn(a_times_fn)
258 , linsol_ctx(linsol_ctx)
259
260 {
261 Assert(a_times_fn != nullptr, ExcInternalError());
262 Assert(linsol_ctx != nullptr, ExcInternalError());
263 }
264# else
265 template <typename VectorType>
267 ATimesFn a_times_fn)
268 : A_data(A_data)
269 , a_times_fn(a_times_fn)
270
271 {
272 Assert(a_times_fn != nullptr, ExcInternalError());
273 }
274# endif
275
276
277
278 template <typename VectorType>
279 void
281 const VectorType &src) const
282 {
283 auto sun_dst = internal::make_nvector_view(dst
285 ,
286 linsol_ctx
287# endif
288 );
289 auto sun_src = internal::make_nvector_view(src
291 ,
292 linsol_ctx
293# endif
294 );
295 int status = a_times_fn(A_data, sun_src, sun_dst);
296 (void)status;
297 AssertSundialsSolver(status);
298 }
299
300
301
302# if DEAL_II_SUNDIALS_VERSION_GTE(6, 0, 0)
303 template <typename VectorType>
305 void * P_data,
306 PSolveFn p_solve_fn,
307 SUNContext linsol_ctx,
308 double tol)
309 : P_data(P_data)
310 , p_solve_fn(p_solve_fn)
311 , linsol_ctx(linsol_ctx)
312 , tol(tol)
313 {}
314# else
315 template <typename VectorType>
317 void *P_data,
318 PSolveFn p_solve_fn,
319 double tol)
320 : P_data(P_data)
321 , p_solve_fn(p_solve_fn)
322 , tol(tol)
323 {}
324# endif
325
326
327
328 template <typename VectorType>
329 void
331 const VectorType &src) const
332 {
333 // apply identity preconditioner if nothing else specified
334 if (!p_solve_fn)
335 {
336 dst = src;
337 return;
338 }
339
340 auto sun_dst = internal::make_nvector_view(dst
342 ,
343 linsol_ctx
344# endif
345 );
346 auto sun_src = internal::make_nvector_view(src
348 ,
349 linsol_ctx
350# endif
351 );
352 // for custom preconditioners no distinction between left and right
353 // preconditioning is made
354 int status =
355 p_solve_fn(P_data, sun_src, sun_dst, tol, 0 /*precondition_type*/);
356 (void)status;
357 AssertSundialsSolver(status);
358 }
359
360 template struct SundialsOperator<Vector<double>>;
363 template struct SundialsOperator<
365
368 template struct SundialsPreconditioner<
370 template struct SundialsPreconditioner<
372
375 template class internal::LinearSolverWrapper<
377 template class internal::LinearSolverWrapper<
379
380# ifdef DEAL_II_WITH_MPI
381# ifdef DEAL_II_WITH_TRILINOS
384
387
389 template class internal::LinearSolverWrapper<
391# endif // DEAL_II_WITH_TRILINOS
392
393# ifdef DEAL_II_WITH_PETSC
394# ifndef PETSC_USE_COMPLEX
395
398
401
404# endif // PETSC_USE_COMPLEX
405# endif // DEAL_II_WITH_PETSC
406
407# endif // DEAL_II_WITH_MPI
408} // namespace SUNDIALS
410
411#endif
std::unique_ptr< LinearSolverContent< VectorType > > content
LinearSolverWrapper(LinearSolveFunction< VectorType > lsolve, SUNContext linsol_ctx)
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:442
#define DEAL_II_SUNDIALS_VERSION_GTE(major, minor, patch)
Definition: config.h:347
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:443
#define Assert(cond, exc)
Definition: exceptions.h:1473
static ::ExceptionBase & ExcInternalError()
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:509
static ::ExceptionBase & ExcSundialsSolverError(int arg1)
std::function< int(SundialsOperator< VectorType > &op, SundialsPreconditioner< VectorType > &prec, VectorType &x, const VectorType &b, double tol)> LinearSolveFunction
STL namespace.
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
#define AssertSundialsSolver(code)