2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
47 #include "gromacs/utility/smalloc.h"
48 #include "types/commrec.h"
55 #include "mtop_util.h"
56 #include "chargegroup.h"
58 static real box_margin;
60 static real max_dist(rvec *x, real *r, int start, int end)
66 for (i = start; i < end; i++)
68 for (j = i+1; j < end; j++)
70 maxd = max(maxd, sqrt(distance2(x[i], x[j]))+0.5*(r[i]+r[j]));
77 static gmx_bool outside_box_minus_margin2(rvec x, matrix box)
79 return ( (x[XX] < 2*box_margin) || (x[XX] > box[XX][XX]-2*box_margin) ||
80 (x[YY] < 2*box_margin) || (x[YY] > box[YY][YY]-2*box_margin) ||
81 (x[ZZ] < 2*box_margin) || (x[ZZ] > box[ZZ][ZZ]-2*box_margin) );
84 static gmx_bool outside_box_plus_margin(rvec x, matrix box)
86 return ( (x[XX] < -box_margin) || (x[XX] > box[XX][XX]+box_margin) ||
87 (x[YY] < -box_margin) || (x[YY] > box[YY][YY]+box_margin) ||
88 (x[ZZ] < -box_margin) || (x[ZZ] > box[ZZ][ZZ]+box_margin) );
91 static int mark_res(int at, gmx_bool *mark, int natoms, t_atom *atom, int *nmark)
95 resind = atom[at].resind;
96 while ( (at > 0) && (resind == atom[at-1].resind) )
100 while ( (at < natoms) && (resind == atom[at].resind) )
113 static real find_max_real(int n, real radius[])
122 for (i = 1; (i < n); i++)
124 rmax = max(rmax, radius[i]);
130 static void combine_atoms(t_atoms *ap, t_atoms *as,
131 rvec xp[], rvec *vp, rvec xs[], rvec *vs,
132 t_atoms **a_comb, rvec **x_comb, rvec **v_comb)
135 rvec *xc, *vc = NULL;
136 int i, j, natot, res0;
138 /* Total number of atoms */
139 natot = ap->nr+as->nr;
142 init_t_atoms(ac, natot, FALSE);
150 /* Fill the new structures */
151 for (i = j = 0; (i < ap->nr); i++, j++)
153 copy_rvec(xp[i], xc[j]);
156 copy_rvec(vp[i], vc[j]);
158 memcpy(&(ac->atom[j]), &(ap->atom[i]), sizeof(ap->atom[i]));
159 ac->atom[j].type = 0;
162 for (i = 0; (i < as->nr); i++, j++)
164 copy_rvec(xs[i], xc[j]);
167 copy_rvec(vs[i], vc[j]);
169 memcpy(&(ac->atom[j]), &(as->atom[i]), sizeof(as->atom[i]));
170 ac->atom[j].type = 0;
171 ac->atom[j].resind += res0;
174 ac->nres = ac->atom[j-1].resind+1;
175 /* Fill all elements to prevent uninitialized memory */
176 for (i = 0; i < ac->nr; i++)
182 ac->atom[i].type = 0;
183 ac->atom[i].typeB = 0;
184 ac->atom[i].ptype = eptAtom;
193 static t_forcerec *fr = NULL;
195 static void do_nsgrid(FILE *fp, gmx_bool bVerbose,
196 matrix box, rvec x[], t_atoms *atoms, real rlong,
197 const output_env_t oenv)
216 /* Charge group index */
217 snew(cg_index, natoms);
218 for (i = 0; (i < natoms); i++)
223 /* Topology needs charge groups and exclusions */
226 mtop->natoms = natoms;
227 /* Make one moltype that contains the whol system */
229 snew(mtop->moltype, mtop->nmoltype);
230 molt = &mtop->moltype[0];
231 molt->name = mtop->name;
232 molt->atoms = *atoms;
233 stupid_fill_block(&molt->cgs, mtop->natoms, FALSE);
234 stupid_fill_blocka(&molt->excls, natoms);
235 /* Make one molblock for the whole system */
237 snew(mtop->molblock, mtop->nmolblock);
238 mtop->molblock[0].type = 0;
239 mtop->molblock[0].nmol = 1;
240 mtop->molblock[0].natoms_mol = natoms;
241 /* Initialize a single energy group */
242 mtop->groups.grps[egcENER].nr = 1;
243 mtop->groups.ngrpnr[egcENER] = 0;
244 mtop->groups.grpnr[egcENER] = NULL;
246 ffp = &mtop->ffparams;
251 snew(ffp->functype, 1);
252 snew(ffp->iparams, 1);
253 ffp->iparams[0].lj.c6 = 1;
254 ffp->iparams[0].lj.c12 = 1;
256 /* inputrec structure */
258 ir->cutoff_scheme = ecutsGROUP;
259 ir->coulomb_modifier = eintmodNONE;
260 ir->vdw_modifier = eintmodNONE;
261 ir->coulombtype = eelCUT;
262 ir->vdwtype = evdwCUT;
264 ir->ns_type = ensGRID;
265 snew(ir->opts.egp_flags, 1);
267 top = gmx_mtop_generate_local_top(mtop, ir);
269 /* Some nasty shortcuts */
272 /* mdatoms structure */
275 md = init_mdatoms(fp, mtop, FALSE);
276 atoms2md(mtop, ir, 0, NULL, mtop->natoms, md);
279 /* forcerec structure */
286 /* cr->nthreads = 1; */
288 /* ir->rlist = ir->rcoulomb = ir->rvdw = rlong;
289 printf("Neighborsearching with a cut-off of %g\n",rlong);
290 init_forcerec(stdout,fr,ir,top,cr,md,box,FALSE,NULL,NULL,NULL,TRUE);*/
292 fr->hcg = top->cgs.nr;
295 /* Prepare for neighboursearching */
298 /* Init things dependent on parameters */
299 ir->rlistlong = ir->rlist = ir->rcoulomb = ir->rvdw = rlong;
300 /* create free energy data to avoid NULLs */
301 snew(ir->fepvals, 1);
302 printf("Neighborsearching with a cut-off of %g\n", rlong);
303 init_forcerec(stdout, oenv, fr, NULL, ir, mtop, cr, box,
304 NULL, NULL, NULL, NULL, NULL, TRUE, -1);
307 pr_forcerec(debug, fr);
310 /* Calculate new stuff dependent on coords and box */
311 for (m = 0; (m < DIM); m++)
313 box_size[m] = box[m][m];
315 calc_shifts(box, fr->shift_vec);
316 put_charge_groups_in_box(fp, 0, cgs->nr, fr->ePBC, box, cgs, x, fr->cg_cm);
318 /* Do the actual neighboursearching */
319 snew(lambda, efptNR);
321 init_neighbor_list(fp, fr, md->homenr);
322 search_neighbours(fp, fr, box, top,
323 &mtop->groups, cr, &nrnb, md, TRUE, FALSE);
327 dump_nblist(debug, cr, fr, 0);
332 fprintf(stderr, "Successfully made neighbourlist\n");
336 static gmx_bool bXor(gmx_bool b1, gmx_bool b2)
338 return (b1 && !b2) || (b2 && !b1);
341 void add_conf(t_atoms *atoms, rvec **x, rvec **v, real **r, gmx_bool bSrenew,
342 int ePBC, matrix box, gmx_bool bInsert,
343 t_atoms *atoms_solvt, rvec *x_solvt, rvec *v_solvt, real *r_solvt,
344 gmx_bool bVerbose, real rshell, int max_sol, const output_env_t oenv)
348 real max_vdw, *r_prot, *r_all, n2, r2, ib1, ib2;
349 int natoms_prot, natoms_solvt;
350 int i, j, jj, m, j0, j1, jjj, jnres, jnr, inr, iprot, is1, is2;
351 int prev, resnr, nresadd, d, k, ncells, maxincell;
352 int dx0, dx1, dy0, dy1, dz0, dz1;
353 int ntest, nremove, nkeep;
354 rvec dx, xi, xj, xpp, *x_all, *v_all;
355 gmx_bool *remove, *keep;
358 natoms_prot = atoms->nr;
359 natoms_solvt = atoms_solvt->nr;
360 if (natoms_solvt <= 0)
362 fprintf(stderr, "WARNING: Nothing to add\n");
366 if (ePBC == epbcSCREW)
368 gmx_fatal(FARGS, "Sorry, %s pbc is not yet supported", epbc_names[ePBC]);
373 fprintf(stderr, "Calculating Overlap...\n");
376 /* Set margin around box edges to largest solvent dimension.
377 * The maximum distance between atoms in a solvent molecule should
378 * be calculated. At the moment a fudge factor of 3 is used.
381 box_margin = 3*find_max_real(natoms_solvt, r_solvt);
382 max_vdw = max(3*find_max_real(natoms_prot, r_prot), box_margin);
383 fprintf(stderr, "box_margin = %g\n", box_margin);
385 snew(remove, natoms_solvt);
390 for (i = 0; i < atoms_solvt->nr; i++)
392 if (outside_box_plus_margin(x_solvt[i], box) )
394 i = mark_res(i, remove, atoms_solvt->nr, atoms_solvt->atom, &nremove);
397 fprintf(stderr, "Removed %d atoms that were outside the box\n", nremove);
400 /* Define grid stuff */
401 /* Largest VDW radius */
402 snew(r_all, natoms_prot+natoms_solvt);
403 for (i = j = 0; i < natoms_prot; i++, j++)
405 r_all[j] = r_prot[i];
407 for (i = 0; i < natoms_solvt; i++, j++)
409 r_all[j] = r_solvt[i];
413 combine_atoms(atoms, atoms_solvt, *x, v ? *v : NULL, x_solvt, v_solvt,
414 &atoms_all, &x_all, &v_all);
416 /* Do neighboursearching step */
417 do_nsgrid(stdout, bVerbose, box, x_all, atoms_all, max_vdw, oenv);
419 /* check solvent with solute */
420 nlist = &(fr->nblists[0].nlist_sr[eNL_VDW]);
421 fprintf(stderr, "nri = %d, nrj = %d\n", nlist->nri, nlist->nrj);
422 for (bSolSol = 0; (bSolSol <= (bInsert ? 0 : 1)); bSolSol++)
425 fprintf(stderr, "Checking %s-Solvent overlap:",
426 bSolSol ? "Solvent" : "Protein");
427 for (i = 0; (i < nlist->nri && nremove < natoms_solvt); i++)
429 inr = nlist->iinr[i];
430 j0 = nlist->jindex[i];
431 j1 = nlist->jindex[i+1];
432 rvec_add(x_all[inr], fr->shift_vec[nlist->shift[i]], xi);
434 for (j = j0; (j < j1 && nremove < natoms_solvt); j++)
436 jnr = nlist->jjnr[j];
437 copy_rvec(x_all[jnr], xj);
439 /* Check solvent-protein and solvent-solvent */
440 is1 = inr-natoms_prot;
441 is2 = jnr-natoms_prot;
443 /* Check if at least one of the atoms is a solvent that is not yet
444 * listed for removal, and if both are solvent, that they are not in the
448 bXor((is1 >= 0), (is2 >= 0)) && /* One atom is protein */
449 ((is1 < 0) || ((is1 >= 0) && !remove[is1])) &&
450 ((is2 < 0) || ((is2 >= 0) && !remove[is2]))) ||
453 (is1 >= 0) && (!remove[is1]) && /* is1 is solvent */
454 (is2 >= 0) && (!remove[is2]) && /* is2 is solvent */
455 (bInsert || /* when inserting also check inside the box */
456 (outside_box_minus_margin2(x_solvt[is1], box) && /* is1 on edge */
457 outside_box_minus_margin2(x_solvt[is2], box)) /* is2 on edge */
459 (atoms_solvt->atom[is1].resind != /* Not the same residue */
460 atoms_solvt->atom[is2].resind)))
464 rvec_sub(xi, xj, dx);
466 r2 = sqr(r_all[inr]+r_all[jnr]);
471 nremove = natoms_solvt;
472 for (k = 0; k < nremove; k++)
477 /* Need only remove one of the solvents... */
480 (void) mark_res(is2, remove, natoms_solvt, atoms_solvt->atom,
485 (void) mark_res(is1, remove, natoms_solvt, atoms_solvt->atom,
490 fprintf(stderr, "Neither atom is solvent%d %d\n", is1, is2);
498 fprintf(stderr, " tested %d pairs, removed %d atoms.\n", ntest, nremove);
503 for (i = 0; i < natoms_solvt; i++)
505 fprintf(debug, "remove[%5d] = %s\n", i, bool_names[remove[i]]);
509 /* Search again, now with another cut-off */
512 do_nsgrid(stdout, bVerbose, box, x_all, atoms_all, rshell, oenv);
513 nlist = &(fr->nblists[0].nlist_sr[eNL_VDW]);
514 fprintf(stderr, "nri = %d, nrj = %d\n", nlist->nri, nlist->nrj);
516 snew(keep, natoms_solvt);
517 for (i = 0; i < nlist->nri; i++)
519 inr = nlist->iinr[i];
520 j0 = nlist->jindex[i];
521 j1 = nlist->jindex[i+1];
523 for (j = j0; j < j1; j++)
525 jnr = nlist->jjnr[j];
527 /* Check solvent-protein and solvent-solvent */
528 is1 = inr-natoms_prot;
529 is2 = jnr-natoms_prot;
531 /* Check if at least one of the atoms is a solvent that is not yet
532 * listed for removal, and if both are solvent, that they are not in the
535 if (is1 >= 0 && is2 < 0)
537 mark_res(is1, keep, natoms_solvt, atoms_solvt->atom, &nkeep);
539 else if (is1 < 0 && is2 >= 0)
541 mark_res(is2, keep, natoms_solvt, atoms_solvt->atom, &nkeep);
545 fprintf(stderr, "Keeping %d solvent atoms after proximity check\n",
547 for (i = 0; i < natoms_solvt; i++)
549 remove[i] = remove[i] || !keep[i];
553 /* count how many atoms and residues will be added and make space */
557 jnres = atoms_solvt->nres;
563 for (i = 0; ((i < atoms_solvt->nr) &&
564 ((max_sol == 0) || (jnres < max_sol))); i++)
570 (atoms_solvt->atom[i].resind != atoms_solvt->atom[i-1].resind))
579 fprintf(debug, "Will add %d atoms in %d residues\n", j, jnres);
583 /* Flag the remaing solvent atoms to be removed */
584 jjj = atoms_solvt->atom[i-1].resind;
585 for (; (i < atoms_solvt->nr); i++)
587 if (atoms_solvt->atom[i].resind > jjj)
600 srenew(atoms->resinfo, atoms->nres+jnres);
601 srenew(atoms->atomname, atoms->nr+j);
602 srenew(atoms->atom, atoms->nr+j);
603 srenew(*x, atoms->nr+j);
606 srenew(*v, atoms->nr+j);
608 srenew(*r, atoms->nr+j);
611 /* add the selected atoms_solvt to atoms */
614 resnr = atoms->resinfo[atoms->atom[atoms->nr-1].resind].nr;
622 for (i = 0; i < atoms_solvt->nr; i++)
627 atoms_solvt->atom[i].resind != atoms_solvt->atom[prev].resind)
632 atoms->resinfo[atoms->nres-1] =
633 atoms_solvt->resinfo[atoms_solvt->atom[i].resind];
634 atoms->resinfo[atoms->nres-1].nr = resnr;
635 /* calculate shift of the solvent molecule using the first atom */
636 copy_rvec(x_solvt[i], dx);
637 put_atoms_in_box(ePBC, box, 1, &dx);
638 rvec_dec(dx, x_solvt[i]);
640 atoms->atom[atoms->nr] = atoms_solvt->atom[i];
641 atoms->atomname[atoms->nr] = atoms_solvt->atomname[i];
642 rvec_add(x_solvt[i], dx, (*x)[atoms->nr]);
645 copy_rvec(v_solvt[i], (*v)[atoms->nr]);
647 (*r)[atoms->nr] = r_solvt[i];
648 atoms->atom[atoms->nr].resind = atoms->nres-1;
655 srenew(atoms->resinfo, atoms->nres+nresadd);
660 fprintf(stderr, "Added %d molecules\n", nresadd);
664 done_atom(atoms_all);