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