51 session->MatchSolverInfo(
"GJPStabilisation",
"SemiImplicit",
78 session, pField->GetGraph(),
"GJP",
true,
false,
83 session, *(dgfield->GetExp()), dgfield->GetGraph(),
true,
"GJP");
87 const std::shared_ptr<LocalRegions::ExpansionVector> exp =
109 for (
int e = 0; e <
m_dgfield->GetExpSize(); ++e)
113 elmt->NormalTraceDerivFactors(dfactors[0], dfactors[1], dfactors[2]);
115 for (
int n = 0; n < elmt->GetNtraces(); ++n, ++cnt)
118 elmt->TraceNormLen(n, h,
p);
119 ASSERTL0(boost::math::isnan(h) ==
false,
120 "h has a nan value when e = " +
121 boost::lexical_cast<std::string>(e) +
122 " n =" + boost::lexical_cast<std::string>(n));
126 jumpScal = 0.02 * h * h;
130 jumpScal = 0.8 * pow(
p + 1, -4.0) * h * h;
133 int nptrace = elmt->GetTraceNumPoints(n);
134 elmt->GetTraceCoeffMap(n, map);
139 e_tmp = LocTrace[i] + offset_phys, 1);
142 offset_phys += nptrace;
163 dgfield->GetExp(0)->GetBase();
165 dgfield->GetExp(0)->StdDerivBaseOnTraceMat(TraceMat);
167 for (
int n = 1; n < dgfield->GetExpSize(); ++n)
170 dgfield->GetExp(n)->GetBase();
173 for (i = 0; i < base.size(); ++i)
175 if (base[i] != base_sav[i])
181 if (i == base.size())
192 dgfield->GetExp(n)->StdDerivBaseOnTraceMat(TraceMat);
194 base_sav = dgfield->GetExp(n)->GetBase();
209 int nTracePts =
m_dgfield->GetTrace()->GetTotPoints();
214 "storage assumptions presume "
215 "that nLocETraceCoeffs < nphys");
223 int nmax = std::max(ncoeffs, nphys);
232 ASSERTL1(LocElmtTracePhys.size() <= nLocETrace,
233 "expect this vector to be at least of size nLocETrace");
257 m_dgfield->PhysDeriv(inarray, deriv[0], deriv[1], deriv[2]);
262 m_dgfield->GetFwdBwdTracePhys(deriv[n], Fwd, Bwd,
true,
true);
268 GradJumpOnTrace, 1, GradJumpOnTrace, 1);
270 GradJumpOnTraceBwd, 1, GradJumpOnTraceBwd, 1);
277 GradJumpOnTrace, 1, GradJumpOnTrace, 1);
287 Vmath::Vmul(nTracePts, unorm, 1, GradJumpOnTrace, 1, GradJumpOnTrace, 1);
296 Vmath::Vmul(nTracePts, unorm, 1, GradJumpOnTraceBwd, 1,
297 GradJumpOnTraceBwd, 1);
317 Vmath::Vadd(ncoeffs, deriv[0], 1, FilterCoeffs, 1, FilterCoeffs, 1);
320 Vmath::Svtvp(ncoeffs, scale, FilterCoeffs, 1, outarray, 1, outarray, 1);
328 if (graph->ExpansionInfoDefined(
"GJP"))
334 graph->GetExpansionInfo(variable);
340 for (
auto expIt = expInfo.begin(); expIt != expInfo.end(); ++expIt)
342 std::vector<LibUtilities::BasisKey> BKeyVector;
344 for (
int i = 0; i < expIt->second->m_basisKeyVector.size(); ++i)
353 case LibUtilities::eGaussRadauMAlpha1Beta0:
354 case LibUtilities::eGaussRadauMAlpha2Beta0:
365 BKeyVector.push_back(bkeynew);
369 BKeyVector.push_back(bkeyold);
374 (*newInfo)[expIt->first] =
376 expIt->second->m_geomShPtr, BKeyVector);
379 graph->SetExpansionInfo(
"GJP", newInfo);
391 int rows = it.second[i]->GetRows();
392 int cols = it.second[i]->GetColumns();
395 &(it.second[i]->GetPtr())[0], rows, &in[0] + cnt, cols, 0.0,
396 &out[0] + cnt1, rows);
398 cnt += cols * it.first;
399 cnt1 += rows * it.first;
#define ASSERTL0(condition, msg)
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
#define sign(a, b)
return the sign(b)*a
Describes the specification for a Basis.
int GetNumPoints() const
Return points order at which basis is defined.
BasisType GetBasisType() const
Return type of expansion basis.
int GetNumModes() const
Returns the order of the basis.
PointsType GetPointsType() const
Return type of quadrature.
Defines a specification for a set of points.
static std::string RegisterEnumValue(std::string pEnum, std::string pString, int pEnumValue)
Registers an enumeration value.
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
void Apply(const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const Array< OneD, NekDouble > &pUnorm=NullNekDouble1DArray, const NekDouble scale=1.0) const
std::vector< std::pair< int, Array< OneD, DNekMatSharedPtr > > > m_StdDBaseOnTraceMat
void SetUpExpansionInfoMapForGJP(SpatialDomains::MeshGraphSharedPtr graph, std::string variable)
void MultiplyByStdDerivBaseOnTraceMat(int i, Array< OneD, NekDouble > &in, Array< OneD, NekDouble > &out) const
MultiRegions::ExpListSharedPtr m_locElmtTrace
Local Elemental trace expansions.
static std::string GJPStabilisationLookupIds[]
GJPStabilisation(ExpListSharedPtr field)
MultiRegions::LocTraceToTraceMapSharedPtr m_locTraceToTraceMap
LocaTraceToTraceMap.
bool m_useGJPSemiImplicit
Array< OneD, Array< OneD, NekDouble > > m_scalTrace
Scale factor for phys values along trace involving the lcoal normals and tangent geometric factors n.
MultiRegions::ExpListSharedPtr m_dgfield
DG expansion for projection evalaution along trace.
static void Dgemm(const char &transa, const char &transb, const int &m, const int &n, const int &k, const double &alpha, const double *a, const int &lda, const double *b, const int &ldb, const double &beta, double *c, const int &ldc)
BLAS level 3: Matrix-matrix multiply C = A x B where op(A)[m x k], op(B)[k x n], C[m x n] DGEMM perfo...
std::shared_ptr< SessionReader > SessionReaderSharedPtr
@ eGaussRadauMLegendre
1D Gauss-Radau-Legendre quadrature points, pinned at x=-1
std::shared_ptr< Expansion > ExpansionSharedPtr
std::shared_ptr< DisContField > DisContFieldSharedPtr
@ eSemiImplicitGJPStabilisation
@ eExplicitGJPStabilisation
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
std::shared_ptr< ExpansionInfoMap > ExpansionInfoMapShPtr
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
std::map< int, ExpansionInfoShPtr > ExpansionInfoMap
The above copyright notice and this permission notice shall be included.
static Array< OneD, NekDouble > NullNekDouble1DArray
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
void Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
svtvp (scalar times vector plus vector): z = alpha*x + y
void Neg(int n, T *x, const int incx)
Negate x = -x.
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
void Zero(int n, T *x, const int incx)
Zero vector.
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.