type_euler_1D_caf_no_foodie Module

  • Uses:

  • IR_Precision
  • wenoof
module~~type_euler_1d_caf_no_foodie~~UsesGraph module~type_euler_1d_caf_no_foodie type_euler_1D_caf_no_foodie wenoof wenoof wenoof->module~type_euler_1d_caf_no_foodie IR_Precision IR_Precision IR_Precision->module~type_euler_1d_caf_no_foodie
Help

Define Euler 1D (CAF enabled) field without extending FOODIE integrand.

Used By

module~~type_euler_1d_caf_no_foodie~~UsedByGraph module~type_euler_1d_caf_no_foodie type_euler_1D_caf_no_foodie program~integrate_euler_1d_caf integrate_euler_1D_caf module~type_euler_1d_caf_no_foodie->program~integrate_euler_1d_caf
Help


Variables

TypeVisibility AttributesNameInitial
real(kind=R_P), private, allocatable:: U_L(:,:)

Integrand (state) variables, left ghost cells [1:Nc,1:Ng].

real(kind=R_P), private, allocatable:: U_R(:,:)

Integrand (state) variables, right ghost cells [1:Nc,Ni-Ng+1:Ni].

real(kind=R_P), private, allocatable:: remote_U(:,:)[:]

CAF buffer for sharing remote conservative variables.

real(kind=R_P), private, allocatable:: remote_U(:,:)

CAF buffer for sharing remote conservative variables.


Derived Types

type, public :: euler_1D_caf_nf

Euler 1D (CAF enabled) PDEs system field.

Components

TypeVisibility AttributesNameInitial
character(len=:), private, allocatable:: BC_L

Left boundary condition type.

character(len=:), private, allocatable:: BC_R

Right boundary condition type.

real(kind=R_P), private :: Dx =0._R_P

Space step.

integer(kind=I_P), private :: Nc =0

Number of conservative variables, Ns+2.

integer(kind=I_P), private :: Ng =0

Number of ghost cells for boundary conditions handling.

integer(kind=I_P), private :: Ni =0

Space dimension.

integer(kind=I_P), private :: Np =0

Number of primitive variables, Ns+4.

integer(kind=I_P), private :: Ns =0

Number of initial species.

real(kind=R_P), private, allocatable:: U(:,:)

Integrand (state) variables, whole physical domain [1:Nc,1:Ni].

real(kind=R_P), private, allocatable:: cp0(:)

Specific heat cp of initial species [1:Ns].

real(kind=R_P), private, allocatable:: cv0(:)

Specific heat cv of initial species [1:Ns].

integer(kind=I_P), private :: me =0

ID of this_image().

integer(kind=I_P), private :: ord =0

Space accuracy formal order.

integer(kind=I_P), private :: we =0

Number of CAF images used.

type(weno_interpolator_upwind), private :: weno

WENO interpolator.

Type-Bound Procedures

procedure, private, pass(lhs) :: assign_integrand => euler_assign_euler

Euler = Euler.

generic, public :: assignment(=) => assign_integrand

Overloading = assignament.

procedure, private, pass(self) :: conservative2primitive

Convert conservative variables to primitive ones.

procedure, public, pass(self) :: destroy

Destroy field.

procedure, public, pass(self) :: dt => compute_dt

Compute the current time step, by means of CFL condition.

procedure, private, pass(self) :: impose_boundary_conditions

Impose boundary conditions.

procedure, public, pass(self) :: init

Init field.

procedure, public, pass(self) :: output

Extract Euler field.

procedure, private, pass(self) :: primitive2conservative

Convert primitive variables to conservative ones.

procedure, private, pass(self) :: reconstruct_interfaces_states

Reconstruct interfaces states.

procedure, private, pass(self) :: riemann_solver

Solve the Riemann Problem at cell interfaces.

procedure, private, pass(self) :: synchronize

Synchronize CAF images.

procedure, private, pass(self) :: t => dEuler_dt

Time derivative, residuals function.

TVD RK integrator.

Components

TypeVisibility AttributesNameInitial
real(kind=R_P), public, allocatable:: alph(:,:)

\(\alpha\) Butcher's coefficients.

real(kind=R_P), public, allocatable:: beta(:)

\(\beta\) Butcher's coefficients.

real(kind=R_P), public, allocatable:: gamm(:)

\(\gamma\) Butcher's coefficients.

integer(kind=I_P), public :: stages =0

Number of stages.

Type-Bound Procedures

procedure, public, pass(self) :: destroy => destroy_rk

Destroy the integrator.

procedure, public, pass(self) :: init => init_rk

Initialize (create) the integrator.

procedure, public, pass(self) :: integrate => integrate_rk

Integrate integrand field.


Functions

private elemental function E(p, r, u, g) result(energy)

Compute total specific energy (per unit of mass). $$ E = \frac{p}{{\left( {\g - 1} \right)\r }} + \frac{{u^2 }}{2} $$

Arguments

Type IntentOptional AttributesName
real(kind=R_P), intent(in) :: p

Pressure.

real(kind=R_P), intent(in) :: r

Density.

real(kind=R_P), intent(in) :: u

