Skip to content

Commit db9feec

Browse files
arun3688lochel
andauthored
Refactor setTime (#1386)
Co-authored-by: Lennart Ochel <lennart.ochel@outlook.com>
1 parent fed7b0e commit db9feec

5 files changed

Lines changed: 18 additions & 55 deletions

File tree

src/OMSimulatorLib/Component.h

Lines changed: 2 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -101,7 +101,8 @@ namespace oms
101101
virtual oms_status_enu_t setString(const ComRef& cref, const std::string& value) { return logError_NotImplemented; }
102102
virtual oms_status_enu_t setUnit(const ComRef& cref, const std::string& value) { return logError_NotImplemented; }
103103
virtual oms_status_enu_t setRealInputDerivative(const ComRef& cref, const SignalDerivative& der) { return logError_NotImplemented; }
104-
virtual oms_status_enu_t stepUntil(double stopTime) { return oms_status_ok; }
104+
virtual oms_status_enu_t stepUntil(double stopTime) { return logError_NotImplemented; }
105+
virtual oms_status_enu_t setTime(double time) { return logError_NotImplemented; }
105106
virtual oms_status_enu_t newResources(const std::string& ssvFileName, const std::string& ssmFileName, bool externalResources) { return logError_NotImplemented; }
106107
virtual oms_status_enu_t addResources(std::string& filename) { return logError_NotImplemented; }
107108
virtual oms_status_enu_t deleteReferencesInSSD(const std::string& filename) {return logError_NotImplemented;}

src/OMSimulatorLib/ComponentFMUME.cpp

Lines changed: 8 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -1018,6 +1018,14 @@ oms_status_enu_t oms::ComponentFMUME::getInteger(const fmi2ValueReference& vr, i
10181018
return oms_status_ok;
10191019
}
10201020

1021+
oms_status_enu_t oms::ComponentFMUME::setTime(double time)
1022+
{
1023+
fmi2Status fmistatus = fmi2_setTime(fmu, time);
1024+
if (fmi2OK != fmistatus)
1025+
return logError_FMUCall("fmi2_setTime", this);
1026+
return oms_status_ok;
1027+
}
1028+
10211029
oms_status_enu_t oms::ComponentFMUME::getInteger(const ComRef& cref, int& value)
10221030
{
10231031
CallClock callClock(clock);

src/OMSimulatorLib/ComponentFMUME.h

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -90,6 +90,7 @@ namespace oms
9090
oms_status_enu_t setReal(const ComRef& cref, double value);
9191
oms_status_enu_t setString(const ComRef& cref, const std::string& value);
9292
oms_status_enu_t setUnit(const ComRef& cref, const std::string& value);
93+
oms_status_enu_t setTime(double time);
9394

9495
oms_status_enu_t getDirectionalDerivative(const ComRef& unknownCref, const ComRef& knownCref, double& value);
9596
oms_status_enu_t getDirectionalDerivativeHeper(const int unknownIndex, const int knownIndex, const std::vector<int>& dependencyList, double& value);

src/OMSimulatorLib/ComponentTable.h

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -63,6 +63,7 @@ namespace oms
6363
oms_status_enu_t reset();
6464

6565
oms_status_enu_t stepUntil(double stopTime) {time = stopTime; return oms_status_ok;}
66+
oms_status_enu_t setTime(double time) {this->time = time; return oms_status_ok;}
6667

6768
Variable* getVariable(const ComRef& cref) {logError_NotImplemented; return NULL;}
6869
oms_status_enu_t getReal(const ComRef& cref, double& value);

src/OMSimulatorLib/SystemSC.cpp

Lines changed: 6 additions & 54 deletions
Original file line numberDiff line numberDiff line change
@@ -51,9 +51,7 @@ int oms::cvode_rhs(realtype t, N_Vector y, N_Vector ydot, void* user_data)
5151
// update states in FMUs
5252
for (size_t i=0, j=0; i < system->fmus.size(); ++i)
5353
{
54-
// set time
55-
fmistatus = fmi2_setTime(system->fmus[i]->getFMU(), t);
56-
if (fmi2OK != fmistatus) logError_FMUCall("fmi2_setTime", system->fmus[i]);
54+
system->fmus[i]->setTime(t);
5755

5856
if (0 == system->nStates[i])
5957
continue;
@@ -95,9 +93,7 @@ int oms::cvode_roots(realtype t, N_Vector y, realtype *gout, void *user_data)
9593
// update states in FMUs
9694
for (size_t i=0, j=0; i < system->fmus.size(); ++i)
9795
{
98-
// set time
99-
fmistatus = fmi2_setTime(system->fmus[i]->getFMU(), t);
100-
if (fmi2OK != fmistatus) logError_FMUCall("fmi2_setTime", system->fmus[i]);
96+
system->fmus[i]->setTime(t);
10197

10298
if (0 == system->nStates[i])
10399
continue;
@@ -663,18 +659,7 @@ oms_status_enu_t oms::SystemSC::doStepEuler()
663659

664660
// set time
665661
for (const auto& component : getComponents())
666-
{
667-
switch (component.second->getType())
668-
{
669-
case oms_component_fmu:
670-
if (fmi2OK != fmi2_setTime(dynamic_cast<ComponentFMUME*>(component.second)->getFMU(), time))
671-
logError_FMUCall("fmi2_setTime", dynamic_cast<ComponentFMUME*>(component.second));
672-
break;
673-
case oms_component_table:
674-
dynamic_cast<ComponentTable*>(component.second)->stepUntil(time);
675-
break;
676-
}
677-
}
662+
component.second->setTime(time);
678663

679664
for (size_t i = 0; i < fmus.size(); ++i)
680665
{
@@ -713,18 +698,7 @@ oms_status_enu_t oms::SystemSC::doStepEuler()
713698

714699
// set time
715700
for (const auto& component : getComponents())
716-
{
717-
switch (component.second->getType())
718-
{
719-
case oms_component_fmu:
720-
if (fmi2OK != fmi2_setTime(dynamic_cast<ComponentFMUME*>(component.second)->getFMU(), time))
721-
logError_FMUCall("fmi2_setTime", dynamic_cast<ComponentFMUME*>(component.second));
722-
break;
723-
case oms_component_table:
724-
dynamic_cast<ComponentTable*>(component.second)->stepUntil(time);
725-
break;
726-
}
727-
}
701+
component.second->setTime(time);
728702

729703
// Enter event mode and handle discrete state updates for each FMU
730704
for (size_t i = 0; i < fmus.size(); ++i)
@@ -841,18 +815,7 @@ oms_status_enu_t oms::SystemSC::doStepCVODE()
841815

842816
// set time
843817
for (const auto& component : getComponents())
844-
{
845-
switch (component.second->getType())
846-
{
847-
case oms_component_fmu:
848-
if (fmi2OK != fmi2_setTime(dynamic_cast<ComponentFMUME*>(component.second)->getFMU(), time))
849-
logError_FMUCall("fmi2_setTime", dynamic_cast<ComponentFMUME*>(component.second));
850-
break;
851-
case oms_component_table:
852-
dynamic_cast<ComponentTable*>(component.second)->stepUntil(time);
853-
break;
854-
}
855-
}
818+
component.second->setTime(time);
856819

857820
for (size_t i = 0; i < fmus.size(); ++i)
858821
{
@@ -931,18 +894,7 @@ oms_status_enu_t oms::SystemSC::doStepCVODE()
931894

932895
// set time
933896
for (const auto& component : getComponents())
934-
{
935-
switch (component.second->getType())
936-
{
937-
case oms_component_fmu:
938-
if (fmi2OK != fmi2_setTime(dynamic_cast<ComponentFMUME*>(component.second)->getFMU(), time))
939-
logError_FMUCall("fmi2_setTime", dynamic_cast<ComponentFMUME*>(component.second));
940-
break;
941-
case oms_component_table:
942-
dynamic_cast<ComponentTable*>(component.second)->stepUntil(time);
943-
break;
944-
}
945-
}
897+
component.second->setTime(time);
946898

947899
for (size_t i = 0; i < fmus.size(); ++i)
948900
{

0 commit comments

Comments
 (0)