2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2010,2011,2012,2013,2014,2015, 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/essentialdynamics/edsam.h"
43 #include "gromacs/fileio/filenm.h"
44 #include "gromacs/fileio/tpxio.h"
45 #include "gromacs/gmxlib/network.h"
46 #include "gromacs/gmxlib/readinp.h"
47 #include "gromacs/math/vec.h"
48 #include "gromacs/mdtypes/commrec.h"
49 #include "gromacs/mdtypes/md_enums.h"
50 #include "gromacs/pbcutil/pbc.h"
51 #include "gromacs/topology/index.h"
52 #include "gromacs/topology/mtop_util.h"
53 #include "gromacs/topology/topology.h"
54 #include "gromacs/utility/cstringutil.h"
55 #include "gromacs/utility/fatalerror.h"
56 #include "gromacs/utility/futil.h"
57 #include "gromacs/utility/smalloc.h"
59 /* information about scaling center */
61 rvec xmin
; /* smallest coordinates of all embedded molecules */
62 rvec xmax
; /* largest coordinates of all embedded molecules */
63 rvec
*geom_cent
; /* scaling center of each independent molecule to embed */
64 int pieces
; /* number of molecules to embed independently */
65 int *nidx
; /* n atoms for every independent embedded molecule (index in subindex) */
66 int **subindex
; /* atomids for independent molecule *
67 * atoms of piece i run from subindex[i][0] to subindex[i][nidx[i]] */
70 /* variables needed in do_md */
72 int it_xy
; /* number of iterations (steps) used to grow something in the xy-plane */
73 int it_z
; /* same, but for z */
74 real xy_step
; /* stepsize used during resize in xy-plane */
75 real z_step
; /* same, but in z */
76 rvec fac
; /* initial scaling of the molecule to grow into the membrane */
77 rvec
*r_ins
; /* final coordinates of the molecule to grow */
78 pos_ins_t
*pos_ins
; /* scaling center for each piece to embed */
81 /* membrane related variables */
83 char *name
; /* name of index group to embed molecule into (usually membrane) */
84 t_block mem_at
; /* list all atoms in membrane */
85 int nmol
; /* number of membrane molecules overlapping with the molecule to embed */
86 int *mol_id
; /* list of molecules in membrane that overlap with the molecule to embed */
87 real lip_area
; /* average area per lipid in membrane (only correct for homogeneous bilayers)*/
88 real zmin
; /* minimum z coordinate of membrane */
89 real zmax
; /* maximum z coordinate of membrane */
90 real zmed
; /* median z coordinate of membrane */
93 /* Lists all molecules in the membrane that overlap with the molecule to be embedded. *
94 * These will then be removed from the system */
96 int nr
; /* number of molecules to remove */
97 int *mol
; /* list of molecule ids to remove */
98 int *block
; /* id of the molblock that the molecule to remove is part of */
101 /* Get the global molecule id, and the corresponding molecule type and id of the *
102 * molblock from the global atom nr. */
103 static int get_mol_id(int at
, gmx_mtop_t
*mtop
, int *type
, int *block
)
108 gmx_mtop_atomlookup_t alook
;
110 alook
= gmx_mtop_atomlookup_settle_init(mtop
);
111 gmx_mtop_atomnr_to_molblock_ind(alook
, at
, block
, &mol_id
, &atnr_mol
);
112 for (i
= 0; i
< *block
; i
++)
114 mol_id
+= mtop
->molblock
[i
].nmol
;
116 *type
= mtop
->molblock
[*block
].type
;
118 gmx_mtop_atomlookup_destroy(alook
);
123 /* Get the id of the molblock from a global molecule id */
124 static int get_molblock(int mol_id
, int nmblock
, gmx_molblock_t
*mblock
)
129 for (i
= 0; i
< nmblock
; i
++)
131 nmol
+= mblock
[i
].nmol
;
138 gmx_fatal(FARGS
, "mol_id %d larger than total number of molecules %d.\n", mol_id
, nmol
);
143 static int get_tpr_version(const char *infile
)
146 int version
, generation
;
148 read_tpxheader(infile
, &header
, TRUE
, &version
, &generation
);
153 /* Get a list of all the molecule types that are present in a group of atoms. *
154 * Because all interaction within the group to embed are removed on the topology *
155 * level, if the same molecule type is found in another part of the system, these *
156 * would also be affected. Therefore we have to check if the embedded and rest group *
157 * share common molecule types. If so, membed will stop with an error. */
158 static int get_mtype_list(t_block
*at
, gmx_mtop_t
*mtop
, t_block
*tlist
)
161 int type
= 0, block
= 0;
165 snew(tlist
->index
, at
->nr
);
166 for (i
= 0; i
< at
->nr
; i
++)
169 get_mol_id(at
->index
[i
], mtop
, &type
, &block
);
170 for (j
= 0; j
< nr
; j
++)
172 if (tlist
->index
[j
] == type
)
180 tlist
->index
[nr
] = type
;
184 srenew(tlist
->index
, nr
);
188 /* Do the actual check of the molecule types between embedded and rest group */
189 static void check_types(t_block
*ins_at
, t_block
*rest_at
, gmx_mtop_t
*mtop
)
191 t_block
*ins_mtype
, *rest_mtype
;
196 ins_mtype
->nr
= get_mtype_list(ins_at
, mtop
, ins_mtype
);
197 rest_mtype
->nr
= get_mtype_list(rest_at
, mtop
, rest_mtype
);
199 for (i
= 0; i
< ins_mtype
->nr
; i
++)
201 for (j
= 0; j
< rest_mtype
->nr
; j
++)
203 if (ins_mtype
->index
[i
] == rest_mtype
->index
[j
])
205 gmx_fatal(FARGS
, "Moleculetype %s is found both in the group to insert and the rest of the system.\n"
206 "1. Your *.ndx and *.top do not match\n"
207 "2. You are inserting some molecules of type %s (for example xray-solvent), while\n"
208 "the same moleculetype is also used in the rest of the system (solvent box). Because\n"
209 "we need to exclude all interactions between the atoms in the group to\n"
210 "insert, the same moleculetype can not be used in both groups. Change the\n"
211 "moleculetype of the molecules %s in the inserted group. Do not forget to provide\n"
212 "an appropriate *.itp file", *(mtop
->moltype
[rest_mtype
->index
[j
]].name
),
213 *(mtop
->moltype
[rest_mtype
->index
[j
]].name
), *(mtop
->moltype
[rest_mtype
->index
[j
]].name
));
218 sfree(ins_mtype
->index
);
219 sfree(rest_mtype
->index
);
224 static void get_input(const char *membed_input
, real
*xy_fac
, real
*xy_max
, real
*z_fac
, real
*z_max
,
225 int *it_xy
, int *it_z
, real
*probe_rad
, int *low_up_rm
, int *maxwarn
,
226 int *pieces
, gmx_bool
*bALLOW_ASYMMETRY
)
232 wi
= init_warning(TRUE
, 0);
234 inp
= read_inpfile(membed_input
, &ninp
, wi
);
235 ITYPE ("nxy", *it_xy
, 1000);
236 ITYPE ("nz", *it_z
, 0);
237 RTYPE ("xyinit", *xy_fac
, 0.5);
238 RTYPE ("xyend", *xy_max
, 1.0);
239 RTYPE ("zinit", *z_fac
, 1.0);
240 RTYPE ("zend", *z_max
, 1.0);
241 RTYPE ("rad", *probe_rad
, 0.22);
242 ITYPE ("ndiff", *low_up_rm
, 0);
243 ITYPE ("maxwarn", *maxwarn
, 0);
244 ITYPE ("pieces", *pieces
, 1);
245 EETYPE("asymmetry", *bALLOW_ASYMMETRY
, yesno_names
);
247 write_inpfile(membed_input
, ninp
, inp
, FALSE
, wi
);
250 /* Obtain the maximum and minimum coordinates of the group to be embedded */
251 static int init_ins_at(t_block
*ins_at
, t_block
*rest_at
, t_state
*state
, pos_ins_t
*pos_ins
,
252 gmx_groups_t
*groups
, int ins_grp_id
, real xy_max
)
255 real x
, xmin
, xmax
, y
, ymin
, ymax
, z
, zmin
, zmax
;
256 const real min_memthick
= 6.0; /* minimum thickness of the bilayer that will be used to *
257 * determine the overlap between molecule to embed and membrane */
258 const real fac_inp_size
= 1.000001; /* scaling factor to obtain input_size + 0.000001 (comparing reals) */
259 snew(rest_at
->index
, state
->natoms
);
261 xmin
= xmax
= state
->x
[ins_at
->index
[0]][XX
];
262 ymin
= ymax
= state
->x
[ins_at
->index
[0]][YY
];
263 zmin
= zmax
= state
->x
[ins_at
->index
[0]][ZZ
];
265 for (i
= 0; i
< state
->natoms
; i
++)
267 gid
= groups
->grpnr
[egcFREEZE
][i
];
268 if (groups
->grps
[egcFREEZE
].nm_ind
[gid
] == ins_grp_id
)
300 rest_at
->index
[c
] = i
;
306 srenew(rest_at
->index
, c
);
308 if (xy_max
> fac_inp_size
)
310 pos_ins
->xmin
[XX
] = xmin
-((xmax
-xmin
)*xy_max
-(xmax
-xmin
))/2;
311 pos_ins
->xmin
[YY
] = ymin
-((ymax
-ymin
)*xy_max
-(ymax
-ymin
))/2;
313 pos_ins
->xmax
[XX
] = xmax
+((xmax
-xmin
)*xy_max
-(xmax
-xmin
))/2;
314 pos_ins
->xmax
[YY
] = ymax
+((ymax
-ymin
)*xy_max
-(ymax
-ymin
))/2;
318 pos_ins
->xmin
[XX
] = xmin
;
319 pos_ins
->xmin
[YY
] = ymin
;
321 pos_ins
->xmax
[XX
] = xmax
;
322 pos_ins
->xmax
[YY
] = ymax
;
325 if ( (zmax
-zmin
) < min_memthick
)
327 pos_ins
->xmin
[ZZ
] = zmin
+(zmax
-zmin
)/2.0-0.5*min_memthick
;
328 pos_ins
->xmax
[ZZ
] = zmin
+(zmax
-zmin
)/2.0+0.5*min_memthick
;
332 pos_ins
->xmin
[ZZ
] = zmin
;
333 pos_ins
->xmax
[ZZ
] = zmax
;
339 /* Estimate the area of the embedded molecule by projecting all coordinates on a grid in the *
340 * xy-plane and counting the number of occupied grid points */
341 static real
est_prot_area(pos_ins_t
*pos_ins
, rvec
*r
, t_block
*ins_at
, mem_t
*mem_p
)
343 real x
, y
, dx
= 0.15, dy
= 0.15, area
= 0.0;
344 real add
, memmin
, memmax
;
347 /* min and max membrane coordinate are altered to reduce the influence of the *
349 memmin
= mem_p
->zmin
+0.1*(mem_p
->zmax
-mem_p
->zmin
);
350 memmax
= mem_p
->zmax
-0.1*(mem_p
->zmax
-mem_p
->zmin
);
352 for (x
= pos_ins
->xmin
[XX
]; x
< pos_ins
->xmax
[XX
]; x
+= dx
)
354 for (y
= pos_ins
->xmin
[YY
]; y
< pos_ins
->xmax
[YY
]; y
+= dy
)
360 at
= ins_at
->index
[c
];
361 if ( (r
[at
][XX
] >= x
) && (r
[at
][XX
] < x
+dx
) &&
362 (r
[at
][YY
] >= y
) && (r
[at
][YY
] < y
+dy
) &&
363 (r
[at
][ZZ
] > memmin
) && (r
[at
][ZZ
] < memmax
) )
369 while ( (c
< ins_at
->nr
) && (add
< 0.5) );
378 static int init_mem_at(mem_t
*mem_p
, gmx_mtop_t
*mtop
, rvec
*r
, matrix box
, pos_ins_t
*pos_ins
)
380 int i
, j
, at
, mol
, nmol
, nmolbox
, count
;
382 real z
, zmin
, zmax
, mem_area
;
385 int type
= 0, block
= 0;
388 mem_a
= &(mem_p
->mem_at
);
389 snew(mol_id
, mem_a
->nr
);
390 zmin
= pos_ins
->xmax
[ZZ
];
391 zmax
= pos_ins
->xmin
[ZZ
];
392 for (i
= 0; i
< mem_a
->nr
; i
++)
394 at
= mem_a
->index
[i
];
395 if ( (r
[at
][XX
] > pos_ins
->xmin
[XX
]) && (r
[at
][XX
] < pos_ins
->xmax
[XX
]) &&
396 (r
[at
][YY
] > pos_ins
->xmin
[YY
]) && (r
[at
][YY
] < pos_ins
->xmax
[YY
]) &&
397 (r
[at
][ZZ
] > pos_ins
->xmin
[ZZ
]) && (r
[at
][ZZ
] < pos_ins
->xmax
[ZZ
]) )
399 mol
= get_mol_id(at
, mtop
, &type
, &block
);
401 for (j
= 0; j
< nmol
; j
++)
403 if (mol
== mol_id
[j
])
431 srenew(mol_id
, nmol
);
432 mem_p
->mol_id
= mol_id
;
434 if ((zmax
-zmin
) > (box
[ZZ
][ZZ
]-0.5))
436 gmx_fatal(FARGS
, "Something is wrong with your membrane. Max and min z values are %f and %f.\n"
437 "Maybe your membrane is not centered in the box, but located at the box edge in the z-direction,\n"
438 "so that one membrane is distributed over two periodic box images. Another possibility is that\n"
439 "your water layer is not thick enough.\n", zmax
, zmin
);
443 mem_p
->zmed
= (zmax
-zmin
)/2+zmin
;
445 /*number of membrane molecules in protein box*/
446 nmolbox
= count
/mtop
->molblock
[block
].natoms_mol
;
447 /*membrane area within the box defined by the min and max coordinates of the embedded molecule*/
448 mem_area
= (pos_ins
->xmax
[XX
]-pos_ins
->xmin
[XX
])*(pos_ins
->xmax
[YY
]-pos_ins
->xmin
[YY
]);
449 /*rough estimate of area per lipid, assuming there is only one type of lipid in the membrane*/
450 mem_p
->lip_area
= 2.0*mem_area
/(double)nmolbox
;
452 return mem_p
->mem_at
.nr
;
455 static void init_resize(t_block
*ins_at
, rvec
*r_ins
, pos_ins_t
*pos_ins
, mem_t
*mem_p
, rvec
*r
,
456 gmx_bool bALLOW_ASYMMETRY
)
458 int i
, j
, at
, c
, outsidesum
, gctr
= 0;
462 for (i
= 0; i
< pos_ins
->pieces
; i
++)
464 idxsum
+= pos_ins
->nidx
[i
];
467 if (idxsum
!= ins_at
->nr
)
469 gmx_fatal(FARGS
, "Piecewise sum of inserted atoms not same as size of group selected to insert.");
472 snew(pos_ins
->geom_cent
, pos_ins
->pieces
);
473 for (i
= 0; i
< pos_ins
->pieces
; i
++)
477 for (j
= 0; j
< DIM
; j
++)
479 pos_ins
->geom_cent
[i
][j
] = 0;
482 for (j
= 0; j
< pos_ins
->nidx
[i
]; j
++)
484 at
= pos_ins
->subindex
[i
][j
];
485 copy_rvec(r
[at
], r_ins
[gctr
]);
486 if ( (r_ins
[gctr
][ZZ
] < mem_p
->zmax
) && (r_ins
[gctr
][ZZ
] > mem_p
->zmin
) )
488 rvec_inc(pos_ins
->geom_cent
[i
], r_ins
[gctr
]);
500 svmul(1/(double)c
, pos_ins
->geom_cent
[i
], pos_ins
->geom_cent
[i
]);
503 if (!bALLOW_ASYMMETRY
)
505 pos_ins
->geom_cent
[i
][ZZ
] = mem_p
->zmed
;
508 fprintf(stderr
, "Embedding piece %d with center of geometry: %f %f %f\n",
509 i
, pos_ins
->geom_cent
[i
][XX
], pos_ins
->geom_cent
[i
][YY
], pos_ins
->geom_cent
[i
][ZZ
]);
511 fprintf(stderr
, "\n");
514 /* resize performed in the md loop */
515 static void resize(rvec
*r_ins
, rvec
*r
, pos_ins_t
*pos_ins
, rvec fac
)
517 int i
, j
, k
, at
, c
= 0;
518 for (k
= 0; k
< pos_ins
->pieces
; k
++)
520 for (i
= 0; i
< pos_ins
->nidx
[k
]; i
++)
522 at
= pos_ins
->subindex
[k
][i
];
523 for (j
= 0; j
< DIM
; j
++)
525 r
[at
][j
] = pos_ins
->geom_cent
[k
][j
]+fac
[j
]*(r_ins
[c
][j
]-pos_ins
->geom_cent
[k
][j
]);
532 /* generate the list of membrane molecules that overlap with the molecule to be embedded. *
533 * The molecule to be embedded is already reduced in size. */
534 static int gen_rm_list(rm_t
*rm_p
, t_block
*ins_at
, t_block
*rest_at
, t_pbc
*pbc
, gmx_mtop_t
*mtop
,
535 rvec
*r
, mem_t
*mem_p
, pos_ins_t
*pos_ins
, real probe_rad
,
536 int low_up_rm
, gmx_bool bALLOW_ASYMMETRY
)
538 int i
, j
, k
, l
, at
, at2
, mol_id
;
539 int type
= 0, block
= 0;
540 int nrm
, nupper
, nlower
;
541 real r_min_rad
, z_lip
, min_norm
;
547 r_min_rad
= probe_rad
*probe_rad
;
548 snew(rm_p
->mol
, mtop
->mols
.nr
);
549 snew(rm_p
->block
, mtop
->mols
.nr
);
552 for (i
= 0; i
< ins_at
->nr
; i
++)
554 at
= ins_at
->index
[i
];
555 for (j
= 0; j
< rest_at
->nr
; j
++)
557 at2
= rest_at
->index
[j
];
558 pbc_dx(pbc
, r
[at
], r
[at2
], dr
);
560 if (norm2(dr
) < r_min_rad
)
562 mol_id
= get_mol_id(at2
, mtop
, &type
, &block
);
564 for (l
= 0; l
< nrm
; l
++)
566 if (rm_p
->mol
[l
] == mol_id
)
574 rm_p
->mol
[nrm
] = mol_id
;
575 rm_p
->block
[nrm
] = block
;
578 for (l
= 0; l
< mem_p
->nmol
; l
++)
580 if (mol_id
== mem_p
->mol_id
[l
])
582 for (k
= mtop
->mols
.index
[mol_id
]; k
< mtop
->mols
.index
[mol_id
+1]; k
++)
586 z_lip
/= mtop
->molblock
[block
].natoms_mol
;
587 if (z_lip
< mem_p
->zmed
)
602 /*make sure equal number of lipids from upper and lower layer are removed */
603 if ( (nupper
!= nlower
) && (!bALLOW_ASYMMETRY
) )
605 snew(dist
, mem_p
->nmol
);
606 snew(order
, mem_p
->nmol
);
607 for (i
= 0; i
< mem_p
->nmol
; i
++)
609 at
= mtop
->mols
.index
[mem_p
->mol_id
[i
]];
610 pbc_dx(pbc
, r
[at
], pos_ins
->geom_cent
[0], dr
);
611 if (pos_ins
->pieces
> 1)
614 min_norm
= norm2(dr
);
615 for (k
= 1; k
< pos_ins
->pieces
; k
++)
617 pbc_dx(pbc
, r
[at
], pos_ins
->geom_cent
[k
], dr_tmp
);
618 if (norm2(dr_tmp
) < min_norm
)
620 min_norm
= norm2(dr_tmp
);
621 copy_rvec(dr_tmp
, dr
);
625 dist
[i
] = dr
[XX
]*dr
[XX
]+dr
[YY
]*dr
[YY
];
627 while (j
>= 0 && dist
[i
] < dist
[order
[j
]])
629 order
[j
+1] = order
[j
];
636 while (nupper
!= nlower
)
638 mol_id
= mem_p
->mol_id
[order
[i
]];
639 block
= get_molblock(mol_id
, mtop
->nmolblock
, mtop
->molblock
);
641 for (l
= 0; l
< nrm
; l
++)
643 if (rm_p
->mol
[l
] == mol_id
)
652 for (k
= mtop
->mols
.index
[mol_id
]; k
< mtop
->mols
.index
[mol_id
+1]; k
++)
656 z_lip
/= mtop
->molblock
[block
].natoms_mol
;
657 if (nupper
> nlower
&& z_lip
< mem_p
->zmed
)
659 rm_p
->mol
[nrm
] = mol_id
;
660 rm_p
->block
[nrm
] = block
;
664 else if (nupper
< nlower
&& z_lip
> mem_p
->zmed
)
666 rm_p
->mol
[nrm
] = mol_id
;
667 rm_p
->block
[nrm
] = block
;
676 gmx_fatal(FARGS
, "Trying to remove more lipid molecules than there are in the membrane");
684 srenew(rm_p
->mol
, nrm
);
685 srenew(rm_p
->block
, nrm
);
687 return nupper
+nlower
;
690 /*remove all lipids and waters overlapping and update all important structures (e.g. state and mtop)*/
691 static void rm_group(gmx_groups_t
*groups
, gmx_mtop_t
*mtop
, rm_t
*rm_p
, t_state
*state
,
692 t_block
*ins_at
, pos_ins_t
*pos_ins
)
694 int i
, j
, k
, n
, rm
, mol_id
, at
, block
;
696 int *list
, *new_mols
;
697 unsigned char *new_egrp
[egcNR
];
701 snew(list
, state
->natoms
);
703 for (i
= 0; i
< rm_p
->nr
; i
++)
705 mol_id
= rm_p
->mol
[i
];
706 at
= mtop
->mols
.index
[mol_id
];
707 block
= rm_p
->block
[i
];
708 mtop
->molblock
[block
].nmol
--;
709 for (j
= 0; j
< mtop
->molblock
[block
].natoms_mol
; j
++)
714 mtop
->mols
.index
[mol_id
] = -1;
717 mtop
->mols
.nr
-= rm_p
->nr
;
718 mtop
->mols
.nalloc_index
-= rm_p
->nr
;
719 snew(new_mols
, mtop
->mols
.nr
);
720 for (i
= 0; i
< mtop
->mols
.nr
+rm_p
->nr
; i
++)
723 if (mtop
->mols
.index
[i
] != -1)
725 new_mols
[j
] = mtop
->mols
.index
[i
];
729 sfree(mtop
->mols
.index
);
730 mtop
->mols
.index
= new_mols
;
733 state
->nalloc
= state
->natoms
;
734 snew(x_tmp
, state
->nalloc
);
735 snew(v_tmp
, state
->nalloc
);
737 for (i
= 0; i
< egcNR
; i
++)
739 if (groups
->grpnr
[i
] != NULL
)
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
] != NULL
)
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
;
795 for (i
= 0; i
< egcNR
; i
++)
797 if (groups
->grpnr
[i
] != NULL
)
799 sfree(groups
->grpnr
[i
]);
800 groups
->grpnr
[i
] = new_egrp
[i
];
804 /* remove empty molblocks */
806 for (i
= 0; i
< mtop
->nmolblock
; i
++)
808 if (mtop
->molblock
[i
].nmol
== 0)
814 mtop
->molblock
[i
-RMmolblock
] = mtop
->molblock
[i
];
817 mtop
->nmolblock
-= RMmolblock
;
820 /* remove al bonded interactions from mtop for the molecule to be embedded */
821 int rm_bonded(t_block
*ins_at
, gmx_mtop_t
*mtop
)
824 int type
, natom
, nmol
, at
, atom1
= 0, rm_at
= 0;
826 /*this routine lives dangerously by assuming that all molecules of a given type are in order in the structure*/
827 /*this routine does not live as dangerously as it seems. There is namely a check in init_membed to make *
828 * sure that g_membed exits with a warning when there are molecules of the same type not in the *
829 * ins_at index group. MGWolf 050710 */
832 snew(bRM
, mtop
->nmoltype
);
833 for (i
= 0; i
< mtop
->nmoltype
; i
++)
838 for (i
= 0; i
< mtop
->nmolblock
; i
++)
840 /*loop over molecule blocks*/
841 type
= mtop
->molblock
[i
].type
;
842 natom
= mtop
->molblock
[i
].natoms_mol
;
843 nmol
= mtop
->molblock
[i
].nmol
;
845 for (j
= 0; j
< natom
*nmol
&& bRM
[type
] == TRUE
; j
++)
847 /*loop over atoms in the block*/
848 at
= j
+atom1
; /*atom index = block index + offset*/
851 for (m
= 0; (m
< ins_at
->nr
) && (bINS
== FALSE
); m
++)
853 /*loop over atoms in insertion index group to determine if we're inserting one*/
854 if (at
== ins_at
->index
[m
])
861 atom1
+= natom
*nmol
; /*update offset*/
864 rm_at
+= natom
*nmol
; /*increment bonded removal counter by # atoms in block*/
868 for (i
= 0; i
< mtop
->nmoltype
; i
++)
872 for (j
= 0; j
< F_LJ
; j
++)
874 mtop
->moltype
[i
].ilist
[j
].nr
= 0;
877 for (j
= F_POSRES
; j
<= F_VSITEN
; j
++)
879 mtop
->moltype
[i
].ilist
[j
].nr
= 0;
888 /* Write a topology where the number of molecules is correct for the system after embedding */
889 static void top_update(const char *topfile
, rm_t
*rm_p
, gmx_mtop_t
*mtop
)
893 char buf
[STRLEN
], buf2
[STRLEN
], *temp
;
894 int i
, *nmol_rm
, nmol
, line
;
895 char temporary_filename
[STRLEN
];
897 fpin
= gmx_ffopen(topfile
, "r");
898 strncpy(temporary_filename
, "temp.topXXXXXX", STRLEN
);
899 gmx_tmpnam(temporary_filename
);
900 fpout
= gmx_ffopen(temporary_filename
, "w");
902 snew(nmol_rm
, mtop
->nmoltype
);
903 for (i
= 0; i
< rm_p
->nr
; i
++)
905 nmol_rm
[rm_p
->block
[i
]]++;
909 while (fgets(buf
, STRLEN
, fpin
))
915 if ((temp
= strchr(buf2
, '\n')) != NULL
)
923 if ((temp
= strchr(buf2
, '\n')) != NULL
)
928 if (buf2
[strlen(buf2
)-1] == ']')
930 buf2
[strlen(buf2
)-1] = '\0';
933 if (gmx_strcasecmp(buf2
, "molecules") == 0)
938 fprintf(fpout
, "%s", buf
);
940 else if (bMolecules
== 1)
942 for (i
= 0; i
< mtop
->nmolblock
; i
++)
944 nmol
= mtop
->molblock
[i
].nmol
;
945 sprintf(buf
, "%-15s %5d\n", *(mtop
->moltype
[mtop
->molblock
[i
].type
].name
), nmol
);
946 fprintf(fpout
, "%s", buf
);
950 else if (bMolecules
== 2)
956 fprintf(fpout
, "%s", buf
);
961 fprintf(fpout
, "%s", buf
);
966 /* use gmx_ffopen to generate backup of topinout */
967 fpout
= gmx_ffopen(topfile
, "w");
969 rename(temporary_filename
, topfile
);
972 void rescale_membed(int step_rel
, gmx_membed_t
*membed
, rvec
*x
)
974 /* Set new positions for the group to embed */
975 if (step_rel
<= membed
->it_xy
)
977 membed
->fac
[0] += membed
->xy_step
;
978 membed
->fac
[1] += membed
->xy_step
;
980 else if (step_rel
<= (membed
->it_xy
+membed
->it_z
))
982 membed
->fac
[2] += membed
->z_step
;
984 resize(membed
->r_ins
, x
, membed
->pos_ins
, membed
->fac
);
987 /* We would like gn to be const as well, but C doesn't allow this */
988 /* TODO this is utility functionality (search for the index of a
989 string in a collection), so should be refactored and located more
990 centrally. Also, it nearly duplicates the same string in readir.c */
991 static int search_string(const char *s
, int ng
, char *gn
[])
995 for (i
= 0; (i
< ng
); i
++)
997 if (gmx_strcasecmp(s
, gn
[i
]) == 0)
1004 "Group %s selected for embedding was not found in the index file.\n"
1005 "Group names must match either [moleculetype] names or custom index group\n"
1006 "names, in which case you must supply an index file to the '-n' option\n"
1013 gmx_membed_t
*init_membed(FILE *fplog
, int nfile
, const t_filenm fnm
[], gmx_mtop_t
*mtop
,
1014 t_inputrec
*inputrec
, t_state
*state
, t_commrec
*cr
, real
*cpt
)
1016 char *ins
, **gnames
;
1017 int i
, rm_bonded_at
, fr_id
, fr_i
= 0, tmp_id
, warn
= 0;
1018 int ng
, j
, max_lip_rm
, ins_grp_id
, ntype
, lip_rm
, tpr_version
;
1021 t_block
*ins_at
, *rest_at
;
1025 gmx_groups_t
*groups
;
1026 gmx_bool bExcl
= FALSE
;
1029 char **piecename
= NULL
;
1030 gmx_membed_t
*membed
= NULL
;
1032 /* input variables */
1033 const char *membed_input
;
1040 real probe_rad
= 0.22;
1044 gmx_bool bALLOW_ASYMMETRY
= FALSE
;
1046 /* sanity check constants */ /* Issue a warning when: */
1047 const int membed_version
= 58; /* tpr version is smaller */
1048 const real min_probe_rad
= 0.2199999; /* A probe radius for overlap between embedded molecule *
1049 * and rest smaller than this value is probably too small */
1050 const real min_xy_init
= 0.0999999; /* the initial shrinking of the molecule to embed is smaller */
1051 const int min_it_xy
= 1000; /* the number of steps to embed in xy-plane is smaller */
1052 const int min_it_z
= 100; /* the number of steps to embed in z is smaller */
1053 const real prot_vs_box
= 7.5; /* molecule to embed is large (more then prot_vs_box) with respect */
1054 const real box_vs_prot
= 50; /* to the box size (less than box_vs_prot) */
1062 /* get input data out membed file */
1063 membed_input
= opt2fn("-membed", nfile
, fnm
);
1064 get_input(membed_input
, &xy_fac
, &xy_max
, &z_fac
, &z_max
, &it_xy
, &it_z
, &probe_rad
, &low_up_rm
,
1065 &maxwarn
, &pieces
, &bALLOW_ASYMMETRY
);
1067 tpr_version
= get_tpr_version(ftp2fn(efTPR
, nfile
, fnm
));
1068 if (tpr_version
< membed_version
)
1070 gmx_fatal(FARGS
, "Version of *.tpr file to old (%d). "
1071 "Rerun grompp with GROMACS version 4.0.3 or newer.\n", tpr_version
);
1074 if (!EI_DYNAMICS(inputrec
->eI
) )
1076 gmx_input("Change integrator to a dynamics integrator in mdp file (e.g. md or sd).");
1081 gmx_input("Sorry, parallel g_membed is not yet fully functional.");
1086 fprintf(stderr
, "\nSetting -cpt to -1, because embedding cannot be restarted from cpt-files.\n");
1089 groups
= &(mtop
->groups
);
1090 snew(gnames
, groups
->ngrpname
);
1091 for (i
= 0; i
< groups
->ngrpname
; i
++)
1093 gnames
[i
] = *(groups
->grpname
[i
]);
1096 atoms
= gmx_mtop_global_atoms(mtop
);
1098 fprintf(stderr
, "\nSelect a group to embed in the membrane:\n");
1099 get_index(&atoms
, opt2fn_null("-mn", nfile
, fnm
), 1, &(ins_at
->nr
), &(ins_at
->index
), &ins
);
1100 ins_grp_id
= search_string(ins
, groups
->ngrpname
, gnames
);
1101 fprintf(stderr
, "\nSelect a group to embed %s into (e.g. the membrane):\n", ins
);
1102 get_index(&atoms
, opt2fn_null("-mn", nfile
, fnm
), 1, &(mem_p
->mem_at
.nr
), &(mem_p
->mem_at
.index
), &(mem_p
->name
));
1104 pos_ins
->pieces
= pieces
;
1105 snew(pos_ins
->nidx
, pieces
);
1106 snew(pos_ins
->subindex
, pieces
);
1107 snew(piecename
, pieces
);
1110 fprintf(stderr
, "\nSelect pieces to embed:\n");
1111 get_index(&atoms
, opt2fn_null("-mn", nfile
, fnm
), pieces
, pos_ins
->nidx
, pos_ins
->subindex
, piecename
);
1115 /*use whole embedded group*/
1116 snew(pos_ins
->nidx
, 1);
1117 snew(pos_ins
->subindex
, 1);
1118 pos_ins
->nidx
[0] = ins_at
->nr
;
1119 pos_ins
->subindex
[0] = ins_at
->index
;
1122 if (probe_rad
< min_probe_rad
)
1125 fprintf(stderr
, "\nWarning %d:\nA probe radius (-rad) smaller than 0.2 nm can result "
1126 "in overlap between waters and the group to embed, which will result "
1127 "in Lincs errors etc.\n\n", warn
);
1130 if (xy_fac
< min_xy_init
)
1133 fprintf(stderr
, "\nWarning %d:\nThe initial size of %s is probably too smal.\n\n", warn
, ins
);
1136 if (it_xy
< min_it_xy
)
1139 fprintf(stderr
, "\nWarning %d;\nThe number of steps used to grow the xy-coordinates of %s (%d)"
1140 " is probably too small.\nIncrease -nxy or.\n\n", warn
, ins
, it_xy
);
1143 if ( (it_z
< min_it_z
) && ( z_fac
< 0.99999999 || z_fac
> 1.0000001) )
1146 fprintf(stderr
, "\nWarning %d;\nThe number of steps used to grow the z-coordinate of %s (%d)"
1147 " is probably too small.\nIncrease -nz or the maxwarn setting in the membed input file.\n\n", warn
, ins
, it_z
);
1150 if (it_xy
+it_z
> inputrec
->nsteps
)
1153 fprintf(stderr
, "\nWarning %d:\nThe number of growth steps (-nxy + -nz) is larger than the "
1154 "number of steps in the tpr.\n"
1155 "(increase maxwarn in the membed input file to override)\n\n", warn
);
1159 if (inputrec
->opts
.ngfrz
== 1)
1161 gmx_fatal(FARGS
, "You did not specify \"%s\" as a freezegroup.", ins
);
1164 for (i
= 0; i
< inputrec
->opts
.ngfrz
; i
++)
1166 tmp_id
= mtop
->groups
.grps
[egcFREEZE
].nm_ind
[i
];
1167 if (ins_grp_id
== tmp_id
)
1176 gmx_fatal(FARGS
, "\"%s\" not as freezegroup defined in the mdp-file.", ins
);
1179 for (i
= 0; i
< DIM
; i
++)
1181 if (inputrec
->opts
.nFreeze
[fr_i
][i
] != 1)
1183 gmx_fatal(FARGS
, "freeze dimensions for %s are not Y Y Y\n", ins
);
1187 ng
= groups
->grps
[egcENER
].nr
;
1190 gmx_input("No energy groups defined. This is necessary for energy exclusion in the freeze group");
1193 for (i
= 0; i
< ng
; i
++)
1195 for (j
= 0; j
< ng
; j
++)
1197 if (inputrec
->opts
.egp_flags
[ng
*i
+j
] == EGP_EXCL
)
1200 if ( (groups
->grps
[egcENER
].nm_ind
[i
] != ins_grp_id
) ||
1201 (groups
->grps
[egcENER
].nm_ind
[j
] != ins_grp_id
) )
1203 gmx_fatal(FARGS
, "Energy exclusions \"%s\" and \"%s\" do not match the group "
1204 "to embed \"%s\"", *groups
->grpname
[groups
->grps
[egcENER
].nm_ind
[i
]],
1205 *groups
->grpname
[groups
->grps
[egcENER
].nm_ind
[j
]], ins
);
1213 gmx_input("No energy exclusion groups defined. This is necessary for energy exclusion in "
1214 "the freeze group");
1217 /* Obtain the maximum and minimum coordinates of the group to be embedded */
1219 init_ins_at(ins_at
, rest_at
, state
, pos_ins
, groups
, ins_grp_id
, xy_max
);
1220 /* Check that moleculetypes in insertion group are not part of the rest of the system */
1221 check_types(ins_at
, rest_at
, mtop
);
1223 init_mem_at(mem_p
, mtop
, state
->x
, state
->box
, pos_ins
);
1225 prot_area
= est_prot_area(pos_ins
, state
->x
, ins_at
, mem_p
);
1226 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
) )
1229 fprintf(stderr
, "\nWarning %d:\nThe xy-area is very small compared to the area of the protein.\n"
1230 "This might cause pressure problems during the growth phase. Just try with\n"
1231 "current setup and increase 'maxwarn' in your membed settings file, but lower the\n"
1232 "compressibility in the mdp-file or disable pressure coupling if problems occur.\n\n", warn
);
1237 gmx_fatal(FARGS
, "Too many warnings (override by setting maxwarn in the membed input file)\n");
1240 printf("The estimated area of the protein in the membrane is %.3f nm^2\n", prot_area
);
1241 printf("\nThere are %d lipids in the membrane part that overlaps the protein.\n"
1242 "The area per lipid is %.4f nm^2.\n", mem_p
->nmol
, mem_p
->lip_area
);
1244 /* Maximum number of lipids to be removed*/
1245 max_lip_rm
= (int)(2*prot_area
/mem_p
->lip_area
);
1246 printf("Maximum number of lipids that will be removed is %d.\n", max_lip_rm
);
1248 printf("\nWill resize the protein by a factor of %.3f in the xy plane and %.3f in the z direction.\n"
1249 "This resizing will be done with respect to the geometrical center of all protein atoms\n"
1250 "that span the membrane region, i.e. z between %.3f and %.3f\n\n",
1251 xy_fac
, z_fac
, mem_p
->zmin
, mem_p
->zmax
);
1253 /* resize the protein by xy and by z if necessary*/
1254 snew(r_ins
, ins_at
->nr
);
1255 init_resize(ins_at
, r_ins
, pos_ins
, mem_p
, state
->x
, bALLOW_ASYMMETRY
);
1256 membed
->fac
[0] = membed
->fac
[1] = xy_fac
;
1257 membed
->fac
[2] = z_fac
;
1259 membed
->xy_step
= (xy_max
-xy_fac
)/(double)(it_xy
);
1260 membed
->z_step
= (z_max
-z_fac
)/(double)(it_z
-1);
1262 resize(r_ins
, state
->x
, pos_ins
, membed
->fac
);
1264 /* remove overlapping lipids and water from the membrane box*/
1265 /*mark molecules to be removed*/
1267 set_pbc(pbc
, inputrec
->ePBC
, state
->box
);
1270 lip_rm
= gen_rm_list(rm_p
, ins_at
, rest_at
, pbc
, mtop
, state
->x
, mem_p
, pos_ins
,
1271 probe_rad
, low_up_rm
, bALLOW_ASYMMETRY
);
1272 lip_rm
-= low_up_rm
;
1276 for (i
= 0; i
< rm_p
->nr
; i
++)
1278 fprintf(fplog
, "rm mol %d\n", rm_p
->mol
[i
]);
1282 for (i
= 0; i
< mtop
->nmolblock
; i
++)
1285 for (j
= 0; j
< rm_p
->nr
; j
++)
1287 if (rm_p
->block
[j
] == i
)
1292 printf("Will remove %d %s molecules\n", ntype
, *(mtop
->moltype
[mtop
->molblock
[i
].type
].name
));
1295 if (lip_rm
> max_lip_rm
)
1298 fprintf(stderr
, "\nWarning %d:\nTrying to remove a larger lipid area than the estimated "
1299 "protein area\nTry making the -xyinit resize factor smaller or increase "
1300 "maxwarn in the membed input file.\n\n", warn
);
1303 /*remove all lipids and waters overlapping and update all important structures*/
1304 rm_group(groups
, mtop
, rm_p
, state
, ins_at
, pos_ins
);
1306 rm_bonded_at
= rm_bonded(ins_at
, mtop
);
1307 if (rm_bonded_at
!= ins_at
->nr
)
1309 fprintf(stderr
, "Warning: The number of atoms for which the bonded interactions are removed is %d, "
1310 "while %d atoms are embedded. Make sure that the atoms to be embedded are not in the same"
1311 "molecule type as atoms that are not to be embedded.\n", rm_bonded_at
, ins_at
->nr
);
1316 gmx_fatal(FARGS
, "Too many warnings.\nIf you are sure these warnings are harmless,\n"
1317 "you can increase the maxwarn setting in the membed input file.");
1320 if (ftp2bSet(efTOP
, nfile
, fnm
))
1322 top_update(opt2fn("-mp", nfile
, fnm
), rm_p
, mtop
);
1332 membed
->it_xy
= it_xy
;
1333 membed
->it_z
= it_z
;
1334 membed
->pos_ins
= pos_ins
;
1335 membed
->r_ins
= r_ins
;
1341 void free_membed(gmx_membed_t
*membed
)