Options for controlling startup output.
[alexxy/gromacs.git] / src / gromacs / gmxlib / main.c
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
2  *
3  *
4  *                This source code is part of
5  *
6  *                 G   R   O   M   A   C   S
7  *
8  *          GROningen MAchine for Chemical Simulations
9  *
10  *                        VERSION 3.2.0
11  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13  * Copyright (c) 2001-2004, The GROMACS development team,
14  * check out http://www.gromacs.org for more information.
15
16  * This program is free software; you can redistribute it and/or
17  * modify it under the terms of the GNU General Public License
18  * as published by the Free Software Foundation; either version 2
19  * of the License, or (at your option) any later version.
20  *
21  * If you want to redistribute modifications, please consider that
22  * scientific software is very special. Version control is crucial -
23  * bugs must be traceable. We will be happy to consider code for
24  * inclusion in the official distribution, but derived work must not
25  * be called official GROMACS. Details are found in the README & COPYING
26  * files - if they are missing, get the official version at www.gromacs.org.
27  *
28  * To help us fund GROMACS development, we humbly ask that you cite
29  * the papers on the package - you can find them in the top README file.
30  *
31  * For more info, check our website at http://www.gromacs.org
32  *
33  * And Hey:
34  * GROningen Mixture of Alchemy and Childrens' Stories
35  */
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
39 #include "gromacs/utility/gmx_header_config.h"
40
41 #include <stdio.h>
42 #include <stdlib.h>
43 #include <string.h>
44 #include <limits.h>
45 #include <time.h>
46
47 #ifdef HAVE_SYS_TIME_H
48 #include <sys/time.h>
49 #endif
50
51 #include "smalloc.h"
52 #include "gmx_fatal.h"
53 #include "network.h"
54 #include "main.h"
55 #include "macros.h"
56 #include "futil.h"
57 #include "filenm.h"
58 #include "gmxfio.h"
59 #include "string2.h"
60 #include "copyrite.h"
61
62 #ifdef GMX_THREAD_MPI
63 #include "thread_mpi.h"
64 #endif
65
66 /* The source code in this file should be thread-safe.
67          Please keep it that way. */
68
69
70 #ifdef HAVE_UNISTD_H
71 #include <unistd.h>
72 #endif
73
74 #ifdef GMX_NATIVE_WINDOWS
75 #include <process.h>
76 #endif
77
78
79 /* Portable version of ctime_r implemented in src/gmxlib/string2.c, but we do not want it declared in public installed headers */
80 char *
81 gmx_ctime_r(const time_t *clock, char *buf, int n);
82
83
84 #define BUFSIZE 1024
85
86
87 static void par_fn(char *base, int ftp, const t_commrec *cr,
88                    gmx_bool bAppendSimId, gmx_bool bAppendNodeId,
89                    char buf[], int bufsize)
90 {
91     int n;
92
93     if ((size_t)bufsize < (strlen(base)+10))
94     {
95         gmx_mem("Character buffer too small!");
96     }
97
98     /* Copy to buf, and strip extension */
99     strcpy(buf, base);
100     buf[strlen(base) - strlen(ftp2ext(fn2ftp(base))) - 1] = '\0';
101
102     if (bAppendSimId)
103     {
104         sprintf(buf+strlen(buf), "%d", cr->ms->sim);
105     }
106     if (bAppendNodeId)
107     {
108         strcat(buf, "_node");
109         sprintf(buf+strlen(buf), "%d", cr->nodeid);
110     }
111     strcat(buf, ".");
112
113     /* Add extension again */
114     strcat(buf, (ftp == efTPX) ? "tpr" : (ftp == efEDR) ? "edr" : ftp2ext(ftp));
115     if (debug)
116     {
117         fprintf(debug, "node %d par_fn '%s'\n", cr->nodeid, buf);
118         if (fn2ftp(buf) == efLOG)
119         {
120             fprintf(debug, "log\n");
121         }
122     }
123 }
124
125 void check_multi_int(FILE *log, const gmx_multisim_t *ms, int val,
126                      const char *name,
127                      gmx_bool bQuiet)
128 {
129     int     *ibuf, p;
130     gmx_bool bCompatible;
131
132     if (NULL != log && !bQuiet)
133     {
134         fprintf(log, "Multi-checking %s ... ", name);
135     }
136
137     if (ms == NULL)
138     {
139         gmx_fatal(FARGS,
140                   "check_multi_int called with a NULL communication pointer");
141     }
142
143     snew(ibuf, ms->nsim);
144     ibuf[ms->sim] = val;
145     gmx_sumi_sim(ms->nsim, ibuf, ms);
146
147     bCompatible = TRUE;
148     for (p = 1; p < ms->nsim; p++)
149     {
150         bCompatible = bCompatible && (ibuf[p-1] == ibuf[p]);
151     }
152
153     if (bCompatible)
154     {
155         if (NULL != log && !bQuiet)
156         {
157             fprintf(log, "OK\n");
158         }
159     }
160     else
161     {
162         if (NULL != log)
163         {
164             fprintf(log, "\n%s is not equal for all subsystems\n", name);
165             for (p = 0; p < ms->nsim; p++)
166             {
167                 fprintf(log, "  subsystem %d: %d\n", p, ibuf[p]);
168             }
169         }
170         gmx_fatal(FARGS, "The %d subsystems are not compatible\n", ms->nsim);
171     }
172
173     sfree(ibuf);
174 }
175
176 void check_multi_large_int(FILE *log, const gmx_multisim_t *ms,
177                            gmx_large_int_t val, const char *name,
178                            gmx_bool bQuiet)
179 {
180     gmx_large_int_t  *ibuf;
181     int               p;
182     gmx_bool          bCompatible;
183
184     if (NULL != log && !bQuiet)
185     {
186         fprintf(log, "Multi-checking %s ... ", name);
187     }
188
189     if (ms == NULL)
190     {
191         gmx_fatal(FARGS,
192                   "check_multi_int called with a NULL communication pointer");
193     }
194
195     snew(ibuf, ms->nsim);
196     ibuf[ms->sim] = val;
197     gmx_sumli_sim(ms->nsim, ibuf, ms);
198
199     bCompatible = TRUE;
200     for (p = 1; p < ms->nsim; p++)
201     {
202         bCompatible = bCompatible && (ibuf[p-1] == ibuf[p]);
203     }
204
205     if (bCompatible)
206     {
207         if (NULL != log && !bQuiet)
208         {
209             fprintf(log, "OK\n");
210         }
211     }
212     else
213     {
214         if (NULL != log)
215         {
216             fprintf(log, "\n%s is not equal for all subsystems\n", name);
217             for (p = 0; p < ms->nsim; p++)
218             {
219                 char strbuf[255];
220                 /* first make the format string */
221                 snprintf(strbuf, 255, "  subsystem %%d: %s\n",
222                          gmx_large_int_pfmt);
223                 fprintf(log, strbuf, p, ibuf[p]);
224             }
225         }
226         gmx_fatal(FARGS, "The %d subsystems are not compatible\n", ms->nsim);
227     }
228
229     sfree(ibuf);
230 }
231
232
233 char *gmx_gethostname(char *name, size_t len)
234 {
235     if (len < 8)
236     {
237         gmx_incons("gmx_gethostname called with len<8");
238     }
239 #ifdef HAVE_UNISTD_H
240     if (gethostname(name, len-1) != 0)
241     {
242         strncpy(name, "unknown", 8);
243     }
244 #else
245     strncpy(name, "unknown", 8);
246 #endif
247
248     return name;
249 }
250
251
252 void gmx_log_open(const char *lognm, const t_commrec *cr, gmx_bool bMasterOnly,
253                   gmx_bool bAppendFiles, FILE** fplog)
254 {
255     int    len, testlen, pid;
256     char   buf[256], host[256];
257     time_t t;
258     char   timebuf[STRLEN];
259     FILE  *fp = *fplog;
260     char  *tmpnm;
261
262     debug_gmx();
263
264     /* Communicate the filename for logfile */
265     if (cr->nnodes > 1 && !bMasterOnly
266 #ifdef GMX_THREAD_MPI
267         /* With thread MPI the non-master log files are opened later
268          * when the files names are already known on all nodes.
269          */
270         && FALSE
271 #endif
272         )
273     {
274         if (MASTER(cr))
275         {
276             len = strlen(lognm) + 1;
277         }
278         gmx_bcast(sizeof(len), &len, cr);
279         if (!MASTER(cr))
280         {
281             snew(tmpnm, len+8);
282         }
283         else
284         {
285             tmpnm = gmx_strdup(lognm);
286         }
287         gmx_bcast(len*sizeof(*tmpnm), tmpnm, cr);
288     }
289     else
290     {
291         tmpnm = gmx_strdup(lognm);
292     }
293
294     debug_gmx();
295
296     if (!bMasterOnly && !MASTER(cr))
297     {
298         /* Since log always ends with '.log' let's use this info */
299         par_fn(tmpnm, efLOG, cr, FALSE, !bMasterOnly, buf, 255);
300         fp = gmx_fio_fopen(buf, bAppendFiles ? "a+" : "w+" );
301     }
302     else if (!bAppendFiles)
303     {
304         fp = gmx_fio_fopen(tmpnm, bAppendFiles ? "a+" : "w+" );
305     }
306
307     sfree(tmpnm);
308
309     gmx_fatal_set_log_file(fp);
310
311     /* Get some machine parameters */
312     gmx_gethostname(host, 256);
313
314     time(&t);
315
316 #ifndef NO_GETPID
317 #   ifdef GMX_NATIVE_WINDOWS
318     pid = _getpid();
319 #   else
320     pid = getpid();
321 #   endif
322 #else
323     pid = 0;
324 #endif
325
326     if (bAppendFiles)
327     {
328         fprintf(fp,
329                 "\n"
330                 "\n"
331                 "-----------------------------------------------------------\n"
332                 "Restarting from checkpoint, appending to previous log file.\n"
333                 "\n"
334                 );
335     }
336
337     gmx_ctime_r(&t, timebuf, STRLEN);
338
339     fprintf(fp,
340             "Log file opened on %s"
341             "Host: %s  pid: %d  nodeid: %d  nnodes:  %d\n",
342             timebuf, host, pid, cr->nodeid, cr->nnodes);
343     gmx_print_version_info(fp);
344     fprintf(fp, "\n\n");
345
346     fflush(fp);
347     debug_gmx();
348
349     *fplog = fp;
350 }
351
352 void gmx_log_close(FILE *fp)
353 {
354     if (fp)
355     {
356         gmx_fatal_set_log_file(NULL);
357         gmx_fio_fclose(fp);
358     }
359 }
360
361 static void comm_args(const t_commrec *cr, int *argc, char ***argv)
362 {
363     int i, len;
364
365     if (PAR(cr))
366     {
367         gmx_bcast(sizeof(*argc), argc, cr);
368     }
369
370     if (!MASTER(cr))
371     {
372         snew(*argv, *argc+1);
373     }
374     if (debug)
375     {
376         fprintf(debug, "NODEID=%d argc=%d\n", cr->nodeid, *argc);
377     }
378     for (i = 0; (i < *argc); i++)
379     {
380         if (MASTER(cr))
381         {
382             len = strlen((*argv)[i])+1;
383         }
384         gmx_bcast(sizeof(len), &len, cr);
385         if (!MASTER(cr))
386         {
387             snew((*argv)[i], len);
388         }
389         /*gmx_bcast(len*sizeof((*argv)[i][0]),(*argv)[i],cr);*/
390         gmx_bcast(len*sizeof(char), (*argv)[i], cr);
391     }
392     debug_gmx();
393 }
394
395 void init_multisystem(t_commrec *cr, int nsim, char **multidirs,
396                       int nfile, const t_filenm fnm[], gmx_bool bParFn)
397 {
398     gmx_multisim_t *ms;
399     int             nnodes, nnodpersim, sim, i, ftp;
400     char            buf[256];
401 #ifdef GMX_MPI
402     MPI_Group       mpi_group_world;
403 #endif
404     int            *rank;
405
406 #ifndef GMX_MPI
407     if (nsim > 1)
408     {
409         gmx_fatal(FARGS, "This binary is compiled without MPI support, can not do multiple simulations.");
410     }
411 #endif
412
413     nnodes  = cr->nnodes;
414     if (nnodes % nsim != 0)
415     {
416         gmx_fatal(FARGS, "The number of nodes (%d) is not a multiple of the number of simulations (%d)", nnodes, nsim);
417     }
418
419     nnodpersim = nnodes/nsim;
420     sim        = cr->nodeid/nnodpersim;
421
422     if (debug)
423     {
424         fprintf(debug, "We have %d simulations, %d nodes per simulation, local simulation is %d\n", nsim, nnodpersim, sim);
425     }
426
427     snew(ms, 1);
428     cr->ms   = ms;
429     ms->nsim = nsim;
430     ms->sim  = sim;
431 #ifdef GMX_MPI
432     /* Create a communicator for the master nodes */
433     snew(rank, ms->nsim);
434     for (i = 0; i < ms->nsim; i++)
435     {
436         rank[i] = i*nnodpersim;
437     }
438     MPI_Comm_group(MPI_COMM_WORLD, &mpi_group_world);
439     MPI_Group_incl(mpi_group_world, nsim, rank, &ms->mpi_group_masters);
440     sfree(rank);
441     MPI_Comm_create(MPI_COMM_WORLD, ms->mpi_group_masters,
442                     &ms->mpi_comm_masters);
443
444 #if !defined(GMX_THREAD_MPI) && !defined(MPI_IN_PLACE_EXISTS)
445     /* initialize the MPI_IN_PLACE replacement buffers */
446     snew(ms->mpb, 1);
447     ms->mpb->ibuf        = NULL;
448     ms->mpb->libuf       = NULL;
449     ms->mpb->fbuf        = NULL;
450     ms->mpb->dbuf        = NULL;
451     ms->mpb->ibuf_alloc  = 0;
452     ms->mpb->libuf_alloc = 0;
453     ms->mpb->fbuf_alloc  = 0;
454     ms->mpb->dbuf_alloc  = 0;
455 #endif
456
457 #endif
458
459     /* Reduce the intra-simulation communication */
460     cr->sim_nodeid = cr->nodeid % nnodpersim;
461     cr->nnodes     = nnodpersim;
462 #ifdef GMX_MPI
463     MPI_Comm_split(MPI_COMM_WORLD, sim, cr->sim_nodeid, &cr->mpi_comm_mysim);
464     cr->mpi_comm_mygroup = cr->mpi_comm_mysim;
465     cr->nodeid           = cr->sim_nodeid;
466 #endif
467
468     if (debug)
469     {
470         fprintf(debug, "This is simulation %d", cr->ms->sim);
471         if (PAR(cr))
472         {
473             fprintf(debug, ", local number of nodes %d, local nodeid %d",
474                     cr->nnodes, cr->sim_nodeid);
475         }
476         fprintf(debug, "\n\n");
477     }
478
479     if (multidirs)
480     {
481         int ret;
482         if (debug)
483         {
484             fprintf(debug, "Changing to directory %s\n", multidirs[cr->ms->sim]);
485         }
486         gmx_chdir(multidirs[cr->ms->sim]);
487     }
488     else if (bParFn)
489     {
490         /* Patch output and tpx, cpt and rerun input file names */
491         for (i = 0; (i < nfile); i++)
492         {
493             /* Because of possible multiple extensions per type we must look
494              * at the actual file name
495              */
496             if (is_output(&fnm[i]) ||
497                 fnm[i].ftp == efTPX || fnm[i].ftp == efCPT ||
498                 strcmp(fnm[i].opt, "-rerun") == 0)
499             {
500                 ftp = fn2ftp(fnm[i].fns[0]);
501                 par_fn(fnm[i].fns[0], ftp, cr, TRUE, FALSE, buf, 255);
502                 sfree(fnm[i].fns[0]);
503                 fnm[i].fns[0] = gmx_strdup(buf);
504             }
505         }
506     }
507 }
508
509 t_commrec *init_par(int gmx_unused *argc, char ***argv_ptr)
510 {
511     t_commrec    *cr;
512     char        **argv;
513     int           i;
514     gmx_bool      pe = FALSE;
515
516     snew(cr, 1);
517
518     argv = argv_ptr ? *argv_ptr : NULL;
519
520 #if defined GMX_MPI && !defined GMX_THREAD_MPI
521     cr->sim_nodeid = gmx_setup(argc, argv, &cr->nnodes);
522
523     if (!PAR(cr) && (cr->sim_nodeid != 0))
524     {
525         gmx_comm("(!PAR(cr) && (cr->sim_nodeid != 0))");
526     }
527
528     cr->mpi_comm_mysim   = MPI_COMM_WORLD;
529     cr->mpi_comm_mygroup = cr->mpi_comm_mysim;
530 #else
531     /* These should never be accessed */
532     cr->mpi_comm_mysim   = NULL;
533     cr->mpi_comm_mygroup = NULL;
534     cr->nnodes           = 1;
535     cr->sim_nodeid       = 0;
536 #endif
537
538     cr->nodeid = cr->sim_nodeid;
539
540     cr->duty = (DUTY_PP | DUTY_PME);
541
542     /* Communicate arguments if parallel */
543 #ifndef GMX_THREAD_MPI
544     if (PAR(cr))
545     {
546         comm_args(cr, argc, argv_ptr);
547     }
548 #endif /* GMX_THREAD_MPI */
549
550 #ifdef GMX_MPI
551 #if !defined(GMX_THREAD_MPI) && !defined(MPI_IN_PLACE_EXISTS)
552     /* initialize the MPI_IN_PLACE replacement buffers */
553     snew(cr->mpb, 1);
554     cr->mpb->ibuf        = NULL;
555     cr->mpb->libuf       = NULL;
556     cr->mpb->fbuf        = NULL;
557     cr->mpb->dbuf        = NULL;
558     cr->mpb->ibuf_alloc  = 0;
559     cr->mpb->libuf_alloc = 0;
560     cr->mpb->fbuf_alloc  = 0;
561     cr->mpb->dbuf_alloc  = 0;
562 #endif
563 #endif
564
565     return cr;
566 }
567
568 t_commrec *init_par_threads(const t_commrec *cro)
569 {
570 #ifdef GMX_THREAD_MPI
571     int        initialized;
572     t_commrec *cr;
573
574     /* make a thread-specific commrec */
575     snew(cr, 1);
576     /* now copy the whole thing, so settings like the number of PME nodes
577        get propagated. */
578     *cr = *cro;
579
580     /* and we start setting our own thread-specific values for things */
581     MPI_Initialized(&initialized);
582     if (!initialized)
583     {
584         gmx_comm("Initializing threads without comm");
585     }
586     /* once threads will be used together with MPI, we'll
587        fill the cr structure with distinct data here. This might even work: */
588     cr->sim_nodeid = gmx_setup(0, NULL, &cr->nnodes);
589
590     cr->mpi_comm_mysim   = MPI_COMM_WORLD;
591     cr->mpi_comm_mygroup = cr->mpi_comm_mysim;
592     cr->nodeid           = cr->sim_nodeid;
593     cr->duty             = (DUTY_PP | DUTY_PME);
594
595     return cr;
596 #else
597     return NULL;
598 #endif
599 }