Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties 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 LibUtilities::SessionReaderSharedPtr &pSession,
122  const SpatialDomains::MeshGraphSharedPtr &graph2D,
123  const std::string &variable,
124  const bool DeclareCoeffPhysArrays,
125  const bool CheckIfSingularSystem,
126  const Collections::ImplementationType ImpType):
127  DisContField2D(pSession,graph2D,variable,false,DeclareCoeffPhysArrays,ImpType),
128  m_globalMat(MemoryManager<GlobalMatrixMap>::AllocateSharedPtr()),
129  m_globalLinSysManager(
130  boost::bind(&ContField2D::GenGlobalLinSys, this, _1),
131  std::string("GlobalLinSys"))
132  {
137  CheckIfSingularSystem,
138  variable,
141 
142  if (m_session->DefinesCmdLineArgument("verbose"))
143  {
144  m_locToGloMap->PrintStats(std::cout, variable);
145  }
146  }
147 
148 
149  /**
150  * Given a mesh \a graph2D, containing information about the domain and
151  * the spectral/hp element expansion, this constructor fills the list
152  * of local expansions #m_exp with the proper expansions, calculates
153  * the total number of quadrature points \f$\boldsymbol{x}_i\f$ and
154  * local expansion coefficients \f$\hat{u}^e_n\f$ and allocates memory
155  * for the arrays #m_coeffs and #m_phys. Furthermore, it constructs the
156  * mapping array (contained in #m_locToGloMap) for the transformation
157  * between local elemental level and global level, it calculates the
158  * total number global expansion coefficients \f$\hat{u}_n\f$ and
159  * allocates memory for the array #m_coeffs. The constructor also
160  * discretises the boundary conditions, specified by the argument \a
161  * bcs, by expressing them in terms of the coefficient of the expansion
162  * on the boundary.
163  *
164  * @param In Existing ContField2D object used to provide the
165  * local to global mapping information and
166  * global solution type.
167  * @param graph2D A mesh, containing information about the domain
168  * and the spectral/hp element expansion.
169  * @param bcs The boundary conditions.
170  * @param bc_loc
171  */
173  const SpatialDomains::MeshGraphSharedPtr &graph2D,
174  const std::string &variable,
175  bool DeclareCoeffPhysArrays,
176  const bool CheckIfSingularSystem):
177  DisContField2D(In,graph2D,variable,false,DeclareCoeffPhysArrays),
178  m_globalMat (MemoryManager<GlobalMatrixMap>::AllocateSharedPtr()),
179  m_globalLinSysManager(
180  boost::bind(&ContField2D::GenGlobalLinSys, this, _1),
181  std::string("GlobalLinSys"))
182  {
183  if(!SameTypeOfBoundaryConditions(In) || CheckIfSingularSystem)
184  {
189  CheckIfSingularSystem,
190  variable,
193 
194  if (m_session->DefinesCmdLineArgument("verbose"))
195  {
196  m_locToGloMap->PrintStats(std::cout, variable);
197  }
198  }
199  else
200  {
202  }
203  }
204 
205 
206  /**
207  * Initialises the object as a copy of an existing ContField2D object.
208  * @param In Existing ContField2D object.
209  * @param DeclareCoeffPhysArrays bool to declare if \a m_phys
210  * and \a m_coeffs should be declared. Default is true
211  */
212  ContField2D::ContField2D(const ContField2D &In, bool DeclareCoeffPhysArrays):
213  DisContField2D(In,DeclareCoeffPhysArrays),
214  m_locToGloMap(In.m_locToGloMap),
215  m_globalMat(In.m_globalMat),
216  m_globalLinSysManager(In.m_globalLinSysManager)
217  {
218  }
219 
220 
221  /**
222  *
223  */
225  {
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 #m_phys of the ExpList object \a Sin. The resulting global
240  * coefficients \f$\hat{u}_g\f$ are stored in the array #m_coeffs.
241  *
242  * @param Sin An ExpList, containing the discrete evaluation
243  * of \f$f(\boldsymbol{x})\f$ at the quadrature
244  * points in its array #m_phys.
245  */
247  Array<OneD, NekDouble> &outarray,
248  CoeffState coeffstate)
249 
250  {
251  // Inner product of forcing
252  int contNcoeffs = m_locToGloMap->GetNumGlobalCoeffs();
253  Array<OneD,NekDouble> wsp(contNcoeffs);
254  IProductWRTBase(inarray,wsp,eGlobal);
255 
256  // Solve the system
258 
259  if(coeffstate == eGlobal)
260  {
261  GlobalSolve(key,wsp,outarray);
262  }
263  else
264  {
265  Array<OneD,NekDouble> tmp(contNcoeffs,0.0);
266  GlobalSolve(key,wsp,tmp);
267  GlobalToLocal(tmp,outarray);
268  }
269  }
270 
271  /**
272  *
273  */
275  {
276  int gloNcoeffs = m_locToGloMap->GetNumGlobalCoeffs();
277  Array<OneD,NekDouble> tmp1(gloNcoeffs);
278  Array<OneD,NekDouble> tmp2(gloNcoeffs);
279 
280  IProductWRTBase(field,tmp1,eGlobal);
281  MultiplyByInvMassMatrix(tmp1,tmp2,eGlobal);
282  BwdTrans(tmp2,field,eGlobal);
283  }
284 
285 
286  /**
287  * Computes the matrix vector product
288  * @f$ \mathbf{y} = \mathbf{M}^{-1}\mathbf{x} @f$. If \a coeffstate == eGlobal
289  * is set then the elemental system is used directly. If not set, the
290  * global system is assembled, the system is solved, and mapped back to
291  * the local elemental system.
292  *
293  * @param inarray Input vector @f$\mathbf{x}@f$.
294  * @param outarray Output vector @f$\mathbf{y}@f$.
295  * @param coeffState Flag for using global system.
296  */
298  const Array<OneD, const NekDouble> &inarray,
299  Array<OneD, NekDouble> &outarray,
300  CoeffState coeffstate)
301 
302  {
304  int contNcoeffs = m_locToGloMap->GetNumGlobalCoeffs();
305 
306  if(coeffstate == eGlobal)
307  {
308  if(inarray.data() == outarray.data())
309  {
310  Array<OneD, NekDouble> tmp(contNcoeffs,0.0);
311  Vmath::Vcopy(contNcoeffs,inarray,1,tmp,1);
312  GlobalSolve(key,tmp,outarray);
313  }
314  else
315  {
316  GlobalSolve(key,inarray,outarray);
317  }
318  }
319  else
320  {
321  Array<OneD, NekDouble> globaltmp(contNcoeffs,0.0);
322 
323  if(inarray.data() == outarray.data())
324  {
325  Array<OneD,NekDouble> tmp(inarray.num_elements());
326  Vmath::Vcopy(inarray.num_elements(),inarray,1,tmp,1);
327  Assemble(tmp,outarray);
328  }
329  else
330  {
331  Assemble(inarray,outarray);
332  }
333 
334  GlobalSolve(key,outarray,globaltmp);
335  GlobalToLocal(globaltmp,outarray);
336  }
337  }
338 
339 
340  /**
341  * Consider the two dimensional Laplace equation,
342  * \f[\nabla\cdot\left(\boldsymbol{\sigma}\nabla
343  * u(\boldsymbol{x})\right) = f(\boldsymbol{x}),\f] supplemented with
344  * appropriate boundary conditions (which are contained in the data
345  * member #m_bndCondExpansions). In the equation above
346  * \f$\boldsymbol{\sigma}\f$ is the (symmetric positive definite)
347  * diffusion tensor:
348  * \f[ \sigma = \left[ \begin{array}{cc}
349  * \sigma_{00}(\boldsymbol{x},t) & \sigma_{01}(\boldsymbol{x},t) \\
350  * \sigma_{01}(\boldsymbol{x},t) & \sigma_{11}(\boldsymbol{x},t)
351  * \end{array} \right]. \f]
352  * Applying a \f$C^0\f$ continuous Galerkin discretisation, this
353  * equation leads to the following linear system:
354  * \f[\boldsymbol{L}
355  * \boldsymbol{\hat{u}}_g=\boldsymbol{\hat{f}}\f]
356  * where \f$\boldsymbol{L}\f$ is the Laplacian matrix. This function
357  * solves the system above for the global coefficients
358  * \f$\boldsymbol{\hat{u}}\f$ by a call to the function #GlobalSolve.
359  *
360  * The values of the function \f$f(\boldsymbol{x})\f$ evaluated at the
361  * quadrature points \f$\boldsymbol{x}_i\f$ should be contained in the
362  * variable #m_phys of the ExpList object \a Sin. The resulting global
363  * coefficients \f$\boldsymbol{\hat{u}}_g\f$ are stored in the array
364  * #m_coeffs.
365  *
366  * @param Sin An ExpList, containing the discrete evaluation
367  * of the forcing function \f$f(\boldsymbol{x})\f$
368  * at the quadrature points in its array #m_phys.
369  * @param variablecoeffs The (optional) parameter containing the
370  * coefficients evaluated at the quadrature
371  * points. It is an Array of (three) arrays which
372  * stores the laplacian coefficients in the
373  * following way
374  * \f[\mathrm{variablecoeffs} = \left[ \begin{array}{c}
375  * \left[\sigma_{00}(\boldsymbol{x_i},t)\right]_i \\
376  * \left[\sigma_{01}(\boldsymbol{x_i},t)\right]_i \\
377  * \left[\sigma_{11}(\boldsymbol{x_i},t)\right]_i
378  * \end{array}\right]
379  * \f]
380  * If this argument is not passed to the function, the following
381  * equation will be solved:
382  * \f[\nabla^2u(\boldsymbol{x}) = f(\boldsymbol{x}),\f]
383  *
384  * @param time The time-level at which the coefficients are
385  * evaluated
386  */
388  const Array<OneD, const NekDouble> &inarray,
389  Array<OneD, NekDouble> &outarray,
390  const Array<OneD, const NekDouble> &dirForcing,
391  const Array<OneD, Array<OneD,NekDouble> >& variablecoeffs,
392  NekDouble time,
393  CoeffState coeffstate)
394  {
395  // Inner product of forcing
396  int contNcoeffs = m_locToGloMap->GetNumGlobalCoeffs();
397  Array<OneD,NekDouble> wsp(contNcoeffs);
398  IProductWRTBase(inarray,wsp,eGlobal);
399  // Note -1.0 term necessary to invert forcing function to
400  // be consistent with matrix definition
401  Vmath::Neg(m_ncoeffs, wsp, 1);
402 
403  // Forcing function with weak boundary conditions
404  int i,j;
405  int bndcnt=0;
406  for(i = 0; i < m_bndCondExpansions.num_elements(); ++i)
407  {
408  if(m_bndConditions[i]->GetBoundaryConditionType() != SpatialDomains::eDirichlet)
409  {
410  for(j = 0; j < (m_bndCondExpansions[i])->GetNcoeffs(); j++)
411  {
412  wsp[m_locToGloMap
413  ->GetBndCondCoeffsToGlobalCoeffsMap(bndcnt++)]
414  += (m_bndCondExpansions[i]->GetCoeffs())[j];
415  }
416  }
417  else
418  {
419  bndcnt += m_bndCondExpansions[i]->GetNcoeffs();
420  }
421  }
422 
423  StdRegions::VarCoeffMap varcoeffs;
424  varcoeffs[StdRegions::eVarCoeffD00] = variablecoeffs[0];
425  varcoeffs[StdRegions::eVarCoeffD11] = variablecoeffs[3];
426  varcoeffs[StdRegions::eVarCoeffD22] = variablecoeffs[5];
428  factors[StdRegions::eFactorTime] = time;
429 
430  // Solve the system
432  varcoeffs);
433 
434  if(coeffstate == eGlobal)
435  {
436  GlobalSolve(key,wsp,outarray,dirForcing);
437  }
438  else
439  {
440  Array<OneD,NekDouble> tmp(contNcoeffs,0.0);
441  GlobalSolve(key,wsp,tmp,dirForcing);
442  GlobalToLocal(tmp,outarray);
443  }
444  }
445 
446 
447  /**
448  * Constructs the GlobalLinearSysKey for the linear advection operator
449  * with the supplied parameters, and computes the eigenvectors and
450  * eigenvalues of the associated matrix.
451  * @param ax Advection parameter, x.
452  * @param ay Advection parameter, y.
453  * @param Real Computed eigenvalues, real component.
454  * @param Imag Computed eigenvalues, imag component.
455  * @param Evecs Computed eigenvectors.
456  */
458  const NekDouble ay,
461  Array<OneD, NekDouble> &Evecs)
462  {
463  // Solve the system
467  vel[0] = vel_x;
468  vel[1] = vel_y;
469 
470  StdRegions::VarCoeffMap varcoeffs;
474  factors[StdRegions::eFactorTime] = 0.0;
476  factors,varcoeffs);
477 
479  Gmat->EigenSolve(Real,Imag,Evecs);
480  }
481 
482 
483 
484 
485  /**
486  * Given a linear system specified by the key \a key,
487  * \f[\boldsymbol{M}\boldsymbol{\hat{u}}_g=\boldsymbol{\hat{f}},\f]
488  * this function solves this linear system taking into account the
489  * boundary conditions specified in the data member
490  * #m_bndCondExpansions. Therefore, it adds an array
491  * \f$\boldsymbol{\hat{g}}\f$ which represents the non-zero surface
492  * integral resulting from the weak boundary conditions (e.g. Neumann
493  * boundary conditions) to the right hand side, that is,
494  * \f[\boldsymbol{M}\boldsymbol{\hat{u}}_g=\boldsymbol{\hat{f}}+
495  * \boldsymbol{\hat{g}}.\f]
496  * Furthermore, it lifts the known degrees of freedom which are
497  * prescribed by the Dirichlet boundary conditions. As these known
498  * coefficients \f$\boldsymbol{\hat{u}}^{\mathcal{D}}\f$ are numbered
499  * first in the global coefficient array \f$\boldsymbol{\hat{u}}_g\f$,
500  * the linear system can be decomposed as,
501  * \f[\left[\begin{array}{cc}
502  * \boldsymbol{M}^{\mathcal{DD}}&\boldsymbol{M}^{\mathcal{DH}}\\
503  * \boldsymbol{M}^{\mathcal{HD}}&\boldsymbol{M}^{\mathcal{HH}}
504  * \end{array}\right]
505  * \left[\begin{array}{c}
506  * \boldsymbol{\hat{u}}^{\mathcal{D}}\\
507  * \boldsymbol{\hat{u}}^{\mathcal{H}}
508  * \end{array}\right]=
509  * \left[\begin{array}{c}
510  * \boldsymbol{\hat{f}}^{\mathcal{D}}\\
511  * \boldsymbol{\hat{f}}^{\mathcal{H}}
512  * \end{array}\right]+
513  * \left[\begin{array}{c}
514  * \boldsymbol{\hat{g}}^{\mathcal{D}}\\
515  * \boldsymbol{\hat{g}}^{\mathcal{H}}
516  * \end{array}\right]
517  * \f]
518  * which will then be solved for the unknown coefficients
519  * \f$\boldsymbol{\hat{u}}^{\mathcal{H}}\f$ as,
520  * \f[
521  * \boldsymbol{M}^{\mathcal{HH}}\boldsymbol{\hat{u}}^{\mathcal{H}}=
522  * \boldsymbol{\hat{f}}^{\mathcal{H}}+
523  * \boldsymbol{\hat{g}}^{\mathcal{H}}-
524  * \boldsymbol{M}^{\mathcal{HD}}\boldsymbol{\hat{u}}^{\mathcal{D}}\f]
525  *
526  * @param mkey This key uniquely defines the linear system to
527  * be solved.
528  * @param Sin An ExpList, containing the discrete evaluation
529  * of the forcing function \f$f(\boldsymbol{x})\f$
530  * at the quadrature points in its array #m_phys.
531  * @param ScaleForcing An optional parameter with which the forcing
532  * vector \f$\boldsymbol{\hat{f}}\f$ should be
533  * multiplied.
534  * @note inout contains initial guess and final output.
535  */
537  const GlobalLinSysKey &key,
538  const Array<OneD, const NekDouble>& rhs,
539  Array<OneD, NekDouble>& inout,
540  const Array<OneD, const NekDouble>& dirForcing)
541  {
542  int NumDirBcs = m_locToGloMap->GetNumGlobalDirBndCoeffs();
543  int contNcoeffs = m_locToGloMap->GetNumGlobalCoeffs();
544 
545  // STEP 1: SET THE DIRICHLET DOFS TO THE RIGHT VALUE
546  // IN THE SOLUTION ARRAY
548 
549  // STEP 2: CALCULATE THE HOMOGENEOUS COEFFICIENTS
550  if(contNcoeffs - NumDirBcs > 0)
551  {
553  LinSys->Solve(rhs,inout,m_locToGloMap,dirForcing);
554  }
555  }
556 
557 
558  /**
559  * Returns the global matrix associated with the given GlobalMatrixKey.
560  * If the global matrix has not yet been constructed on this field,
561  * it is first constructed using GenGlobalMatrix().
562  * @param mkey Global matrix key.
563  * @returns Assocated global matrix.
564  */
566  const GlobalMatrixKey &mkey)
567  {
569  "To use method must have a AssemblyMap "
570  "attached to key");
571 
572  GlobalMatrixSharedPtr glo_matrix;
573  GlobalMatrixMap::iterator matrixIter = m_globalMat->find(mkey);
574 
575  if(matrixIter == m_globalMat->end())
576  {
577  glo_matrix = GenGlobalMatrix(mkey,m_locToGloMap);
578  (*m_globalMat)[mkey] = glo_matrix;
579  }
580  else
581  {
582  glo_matrix = matrixIter->second;
583  }
584 
585  return glo_matrix;
586  }
587 
588 
589  /**
590  * The function searches the map #m_globalLinSys to see if the
591  * global matrix has been created before. If not, it calls the function
592  * #GenGlobalLinSys to generate the requested global system.
593  *
594  * @param mkey This key uniquely defines the requested
595  * linear system.
596  */
598  const GlobalLinSysKey &mkey)
599  {
600  return m_globalLinSysManager[mkey];
601  }
602 
604  const GlobalLinSysKey &mkey)
605  {
607  "To use method must have a AssemblyMap "
608  "attached to key");
610  }
611 
612 
613  /**
614  *
615  */
617  const Array<OneD, const NekDouble> &inarray,
618  Array<OneD, NekDouble> &outarray,
619  CoeffState coeffstate)
620  {
621  BwdTrans(inarray,outarray,coeffstate);
622  }
623 
624 
625  /**
626  *
627  */
629  const Array<OneD, const NekDouble> &inarray,
630  Array<OneD, NekDouble> &outarray,
631  CoeffState coeffstate)
632  {
633  FwdTrans(inarray,outarray,coeffstate);
634  }
635 
637  {
638  int i,j;
639  int bndcnt=0;
640  int nDir = m_locToGloMap->GetNumGlobalDirBndCoeffs();
641 
642  // STEP 1: SET THE DIRICHLET DOFS TO THE RIGHT VALUE IN THE SOLUTION
643  // ARRAY
644  NekDouble sign;
645  const Array<OneD,const int> &bndMap =
646  m_locToGloMap->GetBndCondCoeffsToGlobalCoeffsMap();
647 
649  m_locToGloMap->GetNumGlobalBndCoeffs(), 0.0);
650 
651  // Fill in Dirichlet coefficients that are to be sent to
652  // other processors. This code block uses a
653  // tuple<int,int.NekDouble> which stores the local id of
654  // coefficent the global id of the data location and the
655  // inverse of the values of the data (arising from
656  // periodic boundary conditiosn)
657  map<int, vector<ExtraDirDof> > &extraDirDofs =
658  m_locToGloMap->GetExtraDirDofs();
659  map<int, vector<ExtraDirDof> >::iterator it;
660  for (it = extraDirDofs.begin(); it != extraDirDofs.end(); ++it)
661  {
662  for (i = 0; i < it->second.size(); ++i)
663  {
664  tmp[it->second.at(i).get<1>()] =
665  m_bndCondExpansions[it->first]->GetCoeffs()[
666  it->second.at(i).get<0>()]*it->second.at(i).get<2>();
667  }
668  }
669  m_locToGloMap->UniversalAssembleBnd(tmp);
670 
671  // Now fill in all other Dirichlet coefficients.
672  for(i = 0; i < m_bndCondExpansions.num_elements(); ++i)
673  {
674  if(m_bndConditions[i]->GetBoundaryConditionType() ==
676  {
677  const Array<OneD,const NekDouble>& coeffs =
678  m_bndCondExpansions[i]->GetCoeffs();
679  for(j = 0; j < (m_bndCondExpansions[i])->GetNcoeffs(); ++j)
680  {
681  sign = m_locToGloMap->GetBndCondCoeffsToGlobalCoeffsSign(
682  bndcnt);
683  tmp[bndMap[bndcnt++]] = sign * coeffs[j];
684  }
685  }
686  else
687  {
688  bndcnt += m_bndCondExpansions[i]->GetNcoeffs();
689  }
690  }
691 
692  Vmath::Vcopy(nDir, tmp, 1, outarray, 1);
693  }
694 
696  {
697  NekDouble sign;
698  int bndcnt = 0;
699  const Array<OneD,const int> &bndMap =
700  m_locToGloMap->GetBndCondCoeffsToGlobalCoeffsMap();
701 
702  Array<OneD, NekDouble> tmp(m_locToGloMap->GetNumGlobalCoeffs());
703  LocalToGlobal(m_coeffs,tmp);
704 
705  // Now fill in all other Dirichlet coefficients.
706  for(int i = 0; i < m_bndCondExpansions.num_elements(); ++i)
707  {
708  Array<OneD, NekDouble>& coeffs = m_bndCondExpansions[i]->UpdateCoeffs();
709 
710  for(int j = 0; j < (m_bndCondExpansions[i])->GetNcoeffs(); ++j)
711  {
712  sign = m_locToGloMap->GetBndCondCoeffsToGlobalCoeffsSign(bndcnt);
713  coeffs[j] = sign * tmp[bndMap[bndcnt++]];
714  }
715  }
716  }
717 
718 
720  {
721  NekDouble sign;
722  int bndcnt = 0;
723  const Array<OneD,const int> &bndMap =
724  m_locToGloMap->GetBndCondCoeffsToGlobalCoeffsMap();
725 
726  Array<OneD, NekDouble> tmp(m_locToGloMap->GetNumGlobalCoeffs());
727  LocalToGlobal(m_coeffs,tmp,false);
728 
729  ASSERTL1(nreg < m_bndCondExpansions.num_elements(),
730  "nreg is out or range since this many boundary "
731  "regions to not exist");
732 
733  // Now fill in all other Dirichlet coefficients.
734  Array<OneD, NekDouble>& coeffs = m_bndCondExpansions[nreg]->UpdateCoeffs();
735 
736  for(int j = 0; j < nreg; ++j)
737  {
738  bndcnt += m_bndCondExpansions[j]->GetNcoeffs();
739  }
740 
741  for(int j = 0; j < (m_bndCondExpansions[nreg])->GetNcoeffs(); ++j)
742  {
743  sign = m_locToGloMap->GetBndCondCoeffsToGlobalCoeffsSign(bndcnt);
744  coeffs[j] = sign * tmp[bndMap[bndcnt++]];
745  }
746  }
747 
748 
749  /**
750  * This operation is evaluated as:
751  * \f{tabbing}
752  * \hspace{1cm} \= Do \= $e=$ $1, N_{\mathrm{el}}$ \\
753  * > > Do \= $i=$ $0,N_m^e-1$ \\
754  * > > > $\boldsymbol{\hat{u}}^{e}[i] = \mbox{sign}[e][i] \cdot
755  * \boldsymbol{\hat{u}}_g[\mbox{map}[e][i]]$ \\
756  * > > continue \\
757  * > continue
758  * \f}
759  * where \a map\f$[e][i]\f$ is the mapping array and \a
760  * sign\f$[e][i]\f$ is an array of similar dimensions ensuring the
761  * correct modal connectivity between the different elements (both
762  * these arrays are contained in the data member #m_locToGloMap). This
763  * operation is equivalent to the scatter operation
764  * \f$\boldsymbol{\hat{u}}_l=\mathcal{A}\boldsymbol{\hat{u}}_g\f$,
765  * where \f$\mathcal{A}\f$ is the
766  * \f$N_{\mathrm{eof}}\times N_{\mathrm{dof}}\f$ permutation matrix.
767  *
768  * @note The array #m_coeffs should be filled with the global
769  * coefficients \f$\boldsymbol{\hat{u}}_g\f$ and that the resulting
770  * local coefficients \f$\boldsymbol{\hat{u}}_l\f$ will be stored in
771  * #m_coeffs.
772  */
774  const Array<OneD, const NekDouble> &inarray,
775  Array<OneD,NekDouble> &outarray)
776  {
777  m_locToGloMap->GlobalToLocal(inarray, outarray);
778  }
779 
780 
782  {
783  m_locToGloMap->GlobalToLocal(m_coeffs,m_coeffs);
784  }
785 
786 
787  /**
788  * This operation is evaluated as:
789  * \f{tabbing}
790  * \hspace{1cm} \= Do \= $e=$ $1, N_{\mathrm{el}}$ \\
791  * > > Do \= $i=$ $0,N_m^e-1$ \\
792  * > > > $\boldsymbol{\hat{u}}_g[\mbox{map}[e][i]] =
793  * \mbox{sign}[e][i] \cdot \boldsymbol{\hat{u}}^{e}[i]$\\
794  * > > continue\\
795  * > continue
796  * \f}
797  * where \a map\f$[e][i]\f$ is the mapping array and \a
798  * sign\f$[e][i]\f$ is an array of similar dimensions ensuring the
799  * correct modal connectivity between the different elements (both
800  * these arrays are contained in the data member #m_locToGloMap). This
801  * operation is equivalent to the gather operation
802  * \f$\boldsymbol{\hat{u}}_g=\mathcal{A}^{-1}\boldsymbol{\hat{u}}_l\f$,
803  * where \f$\mathcal{A}\f$ is the
804  * \f$N_{\mathrm{eof}}\times N_{\mathrm{dof}}\f$ permutation matrix.
805  *
806  * @note The array #m_coeffs should be filled with the local
807  * coefficients \f$\boldsymbol{\hat{u}}_l\f$ and that
808  * the resulting global coefficients
809  * \f$\boldsymbol{\hat{u}}_g\f$ will be stored in
810  * #m_coeffs. Also if useComm is set to false then no
811  * communication call will be made to check if all
812  * values are consistent over processors
813  */
814 
816  const Array<OneD, const NekDouble> &inarray,
817  Array<OneD,NekDouble> &outarray,
818  bool useComm)
819  {
820  m_locToGloMap->LocalToGlobal(inarray, outarray, useComm);
821  }
822 
823 
824  void ContField2D::v_LocalToGlobal(bool useComm)
825 
826  {
827  m_locToGloMap->LocalToGlobal(m_coeffs,m_coeffs, useComm);
828  }
829 
830  /**
831  *
832  */
834  const Array<OneD, const NekDouble> &inarray,
835  Array<OneD, NekDouble> &outarray,
836  CoeffState coeffstate)
837  {
838  MultiplyByInvMassMatrix(inarray,outarray,coeffstate);
839  }
840 
841 
842  /**
843  * Consider the two dimensional Helmholtz equation,
844  * \f[\nabla^2u(\boldsymbol{x})-\lambda u(\boldsymbol{x})
845  * = f(\boldsymbol{x}),\f] supplemented with appropriate boundary
846  * conditions (which are contained in the data member
847  * #m_bndCondExpansions). Applying a \f$C^0\f$ continuous Galerkin
848  * discretisation, this equation leads to the following linear system:
849  * \f[\left(\boldsymbol{L}+\lambda\boldsymbol{M}\right)
850  * \boldsymbol{\hat{u}}_g=\boldsymbol{\hat{f}}\f] where
851  * \f$\boldsymbol{L}\f$ and \f$\boldsymbol{M}\f$ are the Laplacian and
852  * mass matrix respectively. This function solves the system above for
853  * the global coefficients \f$\boldsymbol{\hat{u}}\f$ by a call to the
854  * function #GlobalSolve. It is assumed #m_coeff contains an
855  * initial estimate for the solution.
856  *
857  * The values of the function \f$f(\boldsymbol{x})\f$
858  * evaluated at the quadrature points \f$\boldsymbol{x}_i\f$
859  * should be contained in the variable #m_phys of the ExpList
860  * object \a inarray. The resulting global coefficients
861  * \f$\boldsymbol{\hat{u}}_g\f$ are stored in the array
862  * #m_contCoeffs or #m_coeffs depending on whether
863  * \a coeffstate is eGlobal or eLocal
864  *
865  * @param inarray An ExpList, containing the discrete evaluation
866  * of the forcing function \f$f(\boldsymbol{x})\f$
867  * at the quadrature points in its array #m_phys.
868  * @param factors The parameter \f$\lambda\f$ of the Helmholtz
869  * equation is specified through the factors map
870  */
872  const Array<OneD, const NekDouble> &inarray,
873  Array<OneD, NekDouble> &outarray,
874  const FlagList &flags,
875  const StdRegions::ConstFactorMap &factors,
876  const StdRegions::VarCoeffMap &varcoeff,
877  const Array<OneD, const NekDouble> &dirForcing,
878  const bool PhysSpaceForcing)
879 
880  {
881  //----------------------------------
882  // Setup RHS Inner product
883  //----------------------------------
884  // Inner product of forcing
885  int contNcoeffs = m_locToGloMap->GetNumGlobalCoeffs();
886  Array<OneD,NekDouble> wsp(contNcoeffs);
887  if(PhysSpaceForcing)
888  {
889  IProductWRTBase(inarray,wsp,eGlobal);
890  }
891  else
892  {
893  Assemble(inarray,wsp);
894  }
895  // Note -1.0 term necessary to invert forcing function to
896  // be consistent with matrix definition
897  Vmath::Neg(contNcoeffs, wsp, 1);
898 
899  // Fill weak boundary conditions
900  int i,j;
901  int bndcnt=0;
902  Array<OneD, NekDouble> gamma(contNcoeffs, 0.0);
903 
904  for(i = 0; i < m_bndCondExpansions.num_elements(); ++i)
905  {
906  if(m_bndConditions[i]->GetBoundaryConditionType() != SpatialDomains::eDirichlet)
907  {
908  for(j = 0; j < (m_bndCondExpansions[i])->GetNcoeffs(); j++)
909  {
910  gamma[m_locToGloMap
911  ->GetBndCondCoeffsToGlobalCoeffsMap(bndcnt++)]
912  += (m_bndCondExpansions[i]->GetCoeffs())[j];
913  }
914  }
915  else
916  {
917  bndcnt += m_bndCondExpansions[i]->GetNcoeffs();
918  }
919  }
920 
921  m_locToGloMap->UniversalAssemble(gamma);
922 
923  // Add weak boundary conditions to forcing
924  Vmath::Vadd(contNcoeffs, wsp, 1, gamma, 1, wsp, 1);
925 
927 
928  if(flags.isSet(eUseGlobal))
929  {
930  GlobalSolve(key,wsp,outarray,dirForcing);
931  }
932  else
933  {
934  Array<OneD,NekDouble> tmp(contNcoeffs);
935  LocalToGlobal(outarray,tmp);
936  GlobalSolve(key,wsp,tmp,dirForcing);
937  GlobalToLocal(tmp,outarray);
938  }
939  }
940 
941 
942  /**
943  * This is equivalent to the operation:
944  * \f[\boldsymbol{M\hat{u}}_g\f]
945  * where \f$\boldsymbol{M}\f$ is the global matrix of type specified by
946  * \a mkey. After scattering the global array \a inarray to local
947  * level, this operation is evaluated locally by the function
948  * ExpList#GeneralMatrixOp. The global result is then obtained by a
949  * global assembly procedure.
950  *
951  * @param mkey This key uniquely defines the type matrix
952  * required for the operation.
953  * @param inarray The vector \f$\boldsymbol{\hat{u}}_g\f$ of size
954  * \f$N_{\mathrm{dof}}\f$.
955  * @param outarray The resulting vector of size
956  * \f$N_{\mathrm{dof}}\f$.
957  */
959  const GlobalMatrixKey &gkey,
960  const Array<OneD,const NekDouble> &inarray,
961  Array<OneD, NekDouble> &outarray,
962  CoeffState coeffstate)
963  {
964  if(coeffstate == eGlobal)
965  {
966  bool doGlobalOp = m_globalOptParam->DoGlobalMatOp(
967  gkey.GetMatrixType());
968 
969  if(doGlobalOp)
970  {
972  mat->Multiply(inarray,outarray);
973  m_locToGloMap->UniversalAssemble(outarray);
974  }
975  else
976  {
978  Array<OneD,NekDouble> tmp2(tmp1+m_ncoeffs);
979  GlobalToLocal(inarray,tmp1);
980  GeneralMatrixOp_IterPerExp(gkey,tmp1,tmp2);
981  Assemble(tmp2,outarray);
982  }
983  }
984  else
985  {
986  GeneralMatrixOp_IterPerExp(gkey,inarray,outarray);
987  }
988  }
989 
990  /**
991  * First compute the inner product of forcing function with respect to
992  * base, and then solve the system with the linear advection operator.
993  * @param velocity Array of advection velocities in physical space
994  * @param inarray Forcing function.
995  * @param outarray Result.
996  * @param lambda reaction coefficient
997  * @param coeffstate State of Coefficients, Local or Global
998  * @param dirForcing Dirichlet Forcing.
999  */
1000 
1001  // could combine this with HelmholtzCG.
1003  const Array<OneD, const NekDouble> &inarray,
1004  Array<OneD, NekDouble> &outarray,
1005  const NekDouble lambda,
1006  CoeffState coeffstate,
1007  const Array<OneD, const NekDouble>& dirForcing)
1008  {
1009  // Inner product of forcing
1010  int contNcoeffs = m_locToGloMap->GetNumGlobalCoeffs();
1011  Array<OneD,NekDouble> wsp(contNcoeffs);
1012  IProductWRTBase(inarray,wsp,eGlobal);
1013  // Note -1.0 term necessary to invert forcing function to
1014  // be consistent with matrix definition
1015  Vmath::Neg(contNcoeffs, wsp, 1);
1016 
1017  // Forcing function with weak boundary conditions
1018  int i,j;
1019  int bndcnt=0;
1020  Array<OneD, NekDouble> gamma(contNcoeffs, 0.0);
1021  for(i = 0; i < m_bndCondExpansions.num_elements(); ++i)
1022  {
1023  if(m_bndConditions[i]->GetBoundaryConditionType() != SpatialDomains::eDirichlet)
1024  {
1025  for(j = 0; j < (m_bndCondExpansions[i])->GetNcoeffs(); j++)
1026  {
1027  gamma[m_locToGloMap
1028  ->GetBndCondCoeffsToGlobalCoeffsMap(bndcnt++)]
1029  += (m_bndCondExpansions[i]->GetCoeffs())[j];
1030  }
1031  }
1032  else
1033  {
1034  bndcnt += m_bndCondExpansions[i]->GetNcoeffs();
1035  }
1036  }
1037  m_locToGloMap->UniversalAssemble(gamma);
1038  // Add weak boundary conditions to forcing
1039  Vmath::Vadd(contNcoeffs, wsp, 1, gamma, 1, wsp, 1);
1040 
1041  // Solve the system
1043  factors[StdRegions::eFactorLambda] = lambda;
1044  StdRegions::VarCoeffMap varcoeffs;
1045  varcoeffs[StdRegions::eVarCoeffVelX] = velocity[0];
1046  varcoeffs[StdRegions::eVarCoeffVelY] = velocity[1];
1048 
1049  if(coeffstate == eGlobal)
1050  {
1051  GlobalSolve(key,wsp,outarray,dirForcing);
1052  }
1053  else
1054  {
1055  Array<OneD,NekDouble> tmp(contNcoeffs,0.0);
1056  GlobalSolve(key,wsp,tmp,dirForcing);
1057  GlobalToLocal(tmp,outarray);
1058  }
1059  }
1060 
1061  /**
1062  * First compute the inner product of forcing function with respect to
1063  * base, and then solve the system with the linear advection operator.
1064  * @param velocity Array of advection velocities in physical space
1065  * @param inarray Forcing function.
1066  * @param outarray Result.
1067  * @param lambda reaction coefficient
1068  * @param coeffstate State of Coefficients, Local or Global
1069  * @param dirForcing Dirichlet Forcing.
1070  */
1072  const Array<OneD, Array<OneD, NekDouble> > &velocity,
1073  const Array<OneD, const NekDouble> &inarray,
1074  Array<OneD, NekDouble> &outarray,
1075  const NekDouble lambda,
1076  CoeffState coeffstate,
1077  const Array<OneD, const NekDouble>& dirForcing)
1078  {
1079  // Inner product of forcing
1080  int contNcoeffs = m_locToGloMap->GetNumGlobalCoeffs();
1081  Array<OneD,NekDouble> wsp(contNcoeffs);
1082  IProductWRTBase(inarray,wsp,eGlobal);
1083 
1084  // Solve the system
1086  factors[StdRegions::eFactorLambda] = lambda;
1087  StdRegions::VarCoeffMap varcoeffs;
1088  varcoeffs[StdRegions::eVarCoeffVelX] = velocity[0];
1089  varcoeffs[StdRegions::eVarCoeffVelY] = velocity[1];
1091 
1092  if(coeffstate == eGlobal)
1093  {
1094  GlobalSolve(key,wsp,outarray,dirForcing);
1095  }
1096  else
1097  {
1098  Array<OneD,NekDouble> tmp(contNcoeffs,0.0);
1099  GlobalSolve(key,wsp,tmp,dirForcing);
1100  GlobalToLocal(tmp,outarray);
1101  }
1102  }
1103 
1104 
1105  /**
1106  *
1107  */
1110  {
1111  return GetBndConditions();
1112  }
1113 
1114 
1115  /**
1116  * Reset the GlobalLinSys Manager
1117  */
1119  {
1120  m_globalLinSysManager.ClearManager("GlobalLinSys");
1121  }
1122 
1123  } // end of namespace
1124 } //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:997
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 LocalToGlobal(bool useComm=true)
Gathers the global coefficients from the local coefficients .
Definition: ExpList.h:1959
void GeneralMatrixOp_IterPerExp(const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:929
AssemblyMapCGSharedPtr m_locToGloMap
(A shared pointer to) the object which contains all the required information for the transformation f...
Definition: ContField2D.h:167
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:27
NekOptimize::GlobalOptParamSharedPtr m_globalOptParam
Definition: ExpList.h:1052
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:460
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:1134
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:424
STL namespace.
PeriodicMap m_periodicEdges
A map which identifies pairs of periodic edges.
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:252
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:998
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:177
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.
virtual void v_LocalToGlobal(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool useComm)
Gathers the global coefficients from the local 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:227
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:1275
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:976
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:172
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:969
ContField2D()
The default constructor.
Definition: ContField2D.cpp:88
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:396
virtual void v_GlobalToLocal(void)
void Assemble()
Assembles the global coefficients from the local coefficients .
Definition: ContField2D.h:316
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:378
StdRegions::MatrixType GetMatrixType() const
Return the matrix type.
Describes a matrix with ordering defined by a local to global map.
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs typedef NekMatrix< LhsDataType, StandardMatrixTag >::iterator iterator
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, const bool PhysSpaceForcing)
Solves the two-dimensional Helmholtz equation, subject to the boundary conditions specified...
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.
int GetNcoeffs(void) const
Returns the total number of local degrees of freedom .
Definition: ExpList.h:1485
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:228
boost::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:442
void GlobalToLocal(void)
Scatters from the global coefficients to the local coefficients .
Definition: ExpList.h:1972
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1061
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:299
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)