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