Skip to content
This repository was archived by the owner on May 18, 2019. It is now read-only.

Commit f5cbc6b

Browse files
lochelOpenModelica-Hudson
authored andcommitted
Fix CS-FMU doStep function
- Fix event handling - Fix integration method (euler) Belonging to [master]: - #2795
1 parent dc05202 commit f5cbc6b

File tree

1 file changed

+80
-104
lines changed

1 file changed

+80
-104
lines changed

SimulationRuntime/fmi/export/fmi2/fmu2_model_interface.c

Lines changed: 80 additions & 104 deletions
Original file line numberDiff line numberDiff line change
@@ -70,14 +70,14 @@ fmi2ValueReference vrStatesDerivatives[NUMBER_OF_STATES] = STATESDERIVATIVES;
7070
// ---------------------------------------------------------------------------
7171
const char* stateToString(ModelInstance *comp)
7272
{
73-
switch (comp->state) {
73+
switch (comp->state)
74+
{
7475
case modelInstantiated: return "Instantiated";
7576
case modelInitializationMode: return "Initialization Mode";
7677
case modelEventMode: return "Event Mode";
7778
case modelContinuousTimeMode: return "Continuous-Time Mode";
7879
case modelTerminated: return "Terminated";
7980
case modelError: return "Error";
80-
default: break;
8181
}
8282
return "Unknown";
8383
}
@@ -230,7 +230,7 @@ fmi2Status fmi2EventUpdate(fmi2Component c, fmi2EventInfo* eventInfo)
230230

231231
/* activate sample event */
232232
for(i=0; i<comp->fmuData->modelData->nSamples; ++i) {
233-
if(comp->fmuData->simulationInfo->nextSampleTimes[i] <= comp->fmuData->localData[0]->timeValue) {
233+
if (comp->fmuData->simulationInfo->nextSampleTimes[i] <= comp->fmuData->localData[0]->timeValue) {
234234
comp->fmuData->simulationInfo->samples[i] = 1;
235235
infoStreamPrint(LOG_EVENTS, 0, "[%ld] sample(%g, %g)", comp->fmuData->modelData->samplesInfo[i].index, comp->fmuData->modelData->samplesInfo[i].start, comp->fmuData->modelData->samplesInfo[i].interval);
236236
}
@@ -240,14 +240,14 @@ fmi2Status fmi2EventUpdate(fmi2Component c, fmi2EventInfo* eventInfo)
240240

241241
/* deactivate sample events */
242242
for(i=0; i<comp->fmuData->modelData->nSamples; ++i) {
243-
if(comp->fmuData->simulationInfo->samples[i]) {
243+
if (comp->fmuData->simulationInfo->samples[i]) {
244244
comp->fmuData->simulationInfo->samples[i] = 0;
245245
comp->fmuData->simulationInfo->nextSampleTimes[i] += comp->fmuData->modelData->samplesInfo[i].interval;
246246
}
247247
}
248248

249249
for(i=0; i<comp->fmuData->modelData->nSamples; ++i) {
250-
if((i == 0) || (comp->fmuData->simulationInfo->nextSampleTimes[i] < comp->fmuData->simulationInfo->nextSampleEvent)) {
250+
if ((i == 0) || (comp->fmuData->simulationInfo->nextSampleTimes[i] < comp->fmuData->simulationInfo->nextSampleEvent)) {
251251
comp->fmuData->simulationInfo->nextSampleEvent = comp->fmuData->simulationInfo->nextSampleTimes[i];
252252
}
253253
}
@@ -1323,6 +1323,8 @@ fmi2Status fmi2DoStep(fmi2Component c, fmi2Real currentCommunicationPoint, fmi2R
13231323
fmi2Real t = comp->fmuData->localData[0]->timeValue;
13241324
fmi2Real tNext, tEnd;
13251325
fmi2Boolean enterEventMode = fmi2False, terminateSimulation = fmi2False;
1326+
fmi2Real tCommunication;
1327+
13261328
fmi2EventInfo eventInfo;
13271329
eventInfo.newDiscreteStatesNeeded = fmi2False;
13281330
eventInfo.terminateSimulation = fmi2False;
@@ -1331,138 +1333,112 @@ fmi2Status fmi2DoStep(fmi2Component c, fmi2Real currentCommunicationPoint, fmi2R
13311333
eventInfo.nextEventTimeDefined = fmi2False;
13321334
eventInfo.nextEventTime = -0.0;
13331335

1334-
/* fprintf(stderr, "DoStep %g + %g State: %d\n", currentCommunicationPoint, communicationStepSize, comp->state); */
1336+
if (comp->stopTimeDefined)
1337+
tEnd = comp->stopTime;
1338+
else
1339+
tEnd = currentCommunicationPoint + communicationStepSize;
1340+
tCommunication = currentCommunicationPoint;
13351341

13361342
fmi2EnterEventMode(c);
13371343
fmi2EventIteration(c, &eventInfo);
13381344
fmi2EnterContinuousTimeMode(c);
13391345

1340-
if (NUMBER_OF_STATES > 0)
1346+
while (status == fmi2OK && comp->fmuData->localData[0]->timeValue < tEnd)
13411347
{
1342-
status = fmi2GetDerivatives(c, states_der, NUMBER_OF_STATES);
1343-
if (status != fmi2OK)
1348+
while (tCommunication <= comp->fmuData->localData[0]->timeValue)
13441349
{
1345-
functions->freeMemory(states);
1346-
functions->freeMemory(states_der);
1347-
functions->freeMemory(event_indicators);
1348-
functions->freeMemory(event_indicators_prev);
1349-
return fmi2Error;
1350+
tCommunication += communicationStepSize;
13501351
}
13511352

1352-
status = fmi2GetContinuousStates(c, states, NUMBER_OF_STATES);
1353-
if (status != fmi2OK)
1354-
{
1355-
functions->freeMemory(states);
1356-
functions->freeMemory(states_der);
1357-
functions->freeMemory(event_indicators);
1358-
functions->freeMemory(event_indicators_prev);
1359-
return fmi2Error;
1360-
}
1361-
}
1353+
/* fprintf(stderr, "DoStep %g -> %g State: %s\n", comp->fmuData->localData[0]->timeValue, tNext, stateToString(comp)); */
13621354

1363-
if (NUMBER_OF_EVENT_INDICATORS > 0)
1364-
{
1365-
status = fmi2GetEventIndicators(c, event_indicators_prev, NUMBER_OF_EVENT_INDICATORS);
1366-
if (status != fmi2OK)
1355+
if (NUMBER_OF_STATES > 0)
13671356
{
1368-
functions->freeMemory(states);
1369-
functions->freeMemory(states_der);
1370-
functions->freeMemory(event_indicators);
1371-
functions->freeMemory(event_indicators_prev);
1372-
return fmi2Error;
1373-
}
1374-
}
1375-
1376-
tNext = currentCommunicationPoint + communicationStepSize;
1357+
status = fmi2GetDerivatives(c, states_der, NUMBER_OF_STATES);
1358+
if (status != fmi2OK) {status=fmi2Error; break;}
13771359

1378-
/* adjust tNext step to get tEnd exactly */
1379-
if (comp->stopTimeDefined)
1380-
{
1381-
tEnd = comp->stopTime;
1382-
}
1383-
else
1384-
{
1385-
tEnd = currentCommunicationPoint + 2*communicationStepSize + 1;
1386-
}
1387-
if(tNext > tEnd - communicationStepSize/1e16) {
1388-
tNext = tEnd;
1389-
}
13901360

1391-
/* adjust for time events */
1392-
if (eventInfo.nextEventTimeDefined && (tNext >= eventInfo.nextEventTime)) {
1393-
tNext = eventInfo.nextEventTime;
1394-
time_event = 1;
1395-
}
1361+
status = fmi2GetContinuousStates(c, states, NUMBER_OF_STATES);
1362+
if (status != fmi2OK) {status=fmi2Error; break;}
1363+
}
13961364

1397-
fmi2SetTime(c, tNext);
1365+
if (NUMBER_OF_EVENT_INDICATORS > 0)
1366+
{
1367+
status = fmi2GetEventIndicators(c, event_indicators_prev, NUMBER_OF_EVENT_INDICATORS);
1368+
if (status != fmi2OK) {status=fmi2Error; break;}
1369+
}
13981370

1399-
/* integrate */
1400-
for (i = 0; i < NUMBER_OF_STATES; i++) {
1401-
states[i] = states[i] + communicationStepSize * states_der[i];
1402-
}
1371+
/* adjust tNext step to get tEnd exactly */
1372+
if (tCommunication > tEnd - communicationStepSize/1e16)
1373+
tNext = tEnd;
1374+
else
1375+
tNext = tCommunication;
14031376

1404-
/* set the continuous states */
1405-
if (NUMBER_OF_STATES > 0)
1406-
{
1407-
status = fmi2SetContinuousStates(c, states, NUMBER_OF_STATES);
1408-
if (status != fmi2OK)
1377+
/* adjust for time events */
1378+
if (eventInfo.nextEventTimeDefined && (eventInfo.nextEventTime <= tNext))
14091379
{
1410-
functions->freeMemory(states);
1411-
functions->freeMemory(states_der);
1412-
functions->freeMemory(event_indicators);
1413-
functions->freeMemory(event_indicators_prev);
1414-
return fmi2Error;
1380+
tNext = eventInfo.nextEventTime;
1381+
time_event = 1;
14151382
}
1416-
}
14171383

1418-
/* signal completed integrator step */
1419-
status = fmi2CompletedIntegratorStep(c, fmi2True, &enterEventMode, &terminateSimulation);
1420-
if (status != fmi2OK)
1421-
{
1422-
functions->freeMemory(states);
1423-
functions->freeMemory(states_der);
1424-
functions->freeMemory(event_indicators);
1425-
functions->freeMemory(event_indicators_prev);
1426-
return fmi2Error;
1427-
}
1384+
/* integrate */
1385+
for (i = 0; i < NUMBER_OF_STATES; i++)
1386+
{
1387+
states[i] = states[i] + (tNext - comp->fmuData->localData[0]->timeValue) * states_der[i];
1388+
}
1389+
fmi2SetTime(c, tNext);
14281390

1429-
/* check for events */
1430-
if (NUMBER_OF_EVENT_INDICATORS > 0)
1431-
{
1432-
status = fmi2GetEventIndicators(c, event_indicators, NUMBER_OF_EVENT_INDICATORS);
1433-
if (status != fmi2OK)
1391+
/* set the continuous states */
1392+
if (NUMBER_OF_STATES > 0)
14341393
{
1435-
functions->freeMemory(states);
1436-
functions->freeMemory(states_der);
1437-
functions->freeMemory(event_indicators);
1438-
functions->freeMemory(event_indicators_prev);
1439-
return fmi2Error;
1394+
status = fmi2SetContinuousStates(c, states, NUMBER_OF_STATES);
1395+
if (status != fmi2OK) {status=fmi2Error; break;}
14401396
}
14411397

1442-
for (i = 0; i < NUMBER_OF_EVENT_INDICATORS; i++)
1398+
/* signal completed integrator step */
1399+
status = fmi2CompletedIntegratorStep(c, fmi2True, &enterEventMode, &terminateSimulation);
1400+
if (status != fmi2OK) {status=fmi2Error; break;}
1401+
1402+
/* check for events */
1403+
if (NUMBER_OF_EVENT_INDICATORS > 0)
14431404
{
1444-
if (event_indicators[i]*event_indicators_prev[i] < 0) {
1445-
zc_event = 1;
1446-
break;
1405+
status = fmi2GetEventIndicators(c, event_indicators, NUMBER_OF_EVENT_INDICATORS);
1406+
if (status != fmi2OK) {status=fmi2Error; break;}
1407+
1408+
for (i = 0; i < NUMBER_OF_EVENT_INDICATORS; i++)
1409+
{
1410+
if (event_indicators[i]*event_indicators_prev[i] < 0)
1411+
{
1412+
zc_event = 1;
1413+
break;
1414+
}
14471415
}
14481416
}
14491417

1450-
/* fprintf(stderr, "enterEventMode = %d, zc_event = %d, time_event = %d\n", enterEventMode, zc_event, time_event); */
1418+
if (enterEventMode || zc_event || time_event)
1419+
{
1420+
/* fprintf(stderr, "enterEventMode = %d, zc_event = %d, time_event = %d\n", enterEventMode, zc_event, time_event); */
14511421

1452-
if (enterEventMode || zc_event || time_event) {
14531422
fmi2EnterEventMode(c);
1454-
14551423
fmi2EventIteration(c, &eventInfo);
14561424

1457-
if(eventInfo.valuesOfContinuousStatesChanged)
1458-
fmi2GetContinuousStates(c, states, NUMBER_OF_STATES);
1425+
if (eventInfo.valuesOfContinuousStatesChanged)
1426+
{
1427+
status = fmi2GetContinuousStates(c, states, NUMBER_OF_STATES);
1428+
if (status != fmi2OK) {status=fmi2Error; break;}
1429+
}
14591430

1460-
if( eventInfo.nominalsOfContinuousStatesChanged)
1461-
fmi2GetNominalsOfContinuousStates(c, states, NUMBER_OF_STATES);
1431+
if (eventInfo.nominalsOfContinuousStatesChanged)
1432+
{
1433+
status = fmi2GetNominalsOfContinuousStates(c, states, NUMBER_OF_STATES);
1434+
if (status != fmi2OK) {status=fmi2Error; break;}
1435+
}
14621436

1463-
fmi2GetEventIndicators(c, event_indicators_prev, NUMBER_OF_EVENT_INDICATORS);
1437+
status = fmi2GetEventIndicators(c, event_indicators_prev, NUMBER_OF_EVENT_INDICATORS);
1438+
if (status != fmi2OK) {status=fmi2Error; break;}
14641439

1465-
fmi2EnterContinuousTimeMode(c);
1440+
status = fmi2EnterContinuousTimeMode(c);
1441+
if (status != fmi2OK) {status=fmi2Error; break;}
14661442
}
14671443
}
14681444

@@ -1471,7 +1447,7 @@ fmi2Status fmi2DoStep(fmi2Component c, fmi2Real currentCommunicationPoint, fmi2R
14711447
functions->freeMemory(event_indicators);
14721448
functions->freeMemory(event_indicators_prev);
14731449

1474-
return fmi2OK;
1450+
return status;
14751451
}
14761452

14771453
fmi2Status fmi2CancelStep(fmi2Component c)

0 commit comments

Comments
 (0)