MFEM  v4.6.0
Finite element discretization library
Public Member Functions | List of all members
mfem::MortarAssembler Class Reference

This class implements the serial 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 algorithm allows to perform quadrature in the intersection of elements of two separate and unrelated 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. More...

#include <mortarassembler.hpp>

Public Member Functions

 MortarAssembler (const std::shared_ptr< FiniteElementSpace > &source, const std::shared_ptr< FiniteElementSpace > &destination)
 constructs the object with source and destination spaces More...
 
 ~MortarAssembler ()
 
bool Assemble (std::shared_ptr< SparseMatrix > &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). More...
 
bool Transfer (const GridFunction &src_fun, GridFunction &dest_fun)
 transfer a function from source to destination. if the transfer is to be performed multiple times use Assemble instead More...
 
bool Apply (const GridFunction &src_fun, GridFunction &dest_fun)
 transfer a function from source to destination. It requires that the Update function is called before More...
 
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). More...
 
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. More...
 
void SetVerbose (const bool verbose)
 Expose process details with verbose output. More...
 
void SetAssembleMassAndCouplingTogether (const bool value)
 Control if the Mass matrix is computed together with the coupling operator every time. More...
 
void SetMaxSolverIterations (const int max_solver_iterations)
 Control the maximum numbers of conjugate gradients steps for mass matrix inversion. More...
 

Detailed Description

This class implements the serial 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 algorithm allows to perform quadrature in the intersection of elements of two separate and unrelated 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.

Definition at line 34 of file mortarassembler.hpp.

Constructor & Destructor Documentation

◆ MortarAssembler()

mfem::MortarAssembler::MortarAssembler ( const std::shared_ptr< FiniteElementSpace > &  source,
const std::shared_ptr< FiniteElementSpace > &  destination 
)

constructs the object with source and destination spaces

Parameters
sourcethe source space from where we want to transfer the discrete field
destinationthe source space to where we want to transfer the discrete field

Definition at line 155 of file mortarassembler.cpp.

◆ ~MortarAssembler()

mfem::MortarAssembler::~MortarAssembler ( )
default

Member Function Documentation

◆ AddMortarIntegrator()

void mfem::MortarAssembler::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.

Parameters
integratorthe integrator object

Definition at line 83 of file mortarassembler.cpp.

◆ Apply()

bool mfem::MortarAssembler::Apply ( const GridFunction src_fun,
GridFunction dest_fun 
)

transfer a function from source to destination. It requires that the Update function is called before

Parameters
src_funthe function associated with the source finite element space
[out]dest_funthe function associated with the destination finite element space
Returns
true if the transfer was successful, fail otherwise.

Definition at line 335 of file mortarassembler.cpp.

◆ Assemble()

bool mfem::MortarAssembler::Assemble ( std::shared_ptr< SparseMatrix > &  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).

Parameters
Bthe assembled coupling operator. B can be passed uninitialized.
Returns
true if there was an intersection and the operator has been assembled. False otherwise.

Definition at line 171 of file mortarassembler.cpp.

◆ SetAssembleMassAndCouplingTogether()

void mfem::MortarAssembler::SetAssembleMassAndCouplingTogether ( const bool  value)

Control if the Mass matrix is computed together with the coupling operator every time.

Parameters
valueis 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 73 of file mortarassembler.cpp.

◆ SetMaxSolverIterations()

void mfem::MortarAssembler::SetMaxSolverIterations ( const int  max_solver_iterations)

Control the maximum numbers of conjugate gradients steps for mass matrix inversion.

Parameters
max_solver_iterationsthe maximum number of iterations

Definition at line 78 of file mortarassembler.cpp.

◆ SetVerbose()

void mfem::MortarAssembler::SetVerbose ( const bool  verbose)

Expose process details with verbose output.

Parameters
verboseis set to true for verbose output

Definition at line 89 of file mortarassembler.cpp.

◆ Transfer()

bool mfem::MortarAssembler::Transfer ( const GridFunction src_fun,
GridFunction dest_fun 
)

transfer a function from source to destination. if the transfer is to be performed multiple times use Assemble instead

Parameters
src_funthe function associated with the source finite element space
[out]dest_funthe function associated with the destination finite element space
Returns
true if there was an intersection and the output can be used.

Definition at line 329 of file mortarassembler.cpp.

◆ Update()

bool mfem::MortarAssembler::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).

Returns
true if there was an intersection and the operator has been assembled. False otherwise.

Definition at line 363 of file mortarassembler.cpp.


The documentation for this class was generated from the following files: