Moved copyrite.* to fileio from gmxlib and legacyheaders.
[gromacs.git] / src / gromacs / mdlib / constr.cpp
blob9194b1002eef92277fe4ac91162279d92332ecb5
1 /*
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.
37 #include "gmxpre.h"
39 #include "constr.h"
41 #include <assert.h>
42 #include <stdlib.h>
44 #include <cmath>
46 #include <algorithm>
48 #include "gromacs/domdec/domdec.h"
49 #include "gromacs/essentialdynamics/edsam.h"
50 #include "gromacs/fileio/confio.h"
51 #include "gromacs/fileio/copyrite.h"
52 #include "gromacs/fileio/gmxfio.h"
53 #include "gromacs/fileio/pdbio.h"
54 #include "gromacs/gmxlib/gmx_omp_nthreads.h"
55 #include "gromacs/gmxlib/splitter.h"
56 #include "gromacs/legacyheaders/names.h"
57 #include "gromacs/legacyheaders/nrnb.h"
58 #include "gromacs/legacyheaders/txtdump.h"
59 #include "gromacs/legacyheaders/types/commrec.h"
60 #include "gromacs/math/vec.h"
61 #include "gromacs/mdlib/mdrun.h"
62 #include "gromacs/pbcutil/pbc.h"
63 #include "gromacs/pulling/pull.h"
64 #include "gromacs/topology/block.h"
65 #include "gromacs/topology/invblock.h"
66 #include "gromacs/topology/mtop_util.h"
67 #include "gromacs/utility/exceptions.h"
68 #include "gromacs/utility/fatalerror.h"
69 #include "gromacs/utility/smalloc.h"
71 typedef struct gmx_constr {
72 int ncon_tot; /* The total number of constraints */
73 int nflexcon; /* The number of flexible constraints */
74 int n_at2con_mt; /* The size of at2con = #moltypes */
75 t_blocka *at2con_mt; /* A list of atoms to constraints */
76 int n_at2settle_mt; /* The size of at2settle = #moltypes */
77 int **at2settle_mt; /* A list of atoms to settles */
78 gmx_bool bInterCGsettles;
79 gmx_lincsdata_t lincsd; /* LINCS data */
80 gmx_shakedata_t shaked; /* SHAKE data */
81 gmx_settledata_t settled; /* SETTLE data */
82 int nblocks; /* The number of SHAKE blocks */
83 int *sblock; /* The SHAKE blocks */
84 int sblock_nalloc; /* The allocation size of sblock */
85 real *lagr; /* -2 times the Lagrange multipliers for SHAKE */
86 int lagr_nalloc; /* The allocation size of lagr */
87 int maxwarn; /* The maximum number of warnings */
88 int warncount_lincs;
89 int warncount_settle;
90 gmx_edsam_t ed; /* The essential dynamics data */
92 tensor *vir_r_m_dr_th; /* Thread local working data */
93 int *settle_error; /* Thread local working data */
95 const gmx_mtop_t *warn_mtop; /* Only used for printing warnings */
96 } t_gmx_constr;
98 typedef struct {
99 atom_id iatom[3];
100 atom_id blocknr;
101 } t_sortblock;
103 static int pcomp(const void *p1, const void *p2)
105 int db;
106 atom_id min1, min2, max1, max2;
107 t_sortblock *a1 = (t_sortblock *)p1;
108 t_sortblock *a2 = (t_sortblock *)p2;
110 db = a1->blocknr-a2->blocknr;
112 if (db != 0)
114 return db;
117 min1 = std::min(a1->iatom[1], a1->iatom[2]);
118 max1 = std::max(a1->iatom[1], a1->iatom[2]);
119 min2 = std::min(a2->iatom[1], a2->iatom[2]);
120 max2 = std::max(a2->iatom[1], a2->iatom[2]);
122 if (min1 == min2)
124 return max1-max2;
126 else
128 return min1-min2;
132 int n_flexible_constraints(struct gmx_constr *constr)
134 int nflexcon;
136 if (constr)
138 nflexcon = constr->nflexcon;
140 else
142 nflexcon = 0;
145 return nflexcon;
148 static void clear_constraint_quantity_nonlocal(gmx_domdec_t *dd, rvec *q)
150 int nonlocal_at_start, nonlocal_at_end, at;
152 dd_get_constraint_range(dd, &nonlocal_at_start, &nonlocal_at_end);
154 for (at = nonlocal_at_start; at < nonlocal_at_end; at++)
156 clear_rvec(q[at]);
160 void too_many_constraint_warnings(int eConstrAlg, int warncount)
162 gmx_fatal(FARGS,
163 "Too many %s warnings (%d)\n"
164 "If you know what you are doing you can %s"
165 "set the environment variable GMX_MAXCONSTRWARN to -1,\n"
166 "but normally it is better to fix the problem",
167 (eConstrAlg == econtLINCS) ? "LINCS" : "SETTLE", warncount,
168 (eConstrAlg == econtLINCS) ?
169 "adjust the lincs warning threshold in your mdp file\nor " : "\n");
172 static void write_constr_pdb(const char *fn, const char *title,
173 const gmx_mtop_t *mtop,
174 int start, int homenr, t_commrec *cr,
175 rvec x[], matrix box)
177 char fname[STRLEN];
178 FILE *out;
179 int dd_ac0 = 0, dd_ac1 = 0, i, ii, resnr;
180 gmx_domdec_t *dd;
181 char *anm, *resnm;
183 dd = NULL;
184 if (DOMAINDECOMP(cr))
186 dd = cr->dd;
187 dd_get_constraint_range(dd, &dd_ac0, &dd_ac1);
188 start = 0;
189 homenr = dd_ac1;
192 if (PAR(cr))
194 sprintf(fname, "%s_n%d.pdb", fn, cr->sim_nodeid);
196 else
198 sprintf(fname, "%s.pdb", fn);
201 out = gmx_fio_fopen(fname, "w");
203 fprintf(out, "TITLE %s\n", title);
204 gmx_write_pdb_box(out, -1, box);
205 for (i = start; i < start+homenr; i++)
207 if (dd != NULL)
209 if (i >= dd->nat_home && i < dd_ac0)
211 continue;
213 ii = dd->gatindex[i];
215 else
217 ii = i;
219 gmx_mtop_atominfo_global(mtop, ii, &anm, &resnr, &resnm);
220 gmx_fprintf_pdb_atomline(out, epdbATOM, ii+1, anm, ' ', resnm, ' ', resnr, ' ',
221 10*x[i][XX], 10*x[i][YY], 10*x[i][ZZ], 1.0, 0.0, "");
223 fprintf(out, "TER\n");
225 gmx_fio_fclose(out);
228 static void dump_confs(FILE *fplog, gmx_int64_t step, const gmx_mtop_t *mtop,
229 int start, int homenr, t_commrec *cr,
230 rvec x[], rvec xprime[], matrix box)
232 char buf[256], buf2[22];
234 char *env = getenv("GMX_SUPPRESS_DUMP");
235 if (env)
237 return;
240 sprintf(buf, "step%sb", gmx_step_str(step, buf2));
241 write_constr_pdb(buf, "initial coordinates",
242 mtop, start, homenr, cr, x, box);
243 sprintf(buf, "step%sc", gmx_step_str(step, buf2));
244 write_constr_pdb(buf, "coordinates after constraining",
245 mtop, start, homenr, cr, xprime, box);
246 if (fplog)
248 fprintf(fplog, "Wrote pdb files with previous and current coordinates\n");
250 fprintf(stderr, "Wrote pdb files with previous and current coordinates\n");
253 static void pr_sortblock(FILE *fp, const char *title, int nsb, t_sortblock sb[])
255 int i;
257 fprintf(fp, "%s\n", title);
258 for (i = 0; (i < nsb); i++)
260 fprintf(fp, "i: %5d, iatom: (%5d %5d %5d), blocknr: %5d\n",
261 i, sb[i].iatom[0], sb[i].iatom[1], sb[i].iatom[2],
262 sb[i].blocknr);
266 gmx_bool constrain(FILE *fplog, gmx_bool bLog, gmx_bool bEner,
267 struct gmx_constr *constr,
268 t_idef *idef, t_inputrec *ir,
269 t_commrec *cr,
270 gmx_int64_t step, int delta_step,
271 real step_scaling,
272 t_mdatoms *md,
273 rvec *x, rvec *xprime, rvec *min_proj,
274 gmx_bool bMolPBC, matrix box,
275 real lambda, real *dvdlambda,
276 rvec *v, tensor *vir,
277 t_nrnb *nrnb, int econq)
279 gmx_bool bOK, bDump;
280 int start, homenr;
281 int i, j;
282 int settle_error;
283 tensor vir_r_m_dr;
284 real scaled_delta_t;
285 real invdt, vir_fac = 0, t;
286 t_ilist *settle;
287 int nsettle;
288 t_pbc pbc, *pbc_null;
289 char buf[22];
290 int nth, th;
292 if (econq == econqForceDispl && !EI_ENERGY_MINIMIZATION(ir->eI))
294 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");
297 bOK = TRUE;
298 bDump = FALSE;
300 start = 0;
301 homenr = md->homenr;
303 scaled_delta_t = step_scaling * ir->delta_t;
305 /* Prepare time step for use in constraint implementations, and
306 avoid generating inf when ir->delta_t = 0. */
307 if (ir->delta_t == 0)
309 invdt = 0.0;
311 else
313 invdt = 1.0/scaled_delta_t;
316 if (ir->efep != efepNO && EI_DYNAMICS(ir->eI))
318 /* Set the constraint lengths for the step at which this configuration
319 * is meant to be. The invmasses should not be changed.
321 lambda += delta_step*ir->fepvals->delta_lambda;
324 if (vir != NULL)
326 clear_mat(vir_r_m_dr);
329 where();
331 settle = &idef->il[F_SETTLE];
332 nsettle = settle->nr/(1+NRAL(F_SETTLE));
334 if (nsettle > 0)
336 nth = gmx_omp_nthreads_get(emntSETTLE);
338 else
340 nth = 1;
343 if (nth > 1 && constr->vir_r_m_dr_th == NULL)
345 snew(constr->vir_r_m_dr_th, nth);
346 snew(constr->settle_error, nth);
349 settle_error = -1;
351 /* We do not need full pbc when constraints do not cross charge groups,
352 * i.e. when dd->constraint_comm==NULL.
353 * Note that PBC for constraints is different from PBC for bondeds.
354 * For constraints there is both forward and backward communication.
356 if (ir->ePBC != epbcNONE &&
357 (cr->dd || bMolPBC) && !(cr->dd && cr->dd->constraint_comm == NULL))
359 /* With pbc=screw the screw has been changed to a shift
360 * by the constraint coordinate communication routine,
361 * so that here we can use normal pbc.
363 pbc_null = set_pbc_dd(&pbc, ir->ePBC, cr->dd, FALSE, box);
365 else
367 pbc_null = NULL;
370 /* Communicate the coordinates required for the non-local constraints
371 * for LINCS and/or SETTLE.
373 if (cr->dd)
375 dd_move_x_constraints(cr->dd, box, x, xprime, econq == econqCoord);
377 if (v != NULL)
379 /* We need to initialize the non-local components of v.
380 * We never actually use these values, but we do increment them,
381 * so we should avoid uninitialized variables and overflows.
383 clear_constraint_quantity_nonlocal(cr->dd, v);
387 if (constr->lincsd != NULL)
389 bOK = constrain_lincs(fplog, bLog, bEner, ir, step, constr->lincsd, md, cr,
390 x, xprime, min_proj,
391 box, pbc_null, lambda, dvdlambda,
392 invdt, v, vir != NULL, vir_r_m_dr,
393 econq, nrnb,
394 constr->maxwarn, &constr->warncount_lincs);
395 if (!bOK && constr->maxwarn >= 0)
397 if (fplog != NULL)
399 fprintf(fplog, "Constraint error in algorithm %s at step %s\n",
400 econstr_names[econtLINCS], gmx_step_str(step, buf));
402 bDump = TRUE;
406 if (constr->nblocks > 0)
408 switch (econq)
410 case (econqCoord):
411 bOK = bshakef(fplog, constr->shaked,
412 md->invmass, constr->nblocks, constr->sblock,
413 idef, ir, x, xprime, nrnb,
414 constr->lagr, lambda, dvdlambda,
415 invdt, v, vir != NULL, vir_r_m_dr,
416 constr->maxwarn >= 0, econq);
417 break;
418 case (econqVeloc):
419 bOK = bshakef(fplog, constr->shaked,
420 md->invmass, constr->nblocks, constr->sblock,
421 idef, ir, x, min_proj, nrnb,
422 constr->lagr, lambda, dvdlambda,
423 invdt, NULL, vir != NULL, vir_r_m_dr,
424 constr->maxwarn >= 0, econq);
425 break;
426 default:
427 gmx_fatal(FARGS, "Internal error, SHAKE called for constraining something else than coordinates");
428 break;
431 if (!bOK && constr->maxwarn >= 0)
433 if (fplog != NULL)
435 fprintf(fplog, "Constraint error in algorithm %s at step %s\n",
436 econstr_names[econtSHAKE], gmx_step_str(step, buf));
438 bDump = TRUE;
442 if (nsettle > 0)
444 int calcvir_atom_end;
446 if (vir == NULL)
448 calcvir_atom_end = 0;
450 else
452 calcvir_atom_end = md->homenr;
455 switch (econq)
457 case econqCoord:
458 #pragma omp parallel for num_threads(nth) schedule(static)
459 for (th = 0; th < nth; th++)
463 int start_th, end_th;
465 if (th > 0)
467 clear_mat(constr->vir_r_m_dr_th[th]);
470 start_th = (nsettle* th )/nth;
471 end_th = (nsettle*(th+1))/nth;
472 if (start_th >= 0 && end_th - start_th > 0)
474 csettle(constr->settled,
475 end_th-start_th,
476 settle->iatoms+start_th*(1+NRAL(F_SETTLE)),
477 pbc_null,
478 x[0], xprime[0],
479 invdt, v ? v[0] : NULL, calcvir_atom_end,
480 th == 0 ? vir_r_m_dr : constr->vir_r_m_dr_th[th],
481 th == 0 ? &settle_error : &constr->settle_error[th]);
484 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
486 inc_nrnb(nrnb, eNR_SETTLE, nsettle);
487 if (v != NULL)
489 inc_nrnb(nrnb, eNR_CONSTR_V, nsettle*3);
491 if (vir != NULL)
493 inc_nrnb(nrnb, eNR_CONSTR_VIR, nsettle*3);
495 break;
496 case econqVeloc:
497 case econqDeriv:
498 case econqForce:
499 case econqForceDispl:
500 #pragma omp parallel for num_threads(nth) schedule(static)
501 for (th = 0; th < nth; th++)
505 int start_th, end_th;
507 if (th > 0)
509 clear_mat(constr->vir_r_m_dr_th[th]);
512 start_th = (nsettle* th )/nth;
513 end_th = (nsettle*(th+1))/nth;
515 if (start_th >= 0 && end_th - start_th > 0)
517 settle_proj(constr->settled, econq,
518 end_th-start_th,
519 settle->iatoms+start_th*(1+NRAL(F_SETTLE)),
520 pbc_null,
522 xprime, min_proj, calcvir_atom_end,
523 th == 0 ? vir_r_m_dr : constr->vir_r_m_dr_th[th]);
526 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
528 /* This is an overestimate */
529 inc_nrnb(nrnb, eNR_SETTLE, nsettle);
530 break;
531 case econqDeriv_FlexCon:
532 /* Nothing to do, since the are no flexible constraints in settles */
533 break;
534 default:
535 gmx_incons("Unknown constraint quantity for settle");
539 if (settle->nr > 0)
541 /* Combine virial and error info of the other threads */
542 for (i = 1; i < nth; i++)
544 settle_error = constr->settle_error[i];
546 if (vir != NULL)
548 for (i = 1; i < nth; i++)
550 m_add(vir_r_m_dr, constr->vir_r_m_dr_th[i], vir_r_m_dr);
554 if (econq == econqCoord && settle_error >= 0)
556 bOK = FALSE;
557 if (constr->maxwarn >= 0)
559 char buf[256];
560 sprintf(buf,
561 "\nstep " "%" GMX_PRId64 ": Water molecule starting at atom %d can not be "
562 "settled.\nCheck for bad contacts and/or reduce the timestep if appropriate.\n",
563 step, ddglatnr(cr->dd, settle->iatoms[settle_error*(1+NRAL(F_SETTLE))+1]));
564 if (fplog)
566 fprintf(fplog, "%s", buf);
568 fprintf(stderr, "%s", buf);
569 constr->warncount_settle++;
570 if (constr->warncount_settle > constr->maxwarn)
572 too_many_constraint_warnings(-1, constr->warncount_settle);
574 bDump = TRUE;
579 if (vir != NULL)
581 /* The normal uses of constrain() pass step_scaling = 1.0.
582 * The call to constrain() for SD1 that passes step_scaling =
583 * 0.5 also passes vir = NULL, so cannot reach this
584 * assertion. This assertion should remain until someone knows
585 * that this path works for their intended purpose, and then
586 * they can use scaled_delta_t instead of ir->delta_t
587 * below. */
588 assert(gmx_within_tol(step_scaling, 1.0, GMX_REAL_EPS));
589 switch (econq)
591 case econqCoord:
592 vir_fac = 0.5/(ir->delta_t*ir->delta_t);
593 break;
594 case econqVeloc:
595 vir_fac = 0.5/ir->delta_t;
596 break;
597 case econqForce:
598 case econqForceDispl:
599 vir_fac = 0.5;
600 break;
601 default:
602 gmx_incons("Unsupported constraint quantity for virial");
605 if (EI_VV(ir->eI))
607 vir_fac *= 2; /* only constraining over half the distance here */
609 for (i = 0; i < DIM; i++)
611 for (j = 0; j < DIM; j++)
613 (*vir)[i][j] = vir_fac*vir_r_m_dr[i][j];
618 if (bDump)
620 dump_confs(fplog, step, constr->warn_mtop, start, homenr, cr, x, xprime, box);
623 if (econq == econqCoord)
625 if (ir->bPull && pull_have_constraint(ir->pull_work))
627 if (EI_DYNAMICS(ir->eI))
629 t = ir->init_t + (step + delta_step)*ir->delta_t;
631 else
633 t = ir->init_t;
635 set_pbc(&pbc, ir->ePBC, box);
636 pull_constraint(ir->pull_work, md, &pbc, cr, ir->delta_t, t, x, xprime, v, *vir);
638 if (constr->ed && delta_step > 0)
640 /* apply the essential dynamcs constraints here */
641 do_edsam(ir, step, cr, xprime, v, box, constr->ed);
645 return bOK;
648 real *constr_rmsd_data(struct gmx_constr *constr)
650 if (constr->lincsd)
652 return lincs_rmsd_data(constr->lincsd);
654 else
656 return NULL;
660 real constr_rmsd(struct gmx_constr *constr, gmx_bool bSD2)
662 if (constr->lincsd)
664 return lincs_rmsd(constr->lincsd, bSD2);
666 else
668 return 0;
672 static void make_shake_sblock_serial(struct gmx_constr *constr,
673 t_idef *idef, const t_mdatoms *md)
675 int i, j, m, ncons;
676 int bstart, bnr;
677 t_blocka sblocks;
678 t_sortblock *sb;
679 t_iatom *iatom;
680 atom_id *inv_sblock;
682 /* Since we are processing the local topology,
683 * the F_CONSTRNC ilist has been concatenated to the F_CONSTR ilist.
685 ncons = idef->il[F_CONSTR].nr/3;
687 init_blocka(&sblocks);
688 gen_sblocks(NULL, 0, md->homenr, idef, &sblocks, FALSE);
691 bstart=(idef->nodeid > 0) ? blocks->multinr[idef->nodeid-1] : 0;
692 nblocks=blocks->multinr[idef->nodeid] - bstart;
694 bstart = 0;
695 constr->nblocks = sblocks.nr;
696 if (debug)
698 fprintf(debug, "ncons: %d, bstart: %d, nblocks: %d\n",
699 ncons, bstart, constr->nblocks);
702 /* Calculate block number for each atom */
703 inv_sblock = make_invblocka(&sblocks, md->nr);
705 done_blocka(&sblocks);
707 /* Store the block number in temp array and
708 * sort the constraints in order of the sblock number
709 * and the atom numbers, really sorting a segment of the array!
711 #ifdef DEBUGIDEF
712 pr_idef(fplog, 0, "Before Sort", idef);
713 #endif
714 iatom = idef->il[F_CONSTR].iatoms;
715 snew(sb, ncons);
716 for (i = 0; (i < ncons); i++, iatom += 3)
718 for (m = 0; (m < 3); m++)
720 sb[i].iatom[m] = iatom[m];
722 sb[i].blocknr = inv_sblock[iatom[1]];
725 /* Now sort the blocks */
726 if (debug)
728 pr_sortblock(debug, "Before sorting", ncons, sb);
729 fprintf(debug, "Going to sort constraints\n");
732 qsort(sb, ncons, (size_t)sizeof(*sb), pcomp);
734 if (debug)
736 pr_sortblock(debug, "After sorting", ncons, sb);
739 iatom = idef->il[F_CONSTR].iatoms;
740 for (i = 0; (i < ncons); i++, iatom += 3)
742 for (m = 0; (m < 3); m++)
744 iatom[m] = sb[i].iatom[m];
747 #ifdef DEBUGIDEF
748 pr_idef(fplog, 0, "After Sort", idef);
749 #endif
751 j = 0;
752 snew(constr->sblock, constr->nblocks+1);
753 bnr = -2;
754 for (i = 0; (i < ncons); i++)
756 if (sb[i].blocknr != bnr)
758 bnr = sb[i].blocknr;
759 constr->sblock[j++] = 3*i;
762 /* Last block... */
763 constr->sblock[j++] = 3*ncons;
765 if (j != (constr->nblocks+1))
767 fprintf(stderr, "bstart: %d\n", bstart);
768 fprintf(stderr, "j: %d, nblocks: %d, ncons: %d\n",
769 j, constr->nblocks, ncons);
770 for (i = 0; (i < ncons); i++)
772 fprintf(stderr, "i: %5d sb[i].blocknr: %5d\n", i, sb[i].blocknr);
774 for (j = 0; (j <= constr->nblocks); j++)
776 fprintf(stderr, "sblock[%3d]=%5d\n", j, (int)constr->sblock[j]);
778 gmx_fatal(FARGS, "DEATH HORROR: "
779 "sblocks does not match idef->il[F_CONSTR]");
781 sfree(sb);
782 sfree(inv_sblock);
785 static void make_shake_sblock_dd(struct gmx_constr *constr,
786 const t_ilist *ilcon, const t_block *cgs,
787 const gmx_domdec_t *dd)
789 int ncons, c, cg;
790 t_iatom *iatom;
792 if (dd->ncg_home+1 > constr->sblock_nalloc)
794 constr->sblock_nalloc = over_alloc_dd(dd->ncg_home+1);
795 srenew(constr->sblock, constr->sblock_nalloc);
798 ncons = ilcon->nr/3;
799 iatom = ilcon->iatoms;
800 constr->nblocks = 0;
801 cg = 0;
802 for (c = 0; c < ncons; c++)
804 if (c == 0 || iatom[1] >= cgs->index[cg+1])
806 constr->sblock[constr->nblocks++] = 3*c;
807 while (iatom[1] >= cgs->index[cg+1])
809 cg++;
812 iatom += 3;
814 constr->sblock[constr->nblocks] = 3*ncons;
817 t_blocka make_at2con(int start, int natoms,
818 const t_ilist *ilist, const t_iparams *iparams,
819 gmx_bool bDynamics, int *nflexiblecons)
821 int *count, ncon, con, con_tot, nflexcon, ftype, i, a;
822 t_iatom *ia;
823 t_blocka at2con;
824 gmx_bool bFlexCon;
826 snew(count, natoms);
827 nflexcon = 0;
828 for (ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
830 ncon = ilist[ftype].nr/3;
831 ia = ilist[ftype].iatoms;
832 for (con = 0; con < ncon; con++)
834 bFlexCon = (iparams[ia[0]].constr.dA == 0 &&
835 iparams[ia[0]].constr.dB == 0);
836 if (bFlexCon)
838 nflexcon++;
840 if (bDynamics || !bFlexCon)
842 for (i = 1; i < 3; i++)
844 a = ia[i] - start;
845 count[a]++;
848 ia += 3;
851 *nflexiblecons = nflexcon;
853 at2con.nr = natoms;
854 at2con.nalloc_index = at2con.nr+1;
855 snew(at2con.index, at2con.nalloc_index);
856 at2con.index[0] = 0;
857 for (a = 0; a < natoms; a++)
859 at2con.index[a+1] = at2con.index[a] + count[a];
860 count[a] = 0;
862 at2con.nra = at2con.index[natoms];
863 at2con.nalloc_a = at2con.nra;
864 snew(at2con.a, at2con.nalloc_a);
866 /* The F_CONSTRNC constraints have constraint numbers
867 * that continue after the last F_CONSTR constraint.
869 con_tot = 0;
870 for (ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
872 ncon = ilist[ftype].nr/3;
873 ia = ilist[ftype].iatoms;
874 for (con = 0; con < ncon; con++)
876 bFlexCon = (iparams[ia[0]].constr.dA == 0 &&
877 iparams[ia[0]].constr.dB == 0);
878 if (bDynamics || !bFlexCon)
880 for (i = 1; i < 3; i++)
882 a = ia[i] - start;
883 at2con.a[at2con.index[a]+count[a]++] = con_tot;
886 con_tot++;
887 ia += 3;
891 sfree(count);
893 return at2con;
896 static int *make_at2settle(int natoms, const t_ilist *ilist)
898 int *at2s;
899 int a, stride, s;
901 snew(at2s, natoms);
902 /* Set all to no settle */
903 for (a = 0; a < natoms; a++)
905 at2s[a] = -1;
908 stride = 1 + NRAL(F_SETTLE);
910 for (s = 0; s < ilist->nr; s += stride)
912 at2s[ilist->iatoms[s+1]] = s/stride;
913 at2s[ilist->iatoms[s+2]] = s/stride;
914 at2s[ilist->iatoms[s+3]] = s/stride;
917 return at2s;
920 void set_constraints(struct gmx_constr *constr,
921 gmx_localtop_t *top, const t_inputrec *ir,
922 const t_mdatoms *md, t_commrec *cr)
924 t_idef *idef;
925 int ncons;
926 const t_ilist *settle;
927 int iO, iH;
929 idef = &top->idef;
931 if (constr->ncon_tot > 0)
933 /* We are using the local topology,
934 * so there are only F_CONSTR constraints.
936 ncons = idef->il[F_CONSTR].nr/3;
938 /* With DD we might also need to call LINCS with ncons=0 for
939 * communicating coordinates to other nodes that do have constraints.
941 if (ir->eConstrAlg == econtLINCS)
943 set_lincs(idef, md, EI_DYNAMICS(ir->eI), cr, constr->lincsd);
945 if (ir->eConstrAlg == econtSHAKE)
947 if (cr->dd)
949 make_shake_sblock_dd(constr, &idef->il[F_CONSTR], &top->cgs, cr->dd);
951 else
953 make_shake_sblock_serial(constr, idef, md);
955 if (ncons > constr->lagr_nalloc)
957 constr->lagr_nalloc = over_alloc_dd(ncons);
958 srenew(constr->lagr, constr->lagr_nalloc);
963 if (idef->il[F_SETTLE].nr > 0 && constr->settled == NULL)
965 settle = &idef->il[F_SETTLE];
966 iO = settle->iatoms[1];
967 iH = settle->iatoms[2];
968 constr->settled =
969 settle_init(md->massT[iO], md->massT[iH],
970 md->invmass[iO], md->invmass[iH],
971 idef->iparams[settle->iatoms[0]].settle.doh,
972 idef->iparams[settle->iatoms[0]].settle.dhh);
975 /* Make a selection of the local atoms for essential dynamics */
976 if (constr->ed && cr->dd)
978 dd_make_local_ed_indices(cr->dd, constr->ed);
982 static void constr_recur(const t_blocka *at2con,
983 const t_ilist *ilist, const t_iparams *iparams,
984 gmx_bool bTopB,
985 int at, int depth, int nc, int *path,
986 real r0, real r1, real *r2max,
987 int *count)
989 int ncon1;
990 t_iatom *ia1, *ia2;
991 int c, con, a1;
992 gmx_bool bUse;
993 t_iatom *ia;
994 real len, rn0, rn1;
996 (*count)++;
998 ncon1 = ilist[F_CONSTR].nr/3;
999 ia1 = ilist[F_CONSTR].iatoms;
1000 ia2 = ilist[F_CONSTRNC].iatoms;
1002 /* Loop over all constraints connected to this atom */
1003 for (c = at2con->index[at]; c < at2con->index[at+1]; c++)
1005 con = at2con->a[c];
1006 /* Do not walk over already used constraints */
1007 bUse = TRUE;
1008 for (a1 = 0; a1 < depth; a1++)
1010 if (con == path[a1])
1012 bUse = FALSE;
1015 if (bUse)
1017 ia = constr_iatomptr(ncon1, ia1, ia2, con);
1018 /* Flexible constraints currently have length 0, which is incorrect */
1019 if (!bTopB)
1021 len = iparams[ia[0]].constr.dA;
1023 else
1025 len = iparams[ia[0]].constr.dB;
1027 /* In the worst case the bond directions alternate */
1028 if (nc % 2 == 0)
1030 rn0 = r0 + len;
1031 rn1 = r1;
1033 else
1035 rn0 = r0;
1036 rn1 = r1 + len;
1038 /* Assume angles of 120 degrees between all bonds */
1039 if (rn0*rn0 + rn1*rn1 + rn0*rn1 > *r2max)
1041 *r2max = rn0*rn0 + rn1*rn1 + r0*rn1;
1042 if (debug)
1044 fprintf(debug, "Found longer constraint distance: r0 %5.3f r1 %5.3f rmax %5.3f\n", rn0, rn1, sqrt(*r2max));
1045 for (a1 = 0; a1 < depth; a1++)
1047 fprintf(debug, " %d %5.3f",
1048 path[a1],
1049 iparams[constr_iatomptr(ncon1, ia1, ia2, con)[0]].constr.dA);
1051 fprintf(debug, " %d %5.3f\n", con, len);
1054 /* Limit the number of recursions to 1000*nc,
1055 * so a call does not take more than a second,
1056 * even for highly connected systems.
1058 if (depth + 1 < nc && *count < 1000*nc)
1060 if (ia[1] == at)
1062 a1 = ia[2];
1064 else
1066 a1 = ia[1];
1068 /* Recursion */
1069 path[depth] = con;
1070 constr_recur(at2con, ilist, iparams,
1071 bTopB, a1, depth+1, nc, path, rn0, rn1, r2max, count);
1072 path[depth] = -1;
1078 static real constr_r_max_moltype(const gmx_moltype_t *molt,
1079 const t_iparams *iparams,
1080 const t_inputrec *ir)
1082 int natoms, nflexcon, *path, at, count;
1084 t_blocka at2con;
1085 real r0, r1, r2maxA, r2maxB, rmax, lam0, lam1;
1087 if (molt->ilist[F_CONSTR].nr == 0 &&
1088 molt->ilist[F_CONSTRNC].nr == 0)
1090 return 0;
1093 natoms = molt->atoms.nr;
1095 at2con = make_at2con(0, natoms, molt->ilist, iparams,
1096 EI_DYNAMICS(ir->eI), &nflexcon);
1097 snew(path, 1+ir->nProjOrder);
1098 for (at = 0; at < 1+ir->nProjOrder; at++)
1100 path[at] = -1;
1103 r2maxA = 0;
1104 for (at = 0; at < natoms; at++)
1106 r0 = 0;
1107 r1 = 0;
1109 count = 0;
1110 constr_recur(&at2con, molt->ilist, iparams,
1111 FALSE, at, 0, 1+ir->nProjOrder, path, r0, r1, &r2maxA, &count);
1113 if (ir->efep == efepNO)
1115 rmax = sqrt(r2maxA);
1117 else
1119 r2maxB = 0;
1120 for (at = 0; at < natoms; at++)
1122 r0 = 0;
1123 r1 = 0;
1124 count = 0;
1125 constr_recur(&at2con, molt->ilist, iparams,
1126 TRUE, at, 0, 1+ir->nProjOrder, path, r0, r1, &r2maxB, &count);
1128 lam0 = ir->fepvals->init_lambda;
1129 if (EI_DYNAMICS(ir->eI))
1131 lam0 += ir->init_step*ir->fepvals->delta_lambda;
1133 rmax = (1 - lam0)*sqrt(r2maxA) + lam0*sqrt(r2maxB);
1134 if (EI_DYNAMICS(ir->eI))
1136 lam1 = ir->fepvals->init_lambda + (ir->init_step + ir->nsteps)*ir->fepvals->delta_lambda;
1137 rmax = std::max(rmax, (1 - lam1)*std::sqrt(r2maxA) + lam1*std::sqrt(r2maxB));
1141 done_blocka(&at2con);
1142 sfree(path);
1144 return rmax;
1147 real constr_r_max(FILE *fplog, const gmx_mtop_t *mtop, const t_inputrec *ir)
1149 int mt;
1150 real rmax;
1152 rmax = 0;
1153 for (mt = 0; mt < mtop->nmoltype; mt++)
1155 rmax = std::max(rmax,
1156 constr_r_max_moltype(&mtop->moltype[mt],
1157 mtop->ffparams.iparams, ir));
1160 if (fplog)
1162 fprintf(fplog, "Maximum distance for %d constraints, at 120 deg. angles, all-trans: %.3f nm\n", 1+ir->nProjOrder, rmax);
1165 return rmax;
1168 gmx_constr_t init_constraints(FILE *fplog,
1169 const gmx_mtop_t *mtop, const t_inputrec *ir,
1170 gmx_edsam_t ed, t_state *state,
1171 t_commrec *cr)
1173 int ncon, nset, nmol, settle_type, i, mt, nflexcon;
1174 struct gmx_constr *constr;
1175 char *env;
1176 t_ilist *ilist;
1177 gmx_mtop_ilistloop_t iloop;
1179 ncon =
1180 gmx_mtop_ftype_count(mtop, F_CONSTR) +
1181 gmx_mtop_ftype_count(mtop, F_CONSTRNC);
1182 nset = gmx_mtop_ftype_count(mtop, F_SETTLE);
1184 if (ncon+nset == 0 &&
1185 !(ir->bPull && pull_have_constraint(ir->pull_work)) &&
1186 ed == NULL)
1188 return NULL;
1191 snew(constr, 1);
1193 constr->ncon_tot = ncon;
1194 constr->nflexcon = 0;
1195 if (ncon > 0)
1197 constr->n_at2con_mt = mtop->nmoltype;
1198 snew(constr->at2con_mt, constr->n_at2con_mt);
1199 for (mt = 0; mt < mtop->nmoltype; mt++)
1201 constr->at2con_mt[mt] = make_at2con(0, mtop->moltype[mt].atoms.nr,
1202 mtop->moltype[mt].ilist,
1203 mtop->ffparams.iparams,
1204 EI_DYNAMICS(ir->eI), &nflexcon);
1205 for (i = 0; i < mtop->nmolblock; i++)
1207 if (mtop->molblock[i].type == mt)
1209 constr->nflexcon += mtop->molblock[i].nmol*nflexcon;
1214 if (constr->nflexcon > 0)
1216 if (fplog)
1218 fprintf(fplog, "There are %d flexible constraints\n",
1219 constr->nflexcon);
1220 if (ir->fc_stepsize == 0)
1222 fprintf(fplog, "\n"
1223 "WARNING: step size for flexible constraining = 0\n"
1224 " All flexible constraints will be rigid.\n"
1225 " Will try to keep all flexible constraints at their original length,\n"
1226 " but the lengths may exhibit some drift.\n\n");
1227 constr->nflexcon = 0;
1230 if (constr->nflexcon > 0)
1232 please_cite(fplog, "Hess2002");
1236 if (ir->eConstrAlg == econtLINCS)
1238 constr->lincsd = init_lincs(fplog, mtop,
1239 constr->nflexcon, constr->at2con_mt,
1240 DOMAINDECOMP(cr) && cr->dd->bInterCGcons,
1241 ir->nLincsIter, ir->nProjOrder);
1244 if (ir->eConstrAlg == econtSHAKE)
1246 if (DOMAINDECOMP(cr) && cr->dd->bInterCGcons)
1248 gmx_fatal(FARGS, "SHAKE is not supported with domain decomposition and constraint that cross charge group boundaries, use LINCS");
1250 if (constr->nflexcon)
1252 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");
1254 please_cite(fplog, "Ryckaert77a");
1255 if (ir->bShakeSOR)
1257 please_cite(fplog, "Barth95a");
1260 constr->shaked = shake_init();
1264 if (nset > 0)
1266 please_cite(fplog, "Miyamoto92a");
1268 constr->bInterCGsettles = inter_charge_group_settles(mtop);
1270 /* Check that we have only one settle type */
1271 settle_type = -1;
1272 iloop = gmx_mtop_ilistloop_init(mtop);
1273 while (gmx_mtop_ilistloop_next(iloop, &ilist, &nmol))
1275 for (i = 0; i < ilist[F_SETTLE].nr; i += 4)
1277 if (settle_type == -1)
1279 settle_type = ilist[F_SETTLE].iatoms[i];
1281 else if (ilist[F_SETTLE].iatoms[i] != settle_type)
1283 gmx_fatal(FARGS,
1284 "The [molecules] section of your topology specifies more than one block of\n"
1285 "a [moleculetype] with a [settles] block. Only one such is allowed. If you\n"
1286 "are trying to partition your solvent into different *groups* (e.g. for\n"
1287 "freezing, T-coupling, etc.) then you are using the wrong approach. Index\n"
1288 "files specify groups. Otherwise, you may wish to change the least-used\n"
1289 "block of molecules with SETTLE constraints into 3 normal constraints.");
1294 constr->n_at2settle_mt = mtop->nmoltype;
1295 snew(constr->at2settle_mt, constr->n_at2settle_mt);
1296 for (mt = 0; mt < mtop->nmoltype; mt++)
1298 constr->at2settle_mt[mt] =
1299 make_at2settle(mtop->moltype[mt].atoms.nr,
1300 &mtop->moltype[mt].ilist[F_SETTLE]);
1304 if ((ncon + nset) > 0 && ir->epc == epcMTTK)
1306 gmx_fatal(FARGS, "Constraints are not implemented with MTTK pressure control.");
1309 constr->maxwarn = 999;
1310 env = getenv("GMX_MAXCONSTRWARN");
1311 if (env)
1313 constr->maxwarn = 0;
1314 sscanf(env, "%8d", &constr->maxwarn);
1315 if (fplog)
1317 fprintf(fplog,
1318 "Setting the maximum number of constraint warnings to %d\n",
1319 constr->maxwarn);
1321 if (MASTER(cr))
1323 fprintf(stderr,
1324 "Setting the maximum number of constraint warnings to %d\n",
1325 constr->maxwarn);
1328 if (constr->maxwarn < 0 && fplog)
1330 fprintf(fplog, "maxwarn < 0, will not stop on constraint errors\n");
1332 constr->warncount_lincs = 0;
1333 constr->warncount_settle = 0;
1335 /* Initialize the essential dynamics sampling.
1336 * Put the pointer to the ED struct in constr */
1337 constr->ed = ed;
1338 if (ed != NULL || state->edsamstate.nED > 0)
1340 init_edsam(mtop, ir, cr, ed, state->x, state->box, &state->edsamstate);
1343 constr->warn_mtop = mtop;
1345 return constr;
1348 const t_blocka *atom2constraints_moltype(gmx_constr_t constr)
1350 return constr->at2con_mt;
1353 const int **atom2settle_moltype(gmx_constr_t constr)
1355 return (const int **)constr->at2settle_mt;
1359 gmx_bool inter_charge_group_constraints(const gmx_mtop_t *mtop)
1361 const gmx_moltype_t *molt;
1362 const t_block *cgs;
1363 const t_ilist *il;
1364 int mb;
1365 int *at2cg, cg, a, ftype, i;
1366 gmx_bool bInterCG;
1368 bInterCG = FALSE;
1369 for (mb = 0; mb < mtop->nmolblock && !bInterCG; mb++)
1371 molt = &mtop->moltype[mtop->molblock[mb].type];
1373 if (molt->ilist[F_CONSTR].nr > 0 ||
1374 molt->ilist[F_CONSTRNC].nr > 0 ||
1375 molt->ilist[F_SETTLE].nr > 0)
1377 cgs = &molt->cgs;
1378 snew(at2cg, molt->atoms.nr);
1379 for (cg = 0; cg < cgs->nr; cg++)
1381 for (a = cgs->index[cg]; a < cgs->index[cg+1]; a++)
1383 at2cg[a] = cg;
1387 for (ftype = F_CONSTR; ftype <= F_CONSTRNC; ftype++)
1389 il = &molt->ilist[ftype];
1390 for (i = 0; i < il->nr && !bInterCG; i += 1+NRAL(ftype))
1392 if (at2cg[il->iatoms[i+1]] != at2cg[il->iatoms[i+2]])
1394 bInterCG = TRUE;
1399 sfree(at2cg);
1403 return bInterCG;
1406 gmx_bool inter_charge_group_settles(const gmx_mtop_t *mtop)
1408 const gmx_moltype_t *molt;
1409 const t_block *cgs;
1410 const t_ilist *il;
1411 int mb;
1412 int *at2cg, cg, a, ftype, i;
1413 gmx_bool bInterCG;
1415 bInterCG = FALSE;
1416 for (mb = 0; mb < mtop->nmolblock && !bInterCG; mb++)
1418 molt = &mtop->moltype[mtop->molblock[mb].type];
1420 if (molt->ilist[F_SETTLE].nr > 0)
1422 cgs = &molt->cgs;
1423 snew(at2cg, molt->atoms.nr);
1424 for (cg = 0; cg < cgs->nr; cg++)
1426 for (a = cgs->index[cg]; a < cgs->index[cg+1]; a++)
1428 at2cg[a] = cg;
1432 for (ftype = F_SETTLE; ftype <= F_SETTLE; ftype++)
1434 il = &molt->ilist[ftype];
1435 for (i = 0; i < il->nr && !bInterCG; i += 1+NRAL(F_SETTLE))
1437 if (at2cg[il->iatoms[i+1]] != at2cg[il->iatoms[i+2]] ||
1438 at2cg[il->iatoms[i+1]] != at2cg[il->iatoms[i+3]])
1440 bInterCG = TRUE;
1445 sfree(at2cg);
1449 return bInterCG;
1452 /* helper functions for andersen temperature control, because the
1453 * gmx_constr construct is only defined in constr.c. Return the list
1454 * of blocks (get_sblock) and the number of blocks (get_nblocks). */
1456 extern int *get_sblock(struct gmx_constr *constr)
1458 return constr->sblock;
1461 extern int get_nblocks(struct gmx_constr *constr)
1463 return constr->nblocks;