Compiles in MSVC with cmake
[alexxy/gromacs.git] / src / gmxlib / main.c
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.
14
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.
19  * 
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.
26  * 
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.
29  * 
30  * For more info, check our website at http://www.gromacs.org
31  * 
32  * And Hey:
33  * GROningen Mixture of Alchemy and Childrens' Stories
34  */
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
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"
53
54 #ifdef GMX_THREADS
55 #include "thread_mpi.h"
56 #endif
57
58 /* The source code in this file should be thread-safe. 
59          Please keep it that way. */
60
61
62 #ifdef HAVE_UNISTD_H
63 #include <unistd.h>
64 #endif
65
66 #if ((defined WIN32 || defined _WIN32 || defined WIN64 || defined _WIN64) && !defined __CYGWIN__ && !defined __CYGWIN32__)
67 #include <process.h>
68 #endif
69
70
71 #define BUFSIZE 1024
72
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
80
81 bool gmx_parallel_env(void)
82 {
83     bool ret;
84 #ifdef GMX_THREADS
85     tMPI_Thread_mutex_lock(&parallel_env_mutex);
86 #endif
87     ret=parallel_env_val;
88 #ifdef GMX_THREADS
89     tMPI_Thread_mutex_unlock(&parallel_env_mutex);
90 #endif
91     return ret;
92 }
93
94 static void set_parallel_env(bool val)
95 {
96 #ifdef GMX_THREADS
97     tMPI_Thread_mutex_lock(&parallel_env_mutex);
98 #endif
99     if (!parallel_env_val)
100     {
101         /* we only allow it to be set, not unset */
102         parallel_env_val=val;
103     }
104 #ifdef GMX_THREADS
105     tMPI_Thread_mutex_unlock(&parallel_env_mutex);
106 #endif
107 }
108
109
110
111 static void par_fn(char *base,int ftp,const t_commrec *cr,
112                    bool bUnderScore,
113                    char buf[],int bufsize)
114 {
115   int n;
116   
117   if((size_t)bufsize<(strlen(base)+4))
118      gmx_mem("Character buffer too small!");
119
120   /* Copy to buf, and strip extension */
121   strcpy(buf,base);
122   buf[strlen(base) - strlen(ftp2ext(fn2ftp(base))) - 1] = '\0';
123
124   /* Add node info */
125   if (bUnderScore)
126     strcat(buf,"_");
127   if (MULTISIM(cr) && !bUnderScore) {
128     sprintf(buf+strlen(buf),"%d",cr->ms->sim);
129   } else if (PAR(cr)) {
130     sprintf(buf+strlen(buf),"%d",cr->nodeid);
131   }
132   strcat(buf,".");
133   
134   /* Add extension again */
135   strcat(buf,(ftp == efTPX) ? "tpr" : (ftp == efEDR) ? "edr" : ftp2ext(ftp));
136 }
137
138 void check_multi_int(FILE *log,const gmx_multisim_t *ms,int val,
139                      const char *name)
140 {
141   int  *ibuf,p;
142   bool bCompatible;
143
144   fprintf(log,"Multi-checking %s ... ",name);
145   
146   if (ms == NULL)
147     gmx_fatal(FARGS,
148               "check_multi_int called with a NULL communication pointer");
149
150   snew(ibuf,ms->nsim);
151   ibuf[ms->sim] = val;
152   gmx_sumi_sim(ms->nsim,ibuf,ms);
153   
154   bCompatible = TRUE;
155   for(p=1; p<ms->nsim; p++)
156     bCompatible = bCompatible && (ibuf[p-1] == ibuf[p]);
157   
158   if (bCompatible) 
159     fprintf(log,"OK\n");
160   else {
161     fprintf(log,"\n%s is not equal for all subsystems\n",name);
162     for(p=0; p<ms->nsim; p++)
163       fprintf(log,"  subsystem %d: %d\n",p,ibuf[p]);
164     gmx_fatal(FARGS,"The %d subsystems are not compatible\n",ms->nsim);
165   }
166   
167   sfree(ibuf);
168 }
169
170 FILE *gmx_log_open(const char *lognm,const t_commrec *cr,bool bMasterOnly, 
171                    unsigned long Flags)
172 {
173   int  len,testlen,pid;
174   char buf[256],host[256];
175   time_t t;
176   FILE *fp;
177   char *tmpnm;
178
179   bool bAppend = Flags & MD_APPENDFILES;        
180   
181   debug_gmx();
182   
183   /* Communicate the filename for logfile */
184   if (cr->nnodes > 1 && !bMasterOnly) {
185       if (MASTER(cr))
186           len = strlen(lognm)+1;
187       gmx_bcast(sizeof(len),&len,cr);
188       if (!MASTER(cr))
189           snew(tmpnm,len+8);
190       else
191           tmpnm=strdup(lognm);
192       gmx_bcast(len*sizeof(*tmpnm),tmpnm,cr);
193   }
194   else
195   {
196       tmpnm=strdup(lognm);
197   }
198   
199   debug_gmx();
200
201   if (PAR(cr) && !bMasterOnly) {
202     /* Since log always ends with '.log' let's use this info */
203     par_fn(tmpnm,efLOG,cr,cr->ms!=NULL,buf,255);
204           fp = gmx_fio_fopen(buf, bAppend ? "a" : "w" );
205   } else {
206           fp = gmx_fio_fopen(tmpnm, bAppend ? "a" : "w" );
207   }
208
209   sfree(tmpnm);
210
211   gmx_fatal_set_log_file(fp);
212   
213   /* Get some machine parameters */
214 #ifdef HAVE_UNISTD_H
215   if( gethostname(host,255) != 0)
216     sprintf(host,"unknown");
217 #else
218   sprintf(host,"unknown");
219 #endif  
220
221   time(&t);
222
223 #ifndef NO_GETPID
224 #   if ((defined WIN32 || defined _WIN32 || defined WIN64 || defined _WIN64) && !defined __CYGWIN__ && !defined __CYGWIN32__)
225           pid = _getpid();
226 #   else
227           pid = getpid();
228 #   endif
229 #else
230         pid = 0;
231 #endif
232
233   if(bAppend)
234   {
235           fprintf(fp,
236                           "\n\n"
237                           "-----------------------------------------------------------\n"
238                           "Restarting from checkpoint, appending to previous log file.\n\n"
239                           );
240   }
241         
242   fprintf(fp,
243           "Log file opened on %s"
244           "Host: %s  pid: %d  nodeid: %d  nnodes:  %d\n",
245           ctime(&t),host,pid,cr->nodeid,cr->nnodes);
246
247 #if (defined BUILD_MACHINE && defined BUILD_TIME && defined BUILD_USER) 
248   fprintf(fp,
249           "The Gromacs distribution was built %s by\n"
250           "%s (%s)\n\n\n",BUILD_TIME,BUILD_USER,BUILD_MACHINE);
251 #endif
252
253   fflush(fp);
254   debug_gmx();
255
256   return fp;
257 }
258
259 void gmx_log_close(FILE *fp)
260 {
261   if (fp) {
262     gmx_fatal_set_log_file(NULL);
263     gmx_fio_fclose(fp);
264   }
265 }
266
267 static void comm_args(const t_commrec *cr,int *argc,char ***argv)
268 {
269   int i,len;
270   
271   if ((cr) && PAR(cr))
272     gmx_bcast(sizeof(*argc),argc,cr);
273   
274   if (!MASTER(cr))
275     snew(*argv,*argc+1);
276   fprintf(stderr,"NODEID=%d argc=%d\n",cr->nodeid,*argc);
277   for(i=0; (i<*argc); i++) {
278     if (MASTER(cr))
279       len = strlen((*argv)[i])+1;
280     gmx_bcast(sizeof(len),&len,cr);
281     if (!MASTER(cr))
282       snew((*argv)[i],len);
283     /*gmx_bcast(len*sizeof((*argv)[i][0]),(*argv)[i],cr);*/
284     gmx_bcast(len*sizeof(char),(*argv)[i],cr);
285   }
286   debug_gmx();
287 }
288
289 void init_multisystem(t_commrec *cr,int nsim, int nfile,
290                       const t_filenm fnm[],bool bParFn)
291 {
292   gmx_multisim_t *ms;
293   int  nnodes,nnodpersim,sim,i,ftp;
294   char buf[256];
295 #ifdef GMX_MPI
296   MPI_Group mpi_group_world;
297 #endif  
298   int *rank;
299
300   nnodes  = cr->nnodes;
301   if (nnodes % nsim != 0)
302     gmx_fatal(FARGS,"The number of nodes (%d) is not a multiple of the number of simulations (%d)",nnodes,nsim);
303
304   nnodpersim = nnodes/nsim;
305   sim = cr->nodeid/nnodpersim;
306
307   if (debug)
308     fprintf(debug,"We have %d simulations, %d nodes per simulation, local simulation is %d\n",nsim,nnodpersim,sim);
309
310   snew(ms,1);
311   cr->ms = ms;
312   ms->nsim = nsim;
313   ms->sim  = sim;
314 #ifdef GMX_MPI
315   /* Create a communicator for the master nodes */
316   snew(rank,ms->nsim);
317   for(i=0; i<ms->nsim; i++)
318     rank[i] = i*nnodpersim;
319   MPI_Comm_group(MPI_COMM_WORLD,&mpi_group_world);
320   MPI_Group_incl(mpi_group_world,nsim,rank,&ms->mpi_group_masters);
321   sfree(rank);
322   MPI_Comm_create(MPI_COMM_WORLD,ms->mpi_group_masters,
323                   &ms->mpi_comm_masters);
324 #endif
325
326   /* Reduce the intra-simulation communication */
327   cr->sim_nodeid = cr->nodeid % nnodpersim;
328   cr->nnodes = nnodpersim;
329 #ifdef GMX_MPI
330   MPI_Comm_split(MPI_COMM_WORLD,sim,cr->sim_nodeid,&cr->mpi_comm_mysim);
331   cr->mpi_comm_mygroup = cr->mpi_comm_mysim;
332   cr->nodeid = cr->sim_nodeid;
333 #endif
334
335   if (debug) {
336     fprintf(debug,"This is simulation %d",cr->ms->sim);
337     if (PAR(cr))
338       fprintf(debug,", local number of nodes %d, local nodeid %d",
339               cr->nnodes,cr->sim_nodeid);
340     fprintf(debug,"\n\n");
341   }
342
343   if (bParFn) {
344     /* Patch output and tpx file names (except log which has been done already)
345      */
346     for(i=0; (i<nfile); i++) {
347       /* Because of possible multiple extensions per type we must look 
348        * at the actual file name 
349        */
350       if (is_output(&fnm[i]) ||
351           fnm[i].ftp == efTPX || fnm[i].ftp == efCPT ||
352           strcmp(fnm[i].opt,"-rerun") == 0) {
353         ftp = fn2ftp(fnm[i].fns[0]);
354         par_fn(fnm[i].fns[0],ftp,cr,FALSE,buf,255);
355         sfree(fnm[i].fns[0]);
356         fnm[i].fns[0] = strdup(buf);
357       }
358     }
359   }
360 }
361
362 t_commrec *init_par(int *argc,char ***argv_ptr)
363 {
364     t_commrec *cr;
365     char      **argv;
366     int       i;
367     bool      pe=FALSE;
368
369     snew(cr,1);
370
371     argv = *argv_ptr;
372
373 #ifdef GMX_MPI
374 #if 0
375 #ifdef GMX_THREADS
376     if (tMPI_Get_N(argc, argv_ptr)>1)
377         pe=TRUE;
378     else
379         pe=FALSE;
380 #endif /* GMX_THREADS */
381 #endif
382 #ifdef GMX_LIB_MPI
383     pe = TRUE;
384 #ifdef GMX_CHECK_MPI_ENV
385     /* Do not use MPI calls when env.var. GMX_CHECK_MPI_ENV is not set */
386     if (getenv(GMX_CHECK_MPI_ENV) == NULL)
387         pe = FALSE;
388 #endif /* GMX_CHECK_MPI_ENV */
389 #endif /* GMX_LIB_MPI  */
390     set_parallel_env(pe);
391     if (pe) {
392         cr->sim_nodeid = gmx_setup(argc,argv,&cr->nnodes);
393     } else {
394         cr->nnodes     = 1;
395         cr->sim_nodeid = 0;
396     }
397 #else /* GMX_MPI */
398     pe=FALSE;
399     set_parallel_env(pe);
400     cr->sim_nodeid   = 0;
401     cr->nnodes       = 1;
402 #endif /* GMX_MPI */
403
404     if (!PAR(cr) && (cr->sim_nodeid != 0))
405         gmx_comm("(!PAR(cr) && (cr->sim_nodeid != 0))");
406
407     if (PAR(cr)) {
408 #ifdef GMX_MPI
409         cr->mpi_comm_mysim = MPI_COMM_WORLD;
410         cr->mpi_comm_mygroup = cr->mpi_comm_mysim;
411 #endif /* GMX_MPI */
412     }
413     cr->nodeid = cr->sim_nodeid;
414
415     cr->duty = (DUTY_PP | DUTY_PME);
416
417     /* Communicate arguments if parallel */
418 #ifndef GMX_THREADS
419     if (PAR(cr))
420         comm_args(cr,argc,argv_ptr);
421 #endif /* GMX_THREADS */
422
423     return cr;
424 }
425
426 t_commrec *init_par_threads(t_commrec *cro)
427 {
428 #ifdef GMX_THREADS
429     int initialized;
430     t_commrec *cr;
431
432     snew(cr,1);
433     MPI_Initialized(&initialized);
434     if (!initialized)
435         gmx_comm("Initializing threads without comm");
436     set_parallel_env(TRUE);
437     /* once threads will be getting along with MPI, we'll
438        fill the cr structure with more sensible data here */
439     cr->sim_nodeid = gmx_setup(0,NULL, &cr->nnodes);
440
441     cr->mpi_comm_mysim = MPI_COMM_WORLD;
442     cr->mpi_comm_mygroup = cr->mpi_comm_mysim;
443     cr->nodeid = cr->sim_nodeid;
444     cr->duty = (DUTY_PP | DUTY_PME);
445
446     return cr;
447 #else
448     return NULL;
449 #endif
450 }
451
452 t_commrec *init_cr_nopar(void)
453 {
454   t_commrec *cr;
455
456   snew(cr,1);
457
458   cr->nnodes     = 1; 
459   cr->sim_nodeid = 0;
460   cr->nodeid     = 0;
461   cr->nthreads   = 1;
462   cr->threadid   = 0;
463   cr->duty       = (DUTY_PP | DUTY_PME);
464
465   return cr;
466 }