Nektar++
Public Member Functions | Protected Attributes | Private Member Functions | List of all members
Nektar::Collections::Helmholtz_IterPerExp Class Reference

Helmholtz operator using LocalRegions implementation. More...

Inheritance diagram for Nektar::Collections::Helmholtz_IterPerExp:
[legend]

Public Member Functions

 ~Helmholtz_IterPerExp () final
 
void operator() (const Array< OneD, const NekDouble > &input, Array< OneD, NekDouble > &output, Array< OneD, NekDouble > &output1, Array< OneD, NekDouble > &output2, Array< OneD, NekDouble > &wsp) final
 Perform operation. More...
 
void operator() (int dir, const Array< OneD, const NekDouble > &input, Array< OneD, NekDouble > &output, Array< OneD, NekDouble > &wsp) final
 
virtual void CheckFactors (StdRegions::FactorMap factors, int coll_phys_offset)
 Check the validity of supplied constant factors. More...
 
- Public Member Functions inherited from Nektar::Collections::Operator
 Operator (std::vector< StdRegions::StdExpansionSharedPtr > pCollExp, std::shared_ptr< CoalescedGeomData > GeomData, StdRegions::FactorMap factors)
 Constructor. More...
 
virtual COLLECTIONS_EXPORT ~Operator ()
 
unsigned int GetWspSize ()
 Get the size of the required workspace. More...
 
unsigned int GetNumElmt ()
 Get expansion pointer. More...
 
StdRegions::StdExpansionSharedPtr GetExpSharedPtr ()
 Get expansion pointer. More...
 

Protected Attributes

Array< TwoD, const NekDoublem_derivFac
 
Array< OneD, const NekDoublem_jac
 
int m_dim
 
int m_coordim
 
StdRegions::FactorMap m_factors
 
NekDouble m_lambda
 
bool m_HasVarCoeffDiff
 
Array< OneD, Array< OneD, NekDouble > > m_diff
 
const StdRegions::ConstFactorType m_factorCoeffDef [3][3]
 
- Protected Attributes inherited from Nektar::Collections::Operator
bool m_isDeformed
 
StdRegions::StdExpansionSharedPtr m_stdExp
 
unsigned int m_numElmt
 
unsigned int m_nqe
 
unsigned int m_wspSize
 

Private Member Functions

 Helmholtz_IterPerExp (vector< StdRegions::StdExpansionSharedPtr > pCollExp, CoalescedGeomDataSharedPtr pGeomData, StdRegions::FactorMap factors)
 

Detailed Description

Helmholtz operator using LocalRegions implementation.

Definition at line 178 of file library/Collections/Helmholtz.cpp.

Constructor & Destructor Documentation

◆ ~Helmholtz_IterPerExp()

Nektar::Collections::Helmholtz_IterPerExp::~Helmholtz_IterPerExp ( )
inlinefinal

Definition at line 183 of file library/Collections/Helmholtz.cpp.

184  {
185  }

◆ Helmholtz_IterPerExp()

Nektar::Collections::Helmholtz_IterPerExp::Helmholtz_IterPerExp ( vector< StdRegions::StdExpansionSharedPtr pCollExp,
CoalescedGeomDataSharedPtr  pGeomData,
StdRegions::FactorMap  factors 
)
inlineprivate

Definition at line 512 of file library/Collections/Helmholtz.cpp.

516  : Operator(pCollExp, pGeomData, factors)
517  {
518  LibUtilities::PointsKeyVector PtsKey = m_stdExp->GetPointsKeys();
519  m_dim = PtsKey.size();
520  m_coordim = pCollExp[0]->GetCoordim();
521  int nqtot = m_stdExp->GetTotPoints();
522 
523  m_derivFac = pGeomData->GetDerivFactors(pCollExp);
524  m_jac = pGeomData->GetJac(pCollExp);
525  m_wspSize = (2*m_coordim+1)*nqtot;
526 
527  m_lambda = 1.0;
528  m_HasVarCoeffDiff = false;
530  this->CheckFactors(factors, 0);
531  }
virtual void CheckFactors(StdRegions::FactorMap factors, int coll_phys_offset)
Check the validity of supplied constant factors.
StdRegions::StdExpansionSharedPtr m_stdExp
Definition: Operator.h:167
Operator(std::vector< StdRegions::StdExpansionSharedPtr > pCollExp, std::shared_ptr< CoalescedGeomData > GeomData, StdRegions::FactorMap factors)
Constructor.
Definition: Operator.cpp:41
std::vector< PointsKey > PointsKeyVector
Definition: Points.h:246
static FactorMap NullFactorMap
Definition: StdRegions.hpp:319

References Nektar::StdRegions::NullFactorMap.

Member Function Documentation

◆ CheckFactors()

virtual void Nektar::Collections::Helmholtz_IterPerExp::CheckFactors ( StdRegions::FactorMap  factors,
int  coll_phys_offset 
)
inlinevirtual

Check the validity of supplied constant factors.

Parameters
factorsMap of factors
coll_phys_offsetUnused

Implements Nektar::Collections::Operator.

Definition at line 434 of file library/Collections/Helmholtz.cpp.

436  {
437  boost::ignore_unused(coll_phys_offset);
438 
439  // If match previous factors, nothing to do.
440  if (m_factors == factors)
441  {
442  return;
443  }
444 
445  m_factors = factors;
446 
447  // Check Lambda constant of Helmholtz operator
448  auto x = factors.find(StdRegions::eFactorLambda);
449  ASSERTL1(x != factors.end(),
450  "Constant factor not defined: "
451  + std::string(StdRegions::ConstFactorTypeMap
453  m_lambda = x->second;
454 
455  // If varcoeffs not supplied, nothing else to do.
456  m_HasVarCoeffDiff = false;
457  auto d = factors.find(StdRegions::eFactorCoeffD00);
458  if (d == factors.end())
459  {
460  return;
461  }
462 
463  m_diff = Array<OneD, Array<OneD, NekDouble> >(m_coordim);
464  for(int i = 0; i < m_coordim; ++i)
465  {
466  m_diff[i] = Array<OneD, NekDouble>(m_coordim, 0.0);
467  }
468 
469  for(int i = 0; i < m_coordim; ++i)
470  {
471  d = factors.find(m_factorCoeffDef[i][i]);
472  if (d != factors.end())
473  {
474  m_diff[i][i] = d->second;
475  }
476  else
477  {
478  m_diff[i][i] = 1.0;
479  }
480 
481  for(int j = i+1; j < m_coordim; ++j)
482  {
483  d = factors.find(m_factorCoeffDef[i][j]);
484  if (d != factors.end())
485  {
486  m_diff[i][j] = m_diff[j][i] = d->second;
487  }
488  }
489  }
490  m_HasVarCoeffDiff = true;
491  }
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:250
const StdRegions::ConstFactorType m_factorCoeffDef[3][3]
Array< OneD, Array< OneD, NekDouble > > m_diff
const char *const ConstFactorTypeMap[]
Definition: StdRegions.hpp:296

References ASSERTL1, Nektar::StdRegions::ConstFactorTypeMap, Nektar::StdRegions::eFactorCoeffD00, and Nektar::StdRegions::eFactorLambda.

◆ operator()() [1/2]

void Nektar::Collections::Helmholtz_IterPerExp::operator() ( const Array< OneD, const NekDouble > &  input,
Array< OneD, NekDouble > &  output0,
Array< OneD, NekDouble > &  output1,
Array< OneD, NekDouble > &  output2,
Array< OneD, NekDouble > &  wsp 
)
inlinefinalvirtual

Perform operation.

Implements Nektar::Collections::Operator.

Definition at line 187 of file library/Collections/Helmholtz.cpp.

193  {
194  boost::ignore_unused(output1,output2);
195 
196  const int nCoeffs = m_stdExp->GetNcoeffs();
197  const int nPhys = m_stdExp->GetTotPoints();
198 
199  ASSERTL1(input.size() >= m_numElmt*nCoeffs,
200  "input array size is insufficient");
201  ASSERTL1(output.size() >= m_numElmt*nCoeffs,
202  "output array size is insufficient");
203 
204  Array<OneD, NekDouble> tmpphys, t1;
205  Array<OneD, Array<OneD, NekDouble> > dtmp(3);
206  Array<OneD, Array<OneD, NekDouble> > tmp(3);
207 
208  tmpphys = wsp;
209  for(int i = 1; i < m_coordim+1; ++i)
210  {
211  dtmp[i-1] = wsp + i*nPhys;
212  tmp [i-1] = wsp + (i+m_coordim)*nPhys;
213  }
214 
215  for (int i = 0; i < m_numElmt; ++i)
216  {
217  m_stdExp->BwdTrans(input + i*nCoeffs, tmpphys);
218 
219  // local derivative
220  m_stdExp->PhysDeriv(tmpphys, dtmp[0], dtmp[1], dtmp[2]);
221 
222  // determine mass matrix term
223  if(m_isDeformed)
224  {
225  Vmath::Vmul(nPhys,m_jac+i*nPhys,1,tmpphys,1,tmpphys,1);
226  }
227  else
228  {
229  Vmath::Smul(nPhys,m_jac[i],tmpphys,1,tmpphys,1);
230  }
231 
232  m_stdExp->IProductWRTBase(tmpphys,t1 = output + i*nCoeffs);
233  Vmath::Smul(nCoeffs,m_lambda,output + i*nCoeffs,1,
234  t1 = output+i*nCoeffs,1);
235 
236  if(m_isDeformed)
237  {
238  // calculate full derivative
239  for(int j = 0; j < m_coordim; ++j)
240  {
241  Vmath::Vmul(nPhys,
242  m_derivFac[j*m_dim].origin() + i*nPhys, 1,
243  &dtmp[0][0], 1, &tmp[j][0], 1);
244 
245  for(int k = 1; k < m_dim; ++k)
246  {
247  Vmath::Vvtvp (nPhys, m_derivFac[j*m_dim+k].origin()
248  + i*nPhys, 1, &dtmp[k][0], 1,
249  &tmp[j][0], 1, &tmp[j][0], 1);
250  }
251  }
252 
254  {
255  // calculate dtmp[i] = dx/dxi sum_j diff[0][j] tmp[j]
256  // + dy/dxi sum_j diff[1][j] tmp[j]
257  // + dz/dxi sum_j diff[2][j] tmp[j]
258 
259  // First term
260  Vmath::Smul(nPhys, m_diff[0][0], &tmp[0][0], 1,
261  &tmpphys[0], 1);
262  for(int l = 1; l < m_coordim; ++l)
263  {
264  Vmath::Svtvp(nPhys, m_diff[0][l], &tmp[l][0], 1,
265  &tmpphys[0], 1, &tmpphys[0], 1);
266  }
267 
268  for(int j = 0; j < m_dim; ++j)
269  {
270  Vmath::Vmul(nPhys,
271  m_derivFac[j].origin() + i*nPhys, 1,
272  &tmpphys[0], 1, &dtmp[j][0], 1);
273  }
274 
275  // Second and third terms
276  for(int k = 1; k < m_coordim; ++k)
277  {
278  Vmath::Smul(nPhys,m_diff[k][0], &tmp[0][0], 1,
279  &tmpphys[0], 1);
280  for(int l = 1; l < m_coordim; ++l)
281  {
282  Vmath::Svtvp(nPhys, m_diff[k][l], &tmp[l][0], 1,
283  &tmpphys[0], 1, &tmpphys[0], 1);
284  }
285 
286  for(int j = 0; j < m_dim; ++j)
287  {
288  Vmath::Vvtvp (nPhys,
289  m_derivFac[j +k*m_dim].origin()
290  + i*nPhys, 1,
291  &tmpphys[0], 1,
292  &dtmp[j][0], 1, &dtmp[j][0], 1);
293  }
294  }
295  }
296  else
297  {
298  // calculate dx/dxi tmp[0] + dy/dxi tmp[1]
299  // + dz/dxi tmp[2]
300  for(int j = 0; j < m_dim; ++j)
301  {
302  Vmath::Vmul (nPhys,
303  m_derivFac[j].origin() + i*nPhys, 1,
304  &tmp[0][0], 1, &dtmp[j][0], 1);
305 
306  for(int k = 1; k < m_coordim; ++k)
307  {
308  Vmath::Vvtvp (nPhys,
309  m_derivFac[j +k*m_dim].origin()
310  + i*nPhys, 1, &tmp[k][0], 1,
311  &dtmp[j][0], 1, &dtmp[j][0], 1);
312  }
313  }
314  }
315 
316  // calculate Iproduct WRT Std Deriv
317  for(int j = 0; j < m_dim; ++j)
318  {
319 
320  // multiply by Jacobian
321  Vmath::Vmul(nPhys,m_jac+i*nPhys,1,dtmp[j],1,dtmp[j],1);
322 
323  m_stdExp->IProductWRTDerivBase(j,dtmp[j],tmp[0]);
324  Vmath::Vadd(nCoeffs,tmp[0],1,output+i*nCoeffs,1,
325  t1 = output+i*nCoeffs,1);
326  }
327  }
328  else
329  {
330  // calculate full derivative
331  for(int j = 0; j < m_coordim; ++j)
332  {
333  Vmath::Smul(nPhys, m_derivFac[j*m_dim][i],
334  &dtmp[0][0], 1, &tmp[j][0], 1);
335 
336  for(int k = 1; k < m_dim; ++k)
337  {
338  Vmath::Svtvp (nPhys, m_derivFac[j*m_dim+k][i],
339  &dtmp[k][0], 1,
340  &tmp[j][0], 1,
341  &tmp[j][0], 1);
342  }
343 
344  }
345 
347  {
348  // calculate dtmp[i] = dx/dxi sum_j diff[0][j] tmp[j]
349  // + dy/dxi sum_j diff[1][j] tmp[j]
350  // + dz/dxi sum_j diff[2][j] tmp[j]
351 
352  // First term
353  Vmath::Smul(nPhys, m_diff[0][0], &tmp[0][0], 1,
354  &tmpphys[0], 1);
355  for(int l = 1; l < m_coordim; ++l)
356  {
357  Vmath::Svtvp(nPhys, m_diff[0][l], &tmp[l][0], 1,
358  &tmpphys[0], 1, &tmpphys[0], 1);
359  }
360 
361  for(int j = 0; j < m_dim; ++j)
362  {
363  Vmath::Smul (nPhys, m_derivFac[j][i],
364  &tmpphys[0], 1, &dtmp[j][0], 1);
365  }
366 
367  // Second and third terms
368  for(int k = 1; k < m_coordim; ++k)
369  {
370  Vmath::Smul(nPhys, m_diff[k][0], &tmp[0][0], 1,
371  &tmpphys[0], 1);
372  for(int l = 1; l < m_coordim; ++l)
373  {
374  Vmath::Svtvp(nPhys, m_diff[k][l], &tmp[l][0], 1,
375  &tmpphys[0], 1, &tmpphys[0], 1);
376  }
377 
378  for(int j = 0; j < m_dim; ++j)
379  {
380  Vmath::Svtvp (nPhys, m_derivFac[j +k*m_dim][i],
381  &tmpphys[0], 1, &dtmp[j][0], 1,
382  &dtmp[j][0], 1);
383  }
384  }
385  }
386  else
387  {
388  // calculate dx/dxi tmp[0] + dy/dxi tmp[2]
389  // + dz/dxi tmp[3]
390  for(int j = 0; j < m_dim; ++j)
391  {
392  Vmath::Smul (nPhys,m_derivFac[j][i],
393  &tmp[0][0], 1, &dtmp[j][0],1);
394 
395  for(int k = 1; k < m_coordim; ++k)
396  {
397  Vmath::Svtvp (nPhys, m_derivFac[j +k*m_dim][i],
398  &tmp[k][0], 1, &dtmp[j][0], 1,
399  &dtmp[j][0], 1);
400  }
401  }
402  }
403 
404  // calculate Iproduct WRT Std Deriv
405  for(int j = 0; j < m_dim; ++j)
406  {
407  // multiply by Jacobian
408  Vmath::Smul(nPhys,m_jac[i],dtmp[j],1,dtmp[j],1);
409 
410  m_stdExp->IProductWRTDerivBase(j,dtmp[j],tmp[0]);
411  Vmath::Vadd(nCoeffs,tmp[0],1,output+i*nCoeffs,1,
412  t1 = output+i*nCoeffs,1);
413  }
414  }
415  }
416  }
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:192
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:565
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.cpp:513
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.cpp:322
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.cpp:225

References ASSERTL1, Vmath::Smul(), Vmath::Svtvp(), Vmath::Vadd(), Vmath::Vmul(), and Vmath::Vvtvp().

◆ operator()() [2/2]

void Nektar::Collections::Helmholtz_IterPerExp::operator() ( int  dir,
const Array< OneD, const NekDouble > &  input,
Array< OneD, NekDouble > &  output,
Array< OneD, NekDouble > &  wsp 
)
inlinefinalvirtual

Implements Nektar::Collections::Operator.

Definition at line 418 of file library/Collections/Helmholtz.cpp.

422  {
423  boost::ignore_unused(dir, input, output, wsp);
424  NEKERROR(ErrorUtil::efatal, "Not valid for this operator.");
425  }
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209

References NEKERROR.

Member Data Documentation

◆ m_coordim

int Nektar::Collections::Helmholtz_IterPerExp::m_coordim
protected

Definition at line 498 of file library/Collections/Helmholtz.cpp.

◆ m_derivFac

Array<TwoD, const NekDouble> Nektar::Collections::Helmholtz_IterPerExp::m_derivFac
protected

Definition at line 495 of file library/Collections/Helmholtz.cpp.

◆ m_diff

Array<OneD, Array<OneD, NekDouble> > Nektar::Collections::Helmholtz_IterPerExp::m_diff
protected

Definition at line 502 of file library/Collections/Helmholtz.cpp.

◆ m_dim

int Nektar::Collections::Helmholtz_IterPerExp::m_dim
protected

Definition at line 497 of file library/Collections/Helmholtz.cpp.

◆ m_factorCoeffDef

const StdRegions::ConstFactorType Nektar::Collections::Helmholtz_IterPerExp::m_factorCoeffDef[3][3]
protected

◆ m_factors

StdRegions::FactorMap Nektar::Collections::Helmholtz_IterPerExp::m_factors
protected

Definition at line 499 of file library/Collections/Helmholtz.cpp.

◆ m_HasVarCoeffDiff

bool Nektar::Collections::Helmholtz_IterPerExp::m_HasVarCoeffDiff
protected

Definition at line 501 of file library/Collections/Helmholtz.cpp.

◆ m_jac

Array<OneD, const NekDouble> Nektar::Collections::Helmholtz_IterPerExp::m_jac
protected

Definition at line 496 of file library/Collections/Helmholtz.cpp.

◆ m_lambda

NekDouble Nektar::Collections::Helmholtz_IterPerExp::m_lambda
protected

Definition at line 500 of file library/Collections/Helmholtz.cpp.