3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
9 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
10 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
11 * Copyright (c) 2001-2008, The GROMACS development team,
12 * check out http://www.gromacs.org for more information.
14 * This program is free software; you can redistribute it and/or
15 * modify it under the terms of the GNU General Public License
16 * as published by the Free Software Foundation; either version 2
17 * of the License, or (at your option) any later version.
19 * If you want to redistribute modifications, please consider that
20 * scientific software is very special. Version control is crucial -
21 * bugs must be traceable. We will be happy to consider code for
22 * inclusion in the official distribution, but derived work must not
23 * be called official GROMACS. Details are found in the README & COPYING
24 * files - if they are missing, get the official version at www.gromacs.org.
26 * To help us fund GROMACS development, we humbly ask that you cite
27 * the papers on the package - you can find them in the top README file.
29 * For more info, check our website at http://www.gromacs.org
32 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
48 #include "mtop_util.h"
52 #include "gmx_ga2la.h"
54 #include "groupcoord.h"
55 #include "pull_rotation.h"
60 #include "gromacs/fileio/futil.h"
61 #include "gromacs/fileio/gmxfio.h"
62 #include "gromacs/fileio/trnio.h"
63 #include "gromacs/timing/cyclecounter.h"
64 #include "gromacs/timing/wallcycle.h"
66 static char *RotStr = {"Enforced rotation:"};
68 /* Set the minimum weight for the determination of the slab centers */
69 #define WEIGHT_MIN (10*GMX_FLOAT_MIN)
71 /* Helper structure for sorting positions along rotation vector */
73 real xcproj; /* Projection of xc on the rotation vector */
74 int ind; /* Index of xc */
76 rvec x; /* Position */
77 rvec x_ref; /* Reference position */
81 /* Enforced rotation / flexible: determine the angle of each slab */
82 typedef struct gmx_slabdata
84 int nat; /* Number of atoms belonging to this slab */
85 rvec *x; /* The positions belonging to this slab. In
86 general, this should be all positions of the
87 whole rotation group, but we leave those away
88 that have a small enough weight */
89 rvec *ref; /* Same for reference */
90 real *weight; /* The weight for each atom */
94 /* Helper structure for potential fitting */
95 typedef struct gmx_potfit
97 real *degangle; /* Set of angles for which the potential is
98 calculated. The optimum fit is determined as
99 the angle for with the potential is minimal */
100 real *V; /* Potential for the different angles */
101 matrix *rotmat; /* Rotation matrix corresponding to the angles */
105 /* Enforced rotation data for all groups */
106 typedef struct gmx_enfrot
108 FILE *out_rot; /* Output file for rotation data */
109 FILE *out_torque; /* Output file for torque data */
110 FILE *out_angles; /* Output file for slab angles for flexible type */
111 FILE *out_slabs; /* Output file for slab centers */
112 int bufsize; /* Allocation size of buf */
113 rvec *xbuf; /* Coordinate buffer variable for sorting */
114 real *mbuf; /* Masses buffer variable for sorting */
115 sort_along_vec_t *data; /* Buffer variable needed for position sorting */
116 real *mpi_inbuf; /* MPI buffer */
117 real *mpi_outbuf; /* MPI buffer */
118 int mpi_bufsize; /* Allocation size of in & outbuf */
119 unsigned long Flags; /* mdrun flags */
120 gmx_bool bOut; /* Used to skip first output when appending to
121 * avoid duplicate entries in rotation outfiles */
125 /* Global enforced rotation data for a single rotation group */
126 typedef struct gmx_enfrotgrp
128 real degangle; /* Rotation angle in degrees */
129 matrix rotmat; /* Rotation matrix */
130 atom_id *ind_loc; /* Local rotation indices */
131 int nat_loc; /* Number of local group atoms */
132 int nalloc_loc; /* Allocation size for ind_loc and weight_loc */
134 real V; /* Rotation potential for this rotation group */
135 rvec *f_rot_loc; /* Array to store the forces on the local atoms
136 resulting from enforced rotation potential */
138 /* Collective coordinates for the whole rotation group */
139 real *xc_ref_length; /* Length of each x_rotref vector after x_rotref
140 has been put into origin */
141 int *xc_ref_ind; /* Position of each local atom in the collective
143 rvec xc_center; /* Center of the rotation group positions, may
145 rvec xc_ref_center; /* dito, for the reference positions */
146 rvec *xc; /* Current (collective) positions */
147 ivec *xc_shifts; /* Current (collective) shifts */
148 ivec *xc_eshifts; /* Extra shifts since last DD step */
149 rvec *xc_old; /* Old (collective) positions */
150 rvec *xc_norm; /* Normalized form of the current positions */
151 rvec *xc_ref_sorted; /* Reference positions (sorted in the same order
152 as xc when sorted) */
153 int *xc_sortind; /* Where is a position found after sorting? */
154 real *mc; /* Collective masses */
156 real invmass; /* one over the total mass of the rotation group */
158 real torque_v; /* Torque in the direction of rotation vector */
159 real angle_v; /* Actual angle of the whole rotation group */
160 /* Fixed rotation only */
161 real weight_v; /* Weights for angle determination */
162 rvec *xr_loc; /* Local reference coords, correctly rotated */
163 rvec *x_loc_pbc; /* Local current coords, correct PBC image */
164 real *m_loc; /* Masses of the current local atoms */
166 /* Flexible rotation only */
167 int nslabs_alloc; /* For this many slabs memory is allocated */
168 int slab_first; /* Lowermost slab for that the calculation needs
169 to be performed at a given time step */
170 int slab_last; /* Uppermost slab ... */
171 int slab_first_ref; /* First slab for which ref. center is stored */
172 int slab_last_ref; /* Last ... */
173 int slab_buffer; /* Slab buffer region around reference slabs */
174 int *firstatom; /* First relevant atom for a slab */
175 int *lastatom; /* Last relevant atom for a slab */
176 rvec *slab_center; /* Gaussian-weighted slab center */
177 rvec *slab_center_ref; /* Gaussian-weighted slab center for the
178 reference positions */
179 real *slab_weights; /* Sum of gaussian weights in a slab */
180 real *slab_torque_v; /* Torque T = r x f for each slab. */
181 /* torque_v = m.v = angular momentum in the
183 real max_beta; /* min_gaussian from inputrec->rotgrp is the
184 minimum value the gaussian must have so that
185 the force is actually evaluated max_beta is
186 just another way to put it */
187 real *gn_atom; /* Precalculated gaussians for a single atom */
188 int *gn_slabind; /* Tells to which slab each precalculated gaussian
190 rvec *slab_innersumvec; /* Inner sum of the flexible2 potential per slab;
191 this is precalculated for optimization reasons */
192 t_gmx_slabdata *slab_data; /* Holds atom positions and gaussian weights
193 of atoms belonging to a slab */
195 /* For potential fits with varying angle: */
196 t_gmx_potfit *PotAngleFit; /* Used for fit type 'potential' */
200 /* Activate output of forces for correctness checks */
201 /* #define PRINT_FORCES */
203 #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]);
204 #define PRINT_POT_TAU if (MASTER(cr)) { \
205 fprintf(stderr, "potential = %15.8f\n" "torque = %15.8f\n", erg->V, erg->torque_v); \
208 #define PRINT_FORCE_J
209 #define PRINT_POT_TAU
212 /* Shortcuts for often used queries */
213 #define ISFLEX(rg) ( (rg->eType == erotgFLEX) || (rg->eType == erotgFLEXT) || (rg->eType == erotgFLEX2) || (rg->eType == erotgFLEX2T) )
214 #define ISCOLL(rg) ( (rg->eType == erotgFLEX) || (rg->eType == erotgFLEXT) || (rg->eType == erotgFLEX2) || (rg->eType == erotgFLEX2T) || (rg->eType == erotgRMPF) || (rg->eType == erotgRM2PF) )
217 /* Does any of the rotation groups use slab decomposition? */
218 static gmx_bool HaveFlexibleGroups(t_rot *rot)
224 for (g = 0; g < rot->ngrp; g++)
237 /* Is for any group the fit angle determined by finding the minimum of the
238 * rotation potential? */
239 static gmx_bool HavePotFitGroups(t_rot *rot)
245 for (g = 0; g < rot->ngrp; g++)
248 if (erotgFitPOT == rotg->eFittype)
258 static double** allocate_square_matrix(int dim)
265 for (i = 0; i < dim; i++)
274 static void free_square_matrix(double** mat, int dim)
279 for (i = 0; i < dim; i++)
287 /* Return the angle for which the potential is minimal */
288 static real get_fitangle(t_rotgrp *rotg, gmx_enfrotgrp_t erg)
291 real fitangle = -999.9;
292 real pot_min = GMX_FLOAT_MAX;
296 fit = erg->PotAngleFit;
298 for (i = 0; i < rotg->PotAngle_nstep; i++)
300 if (fit->V[i] < pot_min)
303 fitangle = fit->degangle[i];
311 /* Reduce potential angle fit data for this group at this time step? */
312 static gmx_inline gmx_bool bPotAngle(t_rot *rot, t_rotgrp *rotg, gmx_int64_t step)
314 return ( (erotgFitPOT == rotg->eFittype) && (do_per_step(step, rot->nstsout) || do_per_step(step, rot->nstrout)) );
317 /* Reduce slab torqe data for this group at this time step? */
318 static gmx_inline gmx_bool bSlabTau(t_rot *rot, t_rotgrp *rotg, gmx_int64_t step)
320 return ( (ISFLEX(rotg)) && do_per_step(step, rot->nstsout) );
323 /* Output rotation energy, torques, etc. for each rotation group */
324 static void reduce_output(t_commrec *cr, t_rot *rot, real t, gmx_int64_t step)
326 int g, i, islab, nslabs = 0;
327 int count; /* MPI element counter */
329 gmx_enfrot_t er; /* Pointer to the enforced rotation buffer variables */
330 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
337 /* Fill the MPI buffer with stuff to reduce. If items are added for reduction
338 * here, the MPI buffer size has to be enlarged also in calc_mpi_bufsize() */
342 for (g = 0; g < rot->ngrp; g++)
345 erg = rotg->enfrotgrp;
346 nslabs = erg->slab_last - erg->slab_first + 1;
347 er->mpi_inbuf[count++] = erg->V;
348 er->mpi_inbuf[count++] = erg->torque_v;
349 er->mpi_inbuf[count++] = erg->angle_v;
350 er->mpi_inbuf[count++] = erg->weight_v; /* weights are not needed for flex types, but this is just a single value */
352 if (bPotAngle(rot, rotg, step))
354 for (i = 0; i < rotg->PotAngle_nstep; i++)
356 er->mpi_inbuf[count++] = erg->PotAngleFit->V[i];
359 if (bSlabTau(rot, rotg, step))
361 for (i = 0; i < nslabs; i++)
363 er->mpi_inbuf[count++] = erg->slab_torque_v[i];
367 if (count > er->mpi_bufsize)
369 gmx_fatal(FARGS, "%s MPI buffer overflow, please report this error.", RotStr);
373 MPI_Reduce(er->mpi_inbuf, er->mpi_outbuf, count, GMX_MPI_REAL, MPI_SUM, MASTERRANK(cr), cr->mpi_comm_mygroup);
376 /* Copy back the reduced data from the buffer on the master */
380 for (g = 0; g < rot->ngrp; g++)
383 erg = rotg->enfrotgrp;
384 nslabs = erg->slab_last - erg->slab_first + 1;
385 erg->V = er->mpi_outbuf[count++];
386 erg->torque_v = er->mpi_outbuf[count++];
387 erg->angle_v = er->mpi_outbuf[count++];
388 erg->weight_v = er->mpi_outbuf[count++];
390 if (bPotAngle(rot, rotg, step))
392 for (i = 0; i < rotg->PotAngle_nstep; i++)
394 erg->PotAngleFit->V[i] = er->mpi_outbuf[count++];
397 if (bSlabTau(rot, rotg, step))
399 for (i = 0; i < nslabs; i++)
401 erg->slab_torque_v[i] = er->mpi_outbuf[count++];
411 /* Angle and torque for each rotation group */
412 for (g = 0; g < rot->ngrp; g++)
415 bFlex = ISFLEX(rotg);
417 erg = rotg->enfrotgrp;
419 /* Output to main rotation output file: */
420 if (do_per_step(step, rot->nstrout) )
422 if (erotgFitPOT == rotg->eFittype)
424 fitangle = get_fitangle(rotg, erg);
430 fitangle = erg->angle_v; /* RMSD fit angle */
434 fitangle = (erg->angle_v/erg->weight_v)*180.0*M_1_PI;
437 fprintf(er->out_rot, "%12.4f", fitangle);
438 fprintf(er->out_rot, "%12.3e", erg->torque_v);
439 fprintf(er->out_rot, "%12.3e", erg->V);
442 if (do_per_step(step, rot->nstsout) )
444 /* Output to torque log file: */
447 fprintf(er->out_torque, "%12.3e%6d", t, g);
448 for (i = erg->slab_first; i <= erg->slab_last; i++)
450 islab = i - erg->slab_first; /* slab index */
451 /* Only output if enough weight is in slab */
452 if (erg->slab_weights[islab] > rotg->min_gaussian)
454 fprintf(er->out_torque, "%6d%12.3e", i, erg->slab_torque_v[islab]);
457 fprintf(er->out_torque, "\n");
460 /* Output to angles log file: */
461 if (erotgFitPOT == rotg->eFittype)
463 fprintf(er->out_angles, "%12.3e%6d%12.4f", t, g, erg->degangle);
464 /* Output energies at a set of angles around the reference angle */
465 for (i = 0; i < rotg->PotAngle_nstep; i++)
467 fprintf(er->out_angles, "%12.3e", erg->PotAngleFit->V[i]);
469 fprintf(er->out_angles, "\n");
473 if (do_per_step(step, rot->nstrout) )
475 fprintf(er->out_rot, "\n");
481 /* Add the forces from enforced rotation potential to the local forces.
482 * Should be called after the SR forces have been evaluated */
483 extern real add_rot_forces(t_rot *rot, rvec f[], t_commrec *cr, gmx_int64_t step, real t)
487 gmx_enfrot_t er; /* Pointer to the enforced rotation buffer variables */
488 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
489 real Vrot = 0.0; /* If more than one rotation group is present, Vrot
490 assembles the local parts from all groups */
495 /* Loop over enforced rotation groups (usually 1, though)
496 * Apply the forces from rotation potentials */
497 for (g = 0; g < rot->ngrp; g++)
500 erg = rotg->enfrotgrp;
501 Vrot += erg->V; /* add the local parts from the nodes */
502 for (l = 0; l < erg->nat_loc; l++)
504 /* Get the right index of the local force */
505 ii = erg->ind_loc[l];
507 rvec_inc(f[ii], erg->f_rot_loc[l]);
511 /* Reduce energy,torque, angles etc. to get the sum values (per rotation group)
512 * on the master and output these values to file. */
513 if ( (do_per_step(step, rot->nstrout) || do_per_step(step, rot->nstsout)) && er->bOut)
515 reduce_output(cr, rot, t, step);
518 /* When appending, er->bOut is FALSE the first time to avoid duplicate entries */
527 /* The Gaussian norm is chosen such that the sum of the gaussian functions
528 * over the slabs is approximately 1.0 everywhere */
529 #define GAUSS_NORM 0.569917543430618
532 /* Calculate the maximum beta that leads to a gaussian larger min_gaussian,
533 * also does some checks
535 static double calc_beta_max(real min_gaussian, real slab_dist)
541 /* Actually the next two checks are already made in grompp */
544 gmx_fatal(FARGS, "Slab distance of flexible rotation groups must be >=0 !");
546 if (min_gaussian <= 0)
548 gmx_fatal(FARGS, "Cutoff value for Gaussian must be > 0. (You requested %f)");
551 /* Define the sigma value */
552 sigma = 0.7*slab_dist;
554 /* Calculate the argument for the logarithm and check that the log() result is negative or 0 */
555 arg = min_gaussian/GAUSS_NORM;
558 gmx_fatal(FARGS, "min_gaussian of flexible rotation groups must be <%g", GAUSS_NORM);
561 return sqrt(-2.0*sigma*sigma*log(min_gaussian/GAUSS_NORM));
565 static gmx_inline real calc_beta(rvec curr_x, t_rotgrp *rotg, int n)
567 return iprod(curr_x, rotg->vec) - rotg->slab_dist * n;
571 static gmx_inline real gaussian_weight(rvec curr_x, t_rotgrp *rotg, int n)
573 const real norm = GAUSS_NORM;
577 /* Define the sigma value */
578 sigma = 0.7*rotg->slab_dist;
579 /* Calculate the Gaussian value of slab n for position curr_x */
580 return norm * exp( -0.5 * sqr( calc_beta(curr_x, rotg, n)/sigma ) );
584 /* Returns the weight in a single slab, also calculates the Gaussian- and mass-
585 * weighted sum of positions for that slab */
586 static real get_slab_weight(int j, t_rotgrp *rotg, rvec xc[], real mc[], rvec *x_weighted_sum)
588 rvec curr_x; /* The position of an atom */
589 rvec curr_x_weighted; /* The gaussian-weighted position */
590 real gaussian; /* A single gaussian weight */
591 real wgauss; /* gaussian times current mass */
592 real slabweight = 0.0; /* The sum of weights in the slab */
594 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
597 erg = rotg->enfrotgrp;
598 clear_rvec(*x_weighted_sum);
601 islab = j - erg->slab_first;
603 /* Loop over all atoms in the rotation group */
604 for (i = 0; i < rotg->nat; i++)
606 copy_rvec(xc[i], curr_x);
607 gaussian = gaussian_weight(curr_x, rotg, j);
608 wgauss = gaussian * mc[i];
609 svmul(wgauss, curr_x, curr_x_weighted);
610 rvec_add(*x_weighted_sum, curr_x_weighted, *x_weighted_sum);
611 slabweight += wgauss;
612 } /* END of loop over rotation group atoms */
618 static void get_slab_centers(
619 t_rotgrp *rotg, /* The rotation group information */
620 rvec *xc, /* The rotation group positions; will
621 typically be enfrotgrp->xc, but at first call
622 it is enfrotgrp->xc_ref */
623 real *mc, /* The masses of the rotation group atoms */
624 int g, /* The number of the rotation group */
625 real time, /* Used for output only */
626 FILE *out_slabs, /* For outputting center per slab information */
627 gmx_bool bOutStep, /* Is this an output step? */
628 gmx_bool bReference) /* If this routine is called from
629 init_rot_group we need to store
630 the reference slab centers */
633 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
636 erg = rotg->enfrotgrp;
638 /* Loop over slabs */
639 for (j = erg->slab_first; j <= erg->slab_last; j++)
641 islab = j - erg->slab_first;
642 erg->slab_weights[islab] = get_slab_weight(j, rotg, xc, mc, &erg->slab_center[islab]);
644 /* We can do the calculations ONLY if there is weight in the slab! */
645 if (erg->slab_weights[islab] > WEIGHT_MIN)
647 svmul(1.0/erg->slab_weights[islab], erg->slab_center[islab], erg->slab_center[islab]);
651 /* We need to check this here, since we divide through slab_weights
652 * in the flexible low-level routines! */
653 gmx_fatal(FARGS, "Not enough weight in slab %d. Slab center cannot be determined!", j);
656 /* At first time step: save the centers of the reference structure */
659 copy_rvec(erg->slab_center[islab], erg->slab_center_ref[islab]);
661 } /* END of loop over slabs */
663 /* Output on the master */
664 if ( (NULL != out_slabs) && bOutStep)
666 fprintf(out_slabs, "%12.3e%6d", time, g);
667 for (j = erg->slab_first; j <= erg->slab_last; j++)
669 islab = j - erg->slab_first;
670 fprintf(out_slabs, "%6d%12.3e%12.3e%12.3e",
671 j, erg->slab_center[islab][XX], erg->slab_center[islab][YY], erg->slab_center[islab][ZZ]);
673 fprintf(out_slabs, "\n");
678 static void calc_rotmat(
680 real degangle, /* Angle alpha of rotation at time t in degrees */
681 matrix rotmat) /* Rotation matrix */
683 real radangle; /* Rotation angle in radians */
684 real cosa; /* cosine alpha */
685 real sina; /* sine alpha */
686 real OMcosa; /* 1 - cos(alpha) */
687 real dumxy, dumxz, dumyz; /* save computations */
688 rvec rot_vec; /* Rotate around rot_vec ... */
691 radangle = degangle * M_PI/180.0;
692 copy_rvec(vec, rot_vec );
694 /* Precompute some variables: */
695 cosa = cos(radangle);
696 sina = sin(radangle);
698 dumxy = rot_vec[XX]*rot_vec[YY]*OMcosa;
699 dumxz = rot_vec[XX]*rot_vec[ZZ]*OMcosa;
700 dumyz = rot_vec[YY]*rot_vec[ZZ]*OMcosa;
702 /* Construct the rotation matrix for this rotation group: */
704 rotmat[XX][XX] = cosa + rot_vec[XX]*rot_vec[XX]*OMcosa;
705 rotmat[YY][XX] = dumxy + rot_vec[ZZ]*sina;
706 rotmat[ZZ][XX] = dumxz - rot_vec[YY]*sina;
708 rotmat[XX][YY] = dumxy - rot_vec[ZZ]*sina;
709 rotmat[YY][YY] = cosa + rot_vec[YY]*rot_vec[YY]*OMcosa;
710 rotmat[ZZ][YY] = dumyz + rot_vec[XX]*sina;
712 rotmat[XX][ZZ] = dumxz + rot_vec[YY]*sina;
713 rotmat[YY][ZZ] = dumyz - rot_vec[XX]*sina;
714 rotmat[ZZ][ZZ] = cosa + rot_vec[ZZ]*rot_vec[ZZ]*OMcosa;
719 for (iii = 0; iii < 3; iii++)
721 for (jjj = 0; jjj < 3; jjj++)
723 fprintf(stderr, " %10.8f ", rotmat[iii][jjj]);
725 fprintf(stderr, "\n");
731 /* Calculates torque on the rotation axis tau = position x force */
732 static gmx_inline real torque(
733 rvec rotvec, /* rotation vector; MUST be normalized! */
734 rvec force, /* force */
735 rvec x, /* position of atom on which the force acts */
736 rvec pivot) /* pivot point of rotation axis */
741 /* Subtract offset */
742 rvec_sub(x, pivot, vectmp);
744 /* position x force */
745 cprod(vectmp, force, tau);
747 /* Return the part of the torque which is parallel to the rotation vector */
748 return iprod(tau, rotvec);
752 /* Right-aligned output of value with standard width */
753 static void print_aligned(FILE *fp, char *str)
755 fprintf(fp, "%12s", str);
759 /* Right-aligned output of value with standard short width */
760 static void print_aligned_short(FILE *fp, char *str)
762 fprintf(fp, "%6s", str);
766 static FILE *open_output_file(const char *fn, int steps, const char what[])
771 fp = ffopen(fn, "w");
773 fprintf(fp, "# Output of %s is written in intervals of %d time step%s.\n#\n",
774 what, steps, steps > 1 ? "s" : "");
780 /* Open output file for slab center data. Call on master only */
781 static FILE *open_slab_out(const char *fn, t_rot *rot)
788 if (rot->enfrot->Flags & MD_APPENDFILES)
790 fp = gmx_fio_fopen(fn, "a");
794 fp = open_output_file(fn, rot->nstsout, "gaussian weighted slab centers");
796 for (g = 0; g < rot->ngrp; g++)
801 fprintf(fp, "# Rotation group %d (%s), slab distance %f nm, %s.\n",
802 g, erotg_names[rotg->eType], rotg->slab_dist,
803 rotg->bMassW ? "centers of mass" : "geometrical centers");
807 fprintf(fp, "# Reference centers are listed first (t=-1).\n");
808 fprintf(fp, "# The following columns have the syntax:\n");
810 print_aligned_short(fp, "t");
811 print_aligned_short(fp, "grp");
812 /* Print legend for the first two entries only ... */
813 for (i = 0; i < 2; i++)
815 print_aligned_short(fp, "slab");
816 print_aligned(fp, "X center");
817 print_aligned(fp, "Y center");
818 print_aligned(fp, "Z center");
820 fprintf(fp, " ...\n");
828 /* Adds 'buf' to 'str' */
829 static void add_to_string(char **str, char *buf)
834 len = strlen(*str) + strlen(buf) + 1;
840 static void add_to_string_aligned(char **str, char *buf)
842 char buf_aligned[STRLEN];
844 sprintf(buf_aligned, "%12s", buf);
845 add_to_string(str, buf_aligned);
849 /* Open output file and print some general information about the rotation groups.
850 * Call on master only */
851 static FILE *open_rot_out(const char *fn, t_rot *rot, const output_env_t oenv)
856 const char **setname;
857 char buf[50], buf2[75];
858 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
860 char *LegendStr = NULL;
863 if (rot->enfrot->Flags & MD_APPENDFILES)
865 fp = gmx_fio_fopen(fn, "a");
869 fp = xvgropen(fn, "Rotation angles and energy", "Time (ps)", "angles (degrees) and energies (kJ/mol)", oenv);
870 fprintf(fp, "# Output of enforced rotation data is written in intervals of %d time step%s.\n#\n", rot->nstrout, rot->nstrout > 1 ? "s" : "");
871 fprintf(fp, "# The scalar tau is the torque (kJ/mol) in the direction of the rotation vector v.\n");
872 fprintf(fp, "# To obtain the vectorial torque, multiply tau with the group's rot_vec.\n");
873 fprintf(fp, "# For flexible groups, tau(t,n) from all slabs n have been summed in a single value tau(t) here.\n");
874 fprintf(fp, "# The torques tau(t,n) are found in the rottorque.log (-rt) output file\n");
876 for (g = 0; g < rot->ngrp; g++)
879 erg = rotg->enfrotgrp;
880 bFlex = ISFLEX(rotg);
883 fprintf(fp, "# ROTATION GROUP %d, potential type '%s':\n", g, erotg_names[rotg->eType]);
884 fprintf(fp, "# rot_massw%d %s\n", g, yesno_names[rotg->bMassW]);
885 fprintf(fp, "# rot_vec%d %12.5e %12.5e %12.5e\n", g, rotg->vec[XX], rotg->vec[YY], rotg->vec[ZZ]);
886 fprintf(fp, "# rot_rate%d %12.5e degrees/ps\n", g, rotg->rate);
887 fprintf(fp, "# rot_k%d %12.5e kJ/(mol*nm^2)\n", g, rotg->k);
888 if (rotg->eType == erotgISO || rotg->eType == erotgPM || rotg->eType == erotgRM || rotg->eType == erotgRM2)
890 fprintf(fp, "# rot_pivot%d %12.5e %12.5e %12.5e nm\n", g, rotg->pivot[XX], rotg->pivot[YY], rotg->pivot[ZZ]);
895 fprintf(fp, "# rot_slab_distance%d %f nm\n", g, rotg->slab_dist);
896 fprintf(fp, "# rot_min_gaussian%d %12.5e\n", g, rotg->min_gaussian);
899 /* Output the centers of the rotation groups for the pivot-free potentials */
900 if ((rotg->eType == erotgISOPF) || (rotg->eType == erotgPMPF) || (rotg->eType == erotgRMPF) || (rotg->eType == erotgRM2PF
901 || (rotg->eType == erotgFLEXT) || (rotg->eType == erotgFLEX2T)) )
903 fprintf(fp, "# ref. grp. %d center %12.5e %12.5e %12.5e\n", g,
904 erg->xc_ref_center[XX], erg->xc_ref_center[YY], erg->xc_ref_center[ZZ]);
906 fprintf(fp, "# grp. %d init.center %12.5e %12.5e %12.5e\n", g,
907 erg->xc_center[XX], erg->xc_center[YY], erg->xc_center[ZZ]);
910 if ( (rotg->eType == erotgRM2) || (rotg->eType == erotgFLEX2) || (rotg->eType == erotgFLEX2T) )
912 fprintf(fp, "# rot_eps%d %12.5e nm^2\n", g, rotg->eps);
914 if (erotgFitPOT == rotg->eFittype)
917 fprintf(fp, "# theta_fit%d is determined by first evaluating the potential for %d angles around theta_ref%d.\n",
918 g, rotg->PotAngle_nstep, g);
919 fprintf(fp, "# The fit angle is the one with the smallest potential. It is given as the deviation\n");
920 fprintf(fp, "# from the reference angle, i.e. if theta_ref=X and theta_fit=Y, then the angle with\n");
921 fprintf(fp, "# minimal value of the potential is X+Y. Angular resolution is %g degrees.\n", rotg->PotAngle_step);
925 /* Print a nice legend */
928 sprintf(buf, "# %6s", "time");
929 add_to_string_aligned(&LegendStr, buf);
932 snew(setname, 4*rot->ngrp);
934 for (g = 0; g < rot->ngrp; g++)
937 sprintf(buf, "theta_ref%d", g);
938 add_to_string_aligned(&LegendStr, buf);
940 sprintf(buf2, "%s (degrees)", buf);
941 setname[nsets] = strdup(buf2);
944 for (g = 0; g < rot->ngrp; g++)
947 bFlex = ISFLEX(rotg);
949 /* For flexible axis rotation we use RMSD fitting to determine the
950 * actual angle of the rotation group */
951 if (bFlex || erotgFitPOT == rotg->eFittype)
953 sprintf(buf, "theta_fit%d", g);
957 sprintf(buf, "theta_av%d", g);
959 add_to_string_aligned(&LegendStr, buf);
960 sprintf(buf2, "%s (degrees)", buf);
961 setname[nsets] = strdup(buf2);
964 sprintf(buf, "tau%d", g);
965 add_to_string_aligned(&LegendStr, buf);
966 sprintf(buf2, "%s (kJ/mol)", buf);
967 setname[nsets] = strdup(buf2);
970 sprintf(buf, "energy%d", g);
971 add_to_string_aligned(&LegendStr, buf);
972 sprintf(buf2, "%s (kJ/mol)", buf);
973 setname[nsets] = strdup(buf2);
980 xvgr_legend(fp, nsets, setname, oenv);
984 fprintf(fp, "#\n# Legend for the following data columns:\n");
985 fprintf(fp, "%s\n", LegendStr);
995 /* Call on master only */
996 static FILE *open_angles_out(const char *fn, t_rot *rot)
1001 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
1005 if (rot->enfrot->Flags & MD_APPENDFILES)
1007 fp = gmx_fio_fopen(fn, "a");
1011 /* Open output file and write some information about it's structure: */
1012 fp = open_output_file(fn, rot->nstsout, "rotation group angles");
1013 fprintf(fp, "# All angles given in degrees, time in ps.\n");
1014 for (g = 0; g < rot->ngrp; g++)
1016 rotg = &rot->grp[g];
1017 erg = rotg->enfrotgrp;
1019 /* Output for this group happens only if potential type is flexible or
1020 * if fit type is potential! */
1021 if (ISFLEX(rotg) || (erotgFitPOT == rotg->eFittype) )
1025 sprintf(buf, " slab distance %f nm, ", rotg->slab_dist);
1032 fprintf(fp, "#\n# ROTATION GROUP %d '%s',%s fit type '%s'.\n",
1033 g, erotg_names[rotg->eType], buf, erotg_fitnames[rotg->eFittype]);
1035 /* Special type of fitting using the potential minimum. This is
1036 * done for the whole group only, not for the individual slabs. */
1037 if (erotgFitPOT == rotg->eFittype)
1039 fprintf(fp, "# To obtain theta_fit%d, the potential is evaluated for %d angles around theta_ref%d\n", g, rotg->PotAngle_nstep, g);
1040 fprintf(fp, "# The fit angle in the rotation standard outfile is the one with minimal energy E(theta_fit) [kJ/mol].\n");
1044 fprintf(fp, "# Legend for the group %d data columns:\n", g);
1046 print_aligned_short(fp, "time");
1047 print_aligned_short(fp, "grp");
1048 print_aligned(fp, "theta_ref");
1050 if (erotgFitPOT == rotg->eFittype)
1052 /* Output the set of angles around the reference angle */
1053 for (i = 0; i < rotg->PotAngle_nstep; i++)
1055 sprintf(buf, "E(%g)", erg->PotAngleFit->degangle[i]);
1056 print_aligned(fp, buf);
1061 /* Output fit angle for each slab */
1062 print_aligned_short(fp, "slab");
1063 print_aligned_short(fp, "atoms");
1064 print_aligned(fp, "theta_fit");
1065 print_aligned_short(fp, "slab");
1066 print_aligned_short(fp, "atoms");
1067 print_aligned(fp, "theta_fit");
1068 fprintf(fp, " ...");
1080 /* Open torque output file and write some information about it's structure.
1081 * Call on master only */
1082 static FILE *open_torque_out(const char *fn, t_rot *rot)
1089 if (rot->enfrot->Flags & MD_APPENDFILES)
1091 fp = gmx_fio_fopen(fn, "a");
1095 fp = open_output_file(fn, rot->nstsout, "torques");
1097 for (g = 0; g < rot->ngrp; g++)
1099 rotg = &rot->grp[g];
1102 fprintf(fp, "# Rotation group %d (%s), slab distance %f nm.\n", g, erotg_names[rotg->eType], rotg->slab_dist);
1103 fprintf(fp, "# The scalar tau is the torque (kJ/mol) in the direction of the rotation vector.\n");
1104 fprintf(fp, "# To obtain the vectorial torque, multiply tau with\n");
1105 fprintf(fp, "# rot_vec%d %10.3e %10.3e %10.3e\n", g, rotg->vec[XX], rotg->vec[YY], rotg->vec[ZZ]);
1109 fprintf(fp, "# Legend for the following data columns: (tau=torque for that slab):\n");
1111 print_aligned_short(fp, "t");
1112 print_aligned_short(fp, "grp");
1113 print_aligned_short(fp, "slab");
1114 print_aligned(fp, "tau");
1115 print_aligned_short(fp, "slab");
1116 print_aligned(fp, "tau");
1117 fprintf(fp, " ...\n");
1125 static void swap_val(double* vec, int i, int j)
1127 double tmp = vec[j];
1135 static void swap_col(double **mat, int i, int j)
1137 double tmp[3] = {mat[0][j], mat[1][j], mat[2][j]};
1140 mat[0][j] = mat[0][i];
1141 mat[1][j] = mat[1][i];
1142 mat[2][j] = mat[2][i];
1150 /* Eigenvectors are stored in columns of eigen_vec */
1151 static void diagonalize_symmetric(
1159 jacobi(matrix, 3, eigenval, eigen_vec, &n_rot);
1161 /* sort in ascending order */
1162 if (eigenval[0] > eigenval[1])
1164 swap_val(eigenval, 0, 1);
1165 swap_col(eigen_vec, 0, 1);
1167 if (eigenval[1] > eigenval[2])
1169 swap_val(eigenval, 1, 2);
1170 swap_col(eigen_vec, 1, 2);
1172 if (eigenval[0] > eigenval[1])
1174 swap_val(eigenval, 0, 1);
1175 swap_col(eigen_vec, 0, 1);
1180 static void align_with_z(
1181 rvec* s, /* Structure to align */
1186 rvec zet = {0.0, 0.0, 1.0};
1187 rvec rot_axis = {0.0, 0.0, 0.0};
1188 rvec *rotated_str = NULL;
1194 snew(rotated_str, natoms);
1196 /* Normalize the axis */
1197 ooanorm = 1.0/norm(axis);
1198 svmul(ooanorm, axis, axis);
1200 /* Calculate the angle for the fitting procedure */
1201 cprod(axis, zet, rot_axis);
1202 angle = acos(axis[2]);
1208 /* Calculate the rotation matrix */
1209 calc_rotmat(rot_axis, angle*180.0/M_PI, rotmat);
1211 /* Apply the rotation matrix to s */
1212 for (i = 0; i < natoms; i++)
1214 for (j = 0; j < 3; j++)
1216 for (k = 0; k < 3; k++)
1218 rotated_str[i][j] += rotmat[j][k]*s[i][k];
1223 /* Rewrite the rotated structure to s */
1224 for (i = 0; i < natoms; i++)
1226 for (j = 0; j < 3; j++)
1228 s[i][j] = rotated_str[i][j];
1236 static void calc_correl_matrix(rvec* Xstr, rvec* Ystr, double** Rmat, int natoms)
1241 for (i = 0; i < 3; i++)
1243 for (j = 0; j < 3; j++)
1249 for (i = 0; i < 3; i++)
1251 for (j = 0; j < 3; j++)
1253 for (k = 0; k < natoms; k++)
1255 Rmat[i][j] += Ystr[k][i] * Xstr[k][j];
1262 static void weigh_coords(rvec* str, real* weight, int natoms)
1267 for (i = 0; i < natoms; i++)
1269 for (j = 0; j < 3; j++)
1271 str[i][j] *= sqrt(weight[i]);
1277 static real opt_angle_analytic(
1287 rvec *ref_s_1 = NULL;
1288 rvec *act_s_1 = NULL;
1290 double **Rmat, **RtR, **eigvec;
1292 double V[3][3], WS[3][3];
1293 double rot_matrix[3][3];
1297 /* Do not change the original coordinates */
1298 snew(ref_s_1, natoms);
1299 snew(act_s_1, natoms);
1300 for (i = 0; i < natoms; i++)
1302 copy_rvec(ref_s[i], ref_s_1[i]);
1303 copy_rvec(act_s[i], act_s_1[i]);
1306 /* Translate the structures to the origin */
1307 shift[XX] = -ref_com[XX];
1308 shift[YY] = -ref_com[YY];
1309 shift[ZZ] = -ref_com[ZZ];
1310 translate_x(ref_s_1, natoms, shift);
1312 shift[XX] = -act_com[XX];
1313 shift[YY] = -act_com[YY];
1314 shift[ZZ] = -act_com[ZZ];
1315 translate_x(act_s_1, natoms, shift);
1317 /* Align rotation axis with z */
1318 align_with_z(ref_s_1, natoms, axis);
1319 align_with_z(act_s_1, natoms, axis);
1321 /* Correlation matrix */
1322 Rmat = allocate_square_matrix(3);
1324 for (i = 0; i < natoms; i++)
1326 ref_s_1[i][2] = 0.0;
1327 act_s_1[i][2] = 0.0;
1330 /* Weight positions with sqrt(weight) */
1333 weigh_coords(ref_s_1, weight, natoms);
1334 weigh_coords(act_s_1, weight, natoms);
1337 /* Calculate correlation matrices R=YXt (X=ref_s; Y=act_s) */
1338 calc_correl_matrix(ref_s_1, act_s_1, Rmat, natoms);
1341 RtR = allocate_square_matrix(3);
1342 for (i = 0; i < 3; i++)
1344 for (j = 0; j < 3; j++)
1346 for (k = 0; k < 3; k++)
1348 RtR[i][j] += Rmat[k][i] * Rmat[k][j];
1352 /* Diagonalize RtR */
1354 for (i = 0; i < 3; i++)
1359 diagonalize_symmetric(RtR, eigvec, eigval);
1360 swap_col(eigvec, 0, 1);
1361 swap_col(eigvec, 1, 2);
1362 swap_val(eigval, 0, 1);
1363 swap_val(eigval, 1, 2);
1366 for (i = 0; i < 3; i++)
1368 for (j = 0; j < 3; j++)
1375 for (i = 0; i < 2; i++)
1377 for (j = 0; j < 2; j++)
1379 WS[i][j] = eigvec[i][j] / sqrt(eigval[j]);
1383 for (i = 0; i < 3; i++)
1385 for (j = 0; j < 3; j++)
1387 for (k = 0; k < 3; k++)
1389 V[i][j] += Rmat[i][k]*WS[k][j];
1393 free_square_matrix(Rmat, 3);
1395 /* Calculate optimal rotation matrix */
1396 for (i = 0; i < 3; i++)
1398 for (j = 0; j < 3; j++)
1400 rot_matrix[i][j] = 0.0;
1404 for (i = 0; i < 3; i++)
1406 for (j = 0; j < 3; j++)
1408 for (k = 0; k < 3; k++)
1410 rot_matrix[i][j] += eigvec[i][k]*V[j][k];
1414 rot_matrix[2][2] = 1.0;
1416 /* In some cases abs(rot_matrix[0][0]) can be slighly larger
1417 * than unity due to numerical inacurracies. To be able to calculate
1418 * the acos function, we put these values back in range. */
1419 if (rot_matrix[0][0] > 1.0)
1421 rot_matrix[0][0] = 1.0;
1423 else if (rot_matrix[0][0] < -1.0)
1425 rot_matrix[0][0] = -1.0;
1428 /* Determine the optimal rotation angle: */
1429 opt_angle = (-1.0)*acos(rot_matrix[0][0])*180.0/M_PI;
1430 if (rot_matrix[0][1] < 0.0)
1432 opt_angle = (-1.0)*opt_angle;
1435 /* Give back some memory */
1436 free_square_matrix(RtR, 3);
1439 for (i = 0; i < 3; i++)
1445 return (real) opt_angle;
1449 /* Determine angle of the group by RMSD fit to the reference */
1450 /* Not parallelized, call this routine only on the master */
1451 static real flex_fit_angle(t_rotgrp *rotg)
1454 rvec *fitcoords = NULL;
1455 rvec center; /* Center of positions passed to the fit routine */
1456 real fitangle; /* Angle of the rotation group derived by fitting */
1459 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
1462 erg = rotg->enfrotgrp;
1464 /* Get the center of the rotation group.
1465 * Note, again, erg->xc has been sorted in do_flexible */
1466 get_center(erg->xc, erg->mc_sorted, rotg->nat, center);
1468 /* === Determine the optimal fit angle for the rotation group === */
1469 if (rotg->eFittype == erotgFitNORM)
1471 /* Normalize every position to it's reference length */
1472 for (i = 0; i < rotg->nat; i++)
1474 /* Put the center of the positions into the origin */
1475 rvec_sub(erg->xc[i], center, coord);
1476 /* Determine the scaling factor for the length: */
1477 scal = erg->xc_ref_length[erg->xc_sortind[i]] / norm(coord);
1478 /* Get position, multiply with the scaling factor and save */
1479 svmul(scal, coord, erg->xc_norm[i]);
1481 fitcoords = erg->xc_norm;
1485 fitcoords = erg->xc;
1487 /* From the point of view of the current positions, the reference has rotated
1488 * backwards. Since we output the angle relative to the fixed reference,
1489 * we need the minus sign. */
1490 fitangle = -opt_angle_analytic(erg->xc_ref_sorted, fitcoords, erg->mc_sorted,
1491 rotg->nat, erg->xc_ref_center, center, rotg->vec);
1497 /* Determine actual angle of each slab by RMSD fit to the reference */
1498 /* Not parallelized, call this routine only on the master */
1499 static void flex_fit_angle_perslab(
1506 int i, l, n, islab, ind;
1508 rvec act_center; /* Center of actual positions that are passed to the fit routine */
1509 rvec ref_center; /* Same for the reference positions */
1510 real fitangle; /* Angle of a slab derived from an RMSD fit to
1511 * the reference structure at t=0 */
1513 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
1514 real OOm_av; /* 1/average_mass of a rotation group atom */
1515 real m_rel; /* Relative mass of a rotation group atom */
1518 erg = rotg->enfrotgrp;
1520 /* Average mass of a rotation group atom: */
1521 OOm_av = erg->invmass*rotg->nat;
1523 /**********************************/
1524 /* First collect the data we need */
1525 /**********************************/
1527 /* Collect the data for the individual slabs */
1528 for (n = erg->slab_first; n <= erg->slab_last; n++)
1530 islab = n - erg->slab_first; /* slab index */
1531 sd = &(rotg->enfrotgrp->slab_data[islab]);
1532 sd->nat = erg->lastatom[islab]-erg->firstatom[islab]+1;
1535 /* Loop over the relevant atoms in the slab */
1536 for (l = erg->firstatom[islab]; l <= erg->lastatom[islab]; l++)
1538 /* Current position of this atom: x[ii][XX/YY/ZZ] */
1539 copy_rvec(erg->xc[l], curr_x);
1541 /* The (unrotated) reference position of this atom is copied to ref_x.
1542 * Beware, the xc coords have been sorted in do_flexible */
1543 copy_rvec(erg->xc_ref_sorted[l], ref_x);
1545 /* Save data for doing angular RMSD fit later */
1546 /* Save the current atom position */
1547 copy_rvec(curr_x, sd->x[ind]);
1548 /* Save the corresponding reference position */
1549 copy_rvec(ref_x, sd->ref[ind]);
1551 /* Maybe also mass-weighting was requested. If yes, additionally
1552 * multiply the weights with the relative mass of the atom. If not,
1553 * multiply with unity. */
1554 m_rel = erg->mc_sorted[l]*OOm_av;
1556 /* Save the weight for this atom in this slab */
1557 sd->weight[ind] = gaussian_weight(curr_x, rotg, n) * m_rel;
1559 /* Next atom in this slab */
1564 /******************************/
1565 /* Now do the fit calculation */
1566 /******************************/
1568 fprintf(fp, "%12.3e%6d%12.3f", t, g, degangle);
1570 /* === Now do RMSD fitting for each slab === */
1571 /* We require at least SLAB_MIN_ATOMS in a slab, such that the fit makes sense. */
1572 #define SLAB_MIN_ATOMS 4
1574 for (n = erg->slab_first; n <= erg->slab_last; n++)
1576 islab = n - erg->slab_first; /* slab index */
1577 sd = &(rotg->enfrotgrp->slab_data[islab]);
1578 if (sd->nat >= SLAB_MIN_ATOMS)
1580 /* Get the center of the slabs reference and current positions */
1581 get_center(sd->ref, sd->weight, sd->nat, ref_center);
1582 get_center(sd->x, sd->weight, sd->nat, act_center);
1583 if (rotg->eFittype == erotgFitNORM)
1585 /* Normalize every position to it's reference length
1586 * prior to performing the fit */
1587 for (i = 0; i < sd->nat; i++) /* Center */
1589 rvec_dec(sd->ref[i], ref_center);
1590 rvec_dec(sd->x[i], act_center);
1591 /* Normalize x_i such that it gets the same length as ref_i */
1592 svmul( norm(sd->ref[i])/norm(sd->x[i]), sd->x[i], sd->x[i] );
1594 /* We already subtracted the centers */
1595 clear_rvec(ref_center);
1596 clear_rvec(act_center);
1598 fitangle = -opt_angle_analytic(sd->ref, sd->x, sd->weight, sd->nat,
1599 ref_center, act_center, rotg->vec);
1600 fprintf(fp, "%6d%6d%12.3f", n, sd->nat, fitangle);
1605 #undef SLAB_MIN_ATOMS
1609 /* Shift x with is */
1610 static gmx_inline void shift_single_coord(matrix box, rvec x, const ivec is)
1621 x[XX] += tx*box[XX][XX]+ty*box[YY][XX]+tz*box[ZZ][XX];
1622 x[YY] += ty*box[YY][YY]+tz*box[ZZ][YY];
1623 x[ZZ] += tz*box[ZZ][ZZ];
1627 x[XX] += tx*box[XX][XX];
1628 x[YY] += ty*box[YY][YY];
1629 x[ZZ] += tz*box[ZZ][ZZ];
1634 /* Determine the 'home' slab of this atom which is the
1635 * slab with the highest Gaussian weight of all */
1636 #define round(a) (int)(a+0.5)
1637 static gmx_inline int get_homeslab(
1638 rvec curr_x, /* The position for which the home slab shall be determined */
1639 rvec rotvec, /* The rotation vector */
1640 real slabdist) /* The slab distance */
1645 /* The distance of the atom to the coordinate center (where the
1646 * slab with index 0) is */
1647 dist = iprod(rotvec, curr_x);
1649 return round(dist / slabdist);
1653 /* For a local atom determine the relevant slabs, i.e. slabs in
1654 * which the gaussian is larger than min_gaussian
1656 static int get_single_atom_gaussians(
1663 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
1666 erg = rotg->enfrotgrp;
1668 /* Determine the 'home' slab of this atom: */
1669 homeslab = get_homeslab(curr_x, rotg->vec, rotg->slab_dist);
1671 /* First determine the weight in the atoms home slab: */
1672 g = gaussian_weight(curr_x, rotg, homeslab);
1674 erg->gn_atom[count] = g;
1675 erg->gn_slabind[count] = homeslab;
1679 /* Determine the max slab */
1681 while (g > rotg->min_gaussian)
1684 g = gaussian_weight(curr_x, rotg, slab);
1685 erg->gn_slabind[count] = slab;
1686 erg->gn_atom[count] = g;
1691 /* Determine the min slab */
1696 g = gaussian_weight(curr_x, rotg, slab);
1697 erg->gn_slabind[count] = slab;
1698 erg->gn_atom[count] = g;
1701 while (g > rotg->min_gaussian);
1708 static void flex2_precalc_inner_sum(t_rotgrp *rotg)
1711 rvec xi; /* positions in the i-sum */
1712 rvec xcn, ycn; /* the current and the reference slab centers */
1715 rvec rin; /* Helper variables */
1718 real OOpsii, OOpsiistar;
1719 real sin_rin; /* s_ii.r_ii */
1720 rvec s_in, tmpvec, tmpvec2;
1721 real mi, wi; /* Mass-weighting of the positions */
1723 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
1726 erg = rotg->enfrotgrp;
1727 N_M = rotg->nat * erg->invmass;
1729 /* Loop over all slabs that contain something */
1730 for (n = erg->slab_first; n <= erg->slab_last; n++)
1732 islab = n - erg->slab_first; /* slab index */
1734 /* The current center of this slab is saved in xcn: */
1735 copy_rvec(erg->slab_center[islab], xcn);
1736 /* ... and the reference center in ycn: */
1737 copy_rvec(erg->slab_center_ref[islab+erg->slab_buffer], ycn);
1739 /*** D. Calculate the whole inner sum used for second and third sum */
1740 /* For slab n, we need to loop over all atoms i again. Since we sorted
1741 * the atoms with respect to the rotation vector, we know that it is sufficient
1742 * to calculate from firstatom to lastatom only. All other contributions will
1744 clear_rvec(innersumvec);
1745 for (i = erg->firstatom[islab]; i <= erg->lastatom[islab]; i++)
1747 /* Coordinate xi of this atom */
1748 copy_rvec(erg->xc[i], xi);
1751 gaussian_xi = gaussian_weight(xi, rotg, n);
1752 mi = erg->mc_sorted[i]; /* need the sorted mass here */
1756 copy_rvec(erg->xc_ref_sorted[i], yi0); /* Reference position yi0 */
1757 rvec_sub(yi0, ycn, tmpvec2); /* tmpvec2 = yi0 - ycn */
1758 mvmul(erg->rotmat, tmpvec2, rin); /* rin = Omega.(yi0 - ycn) */
1760 /* Calculate psi_i* and sin */
1761 rvec_sub(xi, xcn, tmpvec2); /* tmpvec2 = xi - xcn */
1762 cprod(rotg->vec, tmpvec2, tmpvec); /* tmpvec = v x (xi - xcn) */
1763 OOpsiistar = norm2(tmpvec)+rotg->eps; /* OOpsii* = 1/psii* = |v x (xi-xcn)|^2 + eps */
1764 OOpsii = norm(tmpvec); /* OOpsii = 1 / psii = |v x (xi - xcn)| */
1766 /* * v x (xi - xcn) */
1767 unitv(tmpvec, s_in); /* sin = ---------------- */
1768 /* |v x (xi - xcn)| */
1770 sin_rin = iprod(s_in, rin); /* sin_rin = sin . rin */
1772 /* Now the whole sum */
1773 fac = OOpsii/OOpsiistar;
1774 svmul(fac, rin, tmpvec);
1775 fac2 = fac*fac*OOpsii;
1776 svmul(fac2*sin_rin, s_in, tmpvec2);
1777 rvec_dec(tmpvec, tmpvec2);
1779 svmul(wi*gaussian_xi*sin_rin, tmpvec, tmpvec2);
1781 rvec_inc(innersumvec, tmpvec2);
1782 } /* now we have the inner sum, used both for sum2 and sum3 */
1784 /* Save it to be used in do_flex2_lowlevel */
1785 copy_rvec(innersumvec, erg->slab_innersumvec[islab]);
1786 } /* END of loop over slabs */
1790 static void flex_precalc_inner_sum(t_rotgrp *rotg)
1793 rvec xi; /* position */
1794 rvec xcn, ycn; /* the current and the reference slab centers */
1795 rvec qin, rin; /* q_i^n and r_i^n */
1798 rvec innersumvec; /* Inner part of sum_n2 */
1799 real gaussian_xi; /* Gaussian weight gn(xi) */
1800 real mi, wi; /* Mass-weighting of the positions */
1803 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
1806 erg = rotg->enfrotgrp;
1807 N_M = rotg->nat * erg->invmass;
1809 /* Loop over all slabs that contain something */
1810 for (n = erg->slab_first; n <= erg->slab_last; n++)
1812 islab = n - erg->slab_first; /* slab index */
1814 /* The current center of this slab is saved in xcn: */
1815 copy_rvec(erg->slab_center[islab], xcn);
1816 /* ... and the reference center in ycn: */
1817 copy_rvec(erg->slab_center_ref[islab+erg->slab_buffer], ycn);
1819 /* For slab n, we need to loop over all atoms i again. Since we sorted
1820 * the atoms with respect to the rotation vector, we know that it is sufficient
1821 * to calculate from firstatom to lastatom only. All other contributions will
1823 clear_rvec(innersumvec);
1824 for (i = erg->firstatom[islab]; i <= erg->lastatom[islab]; i++)
1826 /* Coordinate xi of this atom */
1827 copy_rvec(erg->xc[i], xi);
1830 gaussian_xi = gaussian_weight(xi, rotg, n);
1831 mi = erg->mc_sorted[i]; /* need the sorted mass here */
1834 /* Calculate rin and qin */
1835 rvec_sub(erg->xc_ref_sorted[i], ycn, tmpvec); /* tmpvec = yi0-ycn */
1836 mvmul(erg->rotmat, tmpvec, rin); /* rin = Omega.(yi0 - ycn) */
1837 cprod(rotg->vec, rin, tmpvec); /* tmpvec = v x Omega*(yi0-ycn) */
1839 /* * v x Omega*(yi0-ycn) */
1840 unitv(tmpvec, qin); /* qin = --------------------- */
1841 /* |v x Omega*(yi0-ycn)| */
1844 rvec_sub(xi, xcn, tmpvec); /* tmpvec = xi-xcn */
1845 bin = iprod(qin, tmpvec); /* bin = qin*(xi-xcn) */
1847 svmul(wi*gaussian_xi*bin, qin, tmpvec);
1849 /* Add this contribution to the inner sum: */
1850 rvec_add(innersumvec, tmpvec, innersumvec);
1851 } /* now we have the inner sum vector S^n for this slab */
1852 /* Save it to be used in do_flex_lowlevel */
1853 copy_rvec(innersumvec, erg->slab_innersumvec[islab]);
1858 static real do_flex2_lowlevel(
1860 real sigma, /* The Gaussian width sigma */
1862 gmx_bool bOutstepRot,
1863 gmx_bool bOutstepSlab,
1866 int count, ic, ii, j, m, n, islab, iigrp, ifit;
1867 rvec xj; /* position in the i-sum */
1868 rvec yj0; /* the reference position in the j-sum */
1869 rvec xcn, ycn; /* the current and the reference slab centers */
1870 real V; /* This node's part of the rotation pot. energy */
1871 real gaussian_xj; /* Gaussian weight */
1874 real numerator, fit_numerator;
1875 rvec rjn, fit_rjn; /* Helper variables */
1878 real OOpsij, OOpsijstar;
1879 real OOsigma2; /* 1/(sigma^2) */
1882 rvec sjn, tmpvec, tmpvec2, yj0_ycn;
1883 rvec sum1vec_part, sum1vec, sum2vec_part, sum2vec, sum3vec, sum4vec, innersumvec;
1885 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
1886 real mj, wj; /* Mass-weighting of the positions */
1888 real Wjn; /* g_n(x_j) m_j / Mjn */
1889 gmx_bool bCalcPotFit;
1891 /* To calculate the torque per slab */
1892 rvec slab_force; /* Single force from slab n on one atom */
1893 rvec slab_sum1vec_part;
1894 real slab_sum3part, slab_sum4part;
1895 rvec slab_sum1vec, slab_sum2vec, slab_sum3vec, slab_sum4vec;
1898 erg = rotg->enfrotgrp;
1900 /* Pre-calculate the inner sums, so that we do not have to calculate
1901 * them again for every atom */
1902 flex2_precalc_inner_sum(rotg);
1904 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT == rotg->eFittype);
1906 /********************************************************/
1907 /* Main loop over all local atoms of the rotation group */
1908 /********************************************************/
1909 N_M = rotg->nat * erg->invmass;
1911 OOsigma2 = 1.0 / (sigma*sigma);
1912 for (j = 0; j < erg->nat_loc; j++)
1914 /* Local index of a rotation group atom */
1915 ii = erg->ind_loc[j];
1916 /* Position of this atom in the collective array */
1917 iigrp = erg->xc_ref_ind[j];
1918 /* Mass-weighting */
1919 mj = erg->mc[iigrp]; /* need the unsorted mass here */
1922 /* Current position of this atom: x[ii][XX/YY/ZZ]
1923 * Note that erg->xc_center contains the center of mass in case the flex2-t
1924 * potential was chosen. For the flex2 potential erg->xc_center must be
1926 rvec_sub(x[ii], erg->xc_center, xj);
1928 /* Shift this atom such that it is near its reference */
1929 shift_single_coord(box, xj, erg->xc_shifts[iigrp]);
1931 /* Determine the slabs to loop over, i.e. the ones with contributions
1932 * larger than min_gaussian */
1933 count = get_single_atom_gaussians(xj, rotg);
1935 clear_rvec(sum1vec_part);
1936 clear_rvec(sum2vec_part);
1939 /* Loop over the relevant slabs for this atom */
1940 for (ic = 0; ic < count; ic++)
1942 n = erg->gn_slabind[ic];
1944 /* Get the precomputed Gaussian value of curr_slab for curr_x */
1945 gaussian_xj = erg->gn_atom[ic];
1947 islab = n - erg->slab_first; /* slab index */
1949 /* The (unrotated) reference position of this atom is copied to yj0: */
1950 copy_rvec(rotg->x_ref[iigrp], yj0);
1952 beta = calc_beta(xj, rotg, n);
1954 /* The current center of this slab is saved in xcn: */
1955 copy_rvec(erg->slab_center[islab], xcn);
1956 /* ... and the reference center in ycn: */
1957 copy_rvec(erg->slab_center_ref[islab+erg->slab_buffer], ycn);
1959 rvec_sub(yj0, ycn, yj0_ycn); /* yj0_ycn = yj0 - ycn */
1962 mvmul(erg->rotmat, yj0_ycn, rjn); /* rjn = Omega.(yj0 - ycn) */
1964 /* Subtract the slab center from xj */
1965 rvec_sub(xj, xcn, tmpvec2); /* tmpvec2 = xj - xcn */
1967 /* In rare cases, when an atom position coincides with a slab center
1968 * (tmpvec2 == 0) we cannot compute the vector product for sjn.
1969 * However, since the atom is located directly on the pivot, this
1970 * slab's contribution to the force on that atom will be zero
1971 * anyway. Therefore, we directly move on to the next slab. */
1972 if (0 == norm(tmpvec2) )
1978 cprod(rotg->vec, tmpvec2, tmpvec); /* tmpvec = v x (xj - xcn) */
1980 OOpsijstar = norm2(tmpvec)+rotg->eps; /* OOpsij* = 1/psij* = |v x (xj-xcn)|^2 + eps */
1982 numerator = sqr(iprod(tmpvec, rjn));
1984 /*********************************/
1985 /* Add to the rotation potential */
1986 /*********************************/
1987 V += 0.5*rotg->k*wj*gaussian_xj*numerator/OOpsijstar;
1989 /* If requested, also calculate the potential for a set of angles
1990 * near the current reference angle */
1993 for (ifit = 0; ifit < rotg->PotAngle_nstep; ifit++)
1995 mvmul(erg->PotAngleFit->rotmat[ifit], yj0_ycn, fit_rjn);
1996 fit_numerator = sqr(iprod(tmpvec, fit_rjn));
1997 erg->PotAngleFit->V[ifit] += 0.5*rotg->k*wj*gaussian_xj*fit_numerator/OOpsijstar;
2001 /*************************************/
2002 /* Now calculate the force on atom j */
2003 /*************************************/
2005 OOpsij = norm(tmpvec); /* OOpsij = 1 / psij = |v x (xj - xcn)| */
2007 /* * v x (xj - xcn) */
2008 unitv(tmpvec, sjn); /* sjn = ---------------- */
2009 /* |v x (xj - xcn)| */
2011 sjn_rjn = iprod(sjn, rjn); /* sjn_rjn = sjn . rjn */
2014 /*** A. Calculate the first of the four sum terms: ****************/
2015 fac = OOpsij/OOpsijstar;
2016 svmul(fac, rjn, tmpvec);
2017 fac2 = fac*fac*OOpsij;
2018 svmul(fac2*sjn_rjn, sjn, tmpvec2);
2019 rvec_dec(tmpvec, tmpvec2);
2020 fac2 = wj*gaussian_xj; /* also needed for sum4 */
2021 svmul(fac2*sjn_rjn, tmpvec, slab_sum1vec_part);
2022 /********************/
2023 /*** Add to sum1: ***/
2024 /********************/
2025 rvec_inc(sum1vec_part, slab_sum1vec_part); /* sum1 still needs to vector multiplied with v */
2027 /*** B. Calculate the forth of the four sum terms: ****************/
2028 betasigpsi = beta*OOsigma2*OOpsij; /* this is also needed for sum3 */
2029 /********************/
2030 /*** Add to sum4: ***/
2031 /********************/
2032 slab_sum4part = fac2*betasigpsi*fac*sjn_rjn*sjn_rjn; /* Note that fac is still valid from above */
2033 sum4 += slab_sum4part;
2035 /*** C. Calculate Wjn for second and third sum */
2036 /* Note that we can safely divide by slab_weights since we check in
2037 * get_slab_centers that it is non-zero. */
2038 Wjn = gaussian_xj*mj/erg->slab_weights[islab];
2040 /* We already have precalculated the inner sum for slab n */
2041 copy_rvec(erg->slab_innersumvec[islab], innersumvec);
2043 /* Weigh the inner sum vector with Wjn */
2044 svmul(Wjn, innersumvec, innersumvec);
2046 /*** E. Calculate the second of the four sum terms: */
2047 /********************/
2048 /*** Add to sum2: ***/
2049 /********************/
2050 rvec_inc(sum2vec_part, innersumvec); /* sum2 still needs to be vector crossproduct'ed with v */
2052 /*** F. Calculate the third of the four sum terms: */
2053 slab_sum3part = betasigpsi * iprod(sjn, innersumvec);
2054 sum3 += slab_sum3part; /* still needs to be multiplied with v */
2056 /*** G. Calculate the torque on the local slab's axis: */
2060 cprod(slab_sum1vec_part, rotg->vec, slab_sum1vec);
2062 cprod(innersumvec, rotg->vec, slab_sum2vec);
2064 svmul(slab_sum3part, rotg->vec, slab_sum3vec);
2066 svmul(slab_sum4part, rotg->vec, slab_sum4vec);
2068 /* The force on atom ii from slab n only: */
2069 for (m = 0; m < DIM; m++)
2071 slab_force[m] = rotg->k * (-slab_sum1vec[m] + slab_sum2vec[m] - slab_sum3vec[m] + 0.5*slab_sum4vec[m]);
2074 erg->slab_torque_v[islab] += torque(rotg->vec, slab_force, xj, xcn);
2076 } /* END of loop over slabs */
2078 /* Construct the four individual parts of the vector sum: */
2079 cprod(sum1vec_part, rotg->vec, sum1vec); /* sum1vec = { } x v */
2080 cprod(sum2vec_part, rotg->vec, sum2vec); /* sum2vec = { } x v */
2081 svmul(sum3, rotg->vec, sum3vec); /* sum3vec = { } . v */
2082 svmul(sum4, rotg->vec, sum4vec); /* sum4vec = { } . v */
2084 /* Store the additional force so that it can be added to the force
2085 * array after the normal forces have been evaluated */
2086 for (m = 0; m < DIM; m++)
2088 erg->f_rot_loc[j][m] = rotg->k * (-sum1vec[m] + sum2vec[m] - sum3vec[m] + 0.5*sum4vec[m]);
2092 fprintf(stderr, "sum1: %15.8f %15.8f %15.8f\n", -rotg->k*sum1vec[XX], -rotg->k*sum1vec[YY], -rotg->k*sum1vec[ZZ]);
2093 fprintf(stderr, "sum2: %15.8f %15.8f %15.8f\n", rotg->k*sum2vec[XX], rotg->k*sum2vec[YY], rotg->k*sum2vec[ZZ]);
2094 fprintf(stderr, "sum3: %15.8f %15.8f %15.8f\n", -rotg->k*sum3vec[XX], -rotg->k*sum3vec[YY], -rotg->k*sum3vec[ZZ]);
2095 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]);
2100 } /* END of loop over local atoms */
2106 static real do_flex_lowlevel(
2108 real sigma, /* The Gaussian width sigma */
2110 gmx_bool bOutstepRot,
2111 gmx_bool bOutstepSlab,
2114 int count, ic, ifit, ii, j, m, n, islab, iigrp;
2115 rvec xj, yj0; /* current and reference position */
2116 rvec xcn, ycn; /* the current and the reference slab centers */
2117 rvec yj0_ycn; /* yj0 - ycn */
2118 rvec xj_xcn; /* xj - xcn */
2119 rvec qjn, fit_qjn; /* q_i^n */
2120 rvec sum_n1, sum_n2; /* Two contributions to the rotation force */
2121 rvec innersumvec; /* Inner part of sum_n2 */
2123 rvec force_n; /* Single force from slab n on one atom */
2124 rvec force_n1, force_n2; /* First and second part of force_n */
2125 rvec tmpvec, tmpvec2, tmp_f; /* Helper variables */
2126 real V; /* The rotation potential energy */
2127 real OOsigma2; /* 1/(sigma^2) */
2128 real beta; /* beta_n(xj) */
2129 real bjn, fit_bjn; /* b_j^n */
2130 real gaussian_xj; /* Gaussian weight gn(xj) */
2131 real betan_xj_sigma2;
2132 real mj, wj; /* Mass-weighting of the positions */
2134 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
2135 gmx_bool bCalcPotFit;
2138 erg = rotg->enfrotgrp;
2140 /* Pre-calculate the inner sums, so that we do not have to calculate
2141 * them again for every atom */
2142 flex_precalc_inner_sum(rotg);
2144 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT == rotg->eFittype);
2146 /********************************************************/
2147 /* Main loop over all local atoms of the rotation group */
2148 /********************************************************/
2149 OOsigma2 = 1.0/(sigma*sigma);
2150 N_M = rotg->nat * erg->invmass;
2152 for (j = 0; j < erg->nat_loc; j++)
2154 /* Local index of a rotation group atom */
2155 ii = erg->ind_loc[j];
2156 /* Position of this atom in the collective array */
2157 iigrp = erg->xc_ref_ind[j];
2158 /* Mass-weighting */
2159 mj = erg->mc[iigrp]; /* need the unsorted mass here */
2162 /* Current position of this atom: x[ii][XX/YY/ZZ]
2163 * Note that erg->xc_center contains the center of mass in case the flex-t
2164 * potential was chosen. For the flex potential erg->xc_center must be
2166 rvec_sub(x[ii], erg->xc_center, xj);
2168 /* Shift this atom such that it is near its reference */
2169 shift_single_coord(box, xj, erg->xc_shifts[iigrp]);
2171 /* Determine the slabs to loop over, i.e. the ones with contributions
2172 * larger than min_gaussian */
2173 count = get_single_atom_gaussians(xj, rotg);
2178 /* Loop over the relevant slabs for this atom */
2179 for (ic = 0; ic < count; ic++)
2181 n = erg->gn_slabind[ic];
2183 /* Get the precomputed Gaussian for xj in slab n */
2184 gaussian_xj = erg->gn_atom[ic];
2186 islab = n - erg->slab_first; /* slab index */
2188 /* The (unrotated) reference position of this atom is saved in yj0: */
2189 copy_rvec(rotg->x_ref[iigrp], yj0);
2191 beta = calc_beta(xj, rotg, n);
2193 /* The current center of this slab is saved in xcn: */
2194 copy_rvec(erg->slab_center[islab], xcn);
2195 /* ... and the reference center in ycn: */
2196 copy_rvec(erg->slab_center_ref[islab+erg->slab_buffer], ycn);
2198 rvec_sub(yj0, ycn, yj0_ycn); /* yj0_ycn = yj0 - ycn */
2201 mvmul(erg->rotmat, yj0_ycn, tmpvec2); /* tmpvec2= Omega.(yj0-ycn) */
2203 /* Subtract the slab center from xj */
2204 rvec_sub(xj, xcn, xj_xcn); /* xj_xcn = xj - xcn */
2206 /* In rare cases, when an atom position coincides with a slab center
2207 * (xj_xcn == 0) we cannot compute the vector product for qjn.
2208 * However, since the atom is located directly on the pivot, this
2209 * slab's contribution to the force on that atom will be zero
2210 * anyway. Therefore, we directly move on to the next slab. */
2211 if (0 == norm(xj_xcn) )
2217 cprod(rotg->vec, tmpvec2, tmpvec); /* tmpvec= v x Omega.(yj0-ycn) */
2219 /* * v x Omega.(yj0-ycn) */
2220 unitv(tmpvec, qjn); /* qjn = --------------------- */
2221 /* |v x Omega.(yj0-ycn)| */
2223 bjn = iprod(qjn, xj_xcn); /* bjn = qjn * (xj - xcn) */
2225 /*********************************/
2226 /* Add to the rotation potential */
2227 /*********************************/
2228 V += 0.5*rotg->k*wj*gaussian_xj*sqr(bjn);
2230 /* If requested, also calculate the potential for a set of angles
2231 * near the current reference angle */
2234 for (ifit = 0; ifit < rotg->PotAngle_nstep; ifit++)
2236 /* As above calculate Omega.(yj0-ycn), now for the other angles */
2237 mvmul(erg->PotAngleFit->rotmat[ifit], yj0_ycn, tmpvec2); /* tmpvec2= Omega.(yj0-ycn) */
2238 /* As above calculate qjn */
2239 cprod(rotg->vec, tmpvec2, tmpvec); /* tmpvec= v x Omega.(yj0-ycn) */
2240 /* * v x Omega.(yj0-ycn) */
2241 unitv(tmpvec, fit_qjn); /* fit_qjn = --------------------- */
2242 /* |v x Omega.(yj0-ycn)| */
2243 fit_bjn = iprod(fit_qjn, xj_xcn); /* fit_bjn = fit_qjn * (xj - xcn) */
2244 /* Add to the rotation potential for this angle */
2245 erg->PotAngleFit->V[ifit] += 0.5*rotg->k*wj*gaussian_xj*sqr(fit_bjn);
2249 /****************************************************************/
2250 /* sum_n1 will typically be the main contribution to the force: */
2251 /****************************************************************/
2252 betan_xj_sigma2 = beta*OOsigma2; /* beta_n(xj)/sigma^2 */
2254 /* The next lines calculate
2255 * qjn - (bjn*beta(xj)/(2sigma^2))v */
2256 svmul(bjn*0.5*betan_xj_sigma2, rotg->vec, tmpvec2);
2257 rvec_sub(qjn, tmpvec2, tmpvec);
2259 /* Multiply with gn(xj)*bjn: */
2260 svmul(gaussian_xj*bjn, tmpvec, tmpvec2);
2263 rvec_inc(sum_n1, tmpvec2);
2265 /* We already have precalculated the Sn term for slab n */
2266 copy_rvec(erg->slab_innersumvec[islab], s_n);
2268 svmul(betan_xj_sigma2*iprod(s_n, xj_xcn), rotg->vec, tmpvec); /* tmpvec = ---------- s_n (xj-xcn) */
2271 rvec_sub(s_n, tmpvec, innersumvec);
2273 /* We can safely divide by slab_weights since we check in get_slab_centers
2274 * that it is non-zero. */
2275 svmul(gaussian_xj/erg->slab_weights[islab], innersumvec, innersumvec);
2277 rvec_add(sum_n2, innersumvec, sum_n2);
2279 /* Calculate the torque: */
2282 /* The force on atom ii from slab n only: */
2283 svmul(-rotg->k*wj, tmpvec2, force_n1); /* part 1 */
2284 svmul( rotg->k*mj, innersumvec, force_n2); /* part 2 */
2285 rvec_add(force_n1, force_n2, force_n);
2286 erg->slab_torque_v[islab] += torque(rotg->vec, force_n, xj, xcn);
2288 } /* END of loop over slabs */
2290 /* Put both contributions together: */
2291 svmul(wj, sum_n1, sum_n1);
2292 svmul(mj, sum_n2, sum_n2);
2293 rvec_sub(sum_n2, sum_n1, tmp_f); /* F = -grad V */
2295 /* Store the additional force so that it can be added to the force
2296 * array after the normal forces have been evaluated */
2297 for (m = 0; m < DIM; m++)
2299 erg->f_rot_loc[j][m] = rotg->k*tmp_f[m];
2304 } /* END of loop over local atoms */
2310 static void print_coordinates(t_rotgrp *rotg, rvec x[], matrix box, int step)
2314 static char buf[STRLEN];
2315 static gmx_bool bFirst = 1;
2320 sprintf(buf, "coords%d.txt", cr->nodeid);
2321 fp = fopen(buf, "w");
2325 fprintf(fp, "\nStep %d\n", step);
2326 fprintf(fp, "box: %f %f %f %f %f %f %f %f %f\n",
2327 box[XX][XX], box[XX][YY], box[XX][ZZ],
2328 box[YY][XX], box[YY][YY], box[YY][ZZ],
2329 box[ZZ][XX], box[ZZ][ZZ], box[ZZ][ZZ]);
2330 for (i = 0; i < rotg->nat; i++)
2332 fprintf(fp, "%4d %f %f %f\n", i,
2333 erg->xc[i][XX], erg->xc[i][YY], erg->xc[i][ZZ]);
2341 static int projection_compare(const void *a, const void *b)
2343 sort_along_vec_t *xca, *xcb;
2346 xca = (sort_along_vec_t *)a;
2347 xcb = (sort_along_vec_t *)b;
2349 if (xca->xcproj < xcb->xcproj)
2353 else if (xca->xcproj > xcb->xcproj)
2364 static void sort_collective_coordinates(
2365 t_rotgrp *rotg, /* Rotation group */
2366 sort_along_vec_t *data) /* Buffer for sorting the positions */
2369 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
2372 erg = rotg->enfrotgrp;
2374 /* The projection of the position vector on the rotation vector is
2375 * the relevant value for sorting. Fill the 'data' structure */
2376 for (i = 0; i < rotg->nat; i++)
2378 data[i].xcproj = iprod(erg->xc[i], rotg->vec); /* sort criterium */
2379 data[i].m = erg->mc[i];
2381 copy_rvec(erg->xc[i], data[i].x );
2382 copy_rvec(rotg->x_ref[i], data[i].x_ref);
2384 /* Sort the 'data' structure */
2385 gmx_qsort(data, rotg->nat, sizeof(sort_along_vec_t), projection_compare);
2387 /* Copy back the sorted values */
2388 for (i = 0; i < rotg->nat; i++)
2390 copy_rvec(data[i].x, erg->xc[i] );
2391 copy_rvec(data[i].x_ref, erg->xc_ref_sorted[i]);
2392 erg->mc_sorted[i] = data[i].m;
2393 erg->xc_sortind[i] = data[i].ind;
2398 /* For each slab, get the first and the last index of the sorted atom
2400 static void get_firstlast_atom_per_slab(t_rotgrp *rotg)
2404 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
2407 erg = rotg->enfrotgrp;
2409 /* Find the first atom that needs to enter the calculation for each slab */
2410 n = erg->slab_first; /* slab */
2411 i = 0; /* start with the first atom */
2414 /* Find the first atom that significantly contributes to this slab */
2415 do /* move forward in position until a large enough beta is found */
2417 beta = calc_beta(erg->xc[i], rotg, n);
2420 while ((beta < -erg->max_beta) && (i < rotg->nat));
2422 islab = n - erg->slab_first; /* slab index */
2423 erg->firstatom[islab] = i;
2424 /* Proceed to the next slab */
2427 while (n <= erg->slab_last);
2429 /* Find the last atom for each slab */
2430 n = erg->slab_last; /* start with last slab */
2431 i = rotg->nat-1; /* start with the last atom */
2434 do /* move backward in position until a large enough beta is found */
2436 beta = calc_beta(erg->xc[i], rotg, n);
2439 while ((beta > erg->max_beta) && (i > -1));
2441 islab = n - erg->slab_first; /* slab index */
2442 erg->lastatom[islab] = i;
2443 /* Proceed to the next slab */
2446 while (n >= erg->slab_first);
2450 /* Determine the very first and very last slab that needs to be considered
2451 * For the first slab that needs to be considered, we have to find the smallest
2454 * x_first * v - n*Delta_x <= beta_max
2456 * slab index n, slab distance Delta_x, rotation vector v. For the last slab we
2457 * have to find the largest n that obeys
2459 * x_last * v - n*Delta_x >= -beta_max
2462 static gmx_inline int get_first_slab(
2463 t_rotgrp *rotg, /* The rotation group (inputrec data) */
2464 real max_beta, /* The max_beta value, instead of min_gaussian */
2465 rvec firstatom) /* First atom after sorting along the rotation vector v */
2467 /* Find the first slab for the first atom */
2468 return ceil((iprod(firstatom, rotg->vec) - max_beta)/rotg->slab_dist);
2472 static gmx_inline int get_last_slab(
2473 t_rotgrp *rotg, /* The rotation group (inputrec data) */
2474 real max_beta, /* The max_beta value, instead of min_gaussian */
2475 rvec lastatom) /* Last atom along v */
2477 /* Find the last slab for the last atom */
2478 return floor((iprod(lastatom, rotg->vec) + max_beta)/rotg->slab_dist);
2482 static void get_firstlast_slab_check(
2483 t_rotgrp *rotg, /* The rotation group (inputrec data) */
2484 t_gmx_enfrotgrp *erg, /* The rotation group (data only accessible in this file) */
2485 rvec firstatom, /* First atom after sorting along the rotation vector v */
2486 rvec lastatom) /* Last atom along v */
2488 erg->slab_first = get_first_slab(rotg, erg->max_beta, firstatom);
2489 erg->slab_last = get_last_slab(rotg, erg->max_beta, lastatom);
2491 /* Calculate the slab buffer size, which changes when slab_first changes */
2492 erg->slab_buffer = erg->slab_first - erg->slab_first_ref;
2494 /* Check whether we have reference data to compare against */
2495 if (erg->slab_first < erg->slab_first_ref)
2497 gmx_fatal(FARGS, "%s No reference data for first slab (n=%d), unable to proceed.",
2498 RotStr, erg->slab_first);
2501 /* Check whether we have reference data to compare against */
2502 if (erg->slab_last > erg->slab_last_ref)
2504 gmx_fatal(FARGS, "%s No reference data for last slab (n=%d), unable to proceed.",
2505 RotStr, erg->slab_last);
2510 /* Enforced rotation with a flexible axis */
2511 static void do_flexible(
2513 gmx_enfrot_t enfrot, /* Other rotation data */
2514 t_rotgrp *rotg, /* The rotation group */
2515 int g, /* Group number */
2516 rvec x[], /* The local positions */
2518 double t, /* Time in picoseconds */
2519 gmx_bool bOutstepRot, /* Output to main rotation output file */
2520 gmx_bool bOutstepSlab) /* Output per-slab data */
2523 real sigma; /* The Gaussian width sigma */
2524 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
2527 erg = rotg->enfrotgrp;
2529 /* Define the sigma value */
2530 sigma = 0.7*rotg->slab_dist;
2532 /* Sort the collective coordinates erg->xc along the rotation vector. This is
2533 * an optimization for the inner loop. */
2534 sort_collective_coordinates(rotg, enfrot->data);
2536 /* Determine the first relevant slab for the first atom and the last
2537 * relevant slab for the last atom */
2538 get_firstlast_slab_check(rotg, erg, erg->xc[0], erg->xc[rotg->nat-1]);
2540 /* Determine for each slab depending on the min_gaussian cutoff criterium,
2541 * a first and a last atom index inbetween stuff needs to be calculated */
2542 get_firstlast_atom_per_slab(rotg);
2544 /* Determine the gaussian-weighted center of positions for all slabs */
2545 get_slab_centers(rotg, erg->xc, erg->mc_sorted, g, t, enfrot->out_slabs, bOutstepSlab, FALSE);
2547 /* Clear the torque per slab from last time step: */
2548 nslabs = erg->slab_last - erg->slab_first + 1;
2549 for (l = 0; l < nslabs; l++)
2551 erg->slab_torque_v[l] = 0.0;
2554 /* Call the rotational forces kernel */
2555 if (rotg->eType == erotgFLEX || rotg->eType == erotgFLEXT)
2557 erg->V = do_flex_lowlevel(rotg, sigma, x, bOutstepRot, bOutstepSlab, box);
2559 else if (rotg->eType == erotgFLEX2 || rotg->eType == erotgFLEX2T)
2561 erg->V = do_flex2_lowlevel(rotg, sigma, x, bOutstepRot, bOutstepSlab, box);
2565 gmx_fatal(FARGS, "Unknown flexible rotation type");
2568 /* Determine angle by RMSD fit to the reference - Let's hope this */
2569 /* only happens once in a while, since this is not parallelized! */
2570 if (bMaster && (erotgFitPOT != rotg->eFittype) )
2574 /* Fit angle of the whole rotation group */
2575 erg->angle_v = flex_fit_angle(rotg);
2579 /* Fit angle of each slab */
2580 flex_fit_angle_perslab(g, rotg, t, erg->degangle, enfrot->out_angles);
2584 /* Lump together the torques from all slabs: */
2585 erg->torque_v = 0.0;
2586 for (l = 0; l < nslabs; l++)
2588 erg->torque_v += erg->slab_torque_v[l];
2593 /* Calculate the angle between reference and actual rotation group atom,
2594 * both projected into a plane perpendicular to the rotation vector: */
2595 static void angle(t_rotgrp *rotg,
2599 real *weight) /* atoms near the rotation axis should count less than atoms far away */
2601 rvec xp, xrp; /* current and reference positions projected on a plane perpendicular to pg->vec */
2605 /* Project x_ref and x into a plane through the origin perpendicular to rot_vec: */
2606 /* Project x_ref: xrp = x_ref - (vec * x_ref) * vec */
2607 svmul(iprod(rotg->vec, x_ref), rotg->vec, dum);
2608 rvec_sub(x_ref, dum, xrp);
2609 /* Project x_act: */
2610 svmul(iprod(rotg->vec, x_act), rotg->vec, dum);
2611 rvec_sub(x_act, dum, xp);
2613 /* Retrieve information about which vector precedes. gmx_angle always
2614 * returns a positive angle. */
2615 cprod(xp, xrp, dum); /* if reference precedes, this is pointing into the same direction as vec */
2617 if (iprod(rotg->vec, dum) >= 0)
2619 *alpha = -gmx_angle(xrp, xp);
2623 *alpha = +gmx_angle(xrp, xp);
2626 /* Also return the weight */
2631 /* Project first vector onto a plane perpendicular to the second vector
2633 * Note that v must be of unit length.
2635 static gmx_inline void project_onto_plane(rvec dr, const rvec v)
2640 svmul(iprod(dr, v), v, tmp); /* tmp = (dr.v)v */
2641 rvec_dec(dr, tmp); /* dr = dr - (dr.v)v */
2645 /* Fixed rotation: The rotation reference group rotates around the v axis. */
2646 /* The atoms of the actual rotation group are attached with imaginary */
2647 /* springs to the reference atoms. */
2648 static void do_fixed(
2649 t_rotgrp *rotg, /* The rotation group */
2650 gmx_bool bOutstepRot, /* Output to main rotation output file */
2651 gmx_bool bOutstepSlab) /* Output per-slab data */
2655 rvec tmp_f; /* Force */
2656 real alpha; /* a single angle between an actual and a reference position */
2657 real weight; /* single weight for a single angle */
2658 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
2659 rvec xi_xc; /* xi - xc */
2660 gmx_bool bCalcPotFit;
2663 /* for mass weighting: */
2664 real wi; /* Mass-weighting of the positions */
2666 real k_wi; /* k times wi */
2671 erg = rotg->enfrotgrp;
2672 bProject = (rotg->eType == erotgPM) || (rotg->eType == erotgPMPF);
2673 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT == rotg->eFittype);
2675 N_M = rotg->nat * erg->invmass;
2677 /* Each process calculates the forces on its local atoms */
2678 for (j = 0; j < erg->nat_loc; j++)
2680 /* Calculate (x_i-x_c) resp. (x_i-u) */
2681 rvec_sub(erg->x_loc_pbc[j], erg->xc_center, xi_xc);
2683 /* Calculate Omega*(y_i-y_c)-(x_i-x_c) */
2684 rvec_sub(erg->xr_loc[j], xi_xc, dr);
2688 project_onto_plane(dr, rotg->vec);
2691 /* Mass-weighting */
2692 wi = N_M*erg->m_loc[j];
2694 /* Store the additional force so that it can be added to the force
2695 * array after the normal forces have been evaluated */
2697 for (m = 0; m < DIM; m++)
2699 tmp_f[m] = k_wi*dr[m];
2700 erg->f_rot_loc[j][m] = tmp_f[m];
2701 erg->V += 0.5*k_wi*sqr(dr[m]);
2704 /* If requested, also calculate the potential for a set of angles
2705 * near the current reference angle */
2708 for (ifit = 0; ifit < rotg->PotAngle_nstep; ifit++)
2710 /* Index of this rotation group atom with respect to the whole rotation group */
2711 jj = erg->xc_ref_ind[j];
2713 /* Rotate with the alternative angle. Like rotate_local_reference(),
2714 * just for a single local atom */
2715 mvmul(erg->PotAngleFit->rotmat[ifit], rotg->x_ref[jj], fit_xr_loc); /* fit_xr_loc = Omega*(y_i-y_c) */
2717 /* Calculate Omega*(y_i-y_c)-(x_i-x_c) */
2718 rvec_sub(fit_xr_loc, xi_xc, dr);
2722 project_onto_plane(dr, rotg->vec);
2725 /* Add to the rotation potential for this angle: */
2726 erg->PotAngleFit->V[ifit] += 0.5*k_wi*norm2(dr);
2732 /* Add to the torque of this rotation group */
2733 erg->torque_v += torque(rotg->vec, tmp_f, erg->x_loc_pbc[j], erg->xc_center);
2735 /* Calculate the angle between reference and actual rotation group atom. */
2736 angle(rotg, xi_xc, erg->xr_loc[j], &alpha, &weight); /* angle in rad, weighted */
2737 erg->angle_v += alpha * weight;
2738 erg->weight_v += weight;
2740 /* If you want enforced rotation to contribute to the virial,
2741 * activate the following lines:
2744 Add the rotation contribution to the virial
2745 for(j=0; j<DIM; j++)
2747 vir[j][m] += 0.5*f[ii][j]*dr[m];
2753 } /* end of loop over local rotation group atoms */
2757 /* Calculate the radial motion potential and forces */
2758 static void do_radial_motion(
2759 t_rotgrp *rotg, /* The rotation group */
2760 gmx_bool bOutstepRot, /* Output to main rotation output file */
2761 gmx_bool bOutstepSlab) /* Output per-slab data */
2764 rvec tmp_f; /* Force */
2765 real alpha; /* a single angle between an actual and a reference position */
2766 real weight; /* single weight for a single angle */
2767 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
2768 rvec xj_u; /* xj - u */
2769 rvec tmpvec, fit_tmpvec;
2770 real fac, fac2, sum = 0.0;
2772 gmx_bool bCalcPotFit;
2774 /* For mass weighting: */
2775 real wj; /* Mass-weighting of the positions */
2779 erg = rotg->enfrotgrp;
2780 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT == rotg->eFittype);
2782 N_M = rotg->nat * erg->invmass;
2784 /* Each process calculates the forces on its local atoms */
2785 for (j = 0; j < erg->nat_loc; j++)
2787 /* Calculate (xj-u) */
2788 rvec_sub(erg->x_loc_pbc[j], erg->xc_center, xj_u); /* xj_u = xj-u */
2790 /* Calculate Omega.(yj0-u) */
2791 cprod(rotg->vec, erg->xr_loc[j], tmpvec); /* tmpvec = v x Omega.(yj0-u) */
2793 /* * v x Omega.(yj0-u) */
2794 unitv(tmpvec, pj); /* pj = --------------------- */
2795 /* | v x Omega.(yj0-u) | */
2797 fac = iprod(pj, xj_u); /* fac = pj.(xj-u) */
2800 /* Mass-weighting */
2801 wj = N_M*erg->m_loc[j];
2803 /* Store the additional force so that it can be added to the force
2804 * array after the normal forces have been evaluated */
2805 svmul(-rotg->k*wj*fac, pj, tmp_f);
2806 copy_rvec(tmp_f, erg->f_rot_loc[j]);
2809 /* If requested, also calculate the potential for a set of angles
2810 * near the current reference angle */
2813 for (ifit = 0; ifit < rotg->PotAngle_nstep; ifit++)
2815 /* Index of this rotation group atom with respect to the whole rotation group */
2816 jj = erg->xc_ref_ind[j];
2818 /* Rotate with the alternative angle. Like rotate_local_reference(),
2819 * just for a single local atom */
2820 mvmul(erg->PotAngleFit->rotmat[ifit], rotg->x_ref[jj], fit_tmpvec); /* fit_tmpvec = Omega*(yj0-u) */
2822 /* Calculate Omega.(yj0-u) */
2823 cprod(rotg->vec, fit_tmpvec, tmpvec); /* tmpvec = v x Omega.(yj0-u) */
2824 /* * v x Omega.(yj0-u) */
2825 unitv(tmpvec, pj); /* pj = --------------------- */
2826 /* | v x Omega.(yj0-u) | */
2828 fac = iprod(pj, xj_u); /* fac = pj.(xj-u) */
2831 /* Add to the rotation potential for this angle: */
2832 erg->PotAngleFit->V[ifit] += 0.5*rotg->k*wj*fac2;
2838 /* Add to the torque of this rotation group */
2839 erg->torque_v += torque(rotg->vec, tmp_f, erg->x_loc_pbc[j], erg->xc_center);
2841 /* Calculate the angle between reference and actual rotation group atom. */
2842 angle(rotg, xj_u, erg->xr_loc[j], &alpha, &weight); /* angle in rad, weighted */
2843 erg->angle_v += alpha * weight;
2844 erg->weight_v += weight;
2849 } /* end of loop over local rotation group atoms */
2850 erg->V = 0.5*rotg->k*sum;
2854 /* Calculate the radial motion pivot-free potential and forces */
2855 static void do_radial_motion_pf(
2856 t_rotgrp *rotg, /* The rotation group */
2857 rvec x[], /* The positions */
2858 matrix box, /* The simulation box */
2859 gmx_bool bOutstepRot, /* Output to main rotation output file */
2860 gmx_bool bOutstepSlab) /* Output per-slab data */
2862 int i, ii, iigrp, ifit, j;
2863 rvec xj; /* Current position */
2864 rvec xj_xc; /* xj - xc */
2865 rvec yj0_yc0; /* yj0 - yc0 */
2866 rvec tmp_f; /* Force */
2867 real alpha; /* a single angle between an actual and a reference position */
2868 real weight; /* single weight for a single angle */
2869 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
2870 rvec tmpvec, tmpvec2;
2871 rvec innersumvec; /* Precalculation of the inner sum */
2873 real fac, fac2, V = 0.0;
2875 gmx_bool bCalcPotFit;
2877 /* For mass weighting: */
2878 real mj, wi, wj; /* Mass-weighting of the positions */
2882 erg = rotg->enfrotgrp;
2883 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT == rotg->eFittype);
2885 N_M = rotg->nat * erg->invmass;
2887 /* Get the current center of the rotation group: */
2888 get_center(erg->xc, erg->mc, rotg->nat, erg->xc_center);
2890 /* Precalculate Sum_i [ wi qi.(xi-xc) qi ] which is needed for every single j */
2891 clear_rvec(innersumvec);
2892 for (i = 0; i < rotg->nat; i++)
2894 /* Mass-weighting */
2895 wi = N_M*erg->mc[i];
2897 /* Calculate qi. Note that xc_ref_center has already been subtracted from
2898 * x_ref in init_rot_group.*/
2899 mvmul(erg->rotmat, rotg->x_ref[i], tmpvec); /* tmpvec = Omega.(yi0-yc0) */
2901 cprod(rotg->vec, tmpvec, tmpvec2); /* tmpvec2 = v x Omega.(yi0-yc0) */
2903 /* * v x Omega.(yi0-yc0) */
2904 unitv(tmpvec2, qi); /* qi = ----------------------- */
2905 /* | v x Omega.(yi0-yc0) | */
2907 rvec_sub(erg->xc[i], erg->xc_center, tmpvec); /* tmpvec = xi-xc */
2909 svmul(wi*iprod(qi, tmpvec), qi, tmpvec2);
2911 rvec_inc(innersumvec, tmpvec2);
2913 svmul(rotg->k*erg->invmass, innersumvec, innersumveckM);
2915 /* Each process calculates the forces on its local atoms */
2916 for (j = 0; j < erg->nat_loc; j++)
2918 /* Local index of a rotation group atom */
2919 ii = erg->ind_loc[j];
2920 /* Position of this atom in the collective array */
2921 iigrp = erg->xc_ref_ind[j];
2922 /* Mass-weighting */
2923 mj = erg->mc[iigrp]; /* need the unsorted mass here */
2926 /* Current position of this atom: x[ii][XX/YY/ZZ] */
2927 copy_rvec(x[ii], xj);
2929 /* Shift this atom such that it is near its reference */
2930 shift_single_coord(box, xj, erg->xc_shifts[iigrp]);
2932 /* The (unrotated) reference position is yj0. yc0 has already
2933 * been subtracted in init_rot_group */
2934 copy_rvec(rotg->x_ref[iigrp], yj0_yc0); /* yj0_yc0 = yj0 - yc0 */
2936 /* Calculate Omega.(yj0-yc0) */
2937 mvmul(erg->rotmat, yj0_yc0, tmpvec2); /* tmpvec2 = Omega.(yj0 - yc0) */
2939 cprod(rotg->vec, tmpvec2, tmpvec); /* tmpvec = v x Omega.(yj0-yc0) */
2941 /* * v x Omega.(yj0-yc0) */
2942 unitv(tmpvec, qj); /* qj = ----------------------- */
2943 /* | v x Omega.(yj0-yc0) | */
2945 /* Calculate (xj-xc) */
2946 rvec_sub(xj, erg->xc_center, xj_xc); /* xj_xc = xj-xc */
2948 fac = iprod(qj, xj_xc); /* fac = qj.(xj-xc) */
2951 /* Store the additional force so that it can be added to the force
2952 * array after the normal forces have been evaluated */
2953 svmul(-rotg->k*wj*fac, qj, tmp_f); /* part 1 of force */
2954 svmul(mj, innersumveckM, tmpvec); /* part 2 of force */
2955 rvec_inc(tmp_f, tmpvec);
2956 copy_rvec(tmp_f, erg->f_rot_loc[j]);
2959 /* If requested, also calculate the potential for a set of angles
2960 * near the current reference angle */
2963 for (ifit = 0; ifit < rotg->PotAngle_nstep; ifit++)
2965 /* Rotate with the alternative angle. Like rotate_local_reference(),
2966 * just for a single local atom */
2967 mvmul(erg->PotAngleFit->rotmat[ifit], yj0_yc0, tmpvec2); /* tmpvec2 = Omega*(yj0-yc0) */
2969 /* Calculate Omega.(yj0-u) */
2970 cprod(rotg->vec, tmpvec2, tmpvec); /* tmpvec = v x Omega.(yj0-yc0) */
2971 /* * v x Omega.(yj0-yc0) */
2972 unitv(tmpvec, qj); /* qj = ----------------------- */
2973 /* | v x Omega.(yj0-yc0) | */
2975 fac = iprod(qj, xj_xc); /* fac = qj.(xj-xc) */
2978 /* Add to the rotation potential for this angle: */
2979 erg->PotAngleFit->V[ifit] += 0.5*rotg->k*wj*fac2;
2985 /* Add to the torque of this rotation group */
2986 erg->torque_v += torque(rotg->vec, tmp_f, xj, erg->xc_center);
2988 /* Calculate the angle between reference and actual rotation group atom. */
2989 angle(rotg, xj_xc, yj0_yc0, &alpha, &weight); /* angle in rad, weighted */
2990 erg->angle_v += alpha * weight;
2991 erg->weight_v += weight;
2996 } /* end of loop over local rotation group atoms */
2997 erg->V = 0.5*rotg->k*V;
3001 /* Precalculate the inner sum for the radial motion 2 forces */
3002 static void radial_motion2_precalc_inner_sum(t_rotgrp *rotg, rvec innersumvec)
3005 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
3006 rvec xi_xc; /* xj - xc */
3007 rvec tmpvec, tmpvec2;
3011 rvec v_xi_xc; /* v x (xj - u) */
3012 real psii, psiistar;
3013 real wi; /* Mass-weighting of the positions */
3017 erg = rotg->enfrotgrp;
3018 N_M = rotg->nat * erg->invmass;
3020 /* Loop over the collective set of positions */
3022 for (i = 0; i < rotg->nat; i++)
3024 /* Mass-weighting */
3025 wi = N_M*erg->mc[i];
3027 rvec_sub(erg->xc[i], erg->xc_center, xi_xc); /* xi_xc = xi-xc */
3029 /* Calculate ri. Note that xc_ref_center has already been subtracted from
3030 * x_ref in init_rot_group.*/
3031 mvmul(erg->rotmat, rotg->x_ref[i], ri); /* ri = Omega.(yi0-yc0) */
3033 cprod(rotg->vec, xi_xc, v_xi_xc); /* v_xi_xc = v x (xi-u) */
3035 fac = norm2(v_xi_xc);
3037 psiistar = 1.0/(fac + rotg->eps); /* psiistar = --------------------- */
3038 /* |v x (xi-xc)|^2 + eps */
3040 psii = gmx_invsqrt(fac); /* 1 */
3041 /* psii = ------------- */
3044 svmul(psii, v_xi_xc, si); /* si = psii * (v x (xi-xc) ) */
3046 fac = iprod(v_xi_xc, ri); /* fac = (v x (xi-xc)).ri */
3049 siri = iprod(si, ri); /* siri = si.ri */
3051 svmul(psiistar/psii, ri, tmpvec);
3052 svmul(psiistar*psiistar/(psii*psii*psii) * siri, si, tmpvec2);
3053 rvec_dec(tmpvec, tmpvec2);
3054 cprod(tmpvec, rotg->vec, tmpvec2);
3056 svmul(wi*siri, tmpvec2, tmpvec);
3058 rvec_inc(sumvec, tmpvec);
3060 svmul(rotg->k*erg->invmass, sumvec, innersumvec);
3064 /* Calculate the radial motion 2 potential and forces */
3065 static void do_radial_motion2(
3066 t_rotgrp *rotg, /* The rotation group */
3067 rvec x[], /* The positions */
3068 matrix box, /* The simulation box */
3069 gmx_bool bOutstepRot, /* Output to main rotation output file */
3070 gmx_bool bOutstepSlab) /* Output per-slab data */
3072 int ii, iigrp, ifit, j;
3073 rvec xj; /* Position */
3074 real alpha; /* a single angle between an actual and a reference position */
3075 real weight; /* single weight for a single angle */
3076 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
3077 rvec xj_u; /* xj - u */
3078 rvec yj0_yc0; /* yj0 -yc0 */
3079 rvec tmpvec, tmpvec2;
3080 real fac, fit_fac, fac2, Vpart = 0.0;
3081 rvec rj, fit_rj, sj;
3083 rvec v_xj_u; /* v x (xj - u) */
3084 real psij, psijstar;
3085 real mj, wj; /* For mass-weighting of the positions */
3089 gmx_bool bCalcPotFit;
3092 erg = rotg->enfrotgrp;
3094 bPF = rotg->eType == erotgRM2PF;
3095 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT == rotg->eFittype);
3098 clear_rvec(yj0_yc0); /* Make the compiler happy */
3100 clear_rvec(innersumvec);
3103 /* For the pivot-free variant we have to use the current center of
3104 * mass of the rotation group instead of the pivot u */
3105 get_center(erg->xc, erg->mc, rotg->nat, erg->xc_center);
3107 /* Also, we precalculate the second term of the forces that is identical
3108 * (up to the weight factor mj) for all forces */
3109 radial_motion2_precalc_inner_sum(rotg, innersumvec);
3112 N_M = rotg->nat * erg->invmass;
3114 /* Each process calculates the forces on its local atoms */
3115 for (j = 0; j < erg->nat_loc; j++)
3119 /* Local index of a rotation group atom */
3120 ii = erg->ind_loc[j];
3121 /* Position of this atom in the collective array */
3122 iigrp = erg->xc_ref_ind[j];
3123 /* Mass-weighting */
3124 mj = erg->mc[iigrp];
3126 /* Current position of this atom: x[ii] */
3127 copy_rvec(x[ii], xj);
3129 /* Shift this atom such that it is near its reference */
3130 shift_single_coord(box, xj, erg->xc_shifts[iigrp]);
3132 /* The (unrotated) reference position is yj0. yc0 has already
3133 * been subtracted in init_rot_group */
3134 copy_rvec(rotg->x_ref[iigrp], yj0_yc0); /* yj0_yc0 = yj0 - yc0 */
3136 /* Calculate Omega.(yj0-yc0) */
3137 mvmul(erg->rotmat, yj0_yc0, rj); /* rj = Omega.(yj0-yc0) */
3142 copy_rvec(erg->x_loc_pbc[j], xj);
3143 copy_rvec(erg->xr_loc[j], rj); /* rj = Omega.(yj0-u) */
3145 /* Mass-weighting */
3148 /* Calculate (xj-u) resp. (xj-xc) */
3149 rvec_sub(xj, erg->xc_center, xj_u); /* xj_u = xj-u */
3151 cprod(rotg->vec, xj_u, v_xj_u); /* v_xj_u = v x (xj-u) */
3153 fac = norm2(v_xj_u);
3155 psijstar = 1.0/(fac + rotg->eps); /* psistar = -------------------- */
3156 /* |v x (xj-u)|^2 + eps */
3158 psij = gmx_invsqrt(fac); /* 1 */
3159 /* psij = ------------ */
3162 svmul(psij, v_xj_u, sj); /* sj = psij * (v x (xj-u) ) */
3164 fac = iprod(v_xj_u, rj); /* fac = (v x (xj-u)).rj */
3167 sjrj = iprod(sj, rj); /* sjrj = sj.rj */
3169 svmul(psijstar/psij, rj, tmpvec);
3170 svmul(psijstar*psijstar/(psij*psij*psij) * sjrj, sj, tmpvec2);
3171 rvec_dec(tmpvec, tmpvec2);
3172 cprod(tmpvec, rotg->vec, tmpvec2);
3174 /* Store the additional force so that it can be added to the force
3175 * array after the normal forces have been evaluated */
3176 svmul(-rotg->k*wj*sjrj, tmpvec2, tmpvec);
3177 svmul(mj, innersumvec, tmpvec2); /* This is != 0 only for the pivot-free variant */
3179 rvec_add(tmpvec2, tmpvec, erg->f_rot_loc[j]);
3180 Vpart += wj*psijstar*fac2;
3182 /* If requested, also calculate the potential for a set of angles
3183 * near the current reference angle */
3186 for (ifit = 0; ifit < rotg->PotAngle_nstep; ifit++)
3190 mvmul(erg->PotAngleFit->rotmat[ifit], yj0_yc0, fit_rj); /* fit_rj = Omega.(yj0-yc0) */
3194 /* Position of this atom in the collective array */
3195 iigrp = erg->xc_ref_ind[j];
3196 /* Rotate with the alternative angle. Like rotate_local_reference(),
3197 * just for a single local atom */
3198 mvmul(erg->PotAngleFit->rotmat[ifit], rotg->x_ref[iigrp], fit_rj); /* fit_rj = Omega*(yj0-u) */
3200 fit_fac = iprod(v_xj_u, fit_rj); /* fac = (v x (xj-u)).fit_rj */
3201 /* Add to the rotation potential for this angle: */
3202 erg->PotAngleFit->V[ifit] += 0.5*rotg->k*wj*psijstar*fit_fac*fit_fac;
3208 /* Add to the torque of this rotation group */
3209 erg->torque_v += torque(rotg->vec, erg->f_rot_loc[j], xj, erg->xc_center);
3211 /* Calculate the angle between reference and actual rotation group atom. */
3212 angle(rotg, xj_u, rj, &alpha, &weight); /* angle in rad, weighted */
3213 erg->angle_v += alpha * weight;
3214 erg->weight_v += weight;
3219 } /* end of loop over local rotation group atoms */
3220 erg->V = 0.5*rotg->k*Vpart;
3224 /* Determine the smallest and largest position vector (with respect to the
3225 * rotation vector) for the reference group */
3226 static void get_firstlast_atom_ref(
3231 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
3233 real xcproj; /* The projection of a reference position on the
3235 real minproj, maxproj; /* Smallest and largest projection on v */
3239 erg = rotg->enfrotgrp;
3241 /* Start with some value */
3242 minproj = iprod(rotg->x_ref[0], rotg->vec);
3245 /* This is just to ensure that it still works if all the atoms of the
3246 * reference structure are situated in a plane perpendicular to the rotation
3249 *lastindex = rotg->nat-1;
3251 /* Loop over all atoms of the reference group,
3252 * project them on the rotation vector to find the extremes */
3253 for (i = 0; i < rotg->nat; i++)
3255 xcproj = iprod(rotg->x_ref[i], rotg->vec);
3256 if (xcproj < minproj)
3261 if (xcproj > maxproj)
3270 /* Allocate memory for the slabs */
3271 static void allocate_slabs(
3277 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
3281 erg = rotg->enfrotgrp;
3283 /* More slabs than are defined for the reference are never needed */
3284 nslabs = erg->slab_last_ref - erg->slab_first_ref + 1;
3286 /* Remember how many we allocated */
3287 erg->nslabs_alloc = nslabs;
3289 if ( (NULL != fplog) && bVerbose)
3291 fprintf(fplog, "%s allocating memory to store data for %d slabs (rotation group %d).\n",
3294 snew(erg->slab_center, nslabs);
3295 snew(erg->slab_center_ref, nslabs);
3296 snew(erg->slab_weights, nslabs);
3297 snew(erg->slab_torque_v, nslabs);
3298 snew(erg->slab_data, nslabs);
3299 snew(erg->gn_atom, nslabs);
3300 snew(erg->gn_slabind, nslabs);
3301 snew(erg->slab_innersumvec, nslabs);
3302 for (i = 0; i < nslabs; i++)
3304 snew(erg->slab_data[i].x, rotg->nat);
3305 snew(erg->slab_data[i].ref, rotg->nat);
3306 snew(erg->slab_data[i].weight, rotg->nat);
3308 snew(erg->xc_ref_sorted, rotg->nat);
3309 snew(erg->xc_sortind, rotg->nat);
3310 snew(erg->firstatom, nslabs);
3311 snew(erg->lastatom, nslabs);
3315 /* From the extreme positions of the reference group, determine the first
3316 * and last slab of the reference. We can never have more slabs in the real
3317 * simulation than calculated here for the reference.
3319 static void get_firstlast_slab_ref(t_rotgrp *rotg, real mc[], int ref_firstindex, int ref_lastindex)
3321 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
3326 erg = rotg->enfrotgrp;
3327 first = get_first_slab(rotg, erg->max_beta, rotg->x_ref[ref_firstindex]);
3328 last = get_last_slab( rotg, erg->max_beta, rotg->x_ref[ref_lastindex ]);
3330 while (get_slab_weight(first, rotg, rotg->x_ref, mc, &dummy) > WEIGHT_MIN)
3334 erg->slab_first_ref = first+1;
3335 while (get_slab_weight(last, rotg, rotg->x_ref, mc, &dummy) > WEIGHT_MIN)
3339 erg->slab_last_ref = last-1;
3343 /* Special version of copy_rvec:
3344 * During the copy procedure of xcurr to b, the correct PBC image is chosen
3345 * such that the copied vector ends up near its reference position xref */
3346 static inline void copy_correct_pbc_image(
3347 const rvec xcurr, /* copy vector xcurr ... */
3348 rvec b, /* ... to b ... */
3349 const rvec xref, /* choosing the PBC image such that b ends up near xref */
3358 /* Shortest PBC distance between the atom and its reference */
3359 rvec_sub(xcurr, xref, dx);
3361 /* Determine the shift for this atom */
3363 for (m = npbcdim-1; m >= 0; m--)
3365 while (dx[m] < -0.5*box[m][m])
3367 for (d = 0; d < DIM; d++)
3373 while (dx[m] >= 0.5*box[m][m])
3375 for (d = 0; d < DIM; d++)
3383 /* Apply the shift to the position */
3384 copy_rvec(xcurr, b);
3385 shift_single_coord(box, b, shift);
3389 static void init_rot_group(FILE *fplog, t_commrec *cr, int g, t_rotgrp *rotg,
3390 rvec *x, gmx_mtop_t *mtop, gmx_bool bVerbose, FILE *out_slabs, matrix box,
3391 t_inputrec *ir, gmx_bool bOutputCenters)
3394 rvec coord, xref, *xdum;
3395 gmx_bool bFlex, bColl;
3397 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
3398 int ref_firstindex, ref_lastindex;
3399 gmx_mtop_atomlookup_t alook = NULL;
3400 real mass, totalmass;
3405 /* Do we have a flexible axis? */
3406 bFlex = ISFLEX(rotg);
3407 /* Do we use a global set of coordinates? */
3408 bColl = ISCOLL(rotg);
3410 erg = rotg->enfrotgrp;
3412 /* Allocate space for collective coordinates if needed */
3415 snew(erg->xc, rotg->nat);
3416 snew(erg->xc_shifts, rotg->nat);
3417 snew(erg->xc_eshifts, rotg->nat);
3418 snew(erg->xc_old, rotg->nat);
3420 if (rotg->eFittype == erotgFitNORM)
3422 snew(erg->xc_ref_length, rotg->nat); /* in case fit type NORM is chosen */
3423 snew(erg->xc_norm, rotg->nat);
3428 snew(erg->xr_loc, rotg->nat);
3429 snew(erg->x_loc_pbc, rotg->nat);
3432 snew(erg->f_rot_loc, rotg->nat);
3433 snew(erg->xc_ref_ind, rotg->nat);
3435 /* Make space for the calculation of the potential at other angles (used
3436 * for fitting only) */
3437 if (erotgFitPOT == rotg->eFittype)
3439 snew(erg->PotAngleFit, 1);
3440 snew(erg->PotAngleFit->degangle, rotg->PotAngle_nstep);
3441 snew(erg->PotAngleFit->V, rotg->PotAngle_nstep);
3442 snew(erg->PotAngleFit->rotmat, rotg->PotAngle_nstep);
3444 /* Get the set of angles around the reference angle */
3445 start = -0.5 * (rotg->PotAngle_nstep - 1)*rotg->PotAngle_step;
3446 for (i = 0; i < rotg->PotAngle_nstep; i++)
3448 erg->PotAngleFit->degangle[i] = start + i*rotg->PotAngle_step;
3453 erg->PotAngleFit = NULL;
3456 /* xc_ref_ind needs to be set to identity in the serial case */
3459 for (i = 0; i < rotg->nat; i++)
3461 erg->xc_ref_ind[i] = i;
3465 /* Copy the masses so that the center can be determined. For all types of
3466 * enforced rotation, we store the masses in the erg->mc array. */
3469 alook = gmx_mtop_atomlookup_init(mtop);
3471 snew(erg->mc, rotg->nat);
3474 snew(erg->mc_sorted, rotg->nat);
3478 snew(erg->m_loc, rotg->nat);
3481 for (i = 0; i < rotg->nat; i++)
3485 gmx_mtop_atomnr_to_atom(alook, rotg->ind[i], &atom);
3495 erg->invmass = 1.0/totalmass;
3499 gmx_mtop_atomlookup_destroy(alook);
3502 /* Set xc_ref_center for any rotation potential */
3503 if ((rotg->eType == erotgISO) || (rotg->eType == erotgPM) || (rotg->eType == erotgRM) || (rotg->eType == erotgRM2))
3505 /* Set the pivot point for the fixed, stationary-axis potentials. This
3506 * won't change during the simulation */
3507 copy_rvec(rotg->pivot, erg->xc_ref_center);
3508 copy_rvec(rotg->pivot, erg->xc_center );
3512 /* Center of the reference positions */
3513 get_center(rotg->x_ref, erg->mc, rotg->nat, erg->xc_ref_center);
3515 /* Center of the actual positions */
3518 snew(xdum, rotg->nat);
3519 for (i = 0; i < rotg->nat; i++)
3522 copy_rvec(x[ii], xdum[i]);
3524 get_center(xdum, erg->mc, rotg->nat, erg->xc_center);
3530 gmx_bcast(sizeof(erg->xc_center), erg->xc_center, cr);
3537 /* Save the original (whole) set of positions in xc_old such that at later
3538 * steps the rotation group can always be made whole again. If the simulation is
3539 * restarted, we compute the starting reference positions (given the time)
3540 * and assume that the correct PBC image of each position is the one nearest
3541 * to the current reference */
3544 /* Calculate the rotation matrix for this angle: */
3545 t_start = ir->init_t + ir->init_step*ir->delta_t;
3546 erg->degangle = rotg->rate * t_start;
3547 calc_rotmat(rotg->vec, erg->degangle, erg->rotmat);
3549 for (i = 0; i < rotg->nat; i++)
3553 /* Subtract pivot, rotate, and add pivot again. This will yield the
3554 * reference position for time t */
3555 rvec_sub(rotg->x_ref[i], erg->xc_ref_center, coord);
3556 mvmul(erg->rotmat, coord, xref);
3557 rvec_inc(xref, erg->xc_ref_center);
3559 copy_correct_pbc_image(x[ii], erg->xc_old[i], xref, box, 3);
3565 gmx_bcast(rotg->nat*sizeof(erg->xc_old[0]), erg->xc_old, cr);
3570 if ( (rotg->eType != erotgFLEX) && (rotg->eType != erotgFLEX2) )
3572 /* Put the reference positions into origin: */
3573 for (i = 0; i < rotg->nat; i++)
3575 rvec_dec(rotg->x_ref[i], erg->xc_ref_center);
3579 /* Enforced rotation with flexible axis */
3582 /* Calculate maximum beta value from minimum gaussian (performance opt.) */
3583 erg->max_beta = calc_beta_max(rotg->min_gaussian, rotg->slab_dist);
3585 /* Determine the smallest and largest coordinate with respect to the rotation vector */
3586 get_firstlast_atom_ref(rotg, &ref_firstindex, &ref_lastindex);
3588 /* From the extreme positions of the reference group, determine the first
3589 * and last slab of the reference. */
3590 get_firstlast_slab_ref(rotg, erg->mc, ref_firstindex, ref_lastindex);
3592 /* Allocate memory for the slabs */
3593 allocate_slabs(rotg, fplog, g, bVerbose);
3595 /* Flexible rotation: determine the reference centers for the rest of the simulation */
3596 erg->slab_first = erg->slab_first_ref;
3597 erg->slab_last = erg->slab_last_ref;
3598 get_slab_centers(rotg, rotg->x_ref, erg->mc, g, -1, out_slabs, bOutputCenters, TRUE);
3600 /* Length of each x_rotref vector from center (needed if fit routine NORM is chosen): */
3601 if (rotg->eFittype == erotgFitNORM)
3603 for (i = 0; i < rotg->nat; i++)
3605 rvec_sub(rotg->x_ref[i], erg->xc_ref_center, coord);
3606 erg->xc_ref_length[i] = norm(coord);
3613 extern void dd_make_local_rotation_groups(gmx_domdec_t *dd, t_rot *rot)
3618 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
3622 for (g = 0; g < rot->ngrp; g++)
3624 rotg = &rot->grp[g];
3625 erg = rotg->enfrotgrp;
3628 dd_make_local_group_indices(ga2la, rotg->nat, rotg->ind,
3629 &erg->nat_loc, &erg->ind_loc, &erg->nalloc_loc, erg->xc_ref_ind);
3634 /* Calculate the size of the MPI buffer needed in reduce_output() */
3635 static int calc_mpi_bufsize(t_rot *rot)
3638 int count_group, count_total;
3640 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
3644 for (g = 0; g < rot->ngrp; g++)
3646 rotg = &rot->grp[g];
3647 erg = rotg->enfrotgrp;
3649 /* Count the items that are transferred for this group: */
3650 count_group = 4; /* V, torque, angle, weight */
3652 /* Add the maximum number of slabs for flexible groups */
3655 count_group += erg->slab_last_ref - erg->slab_first_ref + 1;
3658 /* Add space for the potentials at different angles: */
3659 if (erotgFitPOT == rotg->eFittype)
3661 count_group += rotg->PotAngle_nstep;
3664 /* Add to the total number: */
3665 count_total += count_group;
3672 extern void init_rot(FILE *fplog, t_inputrec *ir, int nfile, const t_filenm fnm[],
3673 t_commrec *cr, rvec *x, matrix box, gmx_mtop_t *mtop, const output_env_t oenv,
3674 gmx_bool bVerbose, unsigned long Flags)
3679 int nat_max = 0; /* Size of biggest rotation group */
3680 gmx_enfrot_t er; /* Pointer to the enforced rotation buffer variables */
3681 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
3682 rvec *x_pbc = NULL; /* Space for the pbc-correct atom positions */
3685 if ( (PAR(cr)) && !DOMAINDECOMP(cr) )
3687 gmx_fatal(FARGS, "Enforced rotation is only implemented for domain decomposition!");
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);