Merge release-5-0 into release-5-1
[alexxy/gromacs.git] / src / gromacs / swap / swapcoords.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2013,2014,2015, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 /*! \internal \file
36  * \brief
37  * Implements functions in swapcoords.h.
38  *
39  * \author Carsten Kutzner <ckutzne@gwdg.de>
40  * \ingroup module_swap
41  */
42 #include "gmxpre.h"
43
44 #include "swapcoords.h"
45
46 #include <stdio.h>
47 #include <stdlib.h>
48 #include <string.h>
49 #include <time.h>
50
51 #include "gromacs/fileio/confio.h"
52 #include "gromacs/fileio/xvgr.h"
53 #include "gromacs/legacyheaders/copyrite.h"
54 #include "gromacs/legacyheaders/macros.h"
55 #include "gromacs/legacyheaders/mdrun.h"
56 #include "gromacs/legacyheaders/names.h"
57 #include "gromacs/legacyheaders/network.h"
58 #include "gromacs/legacyheaders/typedefs.h"
59 #include "gromacs/math/vec.h"
60 #include "gromacs/mdlib/groupcoord.h"
61 #include "gromacs/pbcutil/pbc.h"
62 #include "gromacs/timing/wallcycle.h"
63 #include "gromacs/topology/mtop_util.h"
64 #include "gromacs/utility/cstringutil.h"
65 #include "gromacs/utility/smalloc.h"
66
67 static const char *SwS      = {"SWAP:"};                                           /**< For output that comes from the swap module */
68 static const char *SwSEmpty = {"     "};                                           /**< Placeholder for multi-line output */
69 static const char* IonString[eIonNR] = {"anion", "cation" };                       /**< Type of ion, used for verbose output */
70 static const char* IonStr[eIonNR]    = {"-", "+"      };                           /**< Type of ion, used for short output */
71 static const char* CompStr[eCompNR] = {"A", "B" };                                 /**< Compartment name */
72 static const char *SwapStr[eSwapTypesNR+1] = { "", "X-", "Y-", "Z-", NULL};        /**< Name for the swap types. */
73 static const char *DimStr[DIM+1] = { "X", "Y", "Z", NULL};                         /**< Name for the swap dimension. */
74
75 /* eGrpSplit0 and eGrpSplit1 _must_ be neighbors in this list because
76  * we sometimes loop from eGrpSplit0 to eGrpSplit1 */
77 enum {
78     eGrpIons, eGrpSplit0, eGrpSplit1, eGrpSolvent, eGrpNr
79 };                                                                               /**< Group identifier */
80 static const char* GrpString[eGrpNr] = { "ion", "split0", "split1", "solvent" }; /**< Group name */
81
82 /** Keep track of through which channel the ions have passed */
83 enum eChannelHistory {
84     eChHistPassedNone, eChHistPassedCh0, eChHistPassedCh1, eChHistNr
85 };
86 static const char* ChannelString[eChHistNr] = { "none", "channel0", "channel1" };  /**< Name for the channels */
87
88 /*! \brief Domain identifier.
89  *
90  * Keeps track of from which compartment the ions came before passing the
91  * channel.
92  */
93 enum eDomain {
94     eDomainNotset, eDomainA, eDomainB, eDomainNr
95 };
96 static const char* DomainString[eDomainNr] = { "not_assigned", "Domain_A", "Domain_B" }; /**< Name for the domains */
97
98
99
100 /*! \internal \brief
101  * Structure containing compartment-specific data.
102  */
103 typedef struct swap_compartment
104 {
105     int                nat;                   /**< Number of atoms matching the
106                                                  compartment conditions.                         */
107     int                nat_old;               /**< Number of atoms before swapping.              */
108     int                nat_req;               /**< Requested number of atoms.                    */
109     real               nat_av;                /**< Time-averaged number of atoms matching
110                                                    the compartment conditions.                   */
111     int               *nat_past;              /**< Past ion counts for time-averaging.           */
112     int               *ind;                   /**< Indices to coll array of atoms.               */
113     real              *dist;                  /**< Distance of atom to compartment center.       */
114     int                nalloc;                /**< Allocation size for ind array.                */
115     int                inflow_netto;          /**< Net inflow of ions into this compartment.     */
116 } t_compartment;
117
118
119 /*! \internal \brief
120  * This structure contains data needed for each of the groups involved in swapping: ions, water,
121  * and channels.
122  */
123 typedef struct swap_group
124 {
125     int               nat;                    /**< Number of atoms in the group                    */
126     int               apm;                    /**< Number of atoms in each molecule                */
127     atom_id          *ind;                    /**< Global atom indices of the group                */
128     atom_id          *ind_loc;                /**< Local atom indices of the group                 */
129     int               nat_loc;                /**< Number of local group atoms                     */
130     int               nalloc_loc;             /**< Allocation size for ind_loc                     */
131     rvec             *xc;                     /**< Collective array of group atom positions        */
132     ivec             *xc_shifts;              /**< Current (collective) shifts                     */
133     ivec             *xc_eshifts;             /**< Extra shifts since last DD step                 */
134     rvec             *xc_old;                 /**< Old (collective) positions                      */
135     real             *qc;                     /**< Collective array of charges                     */
136     int              *c_ind_loc;              /**< Position of local atoms in the
137                                                    collective array, [0..nat_loc]                  */
138     real             *m;                      /**< Masses (can be omitted)                         */
139     unsigned char    *comp_from;              /**< (Collective) Stores from which compartment this
140                                                    atom has come. This way we keep track of through
141                                                    which channel an ion permeates (only used for
142                                                    the ion group)                                  */
143     unsigned char    *comp_now;               /**< In which compartment this ion is now            */
144     unsigned char    *channel_label;          /**< Which channel was passed at last by this ion?   */
145     rvec              center;                 /**< Center of the group; COM if masses are used     */
146 } t_group;
147
148
149 /*! \internal \brief
150  * Main (private) data structure for the position swapping protocol.
151  */
152 typedef struct t_swap
153 {
154     int               swapdim;                       /**< One of XX, YY, ZZ                               */
155     t_pbc            *pbc;                           /**< Needed to make molecules whole.                 */
156     FILE             *fpout;                         /**< Output file.                                    */
157     t_group           group[eGrpNr];                 /**< Ions, solvent or channels?                      */
158     t_compartment     comp[eCompNR][eIonNR];         /**< Data for a specific compartment and ion type.   */
159     t_compartment     compsol[eCompNR];              /**< Solvent compartments.                           */
160     int               fluxfromAtoB[eChanNR][eIonNR]; /**< Net flux per channels and ion type.             */
161     int               ncyl0ions;                     /**< Number of ions residing in channel 0.           */
162     int               ncyl1ions;                     /**< Same for channel 1.                             */
163     int               cyl0and1;                      /**< Ions assigned to cyl0 and cyl1. Not good.       */
164     int              *fluxleak;                      /**< Pointer to a single int value holding the
165                                                           flux not going through any of the channels.     */
166     real              deltaQ;                        /**< The charge imbalance between the compartments.  */
167 } t_swap;
168
169
170
171 /*! \brief Check whether point is in channel.
172  *
173  * A channel is a cylinder defined by a disc
174  * with radius r around its center c. The thickness of the cylinder is
175  * d_up - d_down.
176  *
177  * \code
178  *               ^  d_up
179  *               |
180  *     r         |
181  *     <---------c--------->
182  *               |
183  *               v  d_down
184  *
185  * \endcode
186  */
187 static gmx_bool is_in_channel(
188         rvec   point,  /* Point under consideration */
189         rvec   center, /* 'Center' of cylinder */
190         real   d_up,   /* Upper extension */
191         real   d_down, /* Lower extensions */
192         real   r_cyl2, /* Cylinder radius squared */
193         t_pbc *pbc,
194         int    normal) /* The membrane normal direction is typically 3, i.e. ZZ, but can be X or Y also */
195 {
196     rvec dr;
197     int  plane1, plane2; /* Directions tangential to membrane */
198
199
200     plane1 = (normal + 1) % 3; /* typically 0, i.e. XX */
201     plane2 = (normal + 2) % 3; /* typically 1, i.e. YY */
202
203     /* Get the distance vector dr between the point and the center of the cylinder */
204     pbc_dx(pbc, point, center, dr); /* This puts center in the origin */
205
206     /* Check vertical direction */
207     if ( (dr[normal] > d_up) || (dr[normal] < -d_down) )
208     {
209         return FALSE;
210     }
211
212     /* Check radial direction */
213     if ( (dr[plane1]*dr[plane1] + dr[plane2]*dr[plane2]) > r_cyl2)
214     {
215         return FALSE;
216     }
217
218     /* All check passed, this point is in the cylinder */
219     return TRUE;
220 }
221
222
223 /*! \brief Prints to swap output file which ions are in which compartment. */
224 static void print_ionlist(
225         t_swap       *s,
226         double        time,
227         const char    comment[])
228 {
229     int            itype, icomp, i, j;
230     t_compartment *comp;
231
232
233     fprintf(s->fpout, "%12.5e", time);
234     for (icomp = 0; icomp < eCompNR; icomp++)
235     {
236         for (itype = 0; itype < eIonNR; itype++)
237         {
238             comp = &(s->comp[icomp][itype]);
239             fprintf(s->fpout, "%7d%7.1f%7d", comp->nat, comp->nat_av-comp->nat_req, comp->inflow_netto);
240         }
241     }
242     fprintf(s->fpout, "%12.3e%12.3e",
243             s->group[eGrpSplit0].center[s->swapdim],
244             s->group[eGrpSplit1].center[s->swapdim]);
245
246     for (i = 0; i < eChanNR; i++)
247     {
248         for (j = 0; j < eIonNR; j++)
249         {
250             fprintf(s->fpout, "%12d", s->fluxfromAtoB[i][j]);
251         }
252     }
253
254     /* Also print the number of ions that leaked from A to B: */
255     fprintf(s->fpout, "%12d", *s->fluxleak);
256
257     fprintf(s->fpout, "%s\n", comment);
258 }
259
260
261 /*! \brief Get the center of a group of nat atoms.
262  *
263  * Since with PBC an atom group might not be whole, use the first atom as the
264  * reference atom and determine the center with respect to this reference.
265  */
266 static void get_molecule_center(
267         rvec   x[],
268         int    nat,
269         real  *weights,
270         rvec   center,
271         t_pbc *pbc)
272 {
273     int  i;
274     rvec weightedPBCimage;
275     real wi, wsum;
276     rvec reference, correctPBCimage, dx;
277
278
279     /* Use the first atom as the reference and put other atoms near that one */
280     /* This does not work for large molecules that span > half of the box! */
281     copy_rvec(x[0], reference);
282
283     /* Calculate either the weighted center or simply the center of geometry */
284     wsum = 0.0;
285     clear_rvec(center);
286     for (i = 0; i < nat; i++)
287     {
288         /* PBC distance between position and reference */
289         pbc_dx(pbc, x[i], reference, dx);
290
291         /* Add PBC distance to reference */
292         rvec_add(reference, dx, correctPBCimage);
293
294         /* Take weight into account */
295         if (NULL == weights)
296         {
297             wi = 1.0;
298         }
299         else
300         {
301             wi = weights[i];
302         }
303         wsum += wi;
304         svmul(wi, correctPBCimage, weightedPBCimage);
305
306         /* Add to center */
307         rvec_inc(center, weightedPBCimage);
308     }
309
310     /* Normalize */
311     svmul(1.0/wsum, center, center);
312 }
313
314
315
316 /*! \brief Return TRUE if ion is found in the compartment.
317  *
318  * Returns TRUE if x is between (w1+gap) and (w2-gap)
319  *
320  * \code
321  *
322  * ||-----------|--+--|----------o----------|--+--|---------------------||
323  *                w1   ?????????????????????  w2
324  *
325  * \endcode
326  */
327 static gmx_bool compartment_contains_atom(
328         real  w1, /* position of wall atom 1 */
329         real  w2, /* position of wall atom 2 */
330         real  gap,
331         real  x,
332         real  l,  /* length of the box, from || to || in the sketch */
333         real *distance_from_center)
334 {
335     real m, l_2;
336
337
338     /* First set the origin in the middle of w1 and w2 */
339     m   = 0.5 * (w1 + w2);
340     w1 -= m;
341     w2 -= m;
342     x  -= m;
343
344     /* Now choose the PBC image of x that is closest to the origin: */
345     l_2 = 0.5*l;
346     while (x  > l_2)
347     {
348         x -= l;
349     }
350     while (x <= -l_2)
351     {
352         x += l;
353     }
354
355     *distance_from_center = (real)fabs(x);
356
357     /* Return TRUE if we now are in area "????" */
358     if ( (x >= (w1+gap)) &&  (x < (w2-gap)) )
359     {
360         return TRUE;
361     }
362     else
363     {
364         return FALSE;
365     }
366 }
367
368
369 /*! \brief Updates the time-averaged number of ions in a compartment. */
370 static void update_time_window(t_compartment *comp, int values, int replace)
371 {
372     real average;
373     int  i;
374
375
376     /* Put in the new value */
377     if (replace >= 0)
378     {
379         comp->nat_past[replace] = comp->nat;
380     }
381
382     /* Compute the new time-average */
383     average = 0.0;
384     for (i = 0; i < values; i++)
385     {
386         average += comp->nat_past[i];
387     }
388     average     /= values;
389     comp->nat_av = average;
390 }
391
392
393 /*! \brief Add atom with collective index ci to the list 'comp'. */
394 static void add_to_list(
395         int            ci,       /* index of this ion in the collective array xc, qc */
396         t_compartment *comp,     /* Compartment to add this atom to */
397         real           distance) /* Shortest distance of this atom to the compartment center */
398 {
399     int nr;
400
401
402     nr = comp->nat;
403
404     if (nr >= comp->nalloc)
405     {
406         comp->nalloc = over_alloc_dd(nr+1);
407         srenew(comp->ind, comp->nalloc);
408         srenew(comp->dist, comp->nalloc);
409     }
410     comp->ind[nr]  = ci;
411     comp->dist[nr] = distance;
412     comp->nat++;
413 }
414
415
416 /*! \brief Determine the compartment boundaries from the channel centers. */
417 static void get_compartment_boundaries(
418         int c,
419         t_swap *s,
420         matrix box,
421         real *left, real *right)
422 {
423     real pos0, pos1;
424     real leftpos, rightpos, leftpos_orig;
425
426
427     if (c >= eCompNR)
428     {
429         gmx_fatal(FARGS, "No compartment %d.", c);
430     }
431
432     pos0 = s->group[eGrpSplit0].center[s->swapdim];
433     pos1 = s->group[eGrpSplit1].center[s->swapdim];
434
435     if (pos0 < pos1)
436     {
437         leftpos  = pos0;
438         rightpos = pos1;
439     }
440     else
441     {
442         leftpos  = pos1;
443         rightpos = pos0;
444     }
445
446     /* This gets us the other compartment: */
447     if (c == eCompB)
448     {
449         leftpos_orig = leftpos;
450         leftpos      = rightpos;
451         rightpos     = leftpos_orig + box[s->swapdim][s->swapdim];
452     }
453
454     *left  = leftpos;
455     *right = rightpos;
456 }
457
458
459 /*! \brief Determine the per-channel ion flux.
460  *
461  * To determine the flux through the individual channels, we
462  * remember the compartment and channel history of each ion. An ion can be
463  * either in channel0 or channel1, or in the remaining volume of compartment
464  * A or B.
465  *
466  * \code
467  *    +-----------------+
468  *    |                 | B
469  *    |                 | B compartment
470  *    ||||||||||0|||||||| bilayer with channel 0
471  *    |                 | A
472  *    |                 | A
473  *    |                 | A compartment
474  *    |                 | A
475  *    |||||1||||||||||||| bilayer with channel 1
476  *    |                 | B
477  *    |                 | B compartment
478  *    +-----------------+
479  *
480  * \endcode
481  */
482 static void detect_flux_per_channel(
483         int             iion,
484         int             comp,
485         int             iontype,
486         rvec            ion_pos,
487         unsigned char  *comp_now,
488         unsigned char  *comp_from,
489         unsigned char  *channel_label,
490         t_swapcoords   *sc,
491         real            cyl0_r2,
492         real            cyl1_r2,
493         gmx_int64_t     step,
494         gmx_bool        bRerun,
495         FILE           *fpout)
496 {
497     gmx_swapcoords_t s;
498     int              sd, chan_nr;
499     gmx_bool         in_cyl0, in_cyl1;
500     char             buf[STRLEN];
501
502
503     s    = sc->si_priv;
504     sd   = s->swapdim;
505
506     /* Check whether ion is inside any of the channels */
507     in_cyl0 = is_in_channel(ion_pos, s->group[eGrpSplit0].center, sc->cyl0u, sc->cyl0l, cyl0_r2, s->pbc, sd);
508     in_cyl1 = is_in_channel(ion_pos, s->group[eGrpSplit1].center, sc->cyl1u, sc->cyl1l, cyl1_r2, s->pbc, sd);
509
510     if (in_cyl0 && in_cyl1)
511     {
512         /* Ion appears to be in both channels. Something is severely wrong! */
513         s->cyl0and1++;
514         *comp_now      = eDomainNotset;
515         *comp_from     = eDomainNotset;
516         *channel_label = eChHistPassedNone;
517     }
518     else if (in_cyl0)
519     {
520         /* Ion is in channel 0 now */
521         *channel_label = eChHistPassedCh0;
522         *comp_now      = eDomainNotset;
523         s->ncyl0ions++;
524     }
525     else if (in_cyl1)
526     {
527         /* Ion is in channel 1 now */
528         *channel_label = eChHistPassedCh1;
529         *comp_now      = eDomainNotset;
530         s->ncyl1ions++;
531     }
532     else
533     {
534         /* Ion is not in any of the channels, so it must be in domain A or B */
535         if (eCompA == comp)
536         {
537             *comp_now = eDomainA;
538         }
539         else
540         {
541             *comp_now = eDomainB;
542         }
543     }
544
545     /* Only take action, if ion is now in domain A or B, and was before
546      * in the other domain!
547      */
548     if (eDomainNotset == *comp_from)
549     {
550         /* Maybe we can set the domain now */
551         *comp_from = *comp_now;               /* Could still be eDomainNotset, though */
552     }
553     else if (  (*comp_now  != eDomainNotset ) /* if in channel */
554                && (*comp_from != *comp_now)  )
555     {
556         /* Obviously the ion changed its domain.
557          * Count this for the channel through which it has passed. */
558         switch (*channel_label)
559         {
560             case eChHistPassedNone:
561                 *s->fluxleak = *s->fluxleak + 1;
562
563                 fprintf(stderr, " %s Warning! Step %s, ion %d (%s) moved from %s to %s\n",
564                         SwS, gmx_step_str(step, buf), iion, IonStr[iontype], DomainString[*comp_from], DomainString[*comp_now]);
565                 if (bRerun)
566                 {
567                     fprintf(stderr, ", possibly due to a swap in the original simulation.\n");
568                 }
569                 else
570                 {
571                     fprintf(stderr, "but did not pass cyl0 or cyl1 as defined in the .mdp file.\n"
572                             "Do you have an ion somewhere within the membrane?\n");
573                     /* Write this info to the CompEL output file: */
574                     fprintf(s->fpout, " # Warning: step %s, ion %d (%s) moved from %s to %s (probably through the membrane)\n",
575                             gmx_step_str(step, buf), iion, IonStr[iontype],
576                             DomainString[*comp_from], DomainString[*comp_now]);
577
578                 }
579                 break;
580             case eChHistPassedCh0:
581             case eChHistPassedCh1:
582                 if (*channel_label == eChHistPassedCh0)
583                 {
584                     chan_nr = 0;
585                 }
586                 else
587                 {
588                     chan_nr = 1;
589                 }
590
591                 if (eDomainA == *comp_from)
592                 {
593                     s->fluxfromAtoB[chan_nr][iontype]++;
594                 }
595                 else
596                 {
597                     s->fluxfromAtoB[chan_nr][iontype]--;
598                 }
599                 fprintf(fpout, "# Atom nr. %d finished passing %s.\n", iion, ChannelString[*channel_label]);
600                 break;
601             default:
602                 gmx_fatal(FARGS, "%s Unknown channel history entry!\n", SwS);
603                 break;
604         }
605
606         /* This ion has moved to the _other_ compartment ... */
607         *comp_from = *comp_now;
608         /* ... and it did not pass any channel yet */
609         *channel_label = eChHistPassedNone;
610     }
611 }
612
613
614 /*! \brief Get the lists of ions for the two compartments */
615 static void compartmentalize_ions(
616         t_commrec      *cr,
617         t_swapcoords   *sc,
618         matrix          box,
619         gmx_int64_t     step,
620         FILE           *fpout,
621         gmx_bool        bRerun)
622 {
623     gmx_swapcoords_t s;
624     int              i, sd, replace;
625     real             left, right;
626     t_group         *iong;
627     real             dist;
628     real             cyl0_r2, cyl1_r2;
629     int              comp, type;
630     int              sum, not_in_comp[eCompNR]; /* consistency check */
631     int              ion_nr_global;
632
633
634     s    = sc->si_priv;
635     iong = &s->group[eGrpIons];
636     sd   = s->swapdim;
637
638     cyl0_r2 = sc->cyl0r * sc->cyl0r;
639     cyl1_r2 = sc->cyl1r * sc->cyl1r;
640
641
642     /* Get us a counter that cycles in the range of [0 ... sc->nAverage[ */
643     replace = (step/sc->nstswap) % sc->nAverage;
644
645     for (comp = eCompA; comp <= eCompB; comp++)
646     {
647         /* Get lists of atoms that match criteria for this compartment */
648         get_compartment_boundaries(comp, sc->si_priv, box, &left, &right);
649
650         /* First clear the ion lists */
651         s->comp[comp][eIonNEG].nat = 0;
652         s->comp[comp][eIonPOS].nat = 0;
653         not_in_comp[comp]          = 0; /* consistency check */
654
655         /* Loop over the IONS */
656         for (i = 0; i < iong->nat; i++)
657         {
658             /* Anion or cation? */
659             type = iong->qc[i] < 0 ? eIonNEG : eIonPOS;
660
661             /* Is this ion in the compartment that we look at? */
662             if (compartment_contains_atom(left, right, 0, iong->xc[i][sd], box[sd][sd], &dist) )
663             {
664                 /* Now put it into the list containing only ions of its type */
665                 add_to_list(i, &s->comp[comp][type], dist);
666
667                 /* Master also checks through which channel each ion has passed */
668                 if (MASTER(cr) && (iong->comp_now != NULL))
669                 {
670                     ion_nr_global = iong->ind[i] + 1; /* PDB index starts at 1 ... */
671                     detect_flux_per_channel(ion_nr_global, comp, type, iong->xc[i],
672                                             &iong->comp_now[i], &iong->comp_from[i], &iong->channel_label[i],
673                                             sc, cyl0_r2, cyl1_r2, step, bRerun, fpout);
674                 }
675             }
676             else
677             {
678                 not_in_comp[comp] += 1;
679             }
680         }
681         /* Correct the time-averaged number of ions in both compartments */
682         update_time_window(&s->comp[comp][eIonNEG], sc->nAverage, replace);
683         update_time_window(&s->comp[comp][eIonPOS], sc->nAverage, replace);
684     }
685
686     /* Flux detection warnings */
687     if (MASTER(cr) )
688     {
689         if (s->cyl0and1 > 0)
690         {
691             fprintf(stderr, "\n"
692                     "%s Warning: %d atoms were detected as being in both channels! Probably your split\n"
693                     "%s          cylinder is way too large, or one compartment has collapsed (step %" GMX_PRId64 ")\n",
694                     SwS, s->cyl0and1, SwS, step);
695
696             fprintf(s->fpout, "Warning: %d atoms were assigned to both channels!\n", s->cyl0and1);
697
698             s->cyl0and1 = 0;
699         }
700     }
701
702
703     /* Consistency checks */
704     if (not_in_comp[eCompA] + not_in_comp[eCompB] != iong->nat)
705     {
706         if (NULL != fpout)
707         {
708             fprintf(fpout, "# Warning: Inconsistency during ion compartmentalization. !inA: %d, !inB: %d, total ions %d\n",
709                     not_in_comp[eCompA], not_in_comp[eCompB], iong->nat);
710             fflush(fpout);
711         }
712         else
713         {
714             fprintf(stderr, "%s rank %d: Inconsistency during ion compartmentalization. !inA: %d, !inB: %d, total ions %d\n",
715                     SwS, cr->nodeid, not_in_comp[eCompA], not_in_comp[eCompB], iong->nat);
716
717         }
718     }
719     sum =   s->comp[eCompA][eIonNEG].nat + s->comp[eCompA][eIonPOS].nat
720         + s->comp[eCompB][eIonNEG].nat + s->comp[eCompB][eIonPOS].nat;
721     if (sum != iong->nat)
722     {
723         if (NULL != fpout)
724         {
725             fprintf(fpout, "# Warning: %d atoms are in the ion group, but altogether %d have been assigned to the compartments.\n",
726                     iong->nat, sum);
727             fflush(fpout);
728         }
729         else
730         {
731             fprintf(stderr, "%s rank %d: %d atoms are in the ion group, but altogether %d have been assigned to the compartments.\n",
732                     SwS, cr->nodeid, iong->nat, sum);
733         }
734     }
735
736
737 }
738
739
740 /*! \brief Set up the compartments and get lists of solvent atoms in each compartment */
741 static void compartmentalize_solvent(
742         t_commrec    *cr,
743         t_swapcoords *sc,
744         matrix        box,
745         FILE         *fpout)
746 {
747     gmx_swapcoords_t s;
748     int              apm, i, j, sd;
749     real             left, right;
750     t_group         *solg;
751     real             dist;
752     int              comp;
753     int              sum, not_in_comp[eCompNR]; /* consistency check */
754
755
756     s    = sc->si_priv;
757     solg = &s->group[eGrpSolvent];
758     sd   = s->swapdim;
759     apm  = solg->apm;
760
761     for (comp = eCompA; comp <= eCompB; comp++)
762     {
763         /* Get lists of atoms that match criteria for this compartment */
764         get_compartment_boundaries(comp, sc->si_priv, box, &left, &right);
765
766         /* First clear the solvent molecule lists */
767         s->compsol[comp].nat = 0;
768         not_in_comp[comp]    = 0; /* consistency check */
769
770         /* Loop over the solvent MOLECULES */
771         for (i = 0; i < sc->nat_sol; i += apm)
772         {
773             if (compartment_contains_atom(left, right, 0, solg->xc[i][sd], box[sd][sd], &dist))
774             {
775                 /* Add the whole molecule to the list */
776                 for (j = 0; j < apm; j++)
777                 {
778                     add_to_list(i+j, &s->compsol[comp], dist);
779                 }
780             }
781             else
782             {
783                 not_in_comp[comp] += apm;
784             }
785         }
786     }
787
788     if (NULL != fpout)
789     {
790         fprintf(fpout, "# Solv. molecules in comp.%s: %d   comp.%s: %d\n",
791                 CompStr[eCompA], s->compsol[eCompA].nat/apm,
792                 CompStr[eCompB], s->compsol[eCompB].nat/apm);
793     }
794
795     /* Consistency checks */
796     if (not_in_comp[eCompA] + not_in_comp[eCompB] != solg->nat)
797     {
798         if (NULL != fpout)
799         {
800             fprintf(fpout, "# Warning: Inconsistency during solvent compartmentalization. !inA: %d, !inB: %d, solvent atoms %d\n",
801                     not_in_comp[eCompA], not_in_comp[eCompB], solg->nat);
802             fflush(fpout);
803         }
804         else
805         {
806             fprintf(stderr, "%s rank %d: Inconsistency during solvent compartmentalization. !inA: %d, !inB: %d, solvent atoms %d\n",
807                     SwS, cr->nodeid, not_in_comp[eCompA], not_in_comp[eCompB], solg->nat);
808         }
809     }
810     sum = s->compsol[eCompA].nat + s->compsol[eCompB].nat;
811     if (sum != solg->nat)
812     {
813         if (NULL != fpout)
814         {
815             fprintf(fpout, "# Warning: %d atoms in solvent group, but %d have been assigned to the compartments.\n",
816                     solg->nat, sum);
817             fflush(fpout);
818         }
819         else
820         {
821             fprintf(stderr, "%s rank %d: %d atoms in solvent group, but %d have been assigned to the compartments.\n",
822                     SwS, cr->nodeid, solg->nat, sum);
823         }
824     }
825 }
826
827
828 /*! \brief Find out how many group atoms are in the compartments initially */
829 static void get_initial_ioncounts(
830         t_inputrec       *ir,
831         rvec              x[],   /* the initial positions */
832         matrix            box,
833         t_commrec        *cr,
834         gmx_bool          bRerun)
835 {
836     t_swapcoords *sc;
837     t_swap       *s;
838     int           i, ii, ind, ic;
839     int           req[2], tot[2];
840
841
842     sc = ir->swap;
843     s  = sc->si_priv;
844
845     /* Copy the initial swap group positions to the collective array so
846      * that we can compartmentalize */
847     for (i = 0; i < sc->nat; i++)
848     {
849         ind = sc->ind[i];
850         copy_rvec(x[ind], s->group[eGrpIons].xc[i]);
851     }
852
853     /* Set up the compartments and get lists of atoms in each compartment */
854     compartmentalize_ions(cr, sc, box, 0, s->fpout, bRerun);
855
856     /* Set initial concentrations if requested */
857     for (ic = 0; ic < eCompNR; ic++)
858     {
859         s->comp[ic][eIonPOS].nat_req = sc->ncations[ic];
860         s->comp[ic][eIonNEG].nat_req = sc->nanions[ic];
861     }
862     for (ic = 0; ic < eCompNR; ic++)
863     {
864         for (ii = 0; ii < eIonNR; ii++)
865         {
866             if (s->comp[ic][ii].nat_req < 0)
867             {
868                 s->comp[ic][ii].nat_req = s->comp[ic][ii].nat;
869             }
870         }
871     }
872
873     /* Check whether the number of requested ions adds up to the total number of ions */
874     for (ii = 0; ii < eIonNR; ii++)
875     {
876         req[ii] = s->comp[eCompA][ii].nat_req + s->comp[eCompB][ii].nat_req;
877         tot[ii] = s->comp[eCompA][ii].nat     + s->comp[eCompB][ii].nat;
878     }
879     if ( (req[eCompA] != tot[eCompA]) || (req[eCompB] != tot[eCompB ]) )
880     {
881         gmx_fatal(FARGS, "Mismatch of the number of ions summed over both compartments.\n"
882                   "You requested a total of %d anions and %d cations,\n"
883                   "but there are a total of %d anions and %d cations in the system.\n",
884                   req[eIonNEG], req[eIonPOS],
885                   tot[eIonNEG], tot[eIonPOS]);
886     }
887
888     /* Initialize time-averaging:
889      * Write initial concentrations to all time bins to start with */
890     for (ic = 0; ic < eCompNR; ic++)
891     {
892         for (ii = 0; ii < eIonNR; ii++)
893         {
894             s->comp[ic][ii].nat_av = s->comp[ic][ii].nat;
895             for (i = 0; i < sc->nAverage; i++)
896             {
897                 s->comp[ic][ii].nat_past[i] = s->comp[ic][ii].nat;
898             }
899         }
900     }
901 }
902
903
904 /*! \brief Copy history of ion counts from checkpoint file.
905  *
906  * When called, the checkpoint file has already been read in. Here we copy
907  * over the values from .cpt file to the swap data structure.
908  */
909 static void get_initial_ioncounts_from_cpt(
910         t_inputrec *ir, swapstate_t *swapstate,
911         t_commrec *cr, gmx_bool bVerbose)
912 {
913     t_swapcoords *sc;
914     t_swap       *s;
915     int           ic, ii, j;
916
917     sc = ir->swap;
918     s  = sc->si_priv;
919
920     if (MASTER(cr))
921     {
922         /* Copy the past values from the checkpoint values that have been read in already */
923         if (bVerbose)
924         {
925             fprintf(stderr, "%s Copying values from checkpoint\n", SwS);
926         }
927
928         for (ic = 0; ic < eCompNR; ic++)
929         {
930             for (ii = 0; ii < eIonNR; ii++)
931             {
932                 s->comp[ic][ii].nat_req      = swapstate->nat_req[ic][ii];
933                 s->comp[ic][ii].inflow_netto = swapstate->inflow_netto[ic][ii];
934
935                 if (bVerbose)
936                 {
937                     fprintf(stderr, "%s ... Influx netto: %d   Requested: %d   Past values: ", SwS,
938                             s->comp[ic][ii].inflow_netto, s->comp[ic][ii].nat_req);
939                 }
940
941                 for (j = 0; j < sc->nAverage; j++)
942                 {
943                     s->comp[ic][ii].nat_past[j] = swapstate->nat_past[ic][ii][j];
944                     if (bVerbose)
945                     {
946                         fprintf(stderr, "%d ", s->comp[ic][ii].nat_past[j]);
947                     }
948                 }
949                 if (bVerbose)
950                 {
951                     fprintf(stderr, "\n");
952                 }
953             }
954         }
955     }
956 }
957
958
959 /*! \brief The master lets all others know about the initial ion counts. */
960 static void bc_initial_concentrations(
961         t_commrec    *cr,
962         t_swapcoords *swap)
963 {
964     int     ic, ii;
965     t_swap *s;
966
967     s = swap->si_priv;
968
969     for (ic = 0; ic < eCompNR; ic++)
970     {
971         for (ii = 0; ii < eIonNR; ii++)
972         {
973             gmx_bcast(sizeof(s->comp[ic][ii].nat_req), &(s->comp[ic][ii].nat_req), cr);
974             gmx_bcast(sizeof(s->comp[ic][ii].nat    ), &(s->comp[ic][ii].nat    ), cr);
975             gmx_bcast( swap->nAverage * sizeof(s->comp[ic][ii].nat_past[0]), s->comp[ic][ii].nat_past, cr);
976         }
977     }
978 }
979
980
981 /*! \brief Ensure that each atom belongs to at most one of the swap groups. */
982 static void check_swap_groups(t_swap *s, int nat, gmx_bool bVerbose)
983 {
984     t_group  *g;
985     int       i, j;
986     atom_id  *nGroup    = NULL; /* This array counts for each atom in the MD system to
987                                    how many swap groups it belongs (should be 0 or 1!) */
988     int       ind       = -1;
989     int       nMultiple = 0;    /* Number of atoms belonging to multiple groups */
990
991
992     if (bVerbose)
993     {
994         fprintf(stderr, "%s Making sure each atom belongs to at most one of the swap groups.\n", SwS);
995     }
996
997     /* Add one to the group count of atoms belonging to a swap group: */
998     snew(nGroup, nat);
999     for (i = 0; i < eGrpNr; i++)
1000     {
1001         g = &s->group[i];
1002         for (j = 0; j < g->nat; j++)
1003         {
1004             /* Get the global index of this atom of this group: */
1005             ind = g->ind[j];
1006             nGroup[ind]++;
1007         }
1008     }
1009     /* Make sure each atom belongs to at most one swap group: */
1010     for (j = 0; j < g->nat; j++)
1011     {
1012         if (nGroup[j] > 1)
1013         {
1014             nMultiple++;
1015         }
1016     }
1017     sfree(nGroup);
1018
1019     if (nMultiple)
1020     {
1021         gmx_fatal(FARGS, "%s Cannot perform swapping since %d atom%s allocated to more than one swap index group.\n"
1022                   "%s Each atom must be allocated to at most one of the split groups, the swap group, or the solvent.\n"
1023                   "%s Check the .mdp file settings regarding the swap index groups or the index groups themselves.\n",
1024                   SwS, nMultiple, (1 == nMultiple) ? " is" : "s are", SwSEmpty, SwSEmpty);
1025     }
1026 }
1027
1028
1029 /*! \brief Get the number of atoms per molecule for this group.
1030  *
1031  * Also ensure that all the molecules in this group have this number of atoms.
1032  */
1033 static int get_group_apm_check(
1034         int                         group,
1035         t_swap                     *s,
1036         gmx_bool                    bVerbose,
1037         const gmx_mtop_atomlookup_t alook,
1038         gmx_mtop_t                 *mtop)
1039 {
1040     int *ind;
1041     int  nat, apm, i;
1042     int  molb, molnr, atnr_mol;
1043
1044
1045     ind = s->group[group].ind;
1046     nat = s->group[group].nat;
1047
1048     /* Determine the number of solvent atoms per solvent molecule from the
1049      * first solvent atom: */
1050     i = 0;
1051     gmx_mtop_atomnr_to_molblock_ind(alook, ind[i], &molb, &molnr, &atnr_mol);
1052     apm = mtop->molblock[molb].natoms_mol;
1053
1054     if (bVerbose)
1055     {
1056         fprintf(stderr, "%s Checking whether all %s molecules consist of %d atom%s\n",
1057                 SwS, GrpString[group], apm, apm > 1 ? "s" : "");
1058     }
1059
1060     /* Check whether this is also true for all other solvent atoms */
1061     for (i = 1; i < nat; i++)
1062     {
1063         gmx_mtop_atomnr_to_molblock_ind(alook, ind[i], &molb, &molnr, &atnr_mol);
1064         if (apm != mtop->molblock[molb].natoms_mol)
1065         {
1066             gmx_fatal(FARGS, "Not all %s group molecules consist of %d atoms.",
1067                       GrpString[group], apm);
1068         }
1069     }
1070
1071     return apm;
1072 }
1073
1074
1075 /*! \brief Print the legend to the swap output file.
1076  *
1077  * Also print the initial ion counts
1078  */
1079 static void print_ionlist_legend(t_inputrec *ir, const output_env_t oenv)
1080 {
1081     const char **legend;
1082     int          ic, count, ii;
1083     char         buf[256];
1084     t_swap      *s;
1085
1086
1087     s = ir->swap->si_priv;
1088
1089     snew(legend, eCompNR*eIonNR*3 + 2 + eChanNR*eIonNR + 1);
1090     for (ic = count = 0; ic < eCompNR; ic++)
1091     {
1092         for (ii = 0; ii < eIonNR; ii++)
1093         {
1094             sprintf(buf, "%s %ss", CompStr[ic], IonString[ii]);
1095             legend[count++] = gmx_strdup(buf);
1096             sprintf(buf, "%s av. mismatch to %d%s",
1097                     CompStr[ic], s->comp[ic][ii].nat_req, IonStr[ii]);
1098             legend[count++] = gmx_strdup(buf);
1099             sprintf(buf, "%s netto %s influx", CompStr[ic], IonString[ii]);
1100             legend[count++] = gmx_strdup(buf);
1101         }
1102     }
1103     sprintf(buf, "%scenter of %s of split group 0", SwapStr[ir->eSwapCoords], (NULL != s->group[eGrpSplit0].m) ? "mass" : "geometry");
1104     legend[count++] = gmx_strdup(buf);
1105     sprintf(buf, "%scenter of %s of split group 1", SwapStr[ir->eSwapCoords], (NULL != s->group[eGrpSplit1].m) ? "mass" : "geometry");
1106     legend[count++] = gmx_strdup(buf);
1107
1108     for (ic = 0; ic < eChanNR; ic++)
1109     {
1110         for (ii = 0; ii < eIonNR; ii++)
1111         {
1112             sprintf(buf, "A->ch%d->B %s permeations", ic, IonString[ii]);
1113             legend[count++] = gmx_strdup(buf);
1114         }
1115     }
1116
1117     sprintf(buf, "leakage");
1118     legend[count++] = gmx_strdup(buf);
1119
1120     xvgr_legend(s->fpout, count, legend, oenv);
1121
1122     fprintf(s->fpout, "# Instantaneous ion counts and time-averaged differences to requested numbers\n");
1123     fprintf(s->fpout, "#   time[ps]   A_an   diff   t_in  A_cat   diff   t_in   B_an   diff   t_in  B_cat   diff   t_in ");
1124     fprintf(s->fpout, "   %s-Split0    %s-Split1", DimStr[s->swapdim], DimStr[s->swapdim]);
1125     fprintf(s->fpout, "  A-ch0-B_an A-ch0-B_cat  A-ch1-B_an A-ch1-B_cat ion_leakage\n");
1126     fflush(s->fpout);
1127 }
1128
1129
1130 /*! \brief Initialize channel ion flux detection routine.
1131  *
1132  * Initialize arrays that keep track of where the ions come from and where
1133  * they go.
1134  */
1135 static void detect_flux_per_channel_init(
1136         t_commrec   *cr,
1137         t_swap      *s,
1138         swapstate_t *swapstate,
1139         gmx_bool     bStartFromCpt)
1140 {
1141     int      i, ic, ii;
1142     t_group *g;
1143
1144
1145     g = &(s->group[eGrpIons]);
1146
1147     /* All these flux detection routines run on the master only */
1148     if (!MASTER(cr))
1149     {
1150         g->comp_now      = NULL;
1151         g->comp_from     = NULL;
1152         g->channel_label = NULL;
1153
1154         return;
1155     }
1156
1157     /******************************************************/
1158     /* Channel and domain history for the individual ions */
1159     /******************************************************/
1160     if (bStartFromCpt) /* set the pointers right */
1161     {
1162         g->comp_from     = swapstate->comp_from;
1163         g->channel_label = swapstate->channel_label;
1164     }
1165     else /* allocate memory */
1166     {
1167         snew(g->comp_from, g->nat);
1168         swapstate->comp_from = g->comp_from;
1169         snew(g->channel_label, g->nat);
1170         swapstate->channel_label = g->channel_label;
1171     }
1172     snew(g->comp_now, g->nat);
1173
1174     /* Initialize the channel and domain history counters */
1175     for (i = 0; i < g->nat; i++)
1176     {
1177         g->comp_now[i] = eDomainNotset;
1178         if (!bStartFromCpt)
1179         {
1180             g->comp_from[i]     = eDomainNotset;
1181             g->channel_label[i] = eChHistPassedNone;
1182         }
1183     }
1184
1185     /************************************/
1186     /* Channel fluxes for both channels */
1187     /************************************/
1188     s->ncyl0ions = 0;
1189     s->ncyl1ions = 0;
1190     s->cyl0and1  = 0;
1191
1192     if (bStartFromCpt)
1193     {
1194         fprintf(stderr, "%s Copying channel fluxes from checkpoint file data\n", SwS);
1195     }
1196
1197     for (ic = 0; ic < eChanNR; ic++)
1198     {
1199         fprintf(stderr, "%s Channel %d flux history: ", SwS, ic);
1200         for (ii = 0; ii < eIonNR; ii++)
1201         {
1202             if (bStartFromCpt)
1203             {
1204                 s->fluxfromAtoB[ic][ii] = swapstate->fluxfromAtoB[ic][ii];
1205             }
1206             else
1207             {
1208                 s->fluxfromAtoB[ic][ii] = 0;
1209             }
1210
1211             fprintf(stderr, "%d %s%s   ", s->fluxfromAtoB[ic][ii], IonString[ii], s->fluxfromAtoB[ic][ii] == 1 ? "" : "s");
1212         }
1213         fprintf(stderr, "\n");
1214     }
1215     if (bStartFromCpt)
1216     {
1217         s->fluxleak = swapstate->fluxleak;
1218     }
1219     else
1220     {
1221         snew(s->fluxleak, 1);
1222         *s->fluxleak = 0;
1223         /* Set pointer for checkpoint writing */
1224         swapstate->fluxleak = s->fluxleak;
1225     }
1226
1227     /* Set pointers for checkpoint writing */
1228     for (ic = 0; ic < eChanNR; ic++)
1229     {
1230         for (ii = 0; ii < eIonNR; ii++)
1231         {
1232             swapstate->fluxfromAtoB_p[ic][ii] = &(s->fluxfromAtoB[ic][ii]);
1233         }
1234     }
1235 }
1236
1237
1238 /*! \brief Outputs the initial structure to PDB file for debugging reasons.
1239  *
1240  * Output the starting structure so that in case of multimeric channels
1241  * the user can check whether we have the correct PBC image for all atoms.
1242  * If this is not correct, the ion counts per channel will be very likely
1243  * wrong.
1244  */
1245 static void outputStartStructureIfWanted(gmx_mtop_t *mtop, rvec *x, int ePBC, matrix box)
1246 {
1247     char *env = getenv("GMX_COMPELDUMP");
1248
1249     if (env != NULL)
1250     {
1251         fprintf(stderr, "\n%s Found env.var. GMX_COMPELDUMP, will output CompEL starting structure made whole.\n"
1252                 "%s In case of multimeric channels, please check whether they have the correct PBC representation.\n",
1253                 SwS, SwSEmpty);
1254
1255         write_sto_conf_mtop("CompELAssumedWholeConfiguration.pdb", *mtop->name, mtop, x, NULL, ePBC, box);
1256     }
1257 }
1258
1259
1260 /*! \brief Initialize the swapstate structure, used for checkpoint writing.
1261  *
1262  * The swapstate struct stores the information we need to make the channels
1263  * whole again after restarts from a checkpoint file. Here we do the following:\n
1264  * a) If we did not start from .cpt, we prepare the struct for proper .cpt writing,\n
1265  * b) if we did start from .cpt, we copy over the last whole structures from .cpt,\n
1266  * c) in any case, for subsequent checkpoint writing, we set the pointers in\n
1267  * swapstate to the x_old arrays, which contain the correct PBC representation of
1268  * multimeric channels at the last time step.
1269  */
1270 static void init_swapstate(
1271         swapstate_t      *swapstate,
1272         t_swapcoords     *sc,
1273         gmx_mtop_t       *mtop,
1274         rvec              x[],      /* the initial positions */
1275         matrix            box,
1276         int               ePBC)
1277 {
1278     int                    i, ig;
1279     rvec                  *x_pbc  = NULL;   /* positions of the whole MD system with molecules made whole */
1280     t_group               *g;
1281     t_swap                *s;
1282
1283
1284     s = sc->si_priv;
1285
1286     /* We always need the last whole positions such that
1287      * in the next time step we can make the channels whole again in PBC */
1288     if (swapstate->bFromCpt)
1289     {
1290         /* Copy the last whole positions of each channel from .cpt */
1291         g = &(s->group[eGrpSplit0]);
1292         for (i = 0; i <  g->nat; i++)
1293         {
1294             copy_rvec(swapstate->xc_old_whole[eChan0][i], g->xc_old[i]);
1295         }
1296         g = &(s->group[eGrpSplit1]);
1297         for (i = 0; i <  g->nat; i++)
1298         {
1299             copy_rvec(swapstate->xc_old_whole[eChan1][i], g->xc_old[i]);
1300         }
1301     }
1302     else
1303     {
1304         /* Extract the initial split group positions. */
1305
1306         /* Remove pbc, make molecule whole. */
1307         snew(x_pbc, mtop->natoms);
1308         m_rveccopy(mtop->natoms, x, x_pbc);
1309
1310         /* This can only make individual molecules whole, not multimers */
1311         do_pbc_mtop(NULL, ePBC, box, mtop, x_pbc);
1312
1313         /* Output the starting structure? */
1314         outputStartStructureIfWanted(mtop, x_pbc, ePBC, box);
1315
1316         /* If this is the first run (i.e. no checkpoint present) we assume
1317          * that the starting positions give us the correct PBC representation */
1318         for (ig = eGrpSplit0; ig <= eGrpSplit1; ig++)
1319         {
1320             g = &(s->group[ig]);
1321             for (i = 0; i < g->nat; i++)
1322             {
1323                 copy_rvec(x_pbc[g->ind[i]], g->xc_old[i]);
1324             }
1325         }
1326         sfree(x_pbc);
1327
1328         /* Prepare swapstate arrays for later checkpoint writing */
1329         swapstate->nat[eChan0] = s->group[eGrpSplit0].nat;
1330         swapstate->nat[eChan1] = s->group[eGrpSplit1].nat;
1331     }
1332
1333     /* For subsequent checkpoint writing, set the swapstate pointers to the xc_old
1334      * arrays that get updated at every swapping step */
1335     swapstate->xc_old_whole_p[eChan0] = &s->group[eGrpSplit0].xc_old;
1336     swapstate->xc_old_whole_p[eChan1] = &s->group[eGrpSplit1].xc_old;
1337 }
1338
1339
1340 extern void init_swapcoords(
1341         FILE              *fplog,
1342         gmx_bool           bVerbose,
1343         t_inputrec        *ir,
1344         const char        *fn,
1345         gmx_mtop_t        *mtop,
1346         rvec               x[],
1347         matrix             box,
1348         swapstate_t       *swapstate,
1349         t_commrec         *cr,
1350         const output_env_t oenv,
1351         unsigned long      Flags)
1352 {
1353     int                    i, ic, ig, ii, j;
1354     t_swapcoords          *sc;
1355     t_swap                *s;
1356     t_atom                *atom;
1357     t_group               *g;
1358     gmx_bool               bAppend, bStartFromCpt, bRerun;
1359     gmx_mtop_atomlookup_t  alook = NULL;
1360     matrix                 boxCopy;
1361
1362
1363     alook = gmx_mtop_atomlookup_init(mtop);
1364
1365     if ( (PAR(cr)) && !DOMAINDECOMP(cr) )
1366     {
1367         gmx_fatal(FARGS, "Position swapping is only implemented for domain decomposition!");
1368     }
1369
1370     bAppend       = Flags & MD_APPENDFILES;
1371     bStartFromCpt = Flags & MD_STARTFROMCPT;
1372     bRerun        = Flags & MD_RERUN;
1373
1374     sc = ir->swap;
1375     snew(sc->si_priv, 1);
1376     s = sc->si_priv;
1377
1378     if (bRerun)
1379     {
1380         if (PAR(cr))
1381         {
1382             gmx_fatal(FARGS, "%s This module does not support reruns in parallel\nPlease request a serial run with -nt 1 / -np 1\n", SwS);
1383         }
1384
1385         fprintf(stderr, "%s Rerun - using every available frame\n", SwS);
1386         sc->nstswap  = 1;
1387         sc->nAverage = 1;  /* averaging makes no sense for reruns */
1388     }
1389
1390     if (MASTER(cr) && !bAppend)
1391     {
1392         fprintf(fplog, "\nInitializing ion/water position exchanges\n");
1393         please_cite(fplog, "Kutzner2011b");
1394     }
1395
1396     switch (ir->eSwapCoords)
1397     {
1398         case eswapX:
1399             s->swapdim = XX;
1400             break;
1401         case eswapY:
1402             s->swapdim = YY;
1403             break;
1404         case eswapZ:
1405             s->swapdim = ZZ;
1406             break;
1407         default:
1408             s->swapdim = -1;
1409             break;
1410     }
1411
1412     /* Copy some data to the group structures for convenience */
1413     /* Number of atoms in the group */
1414     s->group[eGrpIons   ].nat = sc->nat;
1415     s->group[eGrpSplit0 ].nat = sc->nat_split[0];
1416     s->group[eGrpSplit1 ].nat = sc->nat_split[1];
1417     s->group[eGrpSolvent].nat = sc->nat_sol;
1418     /* Pointer to the indices */
1419     s->group[eGrpIons   ].ind = sc->ind;
1420     s->group[eGrpSplit0 ].ind = sc->ind_split[0];
1421     s->group[eGrpSplit1 ].ind = sc->ind_split[1];
1422     s->group[eGrpSolvent].ind = sc->ind_sol;
1423
1424     check_swap_groups(s, mtop->natoms, bVerbose && MASTER(cr));
1425
1426     /* Allocate space for the collective arrays for all groups */
1427     for (ig = 0; ig < eGrpNr; ig++)
1428     {
1429         g = &(s->group[ig]);
1430         snew(g->xc, g->nat);
1431         snew(g->c_ind_loc, g->nat);
1432         /* For the split groups (the channels) we need some extra memory to
1433          * be able to make the molecules whole even if they span more than
1434          * half of the box size. */
1435         if (eGrpSplit0 == ig || eGrpSplit1 == ig)
1436         {
1437             snew(g->xc_shifts, g->nat);
1438             snew(g->xc_eshifts, g->nat);
1439             snew(g->xc_old, g->nat);
1440         }
1441     }
1442
1443     if (MASTER(cr))
1444     {
1445         init_swapstate(swapstate, sc, mtop, x, box, ir->ePBC);
1446     }
1447
1448     /* After init_swapstate we have a set of (old) whole positions for our
1449      * channels. Now transfer that to all nodes */
1450     if (PAR(cr))
1451     {
1452         for (ig = eGrpSplit0; ig <= eGrpSplit1; ig++)
1453         {
1454             g = &(s->group[ig]);
1455             gmx_bcast((g->nat)*sizeof((g->xc_old)[0]), g->xc_old, (cr));
1456         }
1457     }
1458
1459     /* Make sure that all molecules in the ion and solvent groups contain the
1460      * same number of atoms each */
1461     s->group[eGrpIons   ].apm = get_group_apm_check(eGrpIons, s, MASTER(cr) && bVerbose, alook, mtop);
1462     s->group[eGrpSolvent].apm = get_group_apm_check(eGrpSolvent, s, MASTER(cr) && bVerbose, alook, mtop);
1463
1464     /* Save masses where needed */
1465     s->group[eGrpIons   ].m = NULL;
1466     /* We only need enough space to determine a single solvent molecule's
1467      * center at at time */
1468     g = &(s->group[eGrpSolvent]);
1469     snew(g->m, g->apm);
1470
1471     /* Need mass-weighted center of split group? */
1472     for (j = 0, ig = eGrpSplit0; j < eChanNR; ig++, j++)
1473     {
1474         g = &(s->group[ig]);
1475         if (TRUE == sc->massw_split[j])
1476         {
1477             /* Save the split group charges if mass-weighting is requested */
1478             snew(g->m, g->nat);
1479             for (i = 0; i < g->nat; i++)
1480             {
1481                 gmx_mtop_atomnr_to_atom(alook, g->ind[i], &atom);
1482                 g->m[i] = atom->m;
1483             }
1484         }
1485         else
1486         {
1487             g->m = NULL;
1488         }
1489     }
1490
1491     /* Save the ionic charges */
1492     g = &(s->group[eGrpIons]);
1493     snew(g->qc, g->nat);
1494     for (i = 0; i < g->nat; i++)
1495     {
1496         gmx_mtop_atomnr_to_atom(alook, g->ind[i], &atom);
1497         g->qc[i] = atom->q;
1498     }
1499
1500     /* Make a t_pbc struct on all nodes so that the molecules
1501      * chosen for an exchange can be made whole. */
1502     snew(s->pbc, 1);
1503     /* Every node needs to call set_pbc() and therefore every node needs
1504      * to know the box dimensions */
1505     copy_mat(box, boxCopy);
1506     if (PAR(cr))
1507     {
1508         gmx_bcast(sizeof(boxCopy), boxCopy, cr);
1509     }
1510     set_pbc(s->pbc, -1, boxCopy);
1511
1512     if (MASTER(cr))
1513     {
1514         if (bVerbose)
1515         {
1516             fprintf(stderr, "%s Opening output file %s%s\n", SwS, fn, bAppend ? " for appending" : "");
1517         }
1518
1519         s->fpout = gmx_fio_fopen(fn, bAppend ? "a" : "w" );
1520
1521         if (!bAppend)
1522         {
1523             xvgr_header(s->fpout, "Ion counts", "Time (ps)", "counts", exvggtXNY, oenv);
1524
1525             for (ig = 0; ig < eGrpNr; ig++)
1526             {
1527                 g = &(s->group[ig]);
1528                 fprintf(s->fpout, "# %s group contains %d atom%s", GrpString[ig], g->nat, (g->nat > 1) ? "s" : "");
1529                 if (eGrpSolvent == ig || eGrpIons == ig)
1530                 {
1531                     fprintf(s->fpout, " with %d atom%s in each molecule", g->apm, (g->apm > 1) ? "s" : "");
1532                 }
1533                 fprintf(s->fpout, ".\n");
1534             }
1535
1536             fprintf(s->fpout, "#\n# Initial positions of split groups:\n");
1537         }
1538
1539         for (j = 0, ig = eGrpSplit0; j < eChanNR; j++, ig++)
1540         {
1541             g = &(s->group[ig]);
1542             for (i = 0; i < g->nat; i++)
1543             {
1544                 copy_rvec(x[sc->ind_split[j][i]], g->xc[i]);
1545             }
1546             if (eGrpSplit0 == ig || eGrpSplit1 == ig)
1547             {
1548                 /* xc has the correct PBC representation for the two channels, so we do
1549                  * not need to correct for that */
1550                 get_center(g->xc, g->m, g->nat, g->center);
1551             }
1552             else
1553             {
1554                 /* For the water molecules, we need to make the molecules whole */
1555                 get_molecule_center(g->xc, g->nat, g->m, g->center, s->pbc);
1556             }
1557             if (!bAppend)
1558             {
1559                 fprintf(s->fpout, "# %s group %s-center %5f nm\n", GrpString[ig],
1560                         DimStr[s->swapdim], g->center[s->swapdim]);
1561             }
1562         }
1563
1564         if (!bAppend)
1565         {
1566             fprintf(s->fpout, "#\n");
1567             fprintf(s->fpout, "# split0 cylinder radius %f nm, up %f nm, down %f nm\n",
1568                     sc->cyl0r, sc->cyl0u, sc->cyl0l);
1569             fprintf(s->fpout, "# split1 cylinder radius %f nm, up %f nm, down %f nm\n",
1570                     sc->cyl1r, sc->cyl1u, sc->cyl1l);
1571         }
1572
1573         if (!bAppend)
1574         {
1575             fprintf(s->fpout, "#\n");
1576             if (!bRerun)
1577             {
1578                 fprintf(s->fpout, "# Coupling constant (number of swap attempt steps to average over): %d  (translates to %f ps).\n",
1579                         sc->nAverage, sc->nAverage*sc->nstswap*ir->delta_t);
1580                 fprintf(s->fpout, "# Threshold is %f\n", sc->threshold);
1581                 fprintf(s->fpout, "#\n");
1582                 fprintf(s->fpout, "# Remarks about which atoms passed which channel use global atoms numbers starting at one.\n");
1583             }
1584         }
1585     }
1586     else
1587     {
1588         s->fpout = NULL;
1589     }
1590
1591     /* Prepare for parallel or serial run */
1592     if (PAR(cr))
1593     {
1594         for (ig = 0; ig < eGrpNr; ig++)
1595         {
1596             g             = &(s->group[ig]);
1597             g->nat_loc    = 0;
1598             g->nalloc_loc = 0;
1599             g->ind_loc    = NULL;
1600         }
1601     }
1602     else
1603     {
1604         for (ig = 0; ig < eGrpNr; ig++)
1605         {
1606             g          = &(s->group[ig]);
1607             g->nat_loc = g->nat;
1608             g->ind_loc = g->ind;
1609             /* c_ind_loc needs to be set to identity in the serial case */
1610             for (i = 0; i < g->nat; i++)
1611             {
1612                 g->c_ind_loc[i] = i;
1613             }
1614         }
1615     }
1616
1617     /* Allocate memory for the ion counts time window */
1618     for (ic = 0; ic < eCompNR; ic++)
1619     {
1620         for (ii = 0; ii < eIonNR; ii++)
1621         {
1622             snew(s->comp[ic][ii].nat_past, sc->nAverage);
1623         }
1624     }
1625
1626     /* Get the initial ion concentrations and let the other nodes know */
1627     if (MASTER(cr))
1628     {
1629         swapstate->nions = s->group[eGrpIons].nat;
1630
1631         if (bStartFromCpt)
1632         {
1633             get_initial_ioncounts_from_cpt(ir, swapstate, cr, bVerbose);
1634         }
1635         else
1636         {
1637             fprintf(stderr, "%s Determining initial ion counts.\n", SwS);
1638             get_initial_ioncounts(ir, x, box, cr, bRerun);
1639         }
1640
1641         /* Prepare (further) checkpoint writes ... */
1642         if (bStartFromCpt)
1643         {
1644             /* Consistency check */
1645             if (swapstate->nAverage != sc->nAverage)
1646             {
1647                 gmx_fatal(FARGS, "%s Ion count averaging steps mismatch! checkpoint: %d, tpr: %d",
1648                           SwS, swapstate->nAverage, sc->nAverage);
1649             }
1650         }
1651         else
1652         {
1653             swapstate->nAverage = sc->nAverage;
1654         }
1655         fprintf(stderr, "%s Setting pointers for checkpoint writing\n", SwS);
1656         for (ic = 0; ic < eCompNR; ic++)
1657         {
1658             for (ii = 0; ii < eIonNR; ii++)
1659             {
1660                 swapstate->nat_req_p[ic][ii]      = &(s->comp[ic][ii].nat_req);
1661                 swapstate->nat_past_p[ic][ii]     = &(s->comp[ic][ii].nat_past[0]);
1662                 swapstate->inflow_netto_p[ic][ii] = &(s->comp[ic][ii].inflow_netto);
1663             }
1664         }
1665
1666         /* Determine the total charge imbalance */
1667         s->deltaQ =  ( (-1) * s->comp[eCompA][eIonNEG].nat_req + s->comp[eCompA][eIonPOS].nat_req )
1668             - ( (-1) * s->comp[eCompB][eIonNEG].nat_req + s->comp[eCompB][eIonPOS].nat_req );
1669
1670         if (bVerbose)
1671         {
1672             fprintf(stderr, "%s Requested charge imbalance is Q(A) - Q(B) = %gz.\n", SwS, s->deltaQ);
1673         }
1674         if (!bAppend)
1675         {
1676             fprintf(s->fpout, "# Requested charge imbalance is Q(A)-Q(B) = %gz.\n", s->deltaQ);
1677         }
1678     }
1679
1680     if (PAR(cr))
1681     {
1682         bc_initial_concentrations(cr, ir->swap);
1683     }
1684
1685     /* Put the time-averaged number of ions for all compartments */
1686     for (ic = 0; ic < eCompNR; ic++)
1687     {
1688         for (ii = 0; ii < eIonNR; ii++)
1689         {
1690             update_time_window(&(s->comp[ic][ii]), sc->nAverage, -1);
1691         }
1692     }
1693
1694     /* Initialize arrays that keep track of through which channel the ions go */
1695     detect_flux_per_channel_init(cr, s, swapstate, bStartFromCpt);
1696
1697     /* We need to print the legend if we open this file for the first time. */
1698     if (MASTER(cr) && !bAppend)
1699     {
1700         print_ionlist_legend(ir, oenv);
1701     }
1702 }
1703
1704
1705 extern void dd_make_local_swap_groups(gmx_domdec_t *dd, t_swapcoords *sc)
1706 {
1707     t_group *g;
1708     int      ig;
1709
1710
1711     /* Make ion group, split groups and solvent group */
1712     for (ig = 0; ig < eGrpNr; ig++)
1713     {
1714         g = &(sc->si_priv->group[ig]);
1715         dd_make_local_group_indices(dd->ga2la, g->nat, g->ind,
1716                                     &(g->nat_loc), &(g->ind_loc), &(g->nalloc_loc), g->c_ind_loc);
1717     }
1718 }
1719
1720
1721 /*! \brief Do we need to swap ions with water molecules at this step?
1722  *
1723  * From the requested and average ion counts we determine whether a swap is needed
1724  * at this time step.
1725  */
1726 static gmx_bool need_swap(t_swapcoords *sc)
1727 {
1728     t_swap *s;
1729     int     ic, ii;
1730
1731
1732     s = sc->si_priv;
1733     for (ic = 0; ic < eCompNR; ic++)
1734     {
1735         for (ii = 0; ii < eIonNR; ii++)
1736         {
1737             if (s->comp[ic][ii].nat_req - s->comp[ic][ii].nat_av >= sc->threshold)
1738             {
1739                 return TRUE;
1740             }
1741         }
1742     }
1743     return FALSE;
1744 }
1745
1746
1747 /*! \brief Return index of atom that we can use for swapping.
1748  *
1749  * Returns the index of an atom that is far off the compartment boundaries.
1750  * Other atoms of the molecule (if any) will directly follow the returned index
1751  */
1752 static int get_index_of_distant_atom(
1753         t_compartment *comp,
1754         int            apm) /* Atoms per molecule - just return the first atom index of a molecule */
1755 {
1756     int  i, ibest = -1;
1757     real d = GMX_REAL_MAX;
1758
1759
1760     /* comp->nat contains the original number of atoms in this compartment
1761      * prior to doing any swaps. Some of these atoms may already have been
1762      * swapped out, but then they are marked with a distance of GMX_REAL_MAX
1763      */
1764     for (i = 0; i < comp->nat_old; i += apm)
1765     {
1766         if (comp->dist[i] < d)
1767         {
1768             ibest = i;
1769             d     = comp->dist[ibest];
1770         }
1771     }
1772
1773     if (ibest < 0)
1774     {
1775         gmx_fatal(FARGS, "Could not get index of swap atom. Compartment atoms %d before swaps, atoms per molecule %d.",
1776                   comp->nat_old, apm);
1777     }
1778
1779     /* Set the distance of this index to infinity such that it won't get selected again in
1780      * this time step
1781      */
1782     comp->dist[ibest] = GMX_REAL_MAX;
1783
1784     return comp->ind[ibest];
1785 }
1786
1787
1788 /*! \brief Swaps centers of mass and makes molecule whole if broken */
1789 static void translate_positions(
1790         rvec  *x,
1791         int    apm,
1792         rvec   old_com,
1793         rvec   new_com,
1794         t_pbc *pbc)
1795 {
1796     int  i;
1797     rvec reference, dx, correctPBCimage;
1798
1799
1800     /* Use the first atom as the reference for PBC */
1801     copy_rvec(x[0], reference);
1802
1803     for (i = 0; i < apm; i++)
1804     {
1805         /* PBC distance between position and reference */
1806         pbc_dx(pbc, x[i], reference, dx);
1807
1808         /* Add PBC distance to reference */
1809         rvec_add(reference, dx, correctPBCimage);
1810
1811         /* Subtract old_com from correct image and add new_com */
1812         rvec_dec(correctPBCimage, old_com);
1813         rvec_inc(correctPBCimage, new_com);
1814
1815         copy_rvec(correctPBCimage, x[i]);
1816     }
1817 }
1818
1819
1820 /*! \brief Write back the the modified local positions from the collective array to the official positions. */
1821 static void apply_modified_positions(
1822         t_group *g,
1823         rvec     x[])
1824 {
1825     int l, ii, cind;
1826
1827
1828     for (l = 0; l < g->nat_loc; l++)
1829     {
1830         /* Get the right local index to write to */
1831         ii = g->ind_loc[l];
1832         /* Where is the local atom in the collective array? */
1833         cind = g->c_ind_loc[l];
1834
1835         /* Copy the possibly modified position */
1836         copy_rvec(g->xc[cind], x[ii]);
1837     }
1838 }
1839
1840
1841 extern gmx_bool do_swapcoords(
1842         t_commrec        *cr,
1843         gmx_int64_t       step,
1844         double            t,
1845         t_inputrec       *ir,
1846         gmx_wallcycle_t   wcycle,
1847         rvec              x[],
1848         matrix            box,
1849         gmx_mtop_t       *mtop,
1850         gmx_bool          bVerbose,
1851         gmx_bool          bRerun)
1852 {
1853     t_swapcoords         *sc;
1854     t_swap               *s;
1855     int                   j, ii, ic, ig, im, gmax, nswaps;
1856     gmx_bool              bSwap = FALSE;
1857     t_group              *g;
1858     real                  vacancy[eCompNR][eIonNR];
1859     int                   isol, iion;
1860     rvec                  solvent_center, ion_center;
1861     t_atom               *atom;
1862     gmx_mtop_atomlookup_t alook = NULL;
1863
1864
1865     wallcycle_start(wcycle, ewcSWAP);
1866
1867     sc  = ir->swap;
1868     s   = sc->si_priv;
1869
1870     /* Assemble all the positions of the swap group (ig = 0), the split groups
1871      * (ig = 1,2), and possibly the solvent group (ig = 3) */
1872     gmax = eGrpNr;
1873
1874     for (ig = 0; ig < gmax; ig++)
1875     {
1876         g = &(s->group[ig]);
1877         if (eGrpSplit0 == ig || eGrpSplit1 == ig)
1878         {
1879             /* The split groups, i.e. the channels. Here we need  the full
1880              * communicate_group_positions(), so that we can make the molecules
1881              * whole even in cases where they span more than half of the box in
1882              * any dimension */
1883             communicate_group_positions(cr, g->xc, g->xc_shifts, g->xc_eshifts, TRUE,
1884                                         x, g->nat, g->nat_loc, g->ind_loc, g->c_ind_loc, g->xc_old, box);
1885
1886             get_center(g->xc, g->m, g->nat, g->center); /* center of split groups == channels */
1887         }
1888         else
1889         {
1890             /* Swap group (ions), and solvent group. These molecules are small
1891              * and we can always make them whole with a simple distance check.
1892              * Therefore we pass NULL as third argument. */
1893             communicate_group_positions(cr, g->xc, NULL, NULL, FALSE,
1894                                         x, g->nat, g->nat_loc, g->ind_loc, g->c_ind_loc, NULL, NULL);
1895         }
1896     }
1897
1898     /* Set up the compartments and get lists of atoms in each compartment,
1899      * determine how many ions each compartment contains */
1900     compartmentalize_ions(cr, sc, box, step, s->fpout, bRerun);
1901
1902     /* Output how many ions are in the compartments */
1903     if (MASTER(cr))
1904     {
1905         print_ionlist(s, t, "");
1906     }
1907
1908     /* If we are doing a rerun, we are finished here, since we cannot perform
1909      * swaps anyway */
1910     if (bRerun)
1911     {
1912         return FALSE;
1913     }
1914
1915     /* Do we have to perform a swap? */
1916     bSwap = need_swap(sc);
1917     if (bSwap)
1918     {
1919         g = &(s->group[eGrpSolvent]);
1920         communicate_group_positions(cr, g->xc, NULL, NULL, FALSE,
1921                                     x, g->nat, g->nat_loc, g->ind_loc, g->c_ind_loc, NULL, NULL);
1922
1923         compartmentalize_solvent(cr, sc, box, s->fpout);
1924
1925         /* Determine where ions are missing and where ions are too many */
1926         for (ic = 0; ic < eCompNR; ic++)
1927         {
1928             for (ii = 0; ii < eIonNR; ii++)
1929             {
1930                 vacancy[ic][ii] = s->comp[ic][ii].nat_req - s->comp[ic][ii].nat_av;
1931             }
1932         }
1933
1934         /* Remember the original number of ions per compartment */
1935         for (ic = 0; ic < eCompNR; ic++)
1936         {
1937             s->compsol[ic].nat_old = s->compsol[ic].nat;
1938             for (ii = 0; ii < eIonNR; ii++)
1939             {
1940                 s->comp[ic][ii].nat_old = s->comp[ic][ii].nat;
1941             }
1942         }
1943
1944         /* Now actually correct the number of ions */
1945         nswaps = 0;
1946         alook  = gmx_mtop_atomlookup_init(mtop);
1947         for (ic = 0; ic < eCompNR; ic++)
1948         {
1949             for (ii = 0; ii < eIonNR; ii++)
1950             {
1951                 while (vacancy[ic][ii] >= sc->threshold)
1952                 {
1953                     /* Swap in an ion */
1954
1955                     /* Get the xc-index of the first atom of a solvent molecule of this compartment */
1956                     isol = get_index_of_distant_atom(&(s->compsol[ic]), s->group[eGrpSolvent].apm );
1957
1958                     /* Get the xc-index of an ion from the other compartment */
1959                     iion = get_index_of_distant_atom(&(s->comp[(ic+1)%eCompNR][ii]), s->group[eGrpIons].apm );
1960
1961                     /* Get the solvent molecule's center of mass */
1962                     for (im = 0; im < s->group[eGrpSolvent].apm; im++)
1963                     {
1964                         gmx_mtop_atomnr_to_atom(alook, s->group[eGrpSolvent].ind[isol+im], &atom);
1965                         s->group[eGrpSolvent].m[im] = atom->m;
1966                     }
1967                     get_molecule_center(&(s->group[eGrpSolvent].xc[isol]), s->group[eGrpSolvent].apm, s->group[eGrpSolvent].m, solvent_center, s->pbc);
1968                     get_molecule_center(&(s->group[eGrpIons   ].xc[iion]), s->group[eGrpIons   ].apm, NULL, ion_center, s->pbc);
1969
1970                     /* subtract com_solvent and add com_ion */
1971                     translate_positions(&(s->group[eGrpSolvent].xc[isol]), s->group[eGrpSolvent].apm, solvent_center, ion_center, s->pbc);
1972                     /* For the ion, subtract com_ion and add com_solvent */
1973                     translate_positions(&(s->group[eGrpIons   ].xc[iion]), s->group[eGrpIons   ].apm, ion_center, solvent_center, s->pbc);
1974
1975                     vacancy[ic              ][ii]--;
1976                     vacancy[(ic+1) % eCompNR][ii]++;
1977
1978                     /* Keep track of the changes */
1979                     s->comp[ic              ][ii].nat++;
1980                     s->comp[(ic+1) % eCompNR][ii].nat--;
1981                     s->comp[ic              ][ii].inflow_netto++;
1982                     s->comp[(ic+1) % eCompNR][ii].inflow_netto--;
1983                     /* Correct the past time window to still get the right averages from now on */
1984                     s->comp[ic              ][ii].nat_av++;
1985                     s->comp[(ic+1) % eCompNR][ii].nat_av--;
1986                     for (j = 0; j < sc->nAverage; j++)
1987                     {
1988                         s->comp[ic              ][ii].nat_past[j]++;
1989                         s->comp[(ic+1) % eCompNR][ii].nat_past[j]--;
1990                     }
1991                     /* Clear ion history */
1992                     if (MASTER(cr))
1993                     {
1994                         s->group[eGrpIons].channel_label[iion] = eChHistPassedNone;
1995                         s->group[eGrpIons].comp_from[iion]     = eDomainNotset;
1996                     }
1997                     /* That was the swap */
1998                     nswaps++;
1999                 }
2000             }
2001         }
2002         gmx_mtop_atomlookup_destroy(alook);
2003
2004         if (bVerbose)
2005         {
2006             fprintf(stderr, "%s Performed %d swap%s in step %" GMX_PRId64 ".\n", SwS, nswaps, nswaps > 1 ? "s" : "", step);
2007         }
2008         if (s->fpout != NULL)
2009         {
2010             print_ionlist(s, t, "  # after swap");
2011         }
2012
2013         /* Write back the the modified local positions from the collective array to the official coordinates */
2014         apply_modified_positions(&(s->group[eGrpIons   ]), x);
2015         apply_modified_positions(&(s->group[eGrpSolvent]), x);
2016     } /* end of if(bSwap) */
2017
2018     wallcycle_stop(wcycle, ewcSWAP);
2019
2020     return bSwap;
2021 }