Nektar++
Loading...
Searching...
No Matches
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
74
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
138
139 /// Generate expansions for the trace space expansions used in
140 /// DisContField.
143 const Array<OneD, const ExpListSharedPtr> &bndConstraint,
145 &bndCond,
146 const LocalRegions::ExpansionVector &locexp,
148 const LibUtilities::CommSharedPtr &comm,
149 const bool DeclareCoeffPhysArrays = true,
150 const std::string variable = "DefaultVar",
151 const Collections::ImplementationType ImpType =
153
154 /// Generate an trace ExpList from a meshgraph \a graph and session file
157 const LocalRegions::ExpansionVector &locexp,
159 const bool DeclareCoeffPhysArrays, const std::string variable,
160 const Collections::ImplementationType ImpType =
162
163 /// Constructor based on domain information only for 1D &
164 /// 2D boundary conditions
167 const SpatialDomains::CompositeMap &domain,
169 const bool DeclareCoeffPhysArrays = true,
170 const std::string variable = "DefaultVar",
171 bool SetToOneSpaceDimension = false,
173 const Collections::ImplementationType ImpType =
175
176 /// The default destructor.
178
179 /// Returns the total number of local degrees of freedom
180 /// \f$N_{\mathrm{eof}}=\sum_{e=1}^{{N_{\mathrm{el}}}}N^{e}_m\f$.
181 inline int GetNcoeffs(void) const;
182
183 /// Returns the total number of local degrees of freedom
184 /// for element eid
185 MULTI_REGIONS_EXPORT int GetNcoeffs(const int eid) const;
186
187 /// Returns the type of the expansion
189
190 /// Returns the type of the expansion
192
193 /// Evaulates the maximum number of modes in the elemental basis
194 /// order over all elements
195 inline int EvalBasisNumModesMax(void) const;
196
197 /// Returns the vector of the number of modes in the elemental
198 /// basis order over all elements.
200 void) const;
201
202 /// Returns the total number of quadrature points #m_npoints
203 /// \f$=Q_{\mathrm{tot}}\f$.
204 inline int GetTotPoints(void) const;
205
206 /// Returns the total number of quadrature points for eid's element
207 /// \f$=Q_{\mathrm{tot}}\f$.
208 inline int GetTotPoints(const int eid) const;
209
210 /// Returns the total number of quadrature points #m_npoints
211 /// \f$=Q_{\mathrm{tot}}\f$.
212 inline int GetNpoints(void) const;
213
214 /// Returns the total number of qudature points scaled by
215 /// the factor scale on each 1D direction
216 inline int Get1DScaledTotPoints(const NekDouble scale) const;
217
218 /// Sets the wave space to the one of the possible configuration
219 /// true or false
220 inline void SetWaveSpace(const bool wavespace);
221
222 /// Set Modified Basis for the stability analysis
223 inline void SetModifiedBasis(const bool modbasis);
224
225 /// This function returns the third direction expansion condition,
226 /// which can be in wave space (coefficient) or not
227 /// It is stored in the variable m_WaveSpace.
228 inline bool GetWaveSpace(void) const;
229
230 /// Set the \a i th value of \a m_phys to value \a val
231 inline void SetPhys(int i, NekDouble val);
232 /// Fills the array #m_phys
233 inline void SetPhys(const Array<OneD, const NekDouble> &inarray);
234 /// Sets the array #m_phys
235 inline void SetPhysArray(Array<OneD, NekDouble> &inarray);
236 /// This function manually sets whether the array of physical
237 /// values \f$\boldsymbol{u}_l\f$ (implemented as #m_phys) is
238 /// filled or not.
239 inline void SetPhysState(const bool physState);
240 /// This function indicates whether the array of physical values
241 /// \f$\boldsymbol{u}_l\f$ (implemented as #m_phys) is filled or
242 /// not.
243 inline bool GetPhysState(void) const;
244
245 /// multiply the metric jacobi and quadrature weights
247 const Array<OneD, const NekDouble> &inarray,
248 Array<OneD, NekDouble> &outarray);
249 /// Divided by the metric jacobi and quadrature weights
251 const Array<OneD, const NekDouble> &inarray,
252 Array<OneD, NekDouble> &outarray);
253 /// This function calculates the inner product of a function
254 /// \f$f(\boldsymbol{x})\f$ with respect to all \em local
255 /// expansion modes \f$\phi_n^e(\boldsymbol{x})\f$.
256 inline void IProductWRTBase(const Array<OneD, const NekDouble> &inarray,
257 Array<OneD, NekDouble> &outarray);
258 /// This function calculates the inner product of a function
259 /// \f$f(\boldsymbol{x})\f$ with respect to the derivative (in
260 /// direction \param dir) of all \em local expansion modes
261 /// \f$\phi_n^e(\boldsymbol{x})\f$.
263 const int dir, const Array<OneD, const NekDouble> &inarray,
264 Array<OneD, NekDouble> &outarray);
265
267 const Array<OneD, const NekDouble> &direction,
268 const Array<OneD, const NekDouble> &inarray,
269 Array<OneD, NekDouble> &outarray);
270
271 /// This function calculates the inner product of a
272 /// function \f$f(\boldsymbol{x})\f$ with respect to the
273 /// derivative of all \em local expansion modes
274 /// \f$\phi_n^e(\boldsymbol{x})\f$.
276 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
277 Array<OneD, NekDouble> &outarray);
278 /// This function elementally evaluates the forward transformation
279 /// of a function \f$u(\boldsymbol{x})\f$ onto the global
280 /// spectral/hp expansion.
281 inline void FwdTransLocalElmt(const Array<OneD, const NekDouble> &inarray,
282 Array<OneD, NekDouble> &outarray);
283 ///
284 inline void FwdTrans(const Array<OneD, const NekDouble> &inarray,
285 Array<OneD, NekDouble> &outarray);
287 const NekDouble alpha,
288 const NekDouble exponent,
289 const NekDouble cutoff);
290 /// This function elementally mulplies the coefficient space of
291 /// Sin my the elemental inverse of the mass matrix.
293 const Array<OneD, const NekDouble> &inarray,
294 Array<OneD, NekDouble> &outarray);
295 inline void MultiplyByInvMassMatrix(
296 const Array<OneD, const NekDouble> &inarray,
297 Array<OneD, NekDouble> &outarray);
299 const Array<OneD, const NekDouble> &inarray,
300 Array<OneD, NekDouble> &outarray)
301 {
303 BwdTrans(inarray, tmp);
304 IProductWRTBase(tmp, outarray);
305 }
306 /// Smooth a field across elements
307 inline void SmoothField(Array<OneD, NekDouble> &field);
308
309 /// Solve helmholtz problem
311 const Array<OneD, const NekDouble> &inarray,
312 Array<OneD, NekDouble> &outarray,
313 const StdRegions::ConstFactorMap &factors,
315 const MultiRegions::VarFactorsMap &varfactors =
318 const bool PhysSpaceForcing = true);
319
320 /// Solve Advection Diffusion Reaction
322 const Array<OneD, const NekDouble> &inarray,
323 Array<OneD, NekDouble> &outarray,
324 const StdRegions::ConstFactorMap &factors,
326 const MultiRegions::VarFactorsMap &varfactors =
329 const bool PhysSpaceForcing = true);
330
331 /// Solve Advection Diffusion Reaction
333 const Array<OneD, const NekDouble> &inarray,
334 Array<OneD, NekDouble> &outarray,
335 const StdRegions::ConstFactorMap &factors,
337 const MultiRegions::VarFactorsMap &varfactors =
340 const bool PhysSpaceForcing = true);
341 ///
343 const Array<OneD, const NekDouble> &inarray,
344 Array<OneD, NekDouble> &outarray);
345 /// This function elementally evaluates the backward transformation
346 /// of the global spectral/hp element expansion.
347 inline void BwdTrans(const Array<OneD, const NekDouble> &inarray,
348 Array<OneD, NekDouble> &outarray);
349 /// This function calculates the coordinates of all the elemental
350 /// quadrature points \f$\boldsymbol{x}_i\f$.
351 inline void GetCoords(
352 Array<OneD, NekDouble> &coord_0,
355
356 inline void GetCoords(
357 const int eid, Array<OneD, NekDouble> &coord_0,
360
361 // Homogeneous transforms
362 inline void HomogeneousFwdTrans(const int npts,
363 const Array<OneD, const NekDouble> &inarray,
364 Array<OneD, NekDouble> &outarray,
365 bool Shuff = true, bool UnShuff = true);
366 inline void HomogeneousBwdTrans(const int npts,
367 const Array<OneD, const NekDouble> &inarray,
368 Array<OneD, NekDouble> &outarray,
369 bool Shuff = true, bool UnShuff = true);
370 inline void DealiasedProd(const int num_dofs,
371 const Array<OneD, NekDouble> &inarray1,
372 const Array<OneD, NekDouble> &inarray2,
373 Array<OneD, NekDouble> &outarray);
374 inline void DealiasedDotProd(
375 const int num_dofs, const Array<OneD, Array<OneD, NekDouble>> &inarray1,
376 const Array<OneD, Array<OneD, NekDouble>> &inarray2,
377 Array<OneD, Array<OneD, NekDouble>> &outarray);
378 inline void GetBCValues(Array<OneD, NekDouble> &BndVals,
379 const Array<OneD, NekDouble> &TotField, int BndID);
382 Array<OneD, NekDouble> &outarray,
383 int BndID);
384 inline void NormVectorIProductWRTBase(
386 Array<OneD, NekDouble> &outarray);
387 /// Apply geometry information to each expansion.
389 /// Reset geometry information and reset matrices
391 {
392 v_Reset();
393 }
394
395 /// Reset matrices
397
398 // not sure we really need these in ExpList
399 void WriteTecplotHeader(std::ostream &outfile, std::string var = "")
400 {
401 v_WriteTecplotHeader(outfile, var);
402 }
403 void WriteTecplotZone(std::ostream &outfile, int expansion = -1)
404 {
405 v_WriteTecplotZone(outfile, expansion);
406 }
407 void WriteTecplotField(std::ostream &outfile, int expansion = -1)
408 {
409 v_WriteTecplotField(outfile, expansion);
410 }
411 void WriteTecplotConnectivity(std::ostream &outfile, int expansion = -1)
412 {
413 v_WriteTecplotConnectivity(outfile, expansion);
414 }
415 MULTI_REGIONS_EXPORT void WriteVtkHeader(std::ostream &outfile);
416 MULTI_REGIONS_EXPORT void WriteVtkFooter(std::ostream &outfile);
417 void WriteVtkPieceHeader(std::ostream &outfile, int expansion,
418 int istrip = 0)
419 {
420 v_WriteVtkPieceHeader(outfile, expansion, istrip);
421 }
422 MULTI_REGIONS_EXPORT void WriteVtkPieceFooter(std::ostream &outfile,
423 int expansion);
424 void WriteVtkPieceData(std::ostream &outfile, int expansion,
425 std::string var = "v")
426 {
427 v_WriteVtkPieceData(outfile, expansion, var);
428 }
429
430 /// This function returns the dimension of the coordinates of the
431 /// element \a eid.
432 // inline
433 MULTI_REGIONS_EXPORT int GetCoordim(int eid);
434
435 /// Set the \a i th coefficiient in \a m_coeffs to value \a val
436 inline void SetCoeff(int i, NekDouble val);
437 /// Set the \a i th coefficiient in #m_coeffs to value \a val
438 inline void SetCoeffs(int i, NekDouble val);
439 /// Set the #m_coeffs array to inarray
440 inline void SetCoeffsArray(Array<OneD, NekDouble> &inarray);
441 /// This function returns the dimension of the shape of the
442 /// element \a eid.
443 // inline
445 /// This function returns (a reference to) the array
446 /// \f$\boldsymbol{\hat{u}}_l\f$ (implemented as #m_coeffs)
447 /// containing all local expansion coefficients.
448 inline const Array<OneD, const NekDouble> &GetCoeffs() const;
449 /// Impose Dirichlet Boundary Conditions onto Array
451 /// Add Neumann Boundary Condition forcing to Array
452 inline void ImposeNeumannConditions(Array<OneD, NekDouble> &outarray);
453 /// Add Robin Boundary Condition forcing to Array
454 inline void ImposeRobinConditions(Array<OneD, NekDouble> &outarray);
455 /// Fill Bnd Condition expansion from the values stored in expansion
456 inline void FillBndCondFromField(const Array<OneD, NekDouble> coeffs);
457 /// Fill Bnd Condition expansion in nreg from the values
458 /// stored in expansion
459 inline void FillBndCondFromField(const int nreg,
460 const Array<OneD, NekDouble> coeffs);
461 /// Gathers the global coefficients \f$\boldsymbol{\hat{u}}_g\f$
462 /// from the local coefficients \f$\boldsymbol{\hat{u}}_l\f$.
463 // inline
464 MULTI_REGIONS_EXPORT inline void LocalToGlobal(bool useComm = true);
466 const Array<OneD, const NekDouble> &inarray,
467 Array<OneD, NekDouble> &outarray, bool useComm = true);
468 /// Scatters from the global coefficients
469 /// \f$\boldsymbol{\hat{u}}_g\f$ to the local coefficients
470 /// \f$\boldsymbol{\hat{u}}_l\f$.
471 // inline
472 MULTI_REGIONS_EXPORT inline void GlobalToLocal(void);
473 /**
474 * This operation is evaluated as:
475 * \f{tabbing}
476 * \hspace{1cm} \= Do \= $e=$ $1, N_{\mathrm{el}}$ \\
477 * > > Do \= $i=$ $0,N_m^e-1$ \\
478 * > > > $\boldsymbol{\hat{u}}^{e}[i] = \mbox{sign}[e][i] \cdot
479 * \boldsymbol{\hat{u}}_g[\mbox{map}[e][i]]$ \\
480 * > > continue \\
481 * > continue
482 * \f}
483 * where \a map\f$[e][i]\f$ is the mapping array and \a
484 * sign\f$[e][i]\f$ is an array of similar dimensions ensuring the
485 * correct modal connectivity between the different elements (both
486 * these arrays are contained in the data member #m_locToGloMap). This
487 * operation is equivalent to the scatter operation
488 * \f$\boldsymbol{\hat{u}}_l=\mathcal{A}\boldsymbol{\hat{u}}_g\f$,
489 * where \f$\mathcal{A}\f$ is the
490 * \f$N_{\mathrm{eof}}\times N_{\mathrm{dof}}\f$ permutation matrix.
491 *
492 * @param inarray An array of size \f$N_\mathrm{dof}\f$
493 * containing the global degrees of freedom
494 * \f$\boldsymbol{x}_g\f$.
495 * @param outarray The resulting local degrees of freedom
496 * \f$\boldsymbol{x}_l\f$ will be stored in this
497 * array of size \f$N_\mathrm{eof}\f$.
498 */
500 const Array<OneD, const NekDouble> &inarray,
501 Array<OneD, NekDouble> &outarray);
502 /// Get the \a i th value (coefficient) of #m_coeffs
503 inline NekDouble GetCoeff(int i);
504 /// Get the \a i th value (coefficient) of #m_coeffs
505 inline NekDouble GetCoeffs(int i);
506 /// This function returns (a reference to) the array
507 /// \f$\boldsymbol{u}_l\f$ (implemented as #m_phys) containing the
508 /// function \f$u^{\delta}(\boldsymbol{x})\f$ evaluated at the
509 /// quadrature points.
510 // inline
512 /// This function calculates the \f$L_\infty\f$ error of the global
513 /// spectral/hp element approximation.
517 /// This function calculates the \f$L_\infty\f$ error of the global
518 /// This function calculates the \f$L_2\f$ error with
519 /// respect to soln of the global
520 /// spectral/hp element approximation.
522 const Array<OneD, const NekDouble> &inarray,
524 {
525 return v_L2(inarray, soln);
526 }
527 /// Calculates the \f$H^1\f$ error of the global spectral/hp
528 /// element approximation.
532 /// Calculates the \f$H^1\f$ error of the global spectral/hp
533 /// element approximation.
534 /**
535 * The integration is evaluated locally, that is
536 * \f[\int
537 * f(\boldsymbol{x})d\boldsymbol{x}=\sum_{e=1}^{{N_{\mathrm{el}}}}
538 * \left\{\int_{\Omega_e}f(\boldsymbol{x})d\boldsymbol{x}\right\}, \f]
539 * where the integration over the separate elements is done by the
540 * function StdRegions#StdExpansion#Integral, which discretely
541 * evaluates the integral using Gaussian quadrature.
542 *
543 * Note that the array #m_phys should be filled with the values of the
544 * function \f$f(\boldsymbol{x})\f$ at the quadrature points
545 * \f$\boldsymbol{x}_i\f$.
546 *
547 * @return The value of the discretely evaluated integral
548 * \f$\int f(\boldsymbol{x})d\boldsymbol{x}\f$.
549 */
551 {
552 ASSERTL1(m_physState == true, "local physical space is not true ");
553 return Integral(m_phys);
554 }
555 /**
556 * The integration is evaluated locally, that is
557 * \f[\int
558 * f(\boldsymbol{x})d\boldsymbol{x}=\sum_{e=1}^{{N_{\mathrm{el}}}}
559 * \left\{\int_{\Omega_e}f(\boldsymbol{x})d\boldsymbol{x}\right\}, \f]
560 * where the integration over the separate elements is done by the
561 * function StdRegions#StdExpansion#Integral, which discretely
562 * evaluates the integral using Gaussian quadrature.
563 *
564 * @param inarray An array of size \f$Q_{\mathrm{tot}}\f$
565 * containing the values of the function
566 * \f$f(\boldsymbol{x})\f$ at the quadrature
567 * points \f$\boldsymbol{x}_i\f$.
568 * @return The value of the discretely evaluated integral
569 * \f$\int f(\boldsymbol{x})d\boldsymbol{x}\f$.
570 */
572 {
573 return v_Integral(inarray);
574 }
576 {
577 return v_VectorFlux(inarray);
578 }
579 /// This function calculates the energy associated with
580 /// each one of the modesof a 3D homogeneous nD expansion
585
586 /// This function sets the Spectral Vanishing Viscosity
587 /// in homogeneous1D expansion.
592
593 /// This function returns a vector containing the wave
594 /// numbers in z-direction associated
595 /// with the 3D homogenous expansion. Required if a
596 /// parellelisation is applied in the Fourier direction
598 {
599 return v_GetZIDs();
600 }
601
602 /// This function returns the transposition class
603 /// associated with the homogeneous expansion.
608
609 /// This function returns the Width of homogeneous direction
610 /// associated with the homogeneous expansion.
612 {
613 return v_GetHomoLen();
614 }
615
616 /// This function sets the Width of homogeneous direction
617 /// associated with the homogeneous expansion.
618 void SetHomoLen(const NekDouble lhom)
619 {
620 return v_SetHomoLen(lhom);
621 }
622
623 /// This function returns a vector containing the wave
624 /// numbers in y-direction associated
625 /// with the 3D homogenous expansion. Required if a
626 /// parellelisation is applied in the Fourier direction
628 {
629 return v_GetYIDs();
630 }
631
632 /// This function interpolates the physical space points in
633 /// \a inarray to \a outarray using the same points defined in the
634 /// expansion but where the number of points are rescaled
635 /// by \a 1DScale
637 const Array<OneD, NekDouble> &inarray,
638 Array<OneD, NekDouble> &outarray)
639 {
640 v_PhysInterp1DScaled(scale, inarray, outarray);
641 }
642
643 /// This function Galerkin projects the physical space points in
644 /// \a inarray to \a outarray where inarray is assumed to
645 /// be defined in the expansion but where the number of
646 /// points are rescaled by \a 1DScale
648 const Array<OneD, NekDouble> &inarray,
649 Array<OneD, NekDouble> &outarray)
650 {
651 v_PhysGalerkinProjection1DScaled(scale, inarray, outarray);
652 }
653
654 /// This function returns the number of elements in the expansion.
655 inline int GetExpSize(void);
656
657 /// This function returns the number of elements in the
658 /// expansion which may be different for a homogeoenous extended
659 /// expansionp.
660 inline size_t GetNumElmts(void)
661 {
662 return v_GetNumElmts();
663 }
664
665 /// This function returns the vector of elements in the expansion.
666 inline const std::shared_ptr<LocalRegions::ExpansionVector> GetExp() const;
667
668 /// This function returns (a shared pointer to) the local elemental
669 /// expansion of the \f$n^{\mathrm{th}}\f$ element.
670 inline LocalRegions::ExpansionSharedPtr &GetExp(int n) const;
671
672 /// This function returns (a shared pointer to) the local elemental
673 /// expansion of the \f$n^{\mathrm{th}}\f$ element given a global
674 /// geometry ID.
676
677 /// This function returns (a shared pointer to) the local elemental
678 /// expansion containing the arbitrary point given by \a gloCoord.
680 const Array<OneD, const NekDouble> &gloCoord);
681
682 /**
683 * @brief This function returns the index of the local elemental
684 * expansion containing the arbitrary point given by \a gloCoord,
685 * within a distance tolerance of tol.
686 *
687 * If returnNearestElmt is true and no element contains the point,
688 * this function returns the nearest element whose bounding box contains
689 * the point. The bounding box has a 10% margin in each direction.
690 *
691 * @param gloCoord (input) coordinate in physics space
692 * @param locCoords (output) local coordinate xi in the returned element
693 * @param tol distance tolerance to judge if a point lies in an
694 * element
695 * @param returnNearestElmt if true and no element contains this point, the
696 * nearest element whose bounding box contains this
697 * point is returned
698 * @param cachedId an initial guess of the most possible element index
699 * @param maxDistance if returnNearestElmt is set as true, the nearest
700 * element will be returned. But the distance of the
701 * nearest element and this point should be <=
702 * maxDistance.
703 *
704 * @return element index; if no element is found, -1 is returned.
705 */
707 const Array<OneD, const NekDouble> &gloCoord, NekDouble tol = 0.0,
708 bool returnNearestElmt = false, int cachedId = -1,
709 NekDouble maxDistance = 1e6);
710
711 /** This function returns the index and the Local
712 * Cartesian Coordinates \a locCoords of the local
713 * elemental expansion containing the arbitrary point
714 * given by \a gloCoords.
715 **/
717 const Array<OneD, const NekDouble> &gloCoords,
718 Array<OneD, NekDouble> &locCoords, NekDouble tol = 0.0,
719 bool returnNearestElmt = false, int cachedId = -1,
720 NekDouble maxDistance = 1e6);
721
722 /** This function return the expansion field value
723 * at the coordinates given as input.
724 **/
727 const Array<OneD, const NekDouble> &phys);
728
729 /// Get the start offset position for a local contiguous list of
730 /// coeffs correspoinding to element n.
731 inline int GetCoeff_Offset(int n) const;
732
733 /// Get the start offset position for a local contiguous list of
734 /// quadrature points in a full array correspoinding to element n.
735 inline int GetPhys_Offset(int n) const;
736
737 /// This function returns (a reference to) the array
738 /// \f$\boldsymbol{\hat{u}}_l\f$ (implemented as #m_coeffs)
739 /// containing all local expansion coefficients.
741 /// This function returns (a reference to) the array
742 /// \f$\boldsymbol{u}_l\f$ (implemented as #m_phys) containing the
743 /// function \f$u^{\delta}(\boldsymbol{x})\f$ evaluated at the
744 /// quadrature points.
746 inline void PhysDeriv(Direction edir,
747 const Array<OneD, const NekDouble> &inarray,
749 /// This function discretely evaluates the derivative of a function
750 /// \f$f(\boldsymbol{x})\f$ on the domain consisting of all
751 /// elements of the expansion.
752 inline void PhysDeriv(
753 const Array<OneD, const NekDouble> &inarray,
757 inline void PhysDeriv(const int dir,
758 const Array<OneD, const NekDouble> &inarray,
760 inline void CurlCurl(Array<OneD, Array<OneD, NekDouble>> &Vel,
762 inline void PhysDirectionalDeriv(
763 const Array<OneD, const NekDouble> &direction,
764 const Array<OneD, const NekDouble> &inarray,
765 Array<OneD, NekDouble> &outarray);
766 inline void GetMovingFrames(const SpatialDomains::GeomMMF MMFdir,
767 const Array<OneD, const NekDouble> &CircCentre,
768 Array<OneD, Array<OneD, NekDouble>> &outarray);
769 // functions associated with DisContField
772 /// Get the weight value for boundary conditions
774 /// Set the weight value for boundary conditions
775 inline void SetBndCondBwdWeight(const int index, const NekDouble value);
776 inline std::shared_ptr<ExpList> &UpdateBndCondExpansion(int i);
777 inline void Upwind(const Array<OneD, const NekDouble> &Vn,
781 inline void Upwind(const Array<OneD, const Array<OneD, NekDouble>> &Vec,
785 /**
786 * Return a reference to the trace space associated with this
787 * expansion list.
788 */
789 inline std::shared_ptr<ExpList> &GetTrace();
790 inline std::shared_ptr<AssemblyMapDG> &GetTraceMap(void);
791 inline std::shared_ptr<InterfaceMapDG> &GetInterfaceMap(void);
792 inline const Array<OneD, const int> &GetTraceBndMap(void);
793 inline void GetNormals(Array<OneD, Array<OneD, NekDouble>> &normals);
794 /// Get the length of elements in boundary normal direction
796 Array<OneD, NekDouble> &lengthsFwd, Array<OneD, NekDouble> &lengthsBwd);
797 /// Get the weight value for boundary conditions
798 /// for boundary average and jump calculations
800 Array<OneD, NekDouble> &weightJump);
802 Array<OneD, NekDouble> &outarray);
805 Array<OneD, NekDouble> &outarray);
808 inline void GetFwdBwdTracePhys(const Array<OneD, const NekDouble> &field,
811 bool FillBnd = true,
812 bool PutFwdInBwdOnBCs = false,
813 bool DoExchange = true);
814 inline void FillBwdWithBoundCond(const Array<OneD, NekDouble> &Fwd,
816 bool PutFwdInBwdOnBCs = false);
817 /// Add Fwd and Bwd value to field,
818 /// a reverse procedure of GetFwdBwdTracePhys
825 {
826 v_AddTraceQuadPhysToOffDiag(Fwd, Bwd, field);
827 }
830 Array<OneD, NekDouble> &locTraceFwd,
831 Array<OneD, NekDouble> &locTraceBwd);
832 /// Fill Bwd with boundary conditions
833 inline void FillBwdWithBwdWeight(Array<OneD, NekDouble> &weightave,
834 Array<OneD, NekDouble> &weightjmp);
835 /// Copy and fill the Periodic boundaries
836 inline void PeriodicBwdCopy(const Array<OneD, const NekDouble> &Fwd,
838 /// Rotate Bwd trace for rotational periodicity boundaries
839 /// when the flow is perpendicular to the rotation axis
841 /// Rotate Bwd trace derivative for rotational periodicity boundaries
842 /// when the flow is perpendicular to the rotation axis
844 /// Rotate local Bwd trace across a rotational interface
845 /// when the flow is perpendicular to the rotation axis
847 /// Rotate local Bwd trace derivatives across a rotational interface
848 /// when the flow is perpendicular to the rotation axis
850 inline const std::vector<bool> &GetLeftAdjacentFaces(void) const;
851 inline void ExtractTracePhys(Array<OneD, NekDouble> &outarray);
852 inline void ExtractTracePhys(const Array<OneD, const NekDouble> &inarray,
853 Array<OneD, NekDouble> &outarray,
854 bool gridVelocity = false);
859 inline void EvaluateBoundaryConditions(
860 const NekDouble time = 0.0, const std::string varName = "",
863 // Routines for continous matrix solution
864 /// This function calculates the result of the multiplication of a
865 /// matrix of type specified by \a mkey with a vector given by \a
866 /// inarray.
868 const GlobalMatrixKey &gkey,
869 const Array<OneD, const NekDouble> &inarray,
870 Array<OneD, NekDouble> &outarray);
871 inline void SetUpPhysNormals();
872 inline void GetBoundaryToElmtMap(Array<OneD, int> &ElmtID,
873 Array<OneD, int> &EdgeID);
874 virtual void GetBndElmtExpansion(int i, std::shared_ptr<ExpList> &result,
875 const bool DeclareCoeffPhysArrays = true);
876
877 inline void ExtractElmtToBndPhys(int i, const Array<OneD, NekDouble> &elmt,
878 Array<OneD, NekDouble> &boundary);
879
880 inline void ExtractPhysToBndElmt(int i,
882 Array<OneD, NekDouble> &bndElmt);
883
884 inline void ExtractPhysToBnd(int i,
887
888 inline void GetBoundaryNormals(
889 int i, Array<OneD, Array<OneD, NekDouble>> &normals);
890
892 std::vector<LibUtilities::FieldDefinitionsSharedPtr> &fielddef,
893 int NumHomoDir = 0,
896 std::vector<NekDouble> &HomoLen = LibUtilities::NullNekDoubleVector,
897 bool homoStrips = false,
898 std::vector<unsigned int> &HomoSIDs =
900 std::vector<unsigned int> &HomoZIDs =
902 std::vector<unsigned int> &HomoYIDs =
904
905 std::map<int, RobinBCInfoSharedPtr> GetRobinBCInfo()
906 {
907 return v_GetRobinBCInfo();
908 }
909
910 void GetPeriodicEntities(PeriodicMap &periodicVerts,
911 PeriodicMap &periodicEdges,
912 PeriodicMap &periodicFaces = NullPeriodicMap)
913 {
914 v_GetPeriodicEntities(periodicVerts, periodicEdges, periodicFaces);
915 }
916
917 std::vector<LibUtilities::FieldDefinitionsSharedPtr> GetFieldDefinitions()
918 {
919 return v_GetFieldDefinitions();
920 }
921
923 std::vector<LibUtilities::FieldDefinitionsSharedPtr> &fielddef)
924 {
925 v_GetFieldDefinitions(fielddef);
926 }
927
928 /// Append the element data listed in elements
929 /// fielddef->m_ElementIDs onto fielddata
931 std::vector<NekDouble> &fielddata)
932 {
933 v_AppendFieldData(fielddef, fielddata);
934 }
935 /// Append the data in coeffs listed in elements
936 /// fielddef->m_ElementIDs onto fielddata
938 std::vector<NekDouble> &fielddata,
940 {
941 v_AppendFieldData(fielddef, fielddata, coeffs);
942 }
943 /** \brief Extract the data in fielddata into the coeffs
944 * using the basic ExpList Elemental expansions rather
945 * than planes in homogeneous case
946 */
949 std::vector<NekDouble> &fielddata, std::string &field,
950 Array<OneD, NekDouble> &coeffs);
951 /** \brief Extract the data from fromField using
952 * fromExpList the coeffs using the basic ExpList
953 * Elemental expansions rather than planes in homogeneous
954 * case
955 */
957 const std::shared_ptr<ExpList> &fromExpList,
958 const Array<OneD, const NekDouble> &fromCoeffs,
959 Array<OneD, NekDouble> &toCoeffs);
960 // Extract data in fielddata into the m_coeffs_list for the 3D stability
961 // analysis (base flow is 2D)
964 std::vector<NekDouble> &fielddata, std::string &field,
966 std::unordered_map<int, int> zIdToPlane =
967 std::unordered_map<int, int>());
968 // Extract data from file fileName and put coefficents into array coefffs
970 const std::string &fileName, LibUtilities::CommSharedPtr comm,
971 const std::string &varName, Array<OneD, NekDouble> &coeffs);
973 const int ElementID, const NekDouble scalar1, const NekDouble scalar2,
974 Array<OneD, NekDouble> &outarray);
975 /// Returns a shared pointer to the current object.
976 std::shared_ptr<ExpList> GetSharedThisPtr()
977 {
978 return shared_from_this();
979 }
980 /// Returns the session object
981 std::shared_ptr<LibUtilities::SessionReader> GetSession() const
982 {
983 return m_session;
984 }
985 /// Returns the comm object
986 std::shared_ptr<LibUtilities::Comm> GetComm() const
987 {
988 return m_comm;
989 }
994 // Wrapper functions for Homogeneous Expansions
999 std::shared_ptr<ExpList> &GetPlane(int n)
1000 {
1001 return v_GetPlane(n);
1002 }
1006
1008
1010
1012 GlobalLinSys> &
1014
1015 /// Get m_coeffs to elemental value map
1017 GetCoeffsToElmt() const;
1023 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
1024 const int nDirctn, Array<OneD, DNekMatSharedPtr> &mtxPerVar);
1026 const TensorOfArray3D<NekDouble> &inarray,
1029 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
1032 const Array<OneD, const NekDouble> &FwdFlux,
1033 const Array<OneD, const NekDouble> &BwdFlux,
1034 Array<OneD, NekDouble> &outarray)
1035 {
1036 v_AddTraceIntegralToOffDiag(FwdFlux, BwdFlux, outarray);
1037 }
1039 const int dir, const Array<OneD, const DNekMatSharedPtr> ElmtJacQuad,
1040 Array<OneD, DNekMatSharedPtr> ElmtJacCoef);
1042 const Array<OneD, const DNekMatSharedPtr> ElmtJacQuad,
1043 Array<OneD, DNekMatSharedPtr> ElmtJacCoef);
1045 GetLocTraceToTraceMap() const;
1046 // Return the internal vector which identifieds if trace
1047 // is left adjacent definiing which trace the normal
1048 // points otwards from
1049 MULTI_REGIONS_EXPORT std::vector<bool> &GetLeftAdjacentTraces(void);
1050
1051 /// This function returns the map of index inside m_exp to geom id
1052 MULTI_REGIONS_EXPORT inline const std::unordered_map<int, int> &GetElmtToExpId(
1053 void)
1054 {
1055 return m_elmtToExpId;
1056 }
1057
1058 /// This function returns the index inside m_exp for a given geom id
1060 {
1061 auto it = m_elmtToExpId.find(elmtId);
1062 ASSERTL0(it != m_elmtToExpId.end(), "Global geometry ID " +
1063 std::to_string(elmtId) +
1064 " not found in element ID to "
1065 "expansion ID map.")
1066 return it->second;
1067 }
1068
1069 /// This function returns collections
1072 {
1073 return m_collections;
1074 }
1075
1077 {
1078 return m_coll_coeff_offset[n];
1079 }
1080
1082 {
1083 return m_coll_phys_offset[n];
1084 }
1085
1086 MULTI_REGIONS_EXPORT inline const Array<OneD,
1087 const Array<OneD, NekDouble>> &
1089 {
1090 return m_gridVelocity;
1091 }
1092
1093protected:
1094 /// Pointer holder for PulseWaveSolver
1096 /// Expansion type
1098 std::shared_ptr<DNekMat> GenGlobalMatrixFull(
1099 const GlobalLinSysKey &mkey,
1100 const std::shared_ptr<AssemblyMapCG> &locToGloMap);
1101 /// Communicator
1103 /// Session
1105 /// Mesh associated with this expansion list.
1107 /// The total number of local degrees of freedom. #m_ncoeffs
1108 /// \f$=N_{\mathrm{eof}}=\sum_{e=1}^{{N_{\mathrm{el}}}}N^{e}_l\f$
1110 /** The total number of quadrature points. #m_npoints
1111 *\f$=Q_{\mathrm{tot}}=\sum_{e=1}^{{N_{\mathrm{el}}}}N^{e}_Q\f$
1112 **/
1114 /**
1115 * \brief Concatenation of all local expansion coefficients.
1116 *
1117 * The array of length #m_ncoeffs\f$=N_{\mathrm{eof}}\f$ which is
1118 * the concatenation of the local expansion coefficients
1119 * \f$\hat{u}_n^e\f$ over all \f$N_{\mathrm{el}}\f$ elements
1120 * \f[\mathrm{\texttt{m\_coeffs}}=\boldsymbol{\hat{u}}_{l} =
1121 * \underline{\boldsymbol{\hat{u}}}^e = \left [ \begin{array}{c}
1122 * \boldsymbol{\hat{u}}^{1} \ \
1123 * \boldsymbol{\hat{u}}^{2} \ \
1124 * \vdots \ \
1125 * \boldsymbol{\hat{u}}^{{{N_{\mathrm{el}}}}} \end{array} \right ],
1126 * \quad
1127 * \mathrm{where}\quad \boldsymbol{\hat{u}}^{e}[n]=\hat{u}_n^{e}\f]
1128 */
1130 /**
1131 * \brief The global expansion evaluated at the quadrature points
1132 *
1133 * The array of length #m_npoints\f$=Q_{\mathrm{tot}}\f$ containing
1134 * the evaluation of \f$u^{\delta}(\boldsymbol{x})\f$ at the
1135 * quadrature points \f$\boldsymbol{x}_i\f$.
1136 * \f[\mathrm{\texttt{m\_phys}}=\boldsymbol{u}_{l} =
1137 * \underline{\boldsymbol{u}}^e = \left [ \begin{array}{c}
1138 * \boldsymbol{u}^{1} \ \
1139 * \boldsymbol{u}^{2} \ \
1140 * \vdots \ \
1141 * \boldsymbol{u}^{{{N_{\mathrm{el}}}}} \end{array} \right ],\quad
1142 * \mathrm{where}\quad
1143 * \boldsymbol{u}^{e}[i]=u^{\delta}(\boldsymbol{x}_i)\f]
1144 */
1146 /**
1147 * \brief The state of the array #m_phys.
1148 *
1149 * Indicates whether the array #m_phys, created to contain the
1150 * evaluation of \f$u^{\delta}(\boldsymbol{x})\f$ at the quadrature
1151 * points \f$\boldsymbol{x}_i\f$, is filled with these values.
1152 */
1154 /**
1155 * \brief The list of local expansions.
1156 *
1157 * The (shared pointer to the) vector containing (shared pointers
1158 * to) all local expansions. The fact that the local expansions are
1159 * all stored as a (pointer to a) #StdExpansion, the abstract base
1160 * class for all local expansions, allows a general implementation
1161 * where most of the routines for the derived classes are defined
1162 * in the #ExpList base class.
1163 */
1164 std::shared_ptr<LocalRegions::ExpansionVector> m_exp;
1166 /// Vector of bools to act as an initialise on first call flag
1167 std::vector<bool> m_collectionsDoInit;
1168 /// Offset of elemental data into the array #m_coeffs
1169 std::vector<int> m_coll_coeff_offset;
1170 /// Offset of elemental data into the array #m_phys
1171 std::vector<int> m_coll_phys_offset;
1172 /// Offset of elemental data into the array #m_coeffs
1174 /// Offset of elemental data into the array #m_phys
1176 /// m_coeffs to elemental value map
1179 //@todo should this be in ExpList or
1180 // ExpListHomogeneous1D.cpp it's a bool which determine if
1181 // the expansion is in the wave space (coefficient space)
1182 // or not
1184 /// Mapping from geometry ID of element to index inside #m_exp
1185 std::unordered_map<int, int> m_elmtToExpId;
1186 /// Grid velocity at quadrature points
1188 /// This function assembles the block diagonal matrix of local
1189 /// matrices of the type \a mtype.
1193 const Array<OneD, const NekDouble> &inarray,
1194 Array<OneD, NekDouble> &outarray);
1195 /// Generates a global matrix from the given key and map.
1196 std::shared_ptr<GlobalMatrix> GenGlobalMatrix(
1197 const GlobalMatrixKey &mkey,
1198 const std::shared_ptr<AssemblyMapCG> &locToGloMap);
1200 const std::shared_ptr<DNekMat> &Gmat,
1201 Array<OneD, NekDouble> &EigValsReal,
1202 Array<OneD, NekDouble> &EigValsImag,
1204 /// This operation constructs the global linear system of type \a
1205 /// mkey.
1206 std::shared_ptr<GlobalLinSys> GenGlobalLinSys(
1207 const GlobalLinSysKey &mkey,
1208 const std::shared_ptr<AssemblyMapCG> &locToGloMap);
1209 /// Generate a GlobalLinSys from information provided by the key
1210 /// "mkey" and the mapping provided in LocToGloBaseMap.
1211 std::shared_ptr<GlobalLinSys> GenGlobalBndLinSys(
1212 const GlobalLinSysKey &mkey, const AssemblyMapSharedPtr &locToGloMap);
1213 // Virtual prototypes
1214 virtual size_t v_GetNumElmts(void)
1215 {
1216 return (*m_exp).size();
1217 }
1221 virtual void v_SetBndCondBwdWeight(const int index, const NekDouble value);
1222 virtual std::shared_ptr<ExpList> &v_UpdateBndCondExpansion(int i);
1223 virtual void v_Upwind(const Array<OneD, const Array<OneD, NekDouble>> &Vec,
1231 virtual std::shared_ptr<ExpList> &v_GetTrace();
1232 virtual std::shared_ptr<AssemblyMapDG> &v_GetTraceMap();
1233 virtual std::shared_ptr<InterfaceMapDG> &v_GetInterfaceMap();
1235 virtual const std::shared_ptr<LocTraceToTraceMap> &v_GetLocTraceToTraceMap(
1236 void) const;
1237 virtual std::vector<bool> &v_GetLeftAdjacentTraces(void);
1238 /// Populate \a normals with the normals of all expansions.
1241 Array<OneD, NekDouble> &outarray);
1245 Array<OneD, NekDouble> &outarray);
1251 bool FillBnd = true,
1252 bool PutFwdInBwdOnBCs = false,
1253 bool DoExchange = true);
1256 bool PutFwdInBwdOnBCs);
1266 Array<OneD, NekDouble> &locTraceFwd,
1267 Array<OneD, NekDouble> &locTraceBwd);
1269 Array<OneD, NekDouble> &weightjmp);
1273
1276
1278
1279 virtual const std::vector<bool> &v_GetLeftAdjacentFaces(void) const;
1282 Array<OneD, NekDouble> &outarray,
1283 bool gridVelocity = false);
1285 const Array<OneD, const NekDouble> &inarray,
1286 Array<OneD, NekDouble> &outarray);
1287
1289 const Array<OneD, const NekDouble> &inarray,
1290 Array<OneD, NekDouble> &outarray,
1291 const StdRegions::ConstFactorMap &factors,
1292 const StdRegions::VarCoeffMap &varcoeff,
1293 const MultiRegions::VarFactorsMap &varfactors,
1294 const Array<OneD, const NekDouble> &dirForcing,
1295 const bool PhysSpaceForcing);
1296
1298 const Array<OneD, const NekDouble> &inarray,
1299 Array<OneD, NekDouble> &outarray,
1300 const StdRegions::ConstFactorMap &factors,
1301 const StdRegions::VarCoeffMap &varcoeff,
1302 const MultiRegions::VarFactorsMap &varfactors,
1303 const Array<OneD, const NekDouble> &dirForcing,
1304 const bool PhysSpaceForcing);
1305
1307 const Array<OneD, const NekDouble> &inarray,
1308 Array<OneD, NekDouble> &outarray,
1309 const StdRegions::ConstFactorMap &factors,
1310 const StdRegions::VarCoeffMap &varcoeff,
1311 const MultiRegions::VarFactorsMap &varfactors,
1312 const Array<OneD, const NekDouble> &dirForcing,
1313 const bool PhysSpaceForcing);
1314
1315 // wrapper functions about virtual functions
1320 virtual void v_FillBndCondFromField(const int nreg,
1321 const Array<OneD, NekDouble> coeffs);
1322 virtual void v_Reset();
1323 virtual void v_LocalToGlobal(bool UseComm);
1325 Array<OneD, NekDouble> &outarray,
1326 bool UseComm);
1327 virtual void v_GlobalToLocal(void);
1329 Array<OneD, NekDouble> &outarray);
1330 virtual void v_BwdTrans(const Array<OneD, const NekDouble> &inarray,
1331 Array<OneD, NekDouble> &outarray);
1332 virtual void v_FwdTrans(const Array<OneD, const NekDouble> &inarray,
1333 Array<OneD, NekDouble> &outarray);
1335 const Array<OneD, const NekDouble> &inarray,
1336 Array<OneD, NekDouble> &outarray);
1338 const Array<OneD, const NekDouble> &inarray,
1339 Array<OneD, NekDouble> &outarray);
1342 Array<OneD, NekDouble> &outarray);
1343
1344 // Define ExpList::IProductWRTDerivBase as virtual function
1346 const int dir, const Array<OneD, const NekDouble> &inarray,
1347 Array<OneD, NekDouble> &outarray);
1348
1350 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
1351 Array<OneD, NekDouble> &outarray);
1352
1353 virtual void v_GetCoords(
1354 Array<OneD, NekDouble> &coord_0,
1357
1358 virtual void v_GetCoords(const int eid, Array<OneD, NekDouble> &xc0,
1361
1362 virtual void v_PhysDeriv(const Array<OneD, const NekDouble> &inarray,
1363 Array<OneD, NekDouble> &out_d0,
1364 Array<OneD, NekDouble> &out_d1,
1365 Array<OneD, NekDouble> &out_d2);
1366 virtual void v_PhysDeriv(const int dir,
1367 const Array<OneD, const NekDouble> &inarray,
1368 Array<OneD, NekDouble> &out_d);
1369 virtual void v_PhysDeriv(Direction edir,
1370 const Array<OneD, const NekDouble> &inarray,
1371 Array<OneD, NekDouble> &out_d);
1372
1375
1378
1380 const Array<OneD, const NekDouble> &direction,
1381 const Array<OneD, const NekDouble> &inarray,
1382 Array<OneD, NekDouble> &outarray);
1383 virtual void v_GetMovingFrames(
1384 const SpatialDomains::GeomMMF MMFdir,
1385 const Array<OneD, const NekDouble> &CircCentre,
1386 Array<OneD, Array<OneD, NekDouble>> &outarray);
1388 const int npts, const Array<OneD, const NekDouble> &inarray,
1389 Array<OneD, NekDouble> &outarray, bool Shuff = true,
1390 bool UnShuff = true);
1392 const int npts, const Array<OneD, const NekDouble> &inarray,
1393 Array<OneD, NekDouble> &outarray, bool Shuff = true,
1394 bool UnShuff = true);
1395 virtual void v_DealiasedProd(const int num_dofs,
1396 const Array<OneD, NekDouble> &inarray1,
1397 const Array<OneD, NekDouble> &inarray2,
1398 Array<OneD, NekDouble> &outarray);
1400 const int num_dofs, const Array<OneD, Array<OneD, NekDouble>> &inarray1,
1401 const Array<OneD, Array<OneD, NekDouble>> &inarray2,
1402 Array<OneD, Array<OneD, NekDouble>> &outarray);
1404 const Array<OneD, NekDouble> &TotField,
1405 int BndID);
1408 Array<OneD, NekDouble> &outarray,
1409 int BndID);
1412 Array<OneD, NekDouble> &outarray);
1413 virtual void v_SetUpPhysNormals();
1415 Array<OneD, int> &EdgeID);
1416 virtual void v_GetBndElmtExpansion(int i, std::shared_ptr<ExpList> &result,
1417 const bool DeclareCoeffPhysArrays);
1418
1419 virtual void v_ExtractElmtToBndPhys(const int i,
1420 const Array<OneD, NekDouble> &elmt,
1421 Array<OneD, NekDouble> &boundary);
1422
1424 const int i, const Array<OneD, const NekDouble> &phys,
1425 Array<OneD, NekDouble> &bndElmt);
1426
1427 virtual void v_ExtractPhysToBnd(const int i,
1428 const Array<OneD, const NekDouble> &phys,
1430
1432 int i, Array<OneD, Array<OneD, NekDouble>> &normals);
1433
1434 virtual std::vector<LibUtilities::FieldDefinitionsSharedPtr>
1436
1438 std::vector<LibUtilities::FieldDefinitionsSharedPtr> &fielddef);
1439
1440 virtual void v_AppendFieldData(
1442 std::vector<NekDouble> &fielddata);
1443 virtual void v_AppendFieldData(
1445 std::vector<NekDouble> &fielddata, Array<OneD, NekDouble> &coeffs);
1448 std::vector<NekDouble> &fielddata, std::string &field,
1449 Array<OneD, NekDouble> &coeffs,
1450 std::unordered_map<int, int> zIdToPlane);
1452 const std::shared_ptr<ExpList> &fromExpList,
1453 const Array<OneD, const NekDouble> &fromCoeffs,
1454 Array<OneD, NekDouble> &toCoeffs);
1455 virtual void v_WriteTecplotHeader(std::ostream &outfile,
1456 std::string var = "");
1457 virtual void v_WriteTecplotZone(std::ostream &outfile, int expansion);
1458 virtual void v_WriteTecplotField(std::ostream &outfile, int expansion);
1459 virtual void v_WriteTecplotConnectivity(std::ostream &outfile,
1460 int expansion);
1461 virtual void v_WriteVtkPieceData(std::ostream &outfile, int expansion,
1462 std::string var);
1463 virtual void v_WriteVtkPieceHeader(std::ostream &outfile, int expansion,
1464 int istrip);
1465
1467 const Array<OneD, const NekDouble> &phys,
1469
1472 const Array<OneD, Array<OneD, NekDouble>> &inarray);
1473
1475
1478 virtual void v_SetHomoLen(const NekDouble lhom);
1481
1482 // 1D Scaling and projection
1483 virtual void v_PhysInterp1DScaled(const NekDouble scale,
1484 const Array<OneD, NekDouble> &inarray,
1485 Array<OneD, NekDouble> &outarray);
1486
1488 const NekDouble scale, const Array<OneD, NekDouble> &inarray,
1489 Array<OneD, NekDouble> &outarray);
1490
1491 virtual void v_ClearGlobalLinSysManager(void);
1492
1493 virtual int v_GetPoolCount(std::string);
1494
1496
1499
1500 void ExtractFileBCs(const std::string &fileName,
1502 const std::string &varName,
1503 const std::shared_ptr<ExpList> locExpList);
1504
1505 // Utility function for a common case of retrieving a
1506 // BoundaryCondition from a boundary condition collection.
1510 unsigned int index, const std::string &variable);
1511
1514
1517
1519 const NekDouble time = 0.0, const std::string varName = "",
1522
1523 virtual std::map<int, RobinBCInfoSharedPtr> v_GetRobinBCInfo(void);
1524
1525 virtual void v_GetPeriodicEntities(PeriodicMap &periodicVerts,
1526 PeriodicMap &periodicEdges,
1527 PeriodicMap &periodicFaces);
1528
1529 // Homogeneous direction wrapper functions.
1531 {
1533 "This method is not defined or valid for this class type");
1535 }
1536
1537 // wrapper function to set viscosity for Homo1D expansion
1539 [[maybe_unused]] Array<OneD, NekDouble> visc)
1540 {
1542 "This method is not defined or valid for this class type");
1543 }
1544
1545 virtual std::shared_ptr<ExpList> &v_GetPlane(int n);
1546
1548 const Array<OneD, const NekDouble> &FwdFlux,
1549 const Array<OneD, const NekDouble> &BwdFlux,
1550 Array<OneD, NekDouble> &outarray);
1551
1552private:
1553 /// Definition of the total number of degrees of freedom and
1554 /// quadrature points and offsets to access data
1555 void SetupCoeffPhys(bool DeclareCoeffPhysArrays = true,
1556 bool SetupOffsets = true);
1557
1558 /// Define a list of elements using the geometry and basis
1559 /// key information in expmap;
1561};
1562
1563/// An empty ExpList object.
1566
1567// An empty GlobaLinSysManager object
1571
1572// Inline routines follow.
1573
1574/**
1575 * Returns the total number of local degrees of freedom
1576 * \f$N_{\mathrm{eof}}=\sum_{e=1}^{{N_{\mathrm{el}}}}N^{e}_m\f$.
1577 */
1578inline int ExpList::GetNcoeffs() const
1579{
1580 return m_ncoeffs;
1581}
1582
1583inline int ExpList::GetNcoeffs(const int eid) const
1584{
1585 return (*m_exp)[eid]->GetNcoeffs();
1586}
1587
1588/**
1589 * Evaulates the maximum number of modes in the elemental basis
1590 * order over all elements
1591 */
1593{
1594 int returnval = 0;
1595
1596 for (size_t i = 0; i < (*m_exp).size(); ++i)
1597 {
1598 returnval = (std::max)(returnval, (*m_exp)[i]->EvalBasisNumModesMax());
1599 }
1600
1601 return returnval;
1602}
1603
1604/**
1605 *
1606 */
1608{
1609 Array<OneD, int> returnval((*m_exp).size(), 0);
1610
1611 for (size_t i = 0; i < (*m_exp).size(); ++i)
1612 {
1613 returnval[i] =
1614 (std::max)(returnval[i], (*m_exp)[i]->EvalBasisNumModesMax());
1615 }
1616
1617 return returnval;
1618}
1619
1620/**
1621 *
1622 */
1623inline int ExpList::GetTotPoints() const
1624{
1625 return m_npoints;
1626}
1627
1628inline int ExpList::GetTotPoints(const int eid) const
1629{
1630 return (*m_exp)[eid]->GetTotPoints();
1631}
1632
1633inline int ExpList::Get1DScaledTotPoints(const NekDouble scale) const
1634{
1635 int returnval = 0;
1636 size_t cnt;
1637 size_t nbase = (*m_exp)[0]->GetNumBases();
1638
1639 for (size_t i = 0; i < (*m_exp).size(); ++i)
1640 {
1641 int npt0 = (*m_exp)[i]->GetNumPoints(0);
1642 cnt = 1;
1643
1644 for (size_t j = 0; j < nbase; ++j)
1645 {
1646 int npt = (*m_exp)[i]->GetNumPoints(j);
1647 cnt *= (npt0 - npt == 1) ? (size_t)(scale * npt0 - 1)
1648 : (size_t)(scale * npt);
1649 }
1650 returnval += cnt;
1651 }
1652 return returnval;
1653}
1654
1655/**
1656 *
1657 */
1658inline int ExpList::GetNpoints() const
1659{
1660 return m_npoints;
1661}
1662
1663/**
1664 *
1665 */
1666inline void ExpList::SetWaveSpace(const bool wavespace)
1667{
1668 m_WaveSpace = wavespace;
1669}
1670
1671/**
1672 *
1673 */
1674inline bool ExpList::GetWaveSpace() const
1675{
1676 return m_WaveSpace;
1677}
1678
1679/// Set the \a i th value of\a m_phys to value \a val
1680inline void ExpList::SetPhys(int i, NekDouble val)
1681{
1682 m_phys[i] = val;
1683}
1684/**
1685 * This function fills the array \f$\boldsymbol{u}_l\f$, the evaluation
1686 * of the expansion at the quadrature points (implemented as #m_phys),
1687 * with the values of the array \a inarray.
1688 *
1689 * @param inarray The array containing the values where
1690 * #m_phys should be filled with.
1691 */
1693{
1694 ASSERTL0((int)inarray.size() == m_npoints,
1695 "Input array does not have correct number of elements.");
1696 Vmath::Vcopy(m_npoints, &inarray[0], 1, &m_phys[0], 1);
1697 m_physState = true;
1698}
1700{
1701 m_phys = inarray;
1702}
1703/**
1704 * @param physState \a true (=filled) or \a false (=not filled).
1705 */
1706inline void ExpList::SetPhysState(const bool physState)
1707{
1708 m_physState = physState;
1709}
1710/**
1711 * @return physState \a true (=filled) or \a false (=not filled).
1712 */
1713inline bool ExpList::GetPhysState() const
1714{
1715 return m_physState;
1716}
1717/**
1718 *
1719 */
1721 const Array<OneD, const NekDouble> &inarray,
1722 Array<OneD, NekDouble> &outarray)
1723{
1724 v_IProductWRTBase(inarray, outarray);
1725}
1726/**
1727 *
1728 */
1730 const int dir, const Array<OneD, const NekDouble> &inarray,
1731 Array<OneD, NekDouble> &outarray)
1732{
1733 v_IProductWRTDerivBase(dir, inarray, outarray);
1734}
1735
1736/**
1737 *
1738 */
1740 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
1741 Array<OneD, NekDouble> &outarray)
1742{
1743 v_IProductWRTDerivBase(inarray, outarray);
1744}
1745
1746/**
1747 *
1748 */
1750 Array<OneD, NekDouble> &outarray)
1751{
1752 v_FwdTrans(inarray, outarray);
1753}
1754/**
1755 *
1756 */
1758 const Array<OneD, const NekDouble> &inarray,
1759 Array<OneD, NekDouble> &outarray)
1760{
1761 v_FwdTransLocalElmt(inarray, outarray);
1762}
1763/**
1764 *
1765 */
1767 const Array<OneD, const NekDouble> &inarray,
1768 Array<OneD, NekDouble> &outarray)
1769{
1770 v_FwdTransBndConstrained(inarray, outarray);
1771}
1772/**
1773 *
1774 */
1776{
1777 v_SmoothField(field);
1778}
1779/**
1780 *
1781 */
1782/**
1783 * Given the coefficients of an expansion, this function evaluates the
1784 * spectral/hp expansion \f$u^{\delta}(\boldsymbol{x})\f$ at the
1785 * quadrature points \f$\boldsymbol{x}_i\f$.
1786 */
1788 Array<OneD, NekDouble> &outarray)
1789{
1790 v_BwdTrans(inarray, outarray);
1791}
1792/**
1793 *
1794 */
1796 const Array<OneD, const NekDouble> &inarray,
1797 Array<OneD, NekDouble> &outarray)
1798{
1799 v_MultiplyByInvMassMatrix(inarray, outarray);
1800}
1801/**
1802 *
1803 */
1805 const Array<OneD, const NekDouble> &inarray,
1806 Array<OneD, NekDouble> &outarray, const StdRegions::ConstFactorMap &factors,
1807 const StdRegions::VarCoeffMap &varcoeff,
1808 const MultiRegions::VarFactorsMap &varfactors,
1809 const Array<OneD, const NekDouble> &dirForcing, const bool PhysSpaceForcing)
1810{
1811 return v_HelmSolve(inarray, outarray, factors, varcoeff, varfactors,
1812 dirForcing, PhysSpaceForcing);
1813}
1814/**
1815 *
1816 */
1818 const Array<OneD, const NekDouble> &inarray,
1819 Array<OneD, NekDouble> &outarray, const StdRegions::ConstFactorMap &factors,
1820 const StdRegions::VarCoeffMap &varcoeff,
1821 const MultiRegions::VarFactorsMap &varfactors,
1822 const Array<OneD, const NekDouble> &dirForcing, const bool PhysSpaceForcing)
1823{
1825 inarray, outarray, factors, varcoeff, varfactors, dirForcing,
1826 PhysSpaceForcing);
1827}
1828
1830 const Array<OneD, const NekDouble> &inarray,
1831 Array<OneD, NekDouble> &outarray, const StdRegions::ConstFactorMap &factors,
1832 const StdRegions::VarCoeffMap &varcoeff,
1833 const MultiRegions::VarFactorsMap &varfactors,
1834 const Array<OneD, const NekDouble> &dirForcing, const bool PhysSpaceForcing)
1835{
1836 return v_LinearAdvectionReactionSolve(inarray, outarray, factors, varcoeff,
1837 varfactors, dirForcing,
1838 PhysSpaceForcing);
1839}
1840/**
1841 *
1842 */
1844 Array<OneD, NekDouble> &coord_1,
1845 Array<OneD, NekDouble> &coord_2)
1846{
1847 v_GetCoords(coord_0, coord_1, coord_2);
1848}
1849
1850inline void ExpList::GetCoords(const int eid, Array<OneD, NekDouble> &xc0,
1853{
1854 v_GetCoords(eid, xc0, xc1, xc2);
1855}
1856
1857/**
1858 *
1859 */
1861 const SpatialDomains::GeomMMF MMFdir,
1862 const Array<OneD, const NekDouble> &CircCentre,
1863 Array<OneD, Array<OneD, NekDouble>> &outarray)
1864{
1865 v_GetMovingFrames(MMFdir, CircCentre, outarray);
1866}
1867/**
1868 *
1869 */
1871 Array<OneD, NekDouble> &out_d0,
1872 Array<OneD, NekDouble> &out_d1,
1873 Array<OneD, NekDouble> &out_d2)
1874{
1875 v_PhysDeriv(inarray, out_d0, out_d1, out_d2);
1876}
1877/**
1878 *
1879 */
1880inline void ExpList::PhysDeriv(const int dir,
1881 const Array<OneD, const NekDouble> &inarray,
1883{
1884 v_PhysDeriv(dir, inarray, out_d);
1885}
1887 const Array<OneD, const NekDouble> &inarray,
1889{
1890 v_PhysDeriv(edir, inarray, out_d);
1891}
1892/**
1893 *
1894 */
1896 const Array<OneD, const NekDouble> &direction,
1897 const Array<OneD, const NekDouble> &inarray,
1898 Array<OneD, NekDouble> &outarray)
1899{
1900 v_PhysDirectionalDeriv(direction, inarray, outarray);
1901}
1902/**
1903 *
1904 */
1907{
1908 v_CurlCurl(Vel, Q);
1909}
1910/**
1911 *
1912 */
1914 const int npts, const Array<OneD, const NekDouble> &inarray,
1915 Array<OneD, NekDouble> &outarray, bool Shuff, bool UnShuff)
1916{
1917 v_HomogeneousFwdTrans(npts, inarray, outarray, Shuff, UnShuff);
1918}
1919/**
1920 *
1921 */
1923 const int npts, const Array<OneD, const NekDouble> &inarray,
1924 Array<OneD, NekDouble> &outarray, bool Shuff, bool UnShuff)
1925{
1926 v_HomogeneousBwdTrans(npts, inarray, outarray, Shuff, UnShuff);
1927}
1928/**
1929 *
1930 */
1931inline void ExpList::DealiasedProd(const int num_dofs,
1932 const Array<OneD, NekDouble> &inarray1,
1933 const Array<OneD, NekDouble> &inarray2,
1934 Array<OneD, NekDouble> &outarray)
1935{
1936 v_DealiasedProd(num_dofs, inarray1, inarray2, outarray);
1937}
1938/**
1939 *
1940 */
1942 const int num_dofs, const Array<OneD, Array<OneD, NekDouble>> &inarray1,
1943 const Array<OneD, Array<OneD, NekDouble>> &inarray2,
1944 Array<OneD, Array<OneD, NekDouble>> &outarray)
1945{
1946 v_DealiasedDotProd(num_dofs, inarray1, inarray2, outarray);
1947}
1948/**
1949 *
1950 */
1952 const Array<OneD, NekDouble> &TotField,
1953 int BndID)
1954{
1955 v_GetBCValues(BndVals, TotField, BndID);
1956}
1957/**
1958 *
1959 */
1962 Array<OneD, NekDouble> &outarray,
1963 int BndID)
1964{
1965 v_NormVectorIProductWRTBase(V1, V2, outarray, BndID);
1966}
1972/**
1973 * @param eid The index of the element to be checked.
1974 * @return The dimension of the coordinates of the specific element.
1975 */
1976inline int ExpList::GetCoordim(int eid)
1977{
1978 ASSERTL2(eid <= (*m_exp).size(), "eid is larger than number of elements");
1979 return (*m_exp)[eid]->GetCoordim();
1980}
1981/**
1982 * @param eid The index of the element to be checked.
1983 * @return The dimension of the shape of the specific element.
1984 */
1986{
1987 return (*m_exp)[0]->GetShapeDimension();
1988}
1989/**
1990 * @param i The index of m_coeffs to be set
1991 * @param val The value which m_coeffs[i] is to be set to.
1992 */
1993inline void ExpList::SetCoeff(int i, NekDouble val)
1994{
1995 m_coeffs[i] = val;
1996}
1997/**
1998 * @param i The index of #m_coeffs to be set.
1999 * @param val The value which #m_coeffs[i] is to be set to.
2000 */
2001inline void ExpList::SetCoeffs(int i, NekDouble val)
2002{
2003 m_coeffs[i] = val;
2004}
2006{
2007 m_coeffs = inarray;
2008}
2009/**
2010 * As the function returns a constant reference to a
2011 * <em>const Array</em>, it is not possible to modify the
2012 * underlying data of the array #m_coeffs. In order to do
2013 * so, use the function #UpdateCoeffs instead.
2014 *
2015 * @return (A constant reference to) the array #m_coeffs.
2016 */
2018{
2019 return m_coeffs;
2020}
2030{
2031 v_ImposeRobinConditions(outarray);
2032}
2034{
2035 v_FillBndCondFromField(coeffs);
2036}
2037inline void ExpList::FillBndCondFromField(const int nreg,
2038 const Array<OneD, NekDouble> coeffs)
2039{
2040 v_FillBndCondFromField(nreg, coeffs);
2041}
2042inline void ExpList::LocalToGlobal(bool useComm)
2043{
2044 v_LocalToGlobal(useComm);
2045}
2047 Array<OneD, NekDouble> &outarray,
2048 bool useComm)
2049{
2050 v_LocalToGlobal(inarray, outarray, useComm);
2051}
2052inline void ExpList::GlobalToLocal(void)
2053{
2055}
2056/**
2057 * This operation is evaluated as:
2058 * \f{tabbing}
2059 * \hspace{1cm} \= Do \= $e=$ $1, N_{\mathrm{el}}$ \\
2060 * > > Do \= $i=$ $0,N_m^e-1$ \\
2061 * > > > $\boldsymbol{\hat{u}}^{e}[i] = \mbox{sign}[e][i] \cdot
2062 * \boldsymbol{\hat{u}}_g[\mbox{map}[e][i]]$ \\
2063 * > > continue \\
2064 * > continue
2065 * \f}
2066 * where \a map\f$[e][i]\f$ is the mapping array and \a
2067 * sign\f$[e][i]\f$ is an array of similar dimensions ensuring the
2068 * correct modal connectivity between the different elements (both
2069 * these arrays are contained in the data member #m_locToGloMap). This
2070 * operation is equivalent to the scatter operation
2071 * \f$\boldsymbol{\hat{u}}_l=\mathcal{A}\boldsymbol{\hat{u}}_g\f$, where
2072 * \f$\mathcal{A}\f$ is the
2073 * \f$N_{\mathrm{eof}}\times N_{\mathrm{dof}}\f$ permutation matrix.
2074 *
2075 * @param inarray An array of size \f$N_\mathrm{dof}\f$
2076 * containing the global degrees of freedom
2077 * \f$\boldsymbol{x}_g\f$.
2078 * @param outarray The resulting local degrees of freedom
2079 * \f$\boldsymbol{x}_l\f$ will be stored in this
2080 * array of size \f$N_\mathrm{eof}\f$.
2081 */
2083 Array<OneD, NekDouble> &outarray)
2084{
2085 v_GlobalToLocal(inarray, outarray);
2086}
2087/**
2088 * @param i The index of #m_coeffs to be returned
2089 * @return The NekDouble held in #m_coeffs[i].
2090 */
2092{
2093 return m_coeffs[i];
2094}
2095/**
2096 * @param i The index of #m_coeffs to be returned
2097 * @return The NekDouble held in #m_coeffs[i].
2098 */
2100{
2101 return m_coeffs[i];
2102}
2103/**
2104 * As the function returns a constant reference to a
2105 * <em>const Array</em> it is not possible to modify the
2106 * underlying data of the array #m_phys. In order to do
2107 * so, use the function #UpdatePhys instead.
2108 *
2109 * @return (A constant reference to) the array #m_phys.
2110 */
2112{
2113 return m_phys;
2114}
2115/**
2116 * @return \f$N_{\mathrm{el}}\f$, the number of elements in the
2117 * expansion.
2118 */
2119inline int ExpList::GetExpSize(void)
2120{
2121 return (*m_exp).size();
2122}
2123/**
2124 * @param n The index of the element concerned.
2125 *
2126 * @return (A shared pointer to) the local expansion of the
2127 * \f$n^{\mathrm{th}}\f$ element.
2128 */
2130{
2131 return (*m_exp)[n];
2132}
2133/**
2134 * @param n The global id of the element concerned.
2135 *
2136 * @return (A shared pointer to) the local expansion of the
2137 * \f$n^{\mathrm{th}}\f$ element.
2138 */
2140{
2141 auto it = m_elmtToExpId.find(n);
2142 ASSERTL0(it != m_elmtToExpId.end(), "Global geometry ID " +
2143 std::to_string(n) +
2144 " not found in element ID to "
2145 "expansion ID map.")
2146 return (*m_exp)[it->second];
2147}
2148/**
2149 * @return (A const shared pointer to) the local expansion vector #m_exp
2150 */
2151inline const std::shared_ptr<LocalRegions::ExpansionVector> ExpList::GetExp(
2152 void) const
2153{
2154 return m_exp;
2155}
2156/**
2157 *
2158 */
2159inline int ExpList::GetCoeff_Offset(int n) const
2160{
2161 return m_coeff_offset[n];
2162}
2163/**
2164 *
2165 */
2166inline int ExpList::GetPhys_Offset(int n) const
2167{
2168 return m_phys_offset[n];
2169}
2170/**
2171 * If one wants to get hold of the underlying data without modifying
2172 * them, rather use the function #GetCoeffs instead.
2173 *
2174 * @return (A reference to) the array #m_coeffs.
2175 */
2177{
2178 return m_coeffs;
2179}
2180/**
2181 * If one wants to get hold of the underlying data without modifying
2182 * them, rather use the function #GetPhys instead.
2183 *
2184 * @return (A reference to) the array #m_phys.
2185 */
2187{
2188 m_physState = true;
2189 return m_phys;
2190}
2191// functions associated with DisContField
2197/// Get m_coeffs to elemental value map
2212inline void ExpList::SetBndCondBwdWeight(const int index, const NekDouble value)
2213{
2214 v_SetBndCondBwdWeight(index, value);
2215}
2216inline std::shared_ptr<ExpList> &ExpList::UpdateBndCondExpansion(int i)
2217{
2218 return v_UpdateBndCondExpansion(i);
2219}
2221 const Array<OneD, const Array<OneD, NekDouble>> &Vec,
2224{
2225 v_Upwind(Vec, Fwd, Bwd, Upwind);
2226}
2230 Array<OneD, NekDouble> &Upwind)
2231{
2232 v_Upwind(Vn, Fwd, Bwd, Upwind);
2233}
2234inline std::shared_ptr<ExpList> &ExpList::GetTrace()
2235{
2236 return v_GetTrace();
2237}
2238inline std::shared_ptr<AssemblyMapDG> &ExpList::GetTraceMap()
2239{
2240 return v_GetTraceMap();
2241}
2242inline std::shared_ptr<InterfaceMapDG> &ExpList::GetInterfaceMap()
2243{
2244 return v_GetInterfaceMap();
2245}
2247{
2248 return v_GetTraceBndMap();
2249}
2251{
2252 v_GetNormals(normals);
2253}
2255 Array<OneD, NekDouble> &outarray)
2256{
2257 v_AddTraceIntegral(Fn, outarray);
2258}
2262{
2263 v_AddFwdBwdTraceIntegral(Fwd, Bwd, outarray);
2264}
2272 Array<OneD, NekDouble> &Bwd, bool FillBnd, bool PutFwdInBwdOnBCs,
2273 bool DoExchange)
2274{
2275 v_GetFwdBwdTracePhys(field, Fwd, Bwd, FillBnd, PutFwdInBwdOnBCs,
2276 DoExchange);
2277}
2280 bool PutFwdInBwdOnBCs)
2281{
2282 v_FillBwdWithBoundCond(Fwd, Bwd, PutFwdInBwdOnBCs);
2283}
2287{
2288 v_AddTraceQuadPhysToField(Fwd, Bwd, field);
2289}
2293 Array<OneD, NekDouble> &locTraceFwd, Array<OneD, NekDouble> &locTraceBwd)
2294{
2295 v_GetLocTraceFromTracePts(Fwd, Bwd, locTraceFwd, locTraceBwd);
2296}
2298 Array<OneD, NekDouble> &weightjmp)
2299{
2300 v_FillBwdWithBwdWeight(weightave, weightjmp);
2301}
2304{
2305 v_PeriodicBwdCopy(Fwd, Bwd);
2306}
2311
2316
2321
2326
2327inline const std::vector<bool> &ExpList::GetLeftAdjacentFaces(void) const
2328{
2329 return v_GetLeftAdjacentFaces();
2330}
2332{
2333 v_ExtractTracePhys(outarray);
2334}
2336 const Array<OneD, const NekDouble> &inarray,
2337 Array<OneD, NekDouble> &outarray, bool gridVelocity)
2338{
2339 v_ExtractTracePhys(inarray, outarray, gridVelocity);
2340}
2352 const std::string varName,
2353 const NekDouble x2_in,
2354 const NekDouble x3_in)
2355{
2356 v_EvaluateBoundaryConditions(time, varName, x2_in, x3_in);
2357}
2359{
2361}
2363 Array<OneD, int> &EdgeID)
2364{
2365 v_GetBoundaryToElmtMap(ElmtID, EdgeID);
2366}
2368 std::shared_ptr<ExpList> &result,
2369 const bool DeclareCoeffPhysArrays)
2370{
2371 v_GetBndElmtExpansion(i, result, DeclareCoeffPhysArrays);
2372}
2373
2375 const Array<OneD, NekDouble> &elmt,
2376 Array<OneD, NekDouble> &boundary)
2377{
2378 v_ExtractElmtToBndPhys(i, elmt, boundary);
2379}
2380
2382 int i, const Array<OneD, const NekDouble> &phys,
2383 Array<OneD, NekDouble> &bndElmt)
2384{
2385 v_ExtractPhysToBndElmt(i, phys, bndElmt);
2386}
2387
2389 const Array<OneD, const NekDouble> &phys,
2391{
2392 v_ExtractPhysToBnd(i, phys, bnd);
2393}
2394
2396 int i, Array<OneD, Array<OneD, NekDouble>> &normals)
2397{
2398 v_GetBoundaryNormals(i, normals);
2399}
2400
2401inline std::vector<bool> &ExpList::GetLeftAdjacentTraces(void)
2402{
2403 return v_GetLeftAdjacentTraces();
2404}
2405
2407
2408} // namespace Nektar::MultiRegions
2409
2410#endif // EXPLIST_H
#define ASSERTL0(condition, msg)
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
#define MULTI_REGIONS_EXPORT
Base class for all multi-elemental spectral/hp expansions.
Definition ExpList.h:99
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition ExpList.h:1129
Array< OneD, NekDouble > & UpdateCoeffs()
This function returns (a reference to) the array (implemented as m_coeffs) containing all local expa...
Definition ExpList.h:2176
void AddTraceJacToElmtJac(const Array< OneD, const DNekMatSharedPtr > &FwdMat, const Array< OneD, const DNekMatSharedPtr > &BwdMat, Array< OneD, DNekMatSharedPtr > &fieldMat)
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:1720
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)
const Array< OneD, const std::shared_ptr< ExpList > > & GetBndCondExpansions()
Definition ExpList.h:2193
std::shared_ptr< ExpList > & GetTrace()
Definition ExpList.h:2234
virtual std::shared_ptr< ExpList > & v_GetTrace()
void WriteTecplotConnectivity(std::ostream &outfile, int expansion=-1)
Definition ExpList.h:411
virtual void v_PeriodicBwdCopy(const Array< OneD, const NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
int GetNcoeffs(void) const
Returns the total number of local degrees of freedom .
Definition ExpList.h:1578
std::shared_ptr< GlobalLinSys > GenGlobalLinSys(const GlobalLinSysKey &mkey, const std::shared_ptr< AssemblyMapCG > &locToGloMap)
This operation constructs the global linear system of type mkey.
const DNekScalBlkMatSharedPtr & GetBlockMatrix(const GlobalMatrixKey &gkey)
LibUtilities::NekManager< GlobalLinSysKey, GlobalLinSys > & GetGlobalLinSysManager(void)
ExpList(const ExpansionType Type=eNoType)
The default constructor using a type.
virtual void v_ImposeNeumannConditions(Array< OneD, NekDouble > &outarray)
void SetCoeffs(int i, NekDouble val)
Set the i th coefficiient in m_coeffs to value val.
Definition ExpList.h:2001
virtual void v_FillBndCondFromField(const Array< OneD, NekDouble > coeffs)
ExpList(const ExpList &in, const std::vector< unsigned int > &eIDs, const bool DeclareCoeffPhysArrays=true, const Collections::ImplementationType ImpType=Collections::eNoImpType)
Constructor copying only elements defined in eIds.
virtual ~ExpList()
The default destructor.
void CreateCollections(Collections::ImplementationType ImpType=Collections::eNoImpType)
void ExponentialFilter(Array< OneD, NekDouble > &array, const NekDouble alpha, const NekDouble exponent, const NekDouble cutoff)
virtual Array< OneD, const NekDouble > v_HomogeneousEnergy(void)
virtual void v_AddFwdBwdTraceIntegral(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &outarray)
virtual void v_WriteTecplotHeader(std::ostream &outfile, std::string var="")
virtual void v_SetHomoLen(const NekDouble lhom)
void MultiplyByMassMatrix(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition ExpList.h:298
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:2278
virtual void v_RotLocalBwdTrace(Array< OneD, Array< OneD, NekDouble > > &Bwd)
void AddRightIPTPhysDerivBase(const int dir, const Array< OneD, const DNekMatSharedPtr > ElmtJacQuad, Array< OneD, DNekMatSharedPtr > ElmtJacCoef)
void ExtractPhysToBndElmt(int i, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &bndElmt)
Definition ExpList.h:2381
virtual void v_SetUpPhysNormals()
void Reset()
Reset geometry information and reset matrices.
Definition ExpList.h:390
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...
void GetFieldDefinitions(std::vector< LibUtilities::FieldDefinitionsSharedPtr > &fielddef)
Definition ExpList.h:922
BlockMatrixMapShPtr m_blockMat
Definition ExpList.h:1178
virtual void v_UnsetGlobalLinSys(GlobalLinSysKey, bool)
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
multiply the metric jacobi and quadrature weights
void SetPhys(int i, NekDouble val)
Set the i th value of m_phys to value val.
Definition ExpList.h:1680
void ExtractTracePhys(Array< OneD, NekDouble > &outarray)
Definition ExpList.h:2331
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)
LocalRegions::ExpansionSharedPtr & GetExp(const Array< OneD, const NekDouble > &gloCoord)
This function returns (a shared pointer to) the local elemental expansion containing the arbitrary po...
void GetBoundaryToElmtMap(Array< OneD, int > &ElmtID, Array< OneD, int > &EdgeID)
Definition ExpList.h:2362
void PeriodicBwdCopy(const Array< OneD, const NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
Copy and fill the Periodic boundaries.
Definition ExpList.h:2302
void WriteVtkPieceFooter(std::ostream &outfile, int expansion)
virtual LibUtilities::TranspositionSharedPtr v_GetTransposition(void)
std::shared_ptr< DNekMat > GenGlobalMatrixFull(const GlobalLinSysKey &mkey, const std::shared_ptr< AssemblyMapCG > &locToGloMap)
std::shared_ptr< InterfaceMapDG > & GetInterfaceMap(void)
Definition ExpList.h:2242
void GetBwdWeight(Array< OneD, NekDouble > &weightAver, Array< OneD, NekDouble > &weightJump)
Get the weight value for boundary conditions for boundary average and jump calculations.
virtual void v_WriteTecplotField(std::ostream &outfile, int expansion)
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)
virtual std::map< int, RobinBCInfoSharedPtr > v_GetRobinBCInfo(void)
virtual void v_ExtractTracePhys(Array< OneD, NekDouble > &outarray)
virtual void v_WriteVtkPieceData(std::ostream &outfile, int expansion, std::string var)
void InitialiseExpVector(const SpatialDomains::ExpansionInfoMap &expmap)
Define a list of elements using the geometry and basis key information in expmap;.
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...
void RotLocalBwdTrace(Array< OneD, Array< OneD, NekDouble > > &Bwd)
Rotate local Bwd trace across a rotational interface when the flow is perpendicular to the rotation a...
Definition ExpList.h:2317
virtual LibUtilities::BasisSharedPtr v_GetHomogeneousBasis(void)
Definition ExpList.h:1530
virtual void v_ExtractPhysToBndElmt(const int i, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &bndElmt)
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:2017
const Array< OneD, const Array< OneD, NekDouble > > & GetGridVelocity()
Definition ExpList.h:1088
virtual void v_IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void GetBndElmtExpansion(int i, std::shared_ptr< ExpList > &result, const bool DeclareCoeffPhysArrays=true)
Definition ExpList.h:2367
void SetPhysArray(Array< OneD, NekDouble > &inarray)
Sets the array m_phys.
Definition ExpList.h:1699
const Collections::CollectionVector & GetCollections() const
This function returns collections.
Definition ExpList.h:1071
virtual void v_ClearGlobalLinSysManager(void)
void FillBndCondFromField(const Array< OneD, NekDouble > coeffs)
Fill Bnd Condition expansion from the values stored in expansion.
Definition ExpList.h:2033
const Array< OneD, const std::pair< int, int > > & GetCoeffsToElmt() const
Get m_coeffs to elemental value map.
Definition ExpList.h:2199
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:2111
ExpList(const LibUtilities::SessionReaderSharedPtr &pSession, const LocalRegions::ExpansionVector &locexp, const SpatialDomains::MeshGraphSharedPtr &graph, const bool DeclareCoeffPhysArrays, const std::string variable, const Collections::ImplementationType ImpType=Collections::eNoImpType)
Generate an trace ExpList from a meshgraph graph and session file.
std::shared_ptr< ExpList > & UpdateBndCondExpansion(int i)
Definition ExpList.h:2216
void NormVectorIProductWRTBase(Array< OneD, const NekDouble > &V1, Array< OneD, const NekDouble > &V2, Array< OneD, NekDouble > &outarray, int BndID)
Definition ExpList.h:1960
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:2159
virtual void v_SetBndCondBwdWeight(const int index, const NekDouble value)
virtual void v_NormVectorIProductWRTBase(Array< OneD, Array< OneD, NekDouble > > &V, Array< OneD, NekDouble > &outarray)
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition ExpList.h:2119
virtual void v_CurlCurl(Array< OneD, Array< OneD, NekDouble > > &Vel, Array< OneD, Array< OneD, NekDouble > > &Q)
SpatialDomains::EntityHolder1D m_holder
Pointer holder for PulseWaveSolver.
Definition ExpList.h:1095
std::shared_ptr< LibUtilities::SessionReader > GetSession() const
Returns the session object.
Definition ExpList.h:981
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...
void GetMatIpwrtDeriveBase(const TensorOfArray3D< NekDouble > &inarray, Array< OneD, DNekMatSharedPtr > &mtxPerVar)
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)
virtual void v_PhysDeriv(Direction edir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d)
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:1729
void SetCoeffsArray(Array< OneD, NekDouble > &inarray)
Set the m_coeffs array to inarray.
Definition ExpList.h:2005
LibUtilities::TranspositionSharedPtr GetTransposition(void)
This function returns the transposition class associated with the homogeneous expansion.
Definition ExpList.h:604
virtual void v_ImposeDirichletConditions(Array< OneD, NekDouble > &outarray)
void SetBndCondBwdWeight(const int index, const NekDouble value)
Set the weight value for boundary conditions.
Definition ExpList.h:2212
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)
void SetHomoLen(const NekDouble lhom)
This function sets the Width of homogeneous direction associated with the homogeneous expansion.
Definition ExpList.h:618
virtual void v_AppendFieldData(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata)
virtual void v_IProductWRTDerivBase(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_WriteTecplotConnectivity(std::ostream &outfile, int expansion)
virtual std::shared_ptr< ExpList > & v_UpdateBndCondExpansion(int i)
int GetNpoints(void) const
Returns the total number of quadrature points m_npoints .
Definition ExpList.h:1658
SpatialDomains::MeshGraphSharedPtr GetGraph()
Definition ExpList.h:990
virtual const std::shared_ptr< LocTraceToTraceMap > & v_GetLocTraceToTraceMap(void) const
void GetFwdBwdTracePhys(Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
Definition ExpList.h:2265
int GetElmtToExpId(int elmtId)
This function returns the index inside m_exp for a given geom id.
Definition ExpList.h:1059
const Array< OneD, const int > & GetTraceBndMap(void)
Definition ExpList.h:2246
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)
Array< OneD, Array< OneD, NekDouble > > m_gridVelocity
Grid velocity at quadrature points.
Definition ExpList.h:1187
int GetPoolCount(std::string)
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:1757
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:597
void FwdTransBndConstrained(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition ExpList.h:1766
virtual size_t v_GetNumElmts(void)
Definition ExpList.h:1214
virtual std::shared_ptr< AssemblyMapDG > & v_GetTraceMap()
Array< OneD, int > m_coeff_offset
Offset of elemental data into the array m_coeffs.
Definition ExpList.h:1173
int EvalBasisNumModesMax(void) const
Evaulates the maximum number of modes in the elemental basis order over all elements.
Definition ExpList.h:1592
void AddFwdBwdTraceIntegral(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &outarray)
Definition ExpList.h:2259
virtual void v_GetLocTraceFromTracePts(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &locTraceFwd, Array< OneD, NekDouble > &locTraceBwd)
std::shared_ptr< AssemblyMapDG > & GetTraceMap(void)
Definition ExpList.h:2238
void WriteTecplotField(std::ostream &outfile, int expansion=-1)
Definition ExpList.h:407
virtual void v_AddTraceQuadPhysToOffDiag(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &field)
size_t GetNumElmts(void)
This function returns the number of elements in the expansion which may be different for a homogeoeno...
Definition ExpList.h:660
virtual LibUtilities::NekManager< GlobalLinSysKey, GlobalLinSys > & v_GetGlobalLinSysManager(void)
Array< OneD, std::pair< int, int > > m_coeffsToElmt
m_coeffs to elemental value map
Definition ExpList.h:1177
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)
NekDouble PhysEvaluate(const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &phys)
const std::unordered_map< int, int > & GetElmtToExpId(void)
This function returns the map of index inside m_exp to geom id.
Definition ExpList.h:1052
virtual void v_GetPeriodicEntities(PeriodicMap &periodicVerts, PeriodicMap &periodicEdges, PeriodicMap &periodicFaces)
void PhysDirectionalDeriv(const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition ExpList.h:1895
virtual NekDouble v_Integral(const Array< OneD, const NekDouble > &inarray)
virtual GlobalLinSysKey v_LinearAdvectionReactionSolve(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)
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void GetBoundaryNormals(int i, Array< OneD, Array< OneD, NekDouble > > &normals)
Definition ExpList.h:2395
virtual void v_RotLocalBwdDeriveTrace(TensorOfArray3D< NekDouble > &Bwd)
virtual void v_LocalToGlobal(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool UseComm)
void ExtractFileBCs(const std::string &fileName, LibUtilities::CommSharedPtr comm, const std::string &varName, const std::shared_ptr< ExpList > locExpList)
virtual const Array< OneD, const std::shared_ptr< ExpList > > & v_GetBndCondExpansions(void)
void GetMovingFrames(const SpatialDomains::GeomMMF MMFdir, const Array< OneD, const NekDouble > &CircCentre, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition ExpList.h:1860
const LocTraceToTraceMapSharedPtr & GetLocTraceToTraceMap() const
Definition ExpList.h:2204
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:581
void SetHomo1DSpecVanVisc(Array< OneD, NekDouble > visc)
This function sets the Spectral Vanishing Viscosity in homogeneous1D expansion.
Definition ExpList.h:588
virtual void v_HomogeneousFwdTrans(const int npts, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
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...
virtual void v_GlobalToLocal(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
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:1843
void ResetMatrices()
Reset matrices.
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:1607
void PeriodicBwdRot(Array< OneD, Array< OneD, NekDouble > > &Bwd)
Rotate Bwd trace for rotational periodicity boundaries when the flow is perpendicular to the rotation...
Definition ExpList.h:2307
virtual void v_HomogeneousBwdTrans(const int npts, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
virtual void v_GetBndElmtExpansion(int i, std::shared_ptr< ExpList > &result, const bool DeclareCoeffPhysArrays)
virtual std::shared_ptr< InterfaceMapDG > & v_GetInterfaceMap()
void UnsetGlobalLinSys(GlobalLinSysKey, bool)
ExpList(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::CompositeMap &domain, const SpatialDomains::MeshGraphSharedPtr &graph, const bool DeclareCoeffPhysArrays=true, const std::string variable="DefaultVar", bool SetToOneSpaceDimension=false, const LibUtilities::CommSharedPtr comm=LibUtilities::CommSharedPtr(), const Collections::ImplementationType ImpType=Collections::eNoImpType)
Constructor based on domain information only for 1D & 2D boundary conditions.
void DealiasedProd(const int num_dofs, const Array< OneD, NekDouble > &inarray1, const Array< OneD, NekDouble > &inarray2, Array< OneD, NekDouble > &outarray)
Definition ExpList.h:1931
void GetNormals(Array< OneD, Array< OneD, NekDouble > > &normals)
Definition ExpList.h:2250
void ExtractPhysToBnd(int i, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &bnd)
Definition ExpList.h:2388
virtual void v_GlobalToLocal(void)
std::vector< bool > m_collectionsDoInit
Vector of bools to act as an initialise on first call flag.
Definition ExpList.h:1167
NekDouble GetCoeff(int i)
Get the i th value (coefficient) of m_coeffs.
Definition ExpList.h:2091
bool m_physState
The state of the array m_phys.
Definition ExpList.h:1153
virtual void v_SmoothField(Array< OneD, NekDouble > &field)
void ExtractCoeffsFromFile(const std::string &fileName, LibUtilities::CommSharedPtr comm, const std::string &varName, Array< OneD, NekDouble > &coeffs)
NekDouble Integral()
Calculates the error of the global spectral/hp element approximation.
Definition ExpList.h:550
std::shared_ptr< ExpList > & GetPlane(int n)
Definition ExpList.h:999
void AddRightIPTBaseMatrix(const Array< OneD, const DNekMatSharedPtr > ElmtJacQuad, Array< OneD, DNekMatSharedPtr > ElmtJacCoef)
virtual std::vector< bool > & v_GetLeftAdjacentTraces(void)
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:2290
virtual void v_PeriodicBwdRot(Array< OneD, Array< OneD, NekDouble > > &Bwd)
std::shared_ptr< ExpList > GetSharedThisPtr()
Returns a shared pointer to the current object.
Definition ExpList.h:976
int GetShapeDimension()
This function returns the dimension of the shape of the element eid.
Definition ExpList.h:1985
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.
virtual void v_ExtractCoeffsToCoeffs(const std::shared_ptr< ExpList > &fromExpList, const Array< OneD, const NekDouble > &fromCoeffs, Array< OneD, NekDouble > &toCoeffs)
std::vector< int > m_coll_phys_offset
Offset of elemental data into the array m_phys.
Definition ExpList.h:1171
LibUtilities::CommSharedPtr m_comm
Communicator.
Definition ExpList.h:1102
void WriteVtkHeader(std::ostream &outfile)
bool GetPhysState(void) const
This function indicates whether the array of physical values (implemented as m_phys) is filled or no...
Definition ExpList.h:1713
std::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition ExpList.h:1164
virtual NekDouble v_VectorFlux(const Array< OneD, Array< OneD, NekDouble > > &inarray)
ExpList(const ExpList &in, const bool DeclareCoeffPhysArrays=true)
The copy constructor.
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Definition ExpList.h:1109
virtual void v_AddTraceQuadPhysToField(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &field)
void GenerateElementVector(const int ElementID, const NekDouble scalar1, const NekDouble scalar2, Array< OneD, NekDouble > &outarray)
NekDouble VectorFlux(const Array< OneD, Array< OneD, NekDouble > > &inarray)
Definition ExpList.h:575
void AddTraceIntegralToOffDiag(const Array< OneD, const NekDouble > &FwdFlux, const Array< OneD, const NekDouble > &BwdFlux, Array< OneD, NekDouble > &outarray)
Definition ExpList.h:1031
virtual void v_ExtractDataToCoeffs(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, std::string &field, Array< OneD, NekDouble > &coeffs, std::unordered_map< int, int > zIdToPlane)
const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > & GetBndConditions()
Definition ExpList.h:2342
virtual Array< OneD, SpatialDomains::BoundaryConditionShPtr > & v_UpdateBndConditions()
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:2139
virtual void v_DealiasedProd(const int num_dofs, const Array< OneD, NekDouble > &inarray1, const Array< OneD, NekDouble > &inarray2, Array< OneD, NekDouble > &outarray)
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:1706
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:647
virtual void v_WriteVtkPieceHeader(std::ostream &outfile, int expansion, int istrip)
void CurlCurl(Array< OneD, Array< OneD, NekDouble > > &Vel, Array< OneD, Array< OneD, NekDouble > > &Q)
Definition ExpList.h:1905
std::map< int, RobinBCInfoSharedPtr > GetRobinBCInfo()
Definition ExpList.h:905
virtual void v_GetBoundaryNormals(int i, Array< OneD, Array< OneD, NekDouble > > &normals)
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:1804
int Get_coll_coeff_offset(int n) const
Definition ExpList.h:1076
void EvaluateBoundaryConditions(const NekDouble time=0.0, const std::string varName="", const NekDouble=NekConstants::kNekUnsetDouble, const NekDouble=NekConstants::kNekUnsetDouble)
Definition ExpList.h:2351
std::vector< LibUtilities::FieldDefinitionsSharedPtr > GetFieldDefinitions()
Definition ExpList.h:917
virtual void v_PeriodicDeriveBwdRot(TensorOfArray3D< NekDouble > &Bwd)
Array< OneD, NekDouble > & UpdatePhys()
This function returns (a reference to) the array (implemented as m_phys) containing the function ev...
Definition ExpList.h:2186
void ImposeRobinConditions(Array< OneD, NekDouble > &outarray)
Add Robin Boundary Condition forcing to Array.
Definition ExpList.h:2029
virtual const Array< OneD, const int > & v_GetTraceBndMap()
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 ...
void ExtractElmtToBndPhys(int i, const Array< OneD, NekDouble > &elmt, Array< OneD, NekDouble > &boundary)
Definition ExpList.h:2374
virtual void v_GetBoundaryToElmtMap(Array< OneD, int > &ElmtID, Array< OneD, int > &EdgeID)
void HomogeneousBwdTrans(const int npts, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
Definition ExpList.h:1922
void SetWaveSpace(const bool wavespace)
Sets the wave space to the one of the possible configuration true or false.
Definition ExpList.h:1666
const std::shared_ptr< LocalRegions::ExpansionVector > GetExp() const
This function returns the vector of elements in the expansion.
Definition ExpList.h:2151
SpatialDomains::MeshGraphSharedPtr m_graph
Mesh associated with this expansion list.
Definition ExpList.h:1106
virtual void v_AddTraceIntegral(const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
virtual void v_GetCoords(const int eid, Array< OneD, NekDouble > &xc0, Array< OneD, NekDouble > &xc1, Array< OneD, NekDouble > &xc2)
virtual void v_GetBCValues(Array< OneD, NekDouble > &BndVals, const Array< OneD, NekDouble > &TotField, int BndID)
virtual NekDouble v_L2(const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray)
void GlobalToLocal(void)
Scatters from the global coefficients to the local coefficients .
Definition ExpList.h:2052
const Array< OneD, const NekDouble > & GetBndCondBwdWeight()
Get the weight value for boundary conditions.
Definition ExpList.h:2208
void WriteVtkFooter(std::ostream &outfile)
virtual void v_PhysDeriv(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d)
void MultiplyByBlockMatrix(const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_WriteTecplotZone(std::ostream &outfile, int expansion)
virtual void v_FwdTransBndConstrained(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void GetElmtNormalLength(Array< OneD, NekDouble > &lengthsFwd, Array< OneD, NekDouble > &lengthsBwd)
Get the length of elements in boundary normal direction.
virtual void v_ExtractPhysToBnd(const int i, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &bnd)
virtual void v_ExtractTracePhys(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool gridVelocity=false)
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:1633
ExpList(SpatialDomains::PointGeom *geom)
Specialised constructors for 0D Expansions Wrapper around LocalRegion::PointExp - used in PrePacing....
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:1817
virtual void v_AppendFieldData(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, Array< OneD, NekDouble > &coeffs)
virtual void v_FillBndCondFromField(const int nreg, const Array< OneD, NekDouble > coeffs)
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:636
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:521
void PeriodicDeriveBwdRot(TensorOfArray3D< NekDouble > &Bwd)
Rotate Bwd trace derivative for rotational periodicity boundaries when the flow is perpendicular to t...
Definition ExpList.h:2312
virtual void v_FwdTransLocalElmt(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
std::vector< int > m_coll_coeff_offset
Offset of elemental data into the array m_coeffs.
Definition ExpList.h:1169
void IProductWRTDirectionalDerivBase(const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void WriteVtkPieceData(std::ostream &outfile, int expansion, std::string var="v")
Definition ExpList.h:424
void FillBwdWithBwdWeight(Array< OneD, NekDouble > &weightave, Array< OneD, NekDouble > &weightjmp)
Fill Bwd with boundary conditions.
Definition ExpList.h:2297
void SetCoeff(int i, NekDouble val)
Set the i th coefficiient in m_coeffs to value val.
Definition ExpList.h:1993
void GetMatIpwrtDeriveBase(const Array< OneD, const Array< OneD, NekDouble > > &inarray, const int nDirctn, Array< OneD, DNekMatSharedPtr > &mtxPerVar)
void LocalToGlobal(bool useComm=true)
Gathers the global coefficients from the local coefficients .
Definition ExpList.h:2042
void GetPeriodicEntities(PeriodicMap &periodicVerts, PeriodicMap &periodicEdges, PeriodicMap &periodicFaces=NullPeriodicMap)
Definition ExpList.h:910
NekDouble Integral(const Array< OneD, const NekDouble > &inarray)
Definition ExpList.h:571
virtual std::shared_ptr< ExpList > & v_GetPlane(int n)
virtual void v_AddTraceIntegralToOffDiag(const Array< OneD, const NekDouble > &FwdFlux, const Array< OneD, const NekDouble > &BwdFlux, Array< OneD, NekDouble > &outarray)
ExpList(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &graph, const bool DeclareCoeffPhysArrays=true, const std::string &var="DefaultVar", const Collections::ImplementationType ImpType=Collections::eNoImpType)
Generate an ExpList from a meshgraph graph and session file.
Collections::CollectionVector m_collections
Definition ExpList.h:1165
virtual const Array< OneD, const NekDouble > & v_GetBndCondBwdWeight()
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:2284
virtual void v_MultiplyByInvMassMatrix(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual int v_GetPoolCount(std::string)
int GetExpIndex(const Array< OneD, const NekDouble > &gloCoords, Array< OneD, NekDouble > &locCoords, NekDouble tol=0.0, bool returnNearestElmt=false, int cachedId=-1, NekDouble maxDistance=1e6)
virtual void v_PhysDirectionalDeriv(const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void ApplyGeomInfo()
Apply geometry information to each expansion.
void AddTraceIntegral(const Array< OneD, const NekDouble > &Fn, Array< OneD, NekDouble > &outarray)
Definition ExpList.h:2254
virtual void v_GetFieldDefinitions(std::vector< LibUtilities::FieldDefinitionsSharedPtr > &fielddef)
std::unordered_map< int, int > m_elmtToExpId
Mapping from geometry ID of element to index inside m_exp.
Definition ExpList.h:1185
void FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition ExpList.h:1749
virtual void v_ImposeRobinConditions(Array< OneD, NekDouble > &outarray)
virtual void v_GetFwdBwdTracePhys(const Array< OneD, const NekDouble > &field, Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd, bool FillBnd=true, bool PutFwdInBwdOnBCs=false, bool DoExchange=true)
virtual void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_ExtractElmtToBndPhys(const int i, const Array< OneD, NekDouble > &elmt, Array< OneD, NekDouble > &boundary)
void SmoothField(Array< OneD, NekDouble > &field)
Smooth a field across elements.
Definition ExpList.h:1775
virtual void v_GetFwdBwdTracePhys(Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
virtual const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > & v_GetBndConditions()
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:627
virtual std::vector< LibUtilities::FieldDefinitionsSharedPtr > v_GetFieldDefinitions(void)
LibUtilities::SessionReaderSharedPtr m_session
Session.
Definition ExpList.h:1104
static SpatialDomains::BoundaryConditionShPtr GetBoundaryCondition(const SpatialDomains::BoundaryConditionCollection &collection, unsigned int index, const std::string &variable)
virtual void v_FillBwdWithBoundCond(const Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd, bool PutFwdInBwdOnBCs)
GlobalLinSysKey LinearAdvectionReactionSolve(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:1829
virtual void v_GetMovingFrames(const SpatialDomains::GeomMMF MMFdir, const Array< OneD, const NekDouble > &CircCentre, Array< OneD, Array< OneD, NekDouble > > &outarray)
void MultiplyByInvMassMatrix(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition ExpList.h:1795
virtual void v_GetNormals(Array< OneD, Array< OneD, NekDouble > > &normals)
Populate normals with the normals of all expansions.
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:937
virtual Array< OneD, const unsigned int > v_GetZIDs(void)
void DivideByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Divided by the metric jacobi and quadrature weights.
void ImposeDirichletConditions(Array< OneD, NekDouble > &outarray)
Impose Dirichlet Boundary Conditions onto Array.
Definition ExpList.h:2021
virtual void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_Upwind(const Array< OneD, const NekDouble > &Vn, const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &Upwind)
Array< OneD, int > m_phys_offset
Offset of elemental data into the array m_phys.
Definition ExpList.h:1175
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:2227
void RotLocalBwdDeriveTrace(TensorOfArray3D< NekDouble > &Bwd)
Rotate local Bwd trace derivatives across a rotational interface when the flow is perpendicular to th...
Definition ExpList.h:2322
Array< OneD, SpatialDomains::BoundaryConditionShPtr > & UpdateBndConditions()
Definition ExpList.h:2347
virtual void v_GetCoords(Array< OneD, NekDouble > &coord_0, Array< OneD, NekDouble > &coord_1=NullNekDouble1DArray, Array< OneD, NekDouble > &coord_2=NullNekDouble1DArray)
const DNekScalBlkMatSharedPtr GenBlockMatrix(const GlobalMatrixKey &gkey)
This function assembles the block diagonal matrix of local matrices of the type mtype.
void PhysDeriv(Direction edir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d)
Definition ExpList.h:1886
int Get_coll_phys_offset(int n) const
Definition ExpList.h:1081
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 ...
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:1787
virtual NekDouble v_GetHomoLen(void)
NekDouble GetHomoLen(void)
This function returns the Width of homogeneous direction associated with the homogeneous expansion.
Definition ExpList.h:611
virtual void v_LocalToGlobal(bool UseComm)
virtual void v_PhysGalerkinProjection1DScaled(const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
ExpList(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::ExpansionInfoMap &expansions, const bool DeclareCoeffPhysArrays=true, const Collections::ImplementationType ImpType=Collections::eNoImpType)
Sets up a list of local expansions based on an expansion Map.
bool GetWaveSpace(void) const
This function returns the third direction expansion condition, which can be in wave space (coefficien...
Definition ExpList.h:1674
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:2166
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 >())
void SetModifiedBasis(const bool modbasis)
Set Modified Basis for the stability analysis.
ExpansionType m_expType
Expansion type.
Definition ExpList.h:1097
void WriteTecplotHeader(std::ostream &outfile, std::string var="")
Definition ExpList.h:399
virtual void v_PhysInterp1DScaled(const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
const std::vector< bool > & GetLeftAdjacentFaces(void) const
Definition ExpList.h:2327
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:1941
std::shared_ptr< GlobalMatrix > GenGlobalMatrix(const GlobalMatrixKey &mkey, const std::shared_ptr< AssemblyMapCG > &locToGloMap)
Generates a global matrix from the given key and map.
void AppendFieldData(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata)
Append the element data listed in elements fielddef->m_ElementIDs onto fielddata.
Definition ExpList.h:930
std::shared_ptr< LibUtilities::Comm > GetComm() const
Returns the comm object.
Definition ExpList.h:986
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
Definition ExpList.h:1145
std::vector< bool > & GetLeftAdjacentTraces(void)
Definition ExpList.h:2401
ExpansionType GetExpType(void)
Returns the type of the expansion.
void AddTraceQuadPhysToOffDiag(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &field)
Definition ExpList.h:822
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.
void HomogeneousFwdTrans(const int npts, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
Definition ExpList.h:1913
virtual void v_Curl(Array< OneD, Array< OneD, NekDouble > > &Vel, Array< OneD, Array< OneD, NekDouble > > &Q)
void WriteVtkPieceHeader(std::ostream &outfile, int expansion, int istrip=0)
Definition ExpList.h:417
virtual void v_FillBwdWithBwdWeight(Array< OneD, NekDouble > &weightave, Array< OneD, NekDouble > &weightjmp)
virtual void v_NormVectorIProductWRTBase(Array< OneD, const NekDouble > &V1, Array< OneD, const NekDouble > &V2, Array< OneD, NekDouble > &outarray, int BndID)
void ImposeNeumannConditions(Array< OneD, NekDouble > &outarray)
Add Neumann Boundary Condition forcing to Array.
Definition ExpList.h:2025
int GetTotPoints(void) const
Returns the total number of quadrature points m_npoints .
Definition ExpList.h:1623
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
virtual Array< OneD, const unsigned int > v_GetYIDs(void)
LibUtilities::BasisSharedPtr GetHomogeneousBasis(void)
Definition ExpList.h:995
void WriteTecplotZone(std::ostream &outfile, int expansion=-1)
Definition ExpList.h:403
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 const std::vector< bool > & v_GetLeftAdjacentFaces(void) const
virtual void v_SetHomo1DSpecVanVisc(Array< OneD, NekDouble > visc)
Definition ExpList.h:1538
void GetBCValues(Array< OneD, NekDouble > &BndVals, const Array< OneD, NekDouble > &TotField, int BndID)
Definition ExpList.h:1951
int GetCoordim(int eid)
This function returns the dimension of the coordinates of the element eid.
Definition ExpList.h:1976
void GetDiagMatIpwrtBase(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, DNekMatSharedPtr > &mtxPerVar)
A global linear system.
Describes a matrix with ordering defined by a local to global map.
std::vector< Collection > CollectionVector
Definition Collection.h:132
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:1564
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:1565
static LibUtilities::NekManager< GlobalLinSysKey, GlobalLinSys > NullGlobalLinSysManager
Definition ExpList.h:1569
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:2406
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:251
std::map< int, BoundaryConditionMapShPtr > BoundaryConditionCollection
Definition Conditions.h:258
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition MeshGraph.h:217
std::map< int, ExpansionInfoShPtr > ExpansionInfoMap
Definition MeshGraph.h:184
std::map< int, CompositeSharedPtr > CompositeMap
Definition MeshGraph.h:179
std::map< ConstFactorType, NekDouble > ConstFactorMap
static VarCoeffMap NullVarCoeffMap
std::map< StdRegions::VarCoeffType, VarCoeffEntry > VarCoeffMap
std::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
static Array< OneD, NekDouble > NullNekDouble1DArray
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition Vmath.hpp:825