Nektar++
Public Member Functions | Private Member Functions | Private Attributes | Static Private Attributes | List of all members
Nektar::MultiRegions::GJPStabilisation Class Reference

#include <GJPStabilisation.h>

Public Member Functions

 GJPStabilisation (ExpListSharedPtr field)
 
 ~GJPStabilisation ()
 
void Apply (const Array< OneD, NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const Array< OneD, NekDouble > &pUnorm=NullNekDouble1DArray, const NekDouble scale=1.0) const
 
Array< OneD, Array< OneD, NekDouble > > & GetTraceNormals (void)
 
int GetNumTracePts (void) const
 
bool IsSemiImplicit () const
 

Private Member Functions

void SetUpExpansionInfoMapForGJP (SpatialDomains::MeshGraphSharedPtr graph, std::string variable)
 
void MultiplyByStdDerivBaseOnTraceMat (int i, Array< OneD, NekDouble > &in, Array< OneD, NekDouble > &out) const
 

Private Attributes

int m_coordDim
 
int m_traceDim
 
bool m_useGJPSemiImplicit
 
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
 
MultiRegions::ExpListSharedPtr m_dgfield
 DG expansion for projection evalaution along trace. More...
 
MultiRegions::LocTraceToTraceMapSharedPtr m_locTraceToTraceMap
 LocaTraceToTraceMap. More...
 
MultiRegions::ExpListSharedPtr m_locElmtTrace
 Local Elemental trace expansions. More...
 
Array< OneD, Array< OneD, NekDouble > > m_scalTrace
 Scale factor for phys values along trace involving the lcoal normals and tangent geometric factors n. More...
 
std::vector< std::pair< int, Array< OneD, DNekMatSharedPtr > > > m_StdDBaseOnTraceMat
 

Static Private Attributes

static std::string GJPStabilisationLookupIds []
 

Detailed Description

Definition at line 48 of file GJPStabilisation.h.

Constructor & Destructor Documentation

◆ GJPStabilisation()

Nektar::MultiRegions::GJPStabilisation::GJPStabilisation ( ExpListSharedPtr  field)

Definition at line 47 of file GJPStabilisation.cpp.

