Class fmipp.FMUModelExchangeV1

class FMUModelExchangeV1([fmuDirUri, ]modelIdentifier[, loggingOn=False, stopBeforeEvent=False, eventSearchPrecision=1e-4, type=integratorBDF])

This class provides a basic Python wrapper that offers a set of convenient methods for accessing and manipulating FMUs for Model Exchange according to the FMI ME V1.0 standard. Instances of class FMUModelExchangeV1 own the actual FMU instance.

When creating the first instance of this class, the FMU has to be loaded. For this, the URI to the unzipped FMU has to be provided:

fmu = fmipp.FMUModelExchangeV1(fmuDirUri, modelIdentifier, False, 1e-4, integratorBDF)

An FMU can be instantiated many times (provided capability flag canBeInstantiatedOnlyOncePerProcess is False). In case you want to have more than one instance of the same FMU, just create another instance of this class without providing the URI to the unzipped FMU (the FMU will have already been loaded in the background the first time):

fmu = fmipp.FMUModelExchangeV1(modelIdentifier, False, 1e-4, integratorBDF)

The available numerical integrators are listed in the following table:

Stepper

Name

Order

Adaptive

Recommended usecases

integratorEU

Explicit Euler

1

No

Testing

integratorRK

4th order Runge-Kutta

4

No

Testing

integratorABM

Adams-Bashforth-Moulton

8

No

Testing

integratorCK

Cash-Karp

5

Yes

Nonstiff Models

integratorDP

Dormand-Prince

5

Yes

Nonstiff Models

integratorFE

Fehlberg

8

Yes

Nonstiff, smooth Models

integratorBS

Bulirsch Stoer

1-16

Yes

High precision required

integratorRO

Rosenbrock

4

Yes

Stiff Models

integratorBDF

Backward Differentiation Formula

1-5

Yes

Stiff Models

integratorABM2

Adams-Bashforth-Moulton

1-12

Yes

Nonstiff Models with expensive right hand side

Parameters
  • fmuDirUri (str) – URI to the unzipped FMU

  • modelIdentifier (str) – model identifier of the FMU

  • loggingOn (bool) – enable logger for debugging

  • stopBeforeEvent (bool) – if True, integration stops immediately before an event

  • eventSearchPrecision (bool) – numerical search precision for events during integration

  • type (int) – the numerical integration method for solving ODEs

Methods for interacting with the FMU

instantiate(instanceName)

This functions creates a new internal FMU instance. This function must be called successfully, before any of the following functions can be called.

Parameters

instanceName (str) – a unique identifier for a given FMI component instance, used to identify a component within a co-simulation graph model and for logging messages

Return type

int (statusOK, statusWarning, statusDiscard, statusError or statusFatal)

initialize([toleranceDefined=False, tolerance=1e-05])

Informs the FMU instance that the simulation run starts now.

Parameters
  • toleranceDefined (bool) – if True, the model is called with a numerical integration scheme where the step size is controlled by using tolerance

  • tolerance (float) – relative tolerance for error estimation

Return type

int (statusOK, statusWarning, statusDiscard, statusError or statusFatal)

integrate(tend[, deltaT=1e-05])

Integrates the FMU model until time tend or until the first event. The exact behaviour depends on the flag stopBeforeEvent.

Parameters
  • tend (float) – stop time for the integration

  • deltaT (int) – starting step size to be used by the integrator, smaller values lead to more accuracy

Returns

integration stop time

Return type

float

integrateN(tend, nsteps)

Integrates the FMU model until time tend or until the first event. The exact behaviour depends on the flag stopBeforeEvent.

Parameters
  • tend (float) – stop time for the integration

  • nsteps (int) – number of integrator steps to be recommended to the integrator, bigger values lead to more accuracy

Returns

integration stop time

Return type

float

Methods for getting/setting values

See also

Use the helper functions from module fmipp.numeric for retrieving derivatives from the FMU model.

getBooleanValue(name)
Parameters

name (str) – variable name

Return type

bool

getIntegerValue(name)
Parameters

name (str) – variable name

Return type

int

getRealValue(name)
Parameters

name (str) – variable name

Return type

float

getStringValue(name)
Parameters

name (str) – variable name

Return type

str

getLastStatus()
Return type

int (statusOK, statusWarning, statusDiscard, statusError or statusFatal)

getTime()
Return type

float

getType(name)

Get information about the type of a variable.

Parameters

name (str) – variable name

Return type

int (typeReal, typeInteger, typeBoolean, typeString or typeUnknown)

getValueRef(name)

Get the value reference of a variable.

Parameters

name (str) – variable name

Return type

int

setBooleanValue(name, val)
Parameters
  • name (str) – variable name

  • val (bool) – vew value for the variable

Return type

int (statusOK, statusWarning, statusDiscard, statusError or statusFatal)

setIntegerValue(name, val)
Parameters
  • name (str) – variable name

  • val (int) – vew value for the variable

Return type

int (statusOK, statusWarning, statusDiscard, statusError or statusFatal)

setRealValue(name, val)
Parameters
  • name (str) – variable name

  • val (float) – vew value for the variable

Return type

int (statusOK, statusWarning, statusDiscard, statusError or statusFatal)

setStringValue(name, val)
Parameters
  • name (str) – variable name

  • val (str) – vew value for the variable

Return type

int (statusOK, statusWarning, statusDiscard, statusError or statusFatal)

Methods for event handling

checkEvents()

Check if any kind of event has happened.

Return type

bool

checkStateEvent()

Check if a state event has happened.

Return type

bool

checkStepEvent()

Check if a step event has happened.

Return type

bool

checkTimeEvent()

Check if a time event has happened.

Return type

bool

getEventSearchPrecision()

Get the event search precision.

Return type

float

getTimeEvent()

Get the time of the next time event (infinity if no time event is returned by the FMU model).

Return type

float

handleEvents()

Handle events. Just call this function if actually an event has happened.

raiseEvent()

Raise an event, i.e., notify that an (external) event has occured.

resetEventFlags()

Reset all internal flags related to event handling.

stepOverEvent()

If stopBeforeEvent is True, use this function to get the right-sided limit of an event.

Methods for retrieving model description flags

canBeInstantiatedOnlyOncePerProcess()
Return type

bool

canHandleEvents()
Return type

bool

canHandleVariableCommunicationStepSize()
Return type

bool

canInterpolateInputs()
Return type

bool

canNotUseMemoryManagementFunctions()
Return type

bool

canRejectSteps()
Return type

bool

canRunAsynchronuously()
Return type

bool

canSignalEvents()
Return type

bool

maxOutputDerivativeOrder()
Return type

int

nEventInds()
Return type

int

nStates()
Return type

int

nValueRefs()
Return type

int

Miscellaneous methods

setCallbacks(logger, allocateMemory, freeMemory)

Set FMU callback functions.

Return type

int (statusOK, statusWarning, statusDiscard, statusError or statusFatal)

logger(status, category, msg)

Call the FMU’s logger.

Parameters
  • status (int) – logger status

  • category (str) – logger category

  • msg (str) – logger message

sendDebugMessage(msg)

Send a debug message.

Parameters

msg (str) – debug message