2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2010-2018, The GROMACS development team.
5 * Copyright (c) 2019,2020, by the GROMACS development team, led by
6 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7 * and including many others, as listed in the AUTHORS file in the
8 * top-level source directory and at http://www.gromacs.org.
10 * GROMACS is free software; you can redistribute it and/or
11 * modify it under the terms of the GNU Lesser General Public License
12 * as published by the Free Software Foundation; either version 2.1
13 * of the License, or (at your option) any later version.
15 * GROMACS is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 * Lesser General Public License for more details.
20 * You should have received a copy of the GNU Lesser General Public
21 * License along with GROMACS; if not, see
22 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
25 * If you want to redistribute modifications to GROMACS, please
26 * consider that scientific software is very special. Version
27 * control is crucial - bugs must be traceable. We will be happy to
28 * consider code for inclusion in the official distribution, but
29 * derived work must not be called official GROMACS. Details are found
30 * in the README & COPYING files - if they are missing, get the
31 * official version at http://www.gromacs.org.
33 * To help us fund GROMACS development, we humbly ask that you cite
34 * the research papers on the package. Check out http://www.gromacs.org.
43 #include "gromacs/commandline/filenm.h"
44 #include "gromacs/fileio/readinp.h"
45 #include "gromacs/fileio/warninp.h"
46 #include "gromacs/gmxlib/network.h"
47 #include "gromacs/math/vec.h"
48 #include "gromacs/mdtypes/commrec.h"
49 #include "gromacs/mdtypes/inputrec.h"
50 #include "gromacs/mdtypes/md_enums.h"
51 #include "gromacs/mdtypes/state.h"
52 #include "gromacs/pbcutil/pbc.h"
53 #include "gromacs/topology/index.h"
54 #include "gromacs/topology/mtop_lookup.h"
55 #include "gromacs/topology/mtop_util.h"
56 #include "gromacs/topology/topology.h"
57 #include "gromacs/utility/cstringutil.h"
58 #include "gromacs/utility/exceptions.h"
59 #include "gromacs/utility/fatalerror.h"
60 #include "gromacs/utility/filestream.h"
61 #include "gromacs/utility/futil.h"
62 #include "gromacs/utility/smalloc.h"
64 /* information about scaling center */
67 rvec xmin; /* smallest coordinates of all embedded molecules */
68 rvec xmax; /* largest coordinates of all embedded molecules */
69 rvec* geom_cent; /* scaling center of each independent molecule to embed */
70 int pieces; /* number of molecules to embed independently */
71 int* nidx; /* n atoms for every independent embedded molecule (index in subindex) */
72 int** subindex; /* atomids for independent molecule *
73 * atoms of piece i run from subindex[i][0] to subindex[i][nidx[i]] */
76 /* variables needed in do_md */
79 int it_xy; /* number of iterations (steps) used to grow something in the xy-plane */
80 int it_z; /* same, but for z */
81 real xy_step; /* stepsize used during resize in xy-plane */
82 real z_step; /* same, but in z */
83 rvec fac; /* initial scaling of the molecule to grow into the membrane */
84 rvec* r_ins; /* final coordinates of the molecule to grow */
85 pos_ins_t* pos_ins; /* scaling center for each piece to embed */
88 /* membrane related variables */
91 char* name; /* name of index group to embed molecule into (usually membrane) */
92 t_block mem_at; /* list all atoms in membrane */
93 int nmol; /* number of membrane molecules overlapping with the molecule to embed */
94 int* mol_id; /* list of molecules in membrane that overlap with the molecule to embed */
95 real lip_area; /* average area per lipid in membrane (only correct for homogeneous bilayers)*/
96 real zmin; /* minimum z coordinate of membrane */
97 real zmax; /* maximum z coordinate of membrane */
98 real zmed; /* median z coordinate of membrane */
101 /* Lists all molecules in the membrane that overlap with the molecule to be embedded. *
102 * These will then be removed from the system */
105 int nr; /* number of molecules to remove */
106 int* mol; /* list of molecule ids to remove */
107 int* block; /* id of the molblock that the molecule to remove is part of */
110 /* Get the global molecule id, and the corresponding molecule type and id of the *
111 * molblock from the global atom nr. */
112 static int get_mol_id(int at, gmx_mtop_t* mtop, int* type, int* block)
119 mtopGetMolblockIndex(mtop, at, block, &mol_id, &atnr_mol);
120 for (i = 0; i < *block; i++)
122 mol_id += mtop->molblock[i].nmol;
124 *type = mtop->molblock[*block].type;
129 /* Get the id of the molblock from a global molecule id */
130 static int get_molblock(int mol_id, const std::vector<gmx_molblock_t>& mblock)
134 for (size_t i = 0; i < mblock.size(); i++)
136 nmol += mblock[i].nmol;
143 gmx_fatal(FARGS, "mol_id %d larger than total number of molecules %d.\n", mol_id, nmol);
146 /* Get a list of all the molecule types that are present in a group of atoms. *
147 * Because all interaction within the group to embed are removed on the topology *
148 * level, if the same molecule type is found in another part of the system, these *
149 * would also be affected. Therefore we have to check if the embedded and rest group *
150 * share common molecule types. If so, membed will stop with an error. */
151 static int get_mtype_list(t_block* at, gmx_mtop_t* mtop, t_block* tlist)
154 int type = 0, block = 0;
158 snew(tlist->index, at->nr);
159 for (i = 0; i < at->nr; i++)
162 get_mol_id(at->index[i], mtop, &type, &block);
163 for (j = 0; j < nr; j++)
165 if (tlist->index[j] == type)
173 tlist->index[nr] = type;
177 srenew(tlist->index, nr);
181 /* Do the actual check of the molecule types between embedded and rest group */
182 static void check_types(t_block* ins_at, t_block* rest_at, gmx_mtop_t* mtop)
184 t_block *ins_mtype, *rest_mtype;
189 ins_mtype->nr = get_mtype_list(ins_at, mtop, ins_mtype);
190 rest_mtype->nr = get_mtype_list(rest_at, mtop, rest_mtype);
192 for (i = 0; i < ins_mtype->nr; i++)
194 for (j = 0; j < rest_mtype->nr; j++)
196 if (ins_mtype->index[i] == rest_mtype->index[j])
200 "Moleculetype %s is found both in the group to insert and the rest of the "
202 "1. Your *.ndx and *.top do not match\n"
203 "2. You are inserting some molecules of type %s (for example "
204 "xray-solvent), while\n"
205 "the same moleculetype is also used in the rest of the system (solvent "
207 "we need to exclude all interactions between the atoms in the group to\n"
208 "insert, the same moleculetype can not be used in both groups. Change the\n"
209 "moleculetype of the molecules %s in the inserted group. Do not forget to "
211 "an appropriate *.itp file",
212 *(mtop->moltype[rest_mtype->index[j]].name),
213 *(mtop->moltype[rest_mtype->index[j]].name),
214 *(mtop->moltype[rest_mtype->index[j]].name));
219 done_block(ins_mtype);
220 done_block(rest_mtype);
225 static void get_input(const char* membed_input,
236 gmx_bool* bALLOW_ASYMMETRY)
239 std::vector<t_inpfile> inp;
241 wi = init_warning(TRUE, 0);
244 gmx::TextInputFile stream(membed_input);
245 inp = read_inpfile(&stream, membed_input, wi);
248 *it_xy = get_eint(&inp, "nxy", 1000, wi);
249 *it_z = get_eint(&inp, "nz", 0, wi);
250 *xy_fac = get_ereal(&inp, "xyinit", 0.5, wi);
251 *xy_max = get_ereal(&inp, "xyend", 1.0, wi);
252 *z_fac = get_ereal(&inp, "zinit", 1.0, wi);
253 *z_max = get_ereal(&inp, "zend", 1.0, wi);
254 *probe_rad = get_ereal(&inp, "rad", 0.22, wi);
255 *low_up_rm = get_eint(&inp, "ndiff", 0, wi);
256 *maxwarn = get_eint(&inp, "maxwarn", 0, wi);
257 *pieces = get_eint(&inp, "pieces", 1, wi);
258 const char* yesno_names[BOOL_NR + 1] = { "no", "yes", nullptr };
259 *bALLOW_ASYMMETRY = (get_eeenum(&inp, "asymmetry", yesno_names, wi) != 0);
261 check_warning_error(wi, FARGS);
263 gmx::TextOutputFile stream(membed_input);
264 write_inpfile(&stream, membed_input, &inp, FALSE, WriteMdpHeader::yes, wi);
267 done_warning(wi, FARGS);
270 /* Obtain the maximum and minimum coordinates of the group to be embedded */
271 static int init_ins_at(t_block* ins_at,
275 SimulationGroups* groups,
280 real xmin, xmax, ymin, ymax, zmin, zmax;
281 const real min_memthick = 6.0; /* minimum thickness of the bilayer that will be used to *
282 * determine the overlap between molecule to embed and membrane */
283 const real fac_inp_size =
284 1.000001; /* scaling factor to obtain input_size + 0.000001 (comparing reals) */
285 snew(rest_at->index, state->natoms);
286 auto x = makeArrayRef(state->x);
288 xmin = xmax = x[ins_at->index[0]][XX];
289 ymin = ymax = x[ins_at->index[0]][YY];
290 zmin = zmax = x[ins_at->index[0]][ZZ];
292 for (i = 0; i < state->natoms; i++)
294 gid = groups->groupNumbers[SimulationAtomGroupType::Freeze][i];
295 if (groups->groups[SimulationAtomGroupType::Freeze][gid] == ins_grp_id)
297 xmin = std::min(xmin, x[i][XX]);
298 xmax = std::max(xmax, x[i][XX]);
299 ymin = std::min(ymin, x[i][YY]);
300 ymax = std::max(ymax, x[i][YY]);
301 zmin = std::min(zmin, x[i][ZZ]);
302 zmax = std::max(zmax, x[i][ZZ]);
306 rest_at->index[c] = i;
312 srenew(rest_at->index, c);
314 if (xy_max > fac_inp_size)
316 pos_ins->xmin[XX] = xmin - ((xmax - xmin) * xy_max - (xmax - xmin)) / 2;
317 pos_ins->xmin[YY] = ymin - ((ymax - ymin) * xy_max - (ymax - ymin)) / 2;
319 pos_ins->xmax[XX] = xmax + ((xmax - xmin) * xy_max - (xmax - xmin)) / 2;
320 pos_ins->xmax[YY] = ymax + ((ymax - ymin) * xy_max - (ymax - ymin)) / 2;
324 pos_ins->xmin[XX] = xmin;
325 pos_ins->xmin[YY] = ymin;
327 pos_ins->xmax[XX] = xmax;
328 pos_ins->xmax[YY] = ymax;
331 if ((zmax - zmin) < min_memthick)
333 pos_ins->xmin[ZZ] = zmin + (zmax - zmin) / 2.0 - 0.5 * min_memthick;
334 pos_ins->xmax[ZZ] = zmin + (zmax - zmin) / 2.0 + 0.5 * min_memthick;
338 pos_ins->xmin[ZZ] = zmin;
339 pos_ins->xmax[ZZ] = zmax;
345 /* Estimate the area of the embedded molecule by projecting all coordinates on a grid in the *
346 * xy-plane and counting the number of occupied grid points */
347 static real est_prot_area(pos_ins_t* pos_ins, rvec* r, t_block* ins_at, mem_t* mem_p)
349 real x, y, dx = 0.15, dy = 0.15, area = 0.0;
350 real add, memmin, memmax;
353 /* min and max membrane coordinate are altered to reduce the influence of the *
355 memmin = mem_p->zmin + 0.1 * (mem_p->zmax - mem_p->zmin);
356 memmax = mem_p->zmax - 0.1 * (mem_p->zmax - mem_p->zmin);
358 //NOLINTNEXTLINE(clang-analyzer-security.FloatLoopCounter)
359 for (x = pos_ins->xmin[XX]; x < pos_ins->xmax[XX]; x += dx)
361 //NOLINTNEXTLINE(clang-analyzer-security.FloatLoopCounter)
362 for (y = pos_ins->xmin[YY]; y < pos_ins->xmax[YY]; y += dy)
368 at = ins_at->index[c];
369 if ((r[at][XX] >= x) && (r[at][XX] < x + dx) && (r[at][YY] >= y)
370 && (r[at][YY] < y + dy) && (r[at][ZZ] > memmin) && (r[at][ZZ] < memmax))
375 } while ((c < ins_at->nr) && (add < 0.5));
379 area = area * dx * dy;
384 static int init_mem_at(mem_t* mem_p, gmx_mtop_t* mtop, rvec* r, matrix box, pos_ins_t* pos_ins)
386 int i, j, at, mol, nmol, nmolbox, count;
388 real z, zmin, zmax, mem_area;
391 int type = 0, block = 0;
394 mem_a = &(mem_p->mem_at);
395 snew(mol_id, mem_a->nr);
396 zmin = pos_ins->xmax[ZZ];
397 zmax = pos_ins->xmin[ZZ];
398 for (i = 0; i < mem_a->nr; i++)
400 at = mem_a->index[i];
401 if ((r[at][XX] > pos_ins->xmin[XX]) && (r[at][XX] < pos_ins->xmax[XX])
402 && (r[at][YY] > pos_ins->xmin[YY]) && (r[at][YY] < pos_ins->xmax[YY])
403 && (r[at][ZZ] > pos_ins->xmin[ZZ]) && (r[at][ZZ] < pos_ins->xmax[ZZ]))
405 mol = get_mol_id(at, mtop, &type, &block);
407 for (j = 0; j < nmol; j++)
409 if (mol == mol_id[j])
437 srenew(mol_id, nmol);
438 mem_p->mol_id = mol_id;
440 if ((zmax - zmin) > (box[ZZ][ZZ] - 0.5))
443 "Something is wrong with your membrane. Max and min z values are %f and %f.\n"
444 "Maybe your membrane is not centered in the box, but located at the box edge in "
446 "so that one membrane is distributed over two periodic box images. Another "
447 "possibility is that\n"
448 "your water layer is not thick enough.\n",
454 mem_p->zmed = (zmax - zmin) / 2 + zmin;
456 /*number of membrane molecules in protein box*/
457 nmolbox = count / mtop->moltype[mtop->molblock[block].type].atoms.nr;
458 /*membrane area within the box defined by the min and max coordinates of the embedded molecule*/
459 mem_area = (pos_ins->xmax[XX] - pos_ins->xmin[XX]) * (pos_ins->xmax[YY] - pos_ins->xmin[YY]);
460 /*rough estimate of area per lipid, assuming there is only one type of lipid in the membrane*/
461 mem_p->lip_area = 2.0 * mem_area / static_cast<double>(nmolbox);
463 return mem_p->mem_at.nr;
466 static void init_resize(t_block* ins_at, rvec* r_ins, pos_ins_t* pos_ins, mem_t* mem_p, rvec* r, gmx_bool bALLOW_ASYMMETRY)
468 int i, j, at, c, outsidesum, gctr = 0;
472 for (i = 0; i < pos_ins->pieces; i++)
474 idxsum += pos_ins->nidx[i];
477 if (idxsum != ins_at->nr)
480 "Piecewise sum of inserted atoms not same as size of group selected to insert.");
483 snew(pos_ins->geom_cent, pos_ins->pieces);
484 for (i = 0; i < pos_ins->pieces; i++)
488 for (j = 0; j < DIM; j++)
490 pos_ins->geom_cent[i][j] = 0;
493 for (j = 0; j < pos_ins->nidx[i]; j++)
495 at = pos_ins->subindex[i][j];
496 copy_rvec(r[at], r_ins[gctr]);
497 if ((r_ins[gctr][ZZ] < mem_p->zmax) && (r_ins[gctr][ZZ] > mem_p->zmin))
499 rvec_inc(pos_ins->geom_cent[i], r_ins[gctr]);
511 svmul(1 / static_cast<double>(c), pos_ins->geom_cent[i], pos_ins->geom_cent[i]);
514 if (!bALLOW_ASYMMETRY)
516 pos_ins->geom_cent[i][ZZ] = mem_p->zmed;
520 "Embedding piece %d with center of geometry: %f %f %f\n",
522 pos_ins->geom_cent[i][XX],
523 pos_ins->geom_cent[i][YY],
524 pos_ins->geom_cent[i][ZZ]);
526 fprintf(stderr, "\n");
529 /* resize performed in the md loop */
530 static void resize(rvec* r_ins, rvec* r, pos_ins_t* pos_ins, const rvec fac)
532 int i, j, k, at, c = 0;
533 for (k = 0; k < pos_ins->pieces; k++)
535 for (i = 0; i < pos_ins->nidx[k]; i++)
537 at = pos_ins->subindex[k][i];
538 for (j = 0; j < DIM; j++)
540 r[at][j] = pos_ins->geom_cent[k][j] + fac[j] * (r_ins[c][j] - pos_ins->geom_cent[k][j]);
547 /* generate the list of membrane molecules that overlap with the molecule to be embedded. *
548 * The molecule to be embedded is already reduced in size. */
549 static int gen_rm_list(rm_t* rm_p,
559 gmx_bool bALLOW_ASYMMETRY)
561 int i, j, k, l, at, at2, mol_id;
562 int type = 0, block = 0;
563 int nrm, nupper, nlower;
564 real r_min_rad, z_lip, min_norm;
570 r_min_rad = probe_rad * probe_rad;
571 gmx::RangePartitioning molecules = gmx_mtop_molecules(*mtop);
572 snew(rm_p->block, molecules.numBlocks());
573 snew(rm_p->mol, molecules.numBlocks());
576 for (i = 0; i < ins_at->nr; i++)
578 at = ins_at->index[i];
579 for (j = 0; j < rest_at->nr; j++)
581 at2 = rest_at->index[j];
582 pbc_dx(pbc, r[at], r[at2], dr);
584 if (norm2(dr) < r_min_rad)
586 mol_id = get_mol_id(at2, mtop, &type, &block);
588 for (l = 0; l < nrm; l++)
590 if (rm_p->mol[l] == mol_id)
598 rm_p->mol[nrm] = mol_id;
599 rm_p->block[nrm] = block;
602 for (l = 0; l < mem_p->nmol; l++)
604 if (mol_id == mem_p->mol_id[l])
606 for (int k : molecules.block(mol_id))
610 z_lip /= molecules.block(mol_id).size();
611 if (z_lip < mem_p->zmed)
626 /*make sure equal number of lipids from upper and lower layer are removed */
627 if ((nupper != nlower) && (!bALLOW_ASYMMETRY))
629 snew(dist, mem_p->nmol);
630 snew(order, mem_p->nmol);
631 for (i = 0; i < mem_p->nmol; i++)
633 at = molecules.block(mem_p->mol_id[i]).begin();
634 pbc_dx(pbc, r[at], pos_ins->geom_cent[0], dr);
635 if (pos_ins->pieces > 1)
638 min_norm = norm2(dr);
639 for (k = 1; k < pos_ins->pieces; k++)
641 pbc_dx(pbc, r[at], pos_ins->geom_cent[k], dr_tmp);
642 if (norm2(dr_tmp) < min_norm)
644 min_norm = norm2(dr_tmp);
645 copy_rvec(dr_tmp, dr);
649 dist[i] = dr[XX] * dr[XX] + dr[YY] * dr[YY];
651 while (j >= 0 && dist[i] < dist[order[j]])
653 order[j + 1] = order[j];
660 while (nupper != nlower)
662 mol_id = mem_p->mol_id[order[i]];
663 block = get_molblock(mol_id, mtop->molblock);
665 for (l = 0; l < nrm; l++)
667 if (rm_p->mol[l] == mol_id)
676 for (int k : molecules.block(mol_id))
680 z_lip /= molecules.block(mol_id).size();
681 if (nupper > nlower && z_lip < mem_p->zmed)
683 rm_p->mol[nrm] = mol_id;
684 rm_p->block[nrm] = block;
688 else if (nupper < nlower && z_lip > mem_p->zmed)
690 rm_p->mol[nrm] = mol_id;
691 rm_p->block[nrm] = block;
701 "Trying to remove more lipid molecules than there are in the membrane");
709 srenew(rm_p->mol, nrm);
710 srenew(rm_p->block, nrm);
712 return nupper + nlower;
715 /*remove all lipids and waters overlapping and update all important structures (e.g. state and mtop)*/
716 static void rm_group(SimulationGroups* groups,
723 int j, k, n, rm, mol_id, at, block;
726 gmx::EnumerationArray<SimulationAtomGroupType, std::vector<unsigned char>> new_egrp;
730 /* Construct the molecule range information */
731 gmx::RangePartitioning molecules = gmx_mtop_molecules(*mtop);
733 snew(list, state->natoms);
735 for (int i = 0; i < rm_p->nr; i++)
737 mol_id = rm_p->mol[i];
738 at = molecules.block(mol_id).begin();
739 block = rm_p->block[i];
740 mtop->molblock[block].nmol--;
741 for (j = 0; j < mtop->moltype[mtop->molblock[block].type].atoms.nr; j++)
749 /* We cannot change the size of the state datastructures here
750 * because we still access the coordinate arrays for all positions
751 * before removing the molecules we want to remove.
753 const int newStateAtomNumber = state->natoms - n;
754 snew(x_tmp, newStateAtomNumber);
755 snew(v_tmp, newStateAtomNumber);
757 for (auto group : keysOf(groups->groupNumbers))
759 if (!groups->groupNumbers[group].empty())
761 groups->groupNumbers[group].resize(newStateAtomNumber);
762 new_egrp[group].resize(newStateAtomNumber);
766 auto x = makeArrayRef(state->x);
767 auto v = makeArrayRef(state->v);
769 for (int i = 0; i < state->natoms; i++)
772 for (j = 0; j < n; j++)
783 for (auto group : keysOf(groups->groupNumbers))
785 if (!groups->groupNumbers[group].empty())
787 new_egrp[group][i - rm] = groups->groupNumbers[group][i];
790 copy_rvec(x[i], x_tmp[i - rm]);
791 copy_rvec(v[i], v_tmp[i - rm]);
792 for (j = 0; j < ins_at->nr; j++)
794 if (i == ins_at->index[j])
796 ins_at->index[j] = i - rm;
800 for (j = 0; j < pos_ins->pieces; j++)
802 for (k = 0; k < pos_ins->nidx[j]; k++)
804 if (i == pos_ins->subindex[j][k])
806 pos_ins->subindex[j][k] = i - rm;
812 state_change_natoms(state, newStateAtomNumber);
813 for (int i = 0; i < state->natoms; i++)
815 copy_rvec(x_tmp[i], x[i]);
818 for (int i = 0; i < state->natoms; i++)
820 copy_rvec(v_tmp[i], v[i]);
824 for (auto group : keysOf(groups->groupNumbers))
826 if (!groups->groupNumbers[group].empty())
828 groups->groupNumbers[group] = new_egrp[group];
832 /* remove empty molblocks */
834 for (size_t i = 0; i < mtop->molblock.size(); i++)
836 if (mtop->molblock[i].nmol == 0)
842 mtop->molblock[i - RMmolblock] = mtop->molblock[i];
845 mtop->molblock.resize(mtop->molblock.size() - RMmolblock);
848 /* remove al bonded interactions from mtop for the molecule to be embedded */
849 static int rm_bonded(t_block* ins_at, gmx_mtop_t* mtop)
852 int type, natom, nmol, at, atom1 = 0, rm_at = 0;
854 /*this routine lives dangerously by assuming that all molecules of a given type are in order in the structure*/
855 /*this routine does not live as dangerously as it seems. There is namely a check in init_membed
856 * to make * sure that g_membed exits with a warning when there are molecules of the same type
857 * not in the * ins_at index group. MGWolf 050710 */
860 snew(bRM, mtop->moltype.size());
861 for (size_t i = 0; i < mtop->moltype.size(); i++)
866 for (size_t i = 0; i < mtop->molblock.size(); i++)
868 /*loop over molecule blocks*/
869 type = mtop->molblock[i].type;
870 natom = mtop->moltype[type].atoms.nr;
871 nmol = mtop->molblock[i].nmol;
873 for (j = 0; j < natom * nmol && bRM[type]; j++)
875 /*loop over atoms in the block*/
876 at = j + atom1; /*atom index = block index + offset*/
879 for (m = 0; (m < ins_at->nr) && (!bINS); m++)
881 /*loop over atoms in insertion index group to determine if we're inserting one*/
882 if (at == ins_at->index[m])
889 atom1 += natom * nmol; /*update offset*/
892 rm_at += natom * nmol; /*increment bonded removal counter by # atoms in block*/
896 for (size_t i = 0; i < mtop->moltype.size(); i++)
900 for (j = 0; j < F_LJ; j++)
902 mtop->moltype[i].ilist[j].iatoms.clear();
905 for (j = F_POSRES; j <= F_VSITEN; j++)
907 mtop->moltype[i].ilist[j].iatoms.clear();
916 /* Write a topology where the number of molecules is correct for the system after embedding */
917 static void top_update(const char* topfile, rm_t* rm_p, gmx_mtop_t* mtop)
921 char buf[STRLEN], buf2[STRLEN], *temp;
922 int i, *nmol_rm, nmol, line;
923 char temporary_filename[STRLEN];
925 fpin = gmx_ffopen(topfile, "r");
926 strncpy(temporary_filename, "temp.topXXXXXX", STRLEN);
927 gmx_tmpnam(temporary_filename);
928 fpout = gmx_ffopen(temporary_filename, "w");
930 snew(nmol_rm, mtop->moltype.size());
931 for (i = 0; i < rm_p->nr; i++)
933 nmol_rm[rm_p->block[i]]++;
937 while (fgets(buf, STRLEN, fpin))
943 if ((temp = strchr(buf2, '\n')) != nullptr)
951 if ((temp = strchr(buf2, '\n')) != nullptr)
956 if (buf2[strlen(buf2) - 1] == ']')
958 buf2[strlen(buf2) - 1] = '\0';
961 if (gmx_strcasecmp(buf2, "molecules") == 0)
966 fprintf(fpout, "%s", buf);
968 else if (bMolecules == 1)
970 for (size_t i = 0; i < mtop->molblock.size(); i++)
972 nmol = mtop->molblock[i].nmol;
973 sprintf(buf, "%-15s %5d\n", *(mtop->moltype[mtop->molblock[i].type].name), nmol);
974 fprintf(fpout, "%s", buf);
978 else if (bMolecules == 2)
984 fprintf(fpout, "%s", buf);
989 fprintf(fpout, "%s", buf);
994 /* use gmx_ffopen to generate backup of topinout */
995 fpout = gmx_ffopen(topfile, "w");
997 rename(temporary_filename, topfile);
1000 void rescale_membed(int step_rel, gmx_membed_t* membed, rvec* x)
1002 /* Set new positions for the group to embed */
1003 if (step_rel <= membed->it_xy)
1005 membed->fac[0] += membed->xy_step;
1006 membed->fac[1] += membed->xy_step;
1008 else if (step_rel <= (membed->it_xy + membed->it_z))
1010 membed->fac[2] += membed->z_step;
1012 resize(membed->r_ins, x, membed->pos_ins, membed->fac);
1015 gmx_membed_t* init_membed(FILE* fplog,
1017 const t_filenm fnm[],
1019 t_inputrec* inputrec,
1025 int i, rm_bonded_at, fr_id, fr_i = 0, tmp_id, warn = 0;
1026 int ng, j, max_lip_rm, ins_grp_id, ntype, lip_rm;
1028 rvec* r_ins = nullptr;
1029 t_block * ins_at, *rest_at;
1033 SimulationGroups* groups;
1034 gmx_bool bExcl = FALSE;
1037 char** piecename = nullptr;
1038 gmx_membed_t* membed = nullptr;
1040 /* input variables */
1047 real probe_rad = 0.22;
1051 gmx_bool bALLOW_ASYMMETRY = FALSE;
1053 /* sanity check constants */ /* Issue a warning when: */
1054 const real min_probe_rad = 0.2199999; /* A probe radius for overlap between embedded molecule *
1055 * and rest smaller than this value is probably too small */
1056 const real min_xy_init = 0.0999999; /* the initial shrinking of the molecule to embed is smaller */
1057 const int min_it_xy = 1000; /* the number of steps to embed in xy-plane is smaller */
1058 const int min_it_z = 100; /* the number of steps to embed in z is smaller */
1059 const real prot_vs_box = 7.5; /* molecule to embed is large (more then prot_vs_box) with respect */
1060 const real box_vs_prot = 50; /* to the box size (less than box_vs_prot) */
1069 "Note: it is expected that in future gmx mdrun -membed will not be the "
1070 "way to access this feature, perhaps changing to e.g. gmx membed.");
1071 /* get input data out membed file */
1074 get_input(opt2fn("-membed", nfile, fnm),
1087 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1089 if (!EI_DYNAMICS(inputrec->eI))
1091 gmx_input("Change integrator to a dynamics integrator in mdp file (e.g. md or sd).");
1096 gmx_input("Sorry, parallel membed is not yet fully functional.");
1102 "\nSetting -cpt to -1, because embedding cannot be restarted from "
1106 groups = &(mtop->groups);
1107 std::vector<std::string> gnames;
1108 for (const auto& groupName : groups->groupNames)
1110 gnames.emplace_back(*groupName);
1113 atoms = gmx_mtop_global_atoms(mtop);
1115 fprintf(stderr, "\nSelect a group to embed in the membrane:\n");
1116 get_index(&atoms, opt2fn_null("-mn", nfile, fnm), 1, &(ins_at->nr), &(ins_at->index), &ins);
1118 auto found = std::find_if(gnames.begin(), gnames.end(), [&ins](const auto& name) {
1119 return gmx::equalCaseInsensitive(ins, name);
1122 if (found == gnames.end())
1125 "Group %s selected for embedding was not found in the index file.\n"
1126 "Group names must match either [moleculetype] names or custom index group\n"
1127 "names, in which case you must supply an index file to the '-n' option\n"
1133 ins_grp_id = std::distance(gnames.begin(), found);
1135 fprintf(stderr, "\nSelect a group to embed %s into (e.g. the membrane):\n", ins);
1137 opt2fn_null("-mn", nfile, fnm),
1139 &(mem_p->mem_at.nr),
1140 &(mem_p->mem_at.index),
1143 pos_ins->pieces = pieces;
1144 snew(pos_ins->nidx, pieces);
1145 snew(pos_ins->subindex, pieces);
1146 snew(piecename, pieces);
1149 fprintf(stderr, "\nSelect pieces to embed:\n");
1150 get_index(&atoms, opt2fn_null("-mn", nfile, fnm), pieces, pos_ins->nidx, pos_ins->subindex, piecename);
1154 /*use whole embedded group*/
1155 snew(pos_ins->nidx, 1);
1156 snew(pos_ins->subindex, 1);
1157 pos_ins->nidx[0] = ins_at->nr;
1158 pos_ins->subindex[0] = ins_at->index;
1161 if (probe_rad < min_probe_rad)
1165 "\nWarning %d:\nA probe radius (-rad) smaller than 0.2 nm can result "
1166 "in overlap between waters and the group to embed, which will result "
1167 "in Lincs errors etc.\n\n",
1171 if (xy_fac < min_xy_init)
1174 fprintf(stderr, "\nWarning %d:\nThe initial size of %s is probably too small.\n\n", warn, ins);
1177 if (it_xy < min_it_xy)
1181 "\nWarning %d;\nThe number of steps used to grow the xy-coordinates of %s (%d)"
1182 " is probably too small.\nIncrease -nxy or.\n\n",
1188 if ((it_z < min_it_z) && (z_fac < 0.99999999 || z_fac > 1.0000001))
1192 "\nWarning %d;\nThe number of steps used to grow the z-coordinate of %s (%d)"
1193 " is probably too small.\nIncrease -nz or the maxwarn setting in the membed "
1200 if (it_xy + it_z > inputrec->nsteps)
1204 "\nWarning %d:\nThe number of growth steps (-nxy + -nz) is larger than the "
1205 "number of steps in the tpr.\n"
1206 "(increase maxwarn in the membed input file to override)\n\n",
1211 if (inputrec->opts.ngfrz == 1)
1213 gmx_fatal(FARGS, "You did not specify \"%s\" as a freezegroup.", ins);
1216 for (i = 0; i < inputrec->opts.ngfrz; i++)
1218 tmp_id = mtop->groups.groups[SimulationAtomGroupType::Freeze][i];
1219 if (ins_grp_id == tmp_id)
1228 gmx_fatal(FARGS, "\"%s\" not as freezegroup defined in the mdp-file.", ins);
1231 for (i = 0; i < DIM; i++)
1233 if (inputrec->opts.nFreeze[fr_i][i] != 1)
1235 gmx_fatal(FARGS, "freeze dimensions for %s are not Y Y Y\n", ins);
1239 ng = groups->groups[SimulationAtomGroupType::EnergyOutput].size();
1243 "No energy groups defined. This is necessary for energy exclusion in the "
1247 for (i = 0; i < ng; i++)
1249 for (j = 0; j < ng; j++)
1251 if (inputrec->opts.egp_flags[ng * i + j] == EGP_EXCL)
1254 if ((groups->groups[SimulationAtomGroupType::EnergyOutput][i] != ins_grp_id)
1255 || (groups->groups[SimulationAtomGroupType::EnergyOutput][j] != ins_grp_id))
1258 "Energy exclusions \"%s\" and \"%s\" do not match the group "
1260 *groups->groupNames[groups->groups[SimulationAtomGroupType::EnergyOutput][i]],
1261 *groups->groupNames[groups->groups[SimulationAtomGroupType::EnergyOutput][j]],
1271 "No energy exclusion groups defined. This is necessary for energy exclusion in "
1272 "the freeze group");
1275 /* Obtain the maximum and minimum coordinates of the group to be embedded */
1277 init_ins_at(ins_at, rest_at, state, pos_ins, groups, ins_grp_id, xy_max);
1278 /* Check that moleculetypes in insertion group are not part of the rest of the system */
1279 check_types(ins_at, rest_at, mtop);
1281 init_mem_at(mem_p, mtop, state->x.rvec_array(), state->box, pos_ins);
1283 prot_area = est_prot_area(pos_ins, state->x.rvec_array(), ins_at, mem_p);
1284 if ((prot_area > prot_vs_box)
1285 && ((state->box[XX][XX] * state->box[YY][YY] - state->box[XX][YY] * state->box[YY][XX])
1290 "\nWarning %d:\nThe xy-area is very small compared to the area of the "
1292 "This might cause pressure problems during the growth phase. Just try with\n"
1293 "current setup and increase 'maxwarn' in your membed settings file, but lower "
1295 "compressibility in the mdp-file or disable pressure coupling if problems "
1303 "Too many warnings (override by setting maxwarn in the membed input file)\n");
1306 printf("The estimated area of the protein in the membrane is %.3f nm^2\n", prot_area);
1307 printf("\nThere are %d lipids in the membrane part that overlaps the protein.\n"
1308 "The area per lipid is %.4f nm^2.\n",
1312 /* Maximum number of lipids to be removed*/
1313 max_lip_rm = static_cast<int>(2 * prot_area / mem_p->lip_area);
1314 printf("Maximum number of lipids that will be removed is %d.\n", max_lip_rm);
1316 printf("\nWill resize the protein by a factor of %.3f in the xy plane and %.3f in the z "
1318 "This resizing will be done with respect to the geometrical center of all protein "
1320 "that span the membrane region, i.e. z between %.3f and %.3f\n\n",
1326 /* resize the protein by xy and by z if necessary*/
1327 snew(r_ins, ins_at->nr);
1328 init_resize(ins_at, r_ins, pos_ins, mem_p, state->x.rvec_array(), bALLOW_ASYMMETRY);
1329 membed->fac[0] = membed->fac[1] = xy_fac;
1330 membed->fac[2] = z_fac;
1332 membed->xy_step = (xy_max - xy_fac) / static_cast<double>(it_xy);
1333 membed->z_step = (z_max - z_fac) / static_cast<double>(it_z - 1);
1335 resize(r_ins, state->x.rvec_array(), pos_ins, membed->fac);
1337 /* remove overlapping lipids and water from the membrane box*/
1338 /*mark molecules to be removed*/
1340 set_pbc(pbc, inputrec->pbcType, state->box);
1343 lip_rm = gen_rm_list(
1344 rm_p, ins_at, rest_at, pbc, mtop, state->x.rvec_array(), mem_p, pos_ins, probe_rad, low_up_rm, bALLOW_ASYMMETRY);
1345 lip_rm -= low_up_rm;
1349 for (i = 0; i < rm_p->nr; i++)
1351 fprintf(fplog, "rm mol %d\n", rm_p->mol[i]);
1355 for (size_t i = 0; i < mtop->molblock.size(); i++)
1358 for (j = 0; j < rm_p->nr; j++)
1360 if (rm_p->block[j] == static_cast<int>(i))
1365 printf("Will remove %d %s molecules\n", ntype, *(mtop->moltype[mtop->molblock[i].type].name));
1368 if (lip_rm > max_lip_rm)
1372 "\nWarning %d:\nTrying to remove a larger lipid area than the estimated "
1373 "protein area\nTry making the -xyinit resize factor smaller or increase "
1374 "maxwarn in the membed input file.\n\n",
1378 /*remove all lipids and waters overlapping and update all important structures*/
1379 rm_group(groups, mtop, rm_p, state, ins_at, pos_ins);
1381 rm_bonded_at = rm_bonded(ins_at, mtop);
1382 if (rm_bonded_at != ins_at->nr)
1385 "Warning: The number of atoms for which the bonded interactions are removed is "
1387 "while %d atoms are embedded. Make sure that the atoms to be embedded are not "
1389 "molecule type as atoms that are not to be embedded.\n",
1397 "Too many warnings.\nIf you are sure these warnings are harmless,\n"
1398 "you can increase the maxwarn setting in the membed input file.");
1401 // Re-establish the invariants of the derived values within
1405 if (ftp2bSet(efTOP, nfile, fnm))
1407 top_update(opt2fn("-mp", nfile, fnm), rm_p, mtop);
1417 membed->it_xy = it_xy;
1418 membed->it_z = it_z;
1419 membed->pos_ins = pos_ins;
1420 membed->r_ins = r_ins;
1426 void free_membed(gmx_membed_t* membed)