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"
54 #include "gmx_omp_nthreads.h"
56 #include "gromacs/fileio/gmxfio.h"
57 #include "gromacs/utility/gmxomp.h"
60 int b0; /* first constraint for this thread */
61 int b1; /* b1-1 is the last constraint for this thread */
62 int nind; /* number of indices */
63 int *ind; /* constraint index for updating atom data */
64 int nind_r; /* number of indices */
65 int *ind_r; /* constraint index for updating atom data */
66 int ind_nalloc; /* allocation size of ind and ind_r */
67 tensor vir_r_m_dr; /* temporary variable for virial calculation */
70 typedef struct gmx_lincsdata {
71 int ncg; /* the global number of constraints */
72 int ncg_flex; /* the global number of flexible constraints */
73 int ncg_triangle; /* the global number of constraints in triangles */
74 int nIter; /* the number of iterations */
75 int nOrder; /* the order of the matrix expansion */
76 int nc; /* the number of constraints */
77 int nc_alloc; /* the number we allocated memory for */
78 int ncc; /* the number of constraint connections */
79 int ncc_alloc; /* the number we allocated memory for */
80 real matlam; /* the FE lambda value used for filling blc and blmf */
81 real *bllen0; /* the reference distance in topology A */
82 real *ddist; /* the reference distance in top B - the r.d. in top A */
83 int *bla; /* the atom pairs involved in the constraints */
84 real *blc; /* 1/sqrt(invmass1 + invmass2) */
85 real *blc1; /* as blc, but with all masses 1 */
86 int *blnr; /* index into blbnb and blmf */
87 int *blbnb; /* list of constraint connections */
88 int ntriangle; /* the local number of constraints in triangles */
89 int *triangle; /* the list of triangle constraints */
90 int *tri_bits; /* the bits tell if the matrix element should be used */
91 int ncc_triangle; /* the number of constraint connections in triangles */
92 gmx_bool bCommIter; /* communicate before each LINCS interation */
93 real *blmf; /* matrix of mass factors for constraint connections */
94 real *blmf1; /* as blmf, but with all masses 1 */
95 real *bllen; /* the reference bond length */
96 int nth; /* The number of threads doing LINCS */
97 lincs_thread_t *th; /* LINCS thread division */
98 unsigned *atf; /* atom flags for thread parallelization */
99 int atf_nalloc; /* allocation size of atf */
100 /* arrays for temporary storage in the LINCS algorithm */
107 real *mlambda; /* the Lagrange multipliers * -1 */
108 /* storage for the constraint RMS relative deviation output */
112 real *lincs_rmsd_data(struct gmx_lincsdata *lincsd)
114 return lincsd->rmsd_data;
117 real lincs_rmsd(struct gmx_lincsdata *lincsd, gmx_bool bSD2)
119 if (lincsd->rmsd_data[0] > 0)
121 return sqrt(lincsd->rmsd_data[bSD2 ? 2 : 1]/lincsd->rmsd_data[0]);
129 /* Do a set of nrec LINCS matrix multiplications.
130 * This function will return with up to date thread-local
131 * constraint data, without an OpenMP barrier.
133 static void lincs_matrix_expand(const struct gmx_lincsdata *lincsd,
136 real *rhs1, real *rhs2, real *sol)
138 int nrec, rec, b, j, n, nr0, nr1;
140 int ntriangle, tb, bits;
141 const int *blnr = lincsd->blnr, *blbnb = lincsd->blbnb;
142 const int *triangle = lincsd->triangle, *tri_bits = lincsd->tri_bits;
144 ntriangle = lincsd->ntriangle;
145 nrec = lincsd->nOrder;
147 for (rec = 0; rec < nrec; rec++)
150 for (b = b0; b < b1; b++)
153 for (n = blnr[b]; n < blnr[b+1]; n++)
156 mvb = mvb + blcc[n]*rhs1[j];
159 sol[b] = sol[b] + mvb;
164 } /* nrec*(ncons+2*nrtot) flops */
168 /* Perform an extra nrec recursions for only the constraints
169 * involved in rigid triangles.
170 * In this way their accuracy should come close to those of the other
171 * constraints, since traingles of constraints can produce eigenvalues
172 * around 0.7, while the effective eigenvalue for bond constraints
173 * is around 0.4 (and 0.7*0.7=0.5).
175 /* We need to copy the temporary array, since only the elements
176 * for constraints involved in triangles are updated and then
177 * the pointers are swapped. This saving copying the whole arrary.
178 * We need barrier as other threads might still be reading from rhs2.
181 for (b = b0; b < b1; b++)
188 for (rec = 0; rec < nrec; rec++)
190 for (tb = 0; tb < ntriangle; tb++)
197 for (n = nr0; n < nr1; n++)
199 if (bits & (1<<(n-nr0)))
202 mvb = mvb + blcc[n]*rhs1[j];
206 sol[b] = sol[b] + mvb;
212 } /* flops count is missing here */
214 /* We need a barrier here as the calling routine will continue
215 * to operate on the thread-local constraints without barrier.
221 static void lincs_update_atoms_noind(int ncons, const int *bla,
223 const real *fac, rvec *r,
228 real mvb, im1, im2, tmp0, tmp1, tmp2;
232 for (b = 0; b < ncons; b++)
248 } /* 16 ncons flops */
252 for (b = 0; b < ncons; b++)
270 static void lincs_update_atoms_ind(int ncons, const int *ind, const int *bla,
272 const real *fac, rvec *r,
277 real mvb, im1, im2, tmp0, tmp1, tmp2;
281 for (bi = 0; bi < ncons; bi++)
298 } /* 16 ncons flops */
302 for (bi = 0; bi < ncons; bi++)
317 } /* 16 ncons flops */
321 static void lincs_update_atoms(struct gmx_lincsdata *li, int th,
323 const real *fac, rvec *r,
329 /* Single thread, we simply update for all constraints */
330 lincs_update_atoms_noind(li->nc, li->bla, prefac, fac, r, invmass, x);
334 /* Update the atom vector components for our thread local
335 * constraints that only access our local atom range.
336 * This can be done without a barrier.
338 lincs_update_atoms_ind(li->th[th].nind, li->th[th].ind,
339 li->bla, prefac, fac, r, invmass, x);
341 if (li->th[li->nth].nind > 0)
343 /* Update the constraints that operate on atoms
344 * in multiple thread atom blocks on the master thread.
349 lincs_update_atoms_ind(li->th[li->nth].nind,
351 li->bla, prefac, fac, r, invmass, x);
357 /* LINCS projection, works on derivatives of the coordinates */
358 static void do_lincsp(rvec *x, rvec *f, rvec *fp, t_pbc *pbc,
359 struct gmx_lincsdata *lincsd, int th,
361 int econq, real *dvdlambda,
362 gmx_bool bCalcVir, tensor rmdf)
364 int b0, b1, b, i, j, k, n;
365 real tmp0, tmp1, tmp2, im1, im2, mvb, rlen, len, wfac, lam;
367 int *bla, *blnr, *blbnb;
369 real *blc, *blmf, *blcc, *rhs1, *rhs2, *sol;
371 b0 = lincsd->th[th].b0;
372 b1 = lincsd->th[th].b1;
377 blbnb = lincsd->blbnb;
378 if (econq != econqForce)
380 /* Use mass-weighted parameters */
386 /* Use non mass-weighted parameters */
388 blmf = lincsd->blmf1;
390 blcc = lincsd->tmpncc;
395 /* Compute normalized i-j vectors */
398 for (b = b0; b < b1; b++)
400 pbc_dx_aiuc(pbc, x[bla[2*b]], x[bla[2*b+1]], dx);
406 for (b = b0; b < b1; b++)
408 rvec_sub(x[bla[2*b]], x[bla[2*b+1]], dx);
410 } /* 16 ncons flops */
414 for (b = b0; b < b1; b++)
421 for (n = blnr[b]; n < blnr[b+1]; n++)
424 blcc[n] = blmf[n]*(tmp0*r[k][0] + tmp1*r[k][1] + tmp2*r[k][2]);
426 mvb = blc[b]*(tmp0*(f[i][0] - f[j][0]) +
427 tmp1*(f[i][1] - f[j][1]) +
428 tmp2*(f[i][2] - f[j][2]));
433 /* Together: 23*ncons + 6*nrtot flops */
435 lincs_matrix_expand(lincsd, b0, b1, blcc, rhs1, rhs2, sol);
436 /* nrec*(ncons+2*nrtot) flops */
438 if (econq == econqDeriv_FlexCon)
440 /* We only want to constraint the flexible constraints,
441 * so we mask out the normal ones by setting sol to 0.
443 for (b = b0; b < b1; b++)
445 if (!(lincsd->bllen0[b] == 0 && lincsd->ddist[b] == 0))
452 /* We multiply sol by blc, so we can use lincs_update_atoms for OpenMP */
453 for (b = b0; b < b1; b++)
458 /* When constraining forces, we should not use mass weighting,
459 * so we pass invmass=NULL, which results in the use of 1 for all atoms.
461 lincs_update_atoms(lincsd, th, 1.0, sol, r,
462 (econq != econqForce) ? invmass : NULL, fp);
464 if (dvdlambda != NULL)
467 for (b = b0; b < b1; b++)
469 *dvdlambda -= sol[b]*lincsd->ddist[b];
476 /* Constraint virial,
477 * determines sum r_bond x delta f,
478 * where delta f is the constraint correction
479 * of the quantity that is being constrained.
481 for (b = b0; b < b1; b++)
483 mvb = lincsd->bllen[b]*sol[b];
484 for (i = 0; i < DIM; i++)
487 for (j = 0; j < DIM; j++)
489 rmdf[i][j] += tmp1*r[b][j];
492 } /* 23 ncons flops */
496 static void do_lincs(rvec *x, rvec *xp, matrix box, t_pbc *pbc,
497 struct gmx_lincsdata *lincsd, int th,
500 gmx_bool bCalcLambda,
501 real wangle, int *warn,
503 gmx_bool bCalcVir, tensor vir_r_m_dr)
505 int b0, b1, b, i, j, k, n, iter;
506 real tmp0, tmp1, tmp2, im1, im2, mvb, rlen, len, len2, dlen2, wfac;
508 int *bla, *blnr, *blbnb;
510 real *blc, *blmf, *bllen, *blcc, *rhs1, *rhs2, *sol, *blc_sol, *mlambda;
513 b0 = lincsd->th[th].b0;
514 b1 = lincsd->th[th].b1;
519 blbnb = lincsd->blbnb;
522 bllen = lincsd->bllen;
523 blcc = lincsd->tmpncc;
527 blc_sol = lincsd->tmp4;
528 mlambda = lincsd->mlambda;
530 if (DOMAINDECOMP(cr) && cr->dd->constraints)
532 nlocat = dd_constraints_nlocalatoms(cr->dd);
534 else if (PARTDECOMP(cr))
536 nlocat = pd_constraints_nlocalatoms(cr->pd);
545 /* Compute normalized i-j vectors */
546 for (b = b0; b < b1; b++)
548 pbc_dx_aiuc(pbc, x[bla[2*b]], x[bla[2*b+1]], dx);
552 for (b = b0; b < b1; b++)
554 for (n = blnr[b]; n < blnr[b+1]; n++)
556 blcc[n] = blmf[n]*iprod(r[b], r[blbnb[n]]);
558 pbc_dx_aiuc(pbc, xp[bla[2*b]], xp[bla[2*b+1]], dx);
559 mvb = blc[b]*(iprod(r[b], dx) - bllen[b]);
566 /* Compute normalized i-j vectors */
567 for (b = b0; b < b1; b++)
571 tmp0 = x[i][0] - x[j][0];
572 tmp1 = x[i][1] - x[j][1];
573 tmp2 = x[i][2] - x[j][2];
574 rlen = gmx_invsqrt(tmp0*tmp0+tmp1*tmp1+tmp2*tmp2);
578 } /* 16 ncons flops */
581 for (b = b0; b < b1; b++)
589 for (n = blnr[b]; n < blnr[b+1]; n++)
592 blcc[n] = blmf[n]*(tmp0*r[k][0] + tmp1*r[k][1] + tmp2*r[k][2]);
594 mvb = blc[b]*(tmp0*(xp[i][0] - xp[j][0]) +
595 tmp1*(xp[i][1] - xp[j][1]) +
596 tmp2*(xp[i][2] - xp[j][2]) - len);
601 /* Together: 26*ncons + 6*nrtot flops */
604 lincs_matrix_expand(lincsd, b0, b1, blcc, rhs1, rhs2, sol);
605 /* nrec*(ncons+2*nrtot) flops */
607 for (b = b0; b < b1; b++)
609 mlambda[b] = blc[b]*sol[b];
612 /* Update the coordinates */
613 lincs_update_atoms(lincsd, th, 1.0, mlambda, r, invmass, xp);
616 ******** Correction for centripetal effects ********
619 wfac = cos(DEG2RAD*wangle);
622 for (iter = 0; iter < lincsd->nIter; iter++)
624 if ((lincsd->bCommIter && DOMAINDECOMP(cr) && cr->dd->constraints) ||
630 /* Communicate the corrected non-local coordinates */
631 if (DOMAINDECOMP(cr))
633 dd_move_x_constraints(cr->dd, box, xp, NULL);
637 pd_move_x_constraints(cr, xp, NULL);
643 for (b = b0; b < b1; b++)
648 pbc_dx_aiuc(pbc, xp[bla[2*b]], xp[bla[2*b+1]], dx);
652 rvec_sub(xp[bla[2*b]], xp[bla[2*b+1]], dx);
655 dlen2 = 2*len2 - norm2(dx);
656 if (dlen2 < wfac*len2 && (nlocat == NULL || nlocat[b]))
662 mvb = blc[b]*(len - dlen2*gmx_invsqrt(dlen2));
670 } /* 20*ncons flops */
672 lincs_matrix_expand(lincsd, b0, b1, blcc, rhs1, rhs2, sol);
673 /* nrec*(ncons+2*nrtot) flops */
675 for (b = b0; b < b1; b++)
682 /* Update the coordinates */
683 lincs_update_atoms(lincsd, th, 1.0, blc_sol, r, invmass, xp);
685 /* nit*ncons*(37+9*nrec) flops */
689 /* Update the velocities */
690 lincs_update_atoms(lincsd, th, invdt, mlambda, r, invmass, v);
694 if (nlocat != NULL && bCalcLambda)
696 /* In lincs_update_atoms thread might cross-read mlambda */
699 /* Only account for local atoms */
700 for (b = b0; b < b1; b++)
702 mlambda[b] *= 0.5*nlocat[b];
708 /* Constraint virial */
709 for (b = b0; b < b1; b++)
711 tmp0 = -bllen[b]*mlambda[b];
712 for (i = 0; i < DIM; i++)
715 for (j = 0; j < DIM; j++)
717 vir_r_m_dr[i][j] -= tmp1*r[b][j];
720 } /* 22 ncons flops */
724 * 26*ncons + 6*nrtot + nrec*(ncons+2*nrtot)
725 * + nit * (20*ncons + nrec*(ncons+2*nrtot) + 17 ncons)
727 * (26+nrec)*ncons + (6+2*nrec)*nrtot
728 * + nit * ((37+nrec)*ncons + 2*nrec*nrtot)
730 * (63+nrec)*ncons + (6+4*nrec)*nrtot
734 void set_lincs_matrix(struct gmx_lincsdata *li, real *invmass, real lambda)
736 int i, a1, a2, n, k, sign, center;
738 const real invsqrt2 = 0.7071067811865475244;
740 for (i = 0; (i < li->nc); i++)
744 li->blc[i] = gmx_invsqrt(invmass[a1] + invmass[a2]);
745 li->blc1[i] = invsqrt2;
748 /* Construct the coupling coefficient matrix blmf */
750 li->ncc_triangle = 0;
751 for (i = 0; (i < li->nc); i++)
755 for (n = li->blnr[i]; (n < li->blnr[i+1]); n++)
758 if (a1 == li->bla[2*k] || a2 == li->bla[2*k+1])
766 if (a1 == li->bla[2*k] || a1 == li->bla[2*k+1])
776 li->blmf[n] = sign*invmass[center]*li->blc[i]*li->blc[k];
777 li->blmf1[n] = sign*0.5;
778 if (li->ncg_triangle > 0)
780 /* Look for constraint triangles */
781 for (nk = li->blnr[k]; (nk < li->blnr[k+1]); nk++)
784 if (kk != i && kk != k &&
785 (li->bla[2*kk] == end || li->bla[2*kk+1] == end))
787 if (li->ntriangle == 0 ||
788 li->triangle[li->ntriangle-1] < i)
790 /* Add this constraint to the triangle list */
791 li->triangle[li->ntriangle] = i;
792 li->tri_bits[li->ntriangle] = 0;
794 if (li->blnr[i+1] - li->blnr[i] > sizeof(li->tri_bits[0])*8 - 1)
796 gmx_fatal(FARGS, "A constraint is connected to %d constraints, this is more than the %d allowed for constraints participating in triangles",
797 li->blnr[i+1] - li->blnr[i],
798 sizeof(li->tri_bits[0])*8-1);
801 li->tri_bits[li->ntriangle-1] |= (1<<(n-li->blnr[i]));
811 fprintf(debug, "Of the %d constraints %d participate in triangles\n",
812 li->nc, li->ntriangle);
813 fprintf(debug, "There are %d couplings of which %d in triangles\n",
814 li->ncc, li->ncc_triangle);
818 * so we know with which lambda value the masses have been set.
823 static int count_triangle_constraints(t_ilist *ilist, t_blocka *at2con)
826 int c0, a00, a01, n1, c1, a10, a11, ac1, n2, c2, a20, a21;
829 t_iatom *ia1, *ia2, *iap;
831 ncon1 = ilist[F_CONSTR].nr/3;
832 ncon_tot = ncon1 + ilist[F_CONSTRNC].nr/3;
834 ia1 = ilist[F_CONSTR].iatoms;
835 ia2 = ilist[F_CONSTRNC].iatoms;
838 for (c0 = 0; c0 < ncon_tot; c0++)
841 iap = constr_iatomptr(ncon1, ia1, ia2, c0);
844 for (n1 = at2con->index[a01]; n1 < at2con->index[a01+1]; n1++)
849 iap = constr_iatomptr(ncon1, ia1, ia2, c1);
860 for (n2 = at2con->index[ac1]; n2 < at2con->index[ac1+1]; n2++)
863 if (c2 != c0 && c2 != c1)
865 iap = constr_iatomptr(ncon1, ia1, ia2, c2);
868 if (a20 == a00 || a21 == a00)
882 return ncon_triangle;
885 static gmx_bool more_than_two_sequential_constraints(const t_ilist *ilist,
886 const t_blocka *at2con)
888 t_iatom *ia1, *ia2, *iap;
889 int ncon1, ncon_tot, c;
891 gmx_bool bMoreThanTwoSequentialConstraints;
893 ncon1 = ilist[F_CONSTR].nr/3;
894 ncon_tot = ncon1 + ilist[F_CONSTRNC].nr/3;
896 ia1 = ilist[F_CONSTR].iatoms;
897 ia2 = ilist[F_CONSTRNC].iatoms;
899 bMoreThanTwoSequentialConstraints = FALSE;
900 for (c = 0; c < ncon_tot && !bMoreThanTwoSequentialConstraints; c++)
902 iap = constr_iatomptr(ncon1, ia1, ia2, c);
905 /* Check if this constraint has constraints connected at both atoms */
906 if (at2con->index[a1+1] - at2con->index[a1] > 1 &&
907 at2con->index[a2+1] - at2con->index[a2] > 1)
909 bMoreThanTwoSequentialConstraints = TRUE;
913 return bMoreThanTwoSequentialConstraints;
916 static int int_comp(const void *a, const void *b)
918 return (*(int *)a) - (*(int *)b);
921 gmx_lincsdata_t init_lincs(FILE *fplog, gmx_mtop_t *mtop,
922 int nflexcon_global, t_blocka *at2con,
923 gmx_bool bPLINCS, int nIter, int nProjOrder)
925 struct gmx_lincsdata *li;
931 fprintf(fplog, "\nInitializing%s LINear Constraint Solver\n",
932 bPLINCS ? " Parallel" : "");
938 gmx_mtop_ftype_count(mtop, F_CONSTR) +
939 gmx_mtop_ftype_count(mtop, F_CONSTRNC);
940 li->ncg_flex = nflexcon_global;
943 li->nOrder = nProjOrder;
945 li->ncg_triangle = 0;
946 li->bCommIter = FALSE;
947 for (mb = 0; mb < mtop->nmolblock; mb++)
949 molt = &mtop->moltype[mtop->molblock[mb].type];
951 mtop->molblock[mb].nmol*
952 count_triangle_constraints(molt->ilist,
953 &at2con[mtop->molblock[mb].type]);
954 if (bPLINCS && li->bCommIter == FALSE)
956 /* Check if we need to communicate not only before LINCS,
957 * but also before each iteration.
958 * The check for only two sequential constraints is only
959 * useful for the common case of H-bond only constraints.
960 * With more effort we could also make it useful for small
961 * molecules with nr. sequential constraints <= nOrder-1.
963 li->bCommIter = (li->nOrder < 1 || more_than_two_sequential_constraints(molt->ilist, &at2con[mtop->molblock[mb].type]));
966 if (debug && bPLINCS)
968 fprintf(debug, "PLINCS communication before each iteration: %d\n",
972 /* LINCS can run on any number of threads.
973 * Currently the number is fixed for the whole simulation,
974 * but it could be set in set_lincs().
976 li->nth = gmx_omp_nthreads_get(emntLINCS);
983 /* Allocate an extra elements for "thread-overlap" constraints */
984 snew(li->th, li->nth+1);
988 fprintf(debug, "LINCS: using %d threads\n", li->nth);
991 if (bPLINCS || li->ncg_triangle > 0)
993 please_cite(fplog, "Hess2008a");
997 please_cite(fplog, "Hess97a");
1002 fprintf(fplog, "The number of constraints is %d\n", li->ncg);
1005 fprintf(fplog, "There are inter charge-group constraints,\n"
1006 "will communicate selected coordinates each lincs iteration\n");
1008 if (li->ncg_triangle > 0)
1011 "%d constraints are involved in constraint triangles,\n"
1012 "will apply an additional matrix expansion of order %d for couplings\n"
1013 "between constraints inside triangles\n",
1014 li->ncg_triangle, li->nOrder);
1021 /* Sets up the work division over the threads */
1022 static void lincs_thread_setup(struct gmx_lincsdata *li, int natoms)
1024 lincs_thread_t *li_m;
1029 if (natoms > li->atf_nalloc)
1031 li->atf_nalloc = over_alloc_large(natoms);
1032 srenew(li->atf, li->atf_nalloc);
1036 /* Clear the atom flags */
1037 for (a = 0; a < natoms; a++)
1042 for (th = 0; th < li->nth; th++)
1044 lincs_thread_t *li_th;
1047 li_th = &li->th[th];
1049 /* The constraints are divided equally over the threads */
1050 li_th->b0 = (li->nc* th )/li->nth;
1051 li_th->b1 = (li->nc*(th+1))/li->nth;
1053 if (th < sizeof(*atf)*8)
1055 /* For each atom set a flag for constraints from each */
1056 for (b = li_th->b0; b < li_th->b1; b++)
1058 atf[li->bla[b*2] ] |= (1U<<th);
1059 atf[li->bla[b*2+1]] |= (1U<<th);
1064 #pragma omp parallel for num_threads(li->nth) schedule(static)
1065 for (th = 0; th < li->nth; th++)
1067 lincs_thread_t *li_th;
1071 li_th = &li->th[th];
1073 if (li_th->b1 - li_th->b0 > li_th->ind_nalloc)
1075 li_th->ind_nalloc = over_alloc_large(li_th->b1-li_th->b0);
1076 srenew(li_th->ind, li_th->ind_nalloc);
1077 srenew(li_th->ind_r, li_th->ind_nalloc);
1080 if (th < sizeof(*atf)*8)
1082 mask = (1U<<th) - 1U;
1086 for (b = li_th->b0; b < li_th->b1; b++)
1088 /* We let the constraint with the lowest thread index
1089 * operate on atoms with constraints from multiple threads.
1091 if (((atf[li->bla[b*2]] & mask) == 0) &&
1092 ((atf[li->bla[b*2+1]] & mask) == 0))
1094 /* Add the constraint to the local atom update index */
1095 li_th->ind[li_th->nind++] = b;
1099 /* Add the constraint to the rest block */
1100 li_th->ind_r[li_th->nind_r++] = b;
1106 /* We are out of bits, assign all constraints to rest */
1107 for (b = li_th->b0; b < li_th->b1; b++)
1109 li_th->ind_r[li_th->nind_r++] = b;
1114 /* We need to copy all constraints which have not be assigned
1115 * to a thread to a separate list which will be handled by one thread.
1117 li_m = &li->th[li->nth];
1120 for (th = 0; th < li->nth; th++)
1122 lincs_thread_t *li_th;
1125 li_th = &li->th[th];
1127 if (li_m->nind + li_th->nind_r > li_m->ind_nalloc)
1129 li_m->ind_nalloc = over_alloc_large(li_m->nind+li_th->nind_r);
1130 srenew(li_m->ind, li_m->ind_nalloc);
1133 for (b = 0; b < li_th->nind_r; b++)
1135 li_m->ind[li_m->nind++] = li_th->ind_r[b];
1140 fprintf(debug, "LINCS thread %d: %d constraints\n",
1147 fprintf(debug, "LINCS thread r: %d constraints\n",
1153 void set_lincs(t_idef *idef, t_mdatoms *md,
1154 gmx_bool bDynamics, t_commrec *cr,
1155 struct gmx_lincsdata *li)
1157 int start, natoms, nflexcon;
1160 int i, k, ncc_alloc, ni, con, nconnect, concon;
1162 real lenA = 0, lenB;
1167 /* Zero the thread index ranges.
1168 * Otherwise without local constraints we could return with old ranges.
1170 for (i = 0; i < li->nth; i++)
1178 li->th[li->nth].nind = 0;
1181 /* This is the local topology, so there are only F_CONSTR constraints */
1182 if (idef->il[F_CONSTR].nr == 0)
1184 /* There are no constraints,
1185 * we do not need to fill any data structures.
1192 fprintf(debug, "Building the LINCS connectivity\n");
1195 if (DOMAINDECOMP(cr))
1197 if (cr->dd->constraints)
1199 dd_get_constraint_range(cr->dd, &start, &natoms);
1203 natoms = cr->dd->nat_home;
1207 else if (PARTDECOMP(cr))
1209 pd_get_constraint_range(cr->pd, &start, &natoms);
1214 natoms = md->homenr;
1216 at2con = make_at2con(start, natoms, idef->il, idef->iparams, bDynamics,
1220 if (idef->il[F_CONSTR].nr/3 > li->nc_alloc || li->nc_alloc == 0)
1222 li->nc_alloc = over_alloc_dd(idef->il[F_CONSTR].nr/3);
1223 srenew(li->bllen0, li->nc_alloc);
1224 srenew(li->ddist, li->nc_alloc);
1225 srenew(li->bla, 2*li->nc_alloc);
1226 srenew(li->blc, li->nc_alloc);
1227 srenew(li->blc1, li->nc_alloc);
1228 srenew(li->blnr, li->nc_alloc+1);
1229 srenew(li->bllen, li->nc_alloc);
1230 srenew(li->tmpv, li->nc_alloc);
1231 srenew(li->tmp1, li->nc_alloc);
1232 srenew(li->tmp2, li->nc_alloc);
1233 srenew(li->tmp3, li->nc_alloc);
1234 srenew(li->tmp4, li->nc_alloc);
1235 srenew(li->mlambda, li->nc_alloc);
1236 if (li->ncg_triangle > 0)
1238 /* This is allocating too much, but it is difficult to improve */
1239 srenew(li->triangle, li->nc_alloc);
1240 srenew(li->tri_bits, li->nc_alloc);
1244 iatom = idef->il[F_CONSTR].iatoms;
1246 ncc_alloc = li->ncc_alloc;
1249 ni = idef->il[F_CONSTR].nr/3;
1253 li->blnr[con] = nconnect;
1254 for (i = 0; i < ni; i++)
1260 lenA = idef->iparams[type].constr.dA;
1261 lenB = idef->iparams[type].constr.dB;
1262 /* Skip the flexible constraints when not doing dynamics */
1263 if (bDynamics || lenA != 0 || lenB != 0)
1265 li->bllen0[con] = lenA;
1266 li->ddist[con] = lenB - lenA;
1267 /* Set the length to the topology A length */
1268 li->bllen[con] = li->bllen0[con];
1269 li->bla[2*con] = a1;
1270 li->bla[2*con+1] = a2;
1271 /* Construct the constraint connection matrix blbnb */
1272 for (k = at2con.index[a1-start]; k < at2con.index[a1-start+1]; k++)
1274 concon = at2con.a[k];
1277 if (nconnect >= ncc_alloc)
1279 ncc_alloc = over_alloc_small(nconnect+1);
1280 srenew(li->blbnb, ncc_alloc);
1282 li->blbnb[nconnect++] = concon;
1285 for (k = at2con.index[a2-start]; k < at2con.index[a2-start+1]; k++)
1287 concon = at2con.a[k];
1290 if (nconnect+1 > ncc_alloc)
1292 ncc_alloc = over_alloc_small(nconnect+1);
1293 srenew(li->blbnb, ncc_alloc);
1295 li->blbnb[nconnect++] = concon;
1298 li->blnr[con+1] = nconnect;
1302 /* Order the blbnb matrix to optimize memory access */
1303 qsort(&(li->blbnb[li->blnr[con]]), li->blnr[con+1]-li->blnr[con],
1304 sizeof(li->blbnb[0]), int_comp);
1306 /* Increase the constraint count */
1311 done_blocka(&at2con);
1313 /* This is the real number of constraints,
1314 * without dynamics the flexible constraints are not present.
1318 li->ncc = li->blnr[con];
1321 /* Since the matrix is static, we can free some memory */
1322 ncc_alloc = li->ncc;
1323 srenew(li->blbnb, ncc_alloc);
1326 if (ncc_alloc > li->ncc_alloc)
1328 li->ncc_alloc = ncc_alloc;
1329 srenew(li->blmf, li->ncc_alloc);
1330 srenew(li->blmf1, li->ncc_alloc);
1331 srenew(li->tmpncc, li->ncc_alloc);
1336 fprintf(debug, "Number of constraints is %d, couplings %d\n",
1343 li->th[0].b1 = li->nc;
1347 lincs_thread_setup(li, md->nr);
1350 set_lincs_matrix(li, md->invmass, md->lambda);
1353 static void lincs_warning(FILE *fplog,
1354 gmx_domdec_t *dd, rvec *x, rvec *xprime, t_pbc *pbc,
1355 int ncons, int *bla, real *bllen, real wangle,
1356 int maxwarn, int *warncount)
1360 real wfac, d0, d1, cosine;
1363 wfac = cos(DEG2RAD*wangle);
1365 sprintf(buf, "bonds that rotated more than %g degrees:\n"
1366 " atom 1 atom 2 angle previous, current, constraint length\n",
1368 fprintf(stderr, "%s", buf);
1371 fprintf(fplog, "%s", buf);
1374 for (b = 0; b < ncons; b++)
1380 pbc_dx_aiuc(pbc, x[i], x[j], v0);
1381 pbc_dx_aiuc(pbc, xprime[i], xprime[j], v1);
1385 rvec_sub(x[i], x[j], v0);
1386 rvec_sub(xprime[i], xprime[j], v1);
1390 cosine = iprod(v0, v1)/(d0*d1);
1393 sprintf(buf, " %6d %6d %5.1f %8.4f %8.4f %8.4f\n",
1394 ddglatnr(dd, i), ddglatnr(dd, j),
1395 RAD2DEG*acos(cosine), d0, d1, bllen[b]);
1396 fprintf(stderr, "%s", buf);
1399 fprintf(fplog, "%s", buf);
1401 if (!gmx_isfinite(d1))
1403 gmx_fatal(FARGS, "Bond length not finite.");
1409 if (*warncount > maxwarn)
1411 too_many_constraint_warnings(econtLINCS, *warncount);
1415 static void cconerr(gmx_domdec_t *dd,
1416 int ncons, int *bla, real *bllen, rvec *x, t_pbc *pbc,
1417 real *ncons_loc, real *ssd, real *max, int *imax)
1419 real len, d, ma, ssd2, r2;
1420 int *nlocat, count, b, im;
1423 if (dd && dd->constraints)
1425 nlocat = dd_constraints_nlocalatoms(dd);
1436 for (b = 0; b < ncons; b++)
1440 pbc_dx_aiuc(pbc, x[bla[2*b]], x[bla[2*b+1]], dx);
1444 rvec_sub(x[bla[2*b]], x[bla[2*b+1]], dx);
1447 len = r2*gmx_invsqrt(r2);
1448 d = fabs(len/bllen[b]-1);
1449 if (d > ma && (nlocat == NULL || nlocat[b]))
1461 ssd2 += nlocat[b]*d*d;
1466 *ncons_loc = (nlocat ? 0.5 : 1)*count;
1467 *ssd = (nlocat ? 0.5 : 1)*ssd2;
1472 gmx_bool constrain_lincs(FILE *fplog, gmx_bool bLog, gmx_bool bEner,
1475 struct gmx_lincsdata *lincsd, t_mdatoms *md,
1477 rvec *x, rvec *xprime, rvec *min_proj,
1478 matrix box, t_pbc *pbc,
1479 real lambda, real *dvdlambda,
1480 real invdt, rvec *v,
1481 gmx_bool bCalcVir, tensor vir_r_m_dr,
1484 int maxwarn, int *warncount)
1486 char buf[STRLEN], buf2[22], buf3[STRLEN];
1487 int i, warn, p_imax, error;
1488 real ncons_loc, p_ssd, p_max = 0;
1494 if (lincsd->nc == 0 && cr->dd == NULL)
1498 lincsd->rmsd_data[0] = 0;
1499 if (ir->eI == eiSD2 && v == NULL)
1507 lincsd->rmsd_data[i] = 0;
1513 if (econq == econqCoord)
1515 if (ir->efep != efepNO)
1517 if (md->nMassPerturbed && lincsd->matlam != md->lambda)
1519 set_lincs_matrix(lincsd, md->invmass, md->lambda);
1522 for (i = 0; i < lincsd->nc; i++)
1524 lincsd->bllen[i] = lincsd->bllen0[i] + lambda*lincsd->ddist[i];
1528 if (lincsd->ncg_flex)
1530 /* Set the flexible constraint lengths to the old lengths */
1533 for (i = 0; i < lincsd->nc; i++)
1535 if (lincsd->bllen[i] == 0)
1537 pbc_dx_aiuc(pbc, x[lincsd->bla[2*i]], x[lincsd->bla[2*i+1]], dx);
1538 lincsd->bllen[i] = norm(dx);
1544 for (i = 0; i < lincsd->nc; i++)
1546 if (lincsd->bllen[i] == 0)
1549 sqrt(distance2(x[lincsd->bla[2*i]],
1550 x[lincsd->bla[2*i+1]]));
1558 cconerr(cr->dd, lincsd->nc, lincsd->bla, lincsd->bllen, xprime, pbc,
1559 &ncons_loc, &p_ssd, &p_max, &p_imax);
1562 /* This warn var can be updated by multiple threads
1563 * at the same time. But as we only need to detect
1564 * if a warning occured or not, this is not an issue.
1568 /* The OpenMP parallel region of constrain_lincs for coords */
1569 #pragma omp parallel num_threads(lincsd->nth)
1571 int th = gmx_omp_get_thread_num();
1573 clear_mat(lincsd->th[th].vir_r_m_dr);
1575 do_lincs(x, xprime, box, pbc, lincsd, th,
1577 bCalcVir || (ir->efep != efepNO),
1578 ir->LincsWarnAngle, &warn,
1580 th == 0 ? vir_r_m_dr : lincsd->th[th].vir_r_m_dr);
1583 if (ir->efep != efepNO)
1585 real dt_2, dvdl = 0;
1587 dt_2 = 1.0/(ir->delta_t*ir->delta_t);
1588 for (i = 0; (i < lincsd->nc); i++)
1590 dvdl -= lincsd->mlambda[i]*dt_2*lincsd->ddist[i];
1595 if (bLog && fplog && lincsd->nc > 0)
1597 fprintf(fplog, " Rel. Constraint Deviation: RMS MAX between atoms\n");
1598 fprintf(fplog, " Before LINCS %.6f %.6f %6d %6d\n",
1599 sqrt(p_ssd/ncons_loc), p_max,
1600 ddglatnr(cr->dd, lincsd->bla[2*p_imax]),
1601 ddglatnr(cr->dd, lincsd->bla[2*p_imax+1]));
1605 cconerr(cr->dd, lincsd->nc, lincsd->bla, lincsd->bllen, xprime, pbc,
1606 &ncons_loc, &p_ssd, &p_max, &p_imax);
1607 /* Check if we are doing the second part of SD */
1608 if (ir->eI == eiSD2 && v == NULL)
1616 lincsd->rmsd_data[0] = ncons_loc;
1617 lincsd->rmsd_data[i] = p_ssd;
1621 lincsd->rmsd_data[0] = 0;
1622 lincsd->rmsd_data[1] = 0;
1623 lincsd->rmsd_data[2] = 0;
1625 if (bLog && fplog && lincsd->nc > 0)
1628 " After LINCS %.6f %.6f %6d %6d\n\n",
1629 sqrt(p_ssd/ncons_loc), p_max,
1630 ddglatnr(cr->dd, lincsd->bla[2*p_imax]),
1631 ddglatnr(cr->dd, lincsd->bla[2*p_imax+1]));
1638 cconerr(cr->dd, lincsd->nc, lincsd->bla, lincsd->bllen, xprime, pbc,
1639 &ncons_loc, &p_ssd, &p_max, &p_imax);
1642 sprintf(buf3, " in simulation %d", cr->ms->sim);
1648 sprintf(buf, "\nStep %s, time %g (ps) LINCS WARNING%s\n"
1649 "relative constraint deviation after LINCS:\n"
1650 "rms %.6f, max %.6f (between atoms %d and %d)\n",
1651 gmx_step_str(step, buf2), ir->init_t+step*ir->delta_t,
1653 sqrt(p_ssd/ncons_loc), p_max,
1654 ddglatnr(cr->dd, lincsd->bla[2*p_imax]),
1655 ddglatnr(cr->dd, lincsd->bla[2*p_imax+1]));
1658 fprintf(fplog, "%s", buf);
1660 fprintf(stderr, "%s", buf);
1661 lincs_warning(fplog, cr->dd, x, xprime, pbc,
1662 lincsd->nc, lincsd->bla, lincsd->bllen,
1663 ir->LincsWarnAngle, maxwarn, warncount);
1665 bOK = (p_max < 0.5);
1668 if (lincsd->ncg_flex)
1670 for (i = 0; (i < lincsd->nc); i++)
1672 if (lincsd->bllen0[i] == 0 && lincsd->ddist[i] == 0)
1674 lincsd->bllen[i] = 0;
1681 /* The OpenMP parallel region of constrain_lincs for derivatives */
1682 #pragma omp parallel num_threads(lincsd->nth)
1684 int th = gmx_omp_get_thread_num();
1686 do_lincsp(x, xprime, min_proj, pbc, lincsd, th,
1687 md->invmass, econq, ir->efep != efepNO ? dvdlambda : NULL,
1688 bCalcVir, th == 0 ? vir_r_m_dr : lincsd->th[th].vir_r_m_dr);
1692 if (bCalcVir && lincsd->nth > 1)
1694 for (i = 1; i < lincsd->nth; i++)
1696 m_add(vir_r_m_dr, lincsd->th[i].vir_r_m_dr, vir_r_m_dr);
1700 /* count assuming nit=1 */
1701 inc_nrnb(nrnb, eNR_LINCS, lincsd->nc);
1702 inc_nrnb(nrnb, eNR_LINCSMAT, (2+lincsd->nOrder)*lincsd->ncc);
1703 if (lincsd->ntriangle > 0)
1705 inc_nrnb(nrnb, eNR_LINCSMAT, lincsd->nOrder*lincsd->ncc_triangle);
1709 inc_nrnb(nrnb, eNR_CONSTR_V, lincsd->nc*2);
1713 inc_nrnb(nrnb, eNR_CONSTR_VIR, lincsd->nc);