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 *rmdr_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;
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->rmdr_th == NULL)
371 snew(constr->rmdr_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,rmdr,
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,rmdr,constr->maxwarn>=0,econq,&vetavar);
438 bOK = bshakef(fplog,constr->shaked,
439 homenr,md->invmass,constr->nblocks,constr->sblock,
440 idef,ir,x,min_proj,nrnb,
441 constr->lagr,lambda,dvdlambda,
442 invdt,NULL,vir!=NULL,rmdr,constr->maxwarn>=0,econq,&vetavar);
445 gmx_fatal(FARGS,"Internal error, SHAKE called for constraining something else than coordinates");
449 if (!bOK && constr->maxwarn >= 0)
453 fprintf(fplog,"Constraint error in algorithm %s at step %s\n",
454 econstr_names[econtSHAKE],gmx_step_str(step,buf));
462 int calcvir_atom_end;
466 calcvir_atom_end = 0;
470 calcvir_atom_end = md->start + md->homenr;
476 #pragma omp parallel for num_threads(nth) schedule(static)
477 for(th=0; th<nth; th++)
483 clear_mat(constr->rmdr_th[th]);
486 start_th = (nsettle* th )/nth;
487 end_th = (nsettle*(th+1))/nth;
488 if (start_th >= 0 && end_th - start_th > 0)
490 csettle(constr->settled,
492 settle->iatoms+start_th*(1+NRAL(F_SETTLE)),
495 invdt,v?v[0]:NULL,calcvir_atom_end,
496 th == 0 ? rmdr : constr->rmdr_th[th],
497 th == 0 ? &settle_error : &constr->settle_error[th],
501 inc_nrnb(nrnb,eNR_SETTLE,nsettle);
504 inc_nrnb(nrnb,eNR_CONSTR_V,nsettle*3);
508 inc_nrnb(nrnb,eNR_CONSTR_VIR,nsettle*3);
514 case econqForceDispl:
515 #pragma omp parallel for num_threads(nth) schedule(static)
516 for(th=0; th<nth; th++)
522 clear_mat(constr->rmdr_th[th]);
525 start_th = (nsettle* th )/nth;
526 end_th = (nsettle*(th+1))/nth;
528 if (start_th >= 0 && end_th - start_th > 0)
530 settle_proj(fplog,constr->settled,econq,
532 settle->iatoms+start_th*(1+NRAL(F_SETTLE)),
535 xprime,min_proj,calcvir_atom_end,
536 th == 0 ? rmdr : constr->rmdr_th[th],
540 /* This is an overestimate */
541 inc_nrnb(nrnb,eNR_SETTLE,nsettle);
543 case econqDeriv_FlexCon:
544 /* Nothing to do, since the are no flexible constraints in settles */
547 gmx_incons("Unknown constraint quantity for settle");
553 /* Combine virial and error info of the other threads */
556 m_add(rmdr,constr->rmdr_th[i],rmdr);
557 settle_error = constr->settle_error[i];
560 if (econq == econqCoord && settle_error >= 0)
563 if (constr->maxwarn >= 0)
567 "\nstep " gmx_large_int_pfmt ": Water molecule starting at atom %d can not be "
568 "settled.\nCheck for bad contacts and/or reduce the timestep if appropriate.\n",
569 step,ddglatnr(cr->dd,settle->iatoms[settle_error*(1+NRAL(F_SETTLE))+1]));
572 fprintf(fplog,"%s",buf);
574 fprintf(stderr,"%s",buf);
575 constr->warncount_settle++;
576 if (constr->warncount_settle > constr->maxwarn)
578 too_many_constraint_warnings(-1,constr->warncount_settle);
585 free_vetavars(&vetavar);
592 vir_fac = 0.5/(ir->delta_t*ir->delta_t);
595 vir_fac = 0.5/ir->delta_t;
598 case econqForceDispl:
603 gmx_incons("Unsupported constraint quantity for virial");
608 vir_fac *= 2; /* only constraining over half the distance here */
614 (*vir)[i][j] = vir_fac*rmdr[i][j];
621 dump_confs(fplog,step,constr->warn_mtop,start,homenr,cr,x,xprime,box);
624 if (econq == econqCoord)
626 if (ir->ePull == epullCONSTRAINT)
628 if (EI_DYNAMICS(ir->eI))
630 t = ir->init_t + (step + delta_step)*ir->delta_t;
636 set_pbc(&pbc,ir->ePBC,box);
637 pull_constraint(ir->pull,md,&pbc,cr,ir->delta_t,t,x,xprime,v,*vir);
639 if (constr->ed && delta_step > 0)
641 /* apply the essential dynamcs constraints here */
642 do_edsam(ir,step,md,cr,xprime,v,box,constr->ed);
649 real *constr_rmsd_data(struct gmx_constr *constr)
652 return lincs_rmsd_data(constr->lincsd);
657 real constr_rmsd(struct gmx_constr *constr,gmx_bool bSD2)
660 return lincs_rmsd(constr->lincsd,bSD2);
665 static void make_shake_sblock_pd(struct gmx_constr *constr,
666 t_idef *idef,t_mdatoms *md)
675 /* Since we are processing the local topology,
676 * the F_CONSTRNC ilist has been concatenated to the F_CONSTR ilist.
678 ncons = idef->il[F_CONSTR].nr/3;
680 init_blocka(&sblocks);
681 gen_sblocks(NULL,md->start,md->start+md->homenr,idef,&sblocks,FALSE);
684 bstart=(idef->nodeid > 0) ? blocks->multinr[idef->nodeid-1] : 0;
685 nblocks=blocks->multinr[idef->nodeid] - bstart;
688 constr->nblocks = sblocks.nr;
690 fprintf(debug,"ncons: %d, bstart: %d, nblocks: %d\n",
691 ncons,bstart,constr->nblocks);
693 /* Calculate block number for each atom */
694 inv_sblock = make_invblocka(&sblocks,md->nr);
696 done_blocka(&sblocks);
698 /* Store the block number in temp array and
699 * sort the constraints in order of the sblock number
700 * and the atom numbers, really sorting a segment of the array!
703 pr_idef(fplog,0,"Before Sort",idef);
705 iatom=idef->il[F_CONSTR].iatoms;
707 for(i=0; (i<ncons); i++,iatom+=3) {
709 sb[i].iatom[m] = iatom[m];
710 sb[i].blocknr = inv_sblock[iatom[1]];
713 /* Now sort the blocks */
715 pr_sortblock(debug,"Before sorting",ncons,sb);
716 fprintf(debug,"Going to sort constraints\n");
719 qsort(sb,ncons,(size_t)sizeof(*sb),pcomp);
722 pr_sortblock(debug,"After sorting",ncons,sb);
725 iatom=idef->il[F_CONSTR].iatoms;
726 for(i=0; (i<ncons); i++,iatom+=3)
728 iatom[m]=sb[i].iatom[m];
730 pr_idef(fplog,0,"After Sort",idef);
734 snew(constr->sblock,constr->nblocks+1);
736 for(i=0; (i<ncons); i++) {
737 if (sb[i].blocknr != bnr) {
739 constr->sblock[j++]=3*i;
743 constr->sblock[j++] = 3*ncons;
745 if (j != (constr->nblocks+1)) {
746 fprintf(stderr,"bstart: %d\n",bstart);
747 fprintf(stderr,"j: %d, nblocks: %d, ncons: %d\n",
748 j,constr->nblocks,ncons);
749 for(i=0; (i<ncons); i++)
750 fprintf(stderr,"i: %5d sb[i].blocknr: %5u\n",i,sb[i].blocknr);
751 for(j=0; (j<=constr->nblocks); j++)
752 fprintf(stderr,"sblock[%3d]=%5d\n",j,(int)constr->sblock[j]);
753 gmx_fatal(FARGS,"DEATH HORROR: "
754 "sblocks does not match idef->il[F_CONSTR]");
760 static void make_shake_sblock_dd(struct gmx_constr *constr,
761 t_ilist *ilcon,t_block *cgs,
767 if (dd->ncg_home+1 > constr->sblock_nalloc) {
768 constr->sblock_nalloc = over_alloc_dd(dd->ncg_home+1);
769 srenew(constr->sblock,constr->sblock_nalloc);
773 iatom = ilcon->iatoms;
776 for(c=0; c<ncons; c++) {
777 if (c == 0 || iatom[1] >= cgs->index[cg+1]) {
778 constr->sblock[constr->nblocks++] = 3*c;
779 while (iatom[1] >= cgs->index[cg+1])
784 constr->sblock[constr->nblocks] = 3*ncons;
787 t_blocka make_at2con(int start,int natoms,
788 t_ilist *ilist,t_iparams *iparams,
789 gmx_bool bDynamics,int *nflexiblecons)
791 int *count,ncon,con,con_tot,nflexcon,ftype,i,a;
798 for(ftype=F_CONSTR; ftype<=F_CONSTRNC; ftype++) {
799 ncon = ilist[ftype].nr/3;
800 ia = ilist[ftype].iatoms;
801 for(con=0; con<ncon; con++) {
802 bFlexCon = (iparams[ia[0]].constr.dA == 0 &&
803 iparams[ia[0]].constr.dB == 0);
806 if (bDynamics || !bFlexCon) {
815 *nflexiblecons = nflexcon;
818 at2con.nalloc_index = at2con.nr+1;
819 snew(at2con.index,at2con.nalloc_index);
821 for(a=0; a<natoms; a++) {
822 at2con.index[a+1] = at2con.index[a] + count[a];
825 at2con.nra = at2con.index[natoms];
826 at2con.nalloc_a = at2con.nra;
827 snew(at2con.a,at2con.nalloc_a);
829 /* The F_CONSTRNC constraints have constraint numbers
830 * that continue after the last F_CONSTR constraint.
833 for(ftype=F_CONSTR; ftype<=F_CONSTRNC; ftype++) {
834 ncon = ilist[ftype].nr/3;
835 ia = ilist[ftype].iatoms;
836 for(con=0; con<ncon; con++) {
837 bFlexCon = (iparams[ia[0]].constr.dA == 0 &&
838 iparams[ia[0]].constr.dB == 0);
839 if (bDynamics || !bFlexCon) {
842 at2con.a[at2con.index[a]+count[a]++] = con_tot;
855 static int *make_at2settle(int natoms,const t_ilist *ilist)
861 /* Set all to no settle */
862 for(a=0; a<natoms; a++)
867 stride = 1 + NRAL(F_SETTLE);
869 for(s=0; s<ilist->nr; s+=stride)
871 at2s[ilist->iatoms[s+1]] = s/stride;
872 at2s[ilist->iatoms[s+2]] = s/stride;
873 at2s[ilist->iatoms[s+3]] = s/stride;
879 void set_constraints(struct gmx_constr *constr,
880 gmx_localtop_t *top,t_inputrec *ir,
881 t_mdatoms *md,t_commrec *cr)
890 if (constr->ncon_tot > 0)
892 /* We are using the local topology,
893 * so there are only F_CONSTR constraints.
895 ncons = idef->il[F_CONSTR].nr/3;
897 /* With DD we might also need to call LINCS with ncons=0 for
898 * communicating coordinates to other nodes that do have constraints.
900 if (ir->eConstrAlg == econtLINCS)
902 set_lincs(idef,md,EI_DYNAMICS(ir->eI),cr,constr->lincsd);
904 if (ir->eConstrAlg == econtSHAKE)
908 make_shake_sblock_dd(constr,&idef->il[F_CONSTR],&top->cgs,cr->dd);
912 make_shake_sblock_pd(constr,idef,md);
914 if (ncons > constr->lagr_nalloc)
916 constr->lagr_nalloc = over_alloc_dd(ncons);
917 srenew(constr->lagr,constr->lagr_nalloc);
922 if (idef->il[F_SETTLE].nr > 0 && constr->settled == NULL)
924 settle = &idef->il[F_SETTLE];
925 iO = settle->iatoms[1];
926 iH = settle->iatoms[2];
928 settle_init(md->massT[iO],md->massT[iH],
929 md->invmass[iO],md->invmass[iH],
930 idef->iparams[settle->iatoms[0]].settle.doh,
931 idef->iparams[settle->iatoms[0]].settle.dhh);
934 /* Make a selection of the local atoms for essential dynamics */
935 if (constr->ed && cr->dd)
937 dd_make_local_ed_indices(cr->dd,constr->ed);
941 static void constr_recur(t_blocka *at2con,
942 t_ilist *ilist,t_iparams *iparams,gmx_bool bTopB,
943 int at,int depth,int nc,int *path,
944 real r0,real r1,real *r2max,
956 ncon1 = ilist[F_CONSTR].nr/3;
957 ia1 = ilist[F_CONSTR].iatoms;
958 ia2 = ilist[F_CONSTRNC].iatoms;
960 /* Loop over all constraints connected to this atom */
961 for(c=at2con->index[at]; c<at2con->index[at+1]; c++) {
963 /* Do not walk over already used constraints */
965 for(a1=0; a1<depth; a1++) {
970 ia = constr_iatomptr(ncon1,ia1,ia2,con);
971 /* Flexible constraints currently have length 0, which is incorrect */
973 len = iparams[ia[0]].constr.dA;
975 len = iparams[ia[0]].constr.dB;
976 /* In the worst case the bond directions alternate */
984 /* Assume angles of 120 degrees between all bonds */
985 if (rn0*rn0 + rn1*rn1 + rn0*rn1 > *r2max) {
986 *r2max = rn0*rn0 + rn1*rn1 + r0*rn1;
988 fprintf(debug,"Found longer constraint distance: r0 %5.3f r1 %5.3f rmax %5.3f\n", rn0,rn1,sqrt(*r2max));
989 for(a1=0; a1<depth; a1++)
990 fprintf(debug," %d %5.3f",
992 iparams[constr_iatomptr(ncon1,ia1,ia2,con)[0]].constr.dA);
993 fprintf(debug," %d %5.3f\n",con,len);
996 /* Limit the number of recursions to 1000*nc,
997 * so a call does not take more than a second,
998 * even for highly connected systems.
1000 if (depth + 1 < nc && *count < 1000*nc) {
1007 constr_recur(at2con,ilist,iparams,
1008 bTopB,a1,depth+1,nc,path,rn0,rn1,r2max,count);
1015 static real constr_r_max_moltype(FILE *fplog,
1016 gmx_moltype_t *molt,t_iparams *iparams,
1019 int natoms,nflexcon,*path,at,count;
1022 real r0,r1,r2maxA,r2maxB,rmax,lam0,lam1;
1024 if (molt->ilist[F_CONSTR].nr == 0 &&
1025 molt->ilist[F_CONSTRNC].nr == 0) {
1029 natoms = molt->atoms.nr;
1031 at2con = make_at2con(0,natoms,molt->ilist,iparams,
1032 EI_DYNAMICS(ir->eI),&nflexcon);
1033 snew(path,1+ir->nProjOrder);
1034 for(at=0; at<1+ir->nProjOrder; at++)
1038 for(at=0; at<natoms; at++) {
1043 constr_recur(&at2con,molt->ilist,iparams,
1044 FALSE,at,0,1+ir->nProjOrder,path,r0,r1,&r2maxA,&count);
1046 if (ir->efep == efepNO) {
1047 rmax = sqrt(r2maxA);
1050 for(at=0; at<natoms; at++) {
1054 constr_recur(&at2con,molt->ilist,iparams,
1055 TRUE,at,0,1+ir->nProjOrder,path,r0,r1,&r2maxB,&count);
1057 lam0 = ir->fepvals->init_lambda;
1058 if (EI_DYNAMICS(ir->eI))
1059 lam0 += ir->init_step*ir->fepvals->delta_lambda;
1060 rmax = (1 - lam0)*sqrt(r2maxA) + lam0*sqrt(r2maxB);
1061 if (EI_DYNAMICS(ir->eI)) {
1062 lam1 = ir->fepvals->init_lambda + (ir->init_step + ir->nsteps)*ir->fepvals->delta_lambda;
1063 rmax = max(rmax,(1 - lam1)*sqrt(r2maxA) + lam1*sqrt(r2maxB));
1067 done_blocka(&at2con);
1073 real constr_r_max(FILE *fplog,gmx_mtop_t *mtop,t_inputrec *ir)
1079 for(mt=0; mt<mtop->nmoltype; mt++) {
1081 constr_r_max_moltype(fplog,&mtop->moltype[mt],
1082 mtop->ffparams.iparams,ir));
1086 fprintf(fplog,"Maximum distance for %d constraints, at 120 deg. angles, all-trans: %.3f nm\n",1+ir->nProjOrder,rmax);
1091 gmx_constr_t init_constraints(FILE *fplog,
1092 gmx_mtop_t *mtop,t_inputrec *ir,
1093 gmx_edsam_t ed,t_state *state,
1096 int ncon,nset,nmol,settle_type,i,natoms,mt,nflexcon;
1097 struct gmx_constr *constr;
1100 gmx_mtop_ilistloop_t iloop;
1103 gmx_mtop_ftype_count(mtop,F_CONSTR) +
1104 gmx_mtop_ftype_count(mtop,F_CONSTRNC);
1105 nset = gmx_mtop_ftype_count(mtop,F_SETTLE);
1107 if (ncon+nset == 0 && ir->ePull != epullCONSTRAINT && ed == NULL)
1114 constr->ncon_tot = ncon;
1115 constr->nflexcon = 0;
1118 constr->n_at2con_mt = mtop->nmoltype;
1119 snew(constr->at2con_mt,constr->n_at2con_mt);
1120 for(mt=0; mt<mtop->nmoltype; mt++)
1122 constr->at2con_mt[mt] = make_at2con(0,mtop->moltype[mt].atoms.nr,
1123 mtop->moltype[mt].ilist,
1124 mtop->ffparams.iparams,
1125 EI_DYNAMICS(ir->eI),&nflexcon);
1126 for(i=0; i<mtop->nmolblock; i++)
1128 if (mtop->molblock[i].type == mt)
1130 constr->nflexcon += mtop->molblock[i].nmol*nflexcon;
1135 if (constr->nflexcon > 0)
1139 fprintf(fplog,"There are %d flexible constraints\n",
1141 if (ir->fc_stepsize == 0)
1144 "WARNING: step size for flexible constraining = 0\n"
1145 " All flexible constraints will be rigid.\n"
1146 " Will try to keep all flexible constraints at their original length,\n"
1147 " but the lengths may exhibit some drift.\n\n");
1148 constr->nflexcon = 0;
1151 if (constr->nflexcon > 0)
1153 please_cite(fplog,"Hess2002");
1157 if (ir->eConstrAlg == econtLINCS)
1159 constr->lincsd = init_lincs(fplog,mtop,
1160 constr->nflexcon,constr->at2con_mt,
1161 DOMAINDECOMP(cr) && cr->dd->bInterCGcons,
1162 ir->nLincsIter,ir->nProjOrder);
1165 if (ir->eConstrAlg == econtSHAKE) {
1166 if (DOMAINDECOMP(cr) && cr->dd->bInterCGcons)
1168 gmx_fatal(FARGS,"SHAKE is not supported with domain decomposition and constraint that cross charge group boundaries, use LINCS");
1170 if (constr->nflexcon)
1172 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");
1174 please_cite(fplog,"Ryckaert77a");
1177 please_cite(fplog,"Barth95a");
1180 constr->shaked = shake_init();
1185 please_cite(fplog,"Miyamoto92a");
1187 constr->bInterCGsettles = inter_charge_group_settles(mtop);
1189 /* Check that we have only one settle type */
1191 iloop = gmx_mtop_ilistloop_init(mtop);
1192 while (gmx_mtop_ilistloop_next(iloop,&ilist,&nmol))
1194 for (i=0; i<ilist[F_SETTLE].nr; i+=4)
1196 if (settle_type == -1)
1198 settle_type = ilist[F_SETTLE].iatoms[i];
1200 else if (ilist[F_SETTLE].iatoms[i] != settle_type)
1203 "The [molecules] section of your topology specifies more than one block of\n"
1204 "a [moleculetype] with a [settles] block. Only one such is allowed. If you\n"
1205 "are trying to partition your solvent into different *groups* (e.g. for\n"
1206 "freezing, T-coupling, etc.) then you are using the wrong approach. Index\n"
1207 "files specify groups. Otherwise, you may wish to change the least-used\n"
1208 "block of molecules with SETTLE constraints into 3 normal constraints.");
1213 constr->n_at2settle_mt = mtop->nmoltype;
1214 snew(constr->at2settle_mt,constr->n_at2settle_mt);
1215 for(mt=0; mt<mtop->nmoltype; mt++)
1217 constr->at2settle_mt[mt] =
1218 make_at2settle(mtop->moltype[mt].atoms.nr,
1219 &mtop->moltype[mt].ilist[F_SETTLE]);
1223 constr->maxwarn = 999;
1224 env = getenv("GMX_MAXCONSTRWARN");
1227 constr->maxwarn = 0;
1228 sscanf(env,"%d",&constr->maxwarn);
1232 "Setting the maximum number of constraint warnings to %d\n",
1238 "Setting the maximum number of constraint warnings to %d\n",
1242 if (constr->maxwarn < 0 && fplog)
1244 fprintf(fplog,"maxwarn < 0, will not stop on constraint errors\n");
1246 constr->warncount_lincs = 0;
1247 constr->warncount_settle = 0;
1249 /* Initialize the essential dynamics sampling.
1250 * Put the pointer to the ED struct in constr */
1254 init_edsam(mtop,ir,cr,ed,state->x,state->box);
1257 constr->warn_mtop = mtop;
1262 const t_blocka *atom2constraints_moltype(gmx_constr_t constr)
1264 return constr->at2con_mt;
1267 const int **atom2settle_moltype(gmx_constr_t constr)
1269 return (const int **)constr->at2settle_mt;
1273 gmx_bool inter_charge_group_constraints(const gmx_mtop_t *mtop)
1275 const gmx_moltype_t *molt;
1279 int nat,*at2cg,cg,a,ftype,i;
1283 for(mb=0; mb<mtop->nmolblock && !bInterCG; mb++)
1285 molt = &mtop->moltype[mtop->molblock[mb].type];
1287 if (molt->ilist[F_CONSTR].nr > 0 ||
1288 molt->ilist[F_CONSTRNC].nr > 0 ||
1289 molt->ilist[F_SETTLE].nr > 0)
1292 snew(at2cg,molt->atoms.nr);
1293 for(cg=0; cg<cgs->nr; cg++)
1295 for(a=cgs->index[cg]; a<cgs->index[cg+1]; a++)
1299 for(ftype=F_CONSTR; ftype<=F_CONSTRNC; ftype++)
1301 il = &molt->ilist[ftype];
1302 for(i=0; i<il->nr && !bInterCG; i+=1+NRAL(ftype))
1304 if (at2cg[il->iatoms[i+1]] != at2cg[il->iatoms[i+2]])
1318 gmx_bool inter_charge_group_settles(const gmx_mtop_t *mtop)
1320 const gmx_moltype_t *molt;
1324 int nat,*at2cg,cg,a,ftype,i;
1328 for(mb=0; mb<mtop->nmolblock && !bInterCG; mb++)
1330 molt = &mtop->moltype[mtop->molblock[mb].type];
1332 if (molt->ilist[F_SETTLE].nr > 0)
1335 snew(at2cg,molt->atoms.nr);
1336 for(cg=0; cg<cgs->nr; cg++)
1338 for(a=cgs->index[cg]; a<cgs->index[cg+1]; a++)
1342 for(ftype=F_SETTLE; ftype<=F_SETTLE; ftype++)
1344 il = &molt->ilist[ftype];
1345 for(i=0; i<il->nr && !bInterCG; i+=1+NRAL(F_SETTLE))
1347 if (at2cg[il->iatoms[i+1]] != at2cg[il->iatoms[i+2]] ||
1348 at2cg[il->iatoms[i+1]] != at2cg[il->iatoms[i+3]])
1362 /* helper functions for andersen temperature control, because the
1363 * gmx_constr construct is only defined in constr.c. Return the list
1364 * of blocks (get_sblock) and the number of blocks (get_nblocks). */
1366 extern int *get_sblock(struct gmx_constr *constr)
1368 return constr->sblock;
1371 extern int get_nblocks(struct gmx_constr *constr)
1373 return constr->nblocks;