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
35#include <boost/core/ignore_unused.hpp>
36
37#include <MatrixFreeOps/Operator.hpp>
38
43
44using namespace std;
45
46namespace Nektar
47{
48namespace Collections
49{
50
58
59/**
60 * @brief Helmholtz operator using LocalRegions implementation.
61 */
62class Helmholtz_NoCollection final : public Operator
63{
64public:
66
68 {
69 }
70
75 Array<OneD, NekDouble> &wsp) override final
76 {
77 boost::ignore_unused(entry2, entry3, wsp);
78
79 unsigned int nmodes = m_expList[0]->GetNcoeffs();
81
82 for (int n = 0; n < m_numElmt; ++n)
83 {
85 (m_expList)[n]->DetShapeType(),
86 *(m_expList)[n], m_factors);
87 m_expList[n]->GeneralMatrixOp(entry0 + n * nmodes,
88 tmp = entry1 + n * nmodes, mkey);
89 }
90 }
91
92 void operator()(int dir, const Array<OneD, const NekDouble> &input,
94 Array<OneD, NekDouble> &wsp) override final
95 {
96 boost::ignore_unused(dir, input, output, wsp);
97 NEKERROR(ErrorUtil::efatal, "Not valid for this operator.");
98 }
99
101 int coll_phys_offset) override
102 {
103 boost::ignore_unused(coll_phys_offset);
105 }
106
107protected:
108 int m_dim;
110 vector<StdRegions::StdExpansionSharedPtr> m_expList;
112
113private:
114 Helmholtz_NoCollection(vector<StdRegions::StdExpansionSharedPtr> pCollExp,
117 : Operator(pCollExp, pGeomData, factors)
118 {
119 m_expList = pCollExp;
120 m_dim = pCollExp[0]->GetNumBases();
121 m_coordim = pCollExp[0]->GetCoordim();
122
124 }
125};
126
127/// Factory initialisation for the Helmholtz_NoCollection operators
128OperatorKey Helmholtz_NoCollection::m_typeArr[] = {
131 Helmholtz_NoCollection::create, "Helmholtz_NoCollection_Seg"),
134 Helmholtz_NoCollection::create, "Helmholtz_NoCollection_Tri"),
137 Helmholtz_NoCollection::create, "Helmholtz_NoCollection_NodalTri"),
140 Helmholtz_NoCollection::create, "Helmholtz_NoCollection_Quad"),
143 Helmholtz_NoCollection::create, "Helmholtz_NoCollection_Tet"),
146 Helmholtz_NoCollection::create, "Helmholtz_NoCollection_NodalTet"),
149 Helmholtz_NoCollection::create, "Helmholtz_NoCollection_Pyr"),
152 Helmholtz_NoCollection::create, "Helmholtz_NoCollection_Prism"),
155 Helmholtz_NoCollection::create, "Helmholtz_NoCollection_NodalPrism"),
158 Helmholtz_NoCollection::create, "Helmholtz_NoCollection_Hex")};
159
160/**
161 * @brief Helmholtz operator using LocalRegions implementation.
162 */
163class Helmholtz_IterPerExp final : public Operator
164{
165public:
167
169 {
170 }
171
174 Array<OneD, NekDouble> &output1,
175 Array<OneD, NekDouble> &output2,
176 Array<OneD, NekDouble> &wsp) override final
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;
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 }
401
402 void operator()(int dir, const Array<OneD, const NekDouble> &input,
404 Array<OneD, NekDouble> &wsp) override final
405 {
406 boost::ignore_unused(dir, input, output, wsp);
407 NEKERROR(ErrorUtil::efatal, "Not valid for this operator.");
408 }
409
410 /**
411 * @brief Check the validity of supplied constant factors.
412 *
413 * @param factors Map of factors
414 * @param coll_phys_offset Unused
415 */
417 int coll_phys_offset) override
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
447 for (int i = 0; i < m_coordim; ++i)
448 {
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 }
475
476protected:
479 int m_dim;
492
493private:
494 Helmholtz_IterPerExp(vector<StdRegions::StdExpansionSharedPtr> pCollExp,
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 }
513};
514
515/// Factory initialisation for the Helmholtz_IterPerExp operators
516OperatorKey Helmholtz_IterPerExp::m_typeArr[] = {
519 Helmholtz_IterPerExp::create, "Helmholtz_IterPerExp_Seg"),
522 Helmholtz_IterPerExp::create, "Helmholtz_IterPerExp_Tri"),
525 Helmholtz_IterPerExp::create, "Helmholtz_IterPerExp_NodalTri"),
528 Helmholtz_IterPerExp::create, "Helmholtz_IterPerExp_Quad"),
531 Helmholtz_IterPerExp::create, "Helmholtz_IterPerExp_Tet"),
534 Helmholtz_IterPerExp::create, "Helmholtz_IterPerExp_NodalTet"),
537 Helmholtz_IterPerExp::create, "Helmholtz_IterPerExp_Pyr"),
540 Helmholtz_IterPerExp::create, "Helmholtz_IterPerExp_Prism"),
543 Helmholtz_IterPerExp::create, "Helmholtz_IterPerExp_NodalPrism"),
546 Helmholtz_IterPerExp::create, "Helmholtz_IterPerExp_Hex")};
547
548/**
549 * @brief Helmholtz operator using matrix free operators.
550 */
552{
553public:
555
557 {
558 }
559
561 Array<OneD, NekDouble> &output0,
562 Array<OneD, NekDouble> &output1,
563 Array<OneD, NekDouble> &output2,
564 Array<OneD, NekDouble> &wsp) override final
565 {
566 boost::ignore_unused(output1, output2, wsp);
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()(int dir, const Array<OneD, const NekDouble> &input,
583 Array<OneD, NekDouble> &wsp) override final
584 {
585 boost::ignore_unused(dir, input, output, wsp);
586 NEKERROR(ErrorUtil::efatal, "Not valid for this operator.");
587 }
588
589 /**
590 *
591 */
593 int coll_phys_offset) override
594 {
595 boost::ignore_unused(coll_phys_offset);
596
597 if (factors == m_factors)
598 {
599 return;
600 }
601
603
604 // Set lambda for this call
605 auto x = factors.find(StdRegions::eFactorLambda);
606 ASSERTL1(
607 x != factors.end(),
608 "Constant factor not defined: " +
609 std::string(
611 m_oper->SetLambda(x->second);
612
613 // set constant diffusion coefficients
614 bool isConstVarDiff = false;
616 diff[0] = diff[2] = diff[5] = 1.0;
617
618 auto xd00 = factors.find(StdRegions::eFactorCoeffD00);
619 if (xd00 != factors.end() && xd00->second != 1.0)
620 {
621 isConstVarDiff = true;
622 diff[0] = xd00->second;
623 }
624
625 auto xd01 = factors.find(StdRegions::eFactorCoeffD01);
626 if (xd01 != factors.end() && xd01->second != 0.0)
627 {
628 isConstVarDiff = true;
629 diff[1] = xd01->second;
630 }
631
632 auto xd11 = factors.find(StdRegions::eFactorCoeffD11);
633 if (xd11 != factors.end() && xd11->second != 1.0)
634 {
635 isConstVarDiff = true;
636 diff[2] = xd11->second;
637 }
638
639 auto xd02 = factors.find(StdRegions::eFactorCoeffD02);
640 if (xd02 != factors.end() && xd02->second != 0.0)
641 {
642 isConstVarDiff = true;
643 diff[3] = xd02->second;
644 }
645
646 auto xd12 = factors.find(StdRegions::eFactorCoeffD12);
647 if (xd12 != factors.end() && xd12->second != 0.0)
648 {
649 isConstVarDiff = true;
650 diff[4] = xd12->second;
651 }
652
653 auto xd22 = factors.find(StdRegions::eFactorCoeffD22);
654 if (xd22 != factors.end() && xd22->second != 1.0)
655 {
656 isConstVarDiff = true;
657 diff[5] = xd22->second;
658 }
659
660 if (isConstVarDiff)
661 {
662 m_oper->SetConstVarDiffusion(diff);
663 }
664
665 // set random here for fn variable diffusion
666 auto k = factors.find(StdRegions::eFactorTau);
667 if (k != factors.end() && k->second != 0.0)
668 {
669 m_oper->SetVarDiffusion(diff);
670 }
671 }
672
673private:
674 std::shared_ptr<MatrixFree::Helmholtz> m_oper;
675 unsigned int m_nmtot;
677
678 Helmholtz_MatrixFree(vector<StdRegions::StdExpansionSharedPtr> pCollExp,
681 : Operator(pCollExp, pGeomData, factors),
682 MatrixFreeOneInOneOut(pCollExp[0]->GetStdExp()->GetNcoeffs(),
683 pCollExp[0]->GetStdExp()->GetNcoeffs(),
684 pCollExp.size())
685 {
686
687 m_nmtot = m_numElmt * pCollExp[0]->GetStdExp()->GetNcoeffs();
688
689 const auto dim = pCollExp[0]->GetStdExp()->GetShapeDimension();
690
691 // Basis vector.
692 std::vector<LibUtilities::BasisSharedPtr> basis(dim);
693 for (auto i = 0; i < dim; ++i)
694 {
695 basis[i] = pCollExp[0]->GetBasis(i);
696 }
697
698 // Get shape type
699 auto shapeType = pCollExp[0]->GetStdExp()->DetShapeType();
700
701 // Generate operator string and create operator.
702 std::string op_string = "Helmholtz";
703 op_string += MatrixFree::GetOpstring(shapeType, m_isDeformed);
705 op_string, basis, m_nElmtPad);
706
707 // Set Jacobian
708 oper->SetJac(pGeomData->GetJacInterLeave(pCollExp, m_nElmtPad));
709
710 // Store derivative factor
711 oper->SetDF(pGeomData->GetDerivFactorsInterLeave(pCollExp, m_nElmtPad));
712
713 m_oper = std::dynamic_pointer_cast<MatrixFree::Helmholtz>(oper);
714 ASSERTL0(m_oper, "Failed to cast pointer.");
715
716 // Set factors
718 this->CheckFactors(factors, 0);
719 }
720};
721
722/// Factory initialisation for the Helmholtz_MatrixFree operators
723OperatorKey Helmholtz_MatrixFree::m_typeArr[] = {
726 Helmholtz_MatrixFree::create, "Helmholtz_MatrixFree_Quad"),
729 Helmholtz_MatrixFree::create, "Helmholtz_MatrixFree_Tri"),
732 Helmholtz_MatrixFree::create, "Helmholtz_MatrixFree_Hex"),
735 Helmholtz_MatrixFree::create, "Helmholtz_MatrixFree_Prism"),
738 Helmholtz_MatrixFree::create, "Helmholtz_MatrixFree_Pyr"),
741 Helmholtz_MatrixFree::create, "Helmholtz_MatrixFree_Tet"),
742};
743
744} // namespace Collections
745} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
Definition: ErrorUtil.hpp:209
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
#define OPERATOR_CREATE(cname)
Definition: Operator.h:45
Helmholtz operator using LocalRegions implementation.
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.
void operator()(int dir, const Array< OneD, const NekDouble > &input, Array< OneD, NekDouble > &output, Array< OneD, NekDouble > &wsp) override final
const StdRegions::ConstFactorType m_factorCoeffDef[3][3]
Helmholtz_IterPerExp(vector< StdRegions::StdExpansionSharedPtr > pCollExp, CoalescedGeomDataSharedPtr pGeomData, StdRegions::FactorMap factors)
virtual 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()(const Array< OneD, const NekDouble > &input, Array< OneD, NekDouble > &output0, Array< OneD, NekDouble > &output1, Array< OneD, NekDouble > &output2, Array< OneD, NekDouble > &wsp) override final
Perform operation.
std::shared_ptr< MatrixFree::Helmholtz > m_oper
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 the supplied factor map.
Helmholtz_MatrixFree(vector< StdRegions::StdExpansionSharedPtr > pCollExp, CoalescedGeomDataSharedPtr pGeomData, StdRegions::FactorMap factors)
Helmholtz operator using LocalRegions implementation.
void operator()(int dir, const Array< OneD, const NekDouble > &input, Array< OneD, NekDouble > &output, Array< OneD, NekDouble > &wsp) override final
void operator()(const Array< OneD, const NekDouble > &entry0, Array< OneD, NekDouble > &entry1, Array< OneD, NekDouble > &entry2, Array< OneD, NekDouble > &entry3, Array< OneD, NekDouble > &wsp) override final
Perform operation.
Helmholtz_NoCollection(vector< StdRegions::StdExpansionSharedPtr > pCollExp, CoalescedGeomDataSharedPtr pGeomData, StdRegions::FactorMap factors)
virtual void CheckFactors(StdRegions::FactorMap factors, int coll_phys_offset) override
Check the validity of the supplied factor map.
vector< StdRegions::StdExpansionSharedPtr > m_expList
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:119
StdRegions::StdExpansionSharedPtr m_stdExp
Definition: Operator.h:165
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
std::tuple< LibUtilities::ShapeType, OperatorType, ImplementationType, ExpansionIsNodal > OperatorKey
Key for describing an Operator.
Definition: Operator.h:177
std::shared_ptr< CoalescedGeomData > CoalescedGeomDataSharedPtr
OperatorFactory & GetOperatorFactory()
Returns the singleton Operator factory object.
Definition: Operator.cpp:117
std::vector< PointsKey > PointsKeyVector
Definition: Points.h:236
static FactorMap NullFactorMap
Definition: StdRegions.hpp:413
const char *const ConstFactorTypeMap[]
Definition: StdRegions.hpp:391
ConstFactorMap FactorMap
Definition: StdRegions.hpp:412
std::vector< double > d(NPUPPER *NPUPPER)
StdRegions::ConstFactorMap factors
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
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.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
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1191