2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
54 #include "thread_mpi/threads.h"
56 #include "gromacs/fileio/filenm.h"
57 #include "gromacs/fileio/md5.h"
58 #include "gromacs/utility/fatalerror.h"
59 #include "gromacs/utility/futil.h"
60 #include "gromacs/utility/smalloc.h"
62 #include "gmxfio-impl.h"
64 /* This is the new improved and thread safe version of gmxfio. */
68 /* the list of open files is a linked list, with a dummy element at its head;
69 it is initialized when the first file is opened. */
70 static t_fileio *open_files = NULL;
73 /* this mutex locks the open_files structure so that no two threads can
76 For now, we use this as a coarse grained lock on all file
77 insertion/deletion operations because it makes avoiding deadlocks
78 easier, and adds almost no overhead: the only overhead is during
79 opening and closing of files, or during global operations like
80 iterating along all open files. All these cases should be rare
81 during the simulation. */
82 static tMPI_Thread_mutex_t open_file_mutex = TMPI_THREAD_MUTEX_INITIALIZER;
84 /******************************************************************
88 ******************************************************************/
90 static int gmx_fio_int_flush(t_fileio* fio)
102 /* lock the mutex associated with this fio. This needs to be done for every
103 type of access to the fio's elements. */
104 void gmx_fio_lock(t_fileio *fio)
106 tMPI_Lock_lock(&(fio->mtx));
108 /* unlock the mutex associated with this fio. */
109 void gmx_fio_unlock(t_fileio *fio)
111 tMPI_Lock_unlock(&(fio->mtx));
114 /* make a dummy head element, assuming we locked everything. */
115 static void gmx_fio_make_dummy(void)
120 open_files->fp = NULL;
121 open_files->fn = NULL;
122 open_files->next = open_files;
123 open_files->prev = open_files;
124 tMPI_Lock_init(&(open_files->mtx));
134 /***********************************************************************
136 * FILE LIST OPERATIONS
138 ***********************************************************************/
141 /* insert a new t_fileio into the list */
142 static void gmx_fio_insert(t_fileio *fio)
145 /* first lock the big open_files mutex. */
146 tMPI_Thread_mutex_lock(&open_file_mutex);
147 /* now check whether the dummy element has been allocated,
148 and allocate it if it hasn't */
149 gmx_fio_make_dummy();
151 /* and lock the fio we got and the list's head **/
153 gmx_fio_lock(open_files);
154 prev = open_files->prev;
155 /* lock the element after the current one */
156 if (prev != open_files)
161 /* now do the actual insertion: */
162 fio->next = open_files;
163 open_files->prev = fio;
167 /* now unlock all our locks */
168 if (prev != open_files)
170 gmx_fio_unlock(prev);
172 gmx_fio_unlock(open_files);
175 /* now unlock the big open_files mutex. */
176 tMPI_Thread_mutex_unlock(&open_file_mutex);
179 /* remove a t_fileio into the list. We assume the fio is locked, and we leave
181 NOTE: We also assume that the open_file_mutex has been locked */
182 static void gmx_fio_remove(t_fileio *fio)
184 /* lock prev, because we're changing it */
185 gmx_fio_lock(fio->prev);
187 /* now set the prev's pointer */
188 fio->prev->next = fio->next;
189 gmx_fio_unlock(fio->prev);
191 /* with the next ptr, we can simply lock while the original was locked */
192 gmx_fio_lock(fio->next);
193 fio->next->prev = fio->prev;
194 gmx_fio_unlock(fio->next);
196 /* and make sure we point nowhere in particular */
197 fio->next = fio->prev = fio;
201 /* get the first open file, or NULL if there is none.
202 Returns a locked fio. */
203 static t_fileio *gmx_fio_get_first(void)
206 /* first lock the big open_files mutex and the dummy's mutex */
208 /* first lock the big open_files mutex. */
209 tMPI_Thread_mutex_lock(&open_file_mutex);
210 gmx_fio_make_dummy();
212 gmx_fio_lock(open_files);
213 ret = open_files->next;
216 /* check whether there were any to begin with */
217 if (ret == open_files)
219 /* after this, the open_file pointer should never change */
224 gmx_fio_lock(open_files->next);
226 gmx_fio_unlock(open_files);
232 /* get the next open file, or NULL if there is none.
233 Unlocks the previous fio and locks the next one. */
234 static t_fileio *gmx_fio_get_next(t_fileio *fio)
239 /* check if that was the last one */
240 if (fio->next == open_files)
243 tMPI_Thread_mutex_unlock(&open_file_mutex);
254 /* Stop looping through the open_files. Unlocks the global lock. */
255 static void gmx_fio_stop_getting_next(t_fileio *fio)
258 tMPI_Thread_mutex_unlock(&open_file_mutex);
264 /*****************************************************************
268 *****************************************************************/
269 t_fileio *gmx_fio_open(const char *fn, const char *mode)
271 t_fileio *fio = NULL;
273 gmx_bool bRead, bReadWrite;
275 /* sanitize the mode string */
276 if (strncmp(mode, "r+", 2) == 0)
278 strcpy(newmode, "r+");
280 else if (mode[0] == 'r')
282 strcpy(newmode, "r");
284 else if (strncmp(mode, "w+", 2) == 0)
286 strcpy(newmode, "w+");
288 else if (mode[0] == 'w')
290 strcpy(newmode, "w");
292 else if (strncmp(mode, "a+", 2) == 0)
294 strcpy(newmode, "a+");
296 else if (mode[0] == 'a')
298 strcpy(newmode, "a");
302 gmx_fatal(FARGS, "DEATH HORROR in gmx_fio_open, mode is '%s'", mode);
305 /* Check if it should be opened as a binary file */
306 if (!ftp_is_text(fn2ftp(fn)))
308 strcat(newmode, "b");
312 tMPI_Lock_init(&(fio->mtx));
313 bRead = (newmode[0] == 'r' && newmode[1] != '+');
314 bReadWrite = (newmode[1] == '+');
319 if (fn2ftp(fn) == efTNG)
321 gmx_incons("gmx_fio_open may not be used to open TNG files");
323 fio->iFTP = fn2ftp(fn);
324 fio->fn = gmx_strdup(fn);
326 fio->fp = gmx_ffopen(fn, newmode);
327 /* If this file type is in the list of XDR files, open it like that */
328 if (ftp_is_xdr(fio->iFTP))
330 /* determine the XDR direction */
331 if (newmode[0] == 'w' || newmode[0] == 'a')
333 fio->xdrmode = XDR_ENCODE;
337 fio->xdrmode = XDR_DECODE;
340 xdrstdio_create(fio->xdr, fio->fp, fio->xdrmode);
343 /* for appending seek to end of file to make sure ftell gives correct position
344 * important for checkpointing */
345 if (newmode[0] == 'a')
347 gmx_fseek(fio->fp, 0, SEEK_END);
351 fio->bReadWrite = bReadWrite;
352 fio->bDouble = (sizeof(real) == sizeof(double));
354 /* and now insert this file into the list of open files. */
359 static int gmx_fio_close_locked(t_fileio *fio)
363 if (fio->xdr != NULL)
365 xdr_destroy(fio->xdr);
371 rc = gmx_ffclose(fio->fp); /* fclose returns 0 if happy */
378 int gmx_fio_close(t_fileio *fio)
382 /* first lock the big open_files mutex. */
383 /* We don't want two processes operating on the list at the same time */
384 tMPI_Thread_mutex_lock(&open_file_mutex);
387 /* first remove it from the list */
389 rc = gmx_fio_close_locked(fio);
395 tMPI_Thread_mutex_unlock(&open_file_mutex);
400 /* close only fp but keep FIO entry. */
401 int gmx_fio_fp_close(t_fileio *fio)
405 if (fio->xdr == NULL)
407 rc = gmx_ffclose(fio->fp); /* fclose returns 0 if happy */
415 FILE * gmx_fio_fopen(const char *fn, const char *mode)
420 fio = gmx_fio_open(fn, mode);
428 int gmx_fio_fclose(FILE *fp)
433 cur = gmx_fio_get_first();
438 rc = gmx_fio_close_locked(cur);
440 gmx_fio_stop_getting_next(cur);
445 cur = gmx_fio_get_next(cur);
451 /* internal variant of get_file_md5 that operates on a locked file */
452 static int gmx_fio_int_get_file_md5(t_fileio *fio, gmx_off_t offset,
453 unsigned char digest[])
455 /*1MB: large size important to catch almost identical files */
456 #define CPT_CHK_LEN 1048576
460 gmx_off_t seek_offset;
463 seek_offset = offset - CPT_CHK_LEN;
468 read_len = offset - seek_offset;
471 if (fio->fp && fio->bReadWrite)
473 ret = gmx_fseek(fio->fp, seek_offset, SEEK_SET);
476 gmx_fseek(fio->fp, 0, SEEK_END);
479 if (ret) /*either no fp, not readwrite, or fseek not successful */
484 snew(buf, CPT_CHK_LEN);
485 /* the read puts the file position back to offset */
486 if ((gmx_off_t)fread(buf, 1, read_len, fio->fp) != read_len)
488 /* not fatal: md5sum check to prevent overwriting files
489 * works (less safe) without
493 fprintf(stderr, "\nTrying to get md5sum: %s: %s\n", fio->fn,
496 else if (feof(fio->fp))
499 * For long runs that checkpoint frequently but write e.g. logs
500 * infrequently we don't want to issue lots of warnings before we
501 * have written anything to the log.
505 fprintf(stderr, "\nTrying to get md5sum: EOF: %s\n", fio->fn);
512 "\nTrying to get md5sum: Unknown reason for short read: %s\n",
516 gmx_fseek(fio->fp, 0, SEEK_END);
520 gmx_fseek(fio->fp, 0, SEEK_END); /*is already at end, but under windows
521 it gives problems otherwise*/
525 fprintf(debug, "chksum %s readlen %ld\n", fio->fn, (long int)read_len);
530 gmx_md5_init(&state);
531 gmx_md5_append(&state, buf, read_len);
532 gmx_md5_finish(&state, digest);
541 * fio: file to compute md5 for
542 * offset: starting pointer of region to use for md5
543 * digest: return array of md5 sum
545 int gmx_fio_get_file_md5(t_fileio *fio, gmx_off_t offset,
546 unsigned char digest[])
551 ret = gmx_fio_int_get_file_md5(fio, offset, digest);
557 /* The fio_mutex should ALWAYS be locked when this function is called */
558 static int gmx_fio_int_get_file_position(t_fileio *fio, gmx_off_t *offset)
560 /* Flush the file, so we are sure it is written */
561 if (gmx_fio_int_flush(fio))
566 "Cannot write file '%s'; maybe you are out of disk space?",
571 /* We cannot count on XDR being able to write 64-bit integers,
572 so separate into high/low 32-bit values.
573 In case the filesystem has 128-bit offsets we only care
574 about the first 64 bits - we'll have to fix
575 this when exabyte-size output files are common...
577 *offset = gmx_ftell(fio->fp);
582 int gmx_fio_get_output_file_positions(gmx_file_position_t **p_outputfiles,
586 gmx_file_position_t * outputfiles;
591 /* pre-allocate 100 files */
593 snew(outputfiles, nalloc);
595 cur = gmx_fio_get_first();
598 /* Skip the checkpoint files themselves, since they could be open when
599 we call this routine... */
600 if (!cur->bRead && cur->iFTP != efCPT)
602 /* This is an output file currently open for writing, add it */
603 if (nfiles == nalloc)
606 srenew(outputfiles, nalloc);
609 strncpy(outputfiles[nfiles].filename, cur->fn, STRLEN - 1);
611 /* Get the file position */
612 gmx_fio_int_get_file_position(cur, &outputfiles[nfiles].offset);
614 outputfiles[nfiles].chksum_size
615 = gmx_fio_int_get_file_md5(cur,
616 outputfiles[nfiles].offset,
617 outputfiles[nfiles].chksum);
622 cur = gmx_fio_get_next(cur);
625 *p_outputfiles = outputfiles;
631 char *gmx_fio_getname(t_fileio *fio)
641 int gmx_fio_getftp(t_fileio* fio)
652 void gmx_fio_rewind(t_fileio* fio)
658 xdr_destroy(fio->xdr);
660 xdrstdio_create(fio->xdr, fio->fp, fio->xdrmode);
670 int gmx_fio_flush(t_fileio* fio)
675 ret = gmx_fio_int_flush(fio);
683 static int gmx_fio_int_fsync(t_fileio *fio)
689 rc = gmx_fsync(fio->fp);
696 int gmx_fio_fsync(t_fileio *fio)
701 rc = gmx_fio_int_fsync(fio);
709 t_fileio *gmx_fio_all_output_fsync(void)
711 t_fileio *ret = NULL;
714 cur = gmx_fio_get_first();
719 /* if any of them fails, return failure code */
720 int rc = gmx_fio_int_fsync(cur);
726 cur = gmx_fio_get_next(cur);
729 /* in addition, we force these to be written out too, if they're being
730 redirected. We don't check for errors because errors most likely mean
731 that they're not redirected. */
734 #if (defined(HAVE_FSYNC))
735 /* again, fahcore defines HAVE_FSYNC and fsync() */
736 fsync(STDOUT_FILENO);
737 fsync(STDERR_FILENO);
744 gmx_off_t gmx_fio_ftell(t_fileio* fio)
751 ret = gmx_ftell(fio->fp);
757 int gmx_fio_seek(t_fileio* fio, gmx_off_t fpos)
764 rc = gmx_fseek(fio->fp, fpos, SEEK_SET);
775 FILE *gmx_fio_getfp(t_fileio *fio)
788 gmx_bool gmx_fio_getread(t_fileio* fio)
799 int xtc_seek_time(t_fileio *fio, real time, int natoms, gmx_bool bSeekForwardOnly)
804 ret = xdr_xtc_seek_time(time, fio->fp, fio->xdr, natoms, bSeekForwardOnly);