rnatm

Load and structure data describing an atmosphere
git clone git://git.meso-star.fr/rnatm.git
Log | Files | Refs | README | LICENSE

commit 6c8ce86b883bdceaf0aef71f2042ab0581b9a249
parent 312ed2569dfa030597764cfda5fb66ddb177a4d2
Author: Vincent Forest <vincent.forest@meso-star.com>
Date:   Thu, 28 Jul 2022 09:03:28 +0200

Updating the memory layout of the partition

In fact, this commit undoes a previous update: one partition only
manages one cluster of voxels belonging to a grid.

Diffstat:
Msrc/rnatm_octree.c | 157++++++++++++++++++++++++++++++++++---------------------------------------------
Msrc/rnatm_voxel_partition.c | 148++++++++++++-------------------------------------------------------------------
Msrc/rnatm_voxel_partition.h | 8+-------
3 files changed, 90 insertions(+), 223 deletions(-)

diff --git a/src/rnatm_octree.c b/src/rnatm_octree.c @@ -51,17 +51,6 @@ static const struct build_octrees_context BUILD_OCTREES_CONTEXT_NULL = /******************************************************************************* * Helper functions ******************************************************************************/ -static size_t -band_get_quad_pt_count(const size_t iband, void* ctx) -{ - struct sck_band band; - struct sck* sck = ctx; - ASSERT(ctx); - - SCK(get_band(sck, iband, &band)); - return band.quad_pts_count; -} - static FINLINE unsigned round_pow2(const unsigned val) { @@ -148,89 +137,79 @@ update_voxel (struct rnatm* atm, const struct suvm_primitive* tetra, struct partition* part, + const size_t iband, + const size_t iquad_pt, const uint64_t vx_mcode) { - size_t iquad_pt; - size_t iband; - size_t nbands; + struct sck_band band; + struct sck_quad_pt quad_pt; + + float tetra_ks[4]; + float tetra_ka[4]; + float tetra_kext[4]; + float tetra_ks_min, tetra_ks_max; + float tetra_ka_min, tetra_ka_max; + float tetra_kext_min, tetra_kext_max; + + float vx_ka_min, vx_ka_max; + float vx_ks_min, vx_ks_max; + float vx_kext_min, vx_kext_max; + float* vx = NULL; ASSERT(atm && tetra && part); ASSERT(tetra->nvertices == 4); - nbands = sck_get_bands_count(atm->gas.ck); - nbands = 1; /* FIXME For tests only */ - #define REDUX_MIN(V4) MMIN(MMIN((V4)[0], (V4)[1]), MMIN((V4)[2], V4[3])) #define REDUX_MAX(V4) MMAX(MMAX((V4)[0], (V4)[1]), MMAX((V4)[2], V4[3])) - FOR_EACH(iband, 0, nbands) { - struct sck_band band; - size_t nquad_pts; - float tetra_ks[4]; - float tetra_ks_min, tetra_ks_max; - - /* Compute the scattering coefficient range of the tetrahedron */ - SCK(get_band(atm->gas.ck, iband, &band)); - tetra_ks[0] = band.ks_list[tetra->indices[0]]; - tetra_ks[1] = band.ks_list[tetra->indices[1]]; - tetra_ks[2] = band.ks_list[tetra->indices[2]]; - tetra_ks[3] = band.ks_list[tetra->indices[3]]; - tetra_ks_min = REDUX_MIN(tetra_ks); - tetra_ks_max = REDUX_MAX(tetra_ks); - - nquad_pts = band.quad_pts_count; - nquad_pts = 1; /* FIXME For tests only */ - - FOR_EACH(iquad_pt, 0, nquad_pts) { - struct sck_quad_pt quad_pt; - float tetra_ka[4]; - float tetra_kext[4]; - float tetra_ka_min, tetra_ka_max; - float tetra_kext_min, tetra_kext_max; - float vx_ka_min, vx_ka_max; - float vx_ks_min, vx_ks_max; - float vx_kext_min, vx_kext_max; - - /* Compute the absorption coefficient range of the tetrahedron */ - SCK(band_get_quad_pt(&band, iquad_pt, &quad_pt)); - tetra_ka[0] = quad_pt.ka_list[tetra->indices[0]]; - tetra_ka[1] = quad_pt.ka_list[tetra->indices[1]]; - tetra_ka[2] = quad_pt.ka_list[tetra->indices[2]]; - tetra_ka[3] = quad_pt.ka_list[tetra->indices[3]]; - tetra_ka_min = REDUX_MIN(tetra_ka); - tetra_ka_max = REDUX_MAX(tetra_ka); - - /* Compute the extinction coefficient range of the tetrahedron */ - tetra_kext[0] = tetra_ka[0] + tetra_ks[0]; - tetra_kext[1] = tetra_ka[1] + tetra_ks[1]; - tetra_kext[2] = tetra_ka[2] + tetra_ks[2]; - tetra_kext[3] = tetra_ka[3] + tetra_ks[3]; - tetra_kext_min = REDUX_MIN(tetra_kext); - tetra_kext_max = REDUX_MAX(tetra_kext); - - vx = partition_get_voxel(part, iband, iquad_pt, vx_mcode); - - /* Update the range of the radiative coefficients of the voxel */ - vx_ka_min = vx[voxel_idata(RNATM_RADCOEF_Ka, RNATM_SVX_OP_MIN)]; - vx_ka_max = vx[voxel_idata(RNATM_RADCOEF_Ka, RNATM_SVX_OP_MAX)]; - vx_ks_min = vx[voxel_idata(RNATM_RADCOEF_Ks, RNATM_SVX_OP_MIN)]; - vx_ks_max = vx[voxel_idata(RNATM_RADCOEF_Ks, RNATM_SVX_OP_MAX)]; - vx_kext_min = vx[voxel_idata(RNATM_RADCOEF_Kext, RNATM_SVX_OP_MIN)]; - vx_kext_max = vx[voxel_idata(RNATM_RADCOEF_Kext, RNATM_SVX_OP_MAX)]; - vx_ka_min = MMIN(vx_ka_min, tetra_ka_min); - vx_ka_max = MMAX(vx_ka_max, tetra_ka_max); - vx_ks_min = MMIN(vx_ks_min, tetra_ks_min); - vx_ks_max = MMAX(vx_ks_max, tetra_ks_max); - vx_kext_min = MMIN(vx_kext_min, tetra_kext_min); - vx_kext_max = MMAX(vx_kext_max, tetra_kext_max); - vx[voxel_idata(RNATM_RADCOEF_Ka, RNATM_SVX_OP_MIN)] = vx_ka_min; - vx[voxel_idata(RNATM_RADCOEF_Ka, RNATM_SVX_OP_MAX)] = vx_ka_max; - vx[voxel_idata(RNATM_RADCOEF_Ks, RNATM_SVX_OP_MIN)] = vx_ks_min; - vx[voxel_idata(RNATM_RADCOEF_Ks, RNATM_SVX_OP_MAX)] = vx_ks_max; - vx[voxel_idata(RNATM_RADCOEF_Kext, RNATM_SVX_OP_MIN)] = vx_kext_min; - vx[voxel_idata(RNATM_RADCOEF_Kext, RNATM_SVX_OP_MAX)] = vx_kext_max; - } - } + /* Compute the scattering coefficient range of the tetrahedron */ + SCK(get_band(atm->gas.ck, iband, &band)); + tetra_ks[0] = band.ks_list[tetra->indices[0]]; + tetra_ks[1] = band.ks_list[tetra->indices[1]]; + tetra_ks[2] = band.ks_list[tetra->indices[2]]; + tetra_ks[3] = band.ks_list[tetra->indices[3]]; + tetra_ks_min = REDUX_MIN(tetra_ks); + tetra_ks_max = REDUX_MAX(tetra_ks); + + /* Compute the absorption coefficient range of the tetrahedron */ + SCK(band_get_quad_pt(&band, iquad_pt, &quad_pt)); + tetra_ka[0] = quad_pt.ka_list[tetra->indices[0]]; + tetra_ka[1] = quad_pt.ka_list[tetra->indices[1]]; + tetra_ka[2] = quad_pt.ka_list[tetra->indices[2]]; + tetra_ka[3] = quad_pt.ka_list[tetra->indices[3]]; + tetra_ka_min = REDUX_MIN(tetra_ka); + tetra_ka_max = REDUX_MAX(tetra_ka); + + /* Compute the extinction coefficient range of the tetrahedron */ + tetra_kext[0] = tetra_ka[0] + tetra_ks[0]; + tetra_kext[1] = tetra_ka[1] + tetra_ks[1]; + tetra_kext[2] = tetra_ka[2] + tetra_ks[2]; + tetra_kext[3] = tetra_ka[3] + tetra_ks[3]; + tetra_kext_min = REDUX_MIN(tetra_kext); + tetra_kext_max = REDUX_MAX(tetra_kext); + + vx = partition_get_voxel(part, vx_mcode); + + /* Update the range of the radiative coefficients of the voxel */ + vx_ka_min = vx[voxel_idata(RNATM_RADCOEF_Ka, RNATM_SVX_OP_MIN)]; + vx_ka_max = vx[voxel_idata(RNATM_RADCOEF_Ka, RNATM_SVX_OP_MAX)]; + vx_ks_min = vx[voxel_idata(RNATM_RADCOEF_Ks, RNATM_SVX_OP_MIN)]; + vx_ks_max = vx[voxel_idata(RNATM_RADCOEF_Ks, RNATM_SVX_OP_MAX)]; + vx_kext_min = vx[voxel_idata(RNATM_RADCOEF_Kext, RNATM_SVX_OP_MIN)]; + vx_kext_max = vx[voxel_idata(RNATM_RADCOEF_Kext, RNATM_SVX_OP_MAX)]; + vx_ka_min = MMIN(vx_ka_min, tetra_ka_min); + vx_ka_max = MMAX(vx_ka_max, tetra_ka_max); + vx_ks_min = MMIN(vx_ks_min, tetra_ks_min); + vx_ks_max = MMAX(vx_ks_max, tetra_ks_max); + vx_kext_min = MMIN(vx_kext_min, tetra_kext_min); + vx_kext_max = MMAX(vx_kext_max, tetra_kext_max); + vx[voxel_idata(RNATM_RADCOEF_Ka, RNATM_SVX_OP_MIN)] = vx_ka_min; + vx[voxel_idata(RNATM_RADCOEF_Ka, RNATM_SVX_OP_MAX)] = vx_ka_max; + vx[voxel_idata(RNATM_RADCOEF_Ks, RNATM_SVX_OP_MIN)] = vx_ks_min; + vx[voxel_idata(RNATM_RADCOEF_Ks, RNATM_SVX_OP_MAX)] = vx_ks_max; + vx[voxel_idata(RNATM_RADCOEF_Kext, RNATM_SVX_OP_MIN)] = vx_kext_min; + vx[voxel_idata(RNATM_RADCOEF_Kext, RNATM_SVX_OP_MAX)] = vx_kext_max; + #undef REDUX_MIN #undef REDUX_MAX @@ -315,7 +294,8 @@ voxelize_gas intersect = suvm_polyhedron_intersect_aabb(&poly, vx_low, vx_upp); if(intersect == SUVM_INTERSECT_NONE) continue; - res = update_voxel(atm, &tetra, part, mcode[0]); + /* TODO setup the band and quadrature point */ + res = update_voxel(atm, &tetra, part, 0/*iband*/, 0/*iquad*/, mcode[0]); if(res != RES_OK) goto error; } } @@ -511,7 +491,7 @@ vx_get(const size_t xyz[3], const uint64_t mcode, void* dst, void* context) } } - vx = partition_get_voxel(ctx->part, 0, 0, ivx); + vx = partition_get_voxel(ctx->part, ivx); memcpy(dst, vx, NFLOATS_PER_VOXEL * sizeof(float)); } @@ -660,9 +640,6 @@ create_octrees(struct rnatm* atm, const struct rnatm_create_args* args) /* Create the vortex partition pool */ pool_args.npartitions = NPARTITIONS_PER_THREAD * atm->nthreads; - pool_args.nbands = sck_get_bands_count(atm->gas.ck); - pool_args.nquad_pts = band_get_quad_pt_count; - pool_args.context = atm->gas.ck; pool_args.partition_definition = 8; pool_args.allocator = atm->allocator; res = pool_create(&pool_args, &pool); diff --git a/src/rnatm_voxel_partition.c b/src/rnatm_voxel_partition.c @@ -28,26 +28,14 @@ #include <rsys/ref_count.h> struct partition { - size_t definition; - - /* Size of a cluster in bytes. A cluster is a list of voxels for a quadrature - * point of a spectral band */ - size_t cluster_size; - - size_t cluster_nvoxels; /* #voxels in a cluster */ + size_t definition; /* #voxels along the 3 dimensions */ struct list_node node; size_t id; /* Unique identifier of the partition */ - /* Offset in #float toward the voxels of a quadrature point in a spectral - * band. The first dimension corresponds to the spectral bands and the second - * dimension to the quadrature points in the spectral band */ - struct darray_size_t_list cluster_offsets; - - /* List of voxels. Use the LUT cluster_offsets to retrieve the list of voxels - * from a specific quadrature point in a given spectral band. The voxels are - * then sorted according to the morton code of the voxel coordinates */ + /* List of voxels sorted according to the morton code of the voxel coords */ float* voxels; + size_t nvoxels; struct mem_allocator* allocator; ref_T ref; @@ -82,90 +70,21 @@ partition_ref_put static INLINE res_T check_pool_create_args(const struct pool_create_args* args) { - size_t iband; - if(!args || !args->npartitions - || !args->nbands - || !args->nquad_pts || !IS_POW2(args->partition_definition)) { return RES_BAD_ARG; } - FOR_EACH(iband, 0, args->nbands) { - const size_t nquad_pts = args->nquad_pts(iband, args->context); - if(!nquad_pts) return RES_BAD_ARG; - } - return RES_OK; } -static res_T -setup_partition_clusters - (struct partition* partition, - const struct pool_create_args* args) -{ - size_t nclusters = 0; - size_t iband = 0; - size_t cluster_nfloats = 0; - res_T res = RES_OK; - ASSERT(partition && args); - - /* Compute the size of a cluster. Ensure a multiple of 64 bytes, i.e. align - * the cluster of voxels on the size of a cache line */ - partition->definition = args->partition_definition; - partition->cluster_nvoxels = - args->partition_definition - * args->partition_definition - * args->partition_definition; - partition->cluster_size = ALIGN_SIZE - (partition->cluster_nvoxels * NFLOATS_PER_VOXEL * sizeof(float), 64u); - cluster_nfloats = partition->cluster_size / sizeof(float); - - /* Allocate the list of per band cluster offsets */ - res = darray_size_t_list_resize(&partition->cluster_offsets, args->nbands); - if(res != RES_OK) goto error; - - FOR_EACH(iband, 0, args->nbands) { - struct darray_size_t* offsets = NULL; - const size_t nquad_pts = args->nquad_pts(iband, args->context); - size_t iquad_pt = 0; - - /* Allocate the cluster offsets for the quadrature points of the band */ - offsets = darray_size_t_list_data_get(&partition->cluster_offsets) + iband; - res = darray_size_t_resize(offsets, nquad_pts); - if(res != RES_OK) goto error; - - /* Setup the cluster offsets for the quadrature points of the band */ - FOR_EACH(iquad_pt, 0, nquad_pts) { - darray_size_t_data_get(offsets)[iquad_pt] = nclusters * cluster_nfloats; - nclusters += 1; - } - } - - /* Allocate the whole set of clusters */ - partition->voxels = MEM_ALLOC_ALIGNED - (partition->allocator, nclusters*partition->cluster_size, 64); - if(!partition->voxels) { res = RES_MEM_ERR; goto error; } - -exit: - return res; -error: - darray_size_t_list_clear(&partition->cluster_offsets); - if(partition->voxels) { - MEM_RM(partition->allocator, partition->voxels); - partition->voxels = NULL; - } - goto exit; -} - static void release_partition(ref_T* ref) { struct partition* partition = CONTAINER_OF(ref, struct partition, ref); ASSERT(partition && is_list_empty(&partition->node)); if(partition->voxels) MEM_RM(partition->allocator, partition->voxels); - darray_size_t_list_release(&partition->cluster_offsets); MEM_RM(partition->allocator, partition); } @@ -208,6 +127,7 @@ partition_create { struct partition* partition = NULL; struct mem_allocator* allocator = NULL; + size_t partsz = 0; res_T res = RES_OK; ASSERT(check_pool_create_args(args) == RES_OK && out_partition); @@ -221,10 +141,17 @@ partition_create ref_init(&partition->ref); list_init(&partition->node); partition->id = SIZE_MAX; - darray_size_t_list_init(partition->allocator, &partition->cluster_offsets); - - res = setup_partition_clusters(partition, args); - if(res != RES_OK) goto error; + partition->definition = args->partition_definition; + partition->nvoxels = + partition->definition + * partition->definition + * partition->definition; + + /* Allocate the partition voxels. Align the voxels on the size of a cache + * line (i.e. 64 bytes) */ + partsz = partition->nvoxels * NFLOATS_PER_VOXEL * sizeof(float); + partition->voxels = MEM_ALLOC_ALIGNED(partition->allocator, partsz, 64); + if(!partition->voxels) { res = RES_MEM_ERR; goto error; } exit: *out_partition = partition; @@ -263,52 +190,21 @@ partition_get_definition(const struct partition* partition) } float* -partition_get_voxel - (struct partition* part, - const size_t iband, - const size_t iquad_pt, - const size_t ivoxel) +partition_get_voxel(struct partition* part, const size_t ivoxel) { - const struct darray_size_t* band_offsets = NULL; - size_t offset = 0; - ASSERT(part); - ASSERT(iband < darray_size_t_list_size_get(&part->cluster_offsets)); - ASSERT(ivoxel < part->cluster_nvoxels); - - /* Fetch the offset toward the list of voxels of the quadrature point of the - * given spectral band */ - band_offsets = darray_size_t_list_cdata_get(&part->cluster_offsets)+iband; - ASSERT(iquad_pt < darray_size_t_size_get(band_offsets)); - offset = darray_size_t_cdata_get(band_offsets)[iquad_pt]; - - return part->voxels + offset + (ivoxel*NFLOATS_PER_VOXEL); + ASSERT(part && ivoxel < part->nvoxels); + return part->voxels + (ivoxel*NFLOATS_PER_VOXEL); } void partition_clear_voxels(struct partition* partition) { - size_t iband, nbands; + size_t ivoxel = 0; ASSERT(partition); - nbands = darray_size_t_list_size_get(&partition->cluster_offsets); - nbands = 1; /* FIXME For tests */ - - FOR_EACH(iband, 0, nbands) { - const struct darray_size_t* band_offsets; - size_t iquad_pt, nquad_pts; - - band_offsets = darray_size_t_list_cdata_get(&partition->cluster_offsets)+iband; - nquad_pts = darray_size_t_size_get(band_offsets); - nquad_pts = 1; /* FIXME For tests */ - - FOR_EACH(iquad_pt, 0, nquad_pts) { - size_t ivoxel = 0; - - FOR_EACH(ivoxel, 0, partition->cluster_nvoxels) { - float* voxel = partition_get_voxel(partition, iband, iquad_pt, ivoxel); - voxel_clear(voxel); - } - } + FOR_EACH(ivoxel, 0, partition->nvoxels) { + float* voxel = partition_get_voxel(partition, ivoxel); + voxel_clear(voxel); } } diff --git a/src/rnatm_voxel_partition.h b/src/rnatm_voxel_partition.h @@ -28,15 +28,11 @@ struct pool_create_args { size_t npartitions; /* Number of partitions to preallocate */ - - size_t nbands; /* Number of spectral band per partition */ - size_t (*nquad_pts)(const size_t iband, void* ctx); /* #quad pts per band */ - void* context; /* User data sent as the last argument of the functor */ size_t partition_definition; /* #voxels along XYZ. Must be a power of 2 */ struct mem_allocator* allocator; /* NULL <=> default allocator */ }; -#define POOL_CREATE_ARGS_DEFAULT__ {0, 0, NULL, NULL, 32, NULL} +#define POOL_CREATE_ARGS_DEFAULT__ {0, 32, NULL} static const struct pool_create_args POOL_CREATE_ARGS_DEFAULT = POOL_CREATE_ARGS_DEFAULT__; @@ -56,8 +52,6 @@ partition_get_definition extern LOCAL_SYM float* partition_get_voxel (struct partition* partition, - const size_t iband, - const size_t iquad_pt, const size_t ivoxel); extern LOCAL_SYM void