Enforced rotation
[alexxy/gromacs.git] / src / gmxlib / mvdata.c
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
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 /* This file is completely threadsafe - keep it that way! */
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
39
40 #include <sysstuff.h>
41 #include <string.h>
42 #include "typedefs.h"
43 #include "main.h"
44 #include "mvdata.h"
45 #include "network.h"
46 #include "smalloc.h"
47 #include "gmx_fatal.h"
48 #include "symtab.h"
49 #include "vec.h"
50 #include "tgroup.h"
51
52 #define   block_bc(cr,   d) gmx_bcast(     sizeof(d),     &(d),(cr))
53 /* Probably the test for (nr) > 0 in the next macro is only needed
54  * on BlueGene(/L), where IBM's MPI_Bcast will segfault after
55  * dereferencing a null pointer, even when no data is to be transferred. */
56 #define  nblock_bc(cr,nr,d) { if ((nr) > 0) gmx_bcast((nr)*sizeof((d)[0]), (d),(cr)); }
57 #define    snew_bc(cr,d,nr) { if (!MASTER(cr)) snew((d),(nr)); }
58 /* Dirty macro with bAlloc not as an argument */
59 #define nblock_abc(cr,nr,d) { if (bAlloc) snew((d),(nr)); nblock_bc(cr,(nr),(d)); }
60
61 static void bc_string(const t_commrec *cr,t_symtab *symtab,char ***s)
62 {
63   int handle;
64   
65   if (MASTER(cr)) {
66     handle = lookup_symtab(symtab,*s);
67   }
68   block_bc(cr,handle);
69   if (!MASTER(cr)) {
70     *s = get_symtab_handle(symtab,handle);
71   }
72 }
73
74 static void bc_strings(const t_commrec *cr,t_symtab *symtab,int nr,char ****nm)
75 {
76   int  i;
77   int  *handle;
78   char ***NM;
79
80   snew(handle,nr);
81   if (MASTER(cr)) {
82     NM = *nm;
83     for(i=0; (i<nr); i++)
84       handle[i] = lookup_symtab(symtab,NM[i]);
85   }
86   nblock_bc(cr,nr,handle);
87
88   if (!MASTER(cr)) {
89     snew_bc(cr,*nm,nr);
90     NM = *nm;
91     for (i=0; (i<nr); i++) 
92       (*nm)[i] = get_symtab_handle(symtab,handle[i]);
93   }
94   sfree(handle);
95 }
96
97 static void bc_strings_resinfo(const t_commrec *cr,t_symtab *symtab,
98                                int nr,t_resinfo *resinfo)
99 {
100   int  i;
101   int  *handle;
102
103   snew(handle,nr);
104   if (MASTER(cr)) {
105     for(i=0; (i<nr); i++)
106       handle[i] = lookup_symtab(symtab,resinfo[i].name);
107   }
108   nblock_bc(cr,nr,handle);
109
110   if (!MASTER(cr)) {
111     for (i=0; (i<nr); i++) 
112       resinfo[i].name = get_symtab_handle(symtab,handle[i]);
113   }
114   sfree(handle);
115 }
116
117 static void bc_symtab(const t_commrec *cr,t_symtab *symtab)
118 {
119   int i,nr,len;
120   t_symbuf *symbuf;
121
122   block_bc(cr,symtab->nr);
123   nr = symtab->nr;
124   snew_bc(cr,symtab->symbuf,1);
125   symbuf = symtab->symbuf;
126   symbuf->bufsize = nr;
127   snew_bc(cr,symbuf->buf,nr);
128   for (i=0; i<nr; i++) {
129     if (MASTER(cr))
130       len = strlen(symbuf->buf[i]) + 1;
131     block_bc(cr,len);
132     snew_bc(cr,symbuf->buf[i],len);
133     nblock_bc(cr,len,symbuf->buf[i]);
134   }
135 }
136
137 static void bc_block(const t_commrec *cr,t_block *block)
138 {
139   block_bc(cr,block->nr);
140   snew_bc(cr,block->index,block->nr+1);
141   nblock_bc(cr,block->nr+1,block->index);
142 }
143
144 static void bc_blocka(const t_commrec *cr,t_blocka *block)
145 {
146   block_bc(cr,block->nr);
147   snew_bc(cr,block->index,block->nr+1);
148   nblock_bc(cr,block->nr+1,block->index);
149   block_bc(cr,block->nra);
150   if (block->nra) {
151     snew_bc(cr,block->a,block->nra);
152     nblock_bc(cr,block->nra,block->a);
153   }
154 }
155
156 static void bc_grps(const t_commrec *cr,t_grps grps[])
157 {
158   int i;
159   
160   for(i=0; (i<egcNR); i++) {
161     block_bc(cr,grps[i].nr);
162     snew_bc(cr,grps[i].nm_ind,grps[i].nr);
163     nblock_bc(cr,grps[i].nr,grps[i].nm_ind);
164   }
165 }
166
167 static void bc_atoms(const t_commrec *cr,t_symtab *symtab,t_atoms *atoms)
168 {
169   int dummy;
170
171   block_bc(cr,atoms->nr);
172   snew_bc(cr,atoms->atom,atoms->nr);
173   nblock_bc(cr,atoms->nr,atoms->atom);
174   bc_strings(cr,symtab,atoms->nr,&atoms->atomname);
175   block_bc(cr,atoms->nres);
176   snew_bc(cr,atoms->resinfo,atoms->nres);
177   nblock_bc(cr,atoms->nres,atoms->resinfo);
178   bc_strings_resinfo(cr,symtab,atoms->nres,atoms->resinfo);
179   /* QMMM requires atomtypes to be known on all nodes as well */
180   bc_strings(cr,symtab,atoms->nr,&atoms->atomtype);
181   bc_strings(cr,symtab,atoms->nr,&atoms->atomtypeB);
182 }
183
184 static void bc_groups(const t_commrec *cr,t_symtab *symtab,
185                       int natoms,gmx_groups_t *groups)
186 {
187   int dummy;
188   int g,n;
189
190   bc_grps(cr,groups->grps);
191   block_bc(cr,groups->ngrpname);
192   bc_strings(cr,symtab,groups->ngrpname,&groups->grpname);
193   for(g=0; g<egcNR; g++) {
194     if (MASTER(cr)) {
195       if (groups->grpnr[g]) {
196         n = natoms;
197       } else {
198         n = 0;
199       }
200     }
201     block_bc(cr,n);
202     if (n == 0) {
203       groups->grpnr[g] = NULL;
204     } else {
205       snew_bc(cr,groups->grpnr[g],n);
206       nblock_bc(cr,n,groups->grpnr[g]);
207     }
208   }
209   if (debug) fprintf(debug,"after bc_groups\n");
210 }
211
212 void bcast_state_setup(const t_commrec *cr,t_state *state)
213 {
214   block_bc(cr,state->natoms);
215   block_bc(cr,state->ngtc);
216   block_bc(cr,state->nnhpres);
217   block_bc(cr,state->nhchainlength);
218   block_bc(cr,state->nrng);
219   block_bc(cr,state->nrngi);
220   block_bc(cr,state->flags);
221 }
222
223 void bcast_state(const t_commrec *cr,t_state *state,gmx_bool bAlloc)
224 {
225   int i,nnht,nnhtp;
226
227   bcast_state_setup(cr,state);
228
229   nnht = (state->ngtc)*(state->nhchainlength); 
230   nnhtp = (state->nnhpres)*(state->nhchainlength); 
231
232   if (MASTER(cr)) {
233     bAlloc = FALSE;
234   }
235   if (bAlloc) {
236     state->nalloc = state->natoms;
237   }
238   for(i=0; i<estNR; i++) {
239     if (state->flags & (1<<i)) {
240       switch (i) {
241       case estLAMBDA:  block_bc(cr,state->lambda); break;
242       case estBOX:     block_bc(cr,state->box); break;
243       case estBOX_REL: block_bc(cr,state->box_rel); break;
244       case estBOXV:    block_bc(cr,state->boxv); break;
245       case estPRES_PREV: block_bc(cr,state->pres_prev); break;
246       case estSVIR_PREV: block_bc(cr,state->svir_prev); break;
247       case estFVIR_PREV: block_bc(cr,state->fvir_prev); break;
248       case estNH_XI:   nblock_abc(cr,nnht,state->nosehoover_xi); break;
249       case estNH_VXI:  nblock_abc(cr,nnht,state->nosehoover_vxi); break;
250       case estNHPRES_XI:   nblock_abc(cr,nnhtp,state->nhpres_xi); break;
251       case estNHPRES_VXI:  nblock_abc(cr,nnhtp,state->nhpres_vxi); break;
252       case estTC_INT:  nblock_abc(cr,state->ngtc,state->therm_integral); break;
253       case estVETA:    block_bc(cr,state->veta); break;
254       case estVOL0:    block_bc(cr,state->vol0); break;
255       case estX:       nblock_abc(cr,state->natoms,state->x); break;
256       case estV:       nblock_abc(cr,state->natoms,state->v); break;
257       case estSDX:     nblock_abc(cr,state->natoms,state->sd_X); break;
258       case estCGP:     nblock_abc(cr,state->natoms,state->cg_p); break;
259           case estLD_RNG:  if(state->nrngi == 1) nblock_abc(cr,state->nrng,state->ld_rng); break;
260           case estLD_RNGI: if(state->nrngi == 1) nblock_abc(cr,state->nrngi,state->ld_rngi); break;
261       case estDISRE_INITF: block_bc(cr,state->hist.disre_initf); break;
262       case estDISRE_RM3TAV:
263           block_bc(cr,state->hist.ndisrepairs);
264           nblock_abc(cr,state->hist.ndisrepairs,state->hist.disre_rm3tav);
265           break;
266       case estORIRE_INITF: block_bc(cr,state->hist.orire_initf); break;
267       case estORIRE_DTAV:
268           block_bc(cr,state->hist.norire_Dtav);
269           nblock_abc(cr,state->hist.norire_Dtav,state->hist.orire_Dtav);
270           break;
271       default:
272           gmx_fatal(FARGS,
273                     "Communication is not implemented for %s in bcast_state",
274                     est_names[i]);
275       }
276     }
277   }
278 }
279
280 static void bc_ilists(const t_commrec *cr,t_ilist *ilist)
281 {
282   int ftype;
283
284   /* Here we only communicate the non-zero length ilists */
285   if (MASTER(cr)) {
286     for(ftype=0; ftype<F_NRE; ftype++) {
287       if (ilist[ftype].nr > 0) {
288         block_bc(cr,ftype);
289         block_bc(cr,ilist[ftype].nr);
290         nblock_bc(cr,ilist[ftype].nr,ilist[ftype].iatoms);
291       }
292     }
293     ftype = -1;
294     block_bc(cr,ftype);
295   } else {
296     for(ftype=0; ftype<F_NRE; ftype++) {
297       ilist[ftype].nr = 0;
298     }
299     do {
300       block_bc(cr,ftype);
301       if (ftype >= 0) {
302         block_bc(cr,ilist[ftype].nr);
303         snew_bc(cr,ilist[ftype].iatoms,ilist[ftype].nr);
304         nblock_bc(cr,ilist[ftype].nr,ilist[ftype].iatoms);
305       }
306     } while (ftype >= 0);
307   }
308
309   if (debug) fprintf(debug,"after bc_ilists\n");
310 }
311
312 static void bc_idef(const t_commrec *cr,t_idef *idef)
313 {
314   block_bc(cr,idef->ntypes);
315   block_bc(cr,idef->atnr);
316   snew_bc(cr,idef->functype,idef->ntypes);
317   snew_bc(cr,idef->iparams,idef->ntypes);
318   nblock_bc(cr,idef->ntypes,idef->functype);
319   nblock_bc(cr,idef->ntypes,idef->iparams);
320   block_bc(cr,idef->fudgeQQ);
321   bc_ilists(cr,idef->il);
322   block_bc(cr,idef->ilsort);
323 }
324
325 static void bc_cmap(const t_commrec *cr, gmx_cmap_t *cmap_grid)
326 {
327         int i,j,nelem,ngrid;
328         
329         block_bc(cr,cmap_grid->ngrid);
330         block_bc(cr,cmap_grid->grid_spacing);
331         
332         ngrid = cmap_grid->ngrid;
333         nelem = cmap_grid->grid_spacing * cmap_grid->grid_spacing;
334         
335         if(ngrid>0)
336         {
337                 snew_bc(cr,cmap_grid->cmapdata,ngrid);
338                 
339                 for(i=0;i<ngrid;i++)
340                 {
341                         snew_bc(cr,cmap_grid->cmapdata[i].cmap,4*nelem);
342                         nblock_bc(cr,4*nelem,cmap_grid->cmapdata[i].cmap);
343                 }
344         }
345 }
346
347 static void bc_ffparams(const t_commrec *cr,gmx_ffparams_t *ffp)
348 {
349   int i;
350   
351   block_bc(cr,ffp->ntypes);
352   block_bc(cr,ffp->atnr);
353   snew_bc(cr,ffp->functype,ffp->ntypes);
354   snew_bc(cr,ffp->iparams,ffp->ntypes);
355   nblock_bc(cr,ffp->ntypes,ffp->functype);
356   nblock_bc(cr,ffp->ntypes,ffp->iparams);
357   block_bc(cr,ffp->reppow);
358   block_bc(cr,ffp->fudgeQQ);
359   bc_cmap(cr,&ffp->cmap_grid);
360 }
361
362 static void bc_grpopts(const t_commrec *cr,t_grpopts *g)
363 {
364     int i,n;
365     
366     block_bc(cr,g->ngtc);
367     block_bc(cr,g->ngacc);
368     block_bc(cr,g->ngfrz);
369     block_bc(cr,g->ngener);
370     snew_bc(cr,g->nrdf,g->ngtc);
371     snew_bc(cr,g->tau_t,g->ngtc);
372     snew_bc(cr,g->ref_t,g->ngtc);
373     snew_bc(cr,g->acc,g->ngacc);
374     snew_bc(cr,g->nFreeze,g->ngfrz);
375     snew_bc(cr,g->egp_flags,g->ngener*g->ngener);
376     
377     nblock_bc(cr,g->ngtc,g->nrdf);
378     nblock_bc(cr,g->ngtc,g->tau_t);
379     nblock_bc(cr,g->ngtc,g->ref_t);
380     nblock_bc(cr,g->ngacc,g->acc);
381     nblock_bc(cr,g->ngfrz,g->nFreeze);
382     nblock_bc(cr,g->ngener*g->ngener,g->egp_flags);
383     snew_bc(cr,g->annealing,g->ngtc);
384     snew_bc(cr,g->anneal_npoints,g->ngtc);
385     snew_bc(cr,g->anneal_time,g->ngtc);
386     snew_bc(cr,g->anneal_temp,g->ngtc);
387     nblock_bc(cr,g->ngtc,g->annealing);
388     nblock_bc(cr,g->ngtc,g->anneal_npoints);
389     for(i=0;(i<g->ngtc); i++) {
390         n = g->anneal_npoints[i];
391         if (n > 0) {
392           snew_bc(cr,g->anneal_time[i],n);
393           snew_bc(cr,g->anneal_temp[i],n);
394           nblock_bc(cr,n,g->anneal_time[i]);
395           nblock_bc(cr,n,g->anneal_temp[i]);
396         }
397     }
398     
399     /* QMMM stuff, see inputrec */
400     block_bc(cr,g->ngQM);
401     snew_bc(cr,g->QMmethod,g->ngQM);
402     snew_bc(cr,g->QMbasis,g->ngQM);
403     snew_bc(cr,g->QMcharge,g->ngQM);
404     snew_bc(cr,g->QMmult,g->ngQM);
405     snew_bc(cr,g->bSH,g->ngQM);
406     snew_bc(cr,g->CASorbitals,g->ngQM);
407     snew_bc(cr,g->CASelectrons,g->ngQM);
408     snew_bc(cr,g->SAon,g->ngQM);
409     snew_bc(cr,g->SAoff,g->ngQM);
410     snew_bc(cr,g->SAsteps,g->ngQM);
411     
412     if (g->ngQM)
413     {
414         nblock_bc(cr,g->ngQM,g->QMmethod);
415         nblock_bc(cr,g->ngQM,g->QMbasis);
416         nblock_bc(cr,g->ngQM,g->QMcharge);
417         nblock_bc(cr,g->ngQM,g->QMmult);
418         nblock_bc(cr,g->ngQM,g->bSH);
419         nblock_bc(cr,g->ngQM,g->CASorbitals);
420         nblock_bc(cr,g->ngQM,g->CASelectrons);
421         nblock_bc(cr,g->ngQM,g->SAon);
422         nblock_bc(cr,g->ngQM,g->SAoff);
423         nblock_bc(cr,g->ngQM,g->SAsteps);
424         /* end of QMMM stuff */
425     }
426 }
427
428 static void bc_cosines(const t_commrec *cr,t_cosines *cs)
429 {
430   block_bc(cr,cs->n);
431   snew_bc(cr,cs->a,cs->n);
432   snew_bc(cr,cs->phi,cs->n);
433   if (cs->n > 0) {
434     nblock_bc(cr,cs->n,cs->a);
435     nblock_bc(cr,cs->n,cs->phi);
436   }
437 }
438
439 static void bc_pullgrp(const t_commrec *cr,t_pullgrp *pgrp)
440 {
441   block_bc(cr,*pgrp);
442   if (pgrp->nat > 0) {
443     snew_bc(cr,pgrp->ind,pgrp->nat);
444     nblock_bc(cr,pgrp->nat,pgrp->ind);
445   }
446   if (pgrp->nweight > 0) {
447     snew_bc(cr,pgrp->weight,pgrp->nweight);
448     nblock_bc(cr,pgrp->nweight,pgrp->weight);
449   }
450 }
451
452 static void bc_pull(const t_commrec *cr,t_pull *pull)
453 {
454   int g;
455
456   block_bc(cr,*pull);
457   snew_bc(cr,pull->grp,pull->ngrp+1);
458   for(g=0; g<pull->ngrp+1; g++)
459   {
460       bc_pullgrp(cr,&pull->grp[g]);
461   }
462 }
463
464 static void bc_rotgrp(const t_commrec *cr,t_rotgrp *rotg)
465 {
466   block_bc(cr,*rotg);
467   if (rotg->nat > 0) {
468     snew_bc(cr,rotg->ind,rotg->nat);
469     nblock_bc(cr,rotg->nat,rotg->ind);
470     snew_bc(cr,rotg->x_ref,rotg->nat);
471     nblock_bc(cr,rotg->nat,rotg->x_ref);
472   }
473 }
474
475 static void bc_rot(const t_commrec *cr,t_rot *rot)
476 {
477   int g;
478
479   block_bc(cr,*rot);
480   snew_bc(cr,rot->grp,rot->ngrp);
481   for(g=0; g<rot->ngrp; g++)
482     bc_rotgrp(cr,&rot->grp[g]);
483 }
484
485 static void bc_inputrec(const t_commrec *cr,t_inputrec *inputrec)
486 {
487   gmx_bool bAlloc=TRUE;
488   int i;
489   
490   block_bc(cr,*inputrec);
491   snew_bc(cr,inputrec->flambda,inputrec->n_flambda);
492   nblock_bc(cr,inputrec->n_flambda,inputrec->flambda);
493   bc_grpopts(cr,&(inputrec->opts));
494   if (inputrec->ePull != epullNO) {
495     snew_bc(cr,inputrec->pull,1);
496     bc_pull(cr,inputrec->pull);
497   }
498   if (inputrec->bRot) {
499     snew_bc(cr,inputrec->rot,1);
500     bc_rot(cr,inputrec->rot);
501   }
502   for(i=0; (i<DIM); i++) {
503     bc_cosines(cr,&(inputrec->ex[i]));
504     bc_cosines(cr,&(inputrec->et[i]));
505   }
506 }
507
508 static void bc_moltype(const t_commrec *cr,t_symtab *symtab,
509                        gmx_moltype_t *moltype)
510 {
511   bc_string(cr,symtab,&moltype->name);
512   bc_atoms(cr,symtab,&moltype->atoms);
513   if (debug) fprintf(debug,"after bc_atoms\n");
514
515   bc_ilists(cr,moltype->ilist);
516   bc_block(cr,&moltype->cgs);
517   bc_blocka(cr,&moltype->excls);
518 }
519
520 static void bc_molblock(const t_commrec *cr,gmx_molblock_t *molb)
521 {
522   gmx_bool bAlloc=TRUE;
523   
524   block_bc(cr,molb->type);
525   block_bc(cr,molb->nmol);
526   block_bc(cr,molb->natoms_mol);
527   block_bc(cr,molb->nposres_xA);
528   if (molb->nposres_xA > 0) {
529     snew_bc(cr,molb->posres_xA,molb->nposres_xA);
530     nblock_bc(cr,molb->nposres_xA*DIM,molb->posres_xA[0]);
531   }
532   block_bc(cr,molb->nposres_xB);
533   if (molb->nposres_xB > 0) {
534     snew_bc(cr,molb->posres_xB,molb->nposres_xB);
535     nblock_bc(cr,molb->nposres_xB*DIM,molb->posres_xB[0]);
536   }
537   if (debug) fprintf(debug,"after bc_molblock\n");
538 }
539
540 static void bc_atomtypes(const t_commrec *cr, t_atomtypes *atomtypes)
541 {
542   int nr;
543
544   block_bc(cr,atomtypes->nr);
545
546   nr = atomtypes->nr;
547
548   snew_bc(cr,atomtypes->radius,nr);
549   snew_bc(cr,atomtypes->vol,nr);
550   snew_bc(cr,atomtypes->surftens,nr);
551   snew_bc(cr,atomtypes->gb_radius,nr);
552   snew_bc(cr,atomtypes->S_hct,nr);
553
554   nblock_bc(cr,nr,atomtypes->radius);
555   nblock_bc(cr,nr,atomtypes->vol);
556   nblock_bc(cr,nr,atomtypes->surftens);
557   nblock_bc(cr,nr,atomtypes->gb_radius);
558   nblock_bc(cr,nr,atomtypes->S_hct);
559 }
560
561
562 void bcast_ir_mtop(const t_commrec *cr,t_inputrec *inputrec,gmx_mtop_t *mtop)
563 {
564   int i; 
565   if (debug) fprintf(debug,"in bc_data\n");
566   bc_inputrec(cr,inputrec);
567   if (debug) fprintf(debug,"after bc_inputrec\n");
568   bc_symtab(cr,&mtop->symtab);
569   if (debug) fprintf(debug,"after bc_symtab\n");
570   bc_string(cr,&mtop->symtab,&mtop->name);
571   if (debug) fprintf(debug,"after bc_name\n");
572
573   bc_ffparams(cr,&mtop->ffparams);
574
575   block_bc(cr,mtop->nmoltype);
576   snew_bc(cr,mtop->moltype,mtop->nmoltype);
577   for(i=0; i<mtop->nmoltype; i++) {
578     bc_moltype(cr,&mtop->symtab,&mtop->moltype[i]);
579   }
580
581   block_bc(cr,mtop->nmolblock);
582   snew_bc(cr,mtop->molblock,mtop->nmolblock);
583   for(i=0; i<mtop->nmolblock; i++) {
584     bc_molblock(cr,&mtop->molblock[i]);
585   }
586
587   block_bc(cr,mtop->natoms);
588
589   bc_atomtypes(cr,&mtop->atomtypes);
590
591   bc_block(cr,&mtop->mols);
592   bc_groups(cr,&mtop->symtab,mtop->natoms,&mtop->groups);
593 }