Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
ContField1D.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File ContField1D.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 1D domain with boundary conditions
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
38 
39 namespace Nektar
40 {
41  namespace MultiRegions
42  {
43  /**
44  * @class ContField1D
45  * As opposed to the class #ContExpList1D, the class #ContField1D is
46  * able to incorporate the boundary conditions imposed to the problem
47  * to be solved. Therefore, the class is equipped with three additional
48  * data members:
49  * - #m_bndCondExpansions
50  * - #m_bndTypes
51  * - #m_bndCondEquations
52  *
53  * The first data structure, #m_bndCondExpansions,
54  * contains the point Expansion on the boundary, #m_bndTypes
55  * stores information about the type of boundary condition on the
56  * different parts of the boundary while #m_bndCondEquations holds the
57  * equation of the imposed boundary conditions.
58  *
59  * Furthermore, in case of Dirichlet boundary conditions,
60  * this class is capable of lifting a known solution satisfying these
61  * boundary conditions. If we denote the unknown solution by
62  * \f$u^{\mathcal{H}}(\boldsymbol{x})\f$ and the known Dirichlet
63  * boundary conditions by \f$u^{\mathcal{D}}(\boldsymbol{x})\f$, the
64  * expansion then can be decomposed as
65  * \f[ u^{\delta}(\boldsymbol{x}_i)=u^{\mathcal{D}}(\boldsymbol{x}_i)+
66  * u^{\mathcal{H}}(\boldsymbol{x}_i)=\sum_{n=0}^{N^{\mathcal{D}}-1}
67  * \hat{u}_n^{\mathcal{D}}\Phi_n(\boldsymbol{x}_i)+
68  * \sum_{n={N^{\mathcal{D}}}}^{N_{\mathrm{dof}}
69  * -1}\hat{u}_n^{\mathcal{H}}
70  * \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 \f$\boldsymbol{\hat{u}}\f$, that is,
74  * \f[\boldsymbol{\hat{u}}=\left[ \begin{array}{c}
75  * \boldsymbol{\hat{u}}^{\mathcal{D}}\\
76  * \boldsymbol{\hat{u}}^{\mathcal{H}}
77  * \end{array} \right].\f]
78  * Such kind of expansions are also referred to as continuoous fields.
79  * This class should be used when solving 2D problems using a standard
80  * Galerkin approach.
81  */
82 
83  /**
84  * Constructs an empty 1D continuous field.
85  */
88  m_locToGloMap(),
89  m_globalLinSysManager(
90  boost::bind(&ContField1D::GenGlobalLinSys, this, _1),
91  std::string("GlobalLinSys"))
92  {
93  }
94 
95 
96  /**
97  * Given a mesh \a graph1D, containing information about the domain and
98  * the spectral/hp element expansion, this constructor fills the list
99  * of local expansions #m_exp with the proper expansions, calculates
100  * the total number of quadrature points \f$\boldsymbol{x}_i\f$ and
101  * local expansion coefficients \f$\hat{u}^e_n\f$ and allocates
102  * memory for the arrays #m_coeffs and #m_phys. Furthermore, it
103  * constructs the mapping array (contained in #m_locToGloMap) for the
104  * transformation between local elemental level and global level, it
105  * calculates the total number global expansion coefficients
106  * \f$\hat{u}_n\f$.
107  * The constructor also discretises the boundary conditions, specified
108  * by the argument \a bcs, by expressing them in terms of the
109  * coefficient of the expansion on the boundary.
110  *
111  * @param graph1D A 1D mesh, containing information about the
112  * domain and the spectral/hp element expansion.
113  * @param bcs The boundary conditions.
114  * @param variable An optional parameter to indicate for which
115  * variable the field should be constructed.
116  */
118  const SpatialDomains::MeshGraphSharedPtr &graph1D,
119  const std::string &variable):
120  DisContField1D(pSession,graph1D,variable,false),
121  m_locToGloMap(),
122  m_globalLinSysManager(
123  boost::bind(&ContField1D::GenGlobalLinSys, this, _1),
124  std::string("GlobalLinSys"))
125  {
126  SpatialDomains::BoundaryConditions bcs(pSession, graph1D);
127 
132  false,
133  variable,
135  }
136 
137 
138  /**
139  * Constructs a 1D continuous field as a copy of an existing field
140  * including all boundary conditions.
141  * @param In Existing continuous field to duplicate.
142  */
144  DisContField1D(In),
145  m_locToGloMap(In.m_locToGloMap),
146  m_globalLinSysManager(
147  boost::bind(&ContField1D::GenGlobalLinSys, this, _1),
148  std::string("GlobalLinSys"))
149  {
150  }
151 
152  /**
153  * Constructs a 1D continuous field as a copy of an existing explist 1D field
154  * and adding all the boundary conditions.
155  * @param In Existing explist1D field .
156  */
158  DisContField1D(In),
159  m_locToGloMap(),
160  m_globalLinSysManager(
161  boost::bind(&ContField1D::GenGlobalLinSys, this, _1),
162  std::string("GlobalLinSys"))
163  {
165  ::AllocateSharedPtr(pSession, m_ncoeffs, In);
166 
167  }
168 
169  /**
170  *
171  */
173  {
174  }
175 
176 
177  /**
178  * Given a function \f$f(x)\f$ defined at the quadrature
179  * points, this function determines the unknown global coefficients
180  * \f$\boldsymbol{\hat{u}}^{\mathcal{H}}\f$ employing a discrete
181  * Galerkin projection from physical space to coefficient
182  * space. The operation is evaluated by the function #GlobalSolve using
183  * the global mass matrix.
184  *
185  * The values of the function \f$f(x)\f$ evaluated at the
186  * quadrature points \f$x_i\f$ should be contained in the
187  * variable #m_phys of the ExpList object \a Sin. The resulting global
188  * coefficients \f$\hat{u}_g\f$ are stored in the array #m_coeffs.
189  *
190  * @param inarray Discrete \f$f(x)\f$.
191  * @param outarray Storage for result.
192  * @param coeffstate
193  */
195  Array<OneD, NekDouble> &outarray,
196  CoeffState coeffstate)
197  {
198  // Inner product of forcing
200  IProductWRTBase(inarray,wsp,eGlobal);
201 
202  // Solve the system
204 
205  GlobalSolve(key,wsp,outarray);
206  if(coeffstate == eLocal)
207  {
208  GlobalToLocal(outarray,outarray);
209  }
210  }
211 
212 
213  /**
214  * Given the coefficients of an expansion, this function evaluates the
215  * spectral/hp expansion \f$u^{\delta}(x)\f$ at the quadrature
216  * points \f$x_i\f$. This operation is evaluated locally by the
217  * function ExpList#BwdTrans.
218  *
219  * The coefficients of the expansion should be contained in the
220  * variable #m_coeffs of the ExpList object \a In. The resulting
221  * physical values at the quadrature points \f$u^{\delta}(x_i)\f$ are
222  * stored in the array #m_phys.
223  *
224  * @param In An ExpList, containing the local
225  * coefficients \f$\hat{u}_n^e\f$ in its array
226  * #m_coeffs.
227  */
229  const Array<OneD, const NekDouble> &inarray,
230  Array<OneD, NekDouble> &outarray,
231  CoeffState coeffstate)
232  {
233  Array<OneD, NekDouble> tmpinarray;
234  if(coeffstate != eLocal)
235  {
236  tmpinarray = Array<OneD, NekDouble>(inarray);
237  GlobalToLocal(inarray,tmpinarray);
238  }
239  else
240  {
241  tmpinarray = inarray;
242  }
243 
244  BwdTrans_IterPerExp(tmpinarray,outarray);
245  }
246 
247 
248  /**
249  *
250  */
252  const Array<OneD, const NekDouble> &inarray,
253  Array<OneD, NekDouble> &outarray,
254  CoeffState coeffstate)
255  {
257  if(coeffstate == eGlobal)
258  {
259  if(inarray.data() == outarray.data())
260  {
261  Array<OneD, NekDouble> tmp(inarray);
262  GlobalSolve(key,tmp,outarray);
263  }
264  else
265  {
266  GlobalSolve(key,inarray,outarray);
267  }
268  }
269  else
270  {
271  Array<OneD, NekDouble> globaltmp(m_ncoeffs,0.0);
272 
273  if(inarray.data() == outarray.data())
274  {
275  Array<OneD,NekDouble> tmp(inarray);
276  Assemble(tmp,outarray);
277  }
278  else
279  {
280  Assemble(inarray,outarray);
281  }
282 
283  GlobalSolve(key,outarray,globaltmp);
284  GlobalToLocal(globaltmp,outarray);
285  }
286  }
287 
288 
289  /**
290  * Given a linear system specified by the key \a key,
291  * \f[\boldsymbol{M}\boldsymbol{\hat{u}}_g=\boldsymbol{\hat{f}},\f]
292  * this function solves this linear system taking into account the
293  * boundary conditions specified in the data member
294  * #m_bndCondExpansions.
295  * Therefore, it adds an array \f$\boldsymbol{\hat{g}}\f$ which
296  * represents the non-zero surface integral resulting from the weak
297  * boundary conditions (e.g. Neumann boundary conditions) to the right
298  * hand side, that is,
299  * \f[\boldsymbol{M}\boldsymbol{\hat{u}}_g=\boldsymbol{\hat{f}}+
300  * \boldsymbol{\hat{g}}.\f]
301  * Furthermore, it lifts the known degrees of freedom which are
302  * prescribed by the Dirichlet boundary conditions. As these known
303  * coefficients \f$\boldsymbol{\hat{u}}^{\mathcal{D}}\f$ are numbered
304  * first in the global coefficient array \f$\boldsymbol{\hat{u}}_g\f$,
305  * the linear system can be decomposed as,
306  * \f[\left[\begin{array}{cc}
307  * \boldsymbol{M}^{\mathcal{DD}}&\boldsymbol{M}^{\mathcal{DH}}\\
308  * \boldsymbol{M}^{\mathcal{HD}}&\boldsymbol{M}^{\mathcal{HH}}
309  * \end{array}\right]
310  * \left[\begin{array}{c}
311  * \boldsymbol{\hat{u}}^{\mathcal{D}}\\
312  * \boldsymbol{\hat{u}}^{\mathcal{H}}
313  * \end{array}\right]=
314  * \left[\begin{array}{c}
315  * \boldsymbol{\hat{f}}^{\mathcal{D}}\\
316  * \boldsymbol{\hat{f}}^{\mathcal{H}}
317  * \end{array}\right]+
318  * \left[\begin{array}{c}
319  * \boldsymbol{\hat{g}}^{\mathcal{D}}\\
320  * \boldsymbol{\hat{g}}^{\mathcal{H}}
321  * \end{array}\right]
322  * \f]
323  * which will then be solved for the unknown coefficients
324  * \f$\boldsymbol{\hat{u}}^{\mathcal{H}}\f$ as,
325  * \f[
326  * \boldsymbol{M}^{\mathcal{HH}}\boldsymbol{\hat{u}}^{\mathcal{H}}
327  * = \boldsymbol{\hat{f}}^{\mathcal{H}}
328  * +\boldsymbol{\hat{g}}^{\mathcal{H}}
329  * -\boldsymbol{M}^{\mathcal{HD}}\boldsymbol{\hat{u}}^{\mathcal{D}}\f]
330  *
331  * @param key Specifes the linear system to solve.
332  * @param rhs Forcing term \f$\boldsymbol{f}\f$.
333  * @param inout Solution vector \f$\boldsymbol{\hat{u}}\f$.
334  * @param dirForcing .
335  */
337  const Array<OneD, const NekDouble>& rhs,
338  Array<OneD, NekDouble>& inout,
339  const Array<OneD, const NekDouble>& dirForcing)
340  {
341  int NumDirBcs = m_locToGloMap->GetNumGlobalDirBndCoeffs();
342  int contNcoeffs = m_locToGloMap->GetNumGlobalCoeffs();
343 
344  // STEP 1: SET THE DIRICHLET DOFS TO THE RIGHT VALUE
345  // IN THE SOLUTION ARRAY
347 
348  // STEP 2: CALCULATE THE HOMOGENEOUS COEFFICIENTS
349  if(contNcoeffs - NumDirBcs > 0)
350  {
352  LinSys->Solve(rhs,inout,m_locToGloMap,dirForcing);
353  }
354  }
355 
356 
357  /**
358  * The function searches the map #m_globalLinSys to see if the global
359  * matrix has been created before. If not, it calls the function
360  * #GenglobalLinSys to generate the requested global system.
361  *
362  * @param mkey Key specifying the linear system.
363  * @returns Pointer to the required linear system.
364  */
366  const GlobalLinSysKey &mkey)
367  {
368  return m_globalLinSysManager[mkey];
369  }
370 
372  const GlobalLinSysKey &mkey)
373  {
375  "To use method must have a AssemblyMap "
376  "attached to key");
378  }
379 
380 
381  /**
382  * The operation is evaluated locally (i.e. with respect to all local
383  * expansion modes) by the function ExpList#IProductWRTBase. The inner
384  * product with respect to the global expansion modes is than obtained
385  * by a global assembly operation.
386  *
387  * The values of the function \f$f(x)\f$ evaluated at the quadrature
388  * points \f$x_i\f$ should be contained in the variable #m_phys of the
389  * ExpList object \a in. The result is stored in the array
390  * #m_coeffs.
391  *
392  * @param In An ExpList, containing the discrete evaluation
393  * of \f$f(x)\f$ at the quadrature points in its
394  * array #m_phys.
395  */
397  const Array<OneD, const NekDouble> &inarray,
398  Array<OneD, NekDouble> &outarray,
399  CoeffState coeffstate)
400  {
401  if(coeffstate == eGlobal)
402  {
404  IProductWRTBase_IterPerExp(inarray,wsp);
405  Assemble(wsp,outarray);
406  }
407  else
408  {
409  IProductWRTBase_IterPerExp(inarray,outarray);
410  }
411  }
412 
413 
415  const Array<OneD, const NekDouble> &inarray,
416  Array<OneD, NekDouble> &outarray,
417  CoeffState coeffstate)
418  {
419  FwdTrans(inarray,outarray,coeffstate);
420  }
421 
423  const Array<OneD, const NekDouble> &inarray,
424  Array<OneD, NekDouble> &outarray,
425  CoeffState coeffstate)
426  {
427  MultiplyByInvMassMatrix(inarray,outarray,coeffstate);
428  }
429 
431  {
432  for(int i = 0; i < m_bndCondExpansions.num_elements(); ++i)
433  {
434  if(m_bndConditions[i]->GetBoundaryConditionType() == SpatialDomains::eDirichlet)
435  {
436  outarray[m_locToGloMap->GetBndCondCoeffsToGlobalCoeffsMap(i)]
437  = m_bndCondExpansions[i]->GetCoeff(0);
438  }
439  }
440  }
441 
442  /**
443  * This operation is evaluated as:
444  * \f{tabbing}
445  * \hspace{1cm} \= Do \= $e=$ $1, N_{\mathrm{el}}$ \\
446  * > > Do \= $i=$ $0,N_m^e-1$ \\
447  * > > > $\boldsymbol{\hat{u}}_g[\mbox{map}[e][i]] =
448  * \mbox{sign}[e][i] \cdot \boldsymbol{\hat{u}}^{e}[i]$\\
449  * > > continue\\
450  * > continue
451  * \f}
452  * where \a map\f$[e][i]\f$ is the mapping array and \a
453  * sign\f$[e][i]\f$ is an array of similar dimensions ensuring the
454  * correct modal connectivity between the different elements (both
455  * these arrays are contained in the data member #m_locToGloMap). This
456  * operation is equivalent to the gather operation
457  * \f$\boldsymbol{\hat{u}}_g=\mathcal{A}^{-1}\boldsymbol{\hat{u}}_l\f$,
458  * where \f$\mathcal{A}\f$ is the
459  * \f$N_{\mathrm{eof}}\times N_{\mathrm{dof}}\f$ permutation matrix.
460  *
461  */
463  {
464  m_locToGloMap->LocalToGlobal(m_coeffs,m_coeffs);
465  }
466 
467  /**
468  * This operation is evaluated as:
469  * \f{tabbing}
470  * \hspace{1cm} \= Do \= $e=$ $1, N_{\mathrm{el}}$ \\
471  * > > Do \= $i=$ $0,N_m^e-1$ \\
472  * > > > $\boldsymbol{\hat{u}}^{e}[i] = \mbox{sign}[e][i] \cdot
473  * \boldsymbol{\hat{u}}_g[\mbox{map}[e][i]]$ \\
474  * > > continue \\
475  * > continue
476  * \f}
477  * where \a map\f$[e][i]\f$ is the mapping array and
478  * \a sign\f$[e][i]\f$ is an array of similar dimensions ensuring the
479  * correct modal connectivity between the different elements (both
480  * these arrays are contained in the data member #m_locToGloMap). This
481  * operation is equivalent to the scatter operation
482  * \f$\boldsymbol{\hat{u}}_l=\mathcal{A}\boldsymbol{\hat{u}}_g\f$, where
483  * \f$\mathcal{A}\f$ is the
484  * \f$N_{\mathrm{eof}}\times N_{\mathrm{dof}}\f$ permutation matrix.
485  *
486  */
488  {
489  m_locToGloMap->GlobalToLocal(m_coeffs,m_coeffs);
490  }
491 
492  /**
493  * Consider the one dimensional Helmholtz equation,
494  * \f[\frac{d^2u}{dx^2}-\lambda u(x) = f(x),\f]
495  * supplemented with appropriate boundary conditions (which are
496  * contained in the data member #m_bndCondExpansions). Applying a
497  * \f$C^0\f$ continuous Galerkin discretisation, this equation leads to
498  * the following linear system:
499  * \f[\left( \boldsymbol{M}+\lambda\boldsymbol{L}\right)
500  * \boldsymbol{\hat{u}}_g=\boldsymbol{\hat{f}}\f]
501  * where \f$\boldsymbol{M}\f$ and \f$\boldsymbol{L}\f$ are the mass and
502  * Laplacian matrix respectively. This function solves the system above
503  * for the global coefficients \f$\boldsymbol{\hat{u}}\f$ by a call to
504  * the function #GlobalSolve.
505  *
506  * The values of the function \f$f(x)\f$ evaluated at the
507  * quadrature points \f$\boldsymbol{x}_i\f$ should be contained in the
508  * variable #m_phys of the ExpList object \a inarray. The resulting
509  * global coefficients \f$\boldsymbol{\hat{u}}_g\f$ are stored in the
510  * array #m_coeffs.
511  *
512  * @param inarray Input containing forcing function
513  * \f$\boldsymbol{f}\f$ at the quadrature points.
514  * @param outarray Output containing the coefficients
515  * \f$\boldsymbol{u}_g\f$
516  * @param lambda Parameter value.
517  * @param Sigma Coefficients of lambda.
518  * @param varcoeff Variable diffusivity coefficients.
519  * @param coeffstate
520  * @param dirForcing Dirichlet Forcing.
521  */
523  const Array<OneD, const NekDouble> &inarray,
524  Array<OneD, NekDouble> &outarray,
525  const FlagList &flags,
526  const StdRegions::ConstFactorMap &factors,
527  const StdRegions::VarCoeffMap &varcoeff,
528  const Array<OneD, const NekDouble> &dirForcing)
529  {
530  // Inner product of forcing
531  int contNcoeffs = m_locToGloMap->GetNumGlobalCoeffs();
532  Array<OneD,NekDouble> wsp(contNcoeffs);
533  IProductWRTBase(inarray,wsp,eGlobal);
534  // Note -1.0 term necessary to invert forcing function to
535  // be consistent with matrix definition
536  Vmath::Neg(contNcoeffs, wsp, 1);
537 
538  // Forcing function with weak boundary conditions
539  int i;
540  for(i = 0; i < m_bndCondExpansions.num_elements(); ++i)
541  {
542  if(m_bndConditions[i]->GetBoundaryConditionType() != SpatialDomains::eDirichlet)
543  {
544  wsp[m_locToGloMap->GetBndCondCoeffsToGlobalCoeffsMap(i)]
545  += m_bndCondExpansions[i]->GetCoeff(0);
546  }
547  }
548 
549  // Solve the system
551  m_locToGloMap,factors,varcoeff);
552 
553  if(flags.isSet(eUseGlobal))
554  {
555  GlobalSolve(key,wsp,outarray,dirForcing);
556  }
557  else
558  {
559  Array<OneD,NekDouble> tmp(contNcoeffs,0.0);
560  GlobalSolve(key,wsp,tmp,dirForcing);
561  GlobalToLocal(tmp,outarray);
562  }
563  }
564 
567  {
568  return GetBndConditions();
569  }
570 
572  const Array<OneD, const NekDouble> &inarray,
573  Array<OneD, NekDouble> &outarray,
574  CoeffState coeffstate)
575  {
576  BwdTrans(inarray,outarray,coeffstate);
577  }
578 
580  const Array<OneD, const NekDouble> &inarray,
581  Array<OneD, NekDouble> &outarray,
582  CoeffState coeffstate)
583  {
584  IProductWRTBase(inarray,outarray,coeffstate);
585  }
586 
587  /**
588  * This is equivalent to the operation:
589  * \f[\boldsymbol{M\hat{u}}_g\f]
590  * where \f$\boldsymbol{M}\f$ is the global matrix of type specified by
591  * \a mkey. After scattering the global array \a inarray to local
592  * level, this operation is evaluated locally by the function
593  * ExpList#GeneralMatrixOp. The global result is then obtained by a
594  * global assembly procedure.
595  *
596  * @param mkey This key uniquely defines the type matrix
597  * required for the operation.
598  * @param inarray The vector \f$\boldsymbol{\hat{u}}_g\f$ of size
599  * \f$N_{\mathrm{dof}}\f$.
600  * @param outarray The resulting vector of size
601  * \f$N_{\mathrm{dof}}\f$.
602  */
604  const GlobalMatrixKey &gkey,
605  const Array<OneD,const NekDouble> &inarray,
606  Array<OneD, NekDouble> &outarray,
607  CoeffState coeffstate)
608  {
609  if(coeffstate == eGlobal)
610  {
612  Array<OneD,NekDouble> tmp2(tmp1+m_ncoeffs);
613  GlobalToLocal(inarray,tmp1);
614  GeneralMatrixOp_IterPerExp(gkey,tmp1,tmp2);
615  Assemble(tmp2,outarray);
616  }
617  else
618  {
619  GeneralMatrixOp_IterPerExp(gkey,inarray,outarray);
620  }
621  }
622 
623 
624 
625  /**
626  * Reset the GlobalLinSys Manager
627  */
629  {
630  m_globalLinSysManager.ClearManager("GlobalLinSys");
631  }
632 
633  } // end of namespace
634 } //end of namespace
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)
void GeneralMatrixOp_IterPerExp(const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:829
Local coefficients.
virtual void v_MultiplyByInvMassMatrix(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
virtual void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
static boost::shared_ptr< DataType > AllocateSharedPtr()
Allocate a shared pointer from the memory pool.
Array< OneD, SpatialDomains::BoundaryConditionShPtr > m_bndConditions
An array which contains the information about the boundary condition on the different boundary region...
STL namespace.
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:251
virtual void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
Perform a forward transform.
Array< OneD, NekDouble > m_coeffs
Concatenation of all local expansion coefficients.
Definition: ExpList.h:939
void MultiplyByInvMassMatrix(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
Definition: MeshPartition.h:51
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate)
GlobalLinSysSharedPtr GetGlobalLinSys(const GlobalLinSysKey &mkey)
Returns the linear system specified by mkey.
Global coefficients.
bool isSet(const FlagType &key) const
Abstraction of a global continuous one-dimensional spectral/hp element expansion which approximates t...
Definition: ContField1D.h:56
void BwdTrans_IterPerExp(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:1623
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
void IProductWRTBase_IterPerExp(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:1573
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Definition: ExpList.h:917
LibUtilities::NekManager< GlobalLinSysKey, GlobalLinSys > m_globalLinSysManager
A manager which collects all the global linear systems being assembled, such that they should be cons...
Definition: ContField1D.h:164
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 ...
This class is the abstraction of a global discontinuous two- dimensional spectral/hp element expansio...
LibUtilities::SessionReaderSharedPtr m_session
Session.
Definition: ExpList.h:910
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:382
Defines a list of flags.
Describe a linear system.
virtual void v_ImposeDirichletConditions(Array< OneD, NekDouble > &outarray)
Impose the Dirichlet Boundary Conditions on outarray.
Describes a matrix with ordering defined by a local to global map.
const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > & GetBndConditions()
Definition: ContField1D.h:248
void GlobalSolve(const GlobalLinSysKey &key, const Array< OneD, const NekDouble > &rhs, Array< OneD, NekDouble > &inout, const Array< OneD, const NekDouble > &dirForcing=NullNekDouble1DArray)
Solve the linear system specified by the key key.
virtual void v_ClearGlobalLinSysManager(void)
virtual void v_GlobalToLocal(void)
Scatters from the global coefficients to the local coefficients .
ContField1D()
Default constructor.
Definition: ContField1D.cpp:86
Array< OneD, MultiRegions::ExpListSharedPtr > m_bndCondExpansions
Discretised boundary conditions.
bool LocToGloMapIsDefined() const
Returns true if a local to global map is defined.
void Assemble()
Assembles the global coefficients from the local coefficients .
Definition: ContField1D.h:309
void BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
This function performs the backward transformation of the spectral/hp element expansion.
This class is the abstraction of a one-dimensional multi-elemental expansions which is merely a colle...
Definition: ExpList1D.h:61
boost::shared_ptr< GlobalLinSys > GlobalLinSysSharedPtr
Pointer to a GlobalLinSys object.
Definition: GlobalLinSys.h:52
virtual ~ContField1D()
Destructor.
GlobalLinSysSharedPtr GenGlobalLinSys(const GlobalLinSysKey &mkey)
void FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
Perform global forward transformation of a function ,.
PeriodicMap m_periodicVerts
A map which identifies groups of periodic vertices.
virtual void v_LocalToGlobal(void)
Gathers the global coefficients from the local coefficients .
#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
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...
AssemblyMapCGSharedPtr m_locToGloMap
(A shared pointer to) the object which contains all the required information for the transformation f...
Definition: ContField1D.h:149
virtual const Array< OneD, const SpatialDomains::BoundaryConditionShPtr > & v_GetBndConditions()