37 #include <boost/core/ignore_unused.hpp>
45 string UnsteadyAdvectionDiffusion::className
47 "UnsteadyAdvectionDiffusion",
48 UnsteadyAdvectionDiffusion::create);
50 UnsteadyAdvectionDiffusion::UnsteadyAdvectionDiffusion(
71 m_session->MatchSolverInfo(
"Extrapolation",
"SubStepping",
76 std::vector<std::string> vel;
105 m_session->LoadSolverInfo(
"AdvectionType", advName,
"WeakDG");
107 CreateInstance(advName, advName);
110 m_session->LoadSolverInfo(
"UpwindType", riemName,
"Upwind");
121 std::string diffName;
122 m_session->LoadSolverInfo(
"DiffusionType", diffName,
"LDG");
124 CreateInstance(diffName, diffName);
139 m_session->LoadSolverInfo(
"AdvectionType", advName,
142 CreateInstance(advName, advName);
146 if(advName.compare(
"WeakDG") == 0)
149 m_session->LoadSolverInfo(
"UpwindType", riemName,
"Upwind");
162 ASSERTL0(
false,
"Explicit Galerkin diffusion not set up.");
169 ASSERTL0(
false,
"Unsupported projection type.");
183 "Projection must be set to Mixed_CG_Discontinuous for "
221 for (i = 0; i < velfield.num_elements(); ++i)
223 m_fields[0]->ExtractTracePhys(velfield[i], tmp);
249 int nVariables = inarray.num_elements();
256 inarray, outarray, time);
259 for (
int i = 0; i < nVariables; ++i)
267 for (
int i = 0; i < nVariables; ++i)
274 for (
int i = 0; i < nVariables; ++i)
277 &outarray[i][0], 1, &outarray[i][0], 1);
297 int nvariables = inarray.num_elements();
306 for(i = 0; i < nvariables; ++i)
317 for(i = 0; i < nvariables; ++i)
319 m_fields[i]->FwdTrans(inarray[i], coeffs);
320 m_fields[i]->BwdTrans_IterPerExp(coeffs, outarray[i]);
326 ASSERTL0(
false,
"Unknown projection scheme");
346 int nvariables = inarray.num_elements();
363 for (
int n = 0; n < nvariables; ++n)
372 for (
int i = 0; i < nvariables; ++i)
376 inarray[i], 1, F[i], 1);
382 for (
int i = 0; i < nvariables; ++i)
403 "Dimension of flux array and velocity array do not match");
405 const int nq =
m_fields[0]->GetNpoints();
407 for (
int i = 0; i < flux.num_elements(); ++i)
409 for (
int j = 0; j < flux[0].num_elements(); ++j)
431 boost::ignore_unused(inarray);
433 unsigned int nDim = qfield.num_elements();
434 unsigned int nConvectiveFields = qfield[0].num_elements();
435 unsigned int nPts = qfield[0][0].num_elements();
436 for (
unsigned int j = 0; j < nDim; ++j)
438 for (
unsigned int i = 0; i < nConvectiveFields; ++i)
441 viscousTensor[j][i], 1 );
482 static int ncalls = 1;
500 cout <<
"Sub-integrating using "<< nsubsteps
505 for (
int m = 0; m < nint; ++m)
508 fields = integrationSoln->UpdateSolutionVector()[m];
517 for(n = 0; n < nsubsteps; ++n)
524 integrationSoln->SetSolVector(m,fields);
534 int n_element =
m_fields[0]->GetExpSize();
546 for(
int el = 0; el < n_element; ++el)
549 (stdVelocity[el] * cLambda *
550 (ExpOrder[el]-1) * (ExpOrder[el]-1));
580 ASSERTL0(0,
"Integration method not suitable: Options include BackwardEuler or BDFImplicitOrder1");
583 m_intSteps = IntegrationScheme->GetIntegrationSteps();
599 int nVariables = inarray.num_elements();
602 int ncoeffs =
m_fields[0]->GetNcoeffs();
607 for(i = 1; i < nVariables; ++i)
609 WeakAdv[i] = WeakAdv[i-1] + ncoeffs;
616 inarray, outarray, time);
618 for(i = 0; i < nVariables; ++i)
620 m_fields[i]->IProductWRTBase(outarray[i],WeakAdv[i]);
629 for(i = 0; i < nVariables; ++i)
635 m_fields[i]->MultiplyByElmtInvMass(WeakAdv[i], WeakAdv[i]);
638 m_fields[i]->BwdTrans(WeakAdv[i], outarray[i]);
650 boost::ignore_unused(time);
652 ASSERTL1(inarray.num_elements() == outarray.num_elements(),
"Inarray and outarray of different sizes ");
654 for(
int i = 0; i < inarray.num_elements(); ++i)
656 Vmath::Vcopy(inarray[i].num_elements(),inarray[i],1,outarray[i],1);
665 ASSERTL1(physfield.num_elements() == Outarray.num_elements(),
666 "Physfield and outarray are of different dimensions");
671 int nTracePts =
m_fields[0]->GetTrace()->GetNpoints();
685 for(i = 0; i < physfield.num_elements(); ++i)
689 m_fields[i]->GetFwdBwdTracePhys(physfield[i], Fwd, Bwd);
692 m_fields[0]->GetTrace()->Upwind(Vn, Fwd, Bwd, numflux);
695 Vmath::Vsub(nTracePts, numflux, 1, Fwd, 1, Fwd, 1);
696 Vmath::Vsub(nTracePts, numflux, 1, Bwd, 1, Bwd, 1);
703 m_fields[0]->AddFwdBwdTraceIntegral(Fwd,Bwd,Outarray[i]);
712 int n_points_0 =
m_fields[0]->GetExp(0)->GetTotPoints();
713 int n_element =
m_fields[0]->GetExpSize();
714 int nvel = inarray.num_elements();
717 ASSERTL0(nvel >= 2,
"Method not implemented for 1D");
726 for (
int i = 0; i < nvel; ++i)
734 for (
int el = 0; el < n_element; ++el)
736 int n_points =
m_fields[0]->GetExp(el)->GetTotPoints();
737 ptsKeys =
m_fields[0]->GetExp(el)->GetPointsKeys();
740 if(n_points != n_points_0)
742 for (
int j = 0; j < nvel; ++j)
746 n_points_0 = n_points;
750 m_fields[0]->GetExp(el)->GetGeom()->GetMetricInfo()->GetDerivFactors(ptsKeys);
752 if (
m_fields[0]->GetExp(el)->GetGeom()->GetMetricInfo()->GetGtype()
755 for (
int i = 0; i < n_points; i++)
757 stdVelocity[0][i] = gmat[0][i]*inarray[0][i+cnt]
758 + gmat[2][i]*inarray[1][i+cnt];
760 stdVelocity[1][i] = gmat[1][i]*inarray[0][i+cnt]
761 + gmat[3][i]*inarray[1][i+cnt];
766 for (
int i = 0; i < n_points; i++)
768 stdVelocity[0][i] = gmat[0][0]*inarray[0][i+cnt]
769 + gmat[2][0]*inarray[1][i+cnt];
771 stdVelocity[1][i] = gmat[1][0]*inarray[0][i+cnt]
772 + gmat[3][0]*inarray[1][i+cnt];
779 for (
int i = 0; i < n_points; i++)
781 pntVelocity = stdVelocity[0][i]*stdVelocity[0][i]
782 + stdVelocity[1][i]*stdVelocity[1][i];
784 if (pntVelocity>maxV[el])
786 maxV[el] = pntVelocity;
789 maxV[el] = sqrt(maxV[el]);
795 for (
int el = 0; el < n_element; ++el)
798 int n_points =
m_fields[0]->GetExp(el)->GetTotPoints();
799 ptsKeys =
m_fields[0]->GetExp(el)->GetPointsKeys();
802 if(n_points != n_points_0)
804 for (
int j = 0; j < nvel; ++j)
808 n_points_0 = n_points;
812 m_fields[0]->GetExp(el)->GetGeom()->GetMetricInfo()->GetDerivFactors(ptsKeys);
814 if (
m_fields[0]->GetExp(el)->GetGeom()->GetMetricInfo()->GetGtype()
817 for (
int i = 0; i < n_points; i++)
819 stdVelocity[0][i] = gmat[0][i]*inarray[0][i+cnt]
820 + gmat[3][i]*inarray[1][i+cnt]
821 + gmat[6][i]*inarray[2][i+cnt];
823 stdVelocity[1][i] = gmat[1][i]*inarray[0][i+cnt]
824 + gmat[4][i]*inarray[1][i+cnt]
825 + gmat[7][i]*inarray[2][i+cnt];
827 stdVelocity[2][i] = gmat[2][i]*inarray[0][i+cnt]
828 + gmat[5][i]*inarray[1][i+cnt]
829 + gmat[8][i]*inarray[2][i+cnt];
834 for (
int i = 0; i < n_points; i++)
836 stdVelocity[0][i] = gmat[0][0]*inarray[0][i+cnt]
837 + gmat[3][0]*inarray[1][i+cnt]
838 + gmat[6][0]*inarray[2][i+cnt];
840 stdVelocity[1][i] = gmat[1][0]*inarray[0][i+cnt]
841 + gmat[4][0]*inarray[1][i+cnt]
842 + gmat[7][0]*inarray[2][i+cnt];
844 stdVelocity[2][i] = gmat[2][0]*inarray[0][i+cnt]
845 + gmat[5][0]*inarray[1][i+cnt]
846 + gmat[8][0]*inarray[2][i+cnt];
852 for (
int i = 0; i < n_points; i++)
854 pntVelocity = stdVelocity[0][i]*stdVelocity[0][i]
855 + stdVelocity[1][i]*stdVelocity[1][i]
856 + stdVelocity[2][i]*stdVelocity[2][i];
858 if (pntVelocity > maxV[el])
860 maxV[el] = pntVelocity;
864 maxV[el] = sqrt(maxV[el]);
#define ASSERTL0(condition, msg)
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
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.
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
void DefineImplicitSolve(FuncPointerT func, ObjectPointerT obj)
A base class for PDEs which include an advection component.
SolverUtils::AdvectionSharedPtr m_advObject
Advection term.
virtual SOLVER_UTILS_EXPORT void v_InitObject()
Init object for UnsteadySystem class.
int m_spacedim
Spatial dimension (>= expansion dim).
SOLVER_UTILS_EXPORT int GetTraceNpoints()
NekDouble m_timestep
Time step size.
NekDouble m_time
Current time of simulation.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
SOLVER_UTILS_EXPORT int GetExpSize()
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array holding trace normals for DG simulations in the forwards direction.
SOLVER_UTILS_EXPORT int GetNpoints()
SOLVER_UTILS_EXPORT int GetNcoeffs()
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
SOLVER_UTILS_EXPORT void SetBoundaryConditions(NekDouble time)
Evaluates the boundary conditions at the given time.
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction(std::string name, const MultiRegions::ExpListSharedPtr &field=MultiRegions::NullExpListSharedPtr, bool cache=false)
Get a SessionFunction by name.
Base class for unsteady solvers.
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
bool m_explicitDiffusion
Indicates if explicit or implicit treatment of diffusion is used.
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &s)
Print a summary of time stepping parameters.
LibUtilities::TimeIntegrationSolutionSharedPtr m_intSoln
LibUtilities::TimeIntegrationWrapperSharedPtr m_intScheme
Wrapper to the time integration scheme.
bool m_homoInitialFwd
Flag to determine if simulation should start in homogeneous forward transformed state.
Array< OneD, NekDouble > m_traceVn
virtual ~UnsteadyAdvectionDiffusion()
Destructor.
SolverUtils::DiffusionSharedPtr m_diffusion
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Perform the projection.
Array< OneD, NekDouble > & GetNormalVelocity()
Get the normal velocity based on m_velocity.
NekDouble m_cflSafetyFactor
virtual void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Compute the RHS.
Array< OneD, Array< OneD, NekDouble > > m_velocity
virtual void DoImplicitSolve(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, NekDouble time, NekDouble lambda)
Solve implicitly the diffusion term.
NekDouble m_sVVCutoffRatio
void AddAdvectionPenaltyFlux(const Array< OneD, const Array< OneD, NekDouble > > &velfield, const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &Outarray)
Array< OneD, NekDouble > & GetNormalVel(const Array< OneD, const Array< OneD, NekDouble > > &velfield)
Get the normal velocity based on input velfield.
void GetFluxVectorAdv(const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
Evaluate the flux at each solution point for the advection part.
void SubStepAdvection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
void SubStepAdvance(const LibUtilities::TimeIntegrationSolutionSharedPtr &integrationSoln, int nstep, NekDouble time)
virtual bool v_PreIntegrate(int step)
PreIntegration step for substepping.
void SetUpSubSteppingTimeIntegration(int intMethod, const LibUtilities::TimeIntegrationWrapperSharedPtr &IntegrationScheme)
Array< OneD, NekDouble > GetMaxStdVelocity(const Array< OneD, Array< OneD, NekDouble > > inarray)
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
void GetFluxVectorDiff(const Array< OneD, Array< OneD, NekDouble > > &inarray, const Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &qfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &viscousTensor)
Evaluate the flux at each solution point for the diffusion part.
void SubStepProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
virtual void v_InitObject()
Initialise the object.
LibUtilities::TimeIntegrationWrapperSharedPtr m_subStepIntegrationScheme
virtual void v_GenerateSummary(SolverUtils::SummaryList &s)
Print Summary.
NekDouble GetSubstepTimeStep()
LibUtilities::TimeIntegrationSchemeOperators m_subStepIntegrationOps
@ eBDFImplicitOrder2
BDF multi-step scheme of order 2 (implicit)
@ eBackwardEuler
Backward Euler scheme.
@ eBDFImplicitOrder1
BDF multi-step scheme of order 1 (implicit)
std::shared_ptr< TimeIntegrationWrapper > TimeIntegrationWrapperSharedPtr
std::shared_ptr< SessionReader > SessionReaderSharedPtr
TimeIntegrationWrapperFactory & GetTimeIntegrationWrapperFactory()
std::shared_ptr< TimeIntegrationSolution > TimeIntegrationSolutionSharedPtr
std::vector< PointsKey > PointsKeyVector
std::shared_ptr< Comm > CommSharedPtr
Pointer to a Communicator object.
@ eMixed_CG_Discontinuous
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
std::vector< std::pair< std::string, std::string > > SummaryList
DiffusionFactory & GetDiffusionFactory()
EquationSystemFactory & GetEquationSystemFactory()
RiemannSolverFactory & GetRiemannSolverFactory()
@ eDeformed
Geometry is curved or has non-constant factors.
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
std::map< ConstFactorType, NekDouble > ConstFactorMap
static FlagList NullFlagList
An empty flag list.
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.
T Vmin(int n, const T *x, const int incx)
Return the minimum element in x - called vmin to avoid conflict with min.
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*y.
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.