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