2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2006,2007,2008,2009,2010,2012,2013,2014,2015, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
38 * \brief This file implements functions for domdec to use
39 * while managing inter-atomic constraints.
41 * \author Berk Hess <hess@kth.se>
42 * \ingroup module_domdec
47 #include "domdec_vsite.h"
53 #include "gromacs/domdec/domdec.h"
54 #include "gromacs/domdec/domdec_struct.h"
55 #include "gromacs/domdec/ga2la.h"
56 #include "gromacs/gmxlib/gmx_omp_nthreads.h"
57 #include "gromacs/legacyheaders/types/commrec.h"
58 #include "gromacs/math/vec.h"
59 #include "gromacs/pbcutil/ishift.h"
60 #include "gromacs/topology/mtop_util.h"
61 #include "gromacs/utility/fatalerror.h"
62 #include "gromacs/utility/gmxassert.h"
63 #include "gromacs/utility/smalloc.h"
65 #include "domdec_specatomcomm.h"
68 void dd_move_f_vsites(gmx_domdec_t
*dd
, rvec
*f
, rvec
*fshift
)
72 dd_move_f_specat(dd
, dd
->vsite_comm
, f
, fshift
);
76 void dd_clear_f_vsites(gmx_domdec_t
*dd
, rvec
*f
)
82 for (i
= dd
->vsite_comm
->at_start
; i
< dd
->vsite_comm
->at_end
; i
++)
89 void dd_move_x_vsites(gmx_domdec_t
*dd
, matrix box
, rvec
*x
)
93 dd_move_x_specat(dd
, dd
->vsite_comm
, box
, x
, NULL
, FALSE
);
97 int dd_make_local_vsites(gmx_domdec_t
*dd
, int at_start
, t_ilist
*lil
)
99 gmx_domdec_specat_comm_t
*spac
;
101 gmx_hash_t
*ga2la_specat
;
102 int ftype
, nral
, i
, j
, a
;
107 spac
= dd
->vsite_comm
;
108 ireq
= &spac
->ireq
[0];
109 ga2la_specat
= dd
->ga2la_vsite
;
112 /* Loop over all the home vsites */
113 for (ftype
= 0; ftype
< F_NRE
; ftype
++)
115 if (interaction_function
[ftype
].flags
& IF_VSITE
)
119 for (i
= 0; i
< lilf
->nr
; i
+= 1+nral
)
121 iatoms
= lilf
->iatoms
+ i
;
122 /* Check if we have the other atoms */
123 for (j
= 1; j
< 1+nral
; j
++)
127 /* This is not a home atom,
128 * we need to ask our neighbors.
131 /* Check to not ask for the same atom more than once */
132 if (gmx_hash_get_minone(dd
->ga2la_vsite
, a
) == -1)
134 /* Add this non-home atom to the list */
135 if (ireq
->n
+1 > ireq
->nalloc
)
137 ireq
->nalloc
= over_alloc_large(ireq
->n
+1);
138 srenew(ireq
->ind
, ireq
->nalloc
);
140 ireq
->ind
[ireq
->n
++] = a
;
141 /* Temporarily mark with -2,
142 * we get the index later.
144 gmx_hash_set(ga2la_specat
, a
, -2);
152 at_end
= setup_specat_communication(dd
, ireq
, dd
->vsite_comm
, ga2la_specat
,
153 at_start
, 1, "vsite", "");
155 /* Fill in the missing indices */
156 for (ftype
= 0; ftype
< F_NRE
; ftype
++)
158 if (interaction_function
[ftype
].flags
& IF_VSITE
)
162 for (i
= 0; i
< lilf
->nr
; i
+= 1+nral
)
164 iatoms
= lilf
->iatoms
+ i
;
165 for (j
= 1; j
< 1+nral
; j
++)
169 iatoms
[j
] = gmx_hash_get_minone(ga2la_specat
, -iatoms
[j
]-1);
179 void init_domdec_vsites(gmx_domdec_t
*dd
, int n_intercg_vsite
)
183 fprintf(debug
, "Begin init_domdec_vsites\n");
186 /* Use a hash table for the global to local index.
187 * The number of keys is a rough estimate, it will be optimized later.
189 dd
->ga2la_vsite
= gmx_hash_init(std::min(n_intercg_vsite
/20,
190 n_intercg_vsite
/(2*dd
->nnodes
)));
192 dd
->vsite_comm
= specat_comm_init(1);