commit a4ae874974154ec141a123a6196cfc019b672fbd
parent 670c81368805ac57807b5305dae9c83111e6fc94
Author: Vincent Forest <vincent.forest@meso-star.com>
Date: Tue, 10 Mar 2026 14:35:12 +0100
Add support for lines in shtr format to the sln-get tool
This binary is more compact and faster to load.
Diffstat:
2 files changed, 64 insertions(+), 22 deletions(-)
diff --git a/doc/sln-get.1 b/doc/sln-get.1
@@ -25,7 +25,7 @@
.\""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""
.Sh SYNOPSIS
.Nm
-.Op Fl hlmnrv
+.Op Fl hlmnrsv
.Op Fl L Ar nlevels
.Op Fl R Ar nlevels
.Op Fl w Ar wavenumber
@@ -68,7 +68,11 @@ Display short help and exit.
.\""""""""""""""""""""""""""""""""""
.It Fl i Ar lines
List of lines from which the tree was built.
-This list is in HITRAN format.
+This list is in binary format as generated by the
+.Xr shtr 1
+binary, or in plain text HITRAN format, depending on whether the
+.Fl s
+option is set or not, respectively.
.\""""""""""""""""""""""""""""""""""
.It Fl L Ar nlevels
Descend the tree by visiting the left child of the nodes traversed up to
@@ -87,6 +91,19 @@ The current node then becomes this left child.
If the node has no left child, i.e., if it is a leaf, then the current
node is not updated.
.\""""""""""""""""""""""""""""""""""
+.It Fl m
+Prints the polyline of the current node, i.e. the mesh representing the high
+resolution spectrum of all the lines it structures.
+.Pp
+The output data is a list of polyline vertices, in plain text, where each
+line represents the two values of a polyline vertex, separated by a space:
+its wavenumber in cm^-1, and its associated spectrum value.
+.\""""""""""""""""""""""""""""""""""
+.It Fl n
+Displays the description of the current node, i.e., its level relative to
+the root, the number of lines it partitions, and the number of mesh
+vertices used to represent them at the node's hierarchical level.
+.\""""""""""""""""""""""""""""""""""
.It Fl p Ar molparams
Isotopologue metadata from which the tree was built.
The data is in HITRAN format.
@@ -108,18 +125,12 @@ This option behaves in the same way as the
option, but here for the right child of the node, instead of the left
child.
.\""""""""""""""""""""""""""""""""""
-.It Fl m
-Prints the polyline of the current node, i.e. the mesh representing the high
-resolution spectrum of all the lines it structures.
-.Pp
-The output data is a list of polyline vertices, in plain text, where each
-line represents the two values of a polyline vertex, separated by a space:
-its wavenumber in cm^-1, and its associated spectrum value.
-.\""""""""""""""""""""""""""""""""""
-.It Fl n
-Displays the description of the current node, i.e., its level relative to
-the root, the number of lines it partitions, and the number of mesh
-vertices used to represent them at the node's hierarchical level.
+.It Fl s
+Specifies that input lines are formatted according to the binary format
+as written by the
+.Xr shtr 1
+utility, and not according to the HITRAN format.
+This format is more compact, allowing for faster loading of line data.
.\""""""""""""""""""""""""""""""""""
.It Fl v
Make
@@ -196,6 +207,7 @@ sln-get -i lines.par -p molparams.txt -L3 -w 50 tree.sln
.Ed
.\""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""
.Sh SEE ALSO
+.Xr shtr 1 ,
.Xr sln-build 1
.Rs
.%T The HITRAN Database
diff --git a/src/sln_get.c b/src/sln_get.c
@@ -70,8 +70,9 @@ struct args {
/* Miscellaneous */
int quit;
int verbose;
+ int lines_in_shtr_format;
};
-#define ARGS_DEFAULT__ {NULL,NULL,NULL,OUTPUT_TREE_DESCRIPTOR,0,0,0,0,0,0}
+#define ARGS_DEFAULT__ {NULL,NULL,NULL,OUTPUT_TREE_DESCRIPTOR,0,0,0,0,0,0,0}
static const struct args ARGS_DEFAULT = ARGS_DEFAULT__;
struct cmd {
@@ -94,7 +95,7 @@ static void
usage(FILE* stream)
{
fprintf(stream,
-"usage: sln-get [-hlmnrv] [-L nlevels] [-R nlevels] [-w wavenumber] -i lines\n"
+"usage: sln-get [-hlmnrsv] [-L nlevels] [-R nlevels] [-w wavenumber] -i lines\n"
" -p molparams [tree]\n");
}
@@ -142,7 +143,7 @@ args_init(struct args* args, int argc, char** argv)
*args = ARGS_DEFAULT;
- while((opt = getopt(argc, argv, "hi:L:lmnp:R:rvw:")) != -1) {
+ while((opt = getopt(argc, argv, "hi:L:lmnp:R:rsvw:")) != -1) {
switch(opt) {
case 'h':
usage(stdout);
@@ -155,6 +156,8 @@ args_init(struct args* args, int argc, char** argv)
}
break;
case 'l': tree_descent(args, LEFT, 1); break;
+ case 'm': args->output_type = OUTPUT_NODE_MESH; break;
+ case 'n': args->output_type = OUTPUT_NODE_DESCRIPTOR; break;
case 'p': args->molparams = optarg; break;
case 'R':
if((res = cstr_to_uint(optarg, &nlvls)) == RES_OK) {
@@ -162,8 +165,7 @@ args_init(struct args* args, int argc, char** argv)
}
break;
case 'r': tree_descent(args, RIGHT, 1); break;
- case 'm': args->output_type = OUTPUT_NODE_MESH; break;
- case 'n': args->output_type = OUTPUT_NODE_DESCRIPTOR; break;
+ case 's': args->lines_in_shtr_format = 1; break;
case 'v': args->verbose += (args->verbose < 3); break;
case 'w':
args->output_type = OUTPUT_NODE_VALUE;
@@ -220,12 +222,41 @@ cmd_release(struct cmd* cmd)
}
static res_T
+load_lines(struct cmd* cmd, const struct args* args)
+{
+ res_T res = RES_OK;
+ ASSERT(cmd && args);
+
+ if(args->lines_in_shtr_format) {
+ struct shtr_line_list_read_args read_args = SHTR_LINE_LIST_READ_ARGS_NULL;
+
+ /* Loads lines from data serialized by the Star-HITRAN library */
+ read_args.filename = args->lines;
+ res = shtr_line_list_read(cmd->shtr, &read_args, &cmd->lines);
+ if(res != RES_OK) goto error;
+
+ } else {
+ struct shtr_line_list_load_args load_args = SHTR_LINE_LIST_LOAD_ARGS_NULL;
+
+ /* Loads lines from a file in HITRAN format */
+ load_args.filename = args->lines;
+ res = shtr_line_list_load(cmd->shtr, &load_args, &cmd->lines);
+ if(res != RES_OK) goto error;
+ }
+
+exit:
+ return res;
+error:
+ if(cmd->lines) { SHTR(line_list_ref_put(cmd->lines)); cmd->lines = NULL; }
+ goto exit;
+}
+
+static res_T
cmd_init(struct cmd* cmd, const struct args* args)
{
struct sln_device_create_args sln_args = SLN_DEVICE_CREATE_ARGS_DEFAULT;
struct sln_tree_read_args tree_args = SLN_TREE_READ_ARGS_NULL;
struct shtr_create_args shtr_args = SHTR_CREATE_ARGS_DEFAULT;
- struct shtr_line_list_load_args lines_args = SHTR_LINE_LIST_LOAD_ARGS_NULL;
res_T res = RES_OK;
ASSERT(cmd && args);
@@ -239,8 +270,7 @@ cmd_init(struct cmd* cmd, const struct args* args)
res = shtr_isotope_metadata_load(cmd->shtr, args->molparams, &cmd->molparams);
if(res != RES_OK) goto error;
- lines_args.filename = args->lines;
- res = shtr_line_list_load(cmd->shtr, &lines_args, &cmd->lines);
+ res = load_lines(cmd, args);
if(res != RES_OK) goto error;
sln_args.verbose = args->verbose;