Nektar++
Public Member Functions | Static Public Member Functions | Static Public Attributes | List of all members
Nektar::Utilities::ProcessOptiExtract Class Reference

#include <ProcessOptiExtract.h>

Inheritance diagram for Nektar::Utilities::ProcessOptiExtract:
[legend]

Public Member Functions

 ProcessOptiExtract (NekMeshUtils::MeshSharedPtr m)
 
virtual ~ProcessOptiExtract ()
 
virtual void Process ()
 Write mesh to output file. More...
 
- Public Member Functions inherited from Nektar::NekMeshUtils::ProcessModule
NEKMESHUTILS_EXPORT ProcessModule (MeshSharedPtr p_m)
 
- Public Member Functions inherited from Nektar::NekMeshUtils::Module
NEKMESHUTILS_EXPORT Module (MeshSharedPtr p_m)
 
NEKMESHUTILS_EXPORT void RegisterConfig (std::string key, std::string value=std::string())
 Register a configuration option with a module. More...
 
NEKMESHUTILS_EXPORT void PrintConfig ()
 Print out all configuration options for a module. More...
 
NEKMESHUTILS_EXPORT void SetDefaults ()
 Sets default configuration options for those which have not been set. More...
 
NEKMESHUTILS_EXPORT MeshSharedPtr GetMesh ()
 
virtual NEKMESHUTILS_EXPORT void ProcessVertices ()
 Extract element vertices. More...
 
virtual NEKMESHUTILS_EXPORT void ProcessEdges (bool ReprocessEdges=true)
 Extract element edges. More...
 
virtual NEKMESHUTILS_EXPORT void ProcessFaces (bool ReprocessFaces=true)
 Extract element faces. More...
 
virtual NEKMESHUTILS_EXPORT void ProcessElements ()
 Generate element IDs. More...
 
virtual NEKMESHUTILS_EXPORT void ProcessComposites ()
 Generate composites. More...
 
virtual NEKMESHUTILS_EXPORT void ClearElementLinks ()
 

Static Public Member Functions

static std::shared_ptr< Modulecreate (NekMeshUtils::MeshSharedPtr m)
 Creates an instance of this class. More...
 

Static Public Attributes

static NekMeshUtils::ModuleKey className
 

Additional Inherited Members

- Protected Member Functions inherited from Nektar::NekMeshUtils::Module
NEKMESHUTILS_EXPORT void ReorderPrisms (PerMap &perFaces)
 Reorder node IDs so that prisms and tetrahedra are aligned correctly. More...
 
NEKMESHUTILS_EXPORT void PrismLines (int prism, PerMap &perFaces, std::set< int > &prismsDone, std::vector< ElementSharedPtr > &line)
 
- Protected Attributes inherited from Nektar::NekMeshUtils::Module
MeshSharedPtr m_mesh
 Mesh object. More...
 
std::map< std::string, ConfigOptionm_config
 List of configuration values. More...
 

Detailed Description

Definition at line 45 of file ProcessOptiExtract.h.

Constructor & Destructor Documentation

◆ ProcessOptiExtract()

Nektar::Utilities::ProcessOptiExtract::ProcessOptiExtract ( NekMeshUtils::MeshSharedPtr  m)

Definition at line 54 of file ProcessOptiExtract.cpp.

References Nektar::NekMeshUtils::Module::m_config.

54  : ProcessModule(m)
55 {
56  m_config["insert"] =
57  ConfigOption(false, "-1", "Name of mesh file to be combined.");
58 }
Represents a command-line configuration option.
NEKMESHUTILS_EXPORT ProcessModule(MeshSharedPtr p_m)
std::map< std::string, ConfigOption > m_config
List of configuration values.

◆ ~ProcessOptiExtract()

Nektar::Utilities::ProcessOptiExtract::~ProcessOptiExtract ( )
virtual

Definition at line 60 of file ProcessOptiExtract.cpp.

61 {
62 }

Member Function Documentation

◆ create()

