MFEM v4.8.0
Finite element discretization library
Loading...
Searching...
No Matches
ptransfermap.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
13
14#ifdef MFEM_USE_MPI
15
16#include "psubmesh.hpp"
17#include "ptransfermap.hpp"
18#include "submesh_utils.hpp"
19
20using namespace mfem;
21
23 const ParFiniteElementSpace &dst)
24{
27 {
28 ParSubMesh* src_sm = static_cast<ParSubMesh*>(src.GetParMesh());
29 ParSubMesh* dst_sm = static_cast<ParSubMesh*>(dst.GetParMesh());
30
31 // There is no immediate relation and both src and dst come from a
32 // SubMesh, check if they have an equivalent root parent.
33 if (SubMeshUtils::GetRootParent(*src_sm) !=
35 {
36 MFEM_ABORT("Can't find a relation between the two GridFunctions");
37 }
38
40
41 {
42 ParMesh * parent_mesh =
43 const_cast<ParMesh *>(SubMeshUtils::GetRootParent(*src_sm));
44
45 int parent_dim = parent_mesh->Dimension();
46 int src_sm_dim = src_sm->Dimension();
47 int dst_sm_dim = dst_sm->Dimension();
48
49 bool root_fes_reset = false;
50 if (src_sm_dim == parent_dim - 1 && dst_sm_dim == parent_dim - 1)
51 {
52 const FiniteElementCollection *src_fec = src.FEColl();
53 const FiniteElementCollection *dst_fec = dst.FEColl();
54
55 const L2_FECollection *src_l2_fec =
56 dynamic_cast<const L2_FECollection*>(src_fec);
57 const L2_FECollection *dst_l2_fec =
58 dynamic_cast<const L2_FECollection*>(dst_fec);
59
60 if (src_l2_fec != NULL && dst_l2_fec != NULL)
61 {
62 // Source and destination are both lower dimension L2 spaces.
63 // Transfer them as the trace of an RT space if possible.
64
65 int src_mt = src_fec->GetMapType(src_sm_dim);
66 int dst_mt = dst_fec->GetMapType(dst_sm_dim);
67
68 int src_bt = src_l2_fec->GetBasisType();
69 int dst_bt = dst_l2_fec->GetBasisType();
70
71 int src_p = src_fec->GetOrder();
72 int dst_p = dst_fec->GetOrder();
73
74 if (src_mt == FiniteElement::INTEGRAL &&
75 dst_mt == FiniteElement::INTEGRAL &&
76 src_bt == BasisType::GaussLegendre &&
77 dst_bt == BasisType::GaussLegendre &&
78 src_p == dst_p)
79 {
80 // The subspaces are consistent with the trace of an RT space
81 root_fec_.reset(new RT_FECollection(src_p, parent_dim));
82 root_fes_.reset(new ParFiniteElementSpace(
83 const_cast<ParMesh *>(
85 root_fec_.get()));
86 root_fes_reset = true;
87 }
88 }
89 }
90
91 if (!root_fes_reset)
92 {
93 root_fes_.reset(new ParFiniteElementSpace(
94 src,
95 const_cast<ParMesh *>(
97 }
98 }
99
100 src_to_parent.reset(new ParTransferMap(src, *root_fes_));
101 dst_to_parent.reset(new ParTransferMap(dst, *root_fes_));
102 parent_to_dst.reset(new ParTransferMap(*root_fes_, dst));
103
104 z_.SetSpace(root_fes_.get());
105 }
106 else if (ParSubMesh::IsParSubMesh(src.GetParMesh()))
107 {
109 ParSubMesh* src_sm = static_cast<ParSubMesh*>(src.GetParMesh());
111 dst,
112 src_sm->GetFrom(),
113 src_sm->GetParentElementIDMap(),
114 sub_to_parent_map_);
115
116 root_gc_ = &dst.GroupComm();
117 CommunicateIndicesSet(sub_to_parent_map_, dst.GetVSize());
118 sub_fes_ = &src;
119 }
120 else if (ParSubMesh::IsParSubMesh(dst.GetParMesh()))
121 {
123 ParSubMesh* dst_sm = static_cast<ParSubMesh*>(dst.GetParMesh());
125 src,
126 dst_sm->GetFrom(),
127 dst_sm->GetParentElementIDMap(),
128 sub_to_parent_map_);
129 sub_fes_ = &dst;
130 }
131 else
132 {
133 MFEM_ABORT("Trying to do a transfer between GridFunctions but none of them is defined on a SubMesh");
134 }
135}
136
138 const ParGridFunction &dst)
139 : ParTransferMap(*src.ParFESpace(), *dst.ParFESpace())
140{ }
141
143 ParGridFunction &dst) const
144{
145 if (category_ == TransferCategory::ParentToSubMesh)
146 {
147 // dst = S1^T src
148 src.HostRead();
149 dst.HostWrite(); // dst is fully overwritten
150 for (int i = 0; i < sub_to_parent_map_.Size(); i++)
151 {
152 real_t s = 1.0;
153 int j = FiniteElementSpace::DecodeDof(sub_to_parent_map_[i], s);
154 dst(i) = s * src(j);
155 }
156
157 CorrectFaceOrientations(*sub_fes_, src, dst);
158 }
159 else if (category_ == TransferCategory::SubMeshToParent)
160 {
161 // dst = G S1 src
162 // = G z
163 //
164 // G is identity if the partitioning matches
165
166 src.HostRead();
167 dst.HostReadWrite(); // dst is only partially overwritten
168 for (int i = 0; i < sub_to_parent_map_.Size(); i++)
169 {
170 real_t s = 1.0;
171 int j = FiniteElementSpace::DecodeDof(sub_to_parent_map_[i], s);
172 dst(j) = s * src(i);
173 }
174
175 CorrectFaceOrientations(*sub_fes_, src, dst,
176 &sub_to_parent_map_);
177
178 CommunicateSharedVdofs(dst);
179 }
180 else if (category_ == TransferCategory::SubMeshToSubMesh)
181 {
182 dst_to_parent->Transfer(dst, z_);
183 src_to_parent->Transfer(src, z_);
184 parent_to_dst->Transfer(z_, dst);
185 }
186 else
187 {
188 MFEM_ABORT("unknown TransferCategory: " << category_);
189 }
190}
191
192void ParTransferMap::CommunicateIndicesSet(Array<int> &map, int dst_sz)
193{
194 indices_set_local_.SetSize(dst_sz);
195 indices_set_local_ = 0;
196 for (int i = 0; i < map.Size(); i++)
197 {
198 indices_set_local_[(map[i]>=0)?map[i]:(-map[i]-1)] = 1;
199 }
200 indices_set_global_ = indices_set_local_;
201 root_gc_->Reduce(indices_set_global_, GroupCommunicator::Sum);
202 root_gc_->Bcast(indices_set_global_);
203}
204
205void ParTransferMap::CommunicateSharedVdofs(Vector &f) const
206{
207 // f is usually defined on the root vdofs
208
209 const Table &group_ldof = root_gc_->GroupLDofTable();
210
211 // Identify indices that were only set by other ranks and clear the dof.
212 for (int i = 0; i < group_ldof.Size_of_connections(); i++)
213 {
214 const int j = group_ldof.GetJ()[i];
215 if (indices_set_global_[j] != 0 && indices_set_local_[j] == 0)
216 {
217 f(j) = 0.0;
218 }
219 }
220
221 // TODO: do the reduce only on dofs of interest
222 root_gc_->Reduce<real_t>(f.HostReadWrite(), GroupCommunicator::Sum);
223
224 // Indices that were set from this rank or other ranks have been summed up
225 // and therefore need to be "averaged". Note that this results in the exact
226 // value that is desired.
227 for (int i = 0; i < group_ldof.Size_of_connections(); i++)
228 {
229 const int j = group_ldof.GetJ()[i];
230 if (indices_set_global_[j] != 0)
231 {
232 f(j) /= indices_set_global_[j];
233 }
234 }
235
236 // Indices for dofs that are shared between processors need to be divided by
237 // the whole group size that share this dof.
238 for (int gr = 1; gr < group_ldof.Size(); gr++)
239 {
240 for (int i = 0; i < group_ldof.RowSize(gr); i++)
241 {
242 const int j = group_ldof.GetRow(gr)[i];
243 if (indices_set_global_[j] == 0)
244 {
245 f(j) /= root_gc_->GetGroupTopology().GetGroupSize(gr);
246 }
247 }
248 }
249
250 root_gc_->Bcast<real_t>(f.HostReadWrite());
251}
252
253void
254ParTransferMap::CorrectFaceOrientations(const ParFiniteElementSpace &fes,
255 const Vector &src,
256 Vector &dst,
257 const Array<int> *sub_to_parent_map)
258{
259 const FiniteElementCollection * fec = fes.FEColl();
260
261 ParSubMesh * mesh = dynamic_cast<ParSubMesh*>(fes.GetParMesh());
262
263 const Array<int>& parent_face_ori = mesh->GetParentFaceOrientations();
264
265 if (parent_face_ori.Size() == 0) { return; }
266
267 DofTransformation doftrans(fes.GetVDim(), fes.GetOrdering());
268
269 int dim = mesh->Dimension();
270 bool face = (dim == 3);
271
272 Array<int> vdofs;
273 Array<int> Fo(1);
274 Vector face_vector;
275
276 for (int i = 0; i < (face ? mesh->GetNumFaces() : mesh->GetNE()); i++)
277 {
278 if (parent_face_ori[i] == 0) { continue; }
279
280 Geometry::Type geom = face ? mesh->GetFaceGeometry(i) :
281 mesh->GetElementGeometry(i);
282
283 if (!fec->DofTransformationForGeometry(geom)) { continue; }
284 doftrans.SetDofTransformation(*fec->DofTransformationForGeometry(geom));
285
286 Fo[0] = parent_face_ori[i];
287 doftrans.SetFaceOrientations(Fo);
288
289 if (face)
290 {
291 fes.GetFaceVDofs(i, vdofs);
292 }
293 else
294 {
295 fes.GetElementVDofs(i, vdofs);
296 }
297
298 if (sub_to_parent_map)
299 {
300 src.GetSubVector(vdofs, face_vector);
301 doftrans.TransformPrimal(face_vector);
302 }
303 else
304 {
305 dst.GetSubVector(vdofs, face_vector);
306 doftrans.InvTransformPrimal(face_vector);
307 }
308
309 for (int j = 0; j < vdofs.Size(); j++)
310 {
311 real_t s = 1.0;
312 int k = FiniteElementSpace::DecodeDof(vdofs[j], s);
313
314 if (sub_to_parent_map)
315 {
316 real_t sps = 1.0;
317 int spk = FiniteElementSpace::DecodeDof((*sub_to_parent_map)[k],
318 sps);
319 s *= sps;
320 k = spk;
321 }
322
323 dst[k] = s * face_vector[j];
324 }
325 }
326}
327
328#endif // MFEM_USE_MPI
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition array.hpp:758
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
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
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition fespace.hpp:848
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
Table & GroupLDofTable()
Fill-in the returned Table reference to initialize the GroupCommunicator then call Finalize().
void Reduce(T *ldata, void(*Op)(OpData< T >)) const
Reduce within each group where the master is the root.
const GroupTopology & GetGroupTopology()
Get a reference to the associated GroupTopology object.
static void Sum(OpData< T >)
Reduce operation Sum, instantiated for int and double.
void Bcast(T *ldata, int layout) const
Broadcast within each group where the master is the root.
int GetGroupSize(int g) const
Get the number of processors in a group.
Arbitrary order "L2-conforming" discontinuous finite elements.
Definition fe_coll.hpp:346
int GetBasisType() const
Definition fe_coll.hpp:390
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
Abstract parallel finite element space.
Definition pfespace.hpp:29
GroupCommunicator & GroupComm()
Return a reference to the internal GroupCommunicator (on VDofs)
Definition pfespace.hpp:405
ParMesh * GetParMesh() const
Definition pfespace.hpp:338
Class for parallel grid function.
Definition pgridfunc.hpp:50
void SetSpace(FiniteElementSpace *f) override
Associate a new FiniteElementSpace with the ParGridFunction.
Definition pgridfunc.cpp:97
Class for parallel meshes.
Definition pmesh.hpp:34
Subdomain representation of a topological parent in another ParMesh.
Definition psubmesh.hpp:54
static bool IsParSubMesh(const ParMesh *m)
Check if ParMesh m is a ParSubMesh.
Definition psubmesh.hpp:232
const Array< int > & GetParentElementIDMap() const
Get the parent element id map.
Definition psubmesh.hpp:110
const Array< int > & GetParentFaceOrientations() const
Get the relative face orientations.
Definition psubmesh.hpp:150
SubMesh::From GetFrom() const
Get the From indicator.
Definition psubmesh.hpp:100
ParTransferMap represents a mapping of degrees of freedom from a source ParGridFunction to a destinat...
void Transfer(const ParGridFunction &src, ParGridFunction &dst) const
Transfer the source ParGridFunction to the destination ParGridFunction.
ParTransferMap(const ParFiniteElementSpace &src, const ParFiniteElementSpace &dst)
Construct a new ParTransferMap object which transfers degrees of freedom from the source ParFiniteEle...
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition fe_coll.hpp:403
int * GetJ()
Definition table.hpp:122
int RowSize(int i) const
Definition table.hpp:116
void GetRow(int i, Array< int > &row) const
Return row i in array row (the Table must be finalized)
Definition table.cpp:259
int Size() const
Returns the number of TYPE I elements.
Definition table.hpp:100
int Size_of_connections() const
Definition table.hpp:106
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
std::function< real_t(const Vector &)> f(real_t mass_coeff)
Definition lor_mms.hpp:30