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
36 /* This file is completely threadsafe - keep it that way! */
53 #include "mtop_util.h"
55 #include "gmx_omp_nthreads.h"
59 int b0; /* first constraint for this thread */
60 int b1; /* b1-1 is the last constraint for this thread */
61 int nind; /* number of indices */
62 int *ind; /* constraint index for updating atom data */
63 int nind_r; /* number of indices */
64 int *ind_r; /* constraint index for updating atom data */
65 int ind_nalloc; /* allocation size of ind and ind_r */
66 tensor vir_r_m_dr; /* temporary variable for virial calculation */
69 typedef struct gmx_lincsdata {
70 int ncg; /* the global number of constraints */
71 int ncg_flex; /* the global number of flexible constraints */
72 int ncg_triangle; /* the global number of constraints in triangles */
73 int nIter; /* the number of iterations */
74 int nOrder; /* the order of the matrix expansion */
75 int nc; /* the number of constraints */
76 int nc_alloc; /* the number we allocated memory for */
77 int ncc; /* the number of constraint connections */
78 int ncc_alloc; /* the number we allocated memory for */
79 real matlam; /* the FE lambda value used for filling blc and blmf */
80 real *bllen0; /* the reference distance in topology A */
81 real *ddist; /* the reference distance in top B - the r.d. in top A */
82 int *bla; /* the atom pairs involved in the constraints */
83 real *blc; /* 1/sqrt(invmass1 + invmass2) */
84 real *blc1; /* as blc, but with all masses 1 */
85 int *blnr; /* index into blbnb and blmf */
86 int *blbnb; /* list of constraint connections */
87 int ntriangle; /* the local number of constraints in triangles */
88 int *triangle; /* the list of triangle constraints */
89 int *tri_bits; /* the bits tell if the matrix element should be used */
90 int ncc_triangle; /* the number of constraint connections in triangles */
91 gmx_bool bCommIter; /* communicate before each LINCS interation */
92 real *blmf; /* matrix of mass factors for constraint connections */
93 real *blmf1; /* as blmf, but with all masses 1 */
94 real *bllen; /* the reference bond length */
95 int nth; /* The number of threads doing LINCS */
96 lincs_thread_t *th; /* LINCS thread division */
97 unsigned *atf; /* atom flags for thread parallelization */
98 int atf_nalloc; /* allocation size of atf */
99 /* arrays for temporary storage in the LINCS algorithm */
106 real *mlambda; /* the Lagrange multipliers * -1 */
107 /* storage for the constraint RMS relative deviation output */
111 real *lincs_rmsd_data(struct gmx_lincsdata *lincsd)
113 return lincsd->rmsd_data;
116 real lincs_rmsd(struct gmx_lincsdata *lincsd, gmx_bool bSD2)
118 if (lincsd->rmsd_data[0] > 0)
120 return sqrt(lincsd->rmsd_data[bSD2 ? 2 : 1]/lincsd->rmsd_data[0]);
128 /* Do a set of nrec LINCS matrix multiplications.
129 * This function will return with up to date thread-local
130 * constraint data, without an OpenMP barrier.
132 static void lincs_matrix_expand(const struct gmx_lincsdata *lincsd,
135 real *rhs1, real *rhs2, real *sol)
137 int nrec, rec, b, j, n, nr0, nr1;
139 int ntriangle, tb, bits;
140 const int *blnr = lincsd->blnr, *blbnb = lincsd->blbnb;
141 const int *triangle = lincsd->triangle, *tri_bits = lincsd->tri_bits;
143 ntriangle = lincsd->ntriangle;
144 nrec = lincsd->nOrder;
146 for (rec = 0; rec < nrec; rec++)
149 for (b = b0; b < b1; b++)
152 for (n = blnr[b]; n < blnr[b+1]; n++)
155 mvb = mvb + blcc[n]*rhs1[j];
158 sol[b] = sol[b] + mvb;
163 } /* nrec*(ncons+2*nrtot) flops */
167 /* Perform an extra nrec recursions for only the constraints
168 * involved in rigid triangles.
169 * In this way their accuracy should come close to those of the other
170 * constraints, since traingles of constraints can produce eigenvalues
171 * around 0.7, while the effective eigenvalue for bond constraints
172 * is around 0.4 (and 0.7*0.7=0.5).
174 /* We need to copy the temporary array, since only the elements
175 * for constraints involved in triangles are updated and then
176 * the pointers are swapped. This saving copying the whole arrary.
177 * We need barrier as other threads might still be reading from rhs2.
180 for (b = b0; b < b1; b++)
187 for (rec = 0; rec < nrec; rec++)
189 for (tb = 0; tb < ntriangle; tb++)
196 for (n = nr0; n < nr1; n++)
198 if (bits & (1<<(n-nr0)))
201 mvb = mvb + blcc[n]*rhs1[j];
205 sol[b] = sol[b] + mvb;
211 } /* flops count is missing here */
213 /* We need a barrier here as the calling routine will continue
214 * to operate on the thread-local constraints without barrier.
220 static void lincs_update_atoms_noind(int ncons, const int *bla,
222 const real *fac, rvec *r,
227 real mvb, im1, im2, tmp0, tmp1, tmp2;
231 for (b = 0; b < ncons; b++)
247 } /* 16 ncons flops */
251 for (b = 0; b < ncons; b++)
269 static void lincs_update_atoms_ind(int ncons, const int *ind, const int *bla,
271 const real *fac, rvec *r,
276 real mvb, im1, im2, tmp0, tmp1, tmp2;
280 for (bi = 0; bi < ncons; bi++)
297 } /* 16 ncons flops */
301 for (bi = 0; bi < ncons; bi++)
316 } /* 16 ncons flops */
320 static void lincs_update_atoms(struct gmx_lincsdata *li, int th,
322 const real *fac, rvec *r,
328 /* Single thread, we simply update for all constraints */
329 lincs_update_atoms_noind(li->nc, li->bla, prefac, fac, r, invmass, x);
333 /* Update the atom vector components for our thread local
334 * constraints that only access our local atom range.
335 * This can be done without a barrier.
337 lincs_update_atoms_ind(li->th[th].nind, li->th[th].ind,
338 li->bla, prefac, fac, r, invmass, x);
340 if (li->th[li->nth].nind > 0)
342 /* Update the constraints that operate on atoms
343 * in multiple thread atom blocks on the master thread.
348 lincs_update_atoms_ind(li->th[li->nth].nind,
350 li->bla, prefac, fac, r, invmass, x);
356 /* LINCS projection, works on derivatives of the coordinates */
357 static void do_lincsp(rvec *x, rvec *f, rvec *fp, t_pbc *pbc,
358 struct gmx_lincsdata *lincsd, int th,
360 int econq, real *dvdlambda,
361 gmx_bool bCalcVir, tensor rmdf)
363 int b0, b1, b, i, j, k, n;
364 real tmp0, tmp1, tmp2, im1, im2, mvb, rlen, len, wfac, lam;
366 int *bla, *blnr, *blbnb;
368 real *blc, *blmf, *blcc, *rhs1, *rhs2, *sol;
370 b0 = lincsd->th[th].b0;
371 b1 = lincsd->th[th].b1;
376 blbnb = lincsd->blbnb;
377 if (econq != econqForce)
379 /* Use mass-weighted parameters */
385 /* Use non mass-weighted parameters */
387 blmf = lincsd->blmf1;
389 blcc = lincsd->tmpncc;
394 /* Compute normalized i-j vectors */
397 for (b = b0; b < b1; b++)
399 pbc_dx_aiuc(pbc, x[bla[2*b]], x[bla[2*b+1]], dx);
405 for (b = b0; b < b1; b++)
407 rvec_sub(x[bla[2*b]], x[bla[2*b+1]], dx);
409 } /* 16 ncons flops */
413 for (b = b0; b < b1; b++)
420 for (n = blnr[b]; n < blnr[b+1]; n++)
423 blcc[n] = blmf[n]*(tmp0*r[k][0] + tmp1*r[k][1] + tmp2*r[k][2]);
425 mvb = blc[b]*(tmp0*(f[i][0] - f[j][0]) +
426 tmp1*(f[i][1] - f[j][1]) +
427 tmp2*(f[i][2] - f[j][2]));
432 /* Together: 23*ncons + 6*nrtot flops */
434 lincs_matrix_expand(lincsd, b0, b1, blcc, rhs1, rhs2, sol);
435 /* nrec*(ncons+2*nrtot) flops */
437 if (econq == econqDeriv_FlexCon)
439 /* We only want to constraint the flexible constraints,
440 * so we mask out the normal ones by setting sol to 0.
442 for (b = b0; b < b1; b++)
444 if (!(lincsd->bllen0[b] == 0 && lincsd->ddist[b] == 0))
451 /* We multiply sol by blc, so we can use lincs_update_atoms for OpenMP */
452 for (b = b0; b < b1; b++)
457 /* When constraining forces, we should not use mass weighting,
458 * so we pass invmass=NULL, which results in the use of 1 for all atoms.
460 lincs_update_atoms(lincsd, th, 1.0, sol, r,
461 (econq != econqForce) ? invmass : NULL, fp);
463 if (dvdlambda != NULL)
466 for (b = b0; b < b1; b++)
468 *dvdlambda -= sol[b]*lincsd->ddist[b];
475 /* Constraint virial,
476 * determines sum r_bond x delta f,
477 * where delta f is the constraint correction
478 * of the quantity that is being constrained.
480 for (b = b0; b < b1; b++)
482 mvb = lincsd->bllen[b]*sol[b];
483 for (i = 0; i < DIM; i++)
486 for (j = 0; j < DIM; j++)
488 rmdf[i][j] += tmp1*r[b][j];
491 } /* 23 ncons flops */
495 static void do_lincs(rvec *x, rvec *xp, matrix box, t_pbc *pbc,
496 struct gmx_lincsdata *lincsd, int th,
499 gmx_bool bCalcLambda,
500 real wangle, int *warn,
502 gmx_bool bCalcVir, tensor vir_r_m_dr)
504 int b0, b1, b, i, j, k, n, iter;
505 real tmp0, tmp1, tmp2, im1, im2, mvb, rlen, len, len2, dlen2, wfac;
507 int *bla, *blnr, *blbnb;
509 real *blc, *blmf, *bllen, *blcc, *rhs1, *rhs2, *sol, *blc_sol, *mlambda;
512 b0 = lincsd->th[th].b0;
513 b1 = lincsd->th[th].b1;
518 blbnb = lincsd->blbnb;
521 bllen = lincsd->bllen;
522 blcc = lincsd->tmpncc;
526 blc_sol = lincsd->tmp4;
527 mlambda = lincsd->mlambda;
529 if (DOMAINDECOMP(cr) && cr->dd->constraints)
531 nlocat = dd_constraints_nlocalatoms(cr->dd);
533 else if (PARTDECOMP(cr))
535 nlocat = pd_constraints_nlocalatoms(cr->pd);
544 /* Compute normalized i-j vectors */
545 for (b = b0; b < b1; b++)
547 pbc_dx_aiuc(pbc, x[bla[2*b]], x[bla[2*b+1]], dx);
551 for (b = b0; b < b1; b++)
553 for (n = blnr[b]; n < blnr[b+1]; n++)
555 blcc[n] = blmf[n]*iprod(r[b], r[blbnb[n]]);
557 pbc_dx_aiuc(pbc, xp[bla[2*b]], xp[bla[2*b+1]], dx);
558 mvb = blc[b]*(iprod(r[b], dx) - bllen[b]);
565 /* Compute normalized i-j vectors */
566 for (b = b0; b < b1; b++)
570 tmp0 = x[i][0] - x[j][0];
571 tmp1 = x[i][1] - x[j][1];
572 tmp2 = x[i][2] - x[j][2];
573 rlen = gmx_invsqrt(tmp0*tmp0+tmp1*tmp1+tmp2*tmp2);
577 } /* 16 ncons flops */
580 for (b = b0; b < b1; b++)
588 for (n = blnr[b]; n < blnr[b+1]; n++)
591 blcc[n] = blmf[n]*(tmp0*r[k][0] + tmp1*r[k][1] + tmp2*r[k][2]);
593 mvb = blc[b]*(tmp0*(xp[i][0] - xp[j][0]) +
594 tmp1*(xp[i][1] - xp[j][1]) +
595 tmp2*(xp[i][2] - xp[j][2]) - len);
600 /* Together: 26*ncons + 6*nrtot flops */
603 lincs_matrix_expand(lincsd, b0, b1, blcc, rhs1, rhs2, sol);
604 /* nrec*(ncons+2*nrtot) flops */
606 for (b = b0; b < b1; b++)
608 mlambda[b] = blc[b]*sol[b];
611 /* Update the coordinates */
612 lincs_update_atoms(lincsd, th, 1.0, mlambda, r, invmass, xp);
615 ******** Correction for centripetal effects ********
618 wfac = cos(DEG2RAD*wangle);
621 for (iter = 0; iter < lincsd->nIter; iter++)
623 if ((lincsd->bCommIter && DOMAINDECOMP(cr) && cr->dd->constraints) ||
629 /* Communicate the corrected non-local coordinates */
630 if (DOMAINDECOMP(cr))
632 dd_move_x_constraints(cr->dd, box, xp, NULL);
636 pd_move_x_constraints(cr, xp, NULL);
642 for (b = b0; b < b1; b++)
647 pbc_dx_aiuc(pbc, xp[bla[2*b]], xp[bla[2*b+1]], dx);
651 rvec_sub(xp[bla[2*b]], xp[bla[2*b+1]], dx);
654 dlen2 = 2*len2 - norm2(dx);
655 if (dlen2 < wfac*len2 && (nlocat == NULL || nlocat[b]))
661 mvb = blc[b]*(len - dlen2*gmx_invsqrt(dlen2));
669 } /* 20*ncons flops */
671 lincs_matrix_expand(lincsd, b0, b1, blcc, rhs1, rhs2, sol);
672 /* nrec*(ncons+2*nrtot) flops */
674 for (b = b0; b < b1; b++)
681 /* Update the coordinates */
682 lincs_update_atoms(lincsd, th, 1.0, blc_sol, r, invmass, xp);
684 /* nit*ncons*(37+9*nrec) flops */
688 /* Update the velocities */
689 lincs_update_atoms(lincsd, th, invdt, mlambda, r, invmass, v);
693 if (nlocat != NULL && bCalcLambda)
695 /* In lincs_update_atoms thread might cross-read mlambda */
698 /* Only account for local atoms */
699 for (b = b0; b < b1; b++)
701 mlambda[b] *= 0.5*nlocat[b];
707 /* Constraint virial */
708 for (b = b0; b < b1; b++)
710 tmp0 = -bllen[b]*mlambda[b];
711 for (i = 0; i < DIM; i++)
714 for (j = 0; j < DIM; j++)
716 vir_r_m_dr[i][j] -= tmp1*r[b][j];
719 } /* 22 ncons flops */
723 * 26*ncons + 6*nrtot + nrec*(ncons+2*nrtot)
724 * + nit * (20*ncons + nrec*(ncons+2*nrtot) + 17 ncons)
726 * (26+nrec)*ncons + (6+2*nrec)*nrtot
727 * + nit * ((37+nrec)*ncons + 2*nrec*nrtot)
729 * (63+nrec)*ncons + (6+4*nrec)*nrtot
733 void set_lincs_matrix(struct gmx_lincsdata *li, real *invmass, real lambda)
735 int i, a1, a2, n, k, sign, center;
737 const real invsqrt2 = 0.7071067811865475244;
739 for (i = 0; (i < li->nc); i++)
743 li->blc[i] = gmx_invsqrt(invmass[a1] + invmass[a2]);
744 li->blc1[i] = invsqrt2;
747 /* Construct the coupling coefficient matrix blmf */
749 li->ncc_triangle = 0;
750 for (i = 0; (i < li->nc); i++)
754 for (n = li->blnr[i]; (n < li->blnr[i+1]); n++)
757 if (a1 == li->bla[2*k] || a2 == li->bla[2*k+1])
765 if (a1 == li->bla[2*k] || a1 == li->bla[2*k+1])
775 li->blmf[n] = sign*invmass[center]*li->blc[i]*li->blc[k];
776 li->blmf1[n] = sign*0.5;
777 if (li->ncg_triangle > 0)
779 /* Look for constraint triangles */
780 for (nk = li->blnr[k]; (nk < li->blnr[k+1]); nk++)
783 if (kk != i && kk != k &&
784 (li->bla[2*kk] == end || li->bla[2*kk+1] == end))
786 if (li->ntriangle == 0 ||
787 li->triangle[li->ntriangle-1] < i)
789 /* Add this constraint to the triangle list */
790 li->triangle[li->ntriangle] = i;
791 li->tri_bits[li->ntriangle] = 0;
793 if (li->blnr[i+1] - li->blnr[i] > sizeof(li->tri_bits[0])*8 - 1)
795 gmx_fatal(FARGS, "A constraint is connected to %d constraints, this is more than the %d allowed for constraints participating in triangles",
796 li->blnr[i+1] - li->blnr[i],
797 sizeof(li->tri_bits[0])*8-1);
800 li->tri_bits[li->ntriangle-1] |= (1<<(n-li->blnr[i]));
810 fprintf(debug, "Of the %d constraints %d participate in triangles\n",
811 li->nc, li->ntriangle);
812 fprintf(debug, "There are %d couplings of which %d in triangles\n",
813 li->ncc, li->ncc_triangle);
817 * so we know with which lambda value the masses have been set.
822 static int count_triangle_constraints(t_ilist *ilist, t_blocka *at2con)
825 int c0, a00, a01, n1, c1, a10, a11, ac1, n2, c2, a20, a21;
828 t_iatom *ia1, *ia2, *iap;
830 ncon1 = ilist[F_CONSTR].nr/3;
831 ncon_tot = ncon1 + ilist[F_CONSTRNC].nr/3;
833 ia1 = ilist[F_CONSTR].iatoms;
834 ia2 = ilist[F_CONSTRNC].iatoms;
837 for (c0 = 0; c0 < ncon_tot; c0++)
840 iap = constr_iatomptr(ncon1, ia1, ia2, c0);
843 for (n1 = at2con->index[a01]; n1 < at2con->index[a01+1]; n1++)
848 iap = constr_iatomptr(ncon1, ia1, ia2, c1);
859 for (n2 = at2con->index[ac1]; n2 < at2con->index[ac1+1]; n2++)
862 if (c2 != c0 && c2 != c1)
864 iap = constr_iatomptr(ncon1, ia1, ia2, c2);
867 if (a20 == a00 || a21 == a00)
881 return ncon_triangle;
884 static gmx_bool more_than_two_sequential_constraints(const t_ilist *ilist,
885 const t_blocka *at2con)
887 t_iatom *ia1, *ia2, *iap;
888 int ncon1, ncon_tot, c;
890 gmx_bool bMoreThanTwoSequentialConstraints;
892 ncon1 = ilist[F_CONSTR].nr/3;
893 ncon_tot = ncon1 + ilist[F_CONSTRNC].nr/3;
895 ia1 = ilist[F_CONSTR].iatoms;
896 ia2 = ilist[F_CONSTRNC].iatoms;
898 bMoreThanTwoSequentialConstraints = FALSE;
899 for (c = 0; c < ncon_tot && !bMoreThanTwoSequentialConstraints; c++)
901 iap = constr_iatomptr(ncon1, ia1, ia2, c);
904 /* Check if this constraint has constraints connected at both atoms */
905 if (at2con->index[a1+1] - at2con->index[a1] > 1 &&
906 at2con->index[a2+1] - at2con->index[a2] > 1)
908 bMoreThanTwoSequentialConstraints = TRUE;
912 return bMoreThanTwoSequentialConstraints;
915 static int int_comp(const void *a, const void *b)
917 return (*(int *)a) - (*(int *)b);
920 gmx_lincsdata_t init_lincs(FILE *fplog, gmx_mtop_t *mtop,
921 int nflexcon_global, t_blocka *at2con,
922 gmx_bool bPLINCS, int nIter, int nProjOrder)
924 struct gmx_lincsdata *li;
930 fprintf(fplog, "\nInitializing%s LINear Constraint Solver\n",
931 bPLINCS ? " Parallel" : "");
937 gmx_mtop_ftype_count(mtop, F_CONSTR) +
938 gmx_mtop_ftype_count(mtop, F_CONSTRNC);
939 li->ncg_flex = nflexcon_global;
942 li->nOrder = nProjOrder;
944 li->ncg_triangle = 0;
945 li->bCommIter = FALSE;
946 for (mb = 0; mb < mtop->nmolblock; mb++)
948 molt = &mtop->moltype[mtop->molblock[mb].type];
950 mtop->molblock[mb].nmol*
951 count_triangle_constraints(molt->ilist,
952 &at2con[mtop->molblock[mb].type]);
953 if (bPLINCS && li->bCommIter == FALSE)
955 /* Check if we need to communicate not only before LINCS,
956 * but also before each iteration.
957 * The check for only two sequential constraints is only
958 * useful for the common case of H-bond only constraints.
959 * With more effort we could also make it useful for small
960 * molecules with nr. sequential constraints <= nOrder-1.
962 li->bCommIter = (li->nOrder < 1 || more_than_two_sequential_constraints(molt->ilist, &at2con[mtop->molblock[mb].type]));
965 if (debug && bPLINCS)
967 fprintf(debug, "PLINCS communication before each iteration: %d\n",
971 /* LINCS can run on any number of threads.
972 * Currently the number is fixed for the whole simulation,
973 * but it could be set in set_lincs().
975 li->nth = gmx_omp_nthreads_get(emntLINCS);
982 /* Allocate an extra elements for "thread-overlap" constraints */
983 snew(li->th, li->nth+1);
987 fprintf(debug, "LINCS: using %d threads\n", li->nth);
990 if (bPLINCS || li->ncg_triangle > 0)
992 please_cite(fplog, "Hess2008a");
996 please_cite(fplog, "Hess97a");
1001 fprintf(fplog, "The number of constraints is %d\n", li->ncg);
1004 fprintf(fplog, "There are inter charge-group constraints,\n"
1005 "will communicate selected coordinates each lincs iteration\n");
1007 if (li->ncg_triangle > 0)
1010 "%d constraints are involved in constraint triangles,\n"
1011 "will apply an additional matrix expansion of order %d for couplings\n"
1012 "between constraints inside triangles\n",
1013 li->ncg_triangle, li->nOrder);
1020 /* Sets up the work division over the threads */
1021 static void lincs_thread_setup(struct gmx_lincsdata *li, int natoms)
1023 lincs_thread_t *li_m;
1028 if (natoms > li->atf_nalloc)
1030 li->atf_nalloc = over_alloc_large(natoms);
1031 srenew(li->atf, li->atf_nalloc);
1035 /* Clear the atom flags */
1036 for (a = 0; a < natoms; a++)
1041 for (th = 0; th < li->nth; th++)
1043 lincs_thread_t *li_th;
1046 li_th = &li->th[th];
1048 /* The constraints are divided equally over the threads */
1049 li_th->b0 = (li->nc* th )/li->nth;
1050 li_th->b1 = (li->nc*(th+1))/li->nth;
1052 if (th < sizeof(*atf)*8)
1054 /* For each atom set a flag for constraints from each */
1055 for (b = li_th->b0; b < li_th->b1; b++)
1057 atf[li->bla[b*2] ] |= (1U<<th);
1058 atf[li->bla[b*2+1]] |= (1U<<th);
1063 #pragma omp parallel for num_threads(li->nth) schedule(static)
1064 for (th = 0; th < li->nth; th++)
1066 lincs_thread_t *li_th;
1070 li_th = &li->th[th];
1072 if (li_th->b1 - li_th->b0 > li_th->ind_nalloc)
1074 li_th->ind_nalloc = over_alloc_large(li_th->b1-li_th->b0);
1075 srenew(li_th->ind, li_th->ind_nalloc);
1076 srenew(li_th->ind_r, li_th->ind_nalloc);
1079 if (th < sizeof(*atf)*8)
1081 mask = (1U<<th) - 1U;
1085 for (b = li_th->b0; b < li_th->b1; b++)
1087 /* We let the constraint with the lowest thread index
1088 * operate on atoms with constraints from multiple threads.
1090 if (((atf[li->bla[b*2]] & mask) == 0) &&
1091 ((atf[li->bla[b*2+1]] & mask) == 0))
1093 /* Add the constraint to the local atom update index */
1094 li_th->ind[li_th->nind++] = b;
1098 /* Add the constraint to the rest block */
1099 li_th->ind_r[li_th->nind_r++] = b;
1105 /* We are out of bits, assign all constraints to rest */
1106 for (b = li_th->b0; b < li_th->b1; b++)
1108 li_th->ind_r[li_th->nind_r++] = b;
1113 /* We need to copy all constraints which have not be assigned
1114 * to a thread to a separate list which will be handled by one thread.
1116 li_m = &li->th[li->nth];
1119 for (th = 0; th < li->nth; th++)
1121 lincs_thread_t *li_th;
1124 li_th = &li->th[th];
1126 if (li_m->nind + li_th->nind_r > li_m->ind_nalloc)
1128 li_m->ind_nalloc = over_alloc_large(li_m->nind+li_th->nind_r);
1129 srenew(li_m->ind, li_m->ind_nalloc);
1132 for (b = 0; b < li_th->nind_r; b++)
1134 li_m->ind[li_m->nind++] = li_th->ind_r[b];
1139 fprintf(debug, "LINCS thread %d: %d constraints\n",
1146 fprintf(debug, "LINCS thread r: %d constraints\n",
1152 void set_lincs(t_idef *idef, t_mdatoms *md,
1153 gmx_bool bDynamics, t_commrec *cr,
1154 struct gmx_lincsdata *li)
1156 int start, natoms, nflexcon;
1159 int i, k, ncc_alloc, ni, con, nconnect, concon;
1161 real lenA = 0, lenB;
1166 /* Zero the thread index ranges.
1167 * Otherwise without local constraints we could return with old ranges.
1169 for (i = 0; i < li->nth; i++)
1177 li->th[li->nth].nind = 0;
1180 /* This is the local topology, so there are only F_CONSTR constraints */
1181 if (idef->il[F_CONSTR].nr == 0)
1183 /* There are no constraints,
1184 * we do not need to fill any data structures.
1191 fprintf(debug, "Building the LINCS connectivity\n");
1194 if (DOMAINDECOMP(cr))
1196 if (cr->dd->constraints)
1198 dd_get_constraint_range(cr->dd, &start, &natoms);
1202 natoms = cr->dd->nat_home;
1206 else if (PARTDECOMP(cr))
1208 pd_get_constraint_range(cr->pd, &start, &natoms);
1213 natoms = md->homenr;
1215 at2con = make_at2con(start, natoms, idef->il, idef->iparams, bDynamics,
1219 if (idef->il[F_CONSTR].nr/3 > li->nc_alloc || li->nc_alloc == 0)
1221 li->nc_alloc = over_alloc_dd(idef->il[F_CONSTR].nr/3);
1222 srenew(li->bllen0, li->nc_alloc);
1223 srenew(li->ddist, li->nc_alloc);
1224 srenew(li->bla, 2*li->nc_alloc);
1225 srenew(li->blc, li->nc_alloc);
1226 srenew(li->blc1, li->nc_alloc);
1227 srenew(li->blnr, li->nc_alloc+1);
1228 srenew(li->bllen, li->nc_alloc);
1229 srenew(li->tmpv, li->nc_alloc);
1230 srenew(li->tmp1, li->nc_alloc);
1231 srenew(li->tmp2, li->nc_alloc);
1232 srenew(li->tmp3, li->nc_alloc);
1233 srenew(li->tmp4, li->nc_alloc);
1234 srenew(li->mlambda, li->nc_alloc);
1235 if (li->ncg_triangle > 0)
1237 /* This is allocating too much, but it is difficult to improve */
1238 srenew(li->triangle, li->nc_alloc);
1239 srenew(li->tri_bits, li->nc_alloc);
1243 iatom = idef->il[F_CONSTR].iatoms;
1245 ncc_alloc = li->ncc_alloc;
1248 ni = idef->il[F_CONSTR].nr/3;
1252 li->blnr[con] = nconnect;
1253 for (i = 0; i < ni; i++)
1259 lenA = idef->iparams[type].constr.dA;
1260 lenB = idef->iparams[type].constr.dB;
1261 /* Skip the flexible constraints when not doing dynamics */
1262 if (bDynamics || lenA != 0 || lenB != 0)
1264 li->bllen0[con] = lenA;
1265 li->ddist[con] = lenB - lenA;
1266 /* Set the length to the topology A length */
1267 li->bllen[con] = li->bllen0[con];
1268 li->bla[2*con] = a1;
1269 li->bla[2*con+1] = a2;
1270 /* Construct the constraint connection matrix blbnb */
1271 for (k = at2con.index[a1-start]; k < at2con.index[a1-start+1]; k++)
1273 concon = at2con.a[k];
1276 if (nconnect >= ncc_alloc)
1278 ncc_alloc = over_alloc_small(nconnect+1);
1279 srenew(li->blbnb, ncc_alloc);
1281 li->blbnb[nconnect++] = concon;
1284 for (k = at2con.index[a2-start]; k < at2con.index[a2-start+1]; k++)
1286 concon = at2con.a[k];
1289 if (nconnect+1 > ncc_alloc)
1291 ncc_alloc = over_alloc_small(nconnect+1);
1292 srenew(li->blbnb, ncc_alloc);
1294 li->blbnb[nconnect++] = concon;
1297 li->blnr[con+1] = nconnect;
1301 /* Order the blbnb matrix to optimize memory access */
1302 qsort(&(li->blbnb[li->blnr[con]]), li->blnr[con+1]-li->blnr[con],
1303 sizeof(li->blbnb[0]), int_comp);
1305 /* Increase the constraint count */
1310 done_blocka(&at2con);
1312 /* This is the real number of constraints,
1313 * without dynamics the flexible constraints are not present.
1317 li->ncc = li->blnr[con];
1320 /* Since the matrix is static, we can free some memory */
1321 ncc_alloc = li->ncc;
1322 srenew(li->blbnb, ncc_alloc);
1325 if (ncc_alloc > li->ncc_alloc)
1327 li->ncc_alloc = ncc_alloc;
1328 srenew(li->blmf, li->ncc_alloc);
1329 srenew(li->blmf1, li->ncc_alloc);
1330 srenew(li->tmpncc, li->ncc_alloc);
1335 fprintf(debug, "Number of constraints is %d, couplings %d\n",
1342 li->th[0].b1 = li->nc;
1346 lincs_thread_setup(li, md->nr);
1349 set_lincs_matrix(li, md->invmass, md->lambda);
1352 static void lincs_warning(FILE *fplog,
1353 gmx_domdec_t *dd, rvec *x, rvec *xprime, t_pbc *pbc,
1354 int ncons, int *bla, real *bllen, real wangle,
1355 int maxwarn, int *warncount)
1359 real wfac, d0, d1, cosine;
1362 wfac = cos(DEG2RAD*wangle);
1364 sprintf(buf, "bonds that rotated more than %g degrees:\n"
1365 " atom 1 atom 2 angle previous, current, constraint length\n",
1367 fprintf(stderr, "%s", buf);
1370 fprintf(fplog, "%s", buf);
1373 for (b = 0; b < ncons; b++)
1379 pbc_dx_aiuc(pbc, x[i], x[j], v0);
1380 pbc_dx_aiuc(pbc, xprime[i], xprime[j], v1);
1384 rvec_sub(x[i], x[j], v0);
1385 rvec_sub(xprime[i], xprime[j], v1);
1389 cosine = iprod(v0, v1)/(d0*d1);
1392 sprintf(buf, " %6d %6d %5.1f %8.4f %8.4f %8.4f\n",
1393 ddglatnr(dd, i), ddglatnr(dd, j),
1394 RAD2DEG*acos(cosine), d0, d1, bllen[b]);
1395 fprintf(stderr, "%s", buf);
1398 fprintf(fplog, "%s", buf);
1400 if (!gmx_isfinite(d1))
1402 gmx_fatal(FARGS, "Bond length not finite.");
1408 if (*warncount > maxwarn)
1410 too_many_constraint_warnings(econtLINCS, *warncount);
1414 static void cconerr(gmx_domdec_t *dd,
1415 int ncons, int *bla, real *bllen, rvec *x, t_pbc *pbc,
1416 real *ncons_loc, real *ssd, real *max, int *imax)
1418 real len, d, ma, ssd2, r2;
1419 int *nlocat, count, b, im;
1422 if (dd && dd->constraints)
1424 nlocat = dd_constraints_nlocalatoms(dd);
1435 for (b = 0; b < ncons; b++)
1439 pbc_dx_aiuc(pbc, x[bla[2*b]], x[bla[2*b+1]], dx);
1443 rvec_sub(x[bla[2*b]], x[bla[2*b+1]], dx);
1446 len = r2*gmx_invsqrt(r2);
1447 d = fabs(len/bllen[b]-1);
1448 if (d > ma && (nlocat == NULL || nlocat[b]))
1460 ssd2 += nlocat[b]*d*d;
1465 *ncons_loc = (nlocat ? 0.5 : 1)*count;
1466 *ssd = (nlocat ? 0.5 : 1)*ssd2;
1471 gmx_bool constrain_lincs(FILE *fplog, gmx_bool bLog, gmx_bool bEner,
1473 gmx_large_int_t step,
1474 struct gmx_lincsdata *lincsd, t_mdatoms *md,
1476 rvec *x, rvec *xprime, rvec *min_proj,
1477 matrix box, t_pbc *pbc,
1478 real lambda, real *dvdlambda,
1479 real invdt, rvec *v,
1480 gmx_bool bCalcVir, tensor vir_r_m_dr,
1483 int maxwarn, int *warncount)
1485 char buf[STRLEN], buf2[22], buf3[STRLEN];
1486 int i, warn, p_imax, error;
1487 real ncons_loc, p_ssd, p_max = 0;
1493 if (lincsd->nc == 0 && cr->dd == NULL)
1497 lincsd->rmsd_data[0] = 0;
1498 if (ir->eI == eiSD2 && v == NULL)
1506 lincsd->rmsd_data[i] = 0;
1512 if (econq == econqCoord)
1514 if (ir->efep != efepNO)
1516 if (md->nMassPerturbed && lincsd->matlam != md->lambda)
1518 set_lincs_matrix(lincsd, md->invmass, md->lambda);
1521 for (i = 0; i < lincsd->nc; i++)
1523 lincsd->bllen[i] = lincsd->bllen0[i] + lambda*lincsd->ddist[i];
1527 if (lincsd->ncg_flex)
1529 /* Set the flexible constraint lengths to the old lengths */
1532 for (i = 0; i < lincsd->nc; i++)
1534 if (lincsd->bllen[i] == 0)
1536 pbc_dx_aiuc(pbc, x[lincsd->bla[2*i]], x[lincsd->bla[2*i+1]], dx);
1537 lincsd->bllen[i] = norm(dx);
1543 for (i = 0; i < lincsd->nc; i++)
1545 if (lincsd->bllen[i] == 0)
1548 sqrt(distance2(x[lincsd->bla[2*i]],
1549 x[lincsd->bla[2*i+1]]));
1557 cconerr(cr->dd, lincsd->nc, lincsd->bla, lincsd->bllen, xprime, pbc,
1558 &ncons_loc, &p_ssd, &p_max, &p_imax);
1561 /* This warn var can be updated by multiple threads
1562 * at the same time. But as we only need to detect
1563 * if a warning occured or not, this is not an issue.
1567 /* The OpenMP parallel region of constrain_lincs for coords */
1568 #pragma omp parallel num_threads(lincsd->nth)
1570 int th = gmx_omp_get_thread_num();
1572 clear_mat(lincsd->th[th].vir_r_m_dr);
1574 do_lincs(x, xprime, box, pbc, lincsd, th,
1576 bCalcVir || (ir->efep != efepNO),
1577 ir->LincsWarnAngle, &warn,
1579 th == 0 ? vir_r_m_dr : lincsd->th[th].vir_r_m_dr);
1582 if (ir->efep != efepNO)
1584 real dt_2, dvdl = 0;
1586 dt_2 = 1.0/(ir->delta_t*ir->delta_t);
1587 for (i = 0; (i < lincsd->nc); i++)
1589 dvdl -= lincsd->mlambda[i]*dt_2*lincsd->ddist[i];
1594 if (bLog && fplog && lincsd->nc > 0)
1596 fprintf(fplog, " Rel. Constraint Deviation: RMS MAX between atoms\n");
1597 fprintf(fplog, " Before LINCS %.6f %.6f %6d %6d\n",
1598 sqrt(p_ssd/ncons_loc), p_max,
1599 ddglatnr(cr->dd, lincsd->bla[2*p_imax]),
1600 ddglatnr(cr->dd, lincsd->bla[2*p_imax+1]));
1604 cconerr(cr->dd, lincsd->nc, lincsd->bla, lincsd->bllen, xprime, pbc,
1605 &ncons_loc, &p_ssd, &p_max, &p_imax);
1606 /* Check if we are doing the second part of SD */
1607 if (ir->eI == eiSD2 && v == NULL)
1615 lincsd->rmsd_data[0] = ncons_loc;
1616 lincsd->rmsd_data[i] = p_ssd;
1620 lincsd->rmsd_data[0] = 0;
1621 lincsd->rmsd_data[1] = 0;
1622 lincsd->rmsd_data[2] = 0;
1624 if (bLog && fplog && lincsd->nc > 0)
1627 " After LINCS %.6f %.6f %6d %6d\n\n",
1628 sqrt(p_ssd/ncons_loc), p_max,
1629 ddglatnr(cr->dd, lincsd->bla[2*p_imax]),
1630 ddglatnr(cr->dd, lincsd->bla[2*p_imax+1]));
1637 cconerr(cr->dd, lincsd->nc, lincsd->bla, lincsd->bllen, xprime, pbc,
1638 &ncons_loc, &p_ssd, &p_max, &p_imax);
1641 sprintf(buf3, " in simulation %d", cr->ms->sim);
1647 sprintf(buf, "\nStep %s, time %g (ps) LINCS WARNING%s\n"
1648 "relative constraint deviation after LINCS:\n"
1649 "rms %.6f, max %.6f (between atoms %d and %d)\n",
1650 gmx_step_str(step, buf2), ir->init_t+step*ir->delta_t,
1652 sqrt(p_ssd/ncons_loc), p_max,
1653 ddglatnr(cr->dd, lincsd->bla[2*p_imax]),
1654 ddglatnr(cr->dd, lincsd->bla[2*p_imax+1]));
1657 fprintf(fplog, "%s", buf);
1659 fprintf(stderr, "%s", buf);
1660 lincs_warning(fplog, cr->dd, x, xprime, pbc,
1661 lincsd->nc, lincsd->bla, lincsd->bllen,
1662 ir->LincsWarnAngle, maxwarn, warncount);
1664 bOK = (p_max < 0.5);
1667 if (lincsd->ncg_flex)
1669 for (i = 0; (i < lincsd->nc); i++)
1671 if (lincsd->bllen0[i] == 0 && lincsd->ddist[i] == 0)
1673 lincsd->bllen[i] = 0;
1680 /* The OpenMP parallel region of constrain_lincs for derivatives */
1681 #pragma omp parallel num_threads(lincsd->nth)
1683 int th = gmx_omp_get_thread_num();
1685 do_lincsp(x, xprime, min_proj, pbc, lincsd, th,
1686 md->invmass, econq, ir->efep != efepNO ? dvdlambda : NULL,
1687 bCalcVir, th == 0 ? vir_r_m_dr : lincsd->th[th].vir_r_m_dr);
1691 if (bCalcVir && lincsd->nth > 1)
1693 for (i = 1; i < lincsd->nth; i++)
1695 m_add(vir_r_m_dr, lincsd->th[i].vir_r_m_dr, vir_r_m_dr);
1699 /* count assuming nit=1 */
1700 inc_nrnb(nrnb, eNR_LINCS, lincsd->nc);
1701 inc_nrnb(nrnb, eNR_LINCSMAT, (2+lincsd->nOrder)*lincsd->ncc);
1702 if (lincsd->ntriangle > 0)
1704 inc_nrnb(nrnb, eNR_LINCSMAT, lincsd->nOrder*lincsd->ncc_triangle);
1708 inc_nrnb(nrnb, eNR_CONSTR_V, lincsd->nc*2);
1712 inc_nrnb(nrnb, eNR_CONSTR_VIR, lincsd->nc);