Nektar++
StdPointExp.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File $Source: /usr/sci/projects/Nektar/cvs/Nektar++/library/LocalRegions/PointExp.cpp,v $
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: Definition of a Point expansion
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
35 #include <boost/core/ignore_unused.hpp>
36 
37 #include <StdRegions/StdPointExp.h>
38 
39 using namespace std;
40 
41 namespace Nektar
42 {
43  namespace StdRegions
44  {
45 
46  StdPointExp::StdPointExp()
47  {
48  }
49 
50 
51  StdPointExp::StdPointExp(const LibUtilities::BasisKey &Ba):
52  StdExpansion(Ba.GetNumModes(), 1, Ba),
53  StdExpansion0D(Ba.GetNumModes(),Ba)
54  {
55  }
56 
57 
58  /** \brief Copy Constructor */
59 
61  StdExpansion(T),
63  {
64  }
65 
66 
67 
69  {
70  }
71 
72 
74  {
75  return LibUtilities::ePoint;
76  }
77 
78 
80  Array<OneD, NekDouble> &coords_0,
81  Array<OneD, NekDouble> &coords_1,
82  Array<OneD, NekDouble> &coords_2)
83  {
84  boost::ignore_unused(coords_1, coords_2);
85  Blas::Dcopy(GetNumPoints(0),(m_base[0]->GetZ()).get(),
86  1,&coords_0[0],1);
87  }
88 
89 
91  const Array<OneD, const NekDouble>& inarray,
92  Array<OneD, NekDouble> &outarray)
93  {
94  int nquad = m_base[0]->GetNumPoints();
95 
96  if(m_base[0]->Collocation())
97  {
98  Vmath::Vcopy(nquad, inarray, 1, outarray, 1);
99  }
100  else
101  {
102  Blas::Dgemv('N',nquad,m_base[0]->GetNumModes(),1.0,
103  (m_base[0]->GetBdata()).get(),
104  nquad,&inarray[0],1,0.0,&outarray[0],1);
105  }
106  }
107 
108 
110  Array<OneD, NekDouble> &outarray)
111  {
112  if(m_base[0]->Collocation())
113  {
114  Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
115  }
116  else
117  {
118  v_IProductWRTBase(inarray,outarray);
119 
120  // get Mass matrix inverse
121  StdMatrixKey masskey(eInvMass,v_DetShapeType(),*this);
122  DNekMatSharedPtr matsys = GetStdMatrix(masskey);
123 
124  NekVector<NekDouble> in(m_ncoeffs,outarray,eCopy);
126 
127  out = (*matsys)*in;
128  }
129  }
130 
132  const Array<OneD, const NekDouble>& inarray,
133  Array<OneD, NekDouble> &outarray)
134  {
135  if(m_base[0]->Collocation())
136  {
137  Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
138  }
139  else
140  {
141  int nInteriorDofs = m_ncoeffs-2;
142  int offset = 0;
143 
144  switch(m_base[0]->GetBasisType())
145  {
147  {
148  offset = 1;
149  }
150  break;
153  {
154  offset = 2;
155  }
156  break;
157  default:
158  ASSERTL0(false,"This type of FwdTrans is not defined for this shapex type");
159  }
160 
161  fill(outarray.get(), outarray.get()+m_ncoeffs, 0.0 );
162 
163  outarray[GetVertexMap(0)] = inarray[0];
164  outarray[GetVertexMap(1)] = inarray[m_base[0]->GetNumPoints()-1];
165 
166  if(m_ncoeffs>2)
167  {
168  // ideally, we would like to have tmp0 to be replaced by
169  // outarray (currently MassMatrixOp does not allow aliasing)
172 
173  StdMatrixKey masskey(eMass,v_DetShapeType(),*this);
174  MassMatrixOp(outarray,tmp0,masskey);
175  v_IProductWRTBase(inarray,tmp1);
176 
177  Vmath::Vsub(m_ncoeffs, tmp1, 1, tmp0, 1, tmp1, 1);
178 
179  // get Mass matrix inverse (only of interior DOF)
180  DNekMatSharedPtr matsys = (m_stdStaticCondMatrixManager[masskey])->GetBlock(1,1);
181 
182  Blas::Dgemv('N',nInteriorDofs,nInteriorDofs,1.0, &(matsys->GetPtr())[0],
183  nInteriorDofs,tmp1.get()+offset,1,0.0,outarray.get()+offset,1);
184  }
185  }
186 
187  }
188 
189 
191  Array<OneD, NekDouble> &outarray)
192  {
193  v_BwdTrans(inarray, outarray);
194  }
195 
196 //Inner product
198  const Array<OneD, const NekDouble>& inarray,
199  Array<OneD, NekDouble> &outarray,
200  int coll_check)
201  {
202  int nquad = m_base[0]->GetNumPoints();
203  Array<OneD, NekDouble> tmp(nquad);
204  Array<OneD, const NekDouble> z = m_base[0]->GetZ();
205  Array<OneD, const NekDouble> w = m_base[0]->GetW();
206 
207  Vmath::Vmul(nquad, inarray, 1, w, 1, tmp, 1);
208 
209  if(coll_check&&m_base[0]->Collocation())
210  {
211  Vmath::Vcopy(nquad, tmp, 1, outarray, 1);
212  }
213  else
214  {
215  Blas::Dgemv('T',nquad,m_ncoeffs,1.0,base.get(),nquad,
216  &tmp[0],1,0.0,outarray.get(),1);
217  }
218  }
219 
221  Array<OneD, NekDouble> &outarray)
222  {
223  v_IProductWRTBase(m_base[0]->GetBdata(),inarray,outarray,1);
224  }
225 
227  const int dir,
228  const Array<OneD, const NekDouble>& inarray,
229  Array<OneD, NekDouble> & outarray)
230  {
231  boost::ignore_unused(dir);
232  ASSERTL1(dir >= 0 && dir < 1,"input dir is out of range");
233  v_IProductWRTBase(m_base[0]->GetDbdata(),inarray,outarray,1);
234  }
235 
237  const Array<OneD, const NekDouble>& inarray,
238  Array<OneD, NekDouble> &outarray,
239  bool multiplybyweights)
240  {
241  boost::ignore_unused(multiplybyweights);
242  v_IProductWRTBase(m_base[0]->GetBdata(),inarray,outarray,1);
243  }
244 
245 
247  {
248  DNekMatSharedPtr Mat;
249  MatrixType mattype;
250 
251  switch(mattype = mkey.GetMatrixType())
252  {
253  case eFwdTrans:
254  {
256  StdMatrixKey iprodkey(eIProductWRTBase,v_DetShapeType(),*this);
257  DNekMat &Iprod = *GetStdMatrix(iprodkey);
258  StdMatrixKey imasskey(eInvMass,v_DetShapeType(),*this);
259  DNekMat &Imass = *GetStdMatrix(imasskey);
260 
261  (*Mat) = Imass*Iprod;
262 
263  }
264  break;
265  default:
266  {
268  }
269  break;
270  }
271 
272  return Mat;
273  }
274 
275 
277  {
278  return v_GenMatrix(mkey);
279  }
280 
281 
282  }
283 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
virtual LibUtilities::ShapeType v_DetShapeType() const
Definition: StdPointExp.cpp:73
void MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:974
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Definition: StdExpansion.h:228
virtual void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdPointExp.cpp:90
Principle Modified Functions .
Definition: BasisType.h:48
STL namespace.
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Calculates the inner product of a given function f with the different modes of the expansion...
virtual void v_GetCoords(Array< OneD, NekDouble > &coords_0, Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2)
Definition: StdPointExp.cpp:79
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:69
LibUtilities::NekManager< StdMatrixKey, DNekBlkMat, StdMatrixKey::opLess > m_stdStaticCondMatrixManager
MatrixType GetMatrixType() const
Definition: StdMatrixKey.h:81
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:714
virtual void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Transform a given function from physical quadrature space to coefficient space.
virtual void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_FwdTrans_BndConstrained(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
The base class for all shapes.
Definition: StdExpansion.h:68
int GetVertexMap(const int localVertexId, bool useCoeffPacking=false)
Definition: StdExpansion.h:822
Principle Modified Functions .
Definition: BasisType.h:49
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
DNekMatSharedPtr CreateGeneralMatrix(const StdMatrixKey &mkey)
this function generates the mass matrix
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
virtual DNekMatSharedPtr v_CreateStdMatrix(const StdMatrixKey &mkey)
static void Dgemv(const char &trans, const int &m, const int &n, const double &alpha, const double *a, const int &lda, const double *x, const int &incx, const double &beta, double *y, const int &incy)
BLAS level 2: Matrix vector multiply y = A x where A[m x n].
Definition: Blas.hpp:168
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.cpp:346
virtual void v_IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:164
Lagrange for SEM basis .
Definition: BasisType.h:54
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:250
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064
static void Dcopy(const int &n, const double *x, const int &incx, double *y, const int &incy)
BLAS level 1: Copy x to y.
Definition: Blas.hpp:103
Describes the specification for a Basis.
Definition: Basis.h:49
virtual DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey)
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.cpp:186