MFEM v4.7.0
Finite element discretization library
|
This class implements the parallel variational transfer between finite element spaces. Variational transfer has been shown to have better approximation properties than standard interpolation. This facilities can be used for supporting applications which require the handling of non matching meshes. For instance: General multi-physics problems, fluid structure interaction, or even visualization of average quantities within subvolumes. This particular code is also used with LLNL for large scale multilevel Monte Carlo simulations. This algorithm allows to perform quadrature in the intersection of elements of two separate, unrelated, and arbitrarily distributed meshes. It generates quadrature rules in the intersection which allows us to integrate with to machine precision using the mfem::MortarIntegrator interface. See https://doi.org/10.1137/15M1008361 for and in-depth explanation. At this time curved elements are not supported. Convex non-affine elements are partially supported, however, high order (>3) finite element discretizations or nonlinear geometric transformations might generate undesired oscillations. Discontinuous fields in general can only be mapped to order 0 destination fields. For such cases localized versions of the projection will have to be developed. More...
#include <pmortarassembler.hpp>
Public Member Functions | |
ParMortarAssembler (const std::shared_ptr< ParFiniteElementSpace > &source, const std::shared_ptr< ParFiniteElementSpace > &destination) | |
constructs the object with source and destination spaces | |
~ParMortarAssembler () | |
bool | Assemble (std::shared_ptr< HypreParMatrix > &B) |
assembles the coupling matrix B. B : source -> destination If u is a coefficient associated with source and v with destination Then v = M^(-1) * B * u; where M is the mass matrix in destination. Works with L2_FECollection, H1_FECollection and DG_FECollection (experimental with RT_FECollection and ND_FECollection). | |
bool | Transfer (const ParGridFunction &src_fun, ParGridFunction &dest_fun) |
transfer a function from source to destination. if the transfer is to be performed multiple times use Assemble or Update/Apply instead | |
bool | Apply (const ParGridFunction &src_fun, ParGridFunction &dest_fun) |
transfer a function from source to destination. It requires that the Update function is called before | |
bool | Update () |
assembles the various components necessary for the transfer. To be called before calling the Apply function if the mesh geometry changed, after previous call. Works with L2_FECollection, H1_FECollection and DG_FECollection (experimental with RT_FECollection and ND_FECollection). | |
void | AddMortarIntegrator (const std::shared_ptr< MortarIntegrator > &integrator) |
This method must be called before Assemble or Transfer. It will assemble the operator in all intersections found. | |
void | SetVerbose (const bool verbose) |
Expose process details with verbose output. | |
void | SetAssembleMassAndCouplingTogether (const bool value) |
Control if the Mass matrix is computed together with the coupling operator every time. | |
void | SetMaxSolverIterations (const int max_solver_iterations) |
Control the maximum numbers of conjugate gradients steps for mass matrix inversion. | |
void | SetSolver (const std::shared_ptr< IterativeSolver > &solver) |
Changes the solver to be used for solving the mass-matrix linear system. | |
This class implements the parallel variational transfer between finite element spaces. Variational transfer has been shown to have better approximation properties than standard interpolation. This facilities can be used for supporting applications which require the handling of non matching meshes. For instance: General multi-physics problems, fluid structure interaction, or even visualization of average quantities within subvolumes. This particular code is also used with LLNL for large scale multilevel Monte Carlo simulations. This algorithm allows to perform quadrature in the intersection of elements of two separate, unrelated, and arbitrarily distributed meshes. It generates quadrature rules in the intersection which allows us to integrate with to machine precision using the mfem::MortarIntegrator interface. See https://doi.org/10.1137/15M1008361 for and in-depth explanation. At this time curved elements are not supported. Convex non-affine elements are partially supported, however, high order (>3) finite element discretizations or nonlinear geometric transformations might generate undesired oscillations. Discontinuous fields in general can only be mapped to order 0 destination fields. For such cases localized versions of the projection will have to be developed.
Definition at line 49 of file pmortarassembler.hpp.
mfem::ParMortarAssembler::ParMortarAssembler | ( | const std::shared_ptr< ParFiniteElementSpace > & | source, |
const std::shared_ptr< ParFiniteElementSpace > & | destination ) |
constructs the object with source and destination spaces
source | the source space from where we want to transfer the discrete field |
destination | the source space to where we want to transfer the discrete field |
Definition at line 1087 of file pmortarassembler.cpp.
|
default |
void mfem::ParMortarAssembler::AddMortarIntegrator | ( | const std::shared_ptr< MortarIntegrator > & | integrator | ) |
This method must be called before Assemble or Transfer. It will assemble the operator in all intersections found.
integrator | the integrator object |
Definition at line 111 of file pmortarassembler.cpp.
bool mfem::ParMortarAssembler::Apply | ( | const ParGridFunction & | src_fun, |
ParGridFunction & | dest_fun ) |
transfer a function from source to destination. It requires that the Update function is called before
src_fun | the function associated with the source finite element space | |
[out] | dest_fun | the function associated with the destination finite element space |
Definition at line 1132 of file pmortarassembler.cpp.
bool mfem::ParMortarAssembler::Assemble | ( | std::shared_ptr< HypreParMatrix > & | B | ) |
assembles the coupling matrix B. B : source -> destination If u is a coefficient associated with source and v with destination Then v = M^(-1) * B * u; where M is the mass matrix in destination. Works with L2_FECollection, H1_FECollection and DG_FECollection (experimental with RT_FECollection and ND_FECollection).
B | the assembled coupling operator. B can be passed uninitialized. |
Definition at line 1097 of file pmortarassembler.cpp.
void mfem::ParMortarAssembler::SetAssembleMassAndCouplingTogether | ( | const bool | value | ) |
Control if the Mass matrix is computed together with the coupling operator every time.
value | is set to true for computing the mass matrix operator with the coupling operator, false otherwise. The option is true by default, set to false if only the coupling operator is needed. |
Definition at line 99 of file pmortarassembler.cpp.
void mfem::ParMortarAssembler::SetMaxSolverIterations | ( | const int | max_solver_iterations | ) |
Control the maximum numbers of conjugate gradients steps for mass matrix inversion.
max_solver_iterations | the maximum number of iterations |
Definition at line 104 of file pmortarassembler.cpp.
void mfem::ParMortarAssembler::SetSolver | ( | const std::shared_ptr< IterativeSolver > & | solver | ) |
Changes the solver to be used for solving the mass-matrix linear system.
solver | new strategy |
Definition at line 94 of file pmortarassembler.cpp.
void mfem::ParMortarAssembler::SetVerbose | ( | const bool | verbose | ) |
Expose process details with verbose output.
verbose | set to true for verbose output |
Definition at line 117 of file pmortarassembler.cpp.
bool mfem::ParMortarAssembler::Transfer | ( | const ParGridFunction & | src_fun, |
ParGridFunction & | dest_fun ) |
transfer a function from source to destination. if the transfer is to be performed multiple times use Assemble or Update/Apply instead
src_fun | the function associated with the source finite element space | |
[out] | dest_fun | the function associated with the destination finite element space |
Definition at line 1125 of file pmortarassembler.cpp.
bool mfem::ParMortarAssembler::Update | ( | ) |
assembles the various components necessary for the transfer. To be called before calling the Apply function if the mesh geometry changed, after previous call. Works with L2_FECollection, H1_FECollection and DG_FECollection (experimental with RT_FECollection and ND_FECollection).
Definition at line 1185 of file pmortarassembler.cpp.