Nektar++
Expansion1D.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: Expansion1D.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: File for Expansion1D routines
32//
33///////////////////////////////////////////////////////////////////////////////
34
37
38using namespace std;
39
41{
42
44{
45 DNekMatSharedPtr returnval;
46
47 switch (mkey.GetMatrixType())
48 {
50 {
52 "HybridDGHelmholtz matrix not set up "
53 "for non boundary-interior expansions");
54 int i;
55 NekDouble lambdaval =
58 int ncoeffs = GetNcoeffs();
59
60 int coordim = GetCoordim();
61
66 DNekMat LocMat(ncoeffs, ncoeffs);
67
68 returnval =
70 DNekMat &Mat = *returnval;
71
72 Vmath::Zero(ncoeffs * ncoeffs, Mat.GetPtr(), 1);
73
74 for (i = 0; i < coordim; ++i)
75 {
76 DNekScalMat &Dmat = *GetLocMatrix(DerivType[i]);
77
78 Mat = Mat + Dmat * invMass * Transpose(Dmat);
79 }
80
81 // Add end Mass Matrix Contribution
83 Mat = Mat + lambdaval * Mass;
84
86 GetBoundaryMap(bmap);
87
88 // Add tau*F_e using elemental mass matrices
89 for (i = 0; i < 2; ++i)
90 {
91 Mat(bmap[i], bmap[i]) = Mat(bmap[i], bmap[i]) + tau;
92 }
93 }
94 break;
96 {
97 int j, k;
98 int nbndry = NumDGBndryCoeffs();
99 int ncoeffs = GetNcoeffs();
105
106 Array<OneD, NekDouble> lambda(nbndry);
107 DNekVec Lambda(nbndry, lambda, eWrapper);
108 Array<OneD, NekDouble> ulam(ncoeffs);
109 DNekVec Ulam(ncoeffs, ulam, eWrapper);
110 Array<OneD, NekDouble> f(ncoeffs);
111 DNekVec F(ncoeffs, f, eWrapper);
112
113 // declare matrix space
114 returnval =
116 DNekMat &Umat = *returnval;
117
118 // Helmholtz matrix
119 DNekScalMat &invHmat =
121
122 // for each degree of freedom of the lambda space
123 // calculate Umat entry
124 // Generate Lambda to U_lambda matrix
125 for (j = 0; j < nbndry; ++j)
126 {
127 Vmath::Zero(nbndry, &lambda[0], 1);
128 Vmath::Zero(ncoeffs, &f[0], 1);
129 lambda[j] = 1.0;
130
132 lambda, f);
133
134 Ulam = invHmat * F; // generate Ulam from lambda
135
136 // fill column of matrix
137 for (k = 0; k < ncoeffs; ++k)
138 {
139 Umat(k, j) = Ulam[k];
140 }
141 }
142 }
143 break;
147 {
148 int j = 0;
149 int k = 0;
150 int dir = 0;
151 int nbndry = NumDGBndryCoeffs();
152 int ncoeffs = GetNcoeffs();
153
154 Array<OneD, NekDouble> lambda(nbndry);
155 DNekVec Lambda(nbndry, lambda, eWrapper);
156
157 Array<OneD, NekDouble> ulam(ncoeffs);
158 DNekVec Ulam(ncoeffs, ulam, eWrapper);
159 Array<OneD, NekDouble> f(ncoeffs);
160 DNekVec F(ncoeffs, f, eWrapper);
166
167 // declare matrix space
168 returnval =
170 DNekMat &Qmat = *returnval;
171
172 // Lambda to U matrix
173 DNekScalMat &lamToU =
175
176 // Inverse mass matrix
178
179 // Weak Derivative matrix
181 switch (mkey.GetMatrixType())
182 {
184 dir = 0;
186 break;
188 dir = 1;
190 break;
192 dir = 2;
194 break;
195 default:
196 ASSERTL0(false, "Direction not known");
197 break;
198 }
199
200 // for each degree of freedom of the lambda space
201 // calculate Qmat entry
202 // Generate Lambda to Q_lambda matrix
203 for (j = 0; j < nbndry; ++j)
204 {
205 Vmath::Zero(nbndry, &lambda[0], 1);
206 lambda[j] = 1.0;
207
208 // for lambda[j] = 1 this is the solution to ulam
209 for (k = 0; k < ncoeffs; ++k)
210 {
211 Ulam[k] = lamToU(k, j);
212 }
213
214 // -D^T ulam
215 Vmath::Neg(ncoeffs, &ulam[0], 1);
216 F = Transpose(*Dmat) * Ulam;
217
218 // + \tilde{G} \lambda
219 AddNormTraceInt(dir, lambda, f);
220
221 // multiply by inverse mass matrix
222 Ulam = invMass * F;
223
224 // fill column of matrix (Qmat is in column major format)
225 Vmath::Vcopy(ncoeffs, &ulam[0], 1,
226 &(Qmat.GetPtr())[0] + j * ncoeffs, 1);
227 }
228 }
229 break;
231 {
232 int j;
233 int nbndry = NumBndryCoeffs();
234
240
243 GetBoundaryMap(bmap);
244
245 // declare matrix space
246 returnval =
248 DNekMat &BndMat = *returnval;
249
250 // Matrix to map Lambda to U
251 DNekScalMat &LamToU =
253
254 // Matrix to map Lambda to Q
255 DNekScalMat &LamToQ =
257
258 lam[0] = 1.0;
259 lam[1] = 0.0;
260 for (j = 0; j < nbndry; ++j)
261 {
262 BndMat(0, j) =
263 -LamToQ(bmap[0], j) - factors[StdRegions::eFactorTau] *
264 (LamToU(bmap[0], j) - lam[j]);
265 }
266
267 lam[0] = 0.0;
268 lam[1] = 1.0;
269 for (j = 0; j < nbndry; ++j)
270 {
271 BndMat(1, j) =
272 LamToQ(bmap[1], j) - factors[StdRegions::eFactorTau] *
273 (LamToU(bmap[1], j) - lam[j]);
274 }
275 }
276 break;
277 default:
278 ASSERTL0(false,
279 "This matrix type cannot be generated from this class");
280 break;
281 }
282
283 return returnval;
284}
285
286void Expansion1D::AddNormTraceInt([[maybe_unused]] const int dir,
288 Array<OneD, NekDouble> &outarray)
289{
290 int k;
291 int nbndry = NumBndryCoeffs();
292 int nquad = GetNumPoints(0);
293 const Array<OneD, const NekDouble> &Basis = GetBasis(0)->GetBdata();
295
296 GetBoundaryMap(vmap);
297
298 // add G \lambda term (can assume G is diagonal since one
299 // of the basis is zero at boundary otherwise)
300 for (k = 0; k < nbndry; ++k)
301 {
302 outarray[vmap[k]] += (Basis[(vmap[k] + 1) * nquad - 1] *
303 Basis[(vmap[k] + 1) * nquad - 1] -
304 Basis[vmap[k] * nquad] * Basis[vmap[k] * nquad]) *
305 inarray[vmap[k]];
306 }
307}
308
310 const NekDouble tau, const Array<OneD, const NekDouble> &inarray,
311 Array<OneD, NekDouble> &outarray)
312{
313 int i, n;
314 int nbndry = NumBndryCoeffs();
315 int nquad = GetNumPoints(0);
316 int ncoeffs = GetNcoeffs();
317 int coordim = GetCoordim();
319
320 ASSERTL0(&inarray[0] != &outarray[0],
321 "Input and output arrays use the same memory");
322
323 const Array<OneD, const NekDouble> &Basis = GetBasis(0)->GetBdata();
325
326 GetBoundaryMap(vmap);
327
328 // Add F = \tau <phi_i,phi_j> (note phi_i is zero if phi_j is non-zero)
329 for (i = 0; i < nbndry; ++i)
330 {
331 outarray[vmap[i]] += tau * Basis[(vmap[i] + 1) * nquad - 1] *
332 Basis[(vmap[i] + 1) * nquad - 1] *
333 inarray[vmap[i]];
334 outarray[vmap[i]] += tau * Basis[vmap[i] * nquad] *
335 Basis[vmap[i] * nquad] * inarray[vmap[i]];
336 }
337
338 //===============================================================
339 // Add -\sum_i D_i^T M^{-1} G_i + E_i M^{-1} G_i =
340 // \sum_i D_i M^{-1} G_i term
341
345 Array<OneD, NekDouble> tmpcoeff(ncoeffs, 0.0);
346 DNekVec Coeffs(ncoeffs, outarray, eWrapper);
347 DNekVec Tmpcoeff(ncoeffs, tmpcoeff, eWrapper);
348
349 for (n = 0; n < coordim; ++n)
350 {
351 // evaluate M^{-1} G
352 for (i = 0; i < ncoeffs; ++i)
353 {
354 // lower boundary (negative normal)
355 tmpcoeff[i] -= invMass(i, vmap[0]) * Basis[vmap[0] * nquad] *
356 Basis[vmap[0] * nquad] * inarray[vmap[0]];
357
358 // upper boundary (positive normal)
359 tmpcoeff[i] += invMass(i, vmap[1]) *
360 Basis[(vmap[1] + 1) * nquad - 1] *
361 Basis[(vmap[1] + 1) * nquad - 1] * inarray[vmap[1]];
362 }
363
364 DNekScalMat &Dmat = *GetLocMatrix(DerivType[n]);
365 Coeffs = Coeffs + Dmat * Tmpcoeff;
366 }
367}
368
370 const int vert, const Array<OneD, const NekDouble> &primCoeffs,
371 DNekMatSharedPtr &inoutmat)
372{
374 "Robin boundary conditions are only implemented for "
375 "boundary-interior expanisons");
376 ASSERTL1(inoutmat->GetRows() == inoutmat->GetColumns(),
377 "Assuming that input matrix was square");
378
379 // Get local Element mapping for vertex point
380 int map = GetVertexMap(vert);
381
382 // Now need to identify a map which takes the local edge
383 // mass matrix to the matrix stored in inoutmat;
384 // This can currently be deduced from the size of the matrix
385 // - if inoutmat.m_rows() == v_NCoeffs() it is a full
386 // matrix system
387 // - if inoutmat.m_rows() == v_NumBndCoeffs() it is a
388 // boundary CG system
389
390 int rows = inoutmat->GetRows();
391
392 if (rows == GetNcoeffs())
393 {
394 // no need to do anything
395 }
396 else if (rows == NumBndryCoeffs()) // same as NumDGBndryCoeffs()
397 {
398 int i;
400 GetBoundaryMap(bmap);
401
402 for (i = 0; i < 2; ++i)
403 {
404 if (map == bmap[i])
405 {
406 map = i;
407 break;
408 }
409 }
410 ASSERTL1(i != 2, "Did not find number in map");
411 }
412
413 // assumes end points have unit magnitude
414 (*inoutmat)(map, map) += primCoeffs[0];
415}
416
417/**
418 * Given an edge and vector of element coefficients:
419 * - maps those elemental coefficients corresponding to the trace into
420 * an vector.
421 * - update the element coefficients
422 * - multiplies the edge vector by the edge mass matrix
423 * - maps the edge coefficients back onto the elemental coefficients
424 */
426 const int vert, const Array<OneD, const NekDouble> &primCoeffs,
427 const Array<OneD, NekDouble> &incoeffs, Array<OneD, NekDouble> &coeffs)
428{
430 "Not set up for non boundary-interior expansions");
431
432 int map = GetVertexMap(vert);
433 coeffs[map] += primCoeffs[0] * incoeffs[map];
434}
435
437 const Array<OneD, Array<OneD, NekDouble>> &vec)
438{
440 GetLeftAdjacentElementExp()->GetTraceNormal(
442
443 int nq = m_base[0]->GetNumPoints();
445 Vmath::Vmul(nq, &vec[0][0], 1, &normals[0][0], 1, &Fn[0], 1);
446 Vmath::Vvtvp(nq, &vec[1][0], 1, &normals[1][0], 1, &Fn[0], 1, &Fn[0], 1);
447
448 return Integral(Fn);
449}
450
451/** @brief: This method gets all of the factors which are
452 required as part of the Gradient Jump Penalty
453 stabilisation and involves the product of the normal and
454 geometric factors along the element trace.
455*/
458 [[maybe_unused]] Array<OneD, Array<OneD, NekDouble>> &d0factors,
459 [[maybe_unused]] Array<OneD, Array<OneD, NekDouble>> &d1factors)
460{
461 int nquad = GetNumPoints(0);
463 m_metricinfo->GetDerivFactors(GetPointsKeys());
464
465 if (factors.size() <= 2)
466 {
470 }
471
472 // Outwards normal
477
478 if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
479 {
480 factors[0][0] = gmat[0][nquad - 1] * normal_0[0][0];
481 factors[1][0] = gmat[0][0] * normal_1[0][0];
482
483 for (int n = 1; n < normal_0.size(); ++n)
484 {
485 factors[0][0] += gmat[n][0] * normal_0[n][0];
486 factors[1][0] += gmat[n][nquad - 1] * normal_1[n][0];
487 }
488 }
489 else
490 {
491 factors[0][0] = gmat[0][0] * normal_0[0][0];
492 factors[1][0] = gmat[0][0] * normal_1[0][0];
493
494 for (int n = 1; n < normal_0.size(); ++n)
495 {
496 factors[0][0] += gmat[n][0] * normal_0[n][0];
497 factors[1][0] += gmat[n][0] * normal_1[n][0];
498 }
499 }
500}
501
503 [[maybe_unused]] const StdRegions::Orientation orient,
504 Array<OneD, int> &idmap, [[maybe_unused]] const int nq0,
505 [[maybe_unused]] const int nq1)
506{
507 if (idmap.size() != 1)
508 {
509 idmap = Array<OneD, int>(1);
510 }
511
512 idmap[0] = 0;
513}
514
515void Expansion1D::v_TraceNormLen([[maybe_unused]] const int traceid,
516 NekDouble &h, NekDouble &p)
517{
518 h = GetGeom()->GetVertex(1)->dist(*GetGeom()->GetVertex(0));
519 p = m_ncoeffs - 1;
520}
521
522} // namespace Nektar::LocalRegions
#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
Represents a basis of a given type.
Definition: Basis.h:198
void v_AddRobinTraceContribution(const int vert, const Array< OneD, const NekDouble > &primCoeffs, const Array< OneD, NekDouble > &incoeffs, Array< OneD, NekDouble > &coeffs) override
void v_ReOrientTracePhysMap(const StdRegions::Orientation orient, Array< OneD, int > &idmap, const int nq0, const int nq1) override
void AddHDGHelmholtzTraceTerms(const NekDouble tau, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
void v_AddRobinMassMatrix(const int vert, const Array< OneD, const NekDouble > &primCoeffs, DNekMatSharedPtr &inoutmat) override
void AddNormTraceInt(const int dir, Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
NekDouble v_VectorFlux(const Array< OneD, Array< OneD, NekDouble > > &vec) override
void v_TraceNormLen(const int traceid, NekDouble &h, NekDouble &p) override
void v_NormalTraceDerivFactors(Array< OneD, Array< OneD, NekDouble > > &factors, Array< OneD, Array< OneD, NekDouble > > &d0factors, Array< OneD, Array< OneD, NekDouble > > &d1factors) override
: This method gets all of the factors which are required as part of the Gradient Jump Penalty stabili...
DNekMatSharedPtr v_GenMatrix(const StdRegions::StdMatrixKey &mkey) override
Definition: Expansion1D.cpp:43
SpatialDomains::GeometrySharedPtr GetGeom() const
Definition: Expansion.cpp:167
ExpansionSharedPtr GetLeftAdjacentElementExp() const
Definition: Expansion.h:441
int GetLeftAdjacentElementTrace() const
Definition: Expansion.h:454
SpatialDomains::GeomFactorsSharedPtr m_metricinfo
Definition: Expansion.h:274
DNekScalMatSharedPtr GetLocMatrix(const LocalRegions::MatrixKey &mkey)
Definition: Expansion.cpp:84
const NormalVector & GetTraceNormal(const int id)
Definition: Expansion.cpp:251
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
void GetBoundaryMap(Array< OneD, unsigned int > &outarray)
Definition: StdExpansion.h:669
const LibUtilities::BasisSharedPtr & GetBasis(int dir) const
This function gets the shared point to basis in the dir direction.
Definition: StdExpansion.h:112
int GetNcoeffs(void) const
This function returns the total number of coefficients used in the expansion.
Definition: StdExpansion.h:124
const LibUtilities::PointsKeyVector GetPointsKeys() const
int GetVertexMap(const int localVertexId, bool useCoeffPacking=false)
Definition: StdExpansion.h:679
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Definition: StdExpansion.h:218
Array< OneD, LibUtilities::BasisSharedPtr > m_base
NekDouble Integral(const Array< OneD, const NekDouble > &inarray)
This function integrates the specified function over the domain.
Definition: StdExpansion.h:473
MatrixType GetMatrixType() const
Definition: StdMatrixKey.h:83
NekDouble GetConstFactor(const ConstFactorType &factor) const
Definition: StdMatrixKey.h:124
@ eDeformed
Geometry is curved or has non-constant factors.
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:430
StdRegions::ConstFactorMap factors
std::shared_ptr< DNekScalMat > DNekScalMatSharedPtr
NekMatrix< InnerMatrixType, BlockMatrixTag > Transpose(NekMatrix< InnerMatrixType, BlockMatrixTag > &rhs)
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75
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 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 Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825
STL namespace.