Nektar++
StdExpansion2D.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File StdExpansion2D.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: Daughter of StdExpansion. This class contains routine
32 // which are common to 2D expansion. Typically this inolves physiocal
33 // space operations.
34 //
35 ///////////////////////////////////////////////////////////////////////////////
36 
38 
39 #ifdef max
40 #undef max
41 #endif
42 
43 namespace Nektar
44 {
45  namespace StdRegions
46  {
47 
49  {
50  }
51 
53  const LibUtilities::BasisKey &Ba,
54  const LibUtilities::BasisKey &Bb):
55  StdExpansion(numcoeffs,2, Ba, Bb)
56  {
57  }
58 
60  StdExpansion(T)
61  {
62  }
63 
65  {
66  }
67 
68  //----------------------------
69  // Differentiation Methods
70  //----------------------------
72  Array<OneD, NekDouble> &outarray_d0,
73  Array<OneD, NekDouble> &outarray_d1)
74  {
75  int nquad0 = m_base[0]->GetNumPoints();
76  int nquad1 = m_base[1]->GetNumPoints();
77 
78  if (outarray_d0.num_elements() > 0) // calculate du/dx_0
79  {
80  DNekMatSharedPtr D0 = m_base[0]->GetD();
81  if(inarray.data() == outarray_d0.data())
82  {
83  Array<OneD, NekDouble> wsp(nquad0 * nquad1);
84  Vmath::Vcopy(nquad0 * nquad1,inarray.get(),1,wsp.get(),1);
85  Blas::Dgemm('N', 'N', nquad0, nquad1, nquad0, 1.0,
86  &(D0->GetPtr())[0], nquad0, &wsp[0], nquad0, 0.0,
87  &outarray_d0[0], nquad0);
88  }
89  else
90  {
91  Blas::Dgemm('N', 'N', nquad0, nquad1, nquad0, 1.0,
92  &(D0->GetPtr())[0], nquad0, &inarray[0], nquad0, 0.0,
93  &outarray_d0[0], nquad0);
94  }
95  }
96 
97  if (outarray_d1.num_elements() > 0) // calculate du/dx_1
98  {
99  DNekMatSharedPtr D1 = m_base[1]->GetD();
100  if(inarray.data() == outarray_d1.data())
101  {
102  Array<OneD, NekDouble> wsp(nquad0 * nquad1);
103  Vmath::Vcopy(nquad0 * nquad1,inarray.get(),1,wsp.get(),1);
104  Blas:: Dgemm('N', 'T', nquad0, nquad1, nquad1, 1.0, &wsp[0], nquad0,
105  &(D1->GetPtr())[0], nquad1, 0.0, &outarray_d1[0], nquad0);
106  }
107  else
108  {
109  Blas:: Dgemm('N', 'T', nquad0, nquad1, nquad1, 1.0, &inarray[0], nquad0,
110  &(D1->GetPtr())[0], nquad1, 0.0, &outarray_d1[0], nquad0);
111  }
112  }
113  }
114 
116  {
117  Array<OneD, NekDouble> coll(2);
119 
120  ASSERTL2(coords[0] > -1 - NekConstants::kNekZeroTol, "coord[0] < -1");
121  ASSERTL2(coords[0] < 1 + NekConstants::kNekZeroTol, "coord[0] > 1");
122  ASSERTL2(coords[1] > -1 - NekConstants::kNekZeroTol, "coord[1] < -1");
123  ASSERTL2(coords[1] < 1 + NekConstants::kNekZeroTol, "coord[1] > 1");
124 
125  LocCoordToLocCollapsed(coords,coll);
126 
127  I[0] = m_base[0]->GetI(coll);
128  I[1] = m_base[1]->GetI(coll+1);
129 
130  return v_PhysEvaluate(I,physvals);
131  }
132 
135  const Array<OneD, const NekDouble> &physvals)
136  {
137  NekDouble val;
138  int i;
139  int nq0 = m_base[0]->GetNumPoints();
140  int nq1 = m_base[1]->GetNumPoints();
141  Array<OneD, NekDouble> wsp1(nq1);
142 
143  // interpolate first coordinate direction
144  for (i = 0; i < nq1;++i)
145  {
146  wsp1[i] = Blas::Ddot(nq0, &(I[0]->GetPtr())[0], 1,
147  &physvals[i * nq0], 1);
148  }
149 
150  // interpolate in second coordinate direction
151  val = Blas::Ddot(nq1, I[1]->GetPtr(), 1, wsp1, 1);
152 
153  return val;
154  }
155 
156  //////////////////////////////
157  // Integration Methods
158  //////////////////////////////
159 
163  {
164  int i;
165  NekDouble Int = 0.0;
166  int nquad0 = m_base[0]->GetNumPoints();
167  int nquad1 = m_base[1]->GetNumPoints();
168  Array<OneD, NekDouble> tmp(nquad0 * nquad1);
169 
170  // multiply by integration constants
171  for (i = 0; i < nquad1; ++i)
172  {
173  Vmath::Vmul(nquad0, &inarray[0] + i*nquad0, 1, w0.get(),
174  1, &tmp[0] + i*nquad0, 1);
175  }
176 
177  for (i = 0; i < nquad0; ++i)
178  {
179  Vmath::Vmul(nquad1, &tmp[0]+ i, nquad0, w1.get(), 1,
180  &tmp[0] + i, nquad0);
181  }
182  Int = Vmath::Vsum(nquad0 * nquad1, tmp, 1);
183 
184  return Int;
185  }
186 
188  const Array<OneD, const NekDouble>& base0,
189  const Array<OneD, const NekDouble>& base1,
190  const Array<OneD, const NekDouble>& inarray,
191  Array<OneD, NekDouble> &outarray,
193  bool doCheckCollDir0,
194  bool doCheckCollDir1)
195  {
196  v_BwdTrans_SumFacKernel(base0, base1, inarray, outarray, wsp, doCheckCollDir0, doCheckCollDir1);
197  }
198 
200  const Array<OneD, const NekDouble>& base0,
201  const Array<OneD, const NekDouble>& base1,
202  const Array<OneD, const NekDouble>& inarray,
203  Array<OneD, NekDouble> &outarray,
205  bool doCheckCollDir0,
206  bool doCheckCollDir1)
207  {
208  v_IProductWRTBase_SumFacKernel(base0, base1, inarray, outarray, wsp, doCheckCollDir0, doCheckCollDir1);
209  }
210 
212  const Array<OneD, const NekDouble> &inarray,
213  Array<OneD,NekDouble> &outarray,
214  const StdRegions::StdMatrixKey &mkey)
215  {
216  if (mkey.GetNVarCoeff() == 0
218  {
219  using std::max;
220 
221  // This implementation is only valid when there are no
222  // coefficients associated to the Laplacian operator
223  int nquad0 = m_base[0]->GetNumPoints();
224  int nquad1 = m_base[1]->GetNumPoints();
225  int nqtot = nquad0*nquad1;
226  int nmodes0 = m_base[0]->GetNumModes();
227  int nmodes1 = m_base[1]->GetNumModes();
228  int wspsize = max(max(max(nqtot,m_ncoeffs),nquad1*nmodes0),nquad0*nmodes1);
229 
230  const Array<OneD, const NekDouble>& base0 = m_base[0]->GetBdata();
231  const Array<OneD, const NekDouble>& base1 = m_base[1]->GetBdata();
232 
233  // Allocate temporary storage
234  Array<OneD,NekDouble> wsp0(4*wspsize); // size wspsize
235  Array<OneD,NekDouble> wsp1(wsp0+wspsize); // size 3*wspsize
236 
237  if(!(m_base[0]->Collocation() && m_base[1]->Collocation()))
238  {
239  // LAPLACIAN MATRIX OPERATION
240  // wsp0 = u = B * u_hat
241  // wsp1 = du_dxi1 = D_xi1 * wsp0 = D_xi1 * u
242  // wsp2 = du_dxi2 = D_xi2 * wsp0 = D_xi2 * u
243  BwdTrans_SumFacKernel(base0,base1,inarray,wsp0,wsp1,true,true);
244  LaplacianMatrixOp_MatFree_Kernel(wsp0, outarray, wsp1);
245  }
246  else
247  {
248  LaplacianMatrixOp_MatFree_Kernel(inarray, outarray, wsp1);
249  }
250  }
251  else
252  {
254  inarray,outarray,mkey);
255  }
256  }
257 
258 
259 
261  const Array<OneD, const NekDouble> &inarray,
262  Array<OneD,NekDouble> &outarray,
263  const StdRegions::StdMatrixKey &mkey)
264  {
265  if (mkey.GetNVarCoeff() == 0
267  {
268  using std::max;
269 
270  int nquad0 = m_base[0]->GetNumPoints();
271  int nquad1 = m_base[1]->GetNumPoints();
272  int nqtot = nquad0*nquad1;
273  int nmodes0 = m_base[0]->GetNumModes();
274  int nmodes1 = m_base[1]->GetNumModes();
275  int wspsize = max(max(max(nqtot,m_ncoeffs),nquad1*nmodes0),
276  nquad0*nmodes1);
277  NekDouble lambda =
279 
280  const Array<OneD, const NekDouble>& base0 = m_base[0]->GetBdata();
281  const Array<OneD, const NekDouble>& base1 = m_base[1]->GetBdata();
282 
283  // Allocate temporary storage
284  Array<OneD,NekDouble> wsp0(5*wspsize); // size wspsize
285  Array<OneD,NekDouble> wsp1(wsp0 + wspsize); // size wspsize
286  Array<OneD,NekDouble> wsp2(wsp0 + 2*wspsize);// size 3*wspsize
287 
288  if (!(m_base[0]->Collocation() && m_base[1]->Collocation()))
289  {
290  // MASS MATRIX OPERATION
291  // The following is being calculated:
292  // wsp0 = B * u_hat = u
293  // wsp1 = W * wsp0
294  // outarray = B^T * wsp1 = B^T * W * B * u_hat = M * u_hat
295  BwdTrans_SumFacKernel (base0, base1, inarray,
296  wsp0, wsp2,true,true);
297  MultiplyByQuadratureMetric (wsp0, wsp1);
298  IProductWRTBase_SumFacKernel(base0, base1, wsp1, outarray,
299  wsp2, true, true);
300 
301  LaplacianMatrixOp_MatFree_Kernel(wsp0, wsp1, wsp2);
302  }
303  else
304  {
305  MultiplyByQuadratureMetric(inarray,outarray);
306  LaplacianMatrixOp_MatFree_Kernel(inarray, wsp1, wsp2);
307  }
308 
309  // outarray = lambda * outarray + wsp1
310  // = (lambda * M + L ) * u_hat
311  Vmath::Svtvp(m_ncoeffs, lambda, &outarray[0], 1,
312  &wsp1[0], 1, &outarray[0], 1);
313  }
314  else
315  {
317  inarray,outarray,mkey);
318  }
319  }
320 
321  } //end namespace
322 } //end namespace
void HelmholtzMatrixOp_MatFree_GenericImpl(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void LaplacianMatrixOp_MatFree_GenericImpl(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
void BwdTrans_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0=true, bool doCheckCollDir1=true)
void MultiplyByQuadratureMetric(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray)
Definition: StdExpansion.h:945
virtual void v_IProductWRTBase_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1)=0
void Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
svtvp (scalar times vector plus vector): z = alpha*x + y
Definition: Vmath.cpp:488
void LocCoordToLocCollapsed(const Array< OneD, const NekDouble > &xi, Array< OneD, NekDouble > &eta)
Convert local cartesian coordinate xi into local collapsed coordinates eta.
NekDouble Integral(const Array< OneD, const NekDouble > &inarray, const Array< OneD, const NekDouble > &w0, const Array< OneD, const NekDouble > &w1)
virtual void v_LaplacianMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)
std::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:69
static void Dgemm(const char &transa, const char &transb, const int &m, const int &n, const int &k, const double &alpha, const double *a, const int &lda, const double *b, const int &ldb, const double &beta, double *c, const int &ldc)
BLAS level 3: Matrix-matrix multiply C = A x B where A[m x n], B[n x k], C[m x k].
Definition: Blas.hpp:213
void PhysTensorDeriv(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray_d0, Array< OneD, NekDouble > &outarray_d1)
Calculate the 2D derivative in the local tensor/collapsed coordinate at the physical points...
The base class for all shapes.
Definition: StdExpansion.h:68
static const NekDouble kNekZeroTol
virtual NekDouble v_PhysEvaluate(const Array< OneD, const NekDouble > &coords, const Array< OneD, const NekDouble > &physvals)
This function evaluates the expansion at a single (arbitrary) point of the domain.
double NekDouble
void LaplacianMatrixOp_MatFree_Kernel(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp)
NekDouble GetConstFactor(const ConstFactorType &factor) const
Definition: StdMatrixKey.h:121
static double Ddot(const int &n, const double *x, const int &incx, const double *y, const int &incy)
BLAS level 1: output = .
Definition: Blas.hpp:140
bool ConstFactorExists(const ConstFactorType &factor) const
Definition: StdMatrixKey.h:130
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:274
void IProductWRTBase_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0=true, bool doCheckCollDir1=true)
T Vsum(int n, const T *x, const int incx)
Subtract return sum(x)
Definition: Vmath.cpp:740
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
Describes the specification for a Basis.
Definition: Basis.h:49
virtual void v_BwdTrans_SumFacKernel(const Array< OneD, const NekDouble > &base0, const Array< OneD, const NekDouble > &base1, const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, Array< OneD, NekDouble > &wsp, bool doCheckCollDir0, bool doCheckCollDir1)=0
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
virtual void v_HelmholtzMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)