def generalized_rush_larsen_solver( ode, function_name="forward_generalized_rush_larsen", delta=1e-8, params=None, ): """ Return an ODEComponent holding expressions for the generalized Rush-Larsen method Arguments --------- ode : gotran.ODE The ODE for which the jacobian expressions should be computed function_name : str The name of the function which should be generated delta : float Value to safeguard the evaluation of the Rush-Larsen step. params : dict Parameters determining how the code should be generated """ if not ode.is_finalized: error("The ODE is not finalized") return GeneralizedRushLarsen( ode, function_name=function_name, delta=delta, params=params, )
def simplified_implicit_euler_solver( ode, function_name="forward_simplified_implicit_euler", numeric_jacobian=False, params=None, ): """ Return an ODEComponent holding expressions for the simplified implicit Euler method Arguments --------- ode : gotran.ODE The ODE for which the jacobian expressions should be computed function_name : str The name of the function which should be generated numeric_jacobian : bool If True use numeric calculated diagonal jacobian. params : dict Parameters determining how the code should be generated """ if not ode.is_finalized: error("The ODE is not finalized") return SimplifiedImplicitEuler( ode, function_name=function_name, numeric_jacobian=numeric_jacobian, params=params, )
def add_states(self, *args, **kwargs): """ Add a number of states to the current ODEComponent Arguments --------- args : list of tuples A list of tuples with states and init values. Use this to set states if you need them ordered. kwargs : dict A dict with states """ states = list(args) + sorted(kwargs.items()) if len(states) == 0: error("Expected at least one state") for arg in states: if not isinstance(arg, tuple) or len(arg) != 2: error( "excpected tuple with lenght 2 with state name (str) " "and init values as the args argument.", ) state_name, init = arg # Add the states self.add_state(state_name, init)
def rhs_expressions(ode, function_name="rhs", result_name="dy", params=None): """ Return a code component with body expressions for the right hand side Arguments --------- ode : gotran.ODE The finalized ODE function_name : str The name of the function which should be generated result_name : str The name of the variable storing the rhs result params : dict Parameters determining how the code should be generated """ check_arg(ode, ODE) if not ode.is_finalized: error( "Cannot compute right hand side expressions if the ODE is " "not finalized", ) descr = f"Compute the right hand side of the {ode} ODE" return CodeComponent( "RHSComponent", ode, function_name, descr, params=params, **{result_name: ode.state_expressions}, )
def hybrid_generalized_rush_larsen_solver( ode, function_name="forward_hybrid_generalized_rush_larsen", delta=1e-8, params=None, stiff_states=None, ): """ Return an ODEComponent holding expressions for the generalized Rush-Larsen method Arguments --------- ode : gotran.ODE The ODE for which the jacobian expressions should be computed function_name : str The name of the function which should be generated delta : float Value to safeguard the evaluation of the Rush-Larsen step. params : dict Parameters determining how the code should be generated stiff_state_variables : list of str States that are stiff and should be solved with GRL1. The remaining states are solved with explicit euler """ if not ode.is_finalized: error("The ODE is not finalized") return HybridGeneralizedRushLarsen( ode, function_name=function_name, delta=delta, params=params, stiff_state_variables=stiff_states, )
def add_parameters(self, *args, **kwargs): """ Add a number of parameters to the current ODEComponent Arguments --------- args : list of tuples A list of tuples with parameters and init values. Use this to set parameters if you need them ordered. kwargs : dict A dict with parameters """ params = list(args) + sorted(kwargs.items()) if len(params) == 0: error("expected at least one parameter") for arg in params: if not isinstance(arg, tuple) or len(arg) != 2: error( "excpected tuple with lenght 2 with parameter name (str) " "and init values as the args argument.", ) parameter_name, value = arg # Add the parameters self.add_parameter(parameter_name, value)
def add_state_solution(self, state, expr, dependent=None): """ Add a solution expression for a state Arguments --------- state : gotran.State, AppliedUndef The State that is solved expr : sympy.Basic The expression that determines the state dependent : gotran.ODEObject If given the count of this expression will follow as a fractional count based on the count of the dependent object """ state = self._expect_state(state) if f"d{state.name}_dt" in self.ode_objects: error( "Cannot registered a state solution for a state " "that has a state derivative registered.", ) if f"alg_{state.name}_0" in self.ode_objects: error( "Cannot registered a state solution for a state " "that has an algebraic expression registered.", ) # Create a StateSolution in the present component obj = StateSolution(state, expr, dependent) self._register_component_object(obj)
def jacobian_expressions( ode, function_name="compute_jacobian", result_name="jac", params=None, ): """ Return an ODEComponent holding expressions for the jacobian Arguments --------- ode : gotran.ODE The ODE for which the jacobian expressions should be computed function_name : str The name of the function which should be generated result_name : str The name of the variable storing the jacobian result params : dict Parameters determining how the code should be generated """ if not ode.is_finalized: error("The ODE is not finalized") return JacobianComponent( ode, function_name=function_name, result_name=result_name, params=params, )
def add_algebraic(self, state, expr, dependent=None): """ Add an algebraic expression which relates a State with an expression which should equal to 0 Arguments --------- state : gotran.State The State which the algebraic expression should determine expr : sympy.Basic The expression that should equal 0 dependent : gotran.ODEObject If given the count of this expression will follow as a fractional count based on the count of the dependent object """ state = self._expect_state(state) if f"d{state.name}_dt" in self.ode_objects: error( "Cannot registered an algebraic expression for a state " "that has a state derivative registered.", ) if state.is_solved: error( "Cannot registered an algebraic expression for a state " "which is registered solved.", ) # Create an AlgebraicExpression in the present component obj = AlgebraicExpression(state, expr, dependent) self._register_component_object(obj, dependent)
def finalize_component(self): """ Called whenever the component should be finalized """ if self._is_finalized: return # A flag to allow adding of for example Markov model rates self._is_finalizing = True # If component is a Markov model if self.rates: self._finalize_markov_model() if not self.is_locally_complete: incomplete_states = [] for obj in self.ode_objects: if isinstance(obj, State): if obj not in self._local_state_expressions: incomplete_states.append(obj) incomplete_state_names = [s.name for s in incomplete_states] msg = f"Cannot finalize component '{self}'. " msg += f"Missing time derivatives for the following states: {', '.join(incomplete_state_names)}" error(msg) self._is_finalizing = False self._is_finalized = True
def _finalize_markov_model(self): """ Finalize a markov model """ # Derivatives states = self.states derivatives = defaultdict(lambda: sp.sympify(0.0)) rate_check = defaultdict(lambda: 0) # Build rate information and check that each rate is added in a # symetric way used_states = [0] * self.num_states for rate in self.rate_expressions: # Get the states to_state, from_state = rate.states # Add to derivatives of the two states derivatives[from_state] -= rate.sym * from_state.sym derivatives[to_state] += rate.sym * from_state.sym if isinstance(from_state, StateSolution): from_state = from_state.state if isinstance(to_state, StateSolution): to_state = to_state.state # Register rate ind_from = states.index(from_state) ind_to = states.index(to_state) ind_tuple = (min(ind_from, ind_to), max(ind_from, ind_to)) rate_check[ind_tuple] += 1 used_states[ind_from] = 1 used_states[ind_to] = 1 # Check used states if 0 in used_states: error( f"No rate registered for state {states[used_states.find(0)]}") # Check rate symetry for (ind_from, ind_to), times in list(rate_check.items()): if times != 2: error( "Only one rate between the states {0} and {1} was " "registered, expected two.".format( states[ind_from], states[ind_to], ), ) # Add derivatives for state in states: # Skip solved states if not isinstance(state, State) or state.is_solved: continue self.add_derivative(state, state.time.sym, derivatives[state])
def __setitem__(self, states, expr): """ Register rate(s) between states Arguments --------- states : tuple of two states, list of States, tuple of two lists of States If tuple of two states is given a single rate is expected. If one list is passed the rate expression should be a square matrix and the states list determines the order of the row and column of the matrix. If two lists are passed the first determines the states in the row and the second the states in the column of the Matrix expr : sympy.Basic, scalar, sympy.MatrixBase A sympy.Basic and scalars is expected for a single rate between two states. A sympy.Matrix is expected if several rate expressions are given, see explaination in for states argument for how the columns and rows are interpreted. """ if isinstance(expr, sp.Matrix): self._comp()._add_rates(states, expr) else: if not isinstance(states, tuple) or len(states) != 2: error( "Expected a tuple of size 2 with states when " "registering a single rate.", ) # NOTE: the actuall item is set by the component while calling this # function, using _register_single_rate. See below. self._comp()._add_single_rate(states[0], states[1], expr)
def parent(self): """ Return the the parent if it is still alive. Otherwise it will return None """ p = self._parent() if p is None: error("Cannot access parent ODEComponent. It is already " "destroyed.") return p
def model_arguments(**kwargs): """ Defines arguments that can be altered while the ODE is loaded Example ------- In gotran model file: >>> ... >>> model_arguments(include_Na=True) >>> if include_Na: >>> states(Na=1.0) >>> ... When the model gets loaded >>> ... >>> load_ode("model", include_Na=False) >>> ... """ # Check the passed load arguments for key in load_arguments: if key not in kwargs: error(f"Name '{key}' is not a model_argument.") # Update the namespace ns = {} for key, value in list(kwargs.items()): if not isinstance(value, (float, int, str, Param)): error( "expected only 'float', 'int', 'str' or 'Param', as model_arguments, " "got: '{}' for '{}'".format(type(value).__name__, key), ) if key not in load_arguments: ns[key] = value.getvalue() if isinstance(value, Param) else value else: # Try to cast the passed load_arguments with the orginal type if isinstance(value, Param): # Cast value new_value = value.value_type(load_arguments[key]) # Try to set new value value.setvalue(new_value) # Assign actual value of Param ns[key] = value.getvalue() else: ns[key] = type(value)(load_arguments[key]) namespace.update(ns)
def _parse_subtree(self, root, parent=None, first_operand=True): op = self._gettag(root) # If the tag i "apply" pick the operator and continue parsing if op == "apply": children = list(root) op = self._gettag(children[0]) root = children[1:] # If use special method to parse if hasattr(self, "_parse_" + op): return getattr(self, "_parse_" + op)(root, parent) elif op in list(self._operators.keys()): # Build the equation string eq = [] # Check if we need parenthesis use_parent = self.use_parenthesis(op, parent, first_operand) # If unary operator if len(root) == 1: # Special case if operator is "minus" if op == "minus": # If an unary minus is infront of a cn or ci we skip # parenthesize if self._gettag(root[0]) in ["ci", "cn"]: use_parent = False # If an unary minus is infront of a plus we always use parenthesize if (self._gettag(root[0]) == "apply" and self._gettag( list(root[0])[0], ) in ["plus", "minus"]): use_parent = True eq += ["-"] else: # Always use paranthesis for unary operators use_parent = True eq += [self._operators[op]] eq += (["("] * use_parent + self._parse_subtree(root[0], op) + [")"] * use_parent) return eq else: # Binary operator eq += ["("] * use_parent + self._parse_subtree(root[0], op) for operand in root[1:]: eq = ( eq + [self._operators[op]] + self._parse_subtree(operand, op, first_operand=False)) eq = eq + [")"] * use_parent return eq else: error("No support for parsing MathML " + op + " operator.")
def is_dae(self): """ Return True if ODE is a DAE """ if not self.is_complete: error("The ODE is not complete") return any( isinstance(expr, AlgebraicExpression) for expr in self.state_expressions)
def _check_name(self, name): """ Check the name """ assert isinstance(name, str) # Check for underscore in name if len(name) > 0 and name[0] == "_": error(f"No ODEObject names can start with an underscore: '{name}'") return name
def add_component(self, name): """ Add a sub ODEComponent """ comp = ODEComponent(name, self) if name in self.root.all_components: error(f"A component with the name '{name}' already excists.") self.children[comp.name] = comp self.root.all_components[comp.name] = comp self.root._present_component = comp.name return comp
def monitored_expressions( ode, monitored, function_name="monitored_expressions", result_name="monitored", params=None, ): """ Return a code component with body expressions to calculate monitored expressions Arguments --------- ode : gotran.ODE The finalized ODE for which the monitored expression should be computed monitored : tuple, list A tuple/list of strings containing the name of the monitored expressions function_name : str The name of the function which should be generated result_name : str The name of the variable storing the rhs result params : dict Parameters determining how the code should be generated """ check_arg(ode, ODE) if not ode.is_finalized: error( "Cannot compute right hand side expressions if the ODE is " "not finalized", ) check_arg(monitored, (tuple, list), itemtypes=str) monitored_exprs = [] for expr_str in monitored: obj = ode.present_ode_objects.get(expr_str) if not isinstance(obj, Expression): error(f"{expr_str} is not an expression in the {ode} ODE") monitored_exprs.append(obj) descr = f"Computes monitored expressions of the {ode} ODE" return CodeComponent( "MonitoredExpressions", ode, function_name, descr, params=params, **{result_name: monitored_exprs}, )
def add_derivative(self, der_expr, dep_var, expr, dependent=None): """ Add a derivative expression Arguments --------- der_expr : gotran.Expression, gotran.State, sympy.AppliedUndef The Expression or State which is differentiated dep_var : gotran.State, gotran.Time, gotran.Expression, sympy.AppliedUndef, sympy.Symbol The dependent variable expr : sympy.Basic The expression which the differetiation should be equal dependent : gotran.ODEObject If given the count of this expression will follow as a fractional count based on the count of the dependent object """ timer = Timer("Add derivatives") # noqa: F841 if isinstance(der_expr, AppliedUndef): name = sympycode(der_expr) der_expr = self.root.present_ode_objects.get(name) if der_expr is None: error(f"{name} is not registered in this ODE") if isinstance(dep_var, (AppliedUndef, sp.Symbol)): name = sympycode(dep_var) dep_var = self.root.present_ode_objects.get(name) if dep_var is None: error(f"{name} is not registered in this ODE") # Check if der_expr is a State if isinstance(der_expr, State): self._expect_state(der_expr) obj = StateDerivative(der_expr, expr, dependent) else: # Create a DerivativeExpression in the present component obj = DerivativeExpression(der_expr, dep_var, expr, dependent) self._register_component_object(obj, dependent) return obj.sym
def add_indexed_object(self, basename, indices, add_offset=False): """ Add an indexed object using a basename and the indices Arguments --------- basename : str The basename of the indexed expression indices : int, tuple of int The fixed indices identifying the expression add_offset : bool Add offset to indices """ timer = Timer("Add indexed object") # noqa: F841 indices = tuplewrap(indices) # Check that provided indices fit with the registered shape if len(self.shapes[basename]) > len(indices): error( "Shape mismatch between indices {0} and registered " "shape for {1}{2}".format(indices, basename, self.shapes[basename]), ) for dim, (index, shape_ind) in enumerate(zip(indices, self.shapes[basename])): if index >= shape_ind: error( "Indices must be smaller or equal to the shape. Mismatch " "in dim {0}: {1}>={2}".format(dim + 1, index, shape_ind), ) # Create IndexedObject obj = IndexedObject( basename, indices, self.shapes[basename], self._params.array, add_offset, ) self._register_component_object(obj) # Return the sympy version of the object return obj.sym
def explicit_euler_solver(ode, function_name="forward_explicit_euler", params=None): """ Return an ODEComponent holding expressions for the explicit Euler method Arguments --------- ode : gotran.ODE The ODE for which the jacobian expressions should be computed function_name : str The name of the function which should be generated params : dict Parameters determining how the code should be generated """ if not ode.is_finalized: error("The ODE is not finalized") return ExplicitEuler(ode, function_name=function_name, params=params)
def componentwise_derivative(ode, indices, params=None, result_name="dy"): """ Return an ODEComponent holding the expressions for the ith state derivative Arguments --------- ode : gotran.ODE The finalized ODE for which the ith derivative should be computed indices : int, list of ints The index params : dict Parameters determining how the code should be generated result_name : str The name of the result variable """ check_arg(ode, ODE) if not ode.is_finalized: error("Cannot compute component wise derivatives if ODE is " "not finalized") indices = listwrap(indices) check_arg(indices, list, itemtypes=int) check_arg(result_name, str) if len(indices) == 0: error("expected at least on index") registered = [] for index in indices: if index < 0 or index >= ode.num_full_states: error( "Expected the passed indices to be between 0 and the " "number of states in the ode, got {0}.".format(index), ) if index in registered: error(f"Index {index} appeared twice.") registered.append(index) # Get state expression exprs = [ode.state_expressions[index] for index in indices] results = {result_name: exprs} return CodeComponent( f"componentwise_derivatives_{'_'.join(expr.state.name for expr in exprs)}", ode, "", "", params=params, **results, )
def linearized_derivatives( ode, function_name="linear_derivatives", result_names=["linearized", "dy"], only_linear=True, include_rhs=False, nonlinear_last=False, params=None, ): """ Return an ODEComponent holding the linearized derivative expressions Arguments --------- ode : gotran.ODE The ODE for which derivatives should be linearized function_name : str The name of the function which should be generated result_names : str The name of the variable storing the linearized derivatives and the rhs evaluation if that is included. only_linear : bool If True, only linear terms will be linearized include_rhs : bool If True, rhs evaluation will be included in the generated code. nonlinear_last : bool If True the nonlinear expressions are added last after a comment params : dict Parameters determining how the code should be generated """ if not ode.is_finalized: error("The ODE is not finalized") return LinearizedDerivativeComponent( ode, function_name, result_names, only_linear, include_rhs, nonlinear_last, params, )
def mass_matrix(self): """ Return the mass matrix as a sympy.Matrix """ if not self.is_finalized: error("The ODE must be finalized") if not self._mass_matrix: state_exprs = self.state_expressions N = len(state_exprs) self._mass_matrix = sp.Matrix( N, N, lambda i, j: 1 if i == j and isinstance( state_exprs[i], StateDerivative) else 0, ) return self._mass_matrix
def _expect_state(self, state, allow_state_solution=False, only_local_states=False): """ Help function to check an argument which should be expected to be a state """ if allow_state_solution: allowed = (State, StateSolution) else: allowed = (State, ) if isinstance(state, AppliedUndef): name = sympycode(state) state = self.root.present_ode_objects.get(name) if state is None: error(f"{name} is not registered in this ODE") if only_local_states and not (state in self.states or (state in self.intermediates and allow_state_solution)): error(f"{name} is not registered in component {self.name}") check_arg(state, allowed, 0) if isinstance(state, State) and state.is_solved: error( "Cannot registered a state expression for a state " "which is registered solved.", ) return state
def add_solve_state(self, state, expr, dependent=None, **solve_flags): """ Add a solve state expression which tries to find a solution to a state by solving an algebraic expression Arguments --------- state : gotran.State, AppliedUndef The State that is solved expr : sympy.Basic The expression that determines the state dependent : gotran.ODEObject If given the count of this expression will follow as a fractional count based on the count of the dependent object solve_flags : dict Flags that are passed directly to sympy.solve """ from modelparameters.sympytools import symbols_from_expr state = self._expect_state(state) # Check the sympy flags if solve_flags.get("dict") or solve_flags.get("set"): error("Cannot use dict=True or set=True as sympy_flags") # Check that there are no expressions that are dependent on the state for sym in symbols_from_expr(expr, include_derivatives=True): if (state.sym is not sym) and (state.sym in sym.atoms()): error( "{0}, a sub expression of the expression, cannot depend " "on the state for which we try to solve for.".format(sym), ) # Try solve the passed expr try: solved_expr = sp.solve(expr, state.sym) except Exception: error("Could not solve the passed expression") assert isinstance(solved_expr, list) # FIXME: Add possibilities to choose solution? if len(solved_expr) != 1: error("expected only 1 solution") # Unpack the solution solved_expr = solved_expr[0] self.add_state_solution(state, solved_expr, dependent)
def _add_single_rate(self, to_state, from_state, expr): """ Add a single rate expression """ if self.state_expressions: error( "A component cannot have both state expressions (derivative " "and algebraic expressions) and rate expressions.", ) check_arg(expr, scalars + (sp.Basic, ), 2, ODEComponent._add_single_rate) expr = sp.sympify(expr) to_state = self._expect_state( to_state, allow_state_solution=True, only_local_states=True, ) from_state = self._expect_state( from_state, allow_state_solution=True, only_local_states=True, ) if to_state == from_state: error("The two states cannot be the same.") if (to_state.sym, from_state.sym) in self.rates: error( f"Rate between state {from_state} and {to_state} is already registered.", ) # FIXME: It should also not be possible to include other # states in the markov model, right? syms_expr = symbols_from_expr(expr) if to_state.sym in syms_expr or from_state.sym in syms_expr: error( "The rate expression cannot be dependent on the " "states it connects.", ) # Create a RateExpression obj = RateExpression(to_state, from_state, expr) self._register_component_object(obj) # Store the rate sym in the rate dict self.rates._register_single_rate(to_state.sym, from_state.sym, obj.sym)
def _add_rates(self, states, rate_matrix): """ Use a rate matrix to set rates between states Arguments --------- states : list of States, tuple of two lists of States If one list is passed the rates should be a square matrix and the states list determines the order of the row and column of the matrix. If two lists are passed the first determines the states in the row and the second the states in the column of the Matrix rates_matrix : sympy.MatrixBase A sympy.Matrix of the rate expressions between the states given in the states argument """ check_arg(states, (tuple, list), 0, ODEComponent._add_rates) check_arg(rate_matrix, sp.MatrixBase, 1, ODEComponent._add_rates) # If list if isinstance(states, list): states = (states, states) # else tuple elif len(states) != 2 and not all( isinstance(list_of_states, list) for list_of_states in states): error("expected a tuple of 2 lists with states as the " "states argument") # Check index arguments # for list_of_states in states: # print list_of_states, local_states # if not all(state in local_states for state in list_of_states): # error("Expected the states arguments to be States in "\ # "the Markov model") # Check that the length of the state lists corresponds with the shape of # the rate matrix if rate_matrix.shape[0] != len( states[0]) or rate_matrix.shape[1] != len(states[1], ): error("Shape of rates does not match given states") for i, state_i in enumerate(states[0]): for j, state_j in enumerate(states[1]): value = rate_matrix[i, j] # If 0 as rate if (isinstance(value, scalars) and value == 0) or (isinstance(value, sp.Basic) and value.is_zero): continue if state_i == state_j: error("Cannot have a nonzero rate value between the " "same states") # Assign the rate self._add_single_rate(state_i, state_j, value)
def __init__(self, der_expr, dep_var, expr, dependent=None): """ Create a DerivativeExpression Arguments --------- der_expr : Expression, State The Expression or State which is differentiated dep_var : State, Time, Expression The dependent variable expr : sympy.Basic The expression which the differetiation should be equal dependent : ODEObject If given the count of this DerivativeExpression will follow as a fractional count based on the count of the dependent object """ check_arg(der_expr, Expression, 0, DerivativeExpression) check_arg(dep_var, (State, Expression, Time), 1, DerivativeExpression) # Check that the der_expr is dependent on var if dep_var.sym not in der_expr.sym.args: error( "Cannot create a DerivativeExpression as {0} is not " "dependent on {1}".format(der_expr, dep_var), ) der_sym = sp.Derivative(der_expr.sym, dep_var.sym) self._der_expr = der_expr self._dep_var = dep_var super(DerivativeExpression, self).__init__(sympycode(der_sym), expr, dependent) self._sym = sp.Derivative(der_expr.sym, dep_var.sym) self._sym._assumptions["real"] = True self._sym._assumptions["commutative"] = True self._sym._assumptions["imaginary"] = False self._sym._assumptions["hermitian"] = True self._sym._assumptions["complex"] = True