Skip to content

adam_nasto_cpu_cns

ADAM, NASTO CPU Compressible-Navier-Stokes fluidyanmics application library.

Source: src/app/nasto/cpu/adam_nasto_cpu_cns.F90

Dependencies

Contents

Subroutines

compute_conservatives

Compute convervative variables from auxiliary ones.

Attributes: pure

fortran
subroutine compute_conservatives(b, i, j, k, ngc, q_aux, q)

Arguments

NameTypeIntentAttributesDescription
binteger(kind=I4P)inCell indexes.
iinteger(kind=I4P)inCell indexes.
jinteger(kind=I4P)inCell indexes.
kinteger(kind=I4P)inCell indexes.
ngcinteger(kind=I4P)inGhost cells number.
q_auxreal(kind=R8P)inAuxiliary variables.
qreal(kind=R8P)inoutConservative varibales.

Call graph

compute_conservative_fluxes

Compute convervative fluxes from auxiliary variables.

Attributes: pure

fortran
subroutine compute_conservative_fluxes(sir, b, i, j, k, ngc, q_aux, f)

Arguments

NameTypeIntentAttributesDescription
sirreal(kind=R8P)inDirectional (1=x,2=y,3=z) increment.
binteger(kind=I4P)inCell indexes.
iinteger(kind=I4P)inCell indexes.
jinteger(kind=I4P)inCell indexes.
kinteger(kind=I4P)inCell indexes.
ngcinteger(kind=I4P)inGhost cells number.
q_auxreal(kind=R8P)inAuxiliary variables.
freal(kind=R8P)inoutConservative fluxes.

Call graph

compute_max_eigenvalues

Attributes: pure

fortran
subroutine compute_max_eigenvalues(si, sir, weno_s, b, i, j, k, ngc, nv, q_aux, evmax)

Arguments

NameTypeIntentAttributesDescription
siinteger(kind=I4P)inDirectional (1=x,2=y,3=z) increment.
sirreal(kind=R8P)inDirectional (1=x,2=y,3=z) increment.
weno_sinteger(kind=I4P)inWeno stencils number/dimension.
binteger(kind=I4P)inCell indexes.
iinteger(kind=I4P)inCell indexes.
jinteger(kind=I4P)inCell indexes.
kinteger(kind=I4P)inCell indexes.
ngcinteger(kind=I4P)inGhost cells number.
nvinteger(kind=I4P)inNumber of conservative varibales.
q_auxreal(kind=R8P)inAuxiliary variables.
evmaxreal(kind=R8P)inoutMaximum eigenvalues in the big stencil.

Call graph

compute_eigenvectors

Attributes: pure

fortran
subroutine compute_eigenvectors(si, sir, b, i, j, k, ngc, nv, g, q_aux, el, er)

Arguments

NameTypeIntentAttributesDescription
siinteger(kind=I4P)inDirectional (1=x,2=y,3=z) increment.
sirreal(kind=R8P)inDirectional (1=x,2=y,3=z) increment.
binteger(kind=I4P)inCell indexes.
iinteger(kind=I4P)inCell indexes.
jinteger(kind=I4P)inCell indexes.
kinteger(kind=I4P)inCell indexes.
ngcinteger(kind=I4P)inGhost cells number.
nvinteger(kind=I4P)inNumber of conservative varibales.
greal(kind=R8P)inSpecific heats ratio.
q_auxreal(kind=R8P)inAuxiliary variables.
elreal(kind=R8P)inoutLeft and right eigenvectors.
erreal(kind=R8P)inoutLeft and right eigenvectors.

Call graph

compute_q_aux

Compute auxiliary variables.

fortran
subroutine compute_q_aux(ni, nj, nk, ngc, blocks_number, R, cv, g, q, q_aux)

Arguments

NameTypeIntentAttributesDescription
niinteger(kind=I4P)inGrid cells number in I direction.
njinteger(kind=I4P)inGrid cells number in J direction.
nkinteger(kind=I4P)inGrid cells number in K direction.
ngcinteger(kind=I4P)inGhost cells number.
blocks_numberinteger(kind=I4P)inNumber of blocks.
Rreal(kind=R8P)inFluid constant, specific heats difference.
cvreal(kind=R8P)inSpecific heat at constant volume.
greal(kind=R8P)inSpecific heats ratio.
qreal(kind=R8P)inConservative variables.
q_auxreal(kind=R8P)outAuxiliary variables.

Call graph

compute_riemann_exact

Solve the Riemann problem between the state 1 and 4 using exact Rainkine-Hugonoit jump relations.

fortran
subroutine compute_riemann_exact(si, sir, nv, q_aux1, q_aux4, g1, g4, F, lmax, ws)

Arguments

NameTypeIntentAttributesDescription
siinteger(kind=I4P)inDirectional (1=x,2=y,3=z) increment.
sirreal(kind=R8P)inDirectional (1=x,2=y,3=z) increment.
nvinteger(kind=I4P)inNumber of conservative varibales.
q_aux1real(kind=R8P)inLeft state, auxiliary variables.
q_aux4real(kind=R8P)inLeft state, auxiliary variables.
g1real(kind=R8P)inSpecific heats ratio of state 1.
g4real(kind=R8P)inSpecific heats ratio of state 4.
Freal(kind=R8P)inoutResulting fluxes.
lmaxreal(kind=R8P)outoptionalMaximum wave speed estimation.
wsreal(kind=R8P)outoptionalMaximum wave speed estimation.

Call graph

compute_riemann_exact_2

Solve the Riemann problem between the state 1 and 4 using exact Rainkine-Hugonoit jump relations.

fortran
subroutine compute_riemann_exact_2(si, sir, nv, q_aux1, q_aux4, g1, g4, F, lmax)

Arguments

NameTypeIntentAttributesDescription
siinteger(kind=I4P)inDirectional (1=x,2=y,3=z) increment.
sirreal(kind=R8P)inDirectional (1=x,2=y,3=z) increment.
nvinteger(kind=I4P)inNumber of conservative varibales.
q_aux1real(kind=R8P)inLeft state, auxiliary variables.
q_aux4real(kind=R8P)inLeft state, auxiliary variables.
g1real(kind=R8P)inSpecific heats ratio of state 1.
g4real(kind=R8P)inSpecific heats ratio of state 4.
Freal(kind=R8P)inoutResulting fluxes.
lmaxreal(kind=R8P)outoptionalMaximum wave speed estimation.

Call graph

compute_riemann_exact_3

Solve the Riemann problem between the state 1 and 4 using exact Rainkine-Hugonoit jump relations.

fortran
subroutine compute_riemann_exact_3(si, sir, nv, q_aux1, q_aux4, g1, g4, F, lmax)

Arguments

NameTypeIntentAttributesDescription
siinteger(kind=I4P)inDirectional (1=x,2=y,3=z) increment.
sirreal(kind=R8P)inDirectional (1=x,2=y,3=z) increment.
nvinteger(kind=I4P)inNumber of conservative varibales.
q_aux1real(kind=R8P)inLeft state, auxiliary variables.
q_aux4real(kind=R8P)inLeft state, auxiliary variables.
g1real(kind=R8P)inSpecific heats ratio of state 1.
g4real(kind=R8P)inSpecific heats ratio of state 4.
Freal(kind=R8P)inoutResulting fluxes.
lmaxreal(kind=R8P)outoptionalMaximum wave speed estimation.

Call graph

compute_riemann_hllc

Solve the Riemann problem between the state 1 and 4 using the HLLC solver.

fortran
subroutine compute_riemann_hllc(si, sir, nv, q_aux1, q_aux4, g1, g4, F, lmax)

Arguments

NameTypeIntentAttributesDescription
siinteger(kind=I4P)inDirectional (1=x,2=y,3=z) increment.
sirreal(kind=R8P)inDirectional (1=x,2=y,3=z) increment.
nvinteger(kind=I4P)inNumber of conservative varibales.
q_aux1real(kind=R8P)inLeft state, auxiliary variables.
q_aux4real(kind=R8P)inLeft state, auxiliary variables.
g1real(kind=R8P)inSpecific heats ratio of state 1.
g4real(kind=R8P)inSpecific heats ratio of state 4.
Freal(kind=R8P)inoutResulting fluxes.
lmaxreal(kind=R8P)outoptionalMaximum wave speed estimation.

Call graph

compute_riemann_llf

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

fortran
subroutine compute_riemann_llf(si, sir, nv, q_aux1, q_aux4, g1, g4, F, lmax)

Arguments

NameTypeIntentAttributesDescription
siinteger(kind=I4P)inDirectional (1=x,2=y,3=z) increment.
sirreal(kind=R8P)inDirectional (1=x,2=y,3=z) increment.
nvinteger(kind=I4P)inNumber of conservative varibales.
q_aux1real(kind=R8P)inLeft state, auxiliary variables.
q_aux4real(kind=R8P)inLeft state, auxiliary variables.
g1real(kind=R8P)inSpecific heats ratio of state 1.
g4real(kind=R8P)inSpecific heats ratio of state 4.
Freal(kind=R8P)inoutResulting fluxes.
lmaxreal(kind=R8P)outoptionalMaximum wave speed estimation.

Call graph

compute_riemann_ts

Solve the Riemann problem between the state 1 and 4 using the TS (two shocks) solver.

fortran
subroutine compute_riemann_ts(si, sir, nv, q_aux1, q_aux4, g1, g4, F, lmax)

Arguments

NameTypeIntentAttributesDescription
siinteger(kind=I4P)inDirectional (1=x,2=y,3=z) increment.
sirreal(kind=R8P)inDirectional (1=x,2=y,3=z) increment.
nvinteger(kind=I4P)inNumber of conservative varibales.
q_aux1real(kind=R8P)inLeft state, auxiliary variables.
q_aux4real(kind=R8P)inLeft state, auxiliary variables.
g1real(kind=R8P)inSpecific heats ratio of state 1.
g4real(kind=R8P)inSpecific heats ratio of state 4.
Freal(kind=R8P)inoutResulting fluxes.
lmaxreal(kind=R8P)outoptionalMaximum wave speed estimation.

Call graph

compute_conservatives_scalar

Compute convervative variables from auxiliary ones, scalar input.

Attributes: pure

fortran
subroutine compute_conservatives_scalar(q_aux, q)

Arguments

NameTypeIntentAttributesDescription
q_auxreal(kind=R8P)inAuxiliary variables.
qreal(kind=R8P)inoutConservative varibales.

Call graph

compute_conservative_fluxes_scalar

Compute convervative fluxes from auxiliary variables, scalar input.

Attributes: pure

fortran
subroutine compute_conservative_fluxes_scalar(sir, q_aux, f)

Arguments

NameTypeIntentAttributesDescription
sirreal(kind=R8P)inDirectional (1=x,2=y,3=z) increment.
q_auxreal(kind=R8P)inAuxiliary variables.
freal(kind=R8P)inoutConservative fluxes.

Call graph

compute_interstates_23_ts

Attributes: elemental

fortran
subroutine compute_interstates_23_ts(p1, r1, u1, a1, g1, p4, r4, u4, a4, g4, u23, p23, r2, r3, S1, S2, S3, S4)

Arguments

NameTypeIntentAttributesDescription
p1real(kind=R8P)in
r1real(kind=R8P)in
u1real(kind=R8P)in
a1real(kind=R8P)in
g1real(kind=R8P)in
p4real(kind=R8P)in
r4real(kind=R8P)in
u4real(kind=R8P)in
a4real(kind=R8P)in
g4real(kind=R8P)in
u23real(kind=R8P)out
p23real(kind=R8P)out
r2real(kind=R8P)out
r3real(kind=R8P)out
S1real(kind=R8P)out
S2real(kind=R8P)out
S3real(kind=R8P)out
S4real(kind=R8P)out

Call graph

compute_interstates_23u

Compute intermediates states knowing the value of speed (u23) of intermediates states.

Attributes: elemental

fortran
subroutine compute_interstates_23u(p1, u1, a1, g1, gm1_1, gp1_1, delta1, eta1, p4, u4, a4, g4, gm1_4, gp1_4, delta4, eta4, u23, r2, p2, a2, r3, p3, a3, S1, S2, S3, S4)

Arguments

NameTypeIntentAttributesDescription
p1real(kind=R8P)inPressure of state 1.
u1real(kind=R8P)inVelocity of state 1.
a1real(kind=R8P)inSpeed of sound of state 1.
g1real(kind=R8P)inSpecific heats ratio of state 1.
gm1_1real(kind=R8P)ing-1, g+1 of state 1.
gp1_1real(kind=R8P)ing-1, g+1 of state 1.
delta1real(kind=R8P)in(g-1)/2, 2*g/(g-1) of state 1.
eta1real(kind=R8P)in(g-1)/2, 2*g/(g-1) of state 1.
p4real(kind=R8P)inPressure of state 4.
u4real(kind=R8P)inVelocity of state 4.
a4real(kind=R8P)inSpeed of sound of state 4.
g4real(kind=R8P)inSpecific heats ratio of state 4.
gm1_4real(kind=R8P)ing-1, g+1 of state 4.
gp1_4real(kind=R8P)ing-1, g+1 of state 4.
delta4real(kind=R8P)in(g-1)/2, 2*g/(g-1) of state 4.
eta4real(kind=R8P)in(g-1)/2, 2*g/(g-1) of state 4.
u23real(kind=R8P)inVelocity of intermediate states.
r2real(kind=R8P)outDensity, pressure and speed of sound of state 2.
p2real(kind=R8P)outDensity, pressure and speed of sound of state 2.
a2real(kind=R8P)outDensity, pressure and speed of sound of state 2.
r3real(kind=R8P)outDensity, pressure and speed of sound of state 3.
p3real(kind=R8P)outDensity, pressure and speed of sound of state 3.
a3real(kind=R8P)outDensity, pressure and speed of sound of state 3.
S1real(kind=R8P)outLeft signal velocities.
S2real(kind=R8P)outLeft signal velocities.
S3real(kind=R8P)outRight signal velocities.
S4real(kind=R8P)outRight signal velocities.

Call graph

compute_roe_average

Compute Roe averaged quantities.

Attributes: pure

fortran
subroutine compute_roe_average(ngc, b, i, j, k, ip, jp, kp, g, q_aux, uu, vv, ww, h, qq, c, ci, b1, b2)

Arguments

NameTypeIntentAttributesDescription
ngcinteger(kind=I4P)inNumber of ghost cells.
binteger(kind=I4P)inLeft/right cells indexes.
iinteger(kind=I4P)inLeft/right cells indexes.
jinteger(kind=I4P)inLeft/right cells indexes.
kinteger(kind=I4P)inLeft/right cells indexes.
ipinteger(kind=I4P)inLeft/right cells indexes.
jpinteger(kind=I4P)inLeft/right cells indexes.
kpinteger(kind=I4P)inLeft/right cells indexes.
greal(kind=R8P)inSpecific heats ratio.
q_auxreal(kind=R8P)inAuxiliary variables.
uureal(kind=R8P)outRoe state average variables.
vvreal(kind=R8P)outRoe state average variables.
wwreal(kind=R8P)outRoe state average variables.
hreal(kind=R8P)outRoe state average variables.
qqreal(kind=R8P)outRoe state average variables.
creal(kind=R8P)outRoe state average variables.
cireal(kind=R8P)outRoe state average variables.
b1real(kind=R8P)outRoe state average variables.
b2real(kind=R8P)outRoe state average variables.

Call graph

compute_rarefaction

Compute an unknown state "x" from a known state "0" when the two states are separated by a rarefaction; it is assumed that the velocity of the unknown state "ux" is known. There is an input variable that indicates if the rarefaction propagates on the "u-a" direction (left) or on the "u+a" one (right).

Attributes: elemental

fortran
subroutine compute_rarefaction(sgn, g, delta, eta, u0, p0, a0, ux, rx, px, ax, s0, sx)

Arguments

NameTypeIntentAttributesDescription
sgnreal(kind=R8P)inSign for distinguishing "left" (-1) from "right" (1) wave.
greal(kind=R8P)inSpecific heats ratio.
deltareal(kind=R8P)in(g-1)/2.
etareal(kind=R8P)in2*g/(g-1).
u0real(kind=R8P)inKnown state (speed, pressure and speed of sound).
p0real(kind=R8P)inKnown state (speed, pressure and speed of sound).
a0real(kind=R8P)inKnown state (speed, pressure and speed of sound).
uxreal(kind=R8P)inKnown speed of unknown state.
rxreal(kind=R8P)outUnknown pressure and density.
pxreal(kind=R8P)outUnknown pressure and density.
axreal(kind=R8P)outUnknown pressure and density.
s0real(kind=R8P)outWave speeds (head and back fronts).
sxreal(kind=R8P)outWave speeds (head and back fronts).

Call graph

compute_shock

Compute an unknown state "x" from a known state "0" when the two states are separated by a shock; it is assumed that the velocity of the unknown state "ux" is known. There is an input variable that indicates if the shock propagates on the "u-a" direction (left) or on the "u+a" one (right).

Attributes: elemental

fortran
subroutine compute_shock(sgn, g, gm1, gp1, u0, p0, a0, ux, rx, px, ax, ss)

Arguments

NameTypeIntentAttributesDescription
sgnreal(kind=R8P)inSign for distinguishing "left" (-1) from "right" (1) wave.
greal(kind=R8P)inSpecific heats ratio.
gm1real(kind=R8P)inGm1 = g - 1 , gp1 = g + 1.
gp1real(kind=R8P)inGm1 = g - 1 , gp1 = g + 1.
u0real(kind=R8P)inKnown state (speed, pressure and speed of sound).
p0real(kind=R8P)inKnown state (speed, pressure and speed of sound).
a0real(kind=R8P)inKnown state (speed, pressure and speed of sound).
uxreal(kind=R8P)inUnknown speed.
rxreal(kind=R8P)outUnknown state (density, pressure and speed of sound).
pxreal(kind=R8P)outUnknown state (density, pressure and speed of sound).
axreal(kind=R8P)outUnknown state (density, pressure and speed of sound).
ssreal(kind=R8P)outShock wave speed.

Call graph

evaluate_waves_pvrs

Evaluate intermediate waves 2 and 3 from the known states 1,4 using the PVRS approximation.

Attributes: pure

fortran
subroutine evaluate_waves_pvrs(uni, q_aux1, q_aux4, g1, g4, S1, S, S4, lmax, p)

Arguments

NameTypeIntentAttributesDescription
uniinteger(kind=I4P)inIndex of normal velocity.
q_aux1real(kind=R8P)inLeft state, auxiliary variables.
q_aux4real(kind=R8P)inLeft state, auxiliary variables.
g1real(kind=R8P)inSpecific heats ratio of state 1.
g4real(kind=R8P)inSpecific heats ratio of state 4.
S1real(kind=R8P)outWaves speed estimations.
Sreal(kind=R8P)outWaves speed estimations.
S4real(kind=R8P)outWaves speed estimations.
lmaxreal(kind=R8P)outoptionalMaximum wave speed estimation.
preal(kind=R8P)outoptionalPressure of the intermediate states.

Call graph