From: Joe Jordan Date: Tue, 23 Mar 2021 09:40:41 +0000 (+0000) Subject: Use ArrayRef in calc_mu X-Git-Url: http://biod.pnpi.spb.ru/gitweb/?a=commitdiff_plain;h=6af73f1a4dd9b17cb0d90bdb959c3489de9f5f3a;p=alexxy%2Fgromacs.git Use ArrayRef in calc_mu --- diff --git a/src/gromacs/mdlib/calcmu.cpp b/src/gromacs/mdlib/calcmu.cpp index 982faf33e0..07b27aad00 100644 --- a/src/gromacs/mdlib/calcmu.cpp +++ b/src/gromacs/mdlib/calcmu.cpp @@ -52,9 +52,9 @@ void calc_mu(int start, int homenr, gmx::ArrayRef x, - const real q[], - const real qB[], - int nChargePerturbed, + gmx::ArrayRef q, + gmx::ArrayRef qB, + bool havePerturbedCharges, dvec mu, dvec mu_B) { @@ -82,7 +82,7 @@ void calc_mu(int start, mu[m] *= gmx::c_enm2Debye; } - if (nChargePerturbed) + if (havePerturbedCharges) { mu_x = mu_y = mu_z = 0.0; #pragma omp parallel for reduction(+: mu_x, mu_y, mu_z) schedule(static) \ diff --git a/src/gromacs/mdlib/calcmu.h b/src/gromacs/mdlib/calcmu.h index 4081afc4ef..23a4fd2907 100644 --- a/src/gromacs/mdlib/calcmu.h +++ b/src/gromacs/mdlib/calcmu.h @@ -4,7 +4,7 @@ * Copyright (c) 1991-2000, University of Groningen, The Netherlands. * Copyright (c) 2001-2004, The GROMACS development team. * Copyright (c) 2010,2014,2015,2017,2018 by the GROMACS development team. - * Copyright (c) 2019,2020, by the GROMACS development team, led by + * Copyright (c) 2019,2020,2021, by the GROMACS development team, led by * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl, * and including many others, as listed in the AUTHORS file in the * top-level source directory and at http://www.gromacs.org. @@ -52,9 +52,9 @@ class ArrayRef; void calc_mu(int start, int homenr, gmx::ArrayRef x, - const real q[], - const real qB[], - int nChargePerturbed, + gmx::ArrayRef q, + gmx::ArrayRef qB, + bool havePerturbedCharges, dvec mu, dvec mu_B); diff --git a/src/gromacs/mdlib/sim_util.cpp b/src/gromacs/mdlib/sim_util.cpp index ec10bb27a0..30985796d9 100644 --- a/src/gromacs/mdlib/sim_util.cpp +++ b/src/gromacs/mdlib/sim_util.cpp @@ -1663,9 +1663,9 @@ void do_force(FILE* fplog, calc_mu(start, mdatoms->homenr, xRef, - mdatoms->chargeA, - mdatoms->chargeB, - mdatoms->nChargePerturbed, + gmx::arrayRefFromArray(mdatoms->chargeA, mdatoms->nr), + gmx::arrayRefFromArray(mdatoms->chargeB, mdatoms->nr), + mdatoms->nChargePerturbed != 0, dipoleData.muStaging[0], dipoleData.muStaging[1]);