2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017,2018,2019, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
47 #include "gromacs/commandline/pargs.h"
48 #include "gromacs/fileio/confio.h"
49 #include "gromacs/fileio/pdbio.h"
50 #include "gromacs/gmxlib/conformation_utilities.h"
51 #include "gromacs/gmxpreprocess/makeexclusiondistances.h"
52 #include "gromacs/math/functions.h"
53 #include "gromacs/math/units.h"
54 #include "gromacs/math/vec.h"
55 #include "gromacs/pbcutil/boxutilities.h"
56 #include "gromacs/pbcutil/pbc.h"
57 #include "gromacs/selection/nbsearch.h"
58 #include "gromacs/topology/atomprop.h"
59 #include "gromacs/topology/atoms.h"
60 #include "gromacs/topology/atomsbuilder.h"
61 #include "gromacs/topology/mtop_util.h"
62 #include "gromacs/topology/topology.h"
63 #include "gromacs/utility/arraysize.h"
64 #include "gromacs/utility/cstringutil.h"
65 #include "gromacs/utility/fatalerror.h"
66 #include "gromacs/utility/futil.h"
67 #include "gromacs/utility/gmxassert.h"
68 #include "gromacs/utility/smalloc.h"
80 static void sort_molecule(t_atoms **atoms_solvt,
85 int atnr, i, j, moltp = 0, nrmoltypes, resi_o, resi_n, resnr;
89 fprintf(stderr, "Sorting configuration\n");
93 /* copy each residue from *atoms to a molecule in *molecule */
96 for (i = 0; i < atoms->nr; i++)
98 if ( (i == 0) || (atoms->atom[i].resind != atoms->atom[i-1].resind) )
100 /* see if this was a molecule type we haven't had yet: */
102 for (j = 0; (j < nrmoltypes) && (moltp == -1); j++)
104 /* moltypes is guaranteed to be allocated because otherwise
105 * nrmoltypes is 0. */
106 if (strcmp(*(atoms->resinfo[atoms->atom[i].resind].name), moltypes[j].name) == 0)
115 srenew(moltypes, nrmoltypes);
116 moltypes[moltp].name = *(atoms->resinfo[atoms->atom[i].resind].name);
118 while ((i+atnr < atoms->nr) &&
119 (atoms->atom[i].resind == atoms->atom[i+atnr].resind))
123 moltypes[moltp].natoms = atnr;
124 moltypes[moltp].nmol = 0;
126 moltypes[moltp].nmol++;
130 fprintf(stderr, "Found %d%s molecule type%s:\n",
131 nrmoltypes, nrmoltypes == 1 ? "" : " different", nrmoltypes == 1 ? "" : "s");
132 for (j = 0; j < nrmoltypes; j++)
136 moltypes[j].res0 = 0;
140 moltypes[j].res0 = moltypes[j-1].res0+moltypes[j-1].nmol;
142 fprintf(stderr, "%7s (%4d atoms): %5d residues\n",
143 moltypes[j].name, moltypes[j].natoms, moltypes[j].nmol);
146 /* if we have only 1 moleculetype, we don't have to sort */
149 /* find out which molecules should go where: */
150 moltypes[0].i = moltypes[0].i0 = 0;
151 for (j = 1; j < nrmoltypes; j++)
155 moltypes[j-1].i0+moltypes[j-1].natoms*moltypes[j-1].nmol;
158 /* now put them there: */
160 init_t_atoms(*newatoms, atoms->nr, FALSE);
161 (*newatoms)->nres = atoms->nres;
162 srenew((*newatoms)->resinfo, atoms->nres);
163 std::vector<RVec> newx(x->size());
164 std::vector<RVec> newv(v->size());
169 for (moltp = 0; moltp < nrmoltypes; moltp++)
172 while (i < atoms->nr)
174 resi_o = atoms->atom[i].resind;
175 if (strcmp(*atoms->resinfo[resi_o].name, moltypes[moltp].name) == 0)
177 /* Copy the residue info */
178 (*newatoms)->resinfo[resi_n] = atoms->resinfo[resi_o];
179 (*newatoms)->resinfo[resi_n].nr = resnr;
180 /* Copy the atom info */
183 (*newatoms)->atom[j] = atoms->atom[i];
184 (*newatoms)->atomname[j] = atoms->atomname[i];
185 (*newatoms)->atom[j].resind = resi_n;
186 copy_rvec((*x)[i], newx[j]);
189 copy_rvec((*v)[i], newv[j]);
194 while (i < atoms->nr && atoms->atom[i].resind == resi_o);
195 /* Increase the new residue counters */
201 /* Skip this residue */
206 while (i < atoms->nr && atoms->atom[i].resind == resi_o);
211 /* put them back into the original arrays and throw away temporary arrays */
213 *atoms_solvt = (*newatoms);
220 static void rm_res_pbc(const t_atoms *atoms, std::vector<RVec> *x, matrix box)
228 for (int n = 0; n < atoms->nr; n++)
230 if (!is_hydrogen(*(atoms->atomname[n])))
233 rvec_inc(xcg, (*x)[n]);
235 if ( (n+1 == atoms->nr) ||
236 (atoms->atom[n+1].resind != atoms->atom[n].resind) )
238 /* if nat==0 we have only hydrogens in the solvent,
239 we take last coordinate as cg */
243 copy_rvec((*x)[n], xcg);
245 svmul(1.0/nat, xcg, xcg);
246 for (int d = 0; d < DIM; d++)
250 for (int i = start; i <= n; i++)
252 (*x)[i][d] += box[d][d];
256 while (xcg[d] >= box[d][d])
258 for (int i = start; i <= n; i++)
260 (*x)[i][d] -= box[d][d];
273 * Generates a solvent configuration of desired size by stacking solvent boxes.
275 * \param[in,out] atoms Solvent atoms.
276 * \param[in,out] x Solvent positions.
277 * \param[in,out] v Solvent velocities (`*v` can be NULL).
278 * \param[in,out] r Solvent exclusion radii.
279 * \param[in] box Initial solvent box.
280 * \param[in] boxTarget Target box size.
282 * The solvent box of desired size is created by stacking the initial box in
283 * the smallest k*l*m array that covers the box, and then removing any residue
284 * where all atoms are outside the target box (with a small margin).
285 * This function does not remove overlap between solvent atoms across the
288 * Note that the input configuration should be in the rectangular unit cell and
289 * have whole residues.
291 static void replicateSolventBox(t_atoms *atoms, std::vector<RVec> *x,
292 std::vector<RVec> *v, std::vector<real> *r,
293 const matrix box, const matrix boxTarget)
295 // Calculate the box multiplication factors.
298 for (int i = 0; i < DIM; ++i)
301 while (n_box[i] * box[i][i] < boxTarget[i][i])
307 fprintf(stderr, "Will generate new solvent configuration of %dx%dx%d boxes\n",
308 n_box[XX], n_box[YY], n_box[ZZ]);
310 // Create arrays for storing the generated system (cannot be done in-place
311 // in case the target box is smaller than the original in one dimension,
314 init_t_atoms(&newAtoms, 0, FALSE);
315 gmx::AtomsBuilder builder(&newAtoms, nullptr);
316 builder.reserve(atoms->nr * nmol, atoms->nres * nmol);
317 std::vector<RVec> newX(atoms->nr * nmol);
318 std::vector<RVec> newV(!v->empty() ? atoms->nr * nmol : 0);
319 std::vector<real> newR(atoms->nr * nmol);
321 const real maxRadius = *std::max_element(r->begin(), r->end());
323 for (int i = 0; i < DIM; ++i)
325 // The code below is only interested about the box diagonal.
326 boxWithMargin[i] = boxTarget[i][i] + 3*maxRadius;
329 for (int ix = 0; ix < n_box[XX]; ++ix)
332 delta[XX] = ix*box[XX][XX];
333 for (int iy = 0; iy < n_box[YY]; ++iy)
335 delta[YY] = iy*box[YY][YY];
336 for (int iz = 0; iz < n_box[ZZ]; ++iz)
338 delta[ZZ] = iz*box[ZZ][ZZ];
339 bool bKeepResidue = false;
340 for (int i = 0; i < atoms->nr; ++i)
342 const int newIndex = builder.currentAtomCount();
343 bool bKeepAtom = true;
344 for (int m = 0; m < DIM; ++m)
346 const real newCoord = delta[m] + (*x)[i][m];
347 bKeepAtom = bKeepAtom && (newCoord < boxWithMargin[m]);
348 newX[newIndex][m] = newCoord;
350 bKeepResidue = bKeepResidue || bKeepAtom;
353 copy_rvec((*v)[i], newV[newIndex]);
355 newR[newIndex] = (*r)[i];
356 builder.addAtom(*atoms, i);
357 if (i == atoms->nr - 1
358 || atoms->atom[i+1].resind != atoms->atom[i].resind)
362 builder.finishResidue(atoms->resinfo[atoms->atom[i].resind]);
366 builder.discardCurrentResidue();
368 // Reset state for the next residue.
369 bKeepResidue = false;
376 sfree(atoms->atomname);
377 sfree(atoms->resinfo);
378 atoms->nr = newAtoms.nr;
379 atoms->nres = newAtoms.nres;
380 atoms->atom = newAtoms.atom;
381 atoms->atomname = newAtoms.atomname;
382 atoms->resinfo = newAtoms.resinfo;
384 newX.resize(atoms->nr);
388 newV.resize(atoms->nr);
391 newR.resize(atoms->nr);
394 fprintf(stderr, "Solvent box contains %d atoms in %d residues\n",
395 atoms->nr, atoms->nres);
399 * Removes overlap of solvent atoms across the edges.
401 * \param[in,out] atoms Solvent atoms.
402 * \param[in,out] x Solvent positions.
403 * \param[in,out] v Solvent velocities (can be empty).
404 * \param[in,out] r Solvent exclusion radii.
405 * \param[in] pbc PBC information.
407 * Solvent residues that lay on the edges that do not touch the origin are
408 * removed if they overlap with other solvent atoms across the PBC.
409 * This is done in this way as the assumption is that the input solvent
410 * configuration is already equilibrated, and so does not contain any
411 * undesirable overlap. The only overlap that should be removed is caused by
412 * cutting the box in half in replicateSolventBox() and leaving a margin of
413 * solvent outside those box edges; these atoms can then overlap with those on
414 * the opposite box edge in a way that is not part of the pre-equilibrated
417 static void removeSolventBoxOverlap(t_atoms *atoms, std::vector<RVec> *x,
418 std::vector<RVec> *v, std::vector<real> *r,
421 gmx::AtomsRemover remover(*atoms);
423 // TODO: We could limit the amount of pairs searched significantly,
424 // since we are only interested in pairs where the positions are on
426 const real maxRadius = *std::max_element(r->begin(), r->end());
427 gmx::AnalysisNeighborhood nb;
428 nb.setCutoff(2*maxRadius);
429 gmx::AnalysisNeighborhoodPositions pos(*x);
430 gmx::AnalysisNeighborhoodSearch search = nb.initSearch(&pbc, pos);
431 gmx::AnalysisNeighborhoodPairSearch pairSearch = search.startPairSearch(pos);
432 gmx::AnalysisNeighborhoodPair pair;
433 while (pairSearch.findNextPair(&pair))
435 const int i1 = pair.refIndex();
436 const int i2 = pair.testIndex();
437 if (remover.isMarked(i2))
439 pairSearch.skipRemainingPairsForTestPosition();
442 if (remover.isMarked(i1) || atoms->atom[i1].resind == atoms->atom[i2].resind)
446 if (pair.distance2() < gmx::square((*r)[i1] + (*r)[i2]))
449 rvec_sub((*x)[i2], (*x)[i1], dx);
450 bool bCandidate1 = false, bCandidate2 = false;
451 // To satisfy Clang static analyzer.
452 GMX_ASSERT(pbc.ndim_ePBC <= DIM, "Too many periodic dimensions");
453 for (int d = 0; d < pbc.ndim_ePBC; ++d)
455 // If the distance in some dimension is larger than the
456 // cutoff, then it means that the distance has been computed
457 // over the PBC. Mark the position with a larger coordinate
458 // for potential removal.
459 if (dx[d] > maxRadius)
463 else if (dx[d] < -maxRadius)
468 // Only mark one of the positions for removal if both were
470 if (bCandidate2 && (!bCandidate1 || i2 > i1))
472 remover.markResidue(*atoms, i2, true);
473 pairSearch.skipRemainingPairsForTestPosition();
475 else if (bCandidate1)
477 remover.markResidue(*atoms, i1, true);
482 remover.removeMarkedElements(x);
485 remover.removeMarkedElements(v);
487 remover.removeMarkedElements(r);
488 const int originalAtomCount = atoms->nr;
489 remover.removeMarkedAtoms(atoms);
490 fprintf(stderr, "Removed %d solvent atoms due to solvent-solvent overlap\n",
491 originalAtomCount - atoms->nr);
495 * Remove all solvent molecules outside a give radius from the solute.
497 * \param[in,out] atoms Solvent atoms.
498 * \param[in,out] x_solvent Solvent positions.
499 * \param[in,out] v_solvent Solvent velocities.
500 * \param[in,out] r Atomic exclusion radii.
501 * \param[in] pbc PBC information.
502 * \param[in] x_solute Solute positions.
503 * \param[in] rshell The radius outside the solute molecule.
505 static void removeSolventOutsideShell(t_atoms *atoms,
506 std::vector<RVec> *x_solvent,
507 std::vector<RVec> *v_solvent,
508 std::vector<real> *r,
510 const std::vector<RVec> &x_solute,
513 gmx::AtomsRemover remover(*atoms);
514 gmx::AnalysisNeighborhood nb;
515 nb.setCutoff(rshell);
516 gmx::AnalysisNeighborhoodPositions posSolute(x_solute);
517 gmx::AnalysisNeighborhoodSearch search = nb.initSearch(&pbc, posSolute);
518 gmx::AnalysisNeighborhoodPositions pos(*x_solvent);
519 gmx::AnalysisNeighborhoodPairSearch pairSearch = search.startPairSearch(pos);
520 gmx::AnalysisNeighborhoodPair pair;
524 // Now put back those within the shell without checking for overlap
525 while (pairSearch.findNextPair(&pair))
527 remover.markResidue(*atoms, pair.testIndex(), false);
528 pairSearch.skipRemainingPairsForTestPosition();
530 remover.removeMarkedElements(x_solvent);
531 if (!v_solvent->empty())
533 remover.removeMarkedElements(v_solvent);
535 remover.removeMarkedElements(r);
536 const int originalAtomCount = atoms->nr;
537 remover.removeMarkedAtoms(atoms);
538 fprintf(stderr, "Removed %d solvent atoms more than %f nm from solute.\n",
539 originalAtomCount - atoms->nr, rshell);
543 * Removes solvent molecules that overlap with the solute.
545 * \param[in,out] atoms Solvent atoms.
546 * \param[in,out] x Solvent positions.
547 * \param[in,out] v Solvent velocities (can be empty).
548 * \param[in,out] r Solvent exclusion radii.
549 * \param[in] pbc PBC information.
550 * \param[in] x_solute Solute positions.
551 * \param[in] r_solute Solute exclusion radii.
553 static void removeSolventOverlappingWithSolute(t_atoms *atoms,
554 std::vector<RVec> *x,
555 std::vector<RVec> *v,
556 std::vector<real> *r,
558 const std::vector<RVec> &x_solute,
559 const std::vector<real> &r_solute)
561 gmx::AtomsRemover remover(*atoms);
562 const real maxRadius1
563 = *std::max_element(r->begin(), r->end());
564 const real maxRadius2
565 = *std::max_element(r_solute.begin(), r_solute.end());
567 // Now check for overlap.
568 gmx::AnalysisNeighborhood nb;
569 gmx::AnalysisNeighborhoodPair pair;
570 nb.setCutoff(maxRadius1 + maxRadius2);
571 gmx::AnalysisNeighborhoodPositions posSolute(x_solute);
572 gmx::AnalysisNeighborhoodSearch search = nb.initSearch(&pbc, posSolute);
573 gmx::AnalysisNeighborhoodPositions pos(*x);
574 gmx::AnalysisNeighborhoodPairSearch pairSearch = search.startPairSearch(pos);
575 while (pairSearch.findNextPair(&pair))
577 if (remover.isMarked(pair.testIndex()))
579 pairSearch.skipRemainingPairsForTestPosition();
582 const real r1 = r_solute[pair.refIndex()];
583 const real r2 = (*r)[pair.testIndex()];
584 const bool bRemove = (pair.distance2() < gmx::square(r1 + r2));
585 remover.markResidue(*atoms, pair.testIndex(), bRemove);
588 remover.removeMarkedElements(x);
591 remover.removeMarkedElements(v);
593 remover.removeMarkedElements(r);
594 const int originalAtomCount = atoms->nr;
595 remover.removeMarkedAtoms(atoms);
596 fprintf(stderr, "Removed %d solvent atoms due to solute-solvent overlap\n",
597 originalAtomCount - atoms->nr);
601 * Removes a given number of solvent residues.
603 * \param[in,out] atoms Solvent atoms.
604 * \param[in,out] x Solvent positions.
605 * \param[in,out] v Solvent velocities (can be empty).
606 * \param[in] numberToRemove Number of residues to remove.
608 * This function is called last in the process of creating the solvent box,
609 * so it does not operate on the exclusion radii, as no code after this needs
612 static void removeExtraSolventMolecules(t_atoms *atoms, std::vector<RVec> *x,
613 std::vector<RVec> *v,
616 gmx::AtomsRemover remover(*atoms);
617 std::random_device rd;
618 std::mt19937 randomNumberGenerator(rd());
619 std::uniform_int_distribution<> randomDistribution(0, atoms->nr - 1);
620 while (numberToRemove > 0)
622 int atomIndex = randomDistribution(randomNumberGenerator);
623 if (!remover.isMarked(atomIndex))
625 remover.markResidue(*atoms, atomIndex, true);
629 remover.removeMarkedElements(x);
632 remover.removeMarkedElements(v);
634 remover.removeMarkedAtoms(atoms);
637 static void add_solv(const char *filename, t_atoms *atoms,
639 std::vector<RVec> *x, std::vector<RVec> *v,
640 int ePBC, matrix box, AtomProperties *aps,
641 real defaultDistance, real scaleFactor,
642 real rshell, int max_sol)
644 gmx_mtop_t topSolvent;
645 std::vector<RVec> xSolvent, vSolvent;
646 matrix boxSolvent = {{ 0 }};
649 fprintf(stderr, "Reading solvent configuration\n");
650 bool bTprFileWasRead;
651 rvec *temporaryX = nullptr, *temporaryV = nullptr;
652 readConfAndTopology(gmx::findLibraryFile(filename).c_str(), &bTprFileWasRead, &topSolvent,
653 &ePBCSolvent, &temporaryX, &temporaryV, boxSolvent);
654 t_atoms *atomsSolvent;
655 snew(atomsSolvent, 1);
656 *atomsSolvent = gmx_mtop_global_atoms(&topSolvent);
657 xSolvent.assign(temporaryX, temporaryX + topSolvent.natoms);
659 vSolvent.assign(temporaryV, temporaryV + topSolvent.natoms);
661 if (gmx::boxIsZero(boxSolvent))
663 gmx_fatal(FARGS, "No box information for solvent in %s, please use a properly formatted file\n",
666 if (0 == atomsSolvent->nr)
668 gmx_fatal(FARGS, "No solvent in %s, please check your input\n", filename);
670 fprintf(stderr, "\n");
672 /* initialise distance arrays for solvent configuration */
673 fprintf(stderr, "Initialising inter-atomic distances...\n");
674 const std::vector<real> exclusionDistances(
675 makeExclusionDistances(atoms, aps, defaultDistance, scaleFactor));
676 std::vector<real> exclusionDistances_solvt(
677 makeExclusionDistances(atomsSolvent, aps, defaultDistance, scaleFactor));
679 /* generate a new solvent configuration */
680 fprintf(stderr, "Generating solvent configuration\n");
682 set_pbc(&pbc, ePBC, box);
683 if (!gmx::boxesAreEqual(boxSolvent, box))
685 if (TRICLINIC(boxSolvent))
687 gmx_fatal(FARGS, "Generating from non-rectangular solvent boxes is currently not supported.\n"
688 "You can try to pass the same box for -cp and -cs.");
690 /* apply pbc for solvent configuration for whole molecules */
691 rm_res_pbc(atomsSolvent, &xSolvent, boxSolvent);
692 replicateSolventBox(atomsSolvent, &xSolvent, &vSolvent, &exclusionDistances_solvt,
694 if (ePBC != epbcNONE)
696 removeSolventBoxOverlap(atomsSolvent, &xSolvent, &vSolvent, &exclusionDistances_solvt, pbc);
703 removeSolventOutsideShell(atomsSolvent, &xSolvent, &vSolvent,
704 &exclusionDistances_solvt, pbc, *x, rshell);
706 removeSolventOverlappingWithSolute(atomsSolvent, &xSolvent, &vSolvent,
707 &exclusionDistances_solvt, pbc, *x,
711 if (max_sol > 0 && atomsSolvent->nres > max_sol)
713 const int numberToRemove = atomsSolvent->nres - max_sol;
714 removeExtraSolventMolecules(atomsSolvent, &xSolvent, &vSolvent, numberToRemove);
717 /* Sort the solvent mixture, not the protein... */
718 t_atoms *newatoms = nullptr;
719 // The sort_molecule function does something creative with the
720 // t_atoms pointers. We need to make sure we neither leak, nor
721 // double-free, so make a shallow pointer that is fine for it to
723 t_atoms *sortedAtomsSolvent = atomsSolvent;
724 sort_molecule(&sortedAtomsSolvent, &newatoms, &xSolvent, &vSolvent);
726 // Merge the two configurations.
727 x->insert(x->end(), xSolvent.begin(), xSolvent.end());
730 v->insert(v->end(), vSolvent.begin(), vSolvent.end());
733 gmx::AtomsBuilder builder(atoms, symtab);
734 builder.mergeAtoms(*sortedAtomsSolvent);
736 fprintf(stderr, "Generated solvent containing %d atoms in %d residues\n",
737 atomsSolvent->nr, atomsSolvent->nres);
746 done_atom(atomsSolvent);
751 static void update_top(t_atoms *atoms,
752 int firstSolventResidueIndex,
759 char buf[STRLEN*2], buf2[STRLEN], *temp;
760 const char *topinout;
767 int nsol = atoms->nres - firstSolventResidueIndex;
770 for (i = 0; (i < atoms->nr); i++)
772 aps->setAtomProperty(epropMass,
773 std::string(*atoms->resinfo[atoms->atom[i].resind].name),
774 std::string(*atoms->atomname[i]), &mm);
780 fprintf(stderr, "Volume : %10g (nm^3)\n", vol);
781 fprintf(stderr, "Density : %10g (g/l)\n",
782 (mtot*1e24)/(AVOGADRO*vol));
783 fprintf(stderr, "Number of solvent molecules: %5d \n\n", nsol);
785 /* open topology file and append sol molecules */
786 topinout = ftp2fn(efTOP, NFILE, fnm);
787 if (ftp2bSet(efTOP, NFILE, fnm) )
789 char temporary_filename[STRLEN];
790 strncpy(temporary_filename, "temp.topXXXXXX", STRLEN);
792 fprintf(stderr, "Processing topology\n");
793 fpin = gmx_ffopen(topinout, "r");
794 fpout = gmx_fopen_temporary(temporary_filename);
797 while (fgets(buf, STRLEN, fpin))
801 if ((temp = strchr(buf2, '\n')) != nullptr)
809 if ((temp = strchr(buf2, '\n')) != nullptr)
814 if (buf2[strlen(buf2)-1] == ']')
816 buf2[strlen(buf2)-1] = '\0';
819 bSystem = (gmx_strcasecmp(buf2, "system") == 0);
822 else if (bSystem && nsol && (buf[0] != ';') )
824 /* if sol present, append "in water" to system name */
826 if (buf2[0] && (!strstr(buf2, " water")) )
828 sprintf(buf, "%s in water\n", buf2);
832 fprintf(fpout, "%s", buf);
836 // Add new solvent molecules to the topology
839 std::string currRes = *atoms->resinfo[firstSolventResidueIndex].name;
842 // Iterate through solvent molecules and increment a count until new resname found
843 for (int i = firstSolventResidueIndex; i < atoms->nres; i++)
845 if ((currRes == *atoms->resinfo[i].name))
851 // Change topology and restart count
852 fprintf(stdout, "Adding line for %d solvent molecules with resname (%s) to "
853 "topology file (%s)\n", resCount, currRes.c_str(), topinout);
854 fprintf(fpout, "%-15s %5d\n", currRes.c_str(), resCount);
855 currRes = *atoms->resinfo[i].name;
859 // One more print needed for last residue type
860 fprintf(stdout, "Adding line for %d solvent molecules with resname (%s) to "
861 "topology file (%s)\n", resCount, currRes.c_str(), topinout);
862 fprintf(fpout, "%-15s %5d\n", currRes.c_str(), resCount);
865 make_backup(topinout);
866 gmx_file_rename(temporary_filename, topinout);
870 int gmx_solvate(int argc, char *argv[])
872 const char *desc[] = {
873 "[THISMODULE] can do one of 2 things:[PAR]",
875 "1) Generate a box of solvent. Specify [TT]-cs[tt] and [TT]-box[tt].",
876 "Or specify [TT]-cs[tt] and [TT]-cp[tt] with a structure file with",
877 "a box, but without atoms.[PAR]",
879 "2) Solvate a solute configuration, e.g. a protein, in a bath of solvent ",
880 "molecules. Specify [TT]-cp[tt] (solute) and [TT]-cs[tt] (solvent). ",
881 "The box specified in the solute coordinate file ([TT]-cp[tt]) is used,",
882 "unless [TT]-box[tt] is set.",
883 "If you want the solute to be centered in the box,",
884 "the program [gmx-editconf] has sophisticated options",
885 "to change the box dimensions and center the solute.",
886 "Solvent molecules are removed from the box where the ",
887 "distance between any atom of the solute molecule(s) and any atom of ",
888 "the solvent molecule is less than the sum of the scaled van der Waals",
889 "radii of both atoms. A database ([TT]vdwradii.dat[tt]) of van der",
890 "Waals radii is read by the program, and the resulting radii scaled",
891 "by [TT]-scale[tt]. If radii are not found in the database, those",
892 "atoms are assigned the (pre-scaled) distance [TT]-radius[tt].",
893 "Note that the usefulness of those radii depends on the atom names,",
894 "and thus varies widely with force field.",
896 "The default solvent is Simple Point Charge water (SPC), with coordinates ",
897 "from [TT]$GMXLIB/spc216.gro[tt]. These coordinates can also be used",
898 "for other 3-site water models, since a short equibilibration will remove",
899 "the small differences between the models.",
900 "Other solvents are also supported, as well as mixed solvents. The",
901 "only restriction to solvent types is that a solvent molecule consists",
902 "of exactly one residue. The residue information in the coordinate",
903 "files is used, and should therefore be more or less consistent.",
904 "In practice this means that two subsequent solvent molecules in the ",
905 "solvent coordinate file should have different residue number.",
906 "The box of solute is built by stacking the coordinates read from",
907 "the coordinate file. This means that these coordinates should be ",
908 "equlibrated in periodic boundary conditions to ensure a good",
909 "alignment of molecules on the stacking interfaces.",
910 "The [TT]-maxsol[tt] option simply adds only the first [TT]-maxsol[tt]",
911 "solvent molecules and leaves out the rest that would have fitted",
912 "into the box. This can create a void that can cause problems later.",
913 "Choose your volume wisely.[PAR]",
915 "Setting [TT]-shell[tt] larger than zero will place a layer of water of",
916 "the specified thickness (nm) around the solute. Hint: it is a good",
917 "idea to put the protein in the center of a box first (using [gmx-editconf]).",
920 "Finally, [THISMODULE] will optionally remove lines from your topology file in ",
921 "which a number of solvent molecules is already added, and adds a ",
922 "line with the total number of solvent molecules in your coordinate file."
925 const char *bugs[] = {
926 "Molecules must be whole in the initial configurations.",
930 gmx_bool bProt, bBox;
931 const char *conf_prot, *confout;
934 { efSTX, "-cp", "protein", ffOPTRD },
935 { efSTX, "-cs", "spc216", ffLIBRD},
936 { efSTO, nullptr, nullptr, ffWRITE},
937 { efTOP, nullptr, nullptr, ffOPTRW},
939 #define NFILE asize(fnm)
941 real defaultDistance = 0.105, r_shell = 0, scaleFactor = 0.57;
942 rvec new_box = {0.0, 0.0, 0.0};
943 gmx_bool bReadV = FALSE;
945 int firstSolventResidueIndex = 0;
946 gmx_output_env_t *oenv;
948 { "-box", FALSE, etRVEC, {new_box},
949 "Box size (in nm)" },
950 { "-radius", FALSE, etREAL, {&defaultDistance},
951 "Default van der Waals distance"},
952 { "-scale", FALSE, etREAL, {&scaleFactor},
953 "Scale factor to multiply Van der Waals radii from the database in share/gromacs/top/vdwradii.dat. The default value of 0.57 yields density close to 1000 g/l for proteins in water." },
954 { "-shell", FALSE, etREAL, {&r_shell},
955 "Thickness of optional water layer around solute" },
956 { "-maxsol", FALSE, etINT, {&max_sol},
957 "Maximum number of solvent molecules to add if they fit in the box. If zero (default) this is ignored" },
958 { "-vel", FALSE, etBOOL, {&bReadV},
959 "Keep velocities from input solute and solvent" },
962 if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa,
963 asize(desc), desc, asize(bugs), bugs, &oenv))
968 const char *solventFileName = opt2fn("-cs", NFILE, fnm);
969 bProt = opt2bSet("-cp", NFILE, fnm);
970 bBox = opt2parg_bSet("-box", asize(pa), pa);
975 gmx_fatal(FARGS, "When no solute (-cp) is specified, "
976 "a box size (-box) must be specified");
981 /* solute configuration data */
983 std::vector<RVec> x, v;
984 matrix box = {{ 0 }};
990 /* Generate a solute configuration */
991 conf_prot = opt2fn("-cp", NFILE, fnm);
992 fprintf(stderr, "Reading solute configuration%s\n",
993 bReadV ? " and velocities" : "");
994 bool bTprFileWasRead;
995 rvec *temporaryX = nullptr, *temporaryV = nullptr;
996 readConfAndTopology(conf_prot, &bTprFileWasRead, &top,
997 &ePBC, &temporaryX, bReadV ? &temporaryV : nullptr, box);
998 *atoms = gmx_mtop_global_atoms(&top);
999 x.assign(temporaryX, temporaryX + top.natoms);
1003 v.assign(temporaryV, temporaryV + top.natoms);
1008 fprintf(stderr, "Note: no velocities found\n");
1012 fprintf(stderr, "Note: no atoms in %s\n", conf_prot);
1017 firstSolventResidueIndex = atoms->nres;
1020 int ePBCForOutput = ePBC;
1023 ePBCForOutput = epbcXYZ;
1025 box[XX][XX] = new_box[XX];
1026 box[YY][YY] = new_box[YY];
1027 box[ZZ][ZZ] = new_box[ZZ];
1031 gmx_fatal(FARGS, "Undefined solute box.\nCreate one with gmx editconf "
1032 "or give explicit -box command line option");
1035 add_solv(solventFileName, atoms, &top.symtab, &x, &v, ePBCForOutput, box,
1036 &aps, defaultDistance, scaleFactor, r_shell, max_sol);
1038 /* write new configuration 1 to file confout */
1039 confout = ftp2fn(efSTO, NFILE, fnm);
1040 fprintf(stderr, "Writing generated configuration to %s\n", confout);
1041 const char *outputTitle = (bProt ? *top.name : "Generated by gmx solvate");
1042 write_sto_conf(confout, outputTitle, atoms, as_rvec_array(x.data()),
1043 !v.empty() ? as_rvec_array(v.data()) : nullptr, ePBCForOutput, box);
1045 /* print size of generated configuration */
1046 fprintf(stderr, "\nOutput configuration contains %d atoms in %d residues\n",
1047 atoms->nr, atoms->nres);
1048 update_top(atoms, firstSolventResidueIndex, box, NFILE, fnm, &aps);
1052 output_env_done(oenv);