Nektar++
ExpList.h
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: ExpList.h
4//
5// For more information, please see: http://www.nektar.info
6//
7// The MIT License
8//
9// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// Permission is hereby granted, free of charge, to any person obtaining a
14// copy of this software and associated documentation files (the "Software"),
15// to deal in the Software without restriction, including without limitation
16// the rights to use, copy, modify, merge, publish, distribute, sublicense,
17// and/or sell copies of the Software, and to permit persons to whom the
18// Software is furnished to do so, subject to the following conditions:
19//
20// The above copyright notice and this permission notice shall be included
21// in all copies or substantial portions of the Software.
22//
23// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29// DEALINGS IN THE SOFTWARE.
30//
31// Description: Expansion list top class definition
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#ifndef NEKTAR_LIBS_MULTIREGIONS_EXPLIST_H
36#define NEKTAR_LIBS_MULTIREGIONS_EXPLIST_H
37
38#include <boost/core/ignore_unused.hpp>
39
54#include <tinyxml.h>
55
56namespace Nektar
57{
58namespace MultiRegions
59{
60
61// Forward declarations
62class ExpList;
63class GlobalLinSys;
64class AssemblyMapDG;
65class AssemblyMapCG;
66class GlobalLinSysKey;
67class GlobalMatrix;
68
70{
75 eN
76};
77
79{
88};
89
91
92/// A map between global matrix keys and their associated block
93/// matrices.
94typedef std::map<GlobalMatrixKey, DNekScalBlkMatSharedPtr> BlockMatrixMap;
95/// A shared pointer to a BlockMatrixMap.
96typedef std::shared_ptr<BlockMatrixMap> BlockMatrixMapShPtr;
97/// Shared pointer to an ExpList object.
98typedef std::shared_ptr<ExpList> ExpListSharedPtr;
99
100/// Base class for all multi-elemental spectral/hp expansions.
101class ExpList : public std::enable_shared_from_this<ExpList>
102{
103public:
104 /// The default constructor using a type
106
107 /// The copy constructor.
109 const bool DeclareCoeffPhysArrays = true);
110
111 /// The copy constructor.
113 const bool DeclareCoeffArrays = true,
114 const bool DeclarePhysArrays = true);
115
116 /// Constructor copying only elements defined in eIds.
118 const std::vector<unsigned int> &eIDs,
119 const bool DeclareCoeffPhysArrays = true,
120 const Collections::ImplementationType ImpType =
122
123 /// Generate an ExpList from a meshgraph \a graph and session file
127 const bool DeclareCoeffPhysArrays = true,
128 const std::string &var = "DefaultVar",
129 const Collections::ImplementationType ImpType =
131
132 /// Sets up a list of local expansions based on an expansion Map
135 const SpatialDomains::ExpansionInfoMap &expansions,
136 const bool DeclareCoeffPhysArrays = true,
137 const Collections::ImplementationType ImpType =
139
140 //---------------------------------------------------------
141 // Specialised constructors in ExpListConstructor.cpp
142 //---------------------------------------------------------
143 /// Specialised constructors for 0D Expansions
144 /// Wrapper around LocalRegion::PointExp - used in PrePacing.cpp
147
148 /// Generate expansions for the trace space expansions used in
149 /// DisContField.
152 const Array<OneD, const ExpListSharedPtr> &bndConstraint,
154 &bndCond,
155 const LocalRegions::ExpansionVector &locexp,
157 const LibUtilities::CommSharedPtr &comm,
158 const bool DeclareCoeffPhysArrays = true,
159 const std::string variable = "DefaultVar",
160 const Collections::ImplementationType ImpType =
162
163 /// Generate an trace ExpList from a meshgraph \a graph and session file
166 const LocalRegions::ExpansionVector &locexp,
168 const bool DeclareCoeffPhysArrays, const std::string variable,
169 const Collections::ImplementationType ImpType =
171
172 /// Constructor based on domain information only for 1D &
173 /// 2D boundary conditions
176 const SpatialDomains::CompositeMap &domain,
178 const bool DeclareCoeffPhysArrays = true,
179 const std::string variable = "DefaultVar",
180 bool SetToOneSpaceDimension = false,
182 const Collections::ImplementationType ImpType =
184
185 /// The default destructor.
187
188 /// Returns the total number of local degrees of freedom
189 /// \f$N_{\mathrm{eof}}=\sum_{e=1}^{{N_{\mathrm{el}}}}N^{e}_m\f$.
190 inline int GetNcoeffs(void) const;
191
192 /// Returns the total number of local degrees of freedom
193 /// for element eid
194 MULTI_REGIONS_EXPORT int GetNcoeffs(const int eid) const;
195
196 /// Returns the type of the expansion
198
199 /// Returns the type of the expansion
201
202 /// Evaulates the maximum number of modes in the elemental basis
203 /// order over all elements
204 inline int EvalBasisNumModesMax(void) const;
205
206 /// Returns the vector of the number of modes in the elemental
207 /// basis order over all elements.
209 void) const;
210
211 /// Returns the total number of quadrature points #m_npoints
212 /// \f$=Q_{\mathrm{tot}}\f$.
213 inline int GetTotPoints(void) const;
214
215 /// Returns the total number of quadrature points for eid's element
216 /// \f$=Q_{\mathrm{tot}}\f$.
217 inline int GetTotPoints(const int eid) const;
218
219 /// Returns the total number of quadrature points #m_npoints
220 /// \f$=Q_{\mathrm{tot}}\f$.
221 inline int GetNpoints(void) const;
222
223 /// Returns the total number of qudature points scaled by
224 /// the factor scale on each 1D direction
225 inline int Get1DScaledTotPoints(const NekDouble scale) const;
226
227 /// Sets the wave space to the one of the possible configuration
228 /// true or false
229 inline void SetWaveSpace(const bool wavespace);
230
231 /// Set Modified Basis for the stability analysis
232 inline void SetModifiedBasis(const bool modbasis);
233
234 /// This function returns the third direction expansion condition,
235 /// which can be in wave space (coefficient) or not
236 /// It is stored in the variable m_WaveSpace.
237 inline bool GetWaveSpace(void) const;
238
239 /// Set the \a i th value of \a m_phys to value \a val
240 inline void SetPhys(int i, NekDouble val);
241 /// Fills the array #m_phys
242 inline void SetPhys(const Array<OneD, const NekDouble> &inarray);
243 /// Sets the array #m_phys
244 inline void SetPhysArray(Array<OneD, NekDouble> &inarray);
245 /// This function manually sets whether the array of physical
246 /// values \f$\boldsymbol{u}_l\f$ (implemented as #m_phys) is
247 /// filled or not.
248 inline void SetPhysState(const bool physState);
249 /// This function indicates whether the array of physical values
250 /// \f$\boldsymbol{u}_l\f$ (implemented as #m_phys) is filled or
251 /// not.
252 inline bool GetPhysState(void) const;
253
254 /// multiply the metric jacobi and quadrature weights
256 const Array<OneD, const NekDouble> &inarray,
257 Array<OneD, NekDouble> &outarray);
258 /// Divided by the metric jacobi and quadrature weights
260 const Array<OneD, const NekDouble> &inarray,
261 Array<OneD, NekDouble> &outarray);
262 /// This function calculates the inner product of a function
263 /// \f$f(\boldsymbol{x})\f$ with respect to all \em local
264 /// expansion modes \f$\phi_n^e(\boldsymbol{x})\f$.
265 inline void IProductWRTBase(const Array<OneD, const NekDouble> &inarray,
266 Array<OneD, NekDouble> &outarray);
267 /// This function calculates the inner product of a function
268 /// \f$f(\boldsymbol{x})\f$ with respect to the derivative (in
269 /// direction \param dir) of all \em local expansion modes
270 /// \f$\phi_n^e(\boldsymbol{x})\f$.
272 const int dir, const Array<OneD, const NekDouble> &inarray,
273 Array<OneD, NekDouble> &outarray);
274
276 const Array<OneD, const NekDouble> &direction,
277 const Array<OneD, const NekDouble> &inarray,
278 Array<OneD, NekDouble> &outarray);
279
280 /// This function calculates the inner product of a
281 /// function \f$f(\boldsymbol{x})\f$ with respect to the
282 /// derivative of all \em local expansion modes
283 /// \f$\phi_n^e(\boldsymbol{x})\f$.
285 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
286 Array<OneD, NekDouble> &outarray);
287 /// This function elementally evaluates the forward transformation
288 /// of a function \f$u(\boldsymbol{x})\f$ onto the global
289 /// spectral/hp expansion.
290 inline void FwdTransLocalElmt(const Array<OneD, const NekDouble> &inarray,
291 Array<OneD, NekDouble> &outarray);
292 ///
293 inline void FwdTrans(const Array<OneD, const NekDouble> &inarray,
294 Array<OneD, NekDouble> &outarray);
296 const NekDouble alpha,
297 const NekDouble exponent,
298 const NekDouble cutoff);
299 /// This function elementally mulplies the coefficient space of
300 /// Sin my the elemental inverse of the mass matrix.
302 const Array<OneD, const NekDouble> &inarray,
303 Array<OneD, NekDouble> &outarray);
304 inline void MultiplyByInvMassMatrix(
305 const Array<OneD, const NekDouble> &inarray,
306 Array<OneD, NekDouble> &outarray);
308 const Array<OneD, const NekDouble> &inarray,
309 Array<OneD, NekDouble> &outarray)
310 {
312 BwdTrans(inarray, tmp);
313 IProductWRTBase(tmp, outarray);
314 }
315 /// Smooth a field across elements
316 inline void SmoothField(Array<OneD, NekDouble> &field);
317
318 /// Solve helmholtz problem
320 const Array<OneD, const NekDouble> &inarray,
321 Array<OneD, NekDouble> &outarray,
324 const MultiRegions::VarFactorsMap &varfactors =
327 const bool PhysSpaceForcing = true);
328
329 /// Solve Advection Diffusion Reaction
331 const Array<OneD, const NekDouble> &inarray,
332 Array<OneD, NekDouble> &outarray,
335 const MultiRegions::VarFactorsMap &varfactors =
338 const bool PhysSpaceForcing = true);
339
340 /// Solve Advection Diffusion Reaction
342 const Array<OneD, Array<OneD, NekDouble>> &velocity,
343 const Array<OneD, const NekDouble> &inarray,
344 Array<OneD, NekDouble> &outarray, const NekDouble lambda,
346 ///
348 const Array<OneD, const NekDouble> &inarray,
349 Array<OneD, NekDouble> &outarray);
350 /// This function elementally evaluates the backward transformation
351 /// of the global spectral/hp element expansion.
352 inline void BwdTrans(const Array<OneD, const NekDouble> &inarray,
353 Array<OneD, NekDouble> &outarray);
354 /// This function calculates the coordinates of all the elemental
355 /// quadrature points \f$\boldsymbol{x}_i\f$.
356 inline void GetCoords(
357 Array<OneD, NekDouble> &coord_0,
360
361 inline void GetCoords(
362 const int eid, Array<OneD, NekDouble> &coord_0,
365
366 // Homogeneous transforms
367 inline void HomogeneousFwdTrans(const int npts,
368 const Array<OneD, const NekDouble> &inarray,
369 Array<OneD, NekDouble> &outarray,
370 bool Shuff = true, bool UnShuff = true);
371 inline void HomogeneousBwdTrans(const int npts,
372 const Array<OneD, const NekDouble> &inarray,
373 Array<OneD, NekDouble> &outarray,
374 bool Shuff = true, bool UnShuff = true);
375 inline void DealiasedProd(const int num_dofs,
376 const Array<OneD, NekDouble> &inarray1,
377 const Array<OneD, NekDouble> &inarray2,
378 Array<OneD, NekDouble> &outarray);
379 inline void DealiasedDotProd(
380 const int num_dofs, const Array<OneD, Array<OneD, NekDouble>> &inarray1,
381 const Array<OneD, Array<OneD, NekDouble>> &inarray2,
382 Array<OneD, Array<OneD, NekDouble>> &outarray);
383 inline void GetBCValues(Array<OneD, NekDouble> &BndVals,
384 const Array<OneD, NekDouble> &TotField, int BndID);
387 Array<OneD, NekDouble> &outarray,
388 int BndID);
389 inline void NormVectorIProductWRTBase(
391 Array<OneD, NekDouble> &outarray);
392 /// Apply geometry information to each expansion.
394 /// Reset geometry information and reset matrices
396 {
397 v_Reset();
398 }
399 // not sure we really need these in ExpList
400 void WriteTecplotHeader(std::ostream &outfile, std::string var = "")
401 {
402 v_WriteTecplotHeader(outfile, var);
403 }
404 void WriteTecplotZone(std::ostream &outfile, int expansion = -1)
405 {
406 v_WriteTecplotZone(outfile, expansion);
407 }
408 void WriteTecplotField(std::ostream &outfile, int expansion = -1)
409 {
410 v_WriteTecplotField(outfile, expansion);
411 }
412 void WriteTecplotConnectivity(std::ostream &outfile, int expansion = -1)
413 {
414 v_WriteTecplotConnectivity(outfile, expansion);
415 }
416 MULTI_REGIONS_EXPORT void WriteVtkHeader(std::ostream &outfile);
417 MULTI_REGIONS_EXPORT void WriteVtkFooter(std::ostream &outfile);
418 void WriteVtkPieceHeader(std::ostream &outfile, int expansion,
419 int istrip = 0)
420 {
421 v_WriteVtkPieceHeader(outfile, expansion, istrip);
422 }
423 MULTI_REGIONS_EXPORT void WriteVtkPieceFooter(std::ostream &outfile,
424 int expansion);
425 void WriteVtkPieceData(std::ostream &outfile, int expansion,
426 std::string var = "v")
427 {
428 v_WriteVtkPieceData(outfile, expansion, var);
429 }
430
431 /// This function returns the dimension of the coordinates of the
432 /// element \a eid.
433 // inline
434 MULTI_REGIONS_EXPORT int GetCoordim(int eid);
435
436 /// Set the \a i th coefficiient in \a m_coeffs to value \a val
437 inline void SetCoeff(int i, NekDouble val);
438 /// Set the \a i th coefficiient in #m_coeffs to value \a val
439 inline void SetCoeffs(int i, NekDouble val);
440 /// Set the #m_coeffs array to inarray
441 inline void SetCoeffsArray(Array<OneD, NekDouble> &inarray);
442 /// This function returns the dimension of the shape of the
443 /// element \a eid.
444 // inline
446 /// This function returns (a reference to) the array
447 /// \f$\boldsymbol{\hat{u}}_l\f$ (implemented as #m_coeffs)
448 /// containing all local expansion coefficients.
449 inline const Array<OneD, const NekDouble> &GetCoeffs() const;
450 /// Impose Dirichlet Boundary Conditions onto Array
452 /// Fill Bnd Condition expansion from the values stored in expansion
453 inline void FillBndCondFromField(const Array<OneD, NekDouble> coeffs);
454 /// Fill Bnd Condition expansion in nreg from the values
455 /// stored in expansion
456 inline void FillBndCondFromField(const int nreg,
457 const Array<OneD, NekDouble> coeffs);
458 /// Gathers the global coefficients \f$\boldsymbol{\hat{u}}_g\f$
459 /// from the local coefficients \f$\boldsymbol{\hat{u}}_l\f$.
460 // inline
461 MULTI_REGIONS_EXPORT inline void LocalToGlobal(bool useComm = true);
463 const Array<OneD, const NekDouble> &inarray,
464 Array<OneD, NekDouble> &outarray, bool useComm = true);
465 /// Scatters from the global coefficients
466 /// \f$\boldsymbol{\hat{u}}_g\f$ to the local coefficients
467 /// \f$\boldsymbol{\hat{u}}_l\f$.
468 // inline
469 MULTI_REGIONS_EXPORT inline void GlobalToLocal(void);
470 /**
471 * This operation is evaluated as:
472 * \f{tabbing}
473 * \hspace{1cm} \= Do \= $e=$ $1, N_{\mathrm{el}}$ \\
474 * > > Do \= $i=$ $0,N_m^e-1$ \\
475 * > > > $\boldsymbol{\hat{u}}^{e}[i] = \mbox{sign}[e][i] \cdot
476 * \boldsymbol{\hat{u}}_g[\mbox{map}[e][i]]$ \\
477 * > > continue \\
478 * > continue
479 * \f}
480 * where \a map\f$[e][i]\f$ is the mapping array and \a
481 * sign\f$[e][i]\f$ is an array of similar dimensions ensuring the
482 * correct modal connectivity between the different elements (both
483 * these arrays are contained in the data member #m_locToGloMap). This
484 * operation is equivalent to the scatter operation
485 * \f$\boldsymbol{\hat{u}}_l=\mathcal{A}\boldsymbol{\hat{u}}_g\f$,
486 * where \f$\mathcal{A}\f$ is the
487 * \f$N_{\mathrm{eof}}\times N_{\mathrm{dof}}\f$ permutation matrix.
488 *
489 * @param inarray An array of size \f$N_\mathrm{dof}\f$
490 * containing the global degrees of freedom
491 * \f$\boldsymbol{x}_g\f$.
492 * @param outarray The resulting local degrees of freedom
493 * \f$\boldsymbol{x}_l\f$ will be stored in this
494 * array of size \f$N_\mathrm{eof}\f$.
495 */
497 const Array<OneD, const NekDouble> &inarray,
498 Array<OneD, NekDouble> &outarray);
499 /// Get the \a i th value (coefficient) of #m_coeffs
500 inline NekDouble GetCoeff(int i);
501 /// Get the \a i th value (coefficient) of #m_coeffs
502 inline NekDouble GetCoeffs(int i);
503 /// This function returns (a reference to) the array
504 /// \f$\boldsymbol{u}_l\f$ (implemented as #m_phys) containing the
505 /// function \f$u^{\delta}(\boldsymbol{x})\f$ evaluated at the
506 /// quadrature points.
507 // inline
509 /// This function calculates the \f$L_\infty\f$ error of the global
510 /// spectral/hp element approximation.
512 Linf(const Array<OneD, const NekDouble> &inarray,
514 /// This function calculates the \f$L_\infty\f$ error of the global
515 /// This function calculates the \f$L_2\f$ error with
516 /// respect to soln of the global
517 /// spectral/hp element approximation.
519 const Array<OneD, const NekDouble> &inarray,
521 {
522 return v_L2(inarray, soln);
523 }
524 /// Calculates the \f$H^1\f$ error of the global spectral/hp
525 /// element approximation.
527 H1(const Array<OneD, const NekDouble> &inarray,
529 /// Calculates the \f$H^1\f$ error of the global spectral/hp
530 /// element approximation.
531 /**
532 * The integration is evaluated locally, that is
533 * \f[\int
534 * f(\boldsymbol{x})d\boldsymbol{x}=\sum_{e=1}^{{N_{\mathrm{el}}}}
535 * \left\{\int_{\Omega_e}f(\boldsymbol{x})d\boldsymbol{x}\right\}, \f]
536 * where the integration over the separate elements is done by the
537 * function StdRegions#StdExpansion#Integral, which discretely
538 * evaluates the integral using Gaussian quadrature.
539 *
540 * Note that the array #m_phys should be filled with the values of the
541 * function \f$f(\boldsymbol{x})\f$ at the quadrature points
542 * \f$\boldsymbol{x}_i\f$.
543 *
544 * @return The value of the discretely evaluated integral
545 * \f$\int f(\boldsymbol{x})d\boldsymbol{x}\f$.
546 */
548 {
549 ASSERTL1(m_physState == true, "local physical space is not true ");
550 return Integral(m_phys);
551 }
552 /**
553 * The integration is evaluated locally, that is
554 * \f[\int
555 * f(\boldsymbol{x})d\boldsymbol{x}=\sum_{e=1}^{{N_{\mathrm{el}}}}
556 * \left\{\int_{\Omega_e}f(\boldsymbol{x})d\boldsymbol{x}\right\}, \f]
557 * where the integration over the separate elements is done by the
558 * function StdRegions#StdExpansion#Integral, which discretely
559 * evaluates the integral using Gaussian quadrature.
560 *
561 * @param inarray An array of size \f$Q_{\mathrm{tot}}\f$
562 * containing the values of the function
563 * \f$f(\boldsymbol{x})\f$ at the quadrature
564 * points \f$\boldsymbol{x}_i\f$.
565 * @return The value of the discretely evaluated integral
566 * \f$\int f(\boldsymbol{x})d\boldsymbol{x}\f$.
567 */
569 {
570 return v_Integral(inarray);
571 }
573 {
574 return v_VectorFlux(inarray);
575 }
576 /// This function calculates the energy associated with
577 /// each one of the modesof a 3D homogeneous nD expansion
579 {
580 return v_HomogeneousEnergy();
581 }
582
583 /// This function sets the Spectral Vanishing Viscosity
584 /// in homogeneous1D expansion.
586 {
588 }
589
590 /// This function returns a vector containing the wave
591 /// numbers in z-direction associated
592 /// with the 3D homogenous expansion. Required if a
593 /// parellelisation is applied in the Fourier direction
595 {
596 return v_GetZIDs();
597 }
598
599 /// This function returns the transposition class
600 /// associated with the homogeneous expansion.
602 {
603 return v_GetTransposition();
604 }
605
606 /// This function returns the Width of homogeneous direction
607 /// associated with the homogeneous expansion.
609 {
610 return v_GetHomoLen();
611 }
612
613 /// This function sets the Width of homogeneous direction
614 /// associated with the homogeneous expansion.
615 void SetHomoLen(const NekDouble lhom)
616 {
617 return v_SetHomoLen(lhom);
618 }
619
620 /// This function returns a vector containing the wave
621 /// numbers in y-direction associated
622 /// with the 3D homogenous expansion. Required if a
623 /// parellelisation is applied in the Fourier direction
625 {
626 return v_GetYIDs();
627 }
628
629 /// This function interpolates the physical space points in
630 /// \a inarray to \a outarray using the same points defined in the
631 /// expansion but where the number of points are rescaled
632 /// by \a 1DScale
634 const Array<OneD, NekDouble> &inarray,
635 Array<OneD, NekDouble> &outarray)
636 {
637 v_PhysInterp1DScaled(scale, inarray, outarray);
638 }
639
640 /// This function Galerkin projects the physical space points in
641 /// \a inarray to \a outarray where inarray is assumed to
642 /// be defined in the expansion but where the number of
643 /// points are rescaled by \a 1DScale
645 const Array<OneD, NekDouble> &inarray,
646 Array<OneD, NekDouble> &outarray)
647 {
648 v_PhysGalerkinProjection1DScaled(scale, inarray, outarray);
649 }
650
651 /// This function returns the number of elements in the expansion.
652 inline int GetExpSize(void);
653
654 /// This function returns the number of elements in the
655 /// expansion which may be different for a homogeoenous extended
656 /// expansionp.
657 inline size_t GetNumElmts(void)
658 {
659 return v_GetNumElmts();
660 }
661
662 /// This function returns the vector of elements in the expansion.
663 inline const std::shared_ptr<LocalRegions::ExpansionVector> GetExp() const;
664
665 /// This function returns (a shared pointer to) the local elemental
666 /// expansion of the \f$n^{\mathrm{th}}\f$ element.
667 inline LocalRegions::ExpansionSharedPtr &GetExp(int n) const;
668
669 /// This function returns (a shared pointer to) the local elemental
670 /// expansion of the \f$n^{\mathrm{th}}\f$ element given a global
671 /// geometry ID.
673
674 /// This function returns (a shared pointer to) the local elemental
675 /// expansion containing the arbitrary point given by \a gloCoord.
677 const Array<OneD, const NekDouble> &gloCoord);
678
679 /**
680 * @brief This function returns the index of the local elemental
681 * expansion containing the arbitrary point given by \a gloCoord,
682 * within a distance tolerance of tol.
683 *
684 * If returnNearestElmt is true and no element contains the point,
685 * this function returns the nearest element whose bounding box contains
686 * the point. The bounding box has a 10% margin in each direction.
687 *
688 * @param gloCoord (input) coordinate in physics space
689 * @param locCoords (output) local coordinate xi in the returned element
690 * @param tol distance tolerance to judge if a point lies in an
691 * element
692 * @param returnNearestElmt if true and no element contains this point, the
693 * nearest element whose bounding box contains this
694 * point is returned
695 * @param cachedId an initial guess of the most possible element index
696 * @param maxDistance if returnNearestElmt is set as true, the nearest
697 * element will be returned. But the distance of the
698 * nearest element and this point should be <=
699 * maxDistance.
700 *
701 * @return element index; if no element is found, -1 is returned.
702 */
704 const Array<OneD, const NekDouble> &gloCoord, NekDouble tol = 0.0,
705 bool returnNearestElmt = false, int cachedId = -1,
706 NekDouble maxDistance = 1e6);
707
708 /** This function returns the index and the Local
709 * Cartesian Coordinates \a locCoords of the local
710 * elemental expansion containing the arbitrary point
711 * given by \a gloCoords.
712 **/
714 const Array<OneD, const NekDouble> &gloCoords,
715 Array<OneD, NekDouble> &locCoords, NekDouble tol = 0.0,
716 bool returnNearestElmt = false, int cachedId = -1,
717 NekDouble maxDistance = 1e6);
718
719 /** This function return the expansion field value
720 * at the coordinates given as input.
721 **/
724 const Array<OneD, const NekDouble> &phys);
725
726 /// Get the start offset position for a local contiguous list of
727 /// coeffs correspoinding to element n.
728 inline int GetCoeff_Offset(int n) const;
729
730 /// Get the start offset position for a local contiguous list of
731 /// quadrature points in a full array correspoinding to element n.
732 inline int GetPhys_Offset(int n) const;
733
734 /// This function returns (a reference to) the array
735 /// \f$\boldsymbol{\hat{u}}_l\f$ (implemented as #m_coeffs)
736 /// containing all local expansion coefficients.
738 /// This function returns (a reference to) the array
739 /// \f$\boldsymbol{u}_l\f$ (implemented as #m_phys) containing the
740 /// function \f$u^{\delta}(\boldsymbol{x})\f$ evaluated at the
741 /// quadrature points.
743 inline void PhysDeriv(Direction edir,
744 const Array<OneD, const NekDouble> &inarray,
746 /// This function discretely evaluates the derivative of a function
747 /// \f$f(\boldsymbol{x})\f$ on the domain consisting of all
748 /// elements of the expansion.
749 inline void PhysDeriv(
750 const Array<OneD, const NekDouble> &inarray,
754 inline void PhysDeriv(const int dir,
755 const Array<OneD, const NekDouble> &inarray,
757
758 inline void Curl(Array<OneD, Array<OneD, NekDouble>> &Vel,
760
761 inline void CurlCurl(Array<OneD, Array<OneD, NekDouble>> &Vel,
763 inline void PhysDirectionalDeriv(
764 const Array<OneD, const NekDouble> &direction,
765 const Array<OneD, const NekDouble> &inarray,
766 Array<OneD, NekDouble> &outarray);
767 inline void GetMovingFrames(const SpatialDomains::GeomMMF MMFdir,
768 const Array<OneD, const NekDouble> &CircCentre,
769 Array<OneD, Array<OneD, NekDouble>> &outarray);
770 // functions associated with DisContField
773 /// Get the weight value for boundary conditions
775 /// Set the weight value for boundary conditions
776 inline void SetBndCondBwdWeight(const int index, const NekDouble value);
777 inline std::shared_ptr<ExpList> &UpdateBndCondExpansion(int i);
778 inline void Upwind(const Array<OneD, const NekDouble> &Vn,
782 inline void Upwind(const Array<OneD, const Array<OneD, NekDouble>> &Vec,
786 /**
787 * Return a reference to the trace space associated with this
788 * expansion list.
789 */
790 inline std::shared_ptr<ExpList> &GetTrace();
791 inline std::shared_ptr<AssemblyMapDG> &GetTraceMap(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);
803 Array<OneD, NekDouble> &outarray);
805 Array<OneD, NekDouble> &outarray);
808 Array<OneD, NekDouble> &outarray);
811 inline void GetFwdBwdTracePhys(const Array<OneD, const NekDouble> &field,
814 bool FillBnd = true,
815 bool PutFwdInBwdOnBCs = false,
816 bool DoExchange = true);
817 inline void FillBwdWithBoundCond(const Array<OneD, NekDouble> &Fwd,
819 bool PutFwdInBwdOnBCs = false);
820 /// Add Fwd and Bwd value to field,
821 /// a reverse procedure of GetFwdBwdTracePhys
828 {
829 v_AddTraceQuadPhysToOffDiag(Fwd, Bwd, field);
830 }
833 Array<OneD, NekDouble> &locTraceFwd,
834 Array<OneD, NekDouble> &locTraceBwd);
835 /// Fill Bwd with boundary conditions
836 inline void FillBwdWithBwdWeight(Array<OneD, NekDouble> &weightave,
837 Array<OneD, NekDouble> &weightjmp);
838 /// Copy and fill the Periodic boundaries
839 inline void PeriodicBwdCopy(const Array<OneD, const NekDouble> &Fwd,
841 inline const std::vector<bool> &GetLeftAdjacentFaces(void) const;
842 inline void ExtractTracePhys(Array<OneD, NekDouble> &outarray);
843 inline void ExtractTracePhys(const Array<OneD, const NekDouble> &inarray,
844 Array<OneD, NekDouble> &outarray);
849 inline void EvaluateBoundaryConditions(
850 const NekDouble time = 0.0, const std::string varName = "",
853 // Routines for continous matrix solution
854 /// This function calculates the result of the multiplication of a
855 /// matrix of type specified by \a mkey with a vector given by \a
856 /// inarray.
858 const GlobalMatrixKey &gkey,
859 const Array<OneD, const NekDouble> &inarray,
860 Array<OneD, NekDouble> &outarray);
861 inline void SetUpPhysNormals();
862 inline void GetBoundaryToElmtMap(Array<OneD, int> &ElmtID,
863 Array<OneD, int> &EdgeID);
864 virtual void GetBndElmtExpansion(int i, std::shared_ptr<ExpList> &result,
865 const bool DeclareCoeffPhysArrays = true);
866
867 inline void ExtractElmtToBndPhys(int i, const Array<OneD, NekDouble> &elmt,
868 Array<OneD, NekDouble> &boundary);
869
870 inline void ExtractPhysToBndElmt(int i,
872 Array<OneD, NekDouble> &bndElmt);
873
874 inline void ExtractPhysToBnd(int i,
877
878 inline void GetBoundaryNormals(
879 int i, Array<OneD, Array<OneD, NekDouble>> &normals);
880
882 std::vector<LibUtilities::FieldDefinitionsSharedPtr> &fielddef,
883 int NumHomoDir = 0,
886 std::vector<NekDouble> &HomoLen = LibUtilities::NullNekDoubleVector,
887 bool homoStrips = false,
888 std::vector<unsigned int> &HomoSIDs =
890 std::vector<unsigned int> &HomoZIDs =
892 std::vector<unsigned int> &HomoYIDs =
894
895 std::map<int, RobinBCInfoSharedPtr> GetRobinBCInfo()
896 {
897 return v_GetRobinBCInfo();
898 }
899
900 void GetPeriodicEntities(PeriodicMap &periodicVerts,
901 PeriodicMap &periodicEdges,
902 PeriodicMap &periodicFaces = NullPeriodicMap)
903 {
904 v_GetPeriodicEntities(periodicVerts, periodicEdges, periodicFaces);
905 }
906
907 std::vector<LibUtilities::FieldDefinitionsSharedPtr> GetFieldDefinitions()
908 {
909 return v_GetFieldDefinitions();
910 }
911
913 std::vector<LibUtilities::FieldDefinitionsSharedPtr> &fielddef)
914 {
915 v_GetFieldDefinitions(fielddef);
916 }
917
918 /// Append the element data listed in elements
919 /// fielddef->m_ElementIDs onto fielddata
921 std::vector<NekDouble> &fielddata)
922 {
923 v_AppendFieldData(fielddef, fielddata);
924 }
925 /// Append the data in coeffs listed in elements
926 /// fielddef->m_ElementIDs onto fielddata
928 std::vector<NekDouble> &fielddata,
930 {
931 v_AppendFieldData(fielddef, fielddata, coeffs);
932 }
933 /** \brief Extract the data in fielddata into the coeffs
934 * using the basic ExpList Elemental expansions rather
935 * than planes in homogeneous case
936 */
939 std::vector<NekDouble> &fielddata, std::string &field,
940 Array<OneD, NekDouble> &coeffs);
941 /** \brief Extract the data from fromField using
942 * fromExpList the coeffs using the basic ExpList
943 * Elemental expansions rather than planes in homogeneous
944 * case
945 */
947 const std::shared_ptr<ExpList> &fromExpList,
948 const Array<OneD, const NekDouble> &fromCoeffs,
949 Array<OneD, NekDouble> &toCoeffs);
950 // Extract data in fielddata into the m_coeffs_list for the 3D stability
951 // analysis (base flow is 2D)
954 std::vector<NekDouble> &fielddata, std::string &field,
956 std::unordered_map<int, int> zIdToPlane =
957 std::unordered_map<int, int>());
958 // Extract data from file fileName and put coefficents into array coefffs
960 const std::string &fileName, LibUtilities::CommSharedPtr comm,
961 const std::string &varName, Array<OneD, NekDouble> &coeffs);
963 const int ElementID, const NekDouble scalar1, const NekDouble scalar2,
964 Array<OneD, NekDouble> &outarray);
965 /// Returns a shared pointer to the current object.
966 std::shared_ptr<ExpList> GetSharedThisPtr()
967 {
968 return shared_from_this();
969 }
970 /// Returns the session object
971 std::shared_ptr<LibUtilities::SessionReader> GetSession() const
972 {
973 return m_session;
974 }
975 /// Returns the comm object
976 std::shared_ptr<LibUtilities::Comm> GetComm() const
977 {
978 return m_comm;
979 }
981 {
982 return m_graph;
983 }
984 // Wrapper functions for Homogeneous Expansions
986 {
987 return v_GetHomogeneousBasis();
988 }
989 std::shared_ptr<ExpList> &GetPlane(int n)
990 {
991 return v_GetPlane(n);
992 }
996
997 MULTI_REGIONS_EXPORT int GetPoolCount(std::string);
998
1000
1003
1004 /// Get m_coeffs to elemental value map
1006 &GetCoeffsToElmt() const;
1012 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
1013 const int nDirctn, Array<OneD, DNekMatSharedPtr> &mtxPerVar);
1015 const TensorOfArray3D<NekDouble> &inarray,
1018 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
1021 const Array<OneD, const NekDouble> &FwdFlux,
1022 const Array<OneD, const NekDouble> &BwdFlux,
1023 Array<OneD, NekDouble> &outarray)
1024 {
1025 v_AddTraceIntegralToOffDiag(FwdFlux, BwdFlux, outarray);
1026 }
1028 const int dir, const Array<OneD, const DNekMatSharedPtr> ElmtJacQuad,
1029 Array<OneD, DNekMatSharedPtr> ElmtJacCoef);
1031 const Array<OneD, const DNekMatSharedPtr> ElmtJacQuad,
1032 Array<OneD, DNekMatSharedPtr> ElmtJacCoef);
1034 GetLocTraceToTraceMap() const;
1035 // Return the internal vector which identifieds if trace
1036 // is left adjacent definiing which trace the normal
1037 // points otwards from
1038 MULTI_REGIONS_EXPORT std::vector<bool> &GetLeftAdjacentTraces(void);
1039
1040 /// This function returns the map of index inside m_exp to geom id
1041 MULTI_REGIONS_EXPORT inline const std::unordered_map<int, int>
1043 {
1044 return m_elmtToExpId;
1045 }
1046
1047 /// This function returns the index inside m_exp for a given geom id
1049 {
1050 auto it = m_elmtToExpId.find(elmtId);
1051 ASSERTL0(it != m_elmtToExpId.end(), "Global geometry ID " +
1052 std::to_string(elmtId) +
1053 " not found in element ID to "
1054 "expansion ID map.")
1055 return it->second;
1056 }
1057
1058protected:
1059 /// Exapnsion type
1061 std::shared_ptr<DNekMat> GenGlobalMatrixFull(
1062 const GlobalLinSysKey &mkey,
1063 const std::shared_ptr<AssemblyMapCG> &locToGloMap);
1064 /// Communicator
1066 /// Session
1068 /// Mesh associated with this expansion list.
1070 /// The total number of local degrees of freedom. #m_ncoeffs
1071 /// \f$=N_{\mathrm{eof}}=\sum_{e=1}^{{N_{\mathrm{el}}}}N^{e}_l\f$
1073 /** The total number of quadrature points. #m_npoints
1074 *\f$=Q_{\mathrm{tot}}=\sum_{e=1}^{{N_{\mathrm{el}}}}N^{e}_Q\f$
1075 **/
1077 /**
1078 * \brief Concatenation of all local expansion coefficients.
1079 *
1080 * The array of length #m_ncoeffs\f$=N_{\mathrm{eof}}\f$ which is
1081 * the concatenation of the local expansion coefficients
1082 * \f$\hat{u}_n^e\f$ over all \f$N_{\mathrm{el}}\f$ elements
1083 * \f[\mathrm{\texttt{m\_coeffs}}=\boldsymbol{\hat{u}}_{l} =
1084 * \underline{\boldsymbol{\hat{u}}}^e = \left [ \begin{array}{c}
1085 * \boldsymbol{\hat{u}}^{1} \ \
1086 * \boldsymbol{\hat{u}}^{2} \ \
1087 * \vdots \ \
1088 * \boldsymbol{\hat{u}}^{{{N_{\mathrm{el}}}}} \end{array} \right ],
1089 * \quad
1090 * \mathrm{where}\quad \boldsymbol{\hat{u}}^{e}[n]=\hat{u}_n^{e}\f]
1091 */
1093 /**
1094 * \brief The global expansion evaluated at the quadrature points
1095 *
1096 * The array of length #m_npoints\f$=Q_{\mathrm{tot}}\f$ containing
1097 * the evaluation of \f$u^{\delta}(\boldsymbol{x})\f$ at the
1098 * quadrature points \f$\boldsymbol{x}_i\f$.
1099 * \f[\mathrm{\texttt{m\_phys}}=\boldsymbol{u}_{l} =
1100 * \underline{\boldsymbol{u}}^e = \left [ \begin{array}{c}
1101 * \boldsymbol{u}^{1} \ \
1102 * \boldsymbol{u}^{2} \ \
1103 * \vdots \ \
1104 * \boldsymbol{u}^{{{N_{\mathrm{el}}}}} \end{array} \right ],\quad
1105 * \mathrm{where}\quad
1106 * \boldsymbol{u}^{e}[i]=u^{\delta}(\boldsymbol{x}_i)\f]
1107 */
1109 /**
1110 * \brief The state of the array #m_phys.
1111 *
1112 * Indicates whether the array #m_phys, created to contain the
1113 * evaluation of \f$u^{\delta}(\boldsymbol{x})\f$ at the quadrature
1114 * points \f$\boldsymbol{x}_i\f$, is filled with these values.
1115 */
1117 /**
1118 * \brief The list of local expansions.
1119 *
1120 * The (shared pointer to the) vector containing (shared pointers
1121 * to) all local expansions. The fact that the local expansions are
1122 * all stored as a (pointer to a) #StdExpansion, the abstract base
1123 * class for all local expansions, allows a general implementation
1124 * where most of the routines for the derived classes are defined
1125 * in the #ExpList base class.
1126 */
1127 std::shared_ptr<LocalRegions::ExpansionVector> m_exp;
1129 /// Vector of bools to act as an initialise on first call flag
1130 std::vector<bool> m_collectionsDoInit;
1131 /// Offset of elemental data into the array #m_coeffs
1132 std::vector<int> m_coll_coeff_offset;
1133 /// Offset of elemental data into the array #m_phys
1134 std::vector<int> m_coll_phys_offset;
1135 /// Offset of elemental data into the array #m_coeffs
1137 /// Offset of elemental data into the array #m_phys
1139 /// m_coeffs to elemental value map
1142 //@todo should this be in ExpList or
1143 // ExpListHomogeneous1D.cpp it's a bool which determine if
1144 // the expansion is in the wave space (coefficient space)
1145 // or not
1147 /// Mapping from geometry ID of element to index inside #m_exp
1148 std::unordered_map<int, int> m_elmtToExpId;
1149 /// This function assembles the block diagonal matrix of local
1150 /// matrices of the type \a mtype.
1153 void MultiplyByBlockMatrix(const GlobalMatrixKey &gkey,
1154 const Array<OneD, const NekDouble> &inarray,
1155 Array<OneD, NekDouble> &outarray);
1156 /// Generates a global matrix from the given key and map.
1157 std::shared_ptr<GlobalMatrix> GenGlobalMatrix(
1158 const GlobalMatrixKey &mkey,
1159 const std::shared_ptr<AssemblyMapCG> &locToGloMap);
1161 const std::shared_ptr<DNekMat> &Gmat,
1162 Array<OneD, NekDouble> &EigValsReal,
1163 Array<OneD, NekDouble> &EigValsImag,
1165 /// This operation constructs the global linear system of type \a
1166 /// mkey.
1167 std::shared_ptr<GlobalLinSys> GenGlobalLinSys(
1168 const GlobalLinSysKey &mkey,
1169 const std::shared_ptr<AssemblyMapCG> &locToGloMap);
1170 /// Generate a GlobalLinSys from information provided by the key
1171 /// "mkey" and the mapping provided in LocToGloBaseMap.
1172 std::shared_ptr<GlobalLinSys> GenGlobalBndLinSys(
1173 const GlobalLinSysKey &mkey, const AssemblyMapSharedPtr &locToGloMap);
1174 // Virtual prototypes
1175 virtual size_t v_GetNumElmts(void)
1176 {
1177 return (*m_exp).size();
1178 }
1182 virtual void v_SetBndCondBwdWeight(const int index, const NekDouble value);
1183 virtual std::shared_ptr<ExpList> &v_UpdateBndCondExpansion(int i);
1184 virtual void v_Upwind(const Array<OneD, const Array<OneD, NekDouble>> &Vec,
1188 virtual void v_Upwind(const Array<OneD, const NekDouble> &Vn,
1192 virtual std::shared_ptr<ExpList> &v_GetTrace();
1193 virtual std::shared_ptr<AssemblyMapDG> &v_GetTraceMap();
1195 virtual const std::shared_ptr<LocTraceToTraceMap> &v_GetLocTraceToTraceMap(
1196 void) const;
1197 virtual std::vector<bool> &v_GetLeftAdjacentTraces(void);
1198 /// Populate \a normals with the normals of all expansions.
1199 virtual void v_GetNormals(Array<OneD, Array<OneD, NekDouble>> &normals);
1200 virtual void v_AddTraceIntegral(const Array<OneD, const NekDouble> &Fx,
1202 Array<OneD, NekDouble> &outarray);
1203 virtual void v_AddTraceIntegral(const Array<OneD, const NekDouble> &Fn,
1204 Array<OneD, NekDouble> &outarray);
1205 virtual void v_AddFwdBwdTraceIntegral(
1208 Array<OneD, NekDouble> &outarray);
1211 virtual void v_GetFwdBwdTracePhys(const Array<OneD, const NekDouble> &field,
1214 bool FillBnd = true,
1215 bool PutFwdInBwdOnBCs = false,
1216 bool DoExchange = true);
1217 virtual void v_FillBwdWithBoundCond(const Array<OneD, NekDouble> &Fwd,
1219 bool PutFwdInBwdOnBCs);
1220 virtual void v_AddTraceQuadPhysToField(
1223 virtual void v_AddTraceQuadPhysToOffDiag(
1226 virtual void v_GetLocTraceFromTracePts(
1229 Array<OneD, NekDouble> &locTraceFwd,
1230 Array<OneD, NekDouble> &locTraceBwd);
1231 virtual void v_FillBwdWithBwdWeight(Array<OneD, NekDouble> &weightave,
1232 Array<OneD, NekDouble> &weightjmp);
1233 virtual void v_PeriodicBwdCopy(const Array<OneD, const NekDouble> &Fwd,
1235 virtual const std::vector<bool> &v_GetLeftAdjacentFaces(void) const;
1236 virtual void v_ExtractTracePhys(Array<OneD, NekDouble> &outarray);
1237 virtual void v_ExtractTracePhys(const Array<OneD, const NekDouble> &inarray,
1238 Array<OneD, NekDouble> &outarray);
1239 virtual void v_MultiplyByInvMassMatrix(
1240 const Array<OneD, const NekDouble> &inarray,
1241 Array<OneD, NekDouble> &outarray);
1242
1244 const Array<OneD, const NekDouble> &inarray,
1245 Array<OneD, NekDouble> &outarray,
1247 const StdRegions::VarCoeffMap &varcoeff,
1248 const MultiRegions::VarFactorsMap &varfactors,
1249 const Array<OneD, const NekDouble> &dirForcing,
1250 const bool PhysSpaceForcing);
1251
1253 const Array<OneD, const NekDouble> &inarray,
1254 Array<OneD, NekDouble> &outarray,
1256 const StdRegions::VarCoeffMap &varcoeff,
1257 const MultiRegions::VarFactorsMap &varfactors,
1258 const Array<OneD, const NekDouble> &dirForcing,
1259 const bool PhysSpaceForcing);
1260
1261 virtual void v_LinearAdvectionReactionSolve(
1262 const Array<OneD, Array<OneD, NekDouble>> &velocity,
1263 const Array<OneD, const NekDouble> &inarray,
1264 Array<OneD, NekDouble> &outarray, const NekDouble lambda,
1266 // wrapper functions about virtual functions
1268 virtual void v_FillBndCondFromField(const Array<OneD, NekDouble> coeffs);
1269 virtual void v_FillBndCondFromField(const int nreg,
1270 const Array<OneD, NekDouble> coeffs);
1271 virtual void v_Reset();
1272 virtual void v_LocalToGlobal(bool UseComm);
1273 virtual void v_LocalToGlobal(const Array<OneD, const NekDouble> &inarray,
1274 Array<OneD, NekDouble> &outarray,
1275 bool UseComm);
1276 virtual void v_GlobalToLocal(void);
1277 virtual void v_GlobalToLocal(const Array<OneD, const NekDouble> &inarray,
1278 Array<OneD, NekDouble> &outarray);
1279 virtual void v_BwdTrans(const Array<OneD, const NekDouble> &inarray,
1280 Array<OneD, NekDouble> &outarray);
1281 virtual void v_FwdTrans(const Array<OneD, const NekDouble> &inarray,
1282 Array<OneD, NekDouble> &outarray);
1283 virtual void v_FwdTransLocalElmt(
1284 const Array<OneD, const NekDouble> &inarray,
1285 Array<OneD, NekDouble> &outarray);
1286 virtual void v_FwdTransBndConstrained(
1287 const Array<OneD, const NekDouble> &inarray,
1288 Array<OneD, NekDouble> &outarray);
1289 virtual void v_SmoothField(Array<OneD, NekDouble> &field);
1290 virtual void v_IProductWRTBase(const Array<OneD, const NekDouble> &inarray,
1291 Array<OneD, NekDouble> &outarray);
1292
1293 // Define ExpList::IProductWRTDerivBase as virtual function
1294 virtual void v_IProductWRTDerivBase(
1295 const int dir, const Array<OneD, const NekDouble> &inarray,
1296 Array<OneD, NekDouble> &outarray);
1297
1298 virtual void v_IProductWRTDerivBase(
1299 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
1300 Array<OneD, NekDouble> &outarray);
1301
1302 virtual void v_GetCoords(
1303 Array<OneD, NekDouble> &coord_0,
1306
1307 virtual void v_GetCoords(const int eid, Array<OneD, NekDouble> &xc0,
1310
1311 virtual void v_PhysDeriv(const Array<OneD, const NekDouble> &inarray,
1312 Array<OneD, NekDouble> &out_d0,
1313 Array<OneD, NekDouble> &out_d1,
1314 Array<OneD, NekDouble> &out_d2);
1315 virtual void v_PhysDeriv(const int dir,
1316 const Array<OneD, const NekDouble> &inarray,
1317 Array<OneD, NekDouble> &out_d);
1318 virtual void v_PhysDeriv(Direction edir,
1319 const Array<OneD, const NekDouble> &inarray,
1320 Array<OneD, NekDouble> &out_d);
1321
1322 virtual void v_Curl(Array<OneD, Array<OneD, NekDouble>> &Vel,
1324
1325 virtual void v_CurlCurl(Array<OneD, Array<OneD, NekDouble>> &Vel,
1327
1328 virtual void v_PhysDirectionalDeriv(
1329 const Array<OneD, const NekDouble> &direction,
1330 const Array<OneD, const NekDouble> &inarray,
1331 Array<OneD, NekDouble> &outarray);
1332 virtual void v_GetMovingFrames(
1333 const SpatialDomains::GeomMMF MMFdir,
1334 const Array<OneD, const NekDouble> &CircCentre,
1335 Array<OneD, Array<OneD, NekDouble>> &outarray);
1336 virtual void v_HomogeneousFwdTrans(
1337 const int npts, const Array<OneD, const NekDouble> &inarray,
1338 Array<OneD, NekDouble> &outarray, bool Shuff = true,
1339 bool UnShuff = true);
1340 virtual void v_HomogeneousBwdTrans(
1341 const int npts, const Array<OneD, const NekDouble> &inarray,
1342 Array<OneD, NekDouble> &outarray, bool Shuff = true,
1343 bool UnShuff = true);
1344 virtual void v_DealiasedProd(const int num_dofs,
1345 const Array<OneD, NekDouble> &inarray1,
1346 const Array<OneD, NekDouble> &inarray2,
1347 Array<OneD, NekDouble> &outarray);
1348 virtual void v_DealiasedDotProd(
1349 const int num_dofs, const Array<OneD, Array<OneD, NekDouble>> &inarray1,
1350 const Array<OneD, Array<OneD, NekDouble>> &inarray2,
1351 Array<OneD, Array<OneD, NekDouble>> &outarray);
1352 virtual void v_GetBCValues(Array<OneD, NekDouble> &BndVals,
1353 const Array<OneD, NekDouble> &TotField,
1354 int BndID);
1357 Array<OneD, NekDouble> &outarray,
1358 int BndID);
1359 virtual void v_NormVectorIProductWRTBase(
1361 Array<OneD, NekDouble> &outarray);
1362 virtual void v_SetUpPhysNormals();
1363 virtual void v_GetBoundaryToElmtMap(Array<OneD, int> &ElmtID,
1364 Array<OneD, int> &EdgeID);
1365 virtual void v_GetBndElmtExpansion(int i, std::shared_ptr<ExpList> &result,
1366 const bool DeclareCoeffPhysArrays);
1367
1368 virtual void v_ExtractElmtToBndPhys(const int i,
1369 const Array<OneD, NekDouble> &elmt,
1370 Array<OneD, NekDouble> &boundary);
1371
1372 virtual void v_ExtractPhysToBndElmt(
1373 const int i, const Array<OneD, const NekDouble> &phys,
1374 Array<OneD, NekDouble> &bndElmt);
1375
1376 virtual void v_ExtractPhysToBnd(const int i,
1377 const Array<OneD, const NekDouble> &phys,
1379
1380 virtual void v_GetBoundaryNormals(
1381 int i, Array<OneD, Array<OneD, NekDouble>> &normals);
1382
1383 virtual std::vector<LibUtilities::FieldDefinitionsSharedPtr>
1385
1386 virtual void v_GetFieldDefinitions(
1387 std::vector<LibUtilities::FieldDefinitionsSharedPtr> &fielddef);
1388
1389 virtual void v_AppendFieldData(
1391 std::vector<NekDouble> &fielddata);
1392 virtual void v_AppendFieldData(
1394 std::vector<NekDouble> &fielddata, Array<OneD, NekDouble> &coeffs);
1395 virtual void v_ExtractDataToCoeffs(
1397 std::vector<NekDouble> &fielddata, std::string &field,
1398 Array<OneD, NekDouble> &coeffs,
1399 std::unordered_map<int, int> zIdToPlane);
1400 virtual void v_ExtractCoeffsToCoeffs(
1401 const std::shared_ptr<ExpList> &fromExpList,
1402 const Array<OneD, const NekDouble> &fromCoeffs,
1403 Array<OneD, NekDouble> &toCoeffs);
1404 virtual void v_WriteTecplotHeader(std::ostream &outfile,
1405 std::string var = "");
1406 virtual void v_WriteTecplotZone(std::ostream &outfile, int expansion);
1407 virtual void v_WriteTecplotField(std::ostream &outfile, int expansion);
1408 virtual void v_WriteTecplotConnectivity(std::ostream &outfile,
1409 int expansion);
1410 virtual void v_WriteVtkPieceData(std::ostream &outfile, int expansion,
1411 std::string var);
1412 virtual void v_WriteVtkPieceHeader(std::ostream &outfile, int expansion,
1413 int istrip);
1414
1415 virtual NekDouble v_L2(
1416 const Array<OneD, const NekDouble> &phys,
1418
1419 virtual NekDouble v_Integral(const Array<OneD, const NekDouble> &inarray);
1420 virtual NekDouble v_VectorFlux(
1421 const Array<OneD, Array<OneD, NekDouble>> &inarray);
1422
1424
1426 virtual NekDouble v_GetHomoLen(void);
1427 virtual void v_SetHomoLen(const NekDouble lhom);
1430
1431 // 1D Scaling and projection
1432 virtual void v_PhysInterp1DScaled(const NekDouble scale,
1433 const Array<OneD, NekDouble> &inarray,
1434 Array<OneD, NekDouble> &outarray);
1435
1437 const NekDouble scale, const Array<OneD, NekDouble> &inarray,
1438 Array<OneD, NekDouble> &outarray);
1439
1440 virtual void v_ClearGlobalLinSysManager(void);
1441
1442 virtual int v_GetPoolCount(std::string);
1443
1444 virtual void v_UnsetGlobalLinSys(GlobalLinSysKey, bool);
1445
1448
1449 void ExtractFileBCs(const std::string &fileName,
1451 const std::string &varName,
1452 const std::shared_ptr<ExpList> locExpList);
1453
1454 // Utility function for a common case of retrieving a
1455 // BoundaryCondition from a boundary condition collection.
1459 unsigned int index, const std::string &variable);
1460
1463
1466
1467 virtual void v_EvaluateBoundaryConditions(
1468 const NekDouble time = 0.0, const std::string varName = "",
1471
1472 virtual std::map<int, RobinBCInfoSharedPtr> v_GetRobinBCInfo(void);
1473
1474 virtual void v_GetPeriodicEntities(PeriodicMap &periodicVerts,
1475 PeriodicMap &periodicEdges,
1476 PeriodicMap &periodicFaces);
1477
1478 // Homogeneous direction wrapper functions.
1480 {
1482 "This method is not defined or valid for this class type");
1484 }
1485
1486 // wrapper function to set viscosity for Homo1D expansion
1488 {
1489 boost::ignore_unused(visc);
1491 "This method is not defined or valid for this class type");
1492 }
1493
1494 virtual std::shared_ptr<ExpList> &v_GetPlane(int n);
1495
1496 virtual void v_AddTraceIntegralToOffDiag(
1497 const Array<OneD, const NekDouble> &FwdFlux,
1498 const Array<OneD, const NekDouble> &BwdFlux,
1499 Array<OneD, NekDouble> &outarray);
1500
1501private:
1502 /// Definition of the total number of degrees of freedom and
1503 /// quadrature points and offsets to access data
1504 void SetupCoeffPhys(bool DeclareCoeffPhysArrays = true,
1505 bool SetupOffsets = true);
1506
1507 /// Define a list of elements using the geometry and basis
1508 /// key information in expmap;
1510};
1511
1512/// An empty ExpList object.
1515
1516// An empty GlobaLinSysManager and GlobalLinSysKey object
1520
1521// Inline routines follow.
1522
1523/**
1524 * Returns the total number of local degrees of freedom
1525 * \f$N_{\mathrm{eof}}=\sum_{e=1}^{{N_{\mathrm{el}}}}N^{e}_m\f$.
1526 */
1527inline int ExpList::GetNcoeffs() const
1528{
1529 return m_ncoeffs;
1530}
1531
1532inline int ExpList::GetNcoeffs(const int eid) const
1533{
1534 return (*m_exp)[eid]->GetNcoeffs();
1535}
1536
1537/**
1538 * Evaulates the maximum number of modes in the elemental basis
1539 * order over all elements
1540 */
1542{
1543 int returnval = 0;
1544
1545 for (size_t i = 0; i < (*m_exp).size(); ++i)
1546 {
1547 returnval = (std::max)(returnval, (*m_exp)[i]->EvalBasisNumModesMax());
1548 }
1549
1550 return returnval;
1551}
1552
1553/**
1554 *
1555 */
1557{
1558 Array<OneD, int> returnval((*m_exp).size(), 0);
1559
1560 for (size_t i = 0; i < (*m_exp).size(); ++i)
1561 {
1562 returnval[i] =
1563 (std::max)(returnval[i], (*m_exp)[i]->EvalBasisNumModesMax());
1564 }
1565
1566 return returnval;
1567}
1568
1569/**
1570 *
1571 */
1572inline int ExpList::GetTotPoints() const
1573{
1574 return m_npoints;
1575}
1576
1577inline int ExpList::GetTotPoints(const int eid) const
1578{
1579 return (*m_exp)[eid]->GetTotPoints();
1580}
1581
1582inline int ExpList::Get1DScaledTotPoints(const NekDouble scale) const
1583{
1584 int returnval = 0;
1585 size_t cnt;
1586 size_t nbase = (*m_exp)[0]->GetNumBases();
1587
1588 for (size_t i = 0; i < (*m_exp).size(); ++i)
1589 {
1590 cnt = 1;
1591 for (size_t j = 0; j < nbase; ++j)
1592 {
1593 cnt *= scale * ((*m_exp)[i]->GetNumPoints(j));
1594 }
1595 returnval += cnt;
1596 }
1597 return returnval;
1598}
1599
1600/**
1601 *
1602 */
1603inline int ExpList::GetNpoints() const
1604{
1605 return m_npoints;
1606}
1607
1608/**
1609 *
1610 */
1611inline void ExpList::SetWaveSpace(const bool wavespace)
1612{
1613 m_WaveSpace = wavespace;
1614}
1615
1616/**
1617 *
1618 */
1619inline bool ExpList::GetWaveSpace() const
1620{
1621 return m_WaveSpace;
1622}
1623
1624/// Set the \a i th value of\a m_phys to value \a val
1625inline void ExpList::SetPhys(int i, NekDouble val)
1626{
1627 m_phys[i] = val;
1628}
1629/**
1630 * This function fills the array \f$\boldsymbol{u}_l\f$, the evaluation
1631 * of the expansion at the quadrature points (implemented as #m_phys),
1632 * with the values of the array \a inarray.
1633 *
1634 * @param inarray The array containing the values where
1635 * #m_phys should be filled with.
1636 */
1638{
1639 ASSERTL0((int)inarray.size() == m_npoints,
1640 "Input array does not have correct number of elements.");
1641 Vmath::Vcopy(m_npoints, &inarray[0], 1, &m_phys[0], 1);
1642 m_physState = true;
1643}
1645{
1646 m_phys = inarray;
1647}
1648/**
1649 * @param physState \a true (=filled) or \a false (=not filled).
1650 */
1651inline void ExpList::SetPhysState(const bool physState)
1652{
1653 m_physState = physState;
1654}
1655/**
1656 * @return physState \a true (=filled) or \a false (=not filled).
1657 */
1658inline bool ExpList::GetPhysState() const
1659{
1660 return m_physState;
1661}
1662/**
1663 *
1664 */
1666 const Array<OneD, const NekDouble> &inarray,
1667 Array<OneD, NekDouble> &outarray)
1668{
1669 v_IProductWRTBase(inarray, outarray);
1670}
1671/**
1672 *
1673 */
1675 const int dir, const Array<OneD, const NekDouble> &inarray,
1676 Array<OneD, NekDouble> &outarray)
1677{
1678 v_IProductWRTDerivBase(dir, inarray, outarray);
1679}
1680
1681/**
1682 *
1683 */
1685 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
1686 Array<OneD, NekDouble> &outarray)
1687{
1688 v_IProductWRTDerivBase(inarray, outarray);
1689}
1690
1691/**
1692 *
1693 */
1695 Array<OneD, NekDouble> &outarray)
1696{
1697 v_FwdTrans(inarray, outarray);
1698}
1699/**
1700 *
1701 */
1703 const Array<OneD, const NekDouble> &inarray,
1704 Array<OneD, NekDouble> &outarray)
1705{
1706 v_FwdTransLocalElmt(inarray, outarray);
1707}
1708/**
1709 *
1710 */
1712 const Array<OneD, const NekDouble> &inarray,
1713 Array<OneD, NekDouble> &outarray)
1714{
1715 v_FwdTransBndConstrained(inarray, outarray);
1716}
1717/**
1718 *
1719 */
1721{
1722 v_SmoothField(field);
1723}
1724/**
1725 *
1726 */
1727/**
1728 * Given the coefficients of an expansion, this function evaluates the
1729 * spectral/hp expansion \f$u^{\delta}(\boldsymbol{x})\f$ at the
1730 * quadrature points \f$\boldsymbol{x}_i\f$.
1731 */
1733 Array<OneD, NekDouble> &outarray)
1734{
1735 v_BwdTrans(inarray, outarray);
1736}
1737/**
1738 *
1739 */
1741 const Array<OneD, const NekDouble> &inarray,
1742 Array<OneD, NekDouble> &outarray)
1743{
1744 v_MultiplyByInvMassMatrix(inarray, outarray);
1745}
1746/**
1747 *
1748 */
1750 const Array<OneD, const NekDouble> &inarray,
1752 const StdRegions::VarCoeffMap &varcoeff,
1753 const MultiRegions::VarFactorsMap &varfactors,
1754 const Array<OneD, const NekDouble> &dirForcing, const bool PhysSpaceForcing)
1755{
1756 return v_HelmSolve(inarray, outarray, factors, varcoeff, varfactors,
1757 dirForcing, PhysSpaceForcing);
1758}
1759/**
1760 *
1761 */
1763 const Array<OneD, const NekDouble> &inarray,
1765 const StdRegions::VarCoeffMap &varcoeff,
1766 const MultiRegions::VarFactorsMap &varfactors,
1767 const Array<OneD, const NekDouble> &dirForcing, const bool PhysSpaceForcing)
1768{
1770 inarray, outarray, factors, varcoeff, varfactors, dirForcing,
1771 PhysSpaceForcing);
1772}
1774 const Array<OneD, Array<OneD, NekDouble>> &velocity,
1775 const Array<OneD, const NekDouble> &inarray,
1776 Array<OneD, NekDouble> &outarray, const NekDouble lambda,
1777 const Array<OneD, const NekDouble> &dirForcing)
1778{
1779 v_LinearAdvectionReactionSolve(velocity, inarray, outarray, lambda,
1780 dirForcing);
1781}
1782/**
1783 *
1784 */
1786 Array<OneD, NekDouble> &coord_1,
1787 Array<OneD, NekDouble> &coord_2)
1788{
1789 v_GetCoords(coord_0, coord_1, coord_2);
1790}
1791
1792inline void ExpList::GetCoords(const int eid, Array<OneD, NekDouble> &xc0,
1795{
1796 v_GetCoords(eid, xc0, xc1, xc2);
1797}
1798
1799/**
1800 *
1801 */
1803 const SpatialDomains::GeomMMF MMFdir,
1804 const Array<OneD, const NekDouble> &CircCentre,
1805 Array<OneD, Array<OneD, NekDouble>> &outarray)
1806{
1807 v_GetMovingFrames(MMFdir, CircCentre, outarray);
1808}
1809/**
1810 *
1811 */
1813 Array<OneD, NekDouble> &out_d0,
1814 Array<OneD, NekDouble> &out_d1,
1815 Array<OneD, NekDouble> &out_d2)
1816{
1817 v_PhysDeriv(inarray, out_d0, out_d1, out_d2);
1818}
1819/**
1820 *
1821 */
1822inline void ExpList::PhysDeriv(const int dir,
1823 const Array<OneD, const NekDouble> &inarray,
1825{
1826 v_PhysDeriv(dir, inarray, out_d);
1827}
1829 const Array<OneD, const NekDouble> &inarray,
1831{
1832 v_PhysDeriv(edir, inarray, out_d);
1833}
1834/**
1835 *
1836 */
1838 const Array<OneD, const NekDouble> &direction,
1839 const Array<OneD, const NekDouble> &inarray,
1840 Array<OneD, NekDouble> &outarray)
1841{
1842 v_PhysDirectionalDeriv(direction, inarray, outarray);
1843}
1844/**
1845 *
1846 */
1849{
1850 v_Curl(Vel, Q);
1851}
1852/**
1853 *
1854 */
1857{
1858 v_CurlCurl(Vel, Q);
1859}
1860/**
1861 *
1862 */
1864 const int npts, const Array<OneD, const NekDouble> &inarray,
1865 Array<OneD, NekDouble> &outarray, bool Shuff, bool UnShuff)
1866{
1867 v_HomogeneousFwdTrans(npts, inarray, outarray, Shuff, UnShuff);
1868}
1869/**
1870 *
1871 */
1873 const int npts, const Array<OneD, const NekDouble> &inarray,
1874 Array<OneD, NekDouble> &outarray, bool Shuff, bool UnShuff)
1875{
1876 v_HomogeneousBwdTrans(npts, inarray, outarray, Shuff, UnShuff);
1877}
1878/**
1879 *
1880 */
1881inline void ExpList::DealiasedProd(const int num_dofs,
1882 const Array<OneD, NekDouble> &inarray1,
1883 const Array<OneD, NekDouble> &inarray2,
1884 Array<OneD, NekDouble> &outarray)
1885{
1886 v_DealiasedProd(num_dofs, inarray1, inarray2, outarray);
1887}
1888/**
1889 *
1890 */
1892 const int num_dofs, const Array<OneD, Array<OneD, NekDouble>> &inarray1,
1893 const Array<OneD, Array<OneD, NekDouble>> &inarray2,
1894 Array<OneD, Array<OneD, NekDouble>> &outarray)
1895{
1896 v_DealiasedDotProd(num_dofs, inarray1, inarray2, outarray);
1897}
1898/**
1899 *
1900 */
1902 const Array<OneD, NekDouble> &TotField,
1903 int BndID)
1904{
1905 v_GetBCValues(BndVals, TotField, BndID);
1906}
1907/**
1908 *
1909 */
1912 Array<OneD, NekDouble> &outarray,
1913 int BndID)
1914{
1915 v_NormVectorIProductWRTBase(V1, V2, outarray, BndID);
1916}
1919{
1920 v_NormVectorIProductWRTBase(V, outarray);
1921}
1922/**
1923 * @param eid The index of the element to be checked.
1924 * @return The dimension of the coordinates of the specific element.
1925 */
1926inline int ExpList::GetCoordim(int eid)
1927{
1928 ASSERTL2(eid <= (*m_exp).size(), "eid is larger than number of elements");
1929 return (*m_exp)[eid]->GetCoordim();
1930}
1931/**
1932 * @param eid The index of the element to be checked.
1933 * @return The dimension of the shape of the specific element.
1934 */
1936{
1937 return (*m_exp)[0]->GetShapeDimension();
1938}
1939/**
1940 * @param i The index of m_coeffs to be set
1941 * @param val The value which m_coeffs[i] is to be set to.
1942 */
1943inline void ExpList::SetCoeff(int i, NekDouble val)
1944{
1945 m_coeffs[i] = val;
1946}
1947/**
1948 * @param i The index of #m_coeffs to be set.
1949 * @param val The value which #m_coeffs[i] is to be set to.
1950 */
1951inline void ExpList::SetCoeffs(int i, NekDouble val)
1952{
1953 m_coeffs[i] = val;
1954}
1956{
1957 m_coeffs = inarray;
1958}
1959/**
1960 * As the function returns a constant reference to a
1961 * <em>const Array</em>, it is not possible to modify the
1962 * underlying data of the array #m_coeffs. In order to do
1963 * so, use the function #UpdateCoeffs instead.
1964 *
1965 * @return (A constant reference to) the array #m_coeffs.
1966 */
1968{
1969 return m_coeffs;
1970}
1972{
1974}
1976{
1977 v_FillBndCondFromField(coeffs);
1978}
1979inline void ExpList::FillBndCondFromField(const int nreg,
1980 const Array<OneD, NekDouble> coeffs)
1981{
1982 v_FillBndCondFromField(nreg, coeffs);
1983}
1984inline void ExpList::LocalToGlobal(bool useComm)
1985{
1986 v_LocalToGlobal(useComm);
1987}
1989 Array<OneD, NekDouble> &outarray,
1990 bool useComm)
1991{
1992 v_LocalToGlobal(inarray, outarray, useComm);
1993}
1994inline void ExpList::GlobalToLocal(void)
1995{
1997}
1998/**
1999 * This operation is evaluated as:
2000 * \f{tabbing}
2001 * \hspace{1cm} \= Do \= $e=$ $1, N_{\mathrm{el}}$ \\
2002 * > > Do \= $i=$ $0,N_m^e-1$ \\
2003 * > > > $\boldsymbol{\hat{u}}^{e}[i] = \mbox{sign}[e][i] \cdot
2004 * \boldsymbol{\hat{u}}_g[\mbox{map}[e][i]]$ \\
2005 * > > continue \\
2006 * > continue
2007 * \f}
2008 * where \a map\f$[e][i]\f$ is the mapping array and \a
2009 * sign\f$[e][i]\f$ is an array of similar dimensions ensuring the
2010 * correct modal connectivity between the different elements (both
2011 * these arrays are contained in the data member #m_locToGloMap). This
2012 * operation is equivalent to the scatter operation
2013 * \f$\boldsymbol{\hat{u}}_l=\mathcal{A}\boldsymbol{\hat{u}}_g\f$, where
2014 * \f$\mathcal{A}\f$ is the
2015 * \f$N_{\mathrm{eof}}\times N_{\mathrm{dof}}\f$ permutation matrix.
2016 *
2017 * @param inarray An array of size \f$N_\mathrm{dof}\f$
2018 * containing the global degrees of freedom
2019 * \f$\boldsymbol{x}_g\f$.
2020 * @param outarray The resulting local degrees of freedom
2021 * \f$\boldsymbol{x}_l\f$ will be stored in this
2022 * array of size \f$N_\mathrm{eof}\f$.
2023 */
2025 Array<OneD, NekDouble> &outarray)
2026{
2027 v_GlobalToLocal(inarray, outarray);
2028}
2029/**
2030 * @param i The index of #m_coeffs to be returned
2031 * @return The NekDouble held in #m_coeffs[i].
2032 */
2034{
2035 return m_coeffs[i];
2036}
2037/**
2038 * @param i The index of #m_coeffs to be returned
2039 * @return The NekDouble held in #m_coeffs[i].
2040 */
2042{
2043 return m_coeffs[i];
2044}
2045/**
2046 * As the function returns a constant reference to a
2047 * <em>const Array</em> it is not possible to modify the
2048 * underlying data of the array #m_phys. In order to do
2049 * so, use the function #UpdatePhys instead.
2050 *
2051 * @return (A constant reference to) the array #m_phys.
2052 */
2054{
2055 return m_phys;
2056}
2057/**
2058 * @return \f$N_{\mathrm{el}}\f$, the number of elements in the
2059 * expansion.
2060 */
2061inline int ExpList::GetExpSize(void)
2062{
2063 return (*m_exp).size();
2064}
2065/**
2066 * @param n The index of the element concerned.
2067 *
2068 * @return (A shared pointer to) the local expansion of the
2069 * \f$n^{\mathrm{th}}\f$ element.
2070 */
2072{
2073 return (*m_exp)[n];
2074}
2075/**
2076 * @param n The global id of the element concerned.
2077 *
2078 * @return (A shared pointer to) the local expansion of the
2079 * \f$n^{\mathrm{th}}\f$ element.
2080 */
2082{
2083 auto it = m_elmtToExpId.find(n);
2084 ASSERTL0(it != m_elmtToExpId.end(), "Global geometry ID " +
2085 std::to_string(n) +
2086 " not found in element ID to "
2087 "expansion ID map.")
2088 return (*m_exp)[it->second];
2089}
2090
2091/**
2092 * @return (A const shared pointer to) the local expansion vector #m_exp
2093 */
2094inline const std::shared_ptr<LocalRegions::ExpansionVector> ExpList::GetExp(
2095 void) const
2096{
2097 return m_exp;
2098}
2099/**
2100 *
2101 */
2102inline int ExpList::GetCoeff_Offset(int n) const
2103{
2104 return m_coeff_offset[n];
2105}
2106/**
2107 *
2108 */
2109inline int ExpList::GetPhys_Offset(int n) const
2110{
2111 return m_phys_offset[n];
2112}
2113/**
2114 * If one wants to get hold of the underlying data without modifying
2115 * them, rather use the function #GetCoeffs instead.
2116 *
2117 * @return (A reference to) the array #m_coeffs.
2118 */
2120{
2121 return m_coeffs;
2122}
2123/**
2124 * If one wants to get hold of the underlying data without modifying
2125 * them, rather use the function #GetPhys instead.
2126 *
2127 * @return (A reference to) the array #m_phys.
2128 */
2130{
2131 m_physState = true;
2132 return m_phys;
2133}
2134// functions associated with DisContField
2137{
2138 return v_GetBndCondExpansions();
2139}
2140/// Get m_coeffs to elemental value map
2143{
2144 return m_coeffsToElmt;
2145}
2148{
2149 return v_GetLocTraceToTraceMap();
2150}
2152{
2153 return v_GetBndCondBwdWeight();
2154}
2155inline void ExpList::SetBndCondBwdWeight(const int index, const NekDouble value)
2156{
2157 v_SetBndCondBwdWeight(index, value);
2158}
2159inline std::shared_ptr<ExpList> &ExpList::UpdateBndCondExpansion(int i)
2160{
2161 return v_UpdateBndCondExpansion(i);
2162}
2164 const Array<OneD, const Array<OneD, NekDouble>> &Vec,
2167{
2168 v_Upwind(Vec, Fwd, Bwd, Upwind);
2169}
2173 Array<OneD, NekDouble> &Upwind)
2174{
2175 v_Upwind(Vn, Fwd, Bwd, Upwind);
2176}
2177inline std::shared_ptr<ExpList> &ExpList::GetTrace()
2178{
2179 return v_GetTrace();
2180}
2181inline std::shared_ptr<AssemblyMapDG> &ExpList::GetTraceMap()
2182{
2183 return v_GetTraceMap();
2184}
2186{
2187 return v_GetTraceBndMap();
2188}
2190{
2191 v_GetNormals(normals);
2192}
2195 Array<OneD, NekDouble> &outarray)
2196{
2197 WARNINGL1(false,
2198 "Deprecated AddTraceIntegral interface, will be removed in "
2199 "the next release");
2200 v_AddTraceIntegral(Fx, Fy, outarray);
2201}
2203 Array<OneD, NekDouble> &outarray)
2204{
2205 v_AddTraceIntegral(Fn, outarray);
2206}
2210{
2211 v_AddFwdBwdTraceIntegral(Fwd, Bwd, outarray);
2212}
2215{
2216 v_GetFwdBwdTracePhys(Fwd, Bwd);
2217}
2220 Array<OneD, NekDouble> &Bwd, bool FillBnd, bool PutFwdInBwdOnBCs,
2221 bool DoExchange)
2222{
2223 v_GetFwdBwdTracePhys(field, Fwd, Bwd, FillBnd, PutFwdInBwdOnBCs,
2224 DoExchange);
2225}
2228 bool PutFwdInBwdOnBCs)
2229{
2230 v_FillBwdWithBoundCond(Fwd, Bwd, PutFwdInBwdOnBCs);
2231}
2235{
2236 v_AddTraceQuadPhysToField(Fwd, Bwd, field);
2237}
2241 Array<OneD, NekDouble> &locTraceFwd, Array<OneD, NekDouble> &locTraceBwd)
2242{
2243 v_GetLocTraceFromTracePts(Fwd, Bwd, locTraceFwd, locTraceBwd);
2244}
2246 Array<OneD, NekDouble> &weightjmp)
2247{
2248 v_FillBwdWithBwdWeight(weightave, weightjmp);
2249}
2252{
2253 v_PeriodicBwdCopy(Fwd, Bwd);
2254}
2255inline const std::vector<bool> &ExpList::GetLeftAdjacentFaces(void) const
2256{
2257 return v_GetLeftAdjacentFaces();
2258}
2260{
2261 v_ExtractTracePhys(outarray);
2262}
2264 const Array<OneD, const NekDouble> &inarray,
2265 Array<OneD, NekDouble> &outarray)
2266{
2267 v_ExtractTracePhys(inarray, outarray);
2268}
2271{
2272 return v_GetBndConditions();
2273}
2276{
2277 return v_UpdateBndConditions();
2278}
2280 const std::string varName,
2281 const NekDouble x2_in,
2282 const NekDouble x3_in)
2283{
2284 v_EvaluateBoundaryConditions(time, varName, x2_in, x3_in);
2285}
2287{
2289}
2291 Array<OneD, int> &EdgeID)
2292{
2293 v_GetBoundaryToElmtMap(ElmtID, EdgeID);
2294}
2296 std::shared_ptr<ExpList> &result,
2297 const bool DeclareCoeffPhysArrays)
2298{
2299 v_GetBndElmtExpansion(i, result, DeclareCoeffPhysArrays);
2300}
2301
2303 const Array<OneD, NekDouble> &elmt,
2304 Array<OneD, NekDouble> &boundary)
2305{
2306 v_ExtractElmtToBndPhys(i, elmt, boundary);
2307}
2308
2310 int i, const Array<OneD, const NekDouble> &phys,
2311 Array<OneD, NekDouble> &bndElmt)
2312{
2313 v_ExtractPhysToBndElmt(i, phys, bndElmt);
2314}
2315
2317 const Array<OneD, const NekDouble> &phys,
2319{
2320 v_ExtractPhysToBnd(i, phys, bnd);
2321}
2322
2324 int i, Array<OneD, Array<OneD, NekDouble>> &normals)
2325{
2326 v_GetBoundaryNormals(i, normals);
2327}
2328
2329inline std::vector<bool> &ExpList::GetLeftAdjacentTraces(void)
2330{
2331 return v_GetLeftAdjacentTraces();
2332}
2333
2335
2336} // namespace MultiRegions
2337} // namespace Nektar
2338
2339#endif // EXPLIST_H
#define WARNINGL1(condition, msg)
Definition: ErrorUtil.hpp:250
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed to...
Definition: ErrorUtil.hpp:272
#define MULTI_REGIONS_EXPORT
Base class for all multi-elemental spectral/hp expansions.
Definition: ExpList.h:102
virtual void v_Reset()
Reset geometry information, metrics, matrix managers and geometry information.
Definition: ExpList.cpp:3071
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:1092
Array< OneD, NekDouble > & UpdateCoeffs()
This function returns (a reference to) the array (implemented as m_coeffs) containing all local expa...
Definition: ExpList.h:2119
void AddTraceJacToElmtJac(const Array< OneD, const DNekMatSharedPtr > &FwdMat, const Array< OneD, const DNekMatSharedPtr > &BwdMat, Array< OneD, DNekMatSharedPtr > &fieldMat)
inverse process of v_GetFwdBwdTracePhys. Given Trace integration of Fwd and Bwd Jacobian,...
Definition: ExpList.cpp:5989
void LinearAdvectionReactionSolve(const Array< OneD, Array< OneD, NekDouble > > &velocity, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const NekDouble lambda, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
Solve Advection Diffusion Reaction.
Definition: ExpList.h:1773
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:1665
const Array< OneD, const std::shared_ptr< ExpList > > & GetBndCondExpansions()
Definition: ExpList.h:2136
std::shared_ptr< ExpList > & GetTrace()
Definition: ExpList.h:2177
void WriteTecplotConnectivity(std::ostream &outfile, int expansion=-1)
Definition: ExpList.h:412
int GetNcoeffs(void) const
Returns the total number of local degrees of freedom .
Definition: ExpList.h:1527
virtual void v_ExtractCoeffsToCoeffs(const std::shared_ptr< ExpList > &fromExpList, const Array< OneD, const NekDouble > &fromCoeffs, Array< OneD, NekDouble > &toCoeffs)
Definition: ExpList.cpp:4073
const DNekScalBlkMatSharedPtr & GetBlockMatrix(const GlobalMatrixKey &gkey)
Definition: ExpList.cpp:2415
virtual void v_AddFwdBwdTraceIntegral(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:4795
LibUtilities::NekManager< GlobalLinSysKey, GlobalLinSys > & GetGlobalLinSysManager(void)
Definition: ExpList.cpp:5685
ExpList(const ExpansionType Type=eNoType)
The default constructor using a type.
Definition: ExpList.cpp:102
void SetCoeffs(int i, NekDouble val)
Set the i th coefficiient in m_coeffs to value val.
Definition: ExpList.h:1951
void CreateCollections(Collections::ImplementationType ImpType=Collections::eNoImpType)
Construct collections of elements containing a single element type and polynomial order from the list...
Definition: ExpList.cpp:5503
virtual void v_GetBoundaryToElmtMap(Array< OneD, int > &ElmtID, Array< OneD, int > &EdgeID)
Definition: ExpList.cpp:5393
void ExponentialFilter(Array< OneD, NekDouble > &array, const NekDouble alpha, const NekDouble exponent, const NekDouble cutoff)
Definition: ExpList.cpp:2149
std::shared_ptr< GlobalMatrix > GenGlobalMatrix(const GlobalMatrixKey &mkey, const std::shared_ptr< AssemblyMapCG > &locToGloMap)
Generates a global matrix from the given key and map.
Definition: ExpList.cpp:2530
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:5118
void MultiplyByMassMatrix(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:307
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:2226
void AddRightIPTPhysDerivBase(const int dir, const Array< OneD, const DNekMatSharedPtr > ElmtJacQuad, Array< OneD, DNekMatSharedPtr > ElmtJacCoef)
Definition: ExpList.cpp:6225
virtual void v_ExtractTracePhys(Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:4883
void ExtractPhysToBndElmt(int i, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &bndElmt)
Definition: ExpList.h:2309
void Reset()
Reset geometry information and reset matrices.
Definition: ExpList.h:395
int GetExpIndex(const Array< OneD, const NekDouble > &gloCoord, NekDouble tol=0.0, bool returnNearestElmt=false, int cachedId=-1, NekDouble maxDistance=1e6)
This function returns the index of the local elemental expansion containing the arbitrary point given...
Definition: ExpList.cpp:2920
virtual GlobalLinSysKey v_HelmSolve(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::ConstFactorMap &factors, const StdRegions::VarCoeffMap &varcoeff, const MultiRegions::VarFactorsMap &varfactors, const Array< OneD, const NekDouble > &dirForcing, const bool PhysSpaceForcing)
Definition: ExpList.cpp:4907
void GetFieldDefinitions(std::vector< LibUtilities::FieldDefinitionsSharedPtr > &fielddef)
Definition: ExpList.h:912
virtual void v_WriteTecplotConnectivity(std::ostream &outfile, int expansion)
Definition: ExpList.cpp:3229
BlockMatrixMapShPtr m_blockMat
Definition: ExpList.h:1141
virtual Array< OneD, const unsigned int > v_GetYIDs(void)
Definition: ExpList.cpp:3693
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
multiply the metric jacobi and quadrature weights
Definition: ExpList.cpp:1693
virtual void v_LinearAdvectionReactionSolve(const Array< OneD, Array< OneD, NekDouble > > &velocity, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const NekDouble lambda, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
Definition: ExpList.cpp:4934
void SetPhys(int i, NekDouble val)
Set the i th value of m_phys to value val.
Definition: ExpList.h:1625
void ExtractTracePhys(Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:2259
void GetBoundaryToElmtMap(Array< OneD, int > &ElmtID, Array< OneD, int > &EdgeID)
Definition: ExpList.h:2290
void PeriodicBwdCopy(const Array< OneD, const NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
Copy and fill the Periodic boundaries.
Definition: ExpList.h:2250
void WriteVtkPieceFooter(std::ostream &outfile, int expansion)
Definition: ExpList.cpp:3490
virtual void v_ExtractPhysToBnd(const int i, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &bnd)
Definition: ExpList.cpp:5312
void GetBwdWeight(Array< OneD, NekDouble > &weightAver, Array< OneD, NekDouble > &weightJump)
Get the weight value for boundary conditions for boundary average and jump calculations.
Definition: ExpList.cpp:4102
virtual void v_AddTraceIntegralToOffDiag(const Array< OneD, const NekDouble > &FwdFlux, const Array< OneD, const NekDouble > &BwdFlux, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:5766
virtual std::vector< bool > & v_GetLeftAdjacentTraces(void)
Definition: ExpList.cpp:4335
void InitialiseExpVector(const SpatialDomains::ExpansionInfoMap &expmap)
Define a list of elements using the geometry and basis key information in expmap;.
Definition: ExpList.cpp:1488
void SetupCoeffPhys(bool DeclareCoeffPhysArrays=true, bool SetupOffsets=true)
Definition of the total number of degrees of freedom and quadrature points and offsets to access data...
Definition: ExpList.cpp:1445
virtual LibUtilities::BasisSharedPtr v_GetHomogeneousBasis(void)
Definition: ExpList.h:1479
static SpatialDomains::BoundaryConditionShPtr GetBoundaryCondition(const SpatialDomains::BoundaryConditionCollection &collection, unsigned int index, const std::string &variable)
Definition: ExpList.cpp:5470
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:1967
virtual void GetBndElmtExpansion(int i, std::shared_ptr< ExpList > &result, const bool DeclareCoeffPhysArrays=true)
Definition: ExpList.h:2295
void SetPhysArray(Array< OneD, NekDouble > &inarray)
Sets the array m_phys.
Definition: ExpList.h:1644
virtual void v_PhysDirectionalDeriv(const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:2119
virtual void v_GetBCValues(Array< OneD, NekDouble > &BndVals, const Array< OneD, NekDouble > &TotField, int BndID)
Definition: ExpList.cpp:4985
void FillBndCondFromField(const Array< OneD, NekDouble > coeffs)
Fill Bnd Condition expansion from the values stored in expansion.
Definition: ExpList.h:1975
const Array< OneD, const std::pair< int, int > > & GetCoeffsToElmt() const
Get m_coeffs to elemental value map.
Definition: ExpList.h:2142
virtual void v_GetLocTraceFromTracePts(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &locTraceFwd, Array< OneD, NekDouble > &locTraceBwd)
Definition: ExpList.cpp:4850
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:2053
std::shared_ptr< ExpList > & UpdateBndCondExpansion(int i)
Definition: ExpList.h:2159
virtual void v_WriteTecplotHeader(std::ostream &outfile, std::string var="")
Definition: ExpList.cpp:3098
virtual std::shared_ptr< ExpList > & v_GetPlane(int n)
Definition: ExpList.cpp:5491
void NormVectorIProductWRTBase(Array< OneD, const NekDouble > &V1, Array< OneD, const NekDouble > &V2, Array< OneD, NekDouble > &outarray, int BndID)
Definition: ExpList.h:1910
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:2102
int GetExpSize(void)
This function returns the number of elements in the expansion.
Definition: ExpList.h:2061
virtual void v_AddTraceQuadPhysToOffDiag(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &field)
Definition: ExpList.cpp:4841
std::shared_ptr< LibUtilities::SessionReader > GetSession() const
Returns the session object.
Definition: ExpList.h:971
virtual void v_ClearGlobalLinSysManager(void)
Definition: ExpList.cpp:3701
virtual void v_UnsetGlobalLinSys(GlobalLinSysKey, bool)
Definition: ExpList.cpp:3714
void GeneralGetFieldDefinitions(std::vector< LibUtilities::FieldDefinitionsSharedPtr > &fielddef, int NumHomoDir=0, Array< OneD, LibUtilities::BasisSharedPtr > &HomoBasis=LibUtilities::NullBasisSharedPtr1DArray, std::vector< NekDouble > &HomoLen=LibUtilities::NullNekDoubleVector, bool homoStrips=false, std::vector< unsigned int > &HomoSIDs=LibUtilities::NullUnsignedIntVector, std::vector< unsigned int > &HomoZIDs=LibUtilities::NullUnsignedIntVector, std::vector< unsigned int > &HomoYIDs=LibUtilities::NullUnsignedIntVector)
Definition: ExpList.cpp:3802
virtual void v_GetPeriodicEntities(PeriodicMap &periodicVerts, PeriodicMap &periodicEdges, PeriodicMap &periodicFaces)
Definition: ExpList.cpp:5461
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:1674
void SetCoeffsArray(Array< OneD, NekDouble > &inarray)
Set the m_coeffs array to inarray.
Definition: ExpList.h:1955
virtual void v_CurlCurl(Array< OneD, Array< OneD, NekDouble > > &Vel, Array< OneD, Array< OneD, NekDouble > > &Q)
Definition: ExpList.cpp:2051
LibUtilities::TranspositionSharedPtr GetTransposition(void)
This function returns the transposition class associated with the homogeneous expansion.
Definition: ExpList.h:601
virtual void v_LocalToGlobal(bool UseComm)
Definition: ExpList.cpp:5072
void SetBndCondBwdWeight(const int index, const NekDouble value)
Set the weight value for boundary conditions.
Definition: ExpList.h:2155
virtual const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > & v_GetBndConditions()
Definition: ExpList.cpp:5418
virtual LibUtilities::NekManager< GlobalLinSysKey, GlobalLinSys > & v_GetGlobalLinSysManager(void)
Definition: ExpList.cpp:3722
void SetHomoLen(const NekDouble lhom)
This function sets the Width of homogeneous direction associated with the homogeneous expansion.
Definition: ExpList.h:615
virtual void v_DealiasedDotProd(const int num_dofs, const Array< OneD, Array< OneD, NekDouble > > &inarray1, const Array< OneD, Array< OneD, NekDouble > > &inarray2, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: ExpList.cpp:4975
int GetNpoints(void) const
Returns the total number of quadrature points m_npoints .
Definition: ExpList.h:1603
SpatialDomains::MeshGraphSharedPtr GetGraph()
Definition: ExpList.h:980
virtual ~ExpList()
The default destructor.
Definition: ExpList.cpp:1661
void GetFwdBwdTracePhys(Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
Definition: ExpList.h:2213
int GetElmtToExpId(int elmtId)
This function returns the index inside m_exp for a given geom id.
Definition: ExpList.h:1048
const Array< OneD, const int > & GetTraceBndMap(void)
Definition: ExpList.h:2185
virtual void v_FwdTransLocalElmt(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:2208
int GetPoolCount(std::string)
Definition: ExpList.cpp:5674
virtual void v_ImposeDirichletConditions(Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:5046
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:1702
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:594
virtual void v_EvaluateBoundaryConditions(const NekDouble time=0.0, const std::string varName="", const NekDouble x2_in=NekConstants::kNekUnsetDouble, const NekDouble x3_in=NekConstants::kNekUnsetDouble)
Definition: ExpList.cpp:5439
void FwdTransBndConstrained(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:1711
virtual size_t v_GetNumElmts(void)
Definition: ExpList.h:1175
virtual const std::vector< bool > & v_GetLeftAdjacentFaces(void) const
Definition: ExpList.cpp:4875
Array< OneD, int > m_coeff_offset
Offset of elemental data into the array m_coeffs.
Definition: ExpList.h:1136
int EvalBasisNumModesMax(void) const
Evaulates the maximum number of modes in the elemental basis order over all elements.
Definition: ExpList.h:1541
void AddFwdBwdTraceIntegral(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:2207
std::shared_ptr< AssemblyMapDG > & GetTraceMap(void)
Definition: ExpList.h:2181
void WriteTecplotField(std::ostream &outfile, int expansion=-1)
Definition: ExpList.h:408
virtual void v_GetMovingFrames(const SpatialDomains::GeomMMF MMFdir, const Array< OneD, const NekDouble > &CircCentre, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: ExpList.cpp:4115
size_t GetNumElmts(void)
This function returns the number of elements in the expansion which may be different for a homogeoeno...
Definition: ExpList.h:657
virtual std::shared_ptr< AssemblyMapDG > & v_GetTraceMap()
Definition: ExpList.cpp:4322
Array< OneD, std::pair< int, int > > m_coeffsToElmt
m_coeffs to elemental value map
Definition: ExpList.h:1140
virtual void v_GetFwdBwdTracePhys(Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
Definition: ExpList.cpp:4804
NekDouble PhysEvaluate(const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &phys)
Definition: ExpList.cpp:3036
const std::unordered_map< int, int > & GetElmtToExpId(void)
This function returns the map of index inside m_exp to geom id.
Definition: ExpList.h:1042
virtual GlobalLinSysKey v_LinearAdvectionDiffusionReactionSolve(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::ConstFactorMap &factors, const StdRegions::VarCoeffMap &varcoeff, const MultiRegions::VarFactorsMap &varfactors, const Array< OneD, const NekDouble > &dirForcing, const bool PhysSpaceForcing)
Definition: ExpList.cpp:4920
virtual NekDouble v_Integral(const Array< OneD, const NekDouble > &inarray)
Definition: ExpList.cpp:3620
virtual void v_ExtractDataToCoeffs(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, std::string &field, Array< OneD, NekDouble > &coeffs, std::unordered_map< int, int > zIdToPlane)
Extract data from raw field data into expansion list.
Definition: ExpList.cpp:3985
void PhysDirectionalDeriv(const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:1837
virtual std::map< int, RobinBCInfoSharedPtr > v_GetRobinBCInfo(void)
Definition: ExpList.cpp:5451
virtual const Array< OneD, const int > & v_GetTraceBndMap()
Definition: ExpList.cpp:4330
virtual void v_GetBoundaryNormals(int i, Array< OneD, Array< OneD, NekDouble > > &normals)
Definition: ExpList.cpp:5349
void GetBoundaryNormals(int i, Array< OneD, Array< OneD, NekDouble > > &normals)
Definition: ExpList.h:2323
virtual Array< OneD, const NekDouble > v_HomogeneousEnergy(void)
Definition: ExpList.cpp:3654
virtual NekDouble v_GetHomoLen(void)
Definition: ExpList.cpp:3670
void ExtractFileBCs(const std::string &fileName, LibUtilities::CommSharedPtr comm, const std::string &varName, const std::shared_ptr< ExpList > locExpList)
void GetMovingFrames(const SpatialDomains::GeomMMF MMFdir, const Array< OneD, const NekDouble > &CircCentre, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: ExpList.h:1802
const LocTraceToTraceMapSharedPtr & GetLocTraceToTraceMap() const
Definition: ExpList.h:2147
virtual Array< OneD, const unsigned int > v_GetZIDs(void)
Definition: ExpList.cpp:3685
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:578
void SetHomo1DSpecVanVisc(Array< OneD, NekDouble > visc)
This function sets the Spectral Vanishing Viscosity in homogeneous1D expansion.
Definition: ExpList.h:585
void ExtractCoeffsToCoeffs(const std::shared_ptr< ExpList > &fromExpList, const Array< OneD, const NekDouble > &fromCoeffs, Array< OneD, NekDouble > &toCoeffs)
Extract the data from fromField using fromExpList the coeffs using the basic ExpList Elemental expans...
Definition: ExpList.cpp:3969
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:1785
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:1556
std::shared_ptr< GlobalLinSys > GenGlobalBndLinSys(const GlobalLinSysKey &mkey, const AssemblyMapSharedPtr &locToGloMap)
Generate a GlobalLinSys from information provided by the key "mkey" and the mapping provided in LocTo...
Definition: ExpList.cpp:2833
virtual void v_SetBndCondBwdWeight(const int index, const NekDouble value)
Definition: ExpList.cpp:4868
void UnsetGlobalLinSys(GlobalLinSysKey, bool)
Definition: ExpList.cpp:5679
void DealiasedProd(const int num_dofs, const Array< OneD, NekDouble > &inarray1, const Array< OneD, NekDouble > &inarray2, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:1881
void GetNormals(Array< OneD, Array< OneD, NekDouble > > &normals)
Definition: ExpList.h:2189
void ExtractPhysToBnd(int i, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &bnd)
Definition: ExpList.h:2316
virtual void v_SmoothField(Array< OneD, NekDouble > &field)
Definition: ExpList.cpp:2240
virtual void v_AppendFieldData(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata)
Definition: ExpList.cpp:3929
std::vector< bool > m_collectionsDoInit
Vector of bools to act as an initialise on first call flag.
Definition: ExpList.h:1130
NekDouble GetCoeff(int i)
Get the i th value (coefficient) of m_coeffs.
Definition: ExpList.h:2033
bool m_physState
The state of the array m_phys.
Definition: ExpList.h:1116
virtual void v_HomogeneousFwdTrans(const int npts, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
Definition: ExpList.cpp:4945
void ExtractCoeffsFromFile(const std::string &fileName, LibUtilities::CommSharedPtr comm, const std::string &varName, Array< OneD, NekDouble > &coeffs)
Definition: ExpList.cpp:3729
NekDouble Integral()
Calculates the error of the global spectral/hp element approximation.
Definition: ExpList.h:547
std::shared_ptr< ExpList > & GetPlane(int n)
Definition: ExpList.h:989
virtual void v_DealiasedProd(const int num_dofs, const Array< OneD, NekDouble > &inarray1, const Array< OneD, NekDouble > &inarray2, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:4965
void AddRightIPTBaseMatrix(const Array< OneD, const DNekMatSharedPtr > ElmtJacQuad, Array< OneD, DNekMatSharedPtr > ElmtJacCoef)
Definition: ExpList.cpp:6282
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:2238
virtual void v_WriteVtkPieceData(std::ostream &outfile, int expansion, std::string var)
Definition: ExpList.cpp:3497
virtual void v_WriteVtkPieceHeader(std::ostream &outfile, int expansion, int istrip)
Definition: ExpList.cpp:3356
std::shared_ptr< ExpList > GetSharedThisPtr()
Returns a shared pointer to the current object.
Definition: ExpList.h:966
int GetShapeDimension()
This function returns the dimension of the shape of the element eid.
Definition: ExpList.h:1935
void ExtractElmtDataToCoeffs(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, std::string &field, Array< OneD, NekDouble > &coeffs)
Extract the data in fielddata into the coeffs using the basic ExpList Elemental expansions rather tha...
NekDouble H1(const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray)
Calculates the error of the global spectral/hp element approximation.
Definition: ExpList.cpp:3784
std::vector< int > m_coll_phys_offset
Offset of elemental data into the array m_phys.
Definition: ExpList.h:1134
LibUtilities::CommSharedPtr m_comm
Communicator.
Definition: ExpList.h:1065
std::shared_ptr< GlobalLinSys > GenGlobalLinSys(const GlobalLinSysKey &mkey, const std::shared_ptr< AssemblyMapCG > &locToGloMap)
This operation constructs the global linear system of type mkey.
Definition: ExpList.cpp:2815
void WriteVtkHeader(std::ostream &outfile)
Definition: ExpList.cpp:3342
bool GetPhysState(void) const
This function indicates whether the array of physical values (implemented as m_phys) is filled or no...
Definition: ExpList.h:1658
std::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:1127
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Definition: ExpList.h:1072
void GenerateElementVector(const int ElementID, const NekDouble scalar1, const NekDouble scalar2, Array< OneD, NekDouble > &outarray)
Generate vector v such that v[i] = scalar1 if i is in the element < ElementID. Otherwise,...
Definition: ExpList.cpp:4159
NekDouble VectorFlux(const Array< OneD, Array< OneD, NekDouble > > &inarray)
Definition: ExpList.h:572
void AddTraceIntegralToOffDiag(const Array< OneD, const NekDouble > &FwdFlux, const Array< OneD, const NekDouble > &BwdFlux, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:1020
std::shared_ptr< DNekMat > GenGlobalMatrixFull(const GlobalLinSysKey &mkey, const std::shared_ptr< AssemblyMapCG > &locToGloMap)
Definition: ExpList.cpp:2663
const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > & GetBndConditions()
Definition: ExpList.h:2270
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:2081
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:1651
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:644
virtual void v_NormVectorIProductWRTBase(Array< OneD, const NekDouble > &V1, Array< OneD, const NekDouble > &V2, Array< OneD, NekDouble > &outarray, int BndID)
Definition: ExpList.cpp:4993
virtual void v_ExtractElmtToBndPhys(const int i, const Array< OneD, NekDouble > &elmt, Array< OneD, NekDouble > &boundary)
Definition: ExpList.cpp:5235
virtual void v_ExtractPhysToBndElmt(const int i, const Array< OneD, const NekDouble > &phys, Array< OneD, NekDouble > &bndElmt)
Definition: ExpList.cpp:5273
void CurlCurl(Array< OneD, Array< OneD, NekDouble > > &Vel, Array< OneD, Array< OneD, NekDouble > > &Q)
Definition: ExpList.h:1855
std::map< int, RobinBCInfoSharedPtr > GetRobinBCInfo()
Definition: ExpList.h:895
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:1749
virtual void v_Upwind(const Array< OneD, const Array< OneD, NekDouble > > &Vec, const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &Upwind)
Definition: ExpList.cpp:4216
void EvaluateBoundaryConditions(const NekDouble time=0.0, const std::string varName="", const NekDouble=NekConstants::kNekUnsetDouble, const NekDouble=NekConstants::kNekUnsetDouble)
Definition: ExpList.h:2279
std::vector< LibUtilities::FieldDefinitionsSharedPtr > GetFieldDefinitions()
Definition: ExpList.h:907
virtual void v_AddTraceIntegral(const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:4778
virtual void v_FillBwdWithBwdWeight(Array< OneD, NekDouble > &weightave, Array< OneD, NekDouble > &weightjmp)
Definition: ExpList.cpp:5401
Array< OneD, NekDouble > & UpdatePhys()
This function returns (a reference to) the array (implemented as m_phys) containing the function ev...
Definition: ExpList.h:2129
virtual void v_Curl(Array< OneD, Array< OneD, NekDouble > > &Vel, Array< OneD, Array< OneD, NekDouble > > &Q)
Definition: ExpList.cpp:2002
virtual void v_PhysGalerkinProjection1DScaled(const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:6339
void GeneralMatrixOp(const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function calculates the result of the multiplication of a matrix of type specified by mkey with ...
Definition: ExpList.cpp:2459
virtual void v_SetUpPhysNormals()
: Set up a normal along the trace elements between two elements at elemental level
Definition: ExpList.cpp:5212
virtual void v_GlobalToLocal(void)
Definition: ExpList.cpp:5087
virtual LibUtilities::TranspositionSharedPtr v_GetTransposition(void)
Definition: ExpList.cpp:3662
void ExtractElmtToBndPhys(int i, const Array< OneD, NekDouble > &elmt, Array< OneD, NekDouble > &boundary)
Definition: ExpList.h:2302
void HomogeneousBwdTrans(const int npts, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
Definition: ExpList.h:1872
void SetWaveSpace(const bool wavespace)
Sets the wave space to the one of the possible configuration true or false.
Definition: ExpList.h:1611
const std::shared_ptr< LocalRegions::ExpansionVector > GetExp() const
This function returns the vector of elements in the expansion.
Definition: ExpList.h:2094
SpatialDomains::MeshGraphSharedPtr m_graph
Mesh associated with this expansion list.
Definition: ExpList.h:1069
virtual void v_MultiplyByInvMassMatrix(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:4898
virtual void v_AddTraceQuadPhysToField(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &field)
Definition: ExpList.cpp:4832
virtual const Array< OneD, const NekDouble > & v_GetBndCondBwdWeight()
Definition: ExpList.cpp:4860
void GlobalToLocal(void)
Scatters from the global coefficients to the local coefficients .
Definition: ExpList.h:1994
virtual NekDouble v_L2(const Array< OneD, const NekDouble > &phys, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray)
Definition: ExpList.cpp:3575
const Array< OneD, const NekDouble > & GetBndCondBwdWeight()
Get the weight value for boundary conditions.
Definition: ExpList.h:2151
void WriteVtkFooter(std::ostream &outfile)
Definition: ExpList.cpp:3350
void MultiplyByBlockMatrix(const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:1673
virtual std::shared_ptr< ExpList > & v_GetTrace()
Definition: ExpList.cpp:4314
void GetElmtNormalLength(Array< OneD, NekDouble > &lengthsFwd, Array< OneD, NekDouble > &lengthsBwd)
Get the length of elements in boundary normal direction.
Definition: ExpList.cpp:4630
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:1582
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:1762
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:633
virtual NekDouble v_VectorFlux(const Array< OneD, Array< OneD, NekDouble > > &inarray)
Definition: ExpList.cpp:3634
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:518
virtual std::vector< LibUtilities::FieldDefinitionsSharedPtr > v_GetFieldDefinitions(void)
Definition: ExpList.cpp:3914
std::vector< int > m_coll_coeff_offset
Offset of elemental data into the array m_coeffs.
Definition: ExpList.h:1132
void IProductWRTDirectionalDerivBase(const Array< OneD, const NekDouble > &direction, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Directional derivative along a given direction.
Definition: ExpList.cpp:1757
void WriteVtkPieceData(std::ostream &outfile, int expansion, std::string var="v")
Definition: ExpList.h:425
void FillBwdWithBwdWeight(Array< OneD, NekDouble > &weightave, Array< OneD, NekDouble > &weightjmp)
Fill Bwd with boundary conditions.
Definition: ExpList.h:2245
virtual void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:2869
virtual Array< OneD, SpatialDomains::BoundaryConditionShPtr > & v_UpdateBndConditions()
Definition: ExpList.cpp:5429
void SetCoeff(int i, NekDouble val)
Set the i th coefficiient in m_coeffs to value val.
Definition: ExpList.h:1943
void GetMatIpwrtDeriveBase(const Array< OneD, const Array< OneD, NekDouble > > &inarray, const int nDirctn, Array< OneD, DNekMatSharedPtr > &mtxPerVar)
Definition: ExpList.cpp:5775
void LocalToGlobal(bool useComm=true)
Gathers the global coefficients from the local coefficients .
Definition: ExpList.h:1984
void GetPeriodicEntities(PeriodicMap &periodicVerts, PeriodicMap &periodicEdges, PeriodicMap &periodicFaces=NullPeriodicMap)
Definition: ExpList.h:900
NekDouble Integral(const Array< OneD, const NekDouble > &inarray)
Definition: ExpList.h:568
virtual void v_GetCoords(Array< OneD, NekDouble > &coord_0, Array< OneD, NekDouble > &coord_1=NullNekDouble1DArray, Array< OneD, NekDouble > &coord_2=NullNekDouble1DArray)
Definition: ExpList.cpp:5152
virtual void v_WriteTecplotField(std::ostream &outfile, int expansion)
Definition: ExpList.cpp:3309
virtual void v_FillBndCondFromField(const Array< OneD, NekDouble > coeffs)
Definition: ExpList.cpp:5055
Collections::CollectionVector m_collections
Definition: ExpList.h:1128
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:2232
virtual int v_GetPoolCount(std::string)
Definition: ExpList.cpp:3707
virtual void v_PhysInterp1DScaled(const NekDouble scale, const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:5690
virtual void v_HomogeneousBwdTrans(const int npts, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
Definition: ExpList.cpp:4955
void ApplyGeomInfo()
Apply geometry information to each expansion.
Definition: ExpList.cpp:3059
virtual void v_PhysDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d0, Array< OneD, NekDouble > &out_d1, Array< OneD, NekDouble > &out_d2)
Definition: ExpList.cpp:1901
virtual void v_FwdTransBndConstrained(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:2217
std::unordered_map< int, int > m_elmtToExpId
Mapping from geometry ID of element to index inside m_exp.
Definition: ExpList.h:1148
void FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:1694
void SmoothField(Array< OneD, NekDouble > &field)
Smooth a field across elements.
Definition: ExpList.h:1720
virtual void v_WriteTecplotZone(std::ostream &outfile, int expansion)
Definition: ExpList.cpp:3136
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:624
LibUtilities::SessionReaderSharedPtr m_session
Session.
Definition: ExpList.h:1067
virtual void v_PeriodicBwdCopy(const Array< OneD, const NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd)
Definition: ExpList.cpp:5408
void MultiplyByInvMassMatrix(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:1740
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:927
void AddTraceIntegral(const Array< OneD, const NekDouble > &Fx, const Array< OneD, const NekDouble > &Fy, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:2193
virtual void v_FillBwdWithBoundCond(const Array< OneD, NekDouble > &Fwd, Array< OneD, NekDouble > &Bwd, bool PutFwdInBwdOnBCs)
Definition: ExpList.cpp:4823
void DivideByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Divided by the metric jacobi and quadrature weights.
Definition: ExpList.cpp:1710
virtual const Array< OneD, const std::shared_ptr< ExpList > > & v_GetBndCondExpansions(void)
Definition: ExpList.cpp:4189
void ImposeDirichletConditions(Array< OneD, NekDouble > &outarray)
Impose Dirichlet Boundary Conditions onto Array.
Definition: ExpList.h:1971
Array< OneD, int > m_phys_offset
Offset of elemental data into the array m_phys.
Definition: ExpList.h:1138
virtual void v_SetHomoLen(const NekDouble lhom)
Definition: ExpList.cpp:3678
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:2170
virtual void v_GetNormals(Array< OneD, Array< OneD, NekDouble > > &normals)
Populate normals with the normals of all expansions.
Definition: ExpList.cpp:4413
Array< OneD, SpatialDomains::BoundaryConditionShPtr > & UpdateBndConditions()
Definition: ExpList.h:2275
const DNekScalBlkMatSharedPtr GenBlockMatrix(const GlobalMatrixKey &gkey)
This function assembles the block diagonal matrix of local matrices of the type mtype.
Definition: ExpList.cpp:2280
void Curl(Array< OneD, Array< OneD, NekDouble > > &Vel, Array< OneD, Array< OneD, NekDouble > > &Q)
Definition: ExpList.h:1847
void PhysDeriv(Direction edir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &out_d)
Definition: ExpList.h:1828
void MultiplyByElmtInvMass(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function elementally mulplies the coefficient space of Sin my the elemental inverse of the mass ...
Definition: ExpList.cpp:2170
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:1732
NekDouble GetHomoLen(void)
This function returns the Width of homogeneous direction associated with the homogeneous expansion.
Definition: ExpList.h:608
virtual void v_GetBndElmtExpansion(int i, std::shared_ptr< ExpList > &result, const bool DeclareCoeffPhysArrays)
Definition: ExpList.cpp:5225
virtual void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:5101
bool GetWaveSpace(void) const
This function returns the third direction expansion condition, which can be in wave space (coefficien...
Definition: ExpList.h:1619
virtual std::shared_ptr< ExpList > & v_UpdateBndCondExpansion(int i)
Definition: ExpList.cpp:4197
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:2109
virtual const std::shared_ptr< LocTraceToTraceMap > & v_GetLocTraceToTraceMap(void) const
Definition: ExpList.cpp:6414
virtual void v_IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:1737
void ExtractDataToCoeffs(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata, std::string &field, Array< OneD, NekDouble > &coeffs, std::unordered_map< int, int > zIdToPlane=std::unordered_map< int, int >())
Extract the data in fielddata into the coeffs.
Definition: ExpList.cpp:3961
void SetModifiedBasis(const bool modbasis)
Set Modified Basis for the stability analysis.
ExpansionType m_expType
Exapnsion type.
Definition: ExpList.h:1060
void WriteTecplotHeader(std::ostream &outfile, std::string var="")
Definition: ExpList.h:400
const std::vector< bool > & GetLeftAdjacentFaces(void) const
Definition: ExpList.h:2255
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:1891
void AppendFieldData(LibUtilities::FieldDefinitionsSharedPtr &fielddef, std::vector< NekDouble > &fielddata)
Append the element data listed in elements fielddef->m_ElementIDs onto fielddata.
Definition: ExpList.h:920
std::shared_ptr< LibUtilities::Comm > GetComm() const
Returns the comm object.
Definition: ExpList.h:976
Array< OneD, NekDouble > m_phys
The global expansion evaluated at the quadrature points.
Definition: ExpList.h:1108
std::vector< bool > & GetLeftAdjacentTraces(void)
Definition: ExpList.h:2329
ExpansionType GetExpType(void)
Returns the type of the expansion.
Definition: ExpList.cpp:1656
void AddTraceQuadPhysToOffDiag(const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &field)
Definition: ExpList.h:825
NekDouble Linf(const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &soln=NullNekDouble1DArray)
This function calculates the error of the global spectral/hp element approximation.
Definition: ExpList.cpp:3537
void HomogeneousFwdTrans(const int npts, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool Shuff=true, bool UnShuff=true)
Definition: ExpList.h:1863
void WriteVtkPieceHeader(std::ostream &outfile, int expansion, int istrip=0)
Definition: ExpList.h:418
int GetTotPoints(void) const
Returns the total number of quadrature points m_npoints .
Definition: ExpList.h:1572
void SetExpType(ExpansionType Type)
Returns the type of the expansion.
LibUtilities::BasisSharedPtr GetHomogeneousBasis(void)
Definition: ExpList.h:985
void WriteTecplotZone(std::ostream &outfile, int expansion=-1)
Definition: ExpList.h:404
ExpList(const LibUtilities::SessionReaderSharedPtr &pSession, const Array< OneD, const ExpListSharedPtr > &bndConstraint, const Array< OneD, const SpatialDomains ::BoundaryConditionShPtr > &bndCond, const LocalRegions::ExpansionVector &locexp, const SpatialDomains::MeshGraphSharedPtr &graph, const LibUtilities::CommSharedPtr &comm, const bool DeclareCoeffPhysArrays=true, const std::string variable="DefaultVar", const Collections::ImplementationType ImpType=Collections::eNoImpType)
Generate expansions for the trace space expansions used in DisContField.
virtual void v_SetHomo1DSpecVanVisc(Array< OneD, NekDouble > visc)
Definition: ExpList.h:1487
void GetBCValues(Array< OneD, NekDouble > &BndVals, const Array< OneD, NekDouble > &TotField, int BndID)
Definition: ExpList.h:1901
int GetCoordim(int eid)
This function returns the dimension of the coordinates of the element eid.
Definition: ExpList.h:1926
void GetDiagMatIpwrtBase(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, DNekMatSharedPtr > &mtxPerVar)
Definition: ExpList.cpp:5918
Describes a matrix with ordering defined by a local to global map.
std::vector< Collection > CollectionVector
Definition: Collection.h:110
static BasisSharedPtr NullBasisSharedPtr
Definition: Basis.h:352
static std::vector< unsigned int > NullUnsignedIntVector
std::shared_ptr< Basis > BasisSharedPtr
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< FieldDefinitions > FieldDefinitionsSharedPtr
Definition: FieldIO.h:186
static std::vector< NekDouble > NullNekDoubleVector
std::shared_ptr< Transposition > TranspositionSharedPtr
static Array< OneD, BasisSharedPtr > NullBasisSharedPtr1DArray
Definition: Basis.h:353
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
Definition: Comm.h:57
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:68
std::vector< ExpansionSharedPtr > ExpansionVector
Definition: Expansion.h:70
static ExpList NullExpList
An empty ExpList object.
Definition: ExpList.h:1513
static PeriodicMap NullPeriodicMap
MultiRegions::Direction const DirCartesianMap[]
Definition: ExpList.h:90
std::shared_ptr< BlockMatrixMap > BlockMatrixMapShPtr
A shared pointer to a BlockMatrixMap.
Definition: ExpList.h:96
static VarFactorsMap NullVarFactorsMap
static ExpListSharedPtr NullExpListSharedPtr
Definition: ExpList.h:1514
static LibUtilities::NekManager< GlobalLinSysKey, GlobalLinSys > NullGlobalLinSysManager
Definition: ExpList.h:1518
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:2334
std::map< StdRegions::ConstFactorType, Array< OneD, NekDouble > > VarFactorsMap
std::shared_ptr< AssemblyMap > AssemblyMapSharedPtr
Definition: AssemblyMap.h:52
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:94
static const NekDouble kNekUnsetDouble
GeomMMF
Principle direction for MMF.
std::shared_ptr< BoundaryConditionBase > BoundaryConditionShPtr
Definition: Conditions.h:212
std::map< int, BoundaryConditionMapShPtr > BoundaryConditionCollection
Definition: Conditions.h:219
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:176
std::shared_ptr< PointGeom > PointGeomSharedPtr
Definition: Geometry.h:60
std::map< int, ExpansionInfoShPtr > ExpansionInfoMap
Definition: MeshGraph.h:143
std::map< int, CompositeSharedPtr > CompositeMap
Definition: MeshGraph.h:138
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:408
static VarCoeffMap NullVarCoeffMap
Definition: StdRegions.hpp:353
std::map< StdRegions::VarCoeffType, VarCoeffEntry > VarCoeffMap
Definition: StdRegions.hpp:352
StdRegions::ConstFactorMap factors
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
std::shared_ptr< DNekScalBlkMat > DNekScalBlkMatSharedPtr
Definition: NekTypeDefs.hpp:79
static Array< OneD, NekDouble > NullNekDouble1DArray
double NekDouble
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1191