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 532 PetscFunctionBegin; 533 PetscOptionsHeadBegin(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. "); 534 PetscCall(PetscOptionsReal("-tao_admm_regularizer_coefficient","regularizer constant","",am->lambda,&am->lambda,NULL)); 535 PetscCall(PetscOptionsReal("-tao_admm_spectral_penalty","Constant for Augmented Lagrangian term.","",am->mu,&am->mu,NULL)); 536 PetscCall(PetscOptionsReal("-tao_admm_relaxation_parameter","x relaxation parameter for Z update.","",am->gamma,&am->gamma,NULL)); 537 PetscCall(PetscOptionsReal("-tao_admm_tolerance_update_factor","ADMM dynamic tolerance update factor.","",am->tol,&am->tol,NULL)); 538 PetscCall(PetscOptionsReal("-tao_admm_spectral_penalty_update_factor","ADMM spectral penalty update curvature safeguard value.","",am->orthval,&am->orthval,NULL)); 539 PetscCall(PetscOptionsReal("-tao_admm_minimum_spectral_penalty","Set ADMM minimum spectral penalty.","",am->mumin,&am->mumin,NULL)); 540 PetscCall(PetscOptionsEnum("-tao_admm_dual_update","Lagrangian dual update policy","TaoADMMUpdateType",TaoADMMUpdateTypes,(PetscEnum)am->update,(PetscEnum*)&am->update,NULL)); 541 PetscCall(PetscOptionsEnum("-tao_admm_regularizer_type","ADMM regularizer update rule","TaoADMMRegularizerType",TaoADMMRegularizerTypes,(PetscEnum)am->regswitch,(PetscEnum*)&am->regswitch,NULL)); 542 PetscOptionsHeadEnd(); 543 PetscCall(TaoSetFromOptions(am->subsolverX)); 544 if (am->regswitch != TAO_ADMM_REGULARIZER_SOFT_THRESH) { 545 PetscCall(TaoSetFromOptions(am->subsolverZ)); 546 } 547 PetscFunctionReturn(0); 548 } 549 550 static PetscErrorCode TaoView_ADMM(Tao tao,PetscViewer viewer) 551 { 552 TAO_ADMM *am = (TAO_ADMM*)tao->data; 553 554 PetscFunctionBegin; 555 PetscCall(PetscViewerASCIIPushTab(viewer)); 556 PetscCall(TaoView(am->subsolverX,viewer)); 557 PetscCall(TaoView(am->subsolverZ,viewer)); 558 PetscCall(PetscViewerASCIIPopTab(viewer)); 559 PetscFunctionReturn(0); 560 } 561 562 static PetscErrorCode TaoSetUp_ADMM(Tao tao) 563 { 564 TAO_ADMM *am = (TAO_ADMM*)tao->data; 565 PetscInt n,N,M; 566 567 PetscFunctionBegin; 568 PetscCall(VecGetLocalSize(tao->solution,&n)); 569 PetscCall(VecGetSize(tao->solution,&N)); 570 /* If Jacobian is given as NULL, it means Jacobian is identity matrix with size of solution vector */ 571 if (!am->JB) { 572 am->zJI = PETSC_TRUE; 573 PetscCall(MatCreateShell(PetscObjectComm((PetscObject)tao),n,n,PETSC_DETERMINE,PETSC_DETERMINE,NULL,&am->JB)); 574 PetscCall(MatShellSetOperation(am->JB,MATOP_MULT,(void (*)(void))JacobianIdentityB)); 575 PetscCall(MatShellSetOperation(am->JB,MATOP_MULT_TRANSPOSE,(void (*)(void))JacobianIdentityB)); 576 am->JBpre = am->JB; 577 } 578 if (!am->JA) { 579 am->xJI = PETSC_TRUE; 580 PetscCall(MatCreateShell(PetscObjectComm((PetscObject)tao),n,n,PETSC_DETERMINE,PETSC_DETERMINE,NULL,&am->JA)); 581 PetscCall(MatShellSetOperation(am->JA,MATOP_MULT,(void (*)(void))JacobianIdentity)); 582 PetscCall(MatShellSetOperation(am->JA,MATOP_MULT_TRANSPOSE,(void (*)(void))JacobianIdentity)); 583 am->JApre = am->JA; 584 } 585 PetscCall(MatCreateVecs(am->JA,NULL,&am->Ax)); 586 if (!tao->gradient) { 587 PetscCall(VecDuplicate(tao->solution,&tao->gradient)); 588 } 589 PetscCall(TaoSetSolution(am->subsolverX, tao->solution)); 590 if (!am->z) { 591 PetscCall(VecDuplicate(tao->solution,&am->z)); 592 PetscCall(VecSet(am->z,0.0)); 593 } 594 PetscCall(TaoSetSolution(am->subsolverZ, am->z)); 595 if (!am->workLeft) { 596 PetscCall(VecDuplicate(tao->solution,&am->workLeft)); 597 } 598 if (!am->Axold) { 599 PetscCall(VecDuplicate(am->Ax,&am->Axold)); 600 } 601 if (!am->workJacobianRight) { 602 PetscCall(VecDuplicate(am->Ax,&am->workJacobianRight)); 603 } 604 if (!am->workJacobianRight2) { 605 PetscCall(VecDuplicate(am->Ax,&am->workJacobianRight2)); 606 } 607 if (!am->Bz) { 608 PetscCall(VecDuplicate(am->Ax,&am->Bz)); 609 } 610 if (!am->Bzold) { 611 PetscCall(VecDuplicate(am->Ax,&am->Bzold)); 612 } 613 if (!am->Bz0) { 614 PetscCall(VecDuplicate(am->Ax,&am->Bz0)); 615 } 616 if (!am->y) { 617 PetscCall(VecDuplicate(am->Ax,&am->y)); 618 PetscCall(VecSet(am->y,0.0)); 619 } 620 if (!am->yold) { 621 PetscCall(VecDuplicate(am->Ax,&am->yold)); 622 PetscCall(VecSet(am->yold,0.0)); 623 } 624 if (!am->y0) { 625 PetscCall(VecDuplicate(am->Ax,&am->y0)); 626 PetscCall(VecSet(am->y0,0.0)); 627 } 628 if (!am->yhat) { 629 PetscCall(VecDuplicate(am->Ax,&am->yhat)); 630 PetscCall(VecSet(am->yhat,0.0)); 631 } 632 if (!am->yhatold) { 633 PetscCall(VecDuplicate(am->Ax,&am->yhatold)); 634 PetscCall(VecSet(am->yhatold,0.0)); 635 } 636 if (!am->residual) { 637 PetscCall(VecDuplicate(am->Ax,&am->residual)); 638 PetscCall(VecSet(am->residual,0.0)); 639 } 640 if (!am->constraint) { 641 am->constraint = NULL; 642 } else { 643 PetscCall(VecGetSize(am->constraint,&M)); 644 PetscCheck(M == N,PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_WRONGSTATE,"Solution vector and constraint vector must be of same size!"); 645 } 646 647 /* Save changed tao tolerance for adaptive tolerance */ 648 if (tao->gatol_changed) { 649 am->gatol_admm = tao->gatol; 650 } 651 if (tao->catol_changed) { 652 am->catol_admm = tao->catol; 653 } 654 655 /*Update spectral and dual elements to X subsolver */ 656 PetscCall(TaoSetObjectiveAndGradient(am->subsolverX, NULL, SubObjGradUpdate, tao)); 657 PetscCall(TaoSetJacobianEqualityRoutine(am->subsolverX,am->JA,am->JApre, am->ops->misfitjac, am->misfitjacobianP)); 658 PetscCall(TaoSetJacobianEqualityRoutine(am->subsolverZ,am->JB,am->JBpre, am->ops->regjac, am->regjacobianP)); 659 PetscFunctionReturn(0); 660 } 661 662 static PetscErrorCode TaoDestroy_ADMM(Tao tao) 663 { 664 TAO_ADMM *am = (TAO_ADMM*)tao->data; 665 666 PetscFunctionBegin; 667 PetscCall(VecDestroy(&am->z)); 668 PetscCall(VecDestroy(&am->Ax)); 669 PetscCall(VecDestroy(&am->Axold)); 670 PetscCall(VecDestroy(&am->Bz)); 671 PetscCall(VecDestroy(&am->Bzold)); 672 PetscCall(VecDestroy(&am->Bz0)); 673 PetscCall(VecDestroy(&am->residual)); 674 PetscCall(VecDestroy(&am->y)); 675 PetscCall(VecDestroy(&am->yold)); 676 PetscCall(VecDestroy(&am->y0)); 677 PetscCall(VecDestroy(&am->yhat)); 678 PetscCall(VecDestroy(&am->yhatold)); 679 PetscCall(VecDestroy(&am->workLeft)); 680 PetscCall(VecDestroy(&am->workJacobianRight)); 681 PetscCall(VecDestroy(&am->workJacobianRight2)); 682 683 PetscCall(MatDestroy(&am->JA)); 684 PetscCall(MatDestroy(&am->JB)); 685 if (!am->xJI) { 686 PetscCall(MatDestroy(&am->JApre)); 687 } 688 if (!am->zJI) { 689 PetscCall(MatDestroy(&am->JBpre)); 690 } 691 if (am->Hx) { 692 PetscCall(MatDestroy(&am->Hx)); 693 PetscCall(MatDestroy(&am->Hxpre)); 694 } 695 if (am->Hz) { 696 PetscCall(MatDestroy(&am->Hz)); 697 PetscCall(MatDestroy(&am->Hzpre)); 698 } 699 PetscCall(MatDestroy(&am->ATA)); 700 PetscCall(MatDestroy(&am->BTB)); 701 PetscCall(TaoDestroy(&am->subsolverX)); 702 PetscCall(TaoDestroy(&am->subsolverZ)); 703 am->parent = NULL; 704 PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMSetRegularizerType_C",NULL)); 705 PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMGetRegularizerType_C",NULL)); 706 PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMSetUpdateType_C",NULL)); 707 PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMGetUpdateType_C",NULL)); 708 PetscCall(PetscFree(tao->data)); 709 PetscFunctionReturn(0); 710 } 711 712 /*MC 713 714 TAOADMM - Alternating direction method of multipliers method fo solving linear problems with 715 constraints. in a min_x f(x) + g(z) s.t. Ax+Bz=c. 716 This algorithm employs two sub Tao solvers, of which type can be specified 717 by the user. User need to provide ObjectiveAndGradient routine, and/or HessianRoutine for both subsolvers. 718 Hessians can be given boolean flag determining whether they change with respect to a input vector. This can be set via 719 TaoADMMSet{Misfit,Regularizer}HessianChangeStatus. 720 Second subsolver does support TAOSHELL. It should be noted that L1-norm is used for objective value for TAOSHELL type. 721 There is option to set regularizer option, and currently soft-threshold is implemented. For spectral penalty update, 722 currently there are basic option and adaptive option. 723 Constraint is set at Ax+Bz=c, and A and B can be set with TaoADMMSet{Misfit,Regularizer}ConstraintJacobian. 724 c can be set with TaoADMMSetConstraintVectorRHS. 725 The user can also provide regularizer weight for second subsolver. 726 727 References: 728 . * - Xu, Zheng and Figueiredo, Mario A. T. and Yuan, Xiaoming and Studer, Christoph and Goldstein, Tom 729 "Adaptive Relaxed ADMM: Convergence Theory and Practical Implementation" 730 The IEEE Conference on Computer Vision and Pattern Recognition (CVPR), July, 2017. 731 732 Options Database Keys: 733 + -tao_admm_regularizer_coefficient - regularizer constant (default 1.e-6) 734 . -tao_admm_spectral_penalty - Constant for Augmented Lagrangian term (default 1.) 735 . -tao_admm_relaxation_parameter - relaxation parameter for Z update (default 1.) 736 . -tao_admm_tolerance_update_factor - ADMM dynamic tolerance update factor (default 1.e-12) 737 . -tao_admm_spectral_penalty_update_factor - ADMM spectral penalty update curvature safeguard value (default 0.2) 738 . -tao_admm_minimum_spectral_penalty - Set ADMM minimum spectral penalty (default 0) 739 . -tao_admm_dual_update - Lagrangian dual update policy ("basic","adaptive","adaptive-relaxed") (default "basic") 740 - -tao_admm_regularizer_type - ADMM regularizer update rule ("user","soft-threshold") (default "soft-threshold") 741 742 Level: beginner 743 744 .seealso: `TaoADMMSetMisfitHessianChangeStatus()`, `TaoADMMSetRegHessianChangeStatus()`, `TaoADMMGetSpectralPenalty()`, 745 `TaoADMMGetMisfitSubsolver()`, `TaoADMMGetRegularizationSubsolver()`, `TaoADMMSetConstraintVectorRHS()`, 746 `TaoADMMSetMinimumSpectralPenalty()`, `TaoADMMSetRegularizerCoefficient()`, 747 `TaoADMMSetRegularizerConstraintJacobian()`, `TaoADMMSetMisfitConstraintJacobian()`, 748 `TaoADMMSetMisfitObjectiveAndGradientRoutine()`, `TaoADMMSetMisfitHessianRoutine()`, 749 `TaoADMMSetRegularizerObjectiveAndGradientRoutine()`, `TaoADMMSetRegularizerHessianRoutine()`, 750 `TaoGetADMMParentTao()`, `TaoADMMGetDualVector()`, `TaoADMMSetRegularizerType()`, 751 `TaoADMMGetRegularizerType()`, `TaoADMMSetUpdateType()`, `TaoADMMGetUpdateType()` 752 M*/ 753 754 PETSC_EXTERN PetscErrorCode TaoCreate_ADMM(Tao tao) 755 { 756 TAO_ADMM *am; 757 758 PetscFunctionBegin; 759 PetscCall(PetscNewLog(tao,&am)); 760 761 tao->ops->destroy = TaoDestroy_ADMM; 762 tao->ops->setup = TaoSetUp_ADMM; 763 tao->ops->setfromoptions = TaoSetFromOptions_ADMM; 764 tao->ops->view = TaoView_ADMM; 765 tao->ops->solve = TaoSolve_ADMM; 766 767 tao->data = (void*)am; 768 am->l1epsilon = 1e-6; 769 am->lambda = 1e-4; 770 am->mu = 1.; 771 am->muold = 0.; 772 am->mueps = PETSC_MACHINE_EPSILON; 773 am->mumin = 0.; 774 am->orthval = 0.2; 775 am->T = 2; 776 am->parent = tao; 777 am->update = TAO_ADMM_UPDATE_BASIC; 778 am->regswitch = TAO_ADMM_REGULARIZER_SOFT_THRESH; 779 am->tol = PETSC_SMALL; 780 am->const_norm = 0; 781 am->resnorm = 0; 782 am->dualres = 0; 783 am->ops->regobjgrad = NULL; 784 am->ops->reghess = NULL; 785 am->gamma = 1; 786 am->regobjgradP = NULL; 787 am->reghessP = NULL; 788 am->gatol_admm = 1e-8; 789 am->catol_admm = 0; 790 am->Hxchange = PETSC_TRUE; 791 am->Hzchange = PETSC_TRUE; 792 am->Hzbool = PETSC_TRUE; 793 am->Hxbool = PETSC_TRUE; 794 795 PetscCall(TaoCreate(PetscObjectComm((PetscObject)tao),&am->subsolverX)); 796 PetscCall(TaoSetOptionsPrefix(am->subsolverX,"misfit_")); 797 PetscCall(PetscObjectIncrementTabLevel((PetscObject)am->subsolverX,(PetscObject)tao,1)); 798 PetscCall(TaoCreate(PetscObjectComm((PetscObject)tao),&am->subsolverZ)); 799 PetscCall(TaoSetOptionsPrefix(am->subsolverZ,"reg_")); 800 PetscCall(PetscObjectIncrementTabLevel((PetscObject)am->subsolverZ,(PetscObject)tao,1)); 801 802 PetscCall(TaoSetType(am->subsolverX,TAONLS)); 803 PetscCall(TaoSetType(am->subsolverZ,TAONLS)); 804 PetscCall(PetscObjectCompose((PetscObject)am->subsolverX,"TaoGetADMMParentTao_ADMM", (PetscObject) tao)); 805 PetscCall(PetscObjectCompose((PetscObject)am->subsolverZ,"TaoGetADMMParentTao_ADMM", (PetscObject) tao)); 806 PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMSetRegularizerType_C",TaoADMMSetRegularizerType_ADMM)); 807 PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMGetRegularizerType_C",TaoADMMGetRegularizerType_ADMM)); 808 PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMSetUpdateType_C",TaoADMMSetUpdateType_ADMM)); 809 PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoADMMGetUpdateType_C",TaoADMMGetUpdateType_ADMM)); 810 PetscFunctionReturn(0); 811 } 812 813 /*@ 814 TaoADMMSetMisfitHessianChangeStatus - Set boolean that determines whether Hessian matrix of misfit subsolver changes with respect to input vector. 815 816 Collective on Tao 817 818 Input Parameters: 819 + tao - the Tao solver context. 820 - b - the Hessian matrix change status boolean, PETSC_FALSE when the Hessian matrix does not change, TRUE otherwise. 821 822 Level: advanced 823 824 .seealso: `TAOADMM` 825 826 @*/ 827 PetscErrorCode TaoADMMSetMisfitHessianChangeStatus(Tao tao, PetscBool b) 828 { 829 TAO_ADMM *am = (TAO_ADMM*)tao->data; 830 831 PetscFunctionBegin; 832 am->Hxchange = b; 833 PetscFunctionReturn(0); 834 } 835 836 /*@ 837 TaoADMMSetRegHessianChangeStatus - Set boolean that determines whether Hessian matrix of regularization subsolver changes with respect to input vector. 838 839 Collective on Tao 840 841 Input Parameters: 842 + tao - the Tao solver context 843 - b - the Hessian matrix change status boolean, PETSC_FALSE when the Hessian matrix does not change, TRUE otherwise. 844 845 Level: advanced 846 847 .seealso: `TAOADMM` 848 849 @*/ 850 PetscErrorCode TaoADMMSetRegHessianChangeStatus(Tao tao, PetscBool b) 851 { 852 TAO_ADMM *am = (TAO_ADMM*)tao->data; 853 854 PetscFunctionBegin; 855 am->Hzchange = b; 856 PetscFunctionReturn(0); 857 } 858 859 /*@ 860 TaoADMMSetSpectralPenalty - Set the spectral penalty (mu) value 861 862 Collective on Tao 863 864 Input Parameters: 865 + tao - the Tao solver context 866 - mu - spectral penalty 867 868 Level: advanced 869 870 .seealso: `TaoADMMSetMinimumSpectralPenalty()`, `TAOADMM` 871 @*/ 872 PetscErrorCode TaoADMMSetSpectralPenalty(Tao tao, PetscReal mu) 873 { 874 TAO_ADMM *am = (TAO_ADMM*)tao->data; 875 876 PetscFunctionBegin; 877 am->mu = mu; 878 PetscFunctionReturn(0); 879 } 880 881 /*@ 882 TaoADMMGetSpectralPenalty - Get the spectral penalty (mu) value 883 884 Collective on Tao 885 886 Input Parameter: 887 . tao - the Tao solver context 888 889 Output Parameter: 890 . mu - spectral penalty 891 892 Level: advanced 893 894 .seealso: `TaoADMMSetMinimumSpectralPenalty()`, `TaoADMMSetSpectralPenalty()`, `TAOADMM` 895 @*/ 896 PetscErrorCode TaoADMMGetSpectralPenalty(Tao tao, PetscReal *mu) 897 { 898 TAO_ADMM *am = (TAO_ADMM*)tao->data; 899 900 PetscFunctionBegin; 901 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 902 PetscValidRealPointer(mu,2); 903 *mu = am->mu; 904 PetscFunctionReturn(0); 905 } 906 907 /*@ 908 TaoADMMGetMisfitSubsolver - Get the pointer to the misfit subsolver inside ADMM 909 910 Collective on Tao 911 912 Input Parameter: 913 . tao - the Tao solver context 914 915 Output Parameter: 916 . misfit - the Tao subsolver context 917 918 Level: advanced 919 920 .seealso: `TAOADMM` 921 922 @*/ 923 PetscErrorCode TaoADMMGetMisfitSubsolver(Tao tao, Tao *misfit) 924 { 925 TAO_ADMM *am = (TAO_ADMM*)tao->data; 926 927 PetscFunctionBegin; 928 *misfit = am->subsolverX; 929 PetscFunctionReturn(0); 930 } 931 932 /*@ 933 TaoADMMGetRegularizationSubsolver - Get the pointer to the regularization subsolver inside ADMM 934 935 Collective on Tao 936 937 Input Parameter: 938 . tao - the Tao solver context 939 940 Output Parameter: 941 . reg - the Tao subsolver context 942 943 Level: advanced 944 945 .seealso: `TAOADMM` 946 947 @*/ 948 PetscErrorCode TaoADMMGetRegularizationSubsolver(Tao tao, Tao *reg) 949 { 950 TAO_ADMM *am = (TAO_ADMM*)tao->data; 951 952 PetscFunctionBegin; 953 *reg = am->subsolverZ; 954 PetscFunctionReturn(0); 955 } 956 957 /*@ 958 TaoADMMSetConstraintVectorRHS - Set the RHS constraint vector for ADMM 959 960 Collective on Tao 961 962 Input Parameters: 963 + tao - the Tao solver context 964 - c - RHS vector 965 966 Level: advanced 967 968 .seealso: `TAOADMM` 969 970 @*/ 971 PetscErrorCode TaoADMMSetConstraintVectorRHS(Tao tao, Vec c) 972 { 973 TAO_ADMM *am = (TAO_ADMM*)tao->data; 974 975 PetscFunctionBegin; 976 am->constraint = c; 977 PetscFunctionReturn(0); 978 } 979 980 /*@ 981 TaoADMMSetMinimumSpectralPenalty - Set the minimum value for the spectral penalty 982 983 Collective on Tao 984 985 Input Parameters: 986 + tao - the Tao solver context 987 - mu - minimum spectral penalty value 988 989 Level: advanced 990 991 .seealso: `TaoADMMGetSpectralPenalty()`, `TAOADMM` 992 @*/ 993 PetscErrorCode TaoADMMSetMinimumSpectralPenalty(Tao tao, PetscReal mu) 994 { 995 TAO_ADMM *am = (TAO_ADMM*)tao->data; 996 997 PetscFunctionBegin; 998 am->mumin= mu; 999 PetscFunctionReturn(0); 1000 } 1001 1002 /*@ 1003 TaoADMMSetRegularizerCoefficient - Set the regularization coefficient lambda for L1 norm regularization case 1004 1005 Collective on Tao 1006 1007 Input Parameters: 1008 + tao - the Tao solver context 1009 - lambda - L1-norm regularizer coefficient 1010 1011 Level: advanced 1012 1013 .seealso: `TaoADMMSetMisfitConstraintJacobian()`, `TaoADMMSetRegularizerConstraintJacobian()`, `TAOADMM` 1014 1015 @*/ 1016 PetscErrorCode TaoADMMSetRegularizerCoefficient(Tao tao, PetscReal lambda) 1017 { 1018 TAO_ADMM *am = (TAO_ADMM*)tao->data; 1019 1020 PetscFunctionBegin; 1021 am->lambda = lambda; 1022 PetscFunctionReturn(0); 1023 } 1024 1025 /*@C 1026 TaoADMMSetMisfitConstraintJacobian - Set the constraint matrix B for the ADMM algorithm. Matrix B constrains the z variable. 1027 1028 Collective on Tao 1029 1030 Input Parameters: 1031 + tao - the Tao solver context 1032 . J - user-created regularizer constraint Jacobian matrix 1033 . Jpre - user-created regularizer Jacobian constraint preconditioner matrix 1034 . func - function pointer for the regularizer constraint Jacobian update function 1035 - ctx - user context for the regularizer Hessian 1036 1037 Level: advanced 1038 1039 .seealso: `TaoADMMSetRegularizerCoefficient()`, `TaoADMMSetRegularizerConstraintJacobian()`, `TAOADMM` 1040 1041 @*/ 1042 PetscErrorCode TaoADMMSetMisfitConstraintJacobian(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx) 1043 { 1044 TAO_ADMM *am = (TAO_ADMM*)tao->data; 1045 1046 PetscFunctionBegin; 1047 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 1048 if (J) { 1049 PetscValidHeaderSpecific(J,MAT_CLASSID,2); 1050 PetscCheckSameComm(tao,1,J,2); 1051 } 1052 if (Jpre) { 1053 PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3); 1054 PetscCheckSameComm(tao,1,Jpre,3); 1055 } 1056 if (ctx) am->misfitjacobianP = ctx; 1057 if (func) am->ops->misfitjac = func; 1058 1059 if (J) { 1060 PetscCall(PetscObjectReference((PetscObject)J)); 1061 PetscCall(MatDestroy(&am->JA)); 1062 am->JA = J; 1063 } 1064 if (Jpre) { 1065 PetscCall(PetscObjectReference((PetscObject)Jpre)); 1066 PetscCall(MatDestroy(&am->JApre)); 1067 am->JApre = Jpre; 1068 } 1069 PetscFunctionReturn(0); 1070 } 1071 1072 /*@C 1073 TaoADMMSetRegularizerConstraintJacobian - Set the constraint matrix B for ADMM algorithm. Matrix B constraints z variable. 1074 1075 Collective on Tao 1076 1077 Input Parameters: 1078 + tao - the Tao solver context 1079 . J - user-created regularizer constraint Jacobian matrix 1080 . Jpre - user-created regularizer Jacobian constraint preconditioner matrix 1081 . func - function pointer for the regularizer constraint Jacobian update function 1082 - ctx - user context for the regularizer Hessian 1083 1084 Level: advanced 1085 1086 .seealso: `TaoADMMSetRegularizerCoefficient()`, `TaoADMMSetMisfitConstraintJacobian()`, `TAOADMM` 1087 1088 @*/ 1089 PetscErrorCode TaoADMMSetRegularizerConstraintJacobian(Tao tao, Mat J, Mat Jpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx) 1090 { 1091 TAO_ADMM *am = (TAO_ADMM*)tao->data; 1092 1093 PetscFunctionBegin; 1094 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 1095 if (J) { 1096 PetscValidHeaderSpecific(J,MAT_CLASSID,2); 1097 PetscCheckSameComm(tao,1,J,2); 1098 } 1099 if (Jpre) { 1100 PetscValidHeaderSpecific(Jpre,MAT_CLASSID,3); 1101 PetscCheckSameComm(tao,1,Jpre,3); 1102 } 1103 if (ctx) am->regjacobianP = ctx; 1104 if (func) am->ops->regjac = func; 1105 1106 if (J) { 1107 PetscCall(PetscObjectReference((PetscObject)J)); 1108 PetscCall(MatDestroy(&am->JB)); 1109 am->JB = J; 1110 } 1111 if (Jpre) { 1112 PetscCall(PetscObjectReference((PetscObject)Jpre)); 1113 PetscCall(MatDestroy(&am->JBpre)); 1114 am->JBpre = Jpre; 1115 } 1116 PetscFunctionReturn(0); 1117 } 1118 1119 /*@C 1120 TaoADMMSetMisfitObjectiveAndGradientRoutine - Sets the user-defined misfit call-back function 1121 1122 Collective on tao 1123 1124 Input Parameters: 1125 + tao - the Tao context 1126 . func - function pointer for the misfit value and gradient evaluation 1127 - ctx - user context for the misfit 1128 1129 Level: advanced 1130 1131 .seealso: `TAOADMM` 1132 1133 @*/ 1134 PetscErrorCode TaoADMMSetMisfitObjectiveAndGradientRoutine(Tao tao, PetscErrorCode (*func)(Tao, Vec, PetscReal*, Vec, void*), void *ctx) 1135 { 1136 TAO_ADMM *am = (TAO_ADMM*)tao->data; 1137 1138 PetscFunctionBegin; 1139 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 1140 am->misfitobjgradP = ctx; 1141 am->ops->misfitobjgrad = func; 1142 PetscFunctionReturn(0); 1143 } 1144 1145 /*@C 1146 TaoADMMSetMisfitHessianRoutine - Sets the user-defined misfit Hessian call-back 1147 function into the algorithm, to be used for subsolverX. 1148 1149 Collective on tao 1150 1151 Input Parameters: 1152 + tao - the Tao context 1153 . H - user-created matrix for the Hessian of the misfit term 1154 . Hpre - user-created matrix for the preconditioner of Hessian of the misfit term 1155 . func - function pointer for the misfit Hessian evaluation 1156 - ctx - user context for the misfit Hessian 1157 1158 Level: advanced 1159 1160 .seealso: `TAOADMM` 1161 1162 @*/ 1163 PetscErrorCode TaoADMMSetMisfitHessianRoutine(Tao tao, Mat H, Mat Hpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx) 1164 { 1165 TAO_ADMM *am = (TAO_ADMM*)tao->data; 1166 1167 PetscFunctionBegin; 1168 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 1169 if (H) { 1170 PetscValidHeaderSpecific(H,MAT_CLASSID,2); 1171 PetscCheckSameComm(tao,1,H,2); 1172 } 1173 if (Hpre) { 1174 PetscValidHeaderSpecific(Hpre,MAT_CLASSID,3); 1175 PetscCheckSameComm(tao,1,Hpre,3); 1176 } 1177 if (ctx) { 1178 am->misfithessP = ctx; 1179 } 1180 if (func) { 1181 am->ops->misfithess = func; 1182 } 1183 if (H) { 1184 PetscCall(PetscObjectReference((PetscObject)H)); 1185 PetscCall(MatDestroy(&am->Hx)); 1186 am->Hx = H; 1187 } 1188 if (Hpre) { 1189 PetscCall(PetscObjectReference((PetscObject)Hpre)); 1190 PetscCall(MatDestroy(&am->Hxpre)); 1191 am->Hxpre = Hpre; 1192 } 1193 PetscFunctionReturn(0); 1194 } 1195 1196 /*@C 1197 TaoADMMSetRegularizerObjectiveAndGradientRoutine - Sets the user-defined regularizer call-back function 1198 1199 Collective on tao 1200 1201 Input Parameters: 1202 + tao - the Tao context 1203 . func - function pointer for the regularizer value and gradient evaluation 1204 - ctx - user context for the regularizer 1205 1206 Level: advanced 1207 1208 .seealso: `TAOADMM` 1209 1210 @*/ 1211 PetscErrorCode TaoADMMSetRegularizerObjectiveAndGradientRoutine(Tao tao, PetscErrorCode (*func)(Tao, Vec, PetscReal*, Vec, void*), void *ctx) 1212 { 1213 TAO_ADMM *am = (TAO_ADMM*)tao->data; 1214 1215 PetscFunctionBegin; 1216 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 1217 am->regobjgradP = ctx; 1218 am->ops->regobjgrad = func; 1219 PetscFunctionReturn(0); 1220 } 1221 1222 /*@C 1223 TaoADMMSetRegularizerHessianRoutine - Sets the user-defined regularizer Hessian call-back 1224 function, to be used for subsolverZ. 1225 1226 Collective on tao 1227 1228 Input Parameters: 1229 + tao - the Tao context 1230 . H - user-created matrix for the Hessian of the regularization term 1231 . Hpre - user-created matrix for the preconditioner of Hessian of the regularization term 1232 . func - function pointer for the regularizer Hessian evaluation 1233 - ctx - user context for the regularizer Hessian 1234 1235 Level: advanced 1236 1237 .seealso: `TAOADMM` 1238 1239 @*/ 1240 PetscErrorCode TaoADMMSetRegularizerHessianRoutine(Tao tao, Mat H, Mat Hpre, PetscErrorCode (*func)(Tao, Vec, Mat, Mat, void*), void *ctx) 1241 { 1242 TAO_ADMM *am = (TAO_ADMM*)tao->data; 1243 1244 PetscFunctionBegin; 1245 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 1246 if (H) { 1247 PetscValidHeaderSpecific(H,MAT_CLASSID,2); 1248 PetscCheckSameComm(tao,1,H,2); 1249 } 1250 if (Hpre) { 1251 PetscValidHeaderSpecific(Hpre,MAT_CLASSID,3); 1252 PetscCheckSameComm(tao,1,Hpre,3); 1253 } 1254 if (ctx) { 1255 am->reghessP = ctx; 1256 } 1257 if (func) { 1258 am->ops->reghess = func; 1259 } 1260 if (H) { 1261 PetscCall(PetscObjectReference((PetscObject)H)); 1262 PetscCall(MatDestroy(&am->Hz)); 1263 am->Hz = H; 1264 } 1265 if (Hpre) { 1266 PetscCall(PetscObjectReference((PetscObject)Hpre)); 1267 PetscCall(MatDestroy(&am->Hzpre)); 1268 am->Hzpre = Hpre; 1269 } 1270 PetscFunctionReturn(0); 1271 } 1272 1273 /*@ 1274 TaoGetADMMParentTao - Gets pointer to parent ADMM tao, used by inner subsolver. 1275 1276 Collective on tao 1277 1278 Input Parameter: 1279 . tao - the Tao context 1280 1281 Output Parameter: 1282 . admm_tao - the parent Tao context 1283 1284 Level: advanced 1285 1286 .seealso: `TAOADMM` 1287 1288 @*/ 1289 PetscErrorCode TaoGetADMMParentTao(Tao tao, Tao *admm_tao) 1290 { 1291 PetscFunctionBegin; 1292 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 1293 PetscCall(PetscObjectQuery((PetscObject)tao,"TaoGetADMMParentTao_ADMM", (PetscObject*) admm_tao)); 1294 PetscFunctionReturn(0); 1295 } 1296 1297 /*@ 1298 TaoADMMGetDualVector - Returns the dual vector associated with the current TAOADMM state 1299 1300 Not Collective 1301 1302 Input Parameter: 1303 . tao - the Tao context 1304 1305 Output Parameter: 1306 . Y - the current solution 1307 1308 Level: intermediate 1309 1310 .seealso: `TAOADMM` 1311 1312 @*/ 1313 PetscErrorCode TaoADMMGetDualVector(Tao tao, Vec *Y) 1314 { 1315 TAO_ADMM *am = (TAO_ADMM*)tao->data; 1316 1317 PetscFunctionBegin; 1318 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 1319 *Y = am->y; 1320 PetscFunctionReturn(0); 1321 } 1322 1323 /*@ 1324 TaoADMMSetRegularizerType - Set regularizer type for ADMM routine 1325 1326 Not Collective 1327 1328 Input Parameters: 1329 + tao - the Tao context 1330 - type - regularizer type 1331 1332 Options Database: 1333 . -tao_admm_regularizer_type <admm_regularizer_user,admm_regularizer_soft_thresh> - select the regularizer 1334 1335 Level: intermediate 1336 1337 .seealso: `TaoADMMGetRegularizerType()`, `TaoADMMRegularizerType`, `TAOADMM` 1338 @*/ 1339 PetscErrorCode TaoADMMSetRegularizerType(Tao tao, TaoADMMRegularizerType type) 1340 { 1341 PetscFunctionBegin; 1342 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 1343 PetscValidLogicalCollectiveEnum(tao,type,2); 1344 PetscTryMethod(tao,"TaoADMMSetRegularizerType_C",(Tao,TaoADMMRegularizerType),(tao,type)); 1345 PetscFunctionReturn(0); 1346 } 1347 1348 /*@ 1349 TaoADMMGetRegularizerType - Gets the type of regularizer routine for ADMM 1350 1351 Not Collective 1352 1353 Input Parameter: 1354 . tao - the Tao context 1355 1356 Output Parameter: 1357 . type - the type of regularizer 1358 1359 Level: intermediate 1360 1361 .seealso: `TaoADMMSetRegularizerType()`, `TaoADMMRegularizerType`, `TAOADMM` 1362 @*/ 1363 PetscErrorCode TaoADMMGetRegularizerType(Tao tao, TaoADMMRegularizerType *type) 1364 { 1365 PetscFunctionBegin; 1366 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 1367 PetscUseMethod(tao,"TaoADMMGetRegularizerType_C",(Tao,TaoADMMRegularizerType*),(tao,type)); 1368 PetscFunctionReturn(0); 1369 } 1370 1371 /*@ 1372 TaoADMMSetUpdateType - Set update routine for ADMM routine 1373 1374 Not Collective 1375 1376 Input Parameters: 1377 + tao - the Tao context 1378 - type - spectral parameter update type 1379 1380 Level: intermediate 1381 1382 .seealso: `TaoADMMGetUpdateType()`, `TaoADMMUpdateType`, `TAOADMM` 1383 @*/ 1384 PetscErrorCode TaoADMMSetUpdateType(Tao tao, TaoADMMUpdateType type) 1385 { 1386 PetscFunctionBegin; 1387 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 1388 PetscValidLogicalCollectiveEnum(tao,type,2); 1389 PetscTryMethod(tao,"TaoADMMSetUpdateType_C",(Tao,TaoADMMUpdateType),(tao,type)); 1390 PetscFunctionReturn(0); 1391 } 1392 1393 /*@ 1394 TaoADMMGetUpdateType - Gets the type of spectral penalty update routine for ADMM 1395 1396 Not Collective 1397 1398 Input Parameter: 1399 . tao - the Tao context 1400 1401 Output Parameter: 1402 . type - the type of spectral penalty update routine 1403 1404 Level: intermediate 1405 1406 .seealso: `TaoADMMSetUpdateType()`, `TaoADMMUpdateType`, `TAOADMM` 1407 @*/ 1408 PetscErrorCode TaoADMMGetUpdateType(Tao tao, TaoADMMUpdateType *type) 1409 { 1410 PetscFunctionBegin; 1411 PetscValidHeaderSpecific(tao,TAO_CLASSID,1); 1412 PetscUseMethod(tao,"TaoADMMGetUpdateType_C",(Tao,TaoADMMUpdateType*),(tao,type)); 1413 PetscFunctionReturn(0); 1414 } 1415