A lightweight, dependency-free HSD (Human-friendly Structured Data) parser for Fortran.
HSD is a data format similar to JSON and YAML, but designed to minimize the effort for humans to read and write it. It was originally developed as the input format for DFTB+.
This library provides:
- Parsing HSD files into a tree structure with full error reporting
- Serialization of data structures back to HSD format
- Include support with cycle detection (
<<+for HSD,<<<for text) - Path-based accessors for convenient data retrieval (
"section/subsection/value") - Schema validation for declarative input validation
- Type introspection to query node types before access
- Tree operations including merge, clone, and visitor pattern traversal
program example
use hsd
implicit none
type(hsd_table) :: root
type(hsd_error_t), allocatable :: error
integer :: max_steps, stat
real(dp) :: temperature
! Load HSD file
call hsd_load("input.hsd", root, error)
if (allocated(error)) then
call error%print()
stop 1
end if
! Access values with path navigation
call hsd_get(root, "Driver/MaxSteps", max_steps, stat)
call hsd_get(root, "Hamiltonian/DFTB/Filling/Fermi/Temperature", temperature, stat)
! Use defaults for optional values
call hsd_get_or(root, "Driver/Timeout", max_steps, default=3600, stat=stat)
! Type introspection
if (hsd_is_table(root, "Hamiltonian/DFTB")) then
print *, "Child count:", hsd_child_count(root, "Hamiltonian/DFTB")
end if
! Modify and save
call hsd_set(root, "Driver/MaxSteps", 200)
call hsd_dump(root, "output.hsd")
end program exampleExample input file:
Driver = ConjugateGradient {
MaxSteps = 100
}
Hamiltonian = DFTB {
SCC = Yes
Filling = Fermi {
Temperature [Kelvin] = 300.0
}
}
- Fortran 2008 compatible compiler (gfortran ≥ 7, ifort ≥ 18, ifx)
- CMake ≥ 3.14
cmake -B build
cmake --build build
ctest --test-dir build # Run tests
cmake --install build --prefix /path/to/install| Option | Default | Description |
|---|---|---|
HSD_ACCEPT_TRUE_FALSE |
ON |
Accept True/False as boolean values |
HSD_BUILD_TESTS |
ON |
Build test suite |
HSD_BUILD_EXAMPLES |
ON |
Build examples |
HSD_BUILD_DOCS |
OFF |
Build API documentation with FORD |
HSD_COVERAGE |
OFF |
Enable code coverage (GCC only) |
fpm build
fpm test
fpm run --example simple_read| Procedure | Description |
|---|---|
hsd_load(file, root, error) |
Parse HSD file into tree |
hsd_load_string(str, root, error) |
Parse HSD string |
hsd_dump(root, file) |
Write tree to file |
hsd_dump_to_string(root, str) |
Serialize tree to string |
| Procedure | Description |
|---|---|
hsd_get(root, path, value, stat) |
Get value at path |
hsd_get_or(root, path, value, default, stat) |
Get with fallback default |
hsd_get_matrix(root, path, matrix, stat) |
Get 2D array |
| Procedure | Description |
|---|---|
hsd_has_child(root, path) |
Check if path exists |
hsd_is_table(root, path) |
Check if node is a table |
hsd_is_value(root, path) |
Check if node is a leaf value |
hsd_child_count(root, path) |
Count children in table |
hsd_get_keys(root, path, keys) |
Get child key names |
hsd_get_attrib(root, path, attrib) |
Get attribute (e.g., unit) |
| Procedure | Description |
|---|---|
hsd_set(root, path, value) |
Set value at path |
hsd_remove_child(root, path) |
Remove child node |
hsd_merge(target, source) |
Merge two trees |
hsd_clone(source, dest) |
Deep copy a tree |
| Procedure | Description |
|---|---|
hsd_require(root, path, value, error) |
Get required value or error |
hsd_validate_range(value, min, max, error) |
Validate numeric range |
hsd_validate_one_of(value, options, error) |
Validate against allowed values |
schema_validate(schema, root, error) |
Validate against schema |
# Comments start with #
TagName {
ChildTag = value
WithUnit [eV/Angstrom] = 0.001
}
# Short form (assigns type to parent)
TagName = TypeName {
...
}
| Type | Examples |
|---|---|
| Strings | name = hello, name = "hello world" |
| Integers | count = 42 |
| Reals | value = 3.14, value = 1.0e-5 |
| Booleans | enabled = Yes, enabled = No, enabled = On/Off |
| Arrays | values = 1 2 3 4 5, values = 1, 2, 3 |
KPoints {
4 0 0
0 4 0
0 0 4
}
# Include and parse another HSD file
<<+ "other.hsd"
# Include raw text content
Data {
<<< "datafile.dat"
}
Circular includes are automatically detected and reported as errors.
All parsing and access operations can return detailed errors:
type(hsd_error_t), allocatable :: error
call hsd_load("config.hsd", root, error)
if (allocated(error)) then
call error%print() ! Prints formatted error with location
! Access details: error%code, error%message, error%line_start, error%hint
end ifError codes include: HSD_STAT_OK, HSD_STAT_SYNTAX_ERROR, HSD_STAT_TYPE_ERROR, HSD_STAT_NOT_FOUND, HSD_STAT_FILE_NOT_FOUND, and more.
Full documentation is available at https://elv3rs.github.io/hsd-fortran:
- API Reference — Auto-generated API docs (FORD)
- User Guide — Tutorials and examples
- HSD Format — Format specification
- Error Handling — Error types and handling
- Thread Safety — Concurrency guidelines
- Fuzz Testing — AFL++ fuzzing guide
Install FORD and build the API documentation:
pip install ford
cmake -B build -DHSD_BUILD_DOCS=ON
cmake --build build --target docs
# Documentation will be in ford_docs/BSD 2-Clause Plus Patent License. See LICENSE for details.
- hsd-python — Python implementation
- DFTB+ — Original user of HSD format