Nektar++
GeomFactors.cpp
Go to the documentation of this file.
1////////////////////////////////////////////////////////////////////////////////
2//
3// File: GeomFactors.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: Geometric factors base class.
32//
33////////////////////////////////////////////////////////////////////////////////
34
37
39{
40/**
41 * @class GeomFactors
42 *
43 * This class stores the various geometric factors associated with a
44 * specific element, necessary for fundamental integration and
45 * differentiation operations as well as edge and surface normals.
46 *
47 * Initially, these algorithms are provided with a mapping from the
48 * reference region element to the physical element. Practically, this
49 * is represented using a corresponding reference region element for
50 * each coordinate component. Note that for straight-sided elements,
51 * these elements will be of linear order. Curved elements are
52 * represented using higher-order coordinate mappings. This geometric
53 * order is in contrast to the order of the spectral/hp expansion order
54 * on the element.
55 *
56 * For application of the chain rule during differentiation we require
57 * the partial derivatives \f[\frac{\partial \xi_i}{\partial \chi_j}\f]
58 * evaluated at the physical points of the expansion basis. We also
59 * construct the inverse metric tensor \f$g^{ij}\f$ which, in the case
60 * of a domain embedded in a higher-dimensional space, supports the
61 * conversion of covariant quantities to contravariant quantities.
62 * When the expansion dimension is equal to the coordinate dimension the
63 * Jacobian of the mapping \f$\chi_j\f$ is a square matrix and
64 * consequently the required terms are the entries of the inverse of the
65 * Jacobian. However, in general this is not the case, so we therefore
66 * implement the construction of these terms following the derivation
67 * in Cantwell, et. al. \cite CaYaKiPeSh13. Given the coordinate maps
68 * \f$\chi_i\f$, this comprises of five steps
69 * -# Compute the terms of the Jacobian
70 * \f$\frac{\partial \chi_i}{\partial \xi_j}\f$.
71 * -# Compute the metric tensor
72 * \f$g_{ij}=\mathbf{t}_{(i)}\cdot\mathbf{t}_{(j)}\f$.
73 * -# Compute the square of the Jacobian determinant
74 * \f$g=|\mathbf{g}|\f$.
75 * -# Compute the inverse metric tensor \f$g^{ij}\f$.
76 * -# Compute the terms \f$\frac{\partial \xi_i}{\partial \chi_j}\f$.
77 */
78
79/**
80 * @param gtype Specified whether the geometry is regular or
81 * deformed.
82 * @param coordim Specifies the dimension of the coordinate
83 * system.
84 * @param Coords Coordinate maps of the element.
85 */
86GeomFactors::GeomFactors(const GeomType gtype, const int coordim,
88 const Array<OneD, Array<OneD, NekDouble>> &coords)
89 : m_type(gtype), m_expDim(xmap->GetShapeDimension()), m_coordDim(coordim),
90 m_valid(true), m_xmap(xmap), m_coords(coords)
91{
93}
94
95/**
96 * @param S An instance of a GeomFactors class from which
97 * to construct a new instance.
98 */
100 : m_type(S.m_type), m_expDim(S.m_expDim), m_coordDim(S.m_coordDim),
101 m_valid(S.m_valid), m_xmap(S.m_xmap), m_coords(S.m_coords)
102{
103}
104
105/**
106 *
107 */
109{
110}
111
112/**
113 * Member data equivalence is tested in the following order: shape type,
114 * expansion dimension, coordinate dimension and coordinates.
115 */
116bool operator==(const GeomFactors &lhs, const GeomFactors &rhs)
117{
118 if (!(lhs.m_type == rhs.m_type))
119 {
120 return false;
121 }
122
123 if (!(lhs.m_expDim == rhs.m_expDim))
124 {
125 return false;
126 }
127
128 if (!(lhs.m_coordDim == rhs.m_coordDim))
129 {
130 return false;
131 }
132
133 const Array<OneD, const NekDouble> jac_lhs =
134 lhs.ComputeJac(lhs.m_xmap->GetPointsKeys());
135 const Array<OneD, const NekDouble> jac_rhs =
136 rhs.ComputeJac(rhs.m_xmap->GetPointsKeys());
137 if (!(jac_lhs == jac_rhs))
138 {
139 return false;
140 }
141
142 return true;
143}
144
145/**
146 * Derivatives are computed at the geometry point distributions and
147 * interpolated to the target point distributions.
148 *
149 * @param tpoints Target point distributions.
150 * @returns Derivative of coordinate map evaluated at
151 * target point distributions.
152 */
154 const LibUtilities::PointsKeyVector &keyTgt) const
155{
156 ASSERTL1(keyTgt.size() == m_expDim,
157 "Dimension of target point distribution does not match "
158 "expansion dimension.");
159
160 int i = 0, j = 0;
161 int nqtot_map = 1;
162 int nqtot_tbasis = 1;
166
167 // Allocate storage and compute number of points
168 for (i = 0; i < m_expDim; ++i)
169 {
170 map_points[i] = m_xmap->GetBasis(i)->GetPointsKey();
171 nqtot_map *= map_points[i].GetNumPoints();
172 nqtot_tbasis *= keyTgt[i].GetNumPoints();
175 }
176
177 // Calculate local derivatives
178 for (i = 0; i < m_coordDim; ++i)
179 {
180 Array<OneD, NekDouble> tmp(nqtot_map);
181 // Transform from coefficient space to physical space
182 m_xmap->BwdTrans(m_coords[i], tmp);
183
184 // Allocate storage and take the derivative (calculated at the
185 // points as specified in 'Coords')
186 for (j = 0; j < m_expDim; ++j)
187 {
188 d_map[j][i] = Array<OneD, NekDouble>(nqtot_map);
189 deriv[j][i] = Array<OneD, NekDouble>(nqtot_tbasis);
190 m_xmap->StdPhysDeriv(j, tmp, d_map[j][i]);
191 }
192 }
193
194 for (i = 0; i < m_coordDim; ++i)
195 {
196 // Interpolate the derivatives:
197 // - from the points as defined in the mapping ('Coords')
198 // - to the points at which we want to know the metrics
199 // ('tbasis')
200 bool same = true;
201 for (j = 0; j < m_expDim; ++j)
202 {
203 same = same && (map_points[j] == keyTgt[j]);
204 }
205 if (same)
206 {
207 for (j = 0; j < m_expDim; ++j)
208 {
209 deriv[j][i] = d_map[j][i];
210 }
211 }
212 else
213 {
214 for (j = 0; j < m_expDim; ++j)
215 {
216 Interp(map_points, d_map[j][i], keyTgt, deriv[j][i]);
217 }
218 }
219 }
220
221 return deriv;
222}
223
224/**
225 * This routine returns an array of values specifying the Jacobian
226 * of the mapping at quadrature points in the element. The array
227 * is either of size 1 in the case of elements having #GeomType
228 * #eRegular, or of size equal to the number of quadrature points for
229 * #eDeformed elements.
230 *
231 * @returns Array containing the Jacobian of the coordinate
232 * mapping at the quadrature points of the element.
233 * @see GeomType
234 */
236 const LibUtilities::PointsKeyVector &keyTgt) const
237{
238 ASSERTL1(keyTgt.size() == m_expDim,
239 "Dimension of target point distribution does not match "
240 "expansion dimension.");
241
242 int i = 0, j = 0, k = 0, l = 0;
243 int ptsTgt = 1;
244
245 if (m_type == eDeformed)
246 {
247 // Allocate storage and compute number of points
248 for (i = 0; i < m_expDim; ++i)
249 {
250 ptsTgt *= keyTgt[i].GetNumPoints();
251 }
252 }
253
254 // Get derivative at geometry points
255 DerivStorage deriv = ComputeDeriv(keyTgt);
256
257 Array<TwoD, NekDouble> tmp(m_expDim * m_expDim, ptsTgt, 0.0);
258 Array<TwoD, NekDouble> gmat(m_expDim * m_expDim, ptsTgt, 0.0);
259 Array<OneD, NekDouble> jac(ptsTgt, 0.0);
260
261 // Compute g_{ij} as t_i \cdot t_j and store in tmp
262 for (i = 0, l = 0; i < m_expDim; ++i)
263 {
264 for (j = 0; j < m_expDim; ++j, ++l)
265 {
266 for (k = 0; k < m_coordDim; ++k)
267 {
268 Vmath::Vvtvp(ptsTgt, &deriv[i][k][0], 1, &deriv[j][k][0], 1,
269 &tmp[l][0], 1, &tmp[l][0], 1);
270 }
271 }
272 }
273
274 Adjoint(tmp, gmat);
275
276 // Compute g = det(g_{ij}) (= Jacobian squared) and store
277 // temporarily in m_jac.
278 for (i = 0; i < m_expDim; ++i)
279 {
280 Vmath::Vvtvp(ptsTgt, &tmp[i][0], 1, &gmat[i * m_expDim][0], 1, &jac[0],
281 1, &jac[0], 1);
282 }
283
284 // Compute the Jacobian = sqrt(g)
285 Vmath::Vsqrt(ptsTgt, &jac[0], 1, &jac[0], 1);
286
287 return jac;
288}
289
290/**
291 * This routine returns a two-dimensional array of values specifying
292 * the inverse metric terms associated with the coordinate mapping of
293 * the corresponding reference region to the physical element. These
294 * terms correspond to the \f$g^{ij}\f$ terms in \cite CaYaKiPeSh13 and,
295 * in the case of an embedded manifold, map covariant quantities to
296 * contravariant quantities. The leading index of the array is the index
297 * of the term in the tensor numbered as
298 * \f[\left(\begin{array}{ccc}
299 * 0 & 1 & 2 \\
300 * 1 & 3 & 4 \\
301 * 2 & 4 & 5
302 * \end{array}\right)\f].
303 * The second dimension is either of size 1 in the case of elements
304 * having #GeomType #eRegular, or of size equal to the number of
305 * quadrature points for #eDeformed elements.
306 *
307 * @see [Wikipedia "Covariance and Contravariance of Vectors"]
308 * (http://en.wikipedia.org/wiki/Covariance_and_contravariance_of_vectors)
309 * @returns Two-dimensional array containing the inverse
310 * metric tensor of the coordinate mapping.
311 */
313 const LibUtilities::PointsKeyVector &keyTgt) const
314{
315 ASSERTL1(keyTgt.size() == m_expDim,
316 "Dimension of target point distribution does not match "
317 "expansion dimension.");
318
319 int i = 0, j = 0, k = 0, l = 0;
320 int ptsTgt = 1;
321
322 if (m_type == eDeformed)
323 {
324 // Allocate storage and compute number of points
325 for (i = 0; i < m_expDim; ++i)
326 {
327 ptsTgt *= keyTgt[i].GetNumPoints();
328 }
329 }
330
331 // Get derivative at geometry points
332 DerivStorage deriv = ComputeDeriv(keyTgt);
333
334 Array<TwoD, NekDouble> tmp(m_expDim * m_expDim, ptsTgt, 0.0);
335 Array<TwoD, NekDouble> gmat(m_expDim * m_expDim, ptsTgt, 0.0);
336 Array<OneD, NekDouble> jac(ptsTgt, 0.0);
337
338 // Compute g_{ij} as t_i \cdot t_j and store in tmp
339 for (i = 0, l = 0; i < m_expDim; ++i)
340 {
341 for (j = 0; j < m_expDim; ++j, ++l)
342 {
343 for (k = 0; k < m_coordDim; ++k)
344 {
345 Vmath::Vvtvp(ptsTgt, &deriv[i][k][0], 1, &deriv[j][k][0], 1,
346 &tmp[l][0], 1, &tmp[l][0], 1);
347 }
348 }
349 }
350
351 Adjoint(tmp, gmat);
352
353 // Compute g = det(g_{ij}) (= Jacobian squared) and store
354 // temporarily in m_jac.
355 for (i = 0; i < m_expDim; ++i)
356 {
357 Vmath::Vvtvp(ptsTgt, &tmp[i][0], 1, &gmat[i * m_expDim][0], 1, &jac[0],
358 1, &jac[0], 1);
359 }
360
361 for (i = 0; i < m_expDim * m_expDim; ++i)
362 {
363 Vmath::Vdiv(ptsTgt, &gmat[i][0], 1, &jac[0], 1, &gmat[i][0], 1);
364 }
365
366 return gmat;
367}
368
369/**
370 * @param keyTgt Target point distributions.
371 * @returns Derivative factors evaluated at the target
372 * point distributions.
373 * A 1D example: /f$ Jac =(\partial x/ \partial \xi) /f$ ;
374 * /f$ factor = 1/Jac = (\partial \xi/ \partial x) /f$
375 */
377 const LibUtilities::PointsKeyVector &keyTgt) const
378{
379 ASSERTL1(keyTgt.size() == m_expDim,
380 "Dimension of target point distribution does not match "
381 "expansion dimension.");
382
383 int i = 0, j = 0, k = 0, l = 0;
384 int ptsTgt = 1;
385
386 if (m_type == eDeformed)
387 {
388 // Allocate storage and compute number of points
389 for (i = 0; i < m_expDim; ++i)
390 {
391 ptsTgt *= keyTgt[i].GetNumPoints();
392 }
393 }
394
395 // Get derivative at geometry points
396 DerivStorage deriv = ComputeDeriv(keyTgt);
397
398 Array<TwoD, NekDouble> tmp(m_expDim * m_expDim, ptsTgt, 0.0);
399 Array<TwoD, NekDouble> gmat(m_expDim * m_expDim, ptsTgt, 0.0);
400 Array<OneD, NekDouble> jac(ptsTgt, 0.0);
402
403 // Compute g_{ij} as t_i \cdot t_j and store in tmp
404 for (i = 0, l = 0; i < m_expDim; ++i)
405 {
406 for (j = 0; j < m_expDim; ++j, ++l)
407 {
408 for (k = 0; k < m_coordDim; ++k)
409 {
410 Vmath::Vvtvp(ptsTgt, &deriv[i][k][0], 1, &deriv[j][k][0], 1,
411 &tmp[l][0], 1, &tmp[l][0], 1);
412 }
413 }
414 }
415
416 Adjoint(tmp, gmat);
417
418 // Compute g = det(g_{ij}) (= Jacobian squared) and store
419 // temporarily in m_jac.
420 for (i = 0; i < m_expDim; ++i)
421 {
422 Vmath::Vvtvp(ptsTgt, &tmp[i][0], 1, &gmat[i * m_expDim][0], 1, &jac[0],
423 1, &jac[0], 1);
424 }
425
426 for (i = 0; i < m_expDim * m_expDim; ++i)
427 {
428 Vmath::Vdiv(ptsTgt, &gmat[i][0], 1, &jac[0], 1, &gmat[i][0], 1);
429 }
430
431 // Compute the Jacobian = sqrt(g)
432 Vmath::Vsqrt(ptsTgt, &jac[0], 1, &jac[0], 1);
433
434 // Compute the derivative factors
435 for (k = 0, l = 0; k < m_coordDim; ++k)
436 {
437 for (j = 0; j < m_expDim; ++j, ++l)
438 {
439 for (i = 0; i < m_expDim; ++i)
440 {
441 Vmath::Vvtvp(ptsTgt, &deriv[i][k][0], 1,
442 &gmat[m_expDim * i + j][0], 1, &factors[l][0], 1,
443 &factors[l][0], 1);
444 }
445 }
446 }
447
448 return factors;
449}
450
452 const LibUtilities::PointsKeyVector &keyTgt,
453 const SpatialDomains::GeomMMF MMFdir,
455 Array<OneD, Array<OneD, NekDouble>> &movingframes)
456{
457 ASSERTL1(keyTgt.size() == m_expDim,
458 "Dimension of target point distribution does not match "
459 "expansion dimension.");
460
461 int i = 0, k = 0;
462 int ptsTgt = 1;
463 int nq = 1;
464
465 for (i = 0; i < m_expDim; ++i)
466 {
467 nq *= keyTgt[i].GetNumPoints();
468 }
469
470 if (m_type == eDeformed)
471 {
472 // Allocate storage and compute number of points
473 for (i = 0; i < m_expDim; ++i)
474 {
475 ptsTgt *= keyTgt[i].GetNumPoints();
476 }
477 }
478
479 // Get derivative at geometry points
480 DerivStorage deriv = ComputeDeriv(keyTgt);
481
482 // number of moving frames is requited to be 3, even for surfaces
483 int MFdim = 3;
484
486
487 // Compute g_{ij} as t_i \cdot t_j and store in tmp
488 for (i = 0; i < MFdim; ++i)
489 {
491 for (k = 0; k < m_coordDim; ++k)
492 {
493 MFtmp[i][k] = Array<OneD, NekDouble>(nq);
494 }
495 }
496
497 // Compute g_{ij} as t_i \cdot t_j and store in tmp
498 for (i = 0; i < MFdim - 1; ++i)
499 {
500 for (k = 0; k < m_coordDim; ++k)
501 {
502 if (m_type == eDeformed)
503 {
504 Vmath::Vcopy(ptsTgt, &deriv[i][k][0], 1, &MFtmp[i][k][0], 1);
505 }
506 else
507 {
508 Vmath::Fill(nq, deriv[i][k][0], MFtmp[i][k], 1);
509 }
510 }
511 }
512
513 // Direction of MF1 is preserved: MF2 is considered in the same
514 // tangent plane as MF1. MF3 is computed by cross product of MF1
515 // and MF2. MF2 is consequently computed as the cross product of
516 // MF3 and MF1.
518 for (int k = 0; k < m_coordDim; k++)
519 {
520 PrincipleDir[k] = Array<OneD, NekDouble>(nq);
521 }
522
523 if (!(MMFdir == eLOCAL))
524 {
525 ComputePrincipleDirection(keyTgt, MMFdir, factors, PrincipleDir);
526 }
527
528 // MF3 = MF1 \times MF2
529 VectorCrossProd(MFtmp[0], MFtmp[1], MFtmp[2]);
530
531 // Normalizing MF3
532 VectorNormalise(MFtmp[2]);
533
534 if (!(MMFdir == eLOCAL))
535 {
536 Array<OneD, NekDouble> temp(nq, 0.0);
537
538 // Reorient MF1 along the PrincipleDir
539 for (i = 0; i < m_coordDim; ++i)
540 {
541 Vmath::Vvtvp(nq, MFtmp[2][i], 1, PrincipleDir[i], 1, temp, 1, temp,
542 1);
543 }
544 Vmath::Neg(nq, temp, 1);
545
546 // u2 = v2 - < u1 , v2 > ( u1 / < u1, u1 > )
547 for (i = 0; i < m_coordDim; ++i)
548 {
549 Vmath::Vvtvp(nq, temp, 1, MFtmp[2][i], 1, PrincipleDir[i], 1,
550 MFtmp[0][i], 1);
551 }
552 }
553
554 // Normalizing MF1
555 VectorNormalise(MFtmp[0]);
556
557 // MF2 = MF3 \times MF1
558 VectorCrossProd(MFtmp[2], MFtmp[0], MFtmp[1]);
559
560 // Normalizing MF2
561 VectorNormalise(MFtmp[1]);
562
563 for (i = 0; i < MFdim; ++i)
564 {
565 for (k = 0; k < m_coordDim; ++k)
566 {
567 Vmath::Vcopy(nq, &MFtmp[i][k][0], 1,
568 &movingframes[i * m_coordDim + k][0], 1);
569 }
570 }
571}
572
573/**
574 * Constructs the Jacobian as per Spencer's book p158 and tests if
575 * negative.
576 */
578{
579 // Jacobian test only makes sense when expdim = coorddim
580 // If one-dimensional then element is valid.
581 if (m_coordDim != m_expDim || m_expDim == 1)
582 {
583 m_valid = true;
584 return;
585 }
586
588 int nqtot = 1;
589 for (int i = 0; i < m_expDim; ++i)
590 {
591 p[i] = m_xmap->GetBasis(i)->GetPointsKey();
592 nqtot *= p[i].GetNumPoints();
593 }
594 int pts = (m_type == eRegular || m_type == eMovingRegular) ? 1 : nqtot;
595 Array<OneD, NekDouble> jac(pts, 0.0);
596
597 DerivStorage deriv = GetDeriv(p);
598
599 switch (m_expDim)
600 {
601 case 2:
602 {
603 Vmath::Vvtvvtm(pts, &deriv[0][0][0], 1, &deriv[1][1][0], 1,
604 &deriv[1][0][0], 1, &deriv[0][1][0], 1, &jac[0], 1);
605 break;
606 }
607 case 3:
608 {
609 Array<OneD, NekDouble> tmp(pts, 0.0);
610
611 Vmath::Vvtvvtm(pts, &deriv[1][1][0], 1, &deriv[2][2][0], 1,
612 &deriv[2][1][0], 1, &deriv[1][2][0], 1, &tmp[0], 1);
613 Vmath::Vvtvp(pts, &deriv[0][0][0], 1, &tmp[0], 1, &jac[0], 1,
614 &jac[0], 1);
615
616 Vmath::Vvtvvtm(pts, &deriv[2][1][0], 1, &deriv[0][2][0], 1,
617 &deriv[0][1][0], 1, &deriv[2][2][0], 1, &tmp[0], 1);
618 Vmath::Vvtvp(pts, &deriv[1][0][0], 1, &tmp[0], 1, &jac[0], 1,
619 &jac[0], 1);
620
621 Vmath::Vvtvvtm(pts, &deriv[0][1][0], 1, &deriv[1][2][0], 1,
622 &deriv[1][1][0], 1, &deriv[0][2][0], 1, &tmp[0], 1);
623 Vmath::Vvtvp(pts, &deriv[2][0][0], 1, &tmp[0], 1, &jac[0], 1,
624 &jac[0], 1);
625
626 break;
627 }
628 }
629
630 if (Vmath::Vmin(pts, &jac[0], 1) < 0)
631 {
632 m_valid = false;
633 }
634}
635
636/**
637 * @param map_points Source data point distribution.
638 * @param src Source data to be interpolated.
639 * @param tpoints Target data point distribution.
640 * @param tgt Target data storage.
641 */
644 const LibUtilities::PointsKeyVector &tgt_points,
645 Array<OneD, NekDouble> &tgt) const
646{
647 ASSERTL1(src_points.size() == tgt_points.size(),
648 "Dimension of target point distribution does not match "
649 "expansion dimension.");
650
651 switch (m_expDim)
652 {
653 case 1:
654 LibUtilities::Interp1D(src_points[0], src, tgt_points[0], tgt);
655 break;
656 case 2:
657 LibUtilities::Interp2D(src_points[0], src_points[1], src,
658 tgt_points[0], tgt_points[1], tgt);
659 break;
660 case 3:
661 LibUtilities::Interp3D(src_points[0], src_points[1], src_points[2],
662 src, tgt_points[0], tgt_points[1],
663 tgt_points[2], tgt);
664 break;
665 }
666}
667
668/**
669 * Input and output arrays are of dimension
670 * (m_expDim*m_expDim) x num_points. The first index of the input and
671 * output arrays are ordered row-by-row.
672 * @param src Input data array.
673 * @param tgt Storage for adjoint matrix data.
674 */
676 Array<TwoD, NekDouble> &tgt) const
677{
678 ASSERTL1(src.size() == tgt.size(),
679 "Source matrix is of different size to destination"
680 "matrix for computing adjoint.");
681
682 int n = src[0].size();
683 switch (m_expDim)
684 {
685 case 1:
686 Vmath::Fill(n, 1.0, &tgt[0][0], 1);
687 break;
688 case 2:
689 Vmath::Vcopy(n, &src[3][0], 1, &tgt[0][0], 1);
690 Vmath::Smul(n, -1.0, &src[1][0], 1, &tgt[1][0], 1);
691 Vmath::Smul(n, -1.0, &src[2][0], 1, &tgt[2][0], 1);
692 Vmath::Vcopy(n, &src[0][0], 1, &tgt[3][0], 1);
693 break;
694 case 3:
695 {
696 int a, b, c, d, e, i, j;
697
698 // Compute g^{ij} by computing Cofactors(g_ij)^T
699 for (i = 0; i < m_expDim; ++i)
700 {
701 for (j = 0; j < m_expDim; ++j)
702 {
703 a = ((i + 1) % m_expDim) * m_expDim + ((j + 1) % m_expDim);
704 b = ((i + 1) % m_expDim) * m_expDim + ((j + 2) % m_expDim);
705 c = ((i + 2) % m_expDim) * m_expDim + ((j + 1) % m_expDim);
706 d = ((i + 2) % m_expDim) * m_expDim + ((j + 2) % m_expDim);
707 e = j * m_expDim + i;
708 Vmath::Vvtvvtm(n, &src[a][0], 1, &src[d][0], 1, &src[b][0],
709 1, &src[c][0], 1, &tgt[e][0], 1);
710 }
711 }
712 break;
713 }
714 }
715}
716
717/**
718 *
719 */
721 const LibUtilities::PointsKeyVector &keyTgt,
722 const SpatialDomains::GeomMMF MMFdir,
725{
726 int nq = output[0].size();
727
729 for (int i = 0; i < m_coordDim; ++i)
730 {
731 output[i] = Array<OneD, NekDouble>(nq, 0.0);
732 }
733
734 // Construction of Connection
735 switch (MMFdir)
736 {
737 // projection to x-axis
738 case eTangentX:
739 {
740 Vmath::Fill(nq, 1.0, output[0], 1);
741 break;
742 }
743 case eTangentY:
744 {
745 Vmath::Fill(nq, 1.0, output[1], 1);
746 break;
747 }
748 case eTangentXY:
749 {
750 Vmath::Fill(nq, sqrt(2.0), output[0], 1);
751 Vmath::Fill(nq, sqrt(2.0), output[1], 1);
752 break;
753 }
754 case eTangentZ:
755 {
756 Vmath::Fill(nq, 1.0, output[2], 1);
757 break;
758 }
759 case eTangentCircular:
760 {
761 // Tangent direction depends on spatial location.
763 for (int k = 0; k < m_coordDim; k++)
764 {
765 x[k] = Array<OneD, NekDouble>(nq);
766 }
767
768 // m_coords are StdExpansions which store the mapping
769 // between the std element and the local element. Bwd
770 // transforming the std element minimum basis gives a
771 // minimum physical basis for geometry. Need to then
772 // interpolate this up to the quadrature basis.
773 int nqtot_map = 1;
775 for (int i = 0; i < m_expDim; ++i)
776 {
777 map_points[i] = m_xmap->GetBasis(i)->GetPointsKey();
778 nqtot_map *= map_points[i].GetNumPoints();
779 }
780 Array<OneD, NekDouble> tmp(nqtot_map);
781 for (int k = 0; k < m_coordDim; k++)
782 {
783 m_xmap->BwdTrans(m_coords[k], tmp);
784 Interp(map_points, tmp, keyTgt, x[k]);
785 }
786
787 // circular around the center of the domain
788 NekDouble radius, xc = 0.0, yc = 0.0, xdis, ydis;
789 NekDouble la, lb;
790
791 ASSERTL1(factors.size() >= 4, "factors is too short.");
792
793 la = factors[0];
794 lb = factors[1];
795 xc = factors[2];
796 yc = factors[3];
797
798 for (int i = 0; i < nq; i++)
799 {
800 xdis = x[0][i] - xc;
801 ydis = x[1][i] - yc;
802 radius = sqrt(xdis * xdis / la / la + ydis * ydis / lb / lb);
803 output[0][i] = ydis / radius;
804 output[1][i] = -1.0 * xdis / radius;
805 }
806 break;
807 }
809 {
810 // Tangent direction depends on spatial location.
812 for (int k = 0; k < m_coordDim; k++)
813 {
814 x[k] = Array<OneD, NekDouble>(nq);
815 }
816
817 int nqtot_map = 1;
819 for (int i = 0; i < m_expDim; ++i)
820 {
821 map_points[i] = m_xmap->GetBasis(i)->GetPointsKey();
822 nqtot_map *= map_points[i].GetNumPoints();
823 }
824 Array<OneD, NekDouble> tmp(nqtot_map);
825 for (int k = 0; k < m_coordDim; k++)
826 {
827 m_xmap->BwdTrans(m_coords[k], tmp);
828 Interp(map_points, tmp, keyTgt, x[k]);
829 }
830
831 // circular around the center of the domain
832 NekDouble xtan, ytan, mag;
833 for (int i = 0; i < nq; i++)
834 {
835 xtan = -1.0 * (x[1][i] * x[1][i] * x[1][i] + x[1][i]);
836 ytan = 2.0 * x[0][i];
837 mag = sqrt(xtan * xtan + ytan * ytan);
838 output[0][i] = xtan / mag;
839 output[1][i] = ytan / mag;
840 }
841 break;
842 }
844 {
845 // Tangent direction depends on spatial location.
847 for (int k = 0; k < m_coordDim; k++)
848 {
849 x[k] = Array<OneD, NekDouble>(nq);
850 }
851
852 int nqtot_map = 1;
854 for (int i = 0; i < m_expDim; ++i)
855 {
856 map_points[i] = m_xmap->GetBasis(i)->GetPointsKey();
857 nqtot_map *= map_points[i].GetNumPoints();
858 }
859 Array<OneD, NekDouble> tmp(nqtot_map);
860 for (int k = 0; k < m_coordDim; k++)
861 {
862 m_xmap->BwdTrans(m_coords[k], tmp);
863 Interp(map_points, tmp, keyTgt, x[k]);
864 }
865
866 // circular around the center of the domain
867 NekDouble xtan, ytan, mag;
868 for (int i = 0; i < nq; i++)
869 {
870 xtan = -2.0 * x[1][i] * x[1][i] * x[1][i] + x[1][i];
871 ytan = sqrt(3.0) * x[0][i];
872 mag = sqrt(xtan * xtan + ytan * ytan);
873 output[0][i] = xtan / mag;
874 output[1][i] = ytan / mag;
875 }
876 break;
877 }
878 default:
879 {
880 break;
881 }
882 }
883}
884
885/**
886 *
887 */
889{
890 int ndim = array.size();
891 ASSERTL0(ndim > 0, "Number of components must be > 0.");
892 for (int i = 1; i < ndim; ++i)
893 {
894 ASSERTL0(array[i].size() == array[0].size(),
895 "Array size mismatch in coordinates.");
896 }
897
898 int nq = array[0].size();
899 Array<OneD, NekDouble> norm(nq, 0.0);
900
901 // Compute the norm of each vector.
902 for (int i = 0; i < ndim; ++i)
903 {
904 Vmath::Vvtvp(nq, array[i], 1, array[i], 1, norm, 1, norm, 1);
905 }
906
907 Vmath::Vsqrt(nq, norm, 1, norm, 1);
908
909 // Normalise the vectors by the norm
910 for (int i = 0; i < ndim; ++i)
911 {
912 Vmath::Vdiv(nq, array[i], 1, norm, 1, array[i], 1);
913 }
914}
915
916/**
917 * @brief Computes the vector cross-product in 3D of \a v1 and \a v2,
918 * storing the result in \a v3.
919 *
920 * @param v1 First input vector.
921 * @param v2 Second input vector.
922 * @param v3 Output vector computed to be orthogonal to
923 * both \a v1 and \a v2.
924 */
926 const Array<OneD, const Array<OneD, NekDouble>> &v1,
927 const Array<OneD, const Array<OneD, NekDouble>> &v2,
929{
930 ASSERTL0(v1.size() == 3, "Input 1 has dimension not equal to 3.");
931 ASSERTL0(v2.size() == 3, "Input 2 has dimension not equal to 3.");
932 ASSERTL0(v3.size() == 3, "Output vector has dimension not equal to 3.");
933
934 int nq = v1[0].size();
935 Array<OneD, NekDouble> temp(nq);
936
937 Vmath::Vmul(nq, v1[2], 1, v2[1], 1, temp, 1);
938 Vmath::Vvtvm(nq, v1[1], 1, v2[2], 1, temp, 1, v3[0], 1);
939
940 Vmath::Vmul(nq, v1[0], 1, v2[2], 1, temp, 1);
941 Vmath::Vvtvm(nq, v1[2], 1, v2[0], 1, temp, 1, v3[1], 1);
942
943 Vmath::Vmul(nq, v1[1], 1, v2[0], 1, temp, 1);
944 Vmath::Vvtvm(nq, v1[0], 1, v2[1], 1, temp, 1, v3[2], 1);
945}
946
947} // namespace Nektar::SpatialDomains
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
Calculation and storage of geometric factors associated with the mapping from StdRegions reference el...
Definition: GeomFactors.h:72
DerivStorage ComputeDeriv(const LibUtilities::PointsKeyVector &keyTgt) const
Array< TwoD, NekDouble > ComputeDerivFactors(const LibUtilities::PointsKeyVector &keyTgt) const
Return the derivative of the reference coordinates with respect to the mapping, .
void Adjoint(const Array< TwoD, const NekDouble > &src, Array< TwoD, NekDouble > &tgt) const
Compute the transpose of the cofactors matrix.
int m_coordDim
Dimension of coordinate system.
Definition: GeomFactors.h:131
void CheckIfValid()
Tests if the element is valid and not self-intersecting.
void Interp(const LibUtilities::PointsKeyVector &src_points, const Array< OneD, const NekDouble > &src, const LibUtilities::PointsKeyVector &tgt_points, Array< OneD, NekDouble > &tgt) const
Perform interpolation of data between two point distributions.
Array< OneD, NekDouble > ComputeJac(const LibUtilities::PointsKeyVector &keyTgt) const
Return the Jacobian of the mapping and cache the result.
int m_expDim
Dimension of expansion.
Definition: GeomFactors.h:129
Array< TwoD, NekDouble > ComputeGmat(const LibUtilities::PointsKeyVector &keyTgt) const
Computes the Laplacian coefficients .
void ComputePrincipleDirection(const LibUtilities::PointsKeyVector &keyTgt, const SpatialDomains::GeomMMF MMFdir, const Array< OneD, const NekDouble > &CircCentre, Array< OneD, Array< OneD, NekDouble > > &output)
bool m_valid
Validity of element (Jacobian positive)
Definition: GeomFactors.h:133
StdRegions::StdExpansionSharedPtr m_xmap
Stores information about the expansion.
Definition: GeomFactors.h:139
GeomFactors(const GeomType gtype, const int coordim, const StdRegions::StdExpansionSharedPtr &xmap, const Array< OneD, Array< OneD, NekDouble > > &coords)
Constructor for GeomFactors class.
Definition: GeomFactors.cpp:86
DerivStorage GetDeriv(const LibUtilities::PointsKeyVector &keyTgt)
Return the derivative of the mapping with respect to the reference coordinates, .
Definition: GeomFactors.h:216
GeomType m_type
Type of geometry (e.g. eRegular, eDeformed, eMovingRegular).
Definition: GeomFactors.h:127
Array< OneD, Array< OneD, NekDouble > > m_coords
Stores coordinates of the geometry.
Definition: GeomFactors.h:141
void ComputeMovingFrames(const LibUtilities::PointsKeyVector &keyTgt, const SpatialDomains::GeomMMF MMFdir, const Array< OneD, const NekDouble > &CircCentre, Array< OneD, Array< OneD, NekDouble > > &movingframes)
void VectorNormalise(Array< OneD, Array< OneD, NekDouble > > &array)
void VectorCrossProd(const Array< OneD, const Array< OneD, NekDouble > > &v1, const Array< OneD, const Array< OneD, NekDouble > > &v2, Array< OneD, Array< OneD, NekDouble > > &v3)
Computes the vector cross-product in 3D of v1 and v2, storing the result in v3.
void Interp1D(const BasisKey &fbasis0, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, Array< OneD, NekDouble > &to)
this function interpolates a 1D function evaluated at the quadrature points of the basis fbasis0 to ...
Definition: Interp.cpp:47
void Interp3D(const BasisKey &fbasis0, const BasisKey &fbasis1, const BasisKey &fbasis2, const Array< OneD, const NekDouble > &from, const BasisKey &tbasis0, const BasisKey &tbasis1, const BasisKey &tbasis2, Array< OneD, NekDouble > &to)
this function interpolates a 3D function evaluated at the quadrature points of the 3D basis,...
Definition: Interp.cpp:162
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:101
std::vector< PointsKey > PointsKeyVector
Definition: Points.h:231
GeomMMF
Principle direction for MMF.
@ eLOCAL
No Principal direction.
@ eTangentIrregular
Circular around the centre of domain.
@ eTangentX
X coordinate direction.
@ eTangentCircular
Circular around the centre of domain.
@ eTangentNonconvex
Circular around the centre of domain.
@ eTangentXY
XY direction.
@ eTangentZ
Z coordinate direction.
@ eTangentY
Y coordinate direction.
bool operator==(const GeomFactors &lhs, const GeomFactors &rhs)
Equivalence test for GeomFactors objects.
GeomType
Indicates the type of element geometry.
@ eRegular
Geometry is straight-sided with constant geometric factors.
@ eMovingRegular
Currently unused.
@ eDeformed
Geometry is curved or has non-constant factors.
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > DerivStorage
Storage type for derivative of mapping.
Definition: GeomFactors.h:66
std::shared_ptr< StdExpansion > StdExpansionSharedPtr
std::vector< double > d(NPUPPER *NPUPPER)
StdRegions::ConstFactorMap factors
double NekDouble
void Vsqrt(int n, const T *x, const int incx, T *y, const int incy)
sqrt y = sqrt(x)
Definition: Vmath.hpp:340
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.hpp:72
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.hpp:292
T Vmin(int n, const T *x, const int incx)
Return the minimum element in x - called vmin to avoid conflict with min.
Definition: Vmath.hpp:725
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.hpp:366
void Vvtvm(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvm (vector times vector minus vector): z = w*x - y
Definition: Vmath.hpp:381
void Vvtvvtm(int n, const T *v, int incv, const T *w, int incw, const T *x, int incx, const T *y, int incy, T *z, int incz)
vvtvvtm (vector times vector minus vector times vector):
Definition: Vmath.hpp:456
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.hpp:100
void Vdiv(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x/y.
Definition: Vmath.hpp:126
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
Definition: Vmath.hpp:54
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294