Nektar++
ContField.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: ContField.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: Field definition for a continuous domain with boundary
32// conditions
33//
34///////////////////////////////////////////////////////////////////////////////
35
38#include <tuple>
39
40using namespace std;
41
43{
44/**
45 * @class ContField
46 * The class #ContField is
47 * able to incorporate the boundary conditions imposed to the problem
48 * to be solved. Therefore, the class is equipped with three additional
49 * data members:
50 * - #m_bndCondExpansions
51 * - #m_bndTypes
52 * - #m_bndCondEquations
53 *
54 * The first data structure, #m_bndCondExpansions, contains the
55 * one-dimensional spectral/hp expansion on the boundary, #m_bndTypes
56 * stores information about the type of boundary condition on the
57 * different parts of the boundary while #m_bndCondEquations holds the
58 * equation of the imposed boundary conditions.
59 *
60 * Furthermore, in case of Dirichlet boundary conditions, this class is
61 * capable of lifting a known solution satisfying these boundary
62 * conditions. If we denote the unknown solution by
63 * \f$u^{\mathcal{H}}(\boldsymbol{x})\f$ and the known Dirichlet
64 * boundary conditions by \f$u^{\mathcal{D}}(\boldsymbol{x})\f$, the
65 * expansion then can be decomposed as
66 * \f[ u^{\delta}(\boldsymbol{x}_i)=u^{\mathcal{D}}(\boldsymbol{x}_i)+
67 * u^{\mathcal{H}}(\boldsymbol{x}_i)=\sum_{n=0}^{N^{\mathcal{D}}-1}
68 * \hat{u}_n^{\mathcal{D}}\Phi_n(\boldsymbol{x}_i)+
69 * \sum_{n={N^{\mathcal{D}}}}^{N_{\mathrm{dof}}-1}
70 * \hat{u}_n^{\mathcal{H}} \Phi_n(\boldsymbol{x}_i).\f]
71 * This lifting is accomplished by ordering the known global degrees of
72 * freedom, prescribed by the Dirichlet boundary conditions, first in
73 * the global array
74 * \f$\boldsymbol{\hat{u}}\f$, that is,
75 * \f[\boldsymbol{\hat{u}}=\left[ \begin{array}{c}
76 * \boldsymbol{\hat{u}}^{\mathcal{D}}\\
77 * \boldsymbol{\hat{u}}^{\mathcal{H}}
78 * \end{array} \right].\f]
79 * Such kind of expansions are also referred to as continuous fields.
80 * This class should be used when solving 2D problems using a standard
81 * Galerkin approach.
82 */
83
84/**
85 *
86 */
88 : DisContField(), m_locToGloMap(), m_globalMat(),
89 m_globalLinSysManager(
90 std::bind(&ContField::GenGlobalLinSys, this, std::placeholders::_1),
91 std::string("GlobalLinSys"))
92{
93}
94
95/**
96 * Given a mesh \a graph, containing information about the domain and
97 * the spectral/hp element expansion, this constructor fills the list
98 * of local expansions #m_exp with the proper expansions, calculates
99 * the total number of quadrature points \f$\boldsymbol{x}_i\f$ and
100 * local expansion coefficients \f$\hat{u}^e_n\f$ and allocates memory
101 * for the arrays #m_coeffs and #m_phys. Furthermore, it constructs the
102 * mapping array (contained in #m_locToGloMap) for the transformation
103 * between local elemental level and global level, it calculates the
104 * total number global expansion coefficients \f$\hat{u}_n\f$ and
105 * allocates memory for the array #m_contCoeffs. The constructor also
106 * discretises the boundary conditions, specified by the argument \a
107 * bcs, by expressing them in terms of the coefficient of the expansion
108 * on the boundary.
109 *
110 * @param graph A mesh, containing information about the domain
111 * and the spectral/hp element expansion.
112 * @param bcs The boundary conditions.
113 * @param variable An optional parameter to indicate for which
114 * variable the field should be constructed.
115 */
118 const std::string &variable,
119 const bool DeclareCoeffPhysArrays,
120 const bool CheckIfSingularSystem,
122 : DisContField(pSession, graph, variable, false, DeclareCoeffPhysArrays,
123 ImpType),
124 m_globalMat(MemoryManager<GlobalMatrixMap>::AllocateSharedPtr()),
125 m_globalLinSysManager(
126 std::bind(&ContField::GenGlobalLinSys, this, std::placeholders::_1),
127 std::string("GlobalLinSys"))
128{
131 CheckIfSingularSystem, variable, m_periodicVerts, m_periodicEdges,
133
134 if (m_session->DefinesCmdLineArgument("verbose"))
135 {
136 m_locToGloMap->PrintStats(std::cout, variable);
137 }
138}
139
140/**
141 * Given a mesh \a graph, containing information about the domain and
142 * the spectral/hp element expansion, this constructor fills the list
143 * of local expansions #m_exp with the proper expansions, calculates
144 * the total number of quadrature points \f$\boldsymbol{x}_i\f$ and
145 * local expansion coefficients \f$\hat{u}^e_n\f$ and allocates memory
146 * for the arrays #m_coeffs and #m_phys. Furthermore, it constructs the
147 * mapping array (contained in #m_locToGloMap) for the transformation
148 * between local elemental level and global level, it calculates the
149 * total number global expansion coefficients \f$\hat{u}_n\f$ and
150 * allocates memory for the array #m_coeffs. The constructor also
151 * discretises the boundary conditions, specified by the argument \a
152 * bcs, by expressing them in terms of the coefficient of the expansion
153 * on the boundary.
154 *
155 * @param In Existing ContField object used to provide the
156 * local to global mapping information and
157 * global solution type.
158 * @param graph A mesh, containing information about the domain
159 * and the spectral/hp element expansion.
160 * @param bcs The boundary conditions.
161 * @param bc_loc
162 */
165 const std::string &variable, bool DeclareCoeffPhysArrays,
166 const bool CheckIfSingularSystem)
167 : DisContField(In, graph, variable, false, DeclareCoeffPhysArrays),
168 m_globalMat(MemoryManager<GlobalMatrixMap>::AllocateSharedPtr()),
169 m_globalLinSysManager(
170 std::bind(&ContField::GenGlobalLinSys, this, std::placeholders::_1),
171 std::string("GlobalLinSys")),
172 m_GJPData(In.m_GJPData)
173{
174 if (!SameTypeOfBoundaryConditions(In) || CheckIfSingularSystem)
175 {
178 CheckIfSingularSystem, variable, m_periodicVerts, m_periodicEdges,
180
181 if (m_session->DefinesCmdLineArgument("verbose"))
182 {
183 m_locToGloMap->PrintStats(std::cout, variable);
184 }
185 }
186 else
187 {
189 }
190}
191
192/**
193 * Initialises the object as a copy of an existing ContField object.
194 * @param In Existing ContField object.
195 * @param DeclareCoeffPhysArrays bool to declare if \a m_phys
196 * and \a m_coeffs should be declared. Default is true
197 */
198ContField::ContField(const ContField &In, bool DeclareCoeffPhysArrays)
199 : DisContField(In, DeclareCoeffPhysArrays), m_locToGloMap(In.m_locToGloMap),
200 m_globalMat(In.m_globalMat),
201 m_globalLinSysManager(In.m_globalLinSysManager), m_GJPData(In.m_GJPData)
202{
203}
204
205/**
206 * Constructs a continuous field as a copy of an existing
207 * explist field and adding all the boundary conditions.
208 *
209 * @param In Existing explist1D field .
210 */
212 const ExpList &In)
213 : DisContField(In), m_locToGloMap(),
214 m_globalLinSysManager(
215 std::bind(&ContField::GenGlobalLinSys, this, std::placeholders::_1),
216 std::string("GlobalLinSys"))
217{
219 pSession, m_ncoeffs, In);
220}
221
222/**
223 *
224 */
226{
227}
228
229/**
230 * Given a function \f$f(\boldsymbol{x})\f$ defined at the quadrature
231 * points, this function determines the unknown global coefficients
232 * \f$\boldsymbol{\hat{u}}^{\mathcal{H}}\f$ employing a discrete
233 * Galerkin projection from physical space to coefficient
234 * space. The operation is evaluated by the function #GlobalSolve using
235 * the global mass matrix.
236 *
237 * The values of the function \f$f(\boldsymbol{x})\f$ evaluated at the
238 * quadrature points \f$\boldsymbol{x}_i\f$ should be contained in the
239 * variable #inarray of the ExpList object \a Sin. The resulting global
240 * coefficients \f$\hat{u}_g\f$ are stored in the array #outarray.
241 */
243 Array<OneD, NekDouble> &outarray)
244
245{
246 // Inner product of forcing
248 IProductWRTBase(inarray, wsp);
249
250 // Solve the system
252
253 GlobalSolve(key, wsp, outarray);
254}
255
256/**
257 *
258 */
260{
261 int Ncoeffs = m_locToGloMap->GetNumLocalCoeffs();
262 Array<OneD, NekDouble> tmp1(Ncoeffs);
263 Array<OneD, NekDouble> tmp2(Ncoeffs);
264
265 IProductWRTBase(field, tmp1);
266 MultiplyByInvMassMatrix(tmp1, tmp2);
267 BwdTrans(tmp2, field);
268}
269
270/**
271 * Computes the matrix vector product
272 * @f$ \mathbf{y} = \mathbf{M}^{-1}\mathbf{x} @f$.
273 *
274 * @param inarray Input vector @f$\mathbf{x}@f$.
275 * @param outarray Output vector @f$\mathbf{y}@f$.
276 */
278 const Array<OneD, const NekDouble> &inarray,
279 Array<OneD, NekDouble> &outarray)
280
281{
283 GlobalSolve(key, inarray, outarray);
284}
285
286/**
287 * Consider the two dimensional Laplace equation,
288 * \f[\nabla\cdot\left(\boldsymbol{\sigma}\nabla
289 * u(\boldsymbol{x})\right) = f(\boldsymbol{x}),\f] supplemented with
290 * appropriate boundary conditions (which are contained in the data
291 * member #m_bndCondExpansions). In the equation above
292 * \f$\boldsymbol{\sigma}\f$ is the (symmetric positive definite)
293 * diffusion tensor:
294 * \f[ \sigma = \left[ \begin{array}{cc}
295 * \sigma_{00}(\boldsymbol{x},t) & \sigma_{01}(\boldsymbol{x},t) \\
296 * \sigma_{01}(\boldsymbol{x},t) & \sigma_{11}(\boldsymbol{x},t)
297 * \end{array} \right]. \f]
298 * Applying a \f$C^0\f$ continuous Galerkin discretisation, this
299 * equation leads to the following linear system:
300 * \f[\boldsymbol{L}
301 * \boldsymbol{\hat{u}}_g=\boldsymbol{\hat{f}}\f]
302 * where \f$\boldsymbol{L}\f$ is the Laplacian matrix. This function
303 * solves the system above for the global coefficients
304 * \f$\boldsymbol{\hat{u}}\f$ by a call to the function #GlobalSolve.
305 *
306 * The values of the function \f$f(\boldsymbol{x})\f$ evaluated at the
307 * quadrature points \f$\boldsymbol{x}_i\f$ should be contained in the
308 * variable #inarray
309 *
310 * @param inarray An Array<OneD, NekDouble> containing the discrete
311 * evaluation of the forcing function
312 * \f$f(\boldsymbol{x})\f$ at the quadrature
313 * points.
314 *
315 * @param outarray An Array<OneD, NekDouble> containing the
316 * coefficients of the solution
317 *
318 * @param variablecoeffs The (optional) parameter containing the
319 * coefficients evaluated at the quadrature
320 * points. It is an Array of (three) arrays which
321 * stores the laplacian coefficients in the
322 * following way
323 * \f[\mathrm{variablecoeffs} = \left[ \begin{array}{c}
324 * \left[\sigma_{00}(\boldsymbol{x_i},t)\right]_i \\
325 * \left[\sigma_{01}(\boldsymbol{x_i},t)\right]_i \\
326 * \left[\sigma_{11}(\boldsymbol{x_i},t)\right]_i
327 * \end{array}\right]
328 * \f]
329 * If this argument is not passed to the function, the following
330 * equation will be solved:
331 * \f[\nabla^2u(\boldsymbol{x}) = f(\boldsymbol{x}),\f]
332 *
333 * @param time The time-level at which the coefficients are
334 * evaluated
335 */
337 const Array<OneD, const NekDouble> &inarray,
338 Array<OneD, NekDouble> &outarray,
339 const Array<OneD, const NekDouble> &dirForcing,
340 const Array<OneD, Array<OneD, NekDouble>> &variablecoeffs, NekDouble time)
341{
342 // Inner product of forcing
344 IProductWRTBase(inarray, wsp);
345
346 // Note -1.0 term necessary to invert forcing function to
347 // be consistent with matrix definition
348 Vmath::Neg(m_ncoeffs, wsp, 1);
349
350 // Forcing function with weak boundary conditions
351 int i, j;
352 int bndcnt = 0;
354 m_locToGloMap->GetBndCondCoeffsToLocalCoeffsSign();
355 const Array<OneD, const int> map =
356 m_locToGloMap->GetBndCondCoeffsToLocalCoeffsMap();
357
358 // Add weak boundary conditions to forcing
359 for (i = 0; i < m_bndCondExpansions.size(); ++i)
360 {
361 if (m_bndConditions[i]->GetBoundaryConditionType() ==
363 m_bndConditions[i]->GetBoundaryConditionType() ==
365 {
366
367 const Array<OneD, const NekDouble> bndcoeff =
369
370 if (m_locToGloMap->GetSignChange())
371 {
372 for (j = 0; j < (m_bndCondExpansions[i])->GetNcoeffs(); j++)
373 {
374 wsp[map[bndcnt + j]] += sign[bndcnt + j] * bndcoeff[j];
375 }
376 }
377 else
378 {
379 for (j = 0; j < (m_bndCondExpansions[i])->GetNcoeffs(); j++)
380 {
381 wsp[map[bndcnt + j]] += bndcoeff[bndcnt + j];
382 }
383 }
384 }
385
386 bndcnt += m_bndCondExpansions[i]->GetNcoeffs();
387 }
388
389 StdRegions::VarCoeffMap varcoeffs;
390 varcoeffs[StdRegions::eVarCoeffD00] = variablecoeffs[0];
391 varcoeffs[StdRegions::eVarCoeffD01] = variablecoeffs[1];
392 varcoeffs[StdRegions::eVarCoeffD11] = variablecoeffs[3];
393 varcoeffs[StdRegions::eVarCoeffD22] = variablecoeffs[5];
396
397 // Solve the system
399 varcoeffs);
400
401 GlobalSolve(key, wsp, outarray, dirForcing);
402}
403
404/**
405 * Constructs the GlobalLinearSysKey for the linear advection operator
406 * with the supplied parameters, and computes the eigenvectors and
407 * eigenvalues of the associated matrix.
408 * @param ax Advection parameter, x.
409 * @param ay Advection parameter, y.
410 * @param Real Computed eigenvalues, real component.
411 * @param Imag Computed eigenvalues, imag component.
412 * @param Evecs Computed eigenvectors.
413 */
418{
419 // Solve the system
423 vel[0] = vel_x;
424 vel[1] = vel_y;
425
426 StdRegions::VarCoeffMap varcoeffs;
427 varcoeffs[StdRegions::eVarCoeffVelX] =
429 varcoeffs[StdRegions::eVarCoeffVelY] =
434 factors, varcoeffs);
435
437 Gmat->EigenSolve(Real, Imag, Evecs);
438}
439
440/**
441 * Given a linear system specified by the key \a key,
442 * \f[\boldsymbol{M}\boldsymbol{\hat{u}}_g=\boldsymbol{\hat{f}},\f]
443 * this function solves this linear system taking into account the
444 * boundary conditions specified in the data member
445 * #m_bndCondExpansions. Therefore, it adds an array
446 * \f$\boldsymbol{\hat{g}}\f$ which represents the non-zero surface
447 * integral resulting from the weak boundary conditions (e.g. Neumann
448 * boundary conditions) to the right hand side, that is,
449 * \f[\boldsymbol{M}\boldsymbol{\hat{u}}_g=\boldsymbol{\hat{f}}+
450 * \boldsymbol{\hat{g}}.\f]
451 * Furthermore, it lifts the known degrees of freedom which are
452 * prescribed by the Dirichlet boundary conditions. As these known
453 * coefficients \f$\boldsymbol{\hat{u}}^{\mathcal{D}}\f$ are numbered
454 * first in the global coefficient array \f$\boldsymbol{\hat{u}}_g\f$,
455 * the linear system can be decomposed as,
456 * \f[\left[\begin{array}{cc}
457 * \boldsymbol{M}^{\mathcal{DD}}&\boldsymbol{M}^{\mathcal{DH}}\\
458 * \boldsymbol{M}^{\mathcal{HD}}&\boldsymbol{M}^{\mathcal{HH}}
459 * \end{array}\right]
460 * \left[\begin{array}{c}
461 * \boldsymbol{\hat{u}}^{\mathcal{D}}\\
462 * \boldsymbol{\hat{u}}^{\mathcal{H}}
463 * \end{array}\right]=
464 * \left[\begin{array}{c}
465 * \boldsymbol{\hat{f}}^{\mathcal{D}}\\
466 * \boldsymbol{\hat{f}}^{\mathcal{H}}
467 * \end{array}\right]+
468 * \left[\begin{array}{c}
469 * \boldsymbol{\hat{g}}^{\mathcal{D}}\\
470 * \boldsymbol{\hat{g}}^{\mathcal{H}}
471 * \end{array}\right]
472 * \f]
473 * which will then be solved for the unknown coefficients
474 * \f$\boldsymbol{\hat{u}}^{\mathcal{H}}\f$ as,
475 * \f[
476 * \boldsymbol{M}^{\mathcal{HH}}\boldsymbol{\hat{u}}^{\mathcal{H}}=
477 * \boldsymbol{\hat{f}}^{\mathcal{H}}+
478 * \boldsymbol{\hat{g}}^{\mathcal{H}}-
479 * \boldsymbol{M}^{\mathcal{HD}}\boldsymbol{\hat{u}}^{\mathcal{D}}\f]
480 *
481 * @param mkey This key uniquely defines the linear system to
482 * be solved.
483 * @param locrhs contains the forcing term in local coefficient space
484 * @note inout contains initial guess and final output in local coeffs.
485 */
487 const Array<OneD, const NekDouble> &locrhs,
489 const Array<OneD, const NekDouble> &dirForcing)
490{
491 int NumDirBcs = m_locToGloMap->GetNumGlobalDirBndCoeffs();
492 int contNcoeffs = m_locToGloMap->GetNumGlobalCoeffs();
493
494 // STEP 1: SET THE DIRICHLET DOFS TO THE RIGHT VALUE
495 // IN THE SOLUTION ARRAY
497
498 // STEP 2: CALCULATE THE HOMOGENEOUS COEFFICIENTS
499 if (contNcoeffs - NumDirBcs > 0)
500 {
502 LinSys->Solve(locrhs, inout, m_locToGloMap, dirForcing);
503 }
504}
505
506/**
507 * Returns the global matrix associated with the given GlobalMatrixKey.
508 * If the global matrix has not yet been constructed on this field,
509 * it is first constructed using GenGlobalMatrix().
510 * @param mkey Global matrix key.
511 * @returns Assocated global matrix.
512 */
514{
516 "To use method must have a AssemblyMap "
517 "attached to key");
518
519 GlobalMatrixSharedPtr glo_matrix;
520 auto matrixIter = m_globalMat->find(mkey);
521
522 if (matrixIter == m_globalMat->end())
523 {
524 glo_matrix = GenGlobalMatrix(mkey, m_locToGloMap);
525 (*m_globalMat)[mkey] = glo_matrix;
526 }
527 else
528 {
529 glo_matrix = matrixIter->second;
530 }
531
532 return glo_matrix;
533}
534
535/**
536 * The function searches the map #m_globalLinSys to see if the
537 * global matrix has been created before. If not, it calls the function
538 * #GenGlobalLinSys to generate the requested global system.
539 *
540 * @param mkey This key uniquely defines the requested
541 * linear system.
542 */
544{
545 return m_globalLinSysManager[mkey];
546}
547
549{
551 "To use method must have a AssemblyMap "
552 "attached to key");
554}
555
557{
558 int i, j;
559 int bndcnt = 0;
560
562 m_locToGloMap->GetBndCondCoeffsToLocalCoeffsSign();
563 const Array<OneD, const int> map =
564 m_locToGloMap->GetBndCondCoeffsToLocalCoeffsMap();
565
566 for (i = 0; i < m_bndCondExpansions.size(); ++i)
567 {
568 if (m_bndConditions[i]->GetBoundaryConditionType() ==
570 {
571
572 const Array<OneD, const NekDouble> bndcoeff =
574
575 if (m_locToGloMap->GetSignChange())
576 {
577 for (j = 0; j < (m_bndCondExpansions[i])->GetNcoeffs(); j++)
578 {
579 outarray[map[bndcnt + j]] = sign[bndcnt + j] * bndcoeff[j];
580 }
581 }
582 else
583 {
584 for (j = 0; j < (m_bndCondExpansions[i])->GetNcoeffs(); j++)
585 {
586 outarray[map[bndcnt + j]] = bndcoeff[j];
587 }
588 }
589 }
590 bndcnt += m_bndCondExpansions[i]->GetNcoeffs();
591 }
592
593 // communicate local Dirichlet coeffs that are just
594 // touching a dirichlet boundary on another partition
595 set<int> &ParallelDirBndSign = m_locToGloMap->GetParallelDirBndSign();
596
597 for (auto &it : ParallelDirBndSign)
598 {
599 outarray[it] *= -1;
600 }
601
602 m_locToGloMap->UniversalAbsMaxBnd(outarray);
603
604 for (auto &it : ParallelDirBndSign)
605 {
606 outarray[it] *= -1;
607 }
608
609 set<ExtraDirDof> &copyLocalDirDofs = m_locToGloMap->GetCopyLocalDirDofs();
610 for (auto &it : copyLocalDirDofs)
611 {
612 outarray[std::get<0>(it)] = outarray[std::get<1>(it)] * std::get<2>(it);
613 }
614}
615
617{
618 int bndcnt = 0;
619
621 m_locToGloMap->GetBndCondCoeffsToLocalCoeffsSign();
622 const Array<OneD, const int> bndmap =
623 m_locToGloMap->GetBndCondCoeffsToLocalCoeffsMap();
624
625 for (int i = 0; i < m_bndCondExpansions.size(); ++i)
626 {
627 Array<OneD, NekDouble> &bcoeffs =
628 m_bndCondExpansions[i]->UpdateCoeffs();
629
630 if (m_locToGloMap->GetSignChange())
631 {
632 for (int j = 0; j < (m_bndCondExpansions[i])->GetNcoeffs(); ++j)
633 {
634 bcoeffs[j] = sign[bndcnt + j] * coeffs[bndmap[bndcnt + j]];
635 }
636 }
637 else
638 {
639 for (int j = 0; j < (m_bndCondExpansions[i])->GetNcoeffs(); ++j)
640 {
641 bcoeffs[j] = coeffs[bndmap[bndcnt + j]];
642 }
643 }
644
645 bndcnt += m_bndCondExpansions[i]->GetNcoeffs();
646 }
647}
648
650 const Array<OneD, NekDouble> coeffs)
651{
652 int bndcnt = 0;
653
654 ASSERTL1(nreg < m_bndCondExpansions.size(),
655 "nreg is out or range since this many boundary "
656 "regions to not exist");
657
659 m_locToGloMap->GetBndCondCoeffsToLocalCoeffsSign();
660 const Array<OneD, const int> bndmap =
661 m_locToGloMap->GetBndCondCoeffsToLocalCoeffsMap();
662
663 // Now fill in all other Dirichlet coefficients.
664 Array<OneD, NekDouble> &bcoeffs = m_bndCondExpansions[nreg]->UpdateCoeffs();
665
666 for (int j = 0; j < nreg; ++j)
667 {
668 bndcnt += m_bndCondExpansions[j]->GetNcoeffs();
669 }
670
671 if (m_locToGloMap->GetSignChange())
672 {
673 for (int j = 0; j < (m_bndCondExpansions[nreg])->GetNcoeffs(); ++j)
674 {
675 bcoeffs[j] = sign[bndcnt + j] * coeffs[bndmap[bndcnt + j]];
676 }
677 }
678 else
679 {
680 for (int j = 0; j < (m_bndCondExpansions[nreg])->GetNcoeffs(); ++j)
681 {
682 bcoeffs[j] = coeffs[bndmap[bndcnt + j]];
683 }
684 }
685}
686
687/**
688 * This operation is evaluated as:
689 * \f{tabbing}
690 * \hspace{1cm} \= Do \= $e=$ $1, N_{\mathrm{el}}$ \\
691 * > > Do \= $i=$ $0,N_m^e-1$ \\
692 * > > > $\boldsymbol{\hat{u}}^{e}[i] = \mbox{sign}[e][i] \cdot
693 * \boldsymbol{\hat{u}}_g[\mbox{map}[e][i]]$ \\
694 * > > continue \\
695 * > continue
696 * \f}
697 * where \a map\f$[e][i]\f$ is the mapping array and \a
698 * sign\f$[e][i]\f$ is an array of similar dimensions ensuring the
699 * correct modal connectivity between the different elements (both
700 * these arrays are contained in the data member #m_locToGloMap). This
701 * operation is equivalent to the scatter operation
702 * \f$\boldsymbol{\hat{u}}_l=\mathcal{A}\boldsymbol{\hat{u}}_g\f$,
703 * where \f$\mathcal{A}\f$ is the
704 * \f$N_{\mathrm{eof}}\times N_{\mathrm{dof}}\f$ permutation matrix.
705 *
706 */
708 Array<OneD, NekDouble> &outarray)
709{
710 m_locToGloMap->GlobalToLocal(inarray, outarray);
711}
712
714{
715 m_locToGloMap->GlobalToLocal(m_coeffs, m_coeffs);
716}
717
718/**
719 * This operation is evaluated as:
720 * \f{tabbing}
721 * \hspace{1cm} \= Do \= $e=$ $1, N_{\mathrm{el}}$ \\
722 * > > Do \= $i=$ $0,N_m^e-1$ \\
723 * > > > $\boldsymbol{\hat{u}}_g[\mbox{map}[e][i]] =
724 * \mbox{sign}[e][i] \cdot \boldsymbol{\hat{u}}^{e}[i]$\\
725 * > > continue\\
726 * > continue
727 * \f}
728 * where \a map\f$[e][i]\f$ is the mapping array and \a
729 * sign\f$[e][i]\f$ is an array of similar dimensions ensuring the
730 * correct modal connectivity between the different elements (both
731 * these arrays are contained in the data member #m_locToGloMap). This
732 * operation is equivalent to the gather operation
733 * \f$\boldsymbol{\hat{u}}_g=\mathcal{A}^{-1}\boldsymbol{\hat{u}}_l\f$,
734 * where \f$\mathcal{A}\f$ is the
735 * \f$N_{\mathrm{eof}}\times N_{\mathrm{dof}}\f$ permutation matrix.
736 *
737 */
738
740 Array<OneD, NekDouble> &outarray, bool useComm)
741{
742 m_locToGloMap->LocalToGlobal(inarray, outarray, useComm);
743}
744
746
747{
748 m_locToGloMap->LocalToGlobal(m_coeffs, m_coeffs, useComm);
749}
750
751/**
752 * Consider the two dimensional Helmholtz equation,
753 * \f[\nabla^2u(\boldsymbol{x})-\lambda u(\boldsymbol{x})
754 * = f(\boldsymbol{x}),\f] supplemented with appropriate boundary
755 * conditions (which are contained in the data member
756 * #m_bndCondExpansions). Applying a \f$C^0\f$ continuous Galerkin
757 * discretisation, this equation leads to the following linear system:
758 * \f[\left(\boldsymbol{L}+\lambda\boldsymbol{M}\right)
759 * \boldsymbol{\hat{u}}_g=\boldsymbol{\hat{f}}\f] where
760 * \f$\boldsymbol{L}\f$ and \f$\boldsymbol{M}\f$ are the Laplacian and
761 * mass matrix respectively. This function solves the system above for
762 * the global coefficients \f$\boldsymbol{\hat{u}}\f$ by a call to the
763 * function #GlobalSolve.
764 *
765 * @param inarray An Array<OneD, NekDouble> , containing the discrete
766 * evaluation of the forcing function
767 * \f$f(\boldsymbol{x})\f$ at the quadrature
768 * points
769 * @param factors The parameter \f$\lambda\f$ of the Helmholtz
770 * equation is specified through the factors map
771 */
773 const Array<OneD, const NekDouble> &inarray,
775 const StdRegions::VarCoeffMap &pvarcoeff,
776 const MultiRegions::VarFactorsMap &varfactors,
777 const Array<OneD, const NekDouble> &dirForcing, const bool PhysSpaceForcing)
778{
779 int i, j;
780
781 //----------------------------------
782 // Setup RHS Inner product
783 //----------------------------------
784 // Inner product of forcing
786 if (PhysSpaceForcing)
787 {
788 IProductWRTBase(inarray, wsp);
789 // Note -1.0 term necessary to invert forcing function to
790 // be consistent with matrix definition
791 Vmath::Neg(m_ncoeffs, wsp, 1);
792 }
793 else
794 {
795 Vmath::Smul(m_ncoeffs, -1.0, inarray, 1, wsp, 1);
796 }
797
798 int bndcnt = 0;
800 m_locToGloMap->GetBndCondCoeffsToLocalCoeffsSign();
801 const Array<OneD, const int> map =
802 m_locToGloMap->GetBndCondCoeffsToLocalCoeffsMap();
803 // Add weak boundary conditions to forcing
804 for (i = 0; i < m_bndCondExpansions.size(); ++i)
805 {
806 if (m_bndConditions[i]->GetBoundaryConditionType() ==
808 m_bndConditions[i]->GetBoundaryConditionType() ==
810 {
811
812 const Array<OneD, const NekDouble> bndcoeff =
814
815 if (m_locToGloMap->GetSignChange())
816 {
817 for (j = 0; j < (m_bndCondExpansions[i])->GetNcoeffs(); j++)
818 {
819 wsp[map[bndcnt + j]] += sign[bndcnt + j] * bndcoeff[j];
820 }
821 }
822 else
823 {
824 for (j = 0; j < (m_bndCondExpansions[i])->GetNcoeffs(); j++)
825 {
826 wsp[map[bndcnt + j]] += bndcoeff[j];
827 }
828 }
829 }
830 bndcnt += m_bndCondExpansions[i]->GetNcoeffs();
831 }
832
834
835 StdRegions::VarCoeffMap varcoeff(pvarcoeff);
837 {
838 // initialize if required
839 if (!m_GJPData)
840 {
843 }
844
845 if (m_GJPData->IsSemiImplicit())
846 {
848 }
849
850 // to set up forcing need initial guess in physical space
852 BwdTrans(outarray, phys);
853 NekDouble scale = -1.0 * factors.find(StdRegions::eFactorGJP)->second;
854
855 m_GJPData->Apply(phys, wsp,
856 pvarcoeff.count(StdRegions::eVarCoeffGJPNormVel)
857 ? pvarcoeff.find(StdRegions::eVarCoeffGJPNormVel)
858 ->second.GetValue()
860 scale);
861
862 varcoeff.erase(StdRegions::eVarCoeffGJPNormVel);
863 }
864
865 GlobalLinSysKey key(mtype, m_locToGloMap, factors, varcoeff, varfactors);
866
867 GlobalSolve(key, wsp, outarray, dirForcing);
868
869 return key;
870}
871
872/**
873 * First compute the inner product of forcing function with respect to
874 * base, and then solve the system with the linear advection operator.
875 * @param velocity Array of advection velocities in physical space
876 * @param inarray Forcing function.
877 * @param outarray Result.
878 * @param lambda reaction coefficient
879 * @param dirForcing Dirichlet Forcing.
880 */
881
882// could combine this with HelmholtzCG.
884 const Array<OneD, const NekDouble> &inarray,
886 const StdRegions::VarCoeffMap &pvarcoeff,
887 const MultiRegions::VarFactorsMap &varfactors,
888 const Array<OneD, const NekDouble> &dirForcing, const bool PhysSpaceForcing)
889{
890 // Inner product of forcing
892 if (PhysSpaceForcing)
893 {
894 IProductWRTBase(inarray, wsp);
895 // Note -1.0 term necessary to invert forcing function to
896 // be consistent with matrix definition
897 Vmath::Neg(m_ncoeffs, wsp, 1);
898 }
899 else
900 {
901 Vmath::Smul(m_ncoeffs, -1.0, inarray, 1, wsp, 1);
902 }
903
904 // Forcing function with weak boundary conditions
905 int i, j;
906 int bndcnt = 0;
908 m_locToGloMap->GetBndCondCoeffsToLocalCoeffsSign();
909 const Array<OneD, const int> map =
910 m_locToGloMap->GetBndCondCoeffsToLocalCoeffsMap();
911 // Add weak boundary conditions to forcing
912 for (i = 0; i < m_bndCondExpansions.size(); ++i)
913 {
914 if (m_bndConditions[i]->GetBoundaryConditionType() ==
916 m_bndConditions[i]->GetBoundaryConditionType() ==
918 {
919
920 const Array<OneD, const NekDouble> bndcoeff =
922
923 if (m_locToGloMap->GetSignChange())
924 {
925 for (j = 0; j < (m_bndCondExpansions[i])->GetNcoeffs(); j++)
926 {
927 wsp[map[bndcnt + j]] += sign[bndcnt + j] * bndcoeff[j];
928 }
929 }
930 else
931 {
932 for (j = 0; j < (m_bndCondExpansions[i])->GetNcoeffs(); j++)
933 {
934 wsp[map[bndcnt + j]] += bndcoeff[bndcnt + j];
935 }
936 }
937 }
938
939 bndcnt += m_bndCondExpansions[i]->GetNcoeffs();
940 }
941
944
945 StdRegions::VarCoeffMap varcoeff(pvarcoeff);
947 {
948 // initialize if required
949 if (!m_GJPData)
950 {
953 }
954
955 if (m_GJPData->IsSemiImplicit())
956 {
958 }
959
960 // to set up forcing need initial guess in physical space
962 BwdTrans(outarray, phys);
963 NekDouble scale = -1.0 * factors.find(StdRegions::eFactorGJP)->second;
964
965 m_GJPData->Apply(phys, wsp,
966 pvarcoeff.count(StdRegions::eVarCoeffGJPNormVel)
967 ? pvarcoeff.find(StdRegions::eVarCoeffGJPNormVel)
968 ->second.GetValue()
970 scale);
971
972 varcoeff.erase(StdRegions::eVarCoeffGJPNormVel);
973 }
974
975 // Solve the system
976 GlobalLinSysKey key(mtype, m_locToGloMap, factors, varcoeff, varfactors);
977
978 GlobalSolve(key, wsp, outarray, dirForcing);
979
980 return key;
981}
982
983/**
984 * First compute the inner product of forcing function with respect to
985 * base, and then solve the system with the linear advection operator.
986 * @param velocity Array of advection velocities in physical space
987 * @param inarray Forcing function.
988 * @param outarray Result.
989 * @param lambda reaction coefficient
990 * @param dirForcing Dirichlet Forcing.
991 */
993 const Array<OneD, Array<OneD, NekDouble>> &velocity,
994 const Array<OneD, const NekDouble> &inarray,
995 Array<OneD, NekDouble> &outarray, const NekDouble lambda,
996 const Array<OneD, const NekDouble> &dirForcing)
997{
998 // Inner product of forcing
1000 IProductWRTBase(inarray, wsp);
1001
1002 // Solve the system
1005 StdRegions::VarCoeffMap varcoeffs;
1006 varcoeffs[StdRegions::eVarCoeffVelX] = velocity[0];
1007 varcoeffs[StdRegions::eVarCoeffVelY] = velocity[1];
1009 factors, varcoeffs);
1010
1011 GlobalSolve(key, wsp, outarray, dirForcing);
1012}
1013
1014/**
1015 * Reset the GlobalLinSys Manager
1016 */
1018{
1019 m_globalLinSysManager.ClearManager("GlobalLinSys");
1020}
1021
1022/**
1023 * Get the pool count for the specified poolName
1024 */
1025int ContField::v_GetPoolCount(std::string poolName)
1026{
1027 return m_globalLinSysManager.PoolCount(poolName);
1028}
1029
1030/**
1031 * Clear all memory for GlobalLinSys
1032 * including StaticCond Blocks and LocalMatrix Blocks.
1033 * Avoids memory leakage if matrices are updated in time
1034 */
1036 bool clearLocalMatrices)
1037{
1038 // Get GlobalLinSys from key
1040
1041 // Loop all expansions
1042 for (int n = 0; n < m_exp->size(); ++n)
1043 {
1044 LinSys->DropStaticCondBlock(n);
1045
1046 if (clearLocalMatrices)
1047 {
1048 LinSys->DropBlock(n);
1049 }
1050 }
1051
1052 m_globalLinSysManager.DeleteObject(key);
1053}
1054
1055} // namespace Nektar::MultiRegions
#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
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
This class is the abstraction of a global continuous two- dimensional spectral/hp element expansion w...
Definition: ContField.h:55
~ContField() override
The default destructor.
Definition: ContField.cpp:225
void v_FillBndCondFromField(const Array< OneD, NekDouble > coeffs) override
Definition: ContField.cpp:616
LibUtilities::NekManager< GlobalLinSysKey, GlobalLinSys > m_globalLinSysManager
A manager which collects all the global linear systems being assembled, such that they should be cons...
Definition: ContField.h:161
void LaplaceSolve(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray, const Array< OneD, Array< OneD, NekDouble > > &variablecoeffs=NullNekDoubleArrayOfArray, NekDouble time=0.0)
Solves the two-dimensional Laplace equation, subject to the boundary conditions specified.
Definition: ContField.cpp:336
GlobalLinSysSharedPtr GenGlobalLinSys(const GlobalLinSysKey &mkey)
Definition: ContField.cpp:548
GlobalLinSysKey v_HelmSolve(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::ConstFactorMap &factors, const StdRegions::VarCoeffMap &varcoeff, const MultiRegions::VarFactorsMap &varfactors, const Array< OneD, const NekDouble > &dirForcing, const bool PhysSpaceForcing) override
Solves the two-dimensional Helmholtz equation, subject to the boundary conditions specified.
Definition: ContField.cpp:772
AssemblyMapCGSharedPtr m_locToGloMap
(A shared pointer to) the object which contains all the required information for the transformation f...
Definition: ContField.h:150
void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Template method virtual forwarder for FwdTrans().
Definition: ContField.cpp:242
int v_GetPoolCount(std::string) override
Definition: ContField.cpp:1025
void v_SmoothField(Array< OneD, NekDouble > &field) override
Template method virtual forwarded for SmoothField().
Definition: ContField.cpp:259
void v_ClearGlobalLinSysManager(void) override
Definition: ContField.cpp:1017
void v_GlobalToLocal(void) override
Definition: ContField.cpp:713
void v_UnsetGlobalLinSys(GlobalLinSysKey, bool) override
Definition: ContField.cpp:1035
void v_LocalToGlobal(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool useComm) override
Gathers the global coefficients from the local coefficients .
Definition: ContField.cpp:739
GlobalMatrixMapShPtr m_globalMat
(A shared pointer to) a list which collects all the global matrices being assembled,...
Definition: ContField.h:155
void v_LinearAdvectionReactionSolve(const Array< OneD, Array< OneD, NekDouble > > &velocity, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const NekDouble lambda, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray) override
Definition: ContField.cpp:992
GlobalLinSysSharedPtr GetGlobalLinSys(const GlobalLinSysKey &mkey)
Returns the linear system specified by the key mkey.
Definition: ContField.cpp:543
void v_ImposeDirichletConditions(Array< OneD, NekDouble > &outarray) override
Impose the Dirichlet Boundary Conditions on outarray.
Definition: ContField.cpp:556
void v_MultiplyByInvMassMatrix(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Template method virtual forwarder for MultiplyByInvMassMatrix().
Definition: ContField.cpp:277
GlobalMatrixSharedPtr GetGlobalMatrix(const GlobalMatrixKey &mkey)
Returns the global matrix specified by mkey.
Definition: ContField.cpp:513
ContField()
The default constructor.
Definition: ContField.cpp:87
void LinearAdvectionEigs(const NekDouble ax, const NekDouble ay, Array< OneD, NekDouble > &Real, Array< OneD, NekDouble > &Imag, Array< OneD, NekDouble > &Evecs=NullNekDouble1DArray)
Compute the eigenvalues of the linear advection operator.
Definition: ContField.cpp:414
GlobalLinSysKey v_LinearAdvectionDiffusionReactionSolve(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::ConstFactorMap &factors, const StdRegions::VarCoeffMap &varcoeff, const MultiRegions::VarFactorsMap &varfactors, const Array< OneD, const NekDouble > &dirForcing, const bool PhysSpaceForcing) override
Definition: ContField.cpp:883
void GlobalSolve(const GlobalLinSysKey &key, const Array< OneD, const NekDouble > &rhs, Array< OneD, NekDouble > &inout, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
Solves the linear system specified by the key key.
Definition: ContField.cpp:486
GJPStabilisationSharedPtr m_GJPData
Data for Gradient Jump Penalisation (GJP) stabilisaiton.
Definition: ContField.h:164
This class is the abstractio n of a global discontinuous two- dimensional spectral/hp element expansi...
Definition: DisContField.h:56
PeriodicMap m_periodicEdges
A map which identifies pairs of periodic edges.
Definition: DisContField.h:185
PeriodicMap m_periodicFaces
A map which identifies pairs of periodic faces.
Definition: DisContField.h:190
bool SameTypeOfBoundaryConditions(const DisContField &In)
Check to see if expansion has the same BCs as In.
Array< OneD, SpatialDomains::BoundaryConditionShPtr > m_bndConditions
An array which contains the information about the boundary condition structure definition on the diff...
Definition: DisContField.h:143
PeriodicMap m_periodicVerts
A map which identifies groups of periodic vertices.
Definition: DisContField.h:180
Array< OneD, MultiRegions::ExpListSharedPtr > m_bndCondExpansions
An object which contains the discretised boundary conditions.
Definition: DisContField.h:156
Base class for all multi-elemental spectral/hp expansions.
Definition: ExpList.h:99
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:1083
void IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function calculates the inner product of a function with respect to all local expansion modes .
Definition: ExpList.h:1650
int GetNcoeffs(void) const
Returns the total number of local degrees of freedom .
Definition: ExpList.h:1512
std::shared_ptr< GlobalMatrix > GenGlobalMatrix(const GlobalMatrixKey &mkey, const std::shared_ptr< AssemblyMapCG > &locToGloMap)
Generates a global matrix from the given key and map.
Definition: ExpList.cpp:2786
const Array< OneD, const NekDouble > & GetCoeffs() const
This function returns (a reference to) the array (implemented as m_coeffs) containing all local expa...
Definition: ExpList.h:1944
std::shared_ptr< ExpList > GetSharedThisPtr()
Returns a shared pointer to the current object.
Definition: ExpList.h:956
std::shared_ptr< GlobalLinSys > GenGlobalLinSys(const GlobalLinSysKey &mkey, const std::shared_ptr< AssemblyMapCG > &locToGloMap)
This operation constructs the global linear system of type mkey.
Definition: ExpList.cpp:3071
std::shared_ptr< LocalRegions::ExpansionVector > m_exp
The list of local expansions.
Definition: ExpList.h:1118
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Definition: ExpList.h:1063
std::shared_ptr< DNekMat > GenGlobalMatrixFull(const GlobalLinSysKey &mkey, const std::shared_ptr< AssemblyMapCG > &locToGloMap)
Definition: ExpList.cpp:2919
LibUtilities::SessionReaderSharedPtr m_session
Session.
Definition: ExpList.h:1058
void MultiplyByInvMassMatrix(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.h:1725
void BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
This function elementally evaluates the backward transformation of the global spectral/hp element exp...
Definition: ExpList.h:1717
Describes a matrix with ordering defined by a local to global map.
bool LocToGloMapIsDefined() const
Returns true if a local to global map is defined.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< GlobalLinSys > GlobalLinSysSharedPtr
Pointer to a GlobalLinSys object.
Definition: GlobalLinSys.h:51
std::shared_ptr< GlobalMatrix > GlobalMatrixSharedPtr
Shared pointer to a GlobalMatrix object.
Definition: GlobalMatrix.h:84
std::map< GlobalMatrixKey, GlobalMatrixSharedPtr > GlobalMatrixMap
Mapping from global matrix keys to global matrices.
Definition: GlobalMatrix.h:86
std::map< StdRegions::ConstFactorType, Array< OneD, NekDouble > > VarFactorsMap
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
@ eLinearAdvectionDiffusionReactionGJP
Definition: StdRegions.hpp:105
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:430
std::map< StdRegions::VarCoeffType, VarCoeffEntry > VarCoeffMap
Definition: StdRegions.hpp:375
StdRegions::ConstFactorMap factors
static Array< OneD, NekDouble > NullNekDouble1DArray
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75
double NekDouble
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.hpp:292
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
STL namespace.