Cleanups in legacyheader/types/commrec.h
[gromacs.git] / src / gromacs / domdec / domdec_vsite.cpp
blob958be5c56fe4302d1c517b157846610162245a6d
1 /*
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.
36 /*! \internal \file
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
45 #include "gmxpre.h"
47 #include "domdec_vsite.h"
49 #include <assert.h>
51 #include <algorithm>
53 #include "gromacs/domdec/domdec.h"
54 #include "gromacs/domdec/ga2la.h"
55 #include "gromacs/gmxlib/gmx_omp_nthreads.h"
56 #include "gromacs/legacyheaders/types/commrec.h"
57 #include "gromacs/math/vec.h"
58 #include "gromacs/pbcutil/ishift.h"
59 #include "gromacs/topology/mtop_util.h"
60 #include "gromacs/utility/fatalerror.h"
61 #include "gromacs/utility/gmxassert.h"
62 #include "gromacs/utility/smalloc.h"
64 #include "domdec_specatomcomm.h"
65 #include "hash.h"
67 void dd_move_f_vsites(gmx_domdec_t *dd, rvec *f, rvec *fshift)
69 if (dd->vsite_comm)
71 dd_move_f_specat(dd, dd->vsite_comm, f, fshift);
75 void dd_clear_f_vsites(gmx_domdec_t *dd, rvec *f)
77 int i;
79 if (dd->vsite_comm)
81 for (i = dd->vsite_comm->at_start; i < dd->vsite_comm->at_end; i++)
83 clear_rvec(f[i]);
88 void dd_move_x_vsites(gmx_domdec_t *dd, matrix box, rvec *x)
90 if (dd->vsite_comm)
92 dd_move_x_specat(dd, dd->vsite_comm, box, x, NULL, FALSE);
96 int dd_make_local_vsites(gmx_domdec_t *dd, int at_start, t_ilist *lil)
98 gmx_domdec_specat_comm_t *spac;
99 ind_req_t *ireq;
100 gmx_hash_t *ga2la_specat;
101 int ftype, nral, i, j, a;
102 t_ilist *lilf;
103 t_iatom *iatoms;
104 int at_end;
106 spac = dd->vsite_comm;
107 ireq = &spac->ireq[0];
108 ga2la_specat = dd->ga2la_vsite;
110 ireq->n = 0;
111 /* Loop over all the home vsites */
112 for (ftype = 0; ftype < F_NRE; ftype++)
114 if (interaction_function[ftype].flags & IF_VSITE)
116 nral = NRAL(ftype);
117 lilf = &lil[ftype];
118 for (i = 0; i < lilf->nr; i += 1+nral)
120 iatoms = lilf->iatoms + i;
121 /* Check if we have the other atoms */
122 for (j = 1; j < 1+nral; j++)
124 if (iatoms[j] < 0)
126 /* This is not a home atom,
127 * we need to ask our neighbors.
129 a = -iatoms[j] - 1;
130 /* Check to not ask for the same atom more than once */
131 if (gmx_hash_get_minone(dd->ga2la_vsite, a) == -1)
133 /* Add this non-home atom to the list */
134 if (ireq->n+1 > ireq->nalloc)
136 ireq->nalloc = over_alloc_large(ireq->n+1);
137 srenew(ireq->ind, ireq->nalloc);
139 ireq->ind[ireq->n++] = a;
140 /* Temporarily mark with -2,
141 * we get the index later.
143 gmx_hash_set(ga2la_specat, a, -2);
151 at_end = setup_specat_communication(dd, ireq, dd->vsite_comm, ga2la_specat,
152 at_start, 1, "vsite", "");
154 /* Fill in the missing indices */
155 for (ftype = 0; ftype < F_NRE; ftype++)
157 if (interaction_function[ftype].flags & IF_VSITE)
159 nral = NRAL(ftype);
160 lilf = &lil[ftype];
161 for (i = 0; i < lilf->nr; i += 1+nral)
163 iatoms = lilf->iatoms + i;
164 for (j = 1; j < 1+nral; j++)
166 if (iatoms[j] < 0)
168 iatoms[j] = gmx_hash_get_minone(ga2la_specat, -iatoms[j]-1);
175 return at_end;
178 void init_domdec_vsites(gmx_domdec_t *dd, int n_intercg_vsite)
180 if (debug)
182 fprintf(debug, "Begin init_domdec_vsites\n");
185 /* Use a hash table for the global to local index.
186 * The number of keys is a rough estimate, it will be optimized later.
188 dd->ga2la_vsite = gmx_hash_init(std::min(n_intercg_vsite/20,
189 n_intercg_vsite/(2*dd->nnodes)));
191 dd->vsite_comm = specat_comm_init(1);