Merge branch release-2016 into release-2018
[gromacs.git] / src / gromacs / swap / swapcoords.cpp
blob801a0574487f0db0f3b6959f4891901332109559
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
35 /*! \internal \file
36 * \brief
37 * Implements functions in swapcoords.h.
39 * \author Carsten Kutzner <ckutzne@gwdg.de>
40 * \ingroup module_swap
42 #include "gmxpre.h"
44 #include "swapcoords.h"
46 #include <stdio.h>
47 #include <stdlib.h>
48 #include <time.h>
50 #include <string>
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/commrec.h"
62 #include "gromacs/mdtypes/inputrec.h"
63 #include "gromacs/mdtypes/md_enums.h"
64 #include "gromacs/mdtypes/observableshistory.h"
65 #include "gromacs/mdtypes/state.h"
66 #include "gromacs/mdtypes/swaphistory.h"
67 #include "gromacs/pbcutil/pbc.h"
68 #include "gromacs/timing/wallcycle.h"
69 #include "gromacs/topology/mtop_lookup.h"
70 #include "gromacs/topology/topology.h"
71 #include "gromacs/utility/cstringutil.h"
72 #include "gromacs/utility/fatalerror.h"
73 #include "gromacs/utility/pleasecite.h"
74 #include "gromacs/utility/smalloc.h"
75 #include "gromacs/utility/snprintf.h"
77 static const char *SwS = {"SWAP:"}; /**< For output that comes from the swap module */
78 static const char *SwSEmpty = {" "}; /**< Placeholder for multi-line output */
79 static const char* CompStr[eCompNR] = {"A", "B" }; /**< Compartment name */
80 static const char *SwapStr[eSwapTypesNR+1] = { "", "X-", "Y-", "Z-", nullptr}; /**< Name for the swap types. */
81 static const char *DimStr[DIM+1] = { "X", "Y", "Z", nullptr}; /**< Name for the swap dimension. */
83 /** Keep track of through which channel the ions have passed */
84 enum eChannelHistory {
85 eChHistPassedNone, eChHistPassedCh0, eChHistPassedCh1, eChHistNr
87 static const char* ChannelString[eChHistNr] = { "none", "channel0", "channel1" }; /**< Name for the channels */
89 /*! \brief Domain identifier.
91 * Keeps track of from which compartment the ions came before passing the
92 * channel.
94 enum eDomain {
95 eDomainNotset, eDomainA, eDomainB, eDomainNr
97 static const char* DomainString[eDomainNr] = { "not_assigned", "Domain_A", "Domain_B" }; /**< Name for the domains */
101 /*! \internal \brief
102 * Structure containing compartment-specific data.
104 typedef struct swap_compartment
106 int nMol; /**< Number of ion or water molecules detected
107 in this compartment. */
108 int nMolBefore; /**< Number of molecules before swapping. */
109 int nMolReq; /**< Requested number of molecules in compartment. */
110 real nMolAv; /**< Time-averaged number of molecules matching
111 the compartment conditions. */
112 int *nMolPast; /**< Past molecule counts for time-averaging. */
113 int *ind; /**< Indices to collective array of atoms. */
114 real *dist; /**< Distance of atom to bulk layer, which is
115 normally the center layer of the compartment */
116 int nalloc; /**< Allocation size for ind array. */
117 int inflow_net; /**< Net inflow of ions into this compartment. */
118 } t_compartment;
121 /*! \internal \brief
122 * This structure contains data needed for the groups involved in swapping:
123 * split group 0, split group 1, solvent group, ion groups.
125 typedef struct swap_group
127 char *molname; /**< Name of the group or ion type */
128 int nat; /**< Number of atoms in the group */
129 int apm; /**< Number of atoms in each molecule */
130 int *ind; /**< Global atom indices of the group (size nat) */
131 int *ind_loc; /**< Local atom indices of the group */
132 int nat_loc; /**< Number of local group atoms */
133 int nalloc_loc; /**< Allocation size for ind_loc */
134 rvec *xc; /**< Collective array of group atom positions (size nat) */
135 ivec *xc_shifts; /**< Current (collective) shifts (size nat) */
136 ivec *xc_eshifts; /**< Extra shifts since last DD step (size nat) */
137 rvec *xc_old; /**< Old (collective) positions (size nat) */
138 real q; /**< Total charge of one molecule of this group */
139 int *c_ind_loc; /**< Position of local atoms in the
140 collective array, [0..nat_loc] */
141 real *m; /**< Masses (can be omitted, size apm) */
142 unsigned char *comp_from; /**< (Collective) Stores from which compartment this
143 molecule has come. This way we keep track of
144 through which channel an ion permeates
145 (size nMol = nat/apm) */
146 unsigned char *comp_now; /**< In which compartment this ion is now (size nMol) */
147 unsigned char *channel_label; /**< Which channel was passed at last by this ion?
148 (size nMol) */
149 rvec center; /**< Center of the group; COM if masses are used */
150 t_compartment comp[eCompNR]; /**< Distribution of particles of this group across
151 the two compartments */
152 real vacancy[eCompNR]; /**< How many molecules need to be swapped in? */
153 int fluxfromAtoB[eChanNR]; /**< Net flux of ions per channel */
154 int nCyl[eChanNR]; /**< Number of ions residing in a channel */
155 int nCylBoth; /**< Ions assigned to cyl0 and cyl1. Not good. */
156 } t_swapgrp;
159 /*! \internal \brief
160 * Main (private) data structure for the position swapping protocol.
162 typedef struct t_swap
164 int swapdim; /**< One of XX, YY, ZZ */
165 t_pbc *pbc; /**< Needed to make molecules whole. */
166 FILE *fpout; /**< Output file. */
167 int ngrp; /**< Number of t_swapgrp groups */
168 t_swapgrp *group; /**< Separate groups for channels, solvent, ions */
169 int fluxleak; /**< Flux not going through any of the channels. */
170 real deltaQ; /**< The charge imbalance between the compartments. */
171 } t_swap;
175 /*! \brief Check whether point is in channel.
177 * A channel is a cylinder defined by a disc
178 * with radius r around its center c. The thickness of the cylinder is
179 * d_up - d_down.
181 * \code
182 * ^ d_up
184 * r |
185 * <---------c--------->
187 * v d_down
189 * \endcode
191 * \param[in] point The position (xyz) under consideration.
192 * \param[in] center The center of the cylinder.
193 * \param[in] d_up The upper extension of the cylinder.
194 * \param[in] d_down The lower extension.
195 * \param[in] r_cyl2 Cylinder radius squared.
196 * \param[in] pbc Structure with info about periodic boundary conditions.
197 * \param[in] normal The membrane normal direction is typically 3, i.e. z, but can be x or y also.
199 * \returns Whether the point is inside the defined cylindric channel.
201 static gmx_bool is_in_channel(
202 rvec point,
203 rvec center,
204 real d_up,
205 real d_down,
206 real r_cyl2,
207 t_pbc *pbc,
208 int normal)
210 rvec dr;
211 int plane1, plane2; /* Directions tangential to membrane */
214 plane1 = (normal + 1) % 3; /* typically 0, i.e. XX */
215 plane2 = (normal + 2) % 3; /* typically 1, i.e. YY */
217 /* Get the distance vector dr between the point and the center of the cylinder */
218 pbc_dx(pbc, point, center, dr); /* This puts center in the origin */
220 /* Check vertical direction */
221 if ( (dr[normal] > d_up) || (dr[normal] < -d_down) )
223 return FALSE;
226 /* Check radial direction */
227 if ( (dr[plane1]*dr[plane1] + dr[plane2]*dr[plane2]) > r_cyl2)
229 return FALSE;
232 /* All check passed, this point is in the cylinder */
233 return TRUE;
237 /*! \brief Prints output to CompEL output file.
239 * Prints to swap output file how many ions are in each compartment,
240 * where the centers of the split groups are, and how many ions of each type
241 * passed the channels.
243 static void print_ionlist(
244 t_swap *s,
245 double time,
246 const char comment[])
248 // Output time
249 fprintf(s->fpout, "%12.5e", time);
251 // Output number of molecules and difference to reference counts for each
252 // compartment and ion type
253 for (int iComp = 0; iComp < eCompNR; iComp++)
255 for (int ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
257 t_compartment *comp = &s->group[ig].comp[iComp];
259 fprintf(s->fpout, "%10d%10.1f%10d", comp->nMol, comp->nMolAv - comp->nMolReq, comp->inflow_net);
263 // Output center of split groups
264 fprintf(s->fpout, "%10g%10g",
265 s->group[eGrpSplit0].center[s->swapdim],
266 s->group[eGrpSplit1].center[s->swapdim]);
268 // Output ion flux for each channel and ion type
269 for (int iChan = 0; iChan < eChanNR; iChan++)
271 for (int ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
273 t_swapgrp *g = &s->group[ig];
274 fprintf(s->fpout, "%10d", g->fluxfromAtoB[iChan]);
278 /* Output the number of molecules that leaked from A to B */
279 fprintf(s->fpout, "%10d", s->fluxleak);
281 fprintf(s->fpout, "%s\n", comment);
285 /*! \brief Get the center of a group of nat atoms.
287 * Since with PBC an atom group might not be whole, use the first atom as the
288 * reference atom and determine the center with respect to this reference.
290 static void get_molecule_center(
291 rvec x[],
292 int nat,
293 real *weights,
294 rvec center,
295 t_pbc *pbc)
297 int i;
298 rvec weightedPBCimage;
299 real wi, wsum;
300 rvec reference, correctPBCimage, dx;
303 /* Use the first atom as the reference and put other atoms near that one */
304 /* This does not work for large molecules that span > half of the box! */
305 copy_rvec(x[0], reference);
307 /* Calculate either the weighted center or simply the center of geometry */
308 wsum = 0.0;
309 clear_rvec(center);
310 for (i = 0; i < nat; i++)
312 /* PBC distance between position and reference */
313 pbc_dx(pbc, x[i], reference, dx);
315 /* Add PBC distance to reference */
316 rvec_add(reference, dx, correctPBCimage);
318 /* Take weight into account */
319 if (nullptr == weights)
321 wi = 1.0;
323 else
325 wi = weights[i];
327 wsum += wi;
328 svmul(wi, correctPBCimage, weightedPBCimage);
330 /* Add to center */
331 rvec_inc(center, weightedPBCimage);
334 /* Normalize */
335 svmul(1.0/wsum, center, center);
340 /*! \brief Return TRUE if position x of ion (or water) is found in the compartment,
341 * i.e. between w1 and w2.
343 * One can define and additional offset "b" if one wants to exchange ions/water
344 * to or from a plane not directly in the middle of w1 and w2. The offset can be
345 * in ]-1.0, ..., +1.0 [.
346 * A bulkOffset of 0.0 means 'no offset', so the swap-layer is directly in the
347 * middle between w1 and w2. Offsets -1.0 < b < 0.0 will yield swaps nearer to w1,
348 * whereas offsets 0.0 < 0 < +1.0 will yield swaps nearer to w2.
350 * \code
352 * ||--------------+-------------|-------------+------------------------||
353 * w1 ? ? ? ? ? ? ? ? ? ? ? w2
354 * ||--------------+-------------|----b--------+------------------------||
355 * -1 0.0 +1
357 * \endcode
359 * \param[in] w1 Position of 'wall' atom 1.
360 * \param[in] w2 Position of 'wall' atom 2.
361 * \param[in] x Position of the ion or the water molecule under consideration.
362 * \param[in] l Length of the box, from || to || in the sketch.
363 * \param[in] bulkOffset Where is the bulk layer "b" to be found between w1 and w2?
364 * \param[out] distance_from_b Distance of x to the bulk layer "b".
366 * \returns TRUE if x is between w1 and w2.
368 * Also computes the distance of x to the compartment center (the layer that is
369 * normally situated in the middle of w1 and w2 that would be considered as having
370 * the bulk concentration of ions).
372 static gmx_bool compartment_contains_atom(
373 real w1,
374 real w2,
375 real x,
376 real l,
377 real bulkOffset,
378 real *distance_from_b)
380 real m, l_2;
381 real width;
384 /* First set the origin in the middle of w1 and w2 */
385 m = 0.5 * (w1 + w2);
386 w1 -= m;
387 w2 -= m;
388 x -= m;
389 width = w2 - w1;
391 /* Now choose the PBC image of x that is closest to the origin: */
392 l_2 = 0.5*l;
393 while (x > l_2)
395 x -= l;
397 while (x <= -l_2)
399 x += l;
402 *distance_from_b = (real)fabs(x - bulkOffset*0.5*width);
404 /* Return TRUE if we now are in area "????" */
405 if ( (x >= w1) && (x < w2) )
407 return TRUE;
409 else
411 return FALSE;
416 /*! \brief Updates the time-averaged number of ions in a compartment. */
417 static void update_time_window(t_compartment *comp, int values, int replace)
419 real average;
420 int i;
423 /* Put in the new value */
424 if (replace >= 0)
426 comp->nMolPast[replace] = comp->nMol;
429 /* Compute the new time-average */
430 average = 0.0;
431 for (i = 0; i < values; i++)
433 average += comp->nMolPast[i];
435 average /= values;
436 comp->nMolAv = average;
440 /*! \brief Add the atom with collective index ci to the atom list in compartment 'comp'.
442 * \param[in] ci Index of this ion in the collective xc array.
443 * \param[inout] comp Compartment to add this atom to.
444 * \param[in] distance Shortest distance of this atom to the bulk layer,
445 * from which ion/water pairs are selected for swapping.
447 static void add_to_list(
448 int ci,
449 t_compartment *comp,
450 real distance)
452 int nr = comp->nMol;
454 if (nr >= comp->nalloc)
456 comp->nalloc = over_alloc_dd(nr+1);
457 srenew(comp->ind, comp->nalloc);
458 srenew(comp->dist, comp->nalloc);
460 comp->ind[nr] = ci;
461 comp->dist[nr] = distance;
462 comp->nMol++;
466 /*! \brief Determine the compartment boundaries from the channel centers. */
467 static void get_compartment_boundaries(
468 int c,
469 t_swap *s,
470 const matrix box,
471 real *left, real *right)
473 real pos0, pos1;
474 real leftpos, rightpos, leftpos_orig;
477 if (c >= eCompNR)
479 gmx_fatal(FARGS, "No compartment %c.", c+'A');
482 pos0 = s->group[eGrpSplit0].center[s->swapdim];
483 pos1 = s->group[eGrpSplit1].center[s->swapdim];
485 if (pos0 < pos1)
487 leftpos = pos0;
488 rightpos = pos1;
490 else
492 leftpos = pos1;
493 rightpos = pos0;
496 /* This gets us the other compartment: */
497 if (c == eCompB)
499 leftpos_orig = leftpos;
500 leftpos = rightpos;
501 rightpos = leftpos_orig + box[s->swapdim][s->swapdim];
504 *left = leftpos;
505 *right = rightpos;
509 /*! \brief Determine the per-channel ion flux.
511 * To determine the flux through the individual channels, we
512 * remember the compartment and channel history of each ion. An ion can be
513 * either in channel0 or channel1, or in the remaining volume of compartment
514 * A or B.
516 * \code
517 * +-----------------+
518 * | | B
519 * | | B compartment
520 * ||||||||||0|||||||| bilayer with channel 0
521 * | | A
522 * | | A
523 * | | A compartment
524 * | | A
525 * |||||1||||||||||||| bilayer with channel 1
526 * | | B
527 * | | B compartment
528 * +-----------------+
530 * \endcode
532 static void detect_flux_per_channel(
533 t_swapgrp *g,
534 int iAtom,
535 int comp,
536 rvec atomPosition,
537 unsigned char *comp_now,
538 unsigned char *comp_from,
539 unsigned char *channel_label,
540 t_swapcoords *sc,
541 real cyl0_r2,
542 real cyl1_r2,
543 gmx_int64_t step,
544 gmx_bool bRerun,
545 FILE *fpout)
547 gmx_swapcoords_t s;
548 int sd, chan_nr;
549 gmx_bool in_cyl0, in_cyl1;
550 char buf[STRLEN];
553 s = sc->si_priv;
554 sd = s->swapdim;
556 /* Check whether ion is inside any of the channels */
557 in_cyl0 = is_in_channel(atomPosition, s->group[eGrpSplit0].center, sc->cyl0u, sc->cyl0l, cyl0_r2, s->pbc, sd);
558 in_cyl1 = is_in_channel(atomPosition, s->group[eGrpSplit1].center, sc->cyl1u, sc->cyl1l, cyl1_r2, s->pbc, sd);
560 if (in_cyl0 && in_cyl1)
562 /* Ion appears to be in both channels. Something is severely wrong! */
563 g->nCylBoth++;
564 *comp_now = eDomainNotset;
565 *comp_from = eDomainNotset;
566 *channel_label = eChHistPassedNone;
568 else if (in_cyl0)
570 /* Ion is in channel 0 now */
571 *channel_label = eChHistPassedCh0;
572 *comp_now = eDomainNotset;
573 g->nCyl[eChan0]++;
575 else if (in_cyl1)
577 /* Ion is in channel 1 now */
578 *channel_label = eChHistPassedCh1;
579 *comp_now = eDomainNotset;
580 g->nCyl[eChan1]++;
582 else
584 /* Ion is not in any of the channels, so it must be in domain A or B */
585 if (eCompA == comp)
587 *comp_now = eDomainA;
589 else
591 *comp_now = eDomainB;
595 /* Only take action, if ion is now in domain A or B, and was before
596 * in the other domain!
598 if (eDomainNotset == *comp_from)
600 /* Maybe we can set the domain now */
601 *comp_from = *comp_now; /* Could still be eDomainNotset, though */
603 else if ( (*comp_now != eDomainNotset ) /* if in channel */
604 && (*comp_from != *comp_now) )
606 /* Obviously the ion changed its domain.
607 * Count this for the channel through which it has passed. */
608 switch (*channel_label)
610 case eChHistPassedNone:
611 ++s->fluxleak;
613 fprintf(stderr, " %s Warning! Step %s, ion %d moved from %s to %s\n",
614 SwS, gmx_step_str(step, buf), iAtom, DomainString[*comp_from], DomainString[*comp_now]);
615 if (bRerun)
617 fprintf(stderr, ", possibly due to a swap in the original simulation.\n");
619 else
621 fprintf(stderr, "but did not pass cyl0 or cyl1 as defined in the .mdp file.\n"
622 "Do you have an ion somewhere within the membrane?\n");
623 /* Write this info to the CompEL output file: */
624 fprintf(s->fpout, " # Warning: step %s, ion %d moved from %s to %s (probably through the membrane)\n",
625 gmx_step_str(step, buf), iAtom,
626 DomainString[*comp_from], DomainString[*comp_now]);
629 break;
630 case eChHistPassedCh0:
631 case eChHistPassedCh1:
632 if (*channel_label == eChHistPassedCh0)
634 chan_nr = 0;
636 else
638 chan_nr = 1;
641 if (eDomainA == *comp_from)
643 g->fluxfromAtoB[chan_nr]++;
645 else
647 g->fluxfromAtoB[chan_nr]--;
649 fprintf(fpout, "# Atom nr. %d finished passing %s.\n", iAtom, ChannelString[*channel_label]);
650 break;
651 default:
652 gmx_fatal(FARGS, "%s Unknown channel history entry for ion type '%s'\n",
653 SwS, g->molname);
654 break;
657 /* This ion has moved to the _other_ compartment ... */
658 *comp_from = *comp_now;
659 /* ... and it did not pass any channel yet */
660 *channel_label = eChHistPassedNone;
665 /*! \brief Determines which ions or solvent molecules are in compartment A and B */
666 static void sortMoleculesIntoCompartments(
667 t_swapgrp *g,
668 t_commrec *cr,
669 t_swapcoords *sc,
670 const matrix box,
671 gmx_int64_t step,
672 FILE *fpout,
673 gmx_bool bRerun,
674 gmx_bool bIsSolvent)
676 gmx_swapcoords_t s = sc->si_priv;
677 int nMolNotInComp[eCompNR]; /* consistency check */
678 real cyl0_r2 = sc->cyl0r * sc->cyl0r;
679 real cyl1_r2 = sc->cyl1r * sc->cyl1r;
681 /* Get us a counter that cycles in the range of [0 ... sc->nAverage[ */
682 int replace = (step/sc->nstswap) % sc->nAverage;
684 for (int comp = eCompA; comp <= eCompB; comp++)
686 real left, right;
688 /* Get lists of atoms that match criteria for this compartment */
689 get_compartment_boundaries(comp, sc->si_priv, box, &left, &right);
691 /* First clear the ion molecule lists */
692 g->comp[comp].nMol = 0;
693 nMolNotInComp[comp] = 0; /* consistency check */
695 /* Loop over the molecules and atoms of this group */
696 for (int iMol = 0, iAtom = 0; iAtom < g->nat; iAtom += g->apm, iMol++)
698 real dist;
699 int sd = s->swapdim;
701 /* Is this first atom of the molecule in the compartment that we look at? */
702 if (compartment_contains_atom(left, right, g->xc[iAtom][sd], box[sd][sd], sc->bulkOffset[comp], &dist) )
704 /* Add the first atom of this molecule to the list of molecules in this compartment */
705 add_to_list(iAtom, &g->comp[comp], dist);
707 /* Master also checks for ion groups through which channel each ion has passed */
708 if (MASTER(cr) && (g->comp_now != nullptr) && !bIsSolvent)
710 int globalAtomNr = g->ind[iAtom] + 1; /* PDB index starts at 1 ... */
711 detect_flux_per_channel(g, globalAtomNr, comp, g->xc[iAtom],
712 &g->comp_now[iMol], &g->comp_from[iMol], &g->channel_label[iMol],
713 sc, cyl0_r2, cyl1_r2, step, bRerun, fpout);
716 else
718 nMolNotInComp[comp]++;
721 /* Correct the time-averaged number of ions in the compartment */
722 if (!bIsSolvent)
724 update_time_window(&g->comp[comp], sc->nAverage, replace);
728 /* Flux detection warnings */
729 if (MASTER(cr) && !bIsSolvent)
731 if (g->nCylBoth > 0)
733 fprintf(stderr, "\n"
734 "%s Warning: %d atoms were detected as being in both channels! Probably your split\n"
735 "%s cylinder is way too large, or one compartment has collapsed (step %" GMX_PRId64 ")\n",
736 SwS, g->nCylBoth, SwS, step);
738 fprintf(s->fpout, "Warning: %d atoms were assigned to both channels!\n", g->nCylBoth);
740 g->nCylBoth = 0;
744 if (bIsSolvent && nullptr != fpout)
746 fprintf(fpout, "# Solv. molecules in comp.%s: %d comp.%s: %d\n",
747 CompStr[eCompA], g->comp[eCompA].nMol,
748 CompStr[eCompB], g->comp[eCompB].nMol);
751 /* Consistency checks */
752 if (nMolNotInComp[eCompA] + nMolNotInComp[eCompB] != g->nat / g->apm)
754 fprintf(stderr, "%s Warning: Inconsistency while assigning '%s' molecules to compartments. !inA: %d, !inB: %d, total molecules %d\n",
755 SwS, g->molname, nMolNotInComp[eCompA], nMolNotInComp[eCompB], g->nat/g->apm);
758 int sum = g->comp[eCompA].nMol + g->comp[eCompB].nMol;
759 if (sum != g->nat/g->apm)
761 fprintf(stderr, "%s Warning: %d molecules are in group '%s', but altogether %d have been assigned to the compartments.\n",
762 SwS, g->nat/g->apm, g->molname, sum);
767 /*! \brief Find out how many group atoms are in the compartments initially */
768 static void get_initial_ioncounts(
769 t_inputrec *ir,
770 const rvec x[], /* the initial positions */
771 const matrix box,
772 t_commrec *cr,
773 gmx_bool bRerun)
775 t_swapcoords *sc;
776 t_swap *s;
777 t_swapgrp *g;
778 int i, ind, ic, ig;
779 int req, tot;
782 sc = ir->swap;
783 s = sc->si_priv;
786 /* Loop over the user-defined (ion) groups */
787 for (ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
789 g = &s->group[ig];
791 /* Copy the initial positions of the atoms in the group
792 * to the collective array so that we can compartmentalize */
793 for (i = 0; i < g->nat; i++)
795 ind = g->ind[i];
796 copy_rvec(x[ind], g->xc[i]);
799 /* Set up the compartments and get lists of atoms in each compartment */
800 sortMoleculesIntoCompartments(g, cr, sc, box, 0, s->fpout, bRerun, FALSE);
802 /* Set initial molecule counts if requested (as signaled by "-1" value) */
803 for (ic = 0; ic < eCompNR; ic++)
805 int requested = sc->grp[ig].nmolReq[ic];
806 if (requested < 0)
808 g->comp[ic].nMolReq = g->comp[ic].nMol;
810 else
812 g->comp[ic].nMolReq = requested;
816 /* Check whether the number of requested molecules adds up to the total number */
817 req = g->comp[eCompA].nMolReq + g->comp[eCompB].nMolReq;
818 tot = g->comp[eCompA].nMol + g->comp[eCompB].nMol;
820 if ( (req != tot) )
822 gmx_fatal(FARGS, "Mismatch of the number of %s ions summed over both compartments.\n"
823 "You requested a total of %d ions (%d in A and %d in B),\n"
824 "but there are a total of %d ions of this type in the system.\n",
825 g->molname, req, g->comp[eCompA].nMolReq,
826 g->comp[eCompB].nMolReq, tot);
829 /* Initialize time-averaging:
830 * Write initial concentrations to all time bins to start with */
831 for (ic = 0; ic < eCompNR; ic++)
833 g->comp[ic].nMolAv = g->comp[ic].nMol;
834 for (i = 0; i < sc->nAverage; i++)
836 g->comp[ic].nMolPast[i] = g->comp[ic].nMol;
843 /*! \brief Copy history of ion counts from checkpoint file.
845 * When called, the checkpoint file has already been read in. Here we copy
846 * over the values from .cpt file to the swap data structure.
848 static void get_initial_ioncounts_from_cpt(
849 t_inputrec *ir, swaphistory_t *swapstate,
850 t_commrec *cr, gmx_bool bVerbose)
852 t_swapcoords *sc;
853 t_swap *s;
854 t_swapgrp *g;
855 swapstateIons_t *gs;
857 sc = ir->swap;
858 s = sc->si_priv;
860 if (MASTER(cr))
862 /* Copy the past values from the checkpoint values that have been read in already */
863 if (bVerbose)
865 fprintf(stderr, "%s Copying values from checkpoint\n", SwS);
868 for (int ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
870 g = &s->group[ig];
871 gs = &swapstate->ionType[ig - eSwapFixedGrpNR];
873 for (int ic = 0; ic < eCompNR; ic++)
875 g->comp[ic].nMolReq = gs->nMolReq[ic];
876 g->comp[ic].inflow_net = gs->inflow_net[ic];
878 if (bVerbose)
880 fprintf(stderr, "%s ... Influx netto: %d Requested: %d Past values: ", SwS,
881 g->comp[ic].inflow_net, g->comp[ic].nMolReq);
884 for (int j = 0; j < sc->nAverage; j++)
886 g->comp[ic].nMolPast[j] = gs->nMolPast[ic][j];
887 if (bVerbose)
889 fprintf(stderr, "%d ", g->comp[ic].nMolPast[j]);
892 if (bVerbose)
894 fprintf(stderr, "\n");
902 /*! \brief The master lets all others know about the initial ion counts. */
903 static void bc_initial_concentrations(
904 t_commrec *cr,
905 t_swapcoords *swap)
907 int ic, ig;
908 t_swap *s;
909 t_swapgrp *g;
912 s = swap->si_priv;
914 for (ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
916 g = &s->group[ig];
918 for (ic = 0; ic < eCompNR; ic++)
920 gmx_bcast(sizeof(g->comp[ic].nMolReq), &(g->comp[ic].nMolReq), cr);
921 gmx_bcast(sizeof(g->comp[ic].nMol ), &(g->comp[ic].nMol ), cr);
922 gmx_bcast( swap->nAverage * sizeof(g->comp[ic].nMolPast[0]), g->comp[ic].nMolPast, cr);
928 /*! \brief Ensure that each atom belongs to at most one of the swap groups. */
929 static void check_swap_groups(t_swap *s, int nat, gmx_bool bVerbose)
931 int *nGroup = nullptr; /* This array counts for each atom in the MD system to
932 how many swap groups it belongs (should be 0 or 1!) */
933 int ind = -1;
934 int nMultiple = 0; /* Number of atoms belonging to multiple groups */
937 if (bVerbose)
939 fprintf(stderr, "%s Making sure each atom belongs to at most one of the swap groups.\n", SwS);
942 /* Add one to the group count of atoms belonging to a swap group: */
943 snew(nGroup, nat);
944 for (int i = 0; i < s->ngrp; i++)
946 t_swapgrp *g = &s->group[i];
947 for (int j = 0; j < g->nat; j++)
949 /* Get the global index of this atom of this group: */
950 ind = g->ind[j];
951 nGroup[ind]++;
954 /* Make sure each atom belongs to at most one of the groups: */
955 for (int i = 0; i < nat; i++)
957 if (nGroup[i] > 1)
959 nMultiple++;
962 sfree(nGroup);
964 if (nMultiple)
966 gmx_fatal(FARGS, "%s Cannot perform swapping since %d atom%s allocated to more than one swap index group.\n"
967 "%s Each atom must be allocated to at most one of the split groups, the swap groups, or the solvent.\n"
968 "%s Check the .mdp file settings regarding the swap index groups or the index groups themselves.\n",
969 SwS, nMultiple, (1 == nMultiple) ? " is" : "s are", SwSEmpty, SwSEmpty);
974 /*! \brief Get the number of atoms per molecule for this group.
976 * Also ensure that all the molecules in this group have this number of atoms.
978 static int get_group_apm_check(
979 int igroup,
980 t_swap *s,
981 gmx_bool bVerbose,
982 gmx_mtop_t *mtop)
984 t_swapgrp *g = &s->group[igroup];
985 int *ind = s->group[igroup].ind;
986 int nat = s->group[igroup].nat;
988 /* Determine the number of solvent atoms per solvent molecule from the
989 * first solvent atom: */
990 int molb = 0;
991 mtopGetMolblockIndex(mtop, ind[0], &molb, nullptr, nullptr);
992 int apm = mtop->molblock[molb].natoms_mol;
994 if (bVerbose)
996 fprintf(stderr, "%s Checking whether all %s molecules consist of %d atom%s\n", SwS,
997 g->molname, apm, apm > 1 ? "s" : "");
1000 /* Check whether this is also true for all other solvent atoms */
1001 for (int i = 1; i < nat; i++)
1003 mtopGetMolblockIndex(mtop, ind[i], &molb, nullptr, nullptr);
1004 if (apm != mtop->molblock[molb].natoms_mol)
1006 gmx_fatal(FARGS, "Not all molecules of swap group %d consist of %d atoms.",
1007 igroup, apm);
1011 //TODO: check whether charges and masses of each molecule are identical!
1012 return apm;
1016 /*! \brief Print the legend to the swap output file.
1018 * Also print the initial values of ion counts and position of split groups.
1020 static void print_ionlist_legend(t_inputrec *ir,
1021 const gmx_output_env_t *oenv)
1023 const char **legend;
1024 int count = 0;
1025 char buf[STRLEN];
1027 t_swap *s = ir->swap->si_priv;
1028 int nIonTypes = ir->swap->ngrp - eSwapFixedGrpNR;
1029 snew(legend, eCompNR*nIonTypes*3 + 2 + eChanNR*nIonTypes + 1);
1031 // Number of molecules and difference to reference counts for each
1032 // compartment and ion type
1033 for (int ic = count = 0; ic < eCompNR; ic++)
1035 for (int ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
1037 t_swapGroup *g = &ir->swap->grp[ig];
1038 real q = s->group[ig].q;
1040 snprintf(buf, STRLEN, "%s %s ions (charge %s%g)", CompStr[ic], g->molname, q > 0 ? "+" : "", q);
1041 legend[count++] = gmx_strdup(buf);
1043 snprintf(buf, STRLEN, "%s av. mismatch to %d %s ions",
1044 CompStr[ic], s->group[ig].comp[ic].nMolReq, g->molname);
1045 legend[count++] = gmx_strdup(buf);
1047 snprintf(buf, STRLEN, "%s net %s ion influx", CompStr[ic], g->molname);
1048 legend[count++] = gmx_strdup(buf);
1052 // Center of split groups
1053 snprintf(buf, STRLEN, "%scenter of %s of split group 0", SwapStr[ir->eSwapCoords], (nullptr != s->group[eGrpSplit0].m) ? "mass" : "geometry");
1054 legend[count++] = gmx_strdup(buf);
1055 snprintf(buf, STRLEN, "%scenter of %s of split group 1", SwapStr[ir->eSwapCoords], (nullptr != s->group[eGrpSplit1].m) ? "mass" : "geometry");
1056 legend[count++] = gmx_strdup(buf);
1058 // Ion flux for each channel and ion type
1059 for (int ic = 0; ic < eChanNR; ic++)
1061 for (int ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
1063 t_swapGroup *g = &ir->swap->grp[ig];
1064 snprintf(buf, STRLEN, "A->ch%d->B %s permeations", ic, g->molname);
1065 legend[count++] = gmx_strdup(buf);
1069 // Number of molecules that leaked from A to B
1070 snprintf(buf, STRLEN, "leakage");
1071 legend[count++] = gmx_strdup(buf);
1073 xvgr_legend(s->fpout, count, legend, oenv);
1075 fprintf(s->fpout, "# Instantaneous ion counts and time-averaged differences to requested numbers\n");
1077 // We add a simple text legend helping to identify the columns with xvgr legend strings
1078 fprintf(s->fpout, "# time (ps)");
1079 for (int i = 0; i < count; i++)
1081 snprintf(buf, STRLEN, "s%d", i);
1082 fprintf(s->fpout, "%10s", buf);
1084 fprintf(s->fpout, "\n");
1085 fflush(s->fpout);
1089 /*! \brief Initialize channel ion flux detection routine.
1091 * Initialize arrays that keep track of where the ions come from and where
1092 * they go.
1094 static void detect_flux_per_channel_init(
1095 t_swap *s,
1096 swaphistory_t *swapstate,
1097 gmx_bool bStartFromCpt)
1099 t_swapgrp *g;
1100 swapstateIons_t *gs;
1102 /* All these flux detection routines run on the master only */
1103 if (swapstate == nullptr)
1105 return;
1108 for (int ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
1110 g = &s->group[ig];
1111 gs = &swapstate->ionType[ig - eSwapFixedGrpNR];
1113 /******************************************************/
1114 /* Channel and domain history for the individual ions */
1115 /******************************************************/
1116 if (bStartFromCpt) /* set the pointers right */
1118 g->comp_from = gs->comp_from;
1119 g->channel_label = gs->channel_label;
1121 else /* allocate memory for molecule counts */
1123 snew(g->comp_from, g->nat/g->apm);
1124 gs->comp_from = g->comp_from;
1125 snew(g->channel_label, g->nat/g->apm);
1126 gs->channel_label = g->channel_label;
1128 snew(g->comp_now, g->nat/g->apm);
1130 /* Initialize the channel and domain history counters */
1131 for (int i = 0; i < g->nat/g->apm; i++)
1133 g->comp_now[i] = eDomainNotset;
1134 if (!bStartFromCpt)
1136 g->comp_from[i] = eDomainNotset;
1137 g->channel_label[i] = eChHistPassedNone;
1141 /************************************/
1142 /* Channel fluxes for both channels */
1143 /************************************/
1144 g->nCyl[eChan0] = 0;
1145 g->nCyl[eChan1] = 0;
1146 g->nCylBoth = 0;
1149 if (bStartFromCpt)
1151 fprintf(stderr, "%s Copying channel fluxes from checkpoint file data\n", SwS);
1155 // Loop over ion types (and both channels)
1156 for (int ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
1158 g = &s->group[ig];
1159 gs = &swapstate->ionType[ig - eSwapFixedGrpNR];
1161 for (int ic = 0; ic < eChanNR; ic++)
1163 fprintf(stderr, "%s Channel %d flux history for ion type %s (charge %g): ", SwS, ic, g->molname, g->q);
1164 if (bStartFromCpt)
1166 g->fluxfromAtoB[ic] = gs->fluxfromAtoB[ic];
1168 else
1170 g->fluxfromAtoB[ic] = 0;
1173 fprintf(stderr, "%d molecule%s",
1174 g->fluxfromAtoB[ic], g->fluxfromAtoB[ic] == 1 ? "" : "s");
1175 fprintf(stderr, "\n");
1179 /* Set pointers for checkpoint writing */
1180 swapstate->fluxleak_p = &s->fluxleak;
1181 for (int ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
1183 g = &s->group[ig];
1184 gs = &swapstate->ionType[ig - eSwapFixedGrpNR];
1186 for (int ic = 0; ic < eChanNR; ic++)
1188 gs->fluxfromAtoB_p[ic] = &g->fluxfromAtoB[ic];
1194 /*! \brief Outputs the initial structure to PDB file for debugging reasons.
1196 * Output the starting structure so that in case of multimeric channels
1197 * the user can check whether we have the correct PBC image for all atoms.
1198 * If this is not correct, the ion counts per channel will be very likely
1199 * wrong.
1201 static void outputStartStructureIfWanted(gmx_mtop_t *mtop, rvec *x, int ePBC, const matrix box)
1203 char *env = getenv("GMX_COMPELDUMP");
1205 if (env != nullptr)
1207 fprintf(stderr, "\n%s Found env.var. GMX_COMPELDUMP, will output CompEL starting structure made whole.\n"
1208 "%s In case of multimeric channels, please check whether they have the correct PBC representation.\n",
1209 SwS, SwSEmpty);
1211 write_sto_conf_mtop("CompELAssumedWholeConfiguration.pdb", *mtop->name, mtop, x, nullptr, ePBC, box);
1216 /*! \brief Initialize the swapstate structure, used for checkpoint writing.
1218 * The swapstate struct stores the information we need to make the channels
1219 * whole again after restarts from a checkpoint file. Here we do the following:
1220 * a) If we did not start from .cpt, we prepare the struct for proper .cpt writing,
1221 * b) if we did start from .cpt, we copy over the last whole structures from .cpt,
1222 * c) in any case, for subsequent checkpoint writing, we set the pointers in
1223 * swapstate to the x_old arrays, which contain the correct PBC representation of
1224 * multimeric channels at the last time step.
1226 static void init_swapstate(
1227 swaphistory_t *swapstate,
1228 t_swapcoords *sc,
1229 gmx_mtop_t *mtop,
1230 const rvec *x, /* the initial positions */
1231 const matrix box,
1232 t_inputrec *ir)
1234 rvec *x_pbc = nullptr; /* positions of the whole MD system with molecules made whole */
1235 t_swapgrp *g;
1236 t_swap *s;
1239 s = sc->si_priv;
1241 /* We always need the last whole positions such that
1242 * in the next time step we can make the channels whole again in PBC */
1243 if (swapstate->bFromCpt)
1245 /* Copy the last whole positions of each channel from .cpt */
1246 g = &(s->group[eGrpSplit0]);
1247 for (int i = 0; i < g->nat; i++)
1249 copy_rvec(swapstate->xc_old_whole[eChan0][i], g->xc_old[i]);
1251 g = &(s->group[eGrpSplit1]);
1252 for (int i = 0; i < g->nat; i++)
1254 copy_rvec(swapstate->xc_old_whole[eChan1][i], g->xc_old[i]);
1257 else
1259 swapstate->eSwapCoords = ir->eSwapCoords;
1261 /* Set the number of ion types and allocate memory for checkpointing */
1262 swapstate->nIonTypes = s->ngrp - eSwapFixedGrpNR;
1263 snew(swapstate->ionType, swapstate->nIonTypes);
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++)
1269 swapstateIons_t *gs = &swapstate->ionType[ii];
1270 gs->nMol = sc->grp[ii + eSwapFixedGrpNR].nat;
1273 /* Extract the initial split group positions. */
1275 /* Remove pbc, make molecule whole. */
1276 snew(x_pbc, mtop->natoms);
1277 copy_rvecn(x, x_pbc, 0, mtop->natoms);
1279 /* This can only make individual molecules whole, not multimers */
1280 do_pbc_mtop(nullptr, ir->ePBC, box, mtop, x_pbc);
1282 /* Output the starting structure? */
1283 outputStartStructureIfWanted(mtop, x_pbc, ir->ePBC, box);
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++)
1289 g = &(s->group[ig]);
1290 for (int i = 0; i < g->nat; i++)
1292 copy_rvec(x_pbc[g->ind[i]], g->xc_old[i]);
1295 sfree(x_pbc);
1297 /* Prepare swapstate arrays for later checkpoint writing */
1298 swapstate->nat[eChan0] = s->group[eGrpSplit0].nat;
1299 swapstate->nat[eChan1] = s->group[eGrpSplit1].nat;
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;
1308 /*! \brief Determine the total charge imbalance resulting from the swap groups */
1309 static real getRequestedChargeImbalance(t_swap *s)
1311 int ig;
1312 real DeltaQ = 0.0;
1313 t_swapgrp *g;
1314 real particle_charge;
1315 real particle_number[eCompNR];
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 );
1320 for (ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
1322 g = &s->group[ig];
1324 particle_charge = g->q;
1325 particle_number[eCompA] = g->comp[eCompA].nMolReq;
1326 particle_number[eCompB] = g->comp[eCompB].nMolReq;
1328 DeltaQ += particle_charge * (particle_number[eCompA] - particle_number[eCompB]);
1331 return DeltaQ;
1335 /*! \brief Sorts anions and cations into two separate groups
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.
1340 static void copyIndicesToGroup(
1341 int *indIons,
1342 int nIons,
1343 t_swapGroup *g,
1344 t_commrec *cr)
1346 g->nat = nIons;
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) )
1353 if (g->nat != (g->nmolReq[eCompA] + g->nmolReq[eCompB]) )
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);
1363 srenew(g->ind, g->nat);
1364 for (int i = 0; i < g->nat; i++)
1366 g->ind[i] = indIons[i];
1371 /*! \brief Converts old .tpr file CompEL contents to new data layout.
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
1384 static void convertOldToNewGroupFormat(
1385 t_swapcoords *sc,
1386 gmx_mtop_t *mtop,
1387 gmx_bool bVerbose,
1388 t_commrec *cr)
1390 t_swapGroup *g = &sc->grp[3];
1392 /* Loop through the atom indices of group #3 (anions) and put all indices
1393 * that belong to cations into the cation group.
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);
1402 int molb = 0;
1403 for (int i = 0; i < g->nat; i++)
1405 const t_atom &atom = mtopGetAtomParameters(mtop, g->ind[i], &molb);
1406 if (atom.q < 0)
1408 // This is an anion, add it to the list of anions
1409 indAnions[nAnions++] = g->ind[i];
1411 else
1413 // This is a cation, add it to the list of cations
1414 indCations[nCations++] = g->ind[i];
1418 if (bVerbose)
1420 fprintf(stdout, "%s Sorted %d ions into separate groups of %d anions and %d cations.\n",
1421 SwS, g->nat, nAnions, nCations);
1425 /* Now we have the correct lists of anions and cations.
1426 * Copy it to the right groups.
1428 copyIndicesToGroup(indAnions, nAnions, g, cr);
1429 g = &sc->grp[4];
1430 copyIndicesToGroup(indCations, nCations, g, cr);
1431 sfree(indAnions);
1432 sfree(indCations);
1434 return;
1438 /*! \brief Returns TRUE if we started from an old .tpr
1440 * Then we need to re-sort anions and cations into separate groups */
1441 static gmx_bool bConvertFromOldTpr(t_swapcoords *sc)
1443 // If the last group has no atoms it means we need to convert!
1444 if ( (sc->ngrp >= 5) && (0 == sc->grp[4].nat) )
1446 return TRUE;
1448 return FALSE;
1452 void init_swapcoords(
1453 FILE *fplog,
1454 t_inputrec *ir,
1455 const char *fn,
1456 gmx_mtop_t *mtop,
1457 const t_state *globalState,
1458 ObservablesHistory *oh,
1459 t_commrec *cr,
1460 const gmx_output_env_t *oenv,
1461 const MdrunOptions &mdrunOptions)
1463 t_swapcoords *sc;
1464 t_swap *s;
1465 t_swapgrp *g;
1466 swapstateIons_t *gs;
1467 gmx_bool bAppend, bStartFromCpt;
1468 swaphistory_t *swapstate = nullptr;
1470 if ( (PAR(cr)) && !DOMAINDECOMP(cr) )
1472 gmx_fatal(FARGS, "Position swapping is only implemented for domain decomposition!");
1475 bAppend = mdrunOptions.continuationOptions.appendFiles;
1476 bStartFromCpt = mdrunOptions.continuationOptions.startedFromCheckpoint;
1478 sc = ir->swap;
1479 snew(sc->si_priv, 1);
1480 s = sc->si_priv;
1482 if (mdrunOptions.rerun)
1484 if (PAR(cr))
1486 gmx_fatal(FARGS, "%s This module does not support reruns in parallel\nPlease request a serial run with -nt 1 / -np 1\n", SwS);
1489 fprintf(stderr, "%s Rerun - using every available frame\n", SwS);
1490 sc->nstswap = 1;
1491 sc->nAverage = 1; /* averaging makes no sense for reruns */
1494 if (MASTER(cr) && !bAppend)
1496 fprintf(fplog, "\nInitializing ion/water position exchanges\n");
1497 please_cite(fplog, "Kutzner2011b");
1500 switch (ir->eSwapCoords)
1502 case eswapX:
1503 s->swapdim = XX;
1504 break;
1505 case eswapY:
1506 s->swapdim = YY;
1507 break;
1508 case eswapZ:
1509 s->swapdim = ZZ;
1510 break;
1511 default:
1512 s->swapdim = -1;
1513 break;
1516 const gmx_bool bVerbose = mdrunOptions.verbose;
1518 // For compatibility with old .tpr files
1519 if (bConvertFromOldTpr(sc) )
1521 convertOldToNewGroupFormat(sc, mtop, bVerbose && MASTER(cr), cr);
1524 /* Copy some data and pointers to the group structures for convenience */
1525 /* Number of atoms in the group */
1526 s->ngrp = sc->ngrp;
1527 snew(s->group, s->ngrp);
1528 for (int i = 0; i < s->ngrp; i++)
1530 s->group[i].nat = sc->grp[i].nat;
1531 s->group[i].ind = sc->grp[i].ind;
1532 s->group[i].molname = sc->grp[i].molname;
1535 /* Check for overlapping atoms */
1536 check_swap_groups(s, mtop->natoms, bVerbose && MASTER(cr));
1538 /* Allocate space for the collective arrays for all groups */
1539 /* For the collective position array */
1540 for (int i = 0; i < s->ngrp; i++)
1542 g = &s->group[i];
1543 snew(g->xc, g->nat);
1544 snew(g->c_ind_loc, g->nat);
1546 /* For the split groups (the channels) we need some extra memory to
1547 * be able to make the molecules whole even if they span more than
1548 * half of the box size. */
1549 if ( (i == eGrpSplit0) || (i == eGrpSplit1) )
1551 snew(g->xc_shifts, g->nat);
1552 snew(g->xc_eshifts, g->nat);
1553 snew(g->xc_old, g->nat);
1557 if (MASTER(cr))
1559 if (oh->swapHistory == nullptr)
1561 oh->swapHistory = std::unique_ptr<swaphistory_t>(new swaphistory_t {});
1563 swapstate = oh->swapHistory.get();
1565 init_swapstate(swapstate, sc, mtop, as_rvec_array(globalState->x.data()), globalState->box, ir);
1568 /* After init_swapstate we have a set of (old) whole positions for our
1569 * channels. Now transfer that to all nodes */
1570 if (PAR(cr))
1572 for (int ig = eGrpSplit0; ig <= eGrpSplit1; ig++)
1574 g = &(s->group[ig]);
1575 gmx_bcast((g->nat)*sizeof((g->xc_old)[0]), g->xc_old, (cr));
1579 /* Make sure that all molecules in the solvent and ion groups contain the
1580 * same number of atoms each */
1581 for (int ig = eGrpSolvent; ig < s->ngrp; ig++)
1583 real charge;
1585 g = &(s->group[ig]);
1586 g->apm = get_group_apm_check(ig, s, MASTER(cr) && bVerbose, mtop);
1588 /* Since all molecules of a group are equal, we only need enough space
1589 * to determine properties of a single molecule at at time */
1590 snew(g->m, g->apm); /* For the center of mass */
1591 charge = 0; /* To determine the charge imbalance */
1592 int molb = 0;
1593 for (int j = 0; j < g->apm; j++)
1595 const t_atom &atom = mtopGetAtomParameters(mtop, g->ind[j], &molb);
1596 g->m[j] = atom.m;
1597 charge += atom.q;
1599 /* Total charge of one molecule of this group: */
1600 g->q = charge;
1604 /* Need mass-weighted center of split group? */
1605 for (int j = eGrpSplit0; j <= eGrpSplit1; j++)
1607 g = &(s->group[j]);
1608 if (TRUE == sc->massw_split[j])
1610 /* Save the split group masses if mass-weighting is requested */
1611 snew(g->m, g->nat);
1612 int molb = 0;
1613 for (int i = 0; i < g->nat; i++)
1615 g->m[i] = mtopGetAtomMass(mtop, g->ind[i], &molb);
1620 /* Make a t_pbc struct on all nodes so that the molecules
1621 * chosen for an exchange can be made whole. */
1622 snew(s->pbc, 1);
1624 if (MASTER(cr))
1626 if (bVerbose)
1628 fprintf(stderr, "%s Opening output file %s%s\n", SwS, fn, bAppend ? " for appending" : "");
1631 s->fpout = gmx_fio_fopen(fn, bAppend ? "a" : "w" );
1633 if (!bAppend)
1635 xvgr_header(s->fpout, "Molecule counts", "Time (ps)", "counts", exvggtXNY, oenv);
1637 for (int ig = 0; ig < s->ngrp; ig++)
1639 g = &(s->group[ig]);
1640 fprintf(s->fpout, "# %s group '%s' contains %d atom%s",
1641 ig < eSwapFixedGrpNR ? eSwapFixedGrp_names[ig] : "Ion",
1642 g->molname, g->nat, (g->nat > 1) ? "s" : "");
1643 if (!(eGrpSplit0 == ig || eGrpSplit1 == ig) )
1645 fprintf(s->fpout, " with %d atom%s in each molecule of charge %g",
1646 g->apm, (g->apm > 1) ? "s" : "", g->q);
1648 fprintf(s->fpout, ".\n");
1651 fprintf(s->fpout, "#\n# Initial positions of split groups:\n");
1654 for (int j = eGrpSplit0; j <= eGrpSplit1; j++)
1656 g = &(s->group[j]);
1657 for (int i = 0; i < g->nat; i++)
1659 copy_rvec(globalState->x[sc->grp[j].ind[i]], g->xc[i]);
1661 /* xc has the correct PBC representation for the two channels, so we do
1662 * not need to correct for that */
1663 get_center(g->xc, g->m, g->nat, g->center);
1664 if (!bAppend)
1666 fprintf(s->fpout, "# %s group %s-center %5f nm\n", eSwapFixedGrp_names[j],
1667 DimStr[s->swapdim], g->center[s->swapdim]);
1671 if (!bAppend)
1673 if ( (0 != sc->bulkOffset[eCompA]) || (0 != sc->bulkOffset[eCompB]) )
1675 fprintf(s->fpout, "#\n");
1676 fprintf(s->fpout, "# You provided an offset for the position of the bulk layer(s).\n");
1677 fprintf(s->fpout, "# That means the layers to/from which ions and water molecules are swapped\n");
1678 fprintf(s->fpout, "# are not midway (= at 0.0) between the compartment-defining layers (at +/- 1.0).\n");
1679 fprintf(s->fpout, "# bulk-offsetA = %g\n", sc->bulkOffset[eCompA]);
1680 fprintf(s->fpout, "# bulk-offsetB = %g\n", sc->bulkOffset[eCompB]);
1683 fprintf(s->fpout, "#\n");
1684 fprintf(s->fpout, "# Split0 cylinder radius %f nm, up %f nm, down %f nm\n",
1685 sc->cyl0r, sc->cyl0u, sc->cyl0l);
1686 fprintf(s->fpout, "# Split1 cylinder radius %f nm, up %f nm, down %f nm\n",
1687 sc->cyl1r, sc->cyl1u, sc->cyl1l);
1689 fprintf(s->fpout, "#\n");
1690 if (!mdrunOptions.rerun)
1692 fprintf(s->fpout, "# Coupling constant (number of swap attempt steps to average over): %d (translates to %f ps).\n",
1693 sc->nAverage, sc->nAverage*sc->nstswap*ir->delta_t);
1694 fprintf(s->fpout, "# Threshold is %f\n", sc->threshold);
1695 fprintf(s->fpout, "#\n");
1696 fprintf(s->fpout, "# Remarks about which atoms passed which channel use global atoms numbers starting at one.\n");
1700 else
1702 s->fpout = nullptr;
1705 /* Prepare for parallel or serial run */
1706 if (PAR(cr))
1708 for (int ig = 0; ig < s->ngrp; ig++)
1710 g = &(s->group[ig]);
1711 g->nat_loc = 0;
1712 g->nalloc_loc = 0;
1713 g->ind_loc = nullptr;
1716 else
1718 for (int ig = 0; ig < s->ngrp; ig++)
1720 g = &(s->group[ig]);
1721 g->nat_loc = g->nat;
1722 g->ind_loc = g->ind;
1723 /* c_ind_loc needs to be set to identity in the serial case */
1724 for (int i = 0; i < g->nat; i++)
1726 g->c_ind_loc[i] = i;
1731 /* Allocate memory to remember the past particle counts for time averaging */
1732 for (int ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
1734 g = &(s->group[ig]);
1735 for (int ic = 0; ic < eCompNR; ic++)
1737 snew(g->comp[ic].nMolPast, sc->nAverage);
1741 /* Get the initial particle concentrations and let the other nodes know */
1742 if (MASTER(cr))
1744 if (bStartFromCpt)
1746 get_initial_ioncounts_from_cpt(ir, swapstate, cr, bVerbose);
1748 else
1750 fprintf(stderr, "%s Determining initial numbers of ions per compartment.\n", SwS);
1751 get_initial_ioncounts(ir, as_rvec_array(globalState->x.data()), globalState->box, cr, mdrunOptions.rerun);
1754 /* Prepare (further) checkpoint writes ... */
1755 if (bStartFromCpt)
1757 /* Consistency check */
1758 if (swapstate->nAverage != sc->nAverage)
1760 gmx_fatal(FARGS, "%s Ion count averaging steps mismatch! checkpoint: %d, tpr: %d",
1761 SwS, swapstate->nAverage, sc->nAverage);
1764 else
1766 swapstate->nAverage = sc->nAverage;
1768 fprintf(stderr, "%s Setting pointers for checkpoint writing\n", SwS);
1769 for (int ic = 0; ic < eCompNR; ic++)
1771 for (int ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
1773 g = &s->group[ig];
1774 gs = &swapstate->ionType[ig - eSwapFixedGrpNR];
1776 gs->nMolReq_p[ic] = &(g->comp[ic].nMolReq);
1777 gs->nMolPast_p[ic] = &(g->comp[ic].nMolPast[0]);
1778 gs->inflow_net_p[ic] = &(g->comp[ic].inflow_net);
1782 /* Determine the total charge imbalance */
1783 s->deltaQ = getRequestedChargeImbalance(s);
1785 if (bVerbose)
1787 fprintf(stderr, "%s Requested charge imbalance is Q(A) - Q(B) = %g e.\n", SwS, s->deltaQ);
1789 if (!bAppend)
1791 fprintf(s->fpout, "# Requested charge imbalance is Q(A)-Q(B) = %g e.\n", s->deltaQ);
1795 if (PAR(cr))
1797 bc_initial_concentrations(cr, ir->swap);
1800 /* Update the time-averaged number of molecules for all groups and compartments */
1801 for (int ig = eSwapFixedGrpNR; ig < sc->ngrp; ig++)
1803 g = &s->group[ig];
1804 for (int ic = 0; ic < eCompNR; ic++)
1806 update_time_window(&g->comp[ic], sc->nAverage, -1);
1810 /* Initialize arrays that keep track of through which channel the ions go */
1811 detect_flux_per_channel_init(s, swapstate, bStartFromCpt);
1813 /* We need to print the legend if we open this file for the first time. */
1814 if (MASTER(cr) && !bAppend)
1816 print_ionlist_legend(ir, oenv);
1821 void finish_swapcoords(t_swapcoords *sc)
1823 if (sc->si_priv->fpout)
1825 // Close the swap output file
1826 gmx_fio_fclose(sc->si_priv->fpout);
1831 void dd_make_local_swap_groups(gmx_domdec_t *dd, t_swapcoords *sc)
1833 t_swapgrp *g;
1834 int ig;
1837 /* Make split groups, solvent group, and user-defined groups of particles
1838 * under control of the swap protocol */
1839 for (ig = 0; ig < sc->ngrp; ig++)
1841 g = &(sc->si_priv->group[ig]);
1842 dd_make_local_group_indices(dd->ga2la, g->nat, g->ind,
1843 &(g->nat_loc), &(g->ind_loc), &(g->nalloc_loc), g->c_ind_loc);
1848 /*! \brief Do we need to swap a molecule in any of the ion groups with a water molecule at this step?
1850 * From the requested and average molecule counts we determine whether a swap is needed
1851 * at this time step.
1853 static gmx_bool need_swap(t_swapcoords *sc)
1855 t_swap *s;
1856 int ic, ig;
1857 t_swapgrp *g;
1859 s = sc->si_priv;
1861 for (ig = eSwapFixedGrpNR; ig < sc->ngrp; ig++)
1863 g = &s->group[ig];
1865 for (ic = 0; ic < eCompNR; ic++)
1867 if (g->comp[ic].nMolReq - g->comp[ic].nMolAv >= sc->threshold)
1869 return TRUE;
1873 return FALSE;
1877 /*! \brief Return the index of an atom or molecule suitable for swapping.
1879 * Returns the index of an atom that is far off the compartment boundaries,
1880 * that is near to the bulk layer to/from which the swaps take place.
1881 * Other atoms of the molecule (if any) will directly follow the returned index.
1883 * \param[in] comp Structure containing compartment-specific data.
1884 * \param[in] molname Name of the molecule.
1886 * \returns Index of the first atom of the molecule chosen for a position exchange.
1888 static int get_index_of_distant_atom(
1889 t_compartment *comp,
1890 const char molname[])
1892 int ibest = -1;
1893 real d = GMX_REAL_MAX;
1896 /* comp->nat contains the original number of atoms in this compartment
1897 * prior to doing any swaps. Some of these atoms may already have been
1898 * swapped out, but then they are marked with a distance of GMX_REAL_MAX
1900 for (int iMol = 0; iMol < comp->nMolBefore; iMol++)
1902 if (comp->dist[iMol] < d)
1904 ibest = iMol;
1905 d = comp->dist[ibest];
1909 if (ibest < 0)
1911 gmx_fatal(FARGS, "Could not get index of %s atom. Compartment contains %d %s molecules before swaps.",
1912 molname, comp->nMolBefore, molname);
1915 /* Set the distance of this index to infinity such that it won't get selected again in
1916 * this time step
1918 comp->dist[ibest] = GMX_REAL_MAX;
1920 return comp->ind[ibest];
1924 /*! \brief Swaps centers of mass and makes molecule whole if broken */
1925 static void translate_positions(
1926 rvec *x,
1927 int apm,
1928 rvec old_com,
1929 rvec new_com,
1930 t_pbc *pbc)
1932 int i;
1933 rvec reference, dx, correctPBCimage;
1936 /* Use the first atom as the reference for PBC */
1937 copy_rvec(x[0], reference);
1939 for (i = 0; i < apm; i++)
1941 /* PBC distance between position and reference */
1942 pbc_dx(pbc, x[i], reference, dx);
1944 /* Add PBC distance to reference */
1945 rvec_add(reference, dx, correctPBCimage);
1947 /* Subtract old_com from correct image and add new_com */
1948 rvec_dec(correctPBCimage, old_com);
1949 rvec_inc(correctPBCimage, new_com);
1951 copy_rvec(correctPBCimage, x[i]);
1956 /*! \brief Write back the the modified local positions from the collective array to the official positions. */
1957 static void apply_modified_positions(
1958 t_swapgrp *g,
1959 rvec x[])
1961 int l, ii, cind;
1964 for (l = 0; l < g->nat_loc; l++)
1966 /* Get the right local index to write to */
1967 ii = g->ind_loc[l];
1968 /* Where is the local atom in the collective array? */
1969 cind = g->c_ind_loc[l];
1971 /* Copy the possibly modified position */
1972 copy_rvec(g->xc[cind], x[ii]);
1977 gmx_bool do_swapcoords(
1978 t_commrec *cr,
1979 gmx_int64_t step,
1980 double t,
1981 t_inputrec *ir,
1982 gmx_wallcycle *wcycle,
1983 rvec x[],
1984 matrix box,
1985 gmx_bool bVerbose,
1986 gmx_bool bRerun)
1988 t_swapcoords *sc;
1989 t_swap *s;
1990 int j, ic, ig, nswaps;
1991 int thisC, otherC; /* Index into this compartment and the other one */
1992 gmx_bool bSwap = FALSE;
1993 t_swapgrp *g, *gsol;
1994 int isol, iion;
1995 rvec com_solvent, com_particle; /* solvent and swap molecule's center of mass */
1998 wallcycle_start(wcycle, ewcSWAP);
2000 sc = ir->swap;
2001 s = sc->si_priv;
2003 set_pbc(s->pbc, ir->ePBC, box);
2005 /* Assemble the positions of the split groups, i.e. the channels.
2006 * Here we also pass a shifts array to communicate_group_positions(), so that it can make
2007 * the molecules whole even in cases where they span more than half of the box in
2008 * any dimension */
2009 for (ig = eGrpSplit0; ig <= eGrpSplit1; ig++)
2011 g = &(s->group[ig]);
2012 communicate_group_positions(cr, g->xc, g->xc_shifts, g->xc_eshifts, TRUE,
2013 x, g->nat, g->nat_loc, g->ind_loc, g->c_ind_loc, g->xc_old, box);
2015 get_center(g->xc, g->m, g->nat, g->center); /* center of split groups == channels */
2018 /* Assemble the positions of the ions (ig = 3, 4, ...). These molecules should
2019 * be small and we can always make them whole with a simple distance check.
2020 * Therefore we pass NULL as third argument. */
2021 for (ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
2023 g = &(s->group[ig]);
2024 communicate_group_positions(cr, g->xc, nullptr, nullptr, FALSE,
2025 x, g->nat, g->nat_loc, g->ind_loc, g->c_ind_loc, nullptr, nullptr);
2027 /* Determine how many ions of this type each compartment contains */
2028 sortMoleculesIntoCompartments(g, cr, sc, box, step, s->fpout, bRerun, FALSE);
2031 /* Output how many ions are in the compartments */
2032 if (MASTER(cr))
2034 print_ionlist(s, t, "");
2037 /* If we are doing a rerun, we are finished here, since we cannot perform
2038 * swaps anyway */
2039 if (bRerun)
2041 return FALSE;
2044 /* Do we have to perform a swap? */
2045 bSwap = need_swap(sc);
2046 if (bSwap)
2048 /* Since we here know that we have to perform ion/water position exchanges,
2049 * we now assemble the solvent positions */
2050 g = &(s->group[eGrpSolvent]);
2051 communicate_group_positions(cr, g->xc, nullptr, nullptr, FALSE,
2052 x, g->nat, g->nat_loc, g->ind_loc, g->c_ind_loc, nullptr, nullptr);
2054 /* Determine how many molecules of solvent each compartment contains */
2055 sortMoleculesIntoCompartments(g, cr, sc, box, step, s->fpout, bRerun, TRUE);
2057 /* Save number of solvent molecules per compartment prior to any swaps */
2058 g->comp[eCompA].nMolBefore = g->comp[eCompA].nMol;
2059 g->comp[eCompB].nMolBefore = g->comp[eCompB].nMol;
2061 for (ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
2063 g = &(s->group[ig]);
2065 for (ic = 0; ic < eCompNR; ic++)
2067 /* Determine in which compartment ions are missing and where they are too many */
2068 g->vacancy[ic] = g->comp[ic].nMolReq - g->comp[ic].nMolAv;
2070 /* Save number of ions per compartment prior to swaps */
2071 g->comp[ic].nMolBefore = g->comp[ic].nMol;
2075 /* Now actually perform the particle exchanges, one swap group after another */
2076 gsol = &s->group[eGrpSolvent];
2077 for (ig = eSwapFixedGrpNR; ig < s->ngrp; ig++)
2079 nswaps = 0;
2080 g = &s->group[ig];
2081 for (thisC = 0; thisC < eCompNR; thisC++)
2083 /* Index to the other compartment */
2084 otherC = (thisC+1) % eCompNR;
2086 while (g->vacancy[thisC] >= sc->threshold)
2088 /* Swap in an ion */
2090 /* Get the xc-index of the first atom of a solvent molecule of this compartment */
2091 isol = get_index_of_distant_atom(&gsol->comp[thisC], gsol->molname);
2093 /* Get the xc-index of a particle from the other compartment */
2094 iion = get_index_of_distant_atom(&g->comp[otherC], g->molname);
2096 get_molecule_center(&gsol->xc[isol], gsol->apm, gsol->m, com_solvent, s->pbc);
2097 get_molecule_center(&g->xc[iion], g->apm, g->m, com_particle, s->pbc);
2099 /* Subtract solvent molecule's center of mass and add swap particle's center of mass */
2100 translate_positions(&gsol->xc[isol], gsol->apm, com_solvent, com_particle, s->pbc);
2101 /* Similarly for the swap particle, subtract com_particle and add com_solvent */
2102 translate_positions(&g->xc[iion], g->apm, com_particle, com_solvent, s->pbc);
2104 /* Keep track of the changes */
2105 g->vacancy[thisC ]--;
2106 g->vacancy[otherC]++;
2107 g->comp [thisC ].nMol++;
2108 g->comp [otherC].nMol--;
2109 g->comp [thisC ].inflow_net++;
2110 g->comp [otherC].inflow_net--;
2111 /* Correct the past time window to still get the right averages from now on */
2112 g->comp [thisC ].nMolAv++;
2113 g->comp [otherC].nMolAv--;
2114 for (j = 0; j < sc->nAverage; j++)
2116 g->comp[thisC ].nMolPast[j]++;
2117 g->comp[otherC].nMolPast[j]--;
2119 /* Clear ion history */
2120 if (MASTER(cr))
2122 int iMol = iion / g->apm;
2123 g->channel_label[iMol] = eChHistPassedNone;
2124 g->comp_from[iMol] = eDomainNotset;
2126 /* That was the swap */
2127 nswaps++;
2131 if (nswaps && bVerbose)
2133 fprintf(stderr, "%s Performed %d swap%s in step %" GMX_PRId64 " for iontype %s.\n",
2134 SwS, nswaps, nswaps > 1 ? "s" : "", step, g->molname);
2138 if (s->fpout != nullptr)
2140 print_ionlist(s, t, " # after swap");
2143 /* For the solvent and user-defined swap groups, each rank writes back its
2144 * (possibly modified) local positions to the official position array. */
2145 for (ig = eGrpSolvent; ig < s->ngrp; ig++)
2147 g = &s->group[ig];
2148 apply_modified_positions(g, x);
2151 } /* end of if(bSwap) */
2153 wallcycle_stop(wcycle, ewcSWAP);
2155 return bSwap;