Uniform "generated by" information in output files.
[alexxy/gromacs.git] / src / gromacs / gmxlib / main.cpp
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 #include "gromacs/utility/exceptions.h"
67 #include "gromacs/utility/programinfo.h"
68
69 /* The source code in this file should be thread-safe.
70          Please keep it that way. */
71
72
73 #ifdef HAVE_UNISTD_H
74 #include <unistd.h>
75 #endif
76
77 #ifdef GMX_NATIVE_WINDOWS
78 #include <process.h>
79 #endif
80
81 #define BUFSIZE 1024
82
83
84 static void par_fn(char *base, int ftp, const t_commrec *cr,
85                    gmx_bool bAppendSimId, gmx_bool bAppendNodeId,
86                    char buf[], int bufsize)
87 {
88     if ((size_t)bufsize < (strlen(base)+10))
89     {
90         gmx_mem("Character buffer too small!");
91     }
92
93     /* Copy to buf, and strip extension */
94     strcpy(buf, base);
95     buf[strlen(base) - strlen(ftp2ext(fn2ftp(base))) - 1] = '\0';
96
97     if (bAppendSimId)
98     {
99         sprintf(buf+strlen(buf), "%d", cr->ms->sim);
100     }
101     if (bAppendNodeId)
102     {
103         strcat(buf, "_node");
104         sprintf(buf+strlen(buf), "%d", cr->nodeid);
105     }
106     strcat(buf, ".");
107
108     /* Add extension again */
109     strcat(buf, (ftp == efTPX) ? "tpr" : (ftp == efEDR) ? "edr" : ftp2ext(ftp));
110     if (debug)
111     {
112         fprintf(debug, "node %d par_fn '%s'\n", cr->nodeid, buf);
113         if (fn2ftp(buf) == efLOG)
114         {
115             fprintf(debug, "log\n");
116         }
117     }
118 }
119
120 void check_multi_int(FILE *log, const gmx_multisim_t *ms, int val,
121                      const char *name,
122                      gmx_bool bQuiet)
123 {
124     int     *ibuf, p;
125     gmx_bool bCompatible;
126
127     if (NULL != log && !bQuiet)
128     {
129         fprintf(log, "Multi-checking %s ... ", name);
130     }
131
132     if (ms == NULL)
133     {
134         gmx_fatal(FARGS,
135                   "check_multi_int called with a NULL communication pointer");
136     }
137
138     snew(ibuf, ms->nsim);
139     ibuf[ms->sim] = val;
140     gmx_sumi_sim(ms->nsim, ibuf, ms);
141
142     bCompatible = TRUE;
143     for (p = 1; p < ms->nsim; p++)
144     {
145         bCompatible = bCompatible && (ibuf[p-1] == ibuf[p]);
146     }
147
148     if (bCompatible)
149     {
150         if (NULL != log && !bQuiet)
151         {
152             fprintf(log, "OK\n");
153         }
154     }
155     else
156     {
157         if (NULL != log)
158         {
159             fprintf(log, "\n%s is not equal for all subsystems\n", name);
160             for (p = 0; p < ms->nsim; p++)
161             {
162                 fprintf(log, "  subsystem %d: %d\n", p, ibuf[p]);
163             }
164         }
165         gmx_fatal(FARGS, "The %d subsystems are not compatible\n", ms->nsim);
166     }
167
168     sfree(ibuf);
169 }
170
171 void check_multi_large_int(FILE *log, const gmx_multisim_t *ms,
172                            gmx_large_int_t val, const char *name,
173                            gmx_bool bQuiet)
174 {
175     gmx_large_int_t  *ibuf;
176     int               p;
177     gmx_bool          bCompatible;
178
179     if (NULL != log && !bQuiet)
180     {
181         fprintf(log, "Multi-checking %s ... ", name);
182     }
183
184     if (ms == NULL)
185     {
186         gmx_fatal(FARGS,
187                   "check_multi_int called with a NULL communication pointer");
188     }
189
190     snew(ibuf, ms->nsim);
191     ibuf[ms->sim] = val;
192     gmx_sumli_sim(ms->nsim, ibuf, ms);
193
194     bCompatible = TRUE;
195     for (p = 1; p < ms->nsim; p++)
196     {
197         bCompatible = bCompatible && (ibuf[p-1] == ibuf[p]);
198     }
199
200     if (bCompatible)
201     {
202         if (NULL != log && !bQuiet)
203         {
204             fprintf(log, "OK\n");
205         }
206     }
207     else
208     {
209         if (NULL != log)
210         {
211             fprintf(log, "\n%s is not equal for all subsystems\n", name);
212             for (p = 0; p < ms->nsim; p++)
213             {
214                 char strbuf[255];
215                 /* first make the format string */
216                 snprintf(strbuf, 255, "  subsystem %%d: %s\n",
217                          gmx_large_int_pfmt);
218                 fprintf(log, strbuf, p, ibuf[p]);
219             }
220         }
221         gmx_fatal(FARGS, "The %d subsystems are not compatible\n", ms->nsim);
222     }
223
224     sfree(ibuf);
225 }
226
227
228 char *gmx_gethostname(char *name, size_t len)
229 {
230     if (len < 8)
231     {
232         gmx_incons("gmx_gethostname called with len<8");
233     }
234 #ifdef HAVE_UNISTD_H
235     if (gethostname(name, len-1) != 0)
236     {
237         strncpy(name, "unknown", 8);
238     }
239 #else
240     strncpy(name, "unknown", 8);
241 #endif
242
243     return name;
244 }
245
246
247 void gmx_log_open(const char *lognm, const t_commrec *cr, gmx_bool bMasterOnly,
248                   gmx_bool bAppendFiles, FILE** fplog)
249 {
250     int    len, pid;
251     char   buf[256], host[256];
252     time_t t;
253     char   timebuf[STRLEN];
254     FILE  *fp = *fplog;
255     char  *tmpnm;
256
257     debug_gmx();
258
259     /* Communicate the filename for logfile */
260     if (cr->nnodes > 1 && !bMasterOnly
261 #ifdef GMX_THREAD_MPI
262         /* With thread MPI the non-master log files are opened later
263          * when the files names are already known on all nodes.
264          */
265         && FALSE
266 #endif
267         )
268     {
269         if (MASTER(cr))
270         {
271             len = strlen(lognm) + 1;
272         }
273         gmx_bcast(sizeof(len), &len, cr);
274         if (!MASTER(cr))
275         {
276             snew(tmpnm, len+8);
277         }
278         else
279         {
280             tmpnm = gmx_strdup(lognm);
281         }
282         gmx_bcast(len*sizeof(*tmpnm), tmpnm, cr);
283     }
284     else
285     {
286         tmpnm = gmx_strdup(lognm);
287     }
288
289     debug_gmx();
290
291     if (!bMasterOnly && !MASTER(cr))
292     {
293         /* Since log always ends with '.log' let's use this info */
294         par_fn(tmpnm, efLOG, cr, FALSE, !bMasterOnly, buf, 255);
295         fp = gmx_fio_fopen(buf, bAppendFiles ? "a+" : "w+" );
296     }
297     else if (!bAppendFiles)
298     {
299         fp = gmx_fio_fopen(tmpnm, bAppendFiles ? "a+" : "w+" );
300     }
301
302     sfree(tmpnm);
303
304     gmx_fatal_set_log_file(fp);
305
306     /* Get some machine parameters */
307     gmx_gethostname(host, 256);
308
309     time(&t);
310
311 #ifndef NO_GETPID
312 #   ifdef GMX_NATIVE_WINDOWS
313     pid = _getpid();
314 #   else
315     pid = getpid();
316 #   endif
317 #else
318     pid = 0;
319 #endif
320
321     if (bAppendFiles)
322     {
323         fprintf(fp,
324                 "\n"
325                 "\n"
326                 "-----------------------------------------------------------\n"
327                 "Restarting from checkpoint, appending to previous log file.\n"
328                 "\n"
329                 );
330     }
331
332     gmx_ctime_r(&t, timebuf, STRLEN);
333
334     fprintf(fp,
335             "Log file opened on %s"
336             "Host: %s  pid: %d  nodeid: %d  nnodes:  %d\n",
337             timebuf, host, pid, cr->nodeid, cr->nnodes);
338     try
339     {
340         gmx::BinaryInformationSettings settings;
341         settings.extendedInfo(true);
342         settings.copyright(!bAppendFiles);
343         gmx::printBinaryInformation(fp, gmx::ProgramInfo::getInstance(), settings);
344     }
345     GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
346     fprintf(fp, "\n\n");
347
348     fflush(fp);
349     debug_gmx();
350
351     *fplog = fp;
352 }
353
354 void gmx_log_close(FILE *fp)
355 {
356     if (fp)
357     {
358         gmx_fatal_set_log_file(NULL);
359         gmx_fio_fclose(fp);
360     }
361 }
362
363 static void comm_args(const t_commrec *cr, int *argc, char ***argv)
364 {
365     int i, len;
366
367     if (PAR(cr))
368     {
369         gmx_bcast(sizeof(*argc), argc, cr);
370     }
371
372     if (!MASTER(cr))
373     {
374         snew(*argv, *argc+1);
375     }
376     if (debug)
377     {
378         fprintf(debug, "NODEID=%d argc=%d\n", cr->nodeid, *argc);
379     }
380     for (i = 0; (i < *argc); i++)
381     {
382         if (MASTER(cr))
383         {
384             len = strlen((*argv)[i])+1;
385         }
386         gmx_bcast(sizeof(len), &len, cr);
387         if (!MASTER(cr))
388         {
389             snew((*argv)[i], len);
390         }
391         /*gmx_bcast(len*sizeof((*argv)[i][0]),(*argv)[i],cr);*/
392         gmx_bcast(len*sizeof(char), (*argv)[i], cr);
393     }
394     debug_gmx();
395 }
396
397 void init_multisystem(t_commrec *cr, int nsim, char **multidirs,
398                       int nfile, const t_filenm fnm[], gmx_bool bParFn)
399 {
400     gmx_multisim_t *ms;
401     int             nnodes, nnodpersim, sim, i, ftp;
402     char            buf[256];
403 #ifdef GMX_MPI
404     MPI_Group       mpi_group_world;
405     int            *rank;
406 #endif
407
408 #ifndef GMX_MPI
409     if (nsim > 1)
410     {
411         gmx_fatal(FARGS, "This binary is compiled without MPI support, can not do multiple simulations.");
412     }
413 #endif
414
415     nnodes  = cr->nnodes;
416     if (nnodes % nsim != 0)
417     {
418         gmx_fatal(FARGS, "The number of nodes (%d) is not a multiple of the number of simulations (%d)", nnodes, nsim);
419     }
420
421     nnodpersim = nnodes/nsim;
422     sim        = cr->nodeid/nnodpersim;
423
424     if (debug)
425     {
426         fprintf(debug, "We have %d simulations, %d nodes per simulation, local simulation is %d\n", nsim, nnodpersim, sim);
427     }
428
429     snew(ms, 1);
430     cr->ms   = ms;
431     ms->nsim = nsim;
432     ms->sim  = sim;
433 #ifdef GMX_MPI
434     /* Create a communicator for the master nodes */
435     snew(rank, ms->nsim);
436     for (i = 0; i < ms->nsim; i++)
437     {
438         rank[i] = i*nnodpersim;
439     }
440     MPI_Comm_group(MPI_COMM_WORLD, &mpi_group_world);
441     MPI_Group_incl(mpi_group_world, nsim, rank, &ms->mpi_group_masters);
442     sfree(rank);
443     MPI_Comm_create(MPI_COMM_WORLD, ms->mpi_group_masters,
444                     &ms->mpi_comm_masters);
445
446 #if !defined(GMX_THREAD_MPI) && !defined(MPI_IN_PLACE_EXISTS)
447     /* initialize the MPI_IN_PLACE replacement buffers */
448     snew(ms->mpb, 1);
449     ms->mpb->ibuf        = NULL;
450     ms->mpb->libuf       = NULL;
451     ms->mpb->fbuf        = NULL;
452     ms->mpb->dbuf        = NULL;
453     ms->mpb->ibuf_alloc  = 0;
454     ms->mpb->libuf_alloc = 0;
455     ms->mpb->fbuf_alloc  = 0;
456     ms->mpb->dbuf_alloc  = 0;
457 #endif
458
459 #endif
460
461     /* Reduce the intra-simulation communication */
462     cr->sim_nodeid = cr->nodeid % nnodpersim;
463     cr->nnodes     = nnodpersim;
464 #ifdef GMX_MPI
465     MPI_Comm_split(MPI_COMM_WORLD, sim, cr->sim_nodeid, &cr->mpi_comm_mysim);
466     cr->mpi_comm_mygroup = cr->mpi_comm_mysim;
467     cr->nodeid           = cr->sim_nodeid;
468 #endif
469
470     if (debug)
471     {
472         fprintf(debug, "This is simulation %d", cr->ms->sim);
473         if (PAR(cr))
474         {
475             fprintf(debug, ", local number of nodes %d, local nodeid %d",
476                     cr->nnodes, cr->sim_nodeid);
477         }
478         fprintf(debug, "\n\n");
479     }
480
481     if (multidirs)
482     {
483         if (debug)
484         {
485             fprintf(debug, "Changing to directory %s\n", multidirs[cr->ms->sim]);
486         }
487         gmx_chdir(multidirs[cr->ms->sim]);
488     }
489     else if (bParFn)
490     {
491         /* Patch output and tpx, cpt and rerun input file names */
492         for (i = 0; (i < nfile); i++)
493         {
494             /* Because of possible multiple extensions per type we must look
495              * at the actual file name
496              */
497             if (is_output(&fnm[i]) ||
498                 fnm[i].ftp == efTPX || fnm[i].ftp == efCPT ||
499                 strcmp(fnm[i].opt, "-rerun") == 0)
500             {
501                 ftp = fn2ftp(fnm[i].fns[0]);
502                 par_fn(fnm[i].fns[0], ftp, cr, TRUE, FALSE, buf, 255);
503                 sfree(fnm[i].fns[0]);
504                 fnm[i].fns[0] = gmx_strdup(buf);
505             }
506         }
507     }
508 }
509
510 t_commrec *init_par(int gmx_unused *argc, char ***argv_ptr)
511 {
512     t_commrec    *cr;
513
514     snew(cr, 1);
515
516 #if defined GMX_MPI && !defined GMX_THREAD_MPI
517     char **argv = argv_ptr ? *argv_ptr : NULL;
518     cr->sim_nodeid = gmx_setup(argc, argv, &cr->nnodes);
519
520     if (!PAR(cr) && (cr->sim_nodeid != 0))
521     {
522         gmx_comm("(!PAR(cr) && (cr->sim_nodeid != 0))");
523     }
524
525     cr->mpi_comm_mysim   = MPI_COMM_WORLD;
526     cr->mpi_comm_mygroup = cr->mpi_comm_mysim;
527 #else
528     /* These should never be accessed */
529     cr->mpi_comm_mysim   = NULL;
530     cr->mpi_comm_mygroup = NULL;
531     cr->nnodes           = 1;
532     cr->sim_nodeid       = 0;
533 #endif
534
535     cr->nodeid = cr->sim_nodeid;
536
537     cr->duty = (DUTY_PP | DUTY_PME);
538
539     /* Communicate arguments if parallel */
540 #ifndef GMX_THREAD_MPI
541     if (PAR(cr))
542     {
543         comm_args(cr, argc, argv_ptr);
544     }
545 #endif /* GMX_THREAD_MPI */
546
547 #ifdef GMX_MPI
548 #if !defined(GMX_THREAD_MPI) && !defined(MPI_IN_PLACE_EXISTS)
549     /* initialize the MPI_IN_PLACE replacement buffers */
550     snew(cr->mpb, 1);
551     cr->mpb->ibuf        = NULL;
552     cr->mpb->libuf       = NULL;
553     cr->mpb->fbuf        = NULL;
554     cr->mpb->dbuf        = NULL;
555     cr->mpb->ibuf_alloc  = 0;
556     cr->mpb->libuf_alloc = 0;
557     cr->mpb->fbuf_alloc  = 0;
558     cr->mpb->dbuf_alloc  = 0;
559 #endif
560 #endif
561
562     return cr;
563 }
564
565 t_commrec *init_par_threads(const t_commrec *cro)
566 {
567 #ifdef GMX_THREAD_MPI
568     int        initialized;
569     t_commrec *cr;
570
571     /* make a thread-specific commrec */
572     snew(cr, 1);
573     /* now copy the whole thing, so settings like the number of PME nodes
574        get propagated. */
575     *cr = *cro;
576
577     /* and we start setting our own thread-specific values for things */
578     MPI_Initialized(&initialized);
579     if (!initialized)
580     {
581         gmx_comm("Initializing threads without comm");
582     }
583     /* once threads will be used together with MPI, we'll
584        fill the cr structure with distinct data here. This might even work: */
585     cr->sim_nodeid = gmx_setup(0, NULL, &cr->nnodes);
586
587     cr->mpi_comm_mysim   = MPI_COMM_WORLD;
588     cr->mpi_comm_mygroup = cr->mpi_comm_mysim;
589     cr->nodeid           = cr->sim_nodeid;
590     cr->duty             = (DUTY_PP | DUTY_PME);
591
592     return cr;
593 #else
594     return NULL;
595 #endif
596 }