My Project
programmer's documentation
Loading...
Searching...
No Matches
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-2019 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
53enum {
54
57 CS_ALE_CDO = 2
58
59};
60
61/*=============================================================================
62 * Global variables
63 *============================================================================*/
64
65extern int cs_glob_ale;
66
67/*============================================================================
68 * Fortran wrapper function definitions
69 *============================================================================*/
70
71/*----------------------------------------------------------------------------
72 * Get pointer to cs_glob_ale
73 *----------------------------------------------------------------------------*/
74
75void
76cs_f_ale_get_pointers(int **iale);
77
78/*============================================================================
79 * Public function prototypes
80 *============================================================================*/
81
82/*----------------------------------------------------------------------------*/
91/*----------------------------------------------------------------------------*/
92
93void
95 cs_real_t *max_vol,
96 cs_real_t *tot_vol);
97
98/*----------------------------------------------------------------------------*/
111/*----------------------------------------------------------------------------*/
112
113void
114cs_ale_project_displacement(const int ale_bc_type[],
115 const cs_real_3_t *meshv,
116 const cs_real_33_t gradm[],
117 const cs_real_3_t *claale,
118 const cs_real_33_t *clbale,
119 const cs_real_t *dt,
120 cs_real_3_t *disp_proj);
121
122/*----------------------------------------------------------------------------*/
129/*----------------------------------------------------------------------------*/
130
131void
132cs_ale_update_mesh(const int itrale,
133 const cs_real_3_t *xyzno0);
134
135/*----------------------------------------------------------------------------*/
142/*----------------------------------------------------------------------------*/
143
144void
145cs_ale_update_bcs(int *ale_bc_type,
146 cs_real_3_t *b_fluid_vel);
147
148/*----------------------------------------------------------------------------*/
159/*----------------------------------------------------------------------------*/
160
161void
162cs_ale_solve_mesh_velocity(const int iterns,
163 const int *impale,
164 const int *ale_bc_type);
165
166/*----------------------------------------------------------------------------*/
170/*----------------------------------------------------------------------------*/
171
172void
173cs_ale_activate(void);
174
175/*----------------------------------------------------------------------------*/
181/*----------------------------------------------------------------------------*/
182
183bool
185
186/*----------------------------------------------------------------------------*/
192/*----------------------------------------------------------------------------*/
193
194void
196
197/*----------------------------------------------------------------------------*/
203/*----------------------------------------------------------------------------*/
204
205void
207
208/*----------------------------------------------------------------------------*/
214/*----------------------------------------------------------------------------*/
215
216void
218
219/*----------------------------------------------------------------------------*/
223/*----------------------------------------------------------------------------*/
224
225void
227
228/*----------------------------------------------------------------------------*/
229
231
232#endif /* __CS_ALE_H__ */
void cs_ale_solve_mesh_velocity(const 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:1220
@ 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:1238
void cs_ale_finalize_setup(cs_domain_t *domain)
Finalize the setup stage for the equation of the mesh velocity.
Definition cs_ale.c:1434
void cs_ale_update_mesh(const int itrale, const cs_real_3_t *xyzno0)
Update mesh in the ALE framework.
Definition cs_ale.c:1121
int cs_glob_ale
Definition cs_ale.c:79
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:1196
void cs_f_ale_get_pointers(int **iale)
Definition cs_ale.c:90
void cs_ale_destroy_all(void)
Free the main structure related to the ALE mesh velocity solving.
Definition cs_ale.c:1460
void cs_ale_init_setup(cs_domain_t *domain)
Setup the equations related to mesh deformation.
Definition cs_ale.c:1293
void cs_ale_setup_boundaries(const cs_domain_t *domain)
Setup the equations solving the mesh velocity.
Definition cs_ale.c:1328
bool cs_ale_is_activated(void)
Test if mesh velocity solving with CDO is activated.
Definition cs_ale.c:1276
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:873
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:904
#define BEGIN_C_DECLS
Definition cs_defs.h:467
double cs_real_t
Floating-point value.
Definition cs_defs.h:302
cs_real_t cs_real_3_t[3]
vector of 3 floating-point values
Definition cs_defs.h:315
#define END_C_DECLS
Definition cs_defs.h:468
cs_real_t cs_real_33_t[3][3]
3x3 matrix of floating-point values
Definition cs_defs.h:321
@ dt
Definition cs_field_pointer.h:65
Structure storing the main features of the computational domain and pointers to the main geometrical ...
Definition cs_domain.h:87