2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2008, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017, by the GROMACS development team, 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.
39 #include "pull_rotation.h"
49 #include "gromacs/commandline/filenm.h"
50 #include "gromacs/domdec/domdec.h"
51 #include "gromacs/domdec/domdec_struct.h"
52 #include "gromacs/domdec/ga2la.h"
53 #include "gromacs/fileio/gmxfio.h"
54 #include "gromacs/fileio/xvgr.h"
55 #include "gromacs/gmxlib/network.h"
56 #include "gromacs/linearalgebra/nrjac.h"
57 #include "gromacs/math/functions.h"
58 #include "gromacs/math/utilities.h"
59 #include "gromacs/math/vec.h"
60 #include "gromacs/mdlib/groupcoord.h"
61 #include "gromacs/mdlib/mdrun.h"
62 #include "gromacs/mdlib/sim_util.h"
63 #include "gromacs/mdtypes/inputrec.h"
64 #include "gromacs/mdtypes/md_enums.h"
65 #include "gromacs/pbcutil/pbc.h"
66 #include "gromacs/timing/cyclecounter.h"
67 #include "gromacs/timing/wallcycle.h"
68 #include "gromacs/topology/mtop_lookup.h"
69 #include "gromacs/utility/fatalerror.h"
70 #include "gromacs/utility/pleasecite.h"
71 #include "gromacs/utility/qsort_threadsafe.h"
72 #include "gromacs/utility/smalloc.h"
74 static char const *RotStr = {"Enforced rotation:"};
76 /* Set the minimum weight for the determination of the slab centers */
77 #define WEIGHT_MIN (10*GMX_FLOAT_MIN)
79 /* Helper structure for sorting positions along rotation vector */
81 real xcproj; /* Projection of xc on the rotation vector */
82 int ind; /* Index of xc */
84 rvec x; /* Position */
85 rvec x_ref; /* Reference position */
89 /* Enforced rotation / flexible: determine the angle of each slab */
90 typedef struct gmx_slabdata
92 int nat; /* Number of atoms belonging to this slab */
93 rvec *x; /* The positions belonging to this slab. In
94 general, this should be all positions of the
95 whole rotation group, but we leave those away
96 that have a small enough weight */
97 rvec *ref; /* Same for reference */
98 real *weight; /* The weight for each atom */
102 /* Helper structure for potential fitting */
103 typedef struct gmx_potfit
105 real *degangle; /* Set of angles for which the potential is
106 calculated. The optimum fit is determined as
107 the angle for with the potential is minimal */
108 real *V; /* Potential for the different angles */
109 matrix *rotmat; /* Rotation matrix corresponding to the angles */
113 /* Enforced rotation data for all groups */
114 typedef struct gmx_enfrot
116 FILE *out_rot; /* Output file for rotation data */
117 FILE *out_torque; /* Output file for torque data */
118 FILE *out_angles; /* Output file for slab angles for flexible type */
119 FILE *out_slabs; /* Output file for slab centers */
120 int bufsize; /* Allocation size of buf */
121 rvec *xbuf; /* Coordinate buffer variable for sorting */
122 real *mbuf; /* Masses buffer variable for sorting */
123 sort_along_vec_t *data; /* Buffer variable needed for position sorting */
124 real *mpi_inbuf; /* MPI buffer */
125 real *mpi_outbuf; /* MPI buffer */
126 int mpi_bufsize; /* Allocation size of in & outbuf */
127 unsigned long Flags; /* mdrun flags */
128 gmx_bool bOut; /* Used to skip first output when appending to
129 * avoid duplicate entries in rotation outfiles */
133 /* Global enforced rotation data for a single rotation group */
134 typedef struct gmx_enfrotgrp
136 real degangle; /* Rotation angle in degrees */
137 matrix rotmat; /* Rotation matrix */
138 int *ind_loc; /* Local rotation indices */
139 int nat_loc; /* Number of local group atoms */
140 int nalloc_loc; /* Allocation size for ind_loc and weight_loc */
142 real V; /* Rotation potential for this rotation group */
143 rvec *f_rot_loc; /* Array to store the forces on the local atoms
144 resulting from enforced rotation potential */
146 /* Collective coordinates for the whole rotation group */
147 real *xc_ref_length; /* Length of each x_rotref vector after x_rotref
148 has been put into origin */
149 int *xc_ref_ind; /* Position of each local atom in the collective
151 rvec xc_center; /* Center of the rotation group positions, may
153 rvec xc_ref_center; /* dito, for the reference positions */
154 rvec *xc; /* Current (collective) positions */
155 ivec *xc_shifts; /* Current (collective) shifts */
156 ivec *xc_eshifts; /* Extra shifts since last DD step */
157 rvec *xc_old; /* Old (collective) positions */
158 rvec *xc_norm; /* Normalized form of the current positions */
159 rvec *xc_ref_sorted; /* Reference positions (sorted in the same order
160 as xc when sorted) */
161 int *xc_sortind; /* Where is a position found after sorting? */
162 real *mc; /* Collective masses */
164 real invmass; /* one over the total mass of the rotation group */
166 real torque_v; /* Torque in the direction of rotation vector */
167 real angle_v; /* Actual angle of the whole rotation group */
168 /* Fixed rotation only */
169 real weight_v; /* Weights for angle determination */
170 rvec *xr_loc; /* Local reference coords, correctly rotated */
171 rvec *x_loc_pbc; /* Local current coords, correct PBC image */
172 real *m_loc; /* Masses of the current local atoms */
174 /* Flexible rotation only */
175 int nslabs_alloc; /* For this many slabs memory is allocated */
176 int slab_first; /* Lowermost slab for that the calculation needs
177 to be performed at a given time step */
178 int slab_last; /* Uppermost slab ... */
179 int slab_first_ref; /* First slab for which ref. center is stored */
180 int slab_last_ref; /* Last ... */
181 int slab_buffer; /* Slab buffer region around reference slabs */
182 int *firstatom; /* First relevant atom for a slab */
183 int *lastatom; /* Last relevant atom for a slab */
184 rvec *slab_center; /* Gaussian-weighted slab center */
185 rvec *slab_center_ref; /* Gaussian-weighted slab center for the
186 reference positions */
187 real *slab_weights; /* Sum of gaussian weights in a slab */
188 real *slab_torque_v; /* Torque T = r x f for each slab. */
189 /* torque_v = m.v = angular momentum in the
191 real max_beta; /* min_gaussian from inputrec->rotgrp is the
192 minimum value the gaussian must have so that
193 the force is actually evaluated max_beta is
194 just another way to put it */
195 real *gn_atom; /* Precalculated gaussians for a single atom */
196 int *gn_slabind; /* Tells to which slab each precalculated gaussian
198 rvec *slab_innersumvec; /* Inner sum of the flexible2 potential per slab;
199 this is precalculated for optimization reasons */
200 t_gmx_slabdata *slab_data; /* Holds atom positions and gaussian weights
201 of atoms belonging to a slab */
203 /* For potential fits with varying angle: */
204 t_gmx_potfit *PotAngleFit; /* Used for fit type 'potential' */
208 /* Activate output of forces for correctness checks */
209 /* #define PRINT_FORCES */
211 #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]);
212 #define PRINT_POT_TAU if (MASTER(cr)) { \
213 fprintf(stderr, "potential = %15.8f\n" "torque = %15.8f\n", erg->V, erg->torque_v); \
216 #define PRINT_FORCE_J
217 #define PRINT_POT_TAU
220 /* Shortcuts for often used queries */
221 #define ISFLEX(rg) ( (rg->eType == erotgFLEX) || (rg->eType == erotgFLEXT) || (rg->eType == erotgFLEX2) || (rg->eType == erotgFLEX2T) )
222 #define ISCOLL(rg) ( (rg->eType == erotgFLEX) || (rg->eType == erotgFLEXT) || (rg->eType == erotgFLEX2) || (rg->eType == erotgFLEX2T) || (rg->eType == erotgRMPF) || (rg->eType == erotgRM2PF) )
225 /* Does any of the rotation groups use slab decomposition? */
226 static gmx_bool HaveFlexibleGroups(t_rot *rot)
232 for (g = 0; g < rot->ngrp; g++)
245 /* Is for any group the fit angle determined by finding the minimum of the
246 * rotation potential? */
247 static gmx_bool HavePotFitGroups(t_rot *rot)
253 for (g = 0; g < rot->ngrp; g++)
256 if (erotgFitPOT == rotg->eFittype)
266 static double** allocate_square_matrix(int dim)
269 double** mat = nullptr;
273 for (i = 0; i < dim; i++)
282 static void free_square_matrix(double** mat, int dim)
287 for (i = 0; i < dim; i++)
295 /* Return the angle for which the potential is minimal */
296 static real get_fitangle(t_rotgrp *rotg, gmx_enfrotgrp_t erg)
299 real fitangle = -999.9;
300 real pot_min = GMX_FLOAT_MAX;
304 fit = erg->PotAngleFit;
306 for (i = 0; i < rotg->PotAngle_nstep; i++)
308 if (fit->V[i] < pot_min)
311 fitangle = fit->degangle[i];
319 /* Reduce potential angle fit data for this group at this time step? */
320 static gmx_inline gmx_bool bPotAngle(t_rot *rot, t_rotgrp *rotg, gmx_int64_t step)
322 return ( (erotgFitPOT == rotg->eFittype) && (do_per_step(step, rot->nstsout) || do_per_step(step, rot->nstrout)) );
325 /* Reduce slab torqe data for this group at this time step? */
326 static gmx_inline gmx_bool bSlabTau(t_rot *rot, t_rotgrp *rotg, gmx_int64_t step)
328 return ( (ISFLEX(rotg)) && do_per_step(step, rot->nstsout) );
331 /* Output rotation energy, torques, etc. for each rotation group */
332 static void reduce_output(t_commrec *cr, t_rot *rot, real t, gmx_int64_t step)
334 int g, i, islab, nslabs = 0;
335 int count; /* MPI element counter */
337 gmx_enfrot_t er; /* Pointer to the enforced rotation buffer variables */
338 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
345 /* Fill the MPI buffer with stuff to reduce. If items are added for reduction
346 * here, the MPI buffer size has to be enlarged also in calc_mpi_bufsize() */
350 for (g = 0; g < rot->ngrp; g++)
353 erg = rotg->enfrotgrp;
354 nslabs = erg->slab_last - erg->slab_first + 1;
355 er->mpi_inbuf[count++] = erg->V;
356 er->mpi_inbuf[count++] = erg->torque_v;
357 er->mpi_inbuf[count++] = erg->angle_v;
358 er->mpi_inbuf[count++] = erg->weight_v; /* weights are not needed for flex types, but this is just a single value */
360 if (bPotAngle(rot, rotg, step))
362 for (i = 0; i < rotg->PotAngle_nstep; i++)
364 er->mpi_inbuf[count++] = erg->PotAngleFit->V[i];
367 if (bSlabTau(rot, rotg, step))
369 for (i = 0; i < nslabs; i++)
371 er->mpi_inbuf[count++] = erg->slab_torque_v[i];
375 if (count > er->mpi_bufsize)
377 gmx_fatal(FARGS, "%s MPI buffer overflow, please report this error.", RotStr);
381 MPI_Reduce(er->mpi_inbuf, er->mpi_outbuf, count, GMX_MPI_REAL, MPI_SUM, MASTERRANK(cr), cr->mpi_comm_mygroup);
384 /* Copy back the reduced data from the buffer on the master */
388 for (g = 0; g < rot->ngrp; g++)
391 erg = rotg->enfrotgrp;
392 nslabs = erg->slab_last - erg->slab_first + 1;
393 erg->V = er->mpi_outbuf[count++];
394 erg->torque_v = er->mpi_outbuf[count++];
395 erg->angle_v = er->mpi_outbuf[count++];
396 erg->weight_v = er->mpi_outbuf[count++];
398 if (bPotAngle(rot, rotg, step))
400 for (i = 0; i < rotg->PotAngle_nstep; i++)
402 erg->PotAngleFit->V[i] = er->mpi_outbuf[count++];
405 if (bSlabTau(rot, rotg, step))
407 for (i = 0; i < nslabs; i++)
409 erg->slab_torque_v[i] = er->mpi_outbuf[count++];
419 /* Angle and torque for each rotation group */
420 for (g = 0; g < rot->ngrp; g++)
423 bFlex = ISFLEX(rotg);
425 erg = rotg->enfrotgrp;
427 /* Output to main rotation output file: */
428 if (do_per_step(step, rot->nstrout) )
430 if (erotgFitPOT == rotg->eFittype)
432 fitangle = get_fitangle(rotg, erg);
438 fitangle = erg->angle_v; /* RMSD fit angle */
442 fitangle = (erg->angle_v/erg->weight_v)*180.0*M_1_PI;
445 fprintf(er->out_rot, "%12.4f", fitangle);
446 fprintf(er->out_rot, "%12.3e", erg->torque_v);
447 fprintf(er->out_rot, "%12.3e", erg->V);
450 if (do_per_step(step, rot->nstsout) )
452 /* Output to torque log file: */
455 fprintf(er->out_torque, "%12.3e%6d", t, g);
456 for (i = erg->slab_first; i <= erg->slab_last; i++)
458 islab = i - erg->slab_first; /* slab index */
459 /* Only output if enough weight is in slab */
460 if (erg->slab_weights[islab] > rotg->min_gaussian)
462 fprintf(er->out_torque, "%6d%12.3e", i, erg->slab_torque_v[islab]);
465 fprintf(er->out_torque, "\n");
468 /* Output to angles log file: */
469 if (erotgFitPOT == rotg->eFittype)
471 fprintf(er->out_angles, "%12.3e%6d%12.4f", t, g, erg->degangle);
472 /* Output energies at a set of angles around the reference angle */
473 for (i = 0; i < rotg->PotAngle_nstep; i++)
475 fprintf(er->out_angles, "%12.3e", erg->PotAngleFit->V[i]);
477 fprintf(er->out_angles, "\n");
481 if (do_per_step(step, rot->nstrout) )
483 fprintf(er->out_rot, "\n");
489 /* Add the forces from enforced rotation potential to the local forces.
490 * Should be called after the SR forces have been evaluated */
491 extern real add_rot_forces(t_rot *rot, rvec f[], t_commrec *cr, gmx_int64_t step, real t)
495 gmx_enfrot_t er; /* Pointer to the enforced rotation buffer variables */
496 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
497 real Vrot = 0.0; /* If more than one rotation group is present, Vrot
498 assembles the local parts from all groups */
503 /* Loop over enforced rotation groups (usually 1, though)
504 * Apply the forces from rotation potentials */
505 for (g = 0; g < rot->ngrp; g++)
508 erg = rotg->enfrotgrp;
509 Vrot += erg->V; /* add the local parts from the nodes */
510 for (l = 0; l < erg->nat_loc; l++)
512 /* Get the right index of the local force */
513 ii = erg->ind_loc[l];
515 rvec_inc(f[ii], erg->f_rot_loc[l]);
519 /* Reduce energy,torque, angles etc. to get the sum values (per rotation group)
520 * on the master and output these values to file. */
521 if ( (do_per_step(step, rot->nstrout) || do_per_step(step, rot->nstsout)) && er->bOut)
523 reduce_output(cr, rot, t, step);
526 /* When appending, er->bOut is FALSE the first time to avoid duplicate entries */
535 /* The Gaussian norm is chosen such that the sum of the gaussian functions
536 * over the slabs is approximately 1.0 everywhere */
537 #define GAUSS_NORM 0.569917543430618
540 /* Calculate the maximum beta that leads to a gaussian larger min_gaussian,
541 * also does some checks
543 static double calc_beta_max(real min_gaussian, real slab_dist)
549 /* Actually the next two checks are already made in grompp */
552 gmx_fatal(FARGS, "Slab distance of flexible rotation groups must be >=0 !");
554 if (min_gaussian <= 0)
556 gmx_fatal(FARGS, "Cutoff value for Gaussian must be > 0. (You requested %f)");
559 /* Define the sigma value */
560 sigma = 0.7*slab_dist;
562 /* Calculate the argument for the logarithm and check that the log() result is negative or 0 */
563 arg = min_gaussian/GAUSS_NORM;
566 gmx_fatal(FARGS, "min_gaussian of flexible rotation groups must be <%g", GAUSS_NORM);
569 return std::sqrt(-2.0*sigma*sigma*log(min_gaussian/GAUSS_NORM));
573 static gmx_inline real calc_beta(rvec curr_x, t_rotgrp *rotg, int n)
575 return iprod(curr_x, rotg->vec) - rotg->slab_dist * n;
579 static gmx_inline real gaussian_weight(rvec curr_x, t_rotgrp *rotg, int n)
581 const real norm = GAUSS_NORM;
585 /* Define the sigma value */
586 sigma = 0.7*rotg->slab_dist;
587 /* Calculate the Gaussian value of slab n for position curr_x */
588 return norm * exp( -0.5 * gmx::square( calc_beta(curr_x, rotg, n)/sigma ) );
592 /* Returns the weight in a single slab, also calculates the Gaussian- and mass-
593 * weighted sum of positions for that slab */
594 static real get_slab_weight(int j, t_rotgrp *rotg, rvec xc[], real mc[], rvec *x_weighted_sum)
596 rvec curr_x; /* The position of an atom */
597 rvec curr_x_weighted; /* The gaussian-weighted position */
598 real gaussian; /* A single gaussian weight */
599 real wgauss; /* gaussian times current mass */
600 real slabweight = 0.0; /* The sum of weights in the slab */
604 clear_rvec(*x_weighted_sum);
606 /* Loop over all atoms in the rotation group */
607 for (i = 0; i < rotg->nat; i++)
609 copy_rvec(xc[i], curr_x);
610 gaussian = gaussian_weight(curr_x, rotg, j);
611 wgauss = gaussian * mc[i];
612 svmul(wgauss, curr_x, curr_x_weighted);
613 rvec_add(*x_weighted_sum, curr_x_weighted, *x_weighted_sum);
614 slabweight += wgauss;
615 } /* END of loop over rotation group atoms */
621 static void get_slab_centers(
622 t_rotgrp *rotg, /* The rotation group information */
623 rvec *xc, /* The rotation group positions; will
624 typically be enfrotgrp->xc, but at first call
625 it is enfrotgrp->xc_ref */
626 real *mc, /* The masses of the rotation group atoms */
627 int g, /* The number of the rotation group */
628 real time, /* Used for output only */
629 FILE *out_slabs, /* For outputting center per slab information */
630 gmx_bool bOutStep, /* Is this an output step? */
631 gmx_bool bReference) /* If this routine is called from
632 init_rot_group we need to store
633 the reference slab centers */
637 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
640 erg = rotg->enfrotgrp;
642 /* Loop over slabs */
643 for (j = erg->slab_first; j <= erg->slab_last; j++)
645 islab = j - erg->slab_first;
646 erg->slab_weights[islab] = get_slab_weight(j, rotg, xc, mc, &erg->slab_center[islab]);
648 /* We can do the calculations ONLY if there is weight in the slab! */
649 if (erg->slab_weights[islab] > WEIGHT_MIN)
651 svmul(1.0/erg->slab_weights[islab], erg->slab_center[islab], erg->slab_center[islab]);
655 /* We need to check this here, since we divide through slab_weights
656 * in the flexible low-level routines! */
657 gmx_fatal(FARGS, "Not enough weight in slab %d. Slab center cannot be determined!", j);
660 /* At first time step: save the centers of the reference structure */
663 copy_rvec(erg->slab_center[islab], erg->slab_center_ref[islab]);
665 } /* END of loop over slabs */
667 /* Output on the master */
668 if ( (nullptr != out_slabs) && bOutStep)
670 fprintf(out_slabs, "%12.3e%6d", time, g);
671 for (j = erg->slab_first; j <= erg->slab_last; j++)
673 islab = j - erg->slab_first;
674 fprintf(out_slabs, "%6d%12.3e%12.3e%12.3e",
675 j, erg->slab_center[islab][XX], erg->slab_center[islab][YY], erg->slab_center[islab][ZZ]);
677 fprintf(out_slabs, "\n");
682 static void calc_rotmat(
684 real degangle, /* Angle alpha of rotation at time t in degrees */
685 matrix rotmat) /* Rotation matrix */
687 real radangle; /* Rotation angle in radians */
688 real cosa; /* cosine alpha */
689 real sina; /* sine alpha */
690 real OMcosa; /* 1 - cos(alpha) */
691 real dumxy, dumxz, dumyz; /* save computations */
692 rvec rot_vec; /* Rotate around rot_vec ... */
695 radangle = degangle * M_PI/180.0;
696 copy_rvec(vec, rot_vec );
698 /* Precompute some variables: */
699 cosa = cos(radangle);
700 sina = sin(radangle);
702 dumxy = rot_vec[XX]*rot_vec[YY]*OMcosa;
703 dumxz = rot_vec[XX]*rot_vec[ZZ]*OMcosa;
704 dumyz = rot_vec[YY]*rot_vec[ZZ]*OMcosa;
706 /* Construct the rotation matrix for this rotation group: */
708 rotmat[XX][XX] = cosa + rot_vec[XX]*rot_vec[XX]*OMcosa;
709 rotmat[YY][XX] = dumxy + rot_vec[ZZ]*sina;
710 rotmat[ZZ][XX] = dumxz - rot_vec[YY]*sina;
712 rotmat[XX][YY] = dumxy - rot_vec[ZZ]*sina;
713 rotmat[YY][YY] = cosa + rot_vec[YY]*rot_vec[YY]*OMcosa;
714 rotmat[ZZ][YY] = dumyz + rot_vec[XX]*sina;
716 rotmat[XX][ZZ] = dumxz + rot_vec[YY]*sina;
717 rotmat[YY][ZZ] = dumyz - rot_vec[XX]*sina;
718 rotmat[ZZ][ZZ] = cosa + rot_vec[ZZ]*rot_vec[ZZ]*OMcosa;
723 for (iii = 0; iii < 3; iii++)
725 for (jjj = 0; jjj < 3; jjj++)
727 fprintf(stderr, " %10.8f ", rotmat[iii][jjj]);
729 fprintf(stderr, "\n");
735 /* Calculates torque on the rotation axis tau = position x force */
736 static gmx_inline real torque(
737 rvec rotvec, /* rotation vector; MUST be normalized! */
738 rvec force, /* force */
739 rvec x, /* position of atom on which the force acts */
740 rvec pivot) /* pivot point of rotation axis */
745 /* Subtract offset */
746 rvec_sub(x, pivot, vectmp);
748 /* position x force */
749 cprod(vectmp, force, tau);
751 /* Return the part of the torque which is parallel to the rotation vector */
752 return iprod(tau, rotvec);
756 /* Right-aligned output of value with standard width */
757 static void print_aligned(FILE *fp, char const *str)
759 fprintf(fp, "%12s", str);
763 /* Right-aligned output of value with standard short width */
764 static void print_aligned_short(FILE *fp, char const *str)
766 fprintf(fp, "%6s", str);
770 static FILE *open_output_file(const char *fn, int steps, const char what[])
775 fp = gmx_ffopen(fn, "w");
777 fprintf(fp, "# Output of %s is written in intervals of %d time step%s.\n#\n",
778 what, steps, steps > 1 ? "s" : "");
784 /* Open output file for slab center data. Call on master only */
785 static FILE *open_slab_out(const char *fn, t_rot *rot)
792 if (rot->enfrot->Flags & MD_APPENDFILES)
794 fp = gmx_fio_fopen(fn, "a");
798 fp = open_output_file(fn, rot->nstsout, "gaussian weighted slab centers");
800 for (g = 0; g < rot->ngrp; g++)
805 fprintf(fp, "# Rotation group %d (%s), slab distance %f nm, %s.\n",
806 g, erotg_names[rotg->eType], rotg->slab_dist,
807 rotg->bMassW ? "centers of mass" : "geometrical centers");
811 fprintf(fp, "# Reference centers are listed first (t=-1).\n");
812 fprintf(fp, "# The following columns have the syntax:\n");
814 print_aligned_short(fp, "t");
815 print_aligned_short(fp, "grp");
816 /* Print legend for the first two entries only ... */
817 for (i = 0; i < 2; i++)
819 print_aligned_short(fp, "slab");
820 print_aligned(fp, "X center");
821 print_aligned(fp, "Y center");
822 print_aligned(fp, "Z center");
824 fprintf(fp, " ...\n");
832 /* Adds 'buf' to 'str' */
833 static void add_to_string(char **str, char *buf)
838 len = strlen(*str) + strlen(buf) + 1;
844 static void add_to_string_aligned(char **str, char *buf)
846 char buf_aligned[STRLEN];
848 sprintf(buf_aligned, "%12s", buf);
849 add_to_string(str, buf_aligned);
853 /* Open output file and print some general information about the rotation groups.
854 * Call on master only */
855 static FILE *open_rot_out(const char *fn, t_rot *rot, const gmx_output_env_t *oenv)
860 const char **setname;
861 char buf[50], buf2[75];
862 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
864 char *LegendStr = nullptr;
867 if (rot->enfrot->Flags & MD_APPENDFILES)
869 fp = gmx_fio_fopen(fn, "a");
873 fp = xvgropen(fn, "Rotation angles and energy", "Time (ps)", "angles (degrees) and energies (kJ/mol)", oenv);
874 fprintf(fp, "# Output of enforced rotation data is written in intervals of %d time step%s.\n#\n", rot->nstrout, rot->nstrout > 1 ? "s" : "");
875 fprintf(fp, "# The scalar tau is the torque (kJ/mol) in the direction of the rotation vector v.\n");
876 fprintf(fp, "# To obtain the vectorial torque, multiply tau with the group's rot-vec.\n");
877 fprintf(fp, "# For flexible groups, tau(t,n) from all slabs n have been summed in a single value tau(t) here.\n");
878 fprintf(fp, "# The torques tau(t,n) are found in the rottorque.log (-rt) output file\n");
880 for (g = 0; g < rot->ngrp; g++)
883 erg = rotg->enfrotgrp;
884 bFlex = ISFLEX(rotg);
887 fprintf(fp, "# ROTATION GROUP %d, potential type '%s':\n", g, erotg_names[rotg->eType]);
888 fprintf(fp, "# rot-massw%d %s\n", g, yesno_names[rotg->bMassW]);
889 fprintf(fp, "# rot-vec%d %12.5e %12.5e %12.5e\n", g, rotg->vec[XX], rotg->vec[YY], rotg->vec[ZZ]);
890 fprintf(fp, "# rot-rate%d %12.5e degrees/ps\n", g, rotg->rate);
891 fprintf(fp, "# rot-k%d %12.5e kJ/(mol*nm^2)\n", g, rotg->k);
892 if (rotg->eType == erotgISO || rotg->eType == erotgPM || rotg->eType == erotgRM || rotg->eType == erotgRM2)
894 fprintf(fp, "# rot-pivot%d %12.5e %12.5e %12.5e nm\n", g, rotg->pivot[XX], rotg->pivot[YY], rotg->pivot[ZZ]);
899 fprintf(fp, "# rot-slab-distance%d %f nm\n", g, rotg->slab_dist);
900 fprintf(fp, "# rot-min-gaussian%d %12.5e\n", g, rotg->min_gaussian);
903 /* Output the centers of the rotation groups for the pivot-free potentials */
904 if ((rotg->eType == erotgISOPF) || (rotg->eType == erotgPMPF) || (rotg->eType == erotgRMPF) || (rotg->eType == erotgRM2PF
905 || (rotg->eType == erotgFLEXT) || (rotg->eType == erotgFLEX2T)) )
907 fprintf(fp, "# ref. grp. %d center %12.5e %12.5e %12.5e\n", g,
908 erg->xc_ref_center[XX], erg->xc_ref_center[YY], erg->xc_ref_center[ZZ]);
910 fprintf(fp, "# grp. %d init.center %12.5e %12.5e %12.5e\n", g,
911 erg->xc_center[XX], erg->xc_center[YY], erg->xc_center[ZZ]);
914 if ( (rotg->eType == erotgRM2) || (rotg->eType == erotgFLEX2) || (rotg->eType == erotgFLEX2T) )
916 fprintf(fp, "# rot-eps%d %12.5e nm^2\n", g, rotg->eps);
918 if (erotgFitPOT == rotg->eFittype)
921 fprintf(fp, "# theta_fit%d is determined by first evaluating the potential for %d angles around theta_ref%d.\n",
922 g, rotg->PotAngle_nstep, g);
923 fprintf(fp, "# The fit angle is the one with the smallest potential. It is given as the deviation\n");
924 fprintf(fp, "# from the reference angle, i.e. if theta_ref=X and theta_fit=Y, then the angle with\n");
925 fprintf(fp, "# minimal value of the potential is X+Y. Angular resolution is %g degrees.\n", rotg->PotAngle_step);
929 /* Print a nice legend */
932 sprintf(buf, "# %6s", "time");
933 add_to_string_aligned(&LegendStr, buf);
936 snew(setname, 4*rot->ngrp);
938 for (g = 0; g < rot->ngrp; g++)
940 sprintf(buf, "theta_ref%d", g);
941 add_to_string_aligned(&LegendStr, buf);
943 sprintf(buf2, "%s (degrees)", buf);
944 setname[nsets] = gmx_strdup(buf2);
947 for (g = 0; g < rot->ngrp; g++)
950 bFlex = ISFLEX(rotg);
952 /* For flexible axis rotation we use RMSD fitting to determine the
953 * actual angle of the rotation group */
954 if (bFlex || erotgFitPOT == rotg->eFittype)
956 sprintf(buf, "theta_fit%d", g);
960 sprintf(buf, "theta_av%d", g);
962 add_to_string_aligned(&LegendStr, buf);
963 sprintf(buf2, "%s (degrees)", buf);
964 setname[nsets] = gmx_strdup(buf2);
967 sprintf(buf, "tau%d", g);
968 add_to_string_aligned(&LegendStr, buf);
969 sprintf(buf2, "%s (kJ/mol)", buf);
970 setname[nsets] = gmx_strdup(buf2);
973 sprintf(buf, "energy%d", g);
974 add_to_string_aligned(&LegendStr, buf);
975 sprintf(buf2, "%s (kJ/mol)", buf);
976 setname[nsets] = gmx_strdup(buf2);
983 xvgr_legend(fp, nsets, setname, oenv);
987 fprintf(fp, "#\n# Legend for the following data columns:\n");
988 fprintf(fp, "%s\n", LegendStr);
998 /* Call on master only */
999 static FILE *open_angles_out(const char *fn, t_rot *rot)
1004 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
1008 if (rot->enfrot->Flags & MD_APPENDFILES)
1010 fp = gmx_fio_fopen(fn, "a");
1014 /* Open output file and write some information about it's structure: */
1015 fp = open_output_file(fn, rot->nstsout, "rotation group angles");
1016 fprintf(fp, "# All angles given in degrees, time in ps.\n");
1017 for (g = 0; g < rot->ngrp; g++)
1019 rotg = &rot->grp[g];
1020 erg = rotg->enfrotgrp;
1022 /* Output for this group happens only if potential type is flexible or
1023 * if fit type is potential! */
1024 if (ISFLEX(rotg) || (erotgFitPOT == rotg->eFittype) )
1028 sprintf(buf, " slab distance %f nm, ", rotg->slab_dist);
1035 fprintf(fp, "#\n# ROTATION GROUP %d '%s',%s fit type '%s'.\n",
1036 g, erotg_names[rotg->eType], buf, erotg_fitnames[rotg->eFittype]);
1038 /* Special type of fitting using the potential minimum. This is
1039 * done for the whole group only, not for the individual slabs. */
1040 if (erotgFitPOT == rotg->eFittype)
1042 fprintf(fp, "# To obtain theta_fit%d, the potential is evaluated for %d angles around theta_ref%d\n", g, rotg->PotAngle_nstep, g);
1043 fprintf(fp, "# The fit angle in the rotation standard outfile is the one with minimal energy E(theta_fit) [kJ/mol].\n");
1047 fprintf(fp, "# Legend for the group %d data columns:\n", g);
1049 print_aligned_short(fp, "time");
1050 print_aligned_short(fp, "grp");
1051 print_aligned(fp, "theta_ref");
1053 if (erotgFitPOT == rotg->eFittype)
1055 /* Output the set of angles around the reference angle */
1056 for (i = 0; i < rotg->PotAngle_nstep; i++)
1058 sprintf(buf, "E(%g)", erg->PotAngleFit->degangle[i]);
1059 print_aligned(fp, buf);
1064 /* Output fit angle for each slab */
1065 print_aligned_short(fp, "slab");
1066 print_aligned_short(fp, "atoms");
1067 print_aligned(fp, "theta_fit");
1068 print_aligned_short(fp, "slab");
1069 print_aligned_short(fp, "atoms");
1070 print_aligned(fp, "theta_fit");
1071 fprintf(fp, " ...");
1083 /* Open torque output file and write some information about it's structure.
1084 * Call on master only */
1085 static FILE *open_torque_out(const char *fn, t_rot *rot)
1092 if (rot->enfrot->Flags & MD_APPENDFILES)
1094 fp = gmx_fio_fopen(fn, "a");
1098 fp = open_output_file(fn, rot->nstsout, "torques");
1100 for (g = 0; g < rot->ngrp; g++)
1102 rotg = &rot->grp[g];
1105 fprintf(fp, "# Rotation group %d (%s), slab distance %f nm.\n", g, erotg_names[rotg->eType], rotg->slab_dist);
1106 fprintf(fp, "# The scalar tau is the torque (kJ/mol) in the direction of the rotation vector.\n");
1107 fprintf(fp, "# To obtain the vectorial torque, multiply tau with\n");
1108 fprintf(fp, "# rot-vec%d %10.3e %10.3e %10.3e\n", g, rotg->vec[XX], rotg->vec[YY], rotg->vec[ZZ]);
1112 fprintf(fp, "# Legend for the following data columns: (tau=torque for that slab):\n");
1114 print_aligned_short(fp, "t");
1115 print_aligned_short(fp, "grp");
1116 print_aligned_short(fp, "slab");
1117 print_aligned(fp, "tau");
1118 print_aligned_short(fp, "slab");
1119 print_aligned(fp, "tau");
1120 fprintf(fp, " ...\n");
1128 static void swap_val(double* vec, int i, int j)
1130 double tmp = vec[j];
1138 static void swap_col(double **mat, int i, int j)
1140 double tmp[3] = {mat[0][j], mat[1][j], mat[2][j]};
1143 mat[0][j] = mat[0][i];
1144 mat[1][j] = mat[1][i];
1145 mat[2][j] = mat[2][i];
1153 /* Eigenvectors are stored in columns of eigen_vec */
1154 static void diagonalize_symmetric(
1162 jacobi(matrix, 3, eigenval, eigen_vec, &n_rot);
1164 /* sort in ascending order */
1165 if (eigenval[0] > eigenval[1])
1167 swap_val(eigenval, 0, 1);
1168 swap_col(eigen_vec, 0, 1);
1170 if (eigenval[1] > eigenval[2])
1172 swap_val(eigenval, 1, 2);
1173 swap_col(eigen_vec, 1, 2);
1175 if (eigenval[0] > eigenval[1])
1177 swap_val(eigenval, 0, 1);
1178 swap_col(eigen_vec, 0, 1);
1183 static void align_with_z(
1184 rvec* s, /* Structure to align */
1189 rvec zet = {0.0, 0.0, 1.0};
1190 rvec rot_axis = {0.0, 0.0, 0.0};
1191 rvec *rotated_str = nullptr;
1197 snew(rotated_str, natoms);
1199 /* Normalize the axis */
1200 ooanorm = 1.0/norm(axis);
1201 svmul(ooanorm, axis, axis);
1203 /* Calculate the angle for the fitting procedure */
1204 cprod(axis, zet, rot_axis);
1205 angle = acos(axis[2]);
1211 /* Calculate the rotation matrix */
1212 calc_rotmat(rot_axis, angle*180.0/M_PI, rotmat);
1214 /* Apply the rotation matrix to s */
1215 for (i = 0; i < natoms; i++)
1217 for (j = 0; j < 3; j++)
1219 for (k = 0; k < 3; k++)
1221 rotated_str[i][j] += rotmat[j][k]*s[i][k];
1226 /* Rewrite the rotated structure to s */
1227 for (i = 0; i < natoms; i++)
1229 for (j = 0; j < 3; j++)
1231 s[i][j] = rotated_str[i][j];
1239 static void calc_correl_matrix(rvec* Xstr, rvec* Ystr, double** Rmat, int natoms)
1244 for (i = 0; i < 3; i++)
1246 for (j = 0; j < 3; j++)
1252 for (i = 0; i < 3; i++)
1254 for (j = 0; j < 3; j++)
1256 for (k = 0; k < natoms; k++)
1258 Rmat[i][j] += Ystr[k][i] * Xstr[k][j];
1265 static void weigh_coords(rvec* str, real* weight, int natoms)
1270 for (i = 0; i < natoms; i++)
1272 for (j = 0; j < 3; j++)
1274 str[i][j] *= std::sqrt(weight[i]);
1280 static real opt_angle_analytic(
1290 rvec *ref_s_1 = nullptr;
1291 rvec *act_s_1 = nullptr;
1293 double **Rmat, **RtR, **eigvec;
1295 double V[3][3], WS[3][3];
1296 double rot_matrix[3][3];
1300 /* Do not change the original coordinates */
1301 snew(ref_s_1, natoms);
1302 snew(act_s_1, natoms);
1303 for (i = 0; i < natoms; i++)
1305 copy_rvec(ref_s[i], ref_s_1[i]);
1306 copy_rvec(act_s[i], act_s_1[i]);
1309 /* Translate the structures to the origin */
1310 shift[XX] = -ref_com[XX];
1311 shift[YY] = -ref_com[YY];
1312 shift[ZZ] = -ref_com[ZZ];
1313 translate_x(ref_s_1, natoms, shift);
1315 shift[XX] = -act_com[XX];
1316 shift[YY] = -act_com[YY];
1317 shift[ZZ] = -act_com[ZZ];
1318 translate_x(act_s_1, natoms, shift);
1320 /* Align rotation axis with z */
1321 align_with_z(ref_s_1, natoms, axis);
1322 align_with_z(act_s_1, natoms, axis);
1324 /* Correlation matrix */
1325 Rmat = allocate_square_matrix(3);
1327 for (i = 0; i < natoms; i++)
1329 ref_s_1[i][2] = 0.0;
1330 act_s_1[i][2] = 0.0;
1333 /* Weight positions with sqrt(weight) */
1334 if (nullptr != weight)
1336 weigh_coords(ref_s_1, weight, natoms);
1337 weigh_coords(act_s_1, weight, natoms);
1340 /* Calculate correlation matrices R=YXt (X=ref_s; Y=act_s) */
1341 calc_correl_matrix(ref_s_1, act_s_1, Rmat, natoms);
1344 RtR = allocate_square_matrix(3);
1345 for (i = 0; i < 3; i++)
1347 for (j = 0; j < 3; j++)
1349 for (k = 0; k < 3; k++)
1351 RtR[i][j] += Rmat[k][i] * Rmat[k][j];
1355 /* Diagonalize RtR */
1357 for (i = 0; i < 3; i++)
1362 diagonalize_symmetric(RtR, eigvec, eigval);
1363 swap_col(eigvec, 0, 1);
1364 swap_col(eigvec, 1, 2);
1365 swap_val(eigval, 0, 1);
1366 swap_val(eigval, 1, 2);
1369 for (i = 0; i < 3; i++)
1371 for (j = 0; j < 3; j++)
1378 for (i = 0; i < 2; i++)
1380 for (j = 0; j < 2; j++)
1382 WS[i][j] = eigvec[i][j] / std::sqrt(eigval[j]);
1386 for (i = 0; i < 3; i++)
1388 for (j = 0; j < 3; j++)
1390 for (k = 0; k < 3; k++)
1392 V[i][j] += Rmat[i][k]*WS[k][j];
1396 free_square_matrix(Rmat, 3);
1398 /* Calculate optimal rotation matrix */
1399 for (i = 0; i < 3; i++)
1401 for (j = 0; j < 3; j++)
1403 rot_matrix[i][j] = 0.0;
1407 for (i = 0; i < 3; i++)
1409 for (j = 0; j < 3; j++)
1411 for (k = 0; k < 3; k++)
1413 rot_matrix[i][j] += eigvec[i][k]*V[j][k];
1417 rot_matrix[2][2] = 1.0;
1419 /* In some cases abs(rot_matrix[0][0]) can be slighly larger
1420 * than unity due to numerical inacurracies. To be able to calculate
1421 * the acos function, we put these values back in range. */
1422 if (rot_matrix[0][0] > 1.0)
1424 rot_matrix[0][0] = 1.0;
1426 else if (rot_matrix[0][0] < -1.0)
1428 rot_matrix[0][0] = -1.0;
1431 /* Determine the optimal rotation angle: */
1432 opt_angle = (-1.0)*acos(rot_matrix[0][0])*180.0/M_PI;
1433 if (rot_matrix[0][1] < 0.0)
1435 opt_angle = (-1.0)*opt_angle;
1438 /* Give back some memory */
1439 free_square_matrix(RtR, 3);
1442 for (i = 0; i < 3; i++)
1448 return (real) opt_angle;
1452 /* Determine angle of the group by RMSD fit to the reference */
1453 /* Not parallelized, call this routine only on the master */
1454 static real flex_fit_angle(t_rotgrp *rotg)
1457 rvec *fitcoords = nullptr;
1458 rvec center; /* Center of positions passed to the fit routine */
1459 real fitangle; /* Angle of the rotation group derived by fitting */
1462 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
1465 erg = rotg->enfrotgrp;
1467 /* Get the center of the rotation group.
1468 * Note, again, erg->xc has been sorted in do_flexible */
1469 get_center(erg->xc, erg->mc_sorted, rotg->nat, center);
1471 /* === Determine the optimal fit angle for the rotation group === */
1472 if (rotg->eFittype == erotgFitNORM)
1474 /* Normalize every position to it's reference length */
1475 for (i = 0; i < rotg->nat; i++)
1477 /* Put the center of the positions into the origin */
1478 rvec_sub(erg->xc[i], center, coord);
1479 /* Determine the scaling factor for the length: */
1480 scal = erg->xc_ref_length[erg->xc_sortind[i]] / norm(coord);
1481 /* Get position, multiply with the scaling factor and save */
1482 svmul(scal, coord, erg->xc_norm[i]);
1484 fitcoords = erg->xc_norm;
1488 fitcoords = erg->xc;
1490 /* From the point of view of the current positions, the reference has rotated
1491 * backwards. Since we output the angle relative to the fixed reference,
1492 * we need the minus sign. */
1493 fitangle = -opt_angle_analytic(erg->xc_ref_sorted, fitcoords, erg->mc_sorted,
1494 rotg->nat, erg->xc_ref_center, center, rotg->vec);
1500 /* Determine actual angle of each slab by RMSD fit to the reference */
1501 /* Not parallelized, call this routine only on the master */
1502 static void flex_fit_angle_perslab(
1509 int i, l, n, islab, ind;
1511 rvec act_center; /* Center of actual positions that are passed to the fit routine */
1512 rvec ref_center; /* Same for the reference positions */
1513 real fitangle; /* Angle of a slab derived from an RMSD fit to
1514 * the reference structure at t=0 */
1516 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
1517 real OOm_av; /* 1/average_mass of a rotation group atom */
1518 real m_rel; /* Relative mass of a rotation group atom */
1521 erg = rotg->enfrotgrp;
1523 /* Average mass of a rotation group atom: */
1524 OOm_av = erg->invmass*rotg->nat;
1526 /**********************************/
1527 /* First collect the data we need */
1528 /**********************************/
1530 /* Collect the data for the individual slabs */
1531 for (n = erg->slab_first; n <= erg->slab_last; n++)
1533 islab = n - erg->slab_first; /* slab index */
1534 sd = &(rotg->enfrotgrp->slab_data[islab]);
1535 sd->nat = erg->lastatom[islab]-erg->firstatom[islab]+1;
1538 /* Loop over the relevant atoms in the slab */
1539 for (l = erg->firstatom[islab]; l <= erg->lastatom[islab]; l++)
1541 /* Current position of this atom: x[ii][XX/YY/ZZ] */
1542 copy_rvec(erg->xc[l], curr_x);
1544 /* The (unrotated) reference position of this atom is copied to ref_x.
1545 * Beware, the xc coords have been sorted in do_flexible */
1546 copy_rvec(erg->xc_ref_sorted[l], ref_x);
1548 /* Save data for doing angular RMSD fit later */
1549 /* Save the current atom position */
1550 copy_rvec(curr_x, sd->x[ind]);
1551 /* Save the corresponding reference position */
1552 copy_rvec(ref_x, sd->ref[ind]);
1554 /* Maybe also mass-weighting was requested. If yes, additionally
1555 * multiply the weights with the relative mass of the atom. If not,
1556 * multiply with unity. */
1557 m_rel = erg->mc_sorted[l]*OOm_av;
1559 /* Save the weight for this atom in this slab */
1560 sd->weight[ind] = gaussian_weight(curr_x, rotg, n) * m_rel;
1562 /* Next atom in this slab */
1567 /******************************/
1568 /* Now do the fit calculation */
1569 /******************************/
1571 fprintf(fp, "%12.3e%6d%12.3f", t, g, degangle);
1573 /* === Now do RMSD fitting for each slab === */
1574 /* We require at least SLAB_MIN_ATOMS in a slab, such that the fit makes sense. */
1575 #define SLAB_MIN_ATOMS 4
1577 for (n = erg->slab_first; n <= erg->slab_last; n++)
1579 islab = n - erg->slab_first; /* slab index */
1580 sd = &(rotg->enfrotgrp->slab_data[islab]);
1581 if (sd->nat >= SLAB_MIN_ATOMS)
1583 /* Get the center of the slabs reference and current positions */
1584 get_center(sd->ref, sd->weight, sd->nat, ref_center);
1585 get_center(sd->x, sd->weight, sd->nat, act_center);
1586 if (rotg->eFittype == erotgFitNORM)
1588 /* Normalize every position to it's reference length
1589 * prior to performing the fit */
1590 for (i = 0; i < sd->nat; i++) /* Center */
1592 rvec_dec(sd->ref[i], ref_center);
1593 rvec_dec(sd->x[i], act_center);
1594 /* Normalize x_i such that it gets the same length as ref_i */
1595 svmul( norm(sd->ref[i])/norm(sd->x[i]), sd->x[i], sd->x[i] );
1597 /* We already subtracted the centers */
1598 clear_rvec(ref_center);
1599 clear_rvec(act_center);
1601 fitangle = -opt_angle_analytic(sd->ref, sd->x, sd->weight, sd->nat,
1602 ref_center, act_center, rotg->vec);
1603 fprintf(fp, "%6d%6d%12.3f", n, sd->nat, fitangle);
1608 #undef SLAB_MIN_ATOMS
1612 /* Shift x with is */
1613 static gmx_inline void shift_single_coord(matrix box, rvec x, const ivec is)
1624 x[XX] += tx*box[XX][XX]+ty*box[YY][XX]+tz*box[ZZ][XX];
1625 x[YY] += ty*box[YY][YY]+tz*box[ZZ][YY];
1626 x[ZZ] += tz*box[ZZ][ZZ];
1630 x[XX] += tx*box[XX][XX];
1631 x[YY] += ty*box[YY][YY];
1632 x[ZZ] += tz*box[ZZ][ZZ];
1637 /* Determine the 'home' slab of this atom which is the
1638 * slab with the highest Gaussian weight of all */
1639 #define round(a) (int)(a+0.5)
1640 static gmx_inline int get_homeslab(
1641 rvec curr_x, /* The position for which the home slab shall be determined */
1642 rvec rotvec, /* The rotation vector */
1643 real slabdist) /* The slab distance */
1648 /* The distance of the atom to the coordinate center (where the
1649 * slab with index 0) is */
1650 dist = iprod(rotvec, curr_x);
1652 return round(dist / slabdist);
1656 /* For a local atom determine the relevant slabs, i.e. slabs in
1657 * which the gaussian is larger than min_gaussian
1659 static int get_single_atom_gaussians(
1666 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
1669 erg = rotg->enfrotgrp;
1671 /* Determine the 'home' slab of this atom: */
1672 homeslab = get_homeslab(curr_x, rotg->vec, rotg->slab_dist);
1674 /* First determine the weight in the atoms home slab: */
1675 g = gaussian_weight(curr_x, rotg, homeslab);
1677 erg->gn_atom[count] = g;
1678 erg->gn_slabind[count] = homeslab;
1682 /* Determine the max slab */
1684 while (g > rotg->min_gaussian)
1687 g = gaussian_weight(curr_x, rotg, slab);
1688 erg->gn_slabind[count] = slab;
1689 erg->gn_atom[count] = g;
1694 /* Determine the min slab */
1699 g = gaussian_weight(curr_x, rotg, slab);
1700 erg->gn_slabind[count] = slab;
1701 erg->gn_atom[count] = g;
1704 while (g > rotg->min_gaussian);
1711 static void flex2_precalc_inner_sum(t_rotgrp *rotg)
1714 rvec xi; /* positions in the i-sum */
1715 rvec xcn, ycn; /* the current and the reference slab centers */
1718 rvec rin; /* Helper variables */
1721 real OOpsii, OOpsiistar;
1722 real sin_rin; /* s_ii.r_ii */
1723 rvec s_in, tmpvec, tmpvec2;
1724 real mi, wi; /* Mass-weighting of the positions */
1726 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
1729 erg = rotg->enfrotgrp;
1730 N_M = rotg->nat * erg->invmass;
1732 /* Loop over all slabs that contain something */
1733 for (n = erg->slab_first; n <= erg->slab_last; n++)
1735 islab = n - erg->slab_first; /* slab index */
1737 /* The current center of this slab is saved in xcn: */
1738 copy_rvec(erg->slab_center[islab], xcn);
1739 /* ... and the reference center in ycn: */
1740 copy_rvec(erg->slab_center_ref[islab+erg->slab_buffer], ycn);
1742 /*** D. Calculate the whole inner sum used for second and third sum */
1743 /* For slab n, we need to loop over all atoms i again. Since we sorted
1744 * the atoms with respect to the rotation vector, we know that it is sufficient
1745 * to calculate from firstatom to lastatom only. All other contributions will
1747 clear_rvec(innersumvec);
1748 for (i = erg->firstatom[islab]; i <= erg->lastatom[islab]; i++)
1750 /* Coordinate xi of this atom */
1751 copy_rvec(erg->xc[i], xi);
1754 gaussian_xi = gaussian_weight(xi, rotg, n);
1755 mi = erg->mc_sorted[i]; /* need the sorted mass here */
1759 copy_rvec(erg->xc_ref_sorted[i], yi0); /* Reference position yi0 */
1760 rvec_sub(yi0, ycn, tmpvec2); /* tmpvec2 = yi0 - ycn */
1761 mvmul(erg->rotmat, tmpvec2, rin); /* rin = Omega.(yi0 - ycn) */
1763 /* Calculate psi_i* and sin */
1764 rvec_sub(xi, xcn, tmpvec2); /* tmpvec2 = xi - xcn */
1765 cprod(rotg->vec, tmpvec2, tmpvec); /* tmpvec = v x (xi - xcn) */
1766 OOpsiistar = norm2(tmpvec)+rotg->eps; /* OOpsii* = 1/psii* = |v x (xi-xcn)|^2 + eps */
1767 OOpsii = norm(tmpvec); /* OOpsii = 1 / psii = |v x (xi - xcn)| */
1769 /* * v x (xi - xcn) */
1770 unitv(tmpvec, s_in); /* sin = ---------------- */
1771 /* |v x (xi - xcn)| */
1773 sin_rin = iprod(s_in, rin); /* sin_rin = sin . rin */
1775 /* Now the whole sum */
1776 fac = OOpsii/OOpsiistar;
1777 svmul(fac, rin, tmpvec);
1778 fac2 = fac*fac*OOpsii;
1779 svmul(fac2*sin_rin, s_in, tmpvec2);
1780 rvec_dec(tmpvec, tmpvec2);
1782 svmul(wi*gaussian_xi*sin_rin, tmpvec, tmpvec2);
1784 rvec_inc(innersumvec, tmpvec2);
1785 } /* now we have the inner sum, used both for sum2 and sum3 */
1787 /* Save it to be used in do_flex2_lowlevel */
1788 copy_rvec(innersumvec, erg->slab_innersumvec[islab]);
1789 } /* END of loop over slabs */
1793 static void flex_precalc_inner_sum(t_rotgrp *rotg)
1796 rvec xi; /* position */
1797 rvec xcn, ycn; /* the current and the reference slab centers */
1798 rvec qin, rin; /* q_i^n and r_i^n */
1801 rvec innersumvec; /* Inner part of sum_n2 */
1802 real gaussian_xi; /* Gaussian weight gn(xi) */
1803 real mi, wi; /* Mass-weighting of the positions */
1806 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
1809 erg = rotg->enfrotgrp;
1810 N_M = rotg->nat * erg->invmass;
1812 /* Loop over all slabs that contain something */
1813 for (n = erg->slab_first; n <= erg->slab_last; n++)
1815 islab = n - erg->slab_first; /* slab index */
1817 /* The current center of this slab is saved in xcn: */
1818 copy_rvec(erg->slab_center[islab], xcn);
1819 /* ... and the reference center in ycn: */
1820 copy_rvec(erg->slab_center_ref[islab+erg->slab_buffer], ycn);
1822 /* For slab n, we need to loop over all atoms i again. Since we sorted
1823 * the atoms with respect to the rotation vector, we know that it is sufficient
1824 * to calculate from firstatom to lastatom only. All other contributions will
1826 clear_rvec(innersumvec);
1827 for (i = erg->firstatom[islab]; i <= erg->lastatom[islab]; i++)
1829 /* Coordinate xi of this atom */
1830 copy_rvec(erg->xc[i], xi);
1833 gaussian_xi = gaussian_weight(xi, rotg, n);
1834 mi = erg->mc_sorted[i]; /* need the sorted mass here */
1837 /* Calculate rin and qin */
1838 rvec_sub(erg->xc_ref_sorted[i], ycn, tmpvec); /* tmpvec = yi0-ycn */
1839 mvmul(erg->rotmat, tmpvec, rin); /* rin = Omega.(yi0 - ycn) */
1840 cprod(rotg->vec, rin, tmpvec); /* tmpvec = v x Omega*(yi0-ycn) */
1842 /* * v x Omega*(yi0-ycn) */
1843 unitv(tmpvec, qin); /* qin = --------------------- */
1844 /* |v x Omega*(yi0-ycn)| */
1847 rvec_sub(xi, xcn, tmpvec); /* tmpvec = xi-xcn */
1848 bin = iprod(qin, tmpvec); /* bin = qin*(xi-xcn) */
1850 svmul(wi*gaussian_xi*bin, qin, tmpvec);
1852 /* Add this contribution to the inner sum: */
1853 rvec_add(innersumvec, tmpvec, innersumvec);
1854 } /* now we have the inner sum vector S^n for this slab */
1855 /* Save it to be used in do_flex_lowlevel */
1856 copy_rvec(innersumvec, erg->slab_innersumvec[islab]);
1861 static real do_flex2_lowlevel(
1863 real sigma, /* The Gaussian width sigma */
1865 gmx_bool bOutstepRot,
1866 gmx_bool bOutstepSlab,
1869 int count, ic, ii, j, m, n, islab, iigrp, ifit;
1870 rvec xj; /* position in the i-sum */
1871 rvec yj0; /* the reference position in the j-sum */
1872 rvec xcn, ycn; /* the current and the reference slab centers */
1873 real V; /* This node's part of the rotation pot. energy */
1874 real gaussian_xj; /* Gaussian weight */
1877 real numerator, fit_numerator;
1878 rvec rjn, fit_rjn; /* Helper variables */
1881 real OOpsij, OOpsijstar;
1882 real OOsigma2; /* 1/(sigma^2) */
1885 rvec sjn, tmpvec, tmpvec2, yj0_ycn;
1886 rvec sum1vec_part, sum1vec, sum2vec_part, sum2vec, sum3vec, sum4vec, innersumvec;
1888 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
1889 real mj, wj; /* Mass-weighting of the positions */
1891 real Wjn; /* g_n(x_j) m_j / Mjn */
1892 gmx_bool bCalcPotFit;
1894 /* To calculate the torque per slab */
1895 rvec slab_force; /* Single force from slab n on one atom */
1896 rvec slab_sum1vec_part;
1897 real slab_sum3part, slab_sum4part;
1898 rvec slab_sum1vec, slab_sum2vec, slab_sum3vec, slab_sum4vec;
1901 erg = rotg->enfrotgrp;
1903 /* Pre-calculate the inner sums, so that we do not have to calculate
1904 * them again for every atom */
1905 flex2_precalc_inner_sum(rotg);
1907 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT == rotg->eFittype);
1909 /********************************************************/
1910 /* Main loop over all local atoms of the rotation group */
1911 /********************************************************/
1912 N_M = rotg->nat * erg->invmass;
1914 OOsigma2 = 1.0 / (sigma*sigma);
1915 for (j = 0; j < erg->nat_loc; j++)
1917 /* Local index of a rotation group atom */
1918 ii = erg->ind_loc[j];
1919 /* Position of this atom in the collective array */
1920 iigrp = erg->xc_ref_ind[j];
1921 /* Mass-weighting */
1922 mj = erg->mc[iigrp]; /* need the unsorted mass here */
1925 /* Current position of this atom: x[ii][XX/YY/ZZ]
1926 * Note that erg->xc_center contains the center of mass in case the flex2-t
1927 * potential was chosen. For the flex2 potential erg->xc_center must be
1929 rvec_sub(x[ii], erg->xc_center, xj);
1931 /* Shift this atom such that it is near its reference */
1932 shift_single_coord(box, xj, erg->xc_shifts[iigrp]);
1934 /* Determine the slabs to loop over, i.e. the ones with contributions
1935 * larger than min_gaussian */
1936 count = get_single_atom_gaussians(xj, rotg);
1938 clear_rvec(sum1vec_part);
1939 clear_rvec(sum2vec_part);
1942 /* Loop over the relevant slabs for this atom */
1943 for (ic = 0; ic < count; ic++)
1945 n = erg->gn_slabind[ic];
1947 /* Get the precomputed Gaussian value of curr_slab for curr_x */
1948 gaussian_xj = erg->gn_atom[ic];
1950 islab = n - erg->slab_first; /* slab index */
1952 /* The (unrotated) reference position of this atom is copied to yj0: */
1953 copy_rvec(rotg->x_ref[iigrp], yj0);
1955 beta = calc_beta(xj, rotg, n);
1957 /* The current center of this slab is saved in xcn: */
1958 copy_rvec(erg->slab_center[islab], xcn);
1959 /* ... and the reference center in ycn: */
1960 copy_rvec(erg->slab_center_ref[islab+erg->slab_buffer], ycn);
1962 rvec_sub(yj0, ycn, yj0_ycn); /* yj0_ycn = yj0 - ycn */
1965 mvmul(erg->rotmat, yj0_ycn, rjn); /* rjn = Omega.(yj0 - ycn) */
1967 /* Subtract the slab center from xj */
1968 rvec_sub(xj, xcn, tmpvec2); /* tmpvec2 = xj - xcn */
1970 /* In rare cases, when an atom position coincides with a slab center
1971 * (tmpvec2 == 0) we cannot compute the vector product for sjn.
1972 * However, since the atom is located directly on the pivot, this
1973 * slab's contribution to the force on that atom will be zero
1974 * anyway. Therefore, we directly move on to the next slab. */
1975 if (gmx_numzero(norm(tmpvec2))) /* 0 == norm(xj - xcn) */
1981 cprod(rotg->vec, tmpvec2, tmpvec); /* tmpvec = v x (xj - xcn) */
1983 OOpsijstar = norm2(tmpvec)+rotg->eps; /* OOpsij* = 1/psij* = |v x (xj-xcn)|^2 + eps */
1985 numerator = gmx::square(iprod(tmpvec, rjn));
1987 /*********************************/
1988 /* Add to the rotation potential */
1989 /*********************************/
1990 V += 0.5*rotg->k*wj*gaussian_xj*numerator/OOpsijstar;
1992 /* If requested, also calculate the potential for a set of angles
1993 * near the current reference angle */
1996 for (ifit = 0; ifit < rotg->PotAngle_nstep; ifit++)
1998 mvmul(erg->PotAngleFit->rotmat[ifit], yj0_ycn, fit_rjn);
1999 fit_numerator = gmx::square(iprod(tmpvec, fit_rjn));
2000 erg->PotAngleFit->V[ifit] += 0.5*rotg->k*wj*gaussian_xj*fit_numerator/OOpsijstar;
2004 /*************************************/
2005 /* Now calculate the force on atom j */
2006 /*************************************/
2008 OOpsij = norm(tmpvec); /* OOpsij = 1 / psij = |v x (xj - xcn)| */
2010 /* * v x (xj - xcn) */
2011 unitv(tmpvec, sjn); /* sjn = ---------------- */
2012 /* |v x (xj - xcn)| */
2014 sjn_rjn = iprod(sjn, rjn); /* sjn_rjn = sjn . rjn */
2017 /*** A. Calculate the first of the four sum terms: ****************/
2018 fac = OOpsij/OOpsijstar;
2019 svmul(fac, rjn, tmpvec);
2020 fac2 = fac*fac*OOpsij;
2021 svmul(fac2*sjn_rjn, sjn, tmpvec2);
2022 rvec_dec(tmpvec, tmpvec2);
2023 fac2 = wj*gaussian_xj; /* also needed for sum4 */
2024 svmul(fac2*sjn_rjn, tmpvec, slab_sum1vec_part);
2025 /********************/
2026 /*** Add to sum1: ***/
2027 /********************/
2028 rvec_inc(sum1vec_part, slab_sum1vec_part); /* sum1 still needs to vector multiplied with v */
2030 /*** B. Calculate the forth of the four sum terms: ****************/
2031 betasigpsi = beta*OOsigma2*OOpsij; /* this is also needed for sum3 */
2032 /********************/
2033 /*** Add to sum4: ***/
2034 /********************/
2035 slab_sum4part = fac2*betasigpsi*fac*sjn_rjn*sjn_rjn; /* Note that fac is still valid from above */
2036 sum4 += slab_sum4part;
2038 /*** C. Calculate Wjn for second and third sum */
2039 /* Note that we can safely divide by slab_weights since we check in
2040 * get_slab_centers that it is non-zero. */
2041 Wjn = gaussian_xj*mj/erg->slab_weights[islab];
2043 /* We already have precalculated the inner sum for slab n */
2044 copy_rvec(erg->slab_innersumvec[islab], innersumvec);
2046 /* Weigh the inner sum vector with Wjn */
2047 svmul(Wjn, innersumvec, innersumvec);
2049 /*** E. Calculate the second of the four sum terms: */
2050 /********************/
2051 /*** Add to sum2: ***/
2052 /********************/
2053 rvec_inc(sum2vec_part, innersumvec); /* sum2 still needs to be vector crossproduct'ed with v */
2055 /*** F. Calculate the third of the four sum terms: */
2056 slab_sum3part = betasigpsi * iprod(sjn, innersumvec);
2057 sum3 += slab_sum3part; /* still needs to be multiplied with v */
2059 /*** G. Calculate the torque on the local slab's axis: */
2063 cprod(slab_sum1vec_part, rotg->vec, slab_sum1vec);
2065 cprod(innersumvec, rotg->vec, slab_sum2vec);
2067 svmul(slab_sum3part, rotg->vec, slab_sum3vec);
2069 svmul(slab_sum4part, rotg->vec, slab_sum4vec);
2071 /* The force on atom ii from slab n only: */
2072 for (m = 0; m < DIM; m++)
2074 slab_force[m] = rotg->k * (-slab_sum1vec[m] + slab_sum2vec[m] - slab_sum3vec[m] + 0.5*slab_sum4vec[m]);
2077 erg->slab_torque_v[islab] += torque(rotg->vec, slab_force, xj, xcn);
2079 } /* END of loop over slabs */
2081 /* Construct the four individual parts of the vector sum: */
2082 cprod(sum1vec_part, rotg->vec, sum1vec); /* sum1vec = { } x v */
2083 cprod(sum2vec_part, rotg->vec, sum2vec); /* sum2vec = { } x v */
2084 svmul(sum3, rotg->vec, sum3vec); /* sum3vec = { } . v */
2085 svmul(sum4, rotg->vec, sum4vec); /* sum4vec = { } . v */
2087 /* Store the additional force so that it can be added to the force
2088 * array after the normal forces have been evaluated */
2089 for (m = 0; m < DIM; m++)
2091 erg->f_rot_loc[j][m] = rotg->k * (-sum1vec[m] + sum2vec[m] - sum3vec[m] + 0.5*sum4vec[m]);
2095 fprintf(stderr, "sum1: %15.8f %15.8f %15.8f\n", -rotg->k*sum1vec[XX], -rotg->k*sum1vec[YY], -rotg->k*sum1vec[ZZ]);
2096 fprintf(stderr, "sum2: %15.8f %15.8f %15.8f\n", rotg->k*sum2vec[XX], rotg->k*sum2vec[YY], rotg->k*sum2vec[ZZ]);
2097 fprintf(stderr, "sum3: %15.8f %15.8f %15.8f\n", -rotg->k*sum3vec[XX], -rotg->k*sum3vec[YY], -rotg->k*sum3vec[ZZ]);
2098 fprintf(stderr, "sum4: %15.8f %15.8f %15.8f\n", 0.5*rotg->k*sum4vec[XX], 0.5*rotg->k*sum4vec[YY], 0.5*rotg->k*sum4vec[ZZ]);
2103 } /* END of loop over local atoms */
2109 static real do_flex_lowlevel(
2111 real sigma, /* The Gaussian width sigma */
2113 gmx_bool bOutstepRot,
2114 gmx_bool bOutstepSlab,
2117 int count, ic, ifit, ii, j, m, n, islab, iigrp;
2118 rvec xj, yj0; /* current and reference position */
2119 rvec xcn, ycn; /* the current and the reference slab centers */
2120 rvec yj0_ycn; /* yj0 - ycn */
2121 rvec xj_xcn; /* xj - xcn */
2122 rvec qjn, fit_qjn; /* q_i^n */
2123 rvec sum_n1, sum_n2; /* Two contributions to the rotation force */
2124 rvec innersumvec; /* Inner part of sum_n2 */
2126 rvec force_n; /* Single force from slab n on one atom */
2127 rvec force_n1, force_n2; /* First and second part of force_n */
2128 rvec tmpvec, tmpvec2, tmp_f; /* Helper variables */
2129 real V; /* The rotation potential energy */
2130 real OOsigma2; /* 1/(sigma^2) */
2131 real beta; /* beta_n(xj) */
2132 real bjn, fit_bjn; /* b_j^n */
2133 real gaussian_xj; /* Gaussian weight gn(xj) */
2134 real betan_xj_sigma2;
2135 real mj, wj; /* Mass-weighting of the positions */
2137 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
2138 gmx_bool bCalcPotFit;
2141 erg = rotg->enfrotgrp;
2143 /* Pre-calculate the inner sums, so that we do not have to calculate
2144 * them again for every atom */
2145 flex_precalc_inner_sum(rotg);
2147 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT == rotg->eFittype);
2149 /********************************************************/
2150 /* Main loop over all local atoms of the rotation group */
2151 /********************************************************/
2152 OOsigma2 = 1.0/(sigma*sigma);
2153 N_M = rotg->nat * erg->invmass;
2155 for (j = 0; j < erg->nat_loc; j++)
2157 /* Local index of a rotation group atom */
2158 ii = erg->ind_loc[j];
2159 /* Position of this atom in the collective array */
2160 iigrp = erg->xc_ref_ind[j];
2161 /* Mass-weighting */
2162 mj = erg->mc[iigrp]; /* need the unsorted mass here */
2165 /* Current position of this atom: x[ii][XX/YY/ZZ]
2166 * Note that erg->xc_center contains the center of mass in case the flex-t
2167 * potential was chosen. For the flex potential erg->xc_center must be
2169 rvec_sub(x[ii], erg->xc_center, xj);
2171 /* Shift this atom such that it is near its reference */
2172 shift_single_coord(box, xj, erg->xc_shifts[iigrp]);
2174 /* Determine the slabs to loop over, i.e. the ones with contributions
2175 * larger than min_gaussian */
2176 count = get_single_atom_gaussians(xj, rotg);
2181 /* Loop over the relevant slabs for this atom */
2182 for (ic = 0; ic < count; ic++)
2184 n = erg->gn_slabind[ic];
2186 /* Get the precomputed Gaussian for xj in slab n */
2187 gaussian_xj = erg->gn_atom[ic];
2189 islab = n - erg->slab_first; /* slab index */
2191 /* The (unrotated) reference position of this atom is saved in yj0: */
2192 copy_rvec(rotg->x_ref[iigrp], yj0);
2194 beta = calc_beta(xj, rotg, n);
2196 /* The current center of this slab is saved in xcn: */
2197 copy_rvec(erg->slab_center[islab], xcn);
2198 /* ... and the reference center in ycn: */
2199 copy_rvec(erg->slab_center_ref[islab+erg->slab_buffer], ycn);
2201 rvec_sub(yj0, ycn, yj0_ycn); /* yj0_ycn = yj0 - ycn */
2203 /* In rare cases, when an atom position coincides with a reference slab
2204 * center (yj0_ycn == 0) we cannot compute the normal vector qjn.
2205 * However, since the atom is located directly on the pivot, this
2206 * slab's contribution to the force on that atom will be zero
2207 * anyway. Therefore, we directly move on to the next slab. */
2208 if (gmx_numzero(norm(yj0_ycn))) /* 0 == norm(yj0 - ycn) */
2214 mvmul(erg->rotmat, yj0_ycn, tmpvec2); /* tmpvec2= Omega.(yj0-ycn) */
2216 /* Subtract the slab center from xj */
2217 rvec_sub(xj, xcn, xj_xcn); /* xj_xcn = xj - xcn */
2220 cprod(rotg->vec, tmpvec2, tmpvec); /* tmpvec= v x Omega.(yj0-ycn) */
2222 /* * v x Omega.(yj0-ycn) */
2223 unitv(tmpvec, qjn); /* qjn = --------------------- */
2224 /* |v x Omega.(yj0-ycn)| */
2226 bjn = iprod(qjn, xj_xcn); /* bjn = qjn * (xj - xcn) */
2228 /*********************************/
2229 /* Add to the rotation potential */
2230 /*********************************/
2231 V += 0.5*rotg->k*wj*gaussian_xj*gmx::square(bjn);
2233 /* If requested, also calculate the potential for a set of angles
2234 * near the current reference angle */
2237 for (ifit = 0; ifit < rotg->PotAngle_nstep; ifit++)
2239 /* As above calculate Omega.(yj0-ycn), now for the other angles */
2240 mvmul(erg->PotAngleFit->rotmat[ifit], yj0_ycn, tmpvec2); /* tmpvec2= Omega.(yj0-ycn) */
2241 /* As above calculate qjn */
2242 cprod(rotg->vec, tmpvec2, tmpvec); /* tmpvec= v x Omega.(yj0-ycn) */
2243 /* * v x Omega.(yj0-ycn) */
2244 unitv(tmpvec, fit_qjn); /* fit_qjn = --------------------- */
2245 /* |v x Omega.(yj0-ycn)| */
2246 fit_bjn = iprod(fit_qjn, xj_xcn); /* fit_bjn = fit_qjn * (xj - xcn) */
2247 /* Add to the rotation potential for this angle */
2248 erg->PotAngleFit->V[ifit] += 0.5*rotg->k*wj*gaussian_xj*gmx::square(fit_bjn);
2252 /****************************************************************/
2253 /* sum_n1 will typically be the main contribution to the force: */
2254 /****************************************************************/
2255 betan_xj_sigma2 = beta*OOsigma2; /* beta_n(xj)/sigma^2 */
2257 /* The next lines calculate
2258 * qjn - (bjn*beta(xj)/(2sigma^2))v */
2259 svmul(bjn*0.5*betan_xj_sigma2, rotg->vec, tmpvec2);
2260 rvec_sub(qjn, tmpvec2, tmpvec);
2262 /* Multiply with gn(xj)*bjn: */
2263 svmul(gaussian_xj*bjn, tmpvec, tmpvec2);
2266 rvec_inc(sum_n1, tmpvec2);
2268 /* We already have precalculated the Sn term for slab n */
2269 copy_rvec(erg->slab_innersumvec[islab], s_n);
2271 svmul(betan_xj_sigma2*iprod(s_n, xj_xcn), rotg->vec, tmpvec); /* tmpvec = ---------- s_n (xj-xcn) */
2274 rvec_sub(s_n, tmpvec, innersumvec);
2276 /* We can safely divide by slab_weights since we check in get_slab_centers
2277 * that it is non-zero. */
2278 svmul(gaussian_xj/erg->slab_weights[islab], innersumvec, innersumvec);
2280 rvec_add(sum_n2, innersumvec, sum_n2);
2282 /* Calculate the torque: */
2285 /* The force on atom ii from slab n only: */
2286 svmul(-rotg->k*wj, tmpvec2, force_n1); /* part 1 */
2287 svmul( rotg->k*mj, innersumvec, force_n2); /* part 2 */
2288 rvec_add(force_n1, force_n2, force_n);
2289 erg->slab_torque_v[islab] += torque(rotg->vec, force_n, xj, xcn);
2291 } /* END of loop over slabs */
2293 /* Put both contributions together: */
2294 svmul(wj, sum_n1, sum_n1);
2295 svmul(mj, sum_n2, sum_n2);
2296 rvec_sub(sum_n2, sum_n1, tmp_f); /* F = -grad V */
2298 /* Store the additional force so that it can be added to the force
2299 * array after the normal forces have been evaluated */
2300 for (m = 0; m < DIM; m++)
2302 erg->f_rot_loc[j][m] = rotg->k*tmp_f[m];
2307 } /* END of loop over local atoms */
2313 static void print_coordinates(t_rotgrp *rotg, rvec x[], matrix box, int step)
2317 static char buf[STRLEN];
2318 static gmx_bool bFirst = 1;
2323 sprintf(buf, "coords%d.txt", cr->nodeid);
2324 fp = fopen(buf, "w");
2328 fprintf(fp, "\nStep %d\n", step);
2329 fprintf(fp, "box: %f %f %f %f %f %f %f %f %f\n",
2330 box[XX][XX], box[XX][YY], box[XX][ZZ],
2331 box[YY][XX], box[YY][YY], box[YY][ZZ],
2332 box[ZZ][XX], box[ZZ][ZZ], box[ZZ][ZZ]);
2333 for (i = 0; i < rotg->nat; i++)
2335 fprintf(fp, "%4d %f %f %f\n", i,
2336 erg->xc[i][XX], erg->xc[i][YY], erg->xc[i][ZZ]);
2344 static int projection_compare(const void *a, const void *b)
2346 sort_along_vec_t *xca, *xcb;
2349 xca = (sort_along_vec_t *)a;
2350 xcb = (sort_along_vec_t *)b;
2352 if (xca->xcproj < xcb->xcproj)
2356 else if (xca->xcproj > xcb->xcproj)
2367 static void sort_collective_coordinates(
2368 t_rotgrp *rotg, /* Rotation group */
2369 sort_along_vec_t *data) /* Buffer for sorting the positions */
2372 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
2375 erg = rotg->enfrotgrp;
2377 /* The projection of the position vector on the rotation vector is
2378 * the relevant value for sorting. Fill the 'data' structure */
2379 for (i = 0; i < rotg->nat; i++)
2381 data[i].xcproj = iprod(erg->xc[i], rotg->vec); /* sort criterium */
2382 data[i].m = erg->mc[i];
2384 copy_rvec(erg->xc[i], data[i].x );
2385 copy_rvec(rotg->x_ref[i], data[i].x_ref);
2387 /* Sort the 'data' structure */
2388 gmx_qsort(data, rotg->nat, sizeof(sort_along_vec_t), projection_compare);
2390 /* Copy back the sorted values */
2391 for (i = 0; i < rotg->nat; i++)
2393 copy_rvec(data[i].x, erg->xc[i] );
2394 copy_rvec(data[i].x_ref, erg->xc_ref_sorted[i]);
2395 erg->mc_sorted[i] = data[i].m;
2396 erg->xc_sortind[i] = data[i].ind;
2401 /* For each slab, get the first and the last index of the sorted atom
2403 static void get_firstlast_atom_per_slab(t_rotgrp *rotg)
2407 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
2410 erg = rotg->enfrotgrp;
2412 /* Find the first atom that needs to enter the calculation for each slab */
2413 n = erg->slab_first; /* slab */
2414 i = 0; /* start with the first atom */
2417 /* Find the first atom that significantly contributes to this slab */
2418 do /* move forward in position until a large enough beta is found */
2420 beta = calc_beta(erg->xc[i], rotg, n);
2423 while ((beta < -erg->max_beta) && (i < rotg->nat));
2425 islab = n - erg->slab_first; /* slab index */
2426 erg->firstatom[islab] = i;
2427 /* Proceed to the next slab */
2430 while (n <= erg->slab_last);
2432 /* Find the last atom for each slab */
2433 n = erg->slab_last; /* start with last slab */
2434 i = rotg->nat-1; /* start with the last atom */
2437 do /* move backward in position until a large enough beta is found */
2439 beta = calc_beta(erg->xc[i], rotg, n);
2442 while ((beta > erg->max_beta) && (i > -1));
2444 islab = n - erg->slab_first; /* slab index */
2445 erg->lastatom[islab] = i;
2446 /* Proceed to the next slab */
2449 while (n >= erg->slab_first);
2453 /* Determine the very first and very last slab that needs to be considered
2454 * For the first slab that needs to be considered, we have to find the smallest
2457 * x_first * v - n*Delta_x <= beta_max
2459 * slab index n, slab distance Delta_x, rotation vector v. For the last slab we
2460 * have to find the largest n that obeys
2462 * x_last * v - n*Delta_x >= -beta_max
2465 static gmx_inline int get_first_slab(
2466 t_rotgrp *rotg, /* The rotation group (inputrec data) */
2467 real max_beta, /* The max_beta value, instead of min_gaussian */
2468 rvec firstatom) /* First atom after sorting along the rotation vector v */
2470 /* Find the first slab for the first atom */
2471 return static_cast<int>(ceil(static_cast<double>((iprod(firstatom, rotg->vec) - max_beta)/rotg->slab_dist)));
2475 static gmx_inline int get_last_slab(
2476 t_rotgrp *rotg, /* The rotation group (inputrec data) */
2477 real max_beta, /* The max_beta value, instead of min_gaussian */
2478 rvec lastatom) /* Last atom along v */
2480 /* Find the last slab for the last atom */
2481 return static_cast<int>(floor(static_cast<double>((iprod(lastatom, rotg->vec) + max_beta)/rotg->slab_dist)));
2485 static void get_firstlast_slab_check(
2486 t_rotgrp *rotg, /* The rotation group (inputrec data) */
2487 t_gmx_enfrotgrp *erg, /* The rotation group (data only accessible in this file) */
2488 rvec firstatom, /* First atom after sorting along the rotation vector v */
2489 rvec lastatom) /* Last atom along v */
2491 erg->slab_first = get_first_slab(rotg, erg->max_beta, firstatom);
2492 erg->slab_last = get_last_slab(rotg, erg->max_beta, lastatom);
2494 /* Calculate the slab buffer size, which changes when slab_first changes */
2495 erg->slab_buffer = erg->slab_first - erg->slab_first_ref;
2497 /* Check whether we have reference data to compare against */
2498 if (erg->slab_first < erg->slab_first_ref)
2500 gmx_fatal(FARGS, "%s No reference data for first slab (n=%d), unable to proceed.",
2501 RotStr, erg->slab_first);
2504 /* Check whether we have reference data to compare against */
2505 if (erg->slab_last > erg->slab_last_ref)
2507 gmx_fatal(FARGS, "%s No reference data for last slab (n=%d), unable to proceed.",
2508 RotStr, erg->slab_last);
2513 /* Enforced rotation with a flexible axis */
2514 static void do_flexible(
2516 gmx_enfrot_t enfrot, /* Other rotation data */
2517 t_rotgrp *rotg, /* The rotation group */
2518 int g, /* Group number */
2519 rvec x[], /* The local positions */
2521 double t, /* Time in picoseconds */
2522 gmx_bool bOutstepRot, /* Output to main rotation output file */
2523 gmx_bool bOutstepSlab) /* Output per-slab data */
2526 real sigma; /* The Gaussian width sigma */
2527 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
2530 erg = rotg->enfrotgrp;
2532 /* Define the sigma value */
2533 sigma = 0.7*rotg->slab_dist;
2535 /* Sort the collective coordinates erg->xc along the rotation vector. This is
2536 * an optimization for the inner loop. */
2537 sort_collective_coordinates(rotg, enfrot->data);
2539 /* Determine the first relevant slab for the first atom and the last
2540 * relevant slab for the last atom */
2541 get_firstlast_slab_check(rotg, erg, erg->xc[0], erg->xc[rotg->nat-1]);
2543 /* Determine for each slab depending on the min_gaussian cutoff criterium,
2544 * a first and a last atom index inbetween stuff needs to be calculated */
2545 get_firstlast_atom_per_slab(rotg);
2547 /* Determine the gaussian-weighted center of positions for all slabs */
2548 get_slab_centers(rotg, erg->xc, erg->mc_sorted, g, t, enfrot->out_slabs, bOutstepSlab, FALSE);
2550 /* Clear the torque per slab from last time step: */
2551 nslabs = erg->slab_last - erg->slab_first + 1;
2552 for (l = 0; l < nslabs; l++)
2554 erg->slab_torque_v[l] = 0.0;
2557 /* Call the rotational forces kernel */
2558 if (rotg->eType == erotgFLEX || rotg->eType == erotgFLEXT)
2560 erg->V = do_flex_lowlevel(rotg, sigma, x, bOutstepRot, bOutstepSlab, box);
2562 else if (rotg->eType == erotgFLEX2 || rotg->eType == erotgFLEX2T)
2564 erg->V = do_flex2_lowlevel(rotg, sigma, x, bOutstepRot, bOutstepSlab, box);
2568 gmx_fatal(FARGS, "Unknown flexible rotation type");
2571 /* Determine angle by RMSD fit to the reference - Let's hope this */
2572 /* only happens once in a while, since this is not parallelized! */
2573 if (bMaster && (erotgFitPOT != rotg->eFittype) )
2577 /* Fit angle of the whole rotation group */
2578 erg->angle_v = flex_fit_angle(rotg);
2582 /* Fit angle of each slab */
2583 flex_fit_angle_perslab(g, rotg, t, erg->degangle, enfrot->out_angles);
2587 /* Lump together the torques from all slabs: */
2588 erg->torque_v = 0.0;
2589 for (l = 0; l < nslabs; l++)
2591 erg->torque_v += erg->slab_torque_v[l];
2596 /* Calculate the angle between reference and actual rotation group atom,
2597 * both projected into a plane perpendicular to the rotation vector: */
2598 static void angle(t_rotgrp *rotg,
2602 real *weight) /* atoms near the rotation axis should count less than atoms far away */
2604 rvec xp, xrp; /* current and reference positions projected on a plane perpendicular to pg->vec */
2608 /* Project x_ref and x into a plane through the origin perpendicular to rot_vec: */
2609 /* Project x_ref: xrp = x_ref - (vec * x_ref) * vec */
2610 svmul(iprod(rotg->vec, x_ref), rotg->vec, dum);
2611 rvec_sub(x_ref, dum, xrp);
2612 /* Project x_act: */
2613 svmul(iprod(rotg->vec, x_act), rotg->vec, dum);
2614 rvec_sub(x_act, dum, xp);
2616 /* Retrieve information about which vector precedes. gmx_angle always
2617 * returns a positive angle. */
2618 cprod(xp, xrp, dum); /* if reference precedes, this is pointing into the same direction as vec */
2620 if (iprod(rotg->vec, dum) >= 0)
2622 *alpha = -gmx_angle(xrp, xp);
2626 *alpha = +gmx_angle(xrp, xp);
2629 /* Also return the weight */
2634 /* Project first vector onto a plane perpendicular to the second vector
2636 * Note that v must be of unit length.
2638 static gmx_inline void project_onto_plane(rvec dr, const rvec v)
2643 svmul(iprod(dr, v), v, tmp); /* tmp = (dr.v)v */
2644 rvec_dec(dr, tmp); /* dr = dr - (dr.v)v */
2648 /* Fixed rotation: The rotation reference group rotates around the v axis. */
2649 /* The atoms of the actual rotation group are attached with imaginary */
2650 /* springs to the reference atoms. */
2651 static void do_fixed(
2652 t_rotgrp *rotg, /* The rotation group */
2653 gmx_bool bOutstepRot, /* Output to main rotation output file */
2654 gmx_bool bOutstepSlab) /* Output per-slab data */
2658 rvec tmp_f; /* Force */
2659 real alpha; /* a single angle between an actual and a reference position */
2660 real weight; /* single weight for a single angle */
2661 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
2662 rvec xi_xc; /* xi - xc */
2663 gmx_bool bCalcPotFit;
2666 /* for mass weighting: */
2667 real wi; /* Mass-weighting of the positions */
2669 real k_wi; /* k times wi */
2674 erg = rotg->enfrotgrp;
2675 bProject = (rotg->eType == erotgPM) || (rotg->eType == erotgPMPF);
2676 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT == rotg->eFittype);
2678 N_M = rotg->nat * erg->invmass;
2680 /* Each process calculates the forces on its local atoms */
2681 for (j = 0; j < erg->nat_loc; j++)
2683 /* Calculate (x_i-x_c) resp. (x_i-u) */
2684 rvec_sub(erg->x_loc_pbc[j], erg->xc_center, xi_xc);
2686 /* Calculate Omega*(y_i-y_c)-(x_i-x_c) */
2687 rvec_sub(erg->xr_loc[j], xi_xc, dr);
2691 project_onto_plane(dr, rotg->vec);
2694 /* Mass-weighting */
2695 wi = N_M*erg->m_loc[j];
2697 /* Store the additional force so that it can be added to the force
2698 * array after the normal forces have been evaluated */
2700 for (m = 0; m < DIM; m++)
2702 tmp_f[m] = k_wi*dr[m];
2703 erg->f_rot_loc[j][m] = tmp_f[m];
2704 erg->V += 0.5*k_wi*gmx::square(dr[m]);
2707 /* If requested, also calculate the potential for a set of angles
2708 * near the current reference angle */
2711 for (ifit = 0; ifit < rotg->PotAngle_nstep; ifit++)
2713 /* Index of this rotation group atom with respect to the whole rotation group */
2714 jj = erg->xc_ref_ind[j];
2716 /* Rotate with the alternative angle. Like rotate_local_reference(),
2717 * just for a single local atom */
2718 mvmul(erg->PotAngleFit->rotmat[ifit], rotg->x_ref[jj], fit_xr_loc); /* fit_xr_loc = Omega*(y_i-y_c) */
2720 /* Calculate Omega*(y_i-y_c)-(x_i-x_c) */
2721 rvec_sub(fit_xr_loc, xi_xc, dr);
2725 project_onto_plane(dr, rotg->vec);
2728 /* Add to the rotation potential for this angle: */
2729 erg->PotAngleFit->V[ifit] += 0.5*k_wi*norm2(dr);
2735 /* Add to the torque of this rotation group */
2736 erg->torque_v += torque(rotg->vec, tmp_f, erg->x_loc_pbc[j], erg->xc_center);
2738 /* Calculate the angle between reference and actual rotation group atom. */
2739 angle(rotg, xi_xc, erg->xr_loc[j], &alpha, &weight); /* angle in rad, weighted */
2740 erg->angle_v += alpha * weight;
2741 erg->weight_v += weight;
2743 /* If you want enforced rotation to contribute to the virial,
2744 * activate the following lines:
2747 Add the rotation contribution to the virial
2748 for(j=0; j<DIM; j++)
2750 vir[j][m] += 0.5*f[ii][j]*dr[m];
2756 } /* end of loop over local rotation group atoms */
2760 /* Calculate the radial motion potential and forces */
2761 static void do_radial_motion(
2762 t_rotgrp *rotg, /* The rotation group */
2763 gmx_bool bOutstepRot, /* Output to main rotation output file */
2764 gmx_bool bOutstepSlab) /* Output per-slab data */
2767 rvec tmp_f; /* Force */
2768 real alpha; /* a single angle between an actual and a reference position */
2769 real weight; /* single weight for a single angle */
2770 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
2771 rvec xj_u; /* xj - u */
2772 rvec tmpvec, fit_tmpvec;
2773 real fac, fac2, sum = 0.0;
2775 gmx_bool bCalcPotFit;
2777 /* For mass weighting: */
2778 real wj; /* Mass-weighting of the positions */
2782 erg = rotg->enfrotgrp;
2783 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT == rotg->eFittype);
2785 N_M = rotg->nat * erg->invmass;
2787 /* Each process calculates the forces on its local atoms */
2788 for (j = 0; j < erg->nat_loc; j++)
2790 /* Calculate (xj-u) */
2791 rvec_sub(erg->x_loc_pbc[j], erg->xc_center, xj_u); /* xj_u = xj-u */
2793 /* Calculate Omega.(yj0-u) */
2794 cprod(rotg->vec, erg->xr_loc[j], tmpvec); /* tmpvec = v x Omega.(yj0-u) */
2796 /* * v x Omega.(yj0-u) */
2797 unitv(tmpvec, pj); /* pj = --------------------- */
2798 /* | v x Omega.(yj0-u) | */
2800 fac = iprod(pj, xj_u); /* fac = pj.(xj-u) */
2803 /* Mass-weighting */
2804 wj = N_M*erg->m_loc[j];
2806 /* Store the additional force so that it can be added to the force
2807 * array after the normal forces have been evaluated */
2808 svmul(-rotg->k*wj*fac, pj, tmp_f);
2809 copy_rvec(tmp_f, erg->f_rot_loc[j]);
2812 /* If requested, also calculate the potential for a set of angles
2813 * near the current reference angle */
2816 for (ifit = 0; ifit < rotg->PotAngle_nstep; ifit++)
2818 /* Index of this rotation group atom with respect to the whole rotation group */
2819 jj = erg->xc_ref_ind[j];
2821 /* Rotate with the alternative angle. Like rotate_local_reference(),
2822 * just for a single local atom */
2823 mvmul(erg->PotAngleFit->rotmat[ifit], rotg->x_ref[jj], fit_tmpvec); /* fit_tmpvec = Omega*(yj0-u) */
2825 /* Calculate Omega.(yj0-u) */
2826 cprod(rotg->vec, fit_tmpvec, tmpvec); /* tmpvec = v x Omega.(yj0-u) */
2827 /* * v x Omega.(yj0-u) */
2828 unitv(tmpvec, pj); /* pj = --------------------- */
2829 /* | v x Omega.(yj0-u) | */
2831 fac = iprod(pj, xj_u); /* fac = pj.(xj-u) */
2834 /* Add to the rotation potential for this angle: */
2835 erg->PotAngleFit->V[ifit] += 0.5*rotg->k*wj*fac2;
2841 /* Add to the torque of this rotation group */
2842 erg->torque_v += torque(rotg->vec, tmp_f, erg->x_loc_pbc[j], erg->xc_center);
2844 /* Calculate the angle between reference and actual rotation group atom. */
2845 angle(rotg, xj_u, erg->xr_loc[j], &alpha, &weight); /* angle in rad, weighted */
2846 erg->angle_v += alpha * weight;
2847 erg->weight_v += weight;
2852 } /* end of loop over local rotation group atoms */
2853 erg->V = 0.5*rotg->k*sum;
2857 /* Calculate the radial motion pivot-free potential and forces */
2858 static void do_radial_motion_pf(
2859 t_rotgrp *rotg, /* The rotation group */
2860 rvec x[], /* The positions */
2861 matrix box, /* The simulation box */
2862 gmx_bool bOutstepRot, /* Output to main rotation output file */
2863 gmx_bool bOutstepSlab) /* Output per-slab data */
2865 int i, ii, iigrp, ifit, j;
2866 rvec xj; /* Current position */
2867 rvec xj_xc; /* xj - xc */
2868 rvec yj0_yc0; /* yj0 - yc0 */
2869 rvec tmp_f; /* Force */
2870 real alpha; /* a single angle between an actual and a reference position */
2871 real weight; /* single weight for a single angle */
2872 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
2873 rvec tmpvec, tmpvec2;
2874 rvec innersumvec; /* Precalculation of the inner sum */
2876 real fac, fac2, V = 0.0;
2878 gmx_bool bCalcPotFit;
2880 /* For mass weighting: */
2881 real mj, wi, wj; /* Mass-weighting of the positions */
2885 erg = rotg->enfrotgrp;
2886 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT == rotg->eFittype);
2888 N_M = rotg->nat * erg->invmass;
2890 /* Get the current center of the rotation group: */
2891 get_center(erg->xc, erg->mc, rotg->nat, erg->xc_center);
2893 /* Precalculate Sum_i [ wi qi.(xi-xc) qi ] which is needed for every single j */
2894 clear_rvec(innersumvec);
2895 for (i = 0; i < rotg->nat; i++)
2897 /* Mass-weighting */
2898 wi = N_M*erg->mc[i];
2900 /* Calculate qi. Note that xc_ref_center has already been subtracted from
2901 * x_ref in init_rot_group.*/
2902 mvmul(erg->rotmat, rotg->x_ref[i], tmpvec); /* tmpvec = Omega.(yi0-yc0) */
2904 cprod(rotg->vec, tmpvec, tmpvec2); /* tmpvec2 = v x Omega.(yi0-yc0) */
2906 /* * v x Omega.(yi0-yc0) */
2907 unitv(tmpvec2, qi); /* qi = ----------------------- */
2908 /* | v x Omega.(yi0-yc0) | */
2910 rvec_sub(erg->xc[i], erg->xc_center, tmpvec); /* tmpvec = xi-xc */
2912 svmul(wi*iprod(qi, tmpvec), qi, tmpvec2);
2914 rvec_inc(innersumvec, tmpvec2);
2916 svmul(rotg->k*erg->invmass, innersumvec, innersumveckM);
2918 /* Each process calculates the forces on its local atoms */
2919 for (j = 0; j < erg->nat_loc; j++)
2921 /* Local index of a rotation group atom */
2922 ii = erg->ind_loc[j];
2923 /* Position of this atom in the collective array */
2924 iigrp = erg->xc_ref_ind[j];
2925 /* Mass-weighting */
2926 mj = erg->mc[iigrp]; /* need the unsorted mass here */
2929 /* Current position of this atom: x[ii][XX/YY/ZZ] */
2930 copy_rvec(x[ii], xj);
2932 /* Shift this atom such that it is near its reference */
2933 shift_single_coord(box, xj, erg->xc_shifts[iigrp]);
2935 /* The (unrotated) reference position is yj0. yc0 has already
2936 * been subtracted in init_rot_group */
2937 copy_rvec(rotg->x_ref[iigrp], yj0_yc0); /* yj0_yc0 = yj0 - yc0 */
2939 /* Calculate Omega.(yj0-yc0) */
2940 mvmul(erg->rotmat, yj0_yc0, tmpvec2); /* tmpvec2 = Omega.(yj0 - yc0) */
2942 cprod(rotg->vec, tmpvec2, tmpvec); /* tmpvec = v x Omega.(yj0-yc0) */
2944 /* * v x Omega.(yj0-yc0) */
2945 unitv(tmpvec, qj); /* qj = ----------------------- */
2946 /* | v x Omega.(yj0-yc0) | */
2948 /* Calculate (xj-xc) */
2949 rvec_sub(xj, erg->xc_center, xj_xc); /* xj_xc = xj-xc */
2951 fac = iprod(qj, xj_xc); /* fac = qj.(xj-xc) */
2954 /* Store the additional force so that it can be added to the force
2955 * array after the normal forces have been evaluated */
2956 svmul(-rotg->k*wj*fac, qj, tmp_f); /* part 1 of force */
2957 svmul(mj, innersumveckM, tmpvec); /* part 2 of force */
2958 rvec_inc(tmp_f, tmpvec);
2959 copy_rvec(tmp_f, erg->f_rot_loc[j]);
2962 /* If requested, also calculate the potential for a set of angles
2963 * near the current reference angle */
2966 for (ifit = 0; ifit < rotg->PotAngle_nstep; ifit++)
2968 /* Rotate with the alternative angle. Like rotate_local_reference(),
2969 * just for a single local atom */
2970 mvmul(erg->PotAngleFit->rotmat[ifit], yj0_yc0, tmpvec2); /* tmpvec2 = Omega*(yj0-yc0) */
2972 /* Calculate Omega.(yj0-u) */
2973 cprod(rotg->vec, tmpvec2, tmpvec); /* tmpvec = v x Omega.(yj0-yc0) */
2974 /* * v x Omega.(yj0-yc0) */
2975 unitv(tmpvec, qj); /* qj = ----------------------- */
2976 /* | v x Omega.(yj0-yc0) | */
2978 fac = iprod(qj, xj_xc); /* fac = qj.(xj-xc) */
2981 /* Add to the rotation potential for this angle: */
2982 erg->PotAngleFit->V[ifit] += 0.5*rotg->k*wj*fac2;
2988 /* Add to the torque of this rotation group */
2989 erg->torque_v += torque(rotg->vec, tmp_f, xj, erg->xc_center);
2991 /* Calculate the angle between reference and actual rotation group atom. */
2992 angle(rotg, xj_xc, yj0_yc0, &alpha, &weight); /* angle in rad, weighted */
2993 erg->angle_v += alpha * weight;
2994 erg->weight_v += weight;
2999 } /* end of loop over local rotation group atoms */
3000 erg->V = 0.5*rotg->k*V;
3004 /* Precalculate the inner sum for the radial motion 2 forces */
3005 static void radial_motion2_precalc_inner_sum(t_rotgrp *rotg, rvec innersumvec)
3008 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
3009 rvec xi_xc; /* xj - xc */
3010 rvec tmpvec, tmpvec2;
3014 rvec v_xi_xc; /* v x (xj - u) */
3015 real psii, psiistar;
3016 real wi; /* Mass-weighting of the positions */
3020 erg = rotg->enfrotgrp;
3021 N_M = rotg->nat * erg->invmass;
3023 /* Loop over the collective set of positions */
3025 for (i = 0; i < rotg->nat; i++)
3027 /* Mass-weighting */
3028 wi = N_M*erg->mc[i];
3030 rvec_sub(erg->xc[i], erg->xc_center, xi_xc); /* xi_xc = xi-xc */
3032 /* Calculate ri. Note that xc_ref_center has already been subtracted from
3033 * x_ref in init_rot_group.*/
3034 mvmul(erg->rotmat, rotg->x_ref[i], ri); /* ri = Omega.(yi0-yc0) */
3036 cprod(rotg->vec, xi_xc, v_xi_xc); /* v_xi_xc = v x (xi-u) */
3038 fac = norm2(v_xi_xc);
3040 psiistar = 1.0/(fac + rotg->eps); /* psiistar = --------------------- */
3041 /* |v x (xi-xc)|^2 + eps */
3043 psii = gmx::invsqrt(fac); /* 1 */
3044 /* psii = ------------- */
3047 svmul(psii, v_xi_xc, si); /* si = psii * (v x (xi-xc) ) */
3049 siri = iprod(si, ri); /* siri = si.ri */
3051 svmul(psiistar/psii, ri, tmpvec);
3052 svmul(psiistar*psiistar/(psii*psii*psii) * siri, si, tmpvec2);
3053 rvec_dec(tmpvec, tmpvec2);
3054 cprod(tmpvec, rotg->vec, tmpvec2);
3056 svmul(wi*siri, tmpvec2, tmpvec);
3058 rvec_inc(sumvec, tmpvec);
3060 svmul(rotg->k*erg->invmass, sumvec, innersumvec);
3064 /* Calculate the radial motion 2 potential and forces */
3065 static void do_radial_motion2(
3066 t_rotgrp *rotg, /* The rotation group */
3067 rvec x[], /* The positions */
3068 matrix box, /* The simulation box */
3069 gmx_bool bOutstepRot, /* Output to main rotation output file */
3070 gmx_bool bOutstepSlab) /* Output per-slab data */
3072 int ii, iigrp, ifit, j;
3073 rvec xj; /* Position */
3074 real alpha; /* a single angle between an actual and a reference position */
3075 real weight; /* single weight for a single angle */
3076 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
3077 rvec xj_u; /* xj - u */
3078 rvec yj0_yc0; /* yj0 -yc0 */
3079 rvec tmpvec, tmpvec2;
3080 real fac, fit_fac, fac2, Vpart = 0.0;
3081 rvec rj, fit_rj, sj;
3083 rvec v_xj_u; /* v x (xj - u) */
3084 real psij, psijstar;
3085 real mj, wj; /* For mass-weighting of the positions */
3089 gmx_bool bCalcPotFit;
3092 erg = rotg->enfrotgrp;
3094 bPF = rotg->eType == erotgRM2PF;
3095 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT == rotg->eFittype);
3098 clear_rvec(yj0_yc0); /* Make the compiler happy */
3100 clear_rvec(innersumvec);
3103 /* For the pivot-free variant we have to use the current center of
3104 * mass of the rotation group instead of the pivot u */
3105 get_center(erg->xc, erg->mc, rotg->nat, erg->xc_center);
3107 /* Also, we precalculate the second term of the forces that is identical
3108 * (up to the weight factor mj) for all forces */
3109 radial_motion2_precalc_inner_sum(rotg, innersumvec);
3112 N_M = rotg->nat * erg->invmass;
3114 /* Each process calculates the forces on its local atoms */
3115 for (j = 0; j < erg->nat_loc; j++)
3119 /* Local index of a rotation group atom */
3120 ii = erg->ind_loc[j];
3121 /* Position of this atom in the collective array */
3122 iigrp = erg->xc_ref_ind[j];
3123 /* Mass-weighting */
3124 mj = erg->mc[iigrp];
3126 /* Current position of this atom: x[ii] */
3127 copy_rvec(x[ii], xj);
3129 /* Shift this atom such that it is near its reference */
3130 shift_single_coord(box, xj, erg->xc_shifts[iigrp]);
3132 /* The (unrotated) reference position is yj0. yc0 has already
3133 * been subtracted in init_rot_group */
3134 copy_rvec(rotg->x_ref[iigrp], yj0_yc0); /* yj0_yc0 = yj0 - yc0 */
3136 /* Calculate Omega.(yj0-yc0) */
3137 mvmul(erg->rotmat, yj0_yc0, rj); /* rj = Omega.(yj0-yc0) */
3142 copy_rvec(erg->x_loc_pbc[j], xj);
3143 copy_rvec(erg->xr_loc[j], rj); /* rj = Omega.(yj0-u) */
3145 /* Mass-weighting */
3148 /* Calculate (xj-u) resp. (xj-xc) */
3149 rvec_sub(xj, erg->xc_center, xj_u); /* xj_u = xj-u */
3151 cprod(rotg->vec, xj_u, v_xj_u); /* v_xj_u = v x (xj-u) */
3153 fac = norm2(v_xj_u);
3155 psijstar = 1.0/(fac + rotg->eps); /* psistar = -------------------- */
3156 /* |v x (xj-u)|^2 + eps */
3158 psij = gmx::invsqrt(fac); /* 1 */
3159 /* psij = ------------ */
3162 svmul(psij, v_xj_u, sj); /* sj = psij * (v x (xj-u) ) */
3164 fac = iprod(v_xj_u, rj); /* fac = (v x (xj-u)).rj */
3167 sjrj = iprod(sj, rj); /* sjrj = sj.rj */
3169 svmul(psijstar/psij, rj, tmpvec);
3170 svmul(psijstar*psijstar/(psij*psij*psij) * sjrj, sj, tmpvec2);
3171 rvec_dec(tmpvec, tmpvec2);
3172 cprod(tmpvec, rotg->vec, tmpvec2);
3174 /* Store the additional force so that it can be added to the force
3175 * array after the normal forces have been evaluated */
3176 svmul(-rotg->k*wj*sjrj, tmpvec2, tmpvec);
3177 svmul(mj, innersumvec, tmpvec2); /* This is != 0 only for the pivot-free variant */
3179 rvec_add(tmpvec2, tmpvec, erg->f_rot_loc[j]);
3180 Vpart += wj*psijstar*fac2;
3182 /* If requested, also calculate the potential for a set of angles
3183 * near the current reference angle */
3186 for (ifit = 0; ifit < rotg->PotAngle_nstep; ifit++)
3190 mvmul(erg->PotAngleFit->rotmat[ifit], yj0_yc0, fit_rj); /* fit_rj = Omega.(yj0-yc0) */
3194 /* Position of this atom in the collective array */
3195 iigrp = erg->xc_ref_ind[j];
3196 /* Rotate with the alternative angle. Like rotate_local_reference(),
3197 * just for a single local atom */
3198 mvmul(erg->PotAngleFit->rotmat[ifit], rotg->x_ref[iigrp], fit_rj); /* fit_rj = Omega*(yj0-u) */
3200 fit_fac = iprod(v_xj_u, fit_rj); /* fac = (v x (xj-u)).fit_rj */
3201 /* Add to the rotation potential for this angle: */
3202 erg->PotAngleFit->V[ifit] += 0.5*rotg->k*wj*psijstar*fit_fac*fit_fac;
3208 /* Add to the torque of this rotation group */
3209 erg->torque_v += torque(rotg->vec, erg->f_rot_loc[j], xj, erg->xc_center);
3211 /* Calculate the angle between reference and actual rotation group atom. */
3212 angle(rotg, xj_u, rj, &alpha, &weight); /* angle in rad, weighted */
3213 erg->angle_v += alpha * weight;
3214 erg->weight_v += weight;
3219 } /* end of loop over local rotation group atoms */
3220 erg->V = 0.5*rotg->k*Vpart;
3224 /* Determine the smallest and largest position vector (with respect to the
3225 * rotation vector) for the reference group */
3226 static void get_firstlast_atom_ref(
3232 real xcproj; /* The projection of a reference position on the
3234 real minproj, maxproj; /* Smallest and largest projection on v */
3238 /* Start with some value */
3239 minproj = iprod(rotg->x_ref[0], rotg->vec);
3242 /* This is just to ensure that it still works if all the atoms of the
3243 * reference structure are situated in a plane perpendicular to the rotation
3246 *lastindex = rotg->nat-1;
3248 /* Loop over all atoms of the reference group,
3249 * project them on the rotation vector to find the extremes */
3250 for (i = 0; i < rotg->nat; i++)
3252 xcproj = iprod(rotg->x_ref[i], rotg->vec);
3253 if (xcproj < minproj)
3258 if (xcproj > maxproj)
3267 /* Allocate memory for the slabs */
3268 static void allocate_slabs(
3274 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
3278 erg = rotg->enfrotgrp;
3280 /* More slabs than are defined for the reference are never needed */
3281 nslabs = erg->slab_last_ref - erg->slab_first_ref + 1;
3283 /* Remember how many we allocated */
3284 erg->nslabs_alloc = nslabs;
3286 if ( (nullptr != fplog) && bVerbose)
3288 fprintf(fplog, "%s allocating memory to store data for %d slabs (rotation group %d).\n",
3291 snew(erg->slab_center, nslabs);
3292 snew(erg->slab_center_ref, nslabs);
3293 snew(erg->slab_weights, nslabs);
3294 snew(erg->slab_torque_v, nslabs);
3295 snew(erg->slab_data, nslabs);
3296 snew(erg->gn_atom, nslabs);
3297 snew(erg->gn_slabind, nslabs);
3298 snew(erg->slab_innersumvec, nslabs);
3299 for (i = 0; i < nslabs; i++)
3301 snew(erg->slab_data[i].x, rotg->nat);
3302 snew(erg->slab_data[i].ref, rotg->nat);
3303 snew(erg->slab_data[i].weight, rotg->nat);
3305 snew(erg->xc_ref_sorted, rotg->nat);
3306 snew(erg->xc_sortind, rotg->nat);
3307 snew(erg->firstatom, nslabs);
3308 snew(erg->lastatom, nslabs);
3312 /* From the extreme positions of the reference group, determine the first
3313 * and last slab of the reference. We can never have more slabs in the real
3314 * simulation than calculated here for the reference.
3316 static void get_firstlast_slab_ref(t_rotgrp *rotg, real mc[], int ref_firstindex, int ref_lastindex)
3318 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
3323 erg = rotg->enfrotgrp;
3324 first = get_first_slab(rotg, erg->max_beta, rotg->x_ref[ref_firstindex]);
3325 last = get_last_slab( rotg, erg->max_beta, rotg->x_ref[ref_lastindex ]);
3327 while (get_slab_weight(first, rotg, rotg->x_ref, mc, &dummy) > WEIGHT_MIN)
3331 erg->slab_first_ref = first+1;
3332 while (get_slab_weight(last, rotg, rotg->x_ref, mc, &dummy) > WEIGHT_MIN)
3336 erg->slab_last_ref = last-1;
3340 /* Special version of copy_rvec:
3341 * During the copy procedure of xcurr to b, the correct PBC image is chosen
3342 * such that the copied vector ends up near its reference position xref */
3343 static gmx_inline void copy_correct_pbc_image(
3344 const rvec xcurr, /* copy vector xcurr ... */
3345 rvec b, /* ... to b ... */
3346 const rvec xref, /* choosing the PBC image such that b ends up near xref */
3355 /* Shortest PBC distance between the atom and its reference */
3356 rvec_sub(xcurr, xref, dx);
3358 /* Determine the shift for this atom */
3360 for (m = npbcdim-1; m >= 0; m--)
3362 while (dx[m] < -0.5*box[m][m])
3364 for (d = 0; d < DIM; d++)
3370 while (dx[m] >= 0.5*box[m][m])
3372 for (d = 0; d < DIM; d++)
3380 /* Apply the shift to the position */
3381 copy_rvec(xcurr, b);
3382 shift_single_coord(box, b, shift);
3386 static void init_rot_group(FILE *fplog, t_commrec *cr, int g, t_rotgrp *rotg,
3387 rvec *x, gmx_mtop_t *mtop, gmx_bool bVerbose, FILE *out_slabs, matrix box,
3388 t_inputrec *ir, gmx_bool bOutputCenters)
3391 rvec coord, xref, *xdum;
3392 gmx_bool bFlex, bColl;
3393 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
3394 int ref_firstindex, ref_lastindex;
3395 real mass, totalmass;
3400 /* Do we have a flexible axis? */
3401 bFlex = ISFLEX(rotg);
3402 /* Do we use a global set of coordinates? */
3403 bColl = ISCOLL(rotg);
3405 erg = rotg->enfrotgrp;
3407 /* Allocate space for collective coordinates if needed */
3410 snew(erg->xc, rotg->nat);
3411 snew(erg->xc_shifts, rotg->nat);
3412 snew(erg->xc_eshifts, rotg->nat);
3413 snew(erg->xc_old, rotg->nat);
3415 if (rotg->eFittype == erotgFitNORM)
3417 snew(erg->xc_ref_length, rotg->nat); /* in case fit type NORM is chosen */
3418 snew(erg->xc_norm, rotg->nat);
3423 snew(erg->xr_loc, rotg->nat);
3424 snew(erg->x_loc_pbc, rotg->nat);
3427 snew(erg->f_rot_loc, rotg->nat);
3428 snew(erg->xc_ref_ind, rotg->nat);
3430 /* Make space for the calculation of the potential at other angles (used
3431 * for fitting only) */
3432 if (erotgFitPOT == rotg->eFittype)
3434 snew(erg->PotAngleFit, 1);
3435 snew(erg->PotAngleFit->degangle, rotg->PotAngle_nstep);
3436 snew(erg->PotAngleFit->V, rotg->PotAngle_nstep);
3437 snew(erg->PotAngleFit->rotmat, rotg->PotAngle_nstep);
3439 /* Get the set of angles around the reference angle */
3440 start = -0.5 * (rotg->PotAngle_nstep - 1)*rotg->PotAngle_step;
3441 for (i = 0; i < rotg->PotAngle_nstep; i++)
3443 erg->PotAngleFit->degangle[i] = start + i*rotg->PotAngle_step;
3448 erg->PotAngleFit = nullptr;
3451 /* xc_ref_ind needs to be set to identity in the serial case */
3454 for (i = 0; i < rotg->nat; i++)
3456 erg->xc_ref_ind[i] = i;
3460 /* Copy the masses so that the center can be determined. For all types of
3461 * enforced rotation, we store the masses in the erg->mc array. */
3462 snew(erg->mc, rotg->nat);
3465 snew(erg->mc_sorted, rotg->nat);
3469 snew(erg->m_loc, rotg->nat);
3473 for (i = 0; i < rotg->nat; i++)
3477 mass = mtopGetAtomMass(mtop, rotg->ind[i], &molb);
3486 erg->invmass = 1.0/totalmass;
3488 /* Set xc_ref_center for any rotation potential */
3489 if ((rotg->eType == erotgISO) || (rotg->eType == erotgPM) || (rotg->eType == erotgRM) || (rotg->eType == erotgRM2))
3491 /* Set the pivot point for the fixed, stationary-axis potentials. This
3492 * won't change during the simulation */
3493 copy_rvec(rotg->pivot, erg->xc_ref_center);
3494 copy_rvec(rotg->pivot, erg->xc_center );
3498 /* Center of the reference positions */
3499 get_center(rotg->x_ref, erg->mc, rotg->nat, erg->xc_ref_center);
3501 /* Center of the actual positions */
3504 snew(xdum, rotg->nat);
3505 for (i = 0; i < rotg->nat; i++)
3508 copy_rvec(x[ii], xdum[i]);
3510 get_center(xdum, erg->mc, rotg->nat, erg->xc_center);
3516 gmx_bcast(sizeof(erg->xc_center), erg->xc_center, cr);
3523 /* Save the original (whole) set of positions in xc_old such that at later
3524 * steps the rotation group can always be made whole again. If the simulation is
3525 * restarted, we compute the starting reference positions (given the time)
3526 * and assume that the correct PBC image of each position is the one nearest
3527 * to the current reference */
3530 /* Calculate the rotation matrix for this angle: */
3531 t_start = ir->init_t + ir->init_step*ir->delta_t;
3532 erg->degangle = rotg->rate * t_start;
3533 calc_rotmat(rotg->vec, erg->degangle, erg->rotmat);
3535 for (i = 0; i < rotg->nat; i++)
3539 /* Subtract pivot, rotate, and add pivot again. This will yield the
3540 * reference position for time t */
3541 rvec_sub(rotg->x_ref[i], erg->xc_ref_center, coord);
3542 mvmul(erg->rotmat, coord, xref);
3543 rvec_inc(xref, erg->xc_ref_center);
3545 copy_correct_pbc_image(x[ii], erg->xc_old[i], xref, box, 3);
3551 gmx_bcast(rotg->nat*sizeof(erg->xc_old[0]), erg->xc_old, cr);
3556 if ( (rotg->eType != erotgFLEX) && (rotg->eType != erotgFLEX2) )
3558 /* Put the reference positions into origin: */
3559 for (i = 0; i < rotg->nat; i++)
3561 rvec_dec(rotg->x_ref[i], erg->xc_ref_center);
3565 /* Enforced rotation with flexible axis */
3568 /* Calculate maximum beta value from minimum gaussian (performance opt.) */
3569 erg->max_beta = calc_beta_max(rotg->min_gaussian, rotg->slab_dist);
3571 /* Determine the smallest and largest coordinate with respect to the rotation vector */
3572 get_firstlast_atom_ref(rotg, &ref_firstindex, &ref_lastindex);
3574 /* From the extreme positions of the reference group, determine the first
3575 * and last slab of the reference. */
3576 get_firstlast_slab_ref(rotg, erg->mc, ref_firstindex, ref_lastindex);
3578 /* Allocate memory for the slabs */
3579 allocate_slabs(rotg, fplog, g, bVerbose);
3581 /* Flexible rotation: determine the reference centers for the rest of the simulation */
3582 erg->slab_first = erg->slab_first_ref;
3583 erg->slab_last = erg->slab_last_ref;
3584 get_slab_centers(rotg, rotg->x_ref, erg->mc, g, -1, out_slabs, bOutputCenters, TRUE);
3586 /* Length of each x_rotref vector from center (needed if fit routine NORM is chosen): */
3587 if (rotg->eFittype == erotgFitNORM)
3589 for (i = 0; i < rotg->nat; i++)
3591 rvec_sub(rotg->x_ref[i], erg->xc_ref_center, coord);
3592 erg->xc_ref_length[i] = norm(coord);
3599 extern void dd_make_local_rotation_groups(gmx_domdec_t *dd, t_rot *rot)
3604 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
3608 for (g = 0; g < rot->ngrp; g++)
3610 rotg = &rot->grp[g];
3611 erg = rotg->enfrotgrp;
3614 dd_make_local_group_indices(ga2la, rotg->nat, rotg->ind,
3615 &erg->nat_loc, &erg->ind_loc, &erg->nalloc_loc, erg->xc_ref_ind);
3620 /* Calculate the size of the MPI buffer needed in reduce_output() */
3621 static int calc_mpi_bufsize(t_rot *rot)
3624 int count_group, count_total;
3626 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
3630 for (g = 0; g < rot->ngrp; g++)
3632 rotg = &rot->grp[g];
3633 erg = rotg->enfrotgrp;
3635 /* Count the items that are transferred for this group: */
3636 count_group = 4; /* V, torque, angle, weight */
3638 /* Add the maximum number of slabs for flexible groups */
3641 count_group += erg->slab_last_ref - erg->slab_first_ref + 1;
3644 /* Add space for the potentials at different angles: */
3645 if (erotgFitPOT == rotg->eFittype)
3647 count_group += rotg->PotAngle_nstep;
3650 /* Add to the total number: */
3651 count_total += count_group;
3658 extern void init_rot(FILE *fplog, t_inputrec *ir, int nfile, const t_filenm fnm[],
3659 t_commrec *cr, rvec *x, matrix box, gmx_mtop_t *mtop, const gmx_output_env_t *oenv,
3660 gmx_bool bVerbose, unsigned long Flags)
3665 int nat_max = 0; /* Size of biggest rotation group */
3666 gmx_enfrot_t er; /* Pointer to the enforced rotation buffer variables */
3667 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
3668 rvec *x_pbc = nullptr; /* Space for the pbc-correct atom positions */
3671 if (MASTER(cr) && bVerbose)
3673 fprintf(stdout, "%s Initializing ...\n", RotStr);
3677 snew(rot->enfrot, 1);
3681 /* When appending, skip first output to avoid duplicate entries in the data files */
3682 if (er->Flags & MD_APPENDFILES)
3691 if (MASTER(cr) && er->bOut)
3693 please_cite(fplog, "Kutzner2011");
3696 /* Output every step for reruns */
3697 if (er->Flags & MD_RERUN)
3699 if (nullptr != fplog)
3701 fprintf(fplog, "%s rerun - will write rotation output every available step.\n", RotStr);
3707 er->out_slabs = nullptr;
3708 if (MASTER(cr) && HaveFlexibleGroups(rot) )
3710 er->out_slabs = open_slab_out(opt2fn("-rs", nfile, fnm), rot);
3715 /* Remove pbc, make molecule whole.
3716 * When ir->bContinuation=TRUE this has already been done, but ok. */
3717 snew(x_pbc, mtop->natoms);
3718 copy_rvecn(x, x_pbc, 0, mtop->natoms);
3719 do_pbc_first_mtop(nullptr, ir->ePBC, box, mtop, x_pbc);
3720 /* All molecules will be whole now, but not necessarily in the home box.
3721 * Additionally, if a rotation group consists of more than one molecule
3722 * (e.g. two strands of DNA), each one of them can end up in a different
3723 * periodic box. This is taken care of in init_rot_group. */
3726 for (g = 0; g < rot->ngrp; g++)
3728 rotg = &rot->grp[g];
3730 if (nullptr != fplog)
3732 fprintf(fplog, "%s group %d type '%s'\n", RotStr, g, erotg_names[rotg->eType]);
3737 /* Allocate space for the rotation group's data: */
3738 snew(rotg->enfrotgrp, 1);
3739 erg = rotg->enfrotgrp;
3741 nat_max = std::max(nat_max, rotg->nat);
3746 erg->nalloc_loc = 0;
3747 erg->ind_loc = nullptr;
3751 erg->nat_loc = rotg->nat;
3752 erg->ind_loc = rotg->ind;
3754 init_rot_group(fplog, cr, g, rotg, x_pbc, mtop, bVerbose, er->out_slabs, box, ir,
3755 !(er->Flags & MD_APPENDFILES) ); /* Do not output the reference centers
3756 * again if we are appending */
3760 /* Allocate space for enforced rotation buffer variables */
3761 er->bufsize = nat_max;
3762 snew(er->data, nat_max);
3763 snew(er->xbuf, nat_max);
3764 snew(er->mbuf, nat_max);
3766 /* Buffers for MPI reducing torques, angles, weights (for each group), and V */
3769 er->mpi_bufsize = calc_mpi_bufsize(rot) + 100; /* larger to catch errors */
3770 snew(er->mpi_inbuf, er->mpi_bufsize);
3771 snew(er->mpi_outbuf, er->mpi_bufsize);
3775 er->mpi_bufsize = 0;
3776 er->mpi_inbuf = nullptr;
3777 er->mpi_outbuf = nullptr;
3780 /* Only do I/O on the MASTER */
3781 er->out_angles = nullptr;
3782 er->out_rot = nullptr;
3783 er->out_torque = nullptr;
3786 er->out_rot = open_rot_out(opt2fn("-ro", nfile, fnm), rot, oenv);
3788 if (rot->nstsout > 0)
3790 if (HaveFlexibleGroups(rot) || HavePotFitGroups(rot) )
3792 er->out_angles = open_angles_out(opt2fn("-ra", nfile, fnm), rot);
3794 if (HaveFlexibleGroups(rot) )
3796 er->out_torque = open_torque_out(opt2fn("-rt", nfile, fnm), rot);
3805 extern void finish_rot(t_rot *rot)
3807 gmx_enfrot_t er; /* Pointer to the enforced rotation buffer variables */
3813 gmx_fio_fclose(er->out_rot);
3817 gmx_fio_fclose(er->out_slabs);
3821 gmx_fio_fclose(er->out_angles);
3825 gmx_fio_fclose(er->out_torque);
3830 /* Rotate the local reference positions and store them in
3831 * erg->xr_loc[0...(nat_loc-1)]
3833 * Note that we already subtracted u or y_c from the reference positions
3834 * in init_rot_group().
3836 static void rotate_local_reference(t_rotgrp *rotg)
3838 gmx_enfrotgrp_t erg;
3842 erg = rotg->enfrotgrp;
3844 for (i = 0; i < erg->nat_loc; i++)
3846 /* Index of this rotation group atom with respect to the whole rotation group */
3847 ii = erg->xc_ref_ind[i];
3849 mvmul(erg->rotmat, rotg->x_ref[ii], erg->xr_loc[i]);
3854 /* Select the PBC representation for each local x position and store that
3855 * for later usage. We assume the right PBC image of an x is the one nearest to
3856 * its rotated reference */
3857 static void choose_pbc_image(rvec x[], t_rotgrp *rotg, matrix box, int npbcdim)
3860 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
3864 erg = rotg->enfrotgrp;
3866 for (i = 0; i < erg->nat_loc; i++)
3868 /* Index of a rotation group atom */
3869 ii = erg->ind_loc[i];
3871 /* Get the correctly rotated reference position. The pivot was already
3872 * subtracted in init_rot_group() from the reference positions. Also,
3873 * the reference positions have already been rotated in
3874 * rotate_local_reference(). For the current reference position we thus
3875 * only need to add the pivot again. */
3876 copy_rvec(erg->xr_loc[i], xref);
3877 rvec_inc(xref, erg->xc_ref_center);
3879 copy_correct_pbc_image(x[ii], erg->x_loc_pbc[i], xref, box, npbcdim);
3884 extern void do_rotation(
3891 gmx_wallcycle_t wcycle,
3897 gmx_bool outstep_slab, outstep_rot;
3899 gmx_enfrot_t er; /* Pointer to the enforced rotation buffer variables */
3900 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
3902 t_gmx_potfit *fit = nullptr; /* For fit type 'potential' determine the fit
3903 angle via the potential minimum */
3905 /* Enforced rotation cycle counting: */
3906 gmx_cycles_t cycles_comp; /* Cycles for the enf. rotation computation
3907 only, does not count communication. This
3908 counter is used for load-balancing */
3917 /* When to output in main rotation output file */
3918 outstep_rot = do_per_step(step, rot->nstrout) && er->bOut;
3919 /* When to output per-slab data */
3920 outstep_slab = do_per_step(step, rot->nstsout) && er->bOut;
3922 /* Output time into rotation output file */
3923 if (outstep_rot && MASTER(cr))
3925 fprintf(er->out_rot, "%12.3e", t);
3928 /**************************************************************************/
3929 /* First do ALL the communication! */
3930 for (g = 0; g < rot->ngrp; g++)
3932 rotg = &rot->grp[g];
3933 erg = rotg->enfrotgrp;
3935 /* Do we use a collective (global) set of coordinates? */
3936 bColl = ISCOLL(rotg);
3938 /* Calculate the rotation matrix for this angle: */
3939 erg->degangle = rotg->rate * t;
3940 calc_rotmat(rotg->vec, erg->degangle, erg->rotmat);
3944 /* Transfer the rotation group's positions such that every node has
3945 * all of them. Every node contributes its local positions x and stores
3946 * it in the collective erg->xc array. */
3947 communicate_group_positions(cr, erg->xc, erg->xc_shifts, erg->xc_eshifts, bNS,
3948 x, rotg->nat, erg->nat_loc, erg->ind_loc, erg->xc_ref_ind, erg->xc_old, box);
3952 /* Fill the local masses array;
3953 * this array changes in DD/neighborsearching steps */
3956 for (i = 0; i < erg->nat_loc; i++)
3958 /* Index of local atom w.r.t. the collective rotation group */
3959 ii = erg->xc_ref_ind[i];
3960 erg->m_loc[i] = erg->mc[ii];
3964 /* Calculate Omega*(y_i-y_c) for the local positions */
3965 rotate_local_reference(rotg);
3967 /* Choose the nearest PBC images of the group atoms with respect
3968 * to the rotated reference positions */
3969 choose_pbc_image(x, rotg, box, 3);
3971 /* Get the center of the rotation group */
3972 if ( (rotg->eType == erotgISOPF) || (rotg->eType == erotgPMPF) )
3974 get_center_comm(cr, erg->x_loc_pbc, erg->m_loc, erg->nat_loc, rotg->nat, erg->xc_center);
3978 } /* End of loop over rotation groups */
3980 /**************************************************************************/
3981 /* Done communicating, we can start to count cycles for the load balancing now ... */
3982 cycles_comp = gmx_cycles_read();
3989 for (g = 0; g < rot->ngrp; g++)
3991 rotg = &rot->grp[g];
3992 erg = rotg->enfrotgrp;
3994 if (outstep_rot && MASTER(cr))
3996 fprintf(er->out_rot, "%12.4f", erg->degangle);
3999 /* Calculate angles and rotation matrices for potential fitting: */
4000 if ( (outstep_rot || outstep_slab) && (erotgFitPOT == rotg->eFittype) )
4002 fit = erg->PotAngleFit;
4003 for (i = 0; i < rotg->PotAngle_nstep; i++)
4005 calc_rotmat(rotg->vec, erg->degangle + fit->degangle[i], fit->rotmat[i]);
4007 /* Clear value from last step */
4008 erg->PotAngleFit->V[i] = 0.0;
4012 /* Clear values from last time step */
4014 erg->torque_v = 0.0;
4016 erg->weight_v = 0.0;
4018 switch (rotg->eType)
4024 do_fixed(rotg, outstep_rot, outstep_slab);
4027 do_radial_motion(rotg, outstep_rot, outstep_slab);
4030 do_radial_motion_pf(rotg, x, box, outstep_rot, outstep_slab);
4034 do_radial_motion2(rotg, x, box, outstep_rot, outstep_slab);
4038 /* Subtract the center of the rotation group from the collective positions array
4039 * Also store the center in erg->xc_center since it needs to be subtracted
4040 * in the low level routines from the local coordinates as well */
4041 get_center(erg->xc, erg->mc, rotg->nat, erg->xc_center);
4042 svmul(-1.0, erg->xc_center, transvec);
4043 translate_x(erg->xc, rotg->nat, transvec);
4044 do_flexible(MASTER(cr), er, rotg, g, x, box, t, outstep_rot, outstep_slab);
4048 /* Do NOT subtract the center of mass in the low level routines! */
4049 clear_rvec(erg->xc_center);
4050 do_flexible(MASTER(cr), er, rotg, g, x, box, t, outstep_rot, outstep_slab);
4053 gmx_fatal(FARGS, "No such rotation potential.");
4061 fprintf(stderr, "%s calculation (step %d) took %g seconds.\n", RotStr, step, MPI_Wtime()-t0);
4065 /* Stop the enforced rotation cycle counter and add the computation-only
4066 * cycles to the force cycles for load balancing */
4067 cycles_comp = gmx_cycles_read() - cycles_comp;
4069 if (DOMAINDECOMP(cr) && wcycle)
4071 dd_cycles_add(cr->dd, cycles_comp, ddCyclF);