Merge release-5-0 into master
[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, 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
391         // Copy and assign disallowed by base.
392 };
393
394 Sasa::Sasa()
395     : TrajectoryAnalysisModule(SasaInfo::name, SasaInfo::shortDescription),
396       solsize_(0.14), ndots_(24), dgsDefault_(0), bIncludeSolute_(true), top_(NULL)
397 {
398     //minarea_ = 0.5;
399     registerAnalysisDataset(&area_, "area");
400     registerAnalysisDataset(&atomArea_, "atomarea");
401     registerAnalysisDataset(&residueArea_, "resarea");
402     registerAnalysisDataset(&dgSolv_, "dgsolv");
403     registerAnalysisDataset(&volume_, "volume");
404 }
405
406 void
407 Sasa::initOptions(Options *options, TrajectoryAnalysisSettings *settings)
408 {
409     static const char *const desc[] = {
410         "[THISMODULE] computes solvent accessible surface areas.",
411         "See Eisenhaber F, Lijnzaad P, Argos P, Sander C, & Scharf M",
412         "(1995) J. Comput. Chem. 16, 273-284 for the algorithm used.",
413         "With [TT]-q[tt], the Connolly surface can be generated as well",
414         "in a [TT].pdb[tt] file where the nodes are represented as atoms",
415         "and the edges connecting the nearest nodes as CONECT records.",
416         "[TT]-odg[tt] allows for estimation of solvation free energies",
417         "from per-atom solvation energies per exposed surface area.[PAR]",
418
419         "The program requires a selection for the surface calculation to be",
420         "specified with [TT]-surface[tt]. This should always consist of all",
421         "non-solvent atoms in the system. The area of this group is always",
422         "calculated. Optionally, [TT]-output[tt] can specify additional",
423         "selections, which should be subsets of the calculation group.",
424         "The solvent-accessible areas for these groups are also extracted",
425         "from the full surface.[PAR]",
426
427         "The average and standard deviation of the area over the trajectory",
428         "can be calculated per residue and atom (options [TT]-or[tt] and",
429         "[TT]-oa[tt]).[PAR]",
430         //"In combination with the latter option an [TT].itp[tt] file can be",
431         //"generated (option [TT]-i[tt])",
432         //"which can be used to restrain surface atoms.[PAR]",
433
434         "With the [TT]-tv[tt] option the total volume and density of the",
435         "molecule can be computed.",
436         "Please consider whether the normal probe radius is appropriate",
437         "in this case or whether you would rather use, e.g., 0. It is good",
438         "to keep in mind that the results for volume and density are very",
439         "approximate. For example, in ice Ih, one can easily fit water molecules in the",
440         "pores which would yield a volume that is too low, and surface area and density",
441         "that are both too high."
442     };
443
444     options->setDescription(concatenateStrings(desc));
445
446     options->addOption(FileNameOption("o").filetype(eftPlot).outputFile().required()
447                            .store(&fnArea_).defaultBasename("area")
448                            .description("Total area as a function of time"));
449     options->addOption(FileNameOption("odg").filetype(eftPlot).outputFile()
450                            .store(&fnDGSolv_).defaultBasename("dgsolv")
451                            .description("Estimated solvation free energy as a function of time"));
452     options->addOption(FileNameOption("or").filetype(eftPlot).outputFile()
453                            .store(&fnResidueArea_).defaultBasename("resarea")
454                            .description("Average area per residue"));
455     options->addOption(FileNameOption("oa").filetype(eftPlot).outputFile()
456                            .store(&fnAtomArea_).defaultBasename("atomarea")
457                            .description("Average area per atom"));
458     options->addOption(FileNameOption("tv").filetype(eftPlot).outputFile()
459                            .store(&fnVolume_).defaultBasename("volume")
460                            .description("Total volume and density as a function of time"));
461     options->addOption(FileNameOption("q").filetype(eftPDB).outputFile()
462                            .store(&fnConnolly_).defaultBasename("connolly")
463                            .description("PDB file for Connolly surface"));
464     //options->addOption(FileNameOption("i").filetype(eftITP).outputFile()
465     //                       .store(&fnRestraints_).defaultBasename("surfat")
466     //                       .description("Topology file for position restraings on surface atoms"));
467
468
469     options->addOption(DoubleOption("probe").store(&solsize_)
470                            .description("Radius of the solvent probe (nm)"));
471     options->addOption(IntegerOption("ndots").store(&ndots_)
472                            .description("Number of dots per sphere, more dots means more accuracy"));
473     //options->addOption(DoubleOption("minarea").store(&minarea_)
474     //                       .description("The minimum area (nm^2) to count an atom as a surface atom when writing a position restraint file (see help)"));
475     options->addOption(BooleanOption("prot").store(&bIncludeSolute_)
476                            .description("Output the protein to the Connolly [TT].pdb[tt] file too"));
477     options->addOption(DoubleOption("dgs").store(&dgsDefault_)
478                            .description("Default value for solvation free energy per area (kJ/mol/nm^2)"));
479
480     // Selections must select atoms for the VdW radii lookup to work.
481     // The calculation group uses dynamicMask() so that the coordinates
482     // match a static array of VdW radii.
483     options->addOption(SelectionOption("surface").store(&surfaceSel_)
484                            .required().onlyAtoms().dynamicMask()
485                            .description("Surface calculation selection"));
486     options->addOption(SelectionOption("output").storeVector(&outputSel_)
487                            .onlyAtoms().multiValue()
488                            .description("Output selection(s)"));
489
490     // Atom names etc. are required for the VdW radii lookup.
491     settings->setFlag(TrajectoryAnalysisSettings::efRequireTop);
492 }
493
494 void
495 Sasa::initAnalysis(const TrajectoryAnalysisSettings &settings,
496                    const TopologyInformation        &top)
497 {
498     const t_atoms &atoms = top.topology()->atoms;
499     top_ = top.topology();
500
501     //bITP   = opt2bSet("-i", nfile, fnm);
502     const bool bResAt =
503         !fnResidueArea_.empty() || !fnAtomArea_.empty(); // || bITP;
504     const bool bDGsol = !fnDGSolv_.empty();
505
506     if (solsize_ < 0)
507     {
508         solsize_ = 1e-3;
509         fprintf(stderr, "Probe size too small, setting it to %g\n", solsize_);
510     }
511     if (ndots_ < 20)
512     {
513         ndots_ = 20;
514         fprintf(stderr, "Ndots too small, setting it to %d\n", ndots_);
515     }
516
517     please_cite(stderr, "Eisenhaber95");
518     //if ((top.ePBC() != epbcXYZ) || (TRICLINIC(fr.box)))
519     //{
520     //    fprintf(stderr, "\n\nWARNING: non-rectangular boxes may give erroneous results or crashes.\n"
521     //            "Analysis based on vacuum simulations (with the possibility of evaporation)\n"
522     //            "will certainly crash the analysis.\n\n");
523     //}
524
525     if (bDGsol)
526     {
527         if (!top.hasFullTopology())
528         {
529             GMX_THROW(InconsistentInputError("Cannot compute Delta G of solvation without a tpr file"));
530         }
531         else
532         {
533             if (strcmp(*(atoms.atomtype[0]), "?") == 0)
534             {
535                 GMX_THROW(InconsistentInputError("Your input tpr file is too old (does not contain atom types). Cannot not compute Delta G of solvation"));
536             }
537             else
538             {
539                 printf("Free energy of solvation predictions:\n");
540                 please_cite(stdout, "Eisenberg86a");
541             }
542         }
543     }
544
545     // Now compute atomic radii including solvent probe size.
546     // Also, fetch solvation free energy coefficients and
547     // compute the residue indices that map the calculation atoms
548     // to the columns of residueArea_.
549     radii_.reserve(surfaceSel_.posCount());
550     if (bDGsol)
551     {
552         dgsFactor_.reserve(surfaceSel_.posCount());
553     }
554
555     const int resCount = surfaceSel_.initOriginalIdsToGroup(top_, INDEX_RES);
556
557     // TODO: Not exception-safe, but nice solution would be to have a C++
558     // atom properties class...
559     gmx_atomprop_t     aps = gmx_atomprop_init();
560
561     ConstArrayRef<int> atomIndices = surfaceSel_.atomIndices();
562     int                ndefault    = 0;
563     for (int i = 0; i < surfaceSel_.posCount(); i++)
564     {
565         const int ii     = atomIndices[i];
566         const int resind = atoms.atom[ii].resind;
567         real      radius;
568         if (!gmx_atomprop_query(aps, epropVDW,
569                                 *(atoms.resinfo[resind].name),
570                                 *(atoms.atomname[ii]), &radius))
571         {
572             ndefault++;
573         }
574         radii_.push_back(radius + solsize_);
575         if (bDGsol)
576         {
577             real dgsFactor;
578             if (!gmx_atomprop_query(aps, epropDGsol,
579                                     *(atoms.resinfo[resind].name),
580                                     *(atoms.atomtype[ii]), &dgsFactor))
581             {
582                 dgsFactor = dgsDefault_;
583             }
584             dgsFactor_.push_back(dgsFactor);
585         }
586     }
587     if (ndefault > 0)
588     {
589         fprintf(stderr, "WARNING: could not find a Van der Waals radius for %d atoms\n", ndefault);
590     }
591     gmx_atomprop_destroy(aps);
592
593     // Pre-compute mapping from the output groups to the calculation group,
594     // and store it in the selection ID map for easy lookup.
595     for (size_t g = 0; g < outputSel_.size(); ++g)
596     {
597         ConstArrayRef<int> outputIndices = outputSel_[g].atomIndices();
598         for (int i = 0, j = 0; i < outputSel_[g].posCount(); ++i)
599         {
600             while (j < surfaceSel_.posCount() && outputIndices[i] > atomIndices[j])
601             {
602                 ++j;
603             }
604             if (j == surfaceSel_.posCount() || outputIndices[i] != atomIndices[j])
605             {
606                 GMX_THROW(InconsistentInputError("Output selection is not a subset of the input selection"));
607             }
608             outputSel_[g].setOriginalId(i, j);
609         }
610     }
611
612     // Initialize all the output data objects and initialize the output plotters.
613
614     area_.setColumnCount(0, 1 + outputSel_.size());
615     {
616         AnalysisDataPlotModulePointer plotm(
617                 new AnalysisDataPlotModule(settings.plotSettings()));
618         plotm->setFileName(fnArea_);
619         plotm->setTitle("Solvent Accessible Surface");
620         plotm->setXAxisIsTime();
621         plotm->setYLabel("Area (nm\\S2\\N)");
622         plotm->appendLegend("Total");
623         for (size_t i = 0; i < outputSel_.size(); ++i)
624         {
625             plotm->appendLegend(outputSel_[i].name());
626         }
627         area_.addModule(plotm);
628     }
629
630     if (bResAt)
631     {
632         atomArea_.setDataSetCount(1 + outputSel_.size());
633         residueArea_.setDataSetCount(1 + outputSel_.size());
634         for (size_t i = 0; i <= outputSel_.size(); ++i)
635         {
636             atomArea_.setColumnCount(i, surfaceSel_.posCount());
637             residueArea_.setColumnCount(i, resCount);
638         }
639         {
640             AnalysisDataAverageModulePointer avem(new AnalysisDataAverageModule);
641             for (int i = 0; i < surfaceSel_.posCount(); ++i)
642             {
643                 avem->setXAxisValue(i, surfaceSel_.position(i).atomIndices()[0] + 1);
644             }
645             atomArea_.addModule(avem);
646             if (!fnAtomArea_.empty())
647             {
648                 AnalysisDataPlotModulePointer plotm(
649                         new AnalysisDataPlotModule(settings.plotSettings()));
650                 plotm->setFileName(fnAtomArea_);
651                 plotm->setTitle("Area per residue over the trajectory");
652                 plotm->setXLabel("Atom");
653                 plotm->setXFormat(8, 0);
654                 plotm->setYLabel("Area (nm\\S2\\N)");
655                 plotm->setErrorsAsSeparateColumn(true);
656                 plotm->appendLegend("Average (nm\\S2\\N)");
657                 plotm->appendLegend("Standard deviation (nm\\S2\\N)");
658                 avem->addModule(plotm);
659             }
660         }
661         {
662             AnalysisDataAverageModulePointer avem(new AnalysisDataAverageModule);
663             for (int i = 0; i < surfaceSel_.posCount(); ++i)
664             {
665                 const int atomIndex     = surfaceSel_.position(i).atomIndices()[0];
666                 const int residueIndex  = atoms.atom[atomIndex].resind;
667                 avem->setXAxisValue(i, atoms.resinfo[residueIndex].nr);
668             }
669             residueArea_.addModule(avem);
670             if (!fnResidueArea_.empty())
671             {
672                 AnalysisDataPlotModulePointer plotm(
673                         new AnalysisDataPlotModule(settings.plotSettings()));
674                 plotm->setFileName(fnResidueArea_);
675                 plotm->setTitle("Area per atom over the trajectory");
676                 plotm->setXLabel("Residue");
677                 plotm->setXFormat(8, 0);
678                 plotm->setYLabel("Area (nm\\S2\\N)");
679                 plotm->setErrorsAsSeparateColumn(true);
680                 plotm->appendLegend("Average (nm\\S2\\N)");
681                 plotm->appendLegend("Standard deviation (nm\\S2\\N)");
682                 avem->addModule(plotm);
683             }
684         }
685     }
686
687     if (!fnDGSolv_.empty())
688     {
689         dgSolv_.setColumnCount(0, 1 + outputSel_.size());
690         AnalysisDataPlotModulePointer plotm(
691                 new AnalysisDataPlotModule(settings.plotSettings()));
692         plotm->setFileName(fnDGSolv_);
693         plotm->setTitle("Free Energy of Solvation");
694         plotm->setXAxisIsTime();
695         plotm->setYLabel("D Gsolv");
696         plotm->appendLegend("Total");
697         for (size_t i = 0; i < outputSel_.size(); ++i)
698         {
699             plotm->appendLegend(outputSel_[i].name());
700         }
701         dgSolv_.addModule(plotm);
702     }
703
704     if (!fnVolume_.empty())
705     {
706         volume_.setColumnCount(0, 2);
707         AnalysisDataPlotModulePointer plotm(
708                 new AnalysisDataPlotModule(settings.plotSettings()));
709         plotm->setFileName(fnVolume_);
710         plotm->setTitle("Volume and Density");
711         plotm->setXAxisIsTime();
712         plotm->appendLegend("Volume (nm\\S3\\N)");
713         plotm->appendLegend("Density (g/l)");
714         volume_.addModule(plotm);
715     }
716 }
717
718 /*! \brief
719  * Temporary memory for use within a single-frame calculation.
720  */
721 class SasaModuleData : public TrajectoryAnalysisModuleData
722 {
723     public:
724         /*! \brief
725          * Reserves memory for the frame-local data.
726          *
727          * `residueCount` will be zero if per-residue data is not being
728          * calculated.
729          */
730         SasaModuleData(TrajectoryAnalysisModule          *module,
731                        const AnalysisDataParallelOptions &opt,
732                        const SelectionCollection         &selections,
733                        int atomCount, int residueCount)
734             : TrajectoryAnalysisModuleData(module, opt, selections)
735         {
736             index_.reserve(atomCount);
737             // If the calculation group is not dynamic, pre-calculate
738             // the index, since it is not going to change.
739             for (int i = 0; i < atomCount; ++i)
740             {
741                 index_.push_back(i);
742             }
743             atomAreas_.resize(atomCount);
744             res_a_.resize(residueCount);
745         }
746
747         virtual void finish() { finishDataHandles(); }
748
749         //! Indices of the calculation selection positions selected for the frame.
750         std::vector<int>        index_;
751         /*! \brief
752          * Atom areas for each calculation selection position for the frame.
753          *
754          * One entry for each position in the calculation group.
755          * Values for atoms not selected are set to zero.
756          */
757         std::vector<real>       atomAreas_;
758         /*! \brief
759          * Working array to accumulate areas for each residue.
760          *
761          * One entry for each distinct residue in the calculation group;
762          * indices are not directly residue numbers or residue indices.
763          *
764          * This vector is empty if residue area calculations are not being
765          * performed.
766          */
767         std::vector<real>       res_a_;
768 };
769
770 TrajectoryAnalysisModuleDataPointer Sasa::startFrames(
771         const AnalysisDataParallelOptions &opt,
772         const SelectionCollection         &selections)
773 {
774     return TrajectoryAnalysisModuleDataPointer(
775             new SasaModuleData(this, opt, selections, surfaceSel_.posCount(),
776                                residueArea_.columnCount(0)));
777 }
778
779 /*! \brief
780  * Helper method to compute the areas for a single selection.
781  *
782  * \param[in]  surfaceSel     The calculation selection.
783  * \param[in]  sel            The selection to compute the areas for (can be
784  *     `surfaceSel` or one of the output selections).
785  * \param[in]  atomAreas      Atom areas for each position in `surfaceSel`.
786  * \param[in]  dgsFactor      Free energy coefficients for each position in
787  *     `surfaceSel`. If empty, free energies are not calculated.
788  * \param[out] totalAreaOut Total area of `sel` (sum of atom areas it selects).
789  * \param[out] dgsolvOut      Solvation free energy.
790  *     Will be zero of `dgsFactor` is empty.
791  * \param      atomAreaHandle Data handle to use for storing atom areas for `sel`.
792  * \param      resAreaHandle  Data handle to use for storing residue areas for `sel`.
793  * \param      resAreaWork    Work array for accumulating the residue areas.
794  *     If empty, atom and residue areas are not calculated.
795  *
796  * `atomAreaHandle` and `resAreaHandle` are not used if `resAreaWork` is empty.
797  */
798 void computeAreas(const Selection &surfaceSel, const Selection &sel,
799                   const std::vector<real> &atomAreas,
800                   const std::vector<real> &dgsFactor,
801                   real *totalAreaOut, real *dgsolvOut,
802                   AnalysisDataHandle atomAreaHandle,
803                   AnalysisDataHandle resAreaHandle,
804                   std::vector<real> *resAreaWork)
805 {
806     const bool bResAt    = !resAreaWork->empty();
807     const bool bDGsolv   = !dgsFactor.empty();
808     real       totalArea = 0;
809     real       dgsolv    = 0;
810
811     if (bResAt)
812     {
813         std::fill(resAreaWork->begin(), resAreaWork->end(),
814                   static_cast<real>(0.0));
815     }
816     for (int i = 0; i < sel.posCount(); ++i)
817     {
818         // Get the index of the atom in the calculation group.
819         // For the output groups, the mapping has been precalculated in
820         // initAnalysis().
821         const int  ii = (sel != surfaceSel ? sel.position(i).mappedId() : i);
822         if (!surfaceSel.position(ii).selected())
823         {
824             // For the calculation group, skip unselected atoms.
825             if (sel == surfaceSel)
826             {
827                 continue;
828             }
829             GMX_THROW(InconsistentInputError("Output selection is not a subset of the surface selection"));
830         }
831         // Get the internal index of the matching residue.
832         // These have been precalculated in initAnalysis().
833         const int  ri       = surfaceSel.position(ii).mappedId();
834         const real atomArea = atomAreas[ii];
835         totalArea += atomArea;
836         if (bResAt)
837         {
838             atomAreaHandle.setPoint(ii, atomArea);
839             (*resAreaWork)[ri] += atomArea;
840         }
841         if (bDGsolv)
842         {
843             dgsolv += atomArea * dgsFactor[ii];
844         }
845     }
846     if (bResAt)
847     {
848         for (size_t i = 0; i < (*resAreaWork).size(); ++i)
849         {
850             resAreaHandle.setPoint(i, (*resAreaWork)[i]);
851         }
852     }
853     *totalAreaOut = totalArea;
854     *dgsolvOut    = dgsolv;
855 }
856
857 void
858 Sasa::analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc,
859                    TrajectoryAnalysisModuleData *pdata)
860 {
861     AnalysisDataHandle   ah         = pdata->dataHandle(area_);
862     AnalysisDataHandle   dgh        = pdata->dataHandle(dgSolv_);
863     AnalysisDataHandle   aah        = pdata->dataHandle(atomArea_);
864     AnalysisDataHandle   rah        = pdata->dataHandle(residueArea_);
865     AnalysisDataHandle   vh         = pdata->dataHandle(volume_);
866     const Selection     &surfaceSel = pdata->parallelSelection(surfaceSel_);
867     const SelectionList &outputSel  = pdata->parallelSelections(outputSel_);
868     SasaModuleData      &frameData  = *static_cast<SasaModuleData *>(pdata);
869
870     const bool           bResAt    = !frameData.res_a_.empty();
871     const bool           bDGsol    = !dgsFactor_.empty();
872     const bool           bConnolly = (frnr == 0 && !fnConnolly_.empty());
873
874     // Update indices of selected atoms in the work array.
875     if (surfaceSel.isDynamic())
876     {
877         frameData.index_.clear();
878         for (int i = 0; i < surfaceSel.posCount(); ++i)
879         {
880             if (surfaceSel.position(i).selected())
881             {
882                 frameData.index_.push_back(i);
883             }
884         }
885     }
886
887     // Determine what needs to be calculated.
888     int                  flag      = 0;
889     if (bResAt || bDGsol || !outputSel.empty())
890     {
891         flag |= FLAG_ATOM_AREA;
892     }
893     if (bConnolly)
894     {
895         flag |= FLAG_DOTS;
896     }
897     if (volume_.columnCount() > 0)
898     {
899         flag |= FLAG_VOLUME;
900     }
901
902     // Do the low-level calculation.
903     // totarea and totvolume receive the values for the calculation group.
904     // area array contains the per-atom areas for the selected positions.
905     // surfacedots contains nsurfacedots entries, and contains the actual
906     // points.
907     real  totarea, totvolume;
908     real *area = NULL, *surfacedots = NULL;
909     int   nsurfacedots;
910     int   retval = nsc_dclm_pbc(surfaceSel.coordinates().data(), &radii_[0],
911                                 frameData.index_.size(), ndots_, flag, &totarea,
912                                 &area, &totvolume, &surfacedots, &nsurfacedots,
913                                 &frameData.index_[0],
914                                 pbc != NULL ? pbc->ePBC : epbcNONE,
915                                 pbc != NULL ? pbc->box : NULL);
916     // Unpack the atomwise areas into the frameData.atomAreas_ array for easier
917     // indexing in the case of dynamic surfaceSel.
918     if (area != NULL)
919     {
920         if (surfaceSel.isDynamic())
921         {
922             std::fill(frameData.atomAreas_.begin(), frameData.atomAreas_.end(),
923                       static_cast<real>(0.0));
924             for (size_t i = 0; i < frameData.index_.size(); ++i)
925             {
926                 frameData.atomAreas_[frameData.index_[i]] = area[i];
927             }
928         }
929         else
930         {
931             std::copy(area, area + surfaceSel.posCount(),
932                       frameData.atomAreas_.begin());
933         }
934         sfree(area);
935     }
936     scoped_guard_sfree dotsGuard(surfacedots);
937     if (retval != 0)
938     {
939         GMX_THROW(InternalError("nsc_dclm_pbc failed"));
940     }
941
942     if (bConnolly)
943     {
944         // This is somewhat nasty, as it modifies the atoms and symtab
945         // structures.  But since it is only used in the first frame, and no
946         // one else uses the topology after initialization, it may just work
947         // even with future parallelization.
948         connolly_plot(fnConnolly_.c_str(),
949                       nsurfacedots, surfacedots, fr.x, &top_->atoms,
950                       &top_->symtab, fr.ePBC, fr.box, bIncludeSolute_);
951     }
952
953     ah.startFrame(frnr, fr.time);
954     if (bResAt)
955     {
956         aah.startFrame(frnr, fr.time);
957         rah.startFrame(frnr, fr.time);
958     }
959     if (bDGsol)
960     {
961         dgh.startFrame(frnr, fr.time);
962     }
963
964     ah.setPoint(0, totarea);
965
966     real totalArea, dgsolv;
967     if (bResAt || bDGsol)
968     {
969         computeAreas(surfaceSel, surfaceSel, frameData.atomAreas_, dgsFactor_,
970                      &totalArea, &dgsolv, aah, rah, &frameData.res_a_);
971         if (bDGsol)
972         {
973             dgh.setPoint(0, dgsolv);
974         }
975     }
976     for (size_t g = 0; g < outputSel.size(); ++g)
977     {
978         if (bResAt)
979         {
980             aah.selectDataSet(g + 1);
981             rah.selectDataSet(g + 1);
982         }
983         computeAreas(surfaceSel, outputSel[g], frameData.atomAreas_, dgsFactor_,
984                      &totalArea, &dgsolv, aah, rah, &frameData.res_a_);
985         ah.setPoint(g + 1, totalArea);
986         if (bDGsol)
987         {
988             dgh.setPoint(g + 1, dgsolv);
989         }
990     }
991
992     ah.finishFrame();
993     if (bResAt)
994     {
995         aah.finishFrame();
996         rah.finishFrame();
997     }
998     if (bDGsol)
999     {
1000         dgh.finishFrame();
1001     }
1002
1003     if (vh.isValid())
1004     {
1005         real totmass = 0;
1006         for (int i = 0; i < surfaceSel.posCount(); ++i)
1007         {
1008             totmass += surfaceSel.position(i).mass();
1009         }
1010         const real density = totmass*AMU/(totvolume*NANO*NANO*NANO);
1011         vh.startFrame(frnr, fr.time);
1012         vh.setPoint(0, totvolume);
1013         vh.setPoint(1, density);
1014         vh.finishFrame();
1015     }
1016 }
1017
1018 void
1019 Sasa::finishAnalysis(int /*nframes*/)
1020 {
1021     //if (bITP)
1022     //{
1023     //    fp3 = ftp2FILE(efITP, nfile, fnm, "w");
1024     //    fprintf(fp3, "[ position_restraints ]\n"
1025     //            "#define FCX 1000\n"
1026     //            "#define FCY 1000\n"
1027     //            "#define FCZ 1000\n"
1028     //            "; Atom  Type  fx   fy   fz\n");
1029     //    for (i = 0; i < nx[0]; i++)
1030     //    {
1031     //        if (atom_area[i] > minarea)
1032     //        {
1033     //            fprintf(fp3, "%5d   1     FCX  FCX  FCZ\n", ii+1);
1034     //        }
1035     //    }
1036     //    ffclose(fp3);
1037     //}
1038 }
1039
1040 void
1041 Sasa::writeOutput()
1042 {
1043 }
1044
1045 //! \}
1046
1047 }       // namespace
1048
1049 const char SasaInfo::name[]             = "sasa";
1050 const char SasaInfo::shortDescription[] =
1051     "Compute solvent accessible surface area";
1052
1053 TrajectoryAnalysisModulePointer SasaInfo::create()
1054 {
1055     return TrajectoryAnalysisModulePointer(new Sasa);
1056 }
1057
1058 } // namespace analysismodules
1059
1060 } // namespace gmx