Nektar++
Loading...
Searching...
No Matches
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.
 

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.
 
MappingFactoryGetMappingFactory ()
 Declaration of the mapping factory singleton.
 

Variables

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

Typedef Documentation

◆ MappingFactory

Declaration of the mapping factory.

Definition at line 63 of file Mapping.h.

Function Documentation

◆ GetMappingFactory()

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

Declaration of the mapping factory singleton.

Definition at line 47 of file Mapping.cpp.

48{
49 static MappingFactory instance;
50 return instance;
51}
Provides a generic Factory class.

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 53 of file Deform.cpp.

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

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::ObjPoolManager< DataType >::AllocateUniquePtr(), 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::SpatialDomains::Geometry::GetEdge(), Nektar::SpatialDomains::Geometry::GetEid(), Nektar::SpatialDomains::Geometry::GetEorient(), Nektar::SpatialDomains::Geometry::GetFace(), Nektar::SpatialDomains::Geometry::GetGlobalID(), Nektar::LibUtilities::BasisKey::GetNumPoints(), Nektar::SpatialDomains::Geometry::GetNumVerts(), Nektar::LibUtilities::BasisKey::GetPointsType(), Nektar::SpatialDomains::Geometry::GetShapeType(), Nektar::SpatialDomains::Geometry::GetVertex(), Nektar::SpatialDomains::Geometry::GetVid(), Nektar::LibUtilities::Interp2D(), Nektar::StdRegions::StdNodalTriExp::ModalToNodal(), Nektar::SpatialDomains::PointGeom::UpdatePosition(), 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