minor fixes in ditribution files
[gromacs/qmmm-gamess-us.git] / src / gmxlib / main.c
blob2a5b318346daa60827162f68d75837524510a694
1 /*
2 *
3 * This source code is part of
4 *
5 * G R O M A C S
6 *
7 * GROningen MAchine for Chemical Simulations
8 *
9 * VERSION 3.2.0
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
32 * And Hey:
33 * GROningen Mixture of Alchemy and Childrens' Stories
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
39 #include <stdio.h>
40 #include <stdlib.h>
41 #include <string.h>
42 #include <limits.h>
43 #include <time.h>
44 #include "smalloc.h"
45 #include "gmx_fatal.h"
46 #include "network.h"
47 #include "main.h"
48 #include "macros.h"
49 #include "futil.h"
50 #include "filenm.h"
51 #include "mdrun.h"
52 #include "gmxfio.h"
54 #ifdef GMX_THREADS
55 #include "thread_mpi.h"
56 #endif
58 /* The source code in this file should be thread-safe.
59 Please keep it that way. */
62 #ifdef HAVE_UNISTD_H
63 #include <unistd.h>
64 #endif
66 #if ((defined WIN32 || defined _WIN32 || defined WIN64 || defined _WIN64) && !defined __CYGWIN__ && !defined __CYGWIN32__)
67 #include <process.h>
68 #endif
71 #define BUFSIZE 1024
73 /* this is not strictly thread-safe, but it's only written to at the beginning
74 of the simulation, once by each thread with the same value. We assume
75 that writing to an int is atomic.*/
76 static bool parallel_env_val;
77 #ifdef GMX_THREADS
78 tMPI_Thread_mutex_t parallel_env_mutex=TMPI_THREAD_MUTEX_INITIALIZER;
79 #endif
82 /* returns 1 when running in a parallel environment, so could also be 1 if
83 mdrun was started with: mpirun -np 1.
85 Use this function only to check whether a parallel environment has
86 been initialized, for example when checking whether gmx_finalize()
87 needs to be called. Use PAR(cr) to check whether the simulation actually
88 has more than one node/thread. */
89 bool gmx_parallel_env_initialized(void)
91 bool ret;
92 #ifdef GMX_THREADS
93 tMPI_Thread_mutex_lock(&parallel_env_mutex);
94 #endif
95 ret=parallel_env_val;
96 #ifdef GMX_THREADS
97 tMPI_Thread_mutex_unlock(&parallel_env_mutex);
98 #endif
99 return ret;
102 static void set_parallel_env(bool val)
104 #ifdef GMX_THREADS
105 tMPI_Thread_mutex_lock(&parallel_env_mutex);
106 #endif
107 if (!parallel_env_val)
109 /* we only allow it to be set, not unset */
110 parallel_env_val=val;
112 #ifdef GMX_THREADS
113 tMPI_Thread_mutex_unlock(&parallel_env_mutex);
114 #endif
119 static void par_fn(char *base,int ftp,const t_commrec *cr,
120 bool bUnderScore,
121 char buf[],int bufsize)
123 int n;
125 if((size_t)bufsize<(strlen(base)+4))
126 gmx_mem("Character buffer too small!");
128 /* Copy to buf, and strip extension */
129 strcpy(buf,base);
130 buf[strlen(base) - strlen(ftp2ext(fn2ftp(base))) - 1] = '\0';
132 /* Add node info */
133 if (bUnderScore)
134 strcat(buf,"_");
135 if (MULTISIM(cr) && !bUnderScore) {
136 sprintf(buf+strlen(buf),"%d",cr->ms->sim);
137 } else if (PAR(cr)) {
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));
146 void check_multi_int(FILE *log,const gmx_multisim_t *ms,int val,
147 const char *name)
149 int *ibuf,p;
150 bool bCompatible;
152 fprintf(log,"Multi-checking %s ... ",name);
154 if (ms == NULL)
155 gmx_fatal(FARGS,
156 "check_multi_int called with a NULL communication pointer");
158 snew(ibuf,ms->nsim);
159 ibuf[ms->sim] = val;
160 gmx_sumi_sim(ms->nsim,ibuf,ms);
162 bCompatible = TRUE;
163 for(p=1; p<ms->nsim; p++)
164 bCompatible = bCompatible && (ibuf[p-1] == ibuf[p]);
166 if (bCompatible)
167 fprintf(log,"OK\n");
168 else {
169 fprintf(log,"\n%s is not equal for all subsystems\n",name);
170 for(p=0; p<ms->nsim; p++)
171 fprintf(log," subsystem %d: %d\n",p,ibuf[p]);
172 gmx_fatal(FARGS,"The %d subsystems are not compatible\n",ms->nsim);
175 sfree(ibuf);
178 void gmx_log_open(const char *lognm,const t_commrec *cr,bool bMasterOnly,
179 unsigned long Flags, FILE** fplog)
181 int len,testlen,pid;
182 char buf[256],host[256];
183 time_t t;
184 FILE *fp=*fplog;
185 char *tmpnm;
187 bool bAppend = Flags & MD_APPENDFILES;
189 debug_gmx();
191 /* Communicate the filename for logfile */
192 if (cr->nnodes > 1 && !bMasterOnly) {
193 if (MASTER(cr))
194 len = strlen(lognm)+1;
195 gmx_bcast(sizeof(len),&len,cr);
196 if (!MASTER(cr))
197 snew(tmpnm,len+8);
198 else
199 tmpnm=strdup(lognm);
200 gmx_bcast(len*sizeof(*tmpnm),tmpnm,cr);
202 else
204 tmpnm=strdup(lognm);
207 debug_gmx();
209 /*for bAppend the log file is opened in checkpoint.c:read_checkpoint
210 * (for locking)
212 if (PAR(cr) && !bMasterOnly && (!bAppend || !MASTER(cr))) {
213 /* Since log always ends with '.log' let's use this info */
214 par_fn(tmpnm,efLOG,cr,cr->ms!=NULL,buf,255);
215 fp = gmx_fio_fopen(buf, bAppend ? "a+" : "w+" );
216 } else if (!bAppend) {
217 fp = gmx_fio_fopen(tmpnm, bAppend ? "a+" : "w+" );
220 sfree(tmpnm);
222 gmx_fatal_set_log_file(fp);
224 /* Get some machine parameters */
225 #ifdef HAVE_UNISTD_H
226 if( gethostname(host,255) != 0)
227 sprintf(host,"unknown");
228 #else
229 sprintf(host,"unknown");
230 #endif
232 time(&t);
234 #ifndef NO_GETPID
235 # if ((defined WIN32 || defined _WIN32 || defined WIN64 || defined _WIN64) && !defined __CYGWIN__ && !defined __CYGWIN32__)
236 pid = _getpid();
237 # else
238 pid = getpid();
239 # endif
240 #else
241 pid = 0;
242 #endif
244 if(bAppend)
246 fprintf(fp,
247 "\n\n"
248 "-----------------------------------------------------------\n"
249 "Restarting from checkpoint, appending to previous log file.\n\n"
253 fprintf(fp,
254 "Log file opened on %s"
255 "Host: %s pid: %d nodeid: %d nnodes: %d\n",
256 ctime(&t),host,pid,cr->nodeid,cr->nnodes);
258 #if (defined BUILD_MACHINE && defined BUILD_TIME && defined BUILD_USER)
259 fprintf(fp,
260 "The Gromacs distribution was built %s by\n"
261 "%s (%s)\n\n\n",BUILD_TIME,BUILD_USER,BUILD_MACHINE);
262 #endif
264 fflush(fp);
265 debug_gmx();
267 *fplog = fp;
270 void gmx_log_close(FILE *fp)
272 if (fp) {
273 gmx_fatal_set_log_file(NULL);
274 gmx_fio_fclose(fp);
278 static void comm_args(const t_commrec *cr,int *argc,char ***argv)
280 int i,len;
282 if ((cr) && PAR(cr))
283 gmx_bcast(sizeof(*argc),argc,cr);
285 if (!MASTER(cr))
286 snew(*argv,*argc+1);
287 fprintf(stderr,"NODEID=%d argc=%d\n",cr->nodeid,*argc);
288 for(i=0; (i<*argc); i++) {
289 if (MASTER(cr))
290 len = strlen((*argv)[i])+1;
291 gmx_bcast(sizeof(len),&len,cr);
292 if (!MASTER(cr))
293 snew((*argv)[i],len);
294 /*gmx_bcast(len*sizeof((*argv)[i][0]),(*argv)[i],cr);*/
295 gmx_bcast(len*sizeof(char),(*argv)[i],cr);
297 debug_gmx();
300 void init_multisystem(t_commrec *cr,int nsim, int nfile,
301 const t_filenm fnm[],bool bParFn)
303 gmx_multisim_t *ms;
304 int nnodes,nnodpersim,sim,i,ftp;
305 char buf[256];
306 #ifdef GMX_MPI
307 MPI_Group mpi_group_world;
308 #endif
309 int *rank;
311 #ifndef GMX_MPI
312 if (nsim > 1) {
313 gmx_fatal(FARGS,"This binary is compiled without MPI support, can not do multiple simulations.");
315 #endif
317 nnodes = cr->nnodes;
318 if (nnodes % nsim != 0)
319 gmx_fatal(FARGS,"The number of nodes (%d) is not a multiple of the number of simulations (%d)",nnodes,nsim);
321 nnodpersim = nnodes/nsim;
322 sim = cr->nodeid/nnodpersim;
324 if (debug)
325 fprintf(debug,"We have %d simulations, %d nodes per simulation, local simulation is %d\n",nsim,nnodpersim,sim);
327 snew(ms,1);
328 cr->ms = ms;
329 ms->nsim = nsim;
330 ms->sim = sim;
331 #ifdef GMX_MPI
332 /* Create a communicator for the master nodes */
333 snew(rank,ms->nsim);
334 for(i=0; i<ms->nsim; i++)
335 rank[i] = i*nnodpersim;
336 MPI_Comm_group(MPI_COMM_WORLD,&mpi_group_world);
337 MPI_Group_incl(mpi_group_world,nsim,rank,&ms->mpi_group_masters);
338 sfree(rank);
339 MPI_Comm_create(MPI_COMM_WORLD,ms->mpi_group_masters,
340 &ms->mpi_comm_masters);
342 #if !defined(GMX_THREADS) && !defined(MPI_IN_PLACE_EXISTS)
343 /* initialize the MPI_IN_PLACE replacement buffers */
344 snew(ms->mpb, 1);
345 ms->mpb->ibuf=NULL;
346 ms->mpb->fbuf=NULL;
347 ms->mpb->dbuf=NULL;
348 ms->mpb->ibuf_alloc=0;
349 ms->mpb->fbuf_alloc=0;
350 ms->mpb->dbuf_alloc=0;
351 #endif
353 #endif
355 /* Reduce the intra-simulation communication */
356 cr->sim_nodeid = cr->nodeid % nnodpersim;
357 cr->nnodes = nnodpersim;
358 #ifdef GMX_MPI
359 MPI_Comm_split(MPI_COMM_WORLD,sim,cr->sim_nodeid,&cr->mpi_comm_mysim);
360 cr->mpi_comm_mygroup = cr->mpi_comm_mysim;
361 cr->nodeid = cr->sim_nodeid;
362 #endif
364 if (debug) {
365 fprintf(debug,"This is simulation %d",cr->ms->sim);
366 if (PAR(cr))
367 fprintf(debug,", local number of nodes %d, local nodeid %d",
368 cr->nnodes,cr->sim_nodeid);
369 fprintf(debug,"\n\n");
372 if (bParFn) {
373 /* Patch output and tpx file names (except log which has been done already)
375 for(i=0; (i<nfile); i++) {
376 /* Because of possible multiple extensions per type we must look
377 * at the actual file name
379 if (is_output(&fnm[i]) ||
380 fnm[i].ftp == efTPX || fnm[i].ftp == efCPT ||
381 strcmp(fnm[i].opt,"-rerun") == 0) {
382 ftp = fn2ftp(fnm[i].fns[0]);
383 par_fn(fnm[i].fns[0],ftp,cr,FALSE,buf,255);
384 sfree(fnm[i].fns[0]);
385 fnm[i].fns[0] = strdup(buf);
391 t_commrec *init_par(int *argc,char ***argv_ptr)
393 t_commrec *cr;
394 char **argv;
395 int i;
396 bool pe=FALSE;
398 snew(cr,1);
400 argv = *argv_ptr;
402 #ifdef GMX_MPI
403 #ifdef GMX_LIB_MPI
404 pe = TRUE;
405 #ifdef GMX_CHECK_MPI_ENV
406 /* Do not use MPI calls when env.var. GMX_CHECK_MPI_ENV is not set */
407 if (getenv(GMX_CHECK_MPI_ENV) == NULL)
408 pe = FALSE;
409 #endif /* GMX_CHECK_MPI_ENV */
410 #endif /* GMX_LIB_MPI */
411 set_parallel_env(pe);
412 if (pe) {
413 cr->sim_nodeid = gmx_setup(argc,argv,&cr->nnodes);
414 } else {
415 cr->nnodes = 1;
416 cr->sim_nodeid = 0;
418 #else /* GMX_MPI */
419 pe=FALSE;
420 set_parallel_env(pe);
421 cr->sim_nodeid = 0;
422 cr->nnodes = 1;
423 #endif /* GMX_MPI */
425 if (!PAR(cr) && (cr->sim_nodeid != 0))
426 gmx_comm("(!PAR(cr) && (cr->sim_nodeid != 0))");
428 if (PAR(cr))
430 #ifdef GMX_MPI
431 cr->mpi_comm_mysim = MPI_COMM_WORLD;
432 cr->mpi_comm_mygroup = cr->mpi_comm_mysim;
433 #endif /* GMX_MPI */
435 cr->nodeid = cr->sim_nodeid;
437 cr->duty = (DUTY_PP | DUTY_PME);
439 /* Communicate arguments if parallel */
440 #ifndef GMX_THREADS
441 if (PAR(cr))
442 comm_args(cr,argc,argv_ptr);
443 #endif /* GMX_THREADS */
445 #ifdef GMX_MPI
446 #if !defined(GMX_THREADS) && !defined(MPI_IN_PLACE_EXISTS)
447 /* initialize the MPI_IN_PLACE replacement buffers */
448 snew(cr->mpb, 1);
449 cr->mpb->ibuf=NULL;
450 cr->mpb->fbuf=NULL;
451 cr->mpb->dbuf=NULL;
452 cr->mpb->ibuf_alloc=0;
453 cr->mpb->fbuf_alloc=0;
454 cr->mpb->dbuf_alloc=0;
455 #endif
456 #endif
458 return cr;
461 t_commrec *init_par_threads(t_commrec *cro)
463 #ifdef GMX_THREADS
464 int initialized;
465 t_commrec *cr;
467 snew(cr,1);
468 /* now copy the whole thing, so settings like the number of PME nodes
469 get propagated. */
470 *cr=*cro;
471 /* and we start setting our own thread-specific values for things */
472 MPI_Initialized(&initialized);
473 if (!initialized)
474 gmx_comm("Initializing threads without comm");
475 set_parallel_env(TRUE);
476 /* once threads will be used together with MPI, we'll
477 fill the cr structure with distinct data here. This might even work: */
478 cr->sim_nodeid = gmx_setup(0,NULL, &cr->nnodes);
479 /* note that we're explicitly using tMPI */
480 tMPI_Comm_size(TMPI_COMM_WORLD, &cr->nthreads);
481 tMPI_Comm_rank(TMPI_COMM_WORLD, &cr->threadid);
483 cr->mpi_comm_mysim = MPI_COMM_WORLD;
484 cr->mpi_comm_mygroup = cr->mpi_comm_mysim;
485 cr->nodeid = cr->sim_nodeid;
486 cr->duty = (DUTY_PP | DUTY_PME);
488 return cr;
489 #else
490 return NULL;
491 #endif
494 void cancel_par_threads(t_commrec *cr)
496 #ifdef GMX_THREADS
497 if (MASTERTHREAD(cr))
499 fprintf(stderr, "\nCANCELING %d THREADS\n", cr->nthreads-1);
500 cr->nnodes = 1;
501 cr->nthreads = 1;
502 cr->sim_nodeid = 0;
503 cr->nodeid = 0;
504 cr->threadid = 0;
505 cr->duty = (DUTY_PP | DUTY_PME);
507 else
509 tMPI_Finalize();
511 #else
512 gmx_incons("cancel_par_threads called without thread support");
513 #endif
516 t_commrec *init_cr_nopar(void)
518 t_commrec *cr;
520 snew(cr,1);
522 cr->nnodes = 1;
523 cr->nthreads = 1;
524 cr->sim_nodeid = 0;
525 cr->nodeid = 0;
526 cr->threadid = 0;
527 cr->duty = (DUTY_PP | DUTY_PME);
529 return cr;