21 #ifdef DEAL_II_WITH_SUNDIALS 26 # ifdef DEAL_II_WITH_TRILINOS 30 # ifdef DEAL_II_WITH_PETSC 37 # include <sundials/sundials_config.h> 38 # if DEAL_II_SUNDIALS_VERSION_GTE(3, 0, 0) 39 # include <kinsol/kinsol_direct.h> 40 # include <sunlinsol/sunlinsol_dense.h> 41 # include <sunmatrix/sunmatrix_dense.h> 43 # include <kinsol/kinsol_dense.h> 57 template <
typename VectorType>
59 t_kinsol_function(N_Vector yy, N_Vector FF,
void *user_data)
61 KINSOL<VectorType> &solver =
62 *
static_cast<KINSOL<VectorType> *
>(user_data);
66 solver.reinit_vector(*src_yy);
69 solver.reinit_vector(*dst_FF);
75 err = solver.residual(*src_yy, *dst_FF);
76 else if (solver.iteration_function)
77 err = solver.iteration_function(*src_yy, *dst_FF);
88 template <
typename VectorType>
90 t_kinsol_setup_jacobian(KINMem kinsol_mem)
92 KINSOL<VectorType> &solver =
93 *
static_cast<KINSOL<VectorType> *
>(kinsol_mem->kin_user_data);
97 solver.reinit_vector(*src_ycur);
100 solver.reinit_vector(*src_fcur);
102 copy(*src_ycur, kinsol_mem->kin_uu);
103 copy(*src_fcur, kinsol_mem->kin_fval);
105 int err = solver.setup_jacobian(*src_ycur, *src_fcur);
111 template <
typename VectorType>
113 t_kinsol_solve_jacobian(KINMem kinsol_mem,
119 KINSOL<VectorType> &solver =
120 *
static_cast<KINSOL<VectorType> *
>(kinsol_mem->kin_user_data);
124 solver.reinit_vector(*src_ycur);
127 solver.reinit_vector(*src_fcur);
129 copy(*src_ycur, kinsol_mem->kin_uu);
130 copy(*src_fcur, kinsol_mem->kin_fval);
133 solver.reinit_vector(*src);
136 solver.reinit_vector(*dst);
140 int err = solver.solve_jacobian_system(*src_ycur, *src_fcur, *src, *dst);
143 *sJpnorm = N_VWL2Norm(b, kinsol_mem->kin_fscale);
144 N_VProd(b, kinsol_mem->kin_fscale, b);
145 N_VProd(b, kinsol_mem->kin_fscale, b);
146 *sFdotJp = N_VDotProd(kinsol_mem->kin_fval, b);
152 template <
typename VectorType>
154 const MPI_Comm mpi_comm)
156 , kinsol_mem(nullptr)
169 template <
typename VectorType>
173 KINFree(&kinsol_mem);
174 # ifdef DEAL_II_WITH_MPI 186 template <
typename VectorType>
190 unsigned int system_size = initial_guess_and_solution.size();
195 # ifdef DEAL_II_WITH_MPI 198 const IndexSet is = initial_guess_and_solution.locally_owned_elements();
199 const unsigned int local_system_size = is.
n_elements();
202 N_VNew_Parallel(
communicator, local_system_size, system_size);
205 N_VConst_Parallel(1.e0,
u_scale);
208 N_VConst_Parallel(1.e0,
f_scale);
215 "Trying to use a serial code with a parallel vector."));
216 solution = N_VNew_Serial(system_size);
217 u_scale = N_VNew_Serial(system_size);
218 N_VConst_Serial(1.e0,
u_scale);
219 f_scale = N_VNew_Serial(system_size);
220 N_VConst_Serial(1.e0,
f_scale);
232 KINFree(&kinsol_mem);
234 kinsol_mem = KINCreate();
236 int status = KINInit(kinsol_mem, t_kinsol_function<VectorType>,
solution);
240 status = KINSetUserData(kinsol_mem, static_cast<void *>(
this));
270 # if DEAL_II_SUNDIALS_VERSION_GTE(3, 0, 0) 271 SUNMatrix J =
nullptr;
272 SUNLinearSolver LS =
nullptr;
277 auto KIN_mem =
static_cast<KINMem
>(
kinsol_mem);
278 KIN_mem->kin_lsolve = t_kinsol_solve_jacobian<VectorType>;
281 KIN_mem->kin_lsetup = t_kinsol_setup_jacobian<VectorType>;
282 # if DEAL_II_SUNDIALS_VERSION_LT(3, 0, 0) 283 KIN_mem->kin_setupNonNull =
true;
289 # if DEAL_II_SUNDIALS_VERSION_GTE(3, 0, 0) 290 J = SUNDenseMatrix(system_size, system_size);
291 LS = SUNDenseLinearSolver(
u_scale, J);
292 status = KINDlsSetLinearSolver(kinsol_mem, LS, J);
294 status = KINDense(kinsol_mem, system_size);
314 # ifdef DEAL_II_WITH_MPI 330 status = KINGetNumNonlinSolvIters(kinsol_mem, &nniters);
333 # if DEAL_II_SUNDIALS_VERSION_GTE(3, 0, 0) 337 KINFree(&kinsol_mem);
339 return static_cast<unsigned int>(nniters);
342 template <
typename VectorType>
354 # ifdef DEAL_II_WITH_MPI 356 # ifdef DEAL_II_WITH_TRILINOS 361 # ifdef DEAL_II_WITH_PETSC 362 # ifndef PETSC_USE_COMPLEX
#define AssertNothrow(cond, exc)
unsigned int maximum_non_linear_iterations
std::function< int(const VectorType &src, VectorType &dst)> residual
void copy(TrilinosWrappers::MPI::Vector &dst, const N_Vector &src)
unsigned int maximum_setup_calls
double function_tolerance
std::function< int(const VectorType &src, VectorType &dst)> iteration_function
void set_functions_to_trigger_an_assert()
#define AssertThrow(cond, exc)
SolutionStrategy strategy
std::function< VectorType &()> get_function_scaling
unsigned int anderson_subspace_size
double maximum_newton_step
#define Assert(cond, exc)
std::function< VectorType &()> get_solution_scaling
static ::ExceptionBase & ExcFunctionNotProvided(std::string arg1)
#define DEAL_II_NAMESPACE_CLOSE
unsigned int maximum_beta_failures
std::function< int(const VectorType &ycur, const VectorType &fcur, const VectorType &rhs, VectorType &dst)> solve_jacobian_system
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
std::function< void(VectorType &)> reinit_vector
KINSOL(const AdditionalData &data=AdditionalData(), const MPI_Comm mpi_comm=MPI_COMM_WORLD)
#define AssertKINSOL(code)
MPI_Comm duplicate_communicator(const MPI_Comm &mpi_communicator)
#define DEAL_II_NAMESPACE_OPEN
std::function< int(const VectorType ¤t_u, const VectorType ¤t_f)> setup_jacobian
unsigned int solve(VectorType &initial_guess_and_solution)
size_type n_elements() const
void copy(const T *begin, const T *end, U *dest)
static ::ExceptionBase & ExcInternalError()