3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
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.
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.
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.
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.
29 * For more info, check our website at http://www.gromacs.org
31 /*! \page sm_insolidangle Selection method: insolidangle
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.
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.
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.
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
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
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
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 * Treat zenith angle bins that are completely covered by the cone (in the
87 * case that the cone is centered close to the pole) as a special case.
88 * -# Using the values calculated above, loop through the azimuthal bins that
89 * are partially or completely covered by the cone and update them.
91 * The total solid angle (for covered fraction calculations) is estimated by
92 * taking the total area of completely covered bins plus
93 * half the area of partially covered bins.
94 * The second one is an approximation, but should give reasonable estimates
95 * for the averages as well as in cases where the bin size is small.
98 * \brief Implementation of the \ref sm_insolidangle "insolidangle"
102 * The implementation could be optimized quite a bit.
104 * \todo Move the covered fraction stuff somewhere else and make it more
105 * generic (along the lines it is handled in selection.h and trajana.h).
120 #include <indexutil.h>
121 #include <position.h>
122 #include <selection.h>
123 #include <selmethod.h>
128 * Internal data structure for the \p insolidangle selection method.
130 * \see \c t_partition
134 /** Left edge of the partition. */
136 /** Bin index corresponding to this partition. */
141 * Internal data structure for the \p insolidangle selection method.
143 * Describes the surface partitioning within one slice along the zenith angle.
144 * The slice from azimuthal angle \p p[i].left to \p p[i+1].left belongs to
149 /** Number of partition items (\p p contains \p n+1 items). */
151 /** Array of partition edges and corresponding bins. */
156 * Internal data structure for the \p insolidangle selection method.
158 * Contains the reference points that partially cover a certain region on the
159 * surface of the unit sphere.
160 * If \p n is -1, the whole region described by the bin is covered.
164 /** Number of points in the array \p x, -1 if whole bin covered. */
166 /** Number of elements allocated for \p x. */
168 /** Array of points that partially cover the bin. */
170 } t_spheresurfacebin;
173 * Data structure for the \p insolidangle selection method.
175 * All angle values are in the units of radians.
179 /** Center of the solid angle. */
180 gmx_ana_pos_t center;
181 /** Positions that span the solid angle. */
185 /** Estimate of the covered fraction. */
188 /** Cutoff for the cosine (equals cos(angcut)). */
190 /** Bin size to be used as the target bin size when constructing the bins. */
193 /** Number of bins in the \p tbin array. */
195 /** Size of one bin in the zenith angle direction. */
197 /** Array of zenith angle slices. */
199 /** Number of elements allocated for the \p bin array. */
201 /** Number of elements used in the \p bin array. */
203 /** Array of individual bins. */
204 t_spheresurfacebin *bin;
205 } t_methoddata_insolidangle;
207 /** Allocates data for the \p insolidangle selection method. */
209 init_data_insolidangle(int npar, gmx_ana_selparam_t *param);
210 /** Initializes the \p insolidangle selection method. */
212 init_insolidangle(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
213 /** Frees the data allocated for the \p insolidangle selection method. */
215 free_data_insolidangle(void *data);
216 /** Initializes the evaluation of the \p insolidangle selection method for a frame. */
218 init_frame_insolidangle(t_topology *top, t_trxframe *fr, t_pbc *pbc, void *data);
219 /** Internal helper function for evaluate_insolidangle(). */
221 accept_insolidangle(rvec x, t_pbc *pbc, void *data);
222 /** Evaluates the \p insolidangle selection method. */
224 evaluate_insolidangle(t_topology *top, t_trxframe *fr, t_pbc *pbc,
225 gmx_ana_pos_t *pos, gmx_ana_selvalue_t *out, void *data);
227 /** Calculates the distance between unit vectors. */
229 sph_distc(rvec x1, rvec x2);
230 /** Does a binary search on a \p t_partition to find a bin for a value. */
232 find_partition_bin(t_partition *p, real value);
233 /** Finds a bin that corresponds to a location on the unit sphere surface. */
235 find_surface_bin(t_methoddata_insolidangle *surf, rvec x);
236 /** Clears/initializes the bins on the unit sphere surface. */
238 clear_surface_points(t_methoddata_insolidangle *surf);
239 /** Frees memory allocated for storing the reference points in the surface bins. */
241 free_surface_points(t_methoddata_insolidangle *surf);
242 /** Adds a reference point to a given bin. */
244 add_surface_point(t_methoddata_insolidangle *surf, int tbin, int pbin, rvec x);
245 /** Marks a bin as completely covered. */
247 mark_surface_covered(t_methoddata_insolidangle *surf, int tbin, int pbin);
248 /** Helper function for store_surface_point() to update a single zenith angle bin. */
250 update_surface_bin(t_methoddata_insolidangle *surf, int tbin,
251 real phi, real pdelta1, real pdelta2, real pdeltamax,
253 /** Adds a single reference point and updates the surface bins. */
255 store_surface_point(t_methoddata_insolidangle *surf, rvec x);
256 /** Optimizes the surface bins for faster searching. */
258 optimize_surface_points(t_methoddata_insolidangle *surf);
259 /** Estimates the area covered by the reference cones. */
261 estimate_covered_fraction(t_methoddata_insolidangle *surf);
262 /** Checks whether a point lies within a solid angle. */
264 is_surface_covered(t_methoddata_insolidangle *surf, rvec x);
266 /** Parameters for the \p insolidangle selection method. */
267 static gmx_ana_selparam_t smparams_insolidangle[] = {
268 {"center", {POS_VALUE, 1, {NULL}}, NULL, SPAR_DYNAMIC},
269 {"span", {POS_VALUE, -1, {NULL}}, NULL, SPAR_DYNAMIC | SPAR_VARNUM},
270 {"cutoff", {REAL_VALUE, 1, {NULL}}, NULL, SPAR_OPTIONAL},
273 /** Help text for the \p insolidangle selection method. */
274 static const char *help_insolidangle[] = {
275 "SELECTING ATOMS IN A SOLID ANGLE[PAR]",
277 "[TT]insolidangle center POS span POS_EXPR [cutoff REAL][tt][PAR]",
279 "This keyword selects atoms that are within [TT]REAL[tt] degrees",
280 "(default=5) of any position in [TT]POS_EXPR[tt] as seen from [TT]POS[tt]",
281 "a position expression that evaluates to a single position), i.e., atoms",
282 "in the solid angle spanned by the positions in [TT]POS_EXPR[tt] and",
283 "centered at [TT]POS[tt].[PAR]"
285 "Technically, the solid angle is constructed as a union of small cones",
286 "whose tip is at [TT]POS[tt] and the axis goes through a point in",
287 "[TT]POS_EXPR[tt]. There is such a cone for each position in",
288 "[TT]POS_EXPR[tt], and point is in the solid angle if it lies within any",
289 "of these cones. The cutoff determines the width of the cones.",
292 /** \internal Selection method data for the \p insolidangle method. */
293 gmx_ana_selmethod_t sm_insolidangle = {
294 "insolidangle", GROUP_VALUE, SMETH_DYNAMIC,
295 asize(smparams_insolidangle), smparams_insolidangle,
296 &init_data_insolidangle,
300 &free_data_insolidangle,
301 &init_frame_insolidangle,
303 &evaluate_insolidangle,
304 {"insolidangle center POS span POS_EXPR [cutoff REAL]",
305 asize(help_insolidangle), help_insolidangle},
309 * \param[in] npar Not used (should be 3).
310 * \param[in,out] param Method parameters (should point to
311 * \ref smparams_insolidangle).
312 * \returns Pointer to the allocated data (\ref t_methoddata_insolidangle).
314 * Allocates memory for a \ref t_methoddata_insolidangle structure and
315 * initializes the parameter as follows:
316 * - \p center defines the value for t_methoddata_insolidangle::center.
317 * - \p span defines the value for t_methoddata_insolidangle::span.
318 * - \p cutoff defines the value for t_methoddata_insolidangle::angcut.
321 init_data_insolidangle(int npar, gmx_ana_selparam_t *param)
323 t_methoddata_insolidangle *data;
327 param[0].val.u.p = &data->center;
328 param[1].val.u.p = &data->span;
329 param[2].val.u.r = &data->angcut;
334 * \param top Not used.
335 * \param npar Not used.
336 * \param param Not used.
337 * \param data Pointer to \ref t_methoddata_insolidangle to initialize.
338 * \returns 0 on success, -1 on failure.
340 * Converts t_methoddata_insolidangle::angcut to radians and allocates
341 * and allocates memory for the bins used during the evaluation.
344 init_insolidangle(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data)
346 t_methoddata_insolidangle *surf = (t_methoddata_insolidangle *)data;
349 if (surf->angcut <= 0)
351 fprintf(stderr, "error: angle cutoff should be > 0");
355 surf->angcut *= DEG2RAD;
357 surf->distccut = -cos(surf->angcut);
358 surf->targetbinsize = surf->angcut / 2;
359 surf->ntbins = (int) (M_PI / surf->targetbinsize);
360 surf->tbinsize = (180.0 / surf->ntbins)*DEG2RAD;
362 snew(surf->tbin, (int)(M_PI/surf->tbinsize) + 1);
364 for (i = 0; i < surf->ntbins; ++i)
366 c = max(sin(surf->tbinsize*i), sin(surf->tbinsize*(i+1)))
367 * M_2PI / surf->targetbinsize + 1;
368 snew(surf->tbin[i].p, c+1);
372 snew(surf->bin, surf->maxbins);
378 * \param data Data to free (should point to a \ref t_methoddata_insolidangle).
380 * Frees the memory allocated for \c t_methoddata_insolidangle::center and
381 * \c t_methoddata_insolidangle::span, as well as the memory for the internal
385 free_data_insolidangle(void *data)
387 t_methoddata_insolidangle *d = (t_methoddata_insolidangle *)data;
392 for (i = 0; i < d->ntbins; ++i)
398 free_surface_points(d);
403 * \param[in] top Not used.
404 * \param[in] fr Current frame.
405 * \param[in] pbc PBC structure.
406 * \param data Should point to a \ref t_methoddata_insolidangle.
407 * \returns 0 on success, a non-zero error code on error.
409 * Creates a lookup structure that enables fast queries of whether a point
410 * is within the solid angle or not.
413 init_frame_insolidangle(t_topology *top, t_trxframe *fr, t_pbc *pbc, void *data)
415 t_methoddata_insolidangle *d = (t_methoddata_insolidangle *)data;
419 free_surface_points(d);
420 clear_surface_points(d);
421 for (i = 0; i < d->span.nr; ++i)
425 pbc_dx(pbc, d->span.x[i], d->center.x[0], dx);
429 rvec_sub(d->span.x[i], d->center.x[0], dx);
432 store_surface_point(d, dx);
434 optimize_surface_points(d);
440 * \param[in] x Test point.
441 * \param[in] pbc PBC data (if NULL, no PBC are used).
442 * \param[in] data Pointer to a \c t_methoddata_insolidangle data structure.
443 * \returns TRUE if \p x is within the solid angle, FALSE otherwise.
446 accept_insolidangle(rvec x, t_pbc *pbc, void *data)
448 t_methoddata_insolidangle *d = (t_methoddata_insolidangle *)data;
453 pbc_dx(pbc, x, d->center.x[0], dx);
457 rvec_sub(x, d->center.x[0], dx);
460 return is_surface_covered(d, dx);
464 * See sel_updatefunc() for description of the parameters.
465 * \p data should point to a \c t_methoddata_insolidangle.
467 * Calculates which atoms in \p g are within the solid angle spanned by
468 * \c t_methoddata_insolidangle::span and centered at
469 * \c t_methoddata_insolidangle::center, and stores the result in \p out->u.g.
472 evaluate_insolidangle(t_topology *top, t_trxframe *fr, t_pbc *pbc,
473 gmx_ana_pos_t *pos, gmx_ana_selvalue_t *out, void *data)
475 t_methoddata_insolidangle *d = (t_methoddata_insolidangle *)data;
479 for (b = 0; b < pos->nr; ++b)
481 if (accept_insolidangle(pos->x[b], pbc, data))
483 gmx_ana_pos_append(NULL, out->u.g, pos, b, 0);
490 * \param[in] sel Selection element to query.
491 * \returns TRUE if the covered fraction can be estimated for \p sel with
492 * _gmx_selelem_estimate_coverfrac(), FALSE otherwise.
495 _gmx_selelem_can_estimate_cover(t_selelem *sel)
501 if (sel->type == SEL_BOOLEAN && sel->u.boolt == BOOL_OR)
510 if (child->type == SEL_EXPRESSION)
512 if (child->u.expr.method->name == sm_insolidangle.name)
514 if (bFound || bDynFound)
520 else if (child->u.expr.method
521 && (child->u.expr.method->flags & SMETH_DYNAMIC))
530 else if (!_gmx_selelem_can_estimate_cover(child))
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).
543 * The return value is undefined if _gmx_selelem_can_estimate_cover() returns
545 * Should be called after gmx_ana_evaluate_selections() has been called for the
549 _gmx_selelem_estimate_coverfrac(t_selelem *sel)
554 if (sel->type == SEL_EXPRESSION && sel->u.expr.method->name == sm_insolidangle.name)
556 t_methoddata_insolidangle *d = (t_methoddata_insolidangle *)sel->u.expr.mdata;
559 d->cfrac = estimate_covered_fraction(d);
563 if (sel->type == SEL_BOOLEAN && sel->u.boolt == BOOL_NOT)
565 cfrac = _gmx_selelem_estimate_coverfrac(sel->child);
573 /* Here, we assume that the selection is simple enough */
577 cfrac = _gmx_selelem_estimate_coverfrac(child);
588 * \param[in] x1 Unit vector 1.
589 * \param[in] x2 Unit vector 2.
590 * \returns Minus the dot product of \p x1 and \p x2.
592 * This function is used internally to calculate the distance between the
593 * unit vectors \p x1 and \p x2 to find out whether \p x2 is within the
594 * cone centered at \p x1. Currently, the cosine of the angle is used
595 * for efficiency, and the minus is there to make it behave like a normal
596 * distance (larger values mean longer distances).
599 sph_distc(rvec x1, rvec x2)
601 return -iprod(x1, x2);
605 * \param[in] p Partition to search.
606 * \param[in] value Value to search for.
607 * \returns The partition index in \p p that contains \p value.
609 * If \p value is outside the range of \p p, the first/last index is returned.
610 * Otherwise, the return value \c i satisfies \c p->p[i].left<=value and
611 * \c p->p[i+1].left>value
614 find_partition_bin(t_partition *p, real value)
616 int pmin, pmax, pbin;
618 /* Binary search the partition */
619 pmin = 0; pmax = p->n;
620 while (pmax > pmin + 1)
622 pbin = pmin + (pmax - pmin) / 2;
623 if (p->p[pbin].left <= value)
637 * \param[in] surf Surface data structure to search.
638 * \param[in] x Unit vector to find.
639 * \returns The bin index that contains \p x.
641 * The return value is an index to the \p surf->bin array.
644 find_surface_bin(t_methoddata_insolidangle *surf, rvec x)
650 phi = atan2(x[YY], x[XX]);
651 tbin = floor(theta / surf->tbinsize);
652 if (tbin >= surf->ntbins)
654 tbin = surf->ntbins - 1;
656 pbin = find_partition_bin(&surf->tbin[tbin], phi);
657 return surf->tbin[tbin].p[pbin].bin;
661 * \param[in,out] surf Surface data structure.
663 * Clears the reference points from the bins and (re)initializes the edges
664 * of the azimuthal bins.
667 clear_surface_points(t_methoddata_insolidangle *surf)
672 for (i = 0; i < surf->ntbins; ++i)
674 c = min(sin(surf->tbinsize*i), sin(surf->tbinsize*(i+1)))
675 * M_2PI / surf->targetbinsize + 1;
681 for (j = 0; j < c; ++j)
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;
688 surf->tbin[i].p[c].left = M_PI + 0.0001;
689 surf->tbin[i].p[c].bin = -1;
694 * \param[in,out] surf Surface data structure.
697 free_surface_points(t_methoddata_insolidangle *surf)
701 for (i = 0; i < surf->nbins; ++i)
705 sfree(surf->bin[i].x);
707 surf->bin[i].n_alloc = 0;
708 surf->bin[i].x = NULL;
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.
719 add_surface_point(t_methoddata_insolidangle *surf, int tbin, int pbin, rvec x)
723 bin = surf->tbin[tbin].p[pbin].bin;
724 /* Return if bin is already completely covered */
725 if (surf->bin[bin].n == -1)
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);
732 /* Add the point to the bin */
733 copy_rvec(x, surf->bin[bin].x[surf->bin[bin].n]);
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.
743 mark_surface_covered(t_methoddata_insolidangle *surf, int tbin, int pbin)
747 bin = surf->tbin[tbin].p[pbin].bin;
748 surf->bin[bin].n = -1;
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).
761 update_surface_bin(t_methoddata_insolidangle *surf, int tbin,
762 real phi, real pdelta1, real pdelta2, real pdeltamax,
765 real pdelta, phi1, phi2;
766 int pbin1, pbin2, pbiniter, pbin;
768 /* Find the edges of the bins affected */
769 pdelta = max(max(pdelta1, pdelta2), pdeltamax);
773 pbin = find_partition_bin(&surf->tbin[tbin], phi1);
778 pbin = find_partition_bin(&surf->tbin[tbin], phi1 + M_2PI);
779 pbin1 = pbin - surf->tbin[tbin].n;
784 pbin2 = find_partition_bin(&surf->tbin[tbin], phi2);
788 pbin2 = find_partition_bin(&surf->tbin[tbin], phi2 - M_2PI);
789 pbin2 += surf->tbin[tbin].n;
792 if (pbin2 - pbin1 > surf->tbin[tbin].n)
794 pbin2 = pbin1 + surf->tbin[tbin].n;
796 /* Find the edges of completely covered region */
797 pdelta = min(pdelta1, pdelta2);
804 /* Loop over all affected bins */
805 for (pbiniter = pbin1; pbiniter != pbin2; ++pbiniter, ++pbin)
807 /* Wrap bin around if end reached */
808 if (pbin == surf->tbin[tbin].n)
814 /* Check if bin is completely covered and update */
815 if (surf->tbin[tbin].p[pbin].left >= phi1
816 && surf->tbin[tbin].p[pbin+1].left <= phi2)
818 mark_surface_covered(surf, tbin, pbin);
822 add_surface_point(surf, tbin, pbin, x);
828 * \param[in,out] surf Surface data structure.
829 * \param[in] x Point to store (should have unit length).
831 * Finds all the bins covered by the cone centered at \p x and calls
832 * update_surface_bin() to update them.
835 store_surface_point(t_methoddata_insolidangle *surf, rvec x)
838 real pdeltamax, tmax;
839 real theta1, theta2, pdelta1, pdelta2;
843 phi = atan2(x[YY], x[XX]);
844 /* Find the maximum extent in the phi direction */
845 if (theta <= surf->angcut)
850 else if (theta >= M_PI - surf->angcut)
857 pdeltamax = asin(sin(surf->angcut) / sin(theta));
858 tmax = acos(cos(theta) / cos(surf->angcut));
860 /* Find the first affected bin */
861 tbin = max(floor((theta - surf->angcut) / surf->tbinsize), 0);
862 theta1 = tbin * surf->tbinsize;
863 if (theta1 < theta - surf->angcut)
871 /* Loop through all affected bins */
872 while (tbin < ceil((theta + surf->angcut) / surf->tbinsize)
873 && tbin < surf->ntbins)
875 /* Calculate the next boundaries */
876 theta2 = (tbin+1) * surf->tbinsize;
877 if (theta2 > theta + surf->angcut)
879 /* The circle is completely outside the cone */
882 else if (theta2 <= -(theta - surf->angcut)
883 || theta2 >= M_2PI - (theta + surf->angcut)
884 || tbin == surf->ntbins - 1)
886 /* The circle is completely inside the cone, or we are in the
887 * 360 degree bin covering the pole. */
892 /* TODO: This formula is numerically unstable if theta is very
893 * close to the pole. In practice, it probably does not matter
894 * much, but it would be nicer to adjust the theta bin boundaries
895 * such that the case above catches this instead of falling through
897 pdelta2 = 2*asin(sqrt(
898 (sqr(sin(surf->angcut/2)) - sqr(sin((theta2-theta)/2))) /
899 (sin(theta) * sin(theta2))));
902 if (tmax >= theta1 && tmax <= theta2)
904 update_surface_bin(surf, tbin, phi, pdelta1, pdelta2, pdeltamax, x);
908 update_surface_bin(surf, tbin, phi, pdelta1, pdelta2, 0, x);
918 * \param[in,out] surf Surface data structure.
920 * Currently, this function does nothing.
923 optimize_surface_points(t_methoddata_insolidangle *surf)
925 /* TODO: Implement */
929 * \param[in] surf Surface data structure.
930 * \returns An estimate for the area covered by the reference points.
933 estimate_covered_fraction(t_methoddata_insolidangle *surf)
936 real cfrac, tfrac, pfrac;
939 for (t = 0; t < surf->ntbins; ++t)
941 tfrac = cos(t * surf->tbinsize) - cos((t+1) * surf->tbinsize);
942 for (p = 0; p < surf->tbin[t].n; ++p)
944 pfrac = surf->tbin[t].p[p+1].left - surf->tbin[t].p[p].left;
945 n = surf->bin[surf->tbin[t].p[p].bin].n;
946 if (n == -1) /* Bin completely covered */
948 cfrac += tfrac * pfrac;
950 else if (n > 0) /* Bin partially covered */
952 cfrac += tfrac * pfrac / 2; /* A rough estimate */
956 return cfrac / (4*M_PI);
960 * \param[in] surf Surface data structure to search.
961 * \param[in] x Unit vector to check.
962 * \returns TRUE if \p x is within the solid angle, FALSE otherwise.
965 is_surface_covered(t_methoddata_insolidangle *surf, rvec x)
969 bin = find_surface_bin(surf, x);
970 /* Check for completely covered bin */
971 if (surf->bin[bin].n == -1)
975 /* Check each point that partially covers the bin */
976 for (i = 0; i < surf->bin[bin].n; ++i)
978 if (sph_distc(x, surf->bin[bin].x[i]) < surf->distccut)