1 #include <../src/tao/constrained/impls/admm/admm.h> /*I "petsctao.h" I*/ 2 #include <petsctao.h> 3 #include <petsc/private/petscimpl.h> 4 5 /* Updates terminating criteria 6 * 7 * 1 ||r_k|| = ||Ax+Bz-c|| =< catol_admm* max{||Ax||,||Bz||,||c||} 8 * 9 * 2. Updates dual residual, d_k 10 * 11 * 3. ||d_k|| = ||mu*A^T*B(z_k-z_{k-1})|| =< gatol_admm * ||A^Ty|| */ 12 13 static PetscBool cited = PETSC_FALSE; 14 static const char citation[] = 15 "@misc{xu2017adaptive,\n" 16 " title={Adaptive Relaxed ADMM: Convergence Theory and Practical Implementation},\n" 17 " author={Zheng Xu and Mario A. T. Figueiredo and Xiaoming Yuan and Christoph Studer and Tom Goldstein},\n" 18 " year={2017},\n" 19 " eprint={1704.02712},\n" 20 " archivePrefix={arXiv},\n" 21 " primaryClass={cs.CV}\n" 22 "} \n"; 23 24 const char *const TaoADMMRegularizerTypes[] = {"REGULARIZER_USER","REGULARIZER_SOFT_THRESH","TaoADMMRegularizerType","TAO_ADMM_",NULL}; 25 const char *const TaoADMMUpdateTypes[] = {"UPDATE_BASIC","UPDATE_ADAPTIVE","UPDATE_ADAPTIVE_RELAXED","TaoADMMUpdateType","TAO_ADMM_",NULL}; 26 const char *const TaoALMMTypes[] = {"CLASSIC","PHR","TaoALMMType","TAO_ALMM_",NULL}; 27 28 static PetscErrorCode TaoADMMToleranceUpdate(Tao tao) 29 { 30 TAO_ADMM *am = (TAO_ADMM*)tao->data; 31 PetscReal Axnorm,Bznorm,ATynorm,temp; 32 Vec tempJR,tempL; 33 Tao mis; 34 35 PetscFunctionBegin; 36 mis = am->subsolverX; 37 tempJR = am->workJacobianRight; 38 tempL = am->workLeft; 39 /* ATy */ 40 PetscCall(TaoComputeJacobianEquality(mis, am->y, mis->jacobian_equality, mis->jacobian_equality_pre)); 41 PetscCall(MatMultTranspose(mis->jacobian_equality,am->y,tempJR)); 42 PetscCall(VecNorm(tempJR,NORM_2,&ATynorm)); 43 /* dualres = mu * ||AT(Bz-Bzold)||_2 */ 44 PetscCall(VecWAXPY(tempJR,-1.,am->Bzold,am->Bz)); 45 PetscCall(MatMultTranspose(mis->jacobian_equality,tempJR,tempL)); 46 PetscCall(VecNorm(tempL,NORM_2,&am->dualres)); 47 am->dualres *= am->mu; 48 49 /* ||Ax||_2, ||Bz||_2 */ 50 PetscCall(VecNorm(am->Ax,NORM_2,&Axnorm)); 51 PetscCall(VecNorm(am->Bz,NORM_2,&Bznorm)); 52 53 /* Set catol to be catol_admm * max{||Ax||,||Bz||,||c||} * 54 * Set gatol to be gatol_admm * ||A^Ty|| * 55 * while cnorm is ||r_k||_2, and gnorm is ||d_k||_2 */ 56 temp = am->catol_admm * PetscMax(Axnorm, (!am->const_norm) ? Bznorm : PetscMax(Bznorm,am->const_norm)); 57 PetscCall(TaoSetConstraintTolerances(tao,temp,PETSC_DEFAULT)); 58 PetscCall(TaoSetTolerances(tao, am->gatol_admm*ATynorm, PETSC_DEFAULT,PETSC_DEFAULT)); 59 PetscFunctionReturn(0); 60 } 61 62 /* Penaly Update for Adaptive ADMM. */ 63 static PetscErrorCode AdaptiveADMMPenaltyUpdate(Tao tao) 64 { 65 TAO_ADMM *am = (TAO_ADMM*)tao->data; 66 PetscReal ydiff_norm, yhatdiff_norm, Axdiff_norm, Bzdiff_norm, Axyhat, Bzy, a_sd, a_mg, a_k, b_sd, b_mg, b_k; 67 PetscBool hflag, gflag; 68 Vec tempJR,tempJR2; 69 70 PetscFunctionBegin; 71 tempJR = am->workJacobianRight; 72 tempJR2 = am->workJacobianRight2; 73 hflag = PETSC_FALSE; 74 gflag = PETSC_FALSE; 75 a_k = -1; 76 b_k = -1; 77 78 PetscCall(VecWAXPY(tempJR,-1.,am->Axold,am->Ax)); 79 PetscCall(VecWAXPY(tempJR2,-1.,am->yhatold,am->yhat)); 80 PetscCall(VecNorm(tempJR,NORM_2,&Axdiff_norm)); 81 PetscCall(VecNorm(tempJR2,NORM_2,&yhatdiff_norm)); 82 PetscCall(VecDot(tempJR,tempJR2,&Axyhat)); 83 84 PetscCall(VecWAXPY(tempJR,-1.,am->Bz0,am->Bz)); 85 PetscCall(VecWAXPY(tempJR2,-1.,am->y,am->y0)); 86 PetscCall(VecNorm(tempJR,NORM_2,&Bzdiff_norm)); 87 PetscCall(VecNorm(tempJR2,NORM_2,&ydiff_norm)); 88 PetscCall(VecDot(tempJR,tempJR2,&Bzy)); 89 90 if (Axyhat > am->orthval*Axdiff_norm*yhatdiff_norm + am->mueps) { 91 hflag = PETSC_TRUE; 92 a_sd = PetscSqr(yhatdiff_norm)/Axyhat; /* alphaSD */ 93 a_mg = Axyhat/PetscSqr(Axdiff_norm); /* alphaMG */ 94 a_k = (a_mg/a_sd) > 0.5 ? a_mg : a_sd - 0.5*a_mg; 95 } 96 if (Bzy > am->orthval*Bzdiff_norm*ydiff_norm + am->mueps) { 97 gflag = PETSC_TRUE; 98 b_sd = PetscSqr(ydiff_norm)/Bzy; /* betaSD */ 99 b_mg = Bzy/PetscSqr(Bzdiff_norm); /* betaMG */ 100 b_k = (b_mg/b_sd) > 0.5 ? b_mg : b_sd - 0.5*b_mg; 101 } 102 am->muold = am->mu; 103 if (gflag && hflag) { 104 am->mu = PetscSqrtReal(a_k*b_k); 105 } else if (hflag) { 106 am->mu = a_k; 107 } else if (gflag) { 108 am->mu = b_k; 109 } 110 if (am->mu > am->muold) { 111 am->mu = am->muold; 112 } 113 if (am->mu < am->mumin) { 114 am->mu = am->mumin; 115 } 116 PetscFunctionReturn(0); 117 } 118 119 static PetscErrorCode TaoADMMSetRegularizerType_ADMM(Tao tao, TaoADMMRegularizerType type) 120 { 121 TAO_ADMM *am = (TAO_ADMM*)tao->data; 122 123 PetscFunctionBegin; 124 am->regswitch = type; 125 PetscFunctionReturn(0); 126 } 127 128 static PetscErrorCode TaoADMMGetRegularizerType_ADMM(Tao tao, TaoADMMRegularizerType *type) 129 { 130 TAO_ADMM *am = (TAO_ADMM*)tao->data; 131 132 PetscFunctionBegin; 133 *type = am->regswitch; 134 PetscFunctionReturn(0); 135 } 136 137 static PetscErrorCode TaoADMMSetUpdateType_ADMM(Tao tao, TaoADMMUpdateType type) 138 { 139 TAO_ADMM *am = (TAO_ADMM*)tao->data; 140 141 PetscFunctionBegin; 142 am->update = type; 143 PetscFunctionReturn(0); 144 } 145 146 static PetscErrorCode TaoADMMGetUpdateType_ADMM(Tao tao, TaoADMMUpdateType *type) 147 { 148 TAO_ADMM *am = (TAO_ADMM*)tao->data; 149 150 PetscFunctionBegin; 151 *type = am->update; 152 PetscFunctionReturn(0); 153 } 154 155 /* This routine updates Jacobians with new x,z vectors, 156 * and then updates Ax and Bz vectors, then computes updated residual vector*/ 157 static PetscErrorCode ADMMUpdateConstraintResidualVector(Tao tao, Vec x, Vec z, Vec Ax, Vec Bz, Vec residual) 158 { 159 TAO_ADMM *am = (TAO_ADMM*)tao->data; 160 Tao mis,reg; 161 162 PetscFunctionBegin; 163 mis = am->subsolverX; 164 reg = am->subsolverZ; 165 PetscCall(TaoComputeJacobianEquality(mis, x, mis->jacobian_equality, mis->jacobian_equality_pre)); 166 PetscCall(MatMult(mis->jacobian_equality,x,Ax)); 167 PetscCall(TaoComputeJacobianEquality(reg, z, reg->jacobian_equality, reg->jacobian_equality_pre)); 168 PetscCall(MatMult(reg->jacobian_equality,z,Bz)); 169 170 PetscCall(VecWAXPY(residual,1.,Bz,Ax)); 171 if (am->constraint != NULL) { 172 PetscCall(VecAXPY(residual,-1.,am->constraint)); 173 } 174 PetscFunctionReturn(0); 175 } 176 177 /* Updates Augmented Lagrangians to given routines * 178 * For subsolverX, routine needs to be ComputeObjectiveAndGraidnet 179 * Separate Objective and Gradient routines are not supported. */ 180 static PetscErrorCode SubObjGradUpdate(Tao tao, Vec x, PetscReal *f, Vec g, void *ptr) 181 { 182 Tao parent = (Tao)ptr; 183 TAO_ADMM *am = (TAO_ADMM*)parent->data; 184 PetscReal temp,temp2; 185 Vec tempJR; 186 187 PetscFunctionBegin; 188 tempJR = am->workJacobianRight; 189 PetscCall(ADMMUpdateConstraintResidualVector(parent, x, am->subsolverZ->solution, am->Ax, am->Bz, am->residual)); 190 PetscCall((*am->ops->misfitobjgrad)(am->subsolverX,x,f,g,am->misfitobjgradP)); 191 192 am->last_misfit_val = *f; 193 /* Objective Add + yT(Ax+Bz-c) + mu/2*||Ax+Bz-c||_2^2 */ 194 PetscCall(VecTDot(am->residual,am->y,&temp)); 195 PetscCall(VecTDot(am->residual,am->residual,&temp2)); 196 *f += temp + (am->mu/2)*temp2; 197 198 /* Gradient. Add + mu*AT(Ax+Bz-c) + yTA*/ 199 PetscCall(MatMultTranspose(tao->jacobian_equality,am->residual,tempJR)); 200 PetscCall(VecAXPY(g,am->mu,tempJR)); 201 PetscCall(MatMultTranspose(tao->jacobian_equality,am->y,tempJR)); 202 PetscCall(VecAXPY(g,1.,tempJR)); 203 PetscFunctionReturn(0); 204 } 205 206 /* Updates Augmented Lagrangians to given routines 207 * For subsolverZ, routine needs to be ComputeObjectiveAndGraidnet 208 * Separate Objective and Gradient routines are not supported. */ 209 static PetscErrorCode RegObjGradUpdate(Tao tao, Vec z, PetscReal *f, Vec g, void *ptr) 210 { 211 Tao parent = (Tao)ptr; 212 TAO_ADMM *am = (TAO_ADMM*)parent->data; 213 PetscReal temp,temp2; 214 Vec tempJR; 215 216 PetscFunctionBegin; 217 tempJR = am->workJacobianRight; 218 PetscCall(ADMMUpdateConstraintResidualVector(parent, am->subsolverX->solution, z, am->Ax, am->Bz, am->residual)); 219 PetscCall((*am->ops->regobjgrad)(am->subsolverZ,z,f,g,am->regobjgradP)); 220 am->last_reg_val= *f; 221 /* Objective Add + yT(Ax+Bz-c) + mu/2*||Ax+Bz-c||_2^2 */ 222 PetscCall(VecTDot(am->residual,am->y,&temp)); 223 PetscCall(VecTDot(am->residual,am->residual,&temp2)); 224 *f += temp + (am->mu/2)*temp2; 225 226 /* Gradient. Add + mu*BT(Ax+Bz-c) + yTB*/ 227 PetscCall(MatMultTranspose(am->subsolverZ->jacobian_equality,am->residual,tempJR)); 228 PetscCall(VecAXPY(g,am->mu,tempJR)); 229 PetscCall(MatMultTranspose(am->subsolverZ->jacobian_equality,am->y,tempJR)); 230 PetscCall(VecAXPY(g,1.,tempJR)); 231 PetscFunctionReturn(0); 232 } 233 234 /* Computes epsilon padded L1 norm lambda*sum(sqrt(x^2+eps^2)-eps */ 235 static PetscErrorCode ADMML1EpsilonNorm(Tao tao, Vec x, PetscReal eps, PetscReal *norm) 236 { 237 TAO_ADMM *am = (TAO_ADMM*)tao->data; 238 PetscInt N; 239 240 PetscFunctionBegin; 241 PetscCall(VecGetSize(am->workLeft,&N)); 242 PetscCall(VecPointwiseMult(am->workLeft,x,x)); 243 PetscCall(VecShift(am->workLeft,am->l1epsilon*am->l1epsilon)); 244 PetscCall(VecSqrtAbs(am->workLeft)); 245 PetscCall(VecSum(am->workLeft,norm)); 246 *norm += N*am->l1epsilon; 247 *norm *= am->lambda; 248 PetscFunctionReturn(0); 249 } 250 251 static PetscErrorCode ADMMInternalHessianUpdate(Mat H, Mat Constraint, PetscBool Identity, void *ptr) 252 { 253 TAO_ADMM *am = (TAO_ADMM*)ptr; 254 255 PetscFunctionBegin; 256 switch (am->update) { 257 case (TAO_ADMM_UPDATE_BASIC): 258 break; 259 case (TAO_ADMM_UPDATE_ADAPTIVE): 260 case (TAO_ADMM_UPDATE_ADAPTIVE_RELAXED): 261 if (H && (am->muold != am->mu)) { 262 if (!Identity) { 263 PetscCall(MatAXPY(H,am->mu-am->muold,Constraint,DIFFERENT_NONZERO_PATTERN)); 264 } else { 265 PetscCall(MatShift(H,am->mu-am->muold)); 266 } 267 } 268 break; 269 } 270 PetscFunctionReturn(0); 271 } 272 273 /* Updates Hessian - adds second derivative of augmented Lagrangian 274 * H \gets H + \rho*ATA 275 * Here, \rho does not change in TAO_ADMM_UPDATE_BASIC - thus no-op 276 * For ADAPTAIVE,ADAPTIVE_RELAXED, 277 * H \gets H + (\rho-\rhoold)*ATA 278 * Here, we assume that A is linear constraint i.e., doesnt change. 279 * Thus, for both ADAPTIVE, and RELAXED, ATA matrix is pre-set (except for A=I (null case)) see TaoSetUp_ADMM */ 280 static PetscErrorCode SubHessianUpdate(Tao tao, Vec x, Mat H, Mat Hpre, void *ptr) 281 { 282 Tao parent = (Tao)ptr; 283 TAO_ADMM *am = (TAO_ADMM*)parent->data; 284 285 PetscFunctionBegin; 286 if (am->Hxchange) { 287 /* Case where Hessian gets updated with respect to x vector input. */ 288 PetscCall((*am->ops->misfithess)(am->subsolverX,x,H,Hpre,am->misfithessP)); 289 PetscCall(ADMMInternalHessianUpdate(am->subsolverX->hessian,am->ATA,am->xJI,am)); 290 } else if (am->Hxbool) { 291 /* Hessian doesn't get updated. H(x) = c */ 292 /* Update Lagrangian only only per TAO call */ 293 PetscCall(ADMMInternalHessianUpdate(am->subsolverX->hessian,am->ATA,am->xJI,am)); 294 am->Hxbool = PETSC_FALSE; 295 } 296 PetscFunctionReturn(0); 297 } 298 299 /* Same as SubHessianUpdate, except for B matrix instead of A matrix */ 300 static PetscErrorCode RegHessianUpdate(Tao tao, Vec z, Mat H, Mat Hpre, void *ptr) 301 { 302 Tao parent = (Tao)ptr; 303 TAO_ADMM *am = (TAO_ADMM*)parent->data; 304 305 PetscFunctionBegin; 306 307 if (am->Hzchange) { 308 /* Case where Hessian gets updated with respect to x vector input. */ 309 PetscCall((*am->ops->reghess)(am->subsolverZ,z,H,Hpre,am->reghessP)); 310 PetscCall(ADMMInternalHessianUpdate(am->subsolverZ->hessian,am->BTB,am->zJI,am)); 311 } else if (am->Hzbool) { 312 /* Hessian doesn't get updated. H(x) = c */ 313 /* Update Lagrangian only only per TAO call */ 314 PetscCall(ADMMInternalHessianUpdate(am->subsolverZ->hessian,am->BTB,am->zJI,am)); 315 am->Hzbool = PETSC_FALSE; 316 } 317 PetscFunctionReturn(0); 318 } 319 320 /* Shell Matrix routine for A matrix. 321 * This gets used when user puts NULL for 322 * TaoSetJacobianEqualityRoutine(tao, NULL,NULL, ...) 323 * Essentially sets A=I*/ 324 static PetscErrorCode JacobianIdentity(Mat mat,Vec in,Vec out) 325 { 326 PetscFunctionBegin; 327 PetscCall(VecCopy(in,out)); 328 PetscFunctionReturn(0); 329 } 330 331 /* Shell Matrix routine for B matrix. 332 * This gets used when user puts NULL for 333 * TaoADMMSetRegularizerConstraintJacobian(tao, NULL,NULL, ...) 334 * Sets B=-I */ 335 static PetscErrorCode JacobianIdentityB(Mat mat,Vec in,Vec out) 336 { 337 PetscFunctionBegin; 338 PetscCall(VecCopy(in,out)); 339 PetscCall(VecScale(out,-1.)); 340 PetscFunctionReturn(0); 341 } 342 343 /* Solve f(x) + g(z) s.t. Ax + Bz = c */ 344 static PetscErrorCode TaoSolve_ADMM(Tao tao) 345 { 346 TAO_ADMM *am = (TAO_ADMM*)tao->data; 347 PetscInt N; 348 PetscReal reg_func; 349 PetscBool is_reg_shell; 350 Vec tempL; 351 352 PetscFunctionBegin; 353 if (am->regswitch != TAO_ADMM_REGULARIZER_SOFT_THRESH) { 354 PetscCheck(am->subsolverX->ops->computejacobianequality,PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_WRONGSTATE,"Must call TaoADMMSetMisfitConstraintJacobian() first"); 355 PetscCheck(am->subsolverZ->ops->computejacobianequality,PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_WRONGSTATE,"Must call TaoADMMSetRegularizerConstraintJacobian() first"); 356 if (am->constraint != NULL) { 357 PetscCall(VecNorm(am->constraint,NORM_2,&am->const_norm)); 358 } 359 } 360 tempL = am->workLeft; 361 PetscCall(VecGetSize(tempL,&N)); 362 363 if (am->Hx && am->ops->misfithess) { 364 PetscCall(TaoSetHessian(am->subsolverX, am->Hx, am->Hx, SubHessianUpdate, tao)); 365 } 366 367 if (!am->zJI) { 368 /* Currently, B is assumed to be a linear system, i.e., not getting updated*/ 369 PetscCall(MatTransposeMatMult(am->JB,am->JB,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&(am->BTB))); 370 } 371 if (!am->xJI) { 372 /* Currently, A is assumed to be a linear system, i.e., not getting updated*/ 373 PetscCall(MatTransposeMatMult(am->subsolverX->jacobian_equality,am->subsolverX->jacobian_equality,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&(am->ATA))); 374 } 375 376 is_reg_shell = PETSC_FALSE; 377 378 PetscCall(PetscObjectTypeCompare((PetscObject)am->subsolverZ, TAOSHELL, &is_reg_shell)); 379 380 if (!is_reg_shell) { 381 switch (am->regswitch) { 382 case (TAO_ADMM_REGULARIZER_USER): 383 PetscCheck(am->ops->regobjgrad,PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_WRONGSTATE,"Must call TaoADMMSetRegularizerObjectiveAndGradientRoutine() first if one wishes to use TAO_ADMM_REGULARIZER_USER with non-TAOSHELL type"); 384 break; 385 case (TAO_ADMM_REGULARIZER_SOFT_THRESH): 386 /* Soft Threshold. */ 387 break; 388 } 389 if (am->ops->regobjgrad) { 390 PetscCall(TaoSetObjectiveAndGradient(am->subsolverZ, NULL, RegObjGradUpdate, tao)); 391 } 392 if (am->Hz && am->ops->reghess) { 393 PetscCall(TaoSetHessian(am->subsolverZ, am->Hz, am->Hzpre, RegHessianUpdate, tao)); 394 } 395 } 396 397 switch (am->update) { 398 case TAO_ADMM_UPDATE_BASIC: 399 if (am->subsolverX->hessian) { 400 /* In basic case, Hessian does not get updated w.r.t. to spectral penalty 401 * Here, when A is set, i.e., am->xJI, add mu*ATA to Hessian*/ 402 if (!am->xJI) { 403 PetscCall(MatAXPY(am->subsolverX->hessian,am->mu,am->ATA,DIFFERENT_NONZERO_PATTERN)); 404 } else { 405 PetscCall(MatShift(am->subsolverX->hessian,am->mu)); 406 } 407 } 408 if (am->subsolverZ->hessian && am->regswitch == TAO_ADMM_REGULARIZER_USER) { 409 if (am->regswitch == TAO_ADMM_REGULARIZER_USER && !am->zJI) { 410 PetscCall(MatAXPY(am->subsolverZ->hessian,am->mu,am->BTB,DIFFERENT_NONZERO_PATTERN)); 411 } else { 412 PetscCall(MatShift(am->subsolverZ->hessian,am->mu)); 413 } 414 } 415 break; 416 case TAO_ADMM_UPDATE_ADAPTIVE: 417 case TAO_ADMM_UPDATE_ADAPTIVE_RELAXED: 418 break; 419 } 420 421 PetscCall(PetscCitationsRegister(citation,&cited)); 422 tao->reason = TAO_CONTINUE_ITERATING; 423 424 while (tao->reason == TAO_CONTINUE_ITERATING) { 425 if (tao->ops->update) { 426 PetscCall((*tao->ops->update)(tao, tao->niter, tao->user_update)); 427 } 428 PetscCall(VecCopy(am->Bz, am->Bzold)); 429 430 /* x update */ 431 PetscCall(TaoSolve(am->subsolverX)); 432 PetscCall(TaoComputeJacobianEquality(am->subsolverX, am->subsolverX->solution, am->subsolverX->jacobian_equality, am->subsolverX->jacobian_equality_pre)); 433 PetscCall(MatMult(am->subsolverX->jacobian_equality, am->subsolverX->solution,am->Ax)); 434 435 am->Hxbool = PETSC_TRUE; 436 437 /* z update */ 438 switch (am->regswitch) { 439 case TAO_ADMM_REGULARIZER_USER: 440 PetscCall(TaoSolve(am->subsolverZ)); 441 break; 442 case TAO_ADMM_REGULARIZER_SOFT_THRESH: 443 /* L1 assumes A,B jacobians are identity nxn matrix */ 444 PetscCall(VecWAXPY(am->workJacobianRight,1/am->mu,am->y,am->Ax)); 445 PetscCall(TaoSoftThreshold(am->workJacobianRight,-am->lambda/am->mu,am->lambda/am->mu,am->subsolverZ->solution)); 446 break; 447 } 448 am->Hzbool = PETSC_TRUE; 449 /* Returns Ax + Bz - c with updated Ax,Bz vectors */ 450 PetscCall(ADMMUpdateConstraintResidualVector(tao, am->subsolverX->solution, am->subsolverZ->solution, am->Ax, am->Bz, am->residual)); 451 /* Dual variable, y += y + mu*(Ax+Bz-c) */ 452 PetscCall(VecWAXPY(am->y, am->mu, am->residual, am->yold)); 453 454 /* stopping tolerance update */ 455 PetscCall(TaoADMMToleranceUpdate(tao)); 456 457 /* Updating Spectral Penalty */ 458 switch (am->update) { 459 case TAO_ADMM_UPDATE_BASIC: 460 am->muold = am->mu; 461 break; 462 case TAO_ADMM_UPDATE_ADAPTIVE: 463 case TAO_ADMM_UPDATE_ADAPTIVE_RELAXED: 464 if (tao->niter == 0) { 465 PetscCall(VecCopy(am->y, am->y0)); 466 PetscCall(VecWAXPY(am->residual, 1., am->Ax, am->Bzold)); 467 if (am->constraint) { 468 PetscCall(VecAXPY(am->residual, -1., am->constraint)); 469 } 470 PetscCall(VecWAXPY(am->yhatold, -am->mu, am->residual, am->yold)); 471 PetscCall(VecCopy(am->Ax, am->Axold)); 472 PetscCall(VecCopy(am->Bz, am->Bz0)); 473 am->muold = am->mu; 474 } else if (tao->niter % am->T == 1) { 475 /* we have compute Bzold in a previous iteration, and we computed Ax above */ 476 PetscCall(VecWAXPY(am->residual, 1., am->Ax, am->Bzold)); 477 if (am->constraint) { 478 PetscCall(VecAXPY(am->residual, -1., am->constraint)); 479 } 480 PetscCall(VecWAXPY(am->yhat, -am->mu, am->residual, am->yold)); 481 PetscCall(AdaptiveADMMPenaltyUpdate(tao)); 482 PetscCall(VecCopy(am->Ax, am->Axold)); 483 PetscCall(VecCopy(am->Bz, am->Bz0)); 484 PetscCall(VecCopy(am->yhat, am->yhatold)); 485 PetscCall(VecCopy(am->y, am->y0)); 486 } else { 487 am->muold = am->mu; 488 } 489 break; 490 default: 491 break; 492 } 493 tao->niter++; 494 495 /* Calculate original function values. misfit part was done in TaoADMMToleranceUpdate*/ 496 switch (am->regswitch) { 497 case TAO_ADMM_REGULARIZER_USER: 498 if (is_reg_shell) { 499 PetscCall(ADMML1EpsilonNorm(tao,am->subsolverZ->solution,am->l1epsilon,®_func)); 500 } else { 501 (*am->ops->regobjgrad)(am->subsolverZ,am->subsolverX->solution,®_func,tempL,am->regobjgradP); 502 } 503 break; 504 case TAO_ADMM_REGULARIZER_SOFT_THRESH: 505 PetscCall(ADMML1EpsilonNorm(tao,am->subsolverZ->solution,am->l1epsilon,®_func)); 506 break; 507 } 508 PetscCall(VecCopy(am->y,am->yold)); 509 PetscCall(ADMMUpdateConstraintResidualVector(tao, am->subsolverX->solution, am->subsolverZ->solution, am->Ax, am->Bz, am->residual)); 510 PetscCall(VecNorm(am->residual,NORM_2,&am->resnorm)); 511 PetscCall(TaoLogConvergenceHistory(tao,am->last_misfit_val + reg_func,am->dualres,am->resnorm,tao->ksp_its)); 512 513 PetscCall(TaoMonitor(tao,tao->niter,am->last_misfit_val + reg_func,am->dualres,am->resnorm,1.0)); 514 PetscCall((*tao->ops->convergencetest)(tao,tao->cnvP)); 515 } 516 /* Update vectors */ 517 PetscCall(VecCopy(am->subsolverX->solution,tao->solution)); 518 PetscCall(VecCopy(am->subsolverX->gradient,tao->gradient)); 519 PetscCall(PetscObjectCompose((PetscObject)am->subsolverX,"TaoGetADMMParentTao_ADMM", NULL)); 520 PetscCall(PetscObjectCompose((PetscObject)am->subsolverZ,"TaoGetADMMParentTao_ADMM", NULL)); 521 PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMSetRegularizerType_C",NULL)); 522 PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMGetRegularizerType_C",NULL)); 523 PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMSetUpdateType_C",NULL)); 524 PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMGetUpdateType_C",NULL)); 525 PetscFunctionReturn(0); 526 } 527 528 static PetscErrorCode TaoSetFromOptions_ADMM(PetscOptionItems *PetscOptionsObject,Tao tao) 529 { 530 TAO_ADMM *am = (TAO_ADMM*)tao->data; 531 PetscErrorCode ierr; 532 533 PetscFunctionBegin; 534 PetscCall(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. ")); 535 PetscCall(PetscOptionsReal("-tao_admm_regularizer_coefficient","regularizer constant","",am->lambda,&am->lambda,NULL)); 536 PetscCall(PetscOptionsReal("-tao_admm_spectral_penalty","Constant for Augmented Lagrangian term.","",am->mu,&am->mu,NULL)); 537 PetscCall(PetscOptionsReal("-tao_admm_relaxation_parameter","x relaxation parameter for Z update.","",am->gamma,&am->gamma,NULL)); 538 PetscCall(PetscOptionsReal("-tao_admm_tolerance_update_factor","ADMM dynamic tolerance update factor.","",am->tol,&am->tol,NULL)); 539 PetscCall(PetscOptionsReal("-tao_admm_spectral_penalty_update_factor","ADMM spectral penalty update curvature safeguard value.","",am->orthval,&am->orthval,NULL)); 540 PetscCall(PetscOptionsReal("-tao_admm_minimum_spectral_penalty","Set ADMM minimum spectral penalty.","",am->mumin,&am->mumin,NULL)); 541 ierr = PetscOptionsEnum("-tao_admm_dual_update","Lagrangian dual update policy","TaoADMMUpdateType", 542 TaoADMMUpdateTypes,(PetscEnum)am->update,(PetscEnum*)&am->update,NULL);PetscCall(ierr); 543 ierr = PetscOptionsEnum("-tao_admm_regularizer_type","ADMM regularizer update rule","TaoADMMRegularizerType", 544 TaoADMMRegularizerTypes,(PetscEnum)am->regswitch,(PetscEnum*)&am->regswitch,NULL);PetscCall(ierr); 545 PetscCall(PetscOptionsTail()); 546 PetscCall(TaoSetFromOptions(am->subsolverX)); 547 if (am->regswitch != TAO_ADMM_REGULARIZER_SOFT_THRESH) { 548 PetscCall(TaoSetFromOptions(am->subsolverZ)); 549 } 550 PetscFunctionReturn(0); 551 } 552 553 static PetscErrorCode TaoView_ADMM(Tao tao,PetscViewer viewer) 554 { 555 TAO_ADMM *am = (TAO_ADMM*)tao->data; 556 557 PetscFunctionBegin; 558 PetscCall(PetscViewerASCIIPushTab(viewer)); 559 PetscCall(TaoView(am->subsolverX,viewer)); 560 PetscCall(TaoView(am->subsolverZ,viewer)); 561 PetscCall(PetscViewerASCIIPopTab(viewer)); 562 PetscFunctionReturn(0); 563 } 564 565 static PetscErrorCode TaoSetUp_ADMM(Tao tao) 566 { 567 TAO_ADMM *am = (TAO_ADMM*)tao->data; 568 PetscInt n,N,M; 569 570 PetscFunctionBegin; 571 PetscCall(VecGetLocalSize(tao->solution,&n)); 572 PetscCall(VecGetSize(tao->solution,&N)); 573 /* If Jacobian is given as NULL, it means Jacobian is identity matrix with size of solution vector */ 574 if (!am->JB) { 575 am->zJI = PETSC_TRUE; 576 PetscCall(MatCreateShell(PetscObjectComm((PetscObject)tao),n,n,PETSC_DETERMINE,PETSC_DETERMINE,NULL,&am->JB)); 577 PetscCall(MatShellSetOperation(am->JB,MATOP_MULT,(void (*)(void))JacobianIdentityB)); 578 PetscCall(MatShellSetOperation(am->JB,MATOP_MULT_TRANSPOSE,(void (*)(void))JacobianIdentityB)); 579 PetscCall(MatShellSetOperation(am->JB,MATOP_TRANSPOSE_MAT_MULT,(void (*)(void))JacobianIdentityB)); 580 am->JBpre = am->JB; 581 } 582 if (!am->JA) { 583 am->xJI = PETSC_TRUE; 584 PetscCall(MatCreateShell(PetscObjectComm((PetscObject)tao),n,n,PETSC_DETERMINE,PETSC_DETERMINE,NULL,&am->JA)); 585 PetscCall(MatShellSetOperation(am->JA,MATOP_MULT,(void (*)(void))JacobianIdentity)); 586 PetscCall(MatShellSetOperation(am->JA,MATOP_MULT_TRANSPOSE,(void (*)(void))JacobianIdentity)); 587 PetscCall(MatShellSetOperation(am->JA,MATOP_TRANSPOSE_MAT_MULT,(void (*)(void))JacobianIdentity)); 588 am->JApre = am->JA; 589 } 590 PetscCall(MatCreateVecs(am->JA,NULL,&am->Ax)); 591 if (!tao->gradient) { 592 PetscCall(VecDuplicate(tao->solution,&tao->gradient)); 593 } 594 PetscCall(TaoSetSolution(am->subsolverX, tao->solution)); 595 if (!am->z) { 596 PetscCall(VecDuplicate(tao->solution,&am->z)); 597 PetscCall(VecSet(am->z,0.0)); 598 } 599 PetscCall(TaoSetSolution(am->subsolverZ, am->z)); 600 if (!am->workLeft) { 601 PetscCall(VecDuplicate(tao->solution,&am->workLeft)); 602 } 603 if (!am->Axold) { 604 PetscCall(VecDuplicate(am->Ax,&am->Axold)); 605 } 606 if (!am->workJacobianRight) { 607 PetscCall(VecDuplicate(am->Ax,&am->workJacobianRight)); 608 } 609 if (!am->workJacobianRight2) { 610 PetscCall(VecDuplicate(am->Ax,&am->workJacobianRight2)); 611 } 612 if (!am->Bz) { 613 PetscCall(VecDuplicate(am->Ax,&am->Bz)); 614 } 615 if (!am->Bzold) { 616 PetscCall(VecDuplicate(am->Ax,&am->Bzold)); 617 } 618 if (!am->Bz0) { 619 PetscCall(VecDuplicate(am->Ax,&am->Bz0)); 620 } 621 if (!am->y) { 622 PetscCall(VecDuplicate(am->Ax,&am->y)); 623 PetscCall(VecSet(am->y,0.0)); 624 } 625 if (!am->yold) { 626 PetscCall(VecDuplicate(am->Ax,&am->yold)); 627 PetscCall(VecSet(am->yold,0.0)); 628 } 629 if (!am->y0) { 630 PetscCall(VecDuplicate(am->Ax,&am->y0)); 631 PetscCall(VecSet(am->y0,0.0)); 632 } 633 if (!am->yhat) { 634 PetscCall(VecDuplicate(am->Ax,&am->yhat)); 635 PetscCall(VecSet(am->yhat,0.0)); 636 } 637 if (!am->yhatold) { 638 PetscCall(VecDuplicate(am->Ax,&am->yhatold)); 639 PetscCall(VecSet(am->yhatold,0.0)); 640 } 641 if (!am->residual) { 642 PetscCall(VecDuplicate(am->Ax,&am->residual)); 643 PetscCall(VecSet(am->residual,0.0)); 644 } 645 if (!am->constraint) { 646 am->constraint = NULL; 647 } else { 648 PetscCall(VecGetSize(am->constraint,&M)); 649 PetscCheck(M == N,PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_WRONGSTATE,"Solution vector and constraint vector must be of same size!"); 650 } 651 652 /* Save changed tao tolerance for adaptive tolerance */ 653 if (tao->gatol_changed) { 654 am->gatol_admm = tao->gatol; 655 } 656 if (tao->catol_changed) { 657 am->catol_admm = tao->catol; 658 } 659 660 /*Update spectral and dual elements to X subsolver */ 661 PetscCall(TaoSetObjectiveAndGradient(am->subsolverX, NULL, SubObjGradUpdate, tao)); 662 PetscCall(TaoSetJacobianEqualityRoutine(am->subsolverX,am->JA,am->JApre, am->ops->misfitjac, am->misfitjacobianP)); 663 PetscCall(TaoSetJacobianEqualityRoutine(am->subsolverZ,am->JB,am->JBpre, am->ops->regjac, am->regjacobianP)); 664 PetscFunctionReturn(0); 665 } 666 667 static PetscErrorCode TaoDestroy_ADMM(Tao tao) 668 { 669 TAO_ADMM *am = (TAO_ADMM*)tao->data; 670 671 PetscFunctionBegin; 672 PetscCall(VecDestroy(&am->z)); 673 PetscCall(VecDestroy(&am->Ax)); 674 PetscCall(VecDestroy(&am->Axold)); 675 PetscCall(VecDestroy(&am->Bz)); 676 PetscCall(VecDestroy(&am->Bzold)); 677 PetscCall(VecDestroy(&am->Bz0)); 678 PetscCall(VecDestroy(&am->residual)); 679 PetscCall(VecDestroy(&am->y)); 680 PetscCall(VecDestroy(&am->yold)); 681 PetscCall(VecDestroy(&am->y0)); 682 PetscCall(VecDestroy(&am->yhat)); 683 PetscCall(VecDestroy(&am->yhatold)); 684 PetscCall(VecDestroy(&am->workLeft)); 685 PetscCall(VecDestroy(&am->workJacobianRight)); 686 PetscCall(VecDestroy(&am->workJacobianRight2)); 687 688 PetscCall(MatDestroy(&am->JA)); 689 PetscCall(MatDestroy(&am->JB)); 690 if (!am->xJI) { 691 PetscCall(MatDestroy(&am->JApre)); 692 } 693 if (!am->zJI) { 694 PetscCall(MatDestroy(&am->JBpre)); 695 } 696 if (am->Hx) { 697 PetscCall(MatDestroy(&am->Hx)); 698 PetscCall(MatDestroy(&am->Hxpre)); 699 } 700 if (am->Hz) { 701 PetscCall(MatDestroy(&am->Hz)); 702 PetscCall(MatDestroy(&am->Hzpre)); 703 } 704 PetscCall(MatDestroy(&am->ATA)); 705 PetscCall(MatDestroy(&am->BTB)); 706 PetscCall(TaoDestroy(&am->subsolverX)); 707 PetscCall(TaoDestroy(&am->subsolverZ)); 708 am->parent = NULL; 709 PetscCall(PetscFree(tao->data)); 710 PetscFunctionReturn(0); 711 } 712 713 /*MC 714 715 TAOADMM - Alternating direction method of multipliers method fo solving linear problems with 716 constraints. in a min_x f(x) + g(z) s.t. Ax+Bz=c. 717 This algorithm employs two sub Tao solvers, of which type can be specified 718 by the user. User need to provide ObjectiveAndGradient routine, and/or HessianRoutine for both subsolvers. 719 Hessians can be given boolean flag determining whether they change with respect to a input vector. This can be set via 720 TaoADMMSet{Misfit,Regularizer}HessianChangeStatus. 721 Second subsolver does support TAOSHELL. It should be noted that L1-norm is used for objective value for TAOSHELL type. 722 There is option to set regularizer option, and currently soft-threshold is implemented. For spectral penalty update, 723 currently there are basic option and adaptive option. 724 Constraint is set at Ax+Bz=c, and A and B can be set with TaoADMMSet{Misfit,Regularizer}ConstraintJacobian. 725 c can be set with TaoADMMSetConstraintVectorRHS. 726 The user can also provide regularizer weight for second subsolver. 727 728 References: 729 . * - Xu, Zheng and Figueiredo, Mario A. T. and Yuan, Xiaoming and Studer, Christoph and Goldstein, Tom 730 "Adaptive Relaxed ADMM: Convergence Theory and Practical Implementation" 731 The IEEE Conference on Computer Vision and Pattern Recognition (CVPR), July, 2017. 732 733 Options Database Keys: 734 + -tao_admm_regularizer_coefficient - regularizer constant (default 1.e-6) 735 . -tao_admm_spectral_penalty - Constant for Augmented Lagrangian term (default 1.) 736 . -tao_admm_relaxation_parameter - relaxation parameter for Z update (default 1.) 737 . -tao_admm_tolerance_update_factor - ADMM dynamic tolerance update factor (default 1.e-12) 738 . -tao_admm_spectral_penalty_update_factor - ADMM spectral penalty update curvature safeguard value (default 0.2) 739 . -tao_admm_minimum_spectral_penalty - Set ADMM minimum spectral penalty (default 0) 740 . -tao_admm_dual_update - Lagrangian dual update policy ("basic","adaptive","adaptive-relaxed") (default "basic") 741 - -tao_admm_regularizer_type - ADMM regularizer update rule ("user","soft-threshold") (default "soft-threshold") 742 743 Level: beginner 744 745 .seealso: TaoADMMSetMisfitHessianChangeStatus(), TaoADMMSetRegHessianChangeStatus(), TaoADMMGetSpectralPenalty(), 746 TaoADMMGetMisfitSubsolver(), TaoADMMGetRegularizationSubsolver(), TaoADMMSetConstraintVectorRHS(), 747 TaoADMMSetMinimumSpectralPenalty(), TaoADMMSetRegularizerCoefficient(), 748 TaoADMMSetRegularizerConstraintJacobian(), TaoADMMSetMisfitConstraintJacobian(), 749 TaoADMMSetMisfitObjectiveAndGradientRoutine(), TaoADMMSetMisfitHessianRoutine(), 750 TaoADMMSetRegularizerObjectiveAndGradientRoutine(), TaoADMMSetRegularizerHessianRoutine(), 751 TaoGetADMMParentTao(), TaoADMMGetDualVector(), TaoADMMSetRegularizerType(), 752 TaoADMMGetRegularizerType(), TaoADMMSetUpdateType(), TaoADMMGetUpdateType() 753 M*/ 754 755 PETSC_EXTERN PetscErrorCode TaoCreate_ADMM(Tao tao) 756 { 757 TAO_ADMM *am; 758 759 PetscFunctionBegin; 760 PetscCall(PetscNewLog(tao,&am)); 761 762 tao->ops->destroy = TaoDestroy_ADMM; 763 tao->ops->setup = TaoSetUp_ADMM; 764 tao->ops->setfromoptions = TaoSetFromOptions_ADMM; 765 tao->ops->view = TaoView_ADMM; 766 tao->ops->solve = TaoSolve_ADMM; 767 768 tao->data = (void*)am; 769 am->l1epsilon = 1e-6; 770 am->lambda = 1e-4; 771 am->mu = 1.; 772 am->muold = 0.; 773 am->mueps = PETSC_MACHINE_EPSILON; 774 am->mumin = 0.; 775 am->orthval = 0.2; 776 am->T = 2; 777 am->parent = tao; 778 am->update = TAO_ADMM_UPDATE_BASIC; 779 am->regswitch = TAO_ADMM_REGULARIZER_SOFT_THRESH; 780 am->tol = PETSC_SMALL; 781 am->const_norm = 0; 782 am->resnorm = 0; 783 am->dualres = 0; 784 am->ops->regobjgrad = NULL; 785 am->ops->reghess = NULL; 786 am->gamma = 1; 787 am->regobjgradP = NULL; 788 am->reghessP = NULL; 789 am->gatol_admm = 1e-8; 790 am->catol_admm = 0; 791 am->Hxchange = PETSC_TRUE; 792 am->Hzchange = PETSC_TRUE; 793 am->Hzbool = PETSC_TRUE; 794 am->Hxbool = PETSC_TRUE; 795 796 PetscCall(TaoCreate(PetscObjectComm((PetscObject)tao),&am->subsolverX)); 797 PetscCall(TaoSetOptionsPrefix(am->subsolverX,"misfit_")); 798 PetscCall(PetscObjectIncrementTabLevel((PetscObject)am->subsolverX,(PetscObject)tao,1)); 799 PetscCall(TaoCreate(PetscObjectComm((PetscObject)tao),&am->subsolverZ)); 800 PetscCall(TaoSetOptionsPrefix(am->subsolverZ,"reg_")); 801 PetscCall(PetscObjectIncrementTabLevel((PetscObject)am->subsolverZ,(PetscObject)tao,1)); 802 803 PetscCall(TaoSetType(am->subsolverX,TAONLS)); 804 PetscCall(TaoSetType(am->subsolverZ,TAONLS)); 805 PetscCall(PetscObjectCompose((PetscObject)am->subsolverX,"TaoGetADMMParentTao_ADMM", (PetscObject) tao)); 806 PetscCall(PetscObjectCompose((PetscObject)am->subsolverZ,"TaoGetADMMParentTao_ADMM", (PetscObject) tao)); 807 PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMSetRegularizerType_C",TaoADMMSetRegularizerType_ADMM)); 808 PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMGetRegularizerType_C",TaoADMMGetRegularizerType_ADMM)); 809 PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMSetUpdateType_C",TaoADMMSetUpdateType_ADMM)); 810 PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMGetUpdateType_C",TaoADMMGetUpdateType_ADMM)); 811 PetscFunctionReturn(0); 812 } 813 814 /*@ 815 TaoADMMSetMisfitHessianChangeStatus - Set boolean that determines whether Hessian matrix of misfit subsolver changes with respect to input vector. 816 817 Collective on Tao 818 819 Input Parameters: 820 + tao - the Tao solver context. 821 - b - the Hessian matrix change status boolean, PETSC_FALSE when the Hessian matrix does not change, TRUE otherwise. 822 823 Level: advanced 824 825 .seealso: TAOADMM 826 827 @*/ 828 PetscErrorCode TaoADMMSetMisfitHessianChangeStatus(Tao tao, PetscBool b) 829 { 830 TAO_ADMM *am = (TAO_ADMM*)tao->data; 831 832 PetscFunctionBegin; 833 am->Hxchange = b; 834 PetscFunctionReturn(0); 835 } 836 837 /*@ 838 TaoADMMSetRegHessianChangeStatus - Set boolean that determines whether Hessian matrix of regularization subsolver changes with respect to input vector. 839 840 Collective on Tao 841 842 Input Parameters: 843 + tao - the Tao solver context 844 - b - the Hessian matrix change status boolean, PETSC_FALSE when the Hessian matrix does not change, TRUE otherwise. 845 846 Level: advanced 847 848 .seealso: TAOADMM 849 850 @*/ 851 PetscErrorCode TaoADMMSetRegHessianChangeStatus(Tao tao, PetscBool b) 852 { 853 TAO_ADMM *am = (TAO_ADMM*)tao->data; 854 855 PetscFunctionBegin; 856 am->Hzchange = b; 857 PetscFunctionReturn(0); 858 } 859 860 /*@ 861 TaoADMMSetSpectralPenalty - Set the spectral penalty (mu) value 862 863 Collective on Tao 864 865 Input Parameters: 866 + tao - the Tao solver context 867 - mu - spectral penalty 868 869 Level: advanced 870 871 .seealso: TaoADMMSetMinimumSpectralPenalty(), TAOADMM 872 @*/ 873 PetscErrorCode TaoADMMSetSpectralPenalty(Tao tao, PetscReal mu) 874 { 875 TAO_ADMM *am = (TAO_ADMM*)tao->data; 876 877 PetscFunctionBegin; 878 am->mu = mu; 879 PetscFunctionReturn(0); 880 } 881 882 /*@ 883 TaoADMMGetSpectralPenalty - Get the spectral penalty (mu) value 884 885 Collective on Tao 886 887 Input Parameter: 888 . tao - the Tao solver context 889 890 Output Parameter: 891 . mu - spectral penalty 892 893 Level: advanced 894 895 .seealso: TaoADMMSetMinimumSpectralPenalty(), TaoADMMSetSpectralPenalty(), TAOADMM 896 @*/ 897 PetscErrorCode TaoADMMGetSpectralPenalty(Tao tao, PetscReal *mu) 898 { 899 TAO_ADMM *am = (TAO_ADMM*)tao->data; 900 901 PetscFunctionBegin; 902 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 903 PetscValidRealPointer(mu,2); 904 *mu = am->mu; 905 PetscFunctionReturn(0); 906 } 907 908 /*@ 909 TaoADMMGetMisfitSubsolver - Get the pointer to the misfit subsolver inside ADMM 910 911 Collective on Tao 912 913 Input Parameter: 914 . tao - the Tao solver context 915 916 Output Parameter: 917 . misfit - the Tao subsolver context 918 919 Level: advanced 920 921 .seealso: TAOADMM 922 923 @*/ 924 PetscErrorCode TaoADMMGetMisfitSubsolver(Tao tao, Tao *misfit) 925 { 926 TAO_ADMM *am = (TAO_ADMM*)tao->data; 927 928 PetscFunctionBegin; 929 *misfit = am->subsolverX; 930 PetscFunctionReturn(0); 931 } 932 933 /*@ 934 TaoADMMGetRegularizationSubsolver - Get the pointer to the regularization subsolver inside ADMM 935 936 Collective on Tao 937 938 Input Parameter: 939 . tao - the Tao solver context 940 941 Output Parameter: 942 . reg - the Tao subsolver context 943 944 Level: advanced 945 946 .seealso: TAOADMM 947 948 @*/ 949 PetscErrorCode TaoADMMGetRegularizationSubsolver(Tao tao, Tao *reg) 950 { 951 TAO_ADMM *am = (TAO_ADMM*)tao->data; 952 953 PetscFunctionBegin; 954 *reg = am->subsolverZ; 955 PetscFunctionReturn(0); 956 } 957 958 /*@ 959 TaoADMMSetConstraintVectorRHS - Set the RHS constraint vector for ADMM 960 961 Collective on Tao 962 963 Input Parameters: 964 + tao - the Tao solver context 965 - c - RHS vector 966 967 Level: advanced 968 969 .seealso: TAOADMM 970 971 @*/ 972 PetscErrorCode TaoADMMSetConstraintVectorRHS(Tao tao, Vec c) 973 { 974 TAO_ADMM *am = (TAO_ADMM*)tao->data; 975 976 PetscFunctionBegin; 977 am->constraint = c; 978 PetscFunctionReturn(0); 979 } 980 981 /*@ 982 TaoADMMSetMinimumSpectralPenalty - Set the minimum value for the spectral penalty 983 984 Collective on Tao 985 986 Input Parameters: 987 + tao - the Tao solver context 988 - mu - minimum spectral penalty value 989 990 Level: advanced 991 992 .seealso: TaoADMMGetSpectralPenalty(), TAOADMM 993 @*/ 994 PetscErrorCode TaoADMMSetMinimumSpectralPenalty(Tao tao, PetscReal mu) 995 { 996 TAO_ADMM *am = (TAO_ADMM*)tao->data; 997 998 PetscFunctionBegin; 999 am->mumin= mu; 1000 PetscFunctionReturn(0); 1001 } 1002 1003 /*@ 1004 TaoADMMSetRegularizerCoefficient - Set the regularization coefficient lambda for L1 norm regularization case 1005 1006 Collective on Tao 1007 1008 Input Parameters: 1009 + tao - the Tao solver context 1010 - lambda - L1-norm regularizer coefficient 1011 1012 Level: advanced 1013 1014 .seealso: TaoADMMSetMisfitConstraintJacobian(), TaoADMMSetRegularizerConstraintJacobian(), TAOADMM 1015 1016 @*/ 1017 PetscErrorCode TaoADMMSetRegularizerCoefficient(Tao tao, PetscReal lambda) 1018 { 1019 TAO_ADMM *am = (TAO_ADMM*)tao->data; 1020 1021 PetscFunctionBegin; 1022 am->lambda = lambda; 1023 PetscFunctionReturn(0); 1024 } 1025 1026 /*@C 1027 TaoADMMSetMisfitConstraintJacobian - Set the constraint matrix B for the ADMM algorithm. Matrix B constrains the z variable. 1028 1029 Collective on Tao 1030 1031 Input Parameters: 1032 + tao - the Tao solver context 1033 . J - user-created regularizer constraint Jacobian matrix 1034 . Jpre - user-created regularizer Jacobian constraint preconditioner matrix 1035 . func - function pointer for the regularizer constraint Jacobian update function 1036 - ctx - user context for the regularizer Hessian 1037 1038 Level: advanced 1039 1040 .seealso: TaoADMMSetRegularizerCoefficient(), TaoADMMSetRegularizerConstraintJacobian(), TAOADMM 1041 1042 @*/ 1043 PetscErrorCode TaoADMMSetMisfitConstraintJacobian(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx) 1044 { 1045 TAO_ADMM *am = (TAO_ADMM*)tao->data; 1046 1047 PetscFunctionBegin; 1048 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 1049 if (J) { 1050 PetscValidHeaderSpecific(J,MAT_CLASSID,2); 1051 PetscCheckSameComm(tao,1,J,2); 1052 } 1053 if (Jpre) { 1054 PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3); 1055 PetscCheckSameComm(tao,1,Jpre,3); 1056 } 1057 if (ctx) am->misfitjacobianP = ctx; 1058 if (func) am->ops->misfitjac = func; 1059 1060 if (J) { 1061 PetscCall(PetscObjectReference((PetscObject)J)); 1062 PetscCall(MatDestroy(&am->JA)); 1063 am->JA = J; 1064 } 1065 if (Jpre) { 1066 PetscCall(PetscObjectReference((PetscObject)Jpre)); 1067 PetscCall(MatDestroy(&am->JApre)); 1068 am->JApre = Jpre; 1069 } 1070 PetscFunctionReturn(0); 1071 } 1072 1073 /*@C 1074 TaoADMMSetRegularizerConstraintJacobian - Set the constraint matrix B for ADMM algorithm. Matrix B constraints z variable. 1075 1076 Collective on Tao 1077 1078 Input Parameters: 1079 + tao - the Tao solver context 1080 . J - user-created regularizer constraint Jacobian matrix 1081 . Jpre - user-created regularizer Jacobian constraint preconditioner matrix 1082 . func - function pointer for the regularizer constraint Jacobian update function 1083 - ctx - user context for the regularizer Hessian 1084 1085 Level: advanced 1086 1087 .seealso: TaoADMMSetRegularizerCoefficient(), TaoADMMSetMisfitConstraintJacobian(), TAOADMM 1088 1089 @*/ 1090 PetscErrorCode TaoADMMSetRegularizerConstraintJacobian(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx) 1091 { 1092 TAO_ADMM *am = (TAO_ADMM*)tao->data; 1093 1094 PetscFunctionBegin; 1095 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 1096 if (J) { 1097 PetscValidHeaderSpecific(J,MAT_CLASSID,2); 1098 PetscCheckSameComm(tao,1,J,2); 1099 } 1100 if (Jpre) { 1101 PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3); 1102 PetscCheckSameComm(tao,1,Jpre,3); 1103 } 1104 if (ctx) am->regjacobianP = ctx; 1105 if (func) am->ops->regjac = func; 1106 1107 if (J) { 1108 PetscCall(PetscObjectReference((PetscObject)J)); 1109 PetscCall(MatDestroy(&am->JB)); 1110 am->JB = J; 1111 } 1112 if (Jpre) { 1113 PetscCall(PetscObjectReference((PetscObject)Jpre)); 1114 PetscCall(MatDestroy(&am->JBpre)); 1115 am->JBpre = Jpre; 1116 } 1117 PetscFunctionReturn(0); 1118 } 1119 1120 /*@C 1121 TaoADMMSetMisfitObjectiveAndGradientRoutine - Sets the user-defined misfit call-back function 1122 1123 Collective on tao 1124 1125 Input Parameters: 1126 + tao - the Tao context 1127 . func - function pointer for the misfit value and gradient evaluation 1128 - ctx - user context for the misfit 1129 1130 Level: advanced 1131 1132 .seealso: TAOADMM 1133 1134 @*/ 1135 PetscErrorCode TaoADMMSetMisfitObjectiveAndGradientRoutine(Tao tao, PetscErrorCode (*func)(Tao, Vec, PetscReal*, Vec, void*), void *ctx) 1136 { 1137 TAO_ADMM *am = (TAO_ADMM*)tao->data; 1138 1139 PetscFunctionBegin; 1140 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 1141 am->misfitobjgradP = ctx; 1142 am->ops->misfitobjgrad = func; 1143 PetscFunctionReturn(0); 1144 } 1145 1146 /*@C 1147 TaoADMMSetMisfitHessianRoutine - Sets the user-defined misfit Hessian call-back 1148 function into the algorithm, to be used for subsolverX. 1149 1150 Collective on tao 1151 1152 Input Parameters: 1153 + tao - the Tao context 1154 . H - user-created matrix for the Hessian of the misfit term 1155 . Hpre - user-created matrix for the preconditioner of Hessian of the misfit term 1156 . func - function pointer for the misfit Hessian evaluation 1157 - ctx - user context for the misfit Hessian 1158 1159 Level: advanced 1160 1161 .seealso: TAOADMM 1162 1163 @*/ 1164 PetscErrorCode TaoADMMSetMisfitHessianRoutine(Tao tao, Mat H, Mat Hpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx) 1165 { 1166 TAO_ADMM *am = (TAO_ADMM*)tao->data; 1167 1168 PetscFunctionBegin; 1169 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 1170 if (H) { 1171 PetscValidHeaderSpecific(H,MAT_CLASSID,2); 1172 PetscCheckSameComm(tao,1,H,2); 1173 } 1174 if (Hpre) { 1175 PetscValidHeaderSpecific(Hpre,MAT_CLASSID,3); 1176 PetscCheckSameComm(tao,1,Hpre,3); 1177 } 1178 if (ctx) { 1179 am->misfithessP = ctx; 1180 } 1181 if (func) { 1182 am->ops->misfithess = func; 1183 } 1184 if (H) { 1185 PetscCall(PetscObjectReference((PetscObject)H)); 1186 PetscCall(MatDestroy(&am->Hx)); 1187 am->Hx = H; 1188 } 1189 if (Hpre) { 1190 PetscCall(PetscObjectReference((PetscObject)Hpre)); 1191 PetscCall(MatDestroy(&am->Hxpre)); 1192 am->Hxpre = Hpre; 1193 } 1194 PetscFunctionReturn(0); 1195 } 1196 1197 /*@C 1198 TaoADMMSetRegularizerObjectiveAndGradientRoutine - Sets the user-defined regularizer call-back function 1199 1200 Collective on tao 1201 1202 Input Parameters: 1203 + tao - the Tao context 1204 . func - function pointer for the regularizer value and gradient evaluation 1205 - ctx - user context for the regularizer 1206 1207 Level: advanced 1208 1209 .seealso: TAOADMM 1210 1211 @*/ 1212 PetscErrorCode TaoADMMSetRegularizerObjectiveAndGradientRoutine(Tao tao, PetscErrorCode (*func)(Tao, Vec, PetscReal*, Vec, void*), void *ctx) 1213 { 1214 TAO_ADMM *am = (TAO_ADMM*)tao->data; 1215 1216 PetscFunctionBegin; 1217 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 1218 am->regobjgradP = ctx; 1219 am->ops->regobjgrad = func; 1220 PetscFunctionReturn(0); 1221 } 1222 1223 /*@C 1224 TaoADMMSetRegularizerHessianRoutine - Sets the user-defined regularizer Hessian call-back 1225 function, to be used for subsolverZ. 1226 1227 Collective on tao 1228 1229 Input Parameters: 1230 + tao - the Tao context 1231 . H - user-created matrix for the Hessian of the regularization term 1232 . Hpre - user-created matrix for the preconditioner of Hessian of the regularization term 1233 . func - function pointer for the regularizer Hessian evaluation 1234 - ctx - user context for the regularizer Hessian 1235 1236 Level: advanced 1237 1238 .seealso: TAOADMM 1239 1240 @*/ 1241 PetscErrorCode TaoADMMSetRegularizerHessianRoutine(Tao tao, Mat H, Mat Hpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx) 1242 { 1243 TAO_ADMM *am = (TAO_ADMM*)tao->data; 1244 1245 PetscFunctionBegin; 1246 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 1247 if (H) { 1248 PetscValidHeaderSpecific(H,MAT_CLASSID,2); 1249 PetscCheckSameComm(tao,1,H,2); 1250 } 1251 if (Hpre) { 1252 PetscValidHeaderSpecific(Hpre,MAT_CLASSID,3); 1253 PetscCheckSameComm(tao,1,Hpre,3); 1254 } 1255 if (ctx) { 1256 am->reghessP = ctx; 1257 } 1258 if (func) { 1259 am->ops->reghess = func; 1260 } 1261 if (H) { 1262 PetscCall(PetscObjectReference((PetscObject)H)); 1263 PetscCall(MatDestroy(&am->Hz)); 1264 am->Hz = H; 1265 } 1266 if (Hpre) { 1267 PetscCall(PetscObjectReference((PetscObject)Hpre)); 1268 PetscCall(MatDestroy(&am->Hzpre)); 1269 am->Hzpre = Hpre; 1270 } 1271 PetscFunctionReturn(0); 1272 } 1273 1274 /*@ 1275 TaoGetADMMParentTao - Gets pointer to parent ADMM tao, used by inner subsolver. 1276 1277 Collective on tao 1278 1279 Input Parameter: 1280 . tao - the Tao context 1281 1282 Output Parameter: 1283 . admm_tao - the parent Tao context 1284 1285 Level: advanced 1286 1287 .seealso: TAOADMM 1288 1289 @*/ 1290 PetscErrorCode TaoGetADMMParentTao(Tao tao, Tao *admm_tao) 1291 { 1292 PetscFunctionBegin; 1293 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 1294 PetscCall(PetscObjectQuery((PetscObject)tao,"TaoGetADMMParentTao_ADMM", (PetscObject*) admm_tao)); 1295 PetscFunctionReturn(0); 1296 } 1297 1298 /*@ 1299 TaoADMMGetDualVector - Returns the dual vector associated with the current TAOADMM state 1300 1301 Not Collective 1302 1303 Input Parameter: 1304 . tao - the Tao context 1305 1306 Output Parameter: 1307 . Y - the current solution 1308 1309 Level: intermediate 1310 1311 .seealso: TAOADMM 1312 1313 @*/ 1314 PetscErrorCode TaoADMMGetDualVector(Tao tao, Vec *Y) 1315 { 1316 TAO_ADMM *am = (TAO_ADMM*)tao->data; 1317 1318 PetscFunctionBegin; 1319 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 1320 *Y = am->y; 1321 PetscFunctionReturn(0); 1322 } 1323 1324 /*@ 1325 TaoADMMSetRegularizerType - Set regularizer type for ADMM routine 1326 1327 Not Collective 1328 1329 Input Parameters: 1330 + tao - the Tao context 1331 - type - regularizer type 1332 1333 Options Database: 1334 . -tao_admm_regularizer_type <admm_regularizer_user,admm_regularizer_soft_thresh> - select the regularizer 1335 1336 Level: intermediate 1337 1338 .seealso: TaoADMMGetRegularizerType(), TaoADMMRegularizerType, TAOADMM 1339 @*/ 1340 PetscErrorCode TaoADMMSetRegularizerType(Tao tao, TaoADMMRegularizerType type) 1341 { 1342 PetscFunctionBegin; 1343 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 1344 PetscValidLogicalCollectiveEnum(tao,type,2); 1345 PetscTryMethod(tao,"TaoADMMSetRegularizerType_C",(Tao,TaoADMMRegularizerType),(tao,type)); 1346 PetscFunctionReturn(0); 1347 } 1348 1349 /*@ 1350 TaoADMMGetRegularizerType - Gets the type of regularizer routine for ADMM 1351 1352 Not Collective 1353 1354 Input Parameter: 1355 . tao - the Tao context 1356 1357 Output Parameter: 1358 . type - the type of regularizer 1359 1360 Level: intermediate 1361 1362 .seealso: TaoADMMSetRegularizerType(), TaoADMMRegularizerType, TAOADMM 1363 @*/ 1364 PetscErrorCode TaoADMMGetRegularizerType(Tao tao, TaoADMMRegularizerType *type) 1365 { 1366 PetscFunctionBegin; 1367 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 1368 PetscUseMethod(tao,"TaoADMMGetRegularizerType_C",(Tao,TaoADMMRegularizerType*),(tao,type)); 1369 PetscFunctionReturn(0); 1370 } 1371 1372 /*@ 1373 TaoADMMSetUpdateType - Set update routine for ADMM routine 1374 1375 Not Collective 1376 1377 Input Parameters: 1378 + tao - the Tao context 1379 - type - spectral parameter update type 1380 1381 Level: intermediate 1382 1383 .seealso: TaoADMMGetUpdateType(), TaoADMMUpdateType, TAOADMM 1384 @*/ 1385 PetscErrorCode TaoADMMSetUpdateType(Tao tao, TaoADMMUpdateType type) 1386 { 1387 PetscFunctionBegin; 1388 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 1389 PetscValidLogicalCollectiveEnum(tao,type,2); 1390 PetscTryMethod(tao,"TaoADMMSetUpdateType_C",(Tao,TaoADMMUpdateType),(tao,type)); 1391 PetscFunctionReturn(0); 1392 } 1393 1394 /*@ 1395 TaoADMMGetUpdateType - Gets the type of spectral penalty update routine for ADMM 1396 1397 Not Collective 1398 1399 Input Parameter: 1400 . tao - the Tao context 1401 1402 Output Parameter: 1403 . type - the type of spectral penalty update routine 1404 1405 Level: intermediate 1406 1407 .seealso: TaoADMMSetUpdateType(), TaoADMMUpdateType, TAOADMM 1408 @*/ 1409 PetscErrorCode TaoADMMGetUpdateType(Tao tao, TaoADMMUpdateType *type) 1410 { 1411 PetscFunctionBegin; 1412 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 1413 PetscUseMethod(tao,"TaoADMMGetUpdateType_C",(Tao,TaoADMMUpdateType*),(tao,type)); 1414 PetscFunctionReturn(0); 1415 } 1416