commit 0dfae75c6705adbb51caba57462236719996e1ee
parent ddd8500078660a9a8f8d11bb5ab9b78c6cd804a8
Author: Vincent Forest <vincent.forest@meso-star.com>
Date: Wed, 6 Jul 2016 11:43:29 +0200
Upd the self hit management and incr the lib verbosity
Diffstat:
1 file changed, 46 insertions(+), 10 deletions(-)
diff --git a/src/sgf_estimator.c b/src/sgf_estimator.c
@@ -13,6 +13,8 @@
* You should have received a copy of the GNU General Public License
* along with this program. If not, see <http://www.gnu.org/licenses/>. */
+#define _POSIX_C_SOURCE 200112L /* nextafterf support */
+
#include "sgf.h"
#include "sgf_device_c.h"
@@ -38,6 +40,10 @@ struct accum {
#define DARRAY_DATA struct accum
#include <rsys/dynamic_array.h>
+/* How many self intersections are tested on the same ray before an error
+ * occurs */
+#define NSELF_HITS_MAX 10
+
/* Estimator of the Gebhart Factor of a given face to all the other ones */
struct sgf_estimator {
struct darray_accum buf; /* Integrated radiative flux of each primitive */
@@ -88,16 +94,18 @@ setup_status
/* Return RES_BAD_OP on numerical issue. i.e. the radiative path is skipped */
static res_T
gebhart_radiative_path
- (struct accum* accums,
+ (struct sgf_device* dev,
+ struct accum* accums,
struct accum* accums_infinity,
struct ssp_rng* rng,
+ const float epsilon,
const size_t primitive_id,
struct sgf_scene_desc* desc)
{
struct s3d_hit hit;
struct s3d_primitive prim_from, prim;
struct s3d_attrib attrib;
-
+ size_t nself_hits = 0;
const double trans_min = 1.e-8; /* Minimum transmissivity threshold */
double proba_reflec_spec;
double transmissivity, emissivity, reflectivity, specularity;
@@ -105,7 +113,6 @@ gebhart_radiative_path
double infinite_radiative_flux = 0.0;
double sum_radiative_flux = 0.f;
#endif
-
float st[2] = {0.f, 0.f}; /* Parametric coordinate of the sampled position */
float vec0[3]; /* Temporary vector */
float normal[3]; /* Geometric normal */
@@ -114,7 +121,7 @@ gebhart_radiative_path
float range[2]; /* Traced ray range */
float u, v;
- ASSERT(accums && accums_infinity && rng && desc);
+ ASSERT(accums && accums_infinity && epsilon > 0.f && rng && desc);
/* Discard faces with no emissivity */
emissivity = desc->get_material_property(desc->material,
@@ -140,14 +147,24 @@ gebhart_radiative_path
transmissivity = 1.0;
for(;;) { /* Here we go */
-
prim_from = prim;
- range[0] = 1.e-6f, range[1] = FLT_MAX;
+ range[0] = epsilon, range[1] = FLT_MAX;
- for(;;) {
+ nself_hits = 0;
+ do { /* Handle auto intersection */
S3D(scene_trace_ray(desc->scene, pos, dir, range, NULL, &hit));
- if(S3D_PRIMITIVE_EQ(&hit.prim, &prim_from)) range[0] += 1.e-6f; /* Self hit */
- else break;
+ if(!S3D_HIT_NONE(&hit)) range[0] = MMAX(range[0], hit.distance);
+ range[0] = nextafterf(range[0], range[1]);
+ } while(S3D_PRIMITIVE_EQ(&hit.prim, &prim_from) /* Self hit */
+ && ++nself_hits < NSELF_HITS_MAX); /* Too many self hits */
+
+ if(nself_hits >= NSELF_HITS_MAX) {
+ if(dev->verbose) {
+ logger_print(dev->logger, LOG_ERROR,
+ "Too many self hit on a given ray. Ray origin: %g, %g, %g\n",
+ SPLIT3(pos));
+ }
+ return RES_BAD_OP;
}
if(S3D_HIT_NONE(&hit)) { /* The ray is outside the volume */
@@ -321,6 +338,8 @@ sgf_integrate
size_t istep;
size_t ispectral_band;
size_t nprims;
+ float epsilon;
+ float lower[3], upper[3];
int mask;
res_T res = RES_OK;
@@ -355,6 +374,22 @@ sgf_integrate
goto error;
}
+ /* Compute an epsilon relatively to the scene size */
+ S3D(scene_get_aabb(desc->scene, lower, upper));
+ epsilon = upper[0] - lower[0];
+ epsilon = MMIN(epsilon, upper[1] - lower[1]);
+ epsilon = MMIN(epsilon, upper[2] - lower[2]);
+ if(epsilon <= 0.f) {
+ if(dev->verbose) {
+ logger_print(dev->logger, LOG_ERROR,
+ "Degenerated scene geometry. lower: %g, %g, %g; upper: %g %g %g\n",
+ SPLIT3(lower), SPLIT3(upper));
+ }
+ res = RES_BAD_ARG;
+ goto error;
+ }
+ epsilon = MMAX(epsilon*1.e-6f, FLT_MIN);
+
/* Allocate and init the Gebhart Factor result buffer */
res = darray_accum_resize(&estimator->buf, nprims*desc->spectral_bands_count);
if(res != RES_OK) {
@@ -380,6 +415,7 @@ sgf_integrate
memset(darray_accum_data_get(&estimator->buf_infinity), 0,
darray_accum_size_get(&estimator->buf_infinity)*sizeof(struct accum));
+
/* Invoke the Gebhart factor integration. */
FOR_EACH(ispectral_band, 0, desc->spectral_bands_count) {
struct accum* accums =
@@ -389,7 +425,7 @@ sgf_integrate
FOR_EACH(istep, 0, steps_count) {
res = gebhart_radiative_path
- (accums, accums_infinity, rng, primitive_id, desc);
+ (dev, accums, accums_infinity, rng, epsilon, primitive_id, desc);
if(res != RES_OK) goto error;
}
}