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