star-line

Structure for accelerating line importance sampling
git clone git://git.meso-star.fr/star-line.git
Log | Files | Refs | README | LICENSE

commit 7a7404001b5b7147e5d778b76e6f8e304195d501
parent c8b444517b909892eb65af5dd8a87b5161bf655d
Author: Vincent Forest <vincent.forest@meso-star.com>
Date:   Wed,  4 May 2022 16:48:57 +0200

Add the sln_node_get_mesh function

Diffstat:
Msrc/sln.h | 6++++++
Msrc/sln_tree.c | 13++++++++++++-
Msrc/test_sln_tree.c | 19+++++++++++++++++++
3 files changed, 37 insertions(+), 1 deletion(-)

diff --git a/src/sln.h b/src/sln.h @@ -205,6 +205,12 @@ sln_node_get_line const size_t iline, const struct shtr_line** line); +SLN_API res_T +sln_node_get_mesh + (const struct sln_tree* tree, + const struct sln_node* node, + struct sln_mesh* mesh); + /******************************************************************************* * Helper functions ******************************************************************************/ diff --git a/src/sln_tree.c b/src/sln_tree.c @@ -279,7 +279,6 @@ setup_molecule_params params[mol->id].concentration = mol->concentration; params[mol->id].cutoff = mol->cutoff; - FOR_EACH(iiso, 0, mol->nisotopes) { const struct sln_isotope* iso = &mol->isotopes[iiso]; ASSERT(iso->id_local < SLN_MAX_ISOTOPES_COUNT); @@ -486,3 +485,15 @@ sln_node_get_line return shtr_lines_view_get_line (tree->lines_view, node->range[0] + iline, line); } + +res_T +sln_node_get_mesh + (const struct sln_tree* tree, + const struct sln_node* node, + struct sln_mesh* mesh) +{ + if(!tree || !node || !mesh) return RES_BAD_ARG; + mesh->vertices = darray_vertex_cdata_get(&tree->vertices) + node->ivertices; + mesh->nvertices = node->nvertices; + return RES_OK; +} diff --git a/src/test_sln_tree.c b/src/test_sln_tree.c @@ -309,6 +309,18 @@ check_tree } static void +dump_line(FILE* stream, const struct sln_mesh* mesh) +{ + size_t i; + CHK(mesh); + FOR_EACH(i, 0, mesh->nvertices) { + fprintf(stream, "%g %g\n", + mesh->vertices[i].wavenumber, + mesh->vertices[i].ka); + } +} + +static void test_tree (struct sln_device* sln, struct shtr_isotope_metadata* metadata, @@ -316,6 +328,7 @@ test_tree { struct sln_tree_create_args tree_args = SLN_TREE_CREATE_ARGS_DEFAULT; struct sln_tree_desc desc = SLN_TREE_DESC_NULL; + struct sln_mesh mesh = SLN_MESH_NULL; struct sln_tree* tree = NULL; const struct sln_node* node = NULL; const struct shtr_line* line = NULL; @@ -368,6 +381,12 @@ test_tree CHK(sln_node_get_line(tree, node, 0, &line) == RES_OK); CHK(sln_node_get_line(tree, sln_tree_get_root(tree), 0, &line) == RES_OK); + CHK(sln_node_get_mesh(NULL, node, &mesh) == RES_BAD_ARG); + CHK(sln_node_get_mesh(tree, NULL, &mesh) == RES_BAD_ARG); + CHK(sln_node_get_mesh(tree, node, NULL) == RES_BAD_ARG); + CHK(sln_node_get_mesh(tree, node, &mesh) == RES_OK); + + dump_line(stdout, &mesh); check_tree(tree, lines_list); CHK(sln_tree_ref_get(NULL) == RES_BAD_ARG);