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-2008, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017 by the GROMACS development team.
7 * Copyright (c) 2018,2019,2020,2021, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
40 #include "pull_rotation.h"
51 #include "gromacs/commandline/filenm.h"
52 #include "gromacs/domdec/dlbtiming.h"
53 #include "gromacs/domdec/domdec_struct.h"
54 #include "gromacs/domdec/ga2la.h"
55 #include "gromacs/domdec/localatomset.h"
56 #include "gromacs/domdec/localatomsetmanager.h"
57 #include "gromacs/fileio/gmxfio.h"
58 #include "gromacs/fileio/xvgr.h"
59 #include "gromacs/gmxlib/network.h"
60 #include "gromacs/linearalgebra/nrjac.h"
61 #include "gromacs/math/functions.h"
62 #include "gromacs/math/utilities.h"
63 #include "gromacs/math/vec.h"
64 #include "gromacs/mdlib/groupcoord.h"
65 #include "gromacs/mdlib/stat.h"
66 #include "gromacs/mdrunutility/handlerestart.h"
67 #include "gromacs/mdtypes/commrec.h"
68 #include "gromacs/mdtypes/inputrec.h"
69 #include "gromacs/mdtypes/md_enums.h"
70 #include "gromacs/mdtypes/mdrunoptions.h"
71 #include "gromacs/mdtypes/state.h"
72 #include "gromacs/pbcutil/pbc.h"
73 #include "gromacs/timing/cyclecounter.h"
74 #include "gromacs/timing/wallcycle.h"
75 #include "gromacs/topology/mtop_lookup.h"
76 #include "gromacs/topology/mtop_util.h"
77 #include "gromacs/utility/basedefinitions.h"
78 #include "gromacs/utility/fatalerror.h"
79 #include "gromacs/utility/pleasecite.h"
80 #include "gromacs/utility/smalloc.h"
82 static char const* RotStr = { "Enforced rotation:" };
84 /* Set the minimum weight for the determination of the slab centers */
85 #define WEIGHT_MIN (10 * GMX_FLOAT_MIN)
87 //! Helper structure for sorting positions along rotation vector
88 struct sort_along_vec_t
90 //! Projection of xc on the rotation vector
98 //! Reference position
103 //! Enforced rotation / flexible: determine the angle of each slab
106 //! Number of atoms belonging to this slab
108 /*! \brief The positions belonging to this slab.
110 * In general, this should be all positions of the whole
111 * rotation group, but we leave those away that have a small
114 //! Same for reference
116 //! The weight for each atom
121 //! Helper structure for potential fitting
124 /*! \brief Set of angles for which the potential is calculated.
126 * The optimum fit is determined as the angle for with the
127 * potential is minimal. */
129 //! Potential for the different angles
131 //! Rotation matrix corresponding to the angles
136 //! Enforced rotation data for a single rotation group
139 //! Input parameters for this group
140 const t_rotgrp* rotg = nullptr;
141 //! Index of this group within the set of groups
143 //! Rotation angle in degrees
147 //! The atoms subject to enforced rotation
148 std::unique_ptr<gmx::LocalAtomSet> atomSet;
150 //! The normalized rotation vector
152 //! Rotation potential for this rotation group
154 //! Array to store the forces on the local atoms resulting from enforced rotation potential
157 /* Collective coordinates for the whole rotation group */
158 //! Length of each x_rotref vector after x_rotref has been put into origin
160 //! Center of the rotation group positions, may be mass weighted
162 //! Center of the rotation group reference positions
164 //! Current (collective) positions
166 //! Current (collective) shifts
168 //! Extra shifts since last DD step
170 //! Old (collective) positions
172 //! Normalized form of the current positions
174 //! Reference positions (sorted in the same order as xc when sorted)
176 //! Where is a position found after sorting?
178 //! Collective masses
180 //! Collective masses sorted
182 //! one over the total mass of the rotation group
185 //! Torque in the direction of rotation vector
187 //! Actual angle of the whole rotation group
189 /* Fixed rotation only */
190 //! Weights for angle determination
192 //! Local reference coords, correctly rotated
194 //! Local current coords, correct PBC image
196 //! Masses of the current local atoms
199 /* Flexible rotation only */
200 //! For this many slabs memory is allocated
202 //! Lowermost slab for that the calculation needs to be performed at a given time step
204 //! Uppermost slab ...
206 //! First slab for which ref. center is stored
210 //! Slab buffer region around reference slabs
212 //! First relevant atom for a slab
214 //! Last relevant atom for a slab
216 //! Gaussian-weighted slab center
218 //! Gaussian-weighted slab center for the reference positions
219 rvec* slab_center_ref;
220 //! Sum of gaussian weights in a slab
222 //! Torque T = r x f for each slab. torque_v = m.v = angular momentum in the direction of v
224 //! min_gaussian from t_rotgrp is the minimum value the gaussian must have so that the force is actually evaluated. max_beta is just another way to put it
226 //! Precalculated gaussians for a single atom
228 //! Tells to which slab each precalculated gaussian belongs
230 //! Inner sum of the flexible2 potential per slab; this is precalculated for optimization reasons
231 rvec* slab_innersumvec;
232 //! Holds atom positions and gaussian weights of atoms belonging to a slab
233 gmx_slabdata* slab_data;
235 /* For potential fits with varying angle: */
236 //! Used for fit type 'potential'
237 gmx_potfit* PotAngleFit;
241 //! Enforced rotation data for all groups
244 //! Input parameters.
245 const t_rot* rot = nullptr;
246 //! Output period for main rotation outfile
248 //! Output period for per-slab data
250 //! Output file for rotation data
251 FILE* out_rot = nullptr;
252 //! Output file for torque data
253 FILE* out_torque = nullptr;
254 //! Output file for slab angles for flexible type
255 FILE* out_angles = nullptr;
256 //! Output file for slab centers
257 FILE* out_slabs = nullptr;
258 //! Allocation size of buf
260 //! Coordinate buffer variable for sorting
261 rvec* xbuf = nullptr;
262 //! Masses buffer variable for sorting
263 real* mbuf = nullptr;
264 //! Buffer variable needed for position sorting
265 sort_along_vec_t* data = nullptr;
267 real* mpi_inbuf = nullptr;
269 real* mpi_outbuf = nullptr;
270 //! Allocation size of in & outbuf
272 //! If true, append output files
273 gmx_bool restartWithAppending = false;
274 //! Used to skip first output when appending to avoid duplicate entries in rotation outfiles
275 gmx_bool bOut = false;
276 //! Stores working data per group
277 std::vector<gmx_enfrotgrp> enfrotgrp;
281 gmx_enfrot::~gmx_enfrot()
285 gmx_fio_fclose(out_rot);
289 gmx_fio_fclose(out_slabs);
293 gmx_fio_fclose(out_angles);
297 gmx_fio_fclose(out_torque);
304 extern template LocalAtomSet LocalAtomSetManager::add<void, void>(ArrayRef<const int> globalAtomIndex);
306 class EnforcedRotation::Impl
309 gmx_enfrot enforcedRotation_;
312 EnforcedRotation::EnforcedRotation() : impl_(new Impl) {}
314 EnforcedRotation::~EnforcedRotation() = default;
316 gmx_enfrot* EnforcedRotation::getLegacyEnfrot()
318 return &impl_->enforcedRotation_;
323 /* Activate output of forces for correctness checks */
324 /* #define PRINT_FORCES */
326 # define PRINT_FORCE_J \
328 "f%d = %15.8f %15.8f %15.8f\n", \
329 erg->xc_ref_ind[j], \
330 erg->f_rot_loc[j][XX], \
331 erg->f_rot_loc[j][YY], \
332 erg->f_rot_loc[j][ZZ]);
333 # define PRINT_POT_TAU \
337 "potential = %15.8f\n" \
338 "torque = %15.8f\n", \
343 # define PRINT_FORCE_J
344 # define PRINT_POT_TAU
347 /* Shortcuts for often used queries */
349 (((rg)->eType == EnforcedRotationGroupType::Flex) || ((rg)->eType == EnforcedRotationGroupType::Flext) \
350 || ((rg)->eType == EnforcedRotationGroupType::Flex2) \
351 || ((rg)->eType == EnforcedRotationGroupType::Flex2t))
353 (((rg)->eType == EnforcedRotationGroupType::Flex) || ((rg)->eType == EnforcedRotationGroupType::Flext) \
354 || ((rg)->eType == EnforcedRotationGroupType::Flex2) \
355 || ((rg)->eType == EnforcedRotationGroupType::Flex2t) \
356 || ((rg)->eType == EnforcedRotationGroupType::Rmpf) \
357 || ((rg)->eType == EnforcedRotationGroupType::Rm2pf))
360 /* Does any of the rotation groups use slab decomposition? */
361 static gmx_bool HaveFlexibleGroups(const t_rot* rot)
363 for (int g = 0; g < rot->ngrp; g++)
365 if (ISFLEX(&rot->grp[g]))
375 /* Is for any group the fit angle determined by finding the minimum of the
376 * rotation potential? */
377 static gmx_bool HavePotFitGroups(const t_rot* rot)
379 for (int g = 0; g < rot->ngrp; g++)
381 if (RotationGroupFitting::Pot == rot->grp[g].eFittype)
391 static double** allocate_square_matrix(int dim)
394 double** mat = nullptr;
398 for (i = 0; i < dim; i++)
407 static void free_square_matrix(double** mat, int dim)
412 for (i = 0; i < dim; i++)
420 /* Return the angle for which the potential is minimal */
421 static real get_fitangle(const gmx_enfrotgrp* erg)
424 real fitangle = -999.9;
425 real pot_min = GMX_FLOAT_MAX;
429 fit = erg->PotAngleFit;
431 for (i = 0; i < erg->rotg->PotAngle_nstep; i++)
433 if (fit->V[i] < pot_min)
436 fitangle = fit->degangle[i];
444 /* Reduce potential angle fit data for this group at this time step? */
445 static inline gmx_bool bPotAngle(const gmx_enfrot* er, const t_rotgrp* rotg, int64_t step)
447 return ((RotationGroupFitting::Pot == rotg->eFittype)
448 && (do_per_step(step, er->nstsout) || do_per_step(step, er->nstrout)));
451 /* Reduce slab torqe data for this group at this time step? */
452 static inline gmx_bool bSlabTau(const gmx_enfrot* er, const t_rotgrp* rotg, int64_t step)
454 return ((ISFLEX(rotg)) && do_per_step(step, er->nstsout));
457 /* Output rotation energy, torques, etc. for each rotation group */
458 static void reduce_output(const t_commrec* cr, gmx_enfrot* er, real t, int64_t step)
460 int i, islab, nslabs = 0;
461 int count; /* MPI element counter */
465 /* Fill the MPI buffer with stuff to reduce. If items are added for reduction
466 * here, the MPI buffer size has to be enlarged also in calc_mpi_bufsize() */
470 for (auto& ergRef : er->enfrotgrp)
472 gmx_enfrotgrp* erg = &ergRef;
473 const t_rotgrp* rotg = erg->rotg;
474 nslabs = erg->slab_last - erg->slab_first + 1;
475 er->mpi_inbuf[count++] = erg->V;
476 er->mpi_inbuf[count++] = erg->torque_v;
477 er->mpi_inbuf[count++] = erg->angle_v;
478 er->mpi_inbuf[count++] =
479 erg->weight_v; /* weights are not needed for flex types, but this is just a single value */
481 if (bPotAngle(er, rotg, step))
483 for (i = 0; i < rotg->PotAngle_nstep; i++)
485 er->mpi_inbuf[count++] = erg->PotAngleFit->V[i];
488 if (bSlabTau(er, rotg, step))
490 for (i = 0; i < nslabs; i++)
492 er->mpi_inbuf[count++] = erg->slab_torque_v[i];
496 if (count > er->mpi_bufsize)
498 gmx_fatal(FARGS, "%s MPI buffer overflow, please report this error.", RotStr);
502 MPI_Reduce(er->mpi_inbuf, er->mpi_outbuf, count, GMX_MPI_REAL, MPI_SUM, MASTERRANK(cr), cr->mpi_comm_mygroup);
505 /* Copy back the reduced data from the buffer on the master */
509 for (auto& ergRef : er->enfrotgrp)
511 gmx_enfrotgrp* erg = &ergRef;
512 const t_rotgrp* rotg = erg->rotg;
513 nslabs = erg->slab_last - erg->slab_first + 1;
514 erg->V = er->mpi_outbuf[count++];
515 erg->torque_v = er->mpi_outbuf[count++];
516 erg->angle_v = er->mpi_outbuf[count++];
517 erg->weight_v = er->mpi_outbuf[count++];
519 if (bPotAngle(er, rotg, step))
521 for (int i = 0; i < rotg->PotAngle_nstep; i++)
523 erg->PotAngleFit->V[i] = er->mpi_outbuf[count++];
526 if (bSlabTau(er, rotg, step))
528 for (int i = 0; i < nslabs; i++)
530 erg->slab_torque_v[i] = er->mpi_outbuf[count++];
540 /* Angle and torque for each rotation group */
541 for (auto& ergRef : er->enfrotgrp)
543 gmx_enfrotgrp* erg = &ergRef;
544 const t_rotgrp* rotg = erg->rotg;
545 bFlex = ISFLEX(rotg);
547 /* Output to main rotation output file: */
548 if (do_per_step(step, er->nstrout))
550 if (RotationGroupFitting::Pot == rotg->eFittype)
552 fitangle = get_fitangle(erg);
558 fitangle = erg->angle_v; /* RMSD fit angle */
562 fitangle = (erg->angle_v / erg->weight_v) * 180.0 * M_1_PI;
565 fprintf(er->out_rot, "%12.4f", fitangle);
566 fprintf(er->out_rot, "%12.3e", erg->torque_v);
567 fprintf(er->out_rot, "%12.3e", erg->V);
570 if (do_per_step(step, er->nstsout))
572 /* Output to torque log file: */
575 fprintf(er->out_torque, "%12.3e%6d", t, erg->groupIndex);
576 for (int i = erg->slab_first; i <= erg->slab_last; i++)
578 islab = i - erg->slab_first; /* slab index */
579 /* Only output if enough weight is in slab */
580 if (erg->slab_weights[islab] > rotg->min_gaussian)
582 fprintf(er->out_torque, "%6d%12.3e", i, erg->slab_torque_v[islab]);
585 fprintf(er->out_torque, "\n");
588 /* Output to angles log file: */
589 if (RotationGroupFitting::Pot == rotg->eFittype)
591 fprintf(er->out_angles, "%12.3e%6d%12.4f", t, erg->groupIndex, erg->degangle);
592 /* Output energies at a set of angles around the reference angle */
593 for (int i = 0; i < rotg->PotAngle_nstep; i++)
595 fprintf(er->out_angles, "%12.3e", erg->PotAngleFit->V[i]);
597 fprintf(er->out_angles, "\n");
601 if (do_per_step(step, er->nstrout))
603 fprintf(er->out_rot, "\n");
609 /* Add the forces from enforced rotation potential to the local forces.
610 * Should be called after the SR forces have been evaluated */
611 real add_rot_forces(gmx_enfrot* er, rvec f[], const t_commrec* cr, int64_t step, real t)
613 real Vrot = 0.0; /* If more than one rotation group is present, Vrot
614 assembles the local parts from all groups */
616 /* Loop over enforced rotation groups (usually 1, though)
617 * Apply the forces from rotation potentials */
618 for (auto& ergRef : er->enfrotgrp)
620 gmx_enfrotgrp* erg = &ergRef;
621 Vrot += erg->V; /* add the local parts from the nodes */
622 const auto& localRotationGroupIndex = erg->atomSet->localIndex();
623 for (gmx::index l = 0; l < localRotationGroupIndex.ssize(); l++)
625 /* Get the right index of the local force */
626 int ii = localRotationGroupIndex[l];
628 rvec_inc(f[ii], erg->f_rot_loc[l]);
632 /* Reduce energy,torque, angles etc. to get the sum values (per rotation group)
633 * on the master and output these values to file. */
634 if ((do_per_step(step, er->nstrout) || do_per_step(step, er->nstsout)) && er->bOut)
636 reduce_output(cr, er, t, step);
639 /* When appending, er->bOut is FALSE the first time to avoid duplicate entries */
648 /* The Gaussian norm is chosen such that the sum of the gaussian functions
649 * over the slabs is approximately 1.0 everywhere */
650 #define GAUSS_NORM 0.569917543430618
653 /* Calculate the maximum beta that leads to a gaussian larger min_gaussian,
654 * also does some checks
656 static double calc_beta_max(real min_gaussian, real slab_dist)
662 /* Actually the next two checks are already made in grompp */
665 gmx_fatal(FARGS, "Slab distance of flexible rotation groups must be >=0 !");
667 if (min_gaussian <= 0)
669 gmx_fatal(FARGS, "Cutoff value for Gaussian must be > 0. (You requested %f)", min_gaussian);
672 /* Define the sigma value */
673 sigma = 0.7 * slab_dist;
675 /* Calculate the argument for the logarithm and check that the log() result is negative or 0 */
676 arg = min_gaussian / GAUSS_NORM;
679 gmx_fatal(FARGS, "min_gaussian of flexible rotation groups must be <%g", GAUSS_NORM);
682 return std::sqrt(-2.0 * sigma * sigma * log(min_gaussian / GAUSS_NORM));
686 static inline real calc_beta(rvec curr_x, const gmx_enfrotgrp* erg, int n)
688 return iprod(curr_x, erg->vec) - erg->rotg->slab_dist * n;
692 static inline real gaussian_weight(rvec curr_x, const gmx_enfrotgrp* erg, int n)
694 const real norm = GAUSS_NORM;
698 /* Define the sigma value */
699 sigma = 0.7 * erg->rotg->slab_dist;
700 /* Calculate the Gaussian value of slab n for position curr_x */
701 return norm * exp(-0.5 * gmx::square(calc_beta(curr_x, erg, n) / sigma));
705 /* Returns the weight in a single slab, also calculates the Gaussian- and mass-
706 * weighted sum of positions for that slab */
707 static real get_slab_weight(int j, const gmx_enfrotgrp* erg, rvec xc[], const real mc[], rvec* x_weighted_sum)
709 rvec curr_x; /* The position of an atom */
710 rvec curr_x_weighted; /* The gaussian-weighted position */
711 real gaussian; /* A single gaussian weight */
712 real wgauss; /* gaussian times current mass */
713 real slabweight = 0.0; /* The sum of weights in the slab */
715 clear_rvec(*x_weighted_sum);
717 /* Loop over all atoms in the rotation group */
718 for (int i = 0; i < erg->rotg->nat; i++)
720 copy_rvec(xc[i], curr_x);
721 gaussian = gaussian_weight(curr_x, erg, j);
722 wgauss = gaussian * mc[i];
723 svmul(wgauss, curr_x, curr_x_weighted);
724 rvec_add(*x_weighted_sum, curr_x_weighted, *x_weighted_sum);
725 slabweight += wgauss;
726 } /* END of loop over rotation group atoms */
732 static void get_slab_centers(gmx_enfrotgrp* erg, /* Enforced rotation group working data */
733 rvec* xc, /* The rotation group positions; will
734 typically be enfrotgrp->xc, but at first call
735 it is enfrotgrp->xc_ref */
736 real* mc, /* The masses of the rotation group atoms */
737 real time, /* Used for output only */
738 FILE* out_slabs, /* For outputting center per slab information */
739 gmx_bool bOutStep, /* Is this an output step? */
740 gmx_bool bReference) /* If this routine is called from
741 init_rot_group we need to store
742 the reference slab centers */
744 /* Loop over slabs */
745 for (int j = erg->slab_first; j <= erg->slab_last; j++)
747 int slabIndex = j - erg->slab_first;
748 erg->slab_weights[slabIndex] = get_slab_weight(j, erg, xc, mc, &erg->slab_center[slabIndex]);
750 /* We can do the calculations ONLY if there is weight in the slab! */
751 if (erg->slab_weights[slabIndex] > WEIGHT_MIN)
753 svmul(1.0 / erg->slab_weights[slabIndex], erg->slab_center[slabIndex], erg->slab_center[slabIndex]);
757 /* We need to check this here, since we divide through slab_weights
758 * in the flexible low-level routines! */
759 gmx_fatal(FARGS, "Not enough weight in slab %d. Slab center cannot be determined!", j);
762 /* At first time step: save the centers of the reference structure */
765 copy_rvec(erg->slab_center[slabIndex], erg->slab_center_ref[slabIndex]);
767 } /* END of loop over slabs */
769 /* Output on the master */
770 if ((nullptr != out_slabs) && bOutStep)
772 fprintf(out_slabs, "%12.3e%6d", time, erg->groupIndex);
773 for (int j = erg->slab_first; j <= erg->slab_last; j++)
775 int slabIndex = j - erg->slab_first;
777 "%6d%12.3e%12.3e%12.3e",
779 erg->slab_center[slabIndex][XX],
780 erg->slab_center[slabIndex][YY],
781 erg->slab_center[slabIndex][ZZ]);
783 fprintf(out_slabs, "\n");
788 static void calc_rotmat(const rvec vec,
789 real degangle, /* Angle alpha of rotation at time t in degrees */
790 matrix rotmat) /* Rotation matrix */
792 real radangle; /* Rotation angle in radians */
793 real cosa; /* cosine alpha */
794 real sina; /* sine alpha */
795 real OMcosa; /* 1 - cos(alpha) */
796 real dumxy, dumxz, dumyz; /* save computations */
797 rvec rot_vec; /* Rotate around rot_vec ... */
800 radangle = degangle * M_PI / 180.0;
801 copy_rvec(vec, rot_vec);
803 /* Precompute some variables: */
804 cosa = cos(radangle);
805 sina = sin(radangle);
807 dumxy = rot_vec[XX] * rot_vec[YY] * OMcosa;
808 dumxz = rot_vec[XX] * rot_vec[ZZ] * OMcosa;
809 dumyz = rot_vec[YY] * rot_vec[ZZ] * OMcosa;
811 /* Construct the rotation matrix for this rotation group: */
813 rotmat[XX][XX] = cosa + rot_vec[XX] * rot_vec[XX] * OMcosa;
814 rotmat[YY][XX] = dumxy + rot_vec[ZZ] * sina;
815 rotmat[ZZ][XX] = dumxz - rot_vec[YY] * sina;
817 rotmat[XX][YY] = dumxy - rot_vec[ZZ] * sina;
818 rotmat[YY][YY] = cosa + rot_vec[YY] * rot_vec[YY] * OMcosa;
819 rotmat[ZZ][YY] = dumyz + rot_vec[XX] * sina;
821 rotmat[XX][ZZ] = dumxz + rot_vec[YY] * sina;
822 rotmat[YY][ZZ] = dumyz - rot_vec[XX] * sina;
823 rotmat[ZZ][ZZ] = cosa + rot_vec[ZZ] * rot_vec[ZZ] * OMcosa;
828 for (iii = 0; iii < 3; iii++)
830 for (jjj = 0; jjj < 3; jjj++)
832 fprintf(stderr, " %10.8f ", rotmat[iii][jjj]);
834 fprintf(stderr, "\n");
840 /* Calculates torque on the rotation axis tau = position x force */
841 static inline real torque(const rvec rotvec, /* rotation vector; MUST be normalized! */
842 rvec force, /* force */
843 rvec x, /* position of atom on which the force acts */
844 rvec pivot) /* pivot point of rotation axis */
849 /* Subtract offset */
850 rvec_sub(x, pivot, vectmp);
852 /* position x force */
853 cprod(vectmp, force, tau);
855 /* Return the part of the torque which is parallel to the rotation vector */
856 return iprod(tau, rotvec);
860 /* Right-aligned output of value with standard width */
861 static void print_aligned(FILE* fp, char const* str)
863 fprintf(fp, "%12s", str);
867 /* Right-aligned output of value with standard short width */
868 static void print_aligned_short(FILE* fp, char const* str)
870 fprintf(fp, "%6s", str);
874 static FILE* open_output_file(const char* fn, int steps, const char what[])
879 fp = gmx_ffopen(fn, "w");
881 fprintf(fp, "# Output of %s is written in intervals of %d time step%s.\n#\n", what, steps, steps > 1 ? "s" : "");
887 /* Open output file for slab center data. Call on master only */
888 static FILE* open_slab_out(const char* fn, gmx_enfrot* er)
892 if (er->restartWithAppending)
894 fp = gmx_fio_fopen(fn, "a");
898 fp = open_output_file(fn, er->nstsout, "gaussian weighted slab centers");
900 for (auto& ergRef : er->enfrotgrp)
902 gmx_enfrotgrp* erg = &ergRef;
903 if (ISFLEX(erg->rotg))
906 "# Rotation group %d (%s), slab distance %f nm, %s.\n",
908 enumValueToString(erg->rotg->eType),
909 erg->rotg->slab_dist,
910 erg->rotg->bMassW ? "centers of mass" : "geometrical centers");
914 fprintf(fp, "# Reference centers are listed first (t=-1).\n");
915 fprintf(fp, "# The following columns have the syntax:\n");
917 print_aligned_short(fp, "t");
918 print_aligned_short(fp, "grp");
919 /* Print legend for the first two entries only ... */
920 for (int i = 0; i < 2; i++)
922 print_aligned_short(fp, "slab");
923 print_aligned(fp, "X center");
924 print_aligned(fp, "Y center");
925 print_aligned(fp, "Z center");
927 fprintf(fp, " ...\n");
935 /* Adds 'buf' to 'str' */
936 static void add_to_string(char** str, char* buf)
941 len = strlen(*str) + strlen(buf) + 1;
947 static void add_to_string_aligned(char** str, char* buf)
949 char buf_aligned[STRLEN];
951 sprintf(buf_aligned, "%12s", buf);
952 add_to_string(str, buf_aligned);
956 /* Open output file and print some general information about the rotation groups.
957 * Call on master only */
958 static FILE* open_rot_out(const char* fn, const gmx_output_env_t* oenv, gmx_enfrot* er)
962 const char** setname;
963 char buf[50], buf2[75];
965 char* LegendStr = nullptr;
966 const t_rot* rot = er->rot;
968 if (er->restartWithAppending)
970 fp = gmx_fio_fopen(fn, "a");
975 "Rotation angles and energy",
977 "angles (degrees) and energies (kJ/mol)",
980 "# Output of enforced rotation data is written in intervals of %d time "
983 er->nstrout > 1 ? "s" : "");
985 "# The scalar tau is the torque (kJ/mol) in the direction of the rotation vector "
987 fprintf(fp, "# To obtain the vectorial torque, multiply tau with the group's rot-vec.\n");
989 "# For flexible groups, tau(t,n) from all slabs n have been summed in a single "
990 "value tau(t) here.\n");
991 fprintf(fp, "# The torques tau(t,n) are found in the rottorque.log (-rt) output file\n");
993 for (int g = 0; g < rot->ngrp; g++)
995 const t_rotgrp* rotg = &rot->grp[g];
996 const gmx_enfrotgrp* erg = &er->enfrotgrp[g];
997 bFlex = ISFLEX(rotg);
1000 fprintf(fp, "# ROTATION GROUP %d, potential type '%s':\n", g, enumValueToString(rotg->eType));
1001 fprintf(fp, "# rot-massw%d %s\n", g, booleanValueToString(rotg->bMassW));
1003 "# rot-vec%d %12.5e %12.5e %12.5e\n",
1008 fprintf(fp, "# rot-rate%d %12.5e degrees/ps\n", g, rotg->rate);
1009 fprintf(fp, "# rot-k%d %12.5e kJ/(mol*nm^2)\n", g, rotg->k);
1010 if (rotg->eType == EnforcedRotationGroupType::Iso || rotg->eType == EnforcedRotationGroupType::Pm
1011 || rotg->eType == EnforcedRotationGroupType::Rm
1012 || rotg->eType == EnforcedRotationGroupType::Rm2)
1015 "# rot-pivot%d %12.5e %12.5e %12.5e nm\n",
1024 fprintf(fp, "# rot-slab-distance%d %f nm\n", g, rotg->slab_dist);
1025 fprintf(fp, "# rot-min-gaussian%d %12.5e\n", g, rotg->min_gaussian);
1028 /* Output the centers of the rotation groups for the pivot-free potentials */
1029 if ((rotg->eType == EnforcedRotationGroupType::Isopf)
1030 || (rotg->eType == EnforcedRotationGroupType::Pmpf)
1031 || (rotg->eType == EnforcedRotationGroupType::Rmpf)
1032 || (rotg->eType == EnforcedRotationGroupType::Rm2pf
1033 || (rotg->eType == EnforcedRotationGroupType::Flext)
1034 || (rotg->eType == EnforcedRotationGroupType::Flex2t)))
1037 "# ref. grp. %d center %12.5e %12.5e %12.5e\n",
1039 erg->xc_ref_center[XX],
1040 erg->xc_ref_center[YY],
1041 erg->xc_ref_center[ZZ]);
1044 "# grp. %d init.center %12.5e %12.5e %12.5e\n",
1048 erg->xc_center[ZZ]);
1051 if ((rotg->eType == EnforcedRotationGroupType::Rm2)
1052 || (rotg->eType == EnforcedRotationGroupType::Flex2)
1053 || (rotg->eType == EnforcedRotationGroupType::Flex2t))
1055 fprintf(fp, "# rot-eps%d %12.5e nm^2\n", g, rotg->eps);
1057 if (RotationGroupFitting::Pot == rotg->eFittype)
1061 "# theta_fit%d is determined by first evaluating the potential for %d "
1062 "angles around theta_ref%d.\n",
1064 rotg->PotAngle_nstep,
1067 "# The fit angle is the one with the smallest potential. It is given as "
1070 "# from the reference angle, i.e. if theta_ref=X and theta_fit=Y, then the "
1073 "# minimal value of the potential is X+Y. Angular resolution is %g "
1075 rotg->PotAngle_step);
1079 /* Print a nice legend */
1081 LegendStr[0] = '\0';
1082 sprintf(buf, "# %6s", "time");
1083 add_to_string_aligned(&LegendStr, buf);
1086 snew(setname, 4 * rot->ngrp);
1088 for (int g = 0; g < rot->ngrp; g++)
1090 sprintf(buf, "theta_ref%d", g);
1091 add_to_string_aligned(&LegendStr, buf);
1093 sprintf(buf2, "%s (degrees)", buf);
1094 setname[nsets] = gmx_strdup(buf2);
1097 for (int g = 0; g < rot->ngrp; g++)
1099 const t_rotgrp* rotg = &rot->grp[g];
1100 bFlex = ISFLEX(rotg);
1102 /* For flexible axis rotation we use RMSD fitting to determine the
1103 * actual angle of the rotation group */
1104 if (bFlex || RotationGroupFitting::Pot == rotg->eFittype)
1106 sprintf(buf, "theta_fit%d", g);
1110 sprintf(buf, "theta_av%d", g);
1112 add_to_string_aligned(&LegendStr, buf);
1113 sprintf(buf2, "%s (degrees)", buf);
1114 setname[nsets] = gmx_strdup(buf2);
1117 sprintf(buf, "tau%d", g);
1118 add_to_string_aligned(&LegendStr, buf);
1119 sprintf(buf2, "%s (kJ/mol)", buf);
1120 setname[nsets] = gmx_strdup(buf2);
1123 sprintf(buf, "energy%d", g);
1124 add_to_string_aligned(&LegendStr, buf);
1125 sprintf(buf2, "%s (kJ/mol)", buf);
1126 setname[nsets] = gmx_strdup(buf2);
1133 xvgr_legend(fp, nsets, setname, oenv);
1137 fprintf(fp, "#\n# Legend for the following data columns:\n");
1138 fprintf(fp, "%s\n", LegendStr);
1148 /* Call on master only */
1149 static FILE* open_angles_out(const char* fn, gmx_enfrot* er)
1153 const t_rot* rot = er->rot;
1155 if (er->restartWithAppending)
1157 fp = gmx_fio_fopen(fn, "a");
1161 /* Open output file and write some information about it's structure: */
1162 fp = open_output_file(fn, er->nstsout, "rotation group angles");
1163 fprintf(fp, "# All angles given in degrees, time in ps.\n");
1164 for (int g = 0; g < rot->ngrp; g++)
1166 const t_rotgrp* rotg = &rot->grp[g];
1167 const gmx_enfrotgrp* erg = &er->enfrotgrp[g];
1169 /* Output for this group happens only if potential type is flexible or
1170 * if fit type is potential! */
1171 if (ISFLEX(rotg) || (RotationGroupFitting::Pot == rotg->eFittype))
1175 sprintf(buf, " slab distance %f nm, ", rotg->slab_dist);
1183 "#\n# ROTATION GROUP %d '%s',%s fit type '%s'.\n",
1185 enumValueToString(rotg->eType),
1187 enumValueToString(rotg->eFittype));
1189 /* Special type of fitting using the potential minimum. This is
1190 * done for the whole group only, not for the individual slabs. */
1191 if (RotationGroupFitting::Pot == rotg->eFittype)
1194 "# To obtain theta_fit%d, the potential is evaluated for %d angles "
1195 "around theta_ref%d\n",
1197 rotg->PotAngle_nstep,
1200 "# The fit angle in the rotation standard outfile is the one with "
1201 "minimal energy E(theta_fit) [kJ/mol].\n");
1205 fprintf(fp, "# Legend for the group %d data columns:\n", g);
1207 print_aligned_short(fp, "time");
1208 print_aligned_short(fp, "grp");
1209 print_aligned(fp, "theta_ref");
1211 if (RotationGroupFitting::Pot == rotg->eFittype)
1213 /* Output the set of angles around the reference angle */
1214 for (int i = 0; i < rotg->PotAngle_nstep; i++)
1216 sprintf(buf, "E(%g)", erg->PotAngleFit->degangle[i]);
1217 print_aligned(fp, buf);
1222 /* Output fit angle for each slab */
1223 print_aligned_short(fp, "slab");
1224 print_aligned_short(fp, "atoms");
1225 print_aligned(fp, "theta_fit");
1226 print_aligned_short(fp, "slab");
1227 print_aligned_short(fp, "atoms");
1228 print_aligned(fp, "theta_fit");
1229 fprintf(fp, " ...");
1241 /* Open torque output file and write some information about it's structure.
1242 * Call on master only */
1243 static FILE* open_torque_out(const char* fn, gmx_enfrot* er)
1246 const t_rot* rot = er->rot;
1248 if (er->restartWithAppending)
1250 fp = gmx_fio_fopen(fn, "a");
1254 fp = open_output_file(fn, er->nstsout, "torques");
1256 for (int g = 0; g < rot->ngrp; g++)
1258 const t_rotgrp* rotg = &rot->grp[g];
1259 const gmx_enfrotgrp* erg = &er->enfrotgrp[g];
1263 "# Rotation group %d (%s), slab distance %f nm.\n",
1265 enumValueToString(rotg->eType),
1268 "# The scalar tau is the torque (kJ/mol) in the direction of the rotation "
1270 fprintf(fp, "# To obtain the vectorial torque, multiply tau with\n");
1272 "# rot-vec%d %10.3e %10.3e %10.3e\n",
1280 fprintf(fp, "# Legend for the following data columns: (tau=torque for that slab):\n");
1282 print_aligned_short(fp, "t");
1283 print_aligned_short(fp, "grp");
1284 print_aligned_short(fp, "slab");
1285 print_aligned(fp, "tau");
1286 print_aligned_short(fp, "slab");
1287 print_aligned(fp, "tau");
1288 fprintf(fp, " ...\n");
1296 static void swap_val(double* vec, int i, int j)
1298 double tmp = vec[j];
1306 static void swap_col(double** mat, int i, int j)
1308 double tmp[3] = { mat[0][j], mat[1][j], mat[2][j] };
1311 mat[0][j] = mat[0][i];
1312 mat[1][j] = mat[1][i];
1313 mat[2][j] = mat[2][i];
1321 /* Eigenvectors are stored in columns of eigen_vec */
1322 static void diagonalize_symmetric(double** matrix, double** eigen_vec, double eigenval[3])
1327 jacobi(matrix, 3, eigenval, eigen_vec, &n_rot);
1329 /* sort in ascending order */
1330 if (eigenval[0] > eigenval[1])
1332 swap_val(eigenval, 0, 1);
1333 swap_col(eigen_vec, 0, 1);
1335 if (eigenval[1] > eigenval[2])
1337 swap_val(eigenval, 1, 2);
1338 swap_col(eigen_vec, 1, 2);
1340 if (eigenval[0] > eigenval[1])
1342 swap_val(eigenval, 0, 1);
1343 swap_col(eigen_vec, 0, 1);
1348 static void align_with_z(rvec* s, /* Structure to align */
1353 rvec zet = { 0.0, 0.0, 1.0 };
1354 rvec rot_axis = { 0.0, 0.0, 0.0 };
1355 rvec* rotated_str = nullptr;
1361 snew(rotated_str, natoms);
1363 /* Normalize the axis */
1364 ooanorm = 1.0 / norm(axis);
1365 svmul(ooanorm, axis, axis);
1367 /* Calculate the angle for the fitting procedure */
1368 cprod(axis, zet, rot_axis);
1369 angle = acos(axis[2]);
1375 /* Calculate the rotation matrix */
1376 calc_rotmat(rot_axis, angle * 180.0 / M_PI, rotmat);
1378 /* Apply the rotation matrix to s */
1379 for (i = 0; i < natoms; i++)
1381 for (j = 0; j < 3; j++)
1383 for (k = 0; k < 3; k++)
1385 rotated_str[i][j] += rotmat[j][k] * s[i][k];
1390 /* Rewrite the rotated structure to s */
1391 for (i = 0; i < natoms; i++)
1393 for (j = 0; j < 3; j++)
1395 s[i][j] = rotated_str[i][j];
1403 static void calc_correl_matrix(rvec* Xstr, rvec* Ystr, double** Rmat, int natoms)
1408 for (i = 0; i < 3; i++)
1410 for (j = 0; j < 3; j++)
1416 for (i = 0; i < 3; i++)
1418 for (j = 0; j < 3; j++)
1420 for (k = 0; k < natoms; k++)
1422 Rmat[i][j] += Ystr[k][i] * Xstr[k][j];
1429 static void weigh_coords(rvec* str, real* weight, int natoms)
1434 for (i = 0; i < natoms; i++)
1436 for (j = 0; j < 3; j++)
1438 str[i][j] *= std::sqrt(weight[i]);
1444 static real opt_angle_analytic(rvec* ref_s,
1453 rvec* ref_s_1 = nullptr;
1454 rvec* act_s_1 = nullptr;
1456 double **Rmat, **RtR, **eigvec;
1458 double V[3][3], WS[3][3];
1459 double rot_matrix[3][3];
1463 /* Do not change the original coordinates */
1464 snew(ref_s_1, natoms);
1465 snew(act_s_1, natoms);
1466 for (i = 0; i < natoms; i++)
1468 copy_rvec(ref_s[i], ref_s_1[i]);
1469 copy_rvec(act_s[i], act_s_1[i]);
1472 /* Translate the structures to the origin */
1473 shift[XX] = -ref_com[XX];
1474 shift[YY] = -ref_com[YY];
1475 shift[ZZ] = -ref_com[ZZ];
1476 translate_x(ref_s_1, natoms, shift);
1478 shift[XX] = -act_com[XX];
1479 shift[YY] = -act_com[YY];
1480 shift[ZZ] = -act_com[ZZ];
1481 translate_x(act_s_1, natoms, shift);
1483 /* Align rotation axis with z */
1484 align_with_z(ref_s_1, natoms, axis);
1485 align_with_z(act_s_1, natoms, axis);
1487 /* Correlation matrix */
1488 Rmat = allocate_square_matrix(3);
1490 for (i = 0; i < natoms; i++)
1492 ref_s_1[i][2] = 0.0;
1493 act_s_1[i][2] = 0.0;
1496 /* Weight positions with sqrt(weight) */
1497 if (nullptr != weight)
1499 weigh_coords(ref_s_1, weight, natoms);
1500 weigh_coords(act_s_1, weight, natoms);
1503 /* Calculate correlation matrices R=YXt (X=ref_s; Y=act_s) */
1504 calc_correl_matrix(ref_s_1, act_s_1, Rmat, natoms);
1507 RtR = allocate_square_matrix(3);
1508 for (i = 0; i < 3; i++)
1510 for (j = 0; j < 3; j++)
1512 for (k = 0; k < 3; k++)
1514 RtR[i][j] += Rmat[k][i] * Rmat[k][j];
1518 /* Diagonalize RtR */
1520 for (i = 0; i < 3; i++)
1525 diagonalize_symmetric(RtR, eigvec, eigval);
1526 swap_col(eigvec, 0, 1);
1527 swap_col(eigvec, 1, 2);
1528 swap_val(eigval, 0, 1);
1529 swap_val(eigval, 1, 2);
1532 for (i = 0; i < 3; i++)
1534 for (j = 0; j < 3; j++)
1541 for (i = 0; i < 2; i++)
1543 for (j = 0; j < 2; j++)
1545 WS[i][j] = eigvec[i][j] / std::sqrt(eigval[j]);
1549 for (i = 0; i < 3; i++)
1551 for (j = 0; j < 3; j++)
1553 for (k = 0; k < 3; k++)
1555 V[i][j] += Rmat[i][k] * WS[k][j];
1559 free_square_matrix(Rmat, 3);
1561 /* Calculate optimal rotation matrix */
1562 for (i = 0; i < 3; i++)
1564 for (j = 0; j < 3; j++)
1566 rot_matrix[i][j] = 0.0;
1570 for (i = 0; i < 3; i++)
1572 for (j = 0; j < 3; j++)
1574 for (k = 0; k < 3; k++)
1576 rot_matrix[i][j] += eigvec[i][k] * V[j][k];
1580 rot_matrix[2][2] = 1.0;
1582 /* In some cases abs(rot_matrix[0][0]) can be slighly larger
1583 * than unity due to numerical inacurracies. To be able to calculate
1584 * the acos function, we put these values back in range. */
1585 if (rot_matrix[0][0] > 1.0)
1587 rot_matrix[0][0] = 1.0;
1589 else if (rot_matrix[0][0] < -1.0)
1591 rot_matrix[0][0] = -1.0;
1594 /* Determine the optimal rotation angle: */
1595 opt_angle = (-1.0) * acos(rot_matrix[0][0]) * 180.0 / M_PI;
1596 if (rot_matrix[0][1] < 0.0)
1598 opt_angle = (-1.0) * opt_angle;
1601 /* Give back some memory */
1602 free_square_matrix(RtR, 3);
1605 for (i = 0; i < 3; i++)
1611 return static_cast<real>(opt_angle);
1615 /* Determine angle of the group by RMSD fit to the reference */
1616 /* Not parallelized, call this routine only on the master */
1617 static real flex_fit_angle(gmx_enfrotgrp* erg)
1619 rvec* fitcoords = nullptr;
1620 rvec center; /* Center of positions passed to the fit routine */
1621 real fitangle; /* Angle of the rotation group derived by fitting */
1625 /* Get the center of the rotation group.
1626 * Note, again, erg->xc has been sorted in do_flexible */
1627 get_center(erg->xc, erg->mc_sorted, erg->rotg->nat, center);
1629 /* === Determine the optimal fit angle for the rotation group === */
1630 if (erg->rotg->eFittype == RotationGroupFitting::Norm)
1632 /* Normalize every position to it's reference length */
1633 for (int i = 0; i < erg->rotg->nat; i++)
1635 /* Put the center of the positions into the origin */
1636 rvec_sub(erg->xc[i], center, coord);
1637 /* Determine the scaling factor for the length: */
1638 scal = erg->xc_ref_length[erg->xc_sortind[i]] / norm(coord);
1639 /* Get position, multiply with the scaling factor and save */
1640 svmul(scal, coord, erg->xc_norm[i]);
1642 fitcoords = erg->xc_norm;
1646 fitcoords = erg->xc;
1648 /* From the point of view of the current positions, the reference has rotated
1649 * backwards. Since we output the angle relative to the fixed reference,
1650 * we need the minus sign. */
1651 fitangle = -opt_angle_analytic(
1652 erg->xc_ref_sorted, fitcoords, erg->mc_sorted, erg->rotg->nat, erg->xc_ref_center, center, erg->vec);
1658 /* Determine actual angle of each slab by RMSD fit to the reference */
1659 /* Not parallelized, call this routine only on the master */
1660 static void flex_fit_angle_perslab(gmx_enfrotgrp* erg, double t, real degangle, FILE* fp)
1663 rvec act_center; /* Center of actual positions that are passed to the fit routine */
1664 rvec ref_center; /* Same for the reference positions */
1665 real fitangle; /* Angle of a slab derived from an RMSD fit to
1666 * the reference structure at t=0 */
1668 real OOm_av; /* 1/average_mass of a rotation group atom */
1669 real m_rel; /* Relative mass of a rotation group atom */
1672 /* Average mass of a rotation group atom: */
1673 OOm_av = erg->invmass * erg->rotg->nat;
1675 /**********************************/
1676 /* First collect the data we need */
1677 /**********************************/
1679 /* Collect the data for the individual slabs */
1680 for (int n = erg->slab_first; n <= erg->slab_last; n++)
1682 int slabIndex = n - erg->slab_first; /* slab index */
1683 sd = &(erg->slab_data[slabIndex]);
1684 sd->nat = erg->lastatom[slabIndex] - erg->firstatom[slabIndex] + 1;
1687 /* Loop over the relevant atoms in the slab */
1688 for (int l = erg->firstatom[slabIndex]; l <= erg->lastatom[slabIndex]; l++)
1690 /* Current position of this atom: x[ii][XX/YY/ZZ] */
1691 copy_rvec(erg->xc[l], curr_x);
1693 /* The (unrotated) reference position of this atom is copied to ref_x.
1694 * Beware, the xc coords have been sorted in do_flexible */
1695 copy_rvec(erg->xc_ref_sorted[l], ref_x);
1697 /* Save data for doing angular RMSD fit later */
1698 /* Save the current atom position */
1699 copy_rvec(curr_x, sd->x[ind]);
1700 /* Save the corresponding reference position */
1701 copy_rvec(ref_x, sd->ref[ind]);
1703 /* Maybe also mass-weighting was requested. If yes, additionally
1704 * multiply the weights with the relative mass of the atom. If not,
1705 * multiply with unity. */
1706 m_rel = erg->mc_sorted[l] * OOm_av;
1708 /* Save the weight for this atom in this slab */
1709 sd->weight[ind] = gaussian_weight(curr_x, erg, n) * m_rel;
1711 /* Next atom in this slab */
1716 /******************************/
1717 /* Now do the fit calculation */
1718 /******************************/
1720 fprintf(fp, "%12.3e%6d%12.3f", t, erg->groupIndex, degangle);
1722 /* === Now do RMSD fitting for each slab === */
1723 /* We require at least SLAB_MIN_ATOMS in a slab, such that the fit makes sense. */
1724 #define SLAB_MIN_ATOMS 4
1726 for (int n = erg->slab_first; n <= erg->slab_last; n++)
1728 int slabIndex = n - erg->slab_first; /* slab index */
1729 sd = &(erg->slab_data[slabIndex]);
1730 if (sd->nat >= SLAB_MIN_ATOMS)
1732 /* Get the center of the slabs reference and current positions */
1733 get_center(sd->ref, sd->weight, sd->nat, ref_center);
1734 get_center(sd->x, sd->weight, sd->nat, act_center);
1735 if (erg->rotg->eFittype == RotationGroupFitting::Norm)
1737 /* Normalize every position to it's reference length
1738 * prior to performing the fit */
1739 for (int i = 0; i < sd->nat; i++) /* Center */
1741 rvec_dec(sd->ref[i], ref_center);
1742 rvec_dec(sd->x[i], act_center);
1743 /* Normalize x_i such that it gets the same length as ref_i */
1744 svmul(norm(sd->ref[i]) / norm(sd->x[i]), sd->x[i], sd->x[i]);
1746 /* We already subtracted the centers */
1747 clear_rvec(ref_center);
1748 clear_rvec(act_center);
1750 fitangle = -opt_angle_analytic(
1751 sd->ref, sd->x, sd->weight, sd->nat, ref_center, act_center, erg->vec);
1752 fprintf(fp, "%6d%6d%12.3f", n, sd->nat, fitangle);
1757 #undef SLAB_MIN_ATOMS
1761 /* Shift x with is */
1762 static inline void shift_single_coord(const matrix box, rvec x, const ivec is)
1773 x[XX] += tx * box[XX][XX] + ty * box[YY][XX] + tz * box[ZZ][XX];
1774 x[YY] += ty * box[YY][YY] + tz * box[ZZ][YY];
1775 x[ZZ] += tz * box[ZZ][ZZ];
1779 x[XX] += tx * box[XX][XX];
1780 x[YY] += ty * box[YY][YY];
1781 x[ZZ] += tz * box[ZZ][ZZ];
1786 /* Determine the 'home' slab of this atom which is the
1787 * slab with the highest Gaussian weight of all */
1788 static inline int get_homeslab(rvec curr_x, /* The position for which the home slab shall be determined */
1789 const rvec rotvec, /* The rotation vector */
1790 real slabdist) /* The slab distance */
1795 /* The distance of the atom to the coordinate center (where the
1796 * slab with index 0) is */
1797 dist = iprod(rotvec, curr_x);
1799 return gmx::roundToInt(dist / slabdist);
1803 /* For a local atom determine the relevant slabs, i.e. slabs in
1804 * which the gaussian is larger than min_gaussian
1806 static int get_single_atom_gaussians(rvec curr_x, gmx_enfrotgrp* erg)
1809 /* Determine the 'home' slab of this atom: */
1810 int homeslab = get_homeslab(curr_x, erg->vec, erg->rotg->slab_dist);
1812 /* First determine the weight in the atoms home slab: */
1813 real g = gaussian_weight(curr_x, erg, homeslab);
1815 erg->gn_atom[count] = g;
1816 erg->gn_slabind[count] = homeslab;
1820 /* Determine the max slab */
1821 int slab = homeslab;
1822 while (g > erg->rotg->min_gaussian)
1825 g = gaussian_weight(curr_x, erg, slab);
1826 erg->gn_slabind[count] = slab;
1827 erg->gn_atom[count] = g;
1832 /* Determine the min slab */
1837 g = gaussian_weight(curr_x, erg, slab);
1838 erg->gn_slabind[count] = slab;
1839 erg->gn_atom[count] = g;
1841 } while (g > erg->rotg->min_gaussian);
1848 static void flex2_precalc_inner_sum(const gmx_enfrotgrp* erg)
1850 rvec xi; /* positions in the i-sum */
1851 rvec xcn, ycn; /* the current and the reference slab centers */
1854 rvec rin; /* Helper variables */
1857 real OOpsii, OOpsiistar;
1858 real sin_rin; /* s_ii.r_ii */
1859 rvec s_in, tmpvec, tmpvec2;
1860 real mi, wi; /* Mass-weighting of the positions */
1864 N_M = erg->rotg->nat * erg->invmass;
1866 /* Loop over all slabs that contain something */
1867 for (int n = erg->slab_first; n <= erg->slab_last; n++)
1869 int slabIndex = n - erg->slab_first; /* slab index */
1871 /* The current center of this slab is saved in xcn: */
1872 copy_rvec(erg->slab_center[slabIndex], xcn);
1873 /* ... and the reference center in ycn: */
1874 copy_rvec(erg->slab_center_ref[slabIndex + erg->slab_buffer], ycn);
1876 /*** D. Calculate the whole inner sum used for second and third sum */
1877 /* For slab n, we need to loop over all atoms i again. Since we sorted
1878 * the atoms with respect to the rotation vector, we know that it is sufficient
1879 * to calculate from firstatom to lastatom only. All other contributions will
1881 clear_rvec(innersumvec);
1882 for (int i = erg->firstatom[slabIndex]; i <= erg->lastatom[slabIndex]; i++)
1884 /* Coordinate xi of this atom */
1885 copy_rvec(erg->xc[i], xi);
1888 gaussian_xi = gaussian_weight(xi, erg, n);
1889 mi = erg->mc_sorted[i]; /* need the sorted mass here */
1893 copy_rvec(erg->xc_ref_sorted[i], yi0); /* Reference position yi0 */
1894 rvec_sub(yi0, ycn, tmpvec2); /* tmpvec2 = yi0 - ycn */
1895 mvmul(erg->rotmat, tmpvec2, rin); /* rin = Omega.(yi0 - ycn) */
1897 /* Calculate psi_i* and sin */
1898 rvec_sub(xi, xcn, tmpvec2); /* tmpvec2 = xi - xcn */
1900 /* In rare cases, when an atom position coincides with a slab center
1901 * (tmpvec2 == 0) we cannot compute the vector product for s_in.
1902 * However, since the atom is located directly on the pivot, this
1903 * slab's contribution to the force on that atom will be zero
1904 * anyway. Therefore, we continue with the next atom. */
1905 if (gmx_numzero(norm(tmpvec2))) /* 0 == norm(xi - xcn) */
1910 cprod(erg->vec, tmpvec2, tmpvec); /* tmpvec = v x (xi - xcn) */
1911 OOpsiistar = norm2(tmpvec) + erg->rotg->eps; /* OOpsii* = 1/psii* = |v x (xi-xcn)|^2 + eps */
1912 OOpsii = norm(tmpvec); /* OOpsii = 1 / psii = |v x (xi - xcn)| */
1914 /* * v x (xi - xcn) */
1915 unitv(tmpvec, s_in); /* sin = ---------------- */
1916 /* |v x (xi - xcn)| */
1918 sin_rin = iprod(s_in, rin); /* sin_rin = sin . rin */
1920 /* Now the whole sum */
1921 fac = OOpsii / OOpsiistar;
1922 svmul(fac, rin, tmpvec);
1923 fac2 = fac * fac * OOpsii;
1924 svmul(fac2 * sin_rin, s_in, tmpvec2);
1925 rvec_dec(tmpvec, tmpvec2);
1927 svmul(wi * gaussian_xi * sin_rin, tmpvec, tmpvec2);
1929 rvec_inc(innersumvec, tmpvec2);
1930 } /* now we have the inner sum, used both for sum2 and sum3 */
1932 /* Save it to be used in do_flex2_lowlevel */
1933 copy_rvec(innersumvec, erg->slab_innersumvec[slabIndex]);
1934 } /* END of loop over slabs */
1938 static void flex_precalc_inner_sum(const gmx_enfrotgrp* erg)
1940 rvec xi; /* position */
1941 rvec xcn, ycn; /* the current and the reference slab centers */
1942 rvec qin, rin; /* q_i^n and r_i^n */
1945 rvec innersumvec; /* Inner part of sum_n2 */
1946 real gaussian_xi; /* Gaussian weight gn(xi) */
1947 real mi, wi; /* Mass-weighting of the positions */
1950 N_M = erg->rotg->nat * erg->invmass;
1952 /* Loop over all slabs that contain something */
1953 for (int n = erg->slab_first; n <= erg->slab_last; n++)
1955 int slabIndex = n - erg->slab_first; /* slab index */
1957 /* The current center of this slab is saved in xcn: */
1958 copy_rvec(erg->slab_center[slabIndex], xcn);
1959 /* ... and the reference center in ycn: */
1960 copy_rvec(erg->slab_center_ref[slabIndex + erg->slab_buffer], ycn);
1962 /* For slab n, we need to loop over all atoms i again. Since we sorted
1963 * the atoms with respect to the rotation vector, we know that it is sufficient
1964 * to calculate from firstatom to lastatom only. All other contributions will
1966 clear_rvec(innersumvec);
1967 for (int i = erg->firstatom[slabIndex]; i <= erg->lastatom[slabIndex]; i++)
1969 /* Coordinate xi of this atom */
1970 copy_rvec(erg->xc[i], xi);
1973 gaussian_xi = gaussian_weight(xi, erg, n);
1974 mi = erg->mc_sorted[i]; /* need the sorted mass here */
1977 /* Calculate rin and qin */
1978 rvec_sub(erg->xc_ref_sorted[i], ycn, tmpvec); /* tmpvec = yi0-ycn */
1980 /* In rare cases, when an atom position coincides with a slab center
1981 * (tmpvec == 0) we cannot compute the vector product for qin.
1982 * However, since the atom is located directly on the pivot, this
1983 * slab's contribution to the force on that atom will be zero
1984 * anyway. Therefore, we continue with the next atom. */
1985 if (gmx_numzero(norm(tmpvec))) /* 0 == norm(yi0 - ycn) */
1990 mvmul(erg->rotmat, tmpvec, rin); /* rin = Omega.(yi0 - ycn) */
1991 cprod(erg->vec, rin, tmpvec); /* tmpvec = v x Omega*(yi0-ycn) */
1993 /* * v x Omega*(yi0-ycn) */
1994 unitv(tmpvec, qin); /* qin = --------------------- */
1995 /* |v x Omega*(yi0-ycn)| */
1998 rvec_sub(xi, xcn, tmpvec); /* tmpvec = xi-xcn */
1999 bin = iprod(qin, tmpvec); /* bin = qin*(xi-xcn) */
2001 svmul(wi * gaussian_xi * bin, qin, tmpvec);
2003 /* Add this contribution to the inner sum: */
2004 rvec_add(innersumvec, tmpvec, innersumvec);
2005 } /* now we have the inner sum vector S^n for this slab */
2006 /* Save it to be used in do_flex_lowlevel */
2007 copy_rvec(innersumvec, erg->slab_innersumvec[slabIndex]);
2012 static real do_flex2_lowlevel(gmx_enfrotgrp* erg,
2013 real sigma, /* The Gaussian width sigma */
2015 gmx_bool bOutstepRot,
2016 gmx_bool bOutstepSlab,
2019 int count, ii, iigrp;
2020 rvec xj; /* position in the i-sum */
2021 rvec yj0; /* the reference position in the j-sum */
2022 rvec xcn, ycn; /* the current and the reference slab centers */
2023 real V; /* This node's part of the rotation pot. energy */
2024 real gaussian_xj; /* Gaussian weight */
2027 real numerator, fit_numerator;
2028 rvec rjn, fit_rjn; /* Helper variables */
2031 real OOpsij, OOpsijstar;
2032 real OOsigma2; /* 1/(sigma^2) */
2035 rvec sjn, tmpvec, tmpvec2, yj0_ycn;
2036 rvec sum1vec_part, sum1vec, sum2vec_part, sum2vec, sum3vec, sum4vec, innersumvec;
2038 real mj, wj; /* Mass-weighting of the positions */
2040 real Wjn; /* g_n(x_j) m_j / Mjn */
2041 gmx_bool bCalcPotFit;
2043 /* To calculate the torque per slab */
2044 rvec slab_force; /* Single force from slab n on one atom */
2045 rvec slab_sum1vec_part;
2046 real slab_sum3part, slab_sum4part;
2047 rvec slab_sum1vec, slab_sum2vec, slab_sum3vec, slab_sum4vec;
2049 /* Pre-calculate the inner sums, so that we do not have to calculate
2050 * them again for every atom */
2051 flex2_precalc_inner_sum(erg);
2053 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (RotationGroupFitting::Pot == erg->rotg->eFittype);
2055 /********************************************************/
2056 /* Main loop over all local atoms of the rotation group */
2057 /********************************************************/
2058 N_M = erg->rotg->nat * erg->invmass;
2060 OOsigma2 = 1.0 / (sigma * sigma);
2061 const auto& localRotationGroupIndex = erg->atomSet->localIndex();
2062 const auto& collectiveRotationGroupIndex = erg->atomSet->collectiveIndex();
2064 for (gmx::index j = 0; j < localRotationGroupIndex.ssize(); j++)
2066 /* Local index of a rotation group atom */
2067 ii = localRotationGroupIndex[j];
2068 /* Position of this atom in the collective array */
2069 iigrp = collectiveRotationGroupIndex[j];
2070 /* Mass-weighting */
2071 mj = erg->mc[iigrp]; /* need the unsorted mass here */
2074 /* Current position of this atom: x[ii][XX/YY/ZZ]
2075 * Note that erg->xc_center contains the center of mass in case the flex2-t
2076 * potential was chosen. For the flex2 potential erg->xc_center must be
2078 rvec_sub(x[ii], erg->xc_center, xj);
2080 /* Shift this atom such that it is near its reference */
2081 shift_single_coord(box, xj, erg->xc_shifts[iigrp]);
2083 /* Determine the slabs to loop over, i.e. the ones with contributions
2084 * larger than min_gaussian */
2085 count = get_single_atom_gaussians(xj, erg);
2087 clear_rvec(sum1vec_part);
2088 clear_rvec(sum2vec_part);
2091 /* Loop over the relevant slabs for this atom */
2092 for (int ic = 0; ic < count; ic++)
2094 int n = erg->gn_slabind[ic];
2096 /* Get the precomputed Gaussian value of curr_slab for curr_x */
2097 gaussian_xj = erg->gn_atom[ic];
2099 int slabIndex = n - erg->slab_first; /* slab index */
2101 /* The (unrotated) reference position of this atom is copied to yj0: */
2102 copy_rvec(erg->rotg->x_ref[iigrp], yj0);
2104 beta = calc_beta(xj, erg, n);
2106 /* The current center of this slab is saved in xcn: */
2107 copy_rvec(erg->slab_center[slabIndex], xcn);
2108 /* ... and the reference center in ycn: */
2109 copy_rvec(erg->slab_center_ref[slabIndex + erg->slab_buffer], ycn);
2111 rvec_sub(yj0, ycn, yj0_ycn); /* yj0_ycn = yj0 - ycn */
2114 mvmul(erg->rotmat, yj0_ycn, rjn); /* rjn = Omega.(yj0 - ycn) */
2116 /* Subtract the slab center from xj */
2117 rvec_sub(xj, xcn, tmpvec2); /* tmpvec2 = xj - xcn */
2119 /* In rare cases, when an atom position coincides with a slab center
2120 * (tmpvec2 == 0) we cannot compute the vector product for sjn.
2121 * However, since the atom is located directly on the pivot, this
2122 * slab's contribution to the force on that atom will be zero
2123 * anyway. Therefore, we directly move on to the next slab. */
2124 if (gmx_numzero(norm(tmpvec2))) /* 0 == norm(xj - xcn) */
2130 cprod(erg->vec, tmpvec2, tmpvec); /* tmpvec = v x (xj - xcn) */
2132 OOpsijstar = norm2(tmpvec) + erg->rotg->eps; /* OOpsij* = 1/psij* = |v x (xj-xcn)|^2 + eps */
2134 numerator = gmx::square(iprod(tmpvec, rjn));
2136 /*********************************/
2137 /* Add to the rotation potential */
2138 /*********************************/
2139 V += 0.5 * erg->rotg->k * wj * gaussian_xj * numerator / OOpsijstar;
2141 /* If requested, also calculate the potential for a set of angles
2142 * near the current reference angle */
2145 for (int ifit = 0; ifit < erg->rotg->PotAngle_nstep; ifit++)
2147 mvmul(erg->PotAngleFit->rotmat[ifit], yj0_ycn, fit_rjn);
2148 fit_numerator = gmx::square(iprod(tmpvec, fit_rjn));
2149 erg->PotAngleFit->V[ifit] +=
2150 0.5 * erg->rotg->k * wj * gaussian_xj * fit_numerator / OOpsijstar;
2154 /*************************************/
2155 /* Now calculate the force on atom j */
2156 /*************************************/
2158 OOpsij = norm(tmpvec); /* OOpsij = 1 / psij = |v x (xj - xcn)| */
2160 /* * v x (xj - xcn) */
2161 unitv(tmpvec, sjn); /* sjn = ---------------- */
2162 /* |v x (xj - xcn)| */
2164 sjn_rjn = iprod(sjn, rjn); /* sjn_rjn = sjn . rjn */
2167 /*** A. Calculate the first of the four sum terms: ****************/
2168 fac = OOpsij / OOpsijstar;
2169 svmul(fac, rjn, tmpvec);
2170 fac2 = fac * fac * OOpsij;
2171 svmul(fac2 * sjn_rjn, sjn, tmpvec2);
2172 rvec_dec(tmpvec, tmpvec2);
2173 fac2 = wj * gaussian_xj; /* also needed for sum4 */
2174 svmul(fac2 * sjn_rjn, tmpvec, slab_sum1vec_part);
2175 /********************/
2176 /*** Add to sum1: ***/
2177 /********************/
2178 rvec_inc(sum1vec_part, slab_sum1vec_part); /* sum1 still needs to vector multiplied with v */
2180 /*** B. Calculate the forth of the four sum terms: ****************/
2181 betasigpsi = beta * OOsigma2 * OOpsij; /* this is also needed for sum3 */
2182 /********************/
2183 /*** Add to sum4: ***/
2184 /********************/
2185 slab_sum4part = fac2 * betasigpsi * fac * sjn_rjn
2186 * sjn_rjn; /* Note that fac is still valid from above */
2187 sum4 += slab_sum4part;
2189 /*** C. Calculate Wjn for second and third sum */
2190 /* Note that we can safely divide by slab_weights since we check in
2191 * get_slab_centers that it is non-zero. */
2192 Wjn = gaussian_xj * mj / erg->slab_weights[slabIndex];
2194 /* We already have precalculated the inner sum for slab n */
2195 copy_rvec(erg->slab_innersumvec[slabIndex], innersumvec);
2197 /* Weigh the inner sum vector with Wjn */
2198 svmul(Wjn, innersumvec, innersumvec);
2200 /*** E. Calculate the second of the four sum terms: */
2201 /********************/
2202 /*** Add to sum2: ***/
2203 /********************/
2204 rvec_inc(sum2vec_part, innersumvec); /* sum2 still needs to be vector crossproduct'ed with v */
2206 /*** F. Calculate the third of the four sum terms: */
2207 slab_sum3part = betasigpsi * iprod(sjn, innersumvec);
2208 sum3 += slab_sum3part; /* still needs to be multiplied with v */
2210 /*** G. Calculate the torque on the local slab's axis: */
2214 cprod(slab_sum1vec_part, erg->vec, slab_sum1vec);
2216 cprod(innersumvec, erg->vec, slab_sum2vec);
2218 svmul(slab_sum3part, erg->vec, slab_sum3vec);
2220 svmul(slab_sum4part, erg->vec, slab_sum4vec);
2222 /* The force on atom ii from slab n only: */
2223 for (int m = 0; m < DIM; m++)
2225 slab_force[m] = erg->rotg->k
2226 * (-slab_sum1vec[m] + slab_sum2vec[m] - slab_sum3vec[m]
2227 + 0.5 * slab_sum4vec[m]);
2230 erg->slab_torque_v[slabIndex] += torque(erg->vec, slab_force, xj, xcn);
2232 } /* END of loop over slabs */
2234 /* Construct the four individual parts of the vector sum: */
2235 cprod(sum1vec_part, erg->vec, sum1vec); /* sum1vec = { } x v */
2236 cprod(sum2vec_part, erg->vec, sum2vec); /* sum2vec = { } x v */
2237 svmul(sum3, erg->vec, sum3vec); /* sum3vec = { } . v */
2238 svmul(sum4, erg->vec, sum4vec); /* sum4vec = { } . v */
2240 /* Store the additional force so that it can be added to the force
2241 * array after the normal forces have been evaluated */
2242 for (int m = 0; m < DIM; m++)
2244 erg->f_rot_loc[j][m] =
2245 erg->rotg->k * (-sum1vec[m] + sum2vec[m] - sum3vec[m] + 0.5 * sum4vec[m]);
2250 "sum1: %15.8f %15.8f %15.8f\n",
2251 -erg->rotg->k * sum1vec[XX],
2252 -erg->rotg->k * sum1vec[YY],
2253 -erg->rotg->k * sum1vec[ZZ]);
2255 "sum2: %15.8f %15.8f %15.8f\n",
2256 erg->rotg->k * sum2vec[XX],
2257 erg->rotg->k * sum2vec[YY],
2258 erg->rotg->k * sum2vec[ZZ]);
2260 "sum3: %15.8f %15.8f %15.8f\n",
2261 -erg->rotg->k * sum3vec[XX],
2262 -erg->rotg->k * sum3vec[YY],
2263 -erg->rotg->k * sum3vec[ZZ]);
2265 "sum4: %15.8f %15.8f %15.8f\n",
2266 0.5 * erg->rotg->k * sum4vec[XX],
2267 0.5 * erg->rotg->k * sum4vec[YY],
2268 0.5 * erg->rotg->k * sum4vec[ZZ]);
2273 } /* END of loop over local atoms */
2279 static real do_flex_lowlevel(gmx_enfrotgrp* erg,
2280 real sigma, /* The Gaussian width sigma */
2282 gmx_bool bOutstepRot,
2283 gmx_bool bOutstepSlab,
2287 rvec xj, yj0; /* current and reference position */
2288 rvec xcn, ycn; /* the current and the reference slab centers */
2289 rvec yj0_ycn; /* yj0 - ycn */
2290 rvec xj_xcn; /* xj - xcn */
2291 rvec qjn, fit_qjn; /* q_i^n */
2292 rvec sum_n1, sum_n2; /* Two contributions to the rotation force */
2293 rvec innersumvec; /* Inner part of sum_n2 */
2295 rvec force_n; /* Single force from slab n on one atom */
2296 rvec force_n1, force_n2; /* First and second part of force_n */
2297 rvec tmpvec, tmpvec2, tmp_f; /* Helper variables */
2298 real V; /* The rotation potential energy */
2299 real OOsigma2; /* 1/(sigma^2) */
2300 real beta; /* beta_n(xj) */
2301 real bjn, fit_bjn; /* b_j^n */
2302 real gaussian_xj; /* Gaussian weight gn(xj) */
2303 real betan_xj_sigma2;
2304 real mj, wj; /* Mass-weighting of the positions */
2306 gmx_bool bCalcPotFit;
2308 /* Pre-calculate the inner sums, so that we do not have to calculate
2309 * them again for every atom */
2310 flex_precalc_inner_sum(erg);
2312 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (RotationGroupFitting::Pot == erg->rotg->eFittype);
2314 /********************************************************/
2315 /* Main loop over all local atoms of the rotation group */
2316 /********************************************************/
2317 OOsigma2 = 1.0 / (sigma * sigma);
2318 N_M = erg->rotg->nat * erg->invmass;
2320 const auto& localRotationGroupIndex = erg->atomSet->localIndex();
2321 const auto& collectiveRotationGroupIndex = erg->atomSet->collectiveIndex();
2323 for (gmx::index j = 0; j < localRotationGroupIndex.ssize(); j++)
2325 /* Local index of a rotation group atom */
2326 int ii = localRotationGroupIndex[j];
2327 /* Position of this atom in the collective array */
2328 iigrp = collectiveRotationGroupIndex[j];
2329 /* Mass-weighting */
2330 mj = erg->mc[iigrp]; /* need the unsorted mass here */
2333 /* Current position of this atom: x[ii][XX/YY/ZZ]
2334 * Note that erg->xc_center contains the center of mass in case the flex-t
2335 * potential was chosen. For the flex potential erg->xc_center must be
2337 rvec_sub(x[ii], erg->xc_center, xj);
2339 /* Shift this atom such that it is near its reference */
2340 shift_single_coord(box, xj, erg->xc_shifts[iigrp]);
2342 /* Determine the slabs to loop over, i.e. the ones with contributions
2343 * larger than min_gaussian */
2344 count = get_single_atom_gaussians(xj, erg);
2349 /* Loop over the relevant slabs for this atom */
2350 for (int ic = 0; ic < count; ic++)
2352 int n = erg->gn_slabind[ic];
2354 /* Get the precomputed Gaussian for xj in slab n */
2355 gaussian_xj = erg->gn_atom[ic];
2357 int slabIndex = n - erg->slab_first; /* slab index */
2359 /* The (unrotated) reference position of this atom is saved in yj0: */
2360 copy_rvec(erg->rotg->x_ref[iigrp], yj0);
2362 beta = calc_beta(xj, erg, n);
2364 /* The current center of this slab is saved in xcn: */
2365 copy_rvec(erg->slab_center[slabIndex], xcn);
2366 /* ... and the reference center in ycn: */
2367 copy_rvec(erg->slab_center_ref[slabIndex + erg->slab_buffer], ycn);
2369 rvec_sub(yj0, ycn, yj0_ycn); /* yj0_ycn = yj0 - ycn */
2371 /* In rare cases, when an atom position coincides with a reference slab
2372 * center (yj0_ycn == 0) we cannot compute the normal vector qjn.
2373 * However, since the atom is located directly on the pivot, this
2374 * slab's contribution to the force on that atom will be zero
2375 * anyway. Therefore, we directly move on to the next slab. */
2376 if (gmx_numzero(norm(yj0_ycn))) /* 0 == norm(yj0 - ycn) */
2382 mvmul(erg->rotmat, yj0_ycn, tmpvec2); /* tmpvec2= Omega.(yj0-ycn) */
2384 /* Subtract the slab center from xj */
2385 rvec_sub(xj, xcn, xj_xcn); /* xj_xcn = xj - xcn */
2388 cprod(erg->vec, tmpvec2, tmpvec); /* tmpvec= v x Omega.(yj0-ycn) */
2390 /* * v x Omega.(yj0-ycn) */
2391 unitv(tmpvec, qjn); /* qjn = --------------------- */
2392 /* |v x Omega.(yj0-ycn)| */
2394 bjn = iprod(qjn, xj_xcn); /* bjn = qjn * (xj - xcn) */
2396 /*********************************/
2397 /* Add to the rotation potential */
2398 /*********************************/
2399 V += 0.5 * erg->rotg->k * wj * gaussian_xj * gmx::square(bjn);
2401 /* If requested, also calculate the potential for a set of angles
2402 * near the current reference angle */
2405 for (int ifit = 0; ifit < erg->rotg->PotAngle_nstep; ifit++)
2407 /* As above calculate Omega.(yj0-ycn), now for the other angles */
2408 mvmul(erg->PotAngleFit->rotmat[ifit], yj0_ycn, tmpvec2); /* tmpvec2= Omega.(yj0-ycn) */
2409 /* As above calculate qjn */
2410 cprod(erg->vec, tmpvec2, tmpvec); /* tmpvec= v x Omega.(yj0-ycn) */
2411 /* * v x Omega.(yj0-ycn) */
2412 unitv(tmpvec, fit_qjn); /* fit_qjn = --------------------- */
2413 /* |v x Omega.(yj0-ycn)| */
2414 fit_bjn = iprod(fit_qjn, xj_xcn); /* fit_bjn = fit_qjn * (xj - xcn) */
2415 /* Add to the rotation potential for this angle */
2416 erg->PotAngleFit->V[ifit] +=
2417 0.5 * erg->rotg->k * wj * gaussian_xj * gmx::square(fit_bjn);
2421 /****************************************************************/
2422 /* sum_n1 will typically be the main contribution to the force: */
2423 /****************************************************************/
2424 betan_xj_sigma2 = beta * OOsigma2; /* beta_n(xj)/sigma^2 */
2426 /* The next lines calculate
2427 * qjn - (bjn*beta(xj)/(2sigma^2))v */
2428 svmul(bjn * 0.5 * betan_xj_sigma2, erg->vec, tmpvec2);
2429 rvec_sub(qjn, tmpvec2, tmpvec);
2431 /* Multiply with gn(xj)*bjn: */
2432 svmul(gaussian_xj * bjn, tmpvec, tmpvec2);
2435 rvec_inc(sum_n1, tmpvec2);
2437 /* We already have precalculated the Sn term for slab n */
2438 copy_rvec(erg->slab_innersumvec[slabIndex], s_n);
2440 svmul(betan_xj_sigma2 * iprod(s_n, xj_xcn), erg->vec, tmpvec); /* tmpvec = ---------- s_n (xj-xcn) */
2443 rvec_sub(s_n, tmpvec, innersumvec);
2445 /* We can safely divide by slab_weights since we check in get_slab_centers
2446 * that it is non-zero. */
2447 svmul(gaussian_xj / erg->slab_weights[slabIndex], innersumvec, innersumvec);
2449 rvec_add(sum_n2, innersumvec, sum_n2);
2451 /* Calculate the torque: */
2454 /* The force on atom ii from slab n only: */
2455 svmul(-erg->rotg->k * wj, tmpvec2, force_n1); /* part 1 */
2456 svmul(erg->rotg->k * mj, innersumvec, force_n2); /* part 2 */
2457 rvec_add(force_n1, force_n2, force_n);
2458 erg->slab_torque_v[slabIndex] += torque(erg->vec, force_n, xj, xcn);
2460 } /* END of loop over slabs */
2462 /* Put both contributions together: */
2463 svmul(wj, sum_n1, sum_n1);
2464 svmul(mj, sum_n2, sum_n2);
2465 rvec_sub(sum_n2, sum_n1, tmp_f); /* F = -grad V */
2467 /* Store the additional force so that it can be added to the force
2468 * array after the normal forces have been evaluated */
2469 for (int m = 0; m < DIM; m++)
2471 erg->f_rot_loc[j][m] = erg->rotg->k * tmp_f[m];
2476 } /* END of loop over local atoms */
2481 static void sort_collective_coordinates(gmx_enfrotgrp* erg,
2482 sort_along_vec_t* data) /* Buffer for sorting the positions */
2484 /* The projection of the position vector on the rotation vector is
2485 * the relevant value for sorting. Fill the 'data' structure */
2486 for (int i = 0; i < erg->rotg->nat; i++)
2488 data[i].xcproj = iprod(erg->xc[i], erg->vec); /* sort criterium */
2489 data[i].m = erg->mc[i];
2491 copy_rvec(erg->xc[i], data[i].x);
2492 copy_rvec(erg->rotg->x_ref[i], data[i].x_ref);
2494 /* Sort the 'data' structure */
2495 std::sort(data, data + erg->rotg->nat, [](const sort_along_vec_t& a, const sort_along_vec_t& b) {
2496 return a.xcproj < b.xcproj;
2499 /* Copy back the sorted values */
2500 for (int i = 0; i < erg->rotg->nat; i++)
2502 copy_rvec(data[i].x, erg->xc[i]);
2503 copy_rvec(data[i].x_ref, erg->xc_ref_sorted[i]);
2504 erg->mc_sorted[i] = data[i].m;
2505 erg->xc_sortind[i] = data[i].ind;
2510 /* For each slab, get the first and the last index of the sorted atom
2512 static void get_firstlast_atom_per_slab(const gmx_enfrotgrp* erg)
2516 /* Find the first atom that needs to enter the calculation for each slab */
2517 int n = erg->slab_first; /* slab */
2518 int i = 0; /* start with the first atom */
2521 /* Find the first atom that significantly contributes to this slab */
2522 do /* move forward in position until a large enough beta is found */
2524 beta = calc_beta(erg->xc[i], erg, n);
2526 } while ((beta < -erg->max_beta) && (i < erg->rotg->nat));
2528 int slabIndex = n - erg->slab_first; /* slab index */
2529 erg->firstatom[slabIndex] = i;
2530 /* Proceed to the next slab */
2532 } while (n <= erg->slab_last);
2534 /* Find the last atom for each slab */
2535 n = erg->slab_last; /* start with last slab */
2536 i = erg->rotg->nat - 1; /* start with the last atom */
2539 do /* move backward in position until a large enough beta is found */
2541 beta = calc_beta(erg->xc[i], erg, n);
2543 } while ((beta > erg->max_beta) && (i > -1));
2545 int slabIndex = n - erg->slab_first; /* slab index */
2546 erg->lastatom[slabIndex] = i;
2547 /* Proceed to the next slab */
2549 } while (n >= erg->slab_first);
2553 /* Determine the very first and very last slab that needs to be considered
2554 * For the first slab that needs to be considered, we have to find the smallest
2557 * x_first * v - n*Delta_x <= beta_max
2559 * slab index n, slab distance Delta_x, rotation vector v. For the last slab we
2560 * have to find the largest n that obeys
2562 * x_last * v - n*Delta_x >= -beta_max
2565 static inline int get_first_slab(const gmx_enfrotgrp* erg,
2566 rvec firstatom) /* First atom after sorting along the rotation vector v */
2568 /* Find the first slab for the first atom */
2569 return static_cast<int>(ceil(
2570 static_cast<double>((iprod(firstatom, erg->vec) - erg->max_beta) / erg->rotg->slab_dist)));
2574 static inline int get_last_slab(const gmx_enfrotgrp* erg, rvec lastatom) /* Last atom along v */
2576 /* Find the last slab for the last atom */
2577 return static_cast<int>(floor(
2578 static_cast<double>((iprod(lastatom, erg->vec) + erg->max_beta) / erg->rotg->slab_dist)));
2582 static void get_firstlast_slab_check(gmx_enfrotgrp* erg, /* The rotation group (data only accessible in this file) */
2583 rvec firstatom, /* First atom after sorting along the rotation vector v */
2584 rvec lastatom) /* Last atom along v */
2586 erg->slab_first = get_first_slab(erg, firstatom);
2587 erg->slab_last = get_last_slab(erg, lastatom);
2589 /* Calculate the slab buffer size, which changes when slab_first changes */
2590 erg->slab_buffer = erg->slab_first - erg->slab_first_ref;
2592 /* Check whether we have reference data to compare against */
2593 if (erg->slab_first < erg->slab_first_ref)
2595 gmx_fatal(FARGS, "%s No reference data for first slab (n=%d), unable to proceed.", RotStr, erg->slab_first);
2598 /* Check whether we have reference data to compare against */
2599 if (erg->slab_last > erg->slab_last_ref)
2601 gmx_fatal(FARGS, "%s No reference data for last slab (n=%d), unable to proceed.", RotStr, erg->slab_last);
2606 /* Enforced rotation with a flexible axis */
2607 static void do_flexible(gmx_bool bMaster,
2608 gmx_enfrot* enfrot, /* Other rotation data */
2610 rvec x[], /* The local positions */
2612 double t, /* Time in picoseconds */
2613 gmx_bool bOutstepRot, /* Output to main rotation output file */
2614 gmx_bool bOutstepSlab) /* Output per-slab data */
2617 real sigma; /* The Gaussian width sigma */
2619 /* Define the sigma value */
2620 sigma = 0.7 * erg->rotg->slab_dist;
2622 /* Sort the collective coordinates erg->xc along the rotation vector. This is
2623 * an optimization for the inner loop. */
2624 sort_collective_coordinates(erg, enfrot->data);
2626 /* Determine the first relevant slab for the first atom and the last
2627 * relevant slab for the last atom */
2628 get_firstlast_slab_check(erg, erg->xc[0], erg->xc[erg->rotg->nat - 1]);
2630 /* Determine for each slab depending on the min_gaussian cutoff criterium,
2631 * a first and a last atom index inbetween stuff needs to be calculated */
2632 get_firstlast_atom_per_slab(erg);
2634 /* Determine the gaussian-weighted center of positions for all slabs */
2635 get_slab_centers(erg, erg->xc, erg->mc_sorted, t, enfrot->out_slabs, bOutstepSlab, FALSE);
2637 /* Clear the torque per slab from last time step: */
2638 nslabs = erg->slab_last - erg->slab_first + 1;
2639 for (int l = 0; l < nslabs; l++)
2641 erg->slab_torque_v[l] = 0.0;
2644 /* Call the rotational forces kernel */
2645 if (erg->rotg->eType == EnforcedRotationGroupType::Flex
2646 || erg->rotg->eType == EnforcedRotationGroupType::Flext)
2648 erg->V = do_flex_lowlevel(erg, sigma, x, bOutstepRot, bOutstepSlab, box);
2650 else if (erg->rotg->eType == EnforcedRotationGroupType::Flex2
2651 || erg->rotg->eType == EnforcedRotationGroupType::Flex2t)
2653 erg->V = do_flex2_lowlevel(erg, sigma, x, bOutstepRot, bOutstepSlab, box);
2657 gmx_fatal(FARGS, "Unknown flexible rotation type");
2660 /* Determine angle by RMSD fit to the reference - Let's hope this */
2661 /* only happens once in a while, since this is not parallelized! */
2662 if (bMaster && (RotationGroupFitting::Pot != erg->rotg->eFittype))
2666 /* Fit angle of the whole rotation group */
2667 erg->angle_v = flex_fit_angle(erg);
2671 /* Fit angle of each slab */
2672 flex_fit_angle_perslab(erg, t, erg->degangle, enfrot->out_angles);
2676 /* Lump together the torques from all slabs: */
2677 erg->torque_v = 0.0;
2678 for (int l = 0; l < nslabs; l++)
2680 erg->torque_v += erg->slab_torque_v[l];
2685 /* Calculate the angle between reference and actual rotation group atom,
2686 * both projected into a plane perpendicular to the rotation vector: */
2687 static void angle(const gmx_enfrotgrp* erg,
2691 real* weight) /* atoms near the rotation axis should count less than atoms far away */
2693 rvec xp, xrp; /* current and reference positions projected on a plane perpendicular to pg->vec */
2697 /* Project x_ref and x into a plane through the origin perpendicular to rot_vec: */
2698 /* Project x_ref: xrp = x_ref - (vec * x_ref) * vec */
2699 svmul(iprod(erg->vec, x_ref), erg->vec, dum);
2700 rvec_sub(x_ref, dum, xrp);
2701 /* Project x_act: */
2702 svmul(iprod(erg->vec, x_act), erg->vec, dum);
2703 rvec_sub(x_act, dum, xp);
2705 /* Retrieve information about which vector precedes. gmx_angle always
2706 * returns a positive angle. */
2707 cprod(xp, xrp, dum); /* if reference precedes, this is pointing into the same direction as vec */
2709 if (iprod(erg->vec, dum) >= 0)
2711 *alpha = -gmx_angle(xrp, xp);
2715 *alpha = +gmx_angle(xrp, xp);
2718 /* Also return the weight */
2723 /* Project first vector onto a plane perpendicular to the second vector
2725 * Note that v must be of unit length.
2727 static inline void project_onto_plane(rvec dr, const rvec v)
2732 svmul(iprod(dr, v), v, tmp); /* tmp = (dr.v)v */
2733 rvec_dec(dr, tmp); /* dr = dr - (dr.v)v */
2737 /* Fixed rotation: The rotation reference group rotates around the v axis. */
2738 /* The atoms of the actual rotation group are attached with imaginary */
2739 /* springs to the reference atoms. */
2740 static void do_fixed(gmx_enfrotgrp* erg,
2741 gmx_bool bOutstepRot, /* Output to main rotation output file */
2742 gmx_bool bOutstepSlab) /* Output per-slab data */
2745 rvec tmp_f; /* Force */
2746 real alpha; /* a single angle between an actual and a reference position */
2747 real weight; /* single weight for a single angle */
2748 rvec xi_xc; /* xi - xc */
2749 gmx_bool bCalcPotFit;
2752 /* for mass weighting: */
2753 real wi; /* Mass-weighting of the positions */
2755 real k_wi; /* k times wi */
2759 bProject = (erg->rotg->eType == EnforcedRotationGroupType::Pm)
2760 || (erg->rotg->eType == EnforcedRotationGroupType::Pmpf);
2761 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (RotationGroupFitting::Pot == erg->rotg->eFittype);
2763 N_M = erg->rotg->nat * erg->invmass;
2764 const auto& collectiveRotationGroupIndex = erg->atomSet->collectiveIndex();
2765 /* Each process calculates the forces on its local atoms */
2766 for (size_t j = 0; j < erg->atomSet->numAtomsLocal(); j++)
2768 /* Calculate (x_i-x_c) resp. (x_i-u) */
2769 rvec_sub(erg->x_loc_pbc[j], erg->xc_center, xi_xc);
2771 /* Calculate Omega*(y_i-y_c)-(x_i-x_c) */
2772 rvec_sub(erg->xr_loc[j], xi_xc, dr);
2776 project_onto_plane(dr, erg->vec);
2779 /* Mass-weighting */
2780 wi = N_M * erg->m_loc[j];
2782 /* Store the additional force so that it can be added to the force
2783 * array after the normal forces have been evaluated */
2784 k_wi = erg->rotg->k * wi;
2785 for (int m = 0; m < DIM; m++)
2787 tmp_f[m] = k_wi * dr[m];
2788 erg->f_rot_loc[j][m] = tmp_f[m];
2789 erg->V += 0.5 * k_wi * gmx::square(dr[m]);
2792 /* If requested, also calculate the potential for a set of angles
2793 * near the current reference angle */
2796 for (int ifit = 0; ifit < erg->rotg->PotAngle_nstep; ifit++)
2798 /* Index of this rotation group atom with respect to the whole rotation group */
2799 int jj = collectiveRotationGroupIndex[j];
2801 /* Rotate with the alternative angle. Like rotate_local_reference(),
2802 * just for a single local atom */
2803 mvmul(erg->PotAngleFit->rotmat[ifit], erg->rotg->x_ref[jj], fit_xr_loc); /* fit_xr_loc = Omega*(y_i-y_c) */
2805 /* Calculate Omega*(y_i-y_c)-(x_i-x_c) */
2806 rvec_sub(fit_xr_loc, xi_xc, dr);
2810 project_onto_plane(dr, erg->vec);
2813 /* Add to the rotation potential for this angle: */
2814 erg->PotAngleFit->V[ifit] += 0.5 * k_wi * norm2(dr);
2820 /* Add to the torque of this rotation group */
2821 erg->torque_v += torque(erg->vec, tmp_f, erg->x_loc_pbc[j], erg->xc_center);
2823 /* Calculate the angle between reference and actual rotation group atom. */
2824 angle(erg, xi_xc, erg->xr_loc[j], &alpha, &weight); /* angle in rad, weighted */
2825 erg->angle_v += alpha * weight;
2826 erg->weight_v += weight;
2828 /* If you want enforced rotation to contribute to the virial,
2829 * activate the following lines:
2832 Add the rotation contribution to the virial
2833 for(j=0; j<DIM; j++)
2835 vir[j][m] += 0.5*f[ii][j]*dr[m];
2841 } /* end of loop over local rotation group atoms */
2845 /* Calculate the radial motion potential and forces */
2846 static void do_radial_motion(gmx_enfrotgrp* erg,
2847 gmx_bool bOutstepRot, /* Output to main rotation output file */
2848 gmx_bool bOutstepSlab) /* Output per-slab data */
2850 rvec tmp_f; /* Force */
2851 real alpha; /* a single angle between an actual and a reference position */
2852 real weight; /* single weight for a single angle */
2853 rvec xj_u; /* xj - u */
2854 rvec tmpvec, fit_tmpvec;
2855 real fac, fac2, sum = 0.0;
2857 gmx_bool bCalcPotFit;
2859 /* For mass weighting: */
2860 real wj; /* Mass-weighting of the positions */
2863 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (RotationGroupFitting::Pot == erg->rotg->eFittype);
2865 N_M = erg->rotg->nat * erg->invmass;
2866 const auto& collectiveRotationGroupIndex = erg->atomSet->collectiveIndex();
2867 /* Each process calculates the forces on its local atoms */
2868 for (size_t j = 0; j < erg->atomSet->numAtomsLocal(); j++)
2870 /* Calculate (xj-u) */
2871 rvec_sub(erg->x_loc_pbc[j], erg->xc_center, xj_u); /* xj_u = xj-u */
2873 /* Calculate Omega.(yj0-u) */
2874 cprod(erg->vec, erg->xr_loc[j], tmpvec); /* tmpvec = v x Omega.(yj0-u) */
2876 /* * v x Omega.(yj0-u) */
2877 unitv(tmpvec, pj); /* pj = --------------------- */
2878 /* | v x Omega.(yj0-u) | */
2880 fac = iprod(pj, xj_u); /* fac = pj.(xj-u) */
2883 /* Mass-weighting */
2884 wj = N_M * erg->m_loc[j];
2886 /* Store the additional force so that it can be added to the force
2887 * array after the normal forces have been evaluated */
2888 svmul(-erg->rotg->k * wj * fac, pj, tmp_f);
2889 copy_rvec(tmp_f, erg->f_rot_loc[j]);
2892 /* If requested, also calculate the potential for a set of angles
2893 * near the current reference angle */
2896 for (int ifit = 0; ifit < erg->rotg->PotAngle_nstep; ifit++)
2898 /* Index of this rotation group atom with respect to the whole rotation group */
2899 int jj = collectiveRotationGroupIndex[j];
2901 /* Rotate with the alternative angle. Like rotate_local_reference(),
2902 * just for a single local atom */
2903 mvmul(erg->PotAngleFit->rotmat[ifit], erg->rotg->x_ref[jj], fit_tmpvec); /* fit_tmpvec = Omega*(yj0-u) */
2905 /* Calculate Omega.(yj0-u) */
2906 cprod(erg->vec, fit_tmpvec, tmpvec); /* tmpvec = v x Omega.(yj0-u) */
2907 /* * v x Omega.(yj0-u) */
2908 unitv(tmpvec, pj); /* pj = --------------------- */
2909 /* | v x Omega.(yj0-u) | */
2911 fac = iprod(pj, xj_u); /* fac = pj.(xj-u) */
2914 /* Add to the rotation potential for this angle: */
2915 erg->PotAngleFit->V[ifit] += 0.5 * erg->rotg->k * wj * fac2;
2921 /* Add to the torque of this rotation group */
2922 erg->torque_v += torque(erg->vec, tmp_f, erg->x_loc_pbc[j], erg->xc_center);
2924 /* Calculate the angle between reference and actual rotation group atom. */
2925 angle(erg, xj_u, erg->xr_loc[j], &alpha, &weight); /* angle in rad, weighted */
2926 erg->angle_v += alpha * weight;
2927 erg->weight_v += weight;
2932 } /* end of loop over local rotation group atoms */
2933 erg->V = 0.5 * erg->rotg->k * sum;
2937 /* Calculate the radial motion pivot-free potential and forces */
2938 static void do_radial_motion_pf(gmx_enfrotgrp* erg,
2939 rvec x[], /* The positions */
2940 const matrix box, /* The simulation box */
2941 gmx_bool bOutstepRot, /* Output to main rotation output file */
2942 gmx_bool bOutstepSlab) /* Output per-slab data */
2944 rvec xj; /* Current position */
2945 rvec xj_xc; /* xj - xc */
2946 rvec yj0_yc0; /* yj0 - yc0 */
2947 rvec tmp_f; /* Force */
2948 real alpha; /* a single angle between an actual and a reference position */
2949 real weight; /* single weight for a single angle */
2950 rvec tmpvec, tmpvec2;
2951 rvec innersumvec; /* Precalculation of the inner sum */
2953 real fac, fac2, V = 0.0;
2955 gmx_bool bCalcPotFit;
2957 /* For mass weighting: */
2958 real mj, wi, wj; /* Mass-weighting of the positions */
2961 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (RotationGroupFitting::Pot == erg->rotg->eFittype);
2963 N_M = erg->rotg->nat * erg->invmass;
2965 /* Get the current center of the rotation group: */
2966 get_center(erg->xc, erg->mc, erg->rotg->nat, erg->xc_center);
2968 /* Precalculate Sum_i [ wi qi.(xi-xc) qi ] which is needed for every single j */
2969 clear_rvec(innersumvec);
2970 for (int i = 0; i < erg->rotg->nat; i++)
2972 /* Mass-weighting */
2973 wi = N_M * erg->mc[i];
2975 /* Calculate qi. Note that xc_ref_center has already been subtracted from
2976 * x_ref in init_rot_group.*/
2977 mvmul(erg->rotmat, erg->rotg->x_ref[i], tmpvec); /* tmpvec = Omega.(yi0-yc0) */
2979 cprod(erg->vec, tmpvec, tmpvec2); /* tmpvec2 = v x Omega.(yi0-yc0) */
2981 /* * v x Omega.(yi0-yc0) */
2982 unitv(tmpvec2, qi); /* qi = ----------------------- */
2983 /* | v x Omega.(yi0-yc0) | */
2985 rvec_sub(erg->xc[i], erg->xc_center, tmpvec); /* tmpvec = xi-xc */
2987 svmul(wi * iprod(qi, tmpvec), qi, tmpvec2);
2989 rvec_inc(innersumvec, tmpvec2);
2991 svmul(erg->rotg->k * erg->invmass, innersumvec, innersumveckM);
2993 /* Each process calculates the forces on its local atoms */
2994 const auto& localRotationGroupIndex = erg->atomSet->localIndex();
2995 const auto& collectiveRotationGroupIndex = erg->atomSet->collectiveIndex();
2996 for (gmx::index j = 0; j < localRotationGroupIndex.ssize(); j++)
2998 /* Local index of a rotation group atom */
2999 int ii = localRotationGroupIndex[j];
3000 /* Position of this atom in the collective array */
3001 int iigrp = collectiveRotationGroupIndex[j];
3002 /* Mass-weighting */
3003 mj = erg->mc[iigrp]; /* need the unsorted mass here */
3006 /* Current position of this atom: x[ii][XX/YY/ZZ] */
3007 copy_rvec(x[ii], xj);
3009 /* Shift this atom such that it is near its reference */
3010 shift_single_coord(box, xj, erg->xc_shifts[iigrp]);
3012 /* The (unrotated) reference position is yj0. yc0 has already
3013 * been subtracted in init_rot_group */
3014 copy_rvec(erg->rotg->x_ref[iigrp], yj0_yc0); /* yj0_yc0 = yj0 - yc0 */
3016 /* Calculate Omega.(yj0-yc0) */
3017 mvmul(erg->rotmat, yj0_yc0, tmpvec2); /* tmpvec2 = Omega.(yj0 - yc0) */
3019 cprod(erg->vec, tmpvec2, tmpvec); /* tmpvec = v x Omega.(yj0-yc0) */
3021 /* * v x Omega.(yj0-yc0) */
3022 unitv(tmpvec, qj); /* qj = ----------------------- */
3023 /* | v x Omega.(yj0-yc0) | */
3025 /* Calculate (xj-xc) */
3026 rvec_sub(xj, erg->xc_center, xj_xc); /* xj_xc = xj-xc */
3028 fac = iprod(qj, xj_xc); /* fac = qj.(xj-xc) */
3031 /* Store the additional force so that it can be added to the force
3032 * array after the normal forces have been evaluated */
3033 svmul(-erg->rotg->k * wj * fac, qj, tmp_f); /* part 1 of force */
3034 svmul(mj, innersumveckM, tmpvec); /* part 2 of force */
3035 rvec_inc(tmp_f, tmpvec);
3036 copy_rvec(tmp_f, erg->f_rot_loc[j]);
3039 /* If requested, also calculate the potential for a set of angles
3040 * near the current reference angle */
3043 for (int ifit = 0; ifit < erg->rotg->PotAngle_nstep; ifit++)
3045 /* Rotate with the alternative angle. Like rotate_local_reference(),
3046 * just for a single local atom */
3047 mvmul(erg->PotAngleFit->rotmat[ifit], yj0_yc0, tmpvec2); /* tmpvec2 = Omega*(yj0-yc0) */
3049 /* Calculate Omega.(yj0-u) */
3050 cprod(erg->vec, tmpvec2, tmpvec); /* tmpvec = v x Omega.(yj0-yc0) */
3051 /* * v x Omega.(yj0-yc0) */
3052 unitv(tmpvec, qj); /* qj = ----------------------- */
3053 /* | v x Omega.(yj0-yc0) | */
3055 fac = iprod(qj, xj_xc); /* fac = qj.(xj-xc) */
3058 /* Add to the rotation potential for this angle: */
3059 erg->PotAngleFit->V[ifit] += 0.5 * erg->rotg->k * wj * fac2;
3065 /* Add to the torque of this rotation group */
3066 erg->torque_v += torque(erg->vec, tmp_f, xj, erg->xc_center);
3068 /* Calculate the angle between reference and actual rotation group atom. */
3069 angle(erg, xj_xc, yj0_yc0, &alpha, &weight); /* angle in rad, weighted */
3070 erg->angle_v += alpha * weight;
3071 erg->weight_v += weight;
3076 } /* end of loop over local rotation group atoms */
3077 erg->V = 0.5 * erg->rotg->k * V;
3081 /* Precalculate the inner sum for the radial motion 2 forces */
3082 static void radial_motion2_precalc_inner_sum(const gmx_enfrotgrp* erg, rvec innersumvec)
3085 rvec xi_xc; /* xj - xc */
3086 rvec tmpvec, tmpvec2;
3090 rvec v_xi_xc; /* v x (xj - u) */
3091 real psii, psiistar;
3092 real wi; /* Mass-weighting of the positions */
3096 N_M = erg->rotg->nat * erg->invmass;
3098 /* Loop over the collective set of positions */
3100 for (i = 0; i < erg->rotg->nat; i++)
3102 /* Mass-weighting */
3103 wi = N_M * erg->mc[i];
3105 rvec_sub(erg->xc[i], erg->xc_center, xi_xc); /* xi_xc = xi-xc */
3107 /* Calculate ri. Note that xc_ref_center has already been subtracted from
3108 * x_ref in init_rot_group.*/
3109 mvmul(erg->rotmat, erg->rotg->x_ref[i], ri); /* ri = Omega.(yi0-yc0) */
3111 cprod(erg->vec, xi_xc, v_xi_xc); /* v_xi_xc = v x (xi-u) */
3113 fac = norm2(v_xi_xc);
3115 psiistar = 1.0 / (fac + erg->rotg->eps); /* psiistar = --------------------- */
3116 /* |v x (xi-xc)|^2 + eps */
3118 psii = gmx::invsqrt(fac); /* 1 */
3119 /* psii = ------------- */
3122 svmul(psii, v_xi_xc, si); /* si = psii * (v x (xi-xc) ) */
3124 siri = iprod(si, ri); /* siri = si.ri */
3126 svmul(psiistar / psii, ri, tmpvec);
3127 svmul(psiistar * psiistar / (psii * psii * psii) * siri, si, tmpvec2);
3128 rvec_dec(tmpvec, tmpvec2);
3129 cprod(tmpvec, erg->vec, tmpvec2);
3131 svmul(wi * siri, tmpvec2, tmpvec);
3133 rvec_inc(sumvec, tmpvec);
3135 svmul(erg->rotg->k * erg->invmass, sumvec, innersumvec);
3139 /* Calculate the radial motion 2 potential and forces */
3140 static void do_radial_motion2(gmx_enfrotgrp* erg,
3141 rvec x[], /* The positions */
3142 const matrix box, /* The simulation box */
3143 gmx_bool bOutstepRot, /* Output to main rotation output file */
3144 gmx_bool bOutstepSlab) /* Output per-slab data */
3146 rvec xj; /* Position */
3147 real alpha; /* a single angle between an actual and a reference position */
3148 real weight; /* single weight for a single angle */
3149 rvec xj_u; /* xj - u */
3150 rvec yj0_yc0; /* yj0 -yc0 */
3151 rvec tmpvec, tmpvec2;
3152 real fac, fit_fac, fac2, Vpart = 0.0;
3153 rvec rj, fit_rj, sj;
3155 rvec v_xj_u; /* v x (xj - u) */
3156 real psij, psijstar;
3157 real mj, wj; /* For mass-weighting of the positions */
3161 gmx_bool bCalcPotFit;
3163 bPF = erg->rotg->eType == EnforcedRotationGroupType::Rm2pf;
3164 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (RotationGroupFitting::Pot == erg->rotg->eFittype);
3166 clear_rvec(yj0_yc0); /* Make the compiler happy */
3168 clear_rvec(innersumvec);
3171 /* For the pivot-free variant we have to use the current center of
3172 * mass of the rotation group instead of the pivot u */
3173 get_center(erg->xc, erg->mc, erg->rotg->nat, erg->xc_center);
3175 /* Also, we precalculate the second term of the forces that is identical
3176 * (up to the weight factor mj) for all forces */
3177 radial_motion2_precalc_inner_sum(erg, innersumvec);
3180 N_M = erg->rotg->nat * erg->invmass;
3182 /* Each process calculates the forces on its local atoms */
3183 const auto& localRotationGroupIndex = erg->atomSet->localIndex();
3184 const auto& collectiveRotationGroupIndex = erg->atomSet->collectiveIndex();
3185 for (gmx::index j = 0; j < localRotationGroupIndex.ssize(); j++)
3189 /* Local index of a rotation group atom */
3190 int ii = localRotationGroupIndex[j];
3191 /* Position of this atom in the collective array */
3192 int iigrp = collectiveRotationGroupIndex[j];
3193 /* Mass-weighting */
3194 mj = erg->mc[iigrp];
3196 /* Current position of this atom: x[ii] */
3197 copy_rvec(x[ii], xj);
3199 /* Shift this atom such that it is near its reference */
3200 shift_single_coord(box, xj, erg->xc_shifts[iigrp]);
3202 /* The (unrotated) reference position is yj0. yc0 has already
3203 * been subtracted in init_rot_group */
3204 copy_rvec(erg->rotg->x_ref[iigrp], yj0_yc0); /* yj0_yc0 = yj0 - yc0 */
3206 /* Calculate Omega.(yj0-yc0) */
3207 mvmul(erg->rotmat, yj0_yc0, rj); /* rj = Omega.(yj0-yc0) */
3212 copy_rvec(erg->x_loc_pbc[j], xj);
3213 copy_rvec(erg->xr_loc[j], rj); /* rj = Omega.(yj0-u) */
3215 /* Mass-weighting */
3218 /* Calculate (xj-u) resp. (xj-xc) */
3219 rvec_sub(xj, erg->xc_center, xj_u); /* xj_u = xj-u */
3221 cprod(erg->vec, xj_u, v_xj_u); /* v_xj_u = v x (xj-u) */
3223 fac = norm2(v_xj_u);
3225 psijstar = 1.0 / (fac + erg->rotg->eps); /* psistar = -------------------- */
3226 /* * |v x (xj-u)|^2 + eps */
3228 psij = gmx::invsqrt(fac); /* 1 */
3229 /* psij = ------------ */
3232 svmul(psij, v_xj_u, sj); /* sj = psij * (v x (xj-u) ) */
3234 fac = iprod(v_xj_u, rj); /* fac = (v x (xj-u)).rj */
3237 sjrj = iprod(sj, rj); /* sjrj = sj.rj */
3239 svmul(psijstar / psij, rj, tmpvec);
3240 svmul(psijstar * psijstar / (psij * psij * psij) * sjrj, sj, tmpvec2);
3241 rvec_dec(tmpvec, tmpvec2);
3242 cprod(tmpvec, erg->vec, tmpvec2);
3244 /* Store the additional force so that it can be added to the force
3245 * array after the normal forces have been evaluated */
3246 svmul(-erg->rotg->k * wj * sjrj, tmpvec2, tmpvec);
3247 svmul(mj, innersumvec, tmpvec2); /* This is != 0 only for the pivot-free variant */
3249 rvec_add(tmpvec2, tmpvec, erg->f_rot_loc[j]);
3250 Vpart += wj * psijstar * fac2;
3252 /* If requested, also calculate the potential for a set of angles
3253 * near the current reference angle */
3256 for (int ifit = 0; ifit < erg->rotg->PotAngle_nstep; ifit++)
3260 mvmul(erg->PotAngleFit->rotmat[ifit], yj0_yc0, fit_rj); /* fit_rj = Omega.(yj0-yc0) */
3264 /* Position of this atom in the collective array */
3265 int iigrp = collectiveRotationGroupIndex[j];
3266 /* Rotate with the alternative angle. Like rotate_local_reference(),
3267 * just for a single local atom */
3268 mvmul(erg->PotAngleFit->rotmat[ifit], erg->rotg->x_ref[iigrp], fit_rj); /* fit_rj = Omega*(yj0-u) */
3270 fit_fac = iprod(v_xj_u, fit_rj); /* fac = (v x (xj-u)).fit_rj */
3271 /* Add to the rotation potential for this angle: */
3272 erg->PotAngleFit->V[ifit] += 0.5 * erg->rotg->k * wj * psijstar * fit_fac * fit_fac;
3278 /* Add to the torque of this rotation group */
3279 erg->torque_v += torque(erg->vec, erg->f_rot_loc[j], xj, erg->xc_center);
3281 /* Calculate the angle between reference and actual rotation group atom. */
3282 angle(erg, xj_u, rj, &alpha, &weight); /* angle in rad, weighted */
3283 erg->angle_v += alpha * weight;
3284 erg->weight_v += weight;
3289 } /* end of loop over local rotation group atoms */
3290 erg->V = 0.5 * erg->rotg->k * Vpart;
3294 /* Determine the smallest and largest position vector (with respect to the
3295 * rotation vector) for the reference group */
3296 static void get_firstlast_atom_ref(const gmx_enfrotgrp* erg, int* firstindex, int* lastindex)
3299 real xcproj; /* The projection of a reference position on the
3301 real minproj, maxproj; /* Smallest and largest projection on v */
3303 /* Start with some value */
3304 minproj = iprod(erg->rotg->x_ref[0], erg->vec);
3307 /* This is just to ensure that it still works if all the atoms of the
3308 * reference structure are situated in a plane perpendicular to the rotation
3311 *lastindex = erg->rotg->nat - 1;
3313 /* Loop over all atoms of the reference group,
3314 * project them on the rotation vector to find the extremes */
3315 for (i = 0; i < erg->rotg->nat; i++)
3317 xcproj = iprod(erg->rotg->x_ref[i], erg->vec);
3318 if (xcproj < minproj)
3323 if (xcproj > maxproj)
3332 /* Allocate memory for the slabs */
3333 static void allocate_slabs(gmx_enfrotgrp* erg, FILE* fplog, gmx_bool bVerbose)
3335 /* More slabs than are defined for the reference are never needed */
3336 int nslabs = erg->slab_last_ref - erg->slab_first_ref + 1;
3338 /* Remember how many we allocated */
3339 erg->nslabs_alloc = nslabs;
3341 if ((nullptr != fplog) && bVerbose)
3344 "%s allocating memory to store data for %d slabs (rotation group %d).\n",
3349 snew(erg->slab_center, nslabs);
3350 snew(erg->slab_center_ref, nslabs);
3351 snew(erg->slab_weights, nslabs);
3352 snew(erg->slab_torque_v, nslabs);
3353 snew(erg->slab_data, nslabs);
3354 snew(erg->gn_atom, nslabs);
3355 snew(erg->gn_slabind, nslabs);
3356 snew(erg->slab_innersumvec, nslabs);
3357 for (int i = 0; i < nslabs; i++)
3359 snew(erg->slab_data[i].x, erg->rotg->nat);
3360 snew(erg->slab_data[i].ref, erg->rotg->nat);
3361 snew(erg->slab_data[i].weight, erg->rotg->nat);
3363 snew(erg->xc_ref_sorted, erg->rotg->nat);
3364 snew(erg->xc_sortind, erg->rotg->nat);
3365 snew(erg->firstatom, nslabs);
3366 snew(erg->lastatom, nslabs);
3370 /* From the extreme positions of the reference group, determine the first
3371 * and last slab of the reference. We can never have more slabs in the real
3372 * simulation than calculated here for the reference.
3374 static void get_firstlast_slab_ref(gmx_enfrotgrp* erg, real mc[], int ref_firstindex, int ref_lastindex)
3378 int first = get_first_slab(erg, erg->rotg->x_ref[ref_firstindex]);
3379 int last = get_last_slab(erg, erg->rotg->x_ref[ref_lastindex]);
3381 while (get_slab_weight(first, erg, erg->rotg->x_ref, mc, &dummy) > WEIGHT_MIN)
3385 erg->slab_first_ref = first + 1;
3386 while (get_slab_weight(last, erg, erg->rotg->x_ref, mc, &dummy) > WEIGHT_MIN)
3390 erg->slab_last_ref = last - 1;
3394 /* Special version of copy_rvec:
3395 * During the copy procedure of xcurr to b, the correct PBC image is chosen
3396 * such that the copied vector ends up near its reference position xref */
3397 static inline void copy_correct_pbc_image(const rvec xcurr, /* copy vector xcurr ... */
3398 rvec b, /* ... to b ... */
3399 const rvec xref, /* choosing the PBC image such that b ends up near xref */
3408 /* Shortest PBC distance between the atom and its reference */
3409 rvec_sub(xcurr, xref, dx);
3411 /* Determine the shift for this atom */
3413 for (m = npbcdim - 1; m >= 0; m--)
3415 while (dx[m] < -0.5 * box[m][m])
3417 for (d = 0; d < DIM; d++)
3423 while (dx[m] >= 0.5 * box[m][m])
3425 for (d = 0; d < DIM; d++)
3433 /* Apply the shift to the position */
3434 copy_rvec(xcurr, b);
3435 shift_single_coord(box, b, shift);
3439 static void init_rot_group(FILE* fplog,
3440 const t_commrec* cr,
3448 gmx_bool bOutputCenters)
3450 rvec coord, xref, *xdum;
3451 gmx_bool bFlex, bColl;
3452 int ref_firstindex, ref_lastindex;
3453 real mass, totalmass;
3456 const t_rotgrp* rotg = erg->rotg;
3459 /* Do we have a flexible axis? */
3460 bFlex = ISFLEX(rotg);
3461 /* Do we use a global set of coordinates? */
3462 bColl = ISCOLL(rotg);
3464 /* Allocate space for collective coordinates if needed */
3467 snew(erg->xc, erg->rotg->nat);
3468 snew(erg->xc_shifts, erg->rotg->nat);
3469 snew(erg->xc_eshifts, erg->rotg->nat);
3470 snew(erg->xc_old, erg->rotg->nat);
3472 if (erg->rotg->eFittype == RotationGroupFitting::Norm)
3474 snew(erg->xc_ref_length, erg->rotg->nat); /* in case fit type NORM is chosen */
3475 snew(erg->xc_norm, erg->rotg->nat);
3480 snew(erg->xr_loc, erg->rotg->nat);
3481 snew(erg->x_loc_pbc, erg->rotg->nat);
3484 copy_rvec(erg->rotg->inputVec, erg->vec);
3485 snew(erg->f_rot_loc, erg->rotg->nat);
3487 /* Make space for the calculation of the potential at other angles (used
3488 * for fitting only) */
3489 if (RotationGroupFitting::Pot == erg->rotg->eFittype)
3491 snew(erg->PotAngleFit, 1);
3492 snew(erg->PotAngleFit->degangle, erg->rotg->PotAngle_nstep);
3493 snew(erg->PotAngleFit->V, erg->rotg->PotAngle_nstep);
3494 snew(erg->PotAngleFit->rotmat, erg->rotg->PotAngle_nstep);
3496 /* Get the set of angles around the reference angle */
3497 start = -0.5 * (erg->rotg->PotAngle_nstep - 1) * erg->rotg->PotAngle_step;
3498 for (int i = 0; i < erg->rotg->PotAngle_nstep; i++)
3500 erg->PotAngleFit->degangle[i] = start + i * erg->rotg->PotAngle_step;
3505 erg->PotAngleFit = nullptr;
3508 /* Copy the masses so that the center can be determined. For all types of
3509 * enforced rotation, we store the masses in the erg->mc array. */
3510 snew(erg->mc, erg->rotg->nat);
3513 snew(erg->mc_sorted, erg->rotg->nat);
3517 snew(erg->m_loc, erg->rotg->nat);
3521 for (int i = 0; i < erg->rotg->nat; i++)
3523 if (erg->rotg->bMassW)
3525 mass = mtopGetAtomMass(mtop, erg->rotg->ind[i], &molb);
3534 erg->invmass = 1.0 / totalmass;
3536 /* Set xc_ref_center for any rotation potential */
3537 if ((erg->rotg->eType == EnforcedRotationGroupType::Iso)
3538 || (erg->rotg->eType == EnforcedRotationGroupType::Pm)
3539 || (erg->rotg->eType == EnforcedRotationGroupType::Rm)
3540 || (erg->rotg->eType == EnforcedRotationGroupType::Rm2))
3542 /* Set the pivot point for the fixed, stationary-axis potentials. This
3543 * won't change during the simulation */
3544 copy_rvec(erg->rotg->pivot, erg->xc_ref_center);
3545 copy_rvec(erg->rotg->pivot, erg->xc_center);
3549 /* Center of the reference positions */
3550 get_center(erg->rotg->x_ref, erg->mc, erg->rotg->nat, erg->xc_ref_center);
3552 /* Center of the actual positions */
3555 snew(xdum, erg->rotg->nat);
3556 for (int i = 0; i < erg->rotg->nat; i++)
3558 int ii = erg->rotg->ind[i];
3559 copy_rvec(x[ii], xdum[i]);
3561 get_center(xdum, erg->mc, erg->rotg->nat, erg->xc_center);
3567 gmx_bcast(sizeof(erg->xc_center), erg->xc_center, cr->mpi_comm_mygroup);
3574 /* Save the original (whole) set of positions in xc_old such that at later
3575 * steps the rotation group can always be made whole again. If the simulation is
3576 * restarted, we compute the starting reference positions (given the time)
3577 * and assume that the correct PBC image of each position is the one nearest
3578 * to the current reference */
3581 /* Calculate the rotation matrix for this angle: */
3582 t_start = ir->init_t + ir->init_step * ir->delta_t;
3583 erg->degangle = erg->rotg->rate * t_start;
3584 calc_rotmat(erg->vec, erg->degangle, erg->rotmat);
3586 for (int i = 0; i < erg->rotg->nat; i++)
3588 int ii = erg->rotg->ind[i];
3590 /* Subtract pivot, rotate, and add pivot again. This will yield the
3591 * reference position for time t */
3592 rvec_sub(erg->rotg->x_ref[i], erg->xc_ref_center, coord);
3593 mvmul(erg->rotmat, coord, xref);
3594 rvec_inc(xref, erg->xc_ref_center);
3596 copy_correct_pbc_image(x[ii], erg->xc_old[i], xref, box, 3);
3602 gmx_bcast(erg->rotg->nat * sizeof(erg->xc_old[0]), erg->xc_old, cr->mpi_comm_mygroup);
3607 if ((erg->rotg->eType != EnforcedRotationGroupType::Flex)
3608 && (erg->rotg->eType != EnforcedRotationGroupType::Flex2))
3610 /* Put the reference positions into origin: */
3611 for (int i = 0; i < erg->rotg->nat; i++)
3613 rvec_dec(erg->rotg->x_ref[i], erg->xc_ref_center);
3617 /* Enforced rotation with flexible axis */
3620 /* Calculate maximum beta value from minimum gaussian (performance opt.) */
3621 erg->max_beta = calc_beta_max(erg->rotg->min_gaussian, erg->rotg->slab_dist);
3623 /* Determine the smallest and largest coordinate with respect to the rotation vector */
3624 get_firstlast_atom_ref(erg, &ref_firstindex, &ref_lastindex);
3626 /* From the extreme positions of the reference group, determine the first
3627 * and last slab of the reference. */
3628 get_firstlast_slab_ref(erg, erg->mc, ref_firstindex, ref_lastindex);
3630 /* Allocate memory for the slabs */
3631 allocate_slabs(erg, fplog, bVerbose);
3633 /* Flexible rotation: determine the reference centers for the rest of the simulation */
3634 erg->slab_first = erg->slab_first_ref;
3635 erg->slab_last = erg->slab_last_ref;
3636 get_slab_centers(erg, erg->rotg->x_ref, erg->mc, -1, out_slabs, bOutputCenters, TRUE);
3638 /* Length of each x_rotref vector from center (needed if fit routine NORM is chosen): */
3639 if (erg->rotg->eFittype == RotationGroupFitting::Norm)
3641 for (int i = 0; i < erg->rotg->nat; i++)
3643 rvec_sub(erg->rotg->x_ref[i], erg->xc_ref_center, coord);
3644 erg->xc_ref_length[i] = norm(coord);
3650 /* Calculate the size of the MPI buffer needed in reduce_output() */
3651 static int calc_mpi_bufsize(const gmx_enfrot* er)
3654 int count_total = 0;
3655 for (int g = 0; g < er->rot->ngrp; g++)
3657 const t_rotgrp* rotg = &er->rot->grp[g];
3658 const gmx_enfrotgrp* erg = &er->enfrotgrp[g];
3660 /* Count the items that are transferred for this group: */
3661 int count_group = 4; /* V, torque, angle, weight */
3663 /* Add the maximum number of slabs for flexible groups */
3666 count_group += erg->slab_last_ref - erg->slab_first_ref + 1;
3669 /* Add space for the potentials at different angles: */
3670 if (RotationGroupFitting::Pot == erg->rotg->eFittype)
3672 count_group += erg->rotg->PotAngle_nstep;
3675 /* Add to the total number: */
3676 count_total += count_group;
3683 std::unique_ptr<gmx::EnforcedRotation> init_rot(FILE* fplog,
3686 const t_filenm fnm[],
3687 const t_commrec* cr,
3688 gmx::LocalAtomSetManager* atomSets,
3689 const t_state* globalState,
3691 const gmx_output_env_t* oenv,
3692 const gmx::MdrunOptions& mdrunOptions,
3693 const gmx::StartingBehavior startingBehavior)
3695 int nat_max = 0; /* Size of biggest rotation group */
3696 rvec* x_pbc = nullptr; /* Space for the pbc-correct atom positions */
3698 if (MASTER(cr) && mdrunOptions.verbose)
3700 fprintf(stdout, "%s Initializing ...\n", RotStr);
3703 auto enforcedRotation = std::make_unique<gmx::EnforcedRotation>();
3704 gmx_enfrot* er = enforcedRotation->getLegacyEnfrot();
3705 // TODO When this module implements IMdpOptions, the ownership will become more clear.
3707 er->restartWithAppending = (startingBehavior == gmx::StartingBehavior::RestartWithAppending);
3709 /* When appending, skip first output to avoid duplicate entries in the data files */
3710 er->bOut = er->restartWithAppending;
3712 if (MASTER(cr) && er->bOut)
3714 please_cite(fplog, "Kutzner2011");
3717 /* Output every step for reruns */
3718 if (mdrunOptions.rerun)
3720 if (nullptr != fplog)
3722 fprintf(fplog, "%s rerun - will write rotation output every available step.\n", RotStr);
3729 er->nstrout = er->rot->nstrout;
3730 er->nstsout = er->rot->nstsout;
3733 er->out_slabs = nullptr;
3734 if (MASTER(cr) && HaveFlexibleGroups(er->rot))
3736 er->out_slabs = open_slab_out(opt2fn("-rs", nfile, fnm), er);
3741 /* Remove pbc, make molecule whole.
3742 * When ir->bContinuation=TRUE this has already been done, but ok. */
3743 snew(x_pbc, mtop->natoms);
3744 copy_rvecn(globalState->x.rvec_array(), x_pbc, 0, mtop->natoms);
3745 do_pbc_first_mtop(nullptr, ir->pbcType, globalState->box, mtop, x_pbc);
3746 /* All molecules will be whole now, but not necessarily in the home box.
3747 * Additionally, if a rotation group consists of more than one molecule
3748 * (e.g. two strands of DNA), each one of them can end up in a different
3749 * periodic box. This is taken care of in init_rot_group. */
3752 /* Allocate space for the per-rotation-group data: */
3753 er->enfrotgrp.resize(er->rot->ngrp);
3755 for (auto& ergRef : er->enfrotgrp)
3757 gmx_enfrotgrp* erg = &ergRef;
3758 erg->rotg = &er->rot->grp[groupIndex];
3759 erg->atomSet = std::make_unique<gmx::LocalAtomSet>(
3760 atomSets->add({ erg->rotg->ind, erg->rotg->ind + erg->rotg->nat }));
3761 erg->groupIndex = groupIndex;
3763 if (nullptr != fplog)
3766 "%s group %d type '%s'\n",
3769 enumValueToString(erg->rotg->eType));
3772 if (erg->rotg->nat > 0)
3774 nat_max = std::max(nat_max, erg->rotg->nat);
3776 init_rot_group(fplog,
3781 mdrunOptions.verbose,
3783 MASTER(cr) ? globalState->box : nullptr,
3785 !er->restartWithAppending); /* Do not output the reference centers
3786 * again if we are appending */
3791 /* Allocate space for enforced rotation buffer variables */
3792 er->bufsize = nat_max;
3793 snew(er->data, nat_max);
3794 snew(er->xbuf, nat_max);
3795 snew(er->mbuf, nat_max);
3797 /* Buffers for MPI reducing torques, angles, weights (for each group), and V */
3800 er->mpi_bufsize = calc_mpi_bufsize(er) + 100; /* larger to catch errors */
3801 snew(er->mpi_inbuf, er->mpi_bufsize);
3802 snew(er->mpi_outbuf, er->mpi_bufsize);
3806 er->mpi_bufsize = 0;
3807 er->mpi_inbuf = nullptr;
3808 er->mpi_outbuf = nullptr;
3811 /* Only do I/O on the MASTER */
3812 er->out_angles = nullptr;
3813 er->out_rot = nullptr;
3814 er->out_torque = nullptr;
3817 er->out_rot = open_rot_out(opt2fn("-ro", nfile, fnm), oenv, er);
3819 if (er->nstsout > 0)
3821 if (HaveFlexibleGroups(er->rot) || HavePotFitGroups(er->rot))
3823 er->out_angles = open_angles_out(opt2fn("-ra", nfile, fnm), er);
3825 if (HaveFlexibleGroups(er->rot))
3827 er->out_torque = open_torque_out(opt2fn("-rt", nfile, fnm), er);
3833 return enforcedRotation;
3836 /* Rotate the local reference positions and store them in
3837 * erg->xr_loc[0...(nat_loc-1)]
3839 * Note that we already subtracted u or y_c from the reference positions
3840 * in init_rot_group().
3842 static void rotate_local_reference(gmx_enfrotgrp* erg)
3844 const auto& collectiveRotationGroupIndex = erg->atomSet->collectiveIndex();
3845 for (size_t i = 0; i < erg->atomSet->numAtomsLocal(); i++)
3847 /* Index of this rotation group atom with respect to the whole rotation group */
3848 int ii = collectiveRotationGroupIndex[i];
3850 mvmul(erg->rotmat, erg->rotg->x_ref[ii], erg->xr_loc[i]);
3855 /* Select the PBC representation for each local x position and store that
3856 * for later usage. We assume the right PBC image of an x is the one nearest to
3857 * its rotated reference */
3858 static void choose_pbc_image(rvec x[], gmx_enfrotgrp* erg, const matrix box, int npbcdim)
3860 const auto& localRotationGroupIndex = erg->atomSet->localIndex();
3861 for (gmx::index i = 0; i < localRotationGroupIndex.ssize(); i++)
3863 /* Index of a rotation group atom */
3864 int ii = localRotationGroupIndex[i];
3866 /* Get the correctly rotated reference position. The pivot was already
3867 * subtracted in init_rot_group() from the reference positions. Also,
3868 * the reference positions have already been rotated in
3869 * rotate_local_reference(). For the current reference position we thus
3870 * only need to add the pivot again. */
3872 copy_rvec(erg->xr_loc[i], xref);
3873 rvec_inc(xref, erg->xc_ref_center);
3875 copy_correct_pbc_image(x[ii], erg->x_loc_pbc[i], xref, box, npbcdim);
3880 void do_rotation(const t_commrec* cr, gmx_enfrot* er, const matrix box, rvec x[], real t, int64_t step, gmx_bool bNS)
3882 gmx_bool outstep_slab, outstep_rot;
3885 gmx_potfit* fit = nullptr; /* For fit type 'potential' determine the fit
3886 angle via the potential minimum */
3892 /* When to output in main rotation output file */
3893 outstep_rot = do_per_step(step, er->nstrout) && er->bOut;
3894 /* When to output per-slab data */
3895 outstep_slab = do_per_step(step, er->nstsout) && er->bOut;
3897 /* Output time into rotation output file */
3898 if (outstep_rot && MASTER(cr))
3900 fprintf(er->out_rot, "%12.3e", t);
3903 /**************************************************************************/
3904 /* First do ALL the communication! */
3905 for (auto& ergRef : er->enfrotgrp)
3907 gmx_enfrotgrp* erg = &ergRef;
3908 const t_rotgrp* rotg = erg->rotg;
3910 /* Do we use a collective (global) set of coordinates? */
3911 bColl = ISCOLL(rotg);
3913 /* Calculate the rotation matrix for this angle: */
3914 erg->degangle = rotg->rate * t;
3915 calc_rotmat(erg->vec, erg->degangle, erg->rotmat);
3919 /* Transfer the rotation group's positions such that every node has
3920 * all of them. Every node contributes its local positions x and stores
3921 * it in the collective erg->xc array. */
3922 communicate_group_positions(cr,
3929 erg->atomSet->numAtomsLocal(),
3930 erg->atomSet->localIndex().data(),
3931 erg->atomSet->collectiveIndex().data(),
3937 /* Fill the local masses array;
3938 * this array changes in DD/neighborsearching steps */
3941 const auto& collectiveRotationGroupIndex = erg->atomSet->collectiveIndex();
3942 for (gmx::index i = 0; i < collectiveRotationGroupIndex.ssize(); i++)
3944 /* Index of local atom w.r.t. the collective rotation group */
3945 int ii = collectiveRotationGroupIndex[i];
3946 erg->m_loc[i] = erg->mc[ii];
3950 /* Calculate Omega*(y_i-y_c) for the local positions */
3951 rotate_local_reference(erg);
3953 /* Choose the nearest PBC images of the group atoms with respect
3954 * to the rotated reference positions */
3955 choose_pbc_image(x, erg, box, 3);
3957 /* Get the center of the rotation group */
3958 if ((rotg->eType == EnforcedRotationGroupType::Isopf)
3959 || (rotg->eType == EnforcedRotationGroupType::Pmpf))
3962 cr, erg->x_loc_pbc, erg->m_loc, erg->atomSet->numAtomsLocal(), rotg->nat, erg->xc_center);
3966 } /* End of loop over rotation groups */
3968 /**************************************************************************/
3969 /* Done communicating, we can start to count cycles for the load balancing now ... */
3970 if (DOMAINDECOMP(cr))
3972 ddReopenBalanceRegionCpu(cr->dd);
3979 for (auto& ergRef : er->enfrotgrp)
3981 gmx_enfrotgrp* erg = &ergRef;
3982 const t_rotgrp* rotg = erg->rotg;
3984 if (outstep_rot && MASTER(cr))
3986 fprintf(er->out_rot, "%12.4f", erg->degangle);
3989 /* Calculate angles and rotation matrices for potential fitting: */
3990 if ((outstep_rot || outstep_slab) && (RotationGroupFitting::Pot == rotg->eFittype))
3992 fit = erg->PotAngleFit;
3993 for (int i = 0; i < rotg->PotAngle_nstep; i++)
3995 calc_rotmat(erg->vec, erg->degangle + fit->degangle[i], fit->rotmat[i]);
3997 /* Clear value from last step */
3998 erg->PotAngleFit->V[i] = 0.0;
4002 /* Clear values from last time step */
4004 erg->torque_v = 0.0;
4006 erg->weight_v = 0.0;
4008 switch (rotg->eType)
4010 case EnforcedRotationGroupType::Iso:
4011 case EnforcedRotationGroupType::Isopf:
4012 case EnforcedRotationGroupType::Pm:
4013 case EnforcedRotationGroupType::Pmpf: do_fixed(erg, outstep_rot, outstep_slab); break;
4014 case EnforcedRotationGroupType::Rm:
4015 do_radial_motion(erg, outstep_rot, outstep_slab);
4017 case EnforcedRotationGroupType::Rmpf:
4018 do_radial_motion_pf(erg, x, box, outstep_rot, outstep_slab);
4020 case EnforcedRotationGroupType::Rm2:
4021 case EnforcedRotationGroupType::Rm2pf:
4022 do_radial_motion2(erg, x, box, outstep_rot, outstep_slab);
4024 case EnforcedRotationGroupType::Flext:
4025 case EnforcedRotationGroupType::Flex2t:
4026 /* Subtract the center of the rotation group from the collective positions array
4027 * Also store the center in erg->xc_center since it needs to be subtracted
4028 * in the low level routines from the local coordinates as well */
4029 get_center(erg->xc, erg->mc, rotg->nat, erg->xc_center);
4030 svmul(-1.0, erg->xc_center, transvec);
4031 translate_x(erg->xc, rotg->nat, transvec);
4032 do_flexible(MASTER(cr), er, erg, x, box, t, outstep_rot, outstep_slab);
4034 case EnforcedRotationGroupType::Flex:
4035 case EnforcedRotationGroupType::Flex2:
4036 /* Do NOT subtract the center of mass in the low level routines! */
4037 clear_rvec(erg->xc_center);
4038 do_flexible(MASTER(cr), er, erg, x, box, t, outstep_rot, outstep_slab);
4040 default: gmx_fatal(FARGS, "No such rotation potential.");
4047 fprintf(stderr, "%s calculation (step %d) took %g seconds.\n", RotStr, step, MPI_Wtime() - t0);