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, Array< OneD, Array< OneD, NekDouble >> &PhysVals, 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 53 of file Mapping.cpp.

54 {
55  static MappingFactory instance;
56  return instance;
57 }
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

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

◆ UpdateGeometry()

GLOBAL_MAPPING_EXPORT void Nektar::GlobalMapping::UpdateGeometry ( SpatialDomains::MeshGraphSharedPtr  graph,
Array< OneD, MultiRegions::ExpListSharedPtr > &  fields,
Array< OneD, Array< OneD, NekDouble >> &  PhysVals,
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 61 of file Deform.cpp.

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

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

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

Variable Documentation

◆ MappingSharedPtr

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