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();
215 TensorOfArray3D<NekDouble> elmtFlux{nDim};
216 TensorOfArray3D<NekDouble> qfield{nDim};
218 for (
int j = 0; j < nDim; ++j)
220 qfield[j] = Array<OneD, Array<OneD, NekDouble>>{nConvectiveFields};
221 elmtFlux[j] = Array<OneD, Array<OneD, NekDouble>>{nConvectiveFields};
222 for (
int i = 0; i < nConvectiveFields; ++i)
224 qfield[j][i] = Array<OneD, NekDouble>{nPts, 0.0};
226 for (
int i = 0; i < nConvectiveFields; ++i)
228 elmtFlux[j][i] = Array<OneD, NekDouble>{nPts, 0.0};
233 Array<OneD, Array<OneD, NekDouble>> vFwd{nConvectiveFields};
234 Array<OneD, Array<OneD, NekDouble>> vBwd{nConvectiveFields};
238 for (
int i = 0; i < nConvectiveFields; ++i)
240 vFwd[i] = Array<OneD, NekDouble>{nTracePts, 0.0};
241 vBwd[i] = Array<OneD, NekDouble>{nTracePts, 0.0};
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)
258 Array<OneD, int> nonZeroIndex;
265 Array<OneD, Array<OneD, NekDouble>> tmpFluxIprdct{nDim};
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)
288 Array<OneD, Array<OneD, NekDouble>> Traceflux{nConvectiveFields};
289 for (
int j = 0; j < nConvectiveFields; ++j)
291 Traceflux[j] = Array<OneD, NekDouble>{nTracePts, 0.0};
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]);
322 const std::size_t nConvectiveFields,
331 int nDim = fields[0]->GetCoordim(0);
332 int nPts = fields[0]->GetTotPoints();
333 int nCoeffs = fields[0]->GetNcoeffs();
334 int nTracePts = fields[0]->GetTrace()->GetTotPoints();
337 for (j = 0; j < nDim; ++j)
340 for (i = 0; i < nConvectiveFields; ++i)
350 for (i = 0; i < nonZeroIndex.size(); ++i)
352 int j = nonZeroIndex[i];
353 for (
int k = 0; k < nDim; ++k)
355 tmpFluxIprdct[k] = elmtFlux[k][j];
357 fields[j]->IProductWRTDerivBase(tmpFluxIprdct, outarray[j]);
362 for (
int j = 0; j < nConvectiveFields; ++j)
371 for (j = 0; j < nDim; ++j)
376 for (i = 0; i < nonZeroIndex.size(); ++i)
378 int j = nonZeroIndex[i];
380 fields[j]->AddTraceIntegral(Traceflux[j], outarray[j]);
381 fields[j]->SetPhysState(
false);
385 elmtFlux, outarray, vFwd, vBwd);
395 boost::ignore_unused(pFwd, pBwd);
398 size_t nDim = fields[0]->GetCoordim(0);
399 for (
int nd = 0; nd < 3; ++nd)
404 size_t nConvectiveFields = inarray.size();
405 for (
int i = 0; i < nConvectiveFields; ++i)
407 for (
int nd = 0; nd < nDim; ++nd)
409 qtmp[nd] = qfield[nd][i];
411 fields[i]->PhysDeriv(inarray[i], qtmp[0], qtmp[1], qtmp[2]);
421 size_t nDim = fields[0]->GetCoordim(0);
442 boost::ignore_unused(VolumeFlux);
445 traceflux3D[0] = TraceFlux;
447 size_t nConvectiveFields = fields.size();
459 const std::size_t nConvectiveFields,
469 size_t nDim = fields[0]->GetCoordim(0);
470 size_t nPts = fields[0]->GetTotPoints();
471 size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
473 for (
int nd = 0; nd < nDim; ++nd)
477 for (
int j = 0; j < nConvectiveFields; ++j)
484 VolumeFlux, traceSymflux, pFwd, pBwd,
488 fields, nonZeroIndex, traceSymflux,
494 const std::size_t nConvectiveFields,
504 size_t nDim = fields[0]->GetCoordim(0);
505 size_t nPts = fields[0]->GetTotPoints();
506 size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
508 for (
int nd = 0; nd < nDim; ++nd)
512 for (
int j = 0; j < nConvectiveFields; ++j)
519 VolumeFlux, traceSymflux, pFwd, pBwd,
523 fields, nonZeroIndex, traceSymflux, outarray);
528 const std::size_t nConvectiveFields,
537 boost::ignore_unused(inarray, qfield, VolumeFlux, pFwd, pBwd);
538 size_t nDim = fields[0]->GetCoordim(0);
541 nonZeroIndex, SymmFlux);
545 const std::size_t nConvectiveFields,
const size_t nDim,
551 size_t nTracePts = solution_jump[nConvectiveFields - 1].size();
556 for (
int nd = 0; nd < nDim; ++nd)
558 for (
int j = 0; j < nonZeroIndexsymm.size(); ++j)
560 int i = nonZeroIndexsymm[j];
562 traceSymflux[nd][i], 1, traceSymflux[nd][i], 1);
568 const std::size_t nConvectiveFields,
const size_t nDim,
const size_t nPts,
569 const size_t nTracePts,
575 boost::ignore_unused(nTracePts);
577 size_t nCoeffs = outarray[nConvectiveFields - 1].size();
580 for (
int i = 0; i < nDim; ++i)
585 for (
int j = 0; j < nonZeroIndex.size(); ++j)
587 nv = nonZeroIndex[j];
589 for (
int nd = 0; nd < nDim; ++nd)
593 tracelist->MultiplyByQuadratureMetric(tracflux[nd][nv],
596 fields[nv]->AddTraceQuadPhysToField(tracflux[nd][nv],
597 tracflux[nd][nv], tmpfield[nd]);
598 fields[nv]->DivideByQuadratureMetric(tmpfield[nd], tmpfield[nd]);
600 fields[nv]->IProductWRTDerivBase(tmpfield, tmpCoeff);
601 Vmath::Vadd(nCoeffs, tmpCoeff, 1, outarray[nv], 1, outarray[nv], 1);
606 const std::size_t nConvectiveFields,
const size_t nDim,
const size_t nPts,
607 const size_t nTracePts,
613 boost::ignore_unused(nTracePts);
615 size_t nCoeffs = outarray[nConvectiveFields - 1].size();
619 for (
int i = 0; i < nDim; ++i)
623 for (
int j = 0; j < nonZeroIndex.size(); ++j)
625 int nv = nonZeroIndex[j];
626 for (
int nd = 0; nd < nDim; ++nd)
630 fields[nv]->AddTraceQuadPhysToField(tracflux[nd][nv],
631 tracflux[nd][nv], tmpfield[nd]);
632 fields[nv]->DivideByQuadratureMetric(tmpfield[nd], tmpfield[nd]);
634 fields[nv]->IProductWRTDerivBase(tmpfield, tmpCoeff);
635 fields[nv]->BwdTrans(tmpCoeff, tmpPhysi);
636 Vmath::Vadd(nPts, tmpPhysi, 1, outarray[nv], 1, outarray[nv], 1);
645 std::shared_ptr<LocalRegions::ExpansionVector> traceExp =
647 size_t ntotTrac = (*traceExp).size();
648 int nTracPnt, noffset;
651 fields[0]->GetLocTraceToTraceMap();
654 locTraceToTraceMap->GetLeftRightAdjacentExpId();
656 locTraceToTraceMap->GetLeftRightAdjacentExpFlag();
658 std::shared_ptr<LocalRegions::ExpansionVector> fieldExp =
667 for (
int ntrace = 0; ntrace < ntotTrac; ++ntrace)
669 noffset = tracelist->GetPhys_Offset(ntrace);
670 nTracPnt = tracelist->GetTotPoints(ntrace);
672 factorFwdBwd[0] = 0.0;
673 factorFwdBwd[1] = 0.0;
675 for (
int nlr = 0; nlr < 2; ++nlr)
677 if (LRAdjflag[nlr][ntrace])
680 for (
int nd = 0; nd < spaceDim; nd++)
683 ->GetExp(LRAdjExpid[nlr][ntrace])
684 ->GetBasisNumModes(nd);
685 numModes = std::max(ntmp, numModes);
687 factorFwdBwd[nlr] = (numModes) * (numModes);
691 for (
int np = 0; np < nTracPnt; ++np)
693 factor[noffset + np] = std::max(factorFwdBwd[0], factorFwdBwd[1]);
702 boost::ignore_unused(fields);
707 const std::size_t nConvectiveFields,
const size_t nPts,
715 std::vector<NekDouble> vFwdTmp(nConvectiveFields),
716 vBwdTmp(nConvectiveFields), averTmp(nConvectiveFields);
717 for (
size_t p = 0;
p < nPts; ++
p)
720 for (
size_t f = 0; f < nConvectiveFields; ++f)
722 vFwdTmp[f] = vFwd[f][
p];
723 vBwdTmp[f] = vBwd[f][
p];
730 for (
size_t f = 0; f < nConvectiveFields; ++f)
732 aver[f][
p] = averTmp[f];
743 for (
size_t f = 0; f < nConvectiveFields; ++f)
745 for (
size_t p = 0;
p < nPts; ++
p)
752 jump[f][
p] = tmpF * Fweight + tmpB * Bweight;
768 size_t nDim = fields[0]->GetCoordim(0);
769 size_t nPts = fields[0]->GetTotPoints();
770 size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
771 size_t nConvectiveFields = fields.size();
773 boost::ignore_unused(inarray);
775 fields[0]->GetTraceMap();
780 for (
int nd = 0; nd < nDim; ++nd)
782 for (
int i = 0; i < nConvectiveFields; ++i)
797 if (fabs(PenaltyFactor2) > 1.0E-12)
807 for (
int nd = 0; nd < nDim; ++nd)
809 for (
int i = 0; i < nConvectiveFields; ++i)
815 fields[i]->GetFwdBwdTracePhys(qfield[nd][i], Fwd, Bwd,
true,
true,
820 for (
size_t p = 0;
p < nTracePts; ++
p)
827 TraceMap->GetAssemblyCommDG()->PerformExchange(
838 ConsVarAveJump(nConvectiveFields, nTracePts, vFwd, vBwd, solution_Aver,
845 for (
size_t p = 0;
p < nTracePts; ++
p)
850 for (
size_t f = 0; f < nConvectiveFields; ++f)
852 NekDouble jumpTmp = solution_jump[f][
p] * PenaltyFactor;
853 for (
size_t d = 0; d < nDim; ++d)
865 traceflux, nonZeroIndexflux,
872 const std::size_t nConvectiveFields,
const size_t nDim,
const size_t nPts,
873 const size_t nTracePts,
const NekDouble PenaltyFactor2,
881 std::vector<NekDouble> tmp(nTracePts);
884 for (
int nd1 = 0; nd1 < nDim; ++nd1)
890 for (
int nd = 0; nd < 3; ++nd)
894 for (
int nd2 = 0; nd2 < nDim; ++nd2)
896 qtmp[nd2] = elmt2ndDerv[nd2];
899 for (
int i = 0; i < nTracePts; ++i)
904 for (
int nd1 = 0; nd1 < nDim; ++nd1)
906 for (
int i = 0; i < nConvectiveFields; ++i)
908 fields[i]->PhysDeriv(qfield[nd1][i], qtmp[0], qtmp[1], qtmp[2]);
910 for (
int nd2 = nd1; nd2 < nDim; ++nd2)
913 fields[i]->GetFwdBwdTracePhys(elmt2ndDerv[nd2], Fwd, Bwd,
true,
915 for (
int p = 0;
p < nTracePts; ++
p)
921 numDerivFwd[nd1][i][
p]);
925 for (
int p = 0;
p < nTracePts; ++
p)
927 numDerivBwd[nd2][i][
p] +=
929 numDerivFwd[nd2][i][
p] =
931 numDerivFwd[nd2][i][
p]);
945 const int nConvectiveFields,
949 int nengy = nConvectiveFields - 1;
953 size_t nBndRegions = fields[nengy]->GetBndCondExpansions().size();
954 for (
size_t j = 0; j < nBndRegions; ++j)
956 if (fields[nengy]->GetBndConditions()[j]->GetBoundaryConditionType() ==
963 fields[nengy]->GetBndCondExpansions()[j]->GetExpSize();
964 for (
size_t e = 0; e < nBndEdges; ++e)
966 size_t nBndEdgePts = fields[nengy]
967 ->GetBndCondExpansions()[j]
971 std::size_t id2 = fields[0]->GetTrace()->GetPhys_Offset(
972 fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt++));
974 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 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 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)
Add symmetric flux to field in coeff space.
DiffusionFluxCons m_FunctorDiffusionfluxConsTrace
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)
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
SOLVER_UTILS_EXPORT 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)
Get the average and jump value of conservative variables on trace.
Array< OneD, NekDouble > m_traceNormDirctnElmtLengthRecip
static DiffusionSharedPtr create(std::string diffType)
Array< OneD, NekDouble > m_tracBwdWeightAver
Array< OneD, NekDouble > m_traceNormDirctnElmtLength
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 GetPenaltyFactorConst(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, Array< OneD, NekDouble > &factor)
Get a constant IP penalty factor.
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)
NekDouble m_IPSymmFluxCoeff
NekDouble m_IP2ndDervCoeff
NekDouble m_IPPenaltyCoeff
virtual void v_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) override
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.
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)
Array< OneD, Array< OneD, NekDouble > > m_traceJump
TensorOfArray3D< NekDouble > m_wspNumDerivBwd
Workspace for CallTraceNumFlux.
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 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.
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.
TensorOfArray3D< NekDouble > m_wspNumDerivFwd
LibUtilities::SessionReaderSharedPtr m_session
virtual void v_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) override
Array< OneD, Array< OneD, NekDouble > > m_traceAver
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
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
virtual void v_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) override
void v_InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields) override
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)
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
Array< OneD, Array< OneD, NekDouble > > m_wspDiff
Workspace for v_Diffusion.
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
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.
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.
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()
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 Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
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)