8.3
general documentation
cs_prototypes.h
Go to the documentation of this file.
1#ifndef __CS_PROTOTYPES_H__
2#define __CS_PROTOTYPES_H__
3
4/*============================================================================
5 * Prototypes for Fortran functions and subroutines callable from C
6 *============================================================================*/
7
8/*
9 This file is part of code_saturne, a general-purpose CFD tool.
10
11 Copyright (C) 1998-2024 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_base.h"
35#include "cs_domain.h"
36#include "cs_field.h"
37#include "cs_mesh.h"
38#include "cs_mesh_quantities.h"
39#include "cs_mesh_bad_cells.h"
40#include "cs_probe.h"
41#include "cs_time_control.h"
42#include "cs_volume_zone.h"
43
44/*----------------------------------------------------------------------------*/
45
47
48/*============================================================================
49 * Macro definitions
50 *============================================================================*/
51
52/*=============================================================================
53 * Fortran function/subroutine prototypes
54 *============================================================================*/
55
56/*----------------------------------------------------------------------------
57 * Main Fortran subroutine
58 *----------------------------------------------------------------------------*/
59
60/*----------------------------------------------------------------------------
61 * Initialize Fortran base common block values
62 *----------------------------------------------------------------------------*/
63
64extern void CS_PROCF (csinit, CSINIT)
65(
66 const int *irgpar, /* <-- MPI Rank in parallel, -1 otherwise */
67 const int *nrgpar /* <-- Number of MPI processes, or 1 */
68);
69
70/*----------------------------------------------------------------------------
71 * Find the nearest cell's center from a node
72 *----------------------------------------------------------------------------*/
73
74extern void CS_PROCF (findpt, FINDPT)
75(
76 const cs_lnum_t *ncelet, /* <-- number of extended (real + ghost) cells */
77 const cs_lnum_t *ncel, /* <-- number of cells */
78 const cs_real_t *xyzcen, /* <-- cell centers */
79 const cs_real_t *xx, /* <-- node coordinate X */
80 const cs_real_t *yy, /* <-- node coordinate Y */
81 const cs_real_t *zz, /* <-- node coordinate Z */
82 cs_lnum_t *node, /* --> node we are looking for, zero if error */
83 int *ndrang /* --> rank of associated process */
84);
85
86/*----------------------------------------------------------------------------
87 * Add field indexes associated with a new non-user solved variable,
88 * with default options
89 *
90 * parameters:
91 * f_id <-- field id
92 *
93 * returns:
94 * variable number for defined field
95 *----------------------------------------------------------------------------*/
96
97int
99
100/*----------------------------------------------------------------------------
101 * Add field indexes associated with a new non-user solved variable,
102 * with default options
103 *
104 * parameters:
105 * f_id <-- field id
106 *
107 * returns:
108 * scalar number for defined field
109 *----------------------------------------------------------------------------*/
110
111int
113
114/*----------------------------------------------------------------------------
115 * Add field indexes associated with a new solved thermal variable,
116 * with default options
117 *
118 * parameters:
119 * f_id <-- field id
120 *----------------------------------------------------------------------------*/
121
122void
124
125/*----------------------------------------------------------------------------
126 * Computes the explicit chemical source term for atmospheric chemistry in
127 * case of a semi-coupled resolution
128 *----------------------------------------------------------------------------*/
129
130void
132 cs_real_t st_exp[],
133 cs_real_t st_imp[]);
134
135/*----------------------------------------------------------------------------*/
140/*----------------------------------------------------------------------------*/
141
142void
144
145/*----------------------------------------------------------------------------*/
150/*----------------------------------------------------------------------------*/
151
152void
154
155/*----------------------------------------------------------------------------*/
164/*----------------------------------------------------------------------------*/
165
168 cs_real_t h);
169
170/*----------------------------------------------------------------------------*/
179/*----------------------------------------------------------------------------*/
180
183 cs_real_t t);
184
185/*----------------------------------------------------------------------------
186 * Return Lagrangian model status.
187 *
188 * parameters:
189 * model_flag --> 0 without Lagrangian, 1 or 2 with Lagrangian
190 * restart_flag --> 1 for Lagrangian restart, 0 otherwise
191 * frozen_flag --> 1 for frozen Eulerian flow, 0 otherwise
192 *----------------------------------------------------------------------------*/
193
194void
195cs_lagr_status(int *model_flag,
196 int *restart_flag,
197 int *frozen_flag);
198
199/*----------------------------------------------------------------------------*/
203/*----------------------------------------------------------------------------*/
204
205void
207
208/*----------------------------------------------------------------------------*/
227/*----------------------------------------------------------------------------*/
228
229void
231
232/*----------------------------------------------------------------------------
233 * Compute source terms for specific physical model scalars.
234 *----------------------------------------------------------------------------*/
235
236void
238 cs_real_t st_imp[],
239 cs_real_t st_exp[]);
240
241/*----------------------------------------------------------------------------*/
262/*----------------------------------------------------------------------------*/
263
264void
266 cs_lnum_t icetsm[],
267 cs_real_t dt[]);
268
269/*============================================================================
270 * User function prototypes
271 *============================================================================*/
272
273/*----------------------------------------------------------------------------
274 * Data Entry of the 1D wall thermal module.
275 *----------------------------------------------------------------------------*/
276
277void
278cs_user_1d_wall_thermal(int iappel);
279
280/*----------------------------------------------------------------------------
281 * Fill in vertical profiles of atmospheric properties prior to solving
282 * 1D radiative transfers.
283 *----------------------------------------------------------------------------*/
284
285void
287 cs_real_t temray[],
288 cs_real_t romray[],
289 cs_real_t qvray[],
290 cs_real_t qlray[],
291 cs_real_t ncray[],
292 cs_real_t aeroso[]);
293
294/*----------------------------------------------------------------------------
295 * Data Entry of the wall condensation module
296 *----------------------------------------------------------------------------*/
297
298void
299cs_user_wall_condensation(int iappel);
300
301/*----------------------------------------------------------------------------
302 * Setup boundary conditions to be applied.
303 *----------------------------------------------------------------------------*/
304
305void
307
308/*----------------------------------------------------------------------------
309 * This function is called at each time step for boundary conditions.
310 *----------------------------------------------------------------------------*/
311
312void
314 int bc_type[]);
315
316/*----------------------------------------------------------------------------
317 * Boundary conditions for (Arbitrary Lagrangian Eulerian).
318 *----------------------------------------------------------------------------*/
319
320void
322 int bc_type[],
323 int ale_bc_type[],
324 int impale[]);
325
326/*----------------------------------------------------------------------------*/
338/*----------------------------------------------------------------------------*/
339
340void
342
343/*----------------------------------------------------------------------------*/
352/*----------------------------------------------------------------------------*/
353
354void
356
357/*----------------------------------------------------------------------------*/
366/*----------------------------------------------------------------------------*/
367
368void
370
371/*----------------------------------------------------------------------------*/
394/*----------------------------------------------------------------------------*/
395
396void
397cs_user_fsi_structure_define(int is_restart,
398 int n_structs,
399 int *plot,
400 cs_time_control_t *plot_time_control,
401 cs_real_t *aexxst,
402 cs_real_t *bexxst,
403 cs_real_t *cfopre,
404 cs_real_t xstr0[][3],
405 cs_real_t vstr0[][3],
406 cs_real_t xstreq[][3]);
407
408/*----------------------------------------------------------------------------*/
423/*----------------------------------------------------------------------------*/
424
425void
426cs_user_fsi_structure_values(int n_structs,
427 const cs_time_step_t *ts,
428 const cs_real_t xstreq[][3],
429 const cs_real_t xstr[][3],
430 const cs_real_t vstr[][3],
431 cs_real_t xmstru[][3][3],
432 cs_real_t xcstru[][3][3],
433 cs_real_t xkstru[][3][3],
434 cs_real_t forstr[][3],
435 cs_real_t dtstr[]);
436
437/*----------------------------------------------------------------------------*/
450/*----------------------------------------------------------------------------*/
451
452void
454 int structure_num[]);
455
456/*----------------------------------------------------------------------------*/
466/*----------------------------------------------------------------------------*/
467
468void
470 cs_real_t cku[][6]);
471
472/*----------------------------------------------------------------------------*/
478/*----------------------------------------------------------------------------*/
479
480void
482
483/*----------------------------------------------------------------------------*/
494/*----------------------------------------------------------------------------*/
495
496void
498
499/*----------------------------------------------------------------------------*/
505/*----------------------------------------------------------------------------*/
506
507void
509
510/*----------------------------------------------------------------------------*/
519/*----------------------------------------------------------------------------*/
520
521void
523
524/*----------------------------------------------------------------------------*/
533/*----------------------------------------------------------------------------*/
534
535void
537
538/*----------------------------------------------------------------------------*/
544/*----------------------------------------------------------------------------*/
545
546void
548
549/*----------------------------------------------------------------------------*/
567/*----------------------------------------------------------------------------*/
568
569void
571 const cs_zone_t *z,
572 bool z_local,
573 const cs_real_t h[],
574 cs_real_t t[]);
575
576/*----------------------------------------------------------------------------*/
594/*----------------------------------------------------------------------------*/
595
596void
598 const cs_zone_t *z,
599 bool z_local,
600 const cs_real_t t[],
601 cs_real_t h[]);
602
603/*----------------------------------------------------------------------------*/
617/*----------------------------------------------------------------------------*/
618
619void
621
622/*----------------------------------------------------------------------------*/
628/*----------------------------------------------------------------------------*/
629
630void
632
633/*----------------------------------------------------------------------------*/
643/*----------------------------------------------------------------------------*/
644
645void
647 int f_id,
648 cs_real_t *st_exp,
649 cs_real_t *st_imp);
650
651/*----------------------------------------------------------------------------*/
661/*----------------------------------------------------------------------------*/
662
663void
665
666/*----------------------------------------------------------------------------
667 * Define mesh joinings.
668 *----------------------------------------------------------------------------*/
669
670void
671cs_user_join(void);
672
673/*----------------------------------------------------------------------------
674 * Define linear solver options.
675 *
676 * This function is called at the setup stage, once user and most model-based
677 * fields are defined.
678 *----------------------------------------------------------------------------*/
679
680void
682
683/*----------------------------------------------------------------------------*/
689/*----------------------------------------------------------------------------*/
690
691void
693
694/*----------------------------------------------------------------------------
695 * Tag bad cells within the mesh based on geometric criteria.
696 *----------------------------------------------------------------------------*/
697
698void
700 cs_mesh_quantities_t *mesh_quantities);
701
702/*----------------------------------------------------------------------------
703 * Define restart mode for mesh input and preprocessing.
704 *----------------------------------------------------------------------------*/
705
706void
708
709/*----------------------------------------------------------------------------
710 * Define mesh files to read and optional associated transformations.
711 *----------------------------------------------------------------------------*/
712
713void
715
716/*----------------------------------------------------------------------------
717 * Modifiy geometry and mesh.
718 *----------------------------------------------------------------------------*/
719
720void
722
723/*----------------------------------------------------------------------------
724 * Insert boundary wall into a mesh.
725 *----------------------------------------------------------------------------*/
726
727void
729
730/*----------------------------------------------------------------------------
731 * Mesh smoothing.
732 *
733 * parameters:
734 * mesh <-> pointer to mesh structure to smoothe
735 *----------------------------------------------------------------------------*/
736
737void
739
740/*----------------------------------------------------------------------------
741 * Enable or disable mesh saving.
742 *
743 * By default, mesh is saved when modified.
744 *
745 * parameters:
746 * mesh <-> pointer to mesh structure
747 *----------------------------------------------------------------------------*/
748
749void
751
752/*----------------------------------------------------------------------------
753 * Set options for cutting of warped faces
754 *
755 * parameters:
756 * mesh <-> pointer to mesh structure to smoothe
757 *----------------------------------------------------------------------------*/
758
759void
761
762/*----------------------------------------------------------------------------*/
770/*----------------------------------------------------------------------------*/
771
772void
774 cs_mesh_quantities_t *mesh_quantities);
775
776/*----------------------------------------------------------------------------*/
780/*----------------------------------------------------------------------------*/
781
782void
784
785/*----------------------------------------------------------------------------
786 * Select physical model options, including user fields.
787 *
788 * This function is called at the earliest stages of the data setup.
789 *----------------------------------------------------------------------------*/
790
791void
792cs_user_model(void);
793
794/*----------------------------------------------------------------------------
795 * Define advanced mesh numbering options.
796 *----------------------------------------------------------------------------*/
797
798void
800
801/*----------------------------------------------------------------------------
802 * Define parallel IO settings.
803 *----------------------------------------------------------------------------*/
804
805void
807
808/*----------------------------------------------------------------------------
809 * Define advanced partitioning options.
810 *----------------------------------------------------------------------------*/
811
812void
814
815/*----------------------------------------------------------------------------
816 * User function for input of radiative transfer module options.
817 *
818 * Deprecated Use cs_user_model instead.
819 *----------------------------------------------------------------------------*/
820
821void
823
824/*----------------------------------------------------------------------------
825 * Define sparse matrix tuning options.
826 *----------------------------------------------------------------------------*/
827
828void
830
831/*----------------------------------------------------------------------------
832 * Define or modify general numerical and physical user parameters.
833 *
834 * At the calling point of this function, most model-related most variables
835 * and other fields have been defined, so specific settings related to those
836 * fields may be set here.
837 *----------------------------------------------------------------------------*/
838
839void
841
842/*-----------------------------------------------------------------------------
843 * User subroutine for input of radiative transfer boundary conditions
844 *----------------------------------------------------------------------------*/
845
846void
848 const int bc_type[],
849 int isothp[],
850 cs_real_t *tmin,
851 cs_real_t *tmax,
852 cs_real_t *tx,
853 const cs_real_t dt[],
854 const cs_real_t thwall[],
855 const cs_real_t qincid[],
856 cs_real_t hfcnvp[],
857 cs_real_t flcnvp[],
858 cs_real_t xlamp[],
859 cs_real_t epap[],
860 cs_real_t epsp[],
861 cs_real_t textp[]);
862
863/*----------------------------------------------------------------------------
864 * Define periodic faces.
865 *----------------------------------------------------------------------------*/
866
867void
869
870/*----------------------------------------------------------------------------
871 * Define post-processing writers.
872 *
873 * The default output format and frequency may be configured, and additional
874 * post-processing writers allowing outputs in different formats or with
875 * different format options and output frequency than the main writer may
876 * be defined.
877 *----------------------------------------------------------------------------*/
878
879void
881
882/*-----------------------------------------------------------------------------
883 * Define monitoring probes and profiles. A profile is seen as a set of probes.
884 *----------------------------------------------------------------------------*/
885
886void
888
889/*----------------------------------------------------------------------------
890 * Define post-processing meshes.
891 *
892 * The main post-processing meshes may be configured, and additional
893 * post-processing meshes may be defined as a subset of the main mesh's
894 * cells or faces (both interior and boundary).
895 *----------------------------------------------------------------------------*/
896
897void
899
900/*----------------------------------------------------------------------------
901 * User function for output of values on a post-processing mesh.
902 *----------------------------------------------------------------------------*/
903
904void
905cs_user_postprocess_values(const char *mesh_name,
906 int mesh_id,
907 int cat_id,
908 cs_probe_set_t *probes,
909 cs_lnum_t n_cells,
910 cs_lnum_t n_i_faces,
911 cs_lnum_t n_b_faces,
912 cs_lnum_t n_vertices,
913 const cs_lnum_t cell_list[],
914 const cs_lnum_t i_face_list[],
915 const cs_lnum_t b_face_list[],
916 const cs_lnum_t vertex_list[],
917 const cs_time_step_t *ts);
918
919/*----------------------------------------------------------------------------
920 * Override default frequency or calculation end based output.
921 *
922 * This allows fine-grained control of activation or deactivation,
923 *
924 * parameters:
925 * nt_max_abs <-- maximum time step number
926 * nt_cur_abs <-- current time step number
927 * t_cur_abs <-- absolute time at the current time step
928 *----------------------------------------------------------------------------*/
929
930void
931cs_user_postprocess_activate(int nt_max_abs,
932 int nt_cur_abs,
933 double t_cur_abs);
934
935/*----------------------------------------------------------------------------
936 * Absorption coefficient for radiative module
937 *----------------------------------------------------------------------------*/
938
939void
940cs_user_rad_transfer_absorption(const int bc_type[],
941 cs_real_t ck[]);
942
943/*----------------------------------------------------------------------------
944 * Compute the net radiation flux
945 *----------------------------------------------------------------------------*/
946
947void
949 const cs_real_t twall[],
950 const cs_real_t qincid[],
951 const cs_real_t xlam[],
952 const cs_real_t epa[],
953 const cs_real_t eps[],
954 const cs_real_t ck[],
955 cs_real_t net_flux[]);
956
957/*----------------------------------------------------------------------------
958 * Set user solver.
959 *----------------------------------------------------------------------------*/
960
961int
963
964/*----------------------------------------------------------------------------
965 * Main call to user solver.
966 *----------------------------------------------------------------------------*/
967
968void
970 const cs_mesh_quantities_t *mesh_quantities);
971
972/*----------------------------------------------------------------------------
973 * Define couplings with other instances of code_saturne.
974 *----------------------------------------------------------------------------*/
975
976void
978
979/*----------------------------------------------------------------------------
980 * Define couplings with SYRTHES code.
981 *----------------------------------------------------------------------------*/
982
983void
985
986/*----------------------------------------------------------------------------*/
996/*----------------------------------------------------------------------------*/
997
998void
1000 const char *syrthes_name,
1001 cs_lnum_t n_elts,
1002 const cs_lnum_t elt_ids[],
1003 cs_real_t h_vol[]);
1004
1005/*----------------------------------------------------------------------------
1006 * Define couplings with CATHARE code.
1007 *----------------------------------------------------------------------------*/
1008
1009void
1011
1012/*----------------------------------------------------------------------------
1013 * Define time moments.
1014 *----------------------------------------------------------------------------*/
1015
1016void
1018
1019
1020/*----------------------------------------------------------------------------
1021 * Define time tables.
1022 *----------------------------------------------------------------------------*/
1023
1024void
1025cs_user_time_table(void);
1026
1027/*----------------------------------------------------------------------------
1028 * Define rotor/stator model.
1029 *----------------------------------------------------------------------------*/
1030
1031void
1033
1034/*----------------------------------------------------------------------------
1035 * Define rotor axes, associated cells, and rotor/stator faces.
1036 *----------------------------------------------------------------------------*/
1037
1038void
1040
1041/*----------------------------------------------------------------------------
1042 * Define rotation velocity of rotor.
1043 *----------------------------------------------------------------------------*/
1044
1045void
1047
1048/*----------------------------------------------------------------------------*/
1052/*----------------------------------------------------------------------------*/
1053
1054void
1055cs_user_zones(void);
1056
1057/*----------------------------------------------------------------------------*/
1061/*----------------------------------------------------------------------------*/
1062
1063void
1065 const cs_mesh_quantities_t *mesh_quantities,
1066 cs_real_t *dt);
1067
1068/*----------------------------------------------------------------------------
1069 * Computation of the relaxation time-scale to equilibrium in the frame of
1070 * the homogeneous two-phase model.
1071 *----------------------------------------------------------------------------*/
1072
1073void
1075 const cs_real_t *alpha_eq,
1076 const cs_real_t *y_eq,
1077 const cs_real_t *z_eq,
1078 const cs_real_t *ei,
1079 const cs_real_t *v,
1080 cs_real_t *relax_tau);
1081
1082/*----------------------------------------------------------------------------*/
1086/*----------------------------------------------------------------------------*/
1087
1088void
1090
1091/*----------------------------------------------------------------------------*/
1095/*----------------------------------------------------------------------------*/
1096
1097void
1099
1100/*----------------------------------------------------------------------------*/
1104/*----------------------------------------------------------------------------*/
1105
1106void
1108
1109/*----------------------------------------------------------------------------*/
1110
1112
1113#endif /* __CS_PROTOTYPES_H__ */
#define BEGIN_C_DECLS
Definition: cs_defs.h:542
double cs_real_t
Floating-point value.
Definition: cs_defs.h:342
#define CS_PROCF(x, y)
Definition: cs_defs.h:576
cs_real_t cs_real_3_t[3]
vector of 3 floating-point values
Definition: cs_defs.h:359
#define END_C_DECLS
Definition: cs_defs.h:543
int cs_lnum_t
local mesh entity id
Definition: cs_defs.h:335
@ t
Definition: cs_field_pointer.h:94
@ xlam
Definition: cs_field_pointer.h:185
@ eps
Definition: cs_field_pointer.h:73
@ h
Definition: cs_field_pointer.h:93
@ epa
Definition: cs_field_pointer.h:186
@ dt
Definition: cs_field_pointer.h:65
struct _cs_probe_set_t cs_probe_set_t
Definition: cs_probe.h:53
void cs_user_postprocess_meshes(void)
Define post-processing meshes.
Definition: cs_user_postprocess.cpp:86
void cs_user_rad_transfer_absorption(const int bc_type[], cs_real_t ck[])
Absorption coefficient for radiative module.
Definition: cs_user_radiative_transfer.cpp:105
void cs_user_mesh_bad_cells_tag(cs_mesh_t *mesh, cs_mesh_quantities_t *mesh_quantities)
Tag bad cells within the mesh based on user-defined geometric criteria.
Definition: cs_user_mesh.cpp:219
void cs_user_mesh_boundary(cs_mesh_t *mesh)
Insert boundaries into a mesh.
Definition: cs_user_mesh.cpp:156
void cs_user_zones(void)
Define volume and surface zones.
Definition: cs_user_zones.cpp:65
void cs_user_model(void)
Select physical model options, including user fields.
Definition: cs_user_parameters.cpp:84
void cs_add_model_thermal_field_indexes(int f_id)
void cs_physical_model_source_terms(int iscal, cs_real_t st_imp[], cs_real_t st_exp[])
void cs_user_join(void)
Define mesh joinings.
Definition: cs_user_mesh.cpp:115
void cs_user_extra_operations_finalize(cs_domain_t *domain)
This function is called at the end of the calculation.
Definition: cs_user_extra_operations.cpp:122
void cs_lagr_status(int *model_flag, int *restart_flag, int *frozen_flag)
void cs_summon_cressman(cs_real_t the_time)
Prepare for the cressman interpolation of the variables.
void cs_user_internal_coupling(void)
Define internal coupling options.
Definition: cs_user_parameters.cpp:162
void cs_user_extra_operations(cs_domain_t *domain)
This function is called at the end of each time step.
Definition: cs_user_extra_operations.cpp:104
void cs_user_syrthes_coupling_volume_h(int coupling_id, const char *syrthes_name, cs_lnum_t n_elts, const cs_lnum_t elt_ids[], cs_real_t h_vol[])
Compute a volume exchange coefficient for SYRTHES couplings.
Definition: cs_user_coupling.cpp:114
void cs_user_syrthes_coupling(void)
Define couplings with SYRTHES code.
Definition: cs_user_coupling.cpp:95
void cs_user_parameters(cs_domain_t *domain)
Define or modify general numerical and physical user parameters.
Definition: cs_user_parameters.cpp:107
void cs_user_boundary_conditions_ale(cs_domain_t *domain, int bc_type[], int ale_bc_type[], int impale[])
User definition of boundary conditions for ALE.
Definition: cs_user_boundary_conditions.cpp:177
void cs_user_physical_properties_td_pressure(cs_real_t *td_p)
User function to define a custom law for the thermodynamic pressure.
Definition: cs_user_physical_properties.cpp:180
void cs_user_physical_properties_turb_viscosity(cs_domain_t *domain)
User modification of the turbulence viscosity.
Definition: cs_user_physical_properties.cpp:157
void cs_user_boundary_conditions_setup(cs_domain_t *domain)
Set boundary conditions to be applied.
Definition: cs_user_boundary_conditions.cpp:77
void cs_user_turbomachinery_rotor(void)
Define rotor axes, associated cells, and rotor/stator faces.
Definition: cs_user_turbomachinery.cpp:90
void cs_field_map_and_init_bcs(void)
Initialize all field BC coefficients.
Definition: cs_field_default.cpp:518
int cs_add_variable_field_indexes(int f_id)
void cs_user_turbomachinery_set_rotation_velocity(void)
Define rotation velocity of rotor.
Definition: cs_user_turbomachinery.cpp:103
void cs_user_postprocess_activate(int nt_max_abs, int nt_cur_abs, double t_cur_abs)
Definition: cs_user_postprocess.cpp:177
void cs_user_linear_solvers(void)
Define linear solver options.
Definition: cs_user_parameters.cpp:130
void cs_user_initial_conditions(cs_domain_t *domain)
This function is called at startup and allows overloading the GUI defined initialization functions.
Definition: cs_user_initialization.cpp:86
void cs_user_periodicity(void)
Define periodic faces.
Definition: cs_user_mesh.cpp:128
void findpt(const cs_lnum_t *ncelet, const cs_lnum_t *ncel, const cs_real_t *xyzcen, const cs_real_t *xx, const cs_real_t *yy, const cs_real_t *zz, cs_lnum_t *node, int *ndrang)
void cs_user_internal_coupling_from_disjoint_meshes(cs_mesh_t *mesh)
Define volumesi from separated meshes as internal coupling zones.
Definition: cs_internal_coupling.cpp:2915
void cs_user_scaling_elec(const cs_mesh_t *mesh, const cs_mesh_quantities_t *mesh_quantities, cs_real_t *dt)
Define scaling parameter for electric model.
Definition: cs_user_electric_scaling.cpp:75
void cs_user_mesh_cartesian_define(void)
Define a cartesian mesh.
Definition: cs_user_mesh.cpp:253
cs_real_t cs_gas_combustion_h_to_t(const cs_real_t x_sp[], cs_real_t h)
Convert an enthalpy to temperature value for gas combustion.
void cs_user_mesh_modify(cs_mesh_t *mesh)
Modify geometry and mesh.
Definition: cs_user_mesh.cpp:171
void cs_mass_flux_prediction(cs_real_t dt[])
Update the convective mass flux before the Navier Stokes equations (prediction and correction steps).
void cs_user_fsi_structure_num(cs_domain_t *domain, int structure_num[])
Define structure numbers for faces associated with internal or external (code_aster) structures.
Definition: cs_user_fluid_structure_interaction.cpp:174
void cs_user_radiative_transfer_bcs(cs_domain_t *domain, const int bc_type[], int isothp[], cs_real_t *tmin, cs_real_t *tmax, cs_real_t *tx, const cs_real_t dt[], const cs_real_t thwall[], const cs_real_t qincid[], cs_real_t hfcnvp[], cs_real_t flcnvp[], cs_real_t xlamp[], cs_real_t epap[], cs_real_t epsp[], cs_real_t textp[])
User definition of radiative transfer boundary conditions.
Definition: cs_user_radiative_transfer_bcs.cpp:124
void cs_user_postprocess_probes(void)
Define monitoring probes and profiles.
Definition: cs_user_postprocess.cpp:101
void cs_user_numbering(void)
Define advanced mesh numbering options.
Definition: cs_user_performance_tuning.cpp:74
void cs_user_fsi_structure_values(int n_structs, const cs_time_step_t *ts, const cs_real_t xstreq[][3], const cs_real_t xstr[][3], const cs_real_t vstr[][3], cs_real_t xmstru[][3][3], cs_real_t xcstru[][3][3], cs_real_t xkstru[][3][3], cs_real_t forstr[][3], cs_real_t dtstr[])
Time-based settings for internal mobile structures.
Definition: cs_user_fluid_structure_interaction.cpp:134
void cs_user_time_moments(void)
Define time moments.
Definition: cs_user_parameters.cpp:147
void cs_user_physical_properties_t_to_h(cs_domain_t *domain, const cs_zone_t *z, bool z_local, const cs_real_t t[], cs_real_t h[])
User definition of temperature to enthalpy conversion.
Definition: cs_user_physical_properties.cpp:134
void cs_atmo_chem_source_terms(int iscal, cs_real_t st_exp[], cs_real_t st_imp[])
void cs_user_wall_condensation(int iappel)
Source terms associated at the boundary faces and the neighboring cells with surface condensation.
Definition: cs_user_wall_condensation.cpp:167
void cs_user_radiative_transfer_parameters(void)
User function for input of radiative transfer module options.
Definition: cs_user_radiative_transfer.cpp:79
int cs_user_solver_set(void)
Set user solver.
Definition: cs_user_solver.cpp:75
void cs_at_source_term_for_inlet(cs_real_3_t st_exp[])
Additional right-hand side source terms for momentum equation in case of free inlet.
void cs_user_mesh_input(void)
Define mesh files to read and optional associated transformations.
Definition: cs_user_mesh.cpp:102
void cs_user_time_table(void)
Define time tables using C API. This function is called at the begin of the simulation only.
Definition: cs_user_time_table.cpp:77
int cs_add_model_field_indexes(int f_id)
void cs_user_1d_wall_thermal(int iappel)
Definition: cs_user_1d_wall_thermal.cpp:74
void cs_user_mesh_restart_mode(void)
Force preprocessing behavior in case of restart.
Definition: cs_user_mesh.cpp:89
void cs_user_mesh_save(cs_mesh_t *mesh)
Enable or disable mesh saving.
Definition: cs_user_mesh.cpp:203
void cs_user_porosity(cs_domain_t *domain)
Compute the porosity (volume factor when the porosity model is activated. (cs_glob_porous_model > 0)...
Definition: cs_user_porosity.cpp:80
void cs_user_source_terms(cs_domain_t *domain, int f_id, cs_real_t *st_exp, cs_real_t *st_imp)
Additional user-defined source terms for variable equations (momentum, scalars, turbulence....
Definition: cs_user_source_terms.cpp:155
void cs_user_atmo_1d_rad_prf(cs_real_t preray[], cs_real_t temray[], cs_real_t romray[], cs_real_t qvray[], cs_real_t qlray[], cs_real_t ncray[], cs_real_t aeroso[])
Fill in vertical profiles of atmospheric properties prior to solving 1D radiative transfers.
Definition: cs_user_atmo.cpp:87
void cs_user_postprocess_values(const char *mesh_name, int mesh_id, int cat_id, cs_probe_set_t *probes, cs_lnum_t n_cells, cs_lnum_t n_i_faces, cs_lnum_t n_b_faces, cs_lnum_t n_vertices, const cs_lnum_t cell_list[], const cs_lnum_t i_face_list[], const cs_lnum_t b_face_list[], const cs_lnum_t vertex_list[], const cs_time_step_t *ts)
User function for output of values on a post-processing mesh.
Definition: cs_user_postprocess.cpp:134
void cs_prediction_mass_flux(cs_lnum_t ncesmp, cs_lnum_t icetsm[], cs_real_t dt[])
Update the convective mass flux before the Navier Stokes equations (prediction and correction steps).
void cs_user_matrix_tuning(void)
Define sparse matrix tuning options.
Definition: cs_user_performance_tuning.cpp:113
void cs_user_paramedmem_define_fields(void)
Define fields to couple with ParaMEDMEM.
Definition: cs_user_paramedmem_coupling.cpp:108
void cs_user_internal_coupling_add_volumes(cs_mesh_t *mesh)
Define volumes as internal coupling zones.
Definition: cs_internal_coupling.cpp:2896
void cs_user_mesh_modify_partial(cs_mesh_t *mesh, cs_mesh_quantities_t *mesh_quantities)
Apply partial modifications to the mesh after the preprocessing and initial postprocessing mesh build...
Definition: cs_user_mesh.cpp:238
void cs_user_physical_properties(cs_domain_t *domain)
This function is called each time step to define physical properties.
Definition: cs_user_physical_properties.cpp:75
void cs_user_head_losses(const cs_zone_t *zone, cs_real_t cku[][6])
Compute GUI-defined head losses for a given volume zone.
Definition: cs_user_head_losses.cpp:86
void cs_user_fsi_structure_define(int is_restart, int n_structs, int *plot, cs_time_control_t *plot_time_control, cs_real_t *aexxst, cs_real_t *bexxst, cs_real_t *cfopre, cs_real_t xstr0[][3], cs_real_t vstr0[][3], cs_real_t xstreq[][3])
Definition of internal mobile structures and corresponding initial conditions (initial displacement a...
Definition: cs_user_fluid_structure_interaction.cpp:92
void cs_user_solver(const cs_mesh_t *mesh, const cs_mesh_quantities_t *mesh_quantities)
Main call to user solver.
Definition: cs_user_solver.cpp:92
void cs_user_postprocess_writers(void)
Define post-processing writers.
Definition: cs_user_postprocess.cpp:69
void cs_user_rad_transfer_net_flux(const int itypfb[], const cs_real_t twall[], const cs_real_t qincid[], const cs_real_t xlam[], const cs_real_t epa[], const cs_real_t eps[], const cs_real_t ck[], cs_real_t net_flux[])
Compute the net radiation flux.
Definition: cs_user_radiative_transfer.cpp:135
void cs_user_boundary_conditions(cs_domain_t *domain, int bc_type[])
User definition of boundary conditions.
Definition: cs_user_boundary_conditions.cpp:123
void cs_user_mesh_smoothe(cs_mesh_t *mesh)
Mesh smoothing.
Definition: cs_user_mesh.cpp:186
void cs_user_physical_properties_h_to_t(cs_domain_t *domain, const cs_zone_t *z, bool z_local, const cs_real_t h[], cs_real_t t[])
User definition of enthalpy to temperature conversion.
Definition: cs_user_physical_properties.cpp:102
void cs_user_initialization(cs_domain_t *domain)
This function is called one time step to initialize problem.
Definition: cs_user_initialization.cpp:107
void cs_user_extra_operations_initialize(cs_domain_t *domain)
Initialize variables.
Definition: cs_user_extra_operations.cpp:86
void cs_user_cathare_coupling(void)
Define couplings with SYRTHES code.
Definition: cs_user_coupling.cpp:138
cs_real_t cs_gas_combustion_t_to_h(const cs_real_t x_sp[], cs_real_t t)
Convert a temperature to enthalpy value for gas combustion.
void csinit(const int *irgpar, const int *nrgpar)
void cs_user_paramedmem_define_couplings(void)
Define ParaMEDMEM coupling(s)
Definition: cs_user_paramedmem_coupling.cpp:82
void cs_user_hgn_thermo_relax_time(const cs_mesh_t *mesh, const cs_real_t *alpha_eq, const cs_real_t *y_eq, const cs_real_t *z_eq, const cs_real_t *ei, const cs_real_t *v, cs_real_t *relax_tau)
Computation of the relaxation time-scale.
Definition: cs_user_hgn.cpp:82
void cs_user_mesh_warping(void)
Set options for cutting of warped faces.
Definition: cs_user_mesh.cpp:141
void cs_user_partition(void)
Define advanced partitioning options.
Definition: cs_user_performance_tuning.cpp:87
void cs_user_parallel_io(void)
Define parallel IO settings.
Definition: cs_user_performance_tuning.cpp:100
void cs_user_paramedmem_define_meshes(void)
Define coupled meshes.
Definition: cs_user_paramedmem_coupling.cpp:95
void cs_user_finalize_setup(cs_domain_t *domain)
Define or modify output user parameters. For CDO schemes, specify the elements such as properties,...
Definition: cs_user_parameters.cpp:180
void cs_user_saturne_coupling(void)
Define couplings with other instances of code_saturne.
Definition: cs_user_coupling.cpp:79
void cs_user_turbomachinery(void)
Define rotor/stator model.
Definition: cs_user_turbomachinery.cpp:77
double precision, dimension(:,:,:), allocatable v
Definition: atimbr.f90:113
integer, dimension(:), pointer, save itypfb
boundary condition type at the boundary face ifac (see cs_user_boundary_conditions)
Definition: pointe.f90:69
integer, save ncelet
number of extended (real + ghost of the 'halo') cells. See Note 1: ghost cells - (halos)
Definition: mesh.f90:46
double precision, dimension(:,:), pointer xyzcen
coordinate of the cell centers
Definition: mesh.f90:68
integer, save ncel
number of real cells in the mesh
Definition: mesh.f90:50
Definition: mesh.f90:26
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
Definition: cs_time_control.h:96
time step descriptor
Definition: cs_time_step.h:64
Definition: cs_zone.h:55