star-line

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

commit 9468a804126acd6ba536ba8894c5532ef957bace
parent e567d47562d1327e2a040e7b0d9f4fb232438fbd
Author: Vincent Forest <vincent.forest@meso-star.com>
Date:   Thu, 12 Mar 2026 11:51:42 +0100

Start the implementation of the sln-slab tool

Its purpose is to perform a spectral computation using Monte Carlo in a
1D slab of a homogeneous semi-transparent medium, primarily for testing
purposes. This will allow the calculated values to be compared with
results from other simulation codes, or verified against a past
reference to ensure there is no regression, particularly in the
construction of acceleration structures.

Currently, only the syntactic analysis of the command arguments is
performed, as well as the loading of the acceleration structure.

Diffstat:
MMakefile | 7+++++--
Asrc/sln_slab.c | 248+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 files changed, 253 insertions(+), 2 deletions(-)

diff --git a/Makefile b/Makefile @@ -82,7 +82,7 @@ libsln.o: $(OBJ) ################################################################################ # Utils ################################################################################ -UTIL_SRC = src/sln_build.c src/sln_get.c +UTIL_SRC = src/sln_build.c src/sln_get.c src/sln_slab.c UTIL_OBJ = $(UTIL_SRC:.c=.o) UTIL_DEP = $(UTIL_SRC:.c=.d) @@ -97,7 +97,7 @@ LDFLAGS_UTIL = $(LDFLAGS_EXE) $(LIBS_UTIL) utils: library $(UTIL_DEP) @$(MAKE) -fMakefile \ $$(for i in $(UTIL_DEP); do printf -- '-f%s\n' "$${i}"; done) \ - sln-build sln-get + sln-build sln-get sln-slab sln-build: config.mk sln-local.pc src/sln_build.o $(LIBNAME) $(CC) $(CFLAGS_UTIL) -o $@ src/sln_build.o $(LDFLAGS_UTIL) @@ -105,6 +105,9 @@ sln-build: config.mk sln-local.pc src/sln_build.o $(LIBNAME) sln-get: config.mk sln-local.pc src/sln_get.o $(LIBNAME) $(CC) $(CFLAGS_UTIL) -o $@ src/sln_get.o $(LDFLAGS_UTIL) +sln-slab: config.mk sln-local.pc src/sln_slab.o $(LIBNAME) + $(CC) $(CFLAGS_UTIL) -o $@ src/sln_slab.o $(LDFLAGS_UTIL) + $(UTIL_DEP): config.mk sln-local.pc @$(CC) $(CFLAGS_UTIL) -MM -MT "$(@:.d=.o) $@" $(@:.d=.c) -MF $@ diff --git a/src/sln_slab.c b/src/sln_slab.c @@ -0,0 +1,248 @@ +/* Copyright (C) 2022, 2026 |Méso|Star> (contact@meso-star.com) + * Copyright (C) 2026 Université de Lorraine + * Copyright (C) 2022 Centre National de la Recherche Scientifique + * Copyright (C) 2022 Université Paul Sabatier + * + * This program is free software: you can redistribute it and/or modify + * it under the terms of the GNU General Public License as published by + * the Free Software Foundation, either version 3 of the License, or + * (at your option) any later version. + * + * This program is distributed in the hope that it will be useful, + * but WITHOUT ANY WARRANTY; without even the implied warranty of + * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the + * GNU General Public License for more details. + * + * You should have received a copy of the GNU General Public License + * along with this program. If not, see <http://www.gnu.org/licenses/>. */ + +#define _POSIX_C_SOURCE 200112L /* getopt */ + +#include "sln.h" + +#include <star/shtr.h> + +#include <rsys/cstr.h> +#include <rsys/mem_allocator.h> + +#include <unistd.h> /* getopt */ + +struct args { + const char* tree; /* Acceleration structure */ + const char* molparams; + const char* lines; + + double thickness; /* Thickness of the slab [m] */ + + unsigned long nrealisations; /* Number of Monte Carlo realisations */ + + /* Miscellaneous */ + unsigned nthreads_hint; /* Hint on the number of threads to use */ + int lines_in_shtr_format; + int verbose; + int quit; +}; +#define ARGS_DEFAULT__ {NULL,NULL,NULL,1,10000,UINT_MAX,0,0,0} +static const struct args ARGS_DEFAULT = ARGS_DEFAULT__; + +struct cmd { + struct args args; + + struct sln_tree* tree; +}; +#define CMD_NULL__ {0} +static const struct cmd CMD_NULL = CMD_NULL__; + +/******************************************************************************* + * Helper functions + ******************************************************************************/ +static void +usage(FILE* stream) +{ + fprintf(stream, +"usage: sln-slab [-hsv] [-n nrealisations] [-T thickness] [-t threads]\n" +" -a accel_struct -m molparams -l lines\n"); +} + +static res_T +args_init(struct args* args, int argc, char** argv) +{ + int opt = 0; + res_T res = RES_OK; + + ASSERT(args); + + *args = ARGS_DEFAULT; + + while((opt = getopt(argc, argv, "a:hl:m:n:sT:t:v")) != -1) { + switch(opt) { + case 'a': args->tree = optarg; break; + case 'h': + usage(stdout); + args->quit = 1; + goto exit; + case 'l': args->lines = optarg; break; + case 'm': args->molparams = optarg; break; + case 'n': res = cstr_to_ulong(optarg, &args->nrealisations); break; + case 's': args->lines_in_shtr_format = 1; break; + case 'T': + res = cstr_to_double(optarg, &args->thickness); + if(res == RES_OK && args->thickness <= 0) res = RES_BAD_ARG; + break; + case 't': + res = cstr_to_uint(optarg, &args->nthreads_hint); + if(res == RES_OK && args->nthreads_hint == 0) res = RES_BAD_ARG; + break; + case 'v': args->verbose += (args->verbose < 3); break; + default: res = RES_BAD_ARG; break; + } + if(res != RES_OK) { + if(optarg) { + fprintf(stderr, "%s: invalid option argument '%s' -- '%c'\n", + argv[0], optarg, opt); + } + goto error; + } + } + + #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", 'm'); + MANDATORY(args->lines, "line list", 'l'); + MANDATORY(args->tree, "acceleration structure", 's'); + #undef MANDATORY + +exit: + return res; +error: + usage(stderr); + goto exit; +} + +static res_T +load_lines + (struct shtr* shtr, + const struct args* args, + struct shtr_line_list** out_lines) +{ + struct shtr_line_list* lines = NULL; + res_T res = RES_OK; + ASSERT(shtr && args && out_lines); + + 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(shtr, &read_args, &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(shtr, &load_args, &lines); + if(res != RES_OK) goto error; + } + +exit: + *out_lines = lines; + return res; +error: + if(lines) { SHTR(line_list_ref_put(lines)); lines = NULL; } + goto exit; +} + + +static void +cmd_release(struct cmd* cmd) +{ + ASSERT(cmd); + if(cmd->tree) SLN(tree_ref_put(cmd->tree)); +} + +static res_T +cmd_init(struct cmd* cmd, const struct args* args) +{ + /* Star Line */ + 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 sln_device* sln = NULL; + + /* Star HITRAN */ + struct shtr_create_args shtr_args = SHTR_CREATE_ARGS_DEFAULT; + struct shtr* shtr = NULL; + struct shtr_isotope_metadata* molparams = NULL; + struct shtr_line_list* lines = NULL; + + /* Miscellaneous */ + res_T res = RES_OK; + + ASSERT(cmd && args); + + *cmd = CMD_NULL; + + shtr_args.verbose = args->verbose; + res = shtr_create(&shtr_args, &shtr); + if(res != RES_OK) goto error; + + res = shtr_isotope_metadata_load(shtr, args->molparams, &molparams); + if(res != RES_OK) goto error; + + res = load_lines(shtr, args, &lines); + if(res != RES_OK) goto error; + + sln_args.verbose = args->verbose; + res = sln_device_create(&sln_args, &sln); + if(res != RES_OK) goto error; + + tree_args.metadata = molparams; + tree_args.lines = lines; + tree_args.filename = args->tree; + res = sln_tree_read(sln, &tree_args, &cmd->tree); + if(res != RES_OK) goto error; + + cmd->args = *args; + +exit: + if(sln) SLN(device_ref_put(sln)); + if(shtr) SHTR(ref_put(shtr)); + if(molparams) SHTR(isotope_metadata_ref_put(molparams)); + if(lines) SHTR(line_list_ref_put(lines)); + return res; +error: + cmd_release(cmd); + *cmd = CMD_NULL; + goto exit; +} + +/******************************************************************************* + * Main function + ******************************************************************************/ +int +main(int argc, char** argv) +{ + struct args args = ARGS_DEFAULT; + struct cmd cmd = CMD_NULL; + int err = 0; + res_T res = RES_OK; + + if((res = args_init(&args, argc, argv)) != RES_OK) goto error; + if(args.quit) goto exit; + + if((res = cmd_init(&cmd, &args)) != RES_OK) goto error; + +exit: + cmd_release(&cmd); + CHK(mem_allocated_size() == 0); + return err; +error: + err = 1; + goto exit; +}