Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
ContField2D.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File ContField2D.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 // License for the specific language governing rights and limitations under
14 // Permission is hereby granted, free of charge, to any person obtaining a
15 // copy of this software and associated documentation files (the "Software"),
16 // to deal in the Software without restriction, including without limitation
17 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
18 // and/or sell copies of the Software, and to permit persons to whom the
19 // Software is furnished to do so, subject to the following conditions:
20 //
21 // The above copyright notice and this permission notice shall be included
22 // in all copies or substantial portions of the Software.
23 //
24 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30 // DEALINGS IN THE SOFTWARE.
31 //
32 // Description: Field definition for 2D domain with boundary conditions
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
38 
39 using namespace std;
40 
41 namespace Nektar
42 {
43  namespace MultiRegions
44  {
45  /**
46  * @class ContField2D
47  * The class #ContField2D is
48  * able to incorporate the boundary conditions imposed to the problem
49  * to be solved. Therefore, the class is equipped with three additional
50  * data members:
51  * - #m_bndCondExpansions
52  * - #m_bndTypes
53  * - #m_bndCondEquations
54  *
55  * The first data structure, #m_bndCondExpansions, contains the
56  * one-dimensional spectral/hp expansion on the boundary, #m_bndTypes
57  * stores information about the type of boundary condition on the
58  * different parts of the boundary while #m_bndCondEquations holds the
59  * equation of the imposed boundary conditions.
60  *
61  * Furthermore, in case of Dirichlet boundary conditions, this class is
62  * capable of lifting a known solution satisfying these boundary
63  * conditions. If we denote the unknown solution by
64  * \f$u^{\mathcal{H}}(\boldsymbol{x})\f$ and the known Dirichlet
65  * boundary conditions by \f$u^{\mathcal{D}}(\boldsymbol{x})\f$, the
66  * expansion then can be decomposed as
67  * \f[ u^{\delta}(\boldsymbol{x}_i)=u^{\mathcal{D}}(\boldsymbol{x}_i)+
68  * u^{\mathcal{H}}(\boldsymbol{x}_i)=\sum_{n=0}^{N^{\mathcal{D}}-1}
69  * \hat{u}_n^{\mathcal{D}}\Phi_n(\boldsymbol{x}_i)+
70  * \sum_{n={N^{\mathcal{D}}}}^{N_{\mathrm{dof}}-1}
71  * \hat{u}_n^{\mathcal{H}} \Phi_n(\boldsymbol{x}_i).\f]
72  * This lifting is accomplished by ordering the known global degrees of
73  * freedom, prescribed by the Dirichlet boundary conditions, first in
74  * the global array
75  * \f$\boldsymbol{\hat{u}}\f$, that is,
76  * \f[\boldsymbol{\hat{u}}=\left[ \begin{array}{c}
77  * \boldsymbol{\hat{u}}^{\mathcal{D}}\\
78  * \boldsymbol{\hat{u}}^{\mathcal{H}}
79  * \end{array} \right].\f]
80  * Such kind of expansions are also referred to as continuous fields.
81  * This class should be used when solving 2D problems using a standard
82  * Galerkin approach.
83  */
84 
85  /**
86  *
87  */
88  ContField2D::ContField2D():
90  m_locToGloMap(),
91  m_globalMat(),
92  m_globalLinSysManager(
93  boost::bind(&ContField2D::GenGlobalLinSys, this, _1),
94  std::string("GlobalLinSys"))
95  {
96  }
97 
98 
99  /**
100  * Given a mesh \a graph2D, containing information about the domain and
101  * the spectral/hp element expansion, this constructor fills the list
102  * of local expansions #m_exp with the proper expansions, calculates
103  * the total number of quadrature points \f$\boldsymbol{x}_i\f$ and
104  * local expansion coefficients \f$\hat{u}^e_n\f$ and allocates memory
105  * for the arrays #m_coeffs and #m_phys. Furthermore, it constructs the
106  * mapping array (contained in #m_locToGloMap) for the transformation
107  * between local elemental level and global level, it calculates the
108  * total number global expansion coefficients \f$\hat{u}_n\f$ and
109  * allocates memory for the array #m_contCoeffs. The constructor also
110  * discretises the boundary conditions, specified by the argument \a
111  * bcs, by expressing them in terms of the coefficient of the expansion
112  * on the boundary.
113  *
114  * @param graph2D A mesh, containing information about the domain
115  * and the spectral/hp element expansion.
116  * @param bcs The boundary conditions.
117  * @param variable An optional parameter to indicate for which
118  * variable the field should be constructed.
119  */
121  const SpatialDomains::MeshGraphSharedPtr &graph2D,
122  const std::string &variable,
123  const bool DeclareCoeffPhysArrays,
124  const bool CheckIfSingularSystem):
125  DisContField2D(pSession,graph2D,variable,false,DeclareCoeffPhysArrays),
126  m_globalMat(MemoryManager<GlobalMatrixMap>::AllocateSharedPtr()),
127  m_globalLinSysManager(
128  boost::bind(&ContField2D::GenGlobalLinSys, this, _1),
129  std::string("GlobalLinSys"))
130  {
135  CheckIfSingularSystem,
136  variable,
139 
140  if (m_session->DefinesCmdLineArgument("verbose"))
141  {
142  m_locToGloMap->PrintStats(std::cout, variable);
143  }
144  }
145 
146 
147  /**
148  * Given a mesh \a graph2D, containing information about the domain and
149  * the spectral/hp element expansion, this constructor fills the list
150  * of local expansions #m_exp with the proper expansions, calculates
151  * the total number of quadrature points \f$\boldsymbol{x}_i\f$ and
152  * local expansion coefficients \f$\hat{u}^e_n\f$ and allocates memory
153  * for the arrays #m_coeffs and #m_phys. Furthermore, it constructs the
154  * mapping array (contained in #m_locToGloMap) for the transformation
155  * between local elemental level and global level, it calculates the
156  * total number global expansion coefficients \f$\hat{u}_n\f$ and
157  * allocates memory for the array #m_coeffs. The constructor also
158  * discretises the boundary conditions, specified by the argument \a
159  * bcs, by expressing them in terms of the coefficient of the expansion
160  * on the boundary.
161  *
162  * @param In Existing ContField2D object used to provide the
163  * local to global mapping information and
164  * global solution type.
165  * @param graph2D A mesh, containing information about the domain
166  * and the spectral/hp element expansion.
167  * @param bcs The boundary conditions.
168  * @param bc_loc
169  */
171  const SpatialDomains::MeshGraphSharedPtr &graph2D,
172  const std::string &variable,
173  bool DeclareCoeffPhysArrays,
174  const bool CheckIfSingularSystem):
175  DisContField2D(In,graph2D,variable,false,DeclareCoeffPhysArrays),
176  m_globalMat (MemoryManager<GlobalMatrixMap>::AllocateSharedPtr()),
177  m_globalLinSysManager(
178  boost::bind(&ContField2D::GenGlobalLinSys, this, _1),
179  std::string("GlobalLinSys"))
180  {
181  if(!SameTypeOfBoundaryConditions(In) || CheckIfSingularSystem)
182  {
187  CheckIfSingularSystem,
188  variable,
191 
192  if (m_session->DefinesCmdLineArgument("verbose"))
193  {
194  m_locToGloMap->PrintStats(std::cout, variable);
195  }
196  }
197  else
198  {
200  }
201  }
202 
203 
204  /**
205  * Initialises the object as a copy of an existing ContField2D object.
206  * @param In Existing ContField2D object.
207  * @param DeclareCoeffPhysArrays bool to declare if \a m_phys
208  * and \a m_coeffs should be declared. Default is true
209  */
210  ContField2D::ContField2D(const ContField2D &In, bool DeclareCoeffPhysArrays):
211  DisContField2D(In,DeclareCoeffPhysArrays),
212  m_locToGloMap(In.m_locToGloMap),
213  m_globalMat(In.m_globalMat),
214  m_globalLinSysManager(In.m_globalLinSysManager)
215  {
216  }
217 
218 
219  /**
220  *
221  */
223  {
224  }
225 
226 
227  /**
228  * Given a function \f$f(\boldsymbol{x})\f$ defined at the quadrature
229  * points, this function determines the unknown global coefficients
230  * \f$\boldsymbol{\hat{u}}^{\mathcal{H}}\f$ employing a discrete
231  * Galerkin projection from physical space to coefficient
232  * space. The operation is evaluated by the function #GlobalSolve using
233  * the global mass matrix.
234  *
235  * The values of the function \f$f(\boldsymbol{x})\f$ evaluated at the
236  * quadrature points \f$\boldsymbol{x}_i\f$ should be contained in the
237  * variable #m_phys of the ExpList object \a Sin. The resulting global
238  * coefficients \f$\hat{u}_g\f$ are stored in the array #m_coeffs.
239  *
240  * @param Sin An ExpList, containing the discrete evaluation
241  * of \f$f(\boldsymbol{x})\f$ at the quadrature
242  * points in its array #m_phys.
243  */
245  Array<OneD, NekDouble> &outarray,
246  CoeffState coeffstate)
247 
248  {
249  // Inner product of forcing
250  int contNcoeffs = m_locToGloMap->GetNumGlobalCoeffs();
251  Array<OneD,NekDouble> wsp(contNcoeffs);
252  IProductWRTBase(inarray,wsp,eGlobal);
253 
254  // Solve the system
256 
257  if(coeffstate == eGlobal)
258  {
259  GlobalSolve(key,wsp,outarray);
260  }
261  else
262  {
263  Array<OneD,NekDouble> tmp(contNcoeffs,0.0);
264  GlobalSolve(key,wsp,tmp);
265  GlobalToLocal(tmp,outarray);
266  }
267  }
268 
269  /**
270  *
271  */
273  {
274  int gloNcoeffs = m_locToGloMap->GetNumGlobalCoeffs();
275  Array<OneD,NekDouble> tmp1(gloNcoeffs);
276  Array<OneD,NekDouble> tmp2(gloNcoeffs);
277 
278  IProductWRTBase(field,tmp1,eGlobal);
279  MultiplyByInvMassMatrix(tmp1,tmp2,eGlobal);
280  BwdTrans(tmp2,field,eGlobal);
281  }
282 
283 
284  /**
285  * Computes the matrix vector product
286  * @f$ \mathbf{y} = \mathbf{M}^{-1}\mathbf{x} @f$. If \a coeffstate == eGlobal
287  * is set then the elemental system is used directly. If not set, the
288  * global system is assembled, the system is solved, and mapped back to
289  * the local elemental system.
290  *
291  * @param inarray Input vector @f$\mathbf{x}@f$.
292  * @param outarray Output vector @f$\mathbf{y}@f$.
293  * @param coeffState Flag for using global system.
294  */
296  const Array<OneD, const NekDouble> &inarray,
297  Array<OneD, NekDouble> &outarray,
298  CoeffState coeffstate)
299 
300  {
302  int contNcoeffs = m_locToGloMap->GetNumGlobalCoeffs();
303 
304  if(coeffstate == eGlobal)
305  {
306  if(inarray.data() == outarray.data())
307  {
308  Array<OneD, NekDouble> tmp(contNcoeffs,0.0);
309  Vmath::Vcopy(contNcoeffs,inarray,1,tmp,1);
310  GlobalSolve(key,tmp,outarray);
311  }
312  else
313  {
314  GlobalSolve(key,inarray,outarray);
315  }
316  }
317  else
318  {
319  Array<OneD, NekDouble> globaltmp(contNcoeffs,0.0);
320 
321  if(inarray.data() == outarray.data())
322  {
323  Array<OneD,NekDouble> tmp(inarray.num_elements());
324  Vmath::Vcopy(inarray.num_elements(),inarray,1,tmp,1);
325  Assemble(tmp,outarray);
326  }
327  else
328  {
329  Assemble(inarray,outarray);
330  }
331 
332  GlobalSolve(key,outarray,globaltmp);
333  GlobalToLocal(globaltmp,outarray);
334  }
335  }
336 
337 
338  /**
339  * Consider the two dimensional Laplace equation,
340  * \f[\nabla\cdot\left(\boldsymbol{\sigma}\nabla
341  * u(\boldsymbol{x})\right) = f(\boldsymbol{x}),\f] supplemented with
342  * appropriate boundary conditions (which are contained in the data
343  * member #m_bndCondExpansions). In the equation above
344  * \f$\boldsymbol{\sigma}\f$ is the (symmetric positive definite)
345  * diffusion tensor:
346  * \f[ \sigma = \left[ \begin{array}{cc}
347  * \sigma_{00}(\boldsymbol{x},t) & \sigma_{01}(\boldsymbol{x},t) \\
348  * \sigma_{01}(\boldsymbol{x},t) & \sigma_{11}(\boldsymbol{x},t)
349  * \end{array} \right]. \f]
350  * Applying a \f$C^0\f$ continuous Galerkin discretisation, this
351  * equation leads to the following linear system:
352  * \f[\boldsymbol{L}
353  * \boldsymbol{\hat{u}}_g=\boldsymbol{\hat{f}}\f]
354  * where \f$\boldsymbol{L}\f$ is the Laplacian matrix. This function
355  * solves the system above for the global coefficients
356  * \f$\boldsymbol{\hat{u}}\f$ by a call to the function #GlobalSolve.
357  *
358  * The values of the function \f$f(\boldsymbol{x})\f$ evaluated at the
359  * quadrature points \f$\boldsymbol{x}_i\f$ should be contained in the
360  * variable #m_phys of the ExpList object \a Sin. The resulting global
361  * coefficients \f$\boldsymbol{\hat{u}}_g\f$ are stored in the array
362  * #m_coeffs.
363  *
364  * @param Sin An ExpList, containing the discrete evaluation
365  * of the forcing function \f$f(\boldsymbol{x})\f$
366  * at the quadrature points in its array #m_phys.
367  * @param variablecoeffs The (optional) parameter containing the
368  * coefficients evaluated at the quadrature
369  * points. It is an Array of (three) arrays which
370  * stores the laplacian coefficients in the
371  * following way
372  * \f[\mathrm{variablecoeffs} = \left[ \begin{array}{c}
373  * \left[\sigma_{00}(\boldsymbol{x_i},t)\right]_i \\
374  * \left[\sigma_{01}(\boldsymbol{x_i},t)\right]_i \\
375  * \left[\sigma_{11}(\boldsymbol{x_i},t)\right]_i
376  * \end{array}\right]
377  * \f]
378  * If this argument is not passed to the function, the following
379  * equation will be solved:
380  * \f[\nabla^2u(\boldsymbol{x}) = f(\boldsymbol{x}),\f]
381  *
382  * @param time The time-level at which the coefficients are
383  * evaluated
384  */
386  const Array<OneD, const NekDouble> &inarray,
387  Array<OneD, NekDouble> &outarray,
388  const Array<OneD, const NekDouble> &dirForcing,
389  const Array<OneD, Array<OneD,NekDouble> >& variablecoeffs,
390  NekDouble time,
391  CoeffState coeffstate)
392  {
393  // Inner product of forcing
394  int contNcoeffs = m_locToGloMap->GetNumGlobalCoeffs();
395  Array<OneD,NekDouble> wsp(contNcoeffs);
396  IProductWRTBase(inarray,wsp,eGlobal);
397  // Note -1.0 term necessary to invert forcing function to
398  // be consistent with matrix definition
399  Vmath::Neg(m_ncoeffs, wsp, 1);
400 
401  // Forcing function with weak boundary conditions
402  int i,j;
403  int bndcnt=0;
404  for(i = 0; i < m_bndCondExpansions.num_elements(); ++i)
405  {
406  if(m_bndConditions[i]->GetBoundaryConditionType() != SpatialDomains::eDirichlet)
407  {
408  for(j = 0; j < (m_bndCondExpansions[i])->GetNcoeffs(); j++)
409  {
410  wsp[m_locToGloMap
411  ->GetBndCondCoeffsToGlobalCoeffsMap(bndcnt++)]
412  += (m_bndCondExpansions[i]->GetCoeffs())[j];
413  }
414  }
415  else
416  {
417  bndcnt += m_bndCondExpansions[i]->GetNcoeffs();
418  }
419  }
420 
421  StdRegions::VarCoeffMap varcoeffs;
422  varcoeffs[StdRegions::eVarCoeffD00] = variablecoeffs[0];
423  varcoeffs[StdRegions::eVarCoeffD11] = variablecoeffs[3];
424  varcoeffs[StdRegions::eVarCoeffD22] = variablecoeffs[5];
426  factors[StdRegions::eFactorTime] = time;
427 
428  // Solve the system
430  varcoeffs);
431 
432  if(coeffstate == eGlobal)
433  {
434  GlobalSolve(key,wsp,outarray,dirForcing);
435  }
436  else
437  {
438  Array<OneD,NekDouble> tmp(contNcoeffs,0.0);
439  GlobalSolve(key,wsp,tmp,dirForcing);
440  GlobalToLocal(tmp,outarray);
441  }
442  }
443 
444 
445  /**
446  * Constructs the GlobalLinearSysKey for the linear advection operator
447  * with the supplied parameters, and computes the eigenvectors and
448  * eigenvalues of the associated matrix.
449  * @param ax Advection parameter, x.
450  * @param ay Advection parameter, y.
451  * @param Real Computed eigenvalues, real component.
452  * @param Imag Computed eigenvalues, imag component.
453  * @param Evecs Computed eigenvectors.
454  */
456  const NekDouble ay,
459  Array<OneD, NekDouble> &Evecs)
460  {
461  // Solve the system
465  vel[0] = vel_x;
466  vel[1] = vel_y;
467 
468  StdRegions::VarCoeffMap varcoeffs;
472  factors[StdRegions::eFactorTime] = 0.0;
474  factors,varcoeffs);
475 
477  Gmat->EigenSolve(Real,Imag,Evecs);
478  }
479 
480 
481 
482 
483  /**
484  * Given a linear system specified by the key \a key,
485  * \f[\boldsymbol{M}\boldsymbol{\hat{u}}_g=\boldsymbol{\hat{f}},\f]
486  * this function solves this linear system taking into account the
487  * boundary conditions specified in the data member
488  * #m_bndCondExpansions. Therefore, it adds an array
489  * \f$\boldsymbol{\hat{g}}\f$ which represents the non-zero surface
490  * integral resulting from the weak boundary conditions (e.g. Neumann
491  * boundary conditions) to the right hand side, that is,
492  * \f[\boldsymbol{M}\boldsymbol{\hat{u}}_g=\boldsymbol{\hat{f}}+
493  * \boldsymbol{\hat{g}}.\f]
494  * Furthermore, it lifts the known degrees of freedom which are
495  * prescribed by the Dirichlet boundary conditions. As these known
496  * coefficients \f$\boldsymbol{\hat{u}}^{\mathcal{D}}\f$ are numbered
497  * first in the global coefficient array \f$\boldsymbol{\hat{u}}_g\f$,
498  * the linear system can be decomposed as,
499  * \f[\left[\begin{array}{cc}
500  * \boldsymbol{M}^{\mathcal{DD}}&\boldsymbol{M}^{\mathcal{DH}}\\
501  * \boldsymbol{M}^{\mathcal{HD}}&\boldsymbol{M}^{\mathcal{HH}}
502  * \end{array}\right]
503  * \left[\begin{array}{c}
504  * \boldsymbol{\hat{u}}^{\mathcal{D}}\\
505  * \boldsymbol{\hat{u}}^{\mathcal{H}}
506  * \end{array}\right]=
507  * \left[\begin{array}{c}
508  * \boldsymbol{\hat{f}}^{\mathcal{D}}\\
509  * \boldsymbol{\hat{f}}^{\mathcal{H}}
510  * \end{array}\right]+
511  * \left[\begin{array}{c}
512  * \boldsymbol{\hat{g}}^{\mathcal{D}}\\
513  * \boldsymbol{\hat{g}}^{\mathcal{H}}
514  * \end{array}\right]
515  * \f]
516  * which will then be solved for the unknown coefficients
517  * \f$\boldsymbol{\hat{u}}^{\mathcal{H}}\f$ as,
518  * \f[
519  * \boldsymbol{M}^{\mathcal{HH}}\boldsymbol{\hat{u}}^{\mathcal{H}}=
520  * \boldsymbol{\hat{f}}^{\mathcal{H}}+
521  * \boldsymbol{\hat{g}}^{\mathcal{H}}-
522  * \boldsymbol{M}^{\mathcal{HD}}\boldsymbol{\hat{u}}^{\mathcal{D}}\f]
523  *
524  * @param mkey This key uniquely defines the linear system to
525  * be solved.
526  * @param Sin An ExpList, containing the discrete evaluation
527  * of the forcing function \f$f(\boldsymbol{x})\f$
528  * at the quadrature points in its array #m_phys.
529  * @param ScaleForcing An optional parameter with which the forcing
530  * vector \f$\boldsymbol{\hat{f}}\f$ should be
531  * multiplied.
532  * @note inout contains initial guess and final output.
533  */
535  const GlobalLinSysKey &key,
536  const Array<OneD, const NekDouble>& rhs,
537  Array<OneD, NekDouble>& inout,
538  const Array<OneD, const NekDouble>& dirForcing)
539  {
540  int NumDirBcs = m_locToGloMap->GetNumGlobalDirBndCoeffs();
541  int contNcoeffs = m_locToGloMap->GetNumGlobalCoeffs();
542 
543  // STEP 1: SET THE DIRICHLET DOFS TO THE RIGHT VALUE
544  // IN THE SOLUTION ARRAY
546 
547  // STEP 2: CALCULATE THE HOMOGENEOUS COEFFICIENTS
548  if(contNcoeffs - NumDirBcs > 0)
549  {
551  LinSys->Solve(rhs,inout,m_locToGloMap,dirForcing);
552  }
553  }
554 
555 
556  /**
557  * Returns the global matrix associated with the given GlobalMatrixKey.
558  * If the global matrix has not yet been constructed on this field,
559  * it is first constructed using GenGlobalMatrix().
560  * @param mkey Global matrix key.
561  * @returns Assocated global matrix.
562  */
564  const GlobalMatrixKey &mkey)
565  {
567  "To use method must have a AssemblyMap "
568  "attached to key");
569 
570  GlobalMatrixSharedPtr glo_matrix;
571  GlobalMatrixMap::iterator matrixIter = m_globalMat->find(mkey);
572 
573  if(matrixIter == m_globalMat->end())
574  {
575  glo_matrix = GenGlobalMatrix(mkey,m_locToGloMap);
576  (*m_globalMat)[mkey] = glo_matrix;
577  }
578  else
579  {
580  glo_matrix = matrixIter->second;
581  }
582 
583  return glo_matrix;
584  }
585 
586 
587  /**
588  * The function searches the map #m_globalLinSys to see if the
589  * global matrix has been created before. If not, it calls the function
590  * #GenGlobalLinSys to generate the requested global system.
591  *
592  * @param mkey This key uniquely defines the requested
593  * linear system.
594  */
596  const GlobalLinSysKey &mkey)
597  {
598  return m_globalLinSysManager[mkey];
599  }
600 
602  const GlobalLinSysKey &mkey)
603  {
605  "To use method must have a AssemblyMap "
606  "attached to key");
608  }
609 
610 
611  /**
612  *
613  */
615  const Array<OneD, const NekDouble> &inarray,
616  Array<OneD, NekDouble> &outarray,
617  CoeffState coeffstate)
618  {
619  BwdTrans(inarray,outarray,coeffstate);
620  }
621 
622 
623  /**
624  *
625  */
627  const Array<OneD, const NekDouble> &inarray,
628  Array<OneD, NekDouble> &outarray,
629  CoeffState coeffstate)
630  {
631  FwdTrans(inarray,outarray,coeffstate);
632  }
633 
635  {
636  int i,j;
637  int bndcnt=0;
638  int nDir = m_locToGloMap->GetNumGlobalDirBndCoeffs();
639 
640  // STEP 1: SET THE DIRICHLET DOFS TO THE RIGHT VALUE IN THE SOLUTION
641  // ARRAY
642  NekDouble sign;
643  const Array<OneD,const int> &bndMap =
644  m_locToGloMap->GetBndCondCoeffsToGlobalCoeffsMap();
645 
647  m_locToGloMap->GetNumGlobalBndCoeffs(), 0.0);
648 
649  // Fill in Dirichlet coefficients that are to be sent to
650  // other processors. This code block uses a
651  // tuple<int,int.NekDouble> which stores the local id of
652  // coefficent the global id of the data location and the
653  // inverse of the values of the data (arising from
654  // periodic boundary conditiosn)
655  map<int, vector<ExtraDirDof> > &extraDirDofs =
656  m_locToGloMap->GetExtraDirDofs();
657  map<int, vector<ExtraDirDof> >::iterator it;
658  for (it = extraDirDofs.begin(); it != extraDirDofs.end(); ++it)
659  {
660  for (i = 0; i < it->second.size(); ++i)
661  {
662  tmp[it->second.at(i).get<1>()] =
663  m_bndCondExpansions[it->first]->GetCoeffs()[
664  it->second.at(i).get<0>()]*it->second.at(i).get<2>();
665  }
666  }
667  m_locToGloMap->UniversalAssembleBnd(tmp);
668 
669  // Now fill in all other Dirichlet coefficients.
670  for(i = 0; i < m_bndCondExpansions.num_elements(); ++i)
671  {
672  if(m_bndConditions[i]->GetBoundaryConditionType() ==
674  {
675  const Array<OneD,const NekDouble>& coeffs =
676  m_bndCondExpansions[i]->GetCoeffs();
677  for(j = 0; j < (m_bndCondExpansions[i])->GetNcoeffs(); ++j)
678  {
679  sign = m_locToGloMap->GetBndCondCoeffsToGlobalCoeffsSign(
680  bndcnt);
681  tmp[bndMap[bndcnt++]] = sign * coeffs[j];
682  }
683  }
684  else
685  {
686  bndcnt += m_bndCondExpansions[i]->GetNcoeffs();
687  }
688  }
689 
690  Vmath::Vcopy(nDir, tmp, 1, outarray, 1);
691  }
692 
694  {
695  NekDouble sign;
696  int bndcnt = 0;
697  const Array<OneD,const int> &bndMap =
698  m_locToGloMap->GetBndCondCoeffsToGlobalCoeffsMap();
699 
700  Array<OneD, NekDouble> tmp(m_locToGloMap->GetNumGlobalCoeffs());
701  LocalToGlobal(m_coeffs,tmp);
702 
703  // Now fill in all other Dirichlet coefficients.
704  for(int i = 0; i < m_bndCondExpansions.num_elements(); ++i)
705  {
706  Array<OneD, NekDouble>& coeffs = m_bndCondExpansions[i]->UpdateCoeffs();
707 
708  for(int j = 0; j < (m_bndCondExpansions[i])->GetNcoeffs(); ++j)
709  {
710  sign = m_locToGloMap->GetBndCondCoeffsToGlobalCoeffsSign(bndcnt);
711  coeffs[j] = sign * tmp[bndMap[bndcnt++]];
712  }
713  }
714  }
715 
716  /**
717  * This operation is evaluated as:
718  * \f{tabbing}
719  * \hspace{1cm} \= Do \= $e=$ $1, N_{\mathrm{el}}$ \\
720  * > > Do \= $i=$ $0,N_m^e-1$ \\
721  * > > > $\boldsymbol{\hat{u}}^{e}[i] = \mbox{sign}[e][i] \cdot
722  * \boldsymbol{\hat{u}}_g[\mbox{map}[e][i]]$ \\
723  * > > continue \\
724  * > continue
725  * \f}
726  * where \a map\f$[e][i]\f$ is the mapping array and \a
727  * sign\f$[e][i]\f$ is an array of similar dimensions ensuring the
728  * correct modal connectivity between the different elements (both
729  * these arrays are contained in the data member #m_locToGloMap). This
730  * operation is equivalent to the scatter operation
731  * \f$\boldsymbol{\hat{u}}_l=\mathcal{A}\boldsymbol{\hat{u}}_g\f$,
732  * where \f$\mathcal{A}\f$ is the
733  * \f$N_{\mathrm{eof}}\times N_{\mathrm{dof}}\f$ permutation matrix.
734  *
735  * @note The array #m_coeffs should be filled with the global
736  * coefficients \f$\boldsymbol{\hat{u}}_g\f$ and that the resulting
737  * local coefficients \f$\boldsymbol{\hat{u}}_l\f$ will be stored in
738  * #m_coeffs.
739  */
741  {
742  m_locToGloMap->GlobalToLocal(m_coeffs,m_coeffs);
743  }
744 
745 
746 
747  /**
748  * This operation is evaluated as:
749  * \f{tabbing}
750  * \hspace{1cm} \= Do \= $e=$ $1, N_{\mathrm{el}}$ \\
751  * > > Do \= $i=$ $0,N_m^e-1$ \\
752  * > > > $\boldsymbol{\hat{u}}_g[\mbox{map}[e][i]] =
753  * \mbox{sign}[e][i] \cdot \boldsymbol{\hat{u}}^{e}[i]$\\
754  * > > continue\\
755  * > continue
756  * \f}
757  * where \a map\f$[e][i]\f$ is the mapping array and \a
758  * sign\f$[e][i]\f$ is an array of similar dimensions ensuring the
759  * correct modal connectivity between the different elements (both
760  * these arrays are contained in the data member #m_locToGloMap). This
761  * operation is equivalent to the gather operation
762  * \f$\boldsymbol{\hat{u}}_g=\mathcal{A}^{-1}\boldsymbol{\hat{u}}_l\f$,
763  * where \f$\mathcal{A}\f$ is the
764  * \f$N_{\mathrm{eof}}\times N_{\mathrm{dof}}\f$ permutation matrix.
765  *
766  * @note The array #m_coeffs should be filled with the local
767  * coefficients \f$\boldsymbol{\hat{u}}_l\f$ and that the
768  * resulting global coefficients \f$\boldsymbol{\hat{u}}_g\f$
769  * will be stored in #m_coeffs.
770  */
772  {
773  m_locToGloMap->LocalToGlobal(m_coeffs,m_coeffs);
774  }
775 
776  /**
777  *
778  */
780  const Array<OneD, const NekDouble> &inarray,
781  Array<OneD, NekDouble> &outarray,
782  CoeffState coeffstate)
783  {
784  MultiplyByInvMassMatrix(inarray,outarray,coeffstate);
785  }
786 
787 
788  /**
789  * Consider the two dimensional Helmholtz equation,
790  * \f[\nabla^2u(\boldsymbol{x})-\lambda u(\boldsymbol{x})
791  * = f(\boldsymbol{x}),\f] supplemented with appropriate boundary
792  * conditions (which are contained in the data member
793  * #m_bndCondExpansions). Applying a \f$C^0\f$ continuous Galerkin
794  * discretisation, this equation leads to the following linear system:
795  * \f[\left(\boldsymbol{L}+\lambda\boldsymbol{M}\right)
796  * \boldsymbol{\hat{u}}_g=\boldsymbol{\hat{f}}\f] where
797  * \f$\boldsymbol{L}\f$ and \f$\boldsymbol{M}\f$ are the Laplacian and
798  * mass matrix respectively. This function solves the system above for
799  * the global coefficients \f$\boldsymbol{\hat{u}}\f$ by a call to the
800  * function #GlobalSolve. It is assumed #m_coeff contains an
801  * initial estimate for the solution.
802  *
803  * The values of the function \f$f(\boldsymbol{x})\f$
804  * evaluated at the quadrature points \f$\boldsymbol{x}_i\f$
805  * should be contained in the variable #m_phys of the ExpList
806  * object \a inarray. The resulting global coefficients
807  * \f$\boldsymbol{\hat{u}}_g\f$ are stored in the array
808  * #m_contCoeffs or #m_coeffs depending on whether
809  * \a coeffstate is eGlobal or eLocal
810  *
811  * @param inarray An ExpList, containing the discrete evaluation
812  * of the forcing function \f$f(\boldsymbol{x})\f$
813  * at the quadrature points in its array #m_phys.
814  * @param factors The parameter \f$\lambda\f$ of the Helmholtz
815  * equation is specified through the factors map
816  */
818  const Array<OneD, const NekDouble> &inarray,
819  Array<OneD, NekDouble> &outarray,
820  const FlagList &flags,
821  const StdRegions::ConstFactorMap &factors,
822  const StdRegions::VarCoeffMap &varcoeff,
823  const Array<OneD, const NekDouble> &dirForcing)
824  {
825  //----------------------------------
826  // Setup RHS Inner product
827  //----------------------------------
828  // Inner product of forcing
829  int contNcoeffs = m_locToGloMap->GetNumGlobalCoeffs();
830  Array<OneD,NekDouble> wsp(contNcoeffs);
831  IProductWRTBase(inarray,wsp,eGlobal);
832  // Note -1.0 term necessary to invert forcing function to
833  // be consistent with matrix definition
834  Vmath::Neg(contNcoeffs, wsp, 1);
835 
836  // Fill weak boundary conditions
837  int i,j;
838  int bndcnt=0;
839  Array<OneD, NekDouble> gamma(contNcoeffs, 0.0);
840 
841  for(i = 0; i < m_bndCondExpansions.num_elements(); ++i)
842  {
843  if(m_bndConditions[i]->GetBoundaryConditionType() != SpatialDomains::eDirichlet)
844  {
845  for(j = 0; j < (m_bndCondExpansions[i])->GetNcoeffs(); j++)
846  {
847  gamma[m_locToGloMap
848  ->GetBndCondCoeffsToGlobalCoeffsMap(bndcnt++)]
849  += (m_bndCondExpansions[i]->GetCoeffs())[j];
850  }
851  }
852  else
853  {
854  bndcnt += m_bndCondExpansions[i]->GetNcoeffs();
855  }
856  }
857 
858  m_locToGloMap->UniversalAssemble(gamma);
859 
860  // Add weak boundary conditions to forcing
861  Vmath::Vadd(contNcoeffs, wsp, 1, gamma, 1, wsp, 1);
862 
864 
865  if(flags.isSet(eUseGlobal))
866  {
867  GlobalSolve(key,wsp,outarray,dirForcing);
868  }
869  else
870  {
871  Array<OneD,NekDouble> tmp(contNcoeffs);
872  LocalToGlobal(outarray,tmp);
873  GlobalSolve(key,wsp,tmp,dirForcing);
874  GlobalToLocal(tmp,outarray);
875  }
876  }
877 
878 
879  /**
880  * This is equivalent to the operation:
881  * \f[\boldsymbol{M\hat{u}}_g\f]
882  * where \f$\boldsymbol{M}\f$ is the global matrix of type specified by
883  * \a mkey. After scattering the global array \a inarray to local
884  * level, this operation is evaluated locally by the function
885  * ExpList#GeneralMatrixOp. The global result is then obtained by a
886  * global assembly procedure.
887  *
888  * @param mkey This key uniquely defines the type matrix
889  * required for the operation.
890  * @param inarray The vector \f$\boldsymbol{\hat{u}}_g\f$ of size
891  * \f$N_{\mathrm{dof}}\f$.
892  * @param outarray The resulting vector of size
893  * \f$N_{\mathrm{dof}}\f$.
894  */
896  const GlobalMatrixKey &gkey,
897  const Array<OneD,const NekDouble> &inarray,
898  Array<OneD, NekDouble> &outarray,
899  CoeffState coeffstate)
900  {
901  if(coeffstate == eGlobal)
902  {
903  bool doGlobalOp = m_globalOptParam->DoGlobalMatOp(
904  gkey.GetMatrixType());
905 
906  if(doGlobalOp)
907  {
909  mat->Multiply(inarray,outarray);
910  m_locToGloMap->UniversalAssemble(outarray);
911  }
912  else
913  {
915  Array<OneD,NekDouble> tmp2(tmp1+m_ncoeffs);
916  GlobalToLocal(inarray,tmp1);
917  GeneralMatrixOp_IterPerExp(gkey,tmp1,tmp2);
918  Assemble(tmp2,outarray);
919  }
920  }
921  else
922  {
923  GeneralMatrixOp_IterPerExp(gkey,inarray,outarray);
924  }
925  }
926 
927  /**
928  * First compute the inner product of forcing function with respect to
929  * base, and then solve the system with the linear advection operator.
930  * @param velocity Array of advection velocities in physical space
931  * @param inarray Forcing function.
932  * @param outarray Result.
933  * @param lambda reaction coefficient
934  * @param coeffstate State of Coefficients, Local or Global
935  * @param dirForcing Dirichlet Forcing.
936  */
937 
938  // could combine this with HelmholtzCG.
940  const Array<OneD, const NekDouble> &inarray,
941  Array<OneD, NekDouble> &outarray,
942  const NekDouble lambda,
943  CoeffState coeffstate,
944  const Array<OneD, const NekDouble>& dirForcing)
945  {
946  // Inner product of forcing
947  int contNcoeffs = m_locToGloMap->GetNumGlobalCoeffs();
948  Array<OneD,NekDouble> wsp(contNcoeffs);
949  IProductWRTBase(inarray,wsp,eGlobal);
950  // Note -1.0 term necessary to invert forcing function to
951  // be consistent with matrix definition
952  Vmath::Neg(contNcoeffs, wsp, 1);
953 
954  // Forcing function with weak boundary conditions
955  int i,j;
956  int bndcnt=0;
957  Array<OneD, NekDouble> gamma(contNcoeffs, 0.0);
958  for(i = 0; i < m_bndCondExpansions.num_elements(); ++i)
959  {
960  if(m_bndConditions[i]->GetBoundaryConditionType() != SpatialDomains::eDirichlet)
961  {
962  for(j = 0; j < (m_bndCondExpansions[i])->GetNcoeffs(); j++)
963  {
964  gamma[m_locToGloMap
965  ->GetBndCondCoeffsToGlobalCoeffsMap(bndcnt++)]
966  += (m_bndCondExpansions[i]->GetCoeffs())[j];
967  }
968  }
969  else
970  {
971  bndcnt += m_bndCondExpansions[i]->GetNcoeffs();
972  }
973  }
974  m_locToGloMap->UniversalAssemble(wsp);
975  // Add weak boundary conditions to forcing
976  Vmath::Vadd(contNcoeffs, wsp, 1, gamma, 1, wsp, 1);
977 
978  // Solve the system
980  factors[StdRegions::eFactorLambda] = lambda;
981  StdRegions::VarCoeffMap varcoeffs;
982  varcoeffs[StdRegions::eVarCoeffVelX] = velocity[0];
983  varcoeffs[StdRegions::eVarCoeffVelY] = velocity[1];
985 
986  if(coeffstate == eGlobal)
987  {
988  GlobalSolve(key,wsp,outarray,dirForcing);
989  }
990  else
991  {
992  Array<OneD,NekDouble> tmp(contNcoeffs,0.0);
993  GlobalSolve(key,wsp,tmp,dirForcing);
994  GlobalToLocal(tmp,outarray);
995  }
996  }
997 
998  /**
999  * First compute the inner product of forcing function with respect to
1000  * base, and then solve the system with the linear advection operator.
1001  * @param velocity Array of advection velocities in physical space
1002  * @param inarray Forcing function.
1003  * @param outarray Result.
1004  * @param lambda reaction coefficient
1005  * @param coeffstate State of Coefficients, Local or Global
1006  * @param dirForcing Dirichlet Forcing.
1007  */
1009  const Array<OneD, const NekDouble> &inarray,
1010  Array<OneD, NekDouble> &outarray,
1011  const NekDouble lambda,
1012  CoeffState coeffstate,
1013  const Array<OneD, const NekDouble>& dirForcing)
1014  {
1015  // Inner product of forcing
1016  int contNcoeffs = m_locToGloMap->GetNumGlobalCoeffs();
1017  Array<OneD,NekDouble> wsp(contNcoeffs);
1018  IProductWRTBase(inarray,wsp,eGlobal);
1019 
1020  // Solve the system
1022  factors[StdRegions::eFactorLambda] = lambda;
1023  StdRegions::VarCoeffMap varcoeffs;
1024  varcoeffs[StdRegions::eVarCoeffVelX] = velocity[0];
1025  varcoeffs[StdRegions::eVarCoeffVelY] = velocity[1];
1027 
1028  if(coeffstate == eGlobal)
1029  {
1030  GlobalSolve(key,wsp,outarray,dirForcing);
1031  }
1032  else
1033  {
1034  Array<OneD,NekDouble> tmp(contNcoeffs,0.0);
1035  GlobalSolve(key,wsp,tmp,dirForcing);
1036  GlobalToLocal(tmp,outarray);
1037  }
1038  }
1039 
1040 
1041  /**
1042  *
1043  */
1046  {
1047  return GetBndConditions();
1048  }
1049 
1050 
1051  /**
1052  * Reset the GlobalLinSys Manager
1053  */
1055  {
1056  m_globalLinSysManager.ClearManager("GlobalLinSys");
1057  }
1058 
1059  } // end of namespace
1060 } //end of namespace
virtual void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
Template method virtual forwarder for FwdTrans().
boost::shared_ptr< GlobalMatrix > GenGlobalMatrix(const GlobalMatrixKey &mkey, const boost::shared_ptr< AssemblyMapCG > &locToGloMap)
Generates a global matrix from the given key and map.
Definition: ExpList.cpp:897
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, CoeffState coeffstate=eLocal)
Solves the two-dimensional Laplace equation, subject to the boundary conditions specified.
void GeneralMatrixOp_IterPerExp(const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:829
AssemblyMapCGSharedPtr m_locToGloMap
(A shared pointer to) the object which contains all the required information for the transformation f...
Definition: ContField2D.h:182
void FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
Performs the global forward transformation of a function , subject to the boundary conditions specifi...
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
#define sign(a, b)
return the sign(b)*a
Definition: Polylib.cpp:22
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
Definition: ExpList.h:1001
void v_LinearAdvectionReactionSolve(const Array< OneD, Array< OneD, NekDouble > > &velocity, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const NekDouble lambda, CoeffState coeffstate=eLocal, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > & GetBndConditions()
Returns the boundary conditions.
Definition: ContField2D.h:536
virtual void v_LocalToGlobal(void)
Gathers the global coefficients from the local coefficients .
virtual void v_ImposeDirichletConditions(Array< OneD, NekDouble > &outarray)
Impose the Dirichlet Boundary Conditions on outarray.
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
std::map< GlobalMatrixKey, GlobalMatrixSharedPtr > GlobalMatrixMap
Mapping from global matrix keys to global matrices.
Definition: GlobalMatrix.h:91
GlobalLinSysSharedPtr GenGlobalLinSys(const GlobalLinSysKey &mkey)
boost::shared_ptr< DNekMat > GenGlobalMatrixFull(const GlobalLinSysKey &mkey, const boost::shared_ptr< AssemblyMapCG > &locToGloMap)
Definition: ExpList.cpp:1034
void BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
Performs the backward transformation of the spectral/hp element expansion.
Definition: ContField2D.h:500
STL namespace.
PeriodicMap m_periodicEdges
A map which identifies pairs of periodic edges.
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:251
PeriodicMap m_periodicVerts
A map which identifies groups of periodic vertices.
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:939
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
LibUtilities::NekManager< GlobalLinSysKey, GlobalLinSys > m_globalLinSysManager
A manager which collects all the global linear systems being assembled, such that they should be cons...
Definition: ContField2D.h:192
boost::shared_ptr< GlobalMatrix > GlobalMatrixSharedPtr
Shared pointer to a GlobalMatrix object.
Definition: GlobalMatrix.h:89
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
Definition: MeshPartition.h:51
Global coefficients.
bool SameTypeOfBoundaryConditions(const DisContField2D &In)
bool isSet(const FlagType &key) const
Array< OneD, MultiRegions::ExpListSharedPtr > m_bndCondExpansions
An object which contains the discretised boundary conditions.
This class is the abstraction of a global continuous two- dimensional spectral/hp element expansion w...
Definition: ContField2D.h:56
std::map< StdRegions::VarCoeffType, Array< OneD, NekDouble > > VarCoeffMap
Definition: StdRegions.hpp:226
boost::shared_ptr< GlobalLinSys > GenGlobalLinSys(const GlobalLinSysKey &mkey, const boost::shared_ptr< AssemblyMapCG > &locToGloMap)
This operation constructs the global linear system of type mkey.
Definition: ExpList.cpp:1175
virtual void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
Template method virtual forwarder for FwdTrans().
virtual void v_SmoothField(Array< OneD, NekDouble > &field)
Template method virtual forwarded for SmoothField().
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Definition: ExpList.h:917
GlobalMatrixMapShPtr m_globalMat
(A shared pointer to) a list which collects all the global matrices being assembled, such that they should be constructed only once.
Definition: ContField2D.h:187
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.
GlobalMatrixSharedPtr GetGlobalMatrix(const GlobalMatrixKey &mkey)
Returns the global matrix specified by mkey.
virtual ~ContField2D()
The default destructor.
virtual void v_GeneralMatrixOp(const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
Calculates the result of the multiplication of a global matrix of type specified by mkey with a vecto...
LibUtilities::SessionReaderSharedPtr m_session
Session.
Definition: ExpList.h:910
ContField2D()
The default constructor.
Definition: ContField2D.cpp:88
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:382
virtual void v_GlobalToLocal(void)
Scatters from the global coefficients to the local coefficients .
void Assemble()
Assembles the global coefficients from the local coefficients .
Definition: ContField2D.h:392
double NekDouble
Defines a list of flags.
Array< OneD, SpatialDomains::BoundaryConditionShPtr > m_bndConditions
An array which contains the information about the boundary condition on the different boundary region...
virtual const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > & v_GetBndConditions()
Template method virtual forwarder for GetBndConditions().
Describe a linear system.
void IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
Calculates the inner product of a function with respect to all global expansion modes ...
Definition: ContField2D.h:454
StdRegions::MatrixType GetMatrixType() const
Return the matrix type.
Describes a matrix with ordering defined by a local to global map.
void LocalToGlobal(void)
Put the coefficients into global ordering using m_coeffs.
Definition: ExpList.h:1853
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator
bool LocToGloMapIsDefined() const
Returns true if a local to global map is defined.
virtual void v_MultiplyByInvMassMatrix(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
Template method virtual forwarder for MultiplyByInvMassMatrix().
void MultiplyByInvMassMatrix(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
Multiply a solution by the inverse mass matrix.
virtual void v_HelmSolve(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const FlagList &flags, const StdRegions::ConstFactorMap &factors, const StdRegions::VarCoeffMap &varcoeff, const Array< OneD, const NekDouble > &dirForcing)
Solves the two-dimensional Helmholtz equation, subject to the boundary conditions specified...
int GetNcoeffs(void) const
Returns the total number of local degrees of freedom .
Definition: ExpList.h:1406
boost::shared_ptr< GlobalLinSys > GlobalLinSysSharedPtr
Pointer to a GlobalLinSys object.
Definition: GlobalLinSys.h:52
GlobalLinSysSharedPtr GetGlobalLinSys(const GlobalLinSysKey &mkey)
Returns the linear system specified by the key mkey.
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.
virtual void v_ClearGlobalLinSysManager(void)
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:218
boost::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:442
void GlobalToLocal(void)
Put the coefficients into local ordering and place in m_coeffs.
Definition: ExpList.h:1858
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
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:285
virtual void v_LinearAdvectionDiffusionReactionSolve(const Array< OneD, Array< OneD, NekDouble > > &velocity, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const NekDouble lambda, CoeffState coeffstate=eLocal, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)