MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
pmortarassembler.hpp
Go to the documentation of this file.
1// Copyright (c) 2010-2024, 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 MFEM_L2P_PAR_MORTAR_ASSEMBLER_HPP
13#define MFEM_L2P_PAR_MORTAR_ASSEMBLER_HPP
14
16
17#ifdef MFEM_USE_MPI
18
19#include <memory>
20#include <vector>
21
22#include "../fem.hpp"
23#include "mortarintegrator.hpp"
24
25namespace mfem
26{
27
28/*!
29 * @brief This class implements the parallel variational transfer between finite
30 * element spaces. Variational transfer has been shown to have better
31 * approximation properties than standard interpolation. This facilities can be
32 * used for supporting applications which require the handling of non matching
33 * meshes. For instance: General multi-physics problems, fluid structure
34 * interaction, or even visualization of average quantities within subvolumes. This
35 * particular code is also used with LLNL for large scale multilevel Monte Carlo
36 * simulations.
37 * This algorithm allows to perform quadrature in the intersection of elements
38 * of two separate, unrelated, and arbitrarily distributed meshes.
39 * It generates quadrature rules in the intersection which allows us to
40 * integrate with to machine precision using the mfem::MortarIntegrator
41 * interface. See https://doi.org/10.1137/15M1008361 for and in-depth
42 * explanation. At this time curved elements are not supported.
43 * Convex non-affine elements are partially supported, however, high order (>3)
44 * finite element discretizations or nonlinear geometric transformations might
45 * generate undesired oscillations. Discontinuous fields in general can only
46 * be mapped to order 0 destination fields.
47 * For such cases localized versions of the projection will have to be developed.
48 */
50{
51public:
52 /*!
53 * @brief constructs the object with source and destination spaces
54 * @param source the source space from where we want to transfer the discrete
55 * field
56 * @param destination the source space to where we want to transfer the
57 * discrete field
58 */
59 ParMortarAssembler(const std::shared_ptr<ParFiniteElementSpace> &source,
60 const std::shared_ptr<ParFiniteElementSpace> &destination);
61
63
64 /*!
65 * @brief assembles the coupling matrix B. B : source -> destination If u is a
66 * coefficient associated with source and v with destination Then v = M^(-1) *
67 * B * u; where M is the mass matrix in destination. Works with
68 * L2_FECollection, H1_FECollection and DG_FECollection (experimental with
69 * RT_FECollection and ND_FECollection).
70 * @param B the assembled coupling operator. B can be passed uninitialized.
71 * @return true if there was an intersection and the operator has been
72 * assembled. False otherwise.
73 */
74 bool Assemble(std::shared_ptr<HypreParMatrix> &B);
75
76 /*!
77 * @brief transfer a function from source to destination. if the transfer is
78 * to be performed multiple times use Assemble or Update/Apply instead
79 * @param src_fun the function associated with the source finite element space
80 * @param[out] dest_fun the function associated with the destination finite
81 * element space
82 * @return true if there was an intersection and the output can be used.
83 */
84 bool Transfer(const ParGridFunction &src_fun, ParGridFunction &dest_fun);
85
86 /*!
87 * @brief transfer a function from source to destination. It requires that
88 * the Update function is called before
89 * @param src_fun the function associated with the source finite element space
90 * @param[out] dest_fun the function associated with the destination finite element space
91 * @return true if the transfer was successful, fail otherwise.
92 */
93 bool Apply(const ParGridFunction &src_fun, ParGridFunction &dest_fun);
94
95 /*!
96 * @brief assembles the various components necessary for the transfer.
97 * To be called before calling the Apply function if the mesh geometry
98 * changed, after previous call. Works with L2_FECollection, H1_FECollection
99 * and DG_FECollection (experimental with RT_FECollection and
100 * ND_FECollection).
101 * @return true if there was an intersection and the operator has been
102 * assembled. False otherwise.
103 */
104 bool Update();
105
106 /*!
107 * @brief This method must be called before Assemble or Transfer.
108 * It will assemble the operator in all intersections found.
109 * @param integrator the integrator object
110 */
111 void AddMortarIntegrator(const std::shared_ptr<MortarIntegrator> &integrator);
112
113 /*!
114 * @brief Expose process details with verbose output
115 * @param verbose set to true for verbose output
116 */
117 void SetVerbose(const bool verbose);
118
119
120 /*!
121 * @brief Control if the Mass matrix is computed together with the coupling
122 * operator every time.
123 * @param value is set to true for computing the mass matrix operator with
124 * the coupling operator, false otherwise. The option is true by default,
125 * set to false if only the coupling operator is needed.
126 */
127 void SetAssembleMassAndCouplingTogether(const bool value);
128
129
130 /*!
131 * @brief Control the maximum numbers of conjugate gradients steps for mass
132 * matrix inversion
133 * @param max_solver_iterations the maximum number of iterations
134 */
135 void SetMaxSolverIterations(const int max_solver_iterations);
136
137 /*!
138 * @brief Changes the solver to be used for solving the mass-matrix linear
139 * system
140 * @param solver new strategy
141 */
142 void SetSolver(const std::shared_ptr<IterativeSolver> &solver);
143
144
145 struct Impl;
146
147 private:
148 std::unique_ptr<Impl> impl_;
149};
150
151} // namespace mfem
152
153#endif // MFEM_USE_MPI
154#endif // MFEM_L2P_PAR_MORTAR_ASSEMBLER_HPP
Class for parallel grid function.
Definition pgridfunc.hpp:33
This class implements the parallel variational transfer between finite element spaces....
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
void SetSolver(const std::shared_ptr< IterativeSolver > &solver)
Changes the solver to be used for solving the mass-matrix linear system.
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...
ParMortarAssembler(const std::shared_ptr< ParFiniteElementSpace > &source, const std::shared_ptr< ParFiniteElementSpace > &destination)
constructs the object with source and destination spaces
void SetAssembleMassAndCouplingTogether(const bool value)
Control if the Mass matrix is computed together with the coupling operator every time.
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 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...
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