cc0bf57893280fb33b50ca465e7c3728b8f05eda
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / addconf.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #include "addconf.h"
38
39 #ifdef HAVE_CONFIG_H
40 #include <config.h>
41 #endif
42
43 #include <stdlib.h>
44 #include <string.h>
45
46 #include "gromacs/math/vec.h"
47 #include "macros.h"
48 #include "types/commrec.h"
49 #include "force.h"
50 #include "names.h"
51 #include "nsgrid.h"
52 #include "mdatoms.h"
53 #include "nrnb.h"
54 #include "ns.h"
55 #include "mtop_util.h"
56 #include "chargegroup.h"
57
58 #include "gromacs/pbcutil/pbc.h"
59 #include "gromacs/topology/block.h"
60 #include "gromacs/utility/smalloc.h"
61
62 static real box_margin;
63
64 static real max_dist(rvec *x, real *r, int start, int end)
65 {
66     real maxd;
67     int  i, j;
68
69     maxd = 0;
70     for (i = start; i < end; i++)
71     {
72         for (j = i+1; j < end; j++)
73         {
74             maxd = max(maxd, sqrt(distance2(x[i], x[j]))+0.5*(r[i]+r[j]));
75         }
76     }
77
78     return 0.5*maxd;
79 }
80
81 static gmx_bool outside_box_minus_margin2(rvec x, matrix box)
82 {
83     return ( (x[XX] < 2*box_margin) || (x[XX] > box[XX][XX]-2*box_margin) ||
84              (x[YY] < 2*box_margin) || (x[YY] > box[YY][YY]-2*box_margin) ||
85              (x[ZZ] < 2*box_margin) || (x[ZZ] > box[ZZ][ZZ]-2*box_margin) );
86 }
87
88 static gmx_bool outside_box_plus_margin(rvec x, matrix box)
89 {
90     return ( (x[XX] < -box_margin) || (x[XX] > box[XX][XX]+box_margin) ||
91              (x[YY] < -box_margin) || (x[YY] > box[YY][YY]+box_margin) ||
92              (x[ZZ] < -box_margin) || (x[ZZ] > box[ZZ][ZZ]+box_margin) );
93 }
94
95 static int mark_res(int at, gmx_bool *mark, int natoms, t_atom *atom, int *nmark)
96 {
97     int resind;
98
99     resind = atom[at].resind;
100     while ( (at > 0) && (resind == atom[at-1].resind) )
101     {
102         at--;
103     }
104     while ( (at < natoms) && (resind == atom[at].resind) )
105     {
106         if (!mark[at])
107         {
108             mark[at] = TRUE;
109             (*nmark)++;
110         }
111         at++;
112     }
113
114     return at;
115 }
116
117 static real find_max_real(int n, real radius[])
118 {
119     int  i;
120     real rmax;
121
122     rmax = 0;
123     if (n > 0)
124     {
125         rmax = radius[0];
126         for (i = 1; (i < n); i++)
127         {
128             rmax = max(rmax, radius[i]);
129         }
130     }
131     return rmax;
132 }
133
134 static void combine_atoms(t_atoms *ap, t_atoms *as,
135                           rvec xp[], rvec *vp, rvec xs[], rvec *vs,
136                           t_atoms **a_comb, rvec **x_comb, rvec **v_comb)
137 {
138     t_atoms *ac;
139     rvec    *xc, *vc = NULL;
140     int      i, j, natot, res0;
141
142     /* Total number of atoms */
143     natot = ap->nr+as->nr;
144
145     snew(ac, 1);
146     init_t_atoms(ac, natot, FALSE);
147
148     snew(xc, natot);
149     if (vp && vs)
150     {
151         snew(vc, natot);
152     }
153
154     /* Fill the new structures */
155     for (i = j = 0; (i < ap->nr); i++, j++)
156     {
157         copy_rvec(xp[i], xc[j]);
158         if (vc)
159         {
160             copy_rvec(vp[i], vc[j]);
161         }
162         memcpy(&(ac->atom[j]), &(ap->atom[i]), sizeof(ap->atom[i]));
163         ac->atom[j].type = 0;
164     }
165     res0 = ap->nres;
166     for (i = 0; (i < as->nr); i++, j++)
167     {
168         copy_rvec(xs[i], xc[j]);
169         if (vc)
170         {
171             copy_rvec(vs[i], vc[j]);
172         }
173         memcpy(&(ac->atom[j]), &(as->atom[i]), sizeof(as->atom[i]));
174         ac->atom[j].type    = 0;
175         ac->atom[j].resind += res0;
176     }
177     ac->nr   = j;
178     ac->nres = ac->atom[j-1].resind+1;
179     /* Fill all elements to prevent uninitialized memory */
180     for (i = 0; i < ac->nr; i++)
181     {
182         ac->atom[i].m     = 1;
183         ac->atom[i].q     = 0;
184         ac->atom[i].mB    = 1;
185         ac->atom[i].qB    = 0;
186         ac->atom[i].type  = 0;
187         ac->atom[i].typeB = 0;
188         ac->atom[i].ptype = eptAtom;
189     }
190
191     /* Return values */
192     *a_comb = ac;
193     *x_comb = xc;
194     *v_comb = vc;
195 }
196
197 static t_forcerec *fr = NULL;
198
199 static void do_nsgrid(FILE *fp, gmx_bool bVerbose,
200                       matrix box, rvec x[], t_atoms *atoms, real rlong,
201                       const output_env_t oenv)
202 {
203     gmx_mtop_t     *mtop;
204     gmx_localtop_t *top;
205     t_mdatoms      *md;
206     t_block        *cgs;
207     t_inputrec     *ir;
208     t_nrnb          nrnb;
209     t_commrec      *cr;
210     int            *cg_index;
211     gmx_moltype_t  *molt;
212     gmx_ffparams_t *ffp;
213     ivec           *nFreeze;
214     int             i, m, natoms;
215     rvec            box_size;
216     real           *lambda, *dvdl;
217
218     natoms = atoms->nr;
219
220     /* Charge group index */
221     snew(cg_index, natoms);
222     for (i = 0; (i < natoms); i++)
223     {
224         cg_index[i] = i;
225     }
226
227     /* Topology needs charge groups and exclusions */
228     snew(mtop, 1);
229     init_mtop(mtop);
230     mtop->natoms = natoms;
231     /* Make one moltype that contains the whol system */
232     mtop->nmoltype = 1;
233     snew(mtop->moltype, mtop->nmoltype);
234     molt        = &mtop->moltype[0];
235     molt->name  = mtop->name;
236     molt->atoms = *atoms;
237     stupid_fill_block(&molt->cgs, mtop->natoms, FALSE);
238     stupid_fill_blocka(&molt->excls, natoms);
239     /* Make one molblock for the whole system */
240     mtop->nmolblock = 1;
241     snew(mtop->molblock, mtop->nmolblock);
242     mtop->molblock[0].type       = 0;
243     mtop->molblock[0].nmol       = 1;
244     mtop->molblock[0].natoms_mol = natoms;
245     /* Initialize a single energy group */
246     mtop->groups.grps[egcENER].nr = 1;
247     mtop->groups.ngrpnr[egcENER]  = 0;
248     mtop->groups.grpnr[egcENER]   = NULL;
249
250     ffp = &mtop->ffparams;
251
252     ffp->ntypes = 1;
253     ffp->atnr   = 1;
254     ffp->reppow = 12;
255     snew(ffp->functype, 1);
256     snew(ffp->iparams, 1);
257     ffp->iparams[0].lj.c6  = 1;
258     ffp->iparams[0].lj.c12 = 1;
259
260     /* inputrec structure */
261     snew(ir, 1);
262     ir->cutoff_scheme    = ecutsGROUP;
263     ir->coulomb_modifier = eintmodNONE;
264     ir->vdw_modifier     = eintmodNONE;
265     ir->coulombtype      = eelCUT;
266     ir->vdwtype          = evdwCUT;
267     ir->ndelta           = 2;
268     ir->ns_type          = ensGRID;
269     snew(ir->opts.egp_flags, 1);
270
271     top = gmx_mtop_generate_local_top(mtop, ir);
272
273     /* Some nasty shortcuts */
274     cgs  = &(top->cgs);
275
276     /* mdatoms structure */
277     snew(nFreeze, 2);
278     snew(md, 1);
279     md = init_mdatoms(fp, mtop, FALSE);
280     atoms2md(mtop, ir, 0, NULL, mtop->natoms, md);
281     sfree(nFreeze);
282
283     /* forcerec structure */
284     if (fr == NULL)
285     {
286         fr = mk_forcerec();
287     }
288     snew(cr, 1);
289     cr->nnodes   = 1;
290     /* cr->nthreads = 1; */
291
292     /*    ir->rlist       = ir->rcoulomb = ir->rvdw = rlong;
293        printf("Neighborsearching with a cut-off of %g\n",rlong);
294        init_forcerec(stdout,fr,ir,top,cr,md,box,FALSE,NULL,NULL,NULL,TRUE);*/
295     fr->cg0     = 0;
296     fr->hcg     = top->cgs.nr;
297     fr->nWatMol = 0;
298
299     /* Prepare for neighboursearching */
300     init_nrnb(&nrnb);
301
302     /* Init things dependent on parameters */
303     ir->rlistlong = ir->rlist = ir->rcoulomb = ir->rvdw = rlong;
304     /* create free energy data to avoid NULLs */
305     snew(ir->fepvals, 1);
306     printf("Neighborsearching with a cut-off of %g\n", rlong);
307     init_forcerec(stdout, oenv, fr, NULL, ir, mtop, cr, box,
308                   NULL, NULL, NULL, NULL, NULL, TRUE, -1);
309     if (debug)
310     {
311         pr_forcerec(debug, fr);
312     }
313
314     /* Calculate new stuff dependent on coords and box */
315     for (m = 0; (m < DIM); m++)
316     {
317         box_size[m] = box[m][m];
318     }
319     calc_shifts(box, fr->shift_vec);
320     put_charge_groups_in_box(fp, 0, cgs->nr, fr->ePBC, box, cgs, x, fr->cg_cm);
321
322     /* Do the actual neighboursearching */
323     snew(lambda, efptNR);
324     snew(dvdl, efptNR);
325     init_neighbor_list(fp, fr, md->homenr);
326     search_neighbours(fp, fr, box, top,
327                       &mtop->groups, cr, &nrnb, md, TRUE, FALSE);
328
329     if (debug)
330     {
331         dump_nblist(debug, cr, fr, 0);
332     }
333
334     if (bVerbose)
335     {
336         fprintf(stderr, "Successfully made neighbourlist\n");
337     }
338 }
339
340 static gmx_bool bXor(gmx_bool b1, gmx_bool b2)
341 {
342     return (b1 && !b2) || (b2 && !b1);
343 }
344
345 void add_conf(t_atoms *atoms, rvec **x, rvec **v, real **r, gmx_bool bSrenew,
346               int ePBC, matrix box, gmx_bool bInsert,
347               t_atoms *atoms_solvt, rvec *x_solvt, rvec *v_solvt, real *r_solvt,
348               gmx_bool bVerbose, real rshell, int max_sol, const output_env_t oenv)
349 {
350     t_nblist       *nlist;
351     t_atoms        *atoms_all;
352     real            max_vdw, *r_prot, *r_all, n2, r2, ib1, ib2;
353     int             natoms_prot, natoms_solvt;
354     int             i, j, jj, m, j0, j1, jjj, jnres, jnr, inr, iprot, is1, is2;
355     int             prev, resnr, nresadd, d, k, ncells, maxincell;
356     int             dx0, dx1, dy0, dy1, dz0, dz1;
357     int             ntest, nremove, nkeep;
358     rvec            dx, xi, xj, xpp, *x_all, *v_all;
359     gmx_bool       *remove, *keep;
360     int             bSolSol;
361
362     natoms_prot  = atoms->nr;
363     natoms_solvt = atoms_solvt->nr;
364     if (natoms_solvt <= 0)
365     {
366         fprintf(stderr, "WARNING: Nothing to add\n");
367         return;
368     }
369
370     if (ePBC == epbcSCREW)
371     {
372         gmx_fatal(FARGS, "Sorry, %s pbc is not yet supported", epbc_names[ePBC]);
373     }
374
375     if (bVerbose)
376     {
377         fprintf(stderr, "Calculating Overlap...\n");
378     }
379
380     /* Set margin around box edges to largest solvent dimension.
381      * The maximum distance between atoms in a solvent molecule should
382      * be calculated. At the moment a fudge factor of 3 is used.
383      */
384     r_prot     = *r;
385     box_margin = 3*find_max_real(natoms_solvt, r_solvt);
386     max_vdw    = max(3*find_max_real(natoms_prot, r_prot), box_margin);
387     fprintf(stderr, "box_margin = %g\n", box_margin);
388
389     snew(remove, natoms_solvt);
390
391     nremove = 0;
392     if (!bInsert)
393     {
394         for (i = 0; i < atoms_solvt->nr; i++)
395         {
396             if (outside_box_plus_margin(x_solvt[i], box) )
397             {
398                 i = mark_res(i, remove, atoms_solvt->nr, atoms_solvt->atom, &nremove);
399             }
400         }
401         fprintf(stderr, "Removed %d atoms that were outside the box\n", nremove);
402     }
403
404     /* Define grid stuff */
405     /* Largest VDW radius */
406     snew(r_all, natoms_prot+natoms_solvt);
407     for (i = j = 0; i < natoms_prot; i++, j++)
408     {
409         r_all[j] = r_prot[i];
410     }
411     for (i = 0; i < natoms_solvt; i++, j++)
412     {
413         r_all[j] = r_solvt[i];
414     }
415
416     /* Combine arrays */
417     combine_atoms(atoms, atoms_solvt, *x, v ? *v : NULL, x_solvt, v_solvt,
418                   &atoms_all, &x_all, &v_all);
419
420     /* Do neighboursearching step */
421     do_nsgrid(stdout, bVerbose, box, x_all, atoms_all, max_vdw, oenv);
422
423     /* check solvent with solute */
424     nlist = &(fr->nblists[0].nlist_sr[eNL_VDW]);
425     fprintf(stderr, "nri = %d, nrj = %d\n", nlist->nri, nlist->nrj);
426     for (bSolSol = 0; (bSolSol <= (bInsert ? 0 : 1)); bSolSol++)
427     {
428         ntest = nremove = 0;
429         fprintf(stderr, "Checking %s-Solvent overlap:",
430                 bSolSol ? "Solvent" : "Protein");
431         for (i = 0; (i < nlist->nri && nremove < natoms_solvt); i++)
432         {
433             inr = nlist->iinr[i];
434             j0  = nlist->jindex[i];
435             j1  = nlist->jindex[i+1];
436             rvec_add(x_all[inr], fr->shift_vec[nlist->shift[i]], xi);
437
438             for (j = j0; (j < j1 && nremove < natoms_solvt); j++)
439             {
440                 jnr = nlist->jjnr[j];
441                 copy_rvec(x_all[jnr], xj);
442
443                 /* Check solvent-protein and solvent-solvent */
444                 is1 = inr-natoms_prot;
445                 is2 = jnr-natoms_prot;
446
447                 /* Check if at least one of the atoms is a solvent that is not yet
448                  * listed for removal, and if both are solvent, that they are not in the
449                  * same residue.
450                  */
451                 if ((!bSolSol &&
452                      bXor((is1 >= 0), (is2 >= 0)) && /* One atom is protein */
453                      ((is1 < 0) || ((is1 >= 0) && !remove[is1])) &&
454                      ((is2 < 0) || ((is2 >= 0) && !remove[is2]))) ||
455
456                     (bSolSol  &&
457                      (is1 >= 0) && (!remove[is1]) &&                   /* is1 is solvent */
458                      (is2 >= 0) && (!remove[is2]) &&                   /* is2 is solvent */
459                      (bInsert ||                                       /* when inserting also check inside the box */
460                       (outside_box_minus_margin2(x_solvt[is1], box) && /* is1 on edge */
461                        outside_box_minus_margin2(x_solvt[is2], box))   /* is2 on edge */
462                      ) &&
463                      (atoms_solvt->atom[is1].resind !=                 /* Not the same residue */
464                       atoms_solvt->atom[is2].resind)))
465                 {
466
467                     ntest++;
468                     rvec_sub(xi, xj, dx);
469                     n2 = norm2(dx);
470                     r2 = sqr(r_all[inr]+r_all[jnr]);
471                     if (n2 < r2)
472                     {
473                         if (bInsert)
474                         {
475                             nremove = natoms_solvt;
476                             for (k = 0; k < nremove; k++)
477                             {
478                                 remove[k] = TRUE;
479                             }
480                         }
481                         /* Need only remove one of the solvents... */
482                         if (is2 >= 0)
483                         {
484                             (void) mark_res(is2, remove, natoms_solvt, atoms_solvt->atom,
485                                             &nremove);
486                         }
487                         else if (is1 >= 0)
488                         {
489                             (void) mark_res(is1, remove, natoms_solvt, atoms_solvt->atom,
490                                             &nremove);
491                         }
492                         else
493                         {
494                             fprintf(stderr, "Neither atom is solvent%d %d\n", is1, is2);
495                         }
496                     }
497                 }
498             }
499         }
500         if (!bInsert)
501         {
502             fprintf(stderr, " tested %d pairs, removed %d atoms.\n", ntest, nremove);
503         }
504     }
505     if (debug)
506     {
507         for (i = 0; i < natoms_solvt; i++)
508         {
509             fprintf(debug, "remove[%5d] = %s\n", i, bool_names[remove[i]]);
510         }
511     }
512
513     /* Search again, now with another cut-off */
514     if (rshell > 0)
515     {
516         do_nsgrid(stdout, bVerbose, box, x_all, atoms_all, rshell, oenv);
517         nlist = &(fr->nblists[0].nlist_sr[eNL_VDW]);
518         fprintf(stderr, "nri = %d, nrj = %d\n", nlist->nri, nlist->nrj);
519         nkeep = 0;
520         snew(keep, natoms_solvt);
521         for (i = 0; i < nlist->nri; i++)
522         {
523             inr = nlist->iinr[i];
524             j0  = nlist->jindex[i];
525             j1  = nlist->jindex[i+1];
526
527             for (j = j0; j < j1; j++)
528             {
529                 jnr = nlist->jjnr[j];
530
531                 /* Check solvent-protein and solvent-solvent */
532                 is1 = inr-natoms_prot;
533                 is2 = jnr-natoms_prot;
534
535                 /* Check if at least one of the atoms is a solvent that is not yet
536                  * listed for removal, and if both are solvent, that they are not in the
537                  * same residue.
538                  */
539                 if (is1 >= 0 && is2 < 0)
540                 {
541                     mark_res(is1, keep, natoms_solvt, atoms_solvt->atom, &nkeep);
542                 }
543                 else if (is1 < 0 && is2 >= 0)
544                 {
545                     mark_res(is2, keep, natoms_solvt, atoms_solvt->atom, &nkeep);
546                 }
547             }
548         }
549         fprintf(stderr, "Keeping %d solvent atoms after proximity check\n",
550                 nkeep);
551         for (i = 0; i < natoms_solvt; i++)
552         {
553             remove[i] = remove[i] || !keep[i];
554         }
555         sfree(keep);
556     }
557     /* count how many atoms and residues will be added and make space */
558     if (bInsert)
559     {
560         j     = atoms_solvt->nr;
561         jnres = atoms_solvt->nres;
562     }
563     else
564     {
565         j     = 0;
566         jnres = 0;
567         for (i = 0; ((i < atoms_solvt->nr) &&
568                      ((max_sol == 0) || (jnres < max_sol))); i++)
569         {
570             if (!remove[i])
571             {
572                 j++;
573                 if ((i == 0) ||
574                     (atoms_solvt->atom[i].resind != atoms_solvt->atom[i-1].resind))
575                 {
576                     jnres++;
577                 }
578             }
579         }
580     }
581     if (debug)
582     {
583         fprintf(debug, "Will add %d atoms in %d residues\n", j, jnres);
584     }
585     if (!bInsert)
586     {
587         /* Flag the remaing solvent atoms to be removed */
588         jjj = atoms_solvt->atom[i-1].resind;
589         for (; (i < atoms_solvt->nr); i++)
590         {
591             if (atoms_solvt->atom[i].resind > jjj)
592             {
593                 remove[i] = TRUE;
594             }
595             else
596             {
597                 j++;
598             }
599         }
600     }
601
602     if (bSrenew)
603     {
604         srenew(atoms->resinfo,  atoms->nres+jnres);
605         srenew(atoms->atomname, atoms->nr+j);
606         srenew(atoms->atom,     atoms->nr+j);
607         srenew(*x,              atoms->nr+j);
608         if (v)
609         {
610             srenew(*v,       atoms->nr+j);
611         }
612         srenew(*r,              atoms->nr+j);
613     }
614
615     /* add the selected atoms_solvt to atoms */
616     if (atoms->nr > 0)
617     {
618         resnr = atoms->resinfo[atoms->atom[atoms->nr-1].resind].nr;
619     }
620     else
621     {
622         resnr = 0;
623     }
624     prev    = -1;
625     nresadd = 0;
626     for (i = 0; i < atoms_solvt->nr; i++)
627     {
628         if (!remove[i])
629         {
630             if (prev == -1 ||
631                 atoms_solvt->atom[i].resind != atoms_solvt->atom[prev].resind)
632             {
633                 nresadd++;
634                 atoms->nres++;
635                 resnr++;
636                 atoms->resinfo[atoms->nres-1] =
637                     atoms_solvt->resinfo[atoms_solvt->atom[i].resind];
638                 atoms->resinfo[atoms->nres-1].nr = resnr;
639                 /* calculate shift of the solvent molecule using the first atom */
640                 copy_rvec(x_solvt[i], dx);
641                 put_atoms_in_box(ePBC, box, 1, &dx);
642                 rvec_dec(dx, x_solvt[i]);
643             }
644             atoms->atom[atoms->nr]     = atoms_solvt->atom[i];
645             atoms->atomname[atoms->nr] = atoms_solvt->atomname[i];
646             rvec_add(x_solvt[i], dx, (*x)[atoms->nr]);
647             if (v)
648             {
649                 copy_rvec(v_solvt[i], (*v)[atoms->nr]);
650             }
651             (*r)[atoms->nr]               = r_solvt[i];
652             atoms->atom[atoms->nr].resind = atoms->nres-1;
653             atoms->nr++;
654             prev = i;
655         }
656     }
657     if (bSrenew)
658     {
659         srenew(atoms->resinfo,  atoms->nres+nresadd);
660     }
661
662     if (bVerbose)
663     {
664         fprintf(stderr, "Added %d molecules\n", nresadd);
665     }
666
667     sfree(remove);
668     done_atom(atoms_all);
669     sfree(x_all);
670     sfree(v_all);
671 }