Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
TetExp.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File TetExp.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:
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
36 #include <LocalRegions/TetExp.h>
37 #include <SpatialDomains/SegGeom.h>
40 
41 namespace Nektar
42 {
43  namespace LocalRegions
44  {
45  /**
46  * @class TetExp
47  * Defines a Tetrahedral local expansion.
48  */
49 
50  /**
51  * \brief Constructor using BasisKey class for quadrature points and
52  * order definition
53  *
54  * @param Ba Basis key for first coordinate.
55  * @param Bb Basis key for second coordinate.
56  * @param Bc Basis key for third coordinate.
57  */
59  const LibUtilities::BasisKey &Bb,
60  const LibUtilities::BasisKey &Bc,
62  ):
63  StdExpansion (LibUtilities::StdTetData::getNumberOfCoefficients(Ba.GetNumModes(),Bb.GetNumModes(),Bc.GetNumModes()),3,Ba,Bb,Bc),
64  StdExpansion3D(LibUtilities::StdTetData::getNumberOfCoefficients(Ba.GetNumModes(),Bb.GetNumModes(),Bc.GetNumModes()),Ba,Bb,Bc),
65  StdRegions::StdTetExp(Ba,Bb,Bc),
66  Expansion (geom),
67  Expansion3D (geom),
68  m_matrixManager(
69  boost::bind(&TetExp::CreateMatrix, this, _1),
70  std::string("TetExpMatrix")),
71  m_staticCondMatrixManager(
72  boost::bind(&TetExp::CreateStaticCondMatrix, this, _1),
73  std::string("TetExpStaticCondMatrix"))
74  {
75  }
76 
77 
78  /**
79  * \brief Copy Constructor
80  */
82  StdExpansion(T),
83  StdExpansion3D(T),
84  StdRegions::StdTetExp(T),
85  Expansion(T),
86  Expansion3D(T),
87  m_matrixManager(T.m_matrixManager),
88  m_staticCondMatrixManager(T.m_staticCondMatrixManager)
89  {
90  }
91 
92  /**
93  * \brief Destructor
94  */
96  {
97  }
98 
99 
100  //-----------------------------
101  // Integration Methods
102  //-----------------------------
103  /**
104  * \brief Integrate the physical point list \a inarray over region
105  *
106  * @param inarray Definition of function to be returned at
107  * quadrature point of expansion.
108  * @returns \f$\int^1_{-1}\int^1_{-1} \int^1_{-1}
109  * u(\eta_1, \eta_2, \eta_3) J[i,j,k] d \eta_1 d \eta_2 d \eta_3 \f$
110  * where \f$inarray[i,j,k] = u(\eta_{1i},\eta_{2j},\eta_{3k})
111  * \f$ and \f$ J[i,j,k] \f$ is the Jacobian evaluated at the quadrature
112  * point.
113  */
115  const Array<OneD, const NekDouble> &inarray)
116  {
117  int nquad0 = m_base[0]->GetNumPoints();
118  int nquad1 = m_base[1]->GetNumPoints();
119  int nquad2 = m_base[2]->GetNumPoints();
120  Array<OneD, const NekDouble> jac = m_metricinfo->GetJac(GetPointsKeys());
121  NekDouble retrunVal;
122  Array<OneD,NekDouble> tmp(nquad0*nquad1*nquad2);
123 
124  // multiply inarray with Jacobian
125  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
126  {
127  Vmath::Vmul(nquad0*nquad1*nquad2,&jac[0],1,
128  (NekDouble*)&inarray[0],1, &tmp[0],1);
129  }
130  else
131  {
132  Vmath::Smul(nquad0*nquad1*nquad2,(NekDouble) jac[0],
133  (NekDouble*)&inarray[0],1,&tmp[0],1);
134  }
135 
136  // call StdTetExp version;
137  retrunVal = StdTetExp::v_Integral(tmp);
138 
139  return retrunVal;
140  }
141 
142 
143  //-----------------------------
144  // Differentiation Methods
145  //-----------------------------
146  /**
147  * \brief Differentiate \a inarray in the three coordinate directions.
148  *
149  * @param inarray Input array of values at quadrature points to
150  * be differentiated.
151  * @param out_d0 Derivative in first coordinate direction.
152  * @param out_d1 Derivative in second coordinate direction.
153  * @param out_d2 Derivative in third coordinate direction.
154  */
156  const Array<OneD, const NekDouble> &inarray,
157  Array<OneD, NekDouble> &out_d0,
158  Array<OneD, NekDouble> &out_d1,
159  Array<OneD, NekDouble> &out_d2)
160  {
161  int TotPts = m_base[0]->GetNumPoints()*m_base[1]->GetNumPoints()*
162  m_base[2]->GetNumPoints();
163 
164  Array<TwoD, const NekDouble> df =
165  m_metricinfo->GetDerivFactors(GetPointsKeys());
166  Array<OneD,NekDouble> Diff0 = Array<OneD,NekDouble>(3*TotPts);
167  Array<OneD,NekDouble> Diff1 = Diff0 + TotPts;
168  Array<OneD,NekDouble> Diff2 = Diff1 + TotPts;
169 
170  StdTetExp::v_PhysDeriv(inarray, Diff0, Diff1, Diff2);
171 
172  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
173  {
174  if(out_d0.num_elements())
175  {
176  Vmath::Vmul (TotPts,&df[0][0],1,&Diff0[0],1, &out_d0[0], 1);
177  Vmath::Vvtvp (TotPts,&df[1][0],1,&Diff1[0],1, &out_d0[0], 1,&out_d0[0],1);
178  Vmath::Vvtvp (TotPts,&df[2][0],1,&Diff2[0],1, &out_d0[0], 1,&out_d0[0],1);
179  }
180 
181  if(out_d1.num_elements())
182  {
183  Vmath::Vmul (TotPts,&df[3][0],1,&Diff0[0],1, &out_d1[0], 1);
184  Vmath::Vvtvp (TotPts,&df[4][0],1,&Diff1[0],1, &out_d1[0], 1,&out_d1[0],1);
185  Vmath::Vvtvp (TotPts,&df[5][0],1,&Diff2[0],1, &out_d1[0], 1,&out_d1[0],1);
186  }
187 
188  if(out_d2.num_elements())
189  {
190  Vmath::Vmul (TotPts,&df[6][0],1,&Diff0[0],1, &out_d2[0], 1);
191  Vmath::Vvtvp (TotPts,&df[7][0],1,&Diff1[0],1, &out_d2[0], 1, &out_d2[0],1);
192  Vmath::Vvtvp (TotPts,&df[8][0],1,&Diff2[0],1, &out_d2[0], 1, &out_d2[0],1);
193  }
194  }
195  else // regular geometry
196  {
197  if(out_d0.num_elements())
198  {
199  Vmath::Smul (TotPts,df[0][0],&Diff0[0],1, &out_d0[0], 1);
200  Blas::Daxpy (TotPts,df[1][0],&Diff1[0],1, &out_d0[0], 1);
201  Blas::Daxpy (TotPts,df[2][0],&Diff2[0],1, &out_d0[0], 1);
202  }
203 
204  if(out_d1.num_elements())
205  {
206  Vmath::Smul (TotPts,df[3][0],&Diff0[0],1, &out_d1[0], 1);
207  Blas::Daxpy (TotPts,df[4][0],&Diff1[0],1, &out_d1[0], 1);
208  Blas::Daxpy (TotPts,df[5][0],&Diff2[0],1, &out_d1[0], 1);
209  }
210 
211  if(out_d2.num_elements())
212  {
213  Vmath::Smul (TotPts,df[6][0],&Diff0[0],1, &out_d2[0], 1);
214  Blas::Daxpy (TotPts,df[7][0],&Diff1[0],1, &out_d2[0], 1);
215  Blas::Daxpy (TotPts,df[8][0],&Diff2[0],1, &out_d2[0], 1);
216  }
217  }
218  }
219 
220 
221  //-----------------------------
222  // Transforms
223  //-----------------------------
224  /**
225  * \brief Forward transform from physical quadrature space stored in
226  * \a inarray and evaluate the expansion coefficients and store
227  * in \a (this)->_coeffs
228  *
229  * @param inarray Array of physical quadrature points to be
230  * transformed.
231  * @param outarray Array of coefficients to update.
232  */
234  const Array<OneD, const NekDouble> & inarray,
235  Array<OneD,NekDouble> &outarray)
236  {
237  if((m_base[0]->Collocation())&&(m_base[1]->Collocation())&&(m_base[2]->Collocation()))
238  {
239  Vmath::Vcopy(GetNcoeffs(),&inarray[0],1,&outarray[0],1);
240  }
241  else
242  {
243  IProductWRTBase(inarray,outarray);
244 
245  // get Mass matrix inverse
247  DetShapeType(),*this);
248  DNekScalMatSharedPtr matsys = m_matrixManager[masskey];
249 
250  // copy inarray in case inarray == outarray
251  DNekVec in (m_ncoeffs,outarray);
252  DNekVec out(m_ncoeffs,outarray,eWrapper);
253 
254  out = (*matsys)*in;
255  }
256  }
257 
258  //-----------------------------
259  // Inner product functions
260  //-----------------------------
261  /**
262  * \brief Calculate the inner product of inarray with respect to the
263  * basis B=m_base0*m_base1*m_base2 and put into outarray:
264  *
265  * \f$ \begin{array}{rcl} I_{pqr} = (\phi_{pqr}, u)_{\delta}
266  * & = & \sum_{i=0}^{nq_0} \sum_{j=0}^{nq_1} \sum_{k=0}^{nq_2}
267  * \psi_{p}^{a} (\eta_{1i}) \psi_{pq}^{b} (\eta_{2j}) \psi_{pqr}^{c}
268  * (\eta_{3k}) w_i w_j w_k u(\eta_{1,i} \eta_{2,j} \eta_{3,k})
269  * J_{i,j,k}\\ & = & \sum_{i=0}^{nq_0} \psi_p^a(\eta_{1,i})
270  * \sum_{j=0}^{nq_1} \psi_{pq}^b(\eta_{2,j}) \sum_{k=0}^{nq_2}
271  * \psi_{pqr}^c u(\eta_{1i},\eta_{2j},\eta_{3k}) J_{i,j,k}
272  * \end{array} \f$ \n
273  * where
274  * \f$ \phi_{pqr} (\xi_1 , \xi_2 , \xi_3)
275  * = \psi_p^a (\eta_1) \psi_{pq}^b (\eta_2) \psi_{pqr}^c (\eta_3) \f$
276  * which can be implemented as \n
277  * \f$f_{pqr} (\xi_{3k})
278  * = \sum_{k=0}^{nq_3} \psi_{pqr}^c u(\eta_{1i},\eta_{2j},\eta_{3k})
279  * J_{i,j,k} = {\bf B_3 U} \f$ \n
280  * \f$ g_{pq} (\xi_{3k})
281  * = \sum_{j=0}^{nq_1} \psi_{pq}^b (\xi_{2j}) f_{pqr} (\xi_{3k})
282  * = {\bf B_2 F} \f$ \n
283  * \f$ (\phi_{pqr}, u)_{\delta}
284  * = \sum_{k=0}^{nq_0} \psi_{p}^a (\xi_{3k}) g_{pq} (\xi_{3k})
285  * = {\bf B_1 G} \f$
286  */
288  const Array<OneD, const NekDouble> &inarray,
289  Array<OneD, NekDouble> &outarray)
290  {
291  v_IProductWRTBase_SumFac(inarray, outarray);
292  }
293 
295  const Array<OneD, const NekDouble> &inarray,
296  Array<OneD, NekDouble> &outarray)
297  {
298  const int nquad0 = m_base[0]->GetNumPoints();
299  const int nquad1 = m_base[1]->GetNumPoints();
300  const int nquad2 = m_base[2]->GetNumPoints();
301  const int order0 = m_base[0]->GetNumModes();
302  const int order1 = m_base[1]->GetNumModes();
303  Array<OneD, NekDouble> wsp(nquad1*nquad2*order0 +
304  nquad2*order0*(order1+1)/2);
305  Array<OneD, NekDouble> tmp(nquad0*nquad1*nquad2);
306 
307  MultiplyByQuadratureMetric(inarray, tmp);
308  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata(),
309  m_base[1]->GetBdata(),
310  m_base[2]->GetBdata(),
311  tmp,outarray,wsp,
312  true,true,true);
313  }
314 
315  /**
316  * @brief Calculates the inner product \f$ I_{pqr} = (u,
317  * \partial_{x_i} \phi_{pqr}) \f$.
318  *
319  * The derivative of the basis functions is performed using the chain
320  * rule in order to incorporate the geometric factors. Assuming that
321  * the basis functions are a tensor product
322  * \f$\phi_{pqr}(\eta_1,\eta_2,\eta_3) =
323  * \phi_1(\eta_1)\phi_2(\eta_2)\phi_3(\eta_3)\f$, this yields the
324  * result
325  *
326  * \f[
327  * I_{pqr} = \sum_{j=1}^3 \left(u, \frac{\partial u}{\partial \eta_j}
328  * \frac{\partial \eta_j}{\partial x_i}\right)
329  * \f]
330  *
331  * In the prismatic element, we must also incorporate a second set of
332  * geometric factors which incorporate the collapsed co-ordinate
333  * system, so that
334  *
335  * \f[ \frac{\partial\eta_j}{\partial x_i} = \sum_{k=1}^3
336  * \frac{\partial\eta_j}{\partial\xi_k}\frac{\partial\xi_k}{\partial
337  * x_i} \f]
338  *
339  * These derivatives can be found on p152 of Sherwin & Karniadakis.
340  *
341  * @param dir Direction in which to take the derivative.
342  * @param inarray The function \f$ u \f$.
343  * @param outarray Value of the inner product.
344  */
346  const int dir,
347  const Array<OneD, const NekDouble> &inarray,
348  Array<OneD, NekDouble> &outarray)
349  {
350  const int nquad0 = m_base[0]->GetNumPoints();
351  const int nquad1 = m_base[1]->GetNumPoints();
352  const int nquad2 = m_base[2]->GetNumPoints();
353  const int order0 = m_base[0]->GetNumModes ();
354  const int order1 = m_base[1]->GetNumModes ();
355  const int nqtot = nquad0*nquad1*nquad2;
356  int i, j;
357 
358  const Array<OneD, const NekDouble> &z0 = m_base[0]->GetZ();
359  const Array<OneD, const NekDouble> &z1 = m_base[1]->GetZ();
360  const Array<OneD, const NekDouble> &z2 = m_base[2]->GetZ();
361 
362  Array<OneD, NekDouble> h0 (nqtot);
363  Array<OneD, NekDouble> h1 (nqtot);
364  Array<OneD, NekDouble> h2 (nqtot);
365  Array<OneD, NekDouble> h3 (nqtot);
366  Array<OneD, NekDouble> tmp1 (nqtot);
367  Array<OneD, NekDouble> tmp2 (nqtot);
368  Array<OneD, NekDouble> tmp3 (nqtot);
369  Array<OneD, NekDouble> tmp4 (nqtot);
370  Array<OneD, NekDouble> tmp5 (nqtot);
371  Array<OneD, NekDouble> tmp6 (m_ncoeffs);
372  Array<OneD, NekDouble> wsp (nquad1*nquad2*order0 +
373  nquad2*order0*(order1+1)/2);
374 
375  const Array<TwoD, const NekDouble>& df =
376  m_metricinfo->GetDerivFactors(GetPointsKeys());
377 
378  MultiplyByQuadratureMetric(inarray,tmp1);
379 
380  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
381  {
382  Vmath::Vmul(nqtot,&df[3*dir][0], 1,tmp1.get(),1,tmp2.get(),1);
383  Vmath::Vmul(nqtot,&df[3*dir+1][0],1,tmp1.get(),1,tmp3.get(),1);
384  Vmath::Vmul(nqtot,&df[3*dir+2][0],1,tmp1.get(),1,tmp4.get(),1);
385  }
386  else
387  {
388  Vmath::Smul(nqtot, df[3*dir ][0],tmp1.get(),1,tmp2.get(), 1);
389  Vmath::Smul(nqtot, df[3*dir+1][0],tmp1.get(),1,tmp3.get(), 1);
390  Vmath::Smul(nqtot, df[3*dir+2][0],tmp1.get(),1,tmp4.get(), 1);
391  }
392 
393  const int nq01 = nquad0*nquad1;
394  const int nq12 = nquad1*nquad2;
395 
396  for(j = 0; j < nquad2; ++j)
397  {
398  for(i = 0; i < nquad1; ++i)
399  {
400  Vmath::Fill(nquad0, 4.0/(1.0-z1[i])/(1.0-z2[j]),
401  &h0[0]+i*nquad0 + j*nq01,1);
402  Vmath::Fill(nquad0, 2.0/(1.0-z1[i])/(1.0-z2[j]),
403  &h1[0]+i*nquad0 + j*nq01,1);
404  Vmath::Fill(nquad0, 2.0/(1.0-z2[j]),
405  &h2[0]+i*nquad0 + j*nq01,1);
406  Vmath::Fill(nquad0, (1.0+z1[i])/(1.0-z2[j]),
407  &h3[0]+i*nquad0 + j*nq01,1);
408  }
409  }
410 
411  for(i = 0; i < nquad0; i++)
412  {
413  Blas::Dscal(nq12, 1+z0[i], &h1[0]+i, nquad0);
414  }
415 
416  // Assemble terms for first IP.
417  Vmath::Vvtvvtp(nqtot, &tmp2[0], 1, &h0[0], 1,
418  &tmp3[0], 1, &h1[0], 1,
419  &tmp5[0], 1);
420  Vmath::Vvtvp (nqtot, &tmp4[0], 1, &h1[0], 1,
421  &tmp5[0], 1, &tmp5[0], 1);
422 
423  IProductWRTBase_SumFacKernel(m_base[0]->GetDbdata(),
424  m_base[1]->GetBdata (),
425  m_base[2]->GetBdata (),
426  tmp5,outarray,wsp,
427  true,true,true);
428 
429  // Assemble terms for second IP.
430  Vmath::Vvtvvtp(nqtot, &tmp3[0], 1, &h2[0], 1,
431  &tmp4[0], 1, &h3[0], 1,
432  &tmp5[0], 1);
433 
434  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata (),
435  m_base[1]->GetDbdata(),
436  m_base[2]->GetBdata (),
437  tmp5,tmp6,wsp,
438  true,true,true);
439  Vmath::Vadd(m_ncoeffs, tmp6, 1, outarray, 1, outarray, 1);
440 
441  // Do third IP.
442  IProductWRTBase_SumFacKernel(m_base[0]->GetBdata (),
443  m_base[1]->GetBdata (),
444  m_base[2]->GetDbdata(),
445  tmp4,tmp6,wsp,
446  true,true,true);
447 
448  // Sum contributions.
449  Vmath::Vadd(m_ncoeffs, tmp6, 1, outarray, 1, outarray, 1);
450  }
451 
452 
453  //-----------------------------
454  // Evaluation functions
455  //-----------------------------
456 
457  /**
458  * Given the local cartesian coordinate \a Lcoord evaluate the
459  * value of physvals at this point by calling through to the
460  * StdExpansion method
461  */
463  const Array<OneD, const NekDouble> &Lcoord,
464  const Array<OneD, const NekDouble> &physvals)
465  {
466  // Evaluate point in local (eta) coordinates.
467  return StdTetExp::v_PhysEvaluate(Lcoord,physvals);
468  }
469 
470  /**
471  * @param coord Physical space coordinate
472  * @returns Evaluation of expansion at given coordinate.
473  */
475  const Array<OneD, const NekDouble> &coord,
476  const Array<OneD, const NekDouble> & physvals)
477  {
478  ASSERTL0(m_geom,"m_geom not defined");
479 
480  Array<OneD,NekDouble> Lcoord = Array<OneD,NekDouble>(3);
481 
482  // Get the local (eta) coordinates of the point
483  m_geom->GetLocCoords(coord,Lcoord);
484 
485  // Evaluate point in local (eta) coordinates.
486  return StdTetExp::v_PhysEvaluate(Lcoord,physvals);
487  }
488 
489  /**
490  * \brief Get the coordinates "coords" at the local coordinates "Lcoords"
491  */
493  const Array<OneD, const NekDouble> &Lcoords,
494  Array<OneD,NekDouble> &coords)
495  {
496  int i;
497 
498  ASSERTL1(Lcoords[0] <= -1.0 && Lcoords[0] >= 1.0 &&
499  Lcoords[1] <= -1.0 && Lcoords[1] >= 1.0 &&
500  Lcoords[2] <= -1.0 && Lcoords[2] >= 1.0,
501  "Local coordinates are not in region [-1,1]");
502 
503  // m_geom->FillGeom(); // TODO: implement FillGeom()
504 
505  for(i = 0; i < m_geom->GetCoordim(); ++i)
506  {
507  coords[i] = m_geom->GetCoord(i,Lcoords);
508  }
509  }
510 
512  Array<OneD, NekDouble> &coords_0,
513  Array<OneD, NekDouble> &coords_1,
514  Array<OneD, NekDouble> &coords_2)
515  {
516  Expansion::v_GetCoords(coords_0, coords_1, coords_2);
517  }
518 
519 
520  //-----------------------------
521  // Helper functions
522  //-----------------------------
523 
524  /**
525  * \brief Return Shape of region, using ShapeType enum list.
526  */
528  {
530  }
531 
533  {
534  return m_geom->GetCoordim();
535  }
536 
537 
538 
540  const NekDouble *data,
541  const std::vector<unsigned int > &nummodes,
542  const int mode_offset,
543  NekDouble * coeffs)
544  {
545  int data_order0 = nummodes[mode_offset];
546  int fillorder0 = min(m_base[0]->GetNumModes(),data_order0);
547  int data_order1 = nummodes[mode_offset+1];
548  int order1 = m_base[1]->GetNumModes();
549  int fillorder1 = min(order1,data_order1);
550  int data_order2 = nummodes[mode_offset+2];
551  int order2 = m_base[2]->GetNumModes();
552  int fillorder2 = min(order2,data_order2);
553 
554  switch(m_base[0]->GetBasisType())
555  {
557  {
558  int i,j;
559  int cnt = 0;
560  int cnt1 = 0;
561 
562  ASSERTL1(m_base[1]->GetBasisType() ==
564  "Extraction routine not set up for this basis");
565  ASSERTL1(m_base[2]->GetBasisType() ==
567  "Extraction routine not set up for this basis");
568 
569  Vmath::Zero(m_ncoeffs,coeffs,1);
570  for(j = 0; j < fillorder0; ++j)
571  {
572  for(i = 0; i < fillorder1-j; ++i)
573  {
574  Vmath::Vcopy(fillorder2-j-i, &data[cnt], 1,
575  &coeffs[cnt1], 1);
576  cnt += data_order2-j-i;
577  cnt1 += order2-j-i;
578  }
579 
580  // count out data for j iteration
581  for(i = fillorder1-j; i < data_order1-j; ++i)
582  {
583  cnt += data_order2-j-i;
584  }
585 
586  for(i = fillorder1-j; i < order1-j; ++i)
587  {
588  cnt1 += order2-j-i;
589  }
590 
591  }
592  }
593  break;
594  default:
595  ASSERTL0(false, "basis is either not set up or not "
596  "hierarchicial");
597  }
598  }
599 
600  /**
601  * \brief Returns the physical values at the quadrature points of a face
602  * Wrapper function to v_GetFacePhysVals
603  */
605  const int face,
606  const StdRegions::StdExpansionSharedPtr &FaceExp,
607  const Array<OneD, const NekDouble> &inarray,
608  Array<OneD, NekDouble> &outarray,
610  {
611  v_GetFacePhysVals(face,FaceExp,inarray,outarray,orient);
612  }
613 
614  /**
615  * \brief Returns the physical values at the quadrature points of a face
616  */
618  const int face,
619  const StdRegions::StdExpansionSharedPtr &FaceExp,
620  const Array<OneD, const NekDouble> &inarray,
621  Array<OneD, NekDouble> &outarray,
623  {
624  int nquad0 = m_base[0]->GetNumPoints();
625  int nquad1 = m_base[1]->GetNumPoints();
626  int nquad2 = m_base[2]->GetNumPoints();
627 
628  Array<OneD,NekDouble> o_tmp (GetFaceNumPoints(face));
629  Array<OneD,NekDouble> o_tmp2(FaceExp->GetTotPoints());
630  Array<OneD,NekDouble> o_tmp3;
631 
632  if (orient == StdRegions::eNoOrientation)
633  {
634  orient = GetForient(face);
635  }
636 
637  switch(face)
638  {
639  case 0:
640  {
641  //Directions A and B positive
642  Vmath::Vcopy(nquad0*nquad1,inarray.get(),1,o_tmp.get(),1);
643  //interpolate
644  LibUtilities::Interp2D(m_base[0]->GetPointsKey(),
645  m_base[1]->GetPointsKey(),
646  o_tmp.get(),
647  FaceExp->GetBasis(0)->GetPointsKey(),
648  FaceExp->GetBasis(1)->GetPointsKey(),
649  o_tmp2.get());
650  break;
651  }
652  case 1:
653  {
654  //Direction A and B positive
655  for (int k=0; k<nquad2; k++)
656  {
657  Vmath::Vcopy(nquad0,inarray.get()+(nquad0*nquad1*k),1,o_tmp.get()+(k*nquad0),1);
658  }
659  //interpolate
660  LibUtilities::Interp2D(m_base[0]->GetPointsKey(),
661  m_base[2]->GetPointsKey(),
662  o_tmp.get(),
663  FaceExp->GetBasis(0)->GetPointsKey(),
664  FaceExp->GetBasis(1)->GetPointsKey(),
665  o_tmp2.get());
666  break;
667  }
668  case 2:
669  {
670  //Directions A and B positive
671  Vmath::Vcopy(nquad1*nquad2,inarray.get()+(nquad0-1),nquad0,o_tmp.get(),1);
672  //interpolate
673  LibUtilities::Interp2D(m_base[1]->GetPointsKey(), m_base[2]->GetPointsKey(), o_tmp.get(),
674  FaceExp->GetBasis(0)->GetPointsKey(),FaceExp->GetBasis(1)->GetPointsKey(),o_tmp2.get());
675  break;
676  }
677  case 3:
678  {
679  //Directions A and B positive
680  Vmath::Vcopy(nquad1*nquad2,inarray.get(),nquad0,o_tmp.get(),1);
681  //interpolate
682  LibUtilities::Interp2D(m_base[1]->GetPointsKey(), m_base[2]->GetPointsKey(), o_tmp.get(),
683  FaceExp->GetBasis(0)->GetPointsKey(),FaceExp->GetBasis(1)->GetPointsKey(),o_tmp2.get());
684  }
685  break;
686  default:
687  ASSERTL0(false,"face value (> 3) is out of range");
688  break;
689  }
690 
691  int nq1 = FaceExp->GetNumPoints(0);
692  int nq2 = FaceExp->GetNumPoints(1);
693 
694  if ((int)orient == 7)
695  {
696  for (int j = 0; j < nq2; ++j)
697  {
698  Vmath::Vcopy(nq1, o_tmp2.get()+((j+1)*nq1-1), -1, outarray.get()+j*nq1, 1);
699  }
700  }
701  else
702  {
703  Vmath::Vcopy(nq1*nq2, o_tmp2.get(), 1, outarray.get(), 1);
704  }
705  }
706 
707 
708  /**
709  * \brief Compute the normal of a triangular face
710  */
711  void TetExp::v_ComputeFaceNormal(const int face)
712  {
713  int i;
714  const SpatialDomains::GeomFactorsSharedPtr &geomFactors =
715  GetGeom()->GetMetricInfo();
717  SpatialDomains::GeomType type = geomFactors->GetGtype();
718  const Array<TwoD, const NekDouble> &df = geomFactors->GetDerivFactors(ptsKeys);
719  const Array<OneD, const NekDouble> &jac = geomFactors->GetJac(ptsKeys);
720 
721 
722  LibUtilities::BasisKey tobasis0 = DetFaceBasisKey(face,0);
723  LibUtilities::BasisKey tobasis1 = DetFaceBasisKey(face,1);
724 
725  // number of face quadrature points
726  int nq_face = tobasis0.GetNumPoints()*tobasis1.GetNumPoints();
727 
728  int vCoordDim = GetCoordim();
729 
730  m_faceNormals[face] = Array<OneD, Array<OneD, NekDouble> >(vCoordDim);
731  Array<OneD, Array<OneD, NekDouble> > &normal = m_faceNormals[face];
732  for (i = 0; i < vCoordDim; ++i)
733  {
734  normal[i] = Array<OneD, NekDouble>(nq_face);
735  }
736 
737  // Regular geometry case
738  if (type == SpatialDomains::eRegular ||
740  {
741  NekDouble fac;
742 
743  // Set up normals
744  switch (face)
745  {
746  case 0:
747  {
748  for (i = 0; i < vCoordDim; ++i)
749  {
750  normal[i][0] = -df[3*i+2][0];
751  }
752 
753  break;
754  }
755  case 1:
756  {
757  for (i = 0; i < vCoordDim; ++i)
758  {
759  normal[i][0] = -df[3*i+1][0];
760  }
761 
762  break;
763  }
764  case 2:
765  {
766  for (i = 0; i < vCoordDim; ++i)
767  {
768  normal[i][0] = df[3*i][0]+df[3*i+1][0]+
769  df[3*i+2][0];
770  }
771 
772  break;
773  }
774  case 3:
775  {
776  for(i = 0; i < vCoordDim; ++i)
777  {
778  normal[i][0] = -df[3*i][0];
779  }
780  break;
781  }
782  default:
783  ASSERTL0(false,"face is out of range (edge < 3)");
784  }
785 
786  // normalise
787  fac = 0.0;
788  for (i = 0; i < vCoordDim; ++i)
789  {
790  fac += normal[i][0]*normal[i][0];
791  }
792  fac = 1.0/sqrt(fac);
793 
794  for (i = 0; i < vCoordDim; ++i)
795  {
796  Vmath::Fill(nq_face,fac*normal[i][0],normal[i],1);
797  }
798  }
799  else
800  {
801  // Set up deformed normals
802  int j, k;
803 
804  int nq0 = ptsKeys[0].GetNumPoints();
805  int nq1 = ptsKeys[1].GetNumPoints();
806  int nq2 = ptsKeys[2].GetNumPoints();
807  int nqtot;
808  int nq01 =nq0*nq1;
809 
810  // number of elemental quad points
811  if (face == 0)
812  {
813  nqtot = nq01;
814  }
815  else if (face == 1)
816  {
817  nqtot = nq0*nq2;
818  }
819  else
820  {
821  nqtot = nq1*nq2;
822  }
823 
824  LibUtilities::PointsKey points0;
825  LibUtilities::PointsKey points1;
826 
827  Array<OneD,NekDouble> normals(vCoordDim*nqtot, 0.0);
828 
829  // Extract Jacobian along face and recover local derivates
830  // (dx/dr) for polynomial interpolation by multiplying m_gmat by
831  // jacobian
832  switch (face)
833  {
834  case 0:
835  {
836  for(j = 0; j < nq01; ++j)
837  {
838  normals[j] = -df[2][j]*jac[j];
839  normals[nqtot+j] = -df[5][j]*jac[j];
840  normals[2*nqtot+j] = -df[8][j]*jac[j];
841  }
842 
843  points0 = ptsKeys[0];
844  points1 = ptsKeys[1];
845  break;
846  }
847 
848  case 1:
849  {
850  for (j = 0; j < nq0; ++j)
851  {
852  for(k = 0; k < nq2; ++k)
853  {
854  int tmp = j+nq01*k;
855  normals[j+k*nq0] =
856  -df[1][tmp]*jac[tmp];
857  normals[nqtot+j+k*nq0] =
858  -df[4][tmp]*jac[tmp];
859  normals[2*nqtot+j+k*nq0] =
860  -df[7][tmp]*jac[tmp];
861  }
862  }
863 
864  points0 = ptsKeys[0];
865  points1 = ptsKeys[2];
866  break;
867  }
868 
869  case 2:
870  {
871  for (j = 0; j < nq1; ++j)
872  {
873  for(k = 0; k < nq2; ++k)
874  {
875  int tmp = nq0-1+nq0*j+nq01*k;
876  normals[j+k*nq1] =
877  (df[0][tmp]+df[1][tmp]+df[2][tmp])*
878  jac[tmp];
879  normals[nqtot+j+k*nq1] =
880  (df[3][tmp]+df[4][tmp]+df[5][tmp])*
881  jac[tmp];
882  normals[2*nqtot+j+k*nq1] =
883  (df[6][tmp]+df[7][tmp]+df[8][tmp])*
884  jac[tmp];
885  }
886  }
887 
888  points0 = ptsKeys[1];
889  points1 = ptsKeys[2];
890  break;
891  }
892 
893  case 3:
894  {
895  for (j = 0; j < nq1; ++j)
896  {
897  for(k = 0; k < nq2; ++k)
898  {
899  int tmp = j*nq0+nq01*k;
900  normals[j+k*nq1] =
901  -df[0][tmp]*jac[tmp];
902  normals[nqtot+j+k*nq1] =
903  -df[3][tmp]*jac[tmp];
904  normals[2*nqtot+j+k*nq1] =
905  -df[6][tmp]*jac[tmp];
906  }
907  }
908 
909  points0 = ptsKeys[1];
910  points1 = ptsKeys[2];
911  break;
912  }
913 
914  default:
915  ASSERTL0(false,"face is out of range (face < 3)");
916  }
917 
918  Array<OneD,NekDouble> work (nq_face, 0.0);
919  // Interpolate Jacobian and invert
920  LibUtilities::Interp2D(points0, points1, jac,
921  tobasis0.GetPointsKey(),
922  tobasis1.GetPointsKey(),
923  work);
924  Vmath::Sdiv(nq_face, 1.0, &work[0], 1, &work[0], 1);
925 
926  // Interpolate normal and multiply by inverse Jacobian.
927  for(i = 0; i < vCoordDim; ++i)
928  {
929  LibUtilities::Interp2D(points0, points1,
930  &normals[i*nqtot],
931  tobasis0.GetPointsKey(),
932  tobasis1.GetPointsKey(),
933  &normal[i][0]);
934  Vmath::Vmul(nq_face,work,1,normal[i],1,normal[i],1);
935  }
936 
937  // Normalise to obtain unit normals.
938  Vmath::Zero(nq_face,work,1);
939  for(i = 0; i < GetCoordim(); ++i)
940  {
941  Vmath::Vvtvp(nq_face,normal[i],1,normal[i],1,work,1,work,1);
942  }
943 
944  Vmath::Vsqrt(nq_face,work,1,work,1);
945  Vmath::Sdiv (nq_face,1.0,work,1,work,1);
946 
947  for(i = 0; i < GetCoordim(); ++i)
948  {
949  Vmath::Vmul(nq_face,normal[i],1,work,1,normal[i],1);
950  }
951  }
952  }
953 
954  //-----------------------------
955  // Operator creation functions
956  //-----------------------------
958  const Array<OneD, const NekDouble> &inarray,
959  Array<OneD,NekDouble> &outarray,
960  const StdRegions::StdMatrixKey &mkey)
961  {
962  TetExp::v_HelmholtzMatrixOp_MatFree(inarray,outarray,mkey);
963  }
964 
965 
967  const Array<OneD, const NekDouble> &inarray,
968  Array<OneD,NekDouble> &outarray,
969  const StdRegions::StdMatrixKey &mkey)
970  {
971  TetExp::v_LaplacianMatrixOp_MatFree(inarray,outarray,mkey);
972  }
973 
975  const int k1,
976  const int k2,
977  const Array<OneD, const NekDouble> &inarray,
978  Array<OneD,NekDouble> &outarray,
979  const StdRegions::StdMatrixKey &mkey)
980  {
981  StdExpansion::LaplacianMatrixOp_MatFree(k1,k2,inarray,outarray,
982  mkey);
983  }
984 
985 
986  //-----------------------------
987  // Matrix creation functions
988  //-----------------------------
990  const StdRegions::StdMatrixKey &mkey)
991  {
992  DNekMatSharedPtr returnval;
993 
994  switch(mkey.GetMatrixType())
995  {
1003  returnval = Expansion3D::v_GenMatrix(mkey);
1004  break;
1005  default:
1006  returnval = StdTetExp::v_GenMatrix(mkey);
1007  }
1008 
1009  return returnval;
1010  }
1011 
1012 
1014  {
1015  DNekScalMatSharedPtr returnval;
1017 
1018  ASSERTL2(m_metricinfo->GetGtype() != SpatialDomains::eNoGeomType,"Geometric information is not set up");
1019 
1020  switch(mkey.GetMatrixType())
1021  {
1022  case StdRegions::eMass:
1023  {
1024  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed ||
1025  mkey.GetNVarCoeff())
1026  {
1027  NekDouble one = 1.0;
1028  DNekMatSharedPtr mat = GenMatrix(mkey);
1029  returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(one,mat);
1030  }
1031  else
1032  {
1033  NekDouble jac = (m_metricinfo->GetJac(ptsKeys))[0];
1034  DNekMatSharedPtr mat = GetStdMatrix(mkey);
1035  returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(jac,mat);
1036  }
1037  }
1038  break;
1039  case StdRegions::eInvMass:
1040  {
1041  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1042  {
1043  NekDouble one = 1.0;
1045  *this);
1046  DNekMatSharedPtr mat = GenMatrix(masskey);
1047  mat->Invert();
1048  returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(one,mat);
1049  }
1050  else
1051  {
1052  NekDouble fac = 1.0/(m_metricinfo->GetJac(ptsKeys))[0];
1053  DNekMatSharedPtr mat = GetStdMatrix(mkey);
1054  returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(fac,mat);
1055  }
1056  }
1057  break;
1061  {
1062  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed ||
1063  mkey.GetNVarCoeff())
1064  {
1065  NekDouble one = 1.0;
1066  DNekMatSharedPtr mat = GenMatrix(mkey);
1067 
1068  returnval = MemoryManager<DNekScalMat>
1069  ::AllocateSharedPtr(one,mat);
1070  }
1071  else
1072  {
1073  NekDouble jac = (m_metricinfo->GetJac(ptsKeys))[0];
1074  Array<TwoD, const NekDouble> df
1075  = m_metricinfo->GetDerivFactors(ptsKeys);
1076  int dir = 0;
1077 
1078  switch(mkey.GetMatrixType())
1079  {
1081  dir = 0;
1082  break;
1084  dir = 1;
1085  break;
1087  dir = 2;
1088  break;
1089  default:
1090  break;
1091  }
1092 
1094  mkey.GetShapeType(), *this);
1096  mkey.GetShapeType(), *this);
1098  mkey.GetShapeType(), *this);
1099 
1100  DNekMat &deriv0 = *GetStdMatrix(deriv0key);
1101  DNekMat &deriv1 = *GetStdMatrix(deriv1key);
1102  DNekMat &deriv2 = *GetStdMatrix(deriv2key);
1103 
1104  int rows = deriv0.GetRows();
1105  int cols = deriv1.GetColumns();
1106 
1108  ::AllocateSharedPtr(rows,cols);
1109  (*WeakDeriv) = df[3*dir][0]*deriv0
1110  + df[3*dir+1][0]*deriv1
1111  + df[3*dir+2][0]*deriv2;
1112 
1113  returnval = MemoryManager<DNekScalMat>
1114  ::AllocateSharedPtr(jac,WeakDeriv);
1115  }
1116  }
1117  break;
1119  {
1120  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed ||
1122  {
1123  NekDouble one = 1.0;
1124  DNekMatSharedPtr mat = GenMatrix(mkey);
1125 
1126  returnval = MemoryManager<DNekScalMat>
1127  ::AllocateSharedPtr(one,mat);
1128  }
1129  else
1130  {
1132  mkey.GetShapeType(), *this);
1134  mkey.GetShapeType(), *this);
1136  mkey.GetShapeType(), *this);
1138  mkey.GetShapeType(), *this);
1140  mkey.GetShapeType(), *this);
1142  mkey.GetShapeType(), *this);
1143 
1144  DNekMat &lap00 = *GetStdMatrix(lap00key);
1145  DNekMat &lap01 = *GetStdMatrix(lap01key);
1146  DNekMat &lap02 = *GetStdMatrix(lap02key);
1147  DNekMat &lap11 = *GetStdMatrix(lap11key);
1148  DNekMat &lap12 = *GetStdMatrix(lap12key);
1149  DNekMat &lap22 = *GetStdMatrix(lap22key);
1150 
1151  NekDouble jac = (m_metricinfo->GetJac(ptsKeys))[0];
1152  Array<TwoD, const NekDouble> gmat
1153  = m_metricinfo->GetGmat(ptsKeys);
1154 
1155  int rows = lap00.GetRows();
1156  int cols = lap00.GetColumns();
1157 
1159  ::AllocateSharedPtr(rows,cols);
1160 
1161  (*lap) = gmat[0][0]*lap00
1162  + gmat[4][0]*lap11
1163  + gmat[8][0]*lap22
1164  + gmat[3][0]*(lap01 + Transpose(lap01))
1165  + gmat[6][0]*(lap02 + Transpose(lap02))
1166  + gmat[7][0]*(lap12 + Transpose(lap12));
1167 
1168  returnval = MemoryManager<DNekScalMat>
1169  ::AllocateSharedPtr(jac,lap);
1170  }
1171  }
1172  break;
1174  {
1176  MatrixKey masskey(StdRegions::eMass, mkey.GetShapeType(), *this);
1177  DNekScalMat &MassMat = *(this->m_matrixManager[masskey]);
1178  MatrixKey lapkey(StdRegions::eLaplacian, mkey.GetShapeType(), *this, mkey.GetConstFactors(), mkey.GetVarCoeffs());
1179  DNekScalMat &LapMat = *(this->m_matrixManager[lapkey]);
1180 
1181  int rows = LapMat.GetRows();
1182  int cols = LapMat.GetColumns();
1183 
1185 
1186  NekDouble one = 1.0;
1187  (*helm) = LapMat + factor*MassMat;
1188 
1189  returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(one, helm);
1190  }
1191  break;
1193  {
1194  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1195  {
1196  NekDouble one = 1.0;
1197  DNekMatSharedPtr mat = GenMatrix(mkey);
1198  returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(one,mat);
1199  }
1200  else
1201  {
1202  NekDouble jac = (m_metricinfo->GetJac(ptsKeys))[0];
1203  DNekMatSharedPtr mat = GetStdMatrix(mkey);
1204  returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(jac,mat);
1205  }
1206  }
1207  break;
1215  {
1216  NekDouble one = 1.0;
1217 
1218  DNekMatSharedPtr mat = GenMatrix(mkey);
1219  returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(one,mat);
1220  }
1221  break;
1223  {
1224  NekDouble one = 1.0;
1225 
1227  DNekMatSharedPtr mat = GenMatrix(hkey);
1228 
1229  mat->Invert();
1230  returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(one,mat);
1231  }
1232  break;
1234  {
1235  NekDouble one = 1.0;
1236  MatrixKey helmkey(StdRegions::eHelmholtz, mkey.GetShapeType(), *this, mkey.GetConstFactors(), mkey.GetVarCoeffs());
1237  DNekScalBlkMatSharedPtr helmStatCond = GetLocStaticCondMatrix(helmkey);
1238  DNekScalMatSharedPtr A =helmStatCond->GetBlock(0,0);
1240 
1242  }
1243  break;
1245  {
1246  NekDouble one = 1.0;
1247  MatrixKey masskey(StdRegions::eMass, mkey.GetShapeType(), *this);
1248  DNekScalBlkMatSharedPtr massStatCond = GetLocStaticCondMatrix(masskey);
1249  DNekScalMatSharedPtr A =massStatCond->GetBlock(0,0);
1251 
1253  }
1254  break;
1255  case StdRegions::ePreconR:
1256  {
1257  NekDouble one = 1.0;
1258  MatrixKey helmkey(StdRegions::eHelmholtz, mkey.GetShapeType(), *this,mkey.GetConstFactors(), mkey.GetVarCoeffs());
1259  DNekScalBlkMatSharedPtr helmStatCond = GetLocStaticCondMatrix(helmkey);
1260  DNekScalMatSharedPtr A =helmStatCond->GetBlock(0,0);
1261 
1262  DNekScalMatSharedPtr Atmp;
1264 
1266  }
1267  break;
1268  case StdRegions::ePreconRT:
1269  {
1270  NekDouble one = 1.0;
1271  MatrixKey helmkey(StdRegions::eHelmholtz, mkey.GetShapeType(), *this,mkey.GetConstFactors(), mkey.GetVarCoeffs());
1272  DNekScalBlkMatSharedPtr helmStatCond = GetLocStaticCondMatrix(helmkey);
1273  DNekScalMatSharedPtr A =helmStatCond->GetBlock(0,0);
1274 
1275  DNekScalMatSharedPtr Atmp;
1277 
1279  }
1280  break;
1282  {
1283  NekDouble one = 1.0;
1284  MatrixKey masskey(StdRegions::eMass, mkey.GetShapeType(), *this);
1285  DNekScalBlkMatSharedPtr StatCond = GetLocStaticCondMatrix(masskey);
1286  DNekScalMatSharedPtr A =StatCond->GetBlock(0,0);
1287 
1288  DNekScalMatSharedPtr Atmp;
1290 
1292  }
1293  break;
1295  {
1296  NekDouble one = 1.0;
1297  MatrixKey masskey(StdRegions::eMass, mkey.GetShapeType(), *this);
1298  DNekScalBlkMatSharedPtr massStatCond = GetLocStaticCondMatrix(masskey);
1299  DNekScalMatSharedPtr A =massStatCond->GetBlock(0,0);
1300 
1301  DNekScalMatSharedPtr Atmp;
1303 
1305  }
1306  break;
1307  default:
1308  {
1309  //ASSERTL0(false, "Missing definition for " + (*StdRegions::MatrixTypeMap[mkey.GetMatrixType()]));
1310  NekDouble one = 1.0;
1311  DNekMatSharedPtr mat = GenMatrix(mkey);
1312 
1313  returnval = MemoryManager<DNekScalMat>::AllocateSharedPtr(one,mat);
1314  }
1315  break;
1316  }
1317 
1318  return returnval;
1319  }
1320 
1321 
1323  const MatrixKey &mkey)
1324  {
1325  DNekScalBlkMatSharedPtr returnval;
1326 
1327  ASSERTL2(m_metricinfo->GetGtype() != SpatialDomains::eNoGeomType,"Geometric information is not set up");
1328 
1329  // set up block matrix system
1330  unsigned int nbdry = NumBndryCoeffs();
1331  unsigned int nint = (unsigned int)(m_ncoeffs - nbdry);
1332  unsigned int exp_size[] = {nbdry, nint};
1333  unsigned int nblks = 2;
1334  returnval = MemoryManager<DNekScalBlkMat>::AllocateSharedPtr(nblks, nblks, exp_size, exp_size);
1335 
1336  NekDouble factor = 1.0;
1337  MatrixStorage AMatStorage = eFULL;
1338 
1339  switch(mkey.GetMatrixType())
1340  {
1342  case StdRegions::eHelmholtz: // special case since Helmholtz not defined in StdRegions
1343  // use Deformed case for both regular and deformed geometries
1344  factor = 1.0;
1345  goto UseLocRegionsMatrix;
1346  break;
1347  default:
1348  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed ||
1349  mkey.GetNVarCoeff())
1350  {
1351  factor = 1.0;
1352  goto UseLocRegionsMatrix;
1353  }
1354  else
1355  {
1356  DNekScalMatSharedPtr mat = GetLocMatrix(mkey);
1357  factor = mat->Scale();
1358  goto UseStdRegionsMatrix;
1359  }
1360  break;
1361  UseStdRegionsMatrix:
1362  {
1363  NekDouble invfactor = 1.0/factor;
1364  NekDouble one = 1.0;
1366  DNekScalMatSharedPtr Atmp;
1367  DNekMatSharedPtr Asubmat;
1368 
1369  //TODO: check below
1370  returnval->SetBlock(0,0,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(factor,Asubmat = mat->GetBlock(0,0)));
1371  returnval->SetBlock(0,1,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(one,Asubmat = mat->GetBlock(0,1)));
1372  returnval->SetBlock(1,0,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(factor,Asubmat = mat->GetBlock(1,0)));
1373  returnval->SetBlock(1,1,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(invfactor,Asubmat = mat->GetBlock(1,1)));
1374  }
1375  break;
1376  UseLocRegionsMatrix:
1377  {
1378  int i,j;
1379  NekDouble invfactor = 1.0/factor;
1380  NekDouble one = 1.0;
1381  DNekScalMat &mat = *GetLocMatrix(mkey);
1382  DNekMatSharedPtr A = MemoryManager<DNekMat>::AllocateSharedPtr(nbdry,nbdry,AMatStorage);
1386 
1387  Array<OneD,unsigned int> bmap(nbdry);
1388  Array<OneD,unsigned int> imap(nint);
1389  GetBoundaryMap(bmap);
1390  GetInteriorMap(imap);
1391 
1392  for(i = 0; i < nbdry; ++i)
1393  {
1394  for(j = 0; j < nbdry; ++j)
1395  {
1396  (*A)(i,j) = mat(bmap[i],bmap[j]);
1397  }
1398 
1399  for(j = 0; j < nint; ++j)
1400  {
1401  (*B)(i,j) = mat(bmap[i],imap[j]);
1402  }
1403  }
1404 
1405  for(i = 0; i < nint; ++i)
1406  {
1407  for(j = 0; j < nbdry; ++j)
1408  {
1409  (*C)(i,j) = mat(imap[i],bmap[j]);
1410  }
1411 
1412  for(j = 0; j < nint; ++j)
1413  {
1414  (*D)(i,j) = mat(imap[i],imap[j]);
1415  }
1416  }
1417 
1418  // Calculate static condensed system
1419  if(nint)
1420  {
1421  D->Invert();
1422  (*B) = (*B)*(*D);
1423  (*A) = (*A) - (*B)*(*C);
1424  }
1425 
1426  DNekScalMatSharedPtr Atmp;
1427 
1428  returnval->SetBlock(0,0,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(factor,A));
1429  returnval->SetBlock(0,1,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(one,B));
1430  returnval->SetBlock(1,0,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(factor,C));
1431  returnval->SetBlock(1,1,Atmp = MemoryManager<DNekScalMat>::AllocateSharedPtr(invfactor,D));
1432 
1433  }
1434  break;
1435  }
1436  return returnval;
1437  }
1438 
1439 
1441  const StdRegions::StdMatrixKey &mkey)
1442  {
1443  LibUtilities::BasisKey bkey0 = m_base[0]->GetBasisKey();
1444  LibUtilities::BasisKey bkey1 = m_base[1]->GetBasisKey();
1445  LibUtilities::BasisKey bkey2 = m_base[2]->GetBasisKey();
1447 
1448  return tmp->GetStdMatrix(mkey);
1449  }
1450 
1452  {
1453  return m_matrixManager[mkey];
1454  }
1455 
1457  {
1458  return m_staticCondMatrixManager[mkey];
1459  }
1460 
1462  {
1464  }
1465 
1467  const Array<OneD, const NekDouble> &inarray,
1468  Array<OneD,NekDouble> &outarray,
1469  const StdRegions::StdMatrixKey &mkey)
1470  {
1471  DNekScalMatSharedPtr mat = GetLocMatrix(mkey);
1472 
1473  if(inarray.get() == outarray.get())
1474  {
1475  Array<OneD,NekDouble> tmp(m_ncoeffs);
1476  Vmath::Vcopy(m_ncoeffs,inarray.get(),1,tmp.get(),1);
1477 
1478  Blas::Dgemv('N',m_ncoeffs,m_ncoeffs,mat->Scale(),(mat->GetOwnedMatrix())->GetPtr().get(),
1479  m_ncoeffs, tmp.get(), 1, 0.0, outarray.get(), 1);
1480  }
1481  else
1482  {
1483  Blas::Dgemv('N',m_ncoeffs,m_ncoeffs,mat->Scale(),(mat->GetOwnedMatrix())->GetPtr().get(),
1484  m_ncoeffs, inarray.get(), 1, 0.0, outarray.get(), 1);
1485  }
1486  }
1487 
1488 
1490  const Array<OneD, const NekDouble> &inarray,
1491  Array<OneD, NekDouble> &outarray,
1492  Array<OneD, NekDouble> &wsp)
1493  {
1494  // This implementation is only valid when there are no
1495  // coefficients associated to the Laplacian operator
1496  if (m_metrics.count(MetricLaplacian00) == 0)
1497  {
1499  }
1500 
1501  int nquad0 = m_base[0]->GetNumPoints();
1502  int nquad1 = m_base[1]->GetNumPoints();
1503  int nquad2 = m_base[2]->GetNumPoints();
1504  int nqtot = nquad0*nquad1*nquad2;
1505 
1506  ASSERTL1(wsp.num_elements() >= 6*nqtot,
1507  "Insufficient workspace size.");
1508  ASSERTL1(m_ncoeffs <= nqtot,
1509  "Workspace not set up for ncoeffs > nqtot");
1510 
1511  const Array<OneD, const NekDouble>& base0 = m_base[0]->GetBdata();
1512  const Array<OneD, const NekDouble>& base1 = m_base[1]->GetBdata();
1513  const Array<OneD, const NekDouble>& base2 = m_base[2]->GetBdata();
1514  const Array<OneD, const NekDouble>& dbase0 = m_base[0]->GetDbdata();
1515  const Array<OneD, const NekDouble>& dbase1 = m_base[1]->GetDbdata();
1516  const Array<OneD, const NekDouble>& dbase2 = m_base[2]->GetDbdata();
1517  const Array<OneD, const NekDouble>& metric00 = m_metrics[MetricLaplacian00];
1518  const Array<OneD, const NekDouble>& metric01 = m_metrics[MetricLaplacian01];
1519  const Array<OneD, const NekDouble>& metric02 = m_metrics[MetricLaplacian02];
1520  const Array<OneD, const NekDouble>& metric11 = m_metrics[MetricLaplacian11];
1521  const Array<OneD, const NekDouble>& metric12 = m_metrics[MetricLaplacian12];
1522  const Array<OneD, const NekDouble>& metric22 = m_metrics[MetricLaplacian22];
1523 
1524  // Allocate temporary storage
1525  Array<OneD,NekDouble> wsp0 (2*nqtot, wsp);
1526  Array<OneD,NekDouble> wsp1 ( nqtot, wsp+1*nqtot);
1527  Array<OneD,NekDouble> wsp2 ( nqtot, wsp+2*nqtot);
1528  Array<OneD,NekDouble> wsp3 ( nqtot, wsp+3*nqtot);
1529  Array<OneD,NekDouble> wsp4 ( nqtot, wsp+4*nqtot);
1530  Array<OneD,NekDouble> wsp5 ( nqtot, wsp+5*nqtot);
1531 
1532  // LAPLACIAN MATRIX OPERATION
1533  // wsp1 = du_dxi1 = D_xi1 * inarray = D_xi1 * u
1534  // wsp2 = du_dxi2 = D_xi2 * inarray = D_xi2 * u
1535  // wsp2 = du_dxi3 = D_xi3 * inarray = D_xi3 * u
1536  StdExpansion3D::PhysTensorDeriv(inarray,wsp0,wsp1,wsp2);
1537 
1538  // wsp0 = k = g0 * wsp1 + g1 * wsp2 = g0 * du_dxi1 + g1 * du_dxi2
1539  // wsp2 = l = g1 * wsp1 + g2 * wsp2 = g0 * du_dxi1 + g1 * du_dxi2
1540  // where g0, g1 and g2 are the metric terms set up in the GeomFactors class
1541  // especially for this purpose
1542  Vmath::Vvtvvtp(nqtot,&metric00[0],1,&wsp0[0],1,&metric01[0],1,&wsp1[0],1,&wsp3[0],1);
1543  Vmath::Vvtvp (nqtot,&metric02[0],1,&wsp2[0],1,&wsp3[0],1,&wsp3[0],1);
1544  Vmath::Vvtvvtp(nqtot,&metric01[0],1,&wsp0[0],1,&metric11[0],1,&wsp1[0],1,&wsp4[0],1);
1545  Vmath::Vvtvp (nqtot,&metric12[0],1,&wsp2[0],1,&wsp4[0],1,&wsp4[0],1);
1546  Vmath::Vvtvvtp(nqtot,&metric02[0],1,&wsp0[0],1,&metric12[0],1,&wsp1[0],1,&wsp5[0],1);
1547  Vmath::Vvtvp (nqtot,&metric22[0],1,&wsp2[0],1,&wsp5[0],1,&wsp5[0],1);
1548 
1549  // outarray = m = (D_xi1 * B)^T * k
1550  // wsp1 = n = (D_xi2 * B)^T * l
1551  IProductWRTBase_SumFacKernel(dbase0,base1,base2,wsp3,outarray,wsp0,false,true,true);
1552  IProductWRTBase_SumFacKernel(base0,dbase1,base2,wsp4,wsp2, wsp0,true,false,true);
1553  Vmath::Vadd(m_ncoeffs,wsp2.get(),1,outarray.get(),1,outarray.get(),1);
1554  IProductWRTBase_SumFacKernel(base0,base1,dbase2,wsp5,wsp2, wsp0,true,true,false);
1555  Vmath::Vadd(m_ncoeffs,wsp2.get(),1,outarray.get(),1,outarray.get(),1);
1556  }
1557 
1558 
1560  {
1561  if (m_metrics.count(MetricQuadrature) == 0)
1562  {
1564  }
1565 
1566  int i, j;
1567  const unsigned int nqtot = GetTotPoints();
1568  const unsigned int dim = 3;
1572  };
1573 
1574  for (unsigned int i = 0; i < dim; ++i)
1575  {
1576  for (unsigned int j = i; j < dim; ++j)
1577  {
1578  m_metrics[m[i][j]] = Array<OneD, NekDouble>(nqtot);
1579  }
1580  }
1581 
1582  // Define shorthand synonyms for m_metrics storage
1583  Array<OneD,NekDouble> g0 (m_metrics[m[0][0]]);
1584  Array<OneD,NekDouble> g1 (m_metrics[m[1][1]]);
1585  Array<OneD,NekDouble> g2 (m_metrics[m[2][2]]);
1586  Array<OneD,NekDouble> g3 (m_metrics[m[0][1]]);
1587  Array<OneD,NekDouble> g4 (m_metrics[m[0][2]]);
1588  Array<OneD,NekDouble> g5 (m_metrics[m[1][2]]);
1589 
1590  // Allocate temporary storage
1591  Array<OneD,NekDouble> alloc(7*nqtot,0.0);
1592  Array<OneD,NekDouble> h0 (alloc); // h0
1593  Array<OneD,NekDouble> h1 (alloc+ 1*nqtot);// h1
1594  Array<OneD,NekDouble> h2 (alloc+ 2*nqtot);// h2
1595  Array<OneD,NekDouble> h3 (alloc+ 3*nqtot);// h3
1596  Array<OneD,NekDouble> wsp4 (alloc+ 4*nqtot);// wsp4
1597  Array<OneD,NekDouble> wsp5 (alloc+ 5*nqtot);// wsp5
1598  Array<OneD,NekDouble> wsp6 (alloc+ 6*nqtot);// wsp6
1599  // Reuse some of the storage as workspace
1600  Array<OneD,NekDouble> wsp7 (alloc); // wsp7
1601  Array<OneD,NekDouble> wsp8 (alloc+ 1*nqtot);// wsp8
1602  Array<OneD,NekDouble> wsp9 (alloc+ 2*nqtot);// wsp9
1603 
1604  const Array<TwoD, const NekDouble>& df =
1605  m_metricinfo->GetDerivFactors(GetPointsKeys());
1606  const Array<OneD, const NekDouble>& z0 = m_base[0]->GetZ();
1607  const Array<OneD, const NekDouble>& z1 = m_base[1]->GetZ();
1608  const Array<OneD, const NekDouble>& z2 = m_base[2]->GetZ();
1609  const unsigned int nquad0 = m_base[0]->GetNumPoints();
1610  const unsigned int nquad1 = m_base[1]->GetNumPoints();
1611  const unsigned int nquad2 = m_base[2]->GetNumPoints();
1612 
1613  for(j = 0; j < nquad2; ++j)
1614  {
1615  for(i = 0; i < nquad1; ++i)
1616  {
1617  Vmath::Fill(nquad0, 4.0/(1.0-z1[i])/(1.0-z2[j]), &h0[0]+i*nquad0 + j*nquad0*nquad1,1);
1618  Vmath::Fill(nquad0, 2.0/(1.0-z1[i])/(1.0-z2[j]), &h1[0]+i*nquad0 + j*nquad0*nquad1,1);
1619  Vmath::Fill(nquad0, 2.0/(1.0-z2[j]), &h2[0]+i*nquad0 + j*nquad0*nquad1,1);
1620  Vmath::Fill(nquad0, (1.0+z1[i])/(1.0-z2[j]), &h3[0]+i*nquad0 + j*nquad0*nquad1,1);
1621  }
1622  }
1623  for(i = 0; i < nquad0; i++)
1624  {
1625  Blas::Dscal(nquad1*nquad2, 1+z0[i], &h1[0]+i, nquad0);
1626  }
1627 
1628  // Step 3. Construct combined metric terms for physical space to
1629  // collapsed coordinate system.
1630  // Order of construction optimised to minimise temporary storage
1631  if(m_metricinfo->GetGtype() == SpatialDomains::eDeformed)
1632  {
1633  // wsp4
1634  Vmath::Vadd(nqtot, &df[1][0], 1, &df[2][0], 1, &wsp4[0], 1);
1635  Vmath::Vvtvvtp(nqtot, &df[0][0], 1, &h0[0], 1, &wsp4[0], 1, &h1[0], 1, &wsp4[0], 1);
1636  // wsp5
1637  Vmath::Vadd(nqtot, &df[4][0], 1, &df[5][0], 1, &wsp5[0], 1);
1638  Vmath::Vvtvvtp(nqtot, &df[3][0], 1, &h0[0], 1, &wsp5[0], 1, &h1[0], 1, &wsp5[0], 1);
1639  // wsp6
1640  Vmath::Vadd(nqtot, &df[7][0], 1, &df[8][0], 1, &wsp6[0], 1);
1641  Vmath::Vvtvvtp(nqtot, &df[6][0], 1, &h0[0], 1, &wsp6[0], 1, &h1[0], 1, &wsp6[0], 1);
1642 
1643  // g0
1644  Vmath::Vvtvvtp(nqtot, &wsp4[0], 1, &wsp4[0], 1, &wsp5[0], 1, &wsp5[0], 1, &g0[0], 1);
1645  Vmath::Vvtvp (nqtot, &wsp6[0], 1, &wsp6[0], 1, &g0[0], 1, &g0[0], 1);
1646 
1647  // g4
1648  Vmath::Vvtvvtp(nqtot, &df[2][0], 1, &wsp4[0], 1, &df[5][0], 1, &wsp5[0], 1, &g4[0], 1);
1649  Vmath::Vvtvp (nqtot, &df[8][0], 1, &wsp6[0], 1, &g4[0], 1, &g4[0], 1);
1650 
1651  // overwrite h0, h1, h2
1652  // wsp7 (h2f1 + h3f2)
1653  Vmath::Vvtvvtp(nqtot, &df[1][0], 1, &h2[0], 1, &df[2][0], 1, &h3[0], 1, &wsp7[0], 1);
1654  // wsp8 (h2f4 + h3f5)
1655  Vmath::Vvtvvtp(nqtot, &df[4][0], 1, &h2[0], 1, &df[5][0], 1, &h3[0], 1, &wsp8[0], 1);
1656  // wsp9 (h2f7 + h3f8)
1657  Vmath::Vvtvvtp(nqtot, &df[7][0], 1, &h2[0], 1, &df[8][0], 1, &h3[0], 1, &wsp9[0], 1);
1658 
1659  // g3
1660  Vmath::Vvtvvtp(nqtot, &wsp4[0], 1, &wsp7[0], 1, &wsp5[0], 1, &wsp8[0], 1, &g3[0], 1);
1661  Vmath::Vvtvp (nqtot, &wsp6[0], 1, &wsp9[0], 1, &g3[0], 1, &g3[0], 1);
1662 
1663  // overwrite wsp4, wsp5, wsp6
1664  // g1
1665  Vmath::Vvtvvtp(nqtot, &wsp7[0], 1, &wsp7[0], 1, &wsp8[0], 1, &wsp8[0], 1, &g1[0], 1);
1666  Vmath::Vvtvp (nqtot, &wsp9[0], 1, &wsp9[0], 1, &g1[0], 1, &g1[0], 1);
1667 
1668  // g5
1669  Vmath::Vvtvvtp(nqtot, &df[2][0], 1, &wsp7[0], 1, &df[5][0], 1, &wsp8[0], 1, &g5[0], 1);
1670  Vmath::Vvtvp (nqtot, &df[8][0], 1, &wsp9[0], 1, &g5[0], 1, &g5[0], 1);
1671 
1672  // g2
1673  Vmath::Vvtvvtp(nqtot, &df[2][0], 1, &df[2][0], 1, &df[5][0], 1, &df[5][0], 1, &g2[0], 1);
1674  Vmath::Vvtvp (nqtot, &df[8][0], 1, &df[8][0], 1, &g2[0], 1, &g2[0], 1);
1675  }
1676  else
1677  {
1678  // wsp4
1679  Vmath::Svtsvtp(nqtot, df[0][0], &h0[0], 1, df[1][0] + df[2][0], &h1[0], 1, &wsp4[0], 1);
1680  // wsp5
1681  Vmath::Svtsvtp(nqtot, df[3][0], &h0[0], 1, df[4][0] + df[5][0], &h1[0], 1, &wsp5[0], 1);
1682  // wsp6
1683  Vmath::Svtsvtp(nqtot, df[6][0], &h0[0], 1, df[7][0] + df[8][0], &h1[0], 1, &wsp6[0], 1);
1684 
1685  // g0
1686  Vmath::Vvtvvtp(nqtot, &wsp4[0], 1, &wsp4[0], 1, &wsp5[0], 1, &wsp5[0], 1, &g0[0], 1);
1687  Vmath::Vvtvp (nqtot, &wsp6[0], 1, &wsp6[0], 1, &g0[0], 1, &g0[0], 1);
1688 
1689  // g4
1690  Vmath::Svtsvtp(nqtot, df[2][0], &wsp4[0], 1, df[5][0], &wsp5[0], 1, &g4[0], 1);
1691  Vmath::Svtvp (nqtot, df[8][0], &wsp6[0], 1, &g4[0], 1, &g4[0], 1);
1692 
1693  // overwrite h0, h1, h2
1694  // wsp7 (h2f1 + h3f2)
1695  Vmath::Svtsvtp(nqtot, df[1][0], &h2[0], 1, df[2][0], &h3[0], 1, &wsp7[0], 1);
1696  // wsp8 (h2f4 + h3f5)
1697  Vmath::Svtsvtp(nqtot, df[4][0], &h2[0], 1, df[5][0], &h3[0], 1, &wsp8[0], 1);
1698  // wsp9 (h2f7 + h3f8)
1699  Vmath::Svtsvtp(nqtot, df[7][0], &h2[0], 1, df[8][0], &h3[0], 1, &wsp9[0], 1);
1700 
1701  // g3
1702  Vmath::Vvtvvtp(nqtot, &wsp4[0], 1, &wsp7[0], 1, &wsp5[0], 1, &wsp8[0], 1, &g3[0], 1);
1703  Vmath::Vvtvp (nqtot, &wsp6[0], 1, &wsp9[0], 1, &g3[0], 1, &g3[0], 1);
1704 
1705  // overwrite wsp4, wsp5, wsp6
1706  // g1
1707  Vmath::Vvtvvtp(nqtot, &wsp7[0], 1, &wsp7[0], 1, &wsp8[0], 1, &wsp8[0], 1, &g1[0], 1);
1708  Vmath::Vvtvp (nqtot, &wsp9[0], 1, &wsp9[0], 1, &g1[0], 1, &g1[0], 1);
1709 
1710  // g5
1711  Vmath::Svtsvtp(nqtot, df[2][0], &wsp7[0], 1, df[5][0], &wsp8[0], 1, &g5[0], 1);
1712  Vmath::Svtvp (nqtot, df[8][0], &wsp9[0], 1, &g5[0], 1, &g5[0], 1);
1713 
1714  // g2
1715  Vmath::Fill(nqtot, df[2][0]*df[2][0] + df[5][0]*df[5][0] + df[8][0]*df[8][0], &g2[0], 1);
1716  }
1717 
1718  for (unsigned int i = 0; i < dim; ++i)
1719  {
1720  for (unsigned int j = i; j < dim; ++j)
1721  {
1723  m_metrics[m[i][j]]);
1724 
1725  }
1726  }
1727 
1728 
1729  }
1730  }//end of namespace
1731 }//end of namespace