commit 670c81368805ac57807b5305dae9c83111e6fc94
parent b39671999e2ef3262b36a2cd26d2a648ff8e9ef6
Author: Vincent Forest <vincent.forest@meso-star.com>
Date: Tue, 10 Mar 2026 13:46:47 +0100
No longer serialize spectroscopic data
This data is not specific to a tree, but shared among different trees.
Storing the list of lines from the HITRAN database and molecular
parameters therefore adds unnecessary fixed storage costs to all trees.
However, the tree is no longer autonomous, and the caller must provide
both the isotopological metadata and the list of lines used to construct
it in order to deserialize it, ensuring that these are indeed the ones
used during its construction.
In any case, an additional check is performed during the deserialization
of the tree in order to verify that the hash of the lines and their
metadata provided by the user correspond to the lines and metadata used
to build the tree, their hase now being serialized with the tree data.
If these signatures are different, an error is returned. This prevents
the use of spectroscopic data that is unsuitable for the tree.
The tests and the sln-get tool have been updated accordingly, as has the
manual page for the latter, which reflects the change in syntax for this
tool in light of this commit. This manual page has also been thoroughly
proofread.
Diffstat:
5 files changed, 170 insertions(+), 69 deletions(-)
diff --git a/doc/sln-get.1 b/doc/sln-get.1
@@ -15,13 +15,13 @@
.\"
.\" You should have received a copy of the GNU General Public License
.\" along with this program. If not, see <http://www.gnu.org/licenses/>.
-.Dd March 6, 2026
+.Dd March 10, 2026
.Dt SLN-GET 1
.Os
.\""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""
.Sh NAME
.Nm sln-get
-.Nd data accessor of a tree used to sample lines by importance
+.Nd data accesser for a tree built to speed up line sampling
.\""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""
.Sh SYNOPSIS
.Nm
@@ -29,18 +29,16 @@
.Op Fl L Ar nlevels
.Op Fl R Ar nlevels
.Op Fl w Ar wavenumber
+.Fl i Ar lines
+.Fl p Ar molparams
.Op Ar tree
.\""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""
.Sh DESCRIPTION
.Nm
-prints on the standard output the data of a
-.Ar tree ,
-built by
-.Xr sln-build 1 ,
-in order to speed up the importance sampling of spectral lines.
-If no structure is defined as an input argument, the
-.Ar tree
-is read on the standard input.
+queries a tree constructed by
+.Xr sln-build 1 .
+If no tree is specified as an input argument, its data is read from
+standard input.
.Pp
By default,
.Nm
@@ -49,13 +47,18 @@ structures, the total number of nodes, or the number of vertices used to
mesh the hierarchical representation of the high-resolution spectrum it
represents.
.Pp
-The data for a node can be printed
-.Pq options Fl m , Fl n , No and Fl w .
-By default, the node for which information is displayed is the root
-node.
+.Nm
+can also display data specific to a particular node, such as its general
+description
+.Pq option Fl n ,
+the mesh of its polyline
+.Pq option Fl m ,
+or the value of the spectrum it represents
+.Pq options Fl w .
+By default, the node queried is the root node.
The caller can select another node by visiting the tree using the
traversal options
-.Pq options Fl L , Fl l , Fl R No and Fl r .
+.Pq Fl L , Fl l , Fl R No and Fl r .
.Pp
The options are as follows:
.Bl -tag -width Ds
@@ -63,6 +66,10 @@ The options are as follows:
.It Fl h
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.
+.\""""""""""""""""""""""""""""""""""
.It Fl L Ar nlevels
Descend the tree by visiting the left child of the nodes traversed up to
.Ar nlevels
@@ -80,6 +87,10 @@ 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 p Ar molparams
+Isotopologue metadata from which the tree was built.
+The data is in HITRAN format.
+.\""""""""""""""""""""""""""""""""""
.It Fl R Ar nlevels
Descend the tree by visiting the right child of the nodes traversed up
to
@@ -98,11 +109,11 @@ option, but here for the right child of the node, instead of the left
child.
.\""""""""""""""""""""""""""""""""""
.It Fl m
-Prints the mesh of the current node, i.e. the mesh representing the high
+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 mesh vertices, in plain text, where each
-line represents the two values of a mesh vertex, separated by a space:
+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
@@ -124,10 +135,10 @@ Calculate the spectrum value of the current node at the given
.Ar wavenumber
in cm^-1.
Both the actual spectrum value, calculated from the lines that the node
-partitions, and the estimated value from its mesh are printed.
+partitions, and the estimated value from its polyline are printed.
The output format is as follows:
.Bd -literal -offset Ds
-"ka(%e) = %e ~ %e\en", wavenumber, ka_node, ka_mesh
+"ka(%e) = %e ~ %e\en", wavenumber, ka_node, ka_polyline
.Ed
.El
.\""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""
@@ -135,22 +146,27 @@ The output format is as follows:
.Ex -std
.\""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""
.Sh EXAMPLES
-Display the description of a tree:
+Display the description of a tree built from the list of lines and
+isotopologue metadata stored in
+.Pa lines.par
+and
+.Pa molparams.txt ,
+respectively:
.Bd -literal -offset Ds
-sln-get tree.sln
+sln-get -i lines.par -p molparams.txt tree.sln
.Ed
.\""""""""""""""""""""""""""""""""""
.Pp
Print the description of the child left of the tree root:
.Bd -literal -offset Ds
-sln-get -ln tree.sln
+sln-get -i lines.par -p molparams.txt -ln tree.sln
.Ed
.\""""""""""""""""""""""""""""""""""
.Pp
Print the description of one of the grandchildren of the root of the
tree, in this case the right child of its left child:
.Bd -literal -offset Ds
-sln-get -lrn tree.sln
+sln-get -i lines.par -p molparams.txt -lrn tree.sln
.Ed
.\""""""""""""""""""""""""""""""""""
.Pp
@@ -162,21 +178,45 @@ then the right child
and finally the left child
.Pq option Fl l
of the next two levels.
-Then output the mesh of the node reached
+Then output the polyline of the node reached
.Pq option Fl m
and save it to the
-.Pa mesh.txt
+.Pa polyline.txt
file:
.Bd -literal -offset Ds
-sln-get -L3 -rl -m tree.sln > mesh.txt
+sln-get -i lines.par -p molparams.txt -L3 -rl -m tree.sln \e
+ > polyline.txt
.Ed
.\""""""""""""""""""""""""""""""""""
.Pp
Print the spectrum value at 50 cm^-1 for one of the great-grandchildren
of the root:
.Bd -literal -offset Ds
-sln-get -L3 -w 50 tree.sln
+sln-get -i lines.par -p molparams.txt -L3 -w 50 tree.sln
.Ed
.\""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""
.Sh SEE ALSO
.Xr sln-build 1
+.Rs
+.%T The HITRAN Database
+.%U https://hitran.org/
+.Re
+.\""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""
+.Sh STANDARDS
+.Rs
+.%A L.S. Rothman et al.
+.%T The HITRAN2012 molecular spectroscopic database
+.%J Journal of Quantitative Spectroscopy & Radiative Transfer
+.%V 130
+.%P pp. 4\(en50
+.%D 2013
+.Re
+.Pp
+.Rs
+.%A L.S. Rothman et al.
+.%T HITEMP, the high-temperature molecular spectroscopic database
+.%J Journal of Quantitative Spectroscopu & Radiative Transfer
+.%V 111
+.%P pp. 2139\(en2150
+.%D 2010
+.Re
diff --git a/src/sln.h b/src/sln.h
@@ -86,7 +86,7 @@ struct sln_molecule {
static const struct sln_molecule SLN_MOLECULE_NULL = SLN_MOLECULE_NULL__;
struct sln_tree_create_args {
- /* Isotope metadata and lise of spectral lines */
+ /* Isotope metadata and list of spectral lines */
struct shtr_isotope_metadata* metadata;
struct shtr_line_list* lines;
@@ -121,9 +121,9 @@ static const struct sln_tree_create_args SLN_TREE_CREATE_ARGS_DEFAULT =
SLN_TREE_CREATE_ARGS_DEFAULT__;
struct sln_tree_read_args {
- /* Handle of the Star-HITRAN library to be used for loading isotopic metadata
- * and line lists */
- struct shtr* shtr;
+ /* Metadata and list of spectral lines from which the tree was constructed */
+ struct shtr_isotope_metadata* metadata;
+ struct shtr_line_list* lines;
/* Name of the file to read or of the provided stream.
* NULL <=> uses a default name for the stream to be read, which must
@@ -131,7 +131,7 @@ struct sln_tree_read_args {
const char* filename; /* Name of the file to read */
FILE* file; /* Stream from where data are read. NULL <=> read from file */
};
-#define SLN_TREE_READ_ARGS_NULL__ {NULL,NULL,NULL}
+#define SLN_TREE_READ_ARGS_NULL__ {NULL,NULL,NULL,NULL}
static const struct sln_tree_read_args SLN_TREE_READ_ARGS_NULL =
SLN_TREE_READ_ARGS_NULL__;
diff --git a/src/sln_get.c b/src/sln_get.c
@@ -48,6 +48,8 @@ enum output_type {
struct args {
const char* tree; /* NULL <=> read from standard input */
+ const char* molparams;
+ const char* lines;
enum output_type output_type;
@@ -69,7 +71,7 @@ struct args {
int quit;
int verbose;
};
-#define ARGS_DEFAULT__ {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}
static const struct args ARGS_DEFAULT = ARGS_DEFAULT__;
struct cmd {
@@ -79,6 +81,8 @@ struct cmd {
struct sln_tree* tree;
struct shtr* shtr;
+ struct shtr_isotope_metadata* molparams;
+ struct shtr_line_list* lines;
};
#define CMD_NULL__ {0}
static const struct cmd CMD_NULL = CMD_NULL__;
@@ -90,7 +94,8 @@ static void
usage(FILE* stream)
{
fprintf(stream,
-"usage: sln-get [-hlmnrv] [-L nlevels] [-R nlevels] [-w wavenumber] [tree]\n");
+"usage: sln-get [-hlmnrv] [-L nlevels] [-R nlevels] [-w wavenumber] -i lines\n"
+" -p molparams [tree]\n");
}
static void
@@ -137,18 +142,20 @@ args_init(struct args* args, int argc, char** argv)
*args = ARGS_DEFAULT;
- while((opt = getopt(argc, argv, "hL:lmnR:rvw:")) != -1) {
+ while((opt = getopt(argc, argv, "hi:L:lmnp:R:rvw:")) != -1) {
switch(opt) {
case 'h':
usage(stdout);
args->quit = 1;
goto exit;
+ case 'i': args->lines = optarg; break;
case 'L':
if((res = cstr_to_uint(optarg, &nlvls)) == RES_OK) {
tree_descent(args, LEFT, nlvls);
}
break;
case 'l': tree_descent(args, LEFT, 1); break;
+ case 'p': args->molparams = optarg; break;
case 'R':
if((res = cstr_to_uint(optarg, &nlvls)) == RES_OK) {
tree_descent(args, RIGHT, nlvls);
@@ -173,6 +180,17 @@ args_init(struct args* args, int argc, char** argv)
}
}
+ #define MANDATORY(Cond, Name, Opt) { \
+ if(!(Cond)) { \
+ fprintf(stderr, "%s: %s missing -- option '-%c'\n", argv[0], (Name), (Opt)); \
+ res = RES_BAD_ARG; \
+ goto error; \
+ } \
+ } (void)0
+ MANDATORY(args->molparams, "molparams", 'p');
+ MANDATORY(args->lines, "line list", 'i');
+ #undef MANDATORY
+
#ifndef NDEBUG
{
const uint64_t descent_mask =
@@ -197,6 +215,8 @@ cmd_release(struct cmd* cmd)
if(cmd->sln) SLN(device_ref_put(cmd->sln));
if(cmd->tree) SLN(tree_ref_put(cmd->tree));
if(cmd->shtr) SHTR(ref_put(cmd->shtr));
+ if(cmd->molparams) SHTR(isotope_metadata_ref_put(cmd->molparams));
+ if(cmd->lines) SHTR(line_list_ref_put(cmd->lines));
}
static res_T
@@ -205,6 +225,7 @@ 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);
@@ -215,11 +236,19 @@ cmd_init(struct cmd* cmd, const struct args* args)
res = shtr_create(&shtr_args, &cmd->shtr);
if(res != RES_OK) goto error;
+ 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);
+ if(res != RES_OK) goto error;
+
sln_args.verbose = args->verbose;
res = sln_device_create(&sln_args, &cmd->sln);
if(res != RES_OK) goto error;
- tree_args.shtr = cmd->shtr;
+ tree_args.metadata = cmd->molparams;
+ tree_args.lines = cmd->lines;
if(args->tree) {
tree_args.file = NULL;
tree_args.filename = args->tree;
diff --git a/src/sln_tree.c b/src/sln_tree.c
@@ -191,7 +191,7 @@ check_sln_tree_create_args
if(!args) return RES_BAD_ARG;
if(!args->metadata) {
- ERROR(sln, "%s: the metadata is missing.\n", caller);
+ ERROR(sln, "%s: the isotope metadata are missing.\n", caller);
return RES_BAD_ARG;
}
@@ -237,9 +237,13 @@ check_sln_tree_read_args
{
if(!args) return RES_BAD_ARG;
- if(!args->shtr) {
- ERROR(sln,
- "%s: the handle to the Star-HITRAN library is missing.\n", caller);
+ if(!args->metadata) {
+ ERROR(sln, "%s: the isotope metadata are missing.\n", caller);
+ return RES_BAD_ARG;
+ }
+
+ if(!args->lines) {
+ ERROR(sln, "%s: the list of lines is missing.\n", caller);
return RES_BAD_ARG;
}
@@ -412,7 +416,11 @@ sln_tree_read
const struct sln_tree_read_args* args,
struct sln_tree** out_tree)
{
- struct shtr_line_list_read_args rlines_args = SHTR_LINE_LIST_READ_ARGS_NULL;
+ hash256_T hash_mdata1;
+ hash256_T hash_mdata2;
+ hash256_T hash_lines1;
+ hash256_T hash_lines2;
+
struct stream stream = STREAM_NULL;
struct sln_tree* tree = NULL;
size_t n = 0;
@@ -438,6 +446,8 @@ sln_tree_read
} else { \
res = RES_UNKNOWN_ERR; \
} \
+ ERROR(sln, "%s: error loading the tree structure -- %s.\n", \
+ stream.name, res_to_cstr(res)); \
goto error; \
} \
} (void)0
@@ -445,11 +455,41 @@ sln_tree_read
if(version != SLN_TREE_VERSION) {
ERROR(sln,
"%s: unexpected tree version %d. Expecting a tree in version %d.\n",
- FUNC_NAME, version, SLN_TREE_VERSION);
+ stream.name, version, SLN_TREE_VERSION);
+ res = RES_BAD_ARG;
+ goto error;
+ }
+
+ res = shtr_isotope_metadata_hash(args->metadata, hash_mdata1);
+ if(res != RES_OK) goto error;
+
+ READ(hash_mdata2, sizeof(hash256_T));
+ if(!hash256_eq(hash_mdata1, hash_mdata2)) {
+ ERROR(sln,
+ "%s: the input isotopic metadata are not those used "
+ "during tree construction.\n", stream.name);
+ res = RES_BAD_ARG;
+ goto error;
+ }
+
+ SHTR(isotope_metadata_ref_get(args->metadata));
+ tree->args.metadata = args->metadata;
+
+ res = shtr_line_list_hash(args->lines, hash_lines1);
+ if(res != RES_OK) goto error;
+
+ READ(hash_lines2, sizeof(hash256_T));
+ if(!hash256_eq(hash_lines1, hash_lines2)) {
+ ERROR(sln,
+ "%s: the input list of lines is not the one used to build the tree.\n",
+ stream.name);
res = RES_BAD_ARG;
goto error;
}
+ SHTR(line_list_ref_get(args->lines));
+ tree->args.lines = args->lines;
+
READ(&n, 1);
if((res = darray_node_resize(&tree->nodes, n)) != RES_OK) goto error;
READ(darray_node_data_get(&tree->nodes), n);
@@ -467,24 +507,12 @@ sln_tree_read
READ(&tree->args.mesh_type, 1);
#undef READ
- res = shtr_isotope_metadata_create_from_stream
- (args->shtr, stream.fp, &tree->args.metadata);
- if(res != RES_OK) goto error;
-
- rlines_args.file = stream.fp;
- res = shtr_line_list_read(args->shtr, &rlines_args, &tree->args.lines);
- if(res != RES_OK) goto error;
-
exit:
stream_release(&stream);
if(out_tree) *out_tree = tree;
return res;
error:
if(tree) { SLN(tree_ref_put(tree)); tree = NULL; }
- if(sln) {
- ERROR(sln, "%s: error loading the tree structure -- %s.\n",
- FUNC_NAME, res_to_cstr(res));
- }
goto exit;
}
@@ -673,15 +701,21 @@ sln_tree_write
(const struct sln_tree* tree,
const struct sln_tree_write_args* args)
{
- struct shtr_line_list_write_args wlines_args = SHTR_LINE_LIST_WRITE_ARGS_NULL;
struct stream stream = STREAM_NULL;
size_t nnodes, nverts;
+ hash256_T hash_mdata;
+ hash256_T hash_lines;
res_T res = RES_OK;
if(!tree) { res = RES_BAD_ARG; goto error; }
res = check_sln_tree_write_args(tree->sln, FUNC_NAME, args);
if(res != RES_OK) goto error;
+ res = shtr_isotope_metadata_hash(tree->args.metadata, hash_mdata);
+ if(res != RES_OK) goto error;
+ res = shtr_line_list_hash(tree->args.lines, hash_lines);
+ if(res != RES_OK) goto error;
+
res = stream_init
(tree->sln, FUNC_NAME, args->filename, args->file, "w", &stream);
if(res != RES_OK) goto error;
@@ -695,6 +729,8 @@ sln_tree_write
} \
} (void)0
WRITE(&SLN_TREE_VERSION, 1);
+ WRITE(hash_mdata, sizeof(hash256_T));
+ WRITE(hash_lines, sizeof(hash256_T));
nnodes = darray_node_size_get(&tree->nodes);
WRITE(&nnodes, 1);
@@ -713,13 +749,6 @@ sln_tree_write
WRITE(&tree->args.mesh_type, 1);
#undef WRITE
- res = shtr_isotope_metadata_write(tree->args.metadata, stream.fp);
- if(res != RES_OK) goto error;
-
- wlines_args.file = stream.fp;
- res = shtr_line_list_write(tree->args.lines, &wlines_args);
- if(res != RES_OK) goto error;
-
exit:
stream_release(&stream);
return res;
diff --git a/src/test_sln_tree.c b/src/test_sln_tree.c
@@ -332,8 +332,7 @@ check_tree_equality
static void
test_tree_serialization
(struct sln_device* sln,
- const struct sln_tree_create_args* tree_args,
- struct shtr* shtr)
+ const struct sln_tree_create_args* tree_args)
{
struct sln_tree_write_args wargs = SLN_TREE_WRITE_ARGS_NULL;
struct sln_tree_read_args rargs = SLN_TREE_READ_ARGS_NULL;
@@ -356,13 +355,17 @@ test_tree_serialization
CHK(sln_tree_write(tree1, &wargs) == RES_OK);
rewind(fp);
- rargs.shtr = shtr;
+ rargs.metadata = tree_args->metadata;
+ rargs.lines = tree_args->lines;
rargs.file = fp;
CHK(sln_tree_read(NULL, &rargs, &tree2) == RES_BAD_ARG);
CHK(sln_tree_read(sln, NULL, &tree2) == RES_BAD_ARG);
- rargs.shtr = NULL;
+ rargs.metadata = NULL;
CHK(sln_tree_read(sln, &rargs, &tree2) == RES_BAD_ARG);
- rargs.shtr = shtr;
+ rargs.metadata = tree_args->metadata;
+ rargs.lines = NULL;
+ CHK(sln_tree_read(sln, &rargs, &tree2) == RES_BAD_ARG);
+ rargs.lines = tree_args->lines;
rargs.file = NULL;
CHK(sln_tree_read(sln, &rargs, &tree2) == RES_BAD_ARG);
rargs.file = fp;
@@ -450,7 +453,7 @@ main(int argc, char** argv)
tree_args.temperature = 296;
test_tree(sln, &tree_args, line_list);
- test_tree_serialization(sln, &tree_args, shtr);
+ test_tree_serialization(sln, &tree_args);
CHK(sln_device_ref_put(sln) == RES_OK);
CHK(shtr_ref_put(shtr) == RES_OK);