12 #ifndef MFEM_L2P_PAR_MORTAR_ASSEMBLER_HPP 13 #define MFEM_L2P_PAR_MORTAR_ASSEMBLER_HPP 15 #include "../../config/config.hpp" 60 const std::shared_ptr<ParFiniteElementSpace> &destination);
74 bool Assemble(std::shared_ptr<HypreParMatrix> &B);
142 void SetSolver(
const std::shared_ptr<IterativeSolver> &solver);
148 std::unique_ptr<Impl> impl_;
153 #endif // MFEM_USE_MPI 154 #endif // MFEM_L2P_PAR_MORTAR_ASSEMBLER_HPP void SetVerbose(const bool verbose)
Expose process details with verbose output.
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...
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.
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...
void AddMortarIntegrator(const std::shared_ptr< MortarIntegrator > &integrator)
This method must be called before Assemble or Transfer. It will assemble the operator in all intersec...
void source(const Vector &x, Vector &f)
ParMortarAssembler(const std::shared_ptr< ParFiniteElementSpace > &source, const std::shared_ptr< ParFiniteElementSpace > &destination)
constructs the object with source and destination spaces
void SetMaxSolverIterations(const int max_solver_iterations)
Control the maximum numbers of conjugate gradients steps for mass matrix inversion.
bool Assemble(std::shared_ptr< HypreParMatrix > &B)
assembles the coupling matrix B. B : source -> destination If u is a coefficient associated with sour...
void SetAssembleMassAndCouplingTogether(const bool value)
Control if the Mass matrix is computed together with the coupling operator every time.
void SetSolver(const std::shared_ptr< IterativeSolver > &solver)
Changes the solver to be used for solving the mass-matrix linear system.
Class for parallel grid function.
bool Update()
assembles the various components necessary for the transfer. To be called before calling the Apply fu...