Unused TestReferenceData no longer does anything.
[alexxy/gromacs.git] / src / gromacs / gmxlib / network.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 <string.h>
40 #include "gmx_fatal.h"
41 #include "main.h"
42 #include "smalloc.h"
43 #include "network.h"
44 #include "copyrite.h"
45 #include "statutil.h"
46 #include "ctype.h"
47 #include "macros.h"
48
49 #ifdef GMX_LIB_MPI
50 #include <mpi.h>
51 #endif
52
53 #ifdef GMX_THREADS
54 #include "tmpi.h"
55 #endif
56
57
58 /* The source code in this file should be thread-safe. 
59       Please keep it that way. */
60
61 gmx_bool gmx_mpi_initialized(void)
62 {
63   int n;
64 #ifndef GMX_MPI
65   return 0;
66 #else
67   MPI_Initialized(&n);
68   
69   return n;
70 #endif
71 }
72
73 int gmx_setup(int *argc,char **argv,int *nnodes)
74 {
75 #ifndef GMX_MPI
76   gmx_call("gmx_setup");
77   return 0;
78 #else
79   char   buf[256];
80   int    resultlen;               /* actual length of node name      */
81   int    i,flag;
82   int  mpi_num_nodes;
83   int  mpi_my_rank;
84   char mpi_hostname[MPI_MAX_PROCESSOR_NAME];
85
86   /* Call the MPI routines */
87 #ifdef GMX_LIB_MPI
88 #ifdef GMX_FAHCORE
89   (void) fah_MPI_Init(argc,&argv);
90 #else
91   (void) MPI_Init(argc,&argv);
92 #endif
93 #endif
94   (void) MPI_Comm_size( MPI_COMM_WORLD, &mpi_num_nodes );
95   (void) MPI_Comm_rank( MPI_COMM_WORLD, &mpi_my_rank );
96   (void) MPI_Get_processor_name( mpi_hostname, &resultlen );
97  
98 #ifdef GMX_LIB_MPI 
99   fprintf(stderr,"NNODES=%d, MYRANK=%d, HOSTNAME=%s\n",
100           mpi_num_nodes,mpi_my_rank,mpi_hostname);
101 #endif
102   
103   *nnodes=mpi_num_nodes;
104   
105   return mpi_my_rank;
106 #endif
107 }
108
109 int  gmx_node_num(void)
110 {
111 #ifndef GMX_MPI
112   return 1;
113 #else
114   int i;
115   (void) MPI_Comm_size(MPI_COMM_WORLD, &i);
116   return i;
117 #endif
118 }
119
120 int gmx_node_rank(void)
121 {
122 #ifndef GMX_MPI
123   return 0;
124 #else
125   int i;
126   (void) MPI_Comm_rank(MPI_COMM_WORLD, &i);
127   return i;
128 #endif
129 }
130
131 void gmx_setup_nodecomm(FILE *fplog,t_commrec *cr)
132 {
133   gmx_nodecomm_t *nc;
134   int  n,rank,resultlen,hostnum,i,j,ng,ni;
135 #ifdef GMX_MPI
136   char mpi_hostname[MPI_MAX_PROCESSOR_NAME],num[MPI_MAX_PROCESSOR_NAME];
137 #endif
138
139   /* Many MPI implementations do not optimize MPI_Allreduce
140    * (and probably also other global communication calls)
141    * for multi-core nodes connected by a network.
142    * We can optimize such communication by using one MPI call
143    * within each node and one between the nodes.
144    * For MVAPICH2 and Intel MPI this reduces the time for
145    * the global_stat communication by 25%
146    * for 2x2-core 3 GHz Woodcrest connected by mixed DDR/SDR Infiniband.
147    * B. Hess, November 2007
148    */
149
150   nc = &cr->nc;
151
152   nc->bUse = FALSE;
153 #ifndef GMX_THREADS
154   if (getenv("GMX_NO_NODECOMM") == NULL) {
155 #ifdef GMX_MPI
156     MPI_Comm_size(cr->mpi_comm_mygroup,&n);
157     MPI_Comm_rank(cr->mpi_comm_mygroup,&rank);
158     MPI_Get_processor_name(mpi_hostname,&resultlen);
159     /* This procedure can only differentiate nodes with host names
160      * that end on unique numbers.
161      */
162     i = 0;
163     j = 0;
164     /* Only parse the host name up to the first dot */
165     while(i < resultlen && mpi_hostname[i] != '.') {
166       if (isdigit(mpi_hostname[i])) {
167         num[j++] = mpi_hostname[i];
168       }
169       i++;
170     }
171     num[j] = '\0';
172     if (j == 0) {
173       hostnum = 0;
174     } else {
175       /* Use only the last 9 decimals, so we don't overflow an int */
176       hostnum = strtol(num + max(0,j-9), NULL, 10); 
177     }
178
179     if (debug) {
180       fprintf(debug,
181               "In gmx_setup_nodecomm: splitting communicator of size %d\n",
182               n);
183       fprintf(debug,"In gmx_setup_nodecomm: hostname '%s', hostnum %d\n",
184               mpi_hostname,hostnum);
185     }
186
187     /* The intra-node communicator, split on node number */
188     MPI_Comm_split(cr->mpi_comm_mygroup,hostnum,rank,&nc->comm_intra);
189     MPI_Comm_rank(nc->comm_intra,&nc->rank_intra);
190     if (debug) {
191       fprintf(debug,"In gmx_setup_nodecomm: node rank %d rank_intra %d\n",
192               rank,nc->rank_intra);
193     }
194     /* The inter-node communicator, split on rank_intra.
195      * We actually only need the one for rank=0,
196      * but it is easier to create them all.
197      */
198     MPI_Comm_split(cr->mpi_comm_mygroup,nc->rank_intra,rank,&nc->comm_inter);
199     /* Check if this really created two step communication */
200     MPI_Comm_size(nc->comm_inter,&ng);
201     MPI_Comm_size(nc->comm_intra,&ni);
202     if (debug) {
203       fprintf(debug,"In gmx_setup_nodecomm: groups %d, my group size %d\n",
204               ng,ni);
205     }
206     if ((ng > 1 && ng < n) || (ni > 1 && ni < n)) {
207       nc->bUse = TRUE;
208       if (fplog)
209         fprintf(fplog,"Using two step summing over %d groups of on average %.1f processes\n\n",ng,(real)n/(real)ng);
210       if (nc->rank_intra > 0)
211         MPI_Comm_free(&nc->comm_inter);
212     } else {
213       /* One group or all processes in a separate group, use normal summing */
214       MPI_Comm_free(&nc->comm_inter);
215       MPI_Comm_free(&nc->comm_intra);
216     }
217 #endif
218   }
219 #endif
220 }
221
222 void gmx_barrier(const t_commrec *cr)
223 {
224 #ifndef GMX_MPI
225   gmx_call("gmx_barrier");
226 #else
227   MPI_Barrier(cr->mpi_comm_mygroup);
228 #endif
229 }
230
231 void gmx_abort(int noderank,int nnodes,int errorno)
232 {
233 #ifndef GMX_MPI
234   gmx_call("gmx_abort");
235 #else
236 #ifdef GMX_THREADS
237   fprintf(stderr,"Halting program %s\n",ShortProgram());
238   thanx(stderr);
239   exit(1);
240 #else
241   if (nnodes > 1)
242   {
243       fprintf(stderr,"Halting parallel program %s on CPU %d out of %d\n",
244               ShortProgram(),noderank,nnodes);
245   }
246   else
247   {
248       fprintf(stderr,"Halting program %s\n",ShortProgram());
249   }
250
251   thanx(stderr);
252   MPI_Abort(MPI_COMM_WORLD,errorno);
253   exit(1);
254 #endif
255 #endif
256 }
257
258 void gmx_bcast(int nbytes,void *b,const t_commrec *cr)
259 {
260 #ifndef GMX_MPI
261   gmx_call("gmx_bast");
262 #else
263   MPI_Bcast(b,nbytes,MPI_BYTE,MASTERRANK(cr),cr->mpi_comm_mygroup);
264 #endif
265 }
266
267 void gmx_bcast_sim(int nbytes,void *b,const t_commrec *cr)
268 {
269 #ifndef GMX_MPI
270   gmx_call("gmx_bast");
271 #else
272   MPI_Bcast(b,nbytes,MPI_BYTE,MASTERRANK(cr),cr->mpi_comm_mysim);
273 #endif
274 }
275
276 void gmx_sumd(int nr,double r[],const t_commrec *cr)
277 {
278 #ifndef GMX_MPI
279     gmx_call("gmx_sumd");
280 #else
281 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREADS)
282     if (cr->nc.bUse) {
283         if (cr->nc.rank_intra == 0)
284         {
285             /* Use two step summing. */
286             MPI_Reduce(MPI_IN_PLACE,r,nr,MPI_DOUBLE,MPI_SUM,0,
287                        cr->nc.comm_intra);
288             /* Sum the roots of the internal (intra) buffers. */
289             MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_DOUBLE,MPI_SUM,
290                           cr->nc.comm_inter);
291         }
292         else
293         {
294             /* This is here because of the silly MPI specification
295                 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
296             MPI_Reduce(r,NULL,nr,MPI_DOUBLE,MPI_SUM,0,cr->nc.comm_intra);
297         }
298         MPI_Bcast(r,nr,MPI_DOUBLE,0,cr->nc.comm_intra);
299     } 
300     else 
301     {
302         MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_DOUBLE,MPI_SUM, 
303                       cr->mpi_comm_mygroup);
304     }
305 #else
306     int i;
307
308     if (nr > cr->mpb->dbuf_alloc) {
309         cr->mpb->dbuf_alloc = nr;
310         srenew(cr->mpb->dbuf,cr->mpb->dbuf_alloc);
311     }
312     if (cr->nc.bUse) {
313         /* Use two step summing */
314         MPI_Allreduce(r,cr->mpb->dbuf,nr,MPI_DOUBLE,MPI_SUM,cr->nc.comm_intra);
315         if (cr->nc.rank_intra == 0) {
316             /* Sum with the buffers reversed */
317             MPI_Allreduce(cr->mpb->dbuf,r,nr,MPI_DOUBLE,MPI_SUM, 
318                           cr->nc.comm_inter);
319         }
320         MPI_Bcast(r,nr,MPI_DOUBLE,0,cr->nc.comm_intra);
321     } else {
322         MPI_Allreduce(r,cr->mpb->dbuf,nr,MPI_DOUBLE,MPI_SUM,
323                       cr->mpi_comm_mygroup);
324         for(i=0; i<nr; i++)
325             r[i] = cr->mpb->dbuf[i];
326     }
327 #endif
328 #endif
329 }
330
331 void gmx_sumf(int nr,float r[],const t_commrec *cr)
332 {
333 #ifndef GMX_MPI
334     gmx_call("gmx_sumf");
335 #else
336 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREADS)
337     if (cr->nc.bUse) {
338         /* Use two step summing.  */
339         if (cr->nc.rank_intra == 0)
340         {
341             MPI_Reduce(MPI_IN_PLACE,r,nr,MPI_FLOAT,MPI_SUM,0,
342                        cr->nc.comm_intra);
343             /* Sum the roots of the internal (intra) buffers */
344             MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_FLOAT,MPI_SUM,
345                           cr->nc.comm_inter);
346         }
347         else
348         {
349             /* This is here because of the silly MPI specification
350                 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
351             MPI_Reduce(r,NULL,nr,MPI_FLOAT,MPI_SUM,0,cr->nc.comm_intra);
352         }
353         MPI_Bcast(r,nr,MPI_FLOAT,0,cr->nc.comm_intra);
354     } 
355     else 
356     {
357         MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_FLOAT,MPI_SUM,cr->mpi_comm_mygroup);
358     }
359 #else
360     int i;
361
362     if (nr > cr->mpb->fbuf_alloc) {
363         cr->mpb->fbuf_alloc = nr;
364         srenew(cr->mpb->fbuf,cr->mpb->fbuf_alloc);
365     }
366     if (cr->nc.bUse) {
367         /* Use two step summing */
368         MPI_Allreduce(r,cr->mpb->fbuf,nr,MPI_FLOAT,MPI_SUM,cr->nc.comm_intra);
369         if (cr->nc.rank_intra == 0) {
370             /* Sum with the buffers reversed */
371             MPI_Allreduce(cr->mpb->fbuf,r,nr,MPI_FLOAT,MPI_SUM, 
372                           cr->nc.comm_inter);
373         }
374         MPI_Bcast(r,nr,MPI_FLOAT,0,cr->nc.comm_intra);
375     } else {
376         MPI_Allreduce(r,cr->mpb->fbuf,nr,MPI_FLOAT,MPI_SUM,
377                       cr->mpi_comm_mygroup);
378         for(i=0; i<nr; i++)
379             r[i] = cr->mpb->fbuf[i];
380     }
381 #endif
382 #endif
383 }
384
385 void gmx_sumi(int nr,int r[],const t_commrec *cr)
386 {
387 #ifndef GMX_MPI
388     gmx_call("gmx_sumi");
389 #else
390 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREADS)
391     if (cr->nc.bUse) {
392         /* Use two step summing */
393         if (cr->nc.rank_intra == 0) 
394         {
395             MPI_Reduce(MPI_IN_PLACE,r,nr,MPI_INT,MPI_SUM,0,cr->nc.comm_intra);
396             /* Sum with the buffers reversed */
397             MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_INT,MPI_SUM,cr->nc.comm_inter);
398         }
399         else
400         {
401             /* This is here because of the silly MPI specification
402                 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
403             MPI_Reduce(r,NULL,nr,MPI_INT,MPI_SUM,0,cr->nc.comm_intra);
404         }
405         MPI_Bcast(r,nr,MPI_INT,0,cr->nc.comm_intra);
406     } 
407     else 
408     {
409         MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_INT,MPI_SUM,cr->mpi_comm_mygroup);
410     }
411 #else
412     int i;
413
414     if (nr > cr->mpb->ibuf_alloc) {
415         cr->mpb->ibuf_alloc = nr;
416         srenew(cr->mpb->ibuf,cr->mpb->ibuf_alloc);
417     }
418     if (cr->nc.bUse) {
419         /* Use two step summing */
420         MPI_Allreduce(r,cr->mpb->ibuf,nr,MPI_INT,MPI_SUM,cr->nc.comm_intra);
421         if (cr->nc.rank_intra == 0) {
422             /* Sum with the buffers reversed */
423             MPI_Allreduce(cr->mpb->ibuf,r,nr,MPI_INT,MPI_SUM,cr->nc.comm_inter);
424         }
425         MPI_Bcast(r,nr,MPI_INT,0,cr->nc.comm_intra);
426     } else {
427         MPI_Allreduce(r,cr->mpb->ibuf,nr,MPI_INT,MPI_SUM,cr->mpi_comm_mygroup);
428         for(i=0; i<nr; i++)
429             r[i] = cr->mpb->ibuf[i];
430     }
431 #endif
432 #endif
433 }
434
435 void gmx_sumli(int nr,gmx_large_int_t r[],const t_commrec *cr)
436 {
437 #ifndef GMX_MPI
438     gmx_call("gmx_sumli");
439 #else
440 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREADS)
441     if (cr->nc.bUse) {
442         /* Use two step summing */
443         if (cr->nc.rank_intra == 0) 
444         {
445             MPI_Reduce(MPI_IN_PLACE,r,nr,GMX_MPI_LARGE_INT,MPI_SUM,0,
446                        cr->nc.comm_intra);
447             /* Sum with the buffers reversed */
448             MPI_Allreduce(MPI_IN_PLACE,r,nr,GMX_MPI_LARGE_INT,MPI_SUM,
449                           cr->nc.comm_inter);
450         }
451         else
452         {
453             /* This is here because of the silly MPI specification
454                 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
455             MPI_Reduce(r,NULL,nr,GMX_MPI_LARGE_INT,MPI_SUM,0,cr->nc.comm_intra);
456         }
457         MPI_Bcast(r,nr,GMX_MPI_LARGE_INT,0,cr->nc.comm_intra);
458     } 
459     else 
460     {
461         MPI_Allreduce(MPI_IN_PLACE,r,nr,GMX_MPI_LARGE_INT,MPI_SUM,cr->mpi_comm_mygroup);
462     }
463 #else
464     int i;
465
466     if (nr > cr->mpb->ibuf_alloc) {
467         cr->mpb->ibuf_alloc = nr;
468         srenew(cr->mpb->ibuf,cr->mpb->ibuf_alloc);
469     }
470     if (cr->nc.bUse) {
471         /* Use two step summing */
472         MPI_Allreduce(r,cr->mpb->ibuf,nr,GMX_MPI_LARGE_INT,MPI_SUM,
473                       cr->nc.comm_intra);
474         if (cr->nc.rank_intra == 0) {
475             /* Sum with the buffers reversed */
476             MPI_Allreduce(cr->mpb->ibuf,r,nr,GMX_MPI_LARGE_INT,MPI_SUM,
477                           cr->nc.comm_inter);
478         }
479         MPI_Bcast(r,nr,GMX_MPI_LARGE_INT,0,cr->nc.comm_intra);
480     } else {
481         MPI_Allreduce(r,cr->mpb->ibuf,nr,GMX_MPI_LARGE_INT,MPI_SUM,
482                       cr->mpi_comm_mygroup);
483         for(i=0; i<nr; i++)
484             r[i] = cr->mpb->ibuf[i];
485     }
486 #endif
487 #endif
488 }
489
490
491
492 #ifdef GMX_MPI
493 void gmx_sumd_comm(int nr,double r[],MPI_Comm mpi_comm)
494 {
495 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREADS)
496     MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_DOUBLE,MPI_SUM,mpi_comm);
497 #else
498     /* this function is only used in code that is not performance critical,
499        (during setup, when comm_rec is not the appropriate communication  
500        structure), so this isn't as bad as it looks. */
501     double *buf;
502     int i;
503
504     snew(buf, nr);
505     MPI_Allreduce(r,buf,nr,MPI_DOUBLE,MPI_SUM,mpi_comm);
506     for(i=0; i<nr; i++)
507         r[i] = buf[i];
508     sfree(buf);
509 #endif
510 }
511 #endif
512
513 #ifdef GMX_MPI
514 void gmx_sumf_comm(int nr,float r[],MPI_Comm mpi_comm)
515 {
516 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREADS)
517     MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_FLOAT,MPI_SUM,mpi_comm);
518 #else
519     /* this function is only used in code that is not performance critical,
520        (during setup, when comm_rec is not the appropriate communication  
521        structure), so this isn't as bad as it looks. */
522     float *buf;
523     int i;
524
525     snew(buf, nr);
526     MPI_Allreduce(r,buf,nr,MPI_FLOAT,MPI_SUM,mpi_comm);
527     for(i=0; i<nr; i++)
528         r[i] = buf[i];
529     sfree(buf);
530 #endif
531 }
532 #endif
533
534 void gmx_sumd_sim(int nr,double r[],const gmx_multisim_t *ms)
535 {
536 #ifndef GMX_MPI
537   gmx_call("gmx_sumd_sim");
538 #else
539   gmx_sumd_comm(nr,r,ms->mpi_comm_masters);
540 #endif
541 }
542
543 void gmx_sumf_sim(int nr,float r[],const gmx_multisim_t *ms)
544 {
545 #ifndef GMX_MPI
546   gmx_call("gmx_sumf_sim");
547 #else
548   gmx_sumf_comm(nr,r,ms->mpi_comm_masters);
549 #endif
550 }
551
552 void gmx_sumi_sim(int nr,int r[], const gmx_multisim_t *ms)
553 {
554 #ifndef GMX_MPI
555     gmx_call("gmx_sumi_sim");
556 #else
557 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREADS)
558     MPI_Allreduce(MPI_IN_PLACE,r,nr,MPI_INT,MPI_SUM,ms->mpi_comm_masters);
559 #else
560     /* this is thread-unsafe, but it will do for now: */
561     int i;
562
563     if (nr > ms->mpb->ibuf_alloc) {
564         ms->mpb->ibuf_alloc = nr;
565         srenew(ms->mpb->ibuf,ms->mpb->ibuf_alloc);
566     }
567     MPI_Allreduce(r,ms->mpb->ibuf,nr,MPI_INT,MPI_SUM,ms->mpi_comm_masters);
568     for(i=0; i<nr; i++)
569         r[i] = ms->mpb->ibuf[i];
570 #endif
571 #endif
572 }
573
574 void gmx_sumli_sim(int nr,gmx_large_int_t r[], const gmx_multisim_t *ms)
575 {
576 #ifndef GMX_MPI
577     gmx_call("gmx_sumli_sim");
578 #else
579 #if defined(MPI_IN_PLACE_EXISTS) || defined(GMX_THREADS)
580     MPI_Allreduce(MPI_IN_PLACE,r,nr,GMX_MPI_LARGE_INT,MPI_SUM,
581                   ms->mpi_comm_masters);
582 #else
583     /* this is thread-unsafe, but it will do for now: */
584     int i;
585
586     if (nr > ms->mpb->ibuf_alloc) {
587         ms->mpb->ibuf_alloc = nr;
588         srenew(ms->mpb->ibuf,ms->mpb->ibuf_alloc);
589     }
590     MPI_Allreduce(r,ms->mpb->ibuf,nr,GMX_MPI_LARGE_INT,MPI_SUM,
591                   ms->mpi_comm_masters);
592     for(i=0; i<nr; i++)
593         r[i] = ms->mpb->ibuf[i];
594 #endif
595 #endif
596 }
597
598
599 void gmx_finalize(void)
600 {
601 #ifndef GMX_MPI
602   gmx_call("gmx_finalize");
603 #else
604   int ret;
605
606   /* just as a check; we don't want to finalize twice */
607   int finalized;
608   MPI_Finalized(&finalized);
609   if (finalized)
610       return;
611
612   /* We sync the processes here to try to avoid problems
613    * with buggy MPI implementations that could cause
614    * unfinished processes to terminate.
615    */
616   MPI_Barrier(MPI_COMM_WORLD);
617
618   /*
619   if (DOMAINDECOMP(cr)) {
620     if (cr->npmenodes > 0 || cr->dd->bCartesian) 
621       MPI_Comm_free(&cr->mpi_comm_mygroup);
622     if (cr->dd->bCartesian)
623       MPI_Comm_free(&cr->mpi_comm_mysim);
624   }
625   */
626
627   /* Apparently certain mpich implementations cause problems
628    * with MPI_Finalize. In that case comment out MPI_Finalize.
629    */
630   if (debug)
631     fprintf(debug,"Will call MPI_Finalize now\n");
632
633   ret = MPI_Finalize();
634   if (debug)
635     fprintf(debug,"Return code from MPI_Finalize = %d\n",ret);
636 #endif
637 }
638