Class FEProjector coordinates data structures for projecting fields in FE models.
More...
#include <ibtk/FEProjector.h>
|
| | FEProjector (libMesh::EquationSystems *equation_systems, const SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > &input_db) |
| | Constructor. More...
|
| |
| | FEProjector (std::shared_ptr< FEData > fe_data, const SAMRAI::tbox::Pointer< SAMRAI::tbox::Database > &input_db) |
| | Alternative constructor that takes in a shared pointer to an FEData object. More...
|
| |
| | FEProjector ()=delete |
| | Deleted default constructor. More...
|
| |
| | FEProjector (const FEProjector &from)=delete |
| | Deleted copy constructor. More...
|
| |
| FEProjector & | operator= (const FEProjector &that)=delete |
| | Deleted assignment operator. More...
|
| |
| | ~FEProjector ()=default |
| | Defaulted destructor. More...
|
| |
| std::pair< libMesh::PetscLinearSolver< double > *, libMesh::PetscMatrix< double > * > | buildL2ProjectionSolver (const std::string &system_name) |
| |
| std::pair< libMesh::PetscLinearSolver< double > *, libMesh::PetscMatrix< double > * > | buildLumpedL2ProjectionSolver (const std::string &system_name) |
| |
| std::pair< libMesh::PetscLinearSolver< double > *, libMesh::PetscMatrix< double > * > | buildStabilizedL2ProjectionSolver (const std::string &system_name, double epsilon) |
| |
| libMesh::PetscVector< double > * | buildDiagonalL2MassMatrix (const std::string &system_name) |
| |
| bool | computeL2Projection (libMesh::PetscVector< double > &U, libMesh::PetscVector< double > &F, const std::string &system_name, bool consistent_mass_matrix=true, bool close_U=true, bool close_F=true, double tol=1.0e-6, unsigned int max_its=100) |
| | Set U to be the L2 projection of F. More...
|
| |
| bool | computeStabilizedL2Projection (libMesh::PetscVector< double > &U, libMesh::PetscVector< double > &F, const std::string &system_name, double epsilon, bool close_U=true, bool close_F=true, double tol=1.0e-6, unsigned int max_its=100) |
| | Set U to be the L2 projection of F with a local projection stabilization term. More...
|
| |
| void | setLoggingEnabled (bool enable_logging=true) |
| | Enable or disable logging. More...
|
| |
| bool | getLoggingEnabled () const |
| | Determine whether logging is enabled or disabled. More...
|
| |
|
| std::shared_ptr< FEData > | d_fe_data |
| |
| std::map< std::string, std::unique_ptr< libMesh::PetscMatrix< double > > > | d_L2_proj_matrix |
| | Data structures for consistent mass matrices and related solvers. More...
|
| |
| std::map< std::string, std::unique_ptr< libMesh::PetscLinearSolver< double > > > | d_L2_proj_solver |
| |
| std::map< std::string, std::unique_ptr< libMesh::PetscMatrix< double > > > | d_lumped_L2_proj_matrix |
| |
| std::map< std::string, std::unique_ptr< libMesh::PetscLinearSolver< double > > > | d_lumped_L2_proj_solver |
| |
| std::map< std::string, std::unique_ptr< libMesh::PetscVector< double > > > | d_diag_L2_proj_matrix |
| |
| std::map< std::string, std::map< double, std::unique_ptr< libMesh::PetscMatrix< double > > > > | d_stab_L2_proj_matrix |
| | Data structures for consistent mass matrices and related solvers with local projection stabilization. More...
|
| |
| std::map< std::string, std::map< double, std::unique_ptr< libMesh::PetscLinearSolver< double > > > > | d_stab_L2_proj_solver |
| |
| std::map< std::string, FischerGuess > | d_initial_guesses |
| |
Parameters read from the input database
-
num_fischer_vectors: Number of previous solution and RHS pairs stored for use in computing the initial guess to each linear system. Defaults to five. Using more vectors will require additional computational work when computing the initial guess (roughly N*N dot products, where N is the number of stored vector pairs) but will decrease the number of solver iterations. The default value lowers the number of solver iterations to, typically, no more than two or three, so its usually the right value.
◆ FEProjector() [1/4]
◆ FEProjector() [2/4]
◆ FEProjector() [3/4]
| IBTK::FEProjector::FEProjector |
( |
| ) |
|
|
delete |
◆ FEProjector() [4/4]
| IBTK::FEProjector::FEProjector |
( |
const FEProjector & |
from | ) |
|
|
delete |
◆ ~FEProjector()
| IBTK::FEProjector::~FEProjector |
( |
| ) |
|
|
default |
◆ operator=()
◆ buildL2ProjectionSolver()
| std::pair<libMesh::PetscLinearSolver<double>*, libMesh::PetscMatrix<double>*> IBTK::FEProjector::buildL2ProjectionSolver |
( |
const std::string & |
system_name | ) |
|
- Returns
- Pointers to a linear solver and sparse matrix corresponding to a L2 projection operator.
◆ buildLumpedL2ProjectionSolver()
| std::pair<libMesh::PetscLinearSolver<double>*, libMesh::PetscMatrix<double>*> IBTK::FEProjector::buildLumpedL2ProjectionSolver |
( |
const std::string & |
system_name | ) |
|
- Returns
- Pointers to a linear solver and sparse matrix corresponding to a L2 projection operator generated with mass lumping. This system will be diagonal if the finite element field does not have any hanging nodes: i.e., the difference between this function and buildDiagonalL2MassMatrix is that this function assembles the mass matrix while imposing constraints.
This mass matrix is constructed using nodal quadrature.
◆ buildStabilizedL2ProjectionSolver()
| std::pair<libMesh::PetscLinearSolver<double>*, libMesh::PetscMatrix<double>*> IBTK::FEProjector::buildStabilizedL2ProjectionSolver |
( |
const std::string & |
system_name, |
|
|
double |
epsilon |
|
) |
| |
- Returns
- Pointers to a linear solver and sparse matrix corresponding to a L2 projection operator with a local projection stabilization term.
This local stabilization approach augments the standard projection with a stabilization of the form (p, q) + epsilon (p - Pi p, q - Pi q) = (f, q). Currently the only supported local projection is to piecewise constant functions.
Ref: https://epubs.siam.org/doi/abs/10.1137/S0036142905444482
◆ buildDiagonalL2MassMatrix()
- Returns
- Pointer to vector representation of diagonal L2 mass matrix. Unlike buildLumpedL2ProjectionSolver this matrix is always diagonal and is always assembled without applying constraints.
This mass matrix is constructed using nodal quadrature.
◆ computeL2Projection()
◆ computeStabilizedL2Projection()
This local stabilization approach augments the standard projection with a stabilization of the form (p, q) + epsilon (p - Pi p, q - Pi q) = (f, q). Currently the only supported local projection is to piecewise constant functions.
Ref: https://epubs.siam.org/doi/abs/10.1137/S0036142905444482
◆ setLoggingEnabled()
| void IBTK::FEProjector::setLoggingEnabled |
( |
bool |
enable_logging = true | ) |
|
◆ getLoggingEnabled()
| bool IBTK::FEProjector::getLoggingEnabled |
( |
| ) |
const |
◆ d_fe_data
| std::shared_ptr<FEData> IBTK::FEProjector::d_fe_data |
|
protected |
FEData object that contains the libMesh data structures.
- Note
- multiple FEDataManager objects may use the same FEData object, usually combined with different hierarchies.
◆ d_L2_proj_matrix
| std::map<std::string, std::unique_ptr<libMesh::PetscMatrix<double> > > IBTK::FEProjector::d_L2_proj_matrix |
|
protected |
◆ d_L2_proj_solver
| std::map<std::string, std::unique_ptr<libMesh::PetscLinearSolver<double> > > IBTK::FEProjector::d_L2_proj_solver |
|
protected |
◆ d_lumped_L2_proj_matrix
| std::map<std::string, std::unique_ptr<libMesh::PetscMatrix<double> > > IBTK::FEProjector::d_lumped_L2_proj_matrix |
|
protected |
Data structures for lumped mass matrices. These are computed in the same way as the normal mass matrix, except the quadrature rule used to compute matrix entries has it's points at the finite element nodes: i.e., in the absence of constraints, the matrix is diagonal.
Here we refer to the unconstrained matrix (which is always diagonal) as proj_matrix_diag and the constrained (should there be constraints) matrix as lumped_L2_proj_matrix.
◆ d_lumped_L2_proj_solver
| std::map<std::string, std::unique_ptr<libMesh::PetscLinearSolver<double> > > IBTK::FEProjector::d_lumped_L2_proj_solver |
|
protected |
◆ d_diag_L2_proj_matrix
◆ d_stab_L2_proj_matrix
| std::map<std::string, std::map<double, std::unique_ptr<libMesh::PetscMatrix<double> > > > IBTK::FEProjector::d_stab_L2_proj_matrix |
|
protected |
◆ d_stab_L2_proj_solver
| std::map<std::string, std::map<double, std::unique_ptr<libMesh::PetscLinearSolver<double> > > > IBTK::FEProjector::d_stab_L2_proj_solver |
|
protected |
◆ d_initial_guesses
| std::map<std::string, FischerGuess> IBTK::FEProjector::d_initial_guesses |
|
protected |
◆ d_linear_solve_system_timers
Pointers for system-specific solver timers.
SAMRAI stores timers in a single unsorted array. To keep timer lookups quick we cache the pointers here.
◆ d_enable_logging
| bool IBTK::FEProjector::d_enable_logging = false |
|
private |
◆ d_num_fischer_vectors
| int IBTK::FEProjector::d_num_fischer_vectors = 5 |
|
private |
The documentation for this class was generated from the following file: