This repository is based upon source code for Herb M. Pickett's SPFIT and SPCAT programs for fitting and simulating rotational spectra. Originally called the CalPGM program suite, it also includes DPFIT and DPCAT, as well as auxiliary programs CALMRG, MOIAM, STARK, TERMVAL, SORTN, REASSIGN, SORTEGY, IAMBAK, and IAMCALC. The original code dates to 1989.
The base version of this code was obtained from the Cologne Database for Molecular Spectroscopy, with the latest modifications by Dr. Holger S. P. Müller (the Pickett_neu.zip file package).
This repository also includes modifications by Kelvin Lee, to allow the partition function routine to calculate the rotational partition function up to 1000 K programmatically without needing to repetitively call SPCAT.
The codebase has been refactored from monolithic C into a modular C++ architecture while preserving exact numerical compatibility with the 2008 reference binaries across a 55-file test suite.
Key changes:
CalculationEngineinterface withSpinvEngineandDpiEngineimplementations (strategy pattern), replacing global state with context structs.CalFitclass encapsulates spfit's iterative Marquardt-Levenberg fitting. I/O separated intoCalFitIO. Entry pointfit_main.cppis a thin wrapper.CalCatclass encapsulates spcat's catalog generation. I/O separated intoCalCatIO. Entry pointcat_main.cppis a thin wrapper.- Both classes take C++ input structs and produce output structs, enabling programmatic use without files.
Command-line usage of all executables (spfit, spcat, etc.) is unchanged from the original C code. Existing input files and workflows work without modification.
See TASKS.md for the full modernization roadmap and status.
make # builds spfit, spcat, calmrg (default targets)
make all # builds all executables
make install # installs to /usr/local/bin
make clean # removes build productsCMake provides IDE integration, dependency tracking, and configurable BLAS selection. Always run from a build subdirectory:
mkdir -p build && cd build
cmake .. -DUSE_SYSTEM_BLAS=OFF
makeDo not run cmake from the project root — it will overwrite the hand-maintained Makefile.
To match the v2008 reference outputs exactly, the build must use the bundled dblas.c fallback BLAS (not OpenBLAS or other optimized BLAS) and disable FMA contraction:
- CMake:
cmake .. -DUSE_SYSTEM_BLAS=OFF(default flags include-ffp-contract=off) - Makefile: leave
BLASLIBundefined (usesdblas.oautomatically);CFLAGSincludes-ffp-contract=off
This is necessary because FMA instructions and optimized BLAS libraries change floating-point accumulation order at the ULP level, which can flip the sign of near-zero Householder reflectors in the least-squares solver. The resulting .var output is physically equivalent but not bit-identical to the reference. See TASKS.md Task 6 for full details.
There is a very useful set of notes and usage information available here.
You will also need to read the included original documentation for SPCAT and SPFIT.
Original documentation for DPFIT and DPCAT are also included.
fit_main.cpp— spfit entry point (thin wrapper aroundCalFit)cat_main.cpp— spcat entry point (thin wrapper aroundCalCat)calmrg.c— catalog merge utility
CalFit.cpp,CalFit_helpers.cpp,CalFit.hpp— fitting logic (Marquardt-Levenberg iteration)CalFitIO.cpp,CalFitIO.hpp— spfit file I/O (.par,.lin,.fit,.var)CalCat.cpp,CalCat_helpers.cpp,CalCat.hpp— catalog generation logicCalCatIO.cpp,CalCatIO.hpp— spcat file I/O (.int,.var)CalculationEngine.hpp— abstract interface for Hamiltonian computation
SpinvEngine.cpp,SpinvEngine.hpp,SpinvContext.hpp— engine for asymmetric/symmetric tops and general moleculesspinv_setup.c,spinv_spin_symmetry.c,spinv_linalg_sort.c,spinv_hamiltonian.c,spinv_utils.c,spinv_internal.h— underlying C implementation (parameterized viaSpinvContext)DpiEngine.cpp,DpiEngine.hpp,DpiContext.hpp— engine for diatomic/linear molecules with nuclear spindpi.c,dpi.h— underlying C implementation (parameterized viaDpiContext)
ulib.c— parameter I/O (getpar,getvar,putvar), error analysis (calerr,prcorr)lsqfit.c,lsqfit.h— least-squares fitting (QR factorization, Marquardt-Levenberg solver)dblas.c— fallback BLAS routines (required for exact numerical reproduction)cnjj.c,cnjj.h— Clebsch-Gordan coefficientsslibgcc.c— system-dependent functionscatutil.c,catutil.h— catalog utility functionssubfit.c— supplementary fitting routinesspinit.c,spinit.h— spin initialization
spinv.md— SPFIT/SPCAT documentationdpi.md— DPFIT/DPCAT documentationTASKS.md— modernization task list and status
These are original Pickett utility programs. They compile and link (make all) but have not been modernized, tested, or verified against the current codebase. They link against the old splib.a static library rather than the new C++ classes.
Data preparation / workflow:
sortn.c— standalone CLI tool for sorting catalog files by frequency (wraps thesortn()function fromsortsub.c, which is also linked directly into spcat)calmrg.c— merges experimental lines into a catalog structurecalbak.c— converts.catformat back to.linfor re-fittingsortegy.c— sorts energy levels by quantum number rulestermval.c— associates spectral lines with upper/lower energy levels from a term listreassign.c— reassigns quantum numbers to lines via a rule file
Internal rotation (IAM) sub-suite:
moiam.c— computes internal rotation structural parameters from molecular coordinatesiamcalc.c— Hamiltonian diagonalization for internal rotation analysisiambak.c— back-transforms fitted IAM parameters between models
Specialized:
stark.c— computes Stark effect coefficients for molecular energy levels
Other:
*.nam— parameter name files for labeling output. Searched in the current directory, then in the directory specified by theSPECNAMEenvironment variable.
The trust region approach to Marquardt-Levenberg parameter adjustment is described in John. E. Dennis
and Robert B. Schnabel, Numerical Methods for Unconstrained Optimization and Non-linear Equations,
Prentice-Hall, 1983. The basic idea is that there is a region over which a linear least squares fit can be
trusted. First the value of each parameter is scaled so that the squares of the derivatives, summed over
all the lines, is unity. Then a simple least squares fit is attempted with a value of zero for the
Marquardt-Levenberg parameter,
This modernization project began as an experiment to test coding assistants handling legacy code and porting to different languages. It has served as a good test of the limits of current AI, particularly given the large code-base size and complexity of monolithic code. Initially, leading LLMs and agents could often give good feedback, but struggled to complete tasks correctly, and frequently got into coding/debugging loops from which they couldn't escape. I often needed to step in and do things correctly or interrupt and redirect. Managing context has been and remains an issue (though it is much easier with now-smaller file sizes). But LLMs and agents have developed impressively over the last year (Spring 2025 to Spring 2026), to the point where they are much more able to complete large and complex tasks, to include debugging well-obscured bugs, nearly without correction. It has been a fruitful experiment. Hopefully the result will also prove useful.