1 /*============================================================================
2  * Writer helper for time-varying 1d data
3  *============================================================================*/
4 
5 /*
6   This file is part of Code_Saturne, a general-purpose CFD tool.
7 
8   Copyright (C) 1998-2021 EDF S.A.
9 
10   This program is free software; you can redistribute it and/or modify it under
11   the terms of the GNU General Public License as published by the Free Software
12   Foundation; either version 2 of the License, or (at your option) any later
13   version.
14 
15   This program is distributed in the hope that it will be useful, but WITHOUT
16   ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
17   FOR A PARTICULAR PURPOSE.  See the GNU General Public License for more
18   details.
19 
20   You should have received a copy of the GNU General Public License along with
21   this program; if not, write to the Free Software Foundation, Inc., 51 Franklin
22   Street, Fifth Floor, Boston, MA 02110-1301, USA.
23 */
24 
25 /*----------------------------------------------------------------------------*/
26 
27 #include "cs_defs.h"
28 
29 /*----------------------------------------------------------------------------
30  * Standard C library headers
31  *----------------------------------------------------------------------------*/
32 
33 #include <ctype.h>
34 #include <errno.h>
35 #include <stdio.h>
36 #include <stdlib.h>
37 #include <string.h>
38 #include <assert.h>
39 
40 /*----------------------------------------------------------------------------
41  * Local headers
42  *----------------------------------------------------------------------------*/
43 
44 #include "bft_mem.h"
45 #include "bft_error.h"
46 
47 #include "cs_base.h"
48 #include "cs_timer.h"
49 #include "cs_time_step.h"
50 
51 /*----------------------------------------------------------------------------
52  *  Header for the current file
53  *----------------------------------------------------------------------------*/
54 
55 #include "cs_time_plot.h"
56 
57 /*----------------------------------------------------------------------------*/
58 
59 BEGIN_C_DECLS
60 
61 /*! \cond DOXYGEN_SHOULD_SKIP_THIS */
62 
63 /*=============================================================================
64  * Local Macro Definitions
65  *============================================================================*/
66 
67 /*=============================================================================
68  * Local Structure Definitions
69  *============================================================================*/
70 
71 /* Basic file output structure */
72 /*-----------------------------*/
73 
74 /* This structure should be non-null only for ranks
75    on which it is active */
76 
77 struct _cs_time_plot_t {
78 
79   char       *plot_name;        /* Associated plot name */
80   char       *file_name;        /* Associated file name */
81 
82   FILE       *f;                /* Associated file */
83 
84   cs_time_plot_format_t   format;  /* Associated format */
85 
86   bool        use_iteration;    /* Use iteration number instead of
87                                  * physical time ? */
88   bool        write_time_value; /* Output time value ? */
89 
90   double      flush_times[2];   /* 0: File flush time interval
91                                  *    (if < 0, no forced flush)
92                                  * 1: last write time
93                                  *    (-2 before first write) */
94   double      buffer_steps[2];  /* 0: Maximum number of time steps in
95                                  *    output buffer (if > 0, file is
96                                  *    not kept open)
97                                  * 1: current number of time steps in
98                                  *    output buffer */
99 
100   size_t      buffer_size;      /* Buffer size if required */
101   size_t      buffer_end;       /* Current buffer end */
102   char       *buffer;           /* Associated buffer if required */
103 
104   struct _cs_time_plot_t  *prev;  /* Previous in flush list */
105   struct _cs_time_plot_t  *next;  /* Next in flush list */
106 
107 };
108 
109 /*============================================================================
110  *  Global variables
111  *============================================================================*/
112 
113 static cs_time_plot_t   *_plots_head = NULL, *_plots_tail = NULL;
114 
115 static size_t            _n_files[2] = {0, 0};
116 static size_t            _n_files_max[2] = {0, 0};
117 static cs_time_plot_t  **_plot_files[2] = {NULL, NULL};
118 
119 static float             _flush_wtime_default = -1;
120 static int               _n_buffer_steps_default = -1;
121 
122 /*============================================================================
123  * Private function definitions
124  *============================================================================*/
125 
126 /*----------------------------------------------------------------------------
127  * Ensure buffer is large enough for upcoming data.
128  *
129  * The buffer is reallocated if necessary.
130  *
131  * parameters:
132  *   p        <-> time plot values file handler
133  *   min_size <-- minimum buffer size
134  *----------------------------------------------------------------------------*/
135 
136 static void
_ensure_buffer_size(cs_time_plot_t * p,size_t min_size)137 _ensure_buffer_size(cs_time_plot_t  *p,
138                     size_t           min_size)
139 {
140   /* Write data to line buffer */
141 
142   if (p->buffer_size < min_size) {
143     p->buffer_size = CS_MAX(1, p->buffer_size);
144     while (p->buffer_size < min_size)
145       p->buffer_size *= 2;
146     BFT_REALLOC(p->buffer, p->buffer_size, char);
147   }
148 }
149 
150 /*----------------------------------------------------------------------------
151  * Write file header for xmgrace/qsplotlib readable .dat files
152  *
153  * parameters:
154  *   p                <-> time plot values file handler
155  *   n_probes         <-- number of probes associated with this variable ?
156  *   probe_list       <-- numbers (1 to n) of probes if filtered, or NULL
157  *   probe_coords     <-- probe coordinates
158  *   probe_names      <-- probe names, or NULL
159  *----------------------------------------------------------------------------*/
160 
161 static void
_write_probe_header_dat(cs_time_plot_t * p,int n_probes,const int * probe_list,const cs_real_t probe_coords[],const char * probe_names[])162 _write_probe_header_dat(cs_time_plot_t    *p,
163                         int                n_probes,
164                         const int         *probe_list,
165                         const cs_real_t    probe_coords[],
166                         const char        *probe_names[])
167 {
168   int i, probe_id;
169   int col_id = 0;
170   FILE *_f = p->f;
171 
172   if (_f != NULL) {
173     fclose(_f);
174     p->f = NULL;
175   }
176 
177   _f = fopen(p->file_name, "w");
178   if (_f == NULL) {
179     bft_error(__FILE__, __LINE__, errno,
180               _("Error opening file: \"%s\""), p->file_name);
181     return;
182   }
183 
184   fprintf(_f, _("# Time varying values for: %s\n"
185                 "#\n"), p->plot_name);
186 
187   if (probe_coords != NULL) {
188     fprintf(_f, _("# Monitoring point coordinates:\n"));
189     for (i = 0; i < n_probes; i++) {
190       probe_id = i;
191       if (probe_list != NULL)
192         probe_id = probe_list[i] - 1;
193       if (probe_names != NULL)
194         fprintf(_f, "# %16s [%14.7e, %14.7e, %14.7e]\n",
195                 probe_names[i],
196                 probe_coords[probe_id*3],
197                 probe_coords[probe_id*3 + 1],
198                 probe_coords[probe_id*3 + 2]);
199       else
200         fprintf(_f, "#   %6i [%14.7e, %14.7e, %14.7e]\n",
201                 probe_id + 1,
202                 probe_coords[probe_id*3],
203                 probe_coords[probe_id*3 + 1],
204                 probe_coords[probe_id*3 + 2]);
205     }
206     fprintf(_f, "#\n");
207   }
208   else if (probe_names != NULL) {
209     fprintf(_f, _("# Monitoring points:\n"));
210     for (i = 0; i < n_probes; i++)
211       fprintf(_f, "# %s\n", probe_names[i]);
212     fprintf(_f, "#\n");
213   }
214 
215   fprintf(_f, _("# Columns:\n"));
216   if (p->use_iteration)
217     fprintf(_f, _("#   %d:     Time step number\n"), col_id++);
218   else
219     fprintf(_f, _("#   %d:     Physical time\n"), col_id++);
220   fprintf(_f, _("#   %d - :  Values at monitoring points\n"), col_id);
221 
222   fprintf(_f,
223           "#\n"
224           "#TITLE: %s\n"
225           "#COLUMN_TITLES: ", p->plot_name);
226   if (p->use_iteration)
227     fprintf(_f, " nt");
228   else
229     fprintf(_f, " t");
230   for (i = 0; i < n_probes; i++) {
231     probe_id = i;
232     if (probe_list != NULL)
233       probe_id = probe_list[i] - 1;
234     if (probe_names != NULL) {
235       if (probe_coords != NULL)
236         fprintf(_f, " | %s [%9.5e, %9.5e, %9.5e]",
237                 probe_names[i],
238                 probe_coords[probe_id*3],
239                 probe_coords[probe_id*3 + 1],
240                 probe_coords[probe_id*3 + 2]);
241       else
242         fprintf(_f, " | %s", probe_names[i]);
243     }
244     else {
245       if (probe_coords != NULL)
246         fprintf(_f, " | %d [%9.5e, %9.5e, %9.5e]",
247                 probe_id + 1,
248                 probe_coords[probe_id*3],
249                 probe_coords[probe_id*3 + 1],
250                 probe_coords[probe_id*3 + 2]);
251       else
252         fprintf(_f, " | %d", probe_id + 1);
253     }
254   }
255   fprintf(_f, "\n");
256 
257   fprintf(_f, "#COLUMN_UNITS: ");
258   if (p->use_iteration)
259     fprintf(_f, " iter");
260   else
261     fprintf(_f, " s");
262   for (probe_id = 0; probe_id < n_probes; probe_id++)
263     fprintf(_f, " -");
264   fprintf(_f, "\n#\n");
265 
266   /* Close file or assign it to handler depending on options */
267 
268   if (p->buffer_steps[0] > 0) {
269     if (fclose(_f) != 0)
270       bft_error(__FILE__, __LINE__, errno,
271                 _("Error closing file: \"%s\""), p->file_name);
272   }
273   else
274     p->f = _f;
275 }
276 
277 /*----------------------------------------------------------------------------
278  * Write probe coordinates for CSV files
279  *
280  * parameters:
281  *   file_prefix      <-- file name prefix
282  *   plot_name        <-- plot name
283  *   n_probes         <-- number of probes associated with this plot
284  *   probe_list       <-- numbers (1 to n) of probes if filtered, or NULL
285  *   probe_coords     <-- probe coordinates
286  *----------------------------------------------------------------------------*/
287 
288 static void
_write_probe_coords_csv(const char * file_prefix,const char * plot_name,int n_probes,const int * probe_list,const cs_real_t probe_coords[])289 _write_probe_coords_csv(const char        *file_prefix,
290                         const char        *plot_name,
291                         int                n_probes,
292                         const int         *probe_list,
293                         const cs_real_t    probe_coords[])
294 {
295   int i, probe_id;
296   char *file_name;
297   FILE *_f;
298 
299   BFT_MALLOC(file_name,
300              strlen(file_prefix) + strlen(plot_name) + strlen("_coords") + 4 + 1,
301              char);
302 
303   if (probe_coords != NULL) {
304     sprintf(file_name, "%s%s%s.csv", file_prefix, plot_name, "_coords");
305 
306     _f = fopen(file_name, "w");
307     if (_f == NULL) {
308       bft_error(__FILE__, __LINE__, errno,
309           _("Error opening file: \"%s\""), file_name);
310       return;
311     }
312 
313     fprintf(_f, "x, y, z\n");
314     for (i = 0; i < n_probes; i++) {
315       probe_id = i;
316       if (probe_list != NULL)
317         probe_id = probe_list[i] - 1;
318       fprintf(_f, "%14.7e, %14.7e, %14.7e\n",
319           probe_coords[probe_id*3],
320           probe_coords[probe_id*3 + 1],
321           probe_coords[probe_id*3 + 2]);
322     }
323 
324     /* Close file or assign it ot handler depending on options */
325 
326     if (fclose(_f) != 0)
327       bft_error(__FILE__, __LINE__, errno,
328           _("Error closing file: \"%s\""), file_name);
329 
330   }
331 
332   BFT_FREE(file_name);
333 }
334 
335 /*----------------------------------------------------------------------------
336  * Write file header for CSV files
337  *
338  * parameters:
339  *   p                <-> time plot values file handler
340  *   n_probes         <-- number of probes associated with this variable ?
341  *   probe_list       <-- numbers (1 to n) of probes if filtered, or NULL
342  *   probe_coords     <-- probe coordinates
343  *   probe_names      <-- probe names, or NULL
344  *----------------------------------------------------------------------------*/
345 
346 static void
_write_probe_header_csv(cs_time_plot_t * p,int n_probes,const int * probe_list,const cs_real_t probe_coords[],const char * probe_names[])347 _write_probe_header_csv(cs_time_plot_t    *p,
348                         int                n_probes,
349                         const int         *probe_list,
350                         const cs_real_t    probe_coords[],
351                         const char        *probe_names[])
352 {
353   int i, probe_id;
354   FILE *_f = p->f;
355 
356   if (_f != NULL) {
357     fclose(_f);
358     p->f = NULL;
359   }
360 
361   _f = fopen(p->file_name, "w");
362   if (_f == NULL) {
363     bft_error(__FILE__, __LINE__, errno,
364               _("Error opening file: \"%s\""), p->file_name);
365     return;
366   }
367 
368   if (p->use_iteration)
369     fprintf(_f, " iteration");
370   else
371     fprintf(_f, "t");
372   for (i = 0; i < n_probes; i++) {
373     probe_id = i;
374     if (probe_list != NULL)
375       probe_id = probe_list[i] - 1;
376     if (probe_coords != NULL) {
377       if (probe_names != NULL)
378         fprintf(_f, ", %s [%9.5e; %9.5e; %9.5e]",
379                 probe_names[i],
380                 probe_coords[probe_id*3],
381                 probe_coords[probe_id*3 + 1],
382                 probe_coords[probe_id*3 + 2]);
383       else
384         fprintf(_f, ", %d [%9.5e; %9.5e; %9.5e]",
385                 probe_id + 1,
386                 probe_coords[probe_id*3],
387                 probe_coords[probe_id*3 + 1],
388                 probe_coords[probe_id*3 + 2]);
389     }
390     else if (probe_names != NULL)
391       fprintf(_f, ", %s", probe_names[i]);
392     else
393       fprintf(_f, ", %d", probe_id + 1);
394   }
395   fprintf(_f, "\n");
396 
397   /* Close file or assign it to handler depending on options */
398 
399   if (p->buffer_steps[0] > 0) {
400     if (fclose(_f) != 0)
401       bft_error(__FILE__, __LINE__, errno,
402                 _("Error closing file: \"%s\""), p->file_name);
403   }
404   else
405     p->f = _f;
406 }
407 
408 /*----------------------------------------------------------------------------
409  * Write file header for xmgrace/qsplotlib readable .dat files
410  *
411  * parameters:
412  *   p                  <-> time plot values file handler
413  *   n_structures       <-- number of structures associated with this plot
414  *   mass_matrixes      <-- mass matrix coefficients (3x3 blocks)
415  *   damping_matrixes   <-- damping matrix coefficients (3x3 blocks)
416  *   stiffness_matrixes <-- stiffness matrix coefficients (3x3 blocks)
417  *----------------------------------------------------------------------------*/
418 
419 static void
_write_struct_header_dat(cs_time_plot_t * p,int n_structures,const cs_real_t mass_matrixes[],const cs_real_t damping_matrixes[],const cs_real_t stiffness_matrixes[])420 _write_struct_header_dat(cs_time_plot_t    *p,
421                          int                n_structures,
422                          const cs_real_t    mass_matrixes[],
423                          const cs_real_t    damping_matrixes[],
424                          const cs_real_t    stiffness_matrixes[])
425 {
426   int i, struct_id;
427   int col_id = 0;
428   FILE *_f = p->f;
429 
430   const int perm_id[9] = {0, 3, 6, 1, 4, 7, 2, 5, 8};
431 
432   if (_f != NULL) {
433     fclose(_f);
434     p->f = NULL;
435   }
436 
437   _f = fopen(p->file_name, "w");
438   if (_f == NULL) {
439     bft_error(__FILE__, __LINE__, errno,
440               _("Error opening file: \"%s\""), p->file_name);
441     return;
442   }
443 
444   fprintf(_f, _("# Time varying values for: %s\n"
445                 "#\n"), p->plot_name);
446 
447   fprintf(_f, _("# Number of structures: %d\n"
448                 "#\n"), n_structures);
449 
450   for (struct_id = 0; struct_id < n_structures; struct_id++) {
451     double m_tmp[9], d_tmp[9], s_tmp[9];
452     for (i = 0; i < 9; i++) {
453       m_tmp[i] = mass_matrixes[perm_id[i] + struct_id*9];
454       d_tmp[i] = damping_matrixes[perm_id[i] + struct_id*9];
455       s_tmp[i] = stiffness_matrixes[perm_id[i] + struct_id*9];
456     }
457     fprintf(_f, _("# Structure: %i\n"
458                   "#\n"), struct_id + 1);
459     fprintf(_f, _("# Mass:       [%14.7e, %14.7e, %14.7e]\n"
460                   "#             [%14.7e, %14.7e, %14.7e]\n"
461                   "#             [%14.7e, %14.7e, %14.7e]\n\n"),
462             m_tmp[0], m_tmp[1], m_tmp[2],
463             m_tmp[3], m_tmp[4], m_tmp[5],
464             m_tmp[6], m_tmp[7], m_tmp[8]);
465     fprintf(_f, _("# Damping:    [%14.7e, %14.7e, %14.7e]\n"
466                   "#             [%14.7e, %14.7e, %14.7e]\n"
467                   "#             [%14.7e, %14.7e, %14.7e]\n\n"),
468             d_tmp[0], d_tmp[1], d_tmp[2],
469             d_tmp[3], d_tmp[4], d_tmp[5],
470             d_tmp[6], d_tmp[7], d_tmp[8]);
471     fprintf(_f, _("# Stiffness:  [%14.7e, %14.7e, %14.7e]\n"
472                   "#             [%14.7e, %14.7e, %14.7e]\n"
473                   "#             [%14.7e, %14.7e, %14.7e]\n\n"),
474             s_tmp[0], s_tmp[1], s_tmp[2],
475             s_tmp[3], s_tmp[4], s_tmp[5],
476             s_tmp[6], s_tmp[7], s_tmp[8]);
477    }
478   fprintf(_f,
479           _("# (when structure characteristics are variable, the values\n"
480             "# above are those at the computation initialization.\n\n"));
481 
482   fprintf(_f, _("# Columns:\n"));
483   if (p->use_iteration)
484     fprintf(_f, _("#   %d:     Time step number\n"), col_id++);
485   else
486     fprintf(_f, _("#   %d:     Physical time\n"), col_id++);
487   fprintf(_f, _("#   %d - :  Values for each structure\n"), col_id);
488 
489 
490   fprintf(_f,
491           "#\n"
492           "#TITLE: %s\n"
493           "#COLUMN_TITLES: ", p->plot_name);
494   if (p->use_iteration)
495     fprintf(_f, " nt");
496   else
497     fprintf(_f, " t");
498   for (i = 0; i < n_structures; i++)
499     fprintf(_f, " | %d", i + 1);
500   fprintf(_f, "\n");
501 
502   fprintf(_f, "#COLUMN_UNITS: ");
503   if (p->use_iteration)
504     fprintf(_f, " iter");
505   else
506     fprintf(_f, " s");
507   for (i = 0; i < n_structures; i++)
508     fprintf(_f, " -");
509   fprintf(_f, "\n#\n");
510 
511   /* Close file or assign it ot handler depending on options */
512 
513   if (p->buffer_steps[0] > 0) {
514     if (fclose(_f) != 0)
515       bft_error(__FILE__, __LINE__, errno,
516                 _("Error closing file: \"%s\""), p->file_name);
517   }
518   else
519     p->f = _f;
520 }
521 
522 /*----------------------------------------------------------------------------
523  * Write file header for CSV files
524  *
525  * parameters:
526  *   p                  <-> time plot values file handler
527  *   n_structures       <-- number of structures associated with this plot
528  *----------------------------------------------------------------------------*/
529 
530 static void
_write_struct_header_csv(cs_time_plot_t * p,int n_structures)531 _write_struct_header_csv(cs_time_plot_t  *p,
532                          int              n_structures)
533 {
534   int i;
535   FILE *_f = p->f;
536 
537   if (_f != NULL) {
538     fclose(_f);
539     p->f = NULL;
540   }
541 
542   _f = fopen(p->file_name, "w");
543   if (_f == NULL) {
544     bft_error(__FILE__, __LINE__, errno,
545               _("Error opening file: \"%s\""), p->file_name);
546     return;
547   }
548 
549   if (p->use_iteration)
550     fprintf(_f, " iteration");
551   else
552     fprintf(_f, "t");
553   for (i = 0; i < n_structures; i++) {
554     fprintf(_f, ",%d", i + 1);
555   }
556   fprintf(_f, "\n");
557 
558   /* Close file or assign it ot handler depending on options */
559 
560   if (p->buffer_steps[0] > 0) {
561     if (fclose(_f) != 0)
562       bft_error(__FILE__, __LINE__, errno,
563                 _("Error closing file: \"%s\""), p->file_name);
564   }
565   else
566     p->f = _f;
567 }
568 
569 /*----------------------------------------------------------------------------
570  * Add a time plot to the global time plots array.
571  *----------------------------------------------------------------------------*/
572 
573 static void
_time_plot_register(cs_time_plot_t * p)574 _time_plot_register(cs_time_plot_t  *p)
575 {
576   p->prev = _plots_tail;
577   p->next = NULL;
578 
579   if (_plots_head == NULL)
580     _plots_head = p;
581   else if (_plots_head->next == NULL)
582     _plots_head->next = p;
583 
584   if (_plots_tail != NULL)
585     _plots_tail->next = p;
586 
587   _plots_tail = p;
588 }
589 
590 /*----------------------------------------------------------------------------
591  * Remove a time plot from the global time plots array.
592  *----------------------------------------------------------------------------*/
593 
594 static void
_time_plot_unregister(cs_time_plot_t * p)595 _time_plot_unregister(cs_time_plot_t  *p)
596 {
597   if (_plots_head == p)
598     _plots_head = p->next;
599   if (_plots_tail == p)
600     _plots_tail = p->prev;
601 
602   if (p->prev != NULL)
603     p->prev->next = p->next;
604   if (p->next != NULL)
605     p->next->prev = p->prev;
606 }
607 
608 /*----------------------------------------------------------------------------
609  * Create time plot writer for a given variable
610  *
611  * parameters:
612  *   plot_name        <-- plot (variable) name
613  *   file_prefix      <-- file name prefix
614  *   format           <-- associated file format
615  *   use_iteration    <-- should we use the iteration number instead of the
616  *                        physical time ?
617  *   flush_wtime      <-- elapsed time interval between file flushes
618  *                        (if < 0, no forced flush)
619  *   n_buffer_steps   <-- number of time steps in output buffer if
620  *                        file is not to be kept open
621  *
622  * returns
623  *   pointer to corresponding probe writer
624  *----------------------------------------------------------------------------*/
625 
626 static cs_time_plot_t *
_plot_file_create(const char * plot_name,const char * file_prefix,cs_time_plot_format_t format,bool use_iteration,double flush_wtime,int n_buffer_steps)627 _plot_file_create(const char             *plot_name,
628                   const char             *file_prefix,
629                   cs_time_plot_format_t   format,
630                   bool                    use_iteration,
631                   double                  flush_wtime,
632                   int                     n_buffer_steps)
633 {
634   size_t i;
635 
636   cs_time_plot_t *p = NULL;
637 
638   BFT_MALLOC(p, 1, cs_time_plot_t);
639   BFT_MALLOC(p->plot_name, strlen(plot_name) + 1, char);
640   BFT_MALLOC(p->file_name,
641              strlen(file_prefix) + strlen(plot_name) + 4 + 1,
642              char);
643 
644   strcpy(p->plot_name, plot_name);
645   switch (format) {
646   case CS_TIME_PLOT_DAT:
647     sprintf(p->file_name, "%s%s.dat", file_prefix, plot_name);
648     break;
649   case CS_TIME_PLOT_CSV:
650     sprintf(p->file_name, "%s%s.csv", file_prefix, plot_name);
651     break;
652   default:
653     break;
654   }
655 
656   for (i = strlen(file_prefix); p->file_name[i] != '\0'; i++) {
657     if (isspace(p->file_name[i]))
658       p->file_name[i] = '_';
659   }
660 
661   p->f = NULL;
662 
663   p->format = format;
664   p->use_iteration = use_iteration;
665 
666   p->flush_times[0] = flush_wtime;
667   p->flush_times[1] = -2;
668 
669   p->buffer_steps[0] = n_buffer_steps;
670   p->buffer_steps[1] = 0;
671 
672   p->buffer_size = 256;
673   p->buffer_end = 0;
674 
675   BFT_MALLOC(p->buffer, p->buffer_size, char);
676 
677   _time_plot_register(p);
678 
679   return p;
680 }
681 
682 /*----------------------------------------------------------------------------
683  * Write buffered values to file if applicable
684  *
685  * parameters:
686  *   p <-> time plot values file handler
687  *----------------------------------------------------------------------------*/
688 
689 static void
_plot_file_check_or_write(cs_time_plot_t * p)690 _plot_file_check_or_write(cs_time_plot_t  *p)
691 {
692   size_t n_written;
693 
694   /* Return immediately if we are buffering and not writing now */
695 
696   if (   p->buffer_steps[0] > 0
697       && p->buffer_steps[1] < p->buffer_steps[0]) {
698     p->buffer_steps[1] += 1;
699       return;
700   }
701 
702   /* Ensure file is open */
703 
704   if (p->f == NULL) {
705     p->f = fopen(p->file_name, "a");
706     if (p->f == NULL) {
707       bft_error(__FILE__, __LINE__, errno,
708                 _("Error re-opening file: \"%s\""), p->file_name);
709       p->buffer_end = 0;
710       return;
711     }
712   }
713 
714   /* Write buffer contents */
715   n_written = fwrite(p->buffer, 1, p->buffer_end, p->f);
716 
717   if (n_written < p->buffer_end)
718     bft_error(__FILE__, __LINE__, ferror(p->f),
719               _("Error writing file: \"%s\""), p->file_name);
720 
721   p->buffer_end = 0;
722 
723   /* Close or flush file depending on options */
724 
725   if (p->buffer_steps[0] > 0) {
726 
727     assert(p->buffer_steps[1] >= p->buffer_steps[0]);
728     if (fclose(p->f) != 0)
729       bft_error(__FILE__, __LINE__, errno,
730                   _("Error closing file: \"%s\""), p->file_name);
731     p->f = NULL;
732     p->buffer_steps[1] = 0;
733 
734   }
735   else {
736     double cur_time = cs_timer_wtime();
737     if (   p->flush_times[0] > 0
738         && (cur_time - p->flush_times[1]) > p->flush_times[0]) {
739       p->flush_times[1] = cur_time;
740       fflush(p->f);
741     }
742   }
743 }
744 
745 /*----------------------------------------------------------------------------
746  * Ensure the array of time plots for the Fortran API is large enough.
747  *
748  * parameters:
749  *   plot_num  <-- associated time plot number (1 to n)
750  *   plot_name <-- plot name
751  *   format    <-- associated file format
752  *----------------------------------------------------------------------------*/
753 
754 static void
_fortran_time_plot_realloc(int plot_num,const char * plot_name,cs_time_plot_format_t format)755 _fortran_time_plot_realloc(int                     plot_num,
756                            const char             *plot_name,
757                            cs_time_plot_format_t   format)
758 
759 {
760   if (plot_num <= 0)
761     bft_error(__FILE__, __LINE__, errno,
762               _("Plot number for \"%s\" must be > 0 and not %d."),
763               plot_name, plot_num);
764 
765   /* Ensure writer helper array is large enough */
766 
767   if (plot_num >= (int)_n_files_max[format]) {
768     int i;
769     int _n_vars_new = 1;
770     while (_n_vars_new < plot_num)
771       _n_vars_new *= 2;
772     BFT_REALLOC(_plot_files[format], _n_vars_new, cs_time_plot_t *);
773     for (i = _n_files_max[format]; i < _n_vars_new; i++)
774       _plot_files[format][i] = NULL;
775     _n_files_max[format] = _n_vars_new;
776   }
777   else if (_plot_files[format][plot_num - 1] != NULL)
778     bft_error(__FILE__, __LINE__, errno,
779               _("Plot number %d is already defined."), plot_num);
780   _n_files[format] += 1;
781 }
782 
783 /*! (DOXYGEN_SHOULD_SKIP_THIS) \endcond */
784 
785 /*============================================================================
786  * Public function definitions for Fortran API
787  *============================================================================*/
788 
789 /*----------------------------------------------------------------------------
790  * Create a writer for time plot structure-type data.
791  *
792  * This subroutine should only be called by one rank for a given data series.
793  *
794  * subroutine tpsini (tplnum, tplnam, tplpre, tplfmt, idtvar,
795  * *****************
796  *                    nprb,   lstprb, xyzprb, lnam,   lpre)
797  *
798  * integer          tplnum      : <-- : number of plot to create (> 0)
799  * character        tplnam      : <-- : name of associated plot
800  * character        tplpre      : <-- : prefix for associated file
801  * integer          tplfmt      : <-- : associated format
802  *                                      (1: dat, 2: csv, 3: both)
803  * integer          idtvar      : <-- : calculation time dependency
804  * integer          nstru       : <-- : number of structures
805  * double precision xmstru      : <-- : mass matrixes
806  * double precision xcstru      : <-- : damping matrixes
807  * double precision xkstru      : <-- : stiffness matrixes
808  * integer          lnam        : <-- : name length
809  * integer          lpre        : <-- : prefix length
810  *----------------------------------------------------------------------------*/
811 
CS_PROCF(tpsini,TPSINI)812 void CS_PROCF (tpsini, TPSINI)
813 (
814  const int       *tplnum,
815  const char      *tplnam,
816  const char      *tplpre,
817  const int       *tplfmt,
818  const int       *idtvar,
819  const int       *nstru,
820  const cs_real_t *xmstru,
821  const cs_real_t *xcstru,
822  const cs_real_t *xkstru,
823  const int       *lnam,
824  const int       *lpre
825  CS_ARGF_SUPP_CHAINE              /*     (possible 'length' arguments added
826                                          by many Fortran compilers) */
827 )
828 {
829   /* Copy Fortran strings to C strings */
830 
831   cs_time_plot_format_t fmt;
832   char *plot_name   = cs_base_string_f_to_c_create(tplnam, *lnam);
833   char *file_prefix = cs_base_string_f_to_c_create(tplpre, *lpre);
834   bool use_iteration = false;
835 
836   if (*idtvar == CS_TIME_STEP_STEADY || *idtvar == CS_TIME_STEP_LOCAL)
837     use_iteration = true;
838 
839   /* Main processing */
840 
841   for (fmt = CS_TIME_PLOT_DAT; fmt <= CS_TIME_PLOT_CSV; fmt++) {
842 
843     int fmt_mask = fmt + 1;
844 
845     if (*tplfmt & fmt_mask) {
846 
847       _fortran_time_plot_realloc(*tplnum, plot_name, fmt);
848 
849       _plot_files[fmt][*tplnum - 1]
850         = cs_time_plot_init_struct(plot_name,
851                                    file_prefix,
852                                    fmt,
853                                    use_iteration,
854                                    _flush_wtime_default,
855                                    _n_buffer_steps_default,
856                                    *nstru,
857                                    xmstru,
858                                    xcstru,
859                                    xkstru);
860     }
861 
862   }
863 
864   /* Free temporary C strings */
865 
866   cs_base_string_f_to_c_free(&plot_name);
867   cs_base_string_f_to_c_free(&file_prefix);
868 }
869 
870 /*----------------------------------------------------------------------------
871  * Finalize a writer for time plot data.
872  *
873  * This subroutine should only be called by one rank for a given data series.
874  *
875  * subroutine tplend (tplnum, tplfmt)
876  * *****************
877  *
878  * integer          tplnum      : <-- : number of plot to finalize (> 0)
879  * integer          tplfmt      : <-- : associated format
880  *                                      (1: dat, 2: csv, 3: both)
881  *----------------------------------------------------------------------------*/
882 
CS_PROCF(tplend,TPLEND)883 void CS_PROCF (tplend, TPLEND)
884 (
885  const int       *tplnum,
886  const int       *tplfmt
887 )
888 {
889   cs_time_plot_format_t fmt;
890   cs_time_plot_t *p = NULL;
891 
892   for (fmt = CS_TIME_PLOT_DAT; fmt <= CS_TIME_PLOT_CSV; fmt++) {
893 
894     int fmt_mask = fmt + 1;
895 
896     if (*tplfmt & fmt_mask) {
897 
898       if (*tplnum <= 0 || *tplnum > (int)_n_files_max[fmt])
899         bft_error(__FILE__, __LINE__, 0,
900                   _("Plot number must be in the interval [1, %d] and not %d."),
901                   (int)(_n_files_max[fmt]), *tplnum);
902 
903       p = _plot_files[fmt][*tplnum - 1];
904 
905       if (p != NULL) {
906         cs_time_plot_finalize(&p);
907         _plot_files[fmt][*tplnum - 1] = NULL;
908         _n_files[fmt] -= 1;
909         if (_n_files[fmt] == 0) {
910           _n_files_max[fmt] = 0;
911           BFT_FREE(_plot_files[fmt]);
912         }
913       }
914     }
915 
916   }
917 
918 }
919 
920 /*----------------------------------------------------------------------------
921  * Write time plot values.
922  *
923  * subroutine tppwri (tplnum, tplfmt, nprb, ntcabs, ttcabs, valprb)
924  * *****************
925  *
926  * integer          tplnum      : <-- : number of associated plot (> 0)
927  * integer          tplfmt      : <-- : associated format
928  *                                      (1: dat, 2: csv, 3: both)
929  * integer          nprb        : <-- : number of probes
930  * integer          ntcabs      : <-- : current time step number
931  * double precision ttcabs      : <-- : current time value
932  * double precision valprb      : <-- : probe values
933  *----------------------------------------------------------------------------*/
934 
CS_PROCF(tplwri,TPLWRI)935 void CS_PROCF (tplwri, TPLWRI)
936 (
937  const int       *tplnum,
938  const int       *tplfmt,
939  const int       *nprb,
940  const int       *ntcabs,
941  const cs_real_t *ttcabs,
942  const cs_real_t *valprb
943 )
944 {
945   cs_time_plot_format_t  fmt;
946 
947   /* Main processing */
948 
949   for (fmt = CS_TIME_PLOT_DAT; fmt <= CS_TIME_PLOT_CSV; fmt++) {
950 
951     int fmt_mask = fmt + 1;
952 
953     if (*tplfmt & fmt_mask) {
954 
955       if (*tplnum > -1 && (size_t)(*tplnum - 1) < _n_files_max[fmt]) {
956         cs_time_plot_t *p = _plot_files[fmt][*tplnum - 1];
957         cs_time_plot_vals_write(p, *ntcabs, *ttcabs, *nprb, valprb);
958       }
959 
960     }
961 
962   }
963 }
964 
965 /*----------------------------------------------------------------------------
966  * Return the number of time plots accessible through the Fortran API
967  *
968  * This subroutine will only return the number of time plots defined by the
969  * local rank
970  *
971  * subroutine tplnbr (ntpl)
972  * *****************
973  *
974  * integer          ntpl        : --> : number of time plots defined
975  *----------------------------------------------------------------------------*/
976 
CS_PROCF(tplnbr,TPLNBR)977 void CS_PROCF (tplnbr, TPLNBR)
978 (
979  int       *ntpl
980 )
981 {
982   int fmt;
983 
984   *ntpl = 0;
985 
986   for (fmt = CS_TIME_PLOT_DAT; fmt <= CS_TIME_PLOT_CSV; fmt++) {
987     if (_n_files_max[fmt] > (size_t)(*ntpl))
988       *ntpl = _n_files_max[fmt];
989   }
990 }
991 
992 /*============================================================================
993  * Public function definitions
994  *============================================================================*/
995 
996 /*----------------------------------------------------------------------------
997  * Initialize a plot file writer for probe-type plots
998  *
999  * This function should only be called by one rank for a given data series.
1000  *
1001  * parameters:
1002  *   plot_name        <-- plot (variable) name
1003  *   file_prefix      <-- file name prefix
1004  *   format           <-- associated file format
1005  *   use_iteration    <-- should we use the iteration number instead of the
1006  *                        physical time ?
1007  *   flush_wtime      <-- elapsed time interval between file flushes
1008  *                        (if < 0, no forced flush)
1009  *   n_buffer_steps   <-- number of time steps in output buffer if
1010  *                        file is not to be kept open
1011  *   n_probes         <-- number of probes associated with this plot
1012  *   probe_list       <-- numbers (1 to n) of probes if filtered, or NULL
1013  *   probe_coords     <-- probe coordinates, or NULL
1014  *   probe_names      <-- probe names, or NULL
1015  *
1016  * returns:
1017  *   pointer to new time plot writer
1018  *----------------------------------------------------------------------------*/
1019 
1020 cs_time_plot_t *
cs_time_plot_init_probe(const char * plot_name,const char * file_prefix,cs_time_plot_format_t format,bool use_iteration,double flush_wtime,int n_buffer_steps,int n_probes,const int * probe_list,const cs_real_t probe_coords[],const char * probe_names[])1021 cs_time_plot_init_probe(const char             *plot_name,
1022                         const char             *file_prefix,
1023                         cs_time_plot_format_t   format,
1024                         bool                    use_iteration,
1025                         double                  flush_wtime,
1026                         int                     n_buffer_steps,
1027                         int                     n_probes,
1028                         const int              *probe_list,
1029                         const cs_real_t         probe_coords[],
1030                         const char             *probe_names[])
1031 {
1032   cs_time_plot_t  *p = _plot_file_create(plot_name,
1033                                          file_prefix,
1034                                          format,
1035                                          use_iteration,
1036                                          flush_wtime,
1037                                          n_buffer_steps);
1038 
1039   switch (format) {
1040   case CS_TIME_PLOT_DAT:
1041     _write_probe_header_dat(p, n_probes, probe_list, probe_coords, probe_names);
1042     break;
1043   case CS_TIME_PLOT_CSV:
1044     _write_probe_coords_csv(file_prefix,
1045                             plot_name,
1046                             n_probes,
1047                             probe_list,
1048                             probe_coords);
1049     _write_probe_header_csv(p, n_probes, probe_list, probe_coords, probe_names);
1050     break;
1051   default:
1052     break;
1053   }
1054 
1055   cs_time_plot_flush(p); /* Flush after creation so plot tools can start
1056                             analyzing file immediately, while computation
1057                             is running */
1058 
1059   return p;
1060 }
1061 
1062 /*----------------------------------------------------------------------------
1063  * Initialize a plot file writer for structure-type plots
1064  *
1065  * This function should only be called by one rank for a given data series.
1066  *
1067  * parameters:
1068  *   plot_name          <-- plot (variable) name
1069  *   file_prefix        <-- file name prefix
1070  *   format             <-- associated file format
1071  *   use_iteration      <-- should we use the iteration number instead of the
1072  *                          physical time ?
1073  *   flush_wtime        <-- elapsed time interval between file flushes
1074  *                          (if < 0, no forced flush)
1075  *   n_buffer_steps     <-- number of time steps in output buffer if
1076  *                          file is not to be kept open
1077  *   n_structures       <-- number of structures associated with this plot
1078  *   mass_matrixes      <-- mass matrix coefficients (3x3 blocks)
1079  *   damping_matrixes   <-- damping matrix coefficients (3x3 blocks)
1080  *   stiffness_matrixes <-- stiffness matrix coefficients (3x3 blocks)
1081  *
1082  * returns:
1083  *   pointer to new time plot writer
1084  *----------------------------------------------------------------------------*/
1085 
1086 cs_time_plot_t *
cs_time_plot_init_struct(const char * plot_name,const char * file_prefix,cs_time_plot_format_t format,bool use_iteration,double flush_wtime,int n_buffer_steps,int n_structures,const cs_real_t mass_matrixes[],const cs_real_t damping_matrixes[],const cs_real_t stiffness_matrixes[])1087 cs_time_plot_init_struct(const char             *plot_name,
1088                          const char             *file_prefix,
1089                          cs_time_plot_format_t   format,
1090                          bool                    use_iteration,
1091                          double                  flush_wtime,
1092                          int                     n_buffer_steps,
1093                          int                     n_structures,
1094                          const cs_real_t         mass_matrixes[],
1095                          const cs_real_t         damping_matrixes[],
1096                          const cs_real_t         stiffness_matrixes[])
1097 {
1098   cs_time_plot_t  *p = _plot_file_create(plot_name,
1099                                          file_prefix,
1100                                          format,
1101                                          use_iteration,
1102                                          flush_wtime,
1103                                          n_buffer_steps);
1104 
1105   switch (format) {
1106   case CS_TIME_PLOT_DAT:
1107     _write_struct_header_dat(p, n_structures,
1108                              mass_matrixes, damping_matrixes, stiffness_matrixes);
1109     break;
1110   case CS_TIME_PLOT_CSV:
1111     _write_struct_header_csv(p, n_structures);
1112   break;
1113   default:
1114     break;
1115   }
1116 
1117   return p;
1118 }
1119 
1120 /*----------------------------------------------------------------------------
1121  * Finalize time plot writer for a given variable
1122  *
1123  * This function should only be called by one rank for a given data series.
1124  *
1125  * parameters:
1126  *   p <-> time plot values file handler
1127  *----------------------------------------------------------------------------*/
1128 
1129 void
cs_time_plot_finalize(cs_time_plot_t ** p)1130 cs_time_plot_finalize(cs_time_plot_t  **p)
1131 {
1132   if (p != NULL) {
1133 
1134     cs_time_plot_t *_p = *p;
1135 
1136     _time_plot_unregister(_p);
1137 
1138     if (_p->buffer_steps[0] > 0)
1139       _p->buffer_steps[1] = _p->buffer_steps[0] + 1;
1140 
1141     _plot_file_check_or_write(_p);
1142 
1143     if (_p->f != NULL) {
1144       if (fclose(_p->f) != 0)
1145         bft_error(__FILE__, __LINE__, errno,
1146                   _("Error closing file: \"%s\""), _p->file_name);
1147     }
1148 
1149     BFT_FREE(_p->buffer);
1150     BFT_FREE(_p->file_name);
1151     BFT_FREE(_p->plot_name);
1152 
1153     BFT_FREE(*p);
1154   }
1155 }
1156 
1157 /*----------------------------------------------------------------------------
1158  * Write time plot values
1159  *
1160  * This function should only be called by one rank for a given data series.
1161  *
1162  * parameters:
1163  *   p      <-- pointer to associated plot structure
1164  *   tn     <-- associated time step number
1165  *   t      <-- associated time value
1166  *   n_vals <-- number of associated time values
1167  *   vals   <-- associated time values
1168  *----------------------------------------------------------------------------*/
1169 
1170 void
cs_time_plot_vals_write(cs_time_plot_t * p,int tn,double t,int n_vals,const cs_real_t vals[])1171 cs_time_plot_vals_write(cs_time_plot_t  *p,
1172                         int              tn,
1173                         double           t,
1174                         int              n_vals,
1175                         const cs_real_t  vals[])
1176 {
1177   int i;
1178 
1179   if (p == NULL)
1180     return;
1181 
1182   /* Write data to line buffer */
1183 
1184   _ensure_buffer_size(p, p->buffer_end + 64);
1185 
1186   switch (p->format) {
1187 
1188   case CS_TIME_PLOT_DAT:
1189 
1190     if (p->use_iteration)
1191       p->buffer_end += sprintf(p->buffer + p->buffer_end, " %8d", tn);
1192     else
1193       p->buffer_end += sprintf(p->buffer + p->buffer_end, " %14.7e", t);
1194 
1195     for (i = 0; i < n_vals; i++) {
1196       _ensure_buffer_size(p, p->buffer_end + 64);
1197       p->buffer_end += sprintf(p->buffer + p->buffer_end, " %14.7e",
1198                                (double)(vals[i]));
1199     }
1200 
1201     p->buffer_end += sprintf(p->buffer + p->buffer_end, "\n");
1202 
1203     break;
1204 
1205   case CS_TIME_PLOT_CSV:
1206 
1207     if (p->use_iteration)
1208       p->buffer_end += sprintf(p->buffer + p->buffer_end, "%8d", tn);
1209     else
1210       p->buffer_end += sprintf(p->buffer + p->buffer_end, "%14.7e", t);
1211 
1212     for (i = 0; i < n_vals; i++) {
1213       _ensure_buffer_size(p, p->buffer_end + 64);
1214       p->buffer_end += sprintf(p->buffer + p->buffer_end, ", %14.7e",
1215                                (double)(vals[i]));
1216     }
1217 
1218     p->buffer_end += sprintf(p->buffer + p->buffer_end, "\n");
1219 
1220     break;
1221 
1222   default:
1223     break;
1224   }
1225 
1226   /* Output buffer if necessary */
1227 
1228   _plot_file_check_or_write(p);
1229 }
1230 
1231 /*----------------------------------------------------------------------------
1232  * Flush buffered values to file if applicable
1233  *
1234  * parameters:
1235  *   p <-> time plot values file handler
1236  *----------------------------------------------------------------------------*/
1237 
1238 void
cs_time_plot_flush(cs_time_plot_t * p)1239 cs_time_plot_flush(cs_time_plot_t  *p)
1240 {
1241   /* Force buffered variant output */
1242 
1243   if (p->buffer_end > 0) {
1244     if (p->buffer_steps[0] > 0)
1245       p->buffer_steps[1] = p->buffer_steps[0];
1246     _plot_file_check_or_write(p);
1247   }
1248 
1249   if (p->f != NULL) {
1250     if (p->flush_times[0] > 0)
1251       p->flush_times[1] = cs_timer_wtime();
1252     fflush(p->f);
1253   }
1254 }
1255 
1256 /*----------------------------------------------------------------------------
1257  * flush all time plots
1258  *----------------------------------------------------------------------------*/
1259 
1260 void
cs_time_plot_flush_all(void)1261 cs_time_plot_flush_all(void)
1262 {
1263   for (cs_time_plot_t *p = _plots_head; p != NULL; p = p->next)
1264     cs_time_plot_flush(p);
1265 }
1266 
1267 /*----------------------------------------------------------------------------
1268  * Set time plot file writer flush behavior defaults.
1269  *
1270  * parameters:
1271  *   flush_wtime     <-- elapsed time interval between file flushes;
1272  *                       if < 0, no forced flush
1273  *   n_buffer_steps  <-- number of time steps in output buffer if
1274  *                       file is not to be kept open
1275  *----------------------------------------------------------------------------*/
1276 
1277 void
cs_time_plot_set_flush_default(float flush_wtime,int n_buffer_steps)1278 cs_time_plot_set_flush_default(float  flush_wtime,
1279                                int    n_buffer_steps)
1280 {
1281   _flush_wtime_default    = flush_wtime;
1282   _n_buffer_steps_default = n_buffer_steps;
1283 }
1284 
1285 /*----------------------------------------------------------------------------
1286  * Return time plot file writer flush behavior defaults.
1287  *
1288  * parameters:
1289  *   flush_wtime     --> elapsed time interval between file flushes;
1290  *                       if < 0, no forced flush (NULL if not queried)
1291  *   n_buffer_steps  <-- number of time steps in output buffer if
1292  *                       file is not to be kept open (NULL if not queried)
1293  *----------------------------------------------------------------------------*/
1294 
1295 void
cs_time_plot_get_flush_default(float * flush_wtime,int * n_buffer_steps)1296 cs_time_plot_get_flush_default(float  *flush_wtime,
1297                                int    *n_buffer_steps)
1298 {
1299   if (flush_wtime != NULL)
1300     *flush_wtime = _flush_wtime_default;
1301 
1302   if (n_buffer_steps != NULL)
1303     *n_buffer_steps = _n_buffer_steps_default;
1304 }
1305 
1306 /*----------------------------------------------------------------------------*/
1307 
1308 END_C_DECLS
1309