9.1
general documentation
cs_convection_diffusion_priv.h
Go to the documentation of this file.
1#ifndef __CS_CONVECTION_DIFFUSION_PRIV_H__
2#define __CS_CONVECTION_DIFFUSION_PRIV_H__
3
4/*============================================================================
5 * Private functions for convection-diffusion operators.
6 *============================================================================*/
7
8/*
9 This file is part of code_saturne, a general-purpose CFD tool.
10
11 Copyright (C) 1998-2025 EDF S.A.
12
13 This program is free software; you can redistribute it and/or modify it under
14 the terms of the GNU General Public License as published by the Free Software
15 Foundation; either version 2 of the License, or (at your option) any later
16 version.
17
18 This program is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
20 FOR A PARTICULAR PURPOSE. See the GNU General Public License for more
21 details.
22
23 You should have received a copy of the GNU General Public License along with
24 this program; if not, write to the Free Software Foundation, Inc., 51 Franklin
25 Street, Fifth Floor, Boston, MA 02110-1301, USA.
26*/
27
28/*----------------------------------------------------------------------------*/
29
30#include "base/cs_defs.h"
31
32/*----------------------------------------------------------------------------
33 * Local headers
34 *----------------------------------------------------------------------------*/
35
36#include "base/cs_math.h"
37#include "base/cs_field.h"
38#include "base/cs_parameters.h" // for BC types
39
41
42/*=============================================================================
43 * Macro definitions
44 *============================================================================*/
45
46/*============================================================================
47 * Type definition
48 *============================================================================*/
49
50/*============================================================================
51 * Global variables
52 *============================================================================*/
53
54/*============================================================================
55 * Public inlined function
56 *============================================================================*/
57
58/*----------------------------------------------------------------------------
59 * Synchronize halos for scalar variables.
60 *
61 * parameters:
62 * m <-- pointer to associated mesh structure
63 * halo_type <-> halo type
64 * pvar <-> variable
65 *----------------------------------------------------------------------------*/
66
67inline static void
69 cs_halo_type_t halo_type,
70 cs_real_t pvar[])
71{
72 if (m->halo != NULL)
73 cs_halo_sync_var(m->halo, halo_type, pvar);
74}
75
76/*----------------------------------------------------------------------------*/
87/*----------------------------------------------------------------------------*/
88
89CS_F_HOST_DEVICE inline static cs_real_t
91 const cs_real_t nvf_p_c,
92 const cs_real_t nvf_r_f,
93 const cs_real_t nvf_r_c)
94{
95 cs_real_t nvf_p_f;
96
97 cs_real_t beta_m, rfc, r1f, r1, r2, r3, b1, b2;
98
99 switch (limiter) {
100 case CS_NVD_GAMMA: /* Gamma scheme */
101 beta_m = 0.1; /* in [0.1, 0.5] */
102 rfc = (nvf_r_f-nvf_r_c)/(1.-nvf_r_c);
103
104 if (nvf_p_c < beta_m) {
105 nvf_p_f = nvf_p_c*(1.+rfc*(1.-nvf_p_c)/beta_m);
106 } else {
107 r1f = (1.-nvf_r_f)/(1.-nvf_r_c);
108
109 nvf_p_f = r1f*nvf_p_c+rfc;
110 }
111
112 break;
113
114 case CS_NVD_SMART: /* SMART scheme */
115 if (nvf_p_c < (nvf_r_c/3.)) {
116 r1 = nvf_r_f*(1.-3.*nvf_r_c+2.*nvf_r_f);
117 r2 = nvf_r_c*(1.-nvf_r_c);
118
119 nvf_p_f = nvf_p_c*r1/r2;
120 } else if (nvf_p_c <= (nvf_r_c*(1.+nvf_r_f-nvf_r_c)/nvf_r_f)) {
121 rfc = (nvf_r_f-nvf_r_c)/(1.-nvf_r_c);
122 r1f = (1.-nvf_r_f)/(1.-nvf_r_c);
123
124 nvf_p_f = nvf_r_f*(r1f*nvf_p_c/nvf_r_c + rfc);
125 } else {
126 nvf_p_f = 1.;
127 }
128
129 break;
130
131 case CS_NVD_CUBISTA: /* CUBISTA scheme */
132 if (nvf_p_c < (3.*nvf_r_c/4.)) {
133 rfc = (nvf_r_f-nvf_r_c)/(1.-nvf_r_c);
134
135 nvf_p_f = nvf_r_f*(1.+rfc/3.)*nvf_p_c/nvf_r_c;
136 } else if (nvf_p_c <= (nvf_r_c*(1.+2.*(nvf_r_f-nvf_r_c))/(2.*nvf_r_f-nvf_r_c))) {
137 rfc = (nvf_r_f-nvf_r_c)/(1.-nvf_r_c);
138 r1f = (1.-nvf_r_f)/(1.-nvf_r_c);
139
140 nvf_p_f = nvf_r_f*(r1f*nvf_p_c/nvf_r_c+rfc);
141 } else {
142 r1f = (1.-nvf_r_f)/(1.-nvf_r_c);
143
144 nvf_p_f = 1.-.5*r1f*(1.-nvf_p_c);
145 }
146
147 break;
148
149 case CS_NVD_SUPERBEE: /* SuperBee scheme */
150 if (nvf_p_c < (nvf_r_c/(2.-nvf_r_c))) {
151 nvf_p_f = (2.*nvf_r_f-nvf_r_c)*nvf_p_c/nvf_r_c;
152 } else if (nvf_p_c < nvf_r_c) {
153 rfc = (nvf_r_f-nvf_r_c)/(1.-nvf_r_c);
154 r1f = (1.-nvf_r_f)/(1.-nvf_r_c);
155
156 nvf_p_f = r1f*nvf_p_c+rfc;
157 } else if (nvf_p_c < (nvf_r_c/nvf_r_f)) {
158 nvf_p_f = nvf_r_f*nvf_p_c/nvf_r_c;
159 } else {
160 nvf_p_f = 1.;
161 }
162
163 break;
164
165 case CS_NVD_MUSCL: /* MUSCL scheme */
166 if (nvf_p_c < (.5*nvf_r_c)) {
167 nvf_p_f = (2.*nvf_r_f-nvf_r_c)*nvf_p_c/nvf_r_c;
168 } else if (nvf_p_c < (1.+nvf_r_c-nvf_r_f)) {
169 nvf_p_f = nvf_p_c+nvf_r_f-nvf_r_c;
170 } else {
171 nvf_p_f = 1.;
172 }
173
174 break;
175
176 case CS_NVD_MINMOD: /* MINMOD scheme */
177 if (nvf_p_c < nvf_r_c) {
178 nvf_p_f = nvf_r_f*nvf_p_c/nvf_r_c;
179 } else {
180 rfc = (nvf_r_f-nvf_r_c)/(1.-nvf_r_c);
181 r1f = (1.-nvf_r_f)/(1.-nvf_r_c);
182
183 nvf_p_f = r1f*nvf_p_c+rfc;
184 }
185
186 break;
187
188 case CS_NVD_CLAM: /* CLAM scheme */
189 r1 = nvf_r_c*nvf_r_c-nvf_r_f;
190 r2 = nvf_r_c*(nvf_r_c-1.);
191 r3 = nvf_r_f-nvf_r_c;
192
193 nvf_p_f = nvf_p_c*(r1+r3*nvf_p_c)/r2;
194
195 break;
196
197 case CS_NVD_STOIC: /* STOIC scheme */
198 b1 = (nvf_r_c-nvf_r_f)*nvf_r_c;
199 b2 = nvf_r_c+nvf_r_f+2.*nvf_r_f*nvf_r_f-4.*nvf_r_f*nvf_r_c;
200
201 if (nvf_p_c < (b1/b2)) {
202 r1 = -nvf_r_f*(1.-3.*nvf_r_c+2.*nvf_r_f);
203 r2 = nvf_r_c*(nvf_r_c-1.);
204
205 nvf_p_f = nvf_p_c*r1/r2;
206 } else if (nvf_p_c < nvf_r_c) {
207 rfc = (nvf_r_f-nvf_r_c)/(1.-nvf_r_c);
208 r1f = (1.-nvf_r_f)/(1.-nvf_r_c);
209
210 nvf_p_f = r1f*nvf_p_c+rfc;
211 } else if (nvf_p_c < (nvf_r_c*(1.+nvf_r_f-nvf_r_c)/nvf_r_f)) {
212 rfc = (nvf_r_f-nvf_r_c)/(1.-nvf_r_c);
213 r1f = (1.-nvf_r_f)/(1.-nvf_r_c);
214
215 nvf_p_f = nvf_r_f*(nvf_p_c*r1f/nvf_r_c+rfc);
216 } else {
217 nvf_p_f = 1.;
218 }
219
220 break;
221
222 case CS_NVD_OSHER: /* OSHER scheme */
223 if (nvf_p_c < (nvf_r_c/nvf_r_f)) {
224 nvf_p_f = nvf_r_f*nvf_p_c/nvf_r_c;
225 } else {
226 nvf_p_f = 1.;
227 }
228
229 break;
230
231 case CS_NVD_WASEB: /* WASEB scheme */
232 r1 = nvf_r_c*nvf_r_f*(nvf_r_f-nvf_r_c);
233 r2 = 2.*nvf_r_c*(1.-nvf_r_c)-nvf_r_f*(1.-nvf_r_f);
234
235 if (nvf_p_c < (r1/r2)) {
236 nvf_p_f = 2.*nvf_p_c;
237 } else if (nvf_p_c <= (nvf_r_c*(1.+nvf_r_f-nvf_r_c)/nvf_r_f)) {
238 rfc = (nvf_r_f-nvf_r_c)/(1.-nvf_r_c);
239 r1f = (1.-nvf_r_f)/(1.-nvf_r_c);
240
241 nvf_p_f = nvf_r_f*(nvf_p_c*r1f/nvf_r_c+rfc);
242 } else {
243 nvf_p_f = 1.;
244 }
245
246 break;
247
248 default: /* Upwinding */
249 nvf_p_f = nvf_p_c;
250 break;
251 }
252
253 return nvf_p_f;
254}
255
256/*----------------------------------------------------------------------------*/
271/*----------------------------------------------------------------------------*/
272
273CS_F_HOST_DEVICE inline static cs_real_t
275 const cs_nreal_t i_face_u_normal[3],
276 cs_real_t nvf_p_c,
277 cs_real_t nvf_r_f,
278 cs_real_t nvf_r_c,
279 const cs_real_t gradv_c[3],
280 const cs_real_t c_courant)
281{
282 assert(limiter >= CS_NVD_VOF_HRIC);
283
284 cs_real_t nvf_p_f;
285 cs_real_t blend, high_order, low_order, ratio;
286
287 /* Compute gradient angle indicator */
288 cs_real_t dotp = cs::abs(cs_math_3_dot_product(gradv_c, i_face_u_normal));
289 cs_real_t sgrad = cs_math_3_norm(gradv_c);
290 cs_real_t denom = sgrad;
291
292 if (limiter == CS_NVD_VOF_HRIC) { /* M-HRIC scheme */
293 /* High order scheme : Bounded Downwind */
294 if (nvf_p_c <= .5) {
295 high_order = 2.*nvf_p_c;
296 } else {
297 high_order = 1.;
298 }
299
300 /* Low order scheme : MUSCL */
302 nvf_p_c,
303 nvf_r_f,
304 nvf_r_c);
305
306 /* Compute the blending factor */
307 if (denom <= (cs_math_epzero*dotp)) {
308 blend = 1.;
309 } else {
310 ratio = dotp/denom;
311 blend = cs_math_fmin(1., pow(ratio, .5));
312 }
313
314 /* Blending */
315 nvf_p_f = blend*high_order + (1.-blend)*low_order;
316
317 /* Extra blending due to the cell Courant number */
318 if (c_courant < .7 && c_courant > .3) {
319 nvf_p_f = nvf_p_f + (nvf_p_f - low_order)*(.7 - c_courant )/.4;
320 } else if (c_courant >= .7) {
321 nvf_p_f = low_order;
322 }
323 } else if (limiter == CS_NVD_VOF_CICSAM) { /* M-CICSAM scheme */
324 /* High order scheme : HYPER-C + SUPERBEE */
325 if (c_courant <= .3) {
326 high_order = cs_math_fmin(1., nvf_p_c/(c_courant+cs_math_epzero));
327 } else if (c_courant <= .6) {
328 high_order = cs_math_fmin(1., nvf_p_c/.3);
329 } else if (c_courant <= .7) {
331 nvf_p_c,
332 nvf_r_f,
333 nvf_r_c);
334 high_order = 10.*( (.7-c_courant)*cs_math_fmin(1., nvf_p_c/.3)
335 + (c_courant-.6)*superbee);
336 }
337 else {
339 nvf_p_c,
340 nvf_r_f,
341 nvf_r_c);
342 }
343
344 /* Low order scheme : MUSCL */
346 nvf_p_c,
347 nvf_r_f,
348 nvf_r_c);
349
350 /* Compute the blending factor */
351 if (denom <= (cs_math_epzero*dotp)) {
352 blend = 1.;
353 } else {
354 ratio = dotp/denom;
355 blend = cs_math_fmin(1., pow(ratio, 2.));
356 }
357
358 /* Blending */
359 nvf_p_f = blend*high_order + (1.-blend)*low_order;
360 } else { /* STACS scheme */
361 /* High order scheme : SUPERBEE */
363 nvf_p_c,
364 nvf_r_f,
365 nvf_r_c);
366
367 /* Low order scheme : STOIC */
369 nvf_p_c,
370 nvf_r_f,
371 nvf_r_c);
372
373 /* Compute the blending factor */
374 if (denom <= (cs_math_epzero*dotp)) {
375 blend = 1.;
376 } else {
377 ratio = dotp/denom;
378 blend = cs_math_fmin(1., pow(ratio, 4.));
379 }
380
381 /* Blending */
382 nvf_p_f = blend*high_order + (1.-blend)*low_order;
383 }
384
385 return nvf_p_f;
386}
387
388/*----------------------------------------------------------------------------*/
404/*----------------------------------------------------------------------------*/
405
406CS_F_HOST_DEVICE inline static void
408 const cs_real_t pj,
409 const cs_real_t distf,
410 const cs_nreal_t i_face_u_normal[3],
411 const cs_real_t gradi[3],
412 const cs_real_t gradj[3],
413 const cs_real_t grdpai[3],
414 const cs_real_t grdpaj[3],
415 const cs_real_t i_massflux,
416 cs_real_t *testij,
417 cs_real_t *tesqck)
418{
419 cs_real_t dcc, ddi, ddj;
420
421 /* Slope test */
422
423 *testij = grdpai[0]*grdpaj[0]
424 + grdpai[1]*grdpaj[1]
425 + grdpai[2]*grdpaj[2];
426
427 if (i_massflux > 0.) {
428 dcc = gradi[0]*i_face_u_normal[0]
429 + gradi[1]*i_face_u_normal[1]
430 + gradi[2]*i_face_u_normal[2];
431 ddi = grdpai[0]*i_face_u_normal[0]
432 + grdpai[1]*i_face_u_normal[1]
433 + grdpai[2]*i_face_u_normal[2];
434 ddj = (pj-pi)/distf;
435 }
436 else {
437 dcc = gradj[0]*i_face_u_normal[0]
438 + gradj[1]*i_face_u_normal[1]
439 + gradj[2]*i_face_u_normal[2];
440 ddi = (pj-pi)/distf;
441 ddj = grdpaj[0]*i_face_u_normal[0]
442 + grdpaj[1]*i_face_u_normal[1]
443 + grdpaj[2]*i_face_u_normal[2];
444 }
445
446 *tesqck = cs_math_sq(dcc) - cs_math_sq(ddi-ddj);
447}
448
449/*----------------------------------------------------------------------------*/
468/*----------------------------------------------------------------------------*/
469
470template <cs_lnum_t stride>
471CS_F_HOST_DEVICE inline static void
473 const cs_real_t pj[stride],
474 const cs_real_t distf,
475 const cs_nreal_t i_face_u_normal[3],
476 const cs_real_t gradi[stride][3],
477 const cs_real_t gradj[stride][3],
478 const cs_real_t gradsti[stride][3],
479 const cs_real_t gradstj[stride][3],
480 const cs_real_t i_massflux,
481 cs_real_t *testij,
482 cs_real_t *tesqck)
483{
484 cs_real_t dcc[stride], ddi[stride], ddj[stride];
485
486 *testij = 0.;
487 *tesqck = 0.;
488
489 /* Slope test */
490
491 for (int i = 0; i < stride; i++) {
492 *testij += gradsti[i][0]*gradstj[i][0]
493 + gradsti[i][1]*gradstj[i][1]
494 + gradsti[i][2]*gradstj[i][2];
495
496 if (i_massflux > 0.) {
497 dcc[i] = gradi[i][0]*i_face_u_normal[0]
498 + gradi[i][1]*i_face_u_normal[1]
499 + gradi[i][2]*i_face_u_normal[2];
500 ddi[i] = gradsti[i][0]*i_face_u_normal[0]
501 + gradsti[i][1]*i_face_u_normal[1]
502 + gradsti[i][2]*i_face_u_normal[2];
503 ddj[i] = (pj[i]-pi[i])/distf;
504 }
505 else {
506 dcc[i] = gradj[i][0]*i_face_u_normal[0]
507 + gradj[i][1]*i_face_u_normal[1]
508 + gradj[i][2]*i_face_u_normal[2];
509 ddi[i] = (pj[i]-pi[i])/distf;
510 ddj[i] = gradstj[i][0]*i_face_u_normal[0]
511 + gradstj[i][1]*i_face_u_normal[1]
512 + gradstj[i][2]*i_face_u_normal[2];
513 }
514
515 *tesqck += cs_math_sq(dcc[i]) - cs_math_sq(ddi[i]-ddj[i]);
516 }
517}
518
519/*----------------------------------------------------------------------------*/
535/*----------------------------------------------------------------------------*/
536
537CS_F_HOST_DEVICE inline static void
539 const cs_rreal_3_t diipf,
540 const cs_rreal_3_t djjpf,
541 const cs_real_3_t gradi,
542 const cs_real_3_t gradj,
543 const cs_real_t pi,
544 const cs_real_t pj,
545 cs_real_t *recoi,
546 cs_real_t *recoj,
547 cs_real_t *pip,
548 cs_real_t *pjp)
549{
550 cs_real_t gradpf[3] = {0.5*(gradi[0] + gradj[0]),
551 0.5*(gradi[1] + gradj[1]),
552 0.5*(gradi[2] + gradj[2])};
553
554 /* reconstruction only if IRCFLP = 1 */
555 *recoi = bldfrp*(cs_math_3_dot_product(gradpf, diipf));
556 *recoj = bldfrp*(cs_math_3_dot_product(gradpf, djjpf));
557 *pip = pi + *recoi;
558 *pjp = pj + *recoj;
559}
560
561/*----------------------------------------------------------------------------*/
580/*----------------------------------------------------------------------------*/
581
582template <cs_lnum_t stride>
583CS_F_HOST_DEVICE inline static void
585 const cs_rreal_t diipf[3],
586 const cs_rreal_t djjpf[3],
587 const cs_real_t gradi[stride][3],
588 const cs_real_t gradj[stride][3],
589 const cs_real_t pi[stride],
590 const cs_real_t pj[stride],
591 cs_real_t recoi[stride],
592 cs_real_t recoj[stride],
593 cs_real_t pip[stride],
594 cs_real_t pjp[stride])
595{
596 cs_real_t dpvf[3];
597
598 /* x-y-z components, p = u, v, w */
599
600 for (int isou = 0; isou < stride; isou++) {
601
602 for (int jsou = 0; jsou < 3; jsou++)
603 dpvf[jsou] = 0.5*( gradi[isou][jsou]
604 + gradj[isou][jsou]);
605
606 /* reconstruction only if IRCFLP = 1 */
607
608 recoi[isou] = bldfrp*(cs_math_3_dot_product(dpvf, diipf));
609 recoj[isou] = bldfrp*(cs_math_3_dot_product(dpvf, djjpf));
610
611 pip[isou] = pi[isou] + recoi[isou];
612 pjp[isou] = pj[isou] + recoj[isou];
613
614 }
615}
616
617/*----------------------------------------------------------------------------*/
633/*----------------------------------------------------------------------------*/
634
635CS_F_HOST_DEVICE inline static void
636cs_i_relax_c_val(const double relaxp,
637 const cs_real_t pia,
638 const cs_real_t pja,
639 const cs_real_t recoi,
640 const cs_real_t recoj,
641 const cs_real_t pi,
642 const cs_real_t pj,
643 cs_real_t *pir,
644 cs_real_t *pjr,
645 cs_real_t *pipr,
646 cs_real_t *pjpr)
647{
648 *pir = pi/relaxp - (1.-relaxp)/relaxp * pia;
649 *pjr = pj/relaxp - (1.-relaxp)/relaxp * pja;
650
651 *pipr = *pir + recoi;
652 *pjpr = *pjr + recoj;
653}
654
655/*----------------------------------------------------------------------------*/
674/*----------------------------------------------------------------------------*/
675
676template <cs_lnum_t stride>
677CS_F_HOST_DEVICE inline static void
679 const cs_real_t pia[stride],
680 const cs_real_t pja[stride],
681 const cs_real_t recoi[stride],
682 const cs_real_t recoj[stride],
683 const cs_real_t pi[stride],
684 const cs_real_t pj[stride],
685 cs_real_t pir[stride],
686 cs_real_t pjr[stride],
687 cs_real_t pipr[stride],
688 cs_real_t pjpr[stride])
689{
690 for (int isou = 0; isou < stride; isou++) {
691 pir[isou] = pi[isou] /relaxp - (1.-relaxp)/relaxp * pia[isou];
692 pjr[isou] = pj[isou] /relaxp - (1.-relaxp)/relaxp * pja[isou];
693
694 pipr[isou] = pir[isou] + recoi[isou];
695 pjpr[isou] = pjr[isou] + recoj[isou];
696 }
697}
698
699/*----------------------------------------------------------------------------*/
706/*----------------------------------------------------------------------------*/
707
708CS_F_HOST_DEVICE inline static void
710 cs_real_t *pf)
711{
712 *pf = p;
713}
714
715/*----------------------------------------------------------------------------*/
725/*----------------------------------------------------------------------------*/
726
727template <cs_lnum_t stride>
728CS_F_HOST_DEVICE inline static void
730 cs_real_t pf[stride])
731{
732 for (int isou = 0; isou < stride; isou++)
733 pf[isou] = p[isou];
734}
735
736/*----------------------------------------------------------------------------*/
745/*----------------------------------------------------------------------------*/
746
747CS_F_HOST_DEVICE inline static void
749 cs_real_t pip,
750 cs_real_t pjp,
751 cs_real_t *pf)
752{
753 *pf = pnd*pip + (1.-pnd)*pjp;
754}
755
756/*----------------------------------------------------------------------------*/
768/*----------------------------------------------------------------------------*/
769
770template <cs_lnum_t stride>
771CS_F_HOST_DEVICE inline static void
773 const cs_real_t pip[stride],
774 const cs_real_t pjp[stride],
775 cs_real_t pf[stride])
776{
777 for (int isou = 0; isou < stride; isou++)
778 pf[isou] = pnd*pip[isou] + (1.-pnd)*pjp[isou];
779}
780
781/*----------------------------------------------------------------------------*/
791/*----------------------------------------------------------------------------*/
792
793CS_F_HOST_DEVICE inline static void
794cs_solu_f_val(const cs_real_t cell_cen[3],
795 const cs_real_t i_face_cog[3],
796 const cs_real_t grad[3],
797 const cs_real_t p,
798 cs_real_t *pf)
799{
800 cs_real_t df[3];
801
802 df[0] = i_face_cog[0] - cell_cen[0];
803 df[1] = i_face_cog[1] - cell_cen[1];
804 df[2] = i_face_cog[2] - cell_cen[2];
805
806 *pf = p + cs_math_3_dot_product(df, grad);
807}
808
809/*----------------------------------------------------------------------------*/
822/*----------------------------------------------------------------------------*/
823
824template <cs_lnum_t stride>
825CS_F_HOST_DEVICE inline static void
827 const cs_real_t i_face_cog[3],
828 const cs_real_t grad[stride][3],
829 const cs_real_t p[stride],
830 cs_real_t pf[stride])
831{
832 cs_real_t df[3];
833
834 for (cs_lnum_t jsou = 0; jsou < 3; jsou++)
835 df[jsou] = i_face_cog[jsou] - cell_cen[jsou];
836
837 for (cs_lnum_t isou = 0; isou < stride; isou++) {
838 pf[isou] = p[isou] + df[0]*grad[isou][0]
839 + df[1]*grad[isou][1]
840 + df[2]*grad[isou][2];
841 }
842}
843
844/*----------------------------------------------------------------------------*/
854/*----------------------------------------------------------------------------*/
855
856CS_F_HOST_DEVICE inline static void
857cs_blend_f_val(const double blencp,
858 const cs_real_t p,
859 cs_real_t *pf)
860{
861 *pf = blencp * (*pf) + (1. - blencp) * p;
862}
863
864/*----------------------------------------------------------------------------*/
877/*----------------------------------------------------------------------------*/
878
879template <cs_lnum_t stride>
880CS_F_HOST_DEVICE inline static void
881cs_blend_f_val_strided(const double blencp,
882 const cs_real_t p[stride],
883 cs_real_t pf[stride])
884{
885 for (int isou = 0; isou < stride; isou++)
886 pf[isou] = blencp*(pf[isou])+(1.-blencp)*p[isou];
887}
888
889/*----------------------------------------------------------------------------*/
910/*----------------------------------------------------------------------------*/
911
912CS_F_HOST_DEVICE inline static void
913cs_i_conv_flux(const int iconvp,
914 const cs_real_t thetap,
915 const int imasac,
916 const cs_real_t pi,
917 const cs_real_t pj,
918 const cs_real_t pifri,
919 const cs_real_t pifrj,
920 const cs_real_t pjfri,
921 const cs_real_t pjfrj,
922 const cs_real_t i_massflux,
923 const cs_real_t xcppi,
924 const cs_real_t xcppj,
925 cs_real_2_t fluxij)
926{
927 cs_real_t flui, fluj;
928
929 flui = 0.5*(i_massflux + fabs(i_massflux));
930 fluj = 0.5*(i_massflux - fabs(i_massflux));
931
932 fluxij[0] += iconvp*xcppi*(thetap*(flui*pifri + fluj*pjfri) - imasac*i_massflux*pi);
933 fluxij[1] += iconvp*xcppj*(thetap*(flui*pifrj + fluj*pjfrj) - imasac*i_massflux*pj);
934}
935
936/*----------------------------------------------------------------------------*/
957/*----------------------------------------------------------------------------*/
958
959template <cs_lnum_t stride>
960CS_F_HOST_DEVICE inline static void
962 cs_real_t thetap,
963 int imasac,
964 const cs_real_t pi[stride],
965 const cs_real_t pj[stride],
966 const cs_real_t pifri[stride],
967 const cs_real_t pifrj[stride],
968 const cs_real_t pjfri[stride],
969 const cs_real_t pjfrj[stride],
970 cs_real_t i_massflux,
971 cs_real_t fluxi[stride],
972 cs_real_t fluxj[stride])
973{
974 cs_real_t flui, fluj;
975
976 flui = 0.5*(i_massflux + fabs(i_massflux));
977 fluj = 0.5*(i_massflux - fabs(i_massflux));
978
979 for (int isou = 0; isou < stride; isou++) {
980 fluxi[isou] += iconvp*( thetap*(flui*pifri[isou] + fluj*pjfri[isou])
981 - imasac*i_massflux*pi[isou]);
982 fluxj[isou] += iconvp*( thetap*(flui*pifrj[isou] + fluj*pjfrj[isou])
983 - imasac*i_massflux*pj[isou]);
984 }
985}
986
987/*----------------------------------------------------------------------------*/
1000/*----------------------------------------------------------------------------*/
1001
1002CS_F_HOST_DEVICE inline static void
1003cs_i_diff_flux(const int idiffp,
1004 const cs_real_t thetap,
1005 const cs_real_t pip,
1006 const cs_real_t pjp,
1007 const cs_real_t pipr,
1008 const cs_real_t pjpr,
1009 const cs_real_t i_visc,
1010 cs_real_t fluxij[2])
1011{
1012 fluxij[0] += idiffp*thetap*i_visc*(pipr -pjp);
1013 fluxij[1] += idiffp*thetap*i_visc*(pip -pjpr);
1014}
1015
1016/*----------------------------------------------------------------------------*/
1033/*----------------------------------------------------------------------------*/
1034
1035template <cs_lnum_t stride>
1036CS_F_HOST_DEVICE inline static void
1038 cs_real_t thetap,
1039 const cs_real_t pip[stride],
1040 const cs_real_t pjp[stride],
1041 const cs_real_t pipr[stride],
1042 const cs_real_t pjpr[stride],
1043 const cs_real_t i_visc,
1044 cs_real_t fluxi[stride],
1045 cs_real_t fluxj[stride])
1046{
1047 for (int isou = 0; isou < stride; isou++) {
1048 fluxi[isou] += idiffp*thetap*i_visc*(pipr[isou] -pjp[isou]);
1049 fluxj[isou] += idiffp*thetap*i_visc*(pip[isou] -pjpr[isou]);
1050 }
1051}
1052
1053/*----------------------------------------------------------------------------*/
1077/*----------------------------------------------------------------------------*/
1078
1079CS_F_HOST_DEVICE inline static void
1081 const cs_real_t relaxp,
1082 const cs_rreal_t diipf[3],
1083 const cs_rreal_t djjpf[3],
1084 const cs_real_t gradi[3],
1085 const cs_real_t gradj[3],
1086 const cs_real_t pi,
1087 const cs_real_t pj,
1088 const cs_real_t pia,
1089 const cs_real_t pja,
1090 cs_real_t *pifri,
1091 cs_real_t *pifrj,
1092 cs_real_t *pjfri,
1093 cs_real_t *pjfrj,
1094 cs_real_t *pip,
1095 cs_real_t *pjp,
1096 cs_real_t *pipr,
1097 cs_real_t *pjpr)
1098{
1099 cs_real_t pir, pjr;
1100 cs_real_t recoi, recoj;
1101
1103 diipf,
1104 djjpf,
1105 gradi,
1106 gradj,
1107 pi,
1108 pj,
1109 &recoi,
1110 &recoj,
1111 pip,
1112 pjp);
1113
1114 cs_i_relax_c_val(relaxp,
1115 pia,
1116 pja,
1117 recoi,
1118 recoj,
1119 pi,
1120 pj,
1121 &pir,
1122 &pjr,
1123 pipr,
1124 pjpr);
1125
1127 pifrj);
1128 cs_upwind_f_val(pir,
1129 pifri);
1130 cs_upwind_f_val(pj,
1131 pjfri);
1132 cs_upwind_f_val(pjr,
1133 pjfrj);
1134}
1135
1136/*----------------------------------------------------------------------------*/
1163/*----------------------------------------------------------------------------*/
1164
1165template <cs_lnum_t stride>
1166CS_F_HOST_DEVICE inline static void
1168 const cs_real_t relaxp,
1169 const cs_rreal_t diipf[3],
1170 const cs_rreal_t djjpf[3],
1171 const cs_real_t gradi[stride][3],
1172 const cs_real_t gradj[stride][3],
1173 const cs_real_t pi[stride],
1174 const cs_real_t pj[stride],
1175 const cs_real_t pia[stride],
1176 const cs_real_t pja[stride],
1177 cs_real_t pifri[stride],
1178 cs_real_t pifrj[stride],
1179 cs_real_t pjfri[stride],
1180 cs_real_t pjfrj[stride],
1181 cs_real_t pip[stride],
1182 cs_real_t pjp[stride],
1183 cs_real_t pipr[stride],
1184 cs_real_t pjpr[stride])
1185{
1186 cs_real_t pir[stride], pjr[stride];
1187 cs_real_t recoi[stride], recoj[stride];
1188
1189 cs_i_compute_quantities_strided<stride>(bldfrp,
1190 diipf,
1191 djjpf,
1192 gradi,
1193 gradj,
1194 pi,
1195 pj,
1196 recoi,
1197 recoj,
1198 pip,
1199 pjp);
1200
1201 cs_i_relax_c_val_strided<stride>(relaxp,
1202 pia,
1203 pja,
1204 recoi,
1205 recoj,
1206 pi,
1207 pj,
1208 pir,
1209 pjr,
1210 pipr,
1211 pjpr);
1212
1213 cs_upwind_f_val_strided<stride>(pi, pifrj);
1214 cs_upwind_f_val_strided<stride>(pir, pifri);
1215 cs_upwind_f_val_strided<stride>(pj, pjfri);
1216 cs_upwind_f_val_strided<stride>(pjr, pjfrj);
1217}
1218
1219/*----------------------------------------------------------------------------*/
1236/*----------------------------------------------------------------------------*/
1237
1238CS_F_HOST_DEVICE inline static void
1240 const cs_rreal_t diipf[3],
1241 const cs_rreal_t djjpf[3],
1242 const cs_real_t gradi[3],
1243 const cs_real_t gradj[3],
1244 const cs_real_t pi,
1245 const cs_real_t pj,
1246 cs_real_t *pif,
1247 cs_real_t *pjf,
1248 cs_real_t *pip,
1249 cs_real_t *pjp)
1250{
1251 cs_real_t recoi, recoj;
1252
1254 diipf,
1255 djjpf,
1256 gradi,
1257 gradj,
1258 pi,
1259 pj,
1260 &recoi,
1261 &recoj,
1262 pip,
1263 pjp);
1264
1265 cs_upwind_f_val(pi, pif);
1266 cs_upwind_f_val(pj, pjf);
1267}
1268
1269/*----------------------------------------------------------------------------*/
1289/*----------------------------------------------------------------------------*/
1290
1291template <cs_lnum_t stride>
1292CS_F_HOST_DEVICE inline static void
1294 const cs_rreal_t diipf[3],
1295 const cs_rreal_t djjpf[3],
1296 const cs_real_t gradi[stride][3],
1297 const cs_real_t gradj[stride][3],
1298 const cs_real_t pi[stride],
1299 const cs_real_t pj[stride],
1300 cs_real_t pif[stride],
1301 cs_real_t pjf[stride],
1302 cs_real_t pip[stride],
1303 cs_real_t pjp[stride])
1304{
1305 cs_real_t recoi[stride], recoj[stride];
1306
1307 cs_i_compute_quantities_strided<stride>(bldfrp,
1308 diipf,
1309 djjpf,
1310 gradi,
1311 gradj,
1312 pi,
1313 pj,
1314 recoi,
1315 recoj,
1316 pip,
1317 pjp);
1318
1319 cs_upwind_f_val_strided<stride>(pi, pif);
1320 cs_upwind_f_val_strided<stride>(pj, pjf);
1321}
1322
1323/*----------------------------------------------------------------------------*/
1356/*----------------------------------------------------------------------------*/
1357
1358CS_F_HOST_DEVICE inline static void
1360 const int ischcp,
1361 const double relaxp,
1362 const double blencp,
1363 const cs_real_t weight,
1364 const cs_real_t cell_ceni[3],
1365 const cs_real_t cell_cenj[3],
1366 const cs_real_t i_face_cog[3],
1367 const cs_rreal_t diipf[3],
1368 const cs_rreal_t djjpf[3],
1369 const cs_real_t gradi[3],
1370 const cs_real_t gradj[3],
1371 const cs_real_t gradupi[3],
1372 const cs_real_t gradupj[3],
1373 const cs_real_t pi,
1374 const cs_real_t pj,
1375 const cs_real_t pia,
1376 const cs_real_t pja,
1377 cs_real_t *pifri,
1378 cs_real_t *pifrj,
1379 cs_real_t *pjfri,
1380 cs_real_t *pjfrj,
1381 cs_real_t *pip,
1382 cs_real_t *pjp,
1383 cs_real_t *pipr,
1384 cs_real_t *pjpr)
1385{
1386 cs_real_t pir, pjr;
1387 cs_real_t recoi, recoj;
1388
1390 diipf,
1391 djjpf,
1392 gradi,
1393 gradj,
1394 pi,
1395 pj,
1396 &recoi,
1397 &recoj,
1398 pip,
1399 pjp);
1400
1401 cs_i_relax_c_val(relaxp,
1402 pia,
1403 pja,
1404 recoi,
1405 recoj,
1406 pi,
1407 pj,
1408 &pir,
1409 &pjr,
1410 pipr,
1411 pjpr);
1412
1413 if (ischcp == 1) {
1414
1415 /* Centered
1416 --------*/
1417
1418 cs_centered_f_val(weight,
1419 *pip,
1420 *pjpr,
1421 pifrj);
1422 cs_centered_f_val(weight,
1423 *pipr,
1424 *pjp,
1425 pifri);
1426 cs_centered_f_val(weight,
1427 *pipr,
1428 *pjp,
1429 pjfri);
1430 cs_centered_f_val(weight,
1431 *pip,
1432 *pjpr,
1433 pjfrj);
1434
1435 } else if (ischcp == 0) {
1436
1437 /* Original SOLU
1438 --------------*/
1439
1440 cs_solu_f_val(cell_ceni,
1441 i_face_cog,
1442 gradi,
1443 pi,
1444 pifrj);
1445 cs_solu_f_val(cell_ceni,
1446 i_face_cog,
1447 gradi,
1448 pir,
1449 pifri);
1450 cs_solu_f_val(cell_cenj,
1451 i_face_cog,
1452 gradj,
1453 pj,
1454 pjfri);
1455 cs_solu_f_val(cell_cenj,
1456 i_face_cog,
1457 gradj,
1458 pjr,
1459 pjfrj);
1460
1461 } else {
1462
1463 /* SOLU
1464 ----*/
1465
1466 cs_solu_f_val(cell_ceni,
1467 i_face_cog,
1468 gradupi,
1469 pi,
1470 pifrj);
1471 cs_solu_f_val(cell_ceni,
1472 i_face_cog,
1473 gradupi,
1474 pir,
1475 pifri);
1476 cs_solu_f_val(cell_cenj,
1477 i_face_cog,
1478 gradupj,
1479 pj,
1480 pjfri);
1481 cs_solu_f_val(cell_cenj,
1482 i_face_cog,
1483 gradupj,
1484 pjr,
1485 pjfrj);
1486
1487 }
1488
1489 /* Blending
1490 --------*/
1491
1492 cs_blend_f_val(blencp,
1493 pi,
1494 pifrj);
1495 cs_blend_f_val(blencp,
1496 pir,
1497 pifri);
1498 cs_blend_f_val(blencp,
1499 pj,
1500 pjfri);
1501 cs_blend_f_val(blencp,
1502 pjr,
1503 pjfrj);
1504}
1505
1506/*----------------------------------------------------------------------------*/
1540/*----------------------------------------------------------------------------*/
1541
1542template <cs_lnum_t stride>
1543CS_F_HOST_DEVICE inline static void
1545 int ischcp,
1546 double relaxp,
1547 double blencp,
1548 cs_real_t weight,
1549 const cs_real_t cell_ceni[3],
1550 const cs_real_t cell_cenj[3],
1551 const cs_real_t i_face_cog[3],
1552 const cs_rreal_t diipf[3],
1553 const cs_rreal_t djjpf[3],
1554 const cs_real_t gradi[stride][3],
1555 const cs_real_t gradj[stride][3],
1556 const cs_real_t pi[stride],
1557 const cs_real_t pj[stride],
1558 const cs_real_t pia[stride],
1559 const cs_real_t pja[stride],
1560 cs_real_t pifri[stride],
1561 cs_real_t pifrj[stride],
1562 cs_real_t pjfri[stride],
1563 cs_real_t pjfrj[stride],
1564 cs_real_t pip[stride],
1565 cs_real_t pjp[stride],
1566 cs_real_t pipr[stride],
1567 cs_real_t pjpr[stride])
1568
1569{
1570 cs_real_t pir[stride], pjr[stride];
1571 cs_real_t recoi[stride], recoj[stride];
1572
1573 cs_i_compute_quantities_strided<stride>(bldfrp,
1574 diipf,
1575 djjpf,
1576 gradi,
1577 gradj,
1578 pi,
1579 pj,
1580 recoi,
1581 recoj,
1582 pip,
1583 pjp);
1584
1585 cs_i_relax_c_val_strided<stride>(relaxp,
1586 pia,
1587 pja,
1588 recoi,
1589 recoj,
1590 pi,
1591 pj,
1592 pir,
1593 pjr,
1594 pipr,
1595 pjpr);
1596
1597 if (ischcp == 1) {
1598
1599 /* Centered
1600 --------*/
1601
1602 cs_centered_f_val_strided<stride>(weight, pip, pjpr, pifrj);
1603 cs_centered_f_val_strided<stride>(weight, pipr, pjp, pifri);
1604 cs_centered_f_val_strided<stride>(weight, pipr, pjp, pjfri);
1605 cs_centered_f_val_strided<stride>(weight, pip, pjpr, pjfrj);
1606
1607 }
1608 else {
1609
1610 /* Second order
1611 ------------*/
1612
1613 cs_solu_f_val_strided<stride>(cell_ceni,
1614 i_face_cog,
1615 gradi,
1616 pi,
1617 pifrj);
1618 cs_solu_f_val_strided<stride>(cell_ceni,
1619 i_face_cog,
1620 gradi,
1621 pir,
1622 pifri);
1623 cs_solu_f_val_strided<stride>(cell_cenj,
1624 i_face_cog,
1625 gradj,
1626 pj,
1627 pjfri);
1628 cs_solu_f_val_strided<stride>(cell_cenj,
1629 i_face_cog,
1630 gradj,
1631 pjr,
1632 pjfrj);
1633
1634 }
1635
1636 /* Blending
1637 --------*/
1638
1639 cs_blend_f_val_strided<stride>(blencp, pi, pifrj);
1640 cs_blend_f_val_strided<stride>(blencp, pir, pifri);
1641 cs_blend_f_val_strided<stride>(blencp, pj, pjfri);
1642 cs_blend_f_val_strided<stride>(blencp, pjr, pjfrj);
1643}
1644
1645/*----------------------------------------------------------------------------*/
1673/*----------------------------------------------------------------------------*/
1674
1675CS_F_HOST_DEVICE inline static void
1677 const int ischcp,
1678 const double blencp,
1679 const cs_real_t weight,
1680 const cs_real_3_t cell_ceni,
1681 const cs_real_3_t cell_cenj,
1682 const cs_real_3_t i_face_cog,
1683 const cs_real_t hybrid_blend_i,
1684 const cs_real_t hybrid_blend_j,
1685 const cs_rreal_3_t diipf,
1686 const cs_rreal_3_t djjpf,
1687 const cs_real_3_t gradi,
1688 const cs_real_3_t gradj,
1689 const cs_real_3_t gradupi,
1690 const cs_real_3_t gradupj,
1691 const cs_real_t pi,
1692 const cs_real_t pj,
1693 cs_real_t *pif,
1694 cs_real_t *pjf,
1695 cs_real_t *pip,
1696 cs_real_t *pjp)
1697{
1698 cs_real_t recoi, recoj;
1699
1701 diipf,
1702 djjpf,
1703 gradi,
1704 gradj,
1705 pi,
1706 pj,
1707 &recoi,
1708 &recoj,
1709 pip,
1710 pjp);
1711
1712
1713 if (ischcp == 1) {
1714
1715 /* Centered
1716 --------*/
1717
1718 cs_centered_f_val(weight,
1719 *pip,
1720 *pjp,
1721 pif);
1722 cs_centered_f_val(weight,
1723 *pip,
1724 *pjp,
1725 pjf);
1726
1727 } else if (ischcp == 0) {
1728
1729 /* Legacy SOLU
1730 -----------*/
1731
1732 cs_solu_f_val(cell_ceni,
1733 i_face_cog,
1734 gradi,
1735 pi,
1736 pif);
1737 cs_solu_f_val(cell_cenj,
1738 i_face_cog,
1739 gradj,
1740 pj,
1741 pjf);
1742
1743 } else if (ischcp == 3) {
1744
1745 /* Centered
1746 --------*/
1747
1748 cs_centered_f_val(weight,
1749 *pip,
1750 *pjp,
1751 pif);
1752 cs_centered_f_val(weight,
1753 *pip,
1754 *pjp,
1755 pjf);
1756
1757 /* Legacy SOLU
1758 -----------*/
1759 cs_real_t pif_up, pjf_up;
1760 cs_real_t hybrid_blend_interp;
1761
1762 cs_solu_f_val(cell_ceni,
1763 i_face_cog,
1764 gradi,
1765 pi,
1766 &pif_up);
1767 cs_solu_f_val(cell_cenj,
1768 i_face_cog,
1769 gradj,
1770 pj,
1771 &pjf_up);
1772
1773 hybrid_blend_interp = fmin(hybrid_blend_i,hybrid_blend_j);
1774 *pif = hybrid_blend_interp*(*pif) + (1. - hybrid_blend_interp)*pif_up;
1775 *pjf = hybrid_blend_interp*(*pjf) + (1. - hybrid_blend_interp)*pjf_up;
1776
1777 } else {
1778
1779 /* SOLU
1780 ----*/
1781
1782 cs_solu_f_val(cell_ceni,
1783 i_face_cog,
1784 gradupi,
1785 pi,
1786 pif);
1787 cs_solu_f_val(cell_cenj,
1788 i_face_cog,
1789 gradupj,
1790 pj,
1791 pjf);
1792
1793 }
1794
1795 /* Blending
1796 --------*/
1797
1798 cs_blend_f_val(blencp,
1799 pi,
1800 pif);
1801 cs_blend_f_val(blencp,
1802 pj,
1803 pjf);
1804}
1805
1806/*----------------------------------------------------------------------------*/
1835/*----------------------------------------------------------------------------*/
1836
1837template <cs_lnum_t stride>
1838CS_F_HOST_DEVICE inline static void
1840 int ischcp,
1841 double blencp,
1842 cs_real_t weight,
1843 const cs_real_t cell_ceni[3],
1844 const cs_real_t cell_cenj[3],
1845 const cs_real_t i_face_cog[3],
1846 const cs_real_t hybrid_blend_i,
1847 const cs_real_t hybrid_blend_j,
1848 const cs_rreal_t diipf[3],
1849 const cs_rreal_t djjpf[3],
1850 const cs_real_t gradi[stride][3],
1851 const cs_real_t gradj[stride][3],
1852 const cs_real_t pi[stride],
1853 const cs_real_t pj[stride],
1854 cs_real_t pif[stride],
1855 cs_real_t pjf[stride],
1856 cs_real_t pip[stride],
1857 cs_real_t pjp[stride])
1858
1859{
1860 cs_real_t recoi[stride], recoj[stride];
1861
1862 cs_i_compute_quantities_strided<stride>(bldfrp,
1863 diipf,
1864 djjpf,
1865 gradi,
1866 gradj,
1867 pi,
1868 pj,
1869 recoi,
1870 recoj,
1871 pip,
1872 pjp);
1873
1874 if (ischcp == 1) {
1875
1876 /* Centered
1877 --------*/
1878
1879 cs_centered_f_val_strided<stride>(weight, pip, pjp, pif);
1880 cs_centered_f_val_strided<stride>(weight, pip, pjp, pjf);
1881
1882 }
1883 else if (ischcp == 3) {
1884
1885 /* Centered
1886 --------*/
1887
1888 cs_centered_f_val_strided<stride>(weight, pip, pjp, pif);
1889 cs_centered_f_val_strided<stride>(weight, pip, pjp, pjf);
1890
1891 /* SOLU
1892 -----*/
1893 cs_real_t pif_up[stride], pjf_up[stride];
1894 cs_real_t hybrid_blend_interp;
1895
1896 cs_solu_f_val_strided<stride>(cell_ceni,
1897 i_face_cog,
1898 gradi,
1899 pi,
1900 pif_up);
1901 cs_solu_f_val_strided<stride>(cell_cenj,
1902 i_face_cog,
1903 gradj,
1904 pj,
1905 pjf_up);
1906
1907 hybrid_blend_interp = fmin(hybrid_blend_i, hybrid_blend_j);
1908 for (int isou = 0; isou < stride; isou++) {
1909 pif[isou] = hybrid_blend_interp *pif[isou]
1910 + (1. - hybrid_blend_interp)*pif_up[isou];
1911 pjf[isou] = hybrid_blend_interp *pjf[isou]
1912 + (1. - hybrid_blend_interp)*pjf_up[isou];
1913 }
1914 }
1915 else {
1916
1917 /* Second order
1918 ------------*/
1919
1920 cs_solu_f_val_strided<stride>(cell_ceni,
1921 i_face_cog,
1922 gradi,
1923 pi,
1924 pif);
1925 cs_solu_f_val_strided<stride>(cell_cenj,
1926 i_face_cog,
1927 gradj,
1928 pj,
1929 pjf);
1930
1931 }
1932
1933 /* Blending
1934 --------*/
1935
1936 cs_blend_f_val_strided<stride>(blencp, pi, pif);
1937 cs_blend_f_val_strided<stride>(blencp, pj, pjf);
1938}
1939
1940/*----------------------------------------------------------------------------*/
1983/*----------------------------------------------------------------------------*/
1984
1985CS_F_HOST_DEVICE inline static void
1986cs_i_cd_steady_slope_test(bool *upwind_switch,
1987 const int iconvp,
1988 const cs_real_t bldfrp,
1989 const int ischcp,
1990 const double relaxp,
1991 const double blencp,
1992 const double blend_st,
1993 const cs_real_t weight,
1994 const cs_real_t i_dist,
1995 const cs_real_t cell_ceni[3],
1996 const cs_real_t cell_cenj[3],
1997 const cs_nreal_t i_face_u_normal[3],
1998 const cs_real_t i_face_cog[3],
1999 const cs_rreal_t diipf[3],
2000 const cs_rreal_t djjpf[3],
2001 const cs_real_t i_massflux,
2002 const cs_real_t gradi[3],
2003 const cs_real_t gradj[3],
2004 const cs_real_t gradupi[3],
2005 const cs_real_t gradupj[3],
2006 const cs_real_t gradsti[3],
2007 const cs_real_t gradstj[3],
2008 const cs_real_t pi,
2009 const cs_real_t pj,
2010 const cs_real_t pia,
2011 const cs_real_t pja,
2012 cs_real_t *pifri,
2013 cs_real_t *pifrj,
2014 cs_real_t *pjfri,
2015 cs_real_t *pjfrj,
2016 cs_real_t *pip,
2017 cs_real_t *pjp,
2018 cs_real_t *pipr,
2019 cs_real_t *pjpr)
2020{
2021 cs_real_t pir, pjr;
2022 cs_real_t recoi, recoj;
2023 cs_real_t testij, tesqck;
2024
2025 *upwind_switch = false;
2026
2028 diipf,
2029 djjpf,
2030 gradi,
2031 gradj,
2032 pi,
2033 pj,
2034 &recoi,
2035 &recoj,
2036 pip,
2037 pjp);
2038
2039 cs_i_relax_c_val(relaxp,
2040 pia,
2041 pja,
2042 recoi,
2043 recoj,
2044 pi,
2045 pj,
2046 &pir,
2047 &pjr,
2048 pipr,
2049 pjpr);
2050
2051 /* Convection slope test is needed only when iconv >0 */
2052 if (iconvp > 0) {
2054 pj,
2055 i_dist,
2056 i_face_u_normal,
2057 gradi,
2058 gradj,
2059 gradsti,
2060 gradstj,
2061 i_massflux,
2062 &testij,
2063 &tesqck);
2064
2065 if (ischcp==1) {
2066
2067 /* Centered
2068 --------*/
2069
2070 cs_centered_f_val(weight,
2071 *pip,
2072 *pjpr,
2073 pifrj);
2074 cs_centered_f_val(weight,
2075 *pipr,
2076 *pjp,
2077 pifri);
2078 cs_centered_f_val(weight,
2079 *pipr,
2080 *pjp,
2081 pjfri);
2082 cs_centered_f_val(weight,
2083 *pip,
2084 *pjpr,
2085 pjfrj);
2086
2087 } else if (ischcp == 0) {
2088
2089 /* Second order
2090 ------------*/
2091
2092 cs_solu_f_val(cell_ceni,
2093 i_face_cog,
2094 gradi,
2095 pi,
2096 pifrj);
2097 cs_solu_f_val(cell_ceni,
2098 i_face_cog,
2099 gradi,
2100 pir,
2101 pifri);
2102 cs_solu_f_val(cell_cenj,
2103 i_face_cog,
2104 gradj,
2105 pj,
2106 pjfri);
2107 cs_solu_f_val(cell_cenj,
2108 i_face_cog,
2109 gradj,
2110 pjr,
2111 pjfrj);
2112
2113 } else {
2114
2115 /* SOLU
2116 -----*/
2117
2118 cs_solu_f_val(cell_ceni,
2119 i_face_cog,
2120 gradupi,
2121 pi,
2122 pifrj);
2123 cs_solu_f_val(cell_ceni,
2124 i_face_cog,
2125 gradupi,
2126 pir,
2127 pifri);
2128 cs_solu_f_val(cell_cenj,
2129 i_face_cog,
2130 gradupj,
2131 pj,
2132 pjfri);
2133 cs_solu_f_val(cell_cenj,
2134 i_face_cog,
2135 gradupj,
2136 pjr,
2137 pjfrj);
2138 }
2139
2140
2141 /* Slope test: percentage of upwind
2142 ----------------------------------*/
2143
2144 if (tesqck <= 0. || testij <= 0.) {
2145
2146 cs_blend_f_val(blend_st,
2147 pi,
2148 pifrj);
2149 cs_blend_f_val(blend_st,
2150 pir,
2151 pifri);
2152 cs_blend_f_val(blend_st,
2153 pj,
2154 pjfri);
2155 cs_blend_f_val(blend_st,
2156 pjr,
2157 pjfrj);
2158
2159 *upwind_switch = true;
2160
2161 }
2162
2163
2164 /* Blending
2165 --------*/
2166
2167 cs_blend_f_val(blencp,
2168 pi,
2169 pifrj);
2170 cs_blend_f_val(blencp,
2171 pir,
2172 pifri);
2173 cs_blend_f_val(blencp,
2174 pj,
2175 pjfri);
2176 cs_blend_f_val(blencp,
2177 pjr,
2178 pjfrj);
2179
2180 /* If iconv=0 p*fr* are useless */
2181 }
2182 else {
2183 cs_upwind_f_val(pi, pifrj);
2184 cs_upwind_f_val(pir, pifri);
2185 cs_upwind_f_val(pj, pjfri);
2186 cs_upwind_f_val(pjr, pjfrj);
2187 }
2188
2189}
2190
2191/*----------------------------------------------------------------------------*/
2235/*----------------------------------------------------------------------------*/
2236
2237template <cs_lnum_t stride>
2238CS_F_HOST_DEVICE inline static void
2240 int iconvp,
2241 cs_real_t bldfrp,
2242 int ischcp,
2243 double relaxp,
2244 double blencp,
2245 double blend_st,
2246 cs_real_t weight,
2247 cs_real_t i_dist,
2248 const cs_real_t cell_ceni[3],
2249 const cs_real_t cell_cenj[3],
2250 const cs_nreal_t i_face_u_normal[3],
2251 const cs_real_t i_face_cog[3],
2252 const cs_rreal_t diipf[3],
2253 const cs_rreal_t djjpf[3],
2254 cs_real_t i_massflux,
2255 const cs_real_t gradi[stride][3],
2256 const cs_real_t gradj[stride][3],
2257 const cs_real_t grdpai[stride][3],
2258 const cs_real_t grdpaj[stride][3],
2259 const cs_real_t pi[stride],
2260 const cs_real_t pj[stride],
2261 const cs_real_t pia[stride],
2262 const cs_real_t pja[stride],
2263 cs_real_t pifri[stride],
2264 cs_real_t pifrj[stride],
2265 cs_real_t pjfri[stride],
2266 cs_real_t pjfrj[stride],
2267 cs_real_t pip[stride],
2268 cs_real_t pjp[stride],
2269 cs_real_t pipr[stride],
2270 cs_real_t pjpr[stride])
2271{
2272 cs_real_t pir[stride], pjr[stride];
2273 cs_real_t recoi[stride], recoj[stride];
2274 cs_real_t testij, tesqck;
2275 int isou;
2276
2277 cs_i_compute_quantities_strided<stride>(bldfrp,
2278 diipf,
2279 djjpf,
2280 gradi,
2281 gradj,
2282 pi,
2283 pj,
2284 recoi,
2285 recoj,
2286 pip,
2287 pjp);
2288
2289 cs_i_relax_c_val_strided<stride>(relaxp,
2290 pia,
2291 pja,
2292 recoi,
2293 recoj,
2294 pi,
2295 pj,
2296 pir,
2297 pjr,
2298 pipr,
2299 pjpr);
2300
2301 /* Convection slope test is needed only when iconv >0 */
2302 if (iconvp > 0) {
2303 cs_slope_test_strided<stride>(pi,
2304 pj,
2305 i_dist,
2306 i_face_u_normal,
2307 gradi,
2308 gradj,
2309 grdpai,
2310 grdpaj,
2311 i_massflux,
2312 &testij,
2313 &tesqck);
2314
2315 for (isou = 0; isou < stride; isou++) {
2316 if (ischcp==1) {
2317
2318 /* Centered
2319 --------*/
2320
2321 cs_centered_f_val(weight, pip[isou], pjpr[isou], &pifrj[isou]);
2322 cs_centered_f_val(weight, pipr[isou], pjp[isou], &pifri[isou]);
2323 cs_centered_f_val(weight, pipr[isou], pjp[isou], &pjfri[isou]);
2324 cs_centered_f_val(weight, pip[isou], pjpr[isou], &pjfrj[isou]);
2325
2326 }
2327 else {
2328
2329 /* Second order
2330 ------------*/
2331
2332 cs_solu_f_val(cell_ceni, i_face_cog, gradi[isou], pi[isou],
2333 &pifrj[isou]);
2334 cs_solu_f_val(cell_ceni, i_face_cog, gradi[isou], pir[isou],
2335 &pifri[isou]);
2336 cs_solu_f_val(cell_cenj, i_face_cog, gradj[isou], pj[isou],
2337 &pjfri[isou]);
2338 cs_solu_f_val(cell_cenj, i_face_cog, gradj[isou], pjr[isou],
2339 &pjfrj[isou]);
2340
2341 }
2342
2343 }
2344
2345 /* Slope test: percentage of upwind
2346 ----------------------------------*/
2347
2348 if (tesqck <= 0. || testij <= 0.) {
2349
2350 cs_blend_f_val_strided<stride>(blend_st, pi, pifrj);
2351 cs_blend_f_val_strided<stride>(blend_st, pir, pifri);
2352 cs_blend_f_val_strided<stride>(blend_st, pj, pjfri);
2353 cs_blend_f_val_strided<stride>(blend_st, pjr, pjfrj);
2354
2355 *upwind_switch = true;
2356
2357 }
2358
2359 /* Blending
2360 --------*/
2361
2362 cs_blend_f_val_strided<stride>(blencp, pi, pifrj);
2363 cs_blend_f_val_strided<stride>(blencp, pir, pifri);
2364 cs_blend_f_val_strided<stride>(blencp, pj, pjfri);
2365 cs_blend_f_val_strided<stride>(blencp, pjr, pjfrj);
2366
2367 /* If iconv=0 p*fr* are useless */
2368 }
2369 else {
2370 for (isou = 0; isou < stride; isou++) {
2371 cs_upwind_f_val(pi[isou], &pifrj[isou]);
2372 cs_upwind_f_val(pir[isou], &pifri[isou]);
2373 cs_upwind_f_val(pj[isou], &pjfri[isou]);
2374 cs_upwind_f_val(pjr[isou], &pjfrj[isou]);
2375 }
2376 }
2377}
2378
2379/*----------------------------------------------------------------------------*/
2415/*----------------------------------------------------------------------------*/
2416
2417CS_F_HOST_DEVICE inline static void
2419 const int iconvp,
2420 const cs_real_t bldfrp,
2421 const int ischcp,
2422 const double blencp,
2423 const double blend_st,
2424 const cs_real_t weight,
2425 const cs_real_t i_dist,
2426 const cs_real_3_t cell_ceni,
2427 const cs_real_3_t cell_cenj,
2428 const cs_nreal_3_t i_face_u_normal,
2429 const cs_real_3_t i_face_cog,
2430 const cs_rreal_3_t diipf,
2431 const cs_rreal_3_t djjpf,
2432 const cs_real_t i_massflux,
2433 const cs_real_3_t gradi,
2434 const cs_real_3_t gradj,
2435 const cs_real_3_t gradupi,
2436 const cs_real_3_t gradupj,
2437 const cs_real_3_t gradsti,
2438 const cs_real_3_t gradstj,
2439 const cs_real_t pi,
2440 const cs_real_t pj,
2441 cs_real_t *pif,
2442 cs_real_t *pjf,
2443 cs_real_t *pip,
2444 cs_real_t *pjp)
2445{
2446 CS_UNUSED(blend_st);
2447
2448 cs_real_t recoi, recoj;
2449 cs_real_t testij, tesqck;
2450
2451 *upwind_switch = false;
2452
2454 diipf,
2455 djjpf,
2456 gradi,
2457 gradj,
2458 pi,
2459 pj,
2460 &recoi,
2461 &recoj,
2462 pip,
2463 pjp);
2464
2465 /* Convection slope test is needed only when iconv >0 */
2466 if (iconvp > 0) {
2468 pj,
2469 i_dist,
2470 i_face_u_normal,
2471 gradi,
2472 gradj,
2473 gradsti,
2474 gradstj,
2475 i_massflux,
2476 &testij,
2477 &tesqck);
2478
2479 if (ischcp==1) {
2480
2481 /* Centered
2482 --------*/
2483
2484 cs_centered_f_val(weight,
2485 *pip,
2486 *pjp,
2487 pif);
2488 cs_centered_f_val(weight,
2489 *pip,
2490 *pjp,
2491 pjf);
2492
2493 } else if (ischcp == 0) {
2494
2495 /* Original SOLU
2496 --------------*/
2497
2498 cs_solu_f_val(cell_ceni,
2499 i_face_cog,
2500 gradi,
2501 pi,
2502 pif);
2503 cs_solu_f_val(cell_cenj,
2504 i_face_cog,
2505 gradj,
2506 pj,
2507 pjf);
2508
2509 } else {
2510
2511 /* SOLU
2512 -----*/
2513
2514 cs_solu_f_val(cell_ceni,
2515 i_face_cog,
2516 gradupi,
2517 pi,
2518 pif);
2519 cs_solu_f_val(cell_cenj,
2520 i_face_cog,
2521 gradupj,
2522 pj,
2523 pjf);
2524
2525 }
2526
2527 /* Slope test: percentage of upwind
2528 -------------------------------- */
2529
2530 if (tesqck<=0. || testij<=0.) {
2531
2532 cs_blend_f_val(blend_st,
2533 pi,
2534 pif);
2535 cs_blend_f_val(blend_st,
2536 pj,
2537 pjf);
2538
2539 *upwind_switch = true;
2540
2541 }
2542
2543 /* Blending
2544 --------*/
2545
2546 cs_blend_f_val(blencp,
2547 pi,
2548 pif);
2549 cs_blend_f_val(blencp,
2550 pj,
2551 pjf);
2552
2553 /* If iconv=0 p*f are useless */
2554 } else {
2556 pif);
2557 cs_upwind_f_val(pj,
2558 pjf);
2559 }
2560
2561}
2562
2563/*----------------------------------------------------------------------------*/
2574/*----------------------------------------------------------------------------*/
2575
2576CS_F_HOST_DEVICE inline static void
2578 const cs_lnum_t jj,
2579 const cs_real_t i_massflux,
2580 cs_lnum_t *ic,
2581 cs_lnum_t *id)
2582{
2583 if (i_massflux >= 0.) {
2584 *ic = ii;
2585 *id = jj;
2586 } else {
2587 *ic = jj;
2588 *id = ii;
2589 }
2590}
2591
2592/*----------------------------------------------------------------------------*/
2613/*----------------------------------------------------------------------------*/
2614
2615CS_F_HOST_DEVICE inline static void
2617 const double beta,
2618 const cs_real_3_t cell_cen_c,
2619 const cs_real_3_t cell_cen_d,
2620 const cs_nreal_3_t i_face_u_normal,
2621 const cs_real_3_t i_face_cog,
2622 const cs_real_3_t gradv_c,
2623 const cs_real_t p_c,
2624 const cs_real_t p_d,
2625 const cs_real_t local_max_c,
2626 const cs_real_t local_min_c,
2627 const cs_real_t courant_c,
2628 cs_real_t *pif,
2629 cs_real_t *pjf)
2630{
2631 /* distance between face center and central cell center */
2632 /* Distance between face center and central cell center */
2633 cs_real_t dist_fc = cs_math_3_distance(cell_cen_c, i_face_cog);
2634
2635 /* Unit vector and distance between central and downwind cells centers */
2636 cs_real_t diff[3] = {cell_cen_d[0] - cell_cen_c[0],
2637 cell_cen_d[1] - cell_cen_c[1],
2638 cell_cen_d[2] - cell_cen_c[2]};
2639
2640 cs_real_t dist_dc = cs_math_3_norm(diff);
2641 cs_real_t invl = 1./dist_dc;
2642
2643 cs_real_t ndc[3] = {invl*diff[0],
2644 invl*diff[1],
2645 invl*diff[2]};
2646
2647 /* Place the upwind point on the line that joins
2648 the two cells on the upwind side and the same
2649 distance as that between the two cells */
2650 const cs_real_t dist_cu = dist_dc;
2651 const cs_real_t dist_du = dist_dc + dist_cu;
2652
2653 /* Compute the property on the upwind assuming a parabolic
2654 variation of the property between the two cells */
2655 const cs_real_t gradc = cs_math_3_dot_product(gradv_c, ndc);
2656
2657 const cs_real_t grad2c = ((p_d - p_c)/dist_dc - gradc)/dist_dc;
2658
2659 cs_real_t p_u = p_c + (grad2c*dist_cu - gradc)*dist_cu;
2660 p_u = cs::max(cs::min(p_u, local_max_c), local_min_c);
2661
2662 /* Compute the normalised distances */
2663 const cs_real_t nvf_r_f = (dist_fc+dist_cu)/dist_du;
2664 const cs_real_t nvf_r_c = dist_cu/dist_du;
2665
2666 /* Check for the bounds of the NVD diagram and compute the face
2667 property according to the selected NVD scheme */
2668 const cs_real_t _small
2669 = cs_math_epzero * (cs::abs(p_u) + cs::abs(p_c) + cs::abs(p_d));
2670
2671 if (cs::abs(p_d-p_u) <= _small) {
2672 *pif = p_c;
2673 *pjf = p_c;
2674 }
2675 else {
2676 const cs_real_t nvf_p_c = (p_c - p_u)/(p_d - p_u);
2677
2678 if (nvf_p_c <= 0. || nvf_p_c >= 1.) {
2679 *pif = p_c;
2680 *pjf = p_c;
2681 }
2682 else {
2683 cs_real_t nvf_p_f;
2684
2685 /* Highly compressive NVD scheme for VOF */
2686 if (limiter >= CS_NVD_VOF_HRIC) {
2687 nvf_p_f = cs_nvd_vof_scheme_scalar(limiter,
2688 i_face_u_normal,
2689 nvf_p_c,
2690 nvf_r_f,
2691 nvf_r_c,
2692 gradv_c,
2693 courant_c);
2694 } else { /* Regular NVD scheme */
2695 nvf_p_f = cs_nvd_scheme_scalar(limiter,
2696 nvf_p_c,
2697 nvf_r_f,
2698 nvf_r_c);
2699 }
2700
2701 *pif = p_u + nvf_p_f*(p_d - p_u);
2702 *pif = cs::max(cs::min(*pif, local_max_c), local_min_c);
2703
2704 cs_blend_f_val(beta,
2705 p_c,
2706 pif);
2707
2708 *pjf = *pif;
2709 }
2710 }
2711}
2712
2713/*----------------------------------------------------------------------------*/
2750/*----------------------------------------------------------------------------*/
2751
2752template <cs_lnum_t stride>
2753CS_F_HOST_DEVICE inline static void
2755 int iconvp,
2756 cs_real_t bldfrp,
2757 int ischcp,
2758 double blencp,
2759 double blend_st,
2760 cs_real_t weight,
2761 cs_real_t i_dist,
2762 const cs_real_t cell_ceni[3],
2763 const cs_real_t cell_cenj[3],
2764 const cs_nreal_t i_face_u_normal[3],
2765 const cs_real_t i_face_cog[3],
2766 const cs_rreal_t diipf[3],
2767 const cs_rreal_t djjpf[3],
2768 cs_real_t i_massflux,
2769 const cs_real_t gradi[stride][3],
2770 const cs_real_t gradj[stride][3],
2771 const cs_real_t grdpai[stride][3],
2772 const cs_real_t grdpaj[stride][3],
2773 const cs_real_t pi[stride],
2774 const cs_real_t pj[stride],
2775 cs_real_t pif[stride],
2776 cs_real_t pjf[stride],
2777 cs_real_t pip[stride],
2778 cs_real_t pjp[stride])
2779{
2780 cs_real_t recoi[stride], recoj[stride];
2781 cs_real_t testij, tesqck;
2782 int isou;
2783
2784 cs_i_compute_quantities_strided<stride>(bldfrp,
2785 diipf,
2786 djjpf,
2787 gradi,
2788 gradj,
2789 pi,
2790 pj,
2791 recoi,
2792 recoj,
2793 pip,
2794 pjp);
2795
2796 /* Convection slope test is needed only when iconv >0 */
2797 if (iconvp > 0) {
2798 cs_slope_test_strided<stride>(pi,
2799 pj,
2800 i_dist,
2801 i_face_u_normal,
2802 gradi,
2803 gradj,
2804 grdpai,
2805 grdpaj,
2806 i_massflux,
2807 &testij,
2808 &tesqck);
2809
2810 for (isou = 0; isou < stride; isou++) {
2811
2812 if (ischcp==1) {
2813
2814 /* Centered
2815 --------*/
2816
2817 cs_centered_f_val(weight,
2818 pip[isou],
2819 pjp[isou],
2820 &pif[isou]);
2821 cs_centered_f_val(weight,
2822 pip[isou],
2823 pjp[isou],
2824 &pjf[isou]);
2825
2826 }
2827 else {
2828
2829 /* Second order
2830 ------------*/
2831
2832 cs_solu_f_val(cell_ceni,
2833 i_face_cog,
2834 gradi[isou],
2835 pi[isou],
2836 &pif[isou]);
2837 cs_solu_f_val(cell_cenj,
2838 i_face_cog,
2839 gradj[isou],
2840 pj[isou],
2841 &pjf[isou]);
2842 }
2843
2844 }
2845
2846 /* Slope test activated: percentage of upwind */
2847 if (tesqck <= 0. || testij <= 0.) {
2848
2849 /* Upwind
2850 --------*/
2851
2852 cs_blend_f_val_strided<stride>(blend_st, pi, pif);
2853 cs_blend_f_val_strided<stride>(blend_st, pj, pjf);
2854
2855 *upwind_switch = true;
2856 }
2857
2858 /* Blending
2859 --------*/
2860
2861 cs_blend_f_val_strided<stride>(blencp, pi, pif);
2862 cs_blend_f_val_strided<stride>(blencp, pj, pjf);
2863
2864 /* If iconv=0 p*fr* are useless */
2865 }
2866 else {
2867
2868 for (isou = 0; isou < stride; isou++) {
2869 cs_upwind_f_val(pi[isou], &pif[isou]);
2870 cs_upwind_f_val(pj[isou], &pjf[isou]);
2871 }
2872 }
2873}
2874
2875/*----------------------------------------------------------------------------*/
2884/*----------------------------------------------------------------------------*/
2885
2886CS_F_HOST_DEVICE inline static void
2888 const cs_real_t gradi[3],
2889 const cs_real_t bldfrp,
2890 cs_real_t *recoi)
2891{
2892 *recoi = bldfrp * ( gradi[0]*diipb[0]
2893 + gradi[1]*diipb[1]
2894 + gradi[2]*diipb[2]);
2895}
2896
2897/*----------------------------------------------------------------------------*/
2909/*----------------------------------------------------------------------------*/
2910
2911template <cs_lnum_t stride>
2912CS_F_HOST_DEVICE inline static void
2914 const cs_real_t gradi[stride][3],
2915 const cs_real_t bldfrp,
2916 cs_real_t recoi[stride])
2917{
2918 for (int isou = 0; isou < stride; isou++) {
2919 recoi[isou] = bldfrp * ( gradi[isou][0]*diipb[0]
2920 + gradi[isou][1]*diipb[1]
2921 + gradi[isou][2]*diipb[2]);
2922 }
2923}
2924
2925/*----------------------------------------------------------------------------*/
2936/*----------------------------------------------------------------------------*/
2937
2938CS_F_HOST_DEVICE inline static void
2939cs_b_relax_c_val(const double relaxp,
2940 const cs_real_t pi,
2941 const cs_real_t pia,
2942 const cs_real_t recoi,
2943 cs_real_t *pir,
2944 cs_real_t *pipr)
2945{
2946 *pir = pi/relaxp - (1.-relaxp)/relaxp*pia;
2947 *pipr = *pir + recoi;
2948}
2949
2950/*----------------------------------------------------------------------------*/
2964/*----------------------------------------------------------------------------*/
2965
2966template <cs_lnum_t stride>
2967CS_F_HOST_DEVICE inline static void
2968cs_b_relax_c_val_strided(const double relaxp,
2969 const cs_real_t pi[stride],
2970 const cs_real_t pia[stride],
2971 const cs_real_t recoi[stride],
2972 cs_real_t pir[stride],
2973 cs_real_t pipr[stride])
2974{
2975 for (int isou = 0; isou < stride; isou++) {
2976 pir[isou] = pi[isou]/relaxp - (1.-relaxp)/relaxp*pia[isou];
2977 pipr[isou] = pir[isou] + recoi[isou];
2978 }
2979}
2980
2981/*----------------------------------------------------------------------------*/
3006/*----------------------------------------------------------------------------*/
3007
3008CS_F_HOST_DEVICE inline static void
3010 cs_real_t thetap,
3011 int imasac,
3012 int inc,
3013 int bc_type,
3014 int icvfli,
3015 cs_real_t pi,
3016 cs_real_t pir,
3017 cs_real_t pipr,
3018 cs_real_t coface,
3019 cs_real_t cofbce,
3020 cs_real_t b_massflux,
3021 cs_real_t xcpp,
3022 cs_real_t val_f,
3023 cs_real_t *flux)
3024{
3025 cs_real_t flui, fluj, pfac;
3026
3027 /* Computed convective flux */
3028
3029 if (icvfli == 0) {
3030
3031 /* Remove decentering for coupled faces */
3032 if (bc_type == CS_COUPLED_FD) {
3033 flui = 0.0;
3034 fluj = b_massflux;
3035 } else {
3036 flui = 0.5*(b_massflux +fabs(b_massflux));
3037 fluj = 0.5*(b_massflux -fabs(b_massflux));
3038 }
3039
3040 *flux += iconvp*xcpp*(thetap*(flui*pir + fluj*val_f) -imasac*( b_massflux*pi));
3041
3042 /* Imposed convective flux */
3043
3044 } else { // TODO: store val_c = inc*coface + cofbce*pipr
3045
3046 pfac = inc*coface + cofbce*pipr;
3047 *flux += iconvp*xcpp*(-imasac*(b_massflux*pi) + thetap*(pfac));
3048
3049 }
3050}
3051
3052/*----------------------------------------------------------------------------*/
3076/*----------------------------------------------------------------------------*/
3077
3078template <cs_lnum_t stride>
3079CS_F_HOST_DEVICE inline static void
3081 cs_real_t thetap,
3082 int imasac,
3083 int inc,
3084 int bc_type,
3085 int icvfli,
3086 const cs_real_t pi[stride],
3087 const cs_real_t pir[stride],
3088 const cs_real_t pipr[stride],
3089 const cs_real_t coface[stride],
3090 const cs_real_t cofbce[stride][stride],
3091 cs_real_t b_massflux,
3092 cs_real_t pfac[stride],
3093 cs_real_t flux[stride])
3094{
3095 cs_real_t flui, fluj;
3096
3097 /* Computed convective flux */
3098
3099 if (icvfli == 0) {
3100
3101 /* Remove decentering for coupled faces */
3102 if (bc_type == CS_COUPLED_FD) {
3103 flui = 0.0;
3104 fluj = b_massflux;
3105 }
3106 else {
3107 flui = 0.5*(b_massflux +fabs(b_massflux));
3108 fluj = 0.5*(b_massflux -fabs(b_massflux));
3109 }
3110
3111 for (int isou = 0; isou < stride; isou++)
3112 flux[isou] += iconvp*( thetap*(flui*pir[isou] + fluj*pfac[isou])
3113 - imasac*b_massflux*pi[isou]);
3114
3115 /* Imposed convective flux */
3116
3117 }
3118 else {
3119
3120 for (int isou = 0; isou < stride; isou++) {
3121 pfac[isou] = inc*coface[isou];
3122 for (int jsou = 0; jsou < stride; jsou++) {
3123 pfac[isou] += cofbce[isou][jsou]*pipr[jsou];
3124 }
3125 flux[isou] += iconvp*( thetap*pfac[isou]
3126 - imasac*b_massflux*pi[isou]);
3127 }
3128
3129 }
3130}
3131
3132/*----------------------------------------------------------------------------*/
3149/*----------------------------------------------------------------------------*/
3150
3151CS_F_HOST_DEVICE inline static void
3152cs_b_upwind_flux(const int iconvp,
3153 const cs_real_t thetap,
3154 const int imasac,
3155 const int bc_type,
3156 const cs_real_t pi,
3157 const cs_real_t pir,
3158 const cs_real_t val_f,
3159 const cs_real_t b_massflux,
3160 const cs_real_t xcpp,
3161 cs_real_t *flux)
3162{
3163 cs_real_t flui, fluj;
3164
3165 /* Remove decentering for coupled faces */
3166 if (bc_type == CS_COUPLED_FD) {
3167 flui = 0.0;
3168 fluj = b_massflux;
3169 } else {
3170 flui = 0.5*(b_massflux +fabs(b_massflux));
3171 fluj = 0.5*(b_massflux -fabs(b_massflux));
3172 }
3173
3174 *flux += iconvp*xcpp*(thetap*(flui*pir + fluj*val_f) -imasac*( b_massflux*pi));
3175}
3176
3177/*----------------------------------------------------------------------------*/
3197/*----------------------------------------------------------------------------*/
3198
3199CS_F_HOST_DEVICE inline static void
3200cs_b_upwind_flux(const int iconvp,
3201 const cs_real_t thetap,
3202 const int imasac,
3203 const int inc,
3204 const int bc_type,
3205 const cs_real_t pi,
3206 const cs_real_t pir,
3207 const cs_real_t pipr,
3208 const cs_real_t coefap,
3209 const cs_real_t coefbp,
3210 const cs_real_t b_massflux,
3211 const cs_real_t xcpp,
3212 cs_real_t *flux)
3213{
3214 cs_real_t flui, fluj, pfac;
3215
3216 /* Remove decentering for coupled faces */
3217 if (bc_type == CS_COUPLED_FD) {
3218 flui = 0.0;
3219 fluj = b_massflux;
3220 } else {
3221 flui = 0.5*(b_massflux +fabs(b_massflux));
3222 fluj = 0.5*(b_massflux -fabs(b_massflux));
3223 }
3224
3225 pfac = inc*coefap + coefbp*pipr;
3226 *flux += iconvp*xcpp*(thetap*(flui*pir + fluj*pfac) -imasac*( b_massflux*pi));
3227}
3228
3229/*----------------------------------------------------------------------------*/
3247/*----------------------------------------------------------------------------*/
3248
3249template <cs_lnum_t stride>
3250CS_F_HOST_DEVICE inline static void
3252 cs_real_t thetap,
3253 int imasac,
3254 int bc_type,
3255 const cs_real_t pi[stride],
3256 const cs_real_t pir[stride],
3257 const cs_real_t b_massflux,
3258 const cs_real_t pfac[stride],
3259 cs_real_t flux[stride])
3260{
3261 cs_real_t flui, fluj;
3262
3263 /* Remove decentering for coupled faces */
3264 if (bc_type == CS_COUPLED_FD) {
3265 flui = 0.0;
3266 fluj = b_massflux;
3267 } else {
3268 flui = 0.5*(b_massflux +fabs(b_massflux));
3269 fluj = 0.5*(b_massflux -fabs(b_massflux));
3270 }
3271
3272 for (int isou = 0; isou < stride; isou++)
3273 flux[isou] += iconvp*( thetap*(flui*pir[isou] + fluj*pfac[isou])
3274 - imasac*b_massflux*pi[isou]);
3275}
3276
3277/*----------------------------------------------------------------------------*/
3287/*----------------------------------------------------------------------------*/
3288
3289CS_F_HOST_DEVICE inline static void
3290cs_b_diff_flux(const int idiffp,
3291 const cs_real_t thetap,
3292 const cs_real_t pfacd,
3293 const cs_real_t b_visc,
3294 cs_real_t *flux)
3295{
3296 *flux += idiffp*thetap*b_visc*pfacd;
3297}
3298
3299/*----------------------------------------------------------------------------*/
3312/*----------------------------------------------------------------------------*/
3313
3314CS_F_HOST_DEVICE inline static void
3315cs_b_diff_flux(const int idiffp,
3316 const cs_real_t thetap,
3317 const int inc,
3318 const cs_real_t pipr,
3319 const cs_real_t cofafp,
3320 const cs_real_t cofbfp,
3321 const cs_real_t b_visc,
3322 cs_real_t *flux)
3323{
3324 cs_real_t pfacd = inc*cofafp + cofbfp*pipr;
3325 *flux += idiffp*thetap*b_visc*pfacd;
3326}
3327
3328/*----------------------------------------------------------------------------*/
3341/*----------------------------------------------------------------------------*/
3342
3343template <cs_lnum_t stride>
3344CS_F_HOST_DEVICE inline static void
3346 cs_real_t thetap,
3347 cs_real_t b_visc,
3348 const cs_real_t pfacd[stride],
3349 cs_real_t flux[stride])
3350{
3351 for (int isou = 0; isou < stride; isou++)
3352 flux[isou] += idiffp*thetap*b_visc*pfacd[isou];
3353}
3354
3355/*----------------------------------------------------------------------------*/
3369/*----------------------------------------------------------------------------*/
3370
3371CS_F_HOST_DEVICE inline static void
3373 const double relaxp,
3374 const cs_rreal_3_t diipb,
3375 const cs_real_3_t gradi,
3376 const cs_real_t pi,
3377 const cs_real_t pia,
3378 cs_real_t *pir,
3379 cs_real_t *pipr)
3380{
3381 cs_real_t recoi;
3382
3384 gradi,
3385 bldfrp,
3386 &recoi);
3387
3388 cs_b_relax_c_val(relaxp,
3389 pi,
3390 pia,
3391 recoi,
3392 pir,
3393 pipr);
3394}
3395
3396/*----------------------------------------------------------------------------*/
3413/*----------------------------------------------------------------------------*/
3414
3415template <cs_lnum_t stride>
3416CS_F_HOST_DEVICE inline static void
3418 double relaxp,
3419 const cs_rreal_t diipb[3],
3420 const cs_real_t gradi[stride][3],
3421 const cs_real_t pi[stride],
3422 const cs_real_t pia[stride],
3423 cs_real_t pir[stride],
3424 cs_real_t pipr[stride])
3425{
3426 cs_real_t recoi[stride];
3427
3428 cs_b_compute_quantities_strided<stride>(diipb,
3429 gradi,
3430 bldfrp,
3431 recoi);
3432
3433 cs_b_relax_c_val_strided<stride>(relaxp,
3434 pi,
3435 pia,
3436 recoi,
3437 pir,
3438 pipr);
3439}
3440
3441/*----------------------------------------------------------------------------*/
3452/*----------------------------------------------------------------------------*/
3453
3454CS_F_HOST_DEVICE inline static void
3456 const cs_rreal_t diipb[3],
3457 const cs_real_t gradi[3],
3458 const cs_real_t pi,
3459 cs_real_t *pip)
3460{
3461 cs_real_t recoi;
3462
3464 gradi,
3465 bldfrp,
3466 &recoi);
3467
3468 *pip = pi + recoi;
3469}
3470
3471/*----------------------------------------------------------------------------*/
3485/*----------------------------------------------------------------------------*/
3486
3487template <cs_lnum_t stride>
3488CS_F_HOST_DEVICE inline static void
3490 const cs_rreal_t diipb[3],
3491 const cs_real_t gradi[stride][3],
3492 const cs_real_t pi[stride],
3493 cs_real_t pip[stride])
3494{
3495 cs_real_t recoi[stride];
3496
3497 cs_b_compute_quantities_strided<stride>(diipb,
3498 gradi,
3499 bldfrp,
3500 recoi);
3501
3502 for (int isou = 0; isou< stride; isou++)
3503 pip[isou] = pi[isou] + recoi[isou];
3504}
3505
3506/*----------------------------------------------------------------------------*/
3517/*----------------------------------------------------------------------------*/
3518
3519CS_F_HOST_DEVICE inline static void
3521 cs_real_t pi,
3522 cs_real_t pj,
3523 cs_real_t b_visc,
3524 cs_real_t *fluxi)
3525{
3526 *fluxi += idiffp*b_visc*(pi - pj);
3527}
3528
3529/*----------------------------------------------------------------------------*/
3543/*----------------------------------------------------------------------------*/
3544
3545template <cs_lnum_t stride>
3546CS_F_HOST_DEVICE inline static void
3548 const cs_real_t pi[stride],
3549 const cs_real_t pj[stride],
3550 cs_real_t b_visc,
3551 cs_real_t fluxi[stride])
3552{
3553 for (int k = 0; k < stride; k++)
3554 fluxi[k] += idiffp*b_visc*(pi[k] - pj[k]);
3555}
3556
3557/*============================================================================
3558 * Semi private functions
3559 *============================================================================*/
3560
3561/*----------------------------------------------------------------------------*
3562 * Add the explicit part of the convection/diffusion terms of a
3563 * standard transport equation of a scalar field \f$ \varia \f$.
3564 *
3565 * Implemented and documented in cs_convection_diffusion_steady.cpp.
3566 *----------------------------------------------------------------------------*/
3567
3568void
3570(
3571 const cs_field_t *f,
3572 const cs_equation_param_t &eqp,
3573 bool icvflb,
3574 int inc,
3575 cs_real_t *restrict pvar,
3576 const cs_real_t *restrict pvara,
3577 const int icvfli[],
3578 const cs_field_bc_coeffs_t *bc_coeffs,
3579 const cs_real_t val_f_g[],
3580 const cs_real_t i_massflux[],
3581 const cs_real_t b_massflux[],
3582 const cs_real_t i_visc[],
3583 const cs_real_t b_visc[],
3584 const cs_real_t xcpp[],
3585 cs_real_t *restrict rhs
3586);
3587
3588/*----------------------------------------------------------------------------
3589 * \brief Add the explicit part of the convection/diffusion terms of a transport
3590 * equation of a vector or tensor field.
3591 *
3592 * Implemented and documented in cs_convection_diffusion_steady.cpp.
3593 *----------------------------------------------------------------------------*/
3594
3595template <cs_lnum_t stride>
3596void
3598(
3599 cs_field_t *f,
3600 const char *var_name,
3601 const cs_equation_param_t &eqp,
3602 int icvflb,
3603 int inc,
3604 cs_real_t (*pvar)[stride],
3605 const cs_real_t (*pvara)[stride],
3606 const int icvfli[],
3607 const cs_field_bc_coeffs_t *bc_coeffs,
3608 const cs_real_t val_f_g[][stride],
3609 const cs_real_t i_massflux[],
3610 const cs_real_t b_massflux[],
3611 const cs_real_t i_visc[],
3612 const cs_real_t b_visc[],
3613 cs_real_t (*restrict grad)[stride][3],
3614 cs_real_t (*restrict rhs)[stride]
3615);
3616
3617/*----------------------------------------------------------------------------
3618 * Update face flux with convection contribution of a standard transport
3619 * equation of a scalar field \f$ \varia \f$.
3620 *
3621 * Implemented and documented in cs_convection_diffusion_steady.cpp.
3622 *----------------------------------------------------------------------------*/
3623
3624void
3626(
3627 const cs_field_t *f,
3628 const cs_equation_param_t eqp,
3629 int icvflb,
3630 int inc,
3631 cs_real_t *restrict pvar,
3632 const cs_real_t *restrict pvara,
3633 const int icvfli[],
3634 const cs_field_bc_coeffs_t *bc_coeffs,
3635 const cs_real_t val_f_g[],
3636 const cs_real_t i_massflux[],
3637 const cs_real_t b_massflux[],
3638 cs_real_t i_conv_flux[][2],
3639 cs_real_t b_conv_flux[]
3640);
3641
3642/*----------------------------------------------------------------------------*/
3643
3644#endif /* __CS_CONVECTION_DIFFUSION_PRIV_H__ */
cs_nvd_type_t
Definition: cs_convection_diffusion.h:59
@ CS_NVD_SUPERBEE
Definition: cs_convection_diffusion.h:64
@ CS_NVD_SMART
Definition: cs_convection_diffusion.h:62
@ CS_NVD_CUBISTA
Definition: cs_convection_diffusion.h:63
@ CS_NVD_STOIC
Definition: cs_convection_diffusion.h:68
@ CS_NVD_WASEB
Definition: cs_convection_diffusion.h:70
@ CS_NVD_CLAM
Definition: cs_convection_diffusion.h:67
@ CS_NVD_OSHER
Definition: cs_convection_diffusion.h:69
@ CS_NVD_VOF_CICSAM
Definition: cs_convection_diffusion.h:72
@ CS_NVD_GAMMA
Definition: cs_convection_diffusion.h:61
@ CS_NVD_MINMOD
Definition: cs_convection_diffusion.h:66
@ CS_NVD_VOF_HRIC
Definition: cs_convection_diffusion.h:71
@ CS_NVD_MUSCL
Definition: cs_convection_diffusion.h:65
static CS_F_HOST_DEVICE void cs_b_diff_flux_coupling(int idiffp, cs_real_t pi, cs_real_t pj, cs_real_t b_visc, cs_real_t *fluxi)
Add diffusive flux to flux at an internal coupling face.
Definition: cs_convection_diffusion_priv.h:3520
void cs_face_convection_steady_scalar(const cs_field_t *f, const cs_equation_param_t eqp, int icvflb, int inc, cs_real_t *restrict pvar, const cs_real_t *restrict pvara, const int icvfli[], const cs_field_bc_coeffs_t *bc_coeffs, const cs_real_t val_f_g[], const cs_real_t i_massflux[], const cs_real_t b_massflux[], cs_real_t i_conv_flux[][2], cs_real_t b_conv_flux[])
static CS_F_HOST_DEVICE void cs_i_diff_flux_strided(int idiffp, cs_real_t thetap, const cs_real_t pip[stride], const cs_real_t pjp[stride], const cs_real_t pipr[stride], const cs_real_t pjpr[stride], const cs_real_t i_visc, cs_real_t fluxi[stride], cs_real_t fluxj[stride])
Add diffusive fluxes to fluxes at face ij.
Definition: cs_convection_diffusion_priv.h:1037
static CS_F_HOST_DEVICE cs_real_t cs_nvd_vof_scheme_scalar(cs_nvd_type_t limiter, const cs_nreal_t i_face_u_normal[3], cs_real_t nvf_p_c, cs_real_t nvf_r_f, cs_real_t nvf_r_c, const cs_real_t gradv_c[3], const cs_real_t c_courant)
Compute the normalised face scalar using the specified NVD scheme for the case of a Volume-of-Fluid (...
Definition: cs_convection_diffusion_priv.h:274
static CS_F_HOST_DEVICE void cs_centered_f_val_strided(double pnd, const cs_real_t pip[stride], const cs_real_t pjp[stride], cs_real_t pf[stride])
Prepare value at face ij by using a centered scheme.
Definition: cs_convection_diffusion_priv.h:772
static CS_F_HOST_DEVICE void cs_i_cd_unsteady_slope_test_strided(bool *upwind_switch, int iconvp, cs_real_t bldfrp, int ischcp, double blencp, double blend_st, cs_real_t weight, cs_real_t i_dist, const cs_real_t cell_ceni[3], const cs_real_t cell_cenj[3], const cs_nreal_t i_face_u_normal[3], const cs_real_t i_face_cog[3], const cs_rreal_t diipf[3], const cs_rreal_t djjpf[3], cs_real_t i_massflux, const cs_real_t gradi[stride][3], const cs_real_t gradj[stride][3], const cs_real_t grdpai[stride][3], const cs_real_t grdpaj[stride][3], const cs_real_t pi[stride], const cs_real_t pj[stride], cs_real_t pif[stride], cs_real_t pjf[stride], cs_real_t pip[stride], cs_real_t pjp[stride])
Handle preparation of internal face values for the fluxes computation in case of a unsteady algorithm...
Definition: cs_convection_diffusion_priv.h:2754
static CS_F_HOST_DEVICE void cs_i_cd_unsteady(const cs_real_t bldfrp, const int ischcp, const double blencp, const cs_real_t weight, const cs_real_3_t cell_ceni, const cs_real_3_t cell_cenj, const cs_real_3_t i_face_cog, const cs_real_t hybrid_blend_i, const cs_real_t hybrid_blend_j, const cs_rreal_3_t diipf, const cs_rreal_3_t djjpf, const cs_real_3_t gradi, const cs_real_3_t gradj, const cs_real_3_t gradupi, const cs_real_3_t gradupj, const cs_real_t pi, const cs_real_t pj, cs_real_t *pif, cs_real_t *pjf, cs_real_t *pip, cs_real_t *pjp)
Handle preparation of internal face values for the fluxes computation in case of a unsteady algorithm...
Definition: cs_convection_diffusion_priv.h:1676
static CS_F_HOST_DEVICE void cs_i_relax_c_val_strided(cs_real_t relaxp, const cs_real_t pia[stride], const cs_real_t pja[stride], const cs_real_t recoi[stride], const cs_real_t recoj[stride], const cs_real_t pi[stride], const cs_real_t pj[stride], cs_real_t pir[stride], cs_real_t pjr[stride], cs_real_t pipr[stride], cs_real_t pjpr[stride])
Compute relaxed values at cell i and j.
Definition: cs_convection_diffusion_priv.h:678
static CS_F_HOST_DEVICE void cs_b_relax_c_val_strided(const double relaxp, const cs_real_t pi[stride], const cs_real_t pia[stride], const cs_real_t recoi[stride], cs_real_t pir[stride], cs_real_t pipr[stride])
Compute relaxed values at boundary cell i.
Definition: cs_convection_diffusion_priv.h:2968
static CS_F_HOST_DEVICE void cs_i_cd_unsteady_upwind_strided(const cs_real_t bldfrp, const cs_rreal_t diipf[3], const cs_rreal_t djjpf[3], const cs_real_t gradi[stride][3], const cs_real_t gradj[stride][3], const cs_real_t pi[stride], const cs_real_t pj[stride], cs_real_t pif[stride], cs_real_t pjf[stride], cs_real_t pip[stride], cs_real_t pjp[stride])
Handle preparation of internal face values for the fluxes computation in case of an unsteady algorith...
Definition: cs_convection_diffusion_priv.h:1293
static CS_F_HOST_DEVICE void cs_solu_f_val(const cs_real_t cell_cen[3], const cs_real_t i_face_cog[3], const cs_real_t grad[3], const cs_real_t p, cs_real_t *pf)
Prepare value at face ij by using a Second Order Linear Upwind scheme.
Definition: cs_convection_diffusion_priv.h:794
static void cs_sync_scalar_halo(const cs_mesh_t *m, cs_halo_type_t halo_type, cs_real_t pvar[])
Definition: cs_convection_diffusion_priv.h:68
static CS_F_HOST_DEVICE void cs_b_cd_steady_strided(cs_real_t bldfrp, double relaxp, const cs_rreal_t diipb[3], const cs_real_t gradi[stride][3], const cs_real_t pi[stride], const cs_real_t pia[stride], cs_real_t pir[stride], cs_real_t pipr[stride])
Handle preparation of boundary face values for the flux computation in case of a steady algorithm.
Definition: cs_convection_diffusion_priv.h:3417
static CS_F_HOST_DEVICE void cs_i_relax_c_val(const double relaxp, const cs_real_t pia, const cs_real_t pja, const cs_real_t recoi, const cs_real_t recoj, const cs_real_t pi, const cs_real_t pj, cs_real_t *pir, cs_real_t *pjr, cs_real_t *pipr, cs_real_t *pjpr)
Compute relaxed values at cell i and j.
Definition: cs_convection_diffusion_priv.h:636
static CS_F_HOST_DEVICE void cs_b_upwind_flux(const int iconvp, const cs_real_t thetap, const int imasac, const int bc_type, const cs_real_t pi, const cs_real_t pir, const cs_real_t val_f, const cs_real_t b_massflux, const cs_real_t xcpp, cs_real_t *flux)
Add convective flux (substracting the mass accumulation from it) to flux at boundary face....
Definition: cs_convection_diffusion_priv.h:3152
static CS_F_HOST_DEVICE void cs_i_cd_steady_upwind_strided(const cs_real_t bldfrp, const cs_real_t relaxp, const cs_rreal_t diipf[3], const cs_rreal_t djjpf[3], const cs_real_t gradi[stride][3], const cs_real_t gradj[stride][3], const cs_real_t pi[stride], const cs_real_t pj[stride], const cs_real_t pia[stride], const cs_real_t pja[stride], cs_real_t pifri[stride], cs_real_t pifrj[stride], cs_real_t pjfri[stride], cs_real_t pjfrj[stride], cs_real_t pip[stride], cs_real_t pjp[stride], cs_real_t pipr[stride], cs_real_t pjpr[stride])
Handle preparation of internal face values for the fluxes computation in case of a steady algorithm a...
Definition: cs_convection_diffusion_priv.h:1167
static CS_F_HOST_DEVICE void cs_i_cd_unsteady_slope_test(bool *upwind_switch, const int iconvp, const cs_real_t bldfrp, const int ischcp, const double blencp, const double blend_st, const cs_real_t weight, const cs_real_t i_dist, const cs_real_3_t cell_ceni, const cs_real_3_t cell_cenj, const cs_nreal_3_t i_face_u_normal, const cs_real_3_t i_face_cog, const cs_rreal_3_t diipf, const cs_rreal_3_t djjpf, const cs_real_t i_massflux, const cs_real_3_t gradi, const cs_real_3_t gradj, const cs_real_3_t gradupi, const cs_real_3_t gradupj, const cs_real_3_t gradsti, const cs_real_3_t gradstj, const cs_real_t pi, const cs_real_t pj, cs_real_t *pif, cs_real_t *pjf, cs_real_t *pip, cs_real_t *pjp)
Handle preparation of internal face values for the fluxes computation in case of a steady algorithm a...
Definition: cs_convection_diffusion_priv.h:2418
static CS_F_HOST_DEVICE void cs_b_upwind_flux_strided(int iconvp, cs_real_t thetap, int imasac, int bc_type, const cs_real_t pi[stride], const cs_real_t pir[stride], const cs_real_t b_massflux, const cs_real_t pfac[stride], cs_real_t flux[stride])
Add convective flux (substracting the mass accumulation from it) to flux at boundary face....
Definition: cs_convection_diffusion_priv.h:3251
static CS_F_HOST_DEVICE void cs_i_cd_unsteady_nvd(const cs_nvd_type_t limiter, const double beta, const cs_real_3_t cell_cen_c, const cs_real_3_t cell_cen_d, const cs_nreal_3_t i_face_u_normal, const cs_real_3_t i_face_cog, const cs_real_3_t gradv_c, const cs_real_t p_c, const cs_real_t p_d, const cs_real_t local_max_c, const cs_real_t local_min_c, const cs_real_t courant_c, cs_real_t *pif, cs_real_t *pjf)
Handle preparation of internal face values for the convection flux computation in case of an unsteady...
Definition: cs_convection_diffusion_priv.h:2616
static CS_F_HOST_DEVICE void cs_b_diff_flux_coupling_strided(int idiffp, const cs_real_t pi[stride], const cs_real_t pj[stride], cs_real_t b_visc, cs_real_t fluxi[stride])
Add diffusive flux to flux at an internal coupling face for a vector.
Definition: cs_convection_diffusion_priv.h:3547
static CS_F_HOST_DEVICE void cs_upwind_f_val_strided(const cs_real_t p[stride], cs_real_t pf[stride])
Prepare value at face ij by using an upwind scheme.
Definition: cs_convection_diffusion_priv.h:729
static CS_F_HOST_DEVICE void cs_b_cd_steady(const cs_real_t bldfrp, const double relaxp, const cs_rreal_3_t diipb, const cs_real_3_t gradi, const cs_real_t pi, const cs_real_t pia, cs_real_t *pir, cs_real_t *pipr)
Handle preparation of boundary face values for the flux computation in case of a steady algorithm.
Definition: cs_convection_diffusion_priv.h:3372
static CS_F_HOST_DEVICE void cs_b_cd_unsteady_strided(cs_real_t bldfrp, const cs_rreal_t diipb[3], const cs_real_t gradi[stride][3], const cs_real_t pi[stride], cs_real_t pip[stride])
Handle preparation of boundary face values for the flux computation in case of a steady algorithm.
Definition: cs_convection_diffusion_priv.h:3489
static CS_F_HOST_DEVICE void cs_b_cd_unsteady(const cs_real_t bldfrp, const cs_rreal_t diipb[3], const cs_real_t gradi[3], const cs_real_t pi, cs_real_t *pip)
Handle preparation of boundary face values for the flux computation in case of an unsteady algorithm.
Definition: cs_convection_diffusion_priv.h:3455
void cs_convection_diffusion_steady_scalar(const cs_field_t *f, const cs_equation_param_t &eqp, bool icvflb, int inc, cs_real_t *restrict pvar, const cs_real_t *restrict pvara, const int icvfli[], const cs_field_bc_coeffs_t *bc_coeffs, const cs_real_t val_f_g[], const cs_real_t i_massflux[], const cs_real_t b_massflux[], const cs_real_t i_visc[], const cs_real_t b_visc[], const cs_real_t xcpp[], cs_real_t *restrict rhs)
static CS_F_HOST_DEVICE void cs_b_diff_flux(const int idiffp, const cs_real_t thetap, const cs_real_t pfacd, const cs_real_t b_visc, cs_real_t *flux)
Add diffusive flux to flux at boundary face.
Definition: cs_convection_diffusion_priv.h:3290
static CS_F_HOST_DEVICE void cs_b_relax_c_val(const double relaxp, const cs_real_t pi, const cs_real_t pia, const cs_real_t recoi, cs_real_t *pir, cs_real_t *pipr)
Compute relaxed values at boundary cell i.
Definition: cs_convection_diffusion_priv.h:2939
static CS_F_HOST_DEVICE void cs_i_compute_quantities(const cs_real_t bldfrp, const cs_rreal_3_t diipf, const cs_rreal_3_t djjpf, const cs_real_3_t gradi, const cs_real_3_t gradj, const cs_real_t pi, const cs_real_t pj, cs_real_t *recoi, cs_real_t *recoj, cs_real_t *pip, cs_real_t *pjp)
Reconstruct values in I' and J'.
Definition: cs_convection_diffusion_priv.h:538
static CS_F_HOST_DEVICE void cs_slope_test(const cs_real_t pi, const cs_real_t pj, const cs_real_t distf, const cs_nreal_t i_face_u_normal[3], const cs_real_t gradi[3], const cs_real_t gradj[3], const cs_real_t grdpai[3], const cs_real_t grdpaj[3], const cs_real_t i_massflux, cs_real_t *testij, cs_real_t *tesqck)
Compute slope test criteria at internal face between cell i and j.
Definition: cs_convection_diffusion_priv.h:407
static CS_F_HOST_DEVICE void cs_centered_f_val(double pnd, cs_real_t pip, cs_real_t pjp, cs_real_t *pf)
Prepare value at face ij by using a centered scheme.
Definition: cs_convection_diffusion_priv.h:748
static CS_F_HOST_DEVICE void cs_b_imposed_conv_flux(int iconvp, cs_real_t thetap, int imasac, int inc, int bc_type, int icvfli, cs_real_t pi, cs_real_t pir, cs_real_t pipr, cs_real_t coface, cs_real_t cofbce, cs_real_t b_massflux, cs_real_t xcpp, cs_real_t val_f, cs_real_t *flux)
Add convective flux (substracting the mass accumulation from it) to flux at boundary face....
Definition: cs_convection_diffusion_priv.h:3009
static CS_F_HOST_DEVICE void cs_i_cd_steady_upwind(const cs_real_t bldfrp, const cs_real_t relaxp, const cs_rreal_t diipf[3], const cs_rreal_t djjpf[3], const cs_real_t gradi[3], const cs_real_t gradj[3], const cs_real_t pi, const cs_real_t pj, const cs_real_t pia, const cs_real_t pja, cs_real_t *pifri, cs_real_t *pifrj, cs_real_t *pjfri, cs_real_t *pjfrj, cs_real_t *pip, cs_real_t *pjp, cs_real_t *pipr, cs_real_t *pjpr)
Handle preparation of internal face values for the fluxes computation in case of a steady algorithm a...
Definition: cs_convection_diffusion_priv.h:1080
static CS_F_HOST_DEVICE void cs_b_imposed_conv_flux_strided(int iconvp, cs_real_t thetap, int imasac, int inc, int bc_type, int icvfli, const cs_real_t pi[stride], const cs_real_t pir[stride], const cs_real_t pipr[stride], const cs_real_t coface[stride], const cs_real_t cofbce[stride][stride], cs_real_t b_massflux, cs_real_t pfac[stride], cs_real_t flux[stride])
Add convective flux (substracting the mass accumulation from it) to flux at boundary face....
Definition: cs_convection_diffusion_priv.h:3080
static CS_F_HOST_DEVICE void cs_i_cd_steady_slope_test_strided(bool *upwind_switch, int iconvp, cs_real_t bldfrp, int ischcp, double relaxp, double blencp, double blend_st, cs_real_t weight, cs_real_t i_dist, const cs_real_t cell_ceni[3], const cs_real_t cell_cenj[3], const cs_nreal_t i_face_u_normal[3], const cs_real_t i_face_cog[3], const cs_rreal_t diipf[3], const cs_rreal_t djjpf[3], cs_real_t i_massflux, const cs_real_t gradi[stride][3], const cs_real_t gradj[stride][3], const cs_real_t grdpai[stride][3], const cs_real_t grdpaj[stride][3], const cs_real_t pi[stride], const cs_real_t pj[stride], const cs_real_t pia[stride], const cs_real_t pja[stride], cs_real_t pifri[stride], cs_real_t pifrj[stride], cs_real_t pjfri[stride], cs_real_t pjfrj[stride], cs_real_t pip[stride], cs_real_t pjp[stride], cs_real_t pipr[stride], cs_real_t pjpr[stride])
Handle preparation of internal face values for the fluxes computation in case of a steady algorithm a...
Definition: cs_convection_diffusion_priv.h:2239
static CS_F_HOST_DEVICE void cs_i_cd_unsteady_upwind(const cs_real_t bldfrp, const cs_rreal_t diipf[3], const cs_rreal_t djjpf[3], const cs_real_t gradi[3], const cs_real_t gradj[3], const cs_real_t pi, const cs_real_t pj, cs_real_t *pif, cs_real_t *pjf, cs_real_t *pip, cs_real_t *pjp)
Handle preparation of internal face values for the fluxes computation in case of an unsteady algorith...
Definition: cs_convection_diffusion_priv.h:1239
static CS_F_HOST_DEVICE void cs_blend_f_val_strided(const double blencp, const cs_real_t p[stride], cs_real_t pf[stride])
Blend face values for a centered or SOLU scheme with face values for an upwind scheme.
Definition: cs_convection_diffusion_priv.h:881
static CS_F_HOST_DEVICE void cs_i_conv_flux_strided(int iconvp, cs_real_t thetap, int imasac, const cs_real_t pi[stride], const cs_real_t pj[stride], const cs_real_t pifri[stride], const cs_real_t pifrj[stride], const cs_real_t pjfri[stride], const cs_real_t pjfrj[stride], cs_real_t i_massflux, cs_real_t fluxi[stride], cs_real_t fluxj[stride])
Add convective fluxes (substracting the mass accumulation from them) to fluxes at face ij.
Definition: cs_convection_diffusion_priv.h:961
static CS_F_HOST_DEVICE void cs_upwind_f_val(const cs_real_t p, cs_real_t *pf)
Prepare value at face ij by using an upwind scheme.
Definition: cs_convection_diffusion_priv.h:709
static CS_F_HOST_DEVICE void cs_i_diff_flux(const int idiffp, const cs_real_t thetap, const cs_real_t pip, const cs_real_t pjp, const cs_real_t pipr, const cs_real_t pjpr, const cs_real_t i_visc, cs_real_t fluxij[2])
Add diffusive fluxes to fluxes at face ij.
Definition: cs_convection_diffusion_priv.h:1003
static CS_F_HOST_DEVICE void cs_i_cd_steady_strided(cs_real_t bldfrp, int ischcp, double relaxp, double blencp, cs_real_t weight, const cs_real_t cell_ceni[3], const cs_real_t cell_cenj[3], const cs_real_t i_face_cog[3], const cs_rreal_t diipf[3], const cs_rreal_t djjpf[3], const cs_real_t gradi[stride][3], const cs_real_t gradj[stride][3], const cs_real_t pi[stride], const cs_real_t pj[stride], const cs_real_t pia[stride], const cs_real_t pja[stride], cs_real_t pifri[stride], cs_real_t pifrj[stride], cs_real_t pjfri[stride], cs_real_t pjfrj[stride], cs_real_t pip[stride], cs_real_t pjp[stride], cs_real_t pipr[stride], cs_real_t pjpr[stride])
Handle preparation of internal face values for the fluxes computation in case of a steady algorithm a...
Definition: cs_convection_diffusion_priv.h:1544
static CS_F_HOST_DEVICE void cs_solu_f_val_strided(const cs_real_t cell_cen[3], const cs_real_t i_face_cog[3], const cs_real_t grad[stride][3], const cs_real_t p[stride], cs_real_t pf[stride])
Prepare value at face ij by using a Second Order Linear Upwind scheme.
Definition: cs_convection_diffusion_priv.h:826
static CS_F_HOST_DEVICE void cs_b_diff_flux_strided(int idiffp, cs_real_t thetap, cs_real_t b_visc, const cs_real_t pfacd[stride], cs_real_t flux[stride])
Add diffusive flux to flux at boundary face.
Definition: cs_convection_diffusion_priv.h:3345
static CS_F_HOST_DEVICE void cs_i_cd_steady(const cs_real_t bldfrp, const int ischcp, const double relaxp, const double blencp, const cs_real_t weight, const cs_real_t cell_ceni[3], const cs_real_t cell_cenj[3], const cs_real_t i_face_cog[3], const cs_rreal_t diipf[3], const cs_rreal_t djjpf[3], const cs_real_t gradi[3], const cs_real_t gradj[3], const cs_real_t gradupi[3], const cs_real_t gradupj[3], const cs_real_t pi, const cs_real_t pj, const cs_real_t pia, const cs_real_t pja, cs_real_t *pifri, cs_real_t *pifrj, cs_real_t *pjfri, cs_real_t *pjfrj, cs_real_t *pip, cs_real_t *pjp, cs_real_t *pipr, cs_real_t *pjpr)
Handle preparation of internal face values for the fluxes computation in case of a steady algorithm a...
Definition: cs_convection_diffusion_priv.h:1359
static CS_F_HOST_DEVICE void cs_b_compute_quantities_strided(const cs_rreal_t diipb[3], const cs_real_t gradi[stride][3], const cs_real_t bldfrp, cs_real_t recoi[stride])
Reconstruct values in I' at boundary cell i.
Definition: cs_convection_diffusion_priv.h:2913
static CS_F_HOST_DEVICE void cs_central_downwind_cells(const cs_lnum_t ii, const cs_lnum_t jj, const cs_real_t i_massflux, cs_lnum_t *ic, cs_lnum_t *id)
Determine the upwind and downwind sides of an internal face and matching cell indices.
Definition: cs_convection_diffusion_priv.h:2577
static CS_F_HOST_DEVICE void cs_i_conv_flux(const int iconvp, const cs_real_t thetap, const int imasac, const cs_real_t pi, const cs_real_t pj, const cs_real_t pifri, const cs_real_t pifrj, const cs_real_t pjfri, const cs_real_t pjfrj, const cs_real_t i_massflux, const cs_real_t xcppi, const cs_real_t xcppj, cs_real_2_t fluxij)
Add convective fluxes (substracting the mass accumulation from them) to fluxes at face ij.
Definition: cs_convection_diffusion_priv.h:913
void cs_convection_diffusion_steady_strided(cs_field_t *f, const char *var_name, const cs_equation_param_t &eqp, int icvflb, int inc, cs_real_t(*pvar)[stride], const cs_real_t(*pvara)[stride], const int icvfli[], const cs_field_bc_coeffs_t *bc_coeffs, const cs_real_t val_f_g[][stride], const cs_real_t i_massflux[], const cs_real_t b_massflux[], const cs_real_t i_visc[], const cs_real_t b_visc[], cs_real_t(*restrict grad)[stride][3], cs_real_t(*restrict rhs)[stride])
static CS_F_HOST_DEVICE void cs_i_compute_quantities_strided(const cs_real_t bldfrp, const cs_rreal_t diipf[3], const cs_rreal_t djjpf[3], const cs_real_t gradi[stride][3], const cs_real_t gradj[stride][3], const cs_real_t pi[stride], const cs_real_t pj[stride], cs_real_t recoi[stride], cs_real_t recoj[stride], cs_real_t pip[stride], cs_real_t pjp[stride])
Reconstruct values in I' and J'.
Definition: cs_convection_diffusion_priv.h:584
static CS_F_HOST_DEVICE cs_real_t cs_nvd_scheme_scalar(const cs_nvd_type_t limiter, const cs_real_t nvf_p_c, const cs_real_t nvf_r_f, const cs_real_t nvf_r_c)
Compute the normalized face scalar using the specified NVD scheme.
Definition: cs_convection_diffusion_priv.h:90
static CS_F_HOST_DEVICE void cs_b_compute_quantities(const cs_rreal_t diipb[3], const cs_real_t gradi[3], const cs_real_t bldfrp, cs_real_t *recoi)
Reconstruct values in I' at boundary cell i.
Definition: cs_convection_diffusion_priv.h:2887
static CS_F_HOST_DEVICE void cs_i_cd_steady_slope_test(bool *upwind_switch, const int iconvp, const cs_real_t bldfrp, const int ischcp, const double relaxp, const double blencp, const double blend_st, const cs_real_t weight, const cs_real_t i_dist, const cs_real_t cell_ceni[3], const cs_real_t cell_cenj[3], const cs_nreal_t i_face_u_normal[3], const cs_real_t i_face_cog[3], const cs_rreal_t diipf[3], const cs_rreal_t djjpf[3], const cs_real_t i_massflux, const cs_real_t gradi[3], const cs_real_t gradj[3], const cs_real_t gradupi[3], const cs_real_t gradupj[3], const cs_real_t gradsti[3], const cs_real_t gradstj[3], const cs_real_t pi, const cs_real_t pj, const cs_real_t pia, const cs_real_t pja, cs_real_t *pifri, cs_real_t *pifrj, cs_real_t *pjfri, cs_real_t *pjfrj, cs_real_t *pip, cs_real_t *pjp, cs_real_t *pipr, cs_real_t *pjpr)
Handle preparation of internal face values for the fluxes computation in case of a steady algorithm a...
Definition: cs_convection_diffusion_priv.h:1986
static CS_F_HOST_DEVICE void cs_blend_f_val(const double blencp, const cs_real_t p, cs_real_t *pf)
Blend face values for a centered or SOLU scheme with face values for an upwind scheme.
Definition: cs_convection_diffusion_priv.h:857
static CS_F_HOST_DEVICE void cs_i_cd_unsteady_strided(cs_real_t bldfrp, int ischcp, double blencp, cs_real_t weight, const cs_real_t cell_ceni[3], const cs_real_t cell_cenj[3], const cs_real_t i_face_cog[3], const cs_real_t hybrid_blend_i, const cs_real_t hybrid_blend_j, const cs_rreal_t diipf[3], const cs_rreal_t djjpf[3], const cs_real_t gradi[stride][3], const cs_real_t gradj[stride][3], const cs_real_t pi[stride], const cs_real_t pj[stride], cs_real_t pif[stride], cs_real_t pjf[stride], cs_real_t pip[stride], cs_real_t pjp[stride])
Handle preparation of internal face values for the fluxes computation in case of an unsteady algorith...
Definition: cs_convection_diffusion_priv.h:1839
static CS_F_HOST_DEVICE void cs_slope_test_strided(const cs_real_t pi[stride], const cs_real_t pj[stride], const cs_real_t distf, const cs_nreal_t i_face_u_normal[3], const cs_real_t gradi[stride][3], const cs_real_t gradj[stride][3], const cs_real_t gradsti[stride][3], const cs_real_t gradstj[stride][3], const cs_real_t i_massflux, cs_real_t *testij, cs_real_t *tesqck)
Compute slope test criteria at internal face between cell i and j.
Definition: cs_convection_diffusion_priv.h:472
#define restrict
Definition: cs_defs.h:158
#define CS_F_HOST_DEVICE
Definition: cs_defs.h:585
double cs_real_t
Floating-point value.
Definition: cs_defs.h:357
cs_rreal_t cs_rreal_3_t[3]
Definition: cs_defs.h:403
cs_nreal_t cs_nreal_3_t[3]
Definition: cs_defs.h:400
cs_real_t cs_real_3_t[3]
vector of 3 floating-point values
Definition: cs_defs.h:374
cs_real_t cs_real_2_t[2]
vector of 2 floating-point values
Definition: cs_defs.h:373
#define CS_UNUSED(x)
Definition: cs_defs.h:543
double cs_rreal_t
Definition: cs_defs.h:363
int cs_lnum_t
local mesh entity id
Definition: cs_defs.h:350
double cs_nreal_t
Definition: cs_defs.h:361
@ p
Definition: cs_field_pointer.h:67
@ k
Definition: cs_field_pointer.h:72
void cs_halo_sync_var(const cs_halo_t *halo, cs_halo_type_t sync_mode, cs_real_t var[])
Definition: cs_halo.cpp:2360
cs_halo_type_t
Definition: cs_halo.h:57
CS_F_HOST_DEVICE cs_real_t cs_math_3_dot_product(const T u[3], const U v[3])
Compute the dot product of two vectors of 3 real values.
Definition: cs_math.h:176
static CS_F_HOST_DEVICE cs_real_t cs_math_3_distance(const cs_real_t xa[3], const cs_real_t xb[3])
Compute the (euclidean) distance between two points xa and xb in a Cartesian coordinate system of dim...
Definition: cs_math.h:855
CS_F_HOST_DEVICE cs_real_t cs_math_3_norm(const T v[3])
Compute the euclidean norm of a vector of dimension 3.
Definition: cs_math.h:198
const cs_real_t cs_math_epzero
static CS_F_HOST_DEVICE cs_real_t cs_math_fmin(cs_real_t x, cs_real_t y)
Compute the min value of two real values.
Definition: cs_math.h:710
static CS_F_HOST_DEVICE cs_real_t cs_math_sq(cs_real_t x)
Compute the square of a real value.
Definition: cs_math.h:771
@ CS_COUPLED_FD
Definition: cs_parameters.h:101
double precision pi
value with 16 digits
Definition: cstnum.f90:48
CS_F_HOST_DEVICE T max(const T a, const T b)
Definition: cs_defs.h:769
CS_F_HOST_DEVICE T min(const T a, const T b)
Definition: cs_defs.h:746
CS_F_HOST_DEVICE T abs(const T a)
Definition: cs_defs.h:724
Set of parameters to handle an unsteady convection-diffusion-reaction equation with term sources.
Definition: cs_equation_param.h:194
Field boundary condition descriptor (for variables)
Definition: cs_field.h:120
Field descriptor.
Definition: cs_field.h:156
Definition: cs_mesh.h:85
cs_halo_t * halo
Definition: cs_mesh.h:156