1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
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
34 * GROningen Mixture of Alchemy and Childrens' Stories
39 #include "gromacs/utility/gmx_header_config.h"
47 #ifdef HAVE_SYS_TIME_H
52 #include "gmx_fatal.h"
63 #include "thread_mpi.h"
66 /* The source code in this file should be thread-safe.
67 Please keep it that way. */
74 #ifdef GMX_NATIVE_WINDOWS
79 /* Portable version of ctime_r implemented in src/gmxlib/string2.c, but we do not want it declared in public installed headers */
81 gmx_ctime_r(const time_t *clock,char *buf, int n);
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)
93 if((size_t)bufsize<(strlen(base)+10))
94 gmx_mem("Character buffer too small!");
96 /* Copy to buf, and strip extension */
98 buf[strlen(base) - strlen(ftp2ext(fn2ftp(base))) - 1] = '\0';
101 sprintf(buf+strlen(buf),"%d",cr->ms->sim);
105 sprintf(buf+strlen(buf),"%d",cr->nodeid);
109 /* Add extension again */
110 strcat(buf,(ftp == efTPX) ? "tpr" : (ftp == efEDR) ? "edr" : ftp2ext(ftp));
113 fprintf(debug, "node %d par_fn '%s'\n",cr->nodeid,buf);
114 if (fn2ftp(buf) == efLOG)
116 fprintf(debug,"log\n");
121 void check_multi_int(FILE *log,const gmx_multisim_t *ms,int val,
125 gmx_bool bCompatible;
128 fprintf(log,"Multi-checking %s ... ",name);
132 "check_multi_int called with a NULL communication pointer");
136 gmx_sumi_sim(ms->nsim,ibuf,ms);
139 for(p=1; p<ms->nsim; p++)
140 bCompatible = bCompatible && (ibuf[p-1] == ibuf[p]);
151 fprintf(log,"\n%s is not equal for all subsystems\n",name);
152 for(p=0; p<ms->nsim; p++)
153 fprintf(log," subsystem %d: %d\n",p,ibuf[p]);
155 gmx_fatal(FARGS,"The %d subsystems are not compatible\n",ms->nsim);
161 void check_multi_large_int(FILE *log,const gmx_multisim_t *ms,
162 gmx_large_int_t val, const char *name)
164 gmx_large_int_t *ibuf;
166 gmx_bool bCompatible;
169 fprintf(log,"Multi-checking %s ... ",name);
173 "check_multi_int called with a NULL communication pointer");
177 gmx_sumli_sim(ms->nsim,ibuf,ms);
180 for(p=1; p<ms->nsim; p++)
181 bCompatible = bCompatible && (ibuf[p-1] == ibuf[p]);
192 fprintf(log,"\n%s is not equal for all subsystems\n",name);
193 for(p=0; p<ms->nsim; p++)
196 /* first make the format string */
197 snprintf(strbuf, 255, " subsystem %%d: %s\n",
199 fprintf(log,strbuf,p,ibuf[p]);
202 gmx_fatal(FARGS,"The %d subsystems are not compatible\n",ms->nsim);
209 char *gmx_gethostname(char *name, size_t len)
213 gmx_incons("gmx_gethostname called with len<8");
216 if (gethostname(name, len-1) != 0)
218 strncpy(name, "unknown",8);
221 strncpy(name, "unknown",8);
228 void gmx_log_open(const char *lognm,const t_commrec *cr,gmx_bool bMasterOnly,
229 gmx_bool bAppendFiles, FILE** fplog)
232 char buf[256],host[256];
234 char timebuf[STRLEN];
240 /* Communicate the filename for logfile */
241 if (cr->nnodes > 1 && !bMasterOnly
242 #ifdef GMX_THREAD_MPI
243 /* With thread MPI the non-master log files are opened later
244 * when the files names are already known on all nodes.
252 len = strlen(lognm) + 1;
254 gmx_bcast(sizeof(len),&len,cr);
261 tmpnm=gmx_strdup(lognm);
263 gmx_bcast(len*sizeof(*tmpnm),tmpnm,cr);
267 tmpnm=gmx_strdup(lognm);
272 if (!bMasterOnly && !MASTER(cr))
274 /* Since log always ends with '.log' let's use this info */
275 par_fn(tmpnm,efLOG,cr,FALSE,!bMasterOnly,buf,255);
276 fp = gmx_fio_fopen(buf, bAppendFiles ? "a+" : "w+" );
278 else if (!bAppendFiles)
280 fp = gmx_fio_fopen(tmpnm, bAppendFiles ? "a+" : "w+" );
285 gmx_fatal_set_log_file(fp);
287 /* Get some machine parameters */
288 gmx_gethostname(host,256);
293 # ifdef GMX_NATIVE_WINDOWS
307 "-----------------------------------------------------------\n"
308 "Restarting from checkpoint, appending to previous log file.\n"
313 gmx_ctime_r(&t,timebuf,STRLEN);
316 "Log file opened on %s"
317 "Host: %s pid: %d nodeid: %d nnodes: %d\n",
318 timebuf,host,pid,cr->nodeid,cr->nnodes);
319 gmx_print_version_info(fp);
328 void gmx_log_close(FILE *fp)
331 gmx_fatal_set_log_file(NULL);
336 static void comm_args(const t_commrec *cr,int *argc,char ***argv)
341 gmx_bcast(sizeof(*argc),argc,cr);
345 fprintf(stderr,"NODEID=%d argc=%d\n",cr->nodeid,*argc);
346 for(i=0; (i<*argc); i++) {
348 len = strlen((*argv)[i])+1;
349 gmx_bcast(sizeof(len),&len,cr);
351 snew((*argv)[i],len);
352 /*gmx_bcast(len*sizeof((*argv)[i][0]),(*argv)[i],cr);*/
353 gmx_bcast(len*sizeof(char),(*argv)[i],cr);
358 void init_multisystem(t_commrec *cr,int nsim, char **multidirs,
359 int nfile, const t_filenm fnm[],gmx_bool bParFn)
362 int nnodes,nnodpersim,sim,i,ftp;
365 MPI_Group mpi_group_world;
372 gmx_fatal(FARGS,"This binary is compiled without MPI support, can not do multiple simulations.");
377 if (nnodes % nsim != 0)
379 gmx_fatal(FARGS,"The number of nodes (%d) is not a multiple of the number of simulations (%d)",nnodes,nsim);
382 nnodpersim = nnodes/nsim;
383 sim = cr->nodeid/nnodpersim;
387 fprintf(debug,"We have %d simulations, %d nodes per simulation, local simulation is %d\n",nsim,nnodpersim,sim);
395 /* Create a communicator for the master nodes */
397 for(i=0; i<ms->nsim; i++)
399 rank[i] = i*nnodpersim;
401 MPI_Comm_group(MPI_COMM_WORLD,&mpi_group_world);
402 MPI_Group_incl(mpi_group_world,nsim,rank,&ms->mpi_group_masters);
404 MPI_Comm_create(MPI_COMM_WORLD,ms->mpi_group_masters,
405 &ms->mpi_comm_masters);
407 #if !defined(GMX_THREAD_MPI) && !defined(MPI_IN_PLACE_EXISTS)
408 /* initialize the MPI_IN_PLACE replacement buffers */
414 ms->mpb->ibuf_alloc=0;
415 ms->mpb->libuf_alloc=0;
416 ms->mpb->fbuf_alloc=0;
417 ms->mpb->dbuf_alloc=0;
422 /* Reduce the intra-simulation communication */
423 cr->sim_nodeid = cr->nodeid % nnodpersim;
424 cr->nnodes = nnodpersim;
426 MPI_Comm_split(MPI_COMM_WORLD,sim,cr->sim_nodeid,&cr->mpi_comm_mysim);
427 cr->mpi_comm_mygroup = cr->mpi_comm_mysim;
428 cr->nodeid = cr->sim_nodeid;
433 fprintf(debug,"This is simulation %d",cr->ms->sim);
436 fprintf(debug,", local number of nodes %d, local nodeid %d",
437 cr->nnodes,cr->sim_nodeid);
439 fprintf(debug,"\n\n");
447 fprintf(debug,"Changing to directory %s\n",multidirs[cr->ms->sim]);
449 gmx_chdir(multidirs[cr->ms->sim]);
453 /* Patch output and tpx, cpt and rerun input file names */
454 for(i=0; (i<nfile); i++)
456 /* Because of possible multiple extensions per type we must look
457 * at the actual file name
459 if (is_output(&fnm[i]) ||
460 fnm[i].ftp == efTPX || fnm[i].ftp == efCPT ||
461 strcmp(fnm[i].opt,"-rerun") == 0)
463 ftp = fn2ftp(fnm[i].fns[0]);
464 par_fn(fnm[i].fns[0],ftp,cr,TRUE,FALSE,buf,255);
465 sfree(fnm[i].fns[0]);
466 fnm[i].fns[0] = gmx_strdup(buf);
472 t_commrec *init_par(int *argc,char ***argv_ptr)
481 argv = argv_ptr ? *argv_ptr : NULL;
483 #if defined GMX_MPI && !defined GMX_THREAD_MPI
484 cr->sim_nodeid = gmx_setup(argc,argv,&cr->nnodes);
486 if (!PAR(cr) && (cr->sim_nodeid != 0))
488 gmx_comm("(!PAR(cr) && (cr->sim_nodeid != 0))");
491 cr->mpi_comm_mysim = MPI_COMM_WORLD;
492 cr->mpi_comm_mygroup = cr->mpi_comm_mysim;
494 /* These should never be accessed */
495 cr->mpi_comm_mysim = NULL;
496 cr->mpi_comm_mygroup = NULL;
501 cr->nodeid = cr->sim_nodeid;
503 cr->duty = (DUTY_PP | DUTY_PME);
505 /* Communicate arguments if parallel */
506 #ifndef GMX_THREAD_MPI
509 comm_args(cr,argc,argv_ptr);
511 #endif /* GMX_THREAD_MPI */
514 #if !defined(GMX_THREAD_MPI) && !defined(MPI_IN_PLACE_EXISTS)
515 /* initialize the MPI_IN_PLACE replacement buffers */
521 cr->mpb->ibuf_alloc=0;
522 cr->mpb->libuf_alloc=0;
523 cr->mpb->fbuf_alloc=0;
524 cr->mpb->dbuf_alloc=0;
531 t_commrec *init_par_threads(const t_commrec *cro)
533 #ifdef GMX_THREAD_MPI
537 /* make a thread-specific commrec */
539 /* now copy the whole thing, so settings like the number of PME nodes
543 /* and we start setting our own thread-specific values for things */
544 MPI_Initialized(&initialized);
547 gmx_comm("Initializing threads without comm");
549 /* once threads will be used together with MPI, we'll
550 fill the cr structure with distinct data here. This might even work: */
551 cr->sim_nodeid = gmx_setup(0,NULL, &cr->nnodes);
553 cr->mpi_comm_mysim = MPI_COMM_WORLD;
554 cr->mpi_comm_mygroup = cr->mpi_comm_mysim;
555 cr->nodeid = cr->sim_nodeid;
556 cr->duty = (DUTY_PP | DUTY_PME);