Code_Saturne
CFD tool
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-2012 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_mesh.h"
36 #include "cs_mesh_quantities.h"
37 #include "cs_mesh_bad_cells.h"
38 
39 /*----------------------------------------------------------------------------*/
40 
42 
43 /*============================================================================
44  * Macro definitions
45  *============================================================================*/
46 
47 /*=============================================================================
48  * Fortran function/subroutine prototypes
49  *============================================================================*/
50 
51 /*----------------------------------------------------------------------------
52  * Main Fortran subroutine
53  *----------------------------------------------------------------------------*/
54 
55 extern void CS_PROCF (caltri, CALTRI)
56 (
57  const cs_int_t *iverif /* <-- activate elementary tests */
58 );
59 
60 /*----------------------------------------------------------------------------
61  * Close log (listing) handled by Fortran: (CLose LIsting)
62  *----------------------------------------------------------------------------*/
63 
64 extern void CS_PROCF (csclli, CSCLLI)
65 (
66  void
67 );
68 
69 /*----------------------------------------------------------------------------
70  * Flush standard output.
71  *----------------------------------------------------------------------------*/
72 
73 extern void CS_PROCF (csflsh, CSFLSH)
74 (
75  void
76 );
77 
78 /*----------------------------------------------------------------------------
79  * Initialize Fortran base common bloc values
80  *----------------------------------------------------------------------------*/
81 
82 extern void CS_PROCF (csinit, CSINIT)
83 (
84  const cs_int_t *irgpar, /* <-- MPI Rank in parallel, -1 otherwise */
85  const cs_int_t *nrgpar, /* <-- Number of MPI processes, or 1 */
86  const cs_int_t *nthpar /* <-- Number of threads */
87 );
88 
89 /*----------------------------------------------------------------------------
90  * Initialize Fortran log (listing) files
91  *----------------------------------------------------------------------------*/
92 
93 extern void CS_PROCF (csopli, CSOPLI)
94 (
95  const cs_int_t *irkpar, /* <-- MPI Rank in parallel, -1 otherwise */
96  const cs_int_t *nrkpar, /* <-- Number of MPI processes, or 1 */
97  const cs_int_t *ilogr0, /* <-- Output of main log (listing (rank 0): */
98  /* 0: non redirected; 1: to 'listing' file */
99  const cs_int_t *ilogrp /* <-- Output of logs for ranks > 0: */
100  /* 0: non redirected; 1: to 'listing_n*' files */
101  /* 2: to '/dev/null' (suppressed) */
102 );
103 
104 /*----------------------------------------------------------------------------
105  * Print a message to standard output.
106  *----------------------------------------------------------------------------*/
107 
108 extern void CS_PROCF (csprnt, CSPRNT)
109 (
110  char *cs_buf_print,
111  cs_int_t *msgsize
112 );
113 
114 /*----------------------------------------------------------------------------
115  * Developper function for output of variables on a post-processing mesh
116  *----------------------------------------------------------------------------*/
117 
118 extern void CS_PROCF (dvvpst, DVVPST)
119 (
120  const cs_int_t *nummai, /* <-- number or post-processing mesh */
121  const cs_int_t *numtyp, /* <-- number or post-processing type
122  * (-1 as volume, -2 as boundary, or nummai) */
123  const cs_int_t *nvar, /* <-- number of variables */
124  const cs_int_t *nscal, /* <-- number of scalars */
125  const cs_int_t *nvlsta, /* <-- number of statistical variables (lagr) */
126  const cs_int_t *nvisbr, /* <-- number of boundary stat. variables (lagr) */
127  const cs_int_t *ncelps, /* <-- number of post-processed cells */
128  const cs_int_t *nfacps, /* <-- number of post processed interior faces */
129  const cs_int_t *nfbrps, /* <-- number of post processed boundary faces */
130  const cs_int_t itypps[3], /* <-- flag (0 or 1) for presence of cells, */
131  const cs_int_t lstcel[], /* <-- list of post-processed cells */
132  const cs_int_t lstfac[], /* <-- list of post-processed interior faces */
133  const cs_int_t lstfbr[], /* <-- list of post-processed boundary faces */
134  const cs_real_t dt[], /* <-- local time step */
135  const cs_real_t rtpa[], /* <-- cell variables at previous time step */
136  const cs_real_t rtp[], /* <-- cell variables */
137  const cs_real_t propce[], /* <-- cell physical properties */
138  const cs_real_t propfa[], /* <-- interior face physical properties */
139  const cs_real_t propfb[], /* <-- boundary face physical properties */
140  const cs_real_t coefa[], /* <-- boundary conditions array */
141  const cs_real_t coefb[], /* <-- boundary conditions array */
142  const cs_real_t statce[], /* <-- cell statistics (Lagrangian) */
143  const cs_real_t stativ[], /* <-- cell variance statistics (Lagrangian) */
144  const cs_real_t statfb[], /* <-- boundary face statistics (Lagrangian) */
145  cs_real_t tracel[], /* --- work array for output cells */
146  cs_real_t trafac[], /* --- work array for output interior faces */
147  cs_real_t trafbr[] /* --- work array for output boundary faces */
148 );
149 
150 /*----------------------------------------------------------------------------
151  * Find the nearest cell's center from a node
152  *----------------------------------------------------------------------------*/
153 
154 extern void CS_PROCF (findpt, FINDPT)
155 (
156  const cs_int_t *ncelet, /* <-- number of extended (real + ghost) cells */
157  const cs_int_t *ncel, /* <-- number of cells */
158  const cs_real_t *xyzcen, /* <-- cell centers */
159  const cs_real_t *xx, /* <-- node coordinate X */
160  const cs_real_t *yy, /* <-- node coordinate Y */
161  const cs_real_t *zz, /* <-- node coordinate Z */
162  cs_int_t *node, /* --> node we are looking for, zero if error */
163  cs_int_t *ndrang /* --> rank of associated process */
164 );
165 
166 /*----------------------------------------------------------------------------
167  * Main Fortran options initialization
168  *----------------------------------------------------------------------------*/
169 
170 extern void CS_PROCF (initi1, INITI1)
171 (
172  const cs_int_t *iverif /* <-- Activate elementary tests */
173 );
174 
175 /*----------------------------------------------------------------------------
176  * Update mesh dimensions in Fortran commons
177  *----------------------------------------------------------------------------*/
178 
179 extern void CS_PROCF (majgeo, MAJGEO)
180 (
181  const cs_int_t *ncel, /* <-- number of cells */
182  const cs_int_t *ncelet, /* <-- number of extended (real + ghost) cells */
183  const cs_int_t *nfac, /* <-- number of interior faces */
184  const cs_int_t *nfabor, /* <-- number of boundary faces */
185  const cs_int_t *nsom, /* <-- number of vertices */
186  const cs_int_t *lndfac, /* <-- interior face -> vertices array size */
187  const cs_int_t *lndfbr, /* <-- boundary face -> vertices array size */
188  const cs_int_t *ncelbr, /* <-- number of boundary cells */
189  const cs_int_t *ncelgb, /* <-- global number of cells */
190  const cs_int_t *nfacgb, /* <-- global number of interior faces */
191  const cs_int_t *nfbrgb, /* <-- global number of boundary faces */
192  const cs_int_t *nsomgb, /* <-- global number of vertices */
193  const cs_int_t *nthrdi, /* <-- max threads per interior faces group */
194  const cs_int_t *nthrdb, /* <-- max threads per boundary faces group */
195  const cs_int_t *ngrpi, /* <-- number of interior face groups */
196  const cs_int_t *ngrpb, /* <-- number of boundary face groups */
197  const cs_int_t *idxfi, /* <-- interior face group/thread start/end ids */
198  const cs_int_t *idxfb, /* <-- boundary face group/thread start/end ids */
199  const cs_int_t ifacel[], /* <-- interior faces -> cells connectivity */
200  const cs_int_t ifabor[], /* <-- boundary faces -> cells connectivity */
201  const cs_int_t ifmfbr[], /* <-- boundary face families */
202  const cs_int_t ifmcel[], /* <-- cell families */
203  const cs_int_t iprfml[], /* <-- list of family properties */
204  const cs_int_t ipnfac[], /* <-- interior faces -> vertices index */
205  const cs_int_t nodfac[], /* <-- interior faces -> vertices connectivity */
206  const cs_int_t ipnfbr[], /* <-- boundary faces -> vertices index */
207  const cs_int_t nodfbr[], /* <-- boundary faces -> vertices connectivity */
208  const cs_int_t icelbr[], /* <-- list of boundary cells */
209  const cs_real_t *volmin, /* <-- minimum control volume */
210  const cs_real_t *volmax, /* <-- maximum control volume */
211  const cs_real_t *voltot, /* <-- total control volume */
212  const cs_real_t xyzcen[], /* <-- cell centers */
213  const cs_real_t surfac[], /* <-- interior face surface vectors */
214  const cs_real_t surfbo[], /* <-- boundary face surface vectors */
215  const cs_real_t cdgfac[], /* <-- interior face centers */
216  const cs_real_t cdgfbo[], /* <-- boundary face centers */
217  const cs_real_t xyznod[], /* <-- vertex coordinates */
218  const cs_real_t volume[], /* <-- cell volumes */
219  const cs_real_t surfan[], /* <-- interior face surfaces */
220  const cs_real_t surfbn[], /* <-- boundary face surfaces */
221  const cs_real_t dist[], /* <-- distance IJ.Nij */
222  const cs_real_t distb[], /* <-- likewise for border faces */
223  const cs_real_t pond[], /* <-- weighting (Aij=pond Ai+(1-pond)Aj */
224  const cs_real_t dijpf[], /* <-- vector I'J' */
225  const cs_real_t diipb[], /* <-- likewise for border faces */
226  const cs_real_t dofij[] /* <-- vector OF at interior faces */
227 );
228 
229 /*----------------------------------------------------------------------------
230  * Free Fortran allocated memory
231  *----------------------------------------------------------------------------*/
232 
233 extern void CS_PROCF (memfin, MEMFIN) (void);
234 
235 /*----------------------------------------------------------------------------
236  * Mesh renumbering for vector processors
237  *----------------------------------------------------------------------------*/
238 
239 extern void CS_PROCF (numvec, NUMVEC)
240 (
241  const cs_int_t *ncelet, /* <-- number of extended (real + ghost) cells */
242  const cs_int_t *ncel, /* <-- number of local cells */
243  const cs_int_t *nfac, /* <-- number of interior faces */
244  const cs_int_t *nfabor, /* <-- number of boundary faces */
245  const cs_int_t *lregis, /* <-- vector registor length */
246  cs_int_t *irveci, /* <-> interior face vectorization indic. */
247  cs_int_t *irvecb, /* <-> boundary face vectorization indic. */
248  cs_int_t ifacel[], /* <-- interior face->cell connectivity */
249  cs_int_t ifabor[], /* <-- boundary face->cell connectivity */
250  cs_int_t inumfi[], /* <-> interior faces renumbering (size: nfac) */
251  cs_int_t inumfb[], /* <-> boundary faces renumbering (size: nfabor) */
252  cs_int_t iworkf[], /* --- work array, size: max(nfac, nfabor) */
253  cs_int_t ismbs[] /* --- work array, size: ncelet */
254 );
255 
256 /*----------------------------------------------------------------------------
257  * Test renumbering for vector processors
258  *----------------------------------------------------------------------------*/
259 
260 extern void CS_PROCF (tstvec, TSTVEC)
261 (
262  const cs_int_t *ncelet, /* <-- number of extended (real + ghost) cells */
263  const cs_int_t *ncel, /* <-- number of local cells */
264  const cs_int_t *nfac, /* <-- number of interior faces */
265  const cs_int_t *nfabor, /* <-- number of boundary faces */
266  const cs_int_t ifacel[], /* <-- interior face->cell connectivity */
267  const cs_int_t ifabor[], /* <-- boundary face->cell connectivity */
268  cs_int_t iworkf[], /* --- work array, size: max(nfac, nfabor) */
269  cs_int_t ismbs[], /* --- work array, size: ncelet */
270  cs_int_t ismbv[], /* --- work array, size: ncelet */
271  cs_real_t rworkf[], /* --- work array, size: max(nfac, nfabor) */
272  cs_real_t rsmbs[], /* --- work array, size: ncelet */
273  cs_real_t rsmbv[] /* --- work array, size: ncelet */
274 );
275 
276 /*----------------------------------------------------------------------------
277  * User override of default frequency or calculation end based output.
278  *
279  * Fortran interface:
280  *
281  * subroutine pstusn (ntmabs, ntcabs, ttcabs)
282  * *****************
283  *
284  * integer ntmabs : <-- : maximum time step number
285  * integer ntcabs : <-- : current time step number
286  * double precision ttcabs : <-- : absolute time at the current time step
287  *----------------------------------------------------------------------------*/
288 
289 void CS_PROCF (pstusn, PSTUSN)
290 (
291  const cs_int_t *ntmabs,
292  const cs_int_t *ntcabs,
293  const cs_real_t *ttcabs
294 );
295 
296 /*----------------------------------------------------------------------------
297  * Indicate if the variable considered is a component of a vector or tensor
298  * in the presence of periodicity of rotation
299  *
300  * Fortran interface:
301  *
302  * subroutine pergra (ivar, ipvar)
303  * *****************
304  *
305  * integer ivar : <-- : variable number
306  * integer idimtr : --> : 0 if ivar does not match a vector or tensor
307  * : : or there is no periodicity of rotation
308  * : : 1 for velocity, 2 for Reynolds stress
309  * : : in case of periodicity of rotation
310  * integer irpvar : --> : -1 if ivar does not match a vector or tensor
311  * : : or there is no periodicity of rotation
312  * : : In presence of periodicity of rotation:
313  * : : 0 for iu, 1 for iv, 2 for iw
314  * : : 0 for ir11, 1 for ir22, 2 for ir33,
315  * : : 3 for ir12, 4 for ir13, 5 for ir23
316  *----------------------------------------------------------------------------*/
317 
318 void CS_PROCF (pergra, PERGRA)
319 (
320  const cs_int_t *ivar,
321  const cs_int_t *idimtr,
322  const cs_int_t *irpvar
323 );
324 
325 /*----------------------------------------------------------------------------
326  * User function for output of variables on a post-processing mesh
327  *----------------------------------------------------------------------------*/
328 
329 void CS_PROCF (usvpst, USVPST)
330 (
331  const cs_int_t *nummai, /* <-- number or post-processing mesh */
332  const cs_int_t *nvar, /* <-- number of variables */
333  const cs_int_t *nscal, /* <-- number of scalars */
334  const cs_int_t *nvlsta, /* <-- number of statistical variables (lagr) */
335  const cs_int_t *ncelps, /* <-- number of post-processed cells */
336  const cs_int_t *nfacps, /* <-- number of post processed interior faces */
337  const cs_int_t *nfbrps, /* <-- number of post processed boundary faces */
338  const cs_int_t itypps[3], /* <-- flag (0 or 1) for presence of cells, */
339  /* interior faces, and boundary faces */
340  const cs_int_t lstcel[], /* <-- list of post-processed cells */
341  const cs_int_t lstfac[], /* <-- list of post-processed interior faces */
342  const cs_int_t lstfbr[], /* <-- list of post-processed boundary faces */
343  const cs_real_t dt[], /* <-- local time step */
344  const cs_real_t rtpa[], /* <-- cell variables at previous time step */
345  const cs_real_t rtp[], /* <-- cell variables */
346  const cs_real_t propce[], /* <-- cell physical properties */
347  const cs_real_t propfa[], /* <-- interior face physical properties */
348  const cs_real_t propfb[], /* <-- boundary face physical properties */
349  const cs_real_t coefa[], /* <-- boundary conditions array */
350  const cs_real_t coefb[], /* <-- boundary conditions array */
351  const cs_real_t statce[], /* <-- cell statistics (Lagrangian) */
352  cs_real_t tracel[], /* --- work array for output cells */
353  cs_real_t trafac[], /* --- work array for output interior faces */
354  cs_real_t trafbr[] /* --- work array for output boundary faces */
355 );
356 
357 /*----------------------------------------------------------------------------
358  * Uniform random number generator
359  *----------------------------------------------------------------------------*/
360 
361 void CS_PROCF (zufall, zufall)
362 (
363  const cs_int_t *n, /* --> size of the vector */
364  const cs_real_t *a /* <-- generated random number vector */
365 );
366 
367 /*----------------------------------------------------------------------------
368  * Gaussian random number generator
369  *----------------------------------------------------------------------------*/
370 
371 void CS_PROCF (normalen, normalen)
372 (
373  const cs_int_t *n, /* --> size of the vector */
374  const cs_real_t *x /* <-- generated random number vector */
375 );
376 
377 /*============================================================================
378  * User function prototypes
379  *============================================================================*/
380 
381 /*----------------------------------------------------------------------------
382  * Define global options for couplings.
383  *
384  * These options allow defining the time step synchronization policy,
385  * as well as a time step multiplier.
386  *----------------------------------------------------------------------------*/
387 
388 void
389 cs_user_coupling(void);
390 
391 /*----------------------------------------------------------------------------
392  * Define mesh joinings.
393  *----------------------------------------------------------------------------*/
394 
395 void
396 cs_user_join(void);
397 
398 /*----------------------------------------------------------------------------
399  * Define mesh files to read and optional associated transformations.
400  *----------------------------------------------------------------------------*/
401 
402 void
403 cs_user_mesh_input(void);
404 
405 /*----------------------------------------------------------------------------
406  * Modifiy geometry and mesh.
407  *----------------------------------------------------------------------------*/
408 
409 void
411 
412 /*----------------------------------------------------------------------------
413  * Insert thin wall into a mesh.
414  *----------------------------------------------------------------------------*/
415 
416 void
418 
419 /*----------------------------------------------------------------------------
420  * Mesh smoothing.
421  *
422  * parameters:
423  * mesh <-> pointer to mesh structure to smoothe
424  *----------------------------------------------------------------------------*/
425 
426 void
428 
429 /*----------------------------------------------------------------------------
430  * Enable or disable mesh saving.
431  *
432  * By default, mesh is saved when modified.
433  *
434  * parameters:
435  * mesh <-> pointer to mesh structure
436  *----------------------------------------------------------------------------*/
437 
438 void
440 
441 /*----------------------------------------------------------------------------
442  * Tag bad cells within the mesh based on geometric criteria.
443  *----------------------------------------------------------------------------*/
444 
445 void
447  cs_mesh_quantities_t *mesh_quantities);
448 
449 /*----------------------------------------------------------------------------
450  * Define advanced partitioning options.
451  *----------------------------------------------------------------------------*/
452 
453 void
454 cs_user_partition(void);
455 
456 /*----------------------------------------------------------------------------
457  * Define periodic faces.
458  *----------------------------------------------------------------------------*/
459 
460 void
461 cs_user_periodicity(void);
462 
463 /*----------------------------------------------------------------------------
464  * Define post-processing writers.
465  *
466  * The default output format and frequency may be configured, and additional
467  * post-processing writers allowing outputs in different formats or with
468  * different format options and output frequency than the main writer may
469  * be defined.
470  *----------------------------------------------------------------------------*/
471 
472 void
474 
475 /*----------------------------------------------------------------------------
476  * Define post-processing meshes.
477  *
478  * The main post-processing meshes may be configured, and additional
479  * post-processing meshes may be defined as a subset of the main mesh's
480  * cells or faces (both interior and boundary).
481  *----------------------------------------------------------------------------*/
482 
483 void
485 
486 /*----------------------------------------------------------------------------
487  * Override default frequency or calculation end based output.
488  *
489  * This allows fine-grained control of activation or deactivation,
490  *
491  * parameters:
492  * nt_max_abs <-- maximum time step number
493  * nt_cur_abs <-- current time step number
494  * t_cur_abs <-- absolute time at the current time step
495  *----------------------------------------------------------------------------*/
496 
497 void
498 cs_user_postprocess_activate(int nt_max_abs,
499  int nt_cur_abs,
500  double t_cur_abs);
501 
502 /*----------------------------------------------------------------------------
503  * Define couplings with other instances of Code_Saturne.
504  *----------------------------------------------------------------------------*/
505 
506 void
508 
509 /*----------------------------------------------------------------------------
510  * Set user solver.
511  *----------------------------------------------------------------------------*/
512 
513 int
514 cs_user_solver_set(void);
515 
516 /*----------------------------------------------------------------------------
517  * Main call to user solver.
518  *----------------------------------------------------------------------------*/
519 
520 void
522  const cs_mesh_quantities_t *mesh_quantities);
523 
524 /*----------------------------------------------------------------------------
525  * Define couplings with SYRTHES code.
526  *----------------------------------------------------------------------------*/
527 
528 void
530 
531 /*----------------------------------------------------------------------------*/
532 
534 
535 #endif /* __CS_PROTOTYPES_H__ */
void cs_user_partition(void)
Definition: cs_user_parallel.c:84
void tstvec(const cs_int_t *ncelet, const cs_int_t *ncel, const cs_int_t *nfac, const cs_int_t *nfabor, const cs_int_t ifacel[], const cs_int_t ifabor[], cs_int_t iworkf[], cs_int_t ismbs[], cs_int_t ismbv[], cs_real_t rworkf[], cs_real_t rsmbs[], cs_real_t rsmbv[])
integer, save ngrpb
Definition: parall.f90:47
integer, save ngrpi
Definition: parall.f90:47
double precision, dimension(:), pointer pond
Definition: mesh.f90:113
double precision, dimension(:,:), pointer xyzcen
Definition: mesh.f90:102
void cs_user_periodicity(void)
Definition: cs_user_mesh.c:292
void csflsh(void)
Definition: csflsh.f90:24
void csopli(const cs_int_t *irkpar, const cs_int_t *nrkpar, const cs_int_t *ilogr0, const cs_int_t *ilogrp)
double precision, dimension(:), pointer surfbn
Definition: mesh.f90:110
integer, dimension(:), pointer nodfbr
Definition: mesh.f90:72
double precision, dimension(:), pointer dist
Definition: mesh.f90:111
integer, save nscal
Definition: dimens.f90:33
void zufall(const cs_int_t *n, const cs_real_t *a)
void cs_user_mesh_thinwall(cs_mesh_t *mesh)
Definition: cs_user_mesh.c:475
integer, save nvlsta
Definition: lagdim.f90:63
#define BEGIN_C_DECLS
Definition: cs_defs.h:365
void cs_user_postprocess_activate(int nt_max_abs, int nt_cur_abs, double t_cur_abs)
Definition: cs_user_postprocess.c:630
void cs_user_postprocess_meshes(void)
Definition: cs_user_postprocess.c:390
double precision, save volmax
Definition: cstphy.f90:141
void cs_user_coupling(void)
Definition: cs_user_coupling.c:75
void cs_user_mesh_modify(cs_mesh_t *mesh)
Definition: cs_user_mesh.c:529
void cs_user_join(void)
Definition: cs_user_mesh.c:152
double precision, dimension(:,:), pointer xyznod
Definition: mesh.f90:107
integer(kind=8), save ncelgb
Definition: parall.f90:58
int cs_int_t
Definition: cs_defs.h:263
integer, dimension(:), pointer ifmcel
Definition: mesh.f90:79
integer, save ntcabs
Definition: optcal.f90:277
integer, dimension(:,:), pointer ifacel
Definition: mesh.f90:61
double precision, dimension(:,:), pointer cdgfac
Definition: mesh.f90:105
double precision, dimension(:,:), pointer diipb
Definition: mesh.f90:115
integer, save nfabor
Definition: mesh.f90:42
void cs_user_mesh_save(cs_mesh_t *mesh)
Definition: cs_user_mesh.c:606
void numvec(const cs_int_t *ncelet, const cs_int_t *ncel, const cs_int_t *nfac, const cs_int_t *nfabor, const cs_int_t *lregis, cs_int_t *irveci, cs_int_t *irvecb, cs_int_t ifacel[], cs_int_t ifabor[], cs_int_t inumfi[], cs_int_t inumfb[], cs_int_t iworkf[], cs_int_t ismbs[])
Definition: cs_mesh.h:62
double precision, dimension(:,:), pointer surfbo
Definition: mesh.f90:104
double precision, dimension(:), pointer volume
Definition: mesh.f90:108
BEGIN_C_DECLS void caltri(const cs_int_t *iverif)
double precision, dimension(:,:), pointer cdgfbo
Definition: mesh.f90:106
void csprnt(char *cs_buf_print, cs_int_t *msgsize)
integer, save nfac
Definition: mesh.f90:42
void cs_user_solver(const cs_mesh_t *mesh, const cs_mesh_quantities_t *mesh_quantities)
Definition: cs_user_solver.c:96
integer, save lndfac
Definition: mesh.f90:51
void initi1(const cs_int_t *iverif)
void pstusn(const cs_int_t *ntmabs, const cs_int_t *ntcabs, const cs_real_t *ttcabs)
Definition: cs_post.c:2322
int cs_user_solver_set(void)
Definition: cs_user_solver.c:85
integer, dimension(:), pointer icelbr
Definition: mesh.f90:84
double precision, save a
Definition: cs_fuel_incl.f90:143
void memfin(void)
Definition: memfin.f90:24
integer, save nvisbr
Definition: lagdim.f90:63
Definition: cs_mesh_quantities.h:51
void usvpst(const cs_int_t *nummai, const cs_int_t *nvar, const cs_int_t *nscal, const cs_int_t *nvlsta, const cs_int_t *ncelps, const cs_int_t *nfacps, const cs_int_t *nfbrps, const cs_int_t itypps[3], const cs_int_t lstcel[], const cs_int_t lstfac[], const cs_int_t lstfbr[], const cs_real_t dt[], const cs_real_t rtpa[], const cs_real_t rtp[], const cs_real_t propce[], const cs_real_t propfa[], const cs_real_t propfb[], const cs_real_t coefa[], const cs_real_t coefb[], const cs_real_t statce[], cs_real_t tracel[], cs_real_t trafac[], cs_real_t trafbr[])
double precision, dimension(:,:), pointer dijpf
Definition: mesh.f90:114
integer, save nthrdb
Definition: parall.f90:47
double precision, dimension(:,:), pointer dofij
Definition: mesh.f90:116
integer, dimension(:), pointer ifabor
Definition: mesh.f90:62
void normalen(const cs_int_t *n, const cs_real_t *x)
integer, dimension(:), pointer ifmfbr
Definition: mesh.f90:78
void cs_user_mesh_input(void)
Definition: cs_user_mesh.c:101
void cs_user_mesh_bad_cells_tag(cs_mesh_t *mesh, cs_mesh_quantities_t *mesh_quantities)
Definition: cs_user_mesh.c:627
integer(kind=8), save nsomgb
Definition: parall.f90:58
void cs_user_postprocess_writers(void)
Definition: cs_user_postprocess.c:264
integer, save ncelet
Definition: mesh.f90:42
void findpt(const cs_int_t *ncelet, const cs_int_t *ncel, const cs_real_t *xyzcen, const cs_real_t *xx, const cs_real_t *yy, const cs_real_t *zz, cs_int_t *node, cs_int_t *ndrang)
double precision, save voltot
Definition: cstphy.f90:141
integer, save ntmabs
Definition: optcal.f90:277
double precision, dimension(:), pointer surfan
Definition: mesh.f90:109
double precision, dimension(:,:), pointer surfac
Definition: mesh.f90:103
void csinit(const cs_int_t *irgpar, const cs_int_t *nrgpar, const cs_int_t *nthpar)
void dvvpst(const cs_int_t *nummai, const cs_int_t *numtyp, const cs_int_t *nvar, const cs_int_t *nscal, const cs_int_t *nvlsta, const cs_int_t *nvisbr, const cs_int_t *ncelps, const cs_int_t *nfacps, const cs_int_t *nfbrps, const cs_int_t itypps[3], const cs_int_t lstcel[], const cs_int_t lstfac[], const cs_int_t lstfbr[], const cs_real_t dt[], const cs_real_t rtpa[], const cs_real_t rtp[], const cs_real_t propce[], const cs_real_t propfa[], const cs_real_t propfb[], const cs_real_t coefa[], const cs_real_t coefb[], const cs_real_t statce[], const cs_real_t stativ[], const cs_real_t statfb[], cs_real_t tracel[], cs_real_t trafac[], cs_real_t trafbr[])
integer, save ncel
Definition: mesh.f90:42
#define END_C_DECLS
Definition: cs_defs.h:366
double precision, save ttcabs
Definition: optcal.f90:278
double cs_real_t
Definition: cs_defs.h:264
integer(kind=8), save nfacgb
Definition: parall.f90:58
void pergra(const cs_int_t *ivar, const cs_int_t *idimtr, const cs_int_t *irpvar)
double precision, save volmin
Definition: cstphy.f90:141
#define CS_PROCF(x, y)
Definition: cs_defs.h:379
integer(kind=8), save nfbrgb
Definition: parall.f90:58
void cs_user_syrthes_coupling(void)
Definition: cs_user_coupling.c:120
integer, dimension(:,:), pointer iprfml
Definition: mesh.f90:80
void csclli(void)
Definition: csclli.f90:24
integer, dimension(:), pointer ipnfbr
Definition: mesh.f90:71
integer, dimension(:), pointer nodfac
Definition: mesh.f90:70
double precision, dimension(:), pointer distb
Definition: mesh.f90:112
integer, save nthrdi
Definition: parall.f90:47
void cs_user_saturne_coupling(void)
Definition: cs_user_coupling.c:208
void cs_user_mesh_smoothe(cs_mesh_t *mesh)
Definition: cs_user_mesh.c:570
Definition: mesh.f90:25
integer, save lndfbr
Definition: mesh.f90:51
integer, save ncelbr
Definition: mesh.f90:46
integer, save nvar
Definition: dimens.f90:33
void majgeo(const cs_int_t *ncel, const cs_int_t *ncelet, const cs_int_t *nfac, const cs_int_t *nfabor, const cs_int_t *nsom, const cs_int_t *lndfac, const cs_int_t *lndfbr, const cs_int_t *ncelbr, const cs_int_t *ncelgb, const cs_int_t *nfacgb, const cs_int_t *nfbrgb, const cs_int_t *nsomgb, const cs_int_t *nthrdi, const cs_int_t *nthrdb, const cs_int_t *ngrpi, const cs_int_t *ngrpb, const cs_int_t *idxfi, const cs_int_t *idxfb, const cs_int_t ifacel[], const cs_int_t ifabor[], const cs_int_t ifmfbr[], const cs_int_t ifmcel[], const cs_int_t iprfml[], const cs_int_t ipnfac[], const cs_int_t nodfac[], const cs_int_t ipnfbr[], const cs_int_t nodfbr[], const cs_int_t icelbr[], const cs_real_t *volmin, const cs_real_t *volmax, const cs_real_t *voltot, const cs_real_t xyzcen[], const cs_real_t surfac[], const cs_real_t surfbo[], const cs_real_t cdgfac[], const cs_real_t cdgfbo[], const cs_real_t xyznod[], const cs_real_t volume[], const cs_real_t surfan[], const cs_real_t surfbn[], const cs_real_t dist[], const cs_real_t distb[], const cs_real_t pond[], const cs_real_t dijpf[], const cs_real_t diipb[], const cs_real_t dofij[])
integer, dimension(:), pointer ipnfac
Definition: mesh.f90:69