48{
49 LibUtilities::SessionReaderSharedPtr session = pField->GetSession();
50
51 session->MatchSolverInfo("GJPStabilisation", "SemiImplicit",
53
54 // Call GetTrace on the initialising field will set up
55 // DG. Store a copoy so that if we make a soft copy of
56 // this class we can re-used this field for operators.
57 pField->GetTrace();
58 m_dgfield = pField;
59
60 m_coordDim = m_dgfield->GetCoordim(0);
61 m_traceDim = m_dgfield->GetShapeDimension() - 1;
62
63 // set up trace normals but would be better if could use
64 // equation system definition
65 m_traceNormals = Array<OneD, Array<OneD, NekDouble>>(m_coordDim);
66 for (int i = 0; i < m_coordDim; ++i)
67 {
69 Array<OneD, NekDouble>(m_dgfield->GetTrace()->GetNpoints());
70 }
71 m_dgfield->GetTrace()->GetNormals(m_traceNormals);
72
73 SetUpExpansionInfoMapForGJP(pField->GetGraph(), session->GetVariable(0));
74
76
78 session, pField->GetGraph(), "GJP", true, false,
79 Collections::eNoImpType, session->GetVariable(0));
80 dgfield->GetLocTraceToTraceMap(m_locTraceToTraceMap);
81
83 session, *(dgfield->GetExp()), dgfield->GetGraph(), true, "GJP");
84
85 m_scalTrace = Array<OneD, Array<OneD, NekDouble>>(m_traceDim + 1);
86
87 const std::shared_ptr<LocalRegions::ExpansionVector> exp =
88 dgfield->GetExp();
89
90 Array<OneD, Array<OneD, NekDouble>> dfactors[3];
91 Array<OneD, Array<OneD, NekDouble>> LocTrace(m_traceDim + 1);
92 Array<OneD, NekDouble> e_tmp;
93
94 for (int i = 0; i < m_traceDim + 1; ++i)
95 {
96 LocTrace[i] =
97 Array<OneD, NekDouble>(m_locTraceToTraceMap->GetNLocTracePts());
98 }
99
100 int cnt = 0;
101 int offset_phys = 0;
102 Array<OneD, Array<OneD, Array<OneD, NekDouble>>> dbasis;
103 Array<OneD, Array<OneD, Array<OneD, unsigned int>>> traceToCoeffMap;
104
105 Array<OneD, unsigned int> map, map1;
106 Array<OneD, int> sign, sign1;
107 NekDouble h, p;
108
109 for (int e = 0; e < m_dgfield->GetExpSize(); ++e)
110 {
111 LocalRegions::ExpansionSharedPtr elmt = (*exp)[e];
112
113 elmt->NormalTraceDerivFactors(dfactors[0], dfactors[1], dfactors[2]);
114
115 for (int n = 0; n < elmt->GetNtraces(); ++n, ++cnt)
116 {
117 NekDouble jumpScal;
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));
123
124 if (p == 1)
125 {
126 jumpScal = 0.02 * h * h;
127 }
128 else
129 {
130 jumpScal = 0.8 * pow(p + 1, -4.0) * h * h;
131 }
132
133 int nptrace = elmt->GetTraceNumPoints(n);
134 elmt->GetTraceCoeffMap(n, map);
135
136 for (int i = 0; i < m_traceDim + 1; ++i)
137 {
138 Vmath::Smul(nptrace, jumpScal, dfactors[i][n], 1,
139 e_tmp = LocTrace[i] + offset_phys, 1);
140 }
141
142 offset_phys += nptrace;
143 }
144 }
145
146 for (int i = 0; i < m_traceDim + 1; ++i)
147 {
148 m_scalTrace[i] = LocTrace[i];
149
150 if (m_traceDim > 0)
151 {
152 // multiply by Jacobian and quadrature points.
153 m_locElmtTrace->MultiplyByQuadratureMetric(m_scalTrace[i],
154 m_scalTrace[i]);
155 }
156 }
157
158 // Assemble list of Matrix Product
159 Array<OneD, DNekMatSharedPtr> TraceMat;
160
161 int nelmt = 1;
162 Array<OneD, const LibUtilities::BasisSharedPtr> base_sav =
163 dgfield->GetExp(0)->GetBase();
164
165 dgfield->GetExp(0)->StdDerivBaseOnTraceMat(TraceMat);
166
167 for (int n = 1; n < dgfield->GetExpSize(); ++n)
168 {
169 const Array<OneD, const LibUtilities::BasisSharedPtr> &base =
170 dgfield->GetExp(n)->GetBase();
171
172 int i;
173 for (i = 0; i < base.size(); ++i)
174 {
175 if (base[i] != base_sav[i])
176 {
177 break;
178 }
179 }
180
181 if (i == base.size())
182 {
183 nelmt++;
184 }
185 else
186 {
187 // save previous block of data.
188 m_StdDBaseOnTraceMat.push_back(
189 std::pair<int, Array<OneD, DNekMatSharedPtr>>(nelmt, TraceMat));
190
191 // start new block
192 dgfield->GetExp(n)->StdDerivBaseOnTraceMat(TraceMat);
193 nelmt = 1;
194 base_sav = dgfield->GetExp(n)->GetBase();
195 }
196 }
197 // save latest block of data.
198 m_StdDBaseOnTraceMat.push_back(
199 std::pair<int, Array<OneD, DNekMatSharedPtr>>(nelmt, TraceMat));
200}
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:49
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
std::vector< std::pair< int, Array< OneD, DNekMatSharedPtr > > > m_StdDBaseOnTraceMat
void SetUpExpansionInfoMapForGJP(SpatialDomains::MeshGraphSharedPtr graph, std::string variable)
MultiRegions::ExpListSharedPtr m_locElmtTrace
Local Elemental trace expansions.
MultiRegions::LocTraceToTraceMapSharedPtr m_locTraceToTraceMap
LocaTraceToTraceMap.
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.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:68
std::shared_ptr< DisContField > DisContFieldSharedPtr
Definition: DisContField.h:341
double NekDouble
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.cpp:245

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), ASSERTL0, Nektar::Collections::eNoImpType, m_coordDim, m_dgfield, m_locElmtTrace, m_locTraceToTraceMap, m_scalTrace, m_StdDBaseOnTraceMat, m_traceDim, m_traceNormals, m_useGJPSemiImplicit, CellMLToNektar.cellml_metadata::p, SetUpExpansionInfoMapForGJP(), sign, and Vmath::Smul().

