1*55a74a43SLisandro Dalcin""" 2*55a74a43SLisandro DalcinThis example demonstrates the use of TAO for Python to solve an 3*55a74a43SLisandro Dalcinunconstrained minimization problem on a single processor. 4*55a74a43SLisandro Dalcin 5*55a74a43SLisandro DalcinWe minimize the extended Rosenbrock function:: 6*55a74a43SLisandro Dalcin 7*55a74a43SLisandro Dalcin sum_{i=0}^{n/2-1} ( alpha*(x_{2i+1}-x_{2i}^2)^2 + (1-x_{2i})^2 ) 8*55a74a43SLisandro Dalcin""" 9*55a74a43SLisandro Dalcin 10*55a74a43SLisandro Dalcintry: range = xrange 11*55a74a43SLisandro Dalcinexcept NameError: pass 12*55a74a43SLisandro Dalcin 13*55a74a43SLisandro Dalcin# the two lines below are only 14*55a74a43SLisandro Dalcin# needed to build options database 15*55a74a43SLisandro Dalcin# from command line arguments 16*55a74a43SLisandro Dalcinimport sys, petsc4py 17*55a74a43SLisandro Dalcinpetsc4py.init(sys.argv) 18*55a74a43SLisandro Dalcin 19*55a74a43SLisandro Dalcinfrom petsc4py import PETSc 20*55a74a43SLisandro Dalcin 21*55a74a43SLisandro Dalcin 22*55a74a43SLisandro Dalcinclass AppCtx(object): 23*55a74a43SLisandro Dalcin 24*55a74a43SLisandro Dalcin """ 25*55a74a43SLisandro Dalcin Extended Rosenbrock function. 26*55a74a43SLisandro Dalcin """ 27*55a74a43SLisandro Dalcin 28*55a74a43SLisandro Dalcin def __init__(self, n=2, alpha=99.0): 29*55a74a43SLisandro Dalcin self.size = int(n) 30*55a74a43SLisandro Dalcin self.alpha = float(alpha) 31*55a74a43SLisandro Dalcin 32*55a74a43SLisandro Dalcin def formObjective(self, tao, x): 33*55a74a43SLisandro Dalcin #print 'AppCtx.formObjective()' 34*55a74a43SLisandro Dalcin alpha = self.alpha 35*55a74a43SLisandro Dalcin nn = self.size // 2 36*55a74a43SLisandro Dalcin ff = 0.0 37*55a74a43SLisandro Dalcin for i in range(nn): 38*55a74a43SLisandro Dalcin t1 = x[2*i+1] - x[2*i] * x[2*i] 39*55a74a43SLisandro Dalcin t2 = 1 - x[2*i]; 40*55a74a43SLisandro Dalcin ff += alpha*t1*t1 + t2*t2; 41*55a74a43SLisandro Dalcin return ff 42*55a74a43SLisandro Dalcin 43*55a74a43SLisandro Dalcin def formGradient(self, tao, x, G): 44*55a74a43SLisandro Dalcin #print 'AppCtx.formGradient()' 45*55a74a43SLisandro Dalcin alpha = self.alpha 46*55a74a43SLisandro Dalcin nn = self.size // 2 47*55a74a43SLisandro Dalcin G.zeroEntries() 48*55a74a43SLisandro Dalcin for i in range(nn): 49*55a74a43SLisandro Dalcin t1 = x[2*i+1] - x[2*i] * x[2*i] 50*55a74a43SLisandro Dalcin t2 = 1 - x[2*i]; 51*55a74a43SLisandro Dalcin G[2*i] = -4*alpha*t1*x[2*i] - 2*t2; 52*55a74a43SLisandro Dalcin G[2*i+1] = 2*alpha*t1; 53*55a74a43SLisandro Dalcin 54*55a74a43SLisandro Dalcin def formObjGrad(self, tao, x, G): 55*55a74a43SLisandro Dalcin #print 'AppCtx.formObjGrad()' 56*55a74a43SLisandro Dalcin alpha = self.alpha 57*55a74a43SLisandro Dalcin nn = self.size // 2 58*55a74a43SLisandro Dalcin ff = 0.0 59*55a74a43SLisandro Dalcin G.zeroEntries() 60*55a74a43SLisandro Dalcin for i in range(nn): 61*55a74a43SLisandro Dalcin t1 = x[2*i+1] - x[2*i] * x[2*i] 62*55a74a43SLisandro Dalcin t2 = 1 - x[2*i]; 63*55a74a43SLisandro Dalcin ff += alpha*t1*t1 + t2*t2; 64*55a74a43SLisandro Dalcin G[2*i] = -4*alpha*t1*x[2*i] - 2*t2; 65*55a74a43SLisandro Dalcin G[2*i+1] = 2*alpha*t1; 66*55a74a43SLisandro Dalcin return ff 67*55a74a43SLisandro Dalcin 68*55a74a43SLisandro Dalcin def formHessian(self, tao, x, H, HP): 69*55a74a43SLisandro Dalcin #print 'AppCtx.formHessian()' 70*55a74a43SLisandro Dalcin alpha = self.alpha 71*55a74a43SLisandro Dalcin nn = self.size // 2 72*55a74a43SLisandro Dalcin idx = [0, 0] 73*55a74a43SLisandro Dalcin v = [[0.0, 0.0], 74*55a74a43SLisandro Dalcin [0.0, 0.0]] 75*55a74a43SLisandro Dalcin H.zeroEntries() 76*55a74a43SLisandro Dalcin for i in range(nn): 77*55a74a43SLisandro Dalcin v[1][1] = 2*alpha 78*55a74a43SLisandro Dalcin v[0][0] = -4*alpha*(x[2*i+1]-3*x[2*i]*x[2*i]) + 2 79*55a74a43SLisandro Dalcin v[1][0] = v[0][1] = -4.0*alpha*x[2*i]; 80*55a74a43SLisandro Dalcin idx[0] = 2*i 81*55a74a43SLisandro Dalcin idx[1] = 2*i+1 82*55a74a43SLisandro Dalcin H[idx,idx] = v 83*55a74a43SLisandro Dalcin H.assemble() 84*55a74a43SLisandro Dalcin 85*55a74a43SLisandro Dalcin# access PETSc options database 86*55a74a43SLisandro DalcinOptDB = PETSc.Options() 87*55a74a43SLisandro Dalcin 88*55a74a43SLisandro Dalcin 89*55a74a43SLisandro Dalcin# create user application context 90*55a74a43SLisandro Dalcin# and configure user parameters 91*55a74a43SLisandro Dalcinuser = AppCtx() 92*55a74a43SLisandro Dalcinuser.size = OptDB.getInt ( 'n', user.size) 93*55a74a43SLisandro Dalcinuser.alpha = OptDB.getReal('alpha', user.alpha) 94*55a74a43SLisandro Dalcin 95*55a74a43SLisandro Dalcin# create solution vector 96*55a74a43SLisandro Dalcinx = PETSc.Vec().create(PETSc.COMM_SELF) 97*55a74a43SLisandro Dalcinx.setSizes(user.size) 98*55a74a43SLisandro Dalcinx.setFromOptions() 99*55a74a43SLisandro Dalcin 100*55a74a43SLisandro Dalcin# create Hessian matrix 101*55a74a43SLisandro DalcinH = PETSc.Mat().create(PETSc.COMM_SELF) 102*55a74a43SLisandro DalcinH.setSizes([user.size, user.size]) 103*55a74a43SLisandro DalcinH.setFromOptions() 104*55a74a43SLisandro DalcinH.setOption(PETSc.Mat.Option.SYMMETRIC, True) 105*55a74a43SLisandro DalcinH.setUp() 106*55a74a43SLisandro Dalcin 107*55a74a43SLisandro Dalcin# pass the following to command line: 108*55a74a43SLisandro Dalcin# $ ... -methods nm,lmvm,nls,ntr,cg,blmvm,tron 109*55a74a43SLisandro Dalcin# to try many methods 110*55a74a43SLisandro Dalcinmethods = OptDB.getString('methods', '') 111*55a74a43SLisandro Dalcinmethods = methods.split(',') 112*55a74a43SLisandro Dalcinfor meth in methods: 113*55a74a43SLisandro Dalcin # create TAO Solver 114*55a74a43SLisandro Dalcin tao = PETSc.TAO().create(PETSc.COMM_SELF) 115*55a74a43SLisandro Dalcin if meth: tao.setType(meth) 116*55a74a43SLisandro Dalcin tao.setFromOptions() 117*55a74a43SLisandro Dalcin # solve the problem 118*55a74a43SLisandro Dalcin tao.setObjectiveGradient(user.formObjGrad) 119*55a74a43SLisandro Dalcin tao.setObjective(user.formObjective) 120*55a74a43SLisandro Dalcin tao.setGradient(user.formGradient) 121*55a74a43SLisandro Dalcin tao.setHessian(user.formHessian, H) 122*55a74a43SLisandro Dalcin #app.getKSP().getPC().setFromOptions() 123*55a74a43SLisandro Dalcin x.set(0) # zero initial guess 124*55a74a43SLisandro Dalcin tao.setSolution(x) 125*55a74a43SLisandro Dalcin tao.solve() 126*55a74a43SLisandro Dalcin tao.destroy() 127*55a74a43SLisandro Dalcin 128*55a74a43SLisandro Dalcin## # this is just for testing 129*55a74a43SLisandro Dalcin## x = app.getSolution() 130*55a74a43SLisandro Dalcin## G = app.getGradient() 131*55a74a43SLisandro Dalcin## H, HP = app.getHessian() 132*55a74a43SLisandro Dalcin 133*55a74a43SLisandro Dalcin## f = tao.computeObjective(x) 134*55a74a43SLisandro Dalcin## tao.computeGradient(x, G) 135*55a74a43SLisandro Dalcin## f = tao.computeObjectiveGradient(x, G) 136*55a74a43SLisandro Dalcin## tao.computeHessian(x, H, HP) 137