My Project
programmer's documentation
cs_boundary_conditions.h
Go to the documentation of this file.
1 #ifndef __CS_BOUNDARY_CONDITIONS_H__
2 #define __CS_BOUNDARY_CONDITIONS_H__
3 
4 /*============================================================================
5  * Boundary condition handling.
6  *============================================================================*/
7 
8 /*
9  This file is part of Code_Saturne, a general-purpose CFD tool.
10 
11  Copyright (C) 1998-2019 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 /*----------------------------------------------------------------------------
31  * Standard C library headers
32  *----------------------------------------------------------------------------*/
33 
34 /*----------------------------------------------------------------------------
35  * Local headers
36  *----------------------------------------------------------------------------*/
37 
38 #include <ple_locator.h>
39 
40 #include "fvm_nodal.h"
41 #include "fvm_writer.h"
42 
43 #include "cs_base.h"
44 #include "cs_field.h"
45 #include "cs_math.h"
46 #include "cs_mesh_location.h"
47 
48 /*----------------------------------------------------------------------------*/
49 
51 
52 /*============================================================================
53  * Macro definitions
54  *============================================================================*/
55 
56 /*============================================================================
57  * Local type definitions
58  *============================================================================*/
59 
60 /*=============================================================================
61  * Global variables
62  *============================================================================*/
63 
66 extern const int *cs_glob_bc_type;
67 
71 extern const int *cs_glob_bc_face_zone;
72 
73 /*============================================================================
74  * Public function prototypes
75  *============================================================================*/
76 
77 /*----------------------------------------------------------------------------
78  * Handling of boundary condition definition errors and associated output.
79  *
80  * This function checks for errors, and simply returns if no errors are
81  * encountered. In case of error, it outputs helpful information so as to
82  * make it easier to locate the matching faces.
83  *
84  * For each boundary face, bc_type defines the boundary condition type.
85  * As a convention here, zero values correspond to undefined types,
86  * positive values to defined types (with no error), and negative values
87  * to defined types with inconsistent or incompatible values, the
88  * absolute value indicating the original boundary condition type.
89  *
90  * An optional label may be used if the error is related to another
91  * attribute than the boundary type, for appropriate error reporting.
92  *
93  * parameters:
94  * bc_flag <-- array of BC type ids
95  * type_name <-- name of attribute in error, or NULL
96  *----------------------------------------------------------------------------*/
97 
98 void
99 cs_boundary_conditions_error(const int *bc_type,
100  const char *type_name);
101 
102 /*----------------------------------------------------------------------------
103  * Locate shifted boundary face coordinates on possibly filtered
104  * cells or boundary faces for later interpolation.
105  *
106  * parameters:
107  * location_type <-- matching values location (CS_MESH_LOCATION_CELLS
108  * or CS_MESH_LOCATION_BOUNDARY_FACES)
109  * n_location_elts <-- number of selected location elements
110  * n_faces <-- number of selected boundary faces
111  * location_elts <-- list of selected location elements (0 to n-1),
112  * or NULL if no indirection is needed
113  * faces <-- list of selected boundary faces (0 to n-1),
114  * or NULL if no indirection is needed
115  * coord_shift <-- array of coordinates shift relative to selected
116  * boundary faces
117  * coord_stride <-- access stride in coord_shift: 0 for uniform
118  * shift, 1 for "per face" shift.
119  * tolerance <-- relative tolerance for point location.
120  *
121  * returns:
122  * associated locator structure
123  *----------------------------------------------------------------------------*/
124 
125 ple_locator_t *
127  cs_lnum_t n_location_elts,
128  cs_lnum_t n_faces,
129  const cs_lnum_t *location_elts,
130  const cs_lnum_t *faces,
131  cs_real_3_t *coord_shift,
132  int coord_stride,
133  double tolerance);
134 
135 /*----------------------------------------------------------------------------
136  * Set mapped boundary conditions for a given field and mapping locator.
137  *
138  * parameters:
139  * field <-- field whose boundary conditions are set
140  * locator <-- associated mapping locator, as returned
141  * by cs_boundary_conditions_map().
142  * location_type <-- matching values location (CS_MESH_LOCATION_CELLS or
143  * CS_MESH_LOCATION_BOUNDARY_FACES)
144  * normalize <-- normalization option:
145  * 0: values are simply mapped
146  * 1: values are mapped, then multiplied
147  * by a constant factor so that their
148  * surface integral on selected faces
149  * is preserved (relative to the
150  * input values)
151  * 2: as 1, but with a boundary-defined
152  * weight, defined by balance_w
153  * 3: as 1, but with a cell-defined
154  * weight, defined by balance_w
155  * interpolate <-- interpolation option:
156  * 0: values are simply based on matching
157  * cell or face center values
158  * 1: values are based on matching cell
159  * or face center values, corrected
160  * by gradient interpolation
161  * n_faces <-- number of selected boundary faces
162  * faces <-- list of selected boundary faces (0 to n-1),
163  * or NULL if no indirection is needed
164  * balance_w <-- optional balance weight, or NULL
165  * nvar <-- number of variables requiring BC's
166  * rcodcl <-> boundary condition values
167  *----------------------------------------------------------------------------*/
168 
169 void
171  ple_locator_t *locator,
172  cs_mesh_location_type_t location_type,
173  int normalize,
174  int interpolate,
175  cs_lnum_t n_faces,
176  const cs_lnum_t *faces,
177  cs_real_t *balance_w,
178  int nvar,
179  cs_real_t rcodcl[]);
180 
181 /*----------------------------------------------------------------------------
182  * Create the boundary conditions face type and face zone arrays
183  *----------------------------------------------------------------------------*/
184 
185 void
187 
188 /*----------------------------------------------------------------------------
189  * Free the boundary conditions face type and face zone arrays
190  *----------------------------------------------------------------------------*/
191 
192 void
194 
195 /*----------------------------------------------------------------------------*/
206 /*----------------------------------------------------------------------------*/
207 
208 inline static void
210  cs_real_t *af,
211  cs_real_t *b,
212  cs_real_t *bf,
213  cs_real_t qimp,
214  cs_real_t hint)
215 {
216  /* Gradient BCs */
217  *a = -qimp/hint;
218  *b = 1.;
219 
220  /* Flux BCs */
221  *af = qimp;
222  *bf = 0.;
223 }
224 
225 /*----------------------------------------------------------------------------*/
238 /*----------------------------------------------------------------------------*/
239 
240 inline static void
242  cs_real_t *af,
243  cs_real_t *b,
244  cs_real_t *bf,
245  cs_real_t pimp,
246  cs_real_t hint,
247  cs_real_t hext)
248 {
249  if (hext < 0.) {
250 
251  /* Gradient BCs */
252  *a = pimp;
253  *b = 0.;
254 
255  /* Flux BCs */
256  *af = -hint*pimp;
257  *bf = hint;
258 
259  }
260  else {
261 
262  /* Gradient BCs */
263  *a = hext*pimp/(hint + hext);
264  *b = hint /(hint + hext);
265 
266  /* Flux BCs */
267  cs_real_t heq = hint*hext/(hint + hext);
268  *af = -heq*pimp;
269  *bf = heq;
270 
271  }
272 }
273 
274 /*----------------------------------------------------------------------------*/
287 /*----------------------------------------------------------------------------*/
288 
289 inline static void
291  cs_real_3_t af,
292  cs_real_33_t b,
293  cs_real_33_t bf,
294  cs_real_3_t pimpv,
295  cs_real_t hint,
296  cs_real_3_t hextv)
297 {
298  for (int isou = 0 ; isou < 3; isou++) {
299  if (fabs(hextv[isou]) > 0.5*cs_math_infinite_r) {
300 
301  /* Gradient BCs */
302  a[isou] = pimpv[isou];
303  for (int jsou = 0 ; jsou < 3 ; jsou++)
304  b[isou][jsou] = 0.;
305 
306  /* Flux BCs */
307  af[isou] = -hint*pimpv[isou];
308  for (int jsou = 0; jsou < 3; jsou++) {
309  if (jsou == isou)
310  bf[isou][jsou] = hint;
311  else
312  bf[isou][jsou] = 0.;
313  }
314  } else {
315 
316  cs_real_t heq = hint*hextv[isou]/(hint + hextv[isou]);
317 
318  /* Gradient BCs */
319  a[isou] = hextv[isou]*pimpv[isou]/(hint + hextv[isou]);
320  for (int jsou = 0 ; jsou < 3 ; jsou++) {
321  if (jsou == isou)
322  b[isou][jsou] = hint/(hint + hextv[isou]);
323  else
324  b[isou][jsou] = 0.;
325  }
326 
327  /* Flux BCs */
328  af[isou] = -heq*pimpv[isou];
329  for (int jsou = 0 ; jsou < 3 ; jsou++) {
330  if (jsou == isou)
331  bf[isou][jsou] = heq;
332  else
333  bf[isou][jsou] = 0.;
334  }
335  }
336  }
337 }
338 
339 /*----------------------------------------------------------------------------*/
352 /*----------------------------------------------------------------------------*/
353 
354 void
356  cs_real_t *cofaf,
357  cs_real_t *coefb,
358  cs_real_t *cofbf,
359  cs_real_t pimp,
360  cs_real_t cfl,
361  cs_real_t hint);
362 
363 /*----------------------------------------------------------------------------*/
377 /*----------------------------------------------------------------------------*/
378 
379 inline static void
381  cs_real_3_t af,
382  cs_real_33_t b,
383  cs_real_33_t bf,
384  cs_real_3_t pimpv,
385  cs_real_6_t hintt,
386  cs_real_3_t hextv)
387 {
388  for (int isou = 0 ; isou < 3 ; isou++) {
389  if (fabs(hextv[isou]) > 0.5*cs_math_infinite_r) {
390 
391  /* Gradient BCs */
392  a[isou] = pimpv[isou];
393  for (int jsou = 0 ; jsou < 3 ; jsou++)
394  b[isou][jsou] = 0.;
395 
396  } else {
397 
398  cs_exit(1);
399 
400  }
401  }
402 
403  /* Flux BCs */
404  cs_math_sym_33_3_product(hintt, pimpv, af);
405  for (int isou = 0 ; isou < 3 ; isou++)
406  af[isou] = -af[isou];
407 
408  bf[0][0] = hintt[0];
409  bf[1][1] = hintt[1];
410  bf[2][2] = hintt[2];
411  bf[0][1] = hintt[3];
412  bf[1][0] = hintt[3];
413  bf[1][2] = hintt[4];
414  bf[2][1] = hintt[4];
415  bf[0][2] = hintt[5];
416  bf[2][0] = hintt[5];
417 }
418 
419 /*----------------------------------------------------------------------------*/
420 
422 
423 #endif /* __CS_BOUNDARY_CONDITIONS_H__ */
cs_mesh_location_type_t
cs_mesh_location_type_t
Definition: cs_mesh_location.h:60
cs_boundary_conditions_error
void cs_boundary_conditions_error(const int *bc_type, const char *type_name)
Handling of boundary condition definition errors and associated output.
Definition: cs_boundary_conditions.c:339
cs_boundary_conditions_map
ple_locator_t * cs_boundary_conditions_map(cs_mesh_location_type_t location_type, cs_lnum_t n_location_elts, cs_lnum_t n_faces, const cs_lnum_t *location_elts, const cs_lnum_t *faces, cs_real_3_t *coord_shift, int coord_stride, double tolerance)
Locate shifted boundary face coordinates on possibly filtered cells or boundary faces for later inter...
Definition: cs_boundary_conditions.c:393
cs_exit
void cs_exit(int status)
Definition: cs_base.c:1401
cs_fuel_incl::a
double precision, save a
Definition: cs_fuel_incl.f90:146
qimp
void const int const int const int const int int int int int int int int int int int int int int int double * qimp
Definition: cs_gui_boundary_conditions.h:64
cs_boundary_conditions_create
void cs_boundary_conditions_create(void)
Definition: cs_boundary_conditions.c:743
cs_glob_bc_type
const int * cs_glob_bc_type
cs_real_3_t
cs_real_t cs_real_3_t[3]
vector of 3 floating-point values
Definition: cs_defs.h:315
cs_fuel_incl::b
double precision, save b
Definition: cs_fuel_incl.f90:146
cs_boundary_conditions_free
void cs_boundary_conditions_free(void)
Definition: cs_boundary_conditions.c:773
END_C_DECLS
#define END_C_DECLS
Definition: cs_defs.h:468
cs_boundary_conditions_mapped_set
void cs_boundary_conditions_mapped_set(const cs_field_t *f, ple_locator_t *locator, cs_mesh_location_type_t location_type, int normalize, int interpolate, cs_lnum_t n_faces, const cs_lnum_t *faces, cs_real_t *balance_w, int nvar, cs_real_t rcodcl[])
Set mapped boundary conditions for a given field and mapping locator.
Definition: cs_boundary_conditions.c:553
cs_real_t
double cs_real_t
Floating-point value.
Definition: cs_defs.h:302
fvm_writer.h
cs_glob_bc_face_zone
const int * cs_glob_bc_face_zone
BEGIN_C_DECLS
#define BEGIN_C_DECLS
Definition: cs_defs.h:467
nvar
void const cs_int_t const cs_int_t * nvar
Definition: cs_prototypes.h:104
cs_boundary_conditions_set_convective_outlet_scalar
void cs_boundary_conditions_set_convective_outlet_scalar(cs_real_t *coefa, cs_real_t *cofaf, cs_real_t *coefb, cs_real_t *cofbf, cs_real_t pimp, cs_real_t cfl, cs_real_t hint)
Set convective oulet boundary condition for a scalar.
Definition: cs_boundary_conditions.c:799
cs_math.h
cs_boundary_conditions_set_dirichlet_scalar
static void cs_boundary_conditions_set_dirichlet_scalar(cs_real_t *a, cs_real_t *af, cs_real_t *b, cs_real_t *bf, cs_real_t pimp, cs_real_t hint, cs_real_t hext)
Set Dirichlet BC for a scalar for a given face.
Definition: cs_boundary_conditions.h:241
cs_math_sym_33_3_product
static void cs_math_sym_33_3_product(const cs_real_t m[6], const cs_real_t v[3], cs_real_t mv[restrict 3])
Compute the product of a symmetric matrix of 3x3 real values by a vector of 3 real values....
Definition: cs_math.h:561
cs_field.h
cs_math_infinite_r
const cs_real_t cs_math_infinite_r
cs_real_6_t
cs_real_t cs_real_6_t[6]
vector of 6 floating-point values
Definition: cs_defs.h:317
cs_lnum_t
int cs_lnum_t
local mesh entity id
Definition: cs_defs.h:298
cs_boundary_conditions_set_neumann_scalar
static void cs_boundary_conditions_set_neumann_scalar(cs_real_t *a, cs_real_t *af, cs_real_t *b, cs_real_t *bf, cs_real_t qimp, cs_real_t hint)
Set Neumann BC for a scalar for a given face.
Definition: cs_boundary_conditions.h:209
cs_real_33_t
cs_real_t cs_real_33_t[3][3]
3x3 matrix of floating-point values
Definition: cs_defs.h:321
cs_boundary_conditions_set_dirichlet_vector
static void cs_boundary_conditions_set_dirichlet_vector(cs_real_3_t a, cs_real_3_t af, cs_real_33_t b, cs_real_33_t bf, cs_real_3_t pimpv, cs_real_t hint, cs_real_3_t hextv)
Set Dirichlet BC for a vector for a given face.
Definition: cs_boundary_conditions.h:290
fvm_nodal.h
cs_mesh_location.h
cs_boundary_conditions_set_dirichlet_vector_aniso
static void cs_boundary_conditions_set_dirichlet_vector_aniso(cs_real_3_t a, cs_real_3_t af, cs_real_33_t b, cs_real_33_t bf, cs_real_3_t pimpv, cs_real_6_t hintt, cs_real_3_t hextv)
Set Dirichlet BC for a vector for a given face with left anisotropic diffusion.
Definition: cs_boundary_conditions.h:380
cs_field_t
Field descriptor.
Definition: cs_field.h:124
cs_base.h
rcodcl
void const int const int const int const int int int int int int int int int int int int int int int double double double double double double double double double double int double * rcodcl
Definition: cs_gui_boundary_conditions.h:64
cs_tagmr::hext
double precision hext
Definition: cs_tagmr.f90:100