Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
ContField1D.h
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File ContField1D.h
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 in one-dimension
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
36 #ifndef NEKTAR_LIBS_MULTIREGIONS_CONTFIELD1D_H
37 #define NEKTAR_LIBS_MULTIREGIONS_CONTFIELD1D_H
38 
39 
48 
49 namespace Nektar
50 {
51  namespace MultiRegions
52  {
53  /// Abstraction of a global continuous one-dimensional spectral/hp
54  /// element expansion which approximates the solution of a set of
55  /// partial differential equations.
57  {
58  public:
59  /// Default constructor.
61 
62  /// Set up global continuous field based on an input mesh and
63  /// boundary conditions.
67  const std::string &variable);
68 
69  /// Copy constructor.
70  MULTI_REGIONS_EXPORT ContField1D(const ContField1D &In);
71 
72  /// Copy constructor.
74  const ExpList1D & In);
75  /// Destructor
77 
78  /// Perform global forward transformation of a function \f$f(x)\f$,
79  // subject to the boundary conditions specified.
80  MULTI_REGIONS_EXPORT void FwdTrans( const Array<OneD, const NekDouble> &inarray,
81  Array<OneD, NekDouble> &outarray,
82  CoeffState coeffstate = eLocal);
83 
84  /// This function performs the backward transformation of the
85  /// spectral/hp element expansion.
86  MULTI_REGIONS_EXPORT void BwdTrans( const Array<OneD, const NekDouble> &inarray,
87  Array<OneD, NekDouble> &outarray,
88  CoeffState coeffstate = eLocal);
89 
90  ///
92  const Array<OneD, const NekDouble> &inarray,
93  Array<OneD, NekDouble> &outarray,
94  CoeffState coeffstate = eLocal);
95 
96  /// Return the boundary conditions expansion.
97  // inline
98  MULTI_REGIONS_EXPORT const Array<OneD,const MultiRegions::ExpListSharedPtr>&
100 
101  // inline
102  MULTI_REGIONS_EXPORT const Array<OneD,const SpatialDomains
104  /// Scatters from the global coefficients
105  /// \f$\boldsymbol{\hat{u}}_g\f$ to the local coefficients
106  /// \f$\boldsymbol{\hat{u}}_l\f$.
107  // inline
108  MULTI_REGIONS_EXPORT void GlobalToLocal( const Array<OneD, const NekDouble> &inarray,
109  Array<OneD,NekDouble> &outarray);
110 
111  /// Gathers the global coefficients \f$\boldsymbol{\hat{u}}_g\f$
112  /// from the local coefficients \f$\boldsymbol{\hat{u}}_l\f$.
113  // inline
115 
116  /// Assembles the global coefficients \f$\boldsymbol{\hat{u}}_g\f$
117  /// from the local coefficients \f$\boldsymbol{\hat{u}}_l\f$.
118  // inline
120 
121  /// Assembles the global coefficients \f$\boldsymbol{\hat{u}}_g\f$
122  /// from the local coefficients \f$\boldsymbol{\hat{u}}_l\f$.
123  // inline
124  MULTI_REGIONS_EXPORT void Assemble(const Array<OneD, const NekDouble> &inarray,
125  Array<OneD,NekDouble> &outarray);
126 
127  /// Returns the map from local to global level.
128  // inline
130 
131  /// Calculates the inner product of a function \f$f(x)\f$ with
132  /// respect to all <em>global</em> expansion modes
133  /// \f$\phi_n^e(x)\f$.
134  MULTI_REGIONS_EXPORT void IProductWRTBase(const Array<OneD, const NekDouble> &inarray,
135  Array<OneD, NekDouble> &outarray,
136  CoeffState coeffstate = eLocal);
137 
138  /// Calculates the result of the multiplication of a global matrix
139  /// of type specified by \a mkey with a vector given by \a inarray.
141  const Array<OneD,const NekDouble> &inarray,
142  Array<OneD, NekDouble> &outarray,
143  CoeffState coeffstate = eLocal);
144 
145  protected:
146  /// (A shared pointer to) the object which contains all the required
147  /// information for the transformation from local to global degrees
148  /// of freedom.
150 
151 
152  /// A enum list declaring how to interpret coeffs,
153  /// i.e. eLocal, eHybrid or eGlobal
155 
156  /// (A shared pointer to) a list which collects all the global
157  /// matrices being assembled, such that they should be constructed
158  /// only once.
160 
161  /// A manager which collects all the global
162  /// linear systems being assembled, such that they should be
163  /// constructed only once.
165 
166  private:
167  /// Returns the linear system specified by \a mkey.
169 
171 
172  /// Solve the linear system specified by the key \a key.
173  void GlobalSolve( const GlobalLinSysKey &key,
174  const Array<OneD, const NekDouble> &rhs,
175  Array<OneD, NekDouble> &inout,
176  const Array<OneD, const NekDouble> &dirForcing
178 
179  /// Perform a forward transform
180  virtual void v_FwdTrans(
181  const Array<OneD, const NekDouble> &inarray,
182  Array<OneD, NekDouble> &outarray,
183  CoeffState coeffstate);
184 
185  virtual void v_MultiplyByInvMassMatrix(
186  const Array<OneD, const NekDouble> &inarray,
187  Array<OneD, NekDouble> &outarray,
188  CoeffState coeffstate);
189 
190  /// Impose the Dirichlet Boundary Conditions on outarray
191  MULTI_REGIONS_EXPORT virtual void v_ImposeDirichletConditions(Array<OneD,NekDouble>& outarray);
192 
193  /// Scatters from the global coefficients
194  /// \f$\boldsymbol{\hat{u}}_g\f$ to the local coefficients
195  /// \f$\boldsymbol{\hat{u}}_l\f$.
196  // inline
197  MULTI_REGIONS_EXPORT virtual void v_GlobalToLocal(void);
198 
199 
200  /// Gathers the global coefficients \f$\boldsymbol{\hat{u}}_g\f$
201  /// from the local coefficients \f$\boldsymbol{\hat{u}}_l\f$.
202  // inline
203  MULTI_REGIONS_EXPORT virtual void v_LocalToGlobal(void);
204 
205  virtual void v_HelmSolve(
206  const Array<OneD, const NekDouble> &inarray,
207  Array<OneD, NekDouble> &outarray,
208  const FlagList &flags,
209  const StdRegions::ConstFactorMap &factors,
210  const StdRegions::VarCoeffMap &varcoeff,
211  const Array<OneD, const NekDouble> &dirForcing);
212 
213  virtual const Array<OneD,const SpatialDomains
215 
216  virtual void v_BwdTrans(
217  const Array<OneD, const NekDouble> &inarray,
218  Array<OneD, NekDouble> &outarray,
219  CoeffState coeffstate);
220 
221  virtual void v_IProductWRTBase(
222  const Array<OneD, const NekDouble> &inarray,
223  Array<OneD, NekDouble> &outarray,
224  CoeffState coeffstate);
225 
226  /// Calculates the result of the multiplication of a global matrix
227  /// of type specified by \a mkey with a vector given by \a inarray.
228  virtual void v_GeneralMatrixOp(
229  const GlobalMatrixKey &gkey,
230  const Array<OneD,const NekDouble> &inarray,
231  Array<OneD, NekDouble> &outarray,
232  CoeffState coeffstate);
233 
234  };
235  typedef boost::shared_ptr<ContField1D> ContField1DSharedPtr;
236 
237  // Inline implementations follow
238 
239  inline const Array<OneD,const MultiRegions::ExpListSharedPtr>&
241  {
242  return m_bndCondExpansions;
243  }
244 
245  inline const Array<OneD,const SpatialDomains::BoundaryConditionShPtr>&
247  {
248  return m_bndConditions;
249  }
250 
251 
252  /**
253  * This operation is evaluated as:
254  * \f{tabbing}
255  * \hspace{1cm} \= Do \= $e=$ $1, N_{\mathrm{el}}$ \\
256  * > > Do \= $i=$ $0,N_m^e-1$ \\
257  * > > > $\boldsymbol{\hat{u}}^{e}[i] = \mbox{sign}[e][i] \cdot
258  * \boldsymbol{\hat{u}}_g[\mbox{map}[e][i]]$ \\
259  * > > continue \\
260  * > continue
261  * \f}
262  * where \a map\f$[e][i]\f$ is the mapping array and \a
263  * sign\f$[e][i]\f$ is an array of similar dimensions ensuring the
264  * correct modal connectivity between the different elements (both
265  * these arrays are contained in the data member #m_locToGloMap). This
266  * operation is equivalent to the scatter operation
267  * \f$\boldsymbol{\hat{u}}_l=\mathcal{A}\boldsymbol{\hat{u}}_g\f$, where
268  * \f$\mathcal{A}\f$ is the
269  * \f$N_{\mathrm{eof}}\times N_{\mathrm{dof}}\f$ permutation matrix.
270  *
271  * @param inarray An array of size \f$N_\mathrm{dof}\f$
272  * containing the global degrees of freedom
273  * \f$\boldsymbol{x}_g\f$.
274  * @param outarray The resulting local degrees of freedom
275  * \f$\boldsymbol{x}_l\f$ will be stored in this
276  * array of size \f$N_\mathrm{eof}\f$.
277  */
279  const Array<OneD, const NekDouble> &inarray,
280  Array<OneD,NekDouble> &outarray)
281  {
282  m_locToGloMap->GlobalToLocal(inarray,outarray);
283  }
284 
285 
286  /**
287  * This operation is evaluated as:
288  * \f{tabbing}
289  * \hspace{1cm} \= Do \= $e=$ $1, N_{\mathrm{el}}$ \\
290  * > > Do \= $i=$ $0,N_m^e-1$ \\
291  * > > > $\boldsymbol{\hat{u}}_g[\mbox{map}[e][i]] =
292  * \boldsymbol{\hat{u}}_g[\mbox{map}[e][i]]+\mbox{sign}[e][i] \cdot
293  * \boldsymbol{\hat{u}}^{e}[i]$\\
294  * > > continue\\
295  * > continue
296  * \f}
297  * where \a map\f$[e][i]\f$ is the mapping array and \a
298  * sign\f$[e][i]\f$ is an array of similar dimensions ensuring the
299  * correct modal connectivity between the different elements (both
300  * these arrays are contained in the data member #m_locToGloMap). This
301  * operation is equivalent to the gather operation
302  * \f$\boldsymbol{\hat{u}}_g=\mathcal{A}^{T}\boldsymbol{\hat{u}}_l\f$,
303  * where \f$\mathcal{A}\f$ is the
304  * \f$N_{\mathrm{eof}}\times N_{\mathrm{dof}}\f$ permutation matrix.
305  *
306  */
307  inline void ContField1D::Assemble()
308  {
309  m_locToGloMap->Assemble(m_coeffs,m_coeffs);
310  }
311 
312  /**
313  * This operation is evaluated as:
314  * \f{tabbing}
315  * \hspace{1cm} \= Do \= $e=$ $1, N_{\mathrm{el}}$ \\
316  * > > Do \= $i=$ $0,N_m^e-1$ \\
317  * > > > $\boldsymbol{\hat{u}}_g[\mbox{map}[e][i]] =
318  * \boldsymbol{\hat{u}}_g[\mbox{map}[e][i]]+\mbox{sign}[e][i] \cdot
319  * \boldsymbol{\hat{u}}^{e}[i]$\\
320  * > > continue\\
321  * > continue
322  * \f}
323  * where \a map\f$[e][i]\f$ is the mapping array and \a
324  * sign\f$[e][i]\f$ is an array of similar dimensions ensuring the
325  * correct modal connectivity between the different elements (both
326  * these arrays are contained in the data member #m_locToGloMap). This
327  * operation is equivalent to the gather operation
328  * \f$\boldsymbol{\hat{u}}_g=\mathcal{A}^{T}\boldsymbol{\hat{u}}_l\f$,
329  * where \f$\mathcal{A}\f$ is the
330  * \f$N_{\mathrm{eof}}\times N_{\mathrm{dof}}\f$ permutation matrix.
331  *
332  * @param inarray An array of size \f$N_\mathrm{eof}\f$
333  * containing the local degrees of freedom
334  * \f$\boldsymbol{x}_l\f$.
335  * @param outarray The resulting global degrees of freedom
336  * \f$\boldsymbol{x}_g\f$ will be stored in this
337  * array of size \f$N_\mathrm{dof}\f$.
338  */
340  const Array<OneD, const NekDouble> &inarray,
341  Array<OneD,NekDouble> &outarray)
342  {
343  m_locToGloMap->Assemble(inarray,outarray);
344  }
345 
346  inline const AssemblyMapCGSharedPtr&
348  {
349  return m_locToGloMap;
350  }
351 
352  } //end of namespace
353 } //end of namespace
354 
355 #endif // MULTIERGIONS_CONTSOLNFIELD1D_H