Framework for printing help for selections.
[gromacs/adressmacs.git] / src / gmxlib / mvdata.c
blob702a1487008a3b2eb069a72471e1aa670e7273ac
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 /* This file is completely threadsafe - keep it that way! */
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
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"
52 #define block_bc(cr, d) gmx_bcast( sizeof(d), &(d),(cr))
53 #define nblock_bc(cr,nr,d) gmx_bcast((nr)*sizeof((d)[0]), (d),(cr))
54 #define snew_bc(cr,d,nr) { if (!MASTER(cr)) snew((d),(nr)); }
55 /* Dirty macro with bAlloc not as an argument */
56 #define nblock_abc(cr,nr,d) { if (bAlloc) snew((d),(nr)); nblock_bc(cr,(nr),(d)); }
58 static void bc_string(const t_commrec *cr,t_symtab *symtab,char ***s)
60 int handle;
62 if (MASTER(cr)) {
63 handle = lookup_symtab(symtab,*s);
65 block_bc(cr,handle);
66 if (!MASTER(cr)) {
67 *s = get_symtab_handle(symtab,handle);
71 static void bc_strings(const t_commrec *cr,t_symtab *symtab,int nr,char ****nm)
73 int i;
74 int *handle;
75 char ***NM;
77 snew(handle,nr);
78 if (MASTER(cr)) {
79 NM = *nm;
80 for(i=0; (i<nr); i++)
81 handle[i] = lookup_symtab(symtab,NM[i]);
83 nblock_bc(cr,nr,handle);
85 if (!MASTER(cr)) {
86 snew_bc(cr,*nm,nr);
87 NM = *nm;
88 for (i=0; (i<nr); i++)
89 (*nm)[i] = get_symtab_handle(symtab,handle[i]);
91 sfree(handle);
94 static void bc_strings_resinfo(const t_commrec *cr,t_symtab *symtab,
95 int nr,t_resinfo *resinfo)
97 int i;
98 int *handle;
100 snew(handle,nr);
101 if (MASTER(cr)) {
102 for(i=0; (i<nr); i++)
103 handle[i] = lookup_symtab(symtab,resinfo[i].name);
105 nblock_bc(cr,nr,handle);
107 if (!MASTER(cr)) {
108 for (i=0; (i<nr); i++)
109 resinfo[i].name = get_symtab_handle(symtab,handle[i]);
111 sfree(handle);
114 static void bc_symtab(const t_commrec *cr,t_symtab *symtab)
116 int i,nr,len;
117 t_symbuf *symbuf;
119 block_bc(cr,symtab->nr);
120 nr = symtab->nr;
121 snew_bc(cr,symtab->symbuf,1);
122 symbuf = symtab->symbuf;
123 symbuf->bufsize = nr;
124 snew_bc(cr,symbuf->buf,nr);
125 for (i=0; i<nr; i++) {
126 if (MASTER(cr))
127 len = strlen(symbuf->buf[i]) + 1;
128 block_bc(cr,len);
129 snew_bc(cr,symbuf->buf[i],len);
130 nblock_bc(cr,len,symbuf->buf[i]);
134 static void bc_block(const t_commrec *cr,t_block *block)
136 block_bc(cr,block->nr);
137 snew_bc(cr,block->index,block->nr+1);
138 nblock_bc(cr,block->nr+1,block->index);
141 static void bc_blocka(const t_commrec *cr,t_blocka *block)
143 block_bc(cr,block->nr);
144 snew_bc(cr,block->index,block->nr+1);
145 nblock_bc(cr,block->nr+1,block->index);
146 block_bc(cr,block->nra);
147 if (block->nra) {
148 snew_bc(cr,block->a,block->nra);
149 nblock_bc(cr,block->nra,block->a);
153 static void bc_grps(const t_commrec *cr,t_grps grps[])
155 int i;
157 for(i=0; (i<egcNR); i++) {
158 block_bc(cr,grps[i].nr);
159 snew_bc(cr,grps[i].nm_ind,grps[i].nr);
160 nblock_bc(cr,grps[i].nr,grps[i].nm_ind);
164 static void bc_atoms(const t_commrec *cr,t_symtab *symtab,t_atoms *atoms)
166 int dummy;
168 block_bc(cr,atoms->nr);
169 snew_bc(cr,atoms->atom,atoms->nr);
170 nblock_bc(cr,atoms->nr,atoms->atom);
171 bc_strings(cr,symtab,atoms->nr,&atoms->atomname);
172 block_bc(cr,atoms->nres);
173 snew_bc(cr,atoms->resinfo,atoms->nres);
174 nblock_bc(cr,atoms->nres,atoms->resinfo);
175 bc_strings_resinfo(cr,symtab,atoms->nres,atoms->resinfo);
176 /* QMMM requires atomtypes to be known on all nodes as well */
177 bc_strings(cr,symtab,atoms->nr,&atoms->atomtype);
178 bc_strings(cr,symtab,atoms->nr,&atoms->atomtypeB);
181 static void bc_groups(const t_commrec *cr,t_symtab *symtab,
182 int natoms,gmx_groups_t *groups)
184 int dummy;
185 int g,n;
187 bc_grps(cr,groups->grps);
188 block_bc(cr,groups->ngrpname);
189 bc_strings(cr,symtab,groups->ngrpname,&groups->grpname);
190 for(g=0; g<egcNR; g++) {
191 if (MASTER(cr)) {
192 if (groups->grpnr[g]) {
193 n = natoms;
194 } else {
195 n = 0;
198 block_bc(cr,n);
199 if (n == 0) {
200 groups->grpnr[g] = NULL;
201 } else {
202 snew_bc(cr,groups->grpnr[g],n);
203 nblock_bc(cr,n,groups->grpnr[g]);
206 if (debug) fprintf(debug,"after bc_groups\n");
209 void bcast_state_setup(const t_commrec *cr,t_state *state)
211 block_bc(cr,state->natoms);
212 block_bc(cr,state->ngtc);
213 block_bc(cr,state->nrng);
214 block_bc(cr,state->nrngi);
215 block_bc(cr,state->flags);
218 void bcast_state(const t_commrec *cr,t_state *state,bool bAlloc)
220 int i;
222 bcast_state_setup(cr,state);
224 if (MASTER(cr)) {
225 bAlloc = FALSE;
227 if (bAlloc) {
228 state->nalloc = state->natoms;
230 for(i=0; i<estNR; i++) {
231 if (state->flags & (1<<i)) {
232 switch (i) {
233 case estLAMBDA: block_bc(cr,state->lambda); break;
234 case estBOX: block_bc(cr,state->box); break;
235 case estBOX_REL: block_bc(cr,state->box_rel); break;
236 case estBOXV: block_bc(cr,state->boxv); break;
237 case estPRES_PREV: block_bc(cr,state->pres_prev); break;
238 case estNH_XI: nblock_abc(cr,state->ngtc,state->nosehoover_xi); break;
239 case estTC_INT: nblock_abc(cr,state->ngtc,state->therm_integral); break;
240 case estX: nblock_abc(cr,state->natoms,state->x); break;
241 case estV: nblock_abc(cr,state->natoms,state->v); break;
242 case estSDX: nblock_abc(cr,state->natoms,state->sd_X); break;
243 case estCGP: nblock_abc(cr,state->natoms,state->cg_p); break;
244 case estLD_RNG: if(state->nrngi == 1) nblock_abc(cr,state->nrng,state->ld_rng); break;
245 case estLD_RNGI: if(state->nrngi == 1) nblock_abc(cr,state->nrngi,state->ld_rngi); break;
246 case estDISRE_INITF: block_bc(cr,state->hist.disre_initf); break;
247 case estDISRE_RM3TAV:
248 block_bc(cr,state->hist.ndisrepairs);
249 nblock_abc(cr,state->hist.ndisrepairs,state->hist.disre_rm3tav);
250 break;
251 case estORIRE_INITF: block_bc(cr,state->hist.orire_initf); break;
252 case estORIRE_DTAV:
253 block_bc(cr,state->hist.norire_Dtav);
254 nblock_abc(cr,state->hist.norire_Dtav,state->hist.orire_Dtav);
255 break;
256 default:
257 gmx_fatal(FARGS,
258 "Communication is not implemented for %s in bcast_state",
259 est_names[i]);
265 static void bc_ilists(const t_commrec *cr,t_ilist *ilist)
267 int ftype;
269 /* Here we only communicate the non-zero length ilists */
270 if (MASTER(cr)) {
271 for(ftype=0; ftype<F_NRE; ftype++) {
272 if (ilist[ftype].nr > 0) {
273 block_bc(cr,ftype);
274 block_bc(cr,ilist[ftype].nr);
275 nblock_bc(cr,ilist[ftype].nr,ilist[ftype].iatoms);
278 ftype = -1;
279 block_bc(cr,ftype);
280 } else {
281 for(ftype=0; ftype<F_NRE; ftype++) {
282 ilist[ftype].nr = 0;
284 do {
285 block_bc(cr,ftype);
286 if (ftype >= 0) {
287 block_bc(cr,ilist[ftype].nr);
288 snew_bc(cr,ilist[ftype].iatoms,ilist[ftype].nr);
289 nblock_bc(cr,ilist[ftype].nr,ilist[ftype].iatoms);
291 } while (ftype >= 0);
294 if (debug) fprintf(debug,"after bc_ilists\n");
297 static void bc_idef(const t_commrec *cr,t_idef *idef)
299 block_bc(cr,idef->ntypes);
300 block_bc(cr,idef->atnr);
301 snew_bc(cr,idef->functype,idef->ntypes);
302 snew_bc(cr,idef->iparams,idef->ntypes);
303 nblock_bc(cr,idef->ntypes,idef->functype);
304 nblock_bc(cr,idef->ntypes,idef->iparams);
305 block_bc(cr,idef->fudgeQQ);
306 bc_ilists(cr,idef->il);
307 block_bc(cr,idef->ilsort);
310 static void bc_ffparams(const t_commrec *cr,gmx_ffparams_t *ffp)
312 int i;
314 block_bc(cr,ffp->ntypes);
315 block_bc(cr,ffp->atnr);
316 snew_bc(cr,ffp->functype,ffp->ntypes);
317 snew_bc(cr,ffp->iparams,ffp->ntypes);
318 nblock_bc(cr,ffp->ntypes,ffp->functype);
319 nblock_bc(cr,ffp->ntypes,ffp->iparams);
320 block_bc(cr,ffp->reppow);
321 block_bc(cr,ffp->fudgeQQ);
324 static void bc_cmap(const t_commrec *cr, gmx_cmap_t *cmap_grid)
326 int i,j,nelem,ngrid;
328 block_bc(cr,cmap_grid->ngrid);
329 block_bc(cr,cmap_grid->grid_spacing);
331 ngrid = cmap_grid->ngrid;
332 nelem = cmap_grid->grid_spacing * cmap_grid->grid_spacing;
334 if(ngrid>0)
336 snew_bc(cr,cmap_grid->cmapdata,ngrid);
338 for(i=0;i<ngrid;i++)
340 snew_bc(cr,cmap_grid->cmapdata[i].cmap,4*nelem);
341 nblock_bc(cr,4*nelem,cmap_grid->cmapdata[i].cmap);
347 static void bc_grpopts(const t_commrec *cr,t_grpopts *g)
349 int i,n;
351 block_bc(cr,g->ngtc);
352 block_bc(cr,g->ngacc);
353 block_bc(cr,g->ngfrz);
354 block_bc(cr,g->ngener);
355 snew_bc(cr,g->nrdf,g->ngtc);
356 snew_bc(cr,g->tau_t,g->ngtc);
357 snew_bc(cr,g->ref_t,g->ngtc);
358 snew_bc(cr,g->acc,g->ngacc);
359 snew_bc(cr,g->nFreeze,g->ngfrz);
360 snew_bc(cr,g->egp_flags,g->ngener*g->ngener);
362 nblock_bc(cr,g->ngtc,g->nrdf);
363 nblock_bc(cr,g->ngtc,g->tau_t);
364 nblock_bc(cr,g->ngtc,g->ref_t);
365 nblock_bc(cr,g->ngacc,g->acc);
366 nblock_bc(cr,g->ngfrz,g->nFreeze);
367 nblock_bc(cr,g->ngener*g->ngener,g->egp_flags);
368 snew_bc(cr,g->annealing,g->ngtc);
369 snew_bc(cr,g->anneal_npoints,g->ngtc);
370 snew_bc(cr,g->anneal_time,g->ngtc);
371 snew_bc(cr,g->anneal_temp,g->ngtc);
372 nblock_bc(cr,g->ngtc,g->annealing);
373 nblock_bc(cr,g->ngtc,g->anneal_npoints);
374 for(i=0;(i<g->ngtc); i++) {
375 n = g->anneal_npoints[i];
376 if (n > 0) {
377 snew_bc(cr,g->anneal_time[i],n);
378 snew_bc(cr,g->anneal_temp[i],n);
379 nblock_bc(cr,n,g->anneal_time[i]);
380 nblock_bc(cr,n,g->anneal_temp[i]);
384 /* QMMM stuff, see inputrec */
385 block_bc(cr,g->ngQM);
386 snew_bc(cr,g->QMmethod,g->ngQM);
387 snew_bc(cr,g->QMbasis,g->ngQM);
388 snew_bc(cr,g->QMcharge,g->ngQM);
389 snew_bc(cr,g->QMmult,g->ngQM);
390 snew_bc(cr,g->bSH,g->ngQM);
391 snew_bc(cr,g->CASorbitals,g->ngQM);
392 snew_bc(cr,g->CASelectrons,g->ngQM);
393 snew_bc(cr,g->SAon,g->ngQM);
394 snew_bc(cr,g->SAoff,g->ngQM);
395 snew_bc(cr,g->SAsteps,g->ngQM);
397 if (g->ngQM)
399 nblock_bc(cr,g->ngQM,g->QMmethod);
400 nblock_bc(cr,g->ngQM,g->QMbasis);
401 nblock_bc(cr,g->ngQM,g->QMcharge);
402 nblock_bc(cr,g->ngQM,g->QMmult);
403 nblock_bc(cr,g->ngQM,g->bSH);
404 nblock_bc(cr,g->ngQM,g->CASorbitals);
405 nblock_bc(cr,g->ngQM,g->CASelectrons);
406 nblock_bc(cr,g->ngQM,g->SAon);
407 nblock_bc(cr,g->ngQM,g->SAoff);
408 nblock_bc(cr,g->ngQM,g->SAsteps);
409 /* end of QMMM stuff */
413 static void bc_cosines(const t_commrec *cr,t_cosines *cs)
415 block_bc(cr,cs->n);
416 snew_bc(cr,cs->a,cs->n);
417 snew_bc(cr,cs->phi,cs->n);
418 if (cs->n > 0) {
419 nblock_bc(cr,cs->n,cs->a);
420 nblock_bc(cr,cs->n,cs->phi);
424 static void bc_pullgrp(const t_commrec *cr,t_pullgrp *pgrp)
426 block_bc(cr,*pgrp);
427 if (pgrp->nat > 0) {
428 snew_bc(cr,pgrp->ind,pgrp->nat);
429 nblock_bc(cr,pgrp->nat,pgrp->ind);
431 if (pgrp->nweight > 0) {
432 snew_bc(cr,pgrp->weight,pgrp->nweight);
433 nblock_bc(cr,pgrp->nweight,pgrp->weight);
437 static void bc_pull(const t_commrec *cr,t_pull *pull)
439 int g;
441 block_bc(cr,*pull);
442 snew_bc(cr,pull->grp,pull->ngrp+1);
443 for(g=0; g<pull->ngrp+1; g++)
444 bc_pullgrp(cr,&pull->grp[g]);
447 static void bc_inputrec(const t_commrec *cr,t_inputrec *inputrec)
449 bool bAlloc=TRUE;
450 int i;
452 block_bc(cr,*inputrec);
453 snew_bc(cr,inputrec->flambda,inputrec->n_flambda);
454 nblock_bc(cr,inputrec->n_flambda,inputrec->flambda);
455 bc_grpopts(cr,&(inputrec->opts));
456 if (inputrec->ePull != epullNO) {
457 snew_bc(cr,inputrec->pull,1);
458 bc_pull(cr,inputrec->pull);
460 for(i=0; (i<DIM); i++) {
461 bc_cosines(cr,&(inputrec->ex[i]));
462 bc_cosines(cr,&(inputrec->et[i]));
466 static void bc_moltype(const t_commrec *cr,t_symtab *symtab,
467 gmx_moltype_t *moltype)
469 bc_string(cr,symtab,&moltype->name);
470 bc_atoms(cr,symtab,&moltype->atoms);
471 if (debug) fprintf(debug,"after bc_atoms\n");
473 bc_ilists(cr,moltype->ilist);
474 bc_block(cr,&moltype->cgs);
475 bc_blocka(cr,&moltype->excls);
478 static void bc_molblock(const t_commrec *cr,gmx_molblock_t *molb)
480 bool bAlloc=TRUE;
482 block_bc(cr,molb->type);
483 block_bc(cr,molb->nmol);
484 block_bc(cr,molb->natoms_mol);
485 block_bc(cr,molb->nposres_xA);
486 if (molb->nposres_xA > 0) {
487 snew_bc(cr,molb->posres_xA,molb->nposres_xA);
488 nblock_bc(cr,molb->nposres_xA*DIM,molb->posres_xA[0]);
490 block_bc(cr,molb->nposres_xB);
491 if (molb->nposres_xB > 0) {
492 snew_bc(cr,molb->posres_xB,molb->nposres_xB);
493 nblock_bc(cr,molb->nposres_xB*DIM,molb->posres_xB[0]);
495 if (debug) fprintf(debug,"after bc_molblock\n");
498 static void bc_atomtypes(const t_commrec *cr, t_atomtypes *atomtypes)
500 int nr;
502 block_bc(cr,atomtypes->nr);
504 nr = atomtypes->nr;
506 snew_bc(cr,atomtypes->radius,nr);
507 snew_bc(cr,atomtypes->vol,nr);
508 snew_bc(cr,atomtypes->surftens,nr);
509 snew_bc(cr,atomtypes->gb_radius,nr);
510 snew_bc(cr,atomtypes->S_hct,nr);
512 nblock_bc(cr,nr,atomtypes->radius);
513 nblock_bc(cr,nr,atomtypes->vol);
514 nblock_bc(cr,nr,atomtypes->surftens);
515 nblock_bc(cr,nr,atomtypes->gb_radius);
516 nblock_bc(cr,nr,atomtypes->S_hct);
520 void bcast_ir_mtop(const t_commrec *cr,t_inputrec *inputrec,gmx_mtop_t *mtop)
522 int i;
523 if (debug) fprintf(debug,"in bc_data\n");
524 bc_inputrec(cr,inputrec);
525 if (debug) fprintf(debug,"after bc_inputrec\n");
526 bc_symtab(cr,&mtop->symtab);
527 if (debug) fprintf(debug,"after bc_symtab\n");
528 bc_string(cr,&mtop->symtab,&mtop->name);
529 if (debug) fprintf(debug,"after bc_name\n");
531 bc_ffparams(cr,&mtop->ffparams);
533 block_bc(cr,mtop->nmoltype);
534 snew_bc(cr,mtop->moltype,mtop->nmoltype);
535 for(i=0; i<mtop->nmoltype; i++) {
536 bc_moltype(cr,&mtop->symtab,&mtop->moltype[i]);
539 block_bc(cr,mtop->nmolblock);
540 snew_bc(cr,mtop->molblock,mtop->nmolblock);
541 for(i=0; i<mtop->nmolblock; i++) {
542 bc_molblock(cr,&mtop->molblock[i]);
545 block_bc(cr,mtop->natoms);
547 bc_atomtypes(cr,&mtop->atomtypes);
549 bc_block(cr,&mtop->mols);
550 bc_groups(cr,&mtop->symtab,mtop->natoms,&mtop->groups);
551 bc_cmap(cr,&mtop->cmap_grid);