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 "mpelogging.h"
59 #include "groupcoord.h"
60 #include "pull_rotation.h"
64 static char *RotStr = {"Enforced rotation:"};
67 /* Set the minimum weight for the determination of the slab centers */
68 #define WEIGHT_MIN (10*GMX_FLOAT_MIN)
70 /* Helper structure for sorting positions along rotation vector */
72 real xcproj; /* Projection of xc on the rotation vector */
73 int ind; /* Index of xc */
75 rvec x; /* Position */
76 rvec x_ref; /* Reference position */
80 /* Enforced rotation / flexible: determine the angle of each slab */
81 typedef struct gmx_slabdata
83 int nat; /* Number of atoms belonging to this slab */
84 rvec *x; /* The positions belonging to this slab. In
85 general, this should be all positions of the
86 whole rotation group, but we leave those away
87 that have a small enough weight */
88 rvec *ref; /* Same for reference */
89 real *weight; /* The weight for each atom */
93 /* Enforced rotation data for all groups */
94 typedef struct gmx_enfrot
96 FILE *out_rot; /* Output file for rotation data */
97 FILE *out_torque; /* Output file for torque data */
98 FILE *out_angles; /* Output file for slab angles for flexible type */
99 FILE *out_slabs; /* Output file for slab centers */
100 int bufsize; /* Allocation size of buf */
101 rvec *xbuf; /* Coordinate buffer variable for sorting */
102 real *mbuf; /* Masses buffer variable for sorting */
103 sort_along_vec_t *data; /* Buffer variable needed for position sorting */
104 real *mpi_inbuf; /* MPI buffer */
105 real *mpi_outbuf; /* MPI buffer */
106 int mpi_bufsize; /* Allocation size of in & outbuf */
107 real Vrot; /* (Local) part of the enf. rotation potential */
111 /* Global enforced rotation data for a single rotation group */
112 typedef struct gmx_enfrotgrp
114 real degangle; /* Rotation angle in degrees */
115 matrix rotmat; /* Rotation matrix */
116 atom_id *ind_loc; /* Local rotation indices */
117 int nat_loc; /* Number of local group atoms */
118 int nalloc_loc; /* Allocation size for ind_loc and weight_loc */
120 real V; /* Rotation potential for this rotation group */
121 rvec *f_rot_loc; /* Array to store the forces on the local atoms
122 resulting from enforced rotation potential */
124 /* Collective coordinates for the whole rotation group */
125 real *xc_ref_length; /* Length of each x_rotref vector after x_rotref
126 has been put into origin */
127 int *xc_ref_ind; /* Position of each local atom in the collective
129 rvec xc_center; /* Center of the rotation group positions, may
131 rvec xc_ref_center; /* dito, for the reference positions */
132 rvec *xc; /* Current (collective) positions */
133 ivec *xc_shifts; /* Current (collective) shifts */
134 ivec *xc_eshifts; /* Extra shifts since last DD step */
135 rvec *xc_old; /* Old (collective) positions */
136 rvec *xc_norm; /* Normalized form of the current positions */
137 rvec *xc_ref_sorted; /* Reference positions (sorted in the same order
138 as xc when sorted) */
139 int *xc_sortind; /* Where is a position found after sorting? */
140 real *mc; /* Collective masses */
142 real invmass; /* one over the total mass of the rotation group */
144 real torque_v; /* Torque in the direction of rotation vector */
145 real angle_v; /* Actual angle of the whole rotation group */
146 /* Fixed rotation only */
147 real weight_v; /* Weights for angle determination */
148 rvec *xr_loc; /* Local reference coords, correctly rotated */
149 rvec *x_loc_pbc; /* Local current coords, correct PBC image */
150 real *m_loc; /* Masses of the current local atoms */
152 /* Flexible rotation only */
153 int nslabs_alloc; /* For this many slabs memory is allocated */
154 int slab_first; /* Lowermost slab for that the calculation needs
155 to be performed at a given time step */
156 int slab_last; /* Uppermost slab ... */
157 int slab_first_ref; /* First slab for which reference COG is stored */
158 int slab_last_ref; /* Last ... */
159 int slab_buffer; /* Slab buffer region around reference slabs */
160 int *firstatom; /* First relevant atom for a slab */
161 int *lastatom; /* Last relevant atom for a slab */
162 rvec *slab_center; /* Gaussian-weighted slab center (COG) */
163 rvec *slab_center_ref; /* Gaussian-weighted slab COG for the
164 reference positions */
165 real *slab_weights; /* Sum of gaussian weights in a slab */
166 real *slab_torque_v; /* Torque T = r x f for each slab. */
167 /* torque_v = m.v = angular momentum in the
169 real max_beta; /* min_gaussian from inputrec->rotgrp is the
170 minimum value the gaussian must have so that
171 the force is actually evaluated max_beta is
172 just another way to put it */
173 real *gn_atom; /* Precalculated gaussians for a single atom */
174 int *gn_slabind; /* Tells to which slab each precalculated gaussian
176 rvec *slab_innersumvec;/* Inner sum of the flexible2 potential per slab;
177 this is precalculated for optimization reasons */
178 t_gmx_slabdata *slab_data; /* Holds atom positions and gaussian weights
179 of atoms belonging to a slab */
183 /* Activate output of forces for correctness checks */
184 /* #define PRINT_FORCES */
186 #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]);
187 #define PRINT_POT_TAU fprintf(stderr,"potential = %15.8f\n" "torque = %15.8f", erg->V, erg->torque_v);
189 #define PRINT_FORCE_J
190 #define PRINT_POT_TAU
193 /* Shortcuts for often used queries */
194 #define ISFLEX(rg) ( (rg->eType==erotgFLEX) || (rg->eType==erotgFLEXT) || (rg->eType==erotgFLEX2) || (rg->eType==erotgFLEX2T) )
195 #define ISCOLL(rg) ( (rg->eType==erotgFLEX) || (rg->eType==erotgFLEXT) || (rg->eType==erotgFLEX2) || (rg->eType==erotgFLEX2T) || (rg->eType==erotgRMPF) || (rg->eType==erotgRM2PF) )
199 static double** allocate_square_matrix(int dim)
213 static void free_square_matrix(double** mat, int dim)
218 for (i=0; i<dim; i++)
224 /* Output rotation energy, torques, etc. for each rotation group */
225 static void reduce_output(t_commrec *cr, t_rot *rot, real t, gmx_large_int_t step)
227 int g,i,islab,nslabs=0;
228 int count; /* MPI element counter */
230 gmx_enfrot_t er; /* Pointer to the enforced rotation buffer variables */
231 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
237 /* Fill the MPI buffer with stuff to reduce: */
241 for (g=0; g < rot->ngrp; g++)
244 erg = rotg->enfrotgrp;
245 nslabs = erg->slab_last - erg->slab_first + 1;
246 er->mpi_inbuf[count++] = erg->V;
247 er->mpi_inbuf[count++] = erg->torque_v;
248 er->mpi_inbuf[count++] = erg->angle_v;
249 er->mpi_inbuf[count++] = erg->weight_v; /* weights are not needed for flex types, but this is just a single value */
252 /* (Re-)allocate memory for MPI buffer: */
253 if (er->mpi_bufsize < count+nslabs)
255 er->mpi_bufsize = count+nslabs;
256 srenew(er->mpi_inbuf , er->mpi_bufsize);
257 srenew(er->mpi_outbuf, er->mpi_bufsize);
259 for (i=0; i<nslabs; i++)
260 er->mpi_inbuf[count++] = erg->slab_torque_v[i];
264 MPI_Reduce(er->mpi_inbuf, er->mpi_outbuf, count, GMX_MPI_REAL, MPI_SUM, MASTERRANK(cr), cr->mpi_comm_mygroup);
266 /* Copy back the reduced data from the buffer on the master */
270 for (g=0; g < rot->ngrp; g++)
273 erg = rotg->enfrotgrp;
274 nslabs = erg->slab_last - erg->slab_first + 1;
275 erg->V = er->mpi_outbuf[count++];
276 erg->torque_v = er->mpi_outbuf[count++];
277 erg->angle_v = er->mpi_outbuf[count++];
278 erg->weight_v = er->mpi_outbuf[count++];
281 for (i=0; i<nslabs; i++)
282 erg->slab_torque_v[i] = er->mpi_outbuf[count++];
291 /* Angle and torque for each rotation group */
292 for (g=0; g < rot->ngrp; g++)
295 bFlex = ISFLEX(rotg);
299 /* Output to main rotation output file: */
300 if ( do_per_step(step, rot->nstrout) )
303 fprintf(er->out_rot, "%12.4f", erg->angle_v); /* RMSD fit angle */
305 fprintf(er->out_rot, "%12.4f", (erg->angle_v/erg->weight_v)*180.0*M_1_PI);
306 fprintf(er->out_rot, "%12.3e", erg->torque_v);
307 fprintf(er->out_rot, "%12.3e", erg->V);
310 /* Output to torque log file: */
311 if ( bFlex && do_per_step(step, rot->nstsout) )
313 fprintf(er->out_torque, "%12.3e%6d", t, g);
314 for (i=erg->slab_first; i<=erg->slab_last; i++)
316 islab = i - erg->slab_first; /* slab index */
317 /* Only output if enough weight is in slab */
318 if (erg->slab_weights[islab] > rotg->min_gaussian)
319 fprintf(er->out_torque, "%6d%12.3e", i, erg->slab_torque_v[islab]);
321 fprintf(er->out_torque , "\n");
324 if ( do_per_step(step, rot->nstrout) )
325 fprintf(er->out_rot, "\n");
330 /* Add the forces from enforced rotation potential to the local forces.
331 * Should be called after the SR forces have been evaluated */
332 extern real add_rot_forces(t_rot *rot, rvec f[], t_commrec *cr, gmx_large_int_t step, real t)
336 gmx_enfrot_t er; /* Pointer to the enforced rotation buffer variables */
337 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
342 GMX_MPE_LOG(ev_add_rot_forces_start);
344 /* Reduce energy,torque, angles etc. to get the sum values (per rotation group)
345 * on the master and output these values to file. */
346 if ( do_per_step(step, rot->nstrout) || do_per_step(step, rot->nstsout) )
347 reduce_output(cr, rot, t, step);
349 /* Total rotation potential is the sum over all rotation groups */
352 /* Loop over enforced rotation groups (usually 1, though)
353 * Apply the forces from rotation potentials */
354 for (g=0; g<rot->ngrp; g++)
359 for (l=0; l<erg->nat_loc; l++)
361 /* Get the right index of the local force */
362 ii = erg->ind_loc[l];
364 rvec_inc(f[ii],erg->f_rot_loc[l]);
368 GMX_MPE_LOG(ev_add_rot_forces_finish);
370 return (MASTER(cr)? er->Vrot : 0.0);
374 /* The Gaussian norm is chosen such that the sum of the gaussian functions
375 * over the slabs is approximately 1.0 everywhere */
376 #define GAUSS_NORM 0.569917543430618
379 /* Calculate the maximum beta that leads to a gaussian larger min_gaussian,
380 * also does some checks
382 static double calc_beta_max(real min_gaussian, real slab_dist)
388 /* Actually the next two checks are already made in grompp */
390 gmx_fatal(FARGS, "Slab distance of flexible rotation groups must be >=0 !");
391 if (min_gaussian <= 0)
392 gmx_fatal(FARGS, "Cutoff value for Gaussian must be > 0. (You requested %f)");
394 /* Define the sigma value */
395 sigma = 0.7*slab_dist;
397 /* Calculate the argument for the logarithm and check that the log() result is negative or 0 */
398 arg = min_gaussian/GAUSS_NORM;
400 gmx_fatal(FARGS, "min_gaussian of flexible rotation groups must be <%g", GAUSS_NORM);
402 return sqrt(-2.0*sigma*sigma*log(min_gaussian/GAUSS_NORM));
406 static inline real calc_beta(rvec curr_x, t_rotgrp *rotg, int n)
408 return iprod(curr_x, rotg->vec) - rotg->slab_dist * n;
412 static inline real gaussian_weight(rvec curr_x, t_rotgrp *rotg, int n)
414 const real norm = GAUSS_NORM;
418 /* Define the sigma value */
419 sigma = 0.7*rotg->slab_dist;
420 /* Calculate the Gaussian value of slab n for position curr_x */
421 return norm * exp( -0.5 * sqr( calc_beta(curr_x, rotg, n)/sigma ) );
425 /* Returns the weight in a single slab, also calculates the Gaussian- and mass-
426 * weighted sum of positions for that slab */
427 static real get_slab_weight(int j, t_rotgrp *rotg, rvec xc[], real mc[], rvec *x_weighted_sum)
429 rvec curr_x; /* The position of an atom */
430 rvec curr_x_weighted; /* The gaussian-weighted position */
431 real gaussian; /* A single gaussian weight */
432 real wgauss; /* gaussian times current mass */
433 real slabweight = 0.0; /* The sum of weights in the slab */
435 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
439 clear_rvec(*x_weighted_sum);
442 islab = j - erg->slab_first;
444 /* Loop over all atoms in the rotation group */
445 for (i=0; i<rotg->nat; i++)
447 copy_rvec(xc[i], curr_x);
448 gaussian = gaussian_weight(curr_x, rotg, j);
449 wgauss = gaussian * mc[i];
450 svmul(wgauss, curr_x, curr_x_weighted);
451 rvec_add(*x_weighted_sum, curr_x_weighted, *x_weighted_sum);
452 slabweight += wgauss;
453 } /* END of loop over rotation group atoms */
459 static void get_slab_centers(
460 t_rotgrp *rotg, /* The rotation group information */
461 rvec *xc, /* The rotation group positions; will
462 typically be enfrotgrp->xc, but at first call
463 it is enfrotgrp->xc_ref */
464 real *mc, /* The masses of the rotation group atoms */
465 t_commrec *cr, /* Communication record */
466 int g, /* The number of the rotation group */
467 real time, /* Used for output only */
468 FILE *out_slabs, /* For outputting center per slab information */
469 gmx_bool bOutStep, /* Is this an output step? */
470 gmx_bool bReference) /* If this routine is called from
471 init_rot_group we need to store
472 the reference slab centers */
475 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
480 /* Loop over slabs */
481 for (j = erg->slab_first; j <= erg->slab_last; j++)
483 islab = j - erg->slab_first;
484 erg->slab_weights[islab] = get_slab_weight(j, rotg, xc, mc, &erg->slab_center[islab]);
486 /* We can do the calculations ONLY if there is weight in the slab! */
487 if (erg->slab_weights[islab] > WEIGHT_MIN)
489 svmul(1.0/erg->slab_weights[islab], erg->slab_center[islab], erg->slab_center[islab]);
493 /* We need to check this here, since we divide through slab_weights
494 * in the flexible low-level routines! */
495 gmx_fatal(FARGS, "Not enough weight in slab %d. Slab center cannot be determined!", j);
498 /* At first time step: save the centers of the reference structure */
500 copy_rvec(erg->slab_center[islab], erg->slab_center_ref[islab]);
501 } /* END of loop over slabs */
503 /* Output on the master */
504 if (MASTER(cr) && bOutStep)
506 fprintf(out_slabs, "%12.3e%6d", time, g);
507 for (j = erg->slab_first; j <= erg->slab_last; j++)
509 islab = j - erg->slab_first;
510 fprintf(out_slabs, "%6d%12.3e%12.3e%12.3e",
511 j,erg->slab_center[islab][XX],erg->slab_center[islab][YY],erg->slab_center[islab][ZZ]);
513 fprintf(out_slabs, "\n");
518 static void calc_rotmat(
520 real degangle, /* Angle alpha of rotation at time t in degrees */
521 matrix rotmat) /* Rotation matrix */
523 real radangle; /* Rotation angle in radians */
524 real cosa; /* cosine alpha */
525 real sina; /* sine alpha */
526 real OMcosa; /* 1 - cos(alpha) */
527 real dumxy, dumxz, dumyz; /* save computations */
528 rvec rot_vec; /* Rotate around rot_vec ... */
531 radangle = degangle * M_PI/180.0;
532 copy_rvec(vec , rot_vec );
534 /* Precompute some variables: */
535 cosa = cos(radangle);
536 sina = sin(radangle);
538 dumxy = rot_vec[XX]*rot_vec[YY]*OMcosa;
539 dumxz = rot_vec[XX]*rot_vec[ZZ]*OMcosa;
540 dumyz = rot_vec[YY]*rot_vec[ZZ]*OMcosa;
542 /* Construct the rotation matrix for this rotation group: */
544 rotmat[XX][XX] = cosa + rot_vec[XX]*rot_vec[XX]*OMcosa;
545 rotmat[YY][XX] = dumxy + rot_vec[ZZ]*sina;
546 rotmat[ZZ][XX] = dumxz - rot_vec[YY]*sina;
548 rotmat[XX][YY] = dumxy - rot_vec[ZZ]*sina;
549 rotmat[YY][YY] = cosa + rot_vec[YY]*rot_vec[YY]*OMcosa;
550 rotmat[ZZ][YY] = dumyz + rot_vec[XX]*sina;
552 rotmat[XX][ZZ] = dumxz + rot_vec[YY]*sina;
553 rotmat[YY][ZZ] = dumyz - rot_vec[XX]*sina;
554 rotmat[ZZ][ZZ] = cosa + rot_vec[ZZ]*rot_vec[ZZ]*OMcosa;
559 for (iii=0; iii<3; iii++) {
560 for (jjj=0; jjj<3; jjj++)
561 fprintf(stderr, " %10.8f ", rotmat[iii][jjj]);
562 fprintf(stderr, "\n");
568 /* Calculates torque on the rotation axis tau = position x force */
569 static inline real torque(
570 rvec rotvec, /* rotation vector; MUST be normalized! */
571 rvec force, /* force */
572 rvec x, /* position of atom on which the force acts */
573 rvec pivot) /* pivot point of rotation axis */
578 /* Subtract offset */
579 rvec_sub(x,pivot,vectmp);
581 /* position x force */
582 cprod(vectmp, force, tau);
584 /* Return the part of the torque which is parallel to the rotation vector */
585 return iprod(tau, rotvec);
589 /* Right-aligned output of value with standard width */
590 static void print_aligned(FILE *fp, char *str)
592 fprintf(fp, "%12s", str);
596 /* Right-aligned output of value with standard short width */
597 static void print_aligned_short(FILE *fp, char *str)
599 fprintf(fp, "%6s", str);
603 static FILE *open_output_file(const char *fn, int steps, const char what[])
608 fp = ffopen(fn, "w");
610 fprintf(fp, "# Output of %s is written in intervals of %d time step%s.\n#\n",
611 what,steps, steps>1 ? "s":"");
617 /* Open output file for slab center data. Call on master only */
618 static FILE *open_slab_out(t_rot *rot, const char *fn)
625 for (g=0; g<rot->ngrp; g++)
631 fp = open_output_file(fn, rot->nstsout, "gaussian weighted slab centers");
632 fprintf(fp, "# Rotation group %d (%s), slab distance %f nm, %s.\n",
633 g, erotg_names[rotg->eType], rotg->slab_dist, rotg->bMassW? "centers of mass":"geometrical centers");
639 fprintf(fp, "# Reference centers are listed first (t=-1).\n");
640 fprintf(fp, "# The following columns have the syntax:\n");
642 print_aligned_short(fp, "t");
643 print_aligned_short(fp, "grp");
644 /* Print ascii legend for the first two entries only ... */
647 print_aligned_short(fp, "slab");
648 print_aligned(fp, "X center");
649 print_aligned(fp, "Y center");
650 print_aligned(fp, "Z center");
652 fprintf(fp, " ...\n");
660 /* Open output file and print some general information about the rotation groups.
661 * Call on master only */
662 static FILE *open_rot_out(const char *fn, t_rot *rot, const output_env_t oenv,
668 const char **setname;
669 char buf[50], buf2[75];
670 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
674 if (Flags & MD_APPENDFILES)
676 fp = gmx_fio_fopen(fn,"a");
680 fp = xvgropen(fn, "Rotation angles and energy", "Time [ps]", "angles [degrees] and energies [kJ/mol]", oenv);
681 fprintf(fp, "# Output of enforced rotation data is written in intervals of %d time step%s.\n#\n", rot->nstrout, rot->nstrout > 1 ? "s":"");
682 fprintf(fp, "# The scalar tau is the torque [kJ/mol] in the direction of the rotation vector v.\n");
683 fprintf(fp, "# To obtain the vectorial torque, multiply tau with the group's rot_vec.\n");
684 fprintf(fp, "# For flexible groups, tau(t,n) from all slabs n have been summed in a single value tau(t) here.\n");
685 fprintf(fp, "# The torques tau(t,n) are found in the rottorque.log (-rt) output file\n");
688 for (g=0; g<rot->ngrp; g++)
692 bFlex = ISFLEX(rotg);
694 fprintf(fp, "# Rotation group %d, potential type '%s':\n" , g, erotg_names[rotg->eType]);
695 fprintf(fp, "# rot_massw%d %s\n" , g, yesno_names[rotg->bMassW]);
696 fprintf(fp, "# rot_vec%d %12.5e %12.5e %12.5e\n" , g, rotg->vec[XX], rotg->vec[YY], rotg->vec[ZZ]);
697 fprintf(fp, "# rot_rate%d %12.5e degrees/ps\n" , g, rotg->rate);
698 fprintf(fp, "# rot_k%d %12.5e kJ/(mol*nm^2)\n" , g, rotg->k);
699 if ( rotg->eType==erotgISO || rotg->eType==erotgPM || rotg->eType==erotgRM || rotg->eType==erotgRM2)
700 fprintf(fp, "# rot_pivot%d %12.5e %12.5e %12.5e nm\n", g, rotg->pivot[XX], rotg->pivot[YY], rotg->pivot[ZZ]);
704 fprintf(fp, "# rot_slab_distance%d %f nm\n", g, rotg->slab_dist);
705 fprintf(fp, "# rot_min_gaussian%d %12.5e\n", g, rotg->min_gaussian);
708 /* Output the centers of the rotation groups for the pivot-free potentials */
709 if ((rotg->eType==erotgISOPF) || (rotg->eType==erotgPMPF) || (rotg->eType==erotgRMPF) || (rotg->eType==erotgRM2PF
710 || (rotg->eType==erotgFLEXT) || (rotg->eType==erotgFLEX2T)) )
712 fprintf(fp, "# %s of ref. grp. %d %12.5e %12.5e %12.5e\n",
713 rotg->bMassW? "COM":"COG", g,
714 erg->xc_ref_center[XX], erg->xc_ref_center[YY], erg->xc_ref_center[ZZ]);
716 fprintf(fp, "# initial %s grp. %d %12.5e %12.5e %12.5e\n",
717 rotg->bMassW? "COM":"COG", g,
718 erg->xc_center[XX], erg->xc_center[YY], erg->xc_center[ZZ]);
721 if ( (rotg->eType == erotgRM2) || (rotg->eType==erotgFLEX2) || (rotg->eType==erotgFLEX2T) )
723 fprintf(fp, "# rot_eps%d %12.5e nm^2\n", g, rotg->eps);
727 fprintf(fp, "#\n# Legend for the following data columns:\n");
729 print_aligned_short(fp, "t");
731 snew(setname, 4*rot->ngrp);
733 for (g=0; g<rot->ngrp; g++)
736 sprintf(buf, "theta_ref%d", g);
737 print_aligned(fp, buf);
738 sprintf(buf2, "%s [degrees]", buf);
739 setname[nsets] = strdup(buf2);
742 for (g=0; g<rot->ngrp; g++)
745 bFlex = ISFLEX(rotg);
747 /* For flexible axis rotation we use RMSD fitting to determine the
748 * actual angle of the rotation group */
750 sprintf(buf, "theta-fit%d", g);
752 sprintf(buf, "theta-av%d", g);
753 print_aligned(fp, buf);
754 sprintf(buf2, "%s [degrees]", buf);
755 setname[nsets] = strdup(buf2);
758 sprintf(buf, "tau%d", g);
759 print_aligned(fp, buf);
760 sprintf(buf2, "%s [kJ/mol]", buf);
761 setname[nsets] = strdup(buf2);
764 sprintf(buf, "energy%d", g);
765 print_aligned(fp, buf);
766 sprintf(buf2, "%s [kJ/mol]", buf);
767 setname[nsets] = strdup(buf2);
770 fprintf(fp, "\n#\n");
773 xvgr_legend(fp, nsets, setname, oenv);
783 /* Call on master only */
784 static FILE *open_angles_out(t_rot *rot, const char *fn)
791 /* Open output file and write some information about it's structure: */
792 fp = open_output_file(fn, rot->nstsout, "rotation group angles");
793 fprintf(fp, "# All angles given in degrees, time in ps.\n");
794 for (g=0; g<rot->ngrp; g++)
799 fprintf(fp, "# Rotation group %d (%s), slab distance %f nm, fit type %s.\n",
800 g, erotg_names[rotg->eType], rotg->slab_dist, erotg_fitnames[rotg->eFittype]);
803 fprintf(fp, "# The following columns will have the syntax:\n");
805 print_aligned_short(fp, "t");
806 print_aligned_short(fp, "grp");
807 print_aligned(fp, "theta_ref");
808 print_aligned_short(fp, "slab");
809 print_aligned_short(fp, "atoms");
810 print_aligned(fp, "theta_fit");
811 print_aligned_short(fp, "slab");
812 print_aligned_short(fp, "atoms");
813 print_aligned(fp, "theta_fit");
814 fprintf(fp, " ...\n");
820 /* Open torque output file and write some information about it's structure.
821 * Call on master only */
822 static FILE *open_torque_out(t_rot *rot, const char *fn)
829 fp = open_output_file(fn, rot->nstsout,"torques");
831 for (g=0; g<rot->ngrp; g++)
836 fprintf(fp, "# Rotation group %d (%s), slab distance %f nm.\n", g, erotg_names[rotg->eType], rotg->slab_dist);
837 fprintf(fp, "# The scalar tau is the torque [kJ/mol] in the direction of the rotation vector.\n");
838 fprintf(fp, "# To obtain the vectorial torque, multiply tau with\n");
839 fprintf(fp, "# rot_vec%d %10.3e %10.3e %10.3e\n", g, rotg->vec[XX], rotg->vec[YY], rotg->vec[ZZ]);
843 fprintf(fp, "# The following columns will have the syntax (tau=torque for that slab):\n");
845 print_aligned_short(fp, "t");
846 print_aligned_short(fp, "grp");
847 print_aligned_short(fp, "slab");
848 print_aligned(fp, "tau");
849 print_aligned_short(fp, "slab");
850 print_aligned(fp, "tau");
851 fprintf(fp, " ...\n");
858 static void swap_val(double* vec, int i, int j)
868 static void swap_col(double **mat, int i, int j)
870 double tmp[3] = {mat[0][j], mat[1][j], mat[2][j]};
883 /* Eigenvectors are stored in columns of eigen_vec */
884 static void diagonalize_symmetric(
892 jacobi(matrix,3,eigenval,eigen_vec,&n_rot);
894 /* sort in ascending order */
895 if (eigenval[0] > eigenval[1])
897 swap_val(eigenval, 0, 1);
898 swap_col(eigen_vec, 0, 1);
900 if (eigenval[1] > eigenval[2])
902 swap_val(eigenval, 1, 2);
903 swap_col(eigen_vec, 1, 2);
905 if (eigenval[0] > eigenval[1])
907 swap_val(eigenval, 0, 1);
908 swap_col(eigen_vec, 0, 1);
913 static void align_with_z(
914 rvec* s, /* Structure to align */
919 rvec zet = {0.0, 0.0, 1.0};
920 rvec rot_axis={0.0, 0.0, 0.0};
921 rvec *rotated_str=NULL;
927 snew(rotated_str, natoms);
929 /* Normalize the axis */
930 ooanorm = 1.0/norm(axis);
931 svmul(ooanorm, axis, axis);
933 /* Calculate the angle for the fitting procedure */
934 cprod(axis, zet, rot_axis);
935 angle = acos(axis[2]);
939 /* Calculate the rotation matrix */
940 calc_rotmat(rot_axis, angle*180.0/M_PI, rotmat);
942 /* Apply the rotation matrix to s */
943 for (i=0; i<natoms; i++)
949 rotated_str[i][j] += rotmat[j][k]*s[i][k];
954 /* Rewrite the rotated structure to s */
955 for(i=0; i<natoms; i++)
959 s[i][j]=rotated_str[i][j];
967 static void calc_correl_matrix(rvec* Xstr, rvec* Ystr, double** Rmat, int natoms)
978 for (k=0; k<natoms; k++)
979 Rmat[i][j] += Ystr[k][i] * Xstr[k][j];
983 static void weigh_coords(rvec* str, real* weight, int natoms)
988 for(i=0; i<natoms; i++)
991 str[i][j] *= sqrt(weight[i]);
996 static real opt_angle_analytic(
1009 double **Rmat, **RtR, **eigvec;
1011 double V[3][3], WS[3][3];
1012 double rot_matrix[3][3];
1016 /* Do not change the original coordinates */
1017 snew(ref_s_1, natoms);
1018 snew(act_s_1, natoms);
1019 for(i=0; i<natoms; i++)
1021 copy_rvec(ref_s[i], ref_s_1[i]);
1022 copy_rvec(act_s[i], act_s_1[i]);
1025 /* Translate the structures to the origin */
1026 shift[XX] = -ref_com[XX];
1027 shift[YY] = -ref_com[YY];
1028 shift[ZZ] = -ref_com[ZZ];
1029 translate_x(ref_s_1, natoms, shift);
1031 shift[XX] = -act_com[XX];
1032 shift[YY] = -act_com[YY];
1033 shift[ZZ] = -act_com[ZZ];
1034 translate_x(act_s_1, natoms, shift);
1036 /* Align rotation axis with z */
1037 align_with_z(ref_s_1, natoms, axis);
1038 align_with_z(act_s_1, natoms, axis);
1040 /* Correlation matrix */
1041 Rmat = allocate_square_matrix(3);
1043 for (i=0; i<natoms; i++)
1049 /* Weight positions with sqrt(weight) */
1052 weigh_coords(ref_s_1, weight, natoms);
1053 weigh_coords(act_s_1, weight, natoms);
1056 /* Calculate correlation matrices R=YXt (X=ref_s; Y=act_s) */
1057 calc_correl_matrix(ref_s_1, act_s_1, Rmat, natoms);
1060 RtR = allocate_square_matrix(3);
1067 RtR[i][j] += Rmat[k][i] * Rmat[k][j];
1071 /* Diagonalize RtR */
1076 diagonalize_symmetric(RtR, eigvec, eigval);
1077 swap_col(eigvec,0,1);
1078 swap_col(eigvec,1,2);
1079 swap_val(eigval,0,1);
1080 swap_val(eigval,1,2);
1094 WS[i][j] = eigvec[i][j] / sqrt(eigval[j]);
1102 V[i][j] += Rmat[i][k]*WS[k][j];
1106 free_square_matrix(Rmat, 3);
1108 /* Calculate optimal rotation matrix */
1111 rot_matrix[i][j] = 0.0;
1118 rot_matrix[i][j] += eigvec[i][k]*V[j][k];
1122 rot_matrix[2][2] = 1.0;
1124 /* In some cases abs(rot_matrix[0][0]) can be slighly larger
1125 * than unity due to numerical inacurracies. To be able to calculate
1126 * the acos function, we put these values back in range. */
1127 if (rot_matrix[0][0] > 1.0)
1129 rot_matrix[0][0] = 1.0;
1131 else if (rot_matrix[0][0] < -1.0)
1133 rot_matrix[0][0] = -1.0;
1136 /* Determine the optimal rotation angle: */
1137 opt_angle = (-1.0)*acos(rot_matrix[0][0])*180.0/M_PI;
1138 if (rot_matrix[0][1] < 0.0)
1139 opt_angle = (-1.0)*opt_angle;
1141 /* Give back some memory */
1142 free_square_matrix(RtR, 3);
1149 return (real) opt_angle;
1153 /* Determine angle of the group by RMSD fit to the reference */
1154 /* Not parallelized, call this routine only on the master */
1155 static real flex_fit_angle(t_rotgrp *rotg)
1158 rvec *fitcoords=NULL;
1159 rvec center; /* Center of positions passed to the fit routine */
1160 real fitangle; /* Angle of the rotation group derived by fitting */
1163 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
1166 erg=rotg->enfrotgrp;
1168 /* Get the center of the rotation group.
1169 * Note, again, erg->xc has been sorted in do_flexible */
1170 get_center(erg->xc, erg->mc_sorted, rotg->nat, center);
1172 /* === Determine the optimal fit angle for the rotation group === */
1173 if (rotg->eFittype == erotgFitNORM)
1175 /* Normalize every position to it's reference length */
1176 for (i=0; i<rotg->nat; i++)
1178 /* Put the center of the positions into the origin */
1179 rvec_sub(erg->xc[i], center, coord);
1180 /* Determine the scaling factor for the length: */
1181 scal = erg->xc_ref_length[erg->xc_sortind[i]] / norm(coord);
1182 /* Get position, multiply with the scaling factor and save */
1183 svmul(scal, coord, erg->xc_norm[i]);
1185 fitcoords = erg->xc_norm;
1189 fitcoords = erg->xc;
1191 /* From the point of view of the current positions, the reference has rotated
1192 * backwards. Since we output the angle relative to the fixed reference,
1193 * we need the minus sign. */
1194 fitangle = -opt_angle_analytic(erg->xc_ref_sorted, fitcoords, erg->mc_sorted,
1195 rotg->nat, erg->xc_ref_center, center, rotg->vec);
1201 /* Determine actual angle of each slab by RMSD fit to the reference */
1202 /* Not parallelized, call this routine only on the master */
1203 static void flex_fit_angle_perslab(
1210 int i,l,n,islab,ind;
1212 rvec act_center; /* Center of actual positions that are passed to the fit routine */
1213 rvec ref_center; /* Same for the reference positions */
1214 real fitangle; /* Angle of a slab derived from an RMSD fit to
1215 * the reference structure at t=0 */
1217 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
1218 real OOm_av; /* 1/average_mass of a rotation group atom */
1219 real m_rel; /* Relative mass of a rotation group atom */
1222 erg=rotg->enfrotgrp;
1224 /* Average mass of a rotation group atom: */
1225 OOm_av = erg->invmass*rotg->nat;
1227 /**********************************/
1228 /* First collect the data we need */
1229 /**********************************/
1231 /* Collect the data for the individual slabs */
1232 for (n = erg->slab_first; n <= erg->slab_last; n++)
1234 islab = n - erg->slab_first; /* slab index */
1235 sd = &(rotg->enfrotgrp->slab_data[islab]);
1236 sd->nat = erg->lastatom[islab]-erg->firstatom[islab]+1;
1239 /* Loop over the relevant atoms in the slab */
1240 for (l=erg->firstatom[islab]; l<=erg->lastatom[islab]; l++)
1242 /* Current position of this atom: x[ii][XX/YY/ZZ] */
1243 copy_rvec(erg->xc[l], curr_x);
1245 /* The (unrotated) reference position of this atom is copied to ref_x.
1246 * Beware, the xc coords have been sorted in do_flexible */
1247 copy_rvec(erg->xc_ref_sorted[l], ref_x);
1249 /* Save data for doing angular RMSD fit later */
1250 /* Save the current atom position */
1251 copy_rvec(curr_x, sd->x[ind]);
1252 /* Save the corresponding reference position */
1253 copy_rvec(ref_x , sd->ref[ind]);
1255 /* Maybe also mass-weighting was requested. If yes, additionally
1256 * multiply the weights with the relative mass of the atom. If not,
1257 * multiply with unity. */
1258 m_rel = erg->mc_sorted[l]*OOm_av;
1260 /* Save the weight for this atom in this slab */
1261 sd->weight[ind] = gaussian_weight(curr_x, rotg, n) * m_rel;
1263 /* Next atom in this slab */
1268 /******************************/
1269 /* Now do the fit calculation */
1270 /******************************/
1272 fprintf(fp, "%12.3e%6d%12.3f", t, g, degangle);
1274 /* === Now do RMSD fitting for each slab === */
1275 /* We require at least SLAB_MIN_ATOMS in a slab, such that the fit makes sense. */
1276 #define SLAB_MIN_ATOMS 4
1278 for (n = erg->slab_first; n <= erg->slab_last; n++)
1280 islab = n - erg->slab_first; /* slab index */
1281 sd = &(rotg->enfrotgrp->slab_data[islab]);
1282 if (sd->nat >= SLAB_MIN_ATOMS)
1284 /* Get the center of the slabs reference and current positions */
1285 get_center(sd->ref, sd->weight, sd->nat, ref_center);
1286 get_center(sd->x , sd->weight, sd->nat, act_center);
1287 if (rotg->eFittype == erotgFitNORM)
1289 /* Normalize every position to it's reference length
1290 * prior to performing the fit */
1291 for (i=0; i<sd->nat;i++) /* Center */
1293 rvec_dec(sd->ref[i], ref_center);
1294 rvec_dec(sd->x[i] , act_center);
1295 /* Normalize x_i such that it gets the same length as ref_i */
1296 svmul( norm(sd->ref[i])/norm(sd->x[i]), sd->x[i], sd->x[i] );
1298 /* We already subtracted the centers */
1299 clear_rvec(ref_center);
1300 clear_rvec(act_center);
1302 fitangle = -opt_angle_analytic(sd->ref, sd->x, sd->weight, sd->nat,
1303 ref_center, act_center, rotg->vec);
1304 fprintf(fp, "%6d%6d%12.3f", n, sd->nat, fitangle);
1309 #undef SLAB_MIN_ATOMS
1313 /* Shift x with is */
1314 static inline void shift_single_coord(matrix box, rvec x, const ivec is)
1325 x[XX] += tx*box[XX][XX]+ty*box[YY][XX]+tz*box[ZZ][XX];
1326 x[YY] += ty*box[YY][YY]+tz*box[ZZ][YY];
1327 x[ZZ] += tz*box[ZZ][ZZ];
1330 x[XX] += tx*box[XX][XX];
1331 x[YY] += ty*box[YY][YY];
1332 x[ZZ] += tz*box[ZZ][ZZ];
1337 /* Determine the 'home' slab of this atom which is the
1338 * slab with the highest Gaussian weight of all */
1339 #define round(a) (int)(a+0.5)
1340 static inline int get_homeslab(
1341 rvec curr_x, /* The position for which the home slab shall be determined */
1342 rvec rotvec, /* The rotation vector */
1343 real slabdist) /* The slab distance */
1348 /* The distance of the atom to the coordinate center (where the
1349 * slab with index 0) is */
1350 dist = iprod(rotvec, curr_x);
1352 return round(dist / slabdist);
1356 /* For a local atom determine the relevant slabs, i.e. slabs in
1357 * which the gaussian is larger than min_gaussian
1359 static int get_single_atom_gaussians(
1367 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
1370 erg=rotg->enfrotgrp;
1372 /* Determine the 'home' slab of this atom: */
1373 homeslab = get_homeslab(curr_x, rotg->vec, rotg->slab_dist);
1375 /* First determine the weight in the atoms home slab: */
1376 g = gaussian_weight(curr_x, rotg, homeslab);
1378 erg->gn_atom[count] = g;
1379 erg->gn_slabind[count] = homeslab;
1383 /* Determine the max slab */
1385 while (g > rotg->min_gaussian)
1388 g = gaussian_weight(curr_x, rotg, slab);
1389 erg->gn_slabind[count]=slab;
1390 erg->gn_atom[count]=g;
1395 /* Determine the max slab */
1400 g = gaussian_weight(curr_x, rotg, slab);
1401 erg->gn_slabind[count]=slab;
1402 erg->gn_atom[count]=g;
1405 while (g > rotg->min_gaussian);
1412 static void flex2_precalc_inner_sum(t_rotgrp *rotg, t_commrec *cr)
1415 rvec xi; /* positions in the i-sum */
1416 rvec xcn, ycn; /* the current and the reference slab centers */
1419 rvec rin; /* Helper variables */
1422 real OOpsii,OOpsiistar;
1423 real sin_rin; /* s_ii.r_ii */
1424 rvec s_in,tmpvec,tmpvec2;
1425 real mi,wi; /* Mass-weighting of the positions */
1427 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
1430 erg=rotg->enfrotgrp;
1431 N_M = rotg->nat * erg->invmass;
1433 /* Loop over all slabs that contain something */
1434 for (n=erg->slab_first; n <= erg->slab_last; n++)
1436 islab = n - erg->slab_first; /* slab index */
1438 /* The current center of this slab is saved in xcn: */
1439 copy_rvec(erg->slab_center[islab], xcn);
1440 /* ... and the reference center in ycn: */
1441 copy_rvec(erg->slab_center_ref[islab+erg->slab_buffer], ycn);
1443 /*** D. Calculate the whole inner sum used for second and third sum */
1444 /* For slab n, we need to loop over all atoms i again. Since we sorted
1445 * the atoms with respect to the rotation vector, we know that it is sufficient
1446 * to calculate from firstatom to lastatom only. All other contributions will
1448 clear_rvec(innersumvec);
1449 for (i = erg->firstatom[islab]; i <= erg->lastatom[islab]; i++)
1451 /* Coordinate xi of this atom */
1452 copy_rvec(erg->xc[i],xi);
1455 gaussian_xi = gaussian_weight(xi,rotg,n);
1456 mi = erg->mc_sorted[i]; /* need the sorted mass here */
1460 copy_rvec(erg->xc_ref_sorted[i],yi0); /* Reference position yi0 */
1461 rvec_sub(yi0, ycn, tmpvec2); /* tmpvec2 = yi0 - ycn */
1462 mvmul(erg->rotmat, tmpvec2, rin); /* rin = Omega.(yi0 - ycn) */
1464 /* Calculate psi_i* and sin */
1465 rvec_sub(xi, xcn, tmpvec2); /* tmpvec2 = xi - xcn */
1466 cprod(rotg->vec, tmpvec2, tmpvec); /* tmpvec = v x (xi - xcn) */
1467 OOpsiistar = norm2(tmpvec)+rotg->eps; /* OOpsii* = 1/psii* = |v x (xi-xcn)|^2 + eps */
1468 OOpsii = norm(tmpvec); /* OOpsii = 1 / psii = |v x (xi - xcn)| */
1470 /* v x (xi - xcn) */
1471 unitv(tmpvec, s_in); /* sin = ---------------- */
1472 /* |v x (xi - xcn)| */
1474 sin_rin=iprod(s_in,rin); /* sin_rin = sin . rin */
1476 /* Now the whole sum */
1477 fac = OOpsii/OOpsiistar;
1478 svmul(fac, rin, tmpvec);
1479 fac2 = fac*fac*OOpsii;
1480 svmul(fac2*sin_rin, s_in, tmpvec2);
1481 rvec_dec(tmpvec, tmpvec2);
1483 svmul(wi*gaussian_xi*sin_rin, tmpvec, tmpvec2);
1485 rvec_inc(innersumvec,tmpvec2);
1486 } /* now we have the inner sum, used both for sum2 and sum3 */
1488 /* Save it to be used in do_flex2_lowlevel */
1489 copy_rvec(innersumvec, erg->slab_innersumvec[islab]);
1490 } /* END of loop over slabs */
1494 static void flex_precalc_inner_sum(t_rotgrp *rotg, t_commrec *cr)
1497 rvec xi; /* position */
1498 rvec xcn, ycn; /* the current and the reference slab centers */
1499 rvec qin,rin; /* q_i^n and r_i^n */
1502 rvec innersumvec; /* Inner part of sum_n2 */
1503 real gaussian_xi; /* Gaussian weight gn(xi) */
1504 real mi,wi; /* Mass-weighting of the positions */
1507 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
1510 erg=rotg->enfrotgrp;
1511 N_M = rotg->nat * erg->invmass;
1513 /* Loop over all slabs that contain something */
1514 for (n=erg->slab_first; n <= erg->slab_last; n++)
1516 islab = n - erg->slab_first; /* slab index */
1518 /* The current center of this slab is saved in xcn: */
1519 copy_rvec(erg->slab_center[islab], xcn);
1520 /* ... and the reference center in ycn: */
1521 copy_rvec(erg->slab_center_ref[islab+erg->slab_buffer], ycn);
1523 /* For slab n, we need to loop over all atoms i again. Since we sorted
1524 * the atoms with respect to the rotation vector, we know that it is sufficient
1525 * to calculate from firstatom to lastatom only. All other contributions will
1527 clear_rvec(innersumvec);
1528 for (i=erg->firstatom[islab]; i<=erg->lastatom[islab]; i++)
1530 /* Coordinate xi of this atom */
1531 copy_rvec(erg->xc[i],xi);
1534 gaussian_xi = gaussian_weight(xi,rotg,n);
1535 mi = erg->mc_sorted[i]; /* need the sorted mass here */
1538 /* Calculate rin and qin */
1539 rvec_sub(erg->xc_ref_sorted[i], ycn, tmpvec); /* tmpvec = yi0-ycn */
1540 mvmul(erg->rotmat, tmpvec, rin); /* rin = Omega.(yi0 - ycn) */
1541 cprod(rotg->vec, rin, tmpvec); /* tmpvec = v x Omega*(yi0-ycn) */
1543 /* v x Omega*(yi0-ycn) */
1544 unitv(tmpvec, qin); /* qin = --------------------- */
1545 /* |v x Omega*(yi0-ycn)| */
1548 rvec_sub(xi, xcn, tmpvec); /* tmpvec = xi-xcn */
1549 bin = iprod(qin, tmpvec); /* bin = qin*(xi-xcn) */
1551 svmul(wi*gaussian_xi*bin, qin, tmpvec);
1553 /* Add this contribution to the inner sum: */
1554 rvec_add(innersumvec, tmpvec, innersumvec);
1555 } /* now we have the inner sum vector S^n for this slab */
1556 /* Save it to be used in do_flex_lowlevel */
1557 copy_rvec(innersumvec, erg->slab_innersumvec[islab]);
1562 static real do_flex2_lowlevel(
1564 real sigma, /* The Gaussian width sigma */
1566 gmx_bool bCalcTorque,
1570 int count,ic,ii,j,m,n,islab,iigrp;
1571 rvec xj; /* position in the i-sum */
1572 rvec yj0; /* the reference position in the j-sum */
1573 rvec xcn, ycn; /* the current and the reference slab centers */
1574 real V; /* This node's part of the rotation pot. energy */
1575 real gaussian_xj; /* Gaussian weight */
1579 rvec rjn; /* Helper variables */
1582 real OOpsij,OOpsijstar;
1583 real OOsigma2; /* 1/(sigma^2) */
1586 rvec sjn,tmpvec,tmpvec2;
1587 rvec sum1vec_part,sum1vec,sum2vec_part,sum2vec,sum3vec,sum4vec,innersumvec;
1589 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
1590 real mj,wj; /* Mass-weighting of the positions */
1592 real Wjn; /* g_n(x_j) m_j / Mjn */
1594 /* To calculate the torque per slab */
1595 rvec slab_force; /* Single force from slab n on one atom */
1596 rvec slab_sum1vec_part;
1597 real slab_sum3part,slab_sum4part;
1598 rvec slab_sum1vec, slab_sum2vec, slab_sum3vec, slab_sum4vec;
1601 erg=rotg->enfrotgrp;
1603 /* Pre-calculate the inner sums, so that we do not have to calculate
1604 * them again for every atom */
1605 flex2_precalc_inner_sum(rotg, cr);
1607 /********************************************************/
1608 /* Main loop over all local atoms of the rotation group */
1609 /********************************************************/
1610 N_M = rotg->nat * erg->invmass;
1612 OOsigma2 = 1.0 / (sigma*sigma);
1613 for (j=0; j<erg->nat_loc; j++)
1615 /* Local index of a rotation group atom */
1616 ii = erg->ind_loc[j];
1617 /* Position of this atom in the collective array */
1618 iigrp = erg->xc_ref_ind[j];
1619 /* Mass-weighting */
1620 mj = erg->mc[iigrp]; /* need the unsorted mass here */
1623 /* Current position of this atom: x[ii][XX/YY/ZZ]
1624 * Note that erg->xc_center contains the center of mass in case the flex2-t
1625 * potential was chosen. For the flex2 potential erg->xc_center must be
1627 rvec_sub(x[ii], erg->xc_center, xj);
1629 /* Shift this atom such that it is near its reference */
1630 shift_single_coord(box, xj, erg->xc_shifts[iigrp]);
1632 /* Determine the slabs to loop over, i.e. the ones with contributions
1633 * larger than min_gaussian */
1634 count = get_single_atom_gaussians(xj, cr, rotg);
1636 clear_rvec(sum1vec_part);
1637 clear_rvec(sum2vec_part);
1640 /* Loop over the relevant slabs for this atom */
1641 for (ic=0; ic < count; ic++)
1643 n = erg->gn_slabind[ic];
1645 /* Get the precomputed Gaussian value of curr_slab for curr_x */
1646 gaussian_xj = erg->gn_atom[ic];
1648 islab = n - erg->slab_first; /* slab index */
1650 /* The (unrotated) reference position of this atom is copied to yj0: */
1651 copy_rvec(rotg->x_ref[iigrp], yj0);
1653 beta = calc_beta(xj, rotg,n);
1655 /* The current center of this slab is saved in xcn: */
1656 copy_rvec(erg->slab_center[islab], xcn);
1657 /* ... and the reference center in ycn: */
1658 copy_rvec(erg->slab_center_ref[islab+erg->slab_buffer], ycn);
1660 rvec_sub(yj0, ycn, tmpvec2); /* tmpvec2 = yj0 - ycn */
1663 mvmul(erg->rotmat, tmpvec2, rjn); /* rjn = Omega.(yj0 - ycn) */
1665 /* Subtract the slab center from xj */
1666 rvec_sub(xj, xcn, tmpvec2); /* tmpvec2 = xj - xcn */
1669 cprod(rotg->vec, tmpvec2, tmpvec); /* tmpvec = v x (xj - xcn) */
1671 OOpsijstar = norm2(tmpvec)+rotg->eps; /* OOpsij* = 1/psij* = |v x (xj-xcn)|^2 + eps */
1673 numerator = sqr(iprod(tmpvec, rjn));
1675 /*********************************/
1676 /* Add to the rotation potential */
1677 /*********************************/
1678 V += 0.5*rotg->k*wj*gaussian_xj*numerator/OOpsijstar;
1681 /*************************************/
1682 /* Now calculate the force on atom j */
1683 /*************************************/
1685 OOpsij = norm(tmpvec); /* OOpsij = 1 / psij = |v x (xj - xcn)| */
1687 /* v x (xj - xcn) */
1688 unitv(tmpvec, sjn); /* sjn = ---------------- */
1689 /* |v x (xj - xcn)| */
1691 sjn_rjn=iprod(sjn,rjn); /* sjn_rjn = sjn . rjn */
1694 /*** A. Calculate the first of the four sum terms: ****************/
1695 fac = OOpsij/OOpsijstar;
1696 svmul(fac, rjn, tmpvec);
1697 fac2 = fac*fac*OOpsij;
1698 svmul(fac2*sjn_rjn, sjn, tmpvec2);
1699 rvec_dec(tmpvec, tmpvec2);
1700 fac2 = wj*gaussian_xj; /* also needed for sum4 */
1701 svmul(fac2*sjn_rjn, tmpvec, slab_sum1vec_part);
1702 /********************/
1703 /*** Add to sum1: ***/
1704 /********************/
1705 rvec_inc(sum1vec_part, slab_sum1vec_part); /* sum1 still needs to vector multiplied with v */
1707 /*** B. Calculate the forth of the four sum terms: ****************/
1708 betasigpsi = beta*OOsigma2*OOpsij; /* this is also needed for sum3 */
1709 /********************/
1710 /*** Add to sum4: ***/
1711 /********************/
1712 slab_sum4part = fac2*betasigpsi*fac*sjn_rjn*sjn_rjn; /* Note that fac is still valid from above */
1713 sum4 += slab_sum4part;
1715 /*** C. Calculate Wjn for second and third sum */
1716 /* Note that we can safely divide by slab_weights since we check in
1717 * get_slab_centers that it is non-zero. */
1718 Wjn = gaussian_xj*mj/erg->slab_weights[islab];
1720 /* We already have precalculated the inner sum for slab n */
1721 copy_rvec(erg->slab_innersumvec[islab], innersumvec);
1723 /* Weigh the inner sum vector with Wjn */
1724 svmul(Wjn, innersumvec, innersumvec);
1726 /*** E. Calculate the second of the four sum terms: */
1727 /********************/
1728 /*** Add to sum2: ***/
1729 /********************/
1730 rvec_inc(sum2vec_part, innersumvec); /* sum2 still needs to be vector crossproduct'ed with v */
1732 /*** F. Calculate the third of the four sum terms: */
1733 slab_sum3part = betasigpsi * iprod(sjn, innersumvec);
1734 sum3 += slab_sum3part; /* still needs to be multiplied with v */
1736 /*** G. Calculate the torque on the local slab's axis: */
1740 cprod(slab_sum1vec_part, rotg->vec, slab_sum1vec);
1742 cprod(innersumvec, rotg->vec, slab_sum2vec);
1744 svmul(slab_sum3part, rotg->vec, slab_sum3vec);
1746 svmul(slab_sum4part, rotg->vec, slab_sum4vec);
1748 /* The force on atom ii from slab n only: */
1749 for (m=0; m<DIM; m++)
1750 slab_force[m] = rotg->k * (-slab_sum1vec[m] + slab_sum2vec[m] - slab_sum3vec[m] + 0.5*slab_sum4vec[m]);
1752 erg->slab_torque_v[islab] += torque(rotg->vec, slab_force, xj, xcn);
1754 } /* END of loop over slabs */
1756 /* Construct the four individual parts of the vector sum: */
1757 cprod(sum1vec_part, rotg->vec, sum1vec); /* sum1vec = { } x v */
1758 cprod(sum2vec_part, rotg->vec, sum2vec); /* sum2vec = { } x v */
1759 svmul(sum3, rotg->vec, sum3vec); /* sum3vec = { } . v */
1760 svmul(sum4, rotg->vec, sum4vec); /* sum4vec = { } . v */
1762 /* Store the additional force so that it can be added to the force
1763 * array after the normal forces have been evaluated */
1764 for (m=0; m<DIM; m++)
1765 erg->f_rot_loc[j][m] = rotg->k * (-sum1vec[m] + sum2vec[m] - sum3vec[m] + 0.5*sum4vec[m]);
1768 fprintf(stderr, "sum1: %15.8f %15.8f %15.8f\n", -rotg->k*sum1vec[XX], -rotg->k*sum1vec[YY], -rotg->k*sum1vec[ZZ]);
1769 fprintf(stderr, "sum2: %15.8f %15.8f %15.8f\n", rotg->k*sum2vec[XX], rotg->k*sum2vec[YY], rotg->k*sum2vec[ZZ]);
1770 fprintf(stderr, "sum3: %15.8f %15.8f %15.8f\n", -rotg->k*sum3vec[XX], -rotg->k*sum3vec[YY], -rotg->k*sum3vec[ZZ]);
1771 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]);
1776 } /* END of loop over local atoms */
1782 static real do_flex_lowlevel(
1784 real sigma, /* The Gaussian width sigma */
1786 gmx_bool bCalcTorque,
1790 int count,ic,ii,j,m,n,islab,iigrp;
1791 rvec xj,yj0; /* current and reference position */
1792 rvec xcn, ycn; /* the current and the reference slab centers */
1793 rvec xj_xcn; /* xj - xcn */
1794 rvec qjn; /* q_i^n */
1795 rvec sum_n1,sum_n2; /* Two contributions to the rotation force */
1796 rvec innersumvec; /* Inner part of sum_n2 */
1798 rvec force_n; /* Single force from slab n on one atom */
1799 rvec tmpvec,tmpvec2,tmp_f; /* Helper variables */
1800 real V; /* The rotation potential energy */
1801 real OOsigma2; /* 1/(sigma^2) */
1802 real beta; /* beta_n(xj) */
1803 real bjn; /* b_j^n */
1804 real gaussian_xj; /* Gaussian weight gn(xj) */
1805 real betan_xj_sigma2;
1806 real mj,wj; /* Mass-weighting of the positions */
1808 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
1811 erg=rotg->enfrotgrp;
1813 /* Pre-calculate the inner sums, so that we do not have to calculate
1814 * them again for every atom */
1815 flex_precalc_inner_sum(rotg, cr);
1817 /********************************************************/
1818 /* Main loop over all local atoms of the rotation group */
1819 /********************************************************/
1820 OOsigma2 = 1.0/(sigma*sigma);
1821 N_M = rotg->nat * erg->invmass;
1823 for (j=0; j<erg->nat_loc; j++)
1825 /* Local index of a rotation group atom */
1826 ii = erg->ind_loc[j];
1827 /* Position of this atom in the collective array */
1828 iigrp = erg->xc_ref_ind[j];
1829 /* Mass-weighting */
1830 mj = erg->mc[iigrp]; /* need the unsorted mass here */
1833 /* Current position of this atom: x[ii][XX/YY/ZZ]
1834 * Note that erg->xc_center contains the center of mass in case the flex-t
1835 * potential was chosen. For the flex potential erg->xc_center must be
1837 rvec_sub(x[ii], erg->xc_center, xj);
1839 /* Shift this atom such that it is near its reference */
1840 shift_single_coord(box, xj, erg->xc_shifts[iigrp]);
1842 /* Determine the slabs to loop over, i.e. the ones with contributions
1843 * larger than min_gaussian */
1844 count = get_single_atom_gaussians(xj, cr, rotg);
1849 /* Loop over the relevant slabs for this atom */
1850 for (ic=0; ic < count; ic++)
1852 n = erg->gn_slabind[ic];
1854 /* Get the precomputed Gaussian for xj in slab n */
1855 gaussian_xj = erg->gn_atom[ic];
1857 islab = n - erg->slab_first; /* slab index */
1859 /* The (unrotated) reference position of this atom is saved in yj0: */
1860 copy_rvec(rotg->x_ref[iigrp], yj0);
1862 beta = calc_beta(xj, rotg, n);
1864 /* The current center of this slab is saved in xcn: */
1865 copy_rvec(erg->slab_center[islab], xcn);
1866 /* ... and the reference center in ycn: */
1867 copy_rvec(erg->slab_center_ref[islab+erg->slab_buffer], ycn);
1869 rvec_sub(yj0, ycn, tmpvec); /* tmpvec = yj0 - ycn */
1872 mvmul(erg->rotmat, tmpvec, tmpvec2); /* tmpvec2 = Omega.(yj0-ycn) */
1874 /* Subtract the slab center from xj */
1875 rvec_sub(xj, xcn, xj_xcn); /* xj_xcn = xj - xcn */
1878 cprod(rotg->vec, tmpvec2, tmpvec); /* tmpvec = v x Omega.(xj-xcn) */
1880 /* v x Omega.(xj-xcn) */
1881 unitv(tmpvec,qjn); /* qjn = -------------------- */
1882 /* |v x Omega.(xj-xcn)| */
1884 bjn = iprod(qjn, xj_xcn); /* bjn = qjn * (xj - xcn) */
1886 /*********************************/
1887 /* Add to the rotation potential */
1888 /*********************************/
1889 V += 0.5*rotg->k*wj*gaussian_xj*sqr(bjn);
1891 /****************************************************************/
1892 /* sum_n1 will typically be the main contribution to the force: */
1893 /****************************************************************/
1894 betan_xj_sigma2 = beta*OOsigma2; /* beta_n(xj)/sigma^2 */
1896 /* The next lines calculate
1897 * qjn - (bjn*beta(xj)/(2sigma^2))v */
1898 svmul(bjn*0.5*betan_xj_sigma2, rotg->vec, tmpvec2);
1899 rvec_sub(qjn,tmpvec2,tmpvec);
1901 /* Multiply with gn(xj)*bjn: */
1902 svmul(gaussian_xj*bjn,tmpvec,tmpvec2);
1905 rvec_inc(sum_n1,tmpvec2);
1907 /* We already have precalculated the Sn term for slab n */
1908 copy_rvec(erg->slab_innersumvec[islab], s_n);
1910 svmul(betan_xj_sigma2*iprod(s_n, xj_xcn), rotg->vec, tmpvec); /* tmpvec = ---------- s_n (xj-xcn) */
1913 rvec_sub(s_n, tmpvec, innersumvec);
1915 /* We can safely divide by slab_weights since we check in get_slab_centers
1916 * that it is non-zero. */
1917 svmul(gaussian_xj/erg->slab_weights[islab], innersumvec, innersumvec);
1919 rvec_add(sum_n2, innersumvec, sum_n2);
1921 GMX_MPE_LOG(ev_inner_loop_finish);
1923 /* Calculate the torque: */
1926 /* The force on atom ii from slab n only: */
1927 rvec_sub(innersumvec, tmpvec2, force_n);
1928 svmul(rotg->k, force_n, force_n);
1929 erg->slab_torque_v[islab] += torque(rotg->vec, force_n, xj, xcn);
1931 } /* END of loop over slabs */
1933 /* Put both contributions together: */
1934 svmul(wj, sum_n1, sum_n1);
1935 svmul(mj, sum_n2, sum_n2);
1936 rvec_sub(sum_n2,sum_n1,tmp_f); /* F = -grad V */
1938 /* Store the additional force so that it can be added to the force
1939 * array after the normal forces have been evaluated */
1940 for(m=0; m<DIM; m++)
1941 erg->f_rot_loc[j][m] = rotg->k*tmp_f[m];
1945 } /* END of loop over local atoms */
1951 static void print_coordinates(t_commrec *cr, t_rotgrp *rotg, rvec x[], matrix box, int step)
1955 static char buf[STRLEN];
1956 static gmx_bool bFirst=1;
1961 sprintf(buf, "coords%d.txt", cr->nodeid);
1962 fp = fopen(buf, "w");
1966 fprintf(fp, "\nStep %d\n", step);
1967 fprintf(fp, "box: %f %f %f %f %f %f %f %f %f\n",
1968 box[XX][XX], box[XX][YY], box[XX][ZZ],
1969 box[YY][XX], box[YY][YY], box[YY][ZZ],
1970 box[ZZ][XX], box[ZZ][ZZ], box[ZZ][ZZ]);
1971 for (i=0; i<rotg->nat; i++)
1973 fprintf(fp, "%4d %f %f %f\n", i,
1974 erg->xc[i][XX], erg->xc[i][YY], erg->xc[i][ZZ]);
1982 static int projection_compare(const void *a, const void *b)
1984 sort_along_vec_t *xca, *xcb;
1987 xca = (sort_along_vec_t *)a;
1988 xcb = (sort_along_vec_t *)b;
1990 if (xca->xcproj < xcb->xcproj)
1992 else if (xca->xcproj > xcb->xcproj)
1999 static void sort_collective_coordinates(
2000 t_rotgrp *rotg, /* Rotation group */
2001 sort_along_vec_t *data) /* Buffer for sorting the positions */
2004 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
2007 erg=rotg->enfrotgrp;
2009 /* The projection of the position vector on the rotation vector is
2010 * the relevant value for sorting. Fill the 'data' structure */
2011 for (i=0; i<rotg->nat; i++)
2013 data[i].xcproj = iprod(erg->xc[i], rotg->vec); /* sort criterium */
2014 data[i].m = erg->mc[i];
2016 copy_rvec(erg->xc[i] , data[i].x );
2017 copy_rvec(rotg->x_ref[i], data[i].x_ref);
2019 /* Sort the 'data' structure */
2020 gmx_qsort(data, rotg->nat, sizeof(sort_along_vec_t), projection_compare);
2022 /* Copy back the sorted values */
2023 for (i=0; i<rotg->nat; i++)
2025 copy_rvec(data[i].x , erg->xc[i] );
2026 copy_rvec(data[i].x_ref, erg->xc_ref_sorted[i]);
2027 erg->mc_sorted[i] = data[i].m;
2028 erg->xc_sortind[i] = data[i].ind;
2033 /* For each slab, get the first and the last index of the sorted atom
2035 static void get_firstlast_atom_per_slab(t_rotgrp *rotg, t_commrec *cr)
2039 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
2042 erg=rotg->enfrotgrp;
2044 GMX_MPE_LOG(ev_get_firstlast_start);
2046 /* Find the first atom that needs to enter the calculation for each slab */
2047 n = erg->slab_first; /* slab */
2048 i = 0; /* start with the first atom */
2051 /* Find the first atom that significantly contributes to this slab */
2052 do /* move forward in position until a large enough beta is found */
2054 beta = calc_beta(erg->xc[i], rotg, n);
2056 } while ((beta < -erg->max_beta) && (i < rotg->nat));
2058 islab = n - erg->slab_first; /* slab index */
2059 erg->firstatom[islab] = i;
2060 /* Proceed to the next slab */
2062 } while (n <= erg->slab_last);
2064 /* Find the last atom for each slab */
2065 n = erg->slab_last; /* start with last slab */
2066 i = rotg->nat-1; /* start with the last atom */
2069 do /* move backward in position until a large enough beta is found */
2071 beta = calc_beta(erg->xc[i], rotg, n);
2073 } while ((beta > erg->max_beta) && (i > -1));
2075 islab = n - erg->slab_first; /* slab index */
2076 erg->lastatom[islab] = i;
2077 /* Proceed to the next slab */
2079 } while (n >= erg->slab_first);
2081 GMX_MPE_LOG(ev_get_firstlast_finish);
2085 /* Determine the very first and very last slab that needs to be considered
2086 * For the first slab that needs to be considered, we have to find the smallest
2089 * x_first * v - n*Delta_x <= beta_max
2091 * slab index n, slab distance Delta_x, rotation vector v. For the last slab we
2092 * have to find the largest n that obeys
2094 * x_last * v - n*Delta_x >= -beta_max
2097 static inline int get_first_slab(
2098 t_rotgrp *rotg, /* The rotation group (inputrec data) */
2099 real max_beta, /* The max_beta value, instead of min_gaussian */
2100 rvec firstatom) /* First atom after sorting along the rotation vector v */
2102 /* Find the first slab for the first atom */
2103 return ceil((iprod(firstatom, rotg->vec) - max_beta)/rotg->slab_dist);
2107 static inline int get_last_slab(
2108 t_rotgrp *rotg, /* The rotation group (inputrec data) */
2109 real max_beta, /* The max_beta value, instead of min_gaussian */
2110 rvec lastatom) /* Last atom along v */
2112 /* Find the last slab for the last atom */
2113 return floor((iprod(lastatom, rotg->vec) + max_beta)/rotg->slab_dist);
2117 static void get_firstlast_slab_check(
2118 t_rotgrp *rotg, /* The rotation group (inputrec data) */
2119 t_gmx_enfrotgrp *erg, /* The rotation group (data only accessible in this file) */
2120 rvec firstatom, /* First atom after sorting along the rotation vector v */
2121 rvec lastatom, /* Last atom along v */
2122 int g, /* The rotation group number */
2125 erg->slab_first = get_first_slab(rotg, erg->max_beta, firstatom);
2126 erg->slab_last = get_last_slab(rotg, erg->max_beta, lastatom);
2128 /* Check whether we have reference data to compare against */
2129 if (erg->slab_first < erg->slab_first_ref)
2130 gmx_fatal(FARGS, "%s No reference data for first slab (n=%d), unable to proceed.",
2131 RotStr, erg->slab_first);
2133 /* Check whether we have reference data to compare against */
2134 if (erg->slab_last > erg->slab_last_ref)
2135 gmx_fatal(FARGS, "%s No reference data for last slab (n=%d), unable to proceed.",
2136 RotStr, erg->slab_last);
2140 /* Enforced rotation with a flexible axis */
2141 static void do_flexible(
2143 gmx_enfrot_t enfrot, /* Other rotation data */
2144 t_rotgrp *rotg, /* The rotation group */
2145 int g, /* Group number */
2146 rvec x[], /* The local positions */
2148 double t, /* Time in picoseconds */
2149 gmx_large_int_t step, /* The time step */
2150 gmx_bool bOutstepRot, /* Output to main rotation output file */
2151 gmx_bool bOutstepSlab) /* Output per-slab data */
2154 real sigma; /* The Gaussian width sigma */
2155 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
2158 erg=rotg->enfrotgrp;
2160 /* Define the sigma value */
2161 sigma = 0.7*rotg->slab_dist;
2163 /* Sort the collective coordinates erg->xc along the rotation vector. This is
2164 * an optimization for the inner loop. */
2165 sort_collective_coordinates(rotg, enfrot->data);
2167 /* Determine the first relevant slab for the first atom and the last
2168 * relevant slab for the last atom */
2169 get_firstlast_slab_check(rotg, erg, erg->xc[0], erg->xc[rotg->nat-1], g, cr);
2171 /* Determine for each slab depending on the min_gaussian cutoff criterium,
2172 * a first and a last atom index inbetween stuff needs to be calculated */
2173 get_firstlast_atom_per_slab(rotg, cr);
2175 /* Determine the gaussian-weighted center of positions for all slabs */
2176 get_slab_centers(rotg,erg->xc,erg->mc_sorted,cr,g,t,enfrot->out_slabs,bOutstepSlab,FALSE);
2178 /* Clear the torque per slab from last time step: */
2179 nslabs = erg->slab_last - erg->slab_first + 1;
2180 for (l=0; l<nslabs; l++)
2181 erg->slab_torque_v[l] = 0.0;
2183 /* Call the rotational forces kernel */
2184 GMX_MPE_LOG(ev_flexll_start);
2185 if (rotg->eType == erotgFLEX || rotg->eType == erotgFLEXT)
2186 erg->V = do_flex_lowlevel(rotg, sigma, x, bOutstepRot, box, cr);
2187 else if (rotg->eType == erotgFLEX2 || rotg->eType == erotgFLEX2T)
2188 erg->V = do_flex2_lowlevel(rotg, sigma, x, bOutstepRot, box, cr);
2190 gmx_fatal(FARGS, "Unknown flexible rotation type");
2191 GMX_MPE_LOG(ev_flexll_finish);
2193 /* Determine angle by RMSD fit to the reference - Let's hope this */
2194 /* only happens once in a while, since this is not parallelized! */
2199 /* Fit angle of the whole rotation group */
2200 erg->angle_v = flex_fit_angle(rotg);
2204 /* Fit angle of each slab */
2205 flex_fit_angle_perslab(g, rotg, t, erg->degangle, enfrot->out_angles);
2209 /* Lump together the torques from all slabs: */
2210 erg->torque_v = 0.0;
2211 for (l=0; l<nslabs; l++)
2212 erg->torque_v += erg->slab_torque_v[l];
2218 /* Calculate the angle between reference and actual rotation group atom,
2219 * both projected into a plane perpendicular to the rotation vector: */
2220 static void angle(t_rotgrp *rotg,
2224 real *weight) /* atoms near the rotation axis should count less than atoms far away */
2226 rvec xp, xrp; /* current and reference positions projected on a plane perpendicular to pg->vec */
2230 /* Project x_ref and x into a plane through the origin perpendicular to rot_vec: */
2231 /* Project x_ref: xrp = x_ref - (vec * x_ref) * vec */
2232 svmul(iprod(rotg->vec, x_ref), rotg->vec, dum);
2233 rvec_sub(x_ref, dum, xrp);
2234 /* Project x_act: */
2235 svmul(iprod(rotg->vec, x_act), rotg->vec, dum);
2236 rvec_sub(x_act, dum, xp);
2238 /* Retrieve information about which vector precedes. gmx_angle always
2239 * returns a positive angle. */
2240 cprod(xp, xrp, dum); /* if reference precedes, this is pointing into the same direction as vec */
2242 if (iprod(rotg->vec, dum) >= 0)
2243 *alpha = -gmx_angle(xrp, xp);
2245 *alpha = +gmx_angle(xrp, xp);
2247 /* Also return the weight */
2252 /* Project first vector onto a plane perpendicular to the second vector
2254 * Note that v must be of unit length.
2256 static inline void project_onto_plane(rvec dr, const rvec v)
2261 svmul(iprod(dr,v),v,tmp); /* tmp = (dr.v)v */
2262 rvec_dec(dr, tmp); /* dr = dr - (dr.v)v */
2266 /* Fixed rotation: The rotation reference group rotates around an axis */
2267 /* The atoms of the actual rotation group are attached with imaginary */
2268 /* springs to the reference atoms. */
2269 static void do_fixed(
2271 t_rotgrp *rotg, /* The rotation group */
2272 rvec x[], /* The positions */
2273 matrix box, /* The simulation box */
2274 double t, /* Time in picoseconds */
2275 gmx_large_int_t step, /* The time step */
2280 rvec tmp_f; /* Force */
2281 real alpha; /* a single angle between an actual and a reference position */
2282 real weight; /* single weight for a single angle */
2283 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
2286 /* for mass weighting: */
2287 real wi; /* Mass-weighting of the positions */
2289 real k_wi; /* k times wi */
2294 erg=rotg->enfrotgrp;
2295 bProject = (rotg->eType==erotgPM) || (rotg->eType==erotgPMPF);
2297 N_M = rotg->nat * erg->invmass;
2299 /* Each process calculates the forces on its local atoms */
2300 for (j=0; j<erg->nat_loc; j++)
2302 /* Calculate (x_i-x_c) resp. (x_i-u) */
2303 rvec_sub(erg->x_loc_pbc[j], erg->xc_center, tmpvec);
2305 /* Calculate Omega*(y_i-y_c)-(x_i-x_c) */
2306 rvec_sub(erg->xr_loc[j], tmpvec, dr);
2309 project_onto_plane(dr, rotg->vec);
2311 /* Mass-weighting */
2312 wi = N_M*erg->m_loc[j];
2314 /* Store the additional force so that it can be added to the force
2315 * array after the normal forces have been evaluated */
2317 for (m=0; m<DIM; m++)
2319 tmp_f[m] = k_wi*dr[m];
2320 erg->f_rot_loc[j][m] = tmp_f[m];
2321 erg->V += 0.5*k_wi*sqr(dr[m]);
2326 /* Add to the torque of this rotation group */
2327 erg->torque_v += torque(rotg->vec, tmp_f, erg->x_loc_pbc[j], erg->xc_center);
2329 /* Calculate the angle between reference and actual rotation group atom. */
2330 angle(rotg, tmpvec, erg->xr_loc[j], &alpha, &weight); /* angle in rad, weighted */
2331 erg->angle_v += alpha * weight;
2332 erg->weight_v += weight;
2334 /* If you want enforced rotation to contribute to the virial,
2335 * activate the following lines:
2338 Add the rotation contribution to the virial
2339 for(j=0; j<DIM; j++)
2341 vir[j][m] += 0.5*f[ii][j]*dr[m];
2347 } /* end of loop over local rotation group atoms */
2353 /* Calculate the radial motion potential and forces */
2354 static void do_radial_motion(
2356 t_rotgrp *rotg, /* The rotation group */
2357 rvec x[], /* The positions */
2358 matrix box, /* The simulation box */
2359 double t, /* Time in picoseconds */
2360 gmx_large_int_t step, /* The time step */
2364 rvec tmp_f; /* Force */
2365 real alpha; /* a single angle between an actual and a reference position */
2366 real weight; /* single weight for a single angle */
2367 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
2368 rvec xj_u; /* xj - u */
2370 real fac,fac2,sum=0.0;
2373 /* For mass weighting: */
2374 real wj; /* Mass-weighting of the positions */
2378 erg=rotg->enfrotgrp;
2380 N_M = rotg->nat * erg->invmass;
2382 /* Each process calculates the forces on its local atoms */
2383 for (j=0; j<erg->nat_loc; j++)
2385 /* Calculate (xj-u) */
2386 rvec_sub(erg->x_loc_pbc[j], erg->xc_center, xj_u); /* xj_u = xj-u */
2388 /* Calculate Omega.(yj-u) */
2389 cprod(rotg->vec, erg->xr_loc[j], tmpvec); /* tmpvec = v x Omega.(yj-u) */
2391 /* v x Omega.(yj-u) */
2392 unitv(tmpvec, pj); /* pj = -------------------- */
2393 /* | v x Omega.(yj-u) | */
2395 fac = iprod(pj, xj_u); /* fac = pj.(xj-u) */
2398 /* Mass-weighting */
2399 wj = N_M*erg->m_loc[j];
2401 /* Store the additional force so that it can be added to the force
2402 * array after the normal forces have been evaluated */
2403 svmul(-rotg->k*wj*fac, pj, tmp_f);
2404 copy_rvec(tmp_f, erg->f_rot_loc[j]);
2408 /* Add to the torque of this rotation group */
2409 erg->torque_v += torque(rotg->vec, tmp_f, erg->x_loc_pbc[j], erg->xc_center);
2411 /* Calculate the angle between reference and actual rotation group atom. */
2412 angle(rotg, xj_u, erg->xr_loc[j], &alpha, &weight); /* angle in rad, weighted */
2413 erg->angle_v += alpha * weight;
2414 erg->weight_v += weight;
2419 } /* end of loop over local rotation group atoms */
2420 erg->V = 0.5*rotg->k*sum;
2427 /* Calculate the radial motion pivot-free potential and forces */
2428 static void do_radial_motion_pf(
2430 t_rotgrp *rotg, /* The rotation group */
2431 rvec x[], /* The positions */
2432 matrix box, /* The simulation box */
2433 double t, /* Time in picoseconds */
2434 gmx_large_int_t step, /* The time step */
2438 rvec xj; /* Current position */
2439 rvec xj_xc; /* xj - xc */
2440 rvec yj0_yc0; /* yj0 - yc0 */
2441 rvec tmp_f; /* Force */
2442 real alpha; /* a single angle between an actual and a reference position */
2443 real weight; /* single weight for a single angle */
2444 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
2445 rvec tmpvec, tmpvec2;
2446 rvec innersumvec; /* Precalculation of the inner sum */
2448 real fac,fac2,V=0.0;
2451 /* For mass weighting: */
2452 real mj,wi,wj; /* Mass-weighting of the positions */
2456 erg=rotg->enfrotgrp;
2458 N_M = rotg->nat * erg->invmass;
2460 /* Get the current center of the rotation group: */
2461 get_center(erg->xc, erg->mc, rotg->nat, erg->xc_center);
2463 /* Precalculate Sum_i [ wi qi.(xi-xc) qi ] which is needed for every single j */
2464 clear_rvec(innersumvec);
2465 for (i=0; i < rotg->nat; i++)
2467 /* Mass-weighting */
2468 wi = N_M*erg->mc[i];
2470 /* Calculate qi. Note that xc_ref_center has already been subtracted from
2471 * x_ref in init_rot_group.*/
2472 mvmul(erg->rotmat, rotg->x_ref[i], tmpvec); /* tmpvec = Omega.(yi0-yc0) */
2474 cprod(rotg->vec, tmpvec, tmpvec2); /* tmpvec2 = v x Omega.(yi0-yc0) */
2476 /* v x Omega.(yi0-yc0) */
2477 unitv(tmpvec2, qi); /* qi = ----------------------- */
2478 /* | v x Omega.(yi0-yc0) | */
2480 rvec_sub(erg->xc[i], erg->xc_center, tmpvec); /* tmpvec = xi-xc */
2482 svmul(wi*iprod(qi, tmpvec), qi, tmpvec2);
2484 rvec_inc(innersumvec, tmpvec2);
2486 svmul(rotg->k*erg->invmass, innersumvec, innersumveckM);
2488 /* Each process calculates the forces on its local atoms */
2489 for (j=0; j<erg->nat_loc; j++)
2491 /* Local index of a rotation group atom */
2492 ii = erg->ind_loc[j];
2493 /* Position of this atom in the collective array */
2494 iigrp = erg->xc_ref_ind[j];
2495 /* Mass-weighting */
2496 mj = erg->mc[iigrp]; /* need the unsorted mass here */
2499 /* Current position of this atom: x[ii][XX/YY/ZZ] */
2500 copy_rvec(x[ii], xj);
2502 /* Shift this atom such that it is near its reference */
2503 shift_single_coord(box, xj, erg->xc_shifts[iigrp]);
2505 /* The (unrotated) reference position is yj0. yc0 has already
2506 * been subtracted in init_rot_group */
2507 copy_rvec(rotg->x_ref[iigrp], yj0_yc0); /* yj0_yc0 = yj0 - yc0 */
2509 /* Calculate Omega.(yj0-yc0) */
2510 mvmul(erg->rotmat, yj0_yc0, tmpvec2); /* tmpvec2 = Omega.(yj0 - yc0) */
2512 cprod(rotg->vec, tmpvec2, tmpvec); /* tmpvec = v x Omega.(yj0-yc0) */
2514 /* v x Omega.(yj0-yc0) */
2515 unitv(tmpvec, qj); /* qj = ----------------------- */
2516 /* | v x Omega.(yj0-yc0) | */
2518 /* Calculate (xj-xc) */
2519 rvec_sub(xj, erg->xc_center, xj_xc); /* xj_xc = xj-xc */
2521 fac = iprod(qj, xj_xc); /* fac = qj.(xj-xc) */
2524 /* Store the additional force so that it can be added to the force
2525 * array after the normal forces have been evaluated */
2526 svmul(-rotg->k*wj*fac, qj, tmp_f); /* part 1 of force */
2527 svmul(mj, innersumveckM, tmpvec); /* part 2 of force */
2528 rvec_inc(tmp_f, tmpvec);
2529 copy_rvec(tmp_f, erg->f_rot_loc[j]);
2533 /* Add to the torque of this rotation group */
2534 erg->torque_v += torque(rotg->vec, tmp_f, xj, erg->xc_center);
2536 /* Calculate the angle between reference and actual rotation group atom. */
2537 angle(rotg, xj_xc, yj0_yc0, &alpha, &weight); /* angle in rad, weighted */
2538 erg->angle_v += alpha * weight;
2539 erg->weight_v += weight;
2544 } /* end of loop over local rotation group atoms */
2545 erg->V = 0.5*rotg->k*V;
2551 /* Precalculate the inner sum for the radial motion 2 forces */
2552 static void radial_motion2_precalc_inner_sum(t_rotgrp *rotg, rvec innersumvec)
2555 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
2556 rvec xi_xc; /* xj - xc */
2557 rvec tmpvec,tmpvec2;
2561 rvec v_xi_xc; /* v x (xj - u) */
2563 real wi; /* Mass-weighting of the positions */
2567 erg=rotg->enfrotgrp;
2568 N_M = rotg->nat * erg->invmass;
2570 /* Loop over the collective set of positions */
2572 for (i=0; i<rotg->nat; i++)
2574 /* Mass-weighting */
2575 wi = N_M*erg->mc[i];
2577 rvec_sub(erg->xc[i], erg->xc_center, xi_xc); /* xi_xc = xi-xc */
2579 /* Calculate ri. Note that xc_ref_center has already been subtracted from
2580 * x_ref in init_rot_group.*/
2581 mvmul(erg->rotmat, rotg->x_ref[i], ri); /* ri = Omega.(yi0-yc0) */
2583 cprod(rotg->vec, xi_xc, v_xi_xc); /* v_xi_xc = v x (xi-u) */
2585 fac = norm2(v_xi_xc);
2587 psiistar = 1.0/(fac + rotg->eps); /* psiistar = --------------------- */
2588 /* |v x (xi-xc)|^2 + eps */
2590 psii = gmx_invsqrt(fac); /* 1 */
2591 /* psii = ------------- */
2594 svmul(psii, v_xi_xc, si); /* si = psii * (v x (xi-xc) ) */
2596 fac = iprod(v_xi_xc, ri); /* fac = (v x (xi-xc)).ri */
2599 siri = iprod(si, ri); /* siri = si.ri */
2601 svmul(psiistar/psii, ri, tmpvec);
2602 svmul(psiistar*psiistar/(psii*psii*psii) * siri, si, tmpvec2);
2603 rvec_dec(tmpvec, tmpvec2);
2604 cprod(tmpvec, rotg->vec, tmpvec2);
2606 svmul(wi*siri, tmpvec2, tmpvec);
2608 rvec_inc(sumvec, tmpvec);
2610 svmul(rotg->k*erg->invmass, sumvec, innersumvec);
2614 /* Calculate the radial motion 2 potential and forces */
2615 static void do_radial_motion2(
2617 t_rotgrp *rotg, /* The rotation group */
2618 rvec x[], /* The positions */
2619 matrix box, /* The simulation box */
2620 double t, /* Time in picoseconds */
2621 gmx_large_int_t step, /* The time step */
2625 rvec xj; /* Position */
2626 real alpha; /* a single angle between an actual and a reference position */
2627 real weight; /* single weight for a single angle */
2628 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
2629 rvec xj_u; /* xj - u */
2630 rvec tmpvec,tmpvec2;
2631 real fac,fac2,Vpart=0.0;
2634 rvec v_xj_u; /* v x (xj - u) */
2636 real mj,wj; /* For mass-weighting of the positions */
2642 erg=rotg->enfrotgrp;
2644 bPF = rotg->eType==erotgRM2PF;
2645 clear_rvec(innersumvec);
2648 /* For the pivot-free variant we have to use the current center of
2649 * mass of the rotation group instead of the pivot u */
2650 get_center(erg->xc, erg->mc, rotg->nat, erg->xc_center);
2652 /* Also, we precalculate the second term of the forces that is identical
2653 * (up to the weight factor mj) for all forces */
2654 radial_motion2_precalc_inner_sum(rotg,innersumvec);
2657 N_M = rotg->nat * erg->invmass;
2659 /* Each process calculates the forces on its local atoms */
2660 for (j=0; j<erg->nat_loc; j++)
2664 /* Local index of a rotation group atom */
2665 ii = erg->ind_loc[j];
2666 /* Position of this atom in the collective array */
2667 iigrp = erg->xc_ref_ind[j];
2668 /* Mass-weighting */
2669 mj = erg->mc[iigrp];
2671 /* Current position of this atom: x[ii] */
2672 copy_rvec(x[ii], xj);
2674 /* Shift this atom such that it is near its reference */
2675 shift_single_coord(box, xj, erg->xc_shifts[iigrp]);
2677 /* The (unrotated) reference position is yj0. yc0 has already
2678 * been subtracted in init_rot_group */
2679 copy_rvec(rotg->x_ref[iigrp], tmpvec); /* tmpvec = yj0 - yc0 */
2681 /* Calculate Omega.(yj0-yc0) */
2682 mvmul(erg->rotmat, tmpvec, rj); /* rj = Omega.(yj0-yc0) */
2687 copy_rvec(erg->x_loc_pbc[j], xj);
2688 copy_rvec(erg->xr_loc[j], rj); /* rj = Omega.(yj0-u) */
2690 /* Mass-weighting */
2693 /* Calculate (xj-u) resp. (xj-xc) */
2694 rvec_sub(xj, erg->xc_center, xj_u); /* xj_u = xj-u */
2696 cprod(rotg->vec, xj_u, v_xj_u); /* v_xj_u = v x (xj-u) */
2698 fac = norm2(v_xj_u);
2700 psijstar = 1.0/(fac + rotg->eps); /* psistar = -------------------- */
2701 /* |v x (xj-u)|^2 + eps */
2703 psij = gmx_invsqrt(fac); /* 1 */
2704 /* psij = ------------ */
2707 svmul(psij, v_xj_u, sj); /* sj = psij * (v x (xj-u) ) */
2709 fac = iprod(v_xj_u, rj); /* fac = (v x (xj-u)).rj */
2712 sjrj = iprod(sj, rj); /* sjrj = sj.rj */
2714 svmul(psijstar/psij, rj, tmpvec);
2715 svmul(psijstar*psijstar/(psij*psij*psij) * sjrj, sj, tmpvec2);
2716 rvec_dec(tmpvec, tmpvec2);
2717 cprod(tmpvec, rotg->vec, tmpvec2);
2719 /* Store the additional force so that it can be added to the force
2720 * array after the normal forces have been evaluated */
2721 svmul(-rotg->k*wj*sjrj, tmpvec2, tmpvec);
2722 svmul(mj, innersumvec, tmpvec2); /* This is != 0 only for the pivot-free variant */
2724 rvec_add(tmpvec2, tmpvec, erg->f_rot_loc[j]);
2725 Vpart += wj*psijstar*fac2;
2728 /* Add to the torque of this rotation group */
2729 erg->torque_v += torque(rotg->vec, erg->f_rot_loc[j], xj, erg->xc_center);
2731 /* Calculate the angle between reference and actual rotation group atom. */
2732 angle(rotg, xj_u, rj, &alpha, &weight); /* angle in rad, weighted */
2733 erg->angle_v += alpha * weight;
2734 erg->weight_v += weight;
2739 } /* end of loop over local rotation group atoms */
2740 erg->V = 0.5*rotg->k*Vpart;
2746 /* Determine the smallest and largest position vector (with respect to the
2747 * rotation vector) for the reference group */
2748 static void get_firstlast_atom_ref(
2753 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
2755 real xcproj; /* The projection of a reference position on the
2757 real minproj, maxproj; /* Smallest and largest projection on v */
2761 erg=rotg->enfrotgrp;
2763 /* Start with some value */
2764 minproj = iprod(rotg->x_ref[0], rotg->vec);
2767 /* This is just to ensure that it still works if all the atoms of the
2768 * reference structure are situated in a plane perpendicular to the rotation
2771 *lastindex = rotg->nat-1;
2773 /* Loop over all atoms of the reference group,
2774 * project them on the rotation vector to find the extremes */
2775 for (i=0; i<rotg->nat; i++)
2777 xcproj = iprod(rotg->x_ref[i], rotg->vec);
2778 if (xcproj < minproj)
2783 if (xcproj > maxproj)
2792 /* Allocate memory for the slabs */
2793 static void allocate_slabs(
2800 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
2804 erg=rotg->enfrotgrp;
2806 /* More slabs than are defined for the reference are never needed */
2807 nslabs = erg->slab_last_ref - erg->slab_first_ref + 1;
2809 /* Remember how many we allocated */
2810 erg->nslabs_alloc = nslabs;
2812 if (MASTER(cr) && bVerbose)
2813 fprintf(fplog, "%s allocating memory to store data for %d slabs (rotation group %d).\n",
2815 snew(erg->slab_center , nslabs);
2816 snew(erg->slab_center_ref , nslabs);
2817 snew(erg->slab_weights , nslabs);
2818 snew(erg->slab_torque_v , nslabs);
2819 snew(erg->slab_data , nslabs);
2820 snew(erg->gn_atom , nslabs);
2821 snew(erg->gn_slabind , nslabs);
2822 snew(erg->slab_innersumvec, nslabs);
2823 for (i=0; i<nslabs; i++)
2825 snew(erg->slab_data[i].x , rotg->nat);
2826 snew(erg->slab_data[i].ref , rotg->nat);
2827 snew(erg->slab_data[i].weight, rotg->nat);
2829 snew(erg->xc_ref_sorted, rotg->nat);
2830 snew(erg->xc_sortind , rotg->nat);
2831 snew(erg->firstatom , nslabs);
2832 snew(erg->lastatom , nslabs);
2836 /* From the extreme coordinates of the reference group, determine the first
2837 * and last slab of the reference. We can never have more slabs in the real
2838 * simulation than calculated here for the reference.
2840 static void get_firstlast_slab_ref(t_rotgrp *rotg, real mc[], int ref_firstindex, int ref_lastindex)
2842 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
2843 int first,last,firststart;
2847 erg=rotg->enfrotgrp;
2848 first = get_first_slab(rotg, erg->max_beta, rotg->x_ref[ref_firstindex]);
2849 last = get_last_slab( rotg, erg->max_beta, rotg->x_ref[ref_lastindex ]);
2852 while (get_slab_weight(first, rotg, rotg->x_ref, mc, &dummy) > WEIGHT_MIN)
2856 erg->slab_first_ref = first+1;
2857 while (get_slab_weight(last, rotg, rotg->x_ref, mc, &dummy) > WEIGHT_MIN)
2861 erg->slab_last_ref = last-1;
2863 erg->slab_buffer = firststart - erg->slab_first_ref;
2868 static void init_rot_group(FILE *fplog,t_commrec *cr,int g,t_rotgrp *rotg,
2869 rvec *x,gmx_mtop_t *mtop,gmx_bool bVerbose,FILE *out_slabs)
2873 gmx_bool bFlex,bColl;
2875 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
2876 int ref_firstindex, ref_lastindex;
2877 real mass,totalmass;
2880 /* Do we have a flexible axis? */
2881 bFlex = ISFLEX(rotg);
2882 /* Do we use a global set of coordinates? */
2883 bColl = ISCOLL(rotg);
2885 erg=rotg->enfrotgrp;
2887 /* Allocate space for collective coordinates if needed */
2890 snew(erg->xc , rotg->nat);
2891 snew(erg->xc_shifts , rotg->nat);
2892 snew(erg->xc_eshifts, rotg->nat);
2894 /* Save the original (whole) set of positions such that later the
2895 * molecule can always be made whole again */
2896 snew(erg->xc_old , rotg->nat);
2899 for (i=0; i<rotg->nat; i++)
2902 copy_rvec(x[ii], erg->xc_old[i]);
2907 gmx_bcast(rotg->nat*sizeof(erg->xc_old[0]),erg->xc_old, cr);
2910 if (rotg->eFittype == erotgFitNORM)
2912 snew(erg->xc_ref_length, rotg->nat); /* in case fit type NORM is chosen */
2913 snew(erg->xc_norm , rotg->nat);
2918 snew(erg->xr_loc , rotg->nat);
2919 snew(erg->x_loc_pbc, rotg->nat);
2922 snew(erg->f_rot_loc , rotg->nat);
2923 snew(erg->xc_ref_ind, rotg->nat);
2925 /* xc_ref_ind needs to be set to identity in the serial case */
2927 for (i=0; i<rotg->nat; i++)
2928 erg->xc_ref_ind[i] = i;
2930 /* Copy the masses so that the COM can be determined. For all types of
2931 * enforced rotation, we store the masses in the erg->mc array. */
2932 snew(erg->mc, rotg->nat);
2934 snew(erg->mc_sorted, rotg->nat);
2936 snew(erg->m_loc, rotg->nat);
2938 for (i=0; i<rotg->nat; i++)
2942 gmx_mtop_atomnr_to_atom(mtop,rotg->ind[i],&atom);
2952 erg->invmass = 1.0/totalmass;
2954 /* Set xc_ref_center for any rotation potential */
2955 if ((rotg->eType==erotgISO) || (rotg->eType==erotgPM) || (rotg->eType==erotgRM) || (rotg->eType==erotgRM2))
2957 /* Set the pivot point for the fixed, stationary-axis potentials. This
2958 * won't change during the simulation */
2959 copy_rvec(rotg->pivot, erg->xc_ref_center);
2960 copy_rvec(rotg->pivot, erg->xc_center );
2964 /* Center of the reference positions */
2965 get_center(rotg->x_ref, erg->mc, rotg->nat, erg->xc_ref_center);
2967 /* Center of the actual positions */
2970 snew(xdum, rotg->nat);
2971 for (i=0; i<rotg->nat; i++)
2974 copy_rvec(x[ii], xdum[i]);
2976 get_center(xdum, erg->mc, rotg->nat, erg->xc_center);
2981 gmx_bcast(sizeof(erg->xc_center), erg->xc_center, cr);
2985 if ( (rotg->eType != erotgFLEX) && (rotg->eType != erotgFLEX2) )
2987 /* Put the reference positions into origin: */
2988 for (i=0; i<rotg->nat; i++)
2989 rvec_dec(rotg->x_ref[i], erg->xc_ref_center);
2992 /* Enforced rotation with flexible axis */
2995 /* Calculate maximum beta value from minimum gaussian (performance opt.) */
2996 erg->max_beta = calc_beta_max(rotg->min_gaussian, rotg->slab_dist);
2998 /* Determine the smallest and largest coordinate with respect to the rotation vector */
2999 get_firstlast_atom_ref(rotg, &ref_firstindex, &ref_lastindex);
3001 /* From the extreme coordinates of the reference group, determine the first
3002 * and last slab of the reference. */
3003 get_firstlast_slab_ref(rotg, erg->mc, ref_firstindex, ref_lastindex);
3005 /* Allocate memory for the slabs */
3006 allocate_slabs(rotg, fplog, g, bVerbose, cr);
3008 /* Flexible rotation: determine the reference centers for the rest of the simulation */
3009 erg->slab_first = erg->slab_first_ref;
3010 erg->slab_last = erg->slab_last_ref;
3011 get_slab_centers(rotg,rotg->x_ref,erg->mc,cr,g,-1,out_slabs,TRUE,TRUE);
3013 /* Length of each x_rotref vector from center (needed if fit routine NORM is chosen): */
3014 if (rotg->eFittype == erotgFitNORM)
3016 for (i=0; i<rotg->nat; i++)
3018 rvec_sub(rotg->x_ref[i], erg->xc_ref_center, coord);
3019 erg->xc_ref_length[i] = norm(coord);
3026 extern void dd_make_local_rotation_groups(gmx_domdec_t *dd,t_rot *rot)
3031 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
3035 for(g=0; g<rot->ngrp; g++)
3037 rotg = &rot->grp[g];
3038 erg = rotg->enfrotgrp;
3041 dd_make_local_group_indices(ga2la,rotg->nat,rotg->ind,
3042 &erg->nat_loc,&erg->ind_loc,&erg->nalloc_loc,erg->xc_ref_ind);
3047 extern void init_rot(FILE *fplog,t_inputrec *ir,int nfile,const t_filenm fnm[],
3048 t_commrec *cr, rvec *x, matrix box, gmx_mtop_t *mtop, const output_env_t oenv,
3049 gmx_bool bVerbose, unsigned long Flags)
3054 int nat_max=0; /* Size of biggest rotation group */
3055 gmx_enfrot_t er; /* Pointer to the enforced rotation buffer variables */
3056 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
3057 rvec *x_pbc=NULL; /* Space for the pbc-correct atom positions */
3058 gmx_bool bHaveFlexGroups = FALSE;
3061 if ( (PAR(cr)) && !DOMAINDECOMP(cr) )
3062 gmx_fatal(FARGS, "Enforced rotation is only implemented for domain decomposition!");
3064 if ( MASTER(cr) && bVerbose)
3065 fprintf(stdout, "%s Initializing ...\n", RotStr);
3069 snew(rot->enfrot, 1);
3072 /* Output every step for reruns */
3073 if (Flags & MD_RERUN)
3076 fprintf(fplog, "%s rerun - will write rotation output every available step.\n", RotStr);
3081 er->out_slabs = NULL;
3083 er->out_slabs = open_slab_out(rot, opt2fn("-rs",nfile,fnm));
3087 /* Remove pbc, make molecule whole.
3088 * When ir->bContinuation=TRUE this has already been done, but ok. */
3089 snew(x_pbc,mtop->natoms);
3090 m_rveccopy(mtop->natoms,x,x_pbc);
3091 do_pbc_first_mtop(NULL,ir->ePBC,box,mtop,x_pbc);
3094 for (g=0; g<rot->ngrp; g++)
3096 rotg = &rot->grp[g];
3098 bHaveFlexGroups = TRUE;
3101 fprintf(fplog,"%s group %d type '%s'\n", RotStr, g, erotg_names[rotg->eType]);
3105 /* Allocate space for the rotation group's data: */
3106 snew(rotg->enfrotgrp, 1);
3107 erg = rotg->enfrotgrp;
3109 nat_max=max(nat_max, rotg->nat);
3114 erg->nalloc_loc = 0;
3115 erg->ind_loc = NULL;
3119 erg->nat_loc = rotg->nat;
3120 erg->ind_loc = rotg->ind;
3122 init_rot_group(fplog,cr,g,rotg,x_pbc,mtop,bVerbose,er->out_slabs);
3126 /* Allocate space for enforced rotation buffer variables */
3127 er->bufsize = nat_max;
3128 snew(er->data, nat_max);
3129 snew(er->xbuf, nat_max);
3130 snew(er->mbuf, nat_max);
3132 /* Buffers for MPI reducing torques, angles, weights (for each group), and V */
3133 er->mpi_bufsize = 4*rot->ngrp; /* To start with */
3134 snew(er->mpi_inbuf , er->mpi_bufsize);
3135 snew(er->mpi_outbuf, er->mpi_bufsize);
3137 /* Only do I/O on the MASTER */
3138 er->out_angles = NULL;
3140 er->out_torque = NULL;
3143 er->out_rot = open_rot_out(opt2fn("-ro",nfile,fnm), rot, oenv, Flags);
3144 if (bHaveFlexGroups)
3146 if (rot->nstrout > 0)
3147 er->out_angles = open_angles_out(rot, opt2fn("-ra",nfile,fnm));
3148 if (rot->nstsout > 0)
3149 er->out_torque = open_torque_out(rot, opt2fn("-rt",nfile,fnm));
3156 extern void finish_rot(FILE *fplog,t_rot *rot)
3158 gmx_enfrot_t er; /* Pointer to the enforced rotation buffer variables */
3163 gmx_fio_fclose(er->out_rot);
3165 gmx_fio_fclose(er->out_slabs);
3167 gmx_fio_fclose(er->out_angles);
3169 gmx_fio_fclose(er->out_torque);
3173 /* Rotate the local reference positions and store them in
3174 * erg->xr_loc[0...(nat_loc-1)]
3176 * Note that we already subtracted u or y_c from the reference positions
3177 * in init_rot_group().
3179 static void rotate_local_reference(t_rotgrp *rotg)
3181 gmx_enfrotgrp_t erg;
3185 erg=rotg->enfrotgrp;
3187 for (i=0; i<erg->nat_loc; i++)
3189 /* Index of this rotation group atom with respect to the whole rotation group */
3190 ii = erg->xc_ref_ind[i];
3192 mvmul(erg->rotmat, rotg->x_ref[ii], erg->xr_loc[i]);
3197 /* Select the PBC representation for each local x position and store that
3198 * for later usage. We assume the right PBC image of an x is the one nearest to
3199 * its rotated reference */
3200 static void choose_pbc_image(rvec x[], t_rotgrp *rotg, matrix box, int npbcdim)
3203 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
3208 erg=rotg->enfrotgrp;
3210 for (i=0; i<erg->nat_loc; i++)
3214 /* Index of a rotation group atom */
3215 ii = erg->ind_loc[i];
3217 /* Get the reference position. The pivot (or COM or COG) was already
3218 * subtracted in init_rot_group() from the reference positions. Also,
3219 * the reference positions have already been rotated in
3220 * rotate_local_reference() */
3221 copy_rvec(erg->xr_loc[i], xref);
3223 /* Subtract the (old) center from the current positions
3224 * (just to determine the shifts!) */
3225 rvec_sub(x[ii], erg->xc_center, xcurr);
3227 /* Shortest PBC distance between the atom and its reference */
3228 rvec_sub(xcurr, xref, dx);
3230 /* Determine the shift for this atom */
3231 for(m=npbcdim-1; m>=0; m--)
3233 while (dx[m] < -0.5*box[m][m])
3235 for(d=0; d<DIM; d++)
3239 while (dx[m] >= 0.5*box[m][m])
3241 for(d=0; d<DIM; d++)
3247 /* Apply the shift to the current atom */
3248 copy_rvec(x[ii], erg->x_loc_pbc[i]);
3249 shift_single_coord(box, erg->x_loc_pbc[i], shift);
3254 extern void do_rotation(
3260 gmx_large_int_t step,
3261 gmx_wallcycle_t wcycle,
3267 gmx_bool outstep_slab, outstep_rot;
3268 gmx_bool bFlex,bColl;
3270 gmx_enfrot_t er; /* Pointer to the enforced rotation buffer variables */
3271 gmx_enfrotgrp_t erg; /* Pointer to enforced rotation group data */
3281 /* When to output in main rotation output file */
3282 outstep_rot = do_per_step(step, rot->nstrout);
3283 /* When to output per-slab data */
3284 outstep_slab = do_per_step(step, rot->nstsout);
3286 /* Output time into rotation output file */
3287 if (outstep_rot && MASTER(cr))
3288 fprintf(er->out_rot, "%12.3e",t);
3290 /**************************************************************************/
3291 /* First do ALL the communication! */
3292 for(g=0; g<rot->ngrp; g++)
3294 rotg = &rot->grp[g];
3295 erg=rotg->enfrotgrp;
3297 /* Do we have a flexible axis? */
3298 bFlex = ISFLEX(rotg);
3299 /* Do we use a collective (global) set of coordinates? */
3300 bColl = ISCOLL(rotg);
3302 /* Calculate the rotation matrix for this angle: */
3303 erg->degangle = rotg->rate * t;
3304 calc_rotmat(rotg->vec,erg->degangle,erg->rotmat);
3308 /* Transfer the rotation group's positions such that every node has
3309 * all of them. Every node contributes its local positions x and stores
3310 * it in the collective erg->xc array. */
3311 communicate_group_positions(cr,erg->xc, erg->xc_shifts, erg->xc_eshifts, bNS,
3312 x, rotg->nat, erg->nat_loc, erg->ind_loc, erg->xc_ref_ind, erg->xc_old, box);
3316 /* Fill the local masses array;
3317 * this array changes in DD/neighborsearching steps */
3320 for (i=0; i<erg->nat_loc; i++)
3322 /* Index of local atom w.r.t. the collective rotation group */
3323 ii = erg->xc_ref_ind[i];
3324 erg->m_loc[i] = erg->mc[ii];
3328 /* Calculate Omega*(y_i-y_c) for the local positions */
3329 rotate_local_reference(rotg);
3331 /* Choose the nearest PBC images of the group atoms with respect
3332 * to the rotated reference positions */
3333 choose_pbc_image(x, rotg, box, 3);
3335 /* Get the center of the rotation group */
3336 if ( (rotg->eType==erotgISOPF) || (rotg->eType==erotgPMPF) )
3337 get_center_comm(cr, erg->x_loc_pbc, erg->m_loc, erg->nat_loc, rotg->nat, erg->xc_center);
3340 } /* End of loop over rotation groups */
3342 /**************************************************************************/
3343 /* Done communicating, we can start to count cycles now ... */
3344 wallcycle_start(wcycle, ewcROT);
3345 GMX_MPE_LOG(ev_rotcycles_start);
3351 for(g=0; g<rot->ngrp; g++)
3353 rotg = &rot->grp[g];
3354 erg=rotg->enfrotgrp;
3356 bFlex = ISFLEX(rotg);
3357 bColl = ISCOLL(rotg);
3359 if (outstep_rot && MASTER(cr))
3360 fprintf(er->out_rot, "%12.4f", erg->degangle);
3362 /* Clear values from last time step */
3364 erg->torque_v = 0.0;
3366 erg->weight_v = 0.0;
3374 do_fixed(cr,rotg,x,box,t,step,outstep_rot);
3377 do_radial_motion(cr,rotg,x,box,t,step,outstep_rot);
3380 do_radial_motion_pf(cr,rotg,x,box,t,step,outstep_rot);
3384 do_radial_motion2(cr,rotg,x,box,t,step,outstep_rot);
3388 /* Subtract the center of the rotation group from the collective positions array
3389 * Also store the center in erg->xc_center since it needs to be subtracted
3390 * in the low level routines from the local coordinates as well */
3391 get_center(erg->xc, erg->mc, rotg->nat, erg->xc_center);
3392 svmul(-1.0, erg->xc_center, transvec);
3393 translate_x(erg->xc, rotg->nat, transvec);
3394 do_flexible(cr,er,rotg,g,x,box,t,step,outstep_rot,outstep_slab);
3398 /* Do NOT subtract the center of mass in the low level routines! */
3399 clear_rvec(erg->xc_center);
3400 do_flexible(cr,er,rotg,g,x,box,t,step,outstep_rot,outstep_slab);
3403 gmx_fatal(FARGS, "No such rotation potential.");
3410 fprintf(stderr, "%s calculation (step %d) took %g seconds.\n", RotStr, step, MPI_Wtime()-t0);
3413 /* Stop the cycle counter and add to the force cycles for load balancing */
3414 cycles_rot = wallcycle_stop(wcycle,ewcROT);
3415 if (DOMAINDECOMP(cr) && wcycle)
3416 dd_cycles_add(cr->dd,cycles_rot,ddCyclF);
3417 GMX_MPE_LOG(ev_rotcycles_finish);