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 
46 BEGIN_C_DECLS
47 
48 /*============================================================================
49  * Macro definitions
50  *============================================================================*/
51 
52 /*!
53  * @defgroup cdo_advfield_flags Flags specifying metadata related to an
54  *  advection field.
55  *
56  * @{
57  */
58 
59 /*!  1: Perform the computation and post-processing of the Courant number */
60 #define CS_ADVECTION_FIELD_POST_COURANT (1 << 0)
61 
62 /*! @} */
63 
64 /*============================================================================
65  * Type definitions
66  *============================================================================*/
67 
68 typedef cs_flag_t  cs_advection_field_status_t;
69 
70 /*! \enum cs_advection_field_status_bit_t
71  *  \brief Bit values for specifying the definition/behavior of an advection
72  *         field
73  *
74  * @name Category of advection field
75  * @{
76  *
77  * \var CS_ADVECTION_FIELD_NAVSTO
78  *      Advection field stemming from the velocity in the (Navier-)Stokes system
79  *
80  * \var CS_ADVECTION_FIELD_GWF
81  *      Advection field stemming from the "GroundWater Flows" module. This is the
82  *      Darcean flux.
83  *
84  * \var CS_ADVECTION_FIELD_USER
85  *      User-defined advection field.
86  *
87  * @}
88  * @name Type of definition
89  * @{
90  *
91  * \var CS_ADVECTION_FIELD_TYPE_VELOCITY_VECTOR
92  *      The advection field stands for a velocity.
93  *      This is described by a vector-valued array or function.
94  *
95  * \var CS_ADVECTION_FIELD_TYPE_SCALAR_FLUX
96  *      The advection field stands for a flux.
97  *      This is described by a scalar-valued array or function.
98  *
99  * @}
100  * @name Optional bits
101  * @{
102  *
103  * \var CS_ADVECTION_FIELD_STEADY
104  *      Advection field is steady-state
105  *
106  * \var CS_ADVECTION_FIELD_LEGACY_FV
107  *      Advection field shared with the legacy Finite Volume solver
108  *
109  * \var CS_ADVECTION_FIELD_DEFINE_AT_VERTICES
110  *      Define a field structure related to the advection field at vertices
111  *      Post-processing and log operations are automatically activated.
112  *
113  * \var CS_ADVECTION_FIELD_DEFINE_AT_BOUNDARY_FACES
114  *      Define a field structure related to the advection field at boundary
115  *      faces Post-processing and log operations are automatically activated.
116  *
117  * @}
118  */
119 
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 
132   CS_ADVECTION_FIELD_TYPE_VELOCITY_VECTOR              = 1<<3, /* =   8 */
133   CS_ADVECTION_FIELD_TYPE_SCALAR_FLUX                  = 1<<4, /* =  16 */
134 
135   /* Optional */
136   /* -------- */
137 
138   CS_ADVECTION_FIELD_STEADY                            = 1<<5, /* =  32 */
139   CS_ADVECTION_FIELD_LEGACY_FV                         = 1<<6, /* =  64 */
140   CS_ADVECTION_FIELD_DEFINE_AT_VERTICES                = 1<<7, /* = 128 */
141   CS_ADVECTION_FIELD_DEFINE_AT_BOUNDARY_FACES          = 1<<8  /* = 256 */
142 
143 } cs_advection_field_status_bit_t;
144 
145 /*! \struct cs_advection_field_t
146  *  \brief Main structure to handle an advection field
147  */
148 
149 typedef struct {
150 
151   /*!
152    * \var id
153    * identification number
154    *
155    * \var name
156    * name of the advection field
157    *
158    * \var status
159    * category (user, navsto, gwf...) and type (velocity, flux...) of the
160    * advection field
161    *
162    * \var post_flag
163    * short descriptor dedicated to postprocessing
164    *
165    * \var vtx_field_id
166    * id to retrieve the related cs_field_t structure (-1 if not used)
167    *
168    * \var cell_field_id
169    * id to retrieve the related cs_field_t structure. It's always
170    * defined since it's used during the building of the advection scheme
171    *
172    * \var bdy_field_id
173    * id to retrieve the related cs_field_t structure corresponding to
174    * the value of the normal flux across boundary faces. It's always
175    * defined since it's used for dealing with the boundary conditions.
176    *
177    * \var int_field_id
178    * id to retrieve the related cs_field_t structure corresponding to
179    * the value of the normal flux across interior faces. It's always
180    * defined when the advection field arise from the resolution of the
181    * Navier-Stokes system with the legacy FV solver.
182    *
183    * \var definition
184    * pointer to the generic structure used to define the advection field. We
185    * assume that only one definition is available (i.e. there is a unique
186    * zone related to "cells").
187    *
188    * \var n_bdy_flux_defs
189    * Number of definitions related to the normal flux at the boundary
190    *
191    * \var bdy_flux_defs
192    * Array of pointers to the definitions of the normal flux at the boundary
193    */
194 
195   int                           id;
196   char                *restrict name;
197   cs_advection_field_status_t   status;
198   cs_flag_t                     post_flag;
199 
200   int                           vtx_field_id;
201   int                           cell_field_id;
202   int                           bdy_field_id;
203   int                           int_field_id;
204 
205   /* We assume that there is only one definition associated to an advection
206      field inside the computational domain */
207   cs_xdef_t                    *definition;
208 
209   /* Optional: Definition(s) for the boundary flux */
210   int                           n_bdy_flux_defs;
211   cs_xdef_t                   **bdy_flux_defs;
212   short int                    *bdy_def_ids;
213 
214 } cs_adv_field_t;
215 
216 /*============================================================================
217  * Global variables
218  *============================================================================*/
219 
220 /*============================================================================
221  * Inline public function prototypes
222  *============================================================================*/
223 
224 /*----------------------------------------------------------------------------*/
225 /*!
226  * \brief  Set a new status for the given advection field structure
227  *
228  * \param[in, out] adv      pointer to an advection field structure
229  * \param[in]      status   status flag to add to the current status
230  */
231 /*----------------------------------------------------------------------------*/
232 
233 static inline void
cs_advection_field_set_status(cs_adv_field_t * adv,cs_advection_field_status_t status)234 cs_advection_field_set_status(cs_adv_field_t               *adv,
235                               cs_advection_field_status_t   status)
236 {
237   if (adv == NULL)
238     return;
239 
240   adv->status = status;
241 }
242 
243 /*----------------------------------------------------------------------------*/
244 /*!
245  * \brief  returns true if the advection field is uniform, otherwise false
246  *
247  * \param[in]    adv    pointer to an advection field to test
248  *
249  * \return  true or false
250  */
251 /*----------------------------------------------------------------------------*/
252 
253 static inline bool
cs_advection_field_is_uniform(const cs_adv_field_t * adv)254 cs_advection_field_is_uniform(const cs_adv_field_t   *adv)
255 {
256   if (adv == NULL)
257     return false;
258 
259   if (adv->definition->state & CS_FLAG_STATE_UNIFORM)
260     return true;
261   else
262     return false;
263 }
264 
265 /*----------------------------------------------------------------------------*/
266 /*!
267  * \brief  returns true if the advection field is uniform in each cell
268  *         otherwise false
269  *
270  * \param[in]    adv    pointer to an advection field to test
271  *
272  * \return  true or false
273  */
274 /*----------------------------------------------------------------------------*/
275 
276 static inline bool
cs_advection_field_is_cellwise(const cs_adv_field_t * adv)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 /*----------------------------------------------------------------------------*/
293 /*!
294  * \brief  Retrieve the name of an advection field
295  *
296  * \param[in]    adv    pointer to an advection field structure
297  *
298  * \return  the name of the related advection field
299  */
300 /*----------------------------------------------------------------------------*/
301 
302 static inline const char *
cs_advection_field_get_name(const cs_adv_field_t * adv)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 /*----------------------------------------------------------------------------*/
312 /*!
313  * \brief  Retrieve the type of definition used to set the current advection
314  *         field structure
315  *
316  * \param[in]    adv    pointer to an advection field structure
317  *
318  * \return  the type of definition
319  */
320 /*----------------------------------------------------------------------------*/
321 
322 static inline cs_xdef_type_t
cs_advection_field_get_deftype(const cs_adv_field_t * adv)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 /*----------------------------------------------------------------------------*/
332 /*!
333  * \brief  Get a cs_field_t structure related to an advection field and a mesh
334  *         location
335  *
336  * \param[in]  adv         pointer to a cs_adv_field_t structure
337  * \param[in]  ml_type     type of mesh location (cells or vertices)
338  *
339  * \return a pointer to a cs_field_t structure
340  */
341 /*----------------------------------------------------------------------------*/
342 
343 static inline cs_field_t *
cs_advection_field_get_field(const cs_adv_field_t * adv,cs_mesh_location_type_t ml_type)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 
353   case CS_MESH_LOCATION_CELLS:
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 
360   case CS_MESH_LOCATION_INTERIOR_FACES:
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 
367   case CS_MESH_LOCATION_BOUNDARY_FACES:
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 
374   case CS_MESH_LOCATION_VERTICES:
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 /*----------------------------------------------------------------------------*/
396 /*!
397  * \brief  Set shared pointers to main domain members
398  *
399  * \param[in]  quant       additional mesh quantities struct.
400  * \param[in]  connect     pointer to a cs_cdo_connect_t struct.
401  */
402 /*----------------------------------------------------------------------------*/
403 
404 void
405 cs_advection_field_set_shared_pointers(const cs_cdo_quantities_t  *quant,
406                                        const cs_cdo_connect_t     *connect);
407 
408 /*----------------------------------------------------------------------------*/
409 /*!
410  * \brief  Get the number of allocated cs_adv_field_t structures
411  *
412  * \return the number of advection fields
413  */
414 /*----------------------------------------------------------------------------*/
415 
416 int
417 cs_advection_field_get_n_fields(void);
418 
419 /*----------------------------------------------------------------------------*/
420 /*!
421  * \brief  Search in the array of advection field structures which one has
422  *         the name given in argument
423  *
424  * \param[in]  name        name of the advection field
425  *
426  * \return a pointer to a cs_adv_field_t structure or NULL if not found
427  */
428 /*----------------------------------------------------------------------------*/
429 
430 cs_adv_field_t *
431 cs_advection_field_by_name(const char   *name);
432 
433 /*----------------------------------------------------------------------------*/
434 /*!
435  * \brief  Search in the array of advection field structures which one has
436  *         the id given in argument
437  *
438  * \param[in]  id        identification number
439  *
440  * \return a pointer to a cs_adv_field_t structure or NULL if not found
441  */
442 /*----------------------------------------------------------------------------*/
443 
444 cs_adv_field_t *
445 cs_advection_field_by_id(int      id);
446 
447 /*----------------------------------------------------------------------------*/
448 /*!
449  * \brief  Add and initialize a new user-defined advection field structure
450  *
451  * \param[in]  name        name of the advection field
452  *
453  * \return a pointer to the new allocated cs_adv_field_t structure
454  */
455 /*----------------------------------------------------------------------------*/
456 
457 cs_adv_field_t *
458 cs_advection_field_add_user(const char  *name);
459 
460 /*----------------------------------------------------------------------------*/
461 /*!
462  * \brief  Add and initialize a new advection field structure
463  *
464  * \param[in]  name        name of the advection field
465  * \param[in]  status      status of the advection field to create
466  *
467  * \return a pointer to the new allocated cs_adv_field_t structure
468  */
469 /*----------------------------------------------------------------------------*/
470 
471 cs_adv_field_t *
472 cs_advection_field_add(const char                    *name,
473                        cs_advection_field_status_t    status);
474 
475 /*----------------------------------------------------------------------------*/
476 /*!
477  * \brief  Free all alllocated cs_adv_field_t structures and its related array
478  *
479  * \return a NULL pointer
480  */
481 /*----------------------------------------------------------------------------*/
482 
483 void
484 cs_advection_field_destroy_all(void);
485 
486 /*----------------------------------------------------------------------------*/
487 /*!
488  * \brief  Check if the given advection field has the name ref_name
489  *
490  * \param[in]  adv         pointer to a cs_adv_field_t structure to test
491  * \param[in]  ref_name    name of the advection field to find
492  *
493  * \return true if the name of the advection field is ref_name otherwise false
494  */
495 /*----------------------------------------------------------------------------*/
496 
497 bool
498 cs_advection_field_check_name(const cs_adv_field_t   *adv,
499                               const char             *ref_name);
500 
501 /*----------------------------------------------------------------------------*/
502 /*!
503  * \brief  Print all setup information related to cs_adv_field_t structures
504  */
505 /*----------------------------------------------------------------------------*/
506 
507 void
508 cs_advection_field_log_setup(void);
509 
510 /*----------------------------------------------------------------------------*/
511 /*!
512  * \brief  Set optional post-processings
513  *
514  * \param[in, out]  adv         pointer to a cs_adv_field_t structure
515  * \param[in]       post_flag   flag to set
516  */
517 /*----------------------------------------------------------------------------*/
518 
519 void
520 cs_advection_field_set_postprocess(cs_adv_field_t            *adv,
521                                    cs_flag_t                  post_flag);
522 
523 /*----------------------------------------------------------------------------*/
524 /*!
525  * \brief  Define the value of a cs_adv_field_t structure
526  *
527  * \param[in, out]  adv       pointer to a cs_adv_field_t structure
528  * \param[in]       vector    values to set
529  */
530 /*----------------------------------------------------------------------------*/
531 
532 void
533 cs_advection_field_def_by_value(cs_adv_field_t    *adv,
534                                 cs_real_t          vector[3]);
535 
536 /*----------------------------------------------------------------------------*/
537 /*!
538  * \brief  Define a cs_adv_field_t structure thanks to an analytic function
539  *
540  * \param[in, out]  adv     pointer to a cs_adv_field_t structure
541  * \param[in]       func    pointer to a function
542  * \param[in]       input   NULL or pointer to a structure cast on-the-fly
543  */
544 /*----------------------------------------------------------------------------*/
545 
546 void
547 cs_advection_field_def_by_analytic(cs_adv_field_t        *adv,
548                                    cs_analytic_func_t    *func,
549                                    void                  *input);
550 
551 /*----------------------------------------------------------------------------*/
552 /*!
553  * \brief  Define a cs_adv_field_t structure using a cs_dof_func_t
554  *
555  * \param[in, out] adv       pointer to a cs_adv_field_t structure
556  * \param[in]      loc_flag  location where values are computed
557  * \param[in]      func      pointer to a cs_dof_func_t function
558  * \param[in]      input     NULL or pointer to a structure cast on-the-fly
559  */
560 /*----------------------------------------------------------------------------*/
561 
562 void
563 cs_advection_field_def_by_dof_func(cs_adv_field_t    *adv,
564                                    cs_flag_t          loc_flag,
565                                    cs_dof_func_t     *func,
566                                    void              *input);
567 
568 /*----------------------------------------------------------------------------*/
569 /*!
570  * \brief  Define a cs_adv_field_t structure thanks to an array of values
571  *
572  * \param[in, out]  adv       pointer to a cs_adv_field_t structure
573  * \param[in]       loc       information to know where are located values
574  * \param[in]       array     pointer to an array
575  * \param[in]       is_owner  transfer the lifecycle to the cs_xdef_t structure
576  *                            (true or false)
577  * \param[in]       index     optional pointer to the array index
578  */
579 /*----------------------------------------------------------------------------*/
580 
581 void
582 cs_advection_field_def_by_array(cs_adv_field_t    *adv,
583                                 cs_flag_t          loc,
584                                 cs_real_t         *array,
585                                 bool               is_owner,
586                                 cs_lnum_t         *index);
587 
588 /*----------------------------------------------------------------------------*/
589 /*!
590  * \brief  Define a cs_adv_field_t structure thanks to a field structure
591  *
592  * \param[in, out]  adv       pointer to a cs_adv_field_t structure
593  * \param[in]       field     pointer to a cs_field_t structure
594  */
595 /*----------------------------------------------------------------------------*/
596 
597 void
598 cs_advection_field_def_by_field(cs_adv_field_t    *adv,
599                                 cs_field_t        *field);
600 
601 /*----------------------------------------------------------------------------*/
602 /*!
603  * \brief  Define the value of the boundary normal flux for the given
604  *         cs_adv_field_t structure
605  *
606  * \param[in, out]  adv           pointer to a cs_adv_field_t structure
607  * \param[in]       zname         name of the boundary zone to consider
608  * \param[in]       normal_flux   value to set
609  */
610 /*----------------------------------------------------------------------------*/
611 
612 void
613 cs_advection_field_def_boundary_flux_by_value(cs_adv_field_t    *adv,
614                                               const char        *zname,
615                                               cs_real_t          normal_flux);
616 
617 /*----------------------------------------------------------------------------*/
618 /*!
619  * \brief  Define the value of the boundary normal flux for the given
620  *         cs_adv_field_t structure using an analytic function
621  *
622  * \param[in, out]  adv     pointer to a cs_adv_field_t structure
623  * \param[in]       zname   name of the boundary zone to consider
624  * \param[in]       func    pointer to a function
625  * \param[in]       input   NULL or pointer to a structure cast on-the-fly
626  */
627 /*----------------------------------------------------------------------------*/
628 
629 void
630 cs_advection_field_def_boundary_flux_by_analytic(cs_adv_field_t        *adv,
631                                                  const char            *zname,
632                                                  cs_analytic_func_t    *func,
633                                                  void                  *input);
634 
635 /*----------------------------------------------------------------------------*/
636 /*!
637  * \brief  Define the value of the boundary normal flux for the given
638  *         cs_adv_field_t structure using an array of values
639  *
640  * \param[in, out]  adv       pointer to a cs_adv_field_t structure
641  * \param[in]       zname     name of the boundary zone to consider
642  * \param[in]       loc       information to know where are located values
643  * \param[in]       array     pointer to an array
644  * \param[in]       is_owner  transfer the lifecycle to the cs_xdef_t structure
645  *                            (true or false)
646  * \param[in]       index     optional pointer to the array index
647  */
648 /*----------------------------------------------------------------------------*/
649 
650 void
651 cs_advection_field_def_boundary_flux_by_array(cs_adv_field_t    *adv,
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 /*----------------------------------------------------------------------------*/
659 /*!
660  * \brief  Define the value of the boundary normal flux for the given
661  *         cs_adv_field_t structure using a field structure
662  *
663  * \param[in, out]  adv       pointer to a cs_adv_field_t structure
664  * \param[in]       field     pointer to a cs_field_t structure
665  */
666 /*----------------------------------------------------------------------------*/
667 
668 void
669 cs_advection_field_def_boundary_flux_by_field(cs_adv_field_t    *adv,
670                                               cs_field_t        *field);
671 
672 /*----------------------------------------------------------------------------*/
673 /*!
674  * \brief  Create all needed cs_field_t structures related to an advection
675  *         field
676  */
677 /*----------------------------------------------------------------------------*/
678 
679 void
680 cs_advection_field_create_fields(void);
681 
682 /*----------------------------------------------------------------------------*/
683 /*!
684  * \brief  Last stage of the definition of an advection field based on several
685  *         definitions (i.e. definition by subdomains on the boundary)
686  */
687 /*----------------------------------------------------------------------------*/
688 
689 void
690 cs_advection_field_finalize_setup(void);
691 
692 /*----------------------------------------------------------------------------*/
693 /*!
694  * \brief  Compute the value of the advection field at the cell center
695  *
696  * \param[in]      c_id    id of the current cell
697  * \param[in]      adv     pointer to a cs_adv_field_t structure
698  * \param[in, out] vect    pointer to a cs_nvec3_t structure (meas + unitv)
699  */
700 /*----------------------------------------------------------------------------*/
701 
702 void
703 cs_advection_field_get_cell_vector(cs_lnum_t               c_id,
704                                    const cs_adv_field_t   *adv,
705                                    cs_nvec3_t             *vect);
706 
707 /*----------------------------------------------------------------------------*/
708 /*!
709  * \brief  Compute the vector-valued interpolation of the advection field at
710  *         a given location inside a cell
711  *
712  * \param[in]      adv          pointer to a cs_adv_field_t structure
713  * \param[in]      cm           pointer to a cs_cell_mesh_t structure
714  * \param[in]      xyz          location where to perform the evaluation
715  * \param[in]      time_eval    physical time at which one evaluates the term
716  * \param[in, out] eval         pointer to a cs_nvec3_t
717  */
718 /*----------------------------------------------------------------------------*/
719 
720 void
721 cs_advection_field_cw_eval_at_xyz(const cs_adv_field_t  *adv,
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 /*----------------------------------------------------------------------------*/
728 /*!
729  * \brief  Compute the mean-value of the vector-valued field related to the
730  *         advection field inside each cell
731  *
732  * \param[in]      adv           pointer to a cs_adv_field_t structure
733  * \param[in]      time_eval     physical time at which one evaluates the term
734  * \param[in, out] cell_values   array of values at cell centers
735  */
736 /*----------------------------------------------------------------------------*/
737 
738 void
739 cs_advection_field_in_cells(const cs_adv_field_t  *adv,
740                             cs_real_t              time_eval,
741                             cs_real_t             *cell_values);
742 
743 /*----------------------------------------------------------------------------*/
744 /*!
745  * \brief  Compute the value of the advection field at vertices
746  *
747  * \param[in]      adv          pointer to a cs_adv_field_t structure
748  * \param[in]      time_eval    physical time at which one evaluates the term
749  * \param[in, out] vtx_values   array storing the results
750  */
751 /*----------------------------------------------------------------------------*/
752 
753 void
754 cs_advection_field_at_vertices(const cs_adv_field_t  *adv,
755                                cs_real_t              time_eval,
756                                cs_real_t             *vtx_values);
757 
758 /*----------------------------------------------------------------------------*/
759 /*!
760  * \brief  Compute the value of the normal flux of the advection field
761  *         across the boundary faces
762  *
763  * \param[in]      adv          pointer to a cs_adv_field_t structure
764  * \param[in]      time_eval    physical time at which one evaluates the term
765  * \param[in, out] flx_values   array storing the results
766  */
767 /*----------------------------------------------------------------------------*/
768 
769 void
770 cs_advection_field_across_boundary(const cs_adv_field_t  *adv,
771                                    cs_real_t              time_eval,
772                                    cs_real_t             *flx_values);
773 
774 /*----------------------------------------------------------------------------*/
775 /*!
776  * \brief  Compute the value of the normal flux of the advection field across
777  *         the closure of the dual cell related to each vertex belonging to the
778  *         boundary face f
779  *
780  * \param[in]      cm         pointer to a cs_cell_mesh_t structure
781  * \param[in]      adv        pointer to a cs_adv_field_t structure
782  * \param[in]      f          face id in the cellwise numbering
783  * \param[in]      time_eval  physical time at which one evaluates the term
784  * \param[in, out] fluxes     normal boundary flux for each vertex of the face
785  */
786 /*----------------------------------------------------------------------------*/
787 
788 void
789 cs_advection_field_cw_boundary_f2v_flux(const cs_cell_mesh_t   *cm,
790                                         const cs_adv_field_t   *adv,
791                                         short int               f,
792                                         cs_real_t               time_eval,
793                                         cs_real_t              *fluxes);
794 
795 /*----------------------------------------------------------------------------*/
796 /*!
797  * \brief  Compute the value of the normal flux of the advection field across
798  *         a boundary face f (cellwise version)
799  *
800  * \param[in] time_eval   physical time at which one evaluates the term
801  * \param[in] f           face id in the cellwise numbering
802  * \param[in] cm          pointer to a cs_cell_mesh_t structure
803  * \param[in] adv         pointer to a cs_adv_field_t structure
804  *
805  * \return  the normal boundary flux for the face f
806  */
807 /*----------------------------------------------------------------------------*/
808 
809 cs_real_t
810 cs_advection_field_cw_boundary_face_flux(const cs_real_t          time_eval,
811                                          const short int          f,
812                                          const cs_cell_mesh_t    *cm,
813                                          const cs_adv_field_t    *adv);
814 
815 /*----------------------------------------------------------------------------*/
816 /*!
817  * \brief  Compute the value of the flux of the advection field across the
818  *         the (primal) faces of a cell
819  *
820  * \param[in]      cm         pointer to a cs_cell_mesh_t structure
821  * \param[in]      adv        pointer to a cs_adv_field_t structure
822  * \param[in]      time_eval  physical time at which one evaluates the term
823  * \param[in, out] fluxes     array of values attached to primal faces of a cell
824  */
825 /*----------------------------------------------------------------------------*/
826 
827 void
828 cs_advection_field_cw_face_flux(const cs_cell_mesh_t       *cm,
829                                 const cs_adv_field_t       *adv,
830                                 cs_real_t                   time_eval,
831                                 cs_real_t                  *fluxes);
832 
833 /*----------------------------------------------------------------------------*/
834 /*!
835  * \brief  Compute the value of the flux of the advection field across the
836  *         the dual faces of a cell
837  *
838  * \param[in]      cm         pointer to a cs_cell_mesh_t structure
839  * \param[in]      adv        pointer to a cs_adv_field_t structure
840  * \param[in]      time_eval  physical time at which one evaluates the term
841  * \param[in, out] fluxes     array of values attached to dual faces of a cell
842  */
843 /*----------------------------------------------------------------------------*/
844 
845 void
846 cs_advection_field_cw_dface_flux(const cs_cell_mesh_t     *cm,
847                                  const cs_adv_field_t     *adv,
848                                  cs_real_t                 time_eval,
849                                  cs_real_t                *fluxes);
850 
851 /*----------------------------------------------------------------------------*/
852 /*!
853  * \brief   For each cs_adv_field_t structures, update the values of the
854  *          related field(s)
855  *
856  * \param[in]  t_eval     physical time at which one evaluates the term
857  * \param[in]  cur2prev   true or false
858  */
859 /*----------------------------------------------------------------------------*/
860 
861 void
862 cs_advection_field_update(cs_real_t    t_eval,
863                           bool         cur2prev);
864 
865 /*----------------------------------------------------------------------------*/
866 /*!
867  * \brief   Compute the Peclet number in each cell
868  *
869  * \param[in]      adv        pointer to the advection field struct.
870  * \param[in]      diff       pointer to the diffusion property struct.
871  * \param[in]      t_eval     time at which one evaluates the advection field
872  * \param[in, out] peclet     pointer to an array storing the Peclet number
873  */
874 /*----------------------------------------------------------------------------*/
875 
876 void
877 cs_advection_get_peclet(const cs_adv_field_t     *adv,
878                         const cs_property_t      *diff,
879                         cs_real_t                 t_eval,
880                         cs_real_t                 peclet[]);
881 
882 /*----------------------------------------------------------------------------*/
883 /*!
884  * \brief   Compute the Courant number in each cell
885  *
886  * \param[in]      adv        pointer to the advection field struct.
887  * \param[in]      dt_cur     current time step
888  * \param[in, out] courant    pointer to an array storing the Courant number
889  */
890 /*----------------------------------------------------------------------------*/
891 
892 void
893 cs_advection_get_courant(const cs_adv_field_t     *adv,
894                          cs_real_t                 dt_cur,
895                          cs_real_t                 courant[]);
896 
897 /*----------------------------------------------------------------------------*/
898 /*!
899  * \brief   Compute the divergence of the advection field at vertices
900  *          Useful for CDO Vertex-based schemes
901  *
902  * \param[in]      adv         pointer to the advection field struct.
903  * \param[in]      t_eval      time at which one evaluates the advection field
904  *
905  * \return a pointer to an array storing the result
906  */
907 /*----------------------------------------------------------------------------*/
908 
909 cs_real_t *
910 cs_advection_field_divergence_at_vertices(const cs_adv_field_t     *adv,
911                                           cs_real_t                 t_eval);
912 
913 /*----------------------------------------------------------------------------*/
914 
915 END_C_DECLS
916 
917 #endif /* __CS_ADVECTION_FIELD_H__ */
918