35 #include <boost/core/ignore_unused.hpp>
37 #include <MatrixFreeOps/Operator.hpp>
77 boost::ignore_unused(entry2, entry3, wsp);
79 unsigned int nmodes = m_expList[0]->GetNcoeffs();
82 for (
int n = 0; n < m_numElmt; ++n)
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);
96 boost::ignore_unused(dir, input, output, wsp);
97 NEKERROR(ErrorUtil::efatal,
"Not valid for this operator.");
101 int coll_phys_offset)
override
103 boost::ignore_unused(coll_phys_offset);
117 :
Operator(pCollExp, pGeomData, factors)
119 m_expList = pCollExp;
120 m_dim = pCollExp[0]->GetNumBases();
121 m_coordim = pCollExp[0]->GetCoordim();
128 OperatorKey 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")};
178 boost::ignore_unused(output1, output2);
180 const int nCoeffs = m_stdExp->GetNcoeffs();
181 const int nPhys = m_stdExp->GetTotPoints();
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");
193 for (
int i = 1; i < m_coordim + 1; ++i)
195 dtmp[i - 1] = wsp + i * nPhys;
196 tmp[i - 1] = wsp + (i + m_coordim) * nPhys;
199 for (
int i = 0; i < m_numElmt; ++i)
201 m_stdExp->BwdTrans(input + i * nCoeffs, tmpphys);
204 m_stdExp->PhysDeriv(tmpphys, dtmp[0], dtmp[1], dtmp[2]);
209 Vmath::Vmul(nPhys, m_jac + i * nPhys, 1, tmpphys, 1, tmpphys,
214 Vmath::Smul(nPhys, m_jac[i], tmpphys, 1, tmpphys, 1);
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);
224 for (
int j = 0; j < m_coordim; ++j)
227 m_derivFac[j * m_dim].origin() + i * nPhys, 1,
228 &dtmp[0][0], 1, &tmp[j][0], 1);
230 for (
int k = 1; k < m_dim; ++k)
234 m_derivFac[j * m_dim + k].origin() + i * nPhys, 1,
235 &dtmp[k][0], 1, &tmp[j][0], 1, &tmp[j][0], 1);
239 if (m_HasVarCoeffDiff)
246 Vmath::Smul(nPhys, m_diff[0][0], &tmp[0][0], 1, &tmpphys[0],
248 for (
int l = 1; l < m_coordim; ++l)
251 &tmpphys[0], 1, &tmpphys[0], 1);
254 for (
int j = 0; j < m_dim; ++j)
256 Vmath::Vmul(nPhys, m_derivFac[j].origin() + i * nPhys,
257 1, &tmpphys[0], 1, &dtmp[j][0], 1);
261 for (
int k = 1; k < m_coordim; ++k)
265 for (
int l = 1; l < m_coordim; ++l)
268 &tmpphys[0], 1, &tmpphys[0], 1);
271 for (
int j = 0; j < m_dim; ++j)
274 m_derivFac[j + k * m_dim].origin() +
276 1, &tmpphys[0], 1, &dtmp[j][0], 1,
285 for (
int j = 0; j < m_dim; ++j)
287 Vmath::Vmul(nPhys, m_derivFac[j].origin() + i * nPhys,
288 1, &tmp[0][0], 1, &dtmp[j][0], 1);
290 for (
int k = 1; k < m_coordim; ++k)
293 m_derivFac[j + k * m_dim].origin() +
295 1, &tmp[k][0], 1, &dtmp[j][0], 1,
302 for (
int j = 0; j < m_dim; ++j)
306 Vmath::Vmul(nPhys, m_jac + i * nPhys, 1, dtmp[j], 1,
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);
317 for (
int j = 0; j < m_coordim; ++j)
319 Vmath::Smul(nPhys, m_derivFac[j * m_dim][i], &dtmp[0][0], 1,
322 for (
int k = 1; k < m_dim; ++k)
325 &dtmp[k][0], 1, &tmp[j][0], 1, &tmp[j][0],
330 if (m_HasVarCoeffDiff)
337 Vmath::Smul(nPhys, m_diff[0][0], &tmp[0][0], 1, &tmpphys[0],
339 for (
int l = 1; l < m_coordim; ++l)
342 &tmpphys[0], 1, &tmpphys[0], 1);
345 for (
int j = 0; j < m_dim; ++j)
347 Vmath::Smul(nPhys, m_derivFac[j][i], &tmpphys[0], 1,
352 for (
int k = 1; k < m_coordim; ++k)
356 for (
int l = 1; l < m_coordim; ++l)
359 &tmpphys[0], 1, &tmpphys[0], 1);
362 for (
int j = 0; j < m_dim; ++j)
365 &tmpphys[0], 1, &dtmp[j][0], 1,
374 for (
int j = 0; j < m_dim; ++j)
376 Vmath::Smul(nPhys, m_derivFac[j][i], &tmp[0][0], 1,
379 for (
int k = 1; k < m_coordim; ++k)
382 &tmp[k][0], 1, &dtmp[j][0], 1,
389 for (
int j = 0; j < m_dim; ++j)
392 Vmath::Smul(nPhys, m_jac[i], dtmp[j], 1, dtmp[j], 1);
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);
406 boost::ignore_unused(dir, input, output, wsp);
407 NEKERROR(ErrorUtil::efatal,
"Not valid for this operator.");
417 int coll_phys_offset)
override
419 boost::ignore_unused(coll_phys_offset);
422 if (m_factors == factors)
433 "Constant factor not defined: " +
436 m_lambda = x->second;
439 m_HasVarCoeffDiff =
false;
441 if (d == factors.end())
447 for (
int i = 0; i < m_coordim; ++i)
452 for (
int i = 0; i < m_coordim; ++i)
454 d = factors.find(m_factorCoeffDef[i][i]);
455 if (d != factors.end())
457 m_diff[i][i] = d->second;
464 for (
int j = i + 1; j < m_coordim; ++j)
466 d = factors.find(m_factorCoeffDef[i][j]);
467 if (d != factors.end())
469 m_diff[i][j] = m_diff[j][i] = d->second;
473 m_HasVarCoeffDiff =
true;
497 :
Operator(pCollExp, pGeomData, factors)
500 m_dim = PtsKey.size();
501 m_coordim = pCollExp[0]->GetCoordim();
502 int nqtot = m_stdExp->GetTotPoints();
504 m_derivFac = pGeomData->GetDerivFactors(pCollExp);
505 m_jac = pGeomData->GetJac(pCollExp);
506 m_wspSize = (2 * m_coordim + 1) * nqtot;
509 m_HasVarCoeffDiff =
false;
511 this->CheckFactors(factors, 0);
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")};
566 boost::ignore_unused(output1, output2, wsp);
572 (*m_oper)(m_input, m_output);
577 (*m_oper)(input, output0);
585 boost::ignore_unused(dir, input, output, wsp);
586 NEKERROR(ErrorUtil::efatal,
"Not valid for this operator.");
593 int coll_phys_offset)
override
595 boost::ignore_unused(coll_phys_offset);
597 if (factors == m_factors)
608 "Constant factor not defined: " +
611 m_oper->SetLambda(x->second);
614 bool isConstVarDiff =
false;
616 diff[0] = diff[2] = diff[5] = 1.0;
619 if (xd00 != factors.end() && xd00->second != 1.0)
621 isConstVarDiff =
true;
622 diff[0] = xd00->second;
626 if (xd01 != factors.end() && xd01->second != 0.0)
628 isConstVarDiff =
true;
629 diff[1] = xd01->second;
633 if (xd11 != factors.end() && xd11->second != 1.0)
635 isConstVarDiff =
true;
636 diff[2] = xd11->second;
640 if (xd02 != factors.end() && xd02->second != 0.0)
642 isConstVarDiff =
true;
643 diff[3] = xd02->second;
647 if (xd12 != factors.end() && xd12->second != 0.0)
649 isConstVarDiff =
true;
650 diff[4] = xd12->second;
654 if (xd22 != factors.end() && xd22->second != 1.0)
656 isConstVarDiff =
true;
657 diff[5] = xd22->second;
662 m_oper->SetConstVarDiffusion(diff);
667 if (k != factors.end() && k->second != 0.0)
669 m_oper->SetVarDiffusion(diff);
674 std::shared_ptr<MatrixFree::Helmholtz>
m_oper;
681 :
Operator(pCollExp, pGeomData, factors),
683 pCollExp[0]->GetStdExp()->GetNcoeffs(),
687 m_nmtot = m_numElmt * pCollExp[0]->GetStdExp()->GetNcoeffs();
689 const auto dim = pCollExp[0]->GetStdExp()->GetShapeDimension();
692 std::vector<LibUtilities::BasisSharedPtr> basis(dim);
693 for (
auto i = 0; i < dim; ++i)
695 basis[i] = pCollExp[0]->GetBasis(i);
699 auto shapeType = pCollExp[0]->GetStdExp()->DetShapeType();
702 std::string op_string =
"Helmholtz";
703 op_string += MatrixFree::GetOpstring(shapeType, m_isDeformed);
705 op_string, basis, m_nElmtPad);
708 oper->SetJac(pGeomData->GetJacInterLeave(pCollExp, m_nElmtPad));
711 oper->SetDF(pGeomData->GetDerivFactorsInterLeave(pCollExp, m_nElmtPad));
713 m_oper = std::dynamic_pointer_cast<MatrixFree::Helmholtz>(oper);
714 ASSERTL0(m_oper,
"Failed to cast pointer.");
718 this->CheckFactors(factors, 0);
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"),
#define ASSERTL0(condition, msg)
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
#define OPERATOR_CREATE(cname)
Helmholtz operator using LocalRegions implementation.
StdRegions::FactorMap m_factors
Array< TwoD, const NekDouble > m_derivFac
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
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
Array< OneD, const NekDouble > m_jac
~Helmholtz_IterPerExp() final
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
StdRegions::FactorMap m_factors
void operator()(int dir, const Array< OneD, const NekDouble > &input, Array< OneD, NekDouble > &output, Array< OneD, NekDouble > &wsp) override final
~Helmholtz_MatrixFree() 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.
StdRegions::FactorMap m_factors
~Helmholtz_NoCollection() final
vector< StdRegions::StdExpansionSharedPtr > m_expList
Base class for operators on a collection of elements.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
std::tuple< LibUtilities::ShapeType, OperatorType, ImplementationType, ExpansionIsNodal > OperatorKey
Key for describing an Operator.
std::shared_ptr< CoalescedGeomData > CoalescedGeomDataSharedPtr
OperatorFactory & GetOperatorFactory()
Returns the singleton Operator factory object.
std::vector< PointsKey > PointsKeyVector
static FactorMap NullFactorMap
const char *const ConstFactorTypeMap[]
The above copyright notice and this permission notice shall be included.
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.
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
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
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.
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)