Redefine the default boolean type to gmx_bool.
[alexxy/gromacs.git] / src / gmxlib / splitter.c
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
2  *
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.
15
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.
20  * 
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.
27  * 
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.
30  * 
31  * For more info, check our website at http://www.gromacs.org
32  * 
33  * And Hey:
34  * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
35  */
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
39
40 #include <stdio.h>
41 #include <string.h>
42
43 #include "sysstuff.h"
44 #include "macros.h"
45 #include "smalloc.h"
46 #include "typedefs.h"
47 #include "mshift.h"
48 #include "invblock.h"
49 #include "txtdump.h"
50 #include "math.h"
51 #include "assert.h"
52 #include "gmx_fatal.h"
53 #include "splitter.h"
54
55 typedef struct {
56   int     nr;
57   t_iatom *ia;
58 } t_sf;
59
60 static t_sf *init_sf(int nr)
61 {
62   t_sf *sf;
63   int  i;
64
65   snew(sf,nr);
66   for(i=0; (i<nr); i++) {
67     sf[i].nr=0;
68     sf[i].ia=NULL;
69   }
70
71   return sf;
72 }
73
74 static void done_sf(int nr,t_sf *sf)
75 {
76   int i;
77
78   for(i=0; (i<nr); i++) {
79     sf[i].nr=0;
80     sfree(sf[i].ia);
81     sf[i].ia=NULL;
82   }
83   sfree(sf);
84 }
85
86 static void push_sf(t_sf *sf,int nr,t_iatom ia[])
87 {
88   int i;
89
90   srenew(sf->ia,sf->nr+nr);
91   for(i=0; (i<nr); i++)
92     sf->ia[sf->nr+i]=ia[i];
93   sf->nr+=nr;
94 }
95
96 static int min_nodeid(int nr,atom_id list[],int hid[])
97 {
98   int i,nodeid,minnodeid;
99
100   if (nr <= 0)
101     gmx_incons("Invalid node number");
102   minnodeid=hid[list[0]];
103   for (i=1; (i<nr); i++) 
104     if ((nodeid=hid[list[i]]) < minnodeid) 
105       minnodeid=nodeid;
106
107   return minnodeid;
108 }
109
110
111 static void split_force2(t_inputrec *ir, int nnodes,int hid[],int ftype,t_ilist *ilist,
112                                                  int *multinr,
113                                                  int *constr_min_nodeid,int * constr_max_nodeid,
114                                                  int *left_range, int *right_range)
115 {
116   int     i,j,k,type,nodeid,nratoms,tnr;
117   int     nvsite_constr;
118   t_iatom ai,aj;
119   int     node_low_ai,node_low_aj,node_high_ai,node_high_aj;
120   int     node_low,node_high;
121   int     nodei,nodej;
122   t_sf    *sf;
123   int     nextra;
124         
125   sf = init_sf(nnodes);
126
127   node_high = node_low = 0;
128   nextra    = 0;
129         
130   /* Walk along all the bonded forces, find the appropriate node
131    * to calc it on, and add it to that nodes list.
132    */
133   for (i=0; i<ilist->nr; i+=(1+nratoms)) 
134   {
135     type    = ilist->iatoms[i];
136     nratoms = interaction_function[ftype].nratoms;
137
138     if (ftype == F_CONSTR) 
139         {
140                 ai = ilist->iatoms[i+1];
141                 aj = ilist->iatoms[i+2];
142                 
143                 nodei = hid[ai];
144                 nodej = hid[aj];
145                 nodeid = nodei;
146
147                 if(ir->eConstrAlg == econtLINCS)
148                 {
149                         node_low_ai   = constr_min_nodeid[ai];
150                         node_low_aj   = constr_min_nodeid[aj];
151                         node_high_ai  = constr_max_nodeid[ai];
152                         node_high_aj  = constr_max_nodeid[aj];
153
154                         node_low      = min(node_low_ai,node_low_aj);
155                         node_high     = max(node_high_ai,node_high_aj);
156                         
157             if (node_high-nodei > 1 || nodei-node_low > 1 ||
158                 node_high-nodej > 1 || nodej-node_low > 1 )
159                         {
160                                 gmx_fatal(FARGS,"Constraint dependencies further away than next-neighbor\n"
161                                                   "in particle decomposition. Constraint between atoms %d--%d evaluated\n"
162                                                   "on node %d and %d, but atom %d has connections within %d bonds (lincs_order)\n"
163                                                   "of node %d, and atom %d has connections within %d bonds of node %d.\n"
164                                                   "Reduce the # nodes, lincs_order, or\n"
165                                                   "try domain decomposition.",ai,aj,nodei,nodej,ai,ir->nProjOrder,node_low,aj,ir->nProjOrder,node_high);
166                         }
167                         
168                         if (node_low < nodei || node_low < nodej)
169                         {
170                                 right_range[node_low] = max(right_range[node_low],aj);
171                         }
172                         if (node_high > nodei || node_high > nodej)
173                         {
174                                 left_range[node_high] = min(left_range[node_high],ai);
175                         }
176                 }
177                 else 
178                 {
179                         /* Shake */
180                         if (hid[ilist->iatoms[i+2]] != nodei) 
181                                 gmx_fatal(FARGS,"Shake block crossing node boundaries\n"
182                                                   "constraint between atoms (%d,%d) (try LINCS instead!)",
183                                                   ilist->iatoms[i+1]+1,ilist->iatoms[i+2]+1);
184                 }
185     }
186     else if (ftype == F_SETTLE) 
187         {
188       /* Only the first particle is stored for settles ... */
189       ai=ilist->iatoms[i+1];
190       nodeid=hid[ai];
191       if ((nodeid != hid[ai+1]) ||
192           (nodeid != hid[ai+2]))
193         gmx_fatal(FARGS,"Settle block crossing node boundaries\n"
194                   "constraint between atoms (%d-%d)",ai+1,ai+2+1);
195     }
196     else if(interaction_function[ftype].flags & IF_VSITE) 
197         {
198       /* Virtual sites are special, since we need to pre-communicate
199            * their coordinates to construct vsites before then main
200            * coordinate communication. 
201            * Vsites can have constructing atoms both larger and smaller than themselves.
202            * To minimize communication and book-keeping, each vsite is constructed on
203            * the home node of the atomnr of the vsite.
204            * Since the vsite coordinates too have to be communicated to the next node,
205            * we need to
206            *
207            * 1. Pre-communicate coordinates of constructing atoms
208            * 2. Construct the vsite
209            * 3. Perform main coordinate communication
210            *
211            * Note that this has change from gromacs 4.0 and earlier, where the vsite
212            * was constructed on the home node of the lowest index of any of the constructing
213            * atoms and the vsite itself.
214            */
215                 
216                 if (ftype==F_VSITE2)
217                         nvsite_constr=2;
218                 else if(ftype==F_VSITE4FD || ftype==F_VSITE4FDN)
219                         nvsite_constr=4;
220                 else
221                         nvsite_constr=3;
222                 
223                 /* Vsites are constructed on the home node of the actual site to save communication 
224                  * and simplify the book-keeping.
225                  */
226                 nodeid=hid[ilist->iatoms[i+1]];
227                 
228                 for(k=2;k<nvsite_constr+2;k++) 
229                 {
230                         if(hid[ilist->iatoms[i+k]]<(nodeid-1) ||
231                            hid[ilist->iatoms[i+k]]>(nodeid+1))
232                                 gmx_fatal(FARGS,"Virtual site %d and its constructing"
233                                                   " atoms are not on the same or adjacent\n" 
234                                                   " nodes. This is necessary to avoid a lot\n"
235                                                   " of extra communication. The easiest way"
236                                                   " to ensure this is to place virtual sites\n"
237                                                   " close to the constructing atoms.\n"
238                                                   " Sorry, but you will have to rework your topology!\n",
239                                                   ilist->iatoms[i+1]);
240                 }
241     } 
242         else
243     {
244                 nodeid=min_nodeid(nratoms,&ilist->iatoms[i+1],hid);
245         }
246           
247         if (ftype == F_CONSTR && ir->eConstrAlg == econtLINCS)
248         {
249                 push_sf(&(sf[nodeid]),nratoms+1,&(ilist->iatoms[i]));
250         
251                 if(node_low<nodeid)
252                 {
253                         push_sf(&(sf[node_low]),nratoms+1,&(ilist->iatoms[i]));
254                         nextra+=nratoms+1;
255                 }
256                 if(node_high>nodeid)
257                 {
258                         push_sf(&(sf[node_high]),nratoms+1,&(ilist->iatoms[i]));
259                         nextra+=nratoms+1;
260                 }
261         }
262         else
263         {
264                 push_sf(&(sf[nodeid]),nratoms+1,&(ilist->iatoms[i]));
265         }
266   }
267
268   if(nextra>0)
269   {
270           ilist->nr += nextra;
271           srenew(ilist->iatoms,ilist->nr);
272   }
273         
274   tnr=0;        
275   for(nodeid=0; (nodeid<nnodes); nodeid++) {
276     for (i=0; (i<sf[nodeid].nr); i++) 
277       ilist->iatoms[tnr++]=sf[nodeid].ia[i];
278
279     multinr[nodeid]=(nodeid==0) ? 0 : multinr[nodeid-1];
280     multinr[nodeid]+=sf[nodeid].nr;       
281   }
282         
283         if (tnr != ilist->nr)
284                 gmx_incons("Splitting forces over processors");
285
286   done_sf(nnodes,sf);
287 }
288
289 static int *home_index(int nnodes,t_block *cgs,int *multinr)
290 {
291   /* This routine determines the node id for each particle */
292   int *hid;
293   int nodeid,j0,j1,j,k;
294   
295   snew(hid,cgs->index[cgs->nr]);
296   /* Initiate to -1 to make it possible to check afterwards,
297    * all hid's should be set in the loop below
298    */
299   for(k=0; (k<cgs->index[cgs->nr]); k++)
300     hid[k]=-1;
301     
302   /* loop over nodes */
303   for(nodeid=0; (nodeid<nnodes); nodeid++) {
304     j0 = (nodeid==0) ? 0 : multinr[nodeid-1];
305     j1 = multinr[nodeid];
306     
307     /* j0 and j1 are the boundariesin the index array */
308     for(j=j0; (j<j1); j++) {
309       for(k=cgs->index[j]; (k<cgs->index[j+1]); k++) {
310         hid[k]=nodeid;
311       }
312     }
313   }
314   /* Now verify that all hid's are not -1 */
315   for(k=0; (k<cgs->index[cgs->nr]); k++)
316     if (hid[k] == -1)
317       gmx_fatal(FARGS,"hid[%d] = -1, cgs->nr = %d, natoms = %d",
318                   k,cgs->nr,cgs->index[cgs->nr]);
319   
320   return hid;
321 }
322
323 typedef struct {
324   int atom,ic,is;
325 } t_border;
326
327 void set_bor(t_border *b,int atom,int ic,int is)
328 {
329   if (debug)
330     fprintf(debug,"border @ atom %5d [ ic = %5d,  is = %5d ]\n",atom,ic,is);
331   b->atom = atom;
332   b->ic   = ic;
333   b->is   = is;
334 }
335
336 static gmx_bool is_bor(atom_id ai[],int i)
337 {
338   return ((ai[i] != ai[i-1]) || ((ai[i] == NO_ATID) && (ai[i-1] == NO_ATID)));
339 }
340
341 static t_border *mk_border(FILE *fp,int natom,atom_id *invcgs,
342                            atom_id *invshk,int *nb)
343 {
344   t_border *bor;
345   atom_id  *sbor,*cbor;
346   int      i,j,is,ic,ns,nc,nbor;
347
348   if (debug) {
349     for(i=0; (i<natom); i++) {
350       fprintf(debug,"atom: %6d  cgindex: %6d  shkindex: %6d\n",
351                           i, invcgs[i], invshk[i]);
352     }
353   }
354     
355   snew(sbor,natom+1);
356   snew(cbor,natom+1);
357   ns = nc = 1;
358   for(i=1; (i<natom); i++) {
359     if (is_bor(invcgs,i)) 
360       cbor[nc++] = i;
361     if (is_bor(invshk,i)) 
362       sbor[ns++] = i;
363   }
364   sbor[ns] = 0;
365   cbor[nc] = 0;
366   if (fp)
367   {
368           fprintf(fp,"There are %d charge group borders",nc);
369           if(invshk!=NULL)
370           {
371                   fprintf(fp," and %d shake borders",ns);
372           }
373           fprintf(fp,".\n");
374   }
375   snew(bor,max(nc,ns));
376   ic = is = nbor = 0;
377   while ((ic < nc) || (is < ns)) {
378     if (sbor[is] == cbor[ic]) {
379       set_bor(&(bor[nbor]),cbor[ic],ic,is);
380       nbor++;
381       if (ic < nc) ic++;
382       if (is < ns) is++;
383     }
384     else if (cbor[ic] > sbor[is]) {
385       if (is == ns) {
386         set_bor(&(bor[nbor]),cbor[ic],ic,is);
387         nbor++;
388         if (ic < nc) ic++;
389       }
390       else if (is < ns) 
391         is++;
392     }
393     else if (ic < nc)
394       ic++;
395     else
396       is++;/*gmx_fatal(FARGS,"Can't happen is=%d, ic=%d (%s, %d)",
397              is,ic,__FILE__,__LINE__);*/
398   }
399   if (fp)
400     fprintf(fp,"There are %d total borders\n",nbor);
401
402   if (debug) {
403     fprintf(debug,"There are %d actual bor entries\n",nbor);
404     for(i=0; (i<nbor); i++) 
405       fprintf(debug,"bor[%5d] = atom: %d  ic: %d  is: %d\n",i,
406               bor[i].atom,bor[i].ic,bor[i].is);
407   }
408   
409   *nb = nbor;
410   
411   return bor;
412 }
413
414 static void split_blocks(FILE *fp,t_inputrec *ir, int nnodes,
415                          t_block *cgs,t_blocka *sblock,real capacity[],
416                          int *multinr_cgs)
417 {
418   int      natoms,*maxatom;
419   int      i,ii,ai,b0,b1;
420   int      nodeid,last_shk,nbor;
421   t_border *border;
422   double   tload,tcap;
423   
424   gmx_bool    bSHK;
425   atom_id *shknum,*cgsnum;
426   
427   natoms = cgs->index[cgs->nr];
428         
429   if (NULL != debug) {
430     pr_block(debug,0,"cgs",cgs,TRUE);
431     pr_blocka(debug,0,"sblock",sblock,TRUE);
432     fflush(debug);
433   }
434
435   cgsnum = make_invblock(cgs,natoms+1); 
436   shknum = make_invblocka(sblock,natoms+1);
437   border = mk_border(fp,natoms,cgsnum,shknum,&nbor);
438
439   snew(maxatom,nnodes);
440   tload  = capacity[0]*natoms;
441   tcap   = 1.0;
442   nodeid = 0;
443   /* Start at bor is 1, to force the first block on the first processor */
444   for(i=0; (i<nbor) && (tload < natoms); i++) {
445     if(i<(nbor-1)) 
446       b1=border[i+1].atom;
447     else
448       b1=natoms;
449
450     b0 = border[i].atom;
451     
452     if ((fabs(b0-tload)<fabs(b1-tload))) {
453       /* New nodeid time */
454       multinr_cgs[nodeid] = border[i].ic;
455       maxatom[nodeid]     = b0;
456       tcap -= capacity[nodeid];
457       nodeid++;
458       
459       /* Recompute target load */
460       tload = b0 + (natoms-b0)*capacity[nodeid]/tcap;
461
462     if(debug)   
463           printf("tload: %g tcap: %g  nodeid: %d\n",tload,tcap,nodeid);
464     } 
465   }
466   /* Now the last one... */
467   while (nodeid < nnodes) {
468     multinr_cgs[nodeid] = cgs->nr;
469     /* Store atom number, see above */
470     maxatom[nodeid]     = natoms;
471     nodeid++;
472   }
473   if (nodeid != nnodes) {
474     gmx_fatal(FARGS,"nodeid = %d, nnodes = %d, file %s, line %d",
475                 nodeid,nnodes,__FILE__,__LINE__);
476   }
477
478   for(i=nnodes-1; (i>0); i--)
479     maxatom[i]-=maxatom[i-1];
480
481   if (fp) {
482     fprintf(fp,"Division over nodes in atoms:\n");
483     for(i=0; (i<nnodes); i++)
484       fprintf(fp," %7d",maxatom[i]);
485     fprintf(fp,"\n");
486   }
487         
488   sfree(maxatom);
489   sfree(shknum);
490   sfree(cgsnum);
491   sfree(border);
492 }
493
494 typedef struct {
495   int atom,sid;
496 } t_sid;
497
498 static int sid_comp(const void *a,const void *b)
499 {
500   t_sid *sa,*sb;
501   int   dd;
502   
503   sa=(t_sid *)a;
504   sb=(t_sid *)b;
505   
506   dd=sa->sid-sb->sid;
507   if (dd == 0)
508     return (sa->atom-sb->atom);
509   else
510     return dd;
511 }
512
513 static int mk_grey(int nnodes,egCol egc[],t_graph *g,int *AtomI,
514                    int maxsid,t_sid sid[])
515 {
516   int  j,ng,ai,aj,g0;
517
518   ng=0;
519   ai=*AtomI;
520   
521   g0=g->start;
522   /* Loop over all the bonds */
523   for(j=0; (j<g->nedge[ai]); j++) {
524     aj=g->edge[ai][j]-g0;
525     /* If there is a white one, make it gray and set pbc */
526     if (egc[aj] == egcolWhite) {
527       if (aj < *AtomI)
528         *AtomI = aj;
529       egc[aj] = egcolGrey;
530
531       /* Check whether this one has been set before... */
532       range_check(aj+g0,0,maxsid);
533       range_check(ai+g0,0,maxsid);
534       if (sid[aj+g0].sid != -1) 
535         gmx_fatal(FARGS,"sid[%d]=%d, sid[%d]=%d, file %s, line %d",
536                     ai,sid[ai+g0].sid,aj,sid[aj+g0].sid,__FILE__,__LINE__);
537       else {
538         sid[aj+g0].sid  = sid[ai+g0].sid;
539         sid[aj+g0].atom = aj+g0;
540       }
541       ng++;
542     }
543   }
544   return ng;
545 }
546
547 static int first_colour(int fC,egCol Col,t_graph *g,egCol egc[])
548 /* Return the first node with colour Col starting at fC.
549  * return -1 if none found.
550  */
551 {
552   int i;
553   
554   for(i=fC; (i<g->nnodes); i++)
555     if ((g->nedge[i] > 0) && (egc[i]==Col))
556       return i;
557   
558   return -1;
559 }
560
561 static int mk_sblocks(FILE *fp,t_graph *g,int maxsid,t_sid sid[])
562 {
563   int    ng,nnodes;
564   int    nW,nG,nB;              /* Number of Grey, Black, White */
565   int    fW,fG;                 /* First of each category       */
566   egCol  *egc=NULL;             /* The colour of each node      */
567   int    g0,nblock;
568   
569   if (!g->nbound) 
570     return 0;
571
572   nblock=0;
573   
574   nnodes=g->nnodes;
575   snew(egc,nnodes);
576   
577   g0=g->start;
578   nW=g->nbound;
579   nG=0;
580   nB=0;
581
582   fW=0;
583
584   /* We even have a loop invariant:
585    * nW+nG+nB == g->nbound
586    */
587   
588   if (fp)
589     fprintf(fp,"Walking down the molecule graph to make constraint-blocks\n");
590
591   while (nW > 0) {
592     /* Find the first white, this will allways be a larger
593      * number than before, because no nodes are made white
594      * in the loop
595      */
596     if ((fW=first_colour(fW,egcolWhite,g,egc)) == -1) 
597       gmx_fatal(FARGS,"No WHITE nodes found while nW=%d\n",nW);
598     
599     /* Make the first white node grey, and set the block number */
600     egc[fW]        = egcolGrey;
601     range_check(fW+g0,0,maxsid);
602     sid[fW+g0].sid = nblock++;
603     nG++;
604     nW--;
605
606     /* Initial value for the first grey */
607     fG=fW;
608
609     if (debug)
610       fprintf(debug,"Starting G loop (nW=%d, nG=%d, nB=%d, total %d)\n",
611               nW,nG,nB,nW+nG+nB);
612               
613     while (nG > 0) {
614       if ((fG=first_colour(fG,egcolGrey,g,egc)) == -1)
615         gmx_fatal(FARGS,"No GREY nodes found while nG=%d\n",nG);
616       
617       /* Make the first grey node black */
618       egc[fG]=egcolBlack;
619       nB++;
620       nG--;
621
622       /* Make all the neighbours of this black node grey
623        * and set their block number
624        */
625       ng=mk_grey(nnodes,egc,g,&fG,maxsid,sid);
626       /* ng is the number of white nodes made grey */
627       nG+=ng;
628       nW-=ng;
629     }
630   }
631   sfree(egc);
632
633   if (debug)
634     fprintf(debug,"Found %d shake blocks\n",nblock);
635     
636   return nblock;
637 }
638
639
640 typedef struct {
641   int first,last,sid;
642 } t_merge_sid;
643
644 static int ms_comp(const void *a, const void *b)
645 {
646   t_merge_sid *ma = (t_merge_sid *)a;
647   t_merge_sid *mb = (t_merge_sid *)b;
648   int d;
649   
650   d = ma->first-mb->first;
651   if (d == 0)
652     return ma->last-mb->last;
653   else
654     return d;
655 }
656
657 static int merge_sid(int i0,int at_start,int at_end,int nsid,t_sid sid[],
658                      t_blocka *sblock)
659 {
660   int  i,j,k,n,isid,ndel;
661   t_merge_sid *ms;
662   int  nChanged;
663
664   /* We try to remdy the following problem:
665    * Atom: 1  2  3  4  5 6 7 8 9 10
666    * Sid:  0 -1  1  0 -1 1 1 1 1 1
667    */
668    
669   /* Determine first and last atom in each shake ID */
670   snew(ms,nsid);
671
672   for(k=0; (k<nsid); k++) {
673     ms[k].first = at_end+1;
674     ms[k].last  = -1;
675     ms[k].sid   = k;
676   }
677   for(i=at_start; (i<at_end); i++) {
678     isid = sid[i].sid;
679     range_check(isid,-1,nsid);
680     if (isid >= 0) {
681       ms[isid].first = min(ms[isid].first,sid[i].atom);
682       ms[isid].last  = max(ms[isid].last,sid[i].atom);
683     }
684   }
685   qsort(ms,nsid,sizeof(ms[0]),ms_comp);
686   
687   /* Now merge the overlapping ones */
688   ndel = 0;
689   for(k=0; (k<nsid); ) {
690     for(j=k+1; (j<nsid); ) {
691       if (ms[j].first <= ms[k].last) {
692         ms[k].last  = max(ms[k].last,ms[j].last);
693         ms[k].first = min(ms[k].first,ms[j].first);
694         ms[j].sid = -1;
695         ndel++;
696         j++;
697       }
698       else {
699         k = j;
700         j = k+1;
701       }
702     }
703     if (j == nsid)
704       k++;
705   }
706   for(k=0; (k<nsid); k++) {
707     while ((k < nsid-1) && (ms[k].sid == -1)) {
708       for(j=k+1; (j<nsid); j++) {
709         memcpy(&(ms[j-1]),&(ms[j]),sizeof(ms[0]));
710       }
711       nsid--;
712     }
713   }
714
715   for(k=at_start; (k<at_end); k++) {
716     sid[k].atom = k;
717     sid[k].sid = -1;
718   }
719   sblock->nr  = nsid;
720   sblock->nalloc_index = sblock->nr+1;
721   snew(sblock->index,sblock->nalloc_index);
722   sblock->nra = at_end - at_start;
723   sblock->nalloc_a = sblock->nra;
724   snew(sblock->a,sblock->nalloc_a);
725   sblock->index[0] = 0;
726   for(k=n=0; (k<nsid); k++) {
727     sblock->index[k+1] = sblock->index[k] + ms[k].last - ms[k].first+1;
728     for(j=ms[k].first; (j<=ms[k].last); j++) {
729       range_check(n,0,sblock->nra);
730       sblock->a[n++] = j;
731       range_check(j,0,at_end);
732       if (sid[j].sid == -1)
733         sid[j].sid = k;
734       else 
735         fprintf(stderr,"Double sids (%d, %d) for atom %d\n",sid[j].sid,k,j);
736     }
737   }
738   assert(k == nsid);
739   /* Removed 2007-09-04
740   sblock->index[k+1] = natoms;
741   for(k=0; (k<natoms); k++)
742     if (sid[k].sid == -1)
743       sblock->a[n++] = k;
744   assert(n == natoms);
745   */
746   sblock->nra = n;
747   assert(sblock->index[k] == sblock->nra);
748   sfree(ms);
749   
750   return nsid;
751 }
752
753 void gen_sblocks(FILE *fp,int at_start,int at_end,
754                                  t_idef *idef,t_blocka *sblock,
755                                  gmx_bool bSettle)
756 {
757   t_graph *g;
758   int     i,i0,j,k,istart,n;
759   t_sid   *sid;
760   int     isid,nsid;
761   
762   g=mk_graph(NULL,idef,at_start,at_end,TRUE,bSettle);
763   if (debug)
764     p_graph(debug,"Graaf Dracula",g);
765   snew(sid,at_end);
766   for(i=at_start; (i<at_end); i++) {
767     sid[i].atom =  i;
768     sid[i].sid  = -1;
769   }
770   nsid=mk_sblocks(fp,g,at_end,sid);
771
772   if (!nsid)
773     return;
774     
775   /* Now sort the shake blocks... */
776   qsort(sid+at_start,at_end-at_start,(size_t)sizeof(sid[0]),sid_comp);
777   
778   if (debug) {
779     fprintf(debug,"Sorted shake block\n");
780     for(i=at_start; (i<at_end); i++) 
781       fprintf(debug,"sid[%5d] = atom:%5d sid:%5d\n",i,sid[i].atom,sid[i].sid);
782   }
783   /* Now check how many are NOT -1, i.e. how many have to be shaken */
784   for(i0=at_start; (i0<at_end); i0++)
785     if (sid[i0].sid > -1)
786       break;
787   
788   /* Now we have the sids that have to be shaken. We'll check the min and
789    * max atom numbers and this determines the shake block. DvdS 2007-07-19.
790    * For the purpose of making boundaries all atoms in between need to be
791    * part of the shake block too. There may be cases where blocks overlap
792    * and they will have to be merged.
793    */
794   nsid = merge_sid(i0,at_start,at_end,nsid,sid,sblock);
795   /* Now sort the shake blocks again... */
796   /*qsort(sid,natoms,(size_t)sizeof(sid[0]),sid_comp);*/
797   
798   /* Fill the sblock struct */    
799   /*  sblock->nr  = nsid;
800   sblock->nra = natoms;
801   srenew(sblock->a,sblock->nra);
802   srenew(sblock->index,sblock->nr+1);
803   
804   i    = i0;
805   isid = sid[i].sid;
806   n    = k = 0;
807   sblock->index[n++]=k;
808   while (i < natoms) {
809     istart = sid[i].atom;
810     while ((i<natoms-1) && (sid[i+1].sid == isid))
811     i++;*/
812     /* After while: we found a new block, or are thru with the atoms */
813   /*    for(j=istart; (j<=sid[i].atom); j++,k++)
814       sblock->a[k]=j;
815     sblock->index[n] = k;
816     if (i < natoms-1)
817       n++;
818     if (n > nsid)
819       gmx_fatal(FARGS,"Death Horror: nsid = %d, n= %d",nsid,n);
820     i++;
821     isid = sid[i].sid;
822   }
823   */
824   sfree(sid);
825   /* Due to unknown reason this free generates a problem sometimes */
826   done_graph(g);
827   sfree(g); 
828   if (debug)
829     fprintf(debug,"Done gen_sblocks\n");
830 }
831
832 static t_blocka block2blocka(t_block *block)
833 {
834   t_blocka blocka;
835   int i;
836
837   blocka.nr = block->nr;
838   blocka.nalloc_index = blocka.nr + 1;
839   snew(blocka.index,blocka.nalloc_index);
840   for(i=0; i<=block->nr; i++)
841     blocka.index[i] = block->index[i];
842   blocka.nra = block->index[block->nr];
843   blocka.nalloc_a = blocka.nra;
844   snew(blocka.a,blocka.nalloc_a);
845   for(i=0; i<blocka.nra; i++)
846     blocka.a[i] = i;
847
848   return blocka;
849 }
850
851 typedef struct
852 {
853         int   nconstr;
854         int   index[10];
855 } pd_constraintlist_t;
856         
857
858 static void
859 find_constraint_range_recursive(pd_constraintlist_t *  constraintlist,
860                                                                 int                    thisatom,
861                                                                 int                    depth,
862                                                                 int *                  min_atomid,
863                                                                 int *                  max_atomid)
864 {
865         int i,j;
866         int nconstr;
867         
868         for(i=0;i<constraintlist[thisatom].nconstr;i++)
869         {
870                 j = constraintlist[thisatom].index[i];
871
872                 *min_atomid = (j<*min_atomid) ? j : *min_atomid;
873                 *max_atomid = (j>*max_atomid) ? j : *max_atomid;
874
875                 if(depth>0)
876                 {
877                         find_constraint_range_recursive(constraintlist,j,depth-1,min_atomid,max_atomid);
878                 }
879         }
880 }
881
882 static void
883 pd_determine_constraints_range(t_inputrec *      ir,
884                                                            int               natoms,
885                                                            t_ilist *         ilist,
886                                                            int               hid[],
887                                                            int *             min_nodeid,
888                                                            int *             max_nodeid)
889 {
890         int i,j,k;
891         int nratoms;
892         int depth;
893         int ai,aj;
894         int min_atomid,max_atomid;
895         pd_constraintlist_t *constraintlist;
896         
897         nratoms = interaction_function[F_CONSTR].nratoms;
898         depth   = ir->nProjOrder;
899
900         snew(constraintlist,natoms);
901         
902         /* Make a list of all the connections */
903         for(i=0;i<ilist->nr;i+=nratoms+1)
904         {
905                 ai=ilist->iatoms[i+1];
906                 aj=ilist->iatoms[i+2];
907                 constraintlist[ai].index[constraintlist[ai].nconstr++]=aj;
908                 constraintlist[aj].index[constraintlist[aj].nconstr++]=ai;
909         }
910                            
911         for(i=0;i<natoms;i++)
912         {
913                 min_atomid = i;
914                 max_atomid = i;
915                 
916                 find_constraint_range_recursive(constraintlist,i,depth,&min_atomid,&max_atomid);
917                 
918                 min_nodeid[i] = hid[min_atomid];
919                 max_nodeid[i] = hid[max_atomid];                
920         }
921         sfree(constraintlist);
922 }
923
924
925 void split_top(FILE *fp,int nnodes,gmx_localtop_t *top,t_inputrec *ir,t_block *mols,
926                real *capacity,int *multinr_cgs,int **multinr_nre, int *left_range, int * right_range)
927 {
928   int     natoms,i,j,k,mj,atom,maxatom,sstart,send,bstart,nodeid;
929   t_blocka sblock;
930   int     *homeind;
931   int ftype,nvsite_constr,nra,nrd;
932   t_iatom   *ia;
933   int minhome,ihome,minidx;
934         int *constr_min_nodeid;
935         int *constr_max_nodeid;
936         
937   if (nnodes <= 1)
938     return;
939   
940   natoms = mols->index[mols->nr];
941
942   if (fp)
943     fprintf(fp,"splitting topology...\n");
944   
945 /*#define MOL_BORDER*/
946 /*Removed the above to allow splitting molecules with h-bond constraints
947   over processors. The results in DP are the same. */
948   init_blocka(&sblock);
949   if(ir->eConstrAlg != econtLINCS)
950   {
951 #ifndef MOL_BORDER
952   /* Make a special shake block that includes settles */
953           gen_sblocks(fp,0,natoms,&top->idef,&sblock,TRUE);             
954 #else
955           sblock = block2blocka(mols);
956 #endif  
957   }
958         
959   split_blocks(fp,ir,nnodes,&top->cgs,&sblock,capacity,multinr_cgs);
960         
961   homeind=home_index(nnodes,&top->cgs,multinr_cgs);
962
963   snew(constr_min_nodeid,natoms);
964   snew(constr_max_nodeid,natoms);
965
966   if(top->idef.il[F_CONSTR].nr>0)
967   {
968           pd_determine_constraints_range(ir,natoms,&top->idef.il[F_CONSTR],homeind,constr_min_nodeid,constr_max_nodeid);
969   }
970   else
971   {
972           /* Not 100% necessary, but it is a bad habit to have uninitialized arrays around... */
973           for(i=0;i<natoms;i++)
974           {
975                   constr_min_nodeid[i] = constr_max_nodeid[i] = homeind[i];  
976           }
977   }
978         
979   /* Default limits (no communication) for PD constraints */
980   left_range[0] = 0;
981   for(i=1;i<nnodes;i++)
982   {
983         left_range[i]    = top->cgs.index[multinr_cgs[i-1]];
984         right_range[i-1] = left_range[i]-1;
985   }
986   right_range[nnodes-1] = top->cgs.index[multinr_cgs[nnodes-1]]-1;
987         
988   for(j=0; (j<F_NRE); j++)
989     split_force2(ir, nnodes,homeind,j,&top->idef.il[j],multinr_nre[j],constr_min_nodeid,constr_max_nodeid,
990                                  left_range,right_range);
991
992   sfree(constr_min_nodeid);
993   sfree(constr_max_nodeid);
994         
995   sfree(homeind);
996   done_blocka(&sblock);
997 }
998