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