8.0
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-2023 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 "cs_base.h"
39 #include "cs_domain.h"
40 
41 /*----------------------------------------------------------------------------*/
42 
44 
45 /*============================================================================
46  * Type definitions
47  *============================================================================*/
48 
49 /*----------------------------------------------------------------------------
50  * ALE type
51  *----------------------------------------------------------------------------*/
52 
53 enum {
54 
57  CS_ALE_CDO = 2
58 
59 };
60 
61 /*=============================================================================
62  * Global variables
63  *============================================================================*/
64 
65 extern int cs_glob_ale;
66 
67 extern int cs_glob_ale_n_ini_f;
70 extern int cs_glob_ale_need_init;
73 /*============================================================================
74  * Public function prototypes
75  *============================================================================*/
76 
77 /*----------------------------------------------------------------------------*/
86 /*----------------------------------------------------------------------------*/
87 
88 void
90  cs_real_t *max_vol,
91  cs_real_t *tot_vol);
92 
93 /*----------------------------------------------------------------------------*/
106 /*----------------------------------------------------------------------------*/
107 
108 void
109 cs_ale_project_displacement(const int ale_bc_type[],
110  const cs_real_3_t *meshv,
111  const cs_real_33_t gradm[],
112  const cs_real_3_t *claale,
113  const cs_real_33_t *clbale,
114  const cs_real_t *dt,
115  cs_real_3_t *disp_proj);
116 
117 /*----------------------------------------------------------------------------*/
123 /*----------------------------------------------------------------------------*/
124 
125 void
126 cs_ale_update_mesh(int itrale);
127 
128 /*----------------------------------------------------------------------------*/
135 /*----------------------------------------------------------------------------*/
136 
137 void
138 cs_ale_update_bcs(int *ale_bc_type,
139  cs_real_3_t *b_fluid_vel);
140 
141 /*----------------------------------------------------------------------------*/
152 /*----------------------------------------------------------------------------*/
153 
154 void
155 cs_ale_solve_mesh_velocity(int iterns,
156  const int *impale,
157  const int *ale_bc_type);
158 
159 /*----------------------------------------------------------------------------*/
163 /*----------------------------------------------------------------------------*/
164 
165 void
166 cs_ale_activate(void);
167 
168 /*----------------------------------------------------------------------------*/
174 /*----------------------------------------------------------------------------*/
175 
176 bool
177 cs_ale_is_activated(void);
178 
179 /*----------------------------------------------------------------------------*/
185 /*----------------------------------------------------------------------------*/
186 
187 void
189 
190 /*----------------------------------------------------------------------------*/
196 /*----------------------------------------------------------------------------*/
197 
198 void
199 cs_ale_setup_boundaries(const cs_domain_t *domain);
200 
201 /*----------------------------------------------------------------------------*/
207 /*----------------------------------------------------------------------------*/
208 
209 void
211 
212 /*----------------------------------------------------------------------------*/
216 /*----------------------------------------------------------------------------*/
217 
218 void
219 cs_ale_destroy_all(void);
220 
221 /*----------------------------------------------------------------------------*/
222 
224 
225 #endif /* __CS_ALE_H__ */
@ CS_ALE_LEGACY
Definition: cs_ale.h:56
@ CS_ALE_NONE
Definition: cs_ale.h:55
@ CS_ALE_CDO
Definition: cs_ale.h:57
void cs_ale_activate(void)
Activate the mesh velocity solving with CDO.
Definition: cs_ale.c:1518
void cs_ale_finalize_setup(cs_domain_t *domain)
Finalize the setup stage for the equation of the mesh velocity.
Definition: cs_ale.c:1747
int cs_glob_ale_n_ini_f
int cs_glob_ale
Definition: cs_ale.c:83
void cs_ale_update_mesh(int itrale)
Update mesh in the ALE framework.
Definition: cs_ale.c:1397
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.c:1477
void cs_ale_destroy_all(void)
Free the main structure related to the ALE mesh velocity solving.
Definition: cs_ale.c:1770
void cs_ale_init_setup(cs_domain_t *domain)
Setup the equations related to mesh deformation.
Definition: cs_ale.c:1577
void cs_ale_solve_mesh_velocity(int iterns, const int *impale, const int *ale_bc_type)
Solve a Poisson equation on the mesh velocity in ALE framework.
Definition: cs_ale.c:1500
void cs_ale_setup_boundaries(const cs_domain_t *domain)
Setup the equations solving the mesh velocity.
Definition: cs_ale.c:1641
bool cs_ale_is_activated(void)
Test if mesh velocity solving with CDO is activated.
Definition: cs_ale.c:1560
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.c:1140
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.c:1173
int cs_glob_ale_need_init
#define BEGIN_C_DECLS
Definition: cs_defs.h:509
double cs_real_t
Floating-point value.
Definition: cs_defs.h:319
cs_real_t cs_real_3_t[3]
vector of 3 floating-point values
Definition: cs_defs.h:332
#define END_C_DECLS
Definition: cs_defs.h:510
cs_real_t cs_real_33_t[3][3]
3x3 matrix of floating-point values
Definition: cs_defs.h:341
@ dt
Definition: cs_field_pointer.h:65
integer, dimension(:), allocatable impale
indicator of imposed displacement
Definition: albase.f90:61
Structure storing the main features of the computational domain and pointers to the main geometrical ...
Definition: cs_domain.h:138