Nektar++
library/Collections/Helmholtz.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: Helmholtz.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// Permission is hereby granted, free of charge, to any person obtaining a
14// copy of this software and associated documentation files (the "Software"),
15// to deal in the Software without restriction, including without limitation
16// the rights to use, copy, modify, merge, publish, distribute, sublicense,
17// and/or sell copies of the Software, and to permit persons to whom the
18// Software is furnished to do so, subject to the following conditions:
19//
20// The above copyright notice and this permission notice shall be included
21// in all copies or substantial portions of the Software.
22//
23// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29// DEALINGS IN THE SOFTWARE.
30//
31// Description: Helmholtz operator implementations
32//
33///////////////////////////////////////////////////////////////////////////////
34
39#include <MatrixFreeOps/Operator.hpp>
40
41using namespace std;
42
43namespace Nektar::Collections
44{
45
53
54/**
55 * @brief Helmholtz help class to calculate the size of the collection
56 * that is given as an input and as an output to the Helmholtz Operator. The
57 * size evaluation takes into account that the evaluation of the Helmholtz
58 * operator takes input from the coeff space and gives the output in the coeff
59 * space too.
60 */
61class Helmholtz_Helper : virtual public Operator
62{
63protected:
65 {
66 // expect input to be number of elements by the number of coefficients
67 m_inputSize = m_numElmt * m_stdExp->GetNcoeffs();
68 // expect output to be number of elements by the number of coefficients
69 // computation is from coeff space to coeff space
71 }
72};
73
74/**
75 * @brief Helmholtz operator using LocalRegions implementation.
76 */
78{
79public:
81
82 ~Helmholtz_NoCollection() final = default;
83
84 void operator()(const Array<OneD, const NekDouble> &entry0,
85 Array<OneD, NekDouble> &entry1,
86 [[maybe_unused]] Array<OneD, NekDouble> &entry2,
87 [[maybe_unused]] Array<OneD, NekDouble> &entry3,
88 [[maybe_unused]] Array<OneD, NekDouble> &wsp) final
89 {
90 unsigned int nmodes = m_expList[0]->GetNcoeffs();
92
93 for (int n = 0; n < m_numElmt; ++n)
94 {
96 (m_expList)[n]->DetShapeType(),
97 *(m_expList)[n], m_factors);
98 m_expList[n]->GeneralMatrixOp(entry0 + n * nmodes,
99 tmp = entry1 + n * nmodes, mkey);
100 }
101 }
102
103 void operator()([[maybe_unused]] int dir,
104 [[maybe_unused]] const Array<OneD, const NekDouble> &input,
105 [[maybe_unused]] Array<OneD, NekDouble> &output,
106 [[maybe_unused]] Array<OneD, NekDouble> &wsp) final
107 {
108 NEKERROR(ErrorUtil::efatal, "Not valid for this operator.");
109 }
110
112 [[maybe_unused]] int coll_phys_offset) override
113 {
115 }
116
117protected:
118 int m_dim;
120 vector<StdRegions::StdExpansionSharedPtr> m_expList;
122
123private:
124 Helmholtz_NoCollection(vector<StdRegions::StdExpansionSharedPtr> pCollExp,
127 : Operator(pCollExp, pGeomData, factors), Helmholtz_Helper()
128 {
129 m_expList = pCollExp;
130 m_dim = pCollExp[0]->GetNumBases();
131 m_coordim = pCollExp[0]->GetCoordim();
132
134 }
135};
136
137/// Factory initialisation for the Helmholtz_NoCollection operators
138OperatorKey Helmholtz_NoCollection::m_typeArr[] = {
141 Helmholtz_NoCollection::create, "Helmholtz_NoCollection_Seg"),
144 Helmholtz_NoCollection::create, "Helmholtz_NoCollection_Tri"),
147 Helmholtz_NoCollection::create, "Helmholtz_NoCollection_NodalTri"),
150 Helmholtz_NoCollection::create, "Helmholtz_NoCollection_Quad"),
153 Helmholtz_NoCollection::create, "Helmholtz_NoCollection_Tet"),
156 Helmholtz_NoCollection::create, "Helmholtz_NoCollection_NodalTet"),
159 Helmholtz_NoCollection::create, "Helmholtz_NoCollection_Pyr"),
162 Helmholtz_NoCollection::create, "Helmholtz_NoCollection_Prism"),
165 Helmholtz_NoCollection::create, "Helmholtz_NoCollection_NodalPrism"),
168 Helmholtz_NoCollection::create, "Helmholtz_NoCollection_Hex")};
169
170/**
171 * @brief Helmholtz operator using LocalRegions implementation.
172 */
173class Helmholtz_IterPerExp final : virtual public Operator, Helmholtz_Helper
174{
175public:
177
178 ~Helmholtz_IterPerExp() final = default;
179
180 void operator()(const Array<OneD, const NekDouble> &input,
181 Array<OneD, NekDouble> &output,
182 [[maybe_unused]] Array<OneD, NekDouble> &output1,
183 [[maybe_unused]] Array<OneD, NekDouble> &output2,
184 Array<OneD, NekDouble> &wsp) final
185 {
186 const int nCoeffs = m_stdExp->GetNcoeffs();
187 const int nPhys = m_stdExp->GetTotPoints();
188
189 ASSERTL1(input.size() >= m_numElmt * nCoeffs,
190 "input array size is insufficient");
191 ASSERTL1(output.size() >= m_numElmt * nCoeffs,
192 "output array size is insufficient");
193
194 Array<OneD, NekDouble> tmpphys, t1;
197
198 tmpphys = wsp;
199 for (int i = 1; i < m_coordim + 1; ++i)
200 {
201 dtmp[i - 1] = wsp + i * nPhys;
202 tmp[i - 1] = wsp + (i + m_coordim) * nPhys;
203 }
204
205 for (int i = 0; i < m_numElmt; ++i)
206 {
207 m_stdExp->BwdTrans(input + i * nCoeffs, tmpphys);
208
209 // local derivative
210 m_stdExp->PhysDeriv(tmpphys, dtmp[0], dtmp[1], dtmp[2]);
211
212 // determine mass matrix term
213 if (m_isDeformed)
214 {
215 Vmath::Vmul(nPhys, m_jac + i * nPhys, 1, tmpphys, 1, tmpphys,
216 1);
217 }
218 else
219 {
220 Vmath::Smul(nPhys, m_jac[i], tmpphys, 1, tmpphys, 1);
221 }
222
223 m_stdExp->IProductWRTBase(tmpphys, t1 = output + i * nCoeffs);
224 Vmath::Smul(nCoeffs, m_lambda, output + i * nCoeffs, 1,
225 t1 = output + i * nCoeffs, 1);
226
227 if (m_isDeformed)
228 {
229 // calculate full derivative
230 for (int j = 0; j < m_coordim; ++j)
231 {
232 Vmath::Vmul(nPhys,
233 m_derivFac[j * m_dim].origin() + i * nPhys, 1,
234 &dtmp[0][0], 1, &tmp[j][0], 1);
235
236 for (int k = 1; k < m_dim; ++k)
237 {
239 nPhys,
240 m_derivFac[j * m_dim + k].origin() + i * nPhys, 1,
241 &dtmp[k][0], 1, &tmp[j][0], 1, &tmp[j][0], 1);
242 }
243 }
244
246 {
247 // calculate dtmp[i] = dx/dxi sum_j diff[0][j] tmp[j]
248 // + dy/dxi sum_j diff[1][j] tmp[j]
249 // + dz/dxi sum_j diff[2][j] tmp[j]
250
251 // First term
252 Vmath::Smul(nPhys, m_diff[0][0], &tmp[0][0], 1, &tmpphys[0],
253 1);
254 for (int l = 1; l < m_coordim; ++l)
255 {
256 Vmath::Svtvp(nPhys, m_diff[0][l], &tmp[l][0], 1,
257 &tmpphys[0], 1, &tmpphys[0], 1);
258 }
259
260 for (int j = 0; j < m_dim; ++j)
261 {
262 Vmath::Vmul(nPhys, m_derivFac[j].origin() + i * nPhys,
263 1, &tmpphys[0], 1, &dtmp[j][0], 1);
264 }
265
266 // Second and third terms
267 for (int k = 1; k < m_coordim; ++k)
268 {
269 Vmath::Smul(nPhys, m_diff[k][0], &tmp[0][0], 1,
270 &tmpphys[0], 1);
271 for (int l = 1; l < m_coordim; ++l)
272 {
273 Vmath::Svtvp(nPhys, m_diff[k][l], &tmp[l][0], 1,
274 &tmpphys[0], 1, &tmpphys[0], 1);
275 }
276
277 for (int j = 0; j < m_dim; ++j)
278 {
279 Vmath::Vvtvp(nPhys,
280 m_derivFac[j + k * m_dim].origin() +
281 i * nPhys,
282 1, &tmpphys[0], 1, &dtmp[j][0], 1,
283 &dtmp[j][0], 1);
284 }
285 }
286 }
287 else
288 {
289 // calculate dx/dxi tmp[0] + dy/dxi tmp[1]
290 // + dz/dxi tmp[2]
291 for (int j = 0; j < m_dim; ++j)
292 {
293 Vmath::Vmul(nPhys, m_derivFac[j].origin() + i * nPhys,
294 1, &tmp[0][0], 1, &dtmp[j][0], 1);
295
296 for (int k = 1; k < m_coordim; ++k)
297 {
298 Vmath::Vvtvp(nPhys,
299 m_derivFac[j + k * m_dim].origin() +
300 i * nPhys,
301 1, &tmp[k][0], 1, &dtmp[j][0], 1,
302 &dtmp[j][0], 1);
303 }
304 }
305 }
306
307 // calculate Iproduct WRT Std Deriv
308 for (int j = 0; j < m_dim; ++j)
309 {
310 // multiply by Jacobian
311 Vmath::Vmul(nPhys, m_jac + i * nPhys, 1, dtmp[j], 1,
312 dtmp[j], 1);
313
314 m_stdExp->IProductWRTDerivBase(j, dtmp[j], tmp[0]);
315 Vmath::Vadd(nCoeffs, tmp[0], 1, output + i * nCoeffs, 1,
316 t1 = output + i * nCoeffs, 1);
317 }
318 }
319 else
320 {
321 // calculate full derivative
322 for (int j = 0; j < m_coordim; ++j)
323 {
324 Vmath::Smul(nPhys, m_derivFac[j * m_dim][i], &dtmp[0][0], 1,
325 &tmp[j][0], 1);
326
327 for (int k = 1; k < m_dim; ++k)
328 {
329 Vmath::Svtvp(nPhys, m_derivFac[j * m_dim + k][i],
330 &dtmp[k][0], 1, &tmp[j][0], 1, &tmp[j][0],
331 1);
332 }
333 }
334
336 {
337 // calculate dtmp[i] = dx/dxi sum_j diff[0][j] tmp[j]
338 // + dy/dxi sum_j diff[1][j] tmp[j]
339 // + dz/dxi sum_j diff[2][j] tmp[j]
340
341 // First term
342 Vmath::Smul(nPhys, m_diff[0][0], &tmp[0][0], 1, &tmpphys[0],
343 1);
344 for (int l = 1; l < m_coordim; ++l)
345 {
346 Vmath::Svtvp(nPhys, m_diff[0][l], &tmp[l][0], 1,
347 &tmpphys[0], 1, &tmpphys[0], 1);
348 }
349
350 for (int j = 0; j < m_dim; ++j)
351 {
352 Vmath::Smul(nPhys, m_derivFac[j][i], &tmpphys[0], 1,
353 &dtmp[j][0], 1);
354 }
355
356 // Second and third terms
357 for (int k = 1; k < m_coordim; ++k)
358 {
359 Vmath::Smul(nPhys, m_diff[k][0], &tmp[0][0], 1,
360 &tmpphys[0], 1);
361 for (int l = 1; l < m_coordim; ++l)
362 {
363 Vmath::Svtvp(nPhys, m_diff[k][l], &tmp[l][0], 1,
364 &tmpphys[0], 1, &tmpphys[0], 1);
365 }
366
367 for (int j = 0; j < m_dim; ++j)
368 {
369 Vmath::Svtvp(nPhys, m_derivFac[j + k * m_dim][i],
370 &tmpphys[0], 1, &dtmp[j][0], 1,
371 &dtmp[j][0], 1);
372 }
373 }
374 }
375 else
376 {
377 // calculate dx/dxi tmp[0] + dy/dxi tmp[2]
378 // + dz/dxi tmp[3]
379 for (int j = 0; j < m_dim; ++j)
380 {
381 Vmath::Smul(nPhys, m_derivFac[j][i], &tmp[0][0], 1,
382 &dtmp[j][0], 1);
383
384 for (int k = 1; k < m_coordim; ++k)
385 {
386 Vmath::Svtvp(nPhys, m_derivFac[j + k * m_dim][i],
387 &tmp[k][0], 1, &dtmp[j][0], 1,
388 &dtmp[j][0], 1);
389 }
390 }
391 }
392
393 // calculate Iproduct WRT Std Deriv
394 for (int j = 0; j < m_dim; ++j)
395 {
396 // multiply by Jacobian
397 Vmath::Smul(nPhys, m_jac[i], dtmp[j], 1, dtmp[j], 1);
398
399 m_stdExp->IProductWRTDerivBase(j, dtmp[j], tmp[0]);
400 Vmath::Vadd(nCoeffs, tmp[0], 1, output + i * nCoeffs, 1,
401 t1 = output + i * nCoeffs, 1);
402 }
403 }
404 }
405 }
406
407 void operator()([[maybe_unused]] int dir,
408 [[maybe_unused]] const Array<OneD, const NekDouble> &input,
409 [[maybe_unused]] Array<OneD, NekDouble> &output,
410 [[maybe_unused]] Array<OneD, NekDouble> &wsp) final
411 {
412 NEKERROR(ErrorUtil::efatal, "Not valid for this operator.");
413 }
414
415 /**
416 * @brief Check the validity of supplied constant factors.
417 *
418 * @param factors Map of factors
419 * @param coll_phys_offset Unused
420 */
422 [[maybe_unused]] int coll_phys_offset) override
423 {
424 // If match previous factors, nothing to do.
425 if (m_factors == factors)
426 {
427 return;
428 }
429
431
432 // Check Lambda constant of Helmholtz operator
433 auto x = factors.find(StdRegions::eFactorLambda);
434 ASSERTL1(
435 x != factors.end(),
436 "Constant factor not defined: " +
437 std::string(
439 m_lambda = x->second;
440
441 // If varcoeffs not supplied, nothing else to do.
442 m_HasVarCoeffDiff = false;
444 if (d == factors.end())
445 {
446 return;
447 }
448
450 for (int i = 0; i < m_coordim; ++i)
451 {
453 }
454
455 for (int i = 0; i < m_coordim; ++i)
456 {
457 d = factors.find(m_factorCoeffDef[i][i]);
458 if (d != factors.end())
459 {
460 m_diff[i][i] = d->second;
461 }
462 else
463 {
464 m_diff[i][i] = 1.0;
465 }
466
467 for (int j = i + 1; j < m_coordim; ++j)
468 {
469 d = factors.find(m_factorCoeffDef[i][j]);
470 if (d != factors.end())
471 {
472 m_diff[i][j] = m_diff[j][i] = d->second;
473 }
474 }
475 }
476 m_HasVarCoeffDiff = true;
477 }
478
479protected:
482 int m_dim;
495
496private:
497 Helmholtz_IterPerExp(vector<StdRegions::StdExpansionSharedPtr> pCollExp,
500 : Operator(pCollExp, pGeomData, factors), Helmholtz_Helper()
501 {
502 m_dim = pCollExp[0]->GetShapeDimension();
503 m_coordim = pCollExp[0]->GetCoordim();
504 int nqtot = m_stdExp->GetTotPoints();
505
506 m_derivFac = pGeomData->GetDerivFactors(pCollExp);
507 m_jac = pGeomData->GetJac(pCollExp);
508 m_wspSize = (2 * m_coordim + 1) * nqtot;
509
510 m_lambda = 1.0;
511 m_HasVarCoeffDiff = false;
513 this->CheckFactors(factors, 0);
514 }
515};
516
517/// Factory initialisation for the Helmholtz_IterPerExp operators
518OperatorKey Helmholtz_IterPerExp::m_typeArr[] = {
521 Helmholtz_IterPerExp::create, "Helmholtz_IterPerExp_Seg"),
524 Helmholtz_IterPerExp::create, "Helmholtz_IterPerExp_Tri"),
527 Helmholtz_IterPerExp::create, "Helmholtz_IterPerExp_NodalTri"),
530 Helmholtz_IterPerExp::create, "Helmholtz_IterPerExp_Quad"),
533 Helmholtz_IterPerExp::create, "Helmholtz_IterPerExp_Tet"),
536 Helmholtz_IterPerExp::create, "Helmholtz_IterPerExp_NodalTet"),
539 Helmholtz_IterPerExp::create, "Helmholtz_IterPerExp_Pyr"),
542 Helmholtz_IterPerExp::create, "Helmholtz_IterPerExp_Prism"),
545 Helmholtz_IterPerExp::create, "Helmholtz_IterPerExp_NodalPrism"),
548 Helmholtz_IterPerExp::create, "Helmholtz_IterPerExp_Hex")};
549
550/**
551 * @brief Helmholtz operator using matrix free operators.
552 */
553class Helmholtz_MatrixFree final : virtual public Operator,
556{
557public:
559
560 ~Helmholtz_MatrixFree() final = default;
561
562 void operator()(const Array<OneD, const NekDouble> &input,
563 Array<OneD, NekDouble> &output0,
564 [[maybe_unused]] Array<OneD, NekDouble> &output1,
565 [[maybe_unused]] Array<OneD, NekDouble> &output2,
566 [[maybe_unused]] Array<OneD, NekDouble> &wsp) final
567 {
568 if (m_isPadded)
569 {
570 // copy into padded vector
571 Vmath::Vcopy(m_nmtot, input, 1, m_input, 1);
572 (*m_oper)(m_input, m_output);
573 Vmath::Vcopy(m_nmtot, m_output, 1, output0, 1);
574 }
575 else
576 {
577 (*m_oper)(input, output0);
578 }
579 }
580
581 void operator()([[maybe_unused]] int dir,
582 [[maybe_unused]] const Array<OneD, const NekDouble> &input,
583 [[maybe_unused]] Array<OneD, NekDouble> &output,
584 [[maybe_unused]] Array<OneD, NekDouble> &wsp) final
585 {
586 NEKERROR(ErrorUtil::efatal, "Not valid for this operator.");
587 }
588
589 /**
590 *
591 */
593 [[maybe_unused]] int coll_phys_offset) override
594 {
595 if (factors == m_factors)
596 {
597 return;
598 }
599
601
602 // Set lambda for this call
603 auto x = factors.find(StdRegions::eFactorLambda);
604 ASSERTL1(
605 x != factors.end(),
606 "Constant factor not defined: " +
607 std::string(
609 m_oper->SetLambda(x->second);
610
611 // set constant diffusion coefficients
612 bool isConstVarDiff = false;
614 diff[0] = diff[2] = diff[5] = 1.0;
615
616 auto xd00 = factors.find(StdRegions::eFactorCoeffD00);
617 if (xd00 != factors.end() && xd00->second != 1.0)
618 {
619 isConstVarDiff = true;
620 diff[0] = xd00->second;
621 }
622
623 auto xd01 = factors.find(StdRegions::eFactorCoeffD01);
624 if (xd01 != factors.end() && xd01->second != 0.0)
625 {
626 isConstVarDiff = true;
627 diff[1] = xd01->second;
628 }
629
630 auto xd11 = factors.find(StdRegions::eFactorCoeffD11);
631 if (xd11 != factors.end() && xd11->second != 1.0)
632 {
633 isConstVarDiff = true;
634 diff[2] = xd11->second;
635 }
636
637 auto xd02 = factors.find(StdRegions::eFactorCoeffD02);
638 if (xd02 != factors.end() && xd02->second != 0.0)
639 {
640 isConstVarDiff = true;
641 diff[3] = xd02->second;
642 }
643
644 auto xd12 = factors.find(StdRegions::eFactorCoeffD12);
645 if (xd12 != factors.end() && xd12->second != 0.0)
646 {
647 isConstVarDiff = true;
648 diff[4] = xd12->second;
649 }
650
651 auto xd22 = factors.find(StdRegions::eFactorCoeffD22);
652 if (xd22 != factors.end() && xd22->second != 1.0)
653 {
654 isConstVarDiff = true;
655 diff[5] = xd22->second;
656 }
657
658 if (isConstVarDiff)
659 {
660 m_oper->SetConstVarDiffusion(diff);
661 }
662
663 // set random here for fn variable diffusion
664 auto k = factors.find(StdRegions::eFactorTau);
665 if (k != factors.end() && k->second != 0.0)
666 {
667 m_oper->SetVarDiffusion(diff);
668 }
669 }
670
671private:
672 std::shared_ptr<MatrixFree::Helmholtz> m_oper;
673 unsigned int m_nmtot;
675
676 Helmholtz_MatrixFree(vector<StdRegions::StdExpansionSharedPtr> pCollExp,
679 : Operator(pCollExp, pGeomData, factors),
680 MatrixFreeOneInOneOut(pCollExp[0]->GetStdExp()->GetNcoeffs(),
681 pCollExp[0]->GetStdExp()->GetNcoeffs(),
682 pCollExp.size()),
684 {
685
686 m_nmtot = m_numElmt * pCollExp[0]->GetStdExp()->GetNcoeffs();
687
688 const auto dim = pCollExp[0]->GetStdExp()->GetShapeDimension();
689
690 // Basis vector.
691 std::vector<LibUtilities::BasisSharedPtr> basis(dim);
692 for (auto i = 0; i < dim; ++i)
693 {
694 basis[i] = pCollExp[0]->GetBasis(i);
695 }
696
697 // Get shape type
698 auto shapeType = pCollExp[0]->GetStdExp()->DetShapeType();
699
700 // Generate operator string and create operator.
701 std::string op_string = "Helmholtz";
702 op_string += MatrixFree::GetOpstring(shapeType, m_isDeformed);
704 op_string, basis, m_nElmtPad);
705
706 // Set Jacobian
707 oper->SetJac(pGeomData->GetJacInterLeave(pCollExp, m_nElmtPad));
708
709 // Store derivative factor
710 oper->SetDF(pGeomData->GetDerivFactorsInterLeave(pCollExp, m_nElmtPad));
711
712 m_oper = std::dynamic_pointer_cast<MatrixFree::Helmholtz>(oper);
713 ASSERTL0(m_oper, "Failed to cast pointer.");
714
715 // Set factors
717 this->CheckFactors(factors, 0);
718 }
719};
720
721/// Factory initialisation for the Helmholtz_MatrixFree operators
722OperatorKey Helmholtz_MatrixFree::m_typeArr[] = {
725 Helmholtz_MatrixFree::create, "Helmholtz_MatrixFree_Quad"),
728 Helmholtz_MatrixFree::create, "Helmholtz_MatrixFree_Tri"),
731 Helmholtz_MatrixFree::create, "Helmholtz_MatrixFree_Hex"),
734 Helmholtz_MatrixFree::create, "Helmholtz_MatrixFree_Prism"),
737 Helmholtz_MatrixFree::create, "Helmholtz_MatrixFree_Pyr"),
740 Helmholtz_MatrixFree::create, "Helmholtz_MatrixFree_Tet"),
741};
742
743} // namespace Nektar::Collections
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:202
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
#define OPERATOR_CREATE(cname)
Definition: Operator.h:43
Helmholtz help class to calculate the size of the collection that is given as an input and as an outp...
Helmholtz operator using LocalRegions implementation.
void operator()(int dir, const Array< OneD, const NekDouble > &input, Array< OneD, NekDouble > &output, Array< OneD, NekDouble > &wsp) final
const StdRegions::ConstFactorType m_factorCoeffDef[3][3]
Helmholtz_IterPerExp(vector< StdRegions::StdExpansionSharedPtr > pCollExp, CoalescedGeomDataSharedPtr pGeomData, StdRegions::FactorMap factors)
void CheckFactors(StdRegions::FactorMap factors, int coll_phys_offset) override
Check the validity of supplied constant factors.
Array< OneD, Array< OneD, NekDouble > > m_diff
Helmholtz operator using matrix free operators.
void operator()(int dir, const Array< OneD, const NekDouble > &input, Array< OneD, NekDouble > &output, Array< OneD, NekDouble > &wsp) final
void CheckFactors(StdRegions::FactorMap factors, int coll_phys_offset) override
Check the validity of the supplied factor map.
std::shared_ptr< MatrixFree::Helmholtz > m_oper
Helmholtz_MatrixFree(vector< StdRegions::StdExpansionSharedPtr > pCollExp, CoalescedGeomDataSharedPtr pGeomData, StdRegions::FactorMap factors)
Helmholtz operator using LocalRegions implementation.
Helmholtz_NoCollection(vector< StdRegions::StdExpansionSharedPtr > pCollExp, CoalescedGeomDataSharedPtr pGeomData, StdRegions::FactorMap factors)
void operator()(int dir, const Array< OneD, const NekDouble > &input, Array< OneD, NekDouble > &output, Array< OneD, NekDouble > &wsp) final
vector< StdRegions::StdExpansionSharedPtr > m_expList
void CheckFactors(StdRegions::FactorMap factors, int coll_phys_offset) override
Check the validity of the supplied factor map.
unsigned int m_nElmtPad
size after padding
Array< OneD, NekDouble > m_input
padded input/output vectors
Base class for operators on a collection of elements.
Definition: Operator.h:133
StdRegions::StdExpansionSharedPtr m_stdExp
Definition: Operator.h:188
unsigned int m_numElmt
number of elements that the operator is applied on
Definition: Operator.h:190
unsigned int m_outputSize
number of modes or quadrature points that are taken as output from an operator
Definition: Operator.h:198
unsigned int m_inputSize
number of modes or quadrature points that are passed as input to an operator
Definition: Operator.h:195
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:197
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:143
std::tuple< LibUtilities::ShapeType, OperatorType, ImplementationType, ExpansionIsNodal > OperatorKey
Key for describing an Operator.
Definition: Operator.h:115
std::shared_ptr< CoalescedGeomData > CoalescedGeomDataSharedPtr
OperatorFactory & GetOperatorFactory()
Returns the singleton Operator factory object.
Definition: Operator.cpp:44
const NekPoint< NekDouble > origin
static FactorMap NullFactorMap
Definition: StdRegions.hpp:407
const char *const ConstFactorTypeMap[]
Definition: StdRegions.hpp:385
ConstFactorMap FactorMap
Definition: StdRegions.hpp:406
std::vector< double > d(NPUPPER *NPUPPER)
StdRegions::ConstFactorMap factors
double NekDouble
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.hpp:72
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.hpp:396
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.hpp:366
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.hpp:180
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.hpp:100
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825