1 /*============================================================================
2  * User initialization prior to solving time steps.
3  *============================================================================*/
4 
5 /* VERS */
6 
7 /*
8   This file is part of Code_Saturne, a general-purpose CFD tool.
9 
10   Copyright (C) 1998-2021 EDF S.A.
11 
12   This program is free software; you can redistribute it and/or modify it under
13   the terms of the GNU General Public License as published by the Free Software
14   Foundation; either version 2 of the License, or (at your option) any later
15   version.
16 
17   This program is distributed in the hope that it will be useful, but WITHOUT
18   ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
19   FOR A PARTICULAR PURPOSE.  See the GNU General Public License for more
20   details.
21 
22   You should have received a copy of the GNU General Public License along with
23   this program; if not, write to the Free Software Foundation, Inc., 51 Franklin
24   Street, Fifth Floor, Boston, MA 02110-1301, USA.
25 */
26 
27 /*----------------------------------------------------------------------------*/
28 
29 #include "cs_defs.h"
30 
31 /*----------------------------------------------------------------------------
32  * Standard C library headers
33  *----------------------------------------------------------------------------*/
34 
35 #include <assert.h>
36 #include <math.h>
37 
38 #if defined(HAVE_MPI)
39 #include <mpi.h>
40 #endif
41 
42 /*----------------------------------------------------------------------------
43  * PLE library headers
44  *----------------------------------------------------------------------------*/
45 
46 #include <ple_coupling.h>
47 
48 /*----------------------------------------------------------------------------
49  * Local headers
50  *----------------------------------------------------------------------------*/
51 
52 #include "cs_headers.h"
53 
54 /*----------------------------------------------------------------------------*/
55 
56 BEGIN_C_DECLS
57 
58 /*----------------------------------------------------------------------------*/
59 /*!
60  * \file cs_user_initialization-electric-arcs.c
61  *
62  * \brief Electric arcs fields initialization example.
63  *
64  * See \ref cs_user_initialization for examples.
65  */
66 /*----------------------------------------------------------------------------*/
67 
68 /*============================================================================
69  * User function definitions
70  *============================================================================*/
71 
72 /*----------------------------------------------------------------------------*/
73 /*!
74  * \brief Initialize variables.
75  *
76  * This function is called at beginning of the computation
77  * (restart or not) before the time step loop.
78  *
79  * This is intended to initialize or modify (when restarted)
80  * variable and time step values.
81  *
82  * \param[in, out]  domain   pointer to a cs_domain_t structure
83  */
84 /*----------------------------------------------------------------------------*/
85 
86 void
cs_user_initialization(cs_domain_t * domain)87 cs_user_initialization(cs_domain_t     *domain)
88 {
89   /*! [loc_var_dec] */
90   const cs_lnum_t n_cells = domain->mesh-> n_cells;
91   const cs_real_3_t *cell_cen
92     = (const cs_real_3_t *)domain->mesh_quantities->cell_cen;
93 
94   const cs_real_t hhot = 0.4075;
95   /*! [loc_var_dec] */
96 
97   /*! [init2] */
98   /* Apply only at the true computation start, not on restarts */
99   if (domain->time_step->nt_prev > 0)
100     return;
101 
102   /* Mass fraction = 1 for first gas, 0 for others */
103 
104   int n_gasses = cs_glob_elec_properties->ngaz;
105 
106   if (n_gasses > 1) {
107     cs_array_set_value_real(n_cells, 1, 1.0, CS_FI_(ycoel, 0)->val);
108 
109     for (int sp_id = 1; sp_id < n_gasses-1; sp_id++) {
110       cs_array_set_value_real(n_cells, 1, 0.0, CS_FI_(ycoel, sp_id)->val);
111     }
112   }
113 
114   /* Enthalpy = H(T0) or 0
115 
116      For electric arc,
117      for the whole compution domain enthalpy is set to H(T) using
118      the model conversion function.
119 
120      For Joule jeating by direct conduction,
121      enthalpy is set to zero, and the user will enter his H(T) function
122      tabulation.
123   */
124 
125   /* Enthaly reinitialization */
126 
127   cs_real_t  *cpro_t = CS_F_(t)->val;
128   cs_real_t  *cvar_h = CS_F_(h)->val;
129 
130   if (cs_glob_physical_model_flag[CS_ELECTRIC_ARCS] >= 1) {
131     cs_elec_convert_t_to_h_cells(cpro_t, cvar_h);
132   }
133   else {
134     cs_user_physical_properties_t_to_h(domain,
135                                        cs_volume_zone_by_id(0),
136                                        false,
137                                        cpro_t,
138                                        cvar_h);
139   }
140 
141   /* Enthalpy modification in given area */
142 
143   for (cs_lnum_t i = 0; i < n_cells; i++) {
144     cs_real_t rayo = sqrt(  cell_cen[i][0]*cell_cen[i][0]
145                           + cell_cen[i][1]*cell_cen[i][1]);
146     if (rayo <= 0.8e-3)
147       cvar_h[i] = hhot;
148   }
149 
150   /* Set electric potential to 0 */
151 
152   /* real component */
153   cs_array_set_value_real(n_cells, 1, 0.0, CS_F_(potr)->val);
154 
155   /* imaginary component (for Joule heating by direct conduction) */
156   if (CS_F_(poti) != NULL)
157     cs_array_set_value_real(n_cells, 1, 0.0, CS_F_(poti)->val);
158 
159   /* Vector potential (electric arcs, 3d) */
160   if (CS_F_(potva) != NULL)
161     cs_array_set_value_real(n_cells, 3, 0.0, CS_F_(potva)->val);
162   /*! [init2] */
163 }
164 
165 /*----------------------------------------------------------------------------*/
166 
167 END_C_DECLS
168