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