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