Program Listing for File multiConditionProblem.h

Return to documentation for file (include/parpeamici/multiConditionProblem.h)

#ifndef PARPE_AMICI_MULTI_CONDITION_PROBLEM_H
#define PARPE_AMICI_MULTI_CONDITION_PROBLEM_H

#include <parpecommon/parpeConfig.h>
#include <parpeoptimization/multiStartOptimization.h>
#include <parpeoptimization/optimizationProblem.h>
#include <parpeamici/amiciSimulationRunner.h>
#include <parpeoptimization/minibatchOptimization.h>

#include <amici/amici.h>
#include <amici/serialization.h>

#include <boost/serialization/map.hpp>

#include <gsl/gsl-lite.hpp>

#include <memory>

namespace parpe {

class LoadBalancerMaster;
class OptimizationResultWriter;
class MultiConditionDataProviderHDF5;
class MultiConditionDataProvider;
AmiciSimulationRunner::AmiciResultPackageSimple runAndLogSimulation(amici::Solver const &solver,
        amici::Model &model,
        int conditionIdx,
        int jobId,
        const MultiConditionDataProvider *dataProvider,
        OptimizationResultWriter *resultWriter,
        bool logLineSearch,
        Logger *logger,
        bool sendStates = false);

FunctionEvaluationStatus getModelOutputsAndSigmas(
        MultiConditionDataProvider *dataProvider,
        LoadBalancerMaster *loadBalancer,
        int maxSimulationsPerPackage,
        OptimizationResultWriter *resultWriter,
        bool logLineSearch,
        gsl::span<const double> parameters,
        std::vector<std::vector<double> > &modelOutputs,
        std::vector<std::vector<double> > &modelSigmas,
        Logger *logger, double *cpuTime, bool sendStates);

void messageHandler(MultiConditionDataProvider *dataProvider,
                    OptimizationResultWriter *resultWriter,
                    bool logLineSearch,
                    std::vector<char> &buffer, int jobId, bool sendStates);

class AmiciSummedGradientFunction : public SummedGradientFunction<int> {

public:
    using WorkPackage = AmiciSimulationRunner::AmiciWorkPackageSimple;
    using ResultPackage = AmiciSimulationRunner::AmiciResultPackageSimple;
    using ResultMap = std::map<int, ResultPackage>;

    AmiciSummedGradientFunction(
            MultiConditionDataProvider *dataProvider,
            LoadBalancerMaster *loadBalancer,
            OptimizationResultWriter *resultWriter);

    ~AmiciSummedGradientFunction() override = default;

    FunctionEvaluationStatus evaluate(
            gsl::span<const double> parameters,
            int dataset,
            double &fval,
            gsl::span<double> gradient,
            Logger *logger,
            double *cpuTime) const override;

    FunctionEvaluationStatus evaluate(
            gsl::span<const double> parameters,
            std::vector<int> datasets,
            double &fval,
            gsl::span<double> gradient,
            Logger *logger,
            double *cpuTime) const override;

    [[nodiscard]] int numParameters() const override;

    [[nodiscard]] std::vector<std::string> getParameterIds() const override;

    virtual FunctionEvaluationStatus getModelOutputsAndSigmas(
        gsl::span<double const> parameters,
        std::vector<std::vector<double> > &modelOutputs,
        std::vector<std::vector<double> > &modelSigmas,
        Logger *logger,
        double *cpuTime) const;

    [[nodiscard]] virtual std::vector<std::vector<double>> getAllMeasurements() const;

    virtual void messageHandler(std::vector<char> &buffer, int jobId) const;

    [[nodiscard]] virtual amici::ParameterScaling getParameterScaling(int parameterIndex) const;

    bool sendStates = false;

protected:// for testing
    AmiciSummedGradientFunction() = default;

    virtual int runSimulations(gsl::span<double const> optimizationParameters,
                               double &nllh,
                               gsl::span<double> objectiveFunctionGradient,
                               const std::vector<int> &dataIndices,
                               Logger *logger,
                               double *cpuTime) const;

    int aggregateLikelihood(JobData &data, double &negLogLikelihood,
                            gsl::span<double> negLogLikelihoodGradient,
                            double &simulationTimeInS,
                            gsl::span<const double> optimizationParameters
                            ) const;


    void addSimulationGradientToObjectiveFunctionGradient(
            int conditionIdx,
            gsl::span<const double> simulationGradient,
            gsl::span<double> objectiveFunctionGradient,
            gsl::span<const double> parameters) const;

    void setSensitivityOptions(bool sensiRequired) const;

private:
    // TODO: make owning
    MultiConditionDataProvider *dataProvider = nullptr;
    // Non-owning
    LoadBalancerMaster *loadBalancer = nullptr;
    std::unique_ptr<amici::Model> model;
    std::unique_ptr<amici::Solver> solver;
    std::unique_ptr<amici::Solver> solverOriginal;
    OptimizationResultWriter *resultWriter = nullptr; // TODO: owning?
    bool logLineSearch = false;
    int maxSimulationsPerPackage = 8;
    int maxGradientSimulationsPerPackage = 1;
};



class MultiConditionProblem
        : public MinibatchOptimizationProblem<int>
{

  public:
    MultiConditionProblem() = default;

    explicit MultiConditionProblem(MultiConditionDataProvider *dp);

    MultiConditionProblem(
            MultiConditionDataProvider *dp,
            LoadBalancerMaster *loadBalancer,
            std::unique_ptr<Logger> logger,
            std::unique_ptr<OptimizationResultWriter> resultWriter);

    ~MultiConditionProblem() override = default;

    virtual int earlyStopping();

    MultiConditionDataProvider *getDataProvider();
    OptimizationResultWriter *getResultWriter();

    //    virtual std::unique_ptr<double[]> getInitialParameters(int multiStartIndex) const override;

    void setInitialParameters(const std::vector<double> &startingPoint);
    void setParametersMin(const std::vector<double> &lowerBounds);
    void setParametersMax(std::vector<double> const& upperBounds);

    void fillParametersMin(gsl::span<double> buffer) const override;
    void fillParametersMax(gsl::span<double> buffer) const override;
    void fillInitialParameters(gsl::span<double> buffer) const override;

    std::unique_ptr<OptimizationReporter> getReporter() const override;

    std::vector<int> getTrainingData() const override;

private:
    //TODO std::unique_ptr<OptimizationProblem> validationProblem;

    MultiConditionDataProvider *dataProvider = nullptr;

    std::unique_ptr<OptimizationResultWriter> resultWriter;

    std::vector<double> startingPoint;
    std::vector<double> parametersMin;
    std::vector<double> parametersMax;
};




class MultiConditionProblemMultiStartOptimizationProblem
    : public MultiStartOptimizationProblem {
  public:
    MultiConditionProblemMultiStartOptimizationProblem(
            MultiConditionDataProviderHDF5 *dp,
            OptimizationOptions options,
            OptimizationResultWriter *resultWriter,
            LoadBalancerMaster *loadBalancer,
            std::unique_ptr<Logger> logger);


    int getNumberOfStarts() const override;

    bool restartOnFailure() const override;

    std::unique_ptr<OptimizationProblem> getLocalProblem(
            int multiStartIndex) const override;

private:
    MultiConditionDataProviderHDF5 *data_provider_ = nullptr;
    OptimizationOptions options_;
    OptimizationResultWriter *result_writer_ = nullptr;
    LoadBalancerMaster *load_balancer_ = nullptr;
    std::unique_ptr<Logger> logger_;
};


void saveSimulation(H5::H5File const& file, const std::string &pathStr,
        const std::vector<double> &parameters, double llh,
        gsl::span<const double> gradient, double timeElapsedInSeconds, gsl::span<const double>, gsl::span<const double>, gsl::span<const double>,
        int jobId, int status, const std::string &label);


} // namespace parpe

#endif