Module of velocity vector.

real(kind=R_P), intent(in) :: g

Specific heats ratio \(\frac{{c_p}}{{c_v}}\).

Return Value real(kind=R_P)

Total specific energy (per unit of mass).

private elemental function H(p, r, u, g) result(entalpy)

Compute total specific entalpy (per unit of mass). $$ H = \frac{{\g p}}{{\left( {\g - 1} \right)\r }} + \frac{{u^2 }}{2} $$

Arguments

Type IntentOptional AttributesName
real(kind=R_P), intent(in) :: p

Pressure.

real(kind=R_P), intent(in) :: r

Density.

real(kind=R_P), intent(in) :: u

Module of velocity vector.

real(kind=R_P), intent(in) :: g

Specific heats ratio \(\frac{{c_p}}{{c_v}}\).

Return Value real(kind=R_P)

Total specific entalpy (per unit of mass).

private elemental function a(p, r, g) result(ss)

Compute the speed of sound for an ideal calorically perfect gas.

Arguments

Type IntentOptional AttributesName
real(kind=R_P), intent(in) :: p

Pressure.

real(kind=R_P), intent(in) :: r

Density.

real(kind=R_P), intent(in) :: g

Specific heats ratio \(\frac{{c_p}}{{c_v}}\).

Return Value real(kind=R_P)

Speed of sound.

private pure function compute_dt(self, Nmax, Tmax, t, CFL) result(Dt)

Compute the current time step by means of CFL condition.

Arguments

Type IntentOptional AttributesName
class(euler_1D_caf_nf), intent(in) :: self

Euler field.

integer(kind=I_P), intent(in) :: Nmax

Maximun number of iterates.

real(kind=R_P), intent(in) :: Tmax

Maximum time (ignored if Nmax>0).

real(kind=R_P), intent(in) :: t

Time.

real(kind=R_P), intent(in) :: CFL

CFL value.

Return Value real(kind=R_P)

Time step.

private pure function conservative2primitive(self, conservative) result(primitive)

Convert conservative variables to primitive variables.

Arguments

Type IntentOptional AttributesName
class(euler_1D_caf_nf), intent(in) :: self

Euler field.

real(kind=R_P), intent(in) :: conservative(:)

Conservative variables.

Return Value real(kind=R_P) (1:self%Np)

Primitive variables.

private function dEuler_dt(self) result(dState_dt)

Time derivative of Euler field, the residuals function.

Arguments

Type IntentOptional AttributesName
class(euler_1D_caf_nf), intent(in) :: self

Euler field.

Return Value type(euler_1D_caf_nf)

Euler field time derivative.

private pure function output(self, conservative) result(state)

Output the Euler field state (primitive variables).

Arguments

Type IntentOptional AttributesName
class(euler_1D_caf_nf), intent(in) :: self

Euler field.

logical, intent(in), optional :: conservative

Output conservative variables instead of primitive.

Return Value real(kind=R_P), dimension(:,:), allocatable

Euler state vector.

private elemental function p(r, a, g) result(pressure)

Compute the pressure for an ideal calorically perfect gas.

Arguments

Type IntentOptional AttributesName
real(kind=R_P), intent(in) :: r

Density.

real(kind=R_P), intent(in) :: a

Speed of sound.

real(kind=R_P), intent(in) :: g

Specific heats ratio \(\frac{{c_p}}{{c_v}}\).

Return Value real(kind=R_P)

Pressure.

private pure function primitive2conservative(self, primitive) result(conservative)

Convert primitive variables to conservative variables.

Arguments

Type IntentOptional AttributesName
class(euler_1D_caf_nf), intent(in) :: self

Euler field.

real(kind=R_P), intent(in) :: primitive(:)

Primitive variables.

Return Value real(kind=R_P) (1:self%Nc)

Conservative variables.

private elemental function r(p, a, g) result(density)

Compute the density for an ideal calorically perfect gas.

Arguments

Type IntentOptional AttributesName
real(kind=R_P), intent(in) :: p

Pressure.

real(kind=R_P), intent(in) :: a

Speed of sound.

real(kind=R_P), intent(in) :: g

Specific heats ratio \(\frac{{c_p}}{{c_v}}\).

Return Value real(kind=R_P)

Density.


Subroutines

private pure subroutine compute_inter_states(r1, p1, u1, g1, r4, p4, u4, g4, p, S, S1, S4)

Compute inter states (23*-states) from state1 and state4.

Arguments

Type IntentOptional AttributesName
real(kind=R_P), intent(in) :: r1

Density of state 1.

real(kind=R_P), intent(in) :: p1

Pressure of state 1.

real(kind=R_P), intent(in) :: u1

Velocity of state 1.

real(kind=R_P), intent(in) :: g1

Specific heat ratio of state 1.

real(kind=R_P), intent(in) :: r4

Density of state 4.

real(kind=R_P), intent(in) :: p4

Pressure of state 4.

real(kind=R_P), intent(in) :: u4

Velocity of state 4.

real(kind=R_P), intent(in) :: g4

Specific heat ratio of state 4.

real(kind=R_P), intent(out) :: p

Pressure of the intermediate states.

