Replace t_grps with std::vector
[alexxy/gromacs.git] / src / gromacs / mdlib / membed.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2010,2011,2012,2013,2014,2015,2016,2017,2018,2019, 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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 #include "gmxpre.h"
36
37 #include "membed.h"
38
39 #include <csignal>
40 #include <cstdlib>
41
42 #include "gromacs/commandline/filenm.h"
43 #include "gromacs/fileio/readinp.h"
44 #include "gromacs/fileio/warninp.h"
45 #include "gromacs/gmxlib/network.h"
46 #include "gromacs/math/vec.h"
47 #include "gromacs/mdtypes/commrec.h"
48 #include "gromacs/mdtypes/inputrec.h"
49 #include "gromacs/mdtypes/md_enums.h"
50 #include "gromacs/mdtypes/state.h"
51 #include "gromacs/pbcutil/pbc.h"
52 #include "gromacs/topology/index.h"
53 #include "gromacs/topology/mtop_lookup.h"
54 #include "gromacs/topology/mtop_util.h"
55 #include "gromacs/topology/topology.h"
56 #include "gromacs/utility/cstringutil.h"
57 #include "gromacs/utility/exceptions.h"
58 #include "gromacs/utility/fatalerror.h"
59 #include "gromacs/utility/filestream.h"
60 #include "gromacs/utility/futil.h"
61 #include "gromacs/utility/smalloc.h"
62
63 /* information about scaling center */
64 typedef struct {
65     rvec      xmin;      /* smallest coordinates of all embedded molecules */
66     rvec      xmax;      /* largest coordinates of all embedded molecules */
67     rvec     *geom_cent; /* scaling center of each independent molecule to embed */
68     int       pieces;    /* number of molecules to embed independently */
69     int      *nidx;      /* n atoms for every independent embedded molecule (index in subindex) */
70     int     **subindex;  /* atomids for independent molecule *
71                           * atoms of piece i run from subindex[i][0] to subindex[i][nidx[i]] */
72 } pos_ins_t;
73
74 /* variables needed in do_md */
75 struct gmx_membed_t {
76     int        it_xy;     /* number of iterations (steps) used to grow something in the xy-plane */
77     int        it_z;      /* same, but for z */
78     real       xy_step;   /* stepsize used during resize in xy-plane */
79     real       z_step;    /* same, but in z */
80     rvec       fac;       /* initial scaling of the molecule to grow into the membrane */
81     rvec      *r_ins;     /* final coordinates of the molecule to grow  */
82     pos_ins_t *pos_ins;   /* scaling center for each piece to embed */
83 };
84
85 /* membrane related variables */
86 typedef struct {
87     char      *name;     /* name of index group to embed molecule into (usually membrane) */
88     t_block    mem_at;   /* list all atoms in membrane */
89     int        nmol;     /* number of membrane molecules overlapping with the molecule to embed */
90     int       *mol_id;   /* list of molecules in membrane that overlap with the molecule to embed */
91     real       lip_area; /* average area per lipid in membrane (only correct for homogeneous bilayers)*/
92     real       zmin;     /* minimum z coordinate of membrane */
93     real       zmax;     /* maximum z coordinate of membrane */
94     real       zmed;     /* median z coordinate of membrane */
95 } mem_t;
96
97 /* Lists all molecules in the membrane that overlap with the molecule to be embedded. *
98  * These will then be removed from the system */
99 typedef struct {
100     int   nr;     /* number of molecules to remove */
101     int  *mol;    /* list of molecule ids to remove */
102     int  *block;  /* id of the molblock that the molecule to remove is part of */
103 } rm_t;
104
105 /* Get the global molecule id, and the corresponding molecule type and id of the *
106  * molblock from the global atom nr. */
107 static int get_mol_id(int at, gmx_mtop_t  *mtop, int *type, int *block)
108 {
109     int                   mol_id = 0;
110     int                   i;
111     int                   atnr_mol;
112
113     *block = 0;
114     mtopGetMolblockIndex(mtop, at, block, &mol_id, &atnr_mol);
115     for (i = 0; i < *block; i++)
116     {
117         mol_id += mtop->molblock[i].nmol;
118     }
119     *type = mtop->molblock[*block].type;
120
121     return mol_id;
122 }
123
124 /* Get the id of the molblock from a global molecule id */
125 static int get_molblock(int mol_id, const std::vector<gmx_molblock_t> &mblock)
126 {
127     int nmol = 0;
128
129     for (size_t i = 0; i < mblock.size(); i++)
130     {
131         nmol += mblock[i].nmol;
132         if (mol_id < nmol)
133         {
134             return i;
135         }
136     }
137
138     gmx_fatal(FARGS, "mol_id %d larger than total number of molecules %d.\n", mol_id, nmol);
139 }
140
141 /* Get a list of all the molecule types that are present in a group of atoms. *
142  * Because all interaction within the group to embed are removed on the topology *
143  * level, if the same molecule type is found in another part of the system, these *
144  * would also be affected. Therefore we have to check if the embedded and rest group *
145  * share common molecule types. If so, membed will stop with an error. */
146 static int get_mtype_list(t_block *at, gmx_mtop_t *mtop, t_block *tlist)
147 {
148     int      i, j, nr;
149     int      type = 0, block = 0;
150     gmx_bool bNEW;
151
152     nr = 0;
153     snew(tlist->index, at->nr);
154     for (i = 0; i < at->nr; i++)
155     {
156         bNEW   = TRUE;
157         get_mol_id(at->index[i], mtop, &type, &block);
158         for (j = 0; j < nr; j++)
159         {
160             if (tlist->index[j] == type)
161             {
162                 bNEW = FALSE;
163             }
164         }
165
166         if (bNEW)
167         {
168             tlist->index[nr] = type;
169             nr++;
170         }
171     }
172     srenew(tlist->index, nr);
173     return nr;
174 }
175
176 /* Do the actual check of the molecule types between embedded and rest group */
177 static void check_types(t_block *ins_at, t_block *rest_at, gmx_mtop_t *mtop)
178 {
179     t_block        *ins_mtype, *rest_mtype;
180     int             i, j;
181
182     snew(ins_mtype, 1);
183     snew(rest_mtype, 1);
184     ins_mtype->nr  = get_mtype_list(ins_at, mtop, ins_mtype );
185     rest_mtype->nr = get_mtype_list(rest_at, mtop, rest_mtype);
186
187     for (i = 0; i < ins_mtype->nr; i++)
188     {
189         for (j = 0; j < rest_mtype->nr; j++)
190         {
191             if (ins_mtype->index[i] == rest_mtype->index[j])
192             {
193                 gmx_fatal(FARGS, "Moleculetype %s is found both in the group to insert and the rest of the system.\n"
194                           "1. Your *.ndx and *.top do not match\n"
195                           "2. You are inserting some molecules of type %s (for example xray-solvent), while\n"
196                           "the same moleculetype is also used in the rest of the system (solvent box). Because\n"
197                           "we need to exclude all interactions between the atoms in the group to\n"
198                           "insert, the same moleculetype can not be used in both groups. Change the\n"
199                           "moleculetype of the molecules %s in the inserted group. Do not forget to provide\n"
200                           "an appropriate *.itp file", *(mtop->moltype[rest_mtype->index[j]].name),
201                           *(mtop->moltype[rest_mtype->index[j]].name), *(mtop->moltype[rest_mtype->index[j]].name));
202             }
203         }
204     }
205
206     done_block(ins_mtype);
207     done_block(rest_mtype);
208     sfree(ins_mtype);
209     sfree(rest_mtype);
210 }
211
212 static void get_input(const char *membed_input, real *xy_fac, real *xy_max, real *z_fac, real *z_max,
213                       int *it_xy, int *it_z, real *probe_rad, int *low_up_rm, int *maxwarn,
214                       int *pieces, gmx_bool *bALLOW_ASYMMETRY)
215 {
216     warninp_t               wi;
217     std::vector <t_inpfile> inp;
218
219     wi = init_warning(TRUE, 0);
220
221     {
222         gmx::TextInputFile stream(membed_input);
223         inp = read_inpfile(&stream, membed_input, wi);
224         stream.close();
225     }
226     *it_xy            = get_eint(&inp, "nxy", 1000, wi);
227     *it_z             = get_eint(&inp, "nz", 0, wi);
228     *xy_fac           = get_ereal(&inp, "xyinit", 0.5, wi);
229     *xy_max           = get_ereal(&inp, "xyend", 1.0, wi);
230     *z_fac            = get_ereal(&inp, "zinit", 1.0, wi);
231     *z_max            = get_ereal(&inp, "zend", 1.0, wi);
232     *probe_rad        = get_ereal(&inp, "rad", 0.22, wi);
233     *low_up_rm        = get_eint(&inp, "ndiff", 0, wi);
234     *maxwarn          = get_eint(&inp, "maxwarn", 0, wi);
235     *pieces           = get_eint(&inp, "pieces", 1, wi);
236     const char *yesno_names[BOOL_NR+1] =
237     {
238         "no", "yes", nullptr
239     };
240     *bALLOW_ASYMMETRY = (get_eeenum(&inp, "asymmetry", yesno_names, wi) != 0);
241
242     check_warning_error(wi, FARGS);
243     {
244         gmx::TextOutputFile stream(membed_input);
245         write_inpfile(&stream, membed_input, &inp, FALSE, WriteMdpHeader::yes, wi);
246         stream.close();
247     }
248     done_warning(wi, FARGS);
249 }
250
251 /* Obtain the maximum and minimum coordinates of the group to be embedded */
252 static int init_ins_at(t_block *ins_at, t_block *rest_at, t_state *state, pos_ins_t *pos_ins,
253                        SimulationGroups *groups, int ins_grp_id, real xy_max)
254 {
255     int        i, gid, c = 0;
256     real       xmin, xmax, ymin, ymax, zmin, zmax;
257     const real min_memthick = 6.0;      /* minimum thickness of the bilayer that will be used to *
258                                          * determine the overlap between molecule to embed and membrane */
259     const real fac_inp_size = 1.000001; /* scaling factor to obtain input_size + 0.000001 (comparing reals) */
260     snew(rest_at->index, state->natoms);
261     auto       x = makeArrayRef(state->x);
262
263     xmin = xmax = x[ins_at->index[0]][XX];
264     ymin = ymax = x[ins_at->index[0]][YY];
265     zmin = zmax = x[ins_at->index[0]][ZZ];
266
267     for (i = 0; i < state->natoms; i++)
268     {
269         gid = groups->groupNumbers[SimulationAtomGroupType::Freeze][i];
270         if (groups->groups[SimulationAtomGroupType::Freeze][gid] == ins_grp_id)
271         {
272             xmin = std::min(xmin, x[i][XX]);
273             xmax = std::max(xmax, x[i][XX]);
274             ymin = std::min(ymin, x[i][YY]);
275             ymax = std::max(ymax, x[i][YY]);
276             zmin = std::min(zmin, x[i][ZZ]);
277             zmax = std::max(zmax, x[i][ZZ]);
278         }
279         else
280         {
281             rest_at->index[c] = i;
282             c++;
283         }
284     }
285
286     rest_at->nr = c;
287     srenew(rest_at->index, c);
288
289     if (xy_max > fac_inp_size)
290     {
291         pos_ins->xmin[XX] = xmin-((xmax-xmin)*xy_max-(xmax-xmin))/2;
292         pos_ins->xmin[YY] = ymin-((ymax-ymin)*xy_max-(ymax-ymin))/2;
293
294         pos_ins->xmax[XX] = xmax+((xmax-xmin)*xy_max-(xmax-xmin))/2;
295         pos_ins->xmax[YY] = ymax+((ymax-ymin)*xy_max-(ymax-ymin))/2;
296     }
297     else
298     {
299         pos_ins->xmin[XX] = xmin;
300         pos_ins->xmin[YY] = ymin;
301
302         pos_ins->xmax[XX] = xmax;
303         pos_ins->xmax[YY] = ymax;
304     }
305
306     if ( (zmax-zmin) < min_memthick)
307     {
308         pos_ins->xmin[ZZ] = zmin+(zmax-zmin)/2.0-0.5*min_memthick;
309         pos_ins->xmax[ZZ] = zmin+(zmax-zmin)/2.0+0.5*min_memthick;
310     }
311     else
312     {
313         pos_ins->xmin[ZZ] = zmin;
314         pos_ins->xmax[ZZ] = zmax;
315     }
316
317     return c;
318 }
319
320 /* Estimate the area of the embedded molecule by projecting all coordinates on a grid in the *
321  * xy-plane and counting the number of occupied grid points */
322 static real est_prot_area(pos_ins_t *pos_ins, rvec *r, t_block *ins_at, mem_t *mem_p)
323 {
324     real x, y, dx = 0.15, dy = 0.15, area = 0.0;
325     real add, memmin, memmax;
326     int  c, at;
327
328     /* min and max membrane coordinate are altered to reduce the influence of the *
329      * boundary region */
330     memmin = mem_p->zmin+0.1*(mem_p->zmax-mem_p->zmin);
331     memmax = mem_p->zmax-0.1*(mem_p->zmax-mem_p->zmin);
332
333     //NOLINTNEXTLINE(clang-analyzer-security.FloatLoopCounter)
334     for (x = pos_ins->xmin[XX]; x < pos_ins->xmax[XX]; x += dx)
335     {
336         //NOLINTNEXTLINE(clang-analyzer-security.FloatLoopCounter)
337         for (y = pos_ins->xmin[YY]; y < pos_ins->xmax[YY]; y += dy)
338         {
339             c   = 0;
340             add = 0.0;
341             do
342             {
343                 at = ins_at->index[c];
344                 if ( (r[at][XX] >= x) && (r[at][XX] < x+dx) &&
345                      (r[at][YY] >= y) && (r[at][YY] < y+dy) &&
346                      (r[at][ZZ] > memmin) && (r[at][ZZ] < memmax) )
347                 {
348                     add = 1.0;
349                 }
350                 c++;
351             }
352             while ( (c < ins_at->nr) && (add < 0.5) );
353             area += add;
354         }
355     }
356     area = area*dx*dy;
357
358     return area;
359 }
360
361 static int init_mem_at(mem_t *mem_p, gmx_mtop_t *mtop, rvec *r, matrix box, pos_ins_t *pos_ins)
362 {
363     int      i, j, at, mol, nmol, nmolbox, count;
364     t_block *mem_a;
365     real     z, zmin, zmax, mem_area;
366     gmx_bool bNew;
367     int     *mol_id;
368     int      type = 0, block = 0;
369
370     nmol  = count = 0;
371     mem_a = &(mem_p->mem_at);
372     snew(mol_id, mem_a->nr);
373     zmin = pos_ins->xmax[ZZ];
374     zmax = pos_ins->xmin[ZZ];
375     for (i = 0; i < mem_a->nr; i++)
376     {
377         at = mem_a->index[i];
378         if ( (r[at][XX] > pos_ins->xmin[XX]) && (r[at][XX] < pos_ins->xmax[XX]) &&
379              (r[at][YY] > pos_ins->xmin[YY]) && (r[at][YY] < pos_ins->xmax[YY]) &&
380              (r[at][ZZ] > pos_ins->xmin[ZZ]) && (r[at][ZZ] < pos_ins->xmax[ZZ]) )
381         {
382             mol  = get_mol_id(at, mtop, &type, &block);
383             bNew = TRUE;
384             for (j = 0; j < nmol; j++)
385             {
386                 if (mol == mol_id[j])
387                 {
388                     bNew = FALSE;
389                 }
390             }
391
392             if (bNew)
393             {
394                 mol_id[nmol] = mol;
395                 nmol++;
396             }
397
398             z = r[at][ZZ];
399             if (z < zmin)
400             {
401                 zmin = z;
402             }
403
404             if (z > zmax)
405             {
406                 zmax = z;
407             }
408
409             count++;
410         }
411     }
412
413     mem_p->nmol = nmol;
414     srenew(mol_id, nmol);
415     mem_p->mol_id = mol_id;
416
417     if ((zmax-zmin) > (box[ZZ][ZZ]-0.5))
418     {
419         gmx_fatal(FARGS, "Something is wrong with your membrane. Max and min z values are %f and %f.\n"
420                   "Maybe your membrane is not centered in the box, but located at the box edge in the z-direction,\n"
421                   "so that one membrane is distributed over two periodic box images. Another possibility is that\n"
422                   "your water layer is not thick enough.\n", zmax, zmin);
423     }
424     mem_p->zmin = zmin;
425     mem_p->zmax = zmax;
426     mem_p->zmed = (zmax-zmin)/2+zmin;
427
428     /*number of membrane molecules in protein box*/
429     nmolbox = count/mtop->moltype[mtop->molblock[block].type].atoms.nr;
430     /*membrane area within the box defined by the min and max coordinates of the embedded molecule*/
431     mem_area = (pos_ins->xmax[XX]-pos_ins->xmin[XX])*(pos_ins->xmax[YY]-pos_ins->xmin[YY]);
432     /*rough estimate of area per lipid, assuming there is only one type of lipid in the membrane*/
433     mem_p->lip_area = 2.0*mem_area/static_cast<double>(nmolbox);
434
435     return mem_p->mem_at.nr;
436 }
437
438 static void init_resize(t_block *ins_at, rvec *r_ins, pos_ins_t *pos_ins, mem_t *mem_p, rvec *r,
439                         gmx_bool bALLOW_ASYMMETRY)
440 {
441     int i, j, at, c, outsidesum, gctr = 0;
442     int idxsum = 0;
443
444     /*sanity check*/
445     for (i = 0; i < pos_ins->pieces; i++)
446     {
447         idxsum += pos_ins->nidx[i];
448     }
449
450     if (idxsum != ins_at->nr)
451     {
452         gmx_fatal(FARGS, "Piecewise sum of inserted atoms not same as size of group selected to insert.");
453     }
454
455     snew(pos_ins->geom_cent, pos_ins->pieces);
456     for (i = 0; i < pos_ins->pieces; i++)
457     {
458         c          = 0;
459         outsidesum = 0;
460         for (j = 0; j < DIM; j++)
461         {
462             pos_ins->geom_cent[i][j] = 0;
463         }
464
465         for (j = 0; j < pos_ins->nidx[i]; j++)
466         {
467             at = pos_ins->subindex[i][j];
468             copy_rvec(r[at], r_ins[gctr]);
469             if ( (r_ins[gctr][ZZ] < mem_p->zmax) && (r_ins[gctr][ZZ] > mem_p->zmin) )
470             {
471                 rvec_inc(pos_ins->geom_cent[i], r_ins[gctr]);
472                 c++;
473             }
474             else
475             {
476                 outsidesum++;
477             }
478             gctr++;
479         }
480
481         if (c > 0)
482         {
483             svmul(1/static_cast<double>(c), pos_ins->geom_cent[i], pos_ins->geom_cent[i]);
484         }
485
486         if (!bALLOW_ASYMMETRY)
487         {
488             pos_ins->geom_cent[i][ZZ] = mem_p->zmed;
489         }
490
491         fprintf(stderr, "Embedding piece %d with center of geometry: %f %f %f\n",
492                 i, pos_ins->geom_cent[i][XX], pos_ins->geom_cent[i][YY], pos_ins->geom_cent[i][ZZ]);
493     }
494     fprintf(stderr, "\n");
495 }
496
497 /* resize performed in the md loop */
498 static void resize(rvec *r_ins, rvec *r, pos_ins_t *pos_ins, const rvec fac)
499 {
500     int i, j, k, at, c = 0;
501     for (k = 0; k < pos_ins->pieces; k++)
502     {
503         for (i = 0; i < pos_ins->nidx[k]; i++)
504         {
505             at = pos_ins->subindex[k][i];
506             for (j = 0; j < DIM; j++)
507             {
508                 r[at][j] = pos_ins->geom_cent[k][j]+fac[j]*(r_ins[c][j]-pos_ins->geom_cent[k][j]);
509             }
510             c++;
511         }
512     }
513 }
514
515 /* generate the list of membrane molecules that overlap with the molecule to be embedded. *
516  * The molecule to be embedded is already reduced in size. */
517 static int gen_rm_list(rm_t *rm_p, t_block *ins_at, t_block *rest_at, t_pbc *pbc, gmx_mtop_t *mtop,
518                        rvec *r, mem_t *mem_p, pos_ins_t *pos_ins, real probe_rad,
519                        int low_up_rm, gmx_bool bALLOW_ASYMMETRY)
520 {
521     int      i, j, k, l, at, at2, mol_id;
522     int      type = 0, block = 0;
523     int      nrm, nupper, nlower;
524     real     r_min_rad, z_lip, min_norm;
525     gmx_bool bRM;
526     rvec     dr, dr_tmp;
527     real    *dist;
528     int     *order;
529
530     r_min_rad = probe_rad*probe_rad;
531     gmx::RangePartitioning molecules = gmx_mtop_molecules(*mtop);
532     snew(rm_p->block, molecules.numBlocks());
533     nrm    = nupper = 0;
534     nlower = low_up_rm;
535     for (i = 0; i < ins_at->nr; i++)
536     {
537         at = ins_at->index[i];
538         for (j = 0; j < rest_at->nr; j++)
539         {
540             at2 = rest_at->index[j];
541             pbc_dx(pbc, r[at], r[at2], dr);
542
543             if (norm2(dr) < r_min_rad)
544             {
545                 mol_id = get_mol_id(at2, mtop, &type, &block);
546                 bRM    = TRUE;
547                 for (l = 0; l < nrm; l++)
548                 {
549                     if (rm_p->mol[l] == mol_id)
550                     {
551                         bRM = FALSE;
552                     }
553                 }
554
555                 if (bRM)
556                 {
557                     rm_p->mol[nrm]   = mol_id;
558                     rm_p->block[nrm] = block;
559                     nrm++;
560                     z_lip = 0.0;
561                     for (l = 0; l < mem_p->nmol; l++)
562                     {
563                         if (mol_id == mem_p->mol_id[l])
564                         {
565                             for (int k : molecules.block(mol_id))
566                             {
567                                 z_lip += r[k][ZZ];
568                             }
569                             z_lip /= molecules.block(mol_id).size();
570                             if (z_lip < mem_p->zmed)
571                             {
572                                 nlower++;
573                             }
574                             else
575                             {
576                                 nupper++;
577                             }
578                         }
579                     }
580                 }
581             }
582         }
583     }
584
585     /*make sure equal number of lipids from upper and lower layer are removed */
586     if ( (nupper != nlower) && (!bALLOW_ASYMMETRY) )
587     {
588         snew(dist, mem_p->nmol);
589         snew(order, mem_p->nmol);
590         for (i = 0; i < mem_p->nmol; i++)
591         {
592             at = molecules.block(mem_p->mol_id[i]).begin();
593             pbc_dx(pbc, r[at], pos_ins->geom_cent[0], dr);
594             if (pos_ins->pieces > 1)
595             {
596                 /*minimum dr value*/
597                 min_norm = norm2(dr);
598                 for (k = 1; k < pos_ins->pieces; k++)
599                 {
600                     pbc_dx(pbc, r[at], pos_ins->geom_cent[k], dr_tmp);
601                     if (norm2(dr_tmp) < min_norm)
602                     {
603                         min_norm = norm2(dr_tmp);
604                         copy_rvec(dr_tmp, dr);
605                     }
606                 }
607             }
608             dist[i] = dr[XX]*dr[XX]+dr[YY]*dr[YY];
609             j       = i-1;
610             while (j >= 0 && dist[i] < dist[order[j]])
611             {
612                 order[j+1] = order[j];
613                 j--;
614             }
615             order[j+1] = i;
616         }
617
618         i = 0;
619         while (nupper != nlower)
620         {
621             mol_id = mem_p->mol_id[order[i]];
622             block  = get_molblock(mol_id, mtop->molblock);
623             bRM    = TRUE;
624             for (l = 0; l < nrm; l++)
625             {
626                 if (rm_p->mol[l] == mol_id)
627                 {
628                     bRM = FALSE;
629                 }
630             }
631
632             if (bRM)
633             {
634                 z_lip = 0;
635                 for (int k : molecules.block(mol_id))
636                 {
637                     z_lip += r[k][ZZ];
638                 }
639                 z_lip /= molecules.block(mol_id).size();
640                 if (nupper > nlower && z_lip < mem_p->zmed)
641                 {
642                     rm_p->mol[nrm]   = mol_id;
643                     rm_p->block[nrm] = block;
644                     nrm++;
645                     nlower++;
646                 }
647                 else if (nupper < nlower && z_lip > mem_p->zmed)
648                 {
649                     rm_p->mol[nrm]   = mol_id;
650                     rm_p->block[nrm] = block;
651                     nrm++;
652                     nupper++;
653                 }
654             }
655             i++;
656
657             if (i > mem_p->nmol)
658             {
659                 gmx_fatal(FARGS, "Trying to remove more lipid molecules than there are in the membrane");
660             }
661         }
662         sfree(dist);
663         sfree(order);
664     }
665
666     rm_p->nr = nrm;
667     srenew(rm_p->mol, nrm);
668     srenew(rm_p->block, nrm);
669
670     return nupper+nlower;
671 }
672
673 /*remove all lipids and waters overlapping and update all important structures (e.g. state and mtop)*/
674 static void rm_group(SimulationGroups *groups, gmx_mtop_t *mtop, rm_t *rm_p, t_state *state,
675                      t_block *ins_at, pos_ins_t *pos_ins)
676 {
677     int             j, k, n, rm, mol_id, at, block;
678     rvec           *x_tmp, *v_tmp;
679     int            *list;
680     gmx::EnumerationArray < SimulationAtomGroupType, std::vector < unsigned char>> new_egrp;
681     gmx_bool        bRM;
682     int             RMmolblock;
683
684     /* Construct the molecule range information */
685     gmx::RangePartitioning molecules = gmx_mtop_molecules(*mtop);
686
687     snew(list, state->natoms);
688     n = 0;
689     for (int i = 0; i < rm_p->nr; i++)
690     {
691         mol_id = rm_p->mol[i];
692         at     = molecules.block(mol_id).size();
693         block  = rm_p->block[i];
694         mtop->molblock[block].nmol--;
695         for (j = 0; j < mtop->moltype[mtop->molblock[block].type].atoms.nr; j++)
696         {
697             list[n] = at+j;
698             n++;
699         }
700     }
701
702     mtop->natoms    -= n;
703     state_change_natoms(state, state->natoms - n);
704     snew(x_tmp, state->natoms);
705     snew(v_tmp, state->natoms);
706
707     for (auto group : keysOf(groups->groupNumbers))
708     {
709         if (!groups->groupNumbers[group].empty())
710         {
711             groups->groupNumbers[group].resize(state->natoms);
712             new_egrp[group].resize(state->natoms);
713         }
714     }
715
716     auto x = makeArrayRef(state->x);
717     auto v = makeArrayRef(state->v);
718     rm = 0;
719     for (int i = 0; i < state->natoms+n; i++)
720     {
721         bRM = FALSE;
722         for (j = 0; j < n; j++)
723         {
724             if (i == list[j])
725             {
726                 bRM = TRUE;
727                 rm++;
728             }
729         }
730
731         if (!bRM)
732         {
733             for (auto group : keysOf(groups->groupNumbers))
734             {
735                 if (!groups->groupNumbers[group].empty())
736                 {
737                     new_egrp[group][i-rm] = groups->groupNumbers[group][i];
738                 }
739             }
740             copy_rvec(x[i], x_tmp[i-rm]);
741             copy_rvec(v[i], v_tmp[i-rm]);
742             for (j = 0; j < ins_at->nr; j++)
743             {
744                 if (i == ins_at->index[j])
745                 {
746                     ins_at->index[j] = i-rm;
747                 }
748             }
749
750             for (j = 0; j < pos_ins->pieces; j++)
751             {
752                 for (k = 0; k < pos_ins->nidx[j]; k++)
753                 {
754                     if (i == pos_ins->subindex[j][k])
755                     {
756                         pos_ins->subindex[j][k] = i-rm;
757                     }
758                 }
759             }
760         }
761     }
762     for (int i = 0; i < state->natoms; i++)
763     {
764         copy_rvec(x_tmp[i], x[i]);
765     }
766     sfree(x_tmp);
767     for (int i = 0; i < state->natoms; i++)
768     {
769         copy_rvec(v_tmp[i], v[i]);
770     }
771     sfree(v_tmp);
772
773     for (auto group : keysOf(groups->groupNumbers))
774     {
775         if (!groups->groupNumbers[group].empty())
776         {
777             groups->groupNumbers[group] = new_egrp[group];
778         }
779     }
780
781     /* remove empty molblocks */
782     RMmolblock = 0;
783     for (size_t i = 0; i < mtop->molblock.size(); i++)
784     {
785         if (mtop->molblock[i].nmol == 0)
786         {
787             RMmolblock++;
788         }
789         else
790         {
791             mtop->molblock[i-RMmolblock] = mtop->molblock[i];
792         }
793     }
794     mtop->molblock.resize(mtop->molblock.size() - RMmolblock);
795 }
796
797 /* remove al bonded interactions from mtop for the molecule to be embedded */
798 static int rm_bonded(t_block *ins_at, gmx_mtop_t *mtop)
799 {
800     int       j, m;
801     int       type, natom, nmol, at, atom1 = 0, rm_at = 0;
802     gmx_bool *bRM, bINS;
803     /*this routine lives dangerously by assuming that all molecules of a given type are in order in the structure*/
804     /*this routine does not live as dangerously as it seems. There is namely a check in init_membed to make *
805      * sure that g_membed exits with a warning when there are molecules of the same type not in the *
806      * ins_at index group. MGWolf 050710 */
807
808
809     snew(bRM, mtop->moltype.size());
810     for (size_t i = 0; i < mtop->moltype.size(); i++)
811     {
812         bRM[i] = TRUE;
813     }
814
815     for (size_t i = 0; i < mtop->molblock.size(); i++)
816     {
817         /*loop over molecule blocks*/
818         type         = mtop->molblock[i].type;
819         natom        = mtop->moltype[type].atoms.nr;
820         nmol         = mtop->molblock[i].nmol;
821
822         for (j = 0; j < natom*nmol && bRM[type]; j++)
823         {
824             /*loop over atoms in the block*/
825             at   = j+atom1; /*atom index = block index + offset*/
826             bINS = FALSE;
827
828             for (m = 0; (m < ins_at->nr) && (!bINS); m++)
829             {
830                 /*loop over atoms in insertion index group to determine if we're inserting one*/
831                 if (at == ins_at->index[m])
832                 {
833                     bINS = TRUE;
834                 }
835             }
836             bRM[type] = bINS;
837         }
838         atom1 += natom*nmol; /*update offset*/
839         if (bRM[type])
840         {
841             rm_at += natom*nmol; /*increment bonded removal counter by # atoms in block*/
842         }
843     }
844
845     for (size_t i = 0; i < mtop->moltype.size(); i++)
846     {
847         if (bRM[i])
848         {
849             for (j = 0; j < F_LJ; j++)
850             {
851                 mtop->moltype[i].ilist[j].iatoms.clear();
852             }
853
854             for (j = F_POSRES; j <= F_VSITEN; j++)
855             {
856                 mtop->moltype[i].ilist[j].iatoms.clear();
857             }
858         }
859     }
860     sfree(bRM);
861
862     return rm_at;
863 }
864
865 /* Write a topology where the number of molecules is correct for the system after embedding */
866 static void top_update(const char *topfile, rm_t *rm_p, gmx_mtop_t *mtop)
867 {
868     int        bMolecules         = 0;
869     FILE      *fpin, *fpout;
870     char       buf[STRLEN], buf2[STRLEN], *temp;
871     int        i, *nmol_rm, nmol, line;
872     char       temporary_filename[STRLEN];
873
874     fpin  = gmx_ffopen(topfile, "r");
875     strncpy(temporary_filename, "temp.topXXXXXX", STRLEN);
876     gmx_tmpnam(temporary_filename);
877     fpout = gmx_ffopen(temporary_filename, "w");
878
879     snew(nmol_rm, mtop->moltype.size());
880     for (i = 0; i < rm_p->nr; i++)
881     {
882         nmol_rm[rm_p->block[i]]++;
883     }
884
885     line = 0;
886     while (fgets(buf, STRLEN, fpin))
887     {
888         line++;
889         if (buf[0] != ';')
890         {
891             strcpy(buf2, buf);
892             if ((temp = strchr(buf2, '\n')) != nullptr)
893             {
894                 temp[0] = '\0';
895             }
896             ltrim(buf2);
897             if (buf2[0] == '[')
898             {
899                 buf2[0] = ' ';
900                 if ((temp = strchr(buf2, '\n')) != nullptr)
901                 {
902                     temp[0] = '\0';
903                 }
904                 rtrim(buf2);
905                 if (buf2[strlen(buf2)-1] == ']')
906                 {
907                     buf2[strlen(buf2)-1] = '\0';
908                     ltrim(buf2);
909                     rtrim(buf2);
910                     if (gmx_strcasecmp(buf2, "molecules") == 0)
911                     {
912                         bMolecules = 1;
913                     }
914                 }
915                 fprintf(fpout, "%s", buf);
916             }
917             else if (bMolecules == 1)
918             {
919                 for (size_t i = 0; i < mtop->molblock.size(); i++)
920                 {
921                     nmol = mtop->molblock[i].nmol;
922                     sprintf(buf, "%-15s %5d\n", *(mtop->moltype[mtop->molblock[i].type].name), nmol);
923                     fprintf(fpout, "%s", buf);
924                 }
925                 bMolecules = 2;
926             }
927             else if (bMolecules == 2)
928             {
929                 /* print nothing */
930             }
931             else
932             {
933                 fprintf(fpout, "%s", buf);
934             }
935         }
936         else
937         {
938             fprintf(fpout, "%s", buf);
939         }
940     }
941
942     gmx_ffclose(fpout);
943     /* use gmx_ffopen to generate backup of topinout */
944     fpout = gmx_ffopen(topfile, "w");
945     gmx_ffclose(fpout);
946     rename(temporary_filename, topfile);
947 }
948
949 void rescale_membed(int step_rel, gmx_membed_t *membed, rvec *x)
950 {
951     /* Set new positions for the group to embed */
952     if (step_rel <= membed->it_xy)
953     {
954         membed->fac[0] += membed->xy_step;
955         membed->fac[1] += membed->xy_step;
956     }
957     else if (step_rel <= (membed->it_xy+membed->it_z))
958     {
959         membed->fac[2] += membed->z_step;
960     }
961     resize(membed->r_ins, x, membed->pos_ins, membed->fac);
962 }
963
964 gmx_membed_t *init_membed(FILE *fplog, int nfile, const t_filenm fnm[], gmx_mtop_t *mtop,
965                           t_inputrec *inputrec, t_state *state, t_commrec *cr, real *cpt)
966 {
967     char                            *ins;
968     int                              i, rm_bonded_at, fr_id, fr_i = 0, tmp_id, warn = 0;
969     int                              ng, j, max_lip_rm, ins_grp_id, ntype, lip_rm;
970     real                             prot_area;
971     rvec                            *r_ins = nullptr;
972     t_block                         *ins_at, *rest_at;
973     pos_ins_t                       *pos_ins;
974     mem_t                           *mem_p;
975     rm_t                            *rm_p;
976     SimulationGroups                *groups;
977     gmx_bool                         bExcl = FALSE;
978     t_atoms                          atoms;
979     t_pbc                           *pbc;
980     char                           **piecename = nullptr;
981     gmx_membed_t                    *membed    = nullptr;
982
983     /* input variables */
984     real        xy_fac           = 0.5;
985     real        xy_max           = 1.0;
986     real        z_fac            = 1.0;
987     real        z_max            = 1.0;
988     int         it_xy            = 1000;
989     int         it_z             = 0;
990     real        probe_rad        = 0.22;
991     int         low_up_rm        = 0;
992     int         maxwarn          = 0;
993     int         pieces           = 1;
994     gmx_bool    bALLOW_ASYMMETRY = FALSE;
995
996     /* sanity check constants */         /* Issue a warning when: */
997     const real min_probe_rad  = 0.2199999; /* A probe radius for overlap between embedded molecule *
998                                             * and rest smaller than this value is probably too small */
999     const real min_xy_init    = 0.0999999; /* the initial shrinking of the molecule to embed is smaller */
1000     const int  min_it_xy      = 1000;      /* the number of steps to embed in xy-plane is smaller */
1001     const int  min_it_z       = 100;       /* the number of steps to embed in z is smaller */
1002     const real prot_vs_box    = 7.5;       /* molecule to embed is large (more then prot_vs_box) with respect */
1003     const real box_vs_prot    = 50;        /* to the box size (less than box_vs_prot) */
1004
1005     snew(membed, 1);
1006     snew(ins_at, 1);
1007     snew(pos_ins, 1);
1008
1009     if (MASTER(cr))
1010     {
1011         fprintf(fplog, "Note: it is expected that in future gmx mdrun -membed will not be the "
1012                 "way to access this feature, perhaps changing to e.g. gmx membed.");
1013         /* get input data out membed file */
1014         try
1015         {
1016             get_input(opt2fn("-membed", nfile, fnm),
1017                       &xy_fac, &xy_max, &z_fac, &z_max, &it_xy, &it_z, &probe_rad, &low_up_rm,
1018                       &maxwarn, &pieces, &bALLOW_ASYMMETRY);
1019         }
1020         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1021
1022         if (!EI_DYNAMICS(inputrec->eI) )
1023         {
1024             gmx_input("Change integrator to a dynamics integrator in mdp file (e.g. md or sd).");
1025         }
1026
1027         if (PAR(cr))
1028         {
1029             gmx_input("Sorry, parallel membed is not yet fully functional.");
1030         }
1031
1032         if (*cpt >= 0)
1033         {
1034             fprintf(stderr, "\nSetting -cpt to -1, because embedding cannot be restarted from cpt-files.\n");
1035             *cpt = -1;
1036         }
1037         groups = &(mtop->groups);
1038         std::vector<std::string> gnames;
1039         for (const auto &groupName : groups->groupNames)
1040         {
1041             gnames.emplace_back(*groupName);
1042         }
1043
1044         atoms = gmx_mtop_global_atoms(mtop);
1045         snew(mem_p, 1);
1046         fprintf(stderr, "\nSelect a group to embed in the membrane:\n");
1047         get_index(&atoms, opt2fn_null("-mn", nfile, fnm), 1, &(ins_at->nr), &(ins_at->index), &ins);
1048
1049         auto found = std::find_if(gnames.begin(), gnames.end(),
1050                                   [&ins](const auto &name)
1051                                   { return gmx::equalCaseInsensitive(ins, name); });
1052
1053         if (found == gnames.end())
1054         {
1055             gmx_fatal(FARGS,
1056                       "Group %s selected for embedding was not found in the index file.\n"
1057                       "Group names must match either [moleculetype] names or custom index group\n"
1058                       "names, in which case you must supply an index file to the '-n' option\n"
1059                       "of grompp.", ins);
1060         }
1061         else
1062         {
1063             ins_grp_id = std::distance(gnames.begin(), found);
1064         }
1065         fprintf(stderr, "\nSelect a group to embed %s into (e.g. the membrane):\n", ins);
1066         get_index(&atoms, opt2fn_null("-mn", nfile, fnm), 1, &(mem_p->mem_at.nr), &(mem_p->mem_at.index), &(mem_p->name));
1067
1068         pos_ins->pieces = pieces;
1069         snew(pos_ins->nidx, pieces);
1070         snew(pos_ins->subindex, pieces);
1071         snew(piecename, pieces);
1072         if (pieces > 1)
1073         {
1074             fprintf(stderr, "\nSelect pieces to embed:\n");
1075             get_index(&atoms, opt2fn_null("-mn", nfile, fnm), pieces, pos_ins->nidx, pos_ins->subindex, piecename);
1076         }
1077         else
1078         {
1079             /*use whole embedded group*/
1080             snew(pos_ins->nidx, 1);
1081             snew(pos_ins->subindex, 1);
1082             pos_ins->nidx[0]     = ins_at->nr;
1083             pos_ins->subindex[0] = ins_at->index;
1084         }
1085
1086         if (probe_rad < min_probe_rad)
1087         {
1088             warn++;
1089             fprintf(stderr, "\nWarning %d:\nA probe radius (-rad) smaller than 0.2 nm can result "
1090                     "in overlap between waters and the group to embed, which will result "
1091                     "in Lincs errors etc.\n\n", warn);
1092         }
1093
1094         if (xy_fac < min_xy_init)
1095         {
1096             warn++;
1097             fprintf(stderr, "\nWarning %d:\nThe initial size of %s is probably too small.\n\n", warn, ins);
1098         }
1099
1100         if (it_xy < min_it_xy)
1101         {
1102             warn++;
1103             fprintf(stderr, "\nWarning %d;\nThe number of steps used to grow the xy-coordinates of %s (%d)"
1104                     " is probably too small.\nIncrease -nxy or.\n\n", warn, ins, it_xy);
1105         }
1106
1107         if ( (it_z < min_it_z) && ( z_fac < 0.99999999 || z_fac > 1.0000001) )
1108         {
1109             warn++;
1110             fprintf(stderr, "\nWarning %d;\nThe number of steps used to grow the z-coordinate of %s (%d)"
1111                     " is probably too small.\nIncrease -nz or the maxwarn setting in the membed input file.\n\n", warn, ins, it_z);
1112         }
1113
1114         if (it_xy+it_z > inputrec->nsteps)
1115         {
1116             warn++;
1117             fprintf(stderr, "\nWarning %d:\nThe number of growth steps (-nxy + -nz) is larger than the "
1118                     "number of steps in the tpr.\n"
1119                     "(increase maxwarn in the membed input file to override)\n\n", warn);
1120         }
1121
1122         fr_id = -1;
1123         if (inputrec->opts.ngfrz == 1)
1124         {
1125             gmx_fatal(FARGS, "You did not specify \"%s\" as a freezegroup.", ins);
1126         }
1127
1128         for (i = 0; i < inputrec->opts.ngfrz; i++)
1129         {
1130             tmp_id = mtop->groups.groups[SimulationAtomGroupType::Freeze][i];
1131             if (ins_grp_id == tmp_id)
1132             {
1133                 fr_id = tmp_id;
1134                 fr_i  = i;
1135             }
1136         }
1137
1138         if (fr_id == -1)
1139         {
1140             gmx_fatal(FARGS, "\"%s\" not as freezegroup defined in the mdp-file.", ins);
1141         }
1142
1143         for (i = 0; i < DIM; i++)
1144         {
1145             if (inputrec->opts.nFreeze[fr_i][i] != 1)
1146             {
1147                 gmx_fatal(FARGS, "freeze dimensions for %s are not Y Y Y\n", ins);
1148             }
1149         }
1150
1151         ng = groups->groups[SimulationAtomGroupType::EnergyOutput].size();
1152         if (ng == 1)
1153         {
1154             gmx_input("No energy groups defined. This is necessary for energy exclusion in the freeze group");
1155         }
1156
1157         for (i = 0; i < ng; i++)
1158         {
1159             for (j = 0; j < ng; j++)
1160             {
1161                 if (inputrec->opts.egp_flags[ng*i+j] == EGP_EXCL)
1162                 {
1163                     bExcl = TRUE;
1164                     if ( (groups->groups[SimulationAtomGroupType::EnergyOutput][i] != ins_grp_id) ||
1165                          (groups->groups[SimulationAtomGroupType::EnergyOutput][j] != ins_grp_id) )
1166                     {
1167                         gmx_fatal(FARGS, "Energy exclusions \"%s\" and  \"%s\" do not match the group "
1168                                   "to embed \"%s\"", *groups->groupNames[groups->groups[SimulationAtomGroupType::EnergyOutput][i]],
1169                                   *groups->groupNames[groups->groups[SimulationAtomGroupType::EnergyOutput][j]], ins);
1170                     }
1171                 }
1172             }
1173         }
1174
1175         if (!bExcl)
1176         {
1177             gmx_input("No energy exclusion groups defined. This is necessary for energy exclusion in "
1178                       "the freeze group");
1179         }
1180
1181         /* Obtain the maximum and minimum coordinates of the group to be embedded */
1182         snew(rest_at, 1);
1183         init_ins_at(ins_at, rest_at, state, pos_ins, groups, ins_grp_id, xy_max);
1184         /* Check that moleculetypes in insertion group are not part of the rest of the system */
1185         check_types(ins_at, rest_at, mtop);
1186
1187         init_mem_at(mem_p, mtop, state->x.rvec_array(), state->box, pos_ins);
1188
1189         prot_area = est_prot_area(pos_ins, state->x.rvec_array(), ins_at, mem_p);
1190         if ( (prot_area > prot_vs_box) && ( (state->box[XX][XX]*state->box[YY][YY]-state->box[XX][YY]*state->box[YY][XX]) < box_vs_prot) )
1191         {
1192             warn++;
1193             fprintf(stderr, "\nWarning %d:\nThe xy-area is very small compared to the area of the protein.\n"
1194                     "This might cause pressure problems during the growth phase. Just try with\n"
1195                     "current setup and increase 'maxwarn' in your membed settings file, but lower the\n"
1196                     "compressibility in the mdp-file or disable pressure coupling if problems occur.\n\n", warn);
1197         }
1198
1199         if (warn > maxwarn)
1200         {
1201             gmx_fatal(FARGS, "Too many warnings (override by setting maxwarn in the membed input file)\n");
1202         }
1203
1204         printf("The estimated area of the protein in the membrane is %.3f nm^2\n", prot_area);
1205         printf("\nThere are %d lipids in the membrane part that overlaps the protein.\n"
1206                "The area per lipid is %.4f nm^2.\n", mem_p->nmol, mem_p->lip_area);
1207
1208         /* Maximum number of lipids to be removed*/
1209         max_lip_rm = static_cast<int>(2*prot_area/mem_p->lip_area);
1210         printf("Maximum number of lipids that will be removed is %d.\n", max_lip_rm);
1211
1212         printf("\nWill resize the protein by a factor of %.3f in the xy plane and %.3f in the z direction.\n"
1213                "This resizing will be done with respect to the geometrical center of all protein atoms\n"
1214                "that span the membrane region, i.e. z between %.3f and %.3f\n\n",
1215                xy_fac, z_fac, mem_p->zmin, mem_p->zmax);
1216
1217         /* resize the protein by xy and by z if necessary*/
1218         snew(r_ins, ins_at->nr);
1219         init_resize(ins_at, r_ins, pos_ins, mem_p, state->x.rvec_array(), bALLOW_ASYMMETRY);
1220         membed->fac[0] = membed->fac[1] = xy_fac;
1221         membed->fac[2] = z_fac;
1222
1223         membed->xy_step = (xy_max-xy_fac)/static_cast<double>(it_xy);
1224         membed->z_step  = (z_max-z_fac)/static_cast<double>(it_z-1);
1225
1226         resize(r_ins, state->x.rvec_array(), pos_ins, membed->fac);
1227
1228         /* remove overlapping lipids and water from the membrane box*/
1229         /*mark molecules to be removed*/
1230         snew(pbc, 1);
1231         set_pbc(pbc, inputrec->ePBC, state->box);
1232
1233         snew(rm_p, 1);
1234         lip_rm = gen_rm_list(rm_p, ins_at, rest_at, pbc, mtop, state->x.rvec_array(), mem_p, pos_ins,
1235                              probe_rad, low_up_rm, bALLOW_ASYMMETRY);
1236         lip_rm -= low_up_rm;
1237
1238         if (fplog)
1239         {
1240             for (i = 0; i < rm_p->nr; i++)
1241             {
1242                 fprintf(fplog, "rm mol %d\n", rm_p->mol[i]);
1243             }
1244         }
1245
1246         for (size_t i = 0; i < mtop->molblock.size(); i++)
1247         {
1248             ntype = 0;
1249             for (j = 0; j < rm_p->nr; j++)
1250             {
1251                 if (rm_p->block[j] == static_cast<int>(i))
1252                 {
1253                     ntype++;
1254                 }
1255             }
1256             printf("Will remove %d %s molecules\n", ntype, *(mtop->moltype[mtop->molblock[i].type].name));
1257         }
1258
1259         if (lip_rm > max_lip_rm)
1260         {
1261             warn++;
1262             fprintf(stderr, "\nWarning %d:\nTrying to remove a larger lipid area than the estimated "
1263                     "protein area\nTry making the -xyinit resize factor smaller or increase "
1264                     "maxwarn in the membed input file.\n\n", warn);
1265         }
1266
1267         /*remove all lipids and waters overlapping and update all important structures*/
1268         rm_group(groups, mtop, rm_p, state, ins_at, pos_ins);
1269
1270         rm_bonded_at = rm_bonded(ins_at, mtop);
1271         if (rm_bonded_at != ins_at->nr)
1272         {
1273             fprintf(stderr, "Warning: The number of atoms for which the bonded interactions are removed is %d, "
1274                     "while %d atoms are embedded. Make sure that the atoms to be embedded are not in the same"
1275                     "molecule type as atoms that are not to be embedded.\n", rm_bonded_at, ins_at->nr);
1276         }
1277
1278         if (warn > maxwarn)
1279         {
1280             gmx_fatal(FARGS, "Too many warnings.\nIf you are sure these warnings are harmless,\n"
1281                       "you can increase the maxwarn setting in the membed input file.");
1282         }
1283
1284         // Re-establish the invariants of the derived values within
1285         // mtop.
1286         gmx_mtop_finalize(mtop);
1287
1288         if (ftp2bSet(efTOP, nfile, fnm))
1289         {
1290             top_update(opt2fn("-mp", nfile, fnm), rm_p, mtop);
1291         }
1292
1293         sfree(pbc);
1294         sfree(rest_at);
1295         if (pieces > 1)
1296         {
1297             sfree(piecename);
1298         }
1299
1300         membed->it_xy   = it_xy;
1301         membed->it_z    = it_z;
1302         membed->pos_ins = pos_ins;
1303         membed->r_ins   = r_ins;
1304     }
1305
1306     return membed;
1307 }
1308
1309 void free_membed(gmx_membed_t *membed)
1310 {
1311     sfree(membed);
1312 }