Introduce SimulatorBuilder
[gromacs.git] / src / gromacs / pulling / pull_rotation.cpp
blobd74aa558f2180ce34b6154f8b810b2113b9fd2a0
1 /*
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,2018,2019, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
37 #include "gmxpre.h"
39 #include "pull_rotation.h"
41 #include "config.h"
43 #include <cstdio>
44 #include <cstdlib>
45 #include <cstring>
47 #include <algorithm>
48 #include <memory>
50 #include "gromacs/commandline/filenm.h"
51 #include "gromacs/domdec/dlbtiming.h"
52 #include "gromacs/domdec/domdec_struct.h"
53 #include "gromacs/domdec/ga2la.h"
54 #include "gromacs/domdec/localatomset.h"
55 #include "gromacs/domdec/localatomsetmanager.h"
56 #include "gromacs/fileio/gmxfio.h"
57 #include "gromacs/fileio/xvgr.h"
58 #include "gromacs/gmxlib/network.h"
59 #include "gromacs/linearalgebra/nrjac.h"
60 #include "gromacs/math/functions.h"
61 #include "gromacs/math/utilities.h"
62 #include "gromacs/math/vec.h"
63 #include "gromacs/mdlib/groupcoord.h"
64 #include "gromacs/mdlib/stat.h"
65 #include "gromacs/mdrunutility/handlerestart.h"
66 #include "gromacs/mdtypes/commrec.h"
67 #include "gromacs/mdtypes/inputrec.h"
68 #include "gromacs/mdtypes/md_enums.h"
69 #include "gromacs/mdtypes/mdrunoptions.h"
70 #include "gromacs/mdtypes/state.h"
71 #include "gromacs/pbcutil/pbc.h"
72 #include "gromacs/timing/cyclecounter.h"
73 #include "gromacs/timing/wallcycle.h"
74 #include "gromacs/topology/mtop_lookup.h"
75 #include "gromacs/topology/mtop_util.h"
76 #include "gromacs/utility/basedefinitions.h"
77 #include "gromacs/utility/fatalerror.h"
78 #include "gromacs/utility/pleasecite.h"
79 #include "gromacs/utility/smalloc.h"
81 static char const *RotStr = {"Enforced rotation:"};
83 /* Set the minimum weight for the determination of the slab centers */
84 #define WEIGHT_MIN (10*GMX_FLOAT_MIN)
86 //! Helper structure for sorting positions along rotation vector
87 struct sort_along_vec_t
89 //! Projection of xc on the rotation vector
90 real xcproj;
91 //! Index of xc
92 int ind;
93 //! Mass
94 real m;
95 //! Position
96 rvec x;
97 //! Reference position
98 rvec x_ref;
102 //! Enforced rotation / flexible: determine the angle of each slab
103 struct gmx_slabdata
105 //! Number of atoms belonging to this slab
106 int nat;
107 /*! \brief The positions belonging to this slab.
109 * In general, this should be all positions of the whole
110 * rotation group, but we leave those away that have a small
111 * enough weight. */
112 rvec *x;
113 //! Same for reference
114 rvec *ref;
115 //! The weight for each atom
116 real *weight;
120 //! Helper structure for potential fitting
121 struct gmx_potfit
123 /*! \brief Set of angles for which the potential is calculated.
125 * The optimum fit is determined as the angle for with the
126 * potential is minimal. */
127 real *degangle;
128 //! Potential for the different angles
129 real *V;
130 //! Rotation matrix corresponding to the angles
131 matrix *rotmat;
135 //! Enforced rotation data for a single rotation group
136 struct gmx_enfrotgrp
138 //! Input parameters for this group
139 const t_rotgrp *rotg = nullptr;
140 //! Index of this group within the set of groups
141 int groupIndex;
142 //! Rotation angle in degrees
143 real degangle;
144 //! Rotation matrix
145 matrix rotmat;
146 //! The atoms subject to enforced rotation
147 std::unique_ptr<gmx::LocalAtomSet> atomSet;
149 //! The normalized rotation vector
150 rvec vec;
151 //! Rotation potential for this rotation group
152 real V;
153 //! Array to store the forces on the local atoms resulting from enforced rotation potential
154 rvec *f_rot_loc;
156 /* Collective coordinates for the whole rotation group */
157 //! Length of each x_rotref vector after x_rotref has been put into origin
158 real *xc_ref_length;
159 //! Center of the rotation group positions, may be mass weighted
160 rvec xc_center;
161 //! Center of the rotation group reference positions
162 rvec xc_ref_center;
163 //! Current (collective) positions
164 rvec *xc;
165 //! Current (collective) shifts
166 ivec *xc_shifts;
167 //! Extra shifts since last DD step
168 ivec *xc_eshifts;
169 //! Old (collective) positions
170 rvec *xc_old;
171 //! Normalized form of the current positions
172 rvec *xc_norm;
173 //! Reference positions (sorted in the same order as xc when sorted)
174 rvec *xc_ref_sorted;
175 //! Where is a position found after sorting?
176 int *xc_sortind;
177 //! Collective masses
178 real *mc;
179 //! Collective masses sorted
180 real *mc_sorted;
181 //! one over the total mass of the rotation group
182 real invmass;
184 //! Torque in the direction of rotation vector
185 real torque_v;
186 //! Actual angle of the whole rotation group
187 real angle_v;
188 /* Fixed rotation only */
189 //! Weights for angle determination
190 real weight_v;
191 //! Local reference coords, correctly rotated
192 rvec *xr_loc;
193 //! Local current coords, correct PBC image
194 rvec *x_loc_pbc;
195 //! Masses of the current local atoms
196 real *m_loc;
198 /* Flexible rotation only */
199 //! For this many slabs memory is allocated
200 int nslabs_alloc;
201 //! Lowermost slab for that the calculation needs to be performed at a given time step
202 int slab_first;
203 //! Uppermost slab ...
204 int slab_last;
205 //! First slab for which ref. center is stored
206 int slab_first_ref;
207 //! Last ...
208 int slab_last_ref;
209 //! Slab buffer region around reference slabs
210 int slab_buffer;
211 //! First relevant atom for a slab
212 int *firstatom;
213 //! Last relevant atom for a slab
214 int *lastatom;
215 //! Gaussian-weighted slab center
216 rvec *slab_center;
217 //! Gaussian-weighted slab center for the reference positions
218 rvec *slab_center_ref;
219 //! Sum of gaussian weights in a slab
220 real *slab_weights;
221 //! Torque T = r x f for each slab. torque_v = m.v = angular momentum in the direction of v
222 real *slab_torque_v;
223 //! 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
224 real max_beta;
225 //! Precalculated gaussians for a single atom
226 real *gn_atom;
227 //! Tells to which slab each precalculated gaussian belongs
228 int *gn_slabind;
229 //! Inner sum of the flexible2 potential per slab; this is precalculated for optimization reasons
230 rvec *slab_innersumvec;
231 //! Holds atom positions and gaussian weights of atoms belonging to a slab
232 gmx_slabdata *slab_data;
234 /* For potential fits with varying angle: */
235 //! Used for fit type 'potential'
236 gmx_potfit *PotAngleFit;
240 //! Enforced rotation data for all groups
241 struct gmx_enfrot
243 //! Input parameters.
244 const t_rot *rot = nullptr;
245 //! Output period for main rotation outfile
246 int nstrout;
247 //! Output period for per-slab data
248 int nstsout;
249 //! Output file for rotation data
250 FILE *out_rot = nullptr;
251 //! Output file for torque data
252 FILE *out_torque = nullptr;
253 //! Output file for slab angles for flexible type
254 FILE *out_angles = nullptr;
255 //! Output file for slab centers
256 FILE *out_slabs = nullptr;
257 //! Allocation size of buf
258 int bufsize = 0;
259 //! Coordinate buffer variable for sorting
260 rvec *xbuf = nullptr;
261 //! Masses buffer variable for sorting
262 real *mbuf = nullptr;
263 //! Buffer variable needed for position sorting
264 sort_along_vec_t *data = nullptr;
265 //! MPI buffer
266 real *mpi_inbuf = nullptr;
267 //! MPI buffer
268 real *mpi_outbuf = nullptr;
269 //! Allocation size of in & outbuf
270 int mpi_bufsize = 0;
271 //! If true, append output files
272 gmx_bool restartWithAppending = false;
273 //! Used to skip first output when appending to avoid duplicate entries in rotation outfiles
274 gmx_bool bOut = false;
275 //! Stores working data per group
276 std::vector<gmx_enfrotgrp> enfrotgrp;
277 ~gmx_enfrot();
280 gmx_enfrot::~gmx_enfrot()
282 if (out_rot)
284 gmx_fio_fclose(out_rot);
286 if (out_slabs)
288 gmx_fio_fclose(out_slabs);
290 if (out_angles)
292 gmx_fio_fclose(out_angles);
294 if (out_torque)
296 gmx_fio_fclose(out_torque);
300 namespace gmx
303 class EnforcedRotation::Impl
305 public:
306 gmx_enfrot enforcedRotation_;
309 EnforcedRotation::EnforcedRotation() : impl_(new Impl)
313 EnforcedRotation::~EnforcedRotation() = default;
315 gmx_enfrot *EnforcedRotation::getLegacyEnfrot()
317 return &impl_->enforcedRotation_;
320 } // namespace gmx
322 /* Activate output of forces for correctness checks */
323 /* #define PRINT_FORCES */
324 #ifdef PRINT_FORCES
325 #define PRINT_FORCE_J fprintf(stderr, "f%d = %15.8f %15.8f %15.8f\n", erg->xc_ref_ind[j], erg->f_rot_loc[j][XX], erg->f_rot_loc[j][YY], erg->f_rot_loc[j][ZZ]);
326 #define PRINT_POT_TAU if (MASTER(cr)) { \
327 fprintf(stderr, "potential = %15.8f\n" "torque = %15.8f\n", erg->V, erg->torque_v); \
329 #else
330 #define PRINT_FORCE_J
331 #define PRINT_POT_TAU
332 #endif
334 /* Shortcuts for often used queries */
335 #define ISFLEX(rg) ( ((rg)->eType == erotgFLEX) || ((rg)->eType == erotgFLEXT) || ((rg)->eType == erotgFLEX2) || ((rg)->eType == erotgFLEX2T) )
336 #define ISCOLL(rg) ( ((rg)->eType == erotgFLEX) || ((rg)->eType == erotgFLEXT) || ((rg)->eType == erotgFLEX2) || ((rg)->eType == erotgFLEX2T) || ((rg)->eType == erotgRMPF) || ((rg)->eType == erotgRM2PF) )
339 /* Does any of the rotation groups use slab decomposition? */
340 static gmx_bool HaveFlexibleGroups(const t_rot *rot)
342 for (int g = 0; g < rot->ngrp; g++)
344 if (ISFLEX(&rot->grp[g]))
346 return TRUE;
350 return FALSE;
354 /* Is for any group the fit angle determined by finding the minimum of the
355 * rotation potential? */
356 static gmx_bool HavePotFitGroups(const t_rot *rot)
358 for (int g = 0; g < rot->ngrp; g++)
360 if (erotgFitPOT == rot->grp[g].eFittype)
362 return TRUE;
366 return FALSE;
370 static double** allocate_square_matrix(int dim)
372 int i;
373 double** mat = nullptr;
376 snew(mat, dim);
377 for (i = 0; i < dim; i++)
379 snew(mat[i], dim);
382 return mat;
386 static void free_square_matrix(double** mat, int dim)
388 int i;
391 for (i = 0; i < dim; i++)
393 sfree(mat[i]);
395 sfree(mat);
399 /* Return the angle for which the potential is minimal */
400 static real get_fitangle(const gmx_enfrotgrp *erg)
402 int i;
403 real fitangle = -999.9;
404 real pot_min = GMX_FLOAT_MAX;
405 gmx_potfit *fit;
408 fit = erg->PotAngleFit;
410 for (i = 0; i < erg->rotg->PotAngle_nstep; i++)
412 if (fit->V[i] < pot_min)
414 pot_min = fit->V[i];
415 fitangle = fit->degangle[i];
419 return fitangle;
423 /* Reduce potential angle fit data for this group at this time step? */
424 static inline gmx_bool bPotAngle(const gmx_enfrot *er, const t_rotgrp *rotg, int64_t step)
426 return ( (erotgFitPOT == rotg->eFittype) && (do_per_step(step, er->nstsout) || do_per_step(step, er->nstrout)) );
429 /* Reduce slab torqe data for this group at this time step? */
430 static inline gmx_bool bSlabTau(const gmx_enfrot *er, const t_rotgrp *rotg, int64_t step)
432 return ( (ISFLEX(rotg)) && do_per_step(step, er->nstsout) );
435 /* Output rotation energy, torques, etc. for each rotation group */
436 static void reduce_output(const t_commrec *cr,
437 gmx_enfrot *er, real t, int64_t step)
439 int i, islab, nslabs = 0;
440 int count; /* MPI element counter */
441 real fitangle;
442 gmx_bool bFlex;
444 /* Fill the MPI buffer with stuff to reduce. If items are added for reduction
445 * here, the MPI buffer size has to be enlarged also in calc_mpi_bufsize() */
446 if (PAR(cr))
448 count = 0;
449 for (auto &ergRef : er->enfrotgrp)
451 gmx_enfrotgrp *erg = &ergRef;
452 const t_rotgrp *rotg = erg->rotg;
453 nslabs = erg->slab_last - erg->slab_first + 1;
454 er->mpi_inbuf[count++] = erg->V;
455 er->mpi_inbuf[count++] = erg->torque_v;
456 er->mpi_inbuf[count++] = erg->angle_v;
457 er->mpi_inbuf[count++] = erg->weight_v; /* weights are not needed for flex types, but this is just a single value */
459 if (bPotAngle(er, rotg, step))
461 for (i = 0; i < rotg->PotAngle_nstep; i++)
463 er->mpi_inbuf[count++] = erg->PotAngleFit->V[i];
466 if (bSlabTau(er, rotg, step))
468 for (i = 0; i < nslabs; i++)
470 er->mpi_inbuf[count++] = erg->slab_torque_v[i];
474 if (count > er->mpi_bufsize)
476 gmx_fatal(FARGS, "%s MPI buffer overflow, please report this error.", RotStr);
479 #if GMX_MPI
480 MPI_Reduce(er->mpi_inbuf, er->mpi_outbuf, count, GMX_MPI_REAL, MPI_SUM, MASTERRANK(cr), cr->mpi_comm_mygroup);
481 #endif
483 /* Copy back the reduced data from the buffer on the master */
484 if (MASTER(cr))
486 count = 0;
487 for (auto &ergRef : er->enfrotgrp)
489 gmx_enfrotgrp *erg = &ergRef;
490 const t_rotgrp *rotg = erg->rotg;
491 nslabs = erg->slab_last - erg->slab_first + 1;
492 erg->V = er->mpi_outbuf[count++];
493 erg->torque_v = er->mpi_outbuf[count++];
494 erg->angle_v = er->mpi_outbuf[count++];
495 erg->weight_v = er->mpi_outbuf[count++];
497 if (bPotAngle(er, rotg, step))
499 for (int i = 0; i < rotg->PotAngle_nstep; i++)
501 erg->PotAngleFit->V[i] = er->mpi_outbuf[count++];
504 if (bSlabTau(er, rotg, step))
506 for (int i = 0; i < nslabs; i++)
508 erg->slab_torque_v[i] = er->mpi_outbuf[count++];
515 /* Output */
516 if (MASTER(cr))
518 /* Angle and torque for each rotation group */
519 for (auto &ergRef : er->enfrotgrp)
521 gmx_enfrotgrp *erg = &ergRef;
522 const t_rotgrp *rotg = erg->rotg;
523 bFlex = ISFLEX(rotg);
525 /* Output to main rotation output file: */
526 if (do_per_step(step, er->nstrout) )
528 if (erotgFitPOT == rotg->eFittype)
530 fitangle = get_fitangle(erg);
532 else
534 if (bFlex)
536 fitangle = erg->angle_v; /* RMSD fit angle */
538 else
540 fitangle = (erg->angle_v/erg->weight_v)*180.0*M_1_PI;
543 fprintf(er->out_rot, "%12.4f", fitangle);
544 fprintf(er->out_rot, "%12.3e", erg->torque_v);
545 fprintf(er->out_rot, "%12.3e", erg->V);
548 if (do_per_step(step, er->nstsout) )
550 /* Output to torque log file: */
551 if (bFlex)
553 fprintf(er->out_torque, "%12.3e%6d", t, erg->groupIndex);
554 for (int i = erg->slab_first; i <= erg->slab_last; i++)
556 islab = i - erg->slab_first; /* slab index */
557 /* Only output if enough weight is in slab */
558 if (erg->slab_weights[islab] > rotg->min_gaussian)
560 fprintf(er->out_torque, "%6d%12.3e", i, erg->slab_torque_v[islab]);
563 fprintf(er->out_torque, "\n");
566 /* Output to angles log file: */
567 if (erotgFitPOT == rotg->eFittype)
569 fprintf(er->out_angles, "%12.3e%6d%12.4f", t, erg->groupIndex, erg->degangle);
570 /* Output energies at a set of angles around the reference angle */
571 for (int i = 0; i < rotg->PotAngle_nstep; i++)
573 fprintf(er->out_angles, "%12.3e", erg->PotAngleFit->V[i]);
575 fprintf(er->out_angles, "\n");
579 if (do_per_step(step, er->nstrout) )
581 fprintf(er->out_rot, "\n");
587 /* Add the forces from enforced rotation potential to the local forces.
588 * Should be called after the SR forces have been evaluated */
589 real add_rot_forces(gmx_enfrot *er,
590 rvec f[], const t_commrec *cr, int64_t step, real t)
592 real Vrot = 0.0; /* If more than one rotation group is present, Vrot
593 assembles the local parts from all groups */
595 /* Loop over enforced rotation groups (usually 1, though)
596 * Apply the forces from rotation potentials */
597 for (auto &ergRef : er->enfrotgrp)
599 gmx_enfrotgrp *erg = &ergRef;
600 Vrot += erg->V; /* add the local parts from the nodes */
601 const auto &localRotationGroupIndex = erg->atomSet->localIndex();
602 for (gmx::index l = 0; l < localRotationGroupIndex.ssize(); l++)
604 /* Get the right index of the local force */
605 int ii = localRotationGroupIndex[l];
606 /* Add */
607 rvec_inc(f[ii], erg->f_rot_loc[l]);
611 /* Reduce energy,torque, angles etc. to get the sum values (per rotation group)
612 * on the master and output these values to file. */
613 if ( (do_per_step(step, er->nstrout) || do_per_step(step, er->nstsout)) && er->bOut)
615 reduce_output(cr, er, t, step);
618 /* When appending, er->bOut is FALSE the first time to avoid duplicate entries */
619 er->bOut = TRUE;
621 PRINT_POT_TAU
623 return Vrot;
627 /* The Gaussian norm is chosen such that the sum of the gaussian functions
628 * over the slabs is approximately 1.0 everywhere */
629 #define GAUSS_NORM 0.569917543430618
632 /* Calculate the maximum beta that leads to a gaussian larger min_gaussian,
633 * also does some checks
635 static double calc_beta_max(real min_gaussian, real slab_dist)
637 double sigma;
638 double arg;
641 /* Actually the next two checks are already made in grompp */
642 if (slab_dist <= 0)
644 gmx_fatal(FARGS, "Slab distance of flexible rotation groups must be >=0 !");
646 if (min_gaussian <= 0)
648 gmx_fatal(FARGS, "Cutoff value for Gaussian must be > 0. (You requested %f)", min_gaussian);
651 /* Define the sigma value */
652 sigma = 0.7*slab_dist;
654 /* Calculate the argument for the logarithm and check that the log() result is negative or 0 */
655 arg = min_gaussian/GAUSS_NORM;
656 if (arg > 1.0)
658 gmx_fatal(FARGS, "min_gaussian of flexible rotation groups must be <%g", GAUSS_NORM);
661 return std::sqrt(-2.0*sigma*sigma*log(min_gaussian/GAUSS_NORM));
665 static inline real calc_beta(rvec curr_x, const gmx_enfrotgrp *erg, int n)
667 return iprod(curr_x, erg->vec) - erg->rotg->slab_dist * n;
671 static inline real gaussian_weight(rvec curr_x, const gmx_enfrotgrp *erg, int n)
673 const real norm = GAUSS_NORM;
674 real sigma;
677 /* Define the sigma value */
678 sigma = 0.7*erg->rotg->slab_dist;
679 /* Calculate the Gaussian value of slab n for position curr_x */
680 return norm * exp( -0.5 * gmx::square( calc_beta(curr_x, erg, n)/sigma ) );
684 /* Returns the weight in a single slab, also calculates the Gaussian- and mass-
685 * weighted sum of positions for that slab */
686 static real get_slab_weight(int j, const gmx_enfrotgrp *erg,
687 rvec xc[], const real mc[], rvec *x_weighted_sum)
689 rvec curr_x; /* The position of an atom */
690 rvec curr_x_weighted; /* The gaussian-weighted position */
691 real gaussian; /* A single gaussian weight */
692 real wgauss; /* gaussian times current mass */
693 real slabweight = 0.0; /* The sum of weights in the slab */
695 clear_rvec(*x_weighted_sum);
697 /* Loop over all atoms in the rotation group */
698 for (int i = 0; i < erg->rotg->nat; i++)
700 copy_rvec(xc[i], curr_x);
701 gaussian = gaussian_weight(curr_x, erg, j);
702 wgauss = gaussian * mc[i];
703 svmul(wgauss, curr_x, curr_x_weighted);
704 rvec_add(*x_weighted_sum, curr_x_weighted, *x_weighted_sum);
705 slabweight += wgauss;
706 } /* END of loop over rotation group atoms */
708 return slabweight;
712 static void get_slab_centers(
713 gmx_enfrotgrp *erg, /* Enforced rotation group working data */
714 rvec *xc, /* The rotation group positions; will
715 typically be enfrotgrp->xc, but at first call
716 it is enfrotgrp->xc_ref */
717 real *mc, /* The masses of the rotation group atoms */
718 real time, /* Used for output only */
719 FILE *out_slabs, /* For outputting center per slab information */
720 gmx_bool bOutStep, /* Is this an output step? */
721 gmx_bool bReference) /* If this routine is called from
722 init_rot_group we need to store
723 the reference slab centers */
725 /* Loop over slabs */
726 for (int j = erg->slab_first; j <= erg->slab_last; j++)
728 int slabIndex = j - erg->slab_first;
729 erg->slab_weights[slabIndex] = get_slab_weight(j, erg, xc, mc, &erg->slab_center[slabIndex]);
731 /* We can do the calculations ONLY if there is weight in the slab! */
732 if (erg->slab_weights[slabIndex] > WEIGHT_MIN)
734 svmul(1.0/erg->slab_weights[slabIndex], erg->slab_center[slabIndex], erg->slab_center[slabIndex]);
736 else
738 /* We need to check this here, since we divide through slab_weights
739 * in the flexible low-level routines! */
740 gmx_fatal(FARGS, "Not enough weight in slab %d. Slab center cannot be determined!", j);
743 /* At first time step: save the centers of the reference structure */
744 if (bReference)
746 copy_rvec(erg->slab_center[slabIndex], erg->slab_center_ref[slabIndex]);
748 } /* END of loop over slabs */
750 /* Output on the master */
751 if ( (nullptr != out_slabs) && bOutStep)
753 fprintf(out_slabs, "%12.3e%6d", time, erg->groupIndex);
754 for (int j = erg->slab_first; j <= erg->slab_last; j++)
756 int slabIndex = j - erg->slab_first;
757 fprintf(out_slabs, "%6d%12.3e%12.3e%12.3e",
758 j, erg->slab_center[slabIndex][XX], erg->slab_center[slabIndex][YY], erg->slab_center[slabIndex][ZZ]);
760 fprintf(out_slabs, "\n");
765 static void calc_rotmat(
766 const rvec vec,
767 real degangle, /* Angle alpha of rotation at time t in degrees */
768 matrix rotmat) /* Rotation matrix */
770 real radangle; /* Rotation angle in radians */
771 real cosa; /* cosine alpha */
772 real sina; /* sine alpha */
773 real OMcosa; /* 1 - cos(alpha) */
774 real dumxy, dumxz, dumyz; /* save computations */
775 rvec rot_vec; /* Rotate around rot_vec ... */
778 radangle = degangle * M_PI/180.0;
779 copy_rvec(vec, rot_vec );
781 /* Precompute some variables: */
782 cosa = cos(radangle);
783 sina = sin(radangle);
784 OMcosa = 1.0 - cosa;
785 dumxy = rot_vec[XX]*rot_vec[YY]*OMcosa;
786 dumxz = rot_vec[XX]*rot_vec[ZZ]*OMcosa;
787 dumyz = rot_vec[YY]*rot_vec[ZZ]*OMcosa;
789 /* Construct the rotation matrix for this rotation group: */
790 /* 1st column: */
791 rotmat[XX][XX] = cosa + rot_vec[XX]*rot_vec[XX]*OMcosa;
792 rotmat[YY][XX] = dumxy + rot_vec[ZZ]*sina;
793 rotmat[ZZ][XX] = dumxz - rot_vec[YY]*sina;
794 /* 2nd column: */
795 rotmat[XX][YY] = dumxy - rot_vec[ZZ]*sina;
796 rotmat[YY][YY] = cosa + rot_vec[YY]*rot_vec[YY]*OMcosa;
797 rotmat[ZZ][YY] = dumyz + rot_vec[XX]*sina;
798 /* 3rd column: */
799 rotmat[XX][ZZ] = dumxz + rot_vec[YY]*sina;
800 rotmat[YY][ZZ] = dumyz - rot_vec[XX]*sina;
801 rotmat[ZZ][ZZ] = cosa + rot_vec[ZZ]*rot_vec[ZZ]*OMcosa;
803 #ifdef PRINTMATRIX
804 int iii, jjj;
806 for (iii = 0; iii < 3; iii++)
808 for (jjj = 0; jjj < 3; jjj++)
810 fprintf(stderr, " %10.8f ", rotmat[iii][jjj]);
812 fprintf(stderr, "\n");
814 #endif
818 /* Calculates torque on the rotation axis tau = position x force */
819 static inline real torque(const rvec rotvec, /* rotation vector; MUST be normalized! */
820 rvec force, /* force */
821 rvec x, /* position of atom on which the force acts */
822 rvec pivot) /* pivot point of rotation axis */
824 rvec vectmp, tau;
827 /* Subtract offset */
828 rvec_sub(x, pivot, vectmp);
830 /* position x force */
831 cprod(vectmp, force, tau);
833 /* Return the part of the torque which is parallel to the rotation vector */
834 return iprod(tau, rotvec);
838 /* Right-aligned output of value with standard width */
839 static void print_aligned(FILE *fp, char const *str)
841 fprintf(fp, "%12s", str);
845 /* Right-aligned output of value with standard short width */
846 static void print_aligned_short(FILE *fp, char const *str)
848 fprintf(fp, "%6s", str);
852 static FILE *open_output_file(const char *fn, int steps, const char what[])
854 FILE *fp;
857 fp = gmx_ffopen(fn, "w");
859 fprintf(fp, "# Output of %s is written in intervals of %d time step%s.\n#\n",
860 what, steps, steps > 1 ? "s" : "");
862 return fp;
866 /* Open output file for slab center data. Call on master only */
867 static FILE *open_slab_out(const char *fn,
868 gmx_enfrot *er)
870 FILE *fp;
872 if (er->restartWithAppending)
874 fp = gmx_fio_fopen(fn, "a");
876 else
878 fp = open_output_file(fn, er->nstsout, "gaussian weighted slab centers");
880 for (auto &ergRef : er->enfrotgrp)
882 gmx_enfrotgrp *erg = &ergRef;
883 if (ISFLEX(erg->rotg))
885 fprintf(fp, "# Rotation group %d (%s), slab distance %f nm, %s.\n",
886 erg->groupIndex, erotg_names[erg->rotg->eType], erg->rotg->slab_dist,
887 erg->rotg->bMassW ? "centers of mass" : "geometrical centers");
891 fprintf(fp, "# Reference centers are listed first (t=-1).\n");
892 fprintf(fp, "# The following columns have the syntax:\n");
893 fprintf(fp, "# ");
894 print_aligned_short(fp, "t");
895 print_aligned_short(fp, "grp");
896 /* Print legend for the first two entries only ... */
897 for (int i = 0; i < 2; i++)
899 print_aligned_short(fp, "slab");
900 print_aligned(fp, "X center");
901 print_aligned(fp, "Y center");
902 print_aligned(fp, "Z center");
904 fprintf(fp, " ...\n");
905 fflush(fp);
908 return fp;
912 /* Adds 'buf' to 'str' */
913 static void add_to_string(char **str, char *buf)
915 int len;
918 len = strlen(*str) + strlen(buf) + 1;
919 srenew(*str, len);
920 strcat(*str, buf);
924 static void add_to_string_aligned(char **str, char *buf)
926 char buf_aligned[STRLEN];
928 sprintf(buf_aligned, "%12s", buf);
929 add_to_string(str, buf_aligned);
933 /* Open output file and print some general information about the rotation groups.
934 * Call on master only */
935 static FILE *open_rot_out(const char *fn,
936 const gmx_output_env_t *oenv,
937 gmx_enfrot *er)
939 FILE *fp;
940 int nsets;
941 const char **setname;
942 char buf[50], buf2[75];
943 gmx_bool bFlex;
944 char *LegendStr = nullptr;
945 const t_rot *rot = er->rot;
947 if (er->restartWithAppending)
949 fp = gmx_fio_fopen(fn, "a");
951 else
953 fp = xvgropen(fn, "Rotation angles and energy", "Time (ps)", "angles (degrees) and energies (kJ/mol)", oenv);
954 fprintf(fp, "# Output of enforced rotation data is written in intervals of %d time step%s.\n#\n", er->nstrout, er->nstrout > 1 ? "s" : "");
955 fprintf(fp, "# The scalar tau is the torque (kJ/mol) in the direction of the rotation vector v.\n");
956 fprintf(fp, "# To obtain the vectorial torque, multiply tau with the group's rot-vec.\n");
957 fprintf(fp, "# For flexible groups, tau(t,n) from all slabs n have been summed in a single value tau(t) here.\n");
958 fprintf(fp, "# The torques tau(t,n) are found in the rottorque.log (-rt) output file\n");
960 for (int g = 0; g < rot->ngrp; g++)
962 const t_rotgrp *rotg = &rot->grp[g];
963 const gmx_enfrotgrp *erg = &er->enfrotgrp[g];
964 bFlex = ISFLEX(rotg);
966 fprintf(fp, "#\n");
967 fprintf(fp, "# ROTATION GROUP %d, potential type '%s':\n", g, erotg_names[rotg->eType]);
968 fprintf(fp, "# rot-massw%d %s\n", g, yesno_names[rotg->bMassW]);
969 fprintf(fp, "# rot-vec%d %12.5e %12.5e %12.5e\n", g, erg->vec[XX], erg->vec[YY], erg->vec[ZZ]);
970 fprintf(fp, "# rot-rate%d %12.5e degrees/ps\n", g, rotg->rate);
971 fprintf(fp, "# rot-k%d %12.5e kJ/(mol*nm^2)\n", g, rotg->k);
972 if (rotg->eType == erotgISO || rotg->eType == erotgPM || rotg->eType == erotgRM || rotg->eType == erotgRM2)
974 fprintf(fp, "# rot-pivot%d %12.5e %12.5e %12.5e nm\n", g, rotg->pivot[XX], rotg->pivot[YY], rotg->pivot[ZZ]);
977 if (bFlex)
979 fprintf(fp, "# rot-slab-distance%d %f nm\n", g, rotg->slab_dist);
980 fprintf(fp, "# rot-min-gaussian%d %12.5e\n", g, rotg->min_gaussian);
983 /* Output the centers of the rotation groups for the pivot-free potentials */
984 if ((rotg->eType == erotgISOPF) || (rotg->eType == erotgPMPF) || (rotg->eType == erotgRMPF) || (rotg->eType == erotgRM2PF
985 || (rotg->eType == erotgFLEXT) || (rotg->eType == erotgFLEX2T)) )
987 fprintf(fp, "# ref. grp. %d center %12.5e %12.5e %12.5e\n", g,
988 erg->xc_ref_center[XX], erg->xc_ref_center[YY], erg->xc_ref_center[ZZ]);
990 fprintf(fp, "# grp. %d init.center %12.5e %12.5e %12.5e\n", g,
991 erg->xc_center[XX], erg->xc_center[YY], erg->xc_center[ZZ]);
994 if ( (rotg->eType == erotgRM2) || (rotg->eType == erotgFLEX2) || (rotg->eType == erotgFLEX2T) )
996 fprintf(fp, "# rot-eps%d %12.5e nm^2\n", g, rotg->eps);
998 if (erotgFitPOT == rotg->eFittype)
1000 fprintf(fp, "#\n");
1001 fprintf(fp, "# theta_fit%d is determined by first evaluating the potential for %d angles around theta_ref%d.\n",
1002 g, rotg->PotAngle_nstep, g);
1003 fprintf(fp, "# The fit angle is the one with the smallest potential. It is given as the deviation\n");
1004 fprintf(fp, "# from the reference angle, i.e. if theta_ref=X and theta_fit=Y, then the angle with\n");
1005 fprintf(fp, "# minimal value of the potential is X+Y. Angular resolution is %g degrees.\n", rotg->PotAngle_step);
1009 /* Print a nice legend */
1010 snew(LegendStr, 1);
1011 LegendStr[0] = '\0';
1012 sprintf(buf, "# %6s", "time");
1013 add_to_string_aligned(&LegendStr, buf);
1015 nsets = 0;
1016 snew(setname, 4*rot->ngrp);
1018 for (int g = 0; g < rot->ngrp; g++)
1020 sprintf(buf, "theta_ref%d", g);
1021 add_to_string_aligned(&LegendStr, buf);
1023 sprintf(buf2, "%s (degrees)", buf);
1024 setname[nsets] = gmx_strdup(buf2);
1025 nsets++;
1027 for (int g = 0; g < rot->ngrp; g++)
1029 const t_rotgrp *rotg = &rot->grp[g];
1030 bFlex = ISFLEX(rotg);
1032 /* For flexible axis rotation we use RMSD fitting to determine the
1033 * actual angle of the rotation group */
1034 if (bFlex || erotgFitPOT == rotg->eFittype)
1036 sprintf(buf, "theta_fit%d", g);
1038 else
1040 sprintf(buf, "theta_av%d", g);
1042 add_to_string_aligned(&LegendStr, buf);
1043 sprintf(buf2, "%s (degrees)", buf);
1044 setname[nsets] = gmx_strdup(buf2);
1045 nsets++;
1047 sprintf(buf, "tau%d", g);
1048 add_to_string_aligned(&LegendStr, buf);
1049 sprintf(buf2, "%s (kJ/mol)", buf);
1050 setname[nsets] = gmx_strdup(buf2);
1051 nsets++;
1053 sprintf(buf, "energy%d", g);
1054 add_to_string_aligned(&LegendStr, buf);
1055 sprintf(buf2, "%s (kJ/mol)", buf);
1056 setname[nsets] = gmx_strdup(buf2);
1057 nsets++;
1059 fprintf(fp, "#\n");
1061 if (nsets > 1)
1063 xvgr_legend(fp, nsets, setname, oenv);
1065 sfree(setname);
1067 fprintf(fp, "#\n# Legend for the following data columns:\n");
1068 fprintf(fp, "%s\n", LegendStr);
1069 sfree(LegendStr);
1071 fflush(fp);
1074 return fp;
1078 /* Call on master only */
1079 static FILE *open_angles_out(const char *fn,
1080 gmx_enfrot *er)
1082 FILE *fp;
1083 char buf[100];
1084 const t_rot *rot = er->rot;
1086 if (er->restartWithAppending)
1088 fp = gmx_fio_fopen(fn, "a");
1090 else
1092 /* Open output file and write some information about it's structure: */
1093 fp = open_output_file(fn, er->nstsout, "rotation group angles");
1094 fprintf(fp, "# All angles given in degrees, time in ps.\n");
1095 for (int g = 0; g < rot->ngrp; g++)
1097 const t_rotgrp *rotg = &rot->grp[g];
1098 const gmx_enfrotgrp *erg = &er->enfrotgrp[g];
1100 /* Output for this group happens only if potential type is flexible or
1101 * if fit type is potential! */
1102 if (ISFLEX(rotg) || (erotgFitPOT == rotg->eFittype) )
1104 if (ISFLEX(rotg))
1106 sprintf(buf, " slab distance %f nm, ", rotg->slab_dist);
1108 else
1110 buf[0] = '\0';
1113 fprintf(fp, "#\n# ROTATION GROUP %d '%s',%s fit type '%s'.\n",
1114 g, erotg_names[rotg->eType], buf, erotg_fitnames[rotg->eFittype]);
1116 /* Special type of fitting using the potential minimum. This is
1117 * done for the whole group only, not for the individual slabs. */
1118 if (erotgFitPOT == rotg->eFittype)
1120 fprintf(fp, "# To obtain theta_fit%d, the potential is evaluated for %d angles around theta_ref%d\n", g, rotg->PotAngle_nstep, g);
1121 fprintf(fp, "# The fit angle in the rotation standard outfile is the one with minimal energy E(theta_fit) [kJ/mol].\n");
1122 fprintf(fp, "#\n");
1125 fprintf(fp, "# Legend for the group %d data columns:\n", g);
1126 fprintf(fp, "# ");
1127 print_aligned_short(fp, "time");
1128 print_aligned_short(fp, "grp");
1129 print_aligned(fp, "theta_ref");
1131 if (erotgFitPOT == rotg->eFittype)
1133 /* Output the set of angles around the reference angle */
1134 for (int i = 0; i < rotg->PotAngle_nstep; i++)
1136 sprintf(buf, "E(%g)", erg->PotAngleFit->degangle[i]);
1137 print_aligned(fp, buf);
1140 else
1142 /* Output fit angle for each slab */
1143 print_aligned_short(fp, "slab");
1144 print_aligned_short(fp, "atoms");
1145 print_aligned(fp, "theta_fit");
1146 print_aligned_short(fp, "slab");
1147 print_aligned_short(fp, "atoms");
1148 print_aligned(fp, "theta_fit");
1149 fprintf(fp, " ...");
1151 fprintf(fp, "\n");
1154 fflush(fp);
1157 return fp;
1161 /* Open torque output file and write some information about it's structure.
1162 * Call on master only */
1163 static FILE *open_torque_out(const char *fn,
1164 gmx_enfrot *er)
1166 FILE *fp;
1167 const t_rot *rot = er->rot;
1169 if (er->restartWithAppending)
1171 fp = gmx_fio_fopen(fn, "a");
1173 else
1175 fp = open_output_file(fn, er->nstsout, "torques");
1177 for (int g = 0; g < rot->ngrp; g++)
1179 const t_rotgrp *rotg = &rot->grp[g];
1180 const gmx_enfrotgrp *erg = &er->enfrotgrp[g];
1181 if (ISFLEX(rotg))
1183 fprintf(fp, "# Rotation group %d (%s), slab distance %f nm.\n", g, erotg_names[rotg->eType], rotg->slab_dist);
1184 fprintf(fp, "# The scalar tau is the torque (kJ/mol) in the direction of the rotation vector.\n");
1185 fprintf(fp, "# To obtain the vectorial torque, multiply tau with\n");
1186 fprintf(fp, "# rot-vec%d %10.3e %10.3e %10.3e\n", g, erg->vec[XX], erg->vec[YY], erg->vec[ZZ]);
1187 fprintf(fp, "#\n");
1190 fprintf(fp, "# Legend for the following data columns: (tau=torque for that slab):\n");
1191 fprintf(fp, "# ");
1192 print_aligned_short(fp, "t");
1193 print_aligned_short(fp, "grp");
1194 print_aligned_short(fp, "slab");
1195 print_aligned(fp, "tau");
1196 print_aligned_short(fp, "slab");
1197 print_aligned(fp, "tau");
1198 fprintf(fp, " ...\n");
1199 fflush(fp);
1202 return fp;
1206 static void swap_val(double* vec, int i, int j)
1208 double tmp = vec[j];
1211 vec[j] = vec[i];
1212 vec[i] = tmp;
1216 static void swap_col(double **mat, int i, int j)
1218 double tmp[3] = {mat[0][j], mat[1][j], mat[2][j]};
1221 mat[0][j] = mat[0][i];
1222 mat[1][j] = mat[1][i];
1223 mat[2][j] = mat[2][i];
1225 mat[0][i] = tmp[0];
1226 mat[1][i] = tmp[1];
1227 mat[2][i] = tmp[2];
1231 /* Eigenvectors are stored in columns of eigen_vec */
1232 static void diagonalize_symmetric(
1233 double **matrix,
1234 double **eigen_vec,
1235 double eigenval[3])
1237 int n_rot;
1240 jacobi(matrix, 3, eigenval, eigen_vec, &n_rot);
1242 /* sort in ascending order */
1243 if (eigenval[0] > eigenval[1])
1245 swap_val(eigenval, 0, 1);
1246 swap_col(eigen_vec, 0, 1);
1248 if (eigenval[1] > eigenval[2])
1250 swap_val(eigenval, 1, 2);
1251 swap_col(eigen_vec, 1, 2);
1253 if (eigenval[0] > eigenval[1])
1255 swap_val(eigenval, 0, 1);
1256 swap_col(eigen_vec, 0, 1);
1261 static void align_with_z(
1262 rvec* s, /* Structure to align */
1263 int natoms,
1264 rvec axis)
1266 int i, j, k;
1267 rvec zet = {0.0, 0.0, 1.0};
1268 rvec rot_axis = {0.0, 0.0, 0.0};
1269 rvec *rotated_str = nullptr;
1270 real ooanorm;
1271 real angle;
1272 matrix rotmat;
1275 snew(rotated_str, natoms);
1277 /* Normalize the axis */
1278 ooanorm = 1.0/norm(axis);
1279 svmul(ooanorm, axis, axis);
1281 /* Calculate the angle for the fitting procedure */
1282 cprod(axis, zet, rot_axis);
1283 angle = acos(axis[2]);
1284 if (angle < 0.0)
1286 angle += M_PI;
1289 /* Calculate the rotation matrix */
1290 calc_rotmat(rot_axis, angle*180.0/M_PI, rotmat);
1292 /* Apply the rotation matrix to s */
1293 for (i = 0; i < natoms; i++)
1295 for (j = 0; j < 3; j++)
1297 for (k = 0; k < 3; k++)
1299 rotated_str[i][j] += rotmat[j][k]*s[i][k];
1304 /* Rewrite the rotated structure to s */
1305 for (i = 0; i < natoms; i++)
1307 for (j = 0; j < 3; j++)
1309 s[i][j] = rotated_str[i][j];
1313 sfree(rotated_str);
1317 static void calc_correl_matrix(rvec* Xstr, rvec* Ystr, double** Rmat, int natoms)
1319 int i, j, k;
1322 for (i = 0; i < 3; i++)
1324 for (j = 0; j < 3; j++)
1326 Rmat[i][j] = 0.0;
1330 for (i = 0; i < 3; i++)
1332 for (j = 0; j < 3; j++)
1334 for (k = 0; k < natoms; k++)
1336 Rmat[i][j] += Ystr[k][i] * Xstr[k][j];
1343 static void weigh_coords(rvec* str, real* weight, int natoms)
1345 int i, j;
1348 for (i = 0; i < natoms; i++)
1350 for (j = 0; j < 3; j++)
1352 str[i][j] *= std::sqrt(weight[i]);
1358 static real opt_angle_analytic(
1359 rvec * ref_s,
1360 rvec * act_s,
1361 real * weight,
1362 int natoms,
1363 const rvec ref_com,
1364 const rvec act_com,
1365 rvec axis)
1367 int i, j, k;
1368 rvec *ref_s_1 = nullptr;
1369 rvec *act_s_1 = nullptr;
1370 rvec shift;
1371 double **Rmat, **RtR, **eigvec;
1372 double eigval[3];
1373 double V[3][3], WS[3][3];
1374 double rot_matrix[3][3];
1375 double opt_angle;
1378 /* Do not change the original coordinates */
1379 snew(ref_s_1, natoms);
1380 snew(act_s_1, natoms);
1381 for (i = 0; i < natoms; i++)
1383 copy_rvec(ref_s[i], ref_s_1[i]);
1384 copy_rvec(act_s[i], act_s_1[i]);
1387 /* Translate the structures to the origin */
1388 shift[XX] = -ref_com[XX];
1389 shift[YY] = -ref_com[YY];
1390 shift[ZZ] = -ref_com[ZZ];
1391 translate_x(ref_s_1, natoms, shift);
1393 shift[XX] = -act_com[XX];
1394 shift[YY] = -act_com[YY];
1395 shift[ZZ] = -act_com[ZZ];
1396 translate_x(act_s_1, natoms, shift);
1398 /* Align rotation axis with z */
1399 align_with_z(ref_s_1, natoms, axis);
1400 align_with_z(act_s_1, natoms, axis);
1402 /* Correlation matrix */
1403 Rmat = allocate_square_matrix(3);
1405 for (i = 0; i < natoms; i++)
1407 ref_s_1[i][2] = 0.0;
1408 act_s_1[i][2] = 0.0;
1411 /* Weight positions with sqrt(weight) */
1412 if (nullptr != weight)
1414 weigh_coords(ref_s_1, weight, natoms);
1415 weigh_coords(act_s_1, weight, natoms);
1418 /* Calculate correlation matrices R=YXt (X=ref_s; Y=act_s) */
1419 calc_correl_matrix(ref_s_1, act_s_1, Rmat, natoms);
1421 /* Calculate RtR */
1422 RtR = allocate_square_matrix(3);
1423 for (i = 0; i < 3; i++)
1425 for (j = 0; j < 3; j++)
1427 for (k = 0; k < 3; k++)
1429 RtR[i][j] += Rmat[k][i] * Rmat[k][j];
1433 /* Diagonalize RtR */
1434 snew(eigvec, 3);
1435 for (i = 0; i < 3; i++)
1437 snew(eigvec[i], 3);
1440 diagonalize_symmetric(RtR, eigvec, eigval);
1441 swap_col(eigvec, 0, 1);
1442 swap_col(eigvec, 1, 2);
1443 swap_val(eigval, 0, 1);
1444 swap_val(eigval, 1, 2);
1446 /* Calculate V */
1447 for (i = 0; i < 3; i++)
1449 for (j = 0; j < 3; j++)
1451 V[i][j] = 0.0;
1452 WS[i][j] = 0.0;
1456 for (i = 0; i < 2; i++)
1458 for (j = 0; j < 2; j++)
1460 WS[i][j] = eigvec[i][j] / std::sqrt(eigval[j]);
1464 for (i = 0; i < 3; i++)
1466 for (j = 0; j < 3; j++)
1468 for (k = 0; k < 3; k++)
1470 V[i][j] += Rmat[i][k]*WS[k][j];
1474 free_square_matrix(Rmat, 3);
1476 /* Calculate optimal rotation matrix */
1477 for (i = 0; i < 3; i++)
1479 for (j = 0; j < 3; j++)
1481 rot_matrix[i][j] = 0.0;
1485 for (i = 0; i < 3; i++)
1487 for (j = 0; j < 3; j++)
1489 for (k = 0; k < 3; k++)
1491 rot_matrix[i][j] += eigvec[i][k]*V[j][k];
1495 rot_matrix[2][2] = 1.0;
1497 /* In some cases abs(rot_matrix[0][0]) can be slighly larger
1498 * than unity due to numerical inacurracies. To be able to calculate
1499 * the acos function, we put these values back in range. */
1500 if (rot_matrix[0][0] > 1.0)
1502 rot_matrix[0][0] = 1.0;
1504 else if (rot_matrix[0][0] < -1.0)
1506 rot_matrix[0][0] = -1.0;
1509 /* Determine the optimal rotation angle: */
1510 opt_angle = (-1.0)*acos(rot_matrix[0][0])*180.0/M_PI;
1511 if (rot_matrix[0][1] < 0.0)
1513 opt_angle = (-1.0)*opt_angle;
1516 /* Give back some memory */
1517 free_square_matrix(RtR, 3);
1518 sfree(ref_s_1);
1519 sfree(act_s_1);
1520 for (i = 0; i < 3; i++)
1522 sfree(eigvec[i]);
1524 sfree(eigvec);
1526 return static_cast<real>(opt_angle);
1530 /* Determine angle of the group by RMSD fit to the reference */
1531 /* Not parallelized, call this routine only on the master */
1532 static real flex_fit_angle(gmx_enfrotgrp *erg)
1534 rvec *fitcoords = nullptr;
1535 rvec center; /* Center of positions passed to the fit routine */
1536 real fitangle; /* Angle of the rotation group derived by fitting */
1537 rvec coord;
1538 real scal;
1540 /* Get the center of the rotation group.
1541 * Note, again, erg->xc has been sorted in do_flexible */
1542 get_center(erg->xc, erg->mc_sorted, erg->rotg->nat, center);
1544 /* === Determine the optimal fit angle for the rotation group === */
1545 if (erg->rotg->eFittype == erotgFitNORM)
1547 /* Normalize every position to it's reference length */
1548 for (int i = 0; i < erg->rotg->nat; i++)
1550 /* Put the center of the positions into the origin */
1551 rvec_sub(erg->xc[i], center, coord);
1552 /* Determine the scaling factor for the length: */
1553 scal = erg->xc_ref_length[erg->xc_sortind[i]] / norm(coord);
1554 /* Get position, multiply with the scaling factor and save */
1555 svmul(scal, coord, erg->xc_norm[i]);
1557 fitcoords = erg->xc_norm;
1559 else
1561 fitcoords = erg->xc;
1563 /* From the point of view of the current positions, the reference has rotated
1564 * backwards. Since we output the angle relative to the fixed reference,
1565 * we need the minus sign. */
1566 fitangle = -opt_angle_analytic(erg->xc_ref_sorted, fitcoords, erg->mc_sorted,
1567 erg->rotg->nat, erg->xc_ref_center, center, erg->vec);
1569 return fitangle;
1573 /* Determine actual angle of each slab by RMSD fit to the reference */
1574 /* Not parallelized, call this routine only on the master */
1575 static void flex_fit_angle_perslab(
1576 gmx_enfrotgrp *erg,
1577 double t,
1578 real degangle,
1579 FILE *fp)
1581 rvec curr_x, ref_x;
1582 rvec act_center; /* Center of actual positions that are passed to the fit routine */
1583 rvec ref_center; /* Same for the reference positions */
1584 real fitangle; /* Angle of a slab derived from an RMSD fit to
1585 * the reference structure at t=0 */
1586 gmx_slabdata *sd;
1587 real OOm_av; /* 1/average_mass of a rotation group atom */
1588 real m_rel; /* Relative mass of a rotation group atom */
1591 /* Average mass of a rotation group atom: */
1592 OOm_av = erg->invmass*erg->rotg->nat;
1594 /**********************************/
1595 /* First collect the data we need */
1596 /**********************************/
1598 /* Collect the data for the individual slabs */
1599 for (int n = erg->slab_first; n <= erg->slab_last; n++)
1601 int slabIndex = n - erg->slab_first; /* slab index */
1602 sd = &(erg->slab_data[slabIndex]);
1603 sd->nat = erg->lastatom[slabIndex]-erg->firstatom[slabIndex]+1;
1604 int ind = 0;
1606 /* Loop over the relevant atoms in the slab */
1607 for (int l = erg->firstatom[slabIndex]; l <= erg->lastatom[slabIndex]; l++)
1609 /* Current position of this atom: x[ii][XX/YY/ZZ] */
1610 copy_rvec(erg->xc[l], curr_x);
1612 /* The (unrotated) reference position of this atom is copied to ref_x.
1613 * Beware, the xc coords have been sorted in do_flexible */
1614 copy_rvec(erg->xc_ref_sorted[l], ref_x);
1616 /* Save data for doing angular RMSD fit later */
1617 /* Save the current atom position */
1618 copy_rvec(curr_x, sd->x[ind]);
1619 /* Save the corresponding reference position */
1620 copy_rvec(ref_x, sd->ref[ind]);
1622 /* Maybe also mass-weighting was requested. If yes, additionally
1623 * multiply the weights with the relative mass of the atom. If not,
1624 * multiply with unity. */
1625 m_rel = erg->mc_sorted[l]*OOm_av;
1627 /* Save the weight for this atom in this slab */
1628 sd->weight[ind] = gaussian_weight(curr_x, erg, n) * m_rel;
1630 /* Next atom in this slab */
1631 ind++;
1635 /******************************/
1636 /* Now do the fit calculation */
1637 /******************************/
1639 fprintf(fp, "%12.3e%6d%12.3f", t, erg->groupIndex, degangle);
1641 /* === Now do RMSD fitting for each slab === */
1642 /* We require at least SLAB_MIN_ATOMS in a slab, such that the fit makes sense. */
1643 #define SLAB_MIN_ATOMS 4
1645 for (int n = erg->slab_first; n <= erg->slab_last; n++)
1647 int slabIndex = n - erg->slab_first; /* slab index */
1648 sd = &(erg->slab_data[slabIndex]);
1649 if (sd->nat >= SLAB_MIN_ATOMS)
1651 /* Get the center of the slabs reference and current positions */
1652 get_center(sd->ref, sd->weight, sd->nat, ref_center);
1653 get_center(sd->x, sd->weight, sd->nat, act_center);
1654 if (erg->rotg->eFittype == erotgFitNORM)
1656 /* Normalize every position to it's reference length
1657 * prior to performing the fit */
1658 for (int i = 0; i < sd->nat; i++) /* Center */
1660 rvec_dec(sd->ref[i], ref_center);
1661 rvec_dec(sd->x[i], act_center);
1662 /* Normalize x_i such that it gets the same length as ref_i */
1663 svmul( norm(sd->ref[i])/norm(sd->x[i]), sd->x[i], sd->x[i] );
1665 /* We already subtracted the centers */
1666 clear_rvec(ref_center);
1667 clear_rvec(act_center);
1669 fitangle = -opt_angle_analytic(sd->ref, sd->x, sd->weight, sd->nat,
1670 ref_center, act_center, erg->vec);
1671 fprintf(fp, "%6d%6d%12.3f", n, sd->nat, fitangle);
1674 fprintf(fp, "\n");
1676 #undef SLAB_MIN_ATOMS
1680 /* Shift x with is */
1681 static inline void shift_single_coord(const matrix box, rvec x, const ivec is)
1683 int tx, ty, tz;
1686 tx = is[XX];
1687 ty = is[YY];
1688 tz = is[ZZ];
1690 if (TRICLINIC(box))
1692 x[XX] += tx*box[XX][XX]+ty*box[YY][XX]+tz*box[ZZ][XX];
1693 x[YY] += ty*box[YY][YY]+tz*box[ZZ][YY];
1694 x[ZZ] += tz*box[ZZ][ZZ];
1696 else
1698 x[XX] += tx*box[XX][XX];
1699 x[YY] += ty*box[YY][YY];
1700 x[ZZ] += tz*box[ZZ][ZZ];
1705 /* Determine the 'home' slab of this atom which is the
1706 * slab with the highest Gaussian weight of all */
1707 static inline int get_homeslab(
1708 rvec curr_x, /* The position for which the home slab shall be determined */
1709 const rvec rotvec, /* The rotation vector */
1710 real slabdist) /* The slab distance */
1712 real dist;
1715 /* The distance of the atom to the coordinate center (where the
1716 * slab with index 0) is */
1717 dist = iprod(rotvec, curr_x);
1719 return gmx::roundToInt(dist / slabdist);
1723 /* For a local atom determine the relevant slabs, i.e. slabs in
1724 * which the gaussian is larger than min_gaussian
1726 static int get_single_atom_gaussians(
1727 rvec curr_x,
1728 gmx_enfrotgrp *erg)
1731 /* Determine the 'home' slab of this atom: */
1732 int homeslab = get_homeslab(curr_x, erg->vec, erg->rotg->slab_dist);
1734 /* First determine the weight in the atoms home slab: */
1735 real g = gaussian_weight(curr_x, erg, homeslab);
1736 int count = 0;
1737 erg->gn_atom[count] = g;
1738 erg->gn_slabind[count] = homeslab;
1739 count++;
1742 /* Determine the max slab */
1743 int slab = homeslab;
1744 while (g > erg->rotg->min_gaussian)
1746 slab++;
1747 g = gaussian_weight(curr_x, erg, slab);
1748 erg->gn_slabind[count] = slab;
1749 erg->gn_atom[count] = g;
1750 count++;
1752 count--;
1754 /* Determine the min slab */
1755 slab = homeslab;
1758 slab--;
1759 g = gaussian_weight(curr_x, erg, slab);
1760 erg->gn_slabind[count] = slab;
1761 erg->gn_atom[count] = g;
1762 count++;
1764 while (g > erg->rotg->min_gaussian);
1765 count--;
1767 return count;
1771 static void flex2_precalc_inner_sum(const gmx_enfrotgrp *erg)
1773 rvec xi; /* positions in the i-sum */
1774 rvec xcn, ycn; /* the current and the reference slab centers */
1775 real gaussian_xi;
1776 rvec yi0;
1777 rvec rin; /* Helper variables */
1778 real fac, fac2;
1779 rvec innersumvec;
1780 real OOpsii, OOpsiistar;
1781 real sin_rin; /* s_ii.r_ii */
1782 rvec s_in, tmpvec, tmpvec2;
1783 real mi, wi; /* Mass-weighting of the positions */
1784 real N_M; /* N/M */
1787 N_M = erg->rotg->nat * erg->invmass;
1789 /* Loop over all slabs that contain something */
1790 for (int n = erg->slab_first; n <= erg->slab_last; n++)
1792 int slabIndex = n - erg->slab_first; /* slab index */
1794 /* The current center of this slab is saved in xcn: */
1795 copy_rvec(erg->slab_center[slabIndex], xcn);
1796 /* ... and the reference center in ycn: */
1797 copy_rvec(erg->slab_center_ref[slabIndex+erg->slab_buffer], ycn);
1799 /*** D. Calculate the whole inner sum used for second and third sum */
1800 /* For slab n, we need to loop over all atoms i again. Since we sorted
1801 * the atoms with respect to the rotation vector, we know that it is sufficient
1802 * to calculate from firstatom to lastatom only. All other contributions will
1803 * be very small. */
1804 clear_rvec(innersumvec);
1805 for (int i = erg->firstatom[slabIndex]; i <= erg->lastatom[slabIndex]; i++)
1807 /* Coordinate xi of this atom */
1808 copy_rvec(erg->xc[i], xi);
1810 /* The i-weights */
1811 gaussian_xi = gaussian_weight(xi, erg, n);
1812 mi = erg->mc_sorted[i]; /* need the sorted mass here */
1813 wi = N_M*mi;
1815 /* Calculate rin */
1816 copy_rvec(erg->xc_ref_sorted[i], yi0); /* Reference position yi0 */
1817 rvec_sub(yi0, ycn, tmpvec2); /* tmpvec2 = yi0 - ycn */
1818 mvmul(erg->rotmat, tmpvec2, rin); /* rin = Omega.(yi0 - ycn) */
1820 /* Calculate psi_i* and sin */
1821 rvec_sub(xi, xcn, tmpvec2); /* tmpvec2 = xi - xcn */
1823 /* In rare cases, when an atom position coincides with a slab center
1824 * (tmpvec2 == 0) we cannot compute the vector product for s_in.
1825 * However, since the atom is located directly on the pivot, this
1826 * slab's contribution to the force on that atom will be zero
1827 * anyway. Therefore, we continue with the next atom. */
1828 if (gmx_numzero(norm(tmpvec2))) /* 0 == norm(xi - xcn) */
1830 continue;
1833 cprod(erg->vec, tmpvec2, tmpvec); /* tmpvec = v x (xi - xcn) */
1834 OOpsiistar = norm2(tmpvec)+erg->rotg->eps; /* OOpsii* = 1/psii* = |v x (xi-xcn)|^2 + eps */
1835 OOpsii = norm(tmpvec); /* OOpsii = 1 / psii = |v x (xi - xcn)| */
1837 /* * v x (xi - xcn) */
1838 unitv(tmpvec, s_in); /* sin = ---------------- */
1839 /* |v x (xi - xcn)| */
1841 sin_rin = iprod(s_in, rin); /* sin_rin = sin . rin */
1843 /* Now the whole sum */
1844 fac = OOpsii/OOpsiistar;
1845 svmul(fac, rin, tmpvec);
1846 fac2 = fac*fac*OOpsii;
1847 svmul(fac2*sin_rin, s_in, tmpvec2);
1848 rvec_dec(tmpvec, tmpvec2);
1850 svmul(wi*gaussian_xi*sin_rin, tmpvec, tmpvec2);
1852 rvec_inc(innersumvec, tmpvec2);
1853 } /* now we have the inner sum, used both for sum2 and sum3 */
1855 /* Save it to be used in do_flex2_lowlevel */
1856 copy_rvec(innersumvec, erg->slab_innersumvec[slabIndex]);
1857 } /* END of loop over slabs */
1861 static void flex_precalc_inner_sum(const gmx_enfrotgrp *erg)
1863 rvec xi; /* position */
1864 rvec xcn, ycn; /* the current and the reference slab centers */
1865 rvec qin, rin; /* q_i^n and r_i^n */
1866 real bin;
1867 rvec tmpvec;
1868 rvec innersumvec; /* Inner part of sum_n2 */
1869 real gaussian_xi; /* Gaussian weight gn(xi) */
1870 real mi, wi; /* Mass-weighting of the positions */
1871 real N_M; /* N/M */
1873 N_M = erg->rotg->nat * erg->invmass;
1875 /* Loop over all slabs that contain something */
1876 for (int n = erg->slab_first; n <= erg->slab_last; n++)
1878 int slabIndex = n - erg->slab_first; /* slab index */
1880 /* The current center of this slab is saved in xcn: */
1881 copy_rvec(erg->slab_center[slabIndex], xcn);
1882 /* ... and the reference center in ycn: */
1883 copy_rvec(erg->slab_center_ref[slabIndex+erg->slab_buffer], ycn);
1885 /* For slab n, we need to loop over all atoms i again. Since we sorted
1886 * the atoms with respect to the rotation vector, we know that it is sufficient
1887 * to calculate from firstatom to lastatom only. All other contributions will
1888 * be very small. */
1889 clear_rvec(innersumvec);
1890 for (int i = erg->firstatom[slabIndex]; i <= erg->lastatom[slabIndex]; i++)
1892 /* Coordinate xi of this atom */
1893 copy_rvec(erg->xc[i], xi);
1895 /* The i-weights */
1896 gaussian_xi = gaussian_weight(xi, erg, n);
1897 mi = erg->mc_sorted[i]; /* need the sorted mass here */
1898 wi = N_M*mi;
1900 /* Calculate rin and qin */
1901 rvec_sub(erg->xc_ref_sorted[i], ycn, tmpvec); /* tmpvec = yi0-ycn */
1903 /* In rare cases, when an atom position coincides with a slab center
1904 * (tmpvec == 0) we cannot compute the vector product for qin.
1905 * However, since the atom is located directly on the pivot, this
1906 * slab's contribution to the force on that atom will be zero
1907 * anyway. Therefore, we continue with the next atom. */
1908 if (gmx_numzero(norm(tmpvec))) /* 0 == norm(yi0 - ycn) */
1910 continue;
1913 mvmul(erg->rotmat, tmpvec, rin); /* rin = Omega.(yi0 - ycn) */
1914 cprod(erg->vec, rin, tmpvec); /* tmpvec = v x Omega*(yi0-ycn) */
1916 /* * v x Omega*(yi0-ycn) */
1917 unitv(tmpvec, qin); /* qin = --------------------- */
1918 /* |v x Omega*(yi0-ycn)| */
1920 /* Calculate bin */
1921 rvec_sub(xi, xcn, tmpvec); /* tmpvec = xi-xcn */
1922 bin = iprod(qin, tmpvec); /* bin = qin*(xi-xcn) */
1924 svmul(wi*gaussian_xi*bin, qin, tmpvec);
1926 /* Add this contribution to the inner sum: */
1927 rvec_add(innersumvec, tmpvec, innersumvec);
1928 } /* now we have the inner sum vector S^n for this slab */
1929 /* Save it to be used in do_flex_lowlevel */
1930 copy_rvec(innersumvec, erg->slab_innersumvec[slabIndex]);
1935 static real do_flex2_lowlevel(
1936 gmx_enfrotgrp *erg,
1937 real sigma, /* The Gaussian width sigma */
1938 rvec x[],
1939 gmx_bool bOutstepRot,
1940 gmx_bool bOutstepSlab,
1941 matrix box)
1943 int count, ii, iigrp;
1944 rvec xj; /* position in the i-sum */
1945 rvec yj0; /* the reference position in the j-sum */
1946 rvec xcn, ycn; /* the current and the reference slab centers */
1947 real V; /* This node's part of the rotation pot. energy */
1948 real gaussian_xj; /* Gaussian weight */
1949 real beta;
1951 real numerator, fit_numerator;
1952 rvec rjn, fit_rjn; /* Helper variables */
1953 real fac, fac2;
1955 real OOpsij, OOpsijstar;
1956 real OOsigma2; /* 1/(sigma^2) */
1957 real sjn_rjn;
1958 real betasigpsi;
1959 rvec sjn, tmpvec, tmpvec2, yj0_ycn;
1960 rvec sum1vec_part, sum1vec, sum2vec_part, sum2vec, sum3vec, sum4vec, innersumvec;
1961 real sum3, sum4;
1962 real mj, wj; /* Mass-weighting of the positions */
1963 real N_M; /* N/M */
1964 real Wjn; /* g_n(x_j) m_j / Mjn */
1965 gmx_bool bCalcPotFit;
1967 /* To calculate the torque per slab */
1968 rvec slab_force; /* Single force from slab n on one atom */
1969 rvec slab_sum1vec_part;
1970 real slab_sum3part, slab_sum4part;
1971 rvec slab_sum1vec, slab_sum2vec, slab_sum3vec, slab_sum4vec;
1973 /* Pre-calculate the inner sums, so that we do not have to calculate
1974 * them again for every atom */
1975 flex2_precalc_inner_sum(erg);
1977 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT == erg->rotg->eFittype);
1979 /********************************************************/
1980 /* Main loop over all local atoms of the rotation group */
1981 /********************************************************/
1982 N_M = erg->rotg->nat * erg->invmass;
1983 V = 0.0;
1984 OOsigma2 = 1.0 / (sigma*sigma);
1985 const auto &localRotationGroupIndex = erg->atomSet->localIndex();
1986 const auto &collectiveRotationGroupIndex = erg->atomSet->collectiveIndex();
1988 for (gmx::index j = 0; j < localRotationGroupIndex.ssize(); j++)
1990 /* Local index of a rotation group atom */
1991 ii = localRotationGroupIndex[j];
1992 /* Position of this atom in the collective array */
1993 iigrp = collectiveRotationGroupIndex[j];
1994 /* Mass-weighting */
1995 mj = erg->mc[iigrp]; /* need the unsorted mass here */
1996 wj = N_M*mj;
1998 /* Current position of this atom: x[ii][XX/YY/ZZ]
1999 * Note that erg->xc_center contains the center of mass in case the flex2-t
2000 * potential was chosen. For the flex2 potential erg->xc_center must be
2001 * zero. */
2002 rvec_sub(x[ii], erg->xc_center, xj);
2004 /* Shift this atom such that it is near its reference */
2005 shift_single_coord(box, xj, erg->xc_shifts[iigrp]);
2007 /* Determine the slabs to loop over, i.e. the ones with contributions
2008 * larger than min_gaussian */
2009 count = get_single_atom_gaussians(xj, erg);
2011 clear_rvec(sum1vec_part);
2012 clear_rvec(sum2vec_part);
2013 sum3 = 0.0;
2014 sum4 = 0.0;
2015 /* Loop over the relevant slabs for this atom */
2016 for (int ic = 0; ic < count; ic++)
2018 int n = erg->gn_slabind[ic];
2020 /* Get the precomputed Gaussian value of curr_slab for curr_x */
2021 gaussian_xj = erg->gn_atom[ic];
2023 int slabIndex = n - erg->slab_first; /* slab index */
2025 /* The (unrotated) reference position of this atom is copied to yj0: */
2026 copy_rvec(erg->rotg->x_ref[iigrp], yj0);
2028 beta = calc_beta(xj, erg, n);
2030 /* The current center of this slab is saved in xcn: */
2031 copy_rvec(erg->slab_center[slabIndex], xcn);
2032 /* ... and the reference center in ycn: */
2033 copy_rvec(erg->slab_center_ref[slabIndex+erg->slab_buffer], ycn);
2035 rvec_sub(yj0, ycn, yj0_ycn); /* yj0_ycn = yj0 - ycn */
2037 /* Rotate: */
2038 mvmul(erg->rotmat, yj0_ycn, rjn); /* rjn = Omega.(yj0 - ycn) */
2040 /* Subtract the slab center from xj */
2041 rvec_sub(xj, xcn, tmpvec2); /* tmpvec2 = xj - xcn */
2043 /* In rare cases, when an atom position coincides with a slab center
2044 * (tmpvec2 == 0) we cannot compute the vector product for sjn.
2045 * However, since the atom is located directly on the pivot, this
2046 * slab's contribution to the force on that atom will be zero
2047 * anyway. Therefore, we directly move on to the next slab. */
2048 if (gmx_numzero(norm(tmpvec2))) /* 0 == norm(xj - xcn) */
2050 continue;
2053 /* Calculate sjn */
2054 cprod(erg->vec, tmpvec2, tmpvec); /* tmpvec = v x (xj - xcn) */
2056 OOpsijstar = norm2(tmpvec)+erg->rotg->eps; /* OOpsij* = 1/psij* = |v x (xj-xcn)|^2 + eps */
2058 numerator = gmx::square(iprod(tmpvec, rjn));
2060 /*********************************/
2061 /* Add to the rotation potential */
2062 /*********************************/
2063 V += 0.5*erg->rotg->k*wj*gaussian_xj*numerator/OOpsijstar;
2065 /* If requested, also calculate the potential for a set of angles
2066 * near the current reference angle */
2067 if (bCalcPotFit)
2069 for (int ifit = 0; ifit < erg->rotg->PotAngle_nstep; ifit++)
2071 mvmul(erg->PotAngleFit->rotmat[ifit], yj0_ycn, fit_rjn);
2072 fit_numerator = gmx::square(iprod(tmpvec, fit_rjn));
2073 erg->PotAngleFit->V[ifit] += 0.5*erg->rotg->k*wj*gaussian_xj*fit_numerator/OOpsijstar;
2077 /*************************************/
2078 /* Now calculate the force on atom j */
2079 /*************************************/
2081 OOpsij = norm(tmpvec); /* OOpsij = 1 / psij = |v x (xj - xcn)| */
2083 /* * v x (xj - xcn) */
2084 unitv(tmpvec, sjn); /* sjn = ---------------- */
2085 /* |v x (xj - xcn)| */
2087 sjn_rjn = iprod(sjn, rjn); /* sjn_rjn = sjn . rjn */
2090 /*** A. Calculate the first of the four sum terms: ****************/
2091 fac = OOpsij/OOpsijstar;
2092 svmul(fac, rjn, tmpvec);
2093 fac2 = fac*fac*OOpsij;
2094 svmul(fac2*sjn_rjn, sjn, tmpvec2);
2095 rvec_dec(tmpvec, tmpvec2);
2096 fac2 = wj*gaussian_xj; /* also needed for sum4 */
2097 svmul(fac2*sjn_rjn, tmpvec, slab_sum1vec_part);
2098 /********************/
2099 /*** Add to sum1: ***/
2100 /********************/
2101 rvec_inc(sum1vec_part, slab_sum1vec_part); /* sum1 still needs to vector multiplied with v */
2103 /*** B. Calculate the forth of the four sum terms: ****************/
2104 betasigpsi = beta*OOsigma2*OOpsij; /* this is also needed for sum3 */
2105 /********************/
2106 /*** Add to sum4: ***/
2107 /********************/
2108 slab_sum4part = fac2*betasigpsi*fac*sjn_rjn*sjn_rjn; /* Note that fac is still valid from above */
2109 sum4 += slab_sum4part;
2111 /*** C. Calculate Wjn for second and third sum */
2112 /* Note that we can safely divide by slab_weights since we check in
2113 * get_slab_centers that it is non-zero. */
2114 Wjn = gaussian_xj*mj/erg->slab_weights[slabIndex];
2116 /* We already have precalculated the inner sum for slab n */
2117 copy_rvec(erg->slab_innersumvec[slabIndex], innersumvec);
2119 /* Weigh the inner sum vector with Wjn */
2120 svmul(Wjn, innersumvec, innersumvec);
2122 /*** E. Calculate the second of the four sum terms: */
2123 /********************/
2124 /*** Add to sum2: ***/
2125 /********************/
2126 rvec_inc(sum2vec_part, innersumvec); /* sum2 still needs to be vector crossproduct'ed with v */
2128 /*** F. Calculate the third of the four sum terms: */
2129 slab_sum3part = betasigpsi * iprod(sjn, innersumvec);
2130 sum3 += slab_sum3part; /* still needs to be multiplied with v */
2132 /*** G. Calculate the torque on the local slab's axis: */
2133 if (bOutstepRot)
2135 /* Sum1 */
2136 cprod(slab_sum1vec_part, erg->vec, slab_sum1vec);
2137 /* Sum2 */
2138 cprod(innersumvec, erg->vec, slab_sum2vec);
2139 /* Sum3 */
2140 svmul(slab_sum3part, erg->vec, slab_sum3vec);
2141 /* Sum4 */
2142 svmul(slab_sum4part, erg->vec, slab_sum4vec);
2144 /* The force on atom ii from slab n only: */
2145 for (int m = 0; m < DIM; m++)
2147 slab_force[m] = erg->rotg->k * (-slab_sum1vec[m] + slab_sum2vec[m] - slab_sum3vec[m] + 0.5*slab_sum4vec[m]);
2150 erg->slab_torque_v[slabIndex] += torque(erg->vec, slab_force, xj, xcn);
2152 } /* END of loop over slabs */
2154 /* Construct the four individual parts of the vector sum: */
2155 cprod(sum1vec_part, erg->vec, sum1vec); /* sum1vec = { } x v */
2156 cprod(sum2vec_part, erg->vec, sum2vec); /* sum2vec = { } x v */
2157 svmul(sum3, erg->vec, sum3vec); /* sum3vec = { } . v */
2158 svmul(sum4, erg->vec, sum4vec); /* sum4vec = { } . v */
2160 /* Store the additional force so that it can be added to the force
2161 * array after the normal forces have been evaluated */
2162 for (int m = 0; m < DIM; m++)
2164 erg->f_rot_loc[j][m] = erg->rotg->k * (-sum1vec[m] + sum2vec[m] - sum3vec[m] + 0.5*sum4vec[m]);
2167 #ifdef SUM_PARTS
2168 fprintf(stderr, "sum1: %15.8f %15.8f %15.8f\n", -erg->rotg->k*sum1vec[XX], -erg->rotg->k*sum1vec[YY], -erg->rotg->k*sum1vec[ZZ]);
2169 fprintf(stderr, "sum2: %15.8f %15.8f %15.8f\n", erg->rotg->k*sum2vec[XX], erg->rotg->k*sum2vec[YY], erg->rotg->k*sum2vec[ZZ]);
2170 fprintf(stderr, "sum3: %15.8f %15.8f %15.8f\n", -erg->rotg->k*sum3vec[XX], -erg->rotg->k*sum3vec[YY], -erg->rotg->k*sum3vec[ZZ]);
2171 fprintf(stderr, "sum4: %15.8f %15.8f %15.8f\n", 0.5*erg->rotg->k*sum4vec[XX], 0.5*erg->rotg->k*sum4vec[YY], 0.5*erg->rotg->k*sum4vec[ZZ]);
2172 #endif
2174 PRINT_FORCE_J
2176 } /* END of loop over local atoms */
2178 return V;
2182 static real do_flex_lowlevel(
2183 gmx_enfrotgrp *erg,
2184 real sigma, /* The Gaussian width sigma */
2185 rvec x[],
2186 gmx_bool bOutstepRot,
2187 gmx_bool bOutstepSlab,
2188 matrix box)
2190 int count, iigrp;
2191 rvec xj, yj0; /* current and reference position */
2192 rvec xcn, ycn; /* the current and the reference slab centers */
2193 rvec yj0_ycn; /* yj0 - ycn */
2194 rvec xj_xcn; /* xj - xcn */
2195 rvec qjn, fit_qjn; /* q_i^n */
2196 rvec sum_n1, sum_n2; /* Two contributions to the rotation force */
2197 rvec innersumvec; /* Inner part of sum_n2 */
2198 rvec s_n;
2199 rvec force_n; /* Single force from slab n on one atom */
2200 rvec force_n1, force_n2; /* First and second part of force_n */
2201 rvec tmpvec, tmpvec2, tmp_f; /* Helper variables */
2202 real V; /* The rotation potential energy */
2203 real OOsigma2; /* 1/(sigma^2) */
2204 real beta; /* beta_n(xj) */
2205 real bjn, fit_bjn; /* b_j^n */
2206 real gaussian_xj; /* Gaussian weight gn(xj) */
2207 real betan_xj_sigma2;
2208 real mj, wj; /* Mass-weighting of the positions */
2209 real N_M; /* N/M */
2210 gmx_bool bCalcPotFit;
2212 /* Pre-calculate the inner sums, so that we do not have to calculate
2213 * them again for every atom */
2214 flex_precalc_inner_sum(erg);
2216 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT == erg->rotg->eFittype);
2218 /********************************************************/
2219 /* Main loop over all local atoms of the rotation group */
2220 /********************************************************/
2221 OOsigma2 = 1.0/(sigma*sigma);
2222 N_M = erg->rotg->nat * erg->invmass;
2223 V = 0.0;
2224 const auto &localRotationGroupIndex = erg->atomSet->localIndex();
2225 const auto &collectiveRotationGroupIndex = erg->atomSet->collectiveIndex();
2227 for (gmx::index j = 0; j < localRotationGroupIndex.ssize(); j++)
2229 /* Local index of a rotation group atom */
2230 int ii = localRotationGroupIndex[j];
2231 /* Position of this atom in the collective array */
2232 iigrp = collectiveRotationGroupIndex[j];
2233 /* Mass-weighting */
2234 mj = erg->mc[iigrp]; /* need the unsorted mass here */
2235 wj = N_M*mj;
2237 /* Current position of this atom: x[ii][XX/YY/ZZ]
2238 * Note that erg->xc_center contains the center of mass in case the flex-t
2239 * potential was chosen. For the flex potential erg->xc_center must be
2240 * zero. */
2241 rvec_sub(x[ii], erg->xc_center, xj);
2243 /* Shift this atom such that it is near its reference */
2244 shift_single_coord(box, xj, erg->xc_shifts[iigrp]);
2246 /* Determine the slabs to loop over, i.e. the ones with contributions
2247 * larger than min_gaussian */
2248 count = get_single_atom_gaussians(xj, erg);
2250 clear_rvec(sum_n1);
2251 clear_rvec(sum_n2);
2253 /* Loop over the relevant slabs for this atom */
2254 for (int ic = 0; ic < count; ic++)
2256 int n = erg->gn_slabind[ic];
2258 /* Get the precomputed Gaussian for xj in slab n */
2259 gaussian_xj = erg->gn_atom[ic];
2261 int slabIndex = n - erg->slab_first; /* slab index */
2263 /* The (unrotated) reference position of this atom is saved in yj0: */
2264 copy_rvec(erg->rotg->x_ref[iigrp], yj0);
2266 beta = calc_beta(xj, erg, n);
2268 /* The current center of this slab is saved in xcn: */
2269 copy_rvec(erg->slab_center[slabIndex], xcn);
2270 /* ... and the reference center in ycn: */
2271 copy_rvec(erg->slab_center_ref[slabIndex+erg->slab_buffer], ycn);
2273 rvec_sub(yj0, ycn, yj0_ycn); /* yj0_ycn = yj0 - ycn */
2275 /* In rare cases, when an atom position coincides with a reference slab
2276 * center (yj0_ycn == 0) we cannot compute the normal vector qjn.
2277 * However, since the atom is located directly on the pivot, this
2278 * slab's contribution to the force on that atom will be zero
2279 * anyway. Therefore, we directly move on to the next slab. */
2280 if (gmx_numzero(norm(yj0_ycn))) /* 0 == norm(yj0 - ycn) */
2282 continue;
2285 /* Rotate: */
2286 mvmul(erg->rotmat, yj0_ycn, tmpvec2); /* tmpvec2= Omega.(yj0-ycn) */
2288 /* Subtract the slab center from xj */
2289 rvec_sub(xj, xcn, xj_xcn); /* xj_xcn = xj - xcn */
2291 /* Calculate qjn */
2292 cprod(erg->vec, tmpvec2, tmpvec); /* tmpvec= v x Omega.(yj0-ycn) */
2294 /* * v x Omega.(yj0-ycn) */
2295 unitv(tmpvec, qjn); /* qjn = --------------------- */
2296 /* |v x Omega.(yj0-ycn)| */
2298 bjn = iprod(qjn, xj_xcn); /* bjn = qjn * (xj - xcn) */
2300 /*********************************/
2301 /* Add to the rotation potential */
2302 /*********************************/
2303 V += 0.5*erg->rotg->k*wj*gaussian_xj*gmx::square(bjn);
2305 /* If requested, also calculate the potential for a set of angles
2306 * near the current reference angle */
2307 if (bCalcPotFit)
2309 for (int ifit = 0; ifit < erg->rotg->PotAngle_nstep; ifit++)
2311 /* As above calculate Omega.(yj0-ycn), now for the other angles */
2312 mvmul(erg->PotAngleFit->rotmat[ifit], yj0_ycn, tmpvec2); /* tmpvec2= Omega.(yj0-ycn) */
2313 /* As above calculate qjn */
2314 cprod(erg->vec, tmpvec2, tmpvec); /* tmpvec= v x Omega.(yj0-ycn) */
2315 /* * v x Omega.(yj0-ycn) */
2316 unitv(tmpvec, fit_qjn); /* fit_qjn = --------------------- */
2317 /* |v x Omega.(yj0-ycn)| */
2318 fit_bjn = iprod(fit_qjn, xj_xcn); /* fit_bjn = fit_qjn * (xj - xcn) */
2319 /* Add to the rotation potential for this angle */
2320 erg->PotAngleFit->V[ifit] += 0.5*erg->rotg->k*wj*gaussian_xj*gmx::square(fit_bjn);
2324 /****************************************************************/
2325 /* sum_n1 will typically be the main contribution to the force: */
2326 /****************************************************************/
2327 betan_xj_sigma2 = beta*OOsigma2; /* beta_n(xj)/sigma^2 */
2329 /* The next lines calculate
2330 * qjn - (bjn*beta(xj)/(2sigma^2))v */
2331 svmul(bjn*0.5*betan_xj_sigma2, erg->vec, tmpvec2);
2332 rvec_sub(qjn, tmpvec2, tmpvec);
2334 /* Multiply with gn(xj)*bjn: */
2335 svmul(gaussian_xj*bjn, tmpvec, tmpvec2);
2337 /* Sum over n: */
2338 rvec_inc(sum_n1, tmpvec2);
2340 /* We already have precalculated the Sn term for slab n */
2341 copy_rvec(erg->slab_innersumvec[slabIndex], s_n);
2342 /* * beta_n(xj) */
2343 svmul(betan_xj_sigma2*iprod(s_n, xj_xcn), erg->vec, tmpvec); /* tmpvec = ---------- s_n (xj-xcn) */
2344 /* sigma^2 */
2346 rvec_sub(s_n, tmpvec, innersumvec);
2348 /* We can safely divide by slab_weights since we check in get_slab_centers
2349 * that it is non-zero. */
2350 svmul(gaussian_xj/erg->slab_weights[slabIndex], innersumvec, innersumvec);
2352 rvec_add(sum_n2, innersumvec, sum_n2);
2354 /* Calculate the torque: */
2355 if (bOutstepRot)
2357 /* The force on atom ii from slab n only: */
2358 svmul(-erg->rotg->k*wj, tmpvec2, force_n1); /* part 1 */
2359 svmul( erg->rotg->k*mj, innersumvec, force_n2); /* part 2 */
2360 rvec_add(force_n1, force_n2, force_n);
2361 erg->slab_torque_v[slabIndex] += torque(erg->vec, force_n, xj, xcn);
2363 } /* END of loop over slabs */
2365 /* Put both contributions together: */
2366 svmul(wj, sum_n1, sum_n1);
2367 svmul(mj, sum_n2, sum_n2);
2368 rvec_sub(sum_n2, sum_n1, tmp_f); /* F = -grad V */
2370 /* Store the additional force so that it can be added to the force
2371 * array after the normal forces have been evaluated */
2372 for (int m = 0; m < DIM; m++)
2374 erg->f_rot_loc[j][m] = erg->rotg->k*tmp_f[m];
2377 PRINT_FORCE_J
2379 } /* END of loop over local atoms */
2381 return V;
2384 static void sort_collective_coordinates(
2385 gmx_enfrotgrp *erg,
2386 sort_along_vec_t *data) /* Buffer for sorting the positions */
2388 /* The projection of the position vector on the rotation vector is
2389 * the relevant value for sorting. Fill the 'data' structure */
2390 for (int i = 0; i < erg->rotg->nat; i++)
2392 data[i].xcproj = iprod(erg->xc[i], erg->vec); /* sort criterium */
2393 data[i].m = erg->mc[i];
2394 data[i].ind = i;
2395 copy_rvec(erg->xc[i], data[i].x );
2396 copy_rvec(erg->rotg->x_ref[i], data[i].x_ref);
2398 /* Sort the 'data' structure */
2399 std::sort(data, data+erg->rotg->nat,
2400 [](const sort_along_vec_t &a, const sort_along_vec_t &b)
2402 return a.xcproj < b.xcproj;
2405 /* Copy back the sorted values */
2406 for (int i = 0; i < erg->rotg->nat; i++)
2408 copy_rvec(data[i].x, erg->xc[i] );
2409 copy_rvec(data[i].x_ref, erg->xc_ref_sorted[i]);
2410 erg->mc_sorted[i] = data[i].m;
2411 erg->xc_sortind[i] = data[i].ind;
2416 /* For each slab, get the first and the last index of the sorted atom
2417 * indices */
2418 static void get_firstlast_atom_per_slab(const gmx_enfrotgrp *erg)
2420 real beta;
2422 /* Find the first atom that needs to enter the calculation for each slab */
2423 int n = erg->slab_first; /* slab */
2424 int i = 0; /* start with the first atom */
2427 /* Find the first atom that significantly contributes to this slab */
2428 do /* move forward in position until a large enough beta is found */
2430 beta = calc_beta(erg->xc[i], erg, n);
2431 i++;
2433 while ((beta < -erg->max_beta) && (i < erg->rotg->nat));
2434 i--;
2435 int slabIndex = n - erg->slab_first; /* slab index */
2436 erg->firstatom[slabIndex] = i;
2437 /* Proceed to the next slab */
2438 n++;
2440 while (n <= erg->slab_last);
2442 /* Find the last atom for each slab */
2443 n = erg->slab_last; /* start with last slab */
2444 i = erg->rotg->nat-1; /* start with the last atom */
2447 do /* move backward in position until a large enough beta is found */
2449 beta = calc_beta(erg->xc[i], erg, n);
2450 i--;
2452 while ((beta > erg->max_beta) && (i > -1));
2453 i++;
2454 int slabIndex = n - erg->slab_first; /* slab index */
2455 erg->lastatom[slabIndex] = i;
2456 /* Proceed to the next slab */
2457 n--;
2459 while (n >= erg->slab_first);
2463 /* Determine the very first and very last slab that needs to be considered
2464 * For the first slab that needs to be considered, we have to find the smallest
2465 * n that obeys:
2467 * x_first * v - n*Delta_x <= beta_max
2469 * slab index n, slab distance Delta_x, rotation vector v. For the last slab we
2470 * have to find the largest n that obeys
2472 * x_last * v - n*Delta_x >= -beta_max
2475 static inline int get_first_slab(
2476 const gmx_enfrotgrp *erg,
2477 rvec firstatom) /* First atom after sorting along the rotation vector v */
2479 /* Find the first slab for the first atom */
2480 return static_cast<int>(ceil(static_cast<double>((iprod(firstatom, erg->vec) - erg->max_beta)/erg->rotg->slab_dist)));
2484 static inline int get_last_slab(
2485 const gmx_enfrotgrp *erg,
2486 rvec lastatom) /* Last atom along v */
2488 /* Find the last slab for the last atom */
2489 return static_cast<int>(floor(static_cast<double>((iprod(lastatom, erg->vec) + erg->max_beta)/erg->rotg->slab_dist)));
2493 static void get_firstlast_slab_check(
2494 gmx_enfrotgrp *erg, /* The rotation group (data only accessible in this file) */
2495 rvec firstatom, /* First atom after sorting along the rotation vector v */
2496 rvec lastatom) /* Last atom along v */
2498 erg->slab_first = get_first_slab(erg, firstatom);
2499 erg->slab_last = get_last_slab(erg, lastatom);
2501 /* Calculate the slab buffer size, which changes when slab_first changes */
2502 erg->slab_buffer = erg->slab_first - erg->slab_first_ref;
2504 /* Check whether we have reference data to compare against */
2505 if (erg->slab_first < erg->slab_first_ref)
2507 gmx_fatal(FARGS, "%s No reference data for first slab (n=%d), unable to proceed.",
2508 RotStr, erg->slab_first);
2511 /* Check whether we have reference data to compare against */
2512 if (erg->slab_last > erg->slab_last_ref)
2514 gmx_fatal(FARGS, "%s No reference data for last slab (n=%d), unable to proceed.",
2515 RotStr, erg->slab_last);
2520 /* Enforced rotation with a flexible axis */
2521 static void do_flexible(
2522 gmx_bool bMaster,
2523 gmx_enfrot *enfrot, /* Other rotation data */
2524 gmx_enfrotgrp *erg,
2525 rvec x[], /* The local positions */
2526 matrix box,
2527 double t, /* Time in picoseconds */
2528 gmx_bool bOutstepRot, /* Output to main rotation output file */
2529 gmx_bool bOutstepSlab) /* Output per-slab data */
2531 int nslabs;
2532 real sigma; /* The Gaussian width sigma */
2534 /* Define the sigma value */
2535 sigma = 0.7*erg->rotg->slab_dist;
2537 /* Sort the collective coordinates erg->xc along the rotation vector. This is
2538 * an optimization for the inner loop. */
2539 sort_collective_coordinates(erg, enfrot->data);
2541 /* Determine the first relevant slab for the first atom and the last
2542 * relevant slab for the last atom */
2543 get_firstlast_slab_check(erg, erg->xc[0], erg->xc[erg->rotg->nat-1]);
2545 /* Determine for each slab depending on the min_gaussian cutoff criterium,
2546 * a first and a last atom index inbetween stuff needs to be calculated */
2547 get_firstlast_atom_per_slab(erg);
2549 /* Determine the gaussian-weighted center of positions for all slabs */
2550 get_slab_centers(erg, erg->xc, erg->mc_sorted, t, enfrot->out_slabs, bOutstepSlab, FALSE);
2552 /* Clear the torque per slab from last time step: */
2553 nslabs = erg->slab_last - erg->slab_first + 1;
2554 for (int l = 0; l < nslabs; l++)
2556 erg->slab_torque_v[l] = 0.0;
2559 /* Call the rotational forces kernel */
2560 if (erg->rotg->eType == erotgFLEX || erg->rotg->eType == erotgFLEXT)
2562 erg->V = do_flex_lowlevel(erg, sigma, x, bOutstepRot, bOutstepSlab, box);
2564 else if (erg->rotg->eType == erotgFLEX2 || erg->rotg->eType == erotgFLEX2T)
2566 erg->V = do_flex2_lowlevel(erg, sigma, x, bOutstepRot, bOutstepSlab, box);
2568 else
2570 gmx_fatal(FARGS, "Unknown flexible rotation type");
2573 /* Determine angle by RMSD fit to the reference - Let's hope this */
2574 /* only happens once in a while, since this is not parallelized! */
2575 if (bMaster && (erotgFitPOT != erg->rotg->eFittype) )
2577 if (bOutstepRot)
2579 /* Fit angle of the whole rotation group */
2580 erg->angle_v = flex_fit_angle(erg);
2582 if (bOutstepSlab)
2584 /* Fit angle of each slab */
2585 flex_fit_angle_perslab(erg, t, erg->degangle, enfrot->out_angles);
2589 /* Lump together the torques from all slabs: */
2590 erg->torque_v = 0.0;
2591 for (int l = 0; l < nslabs; l++)
2593 erg->torque_v += erg->slab_torque_v[l];
2598 /* Calculate the angle between reference and actual rotation group atom,
2599 * both projected into a plane perpendicular to the rotation vector: */
2600 static void angle(const gmx_enfrotgrp *erg,
2601 rvec x_act,
2602 rvec x_ref,
2603 real *alpha,
2604 real *weight) /* atoms near the rotation axis should count less than atoms far away */
2606 rvec xp, xrp; /* current and reference positions projected on a plane perpendicular to pg->vec */
2607 rvec dum;
2610 /* Project x_ref and x into a plane through the origin perpendicular to rot_vec: */
2611 /* Project x_ref: xrp = x_ref - (vec * x_ref) * vec */
2612 svmul(iprod(erg->vec, x_ref), erg->vec, dum);
2613 rvec_sub(x_ref, dum, xrp);
2614 /* Project x_act: */
2615 svmul(iprod(erg->vec, x_act), erg->vec, dum);
2616 rvec_sub(x_act, dum, xp);
2618 /* Retrieve information about which vector precedes. gmx_angle always
2619 * returns a positive angle. */
2620 cprod(xp, xrp, dum); /* if reference precedes, this is pointing into the same direction as vec */
2622 if (iprod(erg->vec, dum) >= 0)
2624 *alpha = -gmx_angle(xrp, xp);
2626 else
2628 *alpha = +gmx_angle(xrp, xp);
2631 /* Also return the weight */
2632 *weight = norm(xp);
2636 /* Project first vector onto a plane perpendicular to the second vector
2637 * dr = dr - (dr.v)v
2638 * Note that v must be of unit length.
2640 static inline void project_onto_plane(rvec dr, const rvec v)
2642 rvec tmp;
2645 svmul(iprod(dr, v), v, tmp); /* tmp = (dr.v)v */
2646 rvec_dec(dr, tmp); /* dr = dr - (dr.v)v */
2650 /* Fixed rotation: The rotation reference group rotates around the v axis. */
2651 /* The atoms of the actual rotation group are attached with imaginary */
2652 /* springs to the reference atoms. */
2653 static void do_fixed(
2654 gmx_enfrotgrp *erg,
2655 gmx_bool bOutstepRot, /* Output to main rotation output file */
2656 gmx_bool bOutstepSlab) /* Output per-slab data */
2658 rvec dr;
2659 rvec tmp_f; /* Force */
2660 real alpha; /* a single angle between an actual and a reference position */
2661 real weight; /* single weight for a single angle */
2662 rvec xi_xc; /* xi - xc */
2663 gmx_bool bCalcPotFit;
2664 rvec fit_xr_loc;
2666 /* for mass weighting: */
2667 real wi; /* Mass-weighting of the positions */
2668 real N_M; /* N/M */
2669 real k_wi; /* k times wi */
2671 gmx_bool bProject;
2673 bProject = (erg->rotg->eType == erotgPM) || (erg->rotg->eType == erotgPMPF);
2674 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT == erg->rotg->eFittype);
2676 N_M = erg->rotg->nat * erg->invmass;
2677 const auto &collectiveRotationGroupIndex = erg->atomSet->collectiveIndex();
2678 /* Each process calculates the forces on its local atoms */
2679 for (size_t j = 0; j < erg->atomSet->numAtomsLocal(); j++)
2681 /* Calculate (x_i-x_c) resp. (x_i-u) */
2682 rvec_sub(erg->x_loc_pbc[j], erg->xc_center, xi_xc);
2684 /* Calculate Omega*(y_i-y_c)-(x_i-x_c) */
2685 rvec_sub(erg->xr_loc[j], xi_xc, dr);
2687 if (bProject)
2689 project_onto_plane(dr, erg->vec);
2692 /* Mass-weighting */
2693 wi = N_M*erg->m_loc[j];
2695 /* Store the additional force so that it can be added to the force
2696 * array after the normal forces have been evaluated */
2697 k_wi = erg->rotg->k*wi;
2698 for (int m = 0; m < DIM; m++)
2700 tmp_f[m] = k_wi*dr[m];
2701 erg->f_rot_loc[j][m] = tmp_f[m];
2702 erg->V += 0.5*k_wi*gmx::square(dr[m]);
2705 /* If requested, also calculate the potential for a set of angles
2706 * near the current reference angle */
2707 if (bCalcPotFit)
2709 for (int ifit = 0; ifit < erg->rotg->PotAngle_nstep; ifit++)
2711 /* Index of this rotation group atom with respect to the whole rotation group */
2712 int jj = collectiveRotationGroupIndex[j];
2714 /* Rotate with the alternative angle. Like rotate_local_reference(),
2715 * just for a single local atom */
2716 mvmul(erg->PotAngleFit->rotmat[ifit], erg->rotg->x_ref[jj], fit_xr_loc); /* fit_xr_loc = Omega*(y_i-y_c) */
2718 /* Calculate Omega*(y_i-y_c)-(x_i-x_c) */
2719 rvec_sub(fit_xr_loc, xi_xc, dr);
2721 if (bProject)
2723 project_onto_plane(dr, erg->vec);
2726 /* Add to the rotation potential for this angle: */
2727 erg->PotAngleFit->V[ifit] += 0.5*k_wi*norm2(dr);
2731 if (bOutstepRot)
2733 /* Add to the torque of this rotation group */
2734 erg->torque_v += torque(erg->vec, tmp_f, erg->x_loc_pbc[j], erg->xc_center);
2736 /* Calculate the angle between reference and actual rotation group atom. */
2737 angle(erg, xi_xc, erg->xr_loc[j], &alpha, &weight); /* angle in rad, weighted */
2738 erg->angle_v += alpha * weight;
2739 erg->weight_v += weight;
2741 /* If you want enforced rotation to contribute to the virial,
2742 * activate the following lines:
2743 if (MASTER(cr))
2745 Add the rotation contribution to the virial
2746 for(j=0; j<DIM; j++)
2747 for(m=0;m<DIM;m++)
2748 vir[j][m] += 0.5*f[ii][j]*dr[m];
2752 PRINT_FORCE_J
2754 } /* end of loop over local rotation group atoms */
2758 /* Calculate the radial motion potential and forces */
2759 static void do_radial_motion(
2760 gmx_enfrotgrp *erg,
2761 gmx_bool bOutstepRot, /* Output to main rotation output file */
2762 gmx_bool bOutstepSlab) /* Output per-slab data */
2764 rvec tmp_f; /* Force */
2765 real alpha; /* a single angle between an actual and a reference position */
2766 real weight; /* single weight for a single angle */
2767 rvec xj_u; /* xj - u */
2768 rvec tmpvec, fit_tmpvec;
2769 real fac, fac2, sum = 0.0;
2770 rvec pj;
2771 gmx_bool bCalcPotFit;
2773 /* For mass weighting: */
2774 real wj; /* Mass-weighting of the positions */
2775 real N_M; /* N/M */
2777 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT == erg->rotg->eFittype);
2779 N_M = erg->rotg->nat * erg->invmass;
2780 const auto &collectiveRotationGroupIndex = erg->atomSet->collectiveIndex();
2781 /* Each process calculates the forces on its local atoms */
2782 for (size_t j = 0; j < erg->atomSet->numAtomsLocal(); j++)
2784 /* Calculate (xj-u) */
2785 rvec_sub(erg->x_loc_pbc[j], erg->xc_center, xj_u); /* xj_u = xj-u */
2787 /* Calculate Omega.(yj0-u) */
2788 cprod(erg->vec, erg->xr_loc[j], tmpvec); /* tmpvec = v x Omega.(yj0-u) */
2790 /* * v x Omega.(yj0-u) */
2791 unitv(tmpvec, pj); /* pj = --------------------- */
2792 /* | v x Omega.(yj0-u) | */
2794 fac = iprod(pj, xj_u); /* fac = pj.(xj-u) */
2795 fac2 = fac*fac;
2797 /* Mass-weighting */
2798 wj = N_M*erg->m_loc[j];
2800 /* Store the additional force so that it can be added to the force
2801 * array after the normal forces have been evaluated */
2802 svmul(-erg->rotg->k*wj*fac, pj, tmp_f);
2803 copy_rvec(tmp_f, erg->f_rot_loc[j]);
2804 sum += wj*fac2;
2806 /* If requested, also calculate the potential for a set of angles
2807 * near the current reference angle */
2808 if (bCalcPotFit)
2810 for (int ifit = 0; ifit < erg->rotg->PotAngle_nstep; ifit++)
2812 /* Index of this rotation group atom with respect to the whole rotation group */
2813 int jj = collectiveRotationGroupIndex[j];
2815 /* Rotate with the alternative angle. Like rotate_local_reference(),
2816 * just for a single local atom */
2817 mvmul(erg->PotAngleFit->rotmat[ifit], erg->rotg->x_ref[jj], fit_tmpvec); /* fit_tmpvec = Omega*(yj0-u) */
2819 /* Calculate Omega.(yj0-u) */
2820 cprod(erg->vec, fit_tmpvec, tmpvec); /* tmpvec = v x Omega.(yj0-u) */
2821 /* * v x Omega.(yj0-u) */
2822 unitv(tmpvec, pj); /* pj = --------------------- */
2823 /* | v x Omega.(yj0-u) | */
2825 fac = iprod(pj, xj_u); /* fac = pj.(xj-u) */
2826 fac2 = fac*fac;
2828 /* Add to the rotation potential for this angle: */
2829 erg->PotAngleFit->V[ifit] += 0.5*erg->rotg->k*wj*fac2;
2833 if (bOutstepRot)
2835 /* Add to the torque of this rotation group */
2836 erg->torque_v += torque(erg->vec, tmp_f, erg->x_loc_pbc[j], erg->xc_center);
2838 /* Calculate the angle between reference and actual rotation group atom. */
2839 angle(erg, xj_u, erg->xr_loc[j], &alpha, &weight); /* angle in rad, weighted */
2840 erg->angle_v += alpha * weight;
2841 erg->weight_v += weight;
2844 PRINT_FORCE_J
2846 } /* end of loop over local rotation group atoms */
2847 erg->V = 0.5*erg->rotg->k*sum;
2851 /* Calculate the radial motion pivot-free potential and forces */
2852 static void do_radial_motion_pf(
2853 gmx_enfrotgrp *erg,
2854 rvec x[], /* The positions */
2855 matrix box, /* The simulation box */
2856 gmx_bool bOutstepRot, /* Output to main rotation output file */
2857 gmx_bool bOutstepSlab) /* Output per-slab data */
2859 rvec xj; /* Current position */
2860 rvec xj_xc; /* xj - xc */
2861 rvec yj0_yc0; /* yj0 - yc0 */
2862 rvec tmp_f; /* Force */
2863 real alpha; /* a single angle between an actual and a reference position */
2864 real weight; /* single weight for a single angle */
2865 rvec tmpvec, tmpvec2;
2866 rvec innersumvec; /* Precalculation of the inner sum */
2867 rvec innersumveckM;
2868 real fac, fac2, V = 0.0;
2869 rvec qi, qj;
2870 gmx_bool bCalcPotFit;
2872 /* For mass weighting: */
2873 real mj, wi, wj; /* Mass-weighting of the positions */
2874 real N_M; /* N/M */
2876 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT == erg->rotg->eFittype);
2878 N_M = erg->rotg->nat * erg->invmass;
2880 /* Get the current center of the rotation group: */
2881 get_center(erg->xc, erg->mc, erg->rotg->nat, erg->xc_center);
2883 /* Precalculate Sum_i [ wi qi.(xi-xc) qi ] which is needed for every single j */
2884 clear_rvec(innersumvec);
2885 for (int i = 0; i < erg->rotg->nat; i++)
2887 /* Mass-weighting */
2888 wi = N_M*erg->mc[i];
2890 /* Calculate qi. Note that xc_ref_center has already been subtracted from
2891 * x_ref in init_rot_group.*/
2892 mvmul(erg->rotmat, erg->rotg->x_ref[i], tmpvec); /* tmpvec = Omega.(yi0-yc0) */
2894 cprod(erg->vec, tmpvec, tmpvec2); /* tmpvec2 = v x Omega.(yi0-yc0) */
2896 /* * v x Omega.(yi0-yc0) */
2897 unitv(tmpvec2, qi); /* qi = ----------------------- */
2898 /* | v x Omega.(yi0-yc0) | */
2900 rvec_sub(erg->xc[i], erg->xc_center, tmpvec); /* tmpvec = xi-xc */
2902 svmul(wi*iprod(qi, tmpvec), qi, tmpvec2);
2904 rvec_inc(innersumvec, tmpvec2);
2906 svmul(erg->rotg->k*erg->invmass, innersumvec, innersumveckM);
2908 /* Each process calculates the forces on its local atoms */
2909 const auto &localRotationGroupIndex = erg->atomSet->localIndex();
2910 const auto &collectiveRotationGroupIndex = erg->atomSet->collectiveIndex();
2911 for (gmx::index j = 0; j < localRotationGroupIndex.ssize(); j++)
2913 /* Local index of a rotation group atom */
2914 int ii = localRotationGroupIndex[j];
2915 /* Position of this atom in the collective array */
2916 int iigrp = collectiveRotationGroupIndex[j];
2917 /* Mass-weighting */
2918 mj = erg->mc[iigrp]; /* need the unsorted mass here */
2919 wj = N_M*mj;
2921 /* Current position of this atom: x[ii][XX/YY/ZZ] */
2922 copy_rvec(x[ii], xj);
2924 /* Shift this atom such that it is near its reference */
2925 shift_single_coord(box, xj, erg->xc_shifts[iigrp]);
2927 /* The (unrotated) reference position is yj0. yc0 has already
2928 * been subtracted in init_rot_group */
2929 copy_rvec(erg->rotg->x_ref[iigrp], yj0_yc0); /* yj0_yc0 = yj0 - yc0 */
2931 /* Calculate Omega.(yj0-yc0) */
2932 mvmul(erg->rotmat, yj0_yc0, tmpvec2); /* tmpvec2 = Omega.(yj0 - yc0) */
2934 cprod(erg->vec, tmpvec2, tmpvec); /* tmpvec = v x Omega.(yj0-yc0) */
2936 /* * v x Omega.(yj0-yc0) */
2937 unitv(tmpvec, qj); /* qj = ----------------------- */
2938 /* | v x Omega.(yj0-yc0) | */
2940 /* Calculate (xj-xc) */
2941 rvec_sub(xj, erg->xc_center, xj_xc); /* xj_xc = xj-xc */
2943 fac = iprod(qj, xj_xc); /* fac = qj.(xj-xc) */
2944 fac2 = fac*fac;
2946 /* Store the additional force so that it can be added to the force
2947 * array after the normal forces have been evaluated */
2948 svmul(-erg->rotg->k*wj*fac, qj, tmp_f); /* part 1 of force */
2949 svmul(mj, innersumveckM, tmpvec); /* part 2 of force */
2950 rvec_inc(tmp_f, tmpvec);
2951 copy_rvec(tmp_f, erg->f_rot_loc[j]);
2952 V += wj*fac2;
2954 /* If requested, also calculate the potential for a set of angles
2955 * near the current reference angle */
2956 if (bCalcPotFit)
2958 for (int ifit = 0; ifit < erg->rotg->PotAngle_nstep; ifit++)
2960 /* Rotate with the alternative angle. Like rotate_local_reference(),
2961 * just for a single local atom */
2962 mvmul(erg->PotAngleFit->rotmat[ifit], yj0_yc0, tmpvec2); /* tmpvec2 = Omega*(yj0-yc0) */
2964 /* Calculate Omega.(yj0-u) */
2965 cprod(erg->vec, tmpvec2, tmpvec); /* tmpvec = v x Omega.(yj0-yc0) */
2966 /* * v x Omega.(yj0-yc0) */
2967 unitv(tmpvec, qj); /* qj = ----------------------- */
2968 /* | v x Omega.(yj0-yc0) | */
2970 fac = iprod(qj, xj_xc); /* fac = qj.(xj-xc) */
2971 fac2 = fac*fac;
2973 /* Add to the rotation potential for this angle: */
2974 erg->PotAngleFit->V[ifit] += 0.5*erg->rotg->k*wj*fac2;
2978 if (bOutstepRot)
2980 /* Add to the torque of this rotation group */
2981 erg->torque_v += torque(erg->vec, tmp_f, xj, erg->xc_center);
2983 /* Calculate the angle between reference and actual rotation group atom. */
2984 angle(erg, xj_xc, yj0_yc0, &alpha, &weight); /* angle in rad, weighted */
2985 erg->angle_v += alpha * weight;
2986 erg->weight_v += weight;
2989 PRINT_FORCE_J
2991 } /* end of loop over local rotation group atoms */
2992 erg->V = 0.5*erg->rotg->k*V;
2996 /* Precalculate the inner sum for the radial motion 2 forces */
2997 static void radial_motion2_precalc_inner_sum(const gmx_enfrotgrp *erg,
2998 rvec innersumvec)
3000 int i;
3001 rvec xi_xc; /* xj - xc */
3002 rvec tmpvec, tmpvec2;
3003 real fac;
3004 rvec ri, si;
3005 real siri;
3006 rvec v_xi_xc; /* v x (xj - u) */
3007 real psii, psiistar;
3008 real wi; /* Mass-weighting of the positions */
3009 real N_M; /* N/M */
3010 rvec sumvec;
3012 N_M = erg->rotg->nat * erg->invmass;
3014 /* Loop over the collective set of positions */
3015 clear_rvec(sumvec);
3016 for (i = 0; i < erg->rotg->nat; i++)
3018 /* Mass-weighting */
3019 wi = N_M*erg->mc[i];
3021 rvec_sub(erg->xc[i], erg->xc_center, xi_xc); /* xi_xc = xi-xc */
3023 /* Calculate ri. Note that xc_ref_center has already been subtracted from
3024 * x_ref in init_rot_group.*/
3025 mvmul(erg->rotmat, erg->rotg->x_ref[i], ri); /* ri = Omega.(yi0-yc0) */
3027 cprod(erg->vec, xi_xc, v_xi_xc); /* v_xi_xc = v x (xi-u) */
3029 fac = norm2(v_xi_xc);
3030 /* * 1 */
3031 psiistar = 1.0/(fac + erg->rotg->eps); /* psiistar = --------------------- */
3032 /* |v x (xi-xc)|^2 + eps */
3034 psii = gmx::invsqrt(fac); /* 1 */
3035 /* psii = ------------- */
3036 /* |v x (xi-xc)| */
3038 svmul(psii, v_xi_xc, si); /* si = psii * (v x (xi-xc) ) */
3040 siri = iprod(si, ri); /* siri = si.ri */
3042 svmul(psiistar/psii, ri, tmpvec);
3043 svmul(psiistar*psiistar/(psii*psii*psii) * siri, si, tmpvec2);
3044 rvec_dec(tmpvec, tmpvec2);
3045 cprod(tmpvec, erg->vec, tmpvec2);
3047 svmul(wi*siri, tmpvec2, tmpvec);
3049 rvec_inc(sumvec, tmpvec);
3051 svmul(erg->rotg->k*erg->invmass, sumvec, innersumvec);
3055 /* Calculate the radial motion 2 potential and forces */
3056 static void do_radial_motion2(
3057 gmx_enfrotgrp *erg,
3058 rvec x[], /* The positions */
3059 matrix box, /* The simulation box */
3060 gmx_bool bOutstepRot, /* Output to main rotation output file */
3061 gmx_bool bOutstepSlab) /* Output per-slab data */
3063 rvec xj; /* Position */
3064 real alpha; /* a single angle between an actual and a reference position */
3065 real weight; /* single weight for a single angle */
3066 rvec xj_u; /* xj - u */
3067 rvec yj0_yc0; /* yj0 -yc0 */
3068 rvec tmpvec, tmpvec2;
3069 real fac, fit_fac, fac2, Vpart = 0.0;
3070 rvec rj, fit_rj, sj;
3071 real sjrj;
3072 rvec v_xj_u; /* v x (xj - u) */
3073 real psij, psijstar;
3074 real mj, wj; /* For mass-weighting of the positions */
3075 real N_M; /* N/M */
3076 gmx_bool bPF;
3077 rvec innersumvec;
3078 gmx_bool bCalcPotFit;
3080 bPF = erg->rotg->eType == erotgRM2PF;
3081 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT == erg->rotg->eFittype);
3083 clear_rvec(yj0_yc0); /* Make the compiler happy */
3085 clear_rvec(innersumvec);
3086 if (bPF)
3088 /* For the pivot-free variant we have to use the current center of
3089 * mass of the rotation group instead of the pivot u */
3090 get_center(erg->xc, erg->mc, erg->rotg->nat, erg->xc_center);
3092 /* Also, we precalculate the second term of the forces that is identical
3093 * (up to the weight factor mj) for all forces */
3094 radial_motion2_precalc_inner_sum(erg, innersumvec);
3097 N_M = erg->rotg->nat * erg->invmass;
3099 /* Each process calculates the forces on its local atoms */
3100 const auto &localRotationGroupIndex = erg->atomSet->localIndex();
3101 const auto &collectiveRotationGroupIndex = erg->atomSet->collectiveIndex();
3102 for (gmx::index j = 0; j < localRotationGroupIndex.ssize(); j++)
3104 if (bPF)
3106 /* Local index of a rotation group atom */
3107 int ii = localRotationGroupIndex[j];
3108 /* Position of this atom in the collective array */
3109 int iigrp = collectiveRotationGroupIndex[j];
3110 /* Mass-weighting */
3111 mj = erg->mc[iigrp];
3113 /* Current position of this atom: x[ii] */
3114 copy_rvec(x[ii], xj);
3116 /* Shift this atom such that it is near its reference */
3117 shift_single_coord(box, xj, erg->xc_shifts[iigrp]);
3119 /* The (unrotated) reference position is yj0. yc0 has already
3120 * been subtracted in init_rot_group */
3121 copy_rvec(erg->rotg->x_ref[iigrp], yj0_yc0); /* yj0_yc0 = yj0 - yc0 */
3123 /* Calculate Omega.(yj0-yc0) */
3124 mvmul(erg->rotmat, yj0_yc0, rj); /* rj = Omega.(yj0-yc0) */
3126 else
3128 mj = erg->m_loc[j];
3129 copy_rvec(erg->x_loc_pbc[j], xj);
3130 copy_rvec(erg->xr_loc[j], rj); /* rj = Omega.(yj0-u) */
3132 /* Mass-weighting */
3133 wj = N_M*mj;
3135 /* Calculate (xj-u) resp. (xj-xc) */
3136 rvec_sub(xj, erg->xc_center, xj_u); /* xj_u = xj-u */
3138 cprod(erg->vec, xj_u, v_xj_u); /* v_xj_u = v x (xj-u) */
3140 fac = norm2(v_xj_u);
3141 /* * 1 */
3142 psijstar = 1.0/(fac + erg->rotg->eps); /* psistar = -------------------- */
3143 /* * |v x (xj-u)|^2 + eps */
3145 psij = gmx::invsqrt(fac); /* 1 */
3146 /* psij = ------------ */
3147 /* |v x (xj-u)| */
3149 svmul(psij, v_xj_u, sj); /* sj = psij * (v x (xj-u) ) */
3151 fac = iprod(v_xj_u, rj); /* fac = (v x (xj-u)).rj */
3152 fac2 = fac*fac;
3154 sjrj = iprod(sj, rj); /* sjrj = sj.rj */
3156 svmul(psijstar/psij, rj, tmpvec);
3157 svmul(psijstar*psijstar/(psij*psij*psij) * sjrj, sj, tmpvec2);
3158 rvec_dec(tmpvec, tmpvec2);
3159 cprod(tmpvec, erg->vec, tmpvec2);
3161 /* Store the additional force so that it can be added to the force
3162 * array after the normal forces have been evaluated */
3163 svmul(-erg->rotg->k*wj*sjrj, tmpvec2, tmpvec);
3164 svmul(mj, innersumvec, tmpvec2); /* This is != 0 only for the pivot-free variant */
3166 rvec_add(tmpvec2, tmpvec, erg->f_rot_loc[j]);
3167 Vpart += wj*psijstar*fac2;
3169 /* If requested, also calculate the potential for a set of angles
3170 * near the current reference angle */
3171 if (bCalcPotFit)
3173 for (int ifit = 0; ifit < erg->rotg->PotAngle_nstep; ifit++)
3175 if (bPF)
3177 mvmul(erg->PotAngleFit->rotmat[ifit], yj0_yc0, fit_rj); /* fit_rj = Omega.(yj0-yc0) */
3179 else
3181 /* Position of this atom in the collective array */
3182 int iigrp = collectiveRotationGroupIndex[j];
3183 /* Rotate with the alternative angle. Like rotate_local_reference(),
3184 * just for a single local atom */
3185 mvmul(erg->PotAngleFit->rotmat[ifit], erg->rotg->x_ref[iigrp], fit_rj); /* fit_rj = Omega*(yj0-u) */
3187 fit_fac = iprod(v_xj_u, fit_rj); /* fac = (v x (xj-u)).fit_rj */
3188 /* Add to the rotation potential for this angle: */
3189 erg->PotAngleFit->V[ifit] += 0.5*erg->rotg->k*wj*psijstar*fit_fac*fit_fac;
3193 if (bOutstepRot)
3195 /* Add to the torque of this rotation group */
3196 erg->torque_v += torque(erg->vec, erg->f_rot_loc[j], xj, erg->xc_center);
3198 /* Calculate the angle between reference and actual rotation group atom. */
3199 angle(erg, xj_u, rj, &alpha, &weight); /* angle in rad, weighted */
3200 erg->angle_v += alpha * weight;
3201 erg->weight_v += weight;
3204 PRINT_FORCE_J
3206 } /* end of loop over local rotation group atoms */
3207 erg->V = 0.5*erg->rotg->k*Vpart;
3211 /* Determine the smallest and largest position vector (with respect to the
3212 * rotation vector) for the reference group */
3213 static void get_firstlast_atom_ref(
3214 const gmx_enfrotgrp *erg,
3215 int *firstindex,
3216 int *lastindex)
3218 int i;
3219 real xcproj; /* The projection of a reference position on the
3220 rotation vector */
3221 real minproj, maxproj; /* Smallest and largest projection on v */
3223 /* Start with some value */
3224 minproj = iprod(erg->rotg->x_ref[0], erg->vec);
3225 maxproj = minproj;
3227 /* This is just to ensure that it still works if all the atoms of the
3228 * reference structure are situated in a plane perpendicular to the rotation
3229 * vector */
3230 *firstindex = 0;
3231 *lastindex = erg->rotg->nat-1;
3233 /* Loop over all atoms of the reference group,
3234 * project them on the rotation vector to find the extremes */
3235 for (i = 0; i < erg->rotg->nat; i++)
3237 xcproj = iprod(erg->rotg->x_ref[i], erg->vec);
3238 if (xcproj < minproj)
3240 minproj = xcproj;
3241 *firstindex = i;
3243 if (xcproj > maxproj)
3245 maxproj = xcproj;
3246 *lastindex = i;
3252 /* Allocate memory for the slabs */
3253 static void allocate_slabs(
3254 gmx_enfrotgrp *erg,
3255 FILE *fplog,
3256 gmx_bool bVerbose)
3258 /* More slabs than are defined for the reference are never needed */
3259 int nslabs = erg->slab_last_ref - erg->slab_first_ref + 1;
3261 /* Remember how many we allocated */
3262 erg->nslabs_alloc = nslabs;
3264 if ( (nullptr != fplog) && bVerbose)
3266 fprintf(fplog, "%s allocating memory to store data for %d slabs (rotation group %d).\n",
3267 RotStr, nslabs, erg->groupIndex);
3269 snew(erg->slab_center, nslabs);
3270 snew(erg->slab_center_ref, nslabs);
3271 snew(erg->slab_weights, nslabs);
3272 snew(erg->slab_torque_v, nslabs);
3273 snew(erg->slab_data, nslabs);
3274 snew(erg->gn_atom, nslabs);
3275 snew(erg->gn_slabind, nslabs);
3276 snew(erg->slab_innersumvec, nslabs);
3277 for (int i = 0; i < nslabs; i++)
3279 snew(erg->slab_data[i].x, erg->rotg->nat);
3280 snew(erg->slab_data[i].ref, erg->rotg->nat);
3281 snew(erg->slab_data[i].weight, erg->rotg->nat);
3283 snew(erg->xc_ref_sorted, erg->rotg->nat);
3284 snew(erg->xc_sortind, erg->rotg->nat);
3285 snew(erg->firstatom, nslabs);
3286 snew(erg->lastatom, nslabs);
3290 /* From the extreme positions of the reference group, determine the first
3291 * and last slab of the reference. We can never have more slabs in the real
3292 * simulation than calculated here for the reference.
3294 static void get_firstlast_slab_ref(gmx_enfrotgrp *erg,
3295 real mc[], int ref_firstindex, int ref_lastindex)
3297 rvec dummy;
3299 int first = get_first_slab(erg, erg->rotg->x_ref[ref_firstindex]);
3300 int last = get_last_slab(erg, erg->rotg->x_ref[ref_lastindex ]);
3302 while (get_slab_weight(first, erg, erg->rotg->x_ref, mc, &dummy) > WEIGHT_MIN)
3304 first--;
3306 erg->slab_first_ref = first+1;
3307 while (get_slab_weight(last, erg, erg->rotg->x_ref, mc, &dummy) > WEIGHT_MIN)
3309 last++;
3311 erg->slab_last_ref = last-1;
3315 /* Special version of copy_rvec:
3316 * During the copy procedure of xcurr to b, the correct PBC image is chosen
3317 * such that the copied vector ends up near its reference position xref */
3318 static inline void copy_correct_pbc_image(
3319 const rvec xcurr, /* copy vector xcurr ... */
3320 rvec b, /* ... to b ... */
3321 const rvec xref, /* choosing the PBC image such that b ends up near xref */
3322 const matrix box,
3323 int npbcdim)
3325 rvec dx;
3326 int d, m;
3327 ivec shift;
3330 /* Shortest PBC distance between the atom and its reference */
3331 rvec_sub(xcurr, xref, dx);
3333 /* Determine the shift for this atom */
3334 clear_ivec(shift);
3335 for (m = npbcdim-1; m >= 0; m--)
3337 while (dx[m] < -0.5*box[m][m])
3339 for (d = 0; d < DIM; d++)
3341 dx[d] += box[m][d];
3343 shift[m]++;
3345 while (dx[m] >= 0.5*box[m][m])
3347 for (d = 0; d < DIM; d++)
3349 dx[d] -= box[m][d];
3351 shift[m]--;
3355 /* Apply the shift to the position */
3356 copy_rvec(xcurr, b);
3357 shift_single_coord(box, b, shift);
3361 static void init_rot_group(FILE *fplog, const t_commrec *cr,
3362 gmx_enfrotgrp *erg,
3363 rvec *x, gmx_mtop_t *mtop, gmx_bool bVerbose, FILE *out_slabs, const matrix box,
3364 t_inputrec *ir, gmx_bool bOutputCenters)
3366 rvec coord, xref, *xdum;
3367 gmx_bool bFlex, bColl;
3368 int ref_firstindex, ref_lastindex;
3369 real mass, totalmass;
3370 real start = 0.0;
3371 double t_start;
3372 const t_rotgrp *rotg = erg->rotg;
3375 /* Do we have a flexible axis? */
3376 bFlex = ISFLEX(rotg);
3377 /* Do we use a global set of coordinates? */
3378 bColl = ISCOLL(rotg);
3380 /* Allocate space for collective coordinates if needed */
3381 if (bColl)
3383 snew(erg->xc, erg->rotg->nat);
3384 snew(erg->xc_shifts, erg->rotg->nat);
3385 snew(erg->xc_eshifts, erg->rotg->nat);
3386 snew(erg->xc_old, erg->rotg->nat);
3388 if (erg->rotg->eFittype == erotgFitNORM)
3390 snew(erg->xc_ref_length, erg->rotg->nat); /* in case fit type NORM is chosen */
3391 snew(erg->xc_norm, erg->rotg->nat);
3394 else
3396 snew(erg->xr_loc, erg->rotg->nat);
3397 snew(erg->x_loc_pbc, erg->rotg->nat);
3400 copy_rvec(erg->rotg->inputVec, erg->vec);
3401 snew(erg->f_rot_loc, erg->rotg->nat);
3403 /* Make space for the calculation of the potential at other angles (used
3404 * for fitting only) */
3405 if (erotgFitPOT == erg->rotg->eFittype)
3407 snew(erg->PotAngleFit, 1);
3408 snew(erg->PotAngleFit->degangle, erg->rotg->PotAngle_nstep);
3409 snew(erg->PotAngleFit->V, erg->rotg->PotAngle_nstep);
3410 snew(erg->PotAngleFit->rotmat, erg->rotg->PotAngle_nstep);
3412 /* Get the set of angles around the reference angle */
3413 start = -0.5 * (erg->rotg->PotAngle_nstep - 1)*erg->rotg->PotAngle_step;
3414 for (int i = 0; i < erg->rotg->PotAngle_nstep; i++)
3416 erg->PotAngleFit->degangle[i] = start + i*erg->rotg->PotAngle_step;
3419 else
3421 erg->PotAngleFit = nullptr;
3424 /* Copy the masses so that the center can be determined. For all types of
3425 * enforced rotation, we store the masses in the erg->mc array. */
3426 snew(erg->mc, erg->rotg->nat);
3427 if (bFlex)
3429 snew(erg->mc_sorted, erg->rotg->nat);
3431 if (!bColl)
3433 snew(erg->m_loc, erg->rotg->nat);
3435 totalmass = 0.0;
3436 int molb = 0;
3437 for (int i = 0; i < erg->rotg->nat; i++)
3439 if (erg->rotg->bMassW)
3441 mass = mtopGetAtomMass(mtop, erg->rotg->ind[i], &molb);
3443 else
3445 mass = 1.0;
3447 erg->mc[i] = mass;
3448 totalmass += mass;
3450 erg->invmass = 1.0/totalmass;
3452 /* Set xc_ref_center for any rotation potential */
3453 if ((erg->rotg->eType == erotgISO) || (erg->rotg->eType == erotgPM) || (erg->rotg->eType == erotgRM) || (erg->rotg->eType == erotgRM2))
3455 /* Set the pivot point for the fixed, stationary-axis potentials. This
3456 * won't change during the simulation */
3457 copy_rvec(erg->rotg->pivot, erg->xc_ref_center);
3458 copy_rvec(erg->rotg->pivot, erg->xc_center );
3460 else
3462 /* Center of the reference positions */
3463 get_center(erg->rotg->x_ref, erg->mc, erg->rotg->nat, erg->xc_ref_center);
3465 /* Center of the actual positions */
3466 if (MASTER(cr))
3468 snew(xdum, erg->rotg->nat);
3469 for (int i = 0; i < erg->rotg->nat; i++)
3471 int ii = erg->rotg->ind[i];
3472 copy_rvec(x[ii], xdum[i]);
3474 get_center(xdum, erg->mc, erg->rotg->nat, erg->xc_center);
3475 sfree(xdum);
3477 #if GMX_MPI
3478 if (PAR(cr))
3480 gmx_bcast(sizeof(erg->xc_center), erg->xc_center, cr);
3482 #endif
3485 if (bColl)
3487 /* Save the original (whole) set of positions in xc_old such that at later
3488 * steps the rotation group can always be made whole again. If the simulation is
3489 * restarted, we compute the starting reference positions (given the time)
3490 * and assume that the correct PBC image of each position is the one nearest
3491 * to the current reference */
3492 if (MASTER(cr))
3494 /* Calculate the rotation matrix for this angle: */
3495 t_start = ir->init_t + ir->init_step*ir->delta_t;
3496 erg->degangle = erg->rotg->rate * t_start;
3497 calc_rotmat(erg->vec, erg->degangle, erg->rotmat);
3499 for (int i = 0; i < erg->rotg->nat; i++)
3501 int ii = erg->rotg->ind[i];
3503 /* Subtract pivot, rotate, and add pivot again. This will yield the
3504 * reference position for time t */
3505 rvec_sub(erg->rotg->x_ref[i], erg->xc_ref_center, coord);
3506 mvmul(erg->rotmat, coord, xref);
3507 rvec_inc(xref, erg->xc_ref_center);
3509 copy_correct_pbc_image(x[ii], erg->xc_old[i], xref, box, 3);
3512 #if GMX_MPI
3513 if (PAR(cr))
3515 gmx_bcast(erg->rotg->nat*sizeof(erg->xc_old[0]), erg->xc_old, cr);
3517 #endif
3520 if ( (erg->rotg->eType != erotgFLEX) && (erg->rotg->eType != erotgFLEX2) )
3522 /* Put the reference positions into origin: */
3523 for (int i = 0; i < erg->rotg->nat; i++)
3525 rvec_dec(erg->rotg->x_ref[i], erg->xc_ref_center);
3529 /* Enforced rotation with flexible axis */
3530 if (bFlex)
3532 /* Calculate maximum beta value from minimum gaussian (performance opt.) */
3533 erg->max_beta = calc_beta_max(erg->rotg->min_gaussian, erg->rotg->slab_dist);
3535 /* Determine the smallest and largest coordinate with respect to the rotation vector */
3536 get_firstlast_atom_ref(erg, &ref_firstindex, &ref_lastindex);
3538 /* From the extreme positions of the reference group, determine the first
3539 * and last slab of the reference. */
3540 get_firstlast_slab_ref(erg, erg->mc, ref_firstindex, ref_lastindex);
3542 /* Allocate memory for the slabs */
3543 allocate_slabs(erg, fplog, bVerbose);
3545 /* Flexible rotation: determine the reference centers for the rest of the simulation */
3546 erg->slab_first = erg->slab_first_ref;
3547 erg->slab_last = erg->slab_last_ref;
3548 get_slab_centers(erg, erg->rotg->x_ref, erg->mc, -1, out_slabs, bOutputCenters, TRUE);
3550 /* Length of each x_rotref vector from center (needed if fit routine NORM is chosen): */
3551 if (erg->rotg->eFittype == erotgFitNORM)
3553 for (int i = 0; i < erg->rotg->nat; i++)
3555 rvec_sub(erg->rotg->x_ref[i], erg->xc_ref_center, coord);
3556 erg->xc_ref_length[i] = norm(coord);
3562 /* Calculate the size of the MPI buffer needed in reduce_output() */
3563 static int calc_mpi_bufsize(const gmx_enfrot *er)
3566 int count_total = 0;
3567 for (int g = 0; g < er->rot->ngrp; g++)
3569 const t_rotgrp *rotg = &er->rot->grp[g];
3570 const gmx_enfrotgrp *erg = &er->enfrotgrp[g];
3572 /* Count the items that are transferred for this group: */
3573 int count_group = 4; /* V, torque, angle, weight */
3575 /* Add the maximum number of slabs for flexible groups */
3576 if (ISFLEX(rotg))
3578 count_group += erg->slab_last_ref - erg->slab_first_ref + 1;
3581 /* Add space for the potentials at different angles: */
3582 if (erotgFitPOT == erg->rotg->eFittype)
3584 count_group += erg->rotg->PotAngle_nstep;
3587 /* Add to the total number: */
3588 count_total += count_group;
3591 return count_total;
3595 std::unique_ptr<gmx::EnforcedRotation>
3596 init_rot(FILE *fplog, t_inputrec *ir, int nfile, const t_filenm fnm[],
3597 const t_commrec *cr, gmx::LocalAtomSetManager * atomSets, const t_state *globalState, gmx_mtop_t *mtop, const gmx_output_env_t *oenv,
3598 const gmx::MdrunOptions &mdrunOptions,
3599 const gmx::StartingBehavior startingBehavior)
3601 int nat_max = 0; /* Size of biggest rotation group */
3602 rvec *x_pbc = nullptr; /* Space for the pbc-correct atom positions */
3604 if (MASTER(cr) && mdrunOptions.verbose)
3606 fprintf(stdout, "%s Initializing ...\n", RotStr);
3609 auto enforcedRotation = std::make_unique<gmx::EnforcedRotation>();
3610 gmx_enfrot *er = enforcedRotation->getLegacyEnfrot();
3611 // TODO When this module implements IMdpOptions, the ownership will become more clear.
3612 er->rot = ir->rot;
3613 er->restartWithAppending = (startingBehavior == gmx::StartingBehavior::RestartWithAppending);
3615 /* When appending, skip first output to avoid duplicate entries in the data files */
3616 if (er->restartWithAppending)
3618 er->bOut = FALSE;
3620 else
3622 er->bOut = TRUE;
3625 if (MASTER(cr) && er->bOut)
3627 please_cite(fplog, "Kutzner2011");
3630 /* Output every step for reruns */
3631 if (mdrunOptions.rerun)
3633 if (nullptr != fplog)
3635 fprintf(fplog, "%s rerun - will write rotation output every available step.\n", RotStr);
3637 er->nstrout = 1;
3638 er->nstsout = 1;
3640 else
3642 er->nstrout = er->rot->nstrout;
3643 er->nstsout = er->rot->nstsout;
3646 er->out_slabs = nullptr;
3647 if (MASTER(cr) && HaveFlexibleGroups(er->rot) )
3649 er->out_slabs = open_slab_out(opt2fn("-rs", nfile, fnm), er);
3652 if (MASTER(cr))
3654 /* Remove pbc, make molecule whole.
3655 * When ir->bContinuation=TRUE this has already been done, but ok. */
3656 snew(x_pbc, mtop->natoms);
3657 copy_rvecn(globalState->x.rvec_array(), x_pbc, 0, mtop->natoms);
3658 do_pbc_first_mtop(nullptr, ir->ePBC, globalState->box, mtop, x_pbc);
3659 /* All molecules will be whole now, but not necessarily in the home box.
3660 * Additionally, if a rotation group consists of more than one molecule
3661 * (e.g. two strands of DNA), each one of them can end up in a different
3662 * periodic box. This is taken care of in init_rot_group. */
3665 /* Allocate space for the per-rotation-group data: */
3666 er->enfrotgrp.resize(er->rot->ngrp);
3667 int groupIndex = 0;
3668 for (auto &ergRef : er->enfrotgrp)
3670 gmx_enfrotgrp *erg = &ergRef;
3671 erg->rotg = &er->rot->grp[groupIndex];
3672 erg->atomSet = std::make_unique<gmx::LocalAtomSet>(atomSets->add({erg->rotg->ind, erg->rotg->ind + erg->rotg->nat}));
3673 erg->groupIndex = groupIndex;
3675 if (nullptr != fplog)
3677 fprintf(fplog, "%s group %d type '%s'\n", RotStr, groupIndex, erotg_names[erg->rotg->eType]);
3680 if (erg->rotg->nat > 0)
3682 nat_max = std::max(nat_max, erg->rotg->nat);
3684 init_rot_group(fplog, cr, erg, x_pbc, mtop, mdrunOptions.verbose, er->out_slabs, MASTER(cr) ? globalState->box : nullptr, ir,
3685 !er->restartWithAppending); /* Do not output the reference centers
3686 * again if we are appending */
3688 ++groupIndex;
3691 /* Allocate space for enforced rotation buffer variables */
3692 er->bufsize = nat_max;
3693 snew(er->data, nat_max);
3694 snew(er->xbuf, nat_max);
3695 snew(er->mbuf, nat_max);
3697 /* Buffers for MPI reducing torques, angles, weights (for each group), and V */
3698 if (PAR(cr))
3700 er->mpi_bufsize = calc_mpi_bufsize(er) + 100; /* larger to catch errors */
3701 snew(er->mpi_inbuf, er->mpi_bufsize);
3702 snew(er->mpi_outbuf, er->mpi_bufsize);
3704 else
3706 er->mpi_bufsize = 0;
3707 er->mpi_inbuf = nullptr;
3708 er->mpi_outbuf = nullptr;
3711 /* Only do I/O on the MASTER */
3712 er->out_angles = nullptr;
3713 er->out_rot = nullptr;
3714 er->out_torque = nullptr;
3715 if (MASTER(cr))
3717 er->out_rot = open_rot_out(opt2fn("-ro", nfile, fnm), oenv, er);
3719 if (er->nstsout > 0)
3721 if (HaveFlexibleGroups(er->rot) || HavePotFitGroups(er->rot) )
3723 er->out_angles = open_angles_out(opt2fn("-ra", nfile, fnm), er);
3725 if (HaveFlexibleGroups(er->rot) )
3727 er->out_torque = open_torque_out(opt2fn("-rt", nfile, fnm), er);
3731 sfree(x_pbc);
3733 return enforcedRotation;
3736 /* Rotate the local reference positions and store them in
3737 * erg->xr_loc[0...(nat_loc-1)]
3739 * Note that we already subtracted u or y_c from the reference positions
3740 * in init_rot_group().
3742 static void rotate_local_reference(gmx_enfrotgrp *erg)
3744 const auto &collectiveRotationGroupIndex = erg->atomSet->collectiveIndex();
3745 for (size_t i = 0; i < erg->atomSet->numAtomsLocal(); i++)
3747 /* Index of this rotation group atom with respect to the whole rotation group */
3748 int ii = collectiveRotationGroupIndex[i];
3749 /* Rotate */
3750 mvmul(erg->rotmat, erg->rotg->x_ref[ii], erg->xr_loc[i]);
3755 /* Select the PBC representation for each local x position and store that
3756 * for later usage. We assume the right PBC image of an x is the one nearest to
3757 * its rotated reference */
3758 static void choose_pbc_image(rvec x[],
3759 gmx_enfrotgrp *erg,
3760 matrix box, int npbcdim)
3762 const auto &localRotationGroupIndex = erg->atomSet->localIndex();
3763 for (gmx::index i = 0; i < localRotationGroupIndex.ssize(); i++)
3765 /* Index of a rotation group atom */
3766 int ii = localRotationGroupIndex[i];
3768 /* Get the correctly rotated reference position. The pivot was already
3769 * subtracted in init_rot_group() from the reference positions. Also,
3770 * the reference positions have already been rotated in
3771 * rotate_local_reference(). For the current reference position we thus
3772 * only need to add the pivot again. */
3773 rvec xref;
3774 copy_rvec(erg->xr_loc[i], xref);
3775 rvec_inc(xref, erg->xc_ref_center);
3777 copy_correct_pbc_image(x[ii], erg->x_loc_pbc[i], xref, box, npbcdim);
3782 void do_rotation(const t_commrec *cr,
3783 gmx_enfrot *er,
3784 matrix box,
3785 rvec x[],
3786 real t,
3787 int64_t step,
3788 gmx_bool bNS)
3790 gmx_bool outstep_slab, outstep_rot;
3791 gmx_bool bColl;
3792 rvec transvec;
3793 gmx_potfit *fit = nullptr; /* For fit type 'potential' determine the fit
3794 angle via the potential minimum */
3796 #ifdef TAKETIME
3797 double t0;
3798 #endif
3800 /* When to output in main rotation output file */
3801 outstep_rot = do_per_step(step, er->nstrout) && er->bOut;
3802 /* When to output per-slab data */
3803 outstep_slab = do_per_step(step, er->nstsout) && er->bOut;
3805 /* Output time into rotation output file */
3806 if (outstep_rot && MASTER(cr))
3808 fprintf(er->out_rot, "%12.3e", t);
3811 /**************************************************************************/
3812 /* First do ALL the communication! */
3813 for (auto &ergRef : er->enfrotgrp)
3815 gmx_enfrotgrp *erg = &ergRef;
3816 const t_rotgrp *rotg = erg->rotg;
3818 /* Do we use a collective (global) set of coordinates? */
3819 bColl = ISCOLL(rotg);
3821 /* Calculate the rotation matrix for this angle: */
3822 erg->degangle = rotg->rate * t;
3823 calc_rotmat(erg->vec, erg->degangle, erg->rotmat);
3825 if (bColl)
3827 /* Transfer the rotation group's positions such that every node has
3828 * all of them. Every node contributes its local positions x and stores
3829 * it in the collective erg->xc array. */
3830 communicate_group_positions(cr, erg->xc, erg->xc_shifts, erg->xc_eshifts, bNS,
3831 x, rotg->nat, erg->atomSet->numAtomsLocal(), erg->atomSet->localIndex().data(), erg->atomSet->collectiveIndex().data(), erg->xc_old, box);
3833 else
3835 /* Fill the local masses array;
3836 * this array changes in DD/neighborsearching steps */
3837 if (bNS)
3839 const auto &collectiveRotationGroupIndex = erg->atomSet->collectiveIndex();
3840 for (gmx::index i = 0; i < collectiveRotationGroupIndex.ssize(); i++)
3842 /* Index of local atom w.r.t. the collective rotation group */
3843 int ii = collectiveRotationGroupIndex[i];
3844 erg->m_loc[i] = erg->mc[ii];
3848 /* Calculate Omega*(y_i-y_c) for the local positions */
3849 rotate_local_reference(erg);
3851 /* Choose the nearest PBC images of the group atoms with respect
3852 * to the rotated reference positions */
3853 choose_pbc_image(x, erg, box, 3);
3855 /* Get the center of the rotation group */
3856 if ( (rotg->eType == erotgISOPF) || (rotg->eType == erotgPMPF) )
3858 get_center_comm(cr, erg->x_loc_pbc, erg->m_loc, erg->atomSet->numAtomsLocal(), rotg->nat, erg->xc_center);
3862 } /* End of loop over rotation groups */
3864 /**************************************************************************/
3865 /* Done communicating, we can start to count cycles for the load balancing now ... */
3866 if (DOMAINDECOMP(cr))
3868 ddReopenBalanceRegionCpu(cr->dd);
3871 #ifdef TAKETIME
3872 t0 = MPI_Wtime();
3873 #endif
3875 for (auto &ergRef : er->enfrotgrp)
3877 gmx_enfrotgrp *erg = &ergRef;
3878 const t_rotgrp *rotg = erg->rotg;
3880 if (outstep_rot && MASTER(cr))
3882 fprintf(er->out_rot, "%12.4f", erg->degangle);
3885 /* Calculate angles and rotation matrices for potential fitting: */
3886 if ( (outstep_rot || outstep_slab) && (erotgFitPOT == rotg->eFittype) )
3888 fit = erg->PotAngleFit;
3889 for (int i = 0; i < rotg->PotAngle_nstep; i++)
3891 calc_rotmat(erg->vec, erg->degangle + fit->degangle[i], fit->rotmat[i]);
3893 /* Clear value from last step */
3894 erg->PotAngleFit->V[i] = 0.0;
3898 /* Clear values from last time step */
3899 erg->V = 0.0;
3900 erg->torque_v = 0.0;
3901 erg->angle_v = 0.0;
3902 erg->weight_v = 0.0;
3904 switch (rotg->eType)
3906 case erotgISO:
3907 case erotgISOPF:
3908 case erotgPM:
3909 case erotgPMPF:
3910 do_fixed(erg, outstep_rot, outstep_slab);
3911 break;
3912 case erotgRM:
3913 do_radial_motion(erg, outstep_rot, outstep_slab);
3914 break;
3915 case erotgRMPF:
3916 do_radial_motion_pf(erg, x, box, outstep_rot, outstep_slab);
3917 break;
3918 case erotgRM2:
3919 case erotgRM2PF:
3920 do_radial_motion2(erg, x, box, outstep_rot, outstep_slab);
3921 break;
3922 case erotgFLEXT:
3923 case erotgFLEX2T:
3924 /* Subtract the center of the rotation group from the collective positions array
3925 * Also store the center in erg->xc_center since it needs to be subtracted
3926 * in the low level routines from the local coordinates as well */
3927 get_center(erg->xc, erg->mc, rotg->nat, erg->xc_center);
3928 svmul(-1.0, erg->xc_center, transvec);
3929 translate_x(erg->xc, rotg->nat, transvec);
3930 do_flexible(MASTER(cr), er, erg, x, box, t, outstep_rot, outstep_slab);
3931 break;
3932 case erotgFLEX:
3933 case erotgFLEX2:
3934 /* Do NOT subtract the center of mass in the low level routines! */
3935 clear_rvec(erg->xc_center);
3936 do_flexible(MASTER(cr), er, erg, x, box, t, outstep_rot, outstep_slab);
3937 break;
3938 default:
3939 gmx_fatal(FARGS, "No such rotation potential.");
3943 #ifdef TAKETIME
3944 if (MASTER(cr))
3946 fprintf(stderr, "%s calculation (step %d) took %g seconds.\n", RotStr, step, MPI_Wtime()-t0);
3948 #endif