2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2010,2011,2012,2013,2014,2015,2016,2017, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
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"
63 /* information about scaling center */
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]] */
74 /* variables needed in do_md */
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 */
85 /* membrane related variables */
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 */
97 /* Lists all molecules in the membrane that overlap with the molecule to be embedded. *
98 * These will then be removed from the system */
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 */
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
)
114 mtopGetMolblockIndex(mtop
, at
, block
, &mol_id
, &atnr_mol
);
115 for (i
= 0; i
< *block
; i
++)
117 mol_id
+= mtop
->molblock
[i
].nmol
;
119 *type
= mtop
->molblock
[*block
].type
;
124 /* Get the id of the molblock from a global molecule id */
125 static int get_molblock(int mol_id
, int nmblock
, gmx_molblock_t
*mblock
)
130 for (i
= 0; i
< nmblock
; i
++)
132 nmol
+= mblock
[i
].nmol
;
139 gmx_fatal(FARGS
, "mol_id %d larger than total number of molecules %d.\n", mol_id
, nmol
);
144 /* Get a list of all the molecule types that are present in a group of atoms. *
145 * Because all interaction within the group to embed are removed on the topology *
146 * level, if the same molecule type is found in another part of the system, these *
147 * would also be affected. Therefore we have to check if the embedded and rest group *
148 * share common molecule types. If so, membed will stop with an error. */
149 static int get_mtype_list(t_block
*at
, gmx_mtop_t
*mtop
, t_block
*tlist
)
152 int type
= 0, block
= 0;
156 snew(tlist
->index
, at
->nr
);
157 for (i
= 0; i
< at
->nr
; i
++)
160 get_mol_id(at
->index
[i
], mtop
, &type
, &block
);
161 for (j
= 0; j
< nr
; j
++)
163 if (tlist
->index
[j
] == type
)
171 tlist
->index
[nr
] = type
;
175 srenew(tlist
->index
, nr
);
179 /* Do the actual check of the molecule types between embedded and rest group */
180 static void check_types(t_block
*ins_at
, t_block
*rest_at
, gmx_mtop_t
*mtop
)
182 t_block
*ins_mtype
, *rest_mtype
;
187 ins_mtype
->nr
= get_mtype_list(ins_at
, mtop
, ins_mtype
);
188 rest_mtype
->nr
= get_mtype_list(rest_at
, mtop
, rest_mtype
);
190 for (i
= 0; i
< ins_mtype
->nr
; i
++)
192 for (j
= 0; j
< rest_mtype
->nr
; j
++)
194 if (ins_mtype
->index
[i
] == rest_mtype
->index
[j
])
196 gmx_fatal(FARGS
, "Moleculetype %s is found both in the group to insert and the rest of the system.\n"
197 "1. Your *.ndx and *.top do not match\n"
198 "2. You are inserting some molecules of type %s (for example xray-solvent), while\n"
199 "the same moleculetype is also used in the rest of the system (solvent box). Because\n"
200 "we need to exclude all interactions between the atoms in the group to\n"
201 "insert, the same moleculetype can not be used in both groups. Change the\n"
202 "moleculetype of the molecules %s in the inserted group. Do not forget to provide\n"
203 "an appropriate *.itp file", *(mtop
->moltype
[rest_mtype
->index
[j
]].name
),
204 *(mtop
->moltype
[rest_mtype
->index
[j
]].name
), *(mtop
->moltype
[rest_mtype
->index
[j
]].name
));
209 sfree(ins_mtype
->index
);
210 sfree(rest_mtype
->index
);
215 static void get_input(const char *membed_input
, real
*xy_fac
, real
*xy_max
, real
*z_fac
, real
*z_max
,
216 int *it_xy
, int *it_z
, real
*probe_rad
, int *low_up_rm
, int *maxwarn
,
217 int *pieces
, gmx_bool
*bALLOW_ASYMMETRY
)
223 wi
= init_warning(TRUE
, 0);
226 gmx::TextInputFile
stream(membed_input
);
227 inp
= read_inpfile(&stream
, membed_input
, &ninp
, wi
);
230 ITYPE ("nxy", *it_xy
, 1000);
231 ITYPE ("nz", *it_z
, 0);
232 RTYPE ("xyinit", *xy_fac
, 0.5);
233 RTYPE ("xyend", *xy_max
, 1.0);
234 RTYPE ("zinit", *z_fac
, 1.0);
235 RTYPE ("zend", *z_max
, 1.0);
236 RTYPE ("rad", *probe_rad
, 0.22);
237 ITYPE ("ndiff", *low_up_rm
, 0);
238 ITYPE ("maxwarn", *maxwarn
, 0);
239 ITYPE ("pieces", *pieces
, 1);
240 EETYPE("asymmetry", *bALLOW_ASYMMETRY
, yesno_names
);
242 check_warning_error(wi
, FARGS
);
244 gmx::TextOutputFile
stream(membed_input
);
245 write_inpfile(&stream
, membed_input
, ninp
, inp
, FALSE
, WriteMdpHeader::yes
, wi
);
248 done_warning(wi
, FARGS
);
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 gmx_groups_t
*groups
, int ins_grp_id
, real xy_max
)
256 real x
, xmin
, xmax
, y
, ymin
, ymax
, z
, 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
);
262 xmin
= xmax
= state
->x
[ins_at
->index
[0]][XX
];
263 ymin
= ymax
= state
->x
[ins_at
->index
[0]][YY
];
264 zmin
= zmax
= state
->x
[ins_at
->index
[0]][ZZ
];
266 for (i
= 0; i
< state
->natoms
; i
++)
268 gid
= groups
->grpnr
[egcFREEZE
][i
];
269 if (groups
->grps
[egcFREEZE
].nm_ind
[gid
] == ins_grp_id
)
301 rest_at
->index
[c
] = i
;
307 srenew(rest_at
->index
, c
);
309 if (xy_max
> fac_inp_size
)
311 pos_ins
->xmin
[XX
] = xmin
-((xmax
-xmin
)*xy_max
-(xmax
-xmin
))/2;
312 pos_ins
->xmin
[YY
] = ymin
-((ymax
-ymin
)*xy_max
-(ymax
-ymin
))/2;
314 pos_ins
->xmax
[XX
] = xmax
+((xmax
-xmin
)*xy_max
-(xmax
-xmin
))/2;
315 pos_ins
->xmax
[YY
] = ymax
+((ymax
-ymin
)*xy_max
-(ymax
-ymin
))/2;
319 pos_ins
->xmin
[XX
] = xmin
;
320 pos_ins
->xmin
[YY
] = ymin
;
322 pos_ins
->xmax
[XX
] = xmax
;
323 pos_ins
->xmax
[YY
] = ymax
;
326 if ( (zmax
-zmin
) < min_memthick
)
328 pos_ins
->xmin
[ZZ
] = zmin
+(zmax
-zmin
)/2.0-0.5*min_memthick
;
329 pos_ins
->xmax
[ZZ
] = zmin
+(zmax
-zmin
)/2.0+0.5*min_memthick
;
333 pos_ins
->xmin
[ZZ
] = zmin
;
334 pos_ins
->xmax
[ZZ
] = zmax
;
340 /* Estimate the area of the embedded molecule by projecting all coordinates on a grid in the *
341 * xy-plane and counting the number of occupied grid points */
342 static real
est_prot_area(pos_ins_t
*pos_ins
, rvec
*r
, t_block
*ins_at
, mem_t
*mem_p
)
344 real x
, y
, dx
= 0.15, dy
= 0.15, area
= 0.0;
345 real add
, memmin
, memmax
;
348 /* min and max membrane coordinate are altered to reduce the influence of the *
350 memmin
= mem_p
->zmin
+0.1*(mem_p
->zmax
-mem_p
->zmin
);
351 memmax
= mem_p
->zmax
-0.1*(mem_p
->zmax
-mem_p
->zmin
);
353 for (x
= pos_ins
->xmin
[XX
]; x
< pos_ins
->xmax
[XX
]; x
+= dx
)
355 for (y
= pos_ins
->xmin
[YY
]; y
< pos_ins
->xmax
[YY
]; y
+= dy
)
361 at
= ins_at
->index
[c
];
362 if ( (r
[at
][XX
] >= x
) && (r
[at
][XX
] < x
+dx
) &&
363 (r
[at
][YY
] >= y
) && (r
[at
][YY
] < y
+dy
) &&
364 (r
[at
][ZZ
] > memmin
) && (r
[at
][ZZ
] < memmax
) )
370 while ( (c
< ins_at
->nr
) && (add
< 0.5) );
379 static int init_mem_at(mem_t
*mem_p
, gmx_mtop_t
*mtop
, rvec
*r
, matrix box
, pos_ins_t
*pos_ins
)
381 int i
, j
, at
, mol
, nmol
, nmolbox
, count
;
383 real z
, zmin
, zmax
, mem_area
;
386 int type
= 0, block
= 0;
389 mem_a
= &(mem_p
->mem_at
);
390 snew(mol_id
, mem_a
->nr
);
391 zmin
= pos_ins
->xmax
[ZZ
];
392 zmax
= pos_ins
->xmin
[ZZ
];
393 for (i
= 0; i
< mem_a
->nr
; i
++)
395 at
= mem_a
->index
[i
];
396 if ( (r
[at
][XX
] > pos_ins
->xmin
[XX
]) && (r
[at
][XX
] < pos_ins
->xmax
[XX
]) &&
397 (r
[at
][YY
] > pos_ins
->xmin
[YY
]) && (r
[at
][YY
] < pos_ins
->xmax
[YY
]) &&
398 (r
[at
][ZZ
] > pos_ins
->xmin
[ZZ
]) && (r
[at
][ZZ
] < pos_ins
->xmax
[ZZ
]) )
400 mol
= get_mol_id(at
, mtop
, &type
, &block
);
402 for (j
= 0; j
< nmol
; j
++)
404 if (mol
== mol_id
[j
])
432 srenew(mol_id
, nmol
);
433 mem_p
->mol_id
= mol_id
;
435 if ((zmax
-zmin
) > (box
[ZZ
][ZZ
]-0.5))
437 gmx_fatal(FARGS
, "Something is wrong with your membrane. Max and min z values are %f and %f.\n"
438 "Maybe your membrane is not centered in the box, but located at the box edge in the z-direction,\n"
439 "so that one membrane is distributed over two periodic box images. Another possibility is that\n"
440 "your water layer is not thick enough.\n", zmax
, zmin
);
444 mem_p
->zmed
= (zmax
-zmin
)/2+zmin
;
446 /*number of membrane molecules in protein box*/
447 nmolbox
= count
/mtop
->molblock
[block
].natoms_mol
;
448 /*membrane area within the box defined by the min and max coordinates of the embedded molecule*/
449 mem_area
= (pos_ins
->xmax
[XX
]-pos_ins
->xmin
[XX
])*(pos_ins
->xmax
[YY
]-pos_ins
->xmin
[YY
]);
450 /*rough estimate of area per lipid, assuming there is only one type of lipid in the membrane*/
451 mem_p
->lip_area
= 2.0*mem_area
/(double)nmolbox
;
453 return mem_p
->mem_at
.nr
;
456 static void init_resize(t_block
*ins_at
, rvec
*r_ins
, pos_ins_t
*pos_ins
, mem_t
*mem_p
, rvec
*r
,
457 gmx_bool bALLOW_ASYMMETRY
)
459 int i
, j
, at
, c
, outsidesum
, gctr
= 0;
463 for (i
= 0; i
< pos_ins
->pieces
; i
++)
465 idxsum
+= pos_ins
->nidx
[i
];
468 if (idxsum
!= ins_at
->nr
)
470 gmx_fatal(FARGS
, "Piecewise sum of inserted atoms not same as size of group selected to insert.");
473 snew(pos_ins
->geom_cent
, pos_ins
->pieces
);
474 for (i
= 0; i
< pos_ins
->pieces
; i
++)
478 for (j
= 0; j
< DIM
; j
++)
480 pos_ins
->geom_cent
[i
][j
] = 0;
483 for (j
= 0; j
< pos_ins
->nidx
[i
]; j
++)
485 at
= pos_ins
->subindex
[i
][j
];
486 copy_rvec(r
[at
], r_ins
[gctr
]);
487 if ( (r_ins
[gctr
][ZZ
] < mem_p
->zmax
) && (r_ins
[gctr
][ZZ
] > mem_p
->zmin
) )
489 rvec_inc(pos_ins
->geom_cent
[i
], r_ins
[gctr
]);
501 svmul(1/(double)c
, pos_ins
->geom_cent
[i
], pos_ins
->geom_cent
[i
]);
504 if (!bALLOW_ASYMMETRY
)
506 pos_ins
->geom_cent
[i
][ZZ
] = mem_p
->zmed
;
509 fprintf(stderr
, "Embedding piece %d with center of geometry: %f %f %f\n",
510 i
, pos_ins
->geom_cent
[i
][XX
], pos_ins
->geom_cent
[i
][YY
], pos_ins
->geom_cent
[i
][ZZ
]);
512 fprintf(stderr
, "\n");
515 /* resize performed in the md loop */
516 static void resize(rvec
*r_ins
, rvec
*r
, pos_ins_t
*pos_ins
, rvec fac
)
518 int i
, j
, k
, at
, c
= 0;
519 for (k
= 0; k
< pos_ins
->pieces
; k
++)
521 for (i
= 0; i
< pos_ins
->nidx
[k
]; i
++)
523 at
= pos_ins
->subindex
[k
][i
];
524 for (j
= 0; j
< DIM
; j
++)
526 r
[at
][j
] = pos_ins
->geom_cent
[k
][j
]+fac
[j
]*(r_ins
[c
][j
]-pos_ins
->geom_cent
[k
][j
]);
533 /* generate the list of membrane molecules that overlap with the molecule to be embedded. *
534 * The molecule to be embedded is already reduced in size. */
535 static int gen_rm_list(rm_t
*rm_p
, t_block
*ins_at
, t_block
*rest_at
, t_pbc
*pbc
, gmx_mtop_t
*mtop
,
536 rvec
*r
, mem_t
*mem_p
, pos_ins_t
*pos_ins
, real probe_rad
,
537 int low_up_rm
, gmx_bool bALLOW_ASYMMETRY
)
539 int i
, j
, k
, l
, at
, at2
, mol_id
;
540 int type
= 0, block
= 0;
541 int nrm
, nupper
, nlower
;
542 real r_min_rad
, z_lip
, min_norm
;
548 r_min_rad
= probe_rad
*probe_rad
;
549 snew(rm_p
->mol
, mtop
->mols
.nr
);
550 snew(rm_p
->block
, mtop
->mols
.nr
);
553 for (i
= 0; i
< ins_at
->nr
; i
++)
555 at
= ins_at
->index
[i
];
556 for (j
= 0; j
< rest_at
->nr
; j
++)
558 at2
= rest_at
->index
[j
];
559 pbc_dx(pbc
, r
[at
], r
[at2
], dr
);
561 if (norm2(dr
) < r_min_rad
)
563 mol_id
= get_mol_id(at2
, mtop
, &type
, &block
);
565 for (l
= 0; l
< nrm
; l
++)
567 if (rm_p
->mol
[l
] == mol_id
)
575 rm_p
->mol
[nrm
] = mol_id
;
576 rm_p
->block
[nrm
] = block
;
579 for (l
= 0; l
< mem_p
->nmol
; l
++)
581 if (mol_id
== mem_p
->mol_id
[l
])
583 for (k
= mtop
->mols
.index
[mol_id
]; k
< mtop
->mols
.index
[mol_id
+1]; k
++)
587 z_lip
/= mtop
->molblock
[block
].natoms_mol
;
588 if (z_lip
< mem_p
->zmed
)
603 /*make sure equal number of lipids from upper and lower layer are removed */
604 if ( (nupper
!= nlower
) && (!bALLOW_ASYMMETRY
) )
606 snew(dist
, mem_p
->nmol
);
607 snew(order
, mem_p
->nmol
);
608 for (i
= 0; i
< mem_p
->nmol
; i
++)
610 at
= mtop
->mols
.index
[mem_p
->mol_id
[i
]];
611 pbc_dx(pbc
, r
[at
], pos_ins
->geom_cent
[0], dr
);
612 if (pos_ins
->pieces
> 1)
615 min_norm
= norm2(dr
);
616 for (k
= 1; k
< pos_ins
->pieces
; k
++)
618 pbc_dx(pbc
, r
[at
], pos_ins
->geom_cent
[k
], dr_tmp
);
619 if (norm2(dr_tmp
) < min_norm
)
621 min_norm
= norm2(dr_tmp
);
622 copy_rvec(dr_tmp
, dr
);
626 dist
[i
] = dr
[XX
]*dr
[XX
]+dr
[YY
]*dr
[YY
];
628 while (j
>= 0 && dist
[i
] < dist
[order
[j
]])
630 order
[j
+1] = order
[j
];
637 while (nupper
!= nlower
)
639 mol_id
= mem_p
->mol_id
[order
[i
]];
640 block
= get_molblock(mol_id
, mtop
->nmolblock
, mtop
->molblock
);
642 for (l
= 0; l
< nrm
; l
++)
644 if (rm_p
->mol
[l
] == mol_id
)
653 for (k
= mtop
->mols
.index
[mol_id
]; k
< mtop
->mols
.index
[mol_id
+1]; k
++)
657 z_lip
/= mtop
->molblock
[block
].natoms_mol
;
658 if (nupper
> nlower
&& z_lip
< mem_p
->zmed
)
660 rm_p
->mol
[nrm
] = mol_id
;
661 rm_p
->block
[nrm
] = block
;
665 else if (nupper
< nlower
&& z_lip
> mem_p
->zmed
)
667 rm_p
->mol
[nrm
] = mol_id
;
668 rm_p
->block
[nrm
] = block
;
677 gmx_fatal(FARGS
, "Trying to remove more lipid molecules than there are in the membrane");
685 srenew(rm_p
->mol
, nrm
);
686 srenew(rm_p
->block
, nrm
);
688 return nupper
+nlower
;
691 /*remove all lipids and waters overlapping and update all important structures (e.g. state and mtop)*/
692 static void rm_group(gmx_groups_t
*groups
, gmx_mtop_t
*mtop
, rm_t
*rm_p
, t_state
*state
,
693 t_block
*ins_at
, pos_ins_t
*pos_ins
)
695 int i
, j
, k
, n
, rm
, mol_id
, at
, block
;
697 int *list
, *new_mols
;
698 unsigned char *new_egrp
[egcNR
];
702 snew(list
, state
->natoms
);
704 for (i
= 0; i
< rm_p
->nr
; i
++)
706 mol_id
= rm_p
->mol
[i
];
707 at
= mtop
->mols
.index
[mol_id
];
708 block
= rm_p
->block
[i
];
709 mtop
->molblock
[block
].nmol
--;
710 for (j
= 0; j
< mtop
->molblock
[block
].natoms_mol
; j
++)
715 mtop
->mols
.index
[mol_id
] = -1;
718 mtop
->mols
.nr
-= rm_p
->nr
;
719 mtop
->mols
.nalloc_index
-= rm_p
->nr
;
720 snew(new_mols
, mtop
->mols
.nr
);
721 for (i
= 0; i
< mtop
->mols
.nr
+rm_p
->nr
; i
++)
724 if (mtop
->mols
.index
[i
] != -1)
726 new_mols
[j
] = mtop
->mols
.index
[i
];
730 sfree(mtop
->mols
.index
);
731 mtop
->mols
.index
= new_mols
;
734 snew(x_tmp
, state
->natoms
);
735 snew(v_tmp
, state
->natoms
);
737 for (i
= 0; i
< egcNR
; i
++)
739 if (groups
->grpnr
[i
] != nullptr)
741 groups
->ngrpnr
[i
] = state
->natoms
;
742 snew(new_egrp
[i
], state
->natoms
);
747 for (i
= 0; i
< state
->natoms
+n
; i
++)
750 for (j
= 0; j
< n
; j
++)
761 for (j
= 0; j
< egcNR
; j
++)
763 if (groups
->grpnr
[j
] != nullptr)
765 new_egrp
[j
][i
-rm
] = groups
->grpnr
[j
][i
];
768 copy_rvec(state
->x
[i
], x_tmp
[i
-rm
]);
769 copy_rvec(state
->v
[i
], v_tmp
[i
-rm
]);
770 for (j
= 0; j
< ins_at
->nr
; j
++)
772 if (i
== ins_at
->index
[j
])
774 ins_at
->index
[j
] = i
-rm
;
778 for (j
= 0; j
< pos_ins
->pieces
; j
++)
780 for (k
= 0; k
< pos_ins
->nidx
[j
]; k
++)
782 if (i
== pos_ins
->subindex
[j
][k
])
784 pos_ins
->subindex
[j
][k
] = i
-rm
;
790 for (int i
= 0; i
< state
->natoms
; i
++)
792 copy_rvec(x_tmp
[i
], state
->x
[i
]);
795 for (int i
= 0; i
< state
->natoms
; i
++)
797 copy_rvec(v_tmp
[i
], state
->v
[i
]);
801 for (i
= 0; i
< egcNR
; i
++)
803 if (groups
->grpnr
[i
] != nullptr)
805 sfree(groups
->grpnr
[i
]);
806 groups
->grpnr
[i
] = new_egrp
[i
];
810 /* remove empty molblocks */
812 for (i
= 0; i
< mtop
->nmolblock
; i
++)
814 if (mtop
->molblock
[i
].nmol
== 0)
820 mtop
->molblock
[i
-RMmolblock
] = mtop
->molblock
[i
];
823 mtop
->nmolblock
-= RMmolblock
;
826 /* remove al bonded interactions from mtop for the molecule to be embedded */
827 int rm_bonded(t_block
*ins_at
, gmx_mtop_t
*mtop
)
830 int type
, natom
, nmol
, at
, atom1
= 0, rm_at
= 0;
832 /*this routine lives dangerously by assuming that all molecules of a given type are in order in the structure*/
833 /*this routine does not live as dangerously as it seems. There is namely a check in init_membed to make *
834 * sure that g_membed exits with a warning when there are molecules of the same type not in the *
835 * ins_at index group. MGWolf 050710 */
838 snew(bRM
, mtop
->nmoltype
);
839 for (i
= 0; i
< mtop
->nmoltype
; i
++)
844 for (i
= 0; i
< mtop
->nmolblock
; i
++)
846 /*loop over molecule blocks*/
847 type
= mtop
->molblock
[i
].type
;
848 natom
= mtop
->molblock
[i
].natoms_mol
;
849 nmol
= mtop
->molblock
[i
].nmol
;
851 for (j
= 0; j
< natom
*nmol
&& bRM
[type
] == TRUE
; j
++)
853 /*loop over atoms in the block*/
854 at
= j
+atom1
; /*atom index = block index + offset*/
857 for (m
= 0; (m
< ins_at
->nr
) && (bINS
== FALSE
); m
++)
859 /*loop over atoms in insertion index group to determine if we're inserting one*/
860 if (at
== ins_at
->index
[m
])
867 atom1
+= natom
*nmol
; /*update offset*/
870 rm_at
+= natom
*nmol
; /*increment bonded removal counter by # atoms in block*/
874 for (i
= 0; i
< mtop
->nmoltype
; i
++)
878 for (j
= 0; j
< F_LJ
; j
++)
880 mtop
->moltype
[i
].ilist
[j
].nr
= 0;
883 for (j
= F_POSRES
; j
<= F_VSITEN
; j
++)
885 mtop
->moltype
[i
].ilist
[j
].nr
= 0;
894 /* Write a topology where the number of molecules is correct for the system after embedding */
895 static void top_update(const char *topfile
, rm_t
*rm_p
, gmx_mtop_t
*mtop
)
899 char buf
[STRLEN
], buf2
[STRLEN
], *temp
;
900 int i
, *nmol_rm
, nmol
, line
;
901 char temporary_filename
[STRLEN
];
903 fpin
= gmx_ffopen(topfile
, "r");
904 strncpy(temporary_filename
, "temp.topXXXXXX", STRLEN
);
905 gmx_tmpnam(temporary_filename
);
906 fpout
= gmx_ffopen(temporary_filename
, "w");
908 snew(nmol_rm
, mtop
->nmoltype
);
909 for (i
= 0; i
< rm_p
->nr
; i
++)
911 nmol_rm
[rm_p
->block
[i
]]++;
915 while (fgets(buf
, STRLEN
, fpin
))
921 if ((temp
= strchr(buf2
, '\n')) != nullptr)
929 if ((temp
= strchr(buf2
, '\n')) != nullptr)
934 if (buf2
[strlen(buf2
)-1] == ']')
936 buf2
[strlen(buf2
)-1] = '\0';
939 if (gmx_strcasecmp(buf2
, "molecules") == 0)
944 fprintf(fpout
, "%s", buf
);
946 else if (bMolecules
== 1)
948 for (i
= 0; i
< mtop
->nmolblock
; i
++)
950 nmol
= mtop
->molblock
[i
].nmol
;
951 sprintf(buf
, "%-15s %5d\n", *(mtop
->moltype
[mtop
->molblock
[i
].type
].name
), nmol
);
952 fprintf(fpout
, "%s", buf
);
956 else if (bMolecules
== 2)
962 fprintf(fpout
, "%s", buf
);
967 fprintf(fpout
, "%s", buf
);
972 /* use gmx_ffopen to generate backup of topinout */
973 fpout
= gmx_ffopen(topfile
, "w");
975 rename(temporary_filename
, topfile
);
978 void rescale_membed(int step_rel
, gmx_membed_t
*membed
, rvec
*x
)
980 /* Set new positions for the group to embed */
981 if (step_rel
<= membed
->it_xy
)
983 membed
->fac
[0] += membed
->xy_step
;
984 membed
->fac
[1] += membed
->xy_step
;
986 else if (step_rel
<= (membed
->it_xy
+membed
->it_z
))
988 membed
->fac
[2] += membed
->z_step
;
990 resize(membed
->r_ins
, x
, membed
->pos_ins
, membed
->fac
);
993 /* We would like gn to be const as well, but C doesn't allow this */
994 /* TODO this is utility functionality (search for the index of a
995 string in a collection), so should be refactored and located more
996 centrally. Also, it nearly duplicates the same string in readir.c */
997 static int search_string(const char *s
, int ng
, char *gn
[])
1001 for (i
= 0; (i
< ng
); i
++)
1003 if (gmx_strcasecmp(s
, gn
[i
]) == 0)
1010 "Group %s selected for embedding was not found in the index file.\n"
1011 "Group names must match either [moleculetype] names or custom index group\n"
1012 "names, in which case you must supply an index file to the '-n' option\n"
1019 gmx_membed_t
*init_membed(FILE *fplog
, int nfile
, const t_filenm fnm
[], gmx_mtop_t
*mtop
,
1020 t_inputrec
*inputrec
, t_state
*state
, t_commrec
*cr
, real
*cpt
)
1022 char *ins
, **gnames
;
1023 int i
, rm_bonded_at
, fr_id
, fr_i
= 0, tmp_id
, warn
= 0;
1024 int ng
, j
, max_lip_rm
, ins_grp_id
, ntype
, lip_rm
;
1026 rvec
*r_ins
= nullptr;
1027 t_block
*ins_at
, *rest_at
;
1031 gmx_groups_t
*groups
;
1032 gmx_bool bExcl
= FALSE
;
1035 char **piecename
= nullptr;
1036 gmx_membed_t
*membed
= nullptr;
1038 /* input variables */
1045 real probe_rad
= 0.22;
1049 gmx_bool bALLOW_ASYMMETRY
= FALSE
;
1051 /* sanity check constants */ /* Issue a warning when: */
1052 const real min_probe_rad
= 0.2199999; /* A probe radius for overlap between embedded molecule *
1053 * and rest smaller than this value is probably too small */
1054 const real min_xy_init
= 0.0999999; /* the initial shrinking of the molecule to embed is smaller */
1055 const int min_it_xy
= 1000; /* the number of steps to embed in xy-plane is smaller */
1056 const int min_it_z
= 100; /* the number of steps to embed in z is smaller */
1057 const real prot_vs_box
= 7.5; /* molecule to embed is large (more then prot_vs_box) with respect */
1058 const real box_vs_prot
= 50; /* to the box size (less than box_vs_prot) */
1066 /* get input data out membed file */
1069 get_input(opt2fn("-membed", nfile
, fnm
),
1070 &xy_fac
, &xy_max
, &z_fac
, &z_max
, &it_xy
, &it_z
, &probe_rad
, &low_up_rm
,
1071 &maxwarn
, &pieces
, &bALLOW_ASYMMETRY
);
1073 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
;
1075 if (!EI_DYNAMICS(inputrec
->eI
) )
1077 gmx_input("Change integrator to a dynamics integrator in mdp file (e.g. md or sd).");
1082 gmx_input("Sorry, parallel g_membed is not yet fully functional.");
1087 fprintf(stderr
, "\nSetting -cpt to -1, because embedding cannot be restarted from cpt-files.\n");
1090 groups
= &(mtop
->groups
);
1091 snew(gnames
, groups
->ngrpname
);
1092 for (i
= 0; i
< groups
->ngrpname
; i
++)
1094 gnames
[i
] = *(groups
->grpname
[i
]);
1097 atoms
= gmx_mtop_global_atoms(mtop
);
1099 fprintf(stderr
, "\nSelect a group to embed in the membrane:\n");
1100 get_index(&atoms
, opt2fn_null("-mn", nfile
, fnm
), 1, &(ins_at
->nr
), &(ins_at
->index
), &ins
);
1101 ins_grp_id
= search_string(ins
, groups
->ngrpname
, gnames
);
1102 fprintf(stderr
, "\nSelect a group to embed %s into (e.g. the membrane):\n", ins
);
1103 get_index(&atoms
, opt2fn_null("-mn", nfile
, fnm
), 1, &(mem_p
->mem_at
.nr
), &(mem_p
->mem_at
.index
), &(mem_p
->name
));
1105 pos_ins
->pieces
= pieces
;
1106 snew(pos_ins
->nidx
, pieces
);
1107 snew(pos_ins
->subindex
, pieces
);
1108 snew(piecename
, pieces
);
1111 fprintf(stderr
, "\nSelect pieces to embed:\n");
1112 get_index(&atoms
, opt2fn_null("-mn", nfile
, fnm
), pieces
, pos_ins
->nidx
, pos_ins
->subindex
, piecename
);
1116 /*use whole embedded group*/
1117 snew(pos_ins
->nidx
, 1);
1118 snew(pos_ins
->subindex
, 1);
1119 pos_ins
->nidx
[0] = ins_at
->nr
;
1120 pos_ins
->subindex
[0] = ins_at
->index
;
1123 if (probe_rad
< min_probe_rad
)
1126 fprintf(stderr
, "\nWarning %d:\nA probe radius (-rad) smaller than 0.2 nm can result "
1127 "in overlap between waters and the group to embed, which will result "
1128 "in Lincs errors etc.\n\n", warn
);
1131 if (xy_fac
< min_xy_init
)
1134 fprintf(stderr
, "\nWarning %d:\nThe initial size of %s is probably too small.\n\n", warn
, ins
);
1137 if (it_xy
< min_it_xy
)
1140 fprintf(stderr
, "\nWarning %d;\nThe number of steps used to grow the xy-coordinates of %s (%d)"
1141 " is probably too small.\nIncrease -nxy or.\n\n", warn
, ins
, it_xy
);
1144 if ( (it_z
< min_it_z
) && ( z_fac
< 0.99999999 || z_fac
> 1.0000001) )
1147 fprintf(stderr
, "\nWarning %d;\nThe number of steps used to grow the z-coordinate of %s (%d)"
1148 " is probably too small.\nIncrease -nz or the maxwarn setting in the membed input file.\n\n", warn
, ins
, it_z
);
1151 if (it_xy
+it_z
> inputrec
->nsteps
)
1154 fprintf(stderr
, "\nWarning %d:\nThe number of growth steps (-nxy + -nz) is larger than the "
1155 "number of steps in the tpr.\n"
1156 "(increase maxwarn in the membed input file to override)\n\n", warn
);
1160 if (inputrec
->opts
.ngfrz
== 1)
1162 gmx_fatal(FARGS
, "You did not specify \"%s\" as a freezegroup.", ins
);
1165 for (i
= 0; i
< inputrec
->opts
.ngfrz
; i
++)
1167 tmp_id
= mtop
->groups
.grps
[egcFREEZE
].nm_ind
[i
];
1168 if (ins_grp_id
== tmp_id
)
1177 gmx_fatal(FARGS
, "\"%s\" not as freezegroup defined in the mdp-file.", ins
);
1180 for (i
= 0; i
< DIM
; i
++)
1182 if (inputrec
->opts
.nFreeze
[fr_i
][i
] != 1)
1184 gmx_fatal(FARGS
, "freeze dimensions for %s are not Y Y Y\n", ins
);
1188 ng
= groups
->grps
[egcENER
].nr
;
1191 gmx_input("No energy groups defined. This is necessary for energy exclusion in the freeze group");
1194 for (i
= 0; i
< ng
; i
++)
1196 for (j
= 0; j
< ng
; j
++)
1198 if (inputrec
->opts
.egp_flags
[ng
*i
+j
] == EGP_EXCL
)
1201 if ( (groups
->grps
[egcENER
].nm_ind
[i
] != ins_grp_id
) ||
1202 (groups
->grps
[egcENER
].nm_ind
[j
] != ins_grp_id
) )
1204 gmx_fatal(FARGS
, "Energy exclusions \"%s\" and \"%s\" do not match the group "
1205 "to embed \"%s\"", *groups
->grpname
[groups
->grps
[egcENER
].nm_ind
[i
]],
1206 *groups
->grpname
[groups
->grps
[egcENER
].nm_ind
[j
]], ins
);
1214 gmx_input("No energy exclusion groups defined. This is necessary for energy exclusion in "
1215 "the freeze group");
1218 /* Obtain the maximum and minimum coordinates of the group to be embedded */
1220 init_ins_at(ins_at
, rest_at
, state
, pos_ins
, groups
, ins_grp_id
, xy_max
);
1221 /* Check that moleculetypes in insertion group are not part of the rest of the system */
1222 check_types(ins_at
, rest_at
, mtop
);
1224 init_mem_at(mem_p
, mtop
, as_rvec_array(state
->x
.data()), state
->box
, pos_ins
);
1226 prot_area
= est_prot_area(pos_ins
, as_rvec_array(state
->x
.data()), ins_at
, mem_p
);
1227 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
) )
1230 fprintf(stderr
, "\nWarning %d:\nThe xy-area is very small compared to the area of the protein.\n"
1231 "This might cause pressure problems during the growth phase. Just try with\n"
1232 "current setup and increase 'maxwarn' in your membed settings file, but lower the\n"
1233 "compressibility in the mdp-file or disable pressure coupling if problems occur.\n\n", warn
);
1238 gmx_fatal(FARGS
, "Too many warnings (override by setting maxwarn in the membed input file)\n");
1241 printf("The estimated area of the protein in the membrane is %.3f nm^2\n", prot_area
);
1242 printf("\nThere are %d lipids in the membrane part that overlaps the protein.\n"
1243 "The area per lipid is %.4f nm^2.\n", mem_p
->nmol
, mem_p
->lip_area
);
1245 /* Maximum number of lipids to be removed*/
1246 max_lip_rm
= (int)(2*prot_area
/mem_p
->lip_area
);
1247 printf("Maximum number of lipids that will be removed is %d.\n", max_lip_rm
);
1249 printf("\nWill resize the protein by a factor of %.3f in the xy plane and %.3f in the z direction.\n"
1250 "This resizing will be done with respect to the geometrical center of all protein atoms\n"
1251 "that span the membrane region, i.e. z between %.3f and %.3f\n\n",
1252 xy_fac
, z_fac
, mem_p
->zmin
, mem_p
->zmax
);
1254 /* resize the protein by xy and by z if necessary*/
1255 snew(r_ins
, ins_at
->nr
);
1256 init_resize(ins_at
, r_ins
, pos_ins
, mem_p
, as_rvec_array(state
->x
.data()), bALLOW_ASYMMETRY
);
1257 membed
->fac
[0] = membed
->fac
[1] = xy_fac
;
1258 membed
->fac
[2] = z_fac
;
1260 membed
->xy_step
= (xy_max
-xy_fac
)/(double)(it_xy
);
1261 membed
->z_step
= (z_max
-z_fac
)/(double)(it_z
-1);
1263 resize(r_ins
, as_rvec_array(state
->x
.data()), pos_ins
, membed
->fac
);
1265 /* remove overlapping lipids and water from the membrane box*/
1266 /*mark molecules to be removed*/
1268 set_pbc(pbc
, inputrec
->ePBC
, state
->box
);
1271 lip_rm
= gen_rm_list(rm_p
, ins_at
, rest_at
, pbc
, mtop
, as_rvec_array(state
->x
.data()), mem_p
, pos_ins
,
1272 probe_rad
, low_up_rm
, bALLOW_ASYMMETRY
);
1273 lip_rm
-= low_up_rm
;
1277 for (i
= 0; i
< rm_p
->nr
; i
++)
1279 fprintf(fplog
, "rm mol %d\n", rm_p
->mol
[i
]);
1283 for (i
= 0; i
< mtop
->nmolblock
; i
++)
1286 for (j
= 0; j
< rm_p
->nr
; j
++)
1288 if (rm_p
->block
[j
] == i
)
1293 printf("Will remove %d %s molecules\n", ntype
, *(mtop
->moltype
[mtop
->molblock
[i
].type
].name
));
1296 if (lip_rm
> max_lip_rm
)
1299 fprintf(stderr
, "\nWarning %d:\nTrying to remove a larger lipid area than the estimated "
1300 "protein area\nTry making the -xyinit resize factor smaller or increase "
1301 "maxwarn in the membed input file.\n\n", warn
);
1304 /*remove all lipids and waters overlapping and update all important structures*/
1305 rm_group(groups
, mtop
, rm_p
, state
, ins_at
, pos_ins
);
1307 rm_bonded_at
= rm_bonded(ins_at
, mtop
);
1308 if (rm_bonded_at
!= ins_at
->nr
)
1310 fprintf(stderr
, "Warning: The number of atoms for which the bonded interactions are removed is %d, "
1311 "while %d atoms are embedded. Make sure that the atoms to be embedded are not in the same"
1312 "molecule type as atoms that are not to be embedded.\n", rm_bonded_at
, ins_at
->nr
);
1317 gmx_fatal(FARGS
, "Too many warnings.\nIf you are sure these warnings are harmless,\n"
1318 "you can increase the maxwarn setting in the membed input file.");
1321 if (ftp2bSet(efTOP
, nfile
, fnm
))
1323 top_update(opt2fn("-mp", nfile
, fnm
), rm_p
, mtop
);
1333 membed
->it_xy
= it_xy
;
1334 membed
->it_z
= it_z
;
1335 membed
->pos_ins
= pos_ins
;
1336 membed
->r_ins
= r_ins
;
1342 void free_membed(gmx_membed_t
*membed
)