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