1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
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
34 * GROwing Monsters And Cloning Shrimps
56 #include "mtop_util.h"
59 #include "gmx_omp_nthreads.h"
61 typedef struct gmx_constr {
62 int ncon_tot; /* The total number of constraints */
63 int nflexcon; /* The number of flexible constraints */
64 int n_at2con_mt; /* The size of at2con = #moltypes */
65 t_blocka *at2con_mt; /* A list of atoms to constraints */
66 int n_at2settle_mt; /* The size of at2settle = #moltypes */
67 int **at2settle_mt; /* A list of atoms to settles */
68 gmx_bool bInterCGsettles;
69 gmx_lincsdata_t lincsd; /* LINCS data */
70 gmx_shakedata_t shaked; /* SHAKE data */
71 gmx_settledata_t settled; /* SETTLE data */
72 int nblocks; /* The number of SHAKE blocks */
73 int *sblock; /* The SHAKE blocks */
74 int sblock_nalloc;/* The allocation size of sblock */
75 real *lagr; /* Lagrange multipliers for SHAKE */
76 int lagr_nalloc; /* The allocation size of lagr */
77 int maxwarn; /* The maximum number of warnings */
80 gmx_edsam_t ed; /* The essential dynamics data */
82 tensor *vir_r_m_dr_th;/* Thread local working data */
83 int *settle_error; /* Thread local working data */
85 gmx_mtop_t *warn_mtop; /* Only used for printing warnings */
93 static void *init_vetavars(t_vetavars *vars,
94 gmx_bool constr_deriv,
95 real veta,real vetanew, t_inputrec *ir, gmx_ekindata_t *ekind, gmx_bool bPscal)
100 /* first, set the alpha integrator variable */
101 if ((ir->opts.nrdf[0] > 0) && bPscal)
103 vars->alpha = 1.0 + DIM/((double)ir->opts.nrdf[0]);
107 g = 0.5*veta*ir->delta_t;
108 vars->rscale = exp(g)*series_sinhx(g);
109 g = -0.25*vars->alpha*veta*ir->delta_t;
110 vars->vscale = exp(g)*series_sinhx(g);
111 vars->rvscale = vars->vscale*vars->rscale;
112 vars->veta = vetanew;
116 snew(vars->vscale_nhc,ir->opts.ngtc);
117 if ((ekind==NULL) || (!bPscal))
119 for (i=0;i<ir->opts.ngtc;i++)
121 vars->vscale_nhc[i] = 1;
126 for (i=0;i<ir->opts.ngtc;i++)
128 vars->vscale_nhc[i] = ekind->tcstat[i].vscale_nhc;
134 vars->vscale_nhc = NULL;
140 static void free_vetavars(t_vetavars *vars)
142 if (vars->vscale_nhc != NULL)
144 sfree(vars->vscale_nhc);
148 static int pcomp(const void *p1, const void *p2)
151 atom_id min1,min2,max1,max2;
152 t_sortblock *a1=(t_sortblock *)p1;
153 t_sortblock *a2=(t_sortblock *)p2;
155 db=a1->blocknr-a2->blocknr;
160 min1=min(a1->iatom[1],a1->iatom[2]);
161 max1=max(a1->iatom[1],a1->iatom[2]);
162 min2=min(a2->iatom[1],a2->iatom[2]);
163 max2=max(a2->iatom[1],a2->iatom[2]);
171 int n_flexible_constraints(struct gmx_constr *constr)
176 nflexcon = constr->nflexcon;
183 void too_many_constraint_warnings(int eConstrAlg,int warncount)
185 const char *abort="- aborting to avoid logfile runaway.\n"
186 "This normally happens when your system is not sufficiently equilibrated,"
187 "or if you are changing lambda too fast in free energy simulations.\n";
190 "Too many %s warnings (%d)\n"
191 "If you know what you are doing you can %s"
192 "set the environment variable GMX_MAXCONSTRWARN to -1,\n"
193 "but normally it is better to fix the problem",
194 (eConstrAlg == econtLINCS) ? "LINCS" : "SETTLE",warncount,
195 (eConstrAlg == econtLINCS) ?
196 "adjust the lincs warning threshold in your mdp file\nor " : "\n");
199 static void write_constr_pdb(const char *fn,const char *title,
201 int start,int homenr,t_commrec *cr,
204 char fname[STRLEN],format[STRLEN];
206 int dd_ac0=0,dd_ac1=0,i,ii,resnr;
213 sprintf(fname,"%s_n%d.pdb",fn,cr->sim_nodeid);
214 if (DOMAINDECOMP(cr))
217 dd_get_constraint_range(dd,&dd_ac0,&dd_ac1);
224 sprintf(fname,"%s.pdb",fn);
226 sprintf(format,"%s\n",get_pdbformat());
228 out = gmx_fio_fopen(fname,"w");
230 fprintf(out,"TITLE %s\n",title);
231 gmx_write_pdb_box(out,-1,box);
232 for(i=start; i<start+homenr; i++)
236 if (i >= dd->nat_home && i < dd_ac0)
240 ii = dd->gatindex[i];
246 gmx_mtop_atominfo_global(mtop,ii,&anm,&resnr,&resnm);
247 fprintf(out,format,"ATOM",(ii+1)%100000,
248 anm,resnm,' ',resnr%10000,' ',
249 10*x[i][XX],10*x[i][YY],10*x[i][ZZ]);
251 fprintf(out,"TER\n");
256 static void dump_confs(FILE *fplog,gmx_large_int_t step,gmx_mtop_t *mtop,
257 int start,int homenr,t_commrec *cr,
258 rvec x[],rvec xprime[],matrix box)
260 char buf[256],buf2[22];
262 char *env=getenv("GMX_SUPPRESS_DUMP");
266 sprintf(buf,"step%sb",gmx_step_str(step,buf2));
267 write_constr_pdb(buf,"initial coordinates",
268 mtop,start,homenr,cr,x,box);
269 sprintf(buf,"step%sc",gmx_step_str(step,buf2));
270 write_constr_pdb(buf,"coordinates after constraining",
271 mtop,start,homenr,cr,xprime,box);
274 fprintf(fplog,"Wrote pdb files with previous and current coordinates\n");
276 fprintf(stderr,"Wrote pdb files with previous and current coordinates\n");
279 static void pr_sortblock(FILE *fp,const char *title,int nsb,t_sortblock sb[])
283 fprintf(fp,"%s\n",title);
284 for(i=0; (i<nsb); i++)
285 fprintf(fp,"i: %5d, iatom: (%5d %5d %5d), blocknr: %5d\n",
286 i,sb[i].iatom[0],sb[i].iatom[1],sb[i].iatom[2],
290 gmx_bool constrain(FILE *fplog,gmx_bool bLog,gmx_bool bEner,
291 struct gmx_constr *constr,
292 t_idef *idef,t_inputrec *ir,gmx_ekindata_t *ekind,
294 gmx_large_int_t step,int delta_step,
296 rvec *x,rvec *xprime,rvec *min_proj,
297 gmx_bool bMolPBC,matrix box,
298 real lambda,real *dvdlambda,
300 t_nrnb *nrnb,int econq,gmx_bool bPscal,
301 real veta, real vetanew)
304 int start,homenr,nrend;
306 int ncons,settle_error;
309 real invdt,vir_fac,t;
317 if (econq == econqForceDispl && !EI_ENERGY_MINIMIZATION(ir->eI))
319 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");
327 nrend = start+homenr;
329 /* set constants for pressure control integration */
330 init_vetavars(&vetavar,econq!=econqCoord,
331 veta,vetanew,ir,ekind,bPscal);
333 if (ir->delta_t == 0)
339 invdt = 1/ir->delta_t;
342 if (ir->efep != efepNO && EI_DYNAMICS(ir->eI))
344 /* Set the constraint lengths for the step at which this configuration
345 * is meant to be. The invmasses should not be changed.
347 lambda += delta_step*ir->fepvals->delta_lambda;
352 clear_mat(vir_r_m_dr);
357 settle = &idef->il[F_SETTLE];
358 nsettle = settle->nr/(1+NRAL(F_SETTLE));
362 nth = gmx_omp_nthreads_get(emntSETTLE);
369 if (nth > 1 && constr->vir_r_m_dr_th == NULL)
371 snew(constr->vir_r_m_dr_th,nth);
372 snew(constr->settle_error,nth);
377 /* We do not need full pbc when constraints do not cross charge groups,
378 * i.e. when dd->constraint_comm==NULL.
379 * Note that PBC for constraints is different from PBC for bondeds.
380 * For constraints there is both forward and backward communication.
382 if (ir->ePBC != epbcNONE &&
383 (cr->dd || bMolPBC) && !(cr->dd && cr->dd->constraint_comm==NULL))
385 /* With pbc=screw the screw has been changed to a shift
386 * by the constraint coordinate communication routine,
387 * so that here we can use normal pbc.
389 pbc_null = set_pbc_dd(&pbc,ir->ePBC,cr->dd,FALSE,box);
396 /* Communicate the coordinates required for the non-local constraints
397 * for LINCS and/or SETTLE.
401 dd_move_x_constraints(cr->dd,box,x,xprime);
403 else if (PARTDECOMP(cr))
405 pd_move_x_constraints(cr,x,xprime);
408 if (constr->lincsd != NULL)
410 bOK = constrain_lincs(fplog,bLog,bEner,ir,step,constr->lincsd,md,cr,
412 box,pbc_null,lambda,dvdlambda,
413 invdt,v,vir!=NULL,vir_r_m_dr,
415 constr->maxwarn,&constr->warncount_lincs);
416 if (!bOK && constr->maxwarn >= 0)
420 fprintf(fplog,"Constraint error in algorithm %s at step %s\n",
421 econstr_names[econtLINCS],gmx_step_str(step,buf));
427 if (constr->nblocks > 0)
431 bOK = bshakef(fplog,constr->shaked,
432 homenr,md->invmass,constr->nblocks,constr->sblock,
433 idef,ir,x,xprime,nrnb,
434 constr->lagr,lambda,dvdlambda,
435 invdt,v,vir!=NULL,vir_r_m_dr,
436 constr->maxwarn>=0,econq,&vetavar);
439 bOK = bshakef(fplog,constr->shaked,
440 homenr,md->invmass,constr->nblocks,constr->sblock,
441 idef,ir,x,min_proj,nrnb,
442 constr->lagr,lambda,dvdlambda,
443 invdt,NULL,vir!=NULL,vir_r_m_dr,
444 constr->maxwarn>=0,econq,&vetavar);
447 gmx_fatal(FARGS,"Internal error, SHAKE called for constraining something else than coordinates");
451 if (!bOK && constr->maxwarn >= 0)
455 fprintf(fplog,"Constraint error in algorithm %s at step %s\n",
456 econstr_names[econtSHAKE],gmx_step_str(step,buf));
464 int calcvir_atom_end;
468 calcvir_atom_end = 0;
472 calcvir_atom_end = md->start + md->homenr;
478 #pragma omp parallel for num_threads(nth) schedule(static)
479 for(th=0; th<nth; th++)
485 clear_mat(constr->vir_r_m_dr_th[th]);
488 start_th = (nsettle* th )/nth;
489 end_th = (nsettle*(th+1))/nth;
490 if (start_th >= 0 && end_th - start_th > 0)
492 csettle(constr->settled,
494 settle->iatoms+start_th*(1+NRAL(F_SETTLE)),
497 invdt,v?v[0]:NULL,calcvir_atom_end,
498 th == 0 ? vir_r_m_dr : constr->vir_r_m_dr_th[th],
499 th == 0 ? &settle_error : &constr->settle_error[th],
503 inc_nrnb(nrnb,eNR_SETTLE,nsettle);
506 inc_nrnb(nrnb,eNR_CONSTR_V,nsettle*3);
510 inc_nrnb(nrnb,eNR_CONSTR_VIR,nsettle*3);
516 case econqForceDispl:
517 #pragma omp parallel for num_threads(nth) schedule(static)
518 for(th=0; th<nth; th++)
524 clear_mat(constr->vir_r_m_dr_th[th]);
527 start_th = (nsettle* th )/nth;
528 end_th = (nsettle*(th+1))/nth;
530 if (start_th >= 0 && end_th - start_th > 0)
532 settle_proj(fplog,constr->settled,econq,
534 settle->iatoms+start_th*(1+NRAL(F_SETTLE)),
537 xprime,min_proj,calcvir_atom_end,
538 th == 0 ? vir_r_m_dr : constr->vir_r_m_dr_th[th],
542 /* This is an overestimate */
543 inc_nrnb(nrnb,eNR_SETTLE,nsettle);
545 case econqDeriv_FlexCon:
546 /* Nothing to do, since the are no flexible constraints in settles */
549 gmx_incons("Unknown constraint quantity for settle");
555 /* Combine virial and error info of the other threads */
558 m_add(vir_r_m_dr,constr->vir_r_m_dr_th[i],vir_r_m_dr);
559 settle_error = constr->settle_error[i];
562 if (econq == econqCoord && settle_error >= 0)
565 if (constr->maxwarn >= 0)
569 "\nstep " gmx_large_int_pfmt ": Water molecule starting at atom %d can not be "
570 "settled.\nCheck for bad contacts and/or reduce the timestep if appropriate.\n",
571 step,ddglatnr(cr->dd,settle->iatoms[settle_error*(1+NRAL(F_SETTLE))+1]));
574 fprintf(fplog,"%s",buf);
576 fprintf(stderr,"%s",buf);
577 constr->warncount_settle++;
578 if (constr->warncount_settle > constr->maxwarn)
580 too_many_constraint_warnings(-1,constr->warncount_settle);
587 free_vetavars(&vetavar);
594 vir_fac = 0.5/(ir->delta_t*ir->delta_t);
597 vir_fac = 0.5/ir->delta_t;
600 case econqForceDispl:
605 gmx_incons("Unsupported constraint quantity for virial");
610 vir_fac *= 2; /* only constraining over half the distance here */
616 (*vir)[i][j] = vir_fac*vir_r_m_dr[i][j];
623 dump_confs(fplog,step,constr->warn_mtop,start,homenr,cr,x,xprime,box);
626 if (econq == econqCoord)
628 if (ir->ePull == epullCONSTRAINT)
630 if (EI_DYNAMICS(ir->eI))
632 t = ir->init_t + (step + delta_step)*ir->delta_t;
638 set_pbc(&pbc,ir->ePBC,box);
639 pull_constraint(ir->pull,md,&pbc,cr,ir->delta_t,t,x,xprime,v,*vir);
641 if (constr->ed && delta_step > 0)
643 /* apply the essential dynamcs constraints here */
644 do_edsam(ir,step,cr,xprime,v,box,constr->ed);
651 real *constr_rmsd_data(struct gmx_constr *constr)
654 return lincs_rmsd_data(constr->lincsd);
659 real constr_rmsd(struct gmx_constr *constr,gmx_bool bSD2)
662 return lincs_rmsd(constr->lincsd,bSD2);
667 static void make_shake_sblock_pd(struct gmx_constr *constr,
668 t_idef *idef,t_mdatoms *md)
677 /* Since we are processing the local topology,
678 * the F_CONSTRNC ilist has been concatenated to the F_CONSTR ilist.
680 ncons = idef->il[F_CONSTR].nr/3;
682 init_blocka(&sblocks);
683 gen_sblocks(NULL,md->start,md->start+md->homenr,idef,&sblocks,FALSE);
686 bstart=(idef->nodeid > 0) ? blocks->multinr[idef->nodeid-1] : 0;
687 nblocks=blocks->multinr[idef->nodeid] - bstart;
690 constr->nblocks = sblocks.nr;
692 fprintf(debug,"ncons: %d, bstart: %d, nblocks: %d\n",
693 ncons,bstart,constr->nblocks);
695 /* Calculate block number for each atom */
696 inv_sblock = make_invblocka(&sblocks,md->nr);
698 done_blocka(&sblocks);
700 /* Store the block number in temp array and
701 * sort the constraints in order of the sblock number
702 * and the atom numbers, really sorting a segment of the array!
705 pr_idef(fplog,0,"Before Sort",idef);
707 iatom=idef->il[F_CONSTR].iatoms;
709 for(i=0; (i<ncons); i++,iatom+=3) {
711 sb[i].iatom[m] = iatom[m];
712 sb[i].blocknr = inv_sblock[iatom[1]];
715 /* Now sort the blocks */
717 pr_sortblock(debug,"Before sorting",ncons,sb);
718 fprintf(debug,"Going to sort constraints\n");
721 qsort(sb,ncons,(size_t)sizeof(*sb),pcomp);
724 pr_sortblock(debug,"After sorting",ncons,sb);
727 iatom=idef->il[F_CONSTR].iatoms;
728 for(i=0; (i<ncons); i++,iatom+=3)
730 iatom[m]=sb[i].iatom[m];
732 pr_idef(fplog,0,"After Sort",idef);
736 snew(constr->sblock,constr->nblocks+1);
738 for(i=0; (i<ncons); i++) {
739 if (sb[i].blocknr != bnr) {
741 constr->sblock[j++]=3*i;
745 constr->sblock[j++] = 3*ncons;
747 if (j != (constr->nblocks+1)) {
748 fprintf(stderr,"bstart: %d\n",bstart);
749 fprintf(stderr,"j: %d, nblocks: %d, ncons: %d\n",
750 j,constr->nblocks,ncons);
751 for(i=0; (i<ncons); i++)
752 fprintf(stderr,"i: %5d sb[i].blocknr: %5u\n",i,sb[i].blocknr);
753 for(j=0; (j<=constr->nblocks); j++)
754 fprintf(stderr,"sblock[%3d]=%5d\n",j,(int)constr->sblock[j]);
755 gmx_fatal(FARGS,"DEATH HORROR: "
756 "sblocks does not match idef->il[F_CONSTR]");
762 static void make_shake_sblock_dd(struct gmx_constr *constr,
763 t_ilist *ilcon,t_block *cgs,
769 if (dd->ncg_home+1 > constr->sblock_nalloc) {
770 constr->sblock_nalloc = over_alloc_dd(dd->ncg_home+1);
771 srenew(constr->sblock,constr->sblock_nalloc);
775 iatom = ilcon->iatoms;
778 for(c=0; c<ncons; c++) {
779 if (c == 0 || iatom[1] >= cgs->index[cg+1]) {
780 constr->sblock[constr->nblocks++] = 3*c;
781 while (iatom[1] >= cgs->index[cg+1])
786 constr->sblock[constr->nblocks] = 3*ncons;
789 t_blocka make_at2con(int start,int natoms,
790 t_ilist *ilist,t_iparams *iparams,
791 gmx_bool bDynamics,int *nflexiblecons)
793 int *count,ncon,con,con_tot,nflexcon,ftype,i,a;
800 for(ftype=F_CONSTR; ftype<=F_CONSTRNC; ftype++) {
801 ncon = ilist[ftype].nr/3;
802 ia = ilist[ftype].iatoms;
803 for(con=0; con<ncon; con++) {
804 bFlexCon = (iparams[ia[0]].constr.dA == 0 &&
805 iparams[ia[0]].constr.dB == 0);
808 if (bDynamics || !bFlexCon) {
817 *nflexiblecons = nflexcon;
820 at2con.nalloc_index = at2con.nr+1;
821 snew(at2con.index,at2con.nalloc_index);
823 for(a=0; a<natoms; a++) {
824 at2con.index[a+1] = at2con.index[a] + count[a];
827 at2con.nra = at2con.index[natoms];
828 at2con.nalloc_a = at2con.nra;
829 snew(at2con.a,at2con.nalloc_a);
831 /* The F_CONSTRNC constraints have constraint numbers
832 * that continue after the last F_CONSTR constraint.
835 for(ftype=F_CONSTR; ftype<=F_CONSTRNC; ftype++) {
836 ncon = ilist[ftype].nr/3;
837 ia = ilist[ftype].iatoms;
838 for(con=0; con<ncon; con++) {
839 bFlexCon = (iparams[ia[0]].constr.dA == 0 &&
840 iparams[ia[0]].constr.dB == 0);
841 if (bDynamics || !bFlexCon) {
844 at2con.a[at2con.index[a]+count[a]++] = con_tot;
857 static int *make_at2settle(int natoms,const t_ilist *ilist)
863 /* Set all to no settle */
864 for(a=0; a<natoms; a++)
869 stride = 1 + NRAL(F_SETTLE);
871 for(s=0; s<ilist->nr; s+=stride)
873 at2s[ilist->iatoms[s+1]] = s/stride;
874 at2s[ilist->iatoms[s+2]] = s/stride;
875 at2s[ilist->iatoms[s+3]] = s/stride;
881 void set_constraints(struct gmx_constr *constr,
882 gmx_localtop_t *top,t_inputrec *ir,
883 t_mdatoms *md,t_commrec *cr)
892 if (constr->ncon_tot > 0)
894 /* We are using the local topology,
895 * so there are only F_CONSTR constraints.
897 ncons = idef->il[F_CONSTR].nr/3;
899 /* With DD we might also need to call LINCS with ncons=0 for
900 * communicating coordinates to other nodes that do have constraints.
902 if (ir->eConstrAlg == econtLINCS)
904 set_lincs(idef,md,EI_DYNAMICS(ir->eI),cr,constr->lincsd);
906 if (ir->eConstrAlg == econtSHAKE)
910 make_shake_sblock_dd(constr,&idef->il[F_CONSTR],&top->cgs,cr->dd);
914 make_shake_sblock_pd(constr,idef,md);
916 if (ncons > constr->lagr_nalloc)
918 constr->lagr_nalloc = over_alloc_dd(ncons);
919 srenew(constr->lagr,constr->lagr_nalloc);
924 if (idef->il[F_SETTLE].nr > 0 && constr->settled == NULL)
926 settle = &idef->il[F_SETTLE];
927 iO = settle->iatoms[1];
928 iH = settle->iatoms[2];
930 settle_init(md->massT[iO],md->massT[iH],
931 md->invmass[iO],md->invmass[iH],
932 idef->iparams[settle->iatoms[0]].settle.doh,
933 idef->iparams[settle->iatoms[0]].settle.dhh);
936 /* Make a selection of the local atoms for essential dynamics */
937 if (constr->ed && cr->dd)
939 dd_make_local_ed_indices(cr->dd,constr->ed);
943 static void constr_recur(t_blocka *at2con,
944 t_ilist *ilist,t_iparams *iparams,gmx_bool bTopB,
945 int at,int depth,int nc,int *path,
946 real r0,real r1,real *r2max,
958 ncon1 = ilist[F_CONSTR].nr/3;
959 ia1 = ilist[F_CONSTR].iatoms;
960 ia2 = ilist[F_CONSTRNC].iatoms;
962 /* Loop over all constraints connected to this atom */
963 for(c=at2con->index[at]; c<at2con->index[at+1]; c++) {
965 /* Do not walk over already used constraints */
967 for(a1=0; a1<depth; a1++) {
972 ia = constr_iatomptr(ncon1,ia1,ia2,con);
973 /* Flexible constraints currently have length 0, which is incorrect */
975 len = iparams[ia[0]].constr.dA;
977 len = iparams[ia[0]].constr.dB;
978 /* In the worst case the bond directions alternate */
986 /* Assume angles of 120 degrees between all bonds */
987 if (rn0*rn0 + rn1*rn1 + rn0*rn1 > *r2max) {
988 *r2max = rn0*rn0 + rn1*rn1 + r0*rn1;
990 fprintf(debug,"Found longer constraint distance: r0 %5.3f r1 %5.3f rmax %5.3f\n", rn0,rn1,sqrt(*r2max));
991 for(a1=0; a1<depth; a1++)
992 fprintf(debug," %d %5.3f",
994 iparams[constr_iatomptr(ncon1,ia1,ia2,con)[0]].constr.dA);
995 fprintf(debug," %d %5.3f\n",con,len);
998 /* Limit the number of recursions to 1000*nc,
999 * so a call does not take more than a second,
1000 * even for highly connected systems.
1002 if (depth + 1 < nc && *count < 1000*nc) {
1009 constr_recur(at2con,ilist,iparams,
1010 bTopB,a1,depth+1,nc,path,rn0,rn1,r2max,count);
1017 static real constr_r_max_moltype(FILE *fplog,
1018 gmx_moltype_t *molt,t_iparams *iparams,
1021 int natoms,nflexcon,*path,at,count;
1024 real r0,r1,r2maxA,r2maxB,rmax,lam0,lam1;
1026 if (molt->ilist[F_CONSTR].nr == 0 &&
1027 molt->ilist[F_CONSTRNC].nr == 0) {
1031 natoms = molt->atoms.nr;
1033 at2con = make_at2con(0,natoms,molt->ilist,iparams,
1034 EI_DYNAMICS(ir->eI),&nflexcon);
1035 snew(path,1+ir->nProjOrder);
1036 for(at=0; at<1+ir->nProjOrder; at++)
1040 for(at=0; at<natoms; at++) {
1045 constr_recur(&at2con,molt->ilist,iparams,
1046 FALSE,at,0,1+ir->nProjOrder,path,r0,r1,&r2maxA,&count);
1048 if (ir->efep == efepNO) {
1049 rmax = sqrt(r2maxA);
1052 for(at=0; at<natoms; at++) {
1056 constr_recur(&at2con,molt->ilist,iparams,
1057 TRUE,at,0,1+ir->nProjOrder,path,r0,r1,&r2maxB,&count);
1059 lam0 = ir->fepvals->init_lambda;
1060 if (EI_DYNAMICS(ir->eI))
1061 lam0 += ir->init_step*ir->fepvals->delta_lambda;
1062 rmax = (1 - lam0)*sqrt(r2maxA) + lam0*sqrt(r2maxB);
1063 if (EI_DYNAMICS(ir->eI)) {
1064 lam1 = ir->fepvals->init_lambda + (ir->init_step + ir->nsteps)*ir->fepvals->delta_lambda;
1065 rmax = max(rmax,(1 - lam1)*sqrt(r2maxA) + lam1*sqrt(r2maxB));
1069 done_blocka(&at2con);
1075 real constr_r_max(FILE *fplog,gmx_mtop_t *mtop,t_inputrec *ir)
1081 for(mt=0; mt<mtop->nmoltype; mt++) {
1083 constr_r_max_moltype(fplog,&mtop->moltype[mt],
1084 mtop->ffparams.iparams,ir));
1088 fprintf(fplog,"Maximum distance for %d constraints, at 120 deg. angles, all-trans: %.3f nm\n",1+ir->nProjOrder,rmax);
1093 gmx_constr_t init_constraints(FILE *fplog,
1094 gmx_mtop_t *mtop,t_inputrec *ir,
1095 gmx_edsam_t ed,t_state *state,
1098 int ncon,nset,nmol,settle_type,i,natoms,mt,nflexcon;
1099 struct gmx_constr *constr;
1102 gmx_mtop_ilistloop_t iloop;
1105 gmx_mtop_ftype_count(mtop,F_CONSTR) +
1106 gmx_mtop_ftype_count(mtop,F_CONSTRNC);
1107 nset = gmx_mtop_ftype_count(mtop,F_SETTLE);
1109 if (ncon+nset == 0 && ir->ePull != epullCONSTRAINT && ed == NULL)
1116 constr->ncon_tot = ncon;
1117 constr->nflexcon = 0;
1120 constr->n_at2con_mt = mtop->nmoltype;
1121 snew(constr->at2con_mt,constr->n_at2con_mt);
1122 for(mt=0; mt<mtop->nmoltype; mt++)
1124 constr->at2con_mt[mt] = make_at2con(0,mtop->moltype[mt].atoms.nr,
1125 mtop->moltype[mt].ilist,
1126 mtop->ffparams.iparams,
1127 EI_DYNAMICS(ir->eI),&nflexcon);
1128 for(i=0; i<mtop->nmolblock; i++)
1130 if (mtop->molblock[i].type == mt)
1132 constr->nflexcon += mtop->molblock[i].nmol*nflexcon;
1137 if (constr->nflexcon > 0)
1141 fprintf(fplog,"There are %d flexible constraints\n",
1143 if (ir->fc_stepsize == 0)
1146 "WARNING: step size for flexible constraining = 0\n"
1147 " All flexible constraints will be rigid.\n"
1148 " Will try to keep all flexible constraints at their original length,\n"
1149 " but the lengths may exhibit some drift.\n\n");
1150 constr->nflexcon = 0;
1153 if (constr->nflexcon > 0)
1155 please_cite(fplog,"Hess2002");
1159 if (ir->eConstrAlg == econtLINCS)
1161 constr->lincsd = init_lincs(fplog,mtop,
1162 constr->nflexcon,constr->at2con_mt,
1163 DOMAINDECOMP(cr) && cr->dd->bInterCGcons,
1164 ir->nLincsIter,ir->nProjOrder);
1167 if (ir->eConstrAlg == econtSHAKE) {
1168 if (DOMAINDECOMP(cr) && cr->dd->bInterCGcons)
1170 gmx_fatal(FARGS,"SHAKE is not supported with domain decomposition and constraint that cross charge group boundaries, use LINCS");
1172 if (constr->nflexcon)
1174 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");
1176 please_cite(fplog,"Ryckaert77a");
1179 please_cite(fplog,"Barth95a");
1182 constr->shaked = shake_init();
1187 please_cite(fplog,"Miyamoto92a");
1189 constr->bInterCGsettles = inter_charge_group_settles(mtop);
1191 /* Check that we have only one settle type */
1193 iloop = gmx_mtop_ilistloop_init(mtop);
1194 while (gmx_mtop_ilistloop_next(iloop,&ilist,&nmol))
1196 for (i=0; i<ilist[F_SETTLE].nr; i+=4)
1198 if (settle_type == -1)
1200 settle_type = ilist[F_SETTLE].iatoms[i];
1202 else if (ilist[F_SETTLE].iatoms[i] != settle_type)
1205 "The [molecules] section of your topology specifies more than one block of\n"
1206 "a [moleculetype] with a [settles] block. Only one such is allowed. If you\n"
1207 "are trying to partition your solvent into different *groups* (e.g. for\n"
1208 "freezing, T-coupling, etc.) then you are using the wrong approach. Index\n"
1209 "files specify groups. Otherwise, you may wish to change the least-used\n"
1210 "block of molecules with SETTLE constraints into 3 normal constraints.");
1215 constr->n_at2settle_mt = mtop->nmoltype;
1216 snew(constr->at2settle_mt,constr->n_at2settle_mt);
1217 for(mt=0; mt<mtop->nmoltype; mt++)
1219 constr->at2settle_mt[mt] =
1220 make_at2settle(mtop->moltype[mt].atoms.nr,
1221 &mtop->moltype[mt].ilist[F_SETTLE]);
1225 constr->maxwarn = 999;
1226 env = getenv("GMX_MAXCONSTRWARN");
1229 constr->maxwarn = 0;
1230 sscanf(env,"%d",&constr->maxwarn);
1234 "Setting the maximum number of constraint warnings to %d\n",
1240 "Setting the maximum number of constraint warnings to %d\n",
1244 if (constr->maxwarn < 0 && fplog)
1246 fprintf(fplog,"maxwarn < 0, will not stop on constraint errors\n");
1248 constr->warncount_lincs = 0;
1249 constr->warncount_settle = 0;
1251 /* Initialize the essential dynamics sampling.
1252 * Put the pointer to the ED struct in constr */
1254 if (ed != NULL || state->edsamstate.nED > 0)
1256 init_edsam(mtop,ir,cr,ed,state->x,state->box,&state->edsamstate);
1259 constr->warn_mtop = mtop;
1264 const t_blocka *atom2constraints_moltype(gmx_constr_t constr)
1266 return constr->at2con_mt;
1269 const int **atom2settle_moltype(gmx_constr_t constr)
1271 return (const int **)constr->at2settle_mt;
1275 gmx_bool inter_charge_group_constraints(const gmx_mtop_t *mtop)
1277 const gmx_moltype_t *molt;
1281 int nat,*at2cg,cg,a,ftype,i;
1285 for(mb=0; mb<mtop->nmolblock && !bInterCG; mb++)
1287 molt = &mtop->moltype[mtop->molblock[mb].type];
1289 if (molt->ilist[F_CONSTR].nr > 0 ||
1290 molt->ilist[F_CONSTRNC].nr > 0 ||
1291 molt->ilist[F_SETTLE].nr > 0)
1294 snew(at2cg,molt->atoms.nr);
1295 for(cg=0; cg<cgs->nr; cg++)
1297 for(a=cgs->index[cg]; a<cgs->index[cg+1]; a++)
1301 for(ftype=F_CONSTR; ftype<=F_CONSTRNC; ftype++)
1303 il = &molt->ilist[ftype];
1304 for(i=0; i<il->nr && !bInterCG; i+=1+NRAL(ftype))
1306 if (at2cg[il->iatoms[i+1]] != at2cg[il->iatoms[i+2]])
1320 gmx_bool inter_charge_group_settles(const gmx_mtop_t *mtop)
1322 const gmx_moltype_t *molt;
1326 int nat,*at2cg,cg,a,ftype,i;
1330 for(mb=0; mb<mtop->nmolblock && !bInterCG; mb++)
1332 molt = &mtop->moltype[mtop->molblock[mb].type];
1334 if (molt->ilist[F_SETTLE].nr > 0)
1337 snew(at2cg,molt->atoms.nr);
1338 for(cg=0; cg<cgs->nr; cg++)
1340 for(a=cgs->index[cg]; a<cgs->index[cg+1]; a++)
1344 for(ftype=F_SETTLE; ftype<=F_SETTLE; ftype++)
1346 il = &molt->ilist[ftype];
1347 for(i=0; i<il->nr && !bInterCG; i+=1+NRAL(F_SETTLE))
1349 if (at2cg[il->iatoms[i+1]] != at2cg[il->iatoms[i+2]] ||
1350 at2cg[il->iatoms[i+1]] != at2cg[il->iatoms[i+3]])
1364 /* helper functions for andersen temperature control, because the
1365 * gmx_constr construct is only defined in constr.c. Return the list
1366 * of blocks (get_sblock) and the number of blocks (get_nblocks). */
1368 extern int *get_sblock(struct gmx_constr *constr)
1370 return constr->sblock;
1373 extern int get_nblocks(struct gmx_constr *constr)
1375 return constr->nblocks;