◆ ~GJPStabilisation()

Nektar::MultiRegions::GJPStabilisation::~GJPStabilisation ( )
inline

Definition at line 53 of file GJPStabilisation.h.

53{};

Member Function Documentation

◆ Apply()

void Nektar::MultiRegions::GJPStabilisation::Apply ( const Array< OneD, NekDouble > &  inarray,
Array< OneD, NekDouble > &  outarray,
const Array< OneD, NekDouble > &  pUnorm = NullNekDouble1DArray,
const NekDouble  scale = 1.0 
) const

Definition at line 202 of file GJPStabilisation.cpp.

206{
207 int ncoeffs = m_dgfield->GetNcoeffs();
208 int nphys = m_dgfield->GetNpoints();
209 int nTracePts = m_dgfield->GetTrace()->GetTotPoints();
210 int nLocETrace = m_locElmtTrace->GetTotPoints();
211
212 ASSERTL1(nLocETrace == m_scalTrace[0].size(), "expect these to be similar");
213 ASSERTL1(m_locElmtTrace->GetNcoeffs() <= nphys,
214 "storage assumptions presume "
215 "that nLocETraceCoeffs < nphys");
216
217 Array<OneD, Array<OneD, NekDouble>> deriv(3, NullNekDouble1DArray);
218 for (int i = 0; i < m_coordDim; ++i)
219 {
220 deriv[i] = Array<OneD, NekDouble>(nphys);
221 }
222
223 int nmax = std::max(ncoeffs, nphys);
224 Array<OneD, NekDouble> FilterCoeffs(nmax, 0.0);
225 Array<OneD, NekDouble> GradJumpOnTrace(nTracePts, 0.0);
226 Array<OneD, NekDouble> Fwd(nTracePts), Bwd(nTracePts);
227
228 Array<OneD, NekDouble> wsp(nLocETrace), tmp;
229 Array<OneD, NekDouble> LocElmtTracePhys = m_locElmtTrace->UpdatePhys();
230 Array<OneD, NekDouble> LocElmtTraceCoeffs = m_locElmtTrace->UpdateCoeffs();
231
232 ASSERTL1(LocElmtTracePhys.size() <= nLocETrace,
233 "expect this vector to be at least of size nLocETrace");
234
235 Array<OneD, NekDouble> unorm;
236 if (pUnorm == NullNekDouble1DArray)
237 {
238 unorm = Array<OneD, NekDouble>(nTracePts, 1.0);
239 }
240 else
241 {
242 unorm = pUnorm;
243 }
244
245 Array<OneD, NekDouble> GradJumpOnTraceBwd;
247 {
248 GradJumpOnTraceBwd = Array<OneD, NekDouble>(nTracePts);
249 }
250
252 {
253 Vmath::Zero(nTracePts, GradJumpOnTraceBwd, 1);
254 }
255
256 // calculate derivative
257 m_dgfield->PhysDeriv(inarray, deriv[0], deriv[1], deriv[2]);
258
259 // Evaluate the normal derivative jump on the trace
260 for (int n = 0; n < m_coordDim; ++n)
261 {
262 m_dgfield->GetFwdBwdTracePhys(deriv[n], Fwd, Bwd, true, true);
263
265 {
266 // want to put Fwd vals on bwd trace and vice versa
267 Vmath::Vvtvp(nTracePts, Bwd, 1, m_traceNormals[n], 1,
268 GradJumpOnTrace, 1, GradJumpOnTrace, 1);
269 Vmath::Vvtvp(nTracePts, Fwd, 1, m_traceNormals[n], 1,
270 GradJumpOnTraceBwd, 1, GradJumpOnTraceBwd, 1);
271 }
272 else
273 {
274 // Multiply by normal and add to trace evaluation
275 Vmath::Vsub(nTracePts, Fwd, 1, Bwd, 1, Fwd, 1);
276 Vmath::Vvtvp(nTracePts, Fwd, 1, m_traceNormals[n], 1,
277 GradJumpOnTrace, 1, GradJumpOnTrace, 1);
278 }
279 }
280
282 {
283 // Need to negate Bwd case when using Fwd normal
284 Vmath::Neg(nTracePts, GradJumpOnTrace, 1);
285 }
286
287 Vmath::Vmul(nTracePts, unorm, 1, GradJumpOnTrace, 1, GradJumpOnTrace, 1);
288
289 // Interpolate GradJumpOnTrace to Local elemental traces.
290 m_locTraceToTraceMap->InterpTraceToLocTrace(0, GradJumpOnTrace, wsp);
291 m_locTraceToTraceMap->UnshuffleLocTraces(0, wsp, LocElmtTracePhys);
292
294 {
295 // Vmath::Neg(nTracePts,GradJumpOnTraceBwd,1);
296 Vmath::Vmul(nTracePts, unorm, 1, GradJumpOnTraceBwd, 1,
297 GradJumpOnTraceBwd, 1);
298 m_locTraceToTraceMap->InterpTraceToLocTrace(1, GradJumpOnTraceBwd, wsp);
299 m_locTraceToTraceMap->UnshuffleLocTraces(1, wsp, LocElmtTracePhys);
300 }
301 else
302 {
303 m_locTraceToTraceMap->InterpTraceToLocTrace(1, GradJumpOnTrace, wsp);
304 m_locTraceToTraceMap->UnshuffleLocTraces(1, wsp, LocElmtTracePhys);
305 }
306
307 // Scale jump on trace
308 Vmath::Vmul(nLocETrace, m_scalTrace[0], 1, LocElmtTracePhys, 1, wsp, 1);
309 MultiplyByStdDerivBaseOnTraceMat(0, wsp, FilterCoeffs);
310
311 for (int i = 0; i < m_traceDim; ++i)
312 {
313 // Scale jump on trace
314 Vmath::Vmul(nLocETrace, m_scalTrace[i + 1], 1, LocElmtTracePhys, 1, wsp,
315 1);
316 MultiplyByStdDerivBaseOnTraceMat(i + 1, wsp, deriv[0]);
317 Vmath::Vadd(ncoeffs, deriv[0], 1, FilterCoeffs, 1, FilterCoeffs, 1);
318 }
319
320 Vmath::Svtvp(ncoeffs, scale, FilterCoeffs, 1, outarray, 1, outarray, 1);
321}
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
void MultiplyByStdDerivBaseOnTraceMat(int i, Array< OneD, NekDouble > &in, Array< OneD, NekDouble > &out) const
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.
Definition: Vmath.cpp:207
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
Definition: Vmath.cpp:617
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:513
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:569
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.
Definition: Vmath.cpp:354
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:487
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.
Definition: Vmath.cpp:414

