Header file for serial LU-SGS method. More...
Go to the source code of this file.
Functions | |
void | ns_serial_pre_blusgs (int neles, int nfvars, int nface, double factor, double *fnorm_vol, double *dt, double *diag, double *fjmat) |
Computes Diagonal matrix for LU-SGS method. More... | |
void | rans_serial_pre_blusgs (int neles, int nvars, int nfvars, int nface, double factor, double betast, double *fnorm_vol, double *uptsb, double *dt, double *tdiag, double *tjmat, double *dsrc) |
Computes Diagonal matrix for Block LU-SGS method for RANS equations. More... | |
void | ns_serial_block_lower_sweep (int neles, int nfvars, int nface, int *nei_ele, int *mapping, int *unmapping, double *fnorm_vol, double *rhsb, double *dub, double *diag, double *fjmat) |
Lower sweep of Block LU-SGS method for Navier-Stokes equations. More... | |
void | rans_serial_block_lower_sweep (int neles, int nvars, int nfvars, int nface, int *nei_ele, int *mapping, int *unmapping, double *fnorm_vol, double *rhsb, double *dub, double *tdiag, double *tjmat) |
Lower sweep of Block LU-SGS method for RANS equations. More... | |
void | ns_serial_block_upper_sweep (int neles, int nfvars, int nface, int *nei_ele, int *mapping, int *unmapping, double *fnorm_vol, double *rhsb, double *dub, double *diag, double *fjmat) |
Upper sweep of Block LU-SGS method for Navier-Stokes equations. More... | |
void | rans_serial_block_upper_sweep (int neles, int nvars, int nfvars, int nface, int *nei_ele, int *mapping, int *unmapping, double *fnorm_vol, double *rhsb, double *dub, double *tdiag, double *tjmat) |
Upper sweep of Block LU-SGS method for RANS equations. More... | |
void | serial_update (int neles, int nvars, double *uptsb, double *dub, double *subres) |
Updates solution array. More... | |
Header file for serial LU-SGS method.
Declaration of each function used in LU-SGS method.
Definition in file blusgs.h.
void ns_serial_block_lower_sweep | ( | int | neles, |
int | nfvars, | ||
int | nface, | ||
int * | nei_ele, | ||
int * | mapping, | ||
int * | unmapping, | ||
double * | fnorm_vol, | ||
double * | rhsb, | ||
double * | dub, | ||
double * | diag, | ||
double * | fjmat | ||
) |
Lower sweep of Block LU-SGS method for Navier-Stokes equations.
neles | Number of element cells |
nfvars | Number of flux variables |
nface | Number of faces depends on element type |
nei_ele | Indices for neighbor cells [nface, neles] |
mapping | Reordered index by Reverse Cuthill-McKee algorithm [neles] |
unmapping | Original index before Reverse Cuthill-McKee algorithm [neles] |
fnorm_vol | Surface magnitude/cell volume [nface, neles] |
rhsb | Residual (RHS) array [nfvars, neles] |
dub | (Output) Difference array for update [nvars, neles] |
diag | Diagonal matrix array [nfvars, nfvars, neles] |
fjmat | Outward block operator include flux Jacobian at each cell face [nfvars, nfvars, nface, neles] |
By processing lower sweep, intermediate solution \(\Delta Q^*\) is computed. This function is used for Euler or Navier-Stokes equations, which has the same flux shape.
solution array is stored in dub
array.
fjmat
is NOT identical with ns_serial_pre_blusgs function. Definition at line 154 of file blusgs.c.
void ns_serial_block_upper_sweep | ( | int | neles, |
int | nfvars, | ||
int | nface, | ||
int * | nei_ele, | ||
int * | mapping, | ||
int * | unmapping, | ||
double * | fnorm_vol, | ||
double * | rhsb, | ||
double * | dub, | ||
double * | diag, | ||
double * | fjmat | ||
) |
Upper sweep of Block LU-SGS method for Navier-Stokes equations.
neles | Number of element cells |
nfvars | Number of flux variables |
nface | Number of faces depends on element type |
nei_ele | Indices for neighbor cells [nface, neles] |
mapping | Reordered index by Reverse Cuthill-McKee algorithm [neles] |
unmapping | Original index before Reverse Cuthill-McKee algorithm [neles] |
fnorm_vol | Surface magnitude/cell volume [nface, neles] |
rhsb | Residual (RHS) array [nfvars, neles] |
dub | (Output) Difference array for update [nvars, neles] |
diag | Diagonal matrix array [nfvars, nfvars, neles] |
fjmat | Outward block operator include flux Jacobian at each cell face [nfvars, nfvars, nface, neles] |
By processing upper sweep, next sub-iteration solution \(\Delta Q^(k+1)\) is computed. This function is used for Euler or Navier-Stokes equations, which has the same flux shape.
solution array is stored in dub
array.
fjmat
is NOT identical with ns_serial_pre_blusgs function. Definition at line 278 of file blusgs.c.
void ns_serial_pre_blusgs | ( | int | neles, |
int | nfvars, | ||
int | nface, | ||
double | factor, | ||
double * | fnorm_vol, | ||
double * | dt, | ||
double * | diag, | ||
double * | fjmat | ||
) |
Computes Diagonal matrix for LU-SGS method.
neles | Number of element cells |
nfvars | Number of flux variables |
nface | Number of faces depends on element type |
factor | Multiplication factor for diffusive margin |
fnorm_vol | Surface magnitude/cell volume [nface, neles] |
dt | Time step [neles] |
diag | (Output) Diagonal matrix for LU-SGS method [nfvars, nfvars, neles] |
fjmat | Inward block operator include flux Jacobian at each cell face [nfvars, nfvars, nface, neles] |
This function computes diagonal matrices of the implicit operator. In Block LU-SGS method, implicit operator is approximated with block operator
. Diagonal matrices is composed of block operator matrix, which size is n-by-n. n
is the number of conservative variables in Navier-Stokes equations, or the number of turbulent variables in RANS equations.
Definition at line 39 of file blusgs.c.
void rans_serial_block_lower_sweep | ( | int | neles, |
int | nvars, | ||
int | nfvars, | ||
int | nface, | ||
int * | nei_ele, | ||
int * | mapping, | ||
int * | unmapping, | ||
double * | fnorm_vol, | ||
double * | rhsb, | ||
double * | dub, | ||
double * | tdiag, | ||
double * | tjmat | ||
) |
Lower sweep of Block LU-SGS method for RANS equations.
neles | Number of element cells |
nvars | Number of conservative variables |
nfvars | Number of flux variables |
nface | Number of faces depends on element type |
nei_ele | Indices for neighbor cells [nface, neles] |
mapping | Reordered index by Reverse Cuthill-McKee algorithm [neles] |
unmapping | Original index before Reverse Cuthill-McKee algorithm [neles] |
fnorm_vol | Surface magnitude/cell volume [nface, neles] |
rhsb | Residual (RHS) array [nvars, neles] |
dub | (Output) Difference array for update [nvars, neles] |
tdiag | Turbulent diagonal matrix array [ntvars, ntvars, neles] |
tjmat | Outward block operator include flux Jacobian at each cell face [ntvars, ntvars, nface, neles] |
Lower sweep of Block LU-SGS.
This function is used for RANS equations. solution array is stored in dub
array.
tjmat
is NOT identical with ns_serial_pre_blusgs function. Definition at line 215 of file blusgs.c.
void rans_serial_block_upper_sweep | ( | int | neles, |
int | nvars, | ||
int | nfvars, | ||
int | nface, | ||
int * | nei_ele, | ||
int * | mapping, | ||
int * | unmapping, | ||
double * | fnorm_vol, | ||
double * | rhsb, | ||
double * | dub, | ||
double * | tdiag, | ||
double * | tjmat | ||
) |
Upper sweep of Block LU-SGS method for RANS equations.
neles | Number of element cells |
nvars | Number of conservative variables |
nfvars | Number of flux variables |
nface | Number of faces depends on element type |
nei_ele | Indices for neighbor cells [nface, neles] |
mapping | Reordered index by Reverse Cuthill-McKee algorithm [neles] |
unmapping | Original index before Reverse Cuthill-McKee algorithm [neles] |
fnorm_vol | Surface magnitude/cell volume [nface, neles] |
rhsb | Residual (RHS) array [nvars, neles] |
dub | (Output) Difference array for update [nvars, neles] |
tdiag | Turbulent diagonal matrix array [ntvars, ntvars, neles] |
tjmat | Outward block operator include flux Jacobian at each cell face [ntvars, ntvars, nface, neles] |
Upper sweep of Block LU-SGS.
This function is used for RANS equations. solution array is stored in dub
array.
tjmat
is NOT identical with ns_serial_pre_blusgs function. Definition at line 339 of file blusgs.c.
void rans_serial_pre_blusgs | ( | int | neles, |
int | nvars, | ||
int | nfvars, | ||
int | nface, | ||
double | factor, | ||
double | betast, | ||
double * | fnorm_vol, | ||
double * | uptsb, | ||
double * | dt, | ||
double * | tdiag, | ||
double * | tjmat, | ||
double * | dsrc | ||
) |
Computes Diagonal matrix for Block LU-SGS method for RANS equations.
neles | Number of element cells |
nvars | Number of conservative variables |
nfvars | Number of flux variables |
nface | Number of faces depends on element type |
factor | Multiplication factor for diffusive margin |
betast | Constant \(\beta^*\) for kw-SST turbulent model |
fnorm_vol | Surface magnitude/cell volume [nface, neles] |
uptsb | Current solution array [nvars, neles] |
dt | Time step [neles] |
tdiag | (Output) Turbulent diagonal matrix array [ntvars, ntvars, neles] |
tjmat | Inward block operator include flux Jacobian at each cell face [ntvars, ntvars, nface, neles] |
dsrc | Source term derivatives [nvars, neles] |
This function computes diagonal matrices of the implicit operator of RANS equations. In Block LU-SGS method, implicit operator is approximated with block operator
. Diagonal matrices is composed of block operator matrix, which size is n-by-n. n
is the number of conservative variables in Navier-Stokes equations, or the number of turbulent variables in RANS equations.
Definition at line 90 of file blusgs.c.
void serial_update | ( | int | neles, |
int | nvars, | ||
double * | uptsb, | ||
double * | dub, | ||
double * | subres | ||
) |
Updates solution array.
neles | Number of element cells |
nvars | Number of conservative variables |
uptsb | Solution array |
dub | Result of Block LU-SGS sweeps |
subres | Residual of each sub-iteration |
solution array is updated by adding \(\Delta Q\).
Definition at line 396 of file blusgs.c.