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 page_module_selection_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 * -# Using the values calculated above, loop through the azimuthal bins that
87 * are partially or completely covered by the cone and update them.
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.
97 * Implements the \ref sm_insolidangle "insolidangle" selection method.
100 * The implementation could be optimized quite a bit.
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
107 * \author Teemu Murtola <teemu.murtola@cbr.su.se>
108 * \ingroup module_selection
123 // FIXME: Should really be in the beginning, but causes compilation errors
126 #include "gromacs/fatalerror/exceptions.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"
138 * Internal data structure for the \p insolidangle selection method.
140 * \see \c t_partition
144 /** Left edge of the partition. */
146 /** Bin index corresponding to this partition. */
151 * Internal data structure for the \p insolidangle selection method.
153 * Describes the surface partitioning within one slice along the zenith angle.
154 * The slice from azimuthal angle \p p[i].left to \p p[i+1].left belongs to
159 /** Number of partition items (\p p contains \p n+1 items). */
161 /** Array of partition edges and corresponding bins. */
166 * Internal data structure for the \p insolidangle selection method.
168 * Contains the reference points that partially cover a certain region on the
169 * surface of the unit sphere.
170 * If \p n is -1, the whole region described by the bin is covered.
174 /** Number of points in the array \p x, -1 if whole bin covered. */
176 /** Number of elements allocated for \p x. */
178 /** Array of points that partially cover the bin. */
180 } t_spheresurfacebin;
183 * Data structure for the \p insolidangle selection method.
185 * All angle values are in the units of radians.
189 /** Center of the solid angle. */
190 gmx_ana_pos_t center;
191 /** Positions that span the solid angle. */
195 /** Estimate of the covered fraction. */
198 /** Cutoff for the cosine (equals cos(angcut)). */
200 /** Bin size to be used as the target bin size when constructing the bins. */
203 /** Number of bins in the \p tbin array. */
205 /** Size of one bin in the zenith angle direction. */
207 /** Array of zenith angle slices. */
209 /** Number of elements allocated for the \p bin array. */
211 /** Number of elements used in the \p bin array. */
213 /** Array of individual bins. */
214 t_spheresurfacebin *bin;
215 } t_methoddata_insolidangle;
217 /** Allocates data for the \p insolidangle selection method. */
219 init_data_insolidangle(int npar, gmx_ana_selparam_t *param);
220 /** Initializes the \p insolidangle selection method. */
222 init_insolidangle(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data);
223 /** Frees the data allocated for the \p insolidangle selection method. */
225 free_data_insolidangle(void *data);
226 /** Initializes the evaluation of the \p insolidangle selection method for a frame. */
228 init_frame_insolidangle(t_topology *top, t_trxframe *fr, t_pbc *pbc, void *data);
229 /** Internal helper function for evaluate_insolidangle(). */
231 accept_insolidangle(rvec x, t_pbc *pbc, void *data);
232 /** Evaluates the \p insolidangle selection method. */
234 evaluate_insolidangle(t_topology *top, t_trxframe *fr, t_pbc *pbc,
235 gmx_ana_pos_t *pos, gmx_ana_selvalue_t *out, void *data);
237 /** Calculates the distance between unit vectors. */
239 sph_distc(rvec x1, rvec x2);
240 /** Does a binary search on a \p t_partition to find a bin for a value. */
242 find_partition_bin(t_partition *p, real value);
243 /** Finds a bin that corresponds to a location on the unit sphere surface. */
245 find_surface_bin(t_methoddata_insolidangle *surf, rvec x);
246 /** Clears/initializes the bins on the unit sphere surface. */
248 clear_surface_points(t_methoddata_insolidangle *surf);
249 /** Frees memory allocated for storing the reference points in the surface bins. */
251 free_surface_points(t_methoddata_insolidangle *surf);
252 /** Adds a reference point to a given bin. */
254 add_surface_point(t_methoddata_insolidangle *surf, int tbin, int pbin, rvec x);
255 /** Marks a bin as completely covered. */
257 mark_surface_covered(t_methoddata_insolidangle *surf, int tbin, int pbin);
258 /** Helper function for store_surface_point() to update a single zenith angle bin. */
260 update_surface_bin(t_methoddata_insolidangle *surf, int tbin,
261 real phi, real pdelta1, real pdelta2, real pdeltamax,
263 /** Adds a single reference point and updates the surface bins. */
265 store_surface_point(t_methoddata_insolidangle *surf, rvec x);
266 /** Optimizes the surface bins for faster searching. */
268 optimize_surface_points(t_methoddata_insolidangle *surf);
269 /** Estimates the area covered by the reference cones. */
271 estimate_covered_fraction(t_methoddata_insolidangle *surf);
272 /** Checks whether a point lies within a solid angle. */
274 is_surface_covered(t_methoddata_insolidangle *surf, rvec x);
276 /** Parameters for the \p insolidangle selection method. */
277 static gmx_ana_selparam_t smparams_insolidangle[] = {
278 {"center", {POS_VALUE, 1, {NULL}}, NULL, SPAR_DYNAMIC},
279 {"span", {POS_VALUE, -1, {NULL}}, NULL, SPAR_DYNAMIC | SPAR_VARNUM},
280 {"cutoff", {REAL_VALUE, 1, {NULL}}, NULL, SPAR_OPTIONAL},
283 /** Help text for the \p insolidangle selection method. */
284 static const char *help_insolidangle[] = {
285 "SELECTING ATOMS IN A SOLID ANGLE[PAR]",
287 "[TT]insolidangle center POS span POS_EXPR [cutoff REAL][tt][PAR]",
289 "This keyword selects atoms that are within [TT]REAL[tt] degrees",
290 "(default=5) of any position in [TT]POS_EXPR[tt] as seen from [TT]POS[tt]",
291 "a position expression that evaluates to a single position), i.e., atoms",
292 "in the solid angle spanned by the positions in [TT]POS_EXPR[tt] and",
293 "centered at [TT]POS[tt].[PAR]"
295 "Technically, the solid angle is constructed as a union of small cones",
296 "whose tip is at [TT]POS[tt] and the axis goes through a point in",
297 "[TT]POS_EXPR[tt]. There is such a cone for each position in",
298 "[TT]POS_EXPR[tt], and point is in the solid angle if it lies within any",
299 "of these cones. The cutoff determines the width of the cones.",
302 /** \internal Selection method data for the \p insolidangle method. */
303 gmx_ana_selmethod_t sm_insolidangle = {
304 "insolidangle", GROUP_VALUE, SMETH_DYNAMIC,
305 asize(smparams_insolidangle), smparams_insolidangle,
306 &init_data_insolidangle,
310 &free_data_insolidangle,
311 &init_frame_insolidangle,
313 &evaluate_insolidangle,
314 {"insolidangle center POS span POS_EXPR [cutoff REAL]",
315 asize(help_insolidangle), help_insolidangle},
319 * \param[in] npar Not used (should be 3).
320 * \param[in,out] param Method parameters (should point to
321 * \ref smparams_insolidangle).
322 * \returns Pointer to the allocated data (\ref t_methoddata_insolidangle).
324 * Allocates memory for a \ref t_methoddata_insolidangle structure and
325 * initializes the parameter as follows:
326 * - \p center defines the value for t_methoddata_insolidangle::center.
327 * - \p span defines the value for t_methoddata_insolidangle::span.
328 * - \p cutoff defines the value for t_methoddata_insolidangle::angcut.
331 init_data_insolidangle(int npar, gmx_ana_selparam_t *param)
333 t_methoddata_insolidangle *data;
337 param[0].val.u.p = &data->center;
338 param[1].val.u.p = &data->span;
339 param[2].val.u.r = &data->angcut;
344 * \param top Not used.
345 * \param npar Not used.
346 * \param param Not used.
347 * \param data Pointer to \ref t_methoddata_insolidangle to initialize.
348 * \returns 0 on success, -1 on failure.
350 * Converts t_methoddata_insolidangle::angcut to radians and allocates
351 * and allocates memory for the bins used during the evaluation.
354 init_insolidangle(t_topology *top, int npar, gmx_ana_selparam_t *param, void *data)
356 t_methoddata_insolidangle *surf = (t_methoddata_insolidangle *)data;
359 if (surf->angcut <= 0)
361 GMX_THROW(gmx::InvalidInputError("Angle cutoff should be > 0"));
364 surf->angcut *= DEG2RAD;
366 surf->distccut = -cos(surf->angcut);
367 surf->targetbinsize = surf->angcut / 2;
368 surf->ntbins = (int) (M_PI / surf->targetbinsize);
369 surf->tbinsize = (180.0 / surf->ntbins)*DEG2RAD;
371 snew(surf->tbin, (int)(M_PI/surf->tbinsize) + 1);
373 for (i = 0; i < surf->ntbins; ++i)
375 c = max(sin(surf->tbinsize*i), sin(surf->tbinsize*(i+1)))
376 * M_2PI / surf->targetbinsize + 1;
377 snew(surf->tbin[i].p, c+1);
381 snew(surf->bin, surf->maxbins);
385 * \param data Data to free (should point to a \ref t_methoddata_insolidangle).
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
392 free_data_insolidangle(void *data)
394 t_methoddata_insolidangle *d = (t_methoddata_insolidangle *)data;
399 for (i = 0; i < d->ntbins; ++i)
405 free_surface_points(d);
410 * \param[in] top Not used.
411 * \param[in] fr Current frame.
412 * \param[in] pbc PBC structure.
413 * \param data Should point to a \ref t_methoddata_insolidangle.
415 * Creates a lookup structure that enables fast queries of whether a point
416 * is within the solid angle or not.
419 init_frame_insolidangle(t_topology *top, t_trxframe *fr, t_pbc *pbc, void *data)
421 t_methoddata_insolidangle *d = (t_methoddata_insolidangle *)data;
425 free_surface_points(d);
426 clear_surface_points(d);
427 for (i = 0; i < d->span.nr; ++i)
431 pbc_dx(pbc, d->span.x[i], d->center.x[0], dx);
435 rvec_sub(d->span.x[i], d->center.x[0], dx);
438 store_surface_point(d, dx);
440 optimize_surface_points(d);
445 * \param[in] x Test point.
446 * \param[in] pbc PBC data (if NULL, no PBC are used).
447 * \param[in] data Pointer to a \c t_methoddata_insolidangle data structure.
448 * \returns TRUE if \p x is within the solid angle, FALSE otherwise.
451 accept_insolidangle(rvec x, t_pbc *pbc, void *data)
453 t_methoddata_insolidangle *d = (t_methoddata_insolidangle *)data;
458 pbc_dx(pbc, x, d->center.x[0], dx);
462 rvec_sub(x, d->center.x[0], dx);
465 return is_surface_covered(d, dx);
469 * See sel_updatefunc() for description of the parameters.
470 * \p data should point to a \c t_methoddata_insolidangle.
472 * Calculates which atoms in \p g are within the solid angle spanned by
473 * \c t_methoddata_insolidangle::span and centered at
474 * \c t_methoddata_insolidangle::center, and stores the result in \p out->u.g.
477 evaluate_insolidangle(t_topology *top, t_trxframe *fr, t_pbc *pbc,
478 gmx_ana_pos_t *pos, gmx_ana_selvalue_t *out, void *data)
480 t_methoddata_insolidangle *d = (t_methoddata_insolidangle *)data;
484 for (b = 0; b < pos->nr; ++b)
486 if (accept_insolidangle(pos->x[b], pbc, data))
488 gmx_ana_pos_append(NULL, out->u.g, pos, b, 0);
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.
499 _gmx_selelem_can_estimate_cover(t_selelem *sel)
505 if (sel->type == SEL_BOOLEAN && sel->u.boolt == BOOL_OR)
514 if (child->type == SEL_EXPRESSION)
516 if (child->u.expr.method->name == sm_insolidangle.name)
518 if (bFound || bDynFound)
524 else if (child->u.expr.method
525 && (child->u.expr.method->flags & SMETH_DYNAMIC))
534 else if (!_gmx_selelem_can_estimate_cover(child))
544 * \param[in] sel Selection for which the fraction should be calculated.
545 * \returns Fraction of angles covered by the selection (between zero and one).
547 * The return value is undefined if _gmx_selelem_can_estimate_cover() returns
549 * Should be called after gmx_ana_evaluate_selections() has been called for the
553 _gmx_selelem_estimate_coverfrac(t_selelem *sel)
558 if (sel->type == SEL_EXPRESSION && sel->u.expr.method->name == sm_insolidangle.name)
560 t_methoddata_insolidangle *d = (t_methoddata_insolidangle *)sel->u.expr.mdata;
563 d->cfrac = estimate_covered_fraction(d);
567 if (sel->type == SEL_BOOLEAN && sel->u.boolt == BOOL_NOT)
569 cfrac = _gmx_selelem_estimate_coverfrac(sel->child);
577 /* Here, we assume that the selection is simple enough */
581 cfrac = _gmx_selelem_estimate_coverfrac(child);
592 * \param[in] x1 Unit vector 1.
593 * \param[in] x2 Unit vector 2.
594 * \returns Minus the dot product of \p x1 and \p x2.
596 * This function is used internally to calculate the distance between the
597 * unit vectors \p x1 and \p x2 to find out whether \p x2 is within the
598 * cone centered at \p x1. Currently, the cosine of the angle is used
599 * for efficiency, and the minus is there to make it behave like a normal
600 * distance (larger values mean longer distances).
603 sph_distc(rvec x1, rvec x2)
605 return -iprod(x1, x2);
609 * \param[in] p Partition to search.
610 * \param[in] value Value to search for.
611 * \returns The partition index in \p p that contains \p value.
613 * If \p value is outside the range of \p p, the first/last index is returned.
614 * Otherwise, the return value \c i satisfies \c p->p[i].left<=value and
615 * \c p->p[i+1].left>value
618 find_partition_bin(t_partition *p, real value)
620 int pmin, pmax, pbin;
622 /* Binary search the partition */
623 pmin = 0; pmax = p->n;
624 while (pmax > pmin + 1)
626 pbin = pmin + (pmax - pmin) / 2;
627 if (p->p[pbin].left <= value)
641 * \param[in] surf Surface data structure to search.
642 * \param[in] x Unit vector to find.
643 * \returns The bin index that contains \p x.
645 * The return value is an index to the \p surf->bin array.
648 find_surface_bin(t_methoddata_insolidangle *surf, rvec x)
654 phi = atan2(x[YY], x[XX]);
655 tbin = floor(theta / surf->tbinsize);
656 if (tbin >= surf->ntbins)
658 tbin = surf->ntbins - 1;
660 pbin = find_partition_bin(&surf->tbin[tbin], phi);
661 return surf->tbin[tbin].p[pbin].bin;
665 * \param[in,out] surf Surface data structure.
667 * Clears the reference points from the bins and (re)initializes the edges
668 * of the azimuthal bins.
671 clear_surface_points(t_methoddata_insolidangle *surf)
676 for (i = 0; i < surf->ntbins; ++i)
678 c = min(sin(surf->tbinsize*i), sin(surf->tbinsize*(i+1)))
679 * M_2PI / surf->targetbinsize + 1;
685 for (j = 0; j < c; ++j)
687 surf->tbin[i].p[j].left = -M_PI + j*M_2PI/c - 0.0001;
688 surf->tbin[i].p[j].bin = surf->nbins;
689 surf->bin[surf->nbins].n = 0;
692 surf->tbin[i].p[c].left = M_PI + 0.0001;
693 surf->tbin[i].p[c].bin = -1;
698 * \param[in,out] surf Surface data structure.
701 free_surface_points(t_methoddata_insolidangle *surf)
705 for (i = 0; i < surf->nbins; ++i)
709 sfree(surf->bin[i].x);
711 surf->bin[i].n_alloc = 0;
712 surf->bin[i].x = NULL;
717 * \param[in,out] surf Surface data structure.
718 * \param[in] tbin Bin number in the zenith angle direction.
719 * \param[in] pbin Bin number in the azimuthal angle direction.
720 * \param[in] x Point to store.
723 add_surface_point(t_methoddata_insolidangle *surf, int tbin, int pbin, rvec x)
727 bin = surf->tbin[tbin].p[pbin].bin;
728 /* Return if bin is already completely covered */
729 if (surf->bin[bin].n == -1)
731 /* Allocate more space if necessary */
732 if (surf->bin[bin].n == surf->bin[bin].n_alloc) {
733 surf->bin[bin].n_alloc += 10;
734 srenew(surf->bin[bin].x, surf->bin[bin].n_alloc);
736 /* Add the point to the bin */
737 copy_rvec(x, surf->bin[bin].x[surf->bin[bin].n]);
742 * \param[in,out] surf Surface data structure.
743 * \param[in] tbin Bin number in the zenith angle direction.
744 * \param[in] pbin Bin number in the azimuthal angle direction.
747 mark_surface_covered(t_methoddata_insolidangle *surf, int tbin, int pbin)
751 bin = surf->tbin[tbin].p[pbin].bin;
752 surf->bin[bin].n = -1;
756 * \param[in,out] surf Surface data structure.
757 * \param[in] tbin Bin number in the zenith angle direction.
758 * \param[in] phi Azimuthal angle of \p x.
759 * \param[in] pdelta1 Width of the cone at the lower edge of \p tbin.
760 * \param[in] pdelta2 Width of the cone at the uppper edge of \p tbin.
761 * \param[in] pdeltamax Max. width of the cone inside \p tbin.
762 * \param[in] x Point to store (should have unit length).
765 update_surface_bin(t_methoddata_insolidangle *surf, int tbin,
766 real phi, real pdelta1, real pdelta2, real pdeltamax,
769 real pdelta, phi1, phi2;
770 int pbin1, pbin2, pbin;
772 /* Find the edges of the bins affected */
773 pdelta = max(max(pdelta1, pdelta2), pdeltamax);
784 pbin1 = find_partition_bin(&surf->tbin[tbin], phi1);
785 pbin2 = find_partition_bin(&surf->tbin[tbin], phi2);
786 /* Find the edges of completely covered region */
787 pdelta = min(pdelta1, pdelta2);
794 /* Loop over all affected bins */
798 /* Wrap bin around if end reached */
799 if (pbin == surf->tbin[tbin].n)
805 /* Check if bin is completely covered and update */
806 if (surf->tbin[tbin].p[pbin].left >= phi1
807 && surf->tbin[tbin].p[pbin+1].left <= phi2)
809 mark_surface_covered(surf, tbin, pbin);
813 add_surface_point(surf, tbin, pbin, x);
816 while (pbin++ != pbin2); /* Loop including pbin2 */
820 * \param[in,out] surf Surface data structure.
821 * \param[in] x Point to store (should have unit length).
823 * Finds all the bins covered by the cone centered at \p x and calls
824 * update_surface_bin() to update them.
827 store_surface_point(t_methoddata_insolidangle *surf, rvec x)
830 real pdeltamax, tmax;
831 real theta1, theta2, pdelta1, pdelta2;
835 phi = atan2(x[YY], x[XX]);
836 /* Find the maximum extent in the phi direction */
837 if (theta <= surf->angcut)
842 else if (theta >= M_PI - surf->angcut)
849 pdeltamax = asin(sin(surf->angcut) / sin(theta));
850 tmax = acos(cos(theta) / cos(surf->angcut));
852 /* Find the first affected bin */
853 tbin = max(floor((theta - surf->angcut) / surf->tbinsize), 0.0);
854 theta1 = tbin * surf->tbinsize;
855 if (theta1 < theta - surf->angcut)
863 /* Loop through all affected bins */
864 while (tbin < ceil((theta + surf->angcut) / surf->tbinsize)
865 && tbin < surf->ntbins)
867 /* Calculate the next boundaries */
868 theta2 = (tbin+1) * surf->tbinsize;
869 if (theta2 > theta + surf->angcut)
873 else if (tbin == surf->ntbins - 1)
879 pdelta2 = 2*asin(sqrt(
880 (sqr(sin(surf->angcut/2)) - sqr(sin((theta2-theta)/2))) /
881 (sin(theta) * sin(theta2))));
884 if (tmax >= theta1 && tmax <= theta2)
886 update_surface_bin(surf, tbin, phi, pdelta1, pdelta2, pdeltamax, x);
890 update_surface_bin(surf, tbin, phi, pdelta1, pdelta2, 0, x);
900 * \param[in,out] surf Surface data structure.
902 * Currently, this function does nothing.
905 optimize_surface_points(t_methoddata_insolidangle *surf)
907 /* TODO: Implement */
911 * \param[in] surf Surface data structure.
912 * \returns An estimate for the area covered by the reference points.
915 estimate_covered_fraction(t_methoddata_insolidangle *surf)
918 real cfrac, tfrac, pfrac;
921 for (t = 0; t < surf->ntbins; ++t)
923 tfrac = cos(t * surf->tbinsize) - cos((t+1) * surf->tbinsize);
924 for (p = 0; p < surf->tbin[t].n; ++p)
926 pfrac = surf->tbin[t].p[p+1].left - surf->tbin[t].p[p].left;
927 n = surf->bin[surf->tbin[t].p[p].bin].n;
928 if (n == -1) /* Bin completely covered */
930 cfrac += tfrac * pfrac;
932 else if (n > 0) /* Bin partially covered */
934 cfrac += tfrac * pfrac / 2; /* A rough estimate */
938 return cfrac / (4*M_PI);
942 * \param[in] surf Surface data structure to search.
943 * \param[in] x Unit vector to check.
944 * \returns TRUE if \p x is within the solid angle, FALSE otherwise.
947 is_surface_covered(t_methoddata_insolidangle *surf, rvec x)
951 bin = find_surface_bin(surf, x);
952 /* Check for completely covered bin */
953 if (surf->bin[bin].n == -1)
957 /* Check each point that partially covers the bin */
958 for (i = 0; i < surf->bin[bin].n; ++i)
960 if (sph_distc(x, surf->bin[bin].x[i]) < surf->distccut)