References ASSERTL1, m_coordDim, m_dgfield, m_locElmtTrace, m_locTraceToTraceMap, m_scalTrace, m_traceDim, m_traceNormals, m_useGJPSemiImplicit, MultiplyByStdDerivBaseOnTraceMat(), Vmath::Neg(), Nektar::NullNekDouble1DArray, Vmath::Svtvp(), Vmath::Vadd(), Vmath::Vmul(), Vmath::Vsub(), Vmath::Vvtvp(), and Vmath::Zero().

◆ GetNumTracePts()

int Nektar::MultiRegions::GJPStabilisation::GetNumTracePts ( void  ) const
inline

Definition at line 65 of file GJPStabilisation.h.

66 {
67 return m_dgfield->GetTrace()->GetTotPoints();
68 }

References m_dgfield.

◆ GetTraceNormals()

Array< OneD, Array< OneD, NekDouble > > & Nektar::MultiRegions::GJPStabilisation::GetTraceNormals ( void  )
inline

Definition at line 60 of file GJPStabilisation.h.

61 {
62 return m_traceNormals;
63 }

References m_traceNormals.

◆ IsSemiImplicit()

bool Nektar::MultiRegions::GJPStabilisation::IsSemiImplicit ( ) const
inline

Definition at line 70 of file GJPStabilisation.h.

