Nektar++
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 * Member data equivalence is tested in the following order: shape type,
107 * expansion dimension, coordinate dimension and coordinates.
108 */
109bool operator==(const GeomFactors &lhs, const GeomFactors &rhs)
110{
111 if (!(lhs.m_type == rhs.m_type))
112 {
113 return false;
114 }
115
116 if (!(lhs.m_expDim == rhs.m_expDim))
117 {
118 return false;
119 }
120
121 if (!(lhs.m_coordDim == rhs.m_coordDim))
122 {
123 return false;
124 }
125
126 const Array<OneD, const NekDouble> jac_lhs =
127 lhs.ComputeJac(lhs.m_xmap->GetPointsKeys());
128 const Array<OneD, const NekDouble> jac_rhs =
129 rhs.ComputeJac(rhs.m_xmap->GetPointsKeys());
130 if (!(jac_lhs == jac_rhs))
131 {
132 return false;
133 }
134
135 return true;
136}
137
138/**
139 * Derivatives are computed at the geometry point distributions and
140 * interpolated to the target point distributions.
141 *
142 * @param tpoints Target point distributions.
143 * @returns Derivative of coordinate map evaluated at
144 * target point distributions.
145 */
147 const LibUtilities::PointsKeyVector &keyTgt) const
148{
149 ASSERTL1(keyTgt.size() == m_expDim,
150 "Dimension of target point distribution does not match "
151 "expansion dimension.");
152
153 int i = 0, j = 0;
154 int nqtot_map = 1;
155 int nqtot_tbasis = 1;
159
160 // Allocate storage and compute number of points
161 for (i = 0; i < m_expDim; ++i)
162 {
163 map_points[i] = m_xmap->GetBasis(i)->GetPointsKey();
164 nqtot_map *= map_points[i].GetNumPoints();
165 nqtot_tbasis *= keyTgt[i].GetNumPoints();
168 }
169
170 // Calculate local derivatives
171 for (i = 0; i < m_coordDim; ++i)
172 {
173 Array<OneD, NekDouble> tmp(nqtot_map);
174 // Transform from coefficient space to physical space
175 m_xmap->BwdTrans(m_coords[i], tmp);
176
177 // Allocate storage and take the derivative (calculated at the
178 // points as specified in 'Coords')
179 for (j = 0; j < m_expDim; ++j)
180 {
181 d_map[j][i] = Array<OneD, NekDouble>(nqtot_map);
182 deriv[j][i] = Array<OneD, NekDouble>(nqtot_tbasis);
183 m_xmap->StdPhysDeriv(j, tmp, d_map[j][i]);
184 }
185 }
186
187 for (i = 0; i < m_coordDim; ++i)
188 {
189 // Interpolate the derivatives:
190 // - from the points as defined in the mapping ('Coords')
191 // - to the points at which we want to know the metrics
192 // ('tbasis')
193 bool same = true;
194 for (j = 0; j < m_expDim; ++j)
195 {
196 same = same && (map_points[j] == keyTgt[j]);
197 }
198 if (same)
199 {
200 for (j = 0; j < m_expDim; ++j)
201 {
202 deriv[j][i] = d_map[j][i];
203 }
204 }
205 else
206 {
207 for (j = 0; j < m_expDim; ++j)
208 {
209 Interp(map_points, d_map[j][i], keyTgt, deriv[j][i]);
210 }
211 }
212 }
213
214 return deriv;
215}
216
217/**
218 * This routine returns an array of values specifying the Jacobian
219 * of the mapping at quadrature points in the element. The array
220 * is either of size 1 in the case of elements having #GeomType
221 * #eRegular, or of size equal to the number of quadrature points for
222 * #eDeformed elements.
223 *
224 * @returns Array containing the Jacobian of the coordinate
225 * mapping at the quadrature points of the element.
226 * @see GeomType
227 */
229 const LibUtilities::PointsKeyVector &keyTgt) const
230{
231 ASSERTL1(keyTgt.size() == m_expDim,
232 "Dimension of target point distribution does not match "
233 "expansion dimension.");
234
235 int i = 0, j = 0, k = 0, l = 0;
236 int ptsTgt = 1;
237
238 if (m_type == eDeformed)
239 {
240 // Allocate storage and compute number of points
241 for (i = 0; i < m_expDim; ++i)
242 {
243 ptsTgt *= keyTgt[i].GetNumPoints();
244 }
245 }
246
247 // Get derivative at geometry points
248 DerivStorage deriv = ComputeDeriv(keyTgt);
249
250 Array<TwoD, NekDouble> tmp(m_expDim * m_expDim, ptsTgt, 0.0);
251 Array<TwoD, NekDouble> gmat(m_expDim * m_expDim, ptsTgt, 0.0);
252 Array<OneD, NekDouble> jac(ptsTgt, 0.0);
253
254 // Compute g_{ij} as t_i \cdot t_j and store in tmp
255 for (i = 0, l = 0; i < m_expDim; ++i)
256 {
257 for (j = 0; j < m_expDim; ++j, ++l)
258 {
259 for (k = 0; k < m_coordDim; ++k)
260 {
261 Vmath::Vvtvp(ptsTgt, &deriv[i][k][0], 1, &deriv[j][k][0], 1,
262 &tmp[l][0], 1, &tmp[l][0], 1);
263 }
264 }
265 }
266
267 Adjoint(tmp, gmat);
268
269 // Compute g = det(g_{ij}) (= Jacobian squared) and store
270 // temporarily in m_jac.
271 for (i = 0; i < m_expDim; ++i)
272 {
273 Vmath::Vvtvp(ptsTgt, &tmp[i][0], 1, &gmat[i * m_expDim][0], 1, &jac[0],
274 1, &jac[0], 1);
275 }
276
277 // Compute the Jacobian = sqrt(g)
278 Vmath::Vsqrt(ptsTgt, &jac[0], 1, &jac[0], 1);
279
280 return jac;
281}
282
283/**
284 * This routine returns a two-dimensional array of values specifying
285 * the inverse metric terms associated with the coordinate mapping of
286 * the corresponding reference region to the physical element. These
287 * terms correspond to the \f$g^{ij}\f$ terms in \cite CaYaKiPeSh13 and,
288 * in the case of an embedded manifold, map covariant quantities to
289 * contravariant quantities. The leading index of the array is the index
290 * of the term in the tensor numbered as
291 * \f[\left(\begin{array}{ccc}
292 * 0 & 1 & 2 \\
293 * 1 & 3 & 4 \\
294 * 2 & 4 & 5
295 * \end{array}\right)\f].
296 * The second dimension is either of size 1 in the case of elements
297 * having #GeomType #eRegular, or of size equal to the number of
298 * quadrature points for #eDeformed elements.
299 *
300 * @see [Wikipedia "Covariance and Contravariance of Vectors"]
301 * (http://en.wikipedia.org/wiki/Covariance_and_contravariance_of_vectors)
302 * @returns Two-dimensional array containing the inverse
303 * metric tensor of the coordinate mapping.
304 */
306 const LibUtilities::PointsKeyVector &keyTgt) const
307{
308 ASSERTL1(keyTgt.size() == m_expDim,
309 "Dimension of target point distribution does not match "
310 "expansion dimension.");
311
312 int i = 0, j = 0, k = 0, l = 0;
313 int ptsTgt = 1;
314
315 if (m_type == eDeformed)
316 {
317 // Allocate storage and compute number of points
318 for (i = 0; i < m_expDim; ++i)
319 {
320 ptsTgt *= keyTgt[i].GetNumPoints();
321 }
322 }
323
324 // Get derivative at geometry points
325 DerivStorage deriv = ComputeDeriv(keyTgt);
326
327 Array<TwoD, NekDouble> tmp(m_expDim * m_expDim, ptsTgt, 0.0);
328 Array<TwoD, NekDouble> gmat(m_expDim * m_expDim, ptsTgt, 0.0);
329 Array<OneD, NekDouble> jac(ptsTgt, 0.0);
330
331 // Compute g_{ij} as t_i \cdot t_j and store in tmp
332 for (i = 0, l = 0; i < m_expDim; ++i)
333 {
334 for (j = 0; j < m_expDim; ++j, ++l)
335 {
336 for (k = 0; k < m_coordDim; ++k)
337 {
338 Vmath::Vvtvp(ptsTgt, &deriv[i][k][0], 1, &deriv[j][k][0], 1,
339 &tmp[l][0], 1, &tmp[l][0], 1);
340 }
341 }
342 }
343
344 Adjoint(tmp, gmat);
345
346 // Compute g = det(g_{ij}) (= Jacobian squared) and store
347 // temporarily in m_jac.
348 for (i = 0; i < m_expDim; ++i)
349 {
350 Vmath::Vvtvp(ptsTgt, &tmp[i][0], 1, &gmat[i * m_expDim][0], 1, &jac[0],
351 1, &jac[0], 1);
352 }
353
354 for (i = 0; i < m_expDim * m_expDim; ++i)
355 {
356 Vmath::Vdiv(ptsTgt, &gmat[i][0], 1, &jac[0], 1, &gmat[i][0], 1);
357 }
358
359 return gmat;
360}
361
362/**
363 * @param keyTgt Target point distributions.
364 * @returns Derivative factors evaluated at the target
365 * point distributions.
366 * A 1D example: /f$ Jac =(\partial x/ \partial \xi) /f$ ;
367 * /f$ factor = 1/Jac = (\partial \xi/ \partial x) /f$
368 */
370 const LibUtilities::PointsKeyVector &keyTgt) const
371{
372 ASSERTL1(keyTgt.size() == m_expDim,
373 "Dimension of target point distribution does not match "
374 "expansion dimension.");
375
376 int i = 0, j = 0, k = 0, l = 0;
377 int ptsTgt = 1;
378
379 if (m_type == eDeformed)
380 {
381 // Allocate storage and compute number of points
382 for (i = 0; i < m_expDim; ++i)
383 {
384 ptsTgt *= keyTgt[i].GetNumPoints();
385 }
386 }
387
388 // Get derivative at geometry points
389 DerivStorage deriv = ComputeDeriv(keyTgt);
390
391 Array<TwoD, NekDouble> tmp(m_expDim * m_expDim, ptsTgt, 0.0);
392 Array<TwoD, NekDouble> gmat(m_expDim * m_expDim, ptsTgt, 0.0);
393 Array<OneD, NekDouble> jac(ptsTgt, 0.0);
395
396 // Compute g_{ij} as t_i \cdot t_j and store in tmp
397 for (i = 0, l = 0; i < m_expDim; ++i)
398 {
399 for (j = 0; j < m_expDim; ++j, ++l)
400 {
401 for (k = 0; k < m_coordDim; ++k)
402 {
403 Vmath::Vvtvp(ptsTgt, &deriv[i][k][0], 1, &deriv[j][k][0], 1,
404 &tmp[l][0], 1, &tmp[l][0], 1);
405 }
406 }
407 }
408
409 Adjoint(tmp, gmat);
410
411 // Compute g = det(g_{ij}) (= Jacobian squared) and store
412 // temporarily in m_jac.
413 for (i = 0; i < m_expDim; ++i)
414 {
415 Vmath::Vvtvp(ptsTgt, &tmp[i][0], 1, &gmat[i * m_expDim][0], 1, &jac[0],
416 1, &jac[0], 1);
417 }
418
419 for (i = 0; i < m_expDim * m_expDim; ++i)
420 {
421 Vmath::Vdiv(ptsTgt, &gmat[i][0], 1, &jac[0], 1, &gmat[i][0], 1);
422 }
423
424 // Compute the Jacobian = sqrt(g)
425 Vmath::Vsqrt(ptsTgt, &jac[0], 1, &jac[0], 1);
426
427 // Compute the derivative factors
428 for (k = 0, l = 0; k < m_coordDim; ++k)
429 {
430 for (j = 0; j < m_expDim; ++j, ++l)
431 {
432 for (i = 0; i < m_expDim; ++i)
433 {
434 Vmath::Vvtvp(ptsTgt, &deriv[i][k][0], 1,
435 &gmat[m_expDim * i + j][0], 1, &factors[l][0], 1,
436 &factors[l][0], 1);
437 }
438 }
439 }
440
441 return factors;
442}
443
445 const LibUtilities::PointsKeyVector &keyTgt,
446 const SpatialDomains::GeomMMF MMFdir,
448 Array<OneD, Array<OneD, NekDouble>> &movingframes)
449{
450 ASSERTL1(keyTgt.size() == m_expDim,
451 "Dimension of target point distribution does not match "
452 "expansion dimension.");
453
454 int i = 0, k = 0;
455 int ptsTgt = 1;
456 int nq = 1;
457
458 for (i = 0; i < m_expDim; ++i)
459 {
460 nq *= keyTgt[i].GetNumPoints();
461 }
462
463 if (m_type == eDeformed)
464 {
465 // Allocate storage and compute number of points
466 for (i = 0; i < m_expDim; ++i)
467 {
468 ptsTgt *= keyTgt[i].GetNumPoints();
469 }
470 }
471
472 // Get derivative at geometry points
473 DerivStorage deriv = ComputeDeriv(keyTgt);
474
475 // number of moving frames is requited to be 3, even for surfaces
476 int MFdim = 3;
477
479
480 // Compute g_{ij} as t_i \cdot t_j and store in tmp
481 for (i = 0; i < MFdim; ++i)
482 {
484 for (k = 0; k < m_coordDim; ++k)
485 {
486 MFtmp[i][k] = Array<OneD, NekDouble>(nq);
487 }
488 }
489
490 // Compute g_{ij} as t_i \cdot t_j and store in tmp
491 for (i = 0; i < MFdim - 1; ++i)
492 {
493 for (k = 0; k < m_coordDim; ++k)
494 {
495 if (m_type == eDeformed)
496 {
497 Vmath::Vcopy(ptsTgt, &deriv[i][k][0], 1, &MFtmp[i][k][0], 1);
498 }
499 else
500 {
501 Vmath::Fill(nq, deriv[i][k][0], MFtmp[i][k], 1);
502 }
503 }
504 }
505
506 // Direction of MF1 is preserved: MF2 is considered in the same
507 // tangent plane as MF1. MF3 is computed by cross product of MF1
508 // and MF2. MF2 is consequently computed as the cross product of
509 // MF3 and MF1.
511 for (k = 0; k < m_coordDim; k++)
512 {
513 PrincipleDir[k] = Array<OneD, NekDouble>(nq);
514 }
515
516 if (!(MMFdir == eLOCAL))
517 {
518 ComputePrincipleDirection(keyTgt, MMFdir, factors, PrincipleDir);
519 }
520
521 // MF3 = MF1 \times MF2
522 VectorCrossProd(MFtmp[0], MFtmp[1], MFtmp[2]);
523
524 // Normalizing MF3
525 VectorNormalise(MFtmp[2]);
526
527 if (!(MMFdir == eLOCAL))
528 {
529 Array<OneD, NekDouble> temp(nq, 0.0);
530
531 // Reorient MF1 along the PrincipleDir
532 for (i = 0; i < m_coordDim; ++i)
533 {
534 Vmath::Vvtvp(nq, MFtmp[2][i], 1, PrincipleDir[i], 1, temp, 1, temp,
535 1);
536 }
537 Vmath::Neg(nq, temp, 1);
538
539 // u2 = v2 - < u1 , v2 > ( u1 / < u1, u1 > )
540 for (i = 0; i < m_coordDim; ++i)
541 {
542 Vmath::Vvtvp(nq, temp, 1, MFtmp[2][i], 1, PrincipleDir[i], 1,
543 MFtmp[0][i], 1);
544 }
545 }
546
547 // Normalizing MF1
548 VectorNormalise(MFtmp[0]);
549
550 // MF2 = MF3 \times MF1
551 VectorCrossProd(MFtmp[2], MFtmp[0], MFtmp[1]);
552
553 // Normalizing MF2
554 VectorNormalise(MFtmp[1]);
555
556 for (i = 0; i < MFdim; ++i)
557 {
558 for (k = 0; k < m_coordDim; ++k)
559 {
560 Vmath::Vcopy(nq, &MFtmp[i][k][0], 1,
561 &movingframes[i * m_coordDim + k][0], 1);
562 }
563 }
564}
565
566/**
567 * Constructs the Jacobian as per Spencer's book p158 and tests if
568 * negative.
569 */
571{
572 // Jacobian test only makes sense when expdim = coorddim
573 // If one-dimensional then element is valid.
574 if (m_coordDim != m_expDim || m_expDim == 1)
575 {
576 m_valid = true;
577 return;
578 }
579
581 int nqtot = 1;
582 for (int i = 0; i < m_expDim; ++i)
583 {
584 p[i] = m_xmap->GetBasis(i)->GetPointsKey();
585 nqtot *= p[i].GetNumPoints();
586 }
587 int pts = (m_type == eRegular || m_type == eMovingRegular) ? 1 : nqtot;
588 Array<OneD, NekDouble> jac(pts, 0.0);
589
590 DerivStorage deriv = GetDeriv(p);
591
592 switch (m_expDim)
593 {
594 case 2:
595 {
596 Vmath::Vvtvvtm(pts, &deriv[0][0][0], 1, &deriv[1][1][0], 1,
597 &deriv[1][0][0], 1, &deriv[0][1][0], 1, &jac[0], 1);
598 break;
599 }
600 case 3:
601 {
602 Array<OneD, NekDouble> tmp(pts, 0.0);
603
604 Vmath::Vvtvvtm(pts, &deriv[1][1][0], 1, &deriv[2][2][0], 1,
605 &deriv[2][1][0], 1, &deriv[1][2][0], 1, &tmp[0], 1);
606 Vmath::Vvtvp(pts, &deriv[0][0][0], 1, &tmp[0], 1, &jac[0], 1,
607 &jac[0], 1);
608
609 Vmath::Vvtvvtm(pts, &deriv[2][1][0], 1, &deriv[0][2][0], 1,
610 &deriv[0][1][0], 1, &deriv[2][2][0], 1, &tmp[0], 1);
611 Vmath::Vvtvp(pts, &deriv[1][0][0], 1, &tmp[0], 1, &jac[0], 1,
612 &jac[0], 1);
613
614 Vmath::Vvtvvtm(pts, &deriv[0][1][0], 1, &deriv[1][2][0], 1,
615 &deriv[1][1][0], 1, &deriv[0][2][0], 1, &tmp[0], 1);
616 Vmath::Vvtvp(pts, &deriv[2][0][0], 1, &tmp[0], 1, &jac[0], 1,
617 &jac[0], 1);
618
619 break;
620 }
621 }
622
623 if (Vmath::Vmin(pts, &jac[0], 1) < 0)
624 {
625 m_valid = false;
626 }
627}
628
629/**
630 * @param map_points Source data point distribution.
631 * @param src Source data to be interpolated.
632 * @param tpoints Target data point distribution.
633 * @param tgt Target data storage.
634 */
637 const LibUtilities::PointsKeyVector &tgt_points,
638 Array<OneD, NekDouble> &tgt) const
639{
640 ASSERTL1(src_points.size() == tgt_points.size(),
641 "Dimension of target point distribution does not match "
642 "expansion dimension.");
643
644 switch (m_expDim)
645 {
646 case 1:
647 LibUtilities::Interp1D(src_points[0], src, tgt_points[0], tgt);
648 break;
649 case 2:
650 LibUtilities::Interp2D(src_points[0], src_points[1], src,
651 tgt_points[0], tgt_points[1], tgt);
652 break;
653 case 3:
654 LibUtilities::Interp3D(src_points[0], src_points[1], src_points[2],
655 src, tgt_points[0], tgt_points[1],
656 tgt_points[2], tgt);
657 break;
658 }
659}
660
661/**
662 * Input and output arrays are of dimension
663 * (m_expDim*m_expDim) x num_points. The first index of the input and
664 * output arrays are ordered row-by-row.
665 * @param src Input data array.
666 * @param tgt Storage for adjoint matrix data.
667 */
669 Array<TwoD, NekDouble> &tgt) const
670{
671 ASSERTL1(src.size() == tgt.size(),
672 "Source matrix is of different size to destination"
673 "matrix for computing adjoint.");
674
675 int n = src[0].size();
676 switch (m_expDim)
677 {
678 case 1:
679 Vmath::Fill(n, 1.0, &tgt[0][0], 1);
680 break;
681 case 2:
682 Vmath::Vcopy(n, &src[3][0], 1, &tgt[0][0], 1);
683 Vmath::Smul(n, -1.0, &src[1][0], 1, &tgt[1][0], 1);
684 Vmath::Smul(n, -1.0, &src[2][0], 1, &tgt[2][0], 1);
685 Vmath::Vcopy(n, &src[0][0], 1, &tgt[3][0], 1);
686 break;
687 case 3:
688 {
689 int a, b, c, d, e, i, j;
690
691 // Compute g^{ij} by computing Cofactors(g_ij)^T
692 for (i = 0; i < m_expDim; ++i)
693 {
694 for (j = 0; j < m_expDim; ++j)
695 {
696 a = ((i + 1) % m_expDim) * m_expDim + ((j + 1) % m_expDim);
697 b = ((i + 1) % m_expDim) * m_expDim + ((j + 2) % m_expDim);
698 c = ((i + 2) % m_expDim) * m_expDim + ((j + 1) % m_expDim);
699 d = ((i + 2) % m_expDim) * m_expDim + ((j + 2) % m_expDim);
700 e = j * m_expDim + i;
701 Vmath::Vvtvvtm(n, &src[a][0], 1, &src[d][0], 1, &src[b][0],
702 1, &src[c][0], 1, &tgt[e][0], 1);
703 }
704 }
705 break;
706 }
707 }
708}
709
710/**
711 *
712 */
714 const LibUtilities::PointsKeyVector &keyTgt,
715 const SpatialDomains::GeomMMF MMFdir,
718{
719 int nq = output[0].size();
720
722 for (int i = 0; i < m_coordDim; ++i)
723 {
724 output[i] = Array<OneD, NekDouble>(nq, 0.0);
725 }
726
727 // Construction of Connection
728 switch (MMFdir)
729 {
730 // projection to x-axis
731 case eTangentX:
732 {
733 Vmath::Fill(nq, 1.0, output[0], 1);
734 break;
735 }
736 case eTangentY:
737 {
738 Vmath::Fill(nq, 1.0, output[1], 1);
739 break;
740 }
741 case eTangentXY:
742 {
743 Vmath::Fill(nq, sqrt(2.0), output[0], 1);
744 Vmath::Fill(nq, sqrt(2.0), output[1], 1);
745 break;
746 }
747 case eTangentZ:
748 {
749 Vmath::Fill(nq, 1.0, output[2], 1);
750 break;
751 }
752 case eTangentCircular:
753 {
754 // Tangent direction depends on spatial location.
756 for (int k = 0; k < m_coordDim; k++)
757 {
758 x[k] = Array<OneD, NekDouble>(nq);
759 }
760
761 // m_coords are StdExpansions which store the mapping
762 // between the std element and the local element. Bwd
763 // transforming the std element minimum basis gives a
764 // minimum physical basis for geometry. Need to then
765 // interpolate this up to the quadrature basis.
766 int nqtot_map = 1;
768 for (int i = 0; i < m_expDim; ++i)
769 {
770 map_points[i] = m_xmap->GetBasis(i)->GetPointsKey();
771 nqtot_map *= map_points[i].GetNumPoints();
772 }
773 Array<OneD, NekDouble> tmp(nqtot_map);
774 for (int k = 0; k < m_coordDim; k++)
775 {
776 m_xmap->BwdTrans(m_coords[k], tmp);
777 Interp(map_points, tmp, keyTgt, x[k]);
778 }
779
780 // circular around the center of the domain
781 NekDouble radius, xc = 0.0, yc = 0.0, xdis, ydis;
782 NekDouble la, lb;
783
784 ASSERTL1(factors.size() >= 4, "factors is too short.");
785
786 la = factors[0];
787 lb = factors[1];
788 xc = factors[2];
789 yc = factors[3];
790
791 for (int i = 0; i < nq; i++)
792 {
793 xdis = x[0][i] - xc;
794 ydis = x[1][i] - yc;
795 radius = sqrt(xdis * xdis / la / la + ydis * ydis / lb / lb);
796 output[0][i] = ydis / radius;
797 output[1][i] = -1.0 * xdis / radius;
798 }
799 break;
800 }
802 {
803 // Tangent direction depends on spatial location.
805 for (int k = 0; k < m_coordDim; k++)
806 {
807 x[k] = Array<OneD, NekDouble>(nq);
808 }
809
810 int nqtot_map = 1;
812 for (int i = 0; i < m_expDim; ++i)
813 {
814 map_points[i] = m_xmap->GetBasis(i)->GetPointsKey();
815 nqtot_map *= map_points[i].GetNumPoints();
816 }
817 Array<OneD, NekDouble> tmp(nqtot_map);
818 for (int k = 0; k < m_coordDim; k++)
819 {
820 m_xmap->BwdTrans(m_coords[k], tmp);
821 Interp(map_points, tmp, keyTgt, x[k]);
822 }
823
824 // circular around the center of the domain
825 NekDouble xtan, ytan, mag;
826 for (int i = 0; i < nq; i++)
827 {
828 xtan = -1.0 * (x[1][i] * x[1][i] * x[1][i] + x[1][i]);
829 ytan = 2.0 * x[0][i];
830 mag = sqrt(xtan * xtan + ytan * ytan);
831 output[0][i] = xtan / mag;
832 output[1][i] = ytan / mag;
833 }
834 break;
835 }
837 {
838 // Tangent direction depends on spatial location.
840 for (int k = 0; k < m_coordDim; k++)
841 {
842 x[k] = Array<OneD, NekDouble>(nq);
843 }
844
845 int nqtot_map = 1;
847 for (int i = 0; i < m_expDim; ++i)
848 {
849 map_points[i] = m_xmap->GetBasis(i)->GetPointsKey();
850 nqtot_map *= map_points[i].GetNumPoints();
851 }
852 Array<OneD, NekDouble> tmp(nqtot_map);
853 for (int k = 0; k < m_coordDim; k++)
854 {
855 m_xmap->BwdTrans(m_coords[k], tmp);
856 Interp(map_points, tmp, keyTgt, x[k]);
857 }
858
859 // circular around the center of the domain
860 NekDouble xtan, ytan, mag;
861 for (int i = 0; i < nq; i++)
862 {
863 xtan = -2.0 * x[1][i] * x[1][i] * x[1][i] + x[1][i];
864 ytan = sqrt(3.0) * x[0][i];
865 mag = sqrt(xtan * xtan + ytan * ytan);
866 output[0][i] = xtan / mag;
867 output[1][i] = ytan / mag;
868 }
869 break;
870 }
871 default:
872 {
873 break;
874 }
875 }
876}
877
878/**
879 *
880 */
882{
883 int ndim = array.size();
884 ASSERTL0(ndim > 0, "Number of components must be > 0.");
885 for (int i = 1; i < ndim; ++i)
886 {
887 ASSERTL0(array[i].size() == array[0].size(),
888 "Array size mismatch in coordinates.");
889 }
890
891 int nq = array[0].size();
892 Array<OneD, NekDouble> norm(nq, 0.0);
893
894 // Compute the norm of each vector.
895 for (int i = 0; i < ndim; ++i)
896 {
897 Vmath::Vvtvp(nq, array[i], 1, array[i], 1, norm, 1, norm, 1);
898 }
899
900 Vmath::Vsqrt(nq, norm, 1, norm, 1);
901
902 // Normalise the vectors by the norm
903 for (int i = 0; i < ndim; ++i)
904 {
905 Vmath::Vdiv(nq, array[i], 1, norm, 1, array[i], 1);
906 }
907}
908
909/**
910 * @brief Computes the vector cross-product in 3D of \a v1 and \a v2,
911 * storing the result in \a v3.
912 *
913 * @param v1 First input vector.
914 * @param v2 Second input vector.
915 * @param v3 Output vector computed to be orthogonal to
916 * both \a v1 and \a v2.
917 */
919 const Array<OneD, const Array<OneD, NekDouble>> &v1,
920 const Array<OneD, const Array<OneD, NekDouble>> &v2,
922{
923 ASSERTL0(v1.size() == 3, "Input 1 has dimension not equal to 3.");
924 ASSERTL0(v2.size() == 3, "Input 2 has dimension not equal to 3.");
925 ASSERTL0(v3.size() == 3, "Output vector has dimension not equal to 3.");
926
927 int nq = v1[0].size();
928 Array<OneD, NekDouble> temp(nq);
929
930 Vmath::Vmul(nq, v1[2], 1, v2[1], 1, temp, 1);
931 Vmath::Vvtvm(nq, v1[1], 1, v2[2], 1, temp, 1, v3[0], 1);
932
933 Vmath::Vmul(nq, v1[0], 1, v2[2], 1, temp, 1);
934 Vmath::Vvtvm(nq, v1[2], 1, v2[0], 1, temp, 1, v3[1], 1);
935
936 Vmath::Vmul(nq, v1[1], 1, v2[0], 1, temp, 1);
937 Vmath::Vvtvm(nq, v1[0], 1, v2[1], 1, temp, 1, v3[2], 1);
938}
939
940} // 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:285