petsc_linear_solver.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2017 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 #include "libmesh/libmesh_common.h"
19 
20 #ifdef LIBMESH_HAVE_PETSC
21 
22 
23 // C++ includes
24 #include <string.h>
25 
26 // Local Includes
27 #include "libmesh/dof_map.h"
30 #include "libmesh/shell_matrix.h"
31 #include "libmesh/petsc_matrix.h"
33 #include "libmesh/petsc_vector.h"
34 #include "libmesh/string_to_enum.h"
35 #include "libmesh/system.h"
38 
39 namespace libMesh
40 {
41 
42 extern "C"
43 {
44 #if PETSC_RELEASE_LESS_THAN(3,0,1)
45  PetscErrorCode __libmesh_petsc_preconditioner_setup (void * ctx)
46  {
47  Preconditioner<Number> * preconditioner = static_cast<Preconditioner<Number> *>(ctx);
48 
49  if (!preconditioner->initialized())
50  libmesh_error_msg("Preconditioner not initialized! Make sure you call init() before solve!");
51 
52  preconditioner->setup();
53 
54  return 0;
55  }
56 
57 
58  PetscErrorCode __libmesh_petsc_preconditioner_apply(void * ctx, Vec x, Vec y)
59  {
60  Preconditioner<Number> * preconditioner = static_cast<Preconditioner<Number> *>(ctx);
61 
62  PetscVector<Number> x_vec(x, preconditioner->comm());
63  PetscVector<Number> y_vec(y, preconditioner->comm());
64 
65  preconditioner->apply(x_vec,y_vec);
66 
67  return 0;
68  }
69 #else
71  {
72  void * ctx;
73  PetscErrorCode ierr = PCShellGetContext(pc,&ctx);CHKERRQ(ierr);
74  Preconditioner<Number> * preconditioner = static_cast<Preconditioner<Number> *>(ctx);
75 
76  if (!preconditioner->initialized())
77  libmesh_error_msg("Preconditioner not initialized! Make sure you call init() before solve!");
78 
79  preconditioner->setup();
80 
81  return 0;
82  }
83 
84  PetscErrorCode __libmesh_petsc_preconditioner_apply(PC pc, Vec x, Vec y)
85  {
86  void * ctx;
87  PetscErrorCode ierr = PCShellGetContext(pc,&ctx);CHKERRQ(ierr);
88  Preconditioner<Number> * preconditioner = static_cast<Preconditioner<Number> *>(ctx);
89 
90  PetscVector<Number> x_vec(x, preconditioner->comm());
91  PetscVector<Number> y_vec(y, preconditioner->comm());
92 
93  preconditioner->apply(x_vec,y_vec);
94 
95  return 0;
96  }
97 #endif
98 } // end extern "C"
99 
100 /*----------------------- functions ----------------------------------*/
101 template <typename T>
103 {
104  if (this->initialized())
105  {
106  /* If we were restricted to some subset, this restriction must
107  be removed and the subset index set destroyed. */
108  if (_restrict_solve_to_is != libmesh_nullptr)
109  {
110  PetscErrorCode ierr = LibMeshISDestroy(&_restrict_solve_to_is);
111  LIBMESH_CHKERR(ierr);
112  _restrict_solve_to_is = libmesh_nullptr;
113  }
114 
115  if(_restrict_solve_to_is_complement != libmesh_nullptr)
116  {
117  PetscErrorCode ierr = LibMeshISDestroy(&_restrict_solve_to_is_complement);
118  LIBMESH_CHKERR(ierr);
119  _restrict_solve_to_is_complement = libmesh_nullptr;
120  }
121 
122  this->_is_initialized = false;
123 
124  PetscErrorCode ierr=0;
125 
126  ierr = LibMeshKSPDestroy(&_ksp);
127  LIBMESH_CHKERR(ierr);
128 
129  // Mimic PETSc default solver and preconditioner
130  this->_solver_type = GMRES;
131 
132  if(!this->_preconditioner)
133  {
134  if (this->n_processors() == 1)
135  this->_preconditioner_type = ILU_PRECOND;
136  else
137  this->_preconditioner_type = BLOCK_JACOBI_PRECOND;
138  }
139  }
140 }
141 
142 
143 
144 template <typename T>
145 void PetscLinearSolver<T>::init (const char * name)
146 {
147  // Initialize the data structures if not done so already.
148  if (!this->initialized())
149  {
150  this->_is_initialized = true;
151 
152  PetscErrorCode ierr=0;
153 
154  // Create the linear solver context
155  ierr = KSPCreate (this->comm().get(), &_ksp);
156  LIBMESH_CHKERR(ierr);
157 
158  if (name)
159  {
160  ierr = KSPSetOptionsPrefix(_ksp, name);
161  LIBMESH_CHKERR(ierr);
162  }
163 
164  // Create the preconditioner context
165  ierr = KSPGetPC (_ksp, &_pc);
166  LIBMESH_CHKERR(ierr);
167 
168  // Set user-specified solver and preconditioner types
169  this->set_petsc_solver_type();
170 
171  // If the SolverConfiguration object is provided, use it to set
172  // options during solver initialization.
173  if(this->_solver_configuration)
174  {
175  this->_solver_configuration->set_options_during_init();
176  }
177 
178  // Set the options from user-input
179  // Set runtime options, e.g.,
180  // -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
181  // These options will override those specified above as long as
182  // KSPSetFromOptions() is called _after_ any other customization
183  // routines.
184 
185  ierr = KSPSetFromOptions (_ksp);
186  LIBMESH_CHKERR(ierr);
187 
188  // Not sure if this is necessary, or if it is already handled by KSPSetFromOptions?
189  // NOT NECESSARY!!!!
190  //ierr = PCSetFromOptions (_pc);
191  //LIBMESH_CHKERR(ierr);
192 
193  // Have the Krylov subspace method use our good initial guess
194  // rather than 0, unless the user requested a KSPType of
195  // preonly, which complains if asked to use initial guesses.
196 #if PETSC_VERSION_LESS_THAN(3,0,0) || !PETSC_RELEASE_LESS_THAN(3,4,0)
197  // Pre-3.0 and petsc-dev (as of October 2012) use non-const versions
198  KSPType ksp_type;
199 #else
200  const KSPType ksp_type;
201 #endif
202 
203  ierr = KSPGetType (_ksp, &ksp_type);
204  LIBMESH_CHKERR(ierr);
205 
206  if (strcmp(ksp_type, "preonly"))
207  {
208  ierr = KSPSetInitialGuessNonzero (_ksp, PETSC_TRUE);
209  LIBMESH_CHKERR(ierr);
210  }
211 
212  // Notify PETSc of location to store residual history.
213  // This needs to be called before any solves, since
214  // it sets the residual history length to zero. The default
215  // behavior is for PETSc to allocate (internally) an array
216  // of size 1000 to hold the residual norm history.
217  ierr = KSPSetResidualHistory(_ksp,
218  PETSC_NULL, // pointer to the array which holds the history
219  PETSC_DECIDE, // size of the array holding the history
220  PETSC_TRUE); // Whether or not to reset the history for each solve.
221  LIBMESH_CHKERR(ierr);
222 
223  PetscPreconditioner<T>::set_petsc_preconditioner_type(this->_preconditioner_type,_pc);
224 
225  //If there is a preconditioner object we need to set the internal setup and apply routines
226  if(this->_preconditioner)
227  {
228  this->_preconditioner->init();
229  PCShellSetContext(_pc,(void *)this->_preconditioner);
230  PCShellSetSetUp(_pc,__libmesh_petsc_preconditioner_setup);
231  PCShellSetApply(_pc,__libmesh_petsc_preconditioner_apply);
232  }
233  }
234 }
235 
236 
237 template <typename T>
239  const char * name)
240 {
241  // Initialize the data structures if not done so already.
242  if (!this->initialized())
243  {
244  this->_is_initialized = true;
245 
246  PetscErrorCode ierr=0;
247 
248  // Create the linear solver context
249  ierr = KSPCreate (this->comm().get(), &_ksp);
250  LIBMESH_CHKERR(ierr);
251 
252  if (name)
253  {
254  ierr = KSPSetOptionsPrefix(_ksp, name);
255  LIBMESH_CHKERR(ierr);
256  }
257 
258  //ierr = PCCreate (this->comm().get(), &_pc);
259  // LIBMESH_CHKERR(ierr);
260 
261  // Create the preconditioner context
262  ierr = KSPGetPC (_ksp, &_pc);
263  LIBMESH_CHKERR(ierr);
264 
265  // Set operators. The input matrix works as the preconditioning matrix
266 #if PETSC_RELEASE_LESS_THAN(3,5,0)
267  ierr = KSPSetOperators(_ksp, matrix->mat(), matrix->mat(),DIFFERENT_NONZERO_PATTERN);
268 #else
269  ierr = KSPSetOperators(_ksp, matrix->mat(), matrix->mat());
270 #endif
271  LIBMESH_CHKERR(ierr);
272 
273  // Set user-specified solver and preconditioner types
274  this->set_petsc_solver_type();
275 
276  // If the SolverConfiguration object is provided, use it to set
277  // options during solver initialization.
278  if(this->_solver_configuration)
279  {
280  this->_solver_configuration->set_options_during_init();
281  }
282 
283  // Set the options from user-input
284  // Set runtime options, e.g.,
285  // -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
286  // These options will override those specified above as long as
287  // KSPSetFromOptions() is called _after_ any other customization
288  // routines.
289 
290  ierr = KSPSetFromOptions (_ksp);
291  LIBMESH_CHKERR(ierr);
292 
293  // Not sure if this is necessary, or if it is already handled by KSPSetFromOptions?
294  // NOT NECESSARY!!!!
295  //ierr = PCSetFromOptions (_pc);
296  //LIBMESH_CHKERR(ierr);
297 
298  // Have the Krylov subspace method use our good initial guess
299  // rather than 0, unless the user requested a KSPType of
300  // preonly, which complains if asked to use initial guesses.
301 #if PETSC_VERSION_LESS_THAN(3,0,0) || !PETSC_RELEASE_LESS_THAN(3,4,0)
302  KSPType ksp_type;
303 #else
304  const KSPType ksp_type;
305 #endif
306 
307  ierr = KSPGetType (_ksp, &ksp_type);
308  LIBMESH_CHKERR(ierr);
309 
310  if (strcmp(ksp_type, "preonly"))
311  {
312  ierr = KSPSetInitialGuessNonzero (_ksp, PETSC_TRUE);
313  LIBMESH_CHKERR(ierr);
314  }
315 
316  // Notify PETSc of location to store residual history.
317  // This needs to be called before any solves, since
318  // it sets the residual history length to zero. The default
319  // behavior is for PETSc to allocate (internally) an array
320  // of size 1000 to hold the residual norm history.
321  ierr = KSPSetResidualHistory(_ksp,
322  PETSC_NULL, // pointer to the array which holds the history
323  PETSC_DECIDE, // size of the array holding the history
324  PETSC_TRUE); // Whether or not to reset the history for each solve.
325  LIBMESH_CHKERR(ierr);
326 
327  PetscPreconditioner<T>::set_petsc_preconditioner_type(this->_preconditioner_type,_pc);
328  if(this->_preconditioner)
329  {
330  this->_preconditioner->set_matrix(*matrix);
331  this->_preconditioner->init();
332  PCShellSetContext(_pc,(void *)this->_preconditioner);
333  PCShellSetSetUp(_pc,__libmesh_petsc_preconditioner_setup);
334  PCShellSetApply(_pc,__libmesh_petsc_preconditioner_apply);
335  }
336  }
337 }
338 
339 
340 
341 template <typename T>
342 void
344 {
345  petsc_auto_fieldsplit(this->pc(), sys);
346 }
347 
348 
349 
350 template <typename T>
351 void
352 PetscLinearSolver<T>::restrict_solve_to (const std::vector<unsigned int> * const dofs,
353  const SubsetSolveMode subset_solve_mode)
354 {
355  PetscErrorCode ierr=0;
356 
357  /* The preconditioner (in particular if a default preconditioner)
358  will have to be reset. We call this->clear() to do that. This
359  call will also remove and free any previous subset that may have
360  been set before. */
361  this->clear();
362 
363  _subset_solve_mode = subset_solve_mode;
364 
365  if (dofs != libmesh_nullptr)
366  {
367  PetscInt * petsc_dofs = libmesh_nullptr;
368  ierr = PetscMalloc(dofs->size()*sizeof(PetscInt), &petsc_dofs);
369  LIBMESH_CHKERR(ierr);
370 
371  for (std::size_t i=0; i<dofs->size(); i++)
372  petsc_dofs[i] = (*dofs)[i];
373 
374  ierr = ISCreateLibMesh(this->comm().get(),dofs->size(),petsc_dofs,PETSC_OWN_POINTER,&_restrict_solve_to_is);
375  LIBMESH_CHKERR(ierr);
376  }
377 }
378 
379 
380 
381 template <typename T>
382 std::pair<unsigned int, Real>
384  SparseMatrix<T> & precond_in,
385  NumericVector<T> & solution_in,
386  NumericVector<T> & rhs_in,
387  const double tol,
388  const unsigned int m_its)
389 {
390  LOG_SCOPE("solve()", "PetscLinearSolver");
391 
392  // Make sure the data passed in are really of Petsc types
393  PetscMatrix<T> * matrix = cast_ptr<PetscMatrix<T> *>(&matrix_in);
394  PetscMatrix<T> * precond = cast_ptr<PetscMatrix<T> *>(&precond_in);
395  PetscVector<T> * solution = cast_ptr<PetscVector<T> *>(&solution_in);
396  PetscVector<T> * rhs = cast_ptr<PetscVector<T> *>(&rhs_in);
397 
398  this->init (matrix);
399 
400  PetscErrorCode ierr=0;
401  PetscInt its=0, max_its = static_cast<PetscInt>(m_its);
402  PetscReal final_resid=0.;
403 
404  // Close the matrices and vectors in case this wasn't already done.
405  matrix->close ();
406  precond->close ();
407  solution->close ();
408  rhs->close ();
409 
410  // // If matrix != precond, then this means we have specified a
411  // // special preconditioner, so reset preconditioner type to PCMAT.
412  // if (matrix != precond)
413  // {
414  // this->_preconditioner_type = USER_PRECOND;
415  // this->set_petsc_preconditioner_type ();
416  // }
417 
418  Mat submat = libmesh_nullptr;
419  Mat subprecond = libmesh_nullptr;
420  Vec subrhs = libmesh_nullptr;
421  Vec subsolution = libmesh_nullptr;
422  VecScatter scatter = libmesh_nullptr;
423  UniquePtr<PetscMatrix<Number> > subprecond_matrix;
424 
425  // Set operators. Also restrict rhs and solution vector to
426  // subdomain if neccessary.
427  if (_restrict_solve_to_is != libmesh_nullptr)
428  {
429  PetscInt is_local_size = this->_restrict_solve_to_is_local_size();
430 
431  ierr = VecCreate(this->comm().get(),&subrhs);
432  LIBMESH_CHKERR(ierr);
433  ierr = VecSetSizes(subrhs,is_local_size,PETSC_DECIDE);
434  LIBMESH_CHKERR(ierr);
435  ierr = VecSetFromOptions(subrhs);
436  LIBMESH_CHKERR(ierr);
437 
438  ierr = VecCreate(this->comm().get(),&subsolution);
439  LIBMESH_CHKERR(ierr);
440  ierr = VecSetSizes(subsolution,is_local_size,PETSC_DECIDE);
441  LIBMESH_CHKERR(ierr);
442  ierr = VecSetFromOptions(subsolution);
443  LIBMESH_CHKERR(ierr);
444 
445  ierr = VecScatterCreate(rhs->vec(), _restrict_solve_to_is, subrhs, libmesh_nullptr, &scatter);
446  LIBMESH_CHKERR(ierr);
447 
448  ierr = VecScatterBegin(scatter,rhs->vec(),subrhs,INSERT_VALUES,SCATTER_FORWARD);
449  LIBMESH_CHKERR(ierr);
450  ierr = VecScatterEnd(scatter,rhs->vec(),subrhs,INSERT_VALUES,SCATTER_FORWARD);
451  LIBMESH_CHKERR(ierr);
452 
453  ierr = VecScatterBegin(scatter,solution->vec(),subsolution,INSERT_VALUES,SCATTER_FORWARD);
454  LIBMESH_CHKERR(ierr);
455  ierr = VecScatterEnd(scatter,solution->vec(),subsolution,INSERT_VALUES,SCATTER_FORWARD);
456  LIBMESH_CHKERR(ierr);
457 
458 #if PETSC_VERSION_LESS_THAN(3,1,0)
459  ierr = MatGetSubMatrix(matrix->mat(),
460  _restrict_solve_to_is,_restrict_solve_to_is,
461  PETSC_DECIDE,MAT_INITIAL_MATRIX,&submat);
462  LIBMESH_CHKERR(ierr);
463  ierr = MatGetSubMatrix(precond->mat(),
464  _restrict_solve_to_is,_restrict_solve_to_is,
465  PETSC_DECIDE,MAT_INITIAL_MATRIX,&subprecond);
466  LIBMESH_CHKERR(ierr);
467 #else
468  ierr = MatGetSubMatrix(matrix->mat(),
469  _restrict_solve_to_is,_restrict_solve_to_is,
470  MAT_INITIAL_MATRIX,&submat);
471  LIBMESH_CHKERR(ierr);
472  ierr = MatGetSubMatrix(precond->mat(),
473  _restrict_solve_to_is,_restrict_solve_to_is,
474  MAT_INITIAL_MATRIX,&subprecond);
475  LIBMESH_CHKERR(ierr);
476 #endif
477 
478  /* Since removing columns of the matrix changes the equation
479  system, we will now change the right hand side to compensate
480  for this. Note that this is not necessary if \p SUBSET_ZERO
481  has been selected. */
482  if(_subset_solve_mode!=SUBSET_ZERO)
483  {
484  _create_complement_is(rhs_in);
485  PetscInt is_complement_local_size =
486  cast_int<PetscInt>(rhs_in.local_size()-is_local_size);
487 
488  Vec subvec1 = libmesh_nullptr;
489  Mat submat1 = libmesh_nullptr;
490  VecScatter scatter1 = libmesh_nullptr;
491 
492  ierr = VecCreate(this->comm().get(),&subvec1);
493  LIBMESH_CHKERR(ierr);
494  ierr = VecSetSizes(subvec1,is_complement_local_size,PETSC_DECIDE);
495  LIBMESH_CHKERR(ierr);
496  ierr = VecSetFromOptions(subvec1);
497  LIBMESH_CHKERR(ierr);
498 
499  ierr = VecScatterCreate(rhs->vec(),_restrict_solve_to_is_complement, subvec1, libmesh_nullptr, &scatter1);
500  LIBMESH_CHKERR(ierr);
501 
502  ierr = VecScatterBegin(scatter1,_subset_solve_mode==SUBSET_COPY_RHS ? rhs->vec() : solution->vec(),subvec1,INSERT_VALUES,SCATTER_FORWARD);
503  LIBMESH_CHKERR(ierr);
504  ierr = VecScatterEnd(scatter1,_subset_solve_mode==SUBSET_COPY_RHS ? rhs->vec() : solution->vec(),subvec1,INSERT_VALUES,SCATTER_FORWARD);
505  LIBMESH_CHKERR(ierr);
506 
507  ierr = VecScale(subvec1,-1.0);
508  LIBMESH_CHKERR(ierr);
509 
510 #if PETSC_VERSION_LESS_THAN(3,1,0)
511  ierr = MatGetSubMatrix(matrix->mat(),
512  _restrict_solve_to_is,_restrict_solve_to_is_complement,
513  PETSC_DECIDE,MAT_INITIAL_MATRIX,&submat1);
514  LIBMESH_CHKERR(ierr);
515 #else
516  ierr = MatGetSubMatrix(matrix->mat(),
517  _restrict_solve_to_is,_restrict_solve_to_is_complement,
518  MAT_INITIAL_MATRIX,&submat1);
519  LIBMESH_CHKERR(ierr);
520 #endif
521 
522  ierr = MatMultAdd(submat1,subvec1,subrhs,subrhs);
523  LIBMESH_CHKERR(ierr);
524 
525  ierr = LibMeshVecScatterDestroy(&scatter1);
526  LIBMESH_CHKERR(ierr);
527  ierr = LibMeshVecDestroy(&subvec1);
528  LIBMESH_CHKERR(ierr);
529  ierr = LibMeshMatDestroy(&submat1);
530  LIBMESH_CHKERR(ierr);
531  }
532 #if PETSC_RELEASE_LESS_THAN(3,5,0)
533  ierr = KSPSetOperators(_ksp, submat, subprecond,
534  this->same_preconditioner ? SAME_PRECONDITIONER : DIFFERENT_NONZERO_PATTERN);
535 #else
536  ierr = KSPSetOperators(_ksp, submat, subprecond);
537 
538  PetscBool ksp_reuse_preconditioner = this->same_preconditioner ? PETSC_TRUE : PETSC_FALSE;
539  ierr = KSPSetReusePreconditioner(_ksp, ksp_reuse_preconditioner);
540 #endif
541  LIBMESH_CHKERR(ierr);
542 
543  if (this->_preconditioner)
544  {
545  subprecond_matrix.reset(new PetscMatrix<Number>(subprecond, this->comm()));
546  this->_preconditioner->set_matrix(*subprecond_matrix);
547  this->_preconditioner->init();
548  }
549  }
550  else
551  {
552 #if PETSC_RELEASE_LESS_THAN(3,5,0)
553  ierr = KSPSetOperators(_ksp, matrix->mat(), precond->mat(),
554  this->same_preconditioner ? SAME_PRECONDITIONER : DIFFERENT_NONZERO_PATTERN);
555 #else
556  ierr = KSPSetOperators(_ksp, matrix->mat(), precond->mat());
557 
558  PetscBool ksp_reuse_preconditioner = this->same_preconditioner ? PETSC_TRUE : PETSC_FALSE;
559  ierr = KSPSetReusePreconditioner(_ksp, ksp_reuse_preconditioner);
560 #endif
561  LIBMESH_CHKERR(ierr);
562 
563  if(this->_preconditioner)
564  {
565  this->_preconditioner->set_matrix(matrix_in);
566  this->_preconditioner->init();
567  }
568  }
569 
570  // Set the tolerances for the iterative solver. Use the user-supplied
571  // tolerance for the relative residual & leave the others at default values.
572  ierr = KSPSetTolerances (_ksp, tol, PETSC_DEFAULT,
573  PETSC_DEFAULT, max_its);
574  LIBMESH_CHKERR(ierr);
575 
576  // Allow command line options to override anything set programmatically.
577  ierr = KSPSetFromOptions(_ksp);
578  LIBMESH_CHKERR(ierr);
579 
580  // If the SolverConfiguration object is provided, use it to override
581  // solver options.
582  if(this->_solver_configuration)
583  {
584  this->_solver_configuration->configure_solver();
585  }
586 
587  // Solve the linear system
588  if (_restrict_solve_to_is != libmesh_nullptr)
589  {
590  ierr = KSPSolve (_ksp, subrhs, subsolution);
591  LIBMESH_CHKERR(ierr);
592  }
593  else
594  {
595  ierr = KSPSolve (_ksp, rhs->vec(), solution->vec());
596  LIBMESH_CHKERR(ierr);
597  }
598 
599  // Get the number of iterations required for convergence
600  ierr = KSPGetIterationNumber (_ksp, &its);
601  LIBMESH_CHKERR(ierr);
602 
603  // Get the norm of the final residual to return to the user.
604  ierr = KSPGetResidualNorm (_ksp, &final_resid);
605  LIBMESH_CHKERR(ierr);
606 
607  if (_restrict_solve_to_is != libmesh_nullptr)
608  {
609  switch(_subset_solve_mode)
610  {
611  case SUBSET_ZERO:
612  ierr = VecZeroEntries(solution->vec());
613  LIBMESH_CHKERR(ierr);
614  break;
615 
616  case SUBSET_COPY_RHS:
617  ierr = VecCopy(rhs->vec(),solution->vec());
618  LIBMESH_CHKERR(ierr);
619  break;
620 
621  case SUBSET_DONT_TOUCH:
622  /* Nothing to do here. */
623  break;
624 
625  default:
626  libmesh_error_msg("Invalid subset solve mode = " << _subset_solve_mode);
627  }
628  ierr = VecScatterBegin(scatter,subsolution,solution->vec(),INSERT_VALUES,SCATTER_REVERSE);
629  LIBMESH_CHKERR(ierr);
630  ierr = VecScatterEnd(scatter,subsolution,solution->vec(),INSERT_VALUES,SCATTER_REVERSE);
631  LIBMESH_CHKERR(ierr);
632 
633  ierr = LibMeshVecScatterDestroy(&scatter);
634  LIBMESH_CHKERR(ierr);
635 
636  if (this->_preconditioner)
637  {
638  // Before subprecond_matrix gets cleaned up, we should give
639  // the _preconditioner a different matrix.
640  this->_preconditioner->set_matrix(matrix_in);
641  this->_preconditioner->init();
642  }
643 
644  ierr = LibMeshVecDestroy(&subsolution);
645  LIBMESH_CHKERR(ierr);
646  ierr = LibMeshVecDestroy(&subrhs);
647  LIBMESH_CHKERR(ierr);
648  ierr = LibMeshMatDestroy(&submat);
649  LIBMESH_CHKERR(ierr);
650  ierr = LibMeshMatDestroy(&subprecond);
651  LIBMESH_CHKERR(ierr);
652  }
653 
654  // return the # of its. and the final residual norm.
655  return std::make_pair(its, final_resid);
656 }
657 
658 template <typename T>
659 std::pair<unsigned int, Real>
661  NumericVector<T> & solution_in,
662  NumericVector<T> & rhs_in,
663  const double tol,
664  const unsigned int m_its)
665 {
666  LOG_SCOPE("solve()", "PetscLinearSolver");
667 
668  // Make sure the data passed in are really of Petsc types
669  PetscMatrix<T> * matrix = cast_ptr<PetscMatrix<T> *>(&matrix_in);
670  // Note that the matrix and precond matrix are the same
671  PetscMatrix<T> * precond = cast_ptr<PetscMatrix<T> *>(&matrix_in);
672  PetscVector<T> * solution = cast_ptr<PetscVector<T> *>(&solution_in);
673  PetscVector<T> * rhs = cast_ptr<PetscVector<T> *>(&rhs_in);
674 
675  this->init (matrix);
676 
677  PetscErrorCode ierr=0;
678  PetscInt its=0, max_its = static_cast<PetscInt>(m_its);
679  PetscReal final_resid=0.;
680 
681  // Close the matrices and vectors in case this wasn't already done.
682  matrix->close ();
683  precond->close ();
684  solution->close ();
685  rhs->close ();
686 
687  Mat submat = libmesh_nullptr;
688  Mat subprecond = libmesh_nullptr;
689  Vec subrhs = libmesh_nullptr;
690  Vec subsolution = libmesh_nullptr;
691  VecScatter scatter = libmesh_nullptr;
692  UniquePtr<PetscMatrix<Number> > subprecond_matrix;
693 
694  // Set operators. Also restrict rhs and solution vector to
695  // subdomain if neccessary.
696  if (_restrict_solve_to_is != libmesh_nullptr)
697  {
698  PetscInt is_local_size = this->_restrict_solve_to_is_local_size();
699 
700  ierr = VecCreate(this->comm().get(),&subrhs);
701  LIBMESH_CHKERR(ierr);
702  ierr = VecSetSizes(subrhs,is_local_size,PETSC_DECIDE);
703  LIBMESH_CHKERR(ierr);
704  ierr = VecSetFromOptions(subrhs);
705  LIBMESH_CHKERR(ierr);
706 
707  ierr = VecCreate(this->comm().get(),&subsolution);
708  LIBMESH_CHKERR(ierr);
709  ierr = VecSetSizes(subsolution,is_local_size,PETSC_DECIDE);
710  LIBMESH_CHKERR(ierr);
711  ierr = VecSetFromOptions(subsolution);
712  LIBMESH_CHKERR(ierr);
713 
714  ierr = VecScatterCreate(rhs->vec(),_restrict_solve_to_is, subrhs, libmesh_nullptr, &scatter);
715  LIBMESH_CHKERR(ierr);
716 
717  ierr = VecScatterBegin(scatter,rhs->vec(),subrhs,INSERT_VALUES,SCATTER_FORWARD);
718  LIBMESH_CHKERR(ierr);
719  ierr = VecScatterEnd(scatter,rhs->vec(),subrhs,INSERT_VALUES,SCATTER_FORWARD);
720  LIBMESH_CHKERR(ierr);
721 
722  ierr = VecScatterBegin(scatter,solution->vec(),subsolution,INSERT_VALUES,SCATTER_FORWARD);
723  LIBMESH_CHKERR(ierr);
724  ierr = VecScatterEnd(scatter,solution->vec(),subsolution,INSERT_VALUES,SCATTER_FORWARD);
725  LIBMESH_CHKERR(ierr);
726 
727 #if PETSC_VERSION_LESS_THAN(3,1,0)
728  ierr = MatGetSubMatrix(matrix->mat(),
729  _restrict_solve_to_is,_restrict_solve_to_is,
730  PETSC_DECIDE,MAT_INITIAL_MATRIX,&submat);
731  LIBMESH_CHKERR(ierr);
732  ierr = MatGetSubMatrix(precond->mat(),
733  _restrict_solve_to_is,_restrict_solve_to_is,
734  PETSC_DECIDE,MAT_INITIAL_MATRIX,&subprecond);
735  LIBMESH_CHKERR(ierr);
736 #else
737  ierr = MatGetSubMatrix(matrix->mat(),
738  _restrict_solve_to_is,_restrict_solve_to_is,
739  MAT_INITIAL_MATRIX,&submat);
740  LIBMESH_CHKERR(ierr);
741  ierr = MatGetSubMatrix(precond->mat(),
742  _restrict_solve_to_is,_restrict_solve_to_is,
743  MAT_INITIAL_MATRIX,&subprecond);
744  LIBMESH_CHKERR(ierr);
745 #endif
746 
747  /* Since removing columns of the matrix changes the equation
748  system, we will now change the right hand side to compensate
749  for this. Note that this is not necessary if \p SUBSET_ZERO
750  has been selected. */
751  if(_subset_solve_mode!=SUBSET_ZERO)
752  {
753  _create_complement_is(rhs_in);
754  PetscInt is_complement_local_size =
755  cast_int<PetscInt>(rhs_in.local_size()-is_local_size);
756 
757  Vec subvec1 = libmesh_nullptr;
758  Mat submat1 = libmesh_nullptr;
759  VecScatter scatter1 = libmesh_nullptr;
760 
761  ierr = VecCreate(this->comm().get(),&subvec1);
762  LIBMESH_CHKERR(ierr);
763  ierr = VecSetSizes(subvec1,is_complement_local_size,PETSC_DECIDE);
764  LIBMESH_CHKERR(ierr);
765  ierr = VecSetFromOptions(subvec1);
766  LIBMESH_CHKERR(ierr);
767 
768  ierr = VecScatterCreate(rhs->vec(),_restrict_solve_to_is_complement, subvec1, libmesh_nullptr, &scatter1);
769  LIBMESH_CHKERR(ierr);
770 
771  ierr = VecScatterBegin(scatter1,_subset_solve_mode==SUBSET_COPY_RHS ? rhs->vec() : solution->vec(),subvec1,INSERT_VALUES,SCATTER_FORWARD);
772  LIBMESH_CHKERR(ierr);
773  ierr = VecScatterEnd(scatter1,_subset_solve_mode==SUBSET_COPY_RHS ? rhs->vec() : solution->vec(),subvec1,INSERT_VALUES,SCATTER_FORWARD);
774  LIBMESH_CHKERR(ierr);
775 
776  ierr = VecScale(subvec1,-1.0);
777  LIBMESH_CHKERR(ierr);
778 
779 #if PETSC_VERSION_LESS_THAN(3,1,0)
780  ierr = MatGetSubMatrix(matrix->mat(),
781  _restrict_solve_to_is,_restrict_solve_to_is_complement,
782  PETSC_DECIDE,MAT_INITIAL_MATRIX,&submat1);
783  LIBMESH_CHKERR(ierr);
784 #else
785  ierr = MatGetSubMatrix(matrix->mat(),
786  _restrict_solve_to_is,_restrict_solve_to_is_complement,
787  MAT_INITIAL_MATRIX,&submat1);
788  LIBMESH_CHKERR(ierr);
789 #endif
790 
791  ierr = MatMultAdd(submat1,subvec1,subrhs,subrhs);
792  LIBMESH_CHKERR(ierr);
793 
794  ierr = LibMeshVecScatterDestroy(&scatter1);
795  LIBMESH_CHKERR(ierr);
796  ierr = LibMeshVecDestroy(&subvec1);
797  LIBMESH_CHKERR(ierr);
798  ierr = LibMeshMatDestroy(&submat1);
799  LIBMESH_CHKERR(ierr);
800  }
801 #if PETSC_RELEASE_LESS_THAN(3,5,0)
802  ierr = KSPSetOperators(_ksp, submat, subprecond,
803  this->same_preconditioner ? SAME_PRECONDITIONER : DIFFERENT_NONZERO_PATTERN);
804 #else
805  ierr = KSPSetOperators(_ksp, submat, subprecond);
806 
807  PetscBool ksp_reuse_preconditioner = this->same_preconditioner ? PETSC_TRUE : PETSC_FALSE;
808  ierr = KSPSetReusePreconditioner(_ksp, ksp_reuse_preconditioner);
809 #endif
810  LIBMESH_CHKERR(ierr);
811 
812  if (this->_preconditioner)
813  {
814  subprecond_matrix.reset(new PetscMatrix<Number>(subprecond, this->comm()));
815  this->_preconditioner->set_matrix(*subprecond_matrix);
816  this->_preconditioner->init();
817  }
818  }
819  else
820  {
821 #if PETSC_RELEASE_LESS_THAN(3,5,0)
822  ierr = KSPSetOperators(_ksp, matrix->mat(), precond->mat(),
823  this->same_preconditioner ? SAME_PRECONDITIONER : DIFFERENT_NONZERO_PATTERN);
824 #else
825  ierr = KSPSetOperators(_ksp, matrix->mat(), precond->mat());
826 
827  PetscBool ksp_reuse_preconditioner = this->same_preconditioner ? PETSC_TRUE : PETSC_FALSE;
828  ierr = KSPSetReusePreconditioner(_ksp, ksp_reuse_preconditioner);
829 #endif
830  LIBMESH_CHKERR(ierr);
831 
832  if(this->_preconditioner)
833  {
834  this->_preconditioner->set_matrix(matrix_in);
835  this->_preconditioner->init();
836  }
837  }
838 
839  // Set the tolerances for the iterative solver. Use the user-supplied
840  // tolerance for the relative residual & leave the others at default values.
841  ierr = KSPSetTolerances (_ksp, tol, PETSC_DEFAULT,
842  PETSC_DEFAULT, max_its);
843  LIBMESH_CHKERR(ierr);
844 
845  // Allow command line options to override anything set programmatically.
846  ierr = KSPSetFromOptions(_ksp);
847  LIBMESH_CHKERR(ierr);
848 
849  // Solve the linear system
850  if (_restrict_solve_to_is != libmesh_nullptr)
851  {
852  ierr = KSPSolveTranspose (_ksp, subrhs, subsolution);
853  LIBMESH_CHKERR(ierr);
854  }
855  else
856  {
857  ierr = KSPSolveTranspose (_ksp, rhs->vec(), solution->vec());
858  LIBMESH_CHKERR(ierr);
859  }
860 
861  // Get the number of iterations required for convergence
862  ierr = KSPGetIterationNumber (_ksp, &its);
863  LIBMESH_CHKERR(ierr);
864 
865  // Get the norm of the final residual to return to the user.
866  ierr = KSPGetResidualNorm (_ksp, &final_resid);
867  LIBMESH_CHKERR(ierr);
868 
869  if (_restrict_solve_to_is != libmesh_nullptr)
870  {
871  switch(_subset_solve_mode)
872  {
873  case SUBSET_ZERO:
874  ierr = VecZeroEntries(solution->vec());
875  LIBMESH_CHKERR(ierr);
876  break;
877 
878  case SUBSET_COPY_RHS:
879  ierr = VecCopy(rhs->vec(),solution->vec());
880  LIBMESH_CHKERR(ierr);
881  break;
882 
883  case SUBSET_DONT_TOUCH:
884  /* Nothing to do here. */
885  break;
886 
887  default:
888  libmesh_error_msg("Invalid subset solve mode = " << _subset_solve_mode);
889  }
890  ierr = VecScatterBegin(scatter,subsolution,solution->vec(),INSERT_VALUES,SCATTER_REVERSE);
891  LIBMESH_CHKERR(ierr);
892  ierr = VecScatterEnd(scatter,subsolution,solution->vec(),INSERT_VALUES,SCATTER_REVERSE);
893  LIBMESH_CHKERR(ierr);
894 
895  ierr = LibMeshVecScatterDestroy(&scatter);
896  LIBMESH_CHKERR(ierr);
897 
898  if (this->_preconditioner)
899  {
900  // Before subprecond_matrix gets cleaned up, we should give
901  // the _preconditioner a different matrix.
902  this->_preconditioner->set_matrix(matrix_in);
903  this->_preconditioner->init();
904  }
905 
906  ierr = LibMeshVecDestroy(&subsolution);
907  LIBMESH_CHKERR(ierr);
908  ierr = LibMeshVecDestroy(&subrhs);
909  LIBMESH_CHKERR(ierr);
910  ierr = LibMeshMatDestroy(&submat);
911  LIBMESH_CHKERR(ierr);
912  ierr = LibMeshMatDestroy(&subprecond);
913  LIBMESH_CHKERR(ierr);
914  }
915 
916  // return the # of its. and the final residual norm.
917  return std::make_pair(its, final_resid);
918 }
919 
920 
921 template <typename T>
922 std::pair<unsigned int, Real>
924  NumericVector<T> & solution_in,
925  NumericVector<T> & rhs_in,
926  const double tol,
927  const unsigned int m_its)
928 {
929 
930 #if PETSC_VERSION_LESS_THAN(3,1,0)
931  if (_restrict_solve_to_is != libmesh_nullptr)
932  libmesh_error_msg("The current implementation of subset solves with " \
933  << "shell matrices requires PETSc version 3.1 or above. " \
934  << "Older PETSc version do not support automatic " \
935  << "submatrix generation of shell matrices.");
936 #endif
937 
938  LOG_SCOPE("solve()", "PetscLinearSolver");
939 
940  // Make sure the data passed in are really of Petsc types
941  PetscVector<T> * solution = cast_ptr<PetscVector<T> *>(&solution_in);
942  PetscVector<T> * rhs = cast_ptr<PetscVector<T> *>(&rhs_in);
943 
944  this->init ();
945 
946  PetscErrorCode ierr=0;
947  PetscInt its=0, max_its = static_cast<PetscInt>(m_its);
948  PetscReal final_resid=0.;
949 
950  Mat submat = libmesh_nullptr;
951  Vec subrhs = libmesh_nullptr;
952  Vec subsolution = libmesh_nullptr;
953  VecScatter scatter = libmesh_nullptr;
954 
955  // Close the matrices and vectors in case this wasn't already done.
956  solution->close ();
957  rhs->close ();
958 
959  // Prepare the matrix.
960  Mat mat;
961  ierr = MatCreateShell(this->comm().get(),
962  rhs_in.local_size(),
963  solution_in.local_size(),
964  rhs_in.size(),
965  solution_in.size(),
966  const_cast<void *>(static_cast<const void *>(&shell_matrix)),
967  &mat);
968  /* Note that the const_cast above is only necessary because PETSc
969  does not accept a const void *. Inside the member function
970  _petsc_shell_matrix() below, the pointer is casted back to a
971  const ShellMatrix<T> *. */
972 
973  LIBMESH_CHKERR(ierr);
974  ierr = MatShellSetOperation(mat,MATOP_MULT,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_mult));
975  LIBMESH_CHKERR(ierr);
976  ierr = MatShellSetOperation(mat,MATOP_MULT_ADD,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_mult_add));
977  LIBMESH_CHKERR(ierr);
978  ierr = MatShellSetOperation(mat,MATOP_GET_DIAGONAL,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_get_diagonal));
979  LIBMESH_CHKERR(ierr);
980 
981  // Restrict rhs and solution vectors and set operators. The input
982  // matrix works as the preconditioning matrix.
983  if (_restrict_solve_to_is != libmesh_nullptr)
984  {
985  PetscInt is_local_size = this->_restrict_solve_to_is_local_size();
986 
987  ierr = VecCreate(this->comm().get(),&subrhs);
988  LIBMESH_CHKERR(ierr);
989  ierr = VecSetSizes(subrhs,is_local_size,PETSC_DECIDE);
990  LIBMESH_CHKERR(ierr);
991  ierr = VecSetFromOptions(subrhs);
992  LIBMESH_CHKERR(ierr);
993 
994  ierr = VecCreate(this->comm().get(),&subsolution);
995  LIBMESH_CHKERR(ierr);
996  ierr = VecSetSizes(subsolution,is_local_size,PETSC_DECIDE);
997  LIBMESH_CHKERR(ierr);
998  ierr = VecSetFromOptions(subsolution);
999  LIBMESH_CHKERR(ierr);
1000 
1001  ierr = VecScatterCreate(rhs->vec(),_restrict_solve_to_is, subrhs, libmesh_nullptr, &scatter);
1002  LIBMESH_CHKERR(ierr);
1003 
1004  ierr = VecScatterBegin(scatter,rhs->vec(),subrhs,INSERT_VALUES,SCATTER_FORWARD);
1005  LIBMESH_CHKERR(ierr);
1006  ierr = VecScatterEnd(scatter,rhs->vec(),subrhs,INSERT_VALUES,SCATTER_FORWARD);
1007  LIBMESH_CHKERR(ierr);
1008 
1009  ierr = VecScatterBegin(scatter,solution->vec(),subsolution,INSERT_VALUES,SCATTER_FORWARD);
1010  LIBMESH_CHKERR(ierr);
1011  ierr = VecScatterEnd(scatter,solution->vec(),subsolution,INSERT_VALUES,SCATTER_FORWARD);
1012  LIBMESH_CHKERR(ierr);
1013 
1014 #if PETSC_VERSION_LESS_THAN(3,1,0)
1015  /* This point can't be reached, see above. */
1016  libmesh_assert(false);
1017 #else
1018  ierr = MatGetSubMatrix(mat,
1019  _restrict_solve_to_is,_restrict_solve_to_is,
1020  MAT_INITIAL_MATRIX,&submat);
1021  LIBMESH_CHKERR(ierr);
1022 #endif
1023 
1024  /* Since removing columns of the matrix changes the equation
1025  system, we will now change the right hand side to compensate
1026  for this. Note that this is not necessary if \p SUBSET_ZERO
1027  has been selected. */
1028  if(_subset_solve_mode!=SUBSET_ZERO)
1029  {
1030  _create_complement_is(rhs_in);
1031  PetscInt is_complement_local_size =
1032  cast_int<PetscInt>(rhs_in.local_size()-is_local_size);
1033 
1034  Vec subvec1 = libmesh_nullptr;
1035  Mat submat1 = libmesh_nullptr;
1036  VecScatter scatter1 = libmesh_nullptr;
1037 
1038  ierr = VecCreate(this->comm().get(),&subvec1);
1039  LIBMESH_CHKERR(ierr);
1040  ierr = VecSetSizes(subvec1,is_complement_local_size,PETSC_DECIDE);
1041  LIBMESH_CHKERR(ierr);
1042  ierr = VecSetFromOptions(subvec1);
1043  LIBMESH_CHKERR(ierr);
1044 
1045  ierr = VecScatterCreate(rhs->vec(),_restrict_solve_to_is_complement, subvec1, libmesh_nullptr, &scatter1);
1046  LIBMESH_CHKERR(ierr);
1047 
1048  ierr = VecScatterBegin(scatter1,_subset_solve_mode==SUBSET_COPY_RHS ? rhs->vec() : solution->vec(),subvec1,INSERT_VALUES,SCATTER_FORWARD);
1049  LIBMESH_CHKERR(ierr);
1050  ierr = VecScatterEnd(scatter1,_subset_solve_mode==SUBSET_COPY_RHS ? rhs->vec() : solution->vec(),subvec1,INSERT_VALUES,SCATTER_FORWARD);
1051  LIBMESH_CHKERR(ierr);
1052 
1053  ierr = VecScale(subvec1,-1.0);
1054  LIBMESH_CHKERR(ierr);
1055 
1056 #if PETSC_VERSION_LESS_THAN(3,1,0)
1057  /* This point can't be reached, see above. */
1058  libmesh_assert(false);
1059 #else
1060  ierr = MatGetSubMatrix(mat,
1061  _restrict_solve_to_is,_restrict_solve_to_is_complement,
1062  MAT_INITIAL_MATRIX,&submat1);
1063  LIBMESH_CHKERR(ierr);
1064 #endif
1065 
1066  // The following lines would be correct, but don't work
1067  // correctly in PETSc up to 3.1.0-p5. See discussion in
1068  // petsc-users of Nov 9, 2010.
1069  //
1070  // ierr = MatMultAdd(submat1,subvec1,subrhs,subrhs);
1071  // LIBMESH_CHKERR(ierr);
1072  //
1073  // We workaround by using a temporary vector. Note that the
1074  // fix in PETsc 3.1.0-p6 uses a temporary vector internally,
1075  // so this is no effective performance loss.
1076  Vec subvec2 = libmesh_nullptr;
1077  ierr = VecCreate(this->comm().get(),&subvec2);
1078  LIBMESH_CHKERR(ierr);
1079  ierr = VecSetSizes(subvec2,is_local_size,PETSC_DECIDE);
1080  LIBMESH_CHKERR(ierr);
1081  ierr = VecSetFromOptions(subvec2);
1082  LIBMESH_CHKERR(ierr);
1083  ierr = MatMult(submat1,subvec1,subvec2);
1084  LIBMESH_CHKERR(ierr);
1085  ierr = VecAXPY(subrhs,1.0,subvec2);
1086 
1087  ierr = LibMeshVecScatterDestroy(&scatter1);
1088  LIBMESH_CHKERR(ierr);
1089  ierr = LibMeshVecDestroy(&subvec1);
1090  LIBMESH_CHKERR(ierr);
1091  ierr = LibMeshMatDestroy(&submat1);
1092  LIBMESH_CHKERR(ierr);
1093  }
1094 #if PETSC_RELEASE_LESS_THAN(3,5,0)
1095  ierr = KSPSetOperators(_ksp, submat, submat,
1096  DIFFERENT_NONZERO_PATTERN);
1097 #else
1098  ierr = KSPSetOperators(_ksp, submat, submat);
1099 #endif
1100  LIBMESH_CHKERR(ierr);
1101  }
1102  else
1103  {
1104 #if PETSC_RELEASE_LESS_THAN(3,5,0)
1105  ierr = KSPSetOperators(_ksp, mat, mat,
1106  DIFFERENT_NONZERO_PATTERN);
1107 #else
1108  ierr = KSPSetOperators(_ksp, mat, mat);
1109 #endif
1110  LIBMESH_CHKERR(ierr);
1111  }
1112 
1113  // Set the tolerances for the iterative solver. Use the user-supplied
1114  // tolerance for the relative residual & leave the others at default values.
1115  ierr = KSPSetTolerances (_ksp, tol, PETSC_DEFAULT,
1116  PETSC_DEFAULT, max_its);
1117  LIBMESH_CHKERR(ierr);
1118 
1119  // Allow command line options to override anything set programmatically.
1120  ierr = KSPSetFromOptions(_ksp);
1121  LIBMESH_CHKERR(ierr);
1122 
1123  // Solve the linear system
1124  if (_restrict_solve_to_is != libmesh_nullptr)
1125  {
1126  ierr = KSPSolve (_ksp, subrhs, subsolution);
1127  LIBMESH_CHKERR(ierr);
1128  }
1129  else
1130  {
1131  ierr = KSPSolve (_ksp, rhs->vec(), solution->vec());
1132  LIBMESH_CHKERR(ierr);
1133  }
1134 
1135  // Get the number of iterations required for convergence
1136  ierr = KSPGetIterationNumber (_ksp, &its);
1137  LIBMESH_CHKERR(ierr);
1138 
1139  // Get the norm of the final residual to return to the user.
1140  ierr = KSPGetResidualNorm (_ksp, &final_resid);
1141  LIBMESH_CHKERR(ierr);
1142 
1143  if (_restrict_solve_to_is != libmesh_nullptr)
1144  {
1145  switch(_subset_solve_mode)
1146  {
1147  case SUBSET_ZERO:
1148  ierr = VecZeroEntries(solution->vec());
1149  LIBMESH_CHKERR(ierr);
1150  break;
1151 
1152  case SUBSET_COPY_RHS:
1153  ierr = VecCopy(rhs->vec(),solution->vec());
1154  LIBMESH_CHKERR(ierr);
1155  break;
1156 
1157  case SUBSET_DONT_TOUCH:
1158  /* Nothing to do here. */
1159  break;
1160 
1161  default:
1162  libmesh_error_msg("Invalid subset solve mode = " << _subset_solve_mode);
1163  }
1164  ierr = VecScatterBegin(scatter,subsolution,solution->vec(),INSERT_VALUES,SCATTER_REVERSE);
1165  LIBMESH_CHKERR(ierr);
1166  ierr = VecScatterEnd(scatter,subsolution,solution->vec(),INSERT_VALUES,SCATTER_REVERSE);
1167  LIBMESH_CHKERR(ierr);
1168 
1169  ierr = LibMeshVecScatterDestroy(&scatter);
1170  LIBMESH_CHKERR(ierr);
1171 
1172  ierr = LibMeshVecDestroy(&subsolution);
1173  LIBMESH_CHKERR(ierr);
1174  ierr = LibMeshVecDestroy(&subrhs);
1175  LIBMESH_CHKERR(ierr);
1176  ierr = LibMeshMatDestroy(&submat);
1177  LIBMESH_CHKERR(ierr);
1178  }
1179 
1180  // Destroy the matrix.
1181  ierr = LibMeshMatDestroy(&mat);
1182  LIBMESH_CHKERR(ierr);
1183 
1184  // return the # of its. and the final residual norm.
1185  return std::make_pair(its, final_resid);
1186 }
1187 
1188 
1189 
1190 template <typename T>
1191 std::pair<unsigned int, Real>
1193  const SparseMatrix<T> & precond_matrix,
1194  NumericVector<T> & solution_in,
1195  NumericVector<T> & rhs_in,
1196  const double tol,
1197  const unsigned int m_its)
1198 {
1199 
1200 #if PETSC_VERSION_LESS_THAN(3,1,0)
1201  if (_restrict_solve_to_is != libmesh_nullptr)
1202  libmesh_error_msg("The current implementation of subset solves with " \
1203  << "shell matrices requires PETSc version 3.1 or above. " \
1204  << "Older PETSc version do not support automatic " \
1205  << "submatrix generation of shell matrices.");
1206 #endif
1207 
1208  LOG_SCOPE("solve()", "PetscLinearSolver");
1209 
1210  // Make sure the data passed in are really of Petsc types
1211  const PetscMatrix<T> * precond = cast_ptr<const PetscMatrix<T> *>(&precond_matrix);
1212  PetscVector<T> * solution = cast_ptr<PetscVector<T> *>(&solution_in);
1213  PetscVector<T> * rhs = cast_ptr<PetscVector<T> *>(&rhs_in);
1214 
1215  this->init ();
1216 
1217  PetscErrorCode ierr=0;
1218  PetscInt its=0, max_its = static_cast<PetscInt>(m_its);
1219  PetscReal final_resid=0.;
1220 
1221  Mat submat = libmesh_nullptr;
1222  Mat subprecond = libmesh_nullptr;
1223  Vec subrhs = libmesh_nullptr;
1224  Vec subsolution = libmesh_nullptr;
1225  VecScatter scatter = libmesh_nullptr;
1226  UniquePtr<PetscMatrix<Number> > subprecond_matrix;
1227 
1228  // Close the matrices and vectors in case this wasn't already done.
1229  solution->close ();
1230  rhs->close ();
1231 
1232  // Prepare the matrix.
1233  Mat mat;
1234  ierr = MatCreateShell(this->comm().get(),
1235  rhs_in.local_size(),
1236  solution_in.local_size(),
1237  rhs_in.size(),
1238  solution_in.size(),
1239  const_cast<void *>(static_cast<const void *>(&shell_matrix)),
1240  &mat);
1241  /* Note that the const_cast above is only necessary because PETSc
1242  does not accept a const void *. Inside the member function
1243  _petsc_shell_matrix() below, the pointer is casted back to a
1244  const ShellMatrix<T> *. */
1245 
1246  LIBMESH_CHKERR(ierr);
1247  ierr = MatShellSetOperation(mat,MATOP_MULT,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_mult));
1248  LIBMESH_CHKERR(ierr);
1249  ierr = MatShellSetOperation(mat,MATOP_MULT_ADD,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_mult_add));
1250  LIBMESH_CHKERR(ierr);
1251  ierr = MatShellSetOperation(mat,MATOP_GET_DIAGONAL,reinterpret_cast<void(*)(void)>(_petsc_shell_matrix_get_diagonal));
1252  LIBMESH_CHKERR(ierr);
1253 
1254  // Restrict rhs and solution vectors and set operators. The input
1255  // matrix works as the preconditioning matrix.
1256  if (_restrict_solve_to_is != libmesh_nullptr)
1257  {
1258  PetscInt is_local_size = this->_restrict_solve_to_is_local_size();
1259 
1260  ierr = VecCreate(this->comm().get(),&subrhs);
1261  LIBMESH_CHKERR(ierr);
1262  ierr = VecSetSizes(subrhs,is_local_size,PETSC_DECIDE);
1263  LIBMESH_CHKERR(ierr);
1264  ierr = VecSetFromOptions(subrhs);
1265  LIBMESH_CHKERR(ierr);
1266 
1267  ierr = VecCreate(this->comm().get(),&subsolution);
1268  LIBMESH_CHKERR(ierr);
1269  ierr = VecSetSizes(subsolution,is_local_size,PETSC_DECIDE);
1270  LIBMESH_CHKERR(ierr);
1271  ierr = VecSetFromOptions(subsolution);
1272  LIBMESH_CHKERR(ierr);
1273 
1274  ierr = VecScatterCreate(rhs->vec(),_restrict_solve_to_is, subrhs, libmesh_nullptr, &scatter);
1275  LIBMESH_CHKERR(ierr);
1276 
1277  ierr = VecScatterBegin(scatter,rhs->vec(),subrhs,INSERT_VALUES,SCATTER_FORWARD);
1278  LIBMESH_CHKERR(ierr);
1279  ierr = VecScatterEnd(scatter,rhs->vec(),subrhs,INSERT_VALUES,SCATTER_FORWARD);
1280  LIBMESH_CHKERR(ierr);
1281 
1282  ierr = VecScatterBegin(scatter,solution->vec(),subsolution,INSERT_VALUES,SCATTER_FORWARD);
1283  LIBMESH_CHKERR(ierr);
1284  ierr = VecScatterEnd(scatter,solution->vec(),subsolution,INSERT_VALUES,SCATTER_FORWARD);
1285  LIBMESH_CHKERR(ierr);
1286 
1287 #if PETSC_VERSION_LESS_THAN(3,1,0)
1288  /* This point can't be reached, see above. */
1289  libmesh_assert(false);
1290 #else
1291  ierr = MatGetSubMatrix(mat,
1292  _restrict_solve_to_is,_restrict_solve_to_is,
1293  MAT_INITIAL_MATRIX,&submat);
1294  LIBMESH_CHKERR(ierr);
1295  ierr = MatGetSubMatrix(const_cast<PetscMatrix<T> *>(precond)->mat(),
1296  _restrict_solve_to_is,_restrict_solve_to_is,
1297  MAT_INITIAL_MATRIX,&subprecond);
1298  LIBMESH_CHKERR(ierr);
1299 #endif
1300 
1301  /* Since removing columns of the matrix changes the equation
1302  system, we will now change the right hand side to compensate
1303  for this. Note that this is not necessary if \p SUBSET_ZERO
1304  has been selected. */
1305  if(_subset_solve_mode!=SUBSET_ZERO)
1306  {
1307  _create_complement_is(rhs_in);
1308  PetscInt is_complement_local_size = rhs_in.local_size()-is_local_size;
1309 
1310  Vec subvec1 = libmesh_nullptr;
1311  Mat submat1 = libmesh_nullptr;
1312  VecScatter scatter1 = libmesh_nullptr;
1313 
1314  ierr = VecCreate(this->comm().get(),&subvec1);
1315  LIBMESH_CHKERR(ierr);
1316  ierr = VecSetSizes(subvec1,is_complement_local_size,PETSC_DECIDE);
1317  LIBMESH_CHKERR(ierr);
1318  ierr = VecSetFromOptions(subvec1);
1319  LIBMESH_CHKERR(ierr);
1320 
1321  ierr = VecScatterCreate(rhs->vec(),_restrict_solve_to_is_complement, subvec1, libmesh_nullptr, &scatter1);
1322  LIBMESH_CHKERR(ierr);
1323 
1324  ierr = VecScatterBegin(scatter1,_subset_solve_mode==SUBSET_COPY_RHS ? rhs->vec() : solution->vec(),subvec1,INSERT_VALUES,SCATTER_FORWARD);
1325  LIBMESH_CHKERR(ierr);
1326  ierr = VecScatterEnd(scatter1,_subset_solve_mode==SUBSET_COPY_RHS ? rhs->vec() : solution->vec(),subvec1,INSERT_VALUES,SCATTER_FORWARD);
1327  LIBMESH_CHKERR(ierr);
1328 
1329  ierr = VecScale(subvec1,-1.0);
1330  LIBMESH_CHKERR(ierr);
1331 
1332 #if PETSC_VERSION_LESS_THAN(3,1,0)
1333  /* This point can't be reached, see above. */
1334  libmesh_assert(false);
1335 #else
1336  ierr = MatGetSubMatrix(mat,
1337  _restrict_solve_to_is,_restrict_solve_to_is_complement,
1338  MAT_INITIAL_MATRIX,&submat1);
1339  LIBMESH_CHKERR(ierr);
1340 #endif
1341 
1342  // The following lines would be correct, but don't work
1343  // correctly in PETSc up to 3.1.0-p5. See discussion in
1344  // petsc-users of Nov 9, 2010.
1345  //
1346  // ierr = MatMultAdd(submat1,subvec1,subrhs,subrhs);
1347  // LIBMESH_CHKERR(ierr);
1348  //
1349  // We workaround by using a temporary vector. Note that the
1350  // fix in PETsc 3.1.0-p6 uses a temporary vector internally,
1351  // so this is no effective performance loss.
1352  Vec subvec2 = libmesh_nullptr;
1353  ierr = VecCreate(this->comm().get(),&subvec2);
1354  LIBMESH_CHKERR(ierr);
1355  ierr = VecSetSizes(subvec2,is_local_size,PETSC_DECIDE);
1356  LIBMESH_CHKERR(ierr);
1357  ierr = VecSetFromOptions(subvec2);
1358  LIBMESH_CHKERR(ierr);
1359  ierr = MatMult(submat1,subvec1,subvec2);
1360  LIBMESH_CHKERR(ierr);
1361  ierr = VecAXPY(subrhs,1.0,subvec2);
1362  LIBMESH_CHKERR(ierr);
1363 
1364  ierr = LibMeshVecScatterDestroy(&scatter1);
1365  LIBMESH_CHKERR(ierr);
1366  ierr = LibMeshVecDestroy(&subvec1);
1367  LIBMESH_CHKERR(ierr);
1368  ierr = LibMeshMatDestroy(&submat1);
1369  LIBMESH_CHKERR(ierr);
1370  }
1371 
1372 #if PETSC_RELEASE_LESS_THAN(3,5,0)
1373  ierr = KSPSetOperators(_ksp, submat, subprecond,
1374  DIFFERENT_NONZERO_PATTERN);
1375 #else
1376  ierr = KSPSetOperators(_ksp, submat, subprecond);
1377 #endif
1378  LIBMESH_CHKERR(ierr);
1379 
1380  if (this->_preconditioner)
1381  {
1382  subprecond_matrix.reset(new PetscMatrix<Number>(subprecond, this->comm()));
1383  this->_preconditioner->set_matrix(*subprecond_matrix);
1384  this->_preconditioner->init();
1385  }
1386  }
1387  else
1388  {
1389 #if PETSC_RELEASE_LESS_THAN(3,5,0)
1390  ierr = KSPSetOperators(_ksp, mat, const_cast<PetscMatrix<T> *>(precond)->mat(),
1391  DIFFERENT_NONZERO_PATTERN);
1392 #else
1393  ierr = KSPSetOperators(_ksp, mat, const_cast<PetscMatrix<T> *>(precond)->mat());
1394 #endif
1395  LIBMESH_CHKERR(ierr);
1396 
1397  if(this->_preconditioner)
1398  {
1399  this->_preconditioner->set_matrix(const_cast<SparseMatrix<Number> &>(precond_matrix));
1400  this->_preconditioner->init();
1401  }
1402  }
1403 
1404  // Set the tolerances for the iterative solver. Use the user-supplied
1405  // tolerance for the relative residual & leave the others at default values.
1406  ierr = KSPSetTolerances (_ksp, tol, PETSC_DEFAULT,
1407  PETSC_DEFAULT, max_its);
1408  LIBMESH_CHKERR(ierr);
1409 
1410  // Allow command line options to override anything set programmatically.
1411  ierr = KSPSetFromOptions(_ksp);
1412  LIBMESH_CHKERR(ierr);
1413 
1414  // Solve the linear system
1415  if (_restrict_solve_to_is != libmesh_nullptr)
1416  {
1417  ierr = KSPSolve (_ksp, subrhs, subsolution);
1418  LIBMESH_CHKERR(ierr);
1419  }
1420  else
1421  {
1422  ierr = KSPSolve (_ksp, rhs->vec(), solution->vec());
1423  LIBMESH_CHKERR(ierr);
1424  }
1425 
1426  // Get the number of iterations required for convergence
1427  ierr = KSPGetIterationNumber (_ksp, &its);
1428  LIBMESH_CHKERR(ierr);
1429 
1430  // Get the norm of the final residual to return to the user.
1431  ierr = KSPGetResidualNorm (_ksp, &final_resid);
1432  LIBMESH_CHKERR(ierr);
1433 
1434  if (_restrict_solve_to_is != libmesh_nullptr)
1435  {
1436  switch(_subset_solve_mode)
1437  {
1438  case SUBSET_ZERO:
1439  ierr = VecZeroEntries(solution->vec());
1440  LIBMESH_CHKERR(ierr);
1441  break;
1442 
1443  case SUBSET_COPY_RHS:
1444  ierr = VecCopy(rhs->vec(),solution->vec());
1445  LIBMESH_CHKERR(ierr);
1446  break;
1447 
1448  case SUBSET_DONT_TOUCH:
1449  /* Nothing to do here. */
1450  break;
1451 
1452  default:
1453  libmesh_error_msg("Invalid subset solve mode = " << _subset_solve_mode);
1454  }
1455  ierr = VecScatterBegin(scatter,subsolution,solution->vec(),INSERT_VALUES,SCATTER_REVERSE);
1456  LIBMESH_CHKERR(ierr);
1457  ierr = VecScatterEnd(scatter,subsolution,solution->vec(),INSERT_VALUES,SCATTER_REVERSE);
1458  LIBMESH_CHKERR(ierr);
1459 
1460  ierr = LibMeshVecScatterDestroy(&scatter);
1461  LIBMESH_CHKERR(ierr);
1462 
1463  if (this->_preconditioner)
1464  {
1465  // Before subprecond_matrix gets cleaned up, we should give
1466  // the _preconditioner a different matrix.
1467  this->_preconditioner->set_matrix(const_cast<SparseMatrix<Number> &>(precond_matrix));
1468  this->_preconditioner->init();
1469  }
1470 
1471  ierr = LibMeshVecDestroy(&subsolution);
1472  LIBMESH_CHKERR(ierr);
1473  ierr = LibMeshVecDestroy(&subrhs);
1474  LIBMESH_CHKERR(ierr);
1475  ierr = LibMeshMatDestroy(&submat);
1476  LIBMESH_CHKERR(ierr);
1477  ierr = LibMeshMatDestroy(&subprecond);
1478  LIBMESH_CHKERR(ierr);
1479  }
1480 
1481  // Destroy the matrix.
1482  ierr = LibMeshMatDestroy(&mat);
1483  LIBMESH_CHKERR(ierr);
1484 
1485  // return the # of its. and the final residual norm.
1486  return std::make_pair(its, final_resid);
1487 }
1488 
1489 
1490 
1491 template <typename T>
1492 void PetscLinearSolver<T>::get_residual_history(std::vector<double> & hist)
1493 {
1494  PetscErrorCode ierr = 0;
1495  PetscInt its = 0;
1496 
1497  // Fill the residual history vector with the residual norms
1498  // Note that GetResidualHistory() does not copy any values, it
1499  // simply sets the pointer p. Note that for some Krylov subspace
1500  // methods, the number of residuals returned in the history
1501  // vector may be different from what you are expecting. For
1502  // example, TFQMR returns two residual values per iteration step.
1503  PetscReal * p;
1504  ierr = KSPGetResidualHistory(_ksp, &p, &its);
1505  LIBMESH_CHKERR(ierr);
1506 
1507  // Check for early return
1508  if (its == 0) return;
1509 
1510  // Create space to store the result
1511  hist.resize(its);
1512 
1513  // Copy history into the vector provided by the user.
1514  for (PetscInt i=0; i<its; ++i)
1515  {
1516  hist[i] = *p;
1517  p++;
1518  }
1519 }
1520 
1521 
1522 
1523 
1524 template <typename T>
1526 {
1527  PetscErrorCode ierr = 0;
1528  PetscInt its = 0;
1529 
1530  // Fill the residual history vector with the residual norms
1531  // Note that GetResidualHistory() does not copy any values, it
1532  // simply sets the pointer p. Note that for some Krylov subspace
1533  // methods, the number of residuals returned in the history
1534  // vector may be different from what you are expecting. For
1535  // example, TFQMR returns two residual values per iteration step.
1536  PetscReal * p;
1537  ierr = KSPGetResidualHistory(_ksp, &p, &its);
1538  LIBMESH_CHKERR(ierr);
1539 
1540  // Check no residual history
1541  if (its == 0)
1542  {
1543  libMesh::err << "No iterations have been performed, returning 0." << std::endl;
1544  return 0.;
1545  }
1546 
1547  // Otherwise, return the value pointed to by p.
1548  return *p;
1549 }
1550 
1551 
1552 
1553 
1554 template <typename T>
1556 {
1557  PetscErrorCode ierr = 0;
1558 
1559  switch (this->_solver_type)
1560  {
1561 
1562  case CG:
1563  ierr = KSPSetType (_ksp, const_cast<KSPType>(KSPCG));
1564  LIBMESH_CHKERR(ierr);
1565  return;
1566 
1567  case CR:
1568  ierr = KSPSetType (_ksp, const_cast<KSPType>(KSPCR));
1569  LIBMESH_CHKERR(ierr);
1570  return;
1571 
1572  case CGS:
1573  ierr = KSPSetType (_ksp, const_cast<KSPType>(KSPCGS));
1574  LIBMESH_CHKERR(ierr);
1575  return;
1576 
1577  case BICG:
1578  ierr = KSPSetType (_ksp, const_cast<KSPType>(KSPBICG));
1579  LIBMESH_CHKERR(ierr);
1580  return;
1581 
1582  case TCQMR:
1583  ierr = KSPSetType (_ksp, const_cast<KSPType>(KSPTCQMR));
1584  LIBMESH_CHKERR(ierr);
1585  return;
1586 
1587  case TFQMR:
1588  ierr = KSPSetType (_ksp, const_cast<KSPType>(KSPTFQMR));
1589  LIBMESH_CHKERR(ierr);
1590  return;
1591 
1592  case LSQR:
1593  ierr = KSPSetType (_ksp, const_cast<KSPType>(KSPLSQR));
1594  LIBMESH_CHKERR(ierr);
1595  return;
1596 
1597  case BICGSTAB:
1598  ierr = KSPSetType (_ksp, const_cast<KSPType>(KSPBCGS));
1599  LIBMESH_CHKERR(ierr);
1600  return;
1601 
1602  case MINRES:
1603  ierr = KSPSetType (_ksp, const_cast<KSPType>(KSPMINRES));
1604  LIBMESH_CHKERR(ierr);
1605  return;
1606 
1607  case GMRES:
1608  ierr = KSPSetType (_ksp, const_cast<KSPType>(KSPGMRES));
1609  LIBMESH_CHKERR(ierr);
1610  return;
1611 
1612  case RICHARDSON:
1613  ierr = KSPSetType (_ksp, const_cast<KSPType>(KSPRICHARDSON));
1614  LIBMESH_CHKERR(ierr);
1615  return;
1616 
1617  case CHEBYSHEV:
1618 #if defined(LIBMESH_HAVE_PETSC) && PETSC_VERSION_LESS_THAN(3,3,0)
1619  ierr = KSPSetType (_ksp, const_cast<KSPType>(KSPCHEBYCHEV));
1620  LIBMESH_CHKERR(ierr);
1621  return;
1622 #else
1623  ierr = KSPSetType (_ksp, const_cast<KSPType>(KSPCHEBYSHEV));
1624  LIBMESH_CHKERR(ierr);
1625  return;
1626 #endif
1627 
1628 
1629  default:
1630  libMesh::err << "ERROR: Unsupported PETSC Solver: "
1631  << Utility::enum_to_string(this->_solver_type) << std::endl
1632  << "Continuing with PETSC defaults" << std::endl;
1633  }
1634 }
1635 
1636 
1637 
1638 template <typename T>
1640 {
1641  KSPConvergedReason reason;
1642  KSPGetConvergedReason(_ksp, &reason);
1643 
1644  switch(reason)
1645  {
1646 #if !PETSC_VERSION_LESS_THAN(3,2,0)
1647  case KSP_CONVERGED_RTOL_NORMAL : return CONVERGED_RTOL_NORMAL;
1648  case KSP_CONVERGED_ATOL_NORMAL : return CONVERGED_ATOL_NORMAL;
1649 #endif
1650  case KSP_CONVERGED_RTOL : return CONVERGED_RTOL;
1651  case KSP_CONVERGED_ATOL : return CONVERGED_ATOL;
1652  case KSP_CONVERGED_ITS : return CONVERGED_ITS;
1653  case KSP_CONVERGED_CG_NEG_CURVE : return CONVERGED_CG_NEG_CURVE;
1654  case KSP_CONVERGED_CG_CONSTRAINED : return CONVERGED_CG_CONSTRAINED;
1655  case KSP_CONVERGED_STEP_LENGTH : return CONVERGED_STEP_LENGTH;
1656  case KSP_CONVERGED_HAPPY_BREAKDOWN : return CONVERGED_HAPPY_BREAKDOWN;
1657  case KSP_DIVERGED_NULL : return DIVERGED_NULL;
1658  case KSP_DIVERGED_ITS : return DIVERGED_ITS;
1659  case KSP_DIVERGED_DTOL : return DIVERGED_DTOL;
1660  case KSP_DIVERGED_BREAKDOWN : return DIVERGED_BREAKDOWN;
1661  case KSP_DIVERGED_BREAKDOWN_BICG : return DIVERGED_BREAKDOWN_BICG;
1662  case KSP_DIVERGED_NONSYMMETRIC : return DIVERGED_NONSYMMETRIC;
1663  case KSP_DIVERGED_INDEFINITE_PC : return DIVERGED_INDEFINITE_PC;
1664 #if PETSC_VERSION_LESS_THAN(3,4,0)
1665  case KSP_DIVERGED_NAN : return DIVERGED_NAN;
1666 #else
1667  case KSP_DIVERGED_NANORINF : return DIVERGED_NAN;
1668 #endif
1669  case KSP_DIVERGED_INDEFINITE_MAT : return DIVERGED_INDEFINITE_MAT;
1670  case KSP_CONVERGED_ITERATING : return CONVERGED_ITERATING;
1671 #if !PETSC_VERSION_LESS_THAN(3,7,0)
1672  case KSP_DIVERGED_PCSETUP_FAILED : return DIVERGED_PCSETUP_FAILED;
1673 #endif
1674  default :
1675  libMesh::err << "Unknown convergence flag!" << std::endl;
1676  return UNKNOWN_FLAG;
1677  }
1678 }
1679 
1680 
1681 template <typename T>
1682 PetscErrorCode PetscLinearSolver<T>::_petsc_shell_matrix_mult(Mat mat, Vec arg, Vec dest)
1683 {
1684  /* Get the matrix context. */
1685  PetscErrorCode ierr=0;
1686  void * ctx;
1687  ierr = MatShellGetContext(mat,&ctx);
1688 
1689  /* Get user shell matrix object. */
1690  const ShellMatrix<T> & shell_matrix = *static_cast<const ShellMatrix<T> *>(ctx);
1691  CHKERRABORT(shell_matrix.comm().get(), ierr);
1692 
1693  /* Make \p NumericVector instances around the vectors. */
1694  PetscVector<T> arg_global(arg, shell_matrix.comm());
1695  PetscVector<T> dest_global(dest, shell_matrix.comm());
1696 
1697  /* Call the user function. */
1698  shell_matrix.vector_mult(dest_global,arg_global);
1699 
1700  return ierr;
1701 }
1702 
1703 
1704 
1705 template <typename T>
1706 PetscErrorCode PetscLinearSolver<T>::_petsc_shell_matrix_mult_add(Mat mat, Vec arg, Vec add, Vec dest)
1707 {
1708  /* Get the matrix context. */
1709  PetscErrorCode ierr=0;
1710  void * ctx;
1711  ierr = MatShellGetContext(mat,&ctx);
1712 
1713  /* Get user shell matrix object. */
1714  const ShellMatrix<T> & shell_matrix = *static_cast<const ShellMatrix<T> *>(ctx);
1715  CHKERRABORT(shell_matrix.comm().get(), ierr);
1716 
1717  /* Make \p NumericVector instances around the vectors. */
1718  PetscVector<T> arg_global(arg, shell_matrix.comm());
1719  PetscVector<T> dest_global(dest, shell_matrix.comm());
1720  PetscVector<T> add_global(add, shell_matrix.comm());
1721 
1722  if(add!=arg)
1723  {
1724  arg_global = add_global;
1725  }
1726 
1727  /* Call the user function. */
1728  shell_matrix.vector_mult_add(dest_global,arg_global);
1729 
1730  return ierr;
1731 }
1732 
1733 
1734 
1735 template <typename T>
1737 {
1738  /* Get the matrix context. */
1739  PetscErrorCode ierr=0;
1740  void * ctx;
1741  ierr = MatShellGetContext(mat,&ctx);
1742 
1743  /* Get user shell matrix object. */
1744  const ShellMatrix<T> & shell_matrix = *static_cast<const ShellMatrix<T> *>(ctx);
1745  CHKERRABORT(shell_matrix.comm().get(), ierr);
1746 
1747  /* Make \p NumericVector instances around the vector. */
1748  PetscVector<T> dest_global(dest, shell_matrix.comm());
1749 
1750  /* Call the user function. */
1751  shell_matrix.get_diagonal(dest_global);
1752 
1753  return ierr;
1754 }
1755 
1756 
1757 
1758 //------------------------------------------------------------------
1759 // Explicit instantiations
1760 template class PetscLinearSolver<Number>;
1761 
1762 } // namespace libMesh
1763 
1764 
1765 
1766 #endif // #ifdef LIBMESH_HAVE_PETSC
std::string name(const ElemQuality q)
Definition: elem_quality.C:39
virtual void get_diagonal(NumericVector< T > &dest) const =0
virtual void vector_mult(NumericVector< T > &dest, const NumericVector< T > &arg) const =0
void petsc_auto_fieldsplit(PC my_pc, const System &sys)
NumericVector interface to PETSc Vec.
Definition: petsc_vector.h:64
void get_residual_history(std::vector< double > &hist)
virtual std::pair< unsigned int, Real > solve(SparseMatrix< T > &matrix_in, NumericVector< T > &solution_in, NumericVector< T > &rhs_in, const double tol, const unsigned int m_its) libmesh_override
virtual numeric_index_type size() const =0
Leaves dofs outside the subset unchanged. This is fastest, but also most confusing because it abandon...
virtual void clear() libmesh_override
virtual LinearConvergenceReason get_converged_reason() const libmesh_override
virtual void close() const libmesh_override
Definition: petsc_matrix.C:881
ImplicitSystem & sys
Set dofs outside the subset to zero.
const class libmesh_nullptr_t libmesh_nullptr
static PetscErrorCode _petsc_shell_matrix_mult_add(Mat mat, Vec arg, Vec add, Vec dest)
PetscErrorCode Vec Mat Mat pc
libmesh_assert(j)
std::unique_ptr< T > UniquePtr
Definition: auto_ptr.h:46
virtual void init_names(const System &) libmesh_override
PetscErrorCode Vec x
void init(triangulateio &t)
PetscErrorCode Vec Mat Mat void * ctx
OStreamProxy err(std::cerr)
virtual std::pair< unsigned int, Real > adjoint_solve(SparseMatrix< T > &matrix_in, NumericVector< T > &solution_in, NumericVector< T > &rhs_in, const double tol, const unsigned int m_its) libmesh_override
virtual void init(const char *name=libmesh_nullptr) libmesh_override
PetscErrorCode __libmesh_petsc_preconditioner_setup(void *ctx)
Set dofs outside the subset to the value of the corresponding dofs of the right hand side...
static void set_petsc_preconditioner_type(const PreconditionerType &preconditioner_type, PC &pc)
virtual numeric_index_type local_size() const =0
std::string enum_to_string(const T e)
virtual void close() libmesh_override
static PetscErrorCode _petsc_shell_matrix_mult(Mat mat, Vec arg, Vec dest)
PetscTruth PetscBool
Definition: petsc_macro.h:64
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
PetscErrorCode ierr
const Parallel::Communicator & comm() const
SparseMatrix interface to PETSc Mat.
bool initialized()
Definition: libmesh.C:272
virtual void vector_mult_add(NumericVector< T > &dest, const NumericVector< T > &arg) const =0
PetscErrorCode __libmesh_petsc_preconditioner_apply(void *ctx, Vec x, Vec y)
static PetscErrorCode _petsc_shell_matrix_get_diagonal(Mat mat, Vec dest)
CHKERRQ(ierr)
virtual void restrict_solve_to(const std::vector< unsigned int > *const dofs, const SubsetSolveMode subset_solve_mode=SUBSET_ZERO) libmesh_override
processor_id_type n_processors()
Definition: libmesh_base.h:88