70 size_t nDim = pFields[0]->GetCoordim(0);
71 size_t nVariable = pFields.size();
72 size_t nTracePts = pFields[0]->GetTrace()->GetTotPoints();
75 for (
int i = 0; i < nDim; ++i)
81 for (
int i = 0; i < nVariable; ++i)
96 pFields[0]->GetTrace()->GetElmtNormalLength(lengthFwd, lengthBwd);
100 pFields[0]->GetTraceMap()->GetAssemblyCommDG()->PerformExchange(lengthFwd,
104 pFields[0]->PeriodicBwdCopy(lengthFwd, lengthBwd);
112 for (
int i = 0; i < nTracePts; ++i)
117 lengthBwd[i] = lengthFwd[i];
124 Vmath::Vadd(nTracePts, lengthBwd, 1, lengthFwd, 1, lengthsum, 1);
125 Vmath::Vmul(nTracePts, lengthBwd, 1, lengthFwd, 1, lengthmul, 1);
126 Vmath::Vdiv(nTracePts, lengthsum, 1, lengthmul, 1, lengthFwd, 1);
137 for (
int i = 1; i < nVariable; ++i)
139 pFields[i]->GetBwdWeight(tmpBwdWeight, tmpBwdWeightJump);
142 Vmath::Vabs(nTracePts, tmpBwdWeight, 1, tmpBwdWeight, 1);
144 for (
int j = 0; j < nTracePts; ++j)
146 norm += tmpBwdWeight[j];
149 "different BWD for different variable not coded yet");
153 size_t nCoeffs = pFields[0]->GetNcoeffs();
155 for (
int i = 0; i < nVariable; ++i)
163 for (
int nd = 0; nd < nDim; ++nd)
167 for (
int i = 0; i < nVariable; ++i)
176 const size_t nConvectiveFields,
188 for (
int i = 0; i < nConvectiveFields; ++i)
190 fields[i]->BwdTrans(
m_wspDiff[i], outarray[i]);
197 const size_t nConvectiveFields,
207 size_t nDim = fields[0]->GetCoordim(0);
208 size_t nPts = fields[0]->GetTotPoints();
209 size_t nCoeffs = fields[0]->GetNcoeffs();
210 size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
218 for (
int j = 0; j < nDim; ++j)
222 for (
int i = 0; i < nConvectiveFields; ++i)
226 for (
int i = 0; i < nConvectiveFields; ++i)
238 for (
int i = 0; i < nConvectiveFields; ++i)
243 for (
int i = 0; i < nConvectiveFields; ++i)
245 fields[i]->GetFwdBwdTracePhys(inarray[i], vFwd[i], vBwd[i]);
250 for (
int i = 0; i < nConvectiveFields; ++i)
267 for (
int i = 0; i < nonZeroIndex.size(); ++i)
269 int j = nonZeroIndex[i];
270 for (
int k = 0; k < nDim; ++k)
272 tmpFluxIprdct[k] = elmtFlux[k][j];
274 fields[j]->IProductWRTDerivBase(tmpFluxIprdct, outarray[j]);
282 for (
int j = 0; j < nDim; ++j)
289 for (
int j = 0; j < nConvectiveFields; ++j)
299 for (
int i = 0; i < nonZeroIndex.size(); ++i)
301 int j = nonZeroIndex[i];
303 fields[j]->AddTraceIntegral(Traceflux[j], outarray[j]);
304 fields[j]->SetPhysState(
false);
308 elmtFlux, outarray, vFwd, vBwd);
311 for (
int i = 0; i < nonZeroIndex.size(); ++i)
313 int j = nonZeroIndex[i];
315 fields[j]->MultiplyByElmtInvMass(outarray[j], outarray[j]);
328 boost::ignore_unused(pFwd, pBwd);
331 size_t nDim = fields[0]->GetCoordim(0);
332 for (
int nd = 0; nd < 3; ++nd)
337 size_t nConvectiveFields = inarray.size();
338 for (
int i = 0; i < nConvectiveFields; ++i)
340 for (
int nd = 0; nd < nDim; ++nd)
342 qtmp[nd] = qfield[nd][i];
344 fields[i]->PhysDeriv(inarray[i], qtmp[0], qtmp[1], qtmp[2]);
354 size_t nDim = fields[0]->GetCoordim(0);
375 boost::ignore_unused(VolumeFlux);
378 traceflux3D[0] = TraceFlux;
380 size_t nConvectiveFields = fields.size();
392 const std::size_t nConvectiveFields,
402 size_t nDim = fields[0]->GetCoordim(0);
403 size_t nPts = fields[0]->GetTotPoints();
404 size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
406 for (
int nd = 0; nd < nDim; ++nd)
410 for (
int j = 0; j < nConvectiveFields; ++j)
417 VolumeFlux, traceSymflux, pFwd, pBwd,
421 fields, nonZeroIndex, traceSymflux,
427 const std::size_t nConvectiveFields,
437 size_t nDim = fields[0]->GetCoordim(0);
438 size_t nPts = fields[0]->GetTotPoints();
439 size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
441 for (
int nd = 0; nd < nDim; ++nd)
445 for (
int j = 0; j < nConvectiveFields; ++j)
452 VolumeFlux, traceSymflux, pFwd, pBwd,
456 fields, nonZeroIndex, traceSymflux, outarray);
461 const std::size_t nConvectiveFields,
470 boost::ignore_unused(inarray, qfield, VolumeFlux, pFwd, pBwd);
471 size_t nDim = fields[0]->GetCoordim(0);
474 nonZeroIndex, SymmFlux);
478 const std::size_t nConvectiveFields,
const size_t nDim,
484 size_t nTracePts = solution_jump[nConvectiveFields - 1].size();
489 for (
int nd = 0; nd < nDim; ++nd)
491 for (
int j = 0; j < nonZeroIndexsymm.size(); ++j)
493 int i = nonZeroIndexsymm[j];
495 traceSymflux[nd][i], 1, traceSymflux[nd][i], 1);
501 const std::size_t nConvectiveFields,
const size_t nDim,
const size_t nPts,
502 const size_t nTracePts,
508 boost::ignore_unused(nTracePts);
510 size_t nCoeffs = outarray[nConvectiveFields - 1].size();
513 for (
int i = 0; i < nDim; ++i)
518 for (
int j = 0; j < nonZeroIndex.size(); ++j)
520 nv = nonZeroIndex[j];
522 for (
int nd = 0; nd < nDim; ++nd)
526 tracelist->MultiplyByQuadratureMetric(tracflux[nd][nv],
529 fields[nv]->AddTraceQuadPhysToField(tracflux[nd][nv],
530 tracflux[nd][nv], tmpfield[nd]);
531 fields[nv]->DivideByQuadratureMetric(tmpfield[nd], tmpfield[nd]);
533 fields[nv]->IProductWRTDerivBase(tmpfield, tmpCoeff);
534 Vmath::Vadd(nCoeffs, tmpCoeff, 1, outarray[nv], 1, outarray[nv], 1);
539 const std::size_t nConvectiveFields,
const size_t nDim,
const size_t nPts,
540 const size_t nTracePts,
546 boost::ignore_unused(nTracePts);
548 size_t nCoeffs = outarray[nConvectiveFields - 1].size();
552 for (
int i = 0; i < nDim; ++i)
556 for (
int j = 0; j < nonZeroIndex.size(); ++j)
558 int nv = nonZeroIndex[j];
559 for (
int nd = 0; nd < nDim; ++nd)
563 fields[nv]->AddTraceQuadPhysToField(tracflux[nd][nv],
564 tracflux[nd][nv], tmpfield[nd]);
565 fields[nv]->DivideByQuadratureMetric(tmpfield[nd], tmpfield[nd]);
567 fields[nv]->IProductWRTDerivBase(tmpfield, tmpCoeff);
568 fields[nv]->BwdTrans(tmpCoeff, tmpPhysi);
569 Vmath::Vadd(nPts, tmpPhysi, 1, outarray[nv], 1, outarray[nv], 1);
578 std::shared_ptr<LocalRegions::ExpansionVector> traceExp =
580 size_t ntotTrac = (*traceExp).size();
581 int nTracPnt, noffset;
584 fields[0]->GetLocTraceToTraceMap();
587 locTraceToTraceMap->GetLeftRightAdjacentExpId();
589 locTraceToTraceMap->GetLeftRightAdjacentExpFlag();
591 std::shared_ptr<LocalRegions::ExpansionVector> fieldExp =
600 for (
int ntrace = 0; ntrace < ntotTrac; ++ntrace)
602 noffset = tracelist->GetPhys_Offset(ntrace);
603 nTracPnt = tracelist->GetTotPoints(ntrace);
605 factorFwdBwd[0] = 0.0;
606 factorFwdBwd[1] = 0.0;
608 for (
int nlr = 0; nlr < 2; ++nlr)
610 if (LRAdjflag[nlr][ntrace])
613 for (
int nd = 0; nd < spaceDim; nd++)
616 ->GetExp(LRAdjExpid[nlr][ntrace])
617 ->GetBasisNumModes(nd);
618 numModes = std::max(ntmp, numModes);
620 factorFwdBwd[nlr] = (numModes) * (numModes);
624 for (
int np = 0; np < nTracPnt; ++np)
626 factor[noffset + np] = std::max(factorFwdBwd[0], factorFwdBwd[1]);
632 const std::size_t nConvectiveFields,
const size_t nPts,
638 std::vector<NekDouble> vFwdTmp(nConvectiveFields),
639 vBwdTmp(nConvectiveFields), averTmp(nConvectiveFields);
640 for (
size_t p = 0;
p < nPts; ++
p)
643 for (
size_t f = 0; f < nConvectiveFields; ++f)
645 vFwdTmp[f] = vFwd[f][
p];
646 vBwdTmp[f] = vBwd[f][
p];
653 for (
size_t f = 0; f < nConvectiveFields; ++f)
655 aver[f][
p] = averTmp[f];
666 for (
size_t f = 0; f < nConvectiveFields; ++f)
668 for (
size_t p = 0;
p < nPts; ++
p)
675 jump[f][
p] = tmpF * Fweight + tmpB * Bweight;
691 size_t nDim = fields[0]->GetCoordim(0);
692 size_t nPts = fields[0]->GetTotPoints();
693 size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
694 size_t nConvectiveFields = fields.size();
696 boost::ignore_unused(inarray);
698 fields[0]->GetTraceMap();
703 for (
int nd = 0; nd < nDim; ++nd)
705 for (
int i = 0; i < nConvectiveFields; ++i)
720 if (fabs(PenaltyFactor2) > 1.0E-12)
730 for (
int nd = 0; nd < nDim; ++nd)
732 for (
int i = 0; i < nConvectiveFields; ++i)
738 fields[i]->GetFwdBwdTracePhys(qfield[nd][i], Fwd, Bwd,
true,
true,
743 for (
size_t p = 0;
p < nTracePts; ++
p)
750 TraceMap->GetAssemblyCommDG()->PerformExchange(
761 ConsVarAveJump(nConvectiveFields, nTracePts, vFwd, vBwd, solution_Aver,
768 for (
size_t p = 0;
p < nTracePts; ++
p)
773 for (
size_t f = 0; f < nConvectiveFields; ++f)
775 NekDouble jumpTmp = solution_jump[f][
p] * PenaltyFactor;
776 for (
size_t d = 0;
d < nDim; ++
d)
788 traceflux, nonZeroIndexflux,
795 const std::size_t nConvectiveFields,
const size_t nDim,
const size_t nPts,
796 const size_t nTracePts,
const NekDouble PenaltyFactor2,
804 std::vector<NekDouble> tmp(nTracePts);
807 for (
int nd1 = 0; nd1 < nDim; ++nd1)
813 for (
int nd = 0; nd < 3; ++nd)
817 for (
int nd2 = 0; nd2 < nDim; ++nd2)
819 qtmp[nd2] = elmt2ndDerv[nd2];
822 for (
int i = 0; i < nTracePts; ++i)
827 for (
int nd1 = 0; nd1 < nDim; ++nd1)
829 for (
int i = 0; i < nConvectiveFields; ++i)
831 fields[i]->PhysDeriv(qfield[nd1][i], qtmp[0], qtmp[1], qtmp[2]);
833 for (
int nd2 = nd1; nd2 < nDim; ++nd2)
836 fields[i]->GetFwdBwdTracePhys(elmt2ndDerv[nd2], Fwd, Bwd,
true,
838 for (
int p = 0;
p < nTracePts; ++
p)
844 numDerivFwd[nd1][i][
p]);
848 for (
int p = 0;
p < nTracePts; ++
p)
850 numDerivBwd[nd2][i][
p] +=
852 numDerivFwd[nd2][i][
p] =
854 numDerivFwd[nd2][i][
p]);
868 const int nConvectiveFields,
872 int nengy = nConvectiveFields - 1;
876 size_t nBndRegions = fields[nengy]->GetBndCondExpansions().size();
877 for (
size_t j = 0; j < nBndRegions; ++j)
879 if (fields[nengy]->GetBndConditions()[j]->GetBoundaryConditionType() ==
886 fields[nengy]->GetBndCondExpansions()[j]->GetExpSize();
887 for (
size_t e = 0; e < nBndEdges; ++e)
889 size_t nBndEdgePts = fields[nengy]
890 ->GetBndCondExpansions()[j]
894 std::size_t id2 = fields[0]->GetTrace()->GetPhys_Offset(
895 fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt++));
897 if (fields[0]->GetBndConditions()[j]->GetUserDefined() ==
#define ASSERTL0(condition, msg)
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
void AccumulateRegion(std::string, int iolevel=0)
Accumulate elapsed time for a region.
SOLVER_UTILS_EXPORT void DiffuseCalcDerivative(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfields, const Array< OneD, Array< OneD, NekDouble > > &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble > > &pBwd=NullNekDoubleArrayOfArray)
DiffusionFluxCons m_FunctorDiffusionfluxConsTrace
SOLVER_UTILS_EXPORT void DiffuseVolumeFlux(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfields, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, int > &nonZeroIndex=NullInt1DArray)
Diffusion Volume FLux.
SOLVER_UTILS_EXPORT void DiffuseTraceFlux(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfields, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, Array< OneD, NekDouble > > &TraceFlux, const Array< OneD, Array< OneD, NekDouble > > &pFwd=NullNekDoubleArrayOfArray, const Array< OneD, Array< OneD, NekDouble > > &pBwd=NullNekDoubleArrayOfArray, Array< OneD, int > &nonZeroIndex=NullInt1DArray)
Diffusion term Trace Flux.
DiffusionFluxCons m_FunctorDiffusionfluxCons
SpecialBndTreat m_SpecialBndTreat
DiffusionSymmFluxCons m_FunctorSymmetricfluxCons
Array< OneD, NekDouble > m_traceNormDirctnElmtLengthRecip
virtual void v_DiffuseCoeffs(const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd) override
static DiffusionSharedPtr create(std::string diffType)
Array< OneD, NekDouble > m_tracBwdWeightAver
Array< OneD, NekDouble > m_traceNormDirctnElmtLength
void ApplyFluxBndConds(const int nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, Array< OneD, Array< OneD, NekDouble > > &flux)
aplly Neuman boundary conditions on flux Currently only consider WallAdiabatic
void GetPenaltyFactor(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, Array< OneD, NekDouble > &factor)
Get IP penalty factor based on order.
void ConsVarAve(const size_t nConvectiveFields, const T &Bweight, const std::vector< T > &vFwd, const std::vector< T > &vBwd, std::vector< T > &aver)
Calculate the average of conservative variables on traces.
void AddSecondDerivToTrace(const std::size_t nConvectiveFields, const size_t nDim, const size_t nPts, const size_t nTracePts, const NekDouble PenaltyFactor2, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const TensorOfArray3D< NekDouble > &qfield, TensorOfArray3D< NekDouble > &numDerivFwd, TensorOfArray3D< NekDouble > &numDerivBwd)
Add second derivative term to trace jump (for DDG scheme)
void CalcTraceNumFlux(const NekDouble PenaltyFactor2, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, const TensorOfArray3D< NekDouble > &qfield, const Array< OneD, Array< OneD, NekDouble > > &vFwd, const Array< OneD, Array< OneD, NekDouble > > &vBwd, Array< OneD, int > &nonZeroIndexflux, TensorOfArray3D< NekDouble > &traceflux, Array< OneD, Array< OneD, NekDouble > > &solution_Aver, Array< OneD, Array< OneD, NekDouble > > &solution_jump)
Calculate numerical flux on traces.
NekDouble m_IPSymmFluxCoeff
NekDouble m_IP2ndDervCoeff
NekDouble m_IPPenaltyCoeff
virtual void v_DiffuseCalcDerivative(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfield, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd) override
Diffusion Flux, calculate the physical derivatives.
void AddDiffusionSymmFluxToPhys(const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfield, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, Array< OneD, NekDouble > > &outarray, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd)
Array< OneD, Array< OneD, NekDouble > > m_traceJump
void AddSymmFluxIntegralToCoeff(const std::size_t nConvectiveFields, const size_t nDim, const size_t nPts, const size_t nTracePts, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, const int > &nonZeroIndex, TensorOfArray3D< NekDouble > &tracflux, Array< OneD, Array< OneD, NekDouble > > &outarray)
Add symmetric flux integration to field (in coefficient space)
TensorOfArray3D< NekDouble > m_wspNumDerivBwd
Workspace for CallTraceNumFlux.
TensorOfArray3D< NekDouble > m_wspNumDerivFwd
LibUtilities::SessionReaderSharedPtr m_session
Array< OneD, Array< OneD, NekDouble > > m_traceAver
Array< OneD, NekDouble > m_tracBwdWeightJump
virtual void v_Diffuse(const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd) override
void ConsVarAveJump(const std::size_t nConvectiveFields, const size_t npnts, const Array< OneD, const Array< OneD, NekDouble > > &vFwd, const Array< OneD, const Array< OneD, NekDouble > > &vBwd, Array< OneD, Array< OneD, NekDouble > > &aver, Array< OneD, Array< OneD, NekDouble > > &jump)
virtual void v_InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields) override
virtual void v_DiffuseVolumeFlux(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfields, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, int > &nonZeroIndex) override
Diffusion Volume Flux.
void AddDiffusionSymmFluxToCoeff(const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfield, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, Array< OneD, NekDouble > > &outarray, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd)
void DiffuseTraceSymmFlux(const std::size_t nConvectiveFields, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfield, TensorOfArray3D< NekDouble > &VolumeFlux, TensorOfArray3D< NekDouble > &SymmFlux, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd, Array< OneD, int > &nonZeroIndex)
Calculate symmetric flux on traces interface.
Array< OneD, Array< OneD, NekDouble > > m_wspDiff
Workspace for v_Diffusion.
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
virtual void v_DiffuseTraceFlux(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfields, TensorOfArray3D< NekDouble > &VolumeFlux, Array< OneD, Array< OneD, NekDouble > > &TraceFlux, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd, Array< OneD, int > &nonZeroIndex) override
Diffusion term Trace Flux.
void CalcTraceSymFlux(const std::size_t nConvectiveFields, const size_t nDim, const Array< OneD, Array< OneD, NekDouble > > &solution_Aver, Array< OneD, Array< OneD, NekDouble > > &solution_jump, Array< OneD, int > &nonZeroIndexsymm, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &traceSymflux)
Calculate symmetric flux on traces.
void AddSymmFluxIntegralToPhys(const std::size_t nConvectiveFields, const size_t nDim, const size_t nPts, const size_t nTracePts, const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, const int > &nonZeroIndex, TensorOfArray3D< NekDouble > &tracflux, Array< OneD, Array< OneD, NekDouble > > &outarray)
Add symmetric flux integration to field (in physical space)
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< AssemblyMapDG > AssemblyMapDGSharedPtr
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
std::shared_ptr< LocTraceToTraceMap > LocTraceToTraceMapSharedPtr
static const NekDouble kNekMachineEpsilon
DiffusionFactory & GetDiffusionFactory()
std::vector< double > d(NPUPPER *NPUPPER)
The above copyright notice and this permission notice shall be included.
static Array< OneD, Array< OneD, NekDouble > > NullNekDoubleArrayOfArray
static Array< OneD, NekDouble > NullNekDouble1DArray
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 Vabs(int n, const T *x, const int incx, T *y, const int incy)
vabs: y = |x|
void Neg(int n, T *x, const int incx)
Negate x = -x.
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 Vdiv(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 Zero(int n, T *x, const int incx)
Zero vector.
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
scalarT< T > abs(scalarT< T > in)