Nektar++
Deform.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: Deform.cpp
4 //
5 // For more information, please see: http://www.nektar.info
6 //
7 // The MIT License
8 //
9 // Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10 // Department of Aeronautics, Imperial College London (UK), and Scientific
11 // Computing and Imaging Institute, University of Utah (USA).
12 //
13 // Permission is hereby granted, free of charge, to any person obtaining a
14 // copy of this software and associated documentation files (the "Software"),
15 // to deal in the Software without restriction, including without limitation
16 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
17 // and/or sell copies of the Software, and to permit persons to whom the
18 // Software is furnished to do so, subject to the following conditions:
19 //
20 // The above copyright notice and this permission notice shall be included
21 // in all copies or substantial portions of the Software.
22 //
23 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29 // DEALINGS IN THE SOFTWARE.
30 //
31 // Description: Deformation of mesh from fields.
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
35 #include <string>
36 
37 #include <GlobalMapping/Deform.h>
41 #include <MultiRegions/ExpList.h>
44 #include <StdRegions/StdQuadExp.h>
45 #include <StdRegions/StdSegExp.h>
46 
47 using namespace std;
48 
49 namespace Nektar
50 {
51 namespace GlobalMapping
52 {
53 
54 /**
55  * @brief Update geometry according to displacement that is in current
56  * fields.
57  *
58  * @param graph The MeshGraph of the current geometry.
59  * @param fields The fields containing the displacement.
60  */
63  bool modal)
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();
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] =
91  Array<OneD, NekDouble>(nquad, fields[j]->UpdatePhys() + offset);
92  coord[j] = Array<OneD, NekDouble>(nquad);
93  }
94 
95  // In 2D loop over edges.
96  if (dim == 2)
97  {
98  exp->GetCoords(coord[0], coord[1]);
99 
101  std::dynamic_pointer_cast<SpatialDomains::Geometry2D>(
102  exp->GetGeom());
103 
104  for (j = 0; j < exp->GetGeom()->GetNumEdges(); ++j)
105  {
106  SpatialDomains::Geometry1DSharedPtr edge = geom->GetEdge(j);
107 
108  // This edge has already been processed.
109  if (updatedEdges.find(edge->GetGlobalID()) !=
110  updatedEdges.end())
111  {
112  continue;
113  }
114 
115  // Extract edge displacement.
116  int nEdgePts = exp->GetTraceNumPoints(j);
117  Array<OneD, Array<OneD, NekDouble>> edgePhys(dim);
118  Array<OneD, Array<OneD, NekDouble>> edgeCoord(dim);
119 
120  const LibUtilities::BasisKey B(
121  LibUtilities::eModified_A, nEdgePts,
126 
127  for (k = 0; k < dim; ++k)
128  {
129  edgePhys[k] = Array<OneD, NekDouble>(nEdgePts);
130  edgeCoord[k] = Array<OneD, NekDouble>(nEdgePts);
131  exp->GetTracePhysVals(j, seg, phys[k], edgePhys[k]);
132  exp->GetTracePhysVals(j, seg, coord[k], edgeCoord[k]);
133  }
134 
135  // Update verts
136  for (k = 0; k < 2; ++k)
137  {
138  int id = edge->GetVid(k);
139  if (updatedVerts.find(id) != updatedVerts.end())
140  {
141  continue;
142  }
143 
144  SpatialDomains::PointGeomSharedPtr pt = edge->GetVertex(k);
145 
146  pt->UpdatePosition(
147  (*pt)(0) + edgePhys[0][k * (nEdgePts - 1)],
148  (*pt)(1) + edgePhys[1][k * (nEdgePts - 1)], (*pt)(2));
149 
150  updatedVerts.insert(id);
151  }
152 
153  // Update curve
156  edge->GetGlobalID(),
158 
159  for (k = 0; k < nEdgePts; ++k)
160  {
163  AllocateSharedPtr(dim, edge->GetGlobalID(),
164  edgeCoord[0][k] + edgePhys[0][k],
165  edgeCoord[1][k] + edgePhys[1][k],
166  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->GetNtraces(); ++j)
184  {
185  SpatialDomains::Geometry2DSharedPtr face = geom->GetFace(j);
186 
188  exp->as<LocalRegions::Expansion3D>();
189 
190  // This edge has already been processed.
191  if (updatedFaces.find(face->GetGlobalID()) !=
192  updatedFaces.end())
193  {
194  continue;
195  }
196 
197  // Extract face displacement.
198  LibUtilities::BasisKey B0 = exp->GetTraceBasisKey(j, 0);
199  LibUtilities::BasisKey B1 = exp->GetTraceBasisKey(j, 1);
200  int nq0 = B0.GetNumPoints();
201  int nq1 = B1.GetNumPoints();
202 
203  ASSERTL1(B0.GetPointsType() ==
205  B1.GetPointsType() ==
207  "Deformation requires GLL points in both "
208  "directions on a face.");
209 
211 
213  StdRegions::Orientation orient = exp->GetTraceOrient(j);
214 
215  if (face->GetShapeType() == LibUtilities::eTriangle)
216  {
217  faceexp =
219  B0, B1);
220  }
221  else
222  {
223  faceexp = MemoryManager<
224  StdRegions::StdQuadExp>::AllocateSharedPtr(B0, B1);
225  }
226 
227  for (k = 0; k < dim; ++k)
228  {
229  Array<OneD, NekDouble> tmp(nq0 * nq1);
230  newPos[k] = Array<OneD, NekDouble>(nq0 * nq1);
231  exp3d->GetTracePhysVals(j, faceexp, phys[k], tmp, orient);
232  exp3d->GetTracePhysVals(j, faceexp, coord[k], newPos[k],
233  orient);
234  Vmath::Vadd(nq0 * nq1, tmp, 1, newPos[k], 1, newPos[k], 1);
235  }
236 
237  // Now interpolate face onto a more reasonable set of
238  // points.
239  int nq = max(nq0, nq1);
240  if (!modal)
241  nq--;
242 
243  LibUtilities::PointsKey edgePts(
246 
248 
249  for (k = 0; k < dim; ++k)
250  {
251  intPos[k] = Array<OneD, NekDouble>(nq * nq);
252  LibUtilities::Interp2D(faceexp->GetPointsKeys()[0],
253  faceexp->GetPointsKeys()[1],
254  newPos[k], edgePts, edgePts,
255  intPos[k]);
256  }
257 
258  int edgeOff[2][4][2] = {
259  {{0, 1}, {nq - 1, nq}, {nq * (nq - 1), -nq}, {-1, -1}},
260  {{0, 1},
261  {nq - 1, nq},
262  {nq * nq - 1, -1},
263  {nq * (nq - 1), -nq}}};
264 
265  for (k = 0; k < face->GetNumVerts(); ++k)
266  {
267  // Update verts
268  int id = face->GetVid(k);
269  const int o =
270  face->GetShapeType() - LibUtilities::eTriangle;
271 
272  if (updatedVerts.find(id) == updatedVerts.end())
273  {
275  face->GetVertex(k);
276  pt->UpdatePosition(intPos[0][edgeOff[o][k][0]],
277  intPos[1][edgeOff[o][k][0]],
278  intPos[2][edgeOff[o][k][0]]);
279  updatedVerts.insert(id);
280  }
281 
282  // Update edges
283  id = face->GetEid(k);
284  if (updatedEdges.find(id) == updatedEdges.end())
285  {
287  face->GetEdge(k);
291  edge->GetGlobalID(),
293 
294  const int offset = edgeOff[o][k][0];
295  const int pos = edgeOff[o][k][1];
296 
297  if (face->GetEorient(k) == StdRegions::eBackwards)
298  {
299  for (l = nq - 1; l >= 0; --l)
300  {
301  int m = offset + pos * l;
305  dim, edge->GetGlobalID(),
306  intPos[0][m], intPos[1][m],
307  intPos[2][m]);
308  curve->m_points.push_back(vert);
309  }
310  }
311  else
312  {
313  for (l = 0; l < nq; ++l)
314  {
315  int m = offset + pos * l;
319  dim, edge->GetGlobalID(),
320  intPos[0][m], intPos[1][m],
321  intPos[2][m]);
322  curve->m_points.push_back(vert);
323  }
324  }
325 
326  curvedEdges[edge->GetGlobalID()] = curve;
327  updatedEdges.insert(edge->GetGlobalID());
328  }
329  }
330 
331  // Update face-interior curvature
333  face->GetShapeType() == LibUtilities::eTriangle
336 
339  face->GetGlobalID(), pType);
340 
341  if (face->GetShapeType() == LibUtilities::eTriangle)
342  {
343  // This code is probably pretty crappy. Have to go from
344  // GLL-GLL points -> GLL-Gauss-Radau -> nodal triangle
345  // points.
346  const LibUtilities::BasisKey B0(
350  const LibUtilities::BasisKey B1(
353  nq, LibUtilities::eGaussRadauMAlpha1Beta0));
354  StdRegions::StdNodalTriExp nodalTri(B0, B1, pType);
355  StdRegions::StdTriExp tri(B0, B1);
356 
357  for (k = 0; k < dim; ++k)
358  {
359  Array<OneD, NekDouble> nodal(nq * nq);
360 
362  faceexp->GetBasis(0)->GetBasisKey(),
363  faceexp->GetBasis(1)->GetBasisKey(), newPos[k], B0,
364  B1, nodal);
365 
366  Array<OneD, NekDouble> tmp1(nq * (nq + 1) / 2);
367  Array<OneD, NekDouble> tmp2(nq * (nq + 1) / 2);
368 
369  tri.FwdTrans(nodal, tmp1);
370  nodalTri.ModalToNodal(tmp1, tmp2);
371  newPos[k] = tmp2;
372  }
373 
374  for (l = 0; l < nq * (nq + 1) / 2; ++l)
375  {
378  AllocateSharedPtr(dim, face->GetGlobalID(),
379  newPos[0][l], newPos[1][l],
380  newPos[2][l]);
381  curve->m_points.push_back(vert);
382  }
383  }
384  else
385  {
386  for (l = 0; l < nq * nq; ++l)
387  {
390  AllocateSharedPtr(dim, face->GetGlobalID(),
391  intPos[0][l], intPos[1][l],
392  intPos[2][l]);
393  curve->m_points.push_back(vert);
394  }
395  }
396 
397  curvedFaces[face->GetGlobalID()] = curve;
398  updatedFaces.insert(face->GetGlobalID());
399  }
400  }
401  }
402 
403  // Reset geometry information
404  for (i = 0; i < fields.size(); ++i)
405  {
406  fields[i]->Reset();
407  }
408 }
409 } // namespace GlobalMapping
410 } // namespace Nektar
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
Describes the specification for a Basis.
Definition: Basis.h:50
int GetNumPoints() const
Return points order at which basis is defined.
Definition: Basis.h:130
PointsType GetPointsType() const
Return type of quadrature.
Definition: Basis.h:153
Defines a specification for a set of points.
Definition: Points.h:59
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
void FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Forward transformation from physical space to coefficient space.
void ModalToNodal(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void UpdateGeometry(SpatialDomains::MeshGraphSharedPtr graph, Array< OneD, MultiRegions::ExpListSharedPtr > &fields, bool modal)
Update geometry according to displacement that is in current fields.
Definition: Deform.cpp:61
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< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:172
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
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
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