d649a874a4f396552de2a8aebf8588eba700ecb8
[alexxy/gromacs.git] / src / gromacs / selection / sm_insolidangle.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2009,2010,2011,2012,2013 by the GROMACS development team.
5  * Copyright (c) 2014,2015,2016,2017,2018 by the GROMACS development team.
6  * Copyright (c) 2019,2020, 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
38  * \page page_module_selection_insolidangle Selection method: insolidangle
39  *
40  * This method selects a subset of particles that are located in a solid
41  * angle defined by a center and a set of points.
42  * The solid angle is constructed as a union of small cones whose axis
43  * goes through the center and a point.
44  * So there's such a cone for each position, and a
45  * point is in the solid angle if it lies within any of these cones.
46  * The width of the cones can be adjusted.
47  *
48  * The method is implemented by partitioning the surface of the unit sphere
49  * into bins using the polar coordinates \f$(\theta, \phi)\f$.
50  * The partitioning is always uniform in the zenith angle \f$\theta\f$,
51  * while the partitioning in the azimuthal angle \f$\phi\f$ varies.
52  * For each reference point, the unit vector from the center to the point
53  * is constructed, and it is stored in all the bins that overlap with the
54  * cone defined by the point.
55  * Bins that are completely covered by a single cone are marked as such.
56  * Checking whether a point is in the solid angle is then straightforward
57  * with this data structure: one finds the bin that corresponds to the point,
58  * and checks whether the bin is completely covered. If it is not, one
59  * additionally needs to check whether it is within the specified cutoff of
60  * any of the stored points.
61  *
62  * The above construction gives quite a lot of flexibility for constructing
63  * the bins without modifying the rest of the code.
64  * The current (quite inefficient) implementation is discussed below, but
65  * it should be optimized to get the most out of the code.
66  *
67  * The current way of constructing the bins constructs the boundaries
68  * statically: the bin size in the zenith direction is set to approximately
69  * half the angle cutoff, and the bins in the azimuthal direction have
70  * sizes such that the shortest edge of the bin is approximately equal to
71  * half the angle cutoff (for the regions close to the poles, a single bin
72  * is used).
73  * Each reference point is then added to the bins as follows:
74  *  -# Find the zenith angle range that is spanned by the cone centered at the
75  *     point (this is simple addition/subtraction).
76  *  -# Calculate the maximal span of the cone in the azimuthal direction using
77  *     the formula
78  *     \f[\sin \Delta \phi_{max} = \frac{\sin \alpha}{\sin \theta}\f]
79  *     (a sine formula in spherical coordinates),
80  *     where \f$\alpha\f$ is the width of the cone and \f$\theta\f$ is the
81  *     zenith angle of the cone center.
82  *     Similarly, the zenith angle at which this extent is achieved is
83  *     calculated using
84  *     \f[\cos \theta_{max} = \frac{\cos \theta}{\cos \alpha}\f]
85  *     (Pythagoras's theorem in spherical coordinates).
86  *  -# For each zenith angle bin that is at least partially covered by the
87  *     cone, calculate the span of the cone at the edges using
88  *     \f[\sin^2 \frac{\Delta \phi}{2} = \frac{\sin^2 \frac{\alpha}{2} - \sin^2 \frac{\theta -
89  * \theta'}{2}}{\sin \theta \sin \theta'}\f] (distance in spherical geometry), where \f$\theta'\f$
90  * is the zenith angle of the bin edge. Treat zenith angle bins that are completely covered by the
91  * cone (in the case that the cone is centered close to the pole) as a special case.
92  *  -# Using the values calculated above, loop through the azimuthal bins that
93  *     are partially or completely covered by the cone and update them.
94  *
95  * The total solid angle (for covered fraction calculations) is estimated by
96  * taking the total area of completely covered bins plus
97  * half the area of partially covered bins.
98  * The second one is an approximation, but should give reasonable estimates
99  * for the averages as well as in cases where the bin size is small.
100  */
101 /*! \internal \file
102  * \brief
103  * Implements the \ref sm_insolidangle "insolidangle" selection method.
104  *
105  * \todo
106  * The implementation could be optimized quite a bit.
107  *
108  * \todo
109  * Move the covered fraction stuff somewhere else and make it more generic
110  * (along the lines it is handled in selection.h and trajana.h in the old C
111  * API).
112  *
113  * \author Teemu Murtola <teemu.murtola@gmail.com>
114  * \ingroup module_selection
115  */
116 #include "gmxpre.h"
117
118 #include <cmath>
119
120 #include <algorithm>
121
122 #include "gromacs/math/functions.h"
123 #include "gromacs/math/units.h"
124 #include "gromacs/math/utilities.h"
125 #include "gromacs/math/vec.h"
126 #include "gromacs/pbcutil/pbc.h"
127 #include "gromacs/selection/indexutil.h"
128 #include "gromacs/selection/selection.h"
129 #include "gromacs/utility/arraysize.h"
130 #include "gromacs/utility/exceptions.h"
131 #include "gromacs/utility/smalloc.h"
132
133 #include "position.h"
134 #include "selelem.h"
135 #include "selmethod.h"
136 #include "selmethod_impl.h"
137
138 using std::max;
139 using std::min;
140
141 /*! \internal
142  * \brief
143  * Internal data structure for the \p insolidangle selection method.
144  *
145  * \see \c t_partition
146  *
147  * \ingroup module_selection
148  */
149 typedef struct
150 {
151     /** Left edge of the partition. */
152     real left;
153     /** Bin index corresponding to this partition. */
154     int bin;
155 } t_partition_item;
156
157 /*! \internal
158  * \brief
159  * Internal data structure for the \p insolidangle selection method.
160  *
161  * Describes the surface partitioning within one slice along the zenith angle.
162  * The slice from azimuthal angle \p p[i].left to \p p[i+1].left belongs to
163  * bin \p p[i].bin.
164  *
165  * \ingroup module_selection
166  */
167 typedef struct partition
168 {
169     /** Number of partition items (\p p contains \p n+1 items). */
170     int n;
171     /** Array of partition edges and corresponding bins. */
172     t_partition_item* p;
173 } t_partition;
174
175 /*! \internal
176  * \brief
177  * Internal data structure for the \p insolidangle selection method.
178  *
179  * Contains the reference points that partially cover a certain region on the
180  * surface of the unit sphere.
181  * If \p n is -1, the whole region described by the bin is covered.
182  *
183  * \ingroup module_selection
184  */
185 typedef struct spheresurfacebin
186 {
187     /** Number of points in the array \p x, -1 if whole bin covered. */
188     int n;
189     /** Number of elements allocated for \p x. */
190     int n_alloc;
191     /** Array of points that partially cover the bin. */
192     rvec* x;
193 } t_spheresurfacebin;
194
195 /*! \internal
196  * \brief
197  * Data structure for the \p insolidangle selection method.
198  *
199  * All angle values are in the units of radians.
200  *
201  * \ingroup module_selection
202  */
203 typedef struct methoddata_insolidangle
204 {
205     /** Center of the solid angle. */
206     gmx_ana_pos_t center;
207     /** Positions that span the solid angle. */
208     gmx_ana_pos_t span;
209     /** Cutoff angle. */
210     real angcut;
211     /** Estimate of the covered fraction. */
212     real cfrac;
213
214     /** Cutoff for the cosine (equals cos(angcut)). */
215     real distccut;
216     /** Bin size to be used as the target bin size when constructing the bins. */
217     real targetbinsize;
218
219     /** Number of bins in the \p tbin array. */
220     int ntbins;
221     /** Size of one bin in the zenith angle direction. */
222     real tbinsize;
223     /** Array of zenith angle slices. */
224     t_partition* tbin;
225     /** Number of elements allocated for the \p bin array. */
226     int maxbins;
227     /** Number of elements used in the \p bin array. */
228     int nbins;
229     /** Array of individual bins. */
230     t_spheresurfacebin* bin;
231 } t_methoddata_insolidangle;
232
233 /*! \brief
234  * Allocates data for the \p insolidangle selection method.
235  *
236  * \param[in]     npar  Not used (should be 3).
237  * \param[in,out] param Method parameters (should point to
238  *   \ref smparams_insolidangle).
239  * \returns Pointer to the allocated data (\ref t_methoddata_insolidangle).
240  *
241  * Allocates memory for a \ref t_methoddata_insolidangle structure and
242  * initializes the parameter as follows:
243  *  - \p center defines the value for t_methoddata_insolidangle::center.
244  *  - \p span   defines the value for t_methoddata_insolidangle::span.
245  *  - \p cutoff defines the value for t_methoddata_insolidangle::angcut.
246  */
247 static void* init_data_insolidangle(int npar, gmx_ana_selparam_t* param);
248 /*! \brief
249  * Initializes the \p insolidangle selection method.
250  *
251  * \param   top  Not used.
252  * \param   npar Not used.
253  * \param   param Not used.
254  * \param   data Pointer to \ref t_methoddata_insolidangle to initialize.
255  * \returns 0 on success, -1 on failure.
256  *
257  * Converts t_methoddata_insolidangle::angcut to radians and allocates
258  * and allocates memory for the bins used during the evaluation.
259  */
260 static void init_insolidangle(const gmx_mtop_t* top, int npar, gmx_ana_selparam_t* param, void* data);
261 /** Frees the data allocated for the \p insolidangle selection method. */
262 static void free_data_insolidangle(void* data);
263 /*! \brief
264  * Initializes the evaluation of the \p insolidangle selection method for a frame.
265  *
266  * \param[in]  context Evaluation context.
267  * \param      data    Should point to a \ref t_methoddata_insolidangle.
268  *
269  * Creates a lookup structure that enables fast queries of whether a point
270  * is within the solid angle or not.
271  */
272 static void init_frame_insolidangle(const gmx::SelMethodEvalContext& context, void* data);
273 /** Internal helper function for evaluate_insolidangle(). */
274 static bool accept_insolidangle(rvec x, const t_pbc* pbc, void* data);
275 /** Evaluates the \p insolidangle selection method. */
276 static void evaluate_insolidangle(const gmx::SelMethodEvalContext& context,
277                                   gmx_ana_pos_t*                   pos,
278                                   gmx_ana_selvalue_t*              out,
279                                   void*                            data);
280
281 /** Calculates the distance between unit vectors. */
282 static real sph_distc(rvec x1, rvec x2);
283 /** Does a binary search on a \p t_partition to find a bin for a value. */
284 static int find_partition_bin(t_partition* p, real value);
285 /** Finds a bin that corresponds to a location on the unit sphere surface. */
286 static int find_surface_bin(t_methoddata_insolidangle* surf, rvec x);
287 /** Clears/initializes the bins on the unit sphere surface. */
288 static void clear_surface_points(t_methoddata_insolidangle* surf);
289 /** Frees memory allocated for storing the reference points in the surface bins. */
290 static void free_surface_points(t_methoddata_insolidangle* surf);
291 /** Adds a reference point to a given bin. */
292 static void add_surface_point(t_methoddata_insolidangle* surf, int tbin, int pbin, rvec x);
293 /** Marks a bin as completely covered. */
294 static void mark_surface_covered(t_methoddata_insolidangle* surf, int tbin, int pbin);
295 /** Helper function for store_surface_point() to update a single zenith angle bin. */
296 static void update_surface_bin(t_methoddata_insolidangle* surf,
297                                int                        tbin,
298                                real                       phi,
299                                real                       pdelta1,
300                                real                       pdelta2,
301                                real                       pdeltamax,
302                                rvec                       x);
303 /** Adds a single reference point and updates the surface bins. */
304 static void store_surface_point(t_methoddata_insolidangle* surf, rvec x);
305 /*! \brief
306  * Optimizes the surface bins for faster searching.
307  *
308  * \param[in,out] surf Surface data structure.
309  *
310  * Currently, this function does nothing.
311  */
312 static void optimize_surface_points(t_methoddata_insolidangle* surf);
313 /** Estimates the area covered by the reference cones. */
314 static real estimate_covered_fraction(t_methoddata_insolidangle* surf);
315 /** Checks whether a point lies within a solid angle. */
316 static bool is_surface_covered(t_methoddata_insolidangle* surf, rvec x);
317
318 /** Parameters for the \p insolidangle selection method. */
319 static gmx_ana_selparam_t smparams_insolidangle[] = {
320     { "center", { POS_VALUE, 1, { nullptr } }, nullptr, SPAR_DYNAMIC },
321     { "span", { POS_VALUE, -1, { nullptr } }, nullptr, SPAR_DYNAMIC | SPAR_VARNUM },
322     { "cutoff", { REAL_VALUE, 1, { nullptr } }, nullptr, SPAR_OPTIONAL },
323 };
324
325 /** Help text for the \p insolidangle selection method. */
326 static const char* const help_insolidangle[] = {
327     "::",
328     "",
329     "  insolidangle center POS span POS_EXPR [cutoff REAL]",
330     "",
331     "This keyword selects atoms that are within [TT]REAL[tt] degrees",
332     "(default=5) of any position in [TT]POS_EXPR[tt] as seen from [TT]POS[tt]",
333     "a position expression that evaluates to a single position), i.e., atoms",
334     "in the solid angle spanned by the positions in [TT]POS_EXPR[tt] and",
335     "centered at [TT]POS[tt].[PAR]",
336
337     "Technically, the solid angle is constructed as a union of small cones",
338     "whose tip is at [TT]POS[tt] and the axis goes through a point in",
339     "[TT]POS_EXPR[tt]. There is such a cone for each position in",
340     "[TT]POS_EXPR[tt], and point is in the solid angle if it lies within any",
341     "of these cones. The cutoff determines the width of the cones.",
342 };
343
344 /** Selection method data for the \p insolidangle method. */
345 gmx_ana_selmethod_t sm_insolidangle = {
346     "insolidangle",
347     GROUP_VALUE,
348     SMETH_DYNAMIC,
349     asize(smparams_insolidangle),
350     smparams_insolidangle,
351     &init_data_insolidangle,
352     nullptr,
353     &init_insolidangle,
354     nullptr,
355     &free_data_insolidangle,
356     &init_frame_insolidangle,
357     nullptr,
358     &evaluate_insolidangle,
359     { "insolidangle center POS span POS_EXPR [cutoff REAL]", "Selecting atoms in a solid angle",
360       asize(help_insolidangle), help_insolidangle },
361 };
362
363 static void* init_data_insolidangle(int /* npar */, gmx_ana_selparam_t* param)
364 {
365     t_methoddata_insolidangle* data = new t_methoddata_insolidangle();
366     data->angcut                    = 5.0;
367     data->cfrac                     = 0.0;
368
369     data->distccut      = 0.0;
370     data->targetbinsize = 0.0;
371
372     data->ntbins   = 0;
373     data->tbinsize = 0.0;
374     data->tbin     = nullptr;
375     data->maxbins  = 0;
376     data->nbins    = 0;
377     data->bin      = nullptr;
378
379     param[0].val.u.p = &data->center;
380     param[1].val.u.p = &data->span;
381     param[2].val.u.r = &data->angcut;
382     return data;
383 }
384
385 static void init_insolidangle(const gmx_mtop_t* /* top */,
386                               int /* npar */,
387                               gmx_ana_selparam_t* /* param */,
388                               void* data)
389 {
390     t_methoddata_insolidangle* surf = static_cast<t_methoddata_insolidangle*>(data);
391     int                        i, c;
392
393     if (surf->angcut <= 0)
394     {
395         GMX_THROW(gmx::InvalidInputError("Angle cutoff should be > 0"));
396     }
397
398     surf->angcut *= DEG2RAD;
399
400     surf->distccut      = -std::cos(surf->angcut);
401     surf->targetbinsize = surf->angcut / 2;
402     surf->ntbins        = static_cast<int>(M_PI / surf->targetbinsize);
403     surf->tbinsize      = (180.0 / surf->ntbins) * DEG2RAD;
404
405     snew(surf->tbin, static_cast<int>(M_PI / surf->tbinsize) + 1);
406     surf->maxbins = 0;
407     for (i = 0; i < surf->ntbins; ++i)
408     {
409         c = static_cast<int>(std::max(std::sin(surf->tbinsize * i), std::sin(surf->tbinsize * (i + 1)))
410                              * M_2PI / surf->targetbinsize)
411             + 1;
412         snew(surf->tbin[i].p, c + 1);
413         surf->maxbins += c;
414     }
415     surf->nbins = 0;
416     snew(surf->bin, surf->maxbins);
417 }
418
419 /*!
420  * \param data Data to free (should point to a \ref t_methoddata_insolidangle).
421  *
422  * Frees the memory allocated for \c t_methoddata_insolidangle::center and
423  * \c t_methoddata_insolidangle::span, as well as the memory for the internal
424  * bin structure.
425  */
426 static void free_data_insolidangle(void* data)
427 {
428     t_methoddata_insolidangle* d = static_cast<t_methoddata_insolidangle*>(data);
429     int                        i;
430
431     if (d->tbin)
432     {
433         for (i = 0; i < d->ntbins; ++i)
434         {
435             sfree(d->tbin[i].p);
436         }
437         sfree(d->tbin);
438     }
439     free_surface_points(d);
440     sfree(d->bin);
441     delete d;
442 }
443
444 static void init_frame_insolidangle(const gmx::SelMethodEvalContext& context, void* data)
445 {
446     t_methoddata_insolidangle* d = static_cast<t_methoddata_insolidangle*>(data);
447     rvec                       dx;
448     int                        i;
449
450     free_surface_points(d);
451     clear_surface_points(d);
452     for (i = 0; i < d->span.count(); ++i)
453     {
454         if (context.pbc)
455         {
456             pbc_dx(context.pbc, d->span.x[i], d->center.x[0], dx);
457         }
458         else
459         {
460             rvec_sub(d->span.x[i], d->center.x[0], dx);
461         }
462         unitv(dx, dx);
463         store_surface_point(d, dx);
464     }
465     optimize_surface_points(d);
466     d->cfrac = -1;
467 }
468
469 /*!
470  * \param[in] x    Test point.
471  * \param[in] pbc  PBC data (if NULL, no PBC are used).
472  * \param[in] data Pointer to a \c t_methoddata_insolidangle data structure.
473  * \returns   true if \p x is within the solid angle, false otherwise.
474  */
475 static bool accept_insolidangle(rvec x, const t_pbc* pbc, void* data)
476 {
477     t_methoddata_insolidangle* d = static_cast<t_methoddata_insolidangle*>(data);
478     rvec                       dx;
479
480     if (pbc)
481     {
482         pbc_dx(pbc, x, d->center.x[0], dx);
483     }
484     else
485     {
486         rvec_sub(x, d->center.x[0], dx);
487     }
488     unitv(dx, dx);
489     return is_surface_covered(d, dx);
490 }
491
492 /*!
493  * See sel_updatefunc() for description of the parameters.
494  * \p data should point to a \c t_methoddata_insolidangle.
495  *
496  * Calculates which atoms in \p g are within the solid angle spanned by
497  * \c t_methoddata_insolidangle::span and centered at
498  * \c t_methoddata_insolidangle::center, and stores the result in \p out->u.g.
499  */
500 static void evaluate_insolidangle(const gmx::SelMethodEvalContext& context,
501                                   gmx_ana_pos_t*                   pos,
502                                   gmx_ana_selvalue_t*              out,
503                                   void*                            data)
504 {
505     out->u.g->isize = 0;
506     for (int b = 0; b < pos->count(); ++b)
507     {
508         if (accept_insolidangle(pos->x[b], context.pbc, data))
509         {
510             gmx_ana_pos_add_to_group(out->u.g, pos, b);
511         }
512     }
513 }
514
515 /*!
516  * \param[in] sel Selection element to query.
517  * \returns   true if the covered fraction can be estimated for \p sel with
518  *   _gmx_selelem_estimate_coverfrac(), false otherwise.
519  */
520 bool _gmx_selelem_can_estimate_cover(const gmx::SelectionTreeElement& sel)
521 {
522     if (sel.type == SEL_BOOLEAN && sel.u.boolt == BOOL_OR)
523     {
524         return false;
525     }
526     bool                             bFound    = false;
527     bool                             bDynFound = false;
528     gmx::SelectionTreeElementPointer child     = sel.child;
529     while (child)
530     {
531         if (child->type == SEL_EXPRESSION)
532         {
533             if (child->u.expr.method->name == sm_insolidangle.name)
534             {
535                 if (bFound || bDynFound)
536                 {
537                     return false;
538                 }
539                 bFound = true;
540             }
541             else if (child->u.expr.method && (child->u.expr.method->flags & SMETH_DYNAMIC))
542             {
543                 if (bFound)
544                 {
545                     return false;
546                 }
547                 bDynFound = true;
548             }
549         }
550         else if (!_gmx_selelem_can_estimate_cover(*child))
551         {
552             return false;
553         }
554         child = child->next;
555     }
556     return true;
557 }
558
559 /*!
560  * \param[in] sel Selection for which the fraction should be calculated.
561  * \returns Fraction of angles covered by the selection (between zero and one).
562  *
563  * The return value is undefined if _gmx_selelem_can_estimate_cover() returns
564  * false.
565  * Should be called after gmx_ana_evaluate_selections() has been called for the
566  * frame.
567  */
568 real _gmx_selelem_estimate_coverfrac(const gmx::SelectionTreeElement& sel)
569 {
570     real cfrac;
571
572     if (sel.type == SEL_EXPRESSION && sel.u.expr.method->name == sm_insolidangle.name)
573     {
574         t_methoddata_insolidangle* d = static_cast<t_methoddata_insolidangle*>(sel.u.expr.mdata);
575         if (d->cfrac < 0)
576         {
577             d->cfrac = estimate_covered_fraction(d);
578         }
579         return d->cfrac;
580     }
581     if (sel.type == SEL_BOOLEAN && sel.u.boolt == BOOL_NOT)
582     {
583         cfrac = _gmx_selelem_estimate_coverfrac(*sel.child);
584         if (cfrac < 1.0)
585         {
586             return 1 - cfrac;
587         }
588         return 1;
589     }
590
591     /* Here, we assume that the selection is simple enough */
592     gmx::SelectionTreeElementPointer child = sel.child;
593     while (child)
594     {
595         cfrac = _gmx_selelem_estimate_coverfrac(*child);
596         if (cfrac < 1.0)
597         {
598             return cfrac;
599         }
600         child = child->next;
601     }
602     return 1.0;
603 }
604
605 /*!
606  * \param[in] x1  Unit vector 1.
607  * \param[in] x2  Unit vector 2.
608  * \returns   Minus the dot product of \p x1 and \p x2.
609  *
610  * This function is used internally to calculate the distance between the
611  * unit vectors \p x1 and \p x2 to find out whether \p x2 is within the
612  * cone centered at \p x1. Currently, the cosine of the angle is used
613  * for efficiency, and the minus is there to make it behave like a normal
614  * distance (larger values mean longer distances).
615  */
616 static real sph_distc(rvec x1, rvec x2)
617 {
618     return -iprod(x1, x2);
619 }
620
621 /*!
622  * \param[in] p     Partition to search.
623  * \param[in] value Value to search for.
624  * \returns   The partition index in \p p that contains \p value.
625  *
626  * If \p value is outside the range of \p p, the first/last index is returned.
627  * Otherwise, the return value \c i satisfies \c p->p[i].left<=value and
628  * \c p->p[i+1].left>value
629  */
630 static int find_partition_bin(t_partition* p, real value)
631 {
632     int pmin, pmax, pbin;
633
634     /* Binary search the partition */
635     pmin = 0;
636     pmax = p->n;
637     while (pmax > pmin + 1)
638     {
639         pbin = pmin + (pmax - pmin) / 2;
640         if (p->p[pbin].left <= value)
641         {
642             pmin = pbin;
643         }
644         else
645         {
646             pmax = pbin;
647         }
648     }
649     pbin = pmin;
650     return pbin;
651 }
652
653 /*!
654  * \param[in] surf  Surface data structure to search.
655  * \param[in] x     Unit vector to find.
656  * \returns   The bin index that contains \p x.
657  *
658  * The return value is an index to the \p surf->bin array.
659  */
660 static int find_surface_bin(t_methoddata_insolidangle* surf, rvec x)
661 {
662     real theta, phi;
663     int  tbin, pbin;
664
665     theta = acos(x[ZZ]);
666     phi   = atan2(x[YY], x[XX]);
667     tbin  = static_cast<int>(std::floor(theta / surf->tbinsize));
668     if (tbin >= surf->ntbins)
669     {
670         tbin = surf->ntbins - 1;
671     }
672     pbin = find_partition_bin(&surf->tbin[tbin], phi);
673     return surf->tbin[tbin].p[pbin].bin;
674 }
675
676 /*!
677  * \param[in,out] surf Surface data structure.
678  *
679  * Clears the reference points from the bins and (re)initializes the edges
680  * of the azimuthal bins.
681  */
682 static void clear_surface_points(t_methoddata_insolidangle* surf)
683 {
684     int i, j, c;
685
686     surf->nbins = 0;
687     for (i = 0; i < surf->ntbins; ++i)
688     {
689         c = static_cast<int>(std::min(std::sin(surf->tbinsize * i), std::sin(surf->tbinsize * (i + 1)))
690                              * M_2PI / surf->targetbinsize)
691             + 1;
692         if (c <= 0)
693         {
694             c = 1;
695         }
696         surf->tbin[i].n = c;
697         for (j = 0; j < c; ++j)
698         {
699             surf->tbin[i].p[j].left  = -M_PI + j * M_2PI / c - 0.0001;
700             surf->tbin[i].p[j].bin   = surf->nbins;
701             surf->bin[surf->nbins].n = 0;
702             surf->nbins++;
703         }
704         surf->tbin[i].p[c].left = M_PI + 0.0001;
705         surf->tbin[i].p[c].bin  = -1;
706     }
707 }
708
709 /*!
710  * \param[in,out] surf Surface data structure.
711  */
712 static void free_surface_points(t_methoddata_insolidangle* surf)
713 {
714     int i;
715
716     for (i = 0; i < surf->nbins; ++i)
717     {
718         if (surf->bin[i].x)
719         {
720             sfree(surf->bin[i].x);
721         }
722         surf->bin[i].n_alloc = 0;
723         surf->bin[i].x       = nullptr;
724     }
725 }
726
727 /*!
728  * \param[in,out] surf Surface data structure.
729  * \param[in]     tbin Bin number in the zenith angle direction.
730  * \param[in]     pbin Bin number in the azimuthal angle direction.
731  * \param[in]     x    Point to store.
732  */
733 static void add_surface_point(t_methoddata_insolidangle* surf, int tbin, int pbin, rvec x)
734 {
735     int bin;
736
737     bin = surf->tbin[tbin].p[pbin].bin;
738     /* Return if bin is already completely covered */
739     if (surf->bin[bin].n == -1)
740     {
741         return;
742     }
743     /* Allocate more space if necessary */
744     if (surf->bin[bin].n == surf->bin[bin].n_alloc)
745     {
746         surf->bin[bin].n_alloc += 10;
747         srenew(surf->bin[bin].x, surf->bin[bin].n_alloc);
748     }
749     /* Add the point to the bin */
750     copy_rvec(x, surf->bin[bin].x[surf->bin[bin].n]);
751     ++surf->bin[bin].n;
752 }
753
754 /*!
755  * \param[in,out] surf Surface data structure.
756  * \param[in]     tbin Bin number in the zenith angle direction.
757  * \param[in]     pbin Bin number in the azimuthal angle direction.
758  */
759 static void mark_surface_covered(t_methoddata_insolidangle* surf, int tbin, int pbin)
760 {
761     int bin;
762
763     bin              = surf->tbin[tbin].p[pbin].bin;
764     surf->bin[bin].n = -1;
765 }
766
767 /*!
768  * \param[in,out] surf      Surface data structure.
769  * \param[in]     tbin      Bin number in the zenith angle direction.
770  * \param[in]     phi       Azimuthal angle of \p x.
771  * \param[in]     pdelta1   Width of the cone at the lower edge of \p tbin.
772  * \param[in]     pdelta2   Width of the cone at the uppper edge of \p tbin.
773  * \param[in]     pdeltamax Max. width of the cone inside \p tbin.
774  * \param[in]     x         Point to store (should have unit length).
775  */
776 static void update_surface_bin(t_methoddata_insolidangle* surf,
777                                int                        tbin,
778                                real                       phi,
779                                real                       pdelta1,
780                                real                       pdelta2,
781                                real                       pdeltamax,
782                                rvec                       x)
783 {
784     real pdelta, phi1, phi2;
785     int  pbin1, pbin2, pbiniter, pbin;
786
787     /* Find the edges of the bins affected */
788     pdelta = max(max(pdelta1, pdelta2), pdeltamax);
789     phi1   = phi - pdelta;
790     if (phi1 >= -M_PI)
791     {
792         pbin  = find_partition_bin(&surf->tbin[tbin], phi1);
793         pbin1 = pbin;
794     }
795     else
796     {
797         pbin  = find_partition_bin(&surf->tbin[tbin], phi1 + M_2PI);
798         pbin1 = pbin - surf->tbin[tbin].n;
799     }
800     phi2 = phi + pdelta;
801     if (phi2 <= M_PI)
802     {
803         pbin2 = find_partition_bin(&surf->tbin[tbin], phi2);
804     }
805     else
806     {
807         pbin2 = find_partition_bin(&surf->tbin[tbin], phi2 - M_2PI);
808         pbin2 += surf->tbin[tbin].n;
809     }
810     ++pbin2;
811     if (pbin2 - pbin1 > surf->tbin[tbin].n)
812     {
813         pbin2 = pbin1 + surf->tbin[tbin].n;
814     }
815     /* Find the edges of completely covered region */
816     pdelta = min(pdelta1, pdelta2);
817     phi1   = phi - pdelta;
818     if (phi1 < -M_PI)
819     {
820         phi1 += M_2PI;
821     }
822     phi2 = phi + pdelta;
823     /* Loop over all affected bins */
824     for (pbiniter = pbin1; pbiniter != pbin2; ++pbiniter, ++pbin)
825     {
826         /* Wrap bin around if end reached */
827         if (pbin == surf->tbin[tbin].n)
828         {
829             pbin = 0;
830             phi1 -= M_2PI;
831             phi2 -= M_2PI;
832         }
833         /* Check if bin is completely covered and update */
834         if (surf->tbin[tbin].p[pbin].left >= phi1 && surf->tbin[tbin].p[pbin + 1].left <= phi2)
835         {
836             mark_surface_covered(surf, tbin, pbin);
837         }
838         else
839         {
840             add_surface_point(surf, tbin, pbin, x);
841         }
842     }
843 }
844
845 /*!
846  * \param[in,out] surf Surface data structure.
847  * \param[in]     x    Point to store (should have unit length).
848  *
849  * Finds all the bins covered by the cone centered at \p x and calls
850  * update_surface_bin() to update them.
851  */
852 static void store_surface_point(t_methoddata_insolidangle* surf, rvec x)
853 {
854     real theta, phi;
855     real pdeltamax, tmax;
856     real theta1, theta2, pdelta1, pdelta2;
857     int  tbin;
858
859     theta = acos(x[ZZ]);
860     phi   = atan2(x[YY], x[XX]);
861     /* Find the maximum extent in the phi direction */
862     if (theta <= surf->angcut)
863     {
864         pdeltamax = M_PI;
865         tmax      = 0;
866     }
867     else if (theta >= M_PI - surf->angcut)
868     {
869         pdeltamax = M_PI;
870         tmax      = M_PI;
871     }
872     else
873     {
874         pdeltamax = std::asin(sin(surf->angcut) / sin(theta));
875         tmax      = std::acos(cos(theta) / cos(surf->angcut));
876     }
877     /* Find the first affected bin */
878     tbin   = max(static_cast<int>(std::floor((theta - surf->angcut) / surf->tbinsize)), 0);
879     theta1 = tbin * surf->tbinsize;
880     if (theta1 < theta - surf->angcut)
881     {
882         pdelta1 = 0;
883     }
884     else
885     {
886         pdelta1 = M_PI;
887     }
888     /* Loop through all affected bins */
889     while (tbin < std::ceil((theta + surf->angcut) / surf->tbinsize) && tbin < surf->ntbins)
890     {
891         /* Calculate the next boundaries */
892         theta2 = (tbin + 1) * surf->tbinsize;
893         if (theta2 > theta + surf->angcut)
894         {
895             /* The circle is completely outside the cone */
896             pdelta2 = 0;
897         }
898         else if (theta2 <= -(theta - surf->angcut) || theta2 >= M_2PI - (theta + surf->angcut)
899                  || tbin == surf->ntbins - 1)
900         {
901             /* The circle is completely inside the cone, or we are in the
902              * 360 degree bin covering the pole. */
903             pdelta2 = M_PI;
904         }
905         else
906         {
907             /* TODO: This formula is numerically unstable if theta is very
908              * close to the pole.  In practice, it probably does not matter
909              * much, but it would be nicer to adjust the theta bin boundaries
910              * such that the case above catches this instead of falling through
911              * here. */
912             pdelta2 = 2
913                       * asin(std::sqrt((gmx::square(std::sin(surf->angcut / 2))
914                                         - gmx::square(std::sin((theta2 - theta) / 2)))
915                                        / (sin(theta) * sin(theta2))));
916         }
917         /* Update the bin */
918         if (tmax >= theta1 && tmax <= theta2)
919         {
920             update_surface_bin(surf, tbin, phi, pdelta1, pdelta2, pdeltamax, x);
921         }
922         else
923         {
924             update_surface_bin(surf, tbin, phi, pdelta1, pdelta2, 0, x);
925         }
926         /* Next bin */
927         theta1  = theta2;
928         pdelta1 = pdelta2;
929         ++tbin;
930     }
931 }
932
933 static void optimize_surface_points(t_methoddata_insolidangle* /* surf */)
934 {
935     /* TODO: Implement */
936 }
937
938 /*!
939  * \param[in] surf Surface data structure.
940  * \returns   An estimate for the area covered by the reference points.
941  */
942 static real estimate_covered_fraction(t_methoddata_insolidangle* surf)
943 {
944     int  t, p, n;
945     real cfrac, tfrac, pfrac;
946
947     cfrac = 0.0;
948     for (t = 0; t < surf->ntbins; ++t)
949     {
950         tfrac = std::cos(t * surf->tbinsize) - std::cos((t + 1) * surf->tbinsize);
951         for (p = 0; p < surf->tbin[t].n; ++p)
952         {
953             pfrac = surf->tbin[t].p[p + 1].left - surf->tbin[t].p[p].left;
954             n     = surf->bin[surf->tbin[t].p[p].bin].n;
955             if (n == -1) /* Bin completely covered */
956             {
957                 cfrac += tfrac * pfrac;
958             }
959             else if (n > 0) /* Bin partially covered */
960             {
961                 cfrac += tfrac * pfrac / 2; /* A rough estimate */
962             }
963         }
964     }
965     return cfrac / (4 * M_PI);
966 }
967
968 /*!
969  * \param[in] surf  Surface data structure to search.
970  * \param[in] x     Unit vector to check.
971  * \returns   true if \p x is within the solid angle, false otherwise.
972  */
973 static bool is_surface_covered(t_methoddata_insolidangle* surf, rvec x)
974 {
975     int bin, i;
976
977     bin = find_surface_bin(surf, x);
978     /* Check for completely covered bin */
979     if (surf->bin[bin].n == -1)
980     {
981         return true;
982     }
983     /* Check each point that partially covers the bin */
984     for (i = 0; i < surf->bin[bin].n; ++i)
985     {
986         if (sph_distc(x, surf->bin[bin].x[i]) < surf->distccut)
987         {
988             return true;
989         }
990     }
991     return false;
992 }