2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2013,2014, 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.
37 * Implements functions in swapcoords.h.
39 * \author Carsten Kutzner <ckutzne@gwdg.de>
40 * \ingroup module_swap
51 #include "gromacs/legacyheaders/typedefs.h"
52 #include "gromacs/utility/cstringutil.h"
53 #include "gromacs/utility/smalloc.h"
54 #include "gromacs/mdlib/groupcoord.h"
55 #include "gromacs/topology/mtop_util.h"
56 #include "gromacs/legacyheaders/macros.h"
57 #include "gromacs/math/vec.h"
58 #include "gromacs/legacyheaders/names.h"
59 #include "gromacs/legacyheaders/network.h"
60 #include "gromacs/legacyheaders/mdrun.h"
61 #include "gromacs/fileio/xvgr.h"
62 #include "gromacs/legacyheaders/copyrite.h"
63 #include "gromacs/fileio/confio.h"
64 #include "gromacs/timing/wallcycle.h"
65 #include "swapcoords.h"
67 #include "gromacs/pbcutil/pbc.h"
69 static char *SwS = {"SWAP:"}; /**< For output that comes from the swap module */
70 static char *SwSEmpty = {" "}; /**< Placeholder for multi-line output */
71 static char* IonString[eIonNR] = {"anion", "cation" }; /**< Type of ion, used for verbose output */
72 static char* IonStr[eIonNR] = {"-", "+" }; /**< Type of ion, used for short output */
73 static char* CompStr[eCompNR] = {"A", "B" }; /**< Compartment name */
74 static char *SwapStr[eSwapTypesNR+1] = { "", "X-", "Y-", "Z-", NULL}; /**< Name for the swap types. */
75 static char *DimStr[DIM+1] = { "X", "Y", "Z", NULL}; /**< Name for the swap dimension. */
77 /* eGrpSplit0 and eGrpSplit1 _must_ be neighbors in this list because
78 * we sometimes loop from eGrpSplit0 to eGrpSplit1 */
80 eGrpIons, eGrpSplit0, eGrpSplit1, eGrpSolvent, eGrpNr
81 }; /**< Group identifier */
82 static char* GrpString[eGrpNr] = { "ion", "split0", "split1", "solvent" }; /**< Group name */
84 /** Keep track of through which channel the ions have passed */
85 enum eChannelHistory {
86 eChHistPassedNone, eChHistPassedCh0, eChHistPassedCh1, eChHistNr
88 static char* ChannelString[eChHistNr] = { "none", "channel0", "channel1" }; /**< Name for the channels */
90 /*! \brief Domain identifier.
92 * Keeps track of from which compartment the ions came before passing the
96 eDomainNotset, eDomainA, eDomainB, eDomainNr
98 static char* DomainString[eDomainNr] = { "not_assigned", "Domain_A", "Domain_B" }; /**< Name for the domains */
103 * Structure containing compartment-specific data.
105 typedef struct swap_compartment
107 int nat; /**< Number of atoms matching the
108 compartment conditions. */
109 int nat_old; /**< Number of atoms before swapping. */
110 int nat_req; /**< Requested number of atoms. */
111 real nat_av; /**< Time-averaged number of atoms matching
112 the compartment conditions. */
113 int *nat_past; /**< Past ion counts for time-averaging. */
114 int *ind; /**< Indices to coll array of atoms. */
115 real *dist; /**< Distance of atom to compartment center. */
116 int nalloc; /**< Allocation size for ind array. */
117 int inflow_netto; /**< Net inflow of ions into this compartment. */
122 * This structure contains data needed for each of the groups involved in swapping: ions, water,
125 typedef struct swap_group
127 int nat; /**< Number of atoms in the group */
128 int apm; /**< Number of atoms in each molecule */
129 atom_id *ind; /**< Global atom indices of the group */
130 atom_id *ind_loc; /**< Local atom indices of the group */
131 int nat_loc; /**< Number of local group atoms */
132 int nalloc_loc; /**< Allocation size for ind_loc */
133 rvec *xc; /**< Collective array of group atom positions */
134 ivec *xc_shifts; /**< Current (collective) shifts */
135 ivec *xc_eshifts; /**< Extra shifts since last DD step */
136 rvec *xc_old; /**< Old (collective) positions */
137 real *qc; /**< Collective array of charges */
138 int *c_ind_loc; /**< Position of local atoms in the
139 collective array, [0..nat_loc] */
140 real *m; /**< Masses (can be omitted) */
141 unsigned char *comp_from; /**< (Collective) Stores from which compartment this
142 atom has come. This way we keep track of through
143 which channel an ion permeates (only used for
145 unsigned char *comp_now; /**< In which compartment this ion is now */
146 unsigned char *channel_label; /**< Which channel was passed at last by this ion? */
147 rvec center; /**< Center of the group; COM if masses are used */
152 * Main (private) data structure for the position swapping protocol.
156 int swapdim; /**< One of XX, YY, ZZ */
157 t_pbc *pbc; /**< Needed to make molecules whole. */
158 FILE *fpout; /**< Output file. */
159 t_group group[eGrpNr]; /**< Ions, solvent or channels? */
160 t_compartment comp[eCompNR][eIonNR]; /**< Data for a specific compartment and ion type. */
161 t_compartment compsol[eCompNR]; /**< Solvent compartments. */
162 int fluxfromAtoB[eChanNR][eIonNR]; /**< Net flux per channels and ion type. */
163 int ncyl0ions; /**< Number of ions residing in channel 0. */
164 int ncyl1ions; /**< Same for channel 1. */
165 int cyl0and1; /**< Ions assigned to cyl0 and cyl1. Not good. */
166 int *fluxleak; /**< Pointer to a single int value holding the
167 flux not going through any of the channels. */
168 real deltaQ; /**< The charge imbalance between the compartments. */
173 /*! \brief Check whether point is in channel.
175 * A channel is a cylinder defined by a disc
176 * with radius r around its center c. The thickness of the cylinder is
183 * <---------c--------->
189 static gmx_bool is_in_channel(
190 rvec point, /* Point under consideration */
191 rvec center, /* 'Center' of cylinder */
192 real d_up, /* Upper extension */
193 real d_down, /* Lower extensions */
194 real r_cyl2, /* Cylinder radius squared */
196 int normal) /* The membrane normal direction is typically 3, i.e. ZZ, but can be X or Y also */
199 int plane1, plane2; /* Directions tangential to membrane */
202 plane1 = (normal + 1) % 3; /* typically 0, i.e. XX */
203 plane2 = (normal + 2) % 3; /* typically 1, i.e. YY */
205 /* Get the distance vector dr between the point and the center of the cylinder */
206 pbc_dx(pbc, point, center, dr); /* This puts center in the origin */
208 /* Check vertical direction */
209 if ( (dr[normal] > d_up) || (dr[normal] < -d_down) )
214 /* Check radial direction */
215 if ( (dr[plane1]*dr[plane1] + dr[plane2]*dr[plane2]) > r_cyl2)
220 /* All check passed, this point is in the cylinder */
225 /*! \brief Prints to swap output file which ions are in which compartment. */
226 static void print_ionlist(
231 int itype, icomp, i, j;
235 fprintf(s->fpout, "%12.5e", time);
236 for (icomp = 0; icomp < eCompNR; icomp++)
238 for (itype = 0; itype < eIonNR; itype++)
240 comp = &(s->comp[icomp][itype]);
241 fprintf(s->fpout, "%7d%7.1f%7d", comp->nat, comp->nat_av-comp->nat_req, comp->inflow_netto);
244 fprintf(s->fpout, "%12.3e%12.3e",
245 s->group[eGrpSplit0].center[s->swapdim],
246 s->group[eGrpSplit1].center[s->swapdim]);
248 for (i = 0; i < eChanNR; i++)
250 for (j = 0; j < eIonNR; j++)
252 fprintf(s->fpout, "%12d", s->fluxfromAtoB[i][j]);
256 /* Also print the number of ions that leaked from A to B: */
257 fprintf(s->fpout, "%12d", *s->fluxleak);
259 fprintf(s->fpout, "%s\n", comment);
263 /*! \brief Get the center of a group of nat atoms.
265 * Since with PBC an atom group might not be whole, use the first atom as the
266 * reference atom and determine the center with respect to this reference.
268 static void get_molecule_center(
276 rvec weightedPBCimage;
278 rvec reference, correctPBCimage, dx;
281 /* Use the first atom as the reference and put other atoms near that one */
282 /* This does not work for large molecules that span > half of the box! */
283 copy_rvec(x[0], reference);
285 /* Calculate either the weighted center or simply the center of geometry */
288 for (i = 0; i < nat; i++)
290 /* PBC distance between position and reference */
291 pbc_dx(pbc, x[i], reference, dx);
293 /* Add PBC distance to reference */
294 rvec_add(reference, dx, correctPBCimage);
296 /* Take weight into account */
306 svmul(wi, correctPBCimage, weightedPBCimage);
309 rvec_inc(center, weightedPBCimage);
313 svmul(1.0/wsum, center, center);
318 /*! \brief Return TRUE if ion is found in the compartment.
320 * Returns TRUE if x is between (w1+gap) and (w2-gap)
324 * ||-----------|--+--|----------o----------|--+--|---------------------||
325 * w1 ????????????????????? w2
329 static gmx_bool compartment_contains_atom(
330 real w1, /* position of wall atom 1 */
331 real w2, /* position of wall atom 2 */
334 real l, /* length of the box, from || to || in the sketch */
335 real *distance_from_center)
340 /* First set the origin in the middle of w1 and w2 */
346 /* Now choose the PBC image of x that is closest to the origin: */
357 *distance_from_center = (real)fabs(x);
359 /* Return TRUE if we now are in area "????" */
360 if ( (x >= (w1+gap)) && (x < (w2-gap)) )
371 /*! \brief Updates the time-averaged number of ions in a compartment. */
372 static void update_time_window(t_compartment *comp, int values, int replace)
378 /* Put in the new value */
381 comp->nat_past[replace] = comp->nat;
384 /* Compute the new time-average */
386 for (i = 0; i < values; i++)
388 average += comp->nat_past[i];
391 comp->nat_av = average;
395 /*! \brief Add atom with collective index ci to the list 'comp'. */
396 static void add_to_list(
397 int ci, /* index of this ion in the collective array xc, qc */
398 t_compartment *comp, /* Compartment to add this atom to */
399 real distance) /* Shortest distance of this atom to the compartment center */
406 if (nr >= comp->nalloc)
408 comp->nalloc = over_alloc_dd(nr+1);
409 srenew(comp->ind, comp->nalloc);
410 srenew(comp->dist, comp->nalloc);
413 comp->dist[nr] = distance;
418 /*! \brief Determine the compartment boundaries from the channel centers. */
419 static void get_compartment_boundaries(
423 real *left, real *right)
426 real leftpos, rightpos, leftpos_orig;
431 gmx_fatal(FARGS, "No compartment %d.", c);
434 pos0 = s->group[eGrpSplit0].center[s->swapdim];
435 pos1 = s->group[eGrpSplit1].center[s->swapdim];
448 /* This gets us the other compartment: */
451 leftpos_orig = leftpos;
453 rightpos = leftpos_orig + box[s->swapdim][s->swapdim];
461 /*! \brief Determine the per-channel ion flux.
463 * To determine the flux through the individual channels, we
464 * remember the compartment and channel history of each ion. An ion can be
465 * either in channel0 or channel1, or in the remaining volume of compartment
469 * +-----------------+
472 * ||||||||||0|||||||| bilayer with channel 0
477 * |||||1||||||||||||| bilayer with channel 1
480 * +-----------------+
484 static void detect_flux_per_channel(
489 unsigned char *comp_now,
490 unsigned char *comp_from,
491 unsigned char *channel_label,
501 gmx_bool in_cyl0, in_cyl1;
508 /* Check whether ion is inside any of the channels */
509 in_cyl0 = is_in_channel(ion_pos, s->group[eGrpSplit0].center, sc->cyl0u, sc->cyl0l, cyl0_r2, s->pbc, sd);
510 in_cyl1 = is_in_channel(ion_pos, s->group[eGrpSplit1].center, sc->cyl1u, sc->cyl1l, cyl1_r2, s->pbc, sd);
512 if (in_cyl0 && in_cyl1)
514 /* Ion appears to be in both channels. Something is severely wrong! */
516 *comp_now = eDomainNotset;
517 *comp_from = eDomainNotset;
518 *channel_label = eChHistPassedNone;
522 /* Ion is in channel 0 now */
523 *channel_label = eChHistPassedCh0;
524 *comp_now = eDomainNotset;
529 /* Ion is in channel 1 now */
530 *channel_label = eChHistPassedCh1;
531 *comp_now = eDomainNotset;
536 /* Ion is not in any of the channels, so it must be in domain A or B */
539 *comp_now = eDomainA;
543 *comp_now = eDomainB;
547 /* Only take action, if ion is now in domain A or B, and was before
548 * in the other domain!
550 if (eDomainNotset == *comp_from)
552 /* Maybe we can set the domain now */
553 *comp_from = *comp_now; /* Could still be eDomainNotset, though */
555 else if ( (*comp_now != eDomainNotset ) /* if in channel */
556 && (*comp_from != *comp_now) )
558 /* Obviously the ion changed its domain.
559 * Count this for the channel through which it has passed. */
560 switch (*channel_label)
562 case eChHistPassedNone:
563 *s->fluxleak = *s->fluxleak + 1;
565 fprintf(stderr, " %s Warning! Step %s, ion %d (%s) moved from %s to %s\n",
566 SwS, gmx_step_str(step, buf), iion, IonStr[iontype], DomainString[*comp_from], DomainString[*comp_now]);
569 fprintf(stderr, ", possibly due to a swap in the original simulation.\n");
573 fprintf(stderr, "but did not pass cyl0 or cyl1 as defined in the .mdp file.\n"
574 "Do you have an ion somewhere within the membrane?\n");
575 /* Write this info to the CompEL output file: */
576 fprintf(s->fpout, " # Warning: step %s, ion %d (%s) moved from %s to %s (probably through the membrane)\n",
577 gmx_step_str(step, buf), iion, IonStr[iontype],
578 DomainString[*comp_from], DomainString[*comp_now]);
582 case eChHistPassedCh0:
583 case eChHistPassedCh1:
584 if (*channel_label == eChHistPassedCh0)
593 if (eDomainA == *comp_from)
595 s->fluxfromAtoB[chan_nr][iontype]++;
599 s->fluxfromAtoB[chan_nr][iontype]--;
601 fprintf(fpout, "# Atom nr. %d finished passing %s.\n", iion, ChannelString[*channel_label]);
604 gmx_fatal(FARGS, "%s Unknown channel history entry!\n", SwS);
608 /* This ion has moved to the _other_ compartment ... */
609 *comp_from = *comp_now;
610 /* ... and it did not pass any channel yet */
611 *channel_label = eChHistPassedNone;
616 /*! \brief Get the lists of ions for the two compartments */
617 static void compartmentalize_ions(
630 real cyl0_r2, cyl1_r2;
632 int sum, not_in_comp[eCompNR]; /* consistency check */
637 iong = &s->group[eGrpIons];
640 cyl0_r2 = sc->cyl0r * sc->cyl0r;
641 cyl1_r2 = sc->cyl1r * sc->cyl1r;
644 /* Get us a counter that cycles in the range of [0 ... sc->nAverage[ */
645 replace = (step/sc->nstswap) % sc->nAverage;
647 for (comp = eCompA; comp <= eCompB; comp++)
649 /* Get lists of atoms that match criteria for this compartment */
650 get_compartment_boundaries(comp, sc->si_priv, box, &left, &right);
652 /* First clear the ion lists */
653 s->comp[comp][eIonNEG].nat = 0;
654 s->comp[comp][eIonPOS].nat = 0;
655 not_in_comp[comp] = 0; /* consistency check */
657 /* Loop over the IONS */
658 for (i = 0; i < iong->nat; i++)
660 /* Anion or cation? */
661 type = iong->qc[i] < 0 ? eIonNEG : eIonPOS;
663 /* Is this ion in the compartment that we look at? */
664 if (compartment_contains_atom(left, right, 0, iong->xc[i][sd], box[sd][sd], &dist) )
666 /* Now put it into the list containing only ions of its type */
667 add_to_list(i, &s->comp[comp][type], dist);
669 /* Master also checks through which channel each ion has passed */
670 if (MASTER(cr) && (iong->comp_now != NULL))
672 ion_nr_global = iong->ind[i] + 1; /* PDB index starts at 1 ... */
673 detect_flux_per_channel(ion_nr_global, comp, type, iong->xc[i],
674 &iong->comp_now[i], &iong->comp_from[i], &iong->channel_label[i],
675 sc, cyl0_r2, cyl1_r2, step, bRerun, fpout);
680 not_in_comp[comp] += 1;
683 /* Correct the time-averaged number of ions in both compartments */
684 update_time_window(&s->comp[comp][eIonNEG], sc->nAverage, replace);
685 update_time_window(&s->comp[comp][eIonPOS], sc->nAverage, replace);
688 /* Flux detection warnings */
694 "%s Warning: %d atoms were detected as being in both channels! Probably your split\n"
695 "%s cylinder is way too large, or one compartment has collapsed (step %"GMX_PRId64 ")\n",
696 SwS, s->cyl0and1, SwS, step);
698 fprintf(s->fpout, "Warning: %d atoms were assigned to both channels!\n", s->cyl0and1);
705 /* Consistency checks */
706 if (not_in_comp[eCompA] + not_in_comp[eCompB] != iong->nat)
710 fprintf(fpout, "# Warning: Inconsistency during ion compartmentalization. !inA: %d, !inB: %d, total ions %d\n",
711 not_in_comp[eCompA], not_in_comp[eCompB], iong->nat);
716 fprintf(stderr, "%s rank %d: Inconsistency during ion compartmentalization. !inA: %d, !inB: %d, total ions %d\n",
717 SwS, cr->nodeid, not_in_comp[eCompA], not_in_comp[eCompB], iong->nat);
721 sum = s->comp[eCompA][eIonNEG].nat + s->comp[eCompA][eIonPOS].nat
722 + s->comp[eCompB][eIonNEG].nat + s->comp[eCompB][eIonPOS].nat;
723 if (sum != iong->nat)
727 fprintf(fpout, "# Warning: %d atoms are in the ion group, but altogether %d have been assigned to the compartments.\n",
733 fprintf(stderr, "%s rank %d: %d atoms are in the ion group, but altogether %d have been assigned to the compartments.\n",
734 SwS, cr->nodeid, iong->nat, sum);
742 /*! \brief Set up the compartments and get lists of solvent atoms in each compartment */
743 static void compartmentalize_solvent(
755 int sum, not_in_comp[eCompNR]; /* consistency check */
759 solg = &s->group[eGrpSolvent];
763 for (comp = eCompA; comp <= eCompB; comp++)
765 /* Get lists of atoms that match criteria for this compartment */
766 get_compartment_boundaries(comp, sc->si_priv, box, &left, &right);
768 /* First clear the solvent molecule lists */
769 s->compsol[comp].nat = 0;
770 not_in_comp[comp] = 0; /* consistency check */
772 /* Loop over the solvent MOLECULES */
773 for (i = 0; i < sc->nat_sol; i += apm)
775 if (compartment_contains_atom(left, right, 0, solg->xc[i][sd], box[sd][sd], &dist))
777 /* Add the whole molecule to the list */
778 for (j = 0; j < apm; j++)
780 add_to_list(i+j, &s->compsol[comp], dist);
785 not_in_comp[comp] += apm;
792 fprintf(fpout, "# Solv. molecules in comp.%s: %d comp.%s: %d\n",
793 CompStr[eCompA], s->compsol[eCompA].nat/apm,
794 CompStr[eCompB], s->compsol[eCompB].nat/apm);
797 /* Consistency checks */
798 if (not_in_comp[eCompA] + not_in_comp[eCompB] != solg->nat)
802 fprintf(fpout, "# Warning: Inconsistency during solvent compartmentalization. !inA: %d, !inB: %d, solvent atoms %d\n",
803 not_in_comp[eCompA], not_in_comp[eCompB], solg->nat);
808 fprintf(stderr, "%s rank %d: Inconsistency during solvent compartmentalization. !inA: %d, !inB: %d, solvent atoms %d\n",
809 SwS, cr->nodeid, not_in_comp[eCompA], not_in_comp[eCompB], solg->nat);
812 sum = s->compsol[eCompA].nat + s->compsol[eCompB].nat;
813 if (sum != solg->nat)
817 fprintf(fpout, "# Warning: %d atoms in solvent group, but %d have been assigned to the compartments.\n",
823 fprintf(stderr, "%s rank %d: %d atoms in solvent group, but %d have been assigned to the compartments.\n",
824 SwS, cr->nodeid, solg->nat, sum);
830 /*! \brief Find out how many group atoms are in the compartments initially */
831 static void get_initial_ioncounts(
833 rvec x[], /* the initial positions */
847 /* Copy the initial swap group positions to the collective array so
848 * that we can compartmentalize */
849 for (i = 0; i < sc->nat; i++)
852 copy_rvec(x[ind], s->group[eGrpIons].xc[i]);
855 /* Set up the compartments and get lists of atoms in each compartment */
856 compartmentalize_ions(cr, sc, box, 0, s->fpout, bRerun);
858 /* Set initial concentrations if requested */
859 for (ic = 0; ic < eCompNR; ic++)
861 s->comp[ic][eIonPOS].nat_req = sc->ncations[ic];
862 s->comp[ic][eIonNEG].nat_req = sc->nanions[ic];
864 for (ic = 0; ic < eCompNR; ic++)
866 for (ii = 0; ii < eIonNR; ii++)
868 if (s->comp[ic][ii].nat_req < 0)
870 s->comp[ic][ii].nat_req = s->comp[ic][ii].nat;
875 /* Check whether the number of requested ions adds up to the total number of ions */
876 for (ii = 0; ii < eIonNR; ii++)
878 req[ii] = s->comp[eCompA][ii].nat_req + s->comp[eCompB][ii].nat_req;
879 tot[ii] = s->comp[eCompA][ii].nat + s->comp[eCompB][ii].nat;
881 if ( (req[eCompA] != tot[eCompA]) || (req[eCompB] != tot[eCompB ]) )
883 gmx_fatal(FARGS, "Mismatch of the number of ions summed over both compartments.\n"
884 "You requested a total of %d anions and %d cations,\n"
885 "but there are a total of %d anions and %d cations in the system.\n",
886 req[eIonNEG], req[eIonPOS],
887 tot[eIonNEG], tot[eIonPOS]);
890 /* Initialize time-averaging:
891 * Write initial concentrations to all time bins to start with */
892 for (ic = 0; ic < eCompNR; ic++)
894 for (ii = 0; ii < eIonNR; ii++)
896 s->comp[ic][ii].nat_av = s->comp[ic][ii].nat;
897 for (i = 0; i < sc->nAverage; i++)
899 s->comp[ic][ii].nat_past[i] = s->comp[ic][ii].nat;
906 /*! \brief Copy history of ion counts from checkpoint file.
908 * When called, the checkpoint file has already been read in. Here we copy
909 * over the values from .cpt file to the swap data structure.
911 static void get_initial_ioncounts_from_cpt(
912 t_inputrec *ir, swapstate_t *swapstate,
913 t_commrec *cr, gmx_bool bVerbose)
924 /* Copy the past values from the checkpoint values that have been read in already */
927 fprintf(stderr, "%s Copying values from checkpoint\n", SwS);
930 for (ic = 0; ic < eCompNR; ic++)
932 for (ii = 0; ii < eIonNR; ii++)
934 s->comp[ic][ii].nat_req = swapstate->nat_req[ic][ii];
935 s->comp[ic][ii].inflow_netto = swapstate->inflow_netto[ic][ii];
939 fprintf(stderr, "%s ... Influx netto: %d Requested: %d Past values: ", SwS,
940 s->comp[ic][ii].inflow_netto, s->comp[ic][ii].nat_req);
943 for (j = 0; j < sc->nAverage; j++)
945 s->comp[ic][ii].nat_past[j] = swapstate->nat_past[ic][ii][j];
948 fprintf(stderr, "%d ", s->comp[ic][ii].nat_past[j]);
953 fprintf(stderr, "\n");
961 /*! \brief The master lets all others know about the initial ion counts. */
962 static void bc_initial_concentrations(
971 for (ic = 0; ic < eCompNR; ic++)
973 for (ii = 0; ii < eIonNR; ii++)
975 gmx_bcast(sizeof(s->comp[ic][ii].nat_req), &(s->comp[ic][ii].nat_req), cr);
976 gmx_bcast(sizeof(s->comp[ic][ii].nat ), &(s->comp[ic][ii].nat ), cr);
977 gmx_bcast( swap->nAverage * sizeof(s->comp[ic][ii].nat_past[0]), s->comp[ic][ii].nat_past, cr);
983 /*! \brief Ensure that each atom belongs to at most one of the swap groups. */
984 static void check_swap_groups(t_swap *s, int nat, gmx_bool bVerbose)
988 atom_id *nGroup = NULL; /* This array counts for each atom in the MD system to
989 how many swap groups it belongs (should be 0 or 1!) */
991 int nMultiple = 0; /* Number of atoms belonging to multiple groups */
996 fprintf(stderr, "%s Making sure each atom belongs to at most one of the swap groups.\n", SwS);
999 /* Add one to the group count of atoms belonging to a swap group: */
1001 for (i = 0; i < eGrpNr; i++)
1004 for (j = 0; j < g->nat; j++)
1006 /* Get the global index of this atom of this group: */
1011 /* Make sure each atom belongs to at most one swap group: */
1012 for (j = 0; j < g->nat; j++)
1023 gmx_fatal(FARGS, "%s Cannot perform swapping since %d atom%s allocated to more than one swap index group.\n"
1024 "%s Each atom must be allocated to at most one of the split groups, the swap group, or the solvent.\n"
1025 "%s Check the .mdp file settings regarding the swap index groups or the index groups themselves.\n",
1026 SwS, nMultiple, (1 == nMultiple) ? " is" : "s are", SwSEmpty, SwSEmpty);
1031 /*! \brief Get the number of atoms per molecule for this group.
1033 * Also ensure that all the molecules in this group have this number of atoms.
1035 static int get_group_apm_check(
1039 const gmx_mtop_atomlookup_t alook,
1044 int molb, molnr, atnr_mol;
1047 ind = s->group[group].ind;
1048 nat = s->group[group].nat;
1050 /* Determine the number of solvent atoms per solvent molecule from the
1051 * first solvent atom: */
1053 gmx_mtop_atomnr_to_molblock_ind(alook, ind[i], &molb, &molnr, &atnr_mol);
1054 apm = mtop->molblock[molb].natoms_mol;
1058 fprintf(stderr, "%s Checking whether all %s molecules consist of %d atom%s\n",
1059 SwS, GrpString[group], apm, apm > 1 ? "s" : "");
1062 /* Check whether this is also true for all other solvent atoms */
1063 for (i = 1; i < nat; i++)
1065 gmx_mtop_atomnr_to_molblock_ind(alook, ind[i], &molb, &molnr, &atnr_mol);
1066 if (apm != mtop->molblock[molb].natoms_mol)
1068 gmx_fatal(FARGS, "Not all %s group molecules consist of %d atoms.",
1069 GrpString[group], apm);
1077 /*! \brief Print the legend to the swap output file.
1079 * Also print the initial ion counts
1081 static void print_ionlist_legend(t_inputrec *ir, const output_env_t oenv)
1083 const char **legend;
1089 s = ir->swap->si_priv;
1091 snew(legend, eCompNR*eIonNR*3 + 2 + eChanNR*eIonNR + 1);
1092 for (ic = count = 0; ic < eCompNR; ic++)
1094 for (ii = 0; ii < eIonNR; ii++)
1096 sprintf(buf, "%s %ss", CompStr[ic], IonString[ii]);
1097 legend[count++] = gmx_strdup(buf);
1098 sprintf(buf, "%s av. mismatch to %d%s",
1099 CompStr[ic], s->comp[ic][ii].nat_req, IonStr[ii]);
1100 legend[count++] = gmx_strdup(buf);
1101 sprintf(buf, "%s netto %s influx", CompStr[ic], IonString[ii]);
1102 legend[count++] = gmx_strdup(buf);
1105 sprintf(buf, "%scenter of %s of split group 0", SwapStr[ir->eSwapCoords], (NULL != s->group[eGrpSplit0].m) ? "mass" : "geometry");
1106 legend[count++] = gmx_strdup(buf);
1107 sprintf(buf, "%scenter of %s of split group 1", SwapStr[ir->eSwapCoords], (NULL != s->group[eGrpSplit1].m) ? "mass" : "geometry");
1108 legend[count++] = gmx_strdup(buf);
1110 for (ic = 0; ic < eChanNR; ic++)
1112 for (ii = 0; ii < eIonNR; ii++)
1114 sprintf(buf, "A->ch%d->B %s permeations", ic, IonString[ii]);
1115 legend[count++] = gmx_strdup(buf);
1119 sprintf(buf, "leakage");
1120 legend[count++] = gmx_strdup(buf);
1122 xvgr_legend(s->fpout, count, legend, oenv);
1124 fprintf(s->fpout, "# Instantaneous ion counts and time-averaged differences to requested numbers\n");
1125 fprintf(s->fpout, "# time[ps] A_an diff t_in A_cat diff t_in B_an diff t_in B_cat diff t_in ");
1126 fprintf(s->fpout, " %s-Split0 %s-Split1", DimStr[s->swapdim], DimStr[s->swapdim]);
1127 fprintf(s->fpout, " A-ch0-B_an A-ch0-B_cat A-ch1-B_an A-ch1-B_cat ion_leakage\n");
1132 /*! \brief Initialize channel ion flux detection routine.
1134 * Initialize arrays that keep track of where the ions come from and where
1137 static void detect_flux_per_channel_init(
1140 swapstate_t *swapstate,
1141 gmx_bool bStartFromCpt)
1147 g = &(s->group[eGrpIons]);
1149 /* All these flux detection routines run on the master only */
1153 g->comp_from = NULL;
1154 g->channel_label = NULL;
1159 /******************************************************/
1160 /* Channel and domain history for the individual ions */
1161 /******************************************************/
1162 if (bStartFromCpt) /* set the pointers right */
1164 g->comp_from = swapstate->comp_from;
1165 g->channel_label = swapstate->channel_label;
1167 else /* allocate memory */
1169 snew(g->comp_from, g->nat);
1170 swapstate->comp_from = g->comp_from;
1171 snew(g->channel_label, g->nat);
1172 swapstate->channel_label = g->channel_label;
1174 snew(g->comp_now, g->nat);
1176 /* Initialize the channel and domain history counters */
1177 for (i = 0; i < g->nat; i++)
1179 g->comp_now[i] = eDomainNotset;
1182 g->comp_from[i] = eDomainNotset;
1183 g->channel_label[i] = eChHistPassedNone;
1187 /************************************/
1188 /* Channel fluxes for both channels */
1189 /************************************/
1196 fprintf(stderr, "%s Copying channel fluxes from checkpoint file data\n", SwS);
1199 for (ic = 0; ic < eChanNR; ic++)
1201 fprintf(stderr, "%s Channel %d flux history: ", SwS, ic);
1202 for (ii = 0; ii < eIonNR; ii++)
1206 s->fluxfromAtoB[ic][ii] = swapstate->fluxfromAtoB[ic][ii];
1210 s->fluxfromAtoB[ic][ii] = 0;
1213 fprintf(stderr, "%d %s%s ", s->fluxfromAtoB[ic][ii], IonString[ii], s->fluxfromAtoB[ic][ii] == 1 ? "" : "s");
1215 fprintf(stderr, "\n");
1219 s->fluxleak = swapstate->fluxleak;
1223 snew(s->fluxleak, 1);
1225 /* Set pointer for checkpoint writing */
1226 swapstate->fluxleak = s->fluxleak;
1229 /* Set pointers for checkpoint writing */
1230 for (ic = 0; ic < eChanNR; ic++)
1232 for (ii = 0; ii < eIonNR; ii++)
1234 swapstate->fluxfromAtoB_p[ic][ii] = &(s->fluxfromAtoB[ic][ii]);
1240 /*! \brief Outputs the initial structure to PDB file for debugging reasons.
1242 * Output the starting structure so that in case of multimeric channels
1243 * the user can check whether we have the correct PBC image for all atoms.
1244 * If this is not correct, the ion counts per channel will be very likely
1247 static void outputStartStructureIfWanted(gmx_mtop_t *mtop, rvec *x, int ePBC, matrix box)
1249 char *env = getenv("GMX_COMPELDUMP");
1253 fprintf(stderr, "\n%s Found env.var. GMX_COMPELDUMP, will output CompEL starting structure made whole.\n"
1254 "%s In case of multimeric channels, please check whether they have the correct PBC representation.\n",
1257 write_sto_conf_mtop("CompELAssumedWholeConfiguration.pdb", *mtop->name, mtop, x, NULL, ePBC, box);
1262 /*! \brief Initialize the swapstate structure, used for checkpoint writing.
1264 * The swapstate struct stores the information we need to make the channels
1265 * whole again after restarts from a checkpoint file. Here we do the following:\n
1266 * a) If we did not start from .cpt, we prepare the struct for proper .cpt writing,\n
1267 * b) if we did start from .cpt, we copy over the last whole structures from .cpt,\n
1268 * c) in any case, for subsequent checkpoint writing, we set the pointers in\n
1269 * swapstate to the x_old arrays, which contain the correct PBC representation of
1270 * multimeric channels at the last time step.
1272 static void init_swapstate(
1273 swapstate_t *swapstate,
1276 rvec x[], /* the initial positions */
1281 rvec *x_pbc = NULL; /* positions of the whole MD system with molecules made whole */
1288 /* We always need the last whole positions such that
1289 * in the next time step we can make the channels whole again in PBC */
1290 if (swapstate->bFromCpt)
1292 /* Copy the last whole positions of each channel from .cpt */
1293 g = &(s->group[eGrpSplit0]);
1294 for (i = 0; i < g->nat; i++)
1296 copy_rvec(swapstate->xc_old_whole[eChan0][i], g->xc_old[i]);
1298 g = &(s->group[eGrpSplit1]);
1299 for (i = 0; i < g->nat; i++)
1301 copy_rvec(swapstate->xc_old_whole[eChan1][i], g->xc_old[i]);
1306 /* Extract the initial split group positions. */
1308 /* Remove pbc, make molecule whole. */
1309 snew(x_pbc, mtop->natoms);
1310 m_rveccopy(mtop->natoms, x, x_pbc);
1312 /* This can only make individual molecules whole, not multimers */
1313 do_pbc_mtop(NULL, ePBC, box, mtop, x_pbc);
1315 /* Output the starting structure? */
1316 outputStartStructureIfWanted(mtop, x_pbc, ePBC, box);
1318 /* If this is the first run (i.e. no checkpoint present) we assume
1319 * that the starting positions give us the correct PBC representation */
1320 for (ig = eGrpSplit0; ig <= eGrpSplit1; ig++)
1322 g = &(s->group[ig]);
1323 for (i = 0; i < g->nat; i++)
1325 copy_rvec(x_pbc[g->ind[i]], g->xc_old[i]);
1330 /* Prepare swapstate arrays for later checkpoint writing */
1331 swapstate->nat[eChan0] = s->group[eGrpSplit0].nat;
1332 swapstate->nat[eChan1] = s->group[eGrpSplit1].nat;
1335 /* For subsequent checkpoint writing, set the swapstate pointers to the xc_old
1336 * arrays that get updated at every swapping step */
1337 swapstate->xc_old_whole_p[eChan0] = &s->group[eGrpSplit0].xc_old;
1338 swapstate->xc_old_whole_p[eChan1] = &s->group[eGrpSplit1].xc_old;
1342 extern void init_swapcoords(
1350 swapstate_t *swapstate,
1352 const output_env_t oenv,
1353 unsigned long Flags)
1355 int i, ic, ig, ii, j;
1360 gmx_bool bAppend, bStartFromCpt, bRerun;
1361 gmx_mtop_atomlookup_t alook = NULL;
1364 alook = gmx_mtop_atomlookup_init(mtop);
1366 if ( (PAR(cr)) && !DOMAINDECOMP(cr) )
1368 gmx_fatal(FARGS, "Position swapping is only implemented for domain decomposition!");
1371 bAppend = Flags & MD_APPENDFILES;
1372 bStartFromCpt = Flags & MD_STARTFROMCPT;
1373 bRerun = Flags & MD_RERUN;
1376 snew(sc->si_priv, 1);
1383 gmx_fatal(FARGS, "%s This module does not support reruns in parallel\nPlease request a serial run with -nt 1 / -np 1\n", SwS);
1386 fprintf(stderr, "%s Rerun - using every available frame\n", SwS);
1388 sc->nAverage = 1; /* averaging makes no sense for reruns */
1391 if (MASTER(cr) && !bAppend)
1393 fprintf(fplog, "\nInitializing ion/water position exchanges\n");
1394 please_cite(fplog, "Kutzner2011b");
1397 switch (ir->eSwapCoords)
1413 /* Copy some data to the group structures for convenience */
1414 /* Number of atoms in the group */
1415 s->group[eGrpIons ].nat = sc->nat;
1416 s->group[eGrpSplit0 ].nat = sc->nat_split[0];
1417 s->group[eGrpSplit1 ].nat = sc->nat_split[1];
1418 s->group[eGrpSolvent].nat = sc->nat_sol;
1419 /* Pointer to the indices */
1420 s->group[eGrpIons ].ind = sc->ind;
1421 s->group[eGrpSplit0 ].ind = sc->ind_split[0];
1422 s->group[eGrpSplit1 ].ind = sc->ind_split[1];
1423 s->group[eGrpSolvent].ind = sc->ind_sol;
1425 check_swap_groups(s, mtop->natoms, bVerbose && MASTER(cr));
1427 /* Allocate space for the collective arrays for all groups */
1428 for (ig = 0; ig < eGrpNr; ig++)
1430 g = &(s->group[ig]);
1431 snew(g->xc, g->nat);
1432 snew(g->c_ind_loc, g->nat);
1433 /* For the split groups (the channels) we need some extra memory to
1434 * be able to make the molecules whole even if they span more than
1435 * half of the box size. */
1436 if (eGrpSplit0 == ig || eGrpSplit1 == ig)
1438 snew(g->xc_shifts, g->nat);
1439 snew(g->xc_eshifts, g->nat);
1440 snew(g->xc_old, g->nat);
1446 init_swapstate(swapstate, sc, mtop, x, box, ir->ePBC);
1449 /* After init_swapstate we have a set of (old) whole positions for our
1450 * channels. Now transfer that to all nodes */
1453 for (ig = eGrpSplit0; ig <= eGrpSplit1; ig++)
1455 g = &(s->group[ig]);
1456 gmx_bcast((g->nat)*sizeof((g->xc_old)[0]), g->xc_old, (cr));
1460 /* Make sure that all molecules in the ion and solvent groups contain the
1461 * same number of atoms each */
1462 s->group[eGrpIons ].apm = get_group_apm_check(eGrpIons, s, MASTER(cr) && bVerbose, alook, mtop);
1463 s->group[eGrpSolvent].apm = get_group_apm_check(eGrpSolvent, s, MASTER(cr) && bVerbose, alook, mtop);
1465 /* Save masses where needed */
1466 s->group[eGrpIons ].m = NULL;
1467 /* We only need enough space to determine a single solvent molecule's
1468 * center at at time */
1469 g = &(s->group[eGrpSolvent]);
1472 /* Need mass-weighted center of split group? */
1473 for (j = 0, ig = eGrpSplit0; j < eChanNR; ig++, j++)
1475 g = &(s->group[ig]);
1476 if (TRUE == sc->massw_split[j])
1478 /* Save the split group charges if mass-weighting is requested */
1480 for (i = 0; i < g->nat; i++)
1482 gmx_mtop_atomnr_to_atom(alook, g->ind[i], &atom);
1492 /* Save the ionic charges */
1493 g = &(s->group[eGrpIons]);
1494 snew(g->qc, g->nat);
1495 for (i = 0; i < g->nat; i++)
1497 gmx_mtop_atomnr_to_atom(alook, g->ind[i], &atom);
1502 set_pbc(s->pbc, -1, box);
1509 fprintf(stderr, "%s Opening output file %s%s\n", SwS, fn, bAppend ? " for appending" : "");
1512 s->fpout = gmx_fio_fopen(fn, bAppend ? "a" : "w" );
1516 xvgr_header(s->fpout, "Ion counts", "Time (ps)", "counts", exvggtXNY, oenv);
1518 for (ig = 0; ig < eGrpNr; ig++)
1520 g = &(s->group[ig]);
1521 fprintf(s->fpout, "# %s group contains %d atom%s", GrpString[ig], g->nat, (g->nat > 1) ? "s" : "");
1522 if (eGrpSolvent == ig || eGrpIons == ig)
1524 fprintf(s->fpout, " with %d atom%s in each molecule", g->apm, (g->apm > 1) ? "s" : "");
1526 fprintf(s->fpout, ".\n");
1529 fprintf(s->fpout, "#\n# Initial positions of split groups:\n");
1532 for (j = 0, ig = eGrpSplit0; j < eChanNR; j++, ig++)
1534 g = &(s->group[ig]);
1535 for (i = 0; i < g->nat; i++)
1537 copy_rvec(x[sc->ind_split[j][i]], g->xc[i]);
1539 if (eGrpSplit0 == ig || eGrpSplit1 == ig)
1541 /* xc has the correct PBC representation for the two channels, so we do
1542 * not need to correct for that */
1543 get_center(g->xc, g->m, g->nat, g->center);
1547 /* For the water molecules, we need to make the molecules whole */
1548 get_molecule_center(g->xc, g->nat, g->m, g->center, s->pbc);
1552 fprintf(s->fpout, "# %s group %s-center %5f nm\n", GrpString[ig],
1553 DimStr[s->swapdim], g->center[s->swapdim]);
1559 fprintf(s->fpout, "#\n");
1560 fprintf(s->fpout, "# split0 cylinder radius %f nm, up %f nm, down %f nm\n",
1561 sc->cyl0r, sc->cyl0u, sc->cyl0l);
1562 fprintf(s->fpout, "# split1 cylinder radius %f nm, up %f nm, down %f nm\n",
1563 sc->cyl1r, sc->cyl1u, sc->cyl1l);
1568 fprintf(s->fpout, "#\n");
1571 fprintf(s->fpout, "# Coupling constant (number of swap attempt steps to average over): %d (translates to %f ps).\n",
1572 sc->nAverage, sc->nAverage*sc->nstswap*ir->delta_t);
1573 fprintf(s->fpout, "# Threshold is %f\n", sc->threshold);
1574 fprintf(s->fpout, "#\n");
1575 fprintf(s->fpout, "# Remarks about which atoms passed which channel use global atoms numbers starting at one.\n");
1584 /* Prepare for parallel or serial run */
1587 for (ig = 0; ig < eGrpNr; ig++)
1589 g = &(s->group[ig]);
1597 for (ig = 0; ig < eGrpNr; ig++)
1599 g = &(s->group[ig]);
1600 g->nat_loc = g->nat;
1601 g->ind_loc = g->ind;
1602 /* c_ind_loc needs to be set to identity in the serial case */
1603 for (i = 0; i < g->nat; i++)
1605 g->c_ind_loc[i] = i;
1610 /* Allocate memory for the ion counts time window */
1611 for (ic = 0; ic < eCompNR; ic++)
1613 for (ii = 0; ii < eIonNR; ii++)
1615 snew(s->comp[ic][ii].nat_past, sc->nAverage);
1619 /* Get the initial ion concentrations and let the other nodes know */
1622 swapstate->nions = s->group[eGrpIons].nat;
1626 get_initial_ioncounts_from_cpt(ir, swapstate, cr, bVerbose);
1630 fprintf(stderr, "%s Determining initial ion counts.\n", SwS);
1631 get_initial_ioncounts(ir, x, box, cr, bRerun);
1634 /* Prepare (further) checkpoint writes ... */
1637 /* Consistency check */
1638 if (swapstate->nAverage != sc->nAverage)
1640 gmx_fatal(FARGS, "%s Ion count averaging steps mismatch! checkpoint: %d, tpr: %d",
1641 SwS, swapstate->nAverage, sc->nAverage);
1646 swapstate->nAverage = sc->nAverage;
1648 fprintf(stderr, "%s Setting pointers for checkpoint writing\n", SwS);
1649 for (ic = 0; ic < eCompNR; ic++)
1651 for (ii = 0; ii < eIonNR; ii++)
1653 swapstate->nat_req_p[ic][ii] = &(s->comp[ic][ii].nat_req);
1654 swapstate->nat_past_p[ic][ii] = &(s->comp[ic][ii].nat_past[0]);
1655 swapstate->inflow_netto_p[ic][ii] = &(s->comp[ic][ii].inflow_netto);
1659 /* Determine the total charge imbalance */
1660 s->deltaQ = ( (-1) * s->comp[eCompA][eIonNEG].nat_req + s->comp[eCompA][eIonPOS].nat_req )
1661 - ( (-1) * s->comp[eCompB][eIonNEG].nat_req + s->comp[eCompB][eIonPOS].nat_req );
1665 fprintf(stderr, "%s Requested charge imbalance is Q(A) - Q(B) = %gz.\n", SwS, s->deltaQ);
1669 fprintf(s->fpout, "# Requested charge imbalance is Q(A)-Q(B) = %gz.\n", s->deltaQ);
1675 bc_initial_concentrations(cr, ir->swap);
1678 /* Put the time-averaged number of ions for all compartments */
1679 for (ic = 0; ic < eCompNR; ic++)
1681 for (ii = 0; ii < eIonNR; ii++)
1683 update_time_window(&(s->comp[ic][ii]), sc->nAverage, -1);
1687 /* Initialize arrays that keep track of through which channel the ions go */
1688 detect_flux_per_channel_init(cr, s, swapstate, bStartFromCpt);
1690 /* We need to print the legend if we open this file for the first time. */
1691 if (MASTER(cr) && !bAppend)
1693 print_ionlist_legend(ir, oenv);
1698 extern void dd_make_local_swap_groups(gmx_domdec_t *dd, t_swapcoords *sc)
1704 /* Make ion group, split groups and solvent group */
1705 for (ig = 0; ig < eGrpNr; ig++)
1707 g = &(sc->si_priv->group[ig]);
1708 dd_make_local_group_indices(dd->ga2la, g->nat, g->ind,
1709 &(g->nat_loc), &(g->ind_loc), &(g->nalloc_loc), g->c_ind_loc);
1714 /*! \brief Do we need to swap ions with water molecules at this step?
1716 * From the requested and average ion counts we determine whether a swap is needed
1717 * at this time step.
1719 static gmx_bool need_swap(t_swapcoords *sc)
1726 for (ic = 0; ic < eCompNR; ic++)
1728 for (ii = 0; ii < eIonNR; ii++)
1730 if (s->comp[ic][ii].nat_req - s->comp[ic][ii].nat_av >= sc->threshold)
1740 /*! \brief Return index of atom that we can use for swapping.
1742 * Returns the index of an atom that is far off the compartment boundaries.
1743 * Other atoms of the molecule (if any) will directly follow the returned index
1745 static int get_index_of_distant_atom(
1746 t_compartment *comp,
1747 int apm) /* Atoms per molecule - just return the first atom index of a molecule */
1750 real d = GMX_REAL_MAX;
1753 /* comp->nat contains the original number of atoms in this compartment
1754 * prior to doing any swaps. Some of these atoms may already have been
1755 * swapped out, but then they are marked with a distance of GMX_REAL_MAX
1757 for (i = 0; i < comp->nat_old; i += apm)
1759 if (comp->dist[i] < d)
1762 d = comp->dist[ibest];
1768 gmx_fatal(FARGS, "Could not get index of swap atom. Compartment atoms %d before swaps, atoms per molecule %d.",
1769 comp->nat_old, apm);
1772 /* Set the distance of this index to infinity such that it won't get selected again in
1775 comp->dist[ibest] = GMX_REAL_MAX;
1777 return comp->ind[ibest];
1781 /*! \brief Swaps centers of mass and makes molecule whole if broken */
1782 static void translate_positions(
1790 rvec reference, dx, correctPBCimage;
1793 /* Use the first atom as the reference for PBC */
1794 copy_rvec(x[0], reference);
1796 for (i = 0; i < apm; i++)
1798 /* PBC distance between position and reference */
1799 pbc_dx(pbc, x[i], reference, dx);
1801 /* Add PBC distance to reference */
1802 rvec_add(reference, dx, correctPBCimage);
1804 /* Subtract old_com from correct image and add new_com */
1805 rvec_dec(correctPBCimage, old_com);
1806 rvec_inc(correctPBCimage, new_com);
1808 copy_rvec(correctPBCimage, x[i]);
1813 /*! \brief Write back the the modified local positions from the collective array to the official positions. */
1814 static void apply_modified_positions(
1821 for (l = 0; l < g->nat_loc; l++)
1823 /* Get the right local index to write to */
1825 /* Where is the local atom in the collective array? */
1826 cind = g->c_ind_loc[l];
1828 /* Copy the possibly modified position */
1829 copy_rvec(g->xc[cind], x[ii]);
1834 extern gmx_bool do_swapcoords(
1839 gmx_wallcycle_t wcycle,
1848 int j, ii, ic, ig, im, gmax, nswaps;
1849 gmx_bool bSwap = FALSE;
1851 real vacancy[eCompNR][eIonNR];
1853 rvec solvent_center, ion_center;
1855 gmx_mtop_atomlookup_t alook = NULL;
1858 wallcycle_start(wcycle, ewcSWAP);
1863 /* Assemble all the positions of the swap group (ig = 0), the split groups
1864 * (ig = 1,2), and possibly the solvent group (ig = 3) */
1867 for (ig = 0; ig < gmax; ig++)
1869 g = &(s->group[ig]);
1870 if (eGrpSplit0 == ig || eGrpSplit1 == ig)
1872 /* The split groups, i.e. the channels. Here we need the full
1873 * communicate_group_positions(), so that we can make the molecules
1874 * whole even in cases where they span more than half of the box in
1876 communicate_group_positions(cr, g->xc, g->xc_shifts, g->xc_eshifts, TRUE,
1877 x, g->nat, g->nat_loc, g->ind_loc, g->c_ind_loc, g->xc_old, box);
1879 get_center(g->xc, g->m, g->nat, g->center); /* center of split groups == channels */
1883 /* Swap group (ions), and solvent group. These molecules are small
1884 * and we can always make them whole with a simple distance check.
1885 * Therefore we pass NULL as third argument. */
1886 communicate_group_positions(cr, g->xc, NULL, NULL, FALSE,
1887 x, g->nat, g->nat_loc, g->ind_loc, g->c_ind_loc, NULL, NULL);
1891 /* Set up the compartments and get lists of atoms in each compartment,
1892 * determine how many ions each compartment contains */
1893 compartmentalize_ions(cr, sc, box, step, s->fpout, bRerun);
1895 /* Output how many ions are in the compartments */
1898 print_ionlist(s, t, "");
1901 /* If we are doing a rerun, we are finished here, since we cannot perform
1908 /* Do we have to perform a swap? */
1909 bSwap = need_swap(sc);
1912 g = &(s->group[eGrpSolvent]);
1913 communicate_group_positions(cr, g->xc, NULL, NULL, FALSE,
1914 x, g->nat, g->nat_loc, g->ind_loc, g->c_ind_loc, NULL, NULL);
1916 compartmentalize_solvent(cr, sc, box, s->fpout);
1918 /* Determine where ions are missing and where ions are too many */
1919 for (ic = 0; ic < eCompNR; ic++)
1921 for (ii = 0; ii < eIonNR; ii++)
1923 vacancy[ic][ii] = s->comp[ic][ii].nat_req - s->comp[ic][ii].nat_av;
1927 /* Remember the original number of ions per compartment */
1928 for (ic = 0; ic < eCompNR; ic++)
1930 s->compsol[ic].nat_old = s->compsol[ic].nat;
1931 for (ii = 0; ii < eIonNR; ii++)
1933 s->comp[ic][ii].nat_old = s->comp[ic][ii].nat;
1937 /* Now actually correct the number of ions */
1938 g = &(s->group[eGrpSolvent]);
1940 alook = gmx_mtop_atomlookup_init(mtop);
1941 for (ic = 0; ic < eCompNR; ic++)
1943 for (ii = 0; ii < eIonNR; ii++)
1945 while (vacancy[ic][ii] >= sc->threshold)
1947 /* Swap in an ion */
1949 /* Get the xc-index of the first atom of a solvent molecule of this compartment */
1950 isol = get_index_of_distant_atom(&(s->compsol[ic]), s->group[eGrpSolvent].apm );
1952 /* Get the xc-index of an ion from the other compartment */
1953 iion = get_index_of_distant_atom(&(s->comp[(ic+1)%eCompNR][ii]), s->group[eGrpIons].apm );
1955 /* Get the solvent molecule's center of mass */
1956 for (im = 0; im < s->group[eGrpSolvent].apm; im++)
1958 gmx_mtop_atomnr_to_atom(alook, s->group[eGrpSolvent].ind[isol+im], &atom);
1959 s->group[eGrpSolvent].m[im] = atom->m;
1961 get_molecule_center(&(s->group[eGrpSolvent].xc[isol]), s->group[eGrpSolvent].apm, s->group[eGrpSolvent].m, solvent_center, s->pbc);
1962 get_molecule_center(&(s->group[eGrpIons ].xc[iion]), s->group[eGrpIons ].apm, NULL, ion_center, s->pbc);
1964 /* subtract com_solvent and add com_ion */
1965 translate_positions(&(s->group[eGrpSolvent].xc[isol]), s->group[eGrpSolvent].apm, solvent_center, ion_center, s->pbc);
1966 /* For the ion, subtract com_ion and add com_solvent */
1967 translate_positions(&(s->group[eGrpIons ].xc[iion]), s->group[eGrpIons ].apm, ion_center, solvent_center, s->pbc);
1970 vacancy[(ic+1) % eCompNR][ii]++;
1972 /* Keep track of the changes */
1973 s->comp[ic ][ii].nat++;
1974 s->comp[(ic+1) % eCompNR][ii].nat--;
1975 s->comp[ic ][ii].inflow_netto++;
1976 s->comp[(ic+1) % eCompNR][ii].inflow_netto--;
1977 /* Correct the past time window to still get the right averages from now on */
1978 s->comp[ic ][ii].nat_av++;
1979 s->comp[(ic+1) % eCompNR][ii].nat_av--;
1980 for (j = 0; j < sc->nAverage; j++)
1982 s->comp[ic ][ii].nat_past[j]++;
1983 s->comp[(ic+1) % eCompNR][ii].nat_past[j]--;
1985 /* Clear ion history */
1988 s->group[eGrpIons].channel_label[iion] = eChHistPassedNone;
1989 s->group[eGrpIons].comp_from[iion] = eDomainNotset;
1991 /* That was the swap */
1996 gmx_mtop_atomlookup_destroy(alook);
2000 fprintf(stderr, "%s Performed %d swap%s in step %"GMX_PRId64 ".\n", SwS, nswaps, nswaps > 1 ? "s" : "", step);
2002 if (s->fpout != NULL)
2004 print_ionlist(s, t, " # after swap");
2007 /* Write back the the modified local positions from the collective array to the official coordinates */
2008 apply_modified_positions(&(s->group[eGrpIons ]), x);
2009 apply_modified_positions(&(s->group[eGrpSolvent]), x);
2010 } /* end of if(bSwap) */
2012 wallcycle_stop(wcycle, ewcSWAP);