Revertiing my fake merge to make history consistent
[gromacs/adressmacs.git] / src / gmxlib / main.c
blob5a89e9b4d7de20d47a222e6581a8b9962f8df175
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
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.
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.
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.
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.
31 * For more info, check our website at http://www.gromacs.org
33 * And Hey:
34 * GROningen Mixture of Alchemy and Childrens' Stories
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
40 #include <stdio.h>
41 #include <stdlib.h>
42 #include <string.h>
43 #include <limits.h>
44 #include <time.h>
45 #include "smalloc.h"
46 #include "gmx_fatal.h"
47 #include "network.h"
48 #include "main.h"
49 #include "macros.h"
50 #include "futil.h"
51 #include "filenm.h"
52 #include "mdrun.h"
53 #include "gmxfio.h"
54 #include "string2.h"
56 #ifdef GMX_THREADS
57 #include "thread_mpi.h"
58 #endif
60 /* The source code in this file should be thread-safe.
61 Please keep it that way. */
64 #ifdef HAVE_UNISTD_H
65 #include <unistd.h>
66 #endif
68 #if ((defined WIN32 || defined _WIN32 || defined WIN64 || defined _WIN64) && !defined __CYGWIN__ && !defined __CYGWIN32__)
69 #include <process.h>
70 #endif
73 #define BUFSIZE 1024
75 /* this is not strictly thread-safe, but it's only written to at the beginning
76 of the simulation, once by each thread with the same value. We assume
77 that writing to an int is atomic.*/
78 static bool parallel_env_val;
79 #ifdef GMX_THREADS
80 tMPI_Thread_mutex_t parallel_env_mutex=TMPI_THREAD_MUTEX_INITIALIZER;
81 #endif
84 /* returns 1 when running in a parallel environment, so could also be 1 if
85 mdrun was started with: mpirun -np 1.
87 Use this function only to check whether a parallel environment has
88 been initialized, for example when checking whether gmx_finalize()
89 needs to be called. Use PAR(cr) to check whether the simulation actually
90 has more than one node/thread. */
91 bool gmx_parallel_env_initialized(void)
93 bool ret;
94 #ifdef GMX_THREADS
95 tMPI_Thread_mutex_lock(&parallel_env_mutex);
96 #endif
97 ret=parallel_env_val;
98 #ifdef GMX_THREADS
99 tMPI_Thread_mutex_unlock(&parallel_env_mutex);
100 #endif
101 return ret;
104 static void set_parallel_env(bool val)
106 #ifdef GMX_THREADS
107 tMPI_Thread_mutex_lock(&parallel_env_mutex);
108 #endif
109 if (!parallel_env_val)
111 /* we only allow it to be set, not unset */
112 parallel_env_val=val;
114 #ifdef GMX_THREADS
115 tMPI_Thread_mutex_unlock(&parallel_env_mutex);
116 #endif
120 static void par_fn(char *base,int ftp,const t_commrec *cr,
121 bool bAppendSimId,bool bAppendNodeId,
122 char buf[],int bufsize)
124 int n;
126 if((size_t)bufsize<(strlen(base)+10))
127 gmx_mem("Character buffer too small!");
129 /* Copy to buf, and strip extension */
130 strcpy(buf,base);
131 buf[strlen(base) - strlen(ftp2ext(fn2ftp(base))) - 1] = '\0';
133 if (bAppendSimId) {
134 sprintf(buf+strlen(buf),"%d",cr->ms->sim);
136 if (bAppendNodeId) {
137 strcat(buf,"_node");
138 sprintf(buf+strlen(buf),"%d",cr->nodeid);
140 strcat(buf,".");
142 /* Add extension again */
143 strcat(buf,(ftp == efTPX) ? "tpr" : (ftp == efEDR) ? "edr" : ftp2ext(ftp));
144 if (cr->nodeid == 0) {
145 printf("node %d par_fn '%s'\n",cr->nodeid,buf);
146 if (fn2ftp(buf) == efLOG) {
147 printf("log\n");
152 void check_multi_int(FILE *log,const gmx_multisim_t *ms,int val,
153 const char *name)
155 int *ibuf,p;
156 bool bCompatible;
158 fprintf(log,"Multi-checking %s ... ",name);
160 if (ms == NULL)
161 gmx_fatal(FARGS,
162 "check_multi_int called with a NULL communication pointer");
164 snew(ibuf,ms->nsim);
165 ibuf[ms->sim] = val;
166 gmx_sumi_sim(ms->nsim,ibuf,ms);
168 bCompatible = TRUE;
169 for(p=1; p<ms->nsim; p++)
170 bCompatible = bCompatible && (ibuf[p-1] == ibuf[p]);
172 if (bCompatible)
173 fprintf(log,"OK\n");
174 else {
175 fprintf(log,"\n%s is not equal for all subsystems\n",name);
176 for(p=0; p<ms->nsim; p++)
177 fprintf(log," subsystem %d: %d\n",p,ibuf[p]);
178 gmx_fatal(FARGS,"The %d subsystems are not compatible\n",ms->nsim);
181 sfree(ibuf);
184 void gmx_log_open(const char *lognm,const t_commrec *cr,bool bMasterOnly,
185 unsigned long Flags, FILE** fplog)
187 int len,testlen,pid;
188 char buf[256],host[256];
189 time_t t;
190 char timebuf[STRLEN];
191 FILE *fp=*fplog;
192 char *tmpnm;
194 bool bAppend = Flags & MD_APPENDFILES;
196 debug_gmx();
198 /* Communicate the filename for logfile */
199 if (cr->nnodes > 1 && !bMasterOnly
200 #ifdef GMX_THREADS
201 /* With thread MPI the non-master log files are opened later
202 * when the files names are already known on all nodes.
204 && FALSE
205 #endif
208 if (MASTER(cr))
210 len = strlen(lognm) + 1;
212 gmx_bcast(sizeof(len),&len,cr);
213 if (!MASTER(cr))
215 snew(tmpnm,len+8);
217 else
219 tmpnm=strdup(lognm);
221 gmx_bcast(len*sizeof(*tmpnm),tmpnm,cr);
223 else
225 tmpnm=strdup(lognm);
228 debug_gmx();
230 if (!bMasterOnly)
232 /* Since log always ends with '.log' let's use this info */
233 par_fn(tmpnm,efLOG,cr,FALSE,!bMasterOnly,buf,255);
234 fp = gmx_fio_fopen(buf, bAppend ? "a+" : "w+" );
236 else
238 fp = gmx_fio_fopen(tmpnm, bAppend ? "a+" : "w+" );
241 sfree(tmpnm);
243 gmx_fatal_set_log_file(fp);
245 /* Get some machine parameters */
246 #ifdef HAVE_UNISTD_H
247 if (gethostname(host,255) != 0)
249 sprintf(host,"unknown");
251 #else
252 sprintf(host,"unknown");
253 #endif
255 time(&t);
257 #ifndef NO_GETPID
258 # if ((defined WIN32 || defined _WIN32 || defined WIN64 || defined _WIN64) && !defined __CYGWIN__ && !defined __CYGWIN32__)
259 pid = _getpid();
260 # else
261 pid = getpid();
262 # endif
263 #else
264 pid = 0;
265 #endif
267 if (bAppend)
269 fprintf(fp,
270 "\n"
271 "\n"
272 "-----------------------------------------------------------\n"
273 "Restarting from checkpoint, appending to previous log file.\n"
274 "\n"
278 gmx_ctime_r(&t,timebuf,STRLEN);
280 fprintf(fp,
281 "Log file opened on %s"
282 "Host: %s pid: %d nodeid: %d nnodes: %d\n",
283 timebuf,host,pid,cr->nodeid,cr->nnodes);
285 #if (defined BUILD_MACHINE && defined BUILD_TIME && defined BUILD_USER)
286 fprintf(fp,
287 "The Gromacs distribution was built %s by\n"
288 "%s (%s)\n\n\n",BUILD_TIME,BUILD_USER,BUILD_MACHINE);
289 #endif
291 fflush(fp);
292 debug_gmx();
294 *fplog = fp;
297 void gmx_log_close(FILE *fp)
299 if (fp) {
300 gmx_fatal_set_log_file(NULL);
301 gmx_fio_fclose(fp);
305 static void comm_args(const t_commrec *cr,int *argc,char ***argv)
307 int i,len;
309 if ((cr) && PAR(cr))
310 gmx_bcast(sizeof(*argc),argc,cr);
312 if (!MASTER(cr))
313 snew(*argv,*argc+1);
314 fprintf(stderr,"NODEID=%d argc=%d\n",cr->nodeid,*argc);
315 for(i=0; (i<*argc); i++) {
316 if (MASTER(cr))
317 len = strlen((*argv)[i])+1;
318 gmx_bcast(sizeof(len),&len,cr);
319 if (!MASTER(cr))
320 snew((*argv)[i],len);
321 /*gmx_bcast(len*sizeof((*argv)[i][0]),(*argv)[i],cr);*/
322 gmx_bcast(len*sizeof(char),(*argv)[i],cr);
324 debug_gmx();
327 void init_multisystem(t_commrec *cr,int nsim, int nfile,
328 const t_filenm fnm[],bool bParFn)
330 gmx_multisim_t *ms;
331 int nnodes,nnodpersim,sim,i,ftp;
332 char buf[256];
333 #ifdef GMX_MPI
334 MPI_Group mpi_group_world;
335 #endif
336 int *rank;
338 #ifndef GMX_MPI
339 if (nsim > 1)
341 gmx_fatal(FARGS,"This binary is compiled without MPI support, can not do multiple simulations.");
343 #endif
345 nnodes = cr->nnodes;
346 if (nnodes % nsim != 0)
348 gmx_fatal(FARGS,"The number of nodes (%d) is not a multiple of the number of simulations (%d)",nnodes,nsim);
351 nnodpersim = nnodes/nsim;
352 sim = cr->nodeid/nnodpersim;
354 if (debug)
356 fprintf(debug,"We have %d simulations, %d nodes per simulation, local simulation is %d\n",nsim,nnodpersim,sim);
359 snew(ms,1);
360 cr->ms = ms;
361 ms->nsim = nsim;
362 ms->sim = sim;
363 #ifdef GMX_MPI
364 /* Create a communicator for the master nodes */
365 snew(rank,ms->nsim);
366 for(i=0; i<ms->nsim; i++)
368 rank[i] = i*nnodpersim;
370 MPI_Comm_group(MPI_COMM_WORLD,&mpi_group_world);
371 MPI_Group_incl(mpi_group_world,nsim,rank,&ms->mpi_group_masters);
372 sfree(rank);
373 MPI_Comm_create(MPI_COMM_WORLD,ms->mpi_group_masters,
374 &ms->mpi_comm_masters);
376 #if !defined(GMX_THREADS) && !defined(MPI_IN_PLACE_EXISTS)
377 /* initialize the MPI_IN_PLACE replacement buffers */
378 snew(ms->mpb, 1);
379 ms->mpb->ibuf=NULL;
380 ms->mpb->fbuf=NULL;
381 ms->mpb->dbuf=NULL;
382 ms->mpb->ibuf_alloc=0;
383 ms->mpb->fbuf_alloc=0;
384 ms->mpb->dbuf_alloc=0;
385 #endif
387 #endif
389 /* Reduce the intra-simulation communication */
390 cr->sim_nodeid = cr->nodeid % nnodpersim;
391 cr->nnodes = nnodpersim;
392 #ifdef GMX_MPI
393 MPI_Comm_split(MPI_COMM_WORLD,sim,cr->sim_nodeid,&cr->mpi_comm_mysim);
394 cr->mpi_comm_mygroup = cr->mpi_comm_mysim;
395 cr->nodeid = cr->sim_nodeid;
396 #endif
398 if (debug)
400 fprintf(debug,"This is simulation %d",cr->ms->sim);
401 if (PAR(cr))
403 fprintf(debug,", local number of nodes %d, local nodeid %d",
404 cr->nnodes,cr->sim_nodeid);
406 fprintf(debug,"\n\n");
409 if (bParFn)
411 /* Patch output and tpx, cpt and rerun input file names */
412 for(i=0; (i<nfile); i++)
414 /* Because of possible multiple extensions per type we must look
415 * at the actual file name
417 if (is_output(&fnm[i]) ||
418 fnm[i].ftp == efTPX || fnm[i].ftp == efCPT ||
419 strcmp(fnm[i].opt,"-rerun") == 0)
421 ftp = fn2ftp(fnm[i].fns[0]);
422 par_fn(fnm[i].fns[0],ftp,cr,TRUE,FALSE,buf,255);
423 sfree(fnm[i].fns[0]);
424 fnm[i].fns[0] = strdup(buf);
430 t_commrec *init_par(int *argc,char ***argv_ptr)
432 t_commrec *cr;
433 char **argv;
434 int i;
435 bool pe=FALSE;
437 snew(cr,1);
439 argv = *argv_ptr;
441 #ifdef GMX_MPI
442 #ifdef GMX_LIB_MPI
443 pe = TRUE;
444 #ifdef GMX_CHECK_MPI_ENV
445 /* Do not use MPI calls when env.var. GMX_CHECK_MPI_ENV is not set */
446 if (getenv(GMX_CHECK_MPI_ENV) == NULL)
447 pe = FALSE;
448 #endif /* GMX_CHECK_MPI_ENV */
449 #endif /* GMX_LIB_MPI */
450 set_parallel_env(pe);
451 if (pe) {
452 cr->sim_nodeid = gmx_setup(argc,argv,&cr->nnodes);
453 } else {
454 cr->nnodes = 1;
455 cr->sim_nodeid = 0;
457 #else /* GMX_MPI */
458 pe=FALSE;
459 set_parallel_env(pe);
460 cr->sim_nodeid = 0;
461 cr->nnodes = 1;
462 #endif /* GMX_MPI */
464 if (!PAR(cr) && (cr->sim_nodeid != 0))
465 gmx_comm("(!PAR(cr) && (cr->sim_nodeid != 0))");
467 if (PAR(cr))
469 #ifdef GMX_MPI
470 cr->mpi_comm_mysim = MPI_COMM_WORLD;
471 cr->mpi_comm_mygroup = cr->mpi_comm_mysim;
472 #endif /* GMX_MPI */
474 cr->nodeid = cr->sim_nodeid;
476 cr->duty = (DUTY_PP | DUTY_PME);
478 /* Communicate arguments if parallel */
479 #ifndef GMX_THREADS
480 if (PAR(cr))
481 comm_args(cr,argc,argv_ptr);
482 #endif /* GMX_THREADS */
484 #ifdef GMX_MPI
485 #if !defined(GMX_THREADS) && !defined(MPI_IN_PLACE_EXISTS)
486 /* initialize the MPI_IN_PLACE replacement buffers */
487 snew(cr->mpb, 1);
488 cr->mpb->ibuf=NULL;
489 cr->mpb->fbuf=NULL;
490 cr->mpb->dbuf=NULL;
491 cr->mpb->ibuf_alloc=0;
492 cr->mpb->fbuf_alloc=0;
493 cr->mpb->dbuf_alloc=0;
494 #endif
495 #endif
497 return cr;
500 t_commrec *init_par_threads(const t_commrec *cro)
502 #ifdef GMX_THREADS
503 int initialized;
504 t_commrec *cr;
506 /* make a thread-specific commrec */
507 snew(cr,1);
508 /* now copy the whole thing, so settings like the number of PME nodes
509 get propagated. */
510 *cr=*cro;
512 /* and we start setting our own thread-specific values for things */
513 MPI_Initialized(&initialized);
514 if (!initialized)
515 gmx_comm("Initializing threads without comm");
516 set_parallel_env(TRUE);
517 /* once threads will be used together with MPI, we'll
518 fill the cr structure with distinct data here. This might even work: */
519 cr->sim_nodeid = gmx_setup(0,NULL, &cr->nnodes);
521 cr->mpi_comm_mysim = MPI_COMM_WORLD;
522 cr->mpi_comm_mygroup = cr->mpi_comm_mysim;
523 cr->nodeid = cr->sim_nodeid;
524 cr->duty = (DUTY_PP | DUTY_PME);
526 return cr;
527 #else
528 return NULL;
529 #endif
533 t_commrec *init_cr_nopar(void)
535 t_commrec *cr;
537 snew(cr,1);
539 cr->nnodes = 1;
540 /* cr->nthreads = 1; */
541 cr->sim_nodeid = 0;
542 cr->nodeid = 0;
543 /* cr->threadid = 0; */
544 cr->duty = (DUTY_PP | DUTY_PME);
546 return cr;