8.2
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 extern void CS_PROCF (caltri, CALTRI)
61 (
62  void
63 );
64 
65 /*----------------------------------------------------------------------------
66  * Initialize Fortran base common block values
67  *----------------------------------------------------------------------------*/
68 
69 extern void CS_PROCF (csinit, CSINIT)
70 (
71  const int *irgpar, /* <-- MPI Rank in parallel, -1 otherwise */
72  const int *nrgpar /* <-- Number of MPI processes, or 1 */
73 );
74 
75 /*----------------------------------------------------------------------------
76  * Find the nearest cell's center from a node
77  *----------------------------------------------------------------------------*/
78 
79 extern void CS_PROCF (findpt, FINDPT)
80 (
81  const cs_lnum_t *ncelet, /* <-- number of extended (real + ghost) cells */
82  const cs_lnum_t *ncel, /* <-- number of cells */
83  const cs_real_t *xyzcen, /* <-- cell centers */
84  const cs_real_t *xx, /* <-- node coordinate X */
85  const cs_real_t *yy, /* <-- node coordinate Y */
86  const cs_real_t *zz, /* <-- node coordinate Z */
87  cs_lnum_t *node, /* --> node we are looking for, zero if error */
88  int *ndrang /* --> rank of associated process */
89 );
90 
91 /*----------------------------------------------------------------------------
92  * Main Fortran options initialization
93  *----------------------------------------------------------------------------*/
94 
95 extern void CS_PROCF (initi1, INITI1)
96 (
97  void
98 );
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  * variable number for defined field
109  *----------------------------------------------------------------------------*/
110 
111 int
113 
114 /*----------------------------------------------------------------------------
115  * Add field indexes associated with a new non-user solved variable,
116  * with default options
117  *
118  * parameters:
119  * f_id <-- field id
120  *
121  * returns:
122  * scalar number for defined field
123  *----------------------------------------------------------------------------*/
124 
125 int
127 
128 /*----------------------------------------------------------------------------
129  * Add field indexes associated with a new solved thermal variable,
130  * with default options
131  *
132  * parameters:
133  * f_id <-- field id
134  *----------------------------------------------------------------------------*/
135 
136 void
138 
139 /*----------------------------------------------------------------------------
140  * Computes the explicit chemical source term for atmospheric chemistry in
141  * case of a semi-coupled resolution
142  *----------------------------------------------------------------------------*/
143 
144 void
146  cs_real_t st_exp[],
147  cs_real_t st_imp[]);
148 
149 /*----------------------------------------------------------------------------*/
154 /*----------------------------------------------------------------------------*/
155 
156 void
158 
159 /*----------------------------------------------------------------------------*/
164 /*----------------------------------------------------------------------------*/
165 
166 void
168 
169 /*----------------------------------------------------------------------------*/
192 /*----------------------------------------------------------------------------*/
193 
194 void
196  const cs_real_t voidf[]);
197 
198 /*----------------------------------------------------------------------------*/
208 /*----------------------------------------------------------------------------*/
209 
210 void
212  const cs_lnum_t face_ids[],
213  const cs_real_t t[],
214  cs_real_t h[]);
215 
216 /*----------------------------------------------------------------------------*/
225 /*----------------------------------------------------------------------------*/
226 
227 cs_real_t
229  cs_real_t h);
230 
231 /*----------------------------------------------------------------------------*/
240 /*----------------------------------------------------------------------------*/
241 
242 cs_real_t
244  cs_real_t t);
245 
246 /*----------------------------------------------------------------------------*/
252 /*----------------------------------------------------------------------------*/
253 
254 cs_real_t *
256 
257 /*----------------------------------------------------------------------------*/
263 /*----------------------------------------------------------------------------*/
264 
265 cs_real_t *
267 
268 /*----------------------------------------------------------------------------
269  * Return Lagrangian model status.
270  *
271  * parameters:
272  * model_flag --> 0 without Lagrangian, 1 or 2 with Lagrangian
273  * restart_flag --> 1 for Lagrangian restart, 0 otherwise
274  * frozen_flag --> 1 for frozen Eulerian flow, 0 otherwise
275  *----------------------------------------------------------------------------*/
276 
277 void
278 cs_lagr_status(int *model_flag,
279  int *restart_flag,
280  int *frozen_flag);
281 
282 /*----------------------------------------------------------------------------*/
286 /*----------------------------------------------------------------------------*/
287 
288 void
290 
291 /*----------------------------------------------------------------------------*/
310 /*----------------------------------------------------------------------------*/
311 
312 void
314 
315 /*----------------------------------------------------------------------------
316  * Compute source terms for specific physical model scalars.
317  *----------------------------------------------------------------------------*/
318 
319 void
321  cs_real_t st_imp[],
322  cs_real_t st_exp[]);
323 
324 /*----------------------------------------------------------------------------*/
345 /*----------------------------------------------------------------------------*/
346 
347 void
349  cs_lnum_t icetsm[],
350  cs_real_t dt[]);
351 
352 /*----------------------------------------------------------------------------*/
356 /*----------------------------------------------------------------------------*/
357 
358 void
360 
361 /*----------------------------------------------------------------------------
362  * Exchange of coupling variables between tow instances
363  * of code_saturne thanks to cells.
364  *----------------------------------------------------------------------------*/
365 
366 void
368  cs_real_t st_exp[],
369  cs_real_t st_imp[]);
370 
371 /*----------------------------------------------------------------------------*/
376 /*----------------------------------------------------------------------------*/
377 
378 void
380 
381 /*----------------------------------------------------------------------------*/
387 /*----------------------------------------------------------------------------*/
388 
389 void
391 
392 /*============================================================================
393  * User function prototypes
394  *============================================================================*/
395 
396 /*----------------------------------------------------------------------------
397  * Data Entry of the 1D wall thermal module.
398  *----------------------------------------------------------------------------*/
399 
400 void
401 cs_user_1d_wall_thermal(int iappel,
402  int isuit1);
403 
404 /*----------------------------------------------------------------------------
405  * Data Entry of the wall condensation module
406  *----------------------------------------------------------------------------*/
407 
408 void
410  int nscal,
411  int iappel);
412 
413 /*----------------------------------------------------------------------------
414  * Setup boundary conditions to be applied.
415  *----------------------------------------------------------------------------*/
416 
417 void
419 
420 /*----------------------------------------------------------------------------
421  * This function is called at each time step for boundary conditions.
422  *----------------------------------------------------------------------------*/
423 
424 void
426  int bc_type[]);
427 
428 /*----------------------------------------------------------------------------
429  * Boundary conditions for (Arbitrary Lagrangian Eulerian).
430  *----------------------------------------------------------------------------*/
431 
432 void
434  int bc_type[],
435  int ale_bc_type[],
436  int impale[]);
437 
438 /*----------------------------------------------------------------------------*/
450 /*----------------------------------------------------------------------------*/
451 
452 void
454 
455 /*----------------------------------------------------------------------------*/
464 /*----------------------------------------------------------------------------*/
465 
466 void
468 
469 /*----------------------------------------------------------------------------*/
478 /*----------------------------------------------------------------------------*/
479 
480 void
482 
483 /*----------------------------------------------------------------------------*/
506 /*----------------------------------------------------------------------------*/
507 
508 void
509 cs_user_fsi_structure_define(int is_restart,
510  int n_structs,
511  int *plot,
512  cs_time_control_t *plot_time_control,
513  cs_real_t *aexxst,
514  cs_real_t *bexxst,
515  cs_real_t *cfopre,
516  cs_real_t xstr0[][3],
517  cs_real_t vstr0[][3],
518  cs_real_t xstreq[][3]);
519 
520 /*----------------------------------------------------------------------------*/
535 /*----------------------------------------------------------------------------*/
536 
537 void
538 cs_user_fsi_structure_values(int n_structs,
539  const cs_time_step_t *ts,
540  const cs_real_t xstreq[][3],
541  const cs_real_t xstr[][3],
542  const cs_real_t vstr[][3],
543  cs_real_t xmstru[][3][3],
544  cs_real_t xcstru[][3][3],
545  cs_real_t xkstru[][3][3],
546  cs_real_t forstr[][3],
547  cs_real_t dtstr[]);
548 
549 /*----------------------------------------------------------------------------*/
562 /*----------------------------------------------------------------------------*/
563 
564 void
566  int structure_num[]);
567 
568 /*----------------------------------------------------------------------------*/
578 /*----------------------------------------------------------------------------*/
579 
580 void
581 cs_user_head_losses(const cs_zone_t *zone,
582  cs_real_t cku[][6]);
583 
584 /*----------------------------------------------------------------------------*/
590 /*----------------------------------------------------------------------------*/
591 
592 void
594 
595 /*----------------------------------------------------------------------------*/
606 /*----------------------------------------------------------------------------*/
607 
608 void
610 
611 /*----------------------------------------------------------------------------*/
617 /*----------------------------------------------------------------------------*/
618 
619 void
621 
622 /*----------------------------------------------------------------------------*/
631 /*----------------------------------------------------------------------------*/
632 
633 void
635 
636 /*----------------------------------------------------------------------------*/
645 /*----------------------------------------------------------------------------*/
646 
647 void
649 
650 /*----------------------------------------------------------------------------*/
656 /*----------------------------------------------------------------------------*/
657 
658 void
660 
661 /*----------------------------------------------------------------------------*/
679 /*----------------------------------------------------------------------------*/
680 
681 void
683  const cs_zone_t *z,
684  bool z_local,
685  const cs_real_t h[restrict],
686  cs_real_t t[restrict]);
687 
688 /*----------------------------------------------------------------------------*/
706 /*----------------------------------------------------------------------------*/
707 
708 void
710  const cs_zone_t *z,
711  bool z_local,
712  const cs_real_t t[restrict],
713  cs_real_t h[restrict]);
714 
715 /*----------------------------------------------------------------------------*/
728 /*----------------------------------------------------------------------------*/
729 
730 void
732 
733 /*----------------------------------------------------------------------------*/
739 /*----------------------------------------------------------------------------*/
740 
741 void
743 
744 /*----------------------------------------------------------------------------*/
754 /*----------------------------------------------------------------------------*/
755 
756 void
758  int f_id,
759  cs_real_t *st_exp,
760  cs_real_t *st_imp);
761 
762 /*----------------------------------------------------------------------------*/
772 /*----------------------------------------------------------------------------*/
773 
774 void
776 
777 /*----------------------------------------------------------------------------
778  * Define mesh joinings.
779  *----------------------------------------------------------------------------*/
780 
781 void
782 cs_user_join(void);
783 
784 /*----------------------------------------------------------------------------
785  * Define linear solver options.
786  *
787  * This function is called at the setup stage, once user and most model-based
788  * fields are defined.
789  *----------------------------------------------------------------------------*/
790 
791 void
793 
794 /*----------------------------------------------------------------------------*/
800 /*----------------------------------------------------------------------------*/
801 
802 void
804 
805 /*----------------------------------------------------------------------------
806  * Tag bad cells within the mesh based on geometric criteria.
807  *----------------------------------------------------------------------------*/
808 
809 void
811  cs_mesh_quantities_t *mesh_quantities);
812 
813 /*----------------------------------------------------------------------------
814  * Define restart mode for mesh input and preprocessing.
815  *----------------------------------------------------------------------------*/
816 
817 void
819 
820 /*----------------------------------------------------------------------------
821  * Define mesh files to read and optional associated transformations.
822  *----------------------------------------------------------------------------*/
823 
824 void
825 cs_user_mesh_input(void);
826 
827 /*----------------------------------------------------------------------------
828  * Modifiy geometry and mesh.
829  *----------------------------------------------------------------------------*/
830 
831 void
833 
834 /*----------------------------------------------------------------------------
835  * Insert boundary wall into a mesh.
836  *----------------------------------------------------------------------------*/
837 
838 void
840 
841 /*----------------------------------------------------------------------------
842  * Mesh smoothing.
843  *
844  * parameters:
845  * mesh <-> pointer to mesh structure to smoothe
846  *----------------------------------------------------------------------------*/
847 
848 void
850 
851 /*----------------------------------------------------------------------------
852  * Enable or disable mesh saving.
853  *
854  * By default, mesh is saved when modified.
855  *
856  * parameters:
857  * mesh <-> pointer to mesh structure
858  *----------------------------------------------------------------------------*/
859 
860 void
862 
863 /*----------------------------------------------------------------------------
864  * Set options for cutting of warped faces
865  *
866  * parameters:
867  * mesh <-> pointer to mesh structure to smoothe
868  *----------------------------------------------------------------------------*/
869 
870 void
872 
873 /*----------------------------------------------------------------------------*/
881 /*----------------------------------------------------------------------------*/
882 
883 void
885  cs_mesh_quantities_t *mesh_quantities);
886 
887 /*----------------------------------------------------------------------------*/
891 /*----------------------------------------------------------------------------*/
892 
893 void
895 
896 /*----------------------------------------------------------------------------
897  * Select physical model options, including user fields.
898  *
899  * This function is called at the earliest stages of the data setup.
900  *----------------------------------------------------------------------------*/
901 
902 void
903 cs_user_model(void);
904 
905 /*----------------------------------------------------------------------------
906  * Define advanced mesh numbering options.
907  *----------------------------------------------------------------------------*/
908 
909 void
910 cs_user_numbering(void);
911 
912 /*----------------------------------------------------------------------------
913  * Define parallel IO settings.
914  *----------------------------------------------------------------------------*/
915 
916 void
917 cs_user_parallel_io(void);
918 
919 /*----------------------------------------------------------------------------
920  * Define advanced partitioning options.
921  *----------------------------------------------------------------------------*/
922 
923 void
924 cs_user_partition(void);
925 
926 /*----------------------------------------------------------------------------
927  * User function for input of radiative transfer module options.
928  *
929  * Deprecated Use cs_user_model instead.
930  *----------------------------------------------------------------------------*/
931 
932 void
934 
935 /*----------------------------------------------------------------------------
936  * Define sparse matrix tuning options.
937  *----------------------------------------------------------------------------*/
938 
939 void
941 
942 /*----------------------------------------------------------------------------
943  * Define or modify general numerical and physical user parameters.
944  *
945  * At the calling point of this function, most model-related most variables
946  * and other fields have been defined, so specific settings related to those
947  * fields may be set here.
948  *----------------------------------------------------------------------------*/
949 
950 void
952 
953 /*-----------------------------------------------------------------------------
954  * User subroutine for input of radiative transfer boundary conditions
955  *----------------------------------------------------------------------------*/
956 
957 void
959  const int bc_type[],
960  int isothp[],
961  cs_real_t *tmin,
962  cs_real_t *tmax,
963  cs_real_t *tx,
964  const cs_real_t dt[],
965  const cs_real_t thwall[],
966  const cs_real_t qincid[],
967  cs_real_t hfcnvp[],
968  cs_real_t flcnvp[],
969  cs_real_t xlamp[],
970  cs_real_t epap[],
971  cs_real_t epsp[],
972  cs_real_t textp[]);
973 
974 /*----------------------------------------------------------------------------
975  * Define periodic faces.
976  *----------------------------------------------------------------------------*/
977 
978 void
979 cs_user_periodicity(void);
980 
981 /*----------------------------------------------------------------------------
982  * Define post-processing writers.
983  *
984  * The default output format and frequency may be configured, and additional
985  * post-processing writers allowing outputs in different formats or with
986  * different format options and output frequency than the main writer may
987  * be defined.
988  *----------------------------------------------------------------------------*/
989 
990 void
992 
993 /*-----------------------------------------------------------------------------
994  * Define monitoring probes and profiles. A profile is seen as a set of probes.
995  *----------------------------------------------------------------------------*/
996 
997 void
999 
1000 /*----------------------------------------------------------------------------
1001  * Define post-processing meshes.
1002  *
1003  * The main post-processing meshes may be configured, and additional
1004  * post-processing meshes may be defined as a subset of the main mesh's
1005  * cells or faces (both interior and boundary).
1006  *----------------------------------------------------------------------------*/
1007 
1008 void
1010 
1011 /*----------------------------------------------------------------------------
1012  * User function for output of values on a post-processing mesh.
1013  *----------------------------------------------------------------------------*/
1014 
1015 void
1016 cs_user_postprocess_values(const char *mesh_name,
1017  int mesh_id,
1018  int cat_id,
1019  cs_probe_set_t *probes,
1020  cs_lnum_t n_cells,
1021  cs_lnum_t n_i_faces,
1022  cs_lnum_t n_b_faces,
1023  cs_lnum_t n_vertices,
1024  const cs_lnum_t cell_list[],
1025  const cs_lnum_t i_face_list[],
1026  const cs_lnum_t b_face_list[],
1027  const cs_lnum_t vertex_list[],
1028  const cs_time_step_t *ts);
1029 
1030 /*----------------------------------------------------------------------------
1031  * Override default frequency or calculation end based output.
1032  *
1033  * This allows fine-grained control of activation or deactivation,
1034  *
1035  * parameters:
1036  * nt_max_abs <-- maximum time step number
1037  * nt_cur_abs <-- current time step number
1038  * t_cur_abs <-- absolute time at the current time step
1039  *----------------------------------------------------------------------------*/
1040 
1041 void
1042 cs_user_postprocess_activate(int nt_max_abs,
1043  int nt_cur_abs,
1044  double t_cur_abs);
1045 
1046 /*----------------------------------------------------------------------------
1047  * Absorption coefficient for radiative module
1048  *----------------------------------------------------------------------------*/
1049 
1050 void
1051 cs_user_rad_transfer_absorption(const int bc_type[],
1052  cs_real_t ck[]);
1053 
1054 /*----------------------------------------------------------------------------
1055  * Compute the net radiation flux
1056  *----------------------------------------------------------------------------*/
1057 
1058 void
1060  const cs_real_t twall[],
1061  const cs_real_t qincid[],
1062  const cs_real_t xlam[],
1063  const cs_real_t epa[],
1064  const cs_real_t eps[],
1065  const cs_real_t ck[],
1066  cs_real_t net_flux[]);
1067 
1068 /*----------------------------------------------------------------------------
1069  * Set user solver.
1070  *----------------------------------------------------------------------------*/
1071 
1072 int
1073 cs_user_solver_set(void);
1074 
1075 /*----------------------------------------------------------------------------
1076  * Main call to user solver.
1077  *----------------------------------------------------------------------------*/
1078 
1079 void
1081  const cs_mesh_quantities_t *mesh_quantities);
1082 
1083 /*----------------------------------------------------------------------------
1084  * Define couplings with other instances of code_saturne.
1085  *----------------------------------------------------------------------------*/
1086 
1087 void
1089 
1090 /*----------------------------------------------------------------------------
1091  * Define couplings with SYRTHES code.
1092  *----------------------------------------------------------------------------*/
1093 
1094 void
1096 
1097 /*----------------------------------------------------------------------------*/
1107 /*----------------------------------------------------------------------------*/
1108 
1109 void
1110 cs_user_syrthes_coupling_volume_h(int coupling_id,
1111  const char *syrthes_name,
1112  cs_lnum_t n_elts,
1113  const cs_lnum_t elt_ids[],
1114  cs_real_t h_vol[]);
1115 
1116 /*----------------------------------------------------------------------------
1117  * Define couplings with CATHARE code.
1118  *----------------------------------------------------------------------------*/
1119 
1120 void
1122 
1123 /*----------------------------------------------------------------------------
1124  * Define time moments.
1125  *----------------------------------------------------------------------------*/
1126 
1127 void
1128 cs_user_time_moments(void);
1129 
1130 
1131 /*----------------------------------------------------------------------------
1132  * Define time tables.
1133  *----------------------------------------------------------------------------*/
1134 
1135 void
1136 cs_user_time_table(void);
1137 
1138 /*----------------------------------------------------------------------------
1139  * Define rotor/stator model.
1140  *----------------------------------------------------------------------------*/
1141 
1142 void
1144 
1145 /*----------------------------------------------------------------------------
1146  * Define rotor axes, associated cells, and rotor/stator faces.
1147  *----------------------------------------------------------------------------*/
1148 
1149 void
1151 
1152 /*----------------------------------------------------------------------------
1153  * Define rotation velocity of rotor.
1154  *----------------------------------------------------------------------------*/
1155 
1156 void
1158 
1159 /*----------------------------------------------------------------------------*/
1163 /*----------------------------------------------------------------------------*/
1164 
1165 void
1166 cs_user_zones(void);
1167 
1168 /*----------------------------------------------------------------------------*/
1172 /*----------------------------------------------------------------------------*/
1173 
1174 void
1176  const cs_mesh_quantities_t *mesh_quantities,
1177  cs_real_t *dt);
1178 
1179 /*----------------------------------------------------------------------------
1180  * Computation of the relaxation time-scale to equilibrium in the frame of
1181  * the homogeneous two-phase model.
1182  *----------------------------------------------------------------------------*/
1183 
1184 void
1186  const cs_real_t *alpha_eq,
1187  const cs_real_t *y_eq,
1188  const cs_real_t *z_eq,
1189  const cs_real_t *ei,
1190  const cs_real_t *v,
1191  cs_real_t *relax_tau);
1192 
1193 /*----------------------------------------------------------------------------*/
1197 /*----------------------------------------------------------------------------*/
1198 
1199 void
1201 
1202 /*----------------------------------------------------------------------------*/
1206 /*----------------------------------------------------------------------------*/
1207 
1208 void
1210 
1211 /*----------------------------------------------------------------------------*/
1215 /*----------------------------------------------------------------------------*/
1216 
1217 void
1219 
1220 /*----------------------------------------------------------------------------*/
1221 
1223 
1224 #endif /* __CS_PROTOTYPES_H__ */
#define restrict
Definition: cs_defs.h:141
#define BEGIN_C_DECLS
Definition: cs_defs.h:528
double cs_real_t
Floating-point value.
Definition: cs_defs.h:332
#define CS_PROCF(x, y)
Definition: cs_defs.h:560
cs_real_t cs_real_3_t[3]
vector of 3 floating-point values
Definition: cs_defs.h:347
#define END_C_DECLS
Definition: cs_defs.h:529
int cs_lnum_t
local mesh entity id
Definition: cs_defs.h:325
@ t
Definition: cs_field_pointer.h:92
@ xlam
Definition: cs_field_pointer.h:183
@ eps
Definition: cs_field_pointer.h:71
@ h
Definition: cs_field_pointer.h:91
@ epa
Definition: cs_field_pointer.h:184
@ 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.c: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.c:105
void cs_coal_bt2h(cs_lnum_t n_faces, const cs_lnum_t face_ids[], const cs_real_t t[], cs_real_t h[])
Convert temperature to enthalpy at boundary for coal combustion.
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.c:219
void cs_user_wall_condensation(int nvar, int nscal, int iappel)
Source terms associated at the boundary faces and the neighboring cells with surface condensation.
Definition: cs_user_wall_condensation.c:169
void cs_user_mesh_boundary(cs_mesh_t *mesh)
Insert boundaries into a mesh.
Definition: cs_user_mesh.c:156
void cs_user_zones(void)
Define volume and surface zones.
Definition: cs_user_zones.c:65
void cs_user_model(void)
Select physical model options, including user fields.
Definition: cs_user_parameters.c: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.c: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.c: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.c: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.c: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.c:114
void cs_user_syrthes_coupling(void)
Define couplings with SYRTHES code.
Definition: cs_user_coupling.c:95
void cs_user_parameters(cs_domain_t *domain)
Define or modify general numerical and physical user parameters.
Definition: cs_user_parameters.c: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.c: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.c:179
void cs_user_physical_properties_turb_viscosity(cs_domain_t *domain)
User modification of the turbulence viscosity.
Definition: cs_user_physical_properties.c:157
cs_real_t cs_gas_combustion_t_to_h(const cs_real_t x_sp[restrict], cs_real_t t)
Convert a temperature to enthalpy value for gas combustion.
void cs_user_boundary_conditions_setup(cs_domain_t *domain)
Set boundary conditions to be applied.
Definition: cs_user_boundary_conditions.c:77
void cs_user_turbomachinery_rotor(void)
Define rotor axes, associated cells, and rotor/stator faces.
Definition: cs_user_turbomachinery.c:90
void cs_field_map_and_init_bcs(void)
Initialize all field BC coefficients.
Definition: cs_field_default.c:514
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.c:103
void cs_user_postprocess_activate(int nt_max_abs, int nt_cur_abs, double t_cur_abs)
Definition: cs_user_postprocess.c:177
void cs_user_linear_solvers(void)
Define linear solver options.
Definition: cs_user_parameters.c: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.c:86
void cs_user_periodicity(void)
Define periodic faces.
Definition: cs_user_mesh.c: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 initi1(void)
Definition: initi1.f90:29
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.c:2954
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.c:75
void cs_user_mesh_cartesian_define(void)
Define a cartesian mesh.
Definition: cs_user_mesh.c:253
void cs_fortran_resize_aux_arrays(void)
resize some Fortran auxiliary arrays
void cs_user_mesh_modify(cs_mesh_t *mesh)
Modify geometry and mesh.
Definition: cs_user_mesh.c: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.c: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.c:124
void cs_user_postprocess_probes(void)
Define monitoring probes and profiles.
Definition: cs_user_postprocess.c:101
void cs_user_numbering(void)
Define advanced mesh numbering options.
Definition: cs_user_performance_tuning.c:74
void cs_user_1d_wall_thermal(int iappel, int isuit1)
Definition: cs_user_1d_wall_thermal.c:82
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[restrict], cs_real_t t[restrict])
User definition of enthalpy to temperature conversion.
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.c:134
void cs_user_time_moments(void)
Define time moments.
Definition: cs_user_parameters.c:147
void cs_atmo_chem_source_terms(int iscal, cs_real_t st_exp[], cs_real_t st_imp[])
cs_real_t cs_gas_combustion_h_to_t(const cs_real_t x_sp[restrict], cs_real_t h)
Convert an enthalpy to temperature value for gas combustion.
void cs_user_radiative_transfer_parameters(void)
User function for input of radiative transfer module options.
Definition: cs_user_radiative_transfer.c:79
int cs_user_solver_set(void)
Set user solver.
Definition: cs_user_solver.c: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.c: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.c:77
int cs_add_model_field_indexes(int f_id)
void cs_wall_condensation_volume_exchange_surf_at_cells(cs_real_t *surf)
Return condensing volume structures surface at each cell.
Definition: cs_wall_condensation.c:1077
void cs_user_mesh_restart_mode(void)
Force preprocessing behavior in case of restart.
Definition: cs_user_mesh.c:89
void cs_user_mesh_save(cs_mesh_t *mesh)
Enable or disable mesh saving.
Definition: cs_user_mesh.c: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.c: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.c:155
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.c: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.c:113
void cs_user_paramedmem_define_fields(void)
Define fields to couple with ParaMEDMEM.
Definition: cs_user_paramedmem_coupling.c:108
void cs_user_internal_coupling_add_volumes(cs_mesh_t *mesh)
Define volumes as internal coupling zones.
Definition: cs_internal_coupling.c:2935
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.c: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.c:75
void cs_cavitation_compute_source_term(const cs_real_t pressure[], const cs_real_t voidf[])
Compute the vaporization source term using the Merkle model:
cs_real_t * cs_get_cavitation_gam(void)
Return pointer to cavitation "gamcav" array.
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.c: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.c: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.c:92
cs_real_t * cs_get_cavitation_dgdp_st(void)
Return pointer to cavitation "dgdpca" array.
void cs_user_postprocess_writers(void)
Define post-processing writers.
Definition: cs_user_postprocess.c: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.c:135
void cs_user_boundary_conditions(cs_domain_t *domain, int bc_type[])
User definition of boundary conditions.
Definition: cs_user_boundary_conditions.c:123
void cs_user_mesh_smoothe(cs_mesh_t *mesh)
Mesh smoothing.
Definition: cs_user_mesh.c:186
void cs_user_initialization(cs_domain_t *domain)
This function is called one time step to initialize problem.
Definition: cs_user_initialization.c:109
void cs_user_extra_operations_initialize(cs_domain_t *domain)
Initialize variables.
Definition: cs_user_extra_operations.c:86
void cs_user_cathare_coupling(void)
Define couplings with SYRTHES code.
Definition: cs_user_coupling.c:138
void caltri(void)
Definition: caltri.f90:36
void csinit(const int *irgpar, const int *nrgpar)
void cs_user_paramedmem_define_couplings(void)
Define ParaMEDMEM coupling(s)
Definition: cs_user_paramedmem_coupling.c: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.c:82
void cs_user_mesh_warping(void)
Set options for cutting of warped faces.
Definition: cs_user_mesh.c:141
void cs_sat_coupling_exchange_at_cells(int f_id, cs_real_t st_exp[], cs_real_t st_imp[])
void cs_user_partition(void)
Define advanced partitioning options.
Definition: cs_user_performance_tuning.c:87
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[restrict], cs_real_t h[restrict])
User definition of temperature to enthalpy conversion.
void cs_user_parallel_io(void)
Define parallel IO settings.
Definition: cs_user_performance_tuning.c:100
void cs_turbomachinery_update(void)
Sync turbomachinery module components to global c turbomachinery structure.
void cs_user_paramedmem_define_meshes(void)
Define coupled meshes.
Definition: cs_user_paramedmem_coupling.c: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.c:180
void cs_user_saturne_coupling(void)
Define couplings with other instances of code_saturne.
Definition: cs_user_coupling.c:79
void cs_user_turbomachinery(void)
Define rotor/stator model.
Definition: cs_user_turbomachinery.c:77
double precision, dimension(:,:,:), allocatable v
Definition: atimbr.f90:113
double precision, dimension(:,:,:), allocatable pressure
Definition: atimbr.f90:121
integer, dimension(:), allocatable, target icetsm
number of the ncetsm cells in which a mass injection is imposed. See iicesm and the cs_equation_add_v...
Definition: pointe.f90:176
integer, dimension(:), pointer, save itypfb
boundary condition type at the boundary face ifac (see cs_user_boundary_conditions)
Definition: pointe.f90:89
integer, save nscal
number of solved user scalars effective number of scalars solutions of an advection equation,...
Definition: dimens.f90:55
integer, save nvar
number of solved variables (must be lower than nvarmx)
Definition: dimens.f90:42
integer, save ncelet
number of extended (real + ghost of the 'halo') cells. See Note 1: ghost cells - (halos)
Definition: mesh.f90:48
double precision, dimension(:,:), pointer xyzcen
coordinate of the cell centers
Definition: mesh.f90:70
integer, save ncel
number of real cells in the mesh
Definition: mesh.f90:52
integer, save isuit1
For the 1D wall thermal module, activation (1) or not(0) of the reading of the mesh and of the wall t...
Definition: optcal.f90:174
Definition: mesh.f90:26
Structure storing the main features of the computational domain and pointers to the main geometrical ...
Definition: cs_domain.h:138
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