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
42 #include "gmx_wallcycle.h"
51 #include "mtop_util.h"
55 #include "gmx_ga2la.h"
58 #include "groupcoord.h"
59 #include "pull_rotation.h"
63 static char *RotStr = {"Enforced rotation:"};
66 /* Set the minimum weight for the determination of the slab centers */
67 #define WEIGHT_MIN (10*GMX_FLOAT_MIN)
69 /* Helper structure for sorting positions along rotation vector */
71 real xcproj; /* Projection of xc on the rotation vector */
72 int ind; /* Index of xc */
74 rvec x; /* Position */
75 rvec x_ref; /* Reference position */
79 /* Enforced rotation / flexible: determine the angle of each slab */
80 typedef struct gmx_slabdata
82 int nat; /* Number of atoms belonging to this slab */
83 rvec *x; /* The positions belonging to this slab. In
84 general, this should be all positions of the
85 whole rotation group, but we leave those away
86 that have a small enough weight */
87 rvec *ref; /* Same for reference */
88 real *weight; /* The weight for each atom */
92 /* Helper structure for potential fitting */
93 typedef struct gmx_potfit
95 real *degangle; /* Set of angles for which the potential is
96 calculated. The optimum fit is determined as
97 the angle for with the potential is minimal */
98 real *V; /* Potential for the different angles */
99 matrix *rotmat; /* Rotation matrix corresponding to the angles */
103 /* Enforced rotation data for all groups */
104 typedef struct gmx_enfrot
106 FILE *out_rot; /* Output file for rotation data */
107 FILE *out_torque; /* Output file for torque data */
108 FILE *out_angles; /* Output file for slab angles for flexible type */
109 FILE *out_slabs; /* Output file for slab centers */
110 int bufsize; /* Allocation size of buf */
111 rvec *xbuf; /* Coordinate buffer variable for sorting */
112 real *mbuf; /* Masses buffer variable for sorting */
113 sort_along_vec_t *data; /* Buffer variable needed for position sorting */
114 real *mpi_inbuf; /* MPI buffer */
115 real *mpi_outbuf; /* MPI buffer */
116 int mpi_bufsize; /* Allocation size of in & outbuf */
117 unsigned long Flags; /* mdrun flags */
118 gmx_bool bOut; /* Used to skip first output when appending to
119 * avoid duplicate entries in rotation outfiles */
123 /* Global enforced rotation data for a single rotation group */
124 typedef struct gmx_enfrotgrp
126 real degangle; /* Rotation angle in degrees */
127 matrix rotmat; /* Rotation matrix */
128 atom_id *ind_loc; /* Local rotation indices */
129 int nat_loc; /* Number of local group atoms */
130 int nalloc_loc; /* Allocation size for ind_loc and weight_loc */
132 real V; /* Rotation potential for this rotation group */
133 rvec *f_rot_loc; /* Array to store the forces on the local atoms
134 resulting from enforced rotation potential */
136 /* Collective coordinates for the whole rotation group */
137 real *xc_ref_length; /* Length of each x_rotref vector after x_rotref
138 has been put into origin */
139 int *xc_ref_ind; /* Position of each local atom in the collective
141 rvec xc_center; /* Center of the rotation group positions, may
143 rvec xc_ref_center; /* dito, for the reference positions */
144 rvec *xc; /* Current (collective) positions */
145 ivec *xc_shifts; /* Current (collective) shifts */
146 ivec *xc_eshifts; /* Extra shifts since last DD step */
147 rvec *xc_old; /* Old (collective) positions */
148 rvec *xc_norm; /* Normalized form of the current positions */
149 rvec *xc_ref_sorted; /* Reference positions (sorted in the same order
150 as xc when sorted) */
151 int *xc_sortind; /* Where is a position found after sorting? */
152 real *mc; /* Collective masses */
154 real invmass; /* one over the total mass of the rotation group */
156 real torque_v; /* Torque in the direction of rotation vector */
157 real angle_v; /* Actual angle of the whole rotation group */
158 /* Fixed rotation only */
159 real weight_v; /* Weights for angle determination */
160 rvec *xr_loc; /* Local reference coords, correctly rotated */
161 rvec *x_loc_pbc; /* Local current coords, correct PBC image */
162 real *m_loc; /* Masses of the current local atoms */
164 /* Flexible rotation only */
165 int nslabs_alloc; /* For this many slabs memory is allocated */
166 int slab_first; /* Lowermost slab for that the calculation needs
167 to be performed at a given time step */
168 int slab_last; /* Uppermost slab ... */
169 int slab_first_ref; /* First slab for which ref. center is stored */
170 int slab_last_ref; /* Last ... */
171 int slab_buffer; /* Slab buffer region around reference slabs */
172 int *firstatom; /* First relevant atom for a slab */
173 int *lastatom; /* Last relevant atom for a slab */
174 rvec *slab_center; /* Gaussian-weighted slab center */
175 rvec *slab_center_ref; /* Gaussian-weighted slab center for the
176 reference positions */
177 real *slab_weights; /* Sum of gaussian weights in a slab */
178 real *slab_torque_v; /* Torque T = r x f for each slab. */
179 /* torque_v = m.v = angular momentum in the
181 real max_beta; /* min_gaussian from inputrec->rotgrp is the
182 minimum value the gaussian must have so that
183 the force is actually evaluated max_beta is
184 just another way to put it */
185 real *gn_atom; /* Precalculated gaussians for a single atom */
186 int *gn_slabind; /* Tells to which slab each precalculated gaussian
188 rvec *slab_innersumvec;/* Inner sum of the flexible2 potential per slab;
189 this is precalculated for optimization reasons */
190 t_gmx_slabdata *slab_data; /* Holds atom positions and gaussian weights
191 of atoms belonging to a slab */
193 /* For potential fits with varying angle: */
194 t_gmx_potfit *PotAngleFit; /* Used for fit type 'potential' */
198 /* Activate output of forces for correctness checks */
199 /* #define PRINT_FORCES */
201 #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]);
202 #define PRINT_POT_TAU if (MASTER(cr)) { \
203 fprintf(stderr,"potential = %15.8f\n" "torque = %15.8f\n", erg->V, erg->torque_v); \
206 #define PRINT_FORCE_J
207 #define PRINT_POT_TAU
210 /* Shortcuts for often used queries */
211 #define ISFLEX(rg) ( (rg->eType==erotgFLEX) || (rg->eType==erotgFLEXT) || (rg->eType==erotgFLEX2) || (rg->eType==erotgFLEX2T) )
212 #define ISCOLL(rg) ( (rg->eType==erotgFLEX) || (rg->eType==erotgFLEXT) || (rg->eType==erotgFLEX2) || (rg->eType==erotgFLEX2T) || (rg->eType==erotgRMPF) || (rg->eType==erotgRM2PF) )
215 /* Does any of the rotation groups use slab decomposition? */
216 static gmx_bool HaveFlexibleGroups(t_rot *rot)
222 for (g=0; g<rot->ngrp; g++)
233 /* Is for any group the fit angle determined by finding the minimum of the
234 * rotation potential? */
235 static gmx_bool HavePotFitGroups(t_rot *rot)
241 for (g=0; g<rot->ngrp; g++)
244 if (erotgFitPOT == rotg->eFittype)
252 static double** allocate_square_matrix(int dim)
266 static void free_square_matrix(double** mat, int dim)
271 for (i=0; i<dim; i++)
277 /* Return the angle for which the potential is minimal */
278 static real get_fitangle(t_rotgrp *rotg, gmx_enfrotgrp_t erg)
281 real fitangle = -999.9;
282 real pot_min = GMX_FLOAT_MAX;
286 fit = erg->PotAngleFit;
288 for (i = 0; i < rotg->PotAngle_nstep; i++)
290 if (fit->V[i] < pot_min)
293 fitangle = fit->degangle[i];
301 /* Reduce potential angle fit data for this group at this time step? */
302 static gmx_inline gmx_bool bPotAngle(t_rot *rot, t_rotgrp *rotg, gmx_large_int_t step)
304 return ( (erotgFitPOT==rotg->eFittype) && (do_per_step(step, rot->nstsout) || do_per_step(step, rot->nstrout)) );
307 /* Reduce slab torqe data for this group at this time step? */
308 static gmx_inline gmx_bool bSlabTau(t_rot *rot, t_rotgrp *rotg, gmx_large_int_t step)
310 return ( (ISFLEX(rotg)) && do_per_step(step, rot->nstsout) );
313 /* Output rotation energy, torques, etc. for each rotation group */
314 static void reduce_output(t_commrec *cr, t_rot *rot, real t, gmx_large_int_t step)
316 int g,i,islab,nslabs=0;
317 int count; /* MPI element counter */
319 gmx_enfrot_t er; /* Pointer to the enforced rotation buffer variables */
320 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
327 /* Fill the MPI buffer with stuff to reduce. If items are added for reduction
328 * here, the MPI buffer size has to be enlarged also in calc_mpi_bufsize() */
332 for (g=0; g < rot->ngrp; g++)
335 erg = rotg->enfrotgrp;
336 nslabs = erg->slab_last - erg->slab_first + 1;
337 er->mpi_inbuf[count++] = erg->V;
338 er->mpi_inbuf[count++] = erg->torque_v;
339 er->mpi_inbuf[count++] = erg->angle_v;
340 er->mpi_inbuf[count++] = erg->weight_v; /* weights are not needed for flex types, but this is just a single value */
342 if (bPotAngle(rot, rotg, step))
344 for (i = 0; i < rotg->PotAngle_nstep; i++)
345 er->mpi_inbuf[count++] = erg->PotAngleFit->V[i];
347 if (bSlabTau(rot, rotg, step))
349 for (i=0; i<nslabs; i++)
350 er->mpi_inbuf[count++] = erg->slab_torque_v[i];
353 if (count > er->mpi_bufsize)
354 gmx_fatal(FARGS, "%s MPI buffer overflow, please report this error.", RotStr);
357 MPI_Reduce(er->mpi_inbuf, er->mpi_outbuf, count, GMX_MPI_REAL, MPI_SUM, MASTERRANK(cr), cr->mpi_comm_mygroup);
360 /* Copy back the reduced data from the buffer on the master */
364 for (g=0; g < rot->ngrp; g++)
367 erg = rotg->enfrotgrp;
368 nslabs = erg->slab_last - erg->slab_first + 1;
369 erg->V = er->mpi_outbuf[count++];
370 erg->torque_v = er->mpi_outbuf[count++];
371 erg->angle_v = er->mpi_outbuf[count++];
372 erg->weight_v = er->mpi_outbuf[count++];
374 if (bPotAngle(rot, rotg, step))
376 for (i = 0; i < rotg->PotAngle_nstep; i++)
377 erg->PotAngleFit->V[i] = er->mpi_outbuf[count++];
379 if (bSlabTau(rot, rotg, step))
381 for (i=0; i<nslabs; i++)
382 erg->slab_torque_v[i] = er->mpi_outbuf[count++];
391 /* Angle and torque for each rotation group */
392 for (g=0; g < rot->ngrp; g++)
395 bFlex = ISFLEX(rotg);
399 /* Output to main rotation output file: */
400 if ( do_per_step(step, rot->nstrout) )
402 if (erotgFitPOT == rotg->eFittype)
404 fitangle = get_fitangle(rotg, erg);
409 fitangle = erg->angle_v; /* RMSD fit angle */
411 fitangle = (erg->angle_v/erg->weight_v)*180.0*M_1_PI;
413 fprintf(er->out_rot, "%12.4f", fitangle);
414 fprintf(er->out_rot, "%12.3e", erg->torque_v);
415 fprintf(er->out_rot, "%12.3e", erg->V);
418 if ( do_per_step(step, rot->nstsout) )
420 /* Output to torque log file: */
423 fprintf(er->out_torque, "%12.3e%6d", t, g);
424 for (i=erg->slab_first; i<=erg->slab_last; i++)
426 islab = i - erg->slab_first; /* slab index */
427 /* Only output if enough weight is in slab */
428 if (erg->slab_weights[islab] > rotg->min_gaussian)
429 fprintf(er->out_torque, "%6d%12.3e", i, erg->slab_torque_v[islab]);
431 fprintf(er->out_torque , "\n");
434 /* Output to angles log file: */
435 if (erotgFitPOT == rotg->eFittype)
437 fprintf(er->out_angles, "%12.3e%6d%12.4f", t, g, erg->degangle);
438 /* Output energies at a set of angles around the reference angle */
439 for (i = 0; i < rotg->PotAngle_nstep; i++)
440 fprintf(er->out_angles, "%12.3e", erg->PotAngleFit->V[i]);
441 fprintf(er->out_angles, "\n");
445 if ( do_per_step(step, rot->nstrout) )
446 fprintf(er->out_rot, "\n");
451 /* Add the forces from enforced rotation potential to the local forces.
452 * Should be called after the SR forces have been evaluated */
453 extern real add_rot_forces(t_rot *rot, rvec f[], t_commrec *cr, gmx_large_int_t step, real t)
457 gmx_enfrot_t er; /* Pointer to the enforced rotation buffer variables */
458 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
459 real Vrot = 0.0; /* If more than one rotation group is present, Vrot
460 assembles the local parts from all groups */
465 /* Loop over enforced rotation groups (usually 1, though)
466 * Apply the forces from rotation potentials */
467 for (g=0; g<rot->ngrp; g++)
471 Vrot += erg->V; /* add the local parts from the nodes */
472 for (l=0; l<erg->nat_loc; l++)
474 /* Get the right index of the local force */
475 ii = erg->ind_loc[l];
477 rvec_inc(f[ii],erg->f_rot_loc[l]);
481 /* Reduce energy,torque, angles etc. to get the sum values (per rotation group)
482 * on the master and output these values to file. */
483 if ( (do_per_step(step, rot->nstrout) || do_per_step(step, rot->nstsout)) && er->bOut)
484 reduce_output(cr, rot, t, step);
486 /* When appending, er->bOut is FALSE the first time to avoid duplicate entries */
495 /* The Gaussian norm is chosen such that the sum of the gaussian functions
496 * over the slabs is approximately 1.0 everywhere */
497 #define GAUSS_NORM 0.569917543430618
500 /* Calculate the maximum beta that leads to a gaussian larger min_gaussian,
501 * also does some checks
503 static double calc_beta_max(real min_gaussian, real slab_dist)
509 /* Actually the next two checks are already made in grompp */
511 gmx_fatal(FARGS, "Slab distance of flexible rotation groups must be >=0 !");
512 if (min_gaussian <= 0)
513 gmx_fatal(FARGS, "Cutoff value for Gaussian must be > 0. (You requested %f)");
515 /* Define the sigma value */
516 sigma = 0.7*slab_dist;
518 /* Calculate the argument for the logarithm and check that the log() result is negative or 0 */
519 arg = min_gaussian/GAUSS_NORM;
521 gmx_fatal(FARGS, "min_gaussian of flexible rotation groups must be <%g", GAUSS_NORM);
523 return sqrt(-2.0*sigma*sigma*log(min_gaussian/GAUSS_NORM));
527 static gmx_inline real calc_beta(rvec curr_x, t_rotgrp *rotg, int n)
529 return iprod(curr_x, rotg->vec) - rotg->slab_dist * n;
533 static gmx_inline real gaussian_weight(rvec curr_x, t_rotgrp *rotg, int n)
535 const real norm = GAUSS_NORM;
539 /* Define the sigma value */
540 sigma = 0.7*rotg->slab_dist;
541 /* Calculate the Gaussian value of slab n for position curr_x */
542 return norm * exp( -0.5 * sqr( calc_beta(curr_x, rotg, n)/sigma ) );
546 /* Returns the weight in a single slab, also calculates the Gaussian- and mass-
547 * weighted sum of positions for that slab */
548 static real get_slab_weight(int j, t_rotgrp *rotg, rvec xc[], real mc[], rvec *x_weighted_sum)
550 rvec curr_x; /* The position of an atom */
551 rvec curr_x_weighted; /* The gaussian-weighted position */
552 real gaussian; /* A single gaussian weight */
553 real wgauss; /* gaussian times current mass */
554 real slabweight = 0.0; /* The sum of weights in the slab */
556 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
560 clear_rvec(*x_weighted_sum);
563 islab = j - erg->slab_first;
565 /* Loop over all atoms in the rotation group */
566 for (i=0; i<rotg->nat; i++)
568 copy_rvec(xc[i], curr_x);
569 gaussian = gaussian_weight(curr_x, rotg, j);
570 wgauss = gaussian * mc[i];
571 svmul(wgauss, curr_x, curr_x_weighted);
572 rvec_add(*x_weighted_sum, curr_x_weighted, *x_weighted_sum);
573 slabweight += wgauss;
574 } /* END of loop over rotation group atoms */
580 static void get_slab_centers(
581 t_rotgrp *rotg, /* The rotation group information */
582 rvec *xc, /* The rotation group positions; will
583 typically be enfrotgrp->xc, but at first call
584 it is enfrotgrp->xc_ref */
585 real *mc, /* The masses of the rotation group atoms */
586 t_commrec *cr, /* Communication record */
587 int g, /* The number of the rotation group */
588 real time, /* Used for output only */
589 FILE *out_slabs, /* For outputting center per slab information */
590 gmx_bool bOutStep, /* Is this an output step? */
591 gmx_bool bReference) /* If this routine is called from
592 init_rot_group we need to store
593 the reference slab centers */
596 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
601 /* Loop over slabs */
602 for (j = erg->slab_first; j <= erg->slab_last; j++)
604 islab = j - erg->slab_first;
605 erg->slab_weights[islab] = get_slab_weight(j, rotg, xc, mc, &erg->slab_center[islab]);
607 /* We can do the calculations ONLY if there is weight in the slab! */
608 if (erg->slab_weights[islab] > WEIGHT_MIN)
610 svmul(1.0/erg->slab_weights[islab], erg->slab_center[islab], erg->slab_center[islab]);
614 /* We need to check this here, since we divide through slab_weights
615 * in the flexible low-level routines! */
616 gmx_fatal(FARGS, "Not enough weight in slab %d. Slab center cannot be determined!", j);
619 /* At first time step: save the centers of the reference structure */
621 copy_rvec(erg->slab_center[islab], erg->slab_center_ref[islab]);
622 } /* END of loop over slabs */
624 /* Output on the master */
625 if (MASTER(cr) && bOutStep)
627 fprintf(out_slabs, "%12.3e%6d", time, g);
628 for (j = erg->slab_first; j <= erg->slab_last; j++)
630 islab = j - erg->slab_first;
631 fprintf(out_slabs, "%6d%12.3e%12.3e%12.3e",
632 j,erg->slab_center[islab][XX],erg->slab_center[islab][YY],erg->slab_center[islab][ZZ]);
634 fprintf(out_slabs, "\n");
639 static void calc_rotmat(
641 real degangle, /* Angle alpha of rotation at time t in degrees */
642 matrix rotmat) /* Rotation matrix */
644 real radangle; /* Rotation angle in radians */
645 real cosa; /* cosine alpha */
646 real sina; /* sine alpha */
647 real OMcosa; /* 1 - cos(alpha) */
648 real dumxy, dumxz, dumyz; /* save computations */
649 rvec rot_vec; /* Rotate around rot_vec ... */
652 radangle = degangle * M_PI/180.0;
653 copy_rvec(vec , rot_vec );
655 /* Precompute some variables: */
656 cosa = cos(radangle);
657 sina = sin(radangle);
659 dumxy = rot_vec[XX]*rot_vec[YY]*OMcosa;
660 dumxz = rot_vec[XX]*rot_vec[ZZ]*OMcosa;
661 dumyz = rot_vec[YY]*rot_vec[ZZ]*OMcosa;
663 /* Construct the rotation matrix for this rotation group: */
665 rotmat[XX][XX] = cosa + rot_vec[XX]*rot_vec[XX]*OMcosa;
666 rotmat[YY][XX] = dumxy + rot_vec[ZZ]*sina;
667 rotmat[ZZ][XX] = dumxz - rot_vec[YY]*sina;
669 rotmat[XX][YY] = dumxy - rot_vec[ZZ]*sina;
670 rotmat[YY][YY] = cosa + rot_vec[YY]*rot_vec[YY]*OMcosa;
671 rotmat[ZZ][YY] = dumyz + rot_vec[XX]*sina;
673 rotmat[XX][ZZ] = dumxz + rot_vec[YY]*sina;
674 rotmat[YY][ZZ] = dumyz - rot_vec[XX]*sina;
675 rotmat[ZZ][ZZ] = cosa + rot_vec[ZZ]*rot_vec[ZZ]*OMcosa;
680 for (iii=0; iii<3; iii++) {
681 for (jjj=0; jjj<3; jjj++)
682 fprintf(stderr, " %10.8f ", rotmat[iii][jjj]);
683 fprintf(stderr, "\n");
689 /* Calculates torque on the rotation axis tau = position x force */
690 static gmx_inline real torque(
691 rvec rotvec, /* rotation vector; MUST be normalized! */
692 rvec force, /* force */
693 rvec x, /* position of atom on which the force acts */
694 rvec pivot) /* pivot point of rotation axis */
699 /* Subtract offset */
700 rvec_sub(x,pivot,vectmp);
702 /* position x force */
703 cprod(vectmp, force, tau);
705 /* Return the part of the torque which is parallel to the rotation vector */
706 return iprod(tau, rotvec);
710 /* Right-aligned output of value with standard width */
711 static void print_aligned(FILE *fp, char *str)
713 fprintf(fp, "%12s", str);
717 /* Right-aligned output of value with standard short width */
718 static void print_aligned_short(FILE *fp, char *str)
720 fprintf(fp, "%6s", str);
724 static FILE *open_output_file(const char *fn, int steps, const char what[])
729 fp = ffopen(fn, "w");
731 fprintf(fp, "# Output of %s is written in intervals of %d time step%s.\n#\n",
732 what,steps, steps>1 ? "s":"");
738 /* Open output file for slab center data. Call on master only */
739 static FILE *open_slab_out(const char *fn, t_rot *rot, const output_env_t oenv)
746 if (rot->enfrot->Flags & MD_APPENDFILES)
748 fp = gmx_fio_fopen(fn,"a");
752 fp = open_output_file(fn, rot->nstsout, "gaussian weighted slab centers");
754 for (g=0; g<rot->ngrp; g++)
759 fprintf(fp, "# Rotation group %d (%s), slab distance %f nm, %s.\n",
760 g, erotg_names[rotg->eType], rotg->slab_dist,
761 rotg->bMassW? "centers of mass":"geometrical centers");
765 fprintf(fp, "# Reference centers are listed first (t=-1).\n");
766 fprintf(fp, "# The following columns have the syntax:\n");
768 print_aligned_short(fp, "t");
769 print_aligned_short(fp, "grp");
770 /* Print legend for the first two entries only ... */
773 print_aligned_short(fp, "slab");
774 print_aligned(fp, "X center");
775 print_aligned(fp, "Y center");
776 print_aligned(fp, "Z center");
778 fprintf(fp, " ...\n");
786 /* Adds 'buf' to 'str' */
787 static void add_to_string(char **str, char *buf)
792 len = strlen(*str) + strlen(buf) + 1;
798 static void add_to_string_aligned(char **str, char *buf)
800 char buf_aligned[STRLEN];
802 sprintf(buf_aligned, "%12s", buf);
803 add_to_string(str, buf_aligned);
807 /* Open output file and print some general information about the rotation groups.
808 * Call on master only */
809 static FILE *open_rot_out(const char *fn, t_rot *rot, const output_env_t oenv)
814 const char **setname;
815 char buf[50], buf2[75];
816 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
818 char *LegendStr=NULL;
821 if (rot->enfrot->Flags & MD_APPENDFILES)
823 fp = gmx_fio_fopen(fn,"a");
827 fp = xvgropen(fn, "Rotation angles and energy", "Time [ps]", "angles [degrees] and energies [kJ/mol]", oenv);
828 fprintf(fp, "# Output of enforced rotation data is written in intervals of %d time step%s.\n#\n", rot->nstrout, rot->nstrout > 1 ? "s":"");
829 fprintf(fp, "# The scalar tau is the torque [kJ/mol] in the direction of the rotation vector v.\n");
830 fprintf(fp, "# To obtain the vectorial torque, multiply tau with the group's rot_vec.\n");
831 fprintf(fp, "# For flexible groups, tau(t,n) from all slabs n have been summed in a single value tau(t) here.\n");
832 fprintf(fp, "# The torques tau(t,n) are found in the rottorque.log (-rt) output file\n");
834 for (g=0; g<rot->ngrp; g++)
838 bFlex = ISFLEX(rotg);
841 fprintf(fp, "# ROTATION GROUP %d, potential type '%s':\n" , g, erotg_names[rotg->eType]);
842 fprintf(fp, "# rot_massw%d %s\n" , g, yesno_names[rotg->bMassW]);
843 fprintf(fp, "# rot_vec%d %12.5e %12.5e %12.5e\n" , g, rotg->vec[XX], rotg->vec[YY], rotg->vec[ZZ]);
844 fprintf(fp, "# rot_rate%d %12.5e degrees/ps\n" , g, rotg->rate);
845 fprintf(fp, "# rot_k%d %12.5e kJ/(mol*nm^2)\n" , g, rotg->k);
846 if ( rotg->eType==erotgISO || rotg->eType==erotgPM || rotg->eType==erotgRM || rotg->eType==erotgRM2)
847 fprintf(fp, "# rot_pivot%d %12.5e %12.5e %12.5e nm\n", g, rotg->pivot[XX], rotg->pivot[YY], rotg->pivot[ZZ]);
851 fprintf(fp, "# rot_slab_distance%d %f nm\n", g, rotg->slab_dist);
852 fprintf(fp, "# rot_min_gaussian%d %12.5e\n", g, rotg->min_gaussian);
855 /* Output the centers of the rotation groups for the pivot-free potentials */
856 if ((rotg->eType==erotgISOPF) || (rotg->eType==erotgPMPF) || (rotg->eType==erotgRMPF) || (rotg->eType==erotgRM2PF
857 || (rotg->eType==erotgFLEXT) || (rotg->eType==erotgFLEX2T)) )
859 fprintf(fp, "# ref. grp. %d center %12.5e %12.5e %12.5e\n", g,
860 erg->xc_ref_center[XX], erg->xc_ref_center[YY], erg->xc_ref_center[ZZ]);
862 fprintf(fp, "# grp. %d init.center %12.5e %12.5e %12.5e\n", g,
863 erg->xc_center[XX], erg->xc_center[YY], erg->xc_center[ZZ]);
866 if ( (rotg->eType == erotgRM2) || (rotg->eType==erotgFLEX2) || (rotg->eType==erotgFLEX2T) )
868 fprintf(fp, "# rot_eps%d %12.5e nm^2\n", g, rotg->eps);
870 if (erotgFitPOT == rotg->eFittype)
873 fprintf(fp, "# theta_fit%d is determined by first evaluating the potential for %d angles around theta_ref%d.\n",
874 g, rotg->PotAngle_nstep, g);
875 fprintf(fp, "# The fit angle is the one with the smallest potential. It is given as the deviation\n");
876 fprintf(fp, "# from the reference angle, i.e. if theta_ref=X and theta_fit=Y, then the angle with\n");
877 fprintf(fp, "# minimal value of the potential is X+Y. Angular resolution is %g degrees.\n", rotg->PotAngle_step);
881 /* Print a nice legend */
884 sprintf(buf, "# %6s", "time");
885 add_to_string_aligned(&LegendStr, buf);
888 snew(setname, 4*rot->ngrp);
890 for (g=0; g<rot->ngrp; g++)
893 sprintf(buf, "theta_ref%d", g);
894 add_to_string_aligned(&LegendStr, buf);
896 sprintf(buf2, "%s [degrees]", buf);
897 setname[nsets] = strdup(buf2);
900 for (g=0; g<rot->ngrp; g++)
903 bFlex = ISFLEX(rotg);
905 /* For flexible axis rotation we use RMSD fitting to determine the
906 * actual angle of the rotation group */
907 if (bFlex || erotgFitPOT == rotg->eFittype)
908 sprintf(buf, "theta_fit%d", g);
910 sprintf(buf, "theta_av%d", g);
911 add_to_string_aligned(&LegendStr, buf);
912 sprintf(buf2, "%s [degrees]", buf);
913 setname[nsets] = strdup(buf2);
916 sprintf(buf, "tau%d", g);
917 add_to_string_aligned(&LegendStr, buf);
918 sprintf(buf2, "%s [kJ/mol]", buf);
919 setname[nsets] = strdup(buf2);
922 sprintf(buf, "energy%d", g);
923 add_to_string_aligned(&LegendStr, buf);
924 sprintf(buf2, "%s [kJ/mol]", buf);
925 setname[nsets] = strdup(buf2);
931 xvgr_legend(fp, nsets, setname, oenv);
934 fprintf(fp, "#\n# Legend for the following data columns:\n");
935 fprintf(fp, "%s\n", LegendStr);
945 /* Call on master only */
946 static FILE *open_angles_out(const char *fn, t_rot *rot, const output_env_t oenv)
951 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
955 if (rot->enfrot->Flags & MD_APPENDFILES)
957 fp = gmx_fio_fopen(fn,"a");
961 /* Open output file and write some information about it's structure: */
962 fp = open_output_file(fn, rot->nstsout, "rotation group angles");
963 fprintf(fp, "# All angles given in degrees, time in ps.\n");
964 for (g=0; g<rot->ngrp; g++)
969 /* Output for this group happens only if potential type is flexible or
970 * if fit type is potential! */
971 if ( ISFLEX(rotg) || (erotgFitPOT == rotg->eFittype) )
974 sprintf(buf, " slab distance %f nm, ", rotg->slab_dist);
978 fprintf(fp, "#\n# ROTATION GROUP %d '%s',%s fit type '%s'.\n",
979 g, erotg_names[rotg->eType], buf, erotg_fitnames[rotg->eFittype]);
981 /* Special type of fitting using the potential minimum. This is
982 * done for the whole group only, not for the individual slabs. */
983 if (erotgFitPOT == rotg->eFittype)
985 fprintf(fp, "# To obtain theta_fit%d, the potential is evaluated for %d angles around theta_ref%d\n", g, rotg->PotAngle_nstep, g);
986 fprintf(fp, "# The fit angle in the rotation standard outfile is the one with minimal energy E(theta_fit) [kJ/mol].\n");
990 fprintf(fp, "# Legend for the group %d data columns:\n", g);
992 print_aligned_short(fp, "time");
993 print_aligned_short(fp, "grp");
994 print_aligned(fp, "theta_ref");
996 if (erotgFitPOT == rotg->eFittype)
998 /* Output the set of angles around the reference angle */
999 for (i = 0; i < rotg->PotAngle_nstep; i++)
1001 sprintf(buf, "E(%g)", erg->PotAngleFit->degangle[i]);
1002 print_aligned(fp, buf);
1007 /* Output fit angle for each slab */
1008 print_aligned_short(fp, "slab");
1009 print_aligned_short(fp, "atoms");
1010 print_aligned(fp, "theta_fit");
1011 print_aligned_short(fp, "slab");
1012 print_aligned_short(fp, "atoms");
1013 print_aligned(fp, "theta_fit");
1014 fprintf(fp, " ...");
1026 /* Open torque output file and write some information about it's structure.
1027 * Call on master only */
1028 static FILE *open_torque_out(const char *fn, t_rot *rot, const output_env_t oenv)
1035 if (rot->enfrot->Flags & MD_APPENDFILES)
1037 fp = gmx_fio_fopen(fn,"a");
1041 fp = open_output_file(fn, rot->nstsout,"torques");
1043 for (g=0; g<rot->ngrp; g++)
1045 rotg = &rot->grp[g];
1048 fprintf(fp, "# Rotation group %d (%s), slab distance %f nm.\n", g, erotg_names[rotg->eType], rotg->slab_dist);
1049 fprintf(fp, "# The scalar tau is the torque [kJ/mol] in the direction of the rotation vector.\n");
1050 fprintf(fp, "# To obtain the vectorial torque, multiply tau with\n");
1051 fprintf(fp, "# rot_vec%d %10.3e %10.3e %10.3e\n", g, rotg->vec[XX], rotg->vec[YY], rotg->vec[ZZ]);
1055 fprintf(fp, "# Legend for the following data columns: (tau=torque for that slab):\n");
1057 print_aligned_short(fp, "t");
1058 print_aligned_short(fp, "grp");
1059 print_aligned_short(fp, "slab");
1060 print_aligned(fp, "tau");
1061 print_aligned_short(fp, "slab");
1062 print_aligned(fp, "tau");
1063 fprintf(fp, " ...\n");
1071 static void swap_val(double* vec, int i, int j)
1073 double tmp = vec[j];
1081 static void swap_col(double **mat, int i, int j)
1083 double tmp[3] = {mat[0][j], mat[1][j], mat[2][j]};
1086 mat[0][j]=mat[0][i];
1087 mat[1][j]=mat[1][i];
1088 mat[2][j]=mat[2][i];
1096 /* Eigenvectors are stored in columns of eigen_vec */
1097 static void diagonalize_symmetric(
1105 jacobi(matrix,3,eigenval,eigen_vec,&n_rot);
1107 /* sort in ascending order */
1108 if (eigenval[0] > eigenval[1])
1110 swap_val(eigenval, 0, 1);
1111 swap_col(eigen_vec, 0, 1);
1113 if (eigenval[1] > eigenval[2])
1115 swap_val(eigenval, 1, 2);
1116 swap_col(eigen_vec, 1, 2);
1118 if (eigenval[0] > eigenval[1])
1120 swap_val(eigenval, 0, 1);
1121 swap_col(eigen_vec, 0, 1);
1126 static void align_with_z(
1127 rvec* s, /* Structure to align */
1132 rvec zet = {0.0, 0.0, 1.0};
1133 rvec rot_axis={0.0, 0.0, 0.0};
1134 rvec *rotated_str=NULL;
1140 snew(rotated_str, natoms);
1142 /* Normalize the axis */
1143 ooanorm = 1.0/norm(axis);
1144 svmul(ooanorm, axis, axis);
1146 /* Calculate the angle for the fitting procedure */
1147 cprod(axis, zet, rot_axis);
1148 angle = acos(axis[2]);
1152 /* Calculate the rotation matrix */
1153 calc_rotmat(rot_axis, angle*180.0/M_PI, rotmat);
1155 /* Apply the rotation matrix to s */
1156 for (i=0; i<natoms; i++)
1162 rotated_str[i][j] += rotmat[j][k]*s[i][k];
1167 /* Rewrite the rotated structure to s */
1168 for(i=0; i<natoms; i++)
1172 s[i][j]=rotated_str[i][j];
1180 static void calc_correl_matrix(rvec* Xstr, rvec* Ystr, double** Rmat, int natoms)
1191 for (k=0; k<natoms; k++)
1192 Rmat[i][j] += Ystr[k][i] * Xstr[k][j];
1196 static void weigh_coords(rvec* str, real* weight, int natoms)
1201 for(i=0; i<natoms; i++)
1204 str[i][j] *= sqrt(weight[i]);
1209 static real opt_angle_analytic(
1222 double **Rmat, **RtR, **eigvec;
1224 double V[3][3], WS[3][3];
1225 double rot_matrix[3][3];
1229 /* Do not change the original coordinates */
1230 snew(ref_s_1, natoms);
1231 snew(act_s_1, natoms);
1232 for(i=0; i<natoms; i++)
1234 copy_rvec(ref_s[i], ref_s_1[i]);
1235 copy_rvec(act_s[i], act_s_1[i]);
1238 /* Translate the structures to the origin */
1239 shift[XX] = -ref_com[XX];
1240 shift[YY] = -ref_com[YY];
1241 shift[ZZ] = -ref_com[ZZ];
1242 translate_x(ref_s_1, natoms, shift);
1244 shift[XX] = -act_com[XX];
1245 shift[YY] = -act_com[YY];
1246 shift[ZZ] = -act_com[ZZ];
1247 translate_x(act_s_1, natoms, shift);
1249 /* Align rotation axis with z */
1250 align_with_z(ref_s_1, natoms, axis);
1251 align_with_z(act_s_1, natoms, axis);
1253 /* Correlation matrix */
1254 Rmat = allocate_square_matrix(3);
1256 for (i=0; i<natoms; i++)
1262 /* Weight positions with sqrt(weight) */
1265 weigh_coords(ref_s_1, weight, natoms);
1266 weigh_coords(act_s_1, weight, natoms);
1269 /* Calculate correlation matrices R=YXt (X=ref_s; Y=act_s) */
1270 calc_correl_matrix(ref_s_1, act_s_1, Rmat, natoms);
1273 RtR = allocate_square_matrix(3);
1280 RtR[i][j] += Rmat[k][i] * Rmat[k][j];
1284 /* Diagonalize RtR */
1289 diagonalize_symmetric(RtR, eigvec, eigval);
1290 swap_col(eigvec,0,1);
1291 swap_col(eigvec,1,2);
1292 swap_val(eigval,0,1);
1293 swap_val(eigval,1,2);
1307 WS[i][j] = eigvec[i][j] / sqrt(eigval[j]);
1315 V[i][j] += Rmat[i][k]*WS[k][j];
1319 free_square_matrix(Rmat, 3);
1321 /* Calculate optimal rotation matrix */
1324 rot_matrix[i][j] = 0.0;
1331 rot_matrix[i][j] += eigvec[i][k]*V[j][k];
1335 rot_matrix[2][2] = 1.0;
1337 /* In some cases abs(rot_matrix[0][0]) can be slighly larger
1338 * than unity due to numerical inacurracies. To be able to calculate
1339 * the acos function, we put these values back in range. */
1340 if (rot_matrix[0][0] > 1.0)
1342 rot_matrix[0][0] = 1.0;
1344 else if (rot_matrix[0][0] < -1.0)
1346 rot_matrix[0][0] = -1.0;
1349 /* Determine the optimal rotation angle: */
1350 opt_angle = (-1.0)*acos(rot_matrix[0][0])*180.0/M_PI;
1351 if (rot_matrix[0][1] < 0.0)
1352 opt_angle = (-1.0)*opt_angle;
1354 /* Give back some memory */
1355 free_square_matrix(RtR, 3);
1362 return (real) opt_angle;
1366 /* Determine angle of the group by RMSD fit to the reference */
1367 /* Not parallelized, call this routine only on the master */
1368 static real flex_fit_angle(t_rotgrp *rotg)
1371 rvec *fitcoords=NULL;
1372 rvec center; /* Center of positions passed to the fit routine */
1373 real fitangle; /* Angle of the rotation group derived by fitting */
1376 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
1379 erg=rotg->enfrotgrp;
1381 /* Get the center of the rotation group.
1382 * Note, again, erg->xc has been sorted in do_flexible */
1383 get_center(erg->xc, erg->mc_sorted, rotg->nat, center);
1385 /* === Determine the optimal fit angle for the rotation group === */
1386 if (rotg->eFittype == erotgFitNORM)
1388 /* Normalize every position to it's reference length */
1389 for (i=0; i<rotg->nat; i++)
1391 /* Put the center of the positions into the origin */
1392 rvec_sub(erg->xc[i], center, coord);
1393 /* Determine the scaling factor for the length: */
1394 scal = erg->xc_ref_length[erg->xc_sortind[i]] / norm(coord);
1395 /* Get position, multiply with the scaling factor and save */
1396 svmul(scal, coord, erg->xc_norm[i]);
1398 fitcoords = erg->xc_norm;
1402 fitcoords = erg->xc;
1404 /* From the point of view of the current positions, the reference has rotated
1405 * backwards. Since we output the angle relative to the fixed reference,
1406 * we need the minus sign. */
1407 fitangle = -opt_angle_analytic(erg->xc_ref_sorted, fitcoords, erg->mc_sorted,
1408 rotg->nat, erg->xc_ref_center, center, rotg->vec);
1414 /* Determine actual angle of each slab by RMSD fit to the reference */
1415 /* Not parallelized, call this routine only on the master */
1416 static void flex_fit_angle_perslab(
1423 int i,l,n,islab,ind;
1425 rvec act_center; /* Center of actual positions that are passed to the fit routine */
1426 rvec ref_center; /* Same for the reference positions */
1427 real fitangle; /* Angle of a slab derived from an RMSD fit to
1428 * the reference structure at t=0 */
1430 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
1431 real OOm_av; /* 1/average_mass of a rotation group atom */
1432 real m_rel; /* Relative mass of a rotation group atom */
1435 erg=rotg->enfrotgrp;
1437 /* Average mass of a rotation group atom: */
1438 OOm_av = erg->invmass*rotg->nat;
1440 /**********************************/
1441 /* First collect the data we need */
1442 /**********************************/
1444 /* Collect the data for the individual slabs */
1445 for (n = erg->slab_first; n <= erg->slab_last; n++)
1447 islab = n - erg->slab_first; /* slab index */
1448 sd = &(rotg->enfrotgrp->slab_data[islab]);
1449 sd->nat = erg->lastatom[islab]-erg->firstatom[islab]+1;
1452 /* Loop over the relevant atoms in the slab */
1453 for (l=erg->firstatom[islab]; l<=erg->lastatom[islab]; l++)
1455 /* Current position of this atom: x[ii][XX/YY/ZZ] */
1456 copy_rvec(erg->xc[l], curr_x);
1458 /* The (unrotated) reference position of this atom is copied to ref_x.
1459 * Beware, the xc coords have been sorted in do_flexible */
1460 copy_rvec(erg->xc_ref_sorted[l], ref_x);
1462 /* Save data for doing angular RMSD fit later */
1463 /* Save the current atom position */
1464 copy_rvec(curr_x, sd->x[ind]);
1465 /* Save the corresponding reference position */
1466 copy_rvec(ref_x , sd->ref[ind]);
1468 /* Maybe also mass-weighting was requested. If yes, additionally
1469 * multiply the weights with the relative mass of the atom. If not,
1470 * multiply with unity. */
1471 m_rel = erg->mc_sorted[l]*OOm_av;
1473 /* Save the weight for this atom in this slab */
1474 sd->weight[ind] = gaussian_weight(curr_x, rotg, n) * m_rel;
1476 /* Next atom in this slab */
1481 /******************************/
1482 /* Now do the fit calculation */
1483 /******************************/
1485 fprintf(fp, "%12.3e%6d%12.3f", t, g, degangle);
1487 /* === Now do RMSD fitting for each slab === */
1488 /* We require at least SLAB_MIN_ATOMS in a slab, such that the fit makes sense. */
1489 #define SLAB_MIN_ATOMS 4
1491 for (n = erg->slab_first; n <= erg->slab_last; n++)
1493 islab = n - erg->slab_first; /* slab index */
1494 sd = &(rotg->enfrotgrp->slab_data[islab]);
1495 if (sd->nat >= SLAB_MIN_ATOMS)
1497 /* Get the center of the slabs reference and current positions */
1498 get_center(sd->ref, sd->weight, sd->nat, ref_center);
1499 get_center(sd->x , sd->weight, sd->nat, act_center);
1500 if (rotg->eFittype == erotgFitNORM)
1502 /* Normalize every position to it's reference length
1503 * prior to performing the fit */
1504 for (i=0; i<sd->nat;i++) /* Center */
1506 rvec_dec(sd->ref[i], ref_center);
1507 rvec_dec(sd->x[i] , act_center);
1508 /* Normalize x_i such that it gets the same length as ref_i */
1509 svmul( norm(sd->ref[i])/norm(sd->x[i]), sd->x[i], sd->x[i] );
1511 /* We already subtracted the centers */
1512 clear_rvec(ref_center);
1513 clear_rvec(act_center);
1515 fitangle = -opt_angle_analytic(sd->ref, sd->x, sd->weight, sd->nat,
1516 ref_center, act_center, rotg->vec);
1517 fprintf(fp, "%6d%6d%12.3f", n, sd->nat, fitangle);
1522 #undef SLAB_MIN_ATOMS
1526 /* Shift x with is */
1527 static gmx_inline void shift_single_coord(matrix box, rvec x, const ivec is)
1538 x[XX] += tx*box[XX][XX]+ty*box[YY][XX]+tz*box[ZZ][XX];
1539 x[YY] += ty*box[YY][YY]+tz*box[ZZ][YY];
1540 x[ZZ] += tz*box[ZZ][ZZ];
1543 x[XX] += tx*box[XX][XX];
1544 x[YY] += ty*box[YY][YY];
1545 x[ZZ] += tz*box[ZZ][ZZ];
1550 /* Determine the 'home' slab of this atom which is the
1551 * slab with the highest Gaussian weight of all */
1552 #define round(a) (int)(a+0.5)
1553 static gmx_inline int get_homeslab(
1554 rvec curr_x, /* The position for which the home slab shall be determined */
1555 rvec rotvec, /* The rotation vector */
1556 real slabdist) /* The slab distance */
1561 /* The distance of the atom to the coordinate center (where the
1562 * slab with index 0) is */
1563 dist = iprod(rotvec, curr_x);
1565 return round(dist / slabdist);
1569 /* For a local atom determine the relevant slabs, i.e. slabs in
1570 * which the gaussian is larger than min_gaussian
1572 static int get_single_atom_gaussians(
1580 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
1583 erg=rotg->enfrotgrp;
1585 /* Determine the 'home' slab of this atom: */
1586 homeslab = get_homeslab(curr_x, rotg->vec, rotg->slab_dist);
1588 /* First determine the weight in the atoms home slab: */
1589 g = gaussian_weight(curr_x, rotg, homeslab);
1591 erg->gn_atom[count] = g;
1592 erg->gn_slabind[count] = homeslab;
1596 /* Determine the max slab */
1598 while (g > rotg->min_gaussian)
1601 g = gaussian_weight(curr_x, rotg, slab);
1602 erg->gn_slabind[count]=slab;
1603 erg->gn_atom[count]=g;
1608 /* Determine the max slab */
1613 g = gaussian_weight(curr_x, rotg, slab);
1614 erg->gn_slabind[count]=slab;
1615 erg->gn_atom[count]=g;
1618 while (g > rotg->min_gaussian);
1625 static void flex2_precalc_inner_sum(t_rotgrp *rotg, t_commrec *cr)
1628 rvec xi; /* positions in the i-sum */
1629 rvec xcn, ycn; /* the current and the reference slab centers */
1632 rvec rin; /* Helper variables */
1635 real OOpsii,OOpsiistar;
1636 real sin_rin; /* s_ii.r_ii */
1637 rvec s_in,tmpvec,tmpvec2;
1638 real mi,wi; /* Mass-weighting of the positions */
1640 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
1643 erg=rotg->enfrotgrp;
1644 N_M = rotg->nat * erg->invmass;
1646 /* Loop over all slabs that contain something */
1647 for (n=erg->slab_first; n <= erg->slab_last; n++)
1649 islab = n - erg->slab_first; /* slab index */
1651 /* The current center of this slab is saved in xcn: */
1652 copy_rvec(erg->slab_center[islab], xcn);
1653 /* ... and the reference center in ycn: */
1654 copy_rvec(erg->slab_center_ref[islab+erg->slab_buffer], ycn);
1656 /*** D. Calculate the whole inner sum used for second and third sum */
1657 /* For slab n, we need to loop over all atoms i again. Since we sorted
1658 * the atoms with respect to the rotation vector, we know that it is sufficient
1659 * to calculate from firstatom to lastatom only. All other contributions will
1661 clear_rvec(innersumvec);
1662 for (i = erg->firstatom[islab]; i <= erg->lastatom[islab]; i++)
1664 /* Coordinate xi of this atom */
1665 copy_rvec(erg->xc[i],xi);
1668 gaussian_xi = gaussian_weight(xi,rotg,n);
1669 mi = erg->mc_sorted[i]; /* need the sorted mass here */
1673 copy_rvec(erg->xc_ref_sorted[i],yi0); /* Reference position yi0 */
1674 rvec_sub(yi0, ycn, tmpvec2); /* tmpvec2 = yi0 - ycn */
1675 mvmul(erg->rotmat, tmpvec2, rin); /* rin = Omega.(yi0 - ycn) */
1677 /* Calculate psi_i* and sin */
1678 rvec_sub(xi, xcn, tmpvec2); /* tmpvec2 = xi - xcn */
1679 cprod(rotg->vec, tmpvec2, tmpvec); /* tmpvec = v x (xi - xcn) */
1680 OOpsiistar = norm2(tmpvec)+rotg->eps; /* OOpsii* = 1/psii* = |v x (xi-xcn)|^2 + eps */
1681 OOpsii = norm(tmpvec); /* OOpsii = 1 / psii = |v x (xi - xcn)| */
1683 /* v x (xi - xcn) */
1684 unitv(tmpvec, s_in); /* sin = ---------------- */
1685 /* |v x (xi - xcn)| */
1687 sin_rin=iprod(s_in,rin); /* sin_rin = sin . rin */
1689 /* Now the whole sum */
1690 fac = OOpsii/OOpsiistar;
1691 svmul(fac, rin, tmpvec);
1692 fac2 = fac*fac*OOpsii;
1693 svmul(fac2*sin_rin, s_in, tmpvec2);
1694 rvec_dec(tmpvec, tmpvec2);
1696 svmul(wi*gaussian_xi*sin_rin, tmpvec, tmpvec2);
1698 rvec_inc(innersumvec,tmpvec2);
1699 } /* now we have the inner sum, used both for sum2 and sum3 */
1701 /* Save it to be used in do_flex2_lowlevel */
1702 copy_rvec(innersumvec, erg->slab_innersumvec[islab]);
1703 } /* END of loop over slabs */
1707 static void flex_precalc_inner_sum(t_rotgrp *rotg, t_commrec *cr)
1710 rvec xi; /* position */
1711 rvec xcn, ycn; /* the current and the reference slab centers */
1712 rvec qin,rin; /* q_i^n and r_i^n */
1715 rvec innersumvec; /* Inner part of sum_n2 */
1716 real gaussian_xi; /* Gaussian weight gn(xi) */
1717 real mi,wi; /* Mass-weighting of the positions */
1720 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
1723 erg=rotg->enfrotgrp;
1724 N_M = rotg->nat * erg->invmass;
1726 /* Loop over all slabs that contain something */
1727 for (n=erg->slab_first; n <= erg->slab_last; n++)
1729 islab = n - erg->slab_first; /* slab index */
1731 /* The current center of this slab is saved in xcn: */
1732 copy_rvec(erg->slab_center[islab], xcn);
1733 /* ... and the reference center in ycn: */
1734 copy_rvec(erg->slab_center_ref[islab+erg->slab_buffer], ycn);
1736 /* For slab n, we need to loop over all atoms i again. Since we sorted
1737 * the atoms with respect to the rotation vector, we know that it is sufficient
1738 * to calculate from firstatom to lastatom only. All other contributions will
1740 clear_rvec(innersumvec);
1741 for (i=erg->firstatom[islab]; i<=erg->lastatom[islab]; i++)
1743 /* Coordinate xi of this atom */
1744 copy_rvec(erg->xc[i],xi);
1747 gaussian_xi = gaussian_weight(xi,rotg,n);
1748 mi = erg->mc_sorted[i]; /* need the sorted mass here */
1751 /* Calculate rin and qin */
1752 rvec_sub(erg->xc_ref_sorted[i], ycn, tmpvec); /* tmpvec = yi0-ycn */
1753 mvmul(erg->rotmat, tmpvec, rin); /* rin = Omega.(yi0 - ycn) */
1754 cprod(rotg->vec, rin, tmpvec); /* tmpvec = v x Omega*(yi0-ycn) */
1756 /* v x Omega*(yi0-ycn) */
1757 unitv(tmpvec, qin); /* qin = --------------------- */
1758 /* |v x Omega*(yi0-ycn)| */
1761 rvec_sub(xi, xcn, tmpvec); /* tmpvec = xi-xcn */
1762 bin = iprod(qin, tmpvec); /* bin = qin*(xi-xcn) */
1764 svmul(wi*gaussian_xi*bin, qin, tmpvec);
1766 /* Add this contribution to the inner sum: */
1767 rvec_add(innersumvec, tmpvec, innersumvec);
1768 } /* now we have the inner sum vector S^n for this slab */
1769 /* Save it to be used in do_flex_lowlevel */
1770 copy_rvec(innersumvec, erg->slab_innersumvec[islab]);
1775 static real do_flex2_lowlevel(
1777 real sigma, /* The Gaussian width sigma */
1779 gmx_bool bOutstepRot,
1780 gmx_bool bOutstepSlab,
1784 int count,ic,ii,j,m,n,islab,iigrp,ifit;
1785 rvec xj; /* position in the i-sum */
1786 rvec yj0; /* the reference position in the j-sum */
1787 rvec xcn, ycn; /* the current and the reference slab centers */
1788 real V; /* This node's part of the rotation pot. energy */
1789 real gaussian_xj; /* Gaussian weight */
1792 real numerator,fit_numerator;
1793 rvec rjn,fit_rjn; /* Helper variables */
1796 real OOpsij,OOpsijstar;
1797 real OOsigma2; /* 1/(sigma^2) */
1800 rvec sjn,tmpvec,tmpvec2,yj0_ycn;
1801 rvec sum1vec_part,sum1vec,sum2vec_part,sum2vec,sum3vec,sum4vec,innersumvec;
1803 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
1804 real mj,wj; /* Mass-weighting of the positions */
1806 real Wjn; /* g_n(x_j) m_j / Mjn */
1807 gmx_bool bCalcPotFit;
1809 /* To calculate the torque per slab */
1810 rvec slab_force; /* Single force from slab n on one atom */
1811 rvec slab_sum1vec_part;
1812 real slab_sum3part,slab_sum4part;
1813 rvec slab_sum1vec, slab_sum2vec, slab_sum3vec, slab_sum4vec;
1816 erg=rotg->enfrotgrp;
1818 /* Pre-calculate the inner sums, so that we do not have to calculate
1819 * them again for every atom */
1820 flex2_precalc_inner_sum(rotg, cr);
1822 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT==rotg->eFittype);
1824 /********************************************************/
1825 /* Main loop over all local atoms of the rotation group */
1826 /********************************************************/
1827 N_M = rotg->nat * erg->invmass;
1829 OOsigma2 = 1.0 / (sigma*sigma);
1830 for (j=0; j<erg->nat_loc; j++)
1832 /* Local index of a rotation group atom */
1833 ii = erg->ind_loc[j];
1834 /* Position of this atom in the collective array */
1835 iigrp = erg->xc_ref_ind[j];
1836 /* Mass-weighting */
1837 mj = erg->mc[iigrp]; /* need the unsorted mass here */
1840 /* Current position of this atom: x[ii][XX/YY/ZZ]
1841 * Note that erg->xc_center contains the center of mass in case the flex2-t
1842 * potential was chosen. For the flex2 potential erg->xc_center must be
1844 rvec_sub(x[ii], erg->xc_center, xj);
1846 /* Shift this atom such that it is near its reference */
1847 shift_single_coord(box, xj, erg->xc_shifts[iigrp]);
1849 /* Determine the slabs to loop over, i.e. the ones with contributions
1850 * larger than min_gaussian */
1851 count = get_single_atom_gaussians(xj, cr, rotg);
1853 clear_rvec(sum1vec_part);
1854 clear_rvec(sum2vec_part);
1857 /* Loop over the relevant slabs for this atom */
1858 for (ic=0; ic < count; ic++)
1860 n = erg->gn_slabind[ic];
1862 /* Get the precomputed Gaussian value of curr_slab for curr_x */
1863 gaussian_xj = erg->gn_atom[ic];
1865 islab = n - erg->slab_first; /* slab index */
1867 /* The (unrotated) reference position of this atom is copied to yj0: */
1868 copy_rvec(rotg->x_ref[iigrp], yj0);
1870 beta = calc_beta(xj, rotg,n);
1872 /* The current center of this slab is saved in xcn: */
1873 copy_rvec(erg->slab_center[islab], xcn);
1874 /* ... and the reference center in ycn: */
1875 copy_rvec(erg->slab_center_ref[islab+erg->slab_buffer], ycn);
1877 rvec_sub(yj0, ycn, yj0_ycn); /* yj0_ycn = yj0 - ycn */
1880 mvmul(erg->rotmat, yj0_ycn, rjn); /* rjn = Omega.(yj0 - ycn) */
1882 /* Subtract the slab center from xj */
1883 rvec_sub(xj, xcn, tmpvec2); /* tmpvec2 = xj - xcn */
1886 cprod(rotg->vec, tmpvec2, tmpvec); /* tmpvec = v x (xj - xcn) */
1888 OOpsijstar = norm2(tmpvec)+rotg->eps; /* OOpsij* = 1/psij* = |v x (xj-xcn)|^2 + eps */
1890 numerator = sqr(iprod(tmpvec, rjn));
1892 /*********************************/
1893 /* Add to the rotation potential */
1894 /*********************************/
1895 V += 0.5*rotg->k*wj*gaussian_xj*numerator/OOpsijstar;
1897 /* If requested, also calculate the potential for a set of angles
1898 * near the current reference angle */
1901 for (ifit = 0; ifit < rotg->PotAngle_nstep; ifit++)
1903 mvmul(erg->PotAngleFit->rotmat[ifit], yj0_ycn, fit_rjn);
1904 fit_numerator = sqr(iprod(tmpvec, fit_rjn));
1905 erg->PotAngleFit->V[ifit] += 0.5*rotg->k*wj*gaussian_xj*fit_numerator/OOpsijstar;
1909 /*************************************/
1910 /* Now calculate the force on atom j */
1911 /*************************************/
1913 OOpsij = norm(tmpvec); /* OOpsij = 1 / psij = |v x (xj - xcn)| */
1915 /* v x (xj - xcn) */
1916 unitv(tmpvec, sjn); /* sjn = ---------------- */
1917 /* |v x (xj - xcn)| */
1919 sjn_rjn=iprod(sjn,rjn); /* sjn_rjn = sjn . rjn */
1922 /*** A. Calculate the first of the four sum terms: ****************/
1923 fac = OOpsij/OOpsijstar;
1924 svmul(fac, rjn, tmpvec);
1925 fac2 = fac*fac*OOpsij;
1926 svmul(fac2*sjn_rjn, sjn, tmpvec2);
1927 rvec_dec(tmpvec, tmpvec2);
1928 fac2 = wj*gaussian_xj; /* also needed for sum4 */
1929 svmul(fac2*sjn_rjn, tmpvec, slab_sum1vec_part);
1930 /********************/
1931 /*** Add to sum1: ***/
1932 /********************/
1933 rvec_inc(sum1vec_part, slab_sum1vec_part); /* sum1 still needs to vector multiplied with v */
1935 /*** B. Calculate the forth of the four sum terms: ****************/
1936 betasigpsi = beta*OOsigma2*OOpsij; /* this is also needed for sum3 */
1937 /********************/
1938 /*** Add to sum4: ***/
1939 /********************/
1940 slab_sum4part = fac2*betasigpsi*fac*sjn_rjn*sjn_rjn; /* Note that fac is still valid from above */
1941 sum4 += slab_sum4part;
1943 /*** C. Calculate Wjn for second and third sum */
1944 /* Note that we can safely divide by slab_weights since we check in
1945 * get_slab_centers that it is non-zero. */
1946 Wjn = gaussian_xj*mj/erg->slab_weights[islab];
1948 /* We already have precalculated the inner sum for slab n */
1949 copy_rvec(erg->slab_innersumvec[islab], innersumvec);
1951 /* Weigh the inner sum vector with Wjn */
1952 svmul(Wjn, innersumvec, innersumvec);
1954 /*** E. Calculate the second of the four sum terms: */
1955 /********************/
1956 /*** Add to sum2: ***/
1957 /********************/
1958 rvec_inc(sum2vec_part, innersumvec); /* sum2 still needs to be vector crossproduct'ed with v */
1960 /*** F. Calculate the third of the four sum terms: */
1961 slab_sum3part = betasigpsi * iprod(sjn, innersumvec);
1962 sum3 += slab_sum3part; /* still needs to be multiplied with v */
1964 /*** G. Calculate the torque on the local slab's axis: */
1968 cprod(slab_sum1vec_part, rotg->vec, slab_sum1vec);
1970 cprod(innersumvec, rotg->vec, slab_sum2vec);
1972 svmul(slab_sum3part, rotg->vec, slab_sum3vec);
1974 svmul(slab_sum4part, rotg->vec, slab_sum4vec);
1976 /* The force on atom ii from slab n only: */
1977 for (m=0; m<DIM; m++)
1978 slab_force[m] = rotg->k * (-slab_sum1vec[m] + slab_sum2vec[m] - slab_sum3vec[m] + 0.5*slab_sum4vec[m]);
1980 erg->slab_torque_v[islab] += torque(rotg->vec, slab_force, xj, xcn);
1982 } /* END of loop over slabs */
1984 /* Construct the four individual parts of the vector sum: */
1985 cprod(sum1vec_part, rotg->vec, sum1vec); /* sum1vec = { } x v */
1986 cprod(sum2vec_part, rotg->vec, sum2vec); /* sum2vec = { } x v */
1987 svmul(sum3, rotg->vec, sum3vec); /* sum3vec = { } . v */
1988 svmul(sum4, rotg->vec, sum4vec); /* sum4vec = { } . v */
1990 /* Store the additional force so that it can be added to the force
1991 * array after the normal forces have been evaluated */
1992 for (m=0; m<DIM; m++)
1993 erg->f_rot_loc[j][m] = rotg->k * (-sum1vec[m] + sum2vec[m] - sum3vec[m] + 0.5*sum4vec[m]);
1996 fprintf(stderr, "sum1: %15.8f %15.8f %15.8f\n", -rotg->k*sum1vec[XX], -rotg->k*sum1vec[YY], -rotg->k*sum1vec[ZZ]);
1997 fprintf(stderr, "sum2: %15.8f %15.8f %15.8f\n", rotg->k*sum2vec[XX], rotg->k*sum2vec[YY], rotg->k*sum2vec[ZZ]);
1998 fprintf(stderr, "sum3: %15.8f %15.8f %15.8f\n", -rotg->k*sum3vec[XX], -rotg->k*sum3vec[YY], -rotg->k*sum3vec[ZZ]);
1999 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]);
2004 } /* END of loop over local atoms */
2010 static real do_flex_lowlevel(
2012 real sigma, /* The Gaussian width sigma */
2014 gmx_bool bOutstepRot,
2015 gmx_bool bOutstepSlab,
2019 int count,ic,ifit,ii,j,m,n,islab,iigrp;
2020 rvec xj,yj0; /* current and reference position */
2021 rvec xcn, ycn; /* the current and the reference slab centers */
2022 rvec yj0_ycn; /* yj0 - ycn */
2023 rvec xj_xcn; /* xj - xcn */
2024 rvec qjn,fit_qjn; /* q_i^n */
2025 rvec sum_n1,sum_n2; /* Two contributions to the rotation force */
2026 rvec innersumvec; /* Inner part of sum_n2 */
2028 rvec force_n; /* Single force from slab n on one atom */
2029 rvec force_n1,force_n2; /* First and second part of force_n */
2030 rvec tmpvec,tmpvec2,tmp_f; /* Helper variables */
2031 real V; /* The rotation potential energy */
2032 real OOsigma2; /* 1/(sigma^2) */
2033 real beta; /* beta_n(xj) */
2034 real bjn, fit_bjn; /* b_j^n */
2035 real gaussian_xj; /* Gaussian weight gn(xj) */
2036 real betan_xj_sigma2;
2037 real mj,wj; /* Mass-weighting of the positions */
2039 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
2040 gmx_bool bCalcPotFit;
2043 erg=rotg->enfrotgrp;
2045 /* Pre-calculate the inner sums, so that we do not have to calculate
2046 * them again for every atom */
2047 flex_precalc_inner_sum(rotg, cr);
2049 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT==rotg->eFittype);
2051 /********************************************************/
2052 /* Main loop over all local atoms of the rotation group */
2053 /********************************************************/
2054 OOsigma2 = 1.0/(sigma*sigma);
2055 N_M = rotg->nat * erg->invmass;
2057 for (j=0; j<erg->nat_loc; j++)
2059 /* Local index of a rotation group atom */
2060 ii = erg->ind_loc[j];
2061 /* Position of this atom in the collective array */
2062 iigrp = erg->xc_ref_ind[j];
2063 /* Mass-weighting */
2064 mj = erg->mc[iigrp]; /* need the unsorted mass here */
2067 /* Current position of this atom: x[ii][XX/YY/ZZ]
2068 * Note that erg->xc_center contains the center of mass in case the flex-t
2069 * potential was chosen. For the flex potential erg->xc_center must be
2071 rvec_sub(x[ii], erg->xc_center, xj);
2073 /* Shift this atom such that it is near its reference */
2074 shift_single_coord(box, xj, erg->xc_shifts[iigrp]);
2076 /* Determine the slabs to loop over, i.e. the ones with contributions
2077 * larger than min_gaussian */
2078 count = get_single_atom_gaussians(xj, cr, rotg);
2083 /* Loop over the relevant slabs for this atom */
2084 for (ic=0; ic < count; ic++)
2086 n = erg->gn_slabind[ic];
2088 /* Get the precomputed Gaussian for xj in slab n */
2089 gaussian_xj = erg->gn_atom[ic];
2091 islab = n - erg->slab_first; /* slab index */
2093 /* The (unrotated) reference position of this atom is saved in yj0: */
2094 copy_rvec(rotg->x_ref[iigrp], yj0);
2096 beta = calc_beta(xj, rotg, n);
2098 /* The current center of this slab is saved in xcn: */
2099 copy_rvec(erg->slab_center[islab], xcn);
2100 /* ... and the reference center in ycn: */
2101 copy_rvec(erg->slab_center_ref[islab+erg->slab_buffer], ycn);
2103 rvec_sub(yj0, ycn, yj0_ycn); /* yj0_ycn = yj0 - ycn */
2106 mvmul(erg->rotmat, yj0_ycn, tmpvec2); /* tmpvec2= Omega.(yj0-ycn) */
2108 /* Subtract the slab center from xj */
2109 rvec_sub(xj, xcn, xj_xcn); /* xj_xcn = xj - xcn */
2112 cprod(rotg->vec, tmpvec2, tmpvec); /* tmpvec= v x Omega.(yj0-ycn) */
2114 /* v x Omega.(yj0-ycn) */
2115 unitv(tmpvec,qjn); /* qjn = --------------------- */
2116 /* |v x Omega.(yj0-ycn)| */
2118 bjn = iprod(qjn, xj_xcn); /* bjn = qjn * (xj - xcn) */
2120 /*********************************/
2121 /* Add to the rotation potential */
2122 /*********************************/
2123 V += 0.5*rotg->k*wj*gaussian_xj*sqr(bjn);
2125 /* If requested, also calculate the potential for a set of angles
2126 * near the current reference angle */
2129 for (ifit = 0; ifit < rotg->PotAngle_nstep; ifit++)
2131 /* As above calculate Omega.(yj0-ycn), now for the other angles */
2132 mvmul(erg->PotAngleFit->rotmat[ifit], yj0_ycn, tmpvec2); /* tmpvec2= Omega.(yj0-ycn) */
2133 /* As above calculate qjn */
2134 cprod(rotg->vec, tmpvec2, tmpvec); /* tmpvec= v x Omega.(yj0-ycn) */
2135 /* v x Omega.(yj0-ycn) */
2136 unitv(tmpvec,fit_qjn); /* fit_qjn = --------------------- */
2137 /* |v x Omega.(yj0-ycn)| */
2138 fit_bjn = iprod(fit_qjn, xj_xcn); /* fit_bjn = fit_qjn * (xj - xcn) */
2139 /* Add to the rotation potential for this angle */
2140 erg->PotAngleFit->V[ifit] += 0.5*rotg->k*wj*gaussian_xj*sqr(fit_bjn);
2144 /****************************************************************/
2145 /* sum_n1 will typically be the main contribution to the force: */
2146 /****************************************************************/
2147 betan_xj_sigma2 = beta*OOsigma2; /* beta_n(xj)/sigma^2 */
2149 /* The next lines calculate
2150 * qjn - (bjn*beta(xj)/(2sigma^2))v */
2151 svmul(bjn*0.5*betan_xj_sigma2, rotg->vec, tmpvec2);
2152 rvec_sub(qjn,tmpvec2,tmpvec);
2154 /* Multiply with gn(xj)*bjn: */
2155 svmul(gaussian_xj*bjn,tmpvec,tmpvec2);
2158 rvec_inc(sum_n1,tmpvec2);
2160 /* We already have precalculated the Sn term for slab n */
2161 copy_rvec(erg->slab_innersumvec[islab], s_n);
2163 svmul(betan_xj_sigma2*iprod(s_n, xj_xcn), rotg->vec, tmpvec); /* tmpvec = ---------- s_n (xj-xcn) */
2166 rvec_sub(s_n, tmpvec, innersumvec);
2168 /* We can safely divide by slab_weights since we check in get_slab_centers
2169 * that it is non-zero. */
2170 svmul(gaussian_xj/erg->slab_weights[islab], innersumvec, innersumvec);
2172 rvec_add(sum_n2, innersumvec, sum_n2);
2174 /* Calculate the torque: */
2177 /* The force on atom ii from slab n only: */
2178 svmul(-rotg->k*wj, tmpvec2 , force_n1); /* part 1 */
2179 svmul( rotg->k*mj, innersumvec, force_n2); /* part 2 */
2180 rvec_add(force_n1, force_n2, force_n);
2181 erg->slab_torque_v[islab] += torque(rotg->vec, force_n, xj, xcn);
2183 } /* END of loop over slabs */
2185 /* Put both contributions together: */
2186 svmul(wj, sum_n1, sum_n1);
2187 svmul(mj, sum_n2, sum_n2);
2188 rvec_sub(sum_n2,sum_n1,tmp_f); /* F = -grad V */
2190 /* Store the additional force so that it can be added to the force
2191 * array after the normal forces have been evaluated */
2192 for(m=0; m<DIM; m++)
2193 erg->f_rot_loc[j][m] = rotg->k*tmp_f[m];
2197 } /* END of loop over local atoms */
2203 static void print_coordinates(t_commrec *cr, t_rotgrp *rotg, rvec x[], matrix box, int step)
2207 static char buf[STRLEN];
2208 static gmx_bool bFirst=1;
2213 sprintf(buf, "coords%d.txt", cr->nodeid);
2214 fp = fopen(buf, "w");
2218 fprintf(fp, "\nStep %d\n", step);
2219 fprintf(fp, "box: %f %f %f %f %f %f %f %f %f\n",
2220 box[XX][XX], box[XX][YY], box[XX][ZZ],
2221 box[YY][XX], box[YY][YY], box[YY][ZZ],
2222 box[ZZ][XX], box[ZZ][ZZ], box[ZZ][ZZ]);
2223 for (i=0; i<rotg->nat; i++)
2225 fprintf(fp, "%4d %f %f %f\n", i,
2226 erg->xc[i][XX], erg->xc[i][YY], erg->xc[i][ZZ]);
2234 static int projection_compare(const void *a, const void *b)
2236 sort_along_vec_t *xca, *xcb;
2239 xca = (sort_along_vec_t *)a;
2240 xcb = (sort_along_vec_t *)b;
2242 if (xca->xcproj < xcb->xcproj)
2244 else if (xca->xcproj > xcb->xcproj)
2251 static void sort_collective_coordinates(
2252 t_rotgrp *rotg, /* Rotation group */
2253 sort_along_vec_t *data) /* Buffer for sorting the positions */
2256 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
2259 erg=rotg->enfrotgrp;
2261 /* The projection of the position vector on the rotation vector is
2262 * the relevant value for sorting. Fill the 'data' structure */
2263 for (i=0; i<rotg->nat; i++)
2265 data[i].xcproj = iprod(erg->xc[i], rotg->vec); /* sort criterium */
2266 data[i].m = erg->mc[i];
2268 copy_rvec(erg->xc[i] , data[i].x );
2269 copy_rvec(rotg->x_ref[i], data[i].x_ref);
2271 /* Sort the 'data' structure */
2272 gmx_qsort(data, rotg->nat, sizeof(sort_along_vec_t), projection_compare);
2274 /* Copy back the sorted values */
2275 for (i=0; i<rotg->nat; i++)
2277 copy_rvec(data[i].x , erg->xc[i] );
2278 copy_rvec(data[i].x_ref, erg->xc_ref_sorted[i]);
2279 erg->mc_sorted[i] = data[i].m;
2280 erg->xc_sortind[i] = data[i].ind;
2285 /* For each slab, get the first and the last index of the sorted atom
2287 static void get_firstlast_atom_per_slab(t_rotgrp *rotg, t_commrec *cr)
2291 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
2294 erg=rotg->enfrotgrp;
2296 /* Find the first atom that needs to enter the calculation for each slab */
2297 n = erg->slab_first; /* slab */
2298 i = 0; /* start with the first atom */
2301 /* Find the first atom that significantly contributes to this slab */
2302 do /* move forward in position until a large enough beta is found */
2304 beta = calc_beta(erg->xc[i], rotg, n);
2306 } while ((beta < -erg->max_beta) && (i < rotg->nat));
2308 islab = n - erg->slab_first; /* slab index */
2309 erg->firstatom[islab] = i;
2310 /* Proceed to the next slab */
2312 } while (n <= erg->slab_last);
2314 /* Find the last atom for each slab */
2315 n = erg->slab_last; /* start with last slab */
2316 i = rotg->nat-1; /* start with the last atom */
2319 do /* move backward in position until a large enough beta is found */
2321 beta = calc_beta(erg->xc[i], rotg, n);
2323 } while ((beta > erg->max_beta) && (i > -1));
2325 islab = n - erg->slab_first; /* slab index */
2326 erg->lastatom[islab] = i;
2327 /* Proceed to the next slab */
2329 } while (n >= erg->slab_first);
2333 /* Determine the very first and very last slab that needs to be considered
2334 * For the first slab that needs to be considered, we have to find the smallest
2337 * x_first * v - n*Delta_x <= beta_max
2339 * slab index n, slab distance Delta_x, rotation vector v. For the last slab we
2340 * have to find the largest n that obeys
2342 * x_last * v - n*Delta_x >= -beta_max
2345 static gmx_inline int get_first_slab(
2346 t_rotgrp *rotg, /* The rotation group (inputrec data) */
2347 real max_beta, /* The max_beta value, instead of min_gaussian */
2348 rvec firstatom) /* First atom after sorting along the rotation vector v */
2350 /* Find the first slab for the first atom */
2351 return ceil((iprod(firstatom, rotg->vec) - max_beta)/rotg->slab_dist);
2355 static gmx_inline int get_last_slab(
2356 t_rotgrp *rotg, /* The rotation group (inputrec data) */
2357 real max_beta, /* The max_beta value, instead of min_gaussian */
2358 rvec lastatom) /* Last atom along v */
2360 /* Find the last slab for the last atom */
2361 return floor((iprod(lastatom, rotg->vec) + max_beta)/rotg->slab_dist);
2365 static void get_firstlast_slab_check(
2366 t_rotgrp *rotg, /* The rotation group (inputrec data) */
2367 t_gmx_enfrotgrp *erg, /* The rotation group (data only accessible in this file) */
2368 rvec firstatom, /* First atom after sorting along the rotation vector v */
2369 rvec lastatom, /* Last atom along v */
2370 int g, /* The rotation group number */
2373 erg->slab_first = get_first_slab(rotg, erg->max_beta, firstatom);
2374 erg->slab_last = get_last_slab(rotg, erg->max_beta, lastatom);
2376 /* Check whether we have reference data to compare against */
2377 if (erg->slab_first < erg->slab_first_ref)
2378 gmx_fatal(FARGS, "%s No reference data for first slab (n=%d), unable to proceed.",
2379 RotStr, erg->slab_first);
2381 /* Check whether we have reference data to compare against */
2382 if (erg->slab_last > erg->slab_last_ref)
2383 gmx_fatal(FARGS, "%s No reference data for last slab (n=%d), unable to proceed.",
2384 RotStr, erg->slab_last);
2388 /* Enforced rotation with a flexible axis */
2389 static void do_flexible(
2391 gmx_enfrot_t enfrot, /* Other rotation data */
2392 t_rotgrp *rotg, /* The rotation group */
2393 int g, /* Group number */
2394 rvec x[], /* The local positions */
2396 double t, /* Time in picoseconds */
2397 gmx_large_int_t step, /* The time step */
2398 gmx_bool bOutstepRot, /* Output to main rotation output file */
2399 gmx_bool bOutstepSlab) /* Output per-slab data */
2402 real sigma; /* The Gaussian width sigma */
2403 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
2406 erg=rotg->enfrotgrp;
2408 /* Define the sigma value */
2409 sigma = 0.7*rotg->slab_dist;
2411 /* Sort the collective coordinates erg->xc along the rotation vector. This is
2412 * an optimization for the inner loop. */
2413 sort_collective_coordinates(rotg, enfrot->data);
2415 /* Determine the first relevant slab for the first atom and the last
2416 * relevant slab for the last atom */
2417 get_firstlast_slab_check(rotg, erg, erg->xc[0], erg->xc[rotg->nat-1], g, cr);
2419 /* Determine for each slab depending on the min_gaussian cutoff criterium,
2420 * a first and a last atom index inbetween stuff needs to be calculated */
2421 get_firstlast_atom_per_slab(rotg, cr);
2423 /* Determine the gaussian-weighted center of positions for all slabs */
2424 get_slab_centers(rotg,erg->xc,erg->mc_sorted,cr,g,t,enfrot->out_slabs,bOutstepSlab,FALSE);
2426 /* Clear the torque per slab from last time step: */
2427 nslabs = erg->slab_last - erg->slab_first + 1;
2428 for (l=0; l<nslabs; l++)
2429 erg->slab_torque_v[l] = 0.0;
2431 /* Call the rotational forces kernel */
2432 if (rotg->eType == erotgFLEX || rotg->eType == erotgFLEXT)
2433 erg->V = do_flex_lowlevel(rotg, sigma, x, bOutstepRot, bOutstepSlab, box, cr);
2434 else if (rotg->eType == erotgFLEX2 || rotg->eType == erotgFLEX2T)
2435 erg->V = do_flex2_lowlevel(rotg, sigma, x, bOutstepRot, bOutstepSlab, box, cr);
2437 gmx_fatal(FARGS, "Unknown flexible rotation type");
2439 /* Determine angle by RMSD fit to the reference - Let's hope this */
2440 /* only happens once in a while, since this is not parallelized! */
2441 if (MASTER(cr) && (erotgFitPOT != rotg->eFittype) )
2445 /* Fit angle of the whole rotation group */
2446 erg->angle_v = flex_fit_angle(rotg);
2450 /* Fit angle of each slab */
2451 flex_fit_angle_perslab(g, rotg, t, erg->degangle, enfrot->out_angles);
2455 /* Lump together the torques from all slabs: */
2456 erg->torque_v = 0.0;
2457 for (l=0; l<nslabs; l++)
2458 erg->torque_v += erg->slab_torque_v[l];
2462 /* Calculate the angle between reference and actual rotation group atom,
2463 * both projected into a plane perpendicular to the rotation vector: */
2464 static void angle(t_rotgrp *rotg,
2468 real *weight) /* atoms near the rotation axis should count less than atoms far away */
2470 rvec xp, xrp; /* current and reference positions projected on a plane perpendicular to pg->vec */
2474 /* Project x_ref and x into a plane through the origin perpendicular to rot_vec: */
2475 /* Project x_ref: xrp = x_ref - (vec * x_ref) * vec */
2476 svmul(iprod(rotg->vec, x_ref), rotg->vec, dum);
2477 rvec_sub(x_ref, dum, xrp);
2478 /* Project x_act: */
2479 svmul(iprod(rotg->vec, x_act), rotg->vec, dum);
2480 rvec_sub(x_act, dum, xp);
2482 /* Retrieve information about which vector precedes. gmx_angle always
2483 * returns a positive angle. */
2484 cprod(xp, xrp, dum); /* if reference precedes, this is pointing into the same direction as vec */
2486 if (iprod(rotg->vec, dum) >= 0)
2487 *alpha = -gmx_angle(xrp, xp);
2489 *alpha = +gmx_angle(xrp, xp);
2491 /* Also return the weight */
2496 /* Project first vector onto a plane perpendicular to the second vector
2498 * Note that v must be of unit length.
2500 static gmx_inline void project_onto_plane(rvec dr, const rvec v)
2505 svmul(iprod(dr,v),v,tmp); /* tmp = (dr.v)v */
2506 rvec_dec(dr, tmp); /* dr = dr - (dr.v)v */
2510 /* Fixed rotation: The rotation reference group rotates around the v axis. */
2511 /* The atoms of the actual rotation group are attached with imaginary */
2512 /* springs to the reference atoms. */
2513 static void do_fixed(
2515 t_rotgrp *rotg, /* The rotation group */
2516 rvec x[], /* The positions */
2517 matrix box, /* The simulation box */
2518 double t, /* Time in picoseconds */
2519 gmx_large_int_t step, /* The time step */
2520 gmx_bool bOutstepRot, /* Output to main rotation output file */
2521 gmx_bool bOutstepSlab) /* Output per-slab data */
2525 rvec tmp_f; /* Force */
2526 real alpha; /* a single angle between an actual and a reference position */
2527 real weight; /* single weight for a single angle */
2528 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
2529 rvec xi_xc; /* xi - xc */
2530 gmx_bool bCalcPotFit;
2533 /* for mass weighting: */
2534 real wi; /* Mass-weighting of the positions */
2536 real k_wi; /* k times wi */
2541 erg=rotg->enfrotgrp;
2542 bProject = (rotg->eType==erotgPM) || (rotg->eType==erotgPMPF);
2543 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT==rotg->eFittype);
2545 N_M = rotg->nat * erg->invmass;
2547 /* Each process calculates the forces on its local atoms */
2548 for (j=0; j<erg->nat_loc; j++)
2550 /* Calculate (x_i-x_c) resp. (x_i-u) */
2551 rvec_sub(erg->x_loc_pbc[j], erg->xc_center, xi_xc);
2553 /* Calculate Omega*(y_i-y_c)-(x_i-x_c) */
2554 rvec_sub(erg->xr_loc[j], xi_xc, dr);
2557 project_onto_plane(dr, rotg->vec);
2559 /* Mass-weighting */
2560 wi = N_M*erg->m_loc[j];
2562 /* Store the additional force so that it can be added to the force
2563 * array after the normal forces have been evaluated */
2565 for (m=0; m<DIM; m++)
2567 tmp_f[m] = k_wi*dr[m];
2568 erg->f_rot_loc[j][m] = tmp_f[m];
2569 erg->V += 0.5*k_wi*sqr(dr[m]);
2572 /* If requested, also calculate the potential for a set of angles
2573 * near the current reference angle */
2576 for (ifit = 0; ifit < rotg->PotAngle_nstep; ifit++)
2578 /* Index of this rotation group atom with respect to the whole rotation group */
2579 jj = erg->xc_ref_ind[j];
2581 /* Rotate with the alternative angle. Like rotate_local_reference(),
2582 * just for a single local atom */
2583 mvmul(erg->PotAngleFit->rotmat[ifit], rotg->x_ref[jj], fit_xr_loc); /* fit_xr_loc = Omega*(y_i-y_c) */
2585 /* Calculate Omega*(y_i-y_c)-(x_i-x_c) */
2586 rvec_sub(fit_xr_loc, xi_xc, dr);
2589 project_onto_plane(dr, rotg->vec);
2591 /* Add to the rotation potential for this angle: */
2592 erg->PotAngleFit->V[ifit] += 0.5*k_wi*norm2(dr);
2598 /* Add to the torque of this rotation group */
2599 erg->torque_v += torque(rotg->vec, tmp_f, erg->x_loc_pbc[j], erg->xc_center);
2601 /* Calculate the angle between reference and actual rotation group atom. */
2602 angle(rotg, xi_xc, erg->xr_loc[j], &alpha, &weight); /* angle in rad, weighted */
2603 erg->angle_v += alpha * weight;
2604 erg->weight_v += weight;
2606 /* If you want enforced rotation to contribute to the virial,
2607 * activate the following lines:
2610 Add the rotation contribution to the virial
2611 for(j=0; j<DIM; j++)
2613 vir[j][m] += 0.5*f[ii][j]*dr[m];
2619 } /* end of loop over local rotation group atoms */
2623 /* Calculate the radial motion potential and forces */
2624 static void do_radial_motion(
2626 t_rotgrp *rotg, /* The rotation group */
2627 rvec x[], /* The positions */
2628 matrix box, /* The simulation box */
2629 double t, /* Time in picoseconds */
2630 gmx_large_int_t step, /* The time step */
2631 gmx_bool bOutstepRot, /* Output to main rotation output file */
2632 gmx_bool bOutstepSlab) /* Output per-slab data */
2635 rvec tmp_f; /* Force */
2636 real alpha; /* a single angle between an actual and a reference position */
2637 real weight; /* single weight for a single angle */
2638 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
2639 rvec xj_u; /* xj - u */
2640 rvec tmpvec,fit_tmpvec;
2641 real fac,fac2,sum=0.0;
2643 gmx_bool bCalcPotFit;
2645 /* For mass weighting: */
2646 real wj; /* Mass-weighting of the positions */
2650 erg=rotg->enfrotgrp;
2651 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT==rotg->eFittype);
2653 N_M = rotg->nat * erg->invmass;
2655 /* Each process calculates the forces on its local atoms */
2656 for (j=0; j<erg->nat_loc; j++)
2658 /* Calculate (xj-u) */
2659 rvec_sub(erg->x_loc_pbc[j], erg->xc_center, xj_u); /* xj_u = xj-u */
2661 /* Calculate Omega.(yj0-u) */
2662 cprod(rotg->vec, erg->xr_loc[j], tmpvec); /* tmpvec = v x Omega.(yj0-u) */
2664 /* v x Omega.(yj0-u) */
2665 unitv(tmpvec, pj); /* pj = --------------------- */
2666 /* | v x Omega.(yj0-u) | */
2668 fac = iprod(pj, xj_u); /* fac = pj.(xj-u) */
2671 /* Mass-weighting */
2672 wj = N_M*erg->m_loc[j];
2674 /* Store the additional force so that it can be added to the force
2675 * array after the normal forces have been evaluated */
2676 svmul(-rotg->k*wj*fac, pj, tmp_f);
2677 copy_rvec(tmp_f, erg->f_rot_loc[j]);
2680 /* If requested, also calculate the potential for a set of angles
2681 * near the current reference angle */
2684 for (ifit = 0; ifit < rotg->PotAngle_nstep; ifit++)
2686 /* Index of this rotation group atom with respect to the whole rotation group */
2687 jj = erg->xc_ref_ind[j];
2689 /* Rotate with the alternative angle. Like rotate_local_reference(),
2690 * just for a single local atom */
2691 mvmul(erg->PotAngleFit->rotmat[ifit], rotg->x_ref[jj], fit_tmpvec); /* fit_tmpvec = Omega*(yj0-u) */
2693 /* Calculate Omega.(yj0-u) */
2694 cprod(rotg->vec, fit_tmpvec, tmpvec); /* tmpvec = v x Omega.(yj0-u) */
2695 /* v x Omega.(yj0-u) */
2696 unitv(tmpvec, pj); /* pj = --------------------- */
2697 /* | v x Omega.(yj0-u) | */
2699 fac = iprod(pj, xj_u); /* fac = pj.(xj-u) */
2702 /* Add to the rotation potential for this angle: */
2703 erg->PotAngleFit->V[ifit] += 0.5*rotg->k*wj*fac2;
2709 /* Add to the torque of this rotation group */
2710 erg->torque_v += torque(rotg->vec, tmp_f, erg->x_loc_pbc[j], erg->xc_center);
2712 /* Calculate the angle between reference and actual rotation group atom. */
2713 angle(rotg, xj_u, erg->xr_loc[j], &alpha, &weight); /* angle in rad, weighted */
2714 erg->angle_v += alpha * weight;
2715 erg->weight_v += weight;
2720 } /* end of loop over local rotation group atoms */
2721 erg->V = 0.5*rotg->k*sum;
2725 /* Calculate the radial motion pivot-free potential and forces */
2726 static void do_radial_motion_pf(
2728 t_rotgrp *rotg, /* The rotation group */
2729 rvec x[], /* The positions */
2730 matrix box, /* The simulation box */
2731 double t, /* Time in picoseconds */
2732 gmx_large_int_t step, /* The time step */
2733 gmx_bool bOutstepRot, /* Output to main rotation output file */
2734 gmx_bool bOutstepSlab) /* Output per-slab data */
2736 int i,ii,iigrp,ifit,j;
2737 rvec xj; /* Current position */
2738 rvec xj_xc; /* xj - xc */
2739 rvec yj0_yc0; /* yj0 - yc0 */
2740 rvec tmp_f; /* Force */
2741 real alpha; /* a single angle between an actual and a reference position */
2742 real weight; /* single weight for a single angle */
2743 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
2744 rvec tmpvec, tmpvec2;
2745 rvec innersumvec; /* Precalculation of the inner sum */
2747 real fac,fac2,V=0.0;
2749 gmx_bool bCalcPotFit;
2751 /* For mass weighting: */
2752 real mj,wi,wj; /* Mass-weighting of the positions */
2756 erg=rotg->enfrotgrp;
2757 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT==rotg->eFittype);
2759 N_M = rotg->nat * erg->invmass;
2761 /* Get the current center of the rotation group: */
2762 get_center(erg->xc, erg->mc, rotg->nat, erg->xc_center);
2764 /* Precalculate Sum_i [ wi qi.(xi-xc) qi ] which is needed for every single j */
2765 clear_rvec(innersumvec);
2766 for (i=0; i < rotg->nat; i++)
2768 /* Mass-weighting */
2769 wi = N_M*erg->mc[i];
2771 /* Calculate qi. Note that xc_ref_center has already been subtracted from
2772 * x_ref in init_rot_group.*/
2773 mvmul(erg->rotmat, rotg->x_ref[i], tmpvec); /* tmpvec = Omega.(yi0-yc0) */
2775 cprod(rotg->vec, tmpvec, tmpvec2); /* tmpvec2 = v x Omega.(yi0-yc0) */
2777 /* v x Omega.(yi0-yc0) */
2778 unitv(tmpvec2, qi); /* qi = ----------------------- */
2779 /* | v x Omega.(yi0-yc0) | */
2781 rvec_sub(erg->xc[i], erg->xc_center, tmpvec); /* tmpvec = xi-xc */
2783 svmul(wi*iprod(qi, tmpvec), qi, tmpvec2);
2785 rvec_inc(innersumvec, tmpvec2);
2787 svmul(rotg->k*erg->invmass, innersumvec, innersumveckM);
2789 /* Each process calculates the forces on its local atoms */
2790 for (j=0; j<erg->nat_loc; j++)
2792 /* Local index of a rotation group atom */
2793 ii = erg->ind_loc[j];
2794 /* Position of this atom in the collective array */
2795 iigrp = erg->xc_ref_ind[j];
2796 /* Mass-weighting */
2797 mj = erg->mc[iigrp]; /* need the unsorted mass here */
2800 /* Current position of this atom: x[ii][XX/YY/ZZ] */
2801 copy_rvec(x[ii], xj);
2803 /* Shift this atom such that it is near its reference */
2804 shift_single_coord(box, xj, erg->xc_shifts[iigrp]);
2806 /* The (unrotated) reference position is yj0. yc0 has already
2807 * been subtracted in init_rot_group */
2808 copy_rvec(rotg->x_ref[iigrp], yj0_yc0); /* yj0_yc0 = yj0 - yc0 */
2810 /* Calculate Omega.(yj0-yc0) */
2811 mvmul(erg->rotmat, yj0_yc0, tmpvec2); /* tmpvec2 = Omega.(yj0 - yc0) */
2813 cprod(rotg->vec, tmpvec2, tmpvec); /* tmpvec = v x Omega.(yj0-yc0) */
2815 /* v x Omega.(yj0-yc0) */
2816 unitv(tmpvec, qj); /* qj = ----------------------- */
2817 /* | v x Omega.(yj0-yc0) | */
2819 /* Calculate (xj-xc) */
2820 rvec_sub(xj, erg->xc_center, xj_xc); /* xj_xc = xj-xc */
2822 fac = iprod(qj, xj_xc); /* fac = qj.(xj-xc) */
2825 /* Store the additional force so that it can be added to the force
2826 * array after the normal forces have been evaluated */
2827 svmul(-rotg->k*wj*fac, qj, tmp_f); /* part 1 of force */
2828 svmul(mj, innersumveckM, tmpvec); /* part 2 of force */
2829 rvec_inc(tmp_f, tmpvec);
2830 copy_rvec(tmp_f, erg->f_rot_loc[j]);
2833 /* If requested, also calculate the potential for a set of angles
2834 * near the current reference angle */
2837 for (ifit = 0; ifit < rotg->PotAngle_nstep; ifit++)
2839 /* Rotate with the alternative angle. Like rotate_local_reference(),
2840 * just for a single local atom */
2841 mvmul(erg->PotAngleFit->rotmat[ifit], yj0_yc0, tmpvec2); /* tmpvec2 = Omega*(yj0-yc0) */
2843 /* Calculate Omega.(yj0-u) */
2844 cprod(rotg->vec, tmpvec2, tmpvec); /* tmpvec = v x Omega.(yj0-yc0) */
2845 /* v x Omega.(yj0-yc0) */
2846 unitv(tmpvec, qj); /* qj = ----------------------- */
2847 /* | v x Omega.(yj0-yc0) | */
2849 fac = iprod(qj, xj_xc); /* fac = qj.(xj-xc) */
2852 /* Add to the rotation potential for this angle: */
2853 erg->PotAngleFit->V[ifit] += 0.5*rotg->k*wj*fac2;
2859 /* Add to the torque of this rotation group */
2860 erg->torque_v += torque(rotg->vec, tmp_f, xj, erg->xc_center);
2862 /* Calculate the angle between reference and actual rotation group atom. */
2863 angle(rotg, xj_xc, yj0_yc0, &alpha, &weight); /* angle in rad, weighted */
2864 erg->angle_v += alpha * weight;
2865 erg->weight_v += weight;
2870 } /* end of loop over local rotation group atoms */
2871 erg->V = 0.5*rotg->k*V;
2875 /* Precalculate the inner sum for the radial motion 2 forces */
2876 static void radial_motion2_precalc_inner_sum(t_rotgrp *rotg, rvec innersumvec)
2879 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
2880 rvec xi_xc; /* xj - xc */
2881 rvec tmpvec,tmpvec2;
2885 rvec v_xi_xc; /* v x (xj - u) */
2887 real wi; /* Mass-weighting of the positions */
2891 erg=rotg->enfrotgrp;
2892 N_M = rotg->nat * erg->invmass;
2894 /* Loop over the collective set of positions */
2896 for (i=0; i<rotg->nat; i++)
2898 /* Mass-weighting */
2899 wi = N_M*erg->mc[i];
2901 rvec_sub(erg->xc[i], erg->xc_center, xi_xc); /* xi_xc = xi-xc */
2903 /* Calculate ri. Note that xc_ref_center has already been subtracted from
2904 * x_ref in init_rot_group.*/
2905 mvmul(erg->rotmat, rotg->x_ref[i], ri); /* ri = Omega.(yi0-yc0) */
2907 cprod(rotg->vec, xi_xc, v_xi_xc); /* v_xi_xc = v x (xi-u) */
2909 fac = norm2(v_xi_xc);
2911 psiistar = 1.0/(fac + rotg->eps); /* psiistar = --------------------- */
2912 /* |v x (xi-xc)|^2 + eps */
2914 psii = gmx_invsqrt(fac); /* 1 */
2915 /* psii = ------------- */
2918 svmul(psii, v_xi_xc, si); /* si = psii * (v x (xi-xc) ) */
2920 fac = iprod(v_xi_xc, ri); /* fac = (v x (xi-xc)).ri */
2923 siri = iprod(si, ri); /* siri = si.ri */
2925 svmul(psiistar/psii, ri, tmpvec);
2926 svmul(psiistar*psiistar/(psii*psii*psii) * siri, si, tmpvec2);
2927 rvec_dec(tmpvec, tmpvec2);
2928 cprod(tmpvec, rotg->vec, tmpvec2);
2930 svmul(wi*siri, tmpvec2, tmpvec);
2932 rvec_inc(sumvec, tmpvec);
2934 svmul(rotg->k*erg->invmass, sumvec, innersumvec);
2938 /* Calculate the radial motion 2 potential and forces */
2939 static void do_radial_motion2(
2941 t_rotgrp *rotg, /* The rotation group */
2942 rvec x[], /* The positions */
2943 matrix box, /* The simulation box */
2944 double t, /* Time in picoseconds */
2945 gmx_large_int_t step, /* The time step */
2946 gmx_bool bOutstepRot, /* Output to main rotation output file */
2947 gmx_bool bOutstepSlab) /* Output per-slab data */
2949 int ii,iigrp,ifit,j;
2950 rvec xj; /* Position */
2951 real alpha; /* a single angle between an actual and a reference position */
2952 real weight; /* single weight for a single angle */
2953 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
2954 rvec xj_u; /* xj - u */
2955 rvec yj0_yc0; /* yj0 -yc0 */
2956 rvec tmpvec,tmpvec2;
2957 real fac,fit_fac,fac2,Vpart=0.0;
2960 rvec v_xj_u; /* v x (xj - u) */
2962 real mj,wj; /* For mass-weighting of the positions */
2966 gmx_bool bCalcPotFit;
2969 erg=rotg->enfrotgrp;
2971 bPF = rotg->eType==erotgRM2PF;
2972 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT==rotg->eFittype);
2975 clear_rvec(yj0_yc0); /* Make the compiler happy */
2977 clear_rvec(innersumvec);
2980 /* For the pivot-free variant we have to use the current center of
2981 * mass of the rotation group instead of the pivot u */
2982 get_center(erg->xc, erg->mc, rotg->nat, erg->xc_center);
2984 /* Also, we precalculate the second term of the forces that is identical
2985 * (up to the weight factor mj) for all forces */
2986 radial_motion2_precalc_inner_sum(rotg,innersumvec);
2989 N_M = rotg->nat * erg->invmass;
2991 /* Each process calculates the forces on its local atoms */
2992 for (j=0; j<erg->nat_loc; j++)
2996 /* Local index of a rotation group atom */
2997 ii = erg->ind_loc[j];
2998 /* Position of this atom in the collective array */
2999 iigrp = erg->xc_ref_ind[j];
3000 /* Mass-weighting */
3001 mj = erg->mc[iigrp];
3003 /* Current position of this atom: x[ii] */
3004 copy_rvec(x[ii], xj);
3006 /* Shift this atom such that it is near its reference */
3007 shift_single_coord(box, xj, erg->xc_shifts[iigrp]);
3009 /* The (unrotated) reference position is yj0. yc0 has already
3010 * been subtracted in init_rot_group */
3011 copy_rvec(rotg->x_ref[iigrp], yj0_yc0); /* yj0_yc0 = yj0 - yc0 */
3013 /* Calculate Omega.(yj0-yc0) */
3014 mvmul(erg->rotmat, yj0_yc0, rj); /* rj = Omega.(yj0-yc0) */
3019 copy_rvec(erg->x_loc_pbc[j], xj);
3020 copy_rvec(erg->xr_loc[j], rj); /* rj = Omega.(yj0-u) */
3022 /* Mass-weighting */
3025 /* Calculate (xj-u) resp. (xj-xc) */
3026 rvec_sub(xj, erg->xc_center, xj_u); /* xj_u = xj-u */
3028 cprod(rotg->vec, xj_u, v_xj_u); /* v_xj_u = v x (xj-u) */
3030 fac = norm2(v_xj_u);
3032 psijstar = 1.0/(fac + rotg->eps); /* psistar = -------------------- */
3033 /* |v x (xj-u)|^2 + eps */
3035 psij = gmx_invsqrt(fac); /* 1 */
3036 /* psij = ------------ */
3039 svmul(psij, v_xj_u, sj); /* sj = psij * (v x (xj-u) ) */
3041 fac = iprod(v_xj_u, rj); /* fac = (v x (xj-u)).rj */
3044 sjrj = iprod(sj, rj); /* sjrj = sj.rj */
3046 svmul(psijstar/psij, rj, tmpvec);
3047 svmul(psijstar*psijstar/(psij*psij*psij) * sjrj, sj, tmpvec2);
3048 rvec_dec(tmpvec, tmpvec2);
3049 cprod(tmpvec, rotg->vec, tmpvec2);
3051 /* Store the additional force so that it can be added to the force
3052 * array after the normal forces have been evaluated */
3053 svmul(-rotg->k*wj*sjrj, tmpvec2, tmpvec);
3054 svmul(mj, innersumvec, tmpvec2); /* This is != 0 only for the pivot-free variant */
3056 rvec_add(tmpvec2, tmpvec, erg->f_rot_loc[j]);
3057 Vpart += wj*psijstar*fac2;
3059 /* If requested, also calculate the potential for a set of angles
3060 * near the current reference angle */
3063 for (ifit = 0; ifit < rotg->PotAngle_nstep; ifit++)
3067 mvmul(erg->PotAngleFit->rotmat[ifit], yj0_yc0, fit_rj); /* fit_rj = Omega.(yj0-yc0) */
3071 /* Position of this atom in the collective array */
3072 iigrp = erg->xc_ref_ind[j];
3073 /* Rotate with the alternative angle. Like rotate_local_reference(),
3074 * just for a single local atom */
3075 mvmul(erg->PotAngleFit->rotmat[ifit], rotg->x_ref[iigrp], fit_rj); /* fit_rj = Omega*(yj0-u) */
3077 fit_fac = iprod(v_xj_u, fit_rj); /* fac = (v x (xj-u)).fit_rj */
3078 /* Add to the rotation potential for this angle: */
3079 erg->PotAngleFit->V[ifit] += 0.5*rotg->k*wj*psijstar*fit_fac*fit_fac;
3085 /* Add to the torque of this rotation group */
3086 erg->torque_v += torque(rotg->vec, erg->f_rot_loc[j], xj, erg->xc_center);
3088 /* Calculate the angle between reference and actual rotation group atom. */
3089 angle(rotg, xj_u, rj, &alpha, &weight); /* angle in rad, weighted */
3090 erg->angle_v += alpha * weight;
3091 erg->weight_v += weight;
3096 } /* end of loop over local rotation group atoms */
3097 erg->V = 0.5*rotg->k*Vpart;
3101 /* Determine the smallest and largest position vector (with respect to the
3102 * rotation vector) for the reference group */
3103 static void get_firstlast_atom_ref(
3108 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
3110 real xcproj; /* The projection of a reference position on the
3112 real minproj, maxproj; /* Smallest and largest projection on v */
3116 erg=rotg->enfrotgrp;
3118 /* Start with some value */
3119 minproj = iprod(rotg->x_ref[0], rotg->vec);
3122 /* This is just to ensure that it still works if all the atoms of the
3123 * reference structure are situated in a plane perpendicular to the rotation
3126 *lastindex = rotg->nat-1;
3128 /* Loop over all atoms of the reference group,
3129 * project them on the rotation vector to find the extremes */
3130 for (i=0; i<rotg->nat; i++)
3132 xcproj = iprod(rotg->x_ref[i], rotg->vec);
3133 if (xcproj < minproj)
3138 if (xcproj > maxproj)
3147 /* Allocate memory for the slabs */
3148 static void allocate_slabs(
3155 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
3159 erg=rotg->enfrotgrp;
3161 /* More slabs than are defined for the reference are never needed */
3162 nslabs = erg->slab_last_ref - erg->slab_first_ref + 1;
3164 /* Remember how many we allocated */
3165 erg->nslabs_alloc = nslabs;
3167 if (MASTER(cr) && bVerbose)
3168 fprintf(fplog, "%s allocating memory to store data for %d slabs (rotation group %d).\n",
3170 snew(erg->slab_center , nslabs);
3171 snew(erg->slab_center_ref , nslabs);
3172 snew(erg->slab_weights , nslabs);
3173 snew(erg->slab_torque_v , nslabs);
3174 snew(erg->slab_data , nslabs);
3175 snew(erg->gn_atom , nslabs);
3176 snew(erg->gn_slabind , nslabs);
3177 snew(erg->slab_innersumvec, nslabs);
3178 for (i=0; i<nslabs; i++)
3180 snew(erg->slab_data[i].x , rotg->nat);
3181 snew(erg->slab_data[i].ref , rotg->nat);
3182 snew(erg->slab_data[i].weight, rotg->nat);
3184 snew(erg->xc_ref_sorted, rotg->nat);
3185 snew(erg->xc_sortind , rotg->nat);
3186 snew(erg->firstatom , nslabs);
3187 snew(erg->lastatom , nslabs);
3191 /* From the extreme coordinates of the reference group, determine the first
3192 * and last slab of the reference. We can never have more slabs in the real
3193 * simulation than calculated here for the reference.
3195 static void get_firstlast_slab_ref(t_rotgrp *rotg, real mc[], int ref_firstindex, int ref_lastindex)
3197 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
3198 int first,last,firststart;
3202 erg=rotg->enfrotgrp;
3203 first = get_first_slab(rotg, erg->max_beta, rotg->x_ref[ref_firstindex]);
3204 last = get_last_slab( rotg, erg->max_beta, rotg->x_ref[ref_lastindex ]);
3207 while (get_slab_weight(first, rotg, rotg->x_ref, mc, &dummy) > WEIGHT_MIN)
3211 erg->slab_first_ref = first+1;
3212 while (get_slab_weight(last, rotg, rotg->x_ref, mc, &dummy) > WEIGHT_MIN)
3216 erg->slab_last_ref = last-1;
3218 erg->slab_buffer = firststart - erg->slab_first_ref;
3223 static void init_rot_group(FILE *fplog,t_commrec *cr,int g,t_rotgrp *rotg,
3224 rvec *x,gmx_mtop_t *mtop,gmx_bool bVerbose,FILE *out_slabs, gmx_bool bOutputCenters)
3228 gmx_bool bFlex,bColl;
3230 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
3231 int ref_firstindex, ref_lastindex;
3232 real mass,totalmass;
3236 /* Do we have a flexible axis? */
3237 bFlex = ISFLEX(rotg);
3238 /* Do we use a global set of coordinates? */
3239 bColl = ISCOLL(rotg);
3241 erg=rotg->enfrotgrp;
3243 /* Allocate space for collective coordinates if needed */
3246 snew(erg->xc , rotg->nat);
3247 snew(erg->xc_shifts , rotg->nat);
3248 snew(erg->xc_eshifts, rotg->nat);
3250 /* Save the original (whole) set of positions such that later the
3251 * molecule can always be made whole again */
3252 snew(erg->xc_old , rotg->nat);
3255 for (i=0; i<rotg->nat; i++)
3258 copy_rvec(x[ii], erg->xc_old[i]);
3263 gmx_bcast(rotg->nat*sizeof(erg->xc_old[0]),erg->xc_old, cr);
3266 if (rotg->eFittype == erotgFitNORM)
3268 snew(erg->xc_ref_length, rotg->nat); /* in case fit type NORM is chosen */
3269 snew(erg->xc_norm , rotg->nat);
3274 snew(erg->xr_loc , rotg->nat);
3275 snew(erg->x_loc_pbc, rotg->nat);
3278 snew(erg->f_rot_loc , rotg->nat);
3279 snew(erg->xc_ref_ind, rotg->nat);
3281 /* Make space for the calculation of the potential at other angles (used
3282 * for fitting only) */
3283 if (erotgFitPOT == rotg->eFittype)
3285 snew(erg->PotAngleFit, 1);
3286 snew(erg->PotAngleFit->degangle, rotg->PotAngle_nstep);
3287 snew(erg->PotAngleFit->V , rotg->PotAngle_nstep);
3288 snew(erg->PotAngleFit->rotmat , rotg->PotAngle_nstep);
3290 /* Get the set of angles around the reference angle */
3291 start = -0.5 * (rotg->PotAngle_nstep - 1)*rotg->PotAngle_step;
3292 for (i = 0; i < rotg->PotAngle_nstep; i++)
3293 erg->PotAngleFit->degangle[i] = start + i*rotg->PotAngle_step;
3297 erg->PotAngleFit = NULL;
3300 /* xc_ref_ind needs to be set to identity in the serial case */
3302 for (i=0; i<rotg->nat; i++)
3303 erg->xc_ref_ind[i] = i;
3305 /* Copy the masses so that the center can be determined. For all types of
3306 * enforced rotation, we store the masses in the erg->mc array. */
3307 snew(erg->mc, rotg->nat);
3309 snew(erg->mc_sorted, rotg->nat);
3311 snew(erg->m_loc, rotg->nat);
3313 for (i=0; i<rotg->nat; i++)
3317 gmx_mtop_atomnr_to_atom(mtop,rotg->ind[i],&atom);
3327 erg->invmass = 1.0/totalmass;
3329 /* Set xc_ref_center for any rotation potential */
3330 if ((rotg->eType==erotgISO) || (rotg->eType==erotgPM) || (rotg->eType==erotgRM) || (rotg->eType==erotgRM2))
3332 /* Set the pivot point for the fixed, stationary-axis potentials. This
3333 * won't change during the simulation */
3334 copy_rvec(rotg->pivot, erg->xc_ref_center);
3335 copy_rvec(rotg->pivot, erg->xc_center );
3339 /* Center of the reference positions */
3340 get_center(rotg->x_ref, erg->mc, rotg->nat, erg->xc_ref_center);
3342 /* Center of the actual positions */
3345 snew(xdum, rotg->nat);
3346 for (i=0; i<rotg->nat; i++)
3349 copy_rvec(x[ii], xdum[i]);
3351 get_center(xdum, erg->mc, rotg->nat, erg->xc_center);
3356 gmx_bcast(sizeof(erg->xc_center), erg->xc_center, cr);
3360 if ( (rotg->eType != erotgFLEX) && (rotg->eType != erotgFLEX2) )
3362 /* Put the reference positions into origin: */
3363 for (i=0; i<rotg->nat; i++)
3364 rvec_dec(rotg->x_ref[i], erg->xc_ref_center);
3367 /* Enforced rotation with flexible axis */
3370 /* Calculate maximum beta value from minimum gaussian (performance opt.) */
3371 erg->max_beta = calc_beta_max(rotg->min_gaussian, rotg->slab_dist);
3373 /* Determine the smallest and largest coordinate with respect to the rotation vector */
3374 get_firstlast_atom_ref(rotg, &ref_firstindex, &ref_lastindex);
3376 /* From the extreme coordinates of the reference group, determine the first
3377 * and last slab of the reference. */
3378 get_firstlast_slab_ref(rotg, erg->mc, ref_firstindex, ref_lastindex);
3380 /* Allocate memory for the slabs */
3381 allocate_slabs(rotg, fplog, g, bVerbose, cr);
3383 /* Flexible rotation: determine the reference centers for the rest of the simulation */
3384 erg->slab_first = erg->slab_first_ref;
3385 erg->slab_last = erg->slab_last_ref;
3386 get_slab_centers(rotg,rotg->x_ref,erg->mc,cr,g,-1,out_slabs,bOutputCenters,TRUE);
3388 /* Length of each x_rotref vector from center (needed if fit routine NORM is chosen): */
3389 if (rotg->eFittype == erotgFitNORM)
3391 for (i=0; i<rotg->nat; i++)
3393 rvec_sub(rotg->x_ref[i], erg->xc_ref_center, coord);
3394 erg->xc_ref_length[i] = norm(coord);
3401 extern void dd_make_local_rotation_groups(gmx_domdec_t *dd,t_rot *rot)
3406 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
3410 for(g=0; g<rot->ngrp; g++)
3412 rotg = &rot->grp[g];
3413 erg = rotg->enfrotgrp;
3416 dd_make_local_group_indices(ga2la,rotg->nat,rotg->ind,
3417 &erg->nat_loc,&erg->ind_loc,&erg->nalloc_loc,erg->xc_ref_ind);
3422 /* Calculate the size of the MPI buffer needed in reduce_output() */
3423 static int calc_mpi_bufsize(t_rot *rot)
3426 int count_group, count_total;
3428 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
3432 for (g=0; g<rot->ngrp; g++)
3434 rotg = &rot->grp[g];
3435 erg = rotg->enfrotgrp;
3437 /* Count the items that are transferred for this group: */
3438 count_group = 4; /* V, torque, angle, weight */
3440 /* Add the maximum number of slabs for flexible groups */
3442 count_group += erg->slab_last_ref - erg->slab_first_ref + 1;
3444 /* Add space for the potentials at different angles: */
3445 if (erotgFitPOT == rotg->eFittype)
3446 count_group += rotg->PotAngle_nstep;
3448 /* Add to the total number: */
3449 count_total += count_group;
3456 extern void init_rot(FILE *fplog,t_inputrec *ir,int nfile,const t_filenm fnm[],
3457 t_commrec *cr, rvec *x, matrix box, gmx_mtop_t *mtop, const output_env_t oenv,
3458 gmx_bool bVerbose, unsigned long Flags)
3463 int nat_max=0; /* Size of biggest rotation group */
3464 gmx_enfrot_t er; /* Pointer to the enforced rotation buffer variables */
3465 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
3466 rvec *x_pbc=NULL; /* Space for the pbc-correct atom positions */
3469 if ( (PAR(cr)) && !DOMAINDECOMP(cr) )
3470 gmx_fatal(FARGS, "Enforced rotation is only implemented for domain decomposition!");
3472 if ( MASTER(cr) && bVerbose)
3473 fprintf(stdout, "%s Initializing ...\n", RotStr);
3477 snew(rot->enfrot, 1);
3481 /* When appending, skip first output to avoid duplicate entries in the data files */
3482 if (er->Flags & MD_APPENDFILES)
3487 /* Output every step for reruns */
3488 if (er->Flags & MD_RERUN)
3491 fprintf(fplog, "%s rerun - will write rotation output every available step.\n", RotStr);
3496 er->out_slabs = NULL;
3497 if ( MASTER(cr) && HaveFlexibleGroups(rot) )
3498 er->out_slabs = open_slab_out(opt2fn("-rs",nfile,fnm), rot, oenv);
3502 /* Remove pbc, make molecule whole.
3503 * When ir->bContinuation=TRUE this has already been done, but ok. */
3504 snew(x_pbc,mtop->natoms);
3505 m_rveccopy(mtop->natoms,x,x_pbc);
3506 do_pbc_first_mtop(NULL,ir->ePBC,box,mtop,x_pbc);
3509 for (g=0; g<rot->ngrp; g++)
3511 rotg = &rot->grp[g];
3514 fprintf(fplog,"%s group %d type '%s'\n", RotStr, g, erotg_names[rotg->eType]);
3518 /* Allocate space for the rotation group's data: */
3519 snew(rotg->enfrotgrp, 1);
3520 erg = rotg->enfrotgrp;
3522 nat_max=max(nat_max, rotg->nat);
3527 erg->nalloc_loc = 0;
3528 erg->ind_loc = NULL;
3532 erg->nat_loc = rotg->nat;
3533 erg->ind_loc = rotg->ind;
3535 init_rot_group(fplog,cr,g,rotg,x_pbc,mtop,bVerbose,er->out_slabs,
3536 !(er->Flags & MD_APPENDFILES) ); /* Do not output the reference centers
3537 * again if we are appending */
3541 /* Allocate space for enforced rotation buffer variables */
3542 er->bufsize = nat_max;
3543 snew(er->data, nat_max);
3544 snew(er->xbuf, nat_max);
3545 snew(er->mbuf, nat_max);
3547 /* Buffers for MPI reducing torques, angles, weights (for each group), and V */
3550 er->mpi_bufsize = calc_mpi_bufsize(rot) + 100; /* larger to catch errors */
3551 snew(er->mpi_inbuf , er->mpi_bufsize);
3552 snew(er->mpi_outbuf, er->mpi_bufsize);
3556 er->mpi_bufsize = 0;
3557 er->mpi_inbuf = NULL;
3558 er->mpi_outbuf = NULL;
3561 /* Only do I/O on the MASTER */
3562 er->out_angles = NULL;
3564 er->out_torque = NULL;
3567 er->out_rot = open_rot_out(opt2fn("-ro",nfile,fnm), rot, oenv);
3569 if (rot->nstsout > 0)
3571 if ( HaveFlexibleGroups(rot) || HavePotFitGroups(rot) )
3572 er->out_angles = open_angles_out(opt2fn("-ra",nfile,fnm), rot, oenv);
3573 if ( HaveFlexibleGroups(rot) )
3574 er->out_torque = open_torque_out(opt2fn("-rt",nfile,fnm), rot, oenv);
3582 extern void finish_rot(FILE *fplog,t_rot *rot)
3584 gmx_enfrot_t er; /* Pointer to the enforced rotation buffer variables */
3589 gmx_fio_fclose(er->out_rot);
3591 gmx_fio_fclose(er->out_slabs);
3593 gmx_fio_fclose(er->out_angles);
3595 gmx_fio_fclose(er->out_torque);
3599 /* Rotate the local reference positions and store them in
3600 * erg->xr_loc[0...(nat_loc-1)]
3602 * Note that we already subtracted u or y_c from the reference positions
3603 * in init_rot_group().
3605 static void rotate_local_reference(t_rotgrp *rotg)
3607 gmx_enfrotgrp_t erg;
3611 erg=rotg->enfrotgrp;
3613 for (i=0; i<erg->nat_loc; i++)
3615 /* Index of this rotation group atom with respect to the whole rotation group */
3616 ii = erg->xc_ref_ind[i];
3618 mvmul(erg->rotmat, rotg->x_ref[ii], erg->xr_loc[i]);
3623 /* Select the PBC representation for each local x position and store that
3624 * for later usage. We assume the right PBC image of an x is the one nearest to
3625 * its rotated reference */
3626 static void choose_pbc_image(rvec x[], t_rotgrp *rotg, matrix box, int npbcdim)
3629 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
3634 erg=rotg->enfrotgrp;
3636 for (i=0; i<erg->nat_loc; i++)
3640 /* Index of a rotation group atom */
3641 ii = erg->ind_loc[i];
3643 /* Get the reference position. The pivot was already
3644 * subtracted in init_rot_group() from the reference positions. Also,
3645 * the reference positions have already been rotated in
3646 * rotate_local_reference() */
3647 copy_rvec(erg->xr_loc[i], xref);
3649 /* Subtract the (old) center from the current positions
3650 * (just to determine the shifts!) */
3651 rvec_sub(x[ii], erg->xc_center, xcurr);
3653 /* Shortest PBC distance between the atom and its reference */
3654 rvec_sub(xcurr, xref, dx);
3656 /* Determine the shift for this atom */
3657 for(m=npbcdim-1; m>=0; m--)
3659 while (dx[m] < -0.5*box[m][m])
3661 for(d=0; d<DIM; d++)
3665 while (dx[m] >= 0.5*box[m][m])
3667 for(d=0; d<DIM; d++)
3673 /* Apply the shift to the current atom */
3674 copy_rvec(x[ii], erg->x_loc_pbc[i]);
3675 shift_single_coord(box, erg->x_loc_pbc[i], shift);
3680 extern void do_rotation(
3686 gmx_large_int_t step,
3687 gmx_wallcycle_t wcycle,
3693 gmx_bool outstep_slab, outstep_rot;
3694 gmx_bool bFlex,bColl;
3696 gmx_enfrot_t er; /* Pointer to the enforced rotation buffer variables */
3697 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
3699 t_gmx_potfit *fit=NULL; /* For fit type 'potential' determine the fit
3700 angle via the potential minimum */
3709 /* When to output in main rotation output file */
3710 outstep_rot = do_per_step(step, rot->nstrout) && er->bOut;
3711 /* When to output per-slab data */
3712 outstep_slab = do_per_step(step, rot->nstsout) && er->bOut;
3714 /* Output time into rotation output file */
3715 if (outstep_rot && MASTER(cr))
3716 fprintf(er->out_rot, "%12.3e",t);
3718 /**************************************************************************/
3719 /* First do ALL the communication! */
3720 for(g=0; g<rot->ngrp; g++)
3722 rotg = &rot->grp[g];
3723 erg=rotg->enfrotgrp;
3725 /* Do we have a flexible axis? */
3726 bFlex = ISFLEX(rotg);
3727 /* Do we use a collective (global) set of coordinates? */
3728 bColl = ISCOLL(rotg);
3730 /* Calculate the rotation matrix for this angle: */
3731 erg->degangle = rotg->rate * t;
3732 calc_rotmat(rotg->vec,erg->degangle,erg->rotmat);
3736 /* Transfer the rotation group's positions such that every node has
3737 * all of them. Every node contributes its local positions x and stores
3738 * it in the collective erg->xc array. */
3739 communicate_group_positions(cr,erg->xc, erg->xc_shifts, erg->xc_eshifts, bNS,
3740 x, rotg->nat, erg->nat_loc, erg->ind_loc, erg->xc_ref_ind, erg->xc_old, box);
3744 /* Fill the local masses array;
3745 * this array changes in DD/neighborsearching steps */
3748 for (i=0; i<erg->nat_loc; i++)
3750 /* Index of local atom w.r.t. the collective rotation group */
3751 ii = erg->xc_ref_ind[i];
3752 erg->m_loc[i] = erg->mc[ii];
3756 /* Calculate Omega*(y_i-y_c) for the local positions */
3757 rotate_local_reference(rotg);
3759 /* Choose the nearest PBC images of the group atoms with respect
3760 * to the rotated reference positions */
3761 choose_pbc_image(x, rotg, box, 3);
3763 /* Get the center of the rotation group */
3764 if ( (rotg->eType==erotgISOPF) || (rotg->eType==erotgPMPF) )
3765 get_center_comm(cr, erg->x_loc_pbc, erg->m_loc, erg->nat_loc, rotg->nat, erg->xc_center);
3768 } /* End of loop over rotation groups */
3770 /**************************************************************************/
3771 /* Done communicating, we can start to count cycles now ... */
3772 wallcycle_start(wcycle, ewcROT);
3778 for(g=0; g<rot->ngrp; g++)
3780 rotg = &rot->grp[g];
3781 erg=rotg->enfrotgrp;
3783 bFlex = ISFLEX(rotg);
3784 bColl = ISCOLL(rotg);
3786 if (outstep_rot && MASTER(cr))
3787 fprintf(er->out_rot, "%12.4f", erg->degangle);
3789 /* Calculate angles and rotation matrices for potential fitting: */
3790 if ( (outstep_rot || outstep_slab) && (erotgFitPOT == rotg->eFittype) )
3792 fit = erg->PotAngleFit;
3793 for (i = 0; i < rotg->PotAngle_nstep; i++)
3795 calc_rotmat(rotg->vec, erg->degangle + fit->degangle[i], fit->rotmat[i]);
3797 /* Clear value from last step */
3798 erg->PotAngleFit->V[i] = 0.0;
3802 /* Clear values from last time step */
3804 erg->torque_v = 0.0;
3806 erg->weight_v = 0.0;
3814 do_fixed(cr,rotg,x,box,t,step,outstep_rot,outstep_slab);
3817 do_radial_motion(cr,rotg,x,box,t,step,outstep_rot,outstep_slab);
3820 do_radial_motion_pf(cr,rotg,x,box,t,step,outstep_rot,outstep_slab);
3824 do_radial_motion2(cr,rotg,x,box,t,step,outstep_rot,outstep_slab);
3828 /* Subtract the center of the rotation group from the collective positions array
3829 * Also store the center in erg->xc_center since it needs to be subtracted
3830 * in the low level routines from the local coordinates as well */
3831 get_center(erg->xc, erg->mc, rotg->nat, erg->xc_center);
3832 svmul(-1.0, erg->xc_center, transvec);
3833 translate_x(erg->xc, rotg->nat, transvec);
3834 do_flexible(cr,er,rotg,g,x,box,t,step,outstep_rot,outstep_slab);
3838 /* Do NOT subtract the center of mass in the low level routines! */
3839 clear_rvec(erg->xc_center);
3840 do_flexible(cr,er,rotg,g,x,box,t,step,outstep_rot,outstep_slab);
3843 gmx_fatal(FARGS, "No such rotation potential.");
3850 fprintf(stderr, "%s calculation (step %d) took %g seconds.\n", RotStr, step, MPI_Wtime()-t0);
3853 /* Stop the cycle counter and add to the force cycles for load balancing */
3854 cycles_rot = wallcycle_stop(wcycle,ewcROT);
3855 if (DOMAINDECOMP(cr) && wcycle)
3856 dd_cycles_add(cr->dd,cycles_rot,ddCyclF);