2 * This file is part of the GROMACS molecular simulation package.
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.
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.
36 * \page page_module_selection_insolidangle Selection method: insolidangle
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.
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.
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.
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
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
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
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.
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.
102 * Implements the \ref sm_insolidangle "insolidangle" selection method.
105 * The implementation could be optimized quite a bit.
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
112 * \author Teemu Murtola <teemu.murtola@gmail.com>
113 * \ingroup module_selection
121 #include "gromacs/legacyheaders/macros.h"
122 #include "gromacs/math/units.h"
123 #include "gromacs/math/utilities.h"
124 #include "gromacs/math/vec.h"
125 #include "gromacs/pbcutil/pbc.h"
126 #include "gromacs/selection/indexutil.h"
127 #include "gromacs/selection/position.h"
128 #include "gromacs/selection/selection.h"
129 #include "gromacs/utility/exceptions.h"
130 #include "gromacs/utility/smalloc.h"
133 #include "selmethod.h"
140 * Internal data structure for the \p insolidangle selection method.
142 * \see \c t_partition
144 * \ingroup module_selection
148 /** Left edge of the partition. */
150 /** Bin index corresponding to this partition. */
156 * Internal data structure for the \p insolidangle selection method.
158 * Describes the surface partitioning within one slice along the zenith angle.
159 * The slice from azimuthal angle \p p[i].left to \p p[i+1].left belongs to
162 * \ingroup module_selection
166 /** Number of partition items (\p p contains \p n+1 items). */
168 /** Array of partition edges and corresponding bins. */
174 * Internal data structure for the \p insolidangle selection method.
176 * Contains the reference points that partially cover a certain region on the
177 * surface of the unit sphere.
178 * If \p n is -1, the whole region described by the bin is covered.
180 * \ingroup module_selection
184 /** Number of points in the array \p x, -1 if whole bin covered. */
186 /** Number of elements allocated for \p x. */
188 /** Array of points that partially cover the bin. */
190 } t_spheresurfacebin;
194 * Data structure for the \p insolidangle selection method.
196 * All angle values are in the units of radians.
198 * \ingroup module_selection
202 /** Center of the solid angle. */
203 gmx_ana_pos_t center;
204 /** Positions that span the solid angle. */
208 /** Estimate of the covered fraction. */
211 /** Cutoff for the cosine (equals cos(angcut)). */
213 /** Bin size to be used as the target bin size when constructing the bins. */
216 /** Number of bins in the \p tbin array. */
218 /** Size of one bin in the zenith angle direction. */
220 /** Array of zenith angle slices. */
222 /** Number of elements allocated for the \p bin array. */
224 /** Number of elements used in the \p bin array. */
226 /** Array of individual bins. */
227 t_spheresurfacebin *bin;
228 } t_methoddata_insolidangle;
231 * Allocates data for the \p insolidangle selection method.
233 * \param[in] npar Not used (should be 3).
234 * \param[in,out] param Method parameters (should point to
235 * \ref smparams_insolidangle).
236 * \returns Pointer to the allocated data (\ref t_methoddata_insolidangle).
238 * Allocates memory for a \ref t_methoddata_insolidangle structure and
239 * initializes the parameter as follows:
240 * - \p center defines the value for t_methoddata_insolidangle::center.
241 * - \p span defines the value for t_methoddata_insolidangle::span.
242 * - \p cutoff defines the value for t_methoddata_insolidangle::angcut.
245 init_data_insolidangle(int npar, gmx_ana_selparam_t *param);
247 * Initializes the \p insolidangle selection method.
249 * \param top Not used.
250 * \param npar Not used.
251 * \param param Not used.
252 * \param data Pointer to \ref t_methoddata_insolidangle to initialize.
253 * \returns 0 on success, -1 on failure.
255 * Converts t_methoddata_insolidangle::angcut to radians and allocates
256 * and allocates memory for the bins used during the evaluation.
259 init_insolidangle(t_topology * top, int npar, gmx_ana_selparam_t * param, void *data);
260 /** Frees the data allocated for the \p insolidangle selection method. */
262 free_data_insolidangle(void *data);
264 * Initializes the evaluation of the \p insolidangle selection method for a frame.
266 * \param[in] top Not used.
267 * \param[in] fr Not used.
268 * \param[in] pbc PBC structure.
269 * \param data Should point to a \ref t_methoddata_insolidangle.
271 * Creates a lookup structure that enables fast queries of whether a point
272 * is within the solid angle or not.
275 init_frame_insolidangle(t_topology * top, t_trxframe * fr, t_pbc *pbc, void *data);
276 /** Internal helper function for evaluate_insolidangle(). */
278 accept_insolidangle(rvec x, t_pbc *pbc, void *data);
279 /** Evaluates the \p insolidangle selection method. */
281 evaluate_insolidangle(t_topology * /* top */, t_trxframe * /* fr */, t_pbc *pbc,
282 gmx_ana_pos_t *pos, gmx_ana_selvalue_t *out, void *data);
284 /** Calculates the distance between unit vectors. */
286 sph_distc(rvec x1, rvec x2);
287 /** Does a binary search on a \p t_partition to find a bin for a value. */
289 find_partition_bin(t_partition *p, real value);
290 /** Finds a bin that corresponds to a location on the unit sphere surface. */
292 find_surface_bin(t_methoddata_insolidangle *surf, rvec x);
293 /** Clears/initializes the bins on the unit sphere surface. */
295 clear_surface_points(t_methoddata_insolidangle *surf);
296 /** Frees memory allocated for storing the reference points in the surface bins. */
298 free_surface_points(t_methoddata_insolidangle *surf);
299 /** Adds a reference point to a given bin. */
301 add_surface_point(t_methoddata_insolidangle *surf, int tbin, int pbin, rvec x);
302 /** Marks a bin as completely covered. */
304 mark_surface_covered(t_methoddata_insolidangle *surf, int tbin, int pbin);
305 /** Helper function for store_surface_point() to update a single zenith angle bin. */
307 update_surface_bin(t_methoddata_insolidangle *surf, int tbin,
308 real phi, real pdelta1, real pdelta2, real pdeltamax,
310 /** Adds a single reference point and updates the surface bins. */
312 store_surface_point(t_methoddata_insolidangle *surf, rvec x);
314 * Optimizes the surface bins for faster searching.
316 * \param[in,out] surf Surface data structure.
318 * Currently, this function does nothing.
321 optimize_surface_points(t_methoddata_insolidangle *surf);
322 /** Estimates the area covered by the reference cones. */
324 estimate_covered_fraction(t_methoddata_insolidangle *surf);
325 /** Checks whether a point lies within a solid angle. */
327 is_surface_covered(t_methoddata_insolidangle *surf, rvec x);
329 /** Parameters for the \p insolidangle selection method. */
330 static gmx_ana_selparam_t smparams_insolidangle[] = {
331 {"center", {POS_VALUE, 1, {NULL}}, NULL, SPAR_DYNAMIC},
332 {"span", {POS_VALUE, -1, {NULL}}, NULL, SPAR_DYNAMIC | SPAR_VARNUM},
333 {"cutoff", {REAL_VALUE, 1, {NULL}}, NULL, SPAR_OPTIONAL},
336 /** Help text for the \p insolidangle selection method. */
337 static const char *help_insolidangle[] = {
338 "SELECTING ATOMS IN A SOLID ANGLE[PAR]",
340 "[TT]insolidangle center POS span POS_EXPR [cutoff REAL][tt][PAR]",
342 "This keyword selects atoms that are within [TT]REAL[tt] degrees",
343 "(default=5) of any position in [TT]POS_EXPR[tt] as seen from [TT]POS[tt]",
344 "a position expression that evaluates to a single position), i.e., atoms",
345 "in the solid angle spanned by the positions in [TT]POS_EXPR[tt] and",
346 "centered at [TT]POS[tt].[PAR]"
348 "Technically, the solid angle is constructed as a union of small cones",
349 "whose tip is at [TT]POS[tt] and the axis goes through a point in",
350 "[TT]POS_EXPR[tt]. There is such a cone for each position in",
351 "[TT]POS_EXPR[tt], and point is in the solid angle if it lies within any",
352 "of these cones. The cutoff determines the width of the cones.",
355 /** Selection method data for the \p insolidangle method. */
356 gmx_ana_selmethod_t sm_insolidangle = {
357 "insolidangle", GROUP_VALUE, SMETH_DYNAMIC,
358 asize(smparams_insolidangle), smparams_insolidangle,
359 &init_data_insolidangle,
363 &free_data_insolidangle,
364 &init_frame_insolidangle,
366 &evaluate_insolidangle,
367 {"insolidangle center POS span POS_EXPR [cutoff REAL]",
368 asize(help_insolidangle), help_insolidangle},
372 init_data_insolidangle(int /* npar */, gmx_ana_selparam_t *param)
374 t_methoddata_insolidangle *data = new t_methoddata_insolidangle();
378 data->distccut = 0.0;
379 data->targetbinsize = 0.0;
382 data->tbinsize = 0.0;
388 param[0].val.u.p = &data->center;
389 param[1].val.u.p = &data->span;
390 param[2].val.u.r = &data->angcut;
395 init_insolidangle(t_topology * /* top */, int /* npar */, gmx_ana_selparam_t * /* param */, void *data)
397 t_methoddata_insolidangle *surf = (t_methoddata_insolidangle *)data;
400 if (surf->angcut <= 0)
402 GMX_THROW(gmx::InvalidInputError("Angle cutoff should be > 0"));
405 surf->angcut *= DEG2RAD;
407 surf->distccut = -cos(surf->angcut);
408 surf->targetbinsize = surf->angcut / 2;
409 surf->ntbins = static_cast<int>(M_PI / surf->targetbinsize);
410 surf->tbinsize = (180.0 / surf->ntbins)*DEG2RAD;
412 snew(surf->tbin, static_cast<int>(M_PI / surf->tbinsize) + 1);
414 for (i = 0; i < surf->ntbins; ++i)
416 c = static_cast<int>(max(sin(surf->tbinsize*i),
417 sin(surf->tbinsize*(i+1)))
418 * M_2PI / surf->targetbinsize) + 1;
419 snew(surf->tbin[i].p, c+1);
423 snew(surf->bin, surf->maxbins);
427 * \param data Data to free (should point to a \ref t_methoddata_insolidangle).
429 * Frees the memory allocated for \c t_methoddata_insolidangle::center and
430 * \c t_methoddata_insolidangle::span, as well as the memory for the internal
434 free_data_insolidangle(void *data)
436 t_methoddata_insolidangle *d = (t_methoddata_insolidangle *)data;
441 for (i = 0; i < d->ntbins; ++i)
447 free_surface_points(d);
453 init_frame_insolidangle(t_topology * /* top */, t_trxframe * /* fr */, t_pbc *pbc, void *data)
455 t_methoddata_insolidangle *d = (t_methoddata_insolidangle *)data;
459 free_surface_points(d);
460 clear_surface_points(d);
461 for (i = 0; i < d->span.count(); ++i)
465 pbc_dx(pbc, d->span.x[i], d->center.x[0], dx);
469 rvec_sub(d->span.x[i], d->center.x[0], dx);
472 store_surface_point(d, dx);
474 optimize_surface_points(d);
479 * \param[in] x Test point.
480 * \param[in] pbc PBC data (if NULL, no PBC are used).
481 * \param[in] data Pointer to a \c t_methoddata_insolidangle data structure.
482 * \returns true if \p x is within the solid angle, false otherwise.
485 accept_insolidangle(rvec x, t_pbc *pbc, void *data)
487 t_methoddata_insolidangle *d = (t_methoddata_insolidangle *)data;
492 pbc_dx(pbc, x, d->center.x[0], dx);
496 rvec_sub(x, d->center.x[0], dx);
499 return is_surface_covered(d, dx);
503 * See sel_updatefunc() for description of the parameters.
504 * \p data should point to a \c t_methoddata_insolidangle.
506 * Calculates which atoms in \p g are within the solid angle spanned by
507 * \c t_methoddata_insolidangle::span and centered at
508 * \c t_methoddata_insolidangle::center, and stores the result in \p out->u.g.
511 evaluate_insolidangle(t_topology * /* top */, t_trxframe * /* fr */, t_pbc *pbc,
512 gmx_ana_pos_t *pos, gmx_ana_selvalue_t *out, void *data)
515 for (int b = 0; b < pos->count(); ++b)
517 if (accept_insolidangle(pos->x[b], pbc, data))
519 gmx_ana_pos_add_to_group(out->u.g, pos, b);
525 * \param[in] sel Selection element to query.
526 * \returns true if the covered fraction can be estimated for \p sel with
527 * _gmx_selelem_estimate_coverfrac(), false otherwise.
530 _gmx_selelem_can_estimate_cover(const gmx::SelectionTreeElement &sel)
532 if (sel.type == SEL_BOOLEAN && sel.u.boolt == BOOL_OR)
537 bool bDynFound = false;
538 gmx::SelectionTreeElementPointer child = sel.child;
541 if (child->type == SEL_EXPRESSION)
543 if (child->u.expr.method->name == sm_insolidangle.name)
545 if (bFound || bDynFound)
551 else if (child->u.expr.method
552 && (child->u.expr.method->flags & SMETH_DYNAMIC))
561 else if (!_gmx_selelem_can_estimate_cover(*child))
571 * \param[in] sel Selection for which the fraction should be calculated.
572 * \returns Fraction of angles covered by the selection (between zero and one).
574 * The return value is undefined if _gmx_selelem_can_estimate_cover() returns
576 * Should be called after gmx_ana_evaluate_selections() has been called for the
580 _gmx_selelem_estimate_coverfrac(const gmx::SelectionTreeElement &sel)
584 if (sel.type == SEL_EXPRESSION && sel.u.expr.method->name == sm_insolidangle.name)
586 t_methoddata_insolidangle *d = (t_methoddata_insolidangle *)sel.u.expr.mdata;
589 d->cfrac = estimate_covered_fraction(d);
593 if (sel.type == SEL_BOOLEAN && sel.u.boolt == BOOL_NOT)
595 cfrac = _gmx_selelem_estimate_coverfrac(*sel.child);
603 /* Here, we assume that the selection is simple enough */
604 gmx::SelectionTreeElementPointer child = sel.child;
607 cfrac = _gmx_selelem_estimate_coverfrac(*child);
618 * \param[in] x1 Unit vector 1.
619 * \param[in] x2 Unit vector 2.
620 * \returns Minus the dot product of \p x1 and \p x2.
622 * This function is used internally to calculate the distance between the
623 * unit vectors \p x1 and \p x2 to find out whether \p x2 is within the
624 * cone centered at \p x1. Currently, the cosine of the angle is used
625 * for efficiency, and the minus is there to make it behave like a normal
626 * distance (larger values mean longer distances).
629 sph_distc(rvec x1, rvec x2)
631 return -iprod(x1, x2);
635 * \param[in] p Partition to search.
636 * \param[in] value Value to search for.
637 * \returns The partition index in \p p that contains \p value.
639 * If \p value is outside the range of \p p, the first/last index is returned.
640 * Otherwise, the return value \c i satisfies \c p->p[i].left<=value and
641 * \c p->p[i+1].left>value
644 find_partition_bin(t_partition *p, real value)
646 int pmin, pmax, pbin;
648 /* Binary search the partition */
649 pmin = 0; pmax = p->n;
650 while (pmax > pmin + 1)
652 pbin = pmin + (pmax - pmin) / 2;
653 if (p->p[pbin].left <= value)
667 * \param[in] surf Surface data structure to search.
668 * \param[in] x Unit vector to find.
669 * \returns The bin index that contains \p x.
671 * The return value is an index to the \p surf->bin array.
674 find_surface_bin(t_methoddata_insolidangle *surf, rvec x)
680 phi = atan2(x[YY], x[XX]);
681 tbin = static_cast<int>(floor(theta / surf->tbinsize));
682 if (tbin >= surf->ntbins)
684 tbin = surf->ntbins - 1;
686 pbin = find_partition_bin(&surf->tbin[tbin], phi);
687 return surf->tbin[tbin].p[pbin].bin;
691 * \param[in,out] surf Surface data structure.
693 * Clears the reference points from the bins and (re)initializes the edges
694 * of the azimuthal bins.
697 clear_surface_points(t_methoddata_insolidangle *surf)
702 for (i = 0; i < surf->ntbins; ++i)
704 c = static_cast<int>(min(sin(surf->tbinsize*i),
705 sin(surf->tbinsize*(i+1)))
706 * M_2PI / surf->targetbinsize) + 1;
712 for (j = 0; j < c; ++j)
714 surf->tbin[i].p[j].left = -M_PI + j*M_2PI/c - 0.0001;
715 surf->tbin[i].p[j].bin = surf->nbins;
716 surf->bin[surf->nbins].n = 0;
719 surf->tbin[i].p[c].left = M_PI + 0.0001;
720 surf->tbin[i].p[c].bin = -1;
725 * \param[in,out] surf Surface data structure.
728 free_surface_points(t_methoddata_insolidangle *surf)
732 for (i = 0; i < surf->nbins; ++i)
736 sfree(surf->bin[i].x);
738 surf->bin[i].n_alloc = 0;
739 surf->bin[i].x = NULL;
744 * \param[in,out] surf Surface data structure.
745 * \param[in] tbin Bin number in the zenith angle direction.
746 * \param[in] pbin Bin number in the azimuthal angle direction.
747 * \param[in] x Point to store.
750 add_surface_point(t_methoddata_insolidangle *surf, int tbin, int pbin, rvec x)
754 bin = surf->tbin[tbin].p[pbin].bin;
755 /* Return if bin is already completely covered */
756 if (surf->bin[bin].n == -1)
760 /* Allocate more space if necessary */
761 if (surf->bin[bin].n == surf->bin[bin].n_alloc)
763 surf->bin[bin].n_alloc += 10;
764 srenew(surf->bin[bin].x, surf->bin[bin].n_alloc);
766 /* Add the point to the bin */
767 copy_rvec(x, surf->bin[bin].x[surf->bin[bin].n]);
772 * \param[in,out] surf Surface data structure.
773 * \param[in] tbin Bin number in the zenith angle direction.
774 * \param[in] pbin Bin number in the azimuthal angle direction.
777 mark_surface_covered(t_methoddata_insolidangle *surf, int tbin, int pbin)
781 bin = surf->tbin[tbin].p[pbin].bin;
782 surf->bin[bin].n = -1;
786 * \param[in,out] surf Surface data structure.
787 * \param[in] tbin Bin number in the zenith angle direction.
788 * \param[in] phi Azimuthal angle of \p x.
789 * \param[in] pdelta1 Width of the cone at the lower edge of \p tbin.
790 * \param[in] pdelta2 Width of the cone at the uppper edge of \p tbin.
791 * \param[in] pdeltamax Max. width of the cone inside \p tbin.
792 * \param[in] x Point to store (should have unit length).
795 update_surface_bin(t_methoddata_insolidangle *surf, int tbin,
796 real phi, real pdelta1, real pdelta2, real pdeltamax,
799 real pdelta, phi1, phi2;
800 int pbin1, pbin2, pbiniter, pbin;
802 /* Find the edges of the bins affected */
803 pdelta = max(max(pdelta1, pdelta2), pdeltamax);
807 pbin = find_partition_bin(&surf->tbin[tbin], phi1);
812 pbin = find_partition_bin(&surf->tbin[tbin], phi1 + M_2PI);
813 pbin1 = pbin - surf->tbin[tbin].n;
818 pbin2 = find_partition_bin(&surf->tbin[tbin], phi2);
822 pbin2 = find_partition_bin(&surf->tbin[tbin], phi2 - M_2PI);
823 pbin2 += surf->tbin[tbin].n;
826 if (pbin2 - pbin1 > surf->tbin[tbin].n)
828 pbin2 = pbin1 + surf->tbin[tbin].n;
830 /* Find the edges of completely covered region */
831 pdelta = min(pdelta1, pdelta2);
838 /* Loop over all affected bins */
839 for (pbiniter = pbin1; pbiniter != pbin2; ++pbiniter, ++pbin)
841 /* Wrap bin around if end reached */
842 if (pbin == surf->tbin[tbin].n)
848 /* Check if bin is completely covered and update */
849 if (surf->tbin[tbin].p[pbin].left >= phi1
850 && surf->tbin[tbin].p[pbin+1].left <= phi2)
852 mark_surface_covered(surf, tbin, pbin);
856 add_surface_point(surf, tbin, pbin, x);
862 * \param[in,out] surf Surface data structure.
863 * \param[in] x Point to store (should have unit length).
865 * Finds all the bins covered by the cone centered at \p x and calls
866 * update_surface_bin() to update them.
869 store_surface_point(t_methoddata_insolidangle *surf, rvec x)
872 real pdeltamax, tmax;
873 real theta1, theta2, pdelta1, pdelta2;
877 phi = atan2(x[YY], x[XX]);
878 /* Find the maximum extent in the phi direction */
879 if (theta <= surf->angcut)
884 else if (theta >= M_PI - surf->angcut)
891 pdeltamax = asin(sin(surf->angcut) / sin(theta));
892 tmax = acos(cos(theta) / cos(surf->angcut));
894 /* Find the first affected bin */
895 tbin = max(static_cast<int>(floor((theta - surf->angcut) / surf->tbinsize)), 0);
896 theta1 = tbin * surf->tbinsize;
897 if (theta1 < theta - surf->angcut)
905 /* Loop through all affected bins */
906 while (tbin < ceil((theta + surf->angcut) / surf->tbinsize)
907 && tbin < surf->ntbins)
909 /* Calculate the next boundaries */
910 theta2 = (tbin+1) * surf->tbinsize;
911 if (theta2 > theta + surf->angcut)
913 /* The circle is completely outside the cone */
916 else if (theta2 <= -(theta - surf->angcut)
917 || theta2 >= M_2PI - (theta + surf->angcut)
918 || tbin == surf->ntbins - 1)
920 /* The circle is completely inside the cone, or we are in the
921 * 360 degree bin covering the pole. */
926 /* TODO: This formula is numerically unstable if theta is very
927 * close to the pole. In practice, it probably does not matter
928 * much, but it would be nicer to adjust the theta bin boundaries
929 * such that the case above catches this instead of falling through
931 pdelta2 = 2*asin(sqrt(
932 (sqr(sin(surf->angcut/2)) - sqr(sin((theta2-theta)/2))) /
933 (sin(theta) * sin(theta2))));
936 if (tmax >= theta1 && tmax <= theta2)
938 update_surface_bin(surf, tbin, phi, pdelta1, pdelta2, pdeltamax, x);
942 update_surface_bin(surf, tbin, phi, pdelta1, pdelta2, 0, x);
952 optimize_surface_points(t_methoddata_insolidangle * /* surf */)
954 /* TODO: Implement */
958 * \param[in] surf Surface data structure.
959 * \returns An estimate for the area covered by the reference points.
962 estimate_covered_fraction(t_methoddata_insolidangle *surf)
965 real cfrac, tfrac, pfrac;
968 for (t = 0; t < surf->ntbins; ++t)
970 tfrac = cos(t * surf->tbinsize) - cos((t+1) * surf->tbinsize);
971 for (p = 0; p < surf->tbin[t].n; ++p)
973 pfrac = surf->tbin[t].p[p+1].left - surf->tbin[t].p[p].left;
974 n = surf->bin[surf->tbin[t].p[p].bin].n;
975 if (n == -1) /* Bin completely covered */
977 cfrac += tfrac * pfrac;
979 else if (n > 0) /* Bin partially covered */
981 cfrac += tfrac * pfrac / 2; /* A rough estimate */
985 return cfrac / (4*M_PI);
989 * \param[in] surf Surface data structure to search.
990 * \param[in] x Unit vector to check.
991 * \returns true if \p x is within the solid angle, false otherwise.
994 is_surface_covered(t_methoddata_insolidangle *surf, rvec x)
998 bin = find_surface_bin(surf, x);
999 /* Check for completely covered bin */
1000 if (surf->bin[bin].n == -1)
1004 /* Check each point that partially covers the bin */
1005 for (i = 0; i < surf->bin[bin].n; ++i)
1007 if (sph_distc(x, surf->bin[bin].x[i]) < surf->distccut)