Remove all unnecessary HAVE_CONFIG_H
[gromacs.git] / src / gromacs / gmxlib / main.cpp
blobe58bb0e3881e34cefdcd45376fd78d317de87bf8
1 /*
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.
37 #include "config.h"
39 #include <stdio.h>
40 #include <stdlib.h>
41 #include <string.h>
42 #include <limits.h>
43 #include <time.h>
45 #ifdef HAVE_SYS_TIME_H
46 #include <sys/time.h>
47 #endif
48 #ifdef HAVE_UNISTD_H
49 #include <unistd.h>
50 #endif
51 #ifdef GMX_NATIVE_WINDOWS
52 #include <process.h>
53 #endif
55 #include "types/commrec.h"
56 #include "network.h"
57 #include "main.h"
58 #include "macros.h"
59 #include "gromacs/utility/futil.h"
60 #include "gromacs/fileio/filenm.h"
61 #include "gromacs/fileio/gmxfio.h"
62 #include "copyrite.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. */
75 #define BUFSIZE 1024
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 */
87 strcpy(buf, base);
88 buf[strlen(base) - strlen(ftp2ext(fn2ftp(base))) - 1] = '\0';
90 if (bAppendSimId)
92 sprintf(buf+strlen(buf), "%d", cr->ms->sim);
94 if (bAppendNodeId)
96 strcat(buf, "_rank");
97 sprintf(buf+strlen(buf), "%d", cr->nodeid);
99 strcat(buf, ".");
101 /* Add extension again */
102 strcat(buf, (ftp == efTPX) ? "tpr" : (ftp == efEDR) ? "edr" : ftp2ext(ftp));
103 if (debug)
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,
114 const char *name,
115 gmx_bool bQuiet)
117 int *ibuf, p;
118 gmx_bool bCompatible;
120 if (NULL != log && !bQuiet)
122 fprintf(log, "Multi-checking %s ... ", name);
125 if (ms == NULL)
127 gmx_fatal(FARGS,
128 "check_multi_int called with a NULL communication pointer");
131 snew(ibuf, ms->nsim);
132 ibuf[ms->sim] = val;
133 gmx_sumi_sim(ms->nsim, ibuf, ms);
135 bCompatible = TRUE;
136 for (p = 1; p < ms->nsim; p++)
138 bCompatible = bCompatible && (ibuf[p-1] == ibuf[p]);
141 if (bCompatible)
143 if (NULL != log && !bQuiet)
145 fprintf(log, "OK\n");
148 else
150 if (NULL != log)
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);
161 sfree(ibuf);
164 void check_multi_int64(FILE *log, const gmx_multisim_t *ms,
165 gmx_int64_t val, const char *name,
166 gmx_bool bQuiet)
168 gmx_int64_t *ibuf;
169 int p;
170 gmx_bool bCompatible;
172 if (NULL != log && !bQuiet)
174 fprintf(log, "Multi-checking %s ... ", name);
177 if (ms == NULL)
179 gmx_fatal(FARGS,
180 "check_multi_int called with a NULL communication pointer");
183 snew(ibuf, ms->nsim);
184 ibuf[ms->sim] = val;
185 gmx_sumli_sim(ms->nsim, ibuf, ms);
187 bCompatible = TRUE;
188 for (p = 1; p < ms->nsim; p++)
190 bCompatible = bCompatible && (ibuf[p-1] == ibuf[p]);
193 if (bCompatible)
195 if (NULL != log && !bQuiet)
197 fprintf(log, "OK\n");
200 else
202 if (NULL != log)
204 fprintf(log, "\n%s is not equal for all subsystems\n", name);
205 for (p = 0; p < ms->nsim; p++)
207 char strbuf[255];
208 /* first make the format string */
209 snprintf(strbuf, 255, " subsystem %%d: %s\n",
210 "%" GMX_PRId64);
211 fprintf(log, strbuf, p, ibuf[p]);
214 gmx_fatal(FARGS, "The %d subsystems are not compatible\n", ms->nsim);
217 sfree(ibuf);
221 void gmx_log_open(const char *lognm, const t_commrec *cr, gmx_bool bMasterOnly,
222 gmx_bool bAppendFiles, FILE** fplog)
224 int len, pid;
225 char buf[256], host[256];
226 time_t t;
227 char timebuf[STRLEN];
228 FILE *fp = *fplog;
229 char *tmpnm;
231 debug_gmx();
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.
239 && FALSE
240 #endif
243 if (MASTER(cr))
245 len = strlen(lognm) + 1;
247 gmx_bcast(sizeof(len), &len, cr);
248 if (!MASTER(cr))
250 snew(tmpnm, len+8);
252 else
254 tmpnm = gmx_strdup(lognm);
256 gmx_bcast(len*sizeof(*tmpnm), tmpnm, cr);
258 else
260 tmpnm = gmx_strdup(lognm);
263 debug_gmx();
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+" );
276 sfree(tmpnm);
278 gmx_fatal_set_log_file(fp);
280 /* Get some machine parameters */
281 gmx_gethostname(host, 256);
283 time(&t);
285 #ifndef NO_GETPID
286 # ifdef GMX_NATIVE_WINDOWS
287 pid = _getpid();
288 # else
289 pid = getpid();
290 # endif
291 #else
292 pid = 0;
293 #endif
295 if (bAppendFiles)
297 fprintf(fp,
298 "\n"
299 "\n"
300 "-----------------------------------------------------------\n"
301 "Restarting from checkpoint, appending to previous log file.\n"
302 "\n"
306 gmx_ctime_r(&t, timebuf, STRLEN);
308 fprintf(fp,
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;
320 fprintf(fp, "\n\n");
322 fflush(fp);
323 debug_gmx();
325 *fplog = fp;
328 void gmx_log_close(FILE *fp)
330 if (fp)
332 gmx_fatal_set_log_file(NULL);
333 gmx_fio_fclose(fp);
337 void init_multisystem(t_commrec *cr, int nsim, char **multidirs,
338 int nfile, const t_filenm fnm[], gmx_bool bParFn)
340 gmx_multisim_t *ms;
341 int nnodes, nnodpersim, sim, i, ftp;
342 char buf[256];
343 #ifdef GMX_MPI
344 MPI_Group mpi_group_world;
345 int *rank;
346 #endif
348 #ifndef GMX_MPI
349 if (nsim > 1)
351 gmx_fatal(FARGS, "This binary is compiled without MPI support, can not do multiple simulations.");
353 #endif
355 nnodes = cr->nnodes;
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;
364 if (debug)
366 fprintf(debug, "We have %d simulations, %d ranks per simulation, local simulation is %d\n", nsim, nnodpersim, sim);
369 snew(ms, 1);
370 cr->ms = ms;
371 ms->nsim = nsim;
372 ms->sim = sim;
373 #ifdef GMX_MPI
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);
382 sfree(rank);
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 */
388 snew(ms->mpb, 1);
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;
397 #endif
399 #endif
401 /* Reduce the intra-simulation communication */
402 cr->sim_nodeid = cr->nodeid % nnodpersim;
403 cr->nnodes = nnodpersim;
404 #ifdef GMX_MPI
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;
408 #endif
410 if (debug)
412 fprintf(debug, "This is simulation %d", cr->ms->sim);
413 if (PAR(cr))
415 fprintf(debug, ", local number of ranks %d, local rank ID %d",
416 cr->nnodes, cr->sim_nodeid);
418 fprintf(debug, "\n\n");
421 if (multidirs)
423 if (debug)
425 fprintf(debug, "Changing to directory %s\n", multidirs[cr->ms->sim]);
427 gmx_chdir(multidirs[cr->ms->sim]);
429 else if (bParFn)
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);