Nektar++
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 // Permission is hereby granted, free of charge, to any person obtaining a
14 // copy of this software and associated documentation files (the "Software"),
15 // to deal in the Software without restriction, including without limitation
16 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
17 // and/or sell copies of the Software, and to permit persons to whom the
18 // Software is furnished to do so, subject to the following conditions:
19 //
20 // The above copyright notice and this permission notice shall be included
21 // in all copies or substantial portions of the Software.
22 //
23 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29 // DEALINGS IN THE SOFTWARE.
30 //
31 // Description: Field definition for 1D domain with boundary conditions
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
37 
38 namespace Nektar
39 {
40  namespace MultiRegions
41  {
42  /**
43  * @class ContField1D
44  * As opposed to the class #ContExpList1D, the class #ContField1D is
45  * able to incorporate the boundary conditions imposed to the problem
46  * to be solved. Therefore, the class is equipped with three additional
47  * data members:
48  * - #m_bndCondExpansions
49  * - #m_bndTypes
50  * - #m_bndCondEquations
51  *
52  * The first data structure, #m_bndCondExpansions,
53  * contains the point Expansion on the boundary, #m_bndTypes
54  * stores information about the type of boundary condition on the
55  * different parts of the boundary while #m_bndCondEquations holds the
56  * equation of the imposed boundary conditions.
57  *
58  * Furthermore, in case of Dirichlet boundary conditions,
59  * this class is capable of lifting a known solution satisfying these
60  * boundary conditions. If we denote the unknown solution by
61  * \f$u^{\mathcal{H}}(\boldsymbol{x})\f$ and the known Dirichlet
62  * boundary conditions by \f$u^{\mathcal{D}}(\boldsymbol{x})\f$, the
63  * expansion then can be decomposed as
64  * \f[ u^{\delta}(\boldsymbol{x}_i)=u^{\mathcal{D}}(\boldsymbol{x}_i)+
65  * u^{\mathcal{H}}(\boldsymbol{x}_i)=\sum_{n=0}^{N^{\mathcal{D}}-1}
66  * \hat{u}_n^{\mathcal{D}}\Phi_n(\boldsymbol{x}_i)+
67  * \sum_{n={N^{\mathcal{D}}}}^{N_{\mathrm{dof}}
68  * -1}\hat{u}_n^{\mathcal{H}}
69  * \Phi_n(\boldsymbol{x}_i).\f]
70  * This lifting is accomplished by ordering the known global degrees of
71  * freedom, prescribed by the Dirichlet boundary conditions, first in
72  * the global array \f$\boldsymbol{\hat{u}}\f$, that is,
73  * \f[\boldsymbol{\hat{u}}=\left[ \begin{array}{c}
74  * \boldsymbol{\hat{u}}^{\mathcal{D}}\\
75  * \boldsymbol{\hat{u}}^{\mathcal{H}}
76  * \end{array} \right].\f]
77  * Such kind of expansions are also referred to as continuoous fields.
78  * This class should be used when solving 2D problems using a standard
79  * Galerkin approach.
80  */
81 
82  /**
83  * Constructs an empty 1D continuous field.
84  */
87  m_locToGloMap(),
88  m_globalLinSysManager(
89  std::bind(&ContField1D::GenGlobalLinSys, this, std::placeholders::_1),
90  std::string("GlobalLinSys"))
91  {
92  }
93 
94 
95  /**
96  * Given a mesh \a graph1D, containing information about the domain and
97  * the spectral/hp element expansion, this constructor fills the list
98  * of local expansions #m_exp with the proper expansions, calculates
99  * the total number of quadrature points \f$\boldsymbol{x}_i\f$ and
100  * local expansion coefficients \f$\hat{u}^e_n\f$ and allocates
101  * memory for the arrays #m_coeffs and #m_phys. Furthermore, it
102  * constructs the mapping array (contained in #m_locToGloMap) for the
103  * transformation between local elemental level and global level, it
104  * calculates the total number global expansion coefficients
105  * \f$\hat{u}_n\f$.
106  * The constructor also discretises the boundary conditions, specified
107  * by the argument \a bcs, by expressing them in terms of the
108  * coefficient of the expansion on the boundary.
109  *
110  * @param graph1D A 1D mesh, containing information about the
111  * domain and the spectral/hp element expansion.
112  * @param bcs The boundary conditions.
113  * @param variable An optional parameter to indicate for which
114  * variable the field should be constructed.
115  */
117  const LibUtilities::SessionReaderSharedPtr &pSession,
118  const SpatialDomains::MeshGraphSharedPtr &graph1D,
119  const std::string &variable,
120  const Collections::ImplementationType ImpType):
121  DisContField1D(pSession,graph1D,variable,false,ImpType),
122  m_locToGloMap(),
124  std::bind(&ContField1D::GenGlobalLinSys, this, std::placeholders::_1),
125  std::string("GlobalLinSys"))
126  {
127  SpatialDomains::BoundaryConditions bcs(pSession, graph1D);
128 
133  false,
134  variable,
136  }
137 
138 
139  /**
140  * Constructs a 1D continuous field as a copy of an existing field
141  * including all boundary conditions.
142  * @param In Existing continuous field to duplicate.
143  */
145  DisContField1D(In),
148  std::bind(&ContField1D::GenGlobalLinSys, this, std::placeholders::_1),
149  std::string("GlobalLinSys"))
150  {
151  }
152 
153  /**
154  * Constructs a 1D continuous field as a copy of an existing explist 1D field
155  * and adding all the boundary conditions.
156  * @param In Existing explist1D field .
157  */
159  DisContField1D(In),
160  m_locToGloMap(),
162  std::bind(&ContField1D::GenGlobalLinSys, this, std::placeholders::_1),
163  std::string("GlobalLinSys"))
164  {
166  ::AllocateSharedPtr(pSession, m_ncoeffs, In);
167 
168  }
169 
170  /**
171  *
172  */
174  {
175  }
176 
177 
178  /**
179  * Given a function \f$f(x)\f$ defined at the quadrature
180  * points, this function determines the unknown global coefficients
181  * \f$\boldsymbol{\hat{u}}^{\mathcal{H}}\f$ employing a discrete
182  * Galerkin projection from physical space to coefficient
183  * space. The operation is evaluated by the function #GlobalSolve using
184  * the global mass matrix.
185  *
186  * The values of the function \f$f(x)\f$ evaluated at the
187  * quadrature points \f$x_i\f$ should be contained in the
188  * variable #m_phys of the ExpList object \a Sin. The resulting global
189  * coefficients \f$\hat{u}_g\f$ are stored in the array #m_coeffs.
190  *
191  * @param inarray Discrete \f$f(x)\f$.
192  * @param outarray Storage for result.
193  * @param coeffstate
194  */
196  Array<OneD, NekDouble> &outarray,
197  CoeffState coeffstate)
198  {
199  // Inner product of forcing
201  IProductWRTBase(inarray,wsp,eGlobal);
202 
203  // Solve the system
205 
206  GlobalSolve(key,wsp,outarray);
207  if(coeffstate == eLocal)
208  {
209  GlobalToLocal(outarray,outarray);
210  }
211  }
212 
213 
214  /**
215  * Given the coefficients of an expansion, this function evaluates the
216  * spectral/hp expansion \f$u^{\delta}(x)\f$ at the quadrature
217  * points \f$x_i\f$. This operation is evaluated locally by the
218  * function ExpList#BwdTrans.
219  *
220  * The coefficients of the expansion should be contained in the
221  * variable #m_coeffs of the ExpList object \a In. The resulting
222  * physical values at the quadrature points \f$u^{\delta}(x_i)\f$ are
223  * stored in the array #m_phys.
224  *
225  * @param In An ExpList, containing the local
226  * coefficients \f$\hat{u}_n^e\f$ in its array
227  * #m_coeffs.
228  */
230  const Array<OneD, const NekDouble> &inarray,
231  Array<OneD, NekDouble> &outarray,
232  CoeffState coeffstate)
233  {
234  Array<OneD, NekDouble> tmpinarray;
235  if(coeffstate != eLocal)
236  {
237  tmpinarray = Array<OneD, NekDouble>(inarray);
238  GlobalToLocal(inarray,tmpinarray);
239  }
240  else
241  {
242  tmpinarray = inarray;
243  }
244 
245  BwdTrans_IterPerExp(tmpinarray,outarray);
246  }
247 
248 
249  /**
250  *
251  */
253  const Array<OneD, const NekDouble> &inarray,
254  Array<OneD, NekDouble> &outarray,
255  CoeffState coeffstate)
256  {
258  if(coeffstate == eGlobal)
259  {
260  if(inarray.data() == outarray.data())
261  {
262  Array<OneD, NekDouble> tmp(inarray);
263  GlobalSolve(key,tmp,outarray);
264  }
265  else
266  {
267  GlobalSolve(key,inarray,outarray);
268  }
269  }
270  else
271  {
272  Array<OneD, NekDouble> globaltmp(m_ncoeffs,0.0);
273 
274  if(inarray.data() == outarray.data())
275  {
276  Array<OneD,NekDouble> tmp(inarray);
277  Assemble(tmp,outarray);
278  }
279  else
280  {
281  Assemble(inarray,outarray);
282  }
283 
284  GlobalSolve(key,outarray,globaltmp);
285  GlobalToLocal(globaltmp,outarray);
286  }
287  }
288 
289 
290  /**
291  * Given a linear system specified by the key \a key,
292  * \f[\boldsymbol{M}\boldsymbol{\hat{u}}_g=\boldsymbol{\hat{f}},\f]
293  * this function solves this linear system taking into account the
294  * boundary conditions specified in the data member
295  * #m_bndCondExpansions.
296  * Therefore, it adds an array \f$\boldsymbol{\hat{g}}\f$ which
297  * represents the non-zero surface integral resulting from the weak
298  * boundary conditions (e.g. Neumann boundary conditions) to the right
299  * hand side, that is,
300  * \f[\boldsymbol{M}\boldsymbol{\hat{u}}_g=\boldsymbol{\hat{f}}+
301  * \boldsymbol{\hat{g}}.\f]
302  * Furthermore, it lifts the known degrees of freedom which are
303  * prescribed by the Dirichlet boundary conditions. As these known
304  * coefficients \f$\boldsymbol{\hat{u}}^{\mathcal{D}}\f$ are numbered
305  * first in the global coefficient array \f$\boldsymbol{\hat{u}}_g\f$,
306  * the linear system can be decomposed as,
307  * \f[\left[\begin{array}{cc}
308  * \boldsymbol{M}^{\mathcal{DD}}&\boldsymbol{M}^{\mathcal{DH}}\\
309  * \boldsymbol{M}^{\mathcal{HD}}&\boldsymbol{M}^{\mathcal{HH}}
310  * \end{array}\right]
311  * \left[\begin{array}{c}
312  * \boldsymbol{\hat{u}}^{\mathcal{D}}\\
313  * \boldsymbol{\hat{u}}^{\mathcal{H}}
314  * \end{array}\right]=
315  * \left[\begin{array}{c}
316  * \boldsymbol{\hat{f}}^{\mathcal{D}}\\
317  * \boldsymbol{\hat{f}}^{\mathcal{H}}
318  * \end{array}\right]+
319  * \left[\begin{array}{c}
320  * \boldsymbol{\hat{g}}^{\mathcal{D}}\\
321  * \boldsymbol{\hat{g}}^{\mathcal{H}}
322  * \end{array}\right]
323  * \f]
324  * which will then be solved for the unknown coefficients
325  * \f$\boldsymbol{\hat{u}}^{\mathcal{H}}\f$ as,
326  * \f[
327  * \boldsymbol{M}^{\mathcal{HH}}\boldsymbol{\hat{u}}^{\mathcal{H}}
328  * = \boldsymbol{\hat{f}}^{\mathcal{H}}
329  * +\boldsymbol{\hat{g}}^{\mathcal{H}}
330  * -\boldsymbol{M}^{\mathcal{HD}}\boldsymbol{\hat{u}}^{\mathcal{D}}\f]
331  *
332  * @param key Specifes the linear system to solve.
333  * @param rhs Forcing term \f$\boldsymbol{f}\f$.
334  * @param inout Solution vector \f$\boldsymbol{\hat{u}}\f$.
335  * @param dirForcing .
336  */
339  Array<OneD, NekDouble>& inout,
340  const Array<OneD, const NekDouble>& dirForcing)
341  {
342  int NumDirBcs = m_locToGloMap->GetNumGlobalDirBndCoeffs();
343  int contNcoeffs = m_locToGloMap->GetNumGlobalCoeffs();
344 
345  // STEP 1: SET THE DIRICHLET DOFS TO THE RIGHT VALUE
346  // IN THE SOLUTION ARRAY
348 
349  // STEP 2: CALCULATE THE HOMOGENEOUS COEFFICIENTS
350  if(contNcoeffs - NumDirBcs > 0)
351  {
353  LinSys->Solve(rhs,inout,m_locToGloMap,dirForcing);
354  }
355  }
356 
357 
358  /**
359  * The function searches the map #m_globalLinSys to see if the global
360  * matrix has been created before. If not, it calls the function
361  * #GenglobalLinSys to generate the requested global system.
362  *
363  * @param mkey Key specifying the linear system.
364  * @returns Pointer to the required linear system.
365  */
367  const GlobalLinSysKey &mkey)
368  {
369  return m_globalLinSysManager[mkey];
370  }
371 
373  const GlobalLinSysKey &mkey)
374  {
376  "To use method must have a AssemblyMap "
377  "attached to key");
379  }
380 
381 
382  /**
383  * The operation is evaluated locally (i.e. with respect to all local
384  * expansion modes) by the function ExpList#IProductWRTBase. The inner
385  * product with respect to the global expansion modes is than obtained
386  * by a global assembly operation.
387  *
388  * The values of the function \f$f(x)\f$ evaluated at the quadrature
389  * points \f$x_i\f$ should be contained in the variable #m_phys of the
390  * ExpList object \a in. The result is stored in the array
391  * #m_coeffs.
392  *
393  * @param In An ExpList, containing the discrete evaluation
394  * of \f$f(x)\f$ at the quadrature points in its
395  * array #m_phys.
396  */
398  const Array<OneD, const NekDouble> &inarray,
399  Array<OneD, NekDouble> &outarray,
400  CoeffState coeffstate)
401  {
402  if(coeffstate == eGlobal)
403  {
405  IProductWRTBase_IterPerExp(inarray,wsp);
406  Assemble(wsp,outarray);
407  }
408  else
409  {
410  IProductWRTBase_IterPerExp(inarray,outarray);
411  }
412  }
413 
414 
416  const Array<OneD, const NekDouble> &inarray,
417  Array<OneD, NekDouble> &outarray,
418  CoeffState coeffstate)
419  {
420  FwdTrans(inarray,outarray,coeffstate);
421  }
422 
424  const Array<OneD, const NekDouble> &inarray,
425  Array<OneD, NekDouble> &outarray,
426  CoeffState coeffstate)
427  {
428  MultiplyByInvMassMatrix(inarray,outarray,coeffstate);
429  }
430 
432  {
433  for(int i = 0; i < m_bndCondExpansions.num_elements(); ++i)
434  {
435  if(m_bndConditions[i]->GetBoundaryConditionType() == SpatialDomains::eDirichlet)
436  {
437  outarray[m_locToGloMap->GetBndCondCoeffsToGlobalCoeffsMap(i)]
438  = m_bndCondExpansions[i]->GetCoeff(0);
439  }
440  }
441  }
442 
443  /**
444  * This operation is evaluated as:
445  * \f{tabbing}
446  * \hspace{1cm} \= Do \= $e=$ $1, N_{\mathrm{el}}$ \\
447  * > > Do \= $i=$ $0,N_m^e-1$ \\
448  * > > > $\boldsymbol{\hat{u}}_g[\mbox{map}[e][i]] =
449  * \mbox{sign}[e][i] \cdot \boldsymbol{\hat{u}}^{e}[i]$\\
450  * > > continue\\
451  * > continue
452  * \f}
453  * where \a map\f$[e][i]\f$ is the mapping array and \a
454  * sign\f$[e][i]\f$ is an array of similar dimensions ensuring the
455  * correct modal connectivity between the different elements (both
456  * these arrays are contained in the data member #m_locToGloMap). This
457  * operation is equivalent to the gather operation
458  * \f$\boldsymbol{\hat{u}}_g=\mathcal{A}^{-1}\boldsymbol{\hat{u}}_l\f$,
459  * where \f$\mathcal{A}\f$ is the
460  * \f$N_{\mathrm{eof}}\times N_{\mathrm{dof}}\f$ permutation matrix.
461  *
462  */
464  const Array<OneD, const NekDouble> &inarray,
465  Array<OneD,NekDouble> &outarray, bool useComm)
466  {
467  m_locToGloMap->LocalToGlobal(inarray, outarray, useComm);
468  }
469 
470 
471  void ContField1D::v_LocalToGlobal(bool useComm)
472  {
473  m_locToGloMap->LocalToGlobal(m_coeffs,m_coeffs, useComm);
474  }
475 
476  /**
477  * This operation is evaluated as:
478  * \f{tabbing}
479  * \hspace{1cm} \= Do \= $e=$ $1, N_{\mathrm{el}}$ \\
480  * > > Do \= $i=$ $0,N_m^e-1$ \\
481  * > > > $\boldsymbol{\hat{u}}^{e}[i] = \mbox{sign}[e][i] \cdot
482  * \boldsymbol{\hat{u}}_g[\mbox{map}[e][i]]$ \\
483  * > > continue \\
484  * > continue
485  * \f}
486  * where \a map\f$[e][i]\f$ is the mapping array and
487  * \a sign\f$[e][i]\f$ is an array of similar dimensions ensuring the
488  * correct modal connectivity between the different elements (both
489  * these arrays are contained in the data member #m_locToGloMap). This
490  * operation is equivalent to the scatter operation
491  * \f$\boldsymbol{\hat{u}}_l=\mathcal{A}\boldsymbol{\hat{u}}_g\f$, where
492  * \f$\mathcal{A}\f$ is the
493  * \f$N_{\mathrm{eof}}\times N_{\mathrm{dof}}\f$ permutation matrix.
494  *
495  */
497  const Array<OneD, const NekDouble> &inarray,
498  Array<OneD,NekDouble> &outarray)
499  {
500  m_locToGloMap->GlobalToLocal(inarray, outarray);
501  }
502 
504  {
505  m_locToGloMap->GlobalToLocal(m_coeffs,m_coeffs);
506  }
507 
508  /**
509  * Consider the one dimensional Helmholtz equation,
510  * \f[\frac{d^2u}{dx^2}-\lambda u(x) = f(x),\f]
511  * supplemented with appropriate boundary conditions (which are
512  * contained in the data member #m_bndCondExpansions). Applying a
513  * \f$C^0\f$ continuous Galerkin discretisation, this equation leads to
514  * the following linear system:
515  * \f[\left( \boldsymbol{M}+\lambda\boldsymbol{L}\right)
516  * \boldsymbol{\hat{u}}_g=\boldsymbol{\hat{f}}\f]
517  * where \f$\boldsymbol{M}\f$ and \f$\boldsymbol{L}\f$ are the mass and
518  * Laplacian matrix respectively. This function solves the system above
519  * for the global coefficients \f$\boldsymbol{\hat{u}}\f$ by a call to
520  * the function #GlobalSolve.
521  *
522  * The values of the function \f$f(x)\f$ evaluated at the
523  * quadrature points \f$\boldsymbol{x}_i\f$ should be contained in the
524  * variable #m_phys of the ExpList object \a inarray. The resulting
525  * global coefficients \f$\boldsymbol{\hat{u}}_g\f$ are stored in the
526  * array #m_coeffs.
527  *
528  * @param inarray Input containing forcing function
529  * \f$\boldsymbol{f}\f$ at the quadrature points.
530  * @param outarray Output containing the coefficients
531  * \f$\boldsymbol{u}_g\f$
532  * @param lambda Parameter value.
533  * @param Sigma Coefficients of lambda.
534  * @param varcoeff Variable diffusivity coefficients.
535  * @param coeffstate
536  * @param dirForcing Dirichlet Forcing.
537  */
539  const Array<OneD, const NekDouble> &inarray,
540  Array<OneD, NekDouble> &outarray,
541  const FlagList &flags,
542  const StdRegions::ConstFactorMap &factors,
543  const StdRegions::VarCoeffMap &varcoeff,
544  const MultiRegions::VarFactorsMap &varfactors,
545  const Array<OneD, const NekDouble> &dirForcing,
546  const bool PhysSpaceForcing)
547  {
548  // Inner product of forcing
549  int contNcoeffs = m_locToGloMap->GetNumGlobalCoeffs();
550  Array<OneD,NekDouble> wsp(contNcoeffs);
551  if(PhysSpaceForcing)
552  {
553  IProductWRTBase(inarray,wsp,eGlobal);
554  }
555  else
556  {
557  Assemble(inarray,wsp);
558  }
559  // Note -1.0 term necessary to invert forcing function to
560  // be consistent with matrix definition
561  Vmath::Neg(contNcoeffs, wsp, 1);
562 
563  // Forcing function with weak boundary conditions
564  int i;
565  for(i = 0; i < m_bndCondExpansions.num_elements(); ++i)
566  {
567  if(m_bndConditions[i]->GetBoundaryConditionType() ==
569  m_bndConditions[i]->GetBoundaryConditionType() ==
571  {
572  wsp[m_locToGloMap->GetBndCondCoeffsToGlobalCoeffsMap(i)]
573  += m_bndCondExpansions[i]->GetCoeff(0);
574  }
575  }
576 
577  // Solve the system
579  m_locToGloMap,factors,varcoeff,varfactors);
580 
581  if(flags.isSet(eUseGlobal))
582  {
583  GlobalSolve(key,wsp,outarray,dirForcing);
584  }
585  else
586  {
587  Array<OneD,NekDouble> tmp(contNcoeffs,0.0);
588  GlobalSolve(key,wsp,tmp,dirForcing);
589  GlobalToLocal(tmp,outarray);
590  }
591  }
592 
595  {
596  return GetBndConditions();
597  }
598 
600  const Array<OneD, const NekDouble> &inarray,
601  Array<OneD, NekDouble> &outarray,
602  CoeffState coeffstate)
603  {
604  BwdTrans(inarray,outarray,coeffstate);
605  }
606 
608  const Array<OneD, const NekDouble> &inarray,
609  Array<OneD, NekDouble> &outarray,
610  CoeffState coeffstate)
611  {
612  IProductWRTBase(inarray,outarray,coeffstate);
613  }
614 
615  /**
616  * This is equivalent to the operation:
617  * \f[\boldsymbol{M\hat{u}}_g\f]
618  * where \f$\boldsymbol{M}\f$ is the global matrix of type specified by
619  * \a mkey. After scattering the global array \a inarray to local
620  * level, this operation is evaluated locally by the function
621  * ExpList#GeneralMatrixOp. The global result is then obtained by a
622  * global assembly procedure.
623  *
624  * @param mkey This key uniquely defines the type matrix
625  * required for the operation.
626  * @param inarray The vector \f$\boldsymbol{\hat{u}}_g\f$ of size
627  * \f$N_{\mathrm{dof}}\f$.
628  * @param outarray The resulting vector of size
629  * \f$N_{\mathrm{dof}}\f$.
630  */
632  const GlobalMatrixKey &gkey,
633  const Array<OneD,const NekDouble> &inarray,
634  Array<OneD, NekDouble> &outarray,
635  CoeffState coeffstate)
636  {
637  if(coeffstate == eGlobal)
638  {
640  Array<OneD,NekDouble> tmp2(tmp1+m_ncoeffs);
641  GlobalToLocal(inarray,tmp1);
642  GeneralMatrixOp_IterPerExp(gkey,tmp1,tmp2);
643  Assemble(tmp2,outarray);
644  }
645  else
646  {
647  GeneralMatrixOp_IterPerExp(gkey,inarray,outarray);
648  }
649  }
650 
651 
652 
653  /**
654  * Reset the GlobalLinSys Manager
655  */
657  {
658  m_globalLinSysManager.ClearManager("GlobalLinSys");
659  }
660 
661  } // end of namespace
662 } //end of namespace
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:163
void GeneralMatrixOp_IterPerExp(const GlobalMatrixKey &gkey, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: ExpList.cpp:1019
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)
Array< OneD, SpatialDomains::BoundaryConditionShPtr > m_bndConditions
An array which contains the information about the boundary condition on the different boundary region...
std::shared_ptr< GlobalLinSys > GlobalLinSysSharedPtr
Pointer to a GlobalLinSys object.
Definition: GlobalLinSys.h:50
STL namespace.
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:294
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 MultiRegions::VarFactorsMap &varfactors, const Array< OneD, const NekDouble > &dirForcing, const bool PhysSpaceForcing)
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:1052
void MultiplyByInvMassMatrix(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, CoeffState coeffstate=eLocal)
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.
Abstraction of a global continuous one-dimensional spectral/hp element expansion which approximates t...
Definition: ContField1D.h:55
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:1786
std::map< StdRegions::VarCoeffType, Array< OneD, NekDouble > > VarCoeffMap
Definition: StdRegions.hpp:264
bool isSet(const FlagType &key) const
bool LocToGloMapIsDefined() const
Returns true if a local to global map is defined.
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:1725
int m_ncoeffs
The total number of local degrees of freedom. m_ncoeffs .
Definition: ExpList.h:1030
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:155
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 ...
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
This class is the abstraction of a global discontinuous two- dimensional spectral/hp element expansio...
LibUtilities::SessionReaderSharedPtr m_session
Session.
Definition: ExpList.h:1023
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:399
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:249
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)
ContField1D()
Default constructor.
Definition: ContField1D.cpp:85
Array< OneD, MultiRegions::ExpListSharedPtr > m_bndCondExpansions
Discretised boundary conditions.
void Assemble()
Assembles the global coefficients from the local coefficients .
Definition: ContField1D.h:276
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:58
std::shared_ptr< GlobalLinSys > GenGlobalLinSys(const GlobalLinSysKey &mkey, const std::shared_ptr< AssemblyMapCG > &locToGloMap)
This operation constructs the global linear system of type mkey.
Definition: ExpList.cpp:1362
virtual ~ContField1D()
Destructor.
virtual void v_LocalToGlobal(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool useComm)
Gathers the global coefficients from the local coefficients .
GlobalLinSysSharedPtr GenGlobalLinSys(const GlobalLinSysKey &mkey)
StandardMatrixTag boost::call_traits< LhsDataType >::const_reference rhs
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.
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:250
void GlobalToLocal(void)
Scatters from the global coefficients to the local coefficients .
Definition: ExpList.h:2096
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...
std::shared_ptr< SessionReader > SessionReaderSharedPtr
AssemblyMapCGSharedPtr m_locToGloMap
(A shared pointer to) the object which contains all the required information for the transformation f...
Definition: ContField1D.h:140
virtual const Array< OneD, const SpatialDomains ::BoundaryConditionShPtr > & v_GetBndConditions()
std::map< StdRegions::ConstFactorType, Array< OneD, NekDouble > > VarFactorsMap