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 10d71ae5a4SJacob Faibussowitsch static PetscErrorCode QPIPComputeResidual(TAO_BQPIP *qp, Tao tao) 11d71ae5a4SJacob Faibussowitsch { 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); 313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3274962610SAlp Dener } 3374962610SAlp Dener 34d71ae5a4SJacob Faibussowitsch static PetscErrorCode QPIPSetInitialPoint(TAO_BQPIP *qp, Tao tao) 35d71ae5a4SJacob Faibussowitsch { 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 48*76c63389SBarry Smith PetscCheck(!PetscIsInfOrNanReal(qp->pobj), PETSC_COMM_SELF, PETSC_ERR_USER, "User provided data contains infinity 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 879371c9d4SSatish Balay qp->mu = 0; 889371c9d4SSatish Balay qp->dinfeas = 1.0; 899371c9d4SSatish Balay qp->pinfeas = 1.0; 9074962610SAlp Dener while ((qp->dinfeas + qp->pinfeas) / (qp->m + qp->n) >= qp->mu) { 919566063dSJacob Faibussowitsch PetscCall(VecScale(qp->G, two)); 929566063dSJacob Faibussowitsch PetscCall(VecScale(qp->Z, two)); 939566063dSJacob Faibussowitsch PetscCall(VecScale(qp->S, two)); 949566063dSJacob Faibussowitsch PetscCall(VecScale(qp->T, two)); 9574962610SAlp Dener 969566063dSJacob Faibussowitsch PetscCall(QPIPComputeResidual(qp, tao)); 9774962610SAlp Dener 989566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution, qp->R3)); 999566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->R3, -1.0, qp->G)); 1009566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->R3, -1.0, qp->XL)); 10174962610SAlp Dener 1029566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->solution, qp->R5)); 1039566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->R5, 1.0, qp->T)); 1049566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->R5, -1.0, qp->XU)); 10574962610SAlp Dener 1069566063dSJacob Faibussowitsch PetscCall(VecNorm(qp->R3, NORM_INFINITY, &gap1)); 1079566063dSJacob Faibussowitsch PetscCall(VecNorm(qp->R5, NORM_INFINITY, &gap2)); 10874962610SAlp Dener qp->pinfeas = PetscMax(gap1, gap2); 10974962610SAlp Dener 11074962610SAlp Dener /* Compute the duality gap */ 1119566063dSJacob Faibussowitsch PetscCall(VecDot(qp->G, qp->Z, &gap1)); 1129566063dSJacob Faibussowitsch PetscCall(VecDot(qp->T, qp->S, &gap2)); 11374962610SAlp Dener 11474962610SAlp Dener qp->gap = gap1 + gap2; 11574962610SAlp Dener qp->dobj = qp->pobj - qp->gap; 11674962610SAlp Dener if (qp->m > 0) { 11774962610SAlp Dener qp->mu = qp->gap / (qp->m); 11874962610SAlp Dener } else { 11974962610SAlp Dener qp->mu = 0.0; 12074962610SAlp Dener } 12174962610SAlp Dener qp->rgap = qp->gap / (PetscAbsReal(qp->dobj) + PetscAbsReal(qp->pobj) + 1.0); 12274962610SAlp Dener } 1233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12474962610SAlp Dener } 12574962610SAlp Dener 126d71ae5a4SJacob Faibussowitsch static PetscErrorCode QPIPStepLength(TAO_BQPIP *qp) 127d71ae5a4SJacob Faibussowitsch { 12874962610SAlp Dener PetscReal tstep1, tstep2, tstep3, tstep4, tstep; 12974962610SAlp Dener 13074962610SAlp Dener PetscFunctionBegin; 13174962610SAlp Dener /* Compute stepsize to the boundary */ 1329566063dSJacob Faibussowitsch PetscCall(VecStepMax(qp->G, qp->DG, &tstep1)); 1339566063dSJacob Faibussowitsch PetscCall(VecStepMax(qp->T, qp->DT, &tstep2)); 1349566063dSJacob Faibussowitsch PetscCall(VecStepMax(qp->S, qp->DS, &tstep3)); 1359566063dSJacob Faibussowitsch PetscCall(VecStepMax(qp->Z, qp->DZ, &tstep4)); 13674962610SAlp Dener 13774962610SAlp Dener tstep = PetscMin(tstep1, tstep2); 13874962610SAlp Dener qp->psteplength = PetscMin(0.95 * tstep, 1.0); 13974962610SAlp Dener 14074962610SAlp Dener tstep = PetscMin(tstep3, tstep4); 14174962610SAlp Dener qp->dsteplength = PetscMin(0.95 * tstep, 1.0); 14274962610SAlp Dener 14374962610SAlp Dener qp->psteplength = PetscMin(qp->psteplength, qp->dsteplength); 14474962610SAlp Dener qp->dsteplength = qp->psteplength; 1453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14674962610SAlp Dener } 14774962610SAlp Dener 148d71ae5a4SJacob Faibussowitsch static PetscErrorCode QPIPComputeNormFromCentralPath(TAO_BQPIP *qp, PetscReal *norm) 149d71ae5a4SJacob Faibussowitsch { 15074962610SAlp Dener PetscReal gap[2], mu[2], nmu; 15174962610SAlp Dener 15274962610SAlp Dener PetscFunctionBegin; 1539566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(qp->GZwork, qp->G, qp->Z)); 1549566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(qp->TSwork, qp->T, qp->S)); 1559566063dSJacob Faibussowitsch PetscCall(VecNorm(qp->TSwork, NORM_1, &mu[0])); 1569566063dSJacob Faibussowitsch PetscCall(VecNorm(qp->GZwork, NORM_1, &mu[1])); 15774962610SAlp Dener 15874962610SAlp Dener nmu = -(mu[0] + mu[1]) / qp->m; 15974962610SAlp Dener 1609566063dSJacob Faibussowitsch PetscCall(VecShift(qp->GZwork, nmu)); 1619566063dSJacob Faibussowitsch PetscCall(VecShift(qp->TSwork, nmu)); 16274962610SAlp Dener 1639566063dSJacob Faibussowitsch PetscCall(VecNorm(qp->GZwork, NORM_2, &gap[0])); 1649566063dSJacob Faibussowitsch PetscCall(VecNorm(qp->TSwork, NORM_2, &gap[1])); 16574962610SAlp Dener gap[0] *= gap[0]; 16674962610SAlp Dener gap[1] *= gap[1]; 16774962610SAlp Dener 16874962610SAlp Dener qp->pathnorm = PetscSqrtScalar(gap[0] + gap[1]); 16974962610SAlp Dener *norm = qp->pathnorm; 1703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 17174962610SAlp Dener } 17274962610SAlp Dener 173d71ae5a4SJacob Faibussowitsch static PetscErrorCode QPIPComputeStepDirection(TAO_BQPIP *qp, Tao tao) 174d71ae5a4SJacob Faibussowitsch { 17574962610SAlp Dener PetscFunctionBegin; 17674962610SAlp Dener /* Calculate DG */ 1779566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->stepdirection, qp->DG)); 1789566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->DG, 1.0, qp->R3)); 17974962610SAlp Dener 18074962610SAlp Dener /* Calculate DT */ 1819566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->stepdirection, qp->DT)); 1829566063dSJacob Faibussowitsch PetscCall(VecScale(qp->DT, -1.0)); 1839566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->DT, -1.0, qp->R5)); 18474962610SAlp Dener 18574962610SAlp Dener /* Calculate DZ */ 1869566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->DZ, -1.0, qp->Z)); 1879566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(qp->GZwork, qp->DG, qp->G)); 1889566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(qp->GZwork, qp->GZwork, qp->Z)); 1899566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->DZ, -1.0, qp->GZwork)); 19074962610SAlp Dener 19174962610SAlp Dener /* Calculate DS */ 1929566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->DS, -1.0, qp->S)); 1939566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(qp->TSwork, qp->DT, qp->T)); 1949566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(qp->TSwork, qp->TSwork, qp->S)); 1959566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->DS, -1.0, qp->TSwork)); 1963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 19774962610SAlp Dener } 19874962610SAlp Dener 199d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSetup_BQPIP(Tao tao) 200d71ae5a4SJacob Faibussowitsch { 20174962610SAlp Dener TAO_BQPIP *qp = (TAO_BQPIP *)tao->data; 20274962610SAlp Dener 20374962610SAlp Dener PetscFunctionBegin; 20474962610SAlp Dener /* Set pointers to Data */ 2059566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->solution, &qp->n)); 20674962610SAlp Dener 20774962610SAlp Dener /* Allocate some arrays */ 20848a46eb9SPierre Jolivet if (!tao->gradient) PetscCall(VecDuplicate(tao->solution, &tao->gradient)); 20948a46eb9SPierre Jolivet if (!tao->stepdirection) PetscCall(VecDuplicate(tao->solution, &tao->stepdirection)); 21074962610SAlp Dener 2119566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &qp->Work)); 2129566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &qp->XU)); 2139566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &qp->XL)); 2149566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &qp->HDiag)); 2159566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &qp->DiagAxpy)); 2169566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &qp->RHS)); 2179566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &qp->RHS2)); 2189566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &qp->C)); 21974962610SAlp Dener 2209566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &qp->G)); 2219566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &qp->DG)); 2229566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &qp->S)); 2239566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &qp->Z)); 2249566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &qp->DZ)); 2259566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &qp->GZwork)); 2269566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &qp->R3)); 22774962610SAlp Dener 2289566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &qp->T)); 2299566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &qp->DT)); 2309566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &qp->DS)); 2319566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &qp->TSwork)); 2329566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &qp->R5)); 23374962610SAlp Dener qp->m = 2 * qp->n; 2343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 23574962610SAlp Dener } 23674962610SAlp Dener 237d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSolve_BQPIP(Tao tao) 238d71ae5a4SJacob Faibussowitsch { 23974962610SAlp Dener TAO_BQPIP *qp = (TAO_BQPIP *)tao->data; 24074962610SAlp Dener PetscInt its; 24174962610SAlp Dener PetscReal d1, d2, ksptol, sigmamu; 24274962610SAlp Dener PetscReal gnorm, dstep, pstep, step = 0; 24374962610SAlp Dener PetscReal gap[4]; 24474962610SAlp Dener PetscBool getdiagop; 24574962610SAlp Dener 24674962610SAlp Dener PetscFunctionBegin; 24774962610SAlp Dener qp->dobj = 0.0; 24874962610SAlp Dener qp->pobj = 1.0; 24974962610SAlp Dener qp->gap = 10.0; 25074962610SAlp Dener qp->rgap = 1.0; 25174962610SAlp Dener qp->mu = 1.0; 25274962610SAlp Dener qp->dinfeas = 1.0; 25374962610SAlp Dener qp->psteplength = 0.0; 25474962610SAlp Dener qp->dsteplength = 0.0; 25574962610SAlp Dener 25674962610SAlp Dener /* TODO 25774962610SAlp Dener - Remove fixed variables and treat them correctly 25874962610SAlp Dener - Use index sets for the infinite versus finite bounds 25974962610SAlp Dener - Update remaining code for fixed and free variables 26074962610SAlp Dener - Fix inexact solves for predictor and corrector 26174962610SAlp Dener */ 26274962610SAlp Dener 26374962610SAlp Dener /* Tighten infinite bounds, things break when we don't do this 26474962610SAlp Dener -- see test_bqpip.c 26574962610SAlp Dener */ 2669566063dSJacob Faibussowitsch PetscCall(TaoComputeVariableBounds(tao)); 2679566063dSJacob Faibussowitsch PetscCall(VecSet(qp->XU, 1.0e20)); 2689566063dSJacob Faibussowitsch PetscCall(VecSet(qp->XL, -1.0e20)); 26976be6f4fSStefano Zampini if (tao->XL) PetscCall(VecPointwiseMax(qp->XL, qp->XL, tao->XL)); 27076be6f4fSStefano Zampini if (tao->XU) PetscCall(VecPointwiseMin(qp->XU, qp->XU, tao->XU)); 2719566063dSJacob Faibussowitsch PetscCall(VecMedian(qp->XL, tao->solution, qp->XU, tao->solution)); 27274962610SAlp Dener 27374962610SAlp Dener /* Evaluate gradient and Hessian at zero to get the correct values 27474962610SAlp Dener without contaminating them with numerical artifacts. 27574962610SAlp Dener */ 2769566063dSJacob Faibussowitsch PetscCall(VecSet(qp->Work, 0)); 2779566063dSJacob Faibussowitsch PetscCall(TaoComputeObjectiveAndGradient(tao, qp->Work, &qp->d, qp->C)); 2789566063dSJacob Faibussowitsch PetscCall(TaoComputeHessian(tao, qp->Work, tao->hessian, tao->hessian_pre)); 2799566063dSJacob Faibussowitsch PetscCall(MatHasOperation(tao->hessian, MATOP_GET_DIAGONAL, &getdiagop)); 2801baa6e33SBarry Smith if (getdiagop) PetscCall(MatGetDiagonal(tao->hessian, qp->HDiag)); 28174962610SAlp Dener 28274962610SAlp Dener /* Initialize starting point and residuals */ 2839566063dSJacob Faibussowitsch PetscCall(QPIPSetInitialPoint(qp, tao)); 2849566063dSJacob Faibussowitsch PetscCall(QPIPComputeResidual(qp, tao)); 28574962610SAlp Dener 28674962610SAlp Dener /* Enter main loop */ 28774962610SAlp Dener tao->reason = TAO_CONTINUE_ITERATING; 28874962610SAlp Dener while (1) { 28974962610SAlp Dener /* Check Stopping Condition */ 29074962610SAlp Dener gnorm = PetscSqrtScalar(qp->gap + qp->dinfeas); 2919566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, qp->pobj, gnorm, qp->pinfeas, tao->ksp_its)); 2929566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, qp->pobj, gnorm, qp->pinfeas, step)); 293dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP); 29474962610SAlp Dener if (tao->reason != TAO_CONTINUE_ITERATING) break; 295e1e80dc8SAlp Dener /* Call general purpose update function */ 296dbbe0bcdSBarry Smith PetscTryTypeMethod(tao, update, tao->niter, tao->user_update); 29774962610SAlp Dener tao->niter++; 29874962610SAlp Dener tao->ksp_its = 0; 29974962610SAlp Dener 30074962610SAlp Dener /* 30174962610SAlp Dener Dual Infeasibility Direction should already be in the right 30274962610SAlp Dener hand side from computing the residuals 30374962610SAlp Dener */ 30474962610SAlp Dener 3059566063dSJacob Faibussowitsch PetscCall(QPIPComputeNormFromCentralPath(qp, &d1)); 30674962610SAlp Dener 3079566063dSJacob Faibussowitsch PetscCall(VecSet(qp->DZ, 0.0)); 3089566063dSJacob Faibussowitsch PetscCall(VecSet(qp->DS, 0.0)); 30974962610SAlp Dener 31074962610SAlp Dener /* 31174962610SAlp Dener Compute the Primal Infeasiblitiy RHS and the 31274962610SAlp Dener Diagonal Matrix to be added to H and store in Work 31374962610SAlp Dener */ 3149566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(qp->DiagAxpy, qp->Z, qp->G)); 3159566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(qp->GZwork, qp->DiagAxpy, qp->R3)); 3169566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->RHS, -1.0, qp->GZwork)); 31774962610SAlp Dener 3189566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(qp->TSwork, qp->S, qp->T)); 3199566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->DiagAxpy, 1.0, qp->TSwork)); 3209566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(qp->TSwork, qp->TSwork, qp->R5)); 3219566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->RHS, -1.0, qp->TSwork)); 32274962610SAlp Dener 32374962610SAlp Dener /* Determine the solving tolerance */ 32474962610SAlp Dener ksptol = qp->mu / 10.0; 32574962610SAlp Dener ksptol = PetscMin(ksptol, 0.001); 3269566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(tao->ksp, ksptol, 1e-30, 1e30, PetscMax(10, qp->n))); 32774962610SAlp Dener 32874962610SAlp Dener /* Shift the diagonals of the Hessian matrix */ 3299566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(tao->hessian, qp->DiagAxpy, ADD_VALUES)); 33074962610SAlp Dener if (!getdiagop) { 3319566063dSJacob Faibussowitsch PetscCall(VecCopy(qp->DiagAxpy, qp->HDiag)); 3329566063dSJacob Faibussowitsch PetscCall(VecScale(qp->HDiag, -1.0)); 33374962610SAlp Dener } 3349566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(tao->hessian, MAT_FINAL_ASSEMBLY)); 3359566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(tao->hessian, MAT_FINAL_ASSEMBLY)); 33674962610SAlp Dener 3379566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(tao->ksp, tao->hessian, tao->hessian_pre)); 3389566063dSJacob Faibussowitsch PetscCall(KSPSolve(tao->ksp, qp->RHS, tao->stepdirection)); 3399566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(tao->ksp, &its)); 34074962610SAlp Dener tao->ksp_its += its; 34174962610SAlp Dener tao->ksp_tot_its += its; 34274962610SAlp Dener 34374962610SAlp Dener /* Restore the true diagonal of the Hessian matrix */ 34474962610SAlp Dener if (getdiagop) { 3459566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(tao->hessian, qp->HDiag, INSERT_VALUES)); 34674962610SAlp Dener } else { 3479566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(tao->hessian, qp->HDiag, ADD_VALUES)); 34874962610SAlp Dener } 3499566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(tao->hessian, MAT_FINAL_ASSEMBLY)); 3509566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(tao->hessian, MAT_FINAL_ASSEMBLY)); 3519566063dSJacob Faibussowitsch PetscCall(QPIPComputeStepDirection(qp, tao)); 3529566063dSJacob Faibussowitsch PetscCall(QPIPStepLength(qp)); 35374962610SAlp Dener 35474962610SAlp Dener /* Calculate New Residual R1 in Work vector */ 3559566063dSJacob Faibussowitsch PetscCall(MatMult(tao->hessian, tao->stepdirection, qp->RHS2)); 3569566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->RHS2, 1.0, qp->DS)); 3579566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->RHS2, -1.0, qp->DZ)); 3589566063dSJacob Faibussowitsch PetscCall(VecAYPX(qp->RHS2, qp->dsteplength, tao->gradient)); 35974962610SAlp Dener 3609566063dSJacob Faibussowitsch PetscCall(VecNorm(qp->RHS2, NORM_2, &qp->dinfeas)); 3619566063dSJacob Faibussowitsch PetscCall(VecDot(qp->DZ, qp->DG, gap)); 3629566063dSJacob Faibussowitsch PetscCall(VecDot(qp->DS, qp->DT, gap + 1)); 36374962610SAlp Dener 36474962610SAlp Dener qp->rnorm = (qp->dinfeas + qp->psteplength * qp->pinfeas) / (qp->m + qp->n); 36574962610SAlp Dener pstep = qp->psteplength; 36674962610SAlp Dener step = PetscMin(qp->psteplength, qp->dsteplength); 36774962610SAlp Dener sigmamu = (pstep * pstep * (gap[0] + gap[1]) + (1 - pstep) * qp->gap) / qp->m; 36874962610SAlp Dener 36974962610SAlp Dener if (qp->predcorr && step < 0.9) { 37074962610SAlp Dener if (sigmamu < qp->mu) { 37174962610SAlp Dener sigmamu = sigmamu / qp->mu; 37274962610SAlp Dener sigmamu = sigmamu * sigmamu * sigmamu; 37374962610SAlp Dener } else { 37474962610SAlp Dener sigmamu = 1.0; 37574962610SAlp Dener } 37674962610SAlp Dener sigmamu = sigmamu * qp->mu; 37774962610SAlp Dener 37874962610SAlp Dener /* Compute Corrector Step */ 3799566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(qp->DZ, qp->DG, qp->DZ)); 3809566063dSJacob Faibussowitsch PetscCall(VecScale(qp->DZ, -1.0)); 3819566063dSJacob Faibussowitsch PetscCall(VecShift(qp->DZ, sigmamu)); 3829566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(qp->DZ, qp->DZ, qp->G)); 38374962610SAlp Dener 3849566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(qp->DS, qp->DS, qp->DT)); 3859566063dSJacob Faibussowitsch PetscCall(VecScale(qp->DS, -1.0)); 3869566063dSJacob Faibussowitsch PetscCall(VecShift(qp->DS, sigmamu)); 3879566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(qp->DS, qp->DS, qp->T)); 38874962610SAlp Dener 3899566063dSJacob Faibussowitsch PetscCall(VecCopy(qp->DZ, qp->RHS2)); 3909566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->RHS2, -1.0, qp->DS)); 3919566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->RHS2, 1.0, qp->RHS)); 39274962610SAlp Dener 39374962610SAlp Dener /* Approximately solve the linear system */ 3949566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(tao->hessian, qp->DiagAxpy, ADD_VALUES)); 39574962610SAlp Dener if (!getdiagop) { 3969566063dSJacob Faibussowitsch PetscCall(VecCopy(qp->DiagAxpy, qp->HDiag)); 3979566063dSJacob Faibussowitsch PetscCall(VecScale(qp->HDiag, -1.0)); 39874962610SAlp Dener } 3999566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(tao->hessian, MAT_FINAL_ASSEMBLY)); 4009566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(tao->hessian, MAT_FINAL_ASSEMBLY)); 40174962610SAlp Dener 40274962610SAlp Dener /* Solve using the previous tolerances that were set */ 4039566063dSJacob Faibussowitsch PetscCall(KSPSolve(tao->ksp, qp->RHS2, tao->stepdirection)); 4049566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(tao->ksp, &its)); 40574962610SAlp Dener tao->ksp_its += its; 40674962610SAlp Dener tao->ksp_tot_its += its; 40774962610SAlp Dener 40874962610SAlp Dener if (getdiagop) { 4099566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(tao->hessian, qp->HDiag, INSERT_VALUES)); 41074962610SAlp Dener } else { 4119566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(tao->hessian, qp->HDiag, ADD_VALUES)); 41274962610SAlp Dener } 4139566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(tao->hessian, MAT_FINAL_ASSEMBLY)); 4149566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(tao->hessian, MAT_FINAL_ASSEMBLY)); 4159566063dSJacob Faibussowitsch PetscCall(QPIPComputeStepDirection(qp, tao)); 4169566063dSJacob Faibussowitsch PetscCall(QPIPStepLength(qp)); 41774962610SAlp Dener } /* End Corrector step */ 41874962610SAlp Dener 41974962610SAlp Dener /* Take the step */ 42074962610SAlp Dener dstep = qp->dsteplength; 42174962610SAlp Dener 4229566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->Z, dstep, qp->DZ)); 4239566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->S, dstep, qp->DS)); 4249566063dSJacob Faibussowitsch PetscCall(VecAXPY(tao->solution, dstep, tao->stepdirection)); 4259566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->G, dstep, qp->DG)); 4269566063dSJacob Faibussowitsch PetscCall(VecAXPY(qp->T, dstep, qp->DT)); 42774962610SAlp Dener 42874962610SAlp Dener /* Compute Residuals */ 4299566063dSJacob Faibussowitsch PetscCall(QPIPComputeResidual(qp, tao)); 43074962610SAlp Dener 43174962610SAlp Dener /* Evaluate quadratic function */ 4329566063dSJacob Faibussowitsch PetscCall(MatMult(tao->hessian, tao->solution, qp->Work)); 43374962610SAlp Dener 4349566063dSJacob Faibussowitsch PetscCall(VecDot(tao->solution, qp->Work, &d1)); 4359566063dSJacob Faibussowitsch PetscCall(VecDot(tao->solution, qp->C, &d2)); 4369566063dSJacob Faibussowitsch PetscCall(VecDot(qp->G, qp->Z, gap)); 4379566063dSJacob Faibussowitsch PetscCall(VecDot(qp->T, qp->S, gap + 1)); 43874962610SAlp Dener 43974962610SAlp Dener /* Compute the duality gap */ 44074962610SAlp Dener qp->pobj = d1 / 2.0 + d2 + qp->d; 44174962610SAlp Dener qp->gap = gap[0] + gap[1]; 44274962610SAlp Dener qp->dobj = qp->pobj - qp->gap; 443ad540459SPierre Jolivet if (qp->m > 0) qp->mu = qp->gap / (qp->m); 44474962610SAlp Dener qp->rgap = qp->gap / (PetscAbsReal(qp->dobj) + PetscAbsReal(qp->pobj) + 1.0); 44574962610SAlp Dener } /* END MAIN LOOP */ 4463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 44774962610SAlp Dener } 44874962610SAlp Dener 449d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoView_BQPIP(Tao tao, PetscViewer viewer) 450d71ae5a4SJacob Faibussowitsch { 45174962610SAlp Dener PetscFunctionBegin; 4523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 45374962610SAlp Dener } 45474962610SAlp Dener 455ce78bad3SBarry Smith static PetscErrorCode TaoSetFromOptions_BQPIP(Tao tao, PetscOptionItems PetscOptionsObject) 456d71ae5a4SJacob Faibussowitsch { 45774962610SAlp Dener TAO_BQPIP *qp = (TAO_BQPIP *)tao->data; 45874962610SAlp Dener 45974962610SAlp Dener PetscFunctionBegin; 460d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Interior point method for bound constrained quadratic optimization"); 4619566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-tao_bqpip_predcorr", "Use a predictor-corrector method", "", qp->predcorr, &qp->predcorr, NULL)); 462d0609cedSBarry Smith PetscOptionsHeadEnd(); 4639566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(tao->ksp)); 4643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 46574962610SAlp Dener } 46674962610SAlp Dener 467d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoDestroy_BQPIP(Tao tao) 468d71ae5a4SJacob Faibussowitsch { 46974962610SAlp Dener TAO_BQPIP *qp = (TAO_BQPIP *)tao->data; 47074962610SAlp Dener 47174962610SAlp Dener PetscFunctionBegin; 47274962610SAlp Dener if (tao->setupcalled) { 4739566063dSJacob Faibussowitsch PetscCall(VecDestroy(&qp->G)); 4749566063dSJacob Faibussowitsch PetscCall(VecDestroy(&qp->DG)); 4759566063dSJacob Faibussowitsch PetscCall(VecDestroy(&qp->Z)); 4769566063dSJacob Faibussowitsch PetscCall(VecDestroy(&qp->DZ)); 4779566063dSJacob Faibussowitsch PetscCall(VecDestroy(&qp->GZwork)); 4789566063dSJacob Faibussowitsch PetscCall(VecDestroy(&qp->R3)); 4799566063dSJacob Faibussowitsch PetscCall(VecDestroy(&qp->S)); 4809566063dSJacob Faibussowitsch PetscCall(VecDestroy(&qp->DS)); 4819566063dSJacob Faibussowitsch PetscCall(VecDestroy(&qp->T)); 48274962610SAlp Dener 4839566063dSJacob Faibussowitsch PetscCall(VecDestroy(&qp->DT)); 4849566063dSJacob Faibussowitsch PetscCall(VecDestroy(&qp->TSwork)); 4859566063dSJacob Faibussowitsch PetscCall(VecDestroy(&qp->R5)); 4869566063dSJacob Faibussowitsch PetscCall(VecDestroy(&qp->HDiag)); 4879566063dSJacob Faibussowitsch PetscCall(VecDestroy(&qp->Work)); 4889566063dSJacob Faibussowitsch PetscCall(VecDestroy(&qp->XL)); 4899566063dSJacob Faibussowitsch PetscCall(VecDestroy(&qp->XU)); 4909566063dSJacob Faibussowitsch PetscCall(VecDestroy(&qp->DiagAxpy)); 4919566063dSJacob Faibussowitsch PetscCall(VecDestroy(&qp->RHS)); 4929566063dSJacob Faibussowitsch PetscCall(VecDestroy(&qp->RHS2)); 4939566063dSJacob Faibussowitsch PetscCall(VecDestroy(&qp->C)); 49474962610SAlp Dener } 495a958fbfcSStefano Zampini PetscCall(KSPDestroy(&tao->ksp)); 4969566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->data)); 4973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 49874962610SAlp Dener } 49974962610SAlp Dener 500d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoComputeDual_BQPIP(Tao tao, Vec DXL, Vec DXU) 501d71ae5a4SJacob Faibussowitsch { 50274962610SAlp Dener TAO_BQPIP *qp = (TAO_BQPIP *)tao->data; 50374962610SAlp Dener 50474962610SAlp Dener PetscFunctionBegin; 5059566063dSJacob Faibussowitsch PetscCall(VecCopy(qp->Z, DXL)); 5069566063dSJacob Faibussowitsch PetscCall(VecCopy(qp->S, DXU)); 5079566063dSJacob Faibussowitsch PetscCall(VecScale(DXU, -1.0)); 5083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 50974962610SAlp Dener } 51074962610SAlp Dener 51174962610SAlp Dener /*MC 51274962610SAlp Dener TAOBQPIP - interior-point method for quadratic programs with 51374962610SAlp Dener box constraints. 51474962610SAlp Dener 51574962610SAlp Dener Notes: 51674962610SAlp Dener This algorithms solves quadratic problems only, the Hessian will 51774962610SAlp Dener only be computed once. 51874962610SAlp Dener 51974962610SAlp Dener Options Database Keys: 52074962610SAlp Dener . -tao_bqpip_predcorr - use a predictor/corrector method 52174962610SAlp Dener 52274962610SAlp Dener Level: beginner 52374962610SAlp Dener M*/ 52474962610SAlp Dener 525d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoCreate_BQPIP(Tao tao) 526d71ae5a4SJacob Faibussowitsch { 52774962610SAlp Dener TAO_BQPIP *qp; 52874962610SAlp Dener 52974962610SAlp Dener PetscFunctionBegin; 5304dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&qp)); 53174962610SAlp Dener 53274962610SAlp Dener tao->ops->setup = TaoSetup_BQPIP; 53374962610SAlp Dener tao->ops->solve = TaoSolve_BQPIP; 53474962610SAlp Dener tao->ops->view = TaoView_BQPIP; 53574962610SAlp Dener tao->ops->setfromoptions = TaoSetFromOptions_BQPIP; 53674962610SAlp Dener tao->ops->destroy = TaoDestroy_BQPIP; 53774962610SAlp Dener tao->ops->computedual = TaoComputeDual_BQPIP; 53874962610SAlp Dener 53974962610SAlp Dener /* Override default settings (unless already changed) */ 540606f75f6SBarry Smith PetscCall(TaoParametersInitialize(tao)); 541606f75f6SBarry Smith PetscObjectParameterSetDefault(tao, max_it, 100); 542606f75f6SBarry Smith PetscObjectParameterSetDefault(tao, max_funcs, 500); 543606f75f6SBarry Smith PetscObjectParameterSetDefault(tao, catol, PetscDefined(USE_REAL_SINGLE) ? 1e-6 : 1e-12); 54474962610SAlp Dener 54574962610SAlp Dener /* Initialize pointers and variables */ 54674962610SAlp Dener qp->n = 0; 54774962610SAlp Dener qp->m = 0; 54874962610SAlp Dener 54974962610SAlp Dener qp->predcorr = 1; 55074962610SAlp Dener tao->data = (void *)qp; 55174962610SAlp Dener 5529566063dSJacob Faibussowitsch PetscCall(KSPCreate(((PetscObject)tao)->comm, &tao->ksp)); 5539566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1)); 5549566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix)); 5559566063dSJacob Faibussowitsch PetscCall(KSPSetType(tao->ksp, KSPCG)); 5569566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(tao->ksp, 1e-14, 1e-30, 1e30, PetscMax(10, qp->n))); 5573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 55874962610SAlp Dener } 559