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, 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.
39 #include "gromacs/pbcutil/pbc.h"
46 #include "gromacs/legacyheaders/types/inputrec.h"
47 #include "gromacs/legacyheaders/types/commrec.h"
48 #include "gromacs/legacyheaders/txtdump.h"
49 #include "gromacs/legacyheaders/names.h"
50 #include "gromacs/legacyheaders/macros.h"
51 #include "gromacs/legacyheaders/gmx_omp_nthreads.h"
53 #include "gromacs/math/utilities.h"
54 #include "gromacs/math/vec.h"
55 #include "gromacs/pbcutil/ishift.h"
56 #include "gromacs/pbcutil/mshift.h"
57 #include "gromacs/utility/fatalerror.h"
58 #include "gromacs/utility/smalloc.h"
60 /* Skip 0 so we have more chance of detecting if we forgot to call set_pbc. */
62 epbcdxRECTANGULAR = 1, epbcdxTRICLINIC,
63 epbcdx2D_RECT, epbcdx2D_TRIC,
64 epbcdx1D_RECT, epbcdx1D_TRIC,
65 epbcdxSCREW_RECT, epbcdxSCREW_TRIC,
66 epbcdxNOPBC, epbcdxUNSUPPORTED
69 /* Margin factor for error message and correction if the box is too skewed */
70 #define BOX_MARGIN 1.0010
71 #define BOX_MARGIN_CORRECT 1.0005
73 int ePBC2npbcdim(int ePBC)
79 case epbcXYZ: npbcdim = 3; break;
80 case epbcXY: npbcdim = 2; break;
81 case epbcSCREW: npbcdim = 3; break;
82 case epbcNONE: npbcdim = 0; break;
83 default: gmx_fatal(FARGS, "Unknown ePBC=%d in ePBC2npbcdim", ePBC);
89 int inputrec2nboundeddim(t_inputrec *ir)
91 if (ir->nwall == 2 && ir->ePBC == epbcXY)
97 return ePBC2npbcdim(ir->ePBC);
101 void dump_pbc(FILE *fp, t_pbc *pbc)
105 fprintf(fp, "ePBCDX = %d\n", pbc->ePBCDX);
106 pr_rvecs(fp, 0, "box", pbc->box, DIM);
107 pr_rvecs(fp, 0, "fbox_diag", &pbc->fbox_diag, 1);
108 pr_rvecs(fp, 0, "hbox_diag", &pbc->hbox_diag, 1);
109 pr_rvecs(fp, 0, "mhbox_diag", &pbc->mhbox_diag, 1);
110 rvec_add(pbc->hbox_diag, pbc->mhbox_diag, sum_box);
111 pr_rvecs(fp, 0, "sum of the above two", &sum_box, 1);
112 fprintf(fp, "max_cutoff2 = %g\n", pbc->max_cutoff2);
113 fprintf(fp, "bLimitDistance = %s\n", EBOOL(pbc->bLimitDistance));
114 fprintf(fp, "limit_distance2 = %g\n", pbc->limit_distance2);
115 fprintf(fp, "ntric_vec = %d\n", pbc->ntric_vec);
116 if (pbc->ntric_vec > 0)
118 pr_ivecs(fp, 0, "tric_shift", pbc->tric_shift, pbc->ntric_vec, FALSE);
119 pr_rvecs(fp, 0, "tric_vec", pbc->tric_vec, pbc->ntric_vec);
123 const char *check_box(int ePBC, matrix box)
129 ePBC = guess_ePBC(box);
132 if (ePBC == epbcNONE)
137 if ((box[XX][YY] != 0) || (box[XX][ZZ] != 0) || (box[YY][ZZ] != 0))
139 ptr = "Only triclinic boxes with the first vector parallel to the x-axis and the second vector in the xy-plane are supported.";
141 else if (ePBC == epbcSCREW && (box[YY][XX] != 0 || box[ZZ][XX] != 0))
143 ptr = "The unit cell can not have off-diagonal x-components with screw pbc";
145 else if (fabs(box[YY][XX]) > BOX_MARGIN*0.5*box[XX][XX] ||
147 (fabs(box[ZZ][XX]) > BOX_MARGIN*0.5*box[XX][XX] ||
148 fabs(box[ZZ][YY]) > BOX_MARGIN*0.5*box[YY][YY])))
150 ptr = "Triclinic box is too skewed.";
160 real max_cutoff2(int ePBC, matrix box)
162 real min_hv2, min_ss;
164 /* Physical limitation of the cut-off
165 * by half the length of the shortest box vector.
167 min_hv2 = min(0.25*norm2(box[XX]), 0.25*norm2(box[YY]));
170 min_hv2 = min(min_hv2, 0.25*norm2(box[ZZ]));
173 /* Limitation to the smallest diagonal element due to optimizations:
174 * checking only linear combinations of single box-vectors (2 in x)
175 * in the grid search and pbc_dx is a lot faster
176 * than checking all possible combinations.
180 min_ss = min(box[XX][XX], box[YY][YY]);
184 min_ss = min(box[XX][XX], min(box[YY][YY]-fabs(box[ZZ][YY]), box[ZZ][ZZ]));
187 return min(min_hv2, min_ss*min_ss);
190 /* this one is mostly harmless... */
191 static gmx_bool bWarnedGuess = FALSE;
193 int guess_ePBC(matrix box)
197 if (box[XX][XX] > 0 && box[YY][YY] > 0 && box[ZZ][ZZ] > 0)
201 else if (box[XX][XX] > 0 && box[YY][YY] > 0 && box[ZZ][ZZ] == 0)
205 else if (box[XX][XX] == 0 && box[YY][YY] == 0 && box[ZZ][ZZ] == 0)
213 fprintf(stderr, "WARNING: Unsupported box diagonal %f %f %f, "
214 "will not use periodic boundary conditions\n\n",
215 box[XX][XX], box[YY][YY], box[ZZ][ZZ]);
223 fprintf(debug, "Guessed pbc = %s from the box matrix\n", epbc_names[ePBC]);
229 static int correct_box_elem(FILE *fplog, int step, tensor box, int v, int d)
231 int shift, maxshift = 10;
235 /* correct elem d of vector v with vector d */
236 while (box[v][d] > BOX_MARGIN_CORRECT*0.5*box[d][d])
240 fprintf(fplog, "Step %d: correcting invalid box:\n", step);
241 pr_rvecs(fplog, 0, "old box", box, DIM);
243 rvec_dec(box[v], box[d]);
247 pr_rvecs(fplog, 0, "new box", box, DIM);
249 if (shift <= -maxshift)
252 "Box was shifted at least %d times. Please see log-file.",
256 while (box[v][d] < -BOX_MARGIN_CORRECT*0.5*box[d][d])
260 fprintf(fplog, "Step %d: correcting invalid box:\n", step);
261 pr_rvecs(fplog, 0, "old box", box, DIM);
263 rvec_inc(box[v], box[d]);
267 pr_rvecs(fplog, 0, "new box", box, DIM);
269 if (shift >= maxshift)
272 "Box was shifted at least %d times. Please see log-file.",
280 gmx_bool correct_box(FILE *fplog, int step, tensor box, t_graph *graph)
285 /* check if the box still obeys the restrictions, if not, correct it */
286 zy = correct_box_elem(fplog, step, box, ZZ, YY);
287 zx = correct_box_elem(fplog, step, box, ZZ, XX);
288 yx = correct_box_elem(fplog, step, box, YY, XX);
290 bCorrected = (zy || zx || yx);
292 if (bCorrected && graph)
294 /* correct the graph */
295 for (i = graph->at_start; i < graph->at_end; i++)
297 graph->ishift[i][YY] -= graph->ishift[i][ZZ]*zy;
298 graph->ishift[i][XX] -= graph->ishift[i][ZZ]*zx;
299 graph->ishift[i][XX] -= graph->ishift[i][YY]*yx;
306 int ndof_com(t_inputrec *ir)
317 n = (ir->nwall == 0 ? 3 : 2);
323 gmx_incons("Unknown pbc in calc_nrdf");
329 static void low_set_pbc(t_pbc *pbc, int ePBC, ivec *dd_nc, matrix box)
331 int order[5] = {0, -1, 1, -2, 2};
332 int ii, jj, kk, i, j, k, d, dd, jc, kc, npbcdim, shift;
334 real d2old, d2new, d2new_c;
340 pbc->ndim_ePBC = ePBC2npbcdim(ePBC);
342 copy_mat(box, pbc->box);
343 pbc->bLimitDistance = FALSE;
344 pbc->max_cutoff2 = 0;
347 for (i = 0; (i < DIM); i++)
349 pbc->fbox_diag[i] = box[i][i];
350 pbc->hbox_diag[i] = pbc->fbox_diag[i]*0.5;
351 pbc->mhbox_diag[i] = -pbc->hbox_diag[i];
354 ptr = check_box(ePBC, box);
355 if (ePBC == epbcNONE)
357 pbc->ePBCDX = epbcdxNOPBC;
361 fprintf(stderr, "Warning: %s\n", ptr);
362 pr_rvecs(stderr, 0, " Box", box, DIM);
363 fprintf(stderr, " Can not fix pbc.\n");
364 pbc->ePBCDX = epbcdxUNSUPPORTED;
365 pbc->bLimitDistance = TRUE;
366 pbc->limit_distance2 = 0;
370 if (ePBC == epbcSCREW && dd_nc)
372 /* This combinated should never appear here */
373 gmx_incons("low_set_pbc called with screw pbc and dd_nc != NULL");
377 for (i = 0; i < DIM; i++)
379 if ((dd_nc && (*dd_nc)[i] > 1) || (ePBC == epbcXY && i == ZZ))
392 /* 1D pbc is not an mdp option and it is therefore only used
393 * with single shifts.
395 pbc->ePBCDX = epbcdx1D_RECT;
396 for (i = 0; i < DIM; i++)
403 assert(pbc->dim < DIM);
404 for (i = 0; i < pbc->dim; i++)
406 if (pbc->box[pbc->dim][i] != 0)
408 pbc->ePBCDX = epbcdx1D_TRIC;
413 pbc->ePBCDX = epbcdx2D_RECT;
414 for (i = 0; i < DIM; i++)
421 for (i = 0; i < DIM; i++)
425 for (j = 0; j < i; j++)
427 if (pbc->box[i][j] != 0)
429 pbc->ePBCDX = epbcdx2D_TRIC;
436 if (ePBC != epbcSCREW)
440 pbc->ePBCDX = epbcdxTRICLINIC;
444 pbc->ePBCDX = epbcdxRECTANGULAR;
449 pbc->ePBCDX = (box[ZZ][YY] == 0 ? epbcdxSCREW_RECT : epbcdxSCREW_TRIC);
450 if (pbc->ePBCDX == epbcdxSCREW_TRIC)
453 "Screw pbc is not yet implemented for triclinic boxes.\n"
454 "Can not fix pbc.\n");
455 pbc->ePBCDX = epbcdxUNSUPPORTED;
460 gmx_fatal(FARGS, "Incorrect number of pbc dimensions with DD: %d",
463 pbc->max_cutoff2 = max_cutoff2(ePBC, box);
465 if (pbc->ePBCDX == epbcdxTRICLINIC ||
466 pbc->ePBCDX == epbcdx2D_TRIC ||
467 pbc->ePBCDX == epbcdxSCREW_TRIC)
471 pr_rvecs(debug, 0, "Box", box, DIM);
472 fprintf(debug, "max cutoff %.3f\n", sqrt(pbc->max_cutoff2));
475 /* We will only use single shifts, but we will check a few
476 * more shifts to see if there is a limiting distance
477 * above which we can not be sure of the correct distance.
479 for (kk = 0; kk < 5; kk++)
482 if (!bPBC[ZZ] && k != 0)
486 for (jj = 0; jj < 5; jj++)
489 if (!bPBC[YY] && j != 0)
493 for (ii = 0; ii < 3; ii++)
496 if (!bPBC[XX] && i != 0)
500 /* A shift is only useful when it is trilinic */
501 if (j != 0 || k != 0)
505 for (d = 0; d < DIM; d++)
507 trial[d] = i*box[XX][d] + j*box[YY][d] + k*box[ZZ][d];
508 /* Choose the vector within the brick around 0,0,0 that
509 * will become the shortest due to shift try.
520 pos[d] = min( pbc->hbox_diag[d], -trial[d]);
524 pos[d] = max(-pbc->hbox_diag[d], -trial[d]);
527 d2old += sqr(pos[d]);
528 d2new += sqr(pos[d] + trial[d]);
530 if (BOX_MARGIN*d2new < d2old)
532 if (j < -1 || j > 1 || k < -1 || k > 1)
534 /* Check if there is a single shift vector
535 * that decreases this distance even more.
548 for (d = 0; d < DIM; d++)
550 d2new_c += sqr(pos[d] + trial[d]
551 - jc*box[YY][d] - kc*box[ZZ][d]);
553 if (d2new_c > BOX_MARGIN*d2new)
555 /* Reject this shift vector, as there is no a priori limit
556 * to the number of shifts that decrease distances.
558 if (!pbc->bLimitDistance || d2new < pbc->limit_distance2)
560 pbc->limit_distance2 = d2new;
562 pbc->bLimitDistance = TRUE;
567 /* Check if shifts with one box vector less do better */
569 for (dd = 0; dd < DIM; dd++)
571 shift = (dd == 0 ? i : (dd == 1 ? j : k));
575 for (d = 0; d < DIM; d++)
577 d2new_c += sqr(pos[d] + trial[d] - shift*box[dd][d]);
579 if (d2new_c <= BOX_MARGIN*d2new)
587 /* Accept this shift vector. */
588 if (pbc->ntric_vec >= MAX_NTRICVEC)
590 fprintf(stderr, "\nWARNING: Found more than %d triclinic correction vectors, ignoring some.\n"
591 " There is probably something wrong with your box.\n", MAX_NTRICVEC);
592 pr_rvecs(stderr, 0, " Box", box, DIM);
596 copy_rvec(trial, pbc->tric_vec[pbc->ntric_vec]);
597 pbc->tric_shift[pbc->ntric_vec][XX] = i;
598 pbc->tric_shift[pbc->ntric_vec][YY] = j;
599 pbc->tric_shift[pbc->ntric_vec][ZZ] = k;
606 fprintf(debug, " tricvec %2d = %2d %2d %2d %5.2f %5.2f %5.2f %5.2f %5.2f %5.2f %5.2f %5.2f\n",
607 pbc->ntric_vec, i, j, k,
608 sqrt(d2old), sqrt(d2new),
609 trial[XX], trial[YY], trial[ZZ],
610 pos[XX], pos[YY], pos[ZZ]);
621 void set_pbc(t_pbc *pbc, int ePBC, matrix box)
625 ePBC = guess_ePBC(box);
628 low_set_pbc(pbc, ePBC, NULL, box);
631 t_pbc *set_pbc_dd(t_pbc *pbc, int ePBC,
632 gmx_domdec_t *dd, gmx_bool bSingleDir, matrix box)
643 if (ePBC == epbcSCREW && dd->nc[XX] > 1)
645 /* The rotation has been taken care of during coordinate communication */
649 for (i = 0; i < DIM; i++)
651 if (dd->nc[i] <= (bSingleDir ? 1 : 2))
654 if (!(ePBC == epbcXY && i == ZZ))
668 low_set_pbc(pbc, ePBC, npbcdim < DIM ? &nc2 : NULL, box);
671 return (npbcdim > 0 ? pbc : NULL);
674 void pbc_dx(const t_pbc *pbc, const rvec x1, const rvec x2, rvec dx)
677 rvec dx_start, trial;
681 rvec_sub(x1, x2, dx);
685 case epbcdxRECTANGULAR:
686 for (i = 0; i < DIM; i++)
688 while (dx[i] > pbc->hbox_diag[i])
690 dx[i] -= pbc->fbox_diag[i];
692 while (dx[i] <= pbc->mhbox_diag[i])
694 dx[i] += pbc->fbox_diag[i];
698 case epbcdxTRICLINIC:
699 for (i = DIM-1; i >= 0; i--)
701 while (dx[i] > pbc->hbox_diag[i])
703 for (j = i; j >= 0; j--)
705 dx[j] -= pbc->box[i][j];
708 while (dx[i] <= pbc->mhbox_diag[i])
710 for (j = i; j >= 0; j--)
712 dx[j] += pbc->box[i][j];
716 /* dx is the distance in a rectangular box */
718 if (d2min > pbc->max_cutoff2)
720 copy_rvec(dx, dx_start);
722 /* Now try all possible shifts, when the distance is within max_cutoff
723 * it must be the shortest possible distance.
726 while ((d2min > pbc->max_cutoff2) && (i < pbc->ntric_vec))
728 rvec_add(dx_start, pbc->tric_vec[i], trial);
729 d2trial = norm2(trial);
732 copy_rvec(trial, dx);
740 for (i = 0; i < DIM; i++)
744 while (dx[i] > pbc->hbox_diag[i])
746 dx[i] -= pbc->fbox_diag[i];
748 while (dx[i] <= pbc->mhbox_diag[i])
750 dx[i] += pbc->fbox_diag[i];
757 for (i = DIM-1; i >= 0; i--)
761 while (dx[i] > pbc->hbox_diag[i])
763 for (j = i; j >= 0; j--)
765 dx[j] -= pbc->box[i][j];
768 while (dx[i] <= pbc->mhbox_diag[i])
770 for (j = i; j >= 0; j--)
772 dx[j] += pbc->box[i][j];
775 d2min += dx[i]*dx[i];
778 if (d2min > pbc->max_cutoff2)
780 copy_rvec(dx, dx_start);
782 /* Now try all possible shifts, when the distance is within max_cutoff
783 * it must be the shortest possible distance.
786 while ((d2min > pbc->max_cutoff2) && (i < pbc->ntric_vec))
788 rvec_add(dx_start, pbc->tric_vec[i], trial);
790 for (j = 0; j < DIM; j++)
794 d2trial += trial[j]*trial[j];
799 copy_rvec(trial, dx);
806 case epbcdxSCREW_RECT:
807 /* The shift definition requires x first */
809 while (dx[XX] > pbc->hbox_diag[XX])
811 dx[XX] -= pbc->fbox_diag[XX];
814 while (dx[XX] <= pbc->mhbox_diag[XX])
816 dx[XX] += pbc->fbox_diag[YY];
821 /* Rotate around the x-axis in the middle of the box */
822 dx[YY] = pbc->box[YY][YY] - x1[YY] - x2[YY];
823 dx[ZZ] = pbc->box[ZZ][ZZ] - x1[ZZ] - x2[ZZ];
825 /* Normal pbc for y and z */
826 for (i = YY; i <= ZZ; i++)
828 while (dx[i] > pbc->hbox_diag[i])
830 dx[i] -= pbc->fbox_diag[i];
832 while (dx[i] <= pbc->mhbox_diag[i])
834 dx[i] += pbc->fbox_diag[i];
839 case epbcdxUNSUPPORTED:
842 gmx_fatal(FARGS, "Internal error in pbc_dx, set_pbc has not been called");
847 int pbc_dx_aiuc(const t_pbc *pbc, const rvec x1, const rvec x2, rvec dx)
850 rvec dx_start, trial;
852 ivec ishift, ishift_start;
854 rvec_sub(x1, x2, dx);
859 case epbcdxRECTANGULAR:
860 for (i = 0; i < DIM; i++)
862 if (dx[i] > pbc->hbox_diag[i])
864 dx[i] -= pbc->fbox_diag[i];
867 else if (dx[i] <= pbc->mhbox_diag[i])
869 dx[i] += pbc->fbox_diag[i];
874 case epbcdxTRICLINIC:
875 /* For triclinic boxes the performance difference between
876 * if/else and two while loops is negligible.
877 * However, the while version can cause extreme delays
878 * before a simulation crashes due to large forces which
879 * can cause unlimited displacements.
880 * Also allowing multiple shifts would index fshift beyond bounds.
882 for (i = DIM-1; i >= 1; i--)
884 if (dx[i] > pbc->hbox_diag[i])
886 for (j = i; j >= 0; j--)
888 dx[j] -= pbc->box[i][j];
892 else if (dx[i] <= pbc->mhbox_diag[i])
894 for (j = i; j >= 0; j--)
896 dx[j] += pbc->box[i][j];
901 /* Allow 2 shifts in x */
902 if (dx[XX] > pbc->hbox_diag[XX])
904 dx[XX] -= pbc->fbox_diag[XX];
906 if (dx[XX] > pbc->hbox_diag[XX])
908 dx[XX] -= pbc->fbox_diag[XX];
912 else if (dx[XX] <= pbc->mhbox_diag[XX])
914 dx[XX] += pbc->fbox_diag[XX];
916 if (dx[XX] <= pbc->mhbox_diag[XX])
918 dx[XX] += pbc->fbox_diag[XX];
922 /* dx is the distance in a rectangular box */
924 if (d2min > pbc->max_cutoff2)
926 copy_rvec(dx, dx_start);
927 copy_ivec(ishift, ishift_start);
929 /* Now try all possible shifts, when the distance is within max_cutoff
930 * it must be the shortest possible distance.
933 while ((d2min > pbc->max_cutoff2) && (i < pbc->ntric_vec))
935 rvec_add(dx_start, pbc->tric_vec[i], trial);
936 d2trial = norm2(trial);
939 copy_rvec(trial, dx);
940 ivec_add(ishift_start, pbc->tric_shift[i], ishift);
948 for (i = 0; i < DIM; i++)
952 if (dx[i] > pbc->hbox_diag[i])
954 dx[i] -= pbc->fbox_diag[i];
957 else if (dx[i] <= pbc->mhbox_diag[i])
959 dx[i] += pbc->fbox_diag[i];
967 for (i = DIM-1; i >= 1; i--)
971 if (dx[i] > pbc->hbox_diag[i])
973 for (j = i; j >= 0; j--)
975 dx[j] -= pbc->box[i][j];
979 else if (dx[i] <= pbc->mhbox_diag[i])
981 for (j = i; j >= 0; j--)
983 dx[j] += pbc->box[i][j];
987 d2min += dx[i]*dx[i];
992 /* Allow 2 shifts in x */
993 if (dx[XX] > pbc->hbox_diag[XX])
995 dx[XX] -= pbc->fbox_diag[XX];
997 if (dx[XX] > pbc->hbox_diag[XX])
999 dx[XX] -= pbc->fbox_diag[XX];
1003 else if (dx[XX] <= pbc->mhbox_diag[XX])
1005 dx[XX] += pbc->fbox_diag[XX];
1007 if (dx[XX] <= pbc->mhbox_diag[XX])
1009 dx[XX] += pbc->fbox_diag[XX];
1013 d2min += dx[XX]*dx[XX];
1015 if (d2min > pbc->max_cutoff2)
1017 copy_rvec(dx, dx_start);
1018 copy_ivec(ishift, ishift_start);
1019 /* Now try all possible shifts, when the distance is within max_cutoff
1020 * it must be the shortest possible distance.
1023 while ((d2min > pbc->max_cutoff2) && (i < pbc->ntric_vec))
1025 rvec_add(dx_start, pbc->tric_vec[i], trial);
1027 for (j = 0; j < DIM; j++)
1031 d2trial += trial[j]*trial[j];
1034 if (d2trial < d2min)
1036 copy_rvec(trial, dx);
1037 ivec_add(ishift_start, pbc->tric_shift[i], ishift);
1046 if (dx[i] > pbc->hbox_diag[i])
1048 dx[i] -= pbc->fbox_diag[i];
1051 else if (dx[i] <= pbc->mhbox_diag[i])
1053 dx[i] += pbc->fbox_diag[i];
1059 if (dx[i] > pbc->hbox_diag[i])
1061 rvec_dec(dx, pbc->box[i]);
1064 else if (dx[i] <= pbc->mhbox_diag[i])
1066 rvec_inc(dx, pbc->box[i]);
1070 case epbcdxSCREW_RECT:
1071 /* The shift definition requires x first */
1072 if (dx[XX] > pbc->hbox_diag[XX])
1074 dx[XX] -= pbc->fbox_diag[XX];
1077 else if (dx[XX] <= pbc->mhbox_diag[XX])
1079 dx[XX] += pbc->fbox_diag[XX];
1082 if (ishift[XX] == 1 || ishift[XX] == -1)
1084 /* Rotate around the x-axis in the middle of the box */
1085 dx[YY] = pbc->box[YY][YY] - x1[YY] - x2[YY];
1086 dx[ZZ] = pbc->box[ZZ][ZZ] - x1[ZZ] - x2[ZZ];
1088 /* Normal pbc for y and z */
1089 for (i = YY; i <= ZZ; i++)
1091 if (dx[i] > pbc->hbox_diag[i])
1093 dx[i] -= pbc->fbox_diag[i];
1096 else if (dx[i] <= pbc->mhbox_diag[i])
1098 dx[i] += pbc->fbox_diag[i];
1104 case epbcdxUNSUPPORTED:
1107 gmx_fatal(FARGS, "Internal error in pbc_dx_aiuc, set_pbc_dd or set_pbc has not been called");
1111 is = IVEC2IS(ishift);
1114 range_check_mesg(is, 0, SHIFTS, "PBC shift vector index range check.");
1120 void pbc_dx_d(const t_pbc *pbc, const dvec x1, const dvec x2, dvec dx)
1123 dvec dx_start, trial;
1124 double d2min, d2trial;
1127 dvec_sub(x1, x2, dx);
1129 switch (pbc->ePBCDX)
1131 case epbcdxRECTANGULAR:
1133 for (i = 0; i < DIM; i++)
1137 while (dx[i] > pbc->hbox_diag[i])
1139 dx[i] -= pbc->fbox_diag[i];
1141 while (dx[i] <= pbc->mhbox_diag[i])
1143 dx[i] += pbc->fbox_diag[i];
1148 case epbcdxTRICLINIC:
1151 for (i = DIM-1; i >= 0; i--)
1155 while (dx[i] > pbc->hbox_diag[i])
1157 for (j = i; j >= 0; j--)
1159 dx[j] -= pbc->box[i][j];
1162 while (dx[i] <= pbc->mhbox_diag[i])
1164 for (j = i; j >= 0; j--)
1166 dx[j] += pbc->box[i][j];
1169 d2min += dx[i]*dx[i];
1172 if (d2min > pbc->max_cutoff2)
1174 copy_dvec(dx, dx_start);
1175 /* Now try all possible shifts, when the distance is within max_cutoff
1176 * it must be the shortest possible distance.
1179 while ((d2min > pbc->max_cutoff2) && (i < pbc->ntric_vec))
1181 for (j = 0; j < DIM; j++)
1183 trial[j] = dx_start[j] + pbc->tric_vec[i][j];
1186 for (j = 0; j < DIM; j++)
1190 d2trial += trial[j]*trial[j];
1193 if (d2trial < d2min)
1195 copy_dvec(trial, dx);
1202 case epbcdxSCREW_RECT:
1203 /* The shift definition requires x first */
1205 while (dx[XX] > pbc->hbox_diag[XX])
1207 dx[XX] -= pbc->fbox_diag[XX];
1210 while (dx[XX] <= pbc->mhbox_diag[XX])
1212 dx[XX] += pbc->fbox_diag[YY];
1217 /* Rotate around the x-axis in the middle of the box */
1218 dx[YY] = pbc->box[YY][YY] - x1[YY] - x2[YY];
1219 dx[ZZ] = pbc->box[ZZ][ZZ] - x1[ZZ] - x2[ZZ];
1221 /* Normal pbc for y and z */
1222 for (i = YY; i <= ZZ; i++)
1224 while (dx[i] > pbc->hbox_diag[i])
1226 dx[i] -= pbc->fbox_diag[i];
1228 while (dx[i] <= pbc->mhbox_diag[i])
1230 dx[i] += pbc->fbox_diag[i];
1235 case epbcdxUNSUPPORTED:
1238 gmx_fatal(FARGS, "Internal error in pbc_dx, set_pbc has not been called");
1243 gmx_bool image_rect(ivec xi, ivec xj, ivec box_size, real rlong2, int *shift, real *r2)
1251 for (m = 0; (m < DIM); m++)
1283 gmx_bool image_cylindric(ivec xi, ivec xj, ivec box_size, real rlong2,
1284 int *shift, real *r2)
1292 for (m = 0; (m < DIM); m++)
1328 void calc_shifts(matrix box, rvec shift_vec[])
1330 int k, l, m, d, n, test;
1333 for (m = -D_BOX_Z; m <= D_BOX_Z; m++)
1335 for (l = -D_BOX_Y; l <= D_BOX_Y; l++)
1337 for (k = -D_BOX_X; k <= D_BOX_X; k++, n++)
1339 test = XYZ2IS(k, l, m);
1342 gmx_incons("inconsistent shift numbering");
1344 for (d = 0; d < DIM; d++)
1346 shift_vec[n][d] = k*box[XX][d] + l*box[YY][d] + m*box[ZZ][d];
1353 void calc_box_center(int ecenter, matrix box, rvec box_center)
1357 clear_rvec(box_center);
1361 for (m = 0; (m < DIM); m++)
1363 for (d = 0; d < DIM; d++)
1365 box_center[d] += 0.5*box[m][d];
1370 for (d = 0; d < DIM; d++)
1372 box_center[d] = 0.5*box[d][d];
1378 gmx_fatal(FARGS, "Unsupported value %d for ecenter", ecenter);
1382 void calc_triclinic_images(matrix box, rvec img[])
1386 /* Calculate 3 adjacent images in the xy-plane */
1387 copy_rvec(box[0], img[0]);
1388 copy_rvec(box[1], img[1]);
1391 svmul(-1, img[1], img[1]);
1393 rvec_sub(img[1], img[0], img[2]);
1395 /* Get the next 3 in the xy-plane as mirror images */
1396 for (i = 0; i < 3; i++)
1398 svmul(-1, img[i], img[3+i]);
1401 /* Calculate the first 4 out of xy-plane images */
1402 copy_rvec(box[2], img[6]);
1405 svmul(-1, img[6], img[6]);
1407 for (i = 0; i < 3; i++)
1409 rvec_add(img[6], img[i+1], img[7+i]);
1412 /* Mirror the last 4 from the previous in opposite rotation */
1413 for (i = 0; i < 4; i++)
1415 svmul(-1, img[6 + (2+i) % 4], img[10+i]);
1419 void calc_compact_unitcell_vertices(int ecenter, matrix box, rvec vert[])
1421 rvec img[NTRICIMG], box_center;
1422 int n, i, j, tmp[4], d;
1424 calc_triclinic_images(box, img);
1427 for (i = 2; i <= 5; i += 3)
1440 for (j = 0; j < 4; j++)
1442 for (d = 0; d < DIM; d++)
1444 vert[n][d] = img[i][d]+img[tmp[j]][d]+img[tmp[(j+1)%4]][d];
1449 for (i = 7; i <= 13; i += 6)
1462 for (j = 0; j < 4; j++)
1464 for (d = 0; d < DIM; d++)
1466 vert[n][d] = img[i][d]+img[tmp[j]][d]+img[tmp[(j+1)%4]][d];
1471 for (i = 9; i <= 11; i += 2)
1491 for (j = 0; j < 4; j++)
1493 for (d = 0; d < DIM; d++)
1495 vert[n][d] = img[i][d]+img[tmp[j]][d]+img[tmp[(j+1)%4]][d];
1501 calc_box_center(ecenter, box, box_center);
1502 for (i = 0; i < NCUCVERT; i++)
1504 for (d = 0; d < DIM; d++)
1506 vert[i][d] = vert[i][d]*0.25+box_center[d];
1511 int *compact_unitcell_edges()
1513 /* this is an index in vert[] (see calc_box_vertices) */
1514 /*static int edge[NCUCEDGE*2];*/
1516 static const int hexcon[24] = {
1517 0, 9, 1, 19, 2, 15, 3, 21,
1518 4, 17, 5, 11, 6, 23, 7, 13,
1519 8, 20, 10, 18, 12, 16, 14, 22
1522 gmx_bool bFirst = TRUE;
1524 snew(edge, NCUCEDGE*2);
1529 for (i = 0; i < 6; i++)
1531 for (j = 0; j < 4; j++)
1533 edge[e++] = 4*i + j;
1534 edge[e++] = 4*i + (j+1) % 4;
1537 for (i = 0; i < 12*2; i++)
1539 edge[e++] = hexcon[i];
1548 void put_atoms_in_box_omp(int ePBC, matrix box, int natoms, rvec x[])
1551 nth = gmx_omp_nthreads_get(emntDefault);
1553 #pragma omp parallel for num_threads(nth) schedule(static)
1554 for (t = 0; t < nth; t++)
1558 offset = (natoms*t )/nth;
1559 len = (natoms*(t + 1))/nth - offset;
1560 put_atoms_in_box(ePBC, box, len, x + offset);
1564 void put_atoms_in_box(int ePBC, matrix box, int natoms, rvec x[])
1566 int npbcdim, i, m, d;
1568 if (ePBC == epbcSCREW)
1570 gmx_fatal(FARGS, "Sorry, %s pbc is not yet supported", epbc_names[ePBC]);
1584 for (i = 0; (i < natoms); i++)
1586 for (m = npbcdim-1; m >= 0; m--)
1590 for (d = 0; d <= m; d++)
1592 x[i][d] += box[m][d];
1595 while (x[i][m] >= box[m][m])
1597 for (d = 0; d <= m; d++)
1599 x[i][d] -= box[m][d];
1607 for (i = 0; i < natoms; i++)
1609 for (d = 0; d < npbcdim; d++)
1613 x[i][d] += box[d][d];
1615 while (x[i][d] >= box[d][d])
1617 x[i][d] -= box[d][d];
1624 void put_atoms_in_triclinic_unitcell(int ecenter, matrix box,
1625 int natoms, rvec x[])
1627 rvec box_center, shift_center;
1628 real shm01, shm02, shm12, shift;
1631 calc_box_center(ecenter, box, box_center);
1633 /* The product of matrix shm with a coordinate gives the shift vector
1634 which is required determine the periodic cell position */
1635 shm01 = box[1][0]/box[1][1];
1636 shm02 = (box[1][1]*box[2][0] - box[2][1]*box[1][0])/(box[1][1]*box[2][2]);
1637 shm12 = box[2][1]/box[2][2];
1639 clear_rvec(shift_center);
1640 for (d = 0; d < DIM; d++)
1642 rvec_inc(shift_center, box[d]);
1644 svmul(0.5, shift_center, shift_center);
1645 rvec_sub(box_center, shift_center, shift_center);
1647 shift_center[0] = shm01*shift_center[1] + shm02*shift_center[2];
1648 shift_center[1] = shm12*shift_center[2];
1649 shift_center[2] = 0;
1651 for (i = 0; (i < natoms); i++)
1653 for (m = DIM-1; m >= 0; m--)
1655 shift = shift_center[m];
1658 shift += shm01*x[i][1] + shm02*x[i][2];
1662 shift += shm12*x[i][2];
1664 while (x[i][m]-shift < 0)
1666 for (d = 0; d <= m; d++)
1668 x[i][d] += box[m][d];
1671 while (x[i][m]-shift >= box[m][m])
1673 for (d = 0; d <= m; d++)
1675 x[i][d] -= box[m][d];
1683 put_atoms_in_compact_unitcell(int ePBC, int ecenter, matrix box,
1684 int natoms, rvec x[])
1687 rvec box_center, dx;
1690 set_pbc(&pbc, ePBC, box);
1691 calc_box_center(ecenter, box, box_center);
1692 for (i = 0; i < natoms; i++)
1694 pbc_dx(&pbc, x[i], box_center, dx);
1695 rvec_add(box_center, dx, x[i]);
1698 return pbc.bLimitDistance ?
1699 "WARNING: Could not put all atoms in the compact unitcell\n"