1aaa7dc30SBarry Smith #include <../src/tao/complementarity/impls/ssls/ssls.h> 2a7e14dcfSSatish Balay /* 3a7e14dcfSSatish Balay Context for ASXLS 4a7e14dcfSSatish Balay -- active-set - reduced matrices formed 5a7e14dcfSSatish Balay - inherit properties of original system 6a7e14dcfSSatish Balay -- semismooth (S) - function not differentiable 7a7e14dcfSSatish Balay - merit function continuously differentiable 8a7e14dcfSSatish Balay - Fischer-Burmeister reformulation of complementarity 9a7e14dcfSSatish Balay - Billups composition for two finite bounds 10a7e14dcfSSatish Balay -- infeasible (I) - iterates not guaranteed to remain within bounds 11a7e14dcfSSatish Balay -- feasible (F) - iterates guaranteed to remain within bounds 12a7e14dcfSSatish Balay -- linesearch (LS) - Armijo rule on direction 13a7e14dcfSSatish Balay 14a7e14dcfSSatish Balay Many other reformulations are possible and combinations of 15a7e14dcfSSatish Balay feasible/infeasible and linesearch/trust region are possible. 16a7e14dcfSSatish Balay 17a7e14dcfSSatish Balay Basic theory 18a7e14dcfSSatish Balay Fischer-Burmeister reformulation is semismooth with a continuously 19a7e14dcfSSatish Balay differentiable merit function and strongly semismooth if the F has 20a7e14dcfSSatish Balay lipschitz continuous derivatives. 21a7e14dcfSSatish Balay 22a7e14dcfSSatish Balay Every accumulation point generated by the algorithm is a stationary 23a7e14dcfSSatish Balay point for the merit function. Stationary points of the merit function 24a7e14dcfSSatish Balay are solutions of the complementarity problem if 25a7e14dcfSSatish Balay a. the stationary point has a BD-regular subdifferential, or 26a7e14dcfSSatish Balay b. the Schur complement F'/F'_ff is a P_0-matrix where ff is the 27a7e14dcfSSatish Balay index set corresponding to the free variables. 28a7e14dcfSSatish Balay 29a7e14dcfSSatish Balay If one of the accumulation points has a BD-regular subdifferential then 30a7e14dcfSSatish Balay a. the entire sequence converges to this accumulation point at 31a7e14dcfSSatish Balay a local q-superlinear rate 32a7e14dcfSSatish Balay b. if in addition the reformulation is strongly semismooth near 33a7e14dcfSSatish Balay this accumulation point, then the algorithm converges at a 34a7e14dcfSSatish Balay local q-quadratic rate. 35a7e14dcfSSatish Balay 36a7e14dcfSSatish Balay The theory for the feasible version follows from the feasible descent 37a7e14dcfSSatish Balay algorithm framework. 38a7e14dcfSSatish Balay 39a7e14dcfSSatish Balay References: 40a7e14dcfSSatish Balay Billups, "Algorithms for Complementarity Problems and Generalized 4196a0c994SBarry Smith Equations," Ph.D thesis, University of Wisconsin Madison, 1995. 42a7e14dcfSSatish Balay De Luca, Facchinei, Kanzow, "A Semismooth Equation Approach to the 43a7e14dcfSSatish Balay Solution of Nonlinear Complementarity Problems," Mathematical 4496a0c994SBarry Smith Programming, 75, 1996. 45a7e14dcfSSatish Balay Ferris, Kanzow, Munson, "Feasible Descent Algorithms for Mixed 46a7e14dcfSSatish Balay Complementarity Problems," Mathematical Programming, 86, 4796a0c994SBarry Smith 1999. 4896a0c994SBarry Smith Fischer, "A Special Newton type Optimization Method," Optimization, 4996a0c994SBarry Smith 24, 1992 50a7e14dcfSSatish Balay Munson, Facchinei, Ferris, Fischer, Kanzow, "The Semismooth Algorithm 5196a0c994SBarry Smith for Large Scale Complementarity Problems," Technical Report, 5296a0c994SBarry Smith University of Wisconsin Madison, 1999. 53a7e14dcfSSatish Balay */ 54a7e14dcfSSatish Balay 55a7e14dcfSSatish Balay 56e0877f53SBarry Smith static PetscErrorCode TaoSetUp_ASILS(Tao tao) 57a7e14dcfSSatish Balay { 58a7e14dcfSSatish Balay TAO_SSLS *asls = (TAO_SSLS *)tao->data; 59a7e14dcfSSatish Balay PetscErrorCode ierr; 60a7e14dcfSSatish Balay 61a7e14dcfSSatish Balay PetscFunctionBegin; 62a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr); 63a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr); 64a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&asls->ff);CHKERRQ(ierr); 65a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&asls->dpsi);CHKERRQ(ierr); 66a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&asls->da);CHKERRQ(ierr); 67a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&asls->db);CHKERRQ(ierr); 68a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&asls->t1);CHKERRQ(ierr); 69a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution,&asls->t2);CHKERRQ(ierr); 706c23d075SBarry Smith asls->fixed = NULL; 716c23d075SBarry Smith asls->free = NULL; 726c23d075SBarry Smith asls->J_sub = NULL; 736c23d075SBarry Smith asls->Jpre_sub = NULL; 746c23d075SBarry Smith asls->w = NULL; 756c23d075SBarry Smith asls->r1 = NULL; 766c23d075SBarry Smith asls->r2 = NULL; 776c23d075SBarry Smith asls->r3 = NULL; 786c23d075SBarry Smith asls->dxfree = NULL; 79a7e14dcfSSatish Balay PetscFunctionReturn(0); 80a7e14dcfSSatish Balay } 81a7e14dcfSSatish Balay 82a7e14dcfSSatish Balay static PetscErrorCode Tao_ASLS_FunctionGradient(TaoLineSearch ls, Vec X, PetscReal *fcn, Vec G, void *ptr) 83a7e14dcfSSatish Balay { 84441846f8SBarry Smith Tao tao = (Tao)ptr; 85a7e14dcfSSatish Balay TAO_SSLS *asls = (TAO_SSLS *)tao->data; 86a7e14dcfSSatish Balay PetscErrorCode ierr; 87a7e14dcfSSatish Balay 88a7e14dcfSSatish Balay PetscFunctionBegin; 89a7e14dcfSSatish Balay ierr = TaoComputeConstraints(tao, X, tao->constraints);CHKERRQ(ierr); 90a7e14dcfSSatish Balay ierr = VecFischer(X,tao->constraints,tao->XL,tao->XU,asls->ff);CHKERRQ(ierr); 91a7e14dcfSSatish Balay ierr = VecNorm(asls->ff,NORM_2,&asls->merit);CHKERRQ(ierr); 92a7e14dcfSSatish Balay *fcn = 0.5*asls->merit*asls->merit; 93a7e14dcfSSatish Balay 94ffad9901SBarry Smith ierr = TaoComputeJacobian(tao,tao->solution,tao->jacobian,tao->jacobian_pre);CHKERRQ(ierr); 95235fd6e6SBarry Smith ierr = MatDFischer(tao->jacobian, tao->solution, tao->constraints,tao->XL, tao->XU, asls->t1, asls->t2,asls->da, asls->db);CHKERRQ(ierr); 96a7e14dcfSSatish Balay ierr = VecPointwiseMult(asls->t1, asls->ff, asls->db);CHKERRQ(ierr); 97a7e14dcfSSatish Balay ierr = MatMultTranspose(tao->jacobian,asls->t1,G);CHKERRQ(ierr); 98a7e14dcfSSatish Balay ierr = VecPointwiseMult(asls->t1, asls->ff, asls->da);CHKERRQ(ierr); 99a7e14dcfSSatish Balay ierr = VecAXPY(G,1.0,asls->t1);CHKERRQ(ierr); 100a7e14dcfSSatish Balay PetscFunctionReturn(0); 101a7e14dcfSSatish Balay } 102a7e14dcfSSatish Balay 103441846f8SBarry Smith static PetscErrorCode TaoDestroy_ASILS(Tao tao) 104a7e14dcfSSatish Balay { 105a7e14dcfSSatish Balay TAO_SSLS *ssls = (TAO_SSLS *)tao->data; 106a7e14dcfSSatish Balay PetscErrorCode ierr; 107a7e14dcfSSatish Balay 108a7e14dcfSSatish Balay PetscFunctionBegin; 109a7e14dcfSSatish Balay ierr = VecDestroy(&ssls->ff);CHKERRQ(ierr); 110a7e14dcfSSatish Balay ierr = VecDestroy(&ssls->dpsi);CHKERRQ(ierr); 111a7e14dcfSSatish Balay ierr = VecDestroy(&ssls->da);CHKERRQ(ierr); 112a7e14dcfSSatish Balay ierr = VecDestroy(&ssls->db);CHKERRQ(ierr); 113a7e14dcfSSatish Balay ierr = VecDestroy(&ssls->w);CHKERRQ(ierr); 114a7e14dcfSSatish Balay ierr = VecDestroy(&ssls->t1);CHKERRQ(ierr); 115a7e14dcfSSatish Balay ierr = VecDestroy(&ssls->t2);CHKERRQ(ierr); 116a7e14dcfSSatish Balay ierr = VecDestroy(&ssls->r1);CHKERRQ(ierr); 117a7e14dcfSSatish Balay ierr = VecDestroy(&ssls->r2);CHKERRQ(ierr); 118a7e14dcfSSatish Balay ierr = VecDestroy(&ssls->r3);CHKERRQ(ierr); 119a7e14dcfSSatish Balay ierr = VecDestroy(&ssls->dxfree);CHKERRQ(ierr); 120a7e14dcfSSatish Balay ierr = MatDestroy(&ssls->J_sub);CHKERRQ(ierr); 121a7e14dcfSSatish Balay ierr = MatDestroy(&ssls->Jpre_sub);CHKERRQ(ierr); 122a7e14dcfSSatish Balay ierr = ISDestroy(&ssls->fixed);CHKERRQ(ierr); 123a7e14dcfSSatish Balay ierr = ISDestroy(&ssls->free);CHKERRQ(ierr); 124a7e14dcfSSatish Balay ierr = PetscFree(tao->data);CHKERRQ(ierr); 125a7e14dcfSSatish Balay PetscFunctionReturn(0); 126a7e14dcfSSatish Balay } 12747a47007SBarry Smith 128441846f8SBarry Smith static PetscErrorCode TaoSolve_ASILS(Tao tao) 129a7e14dcfSSatish Balay { 130a7e14dcfSSatish Balay TAO_SSLS *asls = (TAO_SSLS *)tao->data; 131a7e14dcfSSatish Balay PetscReal psi,ndpsi, normd, innerd, t=0; 1328931d482SJason Sarich PetscInt nf; 133a7e14dcfSSatish Balay PetscErrorCode ierr; 134e4cb33bbSBarry Smith TaoConvergedReason reason; 135e4cb33bbSBarry Smith TaoLineSearchConvergedReason ls_reason; 136a7e14dcfSSatish Balay 137a7e14dcfSSatish Balay PetscFunctionBegin; 138a7e14dcfSSatish Balay /* Assume that Setup has been called! 139a7e14dcfSSatish Balay Set the structure for the Jacobian and create a linear solver. */ 140a7e14dcfSSatish Balay 141a7e14dcfSSatish Balay ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr); 142a7e14dcfSSatish Balay ierr = TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch,Tao_ASLS_FunctionGradient,tao);CHKERRQ(ierr); 143a7e14dcfSSatish Balay ierr = TaoLineSearchSetObjectiveRoutine(tao->linesearch,Tao_SSLS_Function,tao);CHKERRQ(ierr); 144a7e14dcfSSatish Balay 145a7e14dcfSSatish Balay /* Calculate the function value and fischer function value at the 146a7e14dcfSSatish Balay current iterate */ 147a7e14dcfSSatish Balay ierr = TaoLineSearchComputeObjectiveAndGradient(tao->linesearch,tao->solution,&psi,asls->dpsi);CHKERRQ(ierr); 148a7e14dcfSSatish Balay ierr = VecNorm(asls->dpsi,NORM_2,&ndpsi);CHKERRQ(ierr); 149a7e14dcfSSatish Balay 150a7e14dcfSSatish Balay while (1) { 151a7e14dcfSSatish Balay /* Check the termination criteria */ 1528931d482SJason Sarich ierr = PetscInfo3(tao,"iter %D, merit: %g, ||dpsi||: %g\n",tao->niter, (double)asls->merit, (double)ndpsi);CHKERRQ(ierr); 1538931d482SJason Sarich ierr = TaoMonitor(tao, tao->niter, asls->merit, ndpsi, 0.0, t, &reason);CHKERRQ(ierr); 154a7e14dcfSSatish Balay if (TAO_CONTINUE_ITERATING != reason) break; 155e6d4cb7fSJason Sarich tao->niter++; 156a7e14dcfSSatish Balay 157a7e14dcfSSatish Balay /* We are going to solve a linear system of equations. We need to 158a7e14dcfSSatish Balay set the tolerances for the solve so that we maintain an asymptotic 159a7e14dcfSSatish Balay rate of convergence that is superlinear. 160a7e14dcfSSatish Balay Note: these tolerances are for the reduced system. We really need 161a7e14dcfSSatish Balay to make sure that the full system satisfies the full-space conditions. 162a7e14dcfSSatish Balay 163a7e14dcfSSatish Balay This rule gives superlinear asymptotic convergence 164a7e14dcfSSatish Balay asls->atol = min(0.5, asls->merit*sqrt(asls->merit)); 165a7e14dcfSSatish Balay asls->rtol = 0.0; 166a7e14dcfSSatish Balay 167a7e14dcfSSatish Balay This rule gives quadratic asymptotic convergence 168a7e14dcfSSatish Balay asls->atol = min(0.5, asls->merit*asls->merit); 169a7e14dcfSSatish Balay asls->rtol = 0.0; 170a7e14dcfSSatish Balay 171a7e14dcfSSatish Balay Calculate a free and fixed set of variables. The fixed set of 172a7e14dcfSSatish Balay variables are those for the d_b is approximately equal to zero. 173a7e14dcfSSatish Balay The definition of approximately changes as we approach the solution 174a7e14dcfSSatish Balay to the problem. 175a7e14dcfSSatish Balay 176a7e14dcfSSatish Balay No one rule is guaranteed to work in all cases. The following 177a7e14dcfSSatish Balay definition is based on the norm of the Jacobian matrix. If the 178a7e14dcfSSatish Balay norm is large, the tolerance becomes smaller. */ 179a7e14dcfSSatish Balay ierr = MatNorm(tao->jacobian,NORM_1,&asls->identifier);CHKERRQ(ierr); 180a7e14dcfSSatish Balay asls->identifier = PetscMin(asls->merit, 1e-2) / (1 + asls->identifier); 181a7e14dcfSSatish Balay 182a7e14dcfSSatish Balay ierr = VecSet(asls->t1,-asls->identifier);CHKERRQ(ierr); 183a7e14dcfSSatish Balay ierr = VecSet(asls->t2, asls->identifier);CHKERRQ(ierr); 184a7e14dcfSSatish Balay 185a7e14dcfSSatish Balay ierr = ISDestroy(&asls->fixed);CHKERRQ(ierr); 186a7e14dcfSSatish Balay ierr = ISDestroy(&asls->free);CHKERRQ(ierr); 187a7e14dcfSSatish Balay ierr = VecWhichBetweenOrEqual(asls->t1, asls->db, asls->t2, &asls->fixed);CHKERRQ(ierr); 1884473680cSBarry Smith ierr = ISComplementVec(asls->fixed,asls->t1, &asls->free);CHKERRQ(ierr); 189a7e14dcfSSatish Balay 190a7e14dcfSSatish Balay ierr = ISGetSize(asls->fixed,&nf);CHKERRQ(ierr); 191335036cbSBarry Smith ierr = PetscInfo1(tao,"Number of fixed variables: %D\n", nf);CHKERRQ(ierr); 192a7e14dcfSSatish Balay 193a7e14dcfSSatish Balay /* We now have our partition. Now calculate the direction in the 194a7e14dcfSSatish Balay fixed variable space. */ 195302440fdSBarry Smith ierr = TaoVecGetSubVec(asls->ff, asls->fixed, tao->subset_type, 0.0, &asls->r1);CHKERRQ(ierr); 196302440fdSBarry Smith ierr = TaoVecGetSubVec(asls->da, asls->fixed, tao->subset_type, 1.0, &asls->r2);CHKERRQ(ierr); 197a7e14dcfSSatish Balay ierr = VecPointwiseDivide(asls->r1,asls->r1,asls->r2);CHKERRQ(ierr); 198a7e14dcfSSatish Balay ierr = VecSet(tao->stepdirection,0.0);CHKERRQ(ierr); 1994473680cSBarry Smith ierr = VecISAXPY(tao->stepdirection, asls->fixed,1.0,asls->r1);CHKERRQ(ierr); 200a7e14dcfSSatish Balay 201a7e14dcfSSatish Balay /* Our direction in the Fixed Variable Set is fixed. Calculate the 202a7e14dcfSSatish Balay information needed for the step in the Free Variable Set. To 203a7e14dcfSSatish Balay do this, we need to know the diagonal perturbation and the 204a7e14dcfSSatish Balay right hand side. */ 205a7e14dcfSSatish Balay 206b98f30f2SJason Sarich ierr = TaoVecGetSubVec(asls->da, asls->free, tao->subset_type, 0.0, &asls->r1);CHKERRQ(ierr); 207b98f30f2SJason Sarich ierr = TaoVecGetSubVec(asls->ff, asls->free, tao->subset_type, 0.0, &asls->r2);CHKERRQ(ierr); 208b98f30f2SJason Sarich ierr = TaoVecGetSubVec(asls->db, asls->free, tao->subset_type, 1.0, &asls->r3);CHKERRQ(ierr); 209a7e14dcfSSatish Balay ierr = VecPointwiseDivide(asls->r1,asls->r1, asls->r3);CHKERRQ(ierr); 210a7e14dcfSSatish Balay ierr = VecPointwiseDivide(asls->r2,asls->r2, asls->r3);CHKERRQ(ierr); 211a7e14dcfSSatish Balay 212a7e14dcfSSatish Balay /* r1 is the diagonal perturbation 213a7e14dcfSSatish Balay r2 is the right hand side 214a7e14dcfSSatish Balay r3 is no longer needed 215a7e14dcfSSatish Balay 216a7e14dcfSSatish Balay Now need to modify r2 for our direction choice in the fixed 217a7e14dcfSSatish Balay variable set: calculate t1 = J*d, take the reduced vector 218a7e14dcfSSatish Balay of t1 and modify r2. */ 219a7e14dcfSSatish Balay 220a7e14dcfSSatish Balay ierr = MatMult(tao->jacobian, tao->stepdirection, asls->t1);CHKERRQ(ierr); 221b98f30f2SJason Sarich ierr = TaoVecGetSubVec(asls->t1,asls->free,tao->subset_type,0.0,&asls->r3);CHKERRQ(ierr); 222a7e14dcfSSatish Balay ierr = VecAXPY(asls->r2, -1.0, asls->r3);CHKERRQ(ierr); 223a7e14dcfSSatish Balay 224a7e14dcfSSatish Balay /* Calculate the reduced problem matrix and the direction */ 22547a47007SBarry Smith if (!asls->w && (tao->subset_type == TAO_SUBSET_MASK || tao->subset_type == TAO_SUBSET_MATRIXFREE)) { 226a7e14dcfSSatish Balay ierr = VecDuplicate(tao->solution, &asls->w);CHKERRQ(ierr); 227a7e14dcfSSatish Balay } 228b98f30f2SJason Sarich ierr = TaoMatGetSubMat(tao->jacobian, asls->free, asls->w, tao->subset_type,&asls->J_sub);CHKERRQ(ierr); 229a7e14dcfSSatish Balay if (tao->jacobian != tao->jacobian_pre) { 230b98f30f2SJason Sarich ierr = TaoMatGetSubMat(tao->jacobian_pre, asls->free, asls->w, tao->subset_type, &asls->Jpre_sub);CHKERRQ(ierr); 231a7e14dcfSSatish Balay } else { 232a7e14dcfSSatish Balay ierr = MatDestroy(&asls->Jpre_sub);CHKERRQ(ierr); 233a7e14dcfSSatish Balay asls->Jpre_sub = asls->J_sub; 234a7e14dcfSSatish Balay ierr = PetscObjectReference((PetscObject)(asls->Jpre_sub));CHKERRQ(ierr); 235a7e14dcfSSatish Balay } 236a7e14dcfSSatish Balay ierr = MatDiagonalSet(asls->J_sub, asls->r1,ADD_VALUES);CHKERRQ(ierr); 237b98f30f2SJason Sarich ierr = TaoVecGetSubVec(tao->stepdirection, asls->free, tao->subset_type, 0.0, &asls->dxfree);CHKERRQ(ierr); 238a7e14dcfSSatish Balay ierr = VecSet(asls->dxfree, 0.0);CHKERRQ(ierr); 239a7e14dcfSSatish Balay 240a7e14dcfSSatish Balay /* Calculate the reduced direction. (Really negative of Newton 241a7e14dcfSSatish Balay direction. Therefore, rest of the code uses -d.) */ 242302440fdSBarry Smith ierr = KSPReset(tao->ksp);CHKERRQ(ierr); 24323ee1639SBarry Smith ierr = KSPSetOperators(tao->ksp, asls->J_sub, asls->Jpre_sub);CHKERRQ(ierr); 244a7e14dcfSSatish Balay ierr = KSPSolve(tao->ksp, asls->r2, asls->dxfree);CHKERRQ(ierr); 245b0026674SJason Sarich ierr = KSPGetIterationNumber(tao->ksp,&tao->ksp_its);CHKERRQ(ierr); 246b0026674SJason Sarich tao->ksp_tot_its+=tao->ksp_its; 247a7e14dcfSSatish Balay 248a7e14dcfSSatish Balay /* Add the direction in the free variables back into the real direction. */ 2494473680cSBarry Smith ierr = VecISAXPY(tao->stepdirection, asls->free, 1.0,asls->dxfree);CHKERRQ(ierr); 250a7e14dcfSSatish Balay 251a7e14dcfSSatish Balay /* Check the real direction for descent and if not, use the negative 252a7e14dcfSSatish Balay gradient direction. */ 253a7e14dcfSSatish Balay ierr = VecNorm(tao->stepdirection, NORM_2, &normd);CHKERRQ(ierr); 254a7e14dcfSSatish Balay ierr = VecDot(tao->stepdirection, asls->dpsi, &innerd);CHKERRQ(ierr); 255a7e14dcfSSatish Balay 2561118d4bcSLisandro Dalcin if (innerd <= asls->delta*PetscPowReal(normd, asls->rho)) { 257335036cbSBarry Smith ierr = PetscInfo1(tao,"Gradient direction: %5.4e.\n", (double)innerd);CHKERRQ(ierr); 2588931d482SJason Sarich ierr = PetscInfo1(tao, "Iteration %D: newton direction not descent\n", tao->niter);CHKERRQ(ierr); 259a7e14dcfSSatish Balay ierr = VecCopy(asls->dpsi, tao->stepdirection);CHKERRQ(ierr); 260a7e14dcfSSatish Balay ierr = VecDot(asls->dpsi, tao->stepdirection, &innerd);CHKERRQ(ierr); 261a7e14dcfSSatish Balay } 262a7e14dcfSSatish Balay 263a7e14dcfSSatish Balay ierr = VecScale(tao->stepdirection, -1.0);CHKERRQ(ierr); 264a7e14dcfSSatish Balay innerd = -innerd; 265a7e14dcfSSatish Balay 266a7e14dcfSSatish Balay /* We now have a correct descent direction. Apply a linesearch to 267a7e14dcfSSatish Balay find the new iterate. */ 268a7e14dcfSSatish Balay ierr = TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0);CHKERRQ(ierr); 26947a47007SBarry Smith ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &psi,asls->dpsi, tao->stepdirection, &t, &ls_reason);CHKERRQ(ierr); 270a7e14dcfSSatish Balay ierr = VecNorm(asls->dpsi, NORM_2, &ndpsi);CHKERRQ(ierr); 271a7e14dcfSSatish Balay } 272a7e14dcfSSatish Balay PetscFunctionReturn(0); 273a7e14dcfSSatish Balay } 274a7e14dcfSSatish Balay 275a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 2761522df2eSJason Sarich /*MC 2771522df2eSJason Sarich TAOASILS - Active-set infeasible linesearch algorithm for solving 2781522df2eSJason Sarich complementarity constraints 2791522df2eSJason Sarich 2801522df2eSJason Sarich Options Database Keys: 2811522df2eSJason Sarich + -tao_ssls_delta - descent test fraction 2821522df2eSJason Sarich - -tao_ssls_rho - descent test power 2831522df2eSJason Sarich 2841eb8069cSJason Sarich Level: beginner 2851522df2eSJason Sarich M*/ 286728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_ASILS(Tao tao) 287a7e14dcfSSatish Balay { 288a7e14dcfSSatish Balay TAO_SSLS *asls; 289a7e14dcfSSatish Balay PetscErrorCode ierr; 2908caf6e8cSBarry Smith const char *armijo_type = TAOLINESEARCHARMIJO; 291a7e14dcfSSatish Balay 292a7e14dcfSSatish Balay PetscFunctionBegin; 2933c9e27cfSGeoffrey Irving ierr = PetscNewLog(tao,&asls);CHKERRQ(ierr); 294a7e14dcfSSatish Balay tao->data = (void*)asls; 295a7e14dcfSSatish Balay tao->ops->solve = TaoSolve_ASILS; 296a7e14dcfSSatish Balay tao->ops->setup = TaoSetUp_ASILS; 297a7e14dcfSSatish Balay tao->ops->view = TaoView_SSLS; 298a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_SSLS; 299a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_ASILS; 300a7e14dcfSSatish Balay tao->subset_type = TAO_SUBSET_SUBVEC; 301a7e14dcfSSatish Balay asls->delta = 1e-10; 302a7e14dcfSSatish Balay asls->rho = 2.1; 3036c23d075SBarry Smith asls->fixed = NULL; 3046c23d075SBarry Smith asls->free = NULL; 3056c23d075SBarry Smith asls->J_sub = NULL; 3066c23d075SBarry Smith asls->Jpre_sub = NULL; 3076c23d075SBarry Smith asls->w = NULL; 3086c23d075SBarry Smith asls->r1 = NULL; 3096c23d075SBarry Smith asls->r2 = NULL; 3106c23d075SBarry Smith asls->r3 = NULL; 3116c23d075SBarry Smith asls->t1 = NULL; 3126c23d075SBarry Smith asls->t2 = NULL; 3136c23d075SBarry Smith asls->dxfree = NULL; 314a7e14dcfSSatish Balay 315a7e14dcfSSatish Balay asls->identifier = 1e-5; 316a7e14dcfSSatish Balay 317a7e14dcfSSatish Balay ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);CHKERRQ(ierr); 318*63b15415SAlp Dener ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr); 319a7e14dcfSSatish Balay ierr = TaoLineSearchSetType(tao->linesearch, armijo_type);CHKERRQ(ierr); 3205d527766SPatrick Farrell ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr); 321a7e14dcfSSatish Balay ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr); 322a7e14dcfSSatish Balay 323a7e14dcfSSatish Balay ierr = KSPCreate(((PetscObject)tao)->comm, &tao->ksp);CHKERRQ(ierr); 324*63b15415SAlp Dener ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);CHKERRQ(ierr); 3255d527766SPatrick Farrell ierr = KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix);CHKERRQ(ierr); 326a7e14dcfSSatish Balay ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr); 3276552cf8aSJason Sarich 3286552cf8aSJason Sarich /* Override default settings (unless already changed) */ 3296552cf8aSJason Sarich if (!tao->max_it_changed) tao->max_it = 2000; 3306552cf8aSJason Sarich if (!tao->max_funcs_changed) tao->max_funcs = 4000; 3316552cf8aSJason Sarich if (!tao->gttol_changed) tao->gttol = 0; 3326552cf8aSJason Sarich if (!tao->grtol_changed) tao->grtol = 0; 3336f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 3346552cf8aSJason Sarich if (!tao->gatol_changed) tao->gatol = 1.0e-6; 3356552cf8aSJason Sarich if (!tao->fmin_changed) tao->fmin = 1.0e-4; 3366f4723b1SBarry Smith #else 3376552cf8aSJason Sarich if (!tao->gatol_changed) tao->gatol = 1.0e-16; 3386552cf8aSJason Sarich if (!tao->fmin_changed) tao->fmin = 1.0e-8; 3396f4723b1SBarry Smith #endif 340a7e14dcfSSatish Balay PetscFunctionReturn(0); 341a7e14dcfSSatish Balay } 342