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.
44 #include "gromacs/legacyheaders/gmx_omp_nthreads.h"
45 #include "gromacs/legacyheaders/macros.h"
46 #include "gromacs/legacyheaders/names.h"
47 #include "gromacs/legacyheaders/txtdump.h"
48 #include "gromacs/legacyheaders/types/commrec.h"
49 #include "gromacs/legacyheaders/types/inputrec.h"
50 #include "gromacs/math/utilities.h"
51 #include "gromacs/math/vec.h"
52 #include "gromacs/pbcutil/ishift.h"
53 #include "gromacs/pbcutil/mshift.h"
54 #include "gromacs/utility/fatalerror.h"
55 #include "gromacs/utility/smalloc.h"
57 /* Skip 0 so we have more chance of detecting if we forgot to call set_pbc. */
59 epbcdxRECTANGULAR = 1, epbcdxTRICLINIC,
60 epbcdx2D_RECT, epbcdx2D_TRIC,
61 epbcdx1D_RECT, epbcdx1D_TRIC,
62 epbcdxSCREW_RECT, epbcdxSCREW_TRIC,
63 epbcdxNOPBC, epbcdxUNSUPPORTED
66 /* Margin factor for error message and correction if the box is too skewed */
67 #define BOX_MARGIN 1.0010
68 #define BOX_MARGIN_CORRECT 1.0005
70 int ePBC2npbcdim(int ePBC)
76 case epbcXYZ: npbcdim = 3; break;
77 case epbcXY: npbcdim = 2; break;
78 case epbcSCREW: npbcdim = 3; break;
79 case epbcNONE: npbcdim = 0; break;
80 default: gmx_fatal(FARGS, "Unknown ePBC=%d in ePBC2npbcdim", ePBC);
86 int inputrec2nboundeddim(t_inputrec *ir)
88 if (ir->nwall == 2 && ir->ePBC == epbcXY)
94 return ePBC2npbcdim(ir->ePBC);
98 void dump_pbc(FILE *fp, t_pbc *pbc)
102 fprintf(fp, "ePBCDX = %d\n", pbc->ePBCDX);
103 pr_rvecs(fp, 0, "box", pbc->box, DIM);
104 pr_rvecs(fp, 0, "fbox_diag", &pbc->fbox_diag, 1);
105 pr_rvecs(fp, 0, "hbox_diag", &pbc->hbox_diag, 1);
106 pr_rvecs(fp, 0, "mhbox_diag", &pbc->mhbox_diag, 1);
107 rvec_add(pbc->hbox_diag, pbc->mhbox_diag, sum_box);
108 pr_rvecs(fp, 0, "sum of the above two", &sum_box, 1);
109 fprintf(fp, "max_cutoff2 = %g\n", pbc->max_cutoff2);
110 fprintf(fp, "bLimitDistance = %s\n", EBOOL(pbc->bLimitDistance));
111 fprintf(fp, "limit_distance2 = %g\n", pbc->limit_distance2);
112 fprintf(fp, "ntric_vec = %d\n", pbc->ntric_vec);
113 if (pbc->ntric_vec > 0)
115 pr_ivecs(fp, 0, "tric_shift", pbc->tric_shift, pbc->ntric_vec, FALSE);
116 pr_rvecs(fp, 0, "tric_vec", pbc->tric_vec, pbc->ntric_vec);
120 const char *check_box(int ePBC, matrix box)
126 ePBC = guess_ePBC(box);
129 if (ePBC == epbcNONE)
134 if ((box[XX][YY] != 0) || (box[XX][ZZ] != 0) || (box[YY][ZZ] != 0))
136 ptr = "Only triclinic boxes with the first vector parallel to the x-axis and the second vector in the xy-plane are supported.";
138 else if (ePBC == epbcSCREW && (box[YY][XX] != 0 || box[ZZ][XX] != 0))
140 ptr = "The unit cell can not have off-diagonal x-components with screw pbc";
142 else if (fabs(box[YY][XX]) > BOX_MARGIN*0.5*box[XX][XX] ||
144 (fabs(box[ZZ][XX]) > BOX_MARGIN*0.5*box[XX][XX] ||
145 fabs(box[ZZ][YY]) > BOX_MARGIN*0.5*box[YY][YY])))
147 ptr = "Triclinic box is too skewed.";
157 real max_cutoff2(int ePBC, matrix box)
159 real min_hv2, min_ss;
161 /* Physical limitation of the cut-off
162 * by half the length of the shortest box vector.
164 min_hv2 = min(0.25*norm2(box[XX]), 0.25*norm2(box[YY]));
167 min_hv2 = min(min_hv2, 0.25*norm2(box[ZZ]));
170 /* Limitation to the smallest diagonal element due to optimizations:
171 * checking only linear combinations of single box-vectors (2 in x)
172 * in the grid search and pbc_dx is a lot faster
173 * than checking all possible combinations.
177 min_ss = min(box[XX][XX], box[YY][YY]);
181 min_ss = min(box[XX][XX], min(box[YY][YY]-fabs(box[ZZ][YY]), box[ZZ][ZZ]));
184 return min(min_hv2, min_ss*min_ss);
187 /* this one is mostly harmless... */
188 static gmx_bool bWarnedGuess = FALSE;
190 int guess_ePBC(matrix box)
194 if (box[XX][XX] > 0 && box[YY][YY] > 0 && box[ZZ][ZZ] > 0)
198 else if (box[XX][XX] > 0 && box[YY][YY] > 0 && box[ZZ][ZZ] == 0)
202 else if (box[XX][XX] == 0 && box[YY][YY] == 0 && box[ZZ][ZZ] == 0)
210 fprintf(stderr, "WARNING: Unsupported box diagonal %f %f %f, "
211 "will not use periodic boundary conditions\n\n",
212 box[XX][XX], box[YY][YY], box[ZZ][ZZ]);
220 fprintf(debug, "Guessed pbc = %s from the box matrix\n", epbc_names[ePBC]);
226 static int correct_box_elem(FILE *fplog, int step, tensor box, int v, int d)
228 int shift, maxshift = 10;
232 /* correct elem d of vector v with vector d */
233 while (box[v][d] > BOX_MARGIN_CORRECT*0.5*box[d][d])
237 fprintf(fplog, "Step %d: correcting invalid box:\n", step);
238 pr_rvecs(fplog, 0, "old box", box, DIM);
240 rvec_dec(box[v], box[d]);
244 pr_rvecs(fplog, 0, "new box", box, DIM);
246 if (shift <= -maxshift)
249 "Box was shifted at least %d times. Please see log-file.",
253 while (box[v][d] < -BOX_MARGIN_CORRECT*0.5*box[d][d])
257 fprintf(fplog, "Step %d: correcting invalid box:\n", step);
258 pr_rvecs(fplog, 0, "old box", box, DIM);
260 rvec_inc(box[v], box[d]);
264 pr_rvecs(fplog, 0, "new box", box, DIM);
266 if (shift >= maxshift)
269 "Box was shifted at least %d times. Please see log-file.",
277 gmx_bool correct_box(FILE *fplog, int step, tensor box, t_graph *graph)
282 /* check if the box still obeys the restrictions, if not, correct it */
283 zy = correct_box_elem(fplog, step, box, ZZ, YY);
284 zx = correct_box_elem(fplog, step, box, ZZ, XX);
285 yx = correct_box_elem(fplog, step, box, YY, XX);
287 bCorrected = (zy || zx || yx);
289 if (bCorrected && graph)
291 /* correct the graph */
292 for (i = graph->at_start; i < graph->at_end; i++)
294 graph->ishift[i][YY] -= graph->ishift[i][ZZ]*zy;
295 graph->ishift[i][XX] -= graph->ishift[i][ZZ]*zx;
296 graph->ishift[i][XX] -= graph->ishift[i][YY]*yx;
303 int ndof_com(t_inputrec *ir)
314 n = (ir->nwall == 0 ? 3 : 2);
320 gmx_incons("Unknown pbc in calc_nrdf");
326 static void low_set_pbc(t_pbc *pbc, int ePBC, ivec *dd_nc, matrix box)
328 int order[5] = {0, -1, 1, -2, 2};
329 int ii, jj, kk, i, j, k, d, dd, jc, kc, npbcdim, shift;
331 real d2old, d2new, d2new_c;
337 pbc->ndim_ePBC = ePBC2npbcdim(ePBC);
339 copy_mat(box, pbc->box);
340 pbc->bLimitDistance = FALSE;
341 pbc->max_cutoff2 = 0;
344 for (i = 0; (i < DIM); i++)
346 pbc->fbox_diag[i] = box[i][i];
347 pbc->hbox_diag[i] = pbc->fbox_diag[i]*0.5;
348 pbc->mhbox_diag[i] = -pbc->hbox_diag[i];
351 ptr = check_box(ePBC, box);
352 if (ePBC == epbcNONE)
354 pbc->ePBCDX = epbcdxNOPBC;
358 fprintf(stderr, "Warning: %s\n", ptr);
359 pr_rvecs(stderr, 0, " Box", box, DIM);
360 fprintf(stderr, " Can not fix pbc.\n");
361 pbc->ePBCDX = epbcdxUNSUPPORTED;
362 pbc->bLimitDistance = TRUE;
363 pbc->limit_distance2 = 0;
367 if (ePBC == epbcSCREW && dd_nc)
369 /* This combinated should never appear here */
370 gmx_incons("low_set_pbc called with screw pbc and dd_nc != NULL");
374 for (i = 0; i < DIM; i++)
376 if ((dd_nc && (*dd_nc)[i] > 1) || (ePBC == epbcXY && i == ZZ))
389 /* 1D pbc is not an mdp option and it is therefore only used
390 * with single shifts.
392 pbc->ePBCDX = epbcdx1D_RECT;
393 for (i = 0; i < DIM; i++)
400 assert(pbc->dim < DIM);
401 for (i = 0; i < pbc->dim; i++)
403 if (pbc->box[pbc->dim][i] != 0)
405 pbc->ePBCDX = epbcdx1D_TRIC;
410 pbc->ePBCDX = epbcdx2D_RECT;
411 for (i = 0; i < DIM; i++)
418 for (i = 0; i < DIM; i++)
422 for (j = 0; j < i; j++)
424 if (pbc->box[i][j] != 0)
426 pbc->ePBCDX = epbcdx2D_TRIC;
433 if (ePBC != epbcSCREW)
437 pbc->ePBCDX = epbcdxTRICLINIC;
441 pbc->ePBCDX = epbcdxRECTANGULAR;
446 pbc->ePBCDX = (box[ZZ][YY] == 0 ? epbcdxSCREW_RECT : epbcdxSCREW_TRIC);
447 if (pbc->ePBCDX == epbcdxSCREW_TRIC)
450 "Screw pbc is not yet implemented for triclinic boxes.\n"
451 "Can not fix pbc.\n");
452 pbc->ePBCDX = epbcdxUNSUPPORTED;
457 gmx_fatal(FARGS, "Incorrect number of pbc dimensions with DD: %d",
460 pbc->max_cutoff2 = max_cutoff2(ePBC, box);
462 if (pbc->ePBCDX == epbcdxTRICLINIC ||
463 pbc->ePBCDX == epbcdx2D_TRIC ||
464 pbc->ePBCDX == epbcdxSCREW_TRIC)
468 pr_rvecs(debug, 0, "Box", box, DIM);
469 fprintf(debug, "max cutoff %.3f\n", sqrt(pbc->max_cutoff2));
472 /* We will only use single shifts, but we will check a few
473 * more shifts to see if there is a limiting distance
474 * above which we can not be sure of the correct distance.
476 for (kk = 0; kk < 5; kk++)
479 if (!bPBC[ZZ] && k != 0)
483 for (jj = 0; jj < 5; jj++)
486 if (!bPBC[YY] && j != 0)
490 for (ii = 0; ii < 3; ii++)
493 if (!bPBC[XX] && i != 0)
497 /* A shift is only useful when it is trilinic */
498 if (j != 0 || k != 0)
502 for (d = 0; d < DIM; d++)
504 trial[d] = i*box[XX][d] + j*box[YY][d] + k*box[ZZ][d];
505 /* Choose the vector within the brick around 0,0,0 that
506 * will become the shortest due to shift try.
517 pos[d] = min( pbc->hbox_diag[d], -trial[d]);
521 pos[d] = max(-pbc->hbox_diag[d], -trial[d]);
524 d2old += sqr(pos[d]);
525 d2new += sqr(pos[d] + trial[d]);
527 if (BOX_MARGIN*d2new < d2old)
529 if (j < -1 || j > 1 || k < -1 || k > 1)
531 /* Check if there is a single shift vector
532 * that decreases this distance even more.
545 for (d = 0; d < DIM; d++)
547 d2new_c += sqr(pos[d] + trial[d]
548 - jc*box[YY][d] - kc*box[ZZ][d]);
550 if (d2new_c > BOX_MARGIN*d2new)
552 /* Reject this shift vector, as there is no a priori limit
553 * to the number of shifts that decrease distances.
555 if (!pbc->bLimitDistance || d2new < pbc->limit_distance2)
557 pbc->limit_distance2 = d2new;
559 pbc->bLimitDistance = TRUE;
564 /* Check if shifts with one box vector less do better */
566 for (dd = 0; dd < DIM; dd++)
568 shift = (dd == 0 ? i : (dd == 1 ? j : k));
572 for (d = 0; d < DIM; d++)
574 d2new_c += sqr(pos[d] + trial[d] - shift*box[dd][d]);
576 if (d2new_c <= BOX_MARGIN*d2new)
584 /* Accept this shift vector. */
585 if (pbc->ntric_vec >= MAX_NTRICVEC)
587 fprintf(stderr, "\nWARNING: Found more than %d triclinic correction vectors, ignoring some.\n"
588 " There is probably something wrong with your box.\n", MAX_NTRICVEC);
589 pr_rvecs(stderr, 0, " Box", box, DIM);
593 copy_rvec(trial, pbc->tric_vec[pbc->ntric_vec]);
594 pbc->tric_shift[pbc->ntric_vec][XX] = i;
595 pbc->tric_shift[pbc->ntric_vec][YY] = j;
596 pbc->tric_shift[pbc->ntric_vec][ZZ] = k;
603 fprintf(debug, " tricvec %2d = %2d %2d %2d %5.2f %5.2f %5.2f %5.2f %5.2f %5.2f %5.2f %5.2f\n",
604 pbc->ntric_vec, i, j, k,
605 sqrt(d2old), sqrt(d2new),
606 trial[XX], trial[YY], trial[ZZ],
607 pos[XX], pos[YY], pos[ZZ]);
618 void set_pbc(t_pbc *pbc, int ePBC, matrix box)
622 ePBC = guess_ePBC(box);
625 low_set_pbc(pbc, ePBC, NULL, box);
628 t_pbc *set_pbc_dd(t_pbc *pbc, int ePBC,
629 gmx_domdec_t *dd, gmx_bool bSingleDir, matrix box)
640 if (ePBC == epbcSCREW && dd->nc[XX] > 1)
642 /* The rotation has been taken care of during coordinate communication */
646 for (i = 0; i < DIM; i++)
648 if (dd->nc[i] <= (bSingleDir ? 1 : 2))
651 if (!(ePBC == epbcXY && i == ZZ))
665 low_set_pbc(pbc, ePBC, npbcdim < DIM ? &nc2 : NULL, box);
668 return (npbcdim > 0 ? pbc : NULL);
671 void pbc_dx(const t_pbc *pbc, const rvec x1, const rvec x2, rvec dx)
674 rvec dx_start, trial;
678 rvec_sub(x1, x2, dx);
682 case epbcdxRECTANGULAR:
683 for (i = 0; i < DIM; i++)
685 while (dx[i] > pbc->hbox_diag[i])
687 dx[i] -= pbc->fbox_diag[i];
689 while (dx[i] <= pbc->mhbox_diag[i])
691 dx[i] += pbc->fbox_diag[i];
695 case epbcdxTRICLINIC:
696 for (i = DIM-1; i >= 0; i--)
698 while (dx[i] > pbc->hbox_diag[i])
700 for (j = i; j >= 0; j--)
702 dx[j] -= pbc->box[i][j];
705 while (dx[i] <= pbc->mhbox_diag[i])
707 for (j = i; j >= 0; j--)
709 dx[j] += pbc->box[i][j];
713 /* dx is the distance in a rectangular box */
715 if (d2min > pbc->max_cutoff2)
717 copy_rvec(dx, dx_start);
719 /* Now try all possible shifts, when the distance is within max_cutoff
720 * it must be the shortest possible distance.
723 while ((d2min > pbc->max_cutoff2) && (i < pbc->ntric_vec))
725 rvec_add(dx_start, pbc->tric_vec[i], trial);
726 d2trial = norm2(trial);
729 copy_rvec(trial, dx);
737 for (i = 0; i < DIM; i++)
741 while (dx[i] > pbc->hbox_diag[i])
743 dx[i] -= pbc->fbox_diag[i];
745 while (dx[i] <= pbc->mhbox_diag[i])
747 dx[i] += pbc->fbox_diag[i];
754 for (i = DIM-1; i >= 0; i--)
758 while (dx[i] > pbc->hbox_diag[i])
760 for (j = i; j >= 0; j--)
762 dx[j] -= pbc->box[i][j];
765 while (dx[i] <= pbc->mhbox_diag[i])
767 for (j = i; j >= 0; j--)
769 dx[j] += pbc->box[i][j];
772 d2min += dx[i]*dx[i];
775 if (d2min > pbc->max_cutoff2)
777 copy_rvec(dx, dx_start);
779 /* Now try all possible shifts, when the distance is within max_cutoff
780 * it must be the shortest possible distance.
783 while ((d2min > pbc->max_cutoff2) && (i < pbc->ntric_vec))
785 rvec_add(dx_start, pbc->tric_vec[i], trial);
787 for (j = 0; j < DIM; j++)
791 d2trial += trial[j]*trial[j];
796 copy_rvec(trial, dx);
803 case epbcdxSCREW_RECT:
804 /* The shift definition requires x first */
806 while (dx[XX] > pbc->hbox_diag[XX])
808 dx[XX] -= pbc->fbox_diag[XX];
811 while (dx[XX] <= pbc->mhbox_diag[XX])
813 dx[XX] += pbc->fbox_diag[YY];
818 /* Rotate around the x-axis in the middle of the box */
819 dx[YY] = pbc->box[YY][YY] - x1[YY] - x2[YY];
820 dx[ZZ] = pbc->box[ZZ][ZZ] - x1[ZZ] - x2[ZZ];
822 /* Normal pbc for y and z */
823 for (i = YY; i <= ZZ; i++)
825 while (dx[i] > pbc->hbox_diag[i])
827 dx[i] -= pbc->fbox_diag[i];
829 while (dx[i] <= pbc->mhbox_diag[i])
831 dx[i] += pbc->fbox_diag[i];
836 case epbcdxUNSUPPORTED:
839 gmx_fatal(FARGS, "Internal error in pbc_dx, set_pbc has not been called");
844 int pbc_dx_aiuc(const t_pbc *pbc, const rvec x1, const rvec x2, rvec dx)
847 rvec dx_start, trial;
849 ivec ishift, ishift_start;
851 rvec_sub(x1, x2, dx);
856 case epbcdxRECTANGULAR:
857 for (i = 0; i < DIM; i++)
859 if (dx[i] > pbc->hbox_diag[i])
861 dx[i] -= pbc->fbox_diag[i];
864 else if (dx[i] <= pbc->mhbox_diag[i])
866 dx[i] += pbc->fbox_diag[i];
871 case epbcdxTRICLINIC:
872 /* For triclinic boxes the performance difference between
873 * if/else and two while loops is negligible.
874 * However, the while version can cause extreme delays
875 * before a simulation crashes due to large forces which
876 * can cause unlimited displacements.
877 * Also allowing multiple shifts would index fshift beyond bounds.
879 for (i = DIM-1; i >= 1; i--)
881 if (dx[i] > pbc->hbox_diag[i])
883 for (j = i; j >= 0; j--)
885 dx[j] -= pbc->box[i][j];
889 else if (dx[i] <= pbc->mhbox_diag[i])
891 for (j = i; j >= 0; j--)
893 dx[j] += pbc->box[i][j];
898 /* Allow 2 shifts in x */
899 if (dx[XX] > pbc->hbox_diag[XX])
901 dx[XX] -= pbc->fbox_diag[XX];
903 if (dx[XX] > pbc->hbox_diag[XX])
905 dx[XX] -= pbc->fbox_diag[XX];
909 else if (dx[XX] <= pbc->mhbox_diag[XX])
911 dx[XX] += pbc->fbox_diag[XX];
913 if (dx[XX] <= pbc->mhbox_diag[XX])
915 dx[XX] += pbc->fbox_diag[XX];
919 /* dx is the distance in a rectangular box */
921 if (d2min > pbc->max_cutoff2)
923 copy_rvec(dx, dx_start);
924 copy_ivec(ishift, ishift_start);
926 /* Now try all possible shifts, when the distance is within max_cutoff
927 * it must be the shortest possible distance.
930 while ((d2min > pbc->max_cutoff2) && (i < pbc->ntric_vec))
932 rvec_add(dx_start, pbc->tric_vec[i], trial);
933 d2trial = norm2(trial);
936 copy_rvec(trial, dx);
937 ivec_add(ishift_start, pbc->tric_shift[i], ishift);
945 for (i = 0; i < DIM; i++)
949 if (dx[i] > pbc->hbox_diag[i])
951 dx[i] -= pbc->fbox_diag[i];
954 else if (dx[i] <= pbc->mhbox_diag[i])
956 dx[i] += pbc->fbox_diag[i];
964 for (i = DIM-1; i >= 1; i--)
968 if (dx[i] > pbc->hbox_diag[i])
970 for (j = i; j >= 0; j--)
972 dx[j] -= pbc->box[i][j];
976 else if (dx[i] <= pbc->mhbox_diag[i])
978 for (j = i; j >= 0; j--)
980 dx[j] += pbc->box[i][j];
984 d2min += dx[i]*dx[i];
989 /* Allow 2 shifts in x */
990 if (dx[XX] > pbc->hbox_diag[XX])
992 dx[XX] -= pbc->fbox_diag[XX];
994 if (dx[XX] > pbc->hbox_diag[XX])
996 dx[XX] -= pbc->fbox_diag[XX];
1000 else if (dx[XX] <= pbc->mhbox_diag[XX])
1002 dx[XX] += pbc->fbox_diag[XX];
1004 if (dx[XX] <= pbc->mhbox_diag[XX])
1006 dx[XX] += pbc->fbox_diag[XX];
1010 d2min += dx[XX]*dx[XX];
1012 if (d2min > pbc->max_cutoff2)
1014 copy_rvec(dx, dx_start);
1015 copy_ivec(ishift, ishift_start);
1016 /* Now try all possible shifts, when the distance is within max_cutoff
1017 * it must be the shortest possible distance.
1020 while ((d2min > pbc->max_cutoff2) && (i < pbc->ntric_vec))
1022 rvec_add(dx_start, pbc->tric_vec[i], trial);
1024 for (j = 0; j < DIM; j++)
1028 d2trial += trial[j]*trial[j];
1031 if (d2trial < d2min)
1033 copy_rvec(trial, dx);
1034 ivec_add(ishift_start, pbc->tric_shift[i], ishift);
1043 if (dx[i] > pbc->hbox_diag[i])
1045 dx[i] -= pbc->fbox_diag[i];
1048 else if (dx[i] <= pbc->mhbox_diag[i])
1050 dx[i] += pbc->fbox_diag[i];
1056 if (dx[i] > pbc->hbox_diag[i])
1058 rvec_dec(dx, pbc->box[i]);
1061 else if (dx[i] <= pbc->mhbox_diag[i])
1063 rvec_inc(dx, pbc->box[i]);
1067 case epbcdxSCREW_RECT:
1068 /* The shift definition requires x first */
1069 if (dx[XX] > pbc->hbox_diag[XX])
1071 dx[XX] -= pbc->fbox_diag[XX];
1074 else if (dx[XX] <= pbc->mhbox_diag[XX])
1076 dx[XX] += pbc->fbox_diag[XX];
1079 if (ishift[XX] == 1 || ishift[XX] == -1)
1081 /* Rotate around the x-axis in the middle of the box */
1082 dx[YY] = pbc->box[YY][YY] - x1[YY] - x2[YY];
1083 dx[ZZ] = pbc->box[ZZ][ZZ] - x1[ZZ] - x2[ZZ];
1085 /* Normal pbc for y and z */
1086 for (i = YY; i <= ZZ; i++)
1088 if (dx[i] > pbc->hbox_diag[i])
1090 dx[i] -= pbc->fbox_diag[i];
1093 else if (dx[i] <= pbc->mhbox_diag[i])
1095 dx[i] += pbc->fbox_diag[i];
1101 case epbcdxUNSUPPORTED:
1104 gmx_fatal(FARGS, "Internal error in pbc_dx_aiuc, set_pbc_dd or set_pbc has not been called");
1108 is = IVEC2IS(ishift);
1111 range_check_mesg(is, 0, SHIFTS, "PBC shift vector index range check.");
1117 void pbc_dx_d(const t_pbc *pbc, const dvec x1, const dvec x2, dvec dx)
1120 dvec dx_start, trial;
1121 double d2min, d2trial;
1124 dvec_sub(x1, x2, dx);
1126 switch (pbc->ePBCDX)
1128 case epbcdxRECTANGULAR:
1130 for (i = 0; i < DIM; i++)
1134 while (dx[i] > pbc->hbox_diag[i])
1136 dx[i] -= pbc->fbox_diag[i];
1138 while (dx[i] <= pbc->mhbox_diag[i])
1140 dx[i] += pbc->fbox_diag[i];
1145 case epbcdxTRICLINIC:
1148 for (i = DIM-1; i >= 0; i--)
1152 while (dx[i] > pbc->hbox_diag[i])
1154 for (j = i; j >= 0; j--)
1156 dx[j] -= pbc->box[i][j];
1159 while (dx[i] <= pbc->mhbox_diag[i])
1161 for (j = i; j >= 0; j--)
1163 dx[j] += pbc->box[i][j];
1166 d2min += dx[i]*dx[i];
1169 if (d2min > pbc->max_cutoff2)
1171 copy_dvec(dx, dx_start);
1172 /* Now try all possible shifts, when the distance is within max_cutoff
1173 * it must be the shortest possible distance.
1176 while ((d2min > pbc->max_cutoff2) && (i < pbc->ntric_vec))
1178 for (j = 0; j < DIM; j++)
1180 trial[j] = dx_start[j] + pbc->tric_vec[i][j];
1183 for (j = 0; j < DIM; j++)
1187 d2trial += trial[j]*trial[j];
1190 if (d2trial < d2min)
1192 copy_dvec(trial, dx);
1199 case epbcdxSCREW_RECT:
1200 /* The shift definition requires x first */
1202 while (dx[XX] > pbc->hbox_diag[XX])
1204 dx[XX] -= pbc->fbox_diag[XX];
1207 while (dx[XX] <= pbc->mhbox_diag[XX])
1209 dx[XX] += pbc->fbox_diag[YY];
1214 /* Rotate around the x-axis in the middle of the box */
1215 dx[YY] = pbc->box[YY][YY] - x1[YY] - x2[YY];
1216 dx[ZZ] = pbc->box[ZZ][ZZ] - x1[ZZ] - x2[ZZ];
1218 /* Normal pbc for y and z */
1219 for (i = YY; i <= ZZ; i++)
1221 while (dx[i] > pbc->hbox_diag[i])
1223 dx[i] -= pbc->fbox_diag[i];
1225 while (dx[i] <= pbc->mhbox_diag[i])
1227 dx[i] += pbc->fbox_diag[i];
1232 case epbcdxUNSUPPORTED:
1235 gmx_fatal(FARGS, "Internal error in pbc_dx, set_pbc has not been called");
1240 gmx_bool image_rect(ivec xi, ivec xj, ivec box_size, real rlong2, int *shift, real *r2)
1248 for (m = 0; (m < DIM); m++)
1280 gmx_bool image_cylindric(ivec xi, ivec xj, ivec box_size, real rlong2,
1281 int *shift, real *r2)
1289 for (m = 0; (m < DIM); m++)
1325 void calc_shifts(matrix box, rvec shift_vec[])
1327 int k, l, m, d, n, test;
1330 for (m = -D_BOX_Z; m <= D_BOX_Z; m++)
1332 for (l = -D_BOX_Y; l <= D_BOX_Y; l++)
1334 for (k = -D_BOX_X; k <= D_BOX_X; k++, n++)
1336 test = XYZ2IS(k, l, m);
1339 gmx_incons("inconsistent shift numbering");
1341 for (d = 0; d < DIM; d++)
1343 shift_vec[n][d] = k*box[XX][d] + l*box[YY][d] + m*box[ZZ][d];
1350 void calc_box_center(int ecenter, matrix box, rvec box_center)
1354 clear_rvec(box_center);
1358 for (m = 0; (m < DIM); m++)
1360 for (d = 0; d < DIM; d++)
1362 box_center[d] += 0.5*box[m][d];
1367 for (d = 0; d < DIM; d++)
1369 box_center[d] = 0.5*box[d][d];
1375 gmx_fatal(FARGS, "Unsupported value %d for ecenter", ecenter);
1379 void calc_triclinic_images(matrix box, rvec img[])
1383 /* Calculate 3 adjacent images in the xy-plane */
1384 copy_rvec(box[0], img[0]);
1385 copy_rvec(box[1], img[1]);
1388 svmul(-1, img[1], img[1]);
1390 rvec_sub(img[1], img[0], img[2]);
1392 /* Get the next 3 in the xy-plane as mirror images */
1393 for (i = 0; i < 3; i++)
1395 svmul(-1, img[i], img[3+i]);
1398 /* Calculate the first 4 out of xy-plane images */
1399 copy_rvec(box[2], img[6]);
1402 svmul(-1, img[6], img[6]);
1404 for (i = 0; i < 3; i++)
1406 rvec_add(img[6], img[i+1], img[7+i]);
1409 /* Mirror the last 4 from the previous in opposite rotation */
1410 for (i = 0; i < 4; i++)
1412 svmul(-1, img[6 + (2+i) % 4], img[10+i]);
1416 void calc_compact_unitcell_vertices(int ecenter, matrix box, rvec vert[])
1418 rvec img[NTRICIMG], box_center;
1419 int n, i, j, tmp[4], d;
1421 calc_triclinic_images(box, img);
1424 for (i = 2; i <= 5; i += 3)
1437 for (j = 0; j < 4; j++)
1439 for (d = 0; d < DIM; d++)
1441 vert[n][d] = img[i][d]+img[tmp[j]][d]+img[tmp[(j+1)%4]][d];
1446 for (i = 7; i <= 13; i += 6)
1459 for (j = 0; j < 4; j++)
1461 for (d = 0; d < DIM; d++)
1463 vert[n][d] = img[i][d]+img[tmp[j]][d]+img[tmp[(j+1)%4]][d];
1468 for (i = 9; i <= 11; i += 2)
1488 for (j = 0; j < 4; j++)
1490 for (d = 0; d < DIM; d++)
1492 vert[n][d] = img[i][d]+img[tmp[j]][d]+img[tmp[(j+1)%4]][d];
1498 calc_box_center(ecenter, box, box_center);
1499 for (i = 0; i < NCUCVERT; i++)
1501 for (d = 0; d < DIM; d++)
1503 vert[i][d] = vert[i][d]*0.25+box_center[d];
1508 int *compact_unitcell_edges()
1510 /* this is an index in vert[] (see calc_box_vertices) */
1511 /*static int edge[NCUCEDGE*2];*/
1513 static const int hexcon[24] = {
1514 0, 9, 1, 19, 2, 15, 3, 21,
1515 4, 17, 5, 11, 6, 23, 7, 13,
1516 8, 20, 10, 18, 12, 16, 14, 22
1519 gmx_bool bFirst = TRUE;
1521 snew(edge, NCUCEDGE*2);
1526 for (i = 0; i < 6; i++)
1528 for (j = 0; j < 4; j++)
1530 edge[e++] = 4*i + j;
1531 edge[e++] = 4*i + (j+1) % 4;
1534 for (i = 0; i < 12*2; i++)
1536 edge[e++] = hexcon[i];
1545 void put_atoms_in_box_omp(int ePBC, matrix box, int natoms, rvec x[])
1548 nth = gmx_omp_nthreads_get(emntDefault);
1550 #pragma omp parallel for num_threads(nth) schedule(static)
1551 for (t = 0; t < nth; t++)
1555 offset = (natoms*t )/nth;
1556 len = (natoms*(t + 1))/nth - offset;
1557 put_atoms_in_box(ePBC, box, len, x + offset);
1561 void put_atoms_in_box(int ePBC, matrix box, int natoms, rvec x[])
1563 int npbcdim, i, m, d;
1565 if (ePBC == epbcSCREW)
1567 gmx_fatal(FARGS, "Sorry, %s pbc is not yet supported", epbc_names[ePBC]);
1581 for (i = 0; (i < natoms); i++)
1583 for (m = npbcdim-1; m >= 0; m--)
1587 for (d = 0; d <= m; d++)
1589 x[i][d] += box[m][d];
1592 while (x[i][m] >= box[m][m])
1594 for (d = 0; d <= m; d++)
1596 x[i][d] -= box[m][d];
1604 for (i = 0; i < natoms; i++)
1606 for (d = 0; d < npbcdim; d++)
1610 x[i][d] += box[d][d];
1612 while (x[i][d] >= box[d][d])
1614 x[i][d] -= box[d][d];
1621 void put_atoms_in_triclinic_unitcell(int ecenter, matrix box,
1622 int natoms, rvec x[])
1624 rvec box_center, shift_center;
1625 real shm01, shm02, shm12, shift;
1628 calc_box_center(ecenter, box, box_center);
1630 /* The product of matrix shm with a coordinate gives the shift vector
1631 which is required determine the periodic cell position */
1632 shm01 = box[1][0]/box[1][1];
1633 shm02 = (box[1][1]*box[2][0] - box[2][1]*box[1][0])/(box[1][1]*box[2][2]);
1634 shm12 = box[2][1]/box[2][2];
1636 clear_rvec(shift_center);
1637 for (d = 0; d < DIM; d++)
1639 rvec_inc(shift_center, box[d]);
1641 svmul(0.5, shift_center, shift_center);
1642 rvec_sub(box_center, shift_center, shift_center);
1644 shift_center[0] = shm01*shift_center[1] + shm02*shift_center[2];
1645 shift_center[1] = shm12*shift_center[2];
1646 shift_center[2] = 0;
1648 for (i = 0; (i < natoms); i++)
1650 for (m = DIM-1; m >= 0; m--)
1652 shift = shift_center[m];
1655 shift += shm01*x[i][1] + shm02*x[i][2];
1659 shift += shm12*x[i][2];
1661 while (x[i][m]-shift < 0)
1663 for (d = 0; d <= m; d++)
1665 x[i][d] += box[m][d];
1668 while (x[i][m]-shift >= box[m][m])
1670 for (d = 0; d <= m; d++)
1672 x[i][d] -= box[m][d];
1680 put_atoms_in_compact_unitcell(int ePBC, int ecenter, matrix box,
1681 int natoms, rvec x[])
1684 rvec box_center, dx;
1687 set_pbc(&pbc, ePBC, box);
1688 calc_box_center(ecenter, box, box_center);
1689 for (i = 0; i < natoms; i++)
1691 pbc_dx(&pbc, x[i], box_center, dx);
1692 rvec_add(box_center, dx, x[i]);
1695 return pbc.bLimitDistance ?
1696 "WARNING: Could not put all atoms in the compact unitcell\n"