https://mooseframework.inl.gov
MFEMCutTransitionSubMesh.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
10 #ifdef MOOSE_MFEM_ENABLED
11 
13 #include "MFEMProblem.h"
14 
16 
19 {
22  params.addClassDescription(
23  "Class to construct an MFEMSubMesh formed from the set of elements that have least one "
24  "vertex on the specified cut plane, that lie on one side of the plane, "
25  "and that are restricted to the set of user-specified subdomains.");
26  params.addRequiredParam<BoundaryName>("cut_boundary",
27  "The boundary associated with the mesh cut.");
28  params.addRequiredParam<BoundaryName>(
29  "transition_subdomain_boundary",
30  "Name to assign boundary of transition subdomain not shared with cut surface.");
31  params.addRequiredParam<SubdomainName>(
32  "transition_subdomain",
33  "The name of the subdomain to be created on the mesh comprised of the set of elements "
34  "adjacent to the cut surface on one side.");
35  params.addRequiredParam<SubdomainName>(
36  "closed_subdomain",
37  "The name of the subdomain attribute to be created comprised of the set of all elements "
38  "of the closed geometry, including the new transition region.");
39  return params;
40 }
41 
43  : MFEMSubMesh(parameters),
44  MFEMBlockRestrictable(parameters, getMFEMProblem().mesh().getMFEMParMesh()),
45  _cut_boundary(getParam<BoundaryName>("cut_boundary")),
46  _cut_submesh(std::make_shared<mfem::ParSubMesh>(mfem::ParSubMesh::CreateFromBoundary(
47  getMFEMProblem().mesh().getMFEMParMesh(),
48  getMesh().bdr_attribute_sets.GetAttributeSet(_cut_boundary)))),
49  _transition_subdomain_boundary(getParam<BoundaryName>("transition_subdomain_boundary")),
50  _transition_subdomain(getParam<SubdomainName>("transition_subdomain")),
51  _closed_subdomain(getParam<SubdomainName>("closed_subdomain")),
52  _cut_normal(3)
53 {
54 }
55 
56 void
58 {
59  labelMesh(getMFEMProblem().mesh().getMFEMParMesh());
60  _submesh = std::make_shared<mfem::ParSubMesh>(mfem::ParSubMesh::CreateFromDomain(
61  getMesh(), getMesh().attribute_sets.GetAttributeSet(_transition_subdomain)));
62  _submesh->attribute_sets.attr_sets = getMesh().attribute_sets.attr_sets;
63  _submesh->bdr_attribute_sets.attr_sets = getMesh().bdr_attribute_sets.attr_sets;
64 }
65 
66 void
67 MFEMCutTransitionSubMesh::labelMesh(mfem::ParMesh & parent_mesh)
68 {
69  int mpi_comm_rank = getMFEMProblem().getProblemData().myid;
70  int mpi_comm_size = getMFEMProblem().getProblemData().num_procs;
71 
72  // First determine face normal based on the first boundary element found on the cut
73  // to use when determining orientation relative to the cut
74  const mfem::Array<int> & parent_cut_element_id_map = _cut_submesh->GetParentElementIDMap();
75  int rank_with_submesh = -1;
76  if (parent_cut_element_id_map.Size() > 0)
77  {
78  int reference_face = parent_cut_element_id_map[0];
79  _cut_normal = findFaceNormal(parent_mesh, parent_mesh.GetBdrElementFaceIndex(reference_face));
80  rank_with_submesh = mpi_comm_rank;
81  }
82  MPI_Allreduce(MPI_IN_PLACE, &rank_with_submesh, 1, MPI_INT, MPI_MAX, getMFEMProblem().getComm());
83  MPI_Bcast(_cut_normal.GetData(),
84  _cut_normal.Size(),
85  MFEM_MPI_REAL_T,
86  rank_with_submesh,
88  // Iterate over all vertices on cut, find elements with those vertices,
89  // and declare them transition elements if they are on the +ve side of the cut
90  mfem::Array<int> transition_els;
91  std::vector<HYPRE_BigInt> global_cut_vert_ids;
92  mfem::Array<HYPRE_BigInt> gi;
93  parent_mesh.GetGlobalVertexIndices(gi);
94  std::unique_ptr<mfem::Table> vert_to_elem(parent_mesh.GetVertexToElementTable());
95  const mfem::Array<int> & cut_to_parent_vertex_id_map = _cut_submesh->GetParentVertexIDMap();
96  for (int i = 0; i < _cut_submesh->GetNV(); ++i)
97  {
98  int cut_vert = cut_to_parent_vertex_id_map[i];
99  global_cut_vert_ids.push_back(gi[cut_vert]);
100  int ne = vert_to_elem->RowSize(cut_vert); // number of elements touching cut vertex
101  const int * els_adj_to_cut = vert_to_elem->GetRow(cut_vert); // elements touching cut vertex
102  for (int i = 0; i < ne; i++)
103  {
104  const int el_adj_to_cut = els_adj_to_cut[i];
105  if (isInDomain(el_adj_to_cut, getSubdomainAttributes(), parent_mesh) &&
106  isPositiveSideOfCut(el_adj_to_cut, cut_vert, parent_mesh))
107  transition_els.Append(el_adj_to_cut);
108  }
109  }
110 
111  // share cut verts coords or global ids across all procs
112  int n_cut_vertices = global_cut_vert_ids.size();
113  std::vector<int> cut_vert_sizes(mpi_comm_size, 0);
114  MPI_Allgather(
115  &n_cut_vertices, 1, MPI_INT, &cut_vert_sizes[0], 1, MPI_INT, getMFEMProblem().getComm());
116  // Make an offset array and total the sizes.
117  std::vector<int> n_vert_offset(mpi_comm_size, 0);
118  for (int i = 1; i < mpi_comm_size; i++)
119  n_vert_offset[i] = n_vert_offset[i - 1] + cut_vert_sizes[i - 1];
120  int global_n_cut_vertices = 0;
121  for (int i = 0; i < mpi_comm_size; i++)
122  global_n_cut_vertices += cut_vert_sizes[i];
123 
124  // Gather the queries to all ranks.
125  std::vector<HYPRE_BigInt> all_cut_verts(global_n_cut_vertices, 0);
126  MPI_Allgatherv(&global_cut_vert_ids[0],
127  n_cut_vertices,
128  HYPRE_MPI_BIG_INT,
129  &all_cut_verts[0],
130  &cut_vert_sizes[0],
131  &n_vert_offset[0],
132  HYPRE_MPI_BIG_INT,
133  getMFEMProblem().getComm());
134 
135  // Detect shared vertices and add corresponding elements
136  for (int g = 1, sv = 0; g < parent_mesh.GetNGroups(); g++)
137  {
138  for (int gv = 0; gv < parent_mesh.GroupNVertices(g); gv++, sv++)
139  {
140  // all elements touching this shared vertex should be updated
141  int cut_vert = parent_mesh.GroupVertex(g, gv);
142  for (std::size_t i = 0; i < all_cut_verts.size(); i += 1)
143  {
144  if (gi[cut_vert] == all_cut_verts[i]) // check if shared vertex is on the cut plane
145  {
146  int ne = vert_to_elem->RowSize(cut_vert); // number of elements touching cut vertex
147  const int * els_adj_to_cut =
148  vert_to_elem->GetRow(cut_vert); // elements touching cut vertex
149  for (int i = 0; i < ne; i++)
150  {
151  const int el_adj_to_cut = els_adj_to_cut[i];
152  if (isInDomain(el_adj_to_cut, getSubdomainAttributes(), parent_mesh) &&
153  isPositiveSideOfCut(el_adj_to_cut, cut_vert, parent_mesh))
154  transition_els.Append(el_adj_to_cut);
155  }
156  }
157  }
158  }
159  }
160 
161  transition_els.Sort();
162  transition_els.Unique();
163 
164  setAttributes(parent_mesh, transition_els);
165 }
166 
167 mfem::Vector
168 MFEMCutTransitionSubMesh::findFaceNormal(const mfem::ParMesh & mesh, const int & face)
169 {
170  mooseAssert(mesh.SpaceDimension() == 3,
171  "MFEMCutTransitionSubMesh only works in 3-dimensional meshes!");
172  mfem::Vector normal;
173  mfem::Array<int> face_verts;
174  std::vector<mfem::Vector> v;
175  mesh.GetFaceVertices(face, face_verts);
176 
177  // First we get the coordinates of 3 vertices on the face
178  for (auto vtx : face_verts)
179  {
180  mfem::Vector vtx_coords(3);
181  for (int j = 0; j < 3; ++j)
182  vtx_coords[j] = mesh.GetVertex(vtx)[j];
183  v.push_back(vtx_coords);
184  }
185 
186  // Now we find the unit vector normal to the face
187  v[0] -= v[1];
188  v[1] -= v[2];
189  v[0].cross3D(v[1], normal);
190  normal /= normal.Norml2();
191  return normal;
192 }
193 
194 bool
196  const int & el_vertex_on_cut,
197  mfem::ParMesh & parent_mesh)
198 {
199  const int sdim = parent_mesh.SpaceDimension();
200  mfem::Vector el_center(3);
201  parent_mesh.GetElementCenter(el, el_center);
202  mfem::Vector vertex_coords(parent_mesh.GetVertex(el_vertex_on_cut), sdim);
203  mfem::Vector relative_center(sdim);
204  for (int j = 0; j < sdim; j++)
205  relative_center[j] = el_center[j] - vertex_coords[j];
206  return _cut_normal * relative_center > 0;
207 }
208 
209 void
210 MFEMCutTransitionSubMesh::setAttributes(mfem::ParMesh & parent_mesh,
211  mfem::Array<int> & transition_els)
212 {
213  // Generate a set of local new attributes for transition region elements
214  const int old_max_attr = parent_mesh.attributes.Max();
215  mfem::Array<int> new_attrs(old_max_attr);
216  new_attrs = -1;
217  for (auto const & transition_el : transition_els)
218  {
219  const int old_attr = parent_mesh.GetAttribute(transition_el);
220  new_attrs[old_attr - 1] = old_max_attr + old_attr;
221  // Set the new attribute IDs for transition region elements
222  parent_mesh.SetAttribute(transition_el, new_attrs[old_attr - 1]);
223  }
224 
225  // Distribute local attribute IDs to other MPI ranks to create set of new
226  // global attribute IDs for the transition region.
227  MPI_Allreduce(
228  MPI_IN_PLACE, new_attrs, old_max_attr, MPI_INT, MPI_MAX, getMFEMProblem().getComm());
229 
230  mfem::AttributeSets & attr_sets = parent_mesh.attribute_sets;
231  mfem::AttributeSets & bdr_attr_sets = parent_mesh.bdr_attribute_sets;
232  // Create attribute set labelling the newly created transition region on one side of the cut
233  attr_sets.CreateAttributeSet(_transition_subdomain);
234  // Create attribute set labelling the entire closed geometry
235  attr_sets.SetAttributeSet(_closed_subdomain, getSubdomainAttributes());
236  // Add the new domain attributes to new attribute sets
237  const std::set<std::string> attr_set_names = attr_sets.GetAttributeSetNames();
238  for (int old_attr = 1; old_attr <= old_max_attr; ++old_attr)
239  {
240  int new_attr = new_attrs[old_attr - 1];
241  // Add attributes only if they're transition region attributes
242  if (new_attr != -1)
243  {
244  attr_sets.AddToAttributeSet(_transition_subdomain, new_attr);
245  for (auto const & attr_set_name : attr_set_names)
246  {
247  const mfem::Array<int> & attr_set = attr_sets.GetAttributeSet(attr_set_name);
248  // If attribute set has the old attribute of the transition region, add the new one
249  if (attr_set.Find(old_attr) != -1)
250  attr_sets.AddToAttributeSet(attr_set_name, new_attr);
251  }
252  }
253  }
254 
255  // Create boundary attribute set labelling the exterior of the newly created
256  // transition region, excluding the cut
257  bdr_attr_sets.SetAttributeSet(_transition_subdomain_boundary,
258  mfem::Array<int>({parent_mesh.bdr_attributes.Max() + 1}));
259 
260  parent_mesh.SetAttributes();
261 }
262 
263 bool
265  const mfem::Array<int> & subdomains,
266  const mfem::ParMesh & mesh)
267 {
268  // element<0 for ghost elements
269  if (element < 0)
270  return true;
271 
272  for (const auto & subdomain : subdomains)
273  if (mesh.GetAttribute(element) == subdomain)
274  return true;
275  return false;
276 }
277 
278 #endif
Modifies the MFEM Mesh to label a subdomain consisting of elements adjacent to an interior surface on...
MFEMProblemData & getProblemData()
Method to get the current MFEMProblemData object storing the current data specifying the FE problem...
Definition: MFEMProblem.h:186
const SubdomainName & _closed_subdomain
MeshBase & mesh
virtual void buildSubMesh() override
Constructs the submesh.
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
const BoundaryName & _transition_subdomain_boundary
const mfem::Array< int > & getSubdomainAttributes()
void addRequiredParam(const std::string &name, const std::string &doc_string)
This method adds a parameter and documentation string to the InputParameters object that will be extr...
void labelMesh(mfem::ParMesh &parent_mesh)
Add attributes to the parent mesh representing the cut transition region.
registerMooseObject("MooseApp", MFEMCutTransitionSubMesh)
static InputParameters validParams()
Definition: MFEMSubMesh.C:16
MFEMCutTransitionSubMesh(const InputParameters &parameters)
const SubdomainName & _transition_subdomain
Base class for construction of a mfem::ParSubMesh object.
Definition: MFEMSubMesh.h:22
bool isPositiveSideOfCut(const int &el, const int &el_vertex_on_cut, mfem::ParMesh &mesh)
Checks whether an element lies on the positive or negative side of the cut plane. ...
bool isInDomain(const int &el, const mfem::Array< int > &subdomains, const mfem::ParMesh &mesh)
Checks whether a given element is within a certain domain or vector of domains.
MPI_Comm getComm()
Return the MPI communicator associated with this FE problem&#39;s mesh.
Definition: MFEMProblem.h:192
std::shared_ptr< mfem::ParSubMesh > _cut_submesh
void setAttributes(mfem::ParMesh &parent_mesh, mfem::Array< int > &transition_els)
Set new attributes for the provided transition region elements.
MFEMProblem & getMFEMProblem()
Returns a reference to the MFEMProblem instance.
static InputParameters validParams()
void addClassDescription(const std::string &doc_string)
This method adds a description of the class that will be displayed in the input file syntax dump...
const mfem::ParMesh & getMesh() const
mfem::Vector findFaceNormal(const mfem::ParMesh &mesh, const int &face)
Finds the normal vector of a face in the mesh from its vertices.
std::shared_ptr< mfem::ParSubMesh > _submesh
Stores the constructed submesh.
Definition: MFEMSubMesh.h:39
static InputParameters validParams()
Base class for construction of an object that is restricted to a subset of subdomains of the problem ...