2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
48 #include "gromacs/domdec/domdec.h"
49 #include "gromacs/domdec/domdec_struct.h"
50 #include "gromacs/essentialdynamics/edsam.h"
51 #include "gromacs/fileio/confio.h"
52 #include "gromacs/fileio/gmxfio.h"
53 #include "gromacs/fileio/pdbio.h"
54 #include "gromacs/gmxlib/nrnb.h"
55 #include "gromacs/math/vec.h"
56 #include "gromacs/mdlib/gmx_omp_nthreads.h"
57 #include "gromacs/mdlib/mdrun.h"
58 #include "gromacs/mdlib/splitter.h"
59 #include "gromacs/mdtypes/commrec.h"
60 #include "gromacs/mdtypes/md_enums.h"
61 #include "gromacs/pbcutil/pbc.h"
62 #include "gromacs/pulling/pull.h"
63 #include "gromacs/topology/block.h"
64 #include "gromacs/topology/invblock.h"
65 #include "gromacs/topology/mtop_util.h"
66 #include "gromacs/utility/exceptions.h"
67 #include "gromacs/utility/fatalerror.h"
68 #include "gromacs/utility/pleasecite.h"
69 #include "gromacs/utility/smalloc.h"
70 #include "gromacs/utility/txtdump.h"
72 typedef struct gmx_constr
{
73 int ncon_tot
; /* The total number of constraints */
74 int nflexcon
; /* The number of flexible constraints */
75 int n_at2con_mt
; /* The size of at2con = #moltypes */
76 t_blocka
*at2con_mt
; /* A list of atoms to constraints */
77 int n_at2settle_mt
; /* The size of at2settle = #moltypes */
78 int **at2settle_mt
; /* A list of atoms to settles */
79 gmx_bool bInterCGsettles
;
80 gmx_lincsdata_t lincsd
; /* LINCS data */
81 gmx_shakedata_t shaked
; /* SHAKE data */
82 gmx_settledata_t settled
; /* SETTLE data */
83 int nblocks
; /* The number of SHAKE blocks */
84 int *sblock
; /* The SHAKE blocks */
85 int sblock_nalloc
; /* The allocation size of sblock */
86 real
*lagr
; /* -2 times the Lagrange multipliers for SHAKE */
87 int lagr_nalloc
; /* The allocation size of lagr */
88 int maxwarn
; /* The maximum number of warnings */
91 gmx_edsam_t ed
; /* The essential dynamics data */
93 tensor
*vir_r_m_dr_th
; /* Thread local working data */
94 int *settle_error
; /* Thread local working data */
96 const gmx_mtop_t
*warn_mtop
; /* Only used for printing warnings */
104 static int pcomp(const void *p1
, const void *p2
)
107 int min1
, min2
, max1
, max2
;
108 t_sortblock
*a1
= (t_sortblock
*)p1
;
109 t_sortblock
*a2
= (t_sortblock
*)p2
;
111 db
= a1
->blocknr
-a2
->blocknr
;
118 min1
= std::min(a1
->iatom
[1], a1
->iatom
[2]);
119 max1
= std::max(a1
->iatom
[1], a1
->iatom
[2]);
120 min2
= std::min(a2
->iatom
[1], a2
->iatom
[2]);
121 max2
= std::max(a2
->iatom
[1], a2
->iatom
[2]);
133 int n_flexible_constraints(struct gmx_constr
*constr
)
139 nflexcon
= constr
->nflexcon
;
149 static void clear_constraint_quantity_nonlocal(gmx_domdec_t
*dd
, rvec
*q
)
151 int nonlocal_at_start
, nonlocal_at_end
, at
;
153 dd_get_constraint_range(dd
, &nonlocal_at_start
, &nonlocal_at_end
);
155 for (at
= nonlocal_at_start
; at
< nonlocal_at_end
; at
++)
161 void too_many_constraint_warnings(int eConstrAlg
, int warncount
)
164 "Too many %s warnings (%d)\n"
165 "If you know what you are doing you can %s"
166 "set the environment variable GMX_MAXCONSTRWARN to -1,\n"
167 "but normally it is better to fix the problem",
168 (eConstrAlg
== econtLINCS
) ? "LINCS" : "SETTLE", warncount
,
169 (eConstrAlg
== econtLINCS
) ?
170 "adjust the lincs warning threshold in your mdp file\nor " : "\n");
173 static void write_constr_pdb(const char *fn
, const char *title
,
174 const gmx_mtop_t
*mtop
,
175 int start
, int homenr
, t_commrec
*cr
,
176 rvec x
[], matrix box
)
180 int dd_ac0
= 0, dd_ac1
= 0, i
, ii
, resnr
;
185 if (DOMAINDECOMP(cr
))
188 dd_get_constraint_range(dd
, &dd_ac0
, &dd_ac1
);
195 sprintf(fname
, "%s_n%d.pdb", fn
, cr
->sim_nodeid
);
199 sprintf(fname
, "%s.pdb", fn
);
202 out
= gmx_fio_fopen(fname
, "w");
204 fprintf(out
, "TITLE %s\n", title
);
205 gmx_write_pdb_box(out
, -1, box
);
206 for (i
= start
; i
< start
+homenr
; i
++)
210 if (i
>= dd
->nat_home
&& i
< dd_ac0
)
214 ii
= dd
->gatindex
[i
];
220 gmx_mtop_atominfo_global(mtop
, ii
, &anm
, &resnr
, &resnm
);
221 gmx_fprintf_pdb_atomline(out
, epdbATOM
, ii
+1, anm
, ' ', resnm
, ' ', resnr
, ' ',
222 10*x
[i
][XX
], 10*x
[i
][YY
], 10*x
[i
][ZZ
], 1.0, 0.0, "");
224 fprintf(out
, "TER\n");
229 static void dump_confs(FILE *fplog
, gmx_int64_t step
, const gmx_mtop_t
*mtop
,
230 int start
, int homenr
, t_commrec
*cr
,
231 rvec x
[], rvec xprime
[], matrix box
)
233 char buf
[256], buf2
[22];
235 char *env
= getenv("GMX_SUPPRESS_DUMP");
241 sprintf(buf
, "step%sb", gmx_step_str(step
, buf2
));
242 write_constr_pdb(buf
, "initial coordinates",
243 mtop
, start
, homenr
, cr
, x
, box
);
244 sprintf(buf
, "step%sc", gmx_step_str(step
, buf2
));
245 write_constr_pdb(buf
, "coordinates after constraining",
246 mtop
, start
, homenr
, cr
, xprime
, box
);
249 fprintf(fplog
, "Wrote pdb files with previous and current coordinates\n");
251 fprintf(stderr
, "Wrote pdb files with previous and current coordinates\n");
254 static void pr_sortblock(FILE *fp
, const char *title
, int nsb
, t_sortblock sb
[])
258 fprintf(fp
, "%s\n", title
);
259 for (i
= 0; (i
< nsb
); i
++)
261 fprintf(fp
, "i: %5d, iatom: (%5d %5d %5d), blocknr: %5d\n",
262 i
, sb
[i
].iatom
[0], sb
[i
].iatom
[1], sb
[i
].iatom
[2],
267 gmx_bool
constrain(FILE *fplog
, gmx_bool bLog
, gmx_bool bEner
,
268 struct gmx_constr
*constr
,
269 t_idef
*idef
, t_inputrec
*ir
,
271 gmx_int64_t step
, int delta_step
,
274 rvec
*x
, rvec
*xprime
, rvec
*min_proj
,
275 gmx_bool bMolPBC
, matrix box
,
276 real lambda
, real
*dvdlambda
,
277 rvec
*v
, tensor
*vir
,
278 t_nrnb
*nrnb
, int econq
)
286 real invdt
, vir_fac
= 0, t
;
289 t_pbc pbc
, *pbc_null
;
293 if (econq
== econqForceDispl
&& !EI_ENERGY_MINIMIZATION(ir
->eI
))
295 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");
304 scaled_delta_t
= step_scaling
* ir
->delta_t
;
306 /* Prepare time step for use in constraint implementations, and
307 avoid generating inf when ir->delta_t = 0. */
308 if (ir
->delta_t
== 0)
314 invdt
= 1.0/scaled_delta_t
;
317 if (ir
->efep
!= efepNO
&& EI_DYNAMICS(ir
->eI
))
319 /* Set the constraint lengths for the step at which this configuration
320 * is meant to be. The invmasses should not be changed.
322 lambda
+= delta_step
*ir
->fepvals
->delta_lambda
;
327 clear_mat(vir_r_m_dr
);
332 settle
= &idef
->il
[F_SETTLE
];
333 nsettle
= settle
->nr
/(1+NRAL(F_SETTLE
));
337 nth
= gmx_omp_nthreads_get(emntSETTLE
);
344 if (nth
> 1 && constr
->vir_r_m_dr_th
== NULL
)
346 snew(constr
->vir_r_m_dr_th
, nth
);
347 snew(constr
->settle_error
, nth
);
352 /* We do not need full pbc when constraints do not cross charge groups,
353 * i.e. when dd->constraint_comm==NULL.
354 * Note that PBC for constraints is different from PBC for bondeds.
355 * For constraints there is both forward and backward communication.
357 if (ir
->ePBC
!= epbcNONE
&&
358 (cr
->dd
|| bMolPBC
) && !(cr
->dd
&& cr
->dd
->constraint_comm
== NULL
))
360 /* With pbc=screw the screw has been changed to a shift
361 * by the constraint coordinate communication routine,
362 * so that here we can use normal pbc.
364 pbc_null
= set_pbc_dd(&pbc
, ir
->ePBC
,
365 DOMAINDECOMP(cr
) ? cr
->dd
->nc
: nullptr,
373 /* Communicate the coordinates required for the non-local constraints
374 * for LINCS and/or SETTLE.
378 dd_move_x_constraints(cr
->dd
, box
, x
, xprime
, econq
== econqCoord
);
382 /* We need to initialize the non-local components of v.
383 * We never actually use these values, but we do increment them,
384 * so we should avoid uninitialized variables and overflows.
386 clear_constraint_quantity_nonlocal(cr
->dd
, v
);
390 if (constr
->lincsd
!= NULL
)
392 bOK
= constrain_lincs(fplog
, bLog
, bEner
, ir
, step
, constr
->lincsd
, md
, cr
,
394 box
, pbc_null
, lambda
, dvdlambda
,
395 invdt
, v
, vir
!= NULL
, vir_r_m_dr
,
397 constr
->maxwarn
, &constr
->warncount_lincs
);
398 if (!bOK
&& constr
->maxwarn
>= 0)
402 fprintf(fplog
, "Constraint error in algorithm %s at step %s\n",
403 econstr_names
[econtLINCS
], gmx_step_str(step
, buf
));
409 if (constr
->nblocks
> 0)
414 bOK
= bshakef(fplog
, constr
->shaked
,
415 md
->invmass
, constr
->nblocks
, constr
->sblock
,
416 idef
, ir
, x
, xprime
, nrnb
,
417 constr
->lagr
, lambda
, dvdlambda
,
418 invdt
, v
, vir
!= NULL
, vir_r_m_dr
,
419 constr
->maxwarn
>= 0, econq
);
422 bOK
= bshakef(fplog
, constr
->shaked
,
423 md
->invmass
, constr
->nblocks
, constr
->sblock
,
424 idef
, ir
, x
, min_proj
, nrnb
,
425 constr
->lagr
, lambda
, dvdlambda
,
426 invdt
, NULL
, vir
!= NULL
, vir_r_m_dr
,
427 constr
->maxwarn
>= 0, econq
);
430 gmx_fatal(FARGS
, "Internal error, SHAKE called for constraining something else than coordinates");
434 if (!bOK
&& constr
->maxwarn
>= 0)
438 fprintf(fplog
, "Constraint error in algorithm %s at step %s\n",
439 econstr_names
[econtSHAKE
], gmx_step_str(step
, buf
));
447 int calcvir_atom_end
;
451 calcvir_atom_end
= 0;
455 calcvir_atom_end
= md
->homenr
;
461 #pragma omp parallel for num_threads(nth) schedule(static)
462 for (th
= 0; th
< nth
; th
++)
466 int start_th
, end_th
;
470 clear_mat(constr
->vir_r_m_dr_th
[th
]);
473 start_th
= (nsettle
* th
)/nth
;
474 end_th
= (nsettle
*(th
+1))/nth
;
475 if (start_th
>= 0 && end_th
- start_th
> 0)
477 csettle(constr
->settled
,
479 settle
->iatoms
+start_th
*(1+NRAL(F_SETTLE
)),
482 invdt
, v
? v
[0] : NULL
, calcvir_atom_end
,
483 th
== 0 ? vir_r_m_dr
: constr
->vir_r_m_dr_th
[th
],
484 th
== 0 ? &settle_error
: &constr
->settle_error
[th
]);
487 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
;
489 inc_nrnb(nrnb
, eNR_SETTLE
, nsettle
);
492 inc_nrnb(nrnb
, eNR_CONSTR_V
, nsettle
*3);
496 inc_nrnb(nrnb
, eNR_CONSTR_VIR
, nsettle
*3);
502 case econqForceDispl
:
503 #pragma omp parallel for num_threads(nth) schedule(static)
504 for (th
= 0; th
< nth
; th
++)
508 int start_th
, end_th
;
512 clear_mat(constr
->vir_r_m_dr_th
[th
]);
515 start_th
= (nsettle
* th
)/nth
;
516 end_th
= (nsettle
*(th
+1))/nth
;
518 if (start_th
>= 0 && end_th
- start_th
> 0)
520 settle_proj(constr
->settled
, econq
,
522 settle
->iatoms
+start_th
*(1+NRAL(F_SETTLE
)),
525 xprime
, min_proj
, calcvir_atom_end
,
526 th
== 0 ? vir_r_m_dr
: constr
->vir_r_m_dr_th
[th
]);
529 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
;
531 /* This is an overestimate */
532 inc_nrnb(nrnb
, eNR_SETTLE
, nsettle
);
534 case econqDeriv_FlexCon
:
535 /* Nothing to do, since the are no flexible constraints in settles */
538 gmx_incons("Unknown constraint quantity for settle");
544 /* Combine virial and error info of the other threads */
545 for (i
= 1; i
< nth
; i
++)
547 settle_error
= constr
->settle_error
[i
];
551 for (i
= 1; i
< nth
; i
++)
553 m_add(vir_r_m_dr
, constr
->vir_r_m_dr_th
[i
], vir_r_m_dr
);
557 if (econq
== econqCoord
&& settle_error
>= 0)
560 if (constr
->maxwarn
>= 0)
564 "\nstep " "%" GMX_PRId64
": Water molecule starting at atom %d can not be "
565 "settled.\nCheck for bad contacts and/or reduce the timestep if appropriate.\n",
566 step
, ddglatnr(cr
->dd
, settle
->iatoms
[settle_error
*(1+NRAL(F_SETTLE
))+1]));
569 fprintf(fplog
, "%s", buf
);
571 fprintf(stderr
, "%s", buf
);
572 constr
->warncount_settle
++;
573 if (constr
->warncount_settle
> constr
->maxwarn
)
575 too_many_constraint_warnings(-1, constr
->warncount_settle
);
584 /* The normal uses of constrain() pass step_scaling = 1.0.
585 * The call to constrain() for SD1 that passes step_scaling =
586 * 0.5 also passes vir = NULL, so cannot reach this
587 * assertion. This assertion should remain until someone knows
588 * that this path works for their intended purpose, and then
589 * they can use scaled_delta_t instead of ir->delta_t
591 assert(gmx_within_tol(step_scaling
, 1.0, GMX_REAL_EPS
));
595 vir_fac
= 0.5/(ir
->delta_t
*ir
->delta_t
);
598 vir_fac
= 0.5/ir
->delta_t
;
601 case econqForceDispl
:
605 gmx_incons("Unsupported constraint quantity for virial");
610 vir_fac
*= 2; /* only constraining over half the distance here */
612 for (i
= 0; i
< DIM
; i
++)
614 for (j
= 0; j
< DIM
; j
++)
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
->bPull
&& pull_have_constraint(ir
->pull_work
))
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_work
, 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
)
655 return lincs_rmsd_data(constr
->lincsd
);
663 real
constr_rmsd(struct gmx_constr
*constr
, gmx_bool bSD2
)
667 return lincs_rmsd(constr
->lincsd
, bSD2
);
675 static void make_shake_sblock_serial(struct gmx_constr
*constr
,
676 t_idef
*idef
, const t_mdatoms
*md
)
685 /* Since we are processing the local topology,
686 * the F_CONSTRNC ilist has been concatenated to the F_CONSTR ilist.
688 ncons
= idef
->il
[F_CONSTR
].nr
/3;
690 init_blocka(&sblocks
);
691 gen_sblocks(NULL
, 0, md
->homenr
, idef
, &sblocks
, FALSE
);
694 bstart=(idef->nodeid > 0) ? blocks->multinr[idef->nodeid-1] : 0;
695 nblocks=blocks->multinr[idef->nodeid] - bstart;
698 constr
->nblocks
= sblocks
.nr
;
701 fprintf(debug
, "ncons: %d, bstart: %d, nblocks: %d\n",
702 ncons
, bstart
, constr
->nblocks
);
705 /* Calculate block number for each atom */
706 inv_sblock
= make_invblocka(&sblocks
, md
->nr
);
708 done_blocka(&sblocks
);
710 /* Store the block number in temp array and
711 * sort the constraints in order of the sblock number
712 * and the atom numbers, really sorting a segment of the array!
715 pr_idef(fplog
, 0, "Before Sort", idef
);
717 iatom
= idef
->il
[F_CONSTR
].iatoms
;
719 for (i
= 0; (i
< ncons
); i
++, iatom
+= 3)
721 for (m
= 0; (m
< 3); m
++)
723 sb
[i
].iatom
[m
] = iatom
[m
];
725 sb
[i
].blocknr
= inv_sblock
[iatom
[1]];
728 /* Now sort the blocks */
731 pr_sortblock(debug
, "Before sorting", ncons
, sb
);
732 fprintf(debug
, "Going to sort constraints\n");
735 qsort(sb
, ncons
, (size_t)sizeof(*sb
), pcomp
);
739 pr_sortblock(debug
, "After sorting", ncons
, sb
);
742 iatom
= idef
->il
[F_CONSTR
].iatoms
;
743 for (i
= 0; (i
< ncons
); i
++, iatom
+= 3)
745 for (m
= 0; (m
< 3); m
++)
747 iatom
[m
] = sb
[i
].iatom
[m
];
751 pr_idef(fplog
, 0, "After Sort", idef
);
755 snew(constr
->sblock
, constr
->nblocks
+1);
757 for (i
= 0; (i
< ncons
); i
++)
759 if (sb
[i
].blocknr
!= bnr
)
762 constr
->sblock
[j
++] = 3*i
;
766 constr
->sblock
[j
++] = 3*ncons
;
768 if (j
!= (constr
->nblocks
+1))
770 fprintf(stderr
, "bstart: %d\n", bstart
);
771 fprintf(stderr
, "j: %d, nblocks: %d, ncons: %d\n",
772 j
, constr
->nblocks
, ncons
);
773 for (i
= 0; (i
< ncons
); i
++)
775 fprintf(stderr
, "i: %5d sb[i].blocknr: %5d\n", i
, sb
[i
].blocknr
);
777 for (j
= 0; (j
<= constr
->nblocks
); j
++)
779 fprintf(stderr
, "sblock[%3d]=%5d\n", j
, (int)constr
->sblock
[j
]);
781 gmx_fatal(FARGS
, "DEATH HORROR: "
782 "sblocks does not match idef->il[F_CONSTR]");
788 static void make_shake_sblock_dd(struct gmx_constr
*constr
,
789 const t_ilist
*ilcon
, const t_block
*cgs
,
790 const gmx_domdec_t
*dd
)
795 if (dd
->ncg_home
+1 > constr
->sblock_nalloc
)
797 constr
->sblock_nalloc
= over_alloc_dd(dd
->ncg_home
+1);
798 srenew(constr
->sblock
, constr
->sblock_nalloc
);
802 iatom
= ilcon
->iatoms
;
805 for (c
= 0; c
< ncons
; c
++)
807 if (c
== 0 || iatom
[1] >= cgs
->index
[cg
+1])
809 constr
->sblock
[constr
->nblocks
++] = 3*c
;
810 while (iatom
[1] >= cgs
->index
[cg
+1])
817 constr
->sblock
[constr
->nblocks
] = 3*ncons
;
820 t_blocka
make_at2con(int start
, int natoms
,
821 const t_ilist
*ilist
, const t_iparams
*iparams
,
822 gmx_bool bDynamics
, int *nflexiblecons
)
824 int *count
, ncon
, con
, con_tot
, nflexcon
, ftype
, i
, a
;
831 for (ftype
= F_CONSTR
; ftype
<= F_CONSTRNC
; ftype
++)
833 ncon
= ilist
[ftype
].nr
/3;
834 ia
= ilist
[ftype
].iatoms
;
835 for (con
= 0; con
< ncon
; con
++)
837 bFlexCon
= (iparams
[ia
[0]].constr
.dA
== 0 &&
838 iparams
[ia
[0]].constr
.dB
== 0);
843 if (bDynamics
|| !bFlexCon
)
845 for (i
= 1; i
< 3; i
++)
854 *nflexiblecons
= nflexcon
;
857 at2con
.nalloc_index
= at2con
.nr
+1;
858 snew(at2con
.index
, at2con
.nalloc_index
);
860 for (a
= 0; a
< natoms
; a
++)
862 at2con
.index
[a
+1] = at2con
.index
[a
] + count
[a
];
865 at2con
.nra
= at2con
.index
[natoms
];
866 at2con
.nalloc_a
= at2con
.nra
;
867 snew(at2con
.a
, at2con
.nalloc_a
);
869 /* The F_CONSTRNC constraints have constraint numbers
870 * that continue after the last F_CONSTR constraint.
873 for (ftype
= F_CONSTR
; ftype
<= F_CONSTRNC
; ftype
++)
875 ncon
= ilist
[ftype
].nr
/3;
876 ia
= ilist
[ftype
].iatoms
;
877 for (con
= 0; con
< ncon
; con
++)
879 bFlexCon
= (iparams
[ia
[0]].constr
.dA
== 0 &&
880 iparams
[ia
[0]].constr
.dB
== 0);
881 if (bDynamics
|| !bFlexCon
)
883 for (i
= 1; i
< 3; i
++)
886 at2con
.a
[at2con
.index
[a
]+count
[a
]++] = con_tot
;
899 static int *make_at2settle(int natoms
, const t_ilist
*ilist
)
905 /* Set all to no settle */
906 for (a
= 0; a
< natoms
; a
++)
911 stride
= 1 + NRAL(F_SETTLE
);
913 for (s
= 0; s
< ilist
->nr
; s
+= stride
)
915 at2s
[ilist
->iatoms
[s
+1]] = s
/stride
;
916 at2s
[ilist
->iatoms
[s
+2]] = s
/stride
;
917 at2s
[ilist
->iatoms
[s
+3]] = s
/stride
;
923 void set_constraints(struct gmx_constr
*constr
,
924 gmx_localtop_t
*top
, const t_inputrec
*ir
,
925 const t_mdatoms
*md
, t_commrec
*cr
)
929 const t_ilist
*settle
;
934 if (constr
->ncon_tot
> 0)
936 /* We are using the local topology,
937 * so there are only F_CONSTR constraints.
939 ncons
= idef
->il
[F_CONSTR
].nr
/3;
941 /* With DD we might also need to call LINCS with ncons=0 for
942 * communicating coordinates to other nodes that do have constraints.
944 if (ir
->eConstrAlg
== econtLINCS
)
946 set_lincs(idef
, md
, EI_DYNAMICS(ir
->eI
), cr
, constr
->lincsd
);
948 if (ir
->eConstrAlg
== econtSHAKE
)
952 make_shake_sblock_dd(constr
, &idef
->il
[F_CONSTR
], &top
->cgs
, cr
->dd
);
956 make_shake_sblock_serial(constr
, idef
, md
);
958 if (ncons
> constr
->lagr_nalloc
)
960 constr
->lagr_nalloc
= over_alloc_dd(ncons
);
961 srenew(constr
->lagr
, constr
->lagr_nalloc
);
966 if (idef
->il
[F_SETTLE
].nr
> 0 && constr
->settled
== NULL
)
968 settle
= &idef
->il
[F_SETTLE
];
969 iO
= settle
->iatoms
[1];
970 iH
= settle
->iatoms
[2];
972 settle_init(md
->massT
[iO
], md
->massT
[iH
],
973 md
->invmass
[iO
], md
->invmass
[iH
],
974 idef
->iparams
[settle
->iatoms
[0]].settle
.doh
,
975 idef
->iparams
[settle
->iatoms
[0]].settle
.dhh
);
978 /* Make a selection of the local atoms for essential dynamics */
979 if (constr
->ed
&& cr
->dd
)
981 dd_make_local_ed_indices(cr
->dd
, constr
->ed
);
985 static void constr_recur(const t_blocka
*at2con
,
986 const t_ilist
*ilist
, const t_iparams
*iparams
,
988 int at
, int depth
, int nc
, int *path
,
989 real r0
, real r1
, real
*r2max
,
1001 ncon1
= ilist
[F_CONSTR
].nr
/3;
1002 ia1
= ilist
[F_CONSTR
].iatoms
;
1003 ia2
= ilist
[F_CONSTRNC
].iatoms
;
1005 /* Loop over all constraints connected to this atom */
1006 for (c
= at2con
->index
[at
]; c
< at2con
->index
[at
+1]; c
++)
1009 /* Do not walk over already used constraints */
1011 for (a1
= 0; a1
< depth
; a1
++)
1013 if (con
== path
[a1
])
1020 ia
= constr_iatomptr(ncon1
, ia1
, ia2
, con
);
1021 /* Flexible constraints currently have length 0, which is incorrect */
1024 len
= iparams
[ia
[0]].constr
.dA
;
1028 len
= iparams
[ia
[0]].constr
.dB
;
1030 /* In the worst case the bond directions alternate */
1041 /* Assume angles of 120 degrees between all bonds */
1042 if (rn0
*rn0
+ rn1
*rn1
+ rn0
*rn1
> *r2max
)
1044 *r2max
= rn0
*rn0
+ rn1
*rn1
+ r0
*rn1
;
1047 fprintf(debug
, "Found longer constraint distance: r0 %5.3f r1 %5.3f rmax %5.3f\n", rn0
, rn1
, sqrt(*r2max
));
1048 for (a1
= 0; a1
< depth
; a1
++)
1050 fprintf(debug
, " %d %5.3f",
1052 iparams
[constr_iatomptr(ncon1
, ia1
, ia2
, con
)[0]].constr
.dA
);
1054 fprintf(debug
, " %d %5.3f\n", con
, len
);
1057 /* Limit the number of recursions to 1000*nc,
1058 * so a call does not take more than a second,
1059 * even for highly connected systems.
1061 if (depth
+ 1 < nc
&& *count
< 1000*nc
)
1073 constr_recur(at2con
, ilist
, iparams
,
1074 bTopB
, a1
, depth
+1, nc
, path
, rn0
, rn1
, r2max
, count
);
1081 static real
constr_r_max_moltype(const gmx_moltype_t
*molt
,
1082 const t_iparams
*iparams
,
1083 const t_inputrec
*ir
)
1085 int natoms
, nflexcon
, *path
, at
, count
;
1088 real r0
, r1
, r2maxA
, r2maxB
, rmax
, lam0
, lam1
;
1090 if (molt
->ilist
[F_CONSTR
].nr
== 0 &&
1091 molt
->ilist
[F_CONSTRNC
].nr
== 0)
1096 natoms
= molt
->atoms
.nr
;
1098 at2con
= make_at2con(0, natoms
, molt
->ilist
, iparams
,
1099 EI_DYNAMICS(ir
->eI
), &nflexcon
);
1100 snew(path
, 1+ir
->nProjOrder
);
1101 for (at
= 0; at
< 1+ir
->nProjOrder
; at
++)
1107 for (at
= 0; at
< natoms
; at
++)
1113 constr_recur(&at2con
, molt
->ilist
, iparams
,
1114 FALSE
, at
, 0, 1+ir
->nProjOrder
, path
, r0
, r1
, &r2maxA
, &count
);
1116 if (ir
->efep
== efepNO
)
1118 rmax
= sqrt(r2maxA
);
1123 for (at
= 0; at
< natoms
; at
++)
1128 constr_recur(&at2con
, molt
->ilist
, iparams
,
1129 TRUE
, at
, 0, 1+ir
->nProjOrder
, path
, r0
, r1
, &r2maxB
, &count
);
1131 lam0
= ir
->fepvals
->init_lambda
;
1132 if (EI_DYNAMICS(ir
->eI
))
1134 lam0
+= ir
->init_step
*ir
->fepvals
->delta_lambda
;
1136 rmax
= (1 - lam0
)*sqrt(r2maxA
) + lam0
*sqrt(r2maxB
);
1137 if (EI_DYNAMICS(ir
->eI
))
1139 lam1
= ir
->fepvals
->init_lambda
+ (ir
->init_step
+ ir
->nsteps
)*ir
->fepvals
->delta_lambda
;
1140 rmax
= std::max(rmax
, (1 - lam1
)*std::sqrt(r2maxA
) + lam1
*std::sqrt(r2maxB
));
1144 done_blocka(&at2con
);
1150 real
constr_r_max(FILE *fplog
, const gmx_mtop_t
*mtop
, const t_inputrec
*ir
)
1156 for (mt
= 0; mt
< mtop
->nmoltype
; mt
++)
1158 rmax
= std::max(rmax
,
1159 constr_r_max_moltype(&mtop
->moltype
[mt
],
1160 mtop
->ffparams
.iparams
, ir
));
1165 fprintf(fplog
, "Maximum distance for %d constraints, at 120 deg. angles, all-trans: %.3f nm\n", 1+ir
->nProjOrder
, rmax
);
1171 gmx_constr_t
init_constraints(FILE *fplog
,
1172 const gmx_mtop_t
*mtop
, const t_inputrec
*ir
,
1173 gmx_edsam_t ed
, t_state
*state
,
1176 int ncon
, nset
, nmol
, settle_type
, i
, mt
, nflexcon
;
1177 struct gmx_constr
*constr
;
1180 gmx_mtop_ilistloop_t iloop
;
1183 gmx_mtop_ftype_count(mtop
, F_CONSTR
) +
1184 gmx_mtop_ftype_count(mtop
, F_CONSTRNC
);
1185 nset
= gmx_mtop_ftype_count(mtop
, F_SETTLE
);
1187 if (ncon
+nset
== 0 &&
1188 !(ir
->bPull
&& pull_have_constraint(ir
->pull_work
)) &&
1196 constr
->ncon_tot
= ncon
;
1197 constr
->nflexcon
= 0;
1200 constr
->n_at2con_mt
= mtop
->nmoltype
;
1201 snew(constr
->at2con_mt
, constr
->n_at2con_mt
);
1202 for (mt
= 0; mt
< mtop
->nmoltype
; mt
++)
1204 constr
->at2con_mt
[mt
] = make_at2con(0, mtop
->moltype
[mt
].atoms
.nr
,
1205 mtop
->moltype
[mt
].ilist
,
1206 mtop
->ffparams
.iparams
,
1207 EI_DYNAMICS(ir
->eI
), &nflexcon
);
1208 for (i
= 0; i
< mtop
->nmolblock
; i
++)
1210 if (mtop
->molblock
[i
].type
== mt
)
1212 constr
->nflexcon
+= mtop
->molblock
[i
].nmol
*nflexcon
;
1217 if (constr
->nflexcon
> 0)
1221 fprintf(fplog
, "There are %d flexible constraints\n",
1223 if (ir
->fc_stepsize
== 0)
1226 "WARNING: step size for flexible constraining = 0\n"
1227 " All flexible constraints will be rigid.\n"
1228 " Will try to keep all flexible constraints at their original length,\n"
1229 " but the lengths may exhibit some drift.\n\n");
1230 constr
->nflexcon
= 0;
1233 if (constr
->nflexcon
> 0)
1235 please_cite(fplog
, "Hess2002");
1239 if (ir
->eConstrAlg
== econtLINCS
)
1241 constr
->lincsd
= init_lincs(fplog
, mtop
,
1242 constr
->nflexcon
, constr
->at2con_mt
,
1243 DOMAINDECOMP(cr
) && cr
->dd
->bInterCGcons
,
1244 ir
->nLincsIter
, ir
->nProjOrder
);
1247 if (ir
->eConstrAlg
== econtSHAKE
)
1249 if (DOMAINDECOMP(cr
) && cr
->dd
->bInterCGcons
)
1251 gmx_fatal(FARGS
, "SHAKE is not supported with domain decomposition and constraint that cross charge group boundaries, use LINCS");
1253 if (constr
->nflexcon
)
1255 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");
1257 please_cite(fplog
, "Ryckaert77a");
1260 please_cite(fplog
, "Barth95a");
1263 constr
->shaked
= shake_init();
1269 please_cite(fplog
, "Miyamoto92a");
1271 constr
->bInterCGsettles
= inter_charge_group_settles(mtop
);
1273 /* Check that we have only one settle type */
1275 iloop
= gmx_mtop_ilistloop_init(mtop
);
1276 while (gmx_mtop_ilistloop_next(iloop
, &ilist
, &nmol
))
1278 for (i
= 0; i
< ilist
[F_SETTLE
].nr
; i
+= 4)
1280 if (settle_type
== -1)
1282 settle_type
= ilist
[F_SETTLE
].iatoms
[i
];
1284 else if (ilist
[F_SETTLE
].iatoms
[i
] != settle_type
)
1287 "The [molecules] section of your topology specifies more than one block of\n"
1288 "a [moleculetype] with a [settles] block. Only one such is allowed. If you\n"
1289 "are trying to partition your solvent into different *groups* (e.g. for\n"
1290 "freezing, T-coupling, etc.) then you are using the wrong approach. Index\n"
1291 "files specify groups. Otherwise, you may wish to change the least-used\n"
1292 "block of molecules with SETTLE constraints into 3 normal constraints.");
1297 constr
->n_at2settle_mt
= mtop
->nmoltype
;
1298 snew(constr
->at2settle_mt
, constr
->n_at2settle_mt
);
1299 for (mt
= 0; mt
< mtop
->nmoltype
; mt
++)
1301 constr
->at2settle_mt
[mt
] =
1302 make_at2settle(mtop
->moltype
[mt
].atoms
.nr
,
1303 &mtop
->moltype
[mt
].ilist
[F_SETTLE
]);
1307 if ((ncon
+ nset
) > 0 && ir
->epc
== epcMTTK
)
1309 gmx_fatal(FARGS
, "Constraints are not implemented with MTTK pressure control.");
1312 constr
->maxwarn
= 999;
1313 env
= getenv("GMX_MAXCONSTRWARN");
1316 constr
->maxwarn
= 0;
1317 sscanf(env
, "%8d", &constr
->maxwarn
);
1321 "Setting the maximum number of constraint warnings to %d\n",
1327 "Setting the maximum number of constraint warnings to %d\n",
1331 if (constr
->maxwarn
< 0 && fplog
)
1333 fprintf(fplog
, "maxwarn < 0, will not stop on constraint errors\n");
1335 constr
->warncount_lincs
= 0;
1336 constr
->warncount_settle
= 0;
1338 /* Initialize the essential dynamics sampling.
1339 * Put the pointer to the ED struct in constr */
1341 if (ed
!= NULL
|| state
->edsamstate
.nED
> 0)
1343 init_edsam(mtop
, ir
, cr
, ed
, state
->x
, state
->box
, &state
->edsamstate
);
1346 constr
->warn_mtop
= mtop
;
1351 const t_blocka
*atom2constraints_moltype(gmx_constr_t constr
)
1353 return constr
->at2con_mt
;
1356 const int **atom2settle_moltype(gmx_constr_t constr
)
1358 return (const int **)constr
->at2settle_mt
;
1362 gmx_bool
inter_charge_group_constraints(const gmx_mtop_t
*mtop
)
1364 const gmx_moltype_t
*molt
;
1368 int *at2cg
, cg
, a
, ftype
, i
;
1372 for (mb
= 0; mb
< mtop
->nmolblock
&& !bInterCG
; mb
++)
1374 molt
= &mtop
->moltype
[mtop
->molblock
[mb
].type
];
1376 if (molt
->ilist
[F_CONSTR
].nr
> 0 ||
1377 molt
->ilist
[F_CONSTRNC
].nr
> 0 ||
1378 molt
->ilist
[F_SETTLE
].nr
> 0)
1381 snew(at2cg
, molt
->atoms
.nr
);
1382 for (cg
= 0; cg
< cgs
->nr
; cg
++)
1384 for (a
= cgs
->index
[cg
]; a
< cgs
->index
[cg
+1]; a
++)
1390 for (ftype
= F_CONSTR
; ftype
<= F_CONSTRNC
; ftype
++)
1392 il
= &molt
->ilist
[ftype
];
1393 for (i
= 0; i
< il
->nr
&& !bInterCG
; i
+= 1+NRAL(ftype
))
1395 if (at2cg
[il
->iatoms
[i
+1]] != at2cg
[il
->iatoms
[i
+2]])
1409 gmx_bool
inter_charge_group_settles(const gmx_mtop_t
*mtop
)
1411 const gmx_moltype_t
*molt
;
1415 int *at2cg
, cg
, a
, ftype
, i
;
1419 for (mb
= 0; mb
< mtop
->nmolblock
&& !bInterCG
; mb
++)
1421 molt
= &mtop
->moltype
[mtop
->molblock
[mb
].type
];
1423 if (molt
->ilist
[F_SETTLE
].nr
> 0)
1426 snew(at2cg
, molt
->atoms
.nr
);
1427 for (cg
= 0; cg
< cgs
->nr
; cg
++)
1429 for (a
= cgs
->index
[cg
]; a
< cgs
->index
[cg
+1]; a
++)
1435 for (ftype
= F_SETTLE
; ftype
<= F_SETTLE
; ftype
++)
1437 il
= &molt
->ilist
[ftype
];
1438 for (i
= 0; i
< il
->nr
&& !bInterCG
; i
+= 1+NRAL(F_SETTLE
))
1440 if (at2cg
[il
->iatoms
[i
+1]] != at2cg
[il
->iatoms
[i
+2]] ||
1441 at2cg
[il
->iatoms
[i
+1]] != at2cg
[il
->iatoms
[i
+3]])
1455 /* helper functions for andersen temperature control, because the
1456 * gmx_constr construct is only defined in constr.c. Return the list
1457 * of blocks (get_sblock) and the number of blocks (get_nblocks). */
1459 extern int *get_sblock(struct gmx_constr
*constr
)
1461 return constr
->sblock
;
1464 extern int get_nblocks(struct gmx_constr
*constr
)
1466 return constr
->nblocks
;