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, 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.
45 #include "gromacs/legacyheaders/domdec.h"
46 #include "gromacs/utility/smalloc.h"
47 #include "gromacs/legacyheaders/network.h"
48 #include "gromacs/pbcutil/pbc.h"
49 #include "gromacs/legacyheaders/mdrun.h"
50 #include "gromacs/legacyheaders/txtdump.h"
51 #include "gromacs/legacyheaders/names.h"
52 #include "gromacs/topology/mtop_util.h"
53 #include "gromacs/legacyheaders/names.h"
54 #include "gromacs/math/vec.h"
55 #include "gromacs/legacyheaders/gmx_ga2la.h"
56 #include "gromacs/fileio/xvgr.h"
57 #include "gromacs/legacyheaders/copyrite.h"
58 #include "gromacs/legacyheaders/macros.h"
60 #include "gromacs/utility/futil.h"
61 #include "gromacs/fileio/gmxfio.h"
62 #include "gromacs/fileio/trnio.h"
63 #include "gromacs/linearalgebra/nrjac.h"
64 #include "gromacs/timing/cyclecounter.h"
65 #include "gromacs/timing/wallcycle.h"
66 #include "gromacs/utility/qsort_threadsafe.h"
67 #include "gromacs/pulling/pull_rotation.h"
68 #include "gromacs/mdlib/groupcoord.h"
69 #include "gromacs/math/utilities.h"
71 static char *RotStr = {"Enforced rotation:"};
73 /* Set the minimum weight for the determination of the slab centers */
74 #define WEIGHT_MIN (10*GMX_FLOAT_MIN)
76 /* Helper structure for sorting positions along rotation vector */
78 real xcproj; /* Projection of xc on the rotation vector */
79 int ind; /* Index of xc */
81 rvec x; /* Position */
82 rvec x_ref; /* Reference position */
86 /* Enforced rotation / flexible: determine the angle of each slab */
87 typedef struct gmx_slabdata
89 int nat; /* Number of atoms belonging to this slab */
90 rvec *x; /* The positions belonging to this slab. In
91 general, this should be all positions of the
92 whole rotation group, but we leave those away
93 that have a small enough weight */
94 rvec *ref; /* Same for reference */
95 real *weight; /* The weight for each atom */
99 /* Helper structure for potential fitting */
100 typedef struct gmx_potfit
102 real *degangle; /* Set of angles for which the potential is
103 calculated. The optimum fit is determined as
104 the angle for with the potential is minimal */
105 real *V; /* Potential for the different angles */
106 matrix *rotmat; /* Rotation matrix corresponding to the angles */
110 /* Enforced rotation data for all groups */
111 typedef struct gmx_enfrot
113 FILE *out_rot; /* Output file for rotation data */
114 FILE *out_torque; /* Output file for torque data */
115 FILE *out_angles; /* Output file for slab angles for flexible type */
116 FILE *out_slabs; /* Output file for slab centers */
117 int bufsize; /* Allocation size of buf */
118 rvec *xbuf; /* Coordinate buffer variable for sorting */
119 real *mbuf; /* Masses buffer variable for sorting */
120 sort_along_vec_t *data; /* Buffer variable needed for position sorting */
121 real *mpi_inbuf; /* MPI buffer */
122 real *mpi_outbuf; /* MPI buffer */
123 int mpi_bufsize; /* Allocation size of in & outbuf */
124 unsigned long Flags; /* mdrun flags */
125 gmx_bool bOut; /* Used to skip first output when appending to
126 * avoid duplicate entries in rotation outfiles */
130 /* Global enforced rotation data for a single rotation group */
131 typedef struct gmx_enfrotgrp
133 real degangle; /* Rotation angle in degrees */
134 matrix rotmat; /* Rotation matrix */
135 atom_id *ind_loc; /* Local rotation indices */
136 int nat_loc; /* Number of local group atoms */
137 int nalloc_loc; /* Allocation size for ind_loc and weight_loc */
139 real V; /* Rotation potential for this rotation group */
140 rvec *f_rot_loc; /* Array to store the forces on the local atoms
141 resulting from enforced rotation potential */
143 /* Collective coordinates for the whole rotation group */
144 real *xc_ref_length; /* Length of each x_rotref vector after x_rotref
145 has been put into origin */
146 int *xc_ref_ind; /* Position of each local atom in the collective
148 rvec xc_center; /* Center of the rotation group positions, may
150 rvec xc_ref_center; /* dito, for the reference positions */
151 rvec *xc; /* Current (collective) positions */
152 ivec *xc_shifts; /* Current (collective) shifts */
153 ivec *xc_eshifts; /* Extra shifts since last DD step */
154 rvec *xc_old; /* Old (collective) positions */
155 rvec *xc_norm; /* Normalized form of the current positions */
156 rvec *xc_ref_sorted; /* Reference positions (sorted in the same order
157 as xc when sorted) */
158 int *xc_sortind; /* Where is a position found after sorting? */
159 real *mc; /* Collective masses */
161 real invmass; /* one over the total mass of the rotation group */
163 real torque_v; /* Torque in the direction of rotation vector */
164 real angle_v; /* Actual angle of the whole rotation group */
165 /* Fixed rotation only */
166 real weight_v; /* Weights for angle determination */
167 rvec *xr_loc; /* Local reference coords, correctly rotated */
168 rvec *x_loc_pbc; /* Local current coords, correct PBC image */
169 real *m_loc; /* Masses of the current local atoms */
171 /* Flexible rotation only */
172 int nslabs_alloc; /* For this many slabs memory is allocated */
173 int slab_first; /* Lowermost slab for that the calculation needs
174 to be performed at a given time step */
175 int slab_last; /* Uppermost slab ... */
176 int slab_first_ref; /* First slab for which ref. center is stored */
177 int slab_last_ref; /* Last ... */
178 int slab_buffer; /* Slab buffer region around reference slabs */
179 int *firstatom; /* First relevant atom for a slab */
180 int *lastatom; /* Last relevant atom for a slab */
181 rvec *slab_center; /* Gaussian-weighted slab center */
182 rvec *slab_center_ref; /* Gaussian-weighted slab center for the
183 reference positions */
184 real *slab_weights; /* Sum of gaussian weights in a slab */
185 real *slab_torque_v; /* Torque T = r x f for each slab. */
186 /* torque_v = m.v = angular momentum in the
188 real max_beta; /* min_gaussian from inputrec->rotgrp is the
189 minimum value the gaussian must have so that
190 the force is actually evaluated max_beta is
191 just another way to put it */
192 real *gn_atom; /* Precalculated gaussians for a single atom */
193 int *gn_slabind; /* Tells to which slab each precalculated gaussian
195 rvec *slab_innersumvec; /* Inner sum of the flexible2 potential per slab;
196 this is precalculated for optimization reasons */
197 t_gmx_slabdata *slab_data; /* Holds atom positions and gaussian weights
198 of atoms belonging to a slab */
200 /* For potential fits with varying angle: */
201 t_gmx_potfit *PotAngleFit; /* Used for fit type 'potential' */
205 /* Activate output of forces for correctness checks */
206 /* #define PRINT_FORCES */
208 #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]);
209 #define PRINT_POT_TAU if (MASTER(cr)) { \
210 fprintf(stderr, "potential = %15.8f\n" "torque = %15.8f\n", erg->V, erg->torque_v); \
213 #define PRINT_FORCE_J
214 #define PRINT_POT_TAU
217 /* Shortcuts for often used queries */
218 #define ISFLEX(rg) ( (rg->eType == erotgFLEX) || (rg->eType == erotgFLEXT) || (rg->eType == erotgFLEX2) || (rg->eType == erotgFLEX2T) )
219 #define ISCOLL(rg) ( (rg->eType == erotgFLEX) || (rg->eType == erotgFLEXT) || (rg->eType == erotgFLEX2) || (rg->eType == erotgFLEX2T) || (rg->eType == erotgRMPF) || (rg->eType == erotgRM2PF) )
222 /* Does any of the rotation groups use slab decomposition? */
223 static gmx_bool HaveFlexibleGroups(t_rot *rot)
229 for (g = 0; g < rot->ngrp; g++)
242 /* Is for any group the fit angle determined by finding the minimum of the
243 * rotation potential? */
244 static gmx_bool HavePotFitGroups(t_rot *rot)
250 for (g = 0; g < rot->ngrp; g++)
253 if (erotgFitPOT == rotg->eFittype)
263 static double** allocate_square_matrix(int dim)
270 for (i = 0; i < dim; i++)
279 static void free_square_matrix(double** mat, int dim)
284 for (i = 0; i < dim; i++)
292 /* Return the angle for which the potential is minimal */
293 static real get_fitangle(t_rotgrp *rotg, gmx_enfrotgrp_t erg)
296 real fitangle = -999.9;
297 real pot_min = GMX_FLOAT_MAX;
301 fit = erg->PotAngleFit;
303 for (i = 0; i < rotg->PotAngle_nstep; i++)
305 if (fit->V[i] < pot_min)
308 fitangle = fit->degangle[i];
316 /* Reduce potential angle fit data for this group at this time step? */
317 static gmx_inline gmx_bool bPotAngle(t_rot *rot, t_rotgrp *rotg, gmx_int64_t step)
319 return ( (erotgFitPOT == rotg->eFittype) && (do_per_step(step, rot->nstsout) || do_per_step(step, rot->nstrout)) );
322 /* Reduce slab torqe data for this group at this time step? */
323 static gmx_inline gmx_bool bSlabTau(t_rot *rot, t_rotgrp *rotg, gmx_int64_t step)
325 return ( (ISFLEX(rotg)) && do_per_step(step, rot->nstsout) );
328 /* Output rotation energy, torques, etc. for each rotation group */
329 static void reduce_output(t_commrec *cr, t_rot *rot, real t, gmx_int64_t step)
331 int g, i, islab, nslabs = 0;
332 int count; /* MPI element counter */
334 gmx_enfrot_t er; /* Pointer to the enforced rotation buffer variables */
335 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
342 /* Fill the MPI buffer with stuff to reduce. If items are added for reduction
343 * here, the MPI buffer size has to be enlarged also in calc_mpi_bufsize() */
347 for (g = 0; g < rot->ngrp; g++)
350 erg = rotg->enfrotgrp;
351 nslabs = erg->slab_last - erg->slab_first + 1;
352 er->mpi_inbuf[count++] = erg->V;
353 er->mpi_inbuf[count++] = erg->torque_v;
354 er->mpi_inbuf[count++] = erg->angle_v;
355 er->mpi_inbuf[count++] = erg->weight_v; /* weights are not needed for flex types, but this is just a single value */
357 if (bPotAngle(rot, rotg, step))
359 for (i = 0; i < rotg->PotAngle_nstep; i++)
361 er->mpi_inbuf[count++] = erg->PotAngleFit->V[i];
364 if (bSlabTau(rot, rotg, step))
366 for (i = 0; i < nslabs; i++)
368 er->mpi_inbuf[count++] = erg->slab_torque_v[i];
372 if (count > er->mpi_bufsize)
374 gmx_fatal(FARGS, "%s MPI buffer overflow, please report this error.", RotStr);
378 MPI_Reduce(er->mpi_inbuf, er->mpi_outbuf, count, GMX_MPI_REAL, MPI_SUM, MASTERRANK(cr), cr->mpi_comm_mygroup);
381 /* Copy back the reduced data from the buffer on the master */
385 for (g = 0; g < rot->ngrp; g++)
388 erg = rotg->enfrotgrp;
389 nslabs = erg->slab_last - erg->slab_first + 1;
390 erg->V = er->mpi_outbuf[count++];
391 erg->torque_v = er->mpi_outbuf[count++];
392 erg->angle_v = er->mpi_outbuf[count++];
393 erg->weight_v = er->mpi_outbuf[count++];
395 if (bPotAngle(rot, rotg, step))
397 for (i = 0; i < rotg->PotAngle_nstep; i++)
399 erg->PotAngleFit->V[i] = er->mpi_outbuf[count++];
402 if (bSlabTau(rot, rotg, step))
404 for (i = 0; i < nslabs; i++)
406 erg->slab_torque_v[i] = er->mpi_outbuf[count++];
416 /* Angle and torque for each rotation group */
417 for (g = 0; g < rot->ngrp; g++)
420 bFlex = ISFLEX(rotg);
422 erg = rotg->enfrotgrp;
424 /* Output to main rotation output file: */
425 if (do_per_step(step, rot->nstrout) )
427 if (erotgFitPOT == rotg->eFittype)
429 fitangle = get_fitangle(rotg, erg);
435 fitangle = erg->angle_v; /* RMSD fit angle */
439 fitangle = (erg->angle_v/erg->weight_v)*180.0*M_1_PI;
442 fprintf(er->out_rot, "%12.4f", fitangle);
443 fprintf(er->out_rot, "%12.3e", erg->torque_v);
444 fprintf(er->out_rot, "%12.3e", erg->V);
447 if (do_per_step(step, rot->nstsout) )
449 /* Output to torque log file: */
452 fprintf(er->out_torque, "%12.3e%6d", t, g);
453 for (i = erg->slab_first; i <= erg->slab_last; i++)
455 islab = i - erg->slab_first; /* slab index */
456 /* Only output if enough weight is in slab */
457 if (erg->slab_weights[islab] > rotg->min_gaussian)
459 fprintf(er->out_torque, "%6d%12.3e", i, erg->slab_torque_v[islab]);
462 fprintf(er->out_torque, "\n");
465 /* Output to angles log file: */
466 if (erotgFitPOT == rotg->eFittype)
468 fprintf(er->out_angles, "%12.3e%6d%12.4f", t, g, erg->degangle);
469 /* Output energies at a set of angles around the reference angle */
470 for (i = 0; i < rotg->PotAngle_nstep; i++)
472 fprintf(er->out_angles, "%12.3e", erg->PotAngleFit->V[i]);
474 fprintf(er->out_angles, "\n");
478 if (do_per_step(step, rot->nstrout) )
480 fprintf(er->out_rot, "\n");
486 /* Add the forces from enforced rotation potential to the local forces.
487 * Should be called after the SR forces have been evaluated */
488 extern real add_rot_forces(t_rot *rot, rvec f[], t_commrec *cr, gmx_int64_t step, real t)
492 gmx_enfrot_t er; /* Pointer to the enforced rotation buffer variables */
493 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
494 real Vrot = 0.0; /* If more than one rotation group is present, Vrot
495 assembles the local parts from all groups */
500 /* Loop over enforced rotation groups (usually 1, though)
501 * Apply the forces from rotation potentials */
502 for (g = 0; g < rot->ngrp; g++)
505 erg = rotg->enfrotgrp;
506 Vrot += erg->V; /* add the local parts from the nodes */
507 for (l = 0; l < erg->nat_loc; l++)
509 /* Get the right index of the local force */
510 ii = erg->ind_loc[l];
512 rvec_inc(f[ii], erg->f_rot_loc[l]);
516 /* Reduce energy,torque, angles etc. to get the sum values (per rotation group)
517 * on the master and output these values to file. */
518 if ( (do_per_step(step, rot->nstrout) || do_per_step(step, rot->nstsout)) && er->bOut)
520 reduce_output(cr, rot, t, step);
523 /* When appending, er->bOut is FALSE the first time to avoid duplicate entries */
532 /* The Gaussian norm is chosen such that the sum of the gaussian functions
533 * over the slabs is approximately 1.0 everywhere */
534 #define GAUSS_NORM 0.569917543430618
537 /* Calculate the maximum beta that leads to a gaussian larger min_gaussian,
538 * also does some checks
540 static double calc_beta_max(real min_gaussian, real slab_dist)
546 /* Actually the next two checks are already made in grompp */
549 gmx_fatal(FARGS, "Slab distance of flexible rotation groups must be >=0 !");
551 if (min_gaussian <= 0)
553 gmx_fatal(FARGS, "Cutoff value for Gaussian must be > 0. (You requested %f)");
556 /* Define the sigma value */
557 sigma = 0.7*slab_dist;
559 /* Calculate the argument for the logarithm and check that the log() result is negative or 0 */
560 arg = min_gaussian/GAUSS_NORM;
563 gmx_fatal(FARGS, "min_gaussian of flexible rotation groups must be <%g", GAUSS_NORM);
566 return sqrt(-2.0*sigma*sigma*log(min_gaussian/GAUSS_NORM));
570 static gmx_inline real calc_beta(rvec curr_x, t_rotgrp *rotg, int n)
572 return iprod(curr_x, rotg->vec) - rotg->slab_dist * n;
576 static gmx_inline real gaussian_weight(rvec curr_x, t_rotgrp *rotg, int n)
578 const real norm = GAUSS_NORM;
582 /* Define the sigma value */
583 sigma = 0.7*rotg->slab_dist;
584 /* Calculate the Gaussian value of slab n for position curr_x */
585 return norm * exp( -0.5 * sqr( calc_beta(curr_x, rotg, n)/sigma ) );
589 /* Returns the weight in a single slab, also calculates the Gaussian- and mass-
590 * weighted sum of positions for that slab */
591 static real get_slab_weight(int j, t_rotgrp *rotg, rvec xc[], real mc[], rvec *x_weighted_sum)
593 rvec curr_x; /* The position of an atom */
594 rvec curr_x_weighted; /* The gaussian-weighted position */
595 real gaussian; /* A single gaussian weight */
596 real wgauss; /* gaussian times current mass */
597 real slabweight = 0.0; /* The sum of weights in the slab */
599 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
602 erg = rotg->enfrotgrp;
603 clear_rvec(*x_weighted_sum);
606 islab = j - erg->slab_first;
608 /* Loop over all atoms in the rotation group */
609 for (i = 0; i < rotg->nat; i++)
611 copy_rvec(xc[i], curr_x);
612 gaussian = gaussian_weight(curr_x, rotg, j);
613 wgauss = gaussian * mc[i];
614 svmul(wgauss, curr_x, curr_x_weighted);
615 rvec_add(*x_weighted_sum, curr_x_weighted, *x_weighted_sum);
616 slabweight += wgauss;
617 } /* END of loop over rotation group atoms */
623 static void get_slab_centers(
624 t_rotgrp *rotg, /* The rotation group information */
625 rvec *xc, /* The rotation group positions; will
626 typically be enfrotgrp->xc, but at first call
627 it is enfrotgrp->xc_ref */
628 real *mc, /* The masses of the rotation group atoms */
629 int g, /* The number of the rotation group */
630 real time, /* Used for output only */
631 FILE *out_slabs, /* For outputting center per slab information */
632 gmx_bool bOutStep, /* Is this an output step? */
633 gmx_bool bReference) /* If this routine is called from
634 init_rot_group we need to store
635 the reference slab centers */
638 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
641 erg = rotg->enfrotgrp;
643 /* Loop over slabs */
644 for (j = erg->slab_first; j <= erg->slab_last; j++)
646 islab = j - erg->slab_first;
647 erg->slab_weights[islab] = get_slab_weight(j, rotg, xc, mc, &erg->slab_center[islab]);
649 /* We can do the calculations ONLY if there is weight in the slab! */
650 if (erg->slab_weights[islab] > WEIGHT_MIN)
652 svmul(1.0/erg->slab_weights[islab], erg->slab_center[islab], erg->slab_center[islab]);
656 /* We need to check this here, since we divide through slab_weights
657 * in the flexible low-level routines! */
658 gmx_fatal(FARGS, "Not enough weight in slab %d. Slab center cannot be determined!", j);
661 /* At first time step: save the centers of the reference structure */
664 copy_rvec(erg->slab_center[islab], erg->slab_center_ref[islab]);
666 } /* END of loop over slabs */
668 /* Output on the master */
669 if ( (NULL != out_slabs) && bOutStep)
671 fprintf(out_slabs, "%12.3e%6d", time, g);
672 for (j = erg->slab_first; j <= erg->slab_last; j++)
674 islab = j - erg->slab_first;
675 fprintf(out_slabs, "%6d%12.3e%12.3e%12.3e",
676 j, erg->slab_center[islab][XX], erg->slab_center[islab][YY], erg->slab_center[islab][ZZ]);
678 fprintf(out_slabs, "\n");
683 static void calc_rotmat(
685 real degangle, /* Angle alpha of rotation at time t in degrees */
686 matrix rotmat) /* Rotation matrix */
688 real radangle; /* Rotation angle in radians */
689 real cosa; /* cosine alpha */
690 real sina; /* sine alpha */
691 real OMcosa; /* 1 - cos(alpha) */
692 real dumxy, dumxz, dumyz; /* save computations */
693 rvec rot_vec; /* Rotate around rot_vec ... */
696 radangle = degangle * M_PI/180.0;
697 copy_rvec(vec, rot_vec );
699 /* Precompute some variables: */
700 cosa = cos(radangle);
701 sina = sin(radangle);
703 dumxy = rot_vec[XX]*rot_vec[YY]*OMcosa;
704 dumxz = rot_vec[XX]*rot_vec[ZZ]*OMcosa;
705 dumyz = rot_vec[YY]*rot_vec[ZZ]*OMcosa;
707 /* Construct the rotation matrix for this rotation group: */
709 rotmat[XX][XX] = cosa + rot_vec[XX]*rot_vec[XX]*OMcosa;
710 rotmat[YY][XX] = dumxy + rot_vec[ZZ]*sina;
711 rotmat[ZZ][XX] = dumxz - rot_vec[YY]*sina;
713 rotmat[XX][YY] = dumxy - rot_vec[ZZ]*sina;
714 rotmat[YY][YY] = cosa + rot_vec[YY]*rot_vec[YY]*OMcosa;
715 rotmat[ZZ][YY] = dumyz + rot_vec[XX]*sina;
717 rotmat[XX][ZZ] = dumxz + rot_vec[YY]*sina;
718 rotmat[YY][ZZ] = dumyz - rot_vec[XX]*sina;
719 rotmat[ZZ][ZZ] = cosa + rot_vec[ZZ]*rot_vec[ZZ]*OMcosa;
724 for (iii = 0; iii < 3; iii++)
726 for (jjj = 0; jjj < 3; jjj++)
728 fprintf(stderr, " %10.8f ", rotmat[iii][jjj]);
730 fprintf(stderr, "\n");
736 /* Calculates torque on the rotation axis tau = position x force */
737 static gmx_inline real torque(
738 rvec rotvec, /* rotation vector; MUST be normalized! */
739 rvec force, /* force */
740 rvec x, /* position of atom on which the force acts */
741 rvec pivot) /* pivot point of rotation axis */
746 /* Subtract offset */
747 rvec_sub(x, pivot, vectmp);
749 /* position x force */
750 cprod(vectmp, force, tau);
752 /* Return the part of the torque which is parallel to the rotation vector */
753 return iprod(tau, rotvec);
757 /* Right-aligned output of value with standard width */
758 static void print_aligned(FILE *fp, char *str)
760 fprintf(fp, "%12s", str);
764 /* Right-aligned output of value with standard short width */
765 static void print_aligned_short(FILE *fp, char *str)
767 fprintf(fp, "%6s", str);
771 static FILE *open_output_file(const char *fn, int steps, const char what[])
776 fp = gmx_ffopen(fn, "w");
778 fprintf(fp, "# Output of %s is written in intervals of %d time step%s.\n#\n",
779 what, steps, steps > 1 ? "s" : "");
785 /* Open output file for slab center data. Call on master only */
786 static FILE *open_slab_out(const char *fn, t_rot *rot)
793 if (rot->enfrot->Flags & MD_APPENDFILES)
795 fp = gmx_fio_fopen(fn, "a");
799 fp = open_output_file(fn, rot->nstsout, "gaussian weighted slab centers");
801 for (g = 0; g < rot->ngrp; g++)
806 fprintf(fp, "# Rotation group %d (%s), slab distance %f nm, %s.\n",
807 g, erotg_names[rotg->eType], rotg->slab_dist,
808 rotg->bMassW ? "centers of mass" : "geometrical centers");
812 fprintf(fp, "# Reference centers are listed first (t=-1).\n");
813 fprintf(fp, "# The following columns have the syntax:\n");
815 print_aligned_short(fp, "t");
816 print_aligned_short(fp, "grp");
817 /* Print legend for the first two entries only ... */
818 for (i = 0; i < 2; i++)
820 print_aligned_short(fp, "slab");
821 print_aligned(fp, "X center");
822 print_aligned(fp, "Y center");
823 print_aligned(fp, "Z center");
825 fprintf(fp, " ...\n");
833 /* Adds 'buf' to 'str' */
834 static void add_to_string(char **str, char *buf)
839 len = strlen(*str) + strlen(buf) + 1;
845 static void add_to_string_aligned(char **str, char *buf)
847 char buf_aligned[STRLEN];
849 sprintf(buf_aligned, "%12s", buf);
850 add_to_string(str, buf_aligned);
854 /* Open output file and print some general information about the rotation groups.
855 * Call on master only */
856 static FILE *open_rot_out(const char *fn, t_rot *rot, const output_env_t oenv)
861 const char **setname;
862 char buf[50], buf2[75];
863 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
865 char *LegendStr = NULL;
868 if (rot->enfrot->Flags & MD_APPENDFILES)
870 fp = gmx_fio_fopen(fn, "a");
874 fp = xvgropen(fn, "Rotation angles and energy", "Time (ps)", "angles (degrees) and energies (kJ/mol)", oenv);
875 fprintf(fp, "# Output of enforced rotation data is written in intervals of %d time step%s.\n#\n", rot->nstrout, rot->nstrout > 1 ? "s" : "");
876 fprintf(fp, "# The scalar tau is the torque (kJ/mol) in the direction of the rotation vector v.\n");
877 fprintf(fp, "# To obtain the vectorial torque, multiply tau with the group's rot_vec.\n");
878 fprintf(fp, "# For flexible groups, tau(t,n) from all slabs n have been summed in a single value tau(t) here.\n");
879 fprintf(fp, "# The torques tau(t,n) are found in the rottorque.log (-rt) output file\n");
881 for (g = 0; g < rot->ngrp; g++)
884 erg = rotg->enfrotgrp;
885 bFlex = ISFLEX(rotg);
888 fprintf(fp, "# ROTATION GROUP %d, potential type '%s':\n", g, erotg_names[rotg->eType]);
889 fprintf(fp, "# rot_massw%d %s\n", g, yesno_names[rotg->bMassW]);
890 fprintf(fp, "# rot_vec%d %12.5e %12.5e %12.5e\n", g, rotg->vec[XX], rotg->vec[YY], rotg->vec[ZZ]);
891 fprintf(fp, "# rot_rate%d %12.5e degrees/ps\n", g, rotg->rate);
892 fprintf(fp, "# rot_k%d %12.5e kJ/(mol*nm^2)\n", g, rotg->k);
893 if (rotg->eType == erotgISO || rotg->eType == erotgPM || rotg->eType == erotgRM || rotg->eType == erotgRM2)
895 fprintf(fp, "# rot_pivot%d %12.5e %12.5e %12.5e nm\n", g, rotg->pivot[XX], rotg->pivot[YY], rotg->pivot[ZZ]);
900 fprintf(fp, "# rot_slab_distance%d %f nm\n", g, rotg->slab_dist);
901 fprintf(fp, "# rot_min_gaussian%d %12.5e\n", g, rotg->min_gaussian);
904 /* Output the centers of the rotation groups for the pivot-free potentials */
905 if ((rotg->eType == erotgISOPF) || (rotg->eType == erotgPMPF) || (rotg->eType == erotgRMPF) || (rotg->eType == erotgRM2PF
906 || (rotg->eType == erotgFLEXT) || (rotg->eType == erotgFLEX2T)) )
908 fprintf(fp, "# ref. grp. %d center %12.5e %12.5e %12.5e\n", g,
909 erg->xc_ref_center[XX], erg->xc_ref_center[YY], erg->xc_ref_center[ZZ]);
911 fprintf(fp, "# grp. %d init.center %12.5e %12.5e %12.5e\n", g,
912 erg->xc_center[XX], erg->xc_center[YY], erg->xc_center[ZZ]);
915 if ( (rotg->eType == erotgRM2) || (rotg->eType == erotgFLEX2) || (rotg->eType == erotgFLEX2T) )
917 fprintf(fp, "# rot_eps%d %12.5e nm^2\n", g, rotg->eps);
919 if (erotgFitPOT == rotg->eFittype)
922 fprintf(fp, "# theta_fit%d is determined by first evaluating the potential for %d angles around theta_ref%d.\n",
923 g, rotg->PotAngle_nstep, g);
924 fprintf(fp, "# The fit angle is the one with the smallest potential. It is given as the deviation\n");
925 fprintf(fp, "# from the reference angle, i.e. if theta_ref=X and theta_fit=Y, then the angle with\n");
926 fprintf(fp, "# minimal value of the potential is X+Y. Angular resolution is %g degrees.\n", rotg->PotAngle_step);
930 /* Print a nice legend */
933 sprintf(buf, "# %6s", "time");
934 add_to_string_aligned(&LegendStr, buf);
937 snew(setname, 4*rot->ngrp);
939 for (g = 0; g < rot->ngrp; g++)
942 sprintf(buf, "theta_ref%d", g);
943 add_to_string_aligned(&LegendStr, buf);
945 sprintf(buf2, "%s (degrees)", buf);
946 setname[nsets] = gmx_strdup(buf2);
949 for (g = 0; g < rot->ngrp; g++)
952 bFlex = ISFLEX(rotg);
954 /* For flexible axis rotation we use RMSD fitting to determine the
955 * actual angle of the rotation group */
956 if (bFlex || erotgFitPOT == rotg->eFittype)
958 sprintf(buf, "theta_fit%d", g);
962 sprintf(buf, "theta_av%d", g);
964 add_to_string_aligned(&LegendStr, buf);
965 sprintf(buf2, "%s (degrees)", buf);
966 setname[nsets] = gmx_strdup(buf2);
969 sprintf(buf, "tau%d", g);
970 add_to_string_aligned(&LegendStr, buf);
971 sprintf(buf2, "%s (kJ/mol)", buf);
972 setname[nsets] = gmx_strdup(buf2);
975 sprintf(buf, "energy%d", g);
976 add_to_string_aligned(&LegendStr, buf);
977 sprintf(buf2, "%s (kJ/mol)", buf);
978 setname[nsets] = gmx_strdup(buf2);
985 xvgr_legend(fp, nsets, setname, oenv);
989 fprintf(fp, "#\n# Legend for the following data columns:\n");
990 fprintf(fp, "%s\n", LegendStr);
1000 /* Call on master only */
1001 static FILE *open_angles_out(const char *fn, t_rot *rot)
1006 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
1010 if (rot->enfrot->Flags & MD_APPENDFILES)
1012 fp = gmx_fio_fopen(fn, "a");
1016 /* Open output file and write some information about it's structure: */
1017 fp = open_output_file(fn, rot->nstsout, "rotation group angles");
1018 fprintf(fp, "# All angles given in degrees, time in ps.\n");
1019 for (g = 0; g < rot->ngrp; g++)
1021 rotg = &rot->grp[g];
1022 erg = rotg->enfrotgrp;
1024 /* Output for this group happens only if potential type is flexible or
1025 * if fit type is potential! */
1026 if (ISFLEX(rotg) || (erotgFitPOT == rotg->eFittype) )
1030 sprintf(buf, " slab distance %f nm, ", rotg->slab_dist);
1037 fprintf(fp, "#\n# ROTATION GROUP %d '%s',%s fit type '%s'.\n",
1038 g, erotg_names[rotg->eType], buf, erotg_fitnames[rotg->eFittype]);
1040 /* Special type of fitting using the potential minimum. This is
1041 * done for the whole group only, not for the individual slabs. */
1042 if (erotgFitPOT == rotg->eFittype)
1044 fprintf(fp, "# To obtain theta_fit%d, the potential is evaluated for %d angles around theta_ref%d\n", g, rotg->PotAngle_nstep, g);
1045 fprintf(fp, "# The fit angle in the rotation standard outfile is the one with minimal energy E(theta_fit) [kJ/mol].\n");
1049 fprintf(fp, "# Legend for the group %d data columns:\n", g);
1051 print_aligned_short(fp, "time");
1052 print_aligned_short(fp, "grp");
1053 print_aligned(fp, "theta_ref");
1055 if (erotgFitPOT == rotg->eFittype)
1057 /* Output the set of angles around the reference angle */
1058 for (i = 0; i < rotg->PotAngle_nstep; i++)
1060 sprintf(buf, "E(%g)", erg->PotAngleFit->degangle[i]);
1061 print_aligned(fp, buf);
1066 /* Output fit angle for each slab */
1067 print_aligned_short(fp, "slab");
1068 print_aligned_short(fp, "atoms");
1069 print_aligned(fp, "theta_fit");
1070 print_aligned_short(fp, "slab");
1071 print_aligned_short(fp, "atoms");
1072 print_aligned(fp, "theta_fit");
1073 fprintf(fp, " ...");
1085 /* Open torque output file and write some information about it's structure.
1086 * Call on master only */
1087 static FILE *open_torque_out(const char *fn, t_rot *rot)
1094 if (rot->enfrot->Flags & MD_APPENDFILES)
1096 fp = gmx_fio_fopen(fn, "a");
1100 fp = open_output_file(fn, rot->nstsout, "torques");
1102 for (g = 0; g < rot->ngrp; g++)
1104 rotg = &rot->grp[g];
1107 fprintf(fp, "# Rotation group %d (%s), slab distance %f nm.\n", g, erotg_names[rotg->eType], rotg->slab_dist);
1108 fprintf(fp, "# The scalar tau is the torque (kJ/mol) in the direction of the rotation vector.\n");
1109 fprintf(fp, "# To obtain the vectorial torque, multiply tau with\n");
1110 fprintf(fp, "# rot_vec%d %10.3e %10.3e %10.3e\n", g, rotg->vec[XX], rotg->vec[YY], rotg->vec[ZZ]);
1114 fprintf(fp, "# Legend for the following data columns: (tau=torque for that slab):\n");
1116 print_aligned_short(fp, "t");
1117 print_aligned_short(fp, "grp");
1118 print_aligned_short(fp, "slab");
1119 print_aligned(fp, "tau");
1120 print_aligned_short(fp, "slab");
1121 print_aligned(fp, "tau");
1122 fprintf(fp, " ...\n");
1130 static void swap_val(double* vec, int i, int j)
1132 double tmp = vec[j];
1140 static void swap_col(double **mat, int i, int j)
1142 double tmp[3] = {mat[0][j], mat[1][j], mat[2][j]};
1145 mat[0][j] = mat[0][i];
1146 mat[1][j] = mat[1][i];
1147 mat[2][j] = mat[2][i];
1155 /* Eigenvectors are stored in columns of eigen_vec */
1156 static void diagonalize_symmetric(
1164 jacobi(matrix, 3, eigenval, eigen_vec, &n_rot);
1166 /* sort in ascending order */
1167 if (eigenval[0] > eigenval[1])
1169 swap_val(eigenval, 0, 1);
1170 swap_col(eigen_vec, 0, 1);
1172 if (eigenval[1] > eigenval[2])
1174 swap_val(eigenval, 1, 2);
1175 swap_col(eigen_vec, 1, 2);
1177 if (eigenval[0] > eigenval[1])
1179 swap_val(eigenval, 0, 1);
1180 swap_col(eigen_vec, 0, 1);
1185 static void align_with_z(
1186 rvec* s, /* Structure to align */
1191 rvec zet = {0.0, 0.0, 1.0};
1192 rvec rot_axis = {0.0, 0.0, 0.0};
1193 rvec *rotated_str = NULL;
1199 snew(rotated_str, natoms);
1201 /* Normalize the axis */
1202 ooanorm = 1.0/norm(axis);
1203 svmul(ooanorm, axis, axis);
1205 /* Calculate the angle for the fitting procedure */
1206 cprod(axis, zet, rot_axis);
1207 angle = acos(axis[2]);
1213 /* Calculate the rotation matrix */
1214 calc_rotmat(rot_axis, angle*180.0/M_PI, rotmat);
1216 /* Apply the rotation matrix to s */
1217 for (i = 0; i < natoms; i++)
1219 for (j = 0; j < 3; j++)
1221 for (k = 0; k < 3; k++)
1223 rotated_str[i][j] += rotmat[j][k]*s[i][k];
1228 /* Rewrite the rotated structure to s */
1229 for (i = 0; i < natoms; i++)
1231 for (j = 0; j < 3; j++)
1233 s[i][j] = rotated_str[i][j];
1241 static void calc_correl_matrix(rvec* Xstr, rvec* Ystr, double** Rmat, int natoms)
1246 for (i = 0; i < 3; i++)
1248 for (j = 0; j < 3; j++)
1254 for (i = 0; i < 3; i++)
1256 for (j = 0; j < 3; j++)
1258 for (k = 0; k < natoms; k++)
1260 Rmat[i][j] += Ystr[k][i] * Xstr[k][j];
1267 static void weigh_coords(rvec* str, real* weight, int natoms)
1272 for (i = 0; i < natoms; i++)
1274 for (j = 0; j < 3; j++)
1276 str[i][j] *= sqrt(weight[i]);
1282 static real opt_angle_analytic(
1292 rvec *ref_s_1 = NULL;
1293 rvec *act_s_1 = NULL;
1295 double **Rmat, **RtR, **eigvec;
1297 double V[3][3], WS[3][3];
1298 double rot_matrix[3][3];
1302 /* Do not change the original coordinates */
1303 snew(ref_s_1, natoms);
1304 snew(act_s_1, natoms);
1305 for (i = 0; i < natoms; i++)
1307 copy_rvec(ref_s[i], ref_s_1[i]);
1308 copy_rvec(act_s[i], act_s_1[i]);
1311 /* Translate the structures to the origin */
1312 shift[XX] = -ref_com[XX];
1313 shift[YY] = -ref_com[YY];
1314 shift[ZZ] = -ref_com[ZZ];
1315 translate_x(ref_s_1, natoms, shift);
1317 shift[XX] = -act_com[XX];
1318 shift[YY] = -act_com[YY];
1319 shift[ZZ] = -act_com[ZZ];
1320 translate_x(act_s_1, natoms, shift);
1322 /* Align rotation axis with z */
1323 align_with_z(ref_s_1, natoms, axis);
1324 align_with_z(act_s_1, natoms, axis);
1326 /* Correlation matrix */
1327 Rmat = allocate_square_matrix(3);
1329 for (i = 0; i < natoms; i++)
1331 ref_s_1[i][2] = 0.0;
1332 act_s_1[i][2] = 0.0;
1335 /* Weight positions with sqrt(weight) */
1338 weigh_coords(ref_s_1, weight, natoms);
1339 weigh_coords(act_s_1, weight, natoms);
1342 /* Calculate correlation matrices R=YXt (X=ref_s; Y=act_s) */
1343 calc_correl_matrix(ref_s_1, act_s_1, Rmat, natoms);
1346 RtR = allocate_square_matrix(3);
1347 for (i = 0; i < 3; i++)
1349 for (j = 0; j < 3; j++)
1351 for (k = 0; k < 3; k++)
1353 RtR[i][j] += Rmat[k][i] * Rmat[k][j];
1357 /* Diagonalize RtR */
1359 for (i = 0; i < 3; i++)
1364 diagonalize_symmetric(RtR, eigvec, eigval);
1365 swap_col(eigvec, 0, 1);
1366 swap_col(eigvec, 1, 2);
1367 swap_val(eigval, 0, 1);
1368 swap_val(eigval, 1, 2);
1371 for (i = 0; i < 3; i++)
1373 for (j = 0; j < 3; j++)
1380 for (i = 0; i < 2; i++)
1382 for (j = 0; j < 2; j++)
1384 WS[i][j] = eigvec[i][j] / sqrt(eigval[j]);
1388 for (i = 0; i < 3; i++)
1390 for (j = 0; j < 3; j++)
1392 for (k = 0; k < 3; k++)
1394 V[i][j] += Rmat[i][k]*WS[k][j];
1398 free_square_matrix(Rmat, 3);
1400 /* Calculate optimal rotation matrix */
1401 for (i = 0; i < 3; i++)
1403 for (j = 0; j < 3; j++)
1405 rot_matrix[i][j] = 0.0;
1409 for (i = 0; i < 3; i++)
1411 for (j = 0; j < 3; j++)
1413 for (k = 0; k < 3; k++)
1415 rot_matrix[i][j] += eigvec[i][k]*V[j][k];
1419 rot_matrix[2][2] = 1.0;
1421 /* In some cases abs(rot_matrix[0][0]) can be slighly larger
1422 * than unity due to numerical inacurracies. To be able to calculate
1423 * the acos function, we put these values back in range. */
1424 if (rot_matrix[0][0] > 1.0)
1426 rot_matrix[0][0] = 1.0;
1428 else if (rot_matrix[0][0] < -1.0)
1430 rot_matrix[0][0] = -1.0;
1433 /* Determine the optimal rotation angle: */
1434 opt_angle = (-1.0)*acos(rot_matrix[0][0])*180.0/M_PI;
1435 if (rot_matrix[0][1] < 0.0)
1437 opt_angle = (-1.0)*opt_angle;
1440 /* Give back some memory */
1441 free_square_matrix(RtR, 3);
1444 for (i = 0; i < 3; i++)
1450 return (real) opt_angle;
1454 /* Determine angle of the group by RMSD fit to the reference */
1455 /* Not parallelized, call this routine only on the master */
1456 static real flex_fit_angle(t_rotgrp *rotg)
1459 rvec *fitcoords = NULL;
1460 rvec center; /* Center of positions passed to the fit routine */
1461 real fitangle; /* Angle of the rotation group derived by fitting */
1464 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
1467 erg = rotg->enfrotgrp;
1469 /* Get the center of the rotation group.
1470 * Note, again, erg->xc has been sorted in do_flexible */
1471 get_center(erg->xc, erg->mc_sorted, rotg->nat, center);
1473 /* === Determine the optimal fit angle for the rotation group === */
1474 if (rotg->eFittype == erotgFitNORM)
1476 /* Normalize every position to it's reference length */
1477 for (i = 0; i < rotg->nat; i++)
1479 /* Put the center of the positions into the origin */
1480 rvec_sub(erg->xc[i], center, coord);
1481 /* Determine the scaling factor for the length: */
1482 scal = erg->xc_ref_length[erg->xc_sortind[i]] / norm(coord);
1483 /* Get position, multiply with the scaling factor and save */
1484 svmul(scal, coord, erg->xc_norm[i]);
1486 fitcoords = erg->xc_norm;
1490 fitcoords = erg->xc;
1492 /* From the point of view of the current positions, the reference has rotated
1493 * backwards. Since we output the angle relative to the fixed reference,
1494 * we need the minus sign. */
1495 fitangle = -opt_angle_analytic(erg->xc_ref_sorted, fitcoords, erg->mc_sorted,
1496 rotg->nat, erg->xc_ref_center, center, rotg->vec);
1502 /* Determine actual angle of each slab by RMSD fit to the reference */
1503 /* Not parallelized, call this routine only on the master */
1504 static void flex_fit_angle_perslab(
1511 int i, l, n, islab, ind;
1513 rvec act_center; /* Center of actual positions that are passed to the fit routine */
1514 rvec ref_center; /* Same for the reference positions */
1515 real fitangle; /* Angle of a slab derived from an RMSD fit to
1516 * the reference structure at t=0 */
1518 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
1519 real OOm_av; /* 1/average_mass of a rotation group atom */
1520 real m_rel; /* Relative mass of a rotation group atom */
1523 erg = rotg->enfrotgrp;
1525 /* Average mass of a rotation group atom: */
1526 OOm_av = erg->invmass*rotg->nat;
1528 /**********************************/
1529 /* First collect the data we need */
1530 /**********************************/
1532 /* Collect the data for the individual slabs */
1533 for (n = erg->slab_first; n <= erg->slab_last; n++)
1535 islab = n - erg->slab_first; /* slab index */
1536 sd = &(rotg->enfrotgrp->slab_data[islab]);
1537 sd->nat = erg->lastatom[islab]-erg->firstatom[islab]+1;
1540 /* Loop over the relevant atoms in the slab */
1541 for (l = erg->firstatom[islab]; l <= erg->lastatom[islab]; l++)
1543 /* Current position of this atom: x[ii][XX/YY/ZZ] */
1544 copy_rvec(erg->xc[l], curr_x);
1546 /* The (unrotated) reference position of this atom is copied to ref_x.
1547 * Beware, the xc coords have been sorted in do_flexible */
1548 copy_rvec(erg->xc_ref_sorted[l], ref_x);
1550 /* Save data for doing angular RMSD fit later */
1551 /* Save the current atom position */
1552 copy_rvec(curr_x, sd->x[ind]);
1553 /* Save the corresponding reference position */
1554 copy_rvec(ref_x, sd->ref[ind]);
1556 /* Maybe also mass-weighting was requested. If yes, additionally
1557 * multiply the weights with the relative mass of the atom. If not,
1558 * multiply with unity. */
1559 m_rel = erg->mc_sorted[l]*OOm_av;
1561 /* Save the weight for this atom in this slab */
1562 sd->weight[ind] = gaussian_weight(curr_x, rotg, n) * m_rel;
1564 /* Next atom in this slab */
1569 /******************************/
1570 /* Now do the fit calculation */
1571 /******************************/
1573 fprintf(fp, "%12.3e%6d%12.3f", t, g, degangle);
1575 /* === Now do RMSD fitting for each slab === */
1576 /* We require at least SLAB_MIN_ATOMS in a slab, such that the fit makes sense. */
1577 #define SLAB_MIN_ATOMS 4
1579 for (n = erg->slab_first; n <= erg->slab_last; n++)
1581 islab = n - erg->slab_first; /* slab index */
1582 sd = &(rotg->enfrotgrp->slab_data[islab]);
1583 if (sd->nat >= SLAB_MIN_ATOMS)
1585 /* Get the center of the slabs reference and current positions */
1586 get_center(sd->ref, sd->weight, sd->nat, ref_center);
1587 get_center(sd->x, sd->weight, sd->nat, act_center);
1588 if (rotg->eFittype == erotgFitNORM)
1590 /* Normalize every position to it's reference length
1591 * prior to performing the fit */
1592 for (i = 0; i < sd->nat; i++) /* Center */
1594 rvec_dec(sd->ref[i], ref_center);
1595 rvec_dec(sd->x[i], act_center);
1596 /* Normalize x_i such that it gets the same length as ref_i */
1597 svmul( norm(sd->ref[i])/norm(sd->x[i]), sd->x[i], sd->x[i] );
1599 /* We already subtracted the centers */
1600 clear_rvec(ref_center);
1601 clear_rvec(act_center);
1603 fitangle = -opt_angle_analytic(sd->ref, sd->x, sd->weight, sd->nat,
1604 ref_center, act_center, rotg->vec);
1605 fprintf(fp, "%6d%6d%12.3f", n, sd->nat, fitangle);
1610 #undef SLAB_MIN_ATOMS
1614 /* Shift x with is */
1615 static gmx_inline void shift_single_coord(matrix box, rvec x, const ivec is)
1626 x[XX] += tx*box[XX][XX]+ty*box[YY][XX]+tz*box[ZZ][XX];
1627 x[YY] += ty*box[YY][YY]+tz*box[ZZ][YY];
1628 x[ZZ] += tz*box[ZZ][ZZ];
1632 x[XX] += tx*box[XX][XX];
1633 x[YY] += ty*box[YY][YY];
1634 x[ZZ] += tz*box[ZZ][ZZ];
1639 /* Determine the 'home' slab of this atom which is the
1640 * slab with the highest Gaussian weight of all */
1641 #define round(a) (int)(a+0.5)
1642 static gmx_inline int get_homeslab(
1643 rvec curr_x, /* The position for which the home slab shall be determined */
1644 rvec rotvec, /* The rotation vector */
1645 real slabdist) /* The slab distance */
1650 /* The distance of the atom to the coordinate center (where the
1651 * slab with index 0) is */
1652 dist = iprod(rotvec, curr_x);
1654 return round(dist / slabdist);
1658 /* For a local atom determine the relevant slabs, i.e. slabs in
1659 * which the gaussian is larger than min_gaussian
1661 static int get_single_atom_gaussians(
1668 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
1671 erg = rotg->enfrotgrp;
1673 /* Determine the 'home' slab of this atom: */
1674 homeslab = get_homeslab(curr_x, rotg->vec, rotg->slab_dist);
1676 /* First determine the weight in the atoms home slab: */
1677 g = gaussian_weight(curr_x, rotg, homeslab);
1679 erg->gn_atom[count] = g;
1680 erg->gn_slabind[count] = homeslab;
1684 /* Determine the max slab */
1686 while (g > rotg->min_gaussian)
1689 g = gaussian_weight(curr_x, rotg, slab);
1690 erg->gn_slabind[count] = slab;
1691 erg->gn_atom[count] = g;
1696 /* Determine the min slab */
1701 g = gaussian_weight(curr_x, rotg, slab);
1702 erg->gn_slabind[count] = slab;
1703 erg->gn_atom[count] = g;
1706 while (g > rotg->min_gaussian);
1713 static void flex2_precalc_inner_sum(t_rotgrp *rotg)
1716 rvec xi; /* positions in the i-sum */
1717 rvec xcn, ycn; /* the current and the reference slab centers */
1720 rvec rin; /* Helper variables */
1723 real OOpsii, OOpsiistar;
1724 real sin_rin; /* s_ii.r_ii */
1725 rvec s_in, tmpvec, tmpvec2;
1726 real mi, wi; /* Mass-weighting of the positions */
1728 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
1731 erg = rotg->enfrotgrp;
1732 N_M = rotg->nat * erg->invmass;
1734 /* Loop over all slabs that contain something */
1735 for (n = erg->slab_first; n <= erg->slab_last; n++)
1737 islab = n - erg->slab_first; /* slab index */
1739 /* The current center of this slab is saved in xcn: */
1740 copy_rvec(erg->slab_center[islab], xcn);
1741 /* ... and the reference center in ycn: */
1742 copy_rvec(erg->slab_center_ref[islab+erg->slab_buffer], ycn);
1744 /*** D. Calculate the whole inner sum used for second and third sum */
1745 /* For slab n, we need to loop over all atoms i again. Since we sorted
1746 * the atoms with respect to the rotation vector, we know that it is sufficient
1747 * to calculate from firstatom to lastatom only. All other contributions will
1749 clear_rvec(innersumvec);
1750 for (i = erg->firstatom[islab]; i <= erg->lastatom[islab]; i++)
1752 /* Coordinate xi of this atom */
1753 copy_rvec(erg->xc[i], xi);
1756 gaussian_xi = gaussian_weight(xi, rotg, n);
1757 mi = erg->mc_sorted[i]; /* need the sorted mass here */
1761 copy_rvec(erg->xc_ref_sorted[i], yi0); /* Reference position yi0 */
1762 rvec_sub(yi0, ycn, tmpvec2); /* tmpvec2 = yi0 - ycn */
1763 mvmul(erg->rotmat, tmpvec2, rin); /* rin = Omega.(yi0 - ycn) */
1765 /* Calculate psi_i* and sin */
1766 rvec_sub(xi, xcn, tmpvec2); /* tmpvec2 = xi - xcn */
1767 cprod(rotg->vec, tmpvec2, tmpvec); /* tmpvec = v x (xi - xcn) */
1768 OOpsiistar = norm2(tmpvec)+rotg->eps; /* OOpsii* = 1/psii* = |v x (xi-xcn)|^2 + eps */
1769 OOpsii = norm(tmpvec); /* OOpsii = 1 / psii = |v x (xi - xcn)| */
1771 /* * v x (xi - xcn) */
1772 unitv(tmpvec, s_in); /* sin = ---------------- */
1773 /* |v x (xi - xcn)| */
1775 sin_rin = iprod(s_in, rin); /* sin_rin = sin . rin */
1777 /* Now the whole sum */
1778 fac = OOpsii/OOpsiistar;
1779 svmul(fac, rin, tmpvec);
1780 fac2 = fac*fac*OOpsii;
1781 svmul(fac2*sin_rin, s_in, tmpvec2);
1782 rvec_dec(tmpvec, tmpvec2);
1784 svmul(wi*gaussian_xi*sin_rin, tmpvec, tmpvec2);
1786 rvec_inc(innersumvec, tmpvec2);
1787 } /* now we have the inner sum, used both for sum2 and sum3 */
1789 /* Save it to be used in do_flex2_lowlevel */
1790 copy_rvec(innersumvec, erg->slab_innersumvec[islab]);
1791 } /* END of loop over slabs */
1795 static void flex_precalc_inner_sum(t_rotgrp *rotg)
1798 rvec xi; /* position */
1799 rvec xcn, ycn; /* the current and the reference slab centers */
1800 rvec qin, rin; /* q_i^n and r_i^n */
1803 rvec innersumvec; /* Inner part of sum_n2 */
1804 real gaussian_xi; /* Gaussian weight gn(xi) */
1805 real mi, wi; /* Mass-weighting of the positions */
1808 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
1811 erg = rotg->enfrotgrp;
1812 N_M = rotg->nat * erg->invmass;
1814 /* Loop over all slabs that contain something */
1815 for (n = erg->slab_first; n <= erg->slab_last; n++)
1817 islab = n - erg->slab_first; /* slab index */
1819 /* The current center of this slab is saved in xcn: */
1820 copy_rvec(erg->slab_center[islab], xcn);
1821 /* ... and the reference center in ycn: */
1822 copy_rvec(erg->slab_center_ref[islab+erg->slab_buffer], ycn);
1824 /* For slab n, we need to loop over all atoms i again. Since we sorted
1825 * the atoms with respect to the rotation vector, we know that it is sufficient
1826 * to calculate from firstatom to lastatom only. All other contributions will
1828 clear_rvec(innersumvec);
1829 for (i = erg->firstatom[islab]; i <= erg->lastatom[islab]; i++)
1831 /* Coordinate xi of this atom */
1832 copy_rvec(erg->xc[i], xi);
1835 gaussian_xi = gaussian_weight(xi, rotg, n);
1836 mi = erg->mc_sorted[i]; /* need the sorted mass here */
1839 /* Calculate rin and qin */
1840 rvec_sub(erg->xc_ref_sorted[i], ycn, tmpvec); /* tmpvec = yi0-ycn */
1841 mvmul(erg->rotmat, tmpvec, rin); /* rin = Omega.(yi0 - ycn) */
1842 cprod(rotg->vec, rin, tmpvec); /* tmpvec = v x Omega*(yi0-ycn) */
1844 /* * v x Omega*(yi0-ycn) */
1845 unitv(tmpvec, qin); /* qin = --------------------- */
1846 /* |v x Omega*(yi0-ycn)| */
1849 rvec_sub(xi, xcn, tmpvec); /* tmpvec = xi-xcn */
1850 bin = iprod(qin, tmpvec); /* bin = qin*(xi-xcn) */
1852 svmul(wi*gaussian_xi*bin, qin, tmpvec);
1854 /* Add this contribution to the inner sum: */
1855 rvec_add(innersumvec, tmpvec, innersumvec);
1856 } /* now we have the inner sum vector S^n for this slab */
1857 /* Save it to be used in do_flex_lowlevel */
1858 copy_rvec(innersumvec, erg->slab_innersumvec[islab]);
1863 static real do_flex2_lowlevel(
1865 real sigma, /* The Gaussian width sigma */
1867 gmx_bool bOutstepRot,
1868 gmx_bool bOutstepSlab,
1871 int count, ic, ii, j, m, n, islab, iigrp, ifit;
1872 rvec xj; /* position in the i-sum */
1873 rvec yj0; /* the reference position in the j-sum */
1874 rvec xcn, ycn; /* the current and the reference slab centers */
1875 real V; /* This node's part of the rotation pot. energy */
1876 real gaussian_xj; /* Gaussian weight */
1879 real numerator, fit_numerator;
1880 rvec rjn, fit_rjn; /* Helper variables */
1883 real OOpsij, OOpsijstar;
1884 real OOsigma2; /* 1/(sigma^2) */
1887 rvec sjn, tmpvec, tmpvec2, yj0_ycn;
1888 rvec sum1vec_part, sum1vec, sum2vec_part, sum2vec, sum3vec, sum4vec, innersumvec;
1890 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
1891 real mj, wj; /* Mass-weighting of the positions */
1893 real Wjn; /* g_n(x_j) m_j / Mjn */
1894 gmx_bool bCalcPotFit;
1896 /* To calculate the torque per slab */
1897 rvec slab_force; /* Single force from slab n on one atom */
1898 rvec slab_sum1vec_part;
1899 real slab_sum3part, slab_sum4part;
1900 rvec slab_sum1vec, slab_sum2vec, slab_sum3vec, slab_sum4vec;
1903 erg = rotg->enfrotgrp;
1905 /* Pre-calculate the inner sums, so that we do not have to calculate
1906 * them again for every atom */
1907 flex2_precalc_inner_sum(rotg);
1909 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT == rotg->eFittype);
1911 /********************************************************/
1912 /* Main loop over all local atoms of the rotation group */
1913 /********************************************************/
1914 N_M = rotg->nat * erg->invmass;
1916 OOsigma2 = 1.0 / (sigma*sigma);
1917 for (j = 0; j < erg->nat_loc; j++)
1919 /* Local index of a rotation group atom */
1920 ii = erg->ind_loc[j];
1921 /* Position of this atom in the collective array */
1922 iigrp = erg->xc_ref_ind[j];
1923 /* Mass-weighting */
1924 mj = erg->mc[iigrp]; /* need the unsorted mass here */
1927 /* Current position of this atom: x[ii][XX/YY/ZZ]
1928 * Note that erg->xc_center contains the center of mass in case the flex2-t
1929 * potential was chosen. For the flex2 potential erg->xc_center must be
1931 rvec_sub(x[ii], erg->xc_center, xj);
1933 /* Shift this atom such that it is near its reference */
1934 shift_single_coord(box, xj, erg->xc_shifts[iigrp]);
1936 /* Determine the slabs to loop over, i.e. the ones with contributions
1937 * larger than min_gaussian */
1938 count = get_single_atom_gaussians(xj, rotg);
1940 clear_rvec(sum1vec_part);
1941 clear_rvec(sum2vec_part);
1944 /* Loop over the relevant slabs for this atom */
1945 for (ic = 0; ic < count; ic++)
1947 n = erg->gn_slabind[ic];
1949 /* Get the precomputed Gaussian value of curr_slab for curr_x */
1950 gaussian_xj = erg->gn_atom[ic];
1952 islab = n - erg->slab_first; /* slab index */
1954 /* The (unrotated) reference position of this atom is copied to yj0: */
1955 copy_rvec(rotg->x_ref[iigrp], yj0);
1957 beta = calc_beta(xj, rotg, n);
1959 /* The current center of this slab is saved in xcn: */
1960 copy_rvec(erg->slab_center[islab], xcn);
1961 /* ... and the reference center in ycn: */
1962 copy_rvec(erg->slab_center_ref[islab+erg->slab_buffer], ycn);
1964 rvec_sub(yj0, ycn, yj0_ycn); /* yj0_ycn = yj0 - ycn */
1967 mvmul(erg->rotmat, yj0_ycn, rjn); /* rjn = Omega.(yj0 - ycn) */
1969 /* Subtract the slab center from xj */
1970 rvec_sub(xj, xcn, tmpvec2); /* tmpvec2 = xj - xcn */
1972 /* In rare cases, when an atom position coincides with a slab center
1973 * (tmpvec2 == 0) we cannot compute the vector product for sjn.
1974 * However, since the atom is located directly on the pivot, this
1975 * slab's contribution to the force on that atom will be zero
1976 * anyway. Therefore, we directly move on to the next slab. */
1977 if (gmx_numzero(norm(tmpvec2))) /* 0 == norm(xj - xcn) */
1983 cprod(rotg->vec, tmpvec2, tmpvec); /* tmpvec = v x (xj - xcn) */
1985 OOpsijstar = norm2(tmpvec)+rotg->eps; /* OOpsij* = 1/psij* = |v x (xj-xcn)|^2 + eps */
1987 numerator = sqr(iprod(tmpvec, rjn));
1989 /*********************************/
1990 /* Add to the rotation potential */
1991 /*********************************/
1992 V += 0.5*rotg->k*wj*gaussian_xj*numerator/OOpsijstar;
1994 /* If requested, also calculate the potential for a set of angles
1995 * near the current reference angle */
1998 for (ifit = 0; ifit < rotg->PotAngle_nstep; ifit++)
2000 mvmul(erg->PotAngleFit->rotmat[ifit], yj0_ycn, fit_rjn);
2001 fit_numerator = sqr(iprod(tmpvec, fit_rjn));
2002 erg->PotAngleFit->V[ifit] += 0.5*rotg->k*wj*gaussian_xj*fit_numerator/OOpsijstar;
2006 /*************************************/
2007 /* Now calculate the force on atom j */
2008 /*************************************/
2010 OOpsij = norm(tmpvec); /* OOpsij = 1 / psij = |v x (xj - xcn)| */
2012 /* * v x (xj - xcn) */
2013 unitv(tmpvec, sjn); /* sjn = ---------------- */
2014 /* |v x (xj - xcn)| */
2016 sjn_rjn = iprod(sjn, rjn); /* sjn_rjn = sjn . rjn */
2019 /*** A. Calculate the first of the four sum terms: ****************/
2020 fac = OOpsij/OOpsijstar;
2021 svmul(fac, rjn, tmpvec);
2022 fac2 = fac*fac*OOpsij;
2023 svmul(fac2*sjn_rjn, sjn, tmpvec2);
2024 rvec_dec(tmpvec, tmpvec2);
2025 fac2 = wj*gaussian_xj; /* also needed for sum4 */
2026 svmul(fac2*sjn_rjn, tmpvec, slab_sum1vec_part);
2027 /********************/
2028 /*** Add to sum1: ***/
2029 /********************/
2030 rvec_inc(sum1vec_part, slab_sum1vec_part); /* sum1 still needs to vector multiplied with v */
2032 /*** B. Calculate the forth of the four sum terms: ****************/
2033 betasigpsi = beta*OOsigma2*OOpsij; /* this is also needed for sum3 */
2034 /********************/
2035 /*** Add to sum4: ***/
2036 /********************/
2037 slab_sum4part = fac2*betasigpsi*fac*sjn_rjn*sjn_rjn; /* Note that fac is still valid from above */
2038 sum4 += slab_sum4part;
2040 /*** C. Calculate Wjn for second and third sum */
2041 /* Note that we can safely divide by slab_weights since we check in
2042 * get_slab_centers that it is non-zero. */
2043 Wjn = gaussian_xj*mj/erg->slab_weights[islab];
2045 /* We already have precalculated the inner sum for slab n */
2046 copy_rvec(erg->slab_innersumvec[islab], innersumvec);
2048 /* Weigh the inner sum vector with Wjn */
2049 svmul(Wjn, innersumvec, innersumvec);
2051 /*** E. Calculate the second of the four sum terms: */
2052 /********************/
2053 /*** Add to sum2: ***/
2054 /********************/
2055 rvec_inc(sum2vec_part, innersumvec); /* sum2 still needs to be vector crossproduct'ed with v */
2057 /*** F. Calculate the third of the four sum terms: */
2058 slab_sum3part = betasigpsi * iprod(sjn, innersumvec);
2059 sum3 += slab_sum3part; /* still needs to be multiplied with v */
2061 /*** G. Calculate the torque on the local slab's axis: */
2065 cprod(slab_sum1vec_part, rotg->vec, slab_sum1vec);
2067 cprod(innersumvec, rotg->vec, slab_sum2vec);
2069 svmul(slab_sum3part, rotg->vec, slab_sum3vec);
2071 svmul(slab_sum4part, rotg->vec, slab_sum4vec);
2073 /* The force on atom ii from slab n only: */
2074 for (m = 0; m < DIM; m++)
2076 slab_force[m] = rotg->k * (-slab_sum1vec[m] + slab_sum2vec[m] - slab_sum3vec[m] + 0.5*slab_sum4vec[m]);
2079 erg->slab_torque_v[islab] += torque(rotg->vec, slab_force, xj, xcn);
2081 } /* END of loop over slabs */
2083 /* Construct the four individual parts of the vector sum: */
2084 cprod(sum1vec_part, rotg->vec, sum1vec); /* sum1vec = { } x v */
2085 cprod(sum2vec_part, rotg->vec, sum2vec); /* sum2vec = { } x v */
2086 svmul(sum3, rotg->vec, sum3vec); /* sum3vec = { } . v */
2087 svmul(sum4, rotg->vec, sum4vec); /* sum4vec = { } . v */
2089 /* Store the additional force so that it can be added to the force
2090 * array after the normal forces have been evaluated */
2091 for (m = 0; m < DIM; m++)
2093 erg->f_rot_loc[j][m] = rotg->k * (-sum1vec[m] + sum2vec[m] - sum3vec[m] + 0.5*sum4vec[m]);
2097 fprintf(stderr, "sum1: %15.8f %15.8f %15.8f\n", -rotg->k*sum1vec[XX], -rotg->k*sum1vec[YY], -rotg->k*sum1vec[ZZ]);
2098 fprintf(stderr, "sum2: %15.8f %15.8f %15.8f\n", rotg->k*sum2vec[XX], rotg->k*sum2vec[YY], rotg->k*sum2vec[ZZ]);
2099 fprintf(stderr, "sum3: %15.8f %15.8f %15.8f\n", -rotg->k*sum3vec[XX], -rotg->k*sum3vec[YY], -rotg->k*sum3vec[ZZ]);
2100 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]);
2105 } /* END of loop over local atoms */
2111 static real do_flex_lowlevel(
2113 real sigma, /* The Gaussian width sigma */
2115 gmx_bool bOutstepRot,
2116 gmx_bool bOutstepSlab,
2119 int count, ic, ifit, ii, j, m, n, islab, iigrp;
2120 rvec xj, yj0; /* current and reference position */
2121 rvec xcn, ycn; /* the current and the reference slab centers */
2122 rvec yj0_ycn; /* yj0 - ycn */
2123 rvec xj_xcn; /* xj - xcn */
2124 rvec qjn, fit_qjn; /* q_i^n */
2125 rvec sum_n1, sum_n2; /* Two contributions to the rotation force */
2126 rvec innersumvec; /* Inner part of sum_n2 */
2128 rvec force_n; /* Single force from slab n on one atom */
2129 rvec force_n1, force_n2; /* First and second part of force_n */
2130 rvec tmpvec, tmpvec2, tmp_f; /* Helper variables */
2131 real V; /* The rotation potential energy */
2132 real OOsigma2; /* 1/(sigma^2) */
2133 real beta; /* beta_n(xj) */
2134 real bjn, fit_bjn; /* b_j^n */
2135 real gaussian_xj; /* Gaussian weight gn(xj) */
2136 real betan_xj_sigma2;
2137 real mj, wj; /* Mass-weighting of the positions */
2139 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
2140 gmx_bool bCalcPotFit;
2143 erg = rotg->enfrotgrp;
2145 /* Pre-calculate the inner sums, so that we do not have to calculate
2146 * them again for every atom */
2147 flex_precalc_inner_sum(rotg);
2149 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT == rotg->eFittype);
2151 /********************************************************/
2152 /* Main loop over all local atoms of the rotation group */
2153 /********************************************************/
2154 OOsigma2 = 1.0/(sigma*sigma);
2155 N_M = rotg->nat * erg->invmass;
2157 for (j = 0; j < erg->nat_loc; j++)
2159 /* Local index of a rotation group atom */
2160 ii = erg->ind_loc[j];
2161 /* Position of this atom in the collective array */
2162 iigrp = erg->xc_ref_ind[j];
2163 /* Mass-weighting */
2164 mj = erg->mc[iigrp]; /* need the unsorted mass here */
2167 /* Current position of this atom: x[ii][XX/YY/ZZ]
2168 * Note that erg->xc_center contains the center of mass in case the flex-t
2169 * potential was chosen. For the flex potential erg->xc_center must be
2171 rvec_sub(x[ii], erg->xc_center, xj);
2173 /* Shift this atom such that it is near its reference */
2174 shift_single_coord(box, xj, erg->xc_shifts[iigrp]);
2176 /* Determine the slabs to loop over, i.e. the ones with contributions
2177 * larger than min_gaussian */
2178 count = get_single_atom_gaussians(xj, rotg);
2183 /* Loop over the relevant slabs for this atom */
2184 for (ic = 0; ic < count; ic++)
2186 n = erg->gn_slabind[ic];
2188 /* Get the precomputed Gaussian for xj in slab n */
2189 gaussian_xj = erg->gn_atom[ic];
2191 islab = n - erg->slab_first; /* slab index */
2193 /* The (unrotated) reference position of this atom is saved in yj0: */
2194 copy_rvec(rotg->x_ref[iigrp], yj0);
2196 beta = calc_beta(xj, rotg, n);
2198 /* The current center of this slab is saved in xcn: */
2199 copy_rvec(erg->slab_center[islab], xcn);
2200 /* ... and the reference center in ycn: */
2201 copy_rvec(erg->slab_center_ref[islab+erg->slab_buffer], ycn);
2203 rvec_sub(yj0, ycn, yj0_ycn); /* yj0_ycn = yj0 - ycn */
2205 /* In rare cases, when an atom position coincides with a reference slab
2206 * center (yj0_ycn == 0) we cannot compute the normal vector qjn.
2207 * However, since the atom is located directly on the pivot, this
2208 * slab's contribution to the force on that atom will be zero
2209 * anyway. Therefore, we directly move on to the next slab. */
2210 if (gmx_numzero(norm(yj0_ycn))) /* 0 == norm(yj0 - ycn) */
2216 mvmul(erg->rotmat, yj0_ycn, tmpvec2); /* tmpvec2= Omega.(yj0-ycn) */
2218 /* Subtract the slab center from xj */
2219 rvec_sub(xj, xcn, xj_xcn); /* xj_xcn = xj - xcn */
2222 cprod(rotg->vec, tmpvec2, tmpvec); /* tmpvec= v x Omega.(yj0-ycn) */
2224 /* * v x Omega.(yj0-ycn) */
2225 unitv(tmpvec, qjn); /* qjn = --------------------- */
2226 /* |v x Omega.(yj0-ycn)| */
2228 bjn = iprod(qjn, xj_xcn); /* bjn = qjn * (xj - xcn) */
2230 /*********************************/
2231 /* Add to the rotation potential */
2232 /*********************************/
2233 V += 0.5*rotg->k*wj*gaussian_xj*sqr(bjn);
2235 /* If requested, also calculate the potential for a set of angles
2236 * near the current reference angle */
2239 for (ifit = 0; ifit < rotg->PotAngle_nstep; ifit++)
2241 /* As above calculate Omega.(yj0-ycn), now for the other angles */
2242 mvmul(erg->PotAngleFit->rotmat[ifit], yj0_ycn, tmpvec2); /* tmpvec2= Omega.(yj0-ycn) */
2243 /* As above calculate qjn */
2244 cprod(rotg->vec, tmpvec2, tmpvec); /* tmpvec= v x Omega.(yj0-ycn) */
2245 /* * v x Omega.(yj0-ycn) */
2246 unitv(tmpvec, fit_qjn); /* fit_qjn = --------------------- */
2247 /* |v x Omega.(yj0-ycn)| */
2248 fit_bjn = iprod(fit_qjn, xj_xcn); /* fit_bjn = fit_qjn * (xj - xcn) */
2249 /* Add to the rotation potential for this angle */
2250 erg->PotAngleFit->V[ifit] += 0.5*rotg->k*wj*gaussian_xj*sqr(fit_bjn);
2254 /****************************************************************/
2255 /* sum_n1 will typically be the main contribution to the force: */
2256 /****************************************************************/
2257 betan_xj_sigma2 = beta*OOsigma2; /* beta_n(xj)/sigma^2 */
2259 /* The next lines calculate
2260 * qjn - (bjn*beta(xj)/(2sigma^2))v */
2261 svmul(bjn*0.5*betan_xj_sigma2, rotg->vec, tmpvec2);
2262 rvec_sub(qjn, tmpvec2, tmpvec);
2264 /* Multiply with gn(xj)*bjn: */
2265 svmul(gaussian_xj*bjn, tmpvec, tmpvec2);
2268 rvec_inc(sum_n1, tmpvec2);
2270 /* We already have precalculated the Sn term for slab n */
2271 copy_rvec(erg->slab_innersumvec[islab], s_n);
2273 svmul(betan_xj_sigma2*iprod(s_n, xj_xcn), rotg->vec, tmpvec); /* tmpvec = ---------- s_n (xj-xcn) */
2276 rvec_sub(s_n, tmpvec, innersumvec);
2278 /* We can safely divide by slab_weights since we check in get_slab_centers
2279 * that it is non-zero. */
2280 svmul(gaussian_xj/erg->slab_weights[islab], innersumvec, innersumvec);
2282 rvec_add(sum_n2, innersumvec, sum_n2);
2284 /* Calculate the torque: */
2287 /* The force on atom ii from slab n only: */
2288 svmul(-rotg->k*wj, tmpvec2, force_n1); /* part 1 */
2289 svmul( rotg->k*mj, innersumvec, force_n2); /* part 2 */
2290 rvec_add(force_n1, force_n2, force_n);
2291 erg->slab_torque_v[islab] += torque(rotg->vec, force_n, xj, xcn);
2293 } /* END of loop over slabs */
2295 /* Put both contributions together: */
2296 svmul(wj, sum_n1, sum_n1);
2297 svmul(mj, sum_n2, sum_n2);
2298 rvec_sub(sum_n2, sum_n1, tmp_f); /* F = -grad V */
2300 /* Store the additional force so that it can be added to the force
2301 * array after the normal forces have been evaluated */
2302 for (m = 0; m < DIM; m++)
2304 erg->f_rot_loc[j][m] = rotg->k*tmp_f[m];
2309 } /* END of loop over local atoms */
2315 static void print_coordinates(t_rotgrp *rotg, rvec x[], matrix box, int step)
2319 static char buf[STRLEN];
2320 static gmx_bool bFirst = 1;
2325 sprintf(buf, "coords%d.txt", cr->nodeid);
2326 fp = fopen(buf, "w");
2330 fprintf(fp, "\nStep %d\n", step);
2331 fprintf(fp, "box: %f %f %f %f %f %f %f %f %f\n",
2332 box[XX][XX], box[XX][YY], box[XX][ZZ],
2333 box[YY][XX], box[YY][YY], box[YY][ZZ],
2334 box[ZZ][XX], box[ZZ][ZZ], box[ZZ][ZZ]);
2335 for (i = 0; i < rotg->nat; i++)
2337 fprintf(fp, "%4d %f %f %f\n", i,
2338 erg->xc[i][XX], erg->xc[i][YY], erg->xc[i][ZZ]);
2346 static int projection_compare(const void *a, const void *b)
2348 sort_along_vec_t *xca, *xcb;
2351 xca = (sort_along_vec_t *)a;
2352 xcb = (sort_along_vec_t *)b;
2354 if (xca->xcproj < xcb->xcproj)
2358 else if (xca->xcproj > xcb->xcproj)
2369 static void sort_collective_coordinates(
2370 t_rotgrp *rotg, /* Rotation group */
2371 sort_along_vec_t *data) /* Buffer for sorting the positions */
2374 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
2377 erg = rotg->enfrotgrp;
2379 /* The projection of the position vector on the rotation vector is
2380 * the relevant value for sorting. Fill the 'data' structure */
2381 for (i = 0; i < rotg->nat; i++)
2383 data[i].xcproj = iprod(erg->xc[i], rotg->vec); /* sort criterium */
2384 data[i].m = erg->mc[i];
2386 copy_rvec(erg->xc[i], data[i].x );
2387 copy_rvec(rotg->x_ref[i], data[i].x_ref);
2389 /* Sort the 'data' structure */
2390 gmx_qsort(data, rotg->nat, sizeof(sort_along_vec_t), projection_compare);
2392 /* Copy back the sorted values */
2393 for (i = 0; i < rotg->nat; i++)
2395 copy_rvec(data[i].x, erg->xc[i] );
2396 copy_rvec(data[i].x_ref, erg->xc_ref_sorted[i]);
2397 erg->mc_sorted[i] = data[i].m;
2398 erg->xc_sortind[i] = data[i].ind;
2403 /* For each slab, get the first and the last index of the sorted atom
2405 static void get_firstlast_atom_per_slab(t_rotgrp *rotg)
2409 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
2412 erg = rotg->enfrotgrp;
2414 /* Find the first atom that needs to enter the calculation for each slab */
2415 n = erg->slab_first; /* slab */
2416 i = 0; /* start with the first atom */
2419 /* Find the first atom that significantly contributes to this slab */
2420 do /* move forward in position until a large enough beta is found */
2422 beta = calc_beta(erg->xc[i], rotg, n);
2425 while ((beta < -erg->max_beta) && (i < rotg->nat));
2427 islab = n - erg->slab_first; /* slab index */
2428 erg->firstatom[islab] = i;
2429 /* Proceed to the next slab */
2432 while (n <= erg->slab_last);
2434 /* Find the last atom for each slab */
2435 n = erg->slab_last; /* start with last slab */
2436 i = rotg->nat-1; /* start with the last atom */
2439 do /* move backward in position until a large enough beta is found */
2441 beta = calc_beta(erg->xc[i], rotg, n);
2444 while ((beta > erg->max_beta) && (i > -1));
2446 islab = n - erg->slab_first; /* slab index */
2447 erg->lastatom[islab] = i;
2448 /* Proceed to the next slab */
2451 while (n >= erg->slab_first);
2455 /* Determine the very first and very last slab that needs to be considered
2456 * For the first slab that needs to be considered, we have to find the smallest
2459 * x_first * v - n*Delta_x <= beta_max
2461 * slab index n, slab distance Delta_x, rotation vector v. For the last slab we
2462 * have to find the largest n that obeys
2464 * x_last * v - n*Delta_x >= -beta_max
2467 static gmx_inline int get_first_slab(
2468 t_rotgrp *rotg, /* The rotation group (inputrec data) */
2469 real max_beta, /* The max_beta value, instead of min_gaussian */
2470 rvec firstatom) /* First atom after sorting along the rotation vector v */
2472 /* Find the first slab for the first atom */
2473 return ceil((iprod(firstatom, rotg->vec) - max_beta)/rotg->slab_dist);
2477 static gmx_inline int get_last_slab(
2478 t_rotgrp *rotg, /* The rotation group (inputrec data) */
2479 real max_beta, /* The max_beta value, instead of min_gaussian */
2480 rvec lastatom) /* Last atom along v */
2482 /* Find the last slab for the last atom */
2483 return floor((iprod(lastatom, rotg->vec) + max_beta)/rotg->slab_dist);
2487 static void get_firstlast_slab_check(
2488 t_rotgrp *rotg, /* The rotation group (inputrec data) */
2489 t_gmx_enfrotgrp *erg, /* The rotation group (data only accessible in this file) */
2490 rvec firstatom, /* First atom after sorting along the rotation vector v */
2491 rvec lastatom) /* Last atom along v */
2493 erg->slab_first = get_first_slab(rotg, erg->max_beta, firstatom);
2494 erg->slab_last = get_last_slab(rotg, erg->max_beta, lastatom);
2496 /* Calculate the slab buffer size, which changes when slab_first changes */
2497 erg->slab_buffer = erg->slab_first - erg->slab_first_ref;
2499 /* Check whether we have reference data to compare against */
2500 if (erg->slab_first < erg->slab_first_ref)
2502 gmx_fatal(FARGS, "%s No reference data for first slab (n=%d), unable to proceed.",
2503 RotStr, erg->slab_first);
2506 /* Check whether we have reference data to compare against */
2507 if (erg->slab_last > erg->slab_last_ref)
2509 gmx_fatal(FARGS, "%s No reference data for last slab (n=%d), unable to proceed.",
2510 RotStr, erg->slab_last);
2515 /* Enforced rotation with a flexible axis */
2516 static void do_flexible(
2518 gmx_enfrot_t enfrot, /* Other rotation data */
2519 t_rotgrp *rotg, /* The rotation group */
2520 int g, /* Group number */
2521 rvec x[], /* The local positions */
2523 double t, /* Time in picoseconds */
2524 gmx_bool bOutstepRot, /* Output to main rotation output file */
2525 gmx_bool bOutstepSlab) /* Output per-slab data */
2528 real sigma; /* The Gaussian width sigma */
2529 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
2532 erg = rotg->enfrotgrp;
2534 /* Define the sigma value */
2535 sigma = 0.7*rotg->slab_dist;
2537 /* Sort the collective coordinates erg->xc along the rotation vector. This is
2538 * an optimization for the inner loop. */
2539 sort_collective_coordinates(rotg, enfrot->data);
2541 /* Determine the first relevant slab for the first atom and the last
2542 * relevant slab for the last atom */
2543 get_firstlast_slab_check(rotg, erg, erg->xc[0], erg->xc[rotg->nat-1]);
2545 /* Determine for each slab depending on the min_gaussian cutoff criterium,
2546 * a first and a last atom index inbetween stuff needs to be calculated */
2547 get_firstlast_atom_per_slab(rotg);
2549 /* Determine the gaussian-weighted center of positions for all slabs */
2550 get_slab_centers(rotg, erg->xc, erg->mc_sorted, g, t, enfrot->out_slabs, bOutstepSlab, FALSE);
2552 /* Clear the torque per slab from last time step: */
2553 nslabs = erg->slab_last - erg->slab_first + 1;
2554 for (l = 0; l < nslabs; l++)
2556 erg->slab_torque_v[l] = 0.0;
2559 /* Call the rotational forces kernel */
2560 if (rotg->eType == erotgFLEX || rotg->eType == erotgFLEXT)
2562 erg->V = do_flex_lowlevel(rotg, sigma, x, bOutstepRot, bOutstepSlab, box);
2564 else if (rotg->eType == erotgFLEX2 || rotg->eType == erotgFLEX2T)
2566 erg->V = do_flex2_lowlevel(rotg, sigma, x, bOutstepRot, bOutstepSlab, box);
2570 gmx_fatal(FARGS, "Unknown flexible rotation type");
2573 /* Determine angle by RMSD fit to the reference - Let's hope this */
2574 /* only happens once in a while, since this is not parallelized! */
2575 if (bMaster && (erotgFitPOT != rotg->eFittype) )
2579 /* Fit angle of the whole rotation group */
2580 erg->angle_v = flex_fit_angle(rotg);
2584 /* Fit angle of each slab */
2585 flex_fit_angle_perslab(g, rotg, t, erg->degangle, enfrot->out_angles);
2589 /* Lump together the torques from all slabs: */
2590 erg->torque_v = 0.0;
2591 for (l = 0; l < nslabs; l++)
2593 erg->torque_v += erg->slab_torque_v[l];
2598 /* Calculate the angle between reference and actual rotation group atom,
2599 * both projected into a plane perpendicular to the rotation vector: */
2600 static void angle(t_rotgrp *rotg,
2604 real *weight) /* atoms near the rotation axis should count less than atoms far away */
2606 rvec xp, xrp; /* current and reference positions projected on a plane perpendicular to pg->vec */
2610 /* Project x_ref and x into a plane through the origin perpendicular to rot_vec: */
2611 /* Project x_ref: xrp = x_ref - (vec * x_ref) * vec */
2612 svmul(iprod(rotg->vec, x_ref), rotg->vec, dum);
2613 rvec_sub(x_ref, dum, xrp);
2614 /* Project x_act: */
2615 svmul(iprod(rotg->vec, x_act), rotg->vec, dum);
2616 rvec_sub(x_act, dum, xp);
2618 /* Retrieve information about which vector precedes. gmx_angle always
2619 * returns a positive angle. */
2620 cprod(xp, xrp, dum); /* if reference precedes, this is pointing into the same direction as vec */
2622 if (iprod(rotg->vec, dum) >= 0)
2624 *alpha = -gmx_angle(xrp, xp);
2628 *alpha = +gmx_angle(xrp, xp);
2631 /* Also return the weight */
2636 /* Project first vector onto a plane perpendicular to the second vector
2638 * Note that v must be of unit length.
2640 static gmx_inline void project_onto_plane(rvec dr, const rvec v)
2645 svmul(iprod(dr, v), v, tmp); /* tmp = (dr.v)v */
2646 rvec_dec(dr, tmp); /* dr = dr - (dr.v)v */
2650 /* Fixed rotation: The rotation reference group rotates around the v axis. */
2651 /* The atoms of the actual rotation group are attached with imaginary */
2652 /* springs to the reference atoms. */
2653 static void do_fixed(
2654 t_rotgrp *rotg, /* The rotation group */
2655 gmx_bool bOutstepRot, /* Output to main rotation output file */
2656 gmx_bool bOutstepSlab) /* Output per-slab data */
2660 rvec tmp_f; /* Force */
2661 real alpha; /* a single angle between an actual and a reference position */
2662 real weight; /* single weight for a single angle */
2663 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
2664 rvec xi_xc; /* xi - xc */
2665 gmx_bool bCalcPotFit;
2668 /* for mass weighting: */
2669 real wi; /* Mass-weighting of the positions */
2671 real k_wi; /* k times wi */
2676 erg = rotg->enfrotgrp;
2677 bProject = (rotg->eType == erotgPM) || (rotg->eType == erotgPMPF);
2678 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT == rotg->eFittype);
2680 N_M = rotg->nat * erg->invmass;
2682 /* Each process calculates the forces on its local atoms */
2683 for (j = 0; j < erg->nat_loc; j++)
2685 /* Calculate (x_i-x_c) resp. (x_i-u) */
2686 rvec_sub(erg->x_loc_pbc[j], erg->xc_center, xi_xc);
2688 /* Calculate Omega*(y_i-y_c)-(x_i-x_c) */
2689 rvec_sub(erg->xr_loc[j], xi_xc, dr);
2693 project_onto_plane(dr, rotg->vec);
2696 /* Mass-weighting */
2697 wi = N_M*erg->m_loc[j];
2699 /* Store the additional force so that it can be added to the force
2700 * array after the normal forces have been evaluated */
2702 for (m = 0; m < DIM; m++)
2704 tmp_f[m] = k_wi*dr[m];
2705 erg->f_rot_loc[j][m] = tmp_f[m];
2706 erg->V += 0.5*k_wi*sqr(dr[m]);
2709 /* If requested, also calculate the potential for a set of angles
2710 * near the current reference angle */
2713 for (ifit = 0; ifit < rotg->PotAngle_nstep; ifit++)
2715 /* Index of this rotation group atom with respect to the whole rotation group */
2716 jj = erg->xc_ref_ind[j];
2718 /* Rotate with the alternative angle. Like rotate_local_reference(),
2719 * just for a single local atom */
2720 mvmul(erg->PotAngleFit->rotmat[ifit], rotg->x_ref[jj], fit_xr_loc); /* fit_xr_loc = Omega*(y_i-y_c) */
2722 /* Calculate Omega*(y_i-y_c)-(x_i-x_c) */
2723 rvec_sub(fit_xr_loc, xi_xc, dr);
2727 project_onto_plane(dr, rotg->vec);
2730 /* Add to the rotation potential for this angle: */
2731 erg->PotAngleFit->V[ifit] += 0.5*k_wi*norm2(dr);
2737 /* Add to the torque of this rotation group */
2738 erg->torque_v += torque(rotg->vec, tmp_f, erg->x_loc_pbc[j], erg->xc_center);
2740 /* Calculate the angle between reference and actual rotation group atom. */
2741 angle(rotg, xi_xc, erg->xr_loc[j], &alpha, &weight); /* angle in rad, weighted */
2742 erg->angle_v += alpha * weight;
2743 erg->weight_v += weight;
2745 /* If you want enforced rotation to contribute to the virial,
2746 * activate the following lines:
2749 Add the rotation contribution to the virial
2750 for(j=0; j<DIM; j++)
2752 vir[j][m] += 0.5*f[ii][j]*dr[m];
2758 } /* end of loop over local rotation group atoms */
2762 /* Calculate the radial motion potential and forces */
2763 static void do_radial_motion(
2764 t_rotgrp *rotg, /* The rotation group */
2765 gmx_bool bOutstepRot, /* Output to main rotation output file */
2766 gmx_bool bOutstepSlab) /* Output per-slab data */
2769 rvec tmp_f; /* Force */
2770 real alpha; /* a single angle between an actual and a reference position */
2771 real weight; /* single weight for a single angle */
2772 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
2773 rvec xj_u; /* xj - u */
2774 rvec tmpvec, fit_tmpvec;
2775 real fac, fac2, sum = 0.0;
2777 gmx_bool bCalcPotFit;
2779 /* For mass weighting: */
2780 real wj; /* Mass-weighting of the positions */
2784 erg = rotg->enfrotgrp;
2785 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT == rotg->eFittype);
2787 N_M = rotg->nat * erg->invmass;
2789 /* Each process calculates the forces on its local atoms */
2790 for (j = 0; j < erg->nat_loc; j++)
2792 /* Calculate (xj-u) */
2793 rvec_sub(erg->x_loc_pbc[j], erg->xc_center, xj_u); /* xj_u = xj-u */
2795 /* Calculate Omega.(yj0-u) */
2796 cprod(rotg->vec, erg->xr_loc[j], tmpvec); /* tmpvec = v x Omega.(yj0-u) */
2798 /* * v x Omega.(yj0-u) */
2799 unitv(tmpvec, pj); /* pj = --------------------- */
2800 /* | v x Omega.(yj0-u) | */
2802 fac = iprod(pj, xj_u); /* fac = pj.(xj-u) */
2805 /* Mass-weighting */
2806 wj = N_M*erg->m_loc[j];
2808 /* Store the additional force so that it can be added to the force
2809 * array after the normal forces have been evaluated */
2810 svmul(-rotg->k*wj*fac, pj, tmp_f);
2811 copy_rvec(tmp_f, erg->f_rot_loc[j]);
2814 /* If requested, also calculate the potential for a set of angles
2815 * near the current reference angle */
2818 for (ifit = 0; ifit < rotg->PotAngle_nstep; ifit++)
2820 /* Index of this rotation group atom with respect to the whole rotation group */
2821 jj = erg->xc_ref_ind[j];
2823 /* Rotate with the alternative angle. Like rotate_local_reference(),
2824 * just for a single local atom */
2825 mvmul(erg->PotAngleFit->rotmat[ifit], rotg->x_ref[jj], fit_tmpvec); /* fit_tmpvec = Omega*(yj0-u) */
2827 /* Calculate Omega.(yj0-u) */
2828 cprod(rotg->vec, fit_tmpvec, tmpvec); /* tmpvec = v x Omega.(yj0-u) */
2829 /* * v x Omega.(yj0-u) */
2830 unitv(tmpvec, pj); /* pj = --------------------- */
2831 /* | v x Omega.(yj0-u) | */
2833 fac = iprod(pj, xj_u); /* fac = pj.(xj-u) */
2836 /* Add to the rotation potential for this angle: */
2837 erg->PotAngleFit->V[ifit] += 0.5*rotg->k*wj*fac2;
2843 /* Add to the torque of this rotation group */
2844 erg->torque_v += torque(rotg->vec, tmp_f, erg->x_loc_pbc[j], erg->xc_center);
2846 /* Calculate the angle between reference and actual rotation group atom. */
2847 angle(rotg, xj_u, erg->xr_loc[j], &alpha, &weight); /* angle in rad, weighted */
2848 erg->angle_v += alpha * weight;
2849 erg->weight_v += weight;
2854 } /* end of loop over local rotation group atoms */
2855 erg->V = 0.5*rotg->k*sum;
2859 /* Calculate the radial motion pivot-free potential and forces */
2860 static void do_radial_motion_pf(
2861 t_rotgrp *rotg, /* The rotation group */
2862 rvec x[], /* The positions */
2863 matrix box, /* The simulation box */
2864 gmx_bool bOutstepRot, /* Output to main rotation output file */
2865 gmx_bool bOutstepSlab) /* Output per-slab data */
2867 int i, ii, iigrp, ifit, j;
2868 rvec xj; /* Current position */
2869 rvec xj_xc; /* xj - xc */
2870 rvec yj0_yc0; /* yj0 - yc0 */
2871 rvec tmp_f; /* Force */
2872 real alpha; /* a single angle between an actual and a reference position */
2873 real weight; /* single weight for a single angle */
2874 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
2875 rvec tmpvec, tmpvec2;
2876 rvec innersumvec; /* Precalculation of the inner sum */
2878 real fac, fac2, V = 0.0;
2880 gmx_bool bCalcPotFit;
2882 /* For mass weighting: */
2883 real mj, wi, wj; /* Mass-weighting of the positions */
2887 erg = rotg->enfrotgrp;
2888 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT == rotg->eFittype);
2890 N_M = rotg->nat * erg->invmass;
2892 /* Get the current center of the rotation group: */
2893 get_center(erg->xc, erg->mc, rotg->nat, erg->xc_center);
2895 /* Precalculate Sum_i [ wi qi.(xi-xc) qi ] which is needed for every single j */
2896 clear_rvec(innersumvec);
2897 for (i = 0; i < rotg->nat; i++)
2899 /* Mass-weighting */
2900 wi = N_M*erg->mc[i];
2902 /* Calculate qi. Note that xc_ref_center has already been subtracted from
2903 * x_ref in init_rot_group.*/
2904 mvmul(erg->rotmat, rotg->x_ref[i], tmpvec); /* tmpvec = Omega.(yi0-yc0) */
2906 cprod(rotg->vec, tmpvec, tmpvec2); /* tmpvec2 = v x Omega.(yi0-yc0) */
2908 /* * v x Omega.(yi0-yc0) */
2909 unitv(tmpvec2, qi); /* qi = ----------------------- */
2910 /* | v x Omega.(yi0-yc0) | */
2912 rvec_sub(erg->xc[i], erg->xc_center, tmpvec); /* tmpvec = xi-xc */
2914 svmul(wi*iprod(qi, tmpvec), qi, tmpvec2);
2916 rvec_inc(innersumvec, tmpvec2);
2918 svmul(rotg->k*erg->invmass, innersumvec, innersumveckM);
2920 /* Each process calculates the forces on its local atoms */
2921 for (j = 0; j < erg->nat_loc; j++)
2923 /* Local index of a rotation group atom */
2924 ii = erg->ind_loc[j];
2925 /* Position of this atom in the collective array */
2926 iigrp = erg->xc_ref_ind[j];
2927 /* Mass-weighting */
2928 mj = erg->mc[iigrp]; /* need the unsorted mass here */
2931 /* Current position of this atom: x[ii][XX/YY/ZZ] */
2932 copy_rvec(x[ii], xj);
2934 /* Shift this atom such that it is near its reference */
2935 shift_single_coord(box, xj, erg->xc_shifts[iigrp]);
2937 /* The (unrotated) reference position is yj0. yc0 has already
2938 * been subtracted in init_rot_group */
2939 copy_rvec(rotg->x_ref[iigrp], yj0_yc0); /* yj0_yc0 = yj0 - yc0 */
2941 /* Calculate Omega.(yj0-yc0) */
2942 mvmul(erg->rotmat, yj0_yc0, tmpvec2); /* tmpvec2 = Omega.(yj0 - yc0) */
2944 cprod(rotg->vec, tmpvec2, tmpvec); /* tmpvec = v x Omega.(yj0-yc0) */
2946 /* * v x Omega.(yj0-yc0) */
2947 unitv(tmpvec, qj); /* qj = ----------------------- */
2948 /* | v x Omega.(yj0-yc0) | */
2950 /* Calculate (xj-xc) */
2951 rvec_sub(xj, erg->xc_center, xj_xc); /* xj_xc = xj-xc */
2953 fac = iprod(qj, xj_xc); /* fac = qj.(xj-xc) */
2956 /* Store the additional force so that it can be added to the force
2957 * array after the normal forces have been evaluated */
2958 svmul(-rotg->k*wj*fac, qj, tmp_f); /* part 1 of force */
2959 svmul(mj, innersumveckM, tmpvec); /* part 2 of force */
2960 rvec_inc(tmp_f, tmpvec);
2961 copy_rvec(tmp_f, erg->f_rot_loc[j]);
2964 /* If requested, also calculate the potential for a set of angles
2965 * near the current reference angle */
2968 for (ifit = 0; ifit < rotg->PotAngle_nstep; ifit++)
2970 /* Rotate with the alternative angle. Like rotate_local_reference(),
2971 * just for a single local atom */
2972 mvmul(erg->PotAngleFit->rotmat[ifit], yj0_yc0, tmpvec2); /* tmpvec2 = Omega*(yj0-yc0) */
2974 /* Calculate Omega.(yj0-u) */
2975 cprod(rotg->vec, tmpvec2, tmpvec); /* tmpvec = v x Omega.(yj0-yc0) */
2976 /* * v x Omega.(yj0-yc0) */
2977 unitv(tmpvec, qj); /* qj = ----------------------- */
2978 /* | v x Omega.(yj0-yc0) | */
2980 fac = iprod(qj, xj_xc); /* fac = qj.(xj-xc) */
2983 /* Add to the rotation potential for this angle: */
2984 erg->PotAngleFit->V[ifit] += 0.5*rotg->k*wj*fac2;
2990 /* Add to the torque of this rotation group */
2991 erg->torque_v += torque(rotg->vec, tmp_f, xj, erg->xc_center);
2993 /* Calculate the angle between reference and actual rotation group atom. */
2994 angle(rotg, xj_xc, yj0_yc0, &alpha, &weight); /* angle in rad, weighted */
2995 erg->angle_v += alpha * weight;
2996 erg->weight_v += weight;
3001 } /* end of loop over local rotation group atoms */
3002 erg->V = 0.5*rotg->k*V;
3006 /* Precalculate the inner sum for the radial motion 2 forces */
3007 static void radial_motion2_precalc_inner_sum(t_rotgrp *rotg, rvec innersumvec)
3010 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
3011 rvec xi_xc; /* xj - xc */
3012 rvec tmpvec, tmpvec2;
3016 rvec v_xi_xc; /* v x (xj - u) */
3017 real psii, psiistar;
3018 real wi; /* Mass-weighting of the positions */
3022 erg = rotg->enfrotgrp;
3023 N_M = rotg->nat * erg->invmass;
3025 /* Loop over the collective set of positions */
3027 for (i = 0; i < rotg->nat; i++)
3029 /* Mass-weighting */
3030 wi = N_M*erg->mc[i];
3032 rvec_sub(erg->xc[i], erg->xc_center, xi_xc); /* xi_xc = xi-xc */
3034 /* Calculate ri. Note that xc_ref_center has already been subtracted from
3035 * x_ref in init_rot_group.*/
3036 mvmul(erg->rotmat, rotg->x_ref[i], ri); /* ri = Omega.(yi0-yc0) */
3038 cprod(rotg->vec, xi_xc, v_xi_xc); /* v_xi_xc = v x (xi-u) */
3040 fac = norm2(v_xi_xc);
3042 psiistar = 1.0/(fac + rotg->eps); /* psiistar = --------------------- */
3043 /* |v x (xi-xc)|^2 + eps */
3045 psii = gmx_invsqrt(fac); /* 1 */
3046 /* psii = ------------- */
3049 svmul(psii, v_xi_xc, si); /* si = psii * (v x (xi-xc) ) */
3051 fac = iprod(v_xi_xc, ri); /* fac = (v x (xi-xc)).ri */
3054 siri = iprod(si, ri); /* siri = si.ri */
3056 svmul(psiistar/psii, ri, tmpvec);
3057 svmul(psiistar*psiistar/(psii*psii*psii) * siri, si, tmpvec2);
3058 rvec_dec(tmpvec, tmpvec2);
3059 cprod(tmpvec, rotg->vec, tmpvec2);
3061 svmul(wi*siri, tmpvec2, tmpvec);
3063 rvec_inc(sumvec, tmpvec);
3065 svmul(rotg->k*erg->invmass, sumvec, innersumvec);
3069 /* Calculate the radial motion 2 potential and forces */
3070 static void do_radial_motion2(
3071 t_rotgrp *rotg, /* The rotation group */
3072 rvec x[], /* The positions */
3073 matrix box, /* The simulation box */
3074 gmx_bool bOutstepRot, /* Output to main rotation output file */
3075 gmx_bool bOutstepSlab) /* Output per-slab data */
3077 int ii, iigrp, ifit, j;
3078 rvec xj; /* Position */
3079 real alpha; /* a single angle between an actual and a reference position */
3080 real weight; /* single weight for a single angle */
3081 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
3082 rvec xj_u; /* xj - u */
3083 rvec yj0_yc0; /* yj0 -yc0 */
3084 rvec tmpvec, tmpvec2;
3085 real fac, fit_fac, fac2, Vpart = 0.0;
3086 rvec rj, fit_rj, sj;
3088 rvec v_xj_u; /* v x (xj - u) */
3089 real psij, psijstar;
3090 real mj, wj; /* For mass-weighting of the positions */
3094 gmx_bool bCalcPotFit;
3097 erg = rotg->enfrotgrp;
3099 bPF = rotg->eType == erotgRM2PF;
3100 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT == rotg->eFittype);
3103 clear_rvec(yj0_yc0); /* Make the compiler happy */
3105 clear_rvec(innersumvec);
3108 /* For the pivot-free variant we have to use the current center of
3109 * mass of the rotation group instead of the pivot u */
3110 get_center(erg->xc, erg->mc, rotg->nat, erg->xc_center);
3112 /* Also, we precalculate the second term of the forces that is identical
3113 * (up to the weight factor mj) for all forces */
3114 radial_motion2_precalc_inner_sum(rotg, innersumvec);
3117 N_M = rotg->nat * erg->invmass;
3119 /* Each process calculates the forces on its local atoms */
3120 for (j = 0; j < erg->nat_loc; j++)
3124 /* Local index of a rotation group atom */
3125 ii = erg->ind_loc[j];
3126 /* Position of this atom in the collective array */
3127 iigrp = erg->xc_ref_ind[j];
3128 /* Mass-weighting */
3129 mj = erg->mc[iigrp];
3131 /* Current position of this atom: x[ii] */
3132 copy_rvec(x[ii], xj);
3134 /* Shift this atom such that it is near its reference */
3135 shift_single_coord(box, xj, erg->xc_shifts[iigrp]);
3137 /* The (unrotated) reference position is yj0. yc0 has already
3138 * been subtracted in init_rot_group */
3139 copy_rvec(rotg->x_ref[iigrp], yj0_yc0); /* yj0_yc0 = yj0 - yc0 */
3141 /* Calculate Omega.(yj0-yc0) */
3142 mvmul(erg->rotmat, yj0_yc0, rj); /* rj = Omega.(yj0-yc0) */
3147 copy_rvec(erg->x_loc_pbc[j], xj);
3148 copy_rvec(erg->xr_loc[j], rj); /* rj = Omega.(yj0-u) */
3150 /* Mass-weighting */
3153 /* Calculate (xj-u) resp. (xj-xc) */
3154 rvec_sub(xj, erg->xc_center, xj_u); /* xj_u = xj-u */
3156 cprod(rotg->vec, xj_u, v_xj_u); /* v_xj_u = v x (xj-u) */
3158 fac = norm2(v_xj_u);
3160 psijstar = 1.0/(fac + rotg->eps); /* psistar = -------------------- */
3161 /* |v x (xj-u)|^2 + eps */
3163 psij = gmx_invsqrt(fac); /* 1 */
3164 /* psij = ------------ */
3167 svmul(psij, v_xj_u, sj); /* sj = psij * (v x (xj-u) ) */
3169 fac = iprod(v_xj_u, rj); /* fac = (v x (xj-u)).rj */
3172 sjrj = iprod(sj, rj); /* sjrj = sj.rj */
3174 svmul(psijstar/psij, rj, tmpvec);
3175 svmul(psijstar*psijstar/(psij*psij*psij) * sjrj, sj, tmpvec2);
3176 rvec_dec(tmpvec, tmpvec2);
3177 cprod(tmpvec, rotg->vec, tmpvec2);
3179 /* Store the additional force so that it can be added to the force
3180 * array after the normal forces have been evaluated */
3181 svmul(-rotg->k*wj*sjrj, tmpvec2, tmpvec);
3182 svmul(mj, innersumvec, tmpvec2); /* This is != 0 only for the pivot-free variant */
3184 rvec_add(tmpvec2, tmpvec, erg->f_rot_loc[j]);
3185 Vpart += wj*psijstar*fac2;
3187 /* If requested, also calculate the potential for a set of angles
3188 * near the current reference angle */
3191 for (ifit = 0; ifit < rotg->PotAngle_nstep; ifit++)
3195 mvmul(erg->PotAngleFit->rotmat[ifit], yj0_yc0, fit_rj); /* fit_rj = Omega.(yj0-yc0) */
3199 /* Position of this atom in the collective array */
3200 iigrp = erg->xc_ref_ind[j];
3201 /* Rotate with the alternative angle. Like rotate_local_reference(),
3202 * just for a single local atom */
3203 mvmul(erg->PotAngleFit->rotmat[ifit], rotg->x_ref[iigrp], fit_rj); /* fit_rj = Omega*(yj0-u) */
3205 fit_fac = iprod(v_xj_u, fit_rj); /* fac = (v x (xj-u)).fit_rj */
3206 /* Add to the rotation potential for this angle: */
3207 erg->PotAngleFit->V[ifit] += 0.5*rotg->k*wj*psijstar*fit_fac*fit_fac;
3213 /* Add to the torque of this rotation group */
3214 erg->torque_v += torque(rotg->vec, erg->f_rot_loc[j], xj, erg->xc_center);
3216 /* Calculate the angle between reference and actual rotation group atom. */
3217 angle(rotg, xj_u, rj, &alpha, &weight); /* angle in rad, weighted */
3218 erg->angle_v += alpha * weight;
3219 erg->weight_v += weight;
3224 } /* end of loop over local rotation group atoms */
3225 erg->V = 0.5*rotg->k*Vpart;
3229 /* Determine the smallest and largest position vector (with respect to the
3230 * rotation vector) for the reference group */
3231 static void get_firstlast_atom_ref(
3236 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
3238 real xcproj; /* The projection of a reference position on the
3240 real minproj, maxproj; /* Smallest and largest projection on v */
3244 erg = rotg->enfrotgrp;
3246 /* Start with some value */
3247 minproj = iprod(rotg->x_ref[0], rotg->vec);
3250 /* This is just to ensure that it still works if all the atoms of the
3251 * reference structure are situated in a plane perpendicular to the rotation
3254 *lastindex = rotg->nat-1;
3256 /* Loop over all atoms of the reference group,
3257 * project them on the rotation vector to find the extremes */
3258 for (i = 0; i < rotg->nat; i++)
3260 xcproj = iprod(rotg->x_ref[i], rotg->vec);
3261 if (xcproj < minproj)
3266 if (xcproj > maxproj)
3275 /* Allocate memory for the slabs */
3276 static void allocate_slabs(
3282 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
3286 erg = rotg->enfrotgrp;
3288 /* More slabs than are defined for the reference are never needed */
3289 nslabs = erg->slab_last_ref - erg->slab_first_ref + 1;
3291 /* Remember how many we allocated */
3292 erg->nslabs_alloc = nslabs;
3294 if ( (NULL != fplog) && bVerbose)
3296 fprintf(fplog, "%s allocating memory to store data for %d slabs (rotation group %d).\n",
3299 snew(erg->slab_center, nslabs);
3300 snew(erg->slab_center_ref, nslabs);
3301 snew(erg->slab_weights, nslabs);
3302 snew(erg->slab_torque_v, nslabs);
3303 snew(erg->slab_data, nslabs);
3304 snew(erg->gn_atom, nslabs);
3305 snew(erg->gn_slabind, nslabs);
3306 snew(erg->slab_innersumvec, nslabs);
3307 for (i = 0; i < nslabs; i++)
3309 snew(erg->slab_data[i].x, rotg->nat);
3310 snew(erg->slab_data[i].ref, rotg->nat);
3311 snew(erg->slab_data[i].weight, rotg->nat);
3313 snew(erg->xc_ref_sorted, rotg->nat);
3314 snew(erg->xc_sortind, rotg->nat);
3315 snew(erg->firstatom, nslabs);
3316 snew(erg->lastatom, nslabs);
3320 /* From the extreme positions of the reference group, determine the first
3321 * and last slab of the reference. We can never have more slabs in the real
3322 * simulation than calculated here for the reference.
3324 static void get_firstlast_slab_ref(t_rotgrp *rotg, real mc[], int ref_firstindex, int ref_lastindex)
3326 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
3331 erg = rotg->enfrotgrp;
3332 first = get_first_slab(rotg, erg->max_beta, rotg->x_ref[ref_firstindex]);
3333 last = get_last_slab( rotg, erg->max_beta, rotg->x_ref[ref_lastindex ]);
3335 while (get_slab_weight(first, rotg, rotg->x_ref, mc, &dummy) > WEIGHT_MIN)
3339 erg->slab_first_ref = first+1;
3340 while (get_slab_weight(last, rotg, rotg->x_ref, mc, &dummy) > WEIGHT_MIN)
3344 erg->slab_last_ref = last-1;
3348 /* Special version of copy_rvec:
3349 * During the copy procedure of xcurr to b, the correct PBC image is chosen
3350 * such that the copied vector ends up near its reference position xref */
3351 static gmx_inline void copy_correct_pbc_image(
3352 const rvec xcurr, /* copy vector xcurr ... */
3353 rvec b, /* ... to b ... */
3354 const rvec xref, /* choosing the PBC image such that b ends up near xref */
3363 /* Shortest PBC distance between the atom and its reference */
3364 rvec_sub(xcurr, xref, dx);
3366 /* Determine the shift for this atom */
3368 for (m = npbcdim-1; m >= 0; m--)
3370 while (dx[m] < -0.5*box[m][m])
3372 for (d = 0; d < DIM; d++)
3378 while (dx[m] >= 0.5*box[m][m])
3380 for (d = 0; d < DIM; d++)
3388 /* Apply the shift to the position */
3389 copy_rvec(xcurr, b);
3390 shift_single_coord(box, b, shift);
3394 static void init_rot_group(FILE *fplog, t_commrec *cr, int g, t_rotgrp *rotg,
3395 rvec *x, gmx_mtop_t *mtop, gmx_bool bVerbose, FILE *out_slabs, matrix box,
3396 t_inputrec *ir, gmx_bool bOutputCenters)
3399 rvec coord, xref, *xdum;
3400 gmx_bool bFlex, bColl;
3402 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
3403 int ref_firstindex, ref_lastindex;
3404 gmx_mtop_atomlookup_t alook = NULL;
3405 real mass, totalmass;
3410 /* Do we have a flexible axis? */
3411 bFlex = ISFLEX(rotg);
3412 /* Do we use a global set of coordinates? */
3413 bColl = ISCOLL(rotg);
3415 erg = rotg->enfrotgrp;
3417 /* Allocate space for collective coordinates if needed */
3420 snew(erg->xc, rotg->nat);
3421 snew(erg->xc_shifts, rotg->nat);
3422 snew(erg->xc_eshifts, rotg->nat);
3423 snew(erg->xc_old, rotg->nat);
3425 if (rotg->eFittype == erotgFitNORM)
3427 snew(erg->xc_ref_length, rotg->nat); /* in case fit type NORM is chosen */
3428 snew(erg->xc_norm, rotg->nat);
3433 snew(erg->xr_loc, rotg->nat);
3434 snew(erg->x_loc_pbc, rotg->nat);
3437 snew(erg->f_rot_loc, rotg->nat);
3438 snew(erg->xc_ref_ind, rotg->nat);
3440 /* Make space for the calculation of the potential at other angles (used
3441 * for fitting only) */
3442 if (erotgFitPOT == rotg->eFittype)
3444 snew(erg->PotAngleFit, 1);
3445 snew(erg->PotAngleFit->degangle, rotg->PotAngle_nstep);
3446 snew(erg->PotAngleFit->V, rotg->PotAngle_nstep);
3447 snew(erg->PotAngleFit->rotmat, rotg->PotAngle_nstep);
3449 /* Get the set of angles around the reference angle */
3450 start = -0.5 * (rotg->PotAngle_nstep - 1)*rotg->PotAngle_step;
3451 for (i = 0; i < rotg->PotAngle_nstep; i++)
3453 erg->PotAngleFit->degangle[i] = start + i*rotg->PotAngle_step;
3458 erg->PotAngleFit = NULL;
3461 /* xc_ref_ind needs to be set to identity in the serial case */
3464 for (i = 0; i < rotg->nat; i++)
3466 erg->xc_ref_ind[i] = i;
3470 /* Copy the masses so that the center can be determined. For all types of
3471 * enforced rotation, we store the masses in the erg->mc array. */
3474 alook = gmx_mtop_atomlookup_init(mtop);
3476 snew(erg->mc, rotg->nat);
3479 snew(erg->mc_sorted, rotg->nat);
3483 snew(erg->m_loc, rotg->nat);
3486 for (i = 0; i < rotg->nat; i++)
3490 gmx_mtop_atomnr_to_atom(alook, rotg->ind[i], &atom);
3500 erg->invmass = 1.0/totalmass;
3504 gmx_mtop_atomlookup_destroy(alook);
3507 /* Set xc_ref_center for any rotation potential */
3508 if ((rotg->eType == erotgISO) || (rotg->eType == erotgPM) || (rotg->eType == erotgRM) || (rotg->eType == erotgRM2))
3510 /* Set the pivot point for the fixed, stationary-axis potentials. This
3511 * won't change during the simulation */
3512 copy_rvec(rotg->pivot, erg->xc_ref_center);
3513 copy_rvec(rotg->pivot, erg->xc_center );
3517 /* Center of the reference positions */
3518 get_center(rotg->x_ref, erg->mc, rotg->nat, erg->xc_ref_center);
3520 /* Center of the actual positions */
3523 snew(xdum, rotg->nat);
3524 for (i = 0; i < rotg->nat; i++)
3527 copy_rvec(x[ii], xdum[i]);
3529 get_center(xdum, erg->mc, rotg->nat, erg->xc_center);
3535 gmx_bcast(sizeof(erg->xc_center), erg->xc_center, cr);
3542 /* Save the original (whole) set of positions in xc_old such that at later
3543 * steps the rotation group can always be made whole again. If the simulation is
3544 * restarted, we compute the starting reference positions (given the time)
3545 * and assume that the correct PBC image of each position is the one nearest
3546 * to the current reference */
3549 /* Calculate the rotation matrix for this angle: */
3550 t_start = ir->init_t + ir->init_step*ir->delta_t;
3551 erg->degangle = rotg->rate * t_start;
3552 calc_rotmat(rotg->vec, erg->degangle, erg->rotmat);
3554 for (i = 0; i < rotg->nat; i++)
3558 /* Subtract pivot, rotate, and add pivot again. This will yield the
3559 * reference position for time t */
3560 rvec_sub(rotg->x_ref[i], erg->xc_ref_center, coord);
3561 mvmul(erg->rotmat, coord, xref);
3562 rvec_inc(xref, erg->xc_ref_center);
3564 copy_correct_pbc_image(x[ii], erg->xc_old[i], xref, box, 3);
3570 gmx_bcast(rotg->nat*sizeof(erg->xc_old[0]), erg->xc_old, cr);
3575 if ( (rotg->eType != erotgFLEX) && (rotg->eType != erotgFLEX2) )
3577 /* Put the reference positions into origin: */
3578 for (i = 0; i < rotg->nat; i++)
3580 rvec_dec(rotg->x_ref[i], erg->xc_ref_center);
3584 /* Enforced rotation with flexible axis */
3587 /* Calculate maximum beta value from minimum gaussian (performance opt.) */
3588 erg->max_beta = calc_beta_max(rotg->min_gaussian, rotg->slab_dist);
3590 /* Determine the smallest and largest coordinate with respect to the rotation vector */
3591 get_firstlast_atom_ref(rotg, &ref_firstindex, &ref_lastindex);
3593 /* From the extreme positions of the reference group, determine the first
3594 * and last slab of the reference. */
3595 get_firstlast_slab_ref(rotg, erg->mc, ref_firstindex, ref_lastindex);
3597 /* Allocate memory for the slabs */
3598 allocate_slabs(rotg, fplog, g, bVerbose);
3600 /* Flexible rotation: determine the reference centers for the rest of the simulation */
3601 erg->slab_first = erg->slab_first_ref;
3602 erg->slab_last = erg->slab_last_ref;
3603 get_slab_centers(rotg, rotg->x_ref, erg->mc, g, -1, out_slabs, bOutputCenters, TRUE);
3605 /* Length of each x_rotref vector from center (needed if fit routine NORM is chosen): */
3606 if (rotg->eFittype == erotgFitNORM)
3608 for (i = 0; i < rotg->nat; i++)
3610 rvec_sub(rotg->x_ref[i], erg->xc_ref_center, coord);
3611 erg->xc_ref_length[i] = norm(coord);
3618 extern void dd_make_local_rotation_groups(gmx_domdec_t *dd, t_rot *rot)
3623 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
3627 for (g = 0; g < rot->ngrp; g++)
3629 rotg = &rot->grp[g];
3630 erg = rotg->enfrotgrp;
3633 dd_make_local_group_indices(ga2la, rotg->nat, rotg->ind,
3634 &erg->nat_loc, &erg->ind_loc, &erg->nalloc_loc, erg->xc_ref_ind);
3639 /* Calculate the size of the MPI buffer needed in reduce_output() */
3640 static int calc_mpi_bufsize(t_rot *rot)
3643 int count_group, count_total;
3645 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
3649 for (g = 0; g < rot->ngrp; g++)
3651 rotg = &rot->grp[g];
3652 erg = rotg->enfrotgrp;
3654 /* Count the items that are transferred for this group: */
3655 count_group = 4; /* V, torque, angle, weight */
3657 /* Add the maximum number of slabs for flexible groups */
3660 count_group += erg->slab_last_ref - erg->slab_first_ref + 1;
3663 /* Add space for the potentials at different angles: */
3664 if (erotgFitPOT == rotg->eFittype)
3666 count_group += rotg->PotAngle_nstep;
3669 /* Add to the total number: */
3670 count_total += count_group;
3677 extern void init_rot(FILE *fplog, t_inputrec *ir, int nfile, const t_filenm fnm[],
3678 t_commrec *cr, rvec *x, matrix box, gmx_mtop_t *mtop, const output_env_t oenv,
3679 gmx_bool bVerbose, unsigned long Flags)
3684 int nat_max = 0; /* Size of biggest rotation group */
3685 gmx_enfrot_t er; /* Pointer to the enforced rotation buffer variables */
3686 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
3687 rvec *x_pbc = NULL; /* Space for the pbc-correct atom positions */
3690 if (MASTER(cr) && bVerbose)
3692 fprintf(stdout, "%s Initializing ...\n", RotStr);
3696 snew(rot->enfrot, 1);
3700 /* When appending, skip first output to avoid duplicate entries in the data files */
3701 if (er->Flags & MD_APPENDFILES)
3710 if (MASTER(cr) && er->bOut)
3712 please_cite(fplog, "Kutzner2011");
3715 /* Output every step for reruns */
3716 if (er->Flags & MD_RERUN)
3720 fprintf(fplog, "%s rerun - will write rotation output every available step.\n", RotStr);
3726 er->out_slabs = NULL;
3727 if (MASTER(cr) && HaveFlexibleGroups(rot) )
3729 er->out_slabs = open_slab_out(opt2fn("-rs", nfile, fnm), rot);
3734 /* Remove pbc, make molecule whole.
3735 * When ir->bContinuation=TRUE this has already been done, but ok. */
3736 snew(x_pbc, mtop->natoms);
3737 m_rveccopy(mtop->natoms, x, x_pbc);
3738 do_pbc_first_mtop(NULL, ir->ePBC, box, mtop, x_pbc);
3739 /* All molecules will be whole now, but not necessarily in the home box.
3740 * Additionally, if a rotation group consists of more than one molecule
3741 * (e.g. two strands of DNA), each one of them can end up in a different
3742 * periodic box. This is taken care of in init_rot_group. */
3745 for (g = 0; g < rot->ngrp; g++)
3747 rotg = &rot->grp[g];
3751 fprintf(fplog, "%s group %d type '%s'\n", RotStr, g, erotg_names[rotg->eType]);
3756 /* Allocate space for the rotation group's data: */
3757 snew(rotg->enfrotgrp, 1);
3758 erg = rotg->enfrotgrp;
3760 nat_max = max(nat_max, rotg->nat);
3765 erg->nalloc_loc = 0;
3766 erg->ind_loc = NULL;
3770 erg->nat_loc = rotg->nat;
3771 erg->ind_loc = rotg->ind;
3773 init_rot_group(fplog, cr, g, rotg, x_pbc, mtop, bVerbose, er->out_slabs, box, ir,
3774 !(er->Flags & MD_APPENDFILES) ); /* Do not output the reference centers
3775 * again if we are appending */
3779 /* Allocate space for enforced rotation buffer variables */
3780 er->bufsize = nat_max;
3781 snew(er->data, nat_max);
3782 snew(er->xbuf, nat_max);
3783 snew(er->mbuf, nat_max);
3785 /* Buffers for MPI reducing torques, angles, weights (for each group), and V */
3788 er->mpi_bufsize = calc_mpi_bufsize(rot) + 100; /* larger to catch errors */
3789 snew(er->mpi_inbuf, er->mpi_bufsize);
3790 snew(er->mpi_outbuf, er->mpi_bufsize);
3794 er->mpi_bufsize = 0;
3795 er->mpi_inbuf = NULL;
3796 er->mpi_outbuf = NULL;
3799 /* Only do I/O on the MASTER */
3800 er->out_angles = NULL;
3802 er->out_torque = NULL;
3805 er->out_rot = open_rot_out(opt2fn("-ro", nfile, fnm), rot, oenv);
3807 if (rot->nstsout > 0)
3809 if (HaveFlexibleGroups(rot) || HavePotFitGroups(rot) )
3811 er->out_angles = open_angles_out(opt2fn("-ra", nfile, fnm), rot);
3813 if (HaveFlexibleGroups(rot) )
3815 er->out_torque = open_torque_out(opt2fn("-rt", nfile, fnm), rot);
3824 extern void finish_rot(t_rot *rot)
3826 gmx_enfrot_t er; /* Pointer to the enforced rotation buffer variables */
3832 gmx_fio_fclose(er->out_rot);
3836 gmx_fio_fclose(er->out_slabs);
3840 gmx_fio_fclose(er->out_angles);
3844 gmx_fio_fclose(er->out_torque);
3849 /* Rotate the local reference positions and store them in
3850 * erg->xr_loc[0...(nat_loc-1)]
3852 * Note that we already subtracted u or y_c from the reference positions
3853 * in init_rot_group().
3855 static void rotate_local_reference(t_rotgrp *rotg)
3857 gmx_enfrotgrp_t erg;
3861 erg = rotg->enfrotgrp;
3863 for (i = 0; i < erg->nat_loc; i++)
3865 /* Index of this rotation group atom with respect to the whole rotation group */
3866 ii = erg->xc_ref_ind[i];
3868 mvmul(erg->rotmat, rotg->x_ref[ii], erg->xr_loc[i]);
3873 /* Select the PBC representation for each local x position and store that
3874 * for later usage. We assume the right PBC image of an x is the one nearest to
3875 * its rotated reference */
3876 static void choose_pbc_image(rvec x[], t_rotgrp *rotg, matrix box, int npbcdim)
3879 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
3883 erg = rotg->enfrotgrp;
3885 for (i = 0; i < erg->nat_loc; i++)
3887 /* Index of a rotation group atom */
3888 ii = erg->ind_loc[i];
3890 /* Get the correctly rotated reference position. The pivot was already
3891 * subtracted in init_rot_group() from the reference positions. Also,
3892 * the reference positions have already been rotated in
3893 * rotate_local_reference(). For the current reference position we thus
3894 * only need to add the pivot again. */
3895 copy_rvec(erg->xr_loc[i], xref);
3896 rvec_inc(xref, erg->xc_ref_center);
3898 copy_correct_pbc_image(x[ii], erg->x_loc_pbc[i], xref, box, npbcdim);
3903 extern void do_rotation(
3910 gmx_wallcycle_t wcycle,
3916 gmx_bool outstep_slab, outstep_rot;
3917 gmx_bool bFlex, bColl;
3918 gmx_enfrot_t er; /* Pointer to the enforced rotation buffer variables */
3919 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
3921 t_gmx_potfit *fit = NULL; /* For fit type 'potential' determine the fit
3922 angle via the potential minimum */
3924 /* Enforced rotation cycle counting: */
3925 gmx_cycles_t cycles_comp; /* Cycles for the enf. rotation computation
3926 only, does not count communication. This
3927 counter is used for load-balancing */
3936 /* When to output in main rotation output file */
3937 outstep_rot = do_per_step(step, rot->nstrout) && er->bOut;
3938 /* When to output per-slab data */
3939 outstep_slab = do_per_step(step, rot->nstsout) && er->bOut;
3941 /* Output time into rotation output file */
3942 if (outstep_rot && MASTER(cr))
3944 fprintf(er->out_rot, "%12.3e", t);
3947 /**************************************************************************/
3948 /* First do ALL the communication! */
3949 for (g = 0; g < rot->ngrp; g++)
3951 rotg = &rot->grp[g];
3952 erg = rotg->enfrotgrp;
3954 /* Do we have a flexible axis? */
3955 bFlex = ISFLEX(rotg);
3956 /* Do we use a collective (global) set of coordinates? */
3957 bColl = ISCOLL(rotg);
3959 /* Calculate the rotation matrix for this angle: */
3960 erg->degangle = rotg->rate * t;
3961 calc_rotmat(rotg->vec, erg->degangle, erg->rotmat);
3965 /* Transfer the rotation group's positions such that every node has
3966 * all of them. Every node contributes its local positions x and stores
3967 * it in the collective erg->xc array. */
3968 communicate_group_positions(cr, erg->xc, erg->xc_shifts, erg->xc_eshifts, bNS,
3969 x, rotg->nat, erg->nat_loc, erg->ind_loc, erg->xc_ref_ind, erg->xc_old, box);
3973 /* Fill the local masses array;
3974 * this array changes in DD/neighborsearching steps */
3977 for (i = 0; i < erg->nat_loc; i++)
3979 /* Index of local atom w.r.t. the collective rotation group */
3980 ii = erg->xc_ref_ind[i];
3981 erg->m_loc[i] = erg->mc[ii];
3985 /* Calculate Omega*(y_i-y_c) for the local positions */
3986 rotate_local_reference(rotg);
3988 /* Choose the nearest PBC images of the group atoms with respect
3989 * to the rotated reference positions */
3990 choose_pbc_image(x, rotg, box, 3);
3992 /* Get the center of the rotation group */
3993 if ( (rotg->eType == erotgISOPF) || (rotg->eType == erotgPMPF) )
3995 get_center_comm(cr, erg->x_loc_pbc, erg->m_loc, erg->nat_loc, rotg->nat, erg->xc_center);
3999 } /* End of loop over rotation groups */
4001 /**************************************************************************/
4002 /* Done communicating, we can start to count cycles for the load balancing now ... */
4003 cycles_comp = gmx_cycles_read();
4010 for (g = 0; g < rot->ngrp; g++)
4012 rotg = &rot->grp[g];
4013 erg = rotg->enfrotgrp;
4015 bFlex = ISFLEX(rotg);
4016 bColl = ISCOLL(rotg);
4018 if (outstep_rot && MASTER(cr))
4020 fprintf(er->out_rot, "%12.4f", erg->degangle);
4023 /* Calculate angles and rotation matrices for potential fitting: */
4024 if ( (outstep_rot || outstep_slab) && (erotgFitPOT == rotg->eFittype) )
4026 fit = erg->PotAngleFit;
4027 for (i = 0; i < rotg->PotAngle_nstep; i++)
4029 calc_rotmat(rotg->vec, erg->degangle + fit->degangle[i], fit->rotmat[i]);
4031 /* Clear value from last step */
4032 erg->PotAngleFit->V[i] = 0.0;
4036 /* Clear values from last time step */
4038 erg->torque_v = 0.0;
4040 erg->weight_v = 0.0;
4042 switch (rotg->eType)
4048 do_fixed(rotg, outstep_rot, outstep_slab);
4051 do_radial_motion(rotg, outstep_rot, outstep_slab);
4054 do_radial_motion_pf(rotg, x, box, outstep_rot, outstep_slab);
4058 do_radial_motion2(rotg, x, box, outstep_rot, outstep_slab);
4062 /* Subtract the center of the rotation group from the collective positions array
4063 * Also store the center in erg->xc_center since it needs to be subtracted
4064 * in the low level routines from the local coordinates as well */
4065 get_center(erg->xc, erg->mc, rotg->nat, erg->xc_center);
4066 svmul(-1.0, erg->xc_center, transvec);
4067 translate_x(erg->xc, rotg->nat, transvec);
4068 do_flexible(MASTER(cr), er, rotg, g, x, box, t, outstep_rot, outstep_slab);
4072 /* Do NOT subtract the center of mass in the low level routines! */
4073 clear_rvec(erg->xc_center);
4074 do_flexible(MASTER(cr), er, rotg, g, x, box, t, outstep_rot, outstep_slab);
4077 gmx_fatal(FARGS, "No such rotation potential.");
4085 fprintf(stderr, "%s calculation (step %d) took %g seconds.\n", RotStr, step, MPI_Wtime()-t0);
4089 /* Stop the enforced rotation cycle counter and add the computation-only
4090 * cycles to the force cycles for load balancing */
4091 cycles_comp = gmx_cycles_read() - cycles_comp;
4093 if (DOMAINDECOMP(cr) && wcycle)
4095 dd_cycles_add(cr->dd, cycles_comp, ddCyclF);