71 {
73 }

References m_useGJPSemiImplicit.

◆ MultiplyByStdDerivBaseOnTraceMat()

void Nektar::MultiRegions::GJPStabilisation::MultiplyByStdDerivBaseOnTraceMat ( int  i,
Array< OneD, NekDouble > &  in,
Array< OneD, NekDouble > &  out 
) const
private

Definition at line 382 of file GJPStabilisation.cpp.

384{
385 // Should probably be vectorised
386
387 int cnt = 0;
388 int cnt1 = 0;
389 for (auto &it : m_StdDBaseOnTraceMat)
390 {
391 int rows = it.second[i]->GetRows();
392 int cols = it.second[i]->GetColumns();
393
394 Blas::Dgemm('N', 'N', rows, it.first, cols, 1.0,
395 &(it.second[i]->GetPtr())[0], rows, &in[0] + cnt, cols, 0.0,
396 &out[0] + cnt1, rows);
397
398 cnt += cols * it.first;
399 cnt1 += rows * it.first;
400 }
401}
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...
Definition: Blas.hpp:385

References Blas::Dgemm(), and m_StdDBaseOnTraceMat.

Referenced by Apply().

◆ SetUpExpansionInfoMapForGJP()

void Nektar::MultiRegions::GJPStabilisation::SetUpExpansionInfoMapForGJP ( SpatialDomains::MeshGraphSharedPtr  graph,
std::string  variable 
)
private

Definition at line 323 of file GJPStabilisation.cpp.

325{
326
327 // check to see if already deifned and if so return
328 if (graph->ExpansionInfoDefined("GJP"))
329 {
330 return;
331 }
332
334 graph->GetExpansionInfo(variable);
335
338
339 // loop over epxansion info
340 for (auto expIt = expInfo.begin(); expIt != expInfo.end(); ++expIt)
341 {
342 std::vector<LibUtilities::BasisKey> BKeyVector;
343
344 for (int i = 0; i < expIt->second->m_basisKeyVector.size(); ++i)
345 {
346 LibUtilities::BasisKey bkeyold = expIt->second->m_basisKeyVector[i];
347
348 // Reset radauM alpha non-zero cases to radauM
349 // Legendre at one order higher
350
351 switch (bkeyold.GetPointsType())
352 {
353 case LibUtilities::eGaussRadauMAlpha1Beta0:
354 case LibUtilities::eGaussRadauMAlpha2Beta0:
355 {
356 int npts = bkeyold.GetNumPoints();
357
358 // const LibUtilities::PointsKey pkey(npts+1,
359 // LibUtilities::eGaussRadauMLegendre); trying
360 // npts to be consistent for tri faces
361 const LibUtilities::PointsKey pkey(
363 LibUtilities::BasisKey bkeynew(bkeyold.GetBasisType(),
364 bkeyold.GetNumModes(), pkey);
365 BKeyVector.push_back(bkeynew);
366 }
367 break;
368 default:
369 BKeyVector.push_back(bkeyold);
370 break;
371 }
372 }
373
374 (*newInfo)[expIt->first] =
376 expIt->second->m_geomShPtr, BKeyVector);
377 }
378
379 graph->SetExpansionInfo("GJP", newInfo);
380}
@ eGaussRadauMLegendre
1D Gauss-Radau-Legendre quadrature points, pinned at x=-1
Definition: PointsType.h:49
std::shared_ptr< ExpansionInfoMap > ExpansionInfoMapShPtr
Definition: MeshGraph.h:145
std::map< int, ExpansionInfoShPtr > ExpansionInfoMap
Definition: MeshGraph.h:143

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::LibUtilities::eGaussRadauMLegendre, Nektar::LibUtilities::BasisKey::GetBasisType(), Nektar::LibUtilities::BasisKey::GetNumModes(), Nektar::LibUtilities::BasisKey::GetNumPoints(), and Nektar::LibUtilities::BasisKey::GetPointsType().

