16285c0a3SHansol Suh #include <../src/tao/constrained/impls/admm/admm.h> /*I "petsctao.h" I*/ 26285c0a3SHansol Suh #include <petsctao.h> 36285c0a3SHansol Suh #include <petsc/private/petscimpl.h> 46285c0a3SHansol Suh 56285c0a3SHansol Suh /* Updates terminating criteria 66285c0a3SHansol Suh * 76285c0a3SHansol Suh * 1 ||r_k|| = ||Ax+Bz-c|| =< catol_admm* max{||Ax||,||Bz||,||c||} 86285c0a3SHansol Suh * 96285c0a3SHansol Suh * 2. Updates dual residual, d_k 106285c0a3SHansol Suh * 116285c0a3SHansol Suh * 3. ||d_k|| = ||mu*A^T*B(z_k-z_{k-1})|| =< gatol_admm * ||A^Ty|| */ 126285c0a3SHansol Suh 136285c0a3SHansol Suh static PetscBool cited = PETSC_FALSE; 146285c0a3SHansol Suh static const char citation[] = 156285c0a3SHansol Suh "@misc{xu2017adaptive,\n" 166285c0a3SHansol Suh " title={Adaptive Relaxed ADMM: Convergence Theory and Practical Implementation},\n" 176285c0a3SHansol Suh " author={Zheng Xu and Mario A. T. Figueiredo and Xiaoming Yuan and Christoph Studer and Tom Goldstein},\n" 186285c0a3SHansol Suh " year={2017},\n" 196285c0a3SHansol Suh " eprint={1704.02712},\n" 206285c0a3SHansol Suh " archivePrefix={arXiv},\n" 216285c0a3SHansol Suh " primaryClass={cs.CV}\n" 226285c0a3SHansol Suh "} \n"; 236285c0a3SHansol Suh 246285c0a3SHansol Suh static PetscErrorCode TaoADMMToleranceUpdate(Tao tao) 256285c0a3SHansol Suh { 266285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 276285c0a3SHansol Suh PetscErrorCode ierr; 286285c0a3SHansol Suh PetscReal Axnorm,Bznorm,ATynorm,temp; 296285c0a3SHansol Suh Vec tempJR,tempL; 306285c0a3SHansol Suh Tao mis; 316285c0a3SHansol Suh 326285c0a3SHansol Suh PetscFunctionBegin; 336285c0a3SHansol Suh mis = am->subsolverX; 346285c0a3SHansol Suh tempJR = am->workJacobianRight; 356285c0a3SHansol Suh tempL = am->workLeft; 366285c0a3SHansol Suh /* ATy */ 376285c0a3SHansol Suh ierr = TaoComputeJacobianEquality(mis, am->y, mis->jacobian_equality, mis->jacobian_equality_pre);CHKERRQ(ierr); 386285c0a3SHansol Suh ierr = MatMultTranspose(mis->jacobian_equality,am->y,tempJR);CHKERRQ(ierr); 396285c0a3SHansol Suh ierr = VecNorm(tempJR,NORM_2,&ATynorm);CHKERRQ(ierr); 406285c0a3SHansol Suh /* dualres = mu * ||AT(Bz-Bzold)||_2 */ 416285c0a3SHansol Suh ierr = VecWAXPY(tempJR,-1.,am->Bzold,am->Bz);CHKERRQ(ierr); 426285c0a3SHansol Suh ierr = MatMultTranspose(mis->jacobian_equality,tempJR,tempL);CHKERRQ(ierr); 436285c0a3SHansol Suh ierr = VecNorm(tempL,NORM_2,&am->dualres);CHKERRQ(ierr); 446285c0a3SHansol Suh am->dualres *= am->mu; 456285c0a3SHansol Suh 466285c0a3SHansol Suh /* ||Ax||_2, ||Bz||_2 */ 476285c0a3SHansol Suh ierr = VecNorm(am->Ax,NORM_2,&Axnorm);CHKERRQ(ierr); 486285c0a3SHansol Suh ierr = VecNorm(am->Bz,NORM_2,&Bznorm);CHKERRQ(ierr); 496285c0a3SHansol Suh 506285c0a3SHansol Suh /* Set catol to be catol_admm * max{||Ax||,||Bz||,||c||} * 516285c0a3SHansol Suh * Set gatol to be gatol_admm * ||A^Ty|| * 526285c0a3SHansol Suh * while cnorm is ||r_k||_2, and gnorm is ||d_k||_2 */ 536285c0a3SHansol Suh temp = am->catol_admm * PetscMax(Axnorm, (!am->const_norm) ? Bznorm : PetscMax(Bznorm,am->const_norm)); 546285c0a3SHansol Suh ierr = TaoSetConstraintTolerances(tao,temp,PETSC_DEFAULT);CHKERRQ(ierr); 556285c0a3SHansol Suh ierr = TaoSetTolerances(tao, am->gatol_admm*ATynorm, PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 566285c0a3SHansol Suh PetscFunctionReturn(0); 576285c0a3SHansol Suh } 586285c0a3SHansol Suh 596285c0a3SHansol Suh /* Penaly Update for Adaptive ADMM. */ 606285c0a3SHansol Suh static PetscErrorCode AdaptiveADMMPenaltyUpdate(Tao tao) 616285c0a3SHansol Suh { 626285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 636285c0a3SHansol Suh PetscErrorCode ierr; 646285c0a3SHansol Suh PetscReal ydiff_norm, yhatdiff_norm, Axdiff_norm, Bzdiff_norm, Axyhat, Bzy, a_sd, a_mg, a_k, b_sd, b_mg, b_k; 656285c0a3SHansol Suh PetscBool hflag, gflag; 666285c0a3SHansol Suh Vec tempJR,tempJR2; 676285c0a3SHansol Suh 686285c0a3SHansol Suh PetscFunctionBegin; 696285c0a3SHansol Suh tempJR = am->workJacobianRight; 706285c0a3SHansol Suh tempJR2 = am->workJacobianRight2; 716285c0a3SHansol Suh hflag = PETSC_FALSE; 726285c0a3SHansol Suh gflag = PETSC_FALSE; 736285c0a3SHansol Suh a_k = -1; 746285c0a3SHansol Suh b_k = -1; 756285c0a3SHansol Suh 766285c0a3SHansol Suh ierr = VecWAXPY(tempJR,-1.,am->Axold,am->Ax);CHKERRQ(ierr); 776285c0a3SHansol Suh ierr = VecWAXPY(tempJR2,-1.,am->yhatold,am->yhat);CHKERRQ(ierr); 786285c0a3SHansol Suh ierr = VecNorm(tempJR,NORM_2,&Axdiff_norm);CHKERRQ(ierr); 796285c0a3SHansol Suh ierr = VecNorm(tempJR2,NORM_2,&yhatdiff_norm);CHKERRQ(ierr); 806285c0a3SHansol Suh ierr = VecDot(tempJR,tempJR2,&Axyhat);CHKERRQ(ierr); 816285c0a3SHansol Suh 826285c0a3SHansol Suh ierr = VecWAXPY(tempJR,-1.,am->Bz0,am->Bz);CHKERRQ(ierr); 836285c0a3SHansol Suh ierr = VecWAXPY(tempJR2,-1.,am->y,am->y0);CHKERRQ(ierr); 846285c0a3SHansol Suh ierr = VecNorm(tempJR,NORM_2,&Bzdiff_norm);CHKERRQ(ierr); 856285c0a3SHansol Suh ierr = VecNorm(tempJR2,NORM_2,&ydiff_norm);CHKERRQ(ierr); 866285c0a3SHansol Suh ierr = VecDot(tempJR,tempJR2,&Bzy);CHKERRQ(ierr); 876285c0a3SHansol Suh 886285c0a3SHansol Suh if (Axyhat > am->orthval*Axdiff_norm*yhatdiff_norm + am->mueps) { 896285c0a3SHansol Suh hflag = PETSC_TRUE; 906285c0a3SHansol Suh a_sd = PetscSqr(yhatdiff_norm)/Axyhat; /* alphaSD */ 916285c0a3SHansol Suh a_mg = Axyhat/PetscSqr(Axdiff_norm); /* alphaMG */ 926285c0a3SHansol Suh a_k = (a_mg/a_sd) > 0.5 ? a_mg : a_sd - 0.5*a_mg; 936285c0a3SHansol Suh } 946285c0a3SHansol Suh if (Bzy > am->orthval*Bzdiff_norm*ydiff_norm + am->mueps) { 956285c0a3SHansol Suh gflag = PETSC_TRUE; 966285c0a3SHansol Suh b_sd = PetscSqr(ydiff_norm)/Bzy; /* betaSD */ 976285c0a3SHansol Suh b_mg = Bzy/PetscSqr(Bzdiff_norm); /* betaMG */ 986285c0a3SHansol Suh b_k = (b_mg/b_sd) > 0.5 ? b_mg : b_sd - 0.5*b_mg; 996285c0a3SHansol Suh } 1006285c0a3SHansol Suh am->muold = am->mu; 1016285c0a3SHansol Suh if (gflag && hflag) { 1026285c0a3SHansol Suh am->mu = PetscSqrtReal(a_k*b_k); 1036285c0a3SHansol Suh } else if (hflag) { 1046285c0a3SHansol Suh am->mu = a_k; 1056285c0a3SHansol Suh } else if (gflag) { 1066285c0a3SHansol Suh am->mu = b_k; 1076285c0a3SHansol Suh } 1086285c0a3SHansol Suh if (am->mu > am->muold) { 1096285c0a3SHansol Suh am->mu = am->muold; 1106285c0a3SHansol Suh } 1116285c0a3SHansol Suh if (am->mu < am->mumin) { 1126285c0a3SHansol Suh am->mu = am->mumin; 1136285c0a3SHansol Suh } 1146285c0a3SHansol Suh PetscFunctionReturn(0); 1156285c0a3SHansol Suh } 1166285c0a3SHansol Suh 1177f5c9be9SBarry Smith static PetscErrorCode TaoADMMSetRegularizerType_ADMM(Tao tao, TaoADMMRegularizerType type) 1186285c0a3SHansol Suh { 1196285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 1206285c0a3SHansol Suh 1216285c0a3SHansol Suh PetscFunctionBegin; 1226285c0a3SHansol Suh am->regswitch = type; 1236285c0a3SHansol Suh PetscFunctionReturn(0); 1246285c0a3SHansol Suh } 1256285c0a3SHansol Suh 1267f5c9be9SBarry Smith static PetscErrorCode TaoADMMGetRegularizerType_ADMM(Tao tao, TaoADMMRegularizerType *type) 1276285c0a3SHansol Suh { 1286285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 1296285c0a3SHansol Suh 1306285c0a3SHansol Suh PetscFunctionBegin; 1316285c0a3SHansol Suh *type = am->regswitch; 1326285c0a3SHansol Suh PetscFunctionReturn(0); 1336285c0a3SHansol Suh } 1346285c0a3SHansol Suh 1357f5c9be9SBarry Smith static PetscErrorCode TaoADMMSetUpdateType_ADMM(Tao tao, TaoADMMUpdateType type) 1366285c0a3SHansol Suh { 1376285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 1386285c0a3SHansol Suh 1396285c0a3SHansol Suh PetscFunctionBegin; 1406285c0a3SHansol Suh am->update = type; 1416285c0a3SHansol Suh PetscFunctionReturn(0); 1426285c0a3SHansol Suh } 1436285c0a3SHansol Suh 1447f5c9be9SBarry Smith static PetscErrorCode TaoADMMGetUpdateType_ADMM(Tao tao, TaoADMMUpdateType *type) 1456285c0a3SHansol Suh { 1466285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 1476285c0a3SHansol Suh 1486285c0a3SHansol Suh PetscFunctionBegin; 1496285c0a3SHansol Suh *type = am->update; 1506285c0a3SHansol Suh PetscFunctionReturn(0); 1516285c0a3SHansol Suh } 1526285c0a3SHansol Suh 1536285c0a3SHansol Suh /* This routine updates Jacobians with new x,z vectors, 1546285c0a3SHansol Suh * and then updates Ax and Bz vectors, then computes updated residual vector*/ 1556285c0a3SHansol Suh static PetscErrorCode ADMMUpdateConstraintResidualVector(Tao tao, Vec x, Vec z, Vec Ax, Vec Bz, Vec residual) 1566285c0a3SHansol Suh { 1576285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 1586285c0a3SHansol Suh Tao mis,reg; 1596285c0a3SHansol Suh PetscErrorCode ierr; 1606285c0a3SHansol Suh 1616285c0a3SHansol Suh PetscFunctionBegin; 1626285c0a3SHansol Suh mis = am->subsolverX; 1636285c0a3SHansol Suh reg = am->subsolverZ; 1646285c0a3SHansol Suh ierr = TaoComputeJacobianEquality(mis, x, mis->jacobian_equality, mis->jacobian_equality_pre);CHKERRQ(ierr); 1656285c0a3SHansol Suh ierr = MatMult(mis->jacobian_equality,x,Ax);CHKERRQ(ierr); 1666285c0a3SHansol Suh ierr = TaoComputeJacobianEquality(reg, z, reg->jacobian_equality, reg->jacobian_equality_pre);CHKERRQ(ierr); 1676285c0a3SHansol Suh ierr = MatMult(reg->jacobian_equality,z,Bz);CHKERRQ(ierr); 1686285c0a3SHansol Suh 1696285c0a3SHansol Suh ierr = VecWAXPY(residual,1.,Bz,Ax);CHKERRQ(ierr); 1706285c0a3SHansol Suh if (am->constraint != NULL) { 1716285c0a3SHansol Suh ierr = VecAXPY(residual,-1.,am->constraint);CHKERRQ(ierr); 1726285c0a3SHansol Suh } 1736285c0a3SHansol Suh PetscFunctionReturn(0); 1746285c0a3SHansol Suh } 1756285c0a3SHansol Suh 1766285c0a3SHansol Suh /* Updates Augmented Lagrangians to given routines * 1776285c0a3SHansol Suh * For subsolverX, routine needs to be ComputeObjectiveAndGraidnet 1786285c0a3SHansol Suh * Separate Objective and Gradient routines are not supported. */ 1796285c0a3SHansol Suh static PetscErrorCode SubObjGradUpdate(Tao tao, Vec x, PetscReal *f, Vec g, void *ptr) 1806285c0a3SHansol Suh { 1816285c0a3SHansol Suh Tao parent = (Tao)ptr; 1826285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)parent->data; 1836285c0a3SHansol Suh PetscErrorCode ierr; 1846285c0a3SHansol Suh PetscReal temp,temp2; 1856285c0a3SHansol Suh Vec tempJR; 1866285c0a3SHansol Suh 1876285c0a3SHansol Suh PetscFunctionBegin; 1886285c0a3SHansol Suh tempJR = am->workJacobianRight; 1896285c0a3SHansol Suh ierr = ADMMUpdateConstraintResidualVector(parent, x, am->subsolverZ->solution, am->Ax, am->Bz, am->residual);CHKERRQ(ierr); 1906285c0a3SHansol Suh ierr = (*am->ops->misfitobjgrad)(am->subsolverX,x,f,g,am->misfitobjgradP);CHKERRQ(ierr); 1916285c0a3SHansol Suh 1926285c0a3SHansol Suh am->last_misfit_val = *f; 1936285c0a3SHansol Suh /* Objective Add + yT(Ax+Bz-c) + mu/2*||Ax+Bz-c||_2^2 */ 1946285c0a3SHansol Suh ierr = VecTDot(am->residual,am->y,&temp);CHKERRQ(ierr); 1956285c0a3SHansol Suh ierr = VecTDot(am->residual,am->residual,&temp2);CHKERRQ(ierr); 1966285c0a3SHansol Suh *f += temp + (am->mu/2)*temp2; 1976285c0a3SHansol Suh 1986285c0a3SHansol Suh /* Gradient. Add + mu*AT(Ax+Bz-c) + yTA*/ 1996285c0a3SHansol Suh ierr = MatMultTranspose(tao->jacobian_equality,am->residual,tempJR);CHKERRQ(ierr); 2006285c0a3SHansol Suh ierr = VecAXPY(g,am->mu,tempJR);CHKERRQ(ierr); 2016285c0a3SHansol Suh ierr = MatMultTranspose(tao->jacobian_equality,am->y,tempJR);CHKERRQ(ierr); 2026285c0a3SHansol Suh ierr = VecAXPY(g,1.,tempJR);CHKERRQ(ierr); 2036285c0a3SHansol Suh PetscFunctionReturn(0); 2046285c0a3SHansol Suh } 2056285c0a3SHansol Suh 2066285c0a3SHansol Suh /* Updates Augmented Lagrangians to given routines 2076285c0a3SHansol Suh * For subsolverZ, routine needs to be ComputeObjectiveAndGraidnet 2086285c0a3SHansol Suh * Separate Objective and Gradient routines are not supported. */ 2096285c0a3SHansol Suh static PetscErrorCode RegObjGradUpdate(Tao tao, Vec z, PetscReal *f, Vec g, void *ptr) 2106285c0a3SHansol Suh { 2116285c0a3SHansol Suh Tao parent = (Tao)ptr; 2126285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)parent->data; 2136285c0a3SHansol Suh PetscErrorCode ierr; 2146285c0a3SHansol Suh PetscReal temp,temp2; 2156285c0a3SHansol Suh Vec tempJR; 2166285c0a3SHansol Suh 2176285c0a3SHansol Suh PetscFunctionBegin; 2186285c0a3SHansol Suh tempJR = am->workJacobianRight; 2196285c0a3SHansol Suh ierr = ADMMUpdateConstraintResidualVector(parent, am->subsolverX->solution, z, am->Ax, am->Bz, am->residual);CHKERRQ(ierr); 2206285c0a3SHansol Suh ierr = (*am->ops->regobjgrad)(am->subsolverZ,z,f,g,am->regobjgradP);CHKERRQ(ierr); 2216285c0a3SHansol Suh am->last_reg_val= *f; 2226285c0a3SHansol Suh /* Objective Add + yT(Ax+Bz-c) + mu/2*||Ax+Bz-c||_2^2 */ 2236285c0a3SHansol Suh ierr = VecTDot(am->residual,am->y,&temp);CHKERRQ(ierr); 2246285c0a3SHansol Suh ierr = VecTDot(am->residual,am->residual,&temp2);CHKERRQ(ierr); 2256285c0a3SHansol Suh *f += temp + (am->mu/2)*temp2; 2266285c0a3SHansol Suh 2276285c0a3SHansol Suh /* Gradient. Add + mu*BT(Ax+Bz-c) + yTB*/ 2286285c0a3SHansol Suh ierr = MatMultTranspose(am->subsolverZ->jacobian_equality,am->residual,tempJR);CHKERRQ(ierr); 2296285c0a3SHansol Suh ierr = VecAXPY(g,am->mu,tempJR);CHKERRQ(ierr); 2306285c0a3SHansol Suh ierr = MatMultTranspose(am->subsolverZ->jacobian_equality,am->y,tempJR);CHKERRQ(ierr); 2316285c0a3SHansol Suh ierr = VecAXPY(g,1.,tempJR);CHKERRQ(ierr); 2326285c0a3SHansol Suh PetscFunctionReturn(0); 2336285c0a3SHansol Suh } 2346285c0a3SHansol Suh 2356285c0a3SHansol Suh /* Computes epsilon padded L1 norm lambda*sum(sqrt(x^2+eps^2)-eps */ 2366285c0a3SHansol Suh static PetscErrorCode ADMML1EpsilonNorm(Tao tao, Vec x, PetscReal eps, PetscReal *norm) 2376285c0a3SHansol Suh { 2386285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 2396285c0a3SHansol Suh PetscInt N; 2406285c0a3SHansol Suh PetscErrorCode ierr; 2416285c0a3SHansol Suh 2426285c0a3SHansol Suh PetscFunctionBegin; 2436285c0a3SHansol Suh ierr = VecGetSize(am->workLeft,&N);CHKERRQ(ierr); 2446285c0a3SHansol Suh ierr = VecPointwiseMult(am->workLeft,x,x);CHKERRQ(ierr); 2456285c0a3SHansol Suh ierr = VecShift(am->workLeft,am->l1epsilon*am->l1epsilon);CHKERRQ(ierr); 2466285c0a3SHansol Suh ierr = VecSqrtAbs(am->workLeft);CHKERRQ(ierr); 2476285c0a3SHansol Suh ierr = VecSum(am->workLeft,norm);CHKERRQ(ierr); 2486285c0a3SHansol Suh *norm += N*am->l1epsilon; 2496285c0a3SHansol Suh *norm *= am->lambda; 2506285c0a3SHansol Suh PetscFunctionReturn(0); 2516285c0a3SHansol Suh } 2526285c0a3SHansol Suh 2536285c0a3SHansol Suh static PetscErrorCode ADMMInternalHessianUpdate(Mat H, Mat Constraint, PetscBool Identity, void *ptr) 2546285c0a3SHansol Suh { 2556285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)ptr; 2566285c0a3SHansol Suh PetscErrorCode ierr; 2576285c0a3SHansol Suh 2586285c0a3SHansol Suh PetscFunctionBegin; 2596285c0a3SHansol Suh switch (am->update) { 2606285c0a3SHansol Suh case (TAO_ADMM_UPDATE_BASIC): 2616285c0a3SHansol Suh break; 2626285c0a3SHansol Suh case (TAO_ADMM_UPDATE_ADAPTIVE): 2636285c0a3SHansol Suh case (TAO_ADMM_UPDATE_ADAPTIVE_RELAXED): 2646285c0a3SHansol Suh if (H && (am->muold != am->mu)) { 2656285c0a3SHansol Suh if (!Identity) { 2666285c0a3SHansol Suh ierr = MatAXPY(H,am->mu-am->muold,Constraint,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 2676285c0a3SHansol Suh } else { 2686285c0a3SHansol Suh ierr = MatShift(H,am->mu-am->muold);CHKERRQ(ierr); 2696285c0a3SHansol Suh } 2706285c0a3SHansol Suh } 2716285c0a3SHansol Suh break; 2726285c0a3SHansol Suh } 2736285c0a3SHansol Suh PetscFunctionReturn(0); 2746285c0a3SHansol Suh } 2756285c0a3SHansol Suh 2766285c0a3SHansol Suh /* Updates Hessian - adds second derivative of augmented Lagrangian 2776285c0a3SHansol Suh * H \gets H + \rho*ATA 2786285c0a3SHansol Suh * Here, \rho does not change in TAO_ADMM_UPDATE_BASIC - thus no-op 2796285c0a3SHansol Suh * For ADAPTAIVE,ADAPTIVE_RELAXED, 2806285c0a3SHansol Suh * H \gets H + (\rho-\rhoold)*ATA 2816285c0a3SHansol Suh * Here, we assume that A is linear constraint i.e., doesnt change. 2826285c0a3SHansol Suh * Thus, for both ADAPTIVE, and RELAXED, ATA matrix is pre-set (except for A=I (null case)) see TaoSetUp_ADMM */ 2836285c0a3SHansol Suh static PetscErrorCode SubHessianUpdate(Tao tao, Vec x, Mat H, Mat Hpre, void *ptr) 2846285c0a3SHansol Suh { 2856285c0a3SHansol Suh Tao parent = (Tao)ptr; 2866285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)parent->data; 2876285c0a3SHansol Suh PetscErrorCode ierr; 2886285c0a3SHansol Suh 2896285c0a3SHansol Suh PetscFunctionBegin; 2906285c0a3SHansol Suh if (am->Hxchange) { 2916285c0a3SHansol Suh /* Case where Hessian gets updated with respect to x vector input. */ 2926285c0a3SHansol Suh ierr = (*am->ops->misfithess)(am->subsolverX,x,H,Hpre,am->misfithessP);CHKERRQ(ierr); 2936285c0a3SHansol Suh ierr = ADMMInternalHessianUpdate(am->subsolverX->hessian,am->ATA,am->xJI,am);CHKERRQ(ierr); 2946285c0a3SHansol Suh } else if (am->Hxbool) { 2956285c0a3SHansol Suh /* Hessian doesn't get updated. H(x) = c */ 2966285c0a3SHansol Suh /* Update Lagrangian only only per TAO call */ 2976285c0a3SHansol Suh ierr = ADMMInternalHessianUpdate(am->subsolverX->hessian,am->ATA,am->xJI,am);CHKERRQ(ierr); 2986285c0a3SHansol Suh am->Hxbool = PETSC_FALSE; 2996285c0a3SHansol Suh } 3006285c0a3SHansol Suh PetscFunctionReturn(0); 3016285c0a3SHansol Suh } 3026285c0a3SHansol Suh 3036285c0a3SHansol Suh /* Same as SubHessianUpdate, except for B matrix instead of A matrix */ 3046285c0a3SHansol Suh static PetscErrorCode RegHessianUpdate(Tao tao, Vec z, Mat H, Mat Hpre, void *ptr) 3056285c0a3SHansol Suh { 3066285c0a3SHansol Suh Tao parent = (Tao)ptr; 3076285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)parent->data; 3086285c0a3SHansol Suh PetscErrorCode ierr; 3096285c0a3SHansol Suh 3106285c0a3SHansol Suh PetscFunctionBegin; 3116285c0a3SHansol Suh 3126285c0a3SHansol Suh if (am->Hzchange) { 3136285c0a3SHansol Suh /* Case where Hessian gets updated with respect to x vector input. */ 3146285c0a3SHansol Suh ierr = (*am->ops->reghess)(am->subsolverZ,z,H,Hpre,am->reghessP);CHKERRQ(ierr); 3156285c0a3SHansol Suh ierr = ADMMInternalHessianUpdate(am->subsolverZ->hessian,am->BTB,am->zJI,am);CHKERRQ(ierr); 3166285c0a3SHansol Suh } else if (am->Hzbool) { 3176285c0a3SHansol Suh /* Hessian doesn't get updated. H(x) = c */ 3186285c0a3SHansol Suh /* Update Lagrangian only only per TAO call */ 3196285c0a3SHansol Suh ierr = ADMMInternalHessianUpdate(am->subsolverZ->hessian,am->BTB,am->zJI,am);CHKERRQ(ierr); 3206285c0a3SHansol Suh am->Hzbool = PETSC_FALSE; 3216285c0a3SHansol Suh } 3226285c0a3SHansol Suh PetscFunctionReturn(0); 3236285c0a3SHansol Suh } 3246285c0a3SHansol Suh 3256285c0a3SHansol Suh /* Shell Matrix routine for A matrix. 3266285c0a3SHansol Suh * This gets used when user puts NULL for 3276285c0a3SHansol Suh * TaoSetJacobianEqualityRoutine(tao, NULL,NULL, ...) 3286285c0a3SHansol Suh * Essentially sets A=I*/ 3296285c0a3SHansol Suh static PetscErrorCode JacobianIdentity(Mat mat,Vec in,Vec out) 3306285c0a3SHansol Suh { 3316285c0a3SHansol Suh PetscErrorCode ierr; 3326285c0a3SHansol Suh 3336285c0a3SHansol Suh PetscFunctionBegin; 3346285c0a3SHansol Suh ierr = VecCopy(in,out);CHKERRQ(ierr); 3356285c0a3SHansol Suh PetscFunctionReturn(0); 3366285c0a3SHansol Suh } 3376285c0a3SHansol Suh 3386285c0a3SHansol Suh /* Shell Matrix routine for B matrix. 3396285c0a3SHansol Suh * This gets used when user puts NULL for 3406285c0a3SHansol Suh * TaoADMMSetRegularizerConstraintJacobian(tao, NULL,NULL, ...) 3416285c0a3SHansol Suh * Sets B=-I */ 3426285c0a3SHansol Suh static PetscErrorCode JacobianIdentityB(Mat mat,Vec in,Vec out) 3436285c0a3SHansol Suh { 3446285c0a3SHansol Suh PetscErrorCode ierr; 3456285c0a3SHansol Suh 3466285c0a3SHansol Suh PetscFunctionBegin; 3476285c0a3SHansol Suh ierr = VecCopy(in,out);CHKERRQ(ierr); 3486285c0a3SHansol Suh ierr = VecScale(out,-1.);CHKERRQ(ierr); 3496285c0a3SHansol Suh PetscFunctionReturn(0); 3506285c0a3SHansol Suh } 3516285c0a3SHansol Suh 3526285c0a3SHansol Suh /* Solve f(x) + g(z) s.t. Ax + Bz = c */ 3536285c0a3SHansol Suh static PetscErrorCode TaoSolve_ADMM(Tao tao) 3546285c0a3SHansol Suh { 3556285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 3566285c0a3SHansol Suh PetscErrorCode ierr; 3576285c0a3SHansol Suh PetscInt N; 3586285c0a3SHansol Suh PetscReal reg_func; 3596285c0a3SHansol Suh PetscBool is_reg_shell; 3606285c0a3SHansol Suh Vec tempL; 3616285c0a3SHansol Suh 3626285c0a3SHansol Suh PetscFunctionBegin; 3636285c0a3SHansol Suh if (am->regswitch != TAO_ADMM_REGULARIZER_SOFT_THRESH) { 3646285c0a3SHansol Suh if (!am->subsolverX->ops->computejacobianequality) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_WRONGSTATE,"Must call TaoADMMSetMisfitConstraintJacobian() first"); 3656285c0a3SHansol Suh if (!am->subsolverZ->ops->computejacobianequality) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_WRONGSTATE,"Must call TaoADMMSetRegularizerConstraintJacobian() first"); 3666285c0a3SHansol Suh if (am->constraint != NULL) { 3676285c0a3SHansol Suh ierr = VecNorm(am->constraint,NORM_2,&am->const_norm);CHKERRQ(ierr); 3686285c0a3SHansol Suh } 3696285c0a3SHansol Suh } 3706285c0a3SHansol Suh tempL = am->workLeft; 3716285c0a3SHansol Suh ierr = VecGetSize(tempL,&N);CHKERRQ(ierr); 3726285c0a3SHansol Suh 3736285c0a3SHansol Suh if (am->Hx && am->ops->misfithess) { 3746285c0a3SHansol Suh ierr = TaoSetHessianRoutine(am->subsolverX, am->Hx, am->Hx, SubHessianUpdate, tao);CHKERRQ(ierr); 3756285c0a3SHansol Suh } 3766285c0a3SHansol Suh 3776285c0a3SHansol Suh if (!am->zJI) { 3786285c0a3SHansol Suh /* Currently, B is assumed to be a linear system, i.e., not getting updated*/ 3796285c0a3SHansol Suh ierr = MatTransposeMatMult(am->JB,am->JB,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&(am->BTB));CHKERRQ(ierr); 3806285c0a3SHansol Suh } 3816285c0a3SHansol Suh if (!am->xJI) { 3826285c0a3SHansol Suh /* Currently, A is assumed to be a linear system, i.e., not getting updated*/ 3836285c0a3SHansol Suh ierr = MatTransposeMatMult(am->subsolverX->jacobian_equality,am->subsolverX->jacobian_equality,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&(am->ATA));CHKERRQ(ierr); 3846285c0a3SHansol Suh } 3856285c0a3SHansol Suh 3866285c0a3SHansol Suh is_reg_shell = PETSC_FALSE; 3876285c0a3SHansol Suh 3886285c0a3SHansol Suh ierr = PetscObjectTypeCompare((PetscObject)am->subsolverZ, TAOSHELL, &is_reg_shell);CHKERRQ(ierr); 3896285c0a3SHansol Suh 3906285c0a3SHansol Suh if (!is_reg_shell) { 3916285c0a3SHansol Suh switch (am->regswitch) { 3926285c0a3SHansol Suh case (TAO_ADMM_REGULARIZER_USER): 3936285c0a3SHansol Suh if (!am->ops->regobjgrad) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_WRONGSTATE,"Must call TaoADMMSetRegularizerObjectiveAndGradientRoutine() first if one wishes to use TAO_ADMM_REGULARIZER_USER with non-TAOSHELL type"); 3946285c0a3SHansol Suh break; 3956285c0a3SHansol Suh case (TAO_ADMM_REGULARIZER_SOFT_THRESH): 3966285c0a3SHansol Suh /* Soft Threshold. */ 3976285c0a3SHansol Suh break; 3986285c0a3SHansol Suh } 3996285c0a3SHansol Suh if (am->ops->regobjgrad) { 4006285c0a3SHansol Suh ierr = TaoSetObjectiveAndGradientRoutine(am->subsolverZ, RegObjGradUpdate, tao);CHKERRQ(ierr); 4016285c0a3SHansol Suh } 4026285c0a3SHansol Suh if (am->Hz && am->ops->reghess) { 4036285c0a3SHansol Suh ierr = TaoSetHessianRoutine(am->subsolverZ, am->Hz, am->Hzpre, RegHessianUpdate, tao);CHKERRQ(ierr); 4046285c0a3SHansol Suh } 4056285c0a3SHansol Suh } 4066285c0a3SHansol Suh 4076285c0a3SHansol Suh switch (am->update) { 4086285c0a3SHansol Suh case TAO_ADMM_UPDATE_BASIC: 4096285c0a3SHansol Suh if (am->subsolverX->hessian) { 4106285c0a3SHansol Suh /* In basic case, Hessian does not get updated w.r.t. to spectral penalty 4116285c0a3SHansol Suh * Here, when A is set, i.e., am->xJI, add mu*ATA to Hessian*/ 4126285c0a3SHansol Suh if (!am->xJI) { 4136285c0a3SHansol Suh ierr = MatAXPY(am->subsolverX->hessian,am->mu,am->ATA,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 4146285c0a3SHansol Suh } else { 4156285c0a3SHansol Suh ierr = MatShift(am->subsolverX->hessian,am->mu);CHKERRQ(ierr); 4166285c0a3SHansol Suh } 4176285c0a3SHansol Suh } 4186285c0a3SHansol Suh if (am->subsolverZ->hessian && am->regswitch == TAO_ADMM_REGULARIZER_USER) { 4196285c0a3SHansol Suh if (am->regswitch == TAO_ADMM_REGULARIZER_USER && !am->zJI) { 4206285c0a3SHansol Suh ierr = MatAXPY(am->subsolverZ->hessian,am->mu,am->BTB,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 4216285c0a3SHansol Suh } else { 4226285c0a3SHansol Suh ierr = MatShift(am->subsolverZ->hessian,am->mu);CHKERRQ(ierr); 4236285c0a3SHansol Suh } 4246285c0a3SHansol Suh } 4256285c0a3SHansol Suh break; 4266285c0a3SHansol Suh case TAO_ADMM_UPDATE_ADAPTIVE: 4276285c0a3SHansol Suh case TAO_ADMM_UPDATE_ADAPTIVE_RELAXED: 4286285c0a3SHansol Suh break; 4296285c0a3SHansol Suh } 4306285c0a3SHansol Suh 4316285c0a3SHansol Suh ierr = PetscCitationsRegister(citation,&cited);CHKERRQ(ierr); 4326285c0a3SHansol Suh tao->reason = TAO_CONTINUE_ITERATING; 4336285c0a3SHansol Suh 4346285c0a3SHansol Suh while (tao->reason == TAO_CONTINUE_ITERATING) { 4356285c0a3SHansol Suh if (tao->ops->update) { 4366285c0a3SHansol Suh ierr = (*tao->ops->update)(tao, tao->niter, tao->user_update);CHKERRQ(ierr); 4376285c0a3SHansol Suh } 4386285c0a3SHansol Suh ierr = VecCopy(am->Bz, am->Bzold);CHKERRQ(ierr); 4396285c0a3SHansol Suh 4406285c0a3SHansol Suh /* x update */ 4416285c0a3SHansol Suh ierr = TaoSolve(am->subsolverX);CHKERRQ(ierr); 4426285c0a3SHansol Suh ierr = TaoComputeJacobianEquality(am->subsolverX, am->subsolverX->solution, am->subsolverX->jacobian_equality, am->subsolverX->jacobian_equality_pre);CHKERRQ(ierr); 4436285c0a3SHansol Suh ierr = MatMult(am->subsolverX->jacobian_equality, am->subsolverX->solution,am->Ax);CHKERRQ(ierr); 4446285c0a3SHansol Suh 4456285c0a3SHansol Suh am->Hxbool = PETSC_TRUE; 4466285c0a3SHansol Suh 4476285c0a3SHansol Suh /* z update */ 4486285c0a3SHansol Suh switch (am->regswitch) { 4496285c0a3SHansol Suh case TAO_ADMM_REGULARIZER_USER: 4506285c0a3SHansol Suh ierr = TaoSolve(am->subsolverZ);CHKERRQ(ierr); 4516285c0a3SHansol Suh break; 4526285c0a3SHansol Suh case TAO_ADMM_REGULARIZER_SOFT_THRESH: 4536285c0a3SHansol Suh /* L1 assumes A,B jacobians are identity nxn matrix */ 4546285c0a3SHansol Suh ierr = VecWAXPY(am->workJacobianRight,1/am->mu,am->y,am->Ax);CHKERRQ(ierr); 4556285c0a3SHansol Suh ierr = TaoSoftThreshold(am->workJacobianRight,-am->lambda/am->mu,am->lambda/am->mu,am->subsolverZ->solution);CHKERRQ(ierr); 4566285c0a3SHansol Suh break; 4576285c0a3SHansol Suh } 4586285c0a3SHansol Suh am->Hzbool = PETSC_TRUE; 4596285c0a3SHansol Suh /* Returns Ax + Bz - c with updated Ax,Bz vectors */ 4606285c0a3SHansol Suh ierr = ADMMUpdateConstraintResidualVector(tao, am->subsolverX->solution, am->subsolverZ->solution, am->Ax, am->Bz, am->residual);CHKERRQ(ierr); 4616285c0a3SHansol Suh /* Dual variable, y += y + mu*(Ax+Bz-c) */ 4626285c0a3SHansol Suh ierr = VecWAXPY(am->y, am->mu, am->residual, am->yold);CHKERRQ(ierr); 4636285c0a3SHansol Suh 4646285c0a3SHansol Suh /* stopping tolerance update */ 4656285c0a3SHansol Suh ierr = TaoADMMToleranceUpdate(tao);CHKERRQ(ierr); 4666285c0a3SHansol Suh 4676285c0a3SHansol Suh /* Updating Spectral Penalty */ 4686285c0a3SHansol Suh switch (am->update) { 4696285c0a3SHansol Suh case TAO_ADMM_UPDATE_BASIC: 4706285c0a3SHansol Suh am->muold = am->mu; 4716285c0a3SHansol Suh break; 4726285c0a3SHansol Suh case TAO_ADMM_UPDATE_ADAPTIVE: 4736285c0a3SHansol Suh case TAO_ADMM_UPDATE_ADAPTIVE_RELAXED: 4746285c0a3SHansol Suh if (tao->niter == 0) { 4756285c0a3SHansol Suh ierr = VecCopy(am->y, am->y0);CHKERRQ(ierr); 4766285c0a3SHansol Suh ierr = VecWAXPY(am->residual, 1., am->Ax, am->Bzold);CHKERRQ(ierr); 4776285c0a3SHansol Suh if (am->constraint) { 4786285c0a3SHansol Suh ierr = VecAXPY(am->residual, -1., am->constraint);CHKERRQ(ierr); 4796285c0a3SHansol Suh } 4806285c0a3SHansol Suh ierr = VecWAXPY(am->yhatold, -am->mu, am->residual, am->yold);CHKERRQ(ierr); 4816285c0a3SHansol Suh ierr = VecCopy(am->Ax, am->Axold);CHKERRQ(ierr); 4826285c0a3SHansol Suh ierr = VecCopy(am->Bz, am->Bz0);CHKERRQ(ierr); 4836285c0a3SHansol Suh am->muold = am->mu; 4846285c0a3SHansol Suh } else if (tao->niter % am->T == 1) { 4856285c0a3SHansol Suh /* we have compute Bzold in a previous iteration, and we computed Ax above */ 4866285c0a3SHansol Suh ierr = VecWAXPY(am->residual, 1., am->Ax, am->Bzold);CHKERRQ(ierr); 4876285c0a3SHansol Suh if (am->constraint) { 4886285c0a3SHansol Suh ierr = VecAXPY(am->residual, -1., am->constraint);CHKERRQ(ierr); 4896285c0a3SHansol Suh } 4906285c0a3SHansol Suh ierr = VecWAXPY(am->yhat, -am->mu, am->residual, am->yold);CHKERRQ(ierr); 4916285c0a3SHansol Suh ierr = AdaptiveADMMPenaltyUpdate(tao);CHKERRQ(ierr); 4926285c0a3SHansol Suh ierr = VecCopy(am->Ax, am->Axold);CHKERRQ(ierr); 4936285c0a3SHansol Suh ierr = VecCopy(am->Bz, am->Bz0);CHKERRQ(ierr); 4946285c0a3SHansol Suh ierr = VecCopy(am->yhat, am->yhatold);CHKERRQ(ierr); 4956285c0a3SHansol Suh ierr = VecCopy(am->y, am->y0);CHKERRQ(ierr); 4966285c0a3SHansol Suh } else { 4976285c0a3SHansol Suh am->muold = am->mu; 4986285c0a3SHansol Suh } 4996285c0a3SHansol Suh break; 5006285c0a3SHansol Suh default: 5016285c0a3SHansol Suh break; 5026285c0a3SHansol Suh } 5036285c0a3SHansol Suh tao->niter++; 5046285c0a3SHansol Suh 5056285c0a3SHansol Suh /* Calculate original function values. misfit part was done in TaoADMMToleranceUpdate*/ 5066285c0a3SHansol Suh switch (am->regswitch) { 5076285c0a3SHansol Suh case TAO_ADMM_REGULARIZER_USER: 5086285c0a3SHansol Suh if (is_reg_shell) { 5096285c0a3SHansol Suh ierr = ADMML1EpsilonNorm(tao,am->subsolverZ->solution,am->l1epsilon,®_func);CHKERRQ(ierr); 5106285c0a3SHansol Suh } else { 5116285c0a3SHansol Suh (*am->ops->regobjgrad)(am->subsolverZ,am->subsolverX->solution,®_func,tempL,am->regobjgradP); 5126285c0a3SHansol Suh } 5136285c0a3SHansol Suh break; 5146285c0a3SHansol Suh case TAO_ADMM_REGULARIZER_SOFT_THRESH: 5156285c0a3SHansol Suh ierr = ADMML1EpsilonNorm(tao,am->subsolverZ->solution,am->l1epsilon,®_func);CHKERRQ(ierr); 5166285c0a3SHansol Suh break; 5176285c0a3SHansol Suh } 5186285c0a3SHansol Suh ierr = VecCopy(am->y,am->yold);CHKERRQ(ierr); 5196285c0a3SHansol Suh ierr = ADMMUpdateConstraintResidualVector(tao, am->subsolverX->solution, am->subsolverZ->solution, am->Ax, am->Bz, am->residual);CHKERRQ(ierr); 5206285c0a3SHansol Suh ierr = VecNorm(am->residual,NORM_2,&am->resnorm);CHKERRQ(ierr); 5216285c0a3SHansol Suh ierr = TaoLogConvergenceHistory(tao,am->last_misfit_val + reg_func,am->dualres,am->resnorm,tao->ksp_its);CHKERRQ(ierr); 5226285c0a3SHansol Suh 5236285c0a3SHansol Suh ierr = TaoMonitor(tao,tao->niter,am->last_misfit_val + reg_func,am->dualres,am->resnorm,1.0);CHKERRQ(ierr); 5246285c0a3SHansol Suh ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr); 5256285c0a3SHansol Suh } 5266285c0a3SHansol Suh /* Update vectors */ 5276285c0a3SHansol Suh ierr = VecCopy(am->subsolverX->solution,tao->solution);CHKERRQ(ierr); 5286285c0a3SHansol Suh ierr = VecCopy(am->subsolverX->gradient,tao->gradient);CHKERRQ(ierr); 5296285c0a3SHansol Suh ierr = PetscObjectCompose((PetscObject)am->subsolverX,"TaoGetADMMParentTao_ADMM", NULL);CHKERRQ(ierr); 5306285c0a3SHansol Suh ierr = PetscObjectCompose((PetscObject)am->subsolverZ,"TaoGetADMMParentTao_ADMM", NULL);CHKERRQ(ierr); 5316285c0a3SHansol Suh ierr = PetscObjectComposeFunction((PetscObject)tao,"TaoADMMSetRegularizerType_C",NULL);CHKERRQ(ierr); 5326285c0a3SHansol Suh ierr = PetscObjectComposeFunction((PetscObject)tao,"TaoADMMGetRegularizerType_C",NULL);CHKERRQ(ierr); 5336285c0a3SHansol Suh ierr = PetscObjectComposeFunction((PetscObject)tao,"TaoADMMSetUpdateType_C",NULL);CHKERRQ(ierr); 5346285c0a3SHansol Suh ierr = PetscObjectComposeFunction((PetscObject)tao,"TaoADMMGetUpdateType_C",NULL);CHKERRQ(ierr); 5356285c0a3SHansol Suh PetscFunctionReturn(0); 5366285c0a3SHansol Suh } 5376285c0a3SHansol Suh 5386285c0a3SHansol Suh static PetscErrorCode TaoSetFromOptions_ADMM(PetscOptionItems *PetscOptionsObject,Tao tao) 5396285c0a3SHansol Suh { 5406285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 5416285c0a3SHansol Suh PetscErrorCode ierr; 5426285c0a3SHansol Suh 5436285c0a3SHansol Suh PetscFunctionBegin; 5446285c0a3SHansol Suh ierr = PetscOptionsHead(PetscOptionsObject,"ADMM problem that solves f(x) in a form of f(x) + g(z) subject to x - z = 0. Norm 1 and 2 are supported. Different subsolver routines can be selected. ");CHKERRQ(ierr); 5456285c0a3SHansol Suh ierr = PetscOptionsReal("-tao_admm_regularizer_coefficient","regularizer constant","",am->lambda,&am->lambda,NULL);CHKERRQ(ierr); 5466285c0a3SHansol Suh ierr = PetscOptionsReal("-tao_admm_spectral_penalty","Constant for Augmented Lagrangian term.","",am->mu,&am->mu,NULL);CHKERRQ(ierr); 5476285c0a3SHansol Suh ierr = PetscOptionsReal("-tao_admm_relaxation_parameter","x relaxation parameter for Z update.","",am->gamma,&am->gamma,NULL);CHKERRQ(ierr); 5486285c0a3SHansol Suh ierr = PetscOptionsReal("-tao_admm_tolerance_update_factor","ADMM dynamic tolerance update factor.","",am->tol,&am->tol,NULL);CHKERRQ(ierr); 5496285c0a3SHansol Suh ierr = PetscOptionsReal("-tao_admm_spectral_penalty_update_factor","ADMM spectral penalty update curvature safeguard value.","",am->orthval,&am->orthval,NULL);CHKERRQ(ierr); 5506285c0a3SHansol Suh ierr = PetscOptionsReal("-tao_admm_minimum_spectral_penalty","Set ADMM minimum spectral penalty.","",am->mumin,&am->mumin,NULL);CHKERRQ(ierr); 5516285c0a3SHansol Suh ierr = PetscOptionsEnum("-tao_admm_dual_update","Lagrangian dual update policy","TaoADMMUpdateType", 5526285c0a3SHansol Suh TaoADMMUpdateTypes,(PetscEnum)am->update,(PetscEnum*)&am->update,NULL);CHKERRQ(ierr); 5536285c0a3SHansol Suh ierr = PetscOptionsEnum("-tao_admm_regularizer_type","ADMM regularizer update rule","TaoADMMRegularizerType", 5546285c0a3SHansol Suh TaoADMMRegularizerTypes,(PetscEnum)am->regswitch,(PetscEnum*)&am->regswitch,NULL);CHKERRQ(ierr); 5556285c0a3SHansol Suh ierr = PetscOptionsTail();CHKERRQ(ierr); 5566285c0a3SHansol Suh ierr = TaoSetFromOptions(am->subsolverX);CHKERRQ(ierr); 5576285c0a3SHansol Suh if (am->regswitch != TAO_ADMM_REGULARIZER_SOFT_THRESH) { 5586285c0a3SHansol Suh ierr = TaoSetFromOptions(am->subsolverZ);CHKERRQ(ierr); 5596285c0a3SHansol Suh } 5606285c0a3SHansol Suh PetscFunctionReturn(0); 5616285c0a3SHansol Suh } 5626285c0a3SHansol Suh 5636285c0a3SHansol Suh static PetscErrorCode TaoView_ADMM(Tao tao,PetscViewer viewer) 5646285c0a3SHansol Suh { 5656285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 5666285c0a3SHansol Suh PetscErrorCode ierr; 5676285c0a3SHansol Suh 5686285c0a3SHansol Suh PetscFunctionBegin; 5696285c0a3SHansol Suh ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 5706285c0a3SHansol Suh ierr = TaoView(am->subsolverX,viewer);CHKERRQ(ierr); 5716285c0a3SHansol Suh ierr = TaoView(am->subsolverZ,viewer);CHKERRQ(ierr); 5726285c0a3SHansol Suh ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 5736285c0a3SHansol Suh PetscFunctionReturn(0); 5746285c0a3SHansol Suh } 5756285c0a3SHansol Suh 5766285c0a3SHansol Suh static PetscErrorCode TaoSetUp_ADMM(Tao tao) 5776285c0a3SHansol Suh { 5786285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 5796285c0a3SHansol Suh PetscErrorCode ierr; 5806285c0a3SHansol Suh PetscInt n,N,M; 5816285c0a3SHansol Suh 5826285c0a3SHansol Suh PetscFunctionBegin; 5836285c0a3SHansol Suh ierr = VecGetLocalSize(tao->solution,&n);CHKERRQ(ierr); 5846285c0a3SHansol Suh ierr = VecGetSize(tao->solution,&N);CHKERRQ(ierr); 5856285c0a3SHansol Suh /* If Jacobian is given as NULL, it means Jacobian is identity matrix with size of solution vector */ 5866285c0a3SHansol Suh if (!am->JB) { 5876285c0a3SHansol Suh am->zJI = PETSC_TRUE; 588*ea13f565SStefano Zampini ierr = MatCreateShell(PetscObjectComm((PetscObject)tao),n,n,PETSC_DETERMINE,PETSC_DETERMINE,NULL,&am->JB);CHKERRQ(ierr); 5896285c0a3SHansol Suh ierr = MatShellSetOperation(am->JB,MATOP_MULT,(void (*)(void))JacobianIdentityB);CHKERRQ(ierr); 5906285c0a3SHansol Suh ierr = MatShellSetOperation(am->JB,MATOP_MULT_TRANSPOSE,(void (*)(void))JacobianIdentityB);CHKERRQ(ierr); 5916285c0a3SHansol Suh ierr = MatShellSetOperation(am->JB,MATOP_TRANSPOSE_MAT_MULT,(void (*)(void))JacobianIdentityB);CHKERRQ(ierr); 5926285c0a3SHansol Suh am->JBpre = am->JB; 5936285c0a3SHansol Suh } 5946285c0a3SHansol Suh if (!am->JA) { 5956285c0a3SHansol Suh am->xJI = PETSC_TRUE; 596*ea13f565SStefano Zampini ierr = MatCreateShell(PetscObjectComm((PetscObject)tao),n,n,PETSC_DETERMINE,PETSC_DETERMINE,NULL,&am->JA);CHKERRQ(ierr); 5976285c0a3SHansol Suh ierr = MatShellSetOperation(am->JA,MATOP_MULT,(void (*)(void))JacobianIdentity);CHKERRQ(ierr); 5986285c0a3SHansol Suh ierr = MatShellSetOperation(am->JA,MATOP_MULT_TRANSPOSE,(void (*)(void))JacobianIdentity);CHKERRQ(ierr); 5996285c0a3SHansol Suh ierr = MatShellSetOperation(am->JA,MATOP_TRANSPOSE_MAT_MULT,(void (*)(void))JacobianIdentity);CHKERRQ(ierr); 6006285c0a3SHansol Suh am->JApre = am->JA; 6016285c0a3SHansol Suh } 6026285c0a3SHansol Suh ierr = MatCreateVecs(am->JA,NULL,&am->Ax);CHKERRQ(ierr); 6036285c0a3SHansol Suh if (!tao->gradient) { 6046285c0a3SHansol Suh ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr); 6056285c0a3SHansol Suh } 6066285c0a3SHansol Suh ierr = TaoSetInitialVector(am->subsolverX, tao->solution);CHKERRQ(ierr); 6076285c0a3SHansol Suh if (!am->z) { 6086285c0a3SHansol Suh ierr = VecDuplicate(tao->solution,&am->z);CHKERRQ(ierr); 6096285c0a3SHansol Suh ierr = VecSet(am->z,0.0);CHKERRQ(ierr); 6106285c0a3SHansol Suh } 6116285c0a3SHansol Suh ierr = TaoSetInitialVector(am->subsolverZ, am->z);CHKERRQ(ierr); 6126285c0a3SHansol Suh if (!am->workLeft) { 6136285c0a3SHansol Suh ierr = VecDuplicate(tao->solution,&am->workLeft);CHKERRQ(ierr); 6146285c0a3SHansol Suh } 6156285c0a3SHansol Suh if (!am->Axold) { 6166285c0a3SHansol Suh ierr = VecDuplicate(am->Ax,&am->Axold);CHKERRQ(ierr); 6176285c0a3SHansol Suh } 6186285c0a3SHansol Suh if (!am->workJacobianRight) { 6196285c0a3SHansol Suh ierr = VecDuplicate(am->Ax,&am->workJacobianRight);CHKERRQ(ierr); 6206285c0a3SHansol Suh } 6216285c0a3SHansol Suh if (!am->workJacobianRight2) { 6226285c0a3SHansol Suh ierr = VecDuplicate(am->Ax,&am->workJacobianRight2);CHKERRQ(ierr); 6236285c0a3SHansol Suh } 6246285c0a3SHansol Suh if (!am->Bz) { 6256285c0a3SHansol Suh ierr = VecDuplicate(am->Ax,&am->Bz);CHKERRQ(ierr); 6266285c0a3SHansol Suh } 6276285c0a3SHansol Suh if (!am->Bzold) { 6286285c0a3SHansol Suh ierr = VecDuplicate(am->Ax,&am->Bzold);CHKERRQ(ierr); 6296285c0a3SHansol Suh } 6306285c0a3SHansol Suh if (!am->Bz0) { 6316285c0a3SHansol Suh ierr = VecDuplicate(am->Ax,&am->Bz0);CHKERRQ(ierr); 6326285c0a3SHansol Suh } 6336285c0a3SHansol Suh if (!am->y) { 6346285c0a3SHansol Suh ierr = VecDuplicate(am->Ax,&am->y);CHKERRQ(ierr); 6356285c0a3SHansol Suh ierr = VecSet(am->y,0.0);CHKERRQ(ierr); 6366285c0a3SHansol Suh } 6376285c0a3SHansol Suh if (!am->yold) { 6386285c0a3SHansol Suh ierr = VecDuplicate(am->Ax,&am->yold);CHKERRQ(ierr); 6396285c0a3SHansol Suh ierr = VecSet(am->yold,0.0);CHKERRQ(ierr); 6406285c0a3SHansol Suh } 6416285c0a3SHansol Suh if (!am->y0) { 6426285c0a3SHansol Suh ierr = VecDuplicate(am->Ax,&am->y0);CHKERRQ(ierr); 6436285c0a3SHansol Suh ierr = VecSet(am->y0,0.0);CHKERRQ(ierr); 6446285c0a3SHansol Suh } 6456285c0a3SHansol Suh if (!am->yhat) { 6466285c0a3SHansol Suh ierr = VecDuplicate(am->Ax,&am->yhat);CHKERRQ(ierr); 6476285c0a3SHansol Suh ierr = VecSet(am->yhat,0.0);CHKERRQ(ierr); 6486285c0a3SHansol Suh } 6496285c0a3SHansol Suh if (!am->yhatold) { 6506285c0a3SHansol Suh ierr = VecDuplicate(am->Ax,&am->yhatold);CHKERRQ(ierr); 6516285c0a3SHansol Suh ierr = VecSet(am->yhatold,0.0);CHKERRQ(ierr); 6526285c0a3SHansol Suh } 6536285c0a3SHansol Suh if (!am->residual) { 6546285c0a3SHansol Suh ierr = VecDuplicate(am->Ax,&am->residual);CHKERRQ(ierr); 6556285c0a3SHansol Suh ierr = VecSet(am->residual,0.0);CHKERRQ(ierr); 6566285c0a3SHansol Suh } 6576285c0a3SHansol Suh if (!am->constraint) { 6586285c0a3SHansol Suh am->constraint = NULL; 6596285c0a3SHansol Suh } else { 6606285c0a3SHansol Suh ierr = VecGetSize(am->constraint,&M);CHKERRQ(ierr); 6616285c0a3SHansol Suh if (M != N) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_WRONGSTATE,"Solution vector and constraint vector must be of same size!");CHKERRQ(ierr); 6626285c0a3SHansol Suh } 6636285c0a3SHansol Suh 6646285c0a3SHansol Suh /* Save changed tao tolerance for adaptive tolerance */ 6656285c0a3SHansol Suh if (tao->gatol_changed) { 6666285c0a3SHansol Suh am->gatol_admm = tao->gatol; 6676285c0a3SHansol Suh } 6686285c0a3SHansol Suh if (tao->catol_changed) { 6696285c0a3SHansol Suh am->catol_admm = tao->catol; 6706285c0a3SHansol Suh } 6716285c0a3SHansol Suh 6726285c0a3SHansol Suh /*Update spectral and dual elements to X subsolver */ 6736285c0a3SHansol Suh ierr = TaoSetObjectiveAndGradientRoutine(am->subsolverX, SubObjGradUpdate, tao);CHKERRQ(ierr); 6746285c0a3SHansol Suh ierr = TaoSetJacobianEqualityRoutine(am->subsolverX,am->JA,am->JApre, am->ops->misfitjac, am->misfitjacobianP);CHKERRQ(ierr); 6756285c0a3SHansol Suh ierr = TaoSetJacobianEqualityRoutine(am->subsolverZ,am->JB,am->JBpre, am->ops->regjac, am->regjacobianP);CHKERRQ(ierr); 6766285c0a3SHansol Suh PetscFunctionReturn(0); 6776285c0a3SHansol Suh } 6786285c0a3SHansol Suh 6796285c0a3SHansol Suh static PetscErrorCode TaoDestroy_ADMM(Tao tao) 6806285c0a3SHansol Suh { 6816285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 6826285c0a3SHansol Suh PetscErrorCode ierr; 6836285c0a3SHansol Suh 6846285c0a3SHansol Suh PetscFunctionBegin; 6856285c0a3SHansol Suh ierr = VecDestroy(&am->z);CHKERRQ(ierr); 6866285c0a3SHansol Suh ierr = VecDestroy(&am->Ax);CHKERRQ(ierr); 6876285c0a3SHansol Suh ierr = VecDestroy(&am->Axold);CHKERRQ(ierr); 6886285c0a3SHansol Suh ierr = VecDestroy(&am->Bz);CHKERRQ(ierr); 6896285c0a3SHansol Suh ierr = VecDestroy(&am->Bzold);CHKERRQ(ierr); 6906285c0a3SHansol Suh ierr = VecDestroy(&am->Bz0);CHKERRQ(ierr); 6916285c0a3SHansol Suh ierr = VecDestroy(&am->residual);CHKERRQ(ierr); 6926285c0a3SHansol Suh ierr = VecDestroy(&am->y);CHKERRQ(ierr); 6936285c0a3SHansol Suh ierr = VecDestroy(&am->yold);CHKERRQ(ierr); 6946285c0a3SHansol Suh ierr = VecDestroy(&am->y0);CHKERRQ(ierr); 6956285c0a3SHansol Suh ierr = VecDestroy(&am->yhat);CHKERRQ(ierr); 6966285c0a3SHansol Suh ierr = VecDestroy(&am->yhatold);CHKERRQ(ierr); 6976285c0a3SHansol Suh ierr = VecDestroy(&am->workLeft);CHKERRQ(ierr); 6986285c0a3SHansol Suh ierr = VecDestroy(&am->workJacobianRight);CHKERRQ(ierr); 6996285c0a3SHansol Suh ierr = VecDestroy(&am->workJacobianRight2);CHKERRQ(ierr); 7006285c0a3SHansol Suh 7016285c0a3SHansol Suh ierr = MatDestroy(&am->JA);CHKERRQ(ierr); 7026285c0a3SHansol Suh ierr = MatDestroy(&am->JB);CHKERRQ(ierr); 7036285c0a3SHansol Suh if (!am->xJI) { 7046285c0a3SHansol Suh ierr = MatDestroy(&am->JApre);CHKERRQ(ierr); 7056285c0a3SHansol Suh } 7066285c0a3SHansol Suh if (!am->zJI) { 7076285c0a3SHansol Suh ierr = MatDestroy(&am->JBpre);CHKERRQ(ierr); 7086285c0a3SHansol Suh } 7096285c0a3SHansol Suh if (am->Hx) { 7106285c0a3SHansol Suh ierr = MatDestroy(&am->Hx);CHKERRQ(ierr); 7116285c0a3SHansol Suh ierr = MatDestroy(&am->Hxpre);CHKERRQ(ierr); 7126285c0a3SHansol Suh } 7136285c0a3SHansol Suh if (am->Hz) { 7146285c0a3SHansol Suh ierr = MatDestroy(&am->Hz);CHKERRQ(ierr); 7156285c0a3SHansol Suh ierr = MatDestroy(&am->Hzpre);CHKERRQ(ierr); 7166285c0a3SHansol Suh } 7176285c0a3SHansol Suh ierr = MatDestroy(&am->ATA);CHKERRQ(ierr); 7186285c0a3SHansol Suh ierr = MatDestroy(&am->BTB);CHKERRQ(ierr); 7196285c0a3SHansol Suh ierr = TaoDestroy(&am->subsolverX);CHKERRQ(ierr); 7206285c0a3SHansol Suh ierr = TaoDestroy(&am->subsolverZ);CHKERRQ(ierr); 7216285c0a3SHansol Suh am->parent = NULL; 7226285c0a3SHansol Suh ierr = PetscFree(tao->data);CHKERRQ(ierr); 7236285c0a3SHansol Suh PetscFunctionReturn(0); 7246285c0a3SHansol Suh } 7256285c0a3SHansol Suh 7266285c0a3SHansol Suh /*MC 7276285c0a3SHansol Suh 7286285c0a3SHansol Suh TAOADMM - Alternating direction method of multipliers method fo solving linear problems with 7296285c0a3SHansol Suh constraints. in a min_x f(x) + g(z) s.t. Ax+Bz=c. 7306285c0a3SHansol Suh This algorithm employs two sub Tao solvers, of which type can be specificed 7316285c0a3SHansol Suh by the user. User need to provide ObjectiveAndGradient routine, and/or HessianRoutine for both subsolvers. 7326285c0a3SHansol Suh Hessians can be given boolean flag determining whether they change with respect to a input vector. This can be set via 7336285c0a3SHansol Suh TaoADMMSet{Misfit,Regularizer}HessianChangeStatus. 7346285c0a3SHansol Suh Second subsolver does support TAOSHELL. It should be noted that L1-norm is used for objective value for TAOSHELL type. 7356285c0a3SHansol Suh There is option to set regularizer option, and currently soft-threshold is implemented. For spectral penalty update, 7366285c0a3SHansol Suh currently there are baisc option and adaptive option. 7376285c0a3SHansol Suh Constraint is set at Ax+Bz=c, and A and B can be set with TaoADMMSet{Misfit,Regularizer}ConstraintJacobian. 7386285c0a3SHansol Suh c can be set with TaoADMMSetConstraintVectorRHS. 7396285c0a3SHansol Suh The user can also provide regularizer weight for second subsolver. 7406285c0a3SHansol Suh 7416285c0a3SHansol Suh References: 7426285c0a3SHansol Suh . 1. - Xu, Zheng and Figueiredo, Mario A. T. and Yuan, Xiaoming and Studer, Christoph and Goldstein, Tom 7436285c0a3SHansol Suh "Adaptive Relaxed ADMM: Convergence Theory and Practical Implementation" 7446285c0a3SHansol Suh The IEEE Conference on Computer Vision and Pattern Recognition (CVPR), July, 2017. 7456285c0a3SHansol Suh 7466285c0a3SHansol Suh Options Database Keys: 7476285c0a3SHansol Suh + -tao_admm_regularizer_coefficient - regularizer constant (default 1.e-6) 7486285c0a3SHansol Suh . -tao_admm_spectral_penalty - Constant for Augmented Lagrangian term (default 1.) 7496285c0a3SHansol Suh . -tao_admm_relaxation_parameter - relaxation parameter for Z update (default 1.) 7506285c0a3SHansol Suh . -tao_admm_tolerance_update_factor - ADMM dynamic tolerance update factor (default 1.e-12) 7516285c0a3SHansol Suh . -tao_admm_spectral_penalty_update_factor - ADMM spectral penalty update curvature safeguard value (default 0.2) 7526285c0a3SHansol Suh . -tao_admm_minimum_spectral_penalty - Set ADMM minimum spectral penalty (default 0) 7536285c0a3SHansol Suh . -tao_admm_dual_update - Lagrangian dual update policy ("basic","adaptive","adaptive-relaxed") (default "basic") 7546285c0a3SHansol Suh - -tao_admm_regularizer_type - ADMM regularizer update rule ("user","soft-threshold") (default "soft-threshold") 7556285c0a3SHansol Suh 7566285c0a3SHansol Suh Level: beginner 7576285c0a3SHansol Suh 7586285c0a3SHansol Suh .seealso: TaoADMMSetMisfitHessianChangeStatus(), TaoADMMSetRegHessianChangeStatus(), TaoADMMGetSpectralPenalty(), 7596285c0a3SHansol Suh TaoADMMGetMisfitSubsolver(), TaoADMMGetRegularizationSubsolver(), TaoADMMSetConstraintVectorRHS(), 7606285c0a3SHansol Suh TaoADMMSetMinimumSpectralPenalty(), TaoADMMSetRegularizerCoefficient(), 7616285c0a3SHansol Suh TaoADMMSetRegularizerConstraintJacobian(), TaoADMMSetMisfitConstraintJacobian(), 7626285c0a3SHansol Suh TaoADMMSetMisfitObjectiveAndGradientRoutine(), TaoADMMSetMisfitHessianRoutine(), 7636285c0a3SHansol Suh TaoADMMSetRegularizerObjectiveAndGradientRoutine(), TaoADMMSetRegularizerHessianRoutine(), 7646285c0a3SHansol Suh TaoGetADMMParentTao(), TaoADMMGetDualVector(), TaoADMMSetRegularizerType(), 7656285c0a3SHansol Suh TaoADMMGetRegularizerType(), TaoADMMSetUpdateType(), TaoADMMGetUpdateType() 7666285c0a3SHansol Suh M*/ 7676285c0a3SHansol Suh 7686285c0a3SHansol Suh PETSC_EXTERN PetscErrorCode TaoCreate_ADMM(Tao tao) 7696285c0a3SHansol Suh { 7706285c0a3SHansol Suh TAO_ADMM *am; 7716285c0a3SHansol Suh PetscErrorCode ierr; 7726285c0a3SHansol Suh 7736285c0a3SHansol Suh PetscFunctionBegin; 7746285c0a3SHansol Suh ierr = PetscNewLog(tao,&am);CHKERRQ(ierr); 7756285c0a3SHansol Suh 7766285c0a3SHansol Suh tao->ops->destroy = TaoDestroy_ADMM; 7776285c0a3SHansol Suh tao->ops->setup = TaoSetUp_ADMM; 7786285c0a3SHansol Suh tao->ops->setfromoptions = TaoSetFromOptions_ADMM; 7796285c0a3SHansol Suh tao->ops->view = TaoView_ADMM; 7806285c0a3SHansol Suh tao->ops->solve = TaoSolve_ADMM; 7816285c0a3SHansol Suh 7826285c0a3SHansol Suh tao->data = (void*)am; 7836285c0a3SHansol Suh am->l1epsilon = 1e-6; 7846285c0a3SHansol Suh am->lambda = 1e-4; 7856285c0a3SHansol Suh am->mu = 1.; 7866285c0a3SHansol Suh am->muold = 0.; 7876285c0a3SHansol Suh am->mueps = PETSC_MACHINE_EPSILON; 7886285c0a3SHansol Suh am->mumin = 0.; 7896285c0a3SHansol Suh am->orthval = 0.2; 7906285c0a3SHansol Suh am->T = 2; 7916285c0a3SHansol Suh am->parent = tao; 7926285c0a3SHansol Suh am->update = TAO_ADMM_UPDATE_BASIC; 7936285c0a3SHansol Suh am->regswitch = TAO_ADMM_REGULARIZER_SOFT_THRESH; 7946285c0a3SHansol Suh am->tol = PETSC_SMALL; 7956285c0a3SHansol Suh am->const_norm = 0; 7966285c0a3SHansol Suh am->resnorm = 0; 7976285c0a3SHansol Suh am->dualres = 0; 79883c8fe1dSLisandro Dalcin am->ops->regobjgrad = NULL; 79983c8fe1dSLisandro Dalcin am->ops->reghess = NULL; 8006285c0a3SHansol Suh am->gamma = 1; 80183c8fe1dSLisandro Dalcin am->regobjgradP = NULL; 80283c8fe1dSLisandro Dalcin am->reghessP = NULL; 8036285c0a3SHansol Suh am->gatol_admm = 1e-8; 8046285c0a3SHansol Suh am->catol_admm = 0; 8056285c0a3SHansol Suh am->Hxchange = PETSC_TRUE; 8066285c0a3SHansol Suh am->Hzchange = PETSC_TRUE; 8076285c0a3SHansol Suh am->Hzbool = PETSC_TRUE; 8086285c0a3SHansol Suh am->Hxbool = PETSC_TRUE; 8096285c0a3SHansol Suh 8106285c0a3SHansol Suh ierr = TaoCreate(PetscObjectComm((PetscObject)tao),&am->subsolverX);CHKERRQ(ierr); 8116285c0a3SHansol Suh ierr = TaoSetOptionsPrefix(am->subsolverX,"misfit_");CHKERRQ(ierr); 8126285c0a3SHansol Suh ierr = PetscObjectIncrementTabLevel((PetscObject)am->subsolverX,(PetscObject)tao,1);CHKERRQ(ierr); 8136285c0a3SHansol Suh ierr = TaoCreate(PetscObjectComm((PetscObject)tao),&am->subsolverZ);CHKERRQ(ierr); 8146285c0a3SHansol Suh ierr = TaoSetOptionsPrefix(am->subsolverZ,"reg_");CHKERRQ(ierr); 8156285c0a3SHansol Suh ierr = PetscObjectIncrementTabLevel((PetscObject)am->subsolverZ,(PetscObject)tao,1);CHKERRQ(ierr); 8166285c0a3SHansol Suh 8176285c0a3SHansol Suh ierr = TaoSetType(am->subsolverX,TAONLS);CHKERRQ(ierr); 8186285c0a3SHansol Suh ierr = TaoSetType(am->subsolverZ,TAONLS);CHKERRQ(ierr); 8196285c0a3SHansol Suh ierr = PetscObjectCompose((PetscObject)am->subsolverX,"TaoGetADMMParentTao_ADMM", (PetscObject) tao);CHKERRQ(ierr); 8206285c0a3SHansol Suh ierr = PetscObjectCompose((PetscObject)am->subsolverZ,"TaoGetADMMParentTao_ADMM", (PetscObject) tao);CHKERRQ(ierr); 8216285c0a3SHansol Suh ierr = PetscObjectComposeFunction((PetscObject)tao,"TaoADMMSetRegularizerType_C",TaoADMMSetRegularizerType_ADMM);CHKERRQ(ierr); 8226285c0a3SHansol Suh ierr = PetscObjectComposeFunction((PetscObject)tao,"TaoADMMGetRegularizerType_C",TaoADMMGetRegularizerType_ADMM);CHKERRQ(ierr); 8236285c0a3SHansol Suh ierr = PetscObjectComposeFunction((PetscObject)tao,"TaoADMMSetUpdateType_C",TaoADMMSetUpdateType_ADMM);CHKERRQ(ierr); 8246285c0a3SHansol Suh ierr = PetscObjectComposeFunction((PetscObject)tao,"TaoADMMGetUpdateType_C",TaoADMMGetUpdateType_ADMM);CHKERRQ(ierr); 8256285c0a3SHansol Suh PetscFunctionReturn(0); 8266285c0a3SHansol Suh } 8276285c0a3SHansol Suh 8286285c0a3SHansol Suh /*@ 8296285c0a3SHansol Suh TaoADMMSetMisfitHessianChangeStatus - Set boolean that determines whether Hessian matrix of misfit subsolver changes with respect to input vector. 8306285c0a3SHansol Suh 8316285c0a3SHansol Suh Collective on Tao 8326285c0a3SHansol Suh 8336285c0a3SHansol Suh Input Parameters: 8347f5c9be9SBarry Smith + tao - the Tao solver context. 8357f5c9be9SBarry Smith - b - the Hessian matrix change status boolean, PETSC_FALSE when the Hessian matrix does not change, TRUE otherwise. 8366285c0a3SHansol Suh 8376285c0a3SHansol Suh Level: advanced 838a5a2f7acSBarry Smith 839a5a2f7acSBarry Smith .seealso: TAOADMM 840a5a2f7acSBarry Smith 8416285c0a3SHansol Suh @*/ 8426285c0a3SHansol Suh PetscErrorCode TaoADMMSetMisfitHessianChangeStatus(Tao tao, PetscBool b) 8436285c0a3SHansol Suh { 8446285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 8456285c0a3SHansol Suh 8466285c0a3SHansol Suh PetscFunctionBegin; 8476285c0a3SHansol Suh am->Hxchange = b; 8486285c0a3SHansol Suh PetscFunctionReturn(0); 8496285c0a3SHansol Suh } 8506285c0a3SHansol Suh 8516285c0a3SHansol Suh /*@ 8526285c0a3SHansol Suh TaoADMMSetRegHessianChangeStatus - Set boolean that determines whether Hessian matrix of regularization subsolver changes with respect to input vector. 8537f5c9be9SBarry Smith 8546285c0a3SHansol Suh 8556285c0a3SHansol Suh Collective on Tao 8566285c0a3SHansol Suh 8576285c0a3SHansol Suh Input Parameters: 8586285c0a3SHansol Suh + tao - the Tao solver context 8597f5c9be9SBarry Smith - b - the Hessian matrix change status boolean, PETSC_FALSE when the Hessian matrix does not change, TRUE otherwise. 8606285c0a3SHansol Suh 8616285c0a3SHansol Suh Level: advanced 862a5a2f7acSBarry Smith 863a5a2f7acSBarry Smith .seealso: TAOADMM 864a5a2f7acSBarry Smith 8656285c0a3SHansol Suh @*/ 8666285c0a3SHansol Suh PetscErrorCode TaoADMMSetRegHessianChangeStatus(Tao tao, PetscBool b) 8676285c0a3SHansol Suh { 8686285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 8696285c0a3SHansol Suh 8706285c0a3SHansol Suh PetscFunctionBegin; 8716285c0a3SHansol Suh am->Hzchange = b; 8726285c0a3SHansol Suh PetscFunctionReturn(0); 8736285c0a3SHansol Suh } 8746285c0a3SHansol Suh 8756285c0a3SHansol Suh /*@ 8766285c0a3SHansol Suh TaoADMMSetSpectralPenalty - Set the spectral penalty (mu) value 8776285c0a3SHansol Suh 8786285c0a3SHansol Suh Collective on Tao 8796285c0a3SHansol Suh 8806285c0a3SHansol Suh Input Parameters: 8816285c0a3SHansol Suh + tao - the Tao solver context 8827f5c9be9SBarry Smith - mu - spectral penalty 8836285c0a3SHansol Suh 8846285c0a3SHansol Suh Level: advanced 8856285c0a3SHansol Suh 886a5a2f7acSBarry Smith .seealso: TaoADMMSetMinimumSpectralPenalty(), TAOADMM 8876285c0a3SHansol Suh @*/ 8886285c0a3SHansol Suh PetscErrorCode TaoADMMSetSpectralPenalty(Tao tao, PetscReal mu) 8896285c0a3SHansol Suh { 8906285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 8916285c0a3SHansol Suh 8926285c0a3SHansol Suh PetscFunctionBegin; 8936285c0a3SHansol Suh am->mu = mu; 8946285c0a3SHansol Suh PetscFunctionReturn(0); 8956285c0a3SHansol Suh } 8966285c0a3SHansol Suh 8976285c0a3SHansol Suh /*@ 8986285c0a3SHansol Suh TaoADMMGetSpectralPenalty - Get the spectral penalty (mu) value 8996285c0a3SHansol Suh 9006285c0a3SHansol Suh Collective on Tao 9016285c0a3SHansol Suh 9027f5c9be9SBarry Smith Input Parameter: 9037f5c9be9SBarry Smith . tao - the Tao solver context 9047f5c9be9SBarry Smith 9057f5c9be9SBarry Smith Output Parameter: 9067f5c9be9SBarry Smith . mu - spectral penalty 9076285c0a3SHansol Suh 9086285c0a3SHansol Suh Level: advanced 9096285c0a3SHansol Suh 910a5a2f7acSBarry Smith .seealso: TaoADMMSetMinimumSpectralPenalty(), TaoADMMSetSpectralPenalty(), TAOADMM 9116285c0a3SHansol Suh @*/ 9126285c0a3SHansol Suh PetscErrorCode TaoADMMGetSpectralPenalty(Tao tao, PetscReal *mu) 9136285c0a3SHansol Suh { 9146285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 9156285c0a3SHansol Suh 9166285c0a3SHansol Suh PetscFunctionBegin; 9176285c0a3SHansol Suh PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 9186285c0a3SHansol Suh PetscValidRealPointer(mu,2); 9196285c0a3SHansol Suh *mu = am->mu; 9206285c0a3SHansol Suh PetscFunctionReturn(0); 9216285c0a3SHansol Suh } 9226285c0a3SHansol Suh 9236285c0a3SHansol Suh /*@ 9246285c0a3SHansol Suh TaoADMMGetMisfitSubsolver - Get the pointer to the misfit subsolver inside ADMM 9256285c0a3SHansol Suh 9266285c0a3SHansol Suh Collective on Tao 9276285c0a3SHansol Suh 9287f5c9be9SBarry Smith Input Parameter: 9297f5c9be9SBarry Smith . tao - the Tao solver context 9307f5c9be9SBarry Smith 9317f5c9be9SBarry Smith Output Parameter: 9327f5c9be9SBarry Smith . misfit - the Tao subsolver context 9336285c0a3SHansol Suh 9346285c0a3SHansol Suh Level: advanced 935a5a2f7acSBarry Smith 936a5a2f7acSBarry Smith .seealso: TAOADMM 937a5a2f7acSBarry Smith 9386285c0a3SHansol Suh @*/ 9396285c0a3SHansol Suh PetscErrorCode TaoADMMGetMisfitSubsolver(Tao tao, Tao *misfit) 9406285c0a3SHansol Suh { 9416285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 9426285c0a3SHansol Suh 9436285c0a3SHansol Suh PetscFunctionBegin; 9446285c0a3SHansol Suh *misfit = am->subsolverX; 9456285c0a3SHansol Suh PetscFunctionReturn(0); 9466285c0a3SHansol Suh } 9476285c0a3SHansol Suh 9486285c0a3SHansol Suh /*@ 9496285c0a3SHansol Suh TaoADMMGetRegularizationSubsolver - Get the pointer to the regularization subsolver inside ADMM 9506285c0a3SHansol Suh 9516285c0a3SHansol Suh Collective on Tao 9526285c0a3SHansol Suh 9537f5c9be9SBarry Smith Input Parameter: 9547f5c9be9SBarry Smith . tao - the Tao solver context 9557f5c9be9SBarry Smith 9567f5c9be9SBarry Smith Output Parameter: 9577f5c9be9SBarry Smith . reg - the Tao subsolver context 9586285c0a3SHansol Suh 9596285c0a3SHansol Suh Level: advanced 960a5a2f7acSBarry Smith 961a5a2f7acSBarry Smith .seealso: TAOADMM 962a5a2f7acSBarry Smith 9636285c0a3SHansol Suh @*/ 9646285c0a3SHansol Suh PetscErrorCode TaoADMMGetRegularizationSubsolver(Tao tao, Tao *reg) 9656285c0a3SHansol Suh { 9666285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 9676285c0a3SHansol Suh 9686285c0a3SHansol Suh PetscFunctionBegin; 9696285c0a3SHansol Suh *reg = am->subsolverZ; 9706285c0a3SHansol Suh PetscFunctionReturn(0); 9716285c0a3SHansol Suh } 9726285c0a3SHansol Suh 9736285c0a3SHansol Suh /*@ 9746285c0a3SHansol Suh TaoADMMSetConstraintVectorRHS - Set the RHS constraint vector for ADMM 9756285c0a3SHansol Suh 9766285c0a3SHansol Suh Collective on Tao 9776285c0a3SHansol Suh 9786285c0a3SHansol Suh Input Parameters: 9796285c0a3SHansol Suh + tao - the Tao solver context 9806285c0a3SHansol Suh - c - RHS vector 9816285c0a3SHansol Suh 9826285c0a3SHansol Suh Level: advanced 983a5a2f7acSBarry Smith 984a5a2f7acSBarry Smith .seealso: TAOADMM 985a5a2f7acSBarry Smith 9866285c0a3SHansol Suh @*/ 9876285c0a3SHansol Suh PetscErrorCode TaoADMMSetConstraintVectorRHS(Tao tao, Vec c) 9886285c0a3SHansol Suh { 9896285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 9906285c0a3SHansol Suh 9916285c0a3SHansol Suh PetscFunctionBegin; 9926285c0a3SHansol Suh am->constraint = c; 9936285c0a3SHansol Suh PetscFunctionReturn(0); 9946285c0a3SHansol Suh } 9956285c0a3SHansol Suh 9966285c0a3SHansol Suh /*@ 9977f5c9be9SBarry Smith TaoADMMSetMinimumSpectralPenalty - Set the minimum value for the spectral penalty 9986285c0a3SHansol Suh 9996285c0a3SHansol Suh Collective on Tao 10006285c0a3SHansol Suh 10016285c0a3SHansol Suh Input Parameters: 10026285c0a3SHansol Suh + tao - the Tao solver context 10036285c0a3SHansol Suh - mu - minimum spectral penalty value 10046285c0a3SHansol Suh 10056285c0a3SHansol Suh Level: advanced 10066285c0a3SHansol Suh 1007a5a2f7acSBarry Smith .seealso: TaoADMMGetSpectralPenalty(), TAOADMM 10086285c0a3SHansol Suh @*/ 10096285c0a3SHansol Suh PetscErrorCode TaoADMMSetMinimumSpectralPenalty(Tao tao, PetscReal mu) 10106285c0a3SHansol Suh { 10116285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 10126285c0a3SHansol Suh 10136285c0a3SHansol Suh PetscFunctionBegin; 10146285c0a3SHansol Suh am->mumin= mu; 10156285c0a3SHansol Suh PetscFunctionReturn(0); 10166285c0a3SHansol Suh } 10176285c0a3SHansol Suh 10186285c0a3SHansol Suh /*@ 10196285c0a3SHansol Suh TaoADMMSetRegularizerCoefficient - Set the regularization coefficient lambda for L1 norm regularization case 10206285c0a3SHansol Suh 10216285c0a3SHansol Suh Collective on Tao 10226285c0a3SHansol Suh 10236285c0a3SHansol Suh Input Parameters: 10246285c0a3SHansol Suh + tao - the Tao solver context 10256285c0a3SHansol Suh - lambda - L1-norm regularizer coefficient 10266285c0a3SHansol Suh 10276285c0a3SHansol Suh Level: advanced 10287f5c9be9SBarry Smith 1029a5a2f7acSBarry Smith .seealso: TaoADMMSetMisfitConstraintJacobian(), TaoADMMSetRegularizerConstraintJacobian(), TAOADMM 10307f5c9be9SBarry Smith 10316285c0a3SHansol Suh @*/ 10326285c0a3SHansol Suh PetscErrorCode TaoADMMSetRegularizerCoefficient(Tao tao, PetscReal lambda) 10336285c0a3SHansol Suh { 10346285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 10356285c0a3SHansol Suh 10366285c0a3SHansol Suh PetscFunctionBegin; 10376285c0a3SHansol Suh am->lambda = lambda; 10386285c0a3SHansol Suh PetscFunctionReturn(0); 10396285c0a3SHansol Suh } 10406285c0a3SHansol Suh 10416285c0a3SHansol Suh /*@C 10427f5c9be9SBarry Smith TaoADMMSetMisfitConstraintJacobian - Set the constraint matrix B for the ADMM algorithm. Matrix B constrains the z variable. 10436285c0a3SHansol Suh 10446285c0a3SHansol Suh Collective on Tao 10456285c0a3SHansol Suh 1046a5a2f7acSBarry Smith Input Parameters: 1047a5a2f7acSBarry Smith + tao - the Tao solver context 10486285c0a3SHansol Suh . J - user-created regularizer constraint Jacobian matrix 10496285c0a3SHansol Suh . Jpre - user-created regularizer Jacobian constraint preconditioner matrix 10506285c0a3SHansol Suh . func - function pointer for the regularizer constraint Jacobian update function 10516285c0a3SHansol Suh - ctx - user context for the regularizer Hessian 10526285c0a3SHansol Suh 10536285c0a3SHansol Suh Level: advanced 10547f5c9be9SBarry Smith 1055a5a2f7acSBarry Smith .seealso: TaoADMMSetRegularizerCoefficient(), TaoADMMSetRegularizerConstraintJacobian(), TAOADMM 10567f5c9be9SBarry Smith 10576285c0a3SHansol Suh @*/ 10586285c0a3SHansol Suh PetscErrorCode TaoADMMSetMisfitConstraintJacobian(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx) 10596285c0a3SHansol Suh { 10606285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 10616285c0a3SHansol Suh PetscErrorCode ierr; 10626285c0a3SHansol Suh 10636285c0a3SHansol Suh PetscFunctionBegin; 10646285c0a3SHansol Suh PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 10656285c0a3SHansol Suh if (J) { 10666285c0a3SHansol Suh PetscValidHeaderSpecific(J,MAT_CLASSID,2); 10676285c0a3SHansol Suh PetscCheckSameComm(tao,1,J,2); 10686285c0a3SHansol Suh } 10696285c0a3SHansol Suh if (Jpre) { 10706285c0a3SHansol Suh PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3); 10716285c0a3SHansol Suh PetscCheckSameComm(tao,1,Jpre,3); 10726285c0a3SHansol Suh } 10736285c0a3SHansol Suh if (ctx) am->misfitjacobianP = ctx; 10746285c0a3SHansol Suh if (func) am->ops->misfitjac = func; 10756285c0a3SHansol Suh 10766285c0a3SHansol Suh if (J) { 10776285c0a3SHansol Suh ierr = PetscObjectReference((PetscObject)J);CHKERRQ(ierr); 10786285c0a3SHansol Suh ierr = MatDestroy(&am->JA);CHKERRQ(ierr); 10796285c0a3SHansol Suh am->JA = J; 10806285c0a3SHansol Suh } 10816285c0a3SHansol Suh if (Jpre) { 10826285c0a3SHansol Suh ierr = PetscObjectReference((PetscObject)Jpre);CHKERRQ(ierr); 10836285c0a3SHansol Suh ierr = MatDestroy(&am->JApre);CHKERRQ(ierr); 10846285c0a3SHansol Suh am->JApre = Jpre; 10856285c0a3SHansol Suh } 10866285c0a3SHansol Suh PetscFunctionReturn(0); 10876285c0a3SHansol Suh } 10886285c0a3SHansol Suh 10896285c0a3SHansol Suh /*@C 10906285c0a3SHansol Suh TaoADMMSetRegularizerConstraintJacobian - Set the constraint matrix B for ADMM algorithm. Matrix B constraints z variable. 10916285c0a3SHansol Suh 10926285c0a3SHansol Suh Collective on Tao 10936285c0a3SHansol Suh 10946285c0a3SHansol Suh Input Parameters: 10956285c0a3SHansol Suh + tao - the Tao solver context 10966285c0a3SHansol Suh . J - user-created regularizer constraint Jacobian matrix 10976285c0a3SHansol Suh . Jpre - user-created regularizer Jacobian constraint preconditioner matrix 10986285c0a3SHansol Suh . func - function pointer for the regularizer constraint Jacobian update function 10996285c0a3SHansol Suh - ctx - user context for the regularizer Hessian 11006285c0a3SHansol Suh 11016285c0a3SHansol Suh Level: advanced 11027f5c9be9SBarry Smith 1103a5a2f7acSBarry Smith .seealso: TaoADMMSetRegularizerCoefficient(), TaoADMMSetMisfitConstraintJacobian(), TAOADMM 11047f5c9be9SBarry Smith 11056285c0a3SHansol Suh @*/ 11066285c0a3SHansol Suh PetscErrorCode TaoADMMSetRegularizerConstraintJacobian(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx) 11076285c0a3SHansol Suh { 11086285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 11096285c0a3SHansol Suh PetscErrorCode ierr; 11106285c0a3SHansol Suh 11116285c0a3SHansol Suh PetscFunctionBegin; 11126285c0a3SHansol Suh PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 11136285c0a3SHansol Suh if (J) { 11146285c0a3SHansol Suh PetscValidHeaderSpecific(J,MAT_CLASSID,2); 11156285c0a3SHansol Suh PetscCheckSameComm(tao,1,J,2); 11166285c0a3SHansol Suh } 11176285c0a3SHansol Suh if (Jpre) { 11186285c0a3SHansol Suh PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3); 11196285c0a3SHansol Suh PetscCheckSameComm(tao,1,Jpre,3); 11206285c0a3SHansol Suh } 11216285c0a3SHansol Suh if (ctx) am->regjacobianP = ctx; 11226285c0a3SHansol Suh if (func) am->ops->regjac = func; 11236285c0a3SHansol Suh 11246285c0a3SHansol Suh if (J) { 11256285c0a3SHansol Suh ierr = PetscObjectReference((PetscObject)J);CHKERRQ(ierr); 11266285c0a3SHansol Suh ierr = MatDestroy(&am->JB);CHKERRQ(ierr); 11276285c0a3SHansol Suh am->JB = J; 11286285c0a3SHansol Suh } 11296285c0a3SHansol Suh if (Jpre) { 11306285c0a3SHansol Suh ierr = PetscObjectReference((PetscObject)Jpre);CHKERRQ(ierr); 11316285c0a3SHansol Suh ierr = MatDestroy(&am->JBpre);CHKERRQ(ierr); 11326285c0a3SHansol Suh am->JBpre = Jpre; 11336285c0a3SHansol Suh } 11346285c0a3SHansol Suh PetscFunctionReturn(0); 11356285c0a3SHansol Suh } 11366285c0a3SHansol Suh 11376285c0a3SHansol Suh /*@C 11387f5c9be9SBarry Smith TaoADMMSetMisfitObjectiveAndGradientRoutine - Sets the user-defined misfit call-back function 11397f5c9be9SBarry Smith 11407f5c9be9SBarry Smith Collective on tao 11416285c0a3SHansol Suh 11426285c0a3SHansol Suh Input Parameters: 11436285c0a3SHansol Suh + tao - the Tao context 11446285c0a3SHansol Suh . func - function pointer for the misfit value and gradient evaluation 11456285c0a3SHansol Suh - ctx - user context for the misfit 11466285c0a3SHansol Suh 11476285c0a3SHansol Suh Level: advanced 11487f5c9be9SBarry Smith 1149a5a2f7acSBarry Smith .seealso: TAOADMM 11507f5c9be9SBarry Smith 11516285c0a3SHansol Suh @*/ 11526285c0a3SHansol Suh PetscErrorCode TaoADMMSetMisfitObjectiveAndGradientRoutine(Tao tao, PetscErrorCode (*func)(Tao, Vec, PetscReal*, Vec, void*), void *ctx) 11536285c0a3SHansol Suh { 11546285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 11556285c0a3SHansol Suh 11566285c0a3SHansol Suh PetscFunctionBegin; 11576285c0a3SHansol Suh PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 11586285c0a3SHansol Suh am->misfitobjgradP = ctx; 11596285c0a3SHansol Suh am->ops->misfitobjgrad = func; 11606285c0a3SHansol Suh PetscFunctionReturn(0); 11616285c0a3SHansol Suh } 11626285c0a3SHansol Suh 11636285c0a3SHansol Suh /*@C 11646285c0a3SHansol Suh TaoADMMSetMisfitHessianRoutine - Sets the user-defined misfit Hessian call-back 11656285c0a3SHansol Suh function into the algorithm, to be used for subsolverX. 11666285c0a3SHansol Suh 11677f5c9be9SBarry Smith Collective on tao 11687f5c9be9SBarry Smith 11696285c0a3SHansol Suh Input Parameters: 11706285c0a3SHansol Suh + tao - the Tao context 11716285c0a3SHansol Suh . H - user-created matrix for the Hessian of the misfit term 11726285c0a3SHansol Suh . Hpre - user-created matrix for the preconditioner of Hessian of the misfit term 11736285c0a3SHansol Suh . func - function pointer for the misfit Hessian evaluation 11746285c0a3SHansol Suh - ctx - user context for the misfit Hessian 11756285c0a3SHansol Suh 11766285c0a3SHansol Suh Level: advanced 1177a5a2f7acSBarry Smith 1178a5a2f7acSBarry Smith .seealso: TAOADMM 1179a5a2f7acSBarry Smith 11806285c0a3SHansol Suh @*/ 11816285c0a3SHansol Suh PetscErrorCode TaoADMMSetMisfitHessianRoutine(Tao tao, Mat H, Mat Hpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx) 11826285c0a3SHansol Suh { 11836285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 11846285c0a3SHansol Suh PetscErrorCode ierr; 11856285c0a3SHansol Suh 11866285c0a3SHansol Suh PetscFunctionBegin; 11876285c0a3SHansol Suh PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 11886285c0a3SHansol Suh if (H) { 11896285c0a3SHansol Suh PetscValidHeaderSpecific(H,MAT_CLASSID,2); 11906285c0a3SHansol Suh PetscCheckSameComm(tao,1,H,2); 11916285c0a3SHansol Suh } 11926285c0a3SHansol Suh if (Hpre) { 11936285c0a3SHansol Suh PetscValidHeaderSpecific(Hpre,MAT_CLASSID,3); 11946285c0a3SHansol Suh PetscCheckSameComm(tao,1,Hpre,3); 11956285c0a3SHansol Suh } 11966285c0a3SHansol Suh if (ctx) { 11976285c0a3SHansol Suh am->misfithessP = ctx; 11986285c0a3SHansol Suh } 11996285c0a3SHansol Suh if (func) { 12006285c0a3SHansol Suh am->ops->misfithess = func; 12016285c0a3SHansol Suh } 12026285c0a3SHansol Suh if (H) { 12036285c0a3SHansol Suh ierr = PetscObjectReference((PetscObject)H);CHKERRQ(ierr); 12046285c0a3SHansol Suh ierr = MatDestroy(&am->Hx);CHKERRQ(ierr); 12056285c0a3SHansol Suh am->Hx = H; 12066285c0a3SHansol Suh } 12076285c0a3SHansol Suh if (Hpre) { 12086285c0a3SHansol Suh ierr = PetscObjectReference((PetscObject)Hpre);CHKERRQ(ierr); 12096285c0a3SHansol Suh ierr = MatDestroy(&am->Hxpre);CHKERRQ(ierr); 12106285c0a3SHansol Suh am->Hxpre = Hpre; 12116285c0a3SHansol Suh } 12126285c0a3SHansol Suh PetscFunctionReturn(0); 12136285c0a3SHansol Suh } 12146285c0a3SHansol Suh 12156285c0a3SHansol Suh /*@C 12167f5c9be9SBarry Smith TaoADMMSetRegularizerObjectiveAndGradientRoutine - Sets the user-defined regularizer call-back function 12177f5c9be9SBarry Smith 12187f5c9be9SBarry Smith Collective on tao 12196285c0a3SHansol Suh 12206285c0a3SHansol Suh Input Parameters: 12216285c0a3SHansol Suh + tao - the Tao context 12226285c0a3SHansol Suh . func - function pointer for the regularizer value and gradient evaluation 12236285c0a3SHansol Suh - ctx - user context for the regularizer 12246285c0a3SHansol Suh 12256285c0a3SHansol Suh Level: advanced 1226a5a2f7acSBarry Smith 1227a5a2f7acSBarry Smith .seealso: TAOADMM 1228a5a2f7acSBarry Smith 12296285c0a3SHansol Suh @*/ 12306285c0a3SHansol Suh PetscErrorCode TaoADMMSetRegularizerObjectiveAndGradientRoutine(Tao tao, PetscErrorCode (*func)(Tao, Vec, PetscReal*, Vec, void*), void *ctx) 12316285c0a3SHansol Suh { 12326285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 12336285c0a3SHansol Suh 12346285c0a3SHansol Suh PetscFunctionBegin; 12356285c0a3SHansol Suh PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 12366285c0a3SHansol Suh am->regobjgradP = ctx; 12376285c0a3SHansol Suh am->ops->regobjgrad = func; 12386285c0a3SHansol Suh PetscFunctionReturn(0); 12396285c0a3SHansol Suh } 12406285c0a3SHansol Suh 12416285c0a3SHansol Suh /*@C 12426285c0a3SHansol Suh TaoADMMSetRegularizerHessianRoutine - Sets the user-defined regularizer Hessian call-back 12437f5c9be9SBarry Smith function, to be used for subsolverZ. 12447f5c9be9SBarry Smith 12457f5c9be9SBarry Smith Collective on tao 12466285c0a3SHansol Suh 12476285c0a3SHansol Suh Input Parameters: 12486285c0a3SHansol Suh + tao - the Tao context 12496285c0a3SHansol Suh . H - user-created matrix for the Hessian of the regularization term 12506285c0a3SHansol Suh . Hpre - user-created matrix for the preconditioner of Hessian of the regularization term 12516285c0a3SHansol Suh . func - function pointer for the regularizer Hessian evaluation 12526285c0a3SHansol Suh - ctx - user context for the regularizer Hessian 12536285c0a3SHansol Suh 12546285c0a3SHansol Suh Level: advanced 1255a5a2f7acSBarry Smith 1256a5a2f7acSBarry Smith .seealso: TAOADMM 1257a5a2f7acSBarry Smith 12586285c0a3SHansol Suh @*/ 12596285c0a3SHansol Suh PetscErrorCode TaoADMMSetRegularizerHessianRoutine(Tao tao, Mat H, Mat Hpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx) 12606285c0a3SHansol Suh { 12616285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 12626285c0a3SHansol Suh PetscErrorCode ierr; 12636285c0a3SHansol Suh 12646285c0a3SHansol Suh PetscFunctionBegin; 12656285c0a3SHansol Suh PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 12666285c0a3SHansol Suh if (H) { 12676285c0a3SHansol Suh PetscValidHeaderSpecific(H,MAT_CLASSID,2); 12686285c0a3SHansol Suh PetscCheckSameComm(tao,1,H,2); 12696285c0a3SHansol Suh } 12706285c0a3SHansol Suh if (Hpre) { 12716285c0a3SHansol Suh PetscValidHeaderSpecific(Hpre,MAT_CLASSID,3); 12726285c0a3SHansol Suh PetscCheckSameComm(tao,1,Hpre,3); 12736285c0a3SHansol Suh } 12746285c0a3SHansol Suh if (ctx) { 12756285c0a3SHansol Suh am->reghessP = ctx; 12766285c0a3SHansol Suh } 12776285c0a3SHansol Suh if (func) { 12786285c0a3SHansol Suh am->ops->reghess = func; 12796285c0a3SHansol Suh } 12806285c0a3SHansol Suh if (H) { 12816285c0a3SHansol Suh ierr = PetscObjectReference((PetscObject)H);CHKERRQ(ierr); 12826285c0a3SHansol Suh ierr = MatDestroy(&am->Hz);CHKERRQ(ierr); 12836285c0a3SHansol Suh am->Hz = H; 12846285c0a3SHansol Suh } 12856285c0a3SHansol Suh if (Hpre) { 12866285c0a3SHansol Suh ierr = PetscObjectReference((PetscObject)Hpre);CHKERRQ(ierr); 12876285c0a3SHansol Suh ierr = MatDestroy(&am->Hzpre);CHKERRQ(ierr); 12886285c0a3SHansol Suh am->Hzpre = Hpre; 12896285c0a3SHansol Suh } 12906285c0a3SHansol Suh PetscFunctionReturn(0); 12916285c0a3SHansol Suh } 12926285c0a3SHansol Suh 12936285c0a3SHansol Suh /*@ 12947f5c9be9SBarry Smith TaoGetADMMParentTao - Gets pointer to parent ADMM tao, used by inner subsolver. 12956285c0a3SHansol Suh 12967f5c9be9SBarry Smith Collective on tao 12977f5c9be9SBarry Smith 12987f5c9be9SBarry Smith Input Parameter: 12997f5c9be9SBarry Smith . tao - the Tao context 13007f5c9be9SBarry Smith 13017f5c9be9SBarry Smith Output Parameter: 13027f5c9be9SBarry Smith . admm_tao - the parent Tao context 13036285c0a3SHansol Suh 13046285c0a3SHansol Suh Level: advanced 1305a5a2f7acSBarry Smith 1306a5a2f7acSBarry Smith .seealso: TAOADMM 1307a5a2f7acSBarry Smith 13086285c0a3SHansol Suh @*/ 13096285c0a3SHansol Suh PetscErrorCode TaoGetADMMParentTao(Tao tao, Tao *admm_tao) 13106285c0a3SHansol Suh { 13116285c0a3SHansol Suh PetscErrorCode ierr; 13126285c0a3SHansol Suh 13136285c0a3SHansol Suh PetscFunctionBegin; 13146285c0a3SHansol Suh PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 13156285c0a3SHansol Suh ierr = PetscObjectQuery((PetscObject)tao,"TaoGetADMMParentTao_ADMM", (PetscObject*) admm_tao);CHKERRQ(ierr); 13166285c0a3SHansol Suh PetscFunctionReturn(0); 13176285c0a3SHansol Suh } 13186285c0a3SHansol Suh 13196285c0a3SHansol Suh /*@ 13207f5c9be9SBarry Smith TaoADMMGetDualVector - Returns the dual vector associated with the current TAOADMM state 13216285c0a3SHansol Suh 13226285c0a3SHansol Suh Not Collective 13236285c0a3SHansol Suh 13246285c0a3SHansol Suh Input Parameter: 13257f5c9be9SBarry Smith . tao - the Tao context 13267f5c9be9SBarry Smith 13277f5c9be9SBarry Smith Output Parameter: 13287f5c9be9SBarry Smith . Y - the current solution 13296285c0a3SHansol Suh 13306285c0a3SHansol Suh Level: intermediate 1331a5a2f7acSBarry Smith 1332a5a2f7acSBarry Smith .seealso: TAOADMM 1333a5a2f7acSBarry Smith 13346285c0a3SHansol Suh @*/ 13356285c0a3SHansol Suh PetscErrorCode TaoADMMGetDualVector(Tao tao, Vec *Y) 13366285c0a3SHansol Suh { 13376285c0a3SHansol Suh TAO_ADMM *am = (TAO_ADMM*)tao->data; 13386285c0a3SHansol Suh 13396285c0a3SHansol Suh PetscFunctionBegin; 13406285c0a3SHansol Suh PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 13416285c0a3SHansol Suh *Y = am->y; 13426285c0a3SHansol Suh PetscFunctionReturn(0); 13436285c0a3SHansol Suh } 13446285c0a3SHansol Suh 13456285c0a3SHansol Suh /*@ 13466285c0a3SHansol Suh TaoADMMSetRegularizerType - Set regularizer type for ADMM routine 13476285c0a3SHansol Suh 13486285c0a3SHansol Suh Not Collective 13496285c0a3SHansol Suh 13506285c0a3SHansol Suh Input Parameter: 13516285c0a3SHansol Suh + tao - the Tao context 13527f5c9be9SBarry Smith - type - regularizer type 13537f5c9be9SBarry Smith 13547f5c9be9SBarry Smith Options Database: 13557f5c9be9SBarry Smith . -tao_admm_regularizer_type <admm_regularizer_user,admm_regularizer_soft_thresh> 13566285c0a3SHansol Suh 13576285c0a3SHansol Suh Level: intermediate 13586285c0a3SHansol Suh 1359a5a2f7acSBarry Smith .seealso: TaoADMMGetRegularizerType(), TaoADMMRegularizerType, TAOADMM 13606285c0a3SHansol Suh @*/ 13616285c0a3SHansol Suh PetscErrorCode TaoADMMSetRegularizerType(Tao tao, TaoADMMRegularizerType type) 13626285c0a3SHansol Suh { 13636285c0a3SHansol Suh PetscErrorCode ierr; 13646285c0a3SHansol Suh 13656285c0a3SHansol Suh PetscFunctionBegin; 13666285c0a3SHansol Suh PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 13676285c0a3SHansol Suh PetscValidLogicalCollectiveEnum(tao,type,2); 13686285c0a3SHansol Suh ierr = PetscTryMethod(tao,"TaoADMMSetRegularizerType_C",(Tao,TaoADMMRegularizerType),(tao,type));CHKERRQ(ierr); 13696285c0a3SHansol Suh PetscFunctionReturn(0); 13706285c0a3SHansol Suh } 13716285c0a3SHansol Suh 13726285c0a3SHansol Suh /*@ 13736285c0a3SHansol Suh TaoADMMGetRegularizerType - Gets the type of regularizer routine for ADMM 13746285c0a3SHansol Suh 13756285c0a3SHansol Suh Not Collective 13766285c0a3SHansol Suh 13776285c0a3SHansol Suh Input Parameter: 13786285c0a3SHansol Suh . tao - the Tao context 13796285c0a3SHansol Suh 13806285c0a3SHansol Suh Output Parameter: 13816285c0a3SHansol Suh . type - the type of regularizer 13826285c0a3SHansol Suh 13836285c0a3SHansol Suh Options Database: 13846285c0a3SHansol Suh . -tao_admm_regularizer_type <admm_regularizer_user,admm_regularizer_soft_thresh> 13856285c0a3SHansol Suh 13866285c0a3SHansol Suh Level: intermediate 13876285c0a3SHansol Suh 1388a5a2f7acSBarry Smith .seealso: TaoADMMSetRegularizerType(), TaoADMMRegularizerType, TAOADMM 13896285c0a3SHansol Suh @*/ 13906285c0a3SHansol Suh PetscErrorCode TaoADMMGetRegularizerType(Tao tao, TaoADMMRegularizerType *type) 13916285c0a3SHansol Suh { 13926285c0a3SHansol Suh PetscErrorCode ierr; 13936285c0a3SHansol Suh 13946285c0a3SHansol Suh PetscFunctionBegin; 13956285c0a3SHansol Suh PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 13966285c0a3SHansol Suh ierr = PetscUseMethod(tao,"TaoADMMGetRegularizerType_C",(Tao,TaoADMMRegularizerType*),(tao,type));CHKERRQ(ierr); 13976285c0a3SHansol Suh PetscFunctionReturn(0); 13986285c0a3SHansol Suh } 13996285c0a3SHansol Suh 14006285c0a3SHansol Suh /*@ 14016285c0a3SHansol Suh TaoADMMSetUpdateType - Set update routine for ADMM routine 14026285c0a3SHansol Suh 14036285c0a3SHansol Suh Not Collective 14046285c0a3SHansol Suh 14056285c0a3SHansol Suh Input Parameter: 14066285c0a3SHansol Suh + tao - the Tao context 14076285c0a3SHansol Suh - type - spectral parameter update type 14086285c0a3SHansol Suh 14096285c0a3SHansol Suh Level: intermediate 14106285c0a3SHansol Suh 1411a5a2f7acSBarry Smith .seealso: TaoADMMGetUpdateType(), TaoADMMUpdateType, TAOADMM 14126285c0a3SHansol Suh @*/ 14136285c0a3SHansol Suh PetscErrorCode TaoADMMSetUpdateType(Tao tao, TaoADMMUpdateType type) 14146285c0a3SHansol Suh { 14156285c0a3SHansol Suh PetscErrorCode ierr; 14166285c0a3SHansol Suh 14176285c0a3SHansol Suh PetscFunctionBegin; 14186285c0a3SHansol Suh PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 14196285c0a3SHansol Suh PetscValidLogicalCollectiveEnum(tao,type,2); 14206285c0a3SHansol Suh ierr = PetscTryMethod(tao,"TaoADMMSetUpdateType_C",(Tao,TaoADMMUpdateType),(tao,type));CHKERRQ(ierr); 14216285c0a3SHansol Suh PetscFunctionReturn(0); 14226285c0a3SHansol Suh } 14236285c0a3SHansol Suh 14246285c0a3SHansol Suh /*@ 14256285c0a3SHansol Suh TaoADMMGetUpdateType - Gets the type of spectral penalty update routine for ADMM 14266285c0a3SHansol Suh 14276285c0a3SHansol Suh Not Collective 14286285c0a3SHansol Suh 14296285c0a3SHansol Suh Input Parameter: 1430f0fc11ceSJed Brown . tao - the Tao context 14316285c0a3SHansol Suh 14326285c0a3SHansol Suh Output Parameter: 14336285c0a3SHansol Suh . type - the type of spectral penalty update routine 14346285c0a3SHansol Suh 14356285c0a3SHansol Suh Level: intermediate 14366285c0a3SHansol Suh 1437a5a2f7acSBarry Smith .seealso: TaoADMMSetUpdateType(), TaoADMMUpdateType, TAOADMM 14386285c0a3SHansol Suh @*/ 14396285c0a3SHansol Suh PetscErrorCode TaoADMMGetUpdateType(Tao tao, TaoADMMUpdateType *type) 14406285c0a3SHansol Suh { 14416285c0a3SHansol Suh PetscErrorCode ierr; 14426285c0a3SHansol Suh 14436285c0a3SHansol Suh PetscFunctionBegin; 14446285c0a3SHansol Suh PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 14456285c0a3SHansol Suh ierr = PetscUseMethod(tao,"TaoADMMGetUpdateType_C",(Tao,TaoADMMUpdateType*),(tao,type));CHKERRQ(ierr); 14466285c0a3SHansol Suh PetscFunctionReturn(0); 14476285c0a3SHansol Suh } 14486285c0a3SHansol Suh 1449