API Documentation

Complete C++ API reference for the Tomocam library.

Core Components

The core components provide fundamental data structures and operations:

  • Array: Multi-dimensional array container

  • PolarGrid: Polar coordinate grid geometry

  • Projection Operators: Forward and backward projection

Reconstruction Algorithms

Iterative reconstruction methods:

  • MBIR: Model-Based Iterative Reconstruction

  • Conjugate Gradient: CG solver

  • Split Bregman: TV regularization

  • qGGMRF: Generalized Gaussian MRF regularization

I/O Operations

File input/output:

  • TIFF I/O: Read and write TIFF files

  • VTI Writer: Export VTK ImageData format

  • Configuration: TOML configuration parsing

Utilities

Helper functions and utilities:

  • FFT: Fast Fourier Transform wrappers

  • NUFFT: Non-uniform FFT operations

  • Timer: Performance profiling

  • Logger: Logging utilities

Class Reference

template<typename T>
class Array
#include <array.h>

Public Functions

inline Array()
inline Array(size_t x, size_t y, size_t z)
inline Array(dims_t d)
Array(const Array<T>&) = delete
Array<T> &operator=(const Array<T>&) = delete
inline Array(Array<T> &&rhs) noexcept
inline Array<T> &operator=(Array<T> &&rhs) noexcept
inline Array<T> clone() const
inline T *begin()
inline const T *begin() const
inline T *end()
inline const T *end() const
inline dims_t dims() const
inline size_t size() const
inline size_t nslices() const
inline size_t nrows() const
inline size_t ncols() const
inline T &operator[](size_t i)
inline const T &operator[](size_t i) const
inline T at(size_t i, size_t j, size_t k) const
inline T &operator[](dims_t i)
inline const T &operator[](dims_t i) const
inline Slice<T> slice(size_t begin, size_t end)
inline const Slice<T> slice(size_t begin, size_t end) const
inline std::span<T> row(size_t i, size_t j)
inline const std::span<T> row(size_t i, size_t j) const
inline void setPads(dims_t p)
inline dims_t pads() const
inline Array<T> operator*(const Array<T> &v)
inline Array<T> &operator*=(T v)
inline Array<T> operator*(T v)
inline Array<T> &operator/=(T v)
inline Array<T> &operator/=(const Array<T> &v)
inline Array<T> operator/(const Array<T> &v) const
inline Array<T> operator/(T scalar) const
inline Array<T> &operator+=(T v)
inline Array<T> &operator+=(const Array<T> &v)
inline Array<T> operator+(const Array<T> &rhs) const
inline Array<T> &operator+=(const Array<T> &rhs) const
inline Array<T> &operator-=(T v)
inline Array<T> &operator-=(const Array<T> &v)
inline Array<T> operator-(const Array<T> &rhs) const

Public Static Functions

static inline Array<T> zeros(const dims_t &d)
static inline Array<T> ones(const dims_t &d)
static inline Array<T> random(const dims_t &d)

Private Functions

inline size_t flatIdx(size_t i, size_t j, size_t k) const

Private Members

dims_t dims_
dims_t pads_
size_t size_
std::unique_ptr<T[]> ptr_
struct dims_t
#include <dtypes.h>

Public Functions

inline dims_t()
inline dims_t(size_t a, size_t b, size_t c)
inline dims_t(const std::array<size_t, 3> &arr)
inline std::tuple<size_t, size_t, size_t> unravel_idx(size_t idx) const
inline void set_x(size_t v)
inline void set_y(size_t v)
inline void set_z(size_t v)
inline size_t size() const
inline size_t flat_idx(size_t i, size_t j, size_t k) const
inline size_t x() const
inline size_t y() const
inline size_t z() const
inline dims_t operator+(const dims_t &v) const
inline dims_t operator-(const dims_t &v) const
inline dims_t operator*(int v) const
inline dims_t operator/(int v) const
inline bool operator==(const dims_t &v) const
inline bool operator!=(const dims_t &v) const

Public Members

size_t n1
size_t n2
size_t n3
template<typename T>
class FinufftPlanCache

Public Functions

FinufftPlanCache() = default
~FinufftPlanCache() = default
FinufftPlanCache(const FinufftPlanCache&) = delete
FinufftPlanCache &operator=(const FinufftPlanCache&) = delete
FinufftPlanCache(FinufftPlanCache&&) = delete
FinufftPlanCache &operator=(FinufftPlanCache&&) = delete
inline FinufftPlanWrapper<T> &get_plan(int type, int dim, std::array<int64_t, 3> n_modes, int iflag)

Private Members

FinufftPlanWrapper<T> type1_plan_
FinufftPlanWrapper<T> type2_plan_
std::once_flag type1_init_flag_
std::once_flag type2_init_flag_
template<typename T>
class FinufftPlanWrapper
#include <finufft_plan.h>

Public Functions

FinufftPlanWrapper() = default
inline void make_plan(int type, int dim, std::array<int64_t, 3> n_modes, int iflag)
inline void set_points(const PolarGrid<T> &pg)
inline int execute(std::complex<T> *cz, std::complex<T> *fz)
inline ~FinufftPlanWrapper()
FinufftPlanWrapper(const FinufftPlanWrapper&) = delete
FinufftPlanWrapper &operator=(const FinufftPlanWrapper&) = delete
inline bool valid() const
inline FinufftPlanWrapper(FinufftPlanWrapper &&other) noexcept
inline FinufftPlanWrapper &operator=(FinufftPlanWrapper &&other) noexcept

Private Types

using Traits = FinufftTraits<T>

Private Members

Traits::plan_type plan
bool initialized = false
template<typename T>
struct FinufftTraits
template<>
struct FinufftTraits<double>
#include <finufft_plan.h>

Public Types

using plan_type = finufft_plan
using complex_type = std::complex<double>

Public Static Functions

static inline int makeplan(int type, int dim, int64_t *n_modes, int iflag, int ntrans, plan_type *plan, finufft_opts *opts)
static inline int setpts(plan_type plan, int64_t npts, double *x, double *y, double *z, int nk, double *s, double *t, double *u)
static inline int execute(plan_type plan, complex_type *cz, complex_type *fz)
static inline void destroy(plan_type plan)
template<>
struct FinufftTraits<float>
#include <finufft_plan.h>

Public Types

using plan_type = finufftf_plan
using complex_type = std::complex<float>

Public Static Functions

static inline int makeplan(int type, int dim, int64_t *n_modes, int iflag, int ntrans, plan_type *plan, finufft_opts *opts)
static inline int setpts(plan_type plan, int64_t npts, float *x, float *y, float *z, int nk, float *s, float *t, float *u)
static inline int execute(plan_type plan, complex_type *cz, complex_type *fz)
static inline void destroy(plan_type plan)
class Logger
#include <logger.h>

Public Functions

inline Logger(LogMode mode = LogMode::STDOUT, const std::string &filename = "")
template<typename ...Args>
inline void log(const std::string &message)
inline LogMode get_mode() const
inline ~Logger()

Private Members

LogMode mode_
std::unique_ptr<std::ofstream> file_stream_
class NumPyRandom
#include <utils.h>

Public Functions

inline NumPyRandom()
template<typename T>
inline T rand()

Private Members

std::mt19937 gen_
struct OutputParams
#include <config.h>

Public Functions

OutputParams() = default
inline OutputParams(const toml::table &config)
inline bool has_format(const std::string &fmt) const

Public Members

std::string filepath
std::vector<std::string> formats
template<typename T>
struct PolarGrid
#include <polar_grid.h>

Public Functions

inline PolarGrid()
inline PolarGrid(const std::vector<T> &angles, size_t nrows, size_t ncols, T gamma)
PolarGrid(const PolarGrid<T>&) = delete
PolarGrid<T> &operator=(const PolarGrid<T>&) = delete
inline PolarGrid(PolarGrid<T> &&other) noexcept
inline PolarGrid<T> &operator=(PolarGrid<T> &&other) noexcept
inline dims_t dims() const
inline size_t size() const
inline const std::vector<T> &angles() const
inline T angle(size_t i) const
inline size_t nprojs() const

Public Members

size_t npts
std::vector<T> theta
Array<T> x
Array<T> y
Array<T> z
template<typename T>
class RampPreconditioner
#include <precond.h>

Public Functions

inline RampPreconditioner(dims_t dims)
inline Array<T> apply(const Array<T> &input) const

Private Members

Array<T> filter_
struct ReconParams
#include <config.h>

Public Functions

ReconParams() = default
inline ReconParams(const toml::table &config)
inline void print(std::ostream &os) const

Public Members

Regularizer regularizer = Regularizer::UNCONSTRAINED
std::array<size_t, 3> recon_dims = {0, 0, 0}
size_t maxIters = 50
size_t innerIters = 3
float sigma = 1000.0f
float p = 1.2f
float lambda = 0.1f
float mu = 10.0f
float tol = 1e-5f
float xtol = 1e-5f
float PAD_FACTOR = 2.f
template<typename T>
class Slice
#include <slice.h>

Public Functions

inline Slice(T *data, dims_t dims)
inline dims_t dims() const
inline size_t size() const
inline T *begin()
inline T *end()
inline const T *begin() const
inline const T *end() const
inline T &operator[](dims_t idx)
inline const T &operator[](dims_t idx) const

Private Functions

inline size_t flatIdx(dims_t idx) const

Private Members

T *data_
dims_t dims_
size_t size_
template<typename T>
class Support
#include <support.h>

Public Functions

inline Support(dims_t dims, dims_t supp_dims)
inline void apply(Array<T> &data) const

Private Members

Array<T> support_
class Timer
#include <timer.h>

Public Functions

Timer() = default
inline void start()
inline void stop()
inline auto us() const
inline double ms() const
inline double seconds() const
inline void reset()

Private Members

std::chrono::duration<double> elapsed_time_
std::chrono::time_point<std::chrono::high_resolution_clock> start_time_
std::chrono::time_point<std::chrono::high_resolution_clock> end_time_
namespace tomocam

Typedefs

template<typename T>
using Dataset_t = std::tuple<Array<T>, std::vector<T>, T>

Enums

enum class Regularizer

Values:

enumerator qGGMRF
enumerator SPLIT_BREGMAN
enumerator UNCONSTRAINED
enum class LogMode

Values:

enumerator SILENT
enumerator STDOUT
enumerator FILE
enum class PadType

Values:

enumerator LEFT
enumerator RIGHT
enumerator SYMMETRIC

Functions

inline toml::table read_toml_file(const std::string &filepath)
template<typename T>
inline std::vector<T> read_angles_file(const std::string &filepath)
template<typename T>
std::vector<Dataset_t<T>> parse_input_datasets(const toml::table &config)
inline void dump_config(const std::string &filepath = "config.toml")
inline std::ostream &operator<<(std::ostream &outs, dims_t d)
template<typename T>
void apply_filter(Array<std::complex<T>> &arr, const std::string &filter_name)
template<typename T>
Array<T> mask_infs_nans(const Array<T> &projs)

Creates a mask identifying finite values in an array.

Returns a masked array where all the NaNs and Infs are set to 0.

Template Parameters:

T – Numeric type of the projection array elements

Parameters:

projs – Input array to check for finite values

Returns:

Array<T> Masked array with NaNs and Infs replaced by 0

template<typename T>
Array<T> mask_support(const Array<T> &projs, const dims_t dims, const std::vector<T> &theta, const T gamma)

Creates a mask for laminography/tomography support region.

Ptychographic projections are “infinite” in the sense that they add successively add more date to the projection view as the object is rotated. This extra data is missing from projections which are more normal to the beam. This function assumes that the object is contained with rectangular support.

Template Parameters:

T – Numeric type for coordinates and angles

Parameters:
  • projs – Input projection array

  • dims – Dimensions of the object being reconstructed (n1, n2, n3)

  • theta – Vector of tomographic rotation angles (in radians)

  • gamma – Laminography tilt angle (in radians)

Returns:

Array<T> with values outside the support replaced by 0

template<typename T>
auto pad1d(const Array<T> &arr, T factor, PadType pad_type) -> Array<T>
template<typename T>
auto crop1d(const Array<T> &arr, size_t crop_size, PadType pad_type) -> Array<T>
template<typename T>
Array<T> pad2d(const Array<T> &arr, T factor, PadType pad_type)
template<typename T>
Array<T> crop2d(const Array<T> &arr, dims_t new_dims, PadType pad_type)
template<typename T>
auto pad3d(const Array<T> &arr, T factor, PadType pad_type) -> Array<T>
template<typename T>
Array<T> crop3d(const Array<T> &arr, dims_t new_dims, PadType pad_type)
template<typename T>
inline std::array<T, 3> beam_dir_vector(T gamma, T alpha)
template<typename T>
Array<T> forward(const std::array<Array<T>, 3> &magnetization, const PolarGrid<T> &pg, const T &gamma)

Performs a forward projection of 3D magnetization data onto a polar grid.

Parameters:
  • magnetization – The input 3D magnetization components as std::array of 3 Arrays.

  • pg – The polar grid defining the projection geometry.

  • gamma – orientation of the polar grid.

Returns:

The projected data as an Array.

template<typename T>
std::array<Array<T>, 3> adjoint(const Array<T> &proj, const PolarGrid<T> &pg, const dims_t &recon_dims, T gamma)

Adjoint of the forward projection operator for polar grid data.

Parameters:
  • proj – The input polar grid data as an Array.

  • pg – The polar grid defining the projection geometry.

  • recon_dims – The dimensions of the output volume.

  • gamma – orientation of the polar grid.

Returns:

The adjoint projected data as std::array of 3 Arrays.

template<typename Float>
std::array<Array<Float>, 3> sysmat(const std::array<Array<Float>, 3> &x, const PolarGrid<Float> &grid, Float gamma)

Composition of backprojection and projection operators defined as a system matrix.

Parameters:
  • x – 3D vector field represented as an array of three components.

  • grid – The polar grid defining the projection geometry.

  • gamma – Sample orientation in plane normal to beam direction.

Returns:

\( A x = R^T R x \)

template<typename Float>
std::array<Array<Float>, 3> gradient(const std::array<Array<Float>, 3> &m, const std::array<Array<Float>, 3> &pT, const PolarGrid<Float> &grid, Float gamma)

Computes the gradient of the objective function for iterative reconstruction.

Parameters:
  • m – Current estimate of the 3-d vector field

  • pT – backprojection data as an Array.

  • grid – The polar grid defining the projection geometry.

  • gamma – Sample orientation in plane normal to beam direction.

Returns:

The gradient as an Array.

template<typename Float>
Float residual(const std::array<Array<Float>, 3> &m, const std::array<Array<Float>, 3> &pT, const PolarGrid<Float> &grid, Float pTp, Float gamma)

Computes the residual between the projected data and the measured data.

Parameters:
  • m – Current estimate of the 3-d vector field

  • pT – backprojected data as split into its components.

  • grid – The polar grid defining the projection geometry.

  • pTp – Precomputed inner product of the measured data.

  • gamma – Sample orientation in plane normal to beam direction.

Returns:

The computed residual as a scalar value of type T.

template<typename Float>
std::array<Array<Float>, 3> MBIR1(const std::vector<std::tuple<Array<Float>, std::vector<Float>, Float>> &datasets, const dims_t &recon_dims, const ReconParams &params)

Performs Model-Based Iterative Reconstruction (MBIR) for vector tomography using unconstrained Conjugate Gradient (CG) optimization.

Template Parameters:

T – Floating-point type (float or double).

Parameters:
  • datasets – A vector of tuples, where each tuple contains: 1). An Array representing projection data 2). A vector of representing the projection angles for the dataset 3). A orientation angle (gamma), a proxy for angle between the sample magnetization and the beam polarization direction.

  • recon_dims – The dimensions of the reconstructed volume as a dims_t object.

  • params – Reconstruction parameters including regularization type, and optimization parameters.

Returns:

A array of three components representing the reconstructed 3D vector field, cropped to the specified dimensions.

template<typename Float>
std::array<Array<Float>, 3> MBIR2(const std::vector<std::tuple<Array<Float>, std::vector<Float>, Float>> &datasets, const dims_t &recon_dims, const ReconParams &params)

Performs Model-Based Iterative Reconstruction (MBIR) for vector tomographic imaging.

Reconstructs 3D volumetric data from set of projection images using a regularized optimization. Using split-Bregman method for optimization.

Template Parameters:

T – Floating-point type (float or double).

Parameters:
  • datasets – A vector of tuples, where each tuple contains: 1). An Array representing projection data 2). A vector of representing the projection angles for the dataset 3). A orientation angle (gamma), a proxy for angle between the sample magnetization and the beam polarization direction.

  • recon_dims – The dimensions of the reconstructed volume as dims_t object.

  • params – Reconstruction parameters including regularization type, and optimization parameters.

Returns:

A array of three components representing the reconstructed 3D vector field, cropped to the specified dimensions.

template<typename Float>
std::array<Array<Float>, 3> MBIR3(const std::vector<std::tuple<Array<Float>, std::vector<Float>, Float>> &datasets, const dims_t &recon_dims, const ReconParams &params)

Performs Model-Based Iterative Reconstruction (MBIR) for vector tomography using NAG-style optimization + qGGMRF regularization.

Template Parameters:

T – Floating-point type (float or double).

Parameters:
  • datasets – A vector of tuples, where each tuple contains: 1). An Array representing projection data 2). A vector of representing the projection angles for the dataset 3). A orientation angle (gamma), a proxy for angle between the sample magnetization and the beam polarization direction.

  • recon_dims – The dimensions of the reconstructed volume as a dims_t object.

  • params – Reconstruction parameters including regularization type, and optimization parameters.

Returns:

A array of three components representing the reconstructed 3D vector field, cropped to the specified dimensions.

Variables

template<typename T> concept Float   = std::is_floating_point<T>::value

Restrict T to floating-point types.

namespace array

Functions

template<typename Real_t>
Array<std::complex<Real_t>> to_complex(const Array<Real_t> &a)
template<typename Real_t>
Array<Real_t> to_real(const Array<std::complex<Real_t>> &a)
template<typename T>
Array<T> abs(const Array<T> &a)
template<typename T>
T max(const Array<T> &a)
template<typename T>
T min(const Array<T> &a)
template<typename T>
T norm2(const Array<T> &a)
template<typename T>
T norm1(const Array<T> &a)
template<typename T>
T dot(const Array<T> &a, const Array<T> &b)
template<typename T>
void resetPads(Array<T> &a)
template<typename T>
Array<T> transpose(const Array<T> &a, std::array<size_t, 3> axes)

Variables

template<typename T> concept Real_t   = std::is_same_v<T, float> | std::is_same_v<T, double>
namespace fft

Functions

template<typename T>
Array<std::complex<T>> fft2(const Array<std::complex<T>> &input)
template<typename T>
Array<std::complex<T>> ifft2(const Array<std::complex<T>> &input)
template<typename T>
Array<std::complex<T>> fft2_r2c(const Array<T> &input)
template<typename T>
Array<T> fft2_c2r(const Array<std::complex<T>> &input, dims_t out_dims)
template<typename T>
Array<std::complex<T>> fft3_r2c(const Array<T> &input)
template<typename T>
Array<T> fft3_c2r(const Array<std::complex<T>> &input, dims_t out_dims)
template<typename T>
std::vector<T> fftfreq(size_t n_modes)
template<typename T>
std::vector<T> rfftfreq(size_t n_modes)
template<typename T>
Array<T> fftshift1(const Array<T> &input)
template<typename T>
Array<T> fftshift2(const Array<T> &input)
template<typename T>
Array<T> fftshift3(const Array<T> &input)
template<typename T>
Array<T> ifftshift1(const Array<T> &input)
template<typename T>
Array<T> ifftshift2(const Array<T> &input)
template<typename T>
Array<T> ifftshift3(const Array<T> &input)
namespace nufft

Functions

template<typename T>
void nufft3d1(const Array<std::complex<T>> &cz, Array<std::complex<T>> &fz, const PolarGrid<T> &pg)
template<typename T>
void nufft3d2(Array<std::complex<T>> &cz, const Array<std::complex<T>> &fz, const PolarGrid<T> &pg)
namespace plans

Variables

template<typename T>
FinufftPlanCache<T> cache
namespace opt

Typedefs

template<typename T>
using VecArray = std::array<Array<T>, 3>
template<typename T>
using Function = std::function<VecArray<T>(const VecArray<T>&)>
template<typename T>
using Residual = std::function<T(const VecArray<T>&)>

Functions

template<typename T>
std::array<Array<T>, 3> grad_u(const Array<T> &u)

Compute the gradient of an array.

Template Parameters:

T – The type of the array elements.

Parameters:

u – The input array.

Returns:

The gradient of the input array.

template<typename T>
Array<T> divergence(const std::array<Array<T>, 3> &u)

Compute the divergence of an array.

Template Parameters:

T – The type of the array elements.

Parameters:

u – The input array.

Returns:

The divergence of the input array.

template<typename T>
Array<T> laplacian(const Array<T> &u)

Compute the Laplacian of an array as divergence(grad_u(u)).

Template Parameters:

T – The type of the array elements.

Parameters:

u – The input array.

Returns:

The Laplacian of the input array.

template<typename T>
std::array<Array<T>, 3> demag(const std::array<Array<T>, 3> &m)
template<typename T>
VecArray<T> operator+(const VecArray<T> &a, const VecArray<T> &b)
template<typename T>
VecArray<T> &operator+=(VecArray<T> &a, const VecArray<T> &b)
template<typename T>
VecArray<T> operator-(const VecArray<T> &a, const VecArray<T> &b)
template<typename T>
VecArray<T> &operator-=(VecArray<T> &a, const VecArray<T> &b)
template<typename T>
VecArray<T> operator*(const VecArray<T> &a, T scalar)
template<typename T>
VecArray<T> split_bregman(const Function<T> &A, const VecArray<T> &yT, const VecArray<T> &x0, T lambda, T mu, size_t outer_max, size_t inner_max, T tol, T xtol)

Split Bregman method for L1-regularized optimization problems.

implement Split Bregman method for L1-regularized optimization problems

Parameters:
  • A – Function representing the R^T R operator for datasets 1, 2, … (i.e. A(u) = [R1^T R1 u, R2^T R2 u, …])

  • yT – 3-component arrays representing the R1^T y1 + R2^T y2 + …

  • x0 – Initial guess for the solution (std::array of 3 components)

  • lambda – Regularization parameter

  • mu – Bregman parameter

  • outer_max – Maximum number of outer iterations

  • inner_max – Maximum number of inner iterations

  • tol – Tolerance for convergence based on residual norm

  • xtol – Tolerance for convergence based on solution change

Returns:

Reconstructed solution vector

template<typename T>
VecArray<T> cgsolver(const Function<T> &A, const VecArray<T> &b, const VecArray<T> &x0, size_t max_iters, T tol, T xtol, T lambda)

Conjugate Gradient method for solving the linear system [A1, A2, …]^T x = [b1, b2, …]^T, where A1, A2 are self-adjoint operators.

Parameters:
  • A(u) – Function represeting the combined (A1 + A2 + … + μ∇^T∇)u

  • b – 3-component array (R^T y1 + * R^T y2 + … + μ∇^T(d- b))

  • x0 – Initial guess for the solution

  • max_iters – Maximum number of iterations

  • tol – Tolerance for convergence based on residual norm

  • xtol – Tolerance for convergence based on solution change (step norm)

  • lambda – Demagnetization constraint weight

Returns:

Approximate solution vector

template<typename T>
std::array<Array<T>, 3> nagopt(const std::function<std::array<Array<T>, 3>(const std::array<Array<T>, 3>&)> &grad, const std::function<T(const std::array<Array<T>, 3>&)> &loss, std::array<Array<T>, 3> &x, size_t max_iters, T lipschitz, T tol, T xtol, size_t max_inner_iters = 20, Logger *logger = nullptr)

Nesterov’s Optimal Gradient Method with Boyd’s momentum term (vector version)

Parameters:
  • grad – Gradient function for vector reconstruction

  • loss – Loss function for vector reconstruction

  • x – Initial solution (array of 3 components)

  • max_iters – Maximum number of iterations

  • lipschitz – Lipschitz constant of the gradient

  • tol – Tolerance for convergence based on loss change

  • xtol – Tolerance for convergence based on solution change

  • max_inner_iters – Maximum number of inner iterations for line search

  • loggerLogger instance for output control

Returns:

Optimized solution (array of 3 components)

template<typename T>
T lipschitz(const Function<T> &grad, const Array<T> &x0, size_t max_iters = 20, T tol = 1e-5)

Estimate the Lipschitz constant of a gradient function using the power iteration method.

Parameters:
  • function – representing the system matrix

  • x0 – reference for the size of the input

  • max_iters – Maximum number of iterations (default: 20)

  • tol – Tolerance for convergence (default: 1e-5)

Returns:

Estimated Lipschitz constant

template<typename T>
void qggmrf(Array<T> &g, const Array<T> &x, T sigma, T p)

Compute the gradient of the q-generalized Gaussian Markov Random Field (qGGMRF) prior.

Parameters:
  • x – Input array

  • dx – Output gradient array (updated in place)

  • sigma – Scale parameter

  • p – Exponent parameter

namespace tiff

Functions

inline Array<float> read(std::string filename)
inline void write(std::string filename, const Array<float> &data)
template<typename T>
inline void write_vectors(std::string filename, const std::array<Array<T>, 3> &data)
template<typename T>
void write3(const std::string &name, const std::array<Array<T>, 3> &data)
namespace vti

Functions

template<typename T>
void write_vectors(const std::string &filename, const std::array<Array<T>, 3> &vectors, const std::string &field_name = "magnetization")
file array.h
#include <algorithm>
#include <array>
#include <complex>
#include <cstdint>
#include <execution>
#include <memory>
#include <random>
#include <span>
#include <tuple>
#include <type_traits>
#include “dtypes.h
#include “slice.h
file array_ops.h
#include <algorithm>
#include <complex>
#include <execution>
#include <format>
#include <iostream>
#include <numeric>
#include <type_traits>
#include “array.h

Defines

ARRAY_OPS_H
file bregman.h
#include <array>
#include “array.h
file config.h
#include <cmath>
#include <filesystem>
#include <format>
#include <fstream>
#include <iostream>
#include <stdexcept>
#include <string>
#include <toml++/toml.h>
#include <tuple>
#include <vector>
#include “array.h
#include “mask.h
#include “tiff.h
file demag.h
#include <array>
#include <omp.h>
#include “array.h
#include “fft.h
#include “fftutils.h
file dtypes.h
#include <cstdint>
#include <ostream>
#include <stdexcept>
#include <tuple>
file fft.h
#include <complex>
#include <fftw3.h>
#include <new>
#include <type_traits>
#include “array.h
#include “dtypes.h
file fftutils.h
#include <algorithm>
#include <cstddef>
#include <cstdint>
#include <cstring>
#include <execution>
#include <stdexcept>
#include <vector>
#include “array.h
file filter.h
#include <complex>
#include <string>
#include “array.h
file finufft_plan.h
#include <complex>
#include <finufft.h>
#include <type_traits>
#include “polar_grid.h
file finufft_plan_cache.h
#include <array>
#include <mutex>
#include <stdexcept>
#include “finufft_plan.h
file logger.h
#include <fstream>
#include <iostream>
#include <memory>
#include <string>
file mask.h
#include <cmath>
#include <execution>
#include <stdexcept>
#include “array.h
#include <algorithm>
file nufft.h
#include <cmath>
#include <complex>
#include <vector>
#include “array.h
#include “dtypes.h
#include “finufft_plan.h
#include “finufft_plan_cache.h
#include “polar_grid.h
file optimize.h
#include <cmath>
#include <iostream>
#include “array.h
#include “array_ops.h
#include “logger.h
file padding.h
#include <algorithm>
#include <cmath>
#include <cstdint>
#include <system_error>
#include “array.h
#include “dtypes.h
file polar_grid.h
#include <cmath>
#include <cstddef>
#include <cstdint>
#include <iostream>
#include <vector>
#include “array.h
file precond.h
#include <algorithm>
#include <execution>
#include <vector>
#include “array.h
#include “array_ops.h
#include “fft.h
file projection.h
#include <array>
#include <cmath>
#include “array.h
#include “dtypes.h
#include “polar_grid.h
file slice.h
#include “dtypes.h
file support.h
#include <algorithm>
#include <execution>
#include “array.h
file tiff.h
#include <cstdint>
#include <iostream>
#include <tiffio.h>
#include “array.h
file timer.h
#include <chrono>
#include <iostream>

Variables

constexpr double US_PER_SECOND = 1.0e+06
constexpr double US_PER_MS = 1000.0
file tomocam.h
#include <vector>
#include “array.h
#include “config.h
#include “dtypes.h
#include “padding.h
#include “polar_grid.h
#include “projection.h
#include “tiff.h
#include “timer.h
#include “write_vti.h
file utils.h
#include <chrono>
#include <iostream>
#include <random>
file write_vti.h
#include <array>
#include <cstring>
#include <filesystem>
#include <format>
#include <fstream>
#include <stdexcept>
#include <string>
#include <vector>
#include “array.h
dir /home/docs/checkouts/readthedocs.org/user_builds/camera-lamino/checkouts/latest/include

Main Classes

template<typename T>
class Array

Public Functions

inline Array()
inline Array(size_t x, size_t y, size_t z)
inline Array(dims_t d)
Array(const Array<T>&) = delete
Array<T> &operator=(const Array<T>&) = delete
inline Array(Array<T> &&rhs) noexcept
inline Array<T> &operator=(Array<T> &&rhs) noexcept
inline Array<T> clone() const
inline T *begin()
inline const T *begin() const
inline T *end()
inline const T *end() const
inline dims_t dims() const
inline size_t size() const
inline size_t nslices() const
inline size_t nrows() const
inline size_t ncols() const
inline T &operator[](size_t i)
inline const T &operator[](size_t i) const
inline T at(size_t i, size_t j, size_t k) const
inline T &operator[](dims_t i)
inline const T &operator[](dims_t i) const
inline Slice<T> slice(size_t begin, size_t end)
inline const Slice<T> slice(size_t begin, size_t end) const
inline std::span<T> row(size_t i, size_t j)
inline const std::span<T> row(size_t i, size_t j) const
inline void setPads(dims_t p)
inline dims_t pads() const
inline Array<T> operator*(const Array<T> &v)
inline Array<T> &operator*=(T v)
inline Array<T> operator*(T v)
inline Array<T> &operator/=(T v)
inline Array<T> &operator/=(const Array<T> &v)
inline Array<T> operator/(const Array<T> &v) const
inline Array<T> operator/(T scalar) const
inline Array<T> &operator+=(T v)
inline Array<T> &operator+=(const Array<T> &v)
inline Array<T> operator+(const Array<T> &rhs) const
inline Array<T> &operator+=(const Array<T> &rhs) const
inline Array<T> &operator-=(T v)
inline Array<T> &operator-=(const Array<T> &v)
inline Array<T> operator-(const Array<T> &rhs) const

Public Static Functions

static inline Array<T> zeros(const dims_t &d)
static inline Array<T> ones(const dims_t &d)
static inline Array<T> random(const dims_t &d)

Warning

doxygenclass: Cannot find class “tomocam::PolarGrid” in doxygen xml output for project “Tomocam” from directory: ./doxygen/xml

Function Reference

Projection Functions

template<typename T>
Array<T> tomocam::forward(const std::array<Array<T>, 3> &magnetization, const PolarGrid<T> &pg, const T &gamma)

Performs a forward projection of 3D magnetization data onto a polar grid.

Parameters:
  • magnetization – The input 3D magnetization components as std::array of 3 Arrays.

  • pg – The polar grid defining the projection geometry.

  • gamma – orientation of the polar grid.

Returns:

The projected data as an Array.

Warning

doxygenfunction: Cannot find function “tomocam::backward” in doxygen xml output for project “Tomocam” from directory: ./doxygen/xml

Reconstruction Functions

Warning

doxygenfunction: Cannot find function “tomocam::MBIR” in doxygen xml output for project “Tomocam” from directory: ./doxygen/xml

Warning

doxygenfunction: Cannot find function “tomocam::conjugate_gradient” in doxygen xml output for project “Tomocam” from directory: ./doxygen/xml

Type Definitions

Warning

doxygentypedef: Cannot find typedef “tomocam::real_t” in doxygen xml output for project “Tomocam” from directory: ./doxygen/xml

Warning

doxygentypedef: Cannot find typedef “tomocam::complex_t” in doxygen xml output for project “Tomocam” from directory: ./doxygen/xml

Usage Example

Basic API usage:

#include "tomocam.h"

int main() {
    using namespace tomocam;

    // Create polar grid
    PolarGrid<float> grid(
        /* n_radial */ 256,
        /* n_angular */ 512,
        /* radius */ 1.0f
    );

    // Load projection data
    auto projections = read_tiff<float>("projections.tiff");
    auto angles = read_angles("angles.txt");

    // Perform MBIR reconstruction
    Array3D<float> recon_dims = {32, 256, 256};
    auto result = MBIR(
        projections,
        angles,
        recon_dims,
        /* max_iter */ 50,
        /* sigma */ 1000.0f,
        /* p */ 1.2f,
        /* lambda */ 0.1f,
        /* gamma */ 0.0f
    );

    // Save result
    write_tiff("reconstruction.tiff", result);

    return 0;
}

Build Configuration

To use Tomocam in your project:

CMake

find_package(Tomocam REQUIRED)
target_link_libraries(your_target PRIVATE tomocam::tomocam)

Compiler Flags

g++ -std=c++20 -I/path/to/tomocam/include \
    -L/path/to/tomocam/lib -ltomocam \
    -fopenmp -ltbb -lfftw3f -lfftw3 \
    your_code.cpp -o your_program

See Also