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)
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_);
195 CommunicateSharedVdofs(dst);
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_);
228 CommunicateSharedVdofs(z_);
230 for (
int i = 0; i < sub2_to_parent_map_.
Size(); i++)
237 CorrectFaceOrientations(*dst.
ParFESpace(), z_, dst);
241 MFEM_ABORT(
"unknown TransferCategory: " << category_);
245void 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_);
258void 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)
318 if (parent_face_ori.
Size() == 0) {
return; }
323 bool face = (
dim == 3);
329 for (
int i = 0; i < (face ? mesh->
GetNumFaces() : mesh->GetNE()); i++)
331 if (parent_face_ori[i] == 0) {
continue; }
339 Fo[0] = parent_face_ori[i];
340 doftrans.SetFaceOrientations(Fo);
351 if (sub_to_parent_map)
354 doftrans.TransformPrimal(face_vector);
359 doftrans.InvTransformPrimal(face_vector);
362 for (
int j = 0; j < vdofs.
Size(); j++)
367 if (sub_to_parent_map)
376 dst[k] =
s * face_vector[j];
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
int Size() const
Return the logical size of the array.
@ GaussLegendre
Open type.
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
int GetOrder() const
Return the order (polynomial degree) of the FE collection, corresponding to the order/degree returned...
virtual const StatelessDofTransformation * DofTransformationForGeometry(Geometry::Type GeomType) const
Returns a DoF transformation object compatible with this basis and geometry type.
int GetMapType(int dim) const
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 ...
Ordering::Type GetOrdering() const
Return the ordering method.
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...
const FiniteElementCollection * FEColl() const
int GetVDim() const
Returns vector dimension.
static int DecodeDof(int dof)
Helper to return the DOF associated with a sign encoded DOF.
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.
int GetNumFaces() const
Return the number of faces (3D), edges (2D) or vertices (1D).
Geometry::Type GetFaceGeometry(int i) const
Return the Geometry::Type associated with face i.
Geometry::Type GetElementGeometry(int i) const
int Dimension() const
Dimension of the reference space used within the elements.
Abstract parallel finite element space.
GroupCommunicator & GroupComm()
Return a reference to the internal GroupCommunicator (on VDofs)
ParMesh * GetParMesh() const
Class for parallel grid function.
ParFiniteElementSpace * ParFESpace() const
Class for parallel meshes.
Subdomain representation of a topological parent in another ParMesh.
static bool IsParSubMesh(const ParMesh *m)
Check if ParMesh m is a ParSubMesh.
const Array< int > & GetParentElementIDMap() const
Get the parent element id map.
const Array< int > & GetParentFaceOrientations() const
Get the relative face orientations.
SubMesh::From GetFrom() const
Get the From indicator.
void Transfer(const ParGridFunction &src, ParGridFunction &dst) const
Transfer the source ParGridFunction to the destination ParGridFunction.
ParTransferMap(const ParGridFunction &src, const ParGridFunction &dst)
Construct a new ParTransferMap object which transfers degrees of freedom from the source ParGridFunct...
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
void GetRow(int i, Array< int > &row) const
Return row i in array row (the Table must be finalized)
int Size() const
Returns the number of TYPE I elements.
int Size_of_connections() const
virtual const real_t * HostRead() const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), false).
int Size() const
Returns the size of the vector.
void SetSize(int s)
Resize the vector to size s.
virtual real_t * HostWrite()
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), false).
virtual real_t * HostReadWrite()
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), false).
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Extract entries listed in dofs to the output Vector elemvect.
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.
RT GetRootParent(const T &m)
Identify the root parent of a given SubMesh.
std::function< real_t(const Vector &)> f(real_t mass_coeff)