Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 // License for the specific language governing rights and limitations under
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: Daughter of StdExpansion. This class contains routine
33 // which are common to 2D expansion. Typically this inolves physiocal
34 // space operations.
35 //
36 ///////////////////////////////////////////////////////////////////////////////
37 
39 
40 #ifdef max
41 #undef max
42 #endif
43 
44 namespace Nektar
45 {
46  namespace StdRegions
47  {
48 
50  {
51  }
52 
54  const LibUtilities::BasisKey &Ba,
55  const LibUtilities::BasisKey &Bb):
56  StdExpansion(numcoeffs,2, Ba, Bb)
57  {
58  }
59 
61  StdExpansion(T)
62  {
63  }
64 
66  {
67  }
68 
69  //----------------------------
70  // Differentiation Methods
71  //----------------------------
73  Array<OneD, NekDouble> &outarray_d0,
74  Array<OneD, NekDouble> &outarray_d1)
75  {
76  int nquad0 = m_base[0]->GetNumPoints();
77  int nquad1 = m_base[1]->GetNumPoints();
78 
79  if (outarray_d0.num_elements() > 0) // calculate du/dx_0
80  {
81  DNekMatSharedPtr D0 = m_base[0]->GetD();
82  if(inarray.data() == outarray_d0.data())
83  {
84  Array<OneD, NekDouble> wsp(nquad0 * nquad1);
85  Vmath::Vcopy(nquad0 * nquad1,inarray.get(),1,wsp.get(),1);
86  Blas::Dgemm('N', 'N', nquad0, nquad1, nquad0, 1.0,
87  &(D0->GetPtr())[0], nquad0, &wsp[0], nquad0, 0.0,
88  &outarray_d0[0], nquad0);
89  }
90  else
91  {
92  Blas::Dgemm('N', 'N', nquad0, nquad1, nquad0, 1.0,
93  &(D0->GetPtr())[0], nquad0, &inarray[0], nquad0, 0.0,
94  &outarray_d0[0], nquad0);
95  }
96  }
97 
98  if (outarray_d1.num_elements() > 0) // calculate du/dx_1
99  {
100  DNekMatSharedPtr D1 = m_base[1]->GetD();
101  if(inarray.data() == outarray_d1.data())
102  {
103  Array<OneD, NekDouble> wsp(nquad0 * nquad1);
104  Vmath::Vcopy(nquad0 * nquad1,inarray.get(),1,wsp.get(),1);
105  Blas:: Dgemm('N', 'T', nquad0, nquad1, nquad1, 1.0, &wsp[0], nquad0,
106  &(D1->GetPtr())[0], nquad1, 0.0, &outarray_d1[0], nquad0);
107  }
108  else
109  {
110  Blas:: Dgemm('N', 'T', nquad0, nquad1, nquad1, 1.0, &inarray[0], nquad0,
111  &(D1->GetPtr())[0], nquad1, 0.0, &outarray_d1[0], nquad0);
112  }
113  }
114  }
115 
117  {
118  Array<OneD, NekDouble> coll(2);
120 
121  ASSERTL2(coords[0] > -1 - NekConstants::kNekZeroTol, "coord[0] < -1");
122  ASSERTL2(coords[0] < 1 + NekConstants::kNekZeroTol, "coord[0] > 1");
123  ASSERTL2(coords[1] > -1 - NekConstants::kNekZeroTol, "coord[1] < -1");
124  ASSERTL2(coords[1] < 1 + NekConstants::kNekZeroTol, "coord[1] > 1");
125 
126  LocCoordToLocCollapsed(coords,coll);
127 
128  I[0] = m_base[0]->GetI(coll);
129  I[1] = m_base[1]->GetI(coll+1);
130 
131  return v_PhysEvaluate(I,physvals);
132  }
133 
136  const Array<OneD, const NekDouble> &physvals)
137  {
138  NekDouble val;
139  int i;
140  int nq0 = m_base[0]->GetNumPoints();
141  int nq1 = m_base[1]->GetNumPoints();
142  Array<OneD, NekDouble> wsp1(nq1);
143 
144  // interpolate first coordinate direction
145  for (i = 0; i < nq1;++i)
146  {
147  wsp1[i] = Blas::Ddot(nq0, &(I[0]->GetPtr())[0], 1,
148  &physvals[i * nq0], 1);
149  }
150 
151  // interpolate in second coordinate direction
152  val = Blas::Ddot(nq1, I[1]->GetPtr(), 1, wsp1, 1);
153 
154  return val;
155  }
156 
157  //////////////////////////////
158  // Integration Methods
159  //////////////////////////////
160 
164  {
165  int i;
166  NekDouble Int = 0.0;
167  int nquad0 = m_base[0]->GetNumPoints();
168  int nquad1 = m_base[1]->GetNumPoints();
169  Array<OneD, NekDouble> tmp(nquad0 * nquad1);
170 
171  // multiply by integration constants
172  for (i = 0; i < nquad1; ++i)
173  {
174  Vmath::Vmul(nquad0, &inarray[0] + i*nquad0, 1, w0.get(),
175  1, &tmp[0] + i*nquad0, 1);
176  }
177 
178  for (i = 0; i < nquad0; ++i)
179  {
180  Vmath::Vmul(nquad1, &tmp[0]+ i, nquad0, w1.get(), 1,
181  &tmp[0] + i, nquad0);
182  }
183  Int = Vmath::Vsum(nquad0 * nquad1, tmp, 1);
184 
185  return Int;
186  }
187 
189  const Array<OneD, const NekDouble>& base0,
190  const Array<OneD, const NekDouble>& base1,
191  const Array<OneD, const NekDouble>& inarray,
192  Array<OneD, NekDouble> &outarray,
194  bool doCheckCollDir0,
195  bool doCheckCollDir1)
196  {
197  v_BwdTrans_SumFacKernel(base0, base1, inarray, outarray, wsp, doCheckCollDir0, doCheckCollDir1);
198  }
199 
201  const Array<OneD, const NekDouble>& base0,
202  const Array<OneD, const NekDouble>& base1,
203  const Array<OneD, const NekDouble>& inarray,
204  Array<OneD, NekDouble> &outarray,
206  bool doCheckCollDir0,
207  bool doCheckCollDir1)
208  {
209  v_IProductWRTBase_SumFacKernel(base0, base1, inarray, outarray, wsp, doCheckCollDir0, doCheckCollDir1);
210  }
211 
213  const Array<OneD, const NekDouble> &inarray,
214  Array<OneD,NekDouble> &outarray,
215  const StdRegions::StdMatrixKey &mkey)
216  {
217  if (mkey.GetNVarCoeff() == 0
219  {
220  using std::max;
221 
222  // This implementation is only valid when there are no
223  // coefficients associated to the Laplacian operator
224  int nquad0 = m_base[0]->GetNumPoints();
225  int nquad1 = m_base[1]->GetNumPoints();
226  int nqtot = nquad0*nquad1;
227  int nmodes0 = m_base[0]->GetNumModes();
228  int nmodes1 = m_base[1]->GetNumModes();
229  int wspsize = max(max(max(nqtot,m_ncoeffs),nquad1*nmodes0),nquad0*nmodes1);
230 
231  const Array<OneD, const NekDouble>& base0 = m_base[0]->GetBdata();
232  const Array<OneD, const NekDouble>& base1 = m_base[1]->GetBdata();
233 
234  // Allocate temporary storage
235  Array<OneD,NekDouble> wsp0(4*wspsize); // size wspsize
236  Array<OneD,NekDouble> wsp1(wsp0+wspsize); // size 3*wspsize
237 
238  if(!(m_base[0]->Collocation() && m_base[1]->Collocation()))
239  {
240  // LAPLACIAN MATRIX OPERATION
241  // wsp0 = u = B * u_hat
242  // wsp1 = du_dxi1 = D_xi1 * wsp0 = D_xi1 * u
243  // wsp2 = du_dxi2 = D_xi2 * wsp0 = D_xi2 * u
244  BwdTrans_SumFacKernel(base0,base1,inarray,wsp0,wsp1,true,true);
245  LaplacianMatrixOp_MatFree_Kernel(wsp0, outarray, wsp1);
246  }
247  else
248  {
249  LaplacianMatrixOp_MatFree_Kernel(inarray, outarray, wsp1);
250  }
251  }
252  else
253  {
255  inarray,outarray,mkey);
256  }
257  }
258 
259 
260 
262  const Array<OneD, const NekDouble> &inarray,
263  Array<OneD,NekDouble> &outarray,
264  const StdRegions::StdMatrixKey &mkey)
265  {
266  if (mkey.GetNVarCoeff() == 0
268  {
269  using std::max;
270 
271  int nquad0 = m_base[0]->GetNumPoints();
272  int nquad1 = m_base[1]->GetNumPoints();
273  int nqtot = nquad0*nquad1;
274  int nmodes0 = m_base[0]->GetNumModes();
275  int nmodes1 = m_base[1]->GetNumModes();
276  int wspsize = max(max(max(nqtot,m_ncoeffs),nquad1*nmodes0),
277  nquad0*nmodes1);
278  NekDouble lambda =
280 
281  const Array<OneD, const NekDouble>& base0 = m_base[0]->GetBdata();
282  const Array<OneD, const NekDouble>& base1 = m_base[1]->GetBdata();
283 
284  // Allocate temporary storage
285  Array<OneD,NekDouble> wsp0(5*wspsize); // size wspsize
286  Array<OneD,NekDouble> wsp1(wsp0 + wspsize); // size wspsize
287  Array<OneD,NekDouble> wsp2(wsp0 + 2*wspsize);// size 3*wspsize
288 
289  if (!(m_base[0]->Collocation() && m_base[1]->Collocation()))
290  {
291  // MASS MATRIX OPERATION
292  // The following is being calculated:
293  // wsp0 = B * u_hat = u
294  // wsp1 = W * wsp0
295  // outarray = B^T * wsp1 = B^T * W * B * u_hat = M * u_hat
296  BwdTrans_SumFacKernel (base0, base1, inarray,
297  wsp0, wsp2,true,true);
298  MultiplyByQuadratureMetric (wsp0, wsp1);
299  IProductWRTBase_SumFacKernel(base0, base1, wsp1, outarray,
300  wsp2, true, true);
301 
302  LaplacianMatrixOp_MatFree_Kernel(wsp0, wsp1, wsp2);
303  }
304  else
305  {
306  MultiplyByQuadratureMetric(inarray,outarray);
307  LaplacianMatrixOp_MatFree_Kernel(inarray, wsp1, wsp2);
308  }
309 
310  // outarray = lambda * outarray + wsp1
311  // = (lambda * M + L ) * u_hat
312  Vmath::Svtvp(m_ncoeffs, lambda, &outarray[0], 1,
313  &wsp1[0], 1, &outarray[0], 1);
314  }
315  else
316  {
318  inarray,outarray,mkey);
319  }
320  }
321 
322  } //end namespace
323 } //end namespace
void HelmholtzMatrixOp_MatFree_GenericImpl(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdMatrixKey &mkey)
NekDouble GetConstFactor(const ConstFactorType &factor) const
Definition: StdMatrixKey.h:122
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:942
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:471
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)
boost::shared_ptr< DNekMat > DNekMatSharedPtr
Definition: NekTypeDefs.hpp:70
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...
bool ConstFactorExists(const ConstFactorType &factor) const
Definition: StdMatrixKey.h:131
The base class for all shapes.
Definition: StdExpansion.h:69
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)
T Ddot(int n, const Array< OneD, const T > &w, const int incw, const Array< OneD, const T > &x, const int incx, const Array< OneD, const int > &y, const int incy)
Definition: VmathArray.hpp:434
#define ASSERTL2(condition, msg)
Assert Level 2 – Debugging which is used FULLDEBUG compilation mode. This level assert is designed t...
Definition: ErrorUtil.hpp:213
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:723
Array< OneD, LibUtilities::BasisSharedPtr > m_base
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1047
Describes the specification for a Basis.
Definition: Basis.h:50
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:169
virtual void v_HelmholtzMatrixOp_MatFree(const Array< OneD, const NekDouble > &inarray, Array< OneD, NekDouble > &outarray, const StdRegions::StdMatrixKey &mkey)