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