static std::shared_ptr<Module> Nektar::Utilities::ProcessOptiExtract::create ( NekMeshUtils::MeshSharedPtr  m)
inlinestatic

Creates an instance of this class.

Definition at line 49 of file ProcessOptiExtract.h.

References Nektar::MemoryManager< DataType >::AllocateSharedPtr().

50  {
52  }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.

◆ Process()

void Nektar::Utilities::ProcessOptiExtract::Process ( )
virtual

Write mesh to output file.

Implements Nektar::NekMeshUtils::Module.

Definition at line 64 of file ProcessOptiExtract.cpp.

References Nektar::NekMeshUtils::Module::ClearElementLinks(), Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::FieldUtils::eInputModule, Nektar::LibUtilities::eTriangle, CellMLToNektar.pycml::extract(), Nektar::NekMeshUtils::GetElementFactory(), Nektar::FieldUtils::GetModuleFactory(), Nektar::NekMeshUtils::Module::m_config, Nektar::NekMeshUtils::Module::m_mesh, CG_Iterations::Mesh, Nektar::NekMeshUtils::Module::ProcessComposites(), and Nektar::NekMeshUtils::Module::ProcessFaces().

65 {
66  if (m_mesh->m_verbose)
67  {
68  cout << "ProcessOptiExtract: ... " << endl;
69  }
70 
71  string ins = m_config["insert"].as<string>();
72 
73  bool extract = boost::iequals(ins, "-1");
74 
75  if (extract)
76  {
77  vector<ElementSharedPtr> el = m_mesh->m_element[m_mesh->m_expDim];
78 
79  m_mesh->m_element[m_mesh->m_expDim].clear();
80  m_mesh->m_element[m_mesh->m_expDim - 1].clear();
81 
82  vector<ElementSharedPtr> invalid;
83 
84  // get invalid elements
85  for (int i = 0; i < el.size(); ++i)
86  {
87  // Create elemental geometry.
89  el[i]->GetGeom(m_mesh->m_spaceDim);
90 
91  // Generate geometric factors.
92  SpatialDomains::GeomFactorsSharedPtr gfac = geom->GetGeomFactors();
93 
94  // Get the Jacobian and, if it is negative, print a warning
95  // message.
96  if (!gfac->IsValid())
97  {
98  invalid.push_back(el[i]);
99  }
100  }
101 
102  std::unordered_set<int> inmesh;
103  vector<ElementSharedPtr> totest;
104 
105  for (int i = 0; i < invalid.size(); i++)
106  {
107  auto t = inmesh.insert(invalid[i]->GetId());
108  if (t.second)
109  {
110  m_mesh->m_element[m_mesh->m_expDim].push_back(invalid[i]);
111  }
112 
113  vector<FaceSharedPtr> f = invalid[i]->GetFaceList();
114  for (int j = 0; j < f.size(); j++)
115  {
116  for (int k = 0; k < f[j]->m_elLink.size(); k++)
117  {
118  if (f[j]->m_elLink[k].first->GetId() == invalid[i]->GetId())
119  {
120  continue;
121  }
122 
123  t = inmesh.insert(f[j]->m_elLink[k].first->GetId());
124  if (t.second)
125  {
126  m_mesh->m_element[m_mesh->m_expDim].push_back(
127  f[j]->m_elLink[k].first);
128  totest.push_back(f[j]->m_elLink[k].first);
129  }
130  }
131  }
132  }
133 
134  for (int i = 0; i < 12; i++)
135  {
136  vector<ElementSharedPtr> tmp = totest;
137  totest.clear();
138  for (int j = 0; j < tmp.size(); j++)
139  {
140  vector<FaceSharedPtr> f = tmp[j]->GetFaceList();
141  for (int k = 0; k < f.size(); k++)
142  {
143  for (int l = 0; l < f[k]->m_elLink.size(); l++)
144  {
145  if (f[k]->m_elLink[l].first->GetId() == tmp[j]->GetId())
146  {
147  continue;
148  }
149 
150  auto t = inmesh.insert(f[k]->m_elLink[l].first->GetId());
151  if (t.second)
152  {
153  m_mesh->m_element[m_mesh->m_expDim].push_back(
154  f[k]->m_elLink[l].first);
155  totest.push_back(f[k]->m_elLink[l].first);
156  }
157  }
158  }
159  }
160  }
161 
163  m_mesh->m_vertexSet.clear();
164  m_mesh->m_edgeSet.clear();
165  m_mesh->m_faceSet.clear();
166 
167  el = m_mesh->m_element[m_mesh->m_expDim];
168 
169  if (m_mesh->m_verbose)
170  {
171  cout << el.size() << " elements in blobs" << endl;
172  }
173 
174  m_mesh->m_faceSet.clear();
175 
176  // re build face links
177  for (int i = 0; i < el.size(); ++i)
178  {
179  for (int j = 0; j < el[i]->GetFaceCount(); ++j)
180  {
181  auto testIns = m_mesh->m_faceSet.insert(el[i]->GetFace(j));
182 
183  if (testIns.second)
184  {
185  (*(testIns.first))
186  ->m_elLink.push_back(
187  pair<ElementSharedPtr, int>(el[i], j));
188  }
189  else
190  {
191  el[i]->SetFace(j, *testIns.first);
192  // Update face to element map.
193  (*(testIns.first))
194  ->m_elLink.push_back(
195  pair<ElementSharedPtr, int>(el[i], j));
196  }
197  }
198  }
199 
200  // build surface composite from faces
201  for (int i = 0; i < el.size(); i++)
202  {
203  vector<FaceSharedPtr> f = el[i]->GetFaceList();
204  for (int j = 0; j < f.size(); j++)
205  {
206  if (f[j]->m_elLink.size() == 1)
207  {
208  // boundary element make new composite
209  ElmtConfig conf(LibUtilities::eTriangle, 1, false, false);
210 
211  vector<int> tags;
212  tags.push_back(1);
215  conf,
216  f[j]->m_vertexList,
217  tags);
218  m_mesh->m_element[m_mesh->m_expDim - 1].push_back(E);
219  }
220  }
221  }
222 
224  for (int i = 0; i < el.size(); ++i)
225  {
226  for (int j = 0; j < el[i]->GetVertexCount(); ++j)
227  {
228  auto testIns = m_mesh->m_vertexSet.insert(el[i]->GetVertex(j));
229 
230  if (!testIns.second)
231  {
232  el[i]->SetVertex(j, *testIns.first);
233  }
234  }
235  }
236 
237  for (int i = 0; i < el.size(); ++i)
238  {
239  for (int j = 0; j < el[i]->GetEdgeCount(); ++j)
240  {
241  EdgeSharedPtr ed = el[i]->GetEdge(j);
242  auto testIns = m_mesh->m_edgeSet.insert(ed);
243 
244  if (testIns.second)
245  {
246  EdgeSharedPtr ed2 = *testIns.first;
247  ed2->m_elLink.push_back(
248  pair<ElementSharedPtr, int>(el[i], j));
249  }
250  else
251  {
252  EdgeSharedPtr e2 = *(testIns.first);
253  el[i]->SetEdge(j, e2);
254  if (e2->m_edgeNodes.size() == 0 &&
255  ed->m_edgeNodes.size() > 0)
256  {
257  e2->m_curveType = ed->m_curveType;
258  e2->m_edgeNodes = ed->m_edgeNodes;
259 
260  // Reverse nodes if appropriate.
261  if (e2->m_n1->m_id != ed->m_n1->m_id)
262  {
263  reverse(e2->m_edgeNodes.begin(),
264  e2->m_edgeNodes.end());
265  }
266  }
267 
268  // Update edge to element map.
269  e2->m_elLink.push_back(
270  pair<ElementSharedPtr, int>(el[i], j));
271  }
272  }
273  }
274  for (int i = 0; i < el.size(); ++i)
275  {
276  for (int j = 0; j < el[i]->GetFaceCount(); ++j)
277  {
278  auto testIns = m_mesh->m_faceSet.insert(el[i]->GetFace(j));
279 
280  if (testIns.second)
281  {
282  (*(testIns.first))
283  ->m_elLink.push_back(
284  pair<ElementSharedPtr, int>(el[i], j));
285  }
286  else
287  {
288  el[i]->SetFace(j, *testIns.first);
289  // Update face to element map.
290  (*(testIns.first))
291  ->m_elLink.push_back(
292  pair<ElementSharedPtr, int>(el[i], j));
293  }
294  }
295  }
296  ProcessFaces(false);
298  }
299  else
300  {
301  // insert other mesh
302  cout << ins << endl;
303  MeshSharedPtr inp_mesh = std::shared_ptr<Mesh>(new Mesh());
305  ModuleKey(eInputModule, "xml"), inp_mesh);
306  mod->RegisterConfig("infile", ins);
307  mod->Process();
308 
309  // need to update the vertices manually then the edges and faces can
310  // be updated using more simple means.
311  map<int, NodeSharedPtr> nmap;
312 
313  for (auto &node : inp_mesh->m_vertexSet)
314  {
315  nmap[node->m_id] = node;
316  }
317  // for all the nodes in the main mesh see if they are in nmap, if so
318  // update the node
319  for (auto &node : m_mesh->m_vertexSet)
320  {
321  if (nmap.count(node->m_id) == 1)
322  {
323  auto s = nmap.find(node->m_id);
324  NodeSharedPtr n = s->second;
325  node->m_x = n->m_x;
326  node->m_y = n->m_y;
327  node->m_z = n->m_z;
328  }
329  }
330 
331  for (auto &eit : inp_mesh->m_edgeSet)
332  {
333  auto et = m_mesh->m_edgeSet.find(eit);
334  if (et != m_mesh->m_edgeSet.end())
335  {
336  (*et)->m_edgeNodes = eit->m_edgeNodes;
337  (*et)->m_curveType = eit->m_curveType;
338  }
339  }
340 
341  for (auto &fit : inp_mesh->m_faceSet)
342  {
343  auto ft = m_mesh->m_faceSet.find(fit);
344  if (ft != m_mesh->m_faceSet.end())
345  {
346  (*ft)->m_faceNodes = fit->m_faceNodes;
347  (*ft)->m_curveType = fit->m_curveType;
348  }
349  }
350  }
351 }
Basic information about an element.
Definition: ElementConfig.h:49
std::shared_ptr< GeomFactors > GeomFactorsSharedPtr
Pointer to a GeomFactors object.
Definition: GeomFactors.h:62
std::shared_ptr< Module > ModuleSharedPtr
std::shared_ptr< Edge > EdgeSharedPtr
Shared pointer to an edge.
Definition: Edge.h:136
std::shared_ptr< Mesh > MeshSharedPtr
Shared pointer to a mesh.
Definition: Mesh.h:156
ElementFactory & GetElementFactory()
Definition: Element.cpp:44
std::shared_ptr< Node > NodeSharedPtr
Definition: CADVert.h:49
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
std::pair< ModuleType, std::string > ModuleKey
virtual NEKMESHUTILS_EXPORT void ProcessFaces(bool ReprocessFaces=true)
Extract element faces.
std::shared_ptr< Element > ElementSharedPtr
Definition: Edge.h:49
std::shared_ptr< Geometry > GeometrySharedPtr
Definition: Geometry.h:53
std::map< std::string, ConfigOption > m_config
List of configuration values.
virtual NEKMESHUTILS_EXPORT void ClearElementLinks()
def extract(self, check_equality=False)
Definition: pycml.py:2657
virtual NEKMESHUTILS_EXPORT void ProcessComposites()
Generate composites.
ModuleFactory & GetModuleFactory()

Member Data Documentation

◆ className

ModuleKey Nektar::Utilities::ProcessOptiExtract::className
static
Initial value:
=
"Pulls out blobs for linear elastic solver.")

Definition at line 53 of file ProcessOptiExtract.h.