LAL  7.5.0.1-b72065a
LogPrintf.c
Go to the documentation of this file.
1 /*
2  * Copyright (C) 2008 Karl Wette
3  * Copyright (C) 2005 Reinhard Prix
4  *
5  * [partially based on the MSG_LOG class in BOINC:
6  * Copyright (C) 2005 University of California]
7  *
8  * This program is free software; you can redistribute it and/or modify
9  * it under the terms of the GNU General Public License as published by
10  * the Free Software Foundation; either version 2 of the License, or
11  * (at your option) any later version.
12  *
13  * This program is distributed in the hope that it will be useful,
14  * but WITHOUT ANY WARRANTY; without even the implied warranty of
15  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16  * GNU General Public License for more details.
17  *
18  * You should have received a copy of the GNU General Public License
19  * along with with program; see the file COPYING. If not, write to the
20  * Free Software Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston,
21  * MA 02110-1301 USA
22  */
23 
24 /* Windows version fixed by Bernd Machenschalk */
25 
26 #include <config.h>
27 
28 /*---------- INCLUDES ----------*/
29 #include <stdio.h>
30 #include <string.h>
31 #include <strings.h>
32 
33 #include <lal/XLALError.h>
34 #include <lal/LALMalloc.h>
35 #include <lal/LALDebugLevel.h>
36 #include <lal/Date.h>
37 
38 #ifdef _MSC_VER
39 #include <Windows.h>
40 #include <process.h>
41 #define getpid _getpid
42 #else
43 #ifdef HAVE_SYS_TIME_H
44 #include <sys/time.h>
45 #endif
46 #ifdef HAVE_SYS_TYPES_H
47 #include <sys/types.h>
48 #endif
49 #ifdef HAVE_UNISTD_H
50 #include <unistd.h>
51 #endif
52 #endif
53 
54 #ifdef HAVE_SYS_RESOURCE_H
55 #include <sys/resource.h>
56 #endif
57 
58 #include <time.h>
59 
60 #ifndef HAVE_LOCALTIME_R
61 #define localtime_r(timep, result) memcpy((result), localtime(timep), sizeof(struct tm))
62 #endif
63 
64 #include <lal/LogPrintf.h>
65 
66 /*---------- internal types ----------*/
67 
68 /*---------- empty initializers ---------- */
69 
70 /*---------- internal variables ----------*/
71 static FILE *logFile = NULL;
72 
73 /*---------- internal prototypes ----------*/
74 static FILE* LogFile(void);
75 
76 static const char * LogTimeToString(double t);
77 
78 static const char *LogFormatLevel( LogLevel_t level );
79 
80 /*==================== FUNCTION DEFINITIONS ====================*/
81 
82 /**
83  * Get log level by examining lalDebugLevel
84  */
86 {
87  if (lalDebugLevel == 0)
88  return LOG_NONE; // If not printing LAL messages, also do not print log messages
89  if (lalDebugLevel & LALINFO)
90  return LOG_DETAIL; // Print LOG_DETAIL messages if LAL_DEBUG_LEVEL contains 'info'
92  return LOG_DEBUG; // Print LOG_DEBUG messages if LAL_DEBUG_LEVEL contains 'warning'
93  return LOG_NORMAL; // Print LOG_CRITICAL and LOG_NORMAL messages by default
94 }
95 
96 /** Set file to print log messages to */
97 void LogSetFile(FILE *fp)
98 {
99  logFile = fp;
100 }
101 
102 /** Decide where to print log messages */
103 static FILE* LogFile(void)
104 {
105  if (logFile)
106  return logFile; // Print all messages to 'logFile' if defined
107  if (LogLevel() < LOG_NORMAL)
108  return stderr; // Error log messages are printed to standard error
109  return stdout; // All other log messages are printed to standard output
110 }
111 
112 /**
113  * Verbatim output of the given log-message if LogLevel() >= level
114  */
115 void
116 LogPrintfVerbatim (LogLevel_t level, const char* format, ...)
117 {
118 
119  if ( LogLevel() < level )
120  return;
121 
122  va_list va;
123  va_start(va, format);
124 
125  /* simply print this to output */
126  vfprintf (LogFile(), format, va );
127  fflush(LogFile());
128 
129  va_end(va);
130 
131 } /* LogPrintfVerbatim() */
132 
133 
134 /**
135  * Output the given log-message, prefixed by a timestamp and level, if LogLevel() >= level
136  */
137 void
138 LogPrintf (LogLevel_t level, const char* format, ...)
139 {
140 
141  if ( LogLevel() < level )
142  return;
143 
144  va_list va;
145  va_start(va, format);
146 
147  fprintf(LogFile(), "%s (%d) [%s]: ", LogGetTimestamp(), getpid(), LogFormatLevel(level) );
148  vfprintf(LogFile(), format, va);
149  fflush(LogFile());
150 
151  va_end(va);
152 
153 } /* LogPrintf() */
154 
155 
156 /* taken from BOINC: return time-string for given unix-time
157  */
158 const char *
159 LogTimeToString ( double t )
160 {
161  static char buf[100];
162  char finer[16];
163  time_t x = (time_t)t;
164  struct tm tm;
165  localtime_r(&x, &tm);
166 
167  int hundreds_of_microseconds=(int)(10000*(t-(int)t));
168 
169  if (hundreds_of_microseconds == 10000) {
170  /* paranoia -- this should never happen! */
171  hundreds_of_microseconds=0;
172  t+=1.0;
173  }
174 
175  strftime(buf, sizeof(buf), "%Y-%m-%d %H:%M:%S", &tm);
176  sprintf(finer, ".%04d", hundreds_of_microseconds);
177  strcat(buf, finer);
178 
179  return buf;
180 
181 } /* LogTimeToString() */
182 
183 ///
184 /// Returns the peak amount of memory (in MB) allocated on the heap so far
185 /// using either lalMallocTotalPeak if memory-debugging is active or getrusage (if available),
186 /// otherwise returns -1 (without error) for "dont know"
187 ///
188 /// \note the reported number always refers to the 'nominal' memory usage *not* including any extra "padding"
189 /// that would have been added by LALMalloc(), which typically doubles memory usage.
190 ///
191 REAL8
193 {
194  // first see if lal's memory-debugging can be used
195  if ( lalDebugLevel & LALMEMPADBIT ) {
196  return lalMallocTotalPeak / (1024.0 * 1024.0); // lalMallocTotalMax counts bytes, doesn't inlude 'padding'
197  }
198 
199  // otherwise try using getrusage
200 #ifdef HAVE_SYS_RESOURCE_H
201  struct rusage usage;
202  XLAL_CHECK_REAL8 ( getrusage ( RUSAGE_SELF, &usage ) == 0, XLAL_ESYS, "call to getrusage() failed with errno = %d\n", errno );
203  REAL8 peakHeapMB = usage.ru_maxrss / 1024.0; // maxrss is in KB
204 #ifdef __APPLE__
205  // on macOS rusage.ru_maxrss is returned in bytes
206  // https://developer.apple.com/legacy/library/documentation/Darwin/Reference/ManPages/man2/getrusage.2.html
207  peakHeapMB /= 1024.0;
208 #endif
209 
210  if ( lalDebugLevel & LALMEMPADBIT ) {
211  peakHeapMB /= 2.0; // try to correct for memory-padding added by LALMalloc(), which seems ~factor of 2
212  }
213  // we're suspicious of a value of '0', which can also indicate that ru_maxrss is unsupported on this platform
214  if ( usage.ru_maxrss > 0 ) {
215  return peakHeapMB;
216  }
217 #endif
218 
219  return -1; // fallback answer: "dont know"
220 
221 } // XLALGetMaxHeapUsageMB()
222 
223 /**
224  * Return time of day (seconds since 1970) as a double.
225  * Taken from BOINC's dtime():
226  *
227  */
228 REAL8
230 {
231 #ifdef _MSC_VER
232  /* Windows version of dtime() from BOINC.
233  Compile switch is MS compiler macro,
234  because I suspect gettimeofday should be present in MinGW */
235 #define EPOCHFILETIME_SEC (11644473600.)
236 #define TEN_MILLION 10000000.
237 
238  LARGE_INTEGER time;
239  FILETIME sysTime;
240  double t;
241  GetSystemTimeAsFileTime(&sysTime);
242  time.LowPart = sysTime.dwLowDateTime;
243  time.HighPart = sysTime.dwHighDateTime; /* Time is in 100 ns units */
244  t = (double)time.QuadPart; /* Convert to 1 s units */
245  t /= TEN_MILLION; /* In seconds */
246  t -= EPOCHFILETIME_SEC; /* Offset to the Epoch time */
247  return t;
248 #else
249 
250  struct timeval tv;
251  if ( gettimeofday ( &tv, NULL ) != 0 ) {
252  XLALPrintError ("%s: call to gettimeofday() failed with errno = %d\n", __func__, errno );
254  }
255 
256  return tv.tv_sec + (tv.tv_usec/1.e6);
257 
258 #endif
259 } /* XLALGetTimeOfDay() */
260 
261 
262 ///
263 /// High-resolution CPU timer (returns result in seconds), aimed for code-timing purposes.
264 /// Attempts to provide the highest time resolution available, while adding as little overhead as possible.
265 ///
266 /// \note Returns 0.0 if no CPU timer is available.
267 ///
268 ///
269 REAL8
271 {
272 #if defined(HAVE_CLOCK_GETTIME) && HAVE_DECL_CLOCK_PROCESS_CPUTIME_ID
273 
274  struct timespec ut;
275  if ( clock_gettime ( CLOCK_PROCESS_CPUTIME_ID, &ut ) != 0 ) {
276  return 0.0;
277  }
278  return ut.tv_sec + ut.tv_nsec * 1.e-9;
279 
280 #else // no CPU timer available
281 
282  return 0.0;
283 
284 #endif
285 } // XLALGetCPUTime()
286 
287 
288 /* returns static timestamps-string for 'now' */
289 const char *
291 {
292  return ( LogTimeToString ( XLALGetTimeOfDay() ) );
293 
294 } /* LogGetTimestamp() */
295 
296 /* returns static string characterizing this log-level */
297 static const char *
299 {
300  static char buf[100];
301 
302  switch ( level )
303  {
304  case LOG_CRITICAL:
305  sprintf (buf, "CRITICAL" );
306  break;
307  case LOG_NORMAL:
308  sprintf (buf, "normal");
309  break;
310  case LOG_DEBUG:
311  sprintf (buf, "debug");
312  break;
313  case LOG_DETAIL:
314  sprintf (buf, "detail");
315  break;
316  default:
317  sprintf (buf, "unknown(%d)", level );
318  break;
319  } /* switch(level) */
320 
321  return buf;
322 
323 } /* LogFormatLevel() */
324 
325 
326 
327 /**
328  * Output gsl_matrix in octave-format, using the given format for the matrix-entries
329  * return -1 on error, 0 if OK.
330  */
331 int
332 XLALfprintfGSLmatrix ( FILE *fp, const char *fmt, const gsl_matrix *gij )
333 {
334  int cols, rows;
335  int i, j;
336 
337  /* check user input */
338  if ( !fp || !fmt )
339  return -1;
340 
341  if ( gij == NULL ) /* treat as valid empty matrix */
342  {
343  fprintf ( fp, "[ ];\n");
344  return 0;
345  }
346 
347  rows = gij->size1;
348  cols = gij->size2;
349 
350  fprintf (fp, " [ \n" );
351  for ( i=0; i < rows; i ++ )
352  {
353  for (j=0; j < cols; j ++ )
354  {
355  fprintf (fp, fmt, gsl_matrix_get ( gij, i, j ) );
356  if ( j < cols - 1 )
357  fprintf (fp, ", ");
358  else
359  fprintf (fp, ";\n");
360  } /* for j < cols */
361  }
362  fprintf (fp, " ];\n" );
363 
364  return 0;
365 
366 } /* XLALprintGSLmatrix() */
367 
368 
369 /**
370  * Output gsl_matrix in octave-format, using the given format for the matrix-entries
371  * return -1 on error, 0 if OK.
372  */
373 int
374 XLALfprintfGSLvector ( FILE *fp, const char *fmt, const gsl_vector *vect )
375 {
376  int rows;
377  int i;
378 
379  /* check user input */
380  if ( !vect || !fp || !fmt )
381  return -1;
382 
383  rows = vect->size;
384 
385  fprintf (fp, " [ " );
386  for ( i=0; i < rows; i ++ )
387  {
388  fprintf (fp, fmt, gsl_vector_get ( vect, i ) );
389  if ( i < rows - 1 )
390  fprintf (fp, ", ");
391  } /* for i < rows */
392 
393  fprintf (fp, " ];\n" );
394 
395  return 0;
396 
397 } /* XLALprintGSLvector() */
398 
399 int
400 XLALfprintfGSLvector_int ( FILE *fp, const char *fmt, const gsl_vector_int *vect )
401 {
402  int rows;
403  int i;
404 
405  /* check user input */
406  if ( !vect || !fp || !fmt )
407  return -1;
408 
409  rows = vect->size;
410 
411  fprintf (fp, " [ " );
412  for ( i=0; i < rows; i ++ )
413  {
414  fprintf (fp, fmt, gsl_vector_int_get ( vect, i ) );
415  if ( i < rows - 1 )
416  fprintf (fp, ", ");
417  } /* for i < rows */
418 
419  fprintf (fp, " ];\n" );
420 
421  return 0;
422 
423 } /* XLALprintGSLvector_int() */
424 
425 
426 /**
427  * Returns input string with line-breaks '\n' removed (replaced by space)
428  * The original string is unmodified. The returned string is allocated here.
429  */
430 char *
431 XLALClearLinebreaks ( const char *str )
432 {
433  char *ret, *tmp;
434 
435  if ( !str )
436  return NULL;
437 
438  if ( (ret = LALMalloc( strlen ( str ) + 1) ) == NULL )
439  return NULL;
440  strcpy ( ret, str );
441 
442  tmp = ret;
443  while ( (tmp = strchr ( tmp, '\n' ) ) )
444  {
445  *tmp = ' ';
446  tmp ++;
447  }
448 
449  return ret;
450 
451 } /* XLALClearLinebreaks() */
452 
453 
454 /** dump given REAL4 time-series into a text-file */
455 int
456 XLALdumpREAL4TimeSeries ( const char *fname, const REAL4TimeSeries *series )
457 {
458  XLAL_CHECK ( fname != NULL, XLAL_EINVAL );
459  XLAL_CHECK ( series != NULL, XLAL_EINVAL );
460 
461  FILE *fp;
462  XLAL_CHECK ( (fp = fopen (fname, "wb")) != NULL, XLAL_ESYS );
463 
464  REAL8 dt = series->deltaT;
465  UINT4 numSamples = series->data->length;
466  char buf[256];
467  for ( UINT4 i = 0; i < numSamples; i++ )
468  {
469  LIGOTimeGPS ti = series->epoch;
470  XLALGPSAdd ( &ti, i * dt );
471  XLALGPSToStr ( buf, &ti );
472  fprintf( fp, "%20s %20.16g\n", buf, series->data->data[i] );
473  }
474  fclose ( fp );
475 
476  return XLAL_SUCCESS;
477 
478 } // XLALdumpREAL4TimeSeries()
479 
480 /** dump given REAL8 time-series into a text-file */
481 int
482 XLALdumpREAL8TimeSeries ( const char *fname, const REAL8TimeSeries *series )
483 {
484  XLAL_CHECK ( fname != NULL, XLAL_EINVAL );
485  XLAL_CHECK ( series != NULL, XLAL_EINVAL );
486 
487  FILE *fp;
488  XLAL_CHECK ( (fp = fopen (fname, "wb")) != NULL, XLAL_ESYS );
489 
490  REAL8 dt = series->deltaT;
491  UINT4 numSamples = series->data->length;
492  char buf[256];
493  for ( UINT4 i = 0; i < numSamples; i++ )
494  {
495  LIGOTimeGPS ti = series->epoch;
496  XLALGPSAdd ( &ti, i * dt );
497  XLALGPSToStr ( buf, &ti );
498  fprintf( fp, "%20s %20.16g\n", buf, series->data->data[i] );
499  }
500  fclose ( fp );
501 
502  return XLAL_SUCCESS;
503 
504 } // XLALdumpREAL8TimeSeries()
505 
506 
507 /** dump given COMPLEX8 time-series into a text-file */
508 int
509 XLALdumpCOMPLEX8TimeSeries ( const char *fname, const COMPLEX8TimeSeries *series )
510 {
511  XLAL_CHECK ( fname != NULL, XLAL_EINVAL );
512  XLAL_CHECK ( series != NULL, XLAL_EINVAL );
513 
514  FILE *fp;
515  XLAL_CHECK ( (fp = fopen (fname, "wb")) != NULL, XLAL_ESYS );
516 
517  REAL8 dt = series->deltaT;
518  UINT4 numSamples = series->data->length;
519  char buf[256];
520  for ( UINT4 i = 0; i < numSamples; i++ )
521  {
522  LIGOTimeGPS ti = series->epoch;
523  XLALGPSAdd ( &ti, i * dt );
524  XLALGPSToStr ( buf, &ti );
525  fprintf( fp, "%20s %20.16g %20.16g\n", buf, crealf(series->data->data[i]), cimagf(series->data->data[i]) );
526  }
527  fclose ( fp );
528 
529  return XLAL_SUCCESS;
530 
531 } // XLALdumpCOMPLEX8TimeSeries()
#define usage(program)
#define LALMalloc(n)
Definition: LALMalloc.h:93
static FILE * LogFile(void)
Decide where to print log messages.
Definition: LogPrintf.c:103
void LogPrintfVerbatim(LogLevel_t level, const char *format,...)
Verbatim output of the given log-message if LogLevel() >= level.
Definition: LogPrintf.c:116
static const char * LogTimeToString(double t)
Definition: LogPrintf.c:159
static FILE * logFile
Definition: LogPrintf.c:71
static const char * LogFormatLevel(LogLevel_t level)
Definition: LogPrintf.c:298
#define localtime_r(timep, result)
Definition: LogPrintf.c:61
void LogPrintf(LogLevel_t level, const char *format,...)
Output the given log-message, prefixed by a timestamp and level, if LogLevel() >= level.
Definition: LogPrintf.c:138
#define fprintf
int XLALPrintError(const char *fmt,...)
Definition: XLALError.c:68
char * XLALGPSToStr(char *, const LIGOTimeGPS *t)
Return a string containing the ASCII base 10 representation of a LIGOTimeGPS.
Definition: StrToGPS.c:274
double REAL8
Double precision real floating-point number (8 bytes).
uint32_t UINT4
Four-byte unsigned integer.
#define lalDebugLevel
Definition: LALDebugLevel.h:58
@ LALINFO
enable info messages
Definition: LALDebugLevel.h:48
@ LALWARNING
enable warning messages
Definition: LALDebugLevel.h:47
@ LALMEMPADBIT
enable memory padding
Definition: LALDebugLevel.h:38
size_t lalMallocTotalPeak
peak amount of memory allocated so far
Definition: LALMalloc.c:39
void LogSetFile(FILE *fp)
Set file to print log messages to.
Definition: LogPrintf.c:97
int XLALdumpREAL8TimeSeries(const char *fname, const REAL8TimeSeries *series)
dump given REAL8 time-series into a text-file
Definition: LogPrintf.c:482
REAL8 XLALGetCPUTime(void)
High-resolution CPU timer (returns result in seconds), aimed for code-timing purposes.
Definition: LogPrintf.c:270
int XLALdumpREAL4TimeSeries(const char *fname, const REAL4TimeSeries *series)
dump given REAL4 time-series into a text-file
Definition: LogPrintf.c:456
REAL8 XLALGetPeakHeapUsageMB(void)
Returns the peak amount of memory (in MB) allocated on the heap so far using either lalMallocTotalPea...
Definition: LogPrintf.c:192
const char * LogGetTimestamp(void)
Definition: LogPrintf.c:290
LogLevel_t
Argument-type for LogPrintf(): determines log-level of this message.
Definition: LogPrintf.h:53
int XLALdumpCOMPLEX8TimeSeries(const char *fname, const COMPLEX8TimeSeries *series)
dump given COMPLEX8 time-series into a text-file
Definition: LogPrintf.c:509
int XLALfprintfGSLmatrix(FILE *fp, const char *fmt, const gsl_matrix *gij)
Output gsl_matrix in octave-format, using the given format for the matrix-entries return -1 on error,...
Definition: LogPrintf.c:332
int XLALfprintfGSLvector_int(FILE *fp, const char *fmt, const gsl_vector_int *vect)
Definition: LogPrintf.c:400
int XLALfprintfGSLvector(FILE *fp, const char *fmt, const gsl_vector *vect)
Output gsl_matrix in octave-format, using the given format for the matrix-entries return -1 on error,...
Definition: LogPrintf.c:374
REAL8 XLALGetTimeOfDay(void)
Return time of day (seconds since 1970) as a double.
Definition: LogPrintf.c:229
LogLevel_t LogLevel(void)
Get log level by examining lalDebugLevel.
Definition: LogPrintf.c:85
char * XLALClearLinebreaks(const char *str)
Returns input string with line-breaks ' ' removed (replaced by space) The original string is unmodifi...
Definition: LogPrintf.c:431
@ LOG_CRITICAL
log-level for critical errors
Definition: LogPrintf.h:55
@ LOG_NONE
internal: don't use
Definition: LogPrintf.h:54
@ LOG_DEBUG
debug log-level
Definition: LogPrintf.h:57
@ LOG_NORMAL
'normal' log-level
Definition: LogPrintf.h:56
@ LOG_DETAIL
detailed log-level
Definition: LogPrintf.h:58
#define XLAL_ERROR_REAL8(...)
Macro to invoke a failure from a XLAL routine returning a REAL8.
Definition: XLALError.h:752
#define XLAL_CHECK(assertion,...)
Macro to test an assertion and invoke a failure if it is not true in a function that returns an integ...
Definition: XLALError.h:810
#define XLAL_CHECK_REAL8(assertion,...)
Macro to test an assertion and invoke a failure if it is not true in a function that returns a REAL8.
Definition: XLALError.h:870
@ XLAL_SUCCESS
Success return value (not an error number)
Definition: XLALError.h:401
@ XLAL_ESYS
System error.
Definition: XLALError.h:442
@ XLAL_EINVAL
Invalid argument.
Definition: XLALError.h:409
LIGOTimeGPS * XLALGPSAdd(LIGOTimeGPS *epoch, REAL8 dt)
Adds a double to a GPS time.
Definition: XLALTime.c:131
Time series of COMPLEX8 data, see DATATYPE-TimeSeries types for more details.
Definition: LALDatatypes.h:590
REAL8 deltaT
The time step between samples of the time series in seconds.
Definition: LALDatatypes.h:593
LIGOTimeGPS epoch
The start time of the time series.
Definition: LALDatatypes.h:592
COMPLEX8Sequence * data
The sequence of sampled data.
Definition: LALDatatypes.h:596
UINT4 length
Number of elements in array.
Definition: LALDatatypes.h:167
COMPLEX8 * data
Pointer to the data array.
Definition: LALDatatypes.h:168
Epoch relative to GPS epoch, see LIGOTimeGPS type for more details.
Definition: LALDatatypes.h:458
Time series of REAL4 data, see DATATYPE-TimeSeries types for more details.
Definition: LALDatatypes.h:570
REAL4Sequence * data
The sequence of sampled data.
Definition: LALDatatypes.h:576
REAL8 deltaT
The time step between samples of the time series in seconds.
Definition: LALDatatypes.h:573
LIGOTimeGPS epoch
The start time of the time series.
Definition: LALDatatypes.h:572
REAL4 * data
Pointer to the data array.
Definition: LALDatatypes.h:150
UINT4 length
Number of elements in array.
Definition: LALDatatypes.h:149
Time series of REAL8 data, see DATATYPE-TimeSeries types for more details.
Definition: LALDatatypes.h:580
REAL8Sequence * data
The sequence of sampled data.
Definition: LALDatatypes.h:586
REAL8 deltaT
The time step between samples of the time series in seconds.
Definition: LALDatatypes.h:583
LIGOTimeGPS epoch
The start time of the time series.
Definition: LALDatatypes.h:582
REAL8 * data
Pointer to the data array.
Definition: LALDatatypes.h:159
UINT4 length
Number of elements in array.
Definition: LALDatatypes.h:158
FILE * fp
Definition: tconvert.c:105