38#include <boost/algorithm/string/predicate.hpp>
63 std::size_t nDim = pFields[0]->GetCoordim(0);
64 std::size_t nTracePts = pFields[0]->GetTrace()->GetTotPoints();
67 for (std::size_t i = 0; i < nDim; ++i)
75 const std::size_t nConvectiveFields,
82 std::size_t nCoeffs = fields[0]->GetNcoeffs();
85 for (std::size_t i = 0; i < nConvectiveFields; ++i)
93 for (std::size_t i = 0; i < nConvectiveFields; ++i)
95 fields[i]->BwdTrans(tmp[i], outarray[i]);
100 const std::size_t nConvectiveFields,
107 if (fields[0]->GetGraph()->GetMovement()->GetMoveFlag())
113 std::size_t nDim = fields[0]->GetCoordim(0);
114 std::size_t nPts = fields[0]->GetTotPoints();
115 std::size_t nCoeffs = fields[0]->GetNcoeffs();
116 std::size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
119 for (std::size_t j = 0; j < nDim; ++j)
122 for (std::size_t i = 0; i < nConvectiveFields; ++i)
129 for (std::size_t i = 0; i < nConvectiveFields; ++i)
138 for (std::size_t j = 0; j < nDim; ++j)
141 for (std::size_t i = 0; i < nConvectiveFields; ++i)
151 for (std::size_t i = 0; i < nConvectiveFields; ++i)
153 for (std::size_t j = 0; j < nDim; ++j)
155 qdbase[j] = viscTensor[j][i];
157 fields[i]->IProductWRTDerivBase(qdbase, outarray[i]);
160 fields[i]->AddTraceIntegral(traceflux[i], outarray[i]);
161 fields[i]->SetPhysState(
false);
164 if (!fields[0]->GetGraph()->GetMovement()->GetMoveFlag())
167 for (std::size_t i = 0; i < nConvectiveFields; ++i)
169 fields[i]->MultiplyByElmtInvMass(outarray[i], outarray[i]);
181 std::size_t nConvectiveFields = fields.size();
182 std::size_t nDim = fields[0]->GetCoordim(0);
183 std::size_t nCoeffs = fields[0]->GetNcoeffs();
184 std::size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
188 for (std::size_t j = 0; j < nDim; ++j)
191 for (std::size_t i = 0; i < nConvectiveFields; ++i)
199 for (std::size_t j = 0; j < nDim; ++j)
201 for (std::size_t i = 0; i < nConvectiveFields; ++i)
203 fields[i]->IProductWRTDerivBase(j, inarray[i], tmp);
205 fields[i]->AddTraceIntegral(flux[j][i], tmp);
206 fields[i]->SetPhysState(
false);
207 fields[i]->MultiplyByElmtInvMass(tmp, tmp);
208 fields[i]->BwdTrans(tmp, qfield[j][i]);
242 std::size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
243 std::size_t nvariables = fields.size();
244 std::size_t nDim = fields[0]->GetCoordim(0);
253 for (std::size_t i = 0; i < nvariables; ++i)
259 fields[i]->GetFwdBwdTracePhys(ufield[i], Fwd, Bwd);
271 if (fields[0]->GetBndCondExpansions().size())
276 for (std::size_t j = 0; j < nDim; ++j)
286 const std::size_t var,
293 std::size_t nBndRegions = fields[var]->GetBndCondExpansions().size();
295 for (std::size_t i = 0; i < nBndRegions; ++i)
297 if (fields[var]->GetBndConditions()[i]->GetBoundaryConditionType() ==
304 std::size_t nBndEdges =
305 fields[var]->GetBndCondExpansions()[i]->GetExpSize();
308 for (std::size_t e = 0; e < nBndEdges; ++e)
310 std::size_t nBndEdgePts = fields[var]
311 ->GetBndCondExpansions()[i]
316 fields[var]->GetBndCondExpansions()[i]->GetPhys_Offset(e);
318 std::size_t id2 = fields[0]->GetTrace()->GetPhys_Offset(
319 fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt++));
323 fields[var]->GetBndConditions()[i]->GetUserDefined(),
326 fields[var]->GetBndConditions()[i]->GetUserDefined(),
329 fields[var]->GetBndConditions()[i]->GetUserDefined(),
332 fields[var]->GetBndConditions()[i]->GetUserDefined(),
335 fields[var]->GetBndConditions()[i]->GetUserDefined(),
338 Vmath::Vcopy(nBndEdgePts, &Fwd[id2], 1, &penaltyflux[id2], 1);
342 ->GetBndConditions()[i]
343 ->GetBoundaryConditionType() ==
348 &(fields[var]->GetBndCondExpansions()[i]->GetPhys())[id1],
349 1, &penaltyflux[id2], 1);
352 else if ((fields[var]->GetBndConditions()[i])
353 ->GetBoundaryConditionType() ==
356 Vmath::Vcopy(nBndEdgePts, &Fwd[id2], 1, &penaltyflux[id2], 1);
372 std::size_t nTracePts = fields[0]->GetTrace()->GetTotPoints();
373 std::size_t nvariables = fields.size();
374 std::size_t nDim = qfield.size();
385 for (std::size_t i = 0; i < nvariables; ++i)
388 fields[i]->GetFwdBwdTracePhys(ufield[i], Fwd, Bwd);
393 for (std::size_t j = 0; j < nDim; ++j)
396 fields[i]->GetFwdBwdTracePhys(qfield[j][i], qFwd, qBwd);
405 Vmath::Vadd(nTracePts, uterm, 1, qfluxtemp, 1, qfluxtemp, 1);
408 if (fields[0]->GetBndCondExpansions().size())
417 Vmath::Vadd(nTracePts, qfluxtemp, 1, qflux[i], 1, qflux[i], 1);
429 const std::size_t var,
const std::size_t dir,
435 std::size_t nBndRegions = fields[var]->GetBndCondExpansions().size();
438 for (std::size_t i = 0; i < nBndRegions; ++i)
440 if (fields[var]->GetBndConditions()[i]->GetBoundaryConditionType() ==
445 std::size_t nBndEdges =
446 fields[var]->GetBndCondExpansions()[i]->GetExpSize();
449 for (std::size_t e = 0; e < nBndEdges; ++e)
451 std::size_t nBndEdgePts = fields[var]
452 ->GetBndCondExpansions()[i]
457 fields[var]->GetBndCondExpansions()[i]->GetPhys_Offset(e);
459 std::size_t id2 = fields[0]->GetTrace()->GetPhys_Offset(
460 fields[0]->GetTraceMap()->GetBndCondIDToGlobalTraceID(cnt++));
464 fields[var]->GetBndConditions()[i]->GetUserDefined(),
467 fields[var]->GetBndConditions()[i]->GetUserDefined(),
470 fields[var]->GetBndConditions()[i]->GetUserDefined(),
473 fields[var]->GetBndConditions()[i]->GetUserDefined(),
476 fields[var]->GetBndConditions()[i]->GetUserDefined(),
484 ->GetBndConditions()[i]
485 ->GetBoundaryConditionType() ==
489 &qFwd[id2], 1, &penaltyflux[id2], 1);
492 else if ((fields[var]->GetBndConditions()[i])
493 ->GetBoundaryConditionType() ==
498 &(fields[var]->GetBndCondExpansions()[i]->GetPhys())[id1],
499 1, &penaltyflux[id2], 1);
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
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 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.
DiffusionFluxVecCB m_fluxVector
void v_Diffuse(const std::size_t nConvective, 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 ApplyScalarBCs(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const std::size_t var, const Array< OneD, const NekDouble > &ufield, const Array< OneD, const NekDouble > &Fwd, const Array< OneD, const NekDouble > &Bwd, Array< OneD, NekDouble > &penaltyflux)
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.
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
static DiffusionSharedPtr create(std::string diffType)
void ApplyVectorBCs(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const std::size_t var, const std::size_t dir, const Array< OneD, const NekDouble > &qfield, const Array< OneD, const NekDouble > &qFwd, const Array< OneD, const NekDouble > &qBwd, Array< OneD, NekDouble > &penaltyflux)
void v_DiffuseCalcDerivative(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &inarray, TensorOfArray3D< NekDouble > &qfields, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd) override
Diffusion Flux, calculate the physical derivatives.
void NumFluxforScalar(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &ufield, TensorOfArray3D< NekDouble > &uflux, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd)
void v_DiffuseCoeffs(const std::size_t nConvective, 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
LibUtilities::SessionReaderSharedPtr m_session
std::string m_shockCaptureType
void v_InitObject(LibUtilities::SessionReaderSharedPtr pSession, Array< OneD, MultiRegions::ExpListSharedPtr > pFields) override
void NumFluxforVector(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, Array< OneD, NekDouble > > &ufield, TensorOfArray3D< NekDouble > &qfield, Array< OneD, Array< OneD, NekDouble > > &qflux)
Build the numerical flux for the 2nd order derivatives todo: add variable coeff and h dependence to p...
NekDouble m_C11
Coefficient of penalty term.
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
DiffusionFactory & GetDiffusionFactory()
static Array< OneD, Array< OneD, NekDouble > > NullNekDoubleArrayOfArray
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 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 Zero(int n, T *x, const int incx)
Zero vector.
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
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.