class PIMixin(SimulationProblem): """ Adds `Delft-FEWS Published Interface <https://publicwiki.deltares.nl/display/FEWSDOC/The+Delft-Fews+Published+Interface>`_ I/O to your simulation problem. During preprocessing, files named ``rtcDataConfig.xml``, ``timeseries_import.xml``, and``rtcParameterConfig.xml`` are read from the ``input`` subfolder. ``rtcDataConfig.xml`` maps tuples of FEWS identifiers, including location and parameter ID, to RTC-Tools time series identifiers. During postprocessing, a file named ``timeseries_export.xml`` is written to the ``output`` subfolder. :cvar pi_binary_timeseries: Whether to use PI binary timeseries format. Default is ``False``. :cvar pi_parameter_config_basenames: List of parameter config file basenames to read. Default is [``rtcParameterConfig``]. :cvar pi_check_for_duplicate_parameters: Check if duplicate parameters are read. Default is ``False``. :cvar pi_validate_timeseries: Check consistency of timeseries. Default is ``True``. """ #: Whether to use PI binary timeseries format pi_binary_timeseries = False #: Location of rtcParameterConfig files pi_parameter_config_basenames = ['rtcParameterConfig'] #: Check consistency of timeseries pi_validate_timeseries = True #: Check for duplicate parameters pi_check_for_duplicate_parameters = True # Default names for timeseries I/O timeseries_import_basename = 'timeseries_import' timeseries_export_basename = 'timeseries_export' def __init__(self, **kwargs): # Check arguments assert ('input_folder' in kwargs) assert ('output_folder' in kwargs) # Save arguments self.__input_folder = kwargs['input_folder'] self.__output_folder = kwargs['output_folder'] # Load rtcDataConfig.xml. We assume this file does not change over the # life time of this object. self.__data_config = rtc.DataConfig(self.__input_folder) # Call parent class first for default behaviour. super().__init__(**kwargs) def pre(self): # Call parent class first for default behaviour. super().pre() # rtcParameterConfig self.__parameter_config = [] try: for pi_parameter_config_basename in self.pi_parameter_config_basenames: self.__parameter_config.append( pi.ParameterConfig(self.__input_folder, pi_parameter_config_basename)) except FileNotFoundError: raise FileNotFoundError("PIMixin: {}.xml not found in {}.".format( pi_parameter_config_basename, self.__input_folder)) # Make a parameters dict for later access self.__parameters = {} for parameter_config in self.__parameter_config: for location_id, model_id, parameter_id, value in parameter_config: try: parameter = self.__data_config.parameter( parameter_id, location_id, model_id) except KeyError: parameter = parameter_id self.__parameters[parameter] = value try: self.__timeseries_import = pi.Timeseries( self.__data_config, self.__input_folder, self.timeseries_import_basename, binary=self.pi_binary_timeseries, pi_validate_times=self.pi_validate_timeseries) except FileNotFoundError: raise FileNotFoundError('PIMixin: {}.xml not found in {}'.format( self.timeseries_import_basename, self.__input_folder)) self.__timeseries_export = pi.Timeseries( self.__data_config, self.__output_folder, self.timeseries_export_basename, binary=self.pi_binary_timeseries, pi_validate_times=False, make_new_file=True) # Convert timeseries timestamps to seconds since t0 for internal use self.__timeseries_import_times = self.__datetime_to_sec( self.__timeseries_import.times) # Timestamp check if self.pi_validate_timeseries: for i in range(len(self.__timeseries_import_times) - 1): if self.__timeseries_import_times[ i] >= self.__timeseries_import_times[i + 1]: raise ValueError( 'PIMixin: Time stamps must be strictly increasing.') # Check if the timeseries are equidistant self.__dt = self.__timeseries_import_times[ 1] - self.__timeseries_import_times[0] if self.pi_validate_timeseries: for i in range(len(self.__timeseries_import_times) - 1): if self.__timeseries_import_times[ i + 1] - self.__timeseries_import_times[i] != self.__dt: raise ValueError( 'PIMixin: Expecting equidistant timeseries, the time step ' 'towards {} is not the same as the time step(s) before. Set ' 'unit to nonequidistant if this is intended.'.format( self.__timeseries_import.times[i + 1])) # Stick timeseries into an AliasDict self.__timeseries_import_dict = AliasDict(self.alias_relation) debug = logger.getEffectiveLevel() == logging.DEBUG for variable, values in self.__timeseries_import.items(): if debug and self.__timeseries_import_dict.get(variable, None) is not None: logger.debug( 'PIMixin: Timeseries {} replaced another aliased timeseries.' .format(variable)) self.__timeseries_import_dict[variable] = values def initialize(self, config_file=None): # Set up experiment self.setup_experiment(0, self.__timeseries_import_times[-1], self.__dt) # Load parameters from parameter config self.__parameter_variables = set(self.get_parameter_variables()) logger.debug("Model parameters are {}".format( self.__parameter_variables)) for parameter, value in self.__parameters.items(): if parameter in self.__parameter_variables: logger.debug("PIMixin: Setting parameter {} = {}".format( parameter, value)) self.set_var(parameter, value) # Load input variable names self.__input_variables = set(self.get_input_variables().keys()) # Set input values for variable in self.__input_variables: value = self.__timeseries_import_dict[variable][ self.__timeseries_import.forecast_index] if np.isfinite(value): self.set_var(variable, value) else: logger.debug( 'PIMixin: Found bad value {} at index [{}] in timeseries aliased to input {}' .format(value, self.__timeseries_import.forecast_index, variable)) logger.debug("Model inputs are {}".format(self.__input_variables)) # Empty output self.__output_variables = self.get_output_variables() n_times = len(self.__timeseries_import_times) self.__output = AliasDict(self.alias_relation) self.__output.update({ variable: np.full(n_times, np.nan) for variable in self.__output_variables }) # Call super, which will also initialize the model itself super().initialize(config_file) # Extract consistent t0 values for variable in self.__output_variables: self.__output[variable][self.__timeseries_import. forecast_index] = self.get_var(variable) def update(self, dt): # Time step if dt < 0: dt = self.__dt # Current time stamp t = self.get_current_time() # Get current time index t_idx = bisect.bisect_left(self.__timeseries_import_times, t + dt) # Set input values for variable in self.__input_variables: value = self.__timeseries_import_dict[variable][t_idx] if np.isfinite(value): self.set_var(variable, value) else: logger.debug( 'PIMixin: Found bad value {} at index [{}] in timeseries aliased to input {}' .format(value, t_idx, variable)) # Call super super().update(dt) # Extract results for variable in self.__output_variables: self.__output[variable][t_idx] = self.get_var(variable) def post(self): # Call parent class first for default behaviour. super().post() # Start of write output # Write the time range for the export file. self.__timeseries_export.times = self.__timeseries_import.times[ self.__timeseries_import.forecast_index:] # Write other time settings self.__timeseries_export.forecast_datetime = self.__timeseries_import.forecast_datetime self.__timeseries_export.dt = self.__timeseries_import.dt self.__timeseries_export.timezone = self.__timeseries_import.timezone # Write the ensemble properties for the export file. self.__timeseries_export.ensemble_size = 1 self.__timeseries_export.contains_ensemble = self.__timeseries_import.contains_ensemble # For all variables that are output variables the values are # extracted from the results. for variable in self.__output_variables: values = self.__output[variable] # Check if ID mapping is present try: self.__data_config.pi_variable_ids(variable) except KeyError: logger.debug( 'PIMixin: variable {} has no mapping defined in rtcDataConfig ' 'so cannot be added to the output file.'.format(variable)) continue # Add series to output file self.__timeseries_export.set( variable, values, unit=self.__timeseries_import.get_unit(variable)) # Write output file to disk self.__timeseries_export.write() def __datetime_to_sec(self, d): # Return the date/timestamps in seconds since t0. if hasattr(d, '__iter__'): return np.array([ (t - self.__timeseries_import.forecast_datetime).total_seconds() for t in d ]) else: return ( d - self.__timeseries_import.forecast_datetime).total_seconds() def __sec_to_datetime(self, s): # Return the date/timestamps in seconds since t0 as datetime objects. if hasattr(s, '__iter__'): return [ self.__timeseries_import.forecast_datetime + timedelta(seconds=t) for t in s ] else: return self.__timeseries_import.forecast_datetime + timedelta( seconds=s) @cached def parameters(self): """ Return a dictionary of parameters, including parameters in PI Parameter Config XML files. :returns: Dictionary of parameters """ # Call parent class first for default values. parameters = super().parameters() # Load parameters from parameter config parameters.update(self.__parameters) if logger.getEffectiveLevel() == logging.DEBUG: for parameter in self.__parameters: logger.debug("CSVMixin: Read parameter {} ".format(parameter)) return parameters @cached def times(self, variable=None): """ Return a list of all the timesteps in seconds. :param variable: Variable name. :returns: A list of all the timesteps in seconds. """ return self.__timeseries_import_times[self.__timeseries_import. forecast_index:] def timeseries_at(self, variable, t): """ Return the value of a time series at the given time. :param variable: Variable name. :param t: Time. :returns: The interpolated value of the time series. :raises: KeyError """ values = self.__timeseries_import_dict[variable] t_idx = bisect.bisect_left(self.__timeseries_import_times, t) if self.__timeseries_import_times[t_idx] == t: return values[t_idx] else: return np.interp1d(t, self.__timeseries_import_times, values) @property def timeseries_import(self): """ :class:`pi.Timeseries` object containing the input data. """ return self.__timeseries_import @property def timeseries_import_times(self): """ List of time stamps for which input data is specified. The time stamps are in seconds since t0, and may be negative. """ return self.__timeseries_import_times @property def timeseries_export(self): """ :class:`pi.Timeseries` object for holding the output data. """ return self.__timeseries_export def set_timeseries(self, variable, values, output=True, check_consistency=True, unit=None): if check_consistency: if len(self.times()) != len(values): raise ValueError( 'PIMixin: Trying to set/append values {} with a different ' 'length than the forecast length. Please make sure the ' 'values cover forecastDate through endDate with timestep {}.' .format(variable, self.__timeseries_import.dt)) if unit is None: unit = self.__timeseries_import.get_unit(variable) if output: try: self.__data_config.pi_variable_ids(variable) except KeyError: logger.debug( 'PIMixin: variable {} has no mapping defined in rtcDataConfig ' 'so cannot be added to the output file.'.format(variable)) else: self.__timeseries_export.set(variable, values, unit=unit) self.__timeseries_import.set(variable, values, unit=unit) self.__timeseries_import_dict[variable] = values def get_timeseries(self, variable): return self.__timeseries_import_dict[variable] def extract_results(self): """ Extracts the results of output :returns: An AliasDict of output variables and results array format. """ return self.__output
class ModelicaMixin(OptimizationProblem): """ Adds a `Modelica <http://www.modelica.org/>`_ model to your optimization problem. During preprocessing, the Modelica files located inside the ``model`` subfolder are loaded. :cvar modelica_library_folders: Folders in which any referenced Modelica libraries are to be found. Default is an empty list. """ # Folders in which the referenced Modelica libraries are found modelica_library_folders = [] def __init__(self, **kwargs): # Check arguments assert ('model_folder' in kwargs) # Log pymoca version logger.debug("Using pymoca {}.".format(pymoca.__version__)) # Transfer model from the Modelica .mo file to CasADi using pymoca if 'model_name' in kwargs: model_name = kwargs['model_name'] else: if hasattr(self, 'model_name'): model_name = self.model_name else: model_name = self.__class__.__name__ self.__pymoca_model = pymoca.backends.casadi.api.transfer_model( kwargs['model_folder'], model_name, self.compiler_options()) # Extract the CasADi MX variables used in the model self.__mx = {} self.__mx['time'] = [self.__pymoca_model.time] self.__mx['states'] = [v.symbol for v in self.__pymoca_model.states] self.__mx['derivatives'] = [ v.symbol for v in self.__pymoca_model.der_states ] self.__mx['algebraics'] = [ v.symbol for v in self.__pymoca_model.alg_states ] self.__mx['parameters'] = [ v.symbol for v in self.__pymoca_model.parameters ] self.__mx['control_inputs'] = [] self.__mx['constant_inputs'] = [] self.__mx['lookup_tables'] = [] # Merge with user-specified delayed feedback delayed_feedback_variables = list( map(lambda delayed_feedback: delayed_feedback[1], self.delayed_feedback())) for v in self.__pymoca_model.inputs: if v.symbol.name() in delayed_feedback_variables: # Delayed feedback variables are local to each ensemble, and # therefore belong to the collection of algebraic variables, # rather than to the control inputs. self.__mx['algebraics'].append(v.symbol) else: if v.symbol.name() in kwargs.get('lookup_tables', []): self.__mx['lookup_tables'].append(v.symbol) elif v.fixed: self.__mx['constant_inputs'].append(v.symbol) else: self.__mx['control_inputs'].append(v.symbol) # Initialize nominals and types # These are not in @cached dictionary properties for backwards compatibility. self.__python_types = AliasDict(self.alias_relation) for v in itertools.chain(self.__pymoca_model.states, self.__pymoca_model.alg_states, self.__pymoca_model.inputs): self.__python_types[v.symbol.name()] = v.python_type # Initialize dae and initial residuals # These are not in @cached dictionary properties so that we need to create the list # of function arguments only once. variable_lists = [ 'states', 'der_states', 'alg_states', 'inputs', 'constants', 'parameters' ] function_arguments = [self.__pymoca_model.time] + [ ca.veccat(*[ v.symbol for v in getattr(self.__pymoca_model, variable_list) ]) for variable_list in variable_lists ] self.__dae_residual = self.__pymoca_model.dae_residual_function( *function_arguments) if self.__dae_residual is None: self.__dae_residual = ca.MX() self.__initial_residual = self.__pymoca_model.initial_residual_function( *function_arguments) if self.__initial_residual is None: self.__initial_residual = ca.MX() # Log variables in debug mode if logger.getEffectiveLevel() == logging.DEBUG: logger.debug("ModelicaMixin: Found states {}".format(', '.join( [var.name() for var in self.__mx['states']]))) logger.debug("ModelicaMixin: Found derivatives {}".format( ', '.join([var.name() for var in self.__mx['derivatives']]))) logger.debug("ModelicaMixin: Found algebraics {}".format(', '.join( [var.name() for var in self.__mx['algebraics']]))) logger.debug("ModelicaMixin: Found control inputs {}".format( ', '.join([var.name() for var in self.__mx['control_inputs']]))) logger.debug("ModelicaMixin: Found constant inputs {}".format( ', '.join([var.name() for var in self.__mx['constant_inputs']]))) logger.debug("ModelicaMixin: Found parameters {}".format(', '.join( [var.name() for var in self.__mx['parameters']]))) # Call parent class first for default behaviour. super().__init__(**kwargs) @cached def compiler_options(self) -> Dict[str, Union[str, bool]]: """ Subclasses can configure the `pymoca <http://github.com/pymoca/pymoca>`_ compiler options here. :returns: A dictionary of pymoca compiler options. See the pymoca documentation for details. """ # Default options compiler_options = {} # Expand vector states to multiple scalar component states. compiler_options['expand_vectors'] = True # Where imported model libraries are located. library_folders = self.modelica_library_folders.copy() for ep in pkg_resources.iter_entry_points( group='rtctools.libraries.modelica'): if ep.name == "library_folder": library_folders.append( pkg_resources.resource_filename(ep.module_name, ep.attrs[0])) compiler_options['library_folders'] = library_folders # Eliminate equations of the type 'var = const'. compiler_options['eliminate_constant_assignments'] = True # Eliminate constant symbols from model, replacing them with the values # specified in the model. compiler_options['replace_constant_values'] = True # Replace any constant expressions into the model. compiler_options['replace_constant_expressions'] = True # Replace any parameter expressions into the model. compiler_options['replace_parameter_expressions'] = True # Eliminate variables starting with underscores. compiler_options['eliminable_variable_expression'] = r'(.*[.]|^)_\w+\Z' # Automatically detect and eliminate alias variables. compiler_options['detect_aliases'] = True # Cache the model on disk compiler_options['cache'] = True # Done return compiler_options def delayed_feedback(self): delayed_feedback = super().delayed_feedback() delayed_feedback.extend([(dfb.origin, dfb.name, dfb.delay) for dfb in self.__pymoca_model.delayed_states ]) return delayed_feedback @property def dae_residual(self): return self.__dae_residual @property def dae_variables(self): return self.__mx @property @cached def output_variables(self): output_variables = [ ca.MX.sym(variable) for variable in self.__pymoca_model.outputs ] output_variables.extend(self.__mx['control_inputs']) return output_variables @cached def parameters(self, ensemble_member): # Call parent class first for default values. parameters = super().parameters(ensemble_member) # Return parameter values from pymoca model parameters.update( {v.symbol.name(): v.value for v in self.__pymoca_model.parameters}) # Done return parameters @cached def constant_inputs(self, ensemble_member): # Call parent class first for default values. constant_inputs = super().constant_inputs(ensemble_member) # Return input values from pymoca model times = self.times() constant_input_names = { sym.name() for sym in self.__mx['constant_inputs'] } for v in self.__pymoca_model.inputs: if v.symbol.name() in constant_input_names: constant_inputs[v.symbol.name()] = Timeseries( times, np.full_like(times, v.value)) if logger.getEffectiveLevel() == logging.DEBUG: logger.debug( "Read constant input {} = {} from Modelica model". format(v.symbol.name(), v.value)) return constant_inputs @cached def initial_state(self, ensemble_member): initial_state = AliasDict(self.alias_relation) # Initial conditions obtained from start attributes. for v in self.__pymoca_model.states: if v.fixed: initial_state[v.symbol.name()] = v.start return initial_state @property def initial_residual(self): return self.__initial_residual @cached def bounds(self): # Call parent class first for default values. bounds = super().bounds() # Parameter values parameters = self.parameters(0) parameter_values = [ parameters.get(param.name(), param) for param in self.__mx['parameters'] ] # Load additional bounds from model for v in itertools.chain(self.__pymoca_model.states, self.__pymoca_model.alg_states, self.__pymoca_model.inputs): sym_name = v.symbol.name() try: (m, M) = bounds[sym_name] except KeyError: if self.__python_types.get(sym_name, float) == bool: (m, M) = (0, 1) else: (m, M) = (-np.inf, np.inf) m_ = ca.MX(v.min) if not m_.is_constant(): [m_] = substitute_in_external([m_], self.__mx['parameters'], parameter_values) if not m_.is_constant(): raise Exception( 'Could not resolve lower bound for variable {}'.format( sym_name)) m_ = float(m_) M_ = ca.MX(v.max) if not M_.is_constant(): [M_] = substitute_in_external([M_], self.__mx['parameters'], parameter_values) if not M_.is_constant(): raise Exception( 'Could not resolve upper bound for variable {}'.format( sym_name)) M_ = float(M_) # We take the intersection of all provided bounds m = max(m, m_) M = min(M, M_) bounds[sym_name] = (m, M) return bounds @cached def seed(self, ensemble_member): # Call parent class first for default values. seed = super().seed(ensemble_member) # Parameter values parameters = self.parameters(0) parameter_values = [ parameters.get(param.name(), param) for param in self.__mx['parameters'] ] # Load seeds for var in itertools.chain(self.__pymoca_model.states, self.__pymoca_model.alg_states): start = ca.MX(var.start) if not var.fixed and not start.is_zero(): # If start contains symbolics, try substituting parameter values if not start.is_constant(): [start] = substitute_in_external([start], self.__mx['parameters'], parameter_values) # If start is constant, seed it. Else, warn. sym_name = var.symbol.name() if start.is_constant(): times = self.times(sym_name) start = var.python_type(var.start) s = Timeseries(times, np.full_like(times, start)) if logger.getEffectiveLevel() == logging.DEBUG: logger.debug( "ModelicaMixin: Seeded variable {} = {}".format( sym_name, start)) seed[sym_name] = s else: logger.error( 'ModelicaMixin: Could not resolve seed value for {}'. format(sym_name)) return seed def variable_is_discrete(self, variable): return self.__python_types.get(variable, float) != float @property @cached def alias_relation(self): # Initialize aliases alias_relation = AliasRelation() for v in itertools.chain(self.__pymoca_model.states, self.__pymoca_model.der_states, self.__pymoca_model.alg_states, self.__pymoca_model.inputs): for alias in v.aliases: alias_relation.add(v.symbol.name(), alias) if logger.getEffectiveLevel() == logging.DEBUG: logger.debug("ModelicaMixin: Aliased {} to {}".format( v.symbol.name(), alias)) return alias_relation @property @cached def __nominals(self): # Make the dict nominal_dict = AliasDict(self.alias_relation) # Grab parameters and their values parameters = self.parameters(0) parameter_values = [ parameters.get(param.name(), param) for param in self.__mx['parameters'] ] # Iterate over nominalizable states for v in itertools.chain(self.__pymoca_model.states, self.__pymoca_model.alg_states, self.__pymoca_model.inputs): sym_name = v.symbol.name() # For type consistency, cast to MX nominal = ca.MX(v.nominal) # If nominal contains parameter symbols, substitute them if not nominal.is_constant(): [nominal] = substitute_in_external([nominal], self.__mx['parameters'], parameter_values) if nominal.is_constant(): # Take absolute value (nominal sign is meaningless- a nominal is a magnitude) nominal = ca.fabs(nominal) # If nominal is 0 or 1, we just use the default (1.0) if nominal.is_zero() or (nominal - 1).is_zero(): continue # Cast to numpy array nominal = nominal.to_DM().full().flatten() try: # Attempt to cast to python scalar before storing nominal_dict[sym_name] = nominal.item() except ValueError: # Nominal is numpy array- store it as such nominal_dict[sym_name] = nominal if logger.getEffectiveLevel() == logging.DEBUG: logger.debug( "ModelicaMixin: Set nominal value for variable {} to {}" .format(sym_name, nominal_dict[sym_name])) else: logger.warning( "ModelicaMixin: Could not set nominal value for {}".format( sym_name)) return nominal_dict def variable_nominal(self, variable): return self.__nominals.get(variable, 1)
class SimulationProblem: """ Implements the `BMI <http://csdms.colorado.edu/wiki/BMI_Description>`_ Interface. Base class for all Simulation problems. Loads the Modelica Model. :cvar modelica_library_folders: Folders containing any referenced Modelica libraries. Default is an empty list. """ # Folders in which the referenced Modelica libraries are found modelica_library_folders = [] def __init__(self, **kwargs): # Check arguments assert ('model_folder' in kwargs) # Log pymoca version logger.debug("Using pymoca {}.".format(pymoca.__version__)) # Transfer model from the Modelica .mo file to CasADi using pymoca if 'model_name' in kwargs: model_name = kwargs['model_name'] else: if hasattr(self, 'model_name'): model_name = self.model_name else: model_name = self.__class__.__name__ # Load model from pymoca backend self.__pymoca_model = pymoca.backends.casadi.api.transfer_model( kwargs['model_folder'], model_name, self.compiler_options()) # Extract the CasADi MX variables used in the model self.__mx = {} self.__mx['time'] = [self.__pymoca_model.time] self.__mx['states'] = [v.symbol for v in self.__pymoca_model.states] self.__mx['derivatives'] = [ v.symbol for v in self.__pymoca_model.der_states ] self.__mx['algebraics'] = [ v.symbol for v in self.__pymoca_model.alg_states ] self.__mx['parameters'] = [ v.symbol for v in self.__pymoca_model.parameters ] self.__mx['constant_inputs'] = [] self.__mx['lookup_tables'] = [] # TODO: implement delayed feedback delayed_feedback_variables = [] for v in self.__pymoca_model.inputs: if v.symbol.name() in delayed_feedback_variables: # Delayed feedback variables are local to each ensemble, and # therefore belong to the collection of algebraic variables, # rather than to the control inputs. self.__mx['algebraics'].append(v.symbol) else: if v.symbol.name() in kwargs.get('lookup_tables', []): self.__mx['lookup_tables'].append(v.symbol) else: # All inputs are constant inputs self.__mx['constant_inputs'].append(v.symbol) # Log variables in debug mode if logger.getEffectiveLevel() == logging.DEBUG: logger.debug("SimulationProblem: Found states {}".format(', '.join( [var.name() for var in self.__mx['states']]))) logger.debug("SimulationProblem: Found derivatives {}".format( ', '.join([var.name() for var in self.__mx['derivatives']]))) logger.debug("SimulationProblem: Found algebraics {}".format( ', '.join([var.name() for var in self.__mx['algebraics']]))) logger.debug("SimulationProblem: Found constant inputs {}".format( ', '.join([var.name() for var in self.__mx['constant_inputs']]))) logger.debug("SimulationProblem: Found parameters {}".format( ', '.join([var.name() for var in self.__mx['parameters']]))) # Initialize an AliasDict for nominals and types self.__nominals = AliasDict(self.alias_relation) self.__python_types = AliasDict(self.alias_relation) for v in itertools.chain(self.__pymoca_model.states, self.__pymoca_model.alg_states, self.__pymoca_model.inputs): sym_name = v.symbol.name() # Store the types in an AliasDict self.__python_types[sym_name] = v.python_type # If the nominal is 0.0 or 1.0 or -1.0, ignore: get_variable_nominal returns a default of 1.0 # TODO: handle nominal vectors (update() will need to load them) if ca.MX(v.nominal).is_zero() or ca.MX( v.nominal - 1).is_zero() or ca.MX(v.nominal + 1).is_zero(): continue else: if ca.MX(v.nominal).size1() != 1: logger.error( 'Vector Nominals not supported yet. ({})'.format( sym_name)) self.__nominals[sym_name] = ca.fabs(v.nominal) if logger.getEffectiveLevel() == logging.DEBUG: logger.debug( "SimulationProblem: Setting nominal value for variable {} to {}" .format(sym_name, self.__nominals[sym_name])) # Initialize DAE and initial residuals variable_lists = [ 'states', 'der_states', 'alg_states', 'inputs', 'constants', 'parameters' ] function_arguments = [self.__pymoca_model.time] + [ ca.veccat(*[ v.symbol for v in getattr(self.__pymoca_model, variable_list) ]) for variable_list in variable_lists ] self.__dae_residual = self.__pymoca_model.dae_residual_function( *function_arguments) self.__initial_residual = self.__pymoca_model.initial_residual_function( *function_arguments) if self.__initial_residual is None: self.__initial_residual = ca.MX() # Construct state vector self.__sym_list = self.__mx['states'] + self.__mx['algebraics'] + self.__mx['derivatives'] + \ self.__mx['time'] + self.__mx['constant_inputs'] + self.__mx['parameters'] self.__state_vector = np.full(len(self.__sym_list), np.nan) # A very handy index self.__states_end_index = len(self.__mx['states']) + \ len(self.__mx['algebraics']) + len(self.__mx['derivatives']) # Construct a dict to look up symbols by name (or iterate over) self.__sym_dict = OrderedDict( ((sym.name(), sym) for sym in self.__sym_list)) # Assemble some symbolics, including those needed for a backwards Euler derivative approximation X = ca.vertcat(*self.__sym_list[:self.__states_end_index]) X_prev = ca.vertcat(*[ ca.MX.sym(sym.name() + '_prev') for sym in self.__sym_list[:self.__states_end_index] ]) dt = ca.MX.sym("delta_t") # Make a list of derivative approximations using backwards Euler formulation derivative_approximation_residuals = [] for index, derivative_state in enumerate(self.__mx['derivatives']): derivative_approximation_residuals.append( derivative_state - (X[index] - X_prev[index]) / dt) # Append residuals for derivative approximations dae_residual = ca.vertcat(self.__dae_residual, *derivative_approximation_residuals) # TODO: implement lookup_tables # Make a list of unscaled symbols and a list of their scaled equivalent unscaled_symbols = [] scaled_symbols = [] for sym_name, nominal in self.__nominals.items(): index = self.__get_state_vector_index(sym_name) # If the symbol is a state, Add the symbol to the lists if index <= self.__states_end_index: unscaled_symbols.append(X[index]) scaled_symbols.append(X[index] * nominal) # Also scale previous states unscaled_symbols.append(X_prev[index]) scaled_symbols.append(X_prev[index] * nominal) # Substitute unscaled terms for scaled terms dae_residual = ca.substitute(dae_residual, ca.vertcat(*unscaled_symbols), ca.vertcat(*scaled_symbols)) if logger.getEffectiveLevel() == logging.DEBUG: logger.debug('SimulationProblem: DAE Residual is ' + str(dae_residual)) if X.size1() != dae_residual.size1(): logger.error( 'Formulation Error: Number of states ({}) does not equal number of equations ({})' .format(X.size1(), dae_residual.size1())) # Construct function parameters parameters = ca.vertcat(dt, X_prev, *self.__sym_list[self.__states_end_index:]) # Construct a function res_vals that returns the numerical residuals of a numerical state self.__res_vals = ca.Function("res_vals", [X, parameters], [dae_residual]) # Use rootfinder() to make a function that takes a step forward in time by trying to zero res_vals() options = {'nlpsol': 'ipopt', 'nlpsol_options': self.solver_options()} self.__do_step = ca.rootfinder("next_state", "nlpsol", self.__res_vals, options) # Call parent class for default behaviour. super().__init__() def initialize(self, config_file=None): """ Initialize state vector with default values :param config_file: Path to an initialization file. """ if config_file: # TODO read start and stop time from config_file and call: # self.setup_experiment(start,stop) # for now, assume that setup_experiment was called beforehand raise NotImplementedError # Set values of parameters defined in the model into the state vector for var in self.__pymoca_model.parameters: # First check to see if parameter is already set (this allows child classes to override model defaults) if np.isfinite(self.get_var(var.symbol.name())): continue # Also test to see if the value is constant if isinstance(var.value, ca.MX) and not var.value.is_constant(): continue # Try to extract the value try: # Extract the value as a python type val = var.python_type(var.value) except ValueError: # var.value is a float NaN being cast to non-float continue else: # If val is finite, we set it if np.isfinite(val): logger.debug( 'SimulationProblem: Setting parameter {} = {}'.format( var.symbol.name(), val)) self.set_var(var.symbol.name(), val) # Assemble initial residuals and set values from start attributes into the state vector constrained_residuals = [] minimized_residuals = [] for var in itertools.chain(self.__pymoca_model.states, self.__pymoca_model.alg_states): var_name = var.symbol.name() var_nominal = self.get_variable_nominal(var_name) # Attempt to cast var.start to python type mx_start = ca.MX(var.start) if mx_start.is_constant(): # cast var.start to python type start_val = var.python_type(mx_start.to_DM()) else: # var.start is a symbolic expression with unknown value start_val = None if start_val == 0.0 and not var.fixed: # To make initialization easier, we allow setting initial states by providing timeseries # with names that match a symbol in the model. We only check for this matching if the start # and fixed attributes were left as default try: start_val = self.initial_state()[var_name] except KeyError: pass else: # An initial state was found- add it to the constrained residuals logger.debug( 'Initialize: Added {} = {} to initial equations (found matching timeseries).' .format(var_name, start_val)) # Set var to be fixed var.fixed = True # Attempt to set start_val in the state vector. Default to zero if unknown. try: self.set_var(var_name, start_val if start_val is not None else 0.0) except KeyError: logger.warning( 'Initialize: {} not found in state vector. Initial value of {} not set.' .format(var_name, start_val)) # Add a residual for the difference between the state and its starting expression start_expr = start_val if start_val is not None else var.start if var.fixed: # require residual = 0 constrained_residuals.append( (var.symbol - start_expr) / var_nominal) else: # minimize residual minimized_residuals.append( (var.symbol - start_expr) / var_nominal) # Default start var for ders is zero for der_var in self.__mx['derivatives']: self.set_var(der_var.name(), 0.0) # Warn for nans in state vector (verify we didn't miss anything) self.__warn_for_nans() # Optionally encourage a steady-state initial condition if getattr(self, 'encourage_steady_state_initial_conditions', False): # add penalty for der(var) != 0.0 for d in self.__mx['derivatives']: logger.debug('Added {} to the minimized residuals.'.format( d.name())) minimized_residuals.append(d) # Make minimized_residuals into a single symbolic object minimized_residual = ca.vertcat(*minimized_residuals) # Assemble symbolics needed to make a function describing the initial condition of the model # We constrain every entry in this MX to zero equality_constraints = ca.vertcat(self.__dae_residual, self.__initial_residual, *constrained_residuals) # The variables that need a mutually consistent initial condition X = ca.vertcat(*self.__sym_list[:self.__states_end_index]) # Make a list of unscaled symbols and a list of their scaled equivalent unscaled_symbols = [] scaled_symbols = [] for sym_name, nominal in self.__nominals.items(): # Add the symbol to the lists symbol = self.__sym_dict[sym_name] unscaled_symbols.append(symbol) scaled_symbols.append(symbol * nominal) # Make the lists symbolic unscaled_symbols = ca.vertcat(*unscaled_symbols) scaled_symbols = ca.vertcat(*scaled_symbols) # Substitute unscaled terms for scaled terms equality_constraints = ca.substitute(equality_constraints, unscaled_symbols, scaled_symbols) minimized_residual = ca.substitute(minimized_residual, unscaled_symbols, scaled_symbols) logger.debug('SimulationProblem: Initial Equations are ' + str(equality_constraints)) logger.debug('SimulationProblem: Minimized Residuals are ' + str(minimized_residual)) # State bounds can be symbolic, written in terms of parameters. After all # parameter values are known, we evaluate the numeric values of bounds. symbolic_bounds = ca.vertcat(*[ ca.horzcat(v.min, v.max) for v in itertools.chain( self.__pymoca_model.states, self.__pymoca_model.alg_states, self.__pymoca_model.der_states) ]) bound_evaluator = ca.Function('bound_evaluator', self.__mx['parameters'], [symbolic_bounds]) # Evaluate bounds using values of parameters n_parameters = len(self.__mx['parameters']) if n_parameters > 0: [evaluated_bounds ] = bound_evaluator.call(self.__state_vector[-n_parameters:]) else: [evaluated_bounds] = bound_evaluator.call([]) # Construct arrays of state bounds (used in the initialize() nlp, but not in __do_step rootfinder) self.__lbx = evaluated_bounds[:, 0] self.__ubx = evaluated_bounds[:, 1] # Constrain model equation residuals to zero lbg = np.zeros(equality_constraints.size1()) ubg = np.zeros(equality_constraints.size1()) # Construct objective function from the input residual objective_function = ca.dot(minimized_residual, minimized_residual) # Construct nlp and solver to find initial state using ipopt parameters = ca.vertcat(*self.__mx['time'], *self.__mx['constant_inputs'], *self.__mx['parameters']) nlp = { 'x': X, 'f': objective_function, 'g': equality_constraints, 'p': parameters } solver = ca.nlpsol('solver', 'ipopt', nlp, self.solver_options()) # Construct guess guess = ca.vertcat( *np.nan_to_num(self.__state_vector[:self.__states_end_index])) # Find initial state initial_state = solver(x0=guess, lbx=self.__lbx, ubx=self.__ubx, lbg=lbg, ubg=ubg, p=self.__state_vector[self.__states_end_index:]) # If unsuccessful, stop. return_status = solver.stats()['return_status'] if return_status not in { 'Solve_Succeeded', 'Solved_To_Acceptable_Level' }: raise Exception( 'Initialization Failed with return status "{}"'.format( return_status)) # Update state vector with initial conditions self.__state_vector[:self.__states_end_index] = initial_state[ 'x'][:self.__states_end_index].T # make a copy of the initialized initial state vector in case we want to run the model again self.__initialized_state_vector = copy.deepcopy(self.__state_vector) # Warn for nans in state vector after initialization self.__warn_for_nans() def pre(self): """ Any preprocessing takes place here. """ pass def post(self): """ Any postprocessing takes place here. """ pass def setup_experiment(self, start, stop, dt): """ Method for subclasses (PIMixin, CSVMixin, or user classes) to set timing information for a simulation run. :param start: Start time for the simulation. :param stop: Final time for the simulation. :param dt: Time step size. """ # Set class vars with start/stop/dt values self.__start = start self.__stop = stop self.__dt = dt # Set time in state vector self.set_var('time', start) def update(self, dt): """ Performs one timestep. The methods ``setup_experiment`` and ``initialize`` must have been called before. :param dt: Time step size. """ if dt < 0: dt = self.__dt logger.debug("Taking a step at {} with size {}".format( self.get_current_time(), dt)) # increment time self.set_var('time', self.get_current_time() + dt) # take a step guess = self.__state_vector[:self.__states_end_index] next_state = self.__do_step(guess, ca.vertcat(dt, *self.__state_vector)) # make sure that the step converged sufficiently largest_res = ca.norm_inf( self.__res_vals(next_state, ca.vertcat(self.__dt, *self.__state_vector))) tol = self.solver_options().get('ipopt.tol', 1.0e-8) if largest_res > tol: logger.warning( 'Simulation may have failed to converge at time {}. Residual value {} is greater than {}' .format(self.get_current_time(), largest_res, tol)) if logger.getEffectiveLevel() == logging.DEBUG: logger.debug('Residual maximum magnitude: {:.2E}'.format( float(largest_res))) # Update state vector self.__state_vector[:self.__states_end_index] = next_state.T def simulate(self): """ Run model from start_time to end_time. """ # Do any preprocessing, which may include changing parameter values on # the model logger.info("Preprocessing") self.pre() # Initialize model logger.info("Initializing") self.initialize() # Perform all timesteps logger.info("Running") while self.get_current_time() < self.get_end_time(): self.update(-1) # Do any postprocessing logger.info("Postprocessing") self.post() def reset(self): """ Reset the FMU. """ self.__state_vector = copy.deepcopy(self.__initialized_state_vector) def get_start_time(self): """ Return start time of experiment. :returns: The start time of the experiment. """ return self.__start def get_end_time(self): """ Return end time of experiment. :returns: The end time of the experiment. """ return self.__stop def get_current_time(self): """ Return current time of simulation. :returns: The current simulation time. """ return self.get_var('time') def get_time_step(self): """ Return simulation timestep. :returns: The simulation timestep. """ return self.__dt def get_var(self, name): """ Return a numpy array from FMU. :param name: Variable name. :returns: The value of the variable. """ # Get the canonical name and sign name, sign = self.alias_relation.canonical_signed(name) # Get the raw value of the canonical var index = self.__get_state_vector_index(name) value = self.__state_vector[index] # Adjust sign if needed if sign < 0: value *= sign # Adjust for nominal value if not default nominal = self.get_variable_nominal(name) if nominal != 1.0: value *= nominal return value def get_var_count(self): """ Return the number of variables in the model. :returns: The number of variables in the model. """ return len(self.get_variables()) def get_var_name(self, i): """ Returns the name of a variable. :param i: Index in ordered dictionary returned by method get_variables. :returns: The name of the variable. """ return list(self.get_variables())[i] def get_var_type(self, name): """ Return type, compatible with numpy. :param name: String variable name. :returns: The numpy-compatible type of the variable. :raises: KeyError """ return self.__python_types[name] def get_var_rank(self, name): """ Not implemented """ raise NotImplementedError def get_var_shape(self, name): """ Not implemented """ raise NotImplementedError def get_variables(self): """ Return all variables (both internal and user defined) :returns: An ordered dictionary of all variables supported by the model. """ return self.__sym_dict @cached def get_state_variables(self): return AliasDict( self.alias_relation, { sym.name(): sym for sym in (self.__mx['states'] + self.__mx['algebraics']) }) @cached def get_parameter_variables(self): return AliasDict(self.alias_relation, {sym.name(): sym for sym in self.__mx['parameters']}) @cached def get_input_variables(self): return AliasDict( self.alias_relation, {sym.name(): sym for sym in self.__mx['constant_inputs']}) @cached def get_output_variables(self): return self.__pymoca_model.outputs @cached def __get_state_vector_index(self, variable): index = next((i for i, sym in enumerate(self.__sym_list) if sym.name() == variable), None) if index is None: raise KeyError(str(variable) + " does not exist!") return index def __warn_for_nans(self): """ Test state vector for missing values and warn """ value_is_nan = np.isnan(self.__state_vector) if any(value_is_nan): for sym, isnan in zip(self.__sym_list, value_is_nan): if isnan: logger.warning('Variable {} has no value.'.format(sym)) def set_var(self, name, value): """ Set the value of the given variable. :param name: Name of variable to set. :param value: Value(s). """ # TODO: sanitize input # Get the canonical name, adjust sign if needed name, sign = self.alias_relation.canonical_signed(name) if sign < 0: value *= sign # Adjust for nominal value if not default nominal = self.get_variable_nominal(name) if nominal != 1.0: value /= nominal # Store value in state vector index = self.__get_state_vector_index(name) self.__state_vector[index] = value def set_var_slice(self, name, start, count, var): """ Not implemented. """ raise NotImplementedError def set_var_index(self, name, index, var): """ Not implemented. """ raise NotImplementedError def inq_compound(self, name): """ Not implemented. """ raise NotImplementedError def inq_compound_field(self, name, index): """ Not implemented. """ raise NotImplementedError def solver_options(self): """ Returns a dictionary of CasADi root_finder() solver options. :returns: A dictionary of CasADi :class:`root_finder` options. See the CasADi documentation for details. """ return {'ipopt.print_level': 0, 'print_time': False} def get_variable_nominal(self, variable): """ Get the value of the nominal attribute of a variable """ return self.__nominals.get(variable, 1.0) def timeseries_at(self, variable, t): """ Get value of timeseries variable at time t: should be overridden by pi or csv mixin """ raise NotImplementedError @cached def initial_state(self) -> AliasDict: """ The initial state. :returns: A dictionary of variable names and initial state (t0) values. """ t0 = self.get_start_time() initial_state_dict = AliasDict(self.alias_relation) for variable in list(self.get_state_variables()) + list( self.get_input_variables()): try: initial_state_dict[variable] = self.timeseries_at(variable, t0) except KeyError: pass except NotImplementedError: pass else: if logger.getEffectiveLevel() == logging.DEBUG: logger.debug("Read intial state for {}".format(variable)) return initial_state_dict @cached def parameters(self): """ Return a dictionary of parameter values extracted from the Modelica model """ # Create AliasDict parameters = AliasDict(self.alias_relation) # Update with model parameters parameters.update( {p.symbol.name(): p.value for p in self.__pymoca_model.parameters}) return parameters @property @cached def alias_relation(self): # Initialize aliases alias_relation = AliasRelation() for v in itertools.chain(self.__pymoca_model.states, self.__pymoca_model.der_states, self.__pymoca_model.alg_states, self.__pymoca_model.inputs): for alias in v.aliases: alias_relation.add(v.symbol.name(), alias) if logger.getEffectiveLevel() == logging.DEBUG: logger.debug("SimulationProblem: Aliased {} to {}".format( v.symbol.name(), alias)) return alias_relation @cached def compiler_options(self): """ Subclasses can configure the `pymoca <http://github.com/pymoca/pymoca>`_ compiler options here. :returns: A dictionary of pymoca compiler options. See the pymoca documentation for details. """ # Default options compiler_options = {} # Expand vector states to multiple scalar component states. compiler_options['expand_vectors'] = True # Where imported model libraries are located. library_folders = self.modelica_library_folders.copy() for ep in pkg_resources.iter_entry_points( group='rtctools.libraries.modelica'): if ep.name == "library_folder": library_folders.append( pkg_resources.resource_filename(ep.module_name, ep.attrs[0])) compiler_options['library_folders'] = library_folders # Eliminate equations of the type 'var = const'. compiler_options['eliminate_constant_assignments'] = True # Eliminate constant symbols from model, replacing them with the values # specified in the model. compiler_options['replace_constant_values'] = True # Replace any constant expressions into the model. compiler_options['replace_constant_expressions'] = True # Replace any parameter expressions into the model. compiler_options['replace_parameter_expressions'] = True # Eliminate variables starting with underscores. compiler_options['eliminable_variable_expression'] = r'(.*[.]|^)_\w+\Z' # Automatically detect and eliminate alias variables. compiler_options['detect_aliases'] = True # Cache the model on disk compiler_options['cache'] = True # Done return compiler_options