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
36
37using namespace std;
38
39namespace Nektar::StdRegions
40{
41
43 : StdExpansion(Ba.GetNumModes(), 1, Ba),
44 StdExpansion0D(Ba.GetNumModes(), Ba)
45{
46}
47
49{
51}
52
54 [[maybe_unused]] Array<OneD, NekDouble> &coords_1,
55 [[maybe_unused]] Array<OneD, NekDouble> &coords_2)
56{
57 Blas::Dcopy(GetNumPoints(0), (m_base[0]->GetZ()).get(), 1, &coords_0[0], 1);
58}
59
61 Array<OneD, NekDouble> &outarray)
62{
63 int nquad = m_base[0]->GetNumPoints();
64
65 if (m_base[0]->Collocation())
66 {
67 Vmath::Vcopy(nquad, inarray, 1, outarray, 1);
68 }
69 else
70 {
71 Blas::Dgemv('N', nquad, m_base[0]->GetNumModes(), 1.0,
72 (m_base[0]->GetBdata()).get(), nquad, &inarray[0], 1, 0.0,
73 &outarray[0], 1);
74 }
75}
76
78 Array<OneD, NekDouble> &outarray)
79{
80 if (m_base[0]->Collocation())
81 {
82 Vmath::Vcopy(m_ncoeffs, inarray, 1, outarray, 1);
83 }
84 else
85 {
86 v_IProductWRTBase(inarray, outarray);
87
88 // get Mass matrix inverse
89 StdMatrixKey masskey(eInvMass, v_DetShapeType(), *this);
90 DNekMatSharedPtr matsys = GetStdMatrix(masskey);
91
94
95 out = (*matsys) * in;
96 }
97}
98
100 Array<OneD, NekDouble> &outarray)
101{
102 v_BwdTrans(inarray, outarray);
103}
104
105// Inner product
107 const Array<OneD, const NekDouble> &inarray,
108 Array<OneD, NekDouble> &outarray,
109 int coll_check)
110{
111 int nquad = m_base[0]->GetNumPoints();
112 Array<OneD, NekDouble> tmp(nquad);
115
116 Vmath::Vmul(nquad, inarray, 1, w, 1, tmp, 1);
117
118 if (coll_check && m_base[0]->Collocation())
119 {
120 Vmath::Vcopy(nquad, tmp, 1, outarray, 1);
121 }
122 else
123 {
124 Blas::Dgemv('T', nquad, m_ncoeffs, 1.0, base.get(), nquad, &tmp[0], 1,
125 0.0, outarray.get(), 1);
126 }
127}
128
130 Array<OneD, NekDouble> &outarray)
131{
132 v_IProductWRTBase(m_base[0]->GetBdata(), inarray, outarray, 1);
133}
134
136 [[maybe_unused]] const int dir, const Array<OneD, const NekDouble> &inarray,
137 Array<OneD, NekDouble> &outarray)
138{
139 ASSERTL1(dir >= 0 && dir < 1, "input dir is out of range");
140 v_IProductWRTBase(m_base[0]->GetDbdata(), inarray, outarray, 1);
141}
142
144 const Array<OneD, const NekDouble> &inarray,
145 Array<OneD, NekDouble> &outarray, [[maybe_unused]] bool multiplybyweights)
146{
147 v_IProductWRTBase(m_base[0]->GetBdata(), inarray, outarray, 1);
148}
149
151{
153 MatrixType mattype;
154
155 switch (mattype = mkey.GetMatrixType())
156 {
157 case eFwdTrans:
158 {
159 Mat =
161 StdMatrixKey iprodkey(eIProductWRTBase, v_DetShapeType(), *this);
162 DNekMat &Iprod = *GetStdMatrix(iprodkey);
163 StdMatrixKey imasskey(eInvMass, v_DetShapeType(), *this);
164 DNekMat &Imass = *GetStdMatrix(imasskey);
165
166 (*Mat) = Imass * Iprod;
167 }
168 break;
169 default:
170 {
172 }
173 break;
174 }
175
176 return Mat;
177}
178
180{
181 return v_GenMatrix(mkey);
182}
183
184} // namespace Nektar::StdRegions
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
Describes the specification for a Basis.
Definition: Basis.h:45
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:65
DNekMatSharedPtr GetStdMatrix(const StdMatrixKey &mkey)
Definition: StdExpansion.h:603
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:218
Array< OneD, LibUtilities::BasisSharedPtr > m_base
MatrixType GetMatrixType() const
Definition: StdMatrixKey.h:83
DNekMatSharedPtr v_CreateStdMatrix(const StdMatrixKey &mkey) override
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.
void v_IProductWRTDerivBase(const int dir, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
DNekMatSharedPtr v_GenMatrix(const StdMatrixKey &mkey) override
void v_BwdTrans_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: StdPointExp.cpp:99
LibUtilities::ShapeType v_DetShapeType() const override
Definition: StdPointExp.cpp:48
void v_BwdTrans(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray) override
Definition: StdPointExp.cpp:60
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:77
void v_GetCoords(Array< OneD, NekDouble > &coords_0, Array< OneD, NekDouble > &coords_1, Array< OneD, NekDouble > &coords_2) override
Definition: StdPointExp.cpp:53
void v_IProductWRTBase_SumFac(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, bool multiplybyweights=true) override
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:211
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:128
std::vector< double > w(NPUPPER)
std::vector< double > z(NPUPPER)
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.hpp:72
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825