174962610SAlp Dener /* 274962610SAlp Dener This file implements a Mehrotra predictor-corrector method for 374962610SAlp Dener bound-constrained quadratic programs. 474962610SAlp Dener 574962610SAlp Dener */ 674962610SAlp Dener 774962610SAlp Dener #include <../src/tao/quadratic/impls/bqpip/bqpipimpl.h> 874962610SAlp Dener #include <petscksp.h> 974962610SAlp Dener 1074962610SAlp Dener static PetscErrorCode QPIPComputeResidual(TAO_BQPIP *qp,Tao tao) 1174962610SAlp Dener { 1274962610SAlp Dener PetscReal dtmp = 1.0 - qp->psteplength; 1374962610SAlp Dener 1474962610SAlp Dener PetscFunctionBegin; 1574962610SAlp Dener /* Compute R3 and R5 */ 1674962610SAlp Dener 179566063dSJacob Faibussowitsch PetscCall(VecScale(qp->R3,dtmp)); 189566063dSJacob Faibussowitsch PetscCall(VecScale(qp->R5,dtmp)); 1974962610SAlp Dener qp->pinfeas=dtmp*qp->pinfeas; 2074962610SAlp Dener 219566063dSJacob Faibussowitsch PetscCall(VecCopy(qp->S,tao->gradient)); 229566063dSJacob Faibussowitsch PetscCall(VecAXPY(tao->gradient,-1.0,qp->Z)); 2374962610SAlp Dener 249566063dSJacob Faibussowitsch PetscCall(MatMult(tao->hessian,tao->solution,qp->RHS)); 259566063dSJacob Faibussowitsch PetscCall(VecScale(qp->RHS,-1.0)); 269566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->RHS,-1.0,qp->C)); 279566063dSJacob Faibussowitsch PetscCall(VecAXPY(tao->gradient,-1.0,qp->RHS)); 2874962610SAlp Dener 299566063dSJacob Faibussowitsch PetscCall(VecNorm(tao->gradient,NORM_1,&qp->dinfeas)); 3074962610SAlp Dener qp->rnorm=(qp->dinfeas+qp->pinfeas)/(qp->m+qp->n); 3174962610SAlp Dener PetscFunctionReturn(0); 3274962610SAlp Dener } 3374962610SAlp Dener 3474962610SAlp Dener static PetscErrorCode QPIPSetInitialPoint(TAO_BQPIP *qp, Tao tao) 3574962610SAlp Dener { 3674962610SAlp Dener PetscReal two=2.0,p01=1; 3774962610SAlp Dener PetscReal gap1,gap2,fff,mu; 3874962610SAlp Dener 3974962610SAlp Dener PetscFunctionBegin; 4074962610SAlp Dener /* Compute function, Gradient R=Hx+b, and Hessian */ 419566063dSJacob Faibussowitsch PetscCall(MatMult(tao->hessian,tao->solution,tao->gradient)); 429566063dSJacob Faibussowitsch PetscCall(VecCopy(qp->C,qp->Work)); 439566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->Work,0.5,tao->gradient)); 449566063dSJacob Faibussowitsch PetscCall(VecAXPY(tao->gradient,1.0,qp->C)); 459566063dSJacob Faibussowitsch PetscCall(VecDot(tao->solution,qp->Work,&fff)); 4674962610SAlp Dener qp->pobj = fff + qp->d; 4774962610SAlp Dener 483c859ba3SBarry Smith PetscCheck(!PetscIsInfOrNanReal(qp->pobj),PETSC_COMM_SELF,PETSC_ERR_USER, "User provided data contains Inf or NaN"); 4974962610SAlp Dener 5074962610SAlp Dener /* Initialize slack vectors */ 5174962610SAlp Dener /* T = XU - X; G = X - XL */ 529566063dSJacob Faibussowitsch PetscCall(VecCopy(qp->XU,qp->T)); 539566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->T,-1.0,tao->solution)); 549566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution,qp->G)); 559566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->G,-1.0,qp->XL)); 5674962610SAlp Dener 579566063dSJacob Faibussowitsch PetscCall(VecSet(qp->GZwork,p01)); 589566063dSJacob Faibussowitsch PetscCall(VecSet(qp->TSwork,p01)); 5974962610SAlp Dener 609566063dSJacob Faibussowitsch PetscCall(VecPointwiseMax(qp->G,qp->G,qp->GZwork)); 619566063dSJacob Faibussowitsch PetscCall(VecPointwiseMax(qp->T,qp->T,qp->TSwork)); 6274962610SAlp Dener 6374962610SAlp Dener /* Initialize Dual Variable Vectors */ 649566063dSJacob Faibussowitsch PetscCall(VecCopy(qp->G,qp->Z)); 659566063dSJacob Faibussowitsch PetscCall(VecReciprocal(qp->Z)); 6674962610SAlp Dener 679566063dSJacob Faibussowitsch PetscCall(VecCopy(qp->T,qp->S)); 689566063dSJacob Faibussowitsch PetscCall(VecReciprocal(qp->S)); 6974962610SAlp Dener 709566063dSJacob Faibussowitsch PetscCall(MatMult(tao->hessian,qp->Work,qp->RHS)); 719566063dSJacob Faibussowitsch PetscCall(VecAbs(qp->RHS)); 729566063dSJacob Faibussowitsch PetscCall(VecSet(qp->Work,p01)); 739566063dSJacob Faibussowitsch PetscCall(VecPointwiseMax(qp->RHS,qp->RHS,qp->Work)); 7474962610SAlp Dener 759566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(qp->RHS,tao->gradient,qp->RHS)); 769566063dSJacob Faibussowitsch PetscCall(VecNorm(qp->RHS,NORM_1,&gap1)); 7774962610SAlp Dener mu = PetscMin(10.0,(gap1+10.0)/qp->m); 7874962610SAlp Dener 799566063dSJacob Faibussowitsch PetscCall(VecScale(qp->S,mu)); 809566063dSJacob Faibussowitsch PetscCall(VecScale(qp->Z,mu)); 8174962610SAlp Dener 829566063dSJacob Faibussowitsch PetscCall(VecSet(qp->TSwork,p01)); 839566063dSJacob Faibussowitsch PetscCall(VecSet(qp->GZwork,p01)); 849566063dSJacob Faibussowitsch PetscCall(VecPointwiseMax(qp->S,qp->S,qp->TSwork)); 859566063dSJacob Faibussowitsch PetscCall(VecPointwiseMax(qp->Z,qp->Z,qp->GZwork)); 8674962610SAlp Dener 8774962610SAlp Dener qp->mu=0;qp->dinfeas=1.0;qp->pinfeas=1.0; 8874962610SAlp Dener while ((qp->dinfeas+qp->pinfeas)/(qp->m+qp->n) >= qp->mu) { 899566063dSJacob Faibussowitsch PetscCall(VecScale(qp->G,two)); 909566063dSJacob Faibussowitsch PetscCall(VecScale(qp->Z,two)); 919566063dSJacob Faibussowitsch PetscCall(VecScale(qp->S,two)); 929566063dSJacob Faibussowitsch PetscCall(VecScale(qp->T,two)); 9374962610SAlp Dener 949566063dSJacob Faibussowitsch PetscCall(QPIPComputeResidual(qp,tao)); 9574962610SAlp Dener 969566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution,qp->R3)); 979566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->R3,-1.0,qp->G)); 989566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->R3,-1.0,qp->XL)); 9974962610SAlp Dener 1009566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution,qp->R5)); 1019566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->R5,1.0,qp->T)); 1029566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->R5,-1.0,qp->XU)); 10374962610SAlp Dener 1049566063dSJacob Faibussowitsch PetscCall(VecNorm(qp->R3,NORM_INFINITY,&gap1)); 1059566063dSJacob Faibussowitsch PetscCall(VecNorm(qp->R5,NORM_INFINITY,&gap2)); 10674962610SAlp Dener qp->pinfeas=PetscMax(gap1,gap2); 10774962610SAlp Dener 10874962610SAlp Dener /* Compute the duality gap */ 1099566063dSJacob Faibussowitsch PetscCall(VecDot(qp->G,qp->Z,&gap1)); 1109566063dSJacob Faibussowitsch PetscCall(VecDot(qp->T,qp->S,&gap2)); 11174962610SAlp Dener 11274962610SAlp Dener qp->gap = gap1+gap2; 11374962610SAlp Dener qp->dobj = qp->pobj - qp->gap; 11474962610SAlp Dener if (qp->m>0) { 11574962610SAlp Dener qp->mu=qp->gap/(qp->m); 11674962610SAlp Dener } else { 11774962610SAlp Dener qp->mu=0.0; 11874962610SAlp Dener } 11974962610SAlp Dener qp->rgap=qp->gap/(PetscAbsReal(qp->dobj) + PetscAbsReal(qp->pobj) + 1.0); 12074962610SAlp Dener } 12174962610SAlp Dener PetscFunctionReturn(0); 12274962610SAlp Dener } 12374962610SAlp Dener 12474962610SAlp Dener static PetscErrorCode QPIPStepLength(TAO_BQPIP *qp) 12574962610SAlp Dener { 12674962610SAlp Dener PetscReal tstep1,tstep2,tstep3,tstep4,tstep; 12774962610SAlp Dener 12874962610SAlp Dener PetscFunctionBegin; 12974962610SAlp Dener /* Compute stepsize to the boundary */ 1309566063dSJacob Faibussowitsch PetscCall(VecStepMax(qp->G,qp->DG,&tstep1)); 1319566063dSJacob Faibussowitsch PetscCall(VecStepMax(qp->T,qp->DT,&tstep2)); 1329566063dSJacob Faibussowitsch PetscCall(VecStepMax(qp->S,qp->DS,&tstep3)); 1339566063dSJacob Faibussowitsch PetscCall(VecStepMax(qp->Z,qp->DZ,&tstep4)); 13474962610SAlp Dener 13574962610SAlp Dener tstep = PetscMin(tstep1,tstep2); 13674962610SAlp Dener qp->psteplength = PetscMin(0.95*tstep,1.0); 13774962610SAlp Dener 13874962610SAlp Dener tstep = PetscMin(tstep3,tstep4); 13974962610SAlp Dener qp->dsteplength = PetscMin(0.95*tstep,1.0); 14074962610SAlp Dener 14174962610SAlp Dener qp->psteplength = PetscMin(qp->psteplength,qp->dsteplength); 14274962610SAlp Dener qp->dsteplength = qp->psteplength; 14374962610SAlp Dener PetscFunctionReturn(0); 14474962610SAlp Dener } 14574962610SAlp Dener 14674962610SAlp Dener static PetscErrorCode QPIPComputeNormFromCentralPath(TAO_BQPIP *qp,PetscReal *norm) 14774962610SAlp Dener { 14874962610SAlp Dener PetscReal gap[2],mu[2],nmu; 14974962610SAlp Dener 15074962610SAlp Dener PetscFunctionBegin; 1519566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(qp->GZwork,qp->G,qp->Z)); 1529566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(qp->TSwork,qp->T,qp->S)); 1539566063dSJacob Faibussowitsch PetscCall(VecNorm(qp->TSwork,NORM_1,&mu[0])); 1549566063dSJacob Faibussowitsch PetscCall(VecNorm(qp->GZwork,NORM_1,&mu[1])); 15574962610SAlp Dener 15674962610SAlp Dener nmu=-(mu[0]+mu[1])/qp->m; 15774962610SAlp Dener 1589566063dSJacob Faibussowitsch PetscCall(VecShift(qp->GZwork,nmu)); 1599566063dSJacob Faibussowitsch PetscCall(VecShift(qp->TSwork,nmu)); 16074962610SAlp Dener 1619566063dSJacob Faibussowitsch PetscCall(VecNorm(qp->GZwork,NORM_2,&gap[0])); 1629566063dSJacob Faibussowitsch PetscCall(VecNorm(qp->TSwork,NORM_2,&gap[1])); 16374962610SAlp Dener gap[0]*=gap[0]; 16474962610SAlp Dener gap[1]*=gap[1]; 16574962610SAlp Dener 16674962610SAlp Dener qp->pathnorm=PetscSqrtScalar(gap[0]+gap[1]); 16774962610SAlp Dener *norm=qp->pathnorm; 16874962610SAlp Dener PetscFunctionReturn(0); 16974962610SAlp Dener } 17074962610SAlp Dener 17174962610SAlp Dener static PetscErrorCode QPIPComputeStepDirection(TAO_BQPIP *qp,Tao tao) 17274962610SAlp Dener { 17374962610SAlp Dener PetscFunctionBegin; 17474962610SAlp Dener /* Calculate DG */ 1759566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->stepdirection,qp->DG)); 1769566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->DG,1.0,qp->R3)); 17774962610SAlp Dener 17874962610SAlp Dener /* Calculate DT */ 1799566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->stepdirection,qp->DT)); 1809566063dSJacob Faibussowitsch PetscCall(VecScale(qp->DT,-1.0)); 1819566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->DT,-1.0,qp->R5)); 18274962610SAlp Dener 18374962610SAlp Dener /* Calculate DZ */ 1849566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->DZ,-1.0,qp->Z)); 1859566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(qp->GZwork,qp->DG,qp->G)); 1869566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(qp->GZwork,qp->GZwork,qp->Z)); 1879566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->DZ,-1.0,qp->GZwork)); 18874962610SAlp Dener 18974962610SAlp Dener /* Calculate DS */ 1909566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->DS,-1.0,qp->S)); 1919566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(qp->TSwork,qp->DT,qp->T)); 1929566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(qp->TSwork,qp->TSwork,qp->S)); 1939566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->DS,-1.0,qp->TSwork)); 19474962610SAlp Dener PetscFunctionReturn(0); 19574962610SAlp Dener } 19674962610SAlp Dener 19774962610SAlp Dener static PetscErrorCode TaoSetup_BQPIP(Tao tao) 19874962610SAlp Dener { 19974962610SAlp Dener TAO_BQPIP *qp =(TAO_BQPIP*)tao->data; 20074962610SAlp Dener 20174962610SAlp Dener PetscFunctionBegin; 20274962610SAlp Dener /* Set pointers to Data */ 2039566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->solution,&qp->n)); 20474962610SAlp Dener 20574962610SAlp Dener /* Allocate some arrays */ 20674962610SAlp Dener if (!tao->gradient) { 2079566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&tao->gradient)); 20874962610SAlp Dener } 20974962610SAlp Dener if (!tao->stepdirection) { 2109566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&tao->stepdirection)); 21174962610SAlp Dener } 21274962610SAlp Dener 2139566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&qp->Work)); 2149566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&qp->XU)); 2159566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&qp->XL)); 2169566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&qp->HDiag)); 2179566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&qp->DiagAxpy)); 2189566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&qp->RHS)); 2199566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&qp->RHS2)); 2209566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&qp->C)); 22174962610SAlp Dener 2229566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&qp->G)); 2239566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&qp->DG)); 2249566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&qp->S)); 2259566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&qp->Z)); 2269566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&qp->DZ)); 2279566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&qp->GZwork)); 2289566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&qp->R3)); 22974962610SAlp Dener 2309566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&qp->T)); 2319566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&qp->DT)); 2329566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&qp->DS)); 2339566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&qp->TSwork)); 2349566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution,&qp->R5)); 23574962610SAlp Dener qp->m=2*qp->n; 23674962610SAlp Dener PetscFunctionReturn(0); 23774962610SAlp Dener } 23874962610SAlp Dener 23974962610SAlp Dener static PetscErrorCode TaoSolve_BQPIP(Tao tao) 24074962610SAlp Dener { 24174962610SAlp Dener TAO_BQPIP *qp = (TAO_BQPIP*)tao->data; 24274962610SAlp Dener PetscInt its; 24374962610SAlp Dener PetscReal d1,d2,ksptol,sigmamu; 24474962610SAlp Dener PetscReal gnorm,dstep,pstep,step=0; 24574962610SAlp Dener PetscReal gap[4]; 24674962610SAlp Dener PetscBool getdiagop; 24774962610SAlp Dener 24874962610SAlp Dener PetscFunctionBegin; 24974962610SAlp Dener qp->dobj = 0.0; 25074962610SAlp Dener qp->pobj = 1.0; 25174962610SAlp Dener qp->gap = 10.0; 25274962610SAlp Dener qp->rgap = 1.0; 25374962610SAlp Dener qp->mu = 1.0; 25474962610SAlp Dener qp->dinfeas = 1.0; 25574962610SAlp Dener qp->psteplength = 0.0; 25674962610SAlp Dener qp->dsteplength = 0.0; 25774962610SAlp Dener 25874962610SAlp Dener /* TODO 25974962610SAlp Dener - Remove fixed variables and treat them correctly 26074962610SAlp Dener - Use index sets for the infinite versus finite bounds 26174962610SAlp Dener - Update remaining code for fixed and free variables 26274962610SAlp Dener - Fix inexact solves for predictor and corrector 26374962610SAlp Dener */ 26474962610SAlp Dener 26574962610SAlp Dener /* Tighten infinite bounds, things break when we don't do this 26674962610SAlp Dener -- see test_bqpip.c 26774962610SAlp Dener */ 2689566063dSJacob Faibussowitsch PetscCall(TaoComputeVariableBounds(tao)); 2699566063dSJacob Faibussowitsch PetscCall(VecSet(qp->XU,1.0e20)); 2709566063dSJacob Faibussowitsch PetscCall(VecSet(qp->XL,-1.0e20)); 271*76be6f4fSStefano Zampini if (tao->XL) PetscCall(VecPointwiseMax(qp->XL,qp->XL,tao->XL)); 272*76be6f4fSStefano Zampini if (tao->XU) PetscCall(VecPointwiseMin(qp->XU,qp->XU,tao->XU)); 2739566063dSJacob Faibussowitsch PetscCall(VecMedian(qp->XL,tao->solution,qp->XU,tao->solution)); 27474962610SAlp Dener 27574962610SAlp Dener /* Evaluate gradient and Hessian at zero to get the correct values 27674962610SAlp Dener without contaminating them with numerical artifacts. 27774962610SAlp Dener */ 2789566063dSJacob Faibussowitsch PetscCall(VecSet(qp->Work,0)); 2799566063dSJacob Faibussowitsch PetscCall(TaoComputeObjectiveAndGradient(tao,qp->Work,&qp->d,qp->C)); 2809566063dSJacob Faibussowitsch PetscCall(TaoComputeHessian(tao,qp->Work,tao->hessian,tao->hessian_pre)); 2819566063dSJacob Faibussowitsch PetscCall(MatHasOperation(tao->hessian,MATOP_GET_DIAGONAL,&getdiagop)); 28274962610SAlp Dener if (getdiagop) { 2839566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(tao->hessian,qp->HDiag)); 28474962610SAlp Dener } 28574962610SAlp Dener 28674962610SAlp Dener /* Initialize starting point and residuals */ 2879566063dSJacob Faibussowitsch PetscCall(QPIPSetInitialPoint(qp,tao)); 2889566063dSJacob Faibussowitsch PetscCall(QPIPComputeResidual(qp,tao)); 28974962610SAlp Dener 29074962610SAlp Dener /* Enter main loop */ 29174962610SAlp Dener tao->reason = TAO_CONTINUE_ITERATING; 29274962610SAlp Dener while (1) { 29374962610SAlp Dener 29474962610SAlp Dener /* Check Stopping Condition */ 29574962610SAlp Dener gnorm = PetscSqrtScalar(qp->gap + qp->dinfeas); 2969566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao,qp->pobj,gnorm,qp->pinfeas,tao->ksp_its)); 2979566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao,tao->niter,qp->pobj,gnorm,qp->pinfeas,step)); 2989566063dSJacob Faibussowitsch PetscCall((*tao->ops->convergencetest)(tao,tao->cnvP)); 29974962610SAlp Dener if (tao->reason != TAO_CONTINUE_ITERATING) break; 300e1e80dc8SAlp Dener /* Call general purpose update function */ 301e1e80dc8SAlp Dener if (tao->ops->update) { 3029566063dSJacob Faibussowitsch PetscCall((*tao->ops->update)(tao, tao->niter, tao->user_update)); 303e1e80dc8SAlp Dener } 30474962610SAlp Dener tao->niter++; 30574962610SAlp Dener tao->ksp_its = 0; 30674962610SAlp Dener 30774962610SAlp Dener /* 30874962610SAlp Dener Dual Infeasibility Direction should already be in the right 30974962610SAlp Dener hand side from computing the residuals 31074962610SAlp Dener */ 31174962610SAlp Dener 3129566063dSJacob Faibussowitsch PetscCall(QPIPComputeNormFromCentralPath(qp,&d1)); 31374962610SAlp Dener 3149566063dSJacob Faibussowitsch PetscCall(VecSet(qp->DZ,0.0)); 3159566063dSJacob Faibussowitsch PetscCall(VecSet(qp->DS,0.0)); 31674962610SAlp Dener 31774962610SAlp Dener /* 31874962610SAlp Dener Compute the Primal Infeasiblitiy RHS and the 31974962610SAlp Dener Diagonal Matrix to be added to H and store in Work 32074962610SAlp Dener */ 3219566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(qp->DiagAxpy,qp->Z,qp->G)); 3229566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(qp->GZwork,qp->DiagAxpy,qp->R3)); 3239566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->RHS,-1.0,qp->GZwork)); 32474962610SAlp Dener 3259566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(qp->TSwork,qp->S,qp->T)); 3269566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->DiagAxpy,1.0,qp->TSwork)); 3279566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(qp->TSwork,qp->TSwork,qp->R5)); 3289566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->RHS,-1.0,qp->TSwork)); 32974962610SAlp Dener 33074962610SAlp Dener /* Determine the solving tolerance */ 33174962610SAlp Dener ksptol = qp->mu/10.0; 33274962610SAlp Dener ksptol = PetscMin(ksptol,0.001); 3339566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(tao->ksp,ksptol,1e-30,1e30,PetscMax(10,qp->n))); 33474962610SAlp Dener 33574962610SAlp Dener /* Shift the diagonals of the Hessian matrix */ 3369566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(tao->hessian,qp->DiagAxpy,ADD_VALUES)); 33774962610SAlp Dener if (!getdiagop) { 3389566063dSJacob Faibussowitsch PetscCall(VecCopy(qp->DiagAxpy,qp->HDiag)); 3399566063dSJacob Faibussowitsch PetscCall(VecScale(qp->HDiag,-1.0)); 34074962610SAlp Dener } 3419566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(tao->hessian,MAT_FINAL_ASSEMBLY)); 3429566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(tao->hessian,MAT_FINAL_ASSEMBLY)); 34374962610SAlp Dener 3449566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(tao->ksp,tao->hessian,tao->hessian_pre)); 3459566063dSJacob Faibussowitsch PetscCall(KSPSolve(tao->ksp,qp->RHS,tao->stepdirection)); 3469566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(tao->ksp,&its)); 34774962610SAlp Dener tao->ksp_its += its; 34874962610SAlp Dener tao->ksp_tot_its += its; 34974962610SAlp Dener 35074962610SAlp Dener /* Restore the true diagonal of the Hessian matrix */ 35174962610SAlp Dener if (getdiagop) { 3529566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(tao->hessian,qp->HDiag,INSERT_VALUES)); 35374962610SAlp Dener } else { 3549566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(tao->hessian,qp->HDiag,ADD_VALUES)); 35574962610SAlp Dener } 3569566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(tao->hessian,MAT_FINAL_ASSEMBLY)); 3579566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(tao->hessian,MAT_FINAL_ASSEMBLY)); 3589566063dSJacob Faibussowitsch PetscCall(QPIPComputeStepDirection(qp,tao)); 3599566063dSJacob Faibussowitsch PetscCall(QPIPStepLength(qp)); 36074962610SAlp Dener 36174962610SAlp Dener /* Calculate New Residual R1 in Work vector */ 3629566063dSJacob Faibussowitsch PetscCall(MatMult(tao->hessian,tao->stepdirection,qp->RHS2)); 3639566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->RHS2,1.0,qp->DS)); 3649566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->RHS2,-1.0,qp->DZ)); 3659566063dSJacob Faibussowitsch PetscCall(VecAYPX(qp->RHS2,qp->dsteplength,tao->gradient)); 36674962610SAlp Dener 3679566063dSJacob Faibussowitsch PetscCall(VecNorm(qp->RHS2,NORM_2,&qp->dinfeas)); 3689566063dSJacob Faibussowitsch PetscCall(VecDot(qp->DZ,qp->DG,gap)); 3699566063dSJacob Faibussowitsch PetscCall(VecDot(qp->DS,qp->DT,gap+1)); 37074962610SAlp Dener 37174962610SAlp Dener qp->rnorm = (qp->dinfeas+qp->psteplength*qp->pinfeas)/(qp->m+qp->n); 37274962610SAlp Dener pstep = qp->psteplength; 37374962610SAlp Dener step = PetscMin(qp->psteplength,qp->dsteplength); 37474962610SAlp Dener sigmamu = (pstep*pstep*(gap[0]+gap[1]) + (1 - pstep)*qp->gap)/qp->m; 37574962610SAlp Dener 37674962610SAlp Dener if (qp->predcorr && step < 0.9) { 37774962610SAlp Dener if (sigmamu < qp->mu) { 37874962610SAlp Dener sigmamu = sigmamu/qp->mu; 37974962610SAlp Dener sigmamu = sigmamu*sigmamu*sigmamu; 38074962610SAlp Dener } else { 38174962610SAlp Dener sigmamu = 1.0; 38274962610SAlp Dener } 38374962610SAlp Dener sigmamu = sigmamu*qp->mu; 38474962610SAlp Dener 38574962610SAlp Dener /* Compute Corrector Step */ 3869566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(qp->DZ,qp->DG,qp->DZ)); 3879566063dSJacob Faibussowitsch PetscCall(VecScale(qp->DZ,-1.0)); 3889566063dSJacob Faibussowitsch PetscCall(VecShift(qp->DZ,sigmamu)); 3899566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(qp->DZ,qp->DZ,qp->G)); 39074962610SAlp Dener 3919566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(qp->DS,qp->DS,qp->DT)); 3929566063dSJacob Faibussowitsch PetscCall(VecScale(qp->DS,-1.0)); 3939566063dSJacob Faibussowitsch PetscCall(VecShift(qp->DS,sigmamu)); 3949566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(qp->DS,qp->DS,qp->T)); 39574962610SAlp Dener 3969566063dSJacob Faibussowitsch PetscCall(VecCopy(qp->DZ,qp->RHS2)); 3979566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->RHS2,-1.0,qp->DS)); 3989566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->RHS2,1.0,qp->RHS)); 39974962610SAlp Dener 40074962610SAlp Dener /* Approximately solve the linear system */ 4019566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(tao->hessian,qp->DiagAxpy,ADD_VALUES)); 40274962610SAlp Dener if (!getdiagop) { 4039566063dSJacob Faibussowitsch PetscCall(VecCopy(qp->DiagAxpy,qp->HDiag)); 4049566063dSJacob Faibussowitsch PetscCall(VecScale(qp->HDiag,-1.0)); 40574962610SAlp Dener } 4069566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(tao->hessian,MAT_FINAL_ASSEMBLY)); 4079566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(tao->hessian,MAT_FINAL_ASSEMBLY)); 40874962610SAlp Dener 40974962610SAlp Dener /* Solve using the previous tolerances that were set */ 4109566063dSJacob Faibussowitsch PetscCall(KSPSolve(tao->ksp,qp->RHS2,tao->stepdirection)); 4119566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(tao->ksp,&its)); 41274962610SAlp Dener tao->ksp_its += its; 41374962610SAlp Dener tao->ksp_tot_its += its; 41474962610SAlp Dener 41574962610SAlp Dener if (getdiagop) { 4169566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(tao->hessian,qp->HDiag,INSERT_VALUES)); 41774962610SAlp Dener } else { 4189566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(tao->hessian,qp->HDiag,ADD_VALUES)); 41974962610SAlp Dener } 4209566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(tao->hessian,MAT_FINAL_ASSEMBLY)); 4219566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(tao->hessian,MAT_FINAL_ASSEMBLY)); 4229566063dSJacob Faibussowitsch PetscCall(QPIPComputeStepDirection(qp,tao)); 4239566063dSJacob Faibussowitsch PetscCall(QPIPStepLength(qp)); 42474962610SAlp Dener } /* End Corrector step */ 42574962610SAlp Dener 42674962610SAlp Dener /* Take the step */ 42774962610SAlp Dener dstep = qp->dsteplength; 42874962610SAlp Dener 4299566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->Z,dstep,qp->DZ)); 4309566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->S,dstep,qp->DS)); 4319566063dSJacob Faibussowitsch PetscCall(VecAXPY(tao->solution,dstep,tao->stepdirection)); 4329566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->G,dstep,qp->DG)); 4339566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->T,dstep,qp->DT)); 43474962610SAlp Dener 43574962610SAlp Dener /* Compute Residuals */ 4369566063dSJacob Faibussowitsch PetscCall(QPIPComputeResidual(qp,tao)); 43774962610SAlp Dener 43874962610SAlp Dener /* Evaluate quadratic function */ 4399566063dSJacob Faibussowitsch PetscCall(MatMult(tao->hessian,tao->solution,qp->Work)); 44074962610SAlp Dener 4419566063dSJacob Faibussowitsch PetscCall(VecDot(tao->solution,qp->Work,&d1)); 4429566063dSJacob Faibussowitsch PetscCall(VecDot(tao->solution,qp->C,&d2)); 4439566063dSJacob Faibussowitsch PetscCall(VecDot(qp->G,qp->Z,gap)); 4449566063dSJacob Faibussowitsch PetscCall(VecDot(qp->T,qp->S,gap+1)); 44574962610SAlp Dener 44674962610SAlp Dener /* Compute the duality gap */ 44774962610SAlp Dener qp->pobj = d1/2.0 + d2+qp->d; 44874962610SAlp Dener qp->gap = gap[0]+gap[1]; 44974962610SAlp Dener qp->dobj = qp->pobj - qp->gap; 45074962610SAlp Dener if (qp->m > 0) { 45174962610SAlp Dener qp->mu = qp->gap/(qp->m); 45274962610SAlp Dener } 45374962610SAlp Dener qp->rgap = qp->gap/(PetscAbsReal(qp->dobj) + PetscAbsReal(qp->pobj) + 1.0); 45474962610SAlp Dener } /* END MAIN LOOP */ 45574962610SAlp Dener PetscFunctionReturn(0); 45674962610SAlp Dener } 45774962610SAlp Dener 45874962610SAlp Dener static PetscErrorCode TaoView_BQPIP(Tao tao,PetscViewer viewer) 45974962610SAlp Dener { 46074962610SAlp Dener PetscFunctionBegin; 46174962610SAlp Dener PetscFunctionReturn(0); 46274962610SAlp Dener } 46374962610SAlp Dener 46474962610SAlp Dener static PetscErrorCode TaoSetFromOptions_BQPIP(PetscOptionItems *PetscOptionsObject,Tao tao) 46574962610SAlp Dener { 46674962610SAlp Dener TAO_BQPIP *qp = (TAO_BQPIP*)tao->data; 46774962610SAlp Dener 46874962610SAlp Dener PetscFunctionBegin; 469d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"Interior point method for bound constrained quadratic optimization"); 4709566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-tao_bqpip_predcorr","Use a predictor-corrector method","",qp->predcorr,&qp->predcorr,NULL)); 471d0609cedSBarry Smith PetscOptionsHeadEnd(); 4729566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(tao->ksp)); 47374962610SAlp Dener PetscFunctionReturn(0); 47474962610SAlp Dener } 47574962610SAlp Dener 47674962610SAlp Dener static PetscErrorCode TaoDestroy_BQPIP(Tao tao) 47774962610SAlp Dener { 47874962610SAlp Dener TAO_BQPIP *qp = (TAO_BQPIP*)tao->data; 47974962610SAlp Dener 48074962610SAlp Dener PetscFunctionBegin; 48174962610SAlp Dener if (tao->setupcalled) { 4829566063dSJacob Faibussowitsch PetscCall(VecDestroy(&qp->G)); 4839566063dSJacob Faibussowitsch PetscCall(VecDestroy(&qp->DG)); 4849566063dSJacob Faibussowitsch PetscCall(VecDestroy(&qp->Z)); 4859566063dSJacob Faibussowitsch PetscCall(VecDestroy(&qp->DZ)); 4869566063dSJacob Faibussowitsch PetscCall(VecDestroy(&qp->GZwork)); 4879566063dSJacob Faibussowitsch PetscCall(VecDestroy(&qp->R3)); 4889566063dSJacob Faibussowitsch PetscCall(VecDestroy(&qp->S)); 4899566063dSJacob Faibussowitsch PetscCall(VecDestroy(&qp->DS)); 4909566063dSJacob Faibussowitsch PetscCall(VecDestroy(&qp->T)); 49174962610SAlp Dener 4929566063dSJacob Faibussowitsch PetscCall(VecDestroy(&qp->DT)); 4939566063dSJacob Faibussowitsch PetscCall(VecDestroy(&qp->TSwork)); 4949566063dSJacob Faibussowitsch PetscCall(VecDestroy(&qp->R5)); 4959566063dSJacob Faibussowitsch PetscCall(VecDestroy(&qp->HDiag)); 4969566063dSJacob Faibussowitsch PetscCall(VecDestroy(&qp->Work)); 4979566063dSJacob Faibussowitsch PetscCall(VecDestroy(&qp->XL)); 4989566063dSJacob Faibussowitsch PetscCall(VecDestroy(&qp->XU)); 4999566063dSJacob Faibussowitsch PetscCall(VecDestroy(&qp->DiagAxpy)); 5009566063dSJacob Faibussowitsch PetscCall(VecDestroy(&qp->RHS)); 5019566063dSJacob Faibussowitsch PetscCall(VecDestroy(&qp->RHS2)); 5029566063dSJacob Faibussowitsch PetscCall(VecDestroy(&qp->C)); 50374962610SAlp Dener } 5049566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->data)); 50574962610SAlp Dener PetscFunctionReturn(0); 50674962610SAlp Dener } 50774962610SAlp Dener 50874962610SAlp Dener static PetscErrorCode TaoComputeDual_BQPIP(Tao tao,Vec DXL,Vec DXU) 50974962610SAlp Dener { 51074962610SAlp Dener TAO_BQPIP *qp = (TAO_BQPIP*)tao->data; 51174962610SAlp Dener 51274962610SAlp Dener PetscFunctionBegin; 5139566063dSJacob Faibussowitsch PetscCall(VecCopy(qp->Z,DXL)); 5149566063dSJacob Faibussowitsch PetscCall(VecCopy(qp->S,DXU)); 5159566063dSJacob Faibussowitsch PetscCall(VecScale(DXU,-1.0)); 51674962610SAlp Dener PetscFunctionReturn(0); 51774962610SAlp Dener } 51874962610SAlp Dener 51974962610SAlp Dener /*MC 52074962610SAlp Dener TAOBQPIP - interior-point method for quadratic programs with 52174962610SAlp Dener box constraints. 52274962610SAlp Dener 52374962610SAlp Dener Notes: 52474962610SAlp Dener This algorithms solves quadratic problems only, the Hessian will 52574962610SAlp Dener only be computed once. 52674962610SAlp Dener 52774962610SAlp Dener Options Database Keys: 52874962610SAlp Dener . -tao_bqpip_predcorr - use a predictor/corrector method 52974962610SAlp Dener 53074962610SAlp Dener Level: beginner 53174962610SAlp Dener M*/ 53274962610SAlp Dener 53374962610SAlp Dener PETSC_EXTERN PetscErrorCode TaoCreate_BQPIP(Tao tao) 53474962610SAlp Dener { 53574962610SAlp Dener TAO_BQPIP *qp; 53674962610SAlp Dener 53774962610SAlp Dener PetscFunctionBegin; 5389566063dSJacob Faibussowitsch PetscCall(PetscNewLog(tao,&qp)); 53974962610SAlp Dener 54074962610SAlp Dener tao->ops->setup = TaoSetup_BQPIP; 54174962610SAlp Dener tao->ops->solve = TaoSolve_BQPIP; 54274962610SAlp Dener tao->ops->view = TaoView_BQPIP; 54374962610SAlp Dener tao->ops->setfromoptions = TaoSetFromOptions_BQPIP; 54474962610SAlp Dener tao->ops->destroy = TaoDestroy_BQPIP; 54574962610SAlp Dener tao->ops->computedual = TaoComputeDual_BQPIP; 54674962610SAlp Dener 54774962610SAlp Dener /* Override default settings (unless already changed) */ 54874962610SAlp Dener if (!tao->max_it_changed) tao->max_it=100; 54974962610SAlp Dener if (!tao->max_funcs_changed) tao->max_funcs = 500; 55074962610SAlp Dener #if defined(PETSC_USE_REAL_SINGLE) 55174962610SAlp Dener if (!tao->catol_changed) tao->catol=1e-6; 55274962610SAlp Dener #else 55374962610SAlp Dener if (!tao->catol_changed) tao->catol=1e-12; 55474962610SAlp Dener #endif 55574962610SAlp Dener 55674962610SAlp Dener /* Initialize pointers and variables */ 55774962610SAlp Dener qp->n = 0; 55874962610SAlp Dener qp->m = 0; 55974962610SAlp Dener 56074962610SAlp Dener qp->predcorr = 1; 56174962610SAlp Dener tao->data = (void*)qp; 56274962610SAlp Dener 5639566063dSJacob Faibussowitsch PetscCall(KSPCreate(((PetscObject)tao)->comm,&tao->ksp)); 5649566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1)); 5659566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(tao->ksp,tao->hdr.prefix)); 5669566063dSJacob Faibussowitsch PetscCall(KSPSetType(tao->ksp,KSPCG)); 5679566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(tao->ksp,1e-14,1e-30,1e30,PetscMax(10,qp->n))); 56874962610SAlp Dener PetscFunctionReturn(0); 56974962610SAlp Dener } 570