Referenced by GJPStabilisation().

Member Data Documentation

◆ GJPStabilisationLookupIds

std::string Nektar::MultiRegions::GJPStabilisation::GJPStabilisationLookupIds
staticprivate
Initial value:
= {
"GJPStabilisation", "Explicit", eExplicitGJPStabilisation),
"GJPStabilisation", "SemiImplicit", eSemiImplicitGJPStabilisation),
}
static std::string RegisterEnumValue(std::string pEnum, std::string pString, int pEnumValue)
Registers an enumeration value.

Definition at line 79 of file GJPStabilisation.h.

◆ m_coordDim

int Nektar::MultiRegions::GJPStabilisation::m_coordDim
private

Definition at line 76 of file GJPStabilisation.h.

Referenced by Apply(), and GJPStabilisation().

◆ m_dgfield

MultiRegions::ExpListSharedPtr Nektar::MultiRegions::GJPStabilisation::m_dgfield
private

DG expansion for projection evalaution along trace.

Definition at line 85 of file GJPStabilisation.h.

Referenced by Apply(), GetNumTracePts(), and GJPStabilisation().

◆ m_locElmtTrace

MultiRegions::ExpListSharedPtr Nektar::MultiRegions::GJPStabilisation::m_locElmtTrace
private

Local Elemental trace expansions.

Definition at line 89 of file GJPStabilisation.h.

Referenced by Apply(), and GJPStabilisation().

◆ m_locTraceToTraceMap

MultiRegions::LocTraceToTraceMapSharedPtr Nektar::MultiRegions::GJPStabilisation::m_locTraceToTraceMap
private

LocaTraceToTraceMap.

Definition at line 87 of file GJPStabilisation.h.

Referenced by Apply(), and GJPStabilisation().

◆ m_scalTrace

Array<OneD, Array<OneD, NekDouble> > Nektar::MultiRegions::GJPStabilisation::m_scalTrace
private

Scale factor for phys values along trace involving the lcoal normals and tangent geometric factors n.

Definition at line 93 of file GJPStabilisation.h.

Referenced by Apply(), and GJPStabilisation().

◆ m_StdDBaseOnTraceMat

std::vector<std::pair<int, Array<OneD, DNekMatSharedPtr> > > Nektar::MultiRegions::GJPStabilisation::m_StdDBaseOnTraceMat
private

Definition at line 96 of file GJPStabilisation.h.

Referenced by GJPStabilisation(), and MultiplyByStdDerivBaseOnTraceMat().

◆ m_traceDim

int Nektar::MultiRegions::GJPStabilisation::m_traceDim
private

Definition at line 77 of file GJPStabilisation.h.

Referenced by Apply(), and GJPStabilisation().

◆ m_traceNormals

Array<OneD, Array<OneD, NekDouble> > Nektar::MultiRegions::GJPStabilisation::m_traceNormals
private

Definition at line 82 of file GJPStabilisation.h.

Referenced by Apply(), GetTraceNormals(), and GJPStabilisation().

◆ m_useGJPSemiImplicit

bool Nektar::MultiRegions::GJPStabilisation::m_useGJPSemiImplicit
private

Definition at line 78 of file GJPStabilisation.h.

Referenced by Apply(), GJPStabilisation(), and IsSemiImplicit().