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 */ 429566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&(ctx->d))); 439566063dSJacob Faibussowitsch PetscCall(VecSetSizes(ctx->d,PETSC_DECIDE,ctx->m)); 449566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(ctx->d)); 459566063dSJacob Faibussowitsch PetscCall(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 */ 589566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &(ctx->F))); 599566063dSJacob Faibussowitsch PetscCall(MatSetSizes(ctx->F,PETSC_DECIDE, PETSC_DECIDE, ctx->m, ctx->n)); 609566063dSJacob Faibussowitsch PetscCall(MatSetType(ctx->F,MATAIJ)); /* TODO: Decide specific SetType other than dummy*/ 619566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(ctx->F, 5, NULL, 5, NULL)); /*TODO: some number other than 5?*/ 629566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(ctx->F, 5, NULL)); 639566063dSJacob Faibussowitsch PetscCall(MatSetUp(ctx->F)); 649566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(ctx->F,&Istart,&Iend)); 659566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("Assembly", &stage)); 669566063dSJacob Faibussowitsch PetscCall(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; 839566063dSJacob Faibussowitsch PetscCall(MatSetValue(ctx->F, Ii, Ii, 4., INSERT_VALUES)); 849566063dSJacob Faibussowitsch PetscCall(MatSetValue(ctx->F, Ii, I_n, -1., INSERT_VALUES)); 859566063dSJacob Faibussowitsch PetscCall(MatSetValue(ctx->F, Ii, I_s, -1., INSERT_VALUES)); 869566063dSJacob Faibussowitsch PetscCall(MatSetValue(ctx->F, Ii, I_e, -1., INSERT_VALUES)); 879566063dSJacob Faibussowitsch PetscCall(MatSetValue(ctx->F, Ii, I_w, -1., INSERT_VALUES)); 88c4762a1bSJed Brown } 899566063dSJacob Faibussowitsch } else PetscCall(MatSetRandom(ctx->F, ctx->rctx)); 909566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(ctx->F, MAT_FINAL_ASSEMBLY)); 919566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(ctx->F, MAT_FINAL_ASSEMBLY)); 929566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 93c4762a1bSJed Brown /* Stencil matrix is symmetric. Setting symmetric flag for ICC/Cholesky preconditioner */ 94c4762a1bSJed Brown if (!(ctx->matops)) { 959566063dSJacob Faibussowitsch PetscCall(MatSetOption(ctx->F,MAT_SYMMETRIC,PETSC_TRUE)); 96c4762a1bSJed Brown } 979566063dSJacob Faibussowitsch PetscCall(MatTransposeMatMult(ctx->F,ctx->F, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &(ctx->W))); 98c4762a1bSJed Brown /* Setup Hessian Workspace in same shape as W */ 999566063dSJacob Faibussowitsch PetscCall(MatDuplicate(ctx->W,MAT_DO_NOT_COPY_VALUES,&(ctx->Hm))); 1009566063dSJacob Faibussowitsch PetscCall(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; 1099566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(ctx->F, &ctx->workLeft[0], &ctx->workRight[0])); 110c4762a1bSJed Brown for (i=1; i<NWORKLEFT; i++) { 1119566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ctx->workLeft[0], &(ctx->workLeft[i]))); 112c4762a1bSJed Brown } 113c4762a1bSJed Brown for (i=1; i<NWORKRIGHT; i++) { 1149566063dSJacob Faibussowitsch PetscCall(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 PetscFunctionBegin; 122c4762a1bSJed Brown ctx->m = 16; 123c4762a1bSJed Brown ctx->n = 16; 124c4762a1bSJed Brown ctx->eps = 1.e-3; 125c4762a1bSJed Brown ctx->abstol = 1.e-4; 126c4762a1bSJed Brown ctx->reltol = 1.e-2; 127c4762a1bSJed Brown ctx->hStart = 1.; 128c4762a1bSJed Brown ctx->hMin = 1.e-3; 129c4762a1bSJed Brown ctx->hFactor = 0.5; 130c4762a1bSJed Brown ctx->alpha = 1.; 131c4762a1bSJed Brown ctx->mu = 1.0; 132c4762a1bSJed Brown ctx->matops = 0; 133c4762a1bSJed Brown ctx->iter = 10; 134c4762a1bSJed Brown ctx->p = NORM_2; 135c4762a1bSJed Brown ctx->taylor = PETSC_TRUE; 136c4762a1bSJed Brown ctx->use_admm = PETSC_FALSE; 137d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Configure separable objection example", "ex4.c"); 1389566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-m", "The row dimension of matrix F", "ex4.c", ctx->m, &(ctx->m), NULL)); 1399566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-n", "The column dimension of matrix F", "ex4.c", ctx->n, &(ctx->n), NULL)); 1409566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-matrix_format","Decide format of F matrix. 0 for stencil, 1 for random", "ex4.c", ctx->matops, &(ctx->matops), NULL)); 1419566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-iter","Iteration number ADMM", "ex4.c", ctx->iter, &(ctx->iter), NULL)); 1429566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-alpha", "The regularization multiplier. 1 default", "ex4.c", ctx->alpha, &(ctx->alpha), NULL)); 1439566063dSJacob Faibussowitsch PetscCall(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)); 1449566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-mu", "The augmented lagrangian multiplier in ADMM", "ex4.c", ctx->mu, &(ctx->mu), NULL)); 1459566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-hStart", "Taylor test starting point. 1 default.", "ex4.c", ctx->hStart, &(ctx->hStart), NULL)); 1469566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-hFactor", "Taylor test multiplier factor. 0.5 default", "ex4.c", ctx->hFactor, &(ctx->hFactor), NULL)); 1479566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-hMin", "Taylor test ending condition. 1.e-3 default", "ex4.c", ctx->hMin, &(ctx->hMin), NULL)); 1489566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-abstol", "Absolute stopping criterion for ADMM", "ex4.c", ctx->abstol, &(ctx->abstol), NULL)); 1499566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-reltol", "Relative stopping criterion for ADMM", "ex4.c", ctx->reltol, &(ctx->reltol), NULL)); 1509566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-taylor","Flag for Taylor test. Default is true.", "ex4.c", ctx->taylor, &(ctx->taylor), NULL)); 1519566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-use_admm","Use the ADMM solver in this example.", "ex4.c", ctx->use_admm, &(ctx->use_admm), NULL)); 1529566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-p","Norm type.", "ex4.c", NormTypes, (PetscEnum)ctx->p, (PetscEnum *) &(ctx->p), NULL)); 153d0609cedSBarry Smith PetscOptionsEnd(); 154c4762a1bSJed Brown /* Creating random ctx */ 1559566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&(ctx->rctx))); 1569566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(ctx->rctx)); 1579566063dSJacob Faibussowitsch PetscCall(CreateMatrix(ctx)); 1589566063dSJacob Faibussowitsch PetscCall(CreateRHS(ctx)); 1599566063dSJacob Faibussowitsch PetscCall(SetupWorkspace(ctx)); 160c4762a1bSJed Brown PetscFunctionReturn(0); 161c4762a1bSJed Brown } 162c4762a1bSJed Brown 163c4762a1bSJed Brown static PetscErrorCode DestroyContext(UserCtx *ctx) 164c4762a1bSJed Brown { 165c4762a1bSJed Brown PetscInt i; 166c4762a1bSJed Brown 167c4762a1bSJed Brown PetscFunctionBegin; 1689566063dSJacob Faibussowitsch PetscCall(MatDestroy(&((*ctx)->F))); 1699566063dSJacob Faibussowitsch PetscCall(MatDestroy(&((*ctx)->W))); 1709566063dSJacob Faibussowitsch PetscCall(MatDestroy(&((*ctx)->Hm))); 1719566063dSJacob Faibussowitsch PetscCall(MatDestroy(&((*ctx)->Hr))); 1729566063dSJacob Faibussowitsch PetscCall(VecDestroy(&((*ctx)->d))); 173c4762a1bSJed Brown for (i=0; i<NWORKLEFT; i++) { 1749566063dSJacob Faibussowitsch PetscCall(VecDestroy(&((*ctx)->workLeft[i]))); 175c4762a1bSJed Brown } 176c4762a1bSJed Brown for (i=0; i<NWORKRIGHT; i++) { 1779566063dSJacob Faibussowitsch PetscCall(VecDestroy(&((*ctx)->workRight[i]))); 178c4762a1bSJed Brown } 1799566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&((*ctx)->rctx))); 1809566063dSJacob Faibussowitsch PetscCall(PetscFree(*ctx)); 181c4762a1bSJed Brown PetscFunctionReturn(0); 182c4762a1bSJed Brown } 183c4762a1bSJed Brown 184c4762a1bSJed Brown /* compute (1/2) * ||F x - d||^2 */ 185c4762a1bSJed Brown static PetscErrorCode ObjectiveMisfit(Tao tao, Vec x, PetscReal *J, void *_ctx) 186c4762a1bSJed Brown { 187c4762a1bSJed Brown UserCtx ctx = (UserCtx) _ctx; 188c4762a1bSJed Brown Vec y; 189c4762a1bSJed Brown 190c4762a1bSJed Brown PetscFunctionBegin; 191c4762a1bSJed Brown y = ctx->workLeft[0]; 1929566063dSJacob Faibussowitsch PetscCall(MatMult(ctx->F, x, y)); 1939566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1., ctx->d)); 1949566063dSJacob Faibussowitsch PetscCall(VecDot(y, y, J)); 195c4762a1bSJed Brown *J *= 0.5; 196c4762a1bSJed Brown PetscFunctionReturn(0); 197c4762a1bSJed Brown } 198c4762a1bSJed Brown 199c4762a1bSJed Brown /* compute V = FTFx - FTd */ 200c4762a1bSJed Brown static PetscErrorCode GradientMisfit(Tao tao, Vec x, Vec V, void *_ctx) 201c4762a1bSJed Brown { 202c4762a1bSJed Brown UserCtx ctx = (UserCtx) _ctx; 203c4762a1bSJed Brown Vec FTFx, FTd; 204c4762a1bSJed Brown 205c4762a1bSJed Brown PetscFunctionBegin; 206c4762a1bSJed Brown /* work1 is A^T Ax, work2 is Ab, W is A^T A*/ 207c4762a1bSJed Brown FTFx = ctx->workRight[0]; 208c4762a1bSJed Brown FTd = ctx->workRight[1]; 2099566063dSJacob Faibussowitsch PetscCall(MatMult(ctx->W,x,FTFx)); 2109566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(ctx->F, ctx->d, FTd)); 2119566063dSJacob Faibussowitsch PetscCall(VecWAXPY(V, -1., FTd, FTFx)); 212c4762a1bSJed Brown PetscFunctionReturn(0); 213c4762a1bSJed Brown } 214c4762a1bSJed Brown 215c4762a1bSJed Brown /* returns FTF */ 216c4762a1bSJed Brown static PetscErrorCode HessianMisfit(Tao tao, Vec x, Mat H, Mat Hpre, void *_ctx) 217c4762a1bSJed Brown { 218c4762a1bSJed Brown UserCtx ctx = (UserCtx) _ctx; 219c4762a1bSJed Brown 220c4762a1bSJed Brown PetscFunctionBegin; 2219566063dSJacob Faibussowitsch if (H != ctx->W) PetscCall(MatCopy(ctx->W, H, DIFFERENT_NONZERO_PATTERN)); 2229566063dSJacob Faibussowitsch if (Hpre != ctx->W) PetscCall(MatCopy(ctx->W, Hpre, DIFFERENT_NONZERO_PATTERN)); 223c4762a1bSJed Brown PetscFunctionReturn(0); 224c4762a1bSJed Brown } 225c4762a1bSJed Brown 226c4762a1bSJed Brown /* computes augment Lagrangian objective (with scaled dual): 227c4762a1bSJed Brown * 0.5 * ||F x - d||^2 + 0.5 * mu ||x - z + u||^2 */ 228c4762a1bSJed Brown static PetscErrorCode ObjectiveMisfitADMM(Tao tao, Vec x, PetscReal *J, void *_ctx) 229c4762a1bSJed Brown { 230c4762a1bSJed Brown UserCtx ctx = (UserCtx) _ctx; 231c4762a1bSJed Brown PetscReal mu, workNorm, misfit; 232c4762a1bSJed Brown Vec z, u, temp; 233c4762a1bSJed Brown 234c4762a1bSJed Brown PetscFunctionBegin; 235c4762a1bSJed Brown mu = ctx->mu; 236c4762a1bSJed Brown z = ctx->workRight[5]; 237c4762a1bSJed Brown u = ctx->workRight[6]; 238c4762a1bSJed Brown temp = ctx->workRight[10]; 239c4762a1bSJed Brown /* misfit = f(x) */ 2409566063dSJacob Faibussowitsch PetscCall(ObjectiveMisfit(tao, x, &misfit, _ctx)); 2419566063dSJacob Faibussowitsch PetscCall(VecCopy(x,temp)); 242c4762a1bSJed Brown /* temp = x - z + u */ 2439566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(temp,-1.,1.,1.,z,u)); 244c4762a1bSJed Brown /* workNorm = ||x - z + u||^2 */ 2459566063dSJacob Faibussowitsch PetscCall(VecDot(temp, temp, &workNorm)); 246c4762a1bSJed Brown /* augment Lagrangian objective (with scaled dual): f(x) + 0.5 * mu ||x -z + u||^2 */ 247c4762a1bSJed Brown *J = misfit + 0.5 * mu * workNorm; 248c4762a1bSJed Brown PetscFunctionReturn(0); 249c4762a1bSJed Brown } 250c4762a1bSJed Brown 251c4762a1bSJed Brown /* computes FTFx - FTd mu*(x - z + u) */ 252c4762a1bSJed Brown static PetscErrorCode GradientMisfitADMM(Tao tao, Vec x, Vec V, void *_ctx) 253c4762a1bSJed Brown { 254c4762a1bSJed Brown UserCtx ctx = (UserCtx) _ctx; 255c4762a1bSJed Brown PetscReal mu; 256c4762a1bSJed Brown Vec z, u, temp; 257c4762a1bSJed Brown 258c4762a1bSJed Brown PetscFunctionBegin; 259c4762a1bSJed Brown mu = ctx->mu; 260c4762a1bSJed Brown z = ctx->workRight[5]; 261c4762a1bSJed Brown u = ctx->workRight[6]; 262c4762a1bSJed Brown temp = ctx->workRight[10]; 2639566063dSJacob Faibussowitsch PetscCall(GradientMisfit(tao, x, V, _ctx)); 2649566063dSJacob Faibussowitsch PetscCall(VecCopy(x, temp)); 265c4762a1bSJed Brown /* temp = x - z + u */ 2669566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(temp,-1.,1.,1.,z,u)); 267c4762a1bSJed Brown /* V = FTFx - FTd mu*(x - z + u) */ 2689566063dSJacob Faibussowitsch PetscCall(VecAXPY(V, mu, temp)); 269c4762a1bSJed Brown PetscFunctionReturn(0); 270c4762a1bSJed Brown } 271c4762a1bSJed Brown 272c4762a1bSJed Brown /* returns FTF + diag(mu) */ 273c4762a1bSJed Brown static PetscErrorCode HessianMisfitADMM(Tao tao, Vec x, Mat H, Mat Hpre, void *_ctx) 274c4762a1bSJed Brown { 275c4762a1bSJed Brown UserCtx ctx = (UserCtx) _ctx; 276c4762a1bSJed Brown 277c4762a1bSJed Brown PetscFunctionBegin; 2789566063dSJacob Faibussowitsch PetscCall(MatCopy(ctx->W, H, DIFFERENT_NONZERO_PATTERN)); 2799566063dSJacob Faibussowitsch PetscCall(MatShift(H, ctx->mu)); 280c4762a1bSJed Brown if (Hpre != H) { 2819566063dSJacob Faibussowitsch PetscCall(MatCopy(H, Hpre, DIFFERENT_NONZERO_PATTERN)); 282c4762a1bSJed Brown } 283c4762a1bSJed Brown PetscFunctionReturn(0); 284c4762a1bSJed Brown } 285c4762a1bSJed Brown 286c4762a1bSJed Brown /* computes || x ||_p (mult by 0.5 in case of NORM_2) */ 287c4762a1bSJed Brown static PetscErrorCode ObjectiveRegularization(Tao tao, Vec x, PetscReal *J, void *_ctx) 288c4762a1bSJed Brown { 289c4762a1bSJed Brown UserCtx ctx = (UserCtx) _ctx; 290c4762a1bSJed Brown PetscReal norm; 291c4762a1bSJed Brown 292c4762a1bSJed Brown PetscFunctionBegin; 293c4762a1bSJed Brown *J = 0; 2949566063dSJacob Faibussowitsch PetscCall(VecNorm (x, ctx->p, &norm)); 295c4762a1bSJed Brown if (ctx->p == NORM_2) norm = 0.5 * norm * norm; 296c4762a1bSJed Brown *J = ctx->alpha * norm; 297c4762a1bSJed Brown PetscFunctionReturn(0); 298c4762a1bSJed Brown } 299c4762a1bSJed Brown 300c4762a1bSJed Brown /* NORM_2 Case: return x 301c4762a1bSJed Brown * NORM_1 Case: x/(|x| + eps) 302c4762a1bSJed Brown * Else: TODO */ 303c4762a1bSJed Brown static PetscErrorCode GradientRegularization(Tao tao, Vec x, Vec V, void *_ctx) 304c4762a1bSJed Brown { 305c4762a1bSJed Brown UserCtx ctx = (UserCtx) _ctx; 306c4762a1bSJed Brown PetscReal eps = ctx->eps; 307c4762a1bSJed Brown 308c4762a1bSJed Brown PetscFunctionBegin; 309c4762a1bSJed Brown if (ctx->p == NORM_2) { 3109566063dSJacob Faibussowitsch PetscCall(VecCopy(x, V)); 311c4762a1bSJed Brown } else if (ctx->p == NORM_1) { 3129566063dSJacob Faibussowitsch PetscCall(VecCopy(x, ctx->workRight[1])); 3139566063dSJacob Faibussowitsch PetscCall(VecAbs(ctx->workRight[1])); 3149566063dSJacob Faibussowitsch PetscCall(VecShift(ctx->workRight[1], eps)); 3159566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(V, x, ctx->workRight[1])); 316c4762a1bSJed Brown } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Example only works for NORM_1 and NORM_2"); 317c4762a1bSJed Brown PetscFunctionReturn(0); 318c4762a1bSJed Brown } 319c4762a1bSJed Brown 320c4762a1bSJed Brown /* NORM_2 Case: returns diag(mu) 321c4762a1bSJed Brown * NORM_1 Case: diag(mu* 1/sqrt(x_i^2 + eps) * (1 - x_i^2/ABS(x_i^2+eps))) */ 322c4762a1bSJed Brown static PetscErrorCode HessianRegularization(Tao tao, Vec x, Mat H, Mat Hpre, void *_ctx) 323c4762a1bSJed Brown { 324c4762a1bSJed Brown UserCtx ctx = (UserCtx) _ctx; 325c4762a1bSJed Brown PetscReal eps = ctx->eps; 326c4762a1bSJed Brown Vec copy1,copy2,copy3; 327c4762a1bSJed Brown 328c4762a1bSJed Brown PetscFunctionBegin; 329c4762a1bSJed Brown if (ctx->p == NORM_2) { 330c4762a1bSJed Brown /* Identity matrix scaled by mu */ 3319566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(H)); 3329566063dSJacob Faibussowitsch PetscCall(MatShift(H,ctx->mu)); 333c4762a1bSJed Brown if (Hpre != H) { 3349566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(Hpre)); 3359566063dSJacob Faibussowitsch PetscCall(MatShift(Hpre,ctx->mu)); 336c4762a1bSJed Brown } 337c4762a1bSJed Brown } else if (ctx->p == NORM_1) { 338c4762a1bSJed Brown /* 1/sqrt(x_i^2 + eps) * (1 - x_i^2/ABS(x_i^2+eps)) */ 339c4762a1bSJed Brown copy1 = ctx->workRight[1]; 340c4762a1bSJed Brown copy2 = ctx->workRight[2]; 341c4762a1bSJed Brown copy3 = ctx->workRight[3]; 342c4762a1bSJed Brown /* copy1 : 1/sqrt(x_i^2 + eps) */ 3439566063dSJacob Faibussowitsch PetscCall(VecCopy(x, copy1)); 3449566063dSJacob Faibussowitsch PetscCall(VecPow(copy1,2)); 3459566063dSJacob Faibussowitsch PetscCall(VecShift(copy1, eps)); 3469566063dSJacob Faibussowitsch PetscCall(VecSqrtAbs(copy1)); 3479566063dSJacob Faibussowitsch PetscCall(VecReciprocal(copy1)); 348c4762a1bSJed Brown /* copy2: x_i^2.*/ 3499566063dSJacob Faibussowitsch PetscCall(VecCopy(x,copy2)); 3509566063dSJacob Faibussowitsch PetscCall(VecPow(copy2,2)); 351c4762a1bSJed Brown /* copy3: abs(x_i^2 + eps) */ 3529566063dSJacob Faibussowitsch PetscCall(VecCopy(x,copy3)); 3539566063dSJacob Faibussowitsch PetscCall(VecPow(copy3,2)); 3549566063dSJacob Faibussowitsch PetscCall(VecShift(copy3, eps)); 3559566063dSJacob Faibussowitsch PetscCall(VecAbs(copy3)); 356c4762a1bSJed Brown /* copy2: 1 - x_i^2/abs(x_i^2 + eps) */ 3579566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(copy2, copy2,copy3)); 3589566063dSJacob Faibussowitsch PetscCall(VecScale(copy2, -1.)); 3599566063dSJacob Faibussowitsch PetscCall(VecShift(copy2, 1.)); 3609566063dSJacob Faibussowitsch PetscCall(VecAXPY(copy1,1.,copy2)); 3619566063dSJacob Faibussowitsch PetscCall(VecScale(copy1, ctx->mu)); 3629566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(H)); 3639566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(H, copy1,INSERT_VALUES)); 364c4762a1bSJed Brown if (Hpre != H) { 3659566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(Hpre)); 3669566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(Hpre, copy1,INSERT_VALUES)); 367c4762a1bSJed Brown } 368c4762a1bSJed Brown } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Example only works for NORM_1 and NORM_2"); 369c4762a1bSJed Brown PetscFunctionReturn(0); 370c4762a1bSJed Brown } 371c4762a1bSJed Brown 372c4762a1bSJed Brown /* NORM_2 Case: 0.5 || x ||_2 + 0.5 * mu * ||x + u - z||^2 373c4762a1bSJed Brown * Else : || x ||_2 + 0.5 * mu * ||x + u - z||^2 */ 374c4762a1bSJed Brown static PetscErrorCode ObjectiveRegularizationADMM(Tao tao, Vec z, PetscReal *J, void *_ctx) 375c4762a1bSJed Brown { 376c4762a1bSJed Brown UserCtx ctx = (UserCtx) _ctx; 377c4762a1bSJed Brown PetscReal mu, workNorm, reg; 378c4762a1bSJed Brown Vec x, u, temp; 379c4762a1bSJed Brown 380c4762a1bSJed Brown PetscFunctionBegin; 381c4762a1bSJed Brown mu = ctx->mu; 382c4762a1bSJed Brown x = ctx->workRight[4]; 383c4762a1bSJed Brown u = ctx->workRight[6]; 384c4762a1bSJed Brown temp = ctx->workRight[10]; 3859566063dSJacob Faibussowitsch PetscCall(ObjectiveRegularization(tao, z, ®, _ctx)); 3869566063dSJacob Faibussowitsch PetscCall(VecCopy(z,temp)); 387c4762a1bSJed Brown /* temp = x + u -z */ 3889566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(temp,1.,1.,-1.,x,u)); 389c4762a1bSJed Brown /* workNorm = ||x + u - z ||^2 */ 3909566063dSJacob Faibussowitsch PetscCall(VecDot(temp, temp, &workNorm)); 391c4762a1bSJed Brown *J = reg + 0.5 * mu * workNorm; 392c4762a1bSJed Brown PetscFunctionReturn(0); 393c4762a1bSJed Brown } 394c4762a1bSJed Brown 395c4762a1bSJed Brown /* NORM_2 Case: x - mu*(x + u - z) 396c4762a1bSJed Brown * NORM_1 Case: x/(|x| + eps) - mu*(x + u - z) 397c4762a1bSJed Brown * Else: TODO */ 398c4762a1bSJed Brown static PetscErrorCode GradientRegularizationADMM(Tao tao, Vec z, Vec V, void *_ctx) 399c4762a1bSJed Brown { 400c4762a1bSJed Brown UserCtx ctx = (UserCtx) _ctx; 401c4762a1bSJed Brown PetscReal mu; 402c4762a1bSJed Brown Vec x, u, temp; 403c4762a1bSJed Brown 404c4762a1bSJed Brown PetscFunctionBegin; 405c4762a1bSJed Brown mu = ctx->mu; 406c4762a1bSJed Brown x = ctx->workRight[4]; 407c4762a1bSJed Brown u = ctx->workRight[6]; 408c4762a1bSJed Brown temp = ctx->workRight[10]; 4099566063dSJacob Faibussowitsch PetscCall(GradientRegularization(tao, z, V, _ctx)); 4109566063dSJacob Faibussowitsch PetscCall(VecCopy(z, temp)); 411c4762a1bSJed Brown /* temp = x + u -z */ 4129566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(temp,1.,1.,-1.,x,u)); 4139566063dSJacob Faibussowitsch PetscCall(VecAXPY(V, -mu, temp)); 414c4762a1bSJed Brown PetscFunctionReturn(0); 415c4762a1bSJed Brown } 416c4762a1bSJed Brown 417c4762a1bSJed Brown /* NORM_2 Case: returns diag(mu) 418c4762a1bSJed Brown * NORM_1 Case: FTF + diag(mu) */ 419c4762a1bSJed Brown static PetscErrorCode HessianRegularizationADMM(Tao tao, Vec x, Mat H, Mat Hpre, void *_ctx) 420c4762a1bSJed Brown { 421c4762a1bSJed Brown UserCtx ctx = (UserCtx) _ctx; 422c4762a1bSJed Brown 423c4762a1bSJed Brown PetscFunctionBegin; 424c4762a1bSJed Brown if (ctx->p == NORM_2) { 425c4762a1bSJed Brown /* Identity matrix scaled by mu */ 4269566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(H)); 4279566063dSJacob Faibussowitsch PetscCall(MatShift(H,ctx->mu)); 428c4762a1bSJed Brown if (Hpre != H) { 4299566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(Hpre)); 4309566063dSJacob Faibussowitsch PetscCall(MatShift(Hpre,ctx->mu)); 431c4762a1bSJed Brown } 432c4762a1bSJed Brown } else if (ctx->p == NORM_1) { 4339566063dSJacob Faibussowitsch PetscCall(HessianMisfit(tao, x, H, Hpre, (void*) ctx)); 4349566063dSJacob Faibussowitsch PetscCall(MatShift(H, ctx->mu)); 4359566063dSJacob Faibussowitsch if (Hpre != H) PetscCall(MatShift(Hpre, ctx->mu)); 436c4762a1bSJed Brown } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Example only works for NORM_1 and NORM_2"); 437c4762a1bSJed Brown PetscFunctionReturn(0); 438c4762a1bSJed Brown } 439c4762a1bSJed Brown 440c4762a1bSJed Brown /* NORM_2 Case : (1/2) * ||F x - d||^2 + 0.5 * || x ||_p 441c4762a1bSJed Brown * NORM_1 Case : (1/2) * ||F x - d||^2 + || x ||_p */ 442c4762a1bSJed Brown static PetscErrorCode ObjectiveComplete(Tao tao, Vec x, PetscReal *J, void *ctx) 443c4762a1bSJed Brown { 444c4762a1bSJed Brown PetscReal Jm, Jr; 445c4762a1bSJed Brown 446c4762a1bSJed Brown PetscFunctionBegin; 4479566063dSJacob Faibussowitsch PetscCall(ObjectiveMisfit(tao, x, &Jm, ctx)); 4489566063dSJacob Faibussowitsch PetscCall(ObjectiveRegularization(tao, x, &Jr, ctx)); 449c4762a1bSJed Brown *J = Jm + Jr; 450c4762a1bSJed Brown PetscFunctionReturn(0); 451c4762a1bSJed Brown } 452c4762a1bSJed Brown 453c4762a1bSJed Brown /* NORM_2 Case: FTFx - FTd + x 454c4762a1bSJed Brown * NORM_1 Case: FTFx - FTd + x/(|x| + eps) */ 455c4762a1bSJed Brown static PetscErrorCode GradientComplete(Tao tao, Vec x, Vec V, void *ctx) 456c4762a1bSJed Brown { 457c4762a1bSJed Brown UserCtx cntx = (UserCtx) ctx; 458c4762a1bSJed Brown 459c4762a1bSJed Brown PetscFunctionBegin; 4609566063dSJacob Faibussowitsch PetscCall(GradientMisfit(tao, x, cntx->workRight[2], ctx)); 4619566063dSJacob Faibussowitsch PetscCall(GradientRegularization(tao, x, cntx->workRight[3], ctx)); 4629566063dSJacob Faibussowitsch PetscCall(VecWAXPY(V,1,cntx->workRight[2],cntx->workRight[3])); 463c4762a1bSJed Brown PetscFunctionReturn(0); 464c4762a1bSJed Brown } 465c4762a1bSJed Brown 466c4762a1bSJed Brown /* NORM_2 Case: diag(mu) + FTF 467c4762a1bSJed Brown * NORM_1 Case: diag(mu* 1/sqrt(x_i^2 + eps) * (1 - x_i^2/ABS(x_i^2+eps))) + FTF */ 468c4762a1bSJed Brown static PetscErrorCode HessianComplete(Tao tao, Vec x, Mat H, Mat Hpre, void *ctx) 469c4762a1bSJed Brown { 470c4762a1bSJed Brown Mat tempH; 471c4762a1bSJed Brown 472c4762a1bSJed Brown PetscFunctionBegin; 4739566063dSJacob Faibussowitsch PetscCall(MatDuplicate(H, MAT_SHARE_NONZERO_PATTERN, &tempH)); 4749566063dSJacob Faibussowitsch PetscCall(HessianMisfit(tao, x, H, H, ctx)); 4759566063dSJacob Faibussowitsch PetscCall(HessianRegularization(tao, x, tempH, tempH, ctx)); 4769566063dSJacob Faibussowitsch PetscCall(MatAXPY(H, 1., tempH, DIFFERENT_NONZERO_PATTERN)); 477c4762a1bSJed Brown if (Hpre != H) { 4789566063dSJacob Faibussowitsch PetscCall(MatCopy(H, Hpre, DIFFERENT_NONZERO_PATTERN)); 479c4762a1bSJed Brown } 4809566063dSJacob Faibussowitsch PetscCall(MatDestroy(&tempH)); 481c4762a1bSJed Brown PetscFunctionReturn(0); 482c4762a1bSJed Brown } 483c4762a1bSJed Brown 484c4762a1bSJed Brown static PetscErrorCode TaoSolveADMM(UserCtx ctx, Vec x) 485c4762a1bSJed Brown { 486c4762a1bSJed Brown PetscInt i; 487c4762a1bSJed Brown PetscReal u_norm, r_norm, s_norm, primal, dual, x_norm, z_norm; 488c4762a1bSJed Brown Tao tao1,tao2; 489c4762a1bSJed Brown Vec xk,z,u,diff,zold,zdiff,temp; 490c4762a1bSJed Brown PetscReal mu; 491c4762a1bSJed Brown 492c4762a1bSJed Brown PetscFunctionBegin; 493c4762a1bSJed Brown xk = ctx->workRight[4]; 494c4762a1bSJed Brown z = ctx->workRight[5]; 495c4762a1bSJed Brown u = ctx->workRight[6]; 496c4762a1bSJed Brown diff = ctx->workRight[7]; 497c4762a1bSJed Brown zold = ctx->workRight[8]; 498c4762a1bSJed Brown zdiff = ctx->workRight[9]; 499c4762a1bSJed Brown temp = ctx->workRight[11]; 500c4762a1bSJed Brown mu = ctx->mu; 5019566063dSJacob Faibussowitsch PetscCall(VecSet(u, 0.)); 5029566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao1)); 5039566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao1,TAONLS)); 5049566063dSJacob Faibussowitsch PetscCall(TaoSetObjective(tao1, ObjectiveMisfitADMM, (void*) ctx)); 5059566063dSJacob Faibussowitsch PetscCall(TaoSetGradient(tao1, NULL, GradientMisfitADMM, (void*) ctx)); 5069566063dSJacob Faibussowitsch PetscCall(TaoSetHessian(tao1, ctx->Hm, ctx->Hm, HessianMisfitADMM, (void*) ctx)); 5079566063dSJacob Faibussowitsch PetscCall(VecSet(xk, 0.)); 5089566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao1, xk)); 5099566063dSJacob Faibussowitsch PetscCall(TaoSetOptionsPrefix(tao1, "misfit_")); 5109566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao1)); 5119566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao2)); 512c4762a1bSJed Brown if (ctx->p == NORM_2) { 5139566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao2,TAONLS)); 5149566063dSJacob Faibussowitsch PetscCall(TaoSetObjective(tao2, ObjectiveRegularizationADMM, (void*) ctx)); 5159566063dSJacob Faibussowitsch PetscCall(TaoSetGradient(tao2, NULL, GradientRegularizationADMM, (void*) ctx)); 5169566063dSJacob Faibussowitsch PetscCall(TaoSetHessian(tao2, ctx->Hr, ctx->Hr, HessianRegularizationADMM, (void*) ctx)); 517c4762a1bSJed Brown } 5189566063dSJacob Faibussowitsch PetscCall(VecSet(z, 0.)); 5199566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao2, z)); 5209566063dSJacob Faibussowitsch PetscCall(TaoSetOptionsPrefix(tao2, "reg_")); 5219566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao2)); 522c4762a1bSJed Brown 523c4762a1bSJed Brown for (i=0; i<ctx->iter; i++) { 5249566063dSJacob Faibussowitsch PetscCall(VecCopy(z,zold)); 5259566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao1)); /* Updates xk */ 526c4762a1bSJed Brown if (ctx->p == NORM_1) { 5279566063dSJacob Faibussowitsch PetscCall(VecWAXPY(temp,1.,xk,u)); 5289566063dSJacob Faibussowitsch PetscCall(TaoSoftThreshold(temp,-ctx->alpha/mu,ctx->alpha/mu,z)); 529c4762a1bSJed Brown } else { 5309566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao2)); /* Update zk */ 531c4762a1bSJed Brown } 532c4762a1bSJed Brown /* u = u + xk -z */ 5339566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(u,1.,-1.,1.,xk,z)); 534c4762a1bSJed Brown /* r_norm : norm(x-z) */ 5359566063dSJacob Faibussowitsch PetscCall(VecWAXPY(diff,-1.,z,xk)); 5369566063dSJacob Faibussowitsch PetscCall(VecNorm(diff,NORM_2,&r_norm)); 537c4762a1bSJed Brown /* s_norm : norm(-mu(z-zold)) */ 5389566063dSJacob Faibussowitsch PetscCall(VecWAXPY(zdiff, -1.,zold,z)); 5399566063dSJacob Faibussowitsch PetscCall(VecNorm(zdiff,NORM_2,&s_norm)); 540c4762a1bSJed Brown s_norm = s_norm * mu; 541c4762a1bSJed Brown /* primal : sqrt(n)*ABSTOL + RELTOL*max(norm(x), norm(-z))*/ 5429566063dSJacob Faibussowitsch PetscCall(VecNorm(xk,NORM_2,&x_norm)); 5439566063dSJacob Faibussowitsch PetscCall(VecNorm(z,NORM_2,&z_norm)); 544c4762a1bSJed Brown primal = PetscSqrtReal(ctx->n)*ctx->abstol + ctx->reltol*PetscMax(x_norm,z_norm); 545c4762a1bSJed Brown /* Duality : sqrt(n)*ABSTOL + RELTOL*norm(mu*u)*/ 5469566063dSJacob Faibussowitsch PetscCall(VecNorm(u,NORM_2,&u_norm)); 547c4762a1bSJed Brown dual = PetscSqrtReal(ctx->n)*ctx->abstol + ctx->reltol*u_norm*mu; 548*63a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)tao1),"Iter %" PetscInt_FMT " : ||x-z||: %g, mu*||z-zold||: %g\n", i, (double) r_norm, (double) s_norm)); 549c4762a1bSJed Brown if (r_norm < primal && s_norm < dual) break; 550c4762a1bSJed Brown } 5519566063dSJacob Faibussowitsch PetscCall(VecCopy(xk, x)); 5529566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao1)); 5539566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao2)); 554c4762a1bSJed Brown PetscFunctionReturn(0); 555c4762a1bSJed Brown } 556c4762a1bSJed Brown 557c4762a1bSJed Brown /* Second order Taylor remainder convergence test */ 558c4762a1bSJed Brown static PetscErrorCode TaylorTest(UserCtx ctx, Tao tao, Vec x, PetscReal *C) 559c4762a1bSJed Brown { 560c4762a1bSJed Brown PetscReal h,J,temp; 561c4762a1bSJed Brown PetscInt i,j; 562c4762a1bSJed Brown PetscInt numValues; 563c4762a1bSJed Brown PetscReal Jx,Jxhat_comp,Jxhat_pred; 564c4762a1bSJed Brown PetscReal *Js, *hs; 565c4762a1bSJed Brown PetscReal gdotdx; 566c4762a1bSJed Brown PetscReal minrate = PETSC_MAX_REAL; 567c4762a1bSJed Brown MPI_Comm comm = PetscObjectComm((PetscObject)x); 568c4762a1bSJed Brown Vec g, dx, xhat; 569c4762a1bSJed Brown 570c4762a1bSJed Brown PetscFunctionBegin; 5719566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &g)); 5729566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &xhat)); 573c4762a1bSJed Brown /* choose a perturbation direction */ 5749566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &dx)); 5759566063dSJacob Faibussowitsch PetscCall(VecSetRandom(dx,ctx->rctx)); 576c4762a1bSJed Brown /* evaluate objective at x: J(x) */ 5779566063dSJacob Faibussowitsch PetscCall(TaoComputeObjective(tao, x, &Jx)); 578c4762a1bSJed Brown /* evaluate gradient at x, save in vector g */ 5799566063dSJacob Faibussowitsch PetscCall(TaoComputeGradient(tao, x, g)); 5809566063dSJacob Faibussowitsch PetscCall(VecDot(g, dx, &gdotdx)); 581c4762a1bSJed Brown 582c4762a1bSJed Brown for (numValues=0, h=ctx->hStart; h>=ctx->hMin; h*=ctx->hFactor) numValues++; 5839566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(numValues, &Js, numValues, &hs)); 584c4762a1bSJed Brown for (i=0, h=ctx->hStart; h>=ctx->hMin; h*=ctx->hFactor, i++) { 5859566063dSJacob Faibussowitsch PetscCall(VecWAXPY(xhat, h, dx, x)); 5869566063dSJacob Faibussowitsch PetscCall(TaoComputeObjective(tao, xhat, &Jxhat_comp)); 587c4762a1bSJed Brown /* J(\hat(x)) \approx J(x) + g^T (xhat - x) = J(x) + h * g^T dx */ 588c4762a1bSJed Brown Jxhat_pred = Jx + h * gdotdx; 589c4762a1bSJed Brown /* Vector to dJdm scalar? Dot?*/ 590c4762a1bSJed Brown J = PetscAbsReal(Jxhat_comp - Jxhat_pred); 5919566063dSJacob Faibussowitsch PetscCall(PetscPrintf (comm, "J(xhat): %g, predicted: %g, diff %g\n", (double) Jxhat_comp,(double) Jxhat_pred, (double) J)); 592c4762a1bSJed Brown Js[i] = J; 593c4762a1bSJed Brown hs[i] = h; 594c4762a1bSJed Brown } 595c4762a1bSJed Brown for (j=1; j<numValues; j++) { 596c4762a1bSJed Brown temp = PetscLogReal(Js[j]/Js[j-1]) / PetscLogReal (hs[j]/hs[j-1]); 597*63a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf (comm, "Convergence rate step %" PetscInt_FMT ": %g\n", j-1, (double) temp)); 598c4762a1bSJed Brown minrate = PetscMin(minrate, temp); 599c4762a1bSJed Brown } 600c4762a1bSJed Brown /* If O is not ~2, then the test is wrong */ 6019566063dSJacob Faibussowitsch PetscCall(PetscFree2(Js, hs)); 602c4762a1bSJed Brown *C = minrate; 6039566063dSJacob Faibussowitsch PetscCall(VecDestroy(&dx)); 6049566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xhat)); 6059566063dSJacob Faibussowitsch PetscCall(VecDestroy(&g)); 606c4762a1bSJed Brown PetscFunctionReturn(0); 607c4762a1bSJed Brown } 608c4762a1bSJed Brown 609c4762a1bSJed Brown int main(int argc, char ** argv) 610c4762a1bSJed Brown { 611c4762a1bSJed Brown UserCtx ctx; 612c4762a1bSJed Brown Tao tao; 613c4762a1bSJed Brown Vec x; 614c4762a1bSJed Brown Mat H; 615c4762a1bSJed Brown 6169566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL,help)); 6179566063dSJacob Faibussowitsch PetscCall(PetscNew(&ctx)); 6189566063dSJacob Faibussowitsch PetscCall(ConfigureContext(ctx)); 619a82e8c82SStefano Zampini /* Define two functions that could pass as objectives to TaoSetObjective(): one 620c4762a1bSJed Brown * for the misfit component, and one for the regularization component */ 621c4762a1bSJed Brown /* ObjectiveMisfit() and ObjectiveRegularization() */ 622c4762a1bSJed Brown 623c4762a1bSJed Brown /* Define a single function that calls both components adds them together: the complete objective, 624c4762a1bSJed Brown * in the absence of a Tao implementation that handles separability */ 625c4762a1bSJed Brown /* ObjectiveComplete() */ 6269566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao)); 6279566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao,TAONM)); 6289566063dSJacob Faibussowitsch PetscCall(TaoSetObjective(tao, ObjectiveComplete, (void*) ctx)); 6299566063dSJacob Faibussowitsch PetscCall(TaoSetGradient(tao, NULL, GradientComplete, (void*) ctx)); 6309566063dSJacob Faibussowitsch PetscCall(MatDuplicate(ctx->W, MAT_SHARE_NONZERO_PATTERN, &H)); 6319566063dSJacob Faibussowitsch PetscCall(TaoSetHessian(tao, H, H, HessianComplete, (void*) ctx)); 6329566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(ctx->F, NULL, &x)); 6339566063dSJacob Faibussowitsch PetscCall(VecSet(x, 0.)); 6349566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, x)); 6359566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao)); 636c4762a1bSJed Brown if (ctx->use_admm) { 6379566063dSJacob Faibussowitsch PetscCall(TaoSolveADMM(ctx,x)); 6389566063dSJacob Faibussowitsch } else PetscCall(TaoSolve(tao)); 639c4762a1bSJed Brown /* examine solution */ 6409566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(x, NULL, "-view_sol")); 641c4762a1bSJed Brown if (ctx->taylor) { 642c4762a1bSJed Brown PetscReal rate; 6439566063dSJacob Faibussowitsch PetscCall(TaylorTest(ctx, tao, x, &rate)); 644c4762a1bSJed Brown } 6459566063dSJacob Faibussowitsch PetscCall(MatDestroy(&H)); 6469566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao)); 6479566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 6489566063dSJacob Faibussowitsch PetscCall(DestroyContext(&ctx)); 6499566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 650b122ec5aSJacob Faibussowitsch return 0; 651c4762a1bSJed Brown } 652c4762a1bSJed Brown 653c4762a1bSJed Brown /*TEST 654c4762a1bSJed Brown 655c4762a1bSJed Brown build: 656c4762a1bSJed Brown requires: !complex 657c4762a1bSJed Brown 658c4762a1bSJed Brown test: 659c4762a1bSJed Brown suffix: 0 660c4762a1bSJed Brown args: 661c4762a1bSJed Brown 662c4762a1bSJed Brown test: 663c4762a1bSJed Brown suffix: l1_1 664c4762a1bSJed Brown args: -p 1 -tao_type lmvm -alpha 1. -epsilon 1.e-7 -m 64 -n 64 -view_sol -matrix_format 1 665c4762a1bSJed Brown 666c4762a1bSJed Brown test: 667c4762a1bSJed Brown suffix: hessian_1 668c5f5e425SStefano Zampini args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 1 -tao_type nls 669c4762a1bSJed Brown 670c4762a1bSJed Brown test: 671c4762a1bSJed Brown suffix: hessian_2 672c5f5e425SStefano Zampini args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 2 -tao_type nls 673c4762a1bSJed Brown 674c4762a1bSJed Brown test: 675c4762a1bSJed Brown suffix: nm_1 676c4762a1bSJed Brown args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 1 -tao_type nm -tao_max_it 50 677c4762a1bSJed Brown 678c4762a1bSJed Brown test: 679c4762a1bSJed Brown suffix: nm_2 680c4762a1bSJed Brown args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 2 -tao_type nm -tao_max_it 50 681c4762a1bSJed Brown 682c4762a1bSJed Brown test: 683c4762a1bSJed Brown suffix: lmvm_1 684c4762a1bSJed Brown args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 1 -tao_type lmvm -tao_max_it 40 685c4762a1bSJed Brown 686c4762a1bSJed Brown test: 687c4762a1bSJed Brown suffix: lmvm_2 688c4762a1bSJed Brown args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 2 -tao_type lmvm -tao_max_it 15 689c4762a1bSJed Brown 690c4762a1bSJed Brown test: 691c4762a1bSJed Brown suffix: soft_threshold_admm_1 692c4762a1bSJed Brown args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 1 -use_admm 693c4762a1bSJed Brown 694c4762a1bSJed Brown test: 695c4762a1bSJed Brown suffix: hessian_admm_1 696c4762a1bSJed Brown args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 1 -use_admm -reg_tao_type nls -misfit_tao_type nls 697c4762a1bSJed Brown 698c4762a1bSJed Brown test: 699c4762a1bSJed Brown suffix: hessian_admm_2 700c4762a1bSJed Brown args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 2 -use_admm -reg_tao_type nls -misfit_tao_type nls 701c4762a1bSJed Brown 702c4762a1bSJed Brown test: 703c4762a1bSJed Brown suffix: nm_admm_1 704c4762a1bSJed Brown args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 1 -use_admm -reg_tao_type nm -misfit_tao_type nm 705c4762a1bSJed Brown 706c4762a1bSJed Brown test: 707c4762a1bSJed Brown suffix: nm_admm_2 708c4762a1bSJed 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 709c4762a1bSJed Brown 710c4762a1bSJed Brown test: 711c4762a1bSJed Brown suffix: lmvm_admm_1 712c4762a1bSJed Brown args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 1 -use_admm -reg_tao_type lmvm -misfit_tao_type lmvm 713c4762a1bSJed Brown 714c4762a1bSJed Brown test: 715c4762a1bSJed Brown suffix: lmvm_admm_2 716c4762a1bSJed Brown args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 2 -use_admm -reg_tao_type lmvm -misfit_tao_type lmvm 717c4762a1bSJed Brown 718c4762a1bSJed Brown TEST*/ 719