36 MFEM_ABORT(
"Can't find a relation between the two GridFunctions");
45 int parent_dim = parent_mesh->
Dimension();
49 bool root_fes_reset =
false;
50 if (src_sm_dim == parent_dim - 1 && dst_sm_dim == parent_dim - 1)
60 if (src_l2_fec != NULL && dst_l2_fec != NULL)
86 root_fes_reset =
true;
117 CommunicateIndicesSet(sub_to_parent_map_, dst.
GetVSize());
133 MFEM_ABORT(
"Trying to do a transfer between GridFunctions but none of them is defined on a SubMesh");
150 for (
int i = 0; i < sub_to_parent_map_.
Size(); i++)
157 CorrectFaceOrientations(*sub_fes_, src, dst);
168 for (
int i = 0; i < sub_to_parent_map_.
Size(); i++)
175 CorrectFaceOrientations(*sub_fes_, src, dst,
176 &sub_to_parent_map_);
178 CommunicateSharedVdofs(dst);
182 dst_to_parent->Transfer(dst, z_);
183 src_to_parent->Transfer(src, z_);
184 parent_to_dst->Transfer(z_, dst);
188 MFEM_ABORT(
"unknown TransferCategory: " << category_);
192void ParTransferMap::CommunicateIndicesSet(
Array<int> &map,
int dst_sz)
194 indices_set_local_.
SetSize(dst_sz);
195 indices_set_local_ = 0;
196 for (
int i = 0; i < map.Size(); i++)
198 indices_set_local_[(map[i]>=0)?map[i]:(-map[i]-1)] = 1;
200 indices_set_global_ = indices_set_local_;
202 root_gc_->
Bcast(indices_set_global_);
205void ParTransferMap::CommunicateSharedVdofs(
Vector &
f)
const
214 const int j = group_ldof.
GetJ()[i];
215 if (indices_set_global_[j] != 0 && indices_set_local_[j] == 0)
229 const int j = group_ldof.
GetJ()[i];
230 if (indices_set_global_[j] != 0)
232 f(j) /= indices_set_global_[j];
238 for (
int gr = 1; gr < group_ldof.
Size(); gr++)
240 for (
int i = 0; i < group_ldof.
RowSize(gr); i++)
242 const int j = group_ldof.
GetRow(gr)[i];
243 if (indices_set_global_[j] == 0)
265 if (parent_face_ori.
Size() == 0) {
return; }
270 bool face = (
dim == 3);
276 for (
int i = 0; i < (face ? mesh->
GetNumFaces() : mesh->GetNE()); i++)
278 if (parent_face_ori[i] == 0) {
continue; }
286 Fo[0] = parent_face_ori[i];
287 doftrans.SetFaceOrientations(Fo);
298 if (sub_to_parent_map)
301 doftrans.TransformPrimal(face_vector);
306 doftrans.InvTransformPrimal(face_vector);
309 for (
int j = 0; j < vdofs.
Size(); j++)
314 if (sub_to_parent_map)
323 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 GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
int GetVDim() const
Returns the vector dimension of the finite element space.
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.
void SetSpace(FiniteElementSpace *f) override
Associate a new FiniteElementSpace with the ParGridFunction.
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.
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.
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).
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.
auto GetRootParent(const T &m) -> decltype(std::declval< T >().GetParent())
Identify the root parent of a given SubMesh.
std::function< real_t(const Vector &)> f(real_t mass_coeff)