Nektar++
Public Member Functions | Private Member Functions | 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
 

Detailed Description

Definition at line 48 of file GJPStabilisation.h.

Constructor & Destructor Documentation

◆ GJPStabilisation()

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

Definition at line 41 of file GJPStabilisation.cpp.

42 {
43  LibUtilities::SessionReaderSharedPtr session = pField->GetSession();
44 
45  session->MatchSolverInfo("GJPStabilisation", "SemiImplicit",
46  m_useGJPSemiImplicit, false);
47 
48  // Call GetTrace on the initialising field will set up
49  // DG. Store a copoy so that if we make a soft copy of
50  // this class we can re-used this field for operators.
51  pField->GetTrace();
52  m_dgfield = pField;
53 
54  m_coordDim = m_dgfield->GetCoordim(0);
55  m_traceDim = m_dgfield->GetShapeDimension() - 1;
56 
57  // set up trace normals but would be better if could use
58  // equation system definition
59  m_traceNormals = Array<OneD, Array<OneD, NekDouble>>(m_coordDim);
60  for (int i = 0; i < m_coordDim; ++i)
61  {
62  m_traceNormals[i] =
63  Array<OneD, NekDouble>(m_dgfield->GetTrace()->GetNpoints());
64  }
65  m_dgfield->GetTrace()->GetNormals(m_traceNormals);
66 
67  SetUpExpansionInfoMapForGJP(pField->GetGraph(), session->GetVariable(0));
68 
70 
72  session, pField->GetGraph(), "GJP", true, false,
73  Collections::eNoImpType, session->GetVariable(0));
74  dgfield->GetLocTraceToTraceMap(m_locTraceToTraceMap);
75 
77  session, *(dgfield->GetExp()), dgfield->GetGraph(), true, "GJP");
78 
79  m_scalTrace = Array<OneD, Array<OneD, NekDouble>>(m_traceDim + 1);
80 
81  const std::shared_ptr<LocalRegions::ExpansionVector> exp =
82  dgfield->GetExp();
83 
84  Array<OneD, Array<OneD, NekDouble>> dfactors[3];
85  Array<OneD, Array<OneD, NekDouble>> LocTrace(m_traceDim + 1);
86  Array<OneD, NekDouble> e_tmp;
87 
88  for (int i = 0; i < m_traceDim + 1; ++i)
89  {
90  LocTrace[i] =
91  Array<OneD, NekDouble>(m_locTraceToTraceMap->GetNLocTracePts());
92  }
93 
94  int cnt = 0;
95  int offset_phys = 0;
96  Array<OneD, Array<OneD, Array<OneD, NekDouble>>> dbasis;
97  Array<OneD, Array<OneD, Array<OneD, unsigned int>>> traceToCoeffMap;
98 
99  Array<OneD, unsigned int> map, map1;
100  Array<OneD, int> sign, sign1;
101  NekDouble h, p;
102 
103  for (int e = 0; e < m_dgfield->GetExpSize(); ++e)
104  {
105  LocalRegions::ExpansionSharedPtr elmt = (*exp)[e];
106 
107  elmt->NormalTraceDerivFactors(dfactors[0], dfactors[1], dfactors[2]);
108 
109  for (int n = 0; n < elmt->GetNtraces(); ++n, ++cnt)
110  {
111  NekDouble jumpScal;
112  elmt->TraceNormLen(n, h, p);
113  ASSERTL0(boost::math::isnan(h) == false,
114  "h has a nan value when e = " +
115  boost::lexical_cast<std::string>(e) +
116  " n =" + boost::lexical_cast<std::string>(n));
117 
118  if (p == 1)
119  {
120  jumpScal = 0.02 * h * h;
121  }
122  else
123  {
124  jumpScal = 0.8 * pow(p + 1, -4.0) * h * h;
125  }
126 
127  int nptrace = elmt->GetTraceNumPoints(n);
128  elmt->GetTraceCoeffMap(n, map);
129 
130  for (int i = 0; i < m_traceDim + 1; ++i)
131  {
132  Vmath::Smul(nptrace, jumpScal, dfactors[i][n], 1,
133  e_tmp = LocTrace[i] + offset_phys, 1);
134  }
135 
136  offset_phys += nptrace;
137  }
138  }
139 
140  for (int i = 0; i < m_traceDim + 1; ++i)
141  {
142  m_scalTrace[i] = LocTrace[i];
143 
144  if (m_traceDim > 0)
145  {
146  // multiply by Jacobian and quadrature points.
147  m_locElmtTrace->MultiplyByQuadratureMetric(m_scalTrace[i],
148  m_scalTrace[i]);
149  }
150  }
151 
152  // Assemble list of Matrix Product
153  Array<OneD, DNekMatSharedPtr> TraceMat;
154 
155  int nelmt = 1;
156  Array<OneD, const LibUtilities::BasisSharedPtr> base_sav =
157  dgfield->GetExp(0)->GetBase();
158 
159  dgfield->GetExp(0)->StdDerivBaseOnTraceMat(TraceMat);
160 
161  for (int n = 1; n < dgfield->GetExpSize(); ++n)
162  {
163  const Array<OneD, const LibUtilities::BasisSharedPtr> &base =
164  dgfield->GetExp(n)->GetBase();
165 
166  int i;
167  for (i = 0; i < base.size(); ++i)
168  {
169  if (base[i] != base_sav[i])
170  {
171  break;
172  }
173  }
174 
175  if (i == base.size())
176  {
177  nelmt++;
178  }
179  else
180  {
181  // save previous block of data.
182  m_StdDBaseOnTraceMat.push_back(
183  std::pair<int, Array<OneD, DNekMatSharedPtr>>(nelmt, TraceMat));
184 
185  // start new block
186  dgfield->GetExp(n)->StdDerivBaseOnTraceMat(TraceMat);
187  nelmt = 1;
188  base_sav = dgfield->GetExp(n)->GetBase();
189  }
190  }
191  // save latest block of data.
192  m_StdDBaseOnTraceMat.push_back(
193  std::pair<int, Array<OneD, DNekMatSharedPtr>>(nelmt, TraceMat));
194 }
#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:248

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 196 of file GJPStabilisation.cpp.

