40 double *fnorm_vol,
double *dt,
double *diag,
double *fjmat)
42 int idx, jdx, kdx, row, col;
43 int matsize = nfvars*nfvars;
47 for (idx=0; idx<neles; idx++) {
49 for (kdx=0; kdx<matsize; kdx++) {
54 for (jdx=0; jdx<nface; jdx++) {
55 fv = fnorm_vol[neles*jdx + idx];
56 for (row=0; row<nfvars; row++) {
57 for (col=0; col<nfvars; col++) {
58 dmat[nfvars*row+col] \
59 += fjmat[idx+neles*jdx+nface*neles*col+nfvars*nface*neles*row]*fv;
65 dti = 1.0/(dt[idx]*factor);
66 for (kdx=0; kdx<nfvars; kdx++) {
67 dmat[(nfvars+1)*kdx] += dti;
74 for (row=0; row<nfvars; row++) {
75 for (col=0; col<nfvars; col++) {
76 diag[idx+neles*col+neles*nfvars*row] = dmat[nfvars*row+col];
91 double *fnorm_vol,
double *uptsb,
double *dt,
double *tdiag,
double *tjmat,
double *dsrc)
93 int idx, jdx, kdx, row, col, err;
94 int ntvars =
nvars - nfvars;
95 int matsize = ntvars*ntvars;
100 for (idx=0; idx<neles; idx++) {
102 for (kdx=0; kdx<matsize; kdx++)
105 for (kdx=0; kdx<
nvars; kdx++) {
106 uf[kdx] = uptsb[idx+neles*kdx];
107 dsrcf[kdx] = dsrc[idx+neles*kdx];
111 for (jdx=0; jdx<nface; jdx++) {
112 fv = fnorm_vol[neles*jdx + idx];
113 for (row=0; row<ntvars; row++) {
114 for (col=0; col<ntvars; col++) {
115 tmat[ntvars*row+col] \
116 += tjmat[idx+neles*jdx+nface*neles*col+ntvars*nface*neles*row]*fv;
124 printf(
"Warning:::Source term Jacobian of RANS equations does not match\n");
128 for (kdx=0; kdx<ntvars; kdx++) {
129 tmat[(ntvars+1)*kdx] += 1.0/(dt[idx]*factor);
136 for (row=0; row<ntvars; row++) {
137 for (col=0; col<ntvars; col++) {
138 tdiag[idx+neles*col+neles*ntvars*row] = tmat[ntvars*row+col];
155 int *nei_ele,
int *mapping,
int *unmapping,
double *fnorm_vol, \
156 double *rhsb,
double *dub,
double *diag,
double *fjmat)
158 int _idx, idx, jdx, kdx, neib;
160 double rhs[nfvars], dmat[nfvars*nfvars];
164 for (_idx=0; _idx<neles; _idx++) {
168 for (kdx=0; kdx<nfvars; kdx++) {
169 rhs[kdx] = rhsb[idx+kdx*neles];
172 for (row=0; row<nfvars; row++) {
173 for (col=0; col<nfvars; col++) {
174 dmat[col+nfvars*row] = diag[idx+neles*col+neles*nfvars*row];
179 for (jdx=0; jdx<nface; jdx++) {
180 neib = nei_ele[idx+neles*jdx];
182 if (unmapping[neib] != _idx) {
183 fv = fnorm_vol[neles*jdx + idx];
185 for (row=0; row<nfvars; row++) {
187 for (col=0; col<nfvars; col++) {
188 val += fjmat[idx+neles*jdx+nface*neles*col+nfvars*nface*neles*row] \
189 * dub[neib+neles*col];
200 for (kdx=0; kdx<nfvars; kdx++) {
201 dub[idx+neles*kdx] = rhs[kdx];
216 int *nei_ele,
int *mapping,
int *unmapping,
double *fnorm_vol, \
217 double *rhsb,
double *dub,
double *tdiag,
double *tjmat)
219 int _idx, idx, jdx, kdx, neib;
221 int ntvars =
nvars - nfvars;
223 double rhs[ntvars], dmat[ntvars*ntvars];
226 for (_idx=0; _idx<neles; _idx++) {
230 for (kdx=0; kdx<ntvars; kdx++) {
231 rhs[kdx] = rhsb[idx+(kdx+nfvars)*neles];
234 for (row=0; row<ntvars; row++) {
235 for (col=0; col<ntvars; col++) {
236 dmat[col+ntvars*row] = tdiag[idx+neles*col+neles*ntvars*row];
241 for (jdx=0; jdx<nface; jdx++) {
242 neib = nei_ele[idx+neles*jdx];
244 if (unmapping[neib] != _idx) {
245 fv = fnorm_vol[idx+neles*jdx];
247 for (row=0; row<ntvars; row++) {
249 for (col=0; col<ntvars; col++) {
250 val += tjmat[idx+neles*jdx+nface*neles*col+ntvars*nface*neles*row] \
251 * dub[neib+neles*(col+nfvars)];
262 for (kdx=0; kdx<ntvars; kdx++) {
263 dub[idx+neles*(kdx+nfvars)] = rhs[kdx];
279 int *nei_ele,
int *mapping,
int *unmapping,
double *fnorm_vol, \
280 double *rhsb,
double *dub,
double *diag,
double *fjmat)
282 int _idx, idx, jdx, kdx, neib;
285 double rhs[nfvars], dmat[nfvars*nfvars];
288 for (_idx=neles-1; _idx>-1; _idx--) {
292 for (kdx=0; kdx<nfvars; kdx++) {
293 rhs[kdx] = rhsb[idx+kdx*neles];
296 for (row=0; row<nfvars; row++) {
297 for (col=0; col<nfvars; col++) {
298 dmat[col+nfvars*row] = diag[idx+neles*col+neles*nfvars*row];
303 for (jdx=0; jdx<nface; jdx++) {
304 neib = nei_ele[idx+neles*jdx];
306 if (unmapping[neib] != _idx) {
307 fv = fnorm_vol[neles*jdx + idx];
309 for (row=0; row<nfvars; row++) {
311 for (col=0; col<nfvars; col++) {
312 val += fjmat[idx+neles*jdx+nface*neles*col+nfvars*nface*neles*row] \
313 * dub[neib+neles*col];
324 for (kdx=0; kdx<nfvars; kdx++) {
325 dub[idx+neles*kdx] = rhs[kdx];
340 int *nei_ele,
int *mapping,
int *unmapping,
double *fnorm_vol, \
341 double *rhsb,
double *dub,
double *tdiag,
double *tjmat)
343 int _idx, idx, jdx, kdx, neib;
345 int ntvars =
nvars - nfvars;
347 double rhs[ntvars], dmat[ntvars*ntvars];
350 for (_idx=neles-1; _idx>-1; _idx--) {
354 for (kdx=0; kdx<ntvars; kdx++) {
355 rhs[kdx] = rhsb[idx+(kdx+nfvars)*neles];
358 for (row=0; row<ntvars; row++) {
359 for (col=0; col<ntvars; col++) {
360 dmat[col+ntvars*row] = tdiag[idx+neles*col+neles*ntvars*row];
365 for (jdx=0; jdx<nface; jdx++) {
366 neib = nei_ele[idx+neles*jdx];
368 if (unmapping[neib] != _idx) {
369 fv = fnorm_vol[neles*jdx + idx];
371 for (row=0; row<ntvars; row++) {
373 for (col=0; col<ntvars; col++) {
374 val += tjmat[idx+neles*jdx+nface*neles*col+ntvars*nface*neles*row] \
375 * dub[neib+neles*(col+nfvars)];
386 for (kdx=0; kdx<ntvars; kdx++) {
387 dub[idx+neles*(kdx+nfvars)] = rhs[kdx];
400 for (idx=0; idx<neles; idx++) {
401 for (kdx=0; kdx<
nvars; kdx++) {
402 uptsb[idx+neles*kdx] += dub[idx+neles*kdx];
405 dub[idx+neles*kdx] = 0.0;
void serial_update(int neles, int nvars, double *uptsb, double *dub, double *subres)
Updates solution array.
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.
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.
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.
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.
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.
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.
Header file for serial LU-SGS method.
int rans_source_jacobian(int nvars, int ntvars, double betast, double *uf, double *tmat, double *dsrc)
Computes source term Jacobian matrix for RANS equations.
Header file for numerical flux funtions.
void lusubst(int n, double *LU, double *b)
Forward/Backward Substitution function.
void ludcmp(int n, double *A)
LU Decomposition function.
Header file for LU Decomposition/Substitution.