Malnutrition Comorbidity Module

This module implements malnutrition disease modeling and its interactions with tuberculosis, including anthropometric measurements, nutritional interventions, and TB-malnutrition connectors.

Malnutrition Disease Model

Malnutrition Disease Model for TB Simulation

This module implements a malnutrition disease model that can be used in tuberculosis simulation studies. The model tracks anthropometric measurements (weight, height) using the LMS (Lambda-Mu-Sigma) method and simulates the effects of nutritional interventions on growth and development.

The model is designed to integrate with the RATIONS trial framework and supports both macronutrient and micronutrient supplementation interventions.

Mathematical Framework: - Uses Cole’s LMS method for growth reference curves - Implements random walk processes for weight percentile evolution - Applies nutritional intervention effects through drift parameters

Data Requirements: - Requires anthropometry.csv file with LMS parameters (L, M, S) by age and sex - Supports WHO growth standards and custom reference data

References: - https://www.ncbi.nlm.nih.gov/pmc/articles/PMC10876842/ - https://www.ncbi.nlm.nih.gov/pmc/articles/PMC9971264/ - https://www.espen.org/files/ESPEN-guidelines-on-definitions-and-terminology-of-clinical-nutrition.pdf

class tbsim.comorbidities.malnutrition.malnutrition.Malnutrition(pars=None, **kwargs)[source]

Bases: Disease

Malnutrition disease model for tuberculosis simulation studies.

This class implements a comprehensive malnutrition model that tracks anthropometric measurements using the LMS (Lambda-Mu-Sigma) method. It simulates the effects of nutritional interventions on growth and development, with support for both macronutrient and micronutrient supplementation.

Mathematical Model: - Weight percentile evolution: dW/dt = μ(t) + σ(t) * ε(t)

where μ(t) = drift from interventions, σ(t) = time-varying noise, ε(t) ~ N(0,1)

  • LMS transformation: X = M * (L*S*Z + 1)^(1/L) for L ≠ 0

    X = M * exp(S*Z) for L = 0

    where X = measurement, M = median, L = skewness, S = coefficient of variation, Z = z-score

Individual States: - receiving_macro (bool): Whether individual receives macronutrient supplementation - receiving_micro (bool): Whether individual receives micronutrient supplementation - height_percentile (float): Height percentile (0.0-1.0), assumed constant - weight_percentile (float): Weight percentile (0.0-1.0), evolves over time - micro (float): Micronutrient status z-score, evolves over time

LMS_data

Anthropometric reference data indexed by sex with columns: Age, Weight_L, Weight_M, Weight_S, Height_L, Height_M, Height_S

Type:

pd.DataFrame

dweight

Normal distribution for weight changes with location and scale functions

Type:

ss.normal

Parameters:
  • beta (float) – Transmission rate (placeholder, not used in malnutrition model)

  • init_prev (float) – Initial prevalence of malnutrition (default: 0.001)

File Dependencies:

anthropometry.csv: Must contain LMS parameters by age and sex for growth calculations

static dweight_loc(self, sim, uids)[source]

Calculate the location parameter (mean drift) for weight change distribution.

This method determines the mean drift in weight percentiles for individuals receiving macronutrient supplementation. The drift is proportional to the current time index to model cumulative intervention effects.

Mathematical Formula:

μ_i = 0 if not receiving macro supplementation μ_i = 1.0 * t_i if receiving macro supplementation where t_i is the current time index

Parameters:
  • sim (ss.Sim) – The simulation object containing current state

  • uids (np.ndarray) – Array of unique identifiers for individuals (int64)

Returns:

Mean drift values for weight change distribution (float64)

Shape: (len(uids),) Values: 0.0 for non-supplemented, 1.0*ti for supplemented

Return type:

np.ndarray

Implementation Details:
  • Creates zero array for all individuals

  • Sets positive drift only for individuals with receiving_macro=True

  • Drift magnitude scales linearly with simulation time

static dweight_scale(self, sim, uids)[source]

Calculate the scale parameter (standard deviation) for weight change distribution.

This method determines the standard deviation for weight changes over time. The scale increases linearly with time to model increasing variability in weight changes as individuals age and nutritional status becomes more variable.

Mathematical Formula:

σ_i = 0.01 * t_i for all individuals where t_i is the current time index

Parameters:
  • sim (ss.Sim) – The simulation object containing current state

  • uids (np.ndarray) – Array of unique identifiers for individuals (int64)

Returns:

Standard deviation values for weight change distribution (float64)

Shape: (len(uids),) Values: All elements equal to 0.01 * current_time_index

Return type:

np.ndarray

Implementation Details:
  • Uses np.full to create array with identical values

  • Scale parameter increases linearly with simulation time

  • Applied uniformly to all individuals regardless of intervention status

weight(uids=None)[source]

Calculate actual weight values (kg) from weight percentiles using LMS method.

Converts weight percentiles (0.0-1.0) to actual weight measurements in kilograms using the LMS transformation with age and sex-specific reference parameters.

Mathematical Formula:

W = M * (L*S*Z + 1)^(1/L) for L ≠ 0 W = M * exp(S*Z) for L = 0 where W = weight (kg), M = median weight, L = skewness, S = coefficient of variation Z = Φ^(-1)(percentile) is the inverse normal CDF of the percentile

Parameters:

uids (np.ndarray, optional) – Specific individuals to calculate weights for (int64) If None, calculates for all alive individuals

Returns:

Actual weight values in kilograms (float64)

Shape: (len(uids),) Range: Typically 2-100 kg depending on age and sex

Return type:

np.ndarray

Implementation Details:
  • Calls self.lms() with metric=’Weight’

  • Uses weight_percentile state variable

  • Interpolates LMS parameters by age and sex from reference data

height(uids=None)[source]

Calculate actual height values (cm) from height percentiles using LMS method.

Converts height percentiles (0.0-1.0) to actual height measurements in centimeters using the LMS transformation with age and sex-specific reference parameters.

Mathematical Formula:

H = M * (L*S*Z + 1)^(1/L) for L ≠ 0 H = M * exp(S*Z) for L = 0 where H = height (cm), M = median height, L = skewness, S = coefficient of variation Z = Φ^(-1)(percentile) is the inverse normal CDF of the percentile

Parameters:

uids (np.ndarray, optional) – Specific individuals to calculate heights for (int64) If None, calculates for all alive individuals

Returns:

Actual height values in centimeters (float64)

Shape: (len(uids),) Range: Typically 50-200 cm depending on age and sex

Return type:

np.ndarray

Implementation Details:
  • Calls self.lms() with metric=’Height’

  • Uses height_percentile state variable

  • Interpolates LMS parameters by age and sex from reference data

lms(percentile, uids=None, metric='Weight')[source]

Calculate anthropometric measurements using the LMS (Lambda-Mu-Sigma) method.

The LMS method is a statistical approach for constructing growth reference curves that accounts for the skewness of anthropometric data. It uses three parameters: lambda (skewness), mu (median), and sigma (coefficient of variation).

This method interpolates LMS parameters by age and sex, then converts percentiles to actual measurements using the inverse LMS transformation.

Mathematical Formula:

For each individual i: age_months_i = age_years_i * 12 L_i = interpolate(age_months_i, age_bins, L_values) M_i = interpolate(age_months_i, age_bins, M_values) S_i = interpolate(age_months_i, age_bins, S_values) Z_i = Φ^(-1)(percentile_i) # inverse normal CDF X_i = M_i * (L_i*S_i*Z_i + 1)^(1/L_i) if L_i ≠ 0 X_i = M_i * exp(S_i*Z_i) if L_i = 0 where X_i is the measurement (weight in kg or height in cm)

Parameters:
  • percentile (np.ndarray) – Percentile values (0.0-1.0) for the measurements (float64)

  • uids (np.ndarray, optional) – Specific individuals to calculate for (int64) If None, uses all alive individuals

  • metric (str) – Type of measurement to calculate. Must be one of: ‘Weight’ (kg), ‘Height’ (cm), ‘Length’ (cm), ‘BMI’ (kg/m²)

Returns:

Actual anthropometric measurements (float64)

Shape: (len(uids),) Units: kg for Weight, cm for Height/Length, kg/m² for BMI

Return type:

np.ndarray

Raises:

AssertionError – If metric is not one of [‘Weight’, ‘Height’, ‘Length’, ‘BMI’]

Implementation Details:
  • Processes males and females separately due to different reference data

  • Uses linear interpolation (np.interp) for age-specific parameters

  • Converts age from years to months for reference data lookup

  • Handles edge case where lambda parameter equals zero

  • Uses scipy.stats.norm().ppf() for percentile to z-score conversion

References

__init__(pars=None, **kwargs)[source]

Initialize the Malnutrition disease model.

Sets up the model parameters, loads anthropometric reference data, and defines the disease states for tracking individual nutritional status.

Initialization Process: 1. Calls parent class constructor (ss.Disease.__init__) 2. Defines model parameters with default values 3. Loads LMS reference data from anthropometry.csv 4. Defines individual state variables 5. Creates weight change distribution function

Parameters:
  • pars (dict, optional) – Dictionary of parameters to override defaults Keys: ‘beta’, ‘init_prev’ Values: float

  • **kwargs – Additional keyword arguments passed to parent class

  • beta (float) – Transmission rate (placeholder, not used in malnutrition model) Default: 1.0

  • init_prev (float) – Initial prevalence of malnutrition (0.0-1.0) Default: 0.001 (0.1%)

State Variables Created:

receiving_macro (bool): Macronutrient supplementation status receiving_micro (bool): Micronutrient supplementation status height_percentile (float): Height percentile (0.0-1.0) weight_percentile (float): Weight percentile (0.0-1.0) micro (float): Micronutrient status z-score

File Dependencies:
anthropometry.csv: Must be in tbsim/data/ directory with columns:

Sex, Age, Weight_L, Weight_M, Weight_S, Height_L, Height_M, Height_S

set_initial_states(sim)[source]

Set initial values for disease states during simulation initialization.

This method is called during simulation initialization to set up the initial nutritional status of individuals. Currently a placeholder for future implementation of correlated weight and height percentiles.

Current Implementation:
  • No action taken (placeholder method)

  • All state variables initialized with defaults from define_states()

Future Implementation Notes:
  • Could implement correlation between weight and height percentiles

  • Potential approach: bivariate normal distribution with correlation ρ

  • Corner corrections needed for percentiles near 0.0 or 1.0

  • Formula: (W_p, H_p) ~ BVN(μ, Σ) where Σ = [[1, ρ], [ρ, 1]]

Parameters:

sim (ss.Sim) – The simulation object containing population data

Returns:

None

Implementation Details:
  • Called once during simulation setup

  • Runs after individual creation but before first time step

  • Could modify height_percentile, weight_percentile, micro states

step()[source]

Execute one time step of the malnutrition model.

This method is called at each simulation time step to update the nutritional status of individuals. It implements random walk processes for weight percentiles and applies nutritional interventions.

Mathematical Model:

For each individual i at time t: ΔW_i(t) ~ N(μ_i(t), σ_i(t)²) W_i(t+1) = W_i(t) + ΔW_i(t) W_i(t+1) = clip(W_i(t+1), 0.025, 0.975) where W_i(t) is weight percentile, μ_i(t) is drift from dweight_loc(), σ_i(t) is scale from dweight_scale(), and clip() ensures valid percentiles

Step Process: 1. Get all alive individual IDs 2. Sample weight changes from normal distribution 3. Update weight percentiles with random walk 4. Clip percentiles to valid range (0.025-0.975) 5. Apply intervention effects through drift parameters

Parameters:

None (uses self.sim for simulation state)

Returns:

None (modifies state variables in-place)

Implementation Details:
  • Uses self.dweight() which calls dweight_loc() and dweight_scale()

  • Clipping prevents percentiles from reaching exact 0.0 or 1.0

  • Random walk process models natural weight variability

  • Intervention effects applied through location parameter

init_results()[source]

Initialize results tracking for the malnutrition model.

Sets up the results structure to track key metrics during simulation. Currently tracks the proportion of people alive at each time step.

Results Structure:
people_alive (float): Proportion of population alive at each time step

Shape: (n_timesteps,) Range: 0.0-1.0 Units: Proportion (dimensionless)

Parameters:

None (uses self.sim for simulation parameters)

Returns:

None (creates self.results object)

Implementation Details:
  • Calls parent class init_results() first

  • Defines results using ss.Result objects

  • Results stored in self.results dictionary

  • Updated each time step in update_results()

update_results()[source]

Update results at each time step.

Records the current state of the simulation for analysis. Currently tracks the proportion of individuals who are alive at each time step.

Mathematical Formula:

people_alive[ti] = count(alive_individuals) / total_population where ti is the current time index, alive_individuals is boolean array indicating survival status, and total_population is n_agents parameter

Parameters:

None (uses self.sim for simulation state)

Returns:

None (updates self.results in-place)

Implementation Details:
  • Called at each simulation time step

  • Uses self.sim.people.alive boolean array

  • Calculates proportion as count(alive)/n_agents

  • Stores result in self.results.people_alive[ti]

  • ti is current time index from self.sim.ti

TB-Malnutrition Connector

TB-Malnutrition Connector Module

This module implements a connector between Tuberculosis (TB) and Malnutrition disease models in the simulation framework. The connector defines how nutritional status affects TB dynamics through various risk ratios and susceptibility modifiers.

The connector implements three main interaction mechanisms: 1. Activation risk ratio: How malnutrition affects TB activation from latent to active 2. Clearance risk ratio: How malnutrition affects TB clearance/recovery rates 3. Relative susceptibility: How malnutrition affects susceptibility to new TB infection

Mathematical Framework: - Risk ratios (RR) modify disease transition rates multiplicatively - Relative susceptibility modifies infection probability for uninfected individuals - BMI-based risk functions use sigmoid transformations of log-linear relationships

Examples:

Basic usage: Create connector with default functions and add to simulation. Custom configurations: Use BMI-based risk ratios, supplementation effects, or combined functions. Analysis: Access risk ratios and susceptibility modifiers after simulation runs.

See the module documentation and method docstrings for detailed usage examples.

References: - Lönnroth et al. studies on BMI and TB risk - Nutritional supplementation trials and their effects on TB outcomes

class tbsim.comorbidities.malnutrition.tb_malnut_cnn.TB_Nutrition_Connector(pars=None, **kwargs)[source]

Bases: Connector

Connector between Tuberculosis and Malnutrition disease models.

This connector implements the bidirectional interactions between TB and malnutrition, where nutritional status affects TB dynamics and TB infection may affect nutritional status. The connector modifies disease transition rates through risk ratios and susceptibility modifiers.

Mathematical Model: - Activation rate modification: λ_act(t) = λ_act_base * RR_activation(t) - Clearance rate modification: λ_clear(t) = λ_clear_base * RR_clearance(t) - Susceptibility modification: P_inf(t) = P_inf_base * rel_sus(t)

where λ = rates, RR = risk ratios, P_inf = infection probability

Interaction Mechanisms: 1. Supplementation effects: Reduced risk ratios for individuals receiving nutritional interventions 2. BMI-based risk: Sigmoid function of BMI following Lönnroth et al. log-linear relationship 3. Micronutrient effects: Increased susceptibility for individuals with low micronutrient status

Parameters:
  • rr_activation_func (callable) – Function to compute activation risk ratios

  • rr_clearance_func (callable) – Function to compute clearance risk ratios

  • relsus_func (callable) – Function to compute relative susceptibility modifiers

sim

Reference to the simulation object

Type:

ss.Sim

pars

Parameter dictionary containing function references

Type:

dict

__init__(pars=None, **kwargs)[source]

Initialize the TB-Malnutrition connector.

Sets up the connector with default risk ratio and susceptibility functions, and configures the interaction parameters between TB and malnutrition models.

Initialization Process: 1. Calls parent class constructor with label ‘TB-Malnutrition’ 2. Defines parameter functions for risk ratios and susceptibility 3. Updates parameters with any provided overrides

Parameters:
  • pars (dict, optional) – Dictionary of parameters to override defaults Keys: ‘rr_activation_func’, ‘rr_clearance_func’, ‘relsus_func’ Values: callable functions

  • **kwargs – Additional keyword arguments passed to parent class

Default Functions:

rr_activation_func: ones_rr (no effect on activation) rr_clearance_func: ones_rr (no effect on clearance) relsus_func: compute_relsus (micronutrient-based susceptibility)

static supplementation_rr(tb, mn, uids, rate_ratio=0.5)[source]

Calculate risk ratios based on nutritional supplementation status.

This function reduces TB activation and clearance rates for individuals receiving both macronutrient and micronutrient supplementation, modeling the protective effects of comprehensive nutritional interventions.

Mathematical Formula:

RR_i = 1.0 if not receiving both macro and micro supplementation RR_i = rate_ratio if receiving both macro and micro supplementation where rate_ratio < 1.0 indicates reduced risk (protective effect)

Parameters:
  • tb (TB) – Tuberculosis disease model object

  • mn (Malnutrition) – Malnutrition disease model object

  • uids (np.ndarray) – Array of individual identifiers (int64)

  • rate_ratio (float) – Risk ratio for supplemented individuals (default: 0.5) Range: 0.0-1.0, where 0.5 = 50% risk reduction

Returns:

Risk ratios for each individual (float64)

Shape: (len(uids),) Values: 1.0 for non-supplemented, rate_ratio for supplemented

Return type:

np.ndarray

Implementation Details:
  • Creates array of ones for all individuals

  • Identifies individuals receiving both macro and micro supplementation

  • Applies rate_ratio only to fully supplemented individuals

  • Uses boolean indexing with logical AND operation

static lonnroth_bmi_rr(tb, mn, uids, scale=2, slope=3, bmi50=25)[source]

Calculate risk ratios based on BMI using Lönnroth et al. relationship.

This function implements a sigmoid transformation of the log-linear relationship between BMI and TB risk described by Lönnroth et al. The function creates a smooth transition around a reference BMI value with configurable steepness.

Mathematical Formula:

BMI_i = 10,000 * weight_i(kg) / height_i(cm)² x_i = -0.05 * (BMI_i - 15) + 2 # Log-linear relationship from Lönnroth et al. x0 = -0.05 * (bmi50 - 15) + 2 # Center point at reference BMI RR_i = scale / (1 + 10^(-slope * (x_i - x0)))

where: - BMI_i is calculated from weight and height measurements - x_i is the log-linear predictor from Lönnroth et al. - x0 centers the sigmoid at the reference BMI - scale controls the maximum risk ratio - slope controls the steepness of the sigmoid transition

Parameters:
  • tb (TB) – Tuberculosis disease model object

  • mn (Malnutrition) – Malnutrition disease model object

  • uids (np.ndarray) – Array of individual identifiers (int64)

  • scale (float) – Maximum risk ratio value (default: 2.0) Range: > 0, typically 1.0-5.0

  • slope (float) – Steepness of sigmoid transition (default: 3.0) Range: > 0, higher values = steeper transition

  • bmi50 (float) – Reference BMI for sigmoid center (default: 25.0 kg/m²) Range: 15-35 kg/m², typical healthy adult range

Returns:

Risk ratios based on BMI (float64)

Shape: (len(uids),) Range: 0.0 to scale, with sigmoid transition around bmi50

Return type:

np.ndarray

Implementation Details:
  • Calculates BMI using weight (kg) and height (cm) from malnutrition model

  • Applies Lönnroth et al. log-linear transformation

  • Centers sigmoid function at specified reference BMI

  • Uses 10-based logarithm for sigmoid calculation

  • Returns risk ratios where higher values indicate increased risk

References

  • Lönnroth et al. studies on BMI and TB risk relationships

  • Log-linear model: log(incidence) = -0.05*(BMI-15) + 2

static ones_rr(tb, mn, uids)[source]

Return neutral risk ratios (no effect on disease dynamics).

This function serves as a neutral baseline that applies no modification to TB activation or clearance rates. It is used as a default function when no nutritional effects on TB dynamics are desired.

Mathematical Formula:

RR_i = 1.0 for all individuals i This means: λ_modified = λ_base * 1.0 = λ_base (no change)

Parameters:
  • tb (TB) – Tuberculosis disease model object

  • mn (Malnutrition) – Malnutrition disease model object

  • uids (np.ndarray) – Array of individual identifiers (int64)

Returns:

Neutral risk ratios of 1.0 for all individuals (float64)

Shape: (len(uids),) Values: All elements equal to 1.0

Return type:

np.ndarray

Implementation Details:
  • Creates array of ones with same shape as uids

  • Uses np.ones_like() for efficient array creation

  • Serves as identity function for risk ratio calculations

static compute_relsus(tb, mn, uids)[source]

Calculate relative susceptibility based on micronutrient status.

This function modifies the susceptibility to new TB infection based on individual micronutrient status. Individuals with low micronutrient levels experience increased susceptibility to TB infection.

Mathematical Formula:

rel_sus_i = 1.0 if micro_i ≥ 0.2 (normal micronutrient status) rel_sus_i = 2.0 if micro_i < 0.2 (low micronutrient status) where micro_i is the micronutrient z-score from malnutrition model

Threshold Logic:
  • micro_i ≥ 0.2: Normal susceptibility (rel_sus_i = 1.0)

  • micro_i < 0.2: Doubled susceptibility (rel_sus_i = 2.0)

  • Threshold of 0.2 represents approximately 42nd percentile of normal distribution

Parameters:
  • tb (TB) – Tuberculosis disease model object

  • mn (Malnutrition) – Malnutrition disease model object

  • uids (np.ndarray) – Array of individual identifiers (int64)

Returns:

Relative susceptibility modifiers (float64)

Shape: (len(uids),) Values: 1.0 for normal micronutrient status, 2.0 for low status

Return type:

np.ndarray

Implementation Details:
  • Accesses micronutrient status from malnutrition model

  • Applies threshold-based logic with 0.2 z-score cutoff

  • Uses boolean indexing for efficient conditional assignment

  • Returns susceptibility multipliers where higher values = increased risk

step()[source]

Execute one time step of TB-Malnutrition interactions.

This method is called at each simulation time step to apply the nutritional effects on TB dynamics. It modifies TB transition rates and susceptibility based on current nutritional status of individuals.

Mathematical Model:

For infected individuals (latent TB): - RR_activation(t) = RR_activation_base * rr_activation_func(t) - RR_clearance(t) = RR_clearance_base * rr_clearance_func(t)

For uninfected individuals: - rel_sus(t) = relsus_func(t)

Step Process: 1. Get references to TB and malnutrition disease models 2. For infected individuals: modify activation and clearance risk ratios 3. For uninfected individuals: update relative susceptibility 4. Apply multiplicative modifications to existing rates

Parameters:

None (uses self.sim for simulation state)

Returns:

None (modifies TB model state variables in-place)

Implementation Details:
  • Accesses disease models through self.sim.diseases dictionary

  • Processes infected and uninfected individuals separately

  • Uses multiplicative updates ( *= ) to combine multiple effects

  • Modifies tb.rr_activation, tb.rr_clearance, and tb.rel_sus arrays

  • Risk ratios start at 1.0 each time step and are modified by connector

Model Overview

The malnutrition comorbidity module provides comprehensive modeling of nutritional status and its bidirectional interactions with tuberculosis:

Anthropometric Modeling
  • Weight and Height: LMS (Lambda-Mu-Sigma) method for growth reference curves

  • Percentile Tracking: Individual weight and height percentiles (0.0-1.0)

  • Micronutrient Status: Z-score based micronutrient tracking

  • Age and Sex Specific: Reference data from WHO growth standards

Nutritional Interventions
  • Macronutrient Supplementation: Affects weight percentile evolution

  • Micronutrient Supplementation: Influences susceptibility to TB

  • RATIONS Trial Integration: Framework for nutritional intervention studies

  • Dynamic Effects: Time-dependent intervention impacts

TB-Malnutrition Interactions
  • Malnutrition modifies TB activation and clearance rates

  • BMI-based risk functions using Lönnroth et al. relationships

  • Supplementation effects on TB dynamics

  • Bidirectional disease interactions

Key Features

LMS Growth Modeling
  • Cole’s LMS method for anthropometric measurements

  • Age and sex-specific reference parameters

  • Percentile to measurement conversion

  • Growth curve interpolation

Nutritional Status Tracking
  • Weight percentile evolution with random walk processes

  • Height percentile (assumed constant)

  • Micronutrient z-score tracking

  • Intervention effect modeling

Intervention Framework
  • Macronutrient supplementation effects on weight

  • Micronutrient supplementation effects on susceptibility

  • Time-dependent intervention impacts

  • Individual-level intervention assignment

TB Integration
  • Risk ratio calculations for TB progression

  • BMI-based risk functions

  • Supplementation protective effects

  • Real-time parameter adjustment

Usage Examples

Basic malnutrition simulation:

from tbsim.comorbidities.malnutrition.malnutrition import Malnutrition
import starsim as ss

sim = ss.Sim(diseases=Malnutrition())
sim.run()

With custom parameters:

from tbsim.comorbidities.malnutrition.malnutrition import Malnutrition

malnutrition = Malnutrition(pars={
    'init_prev': 0.05,     # 5% initial malnutrition prevalence
    'beta': 0.8            # Custom transmission rate
})

sim = ss.Sim(diseases=malnutrition)
sim.run()

TB-malnutrition integration:

from tbsim.comorbidities.malnutrition.tb_malnut_cnn import TB_Nutrition_Connector
from tbsim import TB, Malnutrition

# Add both diseases and connector for interactions
connector = TB_Nutrition_Connector(pars={
    'rr_activation_func': TB_Nutrition_Connector.lonnroth_bmi_rr,
    'rr_clearance_func': TB_Nutrition_Connector.supplementation_rr,
    'relsus_func': TB_Nutrition_Connector.compute_relsus
})

sim = ss.Sim(
    diseases=[TB(), Malnutrition()],
    connectors=connector
)
sim.run()

Accessing anthropometric data:

# Get current measurements
weight_kg = malnutrition.weight()           # Weight in kilograms
height_cm = malnutrition.height()           # Height in centimeters

# Get percentiles
weight_percentiles = malnutrition.weight_percentile
height_percentiles = malnutrition.height_percentile
micronutrient_status = malnutrition.micro

# Get intervention status
receiving_macro = malnutrition.receiving_macro
receiving_micro = malnutrition.receiving_micro

Mathematical Framework

LMS Transformation
  • X = M × (L×S×Z + 1)^(1/L) for L ≠ 0

  • X = M × exp(S×Z) for L = 0

  • Where X = measurement, M = median, L = skewness, S = coefficient of variation

  • Z = Φ^(-1)(percentile) is the inverse normal CDF

Weight Evolution
  • Random walk process: ΔW_i(t) ~ N(μ_i(t), σ_i(t)²)

  • Drift function: μ_i(t) = 1.0 × t for supplemented individuals

  • Scale function: σ_i(t) = 0.01 × t for all individuals

  • Percentile clipping: W_i(t+1) = clip(W_i(t+1), 0.025, 0.975)

TB Risk Modification
  • BMI calculation: BMI = 10,000 × weight(kg) / height(cm)²

  • Lönnroth relationship: log(incidence) = -0.05×(BMI-15) + 2

  • Sigmoid transformation: RR = scale / (1 + 10^(-slope × (x-x0)))

  • Supplementation effects: RR = 1.0 for non-supplemented, rate_ratio for supplemented

Susceptibility Modification
  • Relative susceptibility based on micronutrient status

  • Threshold-based logic: rel_sus = 1.0 if micro ≥ 0.2, 2.0 if micro < 0.2

  • Z-score threshold represents approximately 42nd percentile

Data Requirements

Anthropometry Reference Data
  • File: tbsim/data/anthropometry.csv

  • Columns: Sex, Age, Weight_L, Weight_M, Weight_S, Height_L, Height_M, Height_S

  • Age in months, LMS parameters for weight and height

  • WHO growth standards or custom reference data

Intervention Parameters
  • Macronutrient supplementation timing and duration

  • Micronutrient supplementation protocols

  • RATIONS trial specific parameters

  • Age and eligibility criteria

For detailed information about specific methods and parameters, see the individual class documentation above. All methods include comprehensive mathematical models and implementation details in their docstrings.