real(kind=R_P), intent(out) :: S

Contact discontinuity signal velocity.

real(kind=R_P), intent(out) :: S1

Left fastest signal velocity.

real(kind=R_P), intent(out) :: S4

Right fastest signal velocity.

private subroutine destroy(self)

Destroy field.

Arguments

Type IntentOptional AttributesName
class(euler_1D_caf_nf), intent(inout) :: self

Euler field.

private elemental subroutine destroy_rk(self)

Destoy the integrator.

Arguments

Type IntentOptional AttributesName
class(tvd_runge_kutta_integrator), intent(inout) :: self

Integrator.

private subroutine euler_assign_euler(lhs, rhs)

Assign one Euler field to another.

Arguments

Type IntentOptional AttributesName
class(euler_1D_caf_nf), intent(inout) :: lhs

Left hand side.

class(euler_1D_caf_nf), intent(in) :: rhs

Right hand side.

private subroutine impose_boundary_conditions(self, primitive)

Impose boundary conditions.

Arguments

Type IntentOptional AttributesName
class(euler_1D_caf_nf), intent(in) :: self

Euler field.

real(kind=R_P), intent(inout) :: primitive(1:self%Np,1-self%Ng:self%Ni+self%Ng)

Primitive variables [1:Np,1-Ng:Ni+Ng].

private subroutine init(self, Ni, Ns, Dx, BC_L, BC_R, initial_state, cp0, cv0, me, we, ord)

Init field.

Arguments

Type IntentOptional AttributesName
class(euler_1D_caf_nf), intent(inout) :: self

Euler field.

integer(kind=I_P), intent(in) :: Ni

Space dimension (local image).

integer(kind=I_P), intent(in) :: Ns

Number of initial species.

real(kind=R_P), intent(in) :: Dx

Space step.

character(len=*), intent(in) :: BC_L

Left boundary condition type.

character(len=*), intent(in) :: BC_R

Right boundary condition type.

real(kind=R_P), intent(in) :: initial_state(:,:)

Initial state of primitive variables.

real(kind=R_P), intent(in) :: cp0(:)

Initial specific heat, constant pressure.

real(kind=R_P), intent(in) :: cv0(:)

Initial specific heat, constant volume.

integer(kind=I_P), intent(in) :: me

ID of this_image().

integer(kind=I_P), intent(in) :: we

Number of CAF images used.

integer(kind=I_P), intent(in), optional :: ord

Space accuracy formal order.

private elemental subroutine init_rk(self, stages)

Create the actual RK integrator: initialize the Butcher' table coefficients.

Arguments

Type IntentOptional AttributesName
class(tvd_runge_kutta_integrator), intent(inout) :: self

RK integrator.

integer(kind=I_P), intent(in) :: stages

Number of stages used.

private subroutine integrate_rk(self, U, stage, Dt, t)

Integrate field with explicit TVD (or SSP) Runge-Kutta scheme.

Arguments

Type IntentOptional AttributesName
class(tvd_runge_kutta_integrator), intent(in) :: self

Actual RK integrator.

class(euler_1D_caf_nf), intent(inout) :: U

Field to be integrated.

class(euler_1D_caf_nf), intent(inout) :: stage(1:)

Runge-Kutta stages [1:stages].

real(kind=R_P), intent(in) :: Dt

Time step.

real(kind=R_P), intent(in) :: t

Time.

private subroutine reconstruct_interfaces_states(self, primitive, r_primitive)

Reconstruct the interfaces states (into primitive variables formulation) by the requested order of accuracy.

Arguments

Type IntentOptional AttributesName
class(euler_1D_caf_nf), intent(in) :: self

Euler field.

real(kind=R_P), intent(in) :: primitive(1:self%Np,1-self%Ng:self%Ni+self%Ng)

Primitive variables.

real(kind=R_P), intent(inout) :: r_primitive(1:self%Np,1:2,0:self%Ni+1)

Reconstructed primitive variables.

private pure subroutine riemann_solver(self, p1, r1, u1, g1, p4, r4, u4, g4, F)

Solve the Riemann problem between the state $1$ and $4$ using the (local) Lax Friedrichs (Rusanov) solver.

Arguments

Type IntentOptional AttributesName
class(euler_1D_caf_nf), intent(in) :: self

Euler field.

real(kind=R_P), intent(in) :: p1

Pressure of state 1.

real(kind=R_P), intent(in) :: r1

Density of state 1.

real(kind=R_P), intent(in) :: u1

Velocity of state 1.

real(kind=R_P), intent(in) :: g1

Specific heats ratio of state 1.

real(kind=R_P), intent(in) :: p4

Pressure of state 4.

real(kind=R_P), intent(in) :: r4

Density of state 4.

real(kind=R_P), intent(in) :: u4

Velocity of state 4.

real(kind=R_P), intent(in) :: g4

Specific heats ratio of state 4.

real(kind=R_P), intent(out) :: F(1:self%Nc)

Resulting fluxes.

private subroutine synchronize(self)

Synchronize CAF images.

Arguments

Type IntentOptional AttributesName
class(euler_1D_caf_nf), intent(in) :: self

Euler field.