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