Nektar++
StdPointExp.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: StdPointExp.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: Definition of a Point expansion
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#include <boost/core/ignore_unused.hpp>
36
38
39using namespace std;
40
41namespace Nektar
42{
43namespace StdRegions
44{
45
47{
48}
49
51 : StdExpansion(Ba.GetNumModes(), 1, Ba),
52 StdExpansion0D(Ba.GetNumModes(), Ba)
53{
54}
55
56/** \brief Copy Constructor */
57
60{
61}
62
64{
65}
66
68{
70}
71
73 Array<OneD, NekDouble> &coords_1,
74 Array<OneD, NekDouble> &coords_2)
75{
76 boost::ignore_unused(coords_1, coords_2);
77 Blas::Dcopy(GetNumPoints(0), (m_base[0]->GetZ()).get(), 1, &coords_0[0], 1);
78}
79
81 Array<OneD, NekDouble> &outarray)
82{
83 int nquad = m_base[0]->GetNumPoints();
84
85 if (m_base[0]->Collocation())
86 {
87 Vmath::Vcopy(nquad, inarray, 1, outarray, 1);
88 }
89 else
90 {
91 Blas::Dgemv('N', nquad, m_base[0]->GetNumModes(), 1.0,
92 (m_base[0]->GetBdata()).get(), nquad, &inarray[0], 1, 0.0,
93 &outarray[0], 1);
94 }
95}
96
98 Array<OneD, NekDouble> &outarray)
99{
100 if (m_base[0]->Collocation())
101 {
102 Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
103 }
104 else
105 {
106 v_IProductWRTBase(inarray, outarray);
107
108 // get Mass matrix inverse
109 StdMatrixKey masskey(eInvMass, v_DetShapeType(), *this);
110 DNekMatSharedPtr matsys = GetStdMatrix(masskey);
111
112 NekVector<NekDouble> in(m_ncoeffs, outarray, eCopy);
114
115 out = (*matsys) * in;
116 }
117}
118
120 Array<OneD, NekDouble> &outarray)
121{
122 v_BwdTrans(inarray, outarray);
123}
124
125// Inner product
127 const Array<OneD, const NekDouble> &inarray,
128 Array<OneD, NekDouble> &outarray,
129 int coll_check)
130{
131 int nquad = m_base[0]->GetNumPoints();
132 Array<OneD, NekDouble> tmp(nquad);
135
136 Vmath::Vmul(nquad, inarray, 1, w, 1, tmp, 1);
137
138 if (coll_check && m_base[0]->Collocation())
139 {
140 Vmath::Vcopy(nquad, tmp, 1, outarray, 1);
141 }
142 else
143 {
144 Blas::Dgemv('T', nquad, m_ncoeffs, 1.0, base.get(), nquad, &tmp[0], 1,
145 0.0, outarray.get(), 1);
146 }
147}
148
150 Array<OneD, NekDouble> &outarray)
151{
152 v_IProductWRTBase(m_base[0]->GetBdata(), inarray, outarray, 1);
153}
154
156 const int dir, const Array<OneD, const NekDouble> &inarray,
157 Array<OneD, NekDouble> &outarray)
158{
159 boost::ignore_unused(dir);
160 ASSERTL1(dir >= 0 && dir < 1, "input dir is out of range");
161 v_IProductWRTBase(m_base[0]->GetDbdata(), inarray, outarray, 1);
162}
163
165 const Array<OneD, const NekDouble> &inarray,
166 Array<OneD, NekDouble> &outarray, bool multiplybyweights)
167{
168 boost::ignore_unused(multiplybyweights);
169 v_IProductWRTBase(m_base[0]->GetBdata(), inarray, outarray, 1);
170}
171
173{
175 MatrixType mattype;
176
177 switch (mattype = mkey.GetMatrixType())
178 {
179 case eFwdTrans:
180 {
181 Mat =
183 StdMatrixKey iprodkey(eIProductWRTBase, v_DetShapeType(), *this);
184 DNekMat &Iprod = *GetStdMatrix(iprodkey);
185 StdMatrixKey imasskey(eInvMass, v_DetShapeType(), *this);
186 DNekMat &Imass = *GetStdMatrix(imasskey);
187
188 (*Mat) = Imass * Iprod;
189 }
190 break;
191 default:
192 {
194 }
195 break;
196 }
197
198 return Mat;
199}
200
202{
203 return v_GenMatrix(mkey);
204}
205
206} // namespace StdRegions
207} // namespace Nektar
#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:47
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
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:609
DNekMatSharedPtr CreateGeneralMatrix(const StdMatrixKey &mkey)
this function generates the mass matrix
int GetNumPoints(const int dir) const
This function returns the number of quadrature points in the dir direction.
Definition: StdExpansion.h:224
Array< OneD, LibUtilities::BasisSharedPtr > m_base
MatrixType GetMatrixType() const
Definition: StdMatrixKey.h:87
virtual DNekMatSharedPtr v_CreateStdMatrix(const StdMatrixKey &mkey) override
virtual void v_IProductWRTBase(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Calculates the inner product of a given function f with the different modes of the expansion.
virtual void v_IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
virtual DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey) override
virtual void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
virtual LibUtilities::ShapeType v_DetShapeType() const override
Definition: StdPointExp.cpp:67
virtual void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: StdPointExp.cpp:80
virtual void v_FwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Transform a given function from physical quadrature space to coefficient space.
Definition: StdPointExp.cpp:97
virtual void v_GetCoords(Array< OneD, NekDouble > &coords_0, Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2) override
Definition: StdPointExp.cpp:72
virtual void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true) override
virtual ~StdPointExp() override
Definition: StdPointExp.cpp:63
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 = alpha A x plus beta y where A[m x n].
Definition: Blas.hpp:213
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:130
std::vector< double > w(NPUPPER)
std::vector< double > z(NPUPPER)
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
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:207
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1191