28644ba909b177380e0626b35ef3743a09968b46
[alexxy/gromacs.git] / src / gromacs / gmxlib / main.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #include "gmxpre.h"
38
39 #include "gromacs/legacyheaders/main.h"
40
41 #include "config.h"
42
43 #include <limits.h>
44 #include <stdio.h>
45 #include <stdlib.h>
46 #include <string.h>
47
48 #include "gromacs/fileio/filenm.h"
49 #include "gromacs/fileio/gmxfio.h"
50 #include "gromacs/legacyheaders/copyrite.h"
51 #include "gromacs/legacyheaders/macros.h"
52 #include "gromacs/legacyheaders/network.h"
53 #include "gromacs/legacyheaders/types/commrec.h"
54 #include "gromacs/utility/cstringutil.h"
55 #include "gromacs/utility/exceptions.h"
56 #include "gromacs/utility/fatalerror.h"
57 #include "gromacs/utility/futil.h"
58 #include "gromacs/utility/gmxmpi.h"
59 #include "gromacs/utility/programcontext.h"
60 #include "gromacs/utility/smalloc.h"
61 #include "gromacs/utility/snprintf.h"
62 #include "gromacs/utility/sysinfo.h"
63
64 /* The source code in this file should be thread-safe.
65          Please keep it that way. */
66
67 #define BUFSIZE 1024
68
69 static void par_fn(char *base, int ftp, const t_commrec *cr,
70                    gmx_bool bAppendSimId, gmx_bool bAppendNodeId,
71                    char buf[], int bufsize)
72 {
73     if ((size_t)bufsize < (strlen(base)+10))
74     {
75         gmx_mem("Character buffer too small!");
76     }
77
78     /* Copy to buf, and strip extension */
79     strcpy(buf, base);
80     buf[strlen(base) - strlen(ftp2ext(fn2ftp(base))) - 1] = '\0';
81
82     if (bAppendSimId)
83     {
84         sprintf(buf+strlen(buf), "%d", cr->ms->sim);
85     }
86     if (bAppendNodeId)
87     {
88         strcat(buf, "_rank");
89         sprintf(buf+strlen(buf), "%d", cr->nodeid);
90     }
91     strcat(buf, ".");
92
93     /* Add extension again */
94     strcat(buf, (ftp == efTPR) ? "tpr" : (ftp == efEDR) ? "edr" : ftp2ext(ftp));
95     if (debug)
96     {
97         fprintf(debug, "rank %d par_fn '%s'\n", cr->nodeid, buf);
98         if (fn2ftp(buf) == efLOG)
99         {
100             fprintf(debug, "log\n");
101         }
102     }
103 }
104
105 void check_multi_int(FILE *log, const gmx_multisim_t *ms, int val,
106                      const char *name,
107                      gmx_bool bQuiet)
108 {
109     int     *ibuf, p;
110     gmx_bool bCompatible;
111
112     if (NULL != log && !bQuiet)
113     {
114         fprintf(log, "Multi-checking %s ... ", name);
115     }
116
117     if (ms == NULL)
118     {
119         gmx_fatal(FARGS,
120                   "check_multi_int called with a NULL communication pointer");
121     }
122
123     snew(ibuf, ms->nsim);
124     ibuf[ms->sim] = val;
125     gmx_sumi_sim(ms->nsim, ibuf, ms);
126
127     bCompatible = TRUE;
128     for (p = 1; p < ms->nsim; p++)
129     {
130         bCompatible = bCompatible && (ibuf[p-1] == ibuf[p]);
131     }
132
133     if (bCompatible)
134     {
135         if (NULL != log && !bQuiet)
136         {
137             fprintf(log, "OK\n");
138         }
139     }
140     else
141     {
142         if (NULL != log)
143         {
144             fprintf(log, "\n%s is not equal for all subsystems\n", name);
145             for (p = 0; p < ms->nsim; p++)
146             {
147                 fprintf(log, "  subsystem %d: %d\n", p, ibuf[p]);
148             }
149         }
150         gmx_fatal(FARGS, "The %d subsystems are not compatible\n", ms->nsim);
151     }
152
153     sfree(ibuf);
154 }
155
156 void check_multi_int64(FILE *log, const gmx_multisim_t *ms,
157                        gmx_int64_t val, const char *name,
158                        gmx_bool bQuiet)
159 {
160     gmx_int64_t      *ibuf;
161     int               p;
162     gmx_bool          bCompatible;
163
164     if (NULL != log && !bQuiet)
165     {
166         fprintf(log, "Multi-checking %s ... ", name);
167     }
168
169     if (ms == NULL)
170     {
171         gmx_fatal(FARGS,
172                   "check_multi_int called with a NULL communication pointer");
173     }
174
175     snew(ibuf, ms->nsim);
176     ibuf[ms->sim] = val;
177     gmx_sumli_sim(ms->nsim, ibuf, ms);
178
179     bCompatible = TRUE;
180     for (p = 1; p < ms->nsim; p++)
181     {
182         bCompatible = bCompatible && (ibuf[p-1] == ibuf[p]);
183     }
184
185     if (bCompatible)
186     {
187         if (NULL != log && !bQuiet)
188         {
189             fprintf(log, "OK\n");
190         }
191     }
192     else
193     {
194         if (NULL != log)
195         {
196             fprintf(log, "\n%s is not equal for all subsystems\n", name);
197             for (p = 0; p < ms->nsim; p++)
198             {
199                 char strbuf[255];
200                 /* first make the format string */
201                 snprintf(strbuf, 255, "  subsystem %%d: %s\n",
202                          "%" GMX_PRId64);
203                 fprintf(log, strbuf, p, ibuf[p]);
204             }
205         }
206         gmx_fatal(FARGS, "The %d subsystems are not compatible\n", ms->nsim);
207     }
208
209     sfree(ibuf);
210 }
211
212
213 void gmx_log_open(const char *lognm, const t_commrec *cr,
214                   gmx_bool bAppendFiles, FILE** fplog)
215 {
216     int    pid;
217     char   host[256];
218     char   timebuf[STRLEN];
219     FILE  *fp = *fplog;
220
221     debug_gmx();
222
223     if (!bAppendFiles)
224     {
225         fp = gmx_fio_fopen(lognm, bAppendFiles ? "a+" : "w+" );
226     }
227
228     gmx_fatal_set_log_file(fp);
229
230     /* Get some machine parameters */
231     gmx_gethostname(host, 256);
232     pid = gmx_getpid();
233     gmx_format_current_time(timebuf, STRLEN);
234
235     if (bAppendFiles)
236     {
237         fprintf(fp,
238                 "\n"
239                 "\n"
240                 "-----------------------------------------------------------\n"
241                 "Restarting from checkpoint, appending to previous log file.\n"
242                 "\n"
243                 );
244     }
245
246     fprintf(fp,
247             "Log file opened on %s"
248             "Host: %s  pid: %d  rank ID: %d  number of ranks:  %d\n",
249             timebuf, host, pid, cr->nodeid, cr->nnodes);
250     try
251     {
252         gmx::BinaryInformationSettings settings;
253         settings.extendedInfo(true);
254         settings.copyright(!bAppendFiles);
255         gmx::printBinaryInformation(fp, gmx::getProgramContext(), settings);
256     }
257     GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
258     fprintf(fp, "\n");
259
260     fflush(fp);
261     debug_gmx();
262
263     *fplog = fp;
264 }
265
266 void gmx_log_close(FILE *fp)
267 {
268     if (fp)
269     {
270         gmx_fatal_set_log_file(NULL);
271         gmx_fio_fclose(fp);
272     }
273 }
274
275 void init_multisystem(t_commrec *cr, int nsim, char **multidirs,
276                       int nfile, const t_filenm fnm[], gmx_bool bParFn)
277 {
278     gmx_multisim_t *ms;
279     int             nnodes, nnodpersim, sim, i, ftp;
280     char            buf[256];
281 #ifdef GMX_MPI
282     MPI_Group       mpi_group_world;
283     int            *rank;
284 #endif
285
286 #ifndef GMX_MPI
287     if (nsim > 1)
288     {
289         gmx_fatal(FARGS, "This binary is compiled without MPI support, can not do multiple simulations.");
290     }
291 #endif
292
293     nnodes  = cr->nnodes;
294     if (nnodes % nsim != 0)
295     {
296         gmx_fatal(FARGS, "The number of ranks (%d) is not a multiple of the number of simulations (%d)", nnodes, nsim);
297     }
298
299     nnodpersim = nnodes/nsim;
300     sim        = cr->nodeid/nnodpersim;
301
302     if (debug)
303     {
304         fprintf(debug, "We have %d simulations, %d ranks per simulation, local simulation is %d\n", nsim, nnodpersim, sim);
305     }
306
307     snew(ms, 1);
308     cr->ms   = ms;
309     ms->nsim = nsim;
310     ms->sim  = sim;
311 #ifdef GMX_MPI
312     /* Create a communicator for the master nodes */
313     snew(rank, ms->nsim);
314     for (i = 0; i < ms->nsim; i++)
315     {
316         rank[i] = i*nnodpersim;
317     }
318     MPI_Comm_group(MPI_COMM_WORLD, &mpi_group_world);
319     MPI_Group_incl(mpi_group_world, nsim, rank, &ms->mpi_group_masters);
320     sfree(rank);
321     MPI_Comm_create(MPI_COMM_WORLD, ms->mpi_group_masters,
322                     &ms->mpi_comm_masters);
323
324 #if !defined(MPI_IN_PLACE_EXISTS)
325     /* initialize the MPI_IN_PLACE replacement buffers */
326     snew(ms->mpb, 1);
327     ms->mpb->ibuf        = NULL;
328     ms->mpb->libuf       = NULL;
329     ms->mpb->fbuf        = NULL;
330     ms->mpb->dbuf        = NULL;
331     ms->mpb->ibuf_alloc  = 0;
332     ms->mpb->libuf_alloc = 0;
333     ms->mpb->fbuf_alloc  = 0;
334     ms->mpb->dbuf_alloc  = 0;
335 #endif
336
337 #endif
338
339     /* Reduce the intra-simulation communication */
340     cr->sim_nodeid = cr->nodeid % nnodpersim;
341     cr->nnodes     = nnodpersim;
342 #ifdef GMX_MPI
343     MPI_Comm_split(MPI_COMM_WORLD, sim, cr->sim_nodeid, &cr->mpi_comm_mysim);
344     cr->mpi_comm_mygroup = cr->mpi_comm_mysim;
345     cr->nodeid           = cr->sim_nodeid;
346 #endif
347
348     if (debug)
349     {
350         fprintf(debug, "This is simulation %d", cr->ms->sim);
351         if (PAR(cr))
352         {
353             fprintf(debug, ", local number of ranks %d, local rank ID %d",
354                     cr->nnodes, cr->sim_nodeid);
355         }
356         fprintf(debug, "\n\n");
357     }
358
359     if (multidirs)
360     {
361         if (debug)
362         {
363             fprintf(debug, "Changing to directory %s\n", multidirs[cr->ms->sim]);
364         }
365         gmx_chdir(multidirs[cr->ms->sim]);
366     }
367     else if (bParFn)
368     {
369         /* Patch output and tpx, cpt and rerun input file names */
370         for (i = 0; (i < nfile); i++)
371         {
372             /* Because of possible multiple extensions per type we must look
373              * at the actual file name
374              */
375             if (is_output(&fnm[i]) ||
376                 fnm[i].ftp == efTPR || fnm[i].ftp == efCPT ||
377                 strcmp(fnm[i].opt, "-rerun") == 0)
378             {
379                 ftp = fn2ftp(fnm[i].fns[0]);
380                 par_fn(fnm[i].fns[0], ftp, cr, TRUE, FALSE, buf, 255);
381                 sfree(fnm[i].fns[0]);
382                 fnm[i].fns[0] = gmx_strdup(buf);
383             }
384         }
385     }
386 }