7.1
general documentation
cs_advection_field.h
Go to the documentation of this file.
1 #ifndef __CS_ADVECTION_FIELD_H__
2 #define __CS_ADVECTION_FIELD_H__
3 
4 /*============================================================================
5  * Manage the definition/setting of advection fields
6  *============================================================================*/
7 
8 /*
9  This file is part of Code_Saturne, a general-purpose CFD tool.
10 
11  Copyright (C) 1998-2021 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  * Local headers
32  *----------------------------------------------------------------------------*/
33 
34 #include "cs_cdo_connect.h"
35 #include "cs_cdo_local.h"
36 #include "cs_cdo_quantities.h"
37 #include "cs_field.h"
38 #include "cs_mesh_location.h"
39 #include "cs_param_types.h"
40 #include "cs_property.h"
41 #include "cs_xdef.h"
42 #include "cs_xdef_eval.h"
43 
44 /*----------------------------------------------------------------------------*/
45 
47 
48 /*============================================================================
49  * Macro definitions
50  *============================================================================*/
51 
60 #define CS_ADVECTION_FIELD_POST_COURANT (1 << 0)
61 
64 /*============================================================================
65  * Type definitions
66  *============================================================================*/
67 
69 
120 typedef enum {
121 
122  /* Category of advection field */
123  /* --------------------------- */
124 
125  CS_ADVECTION_FIELD_NAVSTO = 1<<0, /* = 1 */
126  CS_ADVECTION_FIELD_GWF = 1<<1, /* = 2 */
127  CS_ADVECTION_FIELD_USER = 1<<2, /* = 4 */
128 
129  /* Type of advection field */
130  /* ----------------------- */
131 
134 
135  /* Optional */
136  /* -------- */
137 
138  CS_ADVECTION_FIELD_STEADY = 1<<5, /* = 32 */
139  CS_ADVECTION_FIELD_LEGACY_FV = 1<<6, /* = 64 */
142 
144 
149 typedef struct {
150 
195  int id;
196  char *restrict name;
199 
204 
205  /* We assume that there is only one definition associated to an advection
206  field inside the computational domain */
208 
209  /* Optional: Definition(s) for the boundary flux */
212  short int *bdy_def_ids;
213 
215 
216 /*============================================================================
217  * Global variables
218  *============================================================================*/
219 
220 /*============================================================================
221  * Inline public function prototypes
222  *============================================================================*/
223 
224 /*----------------------------------------------------------------------------*/
231 /*----------------------------------------------------------------------------*/
232 
233 static inline void
234 cs_advection_field_set_status(cs_adv_field_t *adv,
236 {
237  if (adv == NULL)
238  return;
239 
240  adv->status = status;
241 }
242 
243 /*----------------------------------------------------------------------------*/
251 /*----------------------------------------------------------------------------*/
252 
253 static inline bool
254 cs_advection_field_is_uniform(const cs_adv_field_t *adv)
255 {
256  if (adv == NULL)
257  return false;
258 
260  return true;
261  else
262  return false;
263 }
264 
265 /*----------------------------------------------------------------------------*/
274 /*----------------------------------------------------------------------------*/
275 
276 static inline bool
277 cs_advection_field_is_cellwise(const cs_adv_field_t *adv)
278 {
279  if (adv == NULL)
280  return false;
281 
282  cs_flag_t state = adv->definition->state;
283 
284  if (state & CS_FLAG_STATE_UNIFORM)
285  return true;
286  if (state & CS_FLAG_STATE_CELLWISE)
287  return true;
288  else
289  return false;
290 }
291 
292 /*----------------------------------------------------------------------------*/
300 /*----------------------------------------------------------------------------*/
301 
302 static inline const char *
303 cs_advection_field_get_name(const cs_adv_field_t *adv)
304 {
305  if (adv == NULL)
306  return NULL;
307 
308  return adv->name;
309 }
310 
311 /*----------------------------------------------------------------------------*/
320 /*----------------------------------------------------------------------------*/
321 
322 static inline cs_xdef_type_t
323 cs_advection_field_get_deftype(const cs_adv_field_t *adv)
324 {
325  if (adv == NULL)
326  return CS_N_XDEF_TYPES;
327 
328  return cs_xdef_get_type(adv->definition);
329 }
330 
331 /*----------------------------------------------------------------------------*/
341 /*----------------------------------------------------------------------------*/
342 
343 static inline cs_field_t *
344 cs_advection_field_get_field(const cs_adv_field_t *adv,
345  cs_mesh_location_type_t ml_type)
346 {
347  cs_field_t *f = NULL;
348  if (adv == NULL)
349  return f;
350 
351  switch (ml_type) {
352 
354  if (adv->cell_field_id > -1)
355  f = cs_field_by_id(adv->cell_field_id);
356  else
357  f = NULL;
358  break;
359 
361  if (adv->int_field_id > -1)
362  f = cs_field_by_id(adv->int_field_id);
363  else
364  f = NULL;
365  break;
366 
368  if (adv->bdy_field_id > -1)
369  f = cs_field_by_id(adv->bdy_field_id);
370  else
371  f = NULL;
372  break;
373 
375  if (adv->vtx_field_id > -1)
376  f = cs_field_by_id(adv->vtx_field_id);
377  else
378  f = NULL;
379  break;
380 
381  default:
382  bft_error(__FILE__, __LINE__, 0,
383  " %s: Invalid mesh location type %d.\n"
384  " Stop retrieving the advection field.\n",
385  __func__, ml_type);
386  }
387 
388  return f;
389 }
390 
391 /*============================================================================
392  * Public function prototypes
393  *============================================================================*/
394 
395 /*----------------------------------------------------------------------------*/
402 /*----------------------------------------------------------------------------*/
403 
404 void
406  const cs_cdo_connect_t *connect);
407 
408 /*----------------------------------------------------------------------------*/
414 /*----------------------------------------------------------------------------*/
415 
416 int
418 
419 /*----------------------------------------------------------------------------*/
428 /*----------------------------------------------------------------------------*/
429 
431 cs_advection_field_by_name(const char *name);
432 
433 /*----------------------------------------------------------------------------*/
442 /*----------------------------------------------------------------------------*/
443 
446 
447 /*----------------------------------------------------------------------------*/
455 /*----------------------------------------------------------------------------*/
456 
458 cs_advection_field_add_user(const char *name);
459 
460 /*----------------------------------------------------------------------------*/
469 /*----------------------------------------------------------------------------*/
470 
472 cs_advection_field_add(const char *name,
474 
475 /*----------------------------------------------------------------------------*/
481 /*----------------------------------------------------------------------------*/
482 
483 void
485 
486 /*----------------------------------------------------------------------------*/
495 /*----------------------------------------------------------------------------*/
496 
497 bool
499  const char *ref_name);
500 
501 /*----------------------------------------------------------------------------*/
505 /*----------------------------------------------------------------------------*/
506 
507 void
509 
510 /*----------------------------------------------------------------------------*/
517 /*----------------------------------------------------------------------------*/
518 
519 void
521  cs_flag_t post_flag);
522 
523 /*----------------------------------------------------------------------------*/
530 /*----------------------------------------------------------------------------*/
531 
532 void
534  cs_real_t vector[3]);
535 
536 /*----------------------------------------------------------------------------*/
544 /*----------------------------------------------------------------------------*/
545 
546 void
548  cs_analytic_func_t *func,
549  void *input);
550 
551 /*----------------------------------------------------------------------------*/
560 /*----------------------------------------------------------------------------*/
561 
562 void
564  cs_flag_t loc_flag,
565  cs_dof_func_t *func,
566  void *input);
567 
568 /*----------------------------------------------------------------------------*/
579 /*----------------------------------------------------------------------------*/
580 
581 void
583  cs_flag_t loc,
584  cs_real_t *array,
585  bool is_owner,
586  cs_lnum_t *index);
587 
588 /*----------------------------------------------------------------------------*/
595 /*----------------------------------------------------------------------------*/
596 
597 void
599  cs_field_t *field);
600 
601 /*----------------------------------------------------------------------------*/
610 /*----------------------------------------------------------------------------*/
611 
612 void
614  const char *zname,
615  cs_real_t normal_flux);
616 
617 /*----------------------------------------------------------------------------*/
627 /*----------------------------------------------------------------------------*/
628 
629 void
631  const char *zname,
632  cs_analytic_func_t *func,
633  void *input);
634 
635 /*----------------------------------------------------------------------------*/
648 /*----------------------------------------------------------------------------*/
649 
650 void
652  const char *zname,
653  cs_flag_t loc,
654  cs_real_t *array,
655  bool is_owner,
656  cs_lnum_t *index);
657 
658 /*----------------------------------------------------------------------------*/
666 /*----------------------------------------------------------------------------*/
667 
668 void
670  cs_field_t *field);
671 
672 /*----------------------------------------------------------------------------*/
677 /*----------------------------------------------------------------------------*/
678 
679 void
681 
682 /*----------------------------------------------------------------------------*/
687 /*----------------------------------------------------------------------------*/
688 
689 void
691 
692 /*----------------------------------------------------------------------------*/
700 /*----------------------------------------------------------------------------*/
701 
702 void
704  const cs_adv_field_t *adv,
705  cs_nvec3_t *vect);
706 
707 /*----------------------------------------------------------------------------*/
718 /*----------------------------------------------------------------------------*/
719 
720 void
722  const cs_cell_mesh_t *cm,
723  const cs_real_3_t xyz,
724  cs_real_t time_eval,
725  cs_nvec3_t *eval);
726 
727 /*----------------------------------------------------------------------------*/
736 /*----------------------------------------------------------------------------*/
737 
738 void
740  cs_real_t time_eval,
741  cs_real_t *cell_values);
742 
743 /*----------------------------------------------------------------------------*/
751 /*----------------------------------------------------------------------------*/
752 
753 void
755  cs_real_t time_eval,
756  cs_real_t *vtx_values);
757 
758 /*----------------------------------------------------------------------------*/
767 /*----------------------------------------------------------------------------*/
768 
769 void
771  cs_real_t time_eval,
772  cs_real_t *flx_values);
773 
774 /*----------------------------------------------------------------------------*/
786 /*----------------------------------------------------------------------------*/
787 
788 void
790  const cs_adv_field_t *adv,
791  short int f,
792  cs_real_t time_eval,
793  cs_real_t *fluxes);
794 
795 /*----------------------------------------------------------------------------*/
807 /*----------------------------------------------------------------------------*/
808 
809 cs_real_t
811  const short int f,
812  const cs_cell_mesh_t *cm,
813  const cs_adv_field_t *adv);
814 
815 /*----------------------------------------------------------------------------*/
825 /*----------------------------------------------------------------------------*/
826 
827 void
829  const cs_adv_field_t *adv,
830  cs_real_t time_eval,
831  cs_real_t *fluxes);
832 
833 /*----------------------------------------------------------------------------*/
843 /*----------------------------------------------------------------------------*/
844 
845 void
847  const cs_adv_field_t *adv,
848  cs_real_t time_eval,
849  cs_real_t *fluxes);
850 
851 /*----------------------------------------------------------------------------*/
859 /*----------------------------------------------------------------------------*/
860 
861 void
863  bool cur2prev);
864 
865 /*----------------------------------------------------------------------------*/
874 /*----------------------------------------------------------------------------*/
875 
876 void
878  const cs_property_t *diff,
879  cs_real_t t_eval,
880  cs_real_t peclet[]);
881 
882 /*----------------------------------------------------------------------------*/
890 /*----------------------------------------------------------------------------*/
891 
892 void
894  cs_real_t dt_cur,
895  cs_real_t courant[]);
896 
897 /*----------------------------------------------------------------------------*/
907 /*----------------------------------------------------------------------------*/
908 
909 cs_real_t *
911  cs_real_t t_eval);
912 
913 /*----------------------------------------------------------------------------*/
914 
916 
917 #endif /* __CS_ADVECTION_FIELD_H__ */
Definition: cs_advection_field.h:127
#define CS_FLAG_STATE_CELLWISE
Definition: cs_flag.h:88
Definition: cs_advection_field.h:138
#define restrict
Definition: cs_defs.h:142
void cs_advection_field_cw_face_flux(const cs_cell_mesh_t *cm, const cs_adv_field_t *adv, cs_real_t time_eval, cs_real_t *fluxes)
Compute the value of the flux of the advection field across the the (primal) faces of a cell...
Definition: cs_advection_field.c:2632
void cs_advection_field_def_by_dof_func(cs_adv_field_t *adv, cs_flag_t loc_flag, cs_dof_func_t *func, void *input)
Define a cs_adv_field_t structure using a cs_dof_func_t.
Definition: cs_advection_field.c:802
#define CS_FLAG_STATE_UNIFORM
Definition: cs_flag.h:87
void cs_advection_field_cw_dface_flux(const cs_cell_mesh_t *cm, const cs_adv_field_t *adv, cs_real_t time_eval, cs_real_t *fluxes)
Compute the value of the flux of the advection field across the the dual faces of a cell...
Definition: cs_advection_field.c:2826
Definition: cs_advection_field.h:149
void cs_advection_field_def_boundary_flux_by_analytic(cs_adv_field_t *adv, const char *zname, cs_analytic_func_t *func, void *input)
Define the value of the boundary normal flux for the given cs_adv_field_t structure using an analytic...
Definition: cs_advection_field.c:967
cs_advection_field_status_t status
Definition: cs_advection_field.h:197
cs_adv_field_t * cs_advection_field_by_id(int id)
Search in the array of advection field structures which one has the id given in argument.
Definition: cs_advection_field.c:428
cs_xdef_t ** bdy_flux_defs
Definition: cs_advection_field.h:211
Definition: cs_mesh_location.h:65
void cs_advection_field_def_by_array(cs_adv_field_t *adv, cs_flag_t loc, cs_real_t *array, bool is_owner, cs_lnum_t *index)
Define a cs_adv_field_t structure thanks to an array of values.
Definition: cs_advection_field.c:841
Field descriptor.
Definition: cs_field.h:125
cs_mesh_location_type_t
Definition: cs_mesh_location.h:60
void cs_advection_field_def_by_value(cs_adv_field_t *adv, cs_real_t vector[3])
Define the value of a cs_adv_field_t structure.
Definition: cs_advection_field.h:125
void cs_advection_field_def_by_analytic(cs_adv_field_t *adv, cs_analytic_func_t *func, void *input)
Define a cs_adv_field_t structure thanks to an analytic function.
Definition: cs_advection_field.c:767
cs_field_t * cs_field_by_id(int id)
Return a pointer to a field based on its id.
Definition: cs_field.c:2314
void cs_advection_get_courant(const cs_adv_field_t *adv, cs_real_t dt_cur, cs_real_t courant[])
Compute the Courant number in each cell.
Definition: cs_advection_field.c:3257
cs_real_t * cs_advection_field_divergence_at_vertices(const cs_adv_field_t *adv, cs_real_t t_eval)
Compute the divergence of the advection field at vertices Useful for CDO Vertex-based schemes...
Definition: cs_advection_field.c:3302
#define BEGIN_C_DECLS
Definition: cs_defs.h:510
void cs_advection_field_cw_eval_at_xyz(const cs_adv_field_t *adv, const cs_cell_mesh_t *cm, const cs_real_3_t xyz, cs_real_t time_eval, cs_nvec3_t *eval)
Compute the vector-valued interpolation of the advection field at a given location inside a cell...
Definition: cs_advection_field.c:1336
void cs_advection_field_finalize_setup(void)
Last stage of the definition of an advection field based on several definitions (i.e. definition by subdomains on the boundary)
Definition: cs_advection_field.c:1214
cs_real_t cs_advection_field_cw_boundary_face_flux(const cs_real_t time_eval, const short int f, const cs_cell_mesh_t *cm, const cs_adv_field_t *adv)
Compute the value of the normal flux of the advection field across a boundary face f (cellwise versio...
Definition: cs_advection_field.c:2158
int n_bdy_flux_defs
Definition: cs_advection_field.h:210
Definition: cs_advection_field.h:141
Set of local quantities and connectivities related to a mesh cell.
Definition: cs_cdo_local.h:202
void cs_advection_field_at_vertices(const cs_adv_field_t *adv, cs_real_t time_eval, cs_real_t *vtx_values)
Compute the value of the advection field at vertices.
Definition: cs_advection_field.c:1693
Definition: cs_cdo_connect.h:79
void cs_advection_field_in_cells(const cs_adv_field_t *adv, cs_real_t time_eval, cs_real_t *cell_values)
Compute the mean-value of the vector-valued field related to the advection field inside each cell...
Definition: cs_advection_field.c:1495
void cs_advection_field_cw_boundary_f2v_flux(const cs_cell_mesh_t *cm, const cs_adv_field_t *adv, short int f, cs_real_t time_eval, cs_real_t *fluxes)
Compute the value of the normal flux of the advection field across the closure of the dual cell relat...
Definition: cs_advection_field.c:2405
cs_flag_t state
Definition: cs_xdef.h:193
Definition: field.f90:27
Definition: cs_mesh_location.h:66
int id
Definition: cs_advection_field.h:195
void bft_error(const char *const file_name, const int line_num, const int sys_error_code, const char *const format,...)
Calls the error handler (set by bft_error_handler_set() or default).
Definition: bft_error.c:193
double cs_real_t
Floating-point value.
Definition: cs_defs.h:322
int cs_advection_field_get_n_fields(void)
Get the number of allocated cs_adv_field_t structures.
Definition: cs_advection_field.c:383
Definition: cs_cdo_quantities.h:129
Definition: cs_defs.h:368
void cs_advection_field_def_boundary_flux_by_field(cs_adv_field_t *adv, cs_field_t *field)
Define the value of the boundary normal flux for the given cs_adv_field_t structure using a field str...
Definition: cs_advection_field.c:1065
cs_flag_t cs_advection_field_status_t
Definition: cs_advection_field.h:68
void cs_advection_field_def_boundary_flux_by_array(cs_adv_field_t *adv, const char *zname, cs_flag_t loc, cs_real_t *array, bool is_owner, cs_lnum_t *index)
Define the value of the boundary normal flux for the given cs_adv_field_t structure using an array of...
Definition: cs_advection_field.c:1011
cs_xdef_t * definition
Definition: cs_advection_field.h:207
void cs_advection_field_set_postprocess(cs_adv_field_t *adv, cs_flag_t post_flag)
Set optional post-processings.
Definition: cs_advection_field.c:718
void() cs_dof_func_t(cs_lnum_t n_elts, const cs_lnum_t *elt_ids, bool dense_output, void *input, cs_real_t *retval)
Generic function pointer for computing a quantity at predefined locations such as degrees of freedom ...
Definition: cs_param_types.h:127
void cs_advection_field_get_cell_vector(cs_lnum_t c_id, const cs_adv_field_t *adv, cs_nvec3_t *vect)
Compute the value of the advection field at the cell center.
Definition: cs_advection_field.c:1305
void cs_advection_field_across_boundary(const cs_adv_field_t *adv, cs_real_t time_eval, cs_real_t *flx_values)
Compute the value of the normal flux of the advection field across the boundary faces.
Definition: cs_advection_field.c:1842
Definition: cs_mesh_location.h:64
short int * bdy_def_ids
Definition: cs_advection_field.h:212
Definition: cs_advection_field.h:140
int int_field_id
Definition: cs_advection_field.h:203
void() cs_analytic_func_t(cs_real_t time, cs_lnum_t n_elts, const cs_lnum_t *elt_ids, const cs_real_t *coords, bool dense_output, void *input, cs_real_t *retval)
Generic function pointer for an evaluation relying on an analytic function.
Definition: cs_param_types.h:100
cs_xdef_type_t
Definition: cs_xdef.h:108
int vtx_field_id
Definition: cs_advection_field.h:200
void cs_advection_field_update(cs_real_t t_eval, bool cur2prev)
For each cs_adv_field_t structures, update the values of the related field(s)
Definition: cs_advection_field.c:3133
Definition: cs_advection_field.h:126
void cs_advection_field_set_shared_pointers(const cs_cdo_quantities_t *quant, const cs_cdo_connect_t *connect)
Set shared pointers to main domain members.
Definition: cs_advection_field.c:366
Definition: cs_mesh_location.h:63
Structure storing medata for defining a quantity in a very flexible way.
Definition: cs_xdef.h:154
cs_real_t cs_real_3_t[3]
vector of 3 floating-point values
Definition: cs_defs.h:335
int bdy_field_id
Definition: cs_advection_field.h:202
cs_adv_field_t * cs_advection_field_by_name(const char *name)
Search in the array of advection field structures which one has the name given in argument...
Definition: cs_advection_field.c:400
cs_adv_field_t * cs_advection_field_add(const char *name, cs_advection_field_status_t status)
Add and initialize a new advection field structure.
Definition: cs_advection_field.c:468
Definition: cs_advection_field.h:139
cs_adv_field_t * cs_advection_field_add_user(const char *name)
Add and initialize a new user-defined advection field structure.
Definition: cs_advection_field.c:451
int cs_lnum_t
local mesh entity id
Definition: cs_defs.h:316
cs_xdef_type_t cs_xdef_get_type(const cs_xdef_t *d)
Retrieve the flag dedicated to the state.
Definition: cs_xdef.c:853
void cs_advection_field_log_setup(void)
Print all setup information related to cs_adv_field_t structures.
Definition: cs_advection_field.c:620
void cs_advection_field_create_fields(void)
Create all needed cs_field_t structures related to an advection field.
Definition: cs_advection_field.c:1104
int cell_field_id
Definition: cs_advection_field.h:201
#define END_C_DECLS
Definition: cs_defs.h:511
void cs_advection_field_destroy_all(void)
Free all alllocated cs_adv_field_t structures and its related array.
Definition: cs_advection_field.c:554
unsigned short int cs_flag_t
Definition: cs_defs.h:324
Definition: cs_advection_field.h:133
void cs_advection_field_def_boundary_flux_by_value(cs_adv_field_t *adv, const char *zname, cs_real_t normal_flux)
Define the value of the boundary normal flux for the given cs_adv_field_t structure.
Definition: cs_advection_field.c:931
cs_advection_field_status_bit_t
Bit values for specifying the definition/behavior of an advection field.
Definition: cs_advection_field.h:120
char *restrict name
Definition: cs_advection_field.h:196
cs_flag_t post_flag
Definition: cs_advection_field.h:198
Definition: cs_xdef.h:120
Structure associated to the definition of a property relying on the cs_xdef_t structure.
Definition: cs_advection_field.h:132
void cs_advection_field_def_by_field(cs_adv_field_t *adv, cs_field_t *field)
Define a cs_adv_field_t structure thanks to a field structure.
Definition: cs_advection_field.c:893
void cs_advection_get_peclet(const cs_adv_field_t *adv, const cs_property_t *diff, cs_real_t t_eval, cs_real_t peclet[])
Compute the Peclet number in each cell.
Definition: cs_advection_field.c:3209
bool cs_advection_field_check_name(const cs_adv_field_t *adv, const char *ref_name)
Check if the given advection field has the name ref_name.
Definition: cs_advection_field.c:594