Merge branch release-5-1
[alexxy/gromacs.git] / src / gromacs / trajectoryanalysis / modules / sasa.cpp
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-2006, The GROMACS development team.
6  * Copyright (c) 2008,2009,2010,2011,2012,2013,2014,2015, 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 /*! \internal \file
38  * \brief
39  * Implements gmx::analysismodules::Sasa.
40  *
41  * \author Teemu Murtola <teemu.murtola@gmail.com> (C++ conversion)
42  * \ingroup module_trajectoryanalysis
43  */
44 #include "gmxpre.h"
45
46 #include "sasa.h"
47
48 #include <algorithm>
49 #include <string>
50 #include <vector>
51
52 #include "gromacs/analysisdata/analysisdata.h"
53 #include "gromacs/analysisdata/modules/average.h"
54 #include "gromacs/analysisdata/modules/plot.h"
55 #include "gromacs/fileio/confio.h"
56 #include "gromacs/fileio/pdbio.h"
57 #include "gromacs/fileio/trx.h"
58 #include "gromacs/legacyheaders/copyrite.h"
59 #include "gromacs/math/units.h"
60 #include "gromacs/math/vec.h"
61 #include "gromacs/options/basicoptions.h"
62 #include "gromacs/options/filenameoption.h"
63 #include "gromacs/options/options.h"
64 #include "gromacs/pbcutil/pbc.h"
65 #include "gromacs/selection/selection.h"
66 #include "gromacs/selection/selectionoption.h"
67 #include "gromacs/topology/atomprop.h"
68 #include "gromacs/topology/symtab.h"
69 #include "gromacs/topology/topology.h"
70 #include "gromacs/trajectoryanalysis/analysismodule.h"
71 #include "gromacs/trajectoryanalysis/analysissettings.h"
72 #include "gromacs/utility/exceptions.h"
73 #include "gromacs/utility/futil.h"
74 #include "gromacs/utility/scoped_cptr.h"
75 #include "gromacs/utility/smalloc.h"
76 #include "gromacs/utility/stringutil.h"
77
78 #include "surfacearea.h"
79
80 namespace gmx
81 {
82
83 namespace analysismodules
84 {
85
86 namespace
87 {
88
89 //! \addtogroup module_trajectoryanalysis
90 //! \{
91
92 //! Tracks information on two nearest neighbors of a single surface dot.
93 struct t_conect
94 {
95     //! Index of the second nearest neighbor dot.
96     atom_id  aa;
97     //! Index of the nearest neighbor dot.
98     atom_id  ab;
99     //! Squared distance to `aa`.
100     real     d2a;
101     //! Squared distance to `ab`.
102     real     d2b;
103 };
104
105 /*! \brief
106  * Updates nearest neighbor information for a surface dot.
107  *
108  * \param[in,out] c  Nearest neighbor information array to update.
109  * \param[in]     i  Index in `c` to update.
110  * \param[in]     j  Index of the other surface dot to add to the array.
111  * \param[in]     d2 Squared distance between `i` and `j`.
112  */
113 void add_rec(t_conect c[], atom_id i, atom_id j, real d2)
114 {
115     if (c[i].aa == NO_ATID)
116     {
117         c[i].aa  = j;
118         c[i].d2a = d2;
119     }
120     else if (c[i].ab == NO_ATID)
121     {
122         c[i].ab  = j;
123         c[i].d2b = d2;
124     }
125     else if (d2 < c[i].d2a)
126     {
127         c[i].aa  = j;
128         c[i].d2a = d2;
129     }
130     else if (d2 < c[i].d2b)
131     {
132         c[i].ab  = j;
133         c[i].d2b = d2;
134     }
135     /* Swap them if necessary: a must be larger than b */
136     if (c[i].d2a < c[i].d2b)
137     {
138         j        = c[i].ab;
139         c[i].ab  = c[i].aa;
140         c[i].aa  = j;
141         d2       = c[i].d2b;
142         c[i].d2b = c[i].d2a;
143         c[i].d2a = d2;
144     }
145 }
146
147 /*! \brief
148  * Adds CONECT records for surface dots.
149  *
150  * \param[in] fn  PDB file to append the CONECT records to.
151  * \param[in] n   Number of dots in `x`.
152  * \param[in] x   Array of surface dot positions.
153  *
154  * Adds a CONECT record that connects each surface dot to its two nearest
155  * neighbors.  The function is copied verbatim from the old gmx_sas.c
156  * implementation.
157  */
158 void do_conect(const char *fn, int n, rvec x[])
159 {
160     FILE     *fp;
161     int       i, j;
162     t_conect *c;
163     rvec      dx;
164     real      d2;
165
166     fprintf(stderr, "Building CONECT records\n");
167     snew(c, n);
168     for (i = 0; (i < n); i++)
169     {
170         c[i].aa = c[i].ab = NO_ATID;
171     }
172
173     for (i = 0; (i < n); i++)
174     {
175         for (j = i+1; (j < n); j++)
176         {
177             rvec_sub(x[i], x[j], dx);
178             d2 = iprod(dx, dx);
179             add_rec(c, i, j, d2);
180             add_rec(c, j, i, d2);
181         }
182     }
183     fp = gmx_ffopen(fn, "a");
184     for (i = 0; (i < n); i++)
185     {
186         if ((c[i].aa == NO_ATID) || (c[i].ab == NO_ATID))
187         {
188             fprintf(stderr, "Warning dot %d has no conections\n", i+1);
189         }
190         fprintf(fp, "CONECT%5d%5d%5d\n", i+1, c[i].aa+1, c[i].ab+1);
191     }
192     gmx_ffclose(fp);
193     sfree(c);
194 }
195
196 /*! \brief
197  * Plots the surface into a PDB file, optionally including the original atoms.
198  */
199 void connolly_plot(const char *fn, int ndots, real dots[], rvec x[], t_atoms *atoms,
200                    t_symtab *symtab, int ePBC, const matrix box, gmx_bool bIncludeSolute)
201 {
202     const char *const  atomnm = "DOT";
203     const char *const  resnm  = "DOT";
204     const char *const  title  = "Connolly Dot Surface Generated by gmx sasa";
205
206     int                i, i0, r0, ii0, k;
207     rvec              *xnew;
208     t_atoms            aaa;
209
210     if (bIncludeSolute)
211     {
212         i0 = atoms->nr;
213         r0 = atoms->nres;
214         srenew(atoms->atom, atoms->nr+ndots);
215         memset(&atoms->atom[i0], 0, sizeof(*atoms->atom)*ndots);
216         srenew(atoms->atomname, atoms->nr+ndots);
217         srenew(atoms->resinfo, r0+1);
218         atoms->atom[i0].resind = r0;
219         t_atoms_set_resinfo(atoms, i0, symtab, resnm, r0+1, ' ', 0, ' ');
220         if (atoms->pdbinfo != NULL)
221         {
222             srenew(atoms->pdbinfo, atoms->nr+ndots);
223         }
224         snew(xnew, atoms->nr+ndots);
225         for (i = 0; (i < atoms->nr); i++)
226         {
227             copy_rvec(x[i], xnew[i]);
228         }
229         for (i = k = 0; (i < ndots); i++)
230         {
231             ii0                        = i0+i;
232             atoms->atomname[ii0]       = put_symtab(symtab, atomnm);
233             atoms->atom[ii0].resind    = r0;
234             xnew[ii0][XX]              = dots[k++];
235             xnew[ii0][YY]              = dots[k++];
236             xnew[ii0][ZZ]              = dots[k++];
237             if (atoms->pdbinfo != NULL)
238             {
239                 atoms->pdbinfo[ii0].type   = epdbATOM;
240                 atoms->pdbinfo[ii0].atomnr = ii0;
241                 atoms->pdbinfo[ii0].bfac   = 0.0;
242                 atoms->pdbinfo[ii0].occup  = 0.0;
243             }
244         }
245         atoms->nr   = i0+ndots;
246         atoms->nres = r0+1;
247         write_sto_conf(fn, title, atoms, xnew, NULL, ePBC, const_cast<rvec *>(box));
248         atoms->nres = r0;
249         atoms->nr   = i0;
250     }
251     else
252     {
253         init_t_atoms(&aaa, ndots, TRUE);
254         aaa.atom[0].resind = 0;
255         t_atoms_set_resinfo(&aaa, 0, symtab, resnm, 1, ' ', 0, ' ');
256         snew(xnew, ndots);
257         for (i = k = 0; (i < ndots); i++)
258         {
259             ii0                     = i;
260             aaa.atomname[ii0]       = put_symtab(symtab, atomnm);
261             aaa.pdbinfo[ii0].type   = epdbATOM;
262             aaa.pdbinfo[ii0].atomnr = ii0;
263             aaa.atom[ii0].resind    = 0;
264             xnew[ii0][XX]           = dots[k++];
265             xnew[ii0][YY]           = dots[k++];
266             xnew[ii0][ZZ]           = dots[k++];
267             aaa.pdbinfo[ii0].bfac   = 0.0;
268             aaa.pdbinfo[ii0].occup  = 0.0;
269         }
270         aaa.nr = ndots;
271         write_sto_conf(fn, title, &aaa, xnew, NULL, ePBC, const_cast<rvec *>(box));
272         do_conect(fn, ndots, xnew);
273         free_t_atoms(&aaa, FALSE);
274     }
275     sfree(xnew);
276 }
277
278 /********************************************************************
279  * Actual analysis module
280  */
281
282 /*! \brief
283  * Implements `gmx sas` trajectory analysis module.
284  */
285 class Sasa : public TrajectoryAnalysisModule
286 {
287     public:
288         Sasa();
289
290         virtual void initOptions(Options                    *options,
291                                  TrajectoryAnalysisSettings *settings);
292         virtual void initAnalysis(const TrajectoryAnalysisSettings &settings,
293                                   const TopologyInformation        &top);
294
295         virtual TrajectoryAnalysisModuleDataPointer startFrames(
296             const AnalysisDataParallelOptions &opt,
297             const SelectionCollection         &selections);
298         virtual void analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc,
299                                   TrajectoryAnalysisModuleData *pdata);
300
301         virtual void finishAnalysis(int nframes);
302         virtual void writeOutput();
303
304     private:
305         /*! \brief
306          * Surface areas as a function of time.
307          *
308          * First column is for the calculation group, and the rest for the
309          * output groups.  This data is always produced.
310          */
311         AnalysisData            area_;
312         /*! \brief
313          * Per-atom surface areas as a function of time.
314          *
315          * Contains one data set for each column in `area_`.
316          * Each column corresponds to a selection position in `surfaceSel_`.
317          * This data is only produced if atom or residue areas have been
318          * requested.
319          */
320         AnalysisData            atomArea_;
321         /*! \brief
322          * Per-residue surface areas as a function of time.
323          *
324          * Contains one data set for each column in `area_`.
325          * Each column corresponds to a distinct residue `surfaceSel_`.
326          * For example, if `surfaceSel_` selects residues 2, 5, and 7, there
327          * will be three columns here.
328          * This data is only produced if atom or residue areas have been
329          * requested.
330          */
331         AnalysisData            residueArea_;
332         /*! \brief
333          * Free energy estimates as a function of time.
334          *
335          * Column layout is the same as for `area_`.
336          * This data is only produced if the output is requested.
337          */
338         AnalysisData            dgSolv_;
339         /*! \brief
340          * Total volume and density of the calculation group as a function of
341          * time.
342          *
343          * The first column is the volume and the second column is the density.
344          * This data is only produced if the output is requested.
345          */
346         AnalysisData            volume_;
347
348         /*! \brief
349          * The selection to calculate the surface for.
350          *
351          * Selection::originalId() and Selection::mappedId() store the mapping
352          * from the positions to the columns of `residueArea_`.
353          * The selection is computed with SelectionOption::dynamicMask(), i.e.,
354          * even in the presence of a dynamic selection, the number of returned
355          * positions is fixed, and SelectionPosition::selected() is used.
356          */
357         Selection               surfaceSel_;
358         /*! \brief
359          * List of optional additional output groups.
360          *
361          * Each of these must be a subset of the `surfaceSel_`.
362          * Selection::originalId() and Selection::mappedId() store the mapping
363          * from the positions to the corresponsing positions in `surfaceSel_`.
364          */
365         SelectionList           outputSel_;
366
367         std::string             fnArea_;
368         std::string             fnAtomArea_;
369         std::string             fnResidueArea_;
370         std::string             fnDGSolv_;
371         std::string             fnVolume_;
372         std::string             fnConnolly_;
373
374         double                  solsize_;
375         int                     ndots_;
376         //double                  minarea_;
377         double                  dgsDefault_;
378         bool                    bIncludeSolute_;
379
380         t_topology             *top_;
381         //! Combined VdW and probe radii for each atom in the calculation group.
382         std::vector<real>       radii_;
383         /*! \brief
384          * Solvation free energy coefficients for each atom in the calculation
385          * group.
386          *
387          * Empty if the free energy output has not been requested.
388          */
389         std::vector<real>       dgsFactor_;
390         //! Calculation algorithm.
391         SurfaceAreaCalculator   calculator_;
392
393         // Copy and assign disallowed by base.
394 };
395
396 Sasa::Sasa()
397     : TrajectoryAnalysisModule(SasaInfo::name, SasaInfo::shortDescription),
398       solsize_(0.14), ndots_(24), dgsDefault_(0), bIncludeSolute_(true), top_(NULL)
399 {
400     //minarea_ = 0.5;
401     registerAnalysisDataset(&area_, "area");
402     registerAnalysisDataset(&atomArea_, "atomarea");
403     registerAnalysisDataset(&residueArea_, "resarea");
404     registerAnalysisDataset(&dgSolv_, "dgsolv");
405     registerAnalysisDataset(&volume_, "volume");
406 }
407
408 void
409 Sasa::initOptions(Options *options, TrajectoryAnalysisSettings *settings)
410 {
411     static const char *const desc[] = {
412         "[THISMODULE] computes solvent accessible surface areas.",
413         "See Eisenhaber F, Lijnzaad P, Argos P, Sander C, & Scharf M",
414         "(1995) J. Comput. Chem. 16, 273-284 for the algorithm used.",
415         "With [TT]-q[tt], the Connolly surface can be generated as well",
416         "in a [REF].pdb[ref] file where the nodes are represented as atoms",
417         "and the edges connecting the nearest nodes as CONECT records.",
418         "[TT]-odg[tt] allows for estimation of solvation free energies",
419         "from per-atom solvation energies per exposed surface area.[PAR]",
420
421         "The program requires a selection for the surface calculation to be",
422         "specified with [TT]-surface[tt]. This should always consist of all",
423         "non-solvent atoms in the system. The area of this group is always",
424         "calculated. Optionally, [TT]-output[tt] can specify additional",
425         "selections, which should be subsets of the calculation group.",
426         "The solvent-accessible areas for these groups are also extracted",
427         "from the full surface.[PAR]",
428
429         "The average and standard deviation of the area over the trajectory",
430         "can be calculated per residue and atom (options [TT]-or[tt] and",
431         "[TT]-oa[tt]).[PAR]",
432         //"In combination with the latter option an [REF].itp[ref] file can be",
433         //"generated (option [TT]-i[tt])",
434         //"which can be used to restrain surface atoms.[PAR]",
435
436         "With the [TT]-tv[tt] option the total volume and density of the",
437         "molecule can be computed. With [TT]-pbc[tt] (the default), you",
438         "must ensure that your molecule/surface group is not split across PBC.",
439         "Otherwise, you will get non-sensical results.",
440         "Please also consider whether the normal probe radius is appropriate",
441         "in this case or whether you would rather use, e.g., 0. It is good",
442         "to keep in mind that the results for volume and density are very",
443         "approximate. For example, in ice Ih, one can easily fit water molecules in the",
444         "pores which would yield a volume that is too low, and surface area and density",
445         "that are both too high."
446     };
447
448     settings->setHelpText(desc);
449
450     options->addOption(FileNameOption("o").filetype(eftPlot).outputFile().required()
451                            .store(&fnArea_).defaultBasename("area")
452                            .description("Total area as a function of time"));
453     options->addOption(FileNameOption("odg").filetype(eftPlot).outputFile()
454                            .store(&fnDGSolv_).defaultBasename("dgsolv")
455                            .description("Estimated solvation free energy as a function of time"));
456     options->addOption(FileNameOption("or").filetype(eftPlot).outputFile()
457                            .store(&fnResidueArea_).defaultBasename("resarea")
458                            .description("Average area per residue"));
459     options->addOption(FileNameOption("oa").filetype(eftPlot).outputFile()
460                            .store(&fnAtomArea_).defaultBasename("atomarea")
461                            .description("Average area per atom"));
462     options->addOption(FileNameOption("tv").filetype(eftPlot).outputFile()
463                            .store(&fnVolume_).defaultBasename("volume")
464                            .description("Total volume and density as a function of time"));
465     options->addOption(FileNameOption("q").filetype(eftPDB).outputFile()
466                            .store(&fnConnolly_).defaultBasename("connolly")
467                            .description("PDB file for Connolly surface"));
468     //options->addOption(FileNameOption("i").filetype(eftITP).outputFile()
469     //                       .store(&fnRestraints_).defaultBasename("surfat")
470     //                       .description("Topology file for position restraings on surface atoms"));
471
472
473     options->addOption(DoubleOption("probe").store(&solsize_)
474                            .description("Radius of the solvent probe (nm)"));
475     options->addOption(IntegerOption("ndots").store(&ndots_)
476                            .description("Number of dots per sphere, more dots means more accuracy"));
477     //options->addOption(DoubleOption("minarea").store(&minarea_)
478     //                       .description("The minimum area (nm^2) to count an atom as a surface atom when writing a position restraint file (see help)"));
479     options->addOption(BooleanOption("prot").store(&bIncludeSolute_)
480                            .description("Output the protein to the Connolly [REF].pdb[ref] file too"));
481     options->addOption(DoubleOption("dgs").store(&dgsDefault_)
482                            .description("Default value for solvation free energy per area (kJ/mol/nm^2)"));
483
484     // Selections must select atoms for the VdW radii lookup to work.
485     // The calculation group uses dynamicMask() so that the coordinates
486     // match a static array of VdW radii.
487     options->addOption(SelectionOption("surface").store(&surfaceSel_)
488                            .required().onlyAtoms().dynamicMask()
489                            .description("Surface calculation selection"));
490     options->addOption(SelectionOption("output").storeVector(&outputSel_)
491                            .onlyAtoms().multiValue()
492                            .description("Output selection(s)"));
493
494     // Atom names etc. are required for the VdW radii lookup.
495     settings->setFlag(TrajectoryAnalysisSettings::efRequireTop);
496 }
497
498 void
499 Sasa::initAnalysis(const TrajectoryAnalysisSettings &settings,
500                    const TopologyInformation        &top)
501 {
502     const t_atoms &atoms = top.topology()->atoms;
503     top_ = top.topology();
504
505     //bITP   = opt2bSet("-i", nfile, fnm);
506     const bool bResAt =
507         !fnResidueArea_.empty() || !fnAtomArea_.empty(); // || bITP;
508     const bool bDGsol = !fnDGSolv_.empty();
509
510     if (solsize_ < 0)
511     {
512         solsize_ = 1e-3;
513         fprintf(stderr, "Probe size too small, setting it to %g\n", solsize_);
514     }
515     if (ndots_ < 20)
516     {
517         ndots_ = 20;
518         fprintf(stderr, "Ndots too small, setting it to %d\n", ndots_);
519     }
520
521     please_cite(stderr, "Eisenhaber95");
522     //if ((top.ePBC() != epbcXYZ) || (TRICLINIC(fr.box)))
523     //{
524     //    fprintf(stderr, "\n\nWARNING: non-rectangular boxes may give erroneous results or crashes.\n"
525     //            "Analysis based on vacuum simulations (with the possibility of evaporation)\n"
526     //            "will certainly crash the analysis.\n\n");
527     //}
528
529     if (bDGsol)
530     {
531         if (!top.hasFullTopology())
532         {
533             GMX_THROW(InconsistentInputError("Cannot compute Delta G of solvation without a tpr file"));
534         }
535         else
536         {
537             if (strcmp(*(atoms.atomtype[0]), "?") == 0)
538             {
539                 GMX_THROW(InconsistentInputError("Your input tpr file is too old (does not contain atom types). Cannot not compute Delta G of solvation"));
540             }
541             else
542             {
543                 printf("Free energy of solvation predictions:\n");
544                 please_cite(stdout, "Eisenberg86a");
545             }
546         }
547     }
548
549     // Now compute atomic radii including solvent probe size.
550     // Also, fetch solvation free energy coefficients and
551     // compute the residue indices that map the calculation atoms
552     // to the columns of residueArea_.
553     radii_.reserve(surfaceSel_.posCount());
554     if (bDGsol)
555     {
556         dgsFactor_.reserve(surfaceSel_.posCount());
557     }
558
559     const int resCount = surfaceSel_.initOriginalIdsToGroup(top_, INDEX_RES);
560
561     // TODO: Not exception-safe, but nice solution would be to have a C++
562     // atom properties class...
563     gmx_atomprop_t     aps = gmx_atomprop_init();
564
565     ConstArrayRef<int> atomIndices = surfaceSel_.atomIndices();
566     int                ndefault    = 0;
567     for (int i = 0; i < surfaceSel_.posCount(); i++)
568     {
569         const int ii     = atomIndices[i];
570         const int resind = atoms.atom[ii].resind;
571         real      radius;
572         if (!gmx_atomprop_query(aps, epropVDW,
573                                 *(atoms.resinfo[resind].name),
574                                 *(atoms.atomname[ii]), &radius))
575         {
576             ndefault++;
577         }
578         radii_.push_back(radius + solsize_);
579         if (bDGsol)
580         {
581             real dgsFactor;
582             if (!gmx_atomprop_query(aps, epropDGsol,
583                                     *(atoms.resinfo[resind].name),
584                                     *(atoms.atomtype[ii]), &dgsFactor))
585             {
586                 dgsFactor = dgsDefault_;
587             }
588             dgsFactor_.push_back(dgsFactor);
589         }
590     }
591     if (ndefault > 0)
592     {
593         fprintf(stderr, "WARNING: could not find a Van der Waals radius for %d atoms\n", ndefault);
594     }
595     gmx_atomprop_destroy(aps);
596
597     // Pre-compute mapping from the output groups to the calculation group,
598     // and store it in the selection ID map for easy lookup.
599     for (size_t g = 0; g < outputSel_.size(); ++g)
600     {
601         ConstArrayRef<int> outputIndices = outputSel_[g].atomIndices();
602         for (int i = 0, j = 0; i < outputSel_[g].posCount(); ++i)
603         {
604             while (j < surfaceSel_.posCount() && outputIndices[i] > atomIndices[j])
605             {
606                 ++j;
607             }
608             if (j == surfaceSel_.posCount() || outputIndices[i] != atomIndices[j])
609             {
610                 GMX_THROW(InconsistentInputError("Output selection is not a subset of the input selection"));
611             }
612             outputSel_[g].setOriginalId(i, j);
613         }
614     }
615
616     calculator_.setDotCount(ndots_);
617     calculator_.setRadii(radii_);
618
619     // Initialize all the output data objects and initialize the output plotters.
620
621     area_.setColumnCount(0, 1 + outputSel_.size());
622     {
623         AnalysisDataPlotModulePointer plotm(
624                 new AnalysisDataPlotModule(settings.plotSettings()));
625         plotm->setFileName(fnArea_);
626         plotm->setTitle("Solvent Accessible Surface");
627         plotm->setXAxisIsTime();
628         plotm->setYLabel("Area (nm\\S2\\N)");
629         plotm->appendLegend("Total");
630         for (size_t i = 0; i < outputSel_.size(); ++i)
631         {
632             plotm->appendLegend(outputSel_[i].name());
633         }
634         area_.addModule(plotm);
635     }
636
637     if (bResAt)
638     {
639         atomArea_.setDataSetCount(1 + outputSel_.size());
640         residueArea_.setDataSetCount(1 + outputSel_.size());
641         for (size_t i = 0; i <= outputSel_.size(); ++i)
642         {
643             atomArea_.setColumnCount(i, surfaceSel_.posCount());
644             residueArea_.setColumnCount(i, resCount);
645         }
646         {
647             AnalysisDataAverageModulePointer avem(new AnalysisDataAverageModule);
648             for (int i = 0; i < surfaceSel_.posCount(); ++i)
649             {
650                 avem->setXAxisValue(i, surfaceSel_.position(i).atomIndices()[0] + 1);
651             }
652             atomArea_.addModule(avem);
653             if (!fnAtomArea_.empty())
654             {
655                 AnalysisDataPlotModulePointer plotm(
656                         new AnalysisDataPlotModule(settings.plotSettings()));
657                 plotm->setFileName(fnAtomArea_);
658                 plotm->setTitle("Area per atom over the trajectory");
659                 plotm->setXLabel("Atom");
660                 plotm->setXFormat(8, 0);
661                 plotm->setYLabel("Area (nm\\S2\\N)");
662                 plotm->setErrorsAsSeparateColumn(true);
663                 plotm->appendLegend("Average (nm\\S2\\N)");
664                 plotm->appendLegend("Standard deviation (nm\\S2\\N)");
665                 avem->addModule(plotm);
666             }
667         }
668         {
669             AnalysisDataAverageModulePointer avem(new AnalysisDataAverageModule);
670             int prevResind = -1;
671             int row        = 0;
672             for (int i = 0; i < surfaceSel_.posCount(); ++i)
673             {
674                 const int atomIndex     = surfaceSel_.position(i).atomIndices()[0];
675                 const int residueIndex  = atoms.atom[atomIndex].resind;
676                 if (residueIndex != prevResind)
677                 {
678                     avem->setXAxisValue(row, atoms.resinfo[residueIndex].nr);
679                     prevResind = residueIndex;
680                     ++row;
681                 }
682             }
683             residueArea_.addModule(avem);
684             if (!fnResidueArea_.empty())
685             {
686                 AnalysisDataPlotModulePointer plotm(
687                         new AnalysisDataPlotModule(settings.plotSettings()));
688                 plotm->setFileName(fnResidueArea_);
689                 plotm->setTitle("Area per residue over the trajectory");
690                 plotm->setXLabel("Residue");
691                 plotm->setXFormat(8, 0);
692                 plotm->setYLabel("Area (nm\\S2\\N)");
693                 plotm->setErrorsAsSeparateColumn(true);
694                 plotm->appendLegend("Average (nm\\S2\\N)");
695                 plotm->appendLegend("Standard deviation (nm\\S2\\N)");
696                 avem->addModule(plotm);
697             }
698         }
699     }
700
701     if (!fnDGSolv_.empty())
702     {
703         dgSolv_.setColumnCount(0, 1 + outputSel_.size());
704         AnalysisDataPlotModulePointer plotm(
705                 new AnalysisDataPlotModule(settings.plotSettings()));
706         plotm->setFileName(fnDGSolv_);
707         plotm->setTitle("Free Energy of Solvation");
708         plotm->setXAxisIsTime();
709         plotm->setYLabel("D Gsolv");
710         plotm->appendLegend("Total");
711         for (size_t i = 0; i < outputSel_.size(); ++i)
712         {
713             plotm->appendLegend(outputSel_[i].name());
714         }
715         dgSolv_.addModule(plotm);
716     }
717
718     if (!fnVolume_.empty())
719     {
720         volume_.setColumnCount(0, 2);
721         AnalysisDataPlotModulePointer plotm(
722                 new AnalysisDataPlotModule(settings.plotSettings()));
723         plotm->setFileName(fnVolume_);
724         plotm->setTitle("Volume and Density");
725         plotm->setXAxisIsTime();
726         plotm->appendLegend("Volume (nm\\S3\\N)");
727         plotm->appendLegend("Density (g/l)");
728         volume_.addModule(plotm);
729     }
730 }
731
732 /*! \brief
733  * Temporary memory for use within a single-frame calculation.
734  */
735 class SasaModuleData : public TrajectoryAnalysisModuleData
736 {
737     public:
738         /*! \brief
739          * Reserves memory for the frame-local data.
740          *
741          * `residueCount` will be zero if per-residue data is not being
742          * calculated.
743          */
744         SasaModuleData(TrajectoryAnalysisModule          *module,
745                        const AnalysisDataParallelOptions &opt,
746                        const SelectionCollection         &selections,
747                        int atomCount, int residueCount)
748             : TrajectoryAnalysisModuleData(module, opt, selections)
749         {
750             index_.reserve(atomCount);
751             // If the calculation group is not dynamic, pre-calculate
752             // the index, since it is not going to change.
753             for (int i = 0; i < atomCount; ++i)
754             {
755                 index_.push_back(i);
756             }
757             atomAreas_.resize(atomCount);
758             res_a_.resize(residueCount);
759         }
760
761         virtual void finish() { finishDataHandles(); }
762
763         //! Indices of the calculation selection positions selected for the frame.
764         std::vector<int>        index_;
765         /*! \brief
766          * Atom areas for each calculation selection position for the frame.
767          *
768          * One entry for each position in the calculation group.
769          * Values for atoms not selected are set to zero.
770          */
771         std::vector<real>       atomAreas_;
772         /*! \brief
773          * Working array to accumulate areas for each residue.
774          *
775          * One entry for each distinct residue in the calculation group;
776          * indices are not directly residue numbers or residue indices.
777          *
778          * This vector is empty if residue area calculations are not being
779          * performed.
780          */
781         std::vector<real>       res_a_;
782 };
783
784 TrajectoryAnalysisModuleDataPointer Sasa::startFrames(
785         const AnalysisDataParallelOptions &opt,
786         const SelectionCollection         &selections)
787 {
788     return TrajectoryAnalysisModuleDataPointer(
789             new SasaModuleData(this, opt, selections, surfaceSel_.posCount(),
790                                residueArea_.columnCount(0)));
791 }
792
793 /*! \brief
794  * Helper method to compute the areas for a single selection.
795  *
796  * \param[in]  surfaceSel     The calculation selection.
797  * \param[in]  sel            The selection to compute the areas for (can be
798  *     `surfaceSel` or one of the output selections).
799  * \param[in]  atomAreas      Atom areas for each position in `surfaceSel`.
800  * \param[in]  dgsFactor      Free energy coefficients for each position in
801  *     `surfaceSel`. If empty, free energies are not calculated.
802  * \param[out] totalAreaOut Total area of `sel` (sum of atom areas it selects).
803  * \param[out] dgsolvOut      Solvation free energy.
804  *     Will be zero of `dgsFactor` is empty.
805  * \param      atomAreaHandle Data handle to use for storing atom areas for `sel`.
806  * \param      resAreaHandle  Data handle to use for storing residue areas for `sel`.
807  * \param      resAreaWork    Work array for accumulating the residue areas.
808  *     If empty, atom and residue areas are not calculated.
809  *
810  * `atomAreaHandle` and `resAreaHandle` are not used if `resAreaWork` is empty.
811  */
812 void computeAreas(const Selection &surfaceSel, const Selection &sel,
813                   const std::vector<real> &atomAreas,
814                   const std::vector<real> &dgsFactor,
815                   real *totalAreaOut, real *dgsolvOut,
816                   AnalysisDataHandle atomAreaHandle,
817                   AnalysisDataHandle resAreaHandle,
818                   std::vector<real> *resAreaWork)
819 {
820     const bool bResAt    = !resAreaWork->empty();
821     const bool bDGsolv   = !dgsFactor.empty();
822     real       totalArea = 0;
823     real       dgsolv    = 0;
824
825     if (bResAt)
826     {
827         std::fill(resAreaWork->begin(), resAreaWork->end(),
828                   static_cast<real>(0.0));
829     }
830     for (int i = 0; i < sel.posCount(); ++i)
831     {
832         // Get the index of the atom in the calculation group.
833         // For the output groups, the mapping has been precalculated in
834         // initAnalysis().
835         const int  ii = (sel != surfaceSel ? sel.position(i).mappedId() : i);
836         if (!surfaceSel.position(ii).selected())
837         {
838             // For the calculation group, skip unselected atoms.
839             if (sel == surfaceSel)
840             {
841                 continue;
842             }
843             GMX_THROW(InconsistentInputError("Output selection is not a subset of the surface selection"));
844         }
845         // Get the internal index of the matching residue.
846         // These have been precalculated in initAnalysis().
847         const int  ri       = surfaceSel.position(ii).mappedId();
848         const real atomArea = atomAreas[ii];
849         totalArea += atomArea;
850         if (bResAt)
851         {
852             atomAreaHandle.setPoint(ii, atomArea);
853             (*resAreaWork)[ri] += atomArea;
854         }
855         if (bDGsolv)
856         {
857             dgsolv += atomArea * dgsFactor[ii];
858         }
859     }
860     if (bResAt)
861     {
862         for (size_t i = 0; i < (*resAreaWork).size(); ++i)
863         {
864             resAreaHandle.setPoint(i, (*resAreaWork)[i]);
865         }
866     }
867     *totalAreaOut = totalArea;
868     *dgsolvOut    = dgsolv;
869 }
870
871 void
872 Sasa::analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc,
873                    TrajectoryAnalysisModuleData *pdata)
874 {
875     AnalysisDataHandle   ah         = pdata->dataHandle(area_);
876     AnalysisDataHandle   dgh        = pdata->dataHandle(dgSolv_);
877     AnalysisDataHandle   aah        = pdata->dataHandle(atomArea_);
878     AnalysisDataHandle   rah        = pdata->dataHandle(residueArea_);
879     AnalysisDataHandle   vh         = pdata->dataHandle(volume_);
880     const Selection     &surfaceSel = pdata->parallelSelection(surfaceSel_);
881     const SelectionList &outputSel  = pdata->parallelSelections(outputSel_);
882     SasaModuleData      &frameData  = *static_cast<SasaModuleData *>(pdata);
883
884     const bool           bResAt    = !frameData.res_a_.empty();
885     const bool           bDGsol    = !dgsFactor_.empty();
886     const bool           bConnolly = (frnr == 0 && !fnConnolly_.empty());
887
888     // Update indices of selected atoms in the work array.
889     if (surfaceSel.isDynamic())
890     {
891         frameData.index_.clear();
892         for (int i = 0; i < surfaceSel.posCount(); ++i)
893         {
894             if (surfaceSel.position(i).selected())
895             {
896                 frameData.index_.push_back(i);
897             }
898         }
899     }
900
901     // Determine what needs to be calculated.
902     int                  flag      = 0;
903     if (bResAt || bDGsol || !outputSel.empty())
904     {
905         flag |= FLAG_ATOM_AREA;
906     }
907     if (bConnolly)
908     {
909         flag |= FLAG_DOTS;
910     }
911     if (volume_.columnCount() > 0)
912     {
913         flag |= FLAG_VOLUME;
914     }
915
916     // Do the low-level calculation.
917     // totarea and totvolume receive the values for the calculation group.
918     // area array contains the per-atom areas for the selected positions.
919     // surfacedots contains nsurfacedots entries, and contains the actual
920     // points.
921     real  totarea, totvolume;
922     real *area = NULL, *surfacedots = NULL;
923     int   nsurfacedots;
924     calculator_.calculate(surfaceSel.coordinates().data(), pbc,
925                           frameData.index_.size(), &frameData.index_[0], flag,
926                           &totarea, &totvolume, &area,
927                           &surfacedots, &nsurfacedots);
928     // Unpack the atomwise areas into the frameData.atomAreas_ array for easier
929     // indexing in the case of dynamic surfaceSel.
930     if (area != NULL)
931     {
932         if (surfaceSel.isDynamic())
933         {
934             std::fill(frameData.atomAreas_.begin(), frameData.atomAreas_.end(),
935                       static_cast<real>(0.0));
936             for (size_t i = 0; i < frameData.index_.size(); ++i)
937             {
938                 frameData.atomAreas_[frameData.index_[i]] = area[i];
939             }
940         }
941         else
942         {
943             std::copy(area, area + surfaceSel.posCount(),
944                       frameData.atomAreas_.begin());
945         }
946         sfree(area);
947     }
948     scoped_guard_sfree dotsGuard(surfacedots);
949
950     if (bConnolly)
951     {
952         // This is somewhat nasty, as it modifies the atoms and symtab
953         // structures.  But since it is only used in the first frame, and no
954         // one else uses the topology after initialization, it may just work
955         // even with future parallelization.
956         connolly_plot(fnConnolly_.c_str(),
957                       nsurfacedots, surfacedots, fr.x, &top_->atoms,
958                       &top_->symtab, fr.ePBC, fr.box, bIncludeSolute_);
959     }
960
961     ah.startFrame(frnr, fr.time);
962     if (bResAt)
963     {
964         aah.startFrame(frnr, fr.time);
965         rah.startFrame(frnr, fr.time);
966     }
967     if (bDGsol)
968     {
969         dgh.startFrame(frnr, fr.time);
970     }
971
972     ah.setPoint(0, totarea);
973
974     real totalArea, dgsolv;
975     if (bResAt || bDGsol)
976     {
977         computeAreas(surfaceSel, surfaceSel, frameData.atomAreas_, dgsFactor_,
978                      &totalArea, &dgsolv, aah, rah, &frameData.res_a_);
979         if (bDGsol)
980         {
981             dgh.setPoint(0, dgsolv);
982         }
983     }
984     for (size_t g = 0; g < outputSel.size(); ++g)
985     {
986         if (bResAt)
987         {
988             aah.selectDataSet(g + 1);
989             rah.selectDataSet(g + 1);
990         }
991         computeAreas(surfaceSel, outputSel[g], frameData.atomAreas_, dgsFactor_,
992                      &totalArea, &dgsolv, aah, rah, &frameData.res_a_);
993         ah.setPoint(g + 1, totalArea);
994         if (bDGsol)
995         {
996             dgh.setPoint(g + 1, dgsolv);
997         }
998     }
999
1000     ah.finishFrame();
1001     if (bResAt)
1002     {
1003         aah.finishFrame();
1004         rah.finishFrame();
1005     }
1006     if (bDGsol)
1007     {
1008         dgh.finishFrame();
1009     }
1010
1011     if (vh.isValid())
1012     {
1013         real totmass = 0;
1014         for (int i = 0; i < surfaceSel.posCount(); ++i)
1015         {
1016             totmass += surfaceSel.position(i).mass();
1017         }
1018         const real density = totmass*AMU/(totvolume*NANO*NANO*NANO);
1019         vh.startFrame(frnr, fr.time);
1020         vh.setPoint(0, totvolume);
1021         vh.setPoint(1, density);
1022         vh.finishFrame();
1023     }
1024 }
1025
1026 void
1027 Sasa::finishAnalysis(int /*nframes*/)
1028 {
1029     //if (bITP)
1030     //{
1031     //    fp3 = ftp2FILE(efITP, nfile, fnm, "w");
1032     //    fprintf(fp3, "[ position_restraints ]\n"
1033     //            "#define FCX 1000\n"
1034     //            "#define FCY 1000\n"
1035     //            "#define FCZ 1000\n"
1036     //            "; Atom  Type  fx   fy   fz\n");
1037     //    for (i = 0; i < nx[0]; i++)
1038     //    {
1039     //        if (atom_area[i] > minarea)
1040     //        {
1041     //            fprintf(fp3, "%5d   1     FCX  FCX  FCZ\n", ii+1);
1042     //        }
1043     //    }
1044     //    ffclose(fp3);
1045     //}
1046 }
1047
1048 void
1049 Sasa::writeOutput()
1050 {
1051 }
1052
1053 //! \}
1054
1055 }       // namespace
1056
1057 const char SasaInfo::name[]             = "sasa";
1058 const char SasaInfo::shortDescription[] =
1059     "Compute solvent accessible surface area";
1060
1061 TrajectoryAnalysisModulePointer SasaInfo::create()
1062 {
1063     return TrajectoryAnalysisModulePointer(new Sasa);
1064 }
1065
1066 } // namespace analysismodules
1067
1068 } // namespace gmx