Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 // License for the specific language governing rights and limitations under
14 // Permission is hereby granted, free of charge, to any person obtaining a
15 // copy of this software and associated documentation files (the "Software"),
16 // to deal in the Software without restriction, including without limitation
17 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
18 // and/or sell copies of the Software, and to permit persons to whom the
19 // Software is furnished to do so, subject to the following conditions:
20 //
21 // The above copyright notice and this permission notice shall be included
22 // in all copies or substantial portions of the Software.
23 //
24 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30 // DEALINGS IN THE SOFTWARE.
31 //
32 // Description: Deformation of mesh from fields.
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
36 #include <string>
37 
42 #include <StdRegions/StdSegExp.h>
43 #include <StdRegions/StdQuadExp.h>
45 #include <MultiRegions/ExpList.h>
46 
47 namespace Nektar {
48 namespace SolverUtils {
49 
50  /**
51  * @brief Update geometry according to displacement that is in current
52  * fields.
53  *
54  * @param graph The MeshGraph of the current geometry.
55  * @param fields The fields containing the displacement.
56  */
60  bool modal)
61  {
62  // Clear existing curvature.
63  SpatialDomains::CurveMap &curvedEdges = graph->GetCurvedEdges();
64  SpatialDomains::CurveMap &curvedFaces = graph->GetCurvedFaces();
65  curvedEdges.clear();
66  curvedFaces.clear();
67 
68  int i, j, k, l, dim;
69 
70  // Sets to hold IDs of updated vertices to avoid duplicating effort.
71  set<int> updatedVerts, updatedEdges, updatedFaces;
72 
73  dim = graph->GetSpaceDimension();
76 
77  for (i = 0; i < fields[0]->GetExpSize(); ++i)
78  {
79  LocalRegions::ExpansionSharedPtr exp = fields[0]->GetExp(i);
80  int offset = fields[0]->GetPhys_Offset(i);
81  int nquad = exp->GetTotPoints();
82 
83  // Extract displacement for this element, allocate storage for
84  // elemental coordinates.
85  for (j = 0; j < dim; ++j)
86  {
87  phys[j] = Array<OneD, NekDouble>(
88  nquad, fields[j]->UpdatePhys() + offset);
89  coord[j] = Array<OneD, NekDouble>(nquad);
90  }
91 
92  // In 2D loop over edges.
93  if (dim == 2)
94  {
95  exp->GetCoords(coord[0], coord[1]);
96 
98  boost::dynamic_pointer_cast<SpatialDomains::Geometry2D>(
99  exp->GetGeom());
100 
101  for (j = 0; j < exp->GetNedges(); ++j)
102  {
104 
105  // This edge has already been processed.
106  if (updatedEdges.find(edge->GetGlobalID()) !=
107  updatedEdges.end())
108  {
109  continue;
110  }
111 
112  // Extract edge displacement.
113  int nEdgePts = exp->GetEdgeNumPoints(j);
114  Array<OneD, Array<OneD, NekDouble> > edgePhys (dim);
115  Array<OneD, Array<OneD, NekDouble> > edgeCoord(dim);
116 
117  const LibUtilities::BasisKey B(
118  LibUtilities::eModified_A, nEdgePts,
122  StdRegions::StdSegExp>::AllocateSharedPtr(B);
123 
124  for (k = 0; k < dim; ++k)
125  {
126  edgePhys [k] = Array<OneD, NekDouble>(nEdgePts);
127  edgeCoord[k] = Array<OneD, NekDouble>(nEdgePts);
128  exp->GetEdgePhysVals(j, seg, phys [k], edgePhys [k]);
129  exp->GetEdgePhysVals(j, seg, coord[k], edgeCoord[k]);
130  }
131 
132  // Update verts
133  for (k = 0; k < 2; ++k)
134  {
135  int id = edge->GetVid(k);
136  if (updatedVerts.find(id) != updatedVerts.end())
137  {
138  continue;
139  }
140 
142  edge->GetVertex(k);
143 
144  pt->UpdatePosition(
145  (*pt)(0) + edgePhys[0][k*(nEdgePts-1)],
146  (*pt)(1) + edgePhys[1][k*(nEdgePts-1)],
147  (*pt)(2));
148 
149  updatedVerts.insert(id);
150  }
151 
152  // Update curve
154  SpatialDomains::Curve>::AllocateSharedPtr(
155  edge->GetGlobalID(),
157 
158  for (k = 0; k < nEdgePts; ++k)
159  {
163  dim, edge->GetGlobalID(),
164  edgeCoord[0][k] + edgePhys[0][k],
165  edgeCoord[1][k] + edgePhys[1][k], 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  boost::dynamic_pointer_cast<SpatialDomains::Geometry3D>(
180  exp->GetGeom());
181 
182  for (j = 0; j < exp->GetNfaces(); ++j)
183  {
185 
186  // This edge has already been processed.
187  if (updatedFaces.find(face->GetGlobalID()) !=
188  updatedFaces.end())
189  {
190  continue;
191  }
192 
193  // Extract face displacement.
194  LibUtilities::BasisKey B0 = exp->DetFaceBasisKey(j,0);
195  LibUtilities::BasisKey B1 = exp->DetFaceBasisKey(j,1);
196  int nq0 = B0.GetNumPoints();
197  int nq1 = B1.GetNumPoints();
198 
201  B1.GetPointsType()
203  "Deformation requires GLL points in both "
204  "directions on a face.");
205 
207 
209  StdRegions::Orientation orient = exp->GetForient(j);
210 
211  if (face->GetShapeType() == LibUtilities::eTriangle)
212  {
214  AllocateSharedPtr(B0, B1);
215  }
216  else
217  {
219  AllocateSharedPtr(B0, B1);
220  }
221 
222  for (k = 0; k < dim; ++k)
223  {
224  Array<OneD, NekDouble> tmp(nq0*nq1);
225  newPos[k] = Array<OneD, NekDouble>(nq0*nq1);
226  exp->GetFacePhysVals(
227  j, faceexp, phys [k], tmp, orient);
228  exp->GetFacePhysVals(
229  j, faceexp, coord[k], newPos[k], orient);
230  Vmath::Vadd(
231  nq0*nq1, tmp, 1, newPos[k], 1, newPos[k], 1);
232  }
233 
234  // Now interpolate face onto a more reasonable set of
235  // points.
236  int nq = max(nq0, nq1);
237  if(!modal)
238  nq--;
239 
240  LibUtilities::PointsKey edgePts(
244 
246 
247  for (k = 0; k < dim; ++k)
248  {
249  intPos[k] = Array<OneD, NekDouble>(nq*nq);
251  faceexp->GetPointsKeys()[0],
252  faceexp->GetPointsKeys()[1],
253  newPos[k], edgePts, edgePts, intPos[k]);
254  }
255 
256  int edgeOff[2][4][2] = {
257  {
258  {0, 1},
259  {nq-1, nq},
260  {nq*(nq-1), -nq},
261  {-1,-1}
262  },
263  {
264  {0, 1},
265  {nq-1, nq},
266  {nq*nq-1, -1},
267  {nq*(nq-1), -nq}
268  }
269  };
270 
271  for (k = 0; k < face->GetNumVerts(); ++k)
272  {
273  // Update verts
274  int id = face->GetVid(k);
275  const int o =
276  face->GetShapeType() - LibUtilities::eTriangle;
277 
278  if (updatedVerts.find(id) == updatedVerts.end())
279  {
281  face->GetVertex(k);
282  pt->UpdatePosition(
283  intPos[0][edgeOff[o][k][0]],
284  intPos[1][edgeOff[o][k][0]],
285  intPos[2][edgeOff[o][k][0]]);
286  updatedVerts.insert(id);
287  }
288 
289  // Update edges
290  id = face->GetEid(k);
291  if (updatedEdges.find(id) == updatedEdges.end())
292  {
294  = face->GetEdge(k);
298  edge->GetGlobalID(),
300 
301  const int offset = edgeOff[o][k][0];
302  const int pos = edgeOff[o][k][1];
303 
304  if (face->GetEorient(k) == StdRegions::eBackwards)
305  {
306  for (l = nq-1; l >= 0; --l)
307  {
308  int m = offset + pos*l;
312  dim, edge->GetGlobalID(),
313  intPos[0][m], intPos[1][m],
314  intPos[2][m]);
315  curve->m_points.push_back(vert);
316  }
317  }
318  else
319  {
320  for (l = 0; l < nq; ++l)
321  {
322  int m = offset + pos*l;
326  dim, edge->GetGlobalID(),
327  intPos[0][m], intPos[1][m],
328  intPos[2][m]);
329  curve->m_points.push_back(vert);
330  }
331  }
332 
333  curvedEdges[edge->GetGlobalID()] = curve;
334  updatedEdges.insert(edge->GetGlobalID());
335  }
336  }
337 
338  // Update face-interior curvature
340  face->GetShapeType() == LibUtilities::eTriangle ?
343 
345  SpatialDomains::Curve>::AllocateSharedPtr(
346  face->GetGlobalID(),
347  pType);
348 
349  if (face->GetShapeType() == LibUtilities::eTriangle)
350  {
351  // This code is probably pretty crappy. Have to go from
352  // GLL-GLL points -> GLL-Gauss-Radau -> nodal triangle
353  // points.
354  const LibUtilities::BasisKey B0(
358  const LibUtilities::BasisKey B1(
362  StdRegions::StdNodalTriExp nodalTri(B0, B1, pType);
363  StdRegions::StdTriExp tri (B0, B1);
364 
365  for (k = 0; k < dim; ++k)
366  {
367  Array<OneD, NekDouble> nodal(nq*nq);
368 
370  faceexp->GetBasis(0)->GetBasisKey(),
371  faceexp->GetBasis(1)->GetBasisKey(),
372  newPos[k], B0, B1, nodal);
373 
374  Array<OneD, NekDouble> tmp1(nq*(nq+1)/2);
375  Array<OneD, NekDouble> tmp2(nq*(nq+1)/2);
376 
377  tri.FwdTrans(nodal, tmp1);
378  nodalTri.ModalToNodal(tmp1, tmp2);
379  newPos[k] = tmp2;
380  }
381 
382  for (l = 0; l < nq*(nq+1)/2; ++l)
383  {
387  dim, face->GetGlobalID(),
388  newPos[0][l], newPos[1][l], newPos[2][l]);
389  curve->m_points.push_back(vert);
390  }
391  }
392  else
393  {
394  for (l = 0; l < nq*nq; ++l)
395  {
399  dim, face->GetGlobalID(),
400  intPos[0][l], intPos[1][l], intPos[2][l]);
401  curve->m_points.push_back(vert);
402  }
403  }
404 
405  curvedFaces[face->GetGlobalID()] = curve;
406  updatedFaces.insert(face->GetGlobalID());
407  }
408  }
409  }
410 
411  // Reset geometry information
412  for (i = 0; i < fields.num_elements(); ++i)
413  {
414  fields[i]->Reset();
415  }
416  }
417 }
418 }
void ModalToNodal(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
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:57
Principle Modified Functions .
Definition: BasisType.h:49
Geometry2DSharedPtr GetFace(int i)
Return face i in this element.
Definition: Geometry3D.cpp:79
PointsType GetPointsType() const
Return type of quadrature.
Definition: Basis.h:151
boost::shared_ptr< Curve > CurveSharedPtr
Definition: Curve.hpp:62
Gauss Radau pinned at x=-1, .
Definition: PointsType.h:57
Principle Orthogonal Functions .
Definition: BasisType.h:47
Class representing a segment element in reference space.
Definition: StdSegExp.h:54
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:116
boost::shared_ptr< StdExpansion2D > StdExpansion2DSharedPtr
int GetNumPoints() const
Return points order at which basis is defined.
Definition: Basis.h:128
boost::shared_ptr< StdExpansion1D > StdExpansion1DSharedPtr
Principle Orthogonal Functions .
Definition: BasisType.h:46
Defines a specification for a set of points.
Definition: Points.h:58
boost::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:68
boost::shared_ptr< Geometry2D > Geometry2DSharedPtr
Definition: Geometry2D.h:59
3D geometry information
Definition: Geometry3D.h:70
2D geometry information
Definition: Geometry2D.h:65
const Geometry1DSharedPtr GetEdge(int i) const
Definition: Geometry2D.cpp:80
boost::shared_ptr< Geometry1D > Geometry1DSharedPtr
Definition: Geometry1D.h:48
boost::unordered_map< int, CurveSharedPtr > CurveMap
Definition: Curve.hpp:63
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191
boost::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:442
Describes the specification for a Basis.
Definition: Basis.h:50
boost::shared_ptr< PointGeom > PointGeomSharedPtr
Definition: Geometry.h:60
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:50
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:285
2D Nodal Electrostatic Points on a Triangle
Definition: PointsType.h:68
boost::shared_ptr< Geometry3D > Geometry3DSharedPtr
Definition: Geometry3D.h:52
void FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function performs the Forward transformation from physical space to coefficient space...