2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2009,2010,2011,2012,2013, 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.
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.
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.
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.
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.
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.
35 /*! \page page_module_selection_insolidangle Selection method: insolidangle
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.
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.
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.
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
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
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
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.
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.
103 * Implements the \ref sm_insolidangle "insolidangle" selection method.
106 * The implementation could be optimized quite a bit.
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
113 * \author Teemu Murtola <teemu.murtola@gmail.com>
114 * \ingroup module_selection
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"
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"
139 * Internal data structure for the \p insolidangle selection method.
141 * \see \c t_partition
145 /** Left edge of the partition. */
147 /** Bin index corresponding to this partition. */
152 * Internal data structure for the \p insolidangle selection method.
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
160 /** Number of partition items (\p p contains \p n+1 items). */
162 /** Array of partition edges and corresponding bins. */
167 * Internal data structure for the \p insolidangle selection method.
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.
175 /** Number of points in the array \p x, -1 if whole bin covered. */
177 /** Number of elements allocated for \p x. */
179 /** Array of points that partially cover the bin. */
181 } t_spheresurfacebin;
184 * Data structure for the \p insolidangle selection method.
186 * All angle values are in the units of radians.
190 /** Center of the solid angle. */
191 gmx_ana_pos_t center;
192 /** Positions that span the solid angle. */
196 /** Estimate of the covered fraction. */
199 /** Cutoff for the cosine (equals cos(angcut)). */
201 /** Bin size to be used as the target bin size when constructing the bins. */
204 /** Number of bins in the \p tbin array. */
206 /** Size of one bin in the zenith angle direction. */
208 /** Array of zenith angle slices. */
210 /** Number of elements allocated for the \p bin array. */
212 /** Number of elements used in the \p bin array. */
214 /** Array of individual bins. */
215 t_spheresurfacebin *bin;
216 } t_methoddata_insolidangle;
218 /** Allocates data for the \p insolidangle selection method. */
220 init_data_insolidangle(int npar, gmx_ana_selparam_t *param);
221 /** Initializes the \p insolidangle selection method. */
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. */
226 free_data_insolidangle(void *data);
227 /** Initializes the evaluation of the \p insolidangle selection method for a frame. */
229 init_frame_insolidangle(t_topology *top, t_trxframe *fr, t_pbc *pbc, void *data);
230 /** Internal helper function for evaluate_insolidangle(). */
232 accept_insolidangle(rvec x, t_pbc *pbc, void *data);
233 /** Evaluates the \p insolidangle selection method. */
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);
238 /** Calculates the distance between unit vectors. */
240 sph_distc(rvec x1, rvec x2);
241 /** Does a binary search on a \p t_partition to find a bin for a value. */
243 find_partition_bin(t_partition *p, real value);
244 /** Finds a bin that corresponds to a location on the unit sphere surface. */
246 find_surface_bin(t_methoddata_insolidangle *surf, rvec x);
247 /** Clears/initializes the bins on the unit sphere surface. */
249 clear_surface_points(t_methoddata_insolidangle *surf);
250 /** Frees memory allocated for storing the reference points in the surface bins. */
252 free_surface_points(t_methoddata_insolidangle *surf);
253 /** Adds a reference point to a given bin. */
255 add_surface_point(t_methoddata_insolidangle *surf, int tbin, int pbin, rvec x);
256 /** Marks a bin as completely covered. */
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. */
261 update_surface_bin(t_methoddata_insolidangle *surf, int tbin,
262 real phi, real pdelta1, real pdelta2, real pdeltamax,
264 /** Adds a single reference point and updates the surface bins. */
266 store_surface_point(t_methoddata_insolidangle *surf, rvec x);
267 /** Optimizes the surface bins for faster searching. */
269 optimize_surface_points(t_methoddata_insolidangle *surf);
270 /** Estimates the area covered by the reference cones. */
272 estimate_covered_fraction(t_methoddata_insolidangle *surf);
273 /** Checks whether a point lies within a solid angle. */
275 is_surface_covered(t_methoddata_insolidangle *surf, rvec x);
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},
284 /** Help text for the \p insolidangle selection method. */
285 static const char *help_insolidangle[] = {
286 "SELECTING ATOMS IN A SOLID ANGLE[PAR]",
288 "[TT]insolidangle center POS span POS_EXPR [cutoff REAL][tt][PAR]",
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]"
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.",
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,
311 &free_data_insolidangle,
312 &init_frame_insolidangle,
314 &evaluate_insolidangle,
315 {"insolidangle center POS span POS_EXPR [cutoff REAL]",
316 asize(help_insolidangle), help_insolidangle},
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).
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.
332 init_data_insolidangle(int npar, gmx_ana_selparam_t *param)
334 t_methoddata_insolidangle *data = new t_methoddata_insolidangle();
338 data->distccut = 0.0;
339 data->targetbinsize = 0.0;
342 data->tbinsize = 0.0;
348 param[0].val.u.p = &data->center;
349 param[1].val.u.p = &data->span;
350 param[2].val.u.r = &data->angcut;
355 * \param top Not used.
356 * \param npar Not used.
357 * \param param Not used.
358 * \param data Pointer to \ref t_methoddata_insolidangle to initialize.
359 * \returns 0 on success, -1 on failure.
361 * Converts t_methoddata_insolidangle::angcut to radians and allocates
362 * and allocates memory for the bins used during the evaluation.
365 init_insolidangle(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data)
367 t_methoddata_insolidangle *surf = (t_methoddata_insolidangle *)data;
370 if (surf->angcut <= 0)
372 GMX_THROW(gmx::InvalidInputError("Angle cutoff should be > 0"));
375 surf->angcut *= DEG2RAD;
377 surf->distccut = -cos(surf->angcut);
378 surf->targetbinsize = surf->angcut / 2;
379 surf->ntbins = static_cast<int>(M_PI / surf->targetbinsize);
380 surf->tbinsize = (180.0 / surf->ntbins)*DEG2RAD;
382 snew(surf->tbin, static_cast<int>(M_PI / surf->tbinsize) + 1);
384 for (i = 0; i < surf->ntbins; ++i)
386 c = static_cast<int>(max(sin(surf->tbinsize*i),
387 sin(surf->tbinsize*(i+1)))
388 * M_2PI / surf->targetbinsize) + 1;
389 snew(surf->tbin[i].p, c+1);
393 snew(surf->bin, surf->maxbins);
397 * \param data Data to free (should point to a \ref t_methoddata_insolidangle).
399 * Frees the memory allocated for \c t_methoddata_insolidangle::center and
400 * \c t_methoddata_insolidangle::span, as well as the memory for the internal
404 free_data_insolidangle(void *data)
406 t_methoddata_insolidangle *d = (t_methoddata_insolidangle *)data;
411 for (i = 0; i < d->ntbins; ++i)
417 free_surface_points(d);
423 * \param[in] top Not used.
424 * \param[in] fr Current frame.
425 * \param[in] pbc PBC structure.
426 * \param data Should point to a \ref t_methoddata_insolidangle.
428 * Creates a lookup structure that enables fast queries of whether a point
429 * is within the solid angle or not.
432 init_frame_insolidangle(t_topology *top, t_trxframe *fr, t_pbc *pbc, void *data)
434 t_methoddata_insolidangle *d = (t_methoddata_insolidangle *)data;
438 free_surface_points(d);
439 clear_surface_points(d);
440 for (i = 0; i < d->span.count(); ++i)
444 pbc_dx(pbc, d->span.x[i], d->center.x[0], dx);
448 rvec_sub(d->span.x[i], d->center.x[0], dx);
451 store_surface_point(d, dx);
453 optimize_surface_points(d);
458 * \param[in] x Test point.
459 * \param[in] pbc PBC data (if NULL, no PBC are used).
460 * \param[in] data Pointer to a \c t_methoddata_insolidangle data structure.
461 * \returns true if \p x is within the solid angle, false otherwise.
464 accept_insolidangle(rvec x, t_pbc *pbc, void *data)
466 t_methoddata_insolidangle *d = (t_methoddata_insolidangle *)data;
471 pbc_dx(pbc, x, d->center.x[0], dx);
475 rvec_sub(x, d->center.x[0], dx);
478 return is_surface_covered(d, dx);
482 * See sel_updatefunc() for description of the parameters.
483 * \p data should point to a \c t_methoddata_insolidangle.
485 * Calculates which atoms in \p g are within the solid angle spanned by
486 * \c t_methoddata_insolidangle::span and centered at
487 * \c t_methoddata_insolidangle::center, and stores the result in \p out->u.g.
490 evaluate_insolidangle(t_topology *top, t_trxframe *fr, t_pbc *pbc,
491 gmx_ana_pos_t *pos, gmx_ana_selvalue_t *out, void *data)
494 for (int b = 0; b < pos->count(); ++b)
496 if (accept_insolidangle(pos->x[b], pbc, data))
498 gmx_ana_pos_add_to_group(out->u.g, pos, b);
504 * \param[in] sel Selection element to query.
505 * \returns true if the covered fraction can be estimated for \p sel with
506 * _gmx_selelem_estimate_coverfrac(), false otherwise.
509 _gmx_selelem_can_estimate_cover(const gmx::SelectionTreeElement &sel)
511 if (sel.type == SEL_BOOLEAN && sel.u.boolt == BOOL_OR)
516 bool bDynFound = false;
517 gmx::SelectionTreeElementPointer child = sel.child;
520 if (child->type == SEL_EXPRESSION)
522 if (child->u.expr.method->name == sm_insolidangle.name)
524 if (bFound || bDynFound)
530 else if (child->u.expr.method
531 && (child->u.expr.method->flags & SMETH_DYNAMIC))
540 else if (!_gmx_selelem_can_estimate_cover(*child))
550 * \param[in] sel Selection for which the fraction should be calculated.
551 * \returns Fraction of angles covered by the selection (between zero and one).
553 * The return value is undefined if _gmx_selelem_can_estimate_cover() returns
555 * Should be called after gmx_ana_evaluate_selections() has been called for the
559 _gmx_selelem_estimate_coverfrac(const gmx::SelectionTreeElement &sel)
563 if (sel.type == SEL_EXPRESSION && sel.u.expr.method->name == sm_insolidangle.name)
565 t_methoddata_insolidangle *d = (t_methoddata_insolidangle *)sel.u.expr.mdata;
568 d->cfrac = estimate_covered_fraction(d);
572 if (sel.type == SEL_BOOLEAN && sel.u.boolt == BOOL_NOT)
574 cfrac = _gmx_selelem_estimate_coverfrac(*sel.child);
582 /* Here, we assume that the selection is simple enough */
583 gmx::SelectionTreeElementPointer child = sel.child;
586 cfrac = _gmx_selelem_estimate_coverfrac(*child);
597 * \param[in] x1 Unit vector 1.
598 * \param[in] x2 Unit vector 2.
599 * \returns Minus the dot product of \p x1 and \p x2.
601 * This function is used internally to calculate the distance between the
602 * unit vectors \p x1 and \p x2 to find out whether \p x2 is within the
603 * cone centered at \p x1. Currently, the cosine of the angle is used
604 * for efficiency, and the minus is there to make it behave like a normal
605 * distance (larger values mean longer distances).
608 sph_distc(rvec x1, rvec x2)
610 return -iprod(x1, x2);
614 * \param[in] p Partition to search.
615 * \param[in] value Value to search for.
616 * \returns The partition index in \p p that contains \p value.
618 * If \p value is outside the range of \p p, the first/last index is returned.
619 * Otherwise, the return value \c i satisfies \c p->p[i].left<=value and
620 * \c p->p[i+1].left>value
623 find_partition_bin(t_partition *p, real value)
625 int pmin, pmax, pbin;
627 /* Binary search the partition */
628 pmin = 0; pmax = p->n;
629 while (pmax > pmin + 1)
631 pbin = pmin + (pmax - pmin) / 2;
632 if (p->p[pbin].left <= value)
646 * \param[in] surf Surface data structure to search.
647 * \param[in] x Unit vector to find.
648 * \returns The bin index that contains \p x.
650 * The return value is an index to the \p surf->bin array.
653 find_surface_bin(t_methoddata_insolidangle *surf, rvec x)
659 phi = atan2(x[YY], x[XX]);
660 tbin = static_cast<int>(floor(theta / surf->tbinsize));
661 if (tbin >= surf->ntbins)
663 tbin = surf->ntbins - 1;
665 pbin = find_partition_bin(&surf->tbin[tbin], phi);
666 return surf->tbin[tbin].p[pbin].bin;
670 * \param[in,out] surf Surface data structure.
672 * Clears the reference points from the bins and (re)initializes the edges
673 * of the azimuthal bins.
676 clear_surface_points(t_methoddata_insolidangle *surf)
681 for (i = 0; i < surf->ntbins; ++i)
683 c = static_cast<int>(min(sin(surf->tbinsize*i),
684 sin(surf->tbinsize*(i+1)))
685 * M_2PI / surf->targetbinsize) + 1;
691 for (j = 0; j < c; ++j)
693 surf->tbin[i].p[j].left = -M_PI + j*M_2PI/c - 0.0001;
694 surf->tbin[i].p[j].bin = surf->nbins;
695 surf->bin[surf->nbins].n = 0;
698 surf->tbin[i].p[c].left = M_PI + 0.0001;
699 surf->tbin[i].p[c].bin = -1;
704 * \param[in,out] surf Surface data structure.
707 free_surface_points(t_methoddata_insolidangle *surf)
711 for (i = 0; i < surf->nbins; ++i)
715 sfree(surf->bin[i].x);
717 surf->bin[i].n_alloc = 0;
718 surf->bin[i].x = NULL;
723 * \param[in,out] surf Surface data structure.
724 * \param[in] tbin Bin number in the zenith angle direction.
725 * \param[in] pbin Bin number in the azimuthal angle direction.
726 * \param[in] x Point to store.
729 add_surface_point(t_methoddata_insolidangle *surf, int tbin, int pbin, rvec x)
733 bin = surf->tbin[tbin].p[pbin].bin;
734 /* Return if bin is already completely covered */
735 if (surf->bin[bin].n == -1)
739 /* Allocate more space if necessary */
740 if (surf->bin[bin].n == surf->bin[bin].n_alloc)
742 surf->bin[bin].n_alloc += 10;
743 srenew(surf->bin[bin].x, surf->bin[bin].n_alloc);
745 /* Add the point to the bin */
746 copy_rvec(x, surf->bin[bin].x[surf->bin[bin].n]);
751 * \param[in,out] surf Surface data structure.
752 * \param[in] tbin Bin number in the zenith angle direction.
753 * \param[in] pbin Bin number in the azimuthal angle direction.
756 mark_surface_covered(t_methoddata_insolidangle *surf, int tbin, int pbin)
760 bin = surf->tbin[tbin].p[pbin].bin;
761 surf->bin[bin].n = -1;
765 * \param[in,out] surf Surface data structure.
766 * \param[in] tbin Bin number in the zenith angle direction.
767 * \param[in] phi Azimuthal angle of \p x.
768 * \param[in] pdelta1 Width of the cone at the lower edge of \p tbin.
769 * \param[in] pdelta2 Width of the cone at the uppper edge of \p tbin.
770 * \param[in] pdeltamax Max. width of the cone inside \p tbin.
771 * \param[in] x Point to store (should have unit length).
774 update_surface_bin(t_methoddata_insolidangle *surf, int tbin,
775 real phi, real pdelta1, real pdelta2, real pdeltamax,
778 real pdelta, phi1, phi2;
779 int pbin1, pbin2, pbiniter, pbin;
781 /* Find the edges of the bins affected */
782 pdelta = max(max(pdelta1, pdelta2), pdeltamax);
786 pbin = find_partition_bin(&surf->tbin[tbin], phi1);
791 pbin = find_partition_bin(&surf->tbin[tbin], phi1 + M_2PI);
792 pbin1 = pbin - surf->tbin[tbin].n;
797 pbin2 = find_partition_bin(&surf->tbin[tbin], phi2);
801 pbin2 = find_partition_bin(&surf->tbin[tbin], phi2 - M_2PI);
802 pbin2 += surf->tbin[tbin].n;
805 if (pbin2 - pbin1 > surf->tbin[tbin].n)
807 pbin2 = pbin1 + surf->tbin[tbin].n;
809 /* Find the edges of completely covered region */
810 pdelta = min(pdelta1, pdelta2);
817 /* Loop over all affected bins */
818 for (pbiniter = pbin1; pbiniter != pbin2; ++pbiniter, ++pbin)
820 /* Wrap bin around if end reached */
821 if (pbin == surf->tbin[tbin].n)
827 /* Check if bin is completely covered and update */
828 if (surf->tbin[tbin].p[pbin].left >= phi1
829 && surf->tbin[tbin].p[pbin+1].left <= phi2)
831 mark_surface_covered(surf, tbin, pbin);
835 add_surface_point(surf, tbin, pbin, x);
841 * \param[in,out] surf Surface data structure.
842 * \param[in] x Point to store (should have unit length).
844 * Finds all the bins covered by the cone centered at \p x and calls
845 * update_surface_bin() to update them.
848 store_surface_point(t_methoddata_insolidangle *surf, rvec x)
851 real pdeltamax, tmax;
852 real theta1, theta2, pdelta1, pdelta2;
856 phi = atan2(x[YY], x[XX]);
857 /* Find the maximum extent in the phi direction */
858 if (theta <= surf->angcut)
863 else if (theta >= M_PI - surf->angcut)
870 pdeltamax = asin(sin(surf->angcut) / sin(theta));
871 tmax = acos(cos(theta) / cos(surf->angcut));
873 /* Find the first affected bin */
874 tbin = max(static_cast<int>(floor((theta - surf->angcut) / surf->tbinsize)), 0);
875 theta1 = tbin * surf->tbinsize;
876 if (theta1 < theta - surf->angcut)
884 /* Loop through all affected bins */
885 while (tbin < ceil((theta + surf->angcut) / surf->tbinsize)
886 && tbin < surf->ntbins)
888 /* Calculate the next boundaries */
889 theta2 = (tbin+1) * surf->tbinsize;
890 if (theta2 > theta + surf->angcut)
892 /* The circle is completely outside the cone */
895 else if (theta2 <= -(theta - surf->angcut)
896 || theta2 >= M_2PI - (theta + surf->angcut)
897 || tbin == surf->ntbins - 1)
899 /* The circle is completely inside the cone, or we are in the
900 * 360 degree bin covering the pole. */
905 /* TODO: This formula is numerically unstable if theta is very
906 * close to the pole. In practice, it probably does not matter
907 * much, but it would be nicer to adjust the theta bin boundaries
908 * such that the case above catches this instead of falling through
910 pdelta2 = 2*asin(sqrt(
911 (sqr(sin(surf->angcut/2)) - sqr(sin((theta2-theta)/2))) /
912 (sin(theta) * sin(theta2))));
915 if (tmax >= theta1 && tmax <= theta2)
917 update_surface_bin(surf, tbin, phi, pdelta1, pdelta2, pdeltamax, x);
921 update_surface_bin(surf, tbin, phi, pdelta1, pdelta2, 0, x);
931 * \param[in,out] surf Surface data structure.
933 * Currently, this function does nothing.
936 optimize_surface_points(t_methoddata_insolidangle *surf)
938 /* TODO: Implement */
942 * \param[in] surf Surface data structure.
943 * \returns An estimate for the area covered by the reference points.
946 estimate_covered_fraction(t_methoddata_insolidangle *surf)
949 real cfrac, tfrac, pfrac;
952 for (t = 0; t < surf->ntbins; ++t)
954 tfrac = cos(t * surf->tbinsize) - cos((t+1) * surf->tbinsize);
955 for (p = 0; p < surf->tbin[t].n; ++p)
957 pfrac = surf->tbin[t].p[p+1].left - surf->tbin[t].p[p].left;
958 n = surf->bin[surf->tbin[t].p[p].bin].n;
959 if (n == -1) /* Bin completely covered */
961 cfrac += tfrac * pfrac;
963 else if (n > 0) /* Bin partially covered */
965 cfrac += tfrac * pfrac / 2; /* A rough estimate */
969 return cfrac / (4*M_PI);
973 * \param[in] surf Surface data structure to search.
974 * \param[in] x Unit vector to check.
975 * \returns true if \p x is within the solid angle, false otherwise.
978 is_surface_covered(t_methoddata_insolidangle *surf, rvec x)
982 bin = find_surface_bin(surf, x);
983 /* Check for completely covered bin */
984 if (surf->bin[bin].n == -1)
988 /* Check each point that partially covers the bin */
989 for (i = 0; i < surf->bin[bin].n; ++i)
991 if (sph_distc(x, surf->bin[bin].x[i]) < surf->distccut)