Per-point calculations for getting flux. More...
#include "decs.h"
Go to the source code of this file.
Macros | |
#define | HLLCOMPUTE(cmin, cmax, U_l, U_r, F_l, F_r) (( (cmax * F_l + cmin * F_r) - (FLUXDISSIPATION)*(cmax * cmin) * (U_r - U_l) ) / (cmax + cmin) ) |
Harten-Lax-van Leer. More... | |
#define | LAXFCOMPUTE(ctop, U_l, U_r, F_l, F_r) (0.5 * ( (F_l + F_r) - (FLUXDISSIPATION)*ctop * (U_r - U_l) ) ) |
Lax-Friedrichs with Rusanov wave speed. More... | |
#define | LFCOMPUTE(crus, U_l, U_r, F_l, F_r) (0.5*( (F_l+F_r) - (FLUXDISSIPATION)*crus*(U_r - U_l))) |
Lax-Friedrichs flux. More... | |
#define | UHALF(crus, U_l, U_r, F_l, F_r) (0.5*((U_l+U_r)-(FLUXDISSIPATION)*(F_r-F_l)/crus)) |
mid-step for 2-step LW-flux More... | |
#define | LWCOMPUTE(crus, U_l, U_r, F_l, Fhalf, F_r) (Fhalf) |
2-step Lax-Wendroff flux More... | |
#define | FORCECOMPUTE(crus, U_l, U_r, F_l, Fhalf, F_r) (0.25*( (F_l+F_r) + 2.0*Fhalf - (FLUXDISSIPATION)*crus*(U_r - U_l))) |
FORCE. More... | |
#define | GFORCECOMPUTE(cg, crus, U_l, U_r, F_l, Fhalf, F_r) ( (LWCOMPUTE(crus,U_l,U_r,F_l,Fhalf,F_r) + cg*LFCOMPUTE(crus,U_l,U_r,F_l,F_r))/(1.0+cg) ) |
GFORCE // Toro & Titarev JCP 2006. More... | |
#define | HLLLAXF1(cmin, cmax, ctop, f_s, U_l, U_r, F_l, F_r) ( HLLCOMPUTE(cmin,cmax,U_l,U_r,F_l,F_r)*f_s + LAXFCOMPUTE(ctop,U_l,U_r,F_l,F_r)*(1.0-f_s) ) |
use more HLL as f_s from 0 to 1, use more LAXF as f_s from 1 to 0 More... | |
#define | MINDPOP (0.0) |
#define | MAXDPOP (0.3) |
#define | USE_CORRECTED_STATES 0 |
e.g. double rarefaction Einfelt test (test2 in Liska & Wendroff 2003) leads to F at interface corresponding to negative pressure. This is supposed to correct for that. Leads to slightly better solution except for near center More... | |
#define | MUSTACOEF (1.0) |
#define | NUMMUSTAITERS 1 |
number of single-cell musta iterations More... | |
#define | DOMULTICELL 0 |
whether to do mutli-cell version of MUSTA More... | |
#define | NUMMULTIMUSTAITERS 2 |
number of multicell musta iterations More... | |
#define | NUMLOCALCELLS 2 |
NUMMULTIMUSTAITERS=2 and NUMLOCALCELLS=1 is simplest scheme multicell much much slower than single cell version number of musta cells. More... | |
#define | MULTIMUSTACOEF (0.9) |
#define | MUSTAFORCE 0 |
#define | MUSTALAXF 1 |
#define | MUSTAHLL 2 |
#define | WHICHFLUX MUSTAHLL |
#define | HLLBOUNDMUSTA 0 |
#define | MUSTAVERSION 0 |
Functions | |
int | flux_compute_general (int i, int j, int k, int dir, struct of_geom *ptrgeom, FTYPE CUf, FTYPE *p_c, FTYPE *p_l, FTYPE *p_r, FTYPE *F, FTYPE *ctopallptr) |
actually compute the flux and ctop for dt calculation More... | |
int | p2SFUevolve (int dir, int isleftright, FTYPE *p, struct of_geom *ptrgeom, struct of_state **ptrstate, FTYPE *F, FTYPE *U) |
P->sate,flux,U. More... | |
int | p2SFUevolveall (int dir, int isleftright, FTYPE *p, struct of_geom *ptrgeom, struct of_state **ptrstate, FTYPE *FPAKE, FTYPE *FEN, FTYPE *FEM, FTYPE *FRAD, FTYPE *UPAKE, FTYPE *UEN, FTYPE *UEM, FTYPE *URAD) |
int | flux_compute (int i, int j, int k, int dir, struct of_geom *geom, FTYPE *cminmax_l, FTYPE *cminmax_r, FTYPE *cminmax, FTYPE ctopmhd, FTYPE *cminmaxrad_l, FTYPE *cminmaxrad_r, FTYPE *cminmaxrad, FTYPE ctoprad, FTYPE CUf, FTYPE *p_l, FTYPE *p_r, FTYPE *U_l, FTYPE *U_r, FTYPE *F_l, FTYPE *F_r, FTYPE *F) |
actually compute flux from given data GODMARK: right now if splitmaem==1 assumes F is linear in FMA and FEM. More... | |
Per-point calculations for getting flux.
ideas: GODMARK: Consider ZIP-average for non-diffusive part of LAXF or HLL flux Chacon (2004) Computer Physics Communications
Definition in file fluxcompute.c.
#define DOMULTICELL 0 |
whether to do mutli-cell version of MUSTA
Definition at line 711 of file fluxcompute.c.
#define FORCECOMPUTE | ( | crus, | |
U_l, | |||
U_r, | |||
F_l, | |||
Fhalf, | |||
F_r | |||
) | (0.25*( (F_l+F_r) + 2.0*Fhalf - (FLUXDISSIPATION)*crus*(U_r - U_l))) |
FORCE.
Definition at line 331 of file fluxcompute.c.
#define GFORCECOMPUTE | ( | cg, | |
crus, | |||
U_l, | |||
U_r, | |||
F_l, | |||
Fhalf, | |||
F_r | |||
) | ( (LWCOMPUTE(crus,U_l,U_r,F_l,Fhalf,F_r) + cg*LFCOMPUTE(crus,U_l,U_r,F_l,F_r))/(1.0+cg) ) |
GFORCE // Toro & Titarev JCP 2006.
Definition at line 334 of file fluxcompute.c.
#define HLLBOUNDMUSTA 0 |
Definition at line 734 of file fluxcompute.c.
#define HLLCOMPUTE | ( | cmin, | |
cmax, | |||
U_l, | |||
U_r, | |||
F_l, | |||
F_r | |||
) | (( (cmax * F_l + cmin * F_r) - (FLUXDISSIPATION)*(cmax * cmin) * (U_r - U_l) ) / (cmax + cmin) ) |
Harten-Lax-van Leer.
Definition at line 313 of file fluxcompute.c.
#define HLLLAXF1 | ( | cmin, | |
cmax, | |||
ctop, | |||
f_s, | |||
U_l, | |||
U_r, | |||
F_l, | |||
F_r | |||
) | ( HLLCOMPUTE(cmin,cmax,U_l,U_r,F_l,F_r)*f_s + LAXFCOMPUTE(ctop,U_l,U_r,F_l,F_r)*(1.0-f_s) ) |
use more HLL as f_s from 0 to 1, use more LAXF as f_s from 1 to 0
Definition at line 338 of file fluxcompute.c.
#define LAXFCOMPUTE | ( | ctop, | |
U_l, | |||
U_r, | |||
F_l, | |||
F_r | |||
) | (0.5 * ( (F_l + F_r) - (FLUXDISSIPATION)*ctop * (U_r - U_l) ) ) |
Lax-Friedrichs with Rusanov wave speed.
Definition at line 316 of file fluxcompute.c.
#define LFCOMPUTE | ( | crus, | |
U_l, | |||
U_r, | |||
F_l, | |||
F_r | |||
) | (0.5*( (F_l+F_r) - (FLUXDISSIPATION)*crus*(U_r - U_l))) |
Lax-Friedrichs flux.
Definition at line 322 of file fluxcompute.c.
#define LWCOMPUTE | ( | crus, | |
U_l, | |||
U_r, | |||
F_l, | |||
Fhalf, | |||
F_r | |||
) | (Fhalf) |
2-step Lax-Wendroff flux
Definition at line 328 of file fluxcompute.c.
#define MAXDPOP (0.3) |
#define MINDPOP (0.0) |
#define MULTIMUSTACOEF (0.9) |
Definition at line 718 of file fluxcompute.c.
#define MUSTACOEF (1.0) |
Definition at line 594 of file fluxcompute.c.
#define MUSTAFORCE 0 |
Definition at line 721 of file fluxcompute.c.
#define MUSTAHLL 2 |
Definition at line 723 of file fluxcompute.c.
#define MUSTALAXF 1 |
Definition at line 722 of file fluxcompute.c.
#define MUSTAVERSION 0 |
#define NUMLOCALCELLS 2 |
NUMMULTIMUSTAITERS=2 and NUMLOCALCELLS=1 is simplest scheme multicell much much slower than single cell version number of musta cells.
Definition at line 717 of file fluxcompute.c.
#define NUMMULTIMUSTAITERS 2 |
number of multicell musta iterations
Definition at line 713 of file fluxcompute.c.
#define NUMMUSTAITERS 1 |
number of single-cell musta iterations
Definition at line 708 of file fluxcompute.c.
#define UHALF | ( | crus, | |
U_l, | |||
U_r, | |||
F_l, | |||
F_r | |||
) | (0.5*((U_l+U_r)-(FLUXDISSIPATION)*(F_r-F_l)/crus)) |
mid-step for 2-step LW-flux
Definition at line 325 of file fluxcompute.c.
#define USE_CORRECTED_STATES 0 |
e.g. double rarefaction Einfelt test (test2 in Liska & Wendroff 2003) leads to F at interface corresponding to negative pressure. This is supposed to correct for that. Leads to slightly better solution except for near center
Definition at line 544 of file fluxcompute.c.
#define WHICHFLUX MUSTAHLL |
Definition at line 725 of file fluxcompute.c.
int flux_compute | ( | int | i, |
int | j, | ||
int | k, | ||
int | dir, | ||
struct of_geom * | geom, | ||
FTYPE * | cminmax_l, | ||
FTYPE * | cminmax_r, | ||
FTYPE * | cminmax, | ||
FTYPE | ctopmhd, | ||
FTYPE * | cminmaxrad_l, | ||
FTYPE * | cminmaxrad_r, | ||
FTYPE * | cminmaxrad, | ||
FTYPE | ctoprad, | ||
FTYPE | CUf, | ||
FTYPE * | p_l, | ||
FTYPE * | p_r, | ||
FTYPE * | U_l, | ||
FTYPE * | U_r, | ||
FTYPE * | F_l, | ||
FTYPE * | F_r, | ||
FTYPE * | F | ||
) |
actually compute flux from given data GODMARK: right now if splitmaem==1 assumes F is linear in FMA and FEM.
In general Riemann solver should return FMA and FEM in parts
Definition at line 346 of file fluxcompute.c.
int flux_compute_general | ( | int | i, |
int | j, | ||
int | k, | ||
int | dir, | ||
struct of_geom * | ptrgeom, | ||
FTYPE | CUf, | ||
FTYPE * | p_c, | ||
FTYPE * | p_l, | ||
FTYPE * | p_r, | ||
FTYPE * | F, | ||
FTYPE * | ctopallptr | ||
) |
actually compute the flux and ctop for dt calculation
actually compute flux from given data
Definition at line 25 of file fluxcompute.c.
int p2SFUevolve | ( | int | dir, |
int | isleftright, | ||
FTYPE * | p, | ||
struct of_geom * | ptrgeom, | ||
struct of_state ** | ptrstate, | ||
FTYPE * | F, | ||
FTYPE * | U | ||
) |
P->sate,flux,U.
Definition at line 199 of file fluxcompute.c.
int p2SFUevolveall | ( | int | dir, |
int | isleftright, | ||
FTYPE * | p, | ||
struct of_geom * | ptrgeom, | ||
struct of_state ** | ptrstate, | ||
FTYPE * | FPAKE, | ||
FTYPE * | FEN, | ||
FTYPE * | FEM, | ||
FTYPE * | FRAD, | ||
FTYPE * | UPAKE, | ||
FTYPE * | UEN, | ||
FTYPE * | UEM, | ||
FTYPE * | URAD | ||
) |
Definition at line 218 of file fluxcompute.c.