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)
100 rc = fflush((FILE *) fio->xdr->x_private);
106 /* lock the mutex associated with this fio. This needs to be done for every
107 type of access to the fio's elements. */
108 void gmx_fio_lock(t_fileio *fio)
110 tMPI_Lock_lock(&(fio->mtx));
112 /* unlock the mutex associated with this fio. */
113 void gmx_fio_unlock(t_fileio *fio)
115 tMPI_Lock_unlock(&(fio->mtx));
118 /* make a dummy head element, assuming we locked everything. */
119 static void gmx_fio_make_dummy(void)
124 open_files->fp = NULL;
125 open_files->fn = NULL;
126 open_files->next = open_files;
127 open_files->prev = open_files;
128 tMPI_Lock_init(&(open_files->mtx));
138 /***********************************************************************
140 * FILE LIST OPERATIONS
142 ***********************************************************************/
145 /* insert a new t_fileio into the list */
146 static void gmx_fio_insert(t_fileio *fio)
149 /* first lock the big open_files mutex. */
150 tMPI_Thread_mutex_lock(&open_file_mutex);
151 /* now check whether the dummy element has been allocated,
152 and allocate it if it hasn't */
153 gmx_fio_make_dummy();
155 /* and lock the fio we got and the list's head **/
157 gmx_fio_lock(open_files);
158 prev = open_files->prev;
159 /* lock the element after the current one */
160 if (prev != open_files)
165 /* now do the actual insertion: */
166 fio->next = open_files;
167 open_files->prev = fio;
171 /* now unlock all our locks */
172 if (prev != open_files)
174 gmx_fio_unlock(prev);
176 gmx_fio_unlock(open_files);
179 /* now unlock the big open_files mutex. */
180 tMPI_Thread_mutex_unlock(&open_file_mutex);
183 /* remove a t_fileio into the list. We assume the fio is locked, and we leave
185 NOTE: We also assume that the open_file_mutex has been locked */
186 static void gmx_fio_remove(t_fileio *fio)
188 /* lock prev, because we're changing it */
189 gmx_fio_lock(fio->prev);
191 /* now set the prev's pointer */
192 fio->prev->next = fio->next;
193 gmx_fio_unlock(fio->prev);
195 /* with the next ptr, we can simply lock while the original was locked */
196 gmx_fio_lock(fio->next);
197 fio->next->prev = fio->prev;
198 gmx_fio_unlock(fio->next);
200 /* and make sure we point nowhere in particular */
201 fio->next = fio->prev = fio;
205 /* get the first open file, or NULL if there is none.
206 Returns a locked fio. */
207 static t_fileio *gmx_fio_get_first(void)
210 /* first lock the big open_files mutex and the dummy's mutex */
212 /* first lock the big open_files mutex. */
213 tMPI_Thread_mutex_lock(&open_file_mutex);
214 gmx_fio_make_dummy();
216 gmx_fio_lock(open_files);
217 ret = open_files->next;
220 /* check whether there were any to begin with */
221 if (ret == open_files)
223 /* after this, the open_file pointer should never change */
228 gmx_fio_lock(open_files->next);
230 gmx_fio_unlock(open_files);
236 /* get the next open file, or NULL if there is none.
237 Unlocks the previous fio and locks the next one. */
238 static t_fileio *gmx_fio_get_next(t_fileio *fio)
243 /* check if that was the last one */
244 if (fio->next == open_files)
247 tMPI_Thread_mutex_unlock(&open_file_mutex);
258 /* Stop looping through the open_files. Unlocks the global lock. */
259 static void gmx_fio_stop_getting_next(t_fileio *fio)
262 tMPI_Thread_mutex_unlock(&open_file_mutex);
268 /*****************************************************************
272 *****************************************************************/
273 t_fileio *gmx_fio_open(const char *fn, const char *mode)
275 t_fileio *fio = NULL;
277 gmx_bool bRead, bReadWrite;
279 /* sanitize the mode string */
280 if (strncmp(mode, "r+", 2) == 0)
282 strcpy(newmode, "r+");
284 else if (mode[0] == 'r')
286 strcpy(newmode, "r");
288 else if (strncmp(mode, "w+", 2) == 0)
290 strcpy(newmode, "w+");
292 else if (mode[0] == 'w')
294 strcpy(newmode, "w");
296 else if (strncmp(mode, "a+", 2) == 0)
298 strcpy(newmode, "a+");
300 else if (mode[0] == 'a')
302 strcpy(newmode, "a");
306 gmx_fatal(FARGS, "DEATH HORROR in gmx_fio_open, mode is '%s'", mode);
309 /* Check if it should be opened as a binary file */
310 if (!ftp_is_text(fn2ftp(fn)))
312 /* Not ascii, add b to file mode */
313 if ((strchr(newmode, 'b') == NULL) && (strchr(newmode, 'B') == NULL))
315 strcat(newmode, "b");
320 tMPI_Lock_init(&(fio->mtx));
321 bRead = (newmode[0] == 'r' && newmode[1] != '+');
322 bReadWrite = (newmode[1] == '+');
327 if (fn2ftp(fn) == efTNG)
329 gmx_incons("gmx_fio_open may not be used to open TNG files");
331 fio->iFTP = fn2ftp(fn);
332 fio->fn = gmx_strdup(fn);
334 /* If this file type is in the list of XDR files, open it like that */
335 if (ftp_is_xdr(fio->iFTP))
337 /* First check whether we have to make a backup,
338 * only for writing, not for read or append.
340 if (newmode[0] == 'w')
343 /* only make backups for normal gromacs */
349 /* Check whether file exists */
356 fio->fp = gmx_ffopen(fn, newmode);
358 /* determine the XDR direction */
359 if (newmode[0] == 'w' || newmode[0] == 'a')
361 fio->xdrmode = XDR_ENCODE;
365 fio->xdrmode = XDR_DECODE;
369 xdrstdio_create(fio->xdr, fio->fp, fio->xdrmode);
373 /* If it is not, open it as a regular file */
374 fio->fp = gmx_ffopen(fn, newmode);
377 /* for appending seek to end of file to make sure ftell gives correct position
378 * important for checkpointing */
379 if (newmode[0] == 'a')
381 gmx_fseek(fio->fp, 0, SEEK_END);
385 fio->bReadWrite = bReadWrite;
386 fio->bDouble = (sizeof(real) == sizeof(double));
390 /* and now insert this file into the list of open files. */
395 static int gmx_fio_close_locked(t_fileio *fio)
401 gmx_fatal(FARGS, "File %s closed twice!\n", fio->fn);
404 if (ftp_is_xdr(fio->iFTP))
406 xdr_destroy(fio->xdr);
412 rc = gmx_ffclose(fio->fp); /* fclose returns 0 if happy */
420 int gmx_fio_close(t_fileio *fio)
424 /* first lock the big open_files mutex. */
425 /* We don't want two processes operating on the list at the same time */
426 tMPI_Thread_mutex_lock(&open_file_mutex);
428 if (fio->iFTP == efTNG)
430 gmx_incons("gmx_fio_close should not be called on a TNG file");
433 /* first remove it from the list */
435 rc = gmx_fio_close_locked(fio);
441 tMPI_Thread_mutex_unlock(&open_file_mutex);
446 /* close only fp but keep FIO entry. */
447 int gmx_fio_fp_close(t_fileio *fio)
451 if (!ftp_is_xdr(fio->iFTP))
453 rc = gmx_ffclose(fio->fp); /* fclose returns 0 if happy */
461 FILE * gmx_fio_fopen(const char *fn, const char *mode)
466 fio = gmx_fio_open(fn, mode);
474 int gmx_fio_fclose(FILE *fp)
479 cur = gmx_fio_get_first();
484 rc = gmx_fio_close_locked(cur);
486 gmx_fio_stop_getting_next(cur);
491 cur = gmx_fio_get_next(cur);
497 /* internal variant of get_file_md5 that operates on a locked file */
498 static int gmx_fio_int_get_file_md5(t_fileio *fio, gmx_off_t offset,
499 unsigned char digest[])
501 /*1MB: large size important to catch almost identical files */
502 #define CPT_CHK_LEN 1048576
506 gmx_off_t seek_offset;
509 seek_offset = offset - CPT_CHK_LEN;
514 read_len = offset - seek_offset;
517 if (fio->fp && fio->bReadWrite)
519 ret = gmx_fseek(fio->fp, seek_offset, SEEK_SET);
522 gmx_fseek(fio->fp, 0, SEEK_END);
525 if (ret) /*either no fp, not readwrite, or fseek not successful */
530 snew(buf, CPT_CHK_LEN);
531 /* the read puts the file position back to offset */
532 if ((gmx_off_t)fread(buf, 1, read_len, fio->fp) != read_len)
534 /* not fatal: md5sum check to prevent overwriting files
535 * works (less safe) without
539 fprintf(stderr, "\nTrying to get md5sum: %s: %s\n", fio->fn,
542 else if (feof(fio->fp))
545 * For long runs that checkpoint frequently but write e.g. logs
546 * infrequently we don't want to issue lots of warnings before we
547 * have written anything to the log.
551 fprintf(stderr, "\nTrying to get md5sum: EOF: %s\n", fio->fn);
558 "\nTrying to get md5sum: Unknown reason for short read: %s\n",
562 gmx_fseek(fio->fp, 0, SEEK_END);
566 gmx_fseek(fio->fp, 0, SEEK_END); /*is already at end, but under windows
567 it gives problems otherwise*/
571 fprintf(debug, "chksum %s readlen %ld\n", fio->fn, (long int)read_len);
576 gmx_md5_init(&state);
577 gmx_md5_append(&state, buf, read_len);
578 gmx_md5_finish(&state, digest);
587 * fio: file to compute md5 for
588 * offset: starting pointer of region to use for md5
589 * digest: return array of md5 sum
591 int gmx_fio_get_file_md5(t_fileio *fio, gmx_off_t offset,
592 unsigned char digest[])
597 ret = gmx_fio_int_get_file_md5(fio, offset, digest);
603 /* The fio_mutex should ALWAYS be locked when this function is called */
604 static int gmx_fio_int_get_file_position(t_fileio *fio, gmx_off_t *offset)
606 /* Flush the file, so we are sure it is written */
607 if (gmx_fio_int_flush(fio))
612 "Cannot write file '%s'; maybe you are out of disk space?",
617 /* We cannot count on XDR being able to write 64-bit integers,
618 so separate into high/low 32-bit values.
619 In case the filesystem has 128-bit offsets we only care
620 about the first 64 bits - we'll have to fix
621 this when exabyte-size output files are common...
623 *offset = gmx_ftell(fio->fp);
628 int gmx_fio_get_output_file_positions(gmx_file_position_t **p_outputfiles,
632 gmx_file_position_t * outputfiles;
637 /* pre-allocate 100 files */
639 snew(outputfiles, nalloc);
641 cur = gmx_fio_get_first();
644 /* Skip the checkpoint files themselves, since they could be open when
645 we call this routine... */
646 if (cur->bOpen && !cur->bRead && cur->iFTP != efCPT)
648 /* This is an output file currently open for writing, add it */
649 if (nfiles == nalloc)
652 srenew(outputfiles, nalloc);
655 strncpy(outputfiles[nfiles].filename, cur->fn, STRLEN - 1);
657 /* Get the file position */
658 gmx_fio_int_get_file_position(cur, &outputfiles[nfiles].offset);
660 outputfiles[nfiles].chksum_size
661 = gmx_fio_int_get_file_md5(cur,
662 outputfiles[nfiles].offset,
663 outputfiles[nfiles].chksum);
668 cur = gmx_fio_get_next(cur);
671 *p_outputfiles = outputfiles;
677 void gmx_fio_setprecision(t_fileio *fio, gmx_bool bDouble)
680 fio->bDouble = bDouble;
684 void gmx_fio_setdebug(t_fileio *fio, gmx_bool bDebug)
687 fio->bDebug = bDebug;
691 char *gmx_fio_getname(t_fileio *fio)
701 int gmx_fio_getftp(t_fileio* fio)
712 void gmx_fio_rewind(t_fileio* fio)
718 xdr_destroy(fio->xdr);
720 xdrstdio_create(fio->xdr, fio->fp, fio->xdrmode);
730 int gmx_fio_flush(t_fileio* fio)
735 ret = gmx_fio_int_flush(fio);
743 static int gmx_fio_int_fsync(t_fileio *fio)
749 rc = gmx_fsync(fio->fp);
751 else if (fio->xdr) /* this should normally not happen */
753 rc = gmx_fsync((FILE*) fio->xdr->x_private);
754 /* ^ is this actually OK? */
761 int gmx_fio_fsync(t_fileio *fio)
766 rc = gmx_fio_int_fsync(fio);
774 t_fileio *gmx_fio_all_output_fsync(void)
776 t_fileio *ret = NULL;
779 cur = gmx_fio_get_first();
782 if (cur->bOpen && !cur->bRead)
784 /* if any of them fails, return failure code */
785 int rc = gmx_fio_int_fsync(cur);
791 cur = gmx_fio_get_next(cur);
794 /* in addition, we force these to be written out too, if they're being
795 redirected. We don't check for errors because errors most likely mean
796 that they're not redirected. */
799 #if (defined(HAVE_FSYNC))
800 /* again, fahcore defines HAVE_FSYNC and fsync() */
801 fsync(STDOUT_FILENO);
802 fsync(STDERR_FILENO);
809 gmx_off_t gmx_fio_ftell(t_fileio* fio)
816 ret = gmx_ftell(fio->fp);
822 int gmx_fio_seek(t_fileio* fio, gmx_off_t fpos)
829 rc = gmx_fseek(fio->fp, fpos, SEEK_SET);
840 FILE *gmx_fio_getfp(t_fileio *fio)
853 gmx_bool gmx_fio_getread(t_fileio* fio)
864 int xtc_seek_time(t_fileio *fio, real time, int natoms, gmx_bool bSeekForwardOnly)
869 ret = xdr_xtc_seek_time(time, fio->fp, fio->xdr, natoms, bSeekForwardOnly);