https://mooseframework.inl.gov
StaticCondensationFieldSplitPreconditioner.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
11 
12 // MOOSE includes
13 #include "PetscSupport.h"
14 #include "NonlinearSystemBase.h"
15 #include "PetscDMMoose.h"
16 #include "Split.h"
17 
18 #include "libmesh/static_condensation.h"
19 #include "libmesh/static_condensation_dof_map.h"
20 #include "libmesh/petsc_linear_solver.h"
21 
22 using namespace libMesh;
23 
25 
28 {
30 }
31 
33  const InputParameters & params)
35 {
36  std::shared_ptr<Split> top_split = _nl.getSplit(_decomposition_split);
37  top_split->setup(_nl, prefix());
38 }
39 
40 const libMesh::DofMapBase &
42 {
43  return scDofMap();
44 }
45 
46 const libMesh::System &
48 {
49  return scDofMap().reduced_system();
50 }
51 
52 std::string
54 {
56 }
57 
58 void
60 {
61  PetscBool ismoose;
62  DM dm = LIBMESH_PETSC_NULLPTR;
63 
64  // Initialize the part of the DM package that's packaged with Moose; in the PETSc source tree this
65  // call would be in DMInitializePackage()
66  LibmeshPetscCall(DMMooseRegisterAll());
67  // Create and set up the DM that will consume the split options and deal with block matrices.
68  auto & petsc_solver = cast_ref<PetscLinearSolver<Number> &>(scSysMat().reduced_system_solver());
69  auto ksp = petsc_solver.ksp();
70  // if there exists a DMMoose object, not to recreate a new one
71  LibmeshPetscCall(KSPGetDM(ksp, &dm));
72  if (dm)
73  {
74  LibmeshPetscCall(PetscObjectTypeCompare((PetscObject)dm, DMMOOSE, &ismoose));
75  if (ismoose)
76  return;
77  }
78  createMooseDM(&dm);
79  LibmeshPetscCall(KSPSetDM(ksp, dm));
80  // We set the operators ourselves. We do not want the DM to generate the operators
81  LibmeshPetscCall(KSPSetDMActive(ksp, PETSC_FALSE));
82  LibmeshPetscCall(DMDestroy(&dm));
83 }
84 
85 KSP
87 {
88  auto & petsc_solver = cast_ref<PetscLinearSolver<Number> &>(scSysMat().reduced_system_solver());
89  return petsc_solver.ksp();
90 }
registerMooseObjectAliased("MooseApp", StaticCondensationFieldSplitPreconditioner, "SCFSP")
NonlinearSystemBase & _nl
The nonlinear system this FSP is associated with (convenience reference)
static InputParameters validParams()
Constructor.
StaticCondensationFieldSplitPreconditioner(const InputParameters &parameters)
PetscErrorCode PetscOptionItems *PetscErrorCode DM dm
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
The following methods are specializations for using the libMesh::Parallel::packed_range_* routines fo...
Implements a field split preconditioner for a statically condensed system.
Implements a preconditioner designed to map onto PETSc&#39;s PCFieldSplit.
std::shared_ptr< Split > getSplit(const std::string &name)
Retrieves a split by name.
virtual const libMesh::System & system() const override
void createMooseDM(DM *dm)
creates the MOOSE data management object
virtual const libMesh::DofMapBase & dofMap() const override
const System & reduced_system() const
LinearSolver< Number > & reduced_system_solver()
PetscErrorCode DMMooseRegisterAll()
virtual void setupDM() override
setup the data management data structure that manages the field split