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, 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.
45 #ifdef HAVE_SYS_TIME_H
51 #ifdef GMX_NATIVE_WINDOWS
55 #include "types/commrec.h"
59 #include "gromacs/utility/futil.h"
60 #include "gromacs/fileio/filenm.h"
61 #include "gromacs/fileio/gmxfio.h"
64 #include "gromacs/utility/basenetwork.h"
65 #include "gromacs/utility/cstringutil.h"
66 #include "gromacs/utility/exceptions.h"
67 #include "gromacs/utility/fatalerror.h"
68 #include "gromacs/utility/gmxmpi.h"
69 #include "gromacs/utility/programcontext.h"
70 #include "gromacs/utility/smalloc.h"
72 /* The source code in this file should be thread-safe.
73 Please keep it that way. */
77 static void par_fn(char *base
, int ftp
, const t_commrec
*cr
,
78 gmx_bool bAppendSimId
, gmx_bool bAppendNodeId
,
79 char buf
[], int bufsize
)
81 if ((size_t)bufsize
< (strlen(base
)+10))
83 gmx_mem("Character buffer too small!");
86 /* Copy to buf, and strip extension */
88 buf
[strlen(base
) - strlen(ftp2ext(fn2ftp(base
))) - 1] = '\0';
92 sprintf(buf
+strlen(buf
), "%d", cr
->ms
->sim
);
97 sprintf(buf
+strlen(buf
), "%d", cr
->nodeid
);
101 /* Add extension again */
102 strcat(buf
, (ftp
== efTPX
) ? "tpr" : (ftp
== efEDR
) ? "edr" : ftp2ext(ftp
));
105 fprintf(debug
, "rank %d par_fn '%s'\n", cr
->nodeid
, buf
);
106 if (fn2ftp(buf
) == efLOG
)
108 fprintf(debug
, "log\n");
113 void check_multi_int(FILE *log
, const gmx_multisim_t
*ms
, int val
,
118 gmx_bool bCompatible
;
120 if (NULL
!= log
&& !bQuiet
)
122 fprintf(log
, "Multi-checking %s ... ", name
);
128 "check_multi_int called with a NULL communication pointer");
131 snew(ibuf
, ms
->nsim
);
133 gmx_sumi_sim(ms
->nsim
, ibuf
, ms
);
136 for (p
= 1; p
< ms
->nsim
; p
++)
138 bCompatible
= bCompatible
&& (ibuf
[p
-1] == ibuf
[p
]);
143 if (NULL
!= log
&& !bQuiet
)
145 fprintf(log
, "OK\n");
152 fprintf(log
, "\n%s is not equal for all subsystems\n", name
);
153 for (p
= 0; p
< ms
->nsim
; p
++)
155 fprintf(log
, " subsystem %d: %d\n", p
, ibuf
[p
]);
158 gmx_fatal(FARGS
, "The %d subsystems are not compatible\n", ms
->nsim
);
164 void check_multi_int64(FILE *log
, const gmx_multisim_t
*ms
,
165 gmx_int64_t val
, const char *name
,
170 gmx_bool bCompatible
;
172 if (NULL
!= log
&& !bQuiet
)
174 fprintf(log
, "Multi-checking %s ... ", name
);
180 "check_multi_int called with a NULL communication pointer");
183 snew(ibuf
, ms
->nsim
);
185 gmx_sumli_sim(ms
->nsim
, ibuf
, ms
);
188 for (p
= 1; p
< ms
->nsim
; p
++)
190 bCompatible
= bCompatible
&& (ibuf
[p
-1] == ibuf
[p
]);
195 if (NULL
!= log
&& !bQuiet
)
197 fprintf(log
, "OK\n");
204 fprintf(log
, "\n%s is not equal for all subsystems\n", name
);
205 for (p
= 0; p
< ms
->nsim
; p
++)
208 /* first make the format string */
209 snprintf(strbuf
, 255, " subsystem %%d: %s\n",
211 fprintf(log
, strbuf
, p
, ibuf
[p
]);
214 gmx_fatal(FARGS
, "The %d subsystems are not compatible\n", ms
->nsim
);
221 void gmx_log_open(const char *lognm
, const t_commrec
*cr
, gmx_bool bMasterOnly
,
222 gmx_bool bAppendFiles
, FILE** fplog
)
225 char buf
[256], host
[256];
227 char timebuf
[STRLEN
];
233 /* Communicate the filename for logfile */
234 if (cr
->nnodes
> 1 && !bMasterOnly
235 #ifdef GMX_THREAD_MPI
236 /* With thread MPI the non-master log files are opened later
237 * when the files names are already known on all nodes.
245 len
= strlen(lognm
) + 1;
247 gmx_bcast(sizeof(len
), &len
, cr
);
254 tmpnm
= gmx_strdup(lognm
);
256 gmx_bcast(len
*sizeof(*tmpnm
), tmpnm
, cr
);
260 tmpnm
= gmx_strdup(lognm
);
265 if (!bMasterOnly
&& !MASTER(cr
))
267 /* Since log always ends with '.log' let's use this info */
268 par_fn(tmpnm
, efLOG
, cr
, FALSE
, !bMasterOnly
, buf
, 255);
269 fp
= gmx_fio_fopen(buf
, bAppendFiles
? "a+" : "w+" );
271 else if (!bAppendFiles
)
273 fp
= gmx_fio_fopen(tmpnm
, bAppendFiles
? "a+" : "w+" );
278 gmx_fatal_set_log_file(fp
);
280 /* Get some machine parameters */
281 gmx_gethostname(host
, 256);
286 # ifdef GMX_NATIVE_WINDOWS
300 "-----------------------------------------------------------\n"
301 "Restarting from checkpoint, appending to previous log file.\n"
306 gmx_ctime_r(&t
, timebuf
, STRLEN
);
309 "Log file opened on %s"
310 "Host: %s pid: %d rank ID: %d number of ranks: %d\n",
311 timebuf
, host
, pid
, cr
->nodeid
, cr
->nnodes
);
314 gmx::BinaryInformationSettings settings
;
315 settings
.extendedInfo(true);
316 settings
.copyright(!bAppendFiles
);
317 gmx::printBinaryInformation(fp
, gmx::getProgramContext(), settings
);
319 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
;
328 void gmx_log_close(FILE *fp
)
332 gmx_fatal_set_log_file(NULL
);
337 void init_multisystem(t_commrec
*cr
, int nsim
, char **multidirs
,
338 int nfile
, const t_filenm fnm
[], gmx_bool bParFn
)
341 int nnodes
, nnodpersim
, sim
, i
, ftp
;
344 MPI_Group mpi_group_world
;
351 gmx_fatal(FARGS
, "This binary is compiled without MPI support, can not do multiple simulations.");
356 if (nnodes
% nsim
!= 0)
358 gmx_fatal(FARGS
, "The number of ranks (%d) is not a multiple of the number of simulations (%d)", nnodes
, nsim
);
361 nnodpersim
= nnodes
/nsim
;
362 sim
= cr
->nodeid
/nnodpersim
;
366 fprintf(debug
, "We have %d simulations, %d ranks per simulation, local simulation is %d\n", nsim
, nnodpersim
, sim
);
374 /* Create a communicator for the master nodes */
375 snew(rank
, ms
->nsim
);
376 for (i
= 0; i
< ms
->nsim
; i
++)
378 rank
[i
] = i
*nnodpersim
;
380 MPI_Comm_group(MPI_COMM_WORLD
, &mpi_group_world
);
381 MPI_Group_incl(mpi_group_world
, nsim
, rank
, &ms
->mpi_group_masters
);
383 MPI_Comm_create(MPI_COMM_WORLD
, ms
->mpi_group_masters
,
384 &ms
->mpi_comm_masters
);
386 #if !defined(MPI_IN_PLACE_EXISTS)
387 /* initialize the MPI_IN_PLACE replacement buffers */
389 ms
->mpb
->ibuf
= NULL
;
390 ms
->mpb
->libuf
= NULL
;
391 ms
->mpb
->fbuf
= NULL
;
392 ms
->mpb
->dbuf
= NULL
;
393 ms
->mpb
->ibuf_alloc
= 0;
394 ms
->mpb
->libuf_alloc
= 0;
395 ms
->mpb
->fbuf_alloc
= 0;
396 ms
->mpb
->dbuf_alloc
= 0;
401 /* Reduce the intra-simulation communication */
402 cr
->sim_nodeid
= cr
->nodeid
% nnodpersim
;
403 cr
->nnodes
= nnodpersim
;
405 MPI_Comm_split(MPI_COMM_WORLD
, sim
, cr
->sim_nodeid
, &cr
->mpi_comm_mysim
);
406 cr
->mpi_comm_mygroup
= cr
->mpi_comm_mysim
;
407 cr
->nodeid
= cr
->sim_nodeid
;
412 fprintf(debug
, "This is simulation %d", cr
->ms
->sim
);
415 fprintf(debug
, ", local number of ranks %d, local rank ID %d",
416 cr
->nnodes
, cr
->sim_nodeid
);
418 fprintf(debug
, "\n\n");
425 fprintf(debug
, "Changing to directory %s\n", multidirs
[cr
->ms
->sim
]);
427 gmx_chdir(multidirs
[cr
->ms
->sim
]);
431 /* Patch output and tpx, cpt and rerun input file names */
432 for (i
= 0; (i
< nfile
); i
++)
434 /* Because of possible multiple extensions per type we must look
435 * at the actual file name
437 if (is_output(&fnm
[i
]) ||
438 fnm
[i
].ftp
== efTPX
|| fnm
[i
].ftp
== efCPT
||
439 strcmp(fnm
[i
].opt
, "-rerun") == 0)
441 ftp
= fn2ftp(fnm
[i
].fns
[0]);
442 par_fn(fnm
[i
].fns
[0], ftp
, cr
, TRUE
, FALSE
, buf
, 255);
443 sfree(fnm
[i
].fns
[0]);
444 fnm
[i
].fns
[0] = gmx_strdup(buf
);