commit 434046e459c3944170c94904bf740175fa250299
parent 3f2af96d963e5122758ed73c650547c1042a8933
Author: Christophe Coustet <christophe.coustet@meso-star.com>
Date: Tue, 12 Mar 2019 15:59:23 +0100
Merge branch 'release_0.3.0'
Diffstat:
9 files changed, 233 insertions(+), 30 deletions(-)
diff --git a/README.md b/README.md
@@ -35,6 +35,12 @@ variable the install directories of its dependencies.
Release notes
-------------
+### Version 0.3
+
+- Add API calls to access to geometry frontiers.
+- Improve documentation in the header file.
+- BugFix: wrong data cleaning on computation canceling.
+
### Version 0.2.2
- BugFix when grouping components into enclosures.
diff --git a/cmake/CMakeLists.txt b/cmake/CMakeLists.txt
@@ -46,8 +46,8 @@ endif()
# Configure and define targets
################################################################################
set(VERSION_MAJOR 0)
-set(VERSION_MINOR 2)
-set(VERSION_PATCH 2)
+set(VERSION_MINOR 3)
+set(VERSION_PATCH 0)
set(VERSION ${VERSION_MAJOR}.${VERSION_MINOR}.${VERSION_PATCH})
set(SENC2D_FILES_SRC
diff --git a/src/senc2d.h b/src/senc2d.h
@@ -143,6 +143,7 @@ senc2d_device_ref_put
* StarEnclosures2D scene. A scene is a collection of segments. Each segment is
* defined with a medium on each side.
******************************************************************************/
+/* Creates an empty scene */
SENC2D_API res_T
senc2d_scene_create
(struct senc2d_device* device,
@@ -154,7 +155,8 @@ senc2d_scene_create
* Segments can be duplicates as long as they constantly define the same
* medium on both sides (or an error will be reported) and are deduplicated.
* When deduplicating segments, the first occurence is kept (with it original
- * global_id). */
+ * global_id). Users can provide their own global ids for segments; these ids
+ * are not used by the library but are returned as-is by some API calls. */
SENC2D_API res_T
senc2d_scene_add_geometry
(struct senc2d_scene* scene,
@@ -167,31 +169,39 @@ senc2d_scene_add_geometry
void(*position)(const unsigned ivert, double pos[2], void* context),
void* context);
+/* Returns a descriptor of the scene that holds the analysis' result. */
SENC2D_API res_T
senc2d_scene_analyze
(struct senc2d_scene* scene,
struct senc2d_descriptor** descriptor);
+/* Returns the convention flags in use with the scene. */
SENC2D_API res_T
senc2d_scene_get_convention
(const struct senc2d_scene* scene,
enum senc2d_convention* convention);
+/* Returns the number of segments in the scene. */
SENC2D_API res_T
senc2d_scene_get_segments_count
(const struct senc2d_scene* scene,
unsigned* count);
+/* Returns the number of unique segments in the scene (remaining
+ * segments after deduplication). */
SENC2D_API res_T
senc2d_scene_get_unique_segments_count
(const struct senc2d_scene* scene,
unsigned* count);
+/* Returns the number of vertices in the scene. */
SENC2D_API res_T
senc2d_scene_get_vertices_count
(const struct senc2d_scene* scene,
unsigned* count);
+/* Returns the number of unique vertices in the scene (remaining
+ * vertices after deduplication). */
SENC2D_API res_T
senc2d_scene_get_unique_vertices_count
(const struct senc2d_scene* scene,
@@ -208,28 +218,36 @@ senc2d_scene_ref_put
/*******************************************************************************
* StarEnclosures2D descriptor. It is an handle toward an analyze result.
******************************************************************************/
+/* Returns the greater medium id found in added geometry. In API calls using a
+ * medium, any value in the [0 max_medium_id[ range is valid. However there can
+ * be unused ids (no geometry refered to this medium id). */
SENC2D_API res_T
senc2d_descriptor_get_max_medium
(const struct senc2d_descriptor* descriptor,
unsigned* max_medium_id);
+/* Returns the number of enclosures. */
SENC2D_API res_T
senc2d_descriptor_get_enclosure_count
(const struct senc2d_descriptor* descriptor,
unsigned* count);
+/* Returns the number of enclosures that have some geometry refering to the
+ * imed_th medium. */
SENC2D_API res_T
senc2d_descriptor_get_enclosure_count_by_medium
(const struct senc2d_descriptor* descriptor,
- const unsigned imed, /* Must be in [0 max_medium_id] */
+ const unsigned imed,
unsigned* count);
+/* Returns the idx_th enclosure. */
SENC2D_API res_T
senc2d_descriptor_get_enclosure
(struct senc2d_descriptor* descriptor,
const unsigned idx,
struct senc2d_enclosure** enclosure);
+/* Returns the idx_th enclosure using the imed_th medium. */
SENC2D_API res_T
senc2d_descriptor_get_enclosure_by_medium
(struct senc2d_descriptor* descriptor,
@@ -237,46 +255,74 @@ senc2d_descriptor_get_enclosure_by_medium
const unsigned idx,
struct senc2d_enclosure** enclosure);
+/* Returns the number of unique segments (no duplicates here) in the whole
+ * geometry. */
SENC2D_API res_T
senc2d_descriptor_get_global_segments_count
(const struct senc2d_descriptor* descriptor,
- unsigned* count); /* Number of unique segments. */
+ unsigned* count);
+/* Returns the number of unique vertices (no duplicates here) in the whole
+ * geometry. */
SENC2D_API res_T
senc2d_descriptor_get_global_vertices_count
(const struct senc2d_descriptor* descriptor,
- unsigned* count); /* Number of unique vertices. */
+ unsigned* count);
+/* Returns the iseg_th global unique segment; the returned indices are global
+ * unique vertex indices. */
SENC2D_API res_T
senc2d_descriptor_get_global_segment
(const struct senc2d_descriptor* descriptor,
const unsigned iseg,
unsigned indices[2]);
+/* Returns the coordinates of the ivert_th global unique vertex. */
SENC2D_API res_T
senc2d_descriptor_get_global_vertex
(const struct senc2d_descriptor* descriptor,
const unsigned ivert,
double coord[2]);
+/* Returns the front and back media ids of the iseg_th global unique segment. */
SENC2D_API res_T
senc2d_descriptor_get_global_segment_media
(const struct senc2d_descriptor* descriptor,
const unsigned iseg,
unsigned media[2]);
+/* Returns the enclosures the iseg_th global unique segment front and back
+ * sides are member of. */
SENC2D_API res_T
senc2d_descriptor_get_global_segment_enclosures
(const struct senc2d_descriptor* descriptor,
const unsigned iseg,
unsigned enclosures[2]);
+/* Returns the global id of the iseg_th global unique segment, either the user
+ * provided one or the default one. */
SENC2D_API res_T
senc2d_descriptor_get_global_segment_global_id
(const struct senc2d_descriptor* descriptor,
const unsigned iseg,
unsigned* gid);
+/* Returns the number of vertices that are frontier vertices:
+ * - that have arity 1 (single segment using the vertex)
+ * - that connect 2 different media */
+SENC2D_API res_T
+senc2d_descriptor_get_frontier_vertices_count
+ (const struct senc2d_descriptor* descriptor,
+ unsigned* count);
+
+/* Returns the iver_th frontier vertex; the returned index is an global unique
+ * vertex index. */
+SENC2D_API res_T
+senc2d_descriptor_get_frontier_vertex
+ (const struct senc2d_descriptor* descriptor,
+ const unsigned iver,
+ unsigned* vrtx_id);
+
SENC2D_API res_T
senc2d_descriptor_ref_get
(struct senc2d_descriptor* descriptor);
@@ -292,39 +338,49 @@ senc2d_descriptor_ref_put
* An enclosure can list the "same" segment twice if both sides are in. In this
* case the 2 occurences of the segment have reversed vertices order and
* unique_segment_count and segment_count differ.
+ * Vertices and segments numbering schemes are specific to each enclosure:
+ * the "same" item appearing in 2 different enclosures has no reason to get the
+ * same index twice or to have the same index in the global numbering scheme.
* By-index API accesses of segments (or properties) visit unique segments
- * for indexes in the [0 unique_segment_count[ range and back-faces of the
- * doubly-listed segments in the [unique_segment_count segment_count[ range.
+ * for indices in the [0 unique_segment_count[ range and back-faces of the
+ * doubly-included segments in the [unique_segment_count segment_count[ range.
******************************************************************************/
+/* Returns the header of an enclosure. */
SENC2D_API res_T
senc2d_enclosure_get_header
(const struct senc2d_enclosure* enclosure,
struct senc2d_enclosure_header* header);
+/* Returns the iseg_th segment of an enclosure. */
SENC2D_API res_T
senc2d_enclosure_get_segment
(const struct senc2d_enclosure* enclosure,
const unsigned iseg,
unsigned indices[2]);
+/* Returns the coordinates of the ivert_th vertex of an enclosure. */
SENC2D_API res_T
senc2d_enclosure_get_vertex
(const struct senc2d_enclosure* enclosure,
const unsigned ivert,
double coord[2]);
+/* Returns the front and back side media ids of the iseg_th segment of an
+ * enclosure. */
SENC2D_API res_T
senc2d_enclosure_get_segment_media
(const struct senc2d_enclosure* enclosure,
const unsigned iseg,
unsigned medium[2]);
+/* Returns the global id of the iseg_th segment of an enclosure. */
SENC2D_API res_T
senc2d_enclosure_get_segment_global_id
(const struct senc2d_enclosure* enclosure,
const unsigned iseg,
unsigned* gid);
+/* Returns the id of the imed_th medium of an enclosure. */
SENC2D_API res_T
senc2d_enclosure_get_medium
(const struct senc2d_enclosure* enclosure,
diff --git a/src/senc2d_descriptor.c b/src/senc2d_descriptor.c
@@ -37,6 +37,7 @@ descriptor_release(ref_T * ref)
darray_segment_enc_release(&desc->segments_enc);
darray_enclosure_release(&desc->enclosures);
darray_enc_ids_array_release(&desc->enc_ids_array_by_medium);
+ darray_vrtx_id_release(&desc->frontiers);
MEM_RM(scn->dev->allocator, desc);
SENC2D(scene_ref_put(scn));
@@ -61,6 +62,7 @@ descriptor_create(struct senc2d_scene* scn)
darray_enc_ids_array_init(scn->dev->allocator,
&desc->enc_ids_array_by_medium);
OK(darray_enc_ids_array_resize(&desc->enc_ids_array_by_medium, scn->nmeds));
+ darray_vrtx_id_init(scn->dev->allocator, &desc->frontiers);
/* Enclosure 0 is always defined for infinite */
OK(darray_enclosure_resize(&desc->enclosures, 1));
desc->enclosures_count = 1;
@@ -276,6 +278,35 @@ senc2d_descriptor_get_global_segment_global_id
}
res_T
+senc2d_descriptor_get_frontier_vertices_count
+ (const struct senc2d_descriptor* desc,
+ unsigned* count)
+{
+ size_t tmp;
+ if(!desc || !count)
+ return RES_BAD_ARG;
+ tmp = darray_vrtx_id_size_get(&desc->frontiers);
+ ASSERT(tmp < UINT_MAX);
+ *count = (unsigned)tmp;
+ return RES_OK;
+}
+
+res_T
+senc2d_descriptor_get_frontier_vertex
+ (const struct senc2d_descriptor* desc,
+ const unsigned iver,
+ unsigned* vrtx_id)
+{
+ vrtx_id_t vrtx;
+ if(!vrtx_id || !desc
+ || iver >= darray_vrtx_id_size_get(&desc->frontiers))
+ return RES_BAD_ARG;
+ vrtx = darray_vrtx_id_cdata_get(&desc->frontiers)[iver];
+ *vrtx_id = (unsigned)vrtx; /* Back to API type */
+ return RES_OK;
+}
+
+res_T
senc2d_descriptor_ref_get(struct senc2d_descriptor* desc)
{
if(!desc) return RES_BAD_ARG;
diff --git a/src/senc2d_descriptor_c.h b/src/senc2d_descriptor_c.h
@@ -94,6 +94,8 @@ struct senc2d_descriptor {
struct darray_enc_ids_array enc_ids_array_by_medium;
seg_id_t segment_count;
vrtx_id_t vertices_count;
+ /* Store frontiers */
+ struct darray_vrtx_id frontiers;
ref_T ref;
};
diff --git a/src/senc2d_scene_analyze.c b/src/senc2d_scene_analyze.c
@@ -141,7 +141,7 @@ extract_connex_components
struct s2d_scene_view** s2d_view,
ATOMIC* component_count,
/* Shared error status.
- * We accept to overwritte an error with a different error */
+ * We accept to overwrite an error with a different error */
res_T* p_res)
{
/* This function is called from an omp parallel block and executed
@@ -322,6 +322,19 @@ extract_connex_components
* associated with neighbour's medium). */
component_canceled = 1;
darray_side_id_clear(&stack);
+ /* Reverse the used flag for sides in cancelled component */
+ sz = darray_side_id_size_get(¤t_component);
+ FOR_EACH(ii, 0, sz) {
+ side_id_t used_side
+ = darray_side_id_cdata_get(¤t_component)[ii];
+ seg_id_t used_trg_id = SEGSIDE_2_SEG(used_side);
+ enum side_flag used_side_flag
+ = SEGSIDE_2_SIDEFLAG(used_side);
+ unsigned char* used = processed + used_trg_id;
+ ASSERT(*used & (unsigned char)used_side_flag);
+ *used &= 0xFF ^ (unsigned char)used_side_flag;
+ }
+
goto canceled;
}
/* Mark neighbour as processed and stack it */
@@ -363,7 +376,7 @@ canceled:
/* Write component membership in the global structure
* No need for sync here as an unique thread write a given side */
- STATIC_ASSERT(sizeof(cc->cc_id) >= 4, Cannot_write_IDs_sync_free);
+ {STATIC_ASSERT(sizeof(cc->cc_id) >= 4, Cannot_write_IDs_sync_free);}
ASSERT(IS_ALIGNED(segments_comp->component, sizeof(cc->cc_id)));
FOR_EACH(ii, 0, sz) {
const side_id_t s = darray_side_id_cdata_get(¤t_component)[ii];
@@ -376,6 +389,7 @@ canceled:
/* Compute the normal at the max_y vertex. */
max_ny = 0;
sz = darray_side_id_size_get(&ids_of_sides_around_max_y_vertex);
+ ASSERT(sz > 0);
FOR_EACH(ii, 0, sz) {
const side_id_t side_id =
darray_side_id_cdata_get(&ids_of_sides_around_max_y_vertex)[ii];
@@ -518,7 +532,7 @@ group_connex_components
struct s2d_scene_view* s2d_view,
ATOMIC* next_enclosure_id,
/* Shared error status.
- * We accept to overwritte an error with a different error */
+ * We accept to overwrite an error with a different error */
res_T* res)
{
/* This function is called from an omp parallel block and executed
@@ -589,8 +603,9 @@ group_connex_components
enum side_id hit_side =
((hit.normal[1] < 0) /* Facing geometrical normal of hit */
== ((desc->scene->convention & SENC2D_CONVENTION_NORMAL_FRONT) != 0))
- /* Warning: Embree 2 convention for geometrical normals is
- * left-handed and star-enclosure uses right-handed convention */
+ /* Warning: following Embree 2 convention for geometrical normals,
+ * the Star2D hit normal is left-handed while star-enclosure2D uses
+ * right-handed convention */
? SIDE_BACK : SIDE_FRONT;
ASSERT(hit.normal[1] != 0);
ASSERT(hit_seg_id < desc->scene->nusegs);
@@ -654,8 +669,9 @@ collect_and_link_neighbours
struct segside* segsides,
struct darray_segment_tmp* segments_tmp_array,
struct darray_neighbourhood* neighbourhood_by_vertex,
+ struct darray_vrtx_id* frontiers,
/* Shared error status.
- * We accept to overwritte an error with a different error */
+ * We accept to overwrite an error with a different error */
res_T* res)
{
/* This function is called from an omp parallel block and executed
@@ -671,7 +687,7 @@ collect_and_link_neighbours
seg_id_t s;
ASSERT(scn && segsides && segments_tmp_array
- && neighbourhood_by_vertex && res);
+ && neighbourhood_by_vertex && frontiers && res);
ASSERT((size_t)scn->nuverts + (size_t)scn->nusegs + 2 <= EDGE_MAX__);
segments_in = darray_segment_in_cdata_get(&scn->segments_in);
@@ -832,17 +848,25 @@ collect_and_link_neighbours
previous = darray_neighbour_cdata_get(neighbourhood) + i - 1;
prev_id = previous->seg_id;
log_err(scn->dev,
- "%s: found 2 overlying segments (%u & %u).\n",
- FUNC_NAME, prev_id, crt_id);
+ "%s: found 2 overlying segments (%lu & %lu).\n", FUNC_NAME,
+ (unsigned long)segments_in[crt_id].global_id,
+ (unsigned long)segments_in[prev_id].global_id);
+
/* TODO */
*res = RES_BAD_OP;
return;
}
a = current->angle;
/* Link sides */
+ ASSERT(p_crt_side->facing_side_id[crt_end] == SIDE_NULL__);
+ ASSERT(p_ccw_side->facing_side_id[ccw_end] == SIDE_NULL__);
p_crt_side->facing_side_id[crt_end] = ccw_side_idx;
p_ccw_side->facing_side_id[ccw_end] = crt_side_idx;
/* Record media */
+ ASSERT(p_crt_side->medium == MEDIUM_NULL__
+ || p_crt_side->medium == segments_in[crt_id].medium[crt_side]);
+ ASSERT(p_ccw_side->medium == MEDIUM_NULL__
+ || p_ccw_side->medium == segments_in[ccw_id].medium[ccw_side]);
p_crt_side->medium = segments_in[crt_id].medium[crt_side];
p_ccw_side->medium = segments_in[ccw_id].medium[ccw_side];
ASSERT(p_crt_side->medium < scn->nmeds);
@@ -852,10 +876,12 @@ collect_and_link_neighbours
* - different media on its sides */
if(neighbour_count == 1
&& p_crt_side->medium != p_ccw_side->medium)
+#pragma omp critical
{
log_warn(scn->dev,
- "%s: found possible hole involving segments %u.\n",
- FUNC_NAME, crt_id);
+ "%s: found possible frontier involving segment %lu.\n",
+ FUNC_NAME, (unsigned long)segments_in[crt_id].global_id);
+ darray_vrtx_id_push_back(frontiers, &v);
}
}
}
@@ -867,8 +893,9 @@ build_result
(struct senc2d_descriptor* desc,
const struct darray_ptr_component_descriptor* connex_components,
const struct darray_segment_comp* segments_comp_array,
+ struct darray_vrtx_id* frontiers,
/* Shared error status.
- * We accept to overwritte an error with a different error */
+ * We accept to overwrite an error with a different error */
res_T* res)
{
/* This function is called from an omp parallel block and executed
@@ -885,7 +912,7 @@ build_result
int64_t sg;
int64_t ee;
- ASSERT(desc && connex_components && segments_comp_array && res);
+ ASSERT(desc && connex_components && segments_comp_array && frontiers && res);
alloc = descriptor_get_allocator(desc);
scn = desc->scene;
@@ -893,7 +920,8 @@ build_result
normals_front = (scn->convention & SENC2D_CONVENTION_NORMAL_FRONT) != 0;
normals_back = (scn->convention & SENC2D_CONVENTION_NORMAL_BACK) != 0;
ASSERT(normals_back != normals_front);
- ASSERT(output_normal_in != ((scn->convention & SENC2D_CONVENTION_NORMAL_OUTSIDE) != 0));
+ ASSERT(output_normal_in !=
+ ((scn->convention & SENC2D_CONVENTION_NORMAL_OUTSIDE) != 0));
ASSERT(darray_ptr_component_descriptor_size_get(connex_components)
<= COMPONENT_MAX__);
cc_descriptors = darray_ptr_component_descriptor_cdata_get(connex_components);
@@ -1058,6 +1086,10 @@ build_result
*res = tmp_res;
} /* No barrier here */
htable_vrtx_id_release(&vtable);
+ /* The first thread here copies frontiers into descriptor */
+#pragma omp single nowait
+ darray_vrtx_id_copy_and_clear(&desc->frontiers, frontiers);
+ /* No barrier here */
}
/*******************************************************************************
@@ -1076,6 +1108,9 @@ senc2d_scene_analyze
* They are refered to by arrays of ids. */
struct darray_ptr_component_descriptor connex_components;
char connex_components_initialized = 0;
+ /* Array of frontiers vertices */
+ struct darray_vrtx_id frontiers;
+ char frontiers_initialized = 0;
/* Store by-segment components */
struct darray_segment_comp segments_comp;
char segments_comp_initialized = 0;
@@ -1104,6 +1139,8 @@ senc2d_scene_analyze
darray_segment_tmp_init(scn->dev->allocator, &segments_tmp);
segments_tmp_initialized = 1;
+ darray_vrtx_id_init(scn->dev->allocator, &frontiers);
+ frontiers_initialized = 1;
OK(darray_segment_tmp_resize(&segments_tmp, scn->nusegs));
segsides
@@ -1112,6 +1149,15 @@ senc2d_scene_analyze
res = RES_MEM_ERR;
goto error;
}
+#ifndef NDEBUG
+ else {
+ /* Initialise trgsides to allow assert code */
+ size_t i;
+ FOR_EACH(i, 0, 2 * scn->nusegs)
+ init_segside(scn->dev->allocator, segsides + i);
+ }
+#endif
+
/* We create a neighbourhood for every single unique vertex,
* regardless the fact it is used by a segment.
@@ -1126,7 +1172,7 @@ senc2d_scene_analyze
{
/* Step 1: build neighbourhoods */
collect_and_link_neighbours(scn, segsides, &segments_tmp,
- &neighbourhood_by_vertex, &res);
+ &neighbourhood_by_vertex, &frontiers, &res);
/* No barrier at the end of step 1: data used in step 1 cannot be
* released / data produced by step 1 cannot be used
* until next sync point */
@@ -1223,7 +1269,8 @@ senc2d_scene_analyze
} /* No barrier here */
/* Step 4: Build result */
- build_result(desc, &connex_components, &segments_comp, &res);
+ build_result(desc, &connex_components, &segments_comp, &frontiers,
+ &res);
/* No barrier at the end of step 4: data used in step 4 cannot be
* released / data produced by step 4 cannot be used
* until next sync point */
@@ -1266,6 +1313,7 @@ exit:
darray_neighbourhood_release(&neighbourhood_by_vertex);
if(segments_tmp_initialized) darray_segment_tmp_release(&segments_tmp);
if(segments_comp_initialized) darray_segment_comp_release(&segments_comp);
+ if(frontiers_initialized) darray_vrtx_id_release(&frontiers);
if(segsides) MEM_RM(scn->dev->allocator, segsides);
if(desc) *out_desc = desc;
diff --git a/src/senc2d_scene_analyze_c.h b/src/senc2d_scene_analyze_c.h
@@ -37,6 +37,15 @@ struct segside {
* front(seg_0), back(seg_0), front(seg_1), back(seg_1), ... */
};
+static FINLINE void
+init_segside(struct mem_allocator* alloc, struct segside* data)
+{
+ int i;
+ ASSERT(data); (void)alloc;
+ FOR_EACH(i, 0, 3) data->facing_side_id[i] = SIDE_NULL__;
+ data->medium = MEDIUM_NULL__;
+}
+
#define DARRAY_NAME side_id
#define DARRAY_DATA side_id_t
#include <rsys/dynamic_array.h>
diff --git a/src/test_senc2d_descriptor.c b/src/test_senc2d_descriptor.c
@@ -157,16 +157,27 @@ main(int argc, char** argv)
CHK(senc2d_descriptor_get_global_segment_media(desc, 0, media) == RES_OK);
CHK(media[0] == ctx.front_media[0] && media[1] == ctx.back_media[1]);
- CHK(senc2d_descriptor_get_global_segment_enclosures(
- NULL, 0, enclosures) == RES_BAD_ARG);
- CHK(senc2d_descriptor_get_global_segment_enclosures(
- NULL, nvertices, enclosures) == RES_BAD_ARG);
+ CHK(senc2d_descriptor_get_global_segment_enclosures(NULL, 0, enclosures)
+ == RES_BAD_ARG);
+ CHK(senc2d_descriptor_get_global_segment_enclosures(NULL, nvertices, enclosures)
+ == RES_BAD_ARG);
CHK(senc2d_descriptor_get_global_segment_enclosures(desc, 0, NULL)
== RES_BAD_ARG);
CHK(senc2d_descriptor_get_global_segment_enclosures(desc, 0, enclosures)
== RES_OK);
CHK(enclosures[0] == 0 && enclosures[1] == 1);
+ CHK(senc2d_descriptor_get_global_segment_global_id(NULL, 0, indices)
+ == RES_BAD_ARG);
+ CHK(senc2d_descriptor_get_global_segment_global_id(NULL, nvertices, indices)
+ == RES_BAD_ARG);
+ CHK(senc2d_descriptor_get_global_segment_global_id(desc, 0, NULL)
+ == RES_BAD_ARG);
+ CHK(senc2d_descriptor_get_global_segment_global_id(desc, 0, indices)
+ == RES_OK);
+ /* No duplicates and no custom id: user id is unique vertex id */
+ CHK(indices[0] == 0);
+
/* Add valid duplicate geometry */
CHK(senc2d_descriptor_ref_put(desc) == RES_OK);
desc = NULL;
@@ -191,10 +202,50 @@ main(int argc, char** argv)
nvertices, get_position, &ctx) == RES_BAD_ARG);
CHK(senc2d_scene_ref_put(scn) == RES_OK);
- CHK(senc2d_device_ref_put(dev) == RES_OK);
if(desc) CHK(senc2d_descriptor_ref_put(desc) == RES_OK);
+ desc = NULL;
CHK(senc2d_enclosure_ref_put(enc) == RES_OK);
+ /* Same square with a hole (last segment is missing) */
+ CHK(senc2d_scene_create(dev,
+ SENC2D_CONVENTION_NORMAL_FRONT | SENC2D_CONVENTION_NORMAL_INSIDE, &scn)
+ == RES_OK);
+
+ CHK(senc2d_scene_add_geometry(scn, nsegments - 1, get_indices, get_media,
+ NULL, nvertices, get_position, &ctx) == RES_OK);
+
+ CHK(senc2d_scene_analyze(scn, &desc) == RES_OK);
+
+ CHK(senc2d_descriptor_get_frontier_vertices_count(NULL, NULL)
+ == RES_BAD_ARG);
+ CHK(senc2d_descriptor_get_frontier_vertices_count(desc, NULL)
+ == RES_BAD_ARG);
+ CHK(senc2d_descriptor_get_frontier_vertices_count(NULL, &count)
+ == RES_BAD_ARG);
+ CHK(senc2d_descriptor_get_frontier_vertices_count(desc, &count)
+ == RES_OK);
+
+ CHK(senc2d_descriptor_get_frontier_vertex(NULL, count, NULL)
+ == RES_BAD_ARG);
+ CHK(senc2d_descriptor_get_frontier_vertex(desc, count, NULL)
+ == RES_BAD_ARG);
+ CHK(senc2d_descriptor_get_frontier_vertex(NULL, 0, NULL)
+ == RES_BAD_ARG);
+ CHK(senc2d_descriptor_get_frontier_vertex(NULL, count, indices)
+ == RES_BAD_ARG);
+ CHK(senc2d_descriptor_get_frontier_vertex(desc, 0, NULL)
+ == RES_BAD_ARG);
+ CHK(senc2d_descriptor_get_frontier_vertex(desc, count, indices)
+ == RES_BAD_ARG);
+ CHK(senc2d_descriptor_get_frontier_vertex(NULL, 0, indices)
+ == RES_BAD_ARG);
+ CHK(senc2d_descriptor_get_frontier_vertex(desc, 0, indices)
+ == RES_OK);
+
+ CHK(senc2d_scene_ref_put(scn) == RES_OK);
+ CHK(senc2d_device_ref_put(dev) == RES_OK);
+ if(desc) CHK(senc2d_descriptor_ref_put(desc) == RES_OK);
+
check_memory_allocator(&allocator);
mem_shutdown_proxy_allocator(&allocator);
CHK(mem_allocated_size() == 0);
diff --git a/src/test_senc2d_scene.c b/src/test_senc2d_scene.c
@@ -215,7 +215,7 @@ main(int argc, char** argv)
CHK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, NULL,
nvertices, get_position, &ctx) == RES_BAD_ARG);
- /* It is OK dd geometry after a failed add */
+ /* It is OK add geometry after a failed add */
ctx.back_media = medium1;
CHK(senc2d_scene_add_geometry(scn, nsegments, get_indices, get_media, NULL,
nvertices, get_position, &ctx) == RES_OK);