Nektar++
GJPStabilisation.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: GJPStabilisation.cpp
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: GJP data
32//
33///////////////////////////////////////////////////////////////////////////////
34
36
38{
41 "GJPStabilisation", "Explicit", eExplicitGJPStabilisation),
43 "GJPStabilisation", "SemiImplicit", eSemiImplicitGJPStabilisation),
44};
46{
47 LibUtilities::SessionReaderSharedPtr session = pField->GetSession();
48
49 session->MatchSolverInfo("GJPStabilisation", "SemiImplicit",
51
52 // Call GetTrace on the initialising field will set up
53 // DG. Store a copoy so that if we make a soft copy of
54 // this class we can re-used this field for operators.
55 pField->GetTrace();
56 m_dgfield = pField;
57
58 m_coordDim = m_dgfield->GetCoordim(0);
59 m_traceDim = m_dgfield->GetShapeDimension() - 1;
60
61 // set up trace normals but would be better if could use
62 // equation system definition
64 for (int i = 0; i < m_coordDim; ++i)
65 {
67 Array<OneD, NekDouble>(m_dgfield->GetTrace()->GetNpoints());
68 }
69 m_dgfield->GetTrace()->GetNormals(m_traceNormals);
70
71 SetUpExpansionInfoMapForGJP(pField->GetGraph(), session->GetVariable(0));
72
74
76 session, pField->GetGraph(), "GJP", true, false,
77 Collections::eNoImpType, session->GetVariable(0));
78 dgfield->GetLocTraceToTraceMap(m_locTraceToTraceMap);
79
81 session, *(dgfield->GetExp()), dgfield->GetGraph(), true, "GJP");
82
84
85 const std::shared_ptr<LocalRegions::ExpansionVector> exp =
86 dgfield->GetExp();
87
91
92 for (int i = 0; i < m_traceDim + 1; ++i)
93 {
94 LocTrace[i] =
96 }
97
98 int cnt = 0;
99 int offset_phys = 0;
102
104 Array<OneD, int> sign, sign1;
105 NekDouble h, p;
106
107 for (int e = 0; e < m_dgfield->GetExpSize(); ++e)
108 {
109 LocalRegions::ExpansionSharedPtr elmt = (*exp)[e];
110
111 elmt->NormalTraceDerivFactors(dfactors[0], dfactors[1], dfactors[2]);
112
113 for (int n = 0; n < elmt->GetNtraces(); ++n, ++cnt)
114 {
115 NekDouble jumpScal;
116 elmt->TraceNormLen(n, h, p);
117 ASSERTL0(boost::math::isnan(h) == false,
118 "h has a nan value when e = " + std::to_string(e) +
119 " n =" + std::to_string(n));
120
121 if (p == 1)
122 {
123 jumpScal = 0.02 * h * h;
124 }
125 else
126 {
127 jumpScal = 0.8 * pow(p + 1, -4.0) * h * h;
128 }
129
130 int nptrace = elmt->GetTraceNumPoints(n);
131 elmt->GetTraceCoeffMap(n, map);
132
133 for (int i = 0; i < m_traceDim + 1; ++i)
134 {
135 Vmath::Smul(nptrace, jumpScal, dfactors[i][n], 1,
136 e_tmp = LocTrace[i] + offset_phys, 1);
137 }
138
139 offset_phys += nptrace;
140 }
141 }
142
143 for (int i = 0; i < m_traceDim + 1; ++i)
144 {
145 m_scalTrace[i] = LocTrace[i];
146
147 if (m_traceDim > 0)
148 {
149 // multiply by Jacobian and quadrature points.
150 m_locElmtTrace->MultiplyByQuadratureMetric(m_scalTrace[i],
151 m_scalTrace[i]);
152 }
153 }
154
155 // Assemble list of Matrix Product
157
158 int nelmt = 1;
160 dgfield->GetExp(0)->GetBase();
161
162 dgfield->GetExp(0)->StdDerivBaseOnTraceMat(TraceMat);
163
164 for (int n = 1; n < dgfield->GetExpSize(); ++n)
165 {
167 dgfield->GetExp(n)->GetBase();
168
169 int i;
170 for (i = 0; i < base.size(); ++i)
171 {
172 if (base[i] != base_sav[i])
173 {
174 break;
175 }
176 }
177
178 if (i == base.size())
179 {
180 nelmt++;
181 }
182 else
183 {
184 // save previous block of data.
185 m_StdDBaseOnTraceMat.push_back(
186 std::pair<int, Array<OneD, DNekMatSharedPtr>>(nelmt, TraceMat));
187
188 // start new block
189 dgfield->GetExp(n)->StdDerivBaseOnTraceMat(TraceMat);
190 nelmt = 1;
191 base_sav = dgfield->GetExp(n)->GetBase();
192 }
193 }
194 // save latest block of data.
195 m_StdDBaseOnTraceMat.push_back(
196 std::pair<int, Array<OneD, DNekMatSharedPtr>>(nelmt, TraceMat));
197}
198
200 Array<OneD, NekDouble> &outarray,
201 const Array<OneD, NekDouble> &pUnorm,
202 NekDouble scale) const
203{
204 int ncoeffs = m_dgfield->GetNcoeffs();
205 int nphys = m_dgfield->GetNpoints();
206 int nTracePts = m_dgfield->GetTrace()->GetTotPoints();
207 int nLocETrace = m_locElmtTrace->GetTotPoints();
208
209 ASSERTL1(nLocETrace == m_scalTrace[0].size(), "expect these to be similar");
210 ASSERTL1(m_locElmtTrace->GetNcoeffs() <= nphys,
211 "storage assumptions presume "
212 "that nLocETraceCoeffs < nphys");
213
215 for (int i = 0; i < m_coordDim; ++i)
216 {
217 deriv[i] = Array<OneD, NekDouble>(nphys);
218 }
219
220 int nmax = std::max(ncoeffs, nphys);
221 Array<OneD, NekDouble> FilterCoeffs(nmax, 0.0);
222 Array<OneD, NekDouble> GradJumpOnTrace(nTracePts, 0.0);
223 Array<OneD, NekDouble> Fwd(nTracePts), Bwd(nTracePts);
224
225 Array<OneD, NekDouble> wsp(nLocETrace), tmp;
226 Array<OneD, NekDouble> LocElmtTracePhys = m_locElmtTrace->UpdatePhys();
227 Array<OneD, NekDouble> LocElmtTraceCoeffs = m_locElmtTrace->UpdateCoeffs();
228
229 ASSERTL1(LocElmtTracePhys.size() <= nLocETrace,
230 "expect this vector to be at least of size nLocETrace");
231
233 if (pUnorm == NullNekDouble1DArray)
234 {
235 unorm = Array<OneD, NekDouble>(nTracePts, 1.0);
236 }
237 else
238 {
239 unorm = pUnorm;
240 }
241
242 Array<OneD, NekDouble> GradJumpOnTraceBwd;
244 {
245 GradJumpOnTraceBwd = Array<OneD, NekDouble>(nTracePts);
246 }
247
249 {
250 Vmath::Zero(nTracePts, GradJumpOnTraceBwd, 1);
251 }
252
253 // calculate derivative
254 m_dgfield->PhysDeriv(inarray, deriv[0], deriv[1], deriv[2]);
255
256 // Evaluate the normal derivative jump on the trace
257 for (int n = 0; n < m_coordDim; ++n)
258 {
259 m_dgfield->GetFwdBwdTracePhys(deriv[n], Fwd, Bwd, true, true);
260
262 {
263 // want to put Fwd vals on bwd trace and vice versa
264 Vmath::Vvtvp(nTracePts, Bwd, 1, m_traceNormals[n], 1,
265 GradJumpOnTrace, 1, GradJumpOnTrace, 1);
266 Vmath::Vvtvp(nTracePts, Fwd, 1, m_traceNormals[n], 1,
267 GradJumpOnTraceBwd, 1, GradJumpOnTraceBwd, 1);
268 }
269 else
270 {
271 // Multiply by normal and add to trace evaluation
272 Vmath::Vsub(nTracePts, Fwd, 1, Bwd, 1, Fwd, 1);
273 Vmath::Vvtvp(nTracePts, Fwd, 1, m_traceNormals[n], 1,
274 GradJumpOnTrace, 1, GradJumpOnTrace, 1);
275 }
276 }
277
279 {
280 // Need to negate Bwd case when using Fwd normal
281 Vmath::Neg(nTracePts, GradJumpOnTrace, 1);
282 }
283
284 Vmath::Vmul(nTracePts, unorm, 1, GradJumpOnTrace, 1, GradJumpOnTrace, 1);
285
286 // Interpolate GradJumpOnTrace to Local elemental traces.
287 m_locTraceToTraceMap->InterpTraceToLocTrace(0, GradJumpOnTrace, wsp);
288 m_locTraceToTraceMap->UnshuffleLocTraces(0, wsp, LocElmtTracePhys);
289
291 {
292 // Vmath::Neg(nTracePts,GradJumpOnTraceBwd,1);
293 Vmath::Vmul(nTracePts, unorm, 1, GradJumpOnTraceBwd, 1,
294 GradJumpOnTraceBwd, 1);
295 m_locTraceToTraceMap->InterpTraceToLocTrace(1, GradJumpOnTraceBwd, wsp);
296 m_locTraceToTraceMap->UnshuffleLocTraces(1, wsp, LocElmtTracePhys);
297 }
298 else
299 {
300 m_locTraceToTraceMap->InterpTraceToLocTrace(1, GradJumpOnTrace, wsp);
301 m_locTraceToTraceMap->UnshuffleLocTraces(1, wsp, LocElmtTracePhys);
302 }
303
304 // Scale jump on trace
305 Vmath::Vmul(nLocETrace, m_scalTrace[0], 1, LocElmtTracePhys, 1, wsp, 1);
306 MultiplyByStdDerivBaseOnTraceMat(0, wsp, FilterCoeffs);
307
308 for (int i = 0; i < m_traceDim; ++i)
309 {
310 // Scale jump on trace
311 Vmath::Vmul(nLocETrace, m_scalTrace[i + 1], 1, LocElmtTracePhys, 1, wsp,
312 1);
313 MultiplyByStdDerivBaseOnTraceMat(i + 1, wsp, deriv[0]);
314 Vmath::Vadd(ncoeffs, deriv[0], 1, FilterCoeffs, 1, FilterCoeffs, 1);
315 }
316
317 Vmath::Svtvp(ncoeffs, scale, FilterCoeffs, 1, outarray, 1, outarray, 1);
318}
319
321 SpatialDomains::MeshGraphSharedPtr graph, std::string variable)
322{
323
324 // check to see if already deifned and if so return
325 if (graph->ExpansionInfoDefined("GJP"))
326 {
327 return;
328 }
329
331 graph->GetExpansionInfo(variable);
332
335
336 // loop over epxansion info
337 for (auto expIt = expInfo.begin(); expIt != expInfo.end(); ++expIt)
338 {
339 std::vector<LibUtilities::BasisKey> BKeyVector;
340
341 for (int i = 0; i < expIt->second->m_basisKeyVector.size(); ++i)
342 {
343 LibUtilities::BasisKey bkeyold = expIt->second->m_basisKeyVector[i];
344
345 // Reset radauM alpha non-zero cases to radauM
346 // Legendre at one order higher
347
348 switch (bkeyold.GetPointsType())
349 {
350 case LibUtilities::eGaussRadauMAlpha1Beta0:
351 case LibUtilities::eGaussRadauMAlpha2Beta0:
352 {
353 int npts = bkeyold.GetNumPoints();
354
355 // const LibUtilities::PointsKey pkey(npts+1,
356 // LibUtilities::eGaussRadauMLegendre); trying
357 // npts to be consistent for tri faces
358 const LibUtilities::PointsKey pkey(
360 LibUtilities::BasisKey bkeynew(bkeyold.GetBasisType(),
361 bkeyold.GetNumModes(), pkey);
362 BKeyVector.push_back(bkeynew);
363 }
364 break;
365 default:
366 BKeyVector.push_back(bkeyold);
367 break;
368 }
369 }
370
371 (*newInfo)[expIt->first] =
373 expIt->second->m_geomShPtr, BKeyVector);
374 }
375
376 graph->SetExpansionInfo("GJP", newInfo);
377}
378
380 int i, Array<OneD, NekDouble> &in, Array<OneD, NekDouble> &out) const
381{
382 // Should probably be vectorised
383
384 int cnt = 0;
385 int cnt1 = 0;
386 for (auto &it : m_StdDBaseOnTraceMat)
387 {
388 int rows = it.second[i]->GetRows();
389 int cols = it.second[i]->GetColumns();
390
391 Blas::Dgemm('N', 'N', rows, it.first, cols, 1.0,
392 &(it.second[i]->GetPtr())[0], rows, &in[0] + cnt, cols, 0.0,
393 &out[0] + cnt1, rows);
394
395 cnt += cols * it.first;
396 cnt1 += rows * it.first;
397 }
398}
399} // namespace Nektar::MultiRegions
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:47
Describes the specification for a Basis.
Definition: Basis.h:45
int GetNumPoints() const
Return points order at which basis is defined.
Definition: Basis.h:120
BasisType GetBasisType() const
Return type of expansion basis.
Definition: Basis.h:131
int GetNumModes() const
Returns the order of the basis.
Definition: Basis.h:74
PointsType GetPointsType() const
Return type of quadrature.
Definition: Basis.h:143
Defines a specification for a set of points.
Definition: Points.h:50
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[]
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.
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:383
std::shared_ptr< SessionReader > SessionReaderSharedPtr
@ eGaussRadauMLegendre
1D Gauss-Radau-Legendre quadrature points, pinned at x=-1
Definition: PointsType.h:47
std::shared_ptr< Expansion > ExpansionSharedPtr
Definition: Expansion.h:66
std::shared_ptr< DisContField > DisContFieldSharedPtr
Definition: DisContField.h:346
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
std::shared_ptr< ExpansionInfoMap > ExpansionInfoMapShPtr
Definition: MeshGraph.h:143
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
std::map< int, ExpansionInfoShPtr > ExpansionInfoMap
Definition: MeshGraph.h:141
static Array< OneD, NekDouble > NullNekDouble1DArray
double NekDouble
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.hpp:72
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.hpp:396
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.hpp:292
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.hpp:366
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.hpp:180
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.hpp:100
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273
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.hpp:220