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: 40606c0280SSatish Balay + * - Billups, "Algorithms for Complementarity Problems and Generalized 4196a0c994SBarry Smith Equations," Ph.D thesis, University of Wisconsin Madison, 1995. 42606c0280SSatish Balay . * - De Luca, Facchinei, Kanzow, "A Semismooth Equation Approach to the 43a7e14dcfSSatish Balay Solution of Nonlinear Complementarity Problems," Mathematical 4496a0c994SBarry Smith Programming, 75, 1996. 45606c0280SSatish Balay . * - Ferris, Kanzow, Munson, "Feasible Descent Algorithms for Mixed 46a7e14dcfSSatish Balay Complementarity Problems," Mathematical Programming, 86, 4796a0c994SBarry Smith 1999. 48606c0280SSatish Balay . * - Fischer, "A Special Newton type Optimization Method," Optimization, 4996a0c994SBarry Smith 24, 1992 50606c0280SSatish 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 55e0877f53SBarry Smith static PetscErrorCode TaoSetUp_ASILS(Tao tao) 56a7e14dcfSSatish Balay { 57a7e14dcfSSatish Balay TAO_SSLS *asls = (TAO_SSLS *)tao->data; 58a7e14dcfSSatish Balay 59a7e14dcfSSatish Balay PetscFunctionBegin; 609566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&tao->gradient)); 619566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&tao->stepdirection)); 629566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&asls->ff)); 639566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&asls->dpsi)); 649566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&asls->da)); 659566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&asls->db)); 669566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&asls->t1)); 679566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&asls->t2)); 686c23d075SBarry Smith asls->fixed = NULL; 696c23d075SBarry Smith asls->free = NULL; 706c23d075SBarry Smith asls->J_sub = NULL; 716c23d075SBarry Smith asls->Jpre_sub = NULL; 726c23d075SBarry Smith asls->w = NULL; 736c23d075SBarry Smith asls->r1 = NULL; 746c23d075SBarry Smith asls->r2 = NULL; 756c23d075SBarry Smith asls->r3 = NULL; 766c23d075SBarry Smith asls->dxfree = NULL; 77a7e14dcfSSatish Balay PetscFunctionReturn(0); 78a7e14dcfSSatish Balay } 79a7e14dcfSSatish Balay 80a7e14dcfSSatish Balay static PetscErrorCode Tao_ASLS_FunctionGradient(TaoLineSearch ls, Vec X, PetscReal *fcn, Vec G, void *ptr) 81a7e14dcfSSatish Balay { 82441846f8SBarry Smith Tao tao = (Tao)ptr; 83a7e14dcfSSatish Balay TAO_SSLS *asls = (TAO_SSLS *)tao->data; 84a7e14dcfSSatish Balay 85a7e14dcfSSatish Balay PetscFunctionBegin; 869566063dSJacob Faibussowitsch PetscCall(TaoComputeConstraints(tao, X, tao->constraints)); 879566063dSJacob Faibussowitsch PetscCall(VecFischer(X,tao->constraints,tao->XL,tao->XU,asls->ff)); 889566063dSJacob Faibussowitsch PetscCall(VecNorm(asls->ff,NORM_2,&asls->merit)); 89a7e14dcfSSatish Balay *fcn = 0.5*asls->merit*asls->merit; 90a7e14dcfSSatish Balay 919566063dSJacob Faibussowitsch PetscCall(TaoComputeJacobian(tao,tao->solution,tao->jacobian,tao->jacobian_pre)); 929566063dSJacob Faibussowitsch PetscCall(MatDFischer(tao->jacobian, tao->solution, tao->constraints,tao->XL, tao->XU, asls->t1, asls->t2,asls->da, asls->db)); 939566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(asls->t1, asls->ff, asls->db)); 949566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(tao->jacobian,asls->t1,G)); 959566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(asls->t1, asls->ff, asls->da)); 969566063dSJacob Faibussowitsch PetscCall(VecAXPY(G,1.0,asls->t1)); 97a7e14dcfSSatish Balay PetscFunctionReturn(0); 98a7e14dcfSSatish Balay } 99a7e14dcfSSatish Balay 100441846f8SBarry Smith static PetscErrorCode TaoDestroy_ASILS(Tao tao) 101a7e14dcfSSatish Balay { 102a7e14dcfSSatish Balay TAO_SSLS *ssls = (TAO_SSLS *)tao->data; 103a7e14dcfSSatish Balay 104a7e14dcfSSatish Balay PetscFunctionBegin; 1059566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->ff)); 1069566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->dpsi)); 1079566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->da)); 1089566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->db)); 1099566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->w)); 1109566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->t1)); 1119566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->t2)); 1129566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->r1)); 1139566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->r2)); 1149566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->r3)); 1159566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ssls->dxfree)); 1169566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ssls->J_sub)); 1179566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ssls->Jpre_sub)); 1189566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ssls->fixed)); 1199566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ssls->free)); 120*a958fbfcSStefano Zampini PetscCall(KSPDestroy(&tao->ksp)); 1219566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->data)); 122a7e14dcfSSatish Balay PetscFunctionReturn(0); 123a7e14dcfSSatish Balay } 12447a47007SBarry Smith 125441846f8SBarry Smith static PetscErrorCode TaoSolve_ASILS(Tao tao) 126a7e14dcfSSatish Balay { 127a7e14dcfSSatish Balay TAO_SSLS *asls = (TAO_SSLS *)tao->data; 128a7e14dcfSSatish Balay PetscReal psi,ndpsi, normd, innerd, t=0; 1298931d482SJason Sarich PetscInt nf; 130e4cb33bbSBarry Smith TaoLineSearchConvergedReason ls_reason; 131a7e14dcfSSatish Balay 132a7e14dcfSSatish Balay PetscFunctionBegin; 133a7e14dcfSSatish Balay /* Assume that Setup has been called! 134a7e14dcfSSatish Balay Set the structure for the Jacobian and create a linear solver. */ 135a7e14dcfSSatish Balay 1369566063dSJacob Faibussowitsch PetscCall(TaoComputeVariableBounds(tao)); 1379566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch,Tao_ASLS_FunctionGradient,tao)); 1389566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetObjectiveRoutine(tao->linesearch,Tao_SSLS_Function,tao)); 139a7e14dcfSSatish Balay 140a7e14dcfSSatish Balay /* Calculate the function value and fischer function value at the 141a7e14dcfSSatish Balay current iterate */ 1429566063dSJacob Faibussowitsch PetscCall(TaoLineSearchComputeObjectiveAndGradient(tao->linesearch,tao->solution,&psi,asls->dpsi)); 1439566063dSJacob Faibussowitsch PetscCall(VecNorm(asls->dpsi,NORM_2,&ndpsi)); 144a7e14dcfSSatish Balay 145763847b4SAlp Dener tao->reason = TAO_CONTINUE_ITERATING; 146a7e14dcfSSatish Balay while (1) { 147a7e14dcfSSatish Balay /* Check the termination criteria */ 14863a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(tao,"iter %" PetscInt_FMT ", merit: %g, ||dpsi||: %g\n",tao->niter, (double)asls->merit, (double)ndpsi)); 1499566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao,asls->merit,ndpsi,0.0,tao->ksp_its)); 1509566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao,tao->niter,asls->merit,ndpsi,0.0,t)); 1519566063dSJacob Faibussowitsch PetscCall((*tao->ops->convergencetest)(tao,tao->cnvP)); 152763847b4SAlp Dener if (TAO_CONTINUE_ITERATING != tao->reason) break; 153e1e80dc8SAlp Dener 154e1e80dc8SAlp Dener /* Call general purpose update function */ 155e1e80dc8SAlp Dener if (tao->ops->update) { 1569566063dSJacob Faibussowitsch PetscCall((*tao->ops->update)(tao, tao->niter, tao->user_update)); 157e1e80dc8SAlp Dener } 158e6d4cb7fSJason Sarich tao->niter++; 159a7e14dcfSSatish Balay 160a7e14dcfSSatish Balay /* We are going to solve a linear system of equations. We need to 161a7e14dcfSSatish Balay set the tolerances for the solve so that we maintain an asymptotic 162a7e14dcfSSatish Balay rate of convergence that is superlinear. 163a7e14dcfSSatish Balay Note: these tolerances are for the reduced system. We really need 164a7e14dcfSSatish Balay to make sure that the full system satisfies the full-space conditions. 165a7e14dcfSSatish Balay 166a7e14dcfSSatish Balay This rule gives superlinear asymptotic convergence 167a7e14dcfSSatish Balay asls->atol = min(0.5, asls->merit*sqrt(asls->merit)); 168a7e14dcfSSatish Balay asls->rtol = 0.0; 169a7e14dcfSSatish Balay 170a7e14dcfSSatish Balay This rule gives quadratic asymptotic convergence 171a7e14dcfSSatish Balay asls->atol = min(0.5, asls->merit*asls->merit); 172a7e14dcfSSatish Balay asls->rtol = 0.0; 173a7e14dcfSSatish Balay 174a7e14dcfSSatish Balay Calculate a free and fixed set of variables. The fixed set of 175a7e14dcfSSatish Balay variables are those for the d_b is approximately equal to zero. 176a7e14dcfSSatish Balay The definition of approximately changes as we approach the solution 177a7e14dcfSSatish Balay to the problem. 178a7e14dcfSSatish Balay 179a7e14dcfSSatish Balay No one rule is guaranteed to work in all cases. The following 180a7e14dcfSSatish Balay definition is based on the norm of the Jacobian matrix. If the 181a7e14dcfSSatish Balay norm is large, the tolerance becomes smaller. */ 1829566063dSJacob Faibussowitsch PetscCall(MatNorm(tao->jacobian,NORM_1,&asls->identifier)); 183a7e14dcfSSatish Balay asls->identifier = PetscMin(asls->merit, 1e-2) / (1 + asls->identifier); 184a7e14dcfSSatish Balay 1859566063dSJacob Faibussowitsch PetscCall(VecSet(asls->t1,-asls->identifier)); 1869566063dSJacob Faibussowitsch PetscCall(VecSet(asls->t2, asls->identifier)); 187a7e14dcfSSatish Balay 1889566063dSJacob Faibussowitsch PetscCall(ISDestroy(&asls->fixed)); 1899566063dSJacob Faibussowitsch PetscCall(ISDestroy(&asls->free)); 1909566063dSJacob Faibussowitsch PetscCall(VecWhichBetweenOrEqual(asls->t1, asls->db, asls->t2, &asls->fixed)); 1919566063dSJacob Faibussowitsch PetscCall(ISComplementVec(asls->fixed,asls->t1, &asls->free)); 192a7e14dcfSSatish Balay 1939566063dSJacob Faibussowitsch PetscCall(ISGetSize(asls->fixed,&nf)); 19463a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(tao,"Number of fixed variables: %" PetscInt_FMT "\n", nf)); 195a7e14dcfSSatish Balay 196a7e14dcfSSatish Balay /* We now have our partition. Now calculate the direction in the 197a7e14dcfSSatish Balay fixed variable space. */ 1989566063dSJacob Faibussowitsch PetscCall(TaoVecGetSubVec(asls->ff, asls->fixed, tao->subset_type, 0.0, &asls->r1)); 1999566063dSJacob Faibussowitsch PetscCall(TaoVecGetSubVec(asls->da, asls->fixed, tao->subset_type, 1.0, &asls->r2)); 2009566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(asls->r1,asls->r1,asls->r2)); 2019566063dSJacob Faibussowitsch PetscCall(VecSet(tao->stepdirection,0.0)); 2029566063dSJacob Faibussowitsch PetscCall(VecISAXPY(tao->stepdirection, asls->fixed,1.0,asls->r1)); 203a7e14dcfSSatish Balay 204a7e14dcfSSatish Balay /* Our direction in the Fixed Variable Set is fixed. Calculate the 205a7e14dcfSSatish Balay information needed for the step in the Free Variable Set. To 206a7e14dcfSSatish Balay do this, we need to know the diagonal perturbation and the 207a7e14dcfSSatish Balay right hand side. */ 208a7e14dcfSSatish Balay 2099566063dSJacob Faibussowitsch PetscCall(TaoVecGetSubVec(asls->da, asls->free, tao->subset_type, 0.0, &asls->r1)); 2109566063dSJacob Faibussowitsch PetscCall(TaoVecGetSubVec(asls->ff, asls->free, tao->subset_type, 0.0, &asls->r2)); 2119566063dSJacob Faibussowitsch PetscCall(TaoVecGetSubVec(asls->db, asls->free, tao->subset_type, 1.0, &asls->r3)); 2129566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(asls->r1,asls->r1, asls->r3)); 2139566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(asls->r2,asls->r2, asls->r3)); 214a7e14dcfSSatish Balay 215a7e14dcfSSatish Balay /* r1 is the diagonal perturbation 216a7e14dcfSSatish Balay r2 is the right hand side 217a7e14dcfSSatish Balay r3 is no longer needed 218a7e14dcfSSatish Balay 219a7e14dcfSSatish Balay Now need to modify r2 for our direction choice in the fixed 220a7e14dcfSSatish Balay variable set: calculate t1 = J*d, take the reduced vector 221a7e14dcfSSatish Balay of t1 and modify r2. */ 222a7e14dcfSSatish Balay 2239566063dSJacob Faibussowitsch PetscCall(MatMult(tao->jacobian, tao->stepdirection, asls->t1)); 2249566063dSJacob Faibussowitsch PetscCall(TaoVecGetSubVec(asls->t1,asls->free,tao->subset_type,0.0,&asls->r3)); 2259566063dSJacob Faibussowitsch PetscCall(VecAXPY(asls->r2, -1.0, asls->r3)); 226a7e14dcfSSatish Balay 227a7e14dcfSSatish Balay /* Calculate the reduced problem matrix and the direction */ 22847a47007SBarry Smith if (!asls->w && (tao->subset_type == TAO_SUBSET_MASK || tao->subset_type == TAO_SUBSET_MATRIXFREE)) { 2299566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &asls->w)); 230a7e14dcfSSatish Balay } 2319566063dSJacob Faibussowitsch PetscCall(TaoMatGetSubMat(tao->jacobian, asls->free, asls->w, tao->subset_type,&asls->J_sub)); 232a7e14dcfSSatish Balay if (tao->jacobian != tao->jacobian_pre) { 2339566063dSJacob Faibussowitsch PetscCall(TaoMatGetSubMat(tao->jacobian_pre, asls->free, asls->w, tao->subset_type, &asls->Jpre_sub)); 234a7e14dcfSSatish Balay } else { 2359566063dSJacob Faibussowitsch PetscCall(MatDestroy(&asls->Jpre_sub)); 236a7e14dcfSSatish Balay asls->Jpre_sub = asls->J_sub; 2379566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)(asls->Jpre_sub))); 238a7e14dcfSSatish Balay } 2399566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(asls->J_sub, asls->r1,ADD_VALUES)); 2409566063dSJacob Faibussowitsch PetscCall(TaoVecGetSubVec(tao->stepdirection, asls->free, tao->subset_type, 0.0, &asls->dxfree)); 2419566063dSJacob Faibussowitsch PetscCall(VecSet(asls->dxfree, 0.0)); 242a7e14dcfSSatish Balay 243a7e14dcfSSatish Balay /* Calculate the reduced direction. (Really negative of Newton 244a7e14dcfSSatish Balay direction. Therefore, rest of the code uses -d.) */ 2459566063dSJacob Faibussowitsch PetscCall(KSPReset(tao->ksp)); 2469566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(tao->ksp, asls->J_sub, asls->Jpre_sub)); 2479566063dSJacob Faibussowitsch PetscCall(KSPSolve(tao->ksp, asls->r2, asls->dxfree)); 2489566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(tao->ksp,&tao->ksp_its)); 249b0026674SJason Sarich tao->ksp_tot_its+=tao->ksp_its; 250a7e14dcfSSatish Balay 251a7e14dcfSSatish Balay /* Add the direction in the free variables back into the real direction. */ 2529566063dSJacob Faibussowitsch PetscCall(VecISAXPY(tao->stepdirection, asls->free, 1.0,asls->dxfree)); 253a7e14dcfSSatish Balay 254a7e14dcfSSatish Balay /* Check the real direction for descent and if not, use the negative 255a7e14dcfSSatish Balay gradient direction. */ 2569566063dSJacob Faibussowitsch PetscCall(VecNorm(tao->stepdirection, NORM_2, &normd)); 2579566063dSJacob Faibussowitsch PetscCall(VecDot(tao->stepdirection, asls->dpsi, &innerd)); 258a7e14dcfSSatish Balay 2591118d4bcSLisandro Dalcin if (innerd <= asls->delta*PetscPowReal(normd, asls->rho)) { 2609566063dSJacob Faibussowitsch PetscCall(PetscInfo(tao,"Gradient direction: %5.4e.\n", (double)innerd)); 26163a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(tao, "Iteration %" PetscInt_FMT ": newton direction not descent\n", tao->niter)); 2629566063dSJacob Faibussowitsch PetscCall(VecCopy(asls->dpsi, tao->stepdirection)); 2639566063dSJacob Faibussowitsch PetscCall(VecDot(asls->dpsi, tao->stepdirection, &innerd)); 264a7e14dcfSSatish Balay } 265a7e14dcfSSatish Balay 2669566063dSJacob Faibussowitsch PetscCall(VecScale(tao->stepdirection, -1.0)); 267a7e14dcfSSatish Balay innerd = -innerd; 268a7e14dcfSSatish Balay 269a7e14dcfSSatish Balay /* We now have a correct descent direction. Apply a linesearch to 270a7e14dcfSSatish Balay find the new iterate. */ 2719566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch, 1.0)); 2729566063dSJacob Faibussowitsch PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &psi,asls->dpsi, tao->stepdirection, &t, &ls_reason)); 2739566063dSJacob Faibussowitsch PetscCall(VecNorm(asls->dpsi, NORM_2, &ndpsi)); 274a7e14dcfSSatish Balay } 275a7e14dcfSSatish Balay PetscFunctionReturn(0); 276a7e14dcfSSatish Balay } 277a7e14dcfSSatish Balay 278a7e14dcfSSatish Balay /* ---------------------------------------------------------- */ 2791522df2eSJason Sarich /*MC 2801522df2eSJason Sarich TAOASILS - Active-set infeasible linesearch algorithm for solving 2811522df2eSJason Sarich complementarity constraints 2821522df2eSJason Sarich 2831522df2eSJason Sarich Options Database Keys: 2841522df2eSJason Sarich + -tao_ssls_delta - descent test fraction 2851522df2eSJason Sarich - -tao_ssls_rho - descent test power 2861522df2eSJason Sarich 2871eb8069cSJason Sarich Level: beginner 2881522df2eSJason Sarich M*/ 289728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_ASILS(Tao tao) 290a7e14dcfSSatish Balay { 291a7e14dcfSSatish Balay TAO_SSLS *asls; 2928caf6e8cSBarry Smith const char *armijo_type = TAOLINESEARCHARMIJO; 293a7e14dcfSSatish Balay 294a7e14dcfSSatish Balay PetscFunctionBegin; 2959566063dSJacob Faibussowitsch PetscCall(PetscNewLog(tao,&asls)); 296a7e14dcfSSatish Balay tao->data = (void*)asls; 297a7e14dcfSSatish Balay tao->ops->solve = TaoSolve_ASILS; 298a7e14dcfSSatish Balay tao->ops->setup = TaoSetUp_ASILS; 299a7e14dcfSSatish Balay tao->ops->view = TaoView_SSLS; 300a7e14dcfSSatish Balay tao->ops->setfromoptions = TaoSetFromOptions_SSLS; 301a7e14dcfSSatish Balay tao->ops->destroy = TaoDestroy_ASILS; 302a7e14dcfSSatish Balay tao->subset_type = TAO_SUBSET_SUBVEC; 303a7e14dcfSSatish Balay asls->delta = 1e-10; 304a7e14dcfSSatish Balay asls->rho = 2.1; 3056c23d075SBarry Smith asls->fixed = NULL; 3066c23d075SBarry Smith asls->free = NULL; 3076c23d075SBarry Smith asls->J_sub = NULL; 3086c23d075SBarry Smith asls->Jpre_sub = NULL; 3096c23d075SBarry Smith asls->w = NULL; 3106c23d075SBarry Smith asls->r1 = NULL; 3116c23d075SBarry Smith asls->r2 = NULL; 3126c23d075SBarry Smith asls->r3 = NULL; 3136c23d075SBarry Smith asls->t1 = NULL; 3146c23d075SBarry Smith asls->t2 = NULL; 3156c23d075SBarry Smith asls->dxfree = NULL; 316a7e14dcfSSatish Balay 317a7e14dcfSSatish Balay asls->identifier = 1e-5; 318a7e14dcfSSatish Balay 3199566063dSJacob Faibussowitsch PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch)); 3209566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1)); 3219566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetType(tao->linesearch, armijo_type)); 3229566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix)); 3239566063dSJacob Faibussowitsch PetscCall(TaoLineSearchSetFromOptions(tao->linesearch)); 324a7e14dcfSSatish Balay 3259566063dSJacob Faibussowitsch PetscCall(KSPCreate(((PetscObject)tao)->comm, &tao->ksp)); 3269566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1)); 3279566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix)); 3289566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(tao->ksp)); 3296552cf8aSJason Sarich 3306552cf8aSJason Sarich /* Override default settings (unless already changed) */ 3316552cf8aSJason Sarich if (!tao->max_it_changed) tao->max_it = 2000; 3326552cf8aSJason Sarich if (!tao->max_funcs_changed) tao->max_funcs = 4000; 3336552cf8aSJason Sarich if (!tao->gttol_changed) tao->gttol = 0; 3346552cf8aSJason Sarich if (!tao->grtol_changed) tao->grtol = 0; 3356f4723b1SBarry Smith #if defined(PETSC_USE_REAL_SINGLE) 3366552cf8aSJason Sarich if (!tao->gatol_changed) tao->gatol = 1.0e-6; 3376552cf8aSJason Sarich if (!tao->fmin_changed) tao->fmin = 1.0e-4; 3386f4723b1SBarry Smith #else 3396552cf8aSJason Sarich if (!tao->gatol_changed) tao->gatol = 1.0e-16; 3406552cf8aSJason Sarich if (!tao->fmin_changed) tao->fmin = 1.0e-8; 3416f4723b1SBarry Smith #endif 342a7e14dcfSSatish Balay PetscFunctionReturn(0); 343a7e14dcfSSatish Balay } 344