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