200 {
201  int ncoeffs = m_dgfield->GetNcoeffs();
202  int nphys = m_dgfield->GetNpoints();
203  int nTracePts = m_dgfield->GetTrace()->GetTotPoints();
204  int nLocETrace = m_locElmtTrace->GetTotPoints();
205 
206  ASSERTL1(nLocETrace == m_scalTrace[0].size(), "expect these to be similar");
207  ASSERTL1(m_locElmtTrace->GetNcoeffs() <= nphys,
208  "storage assumptions presume "
209  "that nLocETraceCoeffs < nphys");
210 
211  Array<OneD, Array<OneD, NekDouble>> deriv(3, NullNekDouble1DArray);
212  for (int i = 0; i < m_coordDim; ++i)
213  {
214  deriv[i] = Array<OneD, NekDouble>(nphys);
215  }
216 
217  int nmax = std::max(ncoeffs, nphys);
218  Array<OneD, NekDouble> FilterCoeffs(nmax, 0.0);
219  Array<OneD, NekDouble> GradJumpOnTrace(nTracePts, 0.0);
220  Array<OneD, NekDouble> Fwd(nTracePts), Bwd(nTracePts);
221 
222  Array<OneD, NekDouble> wsp(nLocETrace), tmp;
223  Array<OneD, NekDouble> LocElmtTracePhys = m_locElmtTrace->UpdatePhys();
224  Array<OneD, NekDouble> LocElmtTraceCoeffs = m_locElmtTrace->UpdateCoeffs();
225 
226  ASSERTL1(LocElmtTracePhys.size() <= nLocETrace,
227  "expect this vector to be at least of size nLocETrace");
228 
229  Array<OneD, NekDouble> unorm;
230  if (pUnorm == NullNekDouble1DArray)
231  {
232  unorm = Array<OneD, NekDouble>(nTracePts, 1.0);
233  }
234  else
235  {
236  unorm = pUnorm;
237  }
238 
239  Array<OneD, NekDouble> GradJumpOnTraceBwd;
241  {
242  GradJumpOnTraceBwd = Array<OneD, NekDouble>(nTracePts);
243  }
244 
246  {
247  Vmath::Zero(nTracePts, GradJumpOnTraceBwd, 1);
248  }
249 
250  // calculate derivative
251  m_dgfield->PhysDeriv(inarray, deriv[0], deriv[1], deriv[2]);
252 
253  // Evaluate the normal derivative jump on the trace
254  for (int n = 0; n < m_coordDim; ++n)
255  {
256  m_dgfield->GetFwdBwdTracePhys(deriv[n], Fwd, Bwd, true, true);
257 
259  {
260  // want to put Fwd vals on bwd trace and vice versa
261  Vmath::Vvtvp(nTracePts, Bwd, 1, m_traceNormals[n], 1,
262  GradJumpOnTrace, 1, GradJumpOnTrace, 1);
263  Vmath::Vvtvp(nTracePts, Fwd, 1, m_traceNormals[n], 1,
264  GradJumpOnTraceBwd, 1, GradJumpOnTraceBwd, 1);
265  }
266  else
267  {
268  // Multiply by normal and add to trace evaluation
269  Vmath::Vsub(nTracePts, Fwd, 1, Bwd, 1, Fwd, 1);
270  Vmath::Vvtvp(nTracePts, Fwd, 1, m_traceNormals[n], 1,
271  GradJumpOnTrace, 1, GradJumpOnTrace, 1);
272  }
273  }
274 
276  {
277  // Need to negate Bwd case when using Fwd normal
278  Vmath::Neg(nTracePts, GradJumpOnTrace, 1);
279  }
280 
281  Vmath::Vmul(nTracePts, unorm, 1, GradJumpOnTrace, 1, GradJumpOnTrace, 1);
282 
283  // Interpolate GradJumpOnTrace to Local elemental traces.
284  m_locTraceToTraceMap->InterpTraceToLocTrace(0, GradJumpOnTrace, wsp);
285  m_locTraceToTraceMap->UnshuffleLocTraces(0, wsp, LocElmtTracePhys);
286 
288  {
289  // Vmath::Neg(nTracePts,GradJumpOnTraceBwd,1);
290  Vmath::Vmul(nTracePts, unorm, 1, GradJumpOnTraceBwd, 1,
291  GradJumpOnTraceBwd, 1);
292  m_locTraceToTraceMap->InterpTraceToLocTrace(1, GradJumpOnTraceBwd, wsp);
293  m_locTraceToTraceMap->UnshuffleLocTraces(1, wsp, LocElmtTracePhys);
294  }
295  else
296  {
297  m_locTraceToTraceMap->InterpTraceToLocTrace(1, GradJumpOnTrace, wsp);
298  m_locTraceToTraceMap->UnshuffleLocTraces(1, wsp, LocElmtTracePhys);
299  }
300 
301  // Scale jump on trace
302  Vmath::Vmul(nLocETrace, m_scalTrace[0], 1, LocElmtTracePhys, 1, wsp, 1);
303  MultiplyByStdDerivBaseOnTraceMat(0, wsp, FilterCoeffs);
304 
305  for (int i = 0; i < m_traceDim; ++i)
306  {
307  // Scale jump on trace
308  Vmath::Vmul(nLocETrace, m_scalTrace[i + 1], 1, LocElmtTracePhys, 1, wsp,
309  1);
310  MultiplyByStdDerivBaseOnTraceMat(i + 1, wsp, deriv[0]);
311  Vmath::Vadd(ncoeffs, deriv[0], 1, FilterCoeffs, 1, FilterCoeffs, 1);
312  }
313 
314  Vmath::Svtvp(ncoeffs, scale, FilterCoeffs, 1, outarray, 1, outarray, 1);
315 }
#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:209
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:622
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:518
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:574
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:359
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:492
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:419

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  {
72  return m_useGJPSemiImplicit;
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 376 of file GJPStabilisation.cpp.

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

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 317 of file GJPStabilisation.cpp.

319 {
320 
321  // check to see if already deifned and if so return
322  if (graph->ExpansionInfoDefined("GJP"))
323  {
324  return;
325  }
326 
327  const SpatialDomains::ExpansionInfoMap expInfo =
328  graph->GetExpansionInfo(variable);
329 
332 
333  // loop over epxansion info
334  for (auto expIt = expInfo.begin(); expIt != expInfo.end(); ++expIt)
335  {
336  std::vector<LibUtilities::BasisKey> BKeyVector;
337 
338  for (int i = 0; i < expIt->second->m_basisKeyVector.size(); ++i)
339  {
340  LibUtilities::BasisKey bkeyold = expIt->second->m_basisKeyVector[i];
341 
342  // Reset radauM alpha non-zero cases to radauM
343  // Legendre at one order higher
344 
345  switch (bkeyold.GetPointsType())
346  {
347  case LibUtilities::eGaussRadauMAlpha1Beta0:
348  case LibUtilities::eGaussRadauMAlpha2Beta0:
349  {
350  int npts = bkeyold.GetNumPoints();
351 
352  // const LibUtilities::PointsKey pkey(npts+1,
353  // LibUtilities::eGaussRadauMLegendre); trying
354  // npts to be consistent for tri faces
355  const LibUtilities::PointsKey pkey(
357  LibUtilities::BasisKey bkeynew(bkeyold.GetBasisType(),
358  bkeyold.GetNumModes(), pkey);
359  BKeyVector.push_back(bkeynew);
360  }
361  break;
362  default:
363  BKeyVector.push_back(bkeyold);
364  break;
365  }
366  }
367 
368  (*newInfo)[expIt->first] =
370  expIt->second->m_geomShPtr, BKeyVector);
371  }
372 
373  graph->SetExpansionInfo("GJP", newInfo);
374 }
@ 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

◆ 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 84 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 88 of file GJPStabilisation.h.

Referenced by Apply(), and GJPStabilisation().

◆ m_locTraceToTraceMap

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

LocaTraceToTraceMap.

Definition at line 86 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 92 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 95 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 81 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().