MFEM v4.9.0
Finite element discretization library
Loading...
Searching...
No Matches
mortarassembler.hpp
Go to the documentation of this file.
1// Copyright (c) 2010-2025, Lawrence Livermore National Security, LLC. Produced
2// at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3// LICENSE and NOTICE for details. LLNL-CODE-806117.
4//
5// This file is part of the MFEM library. For more information and source code
6// availability visit https://mfem.org.
7//
8// MFEM is free software; you can redistribute it and/or modify it under the
9// terms of the BSD-3 license. We welcome feedback and contributions, see file
10// CONTRIBUTING.md for details.
11
12#ifndef MFEML2P_MORTAR_ASSEMBLER_HPP
13#define MFEML2P_MORTAR_ASSEMBLER_HPP
14
16
17#ifdef MFEM_USE_MOONOLITH
18
19#include "../fem.hpp"
20#include "mortarintegrator.hpp"
21#include <memory>
22
23namespace mfem
24{
25/*!
26 * @brief This class implements the serial variational transfer between finite
27 * element spaces. Variational transfer has been shown to have better
28 * approximation properties than standard interpolation. This facilities can be
29 * used for supporting applications which require the handling of non matching
30 * meshes. For instance: General multi-physics problems, fluid structure
31 * interaction, or even visualization of average quantities within subvolumes
32 * This algorithm allows to perform quadrature in the intersection of elements
33 * of two separate and unrelated meshes. It generates quadrature rules in
34 * the intersection which allows us to integrate-with to machine precision using
35 * the mfem::MortarIntegrator interface. See https://doi.org/10.1137/15M1008361
36 * for and in-depth explanation. At this time curved elements are not supported.
37 */
39{
40public:
41 /*!
42 * @brief constructs the object with source and destination spaces
43 * @param source the source space from where we want to transfer the discrete
44 * field
45 * @param destination the source space to where we want to transfer the
46 * discrete field
47 */
48 MortarAssembler(const std::shared_ptr<FiniteElementSpace> &source,
49 const std::shared_ptr<FiniteElementSpace> &destination);
50
52
53 /*!
54 * @brief assembles the coupling matrix B. B : source -> destination If u is a
55 * coefficient associated with source and v with destination Then v = M^(-1) *
56 * B * u; where M is the mass matrix in destination. Works with
57 * L2_FECollection, H1_FECollection and DG_FECollection (experimental with
58 * RT_FECollection and ND_FECollection).
59 * @param B the assembled coupling operator. B can be passed uninitialized.
60 * @return true if there was an intersection and the operator has been
61 * assembled. False otherwise.
62 */
63 bool Assemble(std::shared_ptr<SparseMatrix> &B);
64
65 /*!
66 * @brief transfer a function from source to destination. if the transfer is
67 * to be performed multiple times use Assemble instead
68 * @param src_fun the function associated with the source finite element space
69 * @param[out] dest_fun the function associated with the destination finite
70 * element space
71 * @return true if there was an intersection and the output can be used.
72 */
73 bool Transfer(const GridFunction &src_fun, GridFunction &dest_fun);
74
75 /*!
76 * @brief transfer a function from source to destination. It requires that
77 * the Update function is called before
78 * @param src_fun the function associated with the source finite element space
79 * @param[out] dest_fun the function associated with the destination finite
80 * element space
81 * @return true if the transfer was successful, fail otherwise.
82 */
83 bool Apply(const GridFunction &src_fun, GridFunction &dest_fun);
84
85 /*!
86 * @brief assembles the various components necessary for the transfer.
87 * To be called before calling the Apply function if the mesh geometry
88 * changed, after previous call. Works with L2_FECollection, H1_FECollection
89 * and DG_FECollection (experimental with RT_FECollection and
90 * ND_FECollection).
91 * @return true if there was an intersection and the operator has been
92 * assembled. False otherwise.
93 */
94 bool Update();
95
96 /*!
97 * @brief This method must be called before Assemble or Transfer.
98 * It will assemble the operator in all intersections found.
99 * @param integrator the integrator object
100 */
101 void AddMortarIntegrator(const std::shared_ptr<MortarIntegrator> &integrator);
102
103 /*!
104 * @brief Expose process details with verbose output
105 * @param verbose is set to true for verbose output
106 */
107 void SetVerbose(const bool verbose);
108
109
110 /*!
111 * @brief Control if the Mass matrix is computed together with the coupling
112 * operator every time.
113 * @param value is set to true for computing the mass matrix operator with
114 * the coupling operator, false otherwise. The option is true by default,
115 * set to false if only the coupling operator is needed.
116 */
117 void SetAssembleMassAndCouplingTogether(const bool value);
118
119 /*!
120 * @brief Control the maximum numbers of conjugate gradients steps for mass
121 * matrix inversion
122 * @param max_solver_iterations the maximum number of iterations
123 */
124 void SetMaxSolverIterations(const int max_solver_iterations);
125
126private:
127 struct Impl;
128 std::unique_ptr<Impl> impl_;
129};
130
131} // namespace mfem
132
133#endif // MFEM_USE_MOONOLITH
134#endif // MFEML2P_MORTAR_ASSEMBLER_HPP
Class for grid function - Vector with associated FE space.
Definition gridfunc.hpp:32
This class implements the serial variational transfer between finite element spaces....
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 SetVerbose(const bool verbose)
Expose process details with verbose output.
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...
MortarAssembler(const std::shared_ptr< FiniteElementSpace > &source, const std::shared_ptr< FiniteElementSpace > &destination)
constructs the object with source and destination spaces
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
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...
bool Assemble(std::shared_ptr< SparseMatrix > &B)
assembles the coupling matrix B. B : source -> destination If u is a coefficient associated with sour...
bool Update()
assembles the various components necessary for the transfer. To be called before calling the Apply fu...
void source(const Vector &x, Vector &f)
Definition ex25.cpp:620