1*c4762a1bSJed Brown static char help[] = "Simple example to test separable objective optimizers.\n"; 2*c4762a1bSJed Brown 3*c4762a1bSJed Brown #include <petsc.h> 4*c4762a1bSJed Brown #include <petsctao.h> 5*c4762a1bSJed Brown #include <petscvec.h> 6*c4762a1bSJed Brown #include <petscmath.h> 7*c4762a1bSJed Brown 8*c4762a1bSJed Brown #define NWORKLEFT 4 9*c4762a1bSJed Brown #define NWORKRIGHT 12 10*c4762a1bSJed Brown 11*c4762a1bSJed Brown typedef struct _UserCtx 12*c4762a1bSJed Brown { 13*c4762a1bSJed Brown PetscInt m; /* The row dimension of F */ 14*c4762a1bSJed Brown PetscInt n; /* The column dimension of F */ 15*c4762a1bSJed Brown PetscInt matops; /* Matrix format. 0 for stencil, 1 for random */ 16*c4762a1bSJed Brown PetscInt iter; /* Numer of iterations for ADMM */ 17*c4762a1bSJed Brown PetscReal hStart; /* Starting point for Taylor test */ 18*c4762a1bSJed Brown PetscReal hFactor; /* Taylor test step factor */ 19*c4762a1bSJed Brown PetscReal hMin; /* Taylor test end goal */ 20*c4762a1bSJed Brown PetscReal alpha; /* regularization constant applied to || x ||_p */ 21*c4762a1bSJed Brown PetscReal eps; /* small constant for approximating gradient of || x ||_1 */ 22*c4762a1bSJed Brown PetscReal mu; /* the augmented Lagrangian term in ADMM */ 23*c4762a1bSJed Brown PetscReal abstol; 24*c4762a1bSJed Brown PetscReal reltol; 25*c4762a1bSJed Brown Mat F; /* matrix in least squares component $(1/2) * || F x - d ||_2^2$ */ 26*c4762a1bSJed Brown Mat W; /* Workspace matrix. ATA */ 27*c4762a1bSJed Brown Mat Hm; /* Hessian Misfit*/ 28*c4762a1bSJed Brown Mat Hr; /* Hessian Reg*/ 29*c4762a1bSJed Brown Vec d; /* RHS in least squares component $(1/2) * || F x - d ||_2^2$ */ 30*c4762a1bSJed Brown Vec workLeft[NWORKLEFT]; /* Workspace for temporary vec */ 31*c4762a1bSJed Brown Vec workRight[NWORKRIGHT]; /* Workspace for temporary vec */ 32*c4762a1bSJed Brown NormType p; 33*c4762a1bSJed Brown PetscRandom rctx; 34*c4762a1bSJed Brown PetscBool taylor; /* Flag to determine whether to run Taylor test or not */ 35*c4762a1bSJed Brown PetscBool use_admm; /* Flag to determine whether to run Taylor test or not */ 36*c4762a1bSJed Brown }* UserCtx; 37*c4762a1bSJed Brown 38*c4762a1bSJed Brown static PetscErrorCode CreateRHS(UserCtx ctx) 39*c4762a1bSJed Brown { 40*c4762a1bSJed Brown PetscErrorCode ierr; 41*c4762a1bSJed Brown 42*c4762a1bSJed Brown PetscFunctionBegin; 43*c4762a1bSJed Brown /* build the rhs d in ctx */ 44*c4762a1bSJed Brown ierr = VecCreate(PETSC_COMM_WORLD,&(ctx->d));CHKERRQ(ierr); 45*c4762a1bSJed Brown ierr = VecSetSizes(ctx->d,PETSC_DECIDE,ctx->m);CHKERRQ(ierr); 46*c4762a1bSJed Brown ierr = VecSetFromOptions(ctx->d);CHKERRQ(ierr); 47*c4762a1bSJed Brown ierr = VecSetRandom(ctx->d,ctx->rctx);CHKERRQ(ierr); 48*c4762a1bSJed Brown PetscFunctionReturn(0); 49*c4762a1bSJed Brown } 50*c4762a1bSJed Brown 51*c4762a1bSJed Brown static PetscErrorCode CreateMatrix(UserCtx ctx) 52*c4762a1bSJed Brown { 53*c4762a1bSJed Brown PetscInt Istart,Iend,i,j,Ii,gridN,I_n, I_s, I_e, I_w; 54*c4762a1bSJed Brown #if defined(PETSC_USE_LOG) 55*c4762a1bSJed Brown PetscLogStage stage; 56*c4762a1bSJed Brown #endif 57*c4762a1bSJed Brown PetscErrorCode ierr; 58*c4762a1bSJed Brown 59*c4762a1bSJed Brown PetscFunctionBegin; 60*c4762a1bSJed Brown /* build the matrix F in ctx */ 61*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD, &(ctx->F));CHKERRQ(ierr); 62*c4762a1bSJed Brown ierr = MatSetSizes(ctx->F,PETSC_DECIDE, PETSC_DECIDE, ctx->m, ctx->n);CHKERRQ(ierr); 63*c4762a1bSJed Brown ierr = MatSetType(ctx->F,MATAIJ);CHKERRQ(ierr); /* TODO: Decide specific SetType other than dummy*/ 64*c4762a1bSJed Brown ierr = MatMPIAIJSetPreallocation(ctx->F, 5, NULL, 5, NULL);CHKERRQ(ierr); /*TODO: some number other than 5?*/ 65*c4762a1bSJed Brown ierr = MatSeqAIJSetPreallocation(ctx->F, 5, NULL);CHKERRQ(ierr); 66*c4762a1bSJed Brown ierr = MatSetUp(ctx->F);CHKERRQ(ierr); 67*c4762a1bSJed Brown ierr = MatGetOwnershipRange(ctx->F,&Istart,&Iend);CHKERRQ(ierr); 68*c4762a1bSJed Brown ierr = PetscLogStageRegister("Assembly", &stage);CHKERRQ(ierr); 69*c4762a1bSJed Brown ierr = PetscLogStagePush(stage);CHKERRQ(ierr); 70*c4762a1bSJed Brown 71*c4762a1bSJed Brown /* Set matrix elements in 2-D fiveopoint stencil format. */ 72*c4762a1bSJed Brown if (!(ctx->matops)) { 73*c4762a1bSJed Brown if (ctx->m != ctx->n) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_SIZ, "Stencil matrix must be square"); 74*c4762a1bSJed Brown gridN = (PetscInt) PetscSqrtReal((PetscReal) ctx->m); 75*c4762a1bSJed Brown if (gridN * gridN != ctx->m) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_SIZ, "Number of rows must be square"); 76*c4762a1bSJed Brown for (Ii=Istart; Ii<Iend; Ii++) { 77*c4762a1bSJed Brown i = Ii / gridN; j = Ii % gridN; 78*c4762a1bSJed Brown I_n = i * gridN + j + 1; 79*c4762a1bSJed Brown if (j + 1 >= gridN) I_n = -1; 80*c4762a1bSJed Brown I_s = i * gridN + j - 1; 81*c4762a1bSJed Brown if (j - 1 < 0) I_s = -1; 82*c4762a1bSJed Brown I_e = (i + 1) * gridN + j; 83*c4762a1bSJed Brown if (i + 1 >= gridN) I_e = -1; 84*c4762a1bSJed Brown I_w = (i - 1) * gridN + j; 85*c4762a1bSJed Brown if (i - 1 < 0) I_w = -1; 86*c4762a1bSJed Brown ierr = MatSetValue(ctx->F, Ii, Ii, 4., INSERT_VALUES);CHKERRQ(ierr); 87*c4762a1bSJed Brown ierr = MatSetValue(ctx->F, Ii, I_n, -1., INSERT_VALUES);CHKERRQ(ierr); 88*c4762a1bSJed Brown ierr = MatSetValue(ctx->F, Ii, I_s, -1., INSERT_VALUES);CHKERRQ(ierr); 89*c4762a1bSJed Brown ierr = MatSetValue(ctx->F, Ii, I_e, -1., INSERT_VALUES);CHKERRQ(ierr); 90*c4762a1bSJed Brown ierr = MatSetValue(ctx->F, Ii, I_w, -1., INSERT_VALUES);CHKERRQ(ierr); 91*c4762a1bSJed Brown } 92*c4762a1bSJed Brown } else {ierr = MatSetRandom(ctx->F, ctx->rctx);CHKERRQ(ierr);} 93*c4762a1bSJed Brown ierr = MatAssemblyBegin(ctx->F, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 94*c4762a1bSJed Brown ierr = MatAssemblyEnd(ctx->F, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 95*c4762a1bSJed Brown ierr = PetscLogStagePop();CHKERRQ(ierr); 96*c4762a1bSJed Brown /* Stencil matrix is symmetric. Setting symmetric flag for ICC/Cholesky preconditioner */ 97*c4762a1bSJed Brown if (!(ctx->matops)) { 98*c4762a1bSJed Brown ierr = MatSetOption(ctx->F,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr); 99*c4762a1bSJed Brown } 100*c4762a1bSJed Brown ierr = MatTransposeMatMult(ctx->F,ctx->F, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &(ctx->W));CHKERRQ(ierr); 101*c4762a1bSJed Brown /* Setup Hessian Workspace in same shape as W */ 102*c4762a1bSJed Brown ierr = MatDuplicate(ctx->W,MAT_DO_NOT_COPY_VALUES,&(ctx->Hm));CHKERRQ(ierr); 103*c4762a1bSJed Brown ierr = MatDuplicate(ctx->W,MAT_DO_NOT_COPY_VALUES,&(ctx->Hr));CHKERRQ(ierr); 104*c4762a1bSJed Brown PetscFunctionReturn(0); 105*c4762a1bSJed Brown } 106*c4762a1bSJed Brown 107*c4762a1bSJed Brown static PetscErrorCode SetupWorkspace(UserCtx ctx) 108*c4762a1bSJed Brown { 109*c4762a1bSJed Brown PetscInt i; 110*c4762a1bSJed Brown PetscErrorCode ierr; 111*c4762a1bSJed Brown 112*c4762a1bSJed Brown PetscFunctionBegin; 113*c4762a1bSJed Brown ierr = MatCreateVecs(ctx->F, &ctx->workLeft[0], &ctx->workRight[0]);CHKERRQ(ierr); 114*c4762a1bSJed Brown for (i=1; i<NWORKLEFT; i++) { 115*c4762a1bSJed Brown ierr = VecDuplicate(ctx->workLeft[0], &(ctx->workLeft[i]));CHKERRQ(ierr); 116*c4762a1bSJed Brown } 117*c4762a1bSJed Brown for (i=1; i<NWORKRIGHT; i++) { 118*c4762a1bSJed Brown ierr = VecDuplicate(ctx->workRight[0], &(ctx->workRight[i]));CHKERRQ(ierr); 119*c4762a1bSJed Brown } 120*c4762a1bSJed Brown PetscFunctionReturn(0); 121*c4762a1bSJed Brown } 122*c4762a1bSJed Brown 123*c4762a1bSJed Brown static PetscErrorCode ConfigureContext(UserCtx ctx) 124*c4762a1bSJed Brown { 125*c4762a1bSJed Brown PetscErrorCode ierr; 126*c4762a1bSJed Brown 127*c4762a1bSJed Brown PetscFunctionBegin; 128*c4762a1bSJed Brown ctx->m = 16; 129*c4762a1bSJed Brown ctx->n = 16; 130*c4762a1bSJed Brown ctx->eps = 1.e-3; 131*c4762a1bSJed Brown ctx->abstol = 1.e-4; 132*c4762a1bSJed Brown ctx->reltol = 1.e-2; 133*c4762a1bSJed Brown ctx->hStart = 1.; 134*c4762a1bSJed Brown ctx->hMin = 1.e-3; 135*c4762a1bSJed Brown ctx->hFactor = 0.5; 136*c4762a1bSJed Brown ctx->alpha = 1.; 137*c4762a1bSJed Brown ctx->mu = 1.0; 138*c4762a1bSJed Brown ctx->matops = 0; 139*c4762a1bSJed Brown ctx->iter = 10; 140*c4762a1bSJed Brown ctx->p = NORM_2; 141*c4762a1bSJed Brown ctx->taylor = PETSC_TRUE; 142*c4762a1bSJed Brown ctx->use_admm = PETSC_FALSE; 143*c4762a1bSJed Brown ierr = PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Configure separable objection example", "ex4.c");CHKERRQ(ierr); 144*c4762a1bSJed Brown ierr = PetscOptionsInt("-m", "The row dimension of matrix F", "ex4.c", ctx->m, &(ctx->m), NULL);CHKERRQ(ierr); 145*c4762a1bSJed Brown ierr = PetscOptionsInt("-n", "The column dimension of matrix F", "ex4.c", ctx->n, &(ctx->n), NULL);CHKERRQ(ierr); 146*c4762a1bSJed Brown ierr = PetscOptionsInt("-matrix_format","Decide format of F matrix. 0 for stencil, 1 for random", "ex4.c", ctx->matops, &(ctx->matops), NULL);CHKERRQ(ierr); 147*c4762a1bSJed Brown ierr = PetscOptionsInt("-iter","Iteration number ADMM", "ex4.c", ctx->iter, &(ctx->iter), NULL);CHKERRQ(ierr); 148*c4762a1bSJed Brown ierr = PetscOptionsReal("-alpha", "The regularization multiplier. 1 default", "ex4.c", ctx->alpha, &(ctx->alpha), NULL);CHKERRQ(ierr); 149*c4762a1bSJed Brown ierr = PetscOptionsReal("-epsilon", "The small constant added to |x_i| in the denominator to approximate the gradient of ||x||_1", "ex4.c", ctx->eps, &(ctx->eps), NULL);CHKERRQ(ierr); 150*c4762a1bSJed Brown ierr = PetscOptionsReal("-mu", "The augmented lagrangian multiplier in ADMM", "ex4.c", ctx->mu, &(ctx->mu), NULL);CHKERRQ(ierr); 151*c4762a1bSJed Brown ierr = PetscOptionsReal("-hStart", "Taylor test starting point. 1 default.", "ex4.c", ctx->hStart, &(ctx->hStart), NULL);CHKERRQ(ierr); 152*c4762a1bSJed Brown ierr = PetscOptionsReal("-hFactor", "Taylor test multiplier factor. 0.5 default", "ex4.c", ctx->hFactor, &(ctx->hFactor), NULL);CHKERRQ(ierr); 153*c4762a1bSJed Brown ierr = PetscOptionsReal("-hMin", "Taylor test ending condition. 1.e-3 default", "ex4.c", ctx->hMin, &(ctx->hMin), NULL);CHKERRQ(ierr); 154*c4762a1bSJed Brown ierr = PetscOptionsReal("-abstol", "Absolute stopping criterion for ADMM", "ex4.c", ctx->abstol, &(ctx->abstol), NULL);CHKERRQ(ierr); 155*c4762a1bSJed Brown ierr = PetscOptionsReal("-reltol", "Relative stopping criterion for ADMM", "ex4.c", ctx->reltol, &(ctx->reltol), NULL);CHKERRQ(ierr); 156*c4762a1bSJed Brown ierr = PetscOptionsBool("-taylor","Flag for Taylor test. Default is true.", "ex4.c", ctx->taylor, &(ctx->taylor), NULL);CHKERRQ(ierr); 157*c4762a1bSJed Brown ierr = PetscOptionsBool("-use_admm","Use the ADMM solver in this example.", "ex4.c", ctx->use_admm, &(ctx->use_admm), NULL);CHKERRQ(ierr); 158*c4762a1bSJed Brown ierr = PetscOptionsEnum("-p","Norm type.", "ex4.c", NormTypes, (PetscEnum)ctx->p, (PetscEnum *) &(ctx->p), NULL);CHKERRQ(ierr); 159*c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 160*c4762a1bSJed Brown /* Creating random ctx */ 161*c4762a1bSJed Brown ierr = PetscRandomCreate(PETSC_COMM_WORLD,&(ctx->rctx));CHKERRQ(ierr); 162*c4762a1bSJed Brown ierr = PetscRandomSetFromOptions(ctx->rctx);CHKERRQ(ierr); 163*c4762a1bSJed Brown ierr = CreateMatrix(ctx);CHKERRQ(ierr); 164*c4762a1bSJed Brown ierr = CreateRHS(ctx);CHKERRQ(ierr); 165*c4762a1bSJed Brown ierr = SetupWorkspace(ctx);CHKERRQ(ierr); 166*c4762a1bSJed Brown PetscFunctionReturn(0); 167*c4762a1bSJed Brown } 168*c4762a1bSJed Brown 169*c4762a1bSJed Brown static PetscErrorCode DestroyContext(UserCtx *ctx) 170*c4762a1bSJed Brown { 171*c4762a1bSJed Brown PetscInt i; 172*c4762a1bSJed Brown PetscErrorCode ierr; 173*c4762a1bSJed Brown 174*c4762a1bSJed Brown PetscFunctionBegin; 175*c4762a1bSJed Brown ierr = MatDestroy(&((*ctx)->F));CHKERRQ(ierr); 176*c4762a1bSJed Brown ierr = MatDestroy(&((*ctx)->W));CHKERRQ(ierr); 177*c4762a1bSJed Brown ierr = MatDestroy(&((*ctx)->Hm));CHKERRQ(ierr); 178*c4762a1bSJed Brown ierr = MatDestroy(&((*ctx)->Hr));CHKERRQ(ierr); 179*c4762a1bSJed Brown ierr = VecDestroy(&((*ctx)->d));CHKERRQ(ierr); 180*c4762a1bSJed Brown for (i=0; i<NWORKLEFT; i++) { 181*c4762a1bSJed Brown ierr = VecDestroy(&((*ctx)->workLeft[i]));CHKERRQ(ierr); 182*c4762a1bSJed Brown } 183*c4762a1bSJed Brown for (i=0; i<NWORKRIGHT; i++) { 184*c4762a1bSJed Brown ierr = VecDestroy(&((*ctx)->workRight[i]));CHKERRQ(ierr); 185*c4762a1bSJed Brown } 186*c4762a1bSJed Brown ierr = PetscRandomDestroy(&((*ctx)->rctx));CHKERRQ(ierr); 187*c4762a1bSJed Brown ierr = PetscFree(*ctx);CHKERRQ(ierr); 188*c4762a1bSJed Brown PetscFunctionReturn(0); 189*c4762a1bSJed Brown } 190*c4762a1bSJed Brown 191*c4762a1bSJed Brown /* compute (1/2) * ||F x - d||^2 */ 192*c4762a1bSJed Brown static PetscErrorCode ObjectiveMisfit(Tao tao, Vec x, PetscReal *J, void *_ctx) 193*c4762a1bSJed Brown { 194*c4762a1bSJed Brown UserCtx ctx = (UserCtx) _ctx; 195*c4762a1bSJed Brown PetscErrorCode ierr; 196*c4762a1bSJed Brown Vec y; 197*c4762a1bSJed Brown 198*c4762a1bSJed Brown PetscFunctionBegin; 199*c4762a1bSJed Brown y = ctx->workLeft[0]; 200*c4762a1bSJed Brown ierr = MatMult(ctx->F, x, y);CHKERRQ(ierr); 201*c4762a1bSJed Brown ierr = VecAXPY(y, -1., ctx->d);CHKERRQ(ierr); 202*c4762a1bSJed Brown ierr = VecDot(y, y, J);CHKERRQ(ierr); 203*c4762a1bSJed Brown *J *= 0.5; 204*c4762a1bSJed Brown PetscFunctionReturn(0); 205*c4762a1bSJed Brown } 206*c4762a1bSJed Brown 207*c4762a1bSJed Brown /* compute V = FTFx - FTd */ 208*c4762a1bSJed Brown static PetscErrorCode GradientMisfit(Tao tao, Vec x, Vec V, void *_ctx) 209*c4762a1bSJed Brown { 210*c4762a1bSJed Brown UserCtx ctx = (UserCtx) _ctx; 211*c4762a1bSJed Brown PetscErrorCode ierr; 212*c4762a1bSJed Brown Vec FTFx, FTd; 213*c4762a1bSJed Brown 214*c4762a1bSJed Brown PetscFunctionBegin; 215*c4762a1bSJed Brown /* work1 is A^T Ax, work2 is Ab, W is A^T A*/ 216*c4762a1bSJed Brown FTFx = ctx->workRight[0]; 217*c4762a1bSJed Brown FTd = ctx->workRight[1]; 218*c4762a1bSJed Brown ierr = MatMult(ctx->W,x,FTFx);CHKERRQ(ierr); 219*c4762a1bSJed Brown ierr = MatMultTranspose(ctx->F, ctx->d, FTd);CHKERRQ(ierr); 220*c4762a1bSJed Brown ierr = VecWAXPY(V, -1., FTd, FTFx);CHKERRQ(ierr); 221*c4762a1bSJed Brown PetscFunctionReturn(0); 222*c4762a1bSJed Brown } 223*c4762a1bSJed Brown 224*c4762a1bSJed Brown /* returns FTF */ 225*c4762a1bSJed Brown static PetscErrorCode HessianMisfit(Tao tao, Vec x, Mat H, Mat Hpre, void *_ctx) 226*c4762a1bSJed Brown { 227*c4762a1bSJed Brown UserCtx ctx = (UserCtx) _ctx; 228*c4762a1bSJed Brown PetscErrorCode ierr; 229*c4762a1bSJed Brown 230*c4762a1bSJed Brown PetscFunctionBegin; 231*c4762a1bSJed Brown if (H != ctx->W) {ierr = MatCopy(ctx->W, H, DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);} 232*c4762a1bSJed Brown if (Hpre != ctx->W) {ierr = MatCopy(ctx->W, Hpre, DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);} 233*c4762a1bSJed Brown PetscFunctionReturn(0); 234*c4762a1bSJed Brown } 235*c4762a1bSJed Brown 236*c4762a1bSJed Brown /* computes augment Lagrangian objective (with scaled dual): 237*c4762a1bSJed Brown * 0.5 * ||F x - d||^2 + 0.5 * mu ||x - z + u||^2 */ 238*c4762a1bSJed Brown static PetscErrorCode ObjectiveMisfitADMM(Tao tao, Vec x, PetscReal *J, void *_ctx) 239*c4762a1bSJed Brown { 240*c4762a1bSJed Brown UserCtx ctx = (UserCtx) _ctx; 241*c4762a1bSJed Brown PetscReal mu, workNorm, misfit; 242*c4762a1bSJed Brown Vec z, u, temp; 243*c4762a1bSJed Brown PetscErrorCode ierr; 244*c4762a1bSJed Brown 245*c4762a1bSJed Brown PetscFunctionBegin; 246*c4762a1bSJed Brown mu = ctx->mu; 247*c4762a1bSJed Brown z = ctx->workRight[5]; 248*c4762a1bSJed Brown u = ctx->workRight[6]; 249*c4762a1bSJed Brown temp = ctx->workRight[10]; 250*c4762a1bSJed Brown /* misfit = f(x) */ 251*c4762a1bSJed Brown ierr = ObjectiveMisfit(tao, x, &misfit, _ctx);CHKERRQ(ierr); 252*c4762a1bSJed Brown ierr = VecCopy(x,temp);CHKERRQ(ierr); 253*c4762a1bSJed Brown /* temp = x - z + u */ 254*c4762a1bSJed Brown ierr = VecAXPBYPCZ(temp,-1.,1.,1.,z,u);CHKERRQ(ierr); 255*c4762a1bSJed Brown /* workNorm = ||x - z + u||^2 */ 256*c4762a1bSJed Brown ierr = VecDot(temp, temp, &workNorm);CHKERRQ(ierr); 257*c4762a1bSJed Brown /* augment Lagrangian objective (with scaled dual): f(x) + 0.5 * mu ||x -z + u||^2 */ 258*c4762a1bSJed Brown *J = misfit + 0.5 * mu * workNorm; 259*c4762a1bSJed Brown PetscFunctionReturn(0); 260*c4762a1bSJed Brown } 261*c4762a1bSJed Brown 262*c4762a1bSJed Brown /* computes FTFx - FTd mu*(x - z + u) */ 263*c4762a1bSJed Brown static PetscErrorCode GradientMisfitADMM(Tao tao, Vec x, Vec V, void *_ctx) 264*c4762a1bSJed Brown { 265*c4762a1bSJed Brown UserCtx ctx = (UserCtx) _ctx; 266*c4762a1bSJed Brown PetscReal mu; 267*c4762a1bSJed Brown Vec z, u, temp; 268*c4762a1bSJed Brown PetscErrorCode ierr; 269*c4762a1bSJed Brown 270*c4762a1bSJed Brown PetscFunctionBegin; 271*c4762a1bSJed Brown mu = ctx->mu; 272*c4762a1bSJed Brown z = ctx->workRight[5]; 273*c4762a1bSJed Brown u = ctx->workRight[6]; 274*c4762a1bSJed Brown temp = ctx->workRight[10]; 275*c4762a1bSJed Brown ierr = GradientMisfit(tao, x, V, _ctx);CHKERRQ(ierr); 276*c4762a1bSJed Brown ierr = VecCopy(x, temp);CHKERRQ(ierr); 277*c4762a1bSJed Brown /* temp = x - z + u */ 278*c4762a1bSJed Brown ierr = VecAXPBYPCZ(temp,-1.,1.,1.,z,u);CHKERRQ(ierr); 279*c4762a1bSJed Brown /* V = FTFx - FTd mu*(x - z + u) */ 280*c4762a1bSJed Brown ierr = VecAXPY(V, mu, temp);CHKERRQ(ierr); 281*c4762a1bSJed Brown PetscFunctionReturn(0); 282*c4762a1bSJed Brown } 283*c4762a1bSJed Brown 284*c4762a1bSJed Brown /* returns FTF + diag(mu) */ 285*c4762a1bSJed Brown static PetscErrorCode HessianMisfitADMM(Tao tao, Vec x, Mat H, Mat Hpre, void *_ctx) 286*c4762a1bSJed Brown { 287*c4762a1bSJed Brown UserCtx ctx = (UserCtx) _ctx; 288*c4762a1bSJed Brown PetscErrorCode ierr; 289*c4762a1bSJed Brown 290*c4762a1bSJed Brown PetscFunctionBegin; 291*c4762a1bSJed Brown ierr = MatCopy(ctx->W, H, DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 292*c4762a1bSJed Brown ierr = MatShift(H, ctx->mu);CHKERRQ(ierr); 293*c4762a1bSJed Brown if (Hpre != H) { 294*c4762a1bSJed Brown ierr = MatCopy(H, Hpre, DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 295*c4762a1bSJed Brown } 296*c4762a1bSJed Brown PetscFunctionReturn(0); 297*c4762a1bSJed Brown } 298*c4762a1bSJed Brown 299*c4762a1bSJed Brown /* computes || x ||_p (mult by 0.5 in case of NORM_2) */ 300*c4762a1bSJed Brown static PetscErrorCode ObjectiveRegularization(Tao tao, Vec x, PetscReal *J, void *_ctx) 301*c4762a1bSJed Brown { 302*c4762a1bSJed Brown UserCtx ctx = (UserCtx) _ctx; 303*c4762a1bSJed Brown PetscReal norm; 304*c4762a1bSJed Brown PetscErrorCode ierr; 305*c4762a1bSJed Brown 306*c4762a1bSJed Brown PetscFunctionBegin; 307*c4762a1bSJed Brown *J = 0; 308*c4762a1bSJed Brown ierr = VecNorm (x, ctx->p, &norm);CHKERRQ(ierr); 309*c4762a1bSJed Brown if (ctx->p == NORM_2) norm = 0.5 * norm * norm; 310*c4762a1bSJed Brown *J = ctx->alpha * norm; 311*c4762a1bSJed Brown PetscFunctionReturn(0); 312*c4762a1bSJed Brown } 313*c4762a1bSJed Brown 314*c4762a1bSJed Brown /* NORM_2 Case: return x 315*c4762a1bSJed Brown * NORM_1 Case: x/(|x| + eps) 316*c4762a1bSJed Brown * Else: TODO */ 317*c4762a1bSJed Brown static PetscErrorCode GradientRegularization(Tao tao, Vec x, Vec V, void *_ctx) 318*c4762a1bSJed Brown { 319*c4762a1bSJed Brown UserCtx ctx = (UserCtx) _ctx; 320*c4762a1bSJed Brown PetscErrorCode ierr; 321*c4762a1bSJed Brown PetscReal eps = ctx->eps; 322*c4762a1bSJed Brown 323*c4762a1bSJed Brown PetscFunctionBegin; 324*c4762a1bSJed Brown if (ctx->p == NORM_2) { 325*c4762a1bSJed Brown ierr = VecCopy(x, V);CHKERRQ(ierr); 326*c4762a1bSJed Brown } else if (ctx->p == NORM_1) { 327*c4762a1bSJed Brown ierr = VecCopy(x, ctx->workRight[1]);CHKERRQ(ierr); 328*c4762a1bSJed Brown ierr = VecAbs(ctx->workRight[1]);CHKERRQ(ierr); 329*c4762a1bSJed Brown ierr = VecShift(ctx->workRight[1], eps);CHKERRQ(ierr); 330*c4762a1bSJed Brown ierr = VecPointwiseDivide(V, x, ctx->workRight[1]);CHKERRQ(ierr); 331*c4762a1bSJed Brown } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Example only works for NORM_1 and NORM_2"); 332*c4762a1bSJed Brown PetscFunctionReturn(0); 333*c4762a1bSJed Brown } 334*c4762a1bSJed Brown 335*c4762a1bSJed Brown /* NORM_2 Case: returns diag(mu) 336*c4762a1bSJed Brown * NORM_1 Case: diag(mu* 1/sqrt(x_i^2 + eps) * ( 1 - x_i^2/ABS(x_i^2+eps))) */ 337*c4762a1bSJed Brown static PetscErrorCode HessianRegularization(Tao tao, Vec x, Mat H, Mat Hpre, void *_ctx) 338*c4762a1bSJed Brown { 339*c4762a1bSJed Brown UserCtx ctx = (UserCtx) _ctx; 340*c4762a1bSJed Brown PetscReal eps = ctx->eps; 341*c4762a1bSJed Brown Vec copy1,copy2,copy3; 342*c4762a1bSJed Brown PetscErrorCode ierr; 343*c4762a1bSJed Brown 344*c4762a1bSJed Brown PetscFunctionBegin; 345*c4762a1bSJed Brown if (ctx->p == NORM_2) { 346*c4762a1bSJed Brown /* Identity matrix scaled by mu */ 347*c4762a1bSJed Brown ierr = MatZeroEntries(H);CHKERRQ(ierr); 348*c4762a1bSJed Brown ierr = MatShift(H,ctx->mu);CHKERRQ(ierr); 349*c4762a1bSJed Brown if (Hpre != H) { 350*c4762a1bSJed Brown ierr = MatZeroEntries(Hpre);CHKERRQ(ierr); 351*c4762a1bSJed Brown ierr = MatShift(Hpre,ctx->mu);CHKERRQ(ierr); 352*c4762a1bSJed Brown } 353*c4762a1bSJed Brown } else if (ctx->p == NORM_1) { 354*c4762a1bSJed Brown /* 1/sqrt(x_i^2 + eps) * ( 1 - x_i^2/ABS(x_i^2+eps) ) */ 355*c4762a1bSJed Brown copy1 = ctx->workRight[1]; 356*c4762a1bSJed Brown copy2 = ctx->workRight[2]; 357*c4762a1bSJed Brown copy3 = ctx->workRight[3]; 358*c4762a1bSJed Brown /* copy1 : 1/sqrt(x_i^2 + eps) */ 359*c4762a1bSJed Brown ierr = VecCopy(x, copy1);CHKERRQ(ierr); 360*c4762a1bSJed Brown ierr = VecPow(copy1,2);CHKERRQ(ierr); 361*c4762a1bSJed Brown ierr = VecShift(copy1, eps);CHKERRQ(ierr); 362*c4762a1bSJed Brown ierr = VecSqrtAbs(copy1);CHKERRQ(ierr); 363*c4762a1bSJed Brown ierr = VecReciprocal(copy1);CHKERRQ(ierr); 364*c4762a1bSJed Brown /* copy2: x_i^2.*/ 365*c4762a1bSJed Brown ierr = VecCopy(x,copy2);CHKERRQ(ierr); 366*c4762a1bSJed Brown ierr = VecPow(copy2,2);CHKERRQ(ierr); 367*c4762a1bSJed Brown /* copy3: abs(x_i^2 + eps) */ 368*c4762a1bSJed Brown ierr = VecCopy(x,copy3);CHKERRQ(ierr); 369*c4762a1bSJed Brown ierr = VecPow(copy3,2);CHKERRQ(ierr); 370*c4762a1bSJed Brown ierr = VecShift(copy3, eps);CHKERRQ(ierr); 371*c4762a1bSJed Brown ierr = VecAbs(copy3);CHKERRQ(ierr); 372*c4762a1bSJed Brown /* copy2: 1 - x_i^2/abs(x_i^2 + eps) */ 373*c4762a1bSJed Brown ierr = VecPointwiseDivide(copy2, copy2,copy3);CHKERRQ(ierr); 374*c4762a1bSJed Brown ierr = VecScale(copy2, -1.);CHKERRQ(ierr); 375*c4762a1bSJed Brown ierr = VecShift(copy2, 1.);CHKERRQ(ierr); 376*c4762a1bSJed Brown ierr = VecAXPY(copy1,1.,copy2);CHKERRQ(ierr); 377*c4762a1bSJed Brown ierr = VecScale(copy1, ctx->mu);CHKERRQ(ierr); 378*c4762a1bSJed Brown ierr = MatZeroEntries(H);CHKERRQ(ierr); 379*c4762a1bSJed Brown ierr = MatDiagonalSet(H, copy1,INSERT_VALUES);CHKERRQ(ierr); 380*c4762a1bSJed Brown if (Hpre != H) { 381*c4762a1bSJed Brown ierr = MatZeroEntries(Hpre);CHKERRQ(ierr); 382*c4762a1bSJed Brown ierr = MatDiagonalSet(Hpre, copy1,INSERT_VALUES);CHKERRQ(ierr); 383*c4762a1bSJed Brown } 384*c4762a1bSJed Brown } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Example only works for NORM_1 and NORM_2"); 385*c4762a1bSJed Brown PetscFunctionReturn(0); 386*c4762a1bSJed Brown } 387*c4762a1bSJed Brown 388*c4762a1bSJed Brown /* NORM_2 Case: 0.5 || x ||_2 + 0.5 * mu * ||x + u - z||^2 389*c4762a1bSJed Brown * Else : || x ||_2 + 0.5 * mu * ||x + u - z||^2 */ 390*c4762a1bSJed Brown static PetscErrorCode ObjectiveRegularizationADMM(Tao tao, Vec z, PetscReal *J, void *_ctx) 391*c4762a1bSJed Brown { 392*c4762a1bSJed Brown UserCtx ctx = (UserCtx) _ctx; 393*c4762a1bSJed Brown PetscReal mu, workNorm, reg; 394*c4762a1bSJed Brown Vec x, u, temp; 395*c4762a1bSJed Brown PetscErrorCode ierr; 396*c4762a1bSJed Brown 397*c4762a1bSJed Brown PetscFunctionBegin; 398*c4762a1bSJed Brown mu = ctx->mu; 399*c4762a1bSJed Brown x = ctx->workRight[4]; 400*c4762a1bSJed Brown u = ctx->workRight[6]; 401*c4762a1bSJed Brown temp = ctx->workRight[10]; 402*c4762a1bSJed Brown ierr = ObjectiveRegularization(tao, z, ®, _ctx);CHKERRQ(ierr); 403*c4762a1bSJed Brown ierr = VecCopy(z,temp);CHKERRQ(ierr); 404*c4762a1bSJed Brown /* temp = x + u -z */ 405*c4762a1bSJed Brown ierr = VecAXPBYPCZ(temp,1.,1.,-1.,x,u);CHKERRQ(ierr); 406*c4762a1bSJed Brown /* workNorm = ||x + u - z ||^2 */ 407*c4762a1bSJed Brown ierr = VecDot(temp, temp, &workNorm);CHKERRQ(ierr); 408*c4762a1bSJed Brown *J = reg + 0.5 * mu * workNorm; 409*c4762a1bSJed Brown PetscFunctionReturn(0); 410*c4762a1bSJed Brown } 411*c4762a1bSJed Brown 412*c4762a1bSJed Brown 413*c4762a1bSJed Brown /* NORM_2 Case: x - mu*(x + u - z) 414*c4762a1bSJed Brown * NORM_1 Case: x/(|x| + eps) - mu*(x + u - z) 415*c4762a1bSJed Brown * Else: TODO */ 416*c4762a1bSJed Brown static PetscErrorCode GradientRegularizationADMM(Tao tao, Vec z, Vec V, void *_ctx) 417*c4762a1bSJed Brown { 418*c4762a1bSJed Brown UserCtx ctx = (UserCtx) _ctx; 419*c4762a1bSJed Brown PetscReal mu; 420*c4762a1bSJed Brown Vec x, u, temp; 421*c4762a1bSJed Brown PetscErrorCode ierr; 422*c4762a1bSJed Brown 423*c4762a1bSJed Brown PetscFunctionBegin; 424*c4762a1bSJed Brown mu = ctx->mu; 425*c4762a1bSJed Brown x = ctx->workRight[4]; 426*c4762a1bSJed Brown u = ctx->workRight[6]; 427*c4762a1bSJed Brown temp = ctx->workRight[10]; 428*c4762a1bSJed Brown ierr = GradientRegularization(tao, z, V, _ctx);CHKERRQ(ierr); 429*c4762a1bSJed Brown ierr = VecCopy(z, temp);CHKERRQ(ierr); 430*c4762a1bSJed Brown /* temp = x + u -z */ 431*c4762a1bSJed Brown ierr = VecAXPBYPCZ(temp,1.,1.,-1.,x,u);CHKERRQ(ierr); 432*c4762a1bSJed Brown ierr = VecAXPY(V, -mu, temp);CHKERRQ(ierr); 433*c4762a1bSJed Brown PetscFunctionReturn(0); 434*c4762a1bSJed Brown } 435*c4762a1bSJed Brown 436*c4762a1bSJed Brown /* NORM_2 Case: returns diag(mu) 437*c4762a1bSJed Brown * NORM_1 Case: FTF + diag(mu) */ 438*c4762a1bSJed Brown static PetscErrorCode HessianRegularizationADMM(Tao tao, Vec x, Mat H, Mat Hpre, void *_ctx) 439*c4762a1bSJed Brown { 440*c4762a1bSJed Brown UserCtx ctx = (UserCtx) _ctx; 441*c4762a1bSJed Brown PetscErrorCode ierr; 442*c4762a1bSJed Brown 443*c4762a1bSJed Brown PetscFunctionBegin; 444*c4762a1bSJed Brown if (ctx->p == NORM_2) { 445*c4762a1bSJed Brown /* Identity matrix scaled by mu */ 446*c4762a1bSJed Brown ierr = MatZeroEntries(H);CHKERRQ(ierr); 447*c4762a1bSJed Brown ierr = MatShift(H,ctx->mu);CHKERRQ(ierr); 448*c4762a1bSJed Brown if (Hpre != H) { 449*c4762a1bSJed Brown ierr = MatZeroEntries(Hpre);CHKERRQ(ierr); 450*c4762a1bSJed Brown ierr = MatShift(Hpre,ctx->mu);CHKERRQ(ierr); 451*c4762a1bSJed Brown } 452*c4762a1bSJed Brown } else if (ctx->p == NORM_1) { 453*c4762a1bSJed Brown ierr = HessianMisfit(tao, x, H, Hpre, (void*) ctx);CHKERRQ(ierr); 454*c4762a1bSJed Brown ierr = MatShift(H, ctx->mu);CHKERRQ(ierr); 455*c4762a1bSJed Brown if (Hpre != H) {ierr = MatShift(Hpre, ctx->mu);CHKERRQ(ierr);} 456*c4762a1bSJed Brown } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Example only works for NORM_1 and NORM_2"); 457*c4762a1bSJed Brown PetscFunctionReturn(0); 458*c4762a1bSJed Brown } 459*c4762a1bSJed Brown 460*c4762a1bSJed Brown /* NORM_2 Case : (1/2) * ||F x - d||^2 + 0.5 * || x ||_p 461*c4762a1bSJed Brown * NORM_1 Case : (1/2) * ||F x - d||^2 + || x ||_p */ 462*c4762a1bSJed Brown static PetscErrorCode ObjectiveComplete(Tao tao, Vec x, PetscReal *J, void *ctx) 463*c4762a1bSJed Brown { 464*c4762a1bSJed Brown PetscReal Jm, Jr; 465*c4762a1bSJed Brown PetscErrorCode ierr; 466*c4762a1bSJed Brown 467*c4762a1bSJed Brown PetscFunctionBegin; 468*c4762a1bSJed Brown ierr = ObjectiveMisfit(tao, x, &Jm, ctx);CHKERRQ(ierr); 469*c4762a1bSJed Brown ierr = ObjectiveRegularization(tao, x, &Jr, ctx);CHKERRQ(ierr); 470*c4762a1bSJed Brown *J = Jm + Jr; 471*c4762a1bSJed Brown PetscFunctionReturn(0); 472*c4762a1bSJed Brown } 473*c4762a1bSJed Brown 474*c4762a1bSJed Brown /* NORM_2 Case: FTFx - FTd + x 475*c4762a1bSJed Brown * NORM_1 Case: FTFx - FTd + x/(|x| + eps) */ 476*c4762a1bSJed Brown static PetscErrorCode GradientComplete(Tao tao, Vec x, Vec V, void *ctx) 477*c4762a1bSJed Brown { 478*c4762a1bSJed Brown UserCtx cntx = (UserCtx) ctx; 479*c4762a1bSJed Brown PetscErrorCode ierr; 480*c4762a1bSJed Brown 481*c4762a1bSJed Brown PetscFunctionBegin; 482*c4762a1bSJed Brown ierr = GradientMisfit(tao, x, cntx->workRight[2], ctx);CHKERRQ(ierr); 483*c4762a1bSJed Brown ierr = GradientRegularization(tao, x, cntx->workRight[3], ctx);CHKERRQ(ierr); 484*c4762a1bSJed Brown ierr = VecWAXPY(V,1,cntx->workRight[2],cntx->workRight[3]);CHKERRQ(ierr); 485*c4762a1bSJed Brown PetscFunctionReturn(0); 486*c4762a1bSJed Brown } 487*c4762a1bSJed Brown 488*c4762a1bSJed Brown /* NORM_2 Case: diag(mu) + FTF 489*c4762a1bSJed Brown * NORM_1 Case: diag(mu* 1/sqrt(x_i^2 + eps) * ( 1 - x_i^2/ABS(x_i^2+eps))) + FTF */ 490*c4762a1bSJed Brown static PetscErrorCode HessianComplete(Tao tao, Vec x, Mat H, Mat Hpre, void *ctx) 491*c4762a1bSJed Brown { 492*c4762a1bSJed Brown Mat tempH; 493*c4762a1bSJed Brown PetscErrorCode ierr; 494*c4762a1bSJed Brown 495*c4762a1bSJed Brown PetscFunctionBegin; 496*c4762a1bSJed Brown ierr = MatDuplicate(H, MAT_SHARE_NONZERO_PATTERN, &tempH);CHKERRQ(ierr); 497*c4762a1bSJed Brown ierr = HessianMisfit(tao, x, H, H, ctx);CHKERRQ(ierr); 498*c4762a1bSJed Brown ierr = HessianRegularization(tao, x, tempH, tempH, ctx);CHKERRQ(ierr); 499*c4762a1bSJed Brown ierr = MatAXPY(H, 1., tempH, DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 500*c4762a1bSJed Brown if (Hpre != H) { 501*c4762a1bSJed Brown ierr = MatCopy(H, Hpre, DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 502*c4762a1bSJed Brown } 503*c4762a1bSJed Brown ierr = MatDestroy(&tempH);CHKERRQ(ierr); 504*c4762a1bSJed Brown PetscFunctionReturn(0); 505*c4762a1bSJed Brown } 506*c4762a1bSJed Brown 507*c4762a1bSJed Brown static PetscErrorCode TaoSolveADMM(UserCtx ctx, Vec x) 508*c4762a1bSJed Brown { 509*c4762a1bSJed Brown PetscErrorCode ierr; 510*c4762a1bSJed Brown PetscInt i; 511*c4762a1bSJed Brown PetscReal u_norm, r_norm, s_norm, primal, dual, x_norm, z_norm; 512*c4762a1bSJed Brown Tao tao1,tao2; 513*c4762a1bSJed Brown Vec xk,z,u,diff,zold,zdiff,temp; 514*c4762a1bSJed Brown PetscReal mu; 515*c4762a1bSJed Brown 516*c4762a1bSJed Brown PetscFunctionBegin; 517*c4762a1bSJed Brown xk = ctx->workRight[4]; 518*c4762a1bSJed Brown z = ctx->workRight[5]; 519*c4762a1bSJed Brown u = ctx->workRight[6]; 520*c4762a1bSJed Brown diff = ctx->workRight[7]; 521*c4762a1bSJed Brown zold = ctx->workRight[8]; 522*c4762a1bSJed Brown zdiff = ctx->workRight[9]; 523*c4762a1bSJed Brown temp = ctx->workRight[11]; 524*c4762a1bSJed Brown mu = ctx->mu; 525*c4762a1bSJed Brown ierr = VecSet(u, 0.);CHKERRQ(ierr); 526*c4762a1bSJed Brown ierr = TaoCreate(PETSC_COMM_WORLD, &tao1);CHKERRQ(ierr); 527*c4762a1bSJed Brown ierr = TaoSetType(tao1,TAONLS);CHKERRQ(ierr); 528*c4762a1bSJed Brown ierr = TaoSetObjectiveRoutine(tao1, ObjectiveMisfitADMM, (void*) ctx);CHKERRQ(ierr); 529*c4762a1bSJed Brown ierr = TaoSetGradientRoutine(tao1, GradientMisfitADMM, (void*) ctx);CHKERRQ(ierr); 530*c4762a1bSJed Brown ierr = TaoSetHessianRoutine(tao1, ctx->Hm, ctx->Hm, HessianMisfitADMM, (void*) ctx);CHKERRQ(ierr); 531*c4762a1bSJed Brown ierr = VecSet(xk, 0.);CHKERRQ(ierr); 532*c4762a1bSJed Brown ierr = TaoSetInitialVector(tao1, xk);CHKERRQ(ierr); 533*c4762a1bSJed Brown ierr = TaoSetOptionsPrefix(tao1, "misfit_");CHKERRQ(ierr); 534*c4762a1bSJed Brown ierr = TaoSetFromOptions(tao1);CHKERRQ(ierr); 535*c4762a1bSJed Brown ierr = TaoCreate(PETSC_COMM_WORLD, &tao2);CHKERRQ(ierr); 536*c4762a1bSJed Brown if (ctx->p == NORM_2) { 537*c4762a1bSJed Brown ierr = TaoSetType(tao2,TAONLS);CHKERRQ(ierr); 538*c4762a1bSJed Brown ierr = TaoSetObjectiveRoutine(tao2, ObjectiveRegularizationADMM, (void*) ctx);CHKERRQ(ierr); 539*c4762a1bSJed Brown ierr = TaoSetGradientRoutine(tao2, GradientRegularizationADMM, (void*) ctx);CHKERRQ(ierr); 540*c4762a1bSJed Brown ierr = TaoSetHessianRoutine(tao2, ctx->Hr, ctx->Hr, HessianRegularizationADMM, (void*) ctx);CHKERRQ(ierr); 541*c4762a1bSJed Brown } 542*c4762a1bSJed Brown ierr = VecSet(z, 0.);CHKERRQ(ierr); 543*c4762a1bSJed Brown ierr = TaoSetInitialVector(tao2, z);CHKERRQ(ierr); 544*c4762a1bSJed Brown ierr = TaoSetOptionsPrefix(tao2, "reg_");CHKERRQ(ierr); 545*c4762a1bSJed Brown ierr = TaoSetFromOptions(tao2);CHKERRQ(ierr); 546*c4762a1bSJed Brown 547*c4762a1bSJed Brown for (i=0; i<ctx->iter; i++) { 548*c4762a1bSJed Brown ierr = VecCopy(z,zold);CHKERRQ(ierr); 549*c4762a1bSJed Brown ierr = TaoSolve(tao1);CHKERRQ(ierr); /* Updates xk */ 550*c4762a1bSJed Brown if (ctx->p == NORM_1){ 551*c4762a1bSJed Brown ierr = VecWAXPY(temp,1.,xk,u);CHKERRQ(ierr); 552*c4762a1bSJed Brown ierr = TaoSoftThreshold(temp,-ctx->alpha/mu,ctx->alpha/mu,z);CHKERRQ(ierr); 553*c4762a1bSJed Brown } else { 554*c4762a1bSJed Brown ierr = TaoSolve(tao2);CHKERRQ(ierr); /* Update zk */ 555*c4762a1bSJed Brown } 556*c4762a1bSJed Brown /* u = u + xk -z */ 557*c4762a1bSJed Brown ierr = VecAXPBYPCZ(u,1.,-1.,1.,xk,z);CHKERRQ(ierr); 558*c4762a1bSJed Brown /* r_norm : norm(x-z) */ 559*c4762a1bSJed Brown ierr = VecWAXPY(diff,-1.,z,xk);CHKERRQ(ierr); 560*c4762a1bSJed Brown ierr = VecNorm(diff,NORM_2,&r_norm);CHKERRQ(ierr); 561*c4762a1bSJed Brown /* s_norm : norm(-mu(z-zold)) */ 562*c4762a1bSJed Brown ierr = VecWAXPY(zdiff, -1.,zold,z);CHKERRQ(ierr); 563*c4762a1bSJed Brown ierr = VecNorm(zdiff,NORM_2,&s_norm);CHKERRQ(ierr); 564*c4762a1bSJed Brown s_norm = s_norm * mu; 565*c4762a1bSJed Brown /* primal : sqrt(n)*ABSTOL + RELTOL*max(norm(x), norm(-z))*/ 566*c4762a1bSJed Brown ierr = VecNorm(xk,NORM_2,&x_norm);CHKERRQ(ierr); 567*c4762a1bSJed Brown ierr = VecNorm(z,NORM_2,&z_norm);CHKERRQ(ierr); 568*c4762a1bSJed Brown primal = PetscSqrtReal(ctx->n)*ctx->abstol + ctx->reltol*PetscMax(x_norm,z_norm); 569*c4762a1bSJed Brown /* Duality : sqrt(n)*ABSTOL + RELTOL*norm(mu*u)*/ 570*c4762a1bSJed Brown ierr = VecNorm(u,NORM_2,&u_norm);CHKERRQ(ierr); 571*c4762a1bSJed Brown dual = PetscSqrtReal(ctx->n)*ctx->abstol + ctx->reltol*u_norm*mu; 572*c4762a1bSJed Brown ierr = PetscPrintf(PetscObjectComm((PetscObject)tao1),"Iter %D : ||x-z||: %g, mu*||z-zold||: %g\n", i, (double) r_norm, (double) s_norm);CHKERRQ(ierr); 573*c4762a1bSJed Brown if (r_norm < primal && s_norm < dual) break; 574*c4762a1bSJed Brown } 575*c4762a1bSJed Brown ierr = VecCopy(xk, x);CHKERRQ(ierr); 576*c4762a1bSJed Brown ierr = TaoDestroy(&tao1);CHKERRQ(ierr); 577*c4762a1bSJed Brown ierr = TaoDestroy(&tao2);CHKERRQ(ierr); 578*c4762a1bSJed Brown PetscFunctionReturn(0); 579*c4762a1bSJed Brown } 580*c4762a1bSJed Brown 581*c4762a1bSJed Brown /* Second order Taylor remainder convergence test */ 582*c4762a1bSJed Brown static PetscErrorCode TaylorTest(UserCtx ctx, Tao tao, Vec x, PetscReal *C) 583*c4762a1bSJed Brown { 584*c4762a1bSJed Brown PetscReal h,J,temp; 585*c4762a1bSJed Brown PetscInt i,j; 586*c4762a1bSJed Brown PetscInt numValues; 587*c4762a1bSJed Brown PetscReal Jx,Jxhat_comp,Jxhat_pred; 588*c4762a1bSJed Brown PetscReal *Js, *hs; 589*c4762a1bSJed Brown PetscReal gdotdx; 590*c4762a1bSJed Brown PetscReal minrate = PETSC_MAX_REAL; 591*c4762a1bSJed Brown MPI_Comm comm = PetscObjectComm((PetscObject)x); 592*c4762a1bSJed Brown Vec g, dx, xhat; 593*c4762a1bSJed Brown PetscErrorCode ierr; 594*c4762a1bSJed Brown 595*c4762a1bSJed Brown PetscFunctionBegin; 596*c4762a1bSJed Brown ierr = VecDuplicate(x, &g);CHKERRQ(ierr); 597*c4762a1bSJed Brown ierr = VecDuplicate(x, &xhat);CHKERRQ(ierr); 598*c4762a1bSJed Brown /* choose a perturbation direction */ 599*c4762a1bSJed Brown ierr = VecDuplicate(x, &dx);CHKERRQ(ierr); 600*c4762a1bSJed Brown ierr = VecSetRandom(dx,ctx->rctx);CHKERRQ(ierr); 601*c4762a1bSJed Brown /* evaluate objective at x: J(x) */ 602*c4762a1bSJed Brown ierr = TaoComputeObjective(tao, x, &Jx);CHKERRQ(ierr); 603*c4762a1bSJed Brown /* evaluate gradient at x, save in vector g */ 604*c4762a1bSJed Brown ierr = TaoComputeGradient(tao, x, g);CHKERRQ(ierr); 605*c4762a1bSJed Brown ierr = VecDot(g, dx, &gdotdx);CHKERRQ(ierr); 606*c4762a1bSJed Brown 607*c4762a1bSJed Brown for (numValues=0, h=ctx->hStart; h>=ctx->hMin; h*=ctx->hFactor) numValues++; 608*c4762a1bSJed Brown ierr = PetscCalloc2(numValues, &Js, numValues, &hs);CHKERRQ(ierr); 609*c4762a1bSJed Brown for (i=0, h=ctx->hStart; h>=ctx->hMin; h*=ctx->hFactor, i++) { 610*c4762a1bSJed Brown ierr = VecWAXPY(xhat, h, dx, x);CHKERRQ(ierr); 611*c4762a1bSJed Brown ierr = TaoComputeObjective(tao, xhat, &Jxhat_comp);CHKERRQ(ierr); 612*c4762a1bSJed Brown /* J(\hat(x)) \approx J(x) + g^T (xhat - x) = J(x) + h * g^T dx */ 613*c4762a1bSJed Brown Jxhat_pred = Jx + h * gdotdx; 614*c4762a1bSJed Brown /* Vector to dJdm scalar? Dot?*/ 615*c4762a1bSJed Brown J = PetscAbsReal(Jxhat_comp - Jxhat_pred); 616*c4762a1bSJed Brown ierr = PetscPrintf (comm, "J(xhat): %g, predicted: %g, diff %g\n", (double) Jxhat_comp,(double) Jxhat_pred, (double) J);CHKERRQ(ierr); 617*c4762a1bSJed Brown Js[i] = J; 618*c4762a1bSJed Brown hs[i] = h; 619*c4762a1bSJed Brown } 620*c4762a1bSJed Brown for (j=1; j<numValues; j++) { 621*c4762a1bSJed Brown temp = PetscLogReal(Js[j]/Js[j-1]) / PetscLogReal (hs[j]/hs[j-1]); 622*c4762a1bSJed Brown ierr = PetscPrintf (comm, "Convergence rate step %D: %g\n", j-1, (double) temp);CHKERRQ(ierr); 623*c4762a1bSJed Brown minrate = PetscMin(minrate, temp); 624*c4762a1bSJed Brown } 625*c4762a1bSJed Brown /* If O is not ~2, then the test is wrong */ 626*c4762a1bSJed Brown ierr = PetscFree2(Js, hs);CHKERRQ(ierr); 627*c4762a1bSJed Brown *C = minrate; 628*c4762a1bSJed Brown ierr = VecDestroy(&dx);CHKERRQ(ierr); 629*c4762a1bSJed Brown ierr = VecDestroy(&xhat);CHKERRQ(ierr); 630*c4762a1bSJed Brown ierr = VecDestroy(&g);CHKERRQ(ierr); 631*c4762a1bSJed Brown PetscFunctionReturn(0); 632*c4762a1bSJed Brown } 633*c4762a1bSJed Brown 634*c4762a1bSJed Brown int main(int argc, char ** argv) 635*c4762a1bSJed Brown { 636*c4762a1bSJed Brown UserCtx ctx; 637*c4762a1bSJed Brown Tao tao; 638*c4762a1bSJed Brown Vec x; 639*c4762a1bSJed Brown Mat H; 640*c4762a1bSJed Brown PetscErrorCode ierr; 641*c4762a1bSJed Brown 642*c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 643*c4762a1bSJed Brown ierr = PetscNew(&ctx);CHKERRQ(ierr); 644*c4762a1bSJed Brown ierr = ConfigureContext(ctx);CHKERRQ(ierr); 645*c4762a1bSJed Brown /* Define two functions that could pass as objectives to TaoSetObjectiveRoutine(): one 646*c4762a1bSJed Brown * for the misfit component, and one for the regularization component */ 647*c4762a1bSJed Brown /* ObjectiveMisfit() and ObjectiveRegularization() */ 648*c4762a1bSJed Brown 649*c4762a1bSJed Brown /* Define a single function that calls both components adds them together: the complete objective, 650*c4762a1bSJed Brown * in the absence of a Tao implementation that handles separability */ 651*c4762a1bSJed Brown /* ObjectiveComplete() */ 652*c4762a1bSJed Brown ierr = TaoCreate(PETSC_COMM_WORLD, &tao);CHKERRQ(ierr); 653*c4762a1bSJed Brown ierr = TaoSetType(tao,TAONM);CHKERRQ(ierr); 654*c4762a1bSJed Brown ierr = TaoSetObjectiveRoutine(tao, ObjectiveComplete, (void*) ctx);CHKERRQ(ierr); 655*c4762a1bSJed Brown ierr = TaoSetGradientRoutine(tao, GradientComplete, (void*) ctx);CHKERRQ(ierr); 656*c4762a1bSJed Brown ierr = MatDuplicate(ctx->W, MAT_SHARE_NONZERO_PATTERN, &H);CHKERRQ(ierr); 657*c4762a1bSJed Brown ierr = TaoSetHessianRoutine(tao, H, H, HessianComplete, (void*) ctx);CHKERRQ(ierr); 658*c4762a1bSJed Brown ierr = MatCreateVecs(ctx->F, NULL, &x);CHKERRQ(ierr); 659*c4762a1bSJed Brown ierr = VecSet(x, 0.);CHKERRQ(ierr); 660*c4762a1bSJed Brown ierr = TaoSetInitialVector(tao, x);CHKERRQ(ierr); 661*c4762a1bSJed Brown ierr = TaoSetFromOptions(tao);CHKERRQ(ierr); 662*c4762a1bSJed Brown if (ctx->use_admm) { 663*c4762a1bSJed Brown ierr = TaoSolveADMM(ctx,x); CHKERRQ(ierr); 664*c4762a1bSJed Brown } else {ierr = TaoSolve(tao);CHKERRQ(ierr);} 665*c4762a1bSJed Brown /* examine solution */ 666*c4762a1bSJed Brown ierr = VecViewFromOptions(x, NULL, "-view_sol");CHKERRQ(ierr); 667*c4762a1bSJed Brown if (ctx->taylor) { 668*c4762a1bSJed Brown PetscReal rate; 669*c4762a1bSJed Brown ierr = TaylorTest(ctx, tao, x, &rate);CHKERRQ(ierr); 670*c4762a1bSJed Brown } 671*c4762a1bSJed Brown ierr = MatDestroy(&H);CHKERRQ(ierr); 672*c4762a1bSJed Brown ierr = TaoDestroy(&tao);CHKERRQ(ierr); 673*c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 674*c4762a1bSJed Brown ierr = DestroyContext(&ctx);CHKERRQ(ierr); 675*c4762a1bSJed Brown ierr = PetscFinalize();CHKERRQ(ierr); 676*c4762a1bSJed Brown return ierr; 677*c4762a1bSJed Brown } 678*c4762a1bSJed Brown 679*c4762a1bSJed Brown /*TEST 680*c4762a1bSJed Brown 681*c4762a1bSJed Brown build: 682*c4762a1bSJed Brown requires: !complex 683*c4762a1bSJed Brown 684*c4762a1bSJed Brown test: 685*c4762a1bSJed Brown suffix: 0 686*c4762a1bSJed Brown args: 687*c4762a1bSJed Brown 688*c4762a1bSJed Brown test: 689*c4762a1bSJed Brown suffix: l1_1 690*c4762a1bSJed Brown args: -p 1 -tao_type lmvm -alpha 1. -epsilon 1.e-7 -m 64 -n 64 -view_sol -matrix_format 1 691*c4762a1bSJed Brown 692*c4762a1bSJed Brown test: 693*c4762a1bSJed Brown suffix: hessian_1 694*c4762a1bSJed Brown args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 1 -tao_type nls -tao_nls_ksp_monitor 695*c4762a1bSJed Brown 696*c4762a1bSJed Brown test: 697*c4762a1bSJed Brown suffix: hessian_2 698*c4762a1bSJed Brown args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 2 -tao_type nls -tao_nls_ksp_monitor 699*c4762a1bSJed Brown 700*c4762a1bSJed Brown test: 701*c4762a1bSJed Brown suffix: nm_1 702*c4762a1bSJed Brown args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 1 -tao_type nm -tao_max_it 50 703*c4762a1bSJed Brown 704*c4762a1bSJed Brown test: 705*c4762a1bSJed Brown suffix: nm_2 706*c4762a1bSJed Brown args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 2 -tao_type nm -tao_max_it 50 707*c4762a1bSJed Brown 708*c4762a1bSJed Brown test: 709*c4762a1bSJed Brown suffix: lmvm_1 710*c4762a1bSJed Brown args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 1 -tao_type lmvm -tao_max_it 40 711*c4762a1bSJed Brown 712*c4762a1bSJed Brown test: 713*c4762a1bSJed Brown suffix: lmvm_2 714*c4762a1bSJed Brown args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 2 -tao_type lmvm -tao_max_it 15 715*c4762a1bSJed Brown 716*c4762a1bSJed Brown test: 717*c4762a1bSJed Brown suffix: soft_threshold_admm_1 718*c4762a1bSJed Brown args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 1 -use_admm 719*c4762a1bSJed Brown 720*c4762a1bSJed Brown test: 721*c4762a1bSJed Brown suffix: hessian_admm_1 722*c4762a1bSJed Brown args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 1 -use_admm -reg_tao_type nls -misfit_tao_type nls 723*c4762a1bSJed Brown 724*c4762a1bSJed Brown test: 725*c4762a1bSJed Brown suffix: hessian_admm_2 726*c4762a1bSJed Brown args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 2 -use_admm -reg_tao_type nls -misfit_tao_type nls 727*c4762a1bSJed Brown 728*c4762a1bSJed Brown test: 729*c4762a1bSJed Brown suffix: nm_admm_1 730*c4762a1bSJed Brown args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 1 -use_admm -reg_tao_type nm -misfit_tao_type nm 731*c4762a1bSJed Brown 732*c4762a1bSJed Brown test: 733*c4762a1bSJed Brown suffix: nm_admm_2 734*c4762a1bSJed Brown args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 2 -use_admm -reg_tao_type nm -misfit_tao_type nm -iter 7 735*c4762a1bSJed Brown 736*c4762a1bSJed Brown test: 737*c4762a1bSJed Brown suffix: lmvm_admm_1 738*c4762a1bSJed Brown args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 1 -use_admm -reg_tao_type lmvm -misfit_tao_type lmvm 739*c4762a1bSJed Brown 740*c4762a1bSJed Brown test: 741*c4762a1bSJed Brown suffix: lmvm_admm_2 742*c4762a1bSJed Brown args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 2 -use_admm -reg_tao_type lmvm -misfit_tao_type lmvm 743*c4762a1bSJed Brown 744*c4762a1bSJed Brown TEST*/ 745