12 #include "../../config/config.hpp" 39 MFEM_ABORT(
"Can't find a relation between the two GridFunctions");
48 int parent_dim = parent_mesh->
Dimension();
52 bool root_fes_reset =
false;
53 if (src_sm_dim == parent_dim - 1 && dst_sm_dim == parent_dim - 1)
66 if (src_l2_fec != NULL && dst_l2_fec != NULL)
89 const_cast<ParMesh *>(
92 root_fes_reset =
true;
113 sub1_to_parent_map_);
119 sub2_to_parent_map_);
121 root_gc_ = &root_fes_->GroupComm();
122 CommunicateIndicesSet(sub1_to_parent_map_, root_fes_->GetVSize());
124 z_.
SetSize(root_fes_->GetVSize());
136 sub1_to_parent_map_);
139 CommunicateIndicesSet(sub1_to_parent_map_, dst.
Size());
151 sub1_to_parent_map_);
155 MFEM_ABORT(
"Trying to do a transfer between GridFunctions but none of them is defined on a SubMesh");
167 for (
int i = 0; i < sub1_to_parent_map_.
Size(); i++)
174 CorrectFaceOrientations(*dst.
ParFESpace(), src, dst);
185 for (
int i = 0; i < sub1_to_parent_map_.
Size(); i++)
192 CorrectFaceOrientations(*src.
ParFESpace(), src, dst,
193 &sub1_to_parent_map_);
208 for (
int i = 0; i < sub2_to_parent_map_.
Size(); i++)
215 CorrectFaceOrientations(*dst.
ParFESpace(), dst, z_,
216 &sub2_to_parent_map_);
218 for (
int i = 0; i < sub1_to_parent_map_.
Size(); i++)
225 CorrectFaceOrientations(*src.
ParFESpace(), src, z_,
226 &sub1_to_parent_map_);
230 for (
int i = 0; i < sub2_to_parent_map_.
Size(); i++)
237 CorrectFaceOrientations(*dst.
ParFESpace(), z_, dst);
241 MFEM_ABORT(
"unknown TransferCategory: " << category_);
245 void ParTransferMap::CommunicateIndicesSet(
Array<int> &map,
int dst_sz)
247 indices_set_local_.
SetSize(dst_sz);
248 indices_set_local_ = 0;
249 for (
int i = 0; i < map.
Size(); i++)
251 indices_set_local_[(map[i]>=0)?map[i]:(-map[i]-1)] = 1;
253 indices_set_global_ = indices_set_local_;
255 root_gc_->
Bcast(indices_set_global_);
258 void ParTransferMap::CommunicateSharedVdofs(
Vector &
f)
const 267 const int j = group_ldof.
GetJ()[i];
268 if (indices_set_global_[j] != 0 && indices_set_local_[j] == 0)
282 const int j = group_ldof.
GetJ()[i];
283 if (indices_set_global_[j] != 0)
285 f(j) /= indices_set_global_[j];
291 for (
int gr = 1; gr < group_ldof.
Size(); gr++)
293 for (
int i = 0; i < group_ldof.
RowSize(gr); i++)
295 const int j = group_ldof.
GetRow(gr)[i];
296 if (indices_set_global_[j] == 0)
303 root_gc_->
Bcast<
double>(
f.HostReadWrite());
318 if (parent_face_ori.
Size() == 0) {
return; }
324 bool face = (
dim == 3);
332 if (parent_face_ori[i] == 0) {
continue; }
340 if (doftrans == NULL) {
continue; }
342 vdoftrans.SetDofTransformation(*doftrans);
344 Fo[0] = parent_face_ori[i];
345 vdoftrans.SetFaceOrientations(Fo);
356 if (sub_to_parent_map)
359 vdoftrans.TransformPrimal(face_vector);
364 vdoftrans.InvTransformPrimal(face_vector);
367 for (
int j = 0; j < vdofs.
Size(); j++)
372 if (sub_to_parent_map)
381 dst[k] =
s * face_vector[j];
386 #endif // MFEM_USE_MPI
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.
const Array< int > & GetParentElementIDMap() const
Get the parent element id map.
int Dimension() const
Dimension of the reference space used within the elements.
void Bcast(T *ldata, int layout) const
Broadcast within each group where the master is the root.
int GetNumFaces() const
Return the number of faces (3D), edges (2D) or vertices (1D).
void SetSize(int s)
Resize the vector to size s.
virtual const double * HostRead() const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), false).
static int DecodeDof(int dof)
Helper to return the DOF associated with a sign encoded DOF.
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...
int GetOrder() const
Return the order (polynomial degree) of the FE collection, corresponding to the order/degree returned...
Geometry::Type GetFaceGeometry(int i) const
Return the Geometry::Type associated with face i.
int Size() const
Returns the size of the vector.
virtual double * HostWrite()
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), false).
Abstract parallel finite element space.
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Extract entries listed in dofs to the output Vector elemvect.
Subdomain representation of a topological parent in another ParMesh.
std::function< double(const Vector &)> f(double mass_coeff)
static void Sum(OpData< T >)
Reduce operation Sum, instantiated for int and double.
const FiniteElementCollection * FEColl() const
ParMesh * GetParMesh() const
static bool IsParSubMesh(const ParMesh *m)
Check if ParMesh m is a ParSubMesh.
void GetRow(int i, Array< int > &row) const
Return row i in array row (the Table must be finalized)
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
int GetVDim() const
Returns vector dimension.
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 ...
Geometry::Type GetElementGeometry(int i) const
void Transfer(const ParGridFunction &src, ParGridFunction &dst) const
Transfer the source ParGridFunction to the destination ParGridFunction.
const GroupTopology & GetGroupTopology()
Get a reference to the associated GroupTopology object.
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
ParFiniteElementSpace * ParFESpace() const
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
GroupCommunicator & GroupComm()
Return a reference to the internal GroupCommunicator (on VDofs)
RT GetRootParent(const T &m)
Identify the root parent of a given SubMesh.
Table & GroupLDofTable()
Fill-in the returned Table reference to initialize the GroupCommunicator then call Finalize()...
int GetNE() const
Returns number of elements.
int GetGroupSize(int g) const
Get the number of processors in a group.
virtual int GetMapType(int dim) const
int Size() const
Returns the number of TYPE I elements.
Ordering::Type GetOrdering() const
Return the ordering method.
void Reduce(T *ldata, void(*Op)(OpData< T >)) const
Reduce within each group where the master is the root.
int Size_of_connections() const
int Size() const
Return the logical size of the array.
ParTransferMap(const ParGridFunction &src, const ParGridFunction &dst)
Construct a new ParTransferMap object which transfers degrees of freedom from the source ParGridFunct...
virtual StatelessDofTransformation * DofTransformationForGeometry(Geometry::Type GeomType) const
Returns a DoF transformation object compatible with this basis and geometry type. ...
const Array< int > & GetParentFaceOrientations() const
Get the relative face orientations.
Class for parallel grid function.
Class for parallel meshes.
SubMesh::From GetFrom() const
Get the From indicator.
virtual double * HostReadWrite()
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), false).
Arbitrary order "L2-conforming" discontinuous finite elements.