MFEM v4.8.0
Finite element discretization library
Loading...
Searching...
No Matches
transfermap.cpp
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#include "submesh.hpp"
13#include "transfermap.hpp"
14#include "submesh_utils.hpp"
15
16using namespace mfem;
17
19 const FiniteElementSpace &dst)
20{
21 if (SubMesh::IsSubMesh(src.GetMesh()) &&
23 {
24 SubMesh* src_sm = static_cast<SubMesh*>(src.GetMesh());
25 SubMesh* dst_sm = static_cast<SubMesh*>(dst.GetMesh());
26
27 // There is no immediate relation and both src and dst come from a
28 // SubMesh, check if they have an equivalent root parent.
29 if (SubMeshUtils::GetRootParent(*src_sm) !=
31 {
32 MFEM_ABORT("Can't find a relation between the two GridFunctions");
33 }
34
36
37 {
38 Mesh * parent_mesh =
39 const_cast<Mesh *>(SubMeshUtils::GetRootParent(*src_sm));
40
41 int parent_dim = parent_mesh->Dimension();
42 int src_sm_dim = src_sm->Dimension();
43 int dst_sm_dim = dst_sm->Dimension();
44
45 bool root_fes_reset = false;
46 if (src_sm_dim == parent_dim - 1 && dst_sm_dim == parent_dim - 1)
47 {
48 const FiniteElementCollection *src_fec = src.FEColl();
49 const FiniteElementCollection *dst_fec = dst.FEColl();
50
51 const L2_FECollection *src_l2_fec =
52 dynamic_cast<const L2_FECollection*>(src_fec);
53 const L2_FECollection *dst_l2_fec =
54 dynamic_cast<const L2_FECollection*>(dst_fec);
55
56 if (src_l2_fec != NULL && dst_l2_fec != NULL)
57 {
58 // Source and destination are both lower dimension L2 spaces.
59 // Transfer them as the trace of an RT space if possible.
60
61 int src_mt = src_fec->GetMapType(src_sm_dim);
62 int dst_mt = dst_fec->GetMapType(dst_sm_dim);
63
64 int src_bt = src_l2_fec->GetBasisType();
65 int dst_bt = dst_l2_fec->GetBasisType();
66
67 int src_p = src_fec->GetOrder();
68 int dst_p = dst_fec->GetOrder();
69
70 if (src_mt == FiniteElement::INTEGRAL &&
71 dst_mt == FiniteElement::INTEGRAL &&
72 src_bt == BasisType::GaussLegendre &&
73 dst_bt == BasisType::GaussLegendre &&
74 src_p == dst_p)
75 {
76 // The subspaces are consistent with the trace of an RT space
77 root_fec_.reset(new RT_FECollection(src_p, parent_dim));
78 root_fes_.reset(new FiniteElementSpace(
79 const_cast<Mesh *>(
81 root_fec_.get()));
82 root_fes_reset = true;
83 }
84 }
85 }
86
87 if (!root_fes_reset)
88 {
89 const FiniteElementCollection *src_fec = src.FEColl();
90 const FiniteElementCollection *dst_fec = dst.FEColl();
91 auto *root_fec = src_sm_dim == parent_dim ? src_fec : dst_fec;
92
93 root_fes_.reset(new FiniteElementSpace(
94 const_cast<Mesh *>(
95 SubMeshUtils::GetRootParent(*src_sm)), root_fec));
96 }
97 }
98
99 src_to_parent.reset(new TransferMap(src, *root_fes_));
100 dst_to_parent.reset(new TransferMap(dst, *root_fes_));
101 parent_to_dst.reset(new TransferMap(*root_fes_, dst));
102
103 z_.SetSpace(root_fes_.get());
104 }
105 else if (SubMesh::IsSubMesh(src.GetMesh()))
106 {
108 SubMesh* src_sm = static_cast<SubMesh*>(src.GetMesh());
110 dst,
111 src_sm->GetFrom(),
112 src_sm->GetParentElementIDMap(),
113 sub_to_parent_map_);
114 sub_fes_ = &src;
115 }
116 else if (SubMesh::IsSubMesh(dst.GetMesh()))
117 {
119 SubMesh* dst_sm = static_cast<SubMesh*>(dst.GetMesh());
121 src,
122 dst_sm->GetFrom(),
123 dst_sm->GetParentElementIDMap(),
124 sub_to_parent_map_);
125 sub_fes_ = &dst;
126 }
127 else
128 {
129 MFEM_ABORT("Trying to do a transfer between GridFunctions but none of them is defined on a SubMesh");
130 }
131}
132
134 const GridFunction &dst)
135 : TransferMap(*src.FESpace(), *dst.FESpace())
136{ }
137
139{
140 if (category_ == TransferCategory::ParentToSubMesh)
141 {
142 // dst = S1^T src
143 src.HostRead();
144 dst.HostWrite(); // dst is fully overwritten
145 for (int i = 0; i < sub_to_parent_map_.Size(); i++)
146 {
147 real_t s = 1.0;
148 int j = FiniteElementSpace::DecodeDof(sub_to_parent_map_[i], s);
149 dst(i) = s * src(j);
150 }
151
152 CorrectFaceOrientations(*sub_fes_, src, dst);
153 }
154 else if (category_ == TransferCategory::SubMeshToParent)
155 {
156 // dst = G S1 src
157 // = G z
158 //
159 // G is identity if the partitioning matches
160
161 src.HostRead();
162 dst.HostReadWrite(); // dst is only partially overwritten
163 for (int i = 0; i < sub_to_parent_map_.Size(); i++)
164 {
165 real_t s = 1.0;
166 int j = FiniteElementSpace::DecodeDof(sub_to_parent_map_[i], s);
167 dst(j) = s * src(i);
168 }
169
170 CorrectFaceOrientations(*sub_fes_, src, dst,
171 &sub_to_parent_map_);
172 }
173 else if (category_ == TransferCategory::SubMeshToSubMesh)
174 {
175 dst_to_parent->Transfer(dst, z_);
176 src_to_parent->Transfer(src, z_);
177 parent_to_dst->Transfer(z_, dst);
178 }
179 else
180 {
181 MFEM_ABORT("unknown TransferCategory: " << category_);
182 }
183}
184
185void TransferMap::CorrectFaceOrientations(const FiniteElementSpace &fes,
186 const Vector &src,
187 Vector &dst,
188 const Array<int> *sub_to_parent_map)
189{
190 const FiniteElementCollection * fec = fes.FEColl();
191
192 SubMesh * mesh = dynamic_cast<SubMesh*>(fes.GetMesh());
193
194 const Array<int>& parent_face_ori = mesh->GetParentFaceOrientations();
195
196 if (parent_face_ori.Size() == 0) { return; }
197
198 DofTransformation doftrans(fes.GetVDim(), fes.GetOrdering());
199
200 int dim = mesh->Dimension();
201 bool face = (dim == 3);
202
203 Array<int> vdofs;
204 Array<int> Fo(1);
205 Vector face_vector;
206
207 for (int i = 0; i < (face ? mesh->GetNumFaces() : mesh->GetNE()); i++)
208 {
209 if (parent_face_ori[i] == 0) { continue; }
210
211 Geometry::Type geom = face ? mesh->GetFaceGeometry(i) :
212 mesh->GetElementGeometry(i);
213
214 if (!fec->DofTransformationForGeometry(geom)) { continue; }
215 doftrans.SetDofTransformation(*fec->DofTransformationForGeometry(geom));
216
217 Fo[0] = parent_face_ori[i];
218 doftrans.SetFaceOrientations(Fo);
219
220 if (face)
221 {
222 fes.GetFaceVDofs(i, vdofs);
223 }
224 else
225 {
226 fes.GetElementVDofs(i, vdofs);
227 }
228
229 if (sub_to_parent_map)
230 {
231 src.GetSubVector(vdofs, face_vector);
232 doftrans.TransformPrimal(face_vector);
233 }
234 else
235 {
236 dst.GetSubVector(vdofs, face_vector);
237 doftrans.InvTransformPrimal(face_vector);
238 }
239
240 for (int j = 0; j < vdofs.Size(); j++)
241 {
242 real_t s = 1.0;
243 int k = FiniteElementSpace::DecodeDof(vdofs[j], s);
244
245 if (sub_to_parent_map)
246 {
247 real_t sps = 1.0;
248 int spk = FiniteElementSpace::DecodeDof((*sub_to_parent_map)[k],
249 sps);
250 s *= sps;
251 k = spk;
252 }
253
254 dst[k] = s * face_vector[j];
255 }
256 }
257}
int Size() const
Return the logical size of the array.
Definition array.hpp:147
@ GaussLegendre
Open type.
Definition fe_base.hpp:35
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition fe_coll.hpp:27
int GetOrder() const
Return the order (polynomial degree) of the FE collection, corresponding to the order/degree returned...
Definition fe_coll.hpp:240
virtual const StatelessDofTransformation * DofTransformationForGeometry(Geometry::Type GeomType) const
Returns a DoF transformation object compatible with this basis and geometry type.
Definition fe_coll.hpp:67
int GetMapType(int dim) const
Definition fe_coll.cpp:60
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:244
DofTransformation * GetElementVDofs(int i, Array< int > &vdofs) const
Returns indices of degrees of freedom for the i'th element. The returned indices are offsets into an ...
Definition fespace.cpp:332
Ordering::Type GetOrdering() const
Return the ordering method.
Definition fespace.hpp:876
void GetFaceVDofs(int i, Array< int > &vdofs) const
Returns the indices of the degrees of freedom for the specified face, including the DOFs for the edge...
Definition fespace.cpp:361
const FiniteElementCollection * FEColl() const
Definition fespace.hpp:878
Mesh * GetMesh() const
Returns the mesh.
Definition fespace.hpp:679
int GetVDim() const
Returns the vector dimension of the finite element space.
Definition fespace.hpp:841
static int DecodeDof(int dof)
Helper to return the DOF associated with a sign encoded DOF.
Definition fespace.hpp:1168
Class for grid function - Vector with associated FE space.
Definition gridfunc.hpp:31
virtual void SetSpace(FiniteElementSpace *f)
Associate a new FiniteElementSpace with the GridFunction.
Definition gridfunc.cpp:225
Arbitrary order "L2-conforming" discontinuous finite elements.
Definition fe_coll.hpp:346
int GetBasisType() const
Definition fe_coll.hpp:390
Mesh data type.
Definition mesh.hpp:64
int GetNumFaces() const
Return the number of faces (3D), edges (2D) or vertices (1D).
Definition mesh.cpp:6531
Geometry::Type GetFaceGeometry(int i) const
Return the Geometry::Type associated with face i.
Definition mesh.cpp:1476
Geometry::Type GetElementGeometry(int i) const
Definition mesh.hpp:1434
int Dimension() const
Dimension of the reference space used within the elements.
Definition mesh.hpp:1216
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition fe_coll.hpp:403
Subdomain representation of a topological parent in another Mesh.
Definition submesh.hpp:44
const Array< int > & GetParentElementIDMap() const
Get the parent element id map.
Definition submesh.hpp:106
From GetFrom() const
Get the From indicator.
Definition submesh.hpp:96
const Array< int > & GetParentFaceOrientations() const
Get the relative face orientations.
Definition submesh.hpp:136
static bool IsSubMesh(const Mesh *m)
Check if Mesh m is a SubMesh.
Definition submesh.hpp:221
TransferMap represents a mapping of degrees of freedom from a source GridFunction to a destination Gr...
void Transfer(const GridFunction &src, GridFunction &dst) const
Transfer the source GridFunction to the destination GridFunction.
TransferMap(const FiniteElementSpace &src, const FiniteElementSpace &dst)
Construct a new TransferMap object which transfers degrees of freedom from the source FiniteElementSp...
Vector data type.
Definition vector.hpp:82
virtual const real_t * HostRead() const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), false).
Definition vector.hpp:498
virtual real_t * HostWrite()
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), false).
Definition vector.hpp:506
virtual real_t * HostReadWrite()
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), false).
Definition vector.hpp:514
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Extract entries listed in dofs to the output Vector elemvect.
Definition vector.cpp:653
int dim
Definition ex24.cpp:53
void BuildVdofToVdofMap(const FiniteElementSpace &subfes, const FiniteElementSpace &parentfes, const SubMesh::From &from, const Array< int > &parent_element_ids, Array< int > &vdof_to_vdof_map)
Build the vdof to vdof mapping between two FiniteElementSpace objects.
auto GetRootParent(const T &m) -> decltype(std::declval< T >().GetParent())
Identify the root parent of a given SubMesh.
float real_t
Definition config.hpp:43