Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
PyrExp.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File PyrExp.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: PyrExp routines
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
36 #include <LocalRegions/PyrExp.h>
38 
39 namespace Nektar
40 {
41  namespace LocalRegions
42  {
43 
45  const LibUtilities::BasisKey &Bb,
46  const LibUtilities::BasisKey &Bc,
48  StdExpansion (LibUtilities::StdPyrData::getNumberOfCoefficients(
49  Ba.GetNumModes(),
50  Bb.GetNumModes(),
51  Bc.GetNumModes()),
52  3, Ba, Bb, Bc),
53  StdExpansion3D(LibUtilities::StdPyrData::getNumberOfCoefficients(
54  Ba.GetNumModes(),
55  Bb.GetNumModes(),
56  Bc.GetNumModes()),
57  Ba, Bb, Bc),
58  StdPyrExp (Ba,Bb,Bc),
59  Expansion (geom),
60  Expansion3D (geom),
61  m_matrixManager(
62  boost::bind(&PyrExp::CreateMatrix, this, _1),
63  std::string("PyrExpMatrix")),
64  m_staticCondMatrixManager(
65  boost::bind(&PyrExp::CreateStaticCondMatrix, this, _1),
66  std::string("PyrExpStaticCondMatrix"))
67  {
68  }
69 
71  StdExpansion (T),
72  StdExpansion3D(T),
73  StdPyrExp (T),
74  Expansion (T),
75  Expansion3D (T),
76  m_matrixManager(T.m_matrixManager),
77  m_staticCondMatrixManager(T.m_staticCondMatrixManager)
78  {
79  }
80 
82  {
83  }
84 
85 
86  //----------------------------
87  // Integration Methods
88  //----------------------------
89 
90  /**
91  * \brief Integrate the physical point list \a inarray over pyramidic
92  * region and return the value.
93  *
94  * Inputs:\n
95  *
96  * - \a inarray: definition of function to be returned at quadrature
97  * point of expansion.
98  *
99  * Outputs:\n
100  *
101  * - returns \f$\int^1_{-1}\int^1_{-1}\int^1_{-1} u(\bar \eta_1,
102  * \eta_2, \eta_3) J[i,j,k] d \bar \eta_1 d \eta_2 d \eta_3\f$ \n \f$=
103  * \sum_{i=0}^{Q_1 - 1} \sum_{j=0}^{Q_2 - 1} \sum_{k=0}^{Q_3 - 1}
104  * u(\bar \eta_{1i}^{0,0}, \eta_{2j}^{0,0},\eta_{3k}^{2,0})w_{i}^{0,0}
105  * w_{j}^{0,0} \hat w_{k}^{2,0} \f$ \n where \f$inarray[i,j, k] =
106  * u(\bar \eta_{1i},\eta_{2j}, \eta_{3k}) \f$, \n \f$\hat w_{k}^{2,0}
107  * = \frac {w^{2,0}} {2} \f$ \n and \f$ J[i,j,k] \f$ is the Jacobian
108  * evaluated at the quadrature point.
109  */
110  NekDouble PyrExp::v_Integral(const Array<OneD, const NekDouble> &inarray)
111  {
112  int nquad0 = m_base[0]->GetNumPoints();
113  int nquad1 = m_base[1]->GetNumPoints();
114  int nquad2 = m_base[2]->GetNumPoints();
115  Array<OneD, const NekDouble> jac = m_metricinfo->GetJac(GetPointsKeys());
116  Array<OneD, NekDouble> tmp(nquad0*nquad1*nquad2);
117 
118  // multiply inarray with Jacobian
119  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
120  {
121  Vmath::Vmul(nquad0*nquad1*nquad2,&jac[0],1,(NekDouble*)&inarray[0],1, &tmp[0],1);
122  }
123  else
124  {
125  Vmath::Smul(nquad0*nquad1*nquad2,(NekDouble) jac[0], (NekDouble*)&inarray[0],1,&tmp[0],1);
126  }
127 
128  // call StdPyrExp version;
129  return StdPyrExp::v_Integral(tmp);
130  }
131 
132 
133  //----------------------------
134  // Differentiation Methods
135  //----------------------------
136 
137  void PyrExp::v_PhysDeriv(const Array<OneD, const NekDouble>& inarray,
138  Array<OneD, NekDouble>& out_d0,
139  Array<OneD, NekDouble>& out_d1,
140  Array<OneD, NekDouble>& out_d2)
141  {
142  int nquad0 = m_base[0]->GetNumPoints();
143  int nquad1 = m_base[1]->GetNumPoints();
144  int nquad2 = m_base[2]->GetNumPoints();
145  Array<TwoD, const NekDouble> gmat =
146  m_metricinfo->GetDerivFactors(GetPointsKeys());
147  Array<OneD,NekDouble> diff0(nquad0*nquad1*nquad2);
148  Array<OneD,NekDouble> diff1(nquad0*nquad1*nquad2);
149  Array<OneD,NekDouble> diff2(nquad0*nquad1*nquad2);
150 
151  StdPyrExp::v_PhysDeriv(inarray, diff0, diff1, diff2);
152 
153  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
154  {
155  if(out_d0.num_elements())
156  {
157  Vmath::Vmul (nquad0*nquad1*nquad2,&gmat[0][0],1,&diff0[0],1, &out_d0[0], 1);
158  Vmath::Vvtvp (nquad0*nquad1*nquad2,&gmat[1][0],1,&diff1[0],1, &out_d0[0], 1,&out_d0[0],1);
159  Vmath::Vvtvp (nquad0*nquad1*nquad2,&gmat[2][0],1,&diff2[0],1, &out_d0[0], 1,&out_d0[0],1);
160  }
161 
162  if(out_d1.num_elements())
163  {
164  Vmath::Vmul (nquad0*nquad1*nquad2,&gmat[3][0],1,&diff0[0],1, &out_d1[0], 1);
165  Vmath::Vvtvp (nquad0*nquad1*nquad2,&gmat[4][0],1,&diff1[0],1, &out_d1[0], 1,&out_d1[0],1);
166  Vmath::Vvtvp (nquad0*nquad1*nquad2,&gmat[5][0],1,&diff2[0],1, &out_d1[0], 1,&out_d1[0],1);
167  }
168 
169  if(out_d2.num_elements())
170  {
171  Vmath::Vmul (nquad0*nquad1*nquad2,&gmat[6][0],1,&diff0[0],1, &out_d2[0], 1);
172  Vmath::Vvtvp (nquad0*nquad1*nquad2,&gmat[7][0],1,&diff1[0],1, &out_d2[0], 1, &out_d2[0],1);
173  Vmath::Vvtvp (nquad0*nquad1*nquad2,&gmat[8][0],1,&diff2[0],1, &out_d2[0], 1, &out_d2[0],1);
174  }
175  }
176  else // regular geometry
177  {
178  if(out_d0.num_elements())
179  {
180  Vmath::Smul (nquad0*nquad1*nquad2,gmat[0][0],&diff0[0],1, &out_d0[0], 1);
181  Blas::Daxpy (nquad0*nquad1*nquad2,gmat[1][0],&diff1[0],1, &out_d0[0], 1);
182  Blas::Daxpy (nquad0*nquad1*nquad2,gmat[2][0],&diff2[0],1, &out_d0[0], 1);
183  }
184 
185  if(out_d1.num_elements())
186  {
187  Vmath::Smul (nquad0*nquad1*nquad2,gmat[3][0],&diff0[0],1, &out_d1[0], 1);
188  Blas::Daxpy (nquad0*nquad1*nquad2,gmat[4][0],&diff1[0],1, &out_d1[0], 1);
189  Blas::Daxpy (nquad0*nquad1*nquad2,gmat[5][0],&diff2[0],1, &out_d1[0], 1);
190  }
191 
192  if(out_d2.num_elements())
193  {
194  Vmath::Smul (nquad0*nquad1*nquad2,gmat[6][0],&diff0[0],1, &out_d2[0], 1);
195  Blas::Daxpy (nquad0*nquad1*nquad2,gmat[7][0],&diff1[0],1, &out_d2[0], 1);
196  Blas::Daxpy (nquad0*nquad1*nquad2,gmat[8][0],&diff2[0],1, &out_d2[0], 1);
197  }
198  }
199  }
200 
201 
202  //---------------------------------------
203  // Transforms
204  //---------------------------------------
205 
206  /**
207  * \brief Forward transform from physical quadrature space stored in
208  * \a inarray and evaluate the expansion coefficients and store in \a
209  * (this)->m_coeffs
210  *
211  * Inputs:\n
212  *
213  * - \a inarray: array of physical quadrature points to be transformed
214  *
215  * Outputs:\n
216  *
217  * - (this)->_coeffs: updated array of expansion coefficients.
218  */
219  void PyrExp::v_FwdTrans(const Array<OneD, const NekDouble>& inarray,
220  Array<OneD, NekDouble>& outarray)
221  {
222  if(m_base[0]->Collocation() &&
223  m_base[1]->Collocation() &&
224  m_base[2]->Collocation())
225  {
226  Vmath::Vcopy(GetNcoeffs(),&inarray[0],1,&outarray[0],1);
227  }
228  else
229  {
230  v_IProductWRTBase(inarray,outarray);
231 
232  // get Mass matrix inverse
234  DetShapeType(),*this);
235  DNekScalMatSharedPtr matsys = m_matrixManager[masskey];
236 
237  // copy inarray in case inarray == outarray
238  DNekVec in (m_ncoeffs,outarray);
239  DNekVec out(m_ncoeffs,outarray,eWrapper);
240 
241  out = (*matsys)*in;
242  }
243  }
244 
245 
246  //---------------------------------------
247  // Inner product functions
248  //---------------------------------------
249 
250  /**
251  * \brief Calculate the inner product of inarray with respect to the
252  * basis B=base0*base1*base2 and put into outarray:
253  *
254  * \f$ \begin{array}{rcl} I_{pqr} = (\phi_{pqr}, u)_{\delta} & = &
255  * \sum_{i=0}^{nq_0} \sum_{j=0}^{nq_1} \sum_{k=0}^{nq_2} \psi_{p}^{a}
256  * (\bar \eta_{1i}) \psi_{q}^{a} (\eta_{2j}) \psi_{pqr}^{c}
257  * (\eta_{3k}) w_i w_j w_k u(\bar \eta_{1,i} \eta_{2,j} \eta_{3,k})
258  * J_{i,j,k}\\ & = & \sum_{i=0}^{nq_0} \psi_p^a(\bar \eta_{1,i})
259  * \sum_{j=0}^{nq_1} \psi_{q}^a(\eta_{2,j}) \sum_{k=0}^{nq_2}
260  * \psi_{pqr}^c u(\bar \eta_{1i},\eta_{2j},\eta_{3k}) J_{i,j,k}
261  * \end{array} \f$ \n
262  *
263  * where
264  *
265  * \f$\phi_{pqr} (\xi_1 , \xi_2 , \xi_3) = \psi_p^a (\bar \eta_1)
266  * \psi_{q}^a (\eta_2) \psi_{pqr}^c (\eta_3) \f$ \n
267  *
268  * which can be implemented as \n \f$f_{pqr} (\xi_{3k}) =
269  * \sum_{k=0}^{nq_3} \psi_{pqr}^c u(\bar
270  * \eta_{1i},\eta_{2j},\eta_{3k}) J_{i,j,k} = {\bf B_3 U} \f$ \n \f$
271  * g_{pq} (\xi_{3k}) = \sum_{j=0}^{nq_1} \psi_{q}^a (\xi_{2j}) f_{pqr}
272  * (\xi_{3k}) = {\bf B_2 F} \f$ \n \f$ (\phi_{pqr}, u)_{\delta} =
273  * \sum_{k=0}^{nq_0} \psi_{p}^a (\xi_{3k}) g_{pq} (\xi_{3k}) = {\bf
274  * B_1 G} \f$
275  */
277  const Array<OneD, const NekDouble> &inarray,
278  Array<OneD, NekDouble> &outarray)
279  {
280  int nquad0 = m_base[0]->GetNumPoints();
281  int nquad1 = m_base[1]->GetNumPoints();
282  int nquad2 = m_base[2]->GetNumPoints();
283  Array<OneD, const NekDouble> jac = m_metricinfo->GetJac(GetPointsKeys());
284  Array<OneD, NekDouble> tmp(nquad0*nquad1*nquad2);
285 
286  // multiply inarray with Jacobian
287  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
288  {
289  Vmath::Vmul(nquad0*nquad1*nquad2,&jac[0],1,(NekDouble*)&inarray[0],1,&tmp[0],1);
290  }
291  else
292  {
293  Vmath::Smul(nquad0*nquad1*nquad2,jac[0],(NekDouble*)&inarray[0],1,&tmp[0],1);
294  }
295 
296  StdPyrExp::v_IProductWRTBase(tmp,outarray);
297  }
298 
299 
300  //---------------------------------------
301  // Evaluation functions
302  //---------------------------------------
303 
304  /*
305  * @brief Get the coordinates #coords at the local coordinates
306  * #Lcoords
307  */
308  void PyrExp::v_GetCoord(const Array<OneD, const NekDouble>& Lcoords,
309  Array<OneD, NekDouble>& coords)
310  {
311  int i;
312 
313  ASSERTL1(Lcoords[0] <= -1.0 && Lcoords[0] >= 1.0 &&
314  Lcoords[1] <= -1.0 && Lcoords[1] >= 1.0 &&
315  Lcoords[2] <= -1.0 && Lcoords[2] >= 1.0,
316  "Local coordinates are not in region [-1,1]");
317 
318  // m_geom->FillGeom(); // TODO: implement FillGeom()
319 
320  for(i = 0; i < m_geom->GetCoordim(); ++i)
321  {
322  coords[i] = m_geom->GetCoord(i,Lcoords);
323  }
324  }
325 
327  Array<OneD, NekDouble> &coords_1,
328  Array<OneD, NekDouble> &coords_2,
329  Array<OneD, NekDouble> &coords_3)
330  {
331  Expansion::v_GetCoords(coords_1, coords_2, coords_3);
332  }
333 
334  NekDouble PyrExp::v_PhysEvaluate(const Array<OneD, const NekDouble>& coord,
335  const Array<OneD, const NekDouble>& physvals)
336  {
337  Array<OneD,NekDouble> Lcoord(3);
338 
339  ASSERTL0(m_geom,"m_geom not defined");
340 
341  //TODO: check GetLocCoords()
342  m_geom->GetLocCoords(coord, Lcoord);
343 
344  return StdPyrExp::v_PhysEvaluate(Lcoord, physvals);
345  }
346 
347 
348  //---------------------------------------
349  // Helper functions
350  //---------------------------------------
351 
353  {
354  return m_geom->GetCoordim();
355  }
356 
358  const int face,
359  const StdRegions::StdExpansionSharedPtr &FaceExp,
360  const Array<OneD, const NekDouble> &inarray,
361  Array<OneD, NekDouble> &outarray,
363  {
364  int nq0 = m_base[0]->GetNumPoints();
365  int nq1 = m_base[1]->GetNumPoints();
366  int nq2 = m_base[2]->GetNumPoints();
367 
368  Array<OneD,NekDouble> o_tmp(GetFaceNumPoints(face));
369 
370  if (orient == StdRegions::eNoOrientation)
371  {
372  orient = GetForient(face);
373  }
374 
375  switch(face)
376  {
377  case 0:
379  {
380  //Directions A and B positive
381  Vmath::Vcopy(nq0*nq1,&(inarray[0]),1,&(o_tmp[0]),1);
382  }
383  else if(orient == StdRegions::eDir1BwdDir1_Dir2FwdDir2)
384  {
385  //Direction A negative and B positive
386  for (int j=0; j<nq1; j++)
387  {
388  Vmath::Vcopy(nq0,&(inarray[0])+(nq0-1)+j*nq0,-1,&(o_tmp[0])+(j*nq0),1);
389  }
390  }
391  else if(orient == StdRegions::eDir1FwdDir1_Dir2BwdDir2)
392  {
393  //Direction A positive and B negative
394  for (int j=0; j<nq1; j++)
395  {
396  Vmath::Vcopy(nq0,&(inarray[0])+nq0*(nq1-1-j),1,&(o_tmp[0])+(j*nq0),1);
397  }
398  }
399  else if(orient == StdRegions::eDir1BwdDir1_Dir2BwdDir2)
400  {
401  //Direction A negative and B negative
402  for(int j=0; j<nq1; j++)
403  {
404  Vmath::Vcopy(nq0,&(inarray[0])+(nq0*nq1-1-j*nq0),-1,&(o_tmp[0])+(j*nq0),1);
405  }
406  }
407  else if(orient == StdRegions::eDir1FwdDir2_Dir2FwdDir1)
408  {
409  //Transposed, Direction A and B positive
410  for (int i=0; i<nq0; i++)
411  {
412  Vmath::Vcopy(nq1,&(inarray[0])+i,nq0,&(o_tmp[0])+(i*nq1),1);
413  }
414  }
415  else if(orient == StdRegions::eDir1FwdDir2_Dir2BwdDir1)
416  {
417  //Transposed, Direction A positive and B negative
418  for (int i=0; i<nq0; i++)
419  {
420  Vmath::Vcopy(nq1,&(inarray[0])+(nq0-1-i),nq0,&(o_tmp[0])+(i*nq1),1);
421  }
422  }
423  else if(orient == StdRegions::eDir1BwdDir2_Dir2FwdDir1)
424  {
425  //Transposed, Direction A negative and B positive
426  for (int i=0; i<nq0; i++)
427  {
428  Vmath::Vcopy(nq1,&(inarray[0])+i+nq0*(nq1-1),-nq0,&(o_tmp[0])+(i*nq1),1);
429  }
430  }
431  else if(orient == StdRegions::eDir1BwdDir2_Dir2BwdDir1)
432  {
433  //Transposed, Direction A and B negative
434  for (int i=0; i<nq0; i++)
435  {
436  Vmath::Vcopy(nq1,&(inarray[0])+(nq0*nq1-1-i),-nq0,&(o_tmp[0])+(i*nq1),1);
437  }
438  }
439  LibUtilities::Interp2D(m_base[0]->GetPointsKey(), m_base[1]->GetPointsKey(), o_tmp,
440  FaceExp->GetBasis(0)->GetPointsKey(),FaceExp->GetBasis(1)->GetPointsKey(),outarray);
441  break;
442 
443  case 1:
444  {
445  for (int k = 0; k < nq2; k++)
446  {
447  Vmath::Vcopy(nq0,inarray.get()+nq0*nq1*k,1,outarray.get()+k*nq0,1);
448  }
449  LibUtilities::Interp2D(m_base[0]->GetPointsKey(), m_base[2]->GetPointsKey(), outarray.get(),
450  FaceExp->GetBasis(0)->GetPointsKey(),FaceExp->GetBasis(1)->GetPointsKey(),o_tmp.get());
451  break;
452  }
453 
454  case 2:
455  {
456  Vmath::Vcopy(nq1*nq2,inarray.get()+(nq0-1),nq0,outarray.get(),1);
457  LibUtilities::Interp2D(m_base[1]->GetPointsKey(), m_base[2]->GetPointsKey(), outarray.get(),
458  FaceExp->GetBasis(0)->GetPointsKey(),FaceExp->GetBasis(1)->GetPointsKey(),o_tmp.get());
459  break;
460  }
461 
462  case 3:
463  {
464  for (int k = 0; k < nq2; k++)
465  {
466  Vmath::Vcopy(nq0,inarray.get()+nq0*(nq1-1)+nq0*nq1*k,1,outarray.get()+(k*nq0),1);
467  }
468  LibUtilities::Interp2D(m_base[0]->GetPointsKey(), m_base[2]->GetPointsKey(), outarray.get(),
469  FaceExp->GetBasis(0)->GetPointsKey(),FaceExp->GetBasis(1)->GetPointsKey(),o_tmp.get());
470  }
471 
472  case 4:
473  {
474  Vmath::Vcopy(nq1*nq2,inarray.get(),nq0,outarray.get(),1);
475  LibUtilities::Interp2D(m_base[1]->GetPointsKey(), m_base[2]->GetPointsKey(), outarray.get(),
476  FaceExp->GetBasis(0)->GetPointsKey(),FaceExp->GetBasis(1)->GetPointsKey(),o_tmp.get());
477  break;
478  }
479 
480  default:
481  ASSERTL0(false,"face value (> 4) is out of range");
482  break;
483  }
484 
485  if (face > 0)
486  {
487  int fnq1 = FaceExp->GetNumPoints(0);
488  int fnq2 = FaceExp->GetNumPoints(1);
489 
490  if ((int)orient == 7)
491  {
492  for (int j = 0; j < fnq2; ++j)
493  {
494  Vmath::Vcopy(fnq1, o_tmp.get()+((j+1)*fnq1-1), -1, outarray.get()+j*fnq1, 1);
495  }
496  }
497  else
498  {
499  Vmath::Vcopy(fnq1*fnq2, o_tmp.get(), 1, outarray.get(), 1);
500  }
501  }
502  }
503 
504  void PyrExp::v_ComputeFaceNormal(const int face)
505  {
506  const SpatialDomains::GeomFactorsSharedPtr &geomFactors =
507  GetGeom()->GetMetricInfo();
509  SpatialDomains::GeomType type = geomFactors->GetGtype();
510  const Array<TwoD, const NekDouble> &df = geomFactors->GetDerivFactors(ptsKeys);
511  const Array<OneD, const NekDouble> &jac = geomFactors->GetJac(ptsKeys);
512 
513  LibUtilities::BasisKey tobasis0 = DetFaceBasisKey(face,0);
514  LibUtilities::BasisKey tobasis1 = DetFaceBasisKey(face,1);
515 
516  // Number of quadrature points in face expansion.
517  int nq_face = tobasis0.GetNumPoints()*tobasis1.GetNumPoints();
518 
519  int vCoordDim = GetCoordim();
520  int i;
521 
522  m_faceNormals[face] = Array<OneD, Array<OneD, NekDouble> >(vCoordDim);
523  Array<OneD, Array<OneD, NekDouble> > &normal = m_faceNormals[face];
524  for (i = 0; i < vCoordDim; ++i)
525  {
526  normal[i] = Array<OneD, NekDouble>(nq_face);
527  }
528 
529  // Regular geometry case
530  if (type == SpatialDomains::eRegular ||
532  {
533  NekDouble fac;
534  // Set up normals
535  switch(face)
536  {
537  case 0:
538  {
539  for(i = 0; i < vCoordDim; ++i)
540  {
541  normal[i][0] = -df[3*i+2][0];
542  }
543  break;
544  }
545  case 1:
546  {
547  for(i = 0; i < vCoordDim; ++i)
548  {
549  normal[i][0] = -df[3*i+1][0];
550  }
551  break;
552  }
553  case 2:
554  {
555  for(i = 0; i < vCoordDim; ++i)
556  {
557  normal[i][0] = df[3*i][0]+df[3*i+2][0];
558  }
559  break;
560  }
561  case 3:
562  {
563  for(i = 0; i < vCoordDim; ++i)
564  {
565  normal[i][0] = df[3*i+1][0]+df[3*i+2][0];
566  }
567  break;
568  }
569  case 4:
570  {
571  for(i = 0; i < vCoordDim; ++i)
572  {
573  normal[i][0] = -df[3*i][0];
574  }
575  break;
576  }
577  default:
578  ASSERTL0(false,"face is out of range (face < 4)");
579  }
580 
581  // Normalise resulting vector.
582  fac = 0.0;
583  for(i = 0; i < vCoordDim; ++i)
584  {
585  fac += normal[i][0]*normal[i][0];
586  }
587  fac = 1.0/sqrt(fac);
588  for (i = 0; i < vCoordDim; ++i)
589  {
590  Vmath::Fill(nq_face,fac*normal[i][0],normal[i],1);
591  }
592  }
593  else
594  {
595  // Set up deformed normals.
596  int j, k;
597 
598  int nq0 = ptsKeys[0].GetNumPoints();
599  int nq1 = ptsKeys[1].GetNumPoints();
600  int nq2 = ptsKeys[2].GetNumPoints();
601  int nq01 = nq0*nq1;
602  int nqtot;
603 
604  // Determine number of quadrature points on the face.
605  if (face == 0)
606  {
607  nqtot = nq0*nq1;
608  }
609  else if (face == 1 || face == 3)
610  {
611  nqtot = nq0*nq2;
612  }
613  else
614  {
615  nqtot = nq1*nq2;
616  }
617 
618  LibUtilities::PointsKey points0;
619  LibUtilities::PointsKey points1;
620 
621  Array<OneD, NekDouble> normals(vCoordDim*nqtot,0.0);
622 
623  // Extract Jacobian along face and recover local derivatives
624  // (dx/dr) for polynomial interpolation by multiplying m_gmat by
625  // jacobian
626  switch(face)
627  {
628  case 0:
629  {
630  for(j = 0; j < nq01; ++j)
631  {
632  normals[j] = -df[2][j]*jac[j];
633  normals[nqtot+j] = -df[5][j]*jac[j];
634  normals[2*nqtot+j] = -df[8][j]*jac[j];
635  }
636 
637  points0 = ptsKeys[0];
638  points1 = ptsKeys[1];
639  break;
640  }
641 
642  case 1:
643  {
644  for (j = 0; j < nq0; ++j)
645  {
646  for(k = 0; k < nq2; ++k)
647  {
648  int tmp = j+nq01*k;
649  normals[j+k*nq0] =
650  -df[1][tmp]*jac[tmp];
651  normals[nqtot+j+k*nq0] =
652  -df[4][tmp]*jac[tmp];
653  normals[2*nqtot+j+k*nq0] =
654  -df[7][tmp]*jac[tmp];
655  }
656  }
657 
658  points0 = ptsKeys[0];
659  points1 = ptsKeys[2];
660  break;
661  }
662 
663  case 2:
664  {
665  for (j = 0; j < nq1; ++j)
666  {
667  for(k = 0; k < nq2; ++k)
668  {
669  int tmp = nq0-1+nq0*j+nq01*k;
670  normals[j+k*nq1] =
671  (df[0][tmp]+df[2][tmp])*jac[tmp];
672  normals[nqtot+j+k*nq1] =
673  (df[3][tmp]+df[5][tmp])*jac[tmp];
674  normals[2*nqtot+j+k*nq1] =
675  (df[6][tmp]+df[8][tmp])*jac[tmp];
676  }
677  }
678 
679  points0 = ptsKeys[1];
680  points1 = ptsKeys[2];
681  break;
682  }
683 
684  case 3:
685  {
686  for (j = 0; j < nq0; ++j)
687  {
688  for(k = 0; k < nq2; ++k)
689  {
690  int tmp = nq0*(nq1-1) + j + nq01*k;
691  normals[j+k*nq0] =
692  (df[1][tmp]+df[2][tmp])*jac[tmp];
693  normals[nqtot+j+k*nq0] =
694  (df[4][tmp]+df[5][tmp])*jac[tmp];
695  normals[2*nqtot+j+k*nq0] =
696  (df[7][tmp]+df[8][tmp])*jac[tmp];
697  }
698  }
699 
700  points0 = ptsKeys[0];
701  points1 = ptsKeys[2];
702  break;
703  }
704 
705  case 4:
706  {
707  for (j = 0; j < nq1; ++j)
708  {
709  for(k = 0; k < nq2; ++k)
710  {
711  int tmp = j*nq0+nq01*k;
712  normals[j+k*nq1] =
713  -df[0][tmp]*jac[tmp];
714  normals[nqtot+j+k*nq1] =
715  -df[3][tmp]*jac[tmp];
716  normals[2*nqtot+j+k*nq1] =
717  -df[6][tmp]*jac[tmp];
718  }
719  }
720 
721  points0 = ptsKeys[1];
722  points1 = ptsKeys[2];
723  break;
724  }
725 
726  default:
727  ASSERTL0(false,"face is out of range (face < 4)");
728  }
729 
730  Array<OneD, NekDouble> work (nq_face, 0.0);
731  // Interpolate Jacobian and invert
732  LibUtilities::Interp2D(points0, points1, jac,
733  tobasis0.GetPointsKey(),
734  tobasis1.GetPointsKey(),
735  work);
736  Vmath::Sdiv(nq_face, 1.0, &work[0], 1, &work[0], 1);
737 
738  // Interpolate normal and multiply by inverse Jacobian.
739  for(i = 0; i < vCoordDim; ++i)
740  {
741  LibUtilities::Interp2D(points0, points1,
742  &normals[i*nqtot],
743  tobasis0.GetPointsKey(),
744  tobasis1.GetPointsKey(),
745  &normal[i][0]);
746  Vmath::Vmul(nq_face,work,1,normal[i],1,normal[i],1);
747  }
748 
749  // Normalise to obtain unit normals.
750  Vmath::Zero(nq_face,work,1);
751  for(i = 0; i < GetCoordim(); ++i)
752  {
753  Vmath::Vvtvp(nq_face,normal[i],1,normal[i],1,work,1,work,1);
754  }
755 
756  Vmath::Vsqrt(nq_face,work,1,work,1);
757  Vmath::Sdiv (nq_face,1.0,work,1,work,1);
758 
759  for(i = 0; i < GetCoordim(); ++i)
760  {
761  Vmath::Vmul(nq_face,normal[i],1,work,1,normal[i],1);
762  }
763  }
764  }
765 
766  //---------------------------------------
767  // Matrix creation functions
768  //---------------------------------------
769 
771  {
772  DNekMatSharedPtr returnval;
773 
774  switch(mkey.GetMatrixType())
775  {
782  returnval = Expansion3D::v_GenMatrix(mkey);
783  break;
784  default:
785  returnval = StdPyrExp::v_GenMatrix(mkey);
786  }
787 
788  return returnval;
789  }
790 
792  {
793  LibUtilities::BasisKey bkey0 = m_base[0]->GetBasisKey();
794  LibUtilities::BasisKey bkey1 = m_base[1]->GetBasisKey();
795  LibUtilities::BasisKey bkey2 = m_base[2]->GetBasisKey();
797  MemoryManager<StdPyrExp>::AllocateSharedPtr(bkey0, bkey1, bkey2);
798 
799  return tmp->GetStdMatrix(mkey);
800  }
801 
803  {
804  return m_matrixManager[mkey];
805  }
806 
808  {
809  return m_staticCondMatrixManager[mkey];
810  }
811 
813  {
815  }
816 
818  {
819  DNekScalMatSharedPtr returnval;
821 
822  ASSERTL2(m_metricinfo->GetGtype() != SpatialDomains::eNoGeomType,"Geometric information is not set up");
823 
824  switch(mkey.GetMatrixType())
825  {
826  case StdRegions::eMass:
827  {
828  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
829  {
830  NekDouble one = 1.0;
831  DNekMatSharedPtr mat = GenMatrix(mkey);
833  }
834  else
835  {
836  NekDouble jac = (m_metricinfo->GetJac(ptsKeys))[0];
837  DNekMatSharedPtr mat = GetStdMatrix(mkey);
839  }
840  }
841  break;
843  {
844  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
845  {
846  NekDouble one = 1.0;
848  *this);
849  DNekMatSharedPtr mat = GenMatrix(masskey);
850  mat->Invert();
852  }
853  else
854  {
855  NekDouble fac = 1.0/(m_metricinfo->GetJac(ptsKeys))[0];
856  DNekMatSharedPtr mat = GetStdMatrix(mkey);
858  }
859  }
860  break;
862  {
863  if (m_metricinfo->GetGtype() == SpatialDomains::eDeformed ||
864  mkey.GetNVarCoeff() > 0 ||
865  mkey.ConstFactorExists(
867  {
868  NekDouble one = 1.0;
869  DNekMatSharedPtr mat = GenMatrix(mkey);
870 
872  }
873  else
874  {
876  mkey.GetShapeType(), *this);
878  mkey.GetShapeType(), *this);
880  mkey.GetShapeType(), *this);
882  mkey.GetShapeType(), *this);
884  mkey.GetShapeType(), *this);
886  mkey.GetShapeType(), *this);
887 
888  DNekMat &lap00 = *GetStdMatrix(lap00key);
889  DNekMat &lap01 = *GetStdMatrix(lap01key);
890  DNekMat &lap02 = *GetStdMatrix(lap02key);
891  DNekMat &lap11 = *GetStdMatrix(lap11key);
892  DNekMat &lap12 = *GetStdMatrix(lap12key);
893  DNekMat &lap22 = *GetStdMatrix(lap22key);
894 
895  NekDouble jac = (m_metricinfo->GetJac(ptsKeys))[0];
896  Array<TwoD, const NekDouble> gmat =
897  m_metricinfo->GetGmat(ptsKeys);
898 
899  int rows = lap00.GetRows();
900  int cols = lap00.GetColumns();
901 
903  ::AllocateSharedPtr(rows,cols);
904 
905  (*lap) = gmat[0][0]*lap00
906  + gmat[4][0]*lap11
907  + gmat[8][0]*lap22
908  + gmat[3][0]*(lap01 + Transpose(lap01))
909  + gmat[6][0]*(lap02 + Transpose(lap02))
910  + gmat[7][0]*(lap12 + Transpose(lap12));
911 
912  returnval = MemoryManager<DNekScalMat>
913  ::AllocateSharedPtr(jac, lap);
914  }
915  }
916  break;
918  {
920  MatrixKey masskey(StdRegions::eMass, mkey.GetShapeType(), *this);
921  DNekScalMat &MassMat = *(this->m_matrixManager[masskey]);
922  MatrixKey lapkey(StdRegions::eLaplacian, mkey.GetShapeType(), *this, mkey.GetConstFactors(), mkey.GetVarCoeffs());
923  DNekScalMat &LapMat = *(this->m_matrixManager[lapkey]);
924 
925  int rows = LapMat.GetRows();
926  int cols = LapMat.GetColumns();
927 
929 
930  (*helm) = LapMat + factor*MassMat;
931 
932  returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(1.0, helm);
933  }
934  break;
935  default:
936  NEKERROR(ErrorUtil::efatal, "Matrix creation not defined");
937  break;
938  }
939 
940  return returnval;
941  }
942 
944  {
945  DNekScalBlkMatSharedPtr returnval;
946 
947  ASSERTL2(m_metricinfo->GetGtype() != SpatialDomains::eNoGeomType,"Geometric information is not set up");
948 
949  // set up block matrix system
950  unsigned int nbdry = NumBndryCoeffs();
951  unsigned int nint = (unsigned int)(m_ncoeffs - nbdry);
952  unsigned int exp_size[] = {nbdry, nint};
953  unsigned int nblks = 2;
954  returnval = MemoryManager<DNekScalBlkMat>::AllocateSharedPtr(nblks, nblks, exp_size, exp_size); //Really need a constructor which takes Arrays
955  NekDouble factor = 1.0;
956 
957  switch(mkey.GetMatrixType())
958  {
960  case StdRegions::eHelmholtz: // special case since Helmholtz not defined in StdRegions
961 
962  // use Deformed case for both regular and deformed geometries
963  factor = 1.0;
964  goto UseLocRegionsMatrix;
965  break;
966  default:
967  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
968  {
969  factor = 1.0;
970  goto UseLocRegionsMatrix;
971  }
972  else
973  {
975  factor = mat->Scale();
976  goto UseStdRegionsMatrix;
977  }
978  break;
979  UseStdRegionsMatrix:
980  {
981  NekDouble invfactor = 1.0/factor;
982  NekDouble one = 1.0;
985  DNekMatSharedPtr Asubmat;
986 
987  //TODO: check below
988  returnval->SetBlock(0,0,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(factor,Asubmat = mat->GetBlock(0,0)));
989  returnval->SetBlock(0,1,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(one,Asubmat = mat->GetBlock(0,1)));
990  returnval->SetBlock(1,0,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(factor,Asubmat = mat->GetBlock(1,0)));
991  returnval->SetBlock(1,1,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(invfactor,Asubmat = mat->GetBlock(1,1)));
992  }
993  break;
994  UseLocRegionsMatrix:
995  {
996  int i,j;
997  NekDouble invfactor = 1.0/factor;
998  NekDouble one = 1.0;
999  DNekScalMat &mat = *GetLocMatrix(mkey);
1004 
1005  Array<OneD,unsigned int> bmap(nbdry);
1006  Array<OneD,unsigned int> imap(nint);
1007  GetBoundaryMap(bmap);
1008  GetInteriorMap(imap);
1009 
1010  for(i = 0; i < nbdry; ++i)
1011  {
1012  for(j = 0; j < nbdry; ++j)
1013  {
1014  (*A)(i,j) = mat(bmap[i],bmap[j]);
1015  }
1016 
1017  for(j = 0; j < nint; ++j)
1018  {
1019  (*B)(i,j) = mat(bmap[i],imap[j]);
1020  }
1021  }
1022 
1023  for(i = 0; i < nint; ++i)
1024  {
1025  for(j = 0; j < nbdry; ++j)
1026  {
1027  (*C)(i,j) = mat(imap[i],bmap[j]);
1028  }
1029 
1030  for(j = 0; j < nint; ++j)
1031  {
1032  (*D)(i,j) = mat(imap[i],imap[j]);
1033  }
1034  }
1035 
1036  // Calculate static condensed system
1037  if(nint)
1038  {
1039  D->Invert();
1040  (*B) = (*B)*(*D);
1041  (*A) = (*A) - (*B)*(*C);
1042  }
1043 
1044  DNekScalMatSharedPtr Atmp;
1045 
1046  returnval->SetBlock(0,0,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(factor,A));
1047  returnval->SetBlock(0,1,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(one,B));
1048  returnval->SetBlock(1,0,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(factor,C));
1049  returnval->SetBlock(1,1,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(invfactor,D));
1050 
1051  }
1052  }
1053  return returnval;
1054  }
1055 
1057  {
1058  if (m_metrics.count(MetricQuadrature) == 0)
1059  {
1061  }
1062 
1063  int i, j;
1064  const unsigned int nqtot = GetTotPoints();
1065  const unsigned int dim = 3;
1066  const MetricType m[3][3] = {
1070  };
1071 
1072  for (unsigned int i = 0; i < dim; ++i)
1073  {
1074  for (unsigned int j = i; j < dim; ++j)
1075  {
1076  m_metrics[m[i][j]] = Array<OneD, NekDouble>(nqtot);
1077  }
1078  }
1079 
1080  // Define shorthand synonyms for m_metrics storage
1081  Array<OneD,NekDouble> g0 (m_metrics[m[0][0]]);
1082  Array<OneD,NekDouble> g1 (m_metrics[m[1][1]]);
1083  Array<OneD,NekDouble> g2 (m_metrics[m[2][2]]);
1084  Array<OneD,NekDouble> g3 (m_metrics[m[0][1]]);
1085  Array<OneD,NekDouble> g4 (m_metrics[m[0][2]]);
1086  Array<OneD,NekDouble> g5 (m_metrics[m[1][2]]);
1087 
1088  // Allocate temporary storage
1089  Array<OneD,NekDouble> alloc(9*nqtot,0.0);
1090  Array<OneD,NekDouble> h0 (nqtot, alloc);
1091  Array<OneD,NekDouble> h1 (nqtot, alloc+ 1*nqtot);
1092  Array<OneD,NekDouble> h2 (nqtot, alloc+ 2*nqtot);
1093  Array<OneD,NekDouble> wsp1 (nqtot, alloc+ 3*nqtot);
1094  Array<OneD,NekDouble> wsp2 (nqtot, alloc+ 4*nqtot);
1095  Array<OneD,NekDouble> wsp3 (nqtot, alloc+ 5*nqtot);
1096  Array<OneD,NekDouble> wsp4 (nqtot, alloc+ 6*nqtot);
1097  Array<OneD,NekDouble> wsp5 (nqtot, alloc+ 7*nqtot);
1098  Array<OneD,NekDouble> wsp6 (nqtot, alloc+ 8*nqtot);
1099 
1100  const Array<TwoD, const NekDouble>& df =
1101  m_metricinfo->GetDerivFactors(GetPointsKeys());
1102  const Array<OneD, const NekDouble>& z0 = m_base[0]->GetZ();
1103  const Array<OneD, const NekDouble>& z1 = m_base[1]->GetZ();
1104  const Array<OneD, const NekDouble>& z2 = m_base[2]->GetZ();
1105  const unsigned int nquad0 = m_base[0]->GetNumPoints();
1106  const unsigned int nquad1 = m_base[1]->GetNumPoints();
1107  const unsigned int nquad2 = m_base[2]->GetNumPoints();
1108 
1109  // Populate collapsed coordinate arrays h0, h1 and h2.
1110  for(j = 0; j < nquad2; ++j)
1111  {
1112  for(i = 0; i < nquad1; ++i)
1113  {
1114  Vmath::Fill(nquad0, 2.0/(1.0-z2[j]), &h0[0]+i*nquad0 + j*nquad0*nquad1,1);
1115  Vmath::Fill(nquad0, 1.0/(1.0-z2[j]), &h1[0]+i*nquad0 + j*nquad0*nquad1,1);
1116  Vmath::Fill(nquad0, (1.0+z1[i])/(1.0-z2[j]), &h2[0]+i*nquad0 + j*nquad0*nquad1,1);
1117  }
1118  }
1119  for(i = 0; i < nquad0; i++)
1120  {
1121  Blas::Dscal(nquad1*nquad2, 1+z0[i], &h1[0]+i, nquad0);
1122  }
1123 
1124  // Step 3. Construct combined metric terms for physical space to
1125  // collapsed coordinate system.
1126  // Order of construction optimised to minimise temporary storage
1127  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1128  {
1129  // f_{1k}
1130  Vmath::Vvtvvtp(nqtot, &df[0][0], 1, &h0[0], 1, &df[2][0], 1, &h1[0], 1, &wsp1[0], 1);
1131  Vmath::Vvtvvtp(nqtot, &df[3][0], 1, &h0[0], 1, &df[5][0], 1, &h1[0], 1, &wsp2[0], 1);
1132  Vmath::Vvtvvtp(nqtot, &df[6][0], 1, &h0[0], 1, &df[8][0], 1, &h1[0], 1, &wsp3[0], 1);
1133 
1134  // g0
1135  Vmath::Vvtvvtp(nqtot, &wsp1[0], 1, &wsp1[0], 1, &wsp2[0], 1, &wsp2[0], 1, &g0[0], 1);
1136  Vmath::Vvtvp (nqtot, &wsp3[0], 1, &wsp3[0], 1, &g0[0], 1, &g0[0], 1);
1137 
1138  // g4
1139  Vmath::Vvtvvtp(nqtot, &df[2][0], 1, &wsp1[0], 1, &df[5][0], 1, &wsp2[0], 1, &g4[0], 1);
1140  Vmath::Vvtvp (nqtot, &df[8][0], 1, &wsp3[0], 1, &g4[0], 1, &g4[0], 1);
1141 
1142  // f_{2k}
1143  Vmath::Vvtvvtp(nqtot, &df[1][0], 1, &h0[0], 1, &df[2][0], 1, &h2[0], 1, &wsp4[0], 1);
1144  Vmath::Vvtvvtp(nqtot, &df[4][0], 1, &h0[0], 1, &df[5][0], 1, &h2[0], 1, &wsp5[0], 1);
1145  Vmath::Vvtvvtp(nqtot, &df[7][0], 1, &h0[0], 1, &df[8][0], 1, &h2[0], 1, &wsp6[0], 1);
1146 
1147  // g1
1148  Vmath::Vvtvvtp(nqtot, &wsp4[0], 1, &wsp4[0], 1, &wsp5[0], 1, &wsp5[0], 1, &g1[0], 1);
1149  Vmath::Vvtvp (nqtot, &wsp6[0], 1, &wsp6[0], 1, &g1[0], 1, &g1[0], 1);
1150 
1151  // g3
1152  Vmath::Vvtvvtp(nqtot, &wsp1[0], 1, &wsp4[0], 1, &wsp2[0], 1, &wsp5[0], 1, &g3[0], 1);
1153  Vmath::Vvtvp (nqtot, &wsp3[0], 1, &wsp6[0], 1, &g3[0], 1, &g3[0], 1);
1154 
1155  // g5
1156  Vmath::Vvtvvtp(nqtot, &df[2][0], 1, &wsp4[0], 1, &df[5][0], 1, &wsp5[0], 1, &g5[0], 1);
1157  Vmath::Vvtvp (nqtot, &df[8][0], 1, &wsp6[0], 1, &g5[0], 1, &g5[0], 1);
1158 
1159  // g2
1160  Vmath::Vvtvvtp(nqtot, &df[2][0], 1, &df[2][0], 1, &df[5][0], 1, &df[5][0], 1, &g2[0], 1);
1161  Vmath::Vvtvp (nqtot, &df[8][0], 1, &df[8][0], 1, &g2[0], 1, &g2[0], 1);
1162  }
1163  else
1164  {
1165  // f_{1k}
1166  Vmath::Svtsvtp(nqtot, df[0][0], &h0[0], 1, df[2][0], &h1[0], 1, &wsp1[0], 1);
1167  Vmath::Svtsvtp(nqtot, df[3][0], &h0[0], 1, df[5][0], &h1[0], 1, &wsp2[0], 1);
1168  Vmath::Svtsvtp(nqtot, df[6][0], &h0[0], 1, df[8][0], &h1[0], 1, &wsp3[0], 1);
1169 
1170  // g0
1171  Vmath::Vvtvvtp(nqtot, &wsp1[0], 1, &wsp1[0], 1, &wsp2[0], 1, &wsp2[0], 1, &g0[0], 1);
1172  Vmath::Vvtvp (nqtot, &wsp3[0], 1, &wsp3[0], 1, &g0[0], 1, &g0[0], 1);
1173 
1174  // g4
1175  Vmath::Svtsvtp(nqtot, df[2][0], &wsp1[0], 1, df[5][0], &wsp2[0], 1, &g4[0], 1);
1176  Vmath::Svtvp (nqtot, df[8][0], &wsp3[0], 1, &g4[0], 1, &g4[0], 1);
1177 
1178  // f_{2k}
1179  Vmath::Svtsvtp(nqtot, df[1][0], &h0[0], 1, df[2][0], &h2[0], 1, &wsp4[0], 1);
1180  Vmath::Svtsvtp(nqtot, df[4][0], &h0[0], 1, df[5][0], &h2[0], 1, &wsp5[0], 1);
1181  Vmath::Svtsvtp(nqtot, df[7][0], &h0[0], 1, df[8][0], &h2[0], 1, &wsp6[0], 1);
1182 
1183  // g1
1184  Vmath::Vvtvvtp(nqtot, &wsp4[0], 1, &wsp4[0], 1, &wsp5[0], 1, &wsp5[0], 1, &g1[0], 1);
1185  Vmath::Vvtvp (nqtot, &wsp6[0], 1, &wsp6[0], 1, &g1[0], 1, &g1[0], 1);
1186 
1187  // g3
1188  Vmath::Vvtvvtp(nqtot, &wsp1[0], 1, &wsp4[0], 1, &wsp2[0], 1, &wsp5[0], 1, &g3[0], 1);
1189  Vmath::Vvtvp (nqtot, &wsp3[0], 1, &wsp6[0], 1, &g3[0], 1, &g3[0], 1);
1190 
1191  // g5
1192  Vmath::Svtsvtp(nqtot, df[2][0], &wsp4[0], 1, df[5][0], &wsp5[0], 1, &g5[0], 1);
1193  Vmath::Svtvp (nqtot, df[8][0], &wsp6[0], 1, &g5[0], 1, &g5[0], 1);
1194 
1195  // g2
1196  Vmath::Fill(nqtot, df[2][0]*df[2][0] + df[5][0]*df[5][0] + df[8][0]*df[8][0], &g2[0], 1);
1197  }
1198 
1199  for (unsigned int i = 0; i < dim; ++i)
1200  {
1201  for (unsigned int j = i; j < dim; ++j)
1202  {
1204  m_metrics[m[i][j]]);
1205 
1206  }
1207  }
1208  }
1209 
1211  const Array<OneD, const NekDouble> &inarray,
1212  Array<OneD, NekDouble> &outarray,
1213  Array<OneD, NekDouble> &wsp)
1214  {
1215  // This implementation is only valid when there are no coefficients
1216  // associated to the Laplacian operator
1217  if (m_metrics.count(MetricLaplacian00) == 0)
1218  {
1220  }
1221 
1222  int nquad0 = m_base[0]->GetNumPoints();
1223  int nquad1 = m_base[1]->GetNumPoints();
1224  int nq2 = m_base[2]->GetNumPoints();
1225  int nqtot = nquad0*nquad1*nq2;
1226 
1227  ASSERTL1(wsp.num_elements() >= 6*nqtot,
1228  "Insufficient workspace size.");
1229  ASSERTL1(m_ncoeffs <= nqtot,
1230  "Workspace not set up for ncoeffs > nqtot");
1231 
1232  const Array<OneD, const NekDouble>& base0 = m_base[0]->GetBdata();
1233  const Array<OneD, const NekDouble>& base1 = m_base[1]->GetBdata();
1234  const Array<OneD, const NekDouble>& base2 = m_base[2]->GetBdata();
1235  const Array<OneD, const NekDouble>& dbase0 = m_base[0]->GetDbdata();
1236  const Array<OneD, const NekDouble>& dbase1 = m_base[1]->GetDbdata();
1237  const Array<OneD, const NekDouble>& dbase2 = m_base[2]->GetDbdata();
1238  const Array<OneD, const NekDouble>& metric00 = m_metrics[MetricLaplacian00];
1239  const Array<OneD, const NekDouble>& metric01 = m_metrics[MetricLaplacian01];
1240  const Array<OneD, const NekDouble>& metric02 = m_metrics[MetricLaplacian02];
1241  const Array<OneD, const NekDouble>& metric11 = m_metrics[MetricLaplacian11];
1242  const Array<OneD, const NekDouble>& metric12 = m_metrics[MetricLaplacian12];
1243  const Array<OneD, const NekDouble>& metric22 = m_metrics[MetricLaplacian22];
1244 
1245  // Allocate temporary storage
1246  Array<OneD,NekDouble> wsp0 (2*nqtot, wsp);
1247  Array<OneD,NekDouble> wsp1 ( nqtot, wsp+1*nqtot);
1248  Array<OneD,NekDouble> wsp2 ( nqtot, wsp+2*nqtot);
1249  Array<OneD,NekDouble> wsp3 ( nqtot, wsp+3*nqtot);
1250  Array<OneD,NekDouble> wsp4 ( nqtot, wsp+4*nqtot);
1251  Array<OneD,NekDouble> wsp5 ( nqtot, wsp+5*nqtot);
1252 
1253  // LAPLACIAN MATRIX OPERATION
1254  // wsp1 = du_dxi1 = D_xi1 * inarray = D_xi1 * u
1255  // wsp2 = du_dxi2 = D_xi2 * inarray = D_xi2 * u
1256  // wsp2 = du_dxi3 = D_xi3 * inarray = D_xi3 * u
1257  StdExpansion3D::PhysTensorDeriv(inarray,wsp0,wsp1,wsp2);
1258 
1259  // wsp0 = k = g0 * wsp1 + g1 * wsp2 = g0 * du_dxi1 + g1 * du_dxi2
1260  // wsp2 = l = g1 * wsp1 + g2 * wsp2 = g0 * du_dxi1 + g1 * du_dxi2
1261  // where g0, g1 and g2 are the metric terms set up in the GeomFactors class
1262  // especially for this purpose
1263  Vmath::Vvtvvtp(nqtot,&metric00[0],1,&wsp0[0],1,&metric01[0],1,&wsp1[0],1,&wsp3[0],1);
1264  Vmath::Vvtvp (nqtot,&metric02[0],1,&wsp2[0],1,&wsp3[0],1,&wsp3[0],1);
1265  Vmath::Vvtvvtp(nqtot,&metric01[0],1,&wsp0[0],1,&metric11[0],1,&wsp1[0],1,&wsp4[0],1);
1266  Vmath::Vvtvp (nqtot,&metric12[0],1,&wsp2[0],1,&wsp4[0],1,&wsp4[0],1);
1267  Vmath::Vvtvvtp(nqtot,&metric02[0],1,&wsp0[0],1,&metric12[0],1,&wsp1[0],1,&wsp5[0],1);
1268  Vmath::Vvtvp (nqtot,&metric22[0],1,&wsp2[0],1,&wsp5[0],1,&wsp5[0],1);
1269 
1270  // outarray = m = (D_xi1 * B)^T * k
1271  // wsp1 = n = (D_xi2 * B)^T * l
1272  IProductWRTBase_SumFacKernel(dbase0,base1,base2,wsp3,outarray,wsp0,false,true,true);
1273  IProductWRTBase_SumFacKernel(base0,dbase1,base2,wsp4,wsp2, wsp0,true,false,true);
1274  Vmath::Vadd(m_ncoeffs,wsp2.get(),1,outarray.get(),1,outarray.get(),1);
1275  IProductWRTBase_SumFacKernel(base0,base1,dbase2,wsp5,wsp2, wsp0,true,true,false);
1276  Vmath::Vadd(m_ncoeffs,wsp2.get(),1,outarray.get(),1,outarray.get(),1);
1277  }
1278  }//end of namespace
1279 }//end of namespace