Nektar++
Classes | Typedefs | Functions | Variables
Nektar::GlobalMapping Namespace Reference

Classes

class  Mapping
 Base class for mapping to be applied to the coordinate system. More...
 
class  MappingGeneral
 
class  MappingTranslation
 
class  MappingXofXZ
 
class  MappingXofZ
 
class  MappingXYofXY
 
class  MappingXYofZ
 

Typedefs

typedef LibUtilities::NekFactory< std::string, Mapping, const LibUtilities::SessionReaderSharedPtr &, const Array< OneD, MultiRegions::ExpListSharedPtr > &, const TiXmlElement * > MappingFactory
 Declaration of the mapping factory. More...
 

Functions

void UpdateGeometry (SpatialDomains::MeshGraphSharedPtr graph, Array< OneD, MultiRegions::ExpListSharedPtr > &fields, bool modal)
 Update geometry according to displacement that is in current fields. More...
 
MappingFactoryGetMappingFactory ()
 Declaration of the mapping factory singleton. More...
 

Variables

GLOBAL_MAPPING_EXPORT typedef std::shared_ptr< MappingMappingSharedPtr
 A shared pointer to a Mapping object. More...
 

Typedef Documentation

◆ MappingFactory

Declaration of the mapping factory.

Definition at line 59 of file Mapping.h.

Function Documentation

◆ GetMappingFactory()

GLOBAL_MAPPING_EXPORT MappingFactory & Nektar::GlobalMapping::GetMappingFactory ( )

Declaration of the mapping factory singleton.

Definition at line 52 of file Mapping.cpp.

Referenced by Nektar::GlobalMapping::Mapping::Load().

53 {
54  static MappingFactory instance;
55  return instance;
56 }
LibUtilities::NekFactory< std::string, Mapping, const LibUtilities::SessionReaderSharedPtr &, const Array< OneD, MultiRegions::ExpListSharedPtr > &, const TiXmlElement * > MappingFactory
Declaration of the mapping factory.
Definition: Mapping.h:59

◆ UpdateGeometry()

GLOBAL_MAPPING_EXPORT void Nektar::GlobalMapping::UpdateGeometry ( SpatialDomains::MeshGraphSharedPtr  graph,
Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
bool  modal 
)

Update geometry according to displacement that is in current fields.

Adds a summary item to the summary info list.

Parameters
graphThe MeshGraph of the current geometry.
fieldsThe fields containing the displacement.

Definition at line 58 of file Deform.cpp.

References ASSERTL1, Nektar::StdRegions::eBackwards, Nektar::LibUtilities::eGaussLobattoLegendre, Nektar::LibUtilities::eGaussRadauMAlpha1Beta0, Nektar::LibUtilities::eModified_A, Nektar::LibUtilities::eNodalTriElec, Nektar::LibUtilities::eOrtho_A, Nektar::LibUtilities::eOrtho_B, Nektar::LibUtilities::eTriangle, Nektar::StdRegions::StdExpansion::FwdTrans(), Nektar::SpatialDomains::Geometry::GetEdge(), Nektar::SpatialDomains::Geometry::GetFace(), Nektar::LibUtilities::BasisKey::GetNumPoints(), Nektar::LibUtilities::BasisKey::GetPointsType(), Nektar::LibUtilities::Interp2D(), Nektar::StdRegions::StdNodalTriExp::ModalToNodal(), and Vmath::Vadd().

Referenced by Nektar::FieldUtils::ProcessDeform::Process(), and Nektar::IterativeElasticSystem::v_DoSolve().

62  {
63  // Clear existing curvature.
64  SpatialDomains::CurveMap &curvedEdges = graph->GetCurvedEdges();
65  SpatialDomains::CurveMap &curvedFaces = graph->GetCurvedFaces();
66  curvedEdges.clear();
67  curvedFaces.clear();
68 
69  int i, j, k, l, dim;
70 
71  // Sets to hold IDs of updated vertices to avoid duplicating effort.
72  set<int> updatedVerts, updatedEdges, updatedFaces;
73 
74  dim = graph->GetSpaceDimension();
75  Array<OneD, Array<OneD, NekDouble> > phys (dim);
76  Array<OneD, Array<OneD, NekDouble> > coord(dim);
77 
78  for (i = 0; i < fields[0]->GetExpSize(); ++i)
79  {
80  LocalRegions::ExpansionSharedPtr exp = fields[0]->GetExp(i);
81  int offset = fields[0]->GetPhys_Offset(i);
82  int nquad = exp->GetTotPoints();
83 
84  // Extract displacement for this element, allocate storage for
85  // elemental coordinates.
86  for (j = 0; j < dim; ++j)
87  {
88  phys[j] = Array<OneD, NekDouble>(
89  nquad, fields[j]->UpdatePhys() + offset);
90  coord[j] = Array<OneD, NekDouble>(nquad);
91  }
92 
93  // In 2D loop over edges.
94  if (dim == 2)
95  {
96  exp->GetCoords(coord[0], coord[1]);
97 
99  std::dynamic_pointer_cast<SpatialDomains::Geometry2D>(
100  exp->GetGeom());
101 
102  for (j = 0; j < exp->GetNedges(); ++j)
103  {
104  SpatialDomains::Geometry1DSharedPtr edge = geom->GetEdge(j);
105 
106  // This edge has already been processed.
107  if (updatedEdges.find(edge->GetGlobalID()) !=
108  updatedEdges.end())
109  {
110  continue;
111  }
112 
113  // Extract edge displacement.
114  int nEdgePts = exp->GetEdgeNumPoints(j);
115  Array<OneD, Array<OneD, NekDouble> > edgePhys (dim);
116  Array<OneD, Array<OneD, NekDouble> > edgeCoord(dim);
117 
118  const LibUtilities::BasisKey B(
119  LibUtilities::eModified_A, nEdgePts,
120  LibUtilities::PointsKey(
122  StdRegions::StdExpansion1DSharedPtr seg = MemoryManager<
123  StdRegions::StdSegExp>::AllocateSharedPtr(B);
124 
125  for (k = 0; k < dim; ++k)
126  {
127  edgePhys [k] = Array<OneD, NekDouble>(nEdgePts);
128  edgeCoord[k] = Array<OneD, NekDouble>(nEdgePts);
129  exp->GetEdgePhysVals(j, seg, phys [k], edgePhys [k]);
130  exp->GetEdgePhysVals(j, seg, coord[k], edgeCoord[k]);
131  }
132 
133  // Update verts
134  for (k = 0; k < 2; ++k)
135  {
136  int id = edge->GetVid(k);
137  if (updatedVerts.find(id) != updatedVerts.end())
138  {
139  continue;
140  }
141 
143  edge->GetVertex(k);
144 
145  pt->UpdatePosition(
146  (*pt)(0) + edgePhys[0][k*(nEdgePts-1)],
147  (*pt)(1) + edgePhys[1][k*(nEdgePts-1)],
148  (*pt)(2));
149 
150  updatedVerts.insert(id);
151  }
152 
153  // Update curve
154  SpatialDomains::CurveSharedPtr curve = MemoryManager<
155  SpatialDomains::Curve>::AllocateSharedPtr(
156  edge->GetGlobalID(),
158 
159  for (k = 0; k < nEdgePts; ++k)
160  {
162  MemoryManager<SpatialDomains::PointGeom>
163  ::AllocateSharedPtr(
164  dim, edge->GetGlobalID(),
165  edgeCoord[0][k] + edgePhys[0][k],
166  edgeCoord[1][k] + edgePhys[1][k], 0.0);
167 
168  curve->m_points.push_back(vert);
169  }
170 
171  curvedEdges[edge->GetGlobalID()] = curve;
172  updatedEdges.insert(edge->GetGlobalID());
173  }
174  }
175  else if (dim == 3)
176  {
177  exp->GetCoords(coord[0], coord[1], coord[2]);
178 
180  std::dynamic_pointer_cast<SpatialDomains::Geometry3D>(
181  exp->GetGeom());
182 
183  for (j = 0; j < exp->GetNfaces(); ++j)
184  {
185  SpatialDomains::Geometry2DSharedPtr face = geom->GetFace(j);
186 
187  // This edge has already been processed.
188  if (updatedFaces.find(face->GetGlobalID()) !=
189  updatedFaces.end())
190  {
191  continue;
192  }
193 
194  // Extract face displacement.
195  LibUtilities::BasisKey B0 = exp->DetFaceBasisKey(j,0);
196  LibUtilities::BasisKey B1 = exp->DetFaceBasisKey(j,1);
197  int nq0 = B0.GetNumPoints();
198  int nq1 = B1.GetNumPoints();
199 
200  ASSERTL1(B0.GetPointsType()
202  B1.GetPointsType()
204  "Deformation requires GLL points in both "
205  "directions on a face.");
206 
207  Array<OneD, Array<OneD, NekDouble> > newPos(dim);
208 
210  StdRegions::Orientation orient = exp->GetForient(j);
211 
212  if (face->GetShapeType() == LibUtilities::eTriangle)
213  {
214  faceexp = MemoryManager<StdRegions::StdTriExp>::
215  AllocateSharedPtr(B0, B1);
216  }
217  else
218  {
219  faceexp = MemoryManager<StdRegions::StdQuadExp>::
220  AllocateSharedPtr(B0, B1);
221  }
222 
223  for (k = 0; k < dim; ++k)
224  {
225  Array<OneD, NekDouble> tmp(nq0*nq1);
226  newPos[k] = Array<OneD, NekDouble>(nq0*nq1);
227  exp->GetFacePhysVals(
228  j, faceexp, phys [k], tmp, orient);
229  exp->GetFacePhysVals(
230  j, faceexp, coord[k], newPos[k], orient);
231  Vmath::Vadd(
232  nq0*nq1, tmp, 1, newPos[k], 1, newPos[k], 1);
233  }
234 
235  // Now interpolate face onto a more reasonable set of
236  // points.
237  int nq = max(nq0, nq1);
238  if(!modal)
239  nq--;
240 
241  LibUtilities::PointsKey edgePts(
243  LibUtilities::PointsKey triPts(
245 
246  Array<OneD, Array<OneD, NekDouble> > intPos(dim);
247 
248  for (k = 0; k < dim; ++k)
249  {
250  intPos[k] = Array<OneD, NekDouble>(nq*nq);
252  faceexp->GetPointsKeys()[0],
253  faceexp->GetPointsKeys()[1],
254  newPos[k], edgePts, edgePts, intPos[k]);
255  }
256 
257  int edgeOff[2][4][2] = {
258  {
259  {0, 1},
260  {nq-1, nq},
261  {nq*(nq-1), -nq},
262  {-1,-1}
263  },
264  {
265  {0, 1},
266  {nq-1, nq},
267  {nq*nq-1, -1},
268  {nq*(nq-1), -nq}
269  }
270  };
271 
272  for (k = 0; k < face->GetNumVerts(); ++k)
273  {
274  // Update verts
275  int id = face->GetVid(k);
276  const int o =
277  face->GetShapeType() - LibUtilities::eTriangle;
278 
279  if (updatedVerts.find(id) == updatedVerts.end())
280  {
282  face->GetVertex(k);
283  pt->UpdatePosition(
284  intPos[0][edgeOff[o][k][0]],
285  intPos[1][edgeOff[o][k][0]],
286  intPos[2][edgeOff[o][k][0]]);
287  updatedVerts.insert(id);
288  }
289 
290  // Update edges
291  id = face->GetEid(k);
292  if (updatedEdges.find(id) == updatedEdges.end())
293  {
295  = face->GetEdge(k);
297  MemoryManager<SpatialDomains::Curve>
298  ::AllocateSharedPtr(
299  edge->GetGlobalID(),
301 
302  const int offset = edgeOff[o][k][0];
303  const int pos = edgeOff[o][k][1];
304 
305  if (face->GetEorient(k) == StdRegions::eBackwards)
306  {
307  for (l = nq-1; l >= 0; --l)
308  {
309  int m = offset + pos*l;
311  MemoryManager<SpatialDomains::PointGeom>
312  ::AllocateSharedPtr(
313  dim, edge->GetGlobalID(),
314  intPos[0][m], intPos[1][m],
315  intPos[2][m]);
316  curve->m_points.push_back(vert);
317  }
318  }
319  else
320  {
321  for (l = 0; l < nq; ++l)
322  {
323  int m = offset + pos*l;
325  MemoryManager<SpatialDomains::PointGeom>
326  ::AllocateSharedPtr(
327  dim, edge->GetGlobalID(),
328  intPos[0][m], intPos[1][m],
329  intPos[2][m]);
330  curve->m_points.push_back(vert);
331  }
332  }
333 
334  curvedEdges[edge->GetGlobalID()] = curve;
335  updatedEdges.insert(edge->GetGlobalID());
336  }
337  }
338 
339  // Update face-interior curvature
341  face->GetShapeType() == LibUtilities::eTriangle ?
344 
345  SpatialDomains::CurveSharedPtr curve = MemoryManager<
346  SpatialDomains::Curve>::AllocateSharedPtr(
347  face->GetGlobalID(),
348  pType);
349 
350  if (face->GetShapeType() == LibUtilities::eTriangle)
351  {
352  // This code is probably pretty crappy. Have to go from
353  // GLL-GLL points -> GLL-Gauss-Radau -> nodal triangle
354  // points.
355  const LibUtilities::BasisKey B0(
357  LibUtilities::PointsKey(
359  const LibUtilities::BasisKey B1(
361  LibUtilities::PointsKey(
363  StdRegions::StdNodalTriExp nodalTri(B0, B1, pType);
364  StdRegions::StdTriExp tri (B0, B1);
365 
366  for (k = 0; k < dim; ++k)
367  {
368  Array<OneD, NekDouble> nodal(nq*nq);
369 
371  faceexp->GetBasis(0)->GetBasisKey(),
372  faceexp->GetBasis(1)->GetBasisKey(),
373  newPos[k], B0, B1, nodal);
374 
375  Array<OneD, NekDouble> tmp1(nq*(nq+1)/2);
376  Array<OneD, NekDouble> tmp2(nq*(nq+1)/2);
377 
378  tri.FwdTrans(nodal, tmp1);
379  nodalTri.ModalToNodal(tmp1, tmp2);
380  newPos[k] = tmp2;
381  }
382 
383  for (l = 0; l < nq*(nq+1)/2; ++l)
384  {
386  MemoryManager<SpatialDomains::PointGeom>
387  ::AllocateSharedPtr(
388  dim, face->GetGlobalID(),
389  newPos[0][l], newPos[1][l], newPos[2][l]);
390  curve->m_points.push_back(vert);
391  }
392  }
393  else
394  {
395  for (l = 0; l < nq*nq; ++l)
396  {
398  MemoryManager<SpatialDomains::PointGeom>
399  ::AllocateSharedPtr(
400  dim, face->GetGlobalID(),
401  intPos[0][l], intPos[1][l], intPos[2][l]);
402  curve->m_points.push_back(vert);
403  }
404  }
405 
406  curvedFaces[face->GetGlobalID()] = curve;
407  updatedFaces.insert(face->GetGlobalID());
408  }
409  }
410  }
411 
412  // Reset geometry information
413  for (i = 0; i < fields.num_elements(); ++i)
414  {
415  fields[i]->Reset();
416  }
417  }
std::shared_ptr< Geometry3D > Geometry3DSharedPtr
Definition: Geometry3D.h:52
std::shared_ptr< Geometry2D > Geometry2DSharedPtr
Definition: Geometry.h:65
std::unordered_map< int, CurveSharedPtr > CurveMap
Definition: Curve.hpp:62
Principle Modified Functions .
Definition: BasisType.h:48
std::shared_ptr< StdExpansion2D > StdExpansion2DSharedPtr
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:58
Principle Orthogonal Functions .
Definition: BasisType.h:46
void Interp2D(const BasisKey &fbasis0, const BasisKey &fbasis1, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, Array< OneD, NekDouble > &to)
this function interpolates a 2D function evaluated at the quadrature points of the 2D basis...
Definition: Interp.cpp:115
std::shared_ptr< StdExpansion1D > StdExpansion1DSharedPtr
Principle Orthogonal Functions .
Definition: BasisType.h:45
std::shared_ptr< PointGeom > PointGeomSharedPtr
Definition: Geometry.h:59
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:65
std::shared_ptr< Curve > CurveSharedPtr
Definition: Curve.hpp:61
std::shared_ptr< Geometry1D > Geometry1DSharedPtr
Definition: Geometry.h:63
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:250
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:51
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.cpp:302
2D Nodal Electrostatic Points on a Triangle
Definition: PointsType.h:69

Variable Documentation

◆ MappingSharedPtr

GLOBAL_MAPPING_EXPORT typedef std::shared_ptr<Mapping> Nektar::GlobalMapping::MappingSharedPtr