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

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) override final
 Perform operation. More...
 
void operator() (int dir, const Array< OneD, const NekDouble > &input, Array< OneD, NekDouble > &output, Array< OneD, NekDouble > &wsp) override final
 
virtual void CheckFactors (StdRegions::FactorMap factors, int coll_phys_offset) override
 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 void operator() (const Array< OneD, const NekDouble > &input, Array< OneD, NekDouble > &output0, Array< OneD, NekDouble > &output1, Array< OneD, NekDouble > &output2, Array< OneD, NekDouble > &wsp=NullNekDouble1DArray)=0
 Perform operation. More...
 
virtual COLLECTIONS_EXPORT void operator() (int dir, const Array< OneD, const NekDouble > &input, Array< OneD, NekDouble > &output, Array< OneD, NekDouble > &wsp=NullNekDouble1DArray)=0
 
virtual COLLECTIONS_EXPORT ~Operator ()
 
virtual COLLECTIONS_EXPORT void CheckFactors (StdRegions::FactorMap factors, int coll_phys_offset)=0
 Check the validity of the supplied factor map. More...
 
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 163 of file library/Collections/Helmholtz.cpp.

Constructor & Destructor Documentation

◆ ~Helmholtz_IterPerExp()

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

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

169 {
170 }

◆ Helmholtz_IterPerExp()

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

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

497 : Operator(pCollExp, pGeomData, factors)
498 {
499 LibUtilities::PointsKeyVector PtsKey = m_stdExp->GetPointsKeys();
500 m_dim = PtsKey.size();
501 m_coordim = pCollExp[0]->GetCoordim();
502 int nqtot = m_stdExp->GetTotPoints();
503
504 m_derivFac = pGeomData->GetDerivFactors(pCollExp);
505 m_jac = pGeomData->GetJac(pCollExp);
506 m_wspSize = (2 * m_coordim + 1) * nqtot;
507
508 m_lambda = 1.0;
509 m_HasVarCoeffDiff = false;
511 this->CheckFactors(factors, 0);
512 }
virtual void CheckFactors(StdRegions::FactorMap factors, int coll_phys_offset) override
Check the validity of supplied constant factors.
StdRegions::StdExpansionSharedPtr m_stdExp
Definition: Operator.h:165
Operator(std::vector< StdRegions::StdExpansionSharedPtr > pCollExp, std::shared_ptr< CoalescedGeomData > GeomData, StdRegions::FactorMap factors)
Constructor.
Definition: Operator.cpp:43
std::vector< PointsKey > PointsKeyVector
Definition: Points.h:236
static FactorMap NullFactorMap
Definition: StdRegions.hpp:413
StdRegions::ConstFactorMap factors

References CheckFactors(), m_coordim, m_derivFac, m_dim, m_factors, m_HasVarCoeffDiff, m_jac, m_lambda, Nektar::Collections::Operator::m_stdExp, Nektar::Collections::Operator::m_wspSize, and Nektar::StdRegions::NullFactorMap.

Member Function Documentation

◆ CheckFactors()

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

Check the validity of supplied constant factors.

Parameters
factorsMap of factors
coll_phys_offsetUnused

Implements Nektar::Collections::Operator.

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

418 {
419 boost::ignore_unused(coll_phys_offset);
420
421 // If match previous factors, nothing to do.
422 if (m_factors == factors)
423 {
424 return;
425 }
426
428
429 // Check Lambda constant of Helmholtz operator
430 auto x = factors.find(StdRegions::eFactorLambda);
431 ASSERTL1(
432 x != factors.end(),
433 "Constant factor not defined: " +
434 std::string(
436 m_lambda = x->second;
437
438 // If varcoeffs not supplied, nothing else to do.
439 m_HasVarCoeffDiff = false;
441 if (d == factors.end())
442 {
443 return;
444 }
445
446 m_diff = Array<OneD, Array<OneD, NekDouble>>(m_coordim);
447 for (int i = 0; i < m_coordim; ++i)
448 {
449 m_diff[i] = Array<OneD, NekDouble>(m_coordim, 0.0);
450 }
451
452 for (int i = 0; i < m_coordim; ++i)
453 {
454 d = factors.find(m_factorCoeffDef[i][i]);
455 if (d != factors.end())
456 {
457 m_diff[i][i] = d->second;
458 }
459 else
460 {
461 m_diff[i][i] = 1.0;
462 }
463
464 for (int j = i + 1; j < m_coordim; ++j)
465 {
466 d = factors.find(m_factorCoeffDef[i][j]);
467 if (d != factors.end())
468 {
469 m_diff[i][j] = m_diff[j][i] = d->second;
470 }
471 }
472 }
473 m_HasVarCoeffDiff = true;
474 }
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
const StdRegions::ConstFactorType m_factorCoeffDef[3][3]
Array< OneD, Array< OneD, NekDouble > > m_diff
const char *const ConstFactorTypeMap[]
Definition: StdRegions.hpp:391
std::vector< double > d(NPUPPER *NPUPPER)

References ASSERTL1, Nektar::StdRegions::ConstFactorTypeMap, Nektar::UnitTests::d(), Nektar::StdRegions::eFactorCoeffD00, Nektar::StdRegions::eFactorLambda, Nektar::VarcoeffHashingTest::factors, m_coordim, m_diff, m_factorCoeffDef, m_factors, m_HasVarCoeffDiff, and m_lambda.

Referenced by Helmholtz_IterPerExp().

◆ 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 
)
inlinefinaloverridevirtual

Perform operation.

Implements Nektar::Collections::Operator.

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

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

References ASSERTL1, m_coordim, m_derivFac, m_diff, m_dim, m_HasVarCoeffDiff, Nektar::Collections::Operator::m_isDeformed, m_jac, m_lambda, Nektar::Collections::Operator::m_numElmt, Nektar::Collections::Operator::m_stdExp, 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 
)
inlinefinaloverridevirtual

Implements Nektar::Collections::Operator.

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

405 {
406 boost::ignore_unused(dir, input, output, wsp);
407 NEKERROR(ErrorUtil::efatal, "Not valid for this operator.");
408 }
#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 Nektar::ErrorUtil::efatal, and NEKERROR.

Member Data Documentation

◆ m_coordim

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

◆ m_derivFac

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

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

Referenced by Helmholtz_IterPerExp(), and operator()().

◆ m_diff

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

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

Referenced by CheckFactors(), and operator()().

◆ m_dim

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

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

Referenced by Helmholtz_IterPerExp(), and operator()().

◆ 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 481 of file library/Collections/Helmholtz.cpp.

Referenced by CheckFactors(), and Helmholtz_IterPerExp().

◆ m_HasVarCoeffDiff

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

◆ m_jac

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

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

Referenced by Helmholtz_IterPerExp(), and operator()().

◆ m_lambda

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