Nektar++
ExpList.h
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: ExpList.h
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: Expansion list top class definition
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#ifndef NEKTAR_LIBS_MULTIREGIONS_EXPLIST_H
36#define NEKTAR_LIBS_MULTIREGIONS_EXPLIST_H
37
52#include <tinyxml.h>
53
55{
56
57// Forward declarations
58class ExpList;
59class GlobalLinSys;
60class AssemblyMapDG;
61class AssemblyMapCG;
62class InterfaceMapDG;
63class GlobalLinSysKey;
64class GlobalMatrix;
65
67{
72 eN
73};
74
76{
85};
86
88
89/// A map between global matrix keys and their associated block
90/// matrices.
91typedef std::map<GlobalMatrixKey, DNekScalBlkMatSharedPtr> BlockMatrixMap;
92/// A shared pointer to a BlockMatrixMap.
93typedef std::shared_ptr<BlockMatrixMap> BlockMatrixMapShPtr;
94/// Shared pointer to an ExpList object.
95typedef std::shared_ptr<ExpList> ExpListSharedPtr;
96
97/// Base class for all multi-elemental spectral/hp expansions.
98class ExpList : public std::enable_shared_from_this<ExpList>
99{
100public:
101 /// The default constructor using a type
103
104 /// The copy constructor.
106 const bool DeclareCoeffPhysArrays = true);
107
108 /// Constructor copying only elements defined in eIds.
110 const std::vector<unsigned int> &eIDs,
111 const bool DeclareCoeffPhysArrays = true,
112 const Collections::ImplementationType ImpType =
114
115 /// Generate an ExpList from a meshgraph \a graph and session file
119 const bool DeclareCoeffPhysArrays = true,
120 const std::string &var = "DefaultVar",
121 const Collections::ImplementationType ImpType =
123
124 /// Sets up a list of local expansions based on an expansion Map
127 const SpatialDomains::ExpansionInfoMap &expansions,
128 const bool DeclareCoeffPhysArrays = true,
129 const Collections::ImplementationType ImpType =
131
132 //---------------------------------------------------------
133 // Specialised constructors in ExpListConstructor.cpp
134 //---------------------------------------------------------
135 /// Specialised constructors for 0D Expansions
136 /// Wrapper around LocalRegion::PointExp - used in PrePacing.cpp
139
140 /// Generate expansions for the trace space expansions used in
141 /// DisContField.
144 const Array<OneD, const ExpListSharedPtr> &bndConstraint,
146 &bndCond,
147 const LocalRegions::ExpansionVector &locexp,
149 const LibUtilities::CommSharedPtr &comm,
150 const bool DeclareCoeffPhysArrays = true,
151 const std::string variable = "DefaultVar",
152 const Collections::ImplementationType ImpType =
154
155 /// Generate an trace ExpList from a meshgraph \a graph and session file
158 const LocalRegions::ExpansionVector &locexp,
160 const bool DeclareCoeffPhysArrays, const std::string variable,
161 const Collections::ImplementationType ImpType =
163
164 /// Constructor based on domain information only for 1D &
165 /// 2D boundary conditions
168 const SpatialDomains::CompositeMap &domain,
170 const bool DeclareCoeffPhysArrays = true,
171 const std::string variable = "DefaultVar",
172 bool SetToOneSpaceDimension = false,
174 const Collections::ImplementationType ImpType =
176
177 /// The default destructor.
179
180 /// Returns the total number of local degrees of freedom
181 /// \f$N_{\mathrm{eof}}=\sum_{e=1}^{{N_{\mathrm{el}}}}N^{e}_m\f$.
182 inline int GetNcoeffs(void) const;
183
184 /// Returns the total number of local degrees of freedom
185 /// for element eid
186 MULTI_REGIONS_EXPORT int GetNcoeffs(const int eid) const;
187
188 /// Returns the type of the expansion
190
191 /// Returns the type of the expansion
193
194 /// Evaulates the maximum number of modes in the elemental basis
195 /// order over all elements
196 inline int EvalBasisNumModesMax(void) const;
197
198 /// Returns the vector of the number of modes in the elemental
199 /// basis order over all elements.
201 void) const;
202
203 /// Returns the total number of quadrature points #m_npoints
204 /// \f$=Q_{\mathrm{tot}}\f$.
205 inline int GetTotPoints(void) const;
206
207 /// Returns the total number of quadrature points for eid's element
208 /// \f$=Q_{\mathrm{tot}}\f$.
209 inline int GetTotPoints(const int eid) const;
210
211 /// Returns the total number of quadrature points #m_npoints
212 /// \f$=Q_{\mathrm{tot}}\f$.
213 inline int GetNpoints(void) const;
214
215 /// Returns the total number of qudature points scaled by
216 /// the factor scale on each 1D direction
217 inline int Get1DScaledTotPoints(const NekDouble scale) const;
218
219 /// Sets the wave space to the one of the possible configuration
220 /// true or false
221 inline void SetWaveSpace(const bool wavespace);
222
223 /// Set Modified Basis for the stability analysis
224 inline void SetModifiedBasis(const bool modbasis);
225
226 /// This function returns the third direction expansion condition,
227 /// which can be in wave space (coefficient) or not
228 /// It is stored in the variable m_WaveSpace.
229 inline bool GetWaveSpace(void) const;
230
231 /// Set the \a i th value of \a m_phys to value \a val
232 inline void SetPhys(int i, NekDouble val);
233 /// Fills the array #m_phys
234 inline void SetPhys(const Array<OneD, const NekDouble> &inarray);
235 /// Sets the array #m_phys
236 inline void SetPhysArray(Array<OneD, NekDouble> &inarray);
237 /// This function manually sets whether the array of physical
238 /// values \f$\boldsymbol{u}_l\f$ (implemented as #m_phys) is
239 /// filled or not.
240 inline void SetPhysState(const bool physState);
241 /// This function indicates whether the array of physical values
242 /// \f$\boldsymbol{u}_l\f$ (implemented as #m_phys) is filled or
243 /// not.
244 inline bool GetPhysState(void) const;
245
246 /// multiply the metric jacobi and quadrature weights
248 const Array<OneD, const NekDouble> &inarray,
249 Array<OneD, NekDouble> &outarray);
250 /// Divided by the metric jacobi and quadrature weights
252 const Array<OneD, const NekDouble> &inarray,
253 Array<OneD, NekDouble> &outarray);
254 /// This function calculates the inner product of a function
255 /// \f$f(\boldsymbol{x})\f$ with respect to all \em local
256 /// expansion modes \f$\phi_n^e(\boldsymbol{x})\f$.
257 inline void IProductWRTBase(const Array<OneD, const NekDouble> &inarray,
258 Array<OneD, NekDouble> &outarray);
259 /// This function calculates the inner product of a function
260 /// \f$f(\boldsymbol{x})\f$ with respect to the derivative (in
261 /// direction \param dir) of all \em local expansion modes
262 /// \f$\phi_n^e(\boldsymbol{x})\f$.
264 const int dir, const Array<OneD, const NekDouble> &inarray,
265 Array<OneD, NekDouble> &outarray);
266
268 const Array<OneD, const NekDouble> &direction,
269 const Array<OneD, const NekDouble> &inarray,
270 Array<OneD, NekDouble> &outarray);
271
272 /// This function calculates the inner product of a
273 /// function \f$f(\boldsymbol{x})\f$ with respect to the
274 /// derivative of all \em local expansion modes
275 /// \f$\phi_n^e(\boldsymbol{x})\f$.
277 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
278 Array<OneD, NekDouble> &outarray);
279 /// This function elementally evaluates the forward transformation
280 /// of a function \f$u(\boldsymbol{x})\f$ onto the global
281 /// spectral/hp expansion.
282 inline void FwdTransLocalElmt(const Array<OneD, const NekDouble> &inarray,
283 Array<OneD, NekDouble> &outarray);
284 ///
285 inline void FwdTrans(const Array<OneD, const NekDouble> &inarray,
286 Array<OneD, NekDouble> &outarray);
288 const NekDouble alpha,
289 const NekDouble exponent,
290 const NekDouble cutoff);
291 /// This function elementally mulplies the coefficient space of
292 /// Sin my the elemental inverse of the mass matrix.
294 const Array<OneD, const NekDouble> &inarray,
295 Array<OneD, NekDouble> &outarray);
296 inline void MultiplyByInvMassMatrix(
297 const Array<OneD, const NekDouble> &inarray,
298 Array<OneD, NekDouble> &outarray);
300 const Array<OneD, const NekDouble> &inarray,
301 Array<OneD, NekDouble> &outarray)
302 {
304 BwdTrans(inarray, tmp);
305 IProductWRTBase(tmp, outarray);
306 }
307 /// Smooth a field across elements
309
310 /// Solve helmholtz problem
312 const Array<OneD, const NekDouble> &inarray,
313 Array<OneD, NekDouble> &outarray,
316 const MultiRegions::VarFactorsMap &varfactors =
319 const bool PhysSpaceForcing = true);
320
321 /// Solve Advection Diffusion Reaction
323 const Array<OneD, const NekDouble> &inarray,
324 Array<OneD, NekDouble> &outarray,
327 const MultiRegions::VarFactorsMap &varfactors =
330 const bool PhysSpaceForcing = true);
331
332 /// Solve Advection Diffusion Reaction
334 const Array<OneD, Array<OneD, NekDouble>> &velocity,
335 const Array<OneD, const NekDouble> &inarray,
336 Array<OneD, NekDouble> &outarray, const NekDouble lambda,
338 ///
340 const Array<OneD, const NekDouble> &inarray,
341 Array<OneD, NekDouble> &outarray);
342 /// This function elementally evaluates the backward transformation
343 /// of the global spectral/hp element expansion.
344 inline void BwdTrans(const Array<OneD, const NekDouble> &inarray,
345 Array<OneD, NekDouble> &outarray);
346 /// This function calculates the coordinates of all the elemental
347 /// quadrature points \f$\boldsymbol{x}_i\f$.
348 inline void GetCoords(
349 Array<OneD, NekDouble> &coord_0,
352
353 inline void GetCoords(
354 const int eid, Array<OneD, NekDouble> &coord_0,
357
358 // Homogeneous transforms
359 inline void HomogeneousFwdTrans(const int npts,
360 const Array<OneD, const NekDouble> &inarray,
361 Array<OneD, NekDouble> &outarray,
362 bool Shuff = true, bool UnShuff = true);
363 inline void HomogeneousBwdTrans(const int npts,
364 const Array<OneD, const NekDouble> &inarray,
365 Array<OneD, NekDouble> &outarray,
366 bool Shuff = true, bool UnShuff = true);
367 inline void DealiasedProd(const int num_dofs,
368 const Array<OneD, NekDouble> &inarray1,
369 const Array<OneD, NekDouble> &inarray2,
370 Array<OneD, NekDouble> &outarray);
371 inline void DealiasedDotProd(
372 const int num_dofs, const Array<OneD, Array<OneD, NekDouble>> &inarray1,
373 const Array<OneD, Array<OneD, NekDouble>> &inarray2,
374 Array<OneD, Array<OneD, NekDouble>> &outarray);
375 inline void GetBCValues(Array<OneD, NekDouble> &BndVals,
376 const Array<OneD, NekDouble> &TotField, int BndID);
379 Array<OneD, NekDouble> &outarray,
380 int BndID);
381 inline void NormVectorIProductWRTBase(
383 Array<OneD, NekDouble> &outarray);
384 /// Apply geometry information to each expansion.
386 /// Reset geometry information and reset matrices
388 {
389 v_Reset();
390 }
391
392 /// Reset matrices
394
395 // not sure we really need these in ExpList
396 void WriteTecplotHeader(std::ostream &outfile, std::string var = "")
397 {
398 v_WriteTecplotHeader(outfile, var);
399 }
400 void WriteTecplotZone(std::ostream &outfile, int expansion = -1)
401 {
402 v_WriteTecplotZone(outfile, expansion);
403 }
404 void WriteTecplotField(std::ostream &outfile, int expansion = -1)
405 {
406 v_WriteTecplotField(outfile, expansion);
407 }
408 void WriteTecplotConnectivity(std::ostream &outfile, int expansion = -1)
409 {
410 v_WriteTecplotConnectivity(outfile, expansion);
411 }
412 MULTI_REGIONS_EXPORT void WriteVtkHeader(std::ostream &outfile);
413 MULTI_REGIONS_EXPORT void WriteVtkFooter(std::ostream &outfile);
414 void WriteVtkPieceHeader(std::ostream &outfile, int expansion,
415 int istrip = 0)
416 {
417 v_WriteVtkPieceHeader(outfile, expansion, istrip);
418 }
419 MULTI_REGIONS_EXPORT void WriteVtkPieceFooter(std::ostream &outfile,
420 int expansion);
421 void WriteVtkPieceData(std::ostream &outfile, int expansion,
422 std::string var = "v")
423 {
424 v_WriteVtkPieceData(outfile, expansion, var);
425 }
426
427 /// This function returns the dimension of the coordinates of the
428 /// element \a eid.
429 // inline
430 MULTI_REGIONS_EXPORT int GetCoordim(int eid);
431
432 /// Set the \a i th coefficiient in \a m_coeffs to value \a val
433 inline void SetCoeff(int i, NekDouble val);
434 /// Set the \a i th coefficiient in #m_coeffs to value \a val
435 inline void SetCoeffs(int i, NekDouble val);
436 /// Set the #m_coeffs array to inarray
437 inline void SetCoeffsArray(Array<OneD, NekDouble> &inarray);
438 /// This function returns the dimension of the shape of the
439 /// element \a eid.
440 // inline
442 /// This function returns (a reference to) the array
443 /// \f$\boldsymbol{\hat{u}}_l\f$ (implemented as #m_coeffs)
444 /// containing all local expansion coefficients.
445 inline const Array<OneD, const NekDouble> &GetCoeffs() const;
446 /// Impose Dirichlet Boundary Conditions onto Array
448 /// Fill Bnd Condition expansion from the values stored in expansion
449 inline void FillBndCondFromField(const Array<OneD, NekDouble> coeffs);
450 /// Fill Bnd Condition expansion in nreg from the values
451 /// stored in expansion
452 inline void FillBndCondFromField(const int nreg,
453 const Array<OneD, NekDouble> coeffs);
454 /// Gathers the global coefficients \f$\boldsymbol{\hat{u}}_g\f$
455 /// from the local coefficients \f$\boldsymbol{\hat{u}}_l\f$.
456 // inline
457 MULTI_REGIONS_EXPORT inline void LocalToGlobal(bool useComm = true);
459 const Array<OneD, const NekDouble> &inarray,
460 Array<OneD, NekDouble> &outarray, bool useComm = true);
461 /// Scatters from the global coefficients
462 /// \f$\boldsymbol{\hat{u}}_g\f$ to the local coefficients
463 /// \f$\boldsymbol{\hat{u}}_l\f$.
464 // inline
465 MULTI_REGIONS_EXPORT inline void GlobalToLocal(void);
466 /**
467 * This operation is evaluated as:
468 * \f{tabbing}
469 * \hspace{1cm} \= Do \= $e=$ $1, N_{\mathrm{el}}$ \\
470 * > > Do \= $i=$ $0,N_m^e-1$ \\
471 * > > > $\boldsymbol{\hat{u}}^{e}[i] = \mbox{sign}[e][i] \cdot
472 * \boldsymbol{\hat{u}}_g[\mbox{map}[e][i]]$ \\
473 * > > continue \\
474 * > continue
475 * \f}
476 * where \a map\f$[e][i]\f$ is the mapping array and \a
477 * sign\f$[e][i]\f$ is an array of similar dimensions ensuring the
478 * correct modal connectivity between the different elements (both
479 * these arrays are contained in the data member #m_locToGloMap). This
480 * operation is equivalent to the scatter operation
481 * \f$\boldsymbol{\hat{u}}_l=\mathcal{A}\boldsymbol{\hat{u}}_g\f$,
482 * where \f$\mathcal{A}\f$ is the
483 * \f$N_{\mathrm{eof}}\times N_{\mathrm{dof}}\f$ permutation matrix.
484 *
485 * @param inarray An array of size \f$N_\mathrm{dof}\f$
486 * containing the global degrees of freedom
487 * \f$\boldsymbol{x}_g\f$.
488 * @param outarray The resulting local degrees of freedom
489 * \f$\boldsymbol{x}_l\f$ will be stored in this
490 * array of size \f$N_\mathrm{eof}\f$.
491 */
493 const Array<OneD, const NekDouble> &inarray,
494 Array<OneD, NekDouble> &outarray);
495 /// Get the \a i th value (coefficient) of #m_coeffs
496 inline NekDouble GetCoeff(int i);
497 /// Get the \a i th value (coefficient) of #m_coeffs
498 inline NekDouble GetCoeffs(int i);
499 /// This function returns (a reference to) the array
500 /// \f$\boldsymbol{u}_l\f$ (implemented as #m_phys) containing the
501 /// function \f$u^{\delta}(\boldsymbol{x})\f$ evaluated at the
502 /// quadrature points.
503 // inline
505 /// This function calculates the \f$L_\infty\f$ error of the global
506 /// spectral/hp element approximation.
508 Linf(const Array<OneD, const NekDouble> &inarray,
510 /// This function calculates the \f$L_\infty\f$ error of the global
511 /// This function calculates the \f$L_2\f$ error with
512 /// respect to soln of the global
513 /// spectral/hp element approximation.
515 const Array<OneD, const NekDouble> &inarray,
517 {
518 return v_L2(inarray, soln);
519 }
520 /// Calculates the \f$H^1\f$ error of the global spectral/hp
521 /// element approximation.
523 H1(const Array<OneD, const NekDouble> &inarray,
525 /// Calculates the \f$H^1\f$ error of the global spectral/hp
526 /// element approximation.
527 /**
528 * The integration is evaluated locally, that is
529 * \f[\int
530 * f(\boldsymbol{x})d\boldsymbol{x}=\sum_{e=1}^{{N_{\mathrm{el}}}}
531 * \left\{\int_{\Omega_e}f(\boldsymbol{x})d\boldsymbol{x}\right\}, \f]
532 * where the integration over the separate elements is done by the
533 * function StdRegions#StdExpansion#Integral, which discretely
534 * evaluates the integral using Gaussian quadrature.
535 *
536 * Note that the array #m_phys should be filled with the values of the
537 * function \f$f(\boldsymbol{x})\f$ at the quadrature points
538 * \f$\boldsymbol{x}_i\f$.
539 *
540 * @return The value of the discretely evaluated integral
541 * \f$\int f(\boldsymbol{x})d\boldsymbol{x}\f$.
542 */
544 {
545 ASSERTL1(m_physState == true, "local physical space is not true ");
546 return Integral(m_phys);
547 }
548 /**
549 * The integration is evaluated locally, that is
550 * \f[\int
551 * f(\boldsymbol{x})d\boldsymbol{x}=\sum_{e=1}^{{N_{\mathrm{el}}}}
552 * \left\{\int_{\Omega_e}f(\boldsymbol{x})d\boldsymbol{x}\right\}, \f]
553 * where the integration over the separate elements is done by the
554 * function StdRegions#StdExpansion#Integral, which discretely
555 * evaluates the integral using Gaussian quadrature.
556 *
557 * @param inarray An array of size \f$Q_{\mathrm{tot}}\f$
558 * containing the values of the function
559 * \f$f(\boldsymbol{x})\f$ at the quadrature
560 * points \f$\boldsymbol{x}_i\f$.
561 * @return The value of the discretely evaluated integral
562 * \f$\int f(\boldsymbol{x})d\boldsymbol{x}\f$.
563 */
565 {
566 return v_Integral(inarray);
567 }
569 {
570 return v_VectorFlux(inarray);
571 }
572 /// This function calculates the energy associated with
573 /// each one of the modesof a 3D homogeneous nD expansion
575 {
576 return v_HomogeneousEnergy();
577 }
578
579 /// This function sets the Spectral Vanishing Viscosity
580 /// in homogeneous1D expansion.
582 {
584 }
585
586 /// This function returns a vector containing the wave
587 /// numbers in z-direction associated
588 /// with the 3D homogenous expansion. Required if a
589 /// parellelisation is applied in the Fourier direction
591 {
592 return v_GetZIDs();
593 }
594
595 /// This function returns the transposition class
596 /// associated with the homogeneous expansion.
598 {
599 return v_GetTransposition();
600 }
601
602 /// This function returns the Width of homogeneous direction
603 /// associated with the homogeneous expansion.
605 {
606 return v_GetHomoLen();
607 }
608
609 /// This function sets the Width of homogeneous direction
610 /// associated with the homogeneous expansion.
611 void SetHomoLen(const NekDouble lhom)
612 {
613 return v_SetHomoLen(lhom);
614 }
615
616 /// This function returns a vector containing the wave
617 /// numbers in y-direction associated
618 /// with the 3D homogenous expansion. Required if a
619 /// parellelisation is applied in the Fourier direction
621 {
622 return v_GetYIDs();
623 }
624
625 /// This function interpolates the physical space points in
626 /// \a inarray to \a outarray using the same points defined in the
627 /// expansion but where the number of points are rescaled
628 /// by \a 1DScale
630 const Array<OneD, NekDouble> &inarray,
631 Array<OneD, NekDouble> &outarray)
632 {
633 v_PhysInterp1DScaled(scale, inarray, outarray);
634 }
635
636 /// This function Galerkin projects the physical space points in
637 /// \a inarray to \a outarray where inarray is assumed to
638 /// be defined in the expansion but where the number of
639 /// points are rescaled by \a 1DScale
641 const Array<OneD, NekDouble> &inarray,
642 Array<OneD, NekDouble> &outarray)
643 {
644 v_PhysGalerkinProjection1DScaled(scale, inarray, outarray);
645 }
646
647 /// This function returns the number of elements in the expansion.
648 inline int GetExpSize(void);
649
650 /// This function returns the number of elements in the
651 /// expansion which may be different for a homogeoenous extended
652 /// expansionp.
653 inline size_t GetNumElmts(void)
654 {
655 return v_GetNumElmts();
656 }
657
658 /// This function returns the vector of elements in the expansion.
659 inline const std::shared_ptr<LocalRegions::ExpansionVector> GetExp() const;
660
661 /// This function returns (a shared pointer to) the local elemental
662 /// expansion of the \f$n^{\mathrm{th}}\f$ element.
663 inline LocalRegions::ExpansionSharedPtr &GetExp(int n) const;
664
665 /// This function returns (a shared pointer to) the local elemental
666 /// expansion of the \f$n^{\mathrm{th}}\f$ element given a global
667 /// geometry ID.
669
670 /// This function returns (a shared pointer to) the local elemental
671 /// expansion containing the arbitrary point given by \a gloCoord.
673 const Array<OneD, const NekDouble> &gloCoord);
674
675 /**
676 * @brief This function returns the index of the local elemental
677 * expansion containing the arbitrary point given by \a gloCoord,
678 * within a distance tolerance of tol.
679 *
680 * If returnNearestElmt is true and no element contains the point,
681 * this function returns the nearest element whose bounding box contains
682 * the point. The bounding box has a 10% margin in each direction.
683 *
684 * @param gloCoord (input) coordinate in physics space
685 * @param locCoords (output) local coordinate xi in the returned element
686 * @param tol distance tolerance to judge if a point lies in an
687 * element
688 * @param returnNearestElmt if true and no element contains this point, the
689 * nearest element whose bounding box contains this
690 * point is returned
691 * @param cachedId an initial guess of the most possible element index
692 * @param maxDistance if returnNearestElmt is set as true, the nearest
693 * element will be returned. But the distance of the
694 * nearest element and this point should be <=
695 * maxDistance.
696 *
697 * @return element index; if no element is found, -1 is returned.
698 */
700 const Array<OneD, const NekDouble> &gloCoord, NekDouble tol = 0.0,
701 bool returnNearestElmt = false, int cachedId = -1,
702 NekDouble maxDistance = 1e6);
703
704 /** This function returns the index and the Local
705 * Cartesian Coordinates \a locCoords of the local
706 * elemental expansion containing the arbitrary point
707 * given by \a gloCoords.
708 **/
710 const Array<OneD, const NekDouble> &gloCoords,
711 Array<OneD, NekDouble> &locCoords, NekDouble tol = 0.0,
712 bool returnNearestElmt = false, int cachedId = -1,
713 NekDouble maxDistance = 1e6);
714
715 /** This function return the expansion field value
716 * at the coordinates given as input.
717 **/
720 const Array<OneD, const NekDouble> &phys);
721
722 /// Get the start offset position for a local contiguous list of
723 /// coeffs correspoinding to element n.
724 inline int GetCoeff_Offset(int n) const;
725
726 /// Get the start offset position for a local contiguous list of
727 /// quadrature points in a full array correspoinding to element n.
728 inline int GetPhys_Offset(int n) const;
729
730 /// This function returns (a reference to) the array
731 /// \f$\boldsymbol{\hat{u}}_l\f$ (implemented as #m_coeffs)
732 /// containing all local expansion coefficients.
734 /// This function returns (a reference to) the array
735 /// \f$\boldsymbol{u}_l\f$ (implemented as #m_phys) containing the
736 /// function \f$u^{\delta}(\boldsymbol{x})\f$ evaluated at the
737 /// quadrature points.
739 inline void PhysDeriv(Direction edir,
740 const Array<OneD, const NekDouble> &inarray,
742 /// This function discretely evaluates the derivative of a function
743 /// \f$f(\boldsymbol{x})\f$ on the domain consisting of all
744 /// elements of the expansion.
745 inline void PhysDeriv(
746 const Array<OneD, const NekDouble> &inarray,
750 inline void PhysDeriv(const int dir,
751 const Array<OneD, const NekDouble> &inarray,
753 inline void CurlCurl(Array<OneD, Array<OneD, NekDouble>> &Vel,
755 inline void PhysDirectionalDeriv(
756 const Array<OneD, const NekDouble> &direction,
757 const Array<OneD, const NekDouble> &inarray,
758 Array<OneD, NekDouble> &outarray);
759 inline void GetMovingFrames(const SpatialDomains::GeomMMF MMFdir,
760 const Array<OneD, const NekDouble> &CircCentre,
761 Array<OneD, Array<OneD, NekDouble>> &outarray);
762 // functions associated with DisContField
765 /// Get the weight value for boundary conditions
767 /// Set the weight value for boundary conditions
768 inline void SetBndCondBwdWeight(const int index, const NekDouble value);
769 inline std::shared_ptr<ExpList> &UpdateBndCondExpansion(int i);
770 inline void Upwind(const Array<OneD, const NekDouble> &Vn,
774 inline void Upwind(const Array<OneD, const Array<OneD, NekDouble>> &Vec,
778 /**
779 * Return a reference to the trace space associated with this
780 * expansion list.
781 */
782 inline std::shared_ptr<ExpList> &GetTrace();
783 inline std::shared_ptr<AssemblyMapDG> &GetTraceMap(void);
784 inline std::shared_ptr<InterfaceMapDG> &GetInterfaceMap(void);
785 inline const Array<OneD, const int> &GetTraceBndMap(void);
786 inline void GetNormals(Array<OneD, Array<OneD, NekDouble>> &normals);
787 /// Get the length of elements in boundary normal direction
789 Array<OneD, NekDouble> &lengthsFwd, Array<OneD, NekDouble> &lengthsBwd);
790 /// Get the weight value for boundary conditions
791 /// for boundary average and jump calculations
793 Array<OneD, NekDouble> &weightJump);
795 Array<OneD, NekDouble> &outarray);
798 Array<OneD, NekDouble> &outarray);
804 bool FillBnd = true,
805 bool PutFwdInBwdOnBCs = false,
806 bool DoExchange = true);
807 inline void FillBwdWithBoundCond(const Array<OneD, NekDouble> &Fwd,
809 bool PutFwdInBwdOnBCs = false);
810 /// Add Fwd and Bwd value to field,
811 /// a reverse procedure of GetFwdBwdTracePhys
818 {
820 }
823 Array<OneD, NekDouble> &locTraceFwd,
824 Array<OneD, NekDouble> &locTraceBwd);
825 /// Fill Bwd with boundary conditions
826 inline void FillBwdWithBwdWeight(Array<OneD, NekDouble> &weightave,
827 Array<OneD, NekDouble> &weightjmp);
828 /// Copy and fill the Periodic boundaries
829 inline void PeriodicBwdCopy(const Array<OneD, const NekDouble> &Fwd,
831 inline const std::vector<bool> &GetLeftAdjacentFaces(void) const;
832 inline void ExtractTracePhys(Array<OneD, NekDouble> &outarray);
833 inline void ExtractTracePhys(const Array<OneD, const NekDouble> &inarray,
834 Array<OneD, NekDouble> &outarray);
839 inline void EvaluateBoundaryConditions(
840 const NekDouble time = 0.0, const std::string varName = "",
843 // Routines for continous matrix solution
844 /// This function calculates the result of the multiplication of a
845 /// matrix of type specified by \a mkey with a vector given by \a
846 /// inarray.
848 const GlobalMatrixKey &gkey,
849 const Array<OneD, const NekDouble> &inarray,
850 Array<OneD, NekDouble> &outarray);
851 inline void SetUpPhysNormals();
852 inline void GetBoundaryToElmtMap(Array<OneD, int> &ElmtID,
853 Array<OneD, int> &EdgeID);
854 virtual void GetBndElmtExpansion(int i, std::shared_ptr<ExpList> &result,
855 const bool DeclareCoeffPhysArrays = true);
856
857 inline void ExtractElmtToBndPhys(int i, const Array<OneD, NekDouble> &elmt,
858 Array<OneD, NekDouble> &boundary);
859
860 inline void ExtractPhysToBndElmt(int i,
862 Array<OneD, NekDouble> &bndElmt);
863
864 inline void ExtractPhysToBnd(int i,
867
868 inline void GetBoundaryNormals(
869 int i, Array<OneD, Array<OneD, NekDouble>> &normals);
870
872 std::vector<LibUtilities::FieldDefinitionsSharedPtr> &fielddef,
873 int NumHomoDir = 0,
876 std::vector<NekDouble> &HomoLen = LibUtilities::NullNekDoubleVector,
877 bool homoStrips = false,
878 std::vector<unsigned int> &HomoSIDs =
880 std::vector<unsigned int> &HomoZIDs =
882 std::vector<unsigned int> &HomoYIDs =
884
885 std::map<int, RobinBCInfoSharedPtr> GetRobinBCInfo()
886 {
887 return v_GetRobinBCInfo();
888 }
889
890 void GetPeriodicEntities(PeriodicMap &periodicVerts,
891 PeriodicMap &periodicEdges,
892 PeriodicMap &periodicFaces = NullPeriodicMap)
893 {
894 v_GetPeriodicEntities(periodicVerts, periodicEdges, periodicFaces);
895 }
896
897 std::vector<LibUtilities::FieldDefinitionsSharedPtr> GetFieldDefinitions()
898 {
899 return v_GetFieldDefinitions();
900 }
901
903 std::vector<LibUtilities::FieldDefinitionsSharedPtr> &fielddef)
904 {
905 v_GetFieldDefinitions(fielddef);
906 }
907
908 /// Append the element data listed in elements
909 /// fielddef->m_ElementIDs onto fielddata
911 std::vector<NekDouble> &fielddata)
912 {
913 v_AppendFieldData(fielddef, fielddata);
914 }
915 /// Append the data in coeffs listed in elements
916 /// fielddef->m_ElementIDs onto fielddata
918 std::vector<NekDouble> &fielddata,
920 {
921 v_AppendFieldData(fielddef, fielddata, coeffs);
922 }
923 /** \brief Extract the data in fielddata into the coeffs
924 * using the basic ExpList Elemental expansions rather
925 * than planes in homogeneous case
926 */
929 std::vector<NekDouble> &fielddata, std::string &field,
930 Array<OneD, NekDouble> &coeffs);
931 /** \brief Extract the data from fromField using
932 * fromExpList the coeffs using the basic ExpList
933 * Elemental expansions rather than planes in homogeneous
934 * case
935 */
937 const std::shared_ptr<ExpList> &fromExpList,
938 const Array<OneD, const NekDouble> &fromCoeffs,
939 Array<OneD, NekDouble> &toCoeffs);
940 // Extract data in fielddata into the m_coeffs_list for the 3D stability
941 // analysis (base flow is 2D)
944 std::vector<NekDouble> &fielddata, std::string &field,
946 std::unordered_map<int, int> zIdToPlane =
947 std::unordered_map<int, int>());
948 // Extract data from file fileName and put coefficents into array coefffs
950 const std::string &fileName, LibUtilities::CommSharedPtr comm,
951 const std::string &varName, Array<OneD, NekDouble> &coeffs);
953 const int ElementID, const NekDouble scalar1, const NekDouble scalar2,
954 Array<OneD, NekDouble> &outarray);
955 /// Returns a shared pointer to the current object.
956 std::shared_ptr<ExpList> GetSharedThisPtr()
957 {
958 return shared_from_this();
959 }
960 /// Returns the session object
961 std::shared_ptr<LibUtilities::SessionReader> GetSession() const
962 {
963 return m_session;
964 }
965 /// Returns the comm object
966 std::shared_ptr<LibUtilities::Comm> GetComm() const
967 {
968 return m_comm;
969 }
971 {
972 return m_graph;
973 }
974 // Wrapper functions for Homogeneous Expansions
976 {
977 return v_GetHomogeneousBasis();
978 }
979 std::shared_ptr<ExpList> &GetPlane(int n)
980 {
981 return v_GetPlane(n);
982 }
986
987 MULTI_REGIONS_EXPORT int GetPoolCount(std::string);
988
990
992 GlobalLinSys> &
994
995 /// Get m_coeffs to elemental value map
997 GetCoeffsToElmt() const;
1003 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
1004 const int nDirctn, Array<OneD, DNekMatSharedPtr> &mtxPerVar);
1006 const TensorOfArray3D<NekDouble> &inarray,
1009 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
1012 const Array<OneD, const NekDouble> &FwdFlux,
1013 const Array<OneD, const NekDouble> &BwdFlux,
1014 Array<OneD, NekDouble> &outarray)
1015 {
1016 v_AddTraceIntegralToOffDiag(FwdFlux, BwdFlux, outarray);
1017 }
1019 const int dir, const Array<OneD, const DNekMatSharedPtr> ElmtJacQuad,
1020 Array<OneD, DNekMatSharedPtr> ElmtJacCoef);
1022 const Array<OneD, const DNekMatSharedPtr> ElmtJacQuad,
1023 Array<OneD, DNekMatSharedPtr> ElmtJacCoef);
1025 GetLocTraceToTraceMap() const;
1026 // Return the internal vector which identifieds if trace
1027 // is left adjacent definiing which trace the normal
1028 // points otwards from
1029 MULTI_REGIONS_EXPORT std::vector<bool> &GetLeftAdjacentTraces(void);
1030
1031 /// This function returns the map of index inside m_exp to geom id
1032 MULTI_REGIONS_EXPORT inline const std::unordered_map<int, int> &GetElmtToExpId(
1033 void)
1034 {
1035 return m_elmtToExpId;
1036 }
1037
1038 /// This function returns the index inside m_exp for a given geom id
1040 {
1041 auto it = m_elmtToExpId.find(elmtId);
1042 ASSERTL0(it != m_elmtToExpId.end(), "Global geometry ID " +
1043 std::to_string(elmtId) +
1044 " not found in element ID to "
1045 "expansion ID map.")
1046 return it->second;
1047 }
1048
1049protected:
1050 /// Expansion type
1052 std::shared_ptr<DNekMat> GenGlobalMatrixFull(
1053 const GlobalLinSysKey &mkey,
1054 const std::shared_ptr<AssemblyMapCG> &locToGloMap);
1055 /// Communicator
1057 /// Session
1059 /// Mesh associated with this expansion list.
1061 /// The total number of local degrees of freedom. #m_ncoeffs
1062 /// \f$=N_{\mathrm{eof}}=\sum_{e=1}^{{N_{\mathrm{el}}}}N^{e}_l\f$
1064 /** The total number of quadrature points. #m_npoints
1065 *\f$=Q_{\mathrm{tot}}=\sum_{e=1}^{{N_{\mathrm{el}}}}N^{e}_Q\f$
1066 **/
1068 /**
1069 * \brief Concatenation of all local expansion coefficients.
1070 *
1071 * The array of length #m_ncoeffs\f$=N_{\mathrm{eof}}\f$ which is
1072 * the concatenation of the local expansion coefficients
1073 * \f$\hat{u}_n^e\f$ over all \f$N_{\mathrm{el}}\f$ elements
1074 * \f[\mathrm{\texttt{m\_coeffs}}=\boldsymbol{\hat{u}}_{l} =
1075 * \underline{\boldsymbol{\hat{u}}}^e = \left [ \begin{array}{c}
1076 * \boldsymbol{\hat{u}}^{1} \ \
1077 * \boldsymbol{\hat{u}}^{2} \ \
1078 * \vdots \ \
1079 * \boldsymbol{\hat{u}}^{{{N_{\mathrm{el}}}}} \end{array} \right ],
1080 * \quad
1081 * \mathrm{where}\quad \boldsymbol{\hat{u}}^{e}[n]=\hat{u}_n^{e}\f]
1082 */
1084 /**
1085 * \brief The global expansion evaluated at the quadrature points
1086 *
1087 * The array of length #m_npoints\f$=Q_{\mathrm{tot}}\f$ containing
1088 * the evaluation of \f$u^{\delta}(\boldsymbol{x})\f$ at the
1089 * quadrature points \f$\boldsymbol{x}_i\f$.
1090 * \f[\mathrm{\texttt{m\_phys}}=\boldsymbol{u}_{l} =
1091 * \underline{\boldsymbol{u}}^e = \left [ \begin{array}{c}
1092 * \boldsymbol{u}^{1} \ \
1093 * \boldsymbol{u}^{2} \ \
1094 * \vdots \ \
1095 * \boldsymbol{u}^{{{N_{\mathrm{el}}}}} \end{array} \right ],\quad
1096 * \mathrm{where}\quad
1097 * \boldsymbol{u}^{e}[i]=u^{\delta}(\boldsymbol{x}_i)\f]
1098 */
1100 /**
1101 * \brief The state of the array #m_phys.
1102 *
1103 * Indicates whether the array #m_phys, created to contain the
1104 * evaluation of \f$u^{\delta}(\boldsymbol{x})\f$ at the quadrature
1105 * points \f$\boldsymbol{x}_i\f$, is filled with these values.
1106 */
1108 /**
1109 * \brief The list of local expansions.
1110 *
1111 * The (shared pointer to the) vector containing (shared pointers
1112 * to) all local expansions. The fact that the local expansions are
1113 * all stored as a (pointer to a) #StdExpansion, the abstract base
1114 * class for all local expansions, allows a general implementation
1115 * where most of the routines for the derived classes are defined
1116 * in the #ExpList base class.
1117 */
1118 std::shared_ptr<LocalRegions::ExpansionVector> m_exp;
1120 /// Vector of bools to act as an initialise on first call flag
1121 std::vector<bool> m_collectionsDoInit;
1122 /// Offset of elemental data into the array #m_coeffs
1124 /// Offset of elemental data into the array #m_phys
1126 /// m_coeffs to elemental value map
1129 //@todo should this be in ExpList or
1130 // ExpListHomogeneous1D.cpp it's a bool which determine if
1131 // the expansion is in the wave space (coefficient space)
1132 // or not
1134 /// Mapping from geometry ID of element to index inside #m_exp
1135 std::unordered_map<int, int> m_elmtToExpId;
1136 /// This function assembles the block diagonal matrix of local
1137 /// matrices of the type \a mtype.
1140 void MultiplyByBlockMatrix(const GlobalMatrixKey &gkey,
1141 const Array<OneD, const NekDouble> &inarray,
1142 Array<OneD, NekDouble> &outarray);
1143 /// Generates a global matrix from the given key and map.
1144 std::shared_ptr<GlobalMatrix> GenGlobalMatrix(
1145 const GlobalMatrixKey &mkey,
1146 const std::shared_ptr<AssemblyMapCG> &locToGloMap);
1148 const std::shared_ptr<DNekMat> &Gmat,
1149 Array<OneD, NekDouble> &EigValsReal,
1150 Array<OneD, NekDouble> &EigValsImag,
1152 /// This operation constructs the global linear system of type \a
1153 /// mkey.
1154 std::shared_ptr<GlobalLinSys> GenGlobalLinSys(
1155 const GlobalLinSysKey &mkey,
1156 const std::shared_ptr<AssemblyMapCG> &locToGloMap);
1157 /// Generate a GlobalLinSys from information provided by the key
1158 /// "mkey" and the mapping provided in LocToGloBaseMap.
1159 std::shared_ptr<GlobalLinSys> GenGlobalBndLinSys(
1160 const GlobalLinSysKey &mkey, const AssemblyMapSharedPtr &locToGloMap);
1161 // Virtual prototypes
1162 virtual size_t v_GetNumElmts(void)
1163 {
1164 return (*m_exp).size();
1165 }
1169 virtual void v_SetBndCondBwdWeight(const int index, const NekDouble value);
1170 virtual std::shared_ptr<ExpList> &v_UpdateBndCondExpansion(int i);
1171 virtual void v_Upwind(const Array<OneD, const Array<OneD, NekDouble>> &Vec,
1175 virtual void v_Upwind(const Array<OneD, const NekDouble> &Vn,
1179 virtual std::shared_ptr<ExpList> &v_GetTrace();
1180 virtual std::shared_ptr<AssemblyMapDG> &v_GetTraceMap();
1181 virtual std::shared_ptr<InterfaceMapDG> &v_GetInterfaceMap();
1183 virtual const std::shared_ptr<LocTraceToTraceMap> &v_GetLocTraceToTraceMap(
1184 void) const;
1185 virtual std::vector<bool> &v_GetLeftAdjacentTraces(void);
1186 /// Populate \a normals with the normals of all expansions.
1187 virtual void v_GetNormals(Array<OneD, Array<OneD, NekDouble>> &normals);
1188 virtual void v_AddTraceIntegral(const Array<OneD, const NekDouble> &Fn,
1189 Array<OneD, NekDouble> &outarray);
1190 virtual void v_AddFwdBwdTraceIntegral(
1193 Array<OneD, NekDouble> &outarray);
1199 bool FillBnd = true,
1200 bool PutFwdInBwdOnBCs = false,
1201 bool DoExchange = true);
1202 virtual void v_FillBwdWithBoundCond(const Array<OneD, NekDouble> &Fwd,
1204 bool PutFwdInBwdOnBCs);
1205 virtual void v_AddTraceQuadPhysToField(
1208 virtual void v_AddTraceQuadPhysToOffDiag(
1211 virtual void v_GetLocTraceFromTracePts(
1214 Array<OneD, NekDouble> &locTraceFwd,
1215 Array<OneD, NekDouble> &locTraceBwd);
1216 virtual void v_FillBwdWithBwdWeight(Array<OneD, NekDouble> &weightave,
1217 Array<OneD, NekDouble> &weightjmp);
1218 virtual void v_PeriodicBwdCopy(const Array<OneD, const NekDouble> &Fwd,
1220 virtual const std::vector<bool> &v_GetLeftAdjacentFaces(void) const;
1221 virtual void v_ExtractTracePhys(Array<OneD, NekDouble> &outarray);
1222 virtual void v_ExtractTracePhys(const Array<OneD, const NekDouble> &inarray,
1223 Array<OneD, NekDouble> &outarray);
1224 virtual void v_MultiplyByInvMassMatrix(
1225 const Array<OneD, const NekDouble> &inarray,
1226 Array<OneD, NekDouble> &outarray);
1227
1229 const Array<OneD, const NekDouble> &inarray,
1230 Array<OneD, NekDouble> &outarray,
1232 const StdRegions::VarCoeffMap &varcoeff,
1233 const MultiRegions::VarFactorsMap &varfactors,
1234 const Array<OneD, const NekDouble> &dirForcing,
1235 const bool PhysSpaceForcing);
1236
1238 const Array<OneD, const NekDouble> &inarray,
1239 Array<OneD, NekDouble> &outarray,
1241 const StdRegions::VarCoeffMap &varcoeff,
1242 const MultiRegions::VarFactorsMap &varfactors,
1243 const Array<OneD, const NekDouble> &dirForcing,
1244 const bool PhysSpaceForcing);
1245
1246 virtual void v_LinearAdvectionReactionSolve(
1247 const Array<OneD, Array<OneD, NekDouble>> &velocity,
1248 const Array<OneD, const NekDouble> &inarray,
1249 Array<OneD, NekDouble> &outarray, const NekDouble lambda,
1251 // wrapper functions about virtual functions
1253 virtual void v_FillBndCondFromField(const Array<OneD, NekDouble> coeffs);
1254 virtual void v_FillBndCondFromField(const int nreg,
1255 const Array<OneD, NekDouble> coeffs);
1256 virtual void v_Reset();
1257 virtual void v_LocalToGlobal(bool UseComm);
1258 virtual void v_LocalToGlobal(const Array<OneD, const NekDouble> &inarray,
1259 Array<OneD, NekDouble> &outarray,
1260 bool UseComm);
1261 virtual void v_GlobalToLocal(void);
1262 virtual void v_GlobalToLocal(const Array<OneD, const NekDouble> &inarray,
1263 Array<OneD, NekDouble> &outarray);
1264 virtual void v_BwdTrans(const Array<OneD, const NekDouble> &inarray,
1265 Array<OneD, NekDouble> &outarray);
1266 virtual void v_FwdTrans(const Array<OneD, const NekDouble> &inarray,
1267 Array<OneD, NekDouble> &outarray);
1268 virtual void v_FwdTransLocalElmt(
1269 const Array<OneD, const NekDouble> &inarray,
1270 Array<OneD, NekDouble> &outarray);
1271 virtual void v_FwdTransBndConstrained(
1272 const Array<OneD, const NekDouble> &inarray,
1273 Array<OneD, NekDouble> &outarray);
1275 virtual void v_IProductWRTBase(const Array<OneD, const NekDouble> &inarray,
1276 Array<OneD, NekDouble> &outarray);
1277
1278 // Define ExpList::IProductWRTDerivBase as virtual function
1279 virtual void v_IProductWRTDerivBase(
1280 const int dir, const Array<OneD, const NekDouble> &inarray,
1281 Array<OneD, NekDouble> &outarray);
1282
1283 virtual void v_IProductWRTDerivBase(
1284 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
1285 Array<OneD, NekDouble> &outarray);
1286
1287 virtual void v_GetCoords(
1288 Array<OneD, NekDouble> &coord_0,
1291
1292 virtual void v_GetCoords(const int eid, Array<OneD, NekDouble> &xc0,
1295
1296 virtual void v_PhysDeriv(const Array<OneD, const NekDouble> &inarray,
1297 Array<OneD, NekDouble> &out_d0,
1298 Array<OneD, NekDouble> &out_d1,
1299 Array<OneD, NekDouble> &out_d2);
1300 virtual void v_PhysDeriv(const int dir,
1301 const Array<OneD, const NekDouble> &inarray,
1302 Array<OneD, NekDouble> &out_d);
1303 virtual void v_PhysDeriv(Direction edir,
1304 const Array<OneD, const NekDouble> &inarray,
1305 Array<OneD, NekDouble> &out_d);
1306
1307 virtual void v_Curl(Array<OneD, Array<OneD, NekDouble>> &Vel,
1309
1310 virtual void v_CurlCurl(Array<OneD, Array<OneD, NekDouble>> &Vel,
1312
1313 virtual void v_PhysDirectionalDeriv(
1314 const Array<OneD, const NekDouble> &direction,
1315 const Array<OneD, const NekDouble> &inarray,
1316 Array<OneD, NekDouble> &outarray);
1317 virtual void v_GetMovingFrames(
1318 const SpatialDomains::GeomMMF MMFdir,
1319 const Array<OneD, const NekDouble> &CircCentre,
1320 Array<OneD, Array<OneD, NekDouble>> &outarray);
1321 virtual void v_HomogeneousFwdTrans(
1322 const int npts, const Array<OneD, const NekDouble> &inarray,
1323 Array<OneD, NekDouble> &outarray, bool Shuff = true,
1324 bool UnShuff = true);
1325 virtual void v_HomogeneousBwdTrans(
1326 const int npts, const Array<OneD, const NekDouble> &inarray,
1327 Array<OneD, NekDouble> &outarray, bool Shuff = true,
1328 bool UnShuff = true);
1329 virtual void v_DealiasedProd(const int num_dofs,
1330 const Array<OneD, NekDouble> &inarray1,
1331 const Array<OneD, NekDouble> &inarray2,
1332 Array<OneD, NekDouble> &outarray);
1333 virtual void v_DealiasedDotProd(
1334 const int num_dofs, const Array<OneD, Array<OneD, NekDouble>> &inarray1,
1335 const Array<OneD, Array<OneD, NekDouble>> &inarray2,
1336 Array<OneD, Array<OneD, NekDouble>> &outarray);
1337 virtual void v_GetBCValues(Array<OneD, NekDouble> &BndVals,
1338 const Array<OneD, NekDouble> &TotField,
1339 int BndID);
1342 Array<OneD, NekDouble> &outarray,
1343 int BndID);
1344 virtual void v_NormVectorIProductWRTBase(
1346 Array<OneD, NekDouble> &outarray);
1347 virtual void v_SetUpPhysNormals();
1348 virtual void v_GetBoundaryToElmtMap(Array<OneD, int> &ElmtID,
1349 Array<OneD, int> &EdgeID);
1350 virtual void v_GetBndElmtExpansion(int i, std::shared_ptr<ExpList> &result,
1351 const bool DeclareCoeffPhysArrays);
1352
1353 virtual void v_ExtractElmtToBndPhys(const int i,
1354 const Array<OneD, NekDouble> &elmt,
1355 Array<OneD, NekDouble> &boundary);
1356
1357 virtual void v_ExtractPhysToBndElmt(
1358 const int i, const Array<OneD, const NekDouble> &phys,
1359 Array<OneD, NekDouble> &bndElmt);
1360
1361 virtual void v_ExtractPhysToBnd(const int i,
1362 const Array<OneD, const NekDouble> &phys,
1364
1365 virtual void v_GetBoundaryNormals(
1366 int i, Array<OneD, Array<OneD, NekDouble>> &normals);
1367
1368 virtual std::vector<LibUtilities::FieldDefinitionsSharedPtr>
1370
1371 virtual void v_GetFieldDefinitions(
1372 std::vector<LibUtilities::FieldDefinitionsSharedPtr> &fielddef);
1373
1374 virtual void v_AppendFieldData(
1376 std::vector<NekDouble> &fielddata);
1377 virtual void v_AppendFieldData(
1379 std::vector<NekDouble> &fielddata, Array<OneD, NekDouble> &coeffs);
1380 virtual void v_ExtractDataToCoeffs(
1382 std::vector<NekDouble> &fielddata, std::string &field,
1383 Array<OneD, NekDouble> &coeffs,
1384 std::unordered_map<int, int> zIdToPlane);
1385 virtual void v_ExtractCoeffsToCoeffs(
1386 const std::shared_ptr<ExpList> &fromExpList,
1387 const Array<OneD, const NekDouble> &fromCoeffs,
1388 Array<OneD, NekDouble> &toCoeffs);
1389 virtual void v_WriteTecplotHeader(std::ostream &outfile,
1390 std::string var = "");
1391 virtual void v_WriteTecplotZone(std::ostream &outfile, int expansion);
1392 virtual void v_WriteTecplotField(std::ostream &outfile, int expansion);
1393 virtual void v_WriteTecplotConnectivity(std::ostream &outfile,
1394 int expansion);
1395 virtual void v_WriteVtkPieceData(std::ostream &outfile, int expansion,
1396 std::string var);
1397 virtual void v_WriteVtkPieceHeader(std::ostream &outfile, int expansion,
1398 int istrip);
1399
1400 virtual NekDouble v_L2(
1401 const Array<OneD, const NekDouble> &phys,
1403
1404 virtual NekDouble v_Integral(const Array<OneD, const NekDouble> &inarray);
1405 virtual NekDouble v_VectorFlux(
1406 const Array<OneD, Array<OneD, NekDouble>> &inarray);
1407
1409
1411 virtual NekDouble v_GetHomoLen(void);
1412 virtual void v_SetHomoLen(const NekDouble lhom);
1415
1416 // 1D Scaling and projection
1417 virtual void v_PhysInterp1DScaled(const NekDouble scale,
1418 const Array<OneD, NekDouble> &inarray,
1419 Array<OneD, NekDouble> &outarray);
1420
1422 const NekDouble scale, const Array<OneD, NekDouble> &inarray,
1423 Array<OneD, NekDouble> &outarray);
1424
1425 virtual void v_ClearGlobalLinSysManager(void);
1426
1427 virtual int v_GetPoolCount(std::string);
1428
1429 virtual void v_UnsetGlobalLinSys(GlobalLinSysKey, bool);
1430
1433
1434 void ExtractFileBCs(const std::string &fileName,
1436 const std::string &varName,
1437 const std::shared_ptr<ExpList> locExpList);
1438
1439 // Utility function for a common case of retrieving a
1440 // BoundaryCondition from a boundary condition collection.
1444 unsigned int index, const std::string &variable);
1445
1448
1451
1452 virtual void v_EvaluateBoundaryConditions(
1453 const NekDouble time = 0.0, const std::string varName = "",
1456
1457 virtual std::map<int, RobinBCInfoSharedPtr> v_GetRobinBCInfo(void);
1458
1459 virtual void v_GetPeriodicEntities(PeriodicMap &periodicVerts,
1460 PeriodicMap &periodicEdges,
1461 PeriodicMap &periodicFaces);
1462
1463 // Homogeneous direction wrapper functions.
1465 {
1467 "This method is not defined or valid for this class type");
1469 }
1470
1471 // wrapper function to set viscosity for Homo1D expansion
1473 [[maybe_unused]] Array<OneD, NekDouble> visc)
1474 {
1476 "This method is not defined or valid for this class type");
1477 }
1478
1479 virtual std::shared_ptr<ExpList> &v_GetPlane(int n);
1480
1481 virtual void v_AddTraceIntegralToOffDiag(
1482 const Array<OneD, const NekDouble> &FwdFlux,
1483 const Array<OneD, const NekDouble> &BwdFlux,
1484 Array<OneD, NekDouble> &outarray);
1485
1486private:
1487 /// Definition of the total number of degrees of freedom and
1488 /// quadrature points and offsets to access data
1489 void SetupCoeffPhys(bool DeclareCoeffPhysArrays = true,
1490 bool SetupOffsets = true);
1491
1492 /// Define a list of elements using the geometry and basis
1493 /// key information in expmap;
1495};
1496
1497/// An empty ExpList object.
1500
1501// An empty GlobaLinSysManager object
1505
1506// Inline routines follow.
1507
1508/**
1509 * Returns the total number of local degrees of freedom
1510 * \f$N_{\mathrm{eof}}=\sum_{e=1}^{{N_{\mathrm{el}}}}N^{e}_m\f$.
1511 */
1512inline int ExpList::GetNcoeffs() const
1513{
1514 return m_ncoeffs;
1515}
1516
1517inline int ExpList::GetNcoeffs(const int eid) const
1518{
1519 return (*m_exp)[eid]->GetNcoeffs();
1520}
1521
1522/**
1523 * Evaulates the maximum number of modes in the elemental basis
1524 * order over all elements
1525 */
1527{
1528 int returnval = 0;
1529
1530 for (size_t i = 0; i < (*m_exp).size(); ++i)
1531 {
1532 returnval = (std::max)(returnval, (*m_exp)[i]->EvalBasisNumModesMax());
1533 }
1534
1535 return returnval;
1536}
1537
1538/**
1539 *
1540 */
1542{
1543 Array<OneD, int> returnval((*m_exp).size(), 0);
1544
1545 for (size_t i = 0; i < (*m_exp).size(); ++i)
1546 {
1547 returnval[i] =
1548 (std::max)(returnval[i], (*m_exp)[i]->EvalBasisNumModesMax());
1549 }
1550
1551 return returnval;
1552}
1553
1554/**
1555 *
1556 */
1557inline int ExpList::GetTotPoints() const
1558{
1559 return m_npoints;
1560}
1561
1562inline int ExpList::GetTotPoints(const int eid) const
1563{
1564 return (*m_exp)[eid]->GetTotPoints();
1565}
1566
1567inline int ExpList::Get1DScaledTotPoints(const NekDouble scale) const
1568{
1569 int returnval = 0;
1570 size_t cnt;
1571 size_t nbase = (*m_exp)[0]->GetNumBases();
1572
1573 for (size_t i = 0; i < (*m_exp).size(); ++i)
1574 {
1575 cnt = 1;
1576 for (size_t j = 0; j < nbase; ++j)
1577 {
1578 cnt *= scale * ((*m_exp)[i]->GetNumPoints(j));
1579 }
1580 returnval += cnt;
1581 }
1582 return returnval;
1583}
1584
1585/**
1586 *
1587 */
1588inline int ExpList::GetNpoints() const
1589{
1590 return m_npoints;
1591}
1592
1593/**
1594 *
1595 */
1596inline void ExpList::SetWaveSpace(const bool wavespace)
1597{
1598 m_WaveSpace = wavespace;
1599}
1600
1601/**
1602 *
1603 */
1604inline bool ExpList::GetWaveSpace() const
1605{
1606 return m_WaveSpace;
1607}
1608
1609/// Set the \a i th value of\a m_phys to value \a val
1610inline void ExpList::SetPhys(int i, NekDouble val)
1611{
1612 m_phys[i] = val;
1613}
1614/**
1615 * This function fills the array \f$\boldsymbol{u}_l\f$, the evaluation
1616 * of the expansion at the quadrature points (implemented as #m_phys),
1617 * with the values of the array \a inarray.
1618 *
1619 * @param inarray The array containing the values where
1620 * #m_phys should be filled with.
1621 */
1623{
1624 ASSERTL0((int)inarray.size() == m_npoints,
1625 "Input array does not have correct number of elements.");
1626 Vmath::Vcopy(m_npoints, &inarray[0], 1, &m_phys[0], 1);
1627 m_physState = true;
1628}
1630{
1631 m_phys = inarray;
1632}
1633/**
1634 * @param physState \a true (=filled) or \a false (=not filled).
1635 */
1636inline void ExpList::SetPhysState(const bool physState)
1637{
1638 m_physState = physState;
1639}
1640/**
1641 * @return physState \a true (=filled) or \a false (=not filled).
1642 */
1643inline bool ExpList::GetPhysState() const
1644{
1645 return m_physState;
1646}
1647/**
1648 *
1649 */
1651 const Array<OneD, const NekDouble> &inarray,
1652 Array<OneD, NekDouble> &outarray)
1653{
1654 v_IProductWRTBase(inarray, outarray);
1655}
1656/**
1657 *
1658 */
1660 const int dir, const Array<OneD, const NekDouble> &inarray,
1661 Array<OneD, NekDouble> &outarray)
1662{
1663 v_IProductWRTDerivBase(dir, inarray, outarray);
1664}
1665
1666/**
1667 *
1668 */
1670 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
1671 Array<OneD, NekDouble> &outarray)
1672{
1673 v_IProductWRTDerivBase(inarray, outarray);
1674}
1675
1676/**
1677 *
1678 */
1680 Array<OneD, NekDouble> &outarray)
1681{
1682 v_FwdTrans(inarray, outarray);
1683}
1684/**
1685 *
1686 */
1688 const Array<OneD, const NekDouble> &inarray,
1689 Array<OneD, NekDouble> &outarray)
1690{
1691 v_FwdTransLocalElmt(inarray, outarray);
1692}
1693/**
1694 *
1695 */
1697 const Array<OneD, const NekDouble> &inarray,
1698 Array<OneD, NekDouble> &outarray)
1699{
1700 v_FwdTransBndConstrained(inarray, outarray);
1701}
1702/**
1703 *
1704 */
1706{
1708}
1709/**
1710 *
1711 */
1712/**
1713 * Given the coefficients of an expansion, this function evaluates the
1714 * spectral/hp expansion \f$u^{\delta}(\boldsymbol{x})\f$ at the
1715 * quadrature points \f$\boldsymbol{x}_i\f$.
1716 */
1718 Array<OneD, NekDouble> &outarray)
1719{
1720 v_BwdTrans(inarray, outarray);
1721}
1722/**
1723 *
1724 */
1726 const Array<OneD, const NekDouble> &inarray,
1727 Array<OneD, NekDouble> &outarray)
1728{
1729 v_MultiplyByInvMassMatrix(inarray, outarray);
1730}
1731/**
1732 *
1733 */
1735 const Array<OneD, const NekDouble> &inarray,
1737 const StdRegions::VarCoeffMap &varcoeff,
1738 const MultiRegions::VarFactorsMap &varfactors,
1739 const Array<OneD, const NekDouble> &dirForcing, const bool PhysSpaceForcing)
1740{
1741 return v_HelmSolve(inarray, outarray, factors, varcoeff, varfactors,
1742 dirForcing, PhysSpaceForcing);
1743}
1744/**
1745 *
1746 */
1748 const Array<OneD, const NekDouble> &inarray,
1750 const StdRegions::VarCoeffMap &varcoeff,
1751 const MultiRegions::VarFactorsMap &varfactors,
1752 const Array<OneD, const NekDouble> &dirForcing, const bool PhysSpaceForcing)
1753{
1755 inarray, outarray, factors, varcoeff, varfactors, dirForcing,
1756 PhysSpaceForcing);
1757}
1759 const Array<OneD, Array<OneD, NekDouble>> &velocity,
1760 const Array<OneD, const NekDouble> &inarray,
1761 Array<OneD, NekDouble> &outarray, const NekDouble lambda,
1762 const Array<OneD, const NekDouble> &dirForcing)
1763{
1764 v_LinearAdvectionReactionSolve(velocity, inarray, outarray, lambda,
1765 dirForcing);
1766}
1767/**
1768 *
1769 */
1771 Array<OneD, NekDouble> &coord_1,
1772 Array<OneD, NekDouble> &coord_2)
1773{
1774 v_GetCoords(coord_0, coord_1, coord_2);
1775}
1776
1777inline void ExpList::GetCoords(const int eid, Array<OneD, NekDouble> &xc0,
1780{
1781 v_GetCoords(eid, xc0, xc1, xc2);
1782}
1783
1784/**
1785 *
1786 */
1788 const SpatialDomains::GeomMMF MMFdir,
1789 const Array<OneD, const NekDouble> &CircCentre,
1790 Array<OneD, Array<OneD, NekDouble>> &outarray)
1791{
1792 v_GetMovingFrames(MMFdir, CircCentre, outarray);
1793}
1794/**
1795 *
1796 */
1798 Array<OneD, NekDouble> &out_d0,
1799 Array<OneD, NekDouble> &out_d1,
1800 Array<OneD, NekDouble> &out_d2)
1801{
1802 v_PhysDeriv(inarray, out_d0, out_d1, out_d2);
1803}
1804/**
1805 *
1806 */
1807inline void ExpList::PhysDeriv(const int dir,
1808 const Array<OneD, const NekDouble> &inarray,
1810{
1811 v_PhysDeriv(dir, inarray, out_d);
1812}
1814 const Array<OneD, const NekDouble> &inarray,
1816{
1817 v_PhysDeriv(edir, inarray, out_d);
1818}
1819/**
1820 *
1821 */
1823 const Array<OneD, const NekDouble> &direction,
1824 const Array<OneD, const NekDouble> &inarray,
1825 Array<OneD, NekDouble> &outarray)
1826{
1827 v_PhysDirectionalDeriv(direction, inarray, outarray);
1828}
1829/**
1830 *
1831 */
1834{
1835 v_CurlCurl(Vel, Q);
1836}
1837/**
1838 *
1839 */
1841 const int npts, const Array<OneD, const NekDouble> &inarray,
1842 Array<OneD, NekDouble> &outarray, bool Shuff, bool UnShuff)
1843{
1844 v_HomogeneousFwdTrans(npts, inarray, outarray, Shuff, UnShuff);
1845}
1846/**
1847 *
1848 */
1850 const int npts, const Array<OneD, const NekDouble> &inarray,
1851 Array<OneD, NekDouble> &outarray, bool Shuff, bool UnShuff)
1852{
1853 v_HomogeneousBwdTrans(npts, inarray, outarray, Shuff, UnShuff);
1854}
1855/**
1856 *
1857 */
1858inline void ExpList::DealiasedProd(const int num_dofs,
1859 const Array<OneD, NekDouble> &inarray1,
1860 const Array<OneD, NekDouble> &inarray2,
1861 Array<OneD, NekDouble> &outarray)
1862{
1863 v_DealiasedProd(num_dofs, inarray1, inarray2, outarray);
1864}
1865/**
1866 *
1867 */
1869 const int num_dofs, const Array<OneD, Array<OneD, NekDouble>> &inarray1,
1870 const Array<OneD, Array<OneD, NekDouble>> &inarray2,
1871 Array<OneD, Array<OneD, NekDouble>> &outarray)
1872{
1873 v_DealiasedDotProd(num_dofs, inarray1, inarray2, outarray);
1874}
1875/**
1876 *
1877 */
1879 const Array<OneD, NekDouble> &TotField,
1880 int BndID)
1881{
1882 v_GetBCValues(BndVals, TotField, BndID);
1883}
1884/**
1885 *
1886 */
1889 Array<OneD, NekDouble> &outarray,
1890 int BndID)
1891{
1892 v_NormVectorIProductWRTBase(V1, V2, outarray, BndID);
1893}
1896{
1897 v_NormVectorIProductWRTBase(V, outarray);
1898}
1899/**
1900 * @param eid The index of the element to be checked.
1901 * @return The dimension of the coordinates of the specific element.
1902 */
1903inline int ExpList::GetCoordim(int eid)
1904{
1905 ASSERTL2(eid <= (*m_exp).size(), "eid is larger than number of elements");
1906 return (*m_exp)[eid]->GetCoordim();
1907}
1908/**
1909 * @param eid The index of the element to be checked.
1910 * @return The dimension of the shape of the specific element.
1911 */
1913{
1914 return (*m_exp)[0]->GetShapeDimension();
1915}
1916/**
1917 * @param i The index of m_coeffs to be set
1918 * @param val The value which m_coeffs[i] is to be set to.
1919 */
1920inline void ExpList::SetCoeff(int i, NekDouble val)
1921{
1922 m_coeffs[i] = val;
1923}
1924/**
1925 * @param i The index of #m_coeffs to be set.
1926 * @param val The value which #m_coeffs[i] is to be set to.
1927 */
1928inline void ExpList::SetCoeffs(int i, NekDouble val)
1929{
1930 m_coeffs[i] = val;
1931}
1933{
1934 m_coeffs = inarray;
1935}
1936/**
1937 * As the function returns a constant reference to a
1938 * <em>const Array</em>, it is not possible to modify the
1939 * underlying data of the array #m_coeffs. In order to do
1940 * so, use the function #UpdateCoeffs instead.
1941 *
1942 * @return (A constant reference to) the array #m_coeffs.
1943 */
1945{
1946 return m_coeffs;
1947}
1949{
1951}
1953{
1954 v_FillBndCondFromField(coeffs);
1955}
1956inline void ExpList::FillBndCondFromField(const int nreg,
1957 const Array<OneD, NekDouble> coeffs)
1958{
1959 v_FillBndCondFromField(nreg, coeffs);
1960}
1961inline void ExpList::LocalToGlobal(bool useComm)
1962{
1963 v_LocalToGlobal(useComm);
1964}
1966 Array<OneD, NekDouble> &outarray,
1967 bool useComm)
1968{
1969 v_LocalToGlobal(inarray, outarray, useComm);
1970}
1971inline void ExpList::GlobalToLocal(void)
1972{
1974}
1975/**
1976 * This operation is evaluated as:
1977 * \f{tabbing}
1978 * \hspace{1cm} \= Do \= $e=$ $1, N_{\mathrm{el}}$ \\
1979 * > > Do \= $i=$ $0,N_m^e-1$ \\
1980 * > > > $\boldsymbol{\hat{u}}^{e}[i] = \mbox{sign}[e][i] \cdot
1981 * \boldsymbol{\hat{u}}_g[\mbox{map}[e][i]]$ \\
1982 * > > continue \\
1983 * > continue
1984 * \f}
1985 * where \a map\f$[e][i]\f$ is the mapping array and \a
1986 * sign\f$[e][i]\f$ is an array of similar dimensions ensuring the
1987 * correct modal connectivity between the different elements (both
1988 * these arrays are contained in the data member #m_locToGloMap). This
1989 * operation is equivalent to the scatter operation
1990 * \f$\boldsymbol{\hat{u}}_l=\mathcal{A}\boldsymbol{\hat{u}}_g\f$, where
1991 * \f$\mathcal{A}\f$ is the
1992 * \f$N_{\mathrm{eof}}\times N_{\mathrm{dof}}\f$ permutation matrix.
1993 *
1994 * @param inarray An array of size \f$N_\mathrm{dof}\f$
1995 * containing the global degrees of freedom
1996 * \f$\boldsymbol{x}_g\f$.
1997 * @param outarray The resulting local degrees of freedom
1998 * \f$\boldsymbol{x}_l\f$ will be stored in this
1999 * array of size \f$N_\mathrm{eof}\f$.
2000 */
2002 Array<OneD, NekDouble> &outarray)
2003{
2004 v_GlobalToLocal(inarray, outarray);
2005}
2006/**
2007 * @param i The index of #m_coeffs to be returned
2008 * @return The NekDouble held in #m_coeffs[i].
2009 */
2011{
2012 return m_coeffs[i];
2013}
2014/**
2015 * @param i The index of #m_coeffs to be returned
2016 * @return The NekDouble held in #m_coeffs[i].
2017 */
2019{
2020 return m_coeffs[i];
2021}
2022/**
2023 * As the function returns a constant reference to a
2024 * <em>const Array</em> it is not possible to modify the
2025 * underlying data of the array #m_phys. In order to do
2026 * so, use the function #UpdatePhys instead.
2027 *
2028 * @return (A constant reference to) the array #m_phys.
2029 */
2031{
2032 return m_phys;
2033}
2034/**
2035 * @return \f$N_{\mathrm{el}}\f$, the number of elements in the
2036 * expansion.
2037 */
2038inline int ExpList::GetExpSize(void)
2039{
2040 return (*m_exp).size();
2041}
2042/**
2043 * @param n The index of the element concerned.
2044 *
2045 * @return (A shared pointer to) the local expansion of the
2046 * \f$n^{\mathrm{th}}\f$ element.
2047 */
2049{
2050 return (*m_exp)[n];
2051}
2052/**
2053 * @param n The global id of the element concerned.
2054 *
2055 * @return (A shared pointer to) the local expansion of the
2056 * \f$n^{\mathrm{th}}\f$ element.
2057 */
2059{
2060 auto it = m_elmtToExpId.find(n);
2061 ASSERTL0(it != m_elmtToExpId.end(), "Global geometry ID " +
2062 std::to_string(n) +
2063 " not found in element ID to "
2064 "expansion ID map.")
2065 return (*m_exp)[it->second];
2066}
2067/**
2068 * @return (A const shared pointer to) the local expansion vector #m_exp
2069 */
2070inline const std::shared_ptr<LocalRegions::ExpansionVector> ExpList::GetExp(
2071 void) const
2072{
2073 return m_exp;
2074}
2075/**
2076 *
2077 */
2078inline int ExpList::GetCoeff_Offset(int n) const
2079{
2080 return m_coeff_offset[n];
2081}
2082/**
2083 *
2084 */
2085inline int ExpList::GetPhys_Offset(int n) const
2086{
2087 return m_phys_offset[n];
2088}
2089/**
2090 * If one wants to get hold of the underlying data without modifying
2091 * them, rather use the function #GetCoeffs instead.
2092 *
2093 * @return (A reference to) the array #m_coeffs.
2094 */
2096{
2097 return m_coeffs;
2098}
2099/**
2100 * If one wants to get hold of the underlying data without modifying
2101 * them, rather use the function #GetPhys instead.
2102 *
2103 * @return (A reference to) the array #m_phys.
2104 */
2106{
2107 m_physState = true;
2108 return m_phys;
2109}
2110// functions associated with DisContField
2113{
2114 return v_GetBndCondExpansions();
2115}
2116/// Get m_coeffs to elemental value map
2118 GetCoeffsToElmt() const
2119{
2120 return m_coeffsToElmt;
2121}
2124{
2125 return v_GetLocTraceToTraceMap();
2126}
2128{
2129 return v_GetBndCondBwdWeight();
2130}
2131inline void ExpList::SetBndCondBwdWeight(const int index, const NekDouble value)
2132{
2133 v_SetBndCondBwdWeight(index, value);
2134}
2135inline std::shared_ptr<ExpList> &ExpList::UpdateBndCondExpansion(int i)
2136{
2137 return v_UpdateBndCondExpansion(i);
2138}
2140 const Array<OneD, const Array<OneD, NekDouble>> &Vec,
2143{
2144 v_Upwind(Vec, Fwd, Bwd, Upwind);
2145}
2149 Array<OneD, NekDouble> &Upwind)
2150{
2151 v_Upwind(Vn, Fwd, Bwd, Upwind);
2152}
2153inline std::shared_ptr<ExpList> &ExpList::GetTrace()
2154{
2155 return v_GetTrace();
2156}
2157inline std::shared_ptr<AssemblyMapDG> &ExpList::GetTraceMap()
2158{
2159 return v_GetTraceMap();
2160}
2161inline std::shared_ptr<InterfaceMapDG> &ExpList::GetInterfaceMap()
2162{
2163 return v_GetInterfaceMap();
2164}
2166{
2167 return v_GetTraceBndMap();
2168}
2170{
2171 v_GetNormals(normals);
2172}
2174 Array<OneD, NekDouble> &outarray)
2175{
2176 v_AddTraceIntegral(Fn, outarray);
2177}
2181{
2182 v_AddFwdBwdTraceIntegral(Fwd, Bwd, outarray);
2183}
2186{
2187 v_GetFwdBwdTracePhys(Fwd, Bwd);
2188}
2191 Array<OneD, NekDouble> &Bwd, bool FillBnd, bool PutFwdInBwdOnBCs,
2192 bool DoExchange)
2193{
2194 v_GetFwdBwdTracePhys(field, Fwd, Bwd, FillBnd, PutFwdInBwdOnBCs,
2195 DoExchange);
2196}
2199 bool PutFwdInBwdOnBCs)
2200{
2201 v_FillBwdWithBoundCond(Fwd, Bwd, PutFwdInBwdOnBCs);
2202}
2206{
2208}
2212 Array<OneD, NekDouble> &locTraceFwd, Array<OneD, NekDouble> &locTraceBwd)
2213{
2214 v_GetLocTraceFromTracePts(Fwd, Bwd, locTraceFwd, locTraceBwd);
2215}
2217 Array<OneD, NekDouble> &weightjmp)
2218{
2219 v_FillBwdWithBwdWeight(weightave, weightjmp);
2220}
2223{
2224 v_PeriodicBwdCopy(Fwd, Bwd);
2225}
2226inline const std::vector<bool> &ExpList::GetLeftAdjacentFaces(void) const
2227{
2228 return v_GetLeftAdjacentFaces();
2229}
2231{
2232 v_ExtractTracePhys(outarray);
2233}
2235 const Array<OneD, const NekDouble> &inarray,
2236 Array<OneD, NekDouble> &outarray)
2237{
2238 v_ExtractTracePhys(inarray, outarray);
2239}
2242{
2243 return v_GetBndConditions();
2244}
2247{
2248 return v_UpdateBndConditions();
2249}
2251 const std::string varName,
2252 const NekDouble x2_in,
2253 const NekDouble x3_in)
2254{
2255 v_EvaluateBoundaryConditions(time, varName, x2_in, x3_in);
2256}
2258{
2260}
2262 Array<OneD, int> &EdgeID)
2263{
2264 v_GetBoundaryToElmtMap(ElmtID, EdgeID);
2265}
2267 std::shared_ptr<ExpList> &result,
2268 const bool DeclareCoeffPhysArrays)
2269{
2270 v_GetBndElmtExpansion(i, result, DeclareCoeffPhysArrays);
2271}
2272
2274 const Array<OneD, NekDouble> &elmt,
2275 Array<OneD, NekDouble> &boundary)
2276{
2277 v_ExtractElmtToBndPhys(i, elmt, boundary);
2278}
2279
2281 int i, const Array<OneD, const NekDouble> &phys,
2282 Array<OneD, NekDouble> &bndElmt)
2283{
2284 v_ExtractPhysToBndElmt(i, phys, bndElmt);
2285}
2286
2288 const Array<OneD, const NekDouble> &phys,
2290{
2291 v_ExtractPhysToBnd(i, phys, bnd);
2292}
2293
2295 int i, Array<OneD, Array<OneD, NekDouble>> &normals)
2296{
2297 v_GetBoundaryNormals(i, normals);
2298}
2299
2300inline std::vector<bool> &ExpList::GetLeftAdjacentTraces(void)
2301{
2302 return v_GetLeftAdjacentTraces();
2303}
2304
2306
2307} // namespace Nektar::MultiRegions
2308
2309#endif // EXPLIST_H
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:202
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Definition: ErrorUtil.hpp:265
#define MULTI_REGIONS_EXPORT
Base class for all multi-elemental spectral/hp expansions.
Definition: ExpList.h:99
virtual void v_Reset()
Reset geometry information, metrics, matrix managers and geometry information.
Definition: ExpList.cpp:3335
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:1083
Array< OneD, NekDouble > & UpdateCoeffs()
This function returns (a reference to) the array (implemented as m_coeffs) containing all local expa...
Definition: ExpList.h:2095
void AddTraceJacToElmtJac(const Array< OneD, const DNekMatSharedPtr > &FwdMat, const Array< OneD, const DNekMatSharedPtr > &BwdMat, Array< OneD, DNekMatSharedPtr > &fieldMat)
inverse process of v_GetFwdBwdTracePhys. Given Trace integration of Fwd and Bwd Jacobian,...
Definition: ExpList.cpp:6256
void LinearAdvectionReactionSolve(const Array< OneD, Array< OneD, NekDouble > > &velocity, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const NekDouble lambda, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
Solve Advection Diffusion Reaction.
Definition: ExpList.h:1758
void IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function calculates the inner product of a function with respect to all local expansion modes .
Definition: ExpList.h:1650
const Array< OneD, const std::shared_ptr< ExpList > > & GetBndCondExpansions()
Definition: ExpList.h:2112
std::shared_ptr< ExpList > & GetTrace()
Definition: ExpList.h:2153
void WriteTecplotConnectivity(std::ostream &outfile, int expansion=-1)
Definition: ExpList.h:408
int GetNcoeffs(void) const
Returns the total number of local degrees of freedom .
Definition: ExpList.h:1512
virtual void v_ExtractCoeffsToCoeffs(const std::shared_ptr< ExpList > &fromExpList, const Array< OneD, const NekDouble > &fromCoeffs, Array< OneD, NekDouble > &toCoeffs)
Definition: ExpList.cpp:4352
const DNekScalBlkMatSharedPtr & GetBlockMatrix(const GlobalMatrixKey &gkey)
Definition: ExpList.cpp:2652
virtual void v_AddFwdBwdTraceIntegral(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:5074
LibUtilities::NekManager< GlobalLinSysKey, GlobalLinSys > & GetGlobalLinSysManager(void)
Definition: ExpList.cpp:5979
ExpList(const ExpansionType Type=eNoType)
The default constructor using a type.
Definition: ExpList.cpp:91
void SetCoeffs(int i, NekDouble val)
Set the i th coefficiient in m_coeffs to value val.
Definition: ExpList.h:1928
void CreateCollections(Collections::ImplementationType ImpType=Collections::eNoImpType)
Construct collections of elements containing a single element type and polynomial order from the list...
Definition: ExpList.cpp:5787
virtual void v_GetBoundaryToElmtMap(Array< OneD, int > &ElmtID, Array< OneD, int > &EdgeID)
Definition: ExpList.cpp:5680
void ExponentialFilter(Array< OneD, NekDouble > &array, const NekDouble alpha, const NekDouble exponent, const NekDouble cutoff)
Definition: ExpList.cpp:2387
std::shared_ptr< GlobalMatrix > GenGlobalMatrix(const GlobalMatrixKey &mkey, const std::shared_ptr< AssemblyMapCG > &locToGloMap)
Generates a global matrix from the given key and map.
Definition: ExpList.cpp:2786
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:5395
void MultiplyByMassMatrix(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:299
void GlobalEigenSystem(const std::shared_ptr< DNekMat > &Gmat, Array< OneD, NekDouble > &EigValsReal, Array< OneD, NekDouble > &EigValsImag, Array< OneD, NekDouble > &EigVecs=NullNekDouble1DArray)
void FillBwdWithBoundCond(const Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd, bool PutFwdInBwdOnBCs=false)
Definition: ExpList.h:2197
void AddRightIPTPhysDerivBase(const int dir, const Array< OneD, const DNekMatSharedPtr > ElmtJacQuad, Array< OneD, DNekMatSharedPtr > ElmtJacCoef)
Definition: ExpList.cpp:6492
virtual void v_ExtractTracePhys(Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:5160
void ExtractPhysToBndElmt(int i, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &bndElmt)
Definition: ExpList.h:2280
void Reset()
Reset geometry information and reset matrices.
Definition: ExpList.h:387
int GetExpIndex(const Array< OneD, const NekDouble > &gloCoord, NekDouble tol=0.0, bool returnNearestElmt=false, int cachedId=-1, NekDouble maxDistance=1e6)
This function returns the index of the local elemental expansion containing the arbitrary point given...
Definition: ExpList.cpp:3181
virtual GlobalLinSysKey v_HelmSolve(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::ConstFactorMap &factors, const StdRegions::VarCoeffMap &varcoeff, const MultiRegions::VarFactorsMap &varfactors, const Array< OneD, const NekDouble > &dirForcing, const bool PhysSpaceForcing)
Definition: ExpList.cpp:5183
void GetFieldDefinitions(std::vector< LibUtilities::FieldDefinitionsSharedPtr > &fielddef)
Definition: ExpList.h:902
virtual void v_WriteTecplotConnectivity(std::ostream &outfile, int expansion)
Definition: ExpList.cpp:3511
BlockMatrixMapShPtr m_blockMat
Definition: ExpList.h:1128
virtual Array< OneD, const unsigned int > v_GetYIDs(void)
Definition: ExpList.cpp:3973
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
multiply the metric jacobi and quadrature weights
Definition: ExpList.cpp:1925
virtual void v_LinearAdvectionReactionSolve(const Array< OneD, Array< OneD, NekDouble > > &velocity, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const NekDouble lambda, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
Definition: ExpList.cpp:5210
void SetPhys(int i, NekDouble val)
Set the i th value of m_phys to value val.
Definition: ExpList.h:1610
void ExtractTracePhys(Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:2230
void GetBoundaryToElmtMap(Array< OneD, int > &ElmtID, Array< OneD, int > &EdgeID)
Definition: ExpList.h:2261
void PeriodicBwdCopy(const Array< OneD, const NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
Copy and fill the Periodic boundaries.
Definition: ExpList.h:2221
void WriteVtkPieceFooter(std::ostream &outfile, int expansion)
Definition: ExpList.cpp:3771
virtual void v_ExtractPhysToBnd(const int i, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &bnd)
Definition: ExpList.cpp:5599
std::shared_ptr< InterfaceMapDG > & GetInterfaceMap(void)
Definition: ExpList.h:2161
void GetBwdWeight(Array< OneD, NekDouble > &weightAver, Array< OneD, NekDouble > &weightJump)
Get the weight value for boundary conditions for boundary average and jump calculations.
Definition: ExpList.cpp:4389
virtual void v_AddTraceIntegralToOffDiag(const Array< OneD, const NekDouble > &FwdFlux, const Array< OneD, const NekDouble > &BwdFlux, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:6034
virtual std::vector< bool > & v_GetLeftAdjacentTraces(void)
Definition: ExpList.cpp:4630
void InitialiseExpVector(const SpatialDomains::ExpansionInfoMap &expmap)
Define a list of elements using the geometry and basis key information in expmap;.
Definition: ExpList.cpp:1672
void SetupCoeffPhys(bool DeclareCoeffPhysArrays=true, bool SetupOffsets=true)
Definition of the total number of degrees of freedom and quadrature points and offsets to access data...
Definition: ExpList.cpp:1615
virtual LibUtilities::BasisSharedPtr v_GetHomogeneousBasis(void)
Definition: ExpList.h:1464
static SpatialDomains::BoundaryConditionShPtr GetBoundaryCondition(const SpatialDomains::BoundaryConditionCollection &collection, unsigned int index, const std::string &variable)
Definition: ExpList.cpp:5755
const Array< OneD, const NekDouble > & GetCoeffs() const
This function returns (a reference to) the array (implemented as m_coeffs) containing all local expa...
Definition: ExpList.h:1944
virtual void GetBndElmtExpansion(int i, std::shared_ptr< ExpList > &result, const bool DeclareCoeffPhysArrays=true)
Definition: ExpList.h:2266
void SetPhysArray(Array< OneD, NekDouble > &inarray)
Sets the array m_phys.
Definition: ExpList.h:1629
virtual void v_PhysDirectionalDeriv(const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:2357
virtual void v_GetBCValues(Array< OneD, NekDouble > &BndVals, const Array< OneD, NekDouble > &TotField, int BndID)
Definition: ExpList.cpp:5261
void FillBndCondFromField(const Array< OneD, NekDouble > coeffs)
Fill Bnd Condition expansion from the values stored in expansion.
Definition: ExpList.h:1952
const Array< OneD, const std::pair< int, int > > & GetCoeffsToElmt() const
Get m_coeffs to elemental value map.
Definition: ExpList.h:2118
virtual void v_GetLocTraceFromTracePts(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &locTraceFwd, Array< OneD, NekDouble > &locTraceBwd)
Definition: ExpList.cpp:5127
const Array< OneD, const NekDouble > & GetPhys() const
This function returns (a reference to) the array (implemented as m_phys) containing the function ev...
Definition: ExpList.h:2030
std::shared_ptr< ExpList > & UpdateBndCondExpansion(int i)
Definition: ExpList.h:2135
virtual void v_WriteTecplotHeader(std::ostream &outfile, std::string var="")
Definition: ExpList.cpp:3380
virtual std::shared_ptr< ExpList > & v_GetPlane(int n)
Definition: ExpList.cpp:5776
void NormVectorIProductWRTBase(Array< OneD, const NekDouble > &V1, Array< OneD, const NekDouble > &V2, Array< OneD, NekDouble > &outarray, int BndID)
Definition: ExpList.h:1887
int GetCoeff_Offset(int n) const
Get the start offset position for a local contiguous list of coeffs correspoinding to element n.
Definition: ExpList.h:2078
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:2038
virtual void v_AddTraceQuadPhysToOffDiag(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &field)
Definition: ExpList.cpp:5118
std::shared_ptr< LibUtilities::SessionReader > GetSession() const
Returns the session object.
Definition: ExpList.h:961
virtual void v_ClearGlobalLinSysManager(void)
Definition: ExpList.cpp:3981
virtual void v_UnsetGlobalLinSys(GlobalLinSysKey, bool)
Definition: ExpList.cpp:3993
void GeneralGetFieldDefinitions(std::vector< LibUtilities::FieldDefinitionsSharedPtr > &fielddef, int NumHomoDir=0, Array< OneD, LibUtilities::BasisSharedPtr > &HomoBasis=LibUtilities::NullBasisSharedPtr1DArray, std::vector< NekDouble > &HomoLen=LibUtilities::NullNekDoubleVector, bool homoStrips=false, std::vector< unsigned int > &HomoSIDs=LibUtilities::NullUnsignedIntVector, std::vector< unsigned int > &HomoZIDs=LibUtilities::NullUnsignedIntVector, std::vector< unsigned int > &HomoYIDs=LibUtilities::NullUnsignedIntVector)
Definition: ExpList.cpp:4081
virtual void v_GetPeriodicEntities(PeriodicMap &periodicVerts, PeriodicMap &periodicEdges, PeriodicMap &periodicFaces)
Definition: ExpList.cpp:5747
void IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function calculates the inner product of a function with respect to the derivative (in directio...
Definition: ExpList.h:1659
void SetCoeffsArray(Array< OneD, NekDouble > &inarray)
Set the m_coeffs array to inarray.
Definition: ExpList.h:1932
virtual void v_CurlCurl(Array< OneD, Array< OneD, NekDouble > > &Vel, Array< OneD, Array< OneD, NekDouble > > &Q)
Definition: ExpList.cpp:2289
LibUtilities::TranspositionSharedPtr GetTransposition(void)
This function returns the transposition class associated with the homogeneous expansion.
Definition: ExpList.h:597
virtual void v_LocalToGlobal(bool UseComm)
Definition: ExpList.cpp:5349
void SetBndCondBwdWeight(const int index, const NekDouble value)
Set the weight value for boundary conditions.
Definition: ExpList.h:2131
virtual const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > & v_GetBndConditions()
Definition: ExpList.cpp:5704
virtual LibUtilities::NekManager< GlobalLinSysKey, GlobalLinSys > & v_GetGlobalLinSysManager(void)
Definition: ExpList.cpp:4001
void SetHomoLen(const NekDouble lhom)
This function sets the Width of homogeneous direction associated with the homogeneous expansion.
Definition: ExpList.h:611
virtual void v_DealiasedDotProd(const int num_dofs, const Array< OneD, Array< OneD, NekDouble > > &inarray1, const Array< OneD, Array< OneD, NekDouble > > &inarray2, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: ExpList.cpp:5251
int GetNpoints(void) const
Returns the total number of quadrature points m_npoints .
Definition: ExpList.h:1588
SpatialDomains::MeshGraphSharedPtr GetGraph()
Definition: ExpList.h:970
virtual ~ExpList()
The default destructor.
Definition: ExpList.cpp:1893
void GetFwdBwdTracePhys(Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
Definition: ExpList.h:2184
int GetElmtToExpId(int elmtId)
This function returns the index inside m_exp for a given geom id.
Definition: ExpList.h:1039
const Array< OneD, const int > & GetTraceBndMap(void)
Definition: ExpList.h:2165
virtual void v_FwdTransLocalElmt(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:2446
int GetPoolCount(std::string)
Definition: ExpList.cpp:5968
virtual void v_ImposeDirichletConditions(Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:5323
void FwdTransLocalElmt(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function elementally evaluates the forward transformation of a function onto the global spectra...
Definition: ExpList.h:1687
Array< OneD, const unsigned int > GetZIDs(void)
This function returns a vector containing the wave numbers in z-direction associated with the 3D homo...
Definition: ExpList.h:590
virtual void v_EvaluateBoundaryConditions(const NekDouble time=0.0, const std::string varName="", const NekDouble x2_in=NekConstants::kNekUnsetDouble, const NekDouble x3_in=NekConstants::kNekUnsetDouble)
Definition: ExpList.cpp:5725
void FwdTransBndConstrained(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:1696
virtual size_t v_GetNumElmts(void)
Definition: ExpList.h:1162
virtual const std::vector< bool > & v_GetLeftAdjacentFaces(void) const
Definition: ExpList.cpp:5152
Array< OneD, int > m_coeff_offset
Offset of elemental data into the array m_coeffs.
Definition: ExpList.h:1123
int EvalBasisNumModesMax(void) const
Evaulates the maximum number of modes in the elemental basis order over all elements.
Definition: ExpList.h:1526
void AddFwdBwdTraceIntegral(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:2178
std::shared_ptr< AssemblyMapDG > & GetTraceMap(void)
Definition: ExpList.h:2157
void WriteTecplotField(std::ostream &outfile, int expansion=-1)
Definition: ExpList.h:404
virtual void v_GetMovingFrames(const SpatialDomains::GeomMMF MMFdir, const Array< OneD, const NekDouble > &CircCentre, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: ExpList.cpp:4402
size_t GetNumElmts(void)
This function returns the number of elements in the expansion which may be different for a homogeoeno...
Definition: ExpList.h:653
virtual std::shared_ptr< AssemblyMapDG > & v_GetTraceMap()
Definition: ExpList.cpp:4609
Array< OneD, std::pair< int, int > > m_coeffsToElmt
m_coeffs to elemental value map
Definition: ExpList.h:1127
virtual void v_GetFwdBwdTracePhys(Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
Definition: ExpList.cpp:5083
NekDouble PhysEvaluate(const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &phys)
Definition: ExpList.cpp:3300
const std::unordered_map< int, int > & GetElmtToExpId(void)
This function returns the map of index inside m_exp to geom id.
Definition: ExpList.h:1032
virtual GlobalLinSysKey v_LinearAdvectionDiffusionReactionSolve(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::ConstFactorMap &factors, const StdRegions::VarCoeffMap &varcoeff, const MultiRegions::VarFactorsMap &varfactors, const Array< OneD, const NekDouble > &dirForcing, const bool PhysSpaceForcing)
Definition: ExpList.cpp:5196
virtual NekDouble v_Integral(const Array< OneD, const NekDouble > &inarray)
Definition: ExpList.cpp:3901
virtual void v_ExtractDataToCoeffs(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, std::string &field, Array< OneD, NekDouble > &coeffs, std::unordered_map< int, int > zIdToPlane)
Extract data from raw field data into expansion list.
Definition: ExpList.cpp:4264
void PhysDirectionalDeriv(const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:1822
virtual std::map< int, RobinBCInfoSharedPtr > v_GetRobinBCInfo(void)
Definition: ExpList.cpp:5737
virtual const Array< OneD, const int > & v_GetTraceBndMap()
Definition: ExpList.cpp:4625
virtual void v_GetBoundaryNormals(int i, Array< OneD, Array< OneD, NekDouble > > &normals)
Definition: ExpList.cpp:5636
void GetBoundaryNormals(int i, Array< OneD, Array< OneD, NekDouble > > &normals)
Definition: ExpList.h:2294
virtual Array< OneD, const NekDouble > v_HomogeneousEnergy(void)
Definition: ExpList.cpp:3935
virtual NekDouble v_GetHomoLen(void)
Definition: ExpList.cpp:3951
void ExtractFileBCs(const std::string &fileName, LibUtilities::CommSharedPtr comm, const std::string &varName, const std::shared_ptr< ExpList > locExpList)
void GetMovingFrames(const SpatialDomains::GeomMMF MMFdir, const Array< OneD, const NekDouble > &CircCentre, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: ExpList.h:1787
const LocTraceToTraceMapSharedPtr & GetLocTraceToTraceMap() const
Definition: ExpList.h:2123
virtual Array< OneD, const unsigned int > v_GetZIDs(void)
Definition: ExpList.cpp:3965
Array< OneD, const NekDouble > HomogeneousEnergy(void)
This function calculates the energy associated with each one of the modesof a 3D homogeneous nD expan...
Definition: ExpList.h:574
void SetHomo1DSpecVanVisc(Array< OneD, NekDouble > visc)
This function sets the Spectral Vanishing Viscosity in homogeneous1D expansion.
Definition: ExpList.h:581
void ExtractCoeffsToCoeffs(const std::shared_ptr< ExpList > &fromExpList, const Array< OneD, const NekDouble > &fromCoeffs, Array< OneD, NekDouble > &toCoeffs)
Extract the data from fromField using fromExpList the coeffs using the basic ExpList Elemental expans...
Definition: ExpList.cpp:4248
void GetCoords(Array< OneD, NekDouble > &coord_0, Array< OneD, NekDouble > &coord_1=NullNekDouble1DArray, Array< OneD, NekDouble > &coord_2=NullNekDouble1DArray)
This function calculates the coordinates of all the elemental quadrature points .
Definition: ExpList.h:1770
void ResetMatrices()
Reset matrices.
Definition: ExpList.cpp:3363
const Array< OneD, int > EvalBasisNumModesMaxPerExp(void) const
Returns the vector of the number of modes in the elemental basis order over all elements.
Definition: ExpList.h:1541
std::shared_ptr< GlobalLinSys > GenGlobalBndLinSys(const GlobalLinSysKey &mkey, const AssemblyMapSharedPtr &locToGloMap)
Generate a GlobalLinSys from information provided by the key "mkey" and the mapping provided in LocTo...
Definition: ExpList.cpp:3089
virtual void v_SetBndCondBwdWeight(const int index, const NekDouble value)
Definition: ExpList.cpp:5145
void UnsetGlobalLinSys(GlobalLinSysKey, bool)
Definition: ExpList.cpp:5973
void DealiasedProd(const int num_dofs, const Array< OneD, NekDouble > &inarray1, const Array< OneD, NekDouble > &inarray2, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:1858
void GetNormals(Array< OneD, Array< OneD, NekDouble > > &normals)
Definition: ExpList.h:2169
void ExtractPhysToBnd(int i, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &bnd)
Definition: ExpList.h:2287
virtual void v_SmoothField(Array< OneD, NekDouble > &field)
Definition: ExpList.cpp:2478
virtual void v_AppendFieldData(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata)
Definition: ExpList.cpp:4208
std::vector< bool > m_collectionsDoInit
Vector of bools to act as an initialise on first call flag.
Definition: ExpList.h:1121
NekDouble GetCoeff(int i)
Get the i th value (coefficient) of m_coeffs.
Definition: ExpList.h:2010
bool m_physState
The state of the array m_phys.
Definition: ExpList.h:1107
virtual void v_HomogeneousFwdTrans(const int npts, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
Definition: ExpList.cpp:5221
void ExtractCoeffsFromFile(const std::string &fileName, LibUtilities::CommSharedPtr comm, const std::string &varName, Array< OneD, NekDouble > &coeffs)
Definition: ExpList.cpp:4008
NekDouble Integral()
Calculates the error of the global spectral/hp element approximation.
Definition: ExpList.h:543
std::shared_ptr< ExpList > & GetPlane(int n)
Definition: ExpList.h:979
virtual void v_DealiasedProd(const int num_dofs, const Array< OneD, NekDouble > &inarray1, const Array< OneD, NekDouble > &inarray2, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:5241
void AddRightIPTBaseMatrix(const Array< OneD, const DNekMatSharedPtr > ElmtJacQuad, Array< OneD, DNekMatSharedPtr > ElmtJacCoef)
Definition: ExpList.cpp:6549
void GetLocTraceFromTracePts(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &locTraceFwd, Array< OneD, NekDouble > &locTraceBwd)
Definition: ExpList.h:2209
virtual void v_WriteVtkPieceData(std::ostream &outfile, int expansion, std::string var)
Definition: ExpList.cpp:3778
virtual void v_WriteVtkPieceHeader(std::ostream &outfile, int expansion, int istrip)
Definition: ExpList.cpp:3638
std::shared_ptr< ExpList > GetSharedThisPtr()
Returns a shared pointer to the current object.
Definition: ExpList.h:956
int GetShapeDimension()
This function returns the dimension of the shape of the element eid.
Definition: ExpList.h:1912
void ExtractElmtDataToCoeffs(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, std::string &field, Array< OneD, NekDouble > &coeffs)
Extract the data in fielddata into the coeffs using the basic ExpList Elemental expansions rather tha...
NekDouble H1(const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray)
Calculates the error of the global spectral/hp element approximation.
Definition: ExpList.cpp:4063
LibUtilities::CommSharedPtr m_comm
Communicator.
Definition: ExpList.h:1056
std::shared_ptr< GlobalLinSys > GenGlobalLinSys(const GlobalLinSysKey &mkey, const std::shared_ptr< AssemblyMapCG > &locToGloMap)
This operation constructs the global linear system of type mkey.
Definition: ExpList.cpp:3071
void WriteVtkHeader(std::ostream &outfile)
Definition: ExpList.cpp:3624
bool GetPhysState(void) const
This function indicates whether the array of physical values (implemented as m_phys) is filled or no...
Definition: ExpList.h:1643
std::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:1118
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Definition: ExpList.h:1063
void GenerateElementVector(const int ElementID, const NekDouble scalar1, const NekDouble scalar2, Array< OneD, NekDouble > &outarray)
Generate vector v such that v[i] = scalar1 if i is in the element < ElementID. Otherwise,...
Definition: ExpList.cpp:4446
NekDouble VectorFlux(const Array< OneD, Array< OneD, NekDouble > > &inarray)
Definition: ExpList.h:568
void AddTraceIntegralToOffDiag(const Array< OneD, const NekDouble > &FwdFlux, const Array< OneD, const NekDouble > &BwdFlux, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:1011
std::shared_ptr< DNekMat > GenGlobalMatrixFull(const GlobalLinSysKey &mkey, const std::shared_ptr< AssemblyMapCG > &locToGloMap)
Definition: ExpList.cpp:2919
const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > & GetBndConditions()
Definition: ExpList.h:2241
LocalRegions::ExpansionSharedPtr & GetExpFromGeomId(int n)
This function returns (a shared pointer to) the local elemental expansion of the element given a glo...
Definition: ExpList.h:2058
void SetPhysState(const bool physState)
This function manually sets whether the array of physical values (implemented as m_phys) is filled o...
Definition: ExpList.h:1636
void PhysGalerkinProjection1DScaled(const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function Galerkin projects the physical space points in inarray to outarray where inarray is ass...
Definition: ExpList.h:640
virtual void v_NormVectorIProductWRTBase(Array< OneD, const NekDouble > &V1, Array< OneD, const NekDouble > &V2, Array< OneD, NekDouble > &outarray, int BndID)
Definition: ExpList.cpp:5270
virtual void v_ExtractElmtToBndPhys(const int i, const Array< OneD, NekDouble > &elmt, Array< OneD, NekDouble > &boundary)
Definition: ExpList.cpp:5522
virtual void v_ExtractPhysToBndElmt(const int i, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &bndElmt)
Definition: ExpList.cpp:5560
virtual std::shared_ptr< InterfaceMapDG > & v_GetInterfaceMap()
Definition: ExpList.cpp:4617
void CurlCurl(Array< OneD, Array< OneD, NekDouble > > &Vel, Array< OneD, Array< OneD, NekDouble > > &Q)
Definition: ExpList.h:1832
std::map< int, RobinBCInfoSharedPtr > GetRobinBCInfo()
Definition: ExpList.h:885
GlobalLinSysKey HelmSolve(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::ConstFactorMap &factors, const StdRegions::VarCoeffMap &varcoeff=StdRegions::NullVarCoeffMap, const MultiRegions::VarFactorsMap &varfactors=MultiRegions::NullVarFactorsMap, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray, const bool PhysSpaceForcing=true)
Solve helmholtz problem.
Definition: ExpList.h:1734
virtual void v_Upwind(const Array< OneD, const Array< OneD, NekDouble > > &Vec, const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &Upwind)
Definition: ExpList.cpp:4503
void EvaluateBoundaryConditions(const NekDouble time=0.0, const std::string varName="", const NekDouble=NekConstants::kNekUnsetDouble, const NekDouble=NekConstants::kNekUnsetDouble)
Definition: ExpList.h:2250
std::vector< LibUtilities::FieldDefinitionsSharedPtr > GetFieldDefinitions()
Definition: ExpList.h:897
virtual void v_FillBwdWithBwdWeight(Array< OneD, NekDouble > &weightave, Array< OneD, NekDouble > &weightjmp)
Definition: ExpList.cpp:5687
Array< OneD, NekDouble > & UpdatePhys()
This function returns (a reference to) the array (implemented as m_phys) containing the function ev...
Definition: ExpList.h:2105
virtual void v_Curl(Array< OneD, Array< OneD, NekDouble > > &Vel, Array< OneD, Array< OneD, NekDouble > > &Q)
Definition: ExpList.cpp:2240
virtual void v_PhysGalerkinProjection1DScaled(const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:6606
void GeneralMatrixOp(const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function calculates the result of the multiplication of a matrix of type specified by mkey with ...
Definition: ExpList.cpp:2696
virtual void v_SetUpPhysNormals()
: Set up a normal along the trace elements between two elements at elemental level
Definition: ExpList.cpp:5499
virtual void v_GlobalToLocal(void)
Definition: ExpList.cpp:5364
virtual LibUtilities::TranspositionSharedPtr v_GetTransposition(void)
Definition: ExpList.cpp:3943
void ExtractElmtToBndPhys(int i, const Array< OneD, NekDouble > &elmt, Array< OneD, NekDouble > &boundary)
Definition: ExpList.h:2273
void HomogeneousBwdTrans(const int npts, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
Definition: ExpList.h:1849
void SetWaveSpace(const bool wavespace)
Sets the wave space to the one of the possible configuration true or false.
Definition: ExpList.h:1596
const std::shared_ptr< LocalRegions::ExpansionVector > GetExp() const
This function returns the vector of elements in the expansion.
Definition: ExpList.h:2070
SpatialDomains::MeshGraphSharedPtr m_graph
Mesh associated with this expansion list.
Definition: ExpList.h:1060
virtual void v_MultiplyByInvMassMatrix(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:5175
virtual void v_AddTraceQuadPhysToField(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &field)
Definition: ExpList.cpp:5109
virtual const Array< OneD, const NekDouble > & v_GetBndCondBwdWeight()
Definition: ExpList.cpp:5137
void GlobalToLocal(void)
Scatters from the global coefficients to the local coefficients .
Definition: ExpList.h:1971
virtual NekDouble v_L2(const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray)
Definition: ExpList.cpp:3856
const Array< OneD, const NekDouble > & GetBndCondBwdWeight()
Get the weight value for boundary conditions.
Definition: ExpList.h:2127
void WriteVtkFooter(std::ostream &outfile)
Definition: ExpList.cpp:3632
void MultiplyByBlockMatrix(const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:1905
virtual std::shared_ptr< ExpList > & v_GetTrace()
Definition: ExpList.cpp:4601
void GetElmtNormalLength(Array< OneD, NekDouble > &lengthsFwd, Array< OneD, NekDouble > &lengthsBwd)
Get the length of elements in boundary normal direction.
Definition: ExpList.cpp:4919
int Get1DScaledTotPoints(const NekDouble scale) const
Returns the total number of qudature points scaled by the factor scale on each 1D direction.
Definition: ExpList.h:1567
GlobalLinSysKey LinearAdvectionDiffusionReactionSolve(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::ConstFactorMap &factors, const StdRegions::VarCoeffMap &varcoeff=StdRegions::NullVarCoeffMap, const MultiRegions::VarFactorsMap &varfactors=MultiRegions::NullVarFactorsMap, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray, const bool PhysSpaceForcing=true)
Solve Advection Diffusion Reaction.
Definition: ExpList.h:1747
void PhysInterp1DScaled(const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function interpolates the physical space points in inarray to outarray using the same points def...
Definition: ExpList.h:629
virtual NekDouble v_VectorFlux(const Array< OneD, Array< OneD, NekDouble > > &inarray)
Definition: ExpList.cpp:3915
NekDouble L2(const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray)
This function calculates the error of the global This function calculates the error with respect to...
Definition: ExpList.h:514
virtual std::vector< LibUtilities::FieldDefinitionsSharedPtr > v_GetFieldDefinitions(void)
Definition: ExpList.cpp:4193
void IProductWRTDirectionalDerivBase(const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Directional derivative along a given direction.
Definition: ExpList.cpp:1988
void WriteVtkPieceData(std::ostream &outfile, int expansion, std::string var="v")
Definition: ExpList.h:421
void FillBwdWithBwdWeight(Array< OneD, NekDouble > &weightave, Array< OneD, NekDouble > &weightjmp)
Fill Bwd with boundary conditions.
Definition: ExpList.h:2216
virtual void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:3125
virtual Array< OneD, SpatialDomains::BoundaryConditionShPtr > & v_UpdateBndConditions()
Definition: ExpList.cpp:5715
void SetCoeff(int i, NekDouble val)
Set the i th coefficiient in m_coeffs to value val.
Definition: ExpList.h:1920
void GetMatIpwrtDeriveBase(const Array< OneD, const Array< OneD, NekDouble > > &inarray, const int nDirctn, Array< OneD, DNekMatSharedPtr > &mtxPerVar)
Definition: ExpList.cpp:6042
void LocalToGlobal(bool useComm=true)
Gathers the global coefficients from the local coefficients .
Definition: ExpList.h:1961
void GetPeriodicEntities(PeriodicMap &periodicVerts, PeriodicMap &periodicEdges, PeriodicMap &periodicFaces=NullPeriodicMap)
Definition: ExpList.h:890
NekDouble Integral(const Array< OneD, const NekDouble > &inarray)
Definition: ExpList.h:564
virtual void v_GetCoords(Array< OneD, NekDouble > &coord_0, Array< OneD, NekDouble > &coord_1=NullNekDouble1DArray, Array< OneD, NekDouble > &coord_2=NullNekDouble1DArray)
Definition: ExpList.cpp:5439
virtual void v_WriteTecplotField(std::ostream &outfile, int expansion)
Definition: ExpList.cpp:3591
virtual void v_FillBndCondFromField(const Array< OneD, NekDouble > coeffs)
Definition: ExpList.cpp:5332
Collections::CollectionVector m_collections
Definition: ExpList.h:1119
void AddTraceQuadPhysToField(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &field)
Add Fwd and Bwd value to field, a reverse procedure of GetFwdBwdTracePhys.
Definition: ExpList.h:2203
virtual int v_GetPoolCount(std::string)
Definition: ExpList.cpp:3987
virtual void v_PhysInterp1DScaled(const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:5984
virtual void v_HomogeneousBwdTrans(const int npts, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
Definition: ExpList.cpp:5231
void ApplyGeomInfo()
Apply geometry information to each expansion.
Definition: ExpList.cpp:3323
virtual void v_PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2)
Definition: ExpList.cpp:2142
virtual void v_FwdTransBndConstrained(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:2455
void AddTraceIntegral(const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:2173
std::unordered_map< int, int > m_elmtToExpId
Mapping from geometry ID of element to index inside m_exp.
Definition: ExpList.h:1135
void FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:1679
void SmoothField(Array< OneD, NekDouble > &field)
Smooth a field across elements.
Definition: ExpList.h:1705
virtual void v_WriteTecplotZone(std::ostream &outfile, int expansion)
Definition: ExpList.cpp:3418
Array< OneD, const unsigned int > GetYIDs(void)
This function returns a vector containing the wave numbers in y-direction associated with the 3D homo...
Definition: ExpList.h:620
LibUtilities::SessionReaderSharedPtr m_session
Session.
Definition: ExpList.h:1058
virtual void v_PeriodicBwdCopy(const Array< OneD, const NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
Definition: ExpList.cpp:5694
void MultiplyByInvMassMatrix(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:1725
void AppendFieldData(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, Array< OneD, NekDouble > &coeffs)
Append the data in coeffs listed in elements fielddef->m_ElementIDs onto fielddata.
Definition: ExpList.h:917
virtual void v_FillBwdWithBoundCond(const Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd, bool PutFwdInBwdOnBCs)
Definition: ExpList.cpp:5100
void DivideByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Divided by the metric jacobi and quadrature weights.
Definition: ExpList.cpp:1942
virtual const Array< OneD, const std::shared_ptr< ExpList > > & v_GetBndCondExpansions(void)
Definition: ExpList.cpp:4476
void ImposeDirichletConditions(Array< OneD, NekDouble > &outarray)
Impose Dirichlet Boundary Conditions onto Array.
Definition: ExpList.h:1948
Array< OneD, int > m_phys_offset
Offset of elemental data into the array m_phys.
Definition: ExpList.h:1125
virtual void v_SetHomoLen(const NekDouble lhom)
Definition: ExpList.cpp:3959
void Upwind(const Array< OneD, const NekDouble > &Vn, const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &Upwind)
Definition: ExpList.h:2146
virtual void v_GetNormals(Array< OneD, Array< OneD, NekDouble > > &normals)
Populate normals with the normals of all expansions.
Definition: ExpList.cpp:4710
Array< OneD, SpatialDomains::BoundaryConditionShPtr > & UpdateBndConditions()
Definition: ExpList.h:2246
const DNekScalBlkMatSharedPtr GenBlockMatrix(const GlobalMatrixKey &gkey)
This function assembles the block diagonal matrix of local matrices of the type mtype.
Definition: ExpList.cpp:2517
void PhysDeriv(Direction edir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d)
Definition: ExpList.h:1813
void MultiplyByElmtInvMass(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function elementally mulplies the coefficient space of Sin my the elemental inverse of the mass ...
Definition: ExpList.cpp:2408
void BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function elementally evaluates the backward transformation of the global spectral/hp element exp...
Definition: ExpList.h:1717
NekDouble GetHomoLen(void)
This function returns the Width of homogeneous direction associated with the homogeneous expansion.
Definition: ExpList.h:604
virtual void v_GetBndElmtExpansion(int i, std::shared_ptr< ExpList > &result, const bool DeclareCoeffPhysArrays)
Definition: ExpList.cpp:5512
virtual void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:5378
bool GetWaveSpace(void) const
This function returns the third direction expansion condition, which can be in wave space (coefficien...
Definition: ExpList.h:1604
virtual std::shared_ptr< ExpList > & v_UpdateBndCondExpansion(int i)
Definition: ExpList.cpp:4484
int GetPhys_Offset(int n) const
Get the start offset position for a local contiguous list of quadrature points in a full array corres...
Definition: ExpList.h:2085
virtual const std::shared_ptr< LocTraceToTraceMap > & v_GetLocTraceToTraceMap(void) const
Definition: ExpList.cpp:6681
virtual void v_IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:1969
void ExtractDataToCoeffs(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, std::string &field, Array< OneD, NekDouble > &coeffs, std::unordered_map< int, int > zIdToPlane=std::unordered_map< int, int >())
Extract the data in fielddata into the coeffs.
Definition: ExpList.cpp:4240
void SetModifiedBasis(const bool modbasis)
Set Modified Basis for the stability analysis.
ExpansionType m_expType
Expansion type.
Definition: ExpList.h:1051
void WriteTecplotHeader(std::ostream &outfile, std::string var="")
Definition: ExpList.h:396
const std::vector< bool > & GetLeftAdjacentFaces(void) const
Definition: ExpList.h:2226
void DealiasedDotProd(const int num_dofs, const Array< OneD, Array< OneD, NekDouble > > &inarray1, const Array< OneD, Array< OneD, NekDouble > > &inarray2, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: ExpList.h:1868
void AppendFieldData(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata)
Append the element data listed in elements fielddef->m_ElementIDs onto fielddata.
Definition: ExpList.h:910
std::shared_ptr< LibUtilities::Comm > GetComm() const
Returns the comm object.
Definition: ExpList.h:966
virtual void v_AddTraceIntegral(const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:5066
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
Definition: ExpList.h:1099
std::vector< bool > & GetLeftAdjacentTraces(void)
Definition: ExpList.h:2300
ExpansionType GetExpType(void)
Returns the type of the expansion.
Definition: ExpList.cpp:1888
void AddTraceQuadPhysToOffDiag(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &field)
Definition: ExpList.h:815
NekDouble Linf(const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray)
This function calculates the error of the global spectral/hp element approximation.
Definition: ExpList.cpp:3818
void HomogeneousFwdTrans(const int npts, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
Definition: ExpList.h:1840
void WriteVtkPieceHeader(std::ostream &outfile, int expansion, int istrip=0)
Definition: ExpList.h:414
int GetTotPoints(void) const
Returns the total number of quadrature points m_npoints .
Definition: ExpList.h:1557
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
LibUtilities::BasisSharedPtr GetHomogeneousBasis(void)
Definition: ExpList.h:975
void WriteTecplotZone(std::ostream &outfile, int expansion=-1)
Definition: ExpList.h:400
ExpList(const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, const ExpListSharedPtr > &bndConstraint, const Array< OneD, const SpatialDomains ::BoundaryConditionShPtr > &bndCond, const LocalRegions::ExpansionVector &locexp, const SpatialDomains::MeshGraphSharedPtr &graph, const LibUtilities::CommSharedPtr &comm, const bool DeclareCoeffPhysArrays=true, const std::string variable="DefaultVar", const Collections::ImplementationType ImpType=Collections::eNoImpType)
Generate expansions for the trace space expansions used in DisContField.
virtual void v_SetHomo1DSpecVanVisc(Array< OneD, NekDouble > visc)
Definition: ExpList.h:1472
void GetBCValues(Array< OneD, NekDouble > &BndVals, const Array< OneD, NekDouble > &TotField, int BndID)
Definition: ExpList.h:1878
int GetCoordim(int eid)
This function returns the dimension of the coordinates of the element eid.
Definition: ExpList.h:1903
void GetDiagMatIpwrtBase(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, DNekMatSharedPtr > &mtxPerVar)
Definition: ExpList.cpp:6185
A global linear system.
Definition: GlobalLinSys.h:70
Describes a matrix with ordering defined by a local to global map.
std::vector< Collection > CollectionVector
Definition: Collection.h:128
static BasisSharedPtr NullBasisSharedPtr
Definition: Basis.h:350
static std::vector< unsigned int > NullUnsignedIntVector
std::shared_ptr< Basis > BasisSharedPtr
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< FieldDefinitions > FieldDefinitionsSharedPtr
Definition: FieldIO.h:184
static std::vector< NekDouble > NullNekDoubleVector
std::shared_ptr< Transposition > TranspositionSharedPtr
static Array< OneD, BasisSharedPtr > NullBasisSharedPtr1DArray
Definition: Basis.h:351
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:55
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:66
std::vector< ExpansionSharedPtr > ExpansionVector
Definition: Expansion.h:68
static ExpList NullExpList
An empty ExpList object.
Definition: ExpList.h:1498
static PeriodicMap NullPeriodicMap
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:87
std::shared_ptr< BlockMatrixMap > BlockMatrixMapShPtr
A shared pointer to a BlockMatrixMap.
Definition: ExpList.h:93
static VarFactorsMap NullVarFactorsMap
static ExpListSharedPtr NullExpListSharedPtr
Definition: ExpList.h:1499
static LibUtilities::NekManager< GlobalLinSysKey, GlobalLinSys > NullGlobalLinSysManager
Definition: ExpList.h:1503
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
static GlobalLinSysKey NullGlobalLinSysKey(StdRegions::eNoMatrixType)
std::shared_ptr< LocTraceToTraceMap > LocTraceToTraceMapSharedPtr
static const Array< OneD, ExpListSharedPtr > NullExpListSharedPtrArray
Definition: ExpList.h:2305
std::map< StdRegions::ConstFactorType, Array< OneD, NekDouble > > VarFactorsMap
std::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
Definition: AssemblyMap.h:50
std::map< int, std::vector< PeriodicEntity > > PeriodicMap
std::map< GlobalMatrixKey, DNekScalBlkMatSharedPtr > BlockMatrixMap
A map between global matrix keys and their associated block matrices.
Definition: ExpList.h:91
static const NekDouble kNekUnsetDouble
GeomMMF
Principle direction for MMF.
std::shared_ptr< BoundaryConditionBase > BoundaryConditionShPtr
Definition: Conditions.h:212
std::map< int, BoundaryConditionMapShPtr > BoundaryConditionCollection
Definition: Conditions.h:219
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
std::shared_ptr< PointGeom > PointGeomSharedPtr
Definition: Geometry.h:57
std::map< int, ExpansionInfoShPtr > ExpansionInfoMap
Definition: MeshGraph.h:141
std::map< int, CompositeSharedPtr > CompositeMap
Definition: MeshGraph.h:136
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:430
static VarCoeffMap NullVarCoeffMap
Definition: StdRegions.hpp:376
std::map< StdRegions::VarCoeffType, VarCoeffEntry > VarCoeffMap
Definition: StdRegions.hpp:375
StdRegions::ConstFactorMap factors
std::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:79
static Array< OneD, NekDouble > NullNekDouble1DArray
double NekDouble
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825