def limitCycle(self): """ integrate the solution for one period, remembering each of time points along the way """ self.ts = np.linspace(0, self.y0[-1], self.intoptions['lc_res']) intlc = cs.CVodesIntegrator(self.model) intlc.setOption("abstol" , self.intoptions['lc_abstol']) intlc.setOption("reltol" , self.intoptions['lc_reltol']) intlc.setOption("max_num_steps", self.intoptions['lc_maxnumsteps']) intlc.setOption("tf" , self.y0[-1]) intsim = cs.Simulator(intlc, self.ts) intsim.init() # Input Arguments intsim.setInput(self.y0[:-1], cs.INTEGRATOR_X0) intsim.setInput(self.paramset, cs.INTEGRATOR_P) intsim.evaluate() self.sol = intsim.output().toArray() # create interpolation object self.lc = self.interp_sol(self.ts, self.sol)
def S0Calc(self): ints0 = cs.CVodesIntegrator(self.model) ints0.setOption("abstol", ABSTOL) ints0.setOption("reltol", RELTOL) ints0.setOption("max_num_steps", MAXNUMSTEPS) ints0.setOption("sensitivity_method", SENSMETHOD) ints0.setOption("tf", self.period) ints0.setOption("numeric_jacobian", True) ints0.setOption("number_of_fwd_dir", self.NP) ints0.setOption("fsens_err_con", 1) ints0.setOption("fsens_abstol", RELTOL) ints0.setOption("fsens_reltol", ABSTOL) ints0.init() ints0.setInput(self.y0[:-1], cs.INTEGRATOR_X0) ints0.setInput(self.paramset, cs.INTEGRATOR_P) p0_seed = np.eye(self.NP) for i in range(0, self.NP): ints0.setFwdSeed(p0_seed[i], cs.INTEGRATOR_P, i) ints0.evaluate(self.NP) S0 = np.zeros((self.NEQ, self.NP)) for i in range(0, self.NP): tempS0 = ints0.fwdSens(cs.INTEGRATOR_XF, i).toArray().squeeze() S0[:, i] = tempS0 self.XS = S0
def average(self): """ integrate the solution with quadrature to find the average species concentration. outputs to self.avg """ ffcn_in = self.model.inputsSX() ode = self.model.outputSX() quad = cs.vertcat([ffcn_in[cs.DAE_X], ffcn_in[cs.DAE_X]**2]) quadmodel = cs.SXFunction(ffcn_in, cs.daeOut(ode=ode, quad=quad)) qint = cs.CVodesIntegrator(quadmodel) qint.setOption("abstol" , self.intoptions['lc_abstol']) qint.setOption("reltol" , self.intoptions['lc_reltol']) qint.setOption("max_num_steps" , self.intoptions['lc_maxnumsteps']) qint.setOption("tf",self.y0[-1]) qint.init() qint.setInput(self.y0[:-1], cs.INTEGRATOR_X0) qint.setInput(self.paramset, cs.INTEGRATOR_P) qint.evaluate() quad_out = qint.output(cs.INTEGRATOR_QF).toArray().squeeze() self.avg = quad_out[:self.NEQ]/self.y0[-1] self.rms = np.sqrt(quad_out[self.NEQ:]/self.y0[-1]) self.std = np.sqrt(self.rms**2 - self.avg**2)
def integrate_segment(x_in, p, t_in): """ Returns ts,xs for integration over one finite element """ f_p = cs.CVodesIntegrator(self.model) f_p.setOption("abstol", 1e-8) # tolerance f_p.setOption("reltol", 1e-8) # tolerance f_p.setOption("steps_per_checkpoint", 1000) f_p.setOption("fsens_err_con", True) f_p.setOption("tf", self.NLPdata['TF'] / self.NLPdata['NK']) # final time f_p.init() f_p.setInput(x_in, cs.INTEGRATOR_X0) f_p.setInput(p, cs.INTEGRATOR_P) f_p.evaluate() f_p.reset() def out(t): f_p.integrate(t) return f_p.output().toArray() ts = np.linspace(0, self.NLPdata['TF'] / self.NLPdata['NK'], self.NLPdata['NK']) return (ts + t_in, np.array([out(t) for t in ts]).squeeze())
def _create_arc_integrator(self, trans_duration=3, res=100, pulse_res=20): """ Create integrator and simulator objects for later use. trans_duration is number of periods to simulate perturbed trajectory and reference, res is the resolution (per period) of the output trajectories """ # Use parameterized period so the integration length can be # controlled without re-initializing self.arcint = cs.CVodesIntegrator(self.modlT) self.arcint.setOption('abstol', self.intoptions['int_abstol']) self.arcint.setOption('reltol', self.intoptions['int_reltol']) self.arcint.setOption('tf', 1.) self.arcint.init() # # Simulate the perturbed trajectory for trans_duration. tf = 2*np.pi*trans_duration traj_res = int(res*trans_duration) self.arc_traj_ts = np.linspace(0, tf, num=traj_res, endpoint=True) self.arcsim = cs.Simulator(self.arcint, self.arc_traj_ts/tf) self.arcsim.init() self.pulsesim = cs.Simulator(self.arcint, np.linspace(0, 1., num=pulse_res, endpoint=True)) self.pulsesim.init()
def SecondOrder(self, t): integrator = cs.CVodesIntegrator(self.model) integrator.setOption("abstol", ABSTOL) integrator.setOption("reltol", RELTOL) integrator.setOption("max_num_steps", MAXNUMSTEPS) integrator.setOption("sensitivity_method", SENSMETHOD) integrator.setOption("tf", t) integrator.setOption("numeric_jacobian", False) integrator.init() integrator.setInput(self.y0[:-1], cs.INTEGRATOR_X0) integrator.setInput(self.paramset, cs.INTEGRATOR_P) int2 = integrator.jacobian(cs.INTEGRATOR_P, cs.INTEGRATOR_XF) int2.setOption("numeric_jacobian", True) int2.init() int2.setInput(self.y0[:-1], cs.INTEGRATOR_X0) int2.setInput(self.paramset, cs.INTEGRATOR_P) int3 = int2.jacobian(cs.INTEGRATOR_P, cs.INTEGRATOR_XF) int3.init() int3.setInput(self.y0[:-1], cs.INTEGRATOR_X0) int3.setInput(self.paramset, cs.INTEGRATOR_P) int3.evaluate() tempout = int3.output().toArray().squeeze() out = np.empty((self.NEQ, self.NP, self.NP)) for i in range(0, self.NEQ): start = i * self.NP end = (i + 1) * self.NP out[i] = tempout[start:end] return out, int2.fwdSens().toArray(), int2.output().toArray()
def MonodromyCalc(self): intmono = cs.CVodesIntegrator(self.model) intmono.setOption("abstol", ABSTOL) intmono.setOption("reltol", RELTOL) intmono.setOption("max_num_steps", MAXNUMSTEPS) intmono.setOption("sensitivity_method", SENSMETHOD) intmono.setOption("t0", 0) intmono.setOption("tf", self.period) intmono.setOption("numeric_jacobian", True) intmono.setOption("number_of_fwd_dir", self.NEQ) intmono.setOption("fsens_err_con", 1) intmono.setOption("fsens_abstol", RELTOL) intmono.setOption("fsens_reltol", ABSTOL) intmono.init() intmono.setInput(self.y0[:-1], cs.INTEGRATOR_X0) intmono.setInput(self.paramset, cs.INTEGRATOR_P) x0_seed = np.eye(self.NEQ) for i in range(0, self.NEQ): intmono.setFwdSeed(x0_seed[i], cs.INTEGRATOR_X0, i) intmono.evaluate(self.NEQ) Monodromy = np.zeros((self.NEQ, self.NEQ)) for i in range(0, self.NEQ): tempM = intmono.fwdSens(cs.INTEGRATOR_XF, i).toArray().squeeze() Monodromy[:, i] = tempM self.M = Monodromy
def _single_pulse_comparison_state(self, state, phase, amount, trans_duration): """ Compares a single perturbation of state to a reference trajectory """ base = self.base_class # Use parameterized period so the integration length can be # controlled without re-initializing self.arcint = cs.CVodesIntegrator(base.modlT) self.arcint.setOption('abstol', self.int_opt['int_abstol']) self.arcint.setOption('reltol', self.int_opt['int_reltol']) self.arcint.setOption('tf', 1.) self.arcint.init() # Find y0 at start of pulse tstart = phase * base.y0[-1] / (2 * np.pi) y0 = base.lc(tstart) y0[state] += amount # Simulate the perturbed trajectory for trans_duration. tf = base.y0[-1] * trans_duration ts = np.linspace(0, tf, num=int(100 * trans_duration), endpoint=True) self.arcsim = cs.Simulator(self.arcint, ts / tf) self.arcsim.init() self.arcsim.setInput(y0, cs.INTEGRATOR_X0) self.arcsim.setInput(list(base.paramset) + [tf], cs.INTEGRATOR_P) self.arcsim.evaluate() trajectory = self.arcsim.output().toArray() yend = trajectory[-1] tend = ts[-1] % base.y0[-1] + tstart def resy(t): return np.linalg.norm(yend - base.lc(t % base.y0[-1])) # Minimize resy(t) tvals = np.linspace(0, base.y0[-1], num=25) tguess = tvals[np.array([resy(t) for t in tvals]).argmin()] tmin = sc.optimize.fmin(resy, tguess, disp=0)[0] % base.y0[-1] assert resy(tmin) / base.NEQ < 1E-3, "transient not converged" if tmin > base.y0[-1] / 2: tmin += -base.y0[-1] tdiff = tmin - tend # rescale tdiff from -T/2 to T/2 tdiff = tdiff % base.y0[-1] if tdiff > base.y0[-1] / 2: tdiff += -base.y0[-1] reference = base.lc((ts + tend + tdiff) % base.y0[-1]) # from scipy.integrate import cumtrapz # amp_change = cumtrapz((trajectory - reference).T, x=ts)[:,-1] return trajectory, reference, tdiff
def findPRC(self, res=100, num_cycles=20): """ Function to calculate the phase response curve with specified resolution """ # Make sure the lc object exists if not hasattr(self, 'lc'): self.limitCycle() # Get a state that is not at a local max/min (0 should be at # max) state_ind = 1 while self.dydt(self.y0[:-1])[state_ind] < 1E-5: state_ind += 1 integrator = cs.CVodesIntegrator(self.model) integrator.setOption("abstol", self.intoptions['sensabstol']) integrator.setOption("reltol", self.intoptions['sensreltol']) integrator.setOption("max_num_steps", self.intoptions['sensmaxnumsteps']) integrator.setOption("sensitivity_method", self.intoptions['sensmethod']); integrator.setOption("t0", 0) integrator.setOption("tf", num_cycles*self.y0[-1]) integrator.setOption("numeric_jacobian", True) integrator.setOption("fsens_err_con", 1) integrator.setOption("fsens_abstol", self.intoptions['sensabstol']) integrator.setOption("fsens_reltol", self.intoptions['sensreltol']) integrator.init() seed = np.zeros(self.NEQ) seed[state_ind] = 1. integrator.setInput(self.y0[:-1], cs.INTEGRATOR_X0) integrator.setInput(self.paramset, cs.INTEGRATOR_P) integrator.setAdjSeed(seed, cs.INTEGRATOR_XF) integrator.evaluate(0, 1) adjsens = integrator.adjSens(cs.INTEGRATOR_X0).toArray().flatten() from scipy.integrate import odeint def adj_func(y, t): """ t will increase, trace limit cycle backwards through -t. y is the vector of adjoint sensitivities """ jac = self.dfdy(self.lc((-t)%self.y0[-1])) return y.dot(jac) seed = adjsens self.prc_ts = np.linspace(0, self.y0[-1], res) P = odeint(adj_func, seed, self.prc_ts)[::-1] # Adjoint matrix self.sPRC = self._t_to_phi(P/self.dydt(self.y0[:-1])[state_ind]) dfdp = np.array([self.dfdp(self.lc(t)) for t in self.prc_ts]) # Must rescale f to \hat{f}, inverse of rescaling t self.pPRC = np.array([self.sPRC[i].dot(self._phi_to_t(dfdp[i])) for i in xrange(len(self.sPRC))]) self.rel_pPRC = self.pPRC*np.array(self.paramset) # Create interpolation object for the state phase response curve self.sPRC_interp = self.interp_sol(self.prc_ts, self.sPRC)
def calcdSdp(self): self.calcExtrema() T = np.array([self.Tmax, self.Tmin]).flatten() inds = T.argsort() ints = cs.CVodesIntegrator(self.model) ints.setOption("abstol", ABSTOL) ints.setOption("reltol", RELTOL) ints.setOption("max_num_steps", MAXNUMSTEPS) ints.setOption("sensitivity_method", SENSMETHOD) ints.setOption("tf", self.period) ints.setOption("number_of_fwd_dir", self.NP) ints.setOption("fsens_err_con", 1) ints.setOption("fsens_abstol", RELTOL) ints.setOption("fsens_reltol", ABSTOL) ints.setOption("numeric_jacobian", True) ints.init() ints.setInput(self.y0[:-1], cs.INTEGRATOR_X0) ints.setInput(self.paramset, cs.INTEGRATOR_P) p0_seed = np.eye(self.NP) for i in range(0, self.NP): ints.setFwdSeed(p0_seed[i], cs.INTEGRATOR_P, i) for i in range(0, self.NP): ints.setFwdSeed(self.S0[:, i], cs.INTEGRATOR_X0, i) def outsens(t): ints.setFinalTime(t) ints.evaluate(self.NP) Sout = np.zeros((self.NEQ, self.NP)) for i in range(0, self.NP): tempSout = ints.fwdSens(np.INTEGRATOR_XF, i).toArray().squeeze() Sout[:, i] = tempSout print "%0.1f" % (100. * t / self.period), print "%" return Sout S = np.array([outsens(t) for t in T[inds]]) Splus = np.zeros([self.NEQ, self.NP]) Sminus = np.zeros([self.NEQ, self.NP]) S = S[inds.argsort()] for i in range(0, self.NEQ): Splus[i] = S[i, i, :] Sminus[i] = S[i + self.NEQ, i, :] self.Samp = Splus - Sminus self.relSamp = (self.Samp.T / (self.Ymax - self.Ymin)).T * self.paramset
def solveBVP_scipy(self, root_method='hybr'): """ Use a scipy optimize function to optimize the BVP function """ # Make sure inputs are the correct format paramset = list(self.paramset) # Here we create and initialize the integrator SXFunction self.bvpint = cs.CVodesIntegrator(self.modlT) self.bvpint.setOption('abstol',self.intoptions['bvp_abstol']) self.bvpint.setOption('reltol',self.intoptions['bvp_reltol']) self.bvpint.setOption('tf',1) self.bvpint.setOption('disable_internal_warnings', True) self.bvpint.setOption('fsens_err_con', True) self.bvpint.init() def bvp_minimize_function(x): """ Minimization objective """ # perhaps penalize in try/catch? if np.any(x < 0): return np.ones(3) self.bvpint.setInput(x[:-1], cs.INTEGRATOR_X0) self.bvpint.setInput(paramset + [x[-1]], cs.INTEGRATOR_P) self.bvpint.evaluate() out = x[:-1] - self.bvpint.output().toArray().flatten() out = out.tolist() self.modlT.setInput(x[:-1], cs.DAE_X) self.modlT.setInput(paramset + [x[-1]], 2) self.modlT.evaluate() out += self.modlT.output()[0].toArray()[0].tolist() return np.array(out) from scipy.optimize import root options = {} root_out = root(bvp_minimize_function, self.y0, tol=self.intoptions['bvp_ftol'], method=root_method, options=options) # Check solve success if not root_out.status: raise RuntimeError("bvpsolve: " + root_out.message) # Check output convergence if np.linalg.norm(root_out.qtf) > self.intoptions['bvp_ftol']*1E4: raise RuntimeError("bvpsolve: nonconvergent") # save output to self.y0 self.y0 = root_out.x
def first_order_sensitivity(self): """ Function to calculate the first order period sensitivity matricies using the direct method. See Wilkins et. al. Only calculates initial conditions and period sensitivities. Functions for amplitude sensitivitys remain in sensitivityfns """ self.check_monodromy() monodromy = self.monodromy integrator = cs.CVodesIntegrator(self.model) integrator.setOption("abstol", self.intoptions['sensabstol']) integrator.setOption("reltol", self.intoptions['sensreltol']) integrator.setOption("max_num_steps", self.intoptions['sensmaxnumsteps']) integrator.setOption("sensitivity_method", self.intoptions['sensmethod']); integrator.setOption("t0", 0) integrator.setOption("tf", self.y0[-1]) integrator.setOption("numeric_jacobian", True) integrator.setOption("fsens_err_con", 1) integrator.setOption("fsens_abstol", self.intoptions['sensabstol']) integrator.setOption("fsens_reltol", self.intoptions['sensreltol']) integrator.init() integrator.setInput(self.y0[:-1],cs.INTEGRATOR_X0) integrator.setInput(self.paramset,cs.INTEGRATOR_P) intdyfdp = integrator.jacobian(cs.INTEGRATOR_P, cs.INTEGRATOR_XF) intdyfdp.evaluate() s0 = intdyfdp.output().toArray() self.model.init() self.model.setInput(self.y0[:-1],cs.DAE_X) self.model.setInput(self.paramset,cs.DAE_P) self.model.evaluate() ydot0 = self.model.output().toArray().squeeze() LHS = np.zeros([(self.NEQ + 1), (self.NEQ + 1)]) LHS[:-1,:-1] = monodromy - np.eye(len(monodromy)) LHS[-1,:-1] = self.dfdy(self.y0[:-1])[0] LHS[:-1,-1] = ydot0 RHS = np.zeros([(self.NEQ + 1), self.NP]) RHS[:-1] = -s0 RHS[-1] = self.dfdp(self.y0[:-1])[0] unk = np.linalg.solve(LHS,RHS) self.S0 = unk[:-1] self.dTdp = unk[-1] self.reldTdp = self.dTdp*self.paramset/self.y0[-1]
def findARC_whole(self, res=100, trans=3): """ Calculate entire sARC matrix, which will be faster than calcualting for each parameter """ # Calculate necessary quantities if not hasattr(self, 'avg'): self.average() if not hasattr(self, 'sPRC'): self.findPRC(res) # Set up quadrature integrator self.sarc_int = cs.CVodesIntegrator( self._create_ARC_model(numstates=self.NEQ)) self.sarc_int.setOption("abstol", self.intoptions['sensabstol']) self.sarc_int.setOption("reltol", self.intoptions['sensreltol']) self.sarc_int.setOption("max_num_steps", self.intoptions['sensmaxnumsteps']) self.sarc_int.setOption("t0", 0) self.sarc_int.setOption("tf", trans*self.y0[-1]) self.sarc_int.setOption("numeric_jacobian", True) self.sarc_int.init() self.arc_ts = np.linspace(0, self.y0[-1], res) amp_change = [] for t in self.arc_ts: # Initialize model and sensitivity states x0 = np.zeros(self.NEQ*(self.NEQ + 1)) x0[:self.NEQ] = self.lc(t) x0[self.NEQ:] = np.eye(self.NEQ).flatten() # Add dphi/dt from seed perturbation param = np.zeros(self.NP + self.NEQ) param[:self.NP] = self.paramset param[self.NP:] = self.sPRC_interp(t) # Evaluate model self.sarc_int.setInput(x0, cs.INTEGRATOR_X0) self.sarc_int.setInput(param, cs.INTEGRATOR_P) self.sarc_int.evaluate() out = self.sarc_int.output(cs.INTEGRATOR_QF).toArray() # amp_change += [out] amp_change += [out*2*np.pi/self.y0[-1]] #[time, state_out, state_in] self.sARC = np.array(amp_change) dfdp = np.array([self.dfdp(self.lc(t)) for t in self.arc_ts]) self.pARC = np.array([self.sARC[i].dot(self._phi_to_t(dfdp[i])) for i in xrange(len(self.sARC))]) self.rel_pARC = (np.array(self.paramset) * self.pARC / np.atleast_2d(self.avg).T)
def solveBVP_casadi(self): """ Uses casadi's interface to sundials to solve the boundary value problem using a single-shooting method with automatic differen- tiation. """ # Here we create and initialize the integrator SXFunction self.bvpint = cs.CVodesIntegrator(self.modlT) self.bvpint.setOption('abstol',self.intoptions['bvp_abstol']) self.bvpint.setOption('reltol',self.intoptions['bvp_reltol']) self.bvpint.setOption('tf',1) self.bvpint.setOption('disable_internal_warnings', True) self.bvpint.setOption('fsens_err_con', True) self.bvpint.init() # Vector of unknowns [y0, T] V = cs.msym("V",self.NEQ+1) y0 = V[:-1] T = V[-1] t = cs.msym('t') param = cs.vertcat([self.paramset, T]) yf = self.bvpint.call(cs.integratorIn(x0=y0,p=param))[0] fout = self.modlT.call(cs.daeIn(t=t, x=y0,p=param))[0] obj = (yf - y0)**2 obj.append(fout[0]) F = cs.MXFunction([V],[obj]) F.init() solver = cs.KinsolSolver(F) solver.setOption('abstol',self.intoptions['bvp_ftol']) solver.setOption('ad_mode', "forward") solver.setOption('strategy','linesearch') solver.setOption('numeric_jacobian', True) solver.setOption('exact_jacobian', False) solver.setOption('pretype', 'both') solver.setOption('use_preconditioner', True) solver.setOption('numeric_hessian', True) solver.setOption('constraints', (2,)*(self.NEQ+1)) solver.setOption('verbose', False) solver.setOption('sparse', False) solver.setOption('linear_solver', 'dense') solver.init() solver.output().set(self.y0) solver.solve() self.y0 = solver.output().toArray().squeeze()
def calc_fundamental_matrices(self, res): self.limitCycle() # Some constants ABSTOL = 1e-11 RELTOL = 1e-9 MAXNUMSTEPS = 80000 SENSMETHOD = "staggered" dt = self.y0[-1] / res integrator = cs.CVodesIntegrator(self.model) integrator.setOption("abstol", ABSTOL) integrator.setOption("reltol", RELTOL) integrator.setOption("max_num_steps", MAXNUMSTEPS) integrator.setOption("sensitivity_method", SENSMETHOD) integrator.setOption("t0", 0) integrator.setOption('tf', dt) integrator.setOption("numeric_jacobian", True) integrator.setOption("fsens_err_con", 1) integrator.setOption("fsens_abstol", RELTOL) integrator.setOption("fsens_reltol", ABSTOL) integrator.init() phi_i = [] for i in xrange(res): integrator.setInput(self.lc(i * dt), cs.INTEGRATOR_X0) integratordyfdy0 = integrator.jacobian(cs.INTEGRATOR_X0, cs.INTEGRATOR_XF) integratordyfdy0.setInput(self.lc(i * dt), cs.INTEGRATOR_X0) integratordyfdy0.init() integratordyfdy0.evaluate() phi_i += [integratordyfdy0.output().toArray()] phi_i = np.array(phi_i) # # sim = cs.Simulator(integrator, ts) # sim.init() # sim.setInput(self.y0[:-1],cs.INTEGRATOR_X0) # sim.setInput(self.paramset,cs.INTEGRATOR_P) # simdyfdy0 = sim.jacobian(cs.INTEGRATOR_X0, # cs.INTEGRATOR_XF) # simdyfdy0.evaluate() # phi_i = simdyfdy0.output().toArray().reshape([len(ts), self.NEQ, # self.NEQ]) # y_i = sim.output().toArray() return phi_i
def phase_of_point(self, point, error=False, tol=1E-3): """ Finds the phase at which the distance from the point to the limit cycle is minimized. phi=0 corresponds to the definition of y0, returns the phase and the minimum distance to the limit cycle """ point = np.asarray(point) for i in xrange(100): dist = cs.ssym("dist") x = self.model.inputSX(cs.DAE_X) ode = self.model.outputSX() dist_ode = cs.sumAll(2*(x - point)*ode) cat_x = cs.vertcat([x, dist]) cat_ode = cs.vertcat([ode, dist_ode]) dist_model = cs.SXFunction( cs.daeIn(t=self.model.inputSX(cs.DAE_T), x=cat_x, p=self.model.inputSX(cs.DAE_P)), cs.daeOut(ode=cat_ode)) dist_model.setOption("name","distance model") dist_0 = ((self.y0[:-1] - point)**2).sum() cat_y0 = np.hstack([self.y0[:-1], dist_0, self.y0[-1]]) roots_class = pBase(dist_model, self.paramset, cat_y0) # roots_class.solveBVP() roots_class.roots() phase = self._t_to_phi(roots_class.Tmin[-1]) distance = roots_class.Ymin[-1] if distance < tol: return phase, distance intr = cs.CVodesIntegrator(self.model) intr.setOption("abstol", self.intoptions['transabstol']) intr.setOption("reltol", self.intoptions['transreltol']) intr.setOption("tf", self.y0[-1]) intr.setOption("max_num_steps", self.intoptions['transmaxnumsteps']) intr.setOption("disable_internal_warnings", True) intr.init() intr.setInput(point, cs.INTEGRATOR_X0) intr.setInput(self.paramset, cs.INTEGRATOR_P) intr.evaluate() point = intr.output().toArray().flatten() raise RuntimeError("Point failed to converge to limit cycle")
def create_integrator_cvodes(): # Time t = C.SX("t") # Differential states s = C.SX("s") v = C.SX("v") y = [s, v] # Control u = C.SX("u") # Parameters k = C.SX("k") b = C.SX("b") m = 1.0 # mass # Differential equation - explicit form sdot = v vdot = (u - k * s - b * v) / m rhs = [sdot, vdot] # Input of the ode rhs ffcn_in = C.ODE_NUM_IN * [[]] ffcn_in[C.ODE_T] = [t] ffcn_in[C.ODE_Y] = y ffcn_in[C.ODE_P] = [u, k, b] # ODE right hand side ffcn = C.SXFunction(ffcn_in, [rhs]) ffcn.setOption("name", "ODE right hand side") # Explicit integrator (CVODES) integrator = C.CVodesIntegrator(ffcn) # Set options integrator.setOption("exact_jacobian", True) #integrator.setOption("linear_multistep_method","bdf") # adams or bdf #integrator.setOption("nonlinear_solver_iteration","newton") # newton or functional integrator.setOption("fsens_err_con", True) integrator.setOption("abstol", 1e-9) integrator.setOption("reltol", 1e-7) integrator.setOption("steps_per_checkpoint", 1000) #integrator.setOption("fsens_all_at_once",False) return integrator
def _create_arc_integrator(self, trans_duration=3): """ Create integrator and simulator objects for later use """ # Use parameterized period so the integration length can be # controlled without re-initializing self.arcint = cs.CVodesIntegrator(self.modlT) self.arcint.setOption('abstol', self.intoptions['int_abstol']) self.arcint.setOption('reltol', self.intoptions['int_reltol']) self.arcint.setOption('tf', 1.) self.arcint.init() # # Simulate the perturbed trajectory for trans_duration. tf = self.y0[-1] * trans_duration self.arc_traj_ts = np.linspace(0, tf, num=int(100 * trans_duration), endpoint=True) self.arcsim = cs.Simulator(self.arcint, self.arc_traj_ts / tf) self.arcsim.init()
def calc_adjoint(self, res): ts = np.linspace(5 * self.y0[-1], 6 * self.y0[-1], num=res) integrator = cs.CVodesIntegrator(self.modlT) integrator.setOption("abstol", self.intoptions['sensabstol']) integrator.setOption("reltol", self.intoptions['sensreltol']) integrator.setOption("max_num_steps", self.intoptions['sensmaxnumsteps']) integrator.setOption("sensitivity_method", self.intoptions['sensmethod']) integrator.setOption("t0", 0) integrator.setOption("tf", 1) integrator.setOption("numeric_jacobian", True) integrator.setOption("fsens_err_con", 1) integrator.setOption("fsens_abstol", self.intoptions['sensabstol']) integrator.setOption("fsens_reltol", self.intoptions['sensreltol']) integrator.init() seed = [1] + [0] * (self.NEQ - 1) y0 = self.y0[:-1] integrator.setInput(y0, cs.INTEGRATOR_X0) adjoint_sensitivities = [] ys = [] dfdp = [] param = list(self.paramset) + [ts[0]] for param[-1] in ts: integrator.setInput(param, cs.INTEGRATOR_P) integrator.setAdjSeed(seed, cs.INTEGRATOR_XF) integrator.evaluate(0, 1) adjsens = integrator.adjSens(cs.INTEGRATOR_X0).toArray().flatten() adjoint_sensitivities.append(adjsens) ys.append(integrator.output().toArray().flatten()) dfdp.append(self.dfdp(ys[-1])) Q = np.array(adjoint_sensitivities) dfdp = np.array(dfdp) self.VRC = np.array([Q[i].dot(dfdp[i]) for i in xrange(len(Q))]) self.relVRC = self.VRC * np.array(self.paramset)
def _findARC_seed(self, seeds, res=100, trans=3): # Calculate necessary quantities if not hasattr(self, 'avg'): self.average() if not hasattr(self, 'sPRC'): self.findPRC(res) # Set up quadrature integrator self.sarc_int = cs.CVodesIntegrator(self._create_ARC_model()) self.sarc_int.setOption("abstol", self.intoptions['sensabstol']) self.sarc_int.setOption("reltol", self.intoptions['sensreltol']) self.sarc_int.setOption("max_num_steps", self.intoptions['sensmaxnumsteps']) self.sarc_int.setOption("t0", 0) self.sarc_int.setOption("tf", trans*self.y0[-1]) self.sarc_int.setOption("numeric_jacobian", True) self.sarc_int.init() t_arc = np.linspace(0, self.y0[-1], res) arc = np.array([self._sarc_single_time(t, seed) for t, seed in zip(t_arc, seeds)]).squeeze() return t_arc, arc
def _simulate(self, ts, y0=None, paramset=None): """ Simulate the class, outputing the solution at the times specified by ts. Optional inputs of y0 and paramsets to use ones other than those currently in the class """ if y0 is None: y0 = self.y0[:-1] if paramset is None: paramset = self.paramset int = cs.CVodesIntegrator(self.model) int.setOption("abstol" , self.intoptions['lc_abstol']) int.setOption("reltol" , self.intoptions['lc_reltol']) int.setOption("max_num_steps", self.intoptions['lc_maxnumsteps']) int.setOption("tf" , self.y0[-1]) sim = cs.Simulator(int, ts) sim.init() # Input Arguments sim.setInput(y0, cs.INTEGRATOR_X0) sim.setInput(paramset, cs.INTEGRATOR_P) sim.evaluate() return sim.output().toArray()
def intPastTrans(self,tf=500): """ integrate the solution until self.tf """ self.integrator = cs.CVodesIntegrator(self.model) self.integrator.setOption("abstol", self.intoptions['transabstol']) self.integrator.setOption("reltol", self.intoptions['transreltol']) self.integrator.setOption("tf", tf) self.integrator.setOption("max_num_steps", self.intoptions['transmaxnumsteps']) self.integrator.setOption("disable_internal_warnings", True) self.integrator.init() self.integrator.setInput((self.y0[:-1]), cs.INTEGRATOR_X0) self.integrator.setInput(self.paramset, cs.INTEGRATOR_P) self.integrator.evaluate() self.y0 = \ np.array(self.integrator.output().toArray().squeeze().tolist() + [1]) if np.abs(self.model.output().toArray()).sum() <= 1E-3*self.NEQ: raise RuntimeWarning("intPastTrans: converged to steady state")
def check_monodromy(self): """ Check the stability of the limit cycle by finding the eigenvalues of the monodromy matrix """ integrator = cs.CVodesIntegrator(self.model) integrator.setOption("abstol", self.intoptions['sensabstol']) integrator.setOption("reltol", self.intoptions['sensreltol']) integrator.setOption("max_num_steps", self.intoptions['sensmaxnumsteps']) integrator.setOption("sensitivity_method", self.intoptions['sensmethod']); integrator.setOption("t0", 0) integrator.setOption("tf", self.y0[-1]) integrator.setOption("numeric_jacobian", True) integrator.setOption("fsens_err_con", 1) integrator.setOption("fsens_abstol", self.intoptions['sensabstol']) integrator.setOption("fsens_reltol", self.intoptions['sensreltol']) integrator.init() integrator.setInput(self.y0[:-1], cs.INTEGRATOR_X0) integrator.setInput(self.paramset, cs.INTEGRATOR_P) intdyfdy0 = integrator.jacobian(cs.INTEGRATOR_X0, cs.INTEGRATOR_XF) intdyfdy0.evaluate() monodromy = intdyfdy0.output().toArray() self.monodromy = monodromy # Calculate Floquet Multipliers, check if all (besides n_0 = 1) # are inside unit circle eigs = np.linalg.eigvals(monodromy) self.floquet_multipliers = np.abs(eigs) self.floquet_multipliers.sort() idx = (np.abs(self.floquet_multipliers - 1.0)).argmin() f = self.floquet_multipliers.tolist() f.pop(idx) return np.all(np.array(f) < 1)
def Integrate(self, t, sim=True, res=500, pars=None): """ Integrate the forced oscillator from t = 0 to t. sim=True, return vector of values. """ # Get parameters over the limit cycle if pars is None: pars = np.array([self._get_paramset(tj) for tj in self.tgrid_u]).flatten() apars = pars.reshape((self.nk, self.NP)) y = self.y0 integrator = cs.CVodesIntegrator(self.model) integrator.setOption('abstol', self.int_opt['int_abstol']) integrator.setOption('reltol', self.int_opt['int_reltol']) integrator.setOption('tf', self.h) integrator.init() t_i = 0 if not sim: while t_i < t - self.h + 1E-5: element_index = int((t_i / self.h) % self.nk) integrator.setInput(y, cs.INTEGRATOR_X0) integrator.setInput(apars[element_index], cs.INTEGRATOR_P) integrator.evaluate() y = integrator.output() t_i += self.h if t - t_i > 1E-5: element_index = int((t_i / self.h) % self.nk) integrator.setOption('tf', t - t_i) integrator.init() integrator.setInput(y, cs.INTEGRATOR_X0) integrator.setInput(apars[element_index], cs.INTEGRATOR_P) integrator.evaluate() y = integrator.output() return t, y else: # Allow for specific ts requests in res try: ts = np.linspace(0, t, res) except TypeError: ts = res sol = np.zeros((0, self.NEQ)) while t_i < t - self.h + 1E-5: element_index = int((t_i / self.h) % self.nk) ts_i = ts[np.all([ts < t_i + self.h, ts >= t_i], 0)] - t_i sim_i = cs.Simulator(integrator, np.hstack([0, ts_i, self.h])) sim_i.init() sim_i.setInput(y, cs.INTEGRATOR_X0) sim_i.setInput(apars[element_index], cs.INTEGRATOR_P) sim_i.evaluate() out = sim_i.output().toArray() sol = np.vstack([sol, out[1:-1]]) y = out[-1] t_i += self.h if t - t_i > 1E-5: element_index = int((t_i / self.h) % self.nk) ts_i = ts[np.all([ts < t_i + self.h, ts >= t_i], 0)] - t_i sim_i = cs.Simulator(integrator, np.hstack([0, ts_i, t - t_i])) sim_i.init() sim_i.setInput(y, cs.INTEGRATOR_X0) sim_i.setInput(apars[element_index], cs.INTEGRATOR_P) sim_i.evaluate() out = sim_i.output().toArray() sol = np.vstack([sol, out[1:-1]]) y = out[-1] # Close enough. elif len(sol) == len(ts) - 1: sol = np.vstack([sol, y]) return ts, sol
def setup_collocation(tf, nk, deg=5): """ Function, possibly moved to utilities, which creates spacing and matricies for collocation with desired polynomials """ # Legendre collocation points legendre_points1 = [0, 0.500000] legendre_points2 = [0, 0.211325, 0.788675] legendre_points3 = [0, 0.112702, 0.500000, 0.887298] legendre_points4 = [0, 0.069432, 0.330009, 0.669991, 0.930568] legendre_points5 = [0, 0.046910, 0.230765, 0.500000, 0.769235, 0.953090] legendre_points = [ 0, legendre_points1, legendre_points2, legendre_points3, legendre_points4, legendre_points5 ] # Radau collocation points radau_points1 = [0, 1.000000] radau_points2 = [0, 0.333333, 1.000000] radau_points3 = [0, 0.155051, 0.644949, 1.000000] radau_points4 = [0, 0.088588, 0.409467, 0.787659, 1.000000] radau_points5 = [0, 0.057104, 0.276843, 0.583590, 0.860240, 1.000000] radau_points = [ 0, radau_points1, radau_points2, radau_points3, radau_points4, radau_points5 ] # Type of collocation points LEGENDRE = 0 RADAU = 1 collocation_points = [legendre_points, radau_points] # Radau collocation points cp = RADAU # Size of the finite elements h = float(tf) / nk # Coefficients of the collocation equation C = np.zeros((deg + 1, deg + 1)) # Coefficients of the continuity equation D = np.zeros(deg + 1) E = np.zeros(deg + 1) # Collocation point tau = cs.ssym("tau") # All collocation time points tau_root = collocation_points[cp][deg] T = np.zeros((nk, deg + 1)) for i in range(nk): for j in range(deg + 1): T[i][j] = h * (i + tau_root[j]) tg = np.array(tau_root) * h for k in range(nk): if k == 0: tgrid = tg else: tgrid = np.append(tgrid, tgrid[-1] + tg) tgrid_u = np.linspace(0, tf, nk) # For all collocation points: eq 10.4 or 10.17 in Biegler's book # Construct Lagrange polynomials to get the polynomial basis at the collocation point for j in range(deg + 1): L = 1 for j2 in range(deg + 1): if j2 != j: L *= (tau - tau_root[j2]) / (tau_root[j] - tau_root[j2]) lfcn = cs.SXFunction([tau], [L]) lfcn.init() lode = cs.SXFunction(cs.daeIn(x=cs.SX('x'), t=tau), cs.daeOut(ode=L)) lint = cs.CVodesIntegrator(lode) lint.setOption('tf', 1.0) lint.init() lint.setInput(0, 0) lint.evaluate() E[j] = lint.output() # Evaluate the polynomial at the final time to get the coefficients of the continuity equation lfcn.setInput(1.0) lfcn.evaluate() D[j] = lfcn.output() # Evaluate the time derivative of the polynomial at all collocation points to get the coefficients of the continuity equation for j2 in range(deg + 1): lfcn.setInput(tau_root[j2]) lfcn.setFwdSeed(1.0) lfcn.evaluate(1, 0) C[j][j2] = lfcn.fwdSens() return_dict = { 'tgrid': tgrid, 'tgrid_u': tgrid_u, 'C': C, 'D': D, 'E': E, 'h': h, 'deg': deg, 'tau': tau_root, } return return_dict
sxs = writeADT("SXFunctionOption", "sxFunctionUnsafeSetOption", "SXFunction", f) for line in sxs[0]: print line print "" for line in sxs[1]: print line print "" s = casadi.IpoptSolver(f) nlps = writeADT("NLPSolverOption", "nlpSolverUnsafeSetOption", "NLPSolver", s) for line in nlps[0]: print line for line in nlps[1]: print line s = casadi.IdasIntegrator(f) nlps = writeADT("IdasOption", "idasUnsafeSetOption", "SXFunction", s) for line in nlps[0]: print line for line in nlps[1]: print line s = casadi.CVodesIntegrator(f) nlps = writeADT("CvodesOption", "cvodesUnsafeSetOption", "SXFunction", s) for line in nlps[0]: print line for line in nlps[1]: print line
def _single_pulse_comparison(self, param, phase, amount, pulse_duration, trans_duration): """ Compares a single perturbation to a reference trajectory """ base = self.base_class try: par_ind = self.base_class.pdict[param] except KeyError: par_ind = int(param) # Use parameterized period so the integration length can be # controlled without re-initializing self.arcint = cs.CVodesIntegrator(base.modlT) self.arcint.setOption('abstol', self.int_opt['int_abstol']) self.arcint.setOption('reltol', self.int_opt['int_reltol']) self.arcint.setOption('tf', 1.) self.arcint.init() # Find y0 at start of pulse tstart = phase * base.y0[-1] / (2 * np.pi) tpulse = pulse_duration * base.y0[-1] / (2 * np.pi) y0 = base.lc(tstart) # Integrate trajectory through pulse # Find parameter set for pulse param_init = np.array(self.paramset) param_init[par_ind] *= (1 + amount) self.arcint.setInput(y0, cs.INTEGRATOR_X0) self.arcint.setInput(param_init.tolist() + [tpulse], cs.INTEGRATOR_P) self.arcint.evaluate() yf = np.array(self.arcint.output()) # Simulate the perturbed trajectory for trans_duration. tf = base.y0[-1] * trans_duration ts = np.linspace(0, tf, num=int(100 * trans_duration), endpoint=True) self.arcsim = cs.Simulator(self.arcint, ts / tf) self.arcsim.init() self.arcsim.setInput(yf, cs.INTEGRATOR_X0) self.arcsim.setInput(list(base.paramset) + [tf], cs.INTEGRATOR_P) self.arcsim.evaluate() trajectory = self.arcsim.output().toArray() yend = trajectory[-1] tend = ts[-1] % base.y0[-1] def resy(t): return np.linalg.norm(yend - base.lc(t % base.y0[-1])) # Minimize resy(t) tvals = np.linspace(0, base.y0[-1], num=25) tguess = tvals[np.array([resy(t) for t in tvals]).argmin()] tmin = sc.optimize.fmin(resy, tguess, disp=0)[0] % base.y0[-1] assert resy(tmin) / base.NEQ < 1E-3, "transient not converged" tdiff = tmin - tend reference = base.lc((ts + tdiff) % base.y0[-1]) from scipy.integrate import cumtrapz amp_change = cumtrapz((trajectory - reference).T, x=ts)[:, -1] return amp_change
ode.init() print "creating outputs function" outputNames = outputs.keys() fOutputs = C.SXFunction( [others['xVec'], C.veccat([others['uVec'], others['pVec']])], [outputs[n] for n in outputNames]) fOutputs.setOption('name', 'fOutputs') fOutputs.init() print "creating communicator" communicator = simutils.Communicator(fOutputs, outputNames) print "creating integrator" f = C.CVodesIntegrator(ode) f.setOption("reltol", 1e-5) f.setOption("abstol", 1e-7) f.setOption("t0", 0) f.setOption("tf", 0.02) f.setOption('name', 'integrator') # f.setOption("linear_solver_creator",C.CSparse) # f.setOption("linear_solver","user_defined") # f.setOption("monitor",["res"]) f.init() def advanceState(): js = sim.handleInput() # saving/loading fstButton = 13
def calcdSdp(self, res=50): """ Function to calculate the sensitivity state profile shapes to parameter perturbations, from wilkins. Might want to move this to a new class. """ ts = np.linspace(0, self.y0[-1], res, endpoint=True) integrator = cs.CVodesIntegrator(self.model) integrator.setOption("abstol", self.intoptions['sensabstol']) integrator.setOption("reltol", self.intoptions['sensreltol']) integrator.setOption("max_num_steps", self.intoptions['sensmaxnumsteps']) integrator.setOption("sensitivity_method", self.intoptions['sensmethod']); integrator.setOption("t0", 0) integrator.setOption("tf", self.y0[-1]) integrator.setOption("numeric_jacobian", True) integrator.setOption("number_of_fwd_dir", self.NP) integrator.setOption("fsens_err_con", 1) integrator.setOption("fsens_abstol", self.intoptions['sensabstol']) integrator.setOption("fsens_reltol", self.intoptions['sensreltol']) integrator.init() integrator.setInput(self.y0[:-1],cs.INTEGRATOR_X0) integrator.setInput(self.paramset,cs.INTEGRATOR_P) p0_seed = np.eye(self.NP) for i in range(0,self.NP): integrator.setFwdSeed(p0_seed[i],cs.INTEGRATOR_P,i) for i in range(0,self.NP): integrator.setFwdSeed(self.S0[:,i],cs.INTEGRATOR_X0,i) sim = cs.Simulator(integrator, ts) sim.setOption("number_of_fwd_dir", self.NP) # sim.setOption("fsens_err_con", 1) # sim.setOption("fsens_abstol", self.intoptions['sensabstol']) # sim.setOption("fsens_reltol", self.intoptions['sensreltol']) sim.init() sim.setInput(self.y0[:-1],cs.INTEGRATOR_X0) sim.setInput(self.paramset,cs.INTEGRATOR_P) p0_seed = np.eye(self.NP) for i in range(0,self.NP): sim.setFwdSeed(p0_seed[i],cs.INTEGRATOR_P,i) for i in range(0,self.NP): sim.setFwdSeed(self.S0[:,i],cs.INTEGRATOR_X0,i) sim.evaluate(nfdir=self.NP) # Raw sensitivity matrix S, calculated with initial conditions # S_0 = Z[0]. This matrix is not periodic, and will grow # unbounded with time S = np.array([sim.fwdSens(cs.INTEGRATOR_X0, i).toArray() for i in xrange(self.NP)]) S = S.swapaxes(0,1).swapaxes(1,2) # S[t,y,p] y = sim.output().toArray() # Periodic Z matrix, defined as the state sensitivity with # constant period (Wilkins 2009, page 2710) Z = np.zeros(S.shape) for i in xrange(res): Z[i] = S[i] + (ts[i]/self.y0[-1])*np.outer(self.dydt(y[i]), self.dTdp) self.Z = Z # Periodic W matrix, a projection of the Z matrix. Captures the # change in amplitude of the state trajectory without taking # into account changes in phase. W = np.zeros(S.shape) for i, (y_i, Z_i) in enumerate(zip(y,Z)): W[i] = (np.eye(self.NEQ) - np.outer(y_i, y_i)/np.linalg.norm(y_i)**2).dot(Z_i) self.W = W