Expand description
§BestDose Algorithm
Bayesian dose optimization algorithm that finds optimal dosing regimens to achieve target drug concentrations or cumulative AUC (Area Under the Curve) values.
The BestDose algorithm combines Bayesian posterior estimation with dual optimization to balance patient-specific adaptation and population-level robustness.
§Quick Start
use pmcore::bestdose::{BestDoseProblem, Target, DoseRange};
// Create optimization problem
let problem = BestDoseProblem::new(
&population_theta, // Population support points from NPAG
&population_weights, // Population probabilities
Some(past_data), // Patient history (None = use prior)
target, // Future template with targets
None, // time_offset (None = standard mode)
eq, // PK/PD model
error_models, // Error specifications
DoseRange::new(0.0, 1000.0), // Dose constraints (0-1000 mg)
0.5, // bias_weight: 0=personalized, 1=population
settings, // NPAG settings
Target::Concentration, // Target type
)?;
// Run optimization
let result = problem.optimize()?;
// Extract results
println!("Optimal dose: {:?} mg", result.dose);
println!("Final cost: {}", result.objf);
println!("Method: {}", result.optimization_method); // "posterior" or "uniform"§Algorithm Overview (Three Stages)
┌─────────────────────────────────────────────────────────────────┐
│ STAGE 1: Posterior Density Calculation │
│ │
│ Prior (N points) │
│ ↓ │
│ Step 1.1: NPAGFULL11 - Bayesian Filtering │
│ Calculate P(data|θᵢ) for each support point │
│ Apply Bayes rule: P(θᵢ|data) ∝ P(data|θᵢ) × P(θᵢ) │
│ Filter: Keep points where P(θᵢ|data) > 1e-100 × max │
│ ↓ │
│ Filtered Posterior (M points, typically 5-50) │
│ ↓ │
│ Step 1.2: NPAGFULL - Local Refinement │
│ For each filtered point: │
│ Run full NPAG optimization │
│ Find refined "daughter" point │
│ ↓ │
│ Refined Posterior (M points with NPAGFULL11 weights) │
└─────────────────────────────────────────────────────────────────┘
┌─────────────────────────────────────────────────────────────────┐
│ STAGE 2: Dual Optimization │
│ │
│ Optimization 1: Posterior Weights (Patient-Specific) │
│ Minimize Cost = (1-λ)×Variance + λ×Bias² │
│ Using NPAGFULL11 posterior weights │
│ ↓ │
│ Result 1: (doses₁, cost₁) │
│ │
│ Optimization 2: Uniform Weights (Population) │
│ Minimize Cost = (1-λ)×Variance + λ×Bias² │
│ Using uniform weights (1/M for all points) │
│ ↓ │
│ Result 2: (doses₂, cost₂) │
│ │
│ Select Best: min(cost₁, cost₂) │
└─────────────────────────────────────────────────────────────────┘
┌─────────────────────────────────────────────────────────────────┐
│ STAGE 3: Final Predictions │
│ │
│ Calculate predictions with optimal doses │
│ For AUC targets: Use dense time grid + trapezoidal rule │
│ - AUCFromZero: Cumulative from time 0 │
│ - AUCFromLastDose: Interval from last dose │
│ Return: Optimal doses, cost, predictions, method used │
└─────────────────────────────────────────────────────────────────┘§Mathematical Foundation
§Bayesian Posterior
The posterior density is calculated via Bayes’ rule:
P(θ | data) = P(data | θ) × P(θ) / P(data)Where:
P(θ | data): Posterior (patient-specific parameters)P(data | θ): Likelihood (from error model)P(θ): Prior (from population)P(data): Normalizing constant
§Cost Function
The optimization minimizes a hybrid cost function:
Cost = (1-λ) × Variance + λ × Bias²Variance Term (Patient-Specific Performance):
Variance = Σᵢ P(θᵢ|data) × Σⱼ (target[j] - pred[i,j])²Expected squared error using posterior weights.
Bias Term (Population-Level Performance):
Bias² = Σⱼ (target[j] - E[pred[j]])²
where E[pred[j]] = Σᵢ P(θᵢ) × pred[i,j]Squared deviation from population mean prediction using prior weights.
Bias Weight Parameter (λ):
λ = 0.0: Fully personalized (minimize variance only)λ = 0.5: Balanced hybrid approachλ = 1.0: Population-based (minimize bias from population)
§Examples
§Single Dose Optimization
use pmcore::bestdose::{BestDoseProblem, Target, DoseRange};
use pharmsol::prelude::Subject;
// Define target: 5 mg/L at 24 hours
let target = Subject::builder("patient_001")
.bolus(0.0, 100.0, 0) // Initial dose (will be optimized)
.observation(24.0, 5.0, 0) // Target: 5 mg/L at 24h
.build();
let problem = BestDoseProblem::new(
&population_theta, &population_weights, Some(past), target, None,
eq, error_models,
DoseRange::new(10.0, 500.0), // 10-500 mg allowed
0.3, // Slight population emphasis
settings, Target::Concentration,
)?;
let result = problem.optimize()?;
println!("Optimal dose: {} mg", result.dose[0]);§Multiple Doses with AUC Target
use pmcore::bestdose::{BestDoseProblem, Target, DoseRange};
use pharmsol::prelude::Subject;
// Target: Achieve AUC₂₄ = 400 mg·h/L
let target = Subject::builder("patient_002")
.bolus(0.0, 100.0, 0) // Dose 1 (optimized)
.bolus(12.0, 100.0, 0) // Dose 2 (optimized)
.observation(24.0, 400.0, 0) // Target: AUC₂₄ = 400
.build();
let problem = BestDoseProblem::new(
&population_theta, &population_weights, Some(past), target, None,
eq, error_models,
DoseRange::new(50.0, 300.0),
0.0, // Full personalization
settings, Target::AUCFromZero, // Cumulative AUC target!
)?;
let result = problem.optimize()?;
println!("Dose 1: {} mg at t=0", result.dose[0]);
println!("Dose 2: {} mg at t=12", result.dose[1]);
if let Some(auc) = result.auc_predictions {
println!("Predicted AUC₂₄: {} mg·h/L", auc[0].1);
}§Population-Only Optimization
// No patient history - use population prior directly
let problem = BestDoseProblem::new(
&population_theta, &population_weights,
None, // No past data
target, None, // time_offset
eq, error_models,
DoseRange::new(0.0, 1000.0),
1.0, // Full population weighting
settings,
Target::Concentration,
)?;
let result = problem.optimize()?;
// Returns population-typical dose§Configuration
§Key Parameters
-
bias_weight(λ): Controls personalization level0.0: Minimize patient-specific variance (full personalization)1.0: Minimize deviation from population (robustness)
-
max_cycles: NPAGFULL refinement iterations0: Skip refinement (use filtered points directly)100-500: Typical range for refinement
-
doserange: Dose constraints- Set clinically appropriate bounds for your drug
-
target_type: Optimization targetTarget::Concentration: Direct concentration targetsTarget::AUCFromZero: Cumulative AUC from time 0Target::AUCFromLastDose: Interval AUC from last dose
§Performance Tuning
For faster optimization:
// Reduce refinement cycles
let problem = BestDoseProblem::new(
&population_theta, &population_weights, None, target, None,
eq, error_models,
DoseRange::new(0.0, 1000.0), 0.5,
settings.clone(),
Target::Concentration,
)?;
// For AUC: use coarser time grid
settings.predictions().idelta = 30.0; // 30-minute intervals§See Also
BestDoseProblem: Main entry point for optimizationBestDoseResult: Output structure with optimal dosesTarget: Enum for concentration vs AUC targetsDoseRange: Dose constraint specification
Modules§
- cost
- Cost function calculation for BestDose optimization
- predictions
- Stage 3: Prediction calculations
Structs§
- Best
Dose Problem - The BestDose optimization problem
- Best
Dose Result - Result from BestDose optimization
- Dose
Range - Allowable dose range constraints
Enums§
- Target
- Target type for dose optimization