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 if (MASTER(cr)) { \
188 fprintf(stderr,"potential = %15.8f\n" "torque = %15.8f", erg->V, erg->torque_v); \
191 #define PRINT_FORCE_J
192 #define PRINT_POT_TAU
195 /* Shortcuts for often used queries */
196 #define ISFLEX(rg) ( (rg->eType==erotgFLEX) || (rg->eType==erotgFLEXT) || (rg->eType==erotgFLEX2) || (rg->eType==erotgFLEX2T) )
197 #define ISCOLL(rg) ( (rg->eType==erotgFLEX) || (rg->eType==erotgFLEXT) || (rg->eType==erotgFLEX2) || (rg->eType==erotgFLEX2T) || (rg->eType==erotgRMPF) || (rg->eType==erotgRM2PF) )
201 static double** allocate_square_matrix(int dim
)
215 static void free_square_matrix(double** mat
, int dim
)
220 for (i
=0; i
<dim
; i
++)
226 /* Output rotation energy, torques, etc. for each rotation group */
227 static void reduce_output(t_commrec
*cr
, t_rot
*rot
, real t
, gmx_large_int_t step
)
229 int g
,i
,islab
,nslabs
=0;
230 int count
; /* MPI element counter */
232 gmx_enfrot_t er
; /* Pointer to the enforced rotation buffer variables */
233 gmx_enfrotgrp_t erg
; /* Pointer to enforced rotation group data */
239 /* Fill the MPI buffer with stuff to reduce: */
243 for (g
=0; g
< rot
->ngrp
; g
++)
246 erg
= rotg
->enfrotgrp
;
247 nslabs
= erg
->slab_last
- erg
->slab_first
+ 1;
248 er
->mpi_inbuf
[count
++] = erg
->V
;
249 er
->mpi_inbuf
[count
++] = erg
->torque_v
;
250 er
->mpi_inbuf
[count
++] = erg
->angle_v
;
251 er
->mpi_inbuf
[count
++] = erg
->weight_v
; /* weights are not needed for flex types, but this is just a single value */
254 /* (Re-)allocate memory for MPI buffer: */
255 if (er
->mpi_bufsize
< count
+nslabs
)
257 er
->mpi_bufsize
= count
+nslabs
;
258 srenew(er
->mpi_inbuf
, er
->mpi_bufsize
);
259 srenew(er
->mpi_outbuf
, er
->mpi_bufsize
);
261 for (i
=0; i
<nslabs
; i
++)
262 er
->mpi_inbuf
[count
++] = erg
->slab_torque_v
[i
];
266 MPI_Reduce(er
->mpi_inbuf
, er
->mpi_outbuf
, count
, GMX_MPI_REAL
, MPI_SUM
, MASTERRANK(cr
), cr
->mpi_comm_mygroup
);
268 /* Copy back the reduced data from the buffer on the master */
272 for (g
=0; g
< rot
->ngrp
; g
++)
275 erg
= rotg
->enfrotgrp
;
276 nslabs
= erg
->slab_last
- erg
->slab_first
+ 1;
277 erg
->V
= er
->mpi_outbuf
[count
++];
278 erg
->torque_v
= er
->mpi_outbuf
[count
++];
279 erg
->angle_v
= er
->mpi_outbuf
[count
++];
280 erg
->weight_v
= er
->mpi_outbuf
[count
++];
283 for (i
=0; i
<nslabs
; i
++)
284 erg
->slab_torque_v
[i
] = er
->mpi_outbuf
[count
++];
293 /* Angle and torque for each rotation group */
294 for (g
=0; g
< rot
->ngrp
; g
++)
297 bFlex
= ISFLEX(rotg
);
301 /* Output to main rotation output file: */
302 if ( do_per_step(step
, rot
->nstrout
) )
305 fprintf(er
->out_rot
, "%12.4f", erg
->angle_v
); /* RMSD fit angle */
307 fprintf(er
->out_rot
, "%12.4f", (erg
->angle_v
/erg
->weight_v
)*180.0*M_1_PI
);
308 fprintf(er
->out_rot
, "%12.3e", erg
->torque_v
);
309 fprintf(er
->out_rot
, "%12.3e", erg
->V
);
312 /* Output to torque log file: */
313 if ( bFlex
&& do_per_step(step
, rot
->nstsout
) )
315 fprintf(er
->out_torque
, "%12.3e%6d", t
, g
);
316 for (i
=erg
->slab_first
; i
<=erg
->slab_last
; i
++)
318 islab
= i
- erg
->slab_first
; /* slab index */
319 /* Only output if enough weight is in slab */
320 if (erg
->slab_weights
[islab
] > rotg
->min_gaussian
)
321 fprintf(er
->out_torque
, "%6d%12.3e", i
, erg
->slab_torque_v
[islab
]);
323 fprintf(er
->out_torque
, "\n");
326 if ( do_per_step(step
, rot
->nstrout
) )
327 fprintf(er
->out_rot
, "\n");
332 /* Add the forces from enforced rotation potential to the local forces.
333 * Should be called after the SR forces have been evaluated */
334 extern real
add_rot_forces(t_rot
*rot
, rvec f
[], t_commrec
*cr
, gmx_large_int_t step
, real t
)
338 gmx_enfrot_t er
; /* Pointer to the enforced rotation buffer variables */
339 gmx_enfrotgrp_t erg
; /* Pointer to enforced rotation group data */
344 GMX_MPE_LOG(ev_add_rot_forces_start
);
346 /* Reduce energy,torque, angles etc. to get the sum values (per rotation group)
347 * on the master and output these values to file. */
348 if ( do_per_step(step
, rot
->nstrout
) || do_per_step(step
, rot
->nstsout
) )
349 reduce_output(cr
, rot
, t
, step
);
351 /* Total rotation potential is the sum over all rotation groups */
354 /* Loop over enforced rotation groups (usually 1, though)
355 * Apply the forces from rotation potentials */
356 for (g
=0; g
<rot
->ngrp
; g
++)
361 for (l
=0; l
<erg
->nat_loc
; l
++)
363 /* Get the right index of the local force */
364 ii
= erg
->ind_loc
[l
];
366 rvec_inc(f
[ii
],erg
->f_rot_loc
[l
]);
372 GMX_MPE_LOG(ev_add_rot_forces_finish
);
374 return (MASTER(cr
)? er
->Vrot
: 0.0);
378 /* The Gaussian norm is chosen such that the sum of the gaussian functions
379 * over the slabs is approximately 1.0 everywhere */
380 #define GAUSS_NORM 0.569917543430618
383 /* Calculate the maximum beta that leads to a gaussian larger min_gaussian,
384 * also does some checks
386 static double calc_beta_max(real min_gaussian
, real slab_dist
)
392 /* Actually the next two checks are already made in grompp */
394 gmx_fatal(FARGS
, "Slab distance of flexible rotation groups must be >=0 !");
395 if (min_gaussian
<= 0)
396 gmx_fatal(FARGS
, "Cutoff value for Gaussian must be > 0. (You requested %f)");
398 /* Define the sigma value */
399 sigma
= 0.7*slab_dist
;
401 /* Calculate the argument for the logarithm and check that the log() result is negative or 0 */
402 arg
= min_gaussian
/GAUSS_NORM
;
404 gmx_fatal(FARGS
, "min_gaussian of flexible rotation groups must be <%g", GAUSS_NORM
);
406 return sqrt(-2.0*sigma
*sigma
*log(min_gaussian
/GAUSS_NORM
));
410 static inline real
calc_beta(rvec curr_x
, t_rotgrp
*rotg
, int n
)
412 return iprod(curr_x
, rotg
->vec
) - rotg
->slab_dist
* n
;
416 static inline real
gaussian_weight(rvec curr_x
, t_rotgrp
*rotg
, int n
)
418 const real norm
= GAUSS_NORM
;
422 /* Define the sigma value */
423 sigma
= 0.7*rotg
->slab_dist
;
424 /* Calculate the Gaussian value of slab n for position curr_x */
425 return norm
* exp( -0.5 * sqr( calc_beta(curr_x
, rotg
, n
)/sigma
) );
429 /* Returns the weight in a single slab, also calculates the Gaussian- and mass-
430 * weighted sum of positions for that slab */
431 static real
get_slab_weight(int j
, t_rotgrp
*rotg
, rvec xc
[], real mc
[], rvec
*x_weighted_sum
)
433 rvec curr_x
; /* The position of an atom */
434 rvec curr_x_weighted
; /* The gaussian-weighted position */
435 real gaussian
; /* A single gaussian weight */
436 real wgauss
; /* gaussian times current mass */
437 real slabweight
= 0.0; /* The sum of weights in the slab */
439 gmx_enfrotgrp_t erg
; /* Pointer to enforced rotation group data */
443 clear_rvec(*x_weighted_sum
);
446 islab
= j
- erg
->slab_first
;
448 /* Loop over all atoms in the rotation group */
449 for (i
=0; i
<rotg
->nat
; i
++)
451 copy_rvec(xc
[i
], curr_x
);
452 gaussian
= gaussian_weight(curr_x
, rotg
, j
);
453 wgauss
= gaussian
* mc
[i
];
454 svmul(wgauss
, curr_x
, curr_x_weighted
);
455 rvec_add(*x_weighted_sum
, curr_x_weighted
, *x_weighted_sum
);
456 slabweight
+= wgauss
;
457 } /* END of loop over rotation group atoms */
463 static void get_slab_centers(
464 t_rotgrp
*rotg
, /* The rotation group information */
465 rvec
*xc
, /* The rotation group positions; will
466 typically be enfrotgrp->xc, but at first call
467 it is enfrotgrp->xc_ref */
468 real
*mc
, /* The masses of the rotation group atoms */
469 t_commrec
*cr
, /* Communication record */
470 int g
, /* The number of the rotation group */
471 real time
, /* Used for output only */
472 FILE *out_slabs
, /* For outputting center per slab information */
473 gmx_bool bOutStep
, /* Is this an output step? */
474 gmx_bool bReference
) /* If this routine is called from
475 init_rot_group we need to store
476 the reference slab centers */
479 gmx_enfrotgrp_t erg
; /* Pointer to enforced rotation group data */
484 /* Loop over slabs */
485 for (j
= erg
->slab_first
; j
<= erg
->slab_last
; j
++)
487 islab
= j
- erg
->slab_first
;
488 erg
->slab_weights
[islab
] = get_slab_weight(j
, rotg
, xc
, mc
, &erg
->slab_center
[islab
]);
490 /* We can do the calculations ONLY if there is weight in the slab! */
491 if (erg
->slab_weights
[islab
] > WEIGHT_MIN
)
493 svmul(1.0/erg
->slab_weights
[islab
], erg
->slab_center
[islab
], erg
->slab_center
[islab
]);
497 /* We need to check this here, since we divide through slab_weights
498 * in the flexible low-level routines! */
499 gmx_fatal(FARGS
, "Not enough weight in slab %d. Slab center cannot be determined!", j
);
502 /* At first time step: save the centers of the reference structure */
504 copy_rvec(erg
->slab_center
[islab
], erg
->slab_center_ref
[islab
]);
505 } /* END of loop over slabs */
507 /* Output on the master */
508 if (MASTER(cr
) && bOutStep
)
510 fprintf(out_slabs
, "%12.3e%6d", time
, g
);
511 for (j
= erg
->slab_first
; j
<= erg
->slab_last
; j
++)
513 islab
= j
- erg
->slab_first
;
514 fprintf(out_slabs
, "%6d%12.3e%12.3e%12.3e",
515 j
,erg
->slab_center
[islab
][XX
],erg
->slab_center
[islab
][YY
],erg
->slab_center
[islab
][ZZ
]);
517 fprintf(out_slabs
, "\n");
522 static void calc_rotmat(
524 real degangle
, /* Angle alpha of rotation at time t in degrees */
525 matrix rotmat
) /* Rotation matrix */
527 real radangle
; /* Rotation angle in radians */
528 real cosa
; /* cosine alpha */
529 real sina
; /* sine alpha */
530 real OMcosa
; /* 1 - cos(alpha) */
531 real dumxy
, dumxz
, dumyz
; /* save computations */
532 rvec rot_vec
; /* Rotate around rot_vec ... */
535 radangle
= degangle
* M_PI
/180.0;
536 copy_rvec(vec
, rot_vec
);
538 /* Precompute some variables: */
539 cosa
= cos(radangle
);
540 sina
= sin(radangle
);
542 dumxy
= rot_vec
[XX
]*rot_vec
[YY
]*OMcosa
;
543 dumxz
= rot_vec
[XX
]*rot_vec
[ZZ
]*OMcosa
;
544 dumyz
= rot_vec
[YY
]*rot_vec
[ZZ
]*OMcosa
;
546 /* Construct the rotation matrix for this rotation group: */
548 rotmat
[XX
][XX
] = cosa
+ rot_vec
[XX
]*rot_vec
[XX
]*OMcosa
;
549 rotmat
[YY
][XX
] = dumxy
+ rot_vec
[ZZ
]*sina
;
550 rotmat
[ZZ
][XX
] = dumxz
- rot_vec
[YY
]*sina
;
552 rotmat
[XX
][YY
] = dumxy
- rot_vec
[ZZ
]*sina
;
553 rotmat
[YY
][YY
] = cosa
+ rot_vec
[YY
]*rot_vec
[YY
]*OMcosa
;
554 rotmat
[ZZ
][YY
] = dumyz
+ rot_vec
[XX
]*sina
;
556 rotmat
[XX
][ZZ
] = dumxz
+ rot_vec
[YY
]*sina
;
557 rotmat
[YY
][ZZ
] = dumyz
- rot_vec
[XX
]*sina
;
558 rotmat
[ZZ
][ZZ
] = cosa
+ rot_vec
[ZZ
]*rot_vec
[ZZ
]*OMcosa
;
563 for (iii
=0; iii
<3; iii
++) {
564 for (jjj
=0; jjj
<3; jjj
++)
565 fprintf(stderr
, " %10.8f ", rotmat
[iii
][jjj
]);
566 fprintf(stderr
, "\n");
572 /* Calculates torque on the rotation axis tau = position x force */
573 static inline real
torque(
574 rvec rotvec
, /* rotation vector; MUST be normalized! */
575 rvec force
, /* force */
576 rvec x
, /* position of atom on which the force acts */
577 rvec pivot
) /* pivot point of rotation axis */
582 /* Subtract offset */
583 rvec_sub(x
,pivot
,vectmp
);
585 /* position x force */
586 cprod(vectmp
, force
, tau
);
588 /* Return the part of the torque which is parallel to the rotation vector */
589 return iprod(tau
, rotvec
);
593 /* Right-aligned output of value with standard width */
594 static void print_aligned(FILE *fp
, char *str
)
596 fprintf(fp
, "%12s", str
);
600 /* Right-aligned output of value with standard short width */
601 static void print_aligned_short(FILE *fp
, char *str
)
603 fprintf(fp
, "%6s", str
);
607 static FILE *open_output_file(const char *fn
, int steps
, const char what
[])
612 fp
= ffopen(fn
, "w");
614 fprintf(fp
, "# Output of %s is written in intervals of %d time step%s.\n#\n",
615 what
,steps
, steps
>1 ? "s":"");
621 /* Open output file for slab center data. Call on master only */
622 static FILE *open_slab_out(t_rot
*rot
, const char *fn
)
629 for (g
=0; g
<rot
->ngrp
; g
++)
635 fp
= open_output_file(fn
, rot
->nstsout
, "gaussian weighted slab centers");
636 fprintf(fp
, "# Rotation group %d (%s), slab distance %f nm, %s.\n",
637 g
, erotg_names
[rotg
->eType
], rotg
->slab_dist
, rotg
->bMassW
? "centers of mass":"geometrical centers");
643 fprintf(fp
, "# Reference centers are listed first (t=-1).\n");
644 fprintf(fp
, "# The following columns have the syntax:\n");
646 print_aligned_short(fp
, "t");
647 print_aligned_short(fp
, "grp");
648 /* Print ascii legend for the first two entries only ... */
651 print_aligned_short(fp
, "slab");
652 print_aligned(fp
, "X center");
653 print_aligned(fp
, "Y center");
654 print_aligned(fp
, "Z center");
656 fprintf(fp
, " ...\n");
664 /* Open output file and print some general information about the rotation groups.
665 * Call on master only */
666 static FILE *open_rot_out(const char *fn
, t_rot
*rot
, const output_env_t oenv
,
672 const char **setname
;
673 char buf
[50], buf2
[75];
674 gmx_enfrotgrp_t erg
; /* Pointer to enforced rotation group data */
678 if (Flags
& MD_APPENDFILES
)
680 fp
= gmx_fio_fopen(fn
,"a");
684 fp
= xvgropen(fn
, "Rotation angles and energy", "Time [ps]", "angles [degrees] and energies [kJ/mol]", oenv
);
685 fprintf(fp
, "# Output of enforced rotation data is written in intervals of %d time step%s.\n#\n", rot
->nstrout
, rot
->nstrout
> 1 ? "s":"");
686 fprintf(fp
, "# The scalar tau is the torque [kJ/mol] in the direction of the rotation vector v.\n");
687 fprintf(fp
, "# To obtain the vectorial torque, multiply tau with the group's rot_vec.\n");
688 fprintf(fp
, "# For flexible groups, tau(t,n) from all slabs n have been summed in a single value tau(t) here.\n");
689 fprintf(fp
, "# The torques tau(t,n) are found in the rottorque.log (-rt) output file\n");
692 for (g
=0; g
<rot
->ngrp
; g
++)
696 bFlex
= ISFLEX(rotg
);
698 fprintf(fp
, "# Rotation group %d, potential type '%s':\n" , g
, erotg_names
[rotg
->eType
]);
699 fprintf(fp
, "# rot_massw%d %s\n" , g
, yesno_names
[rotg
->bMassW
]);
700 fprintf(fp
, "# rot_vec%d %12.5e %12.5e %12.5e\n" , g
, rotg
->vec
[XX
], rotg
->vec
[YY
], rotg
->vec
[ZZ
]);
701 fprintf(fp
, "# rot_rate%d %12.5e degrees/ps\n" , g
, rotg
->rate
);
702 fprintf(fp
, "# rot_k%d %12.5e kJ/(mol*nm^2)\n" , g
, rotg
->k
);
703 if ( rotg
->eType
==erotgISO
|| rotg
->eType
==erotgPM
|| rotg
->eType
==erotgRM
|| rotg
->eType
==erotgRM2
)
704 fprintf(fp
, "# rot_pivot%d %12.5e %12.5e %12.5e nm\n", g
, rotg
->pivot
[XX
], rotg
->pivot
[YY
], rotg
->pivot
[ZZ
]);
708 fprintf(fp
, "# rot_slab_distance%d %f nm\n", g
, rotg
->slab_dist
);
709 fprintf(fp
, "# rot_min_gaussian%d %12.5e\n", g
, rotg
->min_gaussian
);
712 /* Output the centers of the rotation groups for the pivot-free potentials */
713 if ((rotg
->eType
==erotgISOPF
) || (rotg
->eType
==erotgPMPF
) || (rotg
->eType
==erotgRMPF
) || (rotg
->eType
==erotgRM2PF
714 || (rotg
->eType
==erotgFLEXT
) || (rotg
->eType
==erotgFLEX2T
)) )
716 fprintf(fp
, "# %s of ref. grp. %d %12.5e %12.5e %12.5e\n",
717 rotg
->bMassW
? "COM":"COG", g
,
718 erg
->xc_ref_center
[XX
], erg
->xc_ref_center
[YY
], erg
->xc_ref_center
[ZZ
]);
720 fprintf(fp
, "# initial %s grp. %d %12.5e %12.5e %12.5e\n",
721 rotg
->bMassW
? "COM":"COG", g
,
722 erg
->xc_center
[XX
], erg
->xc_center
[YY
], erg
->xc_center
[ZZ
]);
725 if ( (rotg
->eType
== erotgRM2
) || (rotg
->eType
==erotgFLEX2
) || (rotg
->eType
==erotgFLEX2T
) )
727 fprintf(fp
, "# rot_eps%d %12.5e nm^2\n", g
, rotg
->eps
);
731 fprintf(fp
, "#\n# Legend for the following data columns:\n");
733 print_aligned_short(fp
, "t");
735 snew(setname
, 4*rot
->ngrp
);
737 for (g
=0; g
<rot
->ngrp
; g
++)
740 sprintf(buf
, "theta_ref%d", g
);
741 print_aligned(fp
, buf
);
742 sprintf(buf2
, "%s [degrees]", buf
);
743 setname
[nsets
] = strdup(buf2
);
746 for (g
=0; g
<rot
->ngrp
; g
++)
749 bFlex
= ISFLEX(rotg
);
751 /* For flexible axis rotation we use RMSD fitting to determine the
752 * actual angle of the rotation group */
754 sprintf(buf
, "theta-fit%d", g
);
756 sprintf(buf
, "theta-av%d", g
);
757 print_aligned(fp
, buf
);
758 sprintf(buf2
, "%s [degrees]", buf
);
759 setname
[nsets
] = strdup(buf2
);
762 sprintf(buf
, "tau%d", g
);
763 print_aligned(fp
, buf
);
764 sprintf(buf2
, "%s [kJ/mol]", buf
);
765 setname
[nsets
] = strdup(buf2
);
768 sprintf(buf
, "energy%d", g
);
769 print_aligned(fp
, buf
);
770 sprintf(buf2
, "%s [kJ/mol]", buf
);
771 setname
[nsets
] = strdup(buf2
);
774 fprintf(fp
, "\n#\n");
777 xvgr_legend(fp
, nsets
, setname
, oenv
);
787 /* Call on master only */
788 static FILE *open_angles_out(t_rot
*rot
, const char *fn
)
795 /* Open output file and write some information about it's structure: */
796 fp
= open_output_file(fn
, rot
->nstsout
, "rotation group angles");
797 fprintf(fp
, "# All angles given in degrees, time in ps.\n");
798 for (g
=0; g
<rot
->ngrp
; g
++)
803 fprintf(fp
, "# Rotation group %d (%s), slab distance %f nm, fit type %s.\n",
804 g
, erotg_names
[rotg
->eType
], rotg
->slab_dist
, erotg_fitnames
[rotg
->eFittype
]);
807 fprintf(fp
, "# The following columns will have the syntax:\n");
809 print_aligned_short(fp
, "t");
810 print_aligned_short(fp
, "grp");
811 print_aligned(fp
, "theta_ref");
812 print_aligned_short(fp
, "slab");
813 print_aligned_short(fp
, "atoms");
814 print_aligned(fp
, "theta_fit");
815 print_aligned_short(fp
, "slab");
816 print_aligned_short(fp
, "atoms");
817 print_aligned(fp
, "theta_fit");
818 fprintf(fp
, " ...\n");
824 /* Open torque output file and write some information about it's structure.
825 * Call on master only */
826 static FILE *open_torque_out(t_rot
*rot
, const char *fn
)
833 fp
= open_output_file(fn
, rot
->nstsout
,"torques");
835 for (g
=0; g
<rot
->ngrp
; g
++)
840 fprintf(fp
, "# Rotation group %d (%s), slab distance %f nm.\n", g
, erotg_names
[rotg
->eType
], rotg
->slab_dist
);
841 fprintf(fp
, "# The scalar tau is the torque [kJ/mol] in the direction of the rotation vector.\n");
842 fprintf(fp
, "# To obtain the vectorial torque, multiply tau with\n");
843 fprintf(fp
, "# rot_vec%d %10.3e %10.3e %10.3e\n", g
, rotg
->vec
[XX
], rotg
->vec
[YY
], rotg
->vec
[ZZ
]);
847 fprintf(fp
, "# The following columns will have the syntax (tau=torque for that slab):\n");
849 print_aligned_short(fp
, "t");
850 print_aligned_short(fp
, "grp");
851 print_aligned_short(fp
, "slab");
852 print_aligned(fp
, "tau");
853 print_aligned_short(fp
, "slab");
854 print_aligned(fp
, "tau");
855 fprintf(fp
, " ...\n");
862 static void swap_val(double* vec
, int i
, int j
)
872 static void swap_col(double **mat
, int i
, int j
)
874 double tmp
[3] = {mat
[0][j
], mat
[1][j
], mat
[2][j
]};
887 /* Eigenvectors are stored in columns of eigen_vec */
888 static void diagonalize_symmetric(
896 jacobi(matrix
,3,eigenval
,eigen_vec
,&n_rot
);
898 /* sort in ascending order */
899 if (eigenval
[0] > eigenval
[1])
901 swap_val(eigenval
, 0, 1);
902 swap_col(eigen_vec
, 0, 1);
904 if (eigenval
[1] > eigenval
[2])
906 swap_val(eigenval
, 1, 2);
907 swap_col(eigen_vec
, 1, 2);
909 if (eigenval
[0] > eigenval
[1])
911 swap_val(eigenval
, 0, 1);
912 swap_col(eigen_vec
, 0, 1);
917 static void align_with_z(
918 rvec
* s
, /* Structure to align */
923 rvec zet
= {0.0, 0.0, 1.0};
924 rvec rot_axis
={0.0, 0.0, 0.0};
925 rvec
*rotated_str
=NULL
;
931 snew(rotated_str
, natoms
);
933 /* Normalize the axis */
934 ooanorm
= 1.0/norm(axis
);
935 svmul(ooanorm
, axis
, axis
);
937 /* Calculate the angle for the fitting procedure */
938 cprod(axis
, zet
, rot_axis
);
939 angle
= acos(axis
[2]);
943 /* Calculate the rotation matrix */
944 calc_rotmat(rot_axis
, angle
*180.0/M_PI
, rotmat
);
946 /* Apply the rotation matrix to s */
947 for (i
=0; i
<natoms
; i
++)
953 rotated_str
[i
][j
] += rotmat
[j
][k
]*s
[i
][k
];
958 /* Rewrite the rotated structure to s */
959 for(i
=0; i
<natoms
; i
++)
963 s
[i
][j
]=rotated_str
[i
][j
];
971 static void calc_correl_matrix(rvec
* Xstr
, rvec
* Ystr
, double** Rmat
, int natoms
)
982 for (k
=0; k
<natoms
; k
++)
983 Rmat
[i
][j
] += Ystr
[k
][i
] * Xstr
[k
][j
];
987 static void weigh_coords(rvec
* str
, real
* weight
, int natoms
)
992 for(i
=0; i
<natoms
; i
++)
995 str
[i
][j
] *= sqrt(weight
[i
]);
1000 static real
opt_angle_analytic(
1013 double **Rmat
, **RtR
, **eigvec
;
1015 double V
[3][3], WS
[3][3];
1016 double rot_matrix
[3][3];
1020 /* Do not change the original coordinates */
1021 snew(ref_s_1
, natoms
);
1022 snew(act_s_1
, natoms
);
1023 for(i
=0; i
<natoms
; i
++)
1025 copy_rvec(ref_s
[i
], ref_s_1
[i
]);
1026 copy_rvec(act_s
[i
], act_s_1
[i
]);
1029 /* Translate the structures to the origin */
1030 shift
[XX
] = -ref_com
[XX
];
1031 shift
[YY
] = -ref_com
[YY
];
1032 shift
[ZZ
] = -ref_com
[ZZ
];
1033 translate_x(ref_s_1
, natoms
, shift
);
1035 shift
[XX
] = -act_com
[XX
];
1036 shift
[YY
] = -act_com
[YY
];
1037 shift
[ZZ
] = -act_com
[ZZ
];
1038 translate_x(act_s_1
, natoms
, shift
);
1040 /* Align rotation axis with z */
1041 align_with_z(ref_s_1
, natoms
, axis
);
1042 align_with_z(act_s_1
, natoms
, axis
);
1044 /* Correlation matrix */
1045 Rmat
= allocate_square_matrix(3);
1047 for (i
=0; i
<natoms
; i
++)
1053 /* Weight positions with sqrt(weight) */
1056 weigh_coords(ref_s_1
, weight
, natoms
);
1057 weigh_coords(act_s_1
, weight
, natoms
);
1060 /* Calculate correlation matrices R=YXt (X=ref_s; Y=act_s) */
1061 calc_correl_matrix(ref_s_1
, act_s_1
, Rmat
, natoms
);
1064 RtR
= allocate_square_matrix(3);
1071 RtR
[i
][j
] += Rmat
[k
][i
] * Rmat
[k
][j
];
1075 /* Diagonalize RtR */
1080 diagonalize_symmetric(RtR
, eigvec
, eigval
);
1081 swap_col(eigvec
,0,1);
1082 swap_col(eigvec
,1,2);
1083 swap_val(eigval
,0,1);
1084 swap_val(eigval
,1,2);
1098 WS
[i
][j
] = eigvec
[i
][j
] / sqrt(eigval
[j
]);
1106 V
[i
][j
] += Rmat
[i
][k
]*WS
[k
][j
];
1110 free_square_matrix(Rmat
, 3);
1112 /* Calculate optimal rotation matrix */
1115 rot_matrix
[i
][j
] = 0.0;
1122 rot_matrix
[i
][j
] += eigvec
[i
][k
]*V
[j
][k
];
1126 rot_matrix
[2][2] = 1.0;
1128 /* In some cases abs(rot_matrix[0][0]) can be slighly larger
1129 * than unity due to numerical inacurracies. To be able to calculate
1130 * the acos function, we put these values back in range. */
1131 if (rot_matrix
[0][0] > 1.0)
1133 rot_matrix
[0][0] = 1.0;
1135 else if (rot_matrix
[0][0] < -1.0)
1137 rot_matrix
[0][0] = -1.0;
1140 /* Determine the optimal rotation angle: */
1141 opt_angle
= (-1.0)*acos(rot_matrix
[0][0])*180.0/M_PI
;
1142 if (rot_matrix
[0][1] < 0.0)
1143 opt_angle
= (-1.0)*opt_angle
;
1145 /* Give back some memory */
1146 free_square_matrix(RtR
, 3);
1153 return (real
) opt_angle
;
1157 /* Determine angle of the group by RMSD fit to the reference */
1158 /* Not parallelized, call this routine only on the master */
1159 static real
flex_fit_angle(t_rotgrp
*rotg
)
1162 rvec
*fitcoords
=NULL
;
1163 rvec center
; /* Center of positions passed to the fit routine */
1164 real fitangle
; /* Angle of the rotation group derived by fitting */
1167 gmx_enfrotgrp_t erg
; /* Pointer to enforced rotation group data */
1170 erg
=rotg
->enfrotgrp
;
1172 /* Get the center of the rotation group.
1173 * Note, again, erg->xc has been sorted in do_flexible */
1174 get_center(erg
->xc
, erg
->mc_sorted
, rotg
->nat
, center
);
1176 /* === Determine the optimal fit angle for the rotation group === */
1177 if (rotg
->eFittype
== erotgFitNORM
)
1179 /* Normalize every position to it's reference length */
1180 for (i
=0; i
<rotg
->nat
; i
++)
1182 /* Put the center of the positions into the origin */
1183 rvec_sub(erg
->xc
[i
], center
, coord
);
1184 /* Determine the scaling factor for the length: */
1185 scal
= erg
->xc_ref_length
[erg
->xc_sortind
[i
]] / norm(coord
);
1186 /* Get position, multiply with the scaling factor and save */
1187 svmul(scal
, coord
, erg
->xc_norm
[i
]);
1189 fitcoords
= erg
->xc_norm
;
1193 fitcoords
= erg
->xc
;
1195 /* From the point of view of the current positions, the reference has rotated
1196 * backwards. Since we output the angle relative to the fixed reference,
1197 * we need the minus sign. */
1198 fitangle
= -opt_angle_analytic(erg
->xc_ref_sorted
, fitcoords
, erg
->mc_sorted
,
1199 rotg
->nat
, erg
->xc_ref_center
, center
, rotg
->vec
);
1205 /* Determine actual angle of each slab by RMSD fit to the reference */
1206 /* Not parallelized, call this routine only on the master */
1207 static void flex_fit_angle_perslab(
1214 int i
,l
,n
,islab
,ind
;
1216 rvec act_center
; /* Center of actual positions that are passed to the fit routine */
1217 rvec ref_center
; /* Same for the reference positions */
1218 real fitangle
; /* Angle of a slab derived from an RMSD fit to
1219 * the reference structure at t=0 */
1221 gmx_enfrotgrp_t erg
; /* Pointer to enforced rotation group data */
1222 real OOm_av
; /* 1/average_mass of a rotation group atom */
1223 real m_rel
; /* Relative mass of a rotation group atom */
1226 erg
=rotg
->enfrotgrp
;
1228 /* Average mass of a rotation group atom: */
1229 OOm_av
= erg
->invmass
*rotg
->nat
;
1231 /**********************************/
1232 /* First collect the data we need */
1233 /**********************************/
1235 /* Collect the data for the individual slabs */
1236 for (n
= erg
->slab_first
; n
<= erg
->slab_last
; n
++)
1238 islab
= n
- erg
->slab_first
; /* slab index */
1239 sd
= &(rotg
->enfrotgrp
->slab_data
[islab
]);
1240 sd
->nat
= erg
->lastatom
[islab
]-erg
->firstatom
[islab
]+1;
1243 /* Loop over the relevant atoms in the slab */
1244 for (l
=erg
->firstatom
[islab
]; l
<=erg
->lastatom
[islab
]; l
++)
1246 /* Current position of this atom: x[ii][XX/YY/ZZ] */
1247 copy_rvec(erg
->xc
[l
], curr_x
);
1249 /* The (unrotated) reference position of this atom is copied to ref_x.
1250 * Beware, the xc coords have been sorted in do_flexible */
1251 copy_rvec(erg
->xc_ref_sorted
[l
], ref_x
);
1253 /* Save data for doing angular RMSD fit later */
1254 /* Save the current atom position */
1255 copy_rvec(curr_x
, sd
->x
[ind
]);
1256 /* Save the corresponding reference position */
1257 copy_rvec(ref_x
, sd
->ref
[ind
]);
1259 /* Maybe also mass-weighting was requested. If yes, additionally
1260 * multiply the weights with the relative mass of the atom. If not,
1261 * multiply with unity. */
1262 m_rel
= erg
->mc_sorted
[l
]*OOm_av
;
1264 /* Save the weight for this atom in this slab */
1265 sd
->weight
[ind
] = gaussian_weight(curr_x
, rotg
, n
) * m_rel
;
1267 /* Next atom in this slab */
1272 /******************************/
1273 /* Now do the fit calculation */
1274 /******************************/
1276 fprintf(fp
, "%12.3e%6d%12.3f", t
, g
, degangle
);
1278 /* === Now do RMSD fitting for each slab === */
1279 /* We require at least SLAB_MIN_ATOMS in a slab, such that the fit makes sense. */
1280 #define SLAB_MIN_ATOMS 4
1282 for (n
= erg
->slab_first
; n
<= erg
->slab_last
; n
++)
1284 islab
= n
- erg
->slab_first
; /* slab index */
1285 sd
= &(rotg
->enfrotgrp
->slab_data
[islab
]);
1286 if (sd
->nat
>= SLAB_MIN_ATOMS
)
1288 /* Get the center of the slabs reference and current positions */
1289 get_center(sd
->ref
, sd
->weight
, sd
->nat
, ref_center
);
1290 get_center(sd
->x
, sd
->weight
, sd
->nat
, act_center
);
1291 if (rotg
->eFittype
== erotgFitNORM
)
1293 /* Normalize every position to it's reference length
1294 * prior to performing the fit */
1295 for (i
=0; i
<sd
->nat
;i
++) /* Center */
1297 rvec_dec(sd
->ref
[i
], ref_center
);
1298 rvec_dec(sd
->x
[i
] , act_center
);
1299 /* Normalize x_i such that it gets the same length as ref_i */
1300 svmul( norm(sd
->ref
[i
])/norm(sd
->x
[i
]), sd
->x
[i
], sd
->x
[i
] );
1302 /* We already subtracted the centers */
1303 clear_rvec(ref_center
);
1304 clear_rvec(act_center
);
1306 fitangle
= -opt_angle_analytic(sd
->ref
, sd
->x
, sd
->weight
, sd
->nat
,
1307 ref_center
, act_center
, rotg
->vec
);
1308 fprintf(fp
, "%6d%6d%12.3f", n
, sd
->nat
, fitangle
);
1313 #undef SLAB_MIN_ATOMS
1317 /* Shift x with is */
1318 static inline void shift_single_coord(matrix box
, rvec x
, const ivec is
)
1329 x
[XX
] += tx
*box
[XX
][XX
]+ty
*box
[YY
][XX
]+tz
*box
[ZZ
][XX
];
1330 x
[YY
] += ty
*box
[YY
][YY
]+tz
*box
[ZZ
][YY
];
1331 x
[ZZ
] += tz
*box
[ZZ
][ZZ
];
1334 x
[XX
] += tx
*box
[XX
][XX
];
1335 x
[YY
] += ty
*box
[YY
][YY
];
1336 x
[ZZ
] += tz
*box
[ZZ
][ZZ
];
1341 /* Determine the 'home' slab of this atom which is the
1342 * slab with the highest Gaussian weight of all */
1343 #define round(a) (int)(a+0.5)
1344 static inline int get_homeslab(
1345 rvec curr_x
, /* The position for which the home slab shall be determined */
1346 rvec rotvec
, /* The rotation vector */
1347 real slabdist
) /* The slab distance */
1352 /* The distance of the atom to the coordinate center (where the
1353 * slab with index 0) is */
1354 dist
= iprod(rotvec
, curr_x
);
1356 return round(dist
/ slabdist
);
1360 /* For a local atom determine the relevant slabs, i.e. slabs in
1361 * which the gaussian is larger than min_gaussian
1363 static int get_single_atom_gaussians(
1371 gmx_enfrotgrp_t erg
; /* Pointer to enforced rotation group data */
1374 erg
=rotg
->enfrotgrp
;
1376 /* Determine the 'home' slab of this atom: */
1377 homeslab
= get_homeslab(curr_x
, rotg
->vec
, rotg
->slab_dist
);
1379 /* First determine the weight in the atoms home slab: */
1380 g
= gaussian_weight(curr_x
, rotg
, homeslab
);
1382 erg
->gn_atom
[count
] = g
;
1383 erg
->gn_slabind
[count
] = homeslab
;
1387 /* Determine the max slab */
1389 while (g
> rotg
->min_gaussian
)
1392 g
= gaussian_weight(curr_x
, rotg
, slab
);
1393 erg
->gn_slabind
[count
]=slab
;
1394 erg
->gn_atom
[count
]=g
;
1399 /* Determine the max slab */
1404 g
= gaussian_weight(curr_x
, rotg
, slab
);
1405 erg
->gn_slabind
[count
]=slab
;
1406 erg
->gn_atom
[count
]=g
;
1409 while (g
> rotg
->min_gaussian
);
1416 static void flex2_precalc_inner_sum(t_rotgrp
*rotg
, t_commrec
*cr
)
1419 rvec xi
; /* positions in the i-sum */
1420 rvec xcn
, ycn
; /* the current and the reference slab centers */
1423 rvec rin
; /* Helper variables */
1426 real OOpsii
,OOpsiistar
;
1427 real sin_rin
; /* s_ii.r_ii */
1428 rvec s_in
,tmpvec
,tmpvec2
;
1429 real mi
,wi
; /* Mass-weighting of the positions */
1431 gmx_enfrotgrp_t erg
; /* Pointer to enforced rotation group data */
1434 erg
=rotg
->enfrotgrp
;
1435 N_M
= rotg
->nat
* erg
->invmass
;
1437 /* Loop over all slabs that contain something */
1438 for (n
=erg
->slab_first
; n
<= erg
->slab_last
; n
++)
1440 islab
= n
- erg
->slab_first
; /* slab index */
1442 /* The current center of this slab is saved in xcn: */
1443 copy_rvec(erg
->slab_center
[islab
], xcn
);
1444 /* ... and the reference center in ycn: */
1445 copy_rvec(erg
->slab_center_ref
[islab
+erg
->slab_buffer
], ycn
);
1447 /*** D. Calculate the whole inner sum used for second and third sum */
1448 /* For slab n, we need to loop over all atoms i again. Since we sorted
1449 * the atoms with respect to the rotation vector, we know that it is sufficient
1450 * to calculate from firstatom to lastatom only. All other contributions will
1452 clear_rvec(innersumvec
);
1453 for (i
= erg
->firstatom
[islab
]; i
<= erg
->lastatom
[islab
]; i
++)
1455 /* Coordinate xi of this atom */
1456 copy_rvec(erg
->xc
[i
],xi
);
1459 gaussian_xi
= gaussian_weight(xi
,rotg
,n
);
1460 mi
= erg
->mc_sorted
[i
]; /* need the sorted mass here */
1464 copy_rvec(erg
->xc_ref_sorted
[i
],yi0
); /* Reference position yi0 */
1465 rvec_sub(yi0
, ycn
, tmpvec2
); /* tmpvec2 = yi0 - ycn */
1466 mvmul(erg
->rotmat
, tmpvec2
, rin
); /* rin = Omega.(yi0 - ycn) */
1468 /* Calculate psi_i* and sin */
1469 rvec_sub(xi
, xcn
, tmpvec2
); /* tmpvec2 = xi - xcn */
1470 cprod(rotg
->vec
, tmpvec2
, tmpvec
); /* tmpvec = v x (xi - xcn) */
1471 OOpsiistar
= norm2(tmpvec
)+rotg
->eps
; /* OOpsii* = 1/psii* = |v x (xi-xcn)|^2 + eps */
1472 OOpsii
= norm(tmpvec
); /* OOpsii = 1 / psii = |v x (xi - xcn)| */
1474 /* v x (xi - xcn) */
1475 unitv(tmpvec
, s_in
); /* sin = ---------------- */
1476 /* |v x (xi - xcn)| */
1478 sin_rin
=iprod(s_in
,rin
); /* sin_rin = sin . rin */
1480 /* Now the whole sum */
1481 fac
= OOpsii
/OOpsiistar
;
1482 svmul(fac
, rin
, tmpvec
);
1483 fac2
= fac
*fac
*OOpsii
;
1484 svmul(fac2
*sin_rin
, s_in
, tmpvec2
);
1485 rvec_dec(tmpvec
, tmpvec2
);
1487 svmul(wi
*gaussian_xi
*sin_rin
, tmpvec
, tmpvec2
);
1489 rvec_inc(innersumvec
,tmpvec2
);
1490 } /* now we have the inner sum, used both for sum2 and sum3 */
1492 /* Save it to be used in do_flex2_lowlevel */
1493 copy_rvec(innersumvec
, erg
->slab_innersumvec
[islab
]);
1494 } /* END of loop over slabs */
1498 static void flex_precalc_inner_sum(t_rotgrp
*rotg
, t_commrec
*cr
)
1501 rvec xi
; /* position */
1502 rvec xcn
, ycn
; /* the current and the reference slab centers */
1503 rvec qin
,rin
; /* q_i^n and r_i^n */
1506 rvec innersumvec
; /* Inner part of sum_n2 */
1507 real gaussian_xi
; /* Gaussian weight gn(xi) */
1508 real mi
,wi
; /* Mass-weighting of the positions */
1511 gmx_enfrotgrp_t erg
; /* Pointer to enforced rotation group data */
1514 erg
=rotg
->enfrotgrp
;
1515 N_M
= rotg
->nat
* erg
->invmass
;
1517 /* Loop over all slabs that contain something */
1518 for (n
=erg
->slab_first
; n
<= erg
->slab_last
; n
++)
1520 islab
= n
- erg
->slab_first
; /* slab index */
1522 /* The current center of this slab is saved in xcn: */
1523 copy_rvec(erg
->slab_center
[islab
], xcn
);
1524 /* ... and the reference center in ycn: */
1525 copy_rvec(erg
->slab_center_ref
[islab
+erg
->slab_buffer
], ycn
);
1527 /* For slab n, we need to loop over all atoms i again. Since we sorted
1528 * the atoms with respect to the rotation vector, we know that it is sufficient
1529 * to calculate from firstatom to lastatom only. All other contributions will
1531 clear_rvec(innersumvec
);
1532 for (i
=erg
->firstatom
[islab
]; i
<=erg
->lastatom
[islab
]; i
++)
1534 /* Coordinate xi of this atom */
1535 copy_rvec(erg
->xc
[i
],xi
);
1538 gaussian_xi
= gaussian_weight(xi
,rotg
,n
);
1539 mi
= erg
->mc_sorted
[i
]; /* need the sorted mass here */
1542 /* Calculate rin and qin */
1543 rvec_sub(erg
->xc_ref_sorted
[i
], ycn
, tmpvec
); /* tmpvec = yi0-ycn */
1544 mvmul(erg
->rotmat
, tmpvec
, rin
); /* rin = Omega.(yi0 - ycn) */
1545 cprod(rotg
->vec
, rin
, tmpvec
); /* tmpvec = v x Omega*(yi0-ycn) */
1547 /* v x Omega*(yi0-ycn) */
1548 unitv(tmpvec
, qin
); /* qin = --------------------- */
1549 /* |v x Omega*(yi0-ycn)| */
1552 rvec_sub(xi
, xcn
, tmpvec
); /* tmpvec = xi-xcn */
1553 bin
= iprod(qin
, tmpvec
); /* bin = qin*(xi-xcn) */
1555 svmul(wi
*gaussian_xi
*bin
, qin
, tmpvec
);
1557 /* Add this contribution to the inner sum: */
1558 rvec_add(innersumvec
, tmpvec
, innersumvec
);
1559 } /* now we have the inner sum vector S^n for this slab */
1560 /* Save it to be used in do_flex_lowlevel */
1561 copy_rvec(innersumvec
, erg
->slab_innersumvec
[islab
]);
1566 static real
do_flex2_lowlevel(
1568 real sigma
, /* The Gaussian width sigma */
1570 gmx_bool bCalcTorque
,
1574 int count
,ic
,ii
,j
,m
,n
,islab
,iigrp
;
1575 rvec xj
; /* position in the i-sum */
1576 rvec yj0
; /* the reference position in the j-sum */
1577 rvec xcn
, ycn
; /* the current and the reference slab centers */
1578 real V
; /* This node's part of the rotation pot. energy */
1579 real gaussian_xj
; /* Gaussian weight */
1583 rvec rjn
; /* Helper variables */
1586 real OOpsij
,OOpsijstar
;
1587 real OOsigma2
; /* 1/(sigma^2) */
1590 rvec sjn
,tmpvec
,tmpvec2
;
1591 rvec sum1vec_part
,sum1vec
,sum2vec_part
,sum2vec
,sum3vec
,sum4vec
,innersumvec
;
1593 gmx_enfrotgrp_t erg
; /* Pointer to enforced rotation group data */
1594 real mj
,wj
; /* Mass-weighting of the positions */
1596 real Wjn
; /* g_n(x_j) m_j / Mjn */
1598 /* To calculate the torque per slab */
1599 rvec slab_force
; /* Single force from slab n on one atom */
1600 rvec slab_sum1vec_part
;
1601 real slab_sum3part
,slab_sum4part
;
1602 rvec slab_sum1vec
, slab_sum2vec
, slab_sum3vec
, slab_sum4vec
;
1605 erg
=rotg
->enfrotgrp
;
1607 /* Pre-calculate the inner sums, so that we do not have to calculate
1608 * them again for every atom */
1609 flex2_precalc_inner_sum(rotg
, cr
);
1611 /********************************************************/
1612 /* Main loop over all local atoms of the rotation group */
1613 /********************************************************/
1614 N_M
= rotg
->nat
* erg
->invmass
;
1616 OOsigma2
= 1.0 / (sigma
*sigma
);
1617 for (j
=0; j
<erg
->nat_loc
; j
++)
1619 /* Local index of a rotation group atom */
1620 ii
= erg
->ind_loc
[j
];
1621 /* Position of this atom in the collective array */
1622 iigrp
= erg
->xc_ref_ind
[j
];
1623 /* Mass-weighting */
1624 mj
= erg
->mc
[iigrp
]; /* need the unsorted mass here */
1627 /* Current position of this atom: x[ii][XX/YY/ZZ]
1628 * Note that erg->xc_center contains the center of mass in case the flex2-t
1629 * potential was chosen. For the flex2 potential erg->xc_center must be
1631 rvec_sub(x
[ii
], erg
->xc_center
, xj
);
1633 /* Shift this atom such that it is near its reference */
1634 shift_single_coord(box
, xj
, erg
->xc_shifts
[iigrp
]);
1636 /* Determine the slabs to loop over, i.e. the ones with contributions
1637 * larger than min_gaussian */
1638 count
= get_single_atom_gaussians(xj
, cr
, rotg
);
1640 clear_rvec(sum1vec_part
);
1641 clear_rvec(sum2vec_part
);
1644 /* Loop over the relevant slabs for this atom */
1645 for (ic
=0; ic
< count
; ic
++)
1647 n
= erg
->gn_slabind
[ic
];
1649 /* Get the precomputed Gaussian value of curr_slab for curr_x */
1650 gaussian_xj
= erg
->gn_atom
[ic
];
1652 islab
= n
- erg
->slab_first
; /* slab index */
1654 /* The (unrotated) reference position of this atom is copied to yj0: */
1655 copy_rvec(rotg
->x_ref
[iigrp
], yj0
);
1657 beta
= calc_beta(xj
, rotg
,n
);
1659 /* The current center of this slab is saved in xcn: */
1660 copy_rvec(erg
->slab_center
[islab
], xcn
);
1661 /* ... and the reference center in ycn: */
1662 copy_rvec(erg
->slab_center_ref
[islab
+erg
->slab_buffer
], ycn
);
1664 rvec_sub(yj0
, ycn
, tmpvec2
); /* tmpvec2 = yj0 - ycn */
1667 mvmul(erg
->rotmat
, tmpvec2
, rjn
); /* rjn = Omega.(yj0 - ycn) */
1669 /* Subtract the slab center from xj */
1670 rvec_sub(xj
, xcn
, tmpvec2
); /* tmpvec2 = xj - xcn */
1673 cprod(rotg
->vec
, tmpvec2
, tmpvec
); /* tmpvec = v x (xj - xcn) */
1675 OOpsijstar
= norm2(tmpvec
)+rotg
->eps
; /* OOpsij* = 1/psij* = |v x (xj-xcn)|^2 + eps */
1677 numerator
= sqr(iprod(tmpvec
, rjn
));
1679 /*********************************/
1680 /* Add to the rotation potential */
1681 /*********************************/
1682 V
+= 0.5*rotg
->k
*wj
*gaussian_xj
*numerator
/OOpsijstar
;
1685 /*************************************/
1686 /* Now calculate the force on atom j */
1687 /*************************************/
1689 OOpsij
= norm(tmpvec
); /* OOpsij = 1 / psij = |v x (xj - xcn)| */
1691 /* v x (xj - xcn) */
1692 unitv(tmpvec
, sjn
); /* sjn = ---------------- */
1693 /* |v x (xj - xcn)| */
1695 sjn_rjn
=iprod(sjn
,rjn
); /* sjn_rjn = sjn . rjn */
1698 /*** A. Calculate the first of the four sum terms: ****************/
1699 fac
= OOpsij
/OOpsijstar
;
1700 svmul(fac
, rjn
, tmpvec
);
1701 fac2
= fac
*fac
*OOpsij
;
1702 svmul(fac2
*sjn_rjn
, sjn
, tmpvec2
);
1703 rvec_dec(tmpvec
, tmpvec2
);
1704 fac2
= wj
*gaussian_xj
; /* also needed for sum4 */
1705 svmul(fac2
*sjn_rjn
, tmpvec
, slab_sum1vec_part
);
1706 /********************/
1707 /*** Add to sum1: ***/
1708 /********************/
1709 rvec_inc(sum1vec_part
, slab_sum1vec_part
); /* sum1 still needs to vector multiplied with v */
1711 /*** B. Calculate the forth of the four sum terms: ****************/
1712 betasigpsi
= beta
*OOsigma2
*OOpsij
; /* this is also needed for sum3 */
1713 /********************/
1714 /*** Add to sum4: ***/
1715 /********************/
1716 slab_sum4part
= fac2
*betasigpsi
*fac
*sjn_rjn
*sjn_rjn
; /* Note that fac is still valid from above */
1717 sum4
+= slab_sum4part
;
1719 /*** C. Calculate Wjn for second and third sum */
1720 /* Note that we can safely divide by slab_weights since we check in
1721 * get_slab_centers that it is non-zero. */
1722 Wjn
= gaussian_xj
*mj
/erg
->slab_weights
[islab
];
1724 /* We already have precalculated the inner sum for slab n */
1725 copy_rvec(erg
->slab_innersumvec
[islab
], innersumvec
);
1727 /* Weigh the inner sum vector with Wjn */
1728 svmul(Wjn
, innersumvec
, innersumvec
);
1730 /*** E. Calculate the second of the four sum terms: */
1731 /********************/
1732 /*** Add to sum2: ***/
1733 /********************/
1734 rvec_inc(sum2vec_part
, innersumvec
); /* sum2 still needs to be vector crossproduct'ed with v */
1736 /*** F. Calculate the third of the four sum terms: */
1737 slab_sum3part
= betasigpsi
* iprod(sjn
, innersumvec
);
1738 sum3
+= slab_sum3part
; /* still needs to be multiplied with v */
1740 /*** G. Calculate the torque on the local slab's axis: */
1744 cprod(slab_sum1vec_part
, rotg
->vec
, slab_sum1vec
);
1746 cprod(innersumvec
, rotg
->vec
, slab_sum2vec
);
1748 svmul(slab_sum3part
, rotg
->vec
, slab_sum3vec
);
1750 svmul(slab_sum4part
, rotg
->vec
, slab_sum4vec
);
1752 /* The force on atom ii from slab n only: */
1753 for (m
=0; m
<DIM
; m
++)
1754 slab_force
[m
] = rotg
->k
* (-slab_sum1vec
[m
] + slab_sum2vec
[m
] - slab_sum3vec
[m
] + 0.5*slab_sum4vec
[m
]);
1756 erg
->slab_torque_v
[islab
] += torque(rotg
->vec
, slab_force
, xj
, xcn
);
1758 } /* END of loop over slabs */
1760 /* Construct the four individual parts of the vector sum: */
1761 cprod(sum1vec_part
, rotg
->vec
, sum1vec
); /* sum1vec = { } x v */
1762 cprod(sum2vec_part
, rotg
->vec
, sum2vec
); /* sum2vec = { } x v */
1763 svmul(sum3
, rotg
->vec
, sum3vec
); /* sum3vec = { } . v */
1764 svmul(sum4
, rotg
->vec
, sum4vec
); /* sum4vec = { } . v */
1766 /* Store the additional force so that it can be added to the force
1767 * array after the normal forces have been evaluated */
1768 for (m
=0; m
<DIM
; m
++)
1769 erg
->f_rot_loc
[j
][m
] = rotg
->k
* (-sum1vec
[m
] + sum2vec
[m
] - sum3vec
[m
] + 0.5*sum4vec
[m
]);
1772 fprintf(stderr
, "sum1: %15.8f %15.8f %15.8f\n", -rotg
->k
*sum1vec
[XX
], -rotg
->k
*sum1vec
[YY
], -rotg
->k
*sum1vec
[ZZ
]);
1773 fprintf(stderr
, "sum2: %15.8f %15.8f %15.8f\n", rotg
->k
*sum2vec
[XX
], rotg
->k
*sum2vec
[YY
], rotg
->k
*sum2vec
[ZZ
]);
1774 fprintf(stderr
, "sum3: %15.8f %15.8f %15.8f\n", -rotg
->k
*sum3vec
[XX
], -rotg
->k
*sum3vec
[YY
], -rotg
->k
*sum3vec
[ZZ
]);
1775 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
]);
1780 } /* END of loop over local atoms */
1786 static real
do_flex_lowlevel(
1788 real sigma
, /* The Gaussian width sigma */
1790 gmx_bool bCalcTorque
,
1794 int count
,ic
,ii
,j
,m
,n
,islab
,iigrp
;
1795 rvec xj
,yj0
; /* current and reference position */
1796 rvec xcn
, ycn
; /* the current and the reference slab centers */
1797 rvec xj_xcn
; /* xj - xcn */
1798 rvec qjn
; /* q_i^n */
1799 rvec sum_n1
,sum_n2
; /* Two contributions to the rotation force */
1800 rvec innersumvec
; /* Inner part of sum_n2 */
1802 rvec force_n
; /* Single force from slab n on one atom */
1803 rvec force_n1
,force_n2
; /* First and second part of force_n */
1804 rvec tmpvec
,tmpvec2
,tmp_f
; /* Helper variables */
1805 real V
; /* The rotation potential energy */
1806 real OOsigma2
; /* 1/(sigma^2) */
1807 real beta
; /* beta_n(xj) */
1808 real bjn
; /* b_j^n */
1809 real gaussian_xj
; /* Gaussian weight gn(xj) */
1810 real betan_xj_sigma2
;
1811 real mj
,wj
; /* Mass-weighting of the positions */
1813 gmx_enfrotgrp_t erg
; /* Pointer to enforced rotation group data */
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 flex_precalc_inner_sum(rotg
, cr
);
1822 /********************************************************/
1823 /* Main loop over all local atoms of the rotation group */
1824 /********************************************************/
1825 OOsigma2
= 1.0/(sigma
*sigma
);
1826 N_M
= rotg
->nat
* erg
->invmass
;
1828 for (j
=0; j
<erg
->nat_loc
; j
++)
1830 /* Local index of a rotation group atom */
1831 ii
= erg
->ind_loc
[j
];
1832 /* Position of this atom in the collective array */
1833 iigrp
= erg
->xc_ref_ind
[j
];
1834 /* Mass-weighting */
1835 mj
= erg
->mc
[iigrp
]; /* need the unsorted mass here */
1838 /* Current position of this atom: x[ii][XX/YY/ZZ]
1839 * Note that erg->xc_center contains the center of mass in case the flex-t
1840 * potential was chosen. For the flex potential erg->xc_center must be
1842 rvec_sub(x
[ii
], erg
->xc_center
, xj
);
1844 /* Shift this atom such that it is near its reference */
1845 shift_single_coord(box
, xj
, erg
->xc_shifts
[iigrp
]);
1847 /* Determine the slabs to loop over, i.e. the ones with contributions
1848 * larger than min_gaussian */
1849 count
= get_single_atom_gaussians(xj
, cr
, rotg
);
1854 /* Loop over the relevant slabs for this atom */
1855 for (ic
=0; ic
< count
; ic
++)
1857 n
= erg
->gn_slabind
[ic
];
1859 /* Get the precomputed Gaussian for xj in slab n */
1860 gaussian_xj
= erg
->gn_atom
[ic
];
1862 islab
= n
- erg
->slab_first
; /* slab index */
1864 /* The (unrotated) reference position of this atom is saved in yj0: */
1865 copy_rvec(rotg
->x_ref
[iigrp
], yj0
);
1867 beta
= calc_beta(xj
, rotg
, n
);
1869 /* The current center of this slab is saved in xcn: */
1870 copy_rvec(erg
->slab_center
[islab
], xcn
);
1871 /* ... and the reference center in ycn: */
1872 copy_rvec(erg
->slab_center_ref
[islab
+erg
->slab_buffer
], ycn
);
1874 rvec_sub(yj0
, ycn
, tmpvec
); /* tmpvec = yj0 - ycn */
1877 mvmul(erg
->rotmat
, tmpvec
, tmpvec2
); /* tmpvec2 = Omega.(yj0-ycn) */
1879 /* Subtract the slab center from xj */
1880 rvec_sub(xj
, xcn
, xj_xcn
); /* xj_xcn = xj - xcn */
1883 cprod(rotg
->vec
, tmpvec2
, tmpvec
); /* tmpvec = v x Omega.(xj-xcn) */
1885 /* v x Omega.(xj-xcn) */
1886 unitv(tmpvec
,qjn
); /* qjn = -------------------- */
1887 /* |v x Omega.(xj-xcn)| */
1889 bjn
= iprod(qjn
, xj_xcn
); /* bjn = qjn * (xj - xcn) */
1891 /*********************************/
1892 /* Add to the rotation potential */
1893 /*********************************/
1894 V
+= 0.5*rotg
->k
*wj
*gaussian_xj
*sqr(bjn
);
1896 /****************************************************************/
1897 /* sum_n1 will typically be the main contribution to the force: */
1898 /****************************************************************/
1899 betan_xj_sigma2
= beta
*OOsigma2
; /* beta_n(xj)/sigma^2 */
1901 /* The next lines calculate
1902 * qjn - (bjn*beta(xj)/(2sigma^2))v */
1903 svmul(bjn
*0.5*betan_xj_sigma2
, rotg
->vec
, tmpvec2
);
1904 rvec_sub(qjn
,tmpvec2
,tmpvec
);
1906 /* Multiply with gn(xj)*bjn: */
1907 svmul(gaussian_xj
*bjn
,tmpvec
,tmpvec2
);
1910 rvec_inc(sum_n1
,tmpvec2
);
1912 /* We already have precalculated the Sn term for slab n */
1913 copy_rvec(erg
->slab_innersumvec
[islab
], s_n
);
1915 svmul(betan_xj_sigma2
*iprod(s_n
, xj_xcn
), rotg
->vec
, tmpvec
); /* tmpvec = ---------- s_n (xj-xcn) */
1918 rvec_sub(s_n
, tmpvec
, innersumvec
);
1920 /* We can safely divide by slab_weights since we check in get_slab_centers
1921 * that it is non-zero. */
1922 svmul(gaussian_xj
/erg
->slab_weights
[islab
], innersumvec
, innersumvec
);
1924 rvec_add(sum_n2
, innersumvec
, sum_n2
);
1926 GMX_MPE_LOG(ev_inner_loop_finish
);
1928 /* Calculate the torque: */
1931 /* The force on atom ii from slab n only: */
1932 svmul(-rotg
->k
*wj
, tmpvec2
, force_n1
); /* part 1 */
1933 svmul( rotg
->k
*mj
, innersumvec
, force_n2
); /* part 2 */
1934 rvec_add(force_n1
, force_n2
, force_n
);
1935 erg
->slab_torque_v
[islab
] += torque(rotg
->vec
, force_n
, xj
, xcn
);
1937 } /* END of loop over slabs */
1939 /* Put both contributions together: */
1940 svmul(wj
, sum_n1
, sum_n1
);
1941 svmul(mj
, sum_n2
, sum_n2
);
1942 rvec_sub(sum_n2
,sum_n1
,tmp_f
); /* F = -grad V */
1944 /* Store the additional force so that it can be added to the force
1945 * array after the normal forces have been evaluated */
1946 for(m
=0; m
<DIM
; m
++)
1947 erg
->f_rot_loc
[j
][m
] = rotg
->k
*tmp_f
[m
];
1951 } /* END of loop over local atoms */
1957 static void print_coordinates(t_commrec
*cr
, t_rotgrp
*rotg
, rvec x
[], matrix box
, int step
)
1961 static char buf
[STRLEN
];
1962 static gmx_bool bFirst
=1;
1967 sprintf(buf
, "coords%d.txt", cr
->nodeid
);
1968 fp
= fopen(buf
, "w");
1972 fprintf(fp
, "\nStep %d\n", step
);
1973 fprintf(fp
, "box: %f %f %f %f %f %f %f %f %f\n",
1974 box
[XX
][XX
], box
[XX
][YY
], box
[XX
][ZZ
],
1975 box
[YY
][XX
], box
[YY
][YY
], box
[YY
][ZZ
],
1976 box
[ZZ
][XX
], box
[ZZ
][ZZ
], box
[ZZ
][ZZ
]);
1977 for (i
=0; i
<rotg
->nat
; i
++)
1979 fprintf(fp
, "%4d %f %f %f\n", i
,
1980 erg
->xc
[i
][XX
], erg
->xc
[i
][YY
], erg
->xc
[i
][ZZ
]);
1988 static int projection_compare(const void *a
, const void *b
)
1990 sort_along_vec_t
*xca
, *xcb
;
1993 xca
= (sort_along_vec_t
*)a
;
1994 xcb
= (sort_along_vec_t
*)b
;
1996 if (xca
->xcproj
< xcb
->xcproj
)
1998 else if (xca
->xcproj
> xcb
->xcproj
)
2005 static void sort_collective_coordinates(
2006 t_rotgrp
*rotg
, /* Rotation group */
2007 sort_along_vec_t
*data
) /* Buffer for sorting the positions */
2010 gmx_enfrotgrp_t erg
; /* Pointer to enforced rotation group data */
2013 erg
=rotg
->enfrotgrp
;
2015 /* The projection of the position vector on the rotation vector is
2016 * the relevant value for sorting. Fill the 'data' structure */
2017 for (i
=0; i
<rotg
->nat
; i
++)
2019 data
[i
].xcproj
= iprod(erg
->xc
[i
], rotg
->vec
); /* sort criterium */
2020 data
[i
].m
= erg
->mc
[i
];
2022 copy_rvec(erg
->xc
[i
] , data
[i
].x
);
2023 copy_rvec(rotg
->x_ref
[i
], data
[i
].x_ref
);
2025 /* Sort the 'data' structure */
2026 gmx_qsort(data
, rotg
->nat
, sizeof(sort_along_vec_t
), projection_compare
);
2028 /* Copy back the sorted values */
2029 for (i
=0; i
<rotg
->nat
; i
++)
2031 copy_rvec(data
[i
].x
, erg
->xc
[i
] );
2032 copy_rvec(data
[i
].x_ref
, erg
->xc_ref_sorted
[i
]);
2033 erg
->mc_sorted
[i
] = data
[i
].m
;
2034 erg
->xc_sortind
[i
] = data
[i
].ind
;
2039 /* For each slab, get the first and the last index of the sorted atom
2041 static void get_firstlast_atom_per_slab(t_rotgrp
*rotg
, t_commrec
*cr
)
2045 gmx_enfrotgrp_t erg
; /* Pointer to enforced rotation group data */
2048 erg
=rotg
->enfrotgrp
;
2050 GMX_MPE_LOG(ev_get_firstlast_start
);
2052 /* Find the first atom that needs to enter the calculation for each slab */
2053 n
= erg
->slab_first
; /* slab */
2054 i
= 0; /* start with the first atom */
2057 /* Find the first atom that significantly contributes to this slab */
2058 do /* move forward in position until a large enough beta is found */
2060 beta
= calc_beta(erg
->xc
[i
], rotg
, n
);
2062 } while ((beta
< -erg
->max_beta
) && (i
< rotg
->nat
));
2064 islab
= n
- erg
->slab_first
; /* slab index */
2065 erg
->firstatom
[islab
] = i
;
2066 /* Proceed to the next slab */
2068 } while (n
<= erg
->slab_last
);
2070 /* Find the last atom for each slab */
2071 n
= erg
->slab_last
; /* start with last slab */
2072 i
= rotg
->nat
-1; /* start with the last atom */
2075 do /* move backward in position until a large enough beta is found */
2077 beta
= calc_beta(erg
->xc
[i
], rotg
, n
);
2079 } while ((beta
> erg
->max_beta
) && (i
> -1));
2081 islab
= n
- erg
->slab_first
; /* slab index */
2082 erg
->lastatom
[islab
] = i
;
2083 /* Proceed to the next slab */
2085 } while (n
>= erg
->slab_first
);
2087 GMX_MPE_LOG(ev_get_firstlast_finish
);
2091 /* Determine the very first and very last slab that needs to be considered
2092 * For the first slab that needs to be considered, we have to find the smallest
2095 * x_first * v - n*Delta_x <= beta_max
2097 * slab index n, slab distance Delta_x, rotation vector v. For the last slab we
2098 * have to find the largest n that obeys
2100 * x_last * v - n*Delta_x >= -beta_max
2103 static inline int get_first_slab(
2104 t_rotgrp
*rotg
, /* The rotation group (inputrec data) */
2105 real max_beta
, /* The max_beta value, instead of min_gaussian */
2106 rvec firstatom
) /* First atom after sorting along the rotation vector v */
2108 /* Find the first slab for the first atom */
2109 return ceil((iprod(firstatom
, rotg
->vec
) - max_beta
)/rotg
->slab_dist
);
2113 static inline int get_last_slab(
2114 t_rotgrp
*rotg
, /* The rotation group (inputrec data) */
2115 real max_beta
, /* The max_beta value, instead of min_gaussian */
2116 rvec lastatom
) /* Last atom along v */
2118 /* Find the last slab for the last atom */
2119 return floor((iprod(lastatom
, rotg
->vec
) + max_beta
)/rotg
->slab_dist
);
2123 static void get_firstlast_slab_check(
2124 t_rotgrp
*rotg
, /* The rotation group (inputrec data) */
2125 t_gmx_enfrotgrp
*erg
, /* The rotation group (data only accessible in this file) */
2126 rvec firstatom
, /* First atom after sorting along the rotation vector v */
2127 rvec lastatom
, /* Last atom along v */
2128 int g
, /* The rotation group number */
2131 erg
->slab_first
= get_first_slab(rotg
, erg
->max_beta
, firstatom
);
2132 erg
->slab_last
= get_last_slab(rotg
, erg
->max_beta
, lastatom
);
2134 /* Check whether we have reference data to compare against */
2135 if (erg
->slab_first
< erg
->slab_first_ref
)
2136 gmx_fatal(FARGS
, "%s No reference data for first slab (n=%d), unable to proceed.",
2137 RotStr
, erg
->slab_first
);
2139 /* Check whether we have reference data to compare against */
2140 if (erg
->slab_last
> erg
->slab_last_ref
)
2141 gmx_fatal(FARGS
, "%s No reference data for last slab (n=%d), unable to proceed.",
2142 RotStr
, erg
->slab_last
);
2146 /* Enforced rotation with a flexible axis */
2147 static void do_flexible(
2149 gmx_enfrot_t enfrot
, /* Other rotation data */
2150 t_rotgrp
*rotg
, /* The rotation group */
2151 int g
, /* Group number */
2152 rvec x
[], /* The local positions */
2154 double t
, /* Time in picoseconds */
2155 gmx_large_int_t step
, /* The time step */
2156 gmx_bool bOutstepRot
, /* Output to main rotation output file */
2157 gmx_bool bOutstepSlab
) /* Output per-slab data */
2160 real sigma
; /* The Gaussian width sigma */
2161 gmx_enfrotgrp_t erg
; /* Pointer to enforced rotation group data */
2164 erg
=rotg
->enfrotgrp
;
2166 /* Define the sigma value */
2167 sigma
= 0.7*rotg
->slab_dist
;
2169 /* Sort the collective coordinates erg->xc along the rotation vector. This is
2170 * an optimization for the inner loop. */
2171 sort_collective_coordinates(rotg
, enfrot
->data
);
2173 /* Determine the first relevant slab for the first atom and the last
2174 * relevant slab for the last atom */
2175 get_firstlast_slab_check(rotg
, erg
, erg
->xc
[0], erg
->xc
[rotg
->nat
-1], g
, cr
);
2177 /* Determine for each slab depending on the min_gaussian cutoff criterium,
2178 * a first and a last atom index inbetween stuff needs to be calculated */
2179 get_firstlast_atom_per_slab(rotg
, cr
);
2181 /* Determine the gaussian-weighted center of positions for all slabs */
2182 get_slab_centers(rotg
,erg
->xc
,erg
->mc_sorted
,cr
,g
,t
,enfrot
->out_slabs
,bOutstepSlab
,FALSE
);
2184 /* Clear the torque per slab from last time step: */
2185 nslabs
= erg
->slab_last
- erg
->slab_first
+ 1;
2186 for (l
=0; l
<nslabs
; l
++)
2187 erg
->slab_torque_v
[l
] = 0.0;
2189 /* Call the rotational forces kernel */
2190 GMX_MPE_LOG(ev_flexll_start
);
2191 if (rotg
->eType
== erotgFLEX
|| rotg
->eType
== erotgFLEXT
)
2192 erg
->V
= do_flex_lowlevel(rotg
, sigma
, x
, bOutstepRot
, box
, cr
);
2193 else if (rotg
->eType
== erotgFLEX2
|| rotg
->eType
== erotgFLEX2T
)
2194 erg
->V
= do_flex2_lowlevel(rotg
, sigma
, x
, bOutstepRot
, box
, cr
);
2196 gmx_fatal(FARGS
, "Unknown flexible rotation type");
2197 GMX_MPE_LOG(ev_flexll_finish
);
2199 /* Determine angle by RMSD fit to the reference - Let's hope this */
2200 /* only happens once in a while, since this is not parallelized! */
2205 /* Fit angle of the whole rotation group */
2206 erg
->angle_v
= flex_fit_angle(rotg
);
2210 /* Fit angle of each slab */
2211 flex_fit_angle_perslab(g
, rotg
, t
, erg
->degangle
, enfrot
->out_angles
);
2215 /* Lump together the torques from all slabs: */
2216 erg
->torque_v
= 0.0;
2217 for (l
=0; l
<nslabs
; l
++)
2218 erg
->torque_v
+= erg
->slab_torque_v
[l
];
2222 /* Calculate the angle between reference and actual rotation group atom,
2223 * both projected into a plane perpendicular to the rotation vector: */
2224 static void angle(t_rotgrp
*rotg
,
2228 real
*weight
) /* atoms near the rotation axis should count less than atoms far away */
2230 rvec xp
, xrp
; /* current and reference positions projected on a plane perpendicular to pg->vec */
2234 /* Project x_ref and x into a plane through the origin perpendicular to rot_vec: */
2235 /* Project x_ref: xrp = x_ref - (vec * x_ref) * vec */
2236 svmul(iprod(rotg
->vec
, x_ref
), rotg
->vec
, dum
);
2237 rvec_sub(x_ref
, dum
, xrp
);
2238 /* Project x_act: */
2239 svmul(iprod(rotg
->vec
, x_act
), rotg
->vec
, dum
);
2240 rvec_sub(x_act
, dum
, xp
);
2242 /* Retrieve information about which vector precedes. gmx_angle always
2243 * returns a positive angle. */
2244 cprod(xp
, xrp
, dum
); /* if reference precedes, this is pointing into the same direction as vec */
2246 if (iprod(rotg
->vec
, dum
) >= 0)
2247 *alpha
= -gmx_angle(xrp
, xp
);
2249 *alpha
= +gmx_angle(xrp
, xp
);
2251 /* Also return the weight */
2256 /* Project first vector onto a plane perpendicular to the second vector
2258 * Note that v must be of unit length.
2260 static inline void project_onto_plane(rvec dr
, const rvec v
)
2265 svmul(iprod(dr
,v
),v
,tmp
); /* tmp = (dr.v)v */
2266 rvec_dec(dr
, tmp
); /* dr = dr - (dr.v)v */
2270 /* Fixed rotation: The rotation reference group rotates around an axis */
2271 /* The atoms of the actual rotation group are attached with imaginary */
2272 /* springs to the reference atoms. */
2273 static void do_fixed(
2275 t_rotgrp
*rotg
, /* The rotation group */
2276 rvec x
[], /* The positions */
2277 matrix box
, /* The simulation box */
2278 double t
, /* Time in picoseconds */
2279 gmx_large_int_t step
, /* The time step */
2284 rvec tmp_f
; /* Force */
2285 real alpha
; /* a single angle between an actual and a reference position */
2286 real weight
; /* single weight for a single angle */
2287 gmx_enfrotgrp_t erg
; /* Pointer to enforced rotation group data */
2290 /* for mass weighting: */
2291 real wi
; /* Mass-weighting of the positions */
2293 real k_wi
; /* k times wi */
2298 erg
=rotg
->enfrotgrp
;
2299 bProject
= (rotg
->eType
==erotgPM
) || (rotg
->eType
==erotgPMPF
);
2301 N_M
= rotg
->nat
* erg
->invmass
;
2303 /* Each process calculates the forces on its local atoms */
2304 for (j
=0; j
<erg
->nat_loc
; j
++)
2306 /* Calculate (x_i-x_c) resp. (x_i-u) */
2307 rvec_sub(erg
->x_loc_pbc
[j
], erg
->xc_center
, tmpvec
);
2309 /* Calculate Omega*(y_i-y_c)-(x_i-x_c) */
2310 rvec_sub(erg
->xr_loc
[j
], tmpvec
, dr
);
2313 project_onto_plane(dr
, rotg
->vec
);
2315 /* Mass-weighting */
2316 wi
= N_M
*erg
->m_loc
[j
];
2318 /* Store the additional force so that it can be added to the force
2319 * array after the normal forces have been evaluated */
2321 for (m
=0; m
<DIM
; m
++)
2323 tmp_f
[m
] = k_wi
*dr
[m
];
2324 erg
->f_rot_loc
[j
][m
] = tmp_f
[m
];
2325 erg
->V
+= 0.5*k_wi
*sqr(dr
[m
]);
2330 /* Add to the torque of this rotation group */
2331 erg
->torque_v
+= torque(rotg
->vec
, tmp_f
, erg
->x_loc_pbc
[j
], erg
->xc_center
);
2333 /* Calculate the angle between reference and actual rotation group atom. */
2334 angle(rotg
, tmpvec
, erg
->xr_loc
[j
], &alpha
, &weight
); /* angle in rad, weighted */
2335 erg
->angle_v
+= alpha
* weight
;
2336 erg
->weight_v
+= weight
;
2338 /* If you want enforced rotation to contribute to the virial,
2339 * activate the following lines:
2342 Add the rotation contribution to the virial
2343 for(j=0; j<DIM; j++)
2345 vir[j][m] += 0.5*f[ii][j]*dr[m];
2351 } /* end of loop over local rotation group atoms */
2355 /* Calculate the radial motion potential and forces */
2356 static void do_radial_motion(
2358 t_rotgrp
*rotg
, /* The rotation group */
2359 rvec x
[], /* The positions */
2360 matrix box
, /* The simulation box */
2361 double t
, /* Time in picoseconds */
2362 gmx_large_int_t step
, /* The time step */
2366 rvec tmp_f
; /* Force */
2367 real alpha
; /* a single angle between an actual and a reference position */
2368 real weight
; /* single weight for a single angle */
2369 gmx_enfrotgrp_t erg
; /* Pointer to enforced rotation group data */
2370 rvec xj_u
; /* xj - u */
2372 real fac
,fac2
,sum
=0.0;
2375 /* For mass weighting: */
2376 real wj
; /* Mass-weighting of the positions */
2380 erg
=rotg
->enfrotgrp
;
2382 N_M
= rotg
->nat
* erg
->invmass
;
2384 /* Each process calculates the forces on its local atoms */
2385 for (j
=0; j
<erg
->nat_loc
; j
++)
2387 /* Calculate (xj-u) */
2388 rvec_sub(erg
->x_loc_pbc
[j
], erg
->xc_center
, xj_u
); /* xj_u = xj-u */
2390 /* Calculate Omega.(yj-u) */
2391 cprod(rotg
->vec
, erg
->xr_loc
[j
], tmpvec
); /* tmpvec = v x Omega.(yj-u) */
2393 /* v x Omega.(yj-u) */
2394 unitv(tmpvec
, pj
); /* pj = -------------------- */
2395 /* | v x Omega.(yj-u) | */
2397 fac
= iprod(pj
, xj_u
); /* fac = pj.(xj-u) */
2400 /* Mass-weighting */
2401 wj
= N_M
*erg
->m_loc
[j
];
2403 /* Store the additional force so that it can be added to the force
2404 * array after the normal forces have been evaluated */
2405 svmul(-rotg
->k
*wj
*fac
, pj
, tmp_f
);
2406 copy_rvec(tmp_f
, erg
->f_rot_loc
[j
]);
2410 /* Add to the torque of this rotation group */
2411 erg
->torque_v
+= torque(rotg
->vec
, tmp_f
, erg
->x_loc_pbc
[j
], erg
->xc_center
);
2413 /* Calculate the angle between reference and actual rotation group atom. */
2414 angle(rotg
, xj_u
, erg
->xr_loc
[j
], &alpha
, &weight
); /* angle in rad, weighted */
2415 erg
->angle_v
+= alpha
* weight
;
2416 erg
->weight_v
+= weight
;
2421 } /* end of loop over local rotation group atoms */
2422 erg
->V
= 0.5*rotg
->k
*sum
;
2426 /* Calculate the radial motion pivot-free potential and forces */
2427 static void do_radial_motion_pf(
2429 t_rotgrp
*rotg
, /* The rotation group */
2430 rvec x
[], /* The positions */
2431 matrix box
, /* The simulation box */
2432 double t
, /* Time in picoseconds */
2433 gmx_large_int_t step
, /* The time step */
2437 rvec xj
; /* Current position */
2438 rvec xj_xc
; /* xj - xc */
2439 rvec yj0_yc0
; /* yj0 - yc0 */
2440 rvec tmp_f
; /* Force */
2441 real alpha
; /* a single angle between an actual and a reference position */
2442 real weight
; /* single weight for a single angle */
2443 gmx_enfrotgrp_t erg
; /* Pointer to enforced rotation group data */
2444 rvec tmpvec
, tmpvec2
;
2445 rvec innersumvec
; /* Precalculation of the inner sum */
2447 real fac
,fac2
,V
=0.0;
2450 /* For mass weighting: */
2451 real mj
,wi
,wj
; /* Mass-weighting of the positions */
2455 erg
=rotg
->enfrotgrp
;
2457 N_M
= rotg
->nat
* erg
->invmass
;
2459 /* Get the current center of the rotation group: */
2460 get_center(erg
->xc
, erg
->mc
, rotg
->nat
, erg
->xc_center
);
2462 /* Precalculate Sum_i [ wi qi.(xi-xc) qi ] which is needed for every single j */
2463 clear_rvec(innersumvec
);
2464 for (i
=0; i
< rotg
->nat
; i
++)
2466 /* Mass-weighting */
2467 wi
= N_M
*erg
->mc
[i
];
2469 /* Calculate qi. Note that xc_ref_center has already been subtracted from
2470 * x_ref in init_rot_group.*/
2471 mvmul(erg
->rotmat
, rotg
->x_ref
[i
], tmpvec
); /* tmpvec = Omega.(yi0-yc0) */
2473 cprod(rotg
->vec
, tmpvec
, tmpvec2
); /* tmpvec2 = v x Omega.(yi0-yc0) */
2475 /* v x Omega.(yi0-yc0) */
2476 unitv(tmpvec2
, qi
); /* qi = ----------------------- */
2477 /* | v x Omega.(yi0-yc0) | */
2479 rvec_sub(erg
->xc
[i
], erg
->xc_center
, tmpvec
); /* tmpvec = xi-xc */
2481 svmul(wi
*iprod(qi
, tmpvec
), qi
, tmpvec2
);
2483 rvec_inc(innersumvec
, tmpvec2
);
2485 svmul(rotg
->k
*erg
->invmass
, innersumvec
, innersumveckM
);
2487 /* Each process calculates the forces on its local atoms */
2488 for (j
=0; j
<erg
->nat_loc
; j
++)
2490 /* Local index of a rotation group atom */
2491 ii
= erg
->ind_loc
[j
];
2492 /* Position of this atom in the collective array */
2493 iigrp
= erg
->xc_ref_ind
[j
];
2494 /* Mass-weighting */
2495 mj
= erg
->mc
[iigrp
]; /* need the unsorted mass here */
2498 /* Current position of this atom: x[ii][XX/YY/ZZ] */
2499 copy_rvec(x
[ii
], xj
);
2501 /* Shift this atom such that it is near its reference */
2502 shift_single_coord(box
, xj
, erg
->xc_shifts
[iigrp
]);
2504 /* The (unrotated) reference position is yj0. yc0 has already
2505 * been subtracted in init_rot_group */
2506 copy_rvec(rotg
->x_ref
[iigrp
], yj0_yc0
); /* yj0_yc0 = yj0 - yc0 */
2508 /* Calculate Omega.(yj0-yc0) */
2509 mvmul(erg
->rotmat
, yj0_yc0
, tmpvec2
); /* tmpvec2 = Omega.(yj0 - yc0) */
2511 cprod(rotg
->vec
, tmpvec2
, tmpvec
); /* tmpvec = v x Omega.(yj0-yc0) */
2513 /* v x Omega.(yj0-yc0) */
2514 unitv(tmpvec
, qj
); /* qj = ----------------------- */
2515 /* | v x Omega.(yj0-yc0) | */
2517 /* Calculate (xj-xc) */
2518 rvec_sub(xj
, erg
->xc_center
, xj_xc
); /* xj_xc = xj-xc */
2520 fac
= iprod(qj
, xj_xc
); /* fac = qj.(xj-xc) */
2523 /* Store the additional force so that it can be added to the force
2524 * array after the normal forces have been evaluated */
2525 svmul(-rotg
->k
*wj
*fac
, qj
, tmp_f
); /* part 1 of force */
2526 svmul(mj
, innersumveckM
, tmpvec
); /* part 2 of force */
2527 rvec_inc(tmp_f
, tmpvec
);
2528 copy_rvec(tmp_f
, erg
->f_rot_loc
[j
]);
2532 /* Add to the torque of this rotation group */
2533 erg
->torque_v
+= torque(rotg
->vec
, tmp_f
, xj
, erg
->xc_center
);
2535 /* Calculate the angle between reference and actual rotation group atom. */
2536 angle(rotg
, xj_xc
, yj0_yc0
, &alpha
, &weight
); /* angle in rad, weighted */
2537 erg
->angle_v
+= alpha
* weight
;
2538 erg
->weight_v
+= weight
;
2543 } /* end of loop over local rotation group atoms */
2544 erg
->V
= 0.5*rotg
->k
*V
;
2548 /* Precalculate the inner sum for the radial motion 2 forces */
2549 static void radial_motion2_precalc_inner_sum(t_rotgrp
*rotg
, rvec innersumvec
)
2552 gmx_enfrotgrp_t erg
; /* Pointer to enforced rotation group data */
2553 rvec xi_xc
; /* xj - xc */
2554 rvec tmpvec
,tmpvec2
;
2558 rvec v_xi_xc
; /* v x (xj - u) */
2560 real wi
; /* Mass-weighting of the positions */
2564 erg
=rotg
->enfrotgrp
;
2565 N_M
= rotg
->nat
* erg
->invmass
;
2567 /* Loop over the collective set of positions */
2569 for (i
=0; i
<rotg
->nat
; i
++)
2571 /* Mass-weighting */
2572 wi
= N_M
*erg
->mc
[i
];
2574 rvec_sub(erg
->xc
[i
], erg
->xc_center
, xi_xc
); /* xi_xc = xi-xc */
2576 /* Calculate ri. Note that xc_ref_center has already been subtracted from
2577 * x_ref in init_rot_group.*/
2578 mvmul(erg
->rotmat
, rotg
->x_ref
[i
], ri
); /* ri = Omega.(yi0-yc0) */
2580 cprod(rotg
->vec
, xi_xc
, v_xi_xc
); /* v_xi_xc = v x (xi-u) */
2582 fac
= norm2(v_xi_xc
);
2584 psiistar
= 1.0/(fac
+ rotg
->eps
); /* psiistar = --------------------- */
2585 /* |v x (xi-xc)|^2 + eps */
2587 psii
= gmx_invsqrt(fac
); /* 1 */
2588 /* psii = ------------- */
2591 svmul(psii
, v_xi_xc
, si
); /* si = psii * (v x (xi-xc) ) */
2593 fac
= iprod(v_xi_xc
, ri
); /* fac = (v x (xi-xc)).ri */
2596 siri
= iprod(si
, ri
); /* siri = si.ri */
2598 svmul(psiistar
/psii
, ri
, tmpvec
);
2599 svmul(psiistar
*psiistar
/(psii
*psii
*psii
) * siri
, si
, tmpvec2
);
2600 rvec_dec(tmpvec
, tmpvec2
);
2601 cprod(tmpvec
, rotg
->vec
, tmpvec2
);
2603 svmul(wi
*siri
, tmpvec2
, tmpvec
);
2605 rvec_inc(sumvec
, tmpvec
);
2607 svmul(rotg
->k
*erg
->invmass
, sumvec
, innersumvec
);
2611 /* Calculate the radial motion 2 potential and forces */
2612 static void do_radial_motion2(
2614 t_rotgrp
*rotg
, /* The rotation group */
2615 rvec x
[], /* The positions */
2616 matrix box
, /* The simulation box */
2617 double t
, /* Time in picoseconds */
2618 gmx_large_int_t step
, /* The time step */
2622 rvec xj
; /* Position */
2623 real alpha
; /* a single angle between an actual and a reference position */
2624 real weight
; /* single weight for a single angle */
2625 gmx_enfrotgrp_t erg
; /* Pointer to enforced rotation group data */
2626 rvec xj_u
; /* xj - u */
2627 rvec tmpvec
,tmpvec2
;
2628 real fac
,fac2
,Vpart
=0.0;
2631 rvec v_xj_u
; /* v x (xj - u) */
2633 real mj
,wj
; /* For mass-weighting of the positions */
2639 erg
=rotg
->enfrotgrp
;
2641 bPF
= rotg
->eType
==erotgRM2PF
;
2642 clear_rvec(innersumvec
);
2645 /* For the pivot-free variant we have to use the current center of
2646 * mass of the rotation group instead of the pivot u */
2647 get_center(erg
->xc
, erg
->mc
, rotg
->nat
, erg
->xc_center
);
2649 /* Also, we precalculate the second term of the forces that is identical
2650 * (up to the weight factor mj) for all forces */
2651 radial_motion2_precalc_inner_sum(rotg
,innersumvec
);
2654 N_M
= rotg
->nat
* erg
->invmass
;
2656 /* Each process calculates the forces on its local atoms */
2657 for (j
=0; j
<erg
->nat_loc
; j
++)
2661 /* Local index of a rotation group atom */
2662 ii
= erg
->ind_loc
[j
];
2663 /* Position of this atom in the collective array */
2664 iigrp
= erg
->xc_ref_ind
[j
];
2665 /* Mass-weighting */
2666 mj
= erg
->mc
[iigrp
];
2668 /* Current position of this atom: x[ii] */
2669 copy_rvec(x
[ii
], xj
);
2671 /* Shift this atom such that it is near its reference */
2672 shift_single_coord(box
, xj
, erg
->xc_shifts
[iigrp
]);
2674 /* The (unrotated) reference position is yj0. yc0 has already
2675 * been subtracted in init_rot_group */
2676 copy_rvec(rotg
->x_ref
[iigrp
], tmpvec
); /* tmpvec = yj0 - yc0 */
2678 /* Calculate Omega.(yj0-yc0) */
2679 mvmul(erg
->rotmat
, tmpvec
, rj
); /* rj = Omega.(yj0-yc0) */
2684 copy_rvec(erg
->x_loc_pbc
[j
], xj
);
2685 copy_rvec(erg
->xr_loc
[j
], rj
); /* rj = Omega.(yj0-u) */
2687 /* Mass-weighting */
2690 /* Calculate (xj-u) resp. (xj-xc) */
2691 rvec_sub(xj
, erg
->xc_center
, xj_u
); /* xj_u = xj-u */
2693 cprod(rotg
->vec
, xj_u
, v_xj_u
); /* v_xj_u = v x (xj-u) */
2695 fac
= norm2(v_xj_u
);
2697 psijstar
= 1.0/(fac
+ rotg
->eps
); /* psistar = -------------------- */
2698 /* |v x (xj-u)|^2 + eps */
2700 psij
= gmx_invsqrt(fac
); /* 1 */
2701 /* psij = ------------ */
2704 svmul(psij
, v_xj_u
, sj
); /* sj = psij * (v x (xj-u) ) */
2706 fac
= iprod(v_xj_u
, rj
); /* fac = (v x (xj-u)).rj */
2709 sjrj
= iprod(sj
, rj
); /* sjrj = sj.rj */
2711 svmul(psijstar
/psij
, rj
, tmpvec
);
2712 svmul(psijstar
*psijstar
/(psij
*psij
*psij
) * sjrj
, sj
, tmpvec2
);
2713 rvec_dec(tmpvec
, tmpvec2
);
2714 cprod(tmpvec
, rotg
->vec
, tmpvec2
);
2716 /* Store the additional force so that it can be added to the force
2717 * array after the normal forces have been evaluated */
2718 svmul(-rotg
->k
*wj
*sjrj
, tmpvec2
, tmpvec
);
2719 svmul(mj
, innersumvec
, tmpvec2
); /* This is != 0 only for the pivot-free variant */
2721 rvec_add(tmpvec2
, tmpvec
, erg
->f_rot_loc
[j
]);
2722 Vpart
+= wj
*psijstar
*fac2
;
2725 /* Add to the torque of this rotation group */
2726 erg
->torque_v
+= torque(rotg
->vec
, erg
->f_rot_loc
[j
], xj
, erg
->xc_center
);
2728 /* Calculate the angle between reference and actual rotation group atom. */
2729 angle(rotg
, xj_u
, rj
, &alpha
, &weight
); /* angle in rad, weighted */
2730 erg
->angle_v
+= alpha
* weight
;
2731 erg
->weight_v
+= weight
;
2736 } /* end of loop over local rotation group atoms */
2737 erg
->V
= 0.5*rotg
->k
*Vpart
;
2741 /* Determine the smallest and largest position vector (with respect to the
2742 * rotation vector) for the reference group */
2743 static void get_firstlast_atom_ref(
2748 gmx_enfrotgrp_t erg
; /* Pointer to enforced rotation group data */
2750 real xcproj
; /* The projection of a reference position on the
2752 real minproj
, maxproj
; /* Smallest and largest projection on v */
2756 erg
=rotg
->enfrotgrp
;
2758 /* Start with some value */
2759 minproj
= iprod(rotg
->x_ref
[0], rotg
->vec
);
2762 /* This is just to ensure that it still works if all the atoms of the
2763 * reference structure are situated in a plane perpendicular to the rotation
2766 *lastindex
= rotg
->nat
-1;
2768 /* Loop over all atoms of the reference group,
2769 * project them on the rotation vector to find the extremes */
2770 for (i
=0; i
<rotg
->nat
; i
++)
2772 xcproj
= iprod(rotg
->x_ref
[i
], rotg
->vec
);
2773 if (xcproj
< minproj
)
2778 if (xcproj
> maxproj
)
2787 /* Allocate memory for the slabs */
2788 static void allocate_slabs(
2795 gmx_enfrotgrp_t erg
; /* Pointer to enforced rotation group data */
2799 erg
=rotg
->enfrotgrp
;
2801 /* More slabs than are defined for the reference are never needed */
2802 nslabs
= erg
->slab_last_ref
- erg
->slab_first_ref
+ 1;
2804 /* Remember how many we allocated */
2805 erg
->nslabs_alloc
= nslabs
;
2807 if (MASTER(cr
) && bVerbose
)
2808 fprintf(fplog
, "%s allocating memory to store data for %d slabs (rotation group %d).\n",
2810 snew(erg
->slab_center
, nslabs
);
2811 snew(erg
->slab_center_ref
, nslabs
);
2812 snew(erg
->slab_weights
, nslabs
);
2813 snew(erg
->slab_torque_v
, nslabs
);
2814 snew(erg
->slab_data
, nslabs
);
2815 snew(erg
->gn_atom
, nslabs
);
2816 snew(erg
->gn_slabind
, nslabs
);
2817 snew(erg
->slab_innersumvec
, nslabs
);
2818 for (i
=0; i
<nslabs
; i
++)
2820 snew(erg
->slab_data
[i
].x
, rotg
->nat
);
2821 snew(erg
->slab_data
[i
].ref
, rotg
->nat
);
2822 snew(erg
->slab_data
[i
].weight
, rotg
->nat
);
2824 snew(erg
->xc_ref_sorted
, rotg
->nat
);
2825 snew(erg
->xc_sortind
, rotg
->nat
);
2826 snew(erg
->firstatom
, nslabs
);
2827 snew(erg
->lastatom
, nslabs
);
2831 /* From the extreme coordinates of the reference group, determine the first
2832 * and last slab of the reference. We can never have more slabs in the real
2833 * simulation than calculated here for the reference.
2835 static void get_firstlast_slab_ref(t_rotgrp
*rotg
, real mc
[], int ref_firstindex
, int ref_lastindex
)
2837 gmx_enfrotgrp_t erg
; /* Pointer to enforced rotation group data */
2838 int first
,last
,firststart
;
2842 erg
=rotg
->enfrotgrp
;
2843 first
= get_first_slab(rotg
, erg
->max_beta
, rotg
->x_ref
[ref_firstindex
]);
2844 last
= get_last_slab( rotg
, erg
->max_beta
, rotg
->x_ref
[ref_lastindex
]);
2847 while (get_slab_weight(first
, rotg
, rotg
->x_ref
, mc
, &dummy
) > WEIGHT_MIN
)
2851 erg
->slab_first_ref
= first
+1;
2852 while (get_slab_weight(last
, rotg
, rotg
->x_ref
, mc
, &dummy
) > WEIGHT_MIN
)
2856 erg
->slab_last_ref
= last
-1;
2858 erg
->slab_buffer
= firststart
- erg
->slab_first_ref
;
2863 static void init_rot_group(FILE *fplog
,t_commrec
*cr
,int g
,t_rotgrp
*rotg
,
2864 rvec
*x
,gmx_mtop_t
*mtop
,gmx_bool bVerbose
,FILE *out_slabs
)
2868 gmx_bool bFlex
,bColl
;
2870 gmx_enfrotgrp_t erg
; /* Pointer to enforced rotation group data */
2871 int ref_firstindex
, ref_lastindex
;
2872 real mass
,totalmass
;
2875 /* Do we have a flexible axis? */
2876 bFlex
= ISFLEX(rotg
);
2877 /* Do we use a global set of coordinates? */
2878 bColl
= ISCOLL(rotg
);
2880 erg
=rotg
->enfrotgrp
;
2882 /* Allocate space for collective coordinates if needed */
2885 snew(erg
->xc
, rotg
->nat
);
2886 snew(erg
->xc_shifts
, rotg
->nat
);
2887 snew(erg
->xc_eshifts
, rotg
->nat
);
2889 /* Save the original (whole) set of positions such that later the
2890 * molecule can always be made whole again */
2891 snew(erg
->xc_old
, rotg
->nat
);
2894 for (i
=0; i
<rotg
->nat
; i
++)
2897 copy_rvec(x
[ii
], erg
->xc_old
[i
]);
2902 gmx_bcast(rotg
->nat
*sizeof(erg
->xc_old
[0]),erg
->xc_old
, cr
);
2905 if (rotg
->eFittype
== erotgFitNORM
)
2907 snew(erg
->xc_ref_length
, rotg
->nat
); /* in case fit type NORM is chosen */
2908 snew(erg
->xc_norm
, rotg
->nat
);
2913 snew(erg
->xr_loc
, rotg
->nat
);
2914 snew(erg
->x_loc_pbc
, rotg
->nat
);
2917 snew(erg
->f_rot_loc
, rotg
->nat
);
2918 snew(erg
->xc_ref_ind
, rotg
->nat
);
2920 /* xc_ref_ind needs to be set to identity in the serial case */
2922 for (i
=0; i
<rotg
->nat
; i
++)
2923 erg
->xc_ref_ind
[i
] = i
;
2925 /* Copy the masses so that the COM can be determined. For all types of
2926 * enforced rotation, we store the masses in the erg->mc array. */
2927 snew(erg
->mc
, rotg
->nat
);
2929 snew(erg
->mc_sorted
, rotg
->nat
);
2931 snew(erg
->m_loc
, rotg
->nat
);
2933 for (i
=0; i
<rotg
->nat
; i
++)
2937 gmx_mtop_atomnr_to_atom(mtop
,rotg
->ind
[i
],&atom
);
2947 erg
->invmass
= 1.0/totalmass
;
2949 /* Set xc_ref_center for any rotation potential */
2950 if ((rotg
->eType
==erotgISO
) || (rotg
->eType
==erotgPM
) || (rotg
->eType
==erotgRM
) || (rotg
->eType
==erotgRM2
))
2952 /* Set the pivot point for the fixed, stationary-axis potentials. This
2953 * won't change during the simulation */
2954 copy_rvec(rotg
->pivot
, erg
->xc_ref_center
);
2955 copy_rvec(rotg
->pivot
, erg
->xc_center
);
2959 /* Center of the reference positions */
2960 get_center(rotg
->x_ref
, erg
->mc
, rotg
->nat
, erg
->xc_ref_center
);
2962 /* Center of the actual positions */
2965 snew(xdum
, rotg
->nat
);
2966 for (i
=0; i
<rotg
->nat
; i
++)
2969 copy_rvec(x
[ii
], xdum
[i
]);
2971 get_center(xdum
, erg
->mc
, rotg
->nat
, erg
->xc_center
);
2976 gmx_bcast(sizeof(erg
->xc_center
), erg
->xc_center
, cr
);
2980 if ( (rotg
->eType
!= erotgFLEX
) && (rotg
->eType
!= erotgFLEX2
) )
2982 /* Put the reference positions into origin: */
2983 for (i
=0; i
<rotg
->nat
; i
++)
2984 rvec_dec(rotg
->x_ref
[i
], erg
->xc_ref_center
);
2987 /* Enforced rotation with flexible axis */
2990 /* Calculate maximum beta value from minimum gaussian (performance opt.) */
2991 erg
->max_beta
= calc_beta_max(rotg
->min_gaussian
, rotg
->slab_dist
);
2993 /* Determine the smallest and largest coordinate with respect to the rotation vector */
2994 get_firstlast_atom_ref(rotg
, &ref_firstindex
, &ref_lastindex
);
2996 /* From the extreme coordinates of the reference group, determine the first
2997 * and last slab of the reference. */
2998 get_firstlast_slab_ref(rotg
, erg
->mc
, ref_firstindex
, ref_lastindex
);
3000 /* Allocate memory for the slabs */
3001 allocate_slabs(rotg
, fplog
, g
, bVerbose
, cr
);
3003 /* Flexible rotation: determine the reference centers for the rest of the simulation */
3004 erg
->slab_first
= erg
->slab_first_ref
;
3005 erg
->slab_last
= erg
->slab_last_ref
;
3006 get_slab_centers(rotg
,rotg
->x_ref
,erg
->mc
,cr
,g
,-1,out_slabs
,TRUE
,TRUE
);
3008 /* Length of each x_rotref vector from center (needed if fit routine NORM is chosen): */
3009 if (rotg
->eFittype
== erotgFitNORM
)
3011 for (i
=0; i
<rotg
->nat
; i
++)
3013 rvec_sub(rotg
->x_ref
[i
], erg
->xc_ref_center
, coord
);
3014 erg
->xc_ref_length
[i
] = norm(coord
);
3021 extern void dd_make_local_rotation_groups(gmx_domdec_t
*dd
,t_rot
*rot
)
3026 gmx_enfrotgrp_t erg
; /* Pointer to enforced rotation group data */
3030 for(g
=0; g
<rot
->ngrp
; g
++)
3032 rotg
= &rot
->grp
[g
];
3033 erg
= rotg
->enfrotgrp
;
3036 dd_make_local_group_indices(ga2la
,rotg
->nat
,rotg
->ind
,
3037 &erg
->nat_loc
,&erg
->ind_loc
,&erg
->nalloc_loc
,erg
->xc_ref_ind
);
3042 extern void init_rot(FILE *fplog
,t_inputrec
*ir
,int nfile
,const t_filenm fnm
[],
3043 t_commrec
*cr
, rvec
*x
, matrix box
, gmx_mtop_t
*mtop
, const output_env_t oenv
,
3044 gmx_bool bVerbose
, unsigned long Flags
)
3049 int nat_max
=0; /* Size of biggest rotation group */
3050 gmx_enfrot_t er
; /* Pointer to the enforced rotation buffer variables */
3051 gmx_enfrotgrp_t erg
; /* Pointer to enforced rotation group data */
3052 rvec
*x_pbc
=NULL
; /* Space for the pbc-correct atom positions */
3053 gmx_bool bHaveFlexGroups
= FALSE
;
3056 if ( (PAR(cr
)) && !DOMAINDECOMP(cr
) )
3057 gmx_fatal(FARGS
, "Enforced rotation is only implemented for domain decomposition!");
3059 if ( MASTER(cr
) && bVerbose
)
3060 fprintf(stdout
, "%s Initializing ...\n", RotStr
);
3064 snew(rot
->enfrot
, 1);
3067 /* Output every step for reruns */
3068 if (Flags
& MD_RERUN
)
3071 fprintf(fplog
, "%s rerun - will write rotation output every available step.\n", RotStr
);
3076 er
->out_slabs
= NULL
;
3078 er
->out_slabs
= open_slab_out(rot
, opt2fn("-rs",nfile
,fnm
));
3082 /* Remove pbc, make molecule whole.
3083 * When ir->bContinuation=TRUE this has already been done, but ok. */
3084 snew(x_pbc
,mtop
->natoms
);
3085 m_rveccopy(mtop
->natoms
,x
,x_pbc
);
3086 do_pbc_first_mtop(NULL
,ir
->ePBC
,box
,mtop
,x_pbc
);
3089 for (g
=0; g
<rot
->ngrp
; g
++)
3091 rotg
= &rot
->grp
[g
];
3093 bHaveFlexGroups
= TRUE
;
3096 fprintf(fplog
,"%s group %d type '%s'\n", RotStr
, g
, erotg_names
[rotg
->eType
]);
3100 /* Allocate space for the rotation group's data: */
3101 snew(rotg
->enfrotgrp
, 1);
3102 erg
= rotg
->enfrotgrp
;
3104 nat_max
=max(nat_max
, rotg
->nat
);
3109 erg
->nalloc_loc
= 0;
3110 erg
->ind_loc
= NULL
;
3114 erg
->nat_loc
= rotg
->nat
;
3115 erg
->ind_loc
= rotg
->ind
;
3117 init_rot_group(fplog
,cr
,g
,rotg
,x_pbc
,mtop
,bVerbose
,er
->out_slabs
);
3121 /* Allocate space for enforced rotation buffer variables */
3122 er
->bufsize
= nat_max
;
3123 snew(er
->data
, nat_max
);
3124 snew(er
->xbuf
, nat_max
);
3125 snew(er
->mbuf
, nat_max
);
3127 /* Buffers for MPI reducing torques, angles, weights (for each group), and V */
3128 er
->mpi_bufsize
= 4*rot
->ngrp
; /* To start with */
3129 snew(er
->mpi_inbuf
, er
->mpi_bufsize
);
3130 snew(er
->mpi_outbuf
, er
->mpi_bufsize
);
3132 /* Only do I/O on the MASTER */
3133 er
->out_angles
= NULL
;
3135 er
->out_torque
= NULL
;
3138 er
->out_rot
= open_rot_out(opt2fn("-ro",nfile
,fnm
), rot
, oenv
, Flags
);
3139 if (bHaveFlexGroups
)
3141 if (rot
->nstrout
> 0)
3142 er
->out_angles
= open_angles_out(rot
, opt2fn("-ra",nfile
,fnm
));
3143 if (rot
->nstsout
> 0)
3144 er
->out_torque
= open_torque_out(rot
, opt2fn("-rt",nfile
,fnm
));
3151 extern void finish_rot(FILE *fplog
,t_rot
*rot
)
3153 gmx_enfrot_t er
; /* Pointer to the enforced rotation buffer variables */
3158 gmx_fio_fclose(er
->out_rot
);
3160 gmx_fio_fclose(er
->out_slabs
);
3162 gmx_fio_fclose(er
->out_angles
);
3164 gmx_fio_fclose(er
->out_torque
);
3168 /* Rotate the local reference positions and store them in
3169 * erg->xr_loc[0...(nat_loc-1)]
3171 * Note that we already subtracted u or y_c from the reference positions
3172 * in init_rot_group().
3174 static void rotate_local_reference(t_rotgrp
*rotg
)
3176 gmx_enfrotgrp_t erg
;
3180 erg
=rotg
->enfrotgrp
;
3182 for (i
=0; i
<erg
->nat_loc
; i
++)
3184 /* Index of this rotation group atom with respect to the whole rotation group */
3185 ii
= erg
->xc_ref_ind
[i
];
3187 mvmul(erg
->rotmat
, rotg
->x_ref
[ii
], erg
->xr_loc
[i
]);
3192 /* Select the PBC representation for each local x position and store that
3193 * for later usage. We assume the right PBC image of an x is the one nearest to
3194 * its rotated reference */
3195 static void choose_pbc_image(rvec x
[], t_rotgrp
*rotg
, matrix box
, int npbcdim
)
3198 gmx_enfrotgrp_t erg
; /* Pointer to enforced rotation group data */
3203 erg
=rotg
->enfrotgrp
;
3205 for (i
=0; i
<erg
->nat_loc
; i
++)
3209 /* Index of a rotation group atom */
3210 ii
= erg
->ind_loc
[i
];
3212 /* Get the reference position. The pivot (or COM or COG) was already
3213 * subtracted in init_rot_group() from the reference positions. Also,
3214 * the reference positions have already been rotated in
3215 * rotate_local_reference() */
3216 copy_rvec(erg
->xr_loc
[i
], xref
);
3218 /* Subtract the (old) center from the current positions
3219 * (just to determine the shifts!) */
3220 rvec_sub(x
[ii
], erg
->xc_center
, xcurr
);
3222 /* Shortest PBC distance between the atom and its reference */
3223 rvec_sub(xcurr
, xref
, dx
);
3225 /* Determine the shift for this atom */
3226 for(m
=npbcdim
-1; m
>=0; m
--)
3228 while (dx
[m
] < -0.5*box
[m
][m
])
3230 for(d
=0; d
<DIM
; d
++)
3234 while (dx
[m
] >= 0.5*box
[m
][m
])
3236 for(d
=0; d
<DIM
; d
++)
3242 /* Apply the shift to the current atom */
3243 copy_rvec(x
[ii
], erg
->x_loc_pbc
[i
]);
3244 shift_single_coord(box
, erg
->x_loc_pbc
[i
], shift
);
3249 extern void do_rotation(
3255 gmx_large_int_t step
,
3256 gmx_wallcycle_t wcycle
,
3262 gmx_bool outstep_slab
, outstep_rot
;
3263 gmx_bool bFlex
,bColl
;
3265 gmx_enfrot_t er
; /* Pointer to the enforced rotation buffer variables */
3266 gmx_enfrotgrp_t erg
; /* Pointer to enforced rotation group data */
3276 /* When to output in main rotation output file */
3277 outstep_rot
= do_per_step(step
, rot
->nstrout
);
3278 /* When to output per-slab data */
3279 outstep_slab
= do_per_step(step
, rot
->nstsout
);
3281 /* Output time into rotation output file */
3282 if (outstep_rot
&& MASTER(cr
))
3283 fprintf(er
->out_rot
, "%12.3e",t
);
3285 /**************************************************************************/
3286 /* First do ALL the communication! */
3287 for(g
=0; g
<rot
->ngrp
; g
++)
3289 rotg
= &rot
->grp
[g
];
3290 erg
=rotg
->enfrotgrp
;
3292 /* Do we have a flexible axis? */
3293 bFlex
= ISFLEX(rotg
);
3294 /* Do we use a collective (global) set of coordinates? */
3295 bColl
= ISCOLL(rotg
);
3297 /* Calculate the rotation matrix for this angle: */
3298 erg
->degangle
= rotg
->rate
* t
;
3299 calc_rotmat(rotg
->vec
,erg
->degangle
,erg
->rotmat
);
3303 /* Transfer the rotation group's positions such that every node has
3304 * all of them. Every node contributes its local positions x and stores
3305 * it in the collective erg->xc array. */
3306 communicate_group_positions(cr
,erg
->xc
, erg
->xc_shifts
, erg
->xc_eshifts
, bNS
,
3307 x
, rotg
->nat
, erg
->nat_loc
, erg
->ind_loc
, erg
->xc_ref_ind
, erg
->xc_old
, box
);
3311 /* Fill the local masses array;
3312 * this array changes in DD/neighborsearching steps */
3315 for (i
=0; i
<erg
->nat_loc
; i
++)
3317 /* Index of local atom w.r.t. the collective rotation group */
3318 ii
= erg
->xc_ref_ind
[i
];
3319 erg
->m_loc
[i
] = erg
->mc
[ii
];
3323 /* Calculate Omega*(y_i-y_c) for the local positions */
3324 rotate_local_reference(rotg
);
3326 /* Choose the nearest PBC images of the group atoms with respect
3327 * to the rotated reference positions */
3328 choose_pbc_image(x
, rotg
, box
, 3);
3330 /* Get the center of the rotation group */
3331 if ( (rotg
->eType
==erotgISOPF
) || (rotg
->eType
==erotgPMPF
) )
3332 get_center_comm(cr
, erg
->x_loc_pbc
, erg
->m_loc
, erg
->nat_loc
, rotg
->nat
, erg
->xc_center
);
3335 } /* End of loop over rotation groups */
3337 /**************************************************************************/
3338 /* Done communicating, we can start to count cycles now ... */
3339 wallcycle_start(wcycle
, ewcROT
);
3340 GMX_MPE_LOG(ev_rotcycles_start
);
3346 for(g
=0; g
<rot
->ngrp
; g
++)
3348 rotg
= &rot
->grp
[g
];
3349 erg
=rotg
->enfrotgrp
;
3351 bFlex
= ISFLEX(rotg
);
3352 bColl
= ISCOLL(rotg
);
3354 if (outstep_rot
&& MASTER(cr
))
3355 fprintf(er
->out_rot
, "%12.4f", erg
->degangle
);
3357 /* Clear values from last time step */
3359 erg
->torque_v
= 0.0;
3361 erg
->weight_v
= 0.0;
3369 do_fixed(cr
,rotg
,x
,box
,t
,step
,outstep_rot
);
3372 do_radial_motion(cr
,rotg
,x
,box
,t
,step
,outstep_rot
);
3375 do_radial_motion_pf(cr
,rotg
,x
,box
,t
,step
,outstep_rot
);
3379 do_radial_motion2(cr
,rotg
,x
,box
,t
,step
,outstep_rot
);
3383 /* Subtract the center of the rotation group from the collective positions array
3384 * Also store the center in erg->xc_center since it needs to be subtracted
3385 * in the low level routines from the local coordinates as well */
3386 get_center(erg
->xc
, erg
->mc
, rotg
->nat
, erg
->xc_center
);
3387 svmul(-1.0, erg
->xc_center
, transvec
);
3388 translate_x(erg
->xc
, rotg
->nat
, transvec
);
3389 do_flexible(cr
,er
,rotg
,g
,x
,box
,t
,step
,outstep_rot
,outstep_slab
);
3393 /* Do NOT subtract the center of mass in the low level routines! */
3394 clear_rvec(erg
->xc_center
);
3395 do_flexible(cr
,er
,rotg
,g
,x
,box
,t
,step
,outstep_rot
,outstep_slab
);
3398 gmx_fatal(FARGS
, "No such rotation potential.");
3405 fprintf(stderr
, "%s calculation (step %d) took %g seconds.\n", RotStr
, step
, MPI_Wtime()-t0
);
3408 /* Stop the cycle counter and add to the force cycles for load balancing */
3409 cycles_rot
= wallcycle_stop(wcycle
,ewcROT
);
3410 if (DOMAINDECOMP(cr
) && wcycle
)
3411 dd_cycles_add(cr
->dd
,cycles_rot
,ddCyclF
);
3412 GMX_MPE_LOG(ev_rotcycles_finish
);