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/utility/cstringutil.h"
52 #include "gromacs/utility/smalloc.h"
53 #include "gromacs/mdlib/groupcoord.h"
54 #include "mtop_util.h"
62 #include "gromacs/fileio/confio.h"
63 #include "gromacs/timing/wallcycle.h"
64 #include "swapcoords.h"
66 static char *SwS = {"SWAP:"}; /**< For output that comes from the swap module */
67 static char *SwSEmpty = {" "}; /**< Placeholder for multi-line output */
68 static char* IonString[eIonNR] = {"anion", "cation" }; /**< Type of ion, used for verbose output */
69 static char* IonStr[eIonNR] = {"-", "+" }; /**< Type of ion, used for short output */
70 static char* CompStr[eCompNR] = {"A", "B" }; /**< Compartment name */
71 static char *SwapStr[eSwapTypesNR+1] = { "", "X-", "Y-", "Z-", NULL}; /**< Name for the swap types. */
72 static char *DimStr[DIM+1] = { "X", "Y", "Z", NULL}; /**< Name for the swap dimension. */
74 /* eGrpSplit0 and eGrpSplit1 _must_ be neighbors in this list because
75 * we sometimes loop from eGrpSplit0 to eGrpSplit1 */
77 eGrpIons, eGrpSplit0, eGrpSplit1, eGrpSolvent, eGrpNr
78 }; /**< Group identifier */
79 static char* GrpString[eGrpNr] = { "ion", "split0", "split1", "solvent" }; /**< Group name */
81 /** Keep track of through which channel the ions have passed */
82 enum eChannelHistory {
83 eChHistPassedNone, eChHistPassedCh0, eChHistPassedCh1, eChHistNr
85 static char* ChannelString[eChHistNr] = { "none", "channel0", "channel1" }; /**< Name for the channels */
87 /*! \brief Domain identifier.
89 * Keeps track of from which compartment the ions came before passing the
93 eDomainNotset, eDomainA, eDomainB, eDomainNr
95 static char* DomainString[eDomainNr] = { "not_assigned", "Domain_A", "Domain_B" }; /**< Name for the domains */
100 * Structure containing compartment-specific data.
102 typedef struct swap_compartment
104 int nat; /**< Number of atoms matching the
105 compartment conditions. */
106 int nat_old; /**< Number of atoms before swapping. */
107 int nat_req; /**< Requested number of atoms. */
108 real nat_av; /**< Time-averaged number of atoms matching
109 the compartment conditions. */
110 int *nat_past; /**< Past ion counts for time-averaging. */
111 int *ind; /**< Indices to coll array of atoms. */
112 real *dist; /**< Distance of atom to compartment center. */
113 int nalloc; /**< Allocation size for ind array. */
114 int inflow_netto; /**< Net inflow of ions into this compartment. */
119 * This structure contains data needed for each of the groups involved in swapping: ions, water,
122 typedef struct swap_group
124 int nat; /**< Number of atoms in the group */
125 int apm; /**< Number of atoms in each molecule */
126 atom_id *ind; /**< Global atom indices of the group */
127 atom_id *ind_loc; /**< Local atom indices of the group */
128 int nat_loc; /**< Number of local group atoms */
129 int nalloc_loc; /**< Allocation size for ind_loc */
130 rvec *xc; /**< Collective array of group atom positions */
131 ivec *xc_shifts; /**< Current (collective) shifts */
132 ivec *xc_eshifts; /**< Extra shifts since last DD step */
133 rvec *xc_old; /**< Old (collective) positions */
134 real *qc; /**< Collective array of charges */
135 int *c_ind_loc; /**< Position of local atoms in the
136 collective array, [0..nat_loc] */
137 real *m; /**< Masses (can be omitted) */
138 unsigned char *comp_from; /**< (Collective) Stores from which compartment this
139 atom has come. This way we keep track of through
140 which channel an ion permeates (only used for
142 unsigned char *comp_now; /**< In which compartment this ion is now */
143 unsigned char *channel_label; /**< Which channel was passed at last by this ion? */
144 rvec center; /**< Center of the group; COM if masses are used */
149 * Main (private) data structure for the position swapping protocol.
153 int swapdim; /**< One of XX, YY, ZZ */
154 t_pbc *pbc; /**< Needed to make molecules whole. */
155 FILE *fpout; /**< Output file. */
156 t_group group[eGrpNr]; /**< Ions, solvent or channels? */
157 t_compartment comp[eCompNR][eIonNR]; /**< Data for a specific compartment and ion type. */
158 t_compartment compsol[eCompNR]; /**< Solvent compartments. */
159 int fluxfromAtoB[eChanNR][eIonNR]; /**< Net flux per channels and ion type. */
160 int ncyl0ions; /**< Number of ions residing in channel 0. */
161 int ncyl1ions; /**< Same for channel 1. */
162 int cyl0and1; /**< Ions assigned to cyl0 and cyl1. Not good. */
163 int *fluxleak; /**< Pointer to a single int value holding the
164 flux not going through any of the channels. */
165 real deltaQ; /**< The charge imbalance between the compartments. */
170 /*! \brief Check whether point is in channel.
172 * A channel is a cylinder defined by a disc
173 * with radius r around its center c. The thickness of the cylinder is
180 * <---------c--------->
186 static gmx_bool is_in_channel(
187 rvec point, /* Point under consideration */
188 rvec center, /* 'Center' of cylinder */
189 real d_up, /* Upper extension */
190 real d_down, /* Lower extensions */
191 real r_cyl2, /* Cylinder radius squared */
193 int normal) /* The membrane normal direction is typically 3, i.e. ZZ, but can be X or Y also */
196 int plane1, plane2; /* Directions tangential to membrane */
199 plane1 = (normal + 1) % 3; /* typically 0, i.e. XX */
200 plane2 = (normal + 2) % 3; /* typically 1, i.e. YY */
202 /* Get the distance vector dr between the point and the center of the cylinder */
203 pbc_dx(pbc, point, center, dr); /* This puts center in the origin */
205 /* Check vertical direction */
206 if ( (dr[normal] > d_up) || (dr[normal] < -d_down) )
211 /* Check radial direction */
212 if ( (dr[plane1]*dr[plane1] + dr[plane2]*dr[plane2]) > r_cyl2)
217 /* All check passed, this point is in the cylinder */
222 /*! \brief Prints to swap output file which ions are in which compartment. */
223 static void print_ionlist(
228 int itype, icomp, i, j;
232 fprintf(s->fpout, "%12.5e", time);
233 for (icomp = 0; icomp < eCompNR; icomp++)
235 for (itype = 0; itype < eIonNR; itype++)
237 comp = &(s->comp[icomp][itype]);
238 fprintf(s->fpout, "%7d%7.1f%7d", comp->nat, comp->nat_av-comp->nat_req, comp->inflow_netto);
241 fprintf(s->fpout, "%12.3e%12.3e",
242 s->group[eGrpSplit0].center[s->swapdim],
243 s->group[eGrpSplit1].center[s->swapdim]);
245 for (i = 0; i < eChanNR; i++)
247 for (j = 0; j < eIonNR; j++)
249 fprintf(s->fpout, "%12d", s->fluxfromAtoB[i][j]);
253 /* Also print the number of ions that leaked from A to B: */
254 fprintf(s->fpout, "%12d", *s->fluxleak);
256 fprintf(s->fpout, "%s\n", comment);
260 /*! \brief Get the center of a group of nat atoms.
262 * Since with PBC an atom group might not be whole, use the first atom as the
263 * reference atom and determine the center with respect to this reference.
265 static void get_molecule_center(
273 rvec weightedPBCimage;
275 rvec reference, correctPBCimage, dx;
278 /* Use the first atom as the reference and put other atoms near that one */
279 /* This does not work for large molecules that span > half of the box! */
280 copy_rvec(x[0], reference);
282 /* Calculate either the weighted center or simply the center of geometry */
285 for (i = 0; i < nat; i++)
287 /* PBC distance between position and reference */
288 pbc_dx(pbc, x[i], reference, dx);
290 /* Add PBC distance to reference */
291 rvec_add(reference, dx, correctPBCimage);
293 /* Take weight into account */
303 svmul(wi, correctPBCimage, weightedPBCimage);
306 rvec_inc(center, weightedPBCimage);
310 svmul(1.0/wsum, center, center);
315 /*! \brief Return TRUE if ion is found in the compartment.
317 * Returns TRUE if x is between (w1+gap) and (w2-gap)
321 * ||-----------|--+--|----------o----------|--+--|---------------------||
322 * w1 ????????????????????? w2
326 static gmx_bool compartment_contains_atom(
327 real w1, /* position of wall atom 1 */
328 real w2, /* position of wall atom 2 */
331 real l, /* length of the box, from || to || in the sketch */
332 real *distance_from_center)
337 /* First set the origin in the middle of w1 and w2 */
343 /* Now choose the PBC image of x that is closest to the origin: */
354 *distance_from_center = (real)fabs(x);
356 /* Return TRUE if we now are in area "????" */
357 if ( (x >= (w1+gap)) && (x < (w2-gap)) )
368 /*! \brief Updates the time-averaged number of ions in a compartment. */
369 static void update_time_window(t_compartment *comp, int values, int replace)
375 /* Put in the new value */
378 comp->nat_past[replace] = comp->nat;
381 /* Compute the new time-average */
383 for (i = 0; i < values; i++)
385 average += comp->nat_past[i];
388 comp->nat_av = average;
392 /*! \brief Add atom with collective index ci to the list 'comp'. */
393 static void add_to_list(
394 int ci, /* index of this ion in the collective array xc, qc */
395 t_compartment *comp, /* Compartment to add this atom to */
396 real distance) /* Shortest distance of this atom to the compartment center */
403 if (nr >= comp->nalloc)
405 comp->nalloc = over_alloc_dd(nr+1);
406 srenew(comp->ind, comp->nalloc);
407 srenew(comp->dist, comp->nalloc);
410 comp->dist[nr] = distance;
415 /*! \brief Determine the compartment boundaries from the channel centers. */
416 static void get_compartment_boundaries(
420 real *left, real *right)
423 real leftpos, rightpos, leftpos_orig;
428 gmx_fatal(FARGS, "No compartment %d.", c);
431 pos0 = s->group[eGrpSplit0].center[s->swapdim];
432 pos1 = s->group[eGrpSplit1].center[s->swapdim];
445 /* This gets us the other compartment: */
448 leftpos_orig = leftpos;
450 rightpos = leftpos_orig + box[s->swapdim][s->swapdim];
458 /*! \brief Determine the per-channel ion flux.
460 * To determine the flux through the individual channels, we
461 * remember the compartment and channel history of each ion. An ion can be
462 * either in channel0 or channel1, or in the remaining volume of compartment
466 * +-----------------+
469 * ||||||||||0|||||||| bilayer with channel 0
474 * |||||1||||||||||||| bilayer with channel 1
477 * +-----------------+
481 static void detect_flux_per_channel(
486 unsigned char *comp_now,
487 unsigned char *comp_from,
488 unsigned char *channel_label,
498 gmx_bool in_cyl0, in_cyl1;
505 /* Check whether ion is inside any of the channels */
506 in_cyl0 = is_in_channel(ion_pos, s->group[eGrpSplit0].center, sc->cyl0u, sc->cyl0l, cyl0_r2, s->pbc, sd);
507 in_cyl1 = is_in_channel(ion_pos, s->group[eGrpSplit1].center, sc->cyl1u, sc->cyl1l, cyl1_r2, s->pbc, sd);
509 if (in_cyl0 && in_cyl1)
511 /* Ion appears to be in both channels. Something is severely wrong! */
513 *comp_now = eDomainNotset;
514 *comp_from = eDomainNotset;
515 *channel_label = eChHistPassedNone;
519 /* Ion is in channel 0 now */
520 *channel_label = eChHistPassedCh0;
521 *comp_now = eDomainNotset;
526 /* Ion is in channel 1 now */
527 *channel_label = eChHistPassedCh1;
528 *comp_now = eDomainNotset;
533 /* Ion is not in any of the channels, so it must be in domain A or B */
536 *comp_now = eDomainA;
540 *comp_now = eDomainB;
544 /* Only take action, if ion is now in domain A or B, and was before
545 * in the other domain!
547 if (eDomainNotset == *comp_from)
549 /* Maybe we can set the domain now */
550 *comp_from = *comp_now; /* Could still be eDomainNotset, though */
552 else if ( (*comp_now != eDomainNotset ) /* if in channel */
553 && (*comp_from != *comp_now) )
555 /* Obviously the ion changed its domain.
556 * Count this for the channel through which it has passed. */
557 switch (*channel_label)
559 case eChHistPassedNone:
560 *s->fluxleak = *s->fluxleak + 1;
562 fprintf(stderr, " %s Warning! Step %s, ion %d (%s) moved from %s to %s\n",
563 SwS, gmx_step_str(step, buf), iion, IonStr[iontype], DomainString[*comp_from], DomainString[*comp_now]);
566 fprintf(stderr, ", possibly due to a swap in the original simulation.\n");
570 fprintf(stderr, "but did not pass cyl0 or cyl1 as defined in the .mdp file.\n"
571 "Do you have an ion somewhere within the membrane?\n");
572 /* Write this info to the CompEL output file: */
573 fprintf(s->fpout, " # Warning: step %s, ion %d (%s) moved from %s to %s (probably through the membrane)\n",
574 gmx_step_str(step, buf), iion, IonStr[iontype],
575 DomainString[*comp_from], DomainString[*comp_now]);
579 case eChHistPassedCh0:
580 case eChHistPassedCh1:
581 if (*channel_label == eChHistPassedCh0)
590 if (eDomainA == *comp_from)
592 s->fluxfromAtoB[chan_nr][iontype]++;
596 s->fluxfromAtoB[chan_nr][iontype]--;
598 fprintf(fpout, "# Atom nr. %d finished passing %s.\n", iion, ChannelString[*channel_label]);
601 gmx_fatal(FARGS, "%s Unknown channel history entry!\n", SwS);
605 /* This ion has moved to the _other_ compartment ... */
606 *comp_from = *comp_now;
607 /* ... and it did not pass any channel yet */
608 *channel_label = eChHistPassedNone;
613 /*! \brief Get the lists of ions for the two compartments */
614 static void compartmentalize_ions(
627 real cyl0_r2, cyl1_r2;
629 int sum, not_in_comp[eCompNR]; /* consistency check */
634 iong = &s->group[eGrpIons];
637 cyl0_r2 = sc->cyl0r * sc->cyl0r;
638 cyl1_r2 = sc->cyl1r * sc->cyl1r;
641 /* Get us a counter that cycles in the range of [0 ... sc->nAverage[ */
642 replace = (step/sc->nstswap) % sc->nAverage;
644 for (comp = eCompA; comp <= eCompB; comp++)
646 /* Get lists of atoms that match criteria for this compartment */
647 get_compartment_boundaries(comp, sc->si_priv, box, &left, &right);
649 /* First clear the ion lists */
650 s->comp[comp][eIonNEG].nat = 0;
651 s->comp[comp][eIonPOS].nat = 0;
652 not_in_comp[comp] = 0; /* consistency check */
654 /* Loop over the IONS */
655 for (i = 0; i < iong->nat; i++)
657 /* Anion or cation? */
658 type = iong->qc[i] < 0 ? eIonNEG : eIonPOS;
660 /* Is this ion in the compartment that we look at? */
661 if (compartment_contains_atom(left, right, 0, iong->xc[i][sd], box[sd][sd], &dist) )
663 /* Now put it into the list containing only ions of its type */
664 add_to_list(i, &s->comp[comp][type], dist);
666 /* Master also checks through which channel each ion has passed */
667 if (MASTER(cr) && (iong->comp_now != NULL))
669 ion_nr_global = iong->ind[i] + 1; /* PDB index starts at 1 ... */
670 detect_flux_per_channel(ion_nr_global, comp, type, iong->xc[i],
671 &iong->comp_now[i], &iong->comp_from[i], &iong->channel_label[i],
672 sc, cyl0_r2, cyl1_r2, step, bRerun, fpout);
677 not_in_comp[comp] += 1;
680 /* Correct the time-averaged number of ions in both compartments */
681 update_time_window(&s->comp[comp][eIonNEG], sc->nAverage, replace);
682 update_time_window(&s->comp[comp][eIonPOS], sc->nAverage, replace);
685 /* Flux detection warnings */
691 "%s Warning: %d atoms were detected as being in both channels! Probably your split\n"
692 "%s cylinder is way too large, or one compartment has collapsed (step %"GMX_PRId64 ")\n",
693 SwS, s->cyl0and1, SwS, step);
695 fprintf(s->fpout, "Warning: %d atoms were assigned to both channels!\n", s->cyl0and1);
702 /* Consistency checks */
703 if (not_in_comp[eCompA] + not_in_comp[eCompB] != iong->nat)
707 fprintf(fpout, "# Warning: Inconsistency during ion compartmentalization. !inA: %d, !inB: %d, total ions %d\n",
708 not_in_comp[eCompA], not_in_comp[eCompB], iong->nat);
713 fprintf(stderr, "%s rank %d: Inconsistency during ion compartmentalization. !inA: %d, !inB: %d, total ions %d\n",
714 SwS, cr->nodeid, not_in_comp[eCompA], not_in_comp[eCompB], iong->nat);
718 sum = s->comp[eCompA][eIonNEG].nat + s->comp[eCompA][eIonPOS].nat
719 + s->comp[eCompB][eIonNEG].nat + s->comp[eCompB][eIonPOS].nat;
720 if (sum != iong->nat)
724 fprintf(fpout, "# Warning: %d atoms are in the ion group, but altogether %d have been assigned to the compartments.\n",
730 fprintf(stderr, "%s rank %d: %d atoms are in the ion group, but altogether %d have been assigned to the compartments.\n",
731 SwS, cr->nodeid, iong->nat, sum);
739 /*! \brief Set up the compartments and get lists of solvent atoms in each compartment */
740 static void compartmentalize_solvent(
752 int sum, not_in_comp[eCompNR]; /* consistency check */
756 solg = &s->group[eGrpSolvent];
760 for (comp = eCompA; comp <= eCompB; comp++)
762 /* Get lists of atoms that match criteria for this compartment */
763 get_compartment_boundaries(comp, sc->si_priv, box, &left, &right);
765 /* First clear the solvent molecule lists */
766 s->compsol[comp].nat = 0;
767 not_in_comp[comp] = 0; /* consistency check */
769 /* Loop over the solvent MOLECULES */
770 for (i = 0; i < sc->nat_sol; i += apm)
772 if (compartment_contains_atom(left, right, 0, solg->xc[i][sd], box[sd][sd], &dist))
774 /* Add the whole molecule to the list */
775 for (j = 0; j < apm; j++)
777 add_to_list(i+j, &s->compsol[comp], dist);
782 not_in_comp[comp] += apm;
789 fprintf(fpout, "# Solv. molecules in comp.%s: %d comp.%s: %d\n",
790 CompStr[eCompA], s->compsol[eCompA].nat/apm,
791 CompStr[eCompB], s->compsol[eCompB].nat/apm);
794 /* Consistency checks */
795 if (not_in_comp[eCompA] + not_in_comp[eCompB] != solg->nat)
799 fprintf(fpout, "# Warning: Inconsistency during solvent compartmentalization. !inA: %d, !inB: %d, solvent atoms %d\n",
800 not_in_comp[eCompA], not_in_comp[eCompB], solg->nat);
805 fprintf(stderr, "%s rank %d: Inconsistency during solvent compartmentalization. !inA: %d, !inB: %d, solvent atoms %d\n",
806 SwS, cr->nodeid, not_in_comp[eCompA], not_in_comp[eCompB], solg->nat);
809 sum = s->compsol[eCompA].nat + s->compsol[eCompB].nat;
810 if (sum != solg->nat)
814 fprintf(fpout, "# Warning: %d atoms in solvent group, but %d have been assigned to the compartments.\n",
820 fprintf(stderr, "%s rank %d: %d atoms in solvent group, but %d have been assigned to the compartments.\n",
821 SwS, cr->nodeid, solg->nat, sum);
827 /*! \brief Find out how many group atoms are in the compartments initially */
828 static void get_initial_ioncounts(
830 rvec x[], /* the initial positions */
844 /* Copy the initial swap group positions to the collective array so
845 * that we can compartmentalize */
846 for (i = 0; i < sc->nat; i++)
849 copy_rvec(x[ind], s->group[eGrpIons].xc[i]);
852 /* Set up the compartments and get lists of atoms in each compartment */
853 compartmentalize_ions(cr, sc, box, 0, s->fpout, bRerun);
855 /* Set initial concentrations if requested */
856 for (ic = 0; ic < eCompNR; ic++)
858 s->comp[ic][eIonPOS].nat_req = sc->ncations[ic];
859 s->comp[ic][eIonNEG].nat_req = sc->nanions[ic];
861 for (ic = 0; ic < eCompNR; ic++)
863 for (ii = 0; ii < eIonNR; ii++)
865 if (s->comp[ic][ii].nat_req < 0)
867 s->comp[ic][ii].nat_req = s->comp[ic][ii].nat;
872 /* Check whether the number of requested ions adds up to the total number of ions */
873 for (ii = 0; ii < eIonNR; ii++)
875 req[ii] = s->comp[eCompA][ii].nat_req + s->comp[eCompB][ii].nat_req;
876 tot[ii] = s->comp[eCompA][ii].nat + s->comp[eCompB][ii].nat;
878 if ( (req[eCompA] != tot[eCompA]) || (req[eCompB] != tot[eCompB ]) )
880 gmx_fatal(FARGS, "Mismatch of the number of ions summed over both compartments.\n"
881 "You requested a total of %d anions and %d cations,\n"
882 "but there are a total of %d anions and %d cations in the system.\n",
883 req[eIonNEG], req[eIonPOS],
884 tot[eIonNEG], tot[eIonPOS]);
887 /* Initialize time-averaging:
888 * Write initial concentrations to all time bins to start with */
889 for (ic = 0; ic < eCompNR; ic++)
891 for (ii = 0; ii < eIonNR; ii++)
893 s->comp[ic][ii].nat_av = s->comp[ic][ii].nat;
894 for (i = 0; i < sc->nAverage; i++)
896 s->comp[ic][ii].nat_past[i] = s->comp[ic][ii].nat;
903 /*! \brief Copy history of ion counts from checkpoint file.
905 * When called, the checkpoint file has already been read in. Here we copy
906 * over the values from .cpt file to the swap data structure.
908 static void get_initial_ioncounts_from_cpt(
909 t_inputrec *ir, swapstate_t *swapstate,
910 t_commrec *cr, gmx_bool bVerbose)
921 /* Copy the past values from the checkpoint values that have been read in already */
924 fprintf(stderr, "%s Copying values from checkpoint\n", SwS);
927 for (ic = 0; ic < eCompNR; ic++)
929 for (ii = 0; ii < eIonNR; ii++)
931 s->comp[ic][ii].nat_req = swapstate->nat_req[ic][ii];
932 s->comp[ic][ii].inflow_netto = swapstate->inflow_netto[ic][ii];
936 fprintf(stderr, "%s ... Influx netto: %d Requested: %d Past values: ", SwS,
937 s->comp[ic][ii].inflow_netto, s->comp[ic][ii].nat_req);
940 for (j = 0; j < sc->nAverage; j++)
942 s->comp[ic][ii].nat_past[j] = swapstate->nat_past[ic][ii][j];
945 fprintf(stderr, "%d ", s->comp[ic][ii].nat_past[j]);
950 fprintf(stderr, "\n");
958 /*! \brief The master lets all others know about the initial ion counts. */
959 static void bc_initial_concentrations(
968 for (ic = 0; ic < eCompNR; ic++)
970 for (ii = 0; ii < eIonNR; ii++)
972 gmx_bcast(sizeof(s->comp[ic][ii].nat_req), &(s->comp[ic][ii].nat_req), cr);
973 gmx_bcast(sizeof(s->comp[ic][ii].nat ), &(s->comp[ic][ii].nat ), cr);
974 gmx_bcast( swap->nAverage * sizeof(s->comp[ic][ii].nat_past[0]), s->comp[ic][ii].nat_past, cr);
980 /*! \brief Ensure that each atom belongs to at most one of the swap groups. */
981 static void check_swap_groups(t_swap *s, int nat, gmx_bool bVerbose)
985 atom_id *nGroup = NULL; /* This array counts for each atom in the MD system to
986 how many swap groups it belongs (should be 0 or 1!) */
988 int nMultiple = 0; /* Number of atoms belonging to multiple groups */
993 fprintf(stderr, "%s Making sure each atom belongs to at most one of the swap groups.\n", SwS);
996 /* Add one to the group count of atoms belonging to a swap group: */
998 for (i = 0; i < eGrpNr; i++)
1001 for (j = 0; j < g->nat; j++)
1003 /* Get the global index of this atom of this group: */
1008 /* Make sure each atom belongs to at most one swap group: */
1009 for (j = 0; j < g->nat; j++)
1020 gmx_fatal(FARGS, "%s Cannot perform swapping since %d atom%s allocated to more than one swap index group.\n"
1021 "%s Each atom must be allocated to at most one of the split groups, the swap group, or the solvent.\n"
1022 "%s Check the .mdp file settings regarding the swap index groups or the index groups themselves.\n",
1023 SwS, nMultiple, (1 == nMultiple) ? " is" : "s are", SwSEmpty, SwSEmpty);
1028 /*! \brief Get the number of atoms per molecule for this group.
1030 * Also ensure that all the molecules in this group have this number of atoms.
1032 static int get_group_apm_check(
1036 const gmx_mtop_atomlookup_t alook,
1041 int molb, molnr, atnr_mol;
1044 ind = s->group[group].ind;
1045 nat = s->group[group].nat;
1047 /* Determine the number of solvent atoms per solvent molecule from the
1048 * first solvent atom: */
1050 gmx_mtop_atomnr_to_molblock_ind(alook, ind[i], &molb, &molnr, &atnr_mol);
1051 apm = mtop->molblock[molb].natoms_mol;
1055 fprintf(stderr, "%s Checking whether all %s molecules consist of %d atom%s\n",
1056 SwS, GrpString[group], apm, apm > 1 ? "s" : "");
1059 /* Check whether this is also true for all other solvent atoms */
1060 for (i = 1; i < nat; i++)
1062 gmx_mtop_atomnr_to_molblock_ind(alook, ind[i], &molb, &molnr, &atnr_mol);
1063 if (apm != mtop->molblock[molb].natoms_mol)
1065 gmx_fatal(FARGS, "Not all %s group molecules consist of %d atoms.",
1066 GrpString[group], apm);
1074 /*! \brief Print the legend to the swap output file.
1076 * Also print the initial ion counts
1078 static void print_ionlist_legend(t_inputrec *ir, const output_env_t oenv)
1080 const char **legend;
1086 s = ir->swap->si_priv;
1088 snew(legend, eCompNR*eIonNR*3 + 2 + eChanNR*eIonNR + 1);
1089 for (ic = count = 0; ic < eCompNR; ic++)
1091 for (ii = 0; ii < eIonNR; ii++)
1093 sprintf(buf, "%s %ss", CompStr[ic], IonString[ii]);
1094 legend[count++] = gmx_strdup(buf);
1095 sprintf(buf, "%s av. mismatch to %d%s",
1096 CompStr[ic], s->comp[ic][ii].nat_req, IonStr[ii]);
1097 legend[count++] = gmx_strdup(buf);
1098 sprintf(buf, "%s netto %s influx", CompStr[ic], IonString[ii]);
1099 legend[count++] = gmx_strdup(buf);
1102 sprintf(buf, "%scenter of %s of split group 0", SwapStr[ir->eSwapCoords], (NULL != s->group[eGrpSplit0].m) ? "mass" : "geometry");
1103 legend[count++] = gmx_strdup(buf);
1104 sprintf(buf, "%scenter of %s of split group 1", SwapStr[ir->eSwapCoords], (NULL != s->group[eGrpSplit1].m) ? "mass" : "geometry");
1105 legend[count++] = gmx_strdup(buf);
1107 for (ic = 0; ic < eChanNR; ic++)
1109 for (ii = 0; ii < eIonNR; ii++)
1111 sprintf(buf, "A->ch%d->B %s permeations", ic, IonString[ii]);
1112 legend[count++] = gmx_strdup(buf);
1116 sprintf(buf, "leakage");
1117 legend[count++] = gmx_strdup(buf);
1119 xvgr_legend(s->fpout, count, legend, oenv);
1121 fprintf(s->fpout, "# Instantaneous ion counts and time-averaged differences to requested numbers\n");
1122 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 ");
1123 fprintf(s->fpout, " %s-Split0 %s-Split1", DimStr[s->swapdim], DimStr[s->swapdim]);
1124 fprintf(s->fpout, " A-ch0-B_an A-ch0-B_cat A-ch1-B_an A-ch1-B_cat ion_leakage\n");
1129 /*! \brief Initialize channel ion flux detection routine.
1131 * Initialize arrays that keep track of where the ions come from and where
1134 static void detect_flux_per_channel_init(
1137 swapstate_t *swapstate,
1138 gmx_bool bStartFromCpt)
1144 g = &(s->group[eGrpIons]);
1146 /* All these flux detection routines run on the master only */
1150 g->comp_from = NULL;
1151 g->channel_label = NULL;
1156 /******************************************************/
1157 /* Channel and domain history for the individual ions */
1158 /******************************************************/
1159 if (bStartFromCpt) /* set the pointers right */
1161 g->comp_from = swapstate->comp_from;
1162 g->channel_label = swapstate->channel_label;
1164 else /* allocate memory */
1166 snew(g->comp_from, g->nat);
1167 swapstate->comp_from = g->comp_from;
1168 snew(g->channel_label, g->nat);
1169 swapstate->channel_label = g->channel_label;
1171 snew(g->comp_now, g->nat);
1173 /* Initialize the channel and domain history counters */
1174 for (i = 0; i < g->nat; i++)
1176 g->comp_now[i] = eDomainNotset;
1179 g->comp_from[i] = eDomainNotset;
1180 g->channel_label[i] = eChHistPassedNone;
1184 /************************************/
1185 /* Channel fluxes for both channels */
1186 /************************************/
1193 fprintf(stderr, "%s Copying channel fluxes from checkpoint file data\n", SwS);
1196 for (ic = 0; ic < eChanNR; ic++)
1198 fprintf(stderr, "%s Channel %d flux history: ", SwS, ic);
1199 for (ii = 0; ii < eIonNR; ii++)
1203 s->fluxfromAtoB[ic][ii] = swapstate->fluxfromAtoB[ic][ii];
1207 s->fluxfromAtoB[ic][ii] = 0;
1210 fprintf(stderr, "%d %s%s ", s->fluxfromAtoB[ic][ii], IonString[ii], s->fluxfromAtoB[ic][ii] == 1 ? "" : "s");
1212 fprintf(stderr, "\n");
1216 s->fluxleak = swapstate->fluxleak;
1220 snew(s->fluxleak, 1);
1222 /* Set pointer for checkpoint writing */
1223 swapstate->fluxleak = s->fluxleak;
1226 /* Set pointers for checkpoint writing */
1227 for (ic = 0; ic < eChanNR; ic++)
1229 for (ii = 0; ii < eIonNR; ii++)
1231 swapstate->fluxfromAtoB_p[ic][ii] = &(s->fluxfromAtoB[ic][ii]);
1237 /*! \brief Outputs the initial structure to PDB file for debugging reasons.
1239 * Output the starting structure so that in case of multimeric channels
1240 * the user can check whether we have the correct PBC image for all atoms.
1241 * If this is not correct, the ion counts per channel will be very likely
1244 static void outputStartStructureIfWanted(gmx_mtop_t *mtop, rvec *x, int ePBC, matrix box)
1246 char *env = getenv("GMX_COMPELDUMP");
1250 fprintf(stderr, "\n%s Found env.var. GMX_COMPELDUMP, will output CompEL starting structure made whole.\n"
1251 "%s In case of multimeric channels, please check whether they have the correct PBC representation.\n",
1254 write_sto_conf_mtop("CompELAssumedWholeConfiguration.pdb", *mtop->name, mtop, x, NULL, ePBC, box);
1259 /*! \brief Initialize the swapstate structure, used for checkpoint writing.
1261 * The swapstate struct stores the information we need to make the channels
1262 * whole again after restarts from a checkpoint file. Here we do the following:\n
1263 * a) If we did not start from .cpt, we prepare the struct for proper .cpt writing,\n
1264 * b) if we did start from .cpt, we copy over the last whole structures from .cpt,\n
1265 * c) in any case, for subsequent checkpoint writing, we set the pointers in\n
1266 * swapstate to the x_old arrays, which contain the correct PBC representation of
1267 * multimeric channels at the last time step.
1269 static void init_swapstate(
1270 swapstate_t *swapstate,
1273 rvec x[], /* the initial positions */
1278 rvec *x_pbc = NULL; /* positions of the whole MD system with molecules made whole */
1285 /* We always need the last whole positions such that
1286 * in the next time step we can make the channels whole again in PBC */
1287 if (swapstate->bFromCpt)
1289 /* Copy the last whole positions of each channel from .cpt */
1290 g = &(s->group[eGrpSplit0]);
1291 for (i = 0; i < g->nat; i++)
1293 copy_rvec(swapstate->xc_old_whole[eChan0][i], g->xc_old[i]);
1295 g = &(s->group[eGrpSplit1]);
1296 for (i = 0; i < g->nat; i++)
1298 copy_rvec(swapstate->xc_old_whole[eChan1][i], g->xc_old[i]);
1303 /* Extract the initial split group positions. */
1305 /* Remove pbc, make molecule whole. */
1306 snew(x_pbc, mtop->natoms);
1307 m_rveccopy(mtop->natoms, x, x_pbc);
1309 /* This can only make individual molecules whole, not multimers */
1310 do_pbc_mtop(NULL, ePBC, box, mtop, x_pbc);
1312 /* Output the starting structure? */
1313 outputStartStructureIfWanted(mtop, x_pbc, ePBC, box);
1315 /* If this is the first run (i.e. no checkpoint present) we assume
1316 * that the starting positions give us the correct PBC representation */
1317 for (ig = eGrpSplit0; ig <= eGrpSplit1; ig++)
1319 g = &(s->group[ig]);
1320 for (i = 0; i < g->nat; i++)
1322 copy_rvec(x_pbc[g->ind[i]], g->xc_old[i]);
1327 /* Prepare swapstate arrays for later checkpoint writing */
1328 swapstate->nat[eChan0] = s->group[eGrpSplit0].nat;
1329 swapstate->nat[eChan1] = s->group[eGrpSplit1].nat;
1332 /* For subsequent checkpoint writing, set the swapstate pointers to the xc_old
1333 * arrays that get updated at every swapping step */
1334 swapstate->xc_old_whole_p[eChan0] = &s->group[eGrpSplit0].xc_old;
1335 swapstate->xc_old_whole_p[eChan1] = &s->group[eGrpSplit1].xc_old;
1339 extern void init_swapcoords(
1347 swapstate_t *swapstate,
1349 const output_env_t oenv,
1350 unsigned long Flags)
1352 int i, ic, ig, ii, j;
1357 gmx_bool bAppend, bStartFromCpt, bRerun;
1358 gmx_mtop_atomlookup_t alook = NULL;
1361 alook = gmx_mtop_atomlookup_init(mtop);
1363 if ( (PAR(cr)) && !DOMAINDECOMP(cr) )
1365 gmx_fatal(FARGS, "Position swapping is only implemented for domain decomposition!");
1368 bAppend = Flags & MD_APPENDFILES;
1369 bStartFromCpt = Flags & MD_STARTFROMCPT;
1370 bRerun = Flags & MD_RERUN;
1373 snew(sc->si_priv, 1);
1380 gmx_fatal(FARGS, "%s This module does not support reruns in parallel\nPlease request a serial run with -nt 1 / -np 1\n", SwS);
1383 fprintf(stderr, "%s Rerun - using every available frame\n", SwS);
1385 sc->nAverage = 1; /* averaging makes no sense for reruns */
1388 if (MASTER(cr) && !bAppend)
1390 fprintf(fplog, "\nInitializing ion/water position exchanges\n");
1391 please_cite(fplog, "Kutzner2011b");
1394 switch (ir->eSwapCoords)
1410 /* Copy some data to the group structures for convenience */
1411 /* Number of atoms in the group */
1412 s->group[eGrpIons ].nat = sc->nat;
1413 s->group[eGrpSplit0 ].nat = sc->nat_split[0];
1414 s->group[eGrpSplit1 ].nat = sc->nat_split[1];
1415 s->group[eGrpSolvent].nat = sc->nat_sol;
1416 /* Pointer to the indices */
1417 s->group[eGrpIons ].ind = sc->ind;
1418 s->group[eGrpSplit0 ].ind = sc->ind_split[0];
1419 s->group[eGrpSplit1 ].ind = sc->ind_split[1];
1420 s->group[eGrpSolvent].ind = sc->ind_sol;
1422 check_swap_groups(s, mtop->natoms, bVerbose && MASTER(cr));
1424 /* Allocate space for the collective arrays for all groups */
1425 for (ig = 0; ig < eGrpNr; ig++)
1427 g = &(s->group[ig]);
1428 snew(g->xc, g->nat);
1429 snew(g->c_ind_loc, g->nat);
1430 /* For the split groups (the channels) we need some extra memory to
1431 * be able to make the molecules whole even if they span more than
1432 * half of the box size. */
1433 if (eGrpSplit0 == ig || eGrpSplit1 == ig)
1435 snew(g->xc_shifts, g->nat);
1436 snew(g->xc_eshifts, g->nat);
1437 snew(g->xc_old, g->nat);
1443 init_swapstate(swapstate, sc, mtop, x, box, ir->ePBC);
1446 /* After init_swapstate we have a set of (old) whole positions for our
1447 * channels. Now transfer that to all nodes */
1450 for (ig = eGrpSplit0; ig <= eGrpSplit1; ig++)
1452 g = &(s->group[ig]);
1453 gmx_bcast((g->nat)*sizeof((g->xc_old)[0]), g->xc_old, (cr));
1457 /* Make sure that all molecules in the ion and solvent groups contain the
1458 * same number of atoms each */
1459 s->group[eGrpIons ].apm = get_group_apm_check(eGrpIons, s, MASTER(cr) && bVerbose, alook, mtop);
1460 s->group[eGrpSolvent].apm = get_group_apm_check(eGrpSolvent, s, MASTER(cr) && bVerbose, alook, mtop);
1462 /* Save masses where needed */
1463 s->group[eGrpIons ].m = NULL;
1464 /* We only need enough space to determine a single solvent molecule's
1465 * center at at time */
1466 g = &(s->group[eGrpSolvent]);
1469 /* Need mass-weighted center of split group? */
1470 for (j = 0, ig = eGrpSplit0; j < eChanNR; ig++, j++)
1472 g = &(s->group[ig]);
1473 if (TRUE == sc->massw_split[j])
1475 /* Save the split group charges if mass-weighting is requested */
1477 for (i = 0; i < g->nat; i++)
1479 gmx_mtop_atomnr_to_atom(alook, g->ind[i], &atom);
1489 /* Save the ionic charges */
1490 g = &(s->group[eGrpIons]);
1491 snew(g->qc, g->nat);
1492 for (i = 0; i < g->nat; i++)
1494 gmx_mtop_atomnr_to_atom(alook, g->ind[i], &atom);
1499 set_pbc(s->pbc, -1, box);
1506 fprintf(stderr, "%s Opening output file %s%s\n", SwS, fn, bAppend ? " for appending" : "");
1509 s->fpout = gmx_fio_fopen(fn, bAppend ? "a" : "w" );
1513 xvgr_header(s->fpout, "Ion counts", "Time (ps)", "counts", exvggtXNY, oenv);
1515 for (ig = 0; ig < eGrpNr; ig++)
1517 g = &(s->group[ig]);
1518 fprintf(s->fpout, "# %s group contains %d atom%s", GrpString[ig], g->nat, (g->nat > 1) ? "s" : "");
1519 if (eGrpSolvent == ig || eGrpIons == ig)
1521 fprintf(s->fpout, " with %d atom%s in each molecule", g->apm, (g->apm > 1) ? "s" : "");
1523 fprintf(s->fpout, ".\n");
1526 fprintf(s->fpout, "#\n# Initial positions of split groups:\n");
1529 for (j = 0, ig = eGrpSplit0; j < eChanNR; j++, ig++)
1531 g = &(s->group[ig]);
1532 for (i = 0; i < g->nat; i++)
1534 copy_rvec(x[sc->ind_split[j][i]], g->xc[i]);
1536 if (eGrpSplit0 == ig || eGrpSplit1 == ig)
1538 /* xc has the correct PBC representation for the two channels, so we do
1539 * not need to correct for that */
1540 get_center(g->xc, g->m, g->nat, g->center);
1544 /* For the water molecules, we need to make the molecules whole */
1545 get_molecule_center(g->xc, g->nat, g->m, g->center, s->pbc);
1549 fprintf(s->fpout, "# %s group %s-center %5f nm\n", GrpString[ig],
1550 DimStr[s->swapdim], g->center[s->swapdim]);
1556 fprintf(s->fpout, "#\n");
1557 fprintf(s->fpout, "# split0 cylinder radius %f nm, up %f nm, down %f nm\n",
1558 sc->cyl0r, sc->cyl0u, sc->cyl0l);
1559 fprintf(s->fpout, "# split1 cylinder radius %f nm, up %f nm, down %f nm\n",
1560 sc->cyl1r, sc->cyl1u, sc->cyl1l);
1565 fprintf(s->fpout, "#\n");
1568 fprintf(s->fpout, "# Coupling constant (number of swap attempt steps to average over): %d (translates to %f ps).\n",
1569 sc->nAverage, sc->nAverage*sc->nstswap*ir->delta_t);
1570 fprintf(s->fpout, "# Threshold is %f\n", sc->threshold);
1571 fprintf(s->fpout, "#\n");
1572 fprintf(s->fpout, "# Remarks about which atoms passed which channel use global atoms numbers starting at one.\n");
1581 /* Prepare for parallel or serial run */
1584 for (ig = 0; ig < eGrpNr; ig++)
1586 g = &(s->group[ig]);
1594 for (ig = 0; ig < eGrpNr; ig++)
1596 g = &(s->group[ig]);
1597 g->nat_loc = g->nat;
1598 g->ind_loc = g->ind;
1599 /* c_ind_loc needs to be set to identity in the serial case */
1600 for (i = 0; i < g->nat; i++)
1602 g->c_ind_loc[i] = i;
1607 /* Allocate memory for the ion counts time window */
1608 for (ic = 0; ic < eCompNR; ic++)
1610 for (ii = 0; ii < eIonNR; ii++)
1612 snew(s->comp[ic][ii].nat_past, sc->nAverage);
1616 /* Get the initial ion concentrations and let the other nodes know */
1619 swapstate->nions = s->group[eGrpIons].nat;
1623 get_initial_ioncounts_from_cpt(ir, swapstate, cr, bVerbose);
1627 fprintf(stderr, "%s Determining initial ion counts.\n", SwS);
1628 get_initial_ioncounts(ir, x, box, cr, bRerun);
1631 /* Prepare (further) checkpoint writes ... */
1634 /* Consistency check */
1635 if (swapstate->nAverage != sc->nAverage)
1637 gmx_fatal(FARGS, "%s Ion count averaging steps mismatch! checkpoint: %d, tpr: %d",
1638 SwS, swapstate->nAverage, sc->nAverage);
1643 swapstate->nAverage = sc->nAverage;
1645 fprintf(stderr, "%s Setting pointers for checkpoint writing\n", SwS);
1646 for (ic = 0; ic < eCompNR; ic++)
1648 for (ii = 0; ii < eIonNR; ii++)
1650 swapstate->nat_req_p[ic][ii] = &(s->comp[ic][ii].nat_req);
1651 swapstate->nat_past_p[ic][ii] = &(s->comp[ic][ii].nat_past[0]);
1652 swapstate->inflow_netto_p[ic][ii] = &(s->comp[ic][ii].inflow_netto);
1656 /* Determine the total charge imbalance */
1657 s->deltaQ = ( (-1) * s->comp[eCompA][eIonNEG].nat_req + s->comp[eCompA][eIonPOS].nat_req )
1658 - ( (-1) * s->comp[eCompB][eIonNEG].nat_req + s->comp[eCompB][eIonPOS].nat_req );
1662 fprintf(stderr, "%s Requested charge imbalance is Q(A) - Q(B) = %gz.\n", SwS, s->deltaQ);
1666 fprintf(s->fpout, "# Requested charge imbalance is Q(A)-Q(B) = %gz.\n", s->deltaQ);
1672 bc_initial_concentrations(cr, ir->swap);
1675 /* Put the time-averaged number of ions for all compartments */
1676 for (ic = 0; ic < eCompNR; ic++)
1678 for (ii = 0; ii < eIonNR; ii++)
1680 update_time_window(&(s->comp[ic][ii]), sc->nAverage, -1);
1684 /* Initialize arrays that keep track of through which channel the ions go */
1685 detect_flux_per_channel_init(cr, s, swapstate, bStartFromCpt);
1687 /* We need to print the legend if we open this file for the first time. */
1688 if (MASTER(cr) && !bAppend)
1690 print_ionlist_legend(ir, oenv);
1695 extern void dd_make_local_swap_groups(gmx_domdec_t *dd, t_swapcoords *sc)
1701 /* Make ion group, split groups and solvent group */
1702 for (ig = 0; ig < eGrpNr; ig++)
1704 g = &(sc->si_priv->group[ig]);
1705 dd_make_local_group_indices(dd->ga2la, g->nat, g->ind,
1706 &(g->nat_loc), &(g->ind_loc), &(g->nalloc_loc), g->c_ind_loc);
1711 /*! \brief Do we need to swap ions with water molecules at this step?
1713 * From the requested and average ion counts we determine whether a swap is needed
1714 * at this time step.
1716 static gmx_bool need_swap(t_swapcoords *sc)
1723 for (ic = 0; ic < eCompNR; ic++)
1725 for (ii = 0; ii < eIonNR; ii++)
1727 if (s->comp[ic][ii].nat_req - s->comp[ic][ii].nat_av >= sc->threshold)
1737 /*! \brief Return index of atom that we can use for swapping.
1739 * Returns the index of an atom that is far off the compartment boundaries.
1740 * Other atoms of the molecule (if any) will directly follow the returned index
1742 static int get_index_of_distant_atom(
1743 t_compartment *comp,
1744 int apm) /* Atoms per molecule - just return the first atom index of a molecule */
1747 real d = GMX_REAL_MAX;
1750 /* comp->nat contains the original number of atoms in this compartment
1751 * prior to doing any swaps. Some of these atoms may already have been
1752 * swapped out, but then they are marked with a distance of GMX_REAL_MAX
1754 for (i = 0; i < comp->nat_old; i += apm)
1756 if (comp->dist[i] < d)
1759 d = comp->dist[ibest];
1765 gmx_fatal(FARGS, "Could not get index of swap atom. Compartment atoms %d before swaps, atoms per molecule %d.",
1766 comp->nat_old, apm);
1769 /* Set the distance of this index to infinity such that it won't get selected again in
1772 comp->dist[ibest] = GMX_REAL_MAX;
1774 return comp->ind[ibest];
1778 /*! \brief Swaps centers of mass and makes molecule whole if broken */
1779 static void translate_positions(
1787 rvec reference, dx, correctPBCimage;
1790 /* Use the first atom as the reference for PBC */
1791 copy_rvec(x[0], reference);
1793 for (i = 0; i < apm; i++)
1795 /* PBC distance between position and reference */
1796 pbc_dx(pbc, x[i], reference, dx);
1798 /* Add PBC distance to reference */
1799 rvec_add(reference, dx, correctPBCimage);
1801 /* Subtract old_com from correct image and add new_com */
1802 rvec_dec(correctPBCimage, old_com);
1803 rvec_inc(correctPBCimage, new_com);
1805 copy_rvec(correctPBCimage, x[i]);
1810 /*! \brief Write back the the modified local positions from the collective array to the official positions. */
1811 static void apply_modified_positions(
1818 for (l = 0; l < g->nat_loc; l++)
1820 /* Get the right local index to write to */
1822 /* Where is the local atom in the collective array? */
1823 cind = g->c_ind_loc[l];
1825 /* Copy the possibly modified position */
1826 copy_rvec(g->xc[cind], x[ii]);
1831 extern gmx_bool do_swapcoords(
1836 gmx_wallcycle_t wcycle,
1845 int j, ii, ic, ig, im, gmax, nswaps;
1846 gmx_bool bSwap = FALSE;
1848 real vacancy[eCompNR][eIonNR];
1850 rvec solvent_center, ion_center;
1852 gmx_mtop_atomlookup_t alook = NULL;
1855 wallcycle_start(wcycle, ewcSWAP);
1860 /* Assemble all the positions of the swap group (ig = 0), the split groups
1861 * (ig = 1,2), and possibly the solvent group (ig = 3) */
1864 for (ig = 0; ig < gmax; ig++)
1866 g = &(s->group[ig]);
1867 if (eGrpSplit0 == ig || eGrpSplit1 == ig)
1869 /* The split groups, i.e. the channels. Here we need the full
1870 * communicate_group_positions(), so that we can make the molecules
1871 * whole even in cases where they span more than half of the box in
1873 communicate_group_positions(cr, g->xc, g->xc_shifts, g->xc_eshifts, TRUE,
1874 x, g->nat, g->nat_loc, g->ind_loc, g->c_ind_loc, g->xc_old, box);
1876 get_center(g->xc, g->m, g->nat, g->center); /* center of split groups == channels */
1880 /* Swap group (ions), and solvent group. These molecules are small
1881 * and we can always make them whole with a simple distance check.
1882 * Therefore we pass NULL as third argument. */
1883 communicate_group_positions(cr, g->xc, NULL, NULL, FALSE,
1884 x, g->nat, g->nat_loc, g->ind_loc, g->c_ind_loc, NULL, NULL);
1888 /* Set up the compartments and get lists of atoms in each compartment,
1889 * determine how many ions each compartment contains */
1890 compartmentalize_ions(cr, sc, box, step, s->fpout, bRerun);
1892 /* Output how many ions are in the compartments */
1895 print_ionlist(s, t, "");
1898 /* If we are doing a rerun, we are finished here, since we cannot perform
1905 /* Do we have to perform a swap? */
1906 bSwap = need_swap(sc);
1909 g = &(s->group[eGrpSolvent]);
1910 communicate_group_positions(cr, g->xc, NULL, NULL, FALSE,
1911 x, g->nat, g->nat_loc, g->ind_loc, g->c_ind_loc, NULL, NULL);
1913 compartmentalize_solvent(cr, sc, box, s->fpout);
1915 /* Determine where ions are missing and where ions are too many */
1916 for (ic = 0; ic < eCompNR; ic++)
1918 for (ii = 0; ii < eIonNR; ii++)
1920 vacancy[ic][ii] = s->comp[ic][ii].nat_req - s->comp[ic][ii].nat_av;
1924 /* Remember the original number of ions per compartment */
1925 for (ic = 0; ic < eCompNR; ic++)
1927 s->compsol[ic].nat_old = s->compsol[ic].nat;
1928 for (ii = 0; ii < eIonNR; ii++)
1930 s->comp[ic][ii].nat_old = s->comp[ic][ii].nat;
1934 /* Now actually correct the number of ions */
1935 g = &(s->group[eGrpSolvent]);
1937 alook = gmx_mtop_atomlookup_init(mtop);
1938 for (ic = 0; ic < eCompNR; ic++)
1940 for (ii = 0; ii < eIonNR; ii++)
1942 while (vacancy[ic][ii] >= sc->threshold)
1944 /* Swap in an ion */
1946 /* Get the xc-index of the first atom of a solvent molecule of this compartment */
1947 isol = get_index_of_distant_atom(&(s->compsol[ic]), s->group[eGrpSolvent].apm );
1949 /* Get the xc-index of an ion from the other compartment */
1950 iion = get_index_of_distant_atom(&(s->comp[(ic+1)%eCompNR][ii]), s->group[eGrpIons].apm );
1952 /* Get the solvent molecule's center of mass */
1953 for (im = 0; im < s->group[eGrpSolvent].apm; im++)
1955 gmx_mtop_atomnr_to_atom(alook, s->group[eGrpSolvent].ind[isol+im], &atom);
1956 s->group[eGrpSolvent].m[im] = atom->m;
1958 get_molecule_center(&(s->group[eGrpSolvent].xc[isol]), s->group[eGrpSolvent].apm, s->group[eGrpSolvent].m, solvent_center, s->pbc);
1959 get_molecule_center(&(s->group[eGrpIons ].xc[iion]), s->group[eGrpIons ].apm, NULL, ion_center, s->pbc);
1961 /* subtract com_solvent and add com_ion */
1962 translate_positions(&(s->group[eGrpSolvent].xc[isol]), s->group[eGrpSolvent].apm, solvent_center, ion_center, s->pbc);
1963 /* For the ion, subtract com_ion and add com_solvent */
1964 translate_positions(&(s->group[eGrpIons ].xc[iion]), s->group[eGrpIons ].apm, ion_center, solvent_center, s->pbc);
1967 vacancy[(ic+1) % eCompNR][ii]++;
1969 /* Keep track of the changes */
1970 s->comp[ic ][ii].nat++;
1971 s->comp[(ic+1) % eCompNR][ii].nat--;
1972 s->comp[ic ][ii].inflow_netto++;
1973 s->comp[(ic+1) % eCompNR][ii].inflow_netto--;
1974 /* Correct the past time window to still get the right averages from now on */
1975 s->comp[ic ][ii].nat_av++;
1976 s->comp[(ic+1) % eCompNR][ii].nat_av--;
1977 for (j = 0; j < sc->nAverage; j++)
1979 s->comp[ic ][ii].nat_past[j]++;
1980 s->comp[(ic+1) % eCompNR][ii].nat_past[j]--;
1982 /* Clear ion history */
1985 s->group[eGrpIons].channel_label[iion] = eChHistPassedNone;
1986 s->group[eGrpIons].comp_from[iion] = eDomainNotset;
1988 /* That was the swap */
1993 gmx_mtop_atomlookup_destroy(alook);
1997 fprintf(stderr, "%s Performed %d swap%s in step %"GMX_PRId64 ".\n", SwS, nswaps, nswaps > 1 ? "s" : "", step);
1999 if (s->fpout != NULL)
2001 print_ionlist(s, t, " # after swap");
2004 /* Write back the the modified local positions from the collective array to the official coordinates */
2005 apply_modified_positions(&(s->group[eGrpIons ]), x);
2006 apply_modified_positions(&(s->group[eGrpSolvent]), x);
2007 } /* end of if(bSwap) */
2009 wallcycle_stop(wcycle, ewcSWAP);