UCFD_SPARSE  v1.1
Documentation
Loading...
Searching...
No Matches
blusgs.c
Go to the documentation of this file.
1
26#include "blusgs.h"
27#include "flux.h"
28#include "inverse.h"
29#include <stdio.h>
30
31
40 UCFD_FLOAT *fnorm_vol, UCFD_FLOAT *dt, UCFD_FLOAT *diag, UCFD_FLOAT *fjmat)
41{
42 UCFD_INT idx, jdx, kdx, row, col;
43 UCFD_FLOAT fv, dti;
44 UCFD_FLOAT dmat[NFVARS][NFVARS];
45
46 for (idx=0; idx<neles; idx++) {
47 // Initialize diagonal matrix
48 for (row=0; row<NFVARS; row++) {
49 for (col=0; col<NFVARS; col++)
50 dmat[row][col] = 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[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[kdx][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[row][col];
77 }
78 }
79 }
80}
81
82
91 UCFD_FLOAT *fnorm_vol, UCFD_FLOAT *uptsb, UCFD_FLOAT *dt,
92 UCFD_FLOAT *tdiag, UCFD_FLOAT *tjmat, UCFD_FLOAT *dsrc)
93{
94 UCFD_INT idx, jdx, kdx, row, col;
95 UCFD_FLOAT fv;
96 UCFD_FLOAT tmat[NTURBVARS][NTURBVARS];
97 UCFD_FLOAT uf[NVARS], dsrcf[NVARS];
98
99 for (idx=0; idx<neles; idx++) {
100 // Initialize diagonal matrix
101 for (row=0; row<NTURBVARS; row++) {
102 for (col=0; col<NTURBVARS; col++)
103 tmat[row][col] = 0.0;
104 }
105
106 for (kdx=0; kdx<NVARS; kdx++) {
107 uf[kdx] = uptsb[idx+neles*kdx];
108 dsrcf[kdx] = dsrc[idx+neles*kdx];
109 }
110
111 // Computes diagonal matrix based on neighbor cells
112 for (jdx=0; jdx<nface; jdx++) {
113 fv = fnorm_vol[neles*jdx + idx];
114 for (row=0; row<NTURBVARS; row++) {
115 for (col=0; col<NTURBVARS; col++) {
116 tmat[row][col] \
117 += tjmat[idx+neles*jdx+nface*neles*col+NTURBVARS*nface*neles*row]*fv;
118 }
119 }
120 }
121
122 // Computes Source term Jacobian
123 if (rans_source_jacobian(uf, tmat, dsrcf) == UCFD_STATUS_NOT_SUPPORTED)
124 printf("Error::Invalid `NTURBVARS` value\n");
125
126 // Complete implicit operator
127 for (kdx=0; kdx<NTURBVARS; kdx++) {
128 tmat[kdx][kdx] += 1.0/(dt[idx]*factor);
129 }
130
131 // LU decomposition for inverse process
132 ludcmp(NTURBVARS, tmat);
133
134 // Allocate temporal matrix to diag array
135 for (row=0; row<NTURBVARS; row++) {
136 for (col=0; col<NTURBVARS; col++) {
137 tdiag[idx+neles*col+neles*NTURBVARS*row] = tmat[row][col];
138 }
139 }
140 }
141}
142
143
154 UCFD_INT *nei_ele, UCFD_FLOAT *fnorm_vol,
155 UCFD_FLOAT *rhsb, UCFD_FLOAT *dub, UCFD_FLOAT *diag, UCFD_FLOAT *fjmat)
156{
157 UCFD_INT idx, jdx, kdx, neib;
158 UCFD_INT row, col;
159 UCFD_FLOAT rhs[NFVARS], dmat[NFVARS][NFVARS];
160 UCFD_FLOAT val, fv;
161
162 // Lower sweep via mapping
163 for (idx=0; idx<neles; idx++) {
164 // Initialize
165 for (kdx=0; kdx<NFVARS; kdx++) {
166 rhs[kdx] = rhsb[idx+kdx*neles];
167 }
168
169 for (row=0; row<NFVARS; row++) {
170 for (col=0; col<NFVARS; col++) {
171 dmat[row][col] = diag[idx+neles*col+neles*NFVARS*row];
172 }
173 }
174
175 // Only for neighbor cells
176 for (jdx=0; jdx<nface; jdx++) {
177 neib = nei_ele[idx+neles*jdx];
178
179 if (neib != idx) {
180 fv = fnorm_vol[neles*jdx + idx];
181 // Matrix-Vector multiplication
182 for (row=0; row<NFVARS; row++) {
183 val = 0.0;
184 for (col=0; col<NFVARS; col++) {
185 val += fjmat[idx+neles*jdx+nface*neles*col+NFVARS*nface*neles*row] \
186 * dub[neib+neles*col];
187 }
188 rhs[row] -= val*fv;
189 }
190 }
191 }
192
193 // Compute inverse of diagonal matrix multiplication
194 lusub(NFVARS, dmat, rhs);
195
196 // Update dub array
197 for (kdx=0; kdx<NFVARS; kdx++) {
198 dub[idx+neles*kdx] = rhs[kdx];
199 }
200 }
201}
202
203
213 UCFD_INT *nei_ele, UCFD_FLOAT *fnorm_vol,
214 UCFD_FLOAT *rhsb, UCFD_FLOAT *dub, UCFD_FLOAT *tdiag, UCFD_FLOAT *tjmat)
215{
216 UCFD_INT idx, jdx, kdx, neib;
217 UCFD_INT row, col;
218 UCFD_FLOAT val, fv;
219 UCFD_FLOAT rhs[NTURBVARS], tmat[NTURBVARS][NTURBVARS];
220
221 // Lower sweep via mapping
222 for (idx=0; idx<neles; idx++) {
223 // Initialize
224 for (kdx=0; kdx<NTURBVARS; kdx++) {
225 rhs[kdx] = rhsb[idx+(kdx+NFVARS)*neles];
226 }
227
228 for (row=0; row<NTURBVARS; row++) {
229 for (col=0; col<NTURBVARS; col++) {
230 tmat[row][col] = tdiag[idx+neles*col+neles*NTURBVARS*row];
231 }
232 }
233
234 // Only for neighbor cells
235 for (jdx=0; jdx<nface; jdx++) {
236 neib = nei_ele[idx+neles*jdx];
237
238 if (neib != idx) {
239 fv = fnorm_vol[idx+neles*jdx];
240 // Matrix-Vector multiplication
241 for (row=0; row<NTURBVARS; row++) {
242 val = 0.0;
243 for (col=0; col<NTURBVARS; col++) {
244 val += tjmat[idx+neles*jdx+nface*neles*col+NTURBVARS*nface*neles*row] \
245 * dub[neib+neles*(col+NFVARS)];
246 }
247 rhs[row] -= val*fv;
248 }
249 }
250 }
251
252 // Compute inverse of diagonal matrix multiplication
253 lusub(NTURBVARS, tmat, rhs);
254
255 // Update dub array
256 for (kdx=0; kdx<NTURBVARS; kdx++) {
257 dub[idx+neles*(kdx+NFVARS)] = rhs[kdx];
258 }
259 }
260}
261
262
273 UCFD_INT *nei_ele, UCFD_FLOAT *fnorm_vol,
274 UCFD_FLOAT *rhsb, UCFD_FLOAT *dub, UCFD_FLOAT *diag, UCFD_FLOAT *fjmat)
275{
276 UCFD_INT idx, jdx, kdx, neib;
277 UCFD_INT row, col;
278 UCFD_FLOAT val, fv;
279 UCFD_FLOAT rhs[NFVARS], dmat[NFVARS][NFVARS];
280
281 // Upper sweep via mapping
282 for (idx=neles-1; idx>-1; idx--) {
283 // Initialize
284 for (kdx=0; kdx<NFVARS; kdx++) {
285 rhs[kdx] = rhsb[idx+kdx*neles];
286 }
287
288 for (row=0; row<NFVARS; row++) {
289 for (col=0; col<NFVARS; col++) {
290 dmat[row][col] = diag[idx+neles*col+neles*NFVARS*row];
291 }
292 }
293
294 // Only for neighbor cells
295 for (jdx=0; jdx<nface; jdx++) {
296 neib = nei_ele[idx+neles*jdx];
297
298 if (neib != idx) {
299 fv = fnorm_vol[neles*jdx + idx];
300 // Matrix-Vector multiplication
301 for (row=0; row<NFVARS; row++) {
302 val = 0.0;
303 for (col=0; col<NFVARS; col++) {
304 val += fjmat[idx+neles*jdx+nface*neles*col+NFVARS*nface*neles*row] \
305 * dub[neib+neles*col];
306 }
307 rhs[row] -= val*fv;
308 }
309 }
310 }
311
312 // Compute inverse of diagonal matrix multiplication
313 lusub(NFVARS, dmat, rhs);
314
315 // Update dub array
316 for (kdx=0; kdx<NFVARS; kdx++) {
317 dub[idx+neles*kdx] = rhs[kdx];
318 }
319 }
320}
321
322
332 UCFD_INT *nei_ele, UCFD_FLOAT *fnorm_vol,
333 UCFD_FLOAT *rhsb, UCFD_FLOAT *dub, UCFD_FLOAT *tdiag, UCFD_FLOAT *tjmat)
334{
335 UCFD_INT idx, jdx, kdx, neib;
336 UCFD_INT row, col;
337 UCFD_FLOAT val, fv;
338 UCFD_FLOAT rhs[NTURBVARS], tmat[NTURBVARS][NTURBVARS];
339
340 // Lower sweep via mapping
341 for (idx=neles-1; idx>-1; idx--) {
342
343 // Initialize
344 for (kdx=0; kdx<NTURBVARS; kdx++) {
345 rhs[kdx] = rhsb[idx+(kdx+NFVARS)*neles];
346 }
347
348 for (row=0; row<NTURBVARS; row++) {
349 for (col=0; col<NTURBVARS; col++) {
350 tmat[row][col] = tdiag[idx+neles*col+neles*NTURBVARS*row];
351 }
352 }
353
354 // Only for neighbor cells
355 for (jdx=0; jdx<nface; jdx++) {
356 neib = nei_ele[idx+neles*jdx];
357
358 if (neib != idx) {
359 fv = fnorm_vol[neles*jdx + idx];
360 // Matrix-Vector multiplication
361 for (row=0; row<NTURBVARS; row++) {
362 val = 0.0;
363 for (col=0; col<NTURBVARS; col++) {
364 val += tjmat[idx+neles*jdx+nface*neles*col+NTURBVARS*nface*neles*row] \
365 * dub[neib+neles*(col+NFVARS)];
366 }
367 rhs[row] -= val*fv;
368 }
369 }
370 }
371
372 // Compute inverse of diagonal matrix multiplication
373 lusub(NTURBVARS, tmat, rhs);
374
375 // Update dub array
376 for (kdx=0; kdx<NTURBVARS; kdx++) {
377 dub[idx+neles*(kdx+NFVARS)] = rhs[kdx];
378 }
379 }
380}
381
382
387{
388 UCFD_INT idx, kdx;
389
390 for (idx=0; idx<neles; idx++) {
391 for (kdx=0; kdx<NFVARS; kdx++) {
392 uptsb[idx+neles*kdx] += dub[idx+neles*kdx];
393
394 // Initialize dub array
395 dub[idx+neles*kdx] = 0.0;
396 }
397 // Initialize sub-residual array
398 subres[idx] = 0.0;
399 }
400}
401
402
407{
408 UCFD_INT idx, kdx;
409
410 for (idx=0; idx<neles; idx++) {
411 for (kdx=0; kdx<NVARS; kdx++) {
412 uptsb[idx+neles*kdx] += dub[idx+neles*kdx];
413
414 // Initialize dub array
415 dub[idx+neles*kdx] = 0.0;
416 }
417 // Initialize sub-residual array
418 subres[idx] = 0.0;
419 }
420}
void blusgs_serial_update(UCFD_INT neles, UCFD_FLOAT *uptsb, UCFD_FLOAT *dub, UCFD_FLOAT *subres)
Updates solution array.
Definition: blusgs.c:406
void ns_serial_pre_blusgs(UCFD_INT neles, UCFD_INT nface, UCFD_FLOAT factor, UCFD_FLOAT *fnorm_vol, UCFD_FLOAT *dt, UCFD_FLOAT *diag, UCFD_FLOAT *fjmat)
Computes Diagonal matrix for LU-SGS method.
Definition: blusgs.c:39
void ns_serial_block_lower_sweep(UCFD_INT neles, UCFD_INT nface, UCFD_INT *nei_ele, UCFD_FLOAT *fnorm_vol, UCFD_FLOAT *rhsb, UCFD_FLOAT *dub, UCFD_FLOAT *diag, UCFD_FLOAT *fjmat)
Lower sweep of Block LU-SGS method for Navier-Stokes equations.
Definition: blusgs.c:153
void rans_serial_block_upper_sweep(UCFD_INT neles, UCFD_INT nface, UCFD_INT *nei_ele, UCFD_FLOAT *fnorm_vol, UCFD_FLOAT *rhsb, UCFD_FLOAT *dub, UCFD_FLOAT *tdiag, UCFD_FLOAT *tjmat)
Upper sweep of Block LU-SGS method for RANS equations.
Definition: blusgs.c:331
void ns_serial_block_upper_sweep(UCFD_INT neles, UCFD_INT nface, UCFD_INT *nei_ele, UCFD_FLOAT *fnorm_vol, UCFD_FLOAT *rhsb, UCFD_FLOAT *dub, UCFD_FLOAT *diag, UCFD_FLOAT *fjmat)
Upper sweep of Block LU-SGS method for Navier-Stokes equations.
Definition: blusgs.c:272
void blusgs_serial_ns_update(UCFD_INT neles, UCFD_FLOAT *uptsb, UCFD_FLOAT *dub, UCFD_FLOAT *subres)
Updates solution array.
Definition: blusgs.c:386
void rans_serial_pre_blusgs(UCFD_INT neles, UCFD_INT nface, UCFD_FLOAT factor, UCFD_FLOAT *fnorm_vol, UCFD_FLOAT *uptsb, UCFD_FLOAT *dt, UCFD_FLOAT *tdiag, UCFD_FLOAT *tjmat, UCFD_FLOAT *dsrc)
Computes Diagonal matrix for Block LU-SGS method for RANS equations.
Definition: blusgs.c:90
void rans_serial_block_lower_sweep(UCFD_INT neles, UCFD_INT nface, UCFD_INT *nei_ele, UCFD_FLOAT *fnorm_vol, UCFD_FLOAT *rhsb, UCFD_FLOAT *dub, UCFD_FLOAT *tdiag, UCFD_FLOAT *tjmat)
Lower sweep of Block LU-SGS method for RANS equations.
Definition: blusgs.c:212
Header file for serial LU-SGS method.
int32_t UCFD_INT
Definition: config.h:20
double UCFD_FLOAT
Definition: config.h:29
ucfd_status_t rans_source_jacobian(UCFD_FLOAT *uf, UCFD_FLOAT tmat[NTURBVARS][NTURBVARS], UCFD_FLOAT *dsrc)
Computes source term Jacobian matrix for RANS equations.
Definition: flux.c:93
Header file for numerical flux funtions.
void lusub(UCFD_FLOAT *LU, UCFD_FLOAT *b)
Forward/Backward substitution function.
Definition: inverse.c:75
void ludcmp(UCFD_FLOAT *A)
LU Decomposition function.
Definition: inverse.c:37
@ UCFD_STATUS_NOT_SUPPORTED
Definition: ucfd_types.h:18