def initialize_model(self): """ This function will initialize the model using the profile obtained from simulating the dynamic model. """ if self._tsim is None: raise DAE_Error( "Tried to initialize the model without simulating it first") tvals = list(self._contset) # Build list of state and algebraic variables # that can be initialized initvars = self._diffvars + self._simalgvars for idx, v in enumerate(initvars): for idx2, i in enumerate(v._args): if type(i) is IndexTemplate: break valinit = np.interp(tvals, self._tsim, self._simsolution[:, idx]) for i, t in enumerate(tvals): vidx = tuple(v._args[0:idx2]) + (t,) + \ tuple(v._args[idx2 + 1:]) v._base[vidx] = valinit[i]
def _simulate_with_scipy(self, initcon, tsim, switchpts, varying_inputs, integrator, integrator_options): scipyint = \ scipy.ode(self._rhsfun).set_integrator(integrator, **integrator_options) scipyint.set_initial_value(initcon, tsim[0]) profile = np.array(initcon) i = 1 while scipyint.successful() and scipyint.t < tsim[-1]: # check if tsim[i-1] is a switching time and update value if tsim[i - 1] in switchpts: for v in self._siminputvars.keys(): if tsim[i - 1] in varying_inputs[v]: p = self._templatemap[self._siminputvars[v]] p.set_value(varying_inputs[v][tsim[i - 1]]) profilestep = scipyint.integrate(tsim[i]) profile = np.vstack([profile, profilestep]) i += 1 if not scipyint.successful(): raise DAE_Error("The Scipy integrator %s did not terminate " "successfully." % integrator) return [tsim, profile]
def convert_pyomo2scipy(expr, templatemap): """Substitute _GetItem nodes in an expression tree. This substition function is used to replace Pyomo _GetItem nodes with mutable Params. Args: templatemap: dictionary mapping _GetItemIndexer objects to mutable params Returns: a new expression tree with all substitutions done """ if not scipy_available: raise DAE_Error("SciPy is not installed. Cannot substitute SciPy " "intrinsic functions.") visitor = Pyomo2Scipy_Visitor(templatemap) return visitor.dfs_postorder_stack(expr)
def substitute_pyomo2casadi(expr, templatemap): """Substitute IndexTemplates in an expression tree. This substition function is used to replace Pyomo intrinsic functions with CasADi functions. Args: expr: a Pyomo expression templatemap: dictionary mapping _GetItemIndexer objects to mutable params Returns: a new expression tree with all substitutions done """ if not casadi_available: raise DAE_Error("CASADI is not installed. Cannot substitute CasADi " "variables and intrinsic functions.") visitor = Substitute_Pyomo2Casadi_Visitor(templatemap) return visitor.dfs_postorder_stack(expr)
def convert_pyomo2casadi(expr): """Convert a Pyomo expression tree to Casadi. This function replaces a Pyomo expression with a CasADi expression. This assumes that the `substitute_pyomo2casadi` function has been called, so the Pyomo expression contains CasADi variables and intrinsic functions. The resulting expression can be used with the CasADi integrator. Args: expr: a Pyomo expression with CasADi variables and intrinsic functions Returns: a CasADi expression tree. """ if not casadi_available: raise DAE_Error("CASADI is not installed. Cannot convert a Pyomo " "expression to a Casadi expression.") visitor = Convert_Pyomo2Casadi_Visitor() return visitor.dfs_postorder_stack(expr)
def _transformBlock(self, block, currentds): self._fe = {} for ds in block.component_objects(ContinuousSet, descend_into=True): if currentds is None or currentds == ds.name: if 'scheme' in ds.get_discretization_info(): raise DAE_Error( "Attempting to discretize ContinuousSet " "'%s' after it has already been discretized. " % ds.name) generate_finite_elements(ds, self._nfe[currentds]) if not ds.get_changed(): if len(ds) - 1 > self._nfe[currentds]: logger.warning( "More finite elements were found in " "ContinuousSet '%s' than the number of " "finite elements specified in apply. The " "larger number of finite elements will be " "used." % ds.name) self._nfe[ds.name] = len(ds) - 1 self._fe[ds.name] = list(ds) generate_colloc_points(ds, self._tau[currentds]) # Adding discretization information to the continuousset # object itself so that it can be accessed outside of the # discretization object disc_info = ds.get_discretization_info() disc_info['nfe'] = self._nfe[ds.name] disc_info['ncp'] = self._ncp[currentds] disc_info['tau_points'] = self._tau[currentds] disc_info['adot'] = self._adot[currentds] disc_info['adotdot'] = self._adotdot[currentds] disc_info['afinal'] = self._afinal[currentds] disc_info['scheme'] = self._scheme_name expand_components(block) for d in block.component_objects(DerivativeVar, descend_into=True): dsets = d.get_continuousset_list() for i in ComponentSet(dsets): if currentds is None or i.name == currentds: oldexpr = d.get_derivative_expression() loc = d.get_state_var()._contset[i] count = dsets.count(i) if count >= 3: raise DAE_Error( "Error discretizing '%s' with respect to '%s'. " "Current implementation only allows for taking the" " first or second derivative with respect to a " "particular ContinuousSet" % (d.name, i.name)) scheme = self._scheme[count - 1] newexpr = create_partial_expression( scheme, oldexpr, i, loc) d.set_derivative_expression(newexpr) if self._scheme_name == 'LAGRANGE-LEGENDRE': # Add continuity equations to DerivativeVar's parent # block add_continuity_equations(d.parent_block(), d, i, loc) # Reclassify DerivativeVar if all indexing ContinuousSets have # been discretized. Add discretization equations to the # DerivativeVar's parent block. if d.is_fully_discretized(): add_discretization_equations(d.parent_block(), d) d.parent_block().reclassify_component_type(d, Var) # Keep track of any reclassified DerivativeVar components so # that the Simulator can easily identify them if the model # is simulated after discretization # TODO: Update the discretization transformations to use # a Block to add things to the model and store discretization # information. Using a list for now because the simulator # does not yet support models containing active Blocks reclassified_list = getattr( block, '_pyomo_dae_reclassified_derivativevars', None) if reclassified_list is None: block._pyomo_dae_reclassified_derivativevars = list() reclassified_list = \ block._pyomo_dae_reclassified_derivativevars reclassified_list.append(d) # Reclassify Integrals if all ContinuousSets have been discretized if block_fully_discretized(block): if block.contains_component(Integral): for i in block.component_objects(Integral, descend_into=True): i.parent_block().reclassify_component_type(i, Expression) # TODO: The following reproduces the old behavior of # "reconstruct()". We should come up with an # implementation that does not rely on manipulating # private attributes i.clear() i._constructed = False i.construct() # If a model contains integrals they are most likely to appear # in the objective function which will need to be reconstructed # after the model is discretized. for k in block.component_objects(Objective, descend_into=True): # TODO: check this, reconstruct might not work # TODO: The following reproduces the old behavior of # "reconstruct()". We should come up with an # implementation that does not rely on manipulating # private attributes k.clear() k._constructed = False k.construct()
def simulate(self, numpoints=None, tstep=None, integrator=None, varying_inputs=None, initcon=None, integrator_options=None): """ Simulate the model. Integrator-specific options may be specified as keyword arguments and will be passed on to the integrator. Parameters ---------- numpoints : int The number of points for the profiles returned by the simulator. Default is 100 tstep : int or float The time step to use in the profiles returned by the simulator. This is not the time step used internally by the integrators. This is an optional parameter that may be specified in place of 'numpoints'. integrator : string The string name of the integrator to use for simulation. The default is 'lsoda' when using Scipy and 'idas' when using CasADi varying_inputs : ``pyomo.environ.Suffix`` A :py:class:`Suffix<pyomo.environ.Suffix>` object containing the piecewise constant profiles to be used for certain time-varying algebraic variables. initcon : list of floats The initial conditions for the the differential variables. This is an optional argument. If not specified then the simulator will use the current value of the differential variables at the lower bound of the ContinuousSet for the initial condition. integrator_options : dict Dictionary containing options that should be passed to the integrator. See the documentation for a specific integrator for a list of valid options. Returns ------- numpy array, numpy array The first return value is a 1D array of time points corresponding to the second return value which is a 2D array of the profiles for the simulated differential and algebraic variables. """ if not numpy_available: raise ValueError("The numpy module is not available. " "Cannot simulate the model.") if integrator_options is None: integrator_options = {} if self._intpackage == 'scipy': # Specify the scipy integrator to use for simulation valid_integrators = ['vode', 'zvode', 'lsoda', 'dopri5', 'dop853'] if integrator is None: integrator = 'lsoda' elif integrator == 'odeint': integrator = 'lsoda' else: # Specify the casadi integrator to use for simulation. # Only a subset of these integrators may be used for # DAE simulation. We defer this check to CasADi. valid_integrators = ['cvodes', 'idas', 'collocation', 'rk'] if integrator is None: integrator = 'idas' if integrator not in valid_integrators: raise DAE_Error("Unrecognized %s integrator \'%s\'. Please select" " an integrator from %s" % (self._intpackage, integrator, valid_integrators)) # Set the time step or the number of points for the lists # returned by the integrator if tstep is not None and \ tstep > (self._contset.last() - self._contset.first()): raise ValueError( "The step size %6.2f is larger than the span of the " "ContinuousSet %s" % (tstep, self._contset.name())) if tstep is not None and numpoints is not None: raise ValueError( "Cannot specify both the step size and the number of " "points for the simulator") if tstep is None and numpoints is None: # Use 100 points by default numpoints = 100 if tstep is None: tsim = np.linspace( self._contset.first(), self._contset.last(), num=numpoints) # Consider adding an option for log spaced time points. Can be # important for simulating stiff systems. # tsim = np.logspace(-4,6, num=100) # np.log10(self._contset.first()),np.log10( # self._contset.last()),num=1000, endpoint=True) else: tsim = np.arange( self._contset.first(), self._contset.last(), tstep) switchpts = [] self._siminputvars = {} self._simalgvars = [] if varying_inputs is not None: if type(varying_inputs) is not Suffix: raise TypeError( "Varying input values must be specified using a " "Suffix. Please refer to the simulator documentation.") for alg in self._algvars: if alg._base in varying_inputs: # Find all the switching points switchpts += varying_inputs[alg._base].keys() # Add to dictionary of siminputvars self._siminputvars[alg._base] = alg else: self._simalgvars.append(alg) if self._intpackage == 'scipy' and len(self._simalgvars) != 0: raise DAE_Error("When simulating with Scipy you must " "provide values for all parameters " "and algebraic variables that are indexed " "by the ContinuoutSet using the " "'varying_inputs' keyword argument. " "Please refer to the simulator documentation " "for more information.") # Get the set of unique points switchpts = list(set(switchpts)) switchpts.sort() # Make sure all the switchpts are within the bounds of # the ContinuousSet if switchpts[0] < self._contset.first() or \ switchpts[-1] > self._contset.last(): raise ValueError("Found a switching point for one or more of " "the time-varying inputs that is not within " "the bounds of the ContinuousSet.") # Update tsim to include input switching points # This numpy function returns the unique, sorted points tsim = np.union1d(tsim, switchpts) else: self._simalgvars = self._algvars # Check if initial conditions were provided, otherwise obtain # them from the current variable values if initcon is not None: if len(initcon) > len(self._diffvars): raise ValueError( "Too many initial conditions were specified. The " "simulator was expecting a list with %i values." % len(self._diffvars)) if len(initcon) < len(self._diffvars): raise ValueError( "Too few initial conditions were specified. The " "simulator was expecting a list with %i values." % len(self._diffvars)) else: initcon = [] for v in self._diffvars: for idx, i in enumerate(v._args): if type(i) is IndexTemplate: break initpoint = self._contset.first() vidx = tuple(v._args[0:idx]) + (initpoint,) + \ tuple(v._args[idx + 1:]) # This line will raise an error if no value was set initcon.append(value(v._base[vidx])) # Call the integrator if self._intpackage == 'scipy': if not scipy_available: raise ValueError("The scipy module is not available. " "Cannot simulate the model.") if is_pypy: raise ValueError("The scipy ODE integrators do not work " "under pypy. Cannot simulate the model.") tsim, profile = self._simulate_with_scipy(initcon, tsim, switchpts, varying_inputs, integrator, integrator_options) else: if len(switchpts) != 0: tsim, profile = \ self._simulate_with_casadi_with_inputs(initcon, tsim, varying_inputs, integrator, integrator_options) else: tsim, profile = \ self._simulate_with_casadi_no_inputs(initcon, tsim, integrator, integrator_options) self._tsim = tsim self._simsolution = profile return [tsim, profile]
def __init__(self, m, package='scipy'): self._intpackage = package if self._intpackage not in ['scipy', 'casadi']: raise DAE_Error( "Unrecognized simulator package %s. Please select from " "%s" % (self._intpackage, ['scipy', 'casadi'])) if self._intpackage == 'scipy': if not scipy_available: # Converting this to a warning so that Simulator initialization # can be tested even when scipy is unavailable logger.warning( "The scipy module is not available. " "You may build the Simulator object but you will not " "be able to run the simulation.") elif is_pypy: logger.warning( "The scipy ODE integrators do not work in pypy. " "You may build the Simulator object but you will not " "be able to run the simulation.") else: if not casadi_available: # Initializing the simulator for use with casadi requires # access to casadi objects. Therefore, we must throw an error # here instead of a warning. raise ValueError("The casadi module is not available. " "Cannot simulate model.") # Check for active Blocks and throw error if any are found if len(list(m.component_data_objects(Block, active=True, descend_into=False))): raise DAE_Error("The Simulator cannot handle hierarchical models " "at the moment.") temp = m.component_map(ContinuousSet) if len(temp) != 1: raise DAE_Error( "Currently the simulator may only be applied to " "Pyomo models with a single ContinuousSet") # Get the ContinuousSet in the model contset = list(temp.values())[0] # Create a index template for the continuous set cstemplate = IndexTemplate(contset) # Ensure that there is at least one derivative in the model derivs = m.component_map(DerivativeVar) derivs = list(derivs.keys()) if hasattr(m, '_pyomo_dae_reclassified_derivativevars'): for d in m._pyomo_dae_reclassified_derivativevars: derivs.append(d.name) if len(derivs) == 0: raise DAE_Error("Cannot simulate a model with no derivatives") templatemap = {} # Map for template substituter rhsdict = {} # Map of derivative to its RHS templated expr derivlist = [] # Ordered list of derivatives alglist = [] # list of templated algebraic equations # Loop over constraints to find differential equations with separable # RHS. Must find a RHS for every derivative var otherwise ERROR. Build # dictionary of DerivativeVar:RHS equation. for con in m.component_objects(Constraint, active=True): # Skip the discretization equations if model is discretized if '_disc_eq' in con.name: continue # Check dimension of the Constraint. Check if the # Constraint is indexed by the continuous set and # determine its order in the indexing sets if con.dim() == 0: continue conindex = con.index_set() if not hasattr(conindex, 'set_tuple'): # Check if the continuous set is the indexing set if conindex is not contset: continue else: csidx = 0 noncsidx = (None,) else: temp = conindex.set_tuple dimsum = 0 csidx = -1 noncsidx = None for s in temp: if s is contset: if csidx != -1: raise DAE_Error( "Cannot simulate the constraint %s because " "it is indexed by duplicate ContinuousSets" % con.name) csidx = dimsum elif noncsidx is None: noncsidx = s else: noncsidx = noncsidx.cross(s) dimsum += s.dimen if csidx == -1: continue # Get the rule used to construct the constraint conrule = con.rule for i in noncsidx: # Insert the index template and call the rule to # create a templated expression if i is None: tempexp = conrule(m, cstemplate) else: if not isinstance(i, tuple): i = (i,) tempidx = i[0:csidx] + (cstemplate,) + i[csidx:] tempexp = conrule(m, *tempidx) # Check to make sure it's an EqualityExpression if not type(tempexp) is EXPR.EqualityExpression: continue # Check to make sure it's a differential equation with # separable RHS args = None # Case 1: m.dxdt[t] = RHS if type(tempexp.arg(0)) is EXPR.GetItemExpression: args = _check_getitemexpression(tempexp, 0) # Case 2: RHS = m.dxdt[t] if args is None: if type(tempexp.arg(1)) is EXPR.GetItemExpression: args = _check_getitemexpression(tempexp, 1) # Case 3: m.p*m.dxdt[t] = RHS if args is None: if type(tempexp.arg(0)) is EXPR.ProductExpression or \ type(tempexp.arg(0)) is EXPR.ReciprocalExpression or \ type(tempexp.arg(0)) is EXPR.DivisionExpression: args = _check_productexpression(tempexp, 0) # Case 4: RHS = m.p*m.dxdt[t] if args is None: if type(tempexp.arg(1)) is EXPR.ProductExpression or \ type(tempexp.arg(1)) is EXPR.ReciprocalExpression or \ type(tempexp.arg(1)) is EXPR.DivisionExpression: args = _check_productexpression(tempexp, 1) # Case 5: m.dxdt[t] + sum(ELSE) = RHS # or CONSTANT + m.dxdt[t] = RHS if args is None: if type(tempexp.arg(0)) is EXPR.SumExpression: args = _check_viewsumexpression(tempexp, 0) # Case 6: RHS = m.dxdt[t] + sum(ELSE) if args is None: if type(tempexp.arg(1)) is EXPR.SumExpression: args = _check_viewsumexpression(tempexp, 1) # Case 7: RHS = m.p*m.dxdt[t] + CONSTANT # This case will be caught by Case 6 if p is immutable. If # p is mutable then this case will not be detected as a # separable differential equation # Case 8: - dxdt[t] = RHS if args is None: if type(tempexp.arg(0)) is EXPR.NegationExpression: args = _check_negationexpression(tempexp, 0) # Case 9: RHS = - dxdt[t] if args is None: if type(tempexp.arg(1)) is EXPR.NegationExpression: args = _check_negationexpression(tempexp, 1) # At this point if args is not None then args[0] contains # the _GetItemExpression for the DerivativeVar and args[1] # contains the RHS expression. If args is None then the # constraint is considered an algebraic equation if args is None: # Constraint is an algebraic equation or unsupported # differential equation if self._intpackage == 'scipy': raise DAE_Error( "Model contains an algebraic equation or " "unrecognized differential equation. Constraint " "'%s' cannot be simulated using Scipy. If you are " "trying to simulate a DAE model you must use " "CasADi as the integration package." % str(con.name)) tempexp = tempexp.arg(0) - tempexp.arg(1) algexp = substitute_pyomo2casadi(tempexp, templatemap) alglist.append(algexp) continue # Add the differential equation to rhsdict and derivlist dv = args[0] RHS = args[1] dvkey = _GetItemIndexer(dv) if dvkey in rhsdict.keys(): raise DAE_Error( "Found multiple RHS expressions for the " "DerivativeVar %s" % str(dvkey)) derivlist.append(dvkey) if self._intpackage == 'casadi': rhsdict[dvkey] = substitute_pyomo2casadi(RHS, templatemap) else: rhsdict[dvkey] = convert_pyomo2scipy(RHS, templatemap) # Check to see if we found a RHS for every DerivativeVar in # the model # FIXME: Not sure how to rework this for multi-index case # allderivs = derivs.keys() # if set(allderivs) != set(derivlist): # missing = list(set(allderivs)-set(derivlist)) # print("WARNING: Could not find a RHS expression for the " # "following DerivativeVar components "+str(missing)) # Create ordered list of differential variables corresponding # to the list of derivatives. diffvars = [] for deriv in derivlist: sv = deriv.base.get_state_var() diffvars.append(_GetItemIndexer(sv[deriv._args])) # Create ordered list of algebraic variables and time-varying # parameters algvars = [] for item in iterkeys(templatemap): if item.base.name in derivs: # Make sure there are no DerivativeVars in the # template map raise DAE_Error( "Cannot simulate a differential equation with " "multiple DerivativeVars") if item not in diffvars: # Finds time varying parameters and algebraic vars algvars.append(item) if self._intpackage == 'scipy': # Function sent to scipy integrator def _rhsfun(t, x): residual = [] cstemplate.set_value(t) for idx, v in enumerate(diffvars): if v in templatemap: templatemap[v].set_value(x[idx]) for d in derivlist: residual.append(rhsdict[d]()) return residual self._rhsfun = _rhsfun # Add any diffvars not added by expression walker to self._templatemap if self._intpackage == 'casadi': for _id in diffvars: if _id not in templatemap: name = "%s[%s]" % ( _id.base.name, ','.join(str(x) for x in _id.args)) templatemap[_id] = casadi.SX.sym(name) self._contset = contset self._cstemplate = cstemplate self._diffvars = diffvars self._derivlist = derivlist self._templatemap = templatemap self._rhsdict = rhsdict self._alglist = alglist self._algvars = algvars self._model = m self._tsim = None self._simsolution = None # The algebraic vars in the most recent simulation self._simalgvars = None # The time-varying inputs in the most recent simulation self._siminputvars = None
def __init__(self, *args, **kwds): if "wrt" in kwds and "withrespectto" in kwds: raise TypeError( "Cannot specify both 'wrt' and 'withrespectto keywords") wrt = kwds.pop('wrt', None) wrt = kwds.pop('withrespectto', wrt) if wrt is None: # Check to be sure Integral is indexed by single # ContinuousSet and take Integral with respect to that # ContinuousSet if len(args) != 1: raise ValueError( "Integral indexed by multiple ContinuousSets. " "The desired ContinuousSet must be specified using the " "keyword argument 'wrt'") wrt = args[0] if type(wrt) is not ContinuousSet: raise ValueError( "Cannot take the integral with respect to '%s'. Must take an " "integral with respect to a ContinuousSet" % wrt) self._wrt = wrt loc = None for i, s in enumerate(args): if s is wrt: loc = i # Check that the wrt ContinuousSet is in the argument list if loc is None: raise ValueError( "The ContinuousSet '%s' was not found in the indexing sets " "of the Integral" % wrt.name) self.loc = loc # Remove the index that the integral is being expanded over arg = args[0:loc] + args[loc + 1:] # Check that if bounds are given bounds = kwds.pop('bounds', None) if bounds is not None: raise DAE_Error( "Setting bounds on integrals has not yet been implemented. " "Integrals may only be taken over an entire ContinuousSet") # Create integral expression and pass to the expression initialization intexp = kwds.pop('expr', None) intexp = kwds.pop('rule', intexp) if intexp is None: raise ValueError( "Must specify an integral expression") def _trap_rule(m, *a): ds = sorted(m.find_component(wrt.local_name)) return sum(0.5 * (ds[i + 1] - ds[i]) * (intexp(m, * (a[0:loc] + (ds[i + 1],) + a[loc:])) + intexp(m, * (a[0:loc] + (ds[i],) + a[loc:]))) for i in range(len(ds) - 1)) kwds['rule'] = _trap_rule kwds.setdefault('ctype', Integral) Expression.__init__(self, *arg, **kwds)
def _transformBlock(self, block, currentds): self._fe = {} for ds in block.component_objects(ContinuousSet): if currentds is None or currentds == ds.name or currentds is ds: generate_finite_elements(ds, self._nfe[currentds]) if not ds.get_changed(): if len(ds) - 1 > self._nfe[currentds]: logger.warn("More finite elements were found in " "ContinuousSet '%s' than the number of " "finite elements specified in apply. The " "larger number of finite elements will be " "used." % ds.name) self._nfe[ds.name] = len(ds) - 1 self._fe[ds.name] = sorted(ds) # Adding discretization information to the ContinuousSet # object itself so that it can be accessed outside of the # discretization object disc_info = ds.get_discretization_info() disc_info['nfe'] = self._nfe[ds.name] disc_info['scheme'] = self._scheme_name + ' Difference' # Maybe check to see if any of the ContinuousSets have been changed, # if they haven't then the model components need not be updated # or even iterated through expand_components(block) for d in block.component_objects(DerivativeVar, descend_into=True): dsets = d.get_continuousset_list() for i in set(dsets): if currentds is None or i.name == currentds: oldexpr = d.get_derivative_expression() loc = d.get_state_var()._contset[i] count = dsets.count(i) if count >= 3: raise DAE_Error( "Error discretizing '%s' with respect to '%s'. " "Current implementation only allows for taking the" " first or second derivative with respect to " "a particular ContinuousSet" % (d.name, i.name)) scheme = self._scheme[count - 1] newexpr = create_partial_expression( scheme, oldexpr, i, loc) d.set_derivative_expression(newexpr) # Reclassify DerivativeVar if all indexing ContinuousSets have # been discretized. Add discretization equations to the # DerivativeVar's parent block. if d.is_fully_discretized(): add_discretization_equations(d.parent_block(), d) d.parent_block().reclassify_component_type(d, Var) # Keep track of any reclassified DerivativeVar components so # that the Simulator can easily identify them if the model # is simulated after discretization # TODO: Update the discretization transformations to use # a Block to add things to the model and store discretization # information. Using a list for now because the simulator # does not yet support models containing active Blocks reclassified_list = getattr( block, '_pyomo_dae_reclassified_derivativevars', None) if reclassified_list is None: block._pyomo_dae_reclassified_derivativevars = list() reclassified_list = \ block._pyomo_dae_reclassified_derivativevars reclassified_list.append(d) # Reclassify Integrals if all ContinuousSets have been discretized if block_fully_discretized(block): if block.contains_component(Integral): for i in block.component_objects(Integral, descend_into=True): i.reconstruct() i.parent_block().reclassify_component_type(i, Expression) # If a model contains integrals they are most likely to # appear in the objective function which will need to be # reconstructed after the model is discretized. for k in block.component_objects(Objective, descend_into=True): # TODO: check this, reconstruct might not work k.reconstruct()