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 * -# 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.
96 * \brief Implementation of the \ref sm_insolidangle "insolidangle"
100 * The implementation could be optimized quite a bit.
102 * \todo Move the covered fraction stuff somewhere else and make it more
103 * generic (along the lines it is handled in selection.h and trajana.h).
118 #include <indexutil.h>
119 #include <position.h>
120 #include <selection.h>
121 #include <selmethod.h>
126 * Internal data structure for the \p insolidangle selection method.
128 * \see \c t_partition
132 /** Left edge of the partition. */
134 /** Bin index corresponding to this partition. */
139 * Internal data structure for the \p insolidangle selection method.
141 * Describes the surface partitioning within one slice along the zenith angle.
142 * The slice from azimuthal angle \p p[i].left to \p p[i+1].left belongs to
147 /** Number of partition items (\p p contains \p n+1 items). */
149 /** Array of partition edges and corresponding bins. */
154 * Internal data structure for the \p insolidangle selection method.
156 * Contains the reference points that partially cover a certain region on the
157 * surface of the unit sphere.
158 * If \p n is -1, the whole region described by the bin is covered.
162 /** Number of points in the array \p x, -1 if whole bin covered. */
164 /** Number of elements allocated for \p x. */
166 /** Array of points that partially cover the bin. */
168 } t_spheresurfacebin
;
171 * Data structure for the \p insolidangle selection method.
173 * All angle values are in the units of radians.
177 /** Center of the solid angle. */
178 gmx_ana_pos_t center
;
179 /** Positions that span the solid angle. */
183 /** Estimate of the covered fraction. */
186 /** Cutoff for the cosine (equals cos(angcut)). */
188 /** Bin size to be used as the target bin size when constructing the bins. */
191 /** Number of bins in the \p tbin array. */
193 /** Size of one bin in the zenith angle direction. */
195 /** Array of zenith angle slices. */
197 /** Number of elements allocated for the \p bin array. */
199 /** Number of elements used in the \p bin array. */
201 /** Array of individual bins. */
202 t_spheresurfacebin
*bin
;
203 } t_methoddata_insolidangle
;
205 /** Allocates data for the \p insolidangle selection method. */
207 init_data_insolidangle(int npar
, gmx_ana_selparam_t
*param
);
208 /** Initializes the \p insolidangle selection method. */
210 init_insolidangle(t_topology
*top
, int npar
, gmx_ana_selparam_t
*param
, void *data
);
211 /** Sets the COM/COG data for the \p insolidangle selection method. */
213 set_comg_insolidangle(gmx_ana_pos_t
*pos
, void *data
);
214 /** Frees the data allocated for the \p insolidangle selection method. */
216 free_data_insolidangle(void *data
);
217 /** Initializes the evaluation of the \p insolidangle selection method for a frame. */
219 init_frame_insolidangle(t_topology
*top
, t_trxframe
*fr
, t_pbc
*pbc
, void *data
);
220 /** Internal helper function for evaluate_insolidangle(). */
222 accept_insolidangle(rvec x
, t_pbc
*pbc
, void *data
);
223 /** Evaluates the \p insolidangle selection method. */
225 evaluate_insolidangle(t_topology
*top
, t_trxframe
*fr
, t_pbc
*pbc
,
226 gmx_ana_pos_t
*pos
, gmx_ana_selvalue_t
*out
, void *data
);
228 /** Calculates the distance between unit vectors. */
230 sph_distc(rvec x1
, rvec x2
);
231 /** Does a binary search on a \p t_partition to find a bin for a value. */
233 find_partition_bin(t_partition
*p
, real value
);
234 /** Finds a bin that corresponds to a location on the unit sphere surface. */
236 find_surface_bin(t_methoddata_insolidangle
*surf
, rvec x
);
237 /** Clears/initializes the bins on the unit sphere surface. */
239 clear_surface_points(t_methoddata_insolidangle
*surf
);
240 /** Frees memory allocated for storing the reference points in the surface bins. */
242 free_surface_points(t_methoddata_insolidangle
*surf
);
243 /** Adds a reference point to a given bin. */
245 add_surface_point(t_methoddata_insolidangle
*surf
, int tbin
, int pbin
, rvec x
);
246 /** Marks a bin as completely covered. */
248 mark_surface_covered(t_methoddata_insolidangle
*surf
, int tbin
, int pbin
);
249 /** Helper function for store_surface_point() to update a single zenith angle bin. */
251 update_surface_bin(t_methoddata_insolidangle
*surf
, int tbin
,
252 real phi
, real pdelta1
, real pdelta2
, real pdeltamax
,
254 /** Adds a single reference point and updates the surface bins. */
256 store_surface_point(t_methoddata_insolidangle
*surf
, rvec x
);
257 /** Optimizes the surface bins for faster searching. */
259 optimize_surface_points(t_methoddata_insolidangle
*surf
);
260 /** Estimates the area covered by the reference cones. */
262 estimate_covered_fraction(t_methoddata_insolidangle
*surf
);
263 /** Checks whether a point lies within a solid angle. */
265 is_surface_covered(t_methoddata_insolidangle
*surf
, rvec x
);
267 /** Parameters for the \p insolidangle selection method. */
268 static gmx_ana_selparam_t smparams_insolidangle
[] = {
269 {"center", {POS_VALUE
, 1, {NULL
}}, NULL
, SPAR_DYNAMIC
},
270 {"span", {POS_VALUE
, -1, {NULL
}}, NULL
, SPAR_DYNAMIC
| SPAR_VARNUM
},
271 {"cutoff", {REAL_VALUE
, 1, {NULL
}}, NULL
, SPAR_OPTIONAL
},
274 /** Help text for the \p insolidangle selection method. */
275 static const char *help_insolidangle
[] = {
276 "SELECTING ATOMS IN A SOLID ANGLE[PAR]",
278 "[TT]insolidangle center POS span POS_EXPR [cutoff REAL][tt][PAR]",
280 "This keyword selects atoms that are within [TT]REAL[tt] degrees",
281 "(default=5) of any position in [TT]POS_EXPR[tt] as seen from [TT]POS[tt]",
282 "a position expression that evaluates to a single position), i.e., atoms",
283 "in the solid angle spanned by the positions in [TT]POS_EXPR[tt] and",
284 "centered at [TT]POS[tt].[PAR]"
286 "Technically, the solid angle is constructed as a union of small cones",
287 "whose tip is at [TT]POS[tt] and the axis goes through a point in",
288 "[TT]POS_EXPR[tt]. There is such a cone for each position in",
289 "[TT]POS_EXPR[tt], and point is in the solid angle if it lies within any",
290 "of these cones. The cutoff determines the width of the cones.",
293 /** \internal Selection method data for the \p insolidangle method. */
294 gmx_ana_selmethod_t sm_insolidangle
= {
295 "insolidangle", GROUP_VALUE
, SMETH_DYNAMIC
,
296 asize(smparams_insolidangle
), smparams_insolidangle
,
297 &init_data_insolidangle
,
301 &free_data_insolidangle
,
302 &init_frame_insolidangle
,
304 &evaluate_insolidangle
,
305 {"insolidangle center POS span POS_EXPR [cutoff REAL]",
306 asize(help_insolidangle
), help_insolidangle
},
310 * \param[in] npar Not used (should be 3).
311 * \param[in,out] param Method parameters (should point to
312 * \ref smparams_insolidangle).
313 * \returns Pointer to the allocated data (\ref t_methoddata_insolidangle).
315 * Allocates memory for a \ref t_methoddata_insolidangle structure and
316 * initializes the parameter as follows:
317 * - \p center defines the value for t_methoddata_insolidangle::center.
318 * - \p span defines the value for t_methoddata_insolidangle::span.
319 * - \p cutoff defines the value for t_methoddata_insolidangle::angcut.
322 init_data_insolidangle(int npar
, gmx_ana_selparam_t
*param
)
324 t_methoddata_insolidangle
*data
;
328 param
[0].val
.u
.p
= &data
->center
;
329 param
[1].val
.u
.p
= &data
->span
;
330 param
[2].val
.u
.r
= &data
->angcut
;
335 * \param top Not used.
336 * \param npar Not used.
337 * \param param Not used.
338 * \param data Pointer to \ref t_methoddata_insolidangle to initialize.
339 * \returns 0 on success, -1 on failure.
341 * Converts t_methoddata_insolidangle::angcut to radians and allocates
342 * and allocates memory for the bins used during the evaluation.
345 init_insolidangle(t_topology
*top
, int npar
, gmx_ana_selparam_t
*param
, void *data
)
347 t_methoddata_insolidangle
*surf
= (t_methoddata_insolidangle
*)data
;
350 if (surf
->angcut
<= 0)
352 fprintf(stderr
, "error: angle cutoff should be > 0");
356 surf
->angcut
*= DEG2RAD
;
358 surf
->distccut
= -cos(surf
->angcut
);
359 surf
->targetbinsize
= surf
->angcut
/ 2;
360 surf
->ntbins
= (int) (M_PI
/ surf
->targetbinsize
);
361 surf
->tbinsize
= (180.0 / surf
->ntbins
)*DEG2RAD
;
363 snew(surf
->tbin
, (int)(M_PI
/surf
->tbinsize
) + 1);
365 for (i
= 0; i
< surf
->ntbins
; ++i
)
367 c
= max(sin(surf
->tbinsize
*i
), sin(surf
->tbinsize
*(i
+1)))
368 * M_2PI
/ surf
->targetbinsize
+ 1;
369 snew(surf
->tbin
[i
].p
, c
+1);
373 snew(surf
->bin
, surf
->maxbins
);
379 * \param data Data to free (should point to a \ref t_methoddata_insolidangle).
381 * Frees the memory allocated for \c t_methoddata_insolidangle::center and
382 * \c t_methoddata_insolidangle::span, as well as the memory for the internal
386 free_data_insolidangle(void *data
)
388 t_methoddata_insolidangle
*d
= (t_methoddata_insolidangle
*)data
;
391 gmx_ana_pos_deinit(&d
->center
);
392 gmx_ana_pos_deinit(&d
->span
);
395 for (i
= 0; i
< d
->ntbins
; ++i
)
401 free_surface_points(d
);
406 * \param[in] top Not used.
407 * \param[in] fr Current frame.
408 * \param[in] pbc PBC structure.
409 * \param data Should point to a \ref t_methoddata_insolidangle.
410 * \returns 0 on success, a non-zero error code on error.
412 * Creates a lookup structure that enables fast queries of whether a point
413 * is within the solid angle or not.
416 init_frame_insolidangle(t_topology
*top
, t_trxframe
*fr
, t_pbc
*pbc
, void *data
)
418 t_methoddata_insolidangle
*d
= (t_methoddata_insolidangle
*)data
;
422 free_surface_points(d
);
423 clear_surface_points(d
);
424 for (i
= 0; i
< d
->span
.nr
; ++i
)
428 pbc_dx(pbc
, d
->span
.x
[i
], d
->center
.x
[0], dx
);
432 rvec_sub(d
->span
.x
[i
], d
->center
.x
[0], dx
);
435 store_surface_point(d
, dx
);
437 optimize_surface_points(d
);
443 * \param[in] x Test point.
444 * \param[in] pbc PBC data (if NULL, no PBC are used).
445 * \param[in] data Pointer to a \c t_methoddata_insolidangle data structure.
446 * \returns TRUE if \p x is within the solid angle, FALSE otherwise.
449 accept_insolidangle(rvec x
, t_pbc
*pbc
, void *data
)
451 t_methoddata_insolidangle
*d
= (t_methoddata_insolidangle
*)data
;
456 pbc_dx(pbc
, x
, d
->center
.x
[0], dx
);
460 rvec_sub(x
, d
->center
.x
[0], dx
);
463 return is_surface_covered(d
, dx
);
467 * See sel_updatefunc() for description of the parameters.
468 * \p data should point to a \c t_methoddata_insolidangle.
470 * Calculates which atoms in \p g are within the solid angle spanned by
471 * \c t_methoddata_insolidangle::span and centered at
472 * \c t_methoddata_insolidangle::center, and stores the result in \p out->u.g.
475 evaluate_insolidangle(t_topology
*top
, t_trxframe
*fr
, t_pbc
*pbc
,
476 gmx_ana_pos_t
*pos
, gmx_ana_selvalue_t
*out
, void *data
)
478 t_methoddata_insolidangle
*d
= (t_methoddata_insolidangle
*)data
;
482 for (b
= 0; b
< pos
->nr
; ++b
)
484 if (accept_insolidangle(pos
->x
[b
], pbc
, data
))
486 gmx_ana_pos_append(NULL
, out
->u
.g
, pos
, b
, 0);
493 * \param[in] sel Selection element to query.
494 * \returns TRUE if the covered fraction can be estimated for \p sel with
495 * _gmx_selelem_estimate_coverfrac(), FALSE otherwise.
498 _gmx_selelem_can_estimate_cover(t_selelem
*sel
)
504 if (sel
->type
== SEL_BOOLEAN
&& sel
->u
.boolt
== BOOL_OR
)
513 if (child
->type
== SEL_EXPRESSION
)
515 if (child
->u
.expr
.method
->name
== sm_insolidangle
.name
)
517 if (bFound
|| bDynFound
)
523 else if (child
->u
.expr
.method
524 && (child
->u
.expr
.method
->flags
& SMETH_DYNAMIC
))
533 else if (!_gmx_selelem_can_estimate_cover(child
))
543 * \param[in] sel Selection for which the fraction should be calculated.
544 * \returns Fraction of angles covered by the selection (between zero and one).
546 * The return value is undefined if _gmx_selelem_can_estimate_cover() returns
548 * Should be called after gmx_ana_evaluate_selections() has been called for the
552 _gmx_selelem_estimate_coverfrac(t_selelem
*sel
)
557 if (sel
->type
== SEL_EXPRESSION
&& sel
->u
.expr
.method
->name
== sm_insolidangle
.name
)
559 t_methoddata_insolidangle
*d
= (t_methoddata_insolidangle
*)sel
->u
.expr
.mdata
;
562 d
->cfrac
= estimate_covered_fraction(d
);
566 if (sel
->type
== SEL_BOOLEAN
&& sel
->u
.boolt
== BOOL_NOT
)
568 cfrac
= _gmx_selelem_estimate_coverfrac(sel
->child
);
576 /* Here, we assume that the selection is simple enough */
580 cfrac
= _gmx_selelem_estimate_coverfrac(child
);
591 * \param[in] x1 Unit vector 1.
592 * \param[in] x2 Unit vector 2.
593 * \returns Minus the dot product of \p x1 and \p x2.
595 * This function is used internally to calculate the distance between the
596 * unit vectors \p x1 and \p x2 to find out whether \p x2 is within the
597 * cone centered at \p x1. Currently, the cosine of the angle is used
598 * for efficiency, and the minus is there to make it behave like a normal
599 * distance (larger values mean longer distances).
602 sph_distc(rvec x1
, rvec x2
)
604 return -iprod(x1
, x2
);
608 * \param[in] p Partition to search.
609 * \param[in] value Value to search for.
610 * \returns The partition index in \p p that contains \p value.
612 * If \p value is outside the range of \p p, the first/last index is returned.
613 * Otherwise, the return value \c i satisfies \c p->p[i].left<=value and
614 * \c p->p[i+1].left>value
617 find_partition_bin(t_partition
*p
, real value
)
619 int pmin
, pmax
, pbin
;
621 /* Binary search the partition */
622 pmin
= 0; pmax
= p
->n
;
623 while (pmax
> pmin
+ 1)
625 pbin
= pmin
+ (pmax
- pmin
) / 2;
626 if (p
->p
[pbin
].left
<= value
)
640 * \param[in] surf Surface data structure to search.
641 * \param[in] x Unit vector to find.
642 * \returns The bin index that contains \p x.
644 * The return value is an index to the \p surf->bin array.
647 find_surface_bin(t_methoddata_insolidangle
*surf
, rvec x
)
653 phi
= atan2(x
[YY
], x
[XX
]);
654 tbin
= floor(theta
/ surf
->tbinsize
);
655 if (tbin
>= surf
->ntbins
)
657 tbin
= surf
->ntbins
- 1;
659 pbin
= find_partition_bin(&surf
->tbin
[tbin
], phi
);
660 return surf
->tbin
[tbin
].p
[pbin
].bin
;
664 * \param[in,out] surf Surface data structure.
666 * Clears the reference points from the bins and (re)initializes the edges
667 * of the azimuthal bins.
670 clear_surface_points(t_methoddata_insolidangle
*surf
)
675 for (i
= 0; i
< surf
->ntbins
; ++i
)
677 c
= min(sin(surf
->tbinsize
*i
), sin(surf
->tbinsize
*(i
+1)))
678 * M_2PI
/ surf
->targetbinsize
+ 1;
684 for (j
= 0; j
< c
; ++j
)
686 surf
->tbin
[i
].p
[j
].left
= -M_PI
+ j
*M_2PI
/c
- 0.0001;
687 surf
->tbin
[i
].p
[j
].bin
= surf
->nbins
;
688 surf
->bin
[surf
->nbins
].n
= 0;
691 surf
->tbin
[i
].p
[c
].left
= M_PI
+ 0.0001;
692 surf
->tbin
[i
].p
[c
].bin
= -1;
697 * \param[in,out] surf Surface data structure.
700 free_surface_points(t_methoddata_insolidangle
*surf
)
704 for (i
= 0; i
< surf
->nbins
; ++i
)
708 sfree(surf
->bin
[i
].x
);
710 surf
->bin
[i
].n_alloc
= 0;
711 surf
->bin
[i
].x
= NULL
;
716 * \param[in,out] surf Surface data structure.
717 * \param[in] tbin Bin number in the zenith angle direction.
718 * \param[in] pbin Bin number in the azimuthal angle direction.
719 * \param[in] x Point to store.
722 add_surface_point(t_methoddata_insolidangle
*surf
, int tbin
, int pbin
, rvec x
)
726 bin
= surf
->tbin
[tbin
].p
[pbin
].bin
;
727 /* Return if bin is already completely covered */
728 if (surf
->bin
[bin
].n
== -1)
730 /* Allocate more space if necessary */
731 if (surf
->bin
[bin
].n
== surf
->bin
[bin
].n_alloc
) {
732 surf
->bin
[bin
].n_alloc
+= 10;
733 srenew(surf
->bin
[bin
].x
, surf
->bin
[bin
].n_alloc
);
735 /* Add the point to the bin */
736 copy_rvec(x
, surf
->bin
[bin
].x
[surf
->bin
[bin
].n
]);
741 * \param[in,out] surf Surface data structure.
742 * \param[in] tbin Bin number in the zenith angle direction.
743 * \param[in] pbin Bin number in the azimuthal angle direction.
746 mark_surface_covered(t_methoddata_insolidangle
*surf
, int tbin
, int pbin
)
750 bin
= surf
->tbin
[tbin
].p
[pbin
].bin
;
751 surf
->bin
[bin
].n
= -1;
755 * \param[in,out] surf Surface data structure.
756 * \param[in] tbin Bin number in the zenith angle direction.
757 * \param[in] phi Azimuthal angle of \p x.
758 * \param[in] pdelta1 Width of the cone at the lower edge of \p tbin.
759 * \param[in] pdelta2 Width of the cone at the uppper edge of \p tbin.
760 * \param[in] pdeltamax Max. width of the cone inside \p tbin.
761 * \param[in] x Point to store (should have unit length).
764 update_surface_bin(t_methoddata_insolidangle
*surf
, int tbin
,
765 real phi
, real pdelta1
, real pdelta2
, real pdeltamax
,
768 real pdelta
, phi1
, phi2
;
769 int pbin1
, pbin2
, pbin
;
771 /* Find the edges of the bins affected */
772 pdelta
= max(max(pdelta1
, pdelta2
), pdeltamax
);
783 pbin1
= find_partition_bin(&surf
->tbin
[tbin
], phi1
);
784 pbin2
= find_partition_bin(&surf
->tbin
[tbin
], phi2
);
785 /* Find the edges of completely covered region */
786 pdelta
= min(pdelta1
, pdelta2
);
793 /* Loop over all affected bins */
797 /* Wrap bin around if end reached */
798 if (pbin
== surf
->tbin
[tbin
].n
)
804 /* Check if bin is completely covered and update */
805 if (surf
->tbin
[tbin
].p
[pbin
].left
>= phi1
806 && surf
->tbin
[tbin
].p
[pbin
+1].left
<= phi2
)
808 mark_surface_covered(surf
, tbin
, pbin
);
812 add_surface_point(surf
, tbin
, pbin
, x
);
815 while (pbin
++ != pbin2
); /* Loop including pbin2 */
819 * \param[in,out] surf Surface data structure.
820 * \param[in] x Point to store (should have unit length).
822 * Finds all the bins covered by the cone centered at \p x and calls
823 * update_surface_bin() to update them.
826 store_surface_point(t_methoddata_insolidangle
*surf
, rvec x
)
829 real pdeltamax
, tmax
;
830 real theta1
, theta2
, pdelta1
, pdelta2
;
834 phi
= atan2(x
[YY
], x
[XX
]);
835 /* Find the maximum extent in the phi direction */
836 if (theta
<= surf
->angcut
)
841 else if (theta
>= M_PI
- surf
->angcut
)
848 pdeltamax
= asin(sin(surf
->angcut
) / sin(theta
));
849 tmax
= acos(cos(theta
) / cos(surf
->angcut
));
851 /* Find the first affected bin */
852 tbin
= max(floor((theta
- surf
->angcut
) / surf
->tbinsize
), 0);
853 theta1
= tbin
* surf
->tbinsize
;
854 if (theta1
< theta
- surf
->angcut
)
862 /* Loop through all affected bins */
863 while (tbin
< ceil((theta
+ surf
->angcut
) / surf
->tbinsize
)
864 && tbin
< surf
->ntbins
)
866 /* Calculate the next boundaries */
867 theta2
= (tbin
+1) * surf
->tbinsize
;
868 if (theta2
> theta
+ surf
->angcut
)
872 else if (tbin
== surf
->ntbins
- 1)
878 pdelta2
= 2*asin(sqrt(
879 (sqr(sin(surf
->angcut
/2)) - sqr(sin((theta2
-theta
)/2))) /
880 (sin(theta
) * sin(theta2
))));
883 if (tmax
>= theta1
&& tmax
<= theta2
)
885 update_surface_bin(surf
, tbin
, phi
, pdelta1
, pdelta2
, pdeltamax
, x
);
889 update_surface_bin(surf
, tbin
, phi
, pdelta1
, pdelta2
, 0, x
);
899 * \param[in,out] surf Surface data structure.
901 * Currently, this function does nothing.
904 optimize_surface_points(t_methoddata_insolidangle
*surf
)
906 /* TODO: Implement */
910 * \param[in] surf Surface data structure.
911 * \returns An estimate for the area covered by the reference points.
914 estimate_covered_fraction(t_methoddata_insolidangle
*surf
)
917 real cfrac
, tfrac
, pfrac
;
920 for (t
= 0; t
< surf
->ntbins
; ++t
)
922 tfrac
= cos(t
* surf
->tbinsize
) - cos((t
+1) * surf
->tbinsize
);
923 for (p
= 0; p
< surf
->tbin
[t
].n
; ++p
)
925 pfrac
= surf
->tbin
[t
].p
[p
+1].left
- surf
->tbin
[t
].p
[p
].left
;
926 n
= surf
->bin
[surf
->tbin
[t
].p
[p
].bin
].n
;
927 if (n
== -1) /* Bin completely covered */
929 cfrac
+= tfrac
* pfrac
;
931 else if (n
> 0) /* Bin partially covered */
933 cfrac
+= tfrac
* pfrac
/ 2; /* A rough estimate */
937 return cfrac
/ (4*M_PI
);
941 * \param[in] surf Surface data structure to search.
942 * \param[in] x Unit vector to check.
943 * \returns TRUE if \p x is within the solid angle, FALSE otherwise.
946 is_surface_covered(t_methoddata_insolidangle
*surf
, rvec x
)
950 bin
= find_surface_bin(surf
, x
);
951 /* Check for completely covered bin */
952 if (surf
->bin
[bin
].n
== -1)
956 /* Check each point that partially covers the bin */
957 for (i
= 0; i
< surf
->bin
[bin
].n
; ++i
)
959 if (sph_distc(x
, surf
->bin
[bin
].x
[i
]) < surf
->distccut
)