Nektar++
StdPointExp.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File $Source:
4 // /usr/sci/projects/Nektar/cvs/Nektar++/library/LocalRegions/PointExp.cpp,v $
5 //
6 // For more information, please see: http://www.nektar.info
7 //
8 // The MIT License
9 //
10 // Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
11 // Department of Aeronautics, Imperial College London (UK), and Scientific
12 // Computing and Imaging Institute, University of Utah (USA).
13 //
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: Definition of a Point expansion
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
36 #include <boost/core/ignore_unused.hpp>
37 
38 #include <StdRegions/StdPointExp.h>
39 
40 using namespace std;
41 
42 namespace Nektar
43 {
44 namespace StdRegions
45 {
46 
47 StdPointExp::StdPointExp()
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 /** \brief Copy Constructor */
58 
61 {
62 }
63 
65 {
66 }
67 
69 {
70  return LibUtilities::ePoint;
71 }
72 
74  Array<OneD, NekDouble> &coords_1,
75  Array<OneD, NekDouble> &coords_2)
76 {
77  boost::ignore_unused(coords_1, coords_2);
78  Blas::Dcopy(GetNumPoints(0), (m_base[0]->GetZ()).get(), 1, &coords_0[0], 1);
79 }
80 
82  Array<OneD, NekDouble> &outarray)
83 {
84  int nquad = m_base[0]->GetNumPoints();
85 
86  if (m_base[0]->Collocation())
87  {
88  Vmath::Vcopy(nquad, inarray, 1, outarray, 1);
89  }
90  else
91  {
92  Blas::Dgemv('N', nquad, m_base[0]->GetNumModes(), 1.0,
93  (m_base[0]->GetBdata()).get(), nquad, &inarray[0], 1, 0.0,
94  &outarray[0], 1);
95  }
96 }
97 
99  Array<OneD, NekDouble> &outarray)
100 {
101  if (m_base[0]->Collocation())
102  {
103  Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
104  }
105  else
106  {
107  v_IProductWRTBase(inarray, outarray);
108 
109  // get Mass matrix inverse
110  StdMatrixKey masskey(eInvMass, v_DetShapeType(), *this);
111  DNekMatSharedPtr matsys = GetStdMatrix(masskey);
112 
113  NekVector<NekDouble> in(m_ncoeffs, outarray, eCopy);
114  NekVector<NekDouble> out(m_ncoeffs, outarray, eWrapper);
115 
116  out = (*matsys) * in;
117  }
118 }
119 
121  const Array<OneD, const NekDouble> &inarray,
122  Array<OneD, NekDouble> &outarray)
123 {
124  if (m_base[0]->Collocation())
125  {
126  Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
127  }
128  else
129  {
130  int nInteriorDofs = m_ncoeffs - 2;
131  int offset = 0;
132 
133  switch (m_base[0]->GetBasisType())
134  {
136  {
137  offset = 1;
138  }
139  break;
142  {
143  offset = 2;
144  }
145  break;
146  default:
147  ASSERTL0(false, "This type of FwdTrans is not defined for this "
148  "shapex type");
149  }
150 
151  fill(outarray.get(), outarray.get() + m_ncoeffs, 0.0);
152 
153  outarray[GetVertexMap(0)] = inarray[0];
154  outarray[GetVertexMap(1)] = inarray[m_base[0]->GetNumPoints() - 1];
155 
156  if (m_ncoeffs > 2)
157  {
158  // ideally, we would like to have tmp0 to be replaced by
159  // outarray (currently MassMatrixOp does not allow aliasing)
162 
163  StdMatrixKey masskey(eMass, v_DetShapeType(), *this);
164  MassMatrixOp(outarray, tmp0, masskey);
165  v_IProductWRTBase(inarray, tmp1);
166 
167  Vmath::Vsub(m_ncoeffs, tmp1, 1, tmp0, 1, tmp1, 1);
168 
169  // get Mass matrix inverse (only of interior DOF)
170  DNekMatSharedPtr matsys =
171  (m_stdStaticCondMatrixManager[masskey])->GetBlock(1, 1);
172 
173  Blas::Dgemv('N', nInteriorDofs, nInteriorDofs, 1.0,
174  &(matsys->GetPtr())[0], nInteriorDofs,
175  tmp1.get() + offset, 1, 0.0, outarray.get() + offset,
176  1);
177  }
178  }
179 }
180 
182  Array<OneD, NekDouble> &outarray)
183 {
184  v_BwdTrans(inarray, outarray);
185 }
186 
187 // Inner product
189  const Array<OneD, const NekDouble> &inarray,
190  Array<OneD, NekDouble> &outarray,
191  int coll_check)
192 {
193  int nquad = m_base[0]->GetNumPoints();
194  Array<OneD, NekDouble> tmp(nquad);
195  Array<OneD, const NekDouble> z = m_base[0]->GetZ();
196  Array<OneD, const NekDouble> w = m_base[0]->GetW();
197 
198  Vmath::Vmul(nquad, inarray, 1, w, 1, tmp, 1);
199 
200  if (coll_check && m_base[0]->Collocation())
201  {
202  Vmath::Vcopy(nquad, tmp, 1, outarray, 1);
203  }
204  else
205  {
206  Blas::Dgemv('T', nquad, m_ncoeffs, 1.0, base.get(), nquad, &tmp[0], 1,
207  0.0, outarray.get(), 1);
208  }
209 }
210 
212  Array<OneD, NekDouble> &outarray)
213 {
214  v_IProductWRTBase(m_base[0]->GetBdata(), inarray, outarray, 1);
215 }
216 
218  const int dir, const Array<OneD, const NekDouble> &inarray,
219  Array<OneD, NekDouble> &outarray)
220 {
221  boost::ignore_unused(dir);
222  ASSERTL1(dir >= 0 && dir < 1, "input dir is out of range");
223  v_IProductWRTBase(m_base[0]->GetDbdata(), inarray, outarray, 1);
224 }
225 
227  const Array<OneD, const NekDouble> &inarray,
228  Array<OneD, NekDouble> &outarray, bool multiplybyweights)
229 {
230  boost::ignore_unused(multiplybyweights);
231  v_IProductWRTBase(m_base[0]->GetBdata(), inarray, outarray, 1);
232 }
233 
235 {
236  DNekMatSharedPtr Mat;
237  MatrixType mattype;
238 
239  switch (mattype = mkey.GetMatrixType())
240  {
241  case eFwdTrans:
242  {
243  Mat =
245  StdMatrixKey iprodkey(eIProductWRTBase, v_DetShapeType(), *this);
246  DNekMat &Iprod = *GetStdMatrix(iprodkey);
247  StdMatrixKey imasskey(eInvMass, v_DetShapeType(), *this);
248  DNekMat &Imass = *GetStdMatrix(imasskey);
249 
250  (*Mat) = Imass * Iprod;
251  }
252  break;
253  default:
254  {
256  }
257  break;
258  }
259 
260  return Mat;
261 }
262 
264 {
265  return v_GenMatrix(mkey);
266 }
267 
268 } // namespace StdRegions
269 } // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
Describes the specification for a Basis.
Definition: Basis.h:50
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
The base class for all shapes.
Definition: StdExpansion.h:71
LibUtilities::BasisType GetBasisType(const int dir) const
This function returns the type of basis used in the dir direction.
Definition: StdExpansion.h:163
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:611
void MassMatrixOp(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
Definition: StdExpansion.h:760
int GetVertexMap(const int localVertexId, bool useCoeffPacking=false)
Definition: StdExpansion.h:687
DNekMatSharedPtr CreateGeneralMatrix(const StdMatrixKey &mkey)
this function generates the mass matrix
LibUtilities::NekManager< StdMatrixKey, DNekBlkMat, StdMatrixKey::opLess > m_stdStaticCondMatrixManager
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Definition: StdExpansion.h:226
Array< OneD, LibUtilities::BasisSharedPtr > m_base
MatrixType GetMatrixType() const
Definition: StdMatrixKey.h:85
virtual void v_IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
virtual void v_GetCoords(Array< OneD, NekDouble > &coords_0, Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2)
Definition: StdPointExp.cpp:73
virtual void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdPointExp.cpp:81
virtual DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey)
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_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)
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.
Definition: StdPointExp.cpp:98
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true)
virtual LibUtilities::ShapeType v_DetShapeType() const
Definition: StdPointExp.cpp:68
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:246
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:147
@ eModified_B
Principle Modified Functions .
Definition: BasisType.h:51
@ eGLL_Lagrange
Lagrange for SEM basis .
Definition: BasisType.h:58
@ eModified_A
Principle Modified Functions .
Definition: BasisType.h:50
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:75
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:209
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255
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:419