OMSimulatorLib
The OMSimulator project is a FMI-based co-simulation environment.
SystemSC.h
Go to the documentation of this file.
1 /*
2  * This file is part of OpenModelica.
3  *
4  * Copyright (c) 1998-CurrentYear, Open Source Modelica Consortium (OSMC),
5  * c/o Linköpings universitet, Department of Computer and Information Science,
6  * SE-58183 Linköping, Sweden.
7  *
8  * All rights reserved.
9  *
10  * THIS PROGRAM IS PROVIDED UNDER THE TERMS OF GPL VERSION 3 LICENSE OR
11  * THIS OSMC PUBLIC LICENSE (OSMC-PL) VERSION 1.2.
12  * ANY USE, REPRODUCTION OR DISTRIBUTION OF THIS PROGRAM CONSTITUTES
13  * RECIPIENT'S ACCEPTANCE OF THE OSMC PUBLIC LICENSE OR THE GPL VERSION 3,
14  * ACCORDING TO RECIPIENTS CHOICE.
15  *
16  * The OpenModelica software and the Open Source Modelica
17  * Consortium (OSMC) Public License (OSMC-PL) are obtained
18  * from OSMC, either from the above address,
19  * from the URLs: http://www.ida.liu.se/projects/OpenModelica or
20  * http://www.openmodelica.org, and in the OpenModelica distribution.
21  * GNU version 3 is obtained from: http://www.gnu.org/copyleft/gpl.html.
22  *
23  * This program is distributed WITHOUT ANY WARRANTY; without
24  * even the implied warranty of MERCHANTABILITY or FITNESS
25  * FOR A PARTICULAR PURPOSE, EXCEPT AS EXPRESSLY SET FORTH
26  * IN THE BY RECIPIENT SELECTED SUBSIDIARY LICENSE CONDITIONS OF OSMC-PL.
27  *
28  * See the full OSMC Public License conditions for more details.
29  *
30  */
31 
32 #ifndef _OMS_SYSTEM_SC_H_
33 #define _OMS_SYSTEM_SC_H_
34 
35 #include "ComRef.h"
36 #include "System.h"
37 #include "OMSimulator/Types.h"
38 
39 #include <cvode/cvode.h> /* prototypes for CVODE fcts., consts. */
40 #include <nvector/nvector_serial.h> /* serial N_Vector types, fcts., macros */
41 #include <sunlinsol/sunlinsol_dense.h> /* Default dense linear solver */
42 
43 namespace oms
44 {
45  class Model;
46  class ComponentFMUME;
47  int cvode_rhs(realtype t, N_Vector y, N_Vector ydot, void* user_data);
48  int cvode_rhs_algebraic(realtype t, N_Vector y, N_Vector ydot, void* user_data);
49  int cvode_roots(realtype t, N_Vector y, realtype *gout, void* user_data);
50 
51  class SystemSC : public System
52  {
53  public:
54  ~SystemSC();
55 
57  oms_status_enu_t exportToSSD_SimulationInformation(pugi::xml_node& node) const;
58  oms_status_enu_t importFromSSD_SimulationInformation(const pugi::xml_node& node, const std::string& sspVersion);
59 
60  oms_status_enu_t instantiate();
61  oms_status_enu_t initialize();
62  oms_status_enu_t terminate();
63  oms_status_enu_t reset();
64  oms_status_enu_t doStep();
65  oms_status_enu_t stepUntil(double stopTime);
66 
67  oms_status_enu_t updateInputs(DirectedGraph& graph);
68 
69  std::string getSolverName() const;
70  oms_status_enu_t setSolverMethod(std::string);
71 
72  oms_status_enu_t setSolver(oms_solver_enu_t solver) {if (solver > oms_solver_sc_min && solver < oms_solver_sc_max) {solverMethod=solver; return oms_status_ok;} return oms_status_error;}
73 
74  private:
75  oms_status_enu_t doStepEuler();
76  oms_status_enu_t doStepCVODE();
77 
78  protected:
80 
81  // stop the compiler generating methods copying the object
82  SystemSC(SystemSC const& copy);
83  SystemSC& operator=(SystemSC const& copy);
84 
85  private:
86  std::vector<ComponentFMUME*> fmus;
87 
88  std::vector<fmi2Boolean> callEventUpdate;
89  std::vector<fmi2Boolean> terminateSimulation;
90  std::vector<size_t> nStates;
91  std::vector<size_t> nEventIndicators;
92 
93  std::vector<double*> states;
94  std::vector<double*> states_der;
95  std::vector<double*> states_nominal;
96  std::vector<double*> event_indicators;
97  std::vector<double*> event_indicators_prev;
98 
99  bool algebraic = false;
100 
102  {
103  };
104 
106  {
107  void *mem;
108  N_Vector y;
109  SUNLinearSolver linSol; /* linear solver object */
110  SUNMatrix J; /* Matrix used by linear solver */
111  N_Vector liny; /* Vector used by linear solver */
112  N_Vector abstol;
113  };
114 
116  {
120 
121  friend int oms::cvode_rhs(realtype t, N_Vector y, N_Vector ydot, void* user_data);
122  friend int oms::cvode_rhs_algebraic(realtype t, N_Vector y, N_Vector ydot, void* user_data);
123  friend int oms::cvode_roots(realtype t, N_Vector y, realtype *gout, void* user_data);
124  };
125 }
126 
127 #endif
ComRef - component reference.
Definition: ComRef.h:47
Definition: DirectedGraph.h:65
Definition: Model.h:53
Definition: SystemSC.h:52
oms_status_enu_t importFromSSD_SimulationInformation(const pugi::xml_node &node, const std::string &sspVersion)
Definition: SystemSC.cpp:193
oms_status_enu_t setSolverMethod(std::string)
Definition: SystemSC.cpp:160
oms_status_enu_t doStepEuler()
Definition: SystemSC.cpp:549
std::vector< double * > event_indicators
Definition: SystemSC.h:96
bool algebraic
Definition: SystemSC.h:99
std::vector< ComponentFMUME * > fmus
Definition: SystemSC.h:86
std::string getSolverName() const
Definition: SystemSC.cpp:147
oms_status_enu_t terminate()
Definition: SystemSC.cpp:422
static System * NewSystem(const oms::ComRef &cref, Model *parentModel, System *parentSystem)
Definition: SystemSC.cpp:129
oms_status_enu_t reset()
Definition: SystemSC.cpp:487
oms_status_enu_t setSolver(oms_solver_enu_t solver)
Definition: SystemSC.h:72
oms_status_enu_t exportToSSD_SimulationInformation(pugi::xml_node &node) const
Definition: SystemSC.cpp:172
SystemSC & operator=(SystemSC const &copy)
not implemented
union oms::SystemSC::SolverData_t solverData
~SystemSC()
Definition: SystemSC.cpp:125
std::vector< double * > event_indicators_prev
Definition: SystemSC.h:97
std::vector< size_t > nStates
Definition: SystemSC.h:90
std::vector< double * > states
Definition: SystemSC.h:93
SystemSC(SystemSC const &copy)
not implemented
oms_status_enu_t doStepCVODE()
Definition: SystemSC.cpp:771
std::vector< fmi2Boolean > callEventUpdate
Definition: SystemSC.h:88
oms_status_enu_t initialize()
Definition: SystemSC.cpp:266
oms_status_enu_t instantiate()
Definition: SystemSC.cpp:220
SystemSC(const ComRef &cref, Model *parentModel, System *parentSystem)
Definition: SystemSC.cpp:120
std::vector< double * > states_nominal
Definition: SystemSC.h:95
oms_status_enu_t updateInputs(DirectedGraph &graph)
Definition: SystemSC.cpp:956
oms_status_enu_t stepUntil(double stopTime)
Definition: SystemSC.cpp:930
std::vector< double * > states_der
Definition: SystemSC.h:94
oms_status_enu_t doStep()
Definition: SystemSC.cpp:534
std::vector< fmi2Boolean > terminateSimulation
Definition: SystemSC.h:89
std::vector< size_t > nEventIndicators
Definition: SystemSC.h:91
Definition: System.h:62
ComRef cref
Definition: System.h:225
oms_solver_enu_t solverMethod
Definition: System.h:214
System * parentSystem
Definition: System.h:228
Model * parentModel
Definition: System.h:227
Definition: AlgLoop.h:45
int cvode_rhs_algebraic(realtype t, N_Vector y, N_Vector ydot, void *user_data)
Definition: SystemSC.cpp:77
int cvode_roots(realtype t, N_Vector y, realtype *gout, void *user_data)
Definition: SystemSC.cpp:89
int cvode_rhs(realtype t, N_Vector y, N_Vector ydot, void *user_data)
Definition: SystemSC.cpp:45
Definition: SystemSC.h:106
SUNMatrix J
Definition: SystemSC.h:110
SUNLinearSolver linSol
Definition: SystemSC.h:109
N_Vector y
Definition: SystemSC.h:108
N_Vector abstol
Definition: SystemSC.h:112
void * mem
Definition: SystemSC.h:107
N_Vector liny
Definition: SystemSC.h:111
Definition: SystemSC.h:102
Definition: SystemSC.h:116
SolverDataEuler_t euler
Definition: SystemSC.h:117
SolverDataCVODE_t cvode
Definition: SystemSC.h:118