UCFD_SPARSE  v1.0
Documentation
Loading...
Searching...
No Matches
blusgs.c
Go to the documentation of this file.
1
27#include <stdio.h>
28#include "blusgs.h"
29#include "flux.h"
30#include "inverse.h"
31
39void ns_serial_pre_blusgs(int neles, int nfvars, int nface, double factor, \
40 double *fnorm_vol, double *dt, double *diag, double *fjmat)
41{
42 int idx, jdx, kdx, row, col; // Element index
43 int matsize = nfvars*nfvars;
44 double fv, dti;
45 double dmat[matsize]; // Diagonal matrix at each cell
46
47 for (idx=0; idx<neles; idx++) {
48 // Initialize diagonal matrix
49 for (kdx=0; kdx<matsize; kdx++) {
50 dmat[kdx] = 0.0;
51 }
52
53 // Computes diagonal matrix based on neighbor cells
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;
60 }
61 }
62 }
63
64 // Complete implicit operator
65 dti = 1.0/(dt[idx]*factor);
66 for (kdx=0; kdx<nfvars; kdx++) {
67 dmat[(nfvars+1)*kdx] += dti;
68 }
69
70 // LU decomposition for inverse process
71 ludcmp(nfvars, dmat);
72
73 // Allocate temporal matrix to diag array
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];
77 }
78 }
79 }
80}
81
82
90void rans_serial_pre_blusgs(int neles, int nvars, int nfvars, int nface, double factor, double betast, \
91 double *fnorm_vol, double *uptsb, double *dt, double *tdiag, double *tjmat, double *dsrc)
92{
93 int idx, jdx, kdx, row, col, err; // Element index
94 int ntvars = nvars - nfvars;
95 int matsize = ntvars*ntvars;
96 double fv;
97 double tmat[matsize]; // Diagonal matrix at each cell
98 double uf[nvars], dsrcf[nvars];
99
100 for (idx=0; idx<neles; idx++) {
101 // Initialize diagonal matrix
102 for (kdx=0; kdx<matsize; kdx++)
103 tmat[kdx] = 0.0;
104
105 for (kdx=0; kdx<nvars; kdx++) {
106 uf[kdx] = uptsb[idx+neles*kdx];
107 dsrcf[kdx] = dsrc[idx+neles*kdx];
108 }
109
110 // Computes diagonal matrix based on neighbor cells
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;
117 }
118 }
119 }
120
121 // Computes Source term Jacobian
122 err = rans_source_jacobian(nvars, ntvars, betast, uf, tmat, dsrcf);
123 if (err == -1) {
124 printf("Warning:::Source term Jacobian of RANS equations does not match\n");
125 }
126
127 // Complete implicit operator
128 for (kdx=0; kdx<ntvars; kdx++) {
129 tmat[(ntvars+1)*kdx] += 1.0/(dt[idx]*factor);
130 }
131
132 // LU decomposition for inverse process
133 ludcmp(ntvars, tmat);
134
135 // Allocate temporal matrix to diag array
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];
139 }
140 }
141 }
142}
143
144
154void ns_serial_block_lower_sweep(int neles, int nfvars, int nface, \
155 int *nei_ele, int *mapping, int *unmapping, double *fnorm_vol, \
156 double *rhsb, double *dub, double *diag, double *fjmat)
157{
158 int _idx, idx, jdx, kdx, neib; // Index variables
159 int row, col;
160 double rhs[nfvars], dmat[nfvars*nfvars];
161 double val, fv;
162
163 // Lower sweep via mapping
164 for (_idx=0; _idx<neles; _idx++) {
165 idx = mapping[_idx];
166
167 // Initialize
168 for (kdx=0; kdx<nfvars; kdx++) {
169 rhs[kdx] = rhsb[idx+kdx*neles];
170 }
171
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];
175 }
176 }
177
178 // Only for neighbor cells
179 for (jdx=0; jdx<nface; jdx++) {
180 neib = nei_ele[idx+neles*jdx];
181
182 if (unmapping[neib] != _idx) {
183 fv = fnorm_vol[neles*jdx + idx];
184 // Matrix-Vector multiplication
185 for (row=0; row<nfvars; row++) {
186 val = 0.0;
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];
190 }
191 rhs[row] -= val*fv;
192 }
193 }
194 }
195
196 // Compute inverse of diagonal matrix multiplication
197 lusubst(nfvars, dmat, rhs);
198
199 // Update dub array
200 for (kdx=0; kdx<nfvars; kdx++) {
201 dub[idx+neles*kdx] = rhs[kdx];
202 }
203 }
204}
205
206
215void rans_serial_block_lower_sweep(int neles, int nvars, int nfvars, int nface, \
216 int *nei_ele, int *mapping, int *unmapping, double *fnorm_vol, \
217 double *rhsb, double *dub, double *tdiag, double *tjmat)
218{
219 int _idx, idx, jdx, kdx, neib; // Index variables
220 int row, col;
221 int ntvars = nvars - nfvars;
222 double val, fv;
223 double rhs[ntvars], dmat[ntvars*ntvars];
224
225 // Lower sweep via mapping
226 for (_idx=0; _idx<neles; _idx++) {
227 idx = mapping[_idx];
228
229 // Initialize
230 for (kdx=0; kdx<ntvars; kdx++) {
231 rhs[kdx] = rhsb[idx+(kdx+nfvars)*neles];
232 }
233
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];
237 }
238 }
239
240 // Only for neighbor cells
241 for (jdx=0; jdx<nface; jdx++) {
242 neib = nei_ele[idx+neles*jdx];
243
244 if (unmapping[neib] != _idx) {
245 fv = fnorm_vol[idx+neles*jdx];
246 // Matrix-Vector multiplication
247 for (row=0; row<ntvars; row++) {
248 val = 0.0;
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)];
252 }
253 rhs[row] -= val*fv;
254 }
255 }
256 }
257
258 // Compute inverse of diagonal matrix multiplication
259 lusubst(ntvars, dmat, rhs);
260
261 // Update dub array
262 for (kdx=0; kdx<ntvars; kdx++) {
263 dub[idx+neles*(kdx+nfvars)] = rhs[kdx];
264 }
265 }
266}
267
268
278void ns_serial_block_upper_sweep(int neles, int nfvars, int nface, \
279 int *nei_ele, int *mapping, int *unmapping, double *fnorm_vol, \
280 double *rhsb, double *dub, double *diag, double *fjmat)
281{
282 int _idx, idx, jdx, kdx, neib; // Index variables
283 int row, col;
284 double val, fv;
285 double rhs[nfvars], dmat[nfvars*nfvars];
286
287 // Upper sweep via mapping
288 for (_idx=neles-1; _idx>-1; _idx--) {
289 idx = mapping[_idx];
290
291 // Initialize
292 for (kdx=0; kdx<nfvars; kdx++) {
293 rhs[kdx] = rhsb[idx+kdx*neles];
294 }
295
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];
299 }
300 }
301
302 // Only for neighbor cells
303 for (jdx=0; jdx<nface; jdx++) {
304 neib = nei_ele[idx+neles*jdx];
305
306 if (unmapping[neib] != _idx) {
307 fv = fnorm_vol[neles*jdx + idx];
308 // Matrix-Vector multiplication
309 for (row=0; row<nfvars; row++) {
310 val = 0.0;
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];
314 }
315 rhs[row] -= val*fv;
316 }
317 }
318 }
319
320 // Compute inverse of diagonal matrix multiplication
321 lusubst(nfvars, dmat, rhs);
322
323 // Update dub array
324 for (kdx=0; kdx<nfvars; kdx++) {
325 dub[idx+neles*kdx] = rhs[kdx];
326 }
327 }
328}
329
330
339void rans_serial_block_upper_sweep(int neles, int nvars, int nfvars, int nface, \
340 int *nei_ele, int *mapping, int *unmapping, double *fnorm_vol, \
341 double *rhsb, double *dub, double *tdiag, double *tjmat)
342{
343 int _idx, idx, jdx, kdx, neib; // Index variables
344 int row, col;
345 int ntvars = nvars - nfvars;
346 double val, fv;
347 double rhs[ntvars], dmat[ntvars*ntvars];
348
349 // Lower sweep via mapping
350 for (_idx=neles-1; _idx>-1; _idx--) {
351 idx = mapping[_idx];
352
353 // Initialize
354 for (kdx=0; kdx<ntvars; kdx++) {
355 rhs[kdx] = rhsb[idx+(kdx+nfvars)*neles];
356 }
357
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];
361 }
362 }
363
364 // Only for neighbor cells
365 for (jdx=0; jdx<nface; jdx++) {
366 neib = nei_ele[idx+neles*jdx];
367
368 if (unmapping[neib] != _idx) {
369 fv = fnorm_vol[neles*jdx + idx];
370 // Matrix-Vector multiplication
371 for (row=0; row<ntvars; row++) {
372 val = 0.0;
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)];
376 }
377 rhs[row] -= val*fv;
378 }
379 }
380 }
381
382 // Compute inverse of diagonal matrix multiplication
383 lusubst(ntvars, dmat, rhs);
384
385 // Update dub array
386 for (kdx=0; kdx<ntvars; kdx++) {
387 dub[idx+neles*(kdx+nfvars)] = rhs[kdx];
388 }
389 }
390}
391
392
396void serial_update(int neles, int nvars, double *uptsb, double *dub, double *subres)
397{
398 int idx, kdx;
399
400 for (idx=0; idx<neles; idx++) {
401 for (kdx=0; kdx<nvars; kdx++) {
402 uptsb[idx+neles*kdx] += dub[idx+neles*kdx];
403
404 // Initialize dub array
405 dub[idx+neles*kdx] = 0.0;
406 }
407 // Initialize sub-residual array
408 subres[idx] = 0.0;
409 }
410}
void serial_update(int neles, int nvars, double *uptsb, double *dub, double *subres)
Updates solution array.
Definition: blusgs.c:396
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.
Definition: blusgs.c:39
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.
Definition: blusgs.c:215
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.
Definition: blusgs.c:90
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.
Definition: blusgs.c:339
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.
Definition: blusgs.c:278
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.
Definition: blusgs.c:154
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.
Definition: flux.c:101
Header file for numerical flux funtions.
void lusubst(int n, double *LU, double *b)
Forward/Backward Substitution function.
Definition: inverse.c:73
void ludcmp(int n, double *A)
LU Decomposition function.
Definition: inverse.c:34
Header file for LU Decomposition/Substitution.
#define nvars
Definition: mpi3d.c:31