9.1
general documentation
cs_ale.h
Go to the documentation of this file.
1#ifndef __CS_ALE_H__
2#define __CS_ALE_H__
3
4/*============================================================================
5 * Functions associated to ALE formulation
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/*----------------------------------------------------------------------------
31 * Standard C library headers
32 *----------------------------------------------------------------------------*/
33
34/*----------------------------------------------------------------------------
35 * Local headers
36 *----------------------------------------------------------------------------*/
37
38#include "base/cs_array.h"
39#include "base/cs_base.h"
40#include "cdo/cs_domain.h"
41#include "base/cs_restart.h"
42
43/*----------------------------------------------------------------------------*/
44
46
47/*============================================================================
48 * Type definitions
49 *============================================================================*/
50
51/*----------------------------------------------------------------------------
52 * ALE type
53 *----------------------------------------------------------------------------*/
54
57typedef enum {
58
61 CS_ALE_CDO = 2
64
67typedef struct {
68
69 int *impale;
70 int *bc_type;
76 /* For neptune ALE implicit coupling */
88
89/*=============================================================================
90 * Global variables
91 *============================================================================*/
92
94
96
97extern int cs_glob_ale_n_ini_f; /* Number of sub-iterations for fluid
98 flow initialization */
99
100extern int cs_glob_ale_need_init; /* Indicate whether an iteration to
101 initialize ALE is required */
102
103/*============================================================================
104 * Public function prototypes
105 *============================================================================*/
106
107/*----------------------------------------------------------------------------*/
108/*
109 * \brief Compute of ALE volumic flux from displacement and deduced volume
110 * at time n+1.
111 *
112 * \param[in, out] domain pointer to a cs_domain_t structure
113 */
114/*----------------------------------------------------------------------------*/
115
116void
118
119/*----------------------------------------------------------------------------*/
120/*
121 * \brief In the ALE framework, update mass flux by adding mesh velocity.
122 *
123 * \param[in] m pointer to associated mesh structure
124 * \param[in] mq pointer to associated mesh quantities structure
125 * \param[in] dt time step at cells
126 * \param[in] crom density at cells
127 * \param[in] brom density at boundary faces
128 * \param[in, out] imasfl interior face mass flux
129 * \param[in, out] bmasfl boundary face mass flux
130 */
131/*----------------------------------------------------------------------------*/
132
133void
135 const cs_mesh_quantities_t *mq,
136 const cs_real_t dt[],
137 const cs_real_t crom[],
138 const cs_real_t brom[],
139 cs_real_t imasfl[],
140 cs_real_t bmasfl[]);
141
142/*----------------------------------------------------------------------------*/
143/*
144 * Compute boundary condition code for ALE
145 *----------------------------------------------------------------------------*/
146
147void
149 const cs_mesh_quantities_t *mq,
150 const bool init,
151 const cs_real_t dt[],
152 const int bc_type[],
153 cs_real_t *rcodcl1_vel);
154
155/*----------------------------------------------------------------------------*/
156/*
157 * Compute boundary condition code for ALE
158 *----------------------------------------------------------------------------*/
159
160void
162 const cs_mesh_quantities_t *mq,
163 const bool init,
164 const cs_real_t dt[],
165 const int bc_type[],
166 cs_real_t *rcodcl1_vel);
167
168/*----------------------------------------------------------------------------*/
169/*
170 * \brief Allocation of ialtyb and impale for the ALE structure.
171 */
172/*----------------------------------------------------------------------------*/
173
174void cs_ale_allocate(void);
175
176/*----------------------------------------------------------------------------*/
177/*
178 * \brief Compute cell and face centers of gravity, cell volumes
179 * and update bad cells.
180 *
181 * \param[out] min_vol Minimum cell volume
182 * \param[out] max_vol Maximum cell volume
183 * \param[out] tot_vol Total cell volume
184 */
185/*----------------------------------------------------------------------------*/
186
187void
189 cs_real_t *max_vol,
190 cs_real_t *tot_vol);
191
192/*----------------------------------------------------------------------------*/
193/*
194 * \brief Project the displacement on mesh vertices (solved on cell center).
195 *
196 * \param[in] ale_bc_type Type of boundary for ALE
197 * \param[in] meshv Mesh velocity
198 * \param[in] gradm Mesh velocity gradient
199 * (du_i/dx_j : gradv[][i][j])
200 * \param[in] claale Boundary conditions A
201 * \param[in] clbale Boundary conditions B
202 * \param[in] dt Time step
203 * \param[out] disp_proj Displacement projected on vertices
204 */
205/*----------------------------------------------------------------------------*/
206
207void
208cs_ale_project_displacement(const int ale_bc_type[],
209 const cs_real_3_t *meshv,
210 const cs_real_33_t gradm[],
211 const cs_real_3_t *claale,
212 const cs_real_33_t *clbale,
213 const cs_real_t *dt,
214 cs_real_3_t *disp_proj);
215
216/*----------------------------------------------------------------------------*/
217/*
218 * \brief Update mesh in the ALE framework.
219 *
220 * \param[in] itrale number of the current ALE iteration
221 */
222/*----------------------------------------------------------------------------*/
223
224void
225cs_ale_update_mesh(int itrale);
226
227/*----------------------------------------------------------------------------*/
228/*
229 * \brief Update ALE BCs for required for the fluid
230 *
231 * \param[out] ale_bc_type type of ALE bcs
232 * \param[out] b_fluid_vel Fluid velocity at boundary faces
233 */
234/*----------------------------------------------------------------------------*/
235
236void
237cs_ale_update_bcs(int *ale_bc_type,
238 cs_real_3_t *b_fluid_vel);
239
240/*----------------------------------------------------------------------------*/
241/*
242 * \brief Solve a Poisson equation on the mesh velocity in ALE framework.
243 *
244 * It also updates the mesh displacement
245 * so that it can be used to update mass fluxes (due to mesh displacement).
246 *
247 * \param[in] iterns Navier-Stokes iteration number
248 */
249/*----------------------------------------------------------------------------*/
250
251void
252cs_ale_solve_mesh_velocity(const int iterns,
253 const cs_real_t b_rho[],
254 const cs_real_t i_mass_flux[],
255 const cs_real_t b_mass_flux[]);
256
257/*----------------------------------------------------------------------------*/
258/*
259 * \brief Activate the mesh velocity solving with CDO
260 */
261/*----------------------------------------------------------------------------*/
262
263void
264cs_ale_activate(void);
265
266/*----------------------------------------------------------------------------*/
267/*
268 * \brief Test if mesh velocity solving with CDO is activated
269 *
270 * \return true ifmesh velocity solving with CDO is requested, false otherwise
271 */
272/*----------------------------------------------------------------------------*/
273
274bool
276
277/*----------------------------------------------------------------------------*/
278/*
279 * \brief Add "property" fields dedicated to the ALE model.
280 */
281/*----------------------------------------------------------------------------*/
282
283void
285
286/*----------------------------------------------------------------------------*/
287/*
288 * \brief Initialize fields volume dedicated to the ALE model.
289 */
290/*----------------------------------------------------------------------------*/
291
292void
294
295/*----------------------------------------------------------------------------*/
296/*
297 * \brief Setup the equations solving the mesh velocity when CDO is activated
298 *
299 * \param[in, out] domain pointer to a cs_domain_t structure
300 */
301/*----------------------------------------------------------------------------*/
302
303void
305
306/*----------------------------------------------------------------------------
307 *
308 * \brief Print the ALE options to setup.log.
309 *
310 *----------------------------------------------------------------------------*/
311
312void
313cs_ale_log_setup(void);
314
315/*----------------------------------------------------------------------------*/
316/*
317 * \brief Setup the equations solving the mesh velocity
318 *
319 * \param[in] domain pointer to a cs_domain_t structure
320 */
321/*----------------------------------------------------------------------------*/
322
323void
325
326/*----------------------------------------------------------------------------*/
327/*
328 * \brief Finalize the setup stage for the equation of the mesh velocity
329 *
330 * \param[in, out] domain pointer to a cs_domain_t structure
331 */
332/*----------------------------------------------------------------------------*/
333
334void
336
337/*----------------------------------------------------------------------------*/
338/*
339 * \brief Free the main structure related to the ALE mesh velocity solving
340 */
341/*----------------------------------------------------------------------------*/
342
343void
345
346/*----------------------------------------------------------------------------*/
347/*
348 * \brief Read ALE data from restart file.
349 *
350 * \param[in, out] r associated restart file pointer
351 */
352/*----------------------------------------------------------------------------*/
353
354void
356
357/*----------------------------------------------------------------------------*/
358/*
359 * \brief Write ALE data from restart file.
360 *
361 * \param[in, out] r associated restart file pointer
362 */
363/*----------------------------------------------------------------------------*/
364
365void
367
368/*----------------------------------------------------------------------------*/
369
371
372#endif /* __CS_ALE_H__ */
void cs_ale_solve_mesh_velocity(const int iterns, const cs_real_t b_rho[], const cs_real_t i_mass_flux[], const cs_real_t b_mass_flux[])
Solve a Poisson equation on the mesh velocity in ALE framework.
Definition: cs_ale.cpp:2287
void cs_mesh_velocity_mass_flux(const cs_mesh_t *m, const cs_mesh_quantities_t *mq, const cs_real_t dt[], const cs_real_t crom[], const cs_real_t brom[], cs_real_t imasfl[], cs_real_t bmasfl[])
Definition: cs_ale.cpp:1450
void cs_ale_log_setup(void)
Definition: cs_ale.cpp:2543
void cs_ale_activate(void)
Activate the mesh velocity solving with CDO.
Definition: cs_ale.cpp:2310
void cs_ale_finalize_setup(cs_domain_t *domain)
Finalize the setup stage for the equation of the mesh velocity.
Definition: cs_ale.cpp:2692
int cs_glob_ale_n_ini_f
Definition: cs_ale.cpp:113
void cs_ale_update_mesh(int itrale)
Update mesh in the ALE framework.
Definition: cs_ale.cpp:2192
void cs_ale_update_bcs(int *ale_bc_type, cs_real_3_t *b_fluid_vel)
Update ALE BCs for required for the fluid.
Definition: cs_ale.cpp:2266
cs_ale_data_t * cs_glob_ale_data
Definition: cs_ale.cpp:110
void cs_boundary_condition_ale_type_nep(const cs_mesh_t *m, const cs_mesh_quantities_t *mq, const bool init, const cs_real_t dt[], const int bc_type[], cs_real_t *rcodcl1_vel)
void cs_ale_add_property_fields(void)
Add "property" fields dedicated to the ALE model.
Definition: cs_ale.cpp:2367
void cs_ale_initialize_volume_fields(void)
Initialize volume fields dedicated to the ALE model.
Definition: cs_ale.cpp:2440
void cs_ale_destroy_all(void)
Free the main structure related to the ALE mesh velocity solving.
Definition: cs_ale.cpp:2715
cs_ale_type_t cs_glob_ale
Definition: cs_ale.cpp:108
void cs_ale_restart_write(cs_restart_t *r)
Write ALE data from restart file.
Definition: cs_ale.cpp:2804
void cs_ale_compute_volume_from_displacement(cs_domain_t *domain)
Compute of ALE volumic flux from displacement and deduced volume at time step n+1.
Definition: cs_ale.cpp:1225
void cs_ale_init_setup(cs_domain_t *domain)
Setup the equations solving the mesh velocity when CDO is activated.
Definition: cs_ale.cpp:2488
void cs_boundary_condition_ale_type(const cs_mesh_t *m, const cs_mesh_quantities_t *mq, const bool init, const cs_real_t dt[], const int bc_type[], cs_real_t *rcodcl1_vel)
Definition: cs_ale.cpp:1596
void cs_ale_setup_boundaries(const cs_domain_t *domain)
Setup the equations solving the mesh velocity.
Definition: cs_ale.cpp:2581
void cs_ale_allocate(void)
Allocation of ialtyb and impale for the ALE structure.
Definition: cs_ale.cpp:1907
bool cs_ale_is_activated(void)
Test if mesh velocity solving with CDO is activated.
Definition: cs_ale.cpp:2352
cs_ale_type_t
Definition: cs_ale.h:57
@ CS_ALE_LEGACY
Definition: cs_ale.h:60
@ CS_ALE_NONE
Definition: cs_ale.h:59
@ CS_ALE_CDO
Definition: cs_ale.h:61
void cs_ale_restart_read(cs_restart_t *r)
Read ALE data from restart file.
Definition: cs_ale.cpp:2740
void cs_ale_update_mesh_quantities(cs_real_t *min_vol, cs_real_t *max_vol, cs_real_t *tot_vol)
Compute cell and face centers of gravity, cell volumes and update bad cells.
Definition: cs_ale.cpp:1951
void cs_ale_project_displacement(const int ale_bc_type[], const cs_real_3_t *meshv, const cs_real_33_t gradm[], const cs_real_3_t *claale, const cs_real_33_t *clbale, const cs_real_t *dt, cs_real_3_t *disp_proj)
Project the displacement on mesh vertices (solved on cell center).
Definition: cs_ale.cpp:1984
int cs_glob_ale_need_init
Definition: cs_ale.cpp:116
#define BEGIN_C_DECLS
Definition: cs_defs.h:554
double cs_real_t
Floating-point value.
Definition: cs_defs.h:357
cs_real_t cs_real_3_t[3]
vector of 3 floating-point values
Definition: cs_defs.h:374
#define END_C_DECLS
Definition: cs_defs.h:555
cs_real_t cs_real_33_t[3][3]
3x3 matrix of floating-point values
Definition: cs_defs.h:383
cs_real_t cs_real_63_t[6][3]
Definition: cs_defs.h:391
@ dt
Definition: cs_field_pointer.h:65
struct _cs_restart_t cs_restart_t
Definition: cs_restart.h:95
Definition: cs_ale.h:67
cs_real_t ** alprof_pre
Definition: cs_ale.h:84
cs_real_t ** i_mass_flux_pre
Definition: cs_ale.h:82
cs_real_33_t ** gradu_pre
Definition: cs_ale.h:77
cs_real_t * i_mass_flux_ale
Definition: cs_ale.h:73
int * impale
Definition: cs_ale.h:69
cs_real_t * b_mass_flux_ale
Definition: cs_ale.h:74
cs_real_t ** b_mass_flux_pre
Definition: cs_ale.h:83
int * bc_type
Definition: cs_ale.h:70
cs_real_t ** alprob_pre
Definition: cs_ale.h:85
cs_real_3_t ** gradvolf_pre
Definition: cs_ale.h:79
cs_real_t ** gamma_pre
Definition: cs_ale.h:86
cs_real_63_t ** gradrij_pre
Definition: cs_ale.h:81
int ale_iteration
Definition: cs_ale.h:71
int implicit_coupling_ite
Definition: cs_ale.h:72
cs_real_3_t ** gradhtot_pre
Definition: cs_ale.h:80
cs_real_3_t * gradp_pre
Definition: cs_ale.h:78
Structure storing the main features of the computational domain and pointers to the main geometrical ...
Definition: cs_domain.h:129
Definition: cs_mesh_quantities.h:92
Definition: cs_mesh.h:85