File: OptimizerMultiStageDPBase.h

package info (click to toggle)
stopt 5.12%2Bdfsg-3
  • links: PTS, VCS
  • area: main
  • in suites: trixie
  • size: 8,860 kB
  • sloc: cpp: 70,456; python: 5,950; makefile: 72; sh: 57
file content (71 lines) | stat: -rw-r--r-- 3,166 bytes parent folder | download | duplicates (2)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
// Copyright (C) 2023 EDF
// All Rights Reserved
// This code is published under the GNU Lesser General Public License (GNU LGPL)
#ifndef OPTIMIZERMULTISTAGEDPBASE_H
#define OPTIMIZERMULTISTAGEDPBASE_H
#include <Eigen/Dense>
#include "StOpt/core/utils/StateWithStocks.h"
#include "StOpt/core/grids/SpaceGrid.h"
#include "StOpt/dp/OptimizerBase.h"
#include "StOpt/dp/SimulatorMultiStageDPBase.h"
#include "StOpt/regression/GridAndRegressedValue.h"
#include "StOpt/regression/ContinuationValue.h"

/** \file OptimizerMultiStageDPBase.h
 *  \brief Define an abstract class for Dynamic Programming problems  where each transition problem problem
 *         is itself solved using  a dynamic programming approach
 *     \author Benoit Clair, Xavier Warin
 */

namespace StOpt
{

/// \class OptimizerMultiStageDPBase OptimizerMultiStageDPBase.h
///  Base class for optimizer for Dynamic Programming with regressions where each each transition problem
///  is itself solved using a deterministic DP
class OptimizerMultiStageDPBase : public OptimizerBase
{

public :

    OptimizerMultiStageDPBase() {}

    virtual ~OptimizerMultiStageDPBase() {}



    /// \brief defines a step in optimization
    /// \param p_grid           grid at arrival step after command
    /// \param p_stock          coordinates of the stock point to treat
    /// \param p_condEsp        continuation values for each regime
    /// \param p_phiIn          for each regime  gives the solution calculated at the previous step ( next time step by Dynamic Programming resolution) : structure of the 2D array ( nb simulation ,nb stocks )
    /// \return   for each regimes (column) gives the solution for each particle (row)
    [[nodiscard]] virtual  Eigen::ArrayXXd stepOptimize(const std::shared_ptr<StOpt::SpaceGrid> &p_grid,
            const Eigen::ArrayXd &p_stock,
            const std::vector< std::shared_ptr<ContinuationValue> > &p_condEsp,
            const std::vector<std::shared_ptr<Eigen::ArrayXXd>> &p_phiIn) const = 0;

    /// \brief defines a step in simulation
    /// Notice that this implementation is not optimal but is convenient if the control is discrete.
    /// By avoiding interpolation in control we avoid non admissible control
    /// Control are recalculated during simulation.
    /// \param p_grid          grid at arrival step after command
    /// \param p_continuation  defines the continuation operator for each regime
    /// \param p_state         defines the state value (modified)
    /// \param p_phiInOut      defines the value functions (modified) : size number of functions  to follow
    virtual void stepSimulate(const std::shared_ptr< StOpt::SpaceGrid>   &p_grid, const std::vector< StOpt::GridAndRegressedValue  > &p_continuation,
                              StOpt::StateWithStocks &p_state,
                              Eigen::Ref<Eigen::ArrayXd> p_phiInOut) const = 0 ;


    // number of regime used in deterministic part
    virtual   int getNbDetRegime() const = 0 ;


    /// \brief get the simulator back
    virtual std::shared_ptr< SimulatorMultiStageDPBase > getSimulator() const = 0;

};
}
#endif /* OPTIMIZERMULTISTAGEDPBASE_H */