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