Reference documentation for deal.II version 9.0.0
kinsol.cc
1 //-----------------------------------------------------------
2 //
3 // Copyright (C) 2017 - 2018 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 at
12 // the top level of the deal.II distribution.
13 //
14 //-----------------------------------------------------------
15 
16 
17 #include <deal.II/sundials/kinsol.h>
18 #include <deal.II/base/config.h>
19 
20 #ifdef DEAL_II_WITH_SUNDIALS
21 
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>
27 #endif
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>
31 #endif
32 #include <deal.II/base/utilities.h>
33 #include <deal.II/sundials/copy.h>
34 
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>
40 #else
41 # include <kinsol/kinsol_dense.h>
42 #endif
43 
44 #include <iostream>
45 #include <iomanip>
46 
47 DEAL_II_NAMESPACE_OPEN
48 
49 namespace SUNDIALS
50 {
51  using namespace internal;
52 
53  namespace
54  {
55  template<typename VectorType>
56  int t_kinsol_function(N_Vector yy,
57  N_Vector FF,
58  void *user_data)
59  {
60  KINSOL<VectorType> &solver = *static_cast<KINSOL<VectorType> *>(user_data);
62 
63  typename VectorMemory<VectorType>::Pointer src_yy(mem);
64  solver.reinit_vector(*src_yy);
65 
66  typename VectorMemory<VectorType>::Pointer dst_FF(mem);
67  solver.reinit_vector(*dst_FF);
68 
69  copy(*src_yy, yy);
70 
71  int err = 0;
72  if (solver.residual)
73  err = solver.residual(*src_yy, *dst_FF);
74  else if (solver.iteration_function)
75  err = solver.iteration_function(*src_yy, *dst_FF);
76  else
77  Assert(false, ExcInternalError());
78 
79  copy(FF, *dst_FF);
80 
81  return err;
82  }
83 
84 
85 
86  template<typename VectorType>
87  int t_kinsol_setup_jacobian(KINMem kinsol_mem)
88  {
89  KINSOL<VectorType> &solver = *static_cast<KINSOL<VectorType> *>(kinsol_mem->kin_user_data);
91 
92  typename VectorMemory<VectorType>::Pointer src_ycur(mem);
93  solver.reinit_vector(*src_ycur);
94 
95  typename VectorMemory<VectorType>::Pointer src_fcur(mem);
96  solver.reinit_vector(*src_fcur);
97 
98  copy(*src_ycur, kinsol_mem->kin_uu);
99  copy(*src_fcur, kinsol_mem->kin_fval);
100 
101  int err = solver.setup_jacobian(*src_ycur, *src_fcur);
102  return err;
103  }
104 
105 
106 
107  template<typename VectorType>
108  int t_kinsol_solve_jacobian(KINMem kinsol_mem,
109  N_Vector x,
110  N_Vector b,
111  realtype *sJpnorm,
112  realtype *sFdotJp)
113  {
114  KINSOL<VectorType> &solver = *static_cast<KINSOL<VectorType> *>(kinsol_mem->kin_user_data);
116 
117  typename VectorMemory<VectorType>::Pointer src_ycur(mem);
118  solver.reinit_vector(*src_ycur);
119 
120  typename VectorMemory<VectorType>::Pointer src_fcur(mem);
121  solver.reinit_vector(*src_fcur);
122 
123  copy(*src_ycur, kinsol_mem->kin_uu);
124  copy(*src_fcur, kinsol_mem->kin_fval);
125 
126  typename VectorMemory<VectorType>::Pointer src(mem);
127  solver.reinit_vector(*src);
128 
129  typename VectorMemory<VectorType>::Pointer dst(mem);
130  solver.reinit_vector(*dst);
131 
132  copy(*src, b);
133 
134  int err = solver.solve_jacobian_system(*src_ycur, *src_fcur,
135  *src,*dst);
136  copy(x, *dst);
137 
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);
142 
143  return err;
144  }
145  }
146 
147  template <typename VectorType>
148  KINSOL<VectorType>::KINSOL(const AdditionalData &data, const MPI_Comm mpi_comm) :
149  data(data),
150  kinsol_mem(nullptr),
151  solution(nullptr),
152  u_scale(nullptr),
153  f_scale(nullptr),
154  communicator(is_serial_vector<VectorType>::value ?
155  MPI_COMM_SELF :
156  Utilities::MPI::duplicate_communicator(mpi_comm))
157  {
159  }
160 
161 
162 
163  template <typename VectorType>
165  {
166  if (kinsol_mem)
167  KINFree(&kinsol_mem);
168 #ifdef DEAL_II_WITH_MPI
170  {
171  const int ierr = MPI_Comm_free(&communicator);
172  (void)ierr;
173  AssertNothrow(ierr == MPI_SUCCESS, ExcMPI(ierr));
174  }
175 #endif
176  }
177 
178 
179 
180  template <typename VectorType>
181  unsigned int KINSOL<VectorType>::solve(VectorType &initial_guess_and_solution)
182  {
183  unsigned int system_size = initial_guess_and_solution.size();
184 
185  // The solution is stored in
186  // solution. Here we take only a
187  // view of it.
188 #ifdef DEAL_II_WITH_MPI
190  {
191  const IndexSet is = initial_guess_and_solution.locally_owned_elements();
192  const unsigned int local_system_size = is.n_elements();
193 
194  solution = N_VNew_Parallel(communicator,
195  local_system_size,
196  system_size);
197 
198  u_scale = N_VNew_Parallel(communicator,
199  local_system_size,
200  system_size);
201  N_VConst_Parallel( 1.e0, u_scale );
202 
203  f_scale = N_VNew_Parallel(communicator,
204  local_system_size,
205  system_size);
206  N_VConst_Parallel( 1.e0, f_scale );
207  }
208  else
209 #endif
210  {
212  ExcInternalError("Trying to use a serial code with a parallel vector."));
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 );
218  }
219 
220  if (get_solution_scaling)
221  copy(u_scale, get_solution_scaling());
222 
223  if (get_function_scaling)
224  copy(f_scale, get_function_scaling());
225 
226  copy(solution, initial_guess_and_solution);
227 
228  if (kinsol_mem)
229  KINFree(&kinsol_mem);
230 
231  kinsol_mem = KINCreate();
232 
233  int status = KINInit(kinsol_mem, t_kinsol_function<VectorType> , solution);
234  (void) status;
235  AssertKINSOL(status);
236 
237  status = KINSetUserData(kinsol_mem, (void *) this);
238  AssertKINSOL(status);
239 
240  status = KINSetNumMaxIters(kinsol_mem, data.maximum_non_linear_iterations);
241  AssertKINSOL(status);
242 
243  status = KINSetFuncNormTol(kinsol_mem, data.function_tolerance);
244  AssertKINSOL(status);
245 
246  status = KINSetScaledStepTol(kinsol_mem, data.step_tolerance);
247  AssertKINSOL(status);
248 
249  status = KINSetMaxSetupCalls(kinsol_mem, data.maximum_setup_calls);
250  AssertKINSOL(status);
251 
252  status = KINSetNoInitSetup(kinsol_mem, (int) data.no_init_setup);
253  AssertKINSOL(status);
254 
255  status = KINSetMaxNewtonStep(kinsol_mem, data.maximum_newton_step);
256  AssertKINSOL(status);
257 
258  status = KINSetMaxBetaFails(kinsol_mem, data.maximum_beta_failures);
259  AssertKINSOL(status);
260 
261  status = KINSetMAA(kinsol_mem, data.anderson_subspace_size);
262  AssertKINSOL(status);
263 
264  status = KINSetRelErrFunc(kinsol_mem, data.dq_relative_error);
265  AssertKINSOL(status);
266 
267 #if DEAL_II_SUNDIALS_VERSION_GTE(3,0,0)
268  SUNMatrix J = nullptr;
269  SUNLinearSolver LS = nullptr;
270 #endif
271 
272  if (solve_jacobian_system)
273  {
274  KINMem KIN_mem = (KINMem) kinsol_mem;
275  KIN_mem->kin_lsolve = t_kinsol_solve_jacobian<VectorType>;
276  if (setup_jacobian)
277  {
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;
281 #endif
282  }
283  }
284  else
285  {
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);
290 #else
291  status = KINDense(kinsol_mem, system_size);
292 #endif
293  AssertKINSOL(status);
294  }
295 
296  if (data.strategy == AdditionalData::newton ||
297  data.strategy == AdditionalData::linesearch)
298  Assert(residual, ExcFunctionNotProvided("residual"));
299 
300  if (data.strategy == AdditionalData::fixed_point ||
301  data.strategy == AdditionalData::picard)
302  Assert(iteration_function, ExcFunctionNotProvided("iteration_function"));
303 
304  // call to KINSol
305  status = KINSol(kinsol_mem, solution, (int) data.strategy, u_scale, f_scale);
306  AssertKINSOL(status);
307 
308  copy(initial_guess_and_solution, solution );
309 
310  // Free the vectors which are no longer used.
311 #ifdef DEAL_II_WITH_MPI
313  {
314  N_VDestroy_Parallel(solution);
315  N_VDestroy_Parallel(u_scale);
316  N_VDestroy_Parallel(f_scale);
317  }
318  else
319 #endif
320  {
321  N_VDestroy_Serial(solution);
322  N_VDestroy_Serial(u_scale);
323  N_VDestroy_Serial(f_scale);
324  }
325 
326  long nniters;
327  status = KINGetNumNonlinSolvIters(kinsol_mem, &nniters);
328  AssertKINSOL(status);
329 
330 #if DEAL_II_SUNDIALS_VERSION_GTE(3,0,0)
331  SUNMatDestroy(J);
332  SUNLinSolFree(LS);
333 #endif
334  KINFree(&kinsol_mem);
335 
336  return (unsigned int) nniters;
337  }
338 
339  template<typename VectorType>
341  {
342 
343  reinit_vector = [](VectorType &)
344  {
345  AssertThrow(false, ExcFunctionNotProvided("reinit_vector"));
346  };
347 
348  }
349 
350  template class KINSOL<Vector<double> >;
351  template class KINSOL<BlockVector<double> >;
352 
353 #ifdef DEAL_II_WITH_MPI
354 
355 #ifdef DEAL_II_WITH_TRILINOS
358 #endif
359 
360 #ifdef DEAL_II_WITH_PETSC
361 #ifndef PETSC_USE_COMPLEX
362  template class KINSOL<PETScWrappers::MPI::Vector>;
364 #endif
365 #endif
366 
367 #endif
368 
369 }
370 
371 DEAL_II_NAMESPACE_CLOSE
372 
373 #endif
#define AssertNothrow(cond, exc)
Definition: exceptions.h:1183
void set_functions_to_trigger_an_assert()
Definition: kinsol.cc:340
#define AssertThrow(cond, exc)
Definition: exceptions.h:1221
#define Assert(cond, exc)
Definition: exceptions.h:1142
KINSOL(const AdditionalData &data=AdditionalData(), const MPI_Comm mpi_comm=MPI_COMM_WORLD)
Definition: kinsol.cc:148
Definition: cuda.h:31
unsigned int solve(VectorType &initial_guess_and_solution)
Definition: kinsol.cc:181
size_type n_elements() const
Definition: index_set.h:1717
static ::ExceptionBase & ExcInternalError()