https://mooseframework.inl.gov
NavierStokesProblem.h
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 
10 #pragma once
11 
12 #include "FEProblem.h"
13 #include "libmesh/libmesh_config.h"
14 #include <petscsnes.h>
15 
16 class NonlinearSystem;
17 
23 {
24 public:
26 
28 
29 #if PETSC_RELEASE_GREATER_EQUALS(3, 20, 0)
30 
34 
39 
43  void clearIndexSets() { _index_sets.clear(); }
44 
45  /*
46  * Given a \p KSP \p node and where we are in the field split tree, given by \p tree_position,
47  * return the next \p KSP object in the tree on the way to the Schur complement \p KSP object.
48  * Each invocation of this method moves through one level of our \p _index_sets data member. This
49  * method will call itself recursively until it reaches the Schur complement \p KSP
50  */
51  KSP findSchurKSP(KSP node, unsigned int tree_position);
52 
57  void setupLSCMatrices(KSP schur_ksp);
58 
62  virtual ~NavierStokesProblem();
63 
64  virtual void initialSetup() override;
65 
66 protected:
70  virtual void initPetscOutputAndSomeSolverSettings() override;
71 
72 private:
75  const bool _commute_lsc;
77  const TagName & _mass_matrix;
79  const TagName & _L_matrix;
81  const bool _have_mass_matrix;
83  const bool _have_L_matrix;
84 
89 
93  const std::vector<unsigned int> & _schur_fs_index;
94 
96  Mat _Q_scale = nullptr;
98  Mat _L = nullptr;
99 
111  std::vector<IS> _index_sets;
112 #endif
113 };
const bool _have_mass_matrix
Whether the user attached a mass matrix.
const std::vector< unsigned int > & _schur_fs_index
The length of this vector should correspond to the number of split nesting levels there are in the fi...
virtual void initialSetup() override
unsigned int TagID
const InputParameters & parameters() const
void clearIndexSets()
Clear the field split index sets.
std::vector< IS > _index_sets
This will end up being the same length as _schur_fs_index.
Mat _L
The Poisson operator.
virtual ~NavierStokesProblem()
Will destroy any matrices we allocated.
KSP findSchurKSP(KSP node, unsigned int tree_position)
virtual TagID getMatrixTagID(const TagName &tag_name) const
NavierStokesProblem(const InputParameters &parameters)
const bool _commute_lsc
Whether to commute operators in the style of Olshanskii.
A problem that handles Schur complement preconditioning of the incompressible Navier-Stokes equations...
const bool _pressure_mass_matrix_as_pre
Whether to directly use the pressure mass matrix to form the Schur complement preconditioner.
virtual void initPetscOutputAndSomeSolverSettings() override
Reinitialize PETSc output for proper linear/nonlinear iteration display.
Mat _Q_scale
The mass matrix used for scaling.
const TagName & _L_matrix
The tag name of the Poisson operator.
TagID LMatrixTagID() const
const TagName & _mass_matrix
The tag name of the mass matrix.
void setupLSCMatrices(KSP schur_ksp)
Setup the Least Squares Commutator (LSC) preconditioner given the Schur complement KSP object...
static InputParameters validParams()
TagID massMatrixTagID() const
const bool _have_L_matrix
Whether the user attached a Poisson operator matrix.