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