star-line

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

commit aa66e8720a2574ae7fdcfcb176ad36a13857b698
parent a8b999b349bd90dac3ff41fc156b8e636634931c
Author: Vincent Forest <vincent.forest@meso-star.com>
Date:   Thu, 12 Feb 2026 10:41:27 +0100

Fix breaks introduced by updates in the shtr API

The profiles of the [de]serialization functions have been updated.

Diffstat:
Msrc/sln_build.c | 18+++---------------
Msrc/sln_tree.c | 22+++++++++++++---------
2 files changed, 16 insertions(+), 24 deletions(-)

diff --git a/src/sln_build.c b/src/sln_build.c @@ -301,29 +301,17 @@ load_lines_hitran(struct cmd* cmd, const struct args* args) static res_T load_lines_shtr(struct cmd* cmd, const struct args* args) { - FILE* fp = NULL; + struct shtr_line_list_read_args rlines_args = SHTR_LINE_LIST_READ_ARGS_NULL; res_T res = RES_OK; ASSERT(cmd && args); - if(args->lines == NULL) { - fp = stdin; - } else { - fp = fopen(args->lines, "r"); - if(!fp) { - fprintf(stderr, "error opening file `%s' -- %s\n", - args->lines, strerror(errno)); - res = RES_IO_ERR; - goto error; - } - } - - res = shtr_line_list_create_from_stream(cmd->shtr, fp, &cmd->lines); + rlines_args.filename = args->lines; + res = shtr_line_list_read(cmd->shtr, &rlines_args, &cmd->lines); if(res != RES_OK) goto error; exit: if(cmd->lines) SHTR(line_list_ref_put(cmd->lines)); - if(fp && fp != stdin) CHK(fclose(fp) == 0); return res; error: goto exit; diff --git a/src/sln_tree.c b/src/sln_tree.c @@ -318,6 +318,7 @@ sln_tree_create_from_stream FILE* stream, struct sln_tree** out_tree) { + struct shtr_line_list_read_args rlines_args = SHTR_LINE_LIST_READ_ARGS_NULL; struct sln_tree* tree = NULL; size_t n = 0; int version = 0; @@ -372,7 +373,8 @@ sln_tree_create_from_stream res = shtr_isotope_metadata_create_from_stream(shtr, stream, &tree->args.metadata); if(res != RES_OK) goto error; - res = shtr_line_list_create_from_stream(shtr, stream, &tree->args.lines); + rlines_args.file = stream; + res = shtr_line_list_read(shtr, &rlines_args, &tree->args.lines); if(res != RES_OK) goto error; exit: @@ -534,7 +536,8 @@ sln_mesh_eval(const struct sln_mesh* mesh, const double wavenumber) res_T sln_tree_write(const struct sln_tree* tree, FILE* stream) { - size_t n; + struct shtr_line_list_write_args wlines_args = SHTR_LINE_LIST_WRITE_ARGS_NULL; + size_t nnodes, nverts; res_T res = RES_OK; if(!tree || !stream) { res = RES_BAD_ARG; goto error; } @@ -547,13 +550,13 @@ sln_tree_write(const struct sln_tree* tree, FILE* stream) } (void)0 WRITE(&SLN_TREE_VERSION, 1); - n = darray_node_size_get(&tree->nodes); - WRITE(&n, 1); - WRITE(darray_node_cdata_get(&tree->nodes), n); + nnodes = darray_node_size_get(&tree->nodes); + WRITE(&nnodes, 1); + WRITE(darray_node_cdata_get(&tree->nodes), nnodes); - n = darray_vertex_size_get(&tree->vertices); - WRITE(&n, 1); - WRITE(darray_vertex_cdata_get(&tree->vertices), n); + nverts = darray_vertex_size_get(&tree->vertices); + WRITE(&nverts, 1); + WRITE(darray_vertex_cdata_get(&tree->vertices), nverts); WRITE(&tree->args.line_profile, 1); WRITE(&tree->args.molecules, 1); @@ -567,7 +570,8 @@ sln_tree_write(const struct sln_tree* tree, FILE* stream) res = shtr_isotope_metadata_write(tree->args.metadata, stream); if(res != RES_OK) goto error; - res = shtr_line_list_write(tree->args.lines, stream); + wlines_args.file = stream; + res = shtr_line_list_write(tree->args.lines, &wlines_args); if(res != RES_OK) goto error; exit: