17 #include <deal.II/sundials/kinsol.h> 18 #include <deal.II/base/config.h> 20 #ifdef DEAL_II_WITH_SUNDIALS 22 #include <deal.II/base/utilities.h> 23 #include <deal.II/lac/block_vector.h> 24 #ifdef DEAL_II_WITH_TRILINOS 25 # include <deal.II/lac/trilinos_parallel_block_vector.h> 26 # include <deal.II/lac/trilinos_vector.h> 28 #ifdef DEAL_II_WITH_PETSC 29 # include <deal.II/lac/petsc_parallel_block_vector.h> 30 # include <deal.II/lac/petsc_parallel_vector.h> 32 #include <deal.II/base/utilities.h> 33 #include <deal.II/sundials/copy.h> 35 #include <sundials/sundials_config.h> 36 #if DEAL_II_SUNDIALS_VERSION_GTE(3,0,0) 37 # include <sunmatrix/sunmatrix_dense.h> 38 # include <sunlinsol/sunlinsol_dense.h> 39 # include <kinsol/kinsol_direct.h> 41 # include <kinsol/kinsol_dense.h> 47 DEAL_II_NAMESPACE_OPEN
55 template<
typename VectorType>
56 int t_kinsol_function(N_Vector yy,
60 KINSOL<VectorType> &solver = *
static_cast<KINSOL<VectorType> *
>(user_data);
64 solver.reinit_vector(*src_yy);
67 solver.reinit_vector(*dst_FF);
73 err = solver.residual(*src_yy, *dst_FF);
74 else if (solver.iteration_function)
75 err = solver.iteration_function(*src_yy, *dst_FF);
86 template<
typename VectorType>
87 int t_kinsol_setup_jacobian(KINMem kinsol_mem)
89 KINSOL<VectorType> &solver = *
static_cast<KINSOL<VectorType> *
>(kinsol_mem->kin_user_data);
93 solver.reinit_vector(*src_ycur);
96 solver.reinit_vector(*src_fcur);
98 copy(*src_ycur, kinsol_mem->kin_uu);
99 copy(*src_fcur, kinsol_mem->kin_fval);
101 int err = solver.setup_jacobian(*src_ycur, *src_fcur);
107 template<
typename VectorType>
108 int t_kinsol_solve_jacobian(KINMem kinsol_mem,
114 KINSOL<VectorType> &solver = *
static_cast<KINSOL<VectorType> *
>(kinsol_mem->kin_user_data);
118 solver.reinit_vector(*src_ycur);
121 solver.reinit_vector(*src_fcur);
123 copy(*src_ycur, kinsol_mem->kin_uu);
124 copy(*src_fcur, kinsol_mem->kin_fval);
127 solver.reinit_vector(*src);
130 solver.reinit_vector(*dst);
134 int err = solver.solve_jacobian_system(*src_ycur, *src_fcur,
138 *sJpnorm = N_VWL2Norm(b, kinsol_mem->kin_fscale);
139 N_VProd(b, kinsol_mem->kin_fscale, b);
140 N_VProd(b, kinsol_mem->kin_fscale, b);
141 *sFdotJp = N_VDotProd(kinsol_mem->kin_fval, b);
147 template <
typename VectorType>
156 Utilities::MPI::duplicate_communicator(mpi_comm))
163 template <
typename VectorType>
167 KINFree(&kinsol_mem);
168 #ifdef DEAL_II_WITH_MPI 171 const int ierr = MPI_Comm_free(&communicator);
180 template <
typename VectorType>
183 unsigned int system_size = initial_guess_and_solution.size();
188 #ifdef DEAL_II_WITH_MPI 191 const IndexSet is = initial_guess_and_solution.locally_owned_elements();
192 const unsigned int local_system_size = is.
n_elements();
194 solution = N_VNew_Parallel(communicator,
198 u_scale = N_VNew_Parallel(communicator,
201 N_VConst_Parallel( 1.e0, u_scale );
203 f_scale = N_VNew_Parallel(communicator,
206 N_VConst_Parallel( 1.e0, f_scale );
213 solution = N_VNew_Serial(system_size);
214 u_scale = N_VNew_Serial(system_size);
215 N_VConst_Serial( 1.e0, u_scale );
216 f_scale = N_VNew_Serial(system_size);
217 N_VConst_Serial( 1.e0, f_scale );
220 if (get_solution_scaling)
221 copy(u_scale, get_solution_scaling());
223 if (get_function_scaling)
224 copy(f_scale, get_function_scaling());
226 copy(solution, initial_guess_and_solution);
229 KINFree(&kinsol_mem);
231 kinsol_mem = KINCreate();
233 int status = KINInit(kinsol_mem, t_kinsol_function<VectorType> , solution);
235 AssertKINSOL(status);
237 status = KINSetUserData(kinsol_mem, (
void *)
this);
238 AssertKINSOL(status);
240 status = KINSetNumMaxIters(kinsol_mem, data.maximum_non_linear_iterations);
241 AssertKINSOL(status);
243 status = KINSetFuncNormTol(kinsol_mem, data.function_tolerance);
244 AssertKINSOL(status);
246 status = KINSetScaledStepTol(kinsol_mem, data.step_tolerance);
247 AssertKINSOL(status);
249 status = KINSetMaxSetupCalls(kinsol_mem, data.maximum_setup_calls);
250 AssertKINSOL(status);
252 status = KINSetNoInitSetup(kinsol_mem, (
int) data.no_init_setup);
253 AssertKINSOL(status);
255 status = KINSetMaxNewtonStep(kinsol_mem, data.maximum_newton_step);
256 AssertKINSOL(status);
258 status = KINSetMaxBetaFails(kinsol_mem, data.maximum_beta_failures);
259 AssertKINSOL(status);
261 status = KINSetMAA(kinsol_mem, data.anderson_subspace_size);
262 AssertKINSOL(status);
264 status = KINSetRelErrFunc(kinsol_mem, data.dq_relative_error);
265 AssertKINSOL(status);
267 #if DEAL_II_SUNDIALS_VERSION_GTE(3,0,0) 268 SUNMatrix J =
nullptr;
269 SUNLinearSolver LS =
nullptr;
272 if (solve_jacobian_system)
274 KINMem KIN_mem = (KINMem) kinsol_mem;
275 KIN_mem->kin_lsolve = t_kinsol_solve_jacobian<VectorType>;
278 KIN_mem->kin_lsetup = t_kinsol_setup_jacobian<VectorType>;
279 #if DEAL_II_SUNDIALS_VERSION_LT(3,0,0) 280 KIN_mem->kin_setupNonNull =
true;
286 #if DEAL_II_SUNDIALS_VERSION_GTE(3,0,0) 287 J = SUNDenseMatrix(system_size, system_size);
288 LS = SUNDenseLinearSolver(u_scale, J);
289 status = KINDlsSetLinearSolver(kinsol_mem, LS, J);
291 status = KINDense(kinsol_mem, system_size);
293 AssertKINSOL(status);
296 if (data.strategy == AdditionalData::newton ||
297 data.strategy == AdditionalData::linesearch)
298 Assert(residual, ExcFunctionNotProvided(
"residual"));
300 if (data.strategy == AdditionalData::fixed_point ||
301 data.strategy == AdditionalData::picard)
302 Assert(iteration_function, ExcFunctionNotProvided(
"iteration_function"));
305 status = KINSol(kinsol_mem, solution, (
int) data.strategy, u_scale, f_scale);
306 AssertKINSOL(status);
308 copy(initial_guess_and_solution, solution );
311 #ifdef DEAL_II_WITH_MPI 314 N_VDestroy_Parallel(solution);
315 N_VDestroy_Parallel(u_scale);
316 N_VDestroy_Parallel(f_scale);
321 N_VDestroy_Serial(solution);
322 N_VDestroy_Serial(u_scale);
323 N_VDestroy_Serial(f_scale);
327 status = KINGetNumNonlinSolvIters(kinsol_mem, &nniters);
328 AssertKINSOL(status);
330 #if DEAL_II_SUNDIALS_VERSION_GTE(3,0,0) 334 KINFree(&kinsol_mem);
336 return (
unsigned int) nniters;
339 template<
typename VectorType>
343 reinit_vector = [](VectorType &)
345 AssertThrow(
false, ExcFunctionNotProvided(
"reinit_vector"));
353 #ifdef DEAL_II_WITH_MPI 355 #ifdef DEAL_II_WITH_TRILINOS 360 #ifdef DEAL_II_WITH_PETSC 361 #ifndef PETSC_USE_COMPLEX 371 DEAL_II_NAMESPACE_CLOSE
#define AssertNothrow(cond, exc)
void set_functions_to_trigger_an_assert()
#define AssertThrow(cond, exc)
#define Assert(cond, exc)
KINSOL(const AdditionalData &data=AdditionalData(), const MPI_Comm mpi_comm=MPI_COMM_WORLD)
unsigned int solve(VectorType &initial_guess_and_solution)
size_type n_elements() const
static ::ExceptionBase & ExcInternalError()