Fixed typo in CMakeLists for mdrun-gpu
[gromacs/qmmm-gamess-us.git] / src / mdlib / constr.c
blobdfa6b001463eed047d1bbcef77cf903155976272
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
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.
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.
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.
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.
31 * For more info, check our website at http://www.gromacs.org
33 * And Hey:
34 * GROwing Monsters And Cloning Shrimps
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
40 #include "confio.h"
41 #include "constr.h"
42 #include "copyrite.h"
43 #include "invblock.h"
44 #include "main.h"
45 #include "mdrun.h"
46 #include "nrnb.h"
47 #include "smalloc.h"
48 #include "vec.h"
49 #include "physics.h"
50 #include "names.h"
51 #include "txtdump.h"
52 #include "domdec.h"
53 #include "pdbio.h"
54 #include "partdec.h"
55 #include "splitter.h"
56 #include "mtop_util.h"
57 #include "gmxfio.h"
59 typedef struct gmx_constr {
60 int ncon_tot; /* The total number of constraints */
61 int nflexcon; /* The number of flexible constraints */
62 int n_at2con_mt; /* The size of at2con = #moltypes */
63 t_blocka *at2con_mt; /* A list of atoms to constraints */
64 gmx_lincsdata_t lincsd; /* LINCS data */
65 gmx_shakedata_t shaked; /* SHAKE data */
66 gmx_settledata_t settled; /* SETTLE data */
67 int nblocks; /* The number of SHAKE blocks */
68 int *sblock; /* The SHAKE blocks */
69 int sblock_nalloc;/* The allocation size of sblock */
70 real *lagr; /* Lagrange multipliers for SHAKE */
71 int lagr_nalloc; /* The allocation size of lagr */
72 int maxwarn; /* The maximum number of warnings */
73 int warncount_lincs;
74 int warncount_settle;
75 gmx_edsam_t ed; /* The essential dynamics data */
77 gmx_mtop_t *warn_mtop; /* Only used for printing warnings */
78 } t_gmx_constr;
80 typedef struct {
81 atom_id iatom[3];
82 atom_id blocknr;
83 } t_sortblock;
85 static t_vetavars *init_vetavars(real veta,real vetanew, t_inputrec *ir, gmx_ekindata_t *ekind, bool bPscal)
87 t_vetavars *vars;
88 double g;
89 int i;
91 snew(vars,1);
92 snew(vars->vscale_nhc,ir->opts.ngtc);
93 /* first, set the alpha integrator variable */
94 if ((ir->opts.nrdf[0] > 0) && bPscal)
96 vars->alpha = 1.0 + DIM/((double)ir->opts.nrdf[0]);
97 } else {
98 vars->alpha = 1.0;
100 g = 0.5*veta*ir->delta_t;
101 vars->rscale = exp(g)*series_sinhx(g);
102 g = -0.25*vars->alpha*veta*ir->delta_t;
103 vars->vscale = exp(g)*series_sinhx(g);
104 vars->rvscale = vars->vscale*vars->rscale;
105 vars->veta = vetanew;
106 if ((ekind==NULL) || (!bPscal))
108 for (i=0;i<ir->opts.ngtc;i++)
110 vars->vscale_nhc[i] = 1;
112 } else {
113 for (i=0;i<ir->opts.ngtc;i++)
115 vars->vscale_nhc[i] = ekind->tcstat[i].vscale_nhc;
118 return vars;
121 static void free_vetavars(t_vetavars *vars)
123 sfree(vars->vscale_nhc);
124 sfree(vars);
127 static int pcomp(const void *p1, const void *p2)
129 int db;
130 atom_id min1,min2,max1,max2;
131 t_sortblock *a1=(t_sortblock *)p1;
132 t_sortblock *a2=(t_sortblock *)p2;
134 db=a1->blocknr-a2->blocknr;
136 if (db != 0)
137 return db;
139 min1=min(a1->iatom[1],a1->iatom[2]);
140 max1=max(a1->iatom[1],a1->iatom[2]);
141 min2=min(a2->iatom[1],a2->iatom[2]);
142 max2=max(a2->iatom[1],a2->iatom[2]);
144 if (min1 == min2)
145 return max1-max2;
146 else
147 return min1-min2;
150 static int icomp(const void *p1, const void *p2)
152 atom_id *a1=(atom_id *)p1;
153 atom_id *a2=(atom_id *)p2;
155 return (*a1)-(*a2);
158 int n_flexible_constraints(struct gmx_constr *constr)
160 int nflexcon;
162 if (constr)
163 nflexcon = constr->nflexcon;
164 else
165 nflexcon = 0;
167 return nflexcon;
170 void too_many_constraint_warnings(int eConstrAlg,int warncount)
172 const char *abort="- aborting to avoid logfile runaway.\n"
173 "This normally happens when your system is not sufficiently equilibrated,"
174 "or if you are changing lambda too fast in free energy simulations.\n";
176 gmx_fatal(FARGS,
177 "Too many %s warnings (%d)\n"
178 "If you know what you are doing you can %s"
179 "set the environment variable GMX_MAXCONSTRWARN to -1,\n"
180 "but normally it is better to fix the problem",
181 (eConstrAlg == econtLINCS) ? "LINCS" : "SETTLE",warncount,
182 (eConstrAlg == econtLINCS) ?
183 "adjust the lincs warning threshold in your mdp file\nor " : "\n");
186 static void write_constr_pdb(const char *fn,const char *title,
187 gmx_mtop_t *mtop,
188 int start,int homenr,t_commrec *cr,
189 rvec x[],matrix box)
191 char fname[STRLEN],format[STRLEN];
192 FILE *out;
193 int dd_ac0=0,dd_ac1=0,i,ii,resnr;
194 gmx_domdec_t *dd;
195 char *anm,*resnm;
197 dd = NULL;
198 if (PAR(cr)) {
199 sprintf(fname,"%s_n%d.pdb",fn,cr->sim_nodeid);
200 if (DOMAINDECOMP(cr)) {
201 dd = cr->dd;
202 dd_get_constraint_range(dd,&dd_ac0,&dd_ac1);
203 start = 0;
204 homenr = dd_ac1;
206 } else {
207 sprintf(fname,"%s.pdb",fn);
209 sprintf(format,"%s\n",pdbformat);
211 out = gmx_fio_fopen(fname,"w");
213 fprintf(out,"TITLE %s\n",title);
214 gmx_write_pdb_box(out,-1,box);
215 for(i=start; i<start+homenr; i++) {
216 if (dd) {
217 if (i >= dd->nat_home && i < dd_ac0)
218 continue;
219 ii = dd->gatindex[i];
220 } else {
221 ii = i;
223 gmx_mtop_atominfo_global(mtop,ii,&anm,&resnr,&resnm);
224 fprintf(out,format,"ATOM",(ii+1)%100000,
225 anm,resnm,' ',(resnr+1)%10000,' ',
226 10*x[i][XX],10*x[i][YY],10*x[i][ZZ]);
228 fprintf(out,"TER\n");
230 gmx_fio_fclose(out);
233 static void dump_confs(FILE *fplog,gmx_large_int_t step,gmx_mtop_t *mtop,
234 int start,int homenr,t_commrec *cr,
235 rvec x[],rvec xprime[],matrix box)
237 char buf[256],buf2[22];
239 char *env=getenv("GMX_SUPPRESS_DUMP");
240 if (env)
241 return;
243 sprintf(buf,"step%sb",gmx_step_str(step,buf2));
244 write_constr_pdb(buf,"initial coordinates",
245 mtop,start,homenr,cr,x,box);
246 sprintf(buf,"step%sc",gmx_step_str(step,buf2));
247 write_constr_pdb(buf,"coordinates after constraining",
248 mtop,start,homenr,cr,xprime,box);
249 if (fplog)
251 fprintf(fplog,"Wrote pdb files with previous and current coordinates\n");
253 fprintf(stderr,"Wrote pdb files with previous and current coordinates\n");
256 static void pr_sortblock(FILE *fp,const char *title,int nsb,t_sortblock sb[])
258 int i;
260 fprintf(fp,"%s\n",title);
261 for(i=0; (i<nsb); i++)
262 fprintf(fp,"i: %5d, iatom: (%5d %5d %5d), blocknr: %5d\n",
263 i,sb[i].iatom[0],sb[i].iatom[1],sb[i].iatom[2],
264 sb[i].blocknr);
267 bool constrain(FILE *fplog,bool bLog,bool bEner,
268 struct gmx_constr *constr,
269 t_idef *idef,t_inputrec *ir,gmx_ekindata_t *ekind,
270 t_commrec *cr,
271 gmx_large_int_t step,int delta_step,
272 t_mdatoms *md,
273 rvec *x,rvec *xprime,rvec *min_proj,matrix box,
274 real lambda,real *dvdlambda,
275 rvec *v,tensor *vir,
276 t_nrnb *nrnb,int econq,bool bPscal,real veta, real vetanew)
278 bool bOK,bDump;
279 int start,homenr,nrend;
280 int i,j,d;
281 int ncons,error;
282 tensor rmdr;
283 rvec *vstor;
284 real invdt,vir_fac,t;
285 t_ilist *settle;
286 int nsettle;
287 t_pbc pbc;
288 char buf[22];
289 t_vetavars *vetavar;
291 if (econq == econqForceDispl && !EI_ENERGY_MINIMIZATION(ir->eI))
293 gmx_incons("constrain called for forces displacements while not doing energy minimization, can not do this while the LINCS and SETTLE constraint connection matrices are mass weighted");
296 bOK = TRUE;
297 bDump = FALSE;
299 start = md->start;
300 homenr = md->homenr;
301 nrend = start+homenr;
303 /* set constants for pressure control integration */
304 vetavar = init_vetavars(veta,vetanew,ir,ekind,bPscal);
306 if (ir->delta_t == 0)
308 invdt = 0;
310 else
312 invdt = 1/ir->delta_t;
315 if (ir->efep != efepNO && EI_DYNAMICS(ir->eI))
317 /* Set the constraint lengths for the step at which this configuration
318 * is meant to be. The invmasses should not be changed.
320 lambda += delta_step*ir->delta_lambda;
323 if (vir != NULL)
325 clear_mat(rmdr);
328 where();
329 if (constr->lincsd)
331 bOK = constrain_lincs(fplog,bLog,bEner,ir,step,constr->lincsd,md,cr,
332 x,xprime,min_proj,box,lambda,dvdlambda,
333 invdt,v,vir!=NULL,rmdr,
334 econq,nrnb,
335 constr->maxwarn,&constr->warncount_lincs);
336 if (!bOK && constr->maxwarn >= 0)
338 if (fplog != NULL)
340 fprintf(fplog,"Constraint error in algorithm %s at step %s\n",
341 econstr_names[econtLINCS],gmx_step_str(step,buf));
343 bDump = TRUE;
347 if (constr->nblocks > 0)
349 switch (econq) {
350 case (econqCoord):
351 bOK = bshakef(fplog,constr->shaked,
352 homenr,md->invmass,constr->nblocks,constr->sblock,
353 idef,ir,box,x,xprime,nrnb,
354 constr->lagr,lambda,dvdlambda,
355 invdt,v,vir!=NULL,rmdr,constr->maxwarn>=0,econq,vetavar);
356 break;
357 case (econqVeloc):
358 bOK = bshakef(fplog,constr->shaked,
359 homenr,md->invmass,constr->nblocks,constr->sblock,
360 idef,ir,box,x,min_proj,nrnb,
361 constr->lagr,lambda,dvdlambda,
362 invdt,NULL,vir!=NULL,rmdr,constr->maxwarn>=0,econq,vetavar);
363 break;
364 default:
365 gmx_fatal(FARGS,"Internal error, SHAKE called for constraining something else than coordinates");
366 break;
369 if (!bOK && constr->maxwarn >= 0)
371 if (fplog != NULL)
373 fprintf(fplog,"Constraint error in algorithm %s at step %s\n",
374 econstr_names[econtSHAKE],gmx_step_str(step,buf));
376 bDump = TRUE;
380 settle = &idef->il[F_SETTLE];
381 if (settle->nr > 0)
383 nsettle = settle->nr/2;
385 switch (econq)
387 case econqCoord:
388 csettle(constr->settled,
389 nsettle,settle->iatoms,x[0],xprime[0],
390 invdt,v[0],vir!=NULL,rmdr,&error,vetavar);
391 inc_nrnb(nrnb,eNR_SETTLE,nsettle);
392 if (v != NULL)
394 inc_nrnb(nrnb,eNR_CONSTR_V,nsettle*3);
396 if (vir != NULL)
398 inc_nrnb(nrnb,eNR_CONSTR_VIR,nsettle*3);
401 bOK = (error < 0);
402 if (!bOK && constr->maxwarn >= 0)
404 char buf[256];
405 sprintf(buf,
406 "\nt = %.3f ps: Water molecule starting at atom %d can not be "
407 "settled.\nCheck for bad contacts and/or reduce the timestep.\n",
408 ir->init_t+step*ir->delta_t,
409 ddglatnr(cr->dd,settle->iatoms[error*2+1]));
410 if (fplog)
412 fprintf(fplog,"%s",buf);
414 fprintf(stderr,"%s",buf);
415 constr->warncount_settle++;
416 if (constr->warncount_settle > constr->maxwarn)
418 too_many_constraint_warnings(-1,constr->warncount_settle);
420 bDump = TRUE;
421 break;
422 case econqVeloc:
423 case econqDeriv:
424 case econqForce:
425 case econqForceDispl:
426 settle_proj(fplog,constr->settled,econq,
427 nsettle,settle->iatoms,x,
428 xprime,min_proj,vir!=NULL,rmdr,vetavar);
429 /* This is an overestimate */
430 inc_nrnb(nrnb,eNR_SETTLE,nsettle);
431 break;
432 case econqDeriv_FlexCon:
433 /* Nothing to do, since the are no flexible constraints in settles */
434 break;
435 default:
436 gmx_incons("Unknown constraint quantity for settle");
441 free_vetavars(vetavar);
443 if (vir != NULL)
445 switch (econq)
447 case econqCoord:
448 vir_fac = 0.5/(ir->delta_t*ir->delta_t);
449 break;
450 case econqVeloc:
451 vir_fac = 0.5/ir->delta_t;
452 break;
453 case econqForce:
454 case econqForceDispl:
455 vir_fac = 0.5;
456 break;
457 default:
458 vir_fac = 0;
459 gmx_incons("Unsupported constraint quantity for virial");
462 if (EI_VV(ir->eI))
464 vir_fac *= 2; /* only constraining over half the distance here */
466 for(i=0; i<DIM; i++)
468 for(j=0; j<DIM; j++)
470 (*vir)[i][j] = vir_fac*rmdr[i][j];
475 if (bDump)
477 dump_confs(fplog,step,constr->warn_mtop,start,homenr,cr,x,xprime,box);
480 if (econq == econqCoord)
482 if (ir->ePull == epullCONSTRAINT)
484 if (EI_DYNAMICS(ir->eI))
486 t = ir->init_t + (step + delta_step)*ir->delta_t;
488 else
490 t = ir->init_t;
492 set_pbc(&pbc,ir->ePBC,box);
493 pull_constraint(ir->pull,md,&pbc,cr,ir->delta_t,t,x,xprime,v,*vir);
495 if (constr->ed && delta_step > 0)
497 /* apply the essential dynamcs constraints here */
498 do_edsam(ir,step,md,cr,xprime,v,box,constr->ed);
502 return bOK;
505 real *constr_rmsd_data(struct gmx_constr *constr)
507 if (constr->lincsd)
508 return lincs_rmsd_data(constr->lincsd);
509 else
510 return NULL;
513 real constr_rmsd(struct gmx_constr *constr,bool bSD2)
515 if (constr->lincsd)
516 return lincs_rmsd(constr->lincsd,bSD2);
517 else
518 return 0;
521 static void make_shake_sblock_pd(struct gmx_constr *constr,
522 t_idef *idef,t_mdatoms *md)
524 int i,j,m,ncons;
525 int bstart,bnr;
526 t_blocka sblocks;
527 t_sortblock *sb;
528 t_iatom *iatom;
529 atom_id *inv_sblock;
531 /* Since we are processing the local topology,
532 * the F_CONSTRNC ilist has been concatenated to the F_CONSTR ilist.
534 ncons = idef->il[F_CONSTR].nr/3;
536 init_blocka(&sblocks);
537 gen_sblocks(NULL,md->start,md->start+md->homenr,idef,&sblocks,FALSE);
540 bstart=(idef->nodeid > 0) ? blocks->multinr[idef->nodeid-1] : 0;
541 nblocks=blocks->multinr[idef->nodeid] - bstart;
543 bstart = 0;
544 constr->nblocks = sblocks.nr;
545 if (debug)
546 fprintf(debug,"ncons: %d, bstart: %d, nblocks: %d\n",
547 ncons,bstart,constr->nblocks);
549 /* Calculate block number for each atom */
550 inv_sblock = make_invblocka(&sblocks,md->nr);
552 done_blocka(&sblocks);
554 /* Store the block number in temp array and
555 * sort the constraints in order of the sblock number
556 * and the atom numbers, really sorting a segment of the array!
558 #ifdef DEBUGIDEF
559 pr_idef(fplog,0,"Before Sort",idef);
560 #endif
561 iatom=idef->il[F_CONSTR].iatoms;
562 snew(sb,ncons);
563 for(i=0; (i<ncons); i++,iatom+=3) {
564 for(m=0; (m<3); m++)
565 sb[i].iatom[m] = iatom[m];
566 sb[i].blocknr = inv_sblock[iatom[1]];
569 /* Now sort the blocks */
570 if (debug) {
571 pr_sortblock(debug,"Before sorting",ncons,sb);
572 fprintf(debug,"Going to sort constraints\n");
575 qsort(sb,ncons,(size_t)sizeof(*sb),pcomp);
577 if (debug) {
578 pr_sortblock(debug,"After sorting",ncons,sb);
581 iatom=idef->il[F_CONSTR].iatoms;
582 for(i=0; (i<ncons); i++,iatom+=3)
583 for(m=0; (m<3); m++)
584 iatom[m]=sb[i].iatom[m];
585 #ifdef DEBUGIDEF
586 pr_idef(fplog,0,"After Sort",idef);
587 #endif
589 j=0;
590 snew(constr->sblock,constr->nblocks+1);
591 bnr=-2;
592 for(i=0; (i<ncons); i++) {
593 if (sb[i].blocknr != bnr) {
594 bnr=sb[i].blocknr;
595 constr->sblock[j++]=3*i;
598 /* Last block... */
599 constr->sblock[j++] = 3*ncons;
601 if (j != (constr->nblocks+1)) {
602 fprintf(stderr,"bstart: %d\n",bstart);
603 fprintf(stderr,"j: %d, nblocks: %d, ncons: %d\n",
604 j,constr->nblocks,ncons);
605 for(i=0; (i<ncons); i++)
606 fprintf(stderr,"i: %5d sb[i].blocknr: %5u\n",i,sb[i].blocknr);
607 for(j=0; (j<=constr->nblocks); j++)
608 fprintf(stderr,"sblock[%3d]=%5d\n",j,(int)constr->sblock[j]);
609 gmx_fatal(FARGS,"DEATH HORROR: "
610 "sblocks does not match idef->il[F_CONSTR]");
612 sfree(sb);
613 sfree(inv_sblock);
616 static void make_shake_sblock_dd(struct gmx_constr *constr,
617 t_ilist *ilcon,t_block *cgs,
618 gmx_domdec_t *dd)
620 int ncons,c,cg;
621 t_iatom *iatom;
623 if (dd->ncg_home+1 > constr->sblock_nalloc) {
624 constr->sblock_nalloc = over_alloc_dd(dd->ncg_home+1);
625 srenew(constr->sblock,constr->sblock_nalloc);
628 ncons = ilcon->nr/3;
629 iatom = ilcon->iatoms;
630 constr->nblocks = 0;
631 cg = 0;
632 for(c=0; c<ncons; c++) {
633 if (c == 0 || iatom[1] >= cgs->index[cg+1]) {
634 constr->sblock[constr->nblocks++] = 3*c;
635 while (iatom[1] >= cgs->index[cg+1])
636 cg++;
638 iatom += 3;
640 constr->sblock[constr->nblocks] = 3*ncons;
643 t_blocka make_at2con(int start,int natoms,
644 t_ilist *ilist,t_iparams *iparams,
645 bool bDynamics,int *nflexiblecons)
647 int *count,ncon,con,con_tot,nflexcon,ftype,i,a;
648 t_iatom *ia;
649 t_blocka at2con;
650 bool bFlexCon;
652 snew(count,natoms);
653 nflexcon = 0;
654 for(ftype=F_CONSTR; ftype<=F_CONSTRNC; ftype++) {
655 ncon = ilist[ftype].nr/3;
656 ia = ilist[ftype].iatoms;
657 for(con=0; con<ncon; con++) {
658 bFlexCon = (iparams[ia[0]].constr.dA == 0 &&
659 iparams[ia[0]].constr.dB == 0);
660 if (bFlexCon)
661 nflexcon++;
662 if (bDynamics || !bFlexCon) {
663 for(i=1; i<3; i++) {
664 a = ia[i] - start;
665 count[a]++;
668 ia += 3;
671 *nflexiblecons = nflexcon;
673 at2con.nr = natoms;
674 at2con.nalloc_index = at2con.nr+1;
675 snew(at2con.index,at2con.nalloc_index);
676 at2con.index[0] = 0;
677 for(a=0; a<natoms; a++) {
678 at2con.index[a+1] = at2con.index[a] + count[a];
679 count[a] = 0;
681 at2con.nra = at2con.index[natoms];
682 at2con.nalloc_a = at2con.nra;
683 snew(at2con.a,at2con.nalloc_a);
685 /* The F_CONSTRNC constraints have constraint numbers
686 * that continue after the last F_CONSTR constraint.
688 con_tot = 0;
689 for(ftype=F_CONSTR; ftype<=F_CONSTRNC; ftype++) {
690 ncon = ilist[ftype].nr/3;
691 ia = ilist[ftype].iatoms;
692 for(con=0; con<ncon; con++) {
693 bFlexCon = (iparams[ia[0]].constr.dA == 0 &&
694 iparams[ia[0]].constr.dB == 0);
695 if (bDynamics || !bFlexCon) {
696 for(i=1; i<3; i++) {
697 a = ia[i] - start;
698 at2con.a[at2con.index[a]+count[a]++] = con_tot;
701 con_tot++;
702 ia += 3;
706 sfree(count);
708 return at2con;
711 void set_constraints(struct gmx_constr *constr,
712 gmx_localtop_t *top,t_inputrec *ir,
713 t_mdatoms *md,t_commrec *cr)
715 t_idef *idef;
716 int ncons;
717 t_ilist *settle;
718 int iO,iH;
720 idef = &top->idef;
722 if (constr->ncon_tot > 0)
724 /* We are using the local topology,
725 * so there are only F_CONSTR constraints.
727 ncons = idef->il[F_CONSTR].nr/3;
729 /* With DD we might also need to call LINCS with ncons=0 for
730 * communicating coordinates to other nodes that do have constraints.
732 if (ir->eConstrAlg == econtLINCS)
734 set_lincs(idef,md,EI_DYNAMICS(ir->eI),cr,constr->lincsd);
736 if (ir->eConstrAlg == econtSHAKE)
738 if (cr->dd)
740 make_shake_sblock_dd(constr,&idef->il[F_CONSTR],&top->cgs,cr->dd);
742 else
744 make_shake_sblock_pd(constr,idef,md);
746 if (ncons > constr->lagr_nalloc)
748 constr->lagr_nalloc = over_alloc_dd(ncons);
749 srenew(constr->lagr,constr->lagr_nalloc);
752 constr->shaked = shake_init();
756 if (idef->il[F_SETTLE].nr > 0 && constr->settled == NULL)
758 settle = &idef->il[F_SETTLE];
759 iO = settle->iatoms[1];
760 iH = settle->iatoms[1]+1;
761 constr->settled =
762 settle_init(md->massT[iO],md->massT[iH],
763 md->invmass[iO],md->invmass[iH],
764 idef->iparams[settle->iatoms[0]].settle.doh,
765 idef->iparams[settle->iatoms[0]].settle.dhh);
768 /* Make a selection of the local atoms for essential dynamics */
769 if (constr->ed && cr->dd)
771 dd_make_local_ed_indices(cr->dd,constr->ed);
775 static void constr_recur(t_blocka *at2con,
776 t_ilist *ilist,t_iparams *iparams,bool bTopB,
777 int at,int depth,int nc,int *path,
778 real r0,real r1,real *r2max,
779 int *count)
781 int ncon1;
782 t_iatom *ia1,*ia2;
783 int c,con,a1;
784 bool bUse;
785 t_iatom *ia;
786 real len,rn0,rn1;
788 (*count)++;
790 ncon1 = ilist[F_CONSTR].nr/3;
791 ia1 = ilist[F_CONSTR].iatoms;
792 ia2 = ilist[F_CONSTRNC].iatoms;
794 /* Loop over all constraints connected to this atom */
795 for(c=at2con->index[at]; c<at2con->index[at+1]; c++) {
796 con = at2con->a[c];
797 /* Do not walk over already used constraints */
798 bUse = TRUE;
799 for(a1=0; a1<depth; a1++) {
800 if (con == path[a1])
801 bUse = FALSE;
803 if (bUse) {
804 ia = constr_iatomptr(ncon1,ia1,ia2,con);
805 /* Flexible constraints currently have length 0, which is incorrect */
806 if (!bTopB)
807 len = iparams[ia[0]].constr.dA;
808 else
809 len = iparams[ia[0]].constr.dB;
810 /* In the worst case the bond directions alternate */
811 if (nc % 2 == 0) {
812 rn0 = r0 + len;
813 rn1 = r1;
814 } else {
815 rn0 = r0;
816 rn1 = r1 + len;
818 /* Assume angles of 120 degrees between all bonds */
819 if (rn0*rn0 + rn1*rn1 + rn0*rn1 > *r2max) {
820 *r2max = rn0*rn0 + rn1*rn1 + r0*rn1;
821 if (debug) {
822 fprintf(debug,"Found longer constraint distance: r0 %5.3f r1 %5.3f rmax %5.3f\n", rn0,rn1,sqrt(*r2max));
823 for(a1=0; a1<depth; a1++)
824 fprintf(debug," %d %5.3f",
825 path[a1],
826 iparams[constr_iatomptr(ncon1,ia1,ia2,con)[0]].constr.dA);
827 fprintf(debug," %d %5.3f\n",con,len);
830 /* Limit the number of recursions to 1000*nc,
831 * so a call does not take more than a second,
832 * even for highly connected systems.
834 if (depth + 1 < nc && *count < 1000*nc) {
835 if (ia[1] == at)
836 a1 = ia[2];
837 else
838 a1 = ia[1];
839 /* Recursion */
840 path[depth] = con;
841 constr_recur(at2con,ilist,iparams,
842 bTopB,a1,depth+1,nc,path,rn0,rn1,r2max,count);
843 path[depth] = -1;
849 static real constr_r_max_moltype(FILE *fplog,
850 gmx_moltype_t *molt,t_iparams *iparams,
851 t_inputrec *ir)
853 int natoms,nflexcon,*path,at,count;
855 t_blocka at2con;
856 real r0,r1,r2maxA,r2maxB,rmax,lam0,lam1;
858 if (molt->ilist[F_CONSTR].nr == 0 &&
859 molt->ilist[F_CONSTRNC].nr == 0) {
860 return 0;
863 natoms = molt->atoms.nr;
865 at2con = make_at2con(0,natoms,molt->ilist,iparams,
866 EI_DYNAMICS(ir->eI),&nflexcon);
867 snew(path,1+ir->nProjOrder);
868 for(at=0; at<1+ir->nProjOrder; at++)
869 path[at] = -1;
871 r2maxA = 0;
872 for(at=0; at<natoms; at++) {
873 r0 = 0;
874 r1 = 0;
876 count = 0;
877 constr_recur(&at2con,molt->ilist,iparams,
878 FALSE,at,0,1+ir->nProjOrder,path,r0,r1,&r2maxA,&count);
880 if (ir->efep == efepNO) {
881 rmax = sqrt(r2maxA);
882 } else {
883 r2maxB = 0;
884 for(at=0; at<natoms; at++) {
885 r0 = 0;
886 r1 = 0;
887 count = 0;
888 constr_recur(&at2con,molt->ilist,iparams,
889 TRUE,at,0,1+ir->nProjOrder,path,r0,r1,&r2maxB,&count);
891 lam0 = ir->init_lambda;
892 if (EI_DYNAMICS(ir->eI))
893 lam0 += ir->init_step*ir->delta_lambda;
894 rmax = (1 - lam0)*sqrt(r2maxA) + lam0*sqrt(r2maxB);
895 if (EI_DYNAMICS(ir->eI)) {
896 lam1 = ir->init_lambda + (ir->init_step + ir->nsteps)*ir->delta_lambda;
897 rmax = max(rmax,(1 - lam1)*sqrt(r2maxA) + lam1*sqrt(r2maxB));
901 done_blocka(&at2con);
902 sfree(path);
904 return rmax;
907 real constr_r_max(FILE *fplog,gmx_mtop_t *mtop,t_inputrec *ir)
909 int mt;
910 real rmax;
912 rmax = 0;
913 for(mt=0; mt<mtop->nmoltype; mt++) {
914 rmax = max(rmax,
915 constr_r_max_moltype(fplog,&mtop->moltype[mt],
916 mtop->ffparams.iparams,ir));
919 if (fplog)
920 fprintf(fplog,"Maximum distance for %d constraints, at 120 deg. angles, all-trans: %.3f nm\n",1+ir->nProjOrder,rmax);
922 return rmax;
925 gmx_constr_t init_constraints(FILE *fplog,
926 gmx_mtop_t *mtop,t_inputrec *ir,
927 gmx_edsam_t ed,t_state *state,
928 t_commrec *cr)
930 int ncon,nset,nmol,settle_type,i,natoms,mt,nflexcon;
931 struct gmx_constr *constr;
932 char *env;
933 t_ilist *ilist;
934 gmx_mtop_ilistloop_t iloop;
936 ncon =
937 gmx_mtop_ftype_count(mtop,F_CONSTR) +
938 gmx_mtop_ftype_count(mtop,F_CONSTRNC);
939 nset = gmx_mtop_ftype_count(mtop,F_SETTLE);
941 if (ncon+nset == 0 && ir->ePull != epullCONSTRAINT && ed == NULL)
943 return NULL;
946 snew(constr,1);
948 constr->ncon_tot = ncon;
949 constr->nflexcon = 0;
950 if (ncon > 0)
952 constr->n_at2con_mt = mtop->nmoltype;
953 snew(constr->at2con_mt,constr->n_at2con_mt);
954 for(mt=0; mt<mtop->nmoltype; mt++)
956 constr->at2con_mt[mt] = make_at2con(0,mtop->moltype[mt].atoms.nr,
957 mtop->moltype[mt].ilist,
958 mtop->ffparams.iparams,
959 EI_DYNAMICS(ir->eI),&nflexcon);
960 for(i=0; i<mtop->nmolblock; i++)
962 if (mtop->molblock[i].type == mt)
964 constr->nflexcon += mtop->molblock[i].nmol*nflexcon;
969 if (constr->nflexcon > 0)
971 if (fplog)
973 fprintf(fplog,"There are %d flexible constraints\n",
974 constr->nflexcon);
975 if (ir->fc_stepsize == 0)
977 fprintf(fplog,"\n"
978 "WARNING: step size for flexible constraining = 0\n"
979 " All flexible constraints will be rigid.\n"
980 " Will try to keep all flexible constraints at their original length,\n"
981 " but the lengths may exhibit some drift.\n\n");
982 constr->nflexcon = 0;
985 if (constr->nflexcon > 0)
987 please_cite(fplog,"Hess2002");
991 if (ir->eConstrAlg == econtLINCS)
993 constr->lincsd = init_lincs(fplog,mtop,
994 constr->nflexcon,constr->at2con_mt,
995 DOMAINDECOMP(cr) && cr->dd->bInterCGcons,
996 ir->nLincsIter,ir->nProjOrder);
999 if (ir->eConstrAlg == econtSHAKE) {
1000 if (DOMAINDECOMP(cr) && cr->dd->bInterCGcons)
1002 gmx_fatal(FARGS,"SHAKE is not supported with domain decomposition and constraint that cross charge group boundaries, use LINCS");
1004 if (constr->nflexcon)
1006 gmx_fatal(FARGS,"For this system also velocities and/or forces need to be constrained, this can not be done with SHAKE, you should select LINCS");
1008 please_cite(fplog,"Ryckaert77a");
1009 if (ir->bShakeSOR)
1011 please_cite(fplog,"Barth95a");
1016 if (nset > 0) {
1017 please_cite(fplog,"Miyamoto92a");
1019 /* Check that we have only one settle type */
1020 settle_type = -1;
1021 iloop = gmx_mtop_ilistloop_init(mtop);
1022 while (gmx_mtop_ilistloop_next(iloop,&ilist,&nmol))
1024 for (i=0; i<ilist[F_SETTLE].nr; i+=2)
1026 if (settle_type == -1)
1028 settle_type = ilist[F_SETTLE].iatoms[i];
1030 else if (ilist[F_SETTLE].iatoms[i] != settle_type)
1032 gmx_fatal(FARGS,"More than one settle type.\n"
1033 "Suggestion: change the least use settle constraints into 3 normal constraints.");
1039 constr->maxwarn = 999;
1040 env = getenv("GMX_MAXCONSTRWARN");
1041 if (env)
1043 constr->maxwarn = 0;
1044 sscanf(env,"%d",&constr->maxwarn);
1045 if (fplog)
1047 fprintf(fplog,
1048 "Setting the maximum number of constraint warnings to %d\n",
1049 constr->maxwarn);
1051 if (MASTER(cr))
1053 fprintf(stderr,
1054 "Setting the maximum number of constraint warnings to %d\n",
1055 constr->maxwarn);
1058 if (constr->maxwarn < 0 && fplog)
1060 fprintf(fplog,"maxwarn < 0, will not stop on constraint errors\n");
1062 constr->warncount_lincs = 0;
1063 constr->warncount_settle = 0;
1065 /* Initialize the essential dynamics sampling.
1066 * Put the pointer to the ED struct in constr */
1067 constr->ed = ed;
1068 if (ed != NULL)
1070 init_edsam(mtop,ir,cr,ed,state->x,state->box);
1073 constr->warn_mtop = mtop;
1075 return constr;
1078 t_blocka *atom2constraints_moltype(gmx_constr_t constr)
1080 return constr->at2con_mt;
1084 bool inter_charge_group_constraints(gmx_mtop_t *mtop)
1086 const gmx_moltype_t *molt;
1087 const t_block *cgs;
1088 const t_ilist *il;
1089 int mb;
1090 int nat,*at2cg,cg,a,ftype,i;
1091 bool bInterCG;
1093 bInterCG = FALSE;
1094 for(mb=0; mb<mtop->nmolblock && !bInterCG; mb++) {
1095 molt = &mtop->moltype[mtop->molblock[mb].type];
1097 if (molt->ilist[F_CONSTR].nr > 0 ||
1098 molt->ilist[F_CONSTRNC].nr > 0) {
1099 cgs = &molt->cgs;
1100 snew(at2cg,molt->atoms.nr);
1101 for(cg=0; cg<cgs->nr; cg++) {
1102 for(a=cgs->index[cg]; a<cgs->index[cg+1]; a++)
1103 at2cg[a] = cg;
1106 for(ftype=F_CONSTR; ftype<=F_CONSTRNC; ftype++) {
1107 il = &molt->ilist[ftype];
1108 for(i=0; i<il->nr && !bInterCG; i+=3) {
1109 if (at2cg[il->iatoms[i+1]] != at2cg[il->iatoms[i+2]])
1110 bInterCG = TRUE;
1113 sfree(at2cg);
1117 return bInterCG;