Merge branch 'release-4-6'
[alexxy/gromacs.git] / src / gromacs / pulling / pull_rotation.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2008, The GROMACS development team.
6  * Copyright (c) 2013,2014, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
40
41 #include <stdio.h>
42 #include <stdlib.h>
43 #include <string.h>
44 #include "domdec.h"
45 #include "smalloc.h"
46 #include "network.h"
47 #include "pbc.h"
48 #include "mdrun.h"
49 #include "txtdump.h"
50 #include "names.h"
51 #include "mtop_util.h"
52 #include "names.h"
53 #include "nrjac.h"
54 #include "vec.h"
55 #include "gmx_ga2la.h"
56 #include "xvgr.h"
57 #include "copyrite.h"
58 #include "macros.h"
59
60 #include "gromacs/fileio/futil.h"
61 #include "gromacs/fileio/gmxfio.h"
62 #include "gromacs/fileio/trnio.h"
63 #include "gromacs/timing/cyclecounter.h"
64 #include "gromacs/timing/wallcycle.h"
65 #include "gromacs/utility/qsort_threadsafe.h"
66 #include "gromacs/pulling/pull_rotation.h"
67 #include "gromacs/mdlib/groupcoord.h"
68
69 static char *RotStr = {"Enforced rotation:"};
70
71 /* Set the minimum weight for the determination of the slab centers */
72 #define WEIGHT_MIN (10*GMX_FLOAT_MIN)
73
74 /* Helper structure for sorting positions along rotation vector             */
75 typedef struct {
76     real xcproj;            /* Projection of xc on the rotation vector        */
77     int  ind;               /* Index of xc                                    */
78     real m;                 /* Mass                                           */
79     rvec x;                 /* Position                                       */
80     rvec x_ref;             /* Reference position                             */
81 } sort_along_vec_t;
82
83
84 /* Enforced rotation / flexible: determine the angle of each slab             */
85 typedef struct gmx_slabdata
86 {
87     int   nat;              /* Number of atoms belonging to this slab         */
88     rvec *x;                /* The positions belonging to this slab. In
89                                general, this should be all positions of the
90                                whole rotation group, but we leave those away
91                                that have a small enough weight                */
92     rvec *ref;              /* Same for reference                             */
93     real *weight;           /* The weight for each atom                       */
94 } t_gmx_slabdata;
95
96
97 /* Helper structure for potential fitting */
98 typedef struct gmx_potfit
99 {
100     real   *degangle;       /* Set of angles for which the potential is
101                                calculated. The optimum fit is determined as
102                                the angle for with the potential is minimal    */
103     real   *V;              /* Potential for the different angles             */
104     matrix *rotmat;         /* Rotation matrix corresponding to the angles    */
105 } t_gmx_potfit;
106
107
108 /* Enforced rotation data for all groups                                      */
109 typedef struct gmx_enfrot
110 {
111     FILE             *out_rot;     /* Output file for rotation data                  */
112     FILE             *out_torque;  /* Output file for torque data                    */
113     FILE             *out_angles;  /* Output file for slab angles for flexible type  */
114     FILE             *out_slabs;   /* Output file for slab centers                   */
115     int               bufsize;     /* Allocation size of buf                         */
116     rvec             *xbuf;        /* Coordinate buffer variable for sorting         */
117     real             *mbuf;        /* Masses buffer variable for sorting             */
118     sort_along_vec_t *data;        /* Buffer variable needed for position sorting    */
119     real             *mpi_inbuf;   /* MPI buffer                                     */
120     real             *mpi_outbuf;  /* MPI buffer                                     */
121     int               mpi_bufsize; /* Allocation size of in & outbuf                 */
122     unsigned long     Flags;       /* mdrun flags                                    */
123     gmx_bool          bOut;        /* Used to skip first output when appending to
124                                     * avoid duplicate entries in rotation outfiles   */
125 } t_gmx_enfrot;
126
127
128 /* Global enforced rotation data for a single rotation group                  */
129 typedef struct gmx_enfrotgrp
130 {
131     real     degangle;      /* Rotation angle in degrees                      */
132     matrix   rotmat;        /* Rotation matrix                                */
133     atom_id *ind_loc;       /* Local rotation indices                         */
134     int      nat_loc;       /* Number of local group atoms                    */
135     int      nalloc_loc;    /* Allocation size for ind_loc and weight_loc     */
136
137     real     V;             /* Rotation potential for this rotation group     */
138     rvec    *f_rot_loc;     /* Array to store the forces on the local atoms
139                                resulting from enforced rotation potential     */
140
141     /* Collective coordinates for the whole rotation group */
142     real  *xc_ref_length;   /* Length of each x_rotref vector after x_rotref
143                                has been put into origin                       */
144     int   *xc_ref_ind;      /* Position of each local atom in the collective
145                                array                                          */
146     rvec   xc_center;       /* Center of the rotation group positions, may
147                                be mass weighted                               */
148     rvec   xc_ref_center;   /* dito, for the reference positions              */
149     rvec  *xc;              /* Current (collective) positions                 */
150     ivec  *xc_shifts;       /* Current (collective) shifts                    */
151     ivec  *xc_eshifts;      /* Extra shifts since last DD step                */
152     rvec  *xc_old;          /* Old (collective) positions                     */
153     rvec  *xc_norm;         /* Normalized form of the current positions       */
154     rvec  *xc_ref_sorted;   /* Reference positions (sorted in the same order
155                                as xc when sorted)                             */
156     int   *xc_sortind;      /* Where is a position found after sorting?       */
157     real  *mc;              /* Collective masses                              */
158     real  *mc_sorted;
159     real   invmass;         /* one over the total mass of the rotation group  */
160
161     real   torque_v;        /* Torque in the direction of rotation vector     */
162     real   angle_v;         /* Actual angle of the whole rotation group       */
163     /* Fixed rotation only */
164     real   weight_v;        /* Weights for angle determination                */
165     rvec  *xr_loc;          /* Local reference coords, correctly rotated      */
166     rvec  *x_loc_pbc;       /* Local current coords, correct PBC image        */
167     real  *m_loc;           /* Masses of the current local atoms              */
168
169     /* Flexible rotation only */
170     int    nslabs_alloc;              /* For this many slabs memory is allocated        */
171     int    slab_first;                /* Lowermost slab for that the calculation needs
172                                          to be performed at a given time step           */
173     int    slab_last;                 /* Uppermost slab ...                             */
174     int    slab_first_ref;            /* First slab for which ref. center is stored     */
175     int    slab_last_ref;             /* Last ...                                       */
176     int    slab_buffer;               /* Slab buffer region around reference slabs      */
177     int   *firstatom;                 /* First relevant atom for a slab                 */
178     int   *lastatom;                  /* Last relevant atom for a slab                  */
179     rvec  *slab_center;               /* Gaussian-weighted slab center                  */
180     rvec  *slab_center_ref;           /* Gaussian-weighted slab center for the
181                                          reference positions                            */
182     real  *slab_weights;              /* Sum of gaussian weights in a slab              */
183     real  *slab_torque_v;             /* Torque T = r x f for each slab.                */
184                                       /* torque_v = m.v = angular momentum in the
185                                          direction of v                                 */
186     real  max_beta;                   /* min_gaussian from inputrec->rotgrp is the
187                                          minimum value the gaussian must have so that
188                                          the force is actually evaluated max_beta is
189                                          just another way to put it                     */
190     real           *gn_atom;          /* Precalculated gaussians for a single atom      */
191     int            *gn_slabind;       /* Tells to which slab each precalculated gaussian
192                                          belongs                                        */
193     rvec           *slab_innersumvec; /* Inner sum of the flexible2 potential per slab;
194                                          this is precalculated for optimization reasons */
195     t_gmx_slabdata *slab_data;        /* Holds atom positions and gaussian weights
196                                          of atoms belonging to a slab                   */
197
198     /* For potential fits with varying angle: */
199     t_gmx_potfit *PotAngleFit;  /* Used for fit type 'potential'              */
200 } t_gmx_enfrotgrp;
201
202
203 /* Activate output of forces for correctness checks */
204 /* #define PRINT_FORCES */
205 #ifdef PRINT_FORCES
206 #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]);
207 #define PRINT_POT_TAU  if (MASTER(cr)) { \
208         fprintf(stderr, "potential = %15.8f\n" "torque    = %15.8f\n", erg->V, erg->torque_v); \
209 }
210 #else
211 #define PRINT_FORCE_J
212 #define PRINT_POT_TAU
213 #endif
214
215 /* Shortcuts for often used queries */
216 #define ISFLEX(rg) ( (rg->eType == erotgFLEX) || (rg->eType == erotgFLEXT) || (rg->eType == erotgFLEX2) || (rg->eType == erotgFLEX2T) )
217 #define ISCOLL(rg) ( (rg->eType == erotgFLEX) || (rg->eType == erotgFLEXT) || (rg->eType == erotgFLEX2) || (rg->eType == erotgFLEX2T) || (rg->eType == erotgRMPF) || (rg->eType == erotgRM2PF) )
218
219
220 /* Does any of the rotation groups use slab decomposition? */
221 static gmx_bool HaveFlexibleGroups(t_rot *rot)
222 {
223     int       g;
224     t_rotgrp *rotg;
225
226
227     for (g = 0; g < rot->ngrp; g++)
228     {
229         rotg = &rot->grp[g];
230         if (ISFLEX(rotg))
231         {
232             return TRUE;
233         }
234     }
235
236     return FALSE;
237 }
238
239
240 /* Is for any group the fit angle determined by finding the minimum of the
241  * rotation potential? */
242 static gmx_bool HavePotFitGroups(t_rot *rot)
243 {
244     int       g;
245     t_rotgrp *rotg;
246
247
248     for (g = 0; g < rot->ngrp; g++)
249     {
250         rotg = &rot->grp[g];
251         if (erotgFitPOT == rotg->eFittype)
252         {
253             return TRUE;
254         }
255     }
256
257     return FALSE;
258 }
259
260
261 static double** allocate_square_matrix(int dim)
262 {
263     int      i;
264     double** mat = NULL;
265
266
267     snew(mat, dim);
268     for (i = 0; i < dim; i++)
269     {
270         snew(mat[i], dim);
271     }
272
273     return mat;
274 }
275
276
277 static void free_square_matrix(double** mat, int dim)
278 {
279     int i;
280
281
282     for (i = 0; i < dim; i++)
283     {
284         sfree(mat[i]);
285     }
286     sfree(mat);
287 }
288
289
290 /* Return the angle for which the potential is minimal */
291 static real get_fitangle(t_rotgrp *rotg, gmx_enfrotgrp_t erg)
292 {
293     int           i;
294     real          fitangle = -999.9;
295     real          pot_min  = GMX_FLOAT_MAX;
296     t_gmx_potfit *fit;
297
298
299     fit = erg->PotAngleFit;
300
301     for (i = 0; i < rotg->PotAngle_nstep; i++)
302     {
303         if (fit->V[i] < pot_min)
304         {
305             pot_min  = fit->V[i];
306             fitangle = fit->degangle[i];
307         }
308     }
309
310     return fitangle;
311 }
312
313
314 /* Reduce potential angle fit data for this group at this time step? */
315 static gmx_inline gmx_bool bPotAngle(t_rot *rot, t_rotgrp *rotg, gmx_int64_t step)
316 {
317     return ( (erotgFitPOT == rotg->eFittype) && (do_per_step(step, rot->nstsout) || do_per_step(step, rot->nstrout)) );
318 }
319
320 /* Reduce slab torqe data for this group at this time step? */
321 static gmx_inline gmx_bool bSlabTau(t_rot *rot, t_rotgrp *rotg, gmx_int64_t step)
322 {
323     return ( (ISFLEX(rotg)) && do_per_step(step, rot->nstsout) );
324 }
325
326 /* Output rotation energy, torques, etc. for each rotation group */
327 static void reduce_output(t_commrec *cr, t_rot *rot, real t, gmx_int64_t step)
328 {
329     int             g, i, islab, nslabs = 0;
330     int             count; /* MPI element counter                               */
331     t_rotgrp       *rotg;
332     gmx_enfrot_t    er;    /* Pointer to the enforced rotation buffer variables */
333     gmx_enfrotgrp_t erg;   /* Pointer to enforced rotation group data           */
334     real            fitangle;
335     gmx_bool        bFlex;
336
337
338     er = rot->enfrot;
339
340     /* Fill the MPI buffer with stuff to reduce. If items are added for reduction
341      * here, the MPI buffer size has to be enlarged also in calc_mpi_bufsize() */
342     if (PAR(cr))
343     {
344         count = 0;
345         for (g = 0; g < rot->ngrp; g++)
346         {
347             rotg                   = &rot->grp[g];
348             erg                    = rotg->enfrotgrp;
349             nslabs                 = erg->slab_last - erg->slab_first + 1;
350             er->mpi_inbuf[count++] = erg->V;
351             er->mpi_inbuf[count++] = erg->torque_v;
352             er->mpi_inbuf[count++] = erg->angle_v;
353             er->mpi_inbuf[count++] = erg->weight_v; /* weights are not needed for flex types, but this is just a single value */
354
355             if (bPotAngle(rot, rotg, step))
356             {
357                 for (i = 0; i < rotg->PotAngle_nstep; i++)
358                 {
359                     er->mpi_inbuf[count++] = erg->PotAngleFit->V[i];
360                 }
361             }
362             if (bSlabTau(rot, rotg, step))
363             {
364                 for (i = 0; i < nslabs; i++)
365                 {
366                     er->mpi_inbuf[count++] = erg->slab_torque_v[i];
367                 }
368             }
369         }
370         if (count > er->mpi_bufsize)
371         {
372             gmx_fatal(FARGS, "%s MPI buffer overflow, please report this error.", RotStr);
373         }
374
375 #ifdef GMX_MPI
376         MPI_Reduce(er->mpi_inbuf, er->mpi_outbuf, count, GMX_MPI_REAL, MPI_SUM, MASTERRANK(cr), cr->mpi_comm_mygroup);
377 #endif
378
379         /* Copy back the reduced data from the buffer on the master */
380         if (MASTER(cr))
381         {
382             count = 0;
383             for (g = 0; g < rot->ngrp; g++)
384             {
385                 rotg          = &rot->grp[g];
386                 erg           = rotg->enfrotgrp;
387                 nslabs        = erg->slab_last - erg->slab_first + 1;
388                 erg->V        = er->mpi_outbuf[count++];
389                 erg->torque_v = er->mpi_outbuf[count++];
390                 erg->angle_v  = er->mpi_outbuf[count++];
391                 erg->weight_v = er->mpi_outbuf[count++];
392
393                 if (bPotAngle(rot, rotg, step))
394                 {
395                     for (i = 0; i < rotg->PotAngle_nstep; i++)
396                     {
397                         erg->PotAngleFit->V[i] = er->mpi_outbuf[count++];
398                     }
399                 }
400                 if (bSlabTau(rot, rotg, step))
401                 {
402                     for (i = 0; i < nslabs; i++)
403                     {
404                         erg->slab_torque_v[i] = er->mpi_outbuf[count++];
405                     }
406                 }
407             }
408         }
409     }
410
411     /* Output */
412     if (MASTER(cr))
413     {
414         /* Angle and torque for each rotation group */
415         for (g = 0; g < rot->ngrp; g++)
416         {
417             rotg  = &rot->grp[g];
418             bFlex = ISFLEX(rotg);
419
420             erg = rotg->enfrotgrp;
421
422             /* Output to main rotation output file: */
423             if (do_per_step(step, rot->nstrout) )
424             {
425                 if (erotgFitPOT == rotg->eFittype)
426                 {
427                     fitangle = get_fitangle(rotg, erg);
428                 }
429                 else
430                 {
431                     if (bFlex)
432                     {
433                         fitangle = erg->angle_v; /* RMSD fit angle */
434                     }
435                     else
436                     {
437                         fitangle = (erg->angle_v/erg->weight_v)*180.0*M_1_PI;
438                     }
439                 }
440                 fprintf(er->out_rot, "%12.4f", fitangle);
441                 fprintf(er->out_rot, "%12.3e", erg->torque_v);
442                 fprintf(er->out_rot, "%12.3e", erg->V);
443             }
444
445             if (do_per_step(step, rot->nstsout) )
446             {
447                 /* Output to torque log file: */
448                 if (bFlex)
449                 {
450                     fprintf(er->out_torque, "%12.3e%6d", t, g);
451                     for (i = erg->slab_first; i <= erg->slab_last; i++)
452                     {
453                         islab = i - erg->slab_first;  /* slab index */
454                         /* Only output if enough weight is in slab */
455                         if (erg->slab_weights[islab] > rotg->min_gaussian)
456                         {
457                             fprintf(er->out_torque, "%6d%12.3e", i, erg->slab_torque_v[islab]);
458                         }
459                     }
460                     fprintf(er->out_torque, "\n");
461                 }
462
463                 /* Output to angles log file: */
464                 if (erotgFitPOT == rotg->eFittype)
465                 {
466                     fprintf(er->out_angles, "%12.3e%6d%12.4f", t, g, erg->degangle);
467                     /* Output energies at a set of angles around the reference angle */
468                     for (i = 0; i < rotg->PotAngle_nstep; i++)
469                     {
470                         fprintf(er->out_angles, "%12.3e", erg->PotAngleFit->V[i]);
471                     }
472                     fprintf(er->out_angles, "\n");
473                 }
474             }
475         }
476         if (do_per_step(step, rot->nstrout) )
477         {
478             fprintf(er->out_rot, "\n");
479         }
480     }
481 }
482
483
484 /* Add the forces from enforced rotation potential to the local forces.
485  * Should be called after the SR forces have been evaluated */
486 extern real add_rot_forces(t_rot *rot, rvec f[], t_commrec *cr, gmx_int64_t step, real t)
487 {
488     int             g, l, ii;
489     t_rotgrp       *rotg;
490     gmx_enfrot_t    er;         /* Pointer to the enforced rotation buffer variables */
491     gmx_enfrotgrp_t erg;        /* Pointer to enforced rotation group data           */
492     real            Vrot = 0.0; /* If more than one rotation group is present, Vrot
493                                    assembles the local parts from all groups         */
494
495
496     er = rot->enfrot;
497
498     /* Loop over enforced rotation groups (usually 1, though)
499      * Apply the forces from rotation potentials */
500     for (g = 0; g < rot->ngrp; g++)
501     {
502         rotg  = &rot->grp[g];
503         erg   = rotg->enfrotgrp;
504         Vrot += erg->V;  /* add the local parts from the nodes */
505         for (l = 0; l < erg->nat_loc; l++)
506         {
507             /* Get the right index of the local force */
508             ii = erg->ind_loc[l];
509             /* Add */
510             rvec_inc(f[ii], erg->f_rot_loc[l]);
511         }
512     }
513
514     /* Reduce energy,torque, angles etc. to get the sum values (per rotation group)
515      * on the master and output these values to file. */
516     if ( (do_per_step(step, rot->nstrout) || do_per_step(step, rot->nstsout)) && er->bOut)
517     {
518         reduce_output(cr, rot, t, step);
519     }
520
521     /* When appending, er->bOut is FALSE the first time to avoid duplicate entries */
522     er->bOut = TRUE;
523
524     PRINT_POT_TAU
525
526     return Vrot;
527 }
528
529
530 /* The Gaussian norm is chosen such that the sum of the gaussian functions
531  * over the slabs is approximately 1.0 everywhere */
532 #define GAUSS_NORM   0.569917543430618
533
534
535 /* Calculate the maximum beta that leads to a gaussian larger min_gaussian,
536  * also does some checks
537  */
538 static double calc_beta_max(real min_gaussian, real slab_dist)
539 {
540     double sigma;
541     double arg;
542
543
544     /* Actually the next two checks are already made in grompp */
545     if (slab_dist <= 0)
546     {
547         gmx_fatal(FARGS, "Slab distance of flexible rotation groups must be >=0 !");
548     }
549     if (min_gaussian <= 0)
550     {
551         gmx_fatal(FARGS, "Cutoff value for Gaussian must be > 0. (You requested %f)");
552     }
553
554     /* Define the sigma value */
555     sigma = 0.7*slab_dist;
556
557     /* Calculate the argument for the logarithm and check that the log() result is negative or 0 */
558     arg = min_gaussian/GAUSS_NORM;
559     if (arg > 1.0)
560     {
561         gmx_fatal(FARGS, "min_gaussian of flexible rotation groups must be <%g", GAUSS_NORM);
562     }
563
564     return sqrt(-2.0*sigma*sigma*log(min_gaussian/GAUSS_NORM));
565 }
566
567
568 static gmx_inline real calc_beta(rvec curr_x, t_rotgrp *rotg, int n)
569 {
570     return iprod(curr_x, rotg->vec) - rotg->slab_dist * n;
571 }
572
573
574 static gmx_inline real gaussian_weight(rvec curr_x, t_rotgrp *rotg, int n)
575 {
576     const real norm = GAUSS_NORM;
577     real       sigma;
578
579
580     /* Define the sigma value */
581     sigma = 0.7*rotg->slab_dist;
582     /* Calculate the Gaussian value of slab n for position curr_x */
583     return norm * exp( -0.5 * sqr( calc_beta(curr_x, rotg, n)/sigma ) );
584 }
585
586
587 /* Returns the weight in a single slab, also calculates the Gaussian- and mass-
588  * weighted sum of positions for that slab */
589 static real get_slab_weight(int j, t_rotgrp *rotg, rvec xc[], real mc[], rvec *x_weighted_sum)
590 {
591     rvec            curr_x;           /* The position of an atom                      */
592     rvec            curr_x_weighted;  /* The gaussian-weighted position               */
593     real            gaussian;         /* A single gaussian weight                     */
594     real            wgauss;           /* gaussian times current mass                  */
595     real            slabweight = 0.0; /* The sum of weights in the slab               */
596     int             i, islab;
597     gmx_enfrotgrp_t erg;              /* Pointer to enforced rotation group data      */
598
599
600     erg = rotg->enfrotgrp;
601     clear_rvec(*x_weighted_sum);
602
603     /* Slab index */
604     islab = j - erg->slab_first;
605
606     /* Loop over all atoms in the rotation group */
607     for (i = 0; i < rotg->nat; i++)
608     {
609         copy_rvec(xc[i], curr_x);
610         gaussian = gaussian_weight(curr_x, rotg, j);
611         wgauss   = gaussian * mc[i];
612         svmul(wgauss, curr_x, curr_x_weighted);
613         rvec_add(*x_weighted_sum, curr_x_weighted, *x_weighted_sum);
614         slabweight += wgauss;
615     }  /* END of loop over rotation group atoms */
616
617     return slabweight;
618 }
619
620
621 static void get_slab_centers(
622         t_rotgrp  *rotg,       /* The rotation group information               */
623         rvec      *xc,         /* The rotation group positions; will
624                                   typically be enfrotgrp->xc, but at first call
625                                   it is enfrotgrp->xc_ref                      */
626         real      *mc,         /* The masses of the rotation group atoms       */
627         int        g,          /* The number of the rotation group             */
628         real       time,       /* Used for output only                         */
629         FILE      *out_slabs,  /* For outputting center per slab information   */
630         gmx_bool   bOutStep,   /* Is this an output step?                      */
631         gmx_bool   bReference) /* If this routine is called from
632                                   init_rot_group we need to store
633                                   the reference slab centers                   */
634 {
635     int             j, islab;
636     gmx_enfrotgrp_t erg;      /* Pointer to enforced rotation group data */
637
638
639     erg = rotg->enfrotgrp;
640
641     /* Loop over slabs */
642     for (j = erg->slab_first; j <= erg->slab_last; j++)
643     {
644         islab                    = j - erg->slab_first;
645         erg->slab_weights[islab] = get_slab_weight(j, rotg, xc, mc, &erg->slab_center[islab]);
646
647         /* We can do the calculations ONLY if there is weight in the slab! */
648         if (erg->slab_weights[islab] > WEIGHT_MIN)
649         {
650             svmul(1.0/erg->slab_weights[islab], erg->slab_center[islab], erg->slab_center[islab]);
651         }
652         else
653         {
654             /* We need to check this here, since we divide through slab_weights
655              * in the flexible low-level routines! */
656             gmx_fatal(FARGS, "Not enough weight in slab %d. Slab center cannot be determined!", j);
657         }
658
659         /* At first time step: save the centers of the reference structure */
660         if (bReference)
661         {
662             copy_rvec(erg->slab_center[islab], erg->slab_center_ref[islab]);
663         }
664     } /* END of loop over slabs */
665
666     /* Output on the master */
667     if ( (NULL != out_slabs) && bOutStep)
668     {
669         fprintf(out_slabs, "%12.3e%6d", time, g);
670         for (j = erg->slab_first; j <= erg->slab_last; j++)
671         {
672             islab = j - erg->slab_first;
673             fprintf(out_slabs, "%6d%12.3e%12.3e%12.3e",
674                     j, erg->slab_center[islab][XX], erg->slab_center[islab][YY], erg->slab_center[islab][ZZ]);
675         }
676         fprintf(out_slabs, "\n");
677     }
678 }
679
680
681 static void calc_rotmat(
682         rvec   vec,
683         real   degangle,      /* Angle alpha of rotation at time t in degrees       */
684         matrix rotmat)        /* Rotation matrix                                    */
685 {
686     real radangle;            /* Rotation angle in radians */
687     real cosa;                /* cosine alpha              */
688     real sina;                /* sine alpha                */
689     real OMcosa;              /* 1 - cos(alpha)            */
690     real dumxy, dumxz, dumyz; /* save computations         */
691     rvec rot_vec;             /* Rotate around rot_vec ... */
692
693
694     radangle = degangle * M_PI/180.0;
695     copy_rvec(vec, rot_vec );
696
697     /* Precompute some variables: */
698     cosa   = cos(radangle);
699     sina   = sin(radangle);
700     OMcosa = 1.0 - cosa;
701     dumxy  = rot_vec[XX]*rot_vec[YY]*OMcosa;
702     dumxz  = rot_vec[XX]*rot_vec[ZZ]*OMcosa;
703     dumyz  = rot_vec[YY]*rot_vec[ZZ]*OMcosa;
704
705     /* Construct the rotation matrix for this rotation group: */
706     /* 1st column: */
707     rotmat[XX][XX] = cosa  + rot_vec[XX]*rot_vec[XX]*OMcosa;
708     rotmat[YY][XX] = dumxy + rot_vec[ZZ]*sina;
709     rotmat[ZZ][XX] = dumxz - rot_vec[YY]*sina;
710     /* 2nd column: */
711     rotmat[XX][YY] = dumxy - rot_vec[ZZ]*sina;
712     rotmat[YY][YY] = cosa  + rot_vec[YY]*rot_vec[YY]*OMcosa;
713     rotmat[ZZ][YY] = dumyz + rot_vec[XX]*sina;
714     /* 3rd column: */
715     rotmat[XX][ZZ] = dumxz + rot_vec[YY]*sina;
716     rotmat[YY][ZZ] = dumyz - rot_vec[XX]*sina;
717     rotmat[ZZ][ZZ] = cosa  + rot_vec[ZZ]*rot_vec[ZZ]*OMcosa;
718
719 #ifdef PRINTMATRIX
720     int iii, jjj;
721
722     for (iii = 0; iii < 3; iii++)
723     {
724         for (jjj = 0; jjj < 3; jjj++)
725         {
726             fprintf(stderr, " %10.8f ",  rotmat[iii][jjj]);
727         }
728         fprintf(stderr, "\n");
729     }
730 #endif
731 }
732
733
734 /* Calculates torque on the rotation axis tau = position x force */
735 static gmx_inline real torque(
736         rvec rotvec,  /* rotation vector; MUST be normalized!                 */
737         rvec force,   /* force                                                */
738         rvec x,       /* position of atom on which the force acts             */
739         rvec pivot)   /* pivot point of rotation axis                         */
740 {
741     rvec vectmp, tau;
742
743
744     /* Subtract offset */
745     rvec_sub(x, pivot, vectmp);
746
747     /* position x force */
748     cprod(vectmp, force, tau);
749
750     /* Return the part of the torque which is parallel to the rotation vector */
751     return iprod(tau, rotvec);
752 }
753
754
755 /* Right-aligned output of value with standard width */
756 static void print_aligned(FILE *fp, char *str)
757 {
758     fprintf(fp, "%12s", str);
759 }
760
761
762 /* Right-aligned output of value with standard short width */
763 static void print_aligned_short(FILE *fp, char *str)
764 {
765     fprintf(fp, "%6s", str);
766 }
767
768
769 static FILE *open_output_file(const char *fn, int steps, const char what[])
770 {
771     FILE *fp;
772
773
774     fp = ffopen(fn, "w");
775
776     fprintf(fp, "# Output of %s is written in intervals of %d time step%s.\n#\n",
777             what, steps, steps > 1 ? "s" : "");
778
779     return fp;
780 }
781
782
783 /* Open output file for slab center data. Call on master only */
784 static FILE *open_slab_out(const char *fn, t_rot *rot)
785 {
786     FILE      *fp;
787     int        g, i;
788     t_rotgrp  *rotg;
789
790
791     if (rot->enfrot->Flags & MD_APPENDFILES)
792     {
793         fp = gmx_fio_fopen(fn, "a");
794     }
795     else
796     {
797         fp = open_output_file(fn, rot->nstsout, "gaussian weighted slab centers");
798
799         for (g = 0; g < rot->ngrp; g++)
800         {
801             rotg = &rot->grp[g];
802             if (ISFLEX(rotg))
803             {
804                 fprintf(fp, "# Rotation group %d (%s), slab distance %f nm, %s.\n",
805                         g, erotg_names[rotg->eType], rotg->slab_dist,
806                         rotg->bMassW ? "centers of mass" : "geometrical centers");
807             }
808         }
809
810         fprintf(fp, "# Reference centers are listed first (t=-1).\n");
811         fprintf(fp, "# The following columns have the syntax:\n");
812         fprintf(fp, "#     ");
813         print_aligned_short(fp, "t");
814         print_aligned_short(fp, "grp");
815         /* Print legend for the first two entries only ... */
816         for (i = 0; i < 2; i++)
817         {
818             print_aligned_short(fp, "slab");
819             print_aligned(fp, "X center");
820             print_aligned(fp, "Y center");
821             print_aligned(fp, "Z center");
822         }
823         fprintf(fp, " ...\n");
824         fflush(fp);
825     }
826
827     return fp;
828 }
829
830
831 /* Adds 'buf' to 'str' */
832 static void add_to_string(char **str, char *buf)
833 {
834     int len;
835
836
837     len = strlen(*str) + strlen(buf) + 1;
838     srenew(*str, len);
839     strcat(*str, buf);
840 }
841
842
843 static void add_to_string_aligned(char **str, char *buf)
844 {
845     char buf_aligned[STRLEN];
846
847     sprintf(buf_aligned, "%12s", buf);
848     add_to_string(str, buf_aligned);
849 }
850
851
852 /* Open output file and print some general information about the rotation groups.
853  * Call on master only */
854 static FILE *open_rot_out(const char *fn, t_rot *rot, const output_env_t oenv)
855 {
856     FILE           *fp;
857     int             g, nsets;
858     t_rotgrp       *rotg;
859     const char    **setname;
860     char            buf[50], buf2[75];
861     gmx_enfrotgrp_t erg;       /* Pointer to enforced rotation group data */
862     gmx_bool        bFlex;
863     char           *LegendStr = NULL;
864
865
866     if (rot->enfrot->Flags & MD_APPENDFILES)
867     {
868         fp = gmx_fio_fopen(fn, "a");
869     }
870     else
871     {
872         fp = xvgropen(fn, "Rotation angles and energy", "Time (ps)", "angles (degrees) and energies (kJ/mol)", oenv);
873         fprintf(fp, "# Output of enforced rotation data is written in intervals of %d time step%s.\n#\n", rot->nstrout, rot->nstrout > 1 ? "s" : "");
874         fprintf(fp, "# The scalar tau is the torque (kJ/mol) in the direction of the rotation vector v.\n");
875         fprintf(fp, "# To obtain the vectorial torque, multiply tau with the group's rot_vec.\n");
876         fprintf(fp, "# For flexible groups, tau(t,n) from all slabs n have been summed in a single value tau(t) here.\n");
877         fprintf(fp, "# The torques tau(t,n) are found in the rottorque.log (-rt) output file\n");
878
879         for (g = 0; g < rot->ngrp; g++)
880         {
881             rotg  = &rot->grp[g];
882             erg   = rotg->enfrotgrp;
883             bFlex = ISFLEX(rotg);
884
885             fprintf(fp, "#\n");
886             fprintf(fp, "# ROTATION GROUP %d, potential type '%s':\n", g, erotg_names[rotg->eType]);
887             fprintf(fp, "# rot_massw%d          %s\n", g, yesno_names[rotg->bMassW]);
888             fprintf(fp, "# rot_vec%d            %12.5e %12.5e %12.5e\n", g, rotg->vec[XX], rotg->vec[YY], rotg->vec[ZZ]);
889             fprintf(fp, "# rot_rate%d           %12.5e degrees/ps\n", g, rotg->rate);
890             fprintf(fp, "# rot_k%d              %12.5e kJ/(mol*nm^2)\n", g, rotg->k);
891             if (rotg->eType == erotgISO || rotg->eType == erotgPM || rotg->eType == erotgRM || rotg->eType == erotgRM2)
892             {
893                 fprintf(fp, "# rot_pivot%d          %12.5e %12.5e %12.5e  nm\n", g, rotg->pivot[XX], rotg->pivot[YY], rotg->pivot[ZZ]);
894             }
895
896             if (bFlex)
897             {
898                 fprintf(fp, "# rot_slab_distance%d   %f nm\n", g, rotg->slab_dist);
899                 fprintf(fp, "# rot_min_gaussian%d   %12.5e\n", g, rotg->min_gaussian);
900             }
901
902             /* Output the centers of the rotation groups for the pivot-free potentials */
903             if ((rotg->eType == erotgISOPF) || (rotg->eType == erotgPMPF) || (rotg->eType == erotgRMPF) || (rotg->eType == erotgRM2PF
904                                                                                                             || (rotg->eType == erotgFLEXT) || (rotg->eType == erotgFLEX2T)) )
905             {
906                 fprintf(fp, "# ref. grp. %d center  %12.5e %12.5e %12.5e\n", g,
907                         erg->xc_ref_center[XX], erg->xc_ref_center[YY], erg->xc_ref_center[ZZ]);
908
909                 fprintf(fp, "# grp. %d init.center  %12.5e %12.5e %12.5e\n", g,
910                         erg->xc_center[XX], erg->xc_center[YY], erg->xc_center[ZZ]);
911             }
912
913             if ( (rotg->eType == erotgRM2) || (rotg->eType == erotgFLEX2) || (rotg->eType == erotgFLEX2T) )
914             {
915                 fprintf(fp, "# rot_eps%d            %12.5e nm^2\n", g, rotg->eps);
916             }
917             if (erotgFitPOT == rotg->eFittype)
918             {
919                 fprintf(fp, "#\n");
920                 fprintf(fp, "# theta_fit%d is determined by first evaluating the potential for %d angles around theta_ref%d.\n",
921                         g, rotg->PotAngle_nstep, g);
922                 fprintf(fp, "# The fit angle is the one with the smallest potential. It is given as the deviation\n");
923                 fprintf(fp, "# from the reference angle, i.e. if theta_ref=X and theta_fit=Y, then the angle with\n");
924                 fprintf(fp, "# minimal value of the potential is X+Y. Angular resolution is %g degrees.\n", rotg->PotAngle_step);
925             }
926         }
927
928         /* Print a nice legend */
929         snew(LegendStr, 1);
930         LegendStr[0] = '\0';
931         sprintf(buf, "#     %6s", "time");
932         add_to_string_aligned(&LegendStr, buf);
933
934         nsets = 0;
935         snew(setname, 4*rot->ngrp);
936
937         for (g = 0; g < rot->ngrp; g++)
938         {
939             rotg = &rot->grp[g];
940             sprintf(buf, "theta_ref%d", g);
941             add_to_string_aligned(&LegendStr, buf);
942
943             sprintf(buf2, "%s (degrees)", buf);
944             setname[nsets] = strdup(buf2);
945             nsets++;
946         }
947         for (g = 0; g < rot->ngrp; g++)
948         {
949             rotg  = &rot->grp[g];
950             bFlex = ISFLEX(rotg);
951
952             /* For flexible axis rotation we use RMSD fitting to determine the
953              * actual angle of the rotation group */
954             if (bFlex || erotgFitPOT == rotg->eFittype)
955             {
956                 sprintf(buf, "theta_fit%d", g);
957             }
958             else
959             {
960                 sprintf(buf, "theta_av%d", g);
961             }
962             add_to_string_aligned(&LegendStr, buf);
963             sprintf(buf2, "%s (degrees)", buf);
964             setname[nsets] = strdup(buf2);
965             nsets++;
966
967             sprintf(buf, "tau%d", g);
968             add_to_string_aligned(&LegendStr, buf);
969             sprintf(buf2, "%s (kJ/mol)", buf);
970             setname[nsets] = strdup(buf2);
971             nsets++;
972
973             sprintf(buf, "energy%d", g);
974             add_to_string_aligned(&LegendStr, buf);
975             sprintf(buf2, "%s (kJ/mol)", buf);
976             setname[nsets] = strdup(buf2);
977             nsets++;
978         }
979         fprintf(fp, "#\n");
980
981         if (nsets > 1)
982         {
983             xvgr_legend(fp, nsets, setname, oenv);
984         }
985         sfree(setname);
986
987         fprintf(fp, "#\n# Legend for the following data columns:\n");
988         fprintf(fp, "%s\n", LegendStr);
989         sfree(LegendStr);
990
991         fflush(fp);
992     }
993
994     return fp;
995 }
996
997
998 /* Call on master only */
999 static FILE *open_angles_out(const char *fn, t_rot *rot)
1000 {
1001     int             g, i;
1002     FILE           *fp;
1003     t_rotgrp       *rotg;
1004     gmx_enfrotgrp_t erg;        /* Pointer to enforced rotation group data */
1005     char            buf[100];
1006
1007
1008     if (rot->enfrot->Flags & MD_APPENDFILES)
1009     {
1010         fp = gmx_fio_fopen(fn, "a");
1011     }
1012     else
1013     {
1014         /* Open output file and write some information about it's structure: */
1015         fp = open_output_file(fn, rot->nstsout, "rotation group angles");
1016         fprintf(fp, "# All angles given in degrees, time in ps.\n");
1017         for (g = 0; g < rot->ngrp; g++)
1018         {
1019             rotg = &rot->grp[g];
1020             erg  = rotg->enfrotgrp;
1021
1022             /* Output for this group happens only if potential type is flexible or
1023              * if fit type is potential! */
1024             if (ISFLEX(rotg) || (erotgFitPOT == rotg->eFittype) )
1025             {
1026                 if (ISFLEX(rotg))
1027                 {
1028                     sprintf(buf, " slab distance %f nm, ", rotg->slab_dist);
1029                 }
1030                 else
1031                 {
1032                     buf[0] = '\0';
1033                 }
1034
1035                 fprintf(fp, "#\n# ROTATION GROUP %d '%s',%s fit type '%s'.\n",
1036                         g, erotg_names[rotg->eType], buf, erotg_fitnames[rotg->eFittype]);
1037
1038                 /* Special type of fitting using the potential minimum. This is
1039                  * done for the whole group only, not for the individual slabs. */
1040                 if (erotgFitPOT == rotg->eFittype)
1041                 {
1042                     fprintf(fp, "#    To obtain theta_fit%d, the potential is evaluated for %d angles around theta_ref%d\n", g, rotg->PotAngle_nstep, g);
1043                     fprintf(fp, "#    The fit angle in the rotation standard outfile is the one with minimal energy E(theta_fit) [kJ/mol].\n");
1044                     fprintf(fp, "#\n");
1045                 }
1046
1047                 fprintf(fp, "# Legend for the group %d data columns:\n", g);
1048                 fprintf(fp, "#     ");
1049                 print_aligned_short(fp, "time");
1050                 print_aligned_short(fp, "grp");
1051                 print_aligned(fp, "theta_ref");
1052
1053                 if (erotgFitPOT == rotg->eFittype)
1054                 {
1055                     /* Output the set of angles around the reference angle */
1056                     for (i = 0; i < rotg->PotAngle_nstep; i++)
1057                     {
1058                         sprintf(buf, "E(%g)", erg->PotAngleFit->degangle[i]);
1059                         print_aligned(fp, buf);
1060                     }
1061                 }
1062                 else
1063                 {
1064                     /* Output fit angle for each slab */
1065                     print_aligned_short(fp, "slab");
1066                     print_aligned_short(fp, "atoms");
1067                     print_aligned(fp, "theta_fit");
1068                     print_aligned_short(fp, "slab");
1069                     print_aligned_short(fp, "atoms");
1070                     print_aligned(fp, "theta_fit");
1071                     fprintf(fp, " ...");
1072                 }
1073                 fprintf(fp, "\n");
1074             }
1075         }
1076         fflush(fp);
1077     }
1078
1079     return fp;
1080 }
1081
1082
1083 /* Open torque output file and write some information about it's structure.
1084  * Call on master only */
1085 static FILE *open_torque_out(const char *fn, t_rot *rot)
1086 {
1087     FILE      *fp;
1088     int        g;
1089     t_rotgrp  *rotg;
1090
1091
1092     if (rot->enfrot->Flags & MD_APPENDFILES)
1093     {
1094         fp = gmx_fio_fopen(fn, "a");
1095     }
1096     else
1097     {
1098         fp = open_output_file(fn, rot->nstsout, "torques");
1099
1100         for (g = 0; g < rot->ngrp; g++)
1101         {
1102             rotg = &rot->grp[g];
1103             if (ISFLEX(rotg))
1104             {
1105                 fprintf(fp, "# Rotation group %d (%s), slab distance %f nm.\n", g, erotg_names[rotg->eType], rotg->slab_dist);
1106                 fprintf(fp, "# The scalar tau is the torque (kJ/mol) in the direction of the rotation vector.\n");
1107                 fprintf(fp, "# To obtain the vectorial torque, multiply tau with\n");
1108                 fprintf(fp, "# rot_vec%d            %10.3e %10.3e %10.3e\n", g, rotg->vec[XX], rotg->vec[YY], rotg->vec[ZZ]);
1109                 fprintf(fp, "#\n");
1110             }
1111         }
1112         fprintf(fp, "# Legend for the following data columns: (tau=torque for that slab):\n");
1113         fprintf(fp, "#     ");
1114         print_aligned_short(fp, "t");
1115         print_aligned_short(fp, "grp");
1116         print_aligned_short(fp, "slab");
1117         print_aligned(fp, "tau");
1118         print_aligned_short(fp, "slab");
1119         print_aligned(fp, "tau");
1120         fprintf(fp, " ...\n");
1121         fflush(fp);
1122     }
1123
1124     return fp;
1125 }
1126
1127
1128 static void swap_val(double* vec, int i, int j)
1129 {
1130     double tmp = vec[j];
1131
1132
1133     vec[j] = vec[i];
1134     vec[i] = tmp;
1135 }
1136
1137
1138 static void swap_col(double **mat, int i, int j)
1139 {
1140     double tmp[3] = {mat[0][j], mat[1][j], mat[2][j]};
1141
1142
1143     mat[0][j] = mat[0][i];
1144     mat[1][j] = mat[1][i];
1145     mat[2][j] = mat[2][i];
1146
1147     mat[0][i] = tmp[0];
1148     mat[1][i] = tmp[1];
1149     mat[2][i] = tmp[2];
1150 }
1151
1152
1153 /* Eigenvectors are stored in columns of eigen_vec */
1154 static void diagonalize_symmetric(
1155         double **matrix,
1156         double **eigen_vec,
1157         double   eigenval[3])
1158 {
1159     int n_rot;
1160
1161
1162     jacobi(matrix, 3, eigenval, eigen_vec, &n_rot);
1163
1164     /* sort in ascending order */
1165     if (eigenval[0] > eigenval[1])
1166     {
1167         swap_val(eigenval, 0, 1);
1168         swap_col(eigen_vec, 0, 1);
1169     }
1170     if (eigenval[1] > eigenval[2])
1171     {
1172         swap_val(eigenval, 1, 2);
1173         swap_col(eigen_vec, 1, 2);
1174     }
1175     if (eigenval[0] > eigenval[1])
1176     {
1177         swap_val(eigenval, 0, 1);
1178         swap_col(eigen_vec, 0, 1);
1179     }
1180 }
1181
1182
1183 static void align_with_z(
1184         rvec* s,           /* Structure to align */
1185         int   natoms,
1186         rvec  axis)
1187 {
1188     int     i, j, k;
1189     rvec    zet         = {0.0, 0.0, 1.0};
1190     rvec    rot_axis    = {0.0, 0.0, 0.0};
1191     rvec   *rotated_str = NULL;
1192     real    ooanorm;
1193     real    angle;
1194     matrix  rotmat;
1195
1196
1197     snew(rotated_str, natoms);
1198
1199     /* Normalize the axis */
1200     ooanorm = 1.0/norm(axis);
1201     svmul(ooanorm, axis, axis);
1202
1203     /* Calculate the angle for the fitting procedure */
1204     cprod(axis, zet, rot_axis);
1205     angle = acos(axis[2]);
1206     if (angle < 0.0)
1207     {
1208         angle += M_PI;
1209     }
1210
1211     /* Calculate the rotation matrix */
1212     calc_rotmat(rot_axis, angle*180.0/M_PI, rotmat);
1213
1214     /* Apply the rotation matrix to s */
1215     for (i = 0; i < natoms; i++)
1216     {
1217         for (j = 0; j < 3; j++)
1218         {
1219             for (k = 0; k < 3; k++)
1220             {
1221                 rotated_str[i][j] += rotmat[j][k]*s[i][k];
1222             }
1223         }
1224     }
1225
1226     /* Rewrite the rotated structure to s */
1227     for (i = 0; i < natoms; i++)
1228     {
1229         for (j = 0; j < 3; j++)
1230         {
1231             s[i][j] = rotated_str[i][j];
1232         }
1233     }
1234
1235     sfree(rotated_str);
1236 }
1237
1238
1239 static void calc_correl_matrix(rvec* Xstr, rvec* Ystr, double** Rmat, int natoms)
1240 {
1241     int i, j, k;
1242
1243
1244     for (i = 0; i < 3; i++)
1245     {
1246         for (j = 0; j < 3; j++)
1247         {
1248             Rmat[i][j] = 0.0;
1249         }
1250     }
1251
1252     for (i = 0; i < 3; i++)
1253     {
1254         for (j = 0; j < 3; j++)
1255         {
1256             for (k = 0; k < natoms; k++)
1257             {
1258                 Rmat[i][j] += Ystr[k][i] * Xstr[k][j];
1259             }
1260         }
1261     }
1262 }
1263
1264
1265 static void weigh_coords(rvec* str, real* weight, int natoms)
1266 {
1267     int i, j;
1268
1269
1270     for (i = 0; i < natoms; i++)
1271     {
1272         for (j = 0; j < 3; j++)
1273         {
1274             str[i][j] *= sqrt(weight[i]);
1275         }
1276     }
1277 }
1278
1279
1280 static real opt_angle_analytic(
1281         rvec* ref_s,
1282         rvec* act_s,
1283         real* weight,
1284         int   natoms,
1285         rvec  ref_com,
1286         rvec  act_com,
1287         rvec  axis)
1288 {
1289     int      i, j, k;
1290     rvec    *ref_s_1 = NULL;
1291     rvec    *act_s_1 = NULL;
1292     rvec     shift;
1293     double **Rmat, **RtR, **eigvec;
1294     double   eigval[3];
1295     double   V[3][3], WS[3][3];
1296     double   rot_matrix[3][3];
1297     double   opt_angle;
1298
1299
1300     /* Do not change the original coordinates */
1301     snew(ref_s_1, natoms);
1302     snew(act_s_1, natoms);
1303     for (i = 0; i < natoms; i++)
1304     {
1305         copy_rvec(ref_s[i], ref_s_1[i]);
1306         copy_rvec(act_s[i], act_s_1[i]);
1307     }
1308
1309     /* Translate the structures to the origin */
1310     shift[XX] = -ref_com[XX];
1311     shift[YY] = -ref_com[YY];
1312     shift[ZZ] = -ref_com[ZZ];
1313     translate_x(ref_s_1, natoms, shift);
1314
1315     shift[XX] = -act_com[XX];
1316     shift[YY] = -act_com[YY];
1317     shift[ZZ] = -act_com[ZZ];
1318     translate_x(act_s_1, natoms, shift);
1319
1320     /* Align rotation axis with z */
1321     align_with_z(ref_s_1, natoms, axis);
1322     align_with_z(act_s_1, natoms, axis);
1323
1324     /* Correlation matrix */
1325     Rmat = allocate_square_matrix(3);
1326
1327     for (i = 0; i < natoms; i++)
1328     {
1329         ref_s_1[i][2] = 0.0;
1330         act_s_1[i][2] = 0.0;
1331     }
1332
1333     /* Weight positions with sqrt(weight) */
1334     if (NULL != weight)
1335     {
1336         weigh_coords(ref_s_1, weight, natoms);
1337         weigh_coords(act_s_1, weight, natoms);
1338     }
1339
1340     /* Calculate correlation matrices R=YXt (X=ref_s; Y=act_s) */
1341     calc_correl_matrix(ref_s_1, act_s_1, Rmat, natoms);
1342
1343     /* Calculate RtR */
1344     RtR = allocate_square_matrix(3);
1345     for (i = 0; i < 3; i++)
1346     {
1347         for (j = 0; j < 3; j++)
1348         {
1349             for (k = 0; k < 3; k++)
1350             {
1351                 RtR[i][j] += Rmat[k][i] * Rmat[k][j];
1352             }
1353         }
1354     }
1355     /* Diagonalize RtR */
1356     snew(eigvec, 3);
1357     for (i = 0; i < 3; i++)
1358     {
1359         snew(eigvec[i], 3);
1360     }
1361
1362     diagonalize_symmetric(RtR, eigvec, eigval);
1363     swap_col(eigvec, 0, 1);
1364     swap_col(eigvec, 1, 2);
1365     swap_val(eigval, 0, 1);
1366     swap_val(eigval, 1, 2);
1367
1368     /* Calculate V */
1369     for (i = 0; i < 3; i++)
1370     {
1371         for (j = 0; j < 3; j++)
1372         {
1373             V[i][j]  = 0.0;
1374             WS[i][j] = 0.0;
1375         }
1376     }
1377
1378     for (i = 0; i < 2; i++)
1379     {
1380         for (j = 0; j < 2; j++)
1381         {
1382             WS[i][j] = eigvec[i][j] / sqrt(eigval[j]);
1383         }
1384     }
1385
1386     for (i = 0; i < 3; i++)
1387     {
1388         for (j = 0; j < 3; j++)
1389         {
1390             for (k = 0; k < 3; k++)
1391             {
1392                 V[i][j] += Rmat[i][k]*WS[k][j];
1393             }
1394         }
1395     }
1396     free_square_matrix(Rmat, 3);
1397
1398     /* Calculate optimal rotation matrix */
1399     for (i = 0; i < 3; i++)
1400     {
1401         for (j = 0; j < 3; j++)
1402         {
1403             rot_matrix[i][j] = 0.0;
1404         }
1405     }
1406
1407     for (i = 0; i < 3; i++)
1408     {
1409         for (j = 0; j < 3; j++)
1410         {
1411             for (k = 0; k < 3; k++)
1412             {
1413                 rot_matrix[i][j] += eigvec[i][k]*V[j][k];
1414             }
1415         }
1416     }
1417     rot_matrix[2][2] = 1.0;
1418
1419     /* In some cases abs(rot_matrix[0][0]) can be slighly larger
1420      * than unity due to numerical inacurracies. To be able to calculate
1421      * the acos function, we put these values back in range. */
1422     if (rot_matrix[0][0] > 1.0)
1423     {
1424         rot_matrix[0][0] = 1.0;
1425     }
1426     else if (rot_matrix[0][0] < -1.0)
1427     {
1428         rot_matrix[0][0] = -1.0;
1429     }
1430
1431     /* Determine the optimal rotation angle: */
1432     opt_angle = (-1.0)*acos(rot_matrix[0][0])*180.0/M_PI;
1433     if (rot_matrix[0][1] < 0.0)
1434     {
1435         opt_angle = (-1.0)*opt_angle;
1436     }
1437
1438     /* Give back some memory */
1439     free_square_matrix(RtR, 3);
1440     sfree(ref_s_1);
1441     sfree(act_s_1);
1442     for (i = 0; i < 3; i++)
1443     {
1444         sfree(eigvec[i]);
1445     }
1446     sfree(eigvec);
1447
1448     return (real) opt_angle;
1449 }
1450
1451
1452 /* Determine angle of the group by RMSD fit to the reference */
1453 /* Not parallelized, call this routine only on the master */
1454 static real flex_fit_angle(t_rotgrp *rotg)
1455 {
1456     int             i;
1457     rvec           *fitcoords = NULL;
1458     rvec            center;     /* Center of positions passed to the fit routine */
1459     real            fitangle;   /* Angle of the rotation group derived by fitting */
1460     rvec            coord;
1461     real            scal;
1462     gmx_enfrotgrp_t erg;        /* Pointer to enforced rotation group data */
1463
1464
1465     erg = rotg->enfrotgrp;
1466
1467     /* Get the center of the rotation group.
1468      * Note, again, erg->xc has been sorted in do_flexible */
1469     get_center(erg->xc, erg->mc_sorted, rotg->nat, center);
1470
1471     /* === Determine the optimal fit angle for the rotation group === */
1472     if (rotg->eFittype == erotgFitNORM)
1473     {
1474         /* Normalize every position to it's reference length */
1475         for (i = 0; i < rotg->nat; i++)
1476         {
1477             /* Put the center of the positions into the origin */
1478             rvec_sub(erg->xc[i], center, coord);
1479             /* Determine the scaling factor for the length: */
1480             scal = erg->xc_ref_length[erg->xc_sortind[i]] / norm(coord);
1481             /* Get position, multiply with the scaling factor and save  */
1482             svmul(scal, coord, erg->xc_norm[i]);
1483         }
1484         fitcoords = erg->xc_norm;
1485     }
1486     else
1487     {
1488         fitcoords = erg->xc;
1489     }
1490     /* From the point of view of the current positions, the reference has rotated
1491      * backwards. Since we output the angle relative to the fixed reference,
1492      * we need the minus sign. */
1493     fitangle = -opt_angle_analytic(erg->xc_ref_sorted, fitcoords, erg->mc_sorted,
1494                                    rotg->nat, erg->xc_ref_center, center, rotg->vec);
1495
1496     return fitangle;
1497 }
1498
1499
1500 /* Determine actual angle of each slab by RMSD fit to the reference */
1501 /* Not parallelized, call this routine only on the master */
1502 static void flex_fit_angle_perslab(
1503         int       g,
1504         t_rotgrp *rotg,
1505         double    t,
1506         real      degangle,
1507         FILE     *fp)
1508 {
1509     int             i, l, n, islab, ind;
1510     rvec            curr_x, ref_x;
1511     rvec            act_center; /* Center of actual positions that are passed to the fit routine */
1512     rvec            ref_center; /* Same for the reference positions */
1513     real            fitangle;   /* Angle of a slab derived from an RMSD fit to
1514                                  * the reference structure at t=0  */
1515     t_gmx_slabdata *sd;
1516     gmx_enfrotgrp_t erg;        /* Pointer to enforced rotation group data */
1517     real            OOm_av;     /* 1/average_mass of a rotation group atom */
1518     real            m_rel;      /* Relative mass of a rotation group atom  */
1519
1520
1521     erg = rotg->enfrotgrp;
1522
1523     /* Average mass of a rotation group atom: */
1524     OOm_av = erg->invmass*rotg->nat;
1525
1526     /**********************************/
1527     /* First collect the data we need */
1528     /**********************************/
1529
1530     /* Collect the data for the individual slabs */
1531     for (n = erg->slab_first; n <= erg->slab_last; n++)
1532     {
1533         islab   = n - erg->slab_first; /* slab index */
1534         sd      = &(rotg->enfrotgrp->slab_data[islab]);
1535         sd->nat = erg->lastatom[islab]-erg->firstatom[islab]+1;
1536         ind     = 0;
1537
1538         /* Loop over the relevant atoms in the slab */
1539         for (l = erg->firstatom[islab]; l <= erg->lastatom[islab]; l++)
1540         {
1541             /* Current position of this atom: x[ii][XX/YY/ZZ] */
1542             copy_rvec(erg->xc[l], curr_x);
1543
1544             /* The (unrotated) reference position of this atom is copied to ref_x.
1545              * Beware, the xc coords have been sorted in do_flexible */
1546             copy_rvec(erg->xc_ref_sorted[l], ref_x);
1547
1548             /* Save data for doing angular RMSD fit later */
1549             /* Save the current atom position */
1550             copy_rvec(curr_x, sd->x[ind]);
1551             /* Save the corresponding reference position */
1552             copy_rvec(ref_x, sd->ref[ind]);
1553
1554             /* Maybe also mass-weighting was requested. If yes, additionally
1555              * multiply the weights with the relative mass of the atom. If not,
1556              * multiply with unity. */
1557             m_rel = erg->mc_sorted[l]*OOm_av;
1558
1559             /* Save the weight for this atom in this slab */
1560             sd->weight[ind] = gaussian_weight(curr_x, rotg, n) * m_rel;
1561
1562             /* Next atom in this slab */
1563             ind++;
1564         }
1565     }
1566
1567     /******************************/
1568     /* Now do the fit calculation */
1569     /******************************/
1570
1571     fprintf(fp, "%12.3e%6d%12.3f", t, g, degangle);
1572
1573     /* === Now do RMSD fitting for each slab === */
1574     /* We require at least SLAB_MIN_ATOMS in a slab, such that the fit makes sense. */
1575 #define SLAB_MIN_ATOMS 4
1576
1577     for (n = erg->slab_first; n <= erg->slab_last; n++)
1578     {
1579         islab = n - erg->slab_first; /* slab index */
1580         sd    = &(rotg->enfrotgrp->slab_data[islab]);
1581         if (sd->nat >= SLAB_MIN_ATOMS)
1582         {
1583             /* Get the center of the slabs reference and current positions */
1584             get_center(sd->ref, sd->weight, sd->nat, ref_center);
1585             get_center(sd->x, sd->weight, sd->nat, act_center);
1586             if (rotg->eFittype == erotgFitNORM)
1587             {
1588                 /* Normalize every position to it's reference length
1589                  * prior to performing the fit */
1590                 for (i = 0; i < sd->nat; i++) /* Center */
1591                 {
1592                     rvec_dec(sd->ref[i], ref_center);
1593                     rvec_dec(sd->x[i], act_center);
1594                     /* Normalize x_i such that it gets the same length as ref_i */
1595                     svmul( norm(sd->ref[i])/norm(sd->x[i]), sd->x[i], sd->x[i] );
1596                 }
1597                 /* We already subtracted the centers */
1598                 clear_rvec(ref_center);
1599                 clear_rvec(act_center);
1600             }
1601             fitangle = -opt_angle_analytic(sd->ref, sd->x, sd->weight, sd->nat,
1602                                            ref_center, act_center, rotg->vec);
1603             fprintf(fp, "%6d%6d%12.3f", n, sd->nat, fitangle);
1604         }
1605     }
1606     fprintf(fp, "\n");
1607
1608 #undef SLAB_MIN_ATOMS
1609 }
1610
1611
1612 /* Shift x with is */
1613 static gmx_inline void shift_single_coord(matrix box, rvec x, const ivec is)
1614 {
1615     int tx, ty, tz;
1616
1617
1618     tx = is[XX];
1619     ty = is[YY];
1620     tz = is[ZZ];
1621
1622     if (TRICLINIC(box))
1623     {
1624         x[XX] += tx*box[XX][XX]+ty*box[YY][XX]+tz*box[ZZ][XX];
1625         x[YY] += ty*box[YY][YY]+tz*box[ZZ][YY];
1626         x[ZZ] += tz*box[ZZ][ZZ];
1627     }
1628     else
1629     {
1630         x[XX] += tx*box[XX][XX];
1631         x[YY] += ty*box[YY][YY];
1632         x[ZZ] += tz*box[ZZ][ZZ];
1633     }
1634 }
1635
1636
1637 /* Determine the 'home' slab of this atom which is the
1638  * slab with the highest Gaussian weight of all */
1639 #define round(a) (int)(a+0.5)
1640 static gmx_inline int get_homeslab(
1641         rvec curr_x,   /* The position for which the home slab shall be determined */
1642         rvec rotvec,   /* The rotation vector */
1643         real slabdist) /* The slab distance */
1644 {
1645     real dist;
1646
1647
1648     /* The distance of the atom to the coordinate center (where the
1649      * slab with index 0) is */
1650     dist = iprod(rotvec, curr_x);
1651
1652     return round(dist / slabdist);
1653 }
1654
1655
1656 /* For a local atom determine the relevant slabs, i.e. slabs in
1657  * which the gaussian is larger than min_gaussian
1658  */
1659 static int get_single_atom_gaussians(
1660         rvec       curr_x,
1661         t_rotgrp  *rotg)
1662 {
1663     int             slab, homeslab;
1664     real            g;
1665     int             count = 0;
1666     gmx_enfrotgrp_t erg;      /* Pointer to enforced rotation group data */
1667
1668
1669     erg = rotg->enfrotgrp;
1670
1671     /* Determine the 'home' slab of this atom: */
1672     homeslab = get_homeslab(curr_x, rotg->vec, rotg->slab_dist);
1673
1674     /* First determine the weight in the atoms home slab: */
1675     g = gaussian_weight(curr_x, rotg, homeslab);
1676
1677     erg->gn_atom[count]    = g;
1678     erg->gn_slabind[count] = homeslab;
1679     count++;
1680
1681
1682     /* Determine the max slab */
1683     slab = homeslab;
1684     while (g > rotg->min_gaussian)
1685     {
1686         slab++;
1687         g = gaussian_weight(curr_x, rotg, slab);
1688         erg->gn_slabind[count] = slab;
1689         erg->gn_atom[count]    = g;
1690         count++;
1691     }
1692     count--;
1693
1694     /* Determine the min slab */
1695     slab = homeslab;
1696     do
1697     {
1698         slab--;
1699         g = gaussian_weight(curr_x, rotg, slab);
1700         erg->gn_slabind[count] = slab;
1701         erg->gn_atom[count]    = g;
1702         count++;
1703     }
1704     while (g > rotg->min_gaussian);
1705     count--;
1706
1707     return count;
1708 }
1709
1710
1711 static void flex2_precalc_inner_sum(t_rotgrp *rotg)
1712 {
1713     int             i, n, islab;
1714     rvec            xi;       /* positions in the i-sum                        */
1715     rvec            xcn, ycn; /* the current and the reference slab centers    */
1716     real            gaussian_xi;
1717     rvec            yi0;
1718     rvec            rin;     /* Helper variables                              */
1719     real            fac, fac2;
1720     rvec            innersumvec;
1721     real            OOpsii, OOpsiistar;
1722     real            sin_rin; /* s_ii.r_ii */
1723     rvec            s_in, tmpvec, tmpvec2;
1724     real            mi, wi;  /* Mass-weighting of the positions                 */
1725     real            N_M;     /* N/M                                             */
1726     gmx_enfrotgrp_t erg;     /* Pointer to enforced rotation group data */
1727
1728
1729     erg = rotg->enfrotgrp;
1730     N_M = rotg->nat * erg->invmass;
1731
1732     /* Loop over all slabs that contain something */
1733     for (n = erg->slab_first; n <= erg->slab_last; n++)
1734     {
1735         islab = n - erg->slab_first; /* slab index */
1736
1737         /* The current center of this slab is saved in xcn: */
1738         copy_rvec(erg->slab_center[islab], xcn);
1739         /* ... and the reference center in ycn: */
1740         copy_rvec(erg->slab_center_ref[islab+erg->slab_buffer], ycn);
1741
1742         /*** D. Calculate the whole inner sum used for second and third sum */
1743         /* For slab n, we need to loop over all atoms i again. Since we sorted
1744          * the atoms with respect to the rotation vector, we know that it is sufficient
1745          * to calculate from firstatom to lastatom only. All other contributions will
1746          * be very small. */
1747         clear_rvec(innersumvec);
1748         for (i = erg->firstatom[islab]; i <= erg->lastatom[islab]; i++)
1749         {
1750             /* Coordinate xi of this atom */
1751             copy_rvec(erg->xc[i], xi);
1752
1753             /* The i-weights */
1754             gaussian_xi = gaussian_weight(xi, rotg, n);
1755             mi          = erg->mc_sorted[i]; /* need the sorted mass here */
1756             wi          = N_M*mi;
1757
1758             /* Calculate rin */
1759             copy_rvec(erg->xc_ref_sorted[i], yi0); /* Reference position yi0   */
1760             rvec_sub(yi0, ycn, tmpvec2);           /* tmpvec2 = yi0 - ycn      */
1761             mvmul(erg->rotmat, tmpvec2, rin);      /* rin = Omega.(yi0 - ycn)  */
1762
1763             /* Calculate psi_i* and sin */
1764             rvec_sub(xi, xcn, tmpvec2);           /* tmpvec2 = xi - xcn       */
1765             cprod(rotg->vec, tmpvec2, tmpvec);    /* tmpvec = v x (xi - xcn)  */
1766             OOpsiistar = norm2(tmpvec)+rotg->eps; /* OOpsii* = 1/psii* = |v x (xi-xcn)|^2 + eps */
1767             OOpsii     = norm(tmpvec);            /* OOpsii = 1 / psii = |v x (xi - xcn)| */
1768
1769             /*                           *         v x (xi - xcn)          */
1770             unitv(tmpvec, s_in);        /*  sin = ----------------         */
1771                                         /*        |v x (xi - xcn)|         */
1772
1773             sin_rin = iprod(s_in, rin); /* sin_rin = sin . rin             */
1774
1775             /* Now the whole sum */
1776             fac = OOpsii/OOpsiistar;
1777             svmul(fac, rin, tmpvec);
1778             fac2 = fac*fac*OOpsii;
1779             svmul(fac2*sin_rin, s_in, tmpvec2);
1780             rvec_dec(tmpvec, tmpvec2);
1781
1782             svmul(wi*gaussian_xi*sin_rin, tmpvec, tmpvec2);
1783
1784             rvec_inc(innersumvec, tmpvec2);
1785         } /* now we have the inner sum, used both for sum2 and sum3 */
1786
1787         /* Save it to be used in do_flex2_lowlevel */
1788         copy_rvec(innersumvec, erg->slab_innersumvec[islab]);
1789     } /* END of loop over slabs */
1790 }
1791
1792
1793 static void flex_precalc_inner_sum(t_rotgrp *rotg)
1794 {
1795     int             i, n, islab;
1796     rvec            xi;       /* position                                      */
1797     rvec            xcn, ycn; /* the current and the reference slab centers    */
1798     rvec            qin, rin; /* q_i^n and r_i^n                               */
1799     real            bin;
1800     rvec            tmpvec;
1801     rvec            innersumvec; /* Inner part of sum_n2                          */
1802     real            gaussian_xi; /* Gaussian weight gn(xi)                        */
1803     real            mi, wi;      /* Mass-weighting of the positions               */
1804     real            N_M;         /* N/M                                           */
1805
1806     gmx_enfrotgrp_t erg;         /* Pointer to enforced rotation group data */
1807
1808
1809     erg = rotg->enfrotgrp;
1810     N_M = rotg->nat * erg->invmass;
1811
1812     /* Loop over all slabs that contain something */
1813     for (n = erg->slab_first; n <= erg->slab_last; n++)
1814     {
1815         islab = n - erg->slab_first; /* slab index */
1816
1817         /* The current center of this slab is saved in xcn: */
1818         copy_rvec(erg->slab_center[islab], xcn);
1819         /* ... and the reference center in ycn: */
1820         copy_rvec(erg->slab_center_ref[islab+erg->slab_buffer], ycn);
1821
1822         /* For slab n, we need to loop over all atoms i again. Since we sorted
1823          * the atoms with respect to the rotation vector, we know that it is sufficient
1824          * to calculate from firstatom to lastatom only. All other contributions will
1825          * be very small. */
1826         clear_rvec(innersumvec);
1827         for (i = erg->firstatom[islab]; i <= erg->lastatom[islab]; i++)
1828         {
1829             /* Coordinate xi of this atom */
1830             copy_rvec(erg->xc[i], xi);
1831
1832             /* The i-weights */
1833             gaussian_xi = gaussian_weight(xi, rotg, n);
1834             mi          = erg->mc_sorted[i]; /* need the sorted mass here */
1835             wi          = N_M*mi;
1836
1837             /* Calculate rin and qin */
1838             rvec_sub(erg->xc_ref_sorted[i], ycn, tmpvec); /* tmpvec = yi0-ycn */
1839             mvmul(erg->rotmat, tmpvec, rin);              /* rin = Omega.(yi0 - ycn)  */
1840             cprod(rotg->vec, rin, tmpvec);                /* tmpvec = v x Omega*(yi0-ycn) */
1841
1842             /*                                *        v x Omega*(yi0-ycn)    */
1843             unitv(tmpvec, qin);              /* qin = ---------------------   */
1844                                              /*       |v x Omega*(yi0-ycn)|   */
1845
1846             /* Calculate bin */
1847             rvec_sub(xi, xcn, tmpvec);            /* tmpvec = xi-xcn          */
1848             bin = iprod(qin, tmpvec);             /* bin  = qin*(xi-xcn)      */
1849
1850             svmul(wi*gaussian_xi*bin, qin, tmpvec);
1851
1852             /* Add this contribution to the inner sum: */
1853             rvec_add(innersumvec, tmpvec, innersumvec);
1854         } /* now we have the inner sum vector S^n for this slab */
1855           /* Save it to be used in do_flex_lowlevel */
1856         copy_rvec(innersumvec, erg->slab_innersumvec[islab]);
1857     }
1858 }
1859
1860
1861 static real do_flex2_lowlevel(
1862         t_rotgrp  *rotg,
1863         real       sigma,   /* The Gaussian width sigma */
1864         rvec       x[],
1865         gmx_bool   bOutstepRot,
1866         gmx_bool   bOutstepSlab,
1867         matrix     box)
1868 {
1869     int             count, ic, ii, j, m, n, islab, iigrp, ifit;
1870     rvec            xj;          /* position in the i-sum                         */
1871     rvec            yj0;         /* the reference position in the j-sum           */
1872     rvec            xcn, ycn;    /* the current and the reference slab centers    */
1873     real            V;           /* This node's part of the rotation pot. energy  */
1874     real            gaussian_xj; /* Gaussian weight                               */
1875     real            beta;
1876
1877     real            numerator, fit_numerator;
1878     rvec            rjn, fit_rjn; /* Helper variables                              */
1879     real            fac, fac2;
1880
1881     real            OOpsij, OOpsijstar;
1882     real            OOsigma2; /* 1/(sigma^2)                                   */
1883     real            sjn_rjn;
1884     real            betasigpsi;
1885     rvec            sjn, tmpvec, tmpvec2, yj0_ycn;
1886     rvec            sum1vec_part, sum1vec, sum2vec_part, sum2vec, sum3vec, sum4vec, innersumvec;
1887     real            sum3, sum4;
1888     gmx_enfrotgrp_t erg;     /* Pointer to enforced rotation group data       */
1889     real            mj, wj;  /* Mass-weighting of the positions               */
1890     real            N_M;     /* N/M                                           */
1891     real            Wjn;     /* g_n(x_j) m_j / Mjn                            */
1892     gmx_bool        bCalcPotFit;
1893
1894     /* To calculate the torque per slab */
1895     rvec slab_force;         /* Single force from slab n on one atom          */
1896     rvec slab_sum1vec_part;
1897     real slab_sum3part, slab_sum4part;
1898     rvec slab_sum1vec, slab_sum2vec, slab_sum3vec, slab_sum4vec;
1899
1900
1901     erg = rotg->enfrotgrp;
1902
1903     /* Pre-calculate the inner sums, so that we do not have to calculate
1904      * them again for every atom */
1905     flex2_precalc_inner_sum(rotg);
1906
1907     bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT == rotg->eFittype);
1908
1909     /********************************************************/
1910     /* Main loop over all local atoms of the rotation group */
1911     /********************************************************/
1912     N_M      = rotg->nat * erg->invmass;
1913     V        = 0.0;
1914     OOsigma2 = 1.0 / (sigma*sigma);
1915     for (j = 0; j < erg->nat_loc; j++)
1916     {
1917         /* Local index of a rotation group atom  */
1918         ii = erg->ind_loc[j];
1919         /* Position of this atom in the collective array */
1920         iigrp = erg->xc_ref_ind[j];
1921         /* Mass-weighting */
1922         mj = erg->mc[iigrp];  /* need the unsorted mass here */
1923         wj = N_M*mj;
1924
1925         /* Current position of this atom: x[ii][XX/YY/ZZ]
1926          * Note that erg->xc_center contains the center of mass in case the flex2-t
1927          * potential was chosen. For the flex2 potential erg->xc_center must be
1928          * zero. */
1929         rvec_sub(x[ii], erg->xc_center, xj);
1930
1931         /* Shift this atom such that it is near its reference */
1932         shift_single_coord(box, xj, erg->xc_shifts[iigrp]);
1933
1934         /* Determine the slabs to loop over, i.e. the ones with contributions
1935          * larger than min_gaussian */
1936         count = get_single_atom_gaussians(xj, rotg);
1937
1938         clear_rvec(sum1vec_part);
1939         clear_rvec(sum2vec_part);
1940         sum3 = 0.0;
1941         sum4 = 0.0;
1942         /* Loop over the relevant slabs for this atom */
1943         for (ic = 0; ic < count; ic++)
1944         {
1945             n = erg->gn_slabind[ic];
1946
1947             /* Get the precomputed Gaussian value of curr_slab for curr_x */
1948             gaussian_xj = erg->gn_atom[ic];
1949
1950             islab = n - erg->slab_first; /* slab index */
1951
1952             /* The (unrotated) reference position of this atom is copied to yj0: */
1953             copy_rvec(rotg->x_ref[iigrp], yj0);
1954
1955             beta = calc_beta(xj, rotg, n);
1956
1957             /* The current center of this slab is saved in xcn: */
1958             copy_rvec(erg->slab_center[islab], xcn);
1959             /* ... and the reference center in ycn: */
1960             copy_rvec(erg->slab_center_ref[islab+erg->slab_buffer], ycn);
1961
1962             rvec_sub(yj0, ycn, yj0_ycn);          /* yj0_ycn = yj0 - ycn      */
1963
1964             /* Rotate: */
1965             mvmul(erg->rotmat, yj0_ycn, rjn);     /* rjn = Omega.(yj0 - ycn)  */
1966
1967             /* Subtract the slab center from xj */
1968             rvec_sub(xj, xcn, tmpvec2);           /* tmpvec2 = xj - xcn       */
1969
1970             /* In rare cases, when an atom position coincides with a slab center
1971              * (tmpvec2 == 0) we cannot compute the vector product for sjn.
1972              * However, since the atom is located directly on the pivot, this
1973              * slab's contribution to the force on that atom will be zero
1974              * anyway. Therefore, we directly move on to the next slab.       */
1975             if (0 == norm(tmpvec2) )
1976             {
1977                 continue;
1978             }
1979
1980             /* Calculate sjn */
1981             cprod(rotg->vec, tmpvec2, tmpvec);    /* tmpvec = v x (xj - xcn)  */
1982
1983             OOpsijstar = norm2(tmpvec)+rotg->eps; /* OOpsij* = 1/psij* = |v x (xj-xcn)|^2 + eps */
1984
1985             numerator = sqr(iprod(tmpvec, rjn));
1986
1987             /*********************************/
1988             /* Add to the rotation potential */
1989             /*********************************/
1990             V += 0.5*rotg->k*wj*gaussian_xj*numerator/OOpsijstar;
1991
1992             /* If requested, also calculate the potential for a set of angles
1993              * near the current reference angle */
1994             if (bCalcPotFit)
1995             {
1996                 for (ifit = 0; ifit < rotg->PotAngle_nstep; ifit++)
1997                 {
1998                     mvmul(erg->PotAngleFit->rotmat[ifit], yj0_ycn, fit_rjn);
1999                     fit_numerator              = sqr(iprod(tmpvec, fit_rjn));
2000                     erg->PotAngleFit->V[ifit] += 0.5*rotg->k*wj*gaussian_xj*fit_numerator/OOpsijstar;
2001                 }
2002             }
2003
2004             /*************************************/
2005             /* Now calculate the force on atom j */
2006             /*************************************/
2007
2008             OOpsij = norm(tmpvec);    /* OOpsij = 1 / psij = |v x (xj - xcn)| */
2009
2010             /*                              *         v x (xj - xcn)          */
2011             unitv(tmpvec, sjn);            /*  sjn = ----------------         */
2012                                            /*        |v x (xj - xcn)|         */
2013
2014             sjn_rjn = iprod(sjn, rjn);     /* sjn_rjn = sjn . rjn             */
2015
2016
2017             /*** A. Calculate the first of the four sum terms: ****************/
2018             fac = OOpsij/OOpsijstar;
2019             svmul(fac, rjn, tmpvec);
2020             fac2 = fac*fac*OOpsij;
2021             svmul(fac2*sjn_rjn, sjn, tmpvec2);
2022             rvec_dec(tmpvec, tmpvec2);
2023             fac2 = wj*gaussian_xj; /* also needed for sum4 */
2024             svmul(fac2*sjn_rjn, tmpvec, slab_sum1vec_part);
2025             /********************/
2026             /*** Add to sum1: ***/
2027             /********************/
2028             rvec_inc(sum1vec_part, slab_sum1vec_part); /* sum1 still needs to vector multiplied with v */
2029
2030             /*** B. Calculate the forth of the four sum terms: ****************/
2031             betasigpsi = beta*OOsigma2*OOpsij; /* this is also needed for sum3 */
2032             /********************/
2033             /*** Add to sum4: ***/
2034             /********************/
2035             slab_sum4part = fac2*betasigpsi*fac*sjn_rjn*sjn_rjn; /* Note that fac is still valid from above */
2036             sum4         += slab_sum4part;
2037
2038             /*** C. Calculate Wjn for second and third sum */
2039             /* Note that we can safely divide by slab_weights since we check in
2040              * get_slab_centers that it is non-zero. */
2041             Wjn = gaussian_xj*mj/erg->slab_weights[islab];
2042
2043             /* We already have precalculated the inner sum for slab n */
2044             copy_rvec(erg->slab_innersumvec[islab], innersumvec);
2045
2046             /* Weigh the inner sum vector with Wjn */
2047             svmul(Wjn, innersumvec, innersumvec);
2048
2049             /*** E. Calculate the second of the four sum terms: */
2050             /********************/
2051             /*** Add to sum2: ***/
2052             /********************/
2053             rvec_inc(sum2vec_part, innersumvec); /* sum2 still needs to be vector crossproduct'ed with v */
2054
2055             /*** F. Calculate the third of the four sum terms: */
2056             slab_sum3part = betasigpsi * iprod(sjn, innersumvec);
2057             sum3         += slab_sum3part; /* still needs to be multiplied with v */
2058
2059             /*** G. Calculate the torque on the local slab's axis: */
2060             if (bOutstepRot)
2061             {
2062                 /* Sum1 */
2063                 cprod(slab_sum1vec_part, rotg->vec, slab_sum1vec);
2064                 /* Sum2 */
2065                 cprod(innersumvec, rotg->vec, slab_sum2vec);
2066                 /* Sum3 */
2067                 svmul(slab_sum3part, rotg->vec, slab_sum3vec);
2068                 /* Sum4 */
2069                 svmul(slab_sum4part, rotg->vec, slab_sum4vec);
2070
2071                 /* The force on atom ii from slab n only: */
2072                 for (m = 0; m < DIM; m++)
2073                 {
2074                     slab_force[m] = rotg->k * (-slab_sum1vec[m] + slab_sum2vec[m] - slab_sum3vec[m] + 0.5*slab_sum4vec[m]);
2075                 }
2076
2077                 erg->slab_torque_v[islab] += torque(rotg->vec, slab_force, xj, xcn);
2078             }
2079         } /* END of loop over slabs */
2080
2081         /* Construct the four individual parts of the vector sum: */
2082         cprod(sum1vec_part, rotg->vec, sum1vec);      /* sum1vec =   { } x v  */
2083         cprod(sum2vec_part, rotg->vec, sum2vec);      /* sum2vec =   { } x v  */
2084         svmul(sum3, rotg->vec, sum3vec);              /* sum3vec =   { } . v  */
2085         svmul(sum4, rotg->vec, sum4vec);              /* sum4vec =   { } . v  */
2086
2087         /* Store the additional force so that it can be added to the force
2088          * array after the normal forces have been evaluated */
2089         for (m = 0; m < DIM; m++)
2090         {
2091             erg->f_rot_loc[j][m] = rotg->k * (-sum1vec[m] + sum2vec[m] - sum3vec[m] + 0.5*sum4vec[m]);
2092         }
2093
2094 #ifdef SUM_PARTS
2095         fprintf(stderr, "sum1: %15.8f %15.8f %15.8f\n",    -rotg->k*sum1vec[XX],    -rotg->k*sum1vec[YY],    -rotg->k*sum1vec[ZZ]);
2096         fprintf(stderr, "sum2: %15.8f %15.8f %15.8f\n",     rotg->k*sum2vec[XX],     rotg->k*sum2vec[YY],     rotg->k*sum2vec[ZZ]);
2097         fprintf(stderr, "sum3: %15.8f %15.8f %15.8f\n",    -rotg->k*sum3vec[XX],    -rotg->k*sum3vec[YY],    -rotg->k*sum3vec[ZZ]);
2098         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]);
2099 #endif
2100
2101         PRINT_FORCE_J
2102
2103     } /* END of loop over local atoms */
2104
2105     return V;
2106 }
2107
2108
2109 static real do_flex_lowlevel(
2110         t_rotgrp *rotg,
2111         real      sigma,     /* The Gaussian width sigma                      */
2112         rvec      x[],
2113         gmx_bool  bOutstepRot,
2114         gmx_bool  bOutstepSlab,
2115         matrix    box)
2116 {
2117     int             count, ic, ifit, ii, j, m, n, islab, iigrp;
2118     rvec            xj, yj0;                /* current and reference position                */
2119     rvec            xcn, ycn;               /* the current and the reference slab centers    */
2120     rvec            yj0_ycn;                /* yj0 - ycn                                     */
2121     rvec            xj_xcn;                 /* xj - xcn                                      */
2122     rvec            qjn, fit_qjn;           /* q_i^n                                         */
2123     rvec            sum_n1, sum_n2;         /* Two contributions to the rotation force       */
2124     rvec            innersumvec;            /* Inner part of sum_n2                          */
2125     rvec            s_n;
2126     rvec            force_n;                /* Single force from slab n on one atom          */
2127     rvec            force_n1, force_n2;     /* First and second part of force_n              */
2128     rvec            tmpvec, tmpvec2, tmp_f; /* Helper variables                              */
2129     real            V;                      /* The rotation potential energy                 */
2130     real            OOsigma2;               /* 1/(sigma^2)                                   */
2131     real            beta;                   /* beta_n(xj)                                    */
2132     real            bjn, fit_bjn;           /* b_j^n                                         */
2133     real            gaussian_xj;            /* Gaussian weight gn(xj)                        */
2134     real            betan_xj_sigma2;
2135     real            mj, wj;                 /* Mass-weighting of the positions               */
2136     real            N_M;                    /* N/M                                           */
2137     gmx_enfrotgrp_t erg;                    /* Pointer to enforced rotation group data       */
2138     gmx_bool        bCalcPotFit;
2139
2140
2141     erg = rotg->enfrotgrp;
2142
2143     /* Pre-calculate the inner sums, so that we do not have to calculate
2144      * them again for every atom */
2145     flex_precalc_inner_sum(rotg);
2146
2147     bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT == rotg->eFittype);
2148
2149     /********************************************************/
2150     /* Main loop over all local atoms of the rotation group */
2151     /********************************************************/
2152     OOsigma2 = 1.0/(sigma*sigma);
2153     N_M      = rotg->nat * erg->invmass;
2154     V        = 0.0;
2155     for (j = 0; j < erg->nat_loc; j++)
2156     {
2157         /* Local index of a rotation group atom  */
2158         ii = erg->ind_loc[j];
2159         /* Position of this atom in the collective array */
2160         iigrp = erg->xc_ref_ind[j];
2161         /* Mass-weighting */
2162         mj = erg->mc[iigrp];  /* need the unsorted mass here */
2163         wj = N_M*mj;
2164
2165         /* Current position of this atom: x[ii][XX/YY/ZZ]
2166          * Note that erg->xc_center contains the center of mass in case the flex-t
2167          * potential was chosen. For the flex potential erg->xc_center must be
2168          * zero. */
2169         rvec_sub(x[ii], erg->xc_center, xj);
2170
2171         /* Shift this atom such that it is near its reference */
2172         shift_single_coord(box, xj, erg->xc_shifts[iigrp]);
2173
2174         /* Determine the slabs to loop over, i.e. the ones with contributions
2175          * larger than min_gaussian */
2176         count = get_single_atom_gaussians(xj, rotg);
2177
2178         clear_rvec(sum_n1);
2179         clear_rvec(sum_n2);
2180
2181         /* Loop over the relevant slabs for this atom */
2182         for (ic = 0; ic < count; ic++)
2183         {
2184             n = erg->gn_slabind[ic];
2185
2186             /* Get the precomputed Gaussian for xj in slab n */
2187             gaussian_xj = erg->gn_atom[ic];
2188
2189             islab = n - erg->slab_first; /* slab index */
2190
2191             /* The (unrotated) reference position of this atom is saved in yj0: */
2192             copy_rvec(rotg->x_ref[iigrp], yj0);
2193
2194             beta = calc_beta(xj, rotg, n);
2195
2196             /* The current center of this slab is saved in xcn: */
2197             copy_rvec(erg->slab_center[islab], xcn);
2198             /* ... and the reference center in ycn: */
2199             copy_rvec(erg->slab_center_ref[islab+erg->slab_buffer], ycn);
2200
2201             rvec_sub(yj0, ycn, yj0_ycn); /* yj0_ycn = yj0 - ycn */
2202
2203             /* Rotate: */
2204             mvmul(erg->rotmat, yj0_ycn, tmpvec2); /* tmpvec2= Omega.(yj0-ycn) */
2205
2206             /* Subtract the slab center from xj */
2207             rvec_sub(xj, xcn, xj_xcn);           /* xj_xcn = xj - xcn         */
2208
2209             /* In rare cases, when an atom position coincides with a slab center
2210              * (xj_xcn == 0) we cannot compute the vector product for qjn.
2211              * However, since the atom is located directly on the pivot, this
2212              * slab's contribution to the force on that atom will be zero
2213              * anyway. Therefore, we directly move on to the next slab.       */
2214             if (0 == norm(xj_xcn) )
2215             {
2216                 continue;
2217             }
2218
2219             /* Calculate qjn */
2220             cprod(rotg->vec, tmpvec2, tmpvec); /* tmpvec= v x Omega.(yj0-ycn) */
2221
2222             /*                         *         v x Omega.(yj0-ycn)    */
2223             unitv(tmpvec, qjn);       /*  qjn = ---------------------   */
2224                                       /*        |v x Omega.(yj0-ycn)|   */
2225
2226             bjn = iprod(qjn, xj_xcn); /* bjn = qjn * (xj - xcn) */
2227
2228             /*********************************/
2229             /* Add to the rotation potential */
2230             /*********************************/
2231             V += 0.5*rotg->k*wj*gaussian_xj*sqr(bjn);
2232
2233             /* If requested, also calculate the potential for a set of angles
2234              * near the current reference angle */
2235             if (bCalcPotFit)
2236             {
2237                 for (ifit = 0; ifit < rotg->PotAngle_nstep; ifit++)
2238                 {
2239                     /* As above calculate Omega.(yj0-ycn), now for the other angles */
2240                     mvmul(erg->PotAngleFit->rotmat[ifit], yj0_ycn, tmpvec2); /* tmpvec2= Omega.(yj0-ycn) */
2241                     /* As above calculate qjn */
2242                     cprod(rotg->vec, tmpvec2, tmpvec);                       /* tmpvec= v x Omega.(yj0-ycn) */
2243                     /*                                                        *             v x Omega.(yj0-ycn)    */
2244                     unitv(tmpvec, fit_qjn);                                  /*  fit_qjn = ---------------------   */
2245                                                                              /*            |v x Omega.(yj0-ycn)|   */
2246                     fit_bjn = iprod(fit_qjn, xj_xcn);                        /* fit_bjn = fit_qjn * (xj - xcn) */
2247                     /* Add to the rotation potential for this angle */
2248                     erg->PotAngleFit->V[ifit] += 0.5*rotg->k*wj*gaussian_xj*sqr(fit_bjn);
2249                 }
2250             }
2251
2252             /****************************************************************/
2253             /* sum_n1 will typically be the main contribution to the force: */
2254             /****************************************************************/
2255             betan_xj_sigma2 = beta*OOsigma2;  /*  beta_n(xj)/sigma^2  */
2256
2257             /* The next lines calculate
2258              *  qjn - (bjn*beta(xj)/(2sigma^2))v  */
2259             svmul(bjn*0.5*betan_xj_sigma2, rotg->vec, tmpvec2);
2260             rvec_sub(qjn, tmpvec2, tmpvec);
2261
2262             /* Multiply with gn(xj)*bjn: */
2263             svmul(gaussian_xj*bjn, tmpvec, tmpvec2);
2264
2265             /* Sum over n: */
2266             rvec_inc(sum_n1, tmpvec2);
2267
2268             /* We already have precalculated the Sn term for slab n */
2269             copy_rvec(erg->slab_innersumvec[islab], s_n);
2270             /*                                                             *          beta_n(xj)              */
2271             svmul(betan_xj_sigma2*iprod(s_n, xj_xcn), rotg->vec, tmpvec); /* tmpvec = ---------- s_n (xj-xcn) */
2272                                                                           /*            sigma^2               */
2273
2274             rvec_sub(s_n, tmpvec, innersumvec);
2275
2276             /* We can safely divide by slab_weights since we check in get_slab_centers
2277              * that it is non-zero. */
2278             svmul(gaussian_xj/erg->slab_weights[islab], innersumvec, innersumvec);
2279
2280             rvec_add(sum_n2, innersumvec, sum_n2);
2281
2282             /* Calculate the torque: */
2283             if (bOutstepRot)
2284             {
2285                 /* The force on atom ii from slab n only: */
2286                 svmul(-rotg->k*wj, tmpvec2, force_n1);     /* part 1 */
2287                 svmul( rotg->k*mj, innersumvec, force_n2); /* part 2 */
2288                 rvec_add(force_n1, force_n2, force_n);
2289                 erg->slab_torque_v[islab] += torque(rotg->vec, force_n, xj, xcn);
2290             }
2291         } /* END of loop over slabs */
2292
2293         /* Put both contributions together: */
2294         svmul(wj, sum_n1, sum_n1);
2295         svmul(mj, sum_n2, sum_n2);
2296         rvec_sub(sum_n2, sum_n1, tmp_f); /* F = -grad V */
2297
2298         /* Store the additional force so that it can be added to the force
2299          * array after the normal forces have been evaluated */
2300         for (m = 0; m < DIM; m++)
2301         {
2302             erg->f_rot_loc[j][m] = rotg->k*tmp_f[m];
2303         }
2304
2305         PRINT_FORCE_J
2306
2307     } /* END of loop over local atoms */
2308
2309     return V;
2310 }
2311
2312 #ifdef PRINT_COORDS
2313 static void print_coordinates(t_rotgrp *rotg, rvec x[], matrix box, int step)
2314 {
2315     int             i;
2316     static FILE    *fp;
2317     static char     buf[STRLEN];
2318     static gmx_bool bFirst = 1;
2319
2320
2321     if (bFirst)
2322     {
2323         sprintf(buf, "coords%d.txt", cr->nodeid);
2324         fp     = fopen(buf, "w");
2325         bFirst = 0;
2326     }
2327
2328     fprintf(fp, "\nStep %d\n", step);
2329     fprintf(fp, "box: %f %f %f %f %f %f %f %f %f\n",
2330             box[XX][XX], box[XX][YY], box[XX][ZZ],
2331             box[YY][XX], box[YY][YY], box[YY][ZZ],
2332             box[ZZ][XX], box[ZZ][ZZ], box[ZZ][ZZ]);
2333     for (i = 0; i < rotg->nat; i++)
2334     {
2335         fprintf(fp, "%4d  %f %f %f\n", i,
2336                 erg->xc[i][XX], erg->xc[i][YY], erg->xc[i][ZZ]);
2337     }
2338     fflush(fp);
2339
2340 }
2341 #endif
2342
2343
2344 static int projection_compare(const void *a, const void *b)
2345 {
2346     sort_along_vec_t *xca, *xcb;
2347
2348
2349     xca = (sort_along_vec_t *)a;
2350     xcb = (sort_along_vec_t *)b;
2351
2352     if (xca->xcproj < xcb->xcproj)
2353     {
2354         return -1;
2355     }
2356     else if (xca->xcproj > xcb->xcproj)
2357     {
2358         return 1;
2359     }
2360     else
2361     {
2362         return 0;
2363     }
2364 }
2365
2366
2367 static void sort_collective_coordinates(
2368         t_rotgrp         *rotg, /* Rotation group */
2369         sort_along_vec_t *data) /* Buffer for sorting the positions */
2370 {
2371     int             i;
2372     gmx_enfrotgrp_t erg;       /* Pointer to enforced rotation group data */
2373
2374
2375     erg = rotg->enfrotgrp;
2376
2377     /* The projection of the position vector on the rotation vector is
2378      * the relevant value for sorting. Fill the 'data' structure */
2379     for (i = 0; i < rotg->nat; i++)
2380     {
2381         data[i].xcproj = iprod(erg->xc[i], rotg->vec);  /* sort criterium */
2382         data[i].m      = erg->mc[i];
2383         data[i].ind    = i;
2384         copy_rvec(erg->xc[i], data[i].x    );
2385         copy_rvec(rotg->x_ref[i], data[i].x_ref);
2386     }
2387     /* Sort the 'data' structure */
2388     gmx_qsort(data, rotg->nat, sizeof(sort_along_vec_t), projection_compare);
2389
2390     /* Copy back the sorted values */
2391     for (i = 0; i < rotg->nat; i++)
2392     {
2393         copy_rvec(data[i].x, erg->xc[i]           );
2394         copy_rvec(data[i].x_ref, erg->xc_ref_sorted[i]);
2395         erg->mc_sorted[i]  = data[i].m;
2396         erg->xc_sortind[i] = data[i].ind;
2397     }
2398 }
2399
2400
2401 /* For each slab, get the first and the last index of the sorted atom
2402  * indices */
2403 static void get_firstlast_atom_per_slab(t_rotgrp *rotg)
2404 {
2405     int             i, islab, n;
2406     real            beta;
2407     gmx_enfrotgrp_t erg;     /* Pointer to enforced rotation group data */
2408
2409
2410     erg = rotg->enfrotgrp;
2411
2412     /* Find the first atom that needs to enter the calculation for each slab */
2413     n = erg->slab_first; /* slab */
2414     i = 0;               /* start with the first atom */
2415     do
2416     {
2417         /* Find the first atom that significantly contributes to this slab */
2418         do /* move forward in position until a large enough beta is found */
2419         {
2420             beta = calc_beta(erg->xc[i], rotg, n);
2421             i++;
2422         }
2423         while ((beta < -erg->max_beta) && (i < rotg->nat));
2424         i--;
2425         islab                 = n - erg->slab_first; /* slab index */
2426         erg->firstatom[islab] = i;
2427         /* Proceed to the next slab */
2428         n++;
2429     }
2430     while (n <= erg->slab_last);
2431
2432     /* Find the last atom for each slab */
2433     n = erg->slab_last; /* start with last slab */
2434     i = rotg->nat-1;    /* start with the last atom */
2435     do
2436     {
2437         do  /* move backward in position until a large enough beta is found */
2438         {
2439             beta = calc_beta(erg->xc[i], rotg, n);
2440             i--;
2441         }
2442         while ((beta > erg->max_beta) && (i > -1));
2443         i++;
2444         islab                = n - erg->slab_first; /* slab index */
2445         erg->lastatom[islab] = i;
2446         /* Proceed to the next slab */
2447         n--;
2448     }
2449     while (n >= erg->slab_first);
2450 }
2451
2452
2453 /* Determine the very first and very last slab that needs to be considered
2454  * For the first slab that needs to be considered, we have to find the smallest
2455  * n that obeys:
2456  *
2457  * x_first * v - n*Delta_x <= beta_max
2458  *
2459  * slab index n, slab distance Delta_x, rotation vector v. For the last slab we
2460  * have to find the largest n that obeys
2461  *
2462  * x_last * v - n*Delta_x >= -beta_max
2463  *
2464  */
2465 static gmx_inline int get_first_slab(
2466         t_rotgrp *rotg,      /* The rotation group (inputrec data) */
2467         real      max_beta,  /* The max_beta value, instead of min_gaussian */
2468         rvec      firstatom) /* First atom after sorting along the rotation vector v */
2469 {
2470     /* Find the first slab for the first atom */
2471     return ceil((iprod(firstatom, rotg->vec) - max_beta)/rotg->slab_dist);
2472 }
2473
2474
2475 static gmx_inline int get_last_slab(
2476         t_rotgrp *rotg,     /* The rotation group (inputrec data) */
2477         real      max_beta, /* The max_beta value, instead of min_gaussian */
2478         rvec      lastatom) /* Last atom along v */
2479 {
2480     /* Find the last slab for the last atom */
2481     return floor((iprod(lastatom, rotg->vec) + max_beta)/rotg->slab_dist);
2482 }
2483
2484
2485 static void get_firstlast_slab_check(
2486         t_rotgrp        *rotg,      /* The rotation group (inputrec data) */
2487         t_gmx_enfrotgrp *erg,       /* The rotation group (data only accessible in this file) */
2488         rvec             firstatom, /* First atom after sorting along the rotation vector v */
2489         rvec             lastatom)  /* Last atom along v */
2490 {
2491     erg->slab_first = get_first_slab(rotg, erg->max_beta, firstatom);
2492     erg->slab_last  = get_last_slab(rotg, erg->max_beta, lastatom);
2493
2494     /* Calculate the slab buffer size, which changes when slab_first changes */
2495     erg->slab_buffer = erg->slab_first - erg->slab_first_ref;
2496
2497     /* Check whether we have reference data to compare against */
2498     if (erg->slab_first < erg->slab_first_ref)
2499     {
2500         gmx_fatal(FARGS, "%s No reference data for first slab (n=%d), unable to proceed.",
2501                   RotStr, erg->slab_first);
2502     }
2503
2504     /* Check whether we have reference data to compare against */
2505     if (erg->slab_last > erg->slab_last_ref)
2506     {
2507         gmx_fatal(FARGS, "%s No reference data for last slab (n=%d), unable to proceed.",
2508                   RotStr, erg->slab_last);
2509     }
2510 }
2511
2512
2513 /* Enforced rotation with a flexible axis */
2514 static void do_flexible(
2515         gmx_bool        bMaster,
2516         gmx_enfrot_t    enfrot,       /* Other rotation data                        */
2517         t_rotgrp       *rotg,         /* The rotation group                         */
2518         int             g,            /* Group number                               */
2519         rvec            x[],          /* The local positions                        */
2520         matrix          box,
2521         double          t,            /* Time in picoseconds                        */
2522         gmx_bool        bOutstepRot,  /* Output to main rotation output file        */
2523         gmx_bool        bOutstepSlab) /* Output per-slab data                       */
2524 {
2525     int             l, nslabs;
2526     real            sigma;    /* The Gaussian width sigma */
2527     gmx_enfrotgrp_t erg;      /* Pointer to enforced rotation group data */
2528
2529
2530     erg = rotg->enfrotgrp;
2531
2532     /* Define the sigma value */
2533     sigma = 0.7*rotg->slab_dist;
2534
2535     /* Sort the collective coordinates erg->xc along the rotation vector. This is
2536      * an optimization for the inner loop. */
2537     sort_collective_coordinates(rotg, enfrot->data);
2538
2539     /* Determine the first relevant slab for the first atom and the last
2540      * relevant slab for the last atom */
2541     get_firstlast_slab_check(rotg, erg, erg->xc[0], erg->xc[rotg->nat-1]);
2542
2543     /* Determine for each slab depending on the min_gaussian cutoff criterium,
2544      * a first and a last atom index inbetween stuff needs to be calculated */
2545     get_firstlast_atom_per_slab(rotg);
2546
2547     /* Determine the gaussian-weighted center of positions for all slabs */
2548     get_slab_centers(rotg, erg->xc, erg->mc_sorted, g, t, enfrot->out_slabs, bOutstepSlab, FALSE);
2549
2550     /* Clear the torque per slab from last time step: */
2551     nslabs = erg->slab_last - erg->slab_first + 1;
2552     for (l = 0; l < nslabs; l++)
2553     {
2554         erg->slab_torque_v[l] = 0.0;
2555     }
2556
2557     /* Call the rotational forces kernel */
2558     if (rotg->eType == erotgFLEX || rotg->eType == erotgFLEXT)
2559     {
2560         erg->V = do_flex_lowlevel(rotg, sigma, x, bOutstepRot, bOutstepSlab, box);
2561     }
2562     else if (rotg->eType == erotgFLEX2 || rotg->eType == erotgFLEX2T)
2563     {
2564         erg->V = do_flex2_lowlevel(rotg, sigma, x, bOutstepRot, bOutstepSlab, box);
2565     }
2566     else
2567     {
2568         gmx_fatal(FARGS, "Unknown flexible rotation type");
2569     }
2570
2571     /* Determine angle by RMSD fit to the reference - Let's hope this */
2572     /* only happens once in a while, since this is not parallelized! */
2573     if (bMaster && (erotgFitPOT != rotg->eFittype) )
2574     {
2575         if (bOutstepRot)
2576         {
2577             /* Fit angle of the whole rotation group */
2578             erg->angle_v = flex_fit_angle(rotg);
2579         }
2580         if (bOutstepSlab)
2581         {
2582             /* Fit angle of each slab */
2583             flex_fit_angle_perslab(g, rotg, t, erg->degangle, enfrot->out_angles);
2584         }
2585     }
2586
2587     /* Lump together the torques from all slabs: */
2588     erg->torque_v = 0.0;
2589     for (l = 0; l < nslabs; l++)
2590     {
2591         erg->torque_v += erg->slab_torque_v[l];
2592     }
2593 }
2594
2595
2596 /* Calculate the angle between reference and actual rotation group atom,
2597  * both projected into a plane perpendicular to the rotation vector: */
2598 static void angle(t_rotgrp *rotg,
2599                   rvec      x_act,
2600                   rvec      x_ref,
2601                   real     *alpha,
2602                   real     *weight) /* atoms near the rotation axis should count less than atoms far away */
2603 {
2604     rvec xp, xrp;                   /* current and reference positions projected on a plane perpendicular to pg->vec */
2605     rvec dum;
2606
2607
2608     /* Project x_ref and x into a plane through the origin perpendicular to rot_vec: */
2609     /* Project x_ref: xrp = x_ref - (vec * x_ref) * vec */
2610     svmul(iprod(rotg->vec, x_ref), rotg->vec, dum);
2611     rvec_sub(x_ref, dum, xrp);
2612     /* Project x_act: */
2613     svmul(iprod(rotg->vec, x_act), rotg->vec, dum);
2614     rvec_sub(x_act, dum, xp);
2615
2616     /* Retrieve information about which vector precedes. gmx_angle always
2617      * returns a positive angle. */
2618     cprod(xp, xrp, dum); /* if reference precedes, this is pointing into the same direction as vec */
2619
2620     if (iprod(rotg->vec, dum) >= 0)
2621     {
2622         *alpha = -gmx_angle(xrp, xp);
2623     }
2624     else
2625     {
2626         *alpha = +gmx_angle(xrp, xp);
2627     }
2628
2629     /* Also return the weight */
2630     *weight = norm(xp);
2631 }
2632
2633
2634 /* Project first vector onto a plane perpendicular to the second vector
2635  * dr = dr - (dr.v)v
2636  * Note that v must be of unit length.
2637  */
2638 static gmx_inline void project_onto_plane(rvec dr, const rvec v)
2639 {
2640     rvec tmp;
2641
2642
2643     svmul(iprod(dr, v), v, tmp); /* tmp = (dr.v)v */
2644     rvec_dec(dr, tmp);           /* dr = dr - (dr.v)v */
2645 }
2646
2647
2648 /* Fixed rotation: The rotation reference group rotates around the v axis. */
2649 /* The atoms of the actual rotation group are attached with imaginary  */
2650 /* springs to the reference atoms.                                     */
2651 static void do_fixed(
2652         t_rotgrp       *rotg,         /* The rotation group                         */
2653         gmx_bool        bOutstepRot,  /* Output to main rotation output file        */
2654         gmx_bool        bOutstepSlab) /* Output per-slab data                       */
2655 {
2656     int             ifit, j, jj, m;
2657     rvec            dr;
2658     rvec            tmp_f;     /* Force */
2659     real            alpha;     /* a single angle between an actual and a reference position */
2660     real            weight;    /* single weight for a single angle */
2661     gmx_enfrotgrp_t erg;       /* Pointer to enforced rotation group data */
2662     rvec            xi_xc;     /* xi - xc */
2663     gmx_bool        bCalcPotFit;
2664     rvec            fit_xr_loc;
2665
2666     /* for mass weighting: */
2667     real      wi;              /* Mass-weighting of the positions */
2668     real      N_M;             /* N/M */
2669     real      k_wi;            /* k times wi */
2670
2671     gmx_bool  bProject;
2672
2673
2674     erg         = rotg->enfrotgrp;
2675     bProject    = (rotg->eType == erotgPM) || (rotg->eType == erotgPMPF);
2676     bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT == rotg->eFittype);
2677
2678     N_M = rotg->nat * erg->invmass;
2679
2680     /* Each process calculates the forces on its local atoms */
2681     for (j = 0; j < erg->nat_loc; j++)
2682     {
2683         /* Calculate (x_i-x_c) resp. (x_i-u) */
2684         rvec_sub(erg->x_loc_pbc[j], erg->xc_center, xi_xc);
2685
2686         /* Calculate Omega*(y_i-y_c)-(x_i-x_c) */
2687         rvec_sub(erg->xr_loc[j], xi_xc, dr);
2688
2689         if (bProject)
2690         {
2691             project_onto_plane(dr, rotg->vec);
2692         }
2693
2694         /* Mass-weighting */
2695         wi = N_M*erg->m_loc[j];
2696
2697         /* Store the additional force so that it can be added to the force
2698          * array after the normal forces have been evaluated */
2699         k_wi = rotg->k*wi;
2700         for (m = 0; m < DIM; m++)
2701         {
2702             tmp_f[m]             = k_wi*dr[m];
2703             erg->f_rot_loc[j][m] = tmp_f[m];
2704             erg->V              += 0.5*k_wi*sqr(dr[m]);
2705         }
2706
2707         /* If requested, also calculate the potential for a set of angles
2708          * near the current reference angle */
2709         if (bCalcPotFit)
2710         {
2711             for (ifit = 0; ifit < rotg->PotAngle_nstep; ifit++)
2712             {
2713                 /* Index of this rotation group atom with respect to the whole rotation group */
2714                 jj = erg->xc_ref_ind[j];
2715
2716                 /* Rotate with the alternative angle. Like rotate_local_reference(),
2717                  * just for a single local atom */
2718                 mvmul(erg->PotAngleFit->rotmat[ifit], rotg->x_ref[jj], fit_xr_loc); /* fit_xr_loc = Omega*(y_i-y_c) */
2719
2720                 /* Calculate Omega*(y_i-y_c)-(x_i-x_c) */
2721                 rvec_sub(fit_xr_loc, xi_xc, dr);
2722
2723                 if (bProject)
2724                 {
2725                     project_onto_plane(dr, rotg->vec);
2726                 }
2727
2728                 /* Add to the rotation potential for this angle: */
2729                 erg->PotAngleFit->V[ifit] += 0.5*k_wi*norm2(dr);
2730             }
2731         }
2732
2733         if (bOutstepRot)
2734         {
2735             /* Add to the torque of this rotation group */
2736             erg->torque_v += torque(rotg->vec, tmp_f, erg->x_loc_pbc[j], erg->xc_center);
2737
2738             /* Calculate the angle between reference and actual rotation group atom. */
2739             angle(rotg, xi_xc, erg->xr_loc[j], &alpha, &weight);  /* angle in rad, weighted */
2740             erg->angle_v  += alpha * weight;
2741             erg->weight_v += weight;
2742         }
2743         /* If you want enforced rotation to contribute to the virial,
2744          * activate the following lines:
2745             if (MASTER(cr))
2746             {
2747                Add the rotation contribution to the virial
2748               for(j=0; j<DIM; j++)
2749                 for(m=0;m<DIM;m++)
2750                   vir[j][m] += 0.5*f[ii][j]*dr[m];
2751             }
2752          */
2753
2754         PRINT_FORCE_J
2755
2756     } /* end of loop over local rotation group atoms */
2757 }
2758
2759
2760 /* Calculate the radial motion potential and forces */
2761 static void do_radial_motion(
2762         t_rotgrp       *rotg,         /* The rotation group                         */
2763         gmx_bool        bOutstepRot,  /* Output to main rotation output file        */
2764         gmx_bool        bOutstepSlab) /* Output per-slab data                       */
2765 {
2766     int             j, jj, ifit;
2767     rvec            tmp_f;     /* Force */
2768     real            alpha;     /* a single angle between an actual and a reference position */
2769     real            weight;    /* single weight for a single angle */
2770     gmx_enfrotgrp_t erg;       /* Pointer to enforced rotation group data */
2771     rvec            xj_u;      /* xj - u */
2772     rvec            tmpvec, fit_tmpvec;
2773     real            fac, fac2, sum = 0.0;
2774     rvec            pj;
2775     gmx_bool        bCalcPotFit;
2776
2777     /* For mass weighting: */
2778     real      wj;              /* Mass-weighting of the positions */
2779     real      N_M;             /* N/M */
2780
2781
2782     erg         = rotg->enfrotgrp;
2783     bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT == rotg->eFittype);
2784
2785     N_M = rotg->nat * erg->invmass;
2786
2787     /* Each process calculates the forces on its local atoms */
2788     for (j = 0; j < erg->nat_loc; j++)
2789     {
2790         /* Calculate (xj-u) */
2791         rvec_sub(erg->x_loc_pbc[j], erg->xc_center, xj_u);  /* xj_u = xj-u */
2792
2793         /* Calculate Omega.(yj0-u) */
2794         cprod(rotg->vec, erg->xr_loc[j], tmpvec);  /* tmpvec = v x Omega.(yj0-u) */
2795
2796         /*                       *         v x Omega.(yj0-u)     */
2797         unitv(tmpvec, pj);      /*  pj = ---------------------   */
2798                                 /*       | v x Omega.(yj0-u) |   */
2799
2800         fac  = iprod(pj, xj_u); /* fac = pj.(xj-u) */
2801         fac2 = fac*fac;
2802
2803         /* Mass-weighting */
2804         wj = N_M*erg->m_loc[j];
2805
2806         /* Store the additional force so that it can be added to the force
2807          * array after the normal forces have been evaluated */
2808         svmul(-rotg->k*wj*fac, pj, tmp_f);
2809         copy_rvec(tmp_f, erg->f_rot_loc[j]);
2810         sum += wj*fac2;
2811
2812         /* If requested, also calculate the potential for a set of angles
2813          * near the current reference angle */
2814         if (bCalcPotFit)
2815         {
2816             for (ifit = 0; ifit < rotg->PotAngle_nstep; ifit++)
2817             {
2818                 /* Index of this rotation group atom with respect to the whole rotation group */
2819                 jj = erg->xc_ref_ind[j];
2820
2821                 /* Rotate with the alternative angle. Like rotate_local_reference(),
2822                  * just for a single local atom */
2823                 mvmul(erg->PotAngleFit->rotmat[ifit], rotg->x_ref[jj], fit_tmpvec); /* fit_tmpvec = Omega*(yj0-u) */
2824
2825                 /* Calculate Omega.(yj0-u) */
2826                 cprod(rotg->vec, fit_tmpvec, tmpvec); /* tmpvec = v x Omega.(yj0-u) */
2827                 /*                                     *         v x Omega.(yj0-u)     */
2828                 unitv(tmpvec, pj);                    /*  pj = ---------------------   */
2829                                                       /*       | v x Omega.(yj0-u) |   */
2830
2831                 fac  = iprod(pj, xj_u);               /* fac = pj.(xj-u) */
2832                 fac2 = fac*fac;
2833
2834                 /* Add to the rotation potential for this angle: */
2835                 erg->PotAngleFit->V[ifit] += 0.5*rotg->k*wj*fac2;
2836             }
2837         }
2838
2839         if (bOutstepRot)
2840         {
2841             /* Add to the torque of this rotation group */
2842             erg->torque_v += torque(rotg->vec, tmp_f, erg->x_loc_pbc[j], erg->xc_center);
2843
2844             /* Calculate the angle between reference and actual rotation group atom. */
2845             angle(rotg, xj_u, erg->xr_loc[j], &alpha, &weight);  /* angle in rad, weighted */
2846             erg->angle_v  += alpha * weight;
2847             erg->weight_v += weight;
2848         }
2849
2850         PRINT_FORCE_J
2851
2852     } /* end of loop over local rotation group atoms */
2853     erg->V = 0.5*rotg->k*sum;
2854 }
2855
2856
2857 /* Calculate the radial motion pivot-free potential and forces */
2858 static void do_radial_motion_pf(
2859         t_rotgrp       *rotg,         /* The rotation group                         */
2860         rvec            x[],          /* The positions                              */
2861         matrix          box,          /* The simulation box                         */
2862         gmx_bool        bOutstepRot,  /* Output to main rotation output file        */
2863         gmx_bool        bOutstepSlab) /* Output per-slab data                       */
2864 {
2865     int             i, ii, iigrp, ifit, j;
2866     rvec            xj;          /* Current position */
2867     rvec            xj_xc;       /* xj  - xc  */
2868     rvec            yj0_yc0;     /* yj0 - yc0 */
2869     rvec            tmp_f;       /* Force */
2870     real            alpha;       /* a single angle between an actual and a reference position */
2871     real            weight;      /* single weight for a single angle */
2872     gmx_enfrotgrp_t erg;         /* Pointer to enforced rotation group data */
2873     rvec            tmpvec, tmpvec2;
2874     rvec            innersumvec; /* Precalculation of the inner sum */
2875     rvec            innersumveckM;
2876     real            fac, fac2, V = 0.0;
2877     rvec            qi, qj;
2878     gmx_bool        bCalcPotFit;
2879
2880     /* For mass weighting: */
2881     real      mj, wi, wj;      /* Mass-weighting of the positions */
2882     real      N_M;             /* N/M */
2883
2884
2885     erg         = rotg->enfrotgrp;
2886     bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT == rotg->eFittype);
2887
2888     N_M = rotg->nat * erg->invmass;
2889
2890     /* Get the current center of the rotation group: */
2891     get_center(erg->xc, erg->mc, rotg->nat, erg->xc_center);
2892
2893     /* Precalculate Sum_i [ wi qi.(xi-xc) qi ] which is needed for every single j */
2894     clear_rvec(innersumvec);
2895     for (i = 0; i < rotg->nat; i++)
2896     {
2897         /* Mass-weighting */
2898         wi = N_M*erg->mc[i];
2899
2900         /* Calculate qi. Note that xc_ref_center has already been subtracted from
2901          * x_ref in init_rot_group.*/
2902         mvmul(erg->rotmat, rotg->x_ref[i], tmpvec); /* tmpvec  = Omega.(yi0-yc0) */
2903
2904         cprod(rotg->vec, tmpvec, tmpvec2);          /* tmpvec2 = v x Omega.(yi0-yc0) */
2905
2906         /*                                             *         v x Omega.(yi0-yc0)     */
2907         unitv(tmpvec2, qi);                           /*  qi = -----------------------   */
2908                                                       /*       | v x Omega.(yi0-yc0) |   */
2909
2910         rvec_sub(erg->xc[i], erg->xc_center, tmpvec); /* tmpvec = xi-xc */
2911
2912         svmul(wi*iprod(qi, tmpvec), qi, tmpvec2);
2913
2914         rvec_inc(innersumvec, tmpvec2);
2915     }
2916     svmul(rotg->k*erg->invmass, innersumvec, innersumveckM);
2917
2918     /* Each process calculates the forces on its local atoms */
2919     for (j = 0; j < erg->nat_loc; j++)
2920     {
2921         /* Local index of a rotation group atom  */
2922         ii = erg->ind_loc[j];
2923         /* Position of this atom in the collective array */
2924         iigrp = erg->xc_ref_ind[j];
2925         /* Mass-weighting */
2926         mj = erg->mc[iigrp];  /* need the unsorted mass here */
2927         wj = N_M*mj;
2928
2929         /* Current position of this atom: x[ii][XX/YY/ZZ] */
2930         copy_rvec(x[ii], xj);
2931
2932         /* Shift this atom such that it is near its reference */
2933         shift_single_coord(box, xj, erg->xc_shifts[iigrp]);
2934
2935         /* The (unrotated) reference position is yj0. yc0 has already
2936          * been subtracted in init_rot_group */
2937         copy_rvec(rotg->x_ref[iigrp], yj0_yc0);   /* yj0_yc0 = yj0 - yc0      */
2938
2939         /* Calculate Omega.(yj0-yc0) */
2940         mvmul(erg->rotmat, yj0_yc0, tmpvec2); /* tmpvec2 = Omega.(yj0 - yc0)  */
2941
2942         cprod(rotg->vec, tmpvec2, tmpvec);    /* tmpvec = v x Omega.(yj0-yc0) */
2943
2944         /*                     *         v x Omega.(yj0-yc0)     */
2945         unitv(tmpvec, qj);    /*  qj = -----------------------   */
2946                               /*       | v x Omega.(yj0-yc0) |   */
2947
2948         /* Calculate (xj-xc) */
2949         rvec_sub(xj, erg->xc_center, xj_xc); /* xj_xc = xj-xc */
2950
2951         fac  = iprod(qj, xj_xc);             /* fac = qj.(xj-xc) */
2952         fac2 = fac*fac;
2953
2954         /* Store the additional force so that it can be added to the force
2955          * array after the normal forces have been evaluated */
2956         svmul(-rotg->k*wj*fac, qj, tmp_f); /* part 1 of force */
2957         svmul(mj, innersumveckM, tmpvec);  /* part 2 of force */
2958         rvec_inc(tmp_f, tmpvec);
2959         copy_rvec(tmp_f, erg->f_rot_loc[j]);
2960         V += wj*fac2;
2961
2962         /* If requested, also calculate the potential for a set of angles
2963          * near the current reference angle */
2964         if (bCalcPotFit)
2965         {
2966             for (ifit = 0; ifit < rotg->PotAngle_nstep; ifit++)
2967             {
2968                 /* Rotate with the alternative angle. Like rotate_local_reference(),
2969                  * just for a single local atom */
2970                 mvmul(erg->PotAngleFit->rotmat[ifit], yj0_yc0, tmpvec2); /* tmpvec2 = Omega*(yj0-yc0) */
2971
2972                 /* Calculate Omega.(yj0-u) */
2973                 cprod(rotg->vec, tmpvec2, tmpvec); /* tmpvec = v x Omega.(yj0-yc0) */
2974                 /*                                  *         v x Omega.(yj0-yc0)     */
2975                 unitv(tmpvec, qj);                 /*  qj = -----------------------   */
2976                                                    /*       | v x Omega.(yj0-yc0) |   */
2977
2978                 fac  = iprod(qj, xj_xc);           /* fac = qj.(xj-xc) */
2979                 fac2 = fac*fac;
2980
2981                 /* Add to the rotation potential for this angle: */
2982                 erg->PotAngleFit->V[ifit] += 0.5*rotg->k*wj*fac2;
2983             }
2984         }
2985
2986         if (bOutstepRot)
2987         {
2988             /* Add to the torque of this rotation group */
2989             erg->torque_v += torque(rotg->vec, tmp_f, xj, erg->xc_center);
2990
2991             /* Calculate the angle between reference and actual rotation group atom. */
2992             angle(rotg, xj_xc, yj0_yc0, &alpha, &weight);  /* angle in rad, weighted */
2993             erg->angle_v  += alpha * weight;
2994             erg->weight_v += weight;
2995         }
2996
2997         PRINT_FORCE_J
2998
2999     } /* end of loop over local rotation group atoms */
3000     erg->V = 0.5*rotg->k*V;
3001 }
3002
3003
3004 /* Precalculate the inner sum for the radial motion 2 forces */
3005 static void radial_motion2_precalc_inner_sum(t_rotgrp  *rotg, rvec innersumvec)
3006 {
3007     int             i;
3008     gmx_enfrotgrp_t erg;       /* Pointer to enforced rotation group data */
3009     rvec            xi_xc;     /* xj - xc */
3010     rvec            tmpvec, tmpvec2;
3011     real            fac, fac2;
3012     rvec            ri, si;
3013     real            siri;
3014     rvec            v_xi_xc;   /* v x (xj - u) */
3015     real            psii, psiistar;
3016     real            wi;        /* Mass-weighting of the positions */
3017     real            N_M;       /* N/M */
3018     rvec            sumvec;
3019
3020     erg = rotg->enfrotgrp;
3021     N_M = rotg->nat * erg->invmass;
3022
3023     /* Loop over the collective set of positions */
3024     clear_rvec(sumvec);
3025     for (i = 0; i < rotg->nat; i++)
3026     {
3027         /* Mass-weighting */
3028         wi = N_M*erg->mc[i];
3029
3030         rvec_sub(erg->xc[i], erg->xc_center, xi_xc); /* xi_xc = xi-xc         */
3031
3032         /* Calculate ri. Note that xc_ref_center has already been subtracted from
3033          * x_ref in init_rot_group.*/
3034         mvmul(erg->rotmat, rotg->x_ref[i], ri);      /* ri  = Omega.(yi0-yc0) */
3035
3036         cprod(rotg->vec, xi_xc, v_xi_xc);            /* v_xi_xc = v x (xi-u)  */
3037
3038         fac = norm2(v_xi_xc);
3039         /*                                 *                      1           */
3040         psiistar = 1.0/(fac + rotg->eps); /* psiistar = --------------------- */
3041                                           /*            |v x (xi-xc)|^2 + eps */
3042
3043         psii = gmx_invsqrt(fac);          /*                 1                */
3044                                           /*  psii    = -------------         */
3045                                           /*            |v x (xi-xc)|         */
3046
3047         svmul(psii, v_xi_xc, si);         /*  si = psii * (v x (xi-xc) )     */
3048
3049         fac  = iprod(v_xi_xc, ri);        /* fac = (v x (xi-xc)).ri */
3050         fac2 = fac*fac;
3051
3052         siri = iprod(si, ri);                       /* siri = si.ri           */
3053
3054         svmul(psiistar/psii, ri, tmpvec);
3055         svmul(psiistar*psiistar/(psii*psii*psii) * siri, si, tmpvec2);
3056         rvec_dec(tmpvec, tmpvec2);
3057         cprod(tmpvec, rotg->vec, tmpvec2);
3058
3059         svmul(wi*siri, tmpvec2, tmpvec);
3060
3061         rvec_inc(sumvec, tmpvec);
3062     }
3063     svmul(rotg->k*erg->invmass, sumvec, innersumvec);
3064 }
3065
3066
3067 /* Calculate the radial motion 2 potential and forces */
3068 static void do_radial_motion2(
3069         t_rotgrp       *rotg,         /* The rotation group                         */
3070         rvec            x[],          /* The positions                              */
3071         matrix          box,          /* The simulation box                         */
3072         gmx_bool        bOutstepRot,  /* Output to main rotation output file        */
3073         gmx_bool        bOutstepSlab) /* Output per-slab data                       */
3074 {
3075     int             ii, iigrp, ifit, j;
3076     rvec            xj;        /* Position */
3077     real            alpha;     /* a single angle between an actual and a reference position */
3078     real            weight;    /* single weight for a single angle */
3079     gmx_enfrotgrp_t erg;       /* Pointer to enforced rotation group data */
3080     rvec            xj_u;      /* xj - u */
3081     rvec            yj0_yc0;   /* yj0 -yc0 */
3082     rvec            tmpvec, tmpvec2;
3083     real            fac, fit_fac, fac2, Vpart = 0.0;
3084     rvec            rj, fit_rj, sj;
3085     real            sjrj;
3086     rvec            v_xj_u;    /* v x (xj - u) */
3087     real            psij, psijstar;
3088     real            mj, wj;    /* For mass-weighting of the positions */
3089     real            N_M;       /* N/M */
3090     gmx_bool        bPF;
3091     rvec            innersumvec;
3092     gmx_bool        bCalcPotFit;
3093
3094
3095     erg = rotg->enfrotgrp;
3096
3097     bPF         = rotg->eType == erotgRM2PF;
3098     bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT == rotg->eFittype);
3099
3100
3101     clear_rvec(yj0_yc0); /* Make the compiler happy */
3102
3103     clear_rvec(innersumvec);
3104     if (bPF)
3105     {
3106         /* For the pivot-free variant we have to use the current center of
3107          * mass of the rotation group instead of the pivot u */
3108         get_center(erg->xc, erg->mc, rotg->nat, erg->xc_center);
3109
3110         /* Also, we precalculate the second term of the forces that is identical
3111          * (up to the weight factor mj) for all forces */
3112         radial_motion2_precalc_inner_sum(rotg, innersumvec);
3113     }
3114
3115     N_M = rotg->nat * erg->invmass;
3116
3117     /* Each process calculates the forces on its local atoms */
3118     for (j = 0; j < erg->nat_loc; j++)
3119     {
3120         if (bPF)
3121         {
3122             /* Local index of a rotation group atom  */
3123             ii = erg->ind_loc[j];
3124             /* Position of this atom in the collective array */
3125             iigrp = erg->xc_ref_ind[j];
3126             /* Mass-weighting */
3127             mj = erg->mc[iigrp];
3128
3129             /* Current position of this atom: x[ii] */
3130             copy_rvec(x[ii], xj);
3131
3132             /* Shift this atom such that it is near its reference */
3133             shift_single_coord(box, xj, erg->xc_shifts[iigrp]);
3134
3135             /* The (unrotated) reference position is yj0. yc0 has already
3136              * been subtracted in init_rot_group */
3137             copy_rvec(rotg->x_ref[iigrp], yj0_yc0);   /* yj0_yc0 = yj0 - yc0  */
3138
3139             /* Calculate Omega.(yj0-yc0) */
3140             mvmul(erg->rotmat, yj0_yc0, rj);         /* rj = Omega.(yj0-yc0)  */
3141         }
3142         else
3143         {
3144             mj = erg->m_loc[j];
3145             copy_rvec(erg->x_loc_pbc[j], xj);
3146             copy_rvec(erg->xr_loc[j], rj);           /* rj = Omega.(yj0-u)    */
3147         }
3148         /* Mass-weighting */
3149         wj = N_M*mj;
3150
3151         /* Calculate (xj-u) resp. (xj-xc) */
3152         rvec_sub(xj, erg->xc_center, xj_u);          /* xj_u = xj-u           */
3153
3154         cprod(rotg->vec, xj_u, v_xj_u);              /* v_xj_u = v x (xj-u)   */
3155
3156         fac = norm2(v_xj_u);
3157         /*                                 *                      1           */
3158         psijstar = 1.0/(fac + rotg->eps); /*  psistar = --------------------  */
3159                                           /*            |v x (xj-u)|^2 + eps  */
3160
3161         psij = gmx_invsqrt(fac);          /*                 1                */
3162                                           /*  psij    = ------------          */
3163                                           /*            |v x (xj-u)|          */
3164
3165         svmul(psij, v_xj_u, sj);          /*  sj = psij * (v x (xj-u) )       */
3166
3167         fac  = iprod(v_xj_u, rj);         /* fac = (v x (xj-u)).rj */
3168         fac2 = fac*fac;
3169
3170         sjrj = iprod(sj, rj);                        /* sjrj = sj.rj          */
3171
3172         svmul(psijstar/psij, rj, tmpvec);
3173         svmul(psijstar*psijstar/(psij*psij*psij) * sjrj, sj, tmpvec2);
3174         rvec_dec(tmpvec, tmpvec2);
3175         cprod(tmpvec, rotg->vec, tmpvec2);
3176
3177         /* Store the additional force so that it can be added to the force
3178          * array after the normal forces have been evaluated */
3179         svmul(-rotg->k*wj*sjrj, tmpvec2, tmpvec);
3180         svmul(mj, innersumvec, tmpvec2);  /* This is != 0 only for the pivot-free variant */
3181
3182         rvec_add(tmpvec2, tmpvec, erg->f_rot_loc[j]);
3183         Vpart += wj*psijstar*fac2;
3184
3185         /* If requested, also calculate the potential for a set of angles
3186          * near the current reference angle */
3187         if (bCalcPotFit)
3188         {
3189             for (ifit = 0; ifit < rotg->PotAngle_nstep; ifit++)
3190             {
3191                 if (bPF)
3192                 {
3193                     mvmul(erg->PotAngleFit->rotmat[ifit], yj0_yc0, fit_rj); /* fit_rj = Omega.(yj0-yc0) */
3194                 }
3195                 else
3196                 {
3197                     /* Position of this atom in the collective array */
3198                     iigrp = erg->xc_ref_ind[j];
3199                     /* Rotate with the alternative angle. Like rotate_local_reference(),
3200                      * just for a single local atom */
3201                     mvmul(erg->PotAngleFit->rotmat[ifit], rotg->x_ref[iigrp], fit_rj); /* fit_rj = Omega*(yj0-u) */
3202                 }
3203                 fit_fac = iprod(v_xj_u, fit_rj);                                       /* fac = (v x (xj-u)).fit_rj */
3204                 /* Add to the rotation potential for this angle: */
3205                 erg->PotAngleFit->V[ifit] += 0.5*rotg->k*wj*psijstar*fit_fac*fit_fac;
3206             }
3207         }
3208
3209         if (bOutstepRot)
3210         {
3211             /* Add to the torque of this rotation group */
3212             erg->torque_v += torque(rotg->vec, erg->f_rot_loc[j], xj, erg->xc_center);
3213
3214             /* Calculate the angle between reference and actual rotation group atom. */
3215             angle(rotg, xj_u, rj, &alpha, &weight);  /* angle in rad, weighted */
3216             erg->angle_v  += alpha * weight;
3217             erg->weight_v += weight;
3218         }
3219
3220         PRINT_FORCE_J
3221
3222     } /* end of loop over local rotation group atoms */
3223     erg->V = 0.5*rotg->k*Vpart;
3224 }
3225
3226
3227 /* Determine the smallest and largest position vector (with respect to the
3228  * rotation vector) for the reference group */
3229 static void get_firstlast_atom_ref(
3230         t_rotgrp  *rotg,
3231         int       *firstindex,
3232         int       *lastindex)
3233 {
3234     gmx_enfrotgrp_t erg;              /* Pointer to enforced rotation group data */
3235     int             i;
3236     real            xcproj;           /* The projection of a reference position on the
3237                                          rotation vector */
3238     real            minproj, maxproj; /* Smallest and largest projection on v */
3239
3240
3241
3242     erg = rotg->enfrotgrp;
3243
3244     /* Start with some value */
3245     minproj = iprod(rotg->x_ref[0], rotg->vec);
3246     maxproj = minproj;
3247
3248     /* This is just to ensure that it still works if all the atoms of the
3249      * reference structure are situated in a plane perpendicular to the rotation
3250      * vector */
3251     *firstindex = 0;
3252     *lastindex  = rotg->nat-1;
3253
3254     /* Loop over all atoms of the reference group,
3255      * project them on the rotation vector to find the extremes */
3256     for (i = 0; i < rotg->nat; i++)
3257     {
3258         xcproj = iprod(rotg->x_ref[i], rotg->vec);
3259         if (xcproj < minproj)
3260         {
3261             minproj     = xcproj;
3262             *firstindex = i;
3263         }
3264         if (xcproj > maxproj)
3265         {
3266             maxproj    = xcproj;
3267             *lastindex = i;
3268         }
3269     }
3270 }
3271
3272
3273 /* Allocate memory for the slabs */
3274 static void allocate_slabs(
3275         t_rotgrp  *rotg,
3276         FILE      *fplog,
3277         int        g,
3278         gmx_bool   bVerbose)
3279 {
3280     gmx_enfrotgrp_t erg;      /* Pointer to enforced rotation group data */
3281     int             i, nslabs;
3282
3283
3284     erg = rotg->enfrotgrp;
3285
3286     /* More slabs than are defined for the reference are never needed */
3287     nslabs = erg->slab_last_ref - erg->slab_first_ref + 1;
3288
3289     /* Remember how many we allocated */
3290     erg->nslabs_alloc = nslabs;
3291
3292     if ( (NULL != fplog) && bVerbose)
3293     {
3294         fprintf(fplog, "%s allocating memory to store data for %d slabs (rotation group %d).\n",
3295                 RotStr, nslabs, g);
3296     }
3297     snew(erg->slab_center, nslabs);
3298     snew(erg->slab_center_ref, nslabs);
3299     snew(erg->slab_weights, nslabs);
3300     snew(erg->slab_torque_v, nslabs);
3301     snew(erg->slab_data, nslabs);
3302     snew(erg->gn_atom, nslabs);
3303     snew(erg->gn_slabind, nslabs);
3304     snew(erg->slab_innersumvec, nslabs);
3305     for (i = 0; i < nslabs; i++)
3306     {
3307         snew(erg->slab_data[i].x, rotg->nat);
3308         snew(erg->slab_data[i].ref, rotg->nat);
3309         snew(erg->slab_data[i].weight, rotg->nat);
3310     }
3311     snew(erg->xc_ref_sorted, rotg->nat);
3312     snew(erg->xc_sortind, rotg->nat);
3313     snew(erg->firstatom, nslabs);
3314     snew(erg->lastatom, nslabs);
3315 }
3316
3317
3318 /* From the extreme positions of the reference group, determine the first
3319  * and last slab of the reference. We can never have more slabs in the real
3320  * simulation than calculated here for the reference.
3321  */
3322 static void get_firstlast_slab_ref(t_rotgrp *rotg, real mc[], int ref_firstindex, int ref_lastindex)
3323 {
3324     gmx_enfrotgrp_t erg;      /* Pointer to enforced rotation group data */
3325     int             first, last;
3326     rvec            dummy;
3327
3328
3329     erg        = rotg->enfrotgrp;
3330     first      = get_first_slab(rotg, erg->max_beta, rotg->x_ref[ref_firstindex]);
3331     last       = get_last_slab( rotg, erg->max_beta, rotg->x_ref[ref_lastindex ]);
3332
3333     while (get_slab_weight(first, rotg, rotg->x_ref, mc, &dummy) > WEIGHT_MIN)
3334     {
3335         first--;
3336     }
3337     erg->slab_first_ref = first+1;
3338     while (get_slab_weight(last, rotg, rotg->x_ref, mc, &dummy) > WEIGHT_MIN)
3339     {
3340         last++;
3341     }
3342     erg->slab_last_ref  = last-1;
3343 }
3344
3345
3346 /* Special version of copy_rvec:
3347  * During the copy procedure of xcurr to b, the correct PBC image is chosen
3348  * such that the copied vector ends up near its reference position xref */
3349 static inline void copy_correct_pbc_image(
3350         const rvec  xcurr,  /* copy vector xcurr ...                */
3351         rvec        b,      /* ... to b ...                         */
3352         const rvec  xref,   /* choosing the PBC image such that b ends up near xref */
3353         matrix      box,
3354         int         npbcdim)
3355 {
3356     rvec  dx;
3357     int   d, m;
3358     ivec  shift;
3359
3360
3361     /* Shortest PBC distance between the atom and its reference */
3362     rvec_sub(xcurr, xref, dx);
3363
3364     /* Determine the shift for this atom */
3365     clear_ivec(shift);
3366     for (m = npbcdim-1; m >= 0; m--)
3367     {
3368         while (dx[m] < -0.5*box[m][m])
3369         {
3370             for (d = 0; d < DIM; d++)
3371             {
3372                 dx[d] += box[m][d];
3373             }
3374             shift[m]++;
3375         }
3376         while (dx[m] >= 0.5*box[m][m])
3377         {
3378             for (d = 0; d < DIM; d++)
3379             {
3380                 dx[d] -= box[m][d];
3381             }
3382             shift[m]--;
3383         }
3384     }
3385
3386     /* Apply the shift to the position */
3387     copy_rvec(xcurr, b);
3388     shift_single_coord(box, b, shift);
3389 }
3390
3391
3392 static void init_rot_group(FILE *fplog, t_commrec *cr, int g, t_rotgrp *rotg,
3393                            rvec *x, gmx_mtop_t *mtop, gmx_bool bVerbose, FILE *out_slabs, matrix box,
3394                            t_inputrec *ir, gmx_bool bOutputCenters)
3395 {
3396     int                   i, ii;
3397     rvec                  coord, xref, *xdum;
3398     gmx_bool              bFlex, bColl;
3399     t_atom               *atom;
3400     gmx_enfrotgrp_t       erg; /* Pointer to enforced rotation group data */
3401     int                   ref_firstindex, ref_lastindex;
3402     gmx_mtop_atomlookup_t alook = NULL;
3403     real                  mass, totalmass;
3404     real                  start = 0.0;
3405     double                t_start;
3406
3407
3408     /* Do we have a flexible axis? */
3409     bFlex = ISFLEX(rotg);
3410     /* Do we use a global set of coordinates? */
3411     bColl = ISCOLL(rotg);
3412
3413     erg = rotg->enfrotgrp;
3414
3415     /* Allocate space for collective coordinates if needed */
3416     if (bColl)
3417     {
3418         snew(erg->xc, rotg->nat);
3419         snew(erg->xc_shifts, rotg->nat);
3420         snew(erg->xc_eshifts, rotg->nat);
3421         snew(erg->xc_old, rotg->nat);
3422
3423         if (rotg->eFittype == erotgFitNORM)
3424         {
3425             snew(erg->xc_ref_length, rotg->nat); /* in case fit type NORM is chosen */
3426             snew(erg->xc_norm, rotg->nat);
3427         }
3428     }
3429     else
3430     {
3431         snew(erg->xr_loc, rotg->nat);
3432         snew(erg->x_loc_pbc, rotg->nat);
3433     }
3434
3435     snew(erg->f_rot_loc, rotg->nat);
3436     snew(erg->xc_ref_ind, rotg->nat);
3437
3438     /* Make space for the calculation of the potential at other angles (used
3439      * for fitting only) */
3440     if (erotgFitPOT == rotg->eFittype)
3441     {
3442         snew(erg->PotAngleFit, 1);
3443         snew(erg->PotAngleFit->degangle, rotg->PotAngle_nstep);
3444         snew(erg->PotAngleFit->V, rotg->PotAngle_nstep);
3445         snew(erg->PotAngleFit->rotmat, rotg->PotAngle_nstep);
3446
3447         /* Get the set of angles around the reference angle */
3448         start = -0.5 * (rotg->PotAngle_nstep - 1)*rotg->PotAngle_step;
3449         for (i = 0; i < rotg->PotAngle_nstep; i++)
3450         {
3451             erg->PotAngleFit->degangle[i] = start + i*rotg->PotAngle_step;
3452         }
3453     }
3454     else
3455     {
3456         erg->PotAngleFit = NULL;
3457     }
3458
3459     /* xc_ref_ind needs to be set to identity in the serial case */
3460     if (!PAR(cr))
3461     {
3462         for (i = 0; i < rotg->nat; i++)
3463         {
3464             erg->xc_ref_ind[i] = i;
3465         }
3466     }
3467
3468     /* Copy the masses so that the center can be determined. For all types of
3469      * enforced rotation, we store the masses in the erg->mc array. */
3470     if (rotg->bMassW)
3471     {
3472         alook = gmx_mtop_atomlookup_init(mtop);
3473     }
3474     snew(erg->mc, rotg->nat);
3475     if (bFlex)
3476     {
3477         snew(erg->mc_sorted, rotg->nat);
3478     }
3479     if (!bColl)
3480     {
3481         snew(erg->m_loc, rotg->nat);
3482     }
3483     totalmass = 0.0;
3484     for (i = 0; i < rotg->nat; i++)
3485     {
3486         if (rotg->bMassW)
3487         {
3488             gmx_mtop_atomnr_to_atom(alook, rotg->ind[i], &atom);
3489             mass = atom->m;
3490         }
3491         else
3492         {
3493             mass = 1.0;
3494         }
3495         erg->mc[i] = mass;
3496         totalmass += mass;
3497     }
3498     erg->invmass = 1.0/totalmass;
3499
3500     if (rotg->bMassW)
3501     {
3502         gmx_mtop_atomlookup_destroy(alook);
3503     }
3504
3505     /* Set xc_ref_center for any rotation potential */
3506     if ((rotg->eType == erotgISO) || (rotg->eType == erotgPM) || (rotg->eType == erotgRM) || (rotg->eType == erotgRM2))
3507     {
3508         /* Set the pivot point for the fixed, stationary-axis potentials. This
3509          * won't change during the simulation */
3510         copy_rvec(rotg->pivot, erg->xc_ref_center);
3511         copy_rvec(rotg->pivot, erg->xc_center    );
3512     }
3513     else
3514     {
3515         /* Center of the reference positions */
3516         get_center(rotg->x_ref, erg->mc, rotg->nat, erg->xc_ref_center);
3517
3518         /* Center of the actual positions */
3519         if (MASTER(cr))
3520         {
3521             snew(xdum, rotg->nat);
3522             for (i = 0; i < rotg->nat; i++)
3523             {
3524                 ii = rotg->ind[i];
3525                 copy_rvec(x[ii], xdum[i]);
3526             }
3527             get_center(xdum, erg->mc, rotg->nat, erg->xc_center);
3528             sfree(xdum);
3529         }
3530 #ifdef GMX_MPI
3531         if (PAR(cr))
3532         {
3533             gmx_bcast(sizeof(erg->xc_center), erg->xc_center, cr);
3534         }
3535 #endif
3536     }
3537
3538     if (bColl)
3539     {
3540         /* Save the original (whole) set of positions in xc_old such that at later
3541          * steps the rotation group can always be made whole again. If the simulation is
3542          * restarted, we compute the starting reference positions (given the time)
3543          * and assume that the correct PBC image of each position is the one nearest
3544          * to the current reference */
3545         if (MASTER(cr))
3546         {
3547             /* Calculate the rotation matrix for this angle: */
3548             t_start       = ir->init_t + ir->init_step*ir->delta_t;
3549             erg->degangle = rotg->rate * t_start;
3550             calc_rotmat(rotg->vec, erg->degangle, erg->rotmat);
3551
3552             for (i = 0; i < rotg->nat; i++)
3553             {
3554                 ii = rotg->ind[i];
3555
3556                 /* Subtract pivot, rotate, and add pivot again. This will yield the
3557                  * reference position for time t */
3558                 rvec_sub(rotg->x_ref[i], erg->xc_ref_center, coord);
3559                 mvmul(erg->rotmat, coord, xref);
3560                 rvec_inc(xref, erg->xc_ref_center);
3561
3562                 copy_correct_pbc_image(x[ii], erg->xc_old[i], xref, box, 3);
3563             }
3564         }
3565 #ifdef GMX_MPI
3566         if (PAR(cr))
3567         {
3568             gmx_bcast(rotg->nat*sizeof(erg->xc_old[0]), erg->xc_old, cr);
3569         }
3570 #endif
3571     }
3572
3573     if ( (rotg->eType != erotgFLEX) && (rotg->eType != erotgFLEX2) )
3574     {
3575         /* Put the reference positions into origin: */
3576         for (i = 0; i < rotg->nat; i++)
3577         {
3578             rvec_dec(rotg->x_ref[i], erg->xc_ref_center);
3579         }
3580     }
3581
3582     /* Enforced rotation with flexible axis */
3583     if (bFlex)
3584     {
3585         /* Calculate maximum beta value from minimum gaussian (performance opt.) */
3586         erg->max_beta = calc_beta_max(rotg->min_gaussian, rotg->slab_dist);
3587
3588         /* Determine the smallest and largest coordinate with respect to the rotation vector */
3589         get_firstlast_atom_ref(rotg, &ref_firstindex, &ref_lastindex);
3590
3591         /* From the extreme positions of the reference group, determine the first
3592          * and last slab of the reference. */
3593         get_firstlast_slab_ref(rotg, erg->mc, ref_firstindex, ref_lastindex);
3594
3595         /* Allocate memory for the slabs */
3596         allocate_slabs(rotg, fplog, g, bVerbose);
3597
3598         /* Flexible rotation: determine the reference centers for the rest of the simulation */
3599         erg->slab_first = erg->slab_first_ref;
3600         erg->slab_last  = erg->slab_last_ref;
3601         get_slab_centers(rotg, rotg->x_ref, erg->mc, g, -1, out_slabs, bOutputCenters, TRUE);
3602
3603         /* Length of each x_rotref vector from center (needed if fit routine NORM is chosen): */
3604         if (rotg->eFittype == erotgFitNORM)
3605         {
3606             for (i = 0; i < rotg->nat; i++)
3607             {
3608                 rvec_sub(rotg->x_ref[i], erg->xc_ref_center, coord);
3609                 erg->xc_ref_length[i] = norm(coord);
3610             }
3611         }
3612     }
3613 }
3614
3615
3616 extern void dd_make_local_rotation_groups(gmx_domdec_t *dd, t_rot *rot)
3617 {
3618     gmx_ga2la_t     ga2la;
3619     int             g;
3620     t_rotgrp       *rotg;
3621     gmx_enfrotgrp_t erg;      /* Pointer to enforced rotation group data */
3622
3623     ga2la = dd->ga2la;
3624
3625     for (g = 0; g < rot->ngrp; g++)
3626     {
3627         rotg = &rot->grp[g];
3628         erg  = rotg->enfrotgrp;
3629
3630
3631         dd_make_local_group_indices(ga2la, rotg->nat, rotg->ind,
3632                                     &erg->nat_loc, &erg->ind_loc, &erg->nalloc_loc, erg->xc_ref_ind);
3633     }
3634 }
3635
3636
3637 /* Calculate the size of the MPI buffer needed in reduce_output() */
3638 static int calc_mpi_bufsize(t_rot *rot)
3639 {
3640     int             g;
3641     int             count_group, count_total;
3642     t_rotgrp       *rotg;
3643     gmx_enfrotgrp_t erg;      /* Pointer to enforced rotation group data */
3644
3645
3646     count_total = 0;
3647     for (g = 0; g < rot->ngrp; g++)
3648     {
3649         rotg = &rot->grp[g];
3650         erg  = rotg->enfrotgrp;
3651
3652         /* Count the items that are transferred for this group: */
3653         count_group = 4; /* V, torque, angle, weight */
3654
3655         /* Add the maximum number of slabs for flexible groups */
3656         if (ISFLEX(rotg))
3657         {
3658             count_group += erg->slab_last_ref - erg->slab_first_ref + 1;
3659         }
3660
3661         /* Add space for the potentials at different angles: */
3662         if (erotgFitPOT == rotg->eFittype)
3663         {
3664             count_group += rotg->PotAngle_nstep;
3665         }
3666
3667         /* Add to the total number: */
3668         count_total += count_group;
3669     }
3670
3671     return count_total;
3672 }
3673
3674
3675 extern void init_rot(FILE *fplog, t_inputrec *ir, int nfile, const t_filenm fnm[],
3676                      t_commrec *cr, rvec *x, matrix box, gmx_mtop_t *mtop, const output_env_t oenv,
3677                      gmx_bool bVerbose, unsigned long Flags)
3678 {
3679     t_rot          *rot;
3680     t_rotgrp       *rotg;
3681     int             g;
3682     int             nat_max = 0;  /* Size of biggest rotation group */
3683     gmx_enfrot_t    er;           /* Pointer to the enforced rotation buffer variables */
3684     gmx_enfrotgrp_t erg;          /* Pointer to enforced rotation group data */
3685     rvec           *x_pbc = NULL; /* Space for the pbc-correct atom positions */
3686
3687
3688     if ( (PAR(cr)) && !DOMAINDECOMP(cr) )
3689     {
3690         gmx_fatal(FARGS, "Enforced rotation is only implemented for domain decomposition!");
3691     }
3692
3693     if (MASTER(cr) && bVerbose)
3694     {
3695         fprintf(stdout, "%s Initializing ...\n", RotStr);
3696     }
3697
3698     rot = ir->rot;
3699     snew(rot->enfrot, 1);
3700     er        = rot->enfrot;
3701     er->Flags = Flags;
3702
3703     /* When appending, skip first output to avoid duplicate entries in the data files */
3704     if (er->Flags & MD_APPENDFILES)
3705     {
3706         er->bOut = FALSE;
3707     }
3708     else
3709     {
3710         er->bOut = TRUE;
3711     }
3712
3713     if (MASTER(cr) && er->bOut)
3714     {
3715         please_cite(fplog, "Kutzner2011");
3716     }
3717
3718     /* Output every step for reruns */
3719     if (er->Flags & MD_RERUN)
3720     {
3721         if (NULL != fplog)
3722         {
3723             fprintf(fplog, "%s rerun - will write rotation output every available step.\n", RotStr);
3724         }
3725         rot->nstrout = 1;
3726         rot->nstsout = 1;
3727     }
3728
3729     er->out_slabs = NULL;
3730     if (MASTER(cr) && HaveFlexibleGroups(rot) )
3731     {
3732         er->out_slabs = open_slab_out(opt2fn("-rs", nfile, fnm), rot);
3733     }
3734
3735     if (MASTER(cr))
3736     {
3737         /* Remove pbc, make molecule whole.
3738          * When ir->bContinuation=TRUE this has already been done, but ok. */
3739         snew(x_pbc, mtop->natoms);
3740         m_rveccopy(mtop->natoms, x, x_pbc);
3741         do_pbc_first_mtop(NULL, ir->ePBC, box, mtop, x_pbc);
3742         /* All molecules will be whole now, but not necessarily in the home box.
3743          * Additionally, if a rotation group consists of more than one molecule
3744          * (e.g. two strands of DNA), each one of them can end up in a different
3745          * periodic box. This is taken care of in init_rot_group.  */
3746     }
3747
3748     for (g = 0; g < rot->ngrp; g++)
3749     {
3750         rotg = &rot->grp[g];
3751
3752         if (NULL != fplog)
3753         {
3754             fprintf(fplog, "%s group %d type '%s'\n", RotStr, g, erotg_names[rotg->eType]);
3755         }
3756
3757         if (rotg->nat > 0)
3758         {
3759             /* Allocate space for the rotation group's data: */
3760             snew(rotg->enfrotgrp, 1);
3761             erg  = rotg->enfrotgrp;
3762
3763             nat_max = max(nat_max, rotg->nat);
3764
3765             if (PAR(cr))
3766             {
3767                 erg->nat_loc    = 0;
3768                 erg->nalloc_loc = 0;
3769                 erg->ind_loc    = NULL;
3770             }
3771             else
3772             {
3773                 erg->nat_loc = rotg->nat;
3774                 erg->ind_loc = rotg->ind;
3775             }
3776             init_rot_group(fplog, cr, g, rotg, x_pbc, mtop, bVerbose, er->out_slabs, box, ir,
3777                            !(er->Flags & MD_APPENDFILES) ); /* Do not output the reference centers
3778                                                              * again if we are appending */
3779         }
3780     }
3781
3782     /* Allocate space for enforced rotation buffer variables */
3783     er->bufsize = nat_max;
3784     snew(er->data, nat_max);
3785     snew(er->xbuf, nat_max);
3786     snew(er->mbuf, nat_max);
3787
3788     /* Buffers for MPI reducing torques, angles, weights (for each group), and V */
3789     if (PAR(cr))
3790     {
3791         er->mpi_bufsize = calc_mpi_bufsize(rot) + 100; /* larger to catch errors */
3792         snew(er->mpi_inbuf, er->mpi_bufsize);
3793         snew(er->mpi_outbuf, er->mpi_bufsize);
3794     }
3795     else
3796     {
3797         er->mpi_bufsize = 0;
3798         er->mpi_inbuf   = NULL;
3799         er->mpi_outbuf  = NULL;
3800     }
3801
3802     /* Only do I/O on the MASTER */
3803     er->out_angles  = NULL;
3804     er->out_rot     = NULL;
3805     er->out_torque  = NULL;
3806     if (MASTER(cr))
3807     {
3808         er->out_rot = open_rot_out(opt2fn("-ro", nfile, fnm), rot, oenv);
3809
3810         if (rot->nstsout > 0)
3811         {
3812             if (HaveFlexibleGroups(rot) || HavePotFitGroups(rot) )
3813             {
3814                 er->out_angles  = open_angles_out(opt2fn("-ra", nfile, fnm), rot);
3815             }
3816             if (HaveFlexibleGroups(rot) )
3817             {
3818                 er->out_torque  = open_torque_out(opt2fn("-rt", nfile, fnm), rot);
3819             }
3820         }
3821
3822         sfree(x_pbc);
3823     }
3824 }
3825
3826
3827 extern void finish_rot(t_rot *rot)
3828 {
3829     gmx_enfrot_t er;        /* Pointer to the enforced rotation buffer variables */
3830
3831
3832     er = rot->enfrot;
3833     if (er->out_rot)
3834     {
3835         gmx_fio_fclose(er->out_rot);
3836     }
3837     if (er->out_slabs)
3838     {
3839         gmx_fio_fclose(er->out_slabs);
3840     }
3841     if (er->out_angles)
3842     {
3843         gmx_fio_fclose(er->out_angles);
3844     }
3845     if (er->out_torque)
3846     {
3847         gmx_fio_fclose(er->out_torque);
3848     }
3849 }
3850
3851
3852 /* Rotate the local reference positions and store them in
3853  * erg->xr_loc[0...(nat_loc-1)]
3854  *
3855  * Note that we already subtracted u or y_c from the reference positions
3856  * in init_rot_group().
3857  */
3858 static void rotate_local_reference(t_rotgrp *rotg)
3859 {
3860     gmx_enfrotgrp_t erg;
3861     int             i, ii;
3862
3863
3864     erg = rotg->enfrotgrp;
3865
3866     for (i = 0; i < erg->nat_loc; i++)
3867     {
3868         /* Index of this rotation group atom with respect to the whole rotation group */
3869         ii = erg->xc_ref_ind[i];
3870         /* Rotate */
3871         mvmul(erg->rotmat, rotg->x_ref[ii], erg->xr_loc[i]);
3872     }
3873 }
3874
3875
3876 /* Select the PBC representation for each local x position and store that
3877  * for later usage. We assume the right PBC image of an x is the one nearest to
3878  * its rotated reference */
3879 static void choose_pbc_image(rvec x[], t_rotgrp *rotg, matrix box, int npbcdim)
3880 {
3881     int             i, ii;
3882     gmx_enfrotgrp_t erg;       /* Pointer to enforced rotation group data */
3883     rvec            xref;
3884
3885
3886     erg = rotg->enfrotgrp;
3887
3888     for (i = 0; i < erg->nat_loc; i++)
3889     {
3890         /* Index of a rotation group atom  */
3891         ii = erg->ind_loc[i];
3892
3893         /* Get the correctly rotated reference position. The pivot was already
3894          * subtracted in init_rot_group() from the reference positions. Also,
3895          * the reference positions have already been rotated in
3896          * rotate_local_reference(). For the current reference position we thus
3897          * only need to add the pivot again. */
3898         copy_rvec(erg->xr_loc[i], xref);
3899         rvec_inc(xref, erg->xc_ref_center);
3900
3901         copy_correct_pbc_image(x[ii], erg->x_loc_pbc[i], xref, box, npbcdim);
3902     }
3903 }
3904
3905
3906 extern void do_rotation(
3907         t_commrec      *cr,
3908         t_inputrec     *ir,
3909         matrix          box,
3910         rvec            x[],
3911         real            t,
3912         gmx_int64_t     step,
3913         gmx_wallcycle_t wcycle,
3914         gmx_bool        bNS)
3915 {
3916     int             g, i, ii;
3917     t_rot          *rot;
3918     t_rotgrp       *rotg;
3919     gmx_bool        outstep_slab, outstep_rot;
3920     gmx_bool        bFlex, bColl;
3921     gmx_enfrot_t    er;         /* Pointer to the enforced rotation buffer variables */
3922     gmx_enfrotgrp_t erg;        /* Pointer to enforced rotation group data           */
3923     rvec            transvec;
3924     t_gmx_potfit   *fit = NULL; /* For fit type 'potential' determine the fit
3925                                    angle via the potential minimum            */
3926
3927     /* Enforced rotation cycle counting: */
3928     gmx_cycles_t cycles_comp;   /* Cycles for the enf. rotation computation
3929                                    only, does not count communication. This
3930                                    counter is used for load-balancing         */
3931
3932 #ifdef TAKETIME
3933     double t0;
3934 #endif
3935
3936     rot = ir->rot;
3937     er  = rot->enfrot;
3938
3939     /* When to output in main rotation output file */
3940     outstep_rot  = do_per_step(step, rot->nstrout) && er->bOut;
3941     /* When to output per-slab data */
3942     outstep_slab = do_per_step(step, rot->nstsout) && er->bOut;
3943
3944     /* Output time into rotation output file */
3945     if (outstep_rot && MASTER(cr))
3946     {
3947         fprintf(er->out_rot, "%12.3e", t);
3948     }
3949
3950     /**************************************************************************/
3951     /* First do ALL the communication! */
3952     for (g = 0; g < rot->ngrp; g++)
3953     {
3954         rotg = &rot->grp[g];
3955         erg  = rotg->enfrotgrp;
3956
3957         /* Do we have a flexible axis? */
3958         bFlex = ISFLEX(rotg);
3959         /* Do we use a collective (global) set of coordinates? */
3960         bColl = ISCOLL(rotg);
3961
3962         /* Calculate the rotation matrix for this angle: */
3963         erg->degangle = rotg->rate * t;
3964         calc_rotmat(rotg->vec, erg->degangle, erg->rotmat);
3965
3966         if (bColl)
3967         {
3968             /* Transfer the rotation group's positions such that every node has
3969              * all of them. Every node contributes its local positions x and stores
3970              * it in the collective erg->xc array. */
3971             communicate_group_positions(cr, erg->xc, erg->xc_shifts, erg->xc_eshifts, bNS,
3972                                         x, rotg->nat, erg->nat_loc, erg->ind_loc, erg->xc_ref_ind, erg->xc_old, box);
3973         }
3974         else
3975         {
3976             /* Fill the local masses array;
3977              * this array changes in DD/neighborsearching steps */
3978             if (bNS)
3979             {
3980                 for (i = 0; i < erg->nat_loc; i++)
3981                 {
3982                     /* Index of local atom w.r.t. the collective rotation group */
3983                     ii            = erg->xc_ref_ind[i];
3984                     erg->m_loc[i] = erg->mc[ii];
3985                 }
3986             }
3987
3988             /* Calculate Omega*(y_i-y_c) for the local positions */
3989             rotate_local_reference(rotg);
3990
3991             /* Choose the nearest PBC images of the group atoms with respect
3992              * to the rotated reference positions */
3993             choose_pbc_image(x, rotg, box, 3);
3994
3995             /* Get the center of the rotation group */
3996             if ( (rotg->eType == erotgISOPF) || (rotg->eType == erotgPMPF) )
3997             {
3998                 get_center_comm(cr, erg->x_loc_pbc, erg->m_loc, erg->nat_loc, rotg->nat, erg->xc_center);
3999             }
4000         }
4001
4002     } /* End of loop over rotation groups */
4003
4004     /**************************************************************************/
4005     /* Done communicating, we can start to count cycles for the load balancing now ... */
4006     cycles_comp = gmx_cycles_read();
4007
4008
4009 #ifdef TAKETIME
4010     t0 = MPI_Wtime();
4011 #endif
4012
4013     for (g = 0; g < rot->ngrp; g++)
4014     {
4015         rotg = &rot->grp[g];
4016         erg  = rotg->enfrotgrp;
4017
4018         bFlex = ISFLEX(rotg);
4019         bColl = ISCOLL(rotg);
4020
4021         if (outstep_rot && MASTER(cr))
4022         {
4023             fprintf(er->out_rot, "%12.4f", erg->degangle);
4024         }
4025
4026         /* Calculate angles and rotation matrices for potential fitting: */
4027         if ( (outstep_rot || outstep_slab) && (erotgFitPOT == rotg->eFittype) )
4028         {
4029             fit = erg->PotAngleFit;
4030             for (i = 0; i < rotg->PotAngle_nstep; i++)
4031             {
4032                 calc_rotmat(rotg->vec, erg->degangle + fit->degangle[i], fit->rotmat[i]);
4033
4034                 /* Clear value from last step */
4035                 erg->PotAngleFit->V[i] = 0.0;
4036             }
4037         }
4038
4039         /* Clear values from last time step */
4040         erg->V        = 0.0;
4041         erg->torque_v = 0.0;
4042         erg->angle_v  = 0.0;
4043         erg->weight_v = 0.0;
4044
4045         switch (rotg->eType)
4046         {
4047             case erotgISO:
4048             case erotgISOPF:
4049             case erotgPM:
4050             case erotgPMPF:
4051                 do_fixed(rotg, outstep_rot, outstep_slab);
4052                 break;
4053             case erotgRM:
4054                 do_radial_motion(rotg, outstep_rot, outstep_slab);
4055                 break;
4056             case erotgRMPF:
4057                 do_radial_motion_pf(rotg, x, box, outstep_rot, outstep_slab);
4058                 break;
4059             case erotgRM2:
4060             case erotgRM2PF:
4061                 do_radial_motion2(rotg, x, box, outstep_rot, outstep_slab);
4062                 break;
4063             case erotgFLEXT:
4064             case erotgFLEX2T:
4065                 /* Subtract the center of the rotation group from the collective positions array
4066                  * Also store the center in erg->xc_center since it needs to be subtracted
4067                  * in the low level routines from the local coordinates as well */
4068                 get_center(erg->xc, erg->mc, rotg->nat, erg->xc_center);
4069                 svmul(-1.0, erg->xc_center, transvec);
4070                 translate_x(erg->xc, rotg->nat, transvec);
4071                 do_flexible(MASTER(cr), er, rotg, g, x, box, t, outstep_rot, outstep_slab);
4072                 break;
4073             case erotgFLEX:
4074             case erotgFLEX2:
4075                 /* Do NOT subtract the center of mass in the low level routines! */
4076                 clear_rvec(erg->xc_center);
4077                 do_flexible(MASTER(cr), er, rotg, g, x, box, t, outstep_rot, outstep_slab);
4078                 break;
4079             default:
4080                 gmx_fatal(FARGS, "No such rotation potential.");
4081                 break;
4082         }
4083     }
4084
4085 #ifdef TAKETIME
4086     if (MASTER(cr))
4087     {
4088         fprintf(stderr, "%s calculation (step %d) took %g seconds.\n", RotStr, step, MPI_Wtime()-t0);
4089     }
4090 #endif
4091
4092     /* Stop the enforced rotation cycle counter and add the computation-only
4093      * cycles to the force cycles for load balancing */
4094     cycles_comp  = gmx_cycles_read() - cycles_comp;
4095
4096     if (DOMAINDECOMP(cr) && wcycle)
4097     {
4098         dd_cycles_add(cr->dd, cycles_comp, ddCyclF);
4099     }
4100 }