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 119371c9d4SSatish Balay typedef struct _UserCtx { 12c4762a1bSJed Brown PetscInt m; /* The row dimension of F */ 13c4762a1bSJed Brown PetscInt n; /* The column dimension of F */ 14c4762a1bSJed Brown PetscInt matops; /* Matrix format. 0 for stencil, 1 for random */ 15*be87f6c0SPierre Jolivet PetscInt iter; /* Number of iterations for ADMM */ 16c4762a1bSJed Brown PetscReal hStart; /* Starting point for Taylor test */ 17c4762a1bSJed Brown PetscReal hFactor; /* Taylor test step factor */ 18c4762a1bSJed Brown PetscReal hMin; /* Taylor test end goal */ 19c4762a1bSJed Brown PetscReal alpha; /* regularization constant applied to || x ||_p */ 20c4762a1bSJed Brown PetscReal eps; /* small constant for approximating gradient of || x ||_1 */ 21c4762a1bSJed Brown PetscReal mu; /* the augmented Lagrangian term in ADMM */ 22c4762a1bSJed Brown PetscReal abstol; 23c4762a1bSJed Brown PetscReal reltol; 24c4762a1bSJed Brown Mat F; /* matrix in least squares component $(1/2) * || F x - d ||_2^2$ */ 25c4762a1bSJed Brown Mat W; /* Workspace matrix. ATA */ 26c4762a1bSJed Brown Mat Hm; /* Hessian Misfit*/ 27c4762a1bSJed Brown Mat Hr; /* Hessian Reg*/ 28c4762a1bSJed Brown Vec d; /* RHS in least squares component $(1/2) * || F x - d ||_2^2$ */ 29c4762a1bSJed Brown Vec workLeft[NWORKLEFT]; /* Workspace for temporary vec */ 30c4762a1bSJed Brown Vec workRight[NWORKRIGHT]; /* Workspace for temporary vec */ 31c4762a1bSJed Brown NormType p; 32c4762a1bSJed Brown PetscRandom rctx; 33c4762a1bSJed Brown PetscBool taylor; /* Flag to determine whether to run Taylor test or not */ 34c4762a1bSJed Brown PetscBool use_admm; /* Flag to determine whether to run Taylor test or not */ 35c4762a1bSJed Brown } *UserCtx; 36c4762a1bSJed Brown 37d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateRHS(UserCtx ctx) 38d71ae5a4SJacob Faibussowitsch { 39c4762a1bSJed Brown PetscFunctionBegin; 40c4762a1bSJed Brown /* build the rhs d in ctx */ 419566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &(ctx->d))); 429566063dSJacob Faibussowitsch PetscCall(VecSetSizes(ctx->d, PETSC_DECIDE, ctx->m)); 439566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(ctx->d)); 449566063dSJacob Faibussowitsch PetscCall(VecSetRandom(ctx->d, ctx->rctx)); 453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 46c4762a1bSJed Brown } 47c4762a1bSJed Brown 48d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMatrix(UserCtx ctx) 49d71ae5a4SJacob Faibussowitsch { 50c4762a1bSJed Brown PetscInt Istart, Iend, i, j, Ii, gridN, I_n, I_s, I_e, I_w; 51c4762a1bSJed Brown #if defined(PETSC_USE_LOG) 52c4762a1bSJed Brown PetscLogStage stage; 53c4762a1bSJed Brown #endif 54c4762a1bSJed Brown 55c4762a1bSJed Brown PetscFunctionBegin; 56c4762a1bSJed Brown /* build the matrix F in ctx */ 579566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &(ctx->F))); 589566063dSJacob Faibussowitsch PetscCall(MatSetSizes(ctx->F, PETSC_DECIDE, PETSC_DECIDE, ctx->m, ctx->n)); 599566063dSJacob Faibussowitsch PetscCall(MatSetType(ctx->F, MATAIJ)); /* TODO: Decide specific SetType other than dummy*/ 609566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(ctx->F, 5, NULL, 5, NULL)); /*TODO: some number other than 5?*/ 619566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(ctx->F, 5, NULL)); 629566063dSJacob Faibussowitsch PetscCall(MatSetUp(ctx->F)); 639566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(ctx->F, &Istart, &Iend)); 649566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("Assembly", &stage)); 659566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage)); 66c4762a1bSJed Brown 673c859ba3SBarry Smith /* Set matrix elements in 2-D five point stencil format. */ 68c4762a1bSJed Brown if (!(ctx->matops)) { 693c859ba3SBarry Smith PetscCheck(ctx->m == ctx->n, PETSC_COMM_WORLD, PETSC_ERR_ARG_SIZ, "Stencil matrix must be square"); 70c4762a1bSJed Brown gridN = (PetscInt)PetscSqrtReal((PetscReal)ctx->m); 713c859ba3SBarry Smith PetscCheck(gridN * gridN == ctx->m, PETSC_COMM_WORLD, PETSC_ERR_ARG_SIZ, "Number of rows must be square"); 72c4762a1bSJed Brown for (Ii = Istart; Ii < Iend; Ii++) { 739371c9d4SSatish Balay i = Ii / gridN; 749371c9d4SSatish Balay 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 */ 9448a46eb9SPierre Jolivet if (!(ctx->matops)) PetscCall(MatSetOption(ctx->F, MAT_SYMMETRIC, PETSC_TRUE)); 959566063dSJacob Faibussowitsch PetscCall(MatTransposeMatMult(ctx->F, ctx->F, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &(ctx->W))); 96c4762a1bSJed Brown /* Setup Hessian Workspace in same shape as W */ 979566063dSJacob Faibussowitsch PetscCall(MatDuplicate(ctx->W, MAT_DO_NOT_COPY_VALUES, &(ctx->Hm))); 989566063dSJacob Faibussowitsch PetscCall(MatDuplicate(ctx->W, MAT_DO_NOT_COPY_VALUES, &(ctx->Hr))); 993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 100c4762a1bSJed Brown } 101c4762a1bSJed Brown 102d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupWorkspace(UserCtx ctx) 103d71ae5a4SJacob Faibussowitsch { 104c4762a1bSJed Brown PetscInt i; 105c4762a1bSJed Brown 106c4762a1bSJed Brown PetscFunctionBegin; 1079566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(ctx->F, &ctx->workLeft[0], &ctx->workRight[0])); 10848a46eb9SPierre Jolivet for (i = 1; i < NWORKLEFT; i++) PetscCall(VecDuplicate(ctx->workLeft[0], &(ctx->workLeft[i]))); 10948a46eb9SPierre Jolivet for (i = 1; i < NWORKRIGHT; i++) PetscCall(VecDuplicate(ctx->workRight[0], &(ctx->workRight[i]))); 1103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 111c4762a1bSJed Brown } 112c4762a1bSJed Brown 113d71ae5a4SJacob Faibussowitsch static PetscErrorCode ConfigureContext(UserCtx ctx) 114d71ae5a4SJacob Faibussowitsch { 115c4762a1bSJed Brown PetscFunctionBegin; 116c4762a1bSJed Brown ctx->m = 16; 117c4762a1bSJed Brown ctx->n = 16; 118c4762a1bSJed Brown ctx->eps = 1.e-3; 119c4762a1bSJed Brown ctx->abstol = 1.e-4; 120c4762a1bSJed Brown ctx->reltol = 1.e-2; 121c4762a1bSJed Brown ctx->hStart = 1.; 122c4762a1bSJed Brown ctx->hMin = 1.e-3; 123c4762a1bSJed Brown ctx->hFactor = 0.5; 124c4762a1bSJed Brown ctx->alpha = 1.; 125c4762a1bSJed Brown ctx->mu = 1.0; 126c4762a1bSJed Brown ctx->matops = 0; 127c4762a1bSJed Brown ctx->iter = 10; 128c4762a1bSJed Brown ctx->p = NORM_2; 129c4762a1bSJed Brown ctx->taylor = PETSC_TRUE; 130c4762a1bSJed Brown ctx->use_admm = PETSC_FALSE; 131d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Configure separable objection example", "ex4.c"); 1329566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-m", "The row dimension of matrix F", "ex4.c", ctx->m, &(ctx->m), NULL)); 1339566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-n", "The column dimension of matrix F", "ex4.c", ctx->n, &(ctx->n), NULL)); 1349566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-matrix_format", "Decide format of F matrix. 0 for stencil, 1 for random", "ex4.c", ctx->matops, &(ctx->matops), NULL)); 1359566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-iter", "Iteration number ADMM", "ex4.c", ctx->iter, &(ctx->iter), NULL)); 1369566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-alpha", "The regularization multiplier. 1 default", "ex4.c", ctx->alpha, &(ctx->alpha), NULL)); 1379566063dSJacob 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)); 1389566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-mu", "The augmented lagrangian multiplier in ADMM", "ex4.c", ctx->mu, &(ctx->mu), NULL)); 1399566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-hStart", "Taylor test starting point. 1 default.", "ex4.c", ctx->hStart, &(ctx->hStart), NULL)); 1409566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-hFactor", "Taylor test multiplier factor. 0.5 default", "ex4.c", ctx->hFactor, &(ctx->hFactor), NULL)); 1419566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-hMin", "Taylor test ending condition. 1.e-3 default", "ex4.c", ctx->hMin, &(ctx->hMin), NULL)); 1429566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-abstol", "Absolute stopping criterion for ADMM", "ex4.c", ctx->abstol, &(ctx->abstol), NULL)); 1439566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-reltol", "Relative stopping criterion for ADMM", "ex4.c", ctx->reltol, &(ctx->reltol), NULL)); 1449566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-taylor", "Flag for Taylor test. Default is true.", "ex4.c", ctx->taylor, &(ctx->taylor), NULL)); 1459566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-use_admm", "Use the ADMM solver in this example.", "ex4.c", ctx->use_admm, &(ctx->use_admm), NULL)); 1469566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-p", "Norm type.", "ex4.c", NormTypes, (PetscEnum)ctx->p, (PetscEnum *)&(ctx->p), NULL)); 147d0609cedSBarry Smith PetscOptionsEnd(); 148c4762a1bSJed Brown /* Creating random ctx */ 1499566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &(ctx->rctx))); 1509566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(ctx->rctx)); 1519566063dSJacob Faibussowitsch PetscCall(CreateMatrix(ctx)); 1529566063dSJacob Faibussowitsch PetscCall(CreateRHS(ctx)); 1539566063dSJacob Faibussowitsch PetscCall(SetupWorkspace(ctx)); 1543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 155c4762a1bSJed Brown } 156c4762a1bSJed Brown 157d71ae5a4SJacob Faibussowitsch static PetscErrorCode DestroyContext(UserCtx *ctx) 158d71ae5a4SJacob Faibussowitsch { 159c4762a1bSJed Brown PetscInt i; 160c4762a1bSJed Brown 161c4762a1bSJed Brown PetscFunctionBegin; 1629566063dSJacob Faibussowitsch PetscCall(MatDestroy(&((*ctx)->F))); 1639566063dSJacob Faibussowitsch PetscCall(MatDestroy(&((*ctx)->W))); 1649566063dSJacob Faibussowitsch PetscCall(MatDestroy(&((*ctx)->Hm))); 1659566063dSJacob Faibussowitsch PetscCall(MatDestroy(&((*ctx)->Hr))); 1669566063dSJacob Faibussowitsch PetscCall(VecDestroy(&((*ctx)->d))); 16748a46eb9SPierre Jolivet for (i = 0; i < NWORKLEFT; i++) PetscCall(VecDestroy(&((*ctx)->workLeft[i]))); 16848a46eb9SPierre Jolivet for (i = 0; i < NWORKRIGHT; i++) PetscCall(VecDestroy(&((*ctx)->workRight[i]))); 1699566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&((*ctx)->rctx))); 1709566063dSJacob Faibussowitsch PetscCall(PetscFree(*ctx)); 1713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 172c4762a1bSJed Brown } 173c4762a1bSJed Brown 174c4762a1bSJed Brown /* compute (1/2) * ||F x - d||^2 */ 175d71ae5a4SJacob Faibussowitsch static PetscErrorCode ObjectiveMisfit(Tao tao, Vec x, PetscReal *J, void *_ctx) 176d71ae5a4SJacob Faibussowitsch { 177c4762a1bSJed Brown UserCtx ctx = (UserCtx)_ctx; 178c4762a1bSJed Brown Vec y; 179c4762a1bSJed Brown 180c4762a1bSJed Brown PetscFunctionBegin; 181c4762a1bSJed Brown y = ctx->workLeft[0]; 1829566063dSJacob Faibussowitsch PetscCall(MatMult(ctx->F, x, y)); 1839566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1., ctx->d)); 1849566063dSJacob Faibussowitsch PetscCall(VecDot(y, y, J)); 185c4762a1bSJed Brown *J *= 0.5; 1863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 187c4762a1bSJed Brown } 188c4762a1bSJed Brown 189c4762a1bSJed Brown /* compute V = FTFx - FTd */ 190d71ae5a4SJacob Faibussowitsch static PetscErrorCode GradientMisfit(Tao tao, Vec x, Vec V, void *_ctx) 191d71ae5a4SJacob Faibussowitsch { 192c4762a1bSJed Brown UserCtx ctx = (UserCtx)_ctx; 193c4762a1bSJed Brown Vec FTFx, FTd; 194c4762a1bSJed Brown 195c4762a1bSJed Brown PetscFunctionBegin; 196c4762a1bSJed Brown /* work1 is A^T Ax, work2 is Ab, W is A^T A*/ 197c4762a1bSJed Brown FTFx = ctx->workRight[0]; 198c4762a1bSJed Brown FTd = ctx->workRight[1]; 1999566063dSJacob Faibussowitsch PetscCall(MatMult(ctx->W, x, FTFx)); 2009566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(ctx->F, ctx->d, FTd)); 2019566063dSJacob Faibussowitsch PetscCall(VecWAXPY(V, -1., FTd, FTFx)); 2023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 203c4762a1bSJed Brown } 204c4762a1bSJed Brown 205c4762a1bSJed Brown /* returns FTF */ 206d71ae5a4SJacob Faibussowitsch static PetscErrorCode HessianMisfit(Tao tao, Vec x, Mat H, Mat Hpre, void *_ctx) 207d71ae5a4SJacob Faibussowitsch { 208c4762a1bSJed Brown UserCtx ctx = (UserCtx)_ctx; 209c4762a1bSJed Brown 210c4762a1bSJed Brown PetscFunctionBegin; 2119566063dSJacob Faibussowitsch if (H != ctx->W) PetscCall(MatCopy(ctx->W, H, DIFFERENT_NONZERO_PATTERN)); 2129566063dSJacob Faibussowitsch if (Hpre != ctx->W) PetscCall(MatCopy(ctx->W, Hpre, DIFFERENT_NONZERO_PATTERN)); 2133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 214c4762a1bSJed Brown } 215c4762a1bSJed Brown 216c4762a1bSJed Brown /* computes augment Lagrangian objective (with scaled dual): 217c4762a1bSJed Brown * 0.5 * ||F x - d||^2 + 0.5 * mu ||x - z + u||^2 */ 218d71ae5a4SJacob Faibussowitsch static PetscErrorCode ObjectiveMisfitADMM(Tao tao, Vec x, PetscReal *J, void *_ctx) 219d71ae5a4SJacob Faibussowitsch { 220c4762a1bSJed Brown UserCtx ctx = (UserCtx)_ctx; 221c4762a1bSJed Brown PetscReal mu, workNorm, misfit; 222c4762a1bSJed Brown Vec z, u, temp; 223c4762a1bSJed Brown 224c4762a1bSJed Brown PetscFunctionBegin; 225c4762a1bSJed Brown mu = ctx->mu; 226c4762a1bSJed Brown z = ctx->workRight[5]; 227c4762a1bSJed Brown u = ctx->workRight[6]; 228c4762a1bSJed Brown temp = ctx->workRight[10]; 229c4762a1bSJed Brown /* misfit = f(x) */ 2309566063dSJacob Faibussowitsch PetscCall(ObjectiveMisfit(tao, x, &misfit, _ctx)); 2319566063dSJacob Faibussowitsch PetscCall(VecCopy(x, temp)); 232c4762a1bSJed Brown /* temp = x - z + u */ 2339566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(temp, -1., 1., 1., z, u)); 234c4762a1bSJed Brown /* workNorm = ||x - z + u||^2 */ 2359566063dSJacob Faibussowitsch PetscCall(VecDot(temp, temp, &workNorm)); 236c4762a1bSJed Brown /* augment Lagrangian objective (with scaled dual): f(x) + 0.5 * mu ||x -z + u||^2 */ 237c4762a1bSJed Brown *J = misfit + 0.5 * mu * workNorm; 2383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 239c4762a1bSJed Brown } 240c4762a1bSJed Brown 241c4762a1bSJed Brown /* computes FTFx - FTd mu*(x - z + u) */ 242d71ae5a4SJacob Faibussowitsch static PetscErrorCode GradientMisfitADMM(Tao tao, Vec x, Vec V, void *_ctx) 243d71ae5a4SJacob Faibussowitsch { 244c4762a1bSJed Brown UserCtx ctx = (UserCtx)_ctx; 245c4762a1bSJed Brown PetscReal mu; 246c4762a1bSJed Brown Vec z, u, temp; 247c4762a1bSJed Brown 248c4762a1bSJed Brown PetscFunctionBegin; 249c4762a1bSJed Brown mu = ctx->mu; 250c4762a1bSJed Brown z = ctx->workRight[5]; 251c4762a1bSJed Brown u = ctx->workRight[6]; 252c4762a1bSJed Brown temp = ctx->workRight[10]; 2539566063dSJacob Faibussowitsch PetscCall(GradientMisfit(tao, x, V, _ctx)); 2549566063dSJacob Faibussowitsch PetscCall(VecCopy(x, temp)); 255c4762a1bSJed Brown /* temp = x - z + u */ 2569566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(temp, -1., 1., 1., z, u)); 257c4762a1bSJed Brown /* V = FTFx - FTd mu*(x - z + u) */ 2589566063dSJacob Faibussowitsch PetscCall(VecAXPY(V, mu, temp)); 2593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 260c4762a1bSJed Brown } 261c4762a1bSJed Brown 262c4762a1bSJed Brown /* returns FTF + diag(mu) */ 263d71ae5a4SJacob Faibussowitsch static PetscErrorCode HessianMisfitADMM(Tao tao, Vec x, Mat H, Mat Hpre, void *_ctx) 264d71ae5a4SJacob Faibussowitsch { 265c4762a1bSJed Brown UserCtx ctx = (UserCtx)_ctx; 266c4762a1bSJed Brown 267c4762a1bSJed Brown PetscFunctionBegin; 2689566063dSJacob Faibussowitsch PetscCall(MatCopy(ctx->W, H, DIFFERENT_NONZERO_PATTERN)); 2699566063dSJacob Faibussowitsch PetscCall(MatShift(H, ctx->mu)); 27048a46eb9SPierre Jolivet if (Hpre != H) PetscCall(MatCopy(H, Hpre, DIFFERENT_NONZERO_PATTERN)); 2713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 272c4762a1bSJed Brown } 273c4762a1bSJed Brown 274c4762a1bSJed Brown /* computes || x ||_p (mult by 0.5 in case of NORM_2) */ 275d71ae5a4SJacob Faibussowitsch static PetscErrorCode ObjectiveRegularization(Tao tao, Vec x, PetscReal *J, void *_ctx) 276d71ae5a4SJacob Faibussowitsch { 277c4762a1bSJed Brown UserCtx ctx = (UserCtx)_ctx; 278c4762a1bSJed Brown PetscReal norm; 279c4762a1bSJed Brown 280c4762a1bSJed Brown PetscFunctionBegin; 281c4762a1bSJed Brown *J = 0; 2829566063dSJacob Faibussowitsch PetscCall(VecNorm(x, ctx->p, &norm)); 283c4762a1bSJed Brown if (ctx->p == NORM_2) norm = 0.5 * norm * norm; 284c4762a1bSJed Brown *J = ctx->alpha * norm; 2853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 286c4762a1bSJed Brown } 287c4762a1bSJed Brown 288c4762a1bSJed Brown /* NORM_2 Case: return x 289c4762a1bSJed Brown * NORM_1 Case: x/(|x| + eps) 290c4762a1bSJed Brown * Else: TODO */ 291d71ae5a4SJacob Faibussowitsch static PetscErrorCode GradientRegularization(Tao tao, Vec x, Vec V, void *_ctx) 292d71ae5a4SJacob Faibussowitsch { 293c4762a1bSJed Brown UserCtx ctx = (UserCtx)_ctx; 294c4762a1bSJed Brown PetscReal eps = ctx->eps; 295c4762a1bSJed Brown 296c4762a1bSJed Brown PetscFunctionBegin; 297c4762a1bSJed Brown if (ctx->p == NORM_2) { 2989566063dSJacob Faibussowitsch PetscCall(VecCopy(x, V)); 299c4762a1bSJed Brown } else if (ctx->p == NORM_1) { 3009566063dSJacob Faibussowitsch PetscCall(VecCopy(x, ctx->workRight[1])); 3019566063dSJacob Faibussowitsch PetscCall(VecAbs(ctx->workRight[1])); 3029566063dSJacob Faibussowitsch PetscCall(VecShift(ctx->workRight[1], eps)); 3039566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(V, x, ctx->workRight[1])); 304c4762a1bSJed Brown } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Example only works for NORM_1 and NORM_2"); 3053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 306c4762a1bSJed Brown } 307c4762a1bSJed Brown 308c4762a1bSJed Brown /* NORM_2 Case: returns diag(mu) 309c4762a1bSJed Brown * NORM_1 Case: diag(mu* 1/sqrt(x_i^2 + eps) * (1 - x_i^2/ABS(x_i^2+eps))) */ 310d71ae5a4SJacob Faibussowitsch static PetscErrorCode HessianRegularization(Tao tao, Vec x, Mat H, Mat Hpre, void *_ctx) 311d71ae5a4SJacob Faibussowitsch { 312c4762a1bSJed Brown UserCtx ctx = (UserCtx)_ctx; 313c4762a1bSJed Brown PetscReal eps = ctx->eps; 314c4762a1bSJed Brown Vec copy1, copy2, copy3; 315c4762a1bSJed Brown 316c4762a1bSJed Brown PetscFunctionBegin; 317c4762a1bSJed Brown if (ctx->p == NORM_2) { 318c4762a1bSJed Brown /* Identity matrix scaled by mu */ 3199566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(H)); 3209566063dSJacob Faibussowitsch PetscCall(MatShift(H, ctx->mu)); 321c4762a1bSJed Brown if (Hpre != H) { 3229566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(Hpre)); 3239566063dSJacob Faibussowitsch PetscCall(MatShift(Hpre, ctx->mu)); 324c4762a1bSJed Brown } 325c4762a1bSJed Brown } else if (ctx->p == NORM_1) { 326c4762a1bSJed Brown /* 1/sqrt(x_i^2 + eps) * (1 - x_i^2/ABS(x_i^2+eps)) */ 327c4762a1bSJed Brown copy1 = ctx->workRight[1]; 328c4762a1bSJed Brown copy2 = ctx->workRight[2]; 329c4762a1bSJed Brown copy3 = ctx->workRight[3]; 330c4762a1bSJed Brown /* copy1 : 1/sqrt(x_i^2 + eps) */ 3319566063dSJacob Faibussowitsch PetscCall(VecCopy(x, copy1)); 3329566063dSJacob Faibussowitsch PetscCall(VecPow(copy1, 2)); 3339566063dSJacob Faibussowitsch PetscCall(VecShift(copy1, eps)); 3349566063dSJacob Faibussowitsch PetscCall(VecSqrtAbs(copy1)); 3359566063dSJacob Faibussowitsch PetscCall(VecReciprocal(copy1)); 336c4762a1bSJed Brown /* copy2: x_i^2.*/ 3379566063dSJacob Faibussowitsch PetscCall(VecCopy(x, copy2)); 3389566063dSJacob Faibussowitsch PetscCall(VecPow(copy2, 2)); 339c4762a1bSJed Brown /* copy3: abs(x_i^2 + eps) */ 3409566063dSJacob Faibussowitsch PetscCall(VecCopy(x, copy3)); 3419566063dSJacob Faibussowitsch PetscCall(VecPow(copy3, 2)); 3429566063dSJacob Faibussowitsch PetscCall(VecShift(copy3, eps)); 3439566063dSJacob Faibussowitsch PetscCall(VecAbs(copy3)); 344c4762a1bSJed Brown /* copy2: 1 - x_i^2/abs(x_i^2 + eps) */ 3459566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(copy2, copy2, copy3)); 3469566063dSJacob Faibussowitsch PetscCall(VecScale(copy2, -1.)); 3479566063dSJacob Faibussowitsch PetscCall(VecShift(copy2, 1.)); 3489566063dSJacob Faibussowitsch PetscCall(VecAXPY(copy1, 1., copy2)); 3499566063dSJacob Faibussowitsch PetscCall(VecScale(copy1, ctx->mu)); 3509566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(H)); 3519566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(H, copy1, INSERT_VALUES)); 352c4762a1bSJed Brown if (Hpre != H) { 3539566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(Hpre)); 3549566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(Hpre, copy1, INSERT_VALUES)); 355c4762a1bSJed Brown } 356c4762a1bSJed Brown } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Example only works for NORM_1 and NORM_2"); 3573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 358c4762a1bSJed Brown } 359c4762a1bSJed Brown 360c4762a1bSJed Brown /* NORM_2 Case: 0.5 || x ||_2 + 0.5 * mu * ||x + u - z||^2 361c4762a1bSJed Brown * Else : || x ||_2 + 0.5 * mu * ||x + u - z||^2 */ 362d71ae5a4SJacob Faibussowitsch static PetscErrorCode ObjectiveRegularizationADMM(Tao tao, Vec z, PetscReal *J, void *_ctx) 363d71ae5a4SJacob Faibussowitsch { 364c4762a1bSJed Brown UserCtx ctx = (UserCtx)_ctx; 365c4762a1bSJed Brown PetscReal mu, workNorm, reg; 366c4762a1bSJed Brown Vec x, u, temp; 367c4762a1bSJed Brown 368c4762a1bSJed Brown PetscFunctionBegin; 369c4762a1bSJed Brown mu = ctx->mu; 370c4762a1bSJed Brown x = ctx->workRight[4]; 371c4762a1bSJed Brown u = ctx->workRight[6]; 372c4762a1bSJed Brown temp = ctx->workRight[10]; 3739566063dSJacob Faibussowitsch PetscCall(ObjectiveRegularization(tao, z, ®, _ctx)); 3749566063dSJacob Faibussowitsch PetscCall(VecCopy(z, temp)); 375c4762a1bSJed Brown /* temp = x + u -z */ 3769566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(temp, 1., 1., -1., x, u)); 377c4762a1bSJed Brown /* workNorm = ||x + u - z ||^2 */ 3789566063dSJacob Faibussowitsch PetscCall(VecDot(temp, temp, &workNorm)); 379c4762a1bSJed Brown *J = reg + 0.5 * mu * workNorm; 3803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 381c4762a1bSJed Brown } 382c4762a1bSJed Brown 383c4762a1bSJed Brown /* NORM_2 Case: x - mu*(x + u - z) 384c4762a1bSJed Brown * NORM_1 Case: x/(|x| + eps) - mu*(x + u - z) 385c4762a1bSJed Brown * Else: TODO */ 386d71ae5a4SJacob Faibussowitsch static PetscErrorCode GradientRegularizationADMM(Tao tao, Vec z, Vec V, void *_ctx) 387d71ae5a4SJacob Faibussowitsch { 388c4762a1bSJed Brown UserCtx ctx = (UserCtx)_ctx; 389c4762a1bSJed Brown PetscReal mu; 390c4762a1bSJed Brown Vec x, u, temp; 391c4762a1bSJed Brown 392c4762a1bSJed Brown PetscFunctionBegin; 393c4762a1bSJed Brown mu = ctx->mu; 394c4762a1bSJed Brown x = ctx->workRight[4]; 395c4762a1bSJed Brown u = ctx->workRight[6]; 396c4762a1bSJed Brown temp = ctx->workRight[10]; 3979566063dSJacob Faibussowitsch PetscCall(GradientRegularization(tao, z, V, _ctx)); 3989566063dSJacob Faibussowitsch PetscCall(VecCopy(z, temp)); 399c4762a1bSJed Brown /* temp = x + u -z */ 4009566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(temp, 1., 1., -1., x, u)); 4019566063dSJacob Faibussowitsch PetscCall(VecAXPY(V, -mu, temp)); 4023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 403c4762a1bSJed Brown } 404c4762a1bSJed Brown 405c4762a1bSJed Brown /* NORM_2 Case: returns diag(mu) 406c4762a1bSJed Brown * NORM_1 Case: FTF + diag(mu) */ 407d71ae5a4SJacob Faibussowitsch static PetscErrorCode HessianRegularizationADMM(Tao tao, Vec x, Mat H, Mat Hpre, void *_ctx) 408d71ae5a4SJacob Faibussowitsch { 409c4762a1bSJed Brown UserCtx ctx = (UserCtx)_ctx; 410c4762a1bSJed Brown 411c4762a1bSJed Brown PetscFunctionBegin; 412c4762a1bSJed Brown if (ctx->p == NORM_2) { 413c4762a1bSJed Brown /* Identity matrix scaled by mu */ 4149566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(H)); 4159566063dSJacob Faibussowitsch PetscCall(MatShift(H, ctx->mu)); 416c4762a1bSJed Brown if (Hpre != H) { 4179566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(Hpre)); 4189566063dSJacob Faibussowitsch PetscCall(MatShift(Hpre, ctx->mu)); 419c4762a1bSJed Brown } 420c4762a1bSJed Brown } else if (ctx->p == NORM_1) { 4219566063dSJacob Faibussowitsch PetscCall(HessianMisfit(tao, x, H, Hpre, (void *)ctx)); 4229566063dSJacob Faibussowitsch PetscCall(MatShift(H, ctx->mu)); 4239566063dSJacob Faibussowitsch if (Hpre != H) PetscCall(MatShift(Hpre, ctx->mu)); 424c4762a1bSJed Brown } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Example only works for NORM_1 and NORM_2"); 4253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 426c4762a1bSJed Brown } 427c4762a1bSJed Brown 428c4762a1bSJed Brown /* NORM_2 Case : (1/2) * ||F x - d||^2 + 0.5 * || x ||_p 429c4762a1bSJed Brown * NORM_1 Case : (1/2) * ||F x - d||^2 + || x ||_p */ 430d71ae5a4SJacob Faibussowitsch static PetscErrorCode ObjectiveComplete(Tao tao, Vec x, PetscReal *J, void *ctx) 431d71ae5a4SJacob Faibussowitsch { 432c4762a1bSJed Brown PetscReal Jm, Jr; 433c4762a1bSJed Brown 434c4762a1bSJed Brown PetscFunctionBegin; 4359566063dSJacob Faibussowitsch PetscCall(ObjectiveMisfit(tao, x, &Jm, ctx)); 4369566063dSJacob Faibussowitsch PetscCall(ObjectiveRegularization(tao, x, &Jr, ctx)); 437c4762a1bSJed Brown *J = Jm + Jr; 4383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 439c4762a1bSJed Brown } 440c4762a1bSJed Brown 441c4762a1bSJed Brown /* NORM_2 Case: FTFx - FTd + x 442c4762a1bSJed Brown * NORM_1 Case: FTFx - FTd + x/(|x| + eps) */ 443d71ae5a4SJacob Faibussowitsch static PetscErrorCode GradientComplete(Tao tao, Vec x, Vec V, void *ctx) 444d71ae5a4SJacob Faibussowitsch { 445c4762a1bSJed Brown UserCtx cntx = (UserCtx)ctx; 446c4762a1bSJed Brown 447c4762a1bSJed Brown PetscFunctionBegin; 4489566063dSJacob Faibussowitsch PetscCall(GradientMisfit(tao, x, cntx->workRight[2], ctx)); 4499566063dSJacob Faibussowitsch PetscCall(GradientRegularization(tao, x, cntx->workRight[3], ctx)); 4509566063dSJacob Faibussowitsch PetscCall(VecWAXPY(V, 1, cntx->workRight[2], cntx->workRight[3])); 4513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 452c4762a1bSJed Brown } 453c4762a1bSJed Brown 454c4762a1bSJed Brown /* NORM_2 Case: diag(mu) + FTF 455c4762a1bSJed Brown * NORM_1 Case: diag(mu* 1/sqrt(x_i^2 + eps) * (1 - x_i^2/ABS(x_i^2+eps))) + FTF */ 456d71ae5a4SJacob Faibussowitsch static PetscErrorCode HessianComplete(Tao tao, Vec x, Mat H, Mat Hpre, void *ctx) 457d71ae5a4SJacob Faibussowitsch { 458c4762a1bSJed Brown Mat tempH; 459c4762a1bSJed Brown 460c4762a1bSJed Brown PetscFunctionBegin; 4619566063dSJacob Faibussowitsch PetscCall(MatDuplicate(H, MAT_SHARE_NONZERO_PATTERN, &tempH)); 4629566063dSJacob Faibussowitsch PetscCall(HessianMisfit(tao, x, H, H, ctx)); 4639566063dSJacob Faibussowitsch PetscCall(HessianRegularization(tao, x, tempH, tempH, ctx)); 4649566063dSJacob Faibussowitsch PetscCall(MatAXPY(H, 1., tempH, DIFFERENT_NONZERO_PATTERN)); 46548a46eb9SPierre Jolivet if (Hpre != H) PetscCall(MatCopy(H, Hpre, DIFFERENT_NONZERO_PATTERN)); 4669566063dSJacob Faibussowitsch PetscCall(MatDestroy(&tempH)); 4673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 468c4762a1bSJed Brown } 469c4762a1bSJed Brown 470d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSolveADMM(UserCtx ctx, Vec x) 471d71ae5a4SJacob Faibussowitsch { 472c4762a1bSJed Brown PetscInt i; 473c4762a1bSJed Brown PetscReal u_norm, r_norm, s_norm, primal, dual, x_norm, z_norm; 474c4762a1bSJed Brown Tao tao1, tao2; 475c4762a1bSJed Brown Vec xk, z, u, diff, zold, zdiff, temp; 476c4762a1bSJed Brown PetscReal mu; 477c4762a1bSJed Brown 478c4762a1bSJed Brown PetscFunctionBegin; 479c4762a1bSJed Brown xk = ctx->workRight[4]; 480c4762a1bSJed Brown z = ctx->workRight[5]; 481c4762a1bSJed Brown u = ctx->workRight[6]; 482c4762a1bSJed Brown diff = ctx->workRight[7]; 483c4762a1bSJed Brown zold = ctx->workRight[8]; 484c4762a1bSJed Brown zdiff = ctx->workRight[9]; 485c4762a1bSJed Brown temp = ctx->workRight[11]; 486c4762a1bSJed Brown mu = ctx->mu; 4879566063dSJacob Faibussowitsch PetscCall(VecSet(u, 0.)); 4889566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao1)); 4899566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao1, TAONLS)); 4909566063dSJacob Faibussowitsch PetscCall(TaoSetObjective(tao1, ObjectiveMisfitADMM, (void *)ctx)); 4919566063dSJacob Faibussowitsch PetscCall(TaoSetGradient(tao1, NULL, GradientMisfitADMM, (void *)ctx)); 4929566063dSJacob Faibussowitsch PetscCall(TaoSetHessian(tao1, ctx->Hm, ctx->Hm, HessianMisfitADMM, (void *)ctx)); 4939566063dSJacob Faibussowitsch PetscCall(VecSet(xk, 0.)); 4949566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao1, xk)); 4959566063dSJacob Faibussowitsch PetscCall(TaoSetOptionsPrefix(tao1, "misfit_")); 4969566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao1)); 4979566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao2)); 498c4762a1bSJed Brown if (ctx->p == NORM_2) { 4999566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao2, TAONLS)); 5009566063dSJacob Faibussowitsch PetscCall(TaoSetObjective(tao2, ObjectiveRegularizationADMM, (void *)ctx)); 5019566063dSJacob Faibussowitsch PetscCall(TaoSetGradient(tao2, NULL, GradientRegularizationADMM, (void *)ctx)); 5029566063dSJacob Faibussowitsch PetscCall(TaoSetHessian(tao2, ctx->Hr, ctx->Hr, HessianRegularizationADMM, (void *)ctx)); 503c4762a1bSJed Brown } 5049566063dSJacob Faibussowitsch PetscCall(VecSet(z, 0.)); 5059566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao2, z)); 5069566063dSJacob Faibussowitsch PetscCall(TaoSetOptionsPrefix(tao2, "reg_")); 5079566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao2)); 508c4762a1bSJed Brown 509c4762a1bSJed Brown for (i = 0; i < ctx->iter; i++) { 5109566063dSJacob Faibussowitsch PetscCall(VecCopy(z, zold)); 5119566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao1)); /* Updates xk */ 512c4762a1bSJed Brown if (ctx->p == NORM_1) { 5139566063dSJacob Faibussowitsch PetscCall(VecWAXPY(temp, 1., xk, u)); 5149566063dSJacob Faibussowitsch PetscCall(TaoSoftThreshold(temp, -ctx->alpha / mu, ctx->alpha / mu, z)); 515c4762a1bSJed Brown } else { 5169566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao2)); /* Update zk */ 517c4762a1bSJed Brown } 518c4762a1bSJed Brown /* u = u + xk -z */ 5199566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(u, 1., -1., 1., xk, z)); 520c4762a1bSJed Brown /* r_norm : norm(x-z) */ 5219566063dSJacob Faibussowitsch PetscCall(VecWAXPY(diff, -1., z, xk)); 5229566063dSJacob Faibussowitsch PetscCall(VecNorm(diff, NORM_2, &r_norm)); 523c4762a1bSJed Brown /* s_norm : norm(-mu(z-zold)) */ 5249566063dSJacob Faibussowitsch PetscCall(VecWAXPY(zdiff, -1., zold, z)); 5259566063dSJacob Faibussowitsch PetscCall(VecNorm(zdiff, NORM_2, &s_norm)); 526c4762a1bSJed Brown s_norm = s_norm * mu; 527c4762a1bSJed Brown /* primal : sqrt(n)*ABSTOL + RELTOL*max(norm(x), norm(-z))*/ 5289566063dSJacob Faibussowitsch PetscCall(VecNorm(xk, NORM_2, &x_norm)); 5299566063dSJacob Faibussowitsch PetscCall(VecNorm(z, NORM_2, &z_norm)); 530c4762a1bSJed Brown primal = PetscSqrtReal(ctx->n) * ctx->abstol + ctx->reltol * PetscMax(x_norm, z_norm); 531c4762a1bSJed Brown /* Duality : sqrt(n)*ABSTOL + RELTOL*norm(mu*u)*/ 5329566063dSJacob Faibussowitsch PetscCall(VecNorm(u, NORM_2, &u_norm)); 533c4762a1bSJed Brown dual = PetscSqrtReal(ctx->n) * ctx->abstol + ctx->reltol * u_norm * mu; 53463a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PetscObjectComm((PetscObject)tao1), "Iter %" PetscInt_FMT " : ||x-z||: %g, mu*||z-zold||: %g\n", i, (double)r_norm, (double)s_norm)); 535c4762a1bSJed Brown if (r_norm < primal && s_norm < dual) break; 536c4762a1bSJed Brown } 5379566063dSJacob Faibussowitsch PetscCall(VecCopy(xk, x)); 5389566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao1)); 5399566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao2)); 5403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 541c4762a1bSJed Brown } 542c4762a1bSJed Brown 543c4762a1bSJed Brown /* Second order Taylor remainder convergence test */ 544d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaylorTest(UserCtx ctx, Tao tao, Vec x, PetscReal *C) 545d71ae5a4SJacob Faibussowitsch { 546c4762a1bSJed Brown PetscReal h, J, temp; 547c4762a1bSJed Brown PetscInt i, j; 548c4762a1bSJed Brown PetscInt numValues; 549c4762a1bSJed Brown PetscReal Jx, Jxhat_comp, Jxhat_pred; 550c4762a1bSJed Brown PetscReal *Js, *hs; 551c4762a1bSJed Brown PetscReal gdotdx; 552c4762a1bSJed Brown PetscReal minrate = PETSC_MAX_REAL; 553c4762a1bSJed Brown MPI_Comm comm = PetscObjectComm((PetscObject)x); 554c4762a1bSJed Brown Vec g, dx, xhat; 555c4762a1bSJed Brown 556c4762a1bSJed Brown PetscFunctionBegin; 5579566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &g)); 5589566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &xhat)); 559c4762a1bSJed Brown /* choose a perturbation direction */ 5609566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &dx)); 5619566063dSJacob Faibussowitsch PetscCall(VecSetRandom(dx, ctx->rctx)); 562c4762a1bSJed Brown /* evaluate objective at x: J(x) */ 5639566063dSJacob Faibussowitsch PetscCall(TaoComputeObjective(tao, x, &Jx)); 564c4762a1bSJed Brown /* evaluate gradient at x, save in vector g */ 5659566063dSJacob Faibussowitsch PetscCall(TaoComputeGradient(tao, x, g)); 5669566063dSJacob Faibussowitsch PetscCall(VecDot(g, dx, &gdotdx)); 567c4762a1bSJed Brown 568c4762a1bSJed Brown for (numValues = 0, h = ctx->hStart; h >= ctx->hMin; h *= ctx->hFactor) numValues++; 5699566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(numValues, &Js, numValues, &hs)); 570c4762a1bSJed Brown for (i = 0, h = ctx->hStart; h >= ctx->hMin; h *= ctx->hFactor, i++) { 5719566063dSJacob Faibussowitsch PetscCall(VecWAXPY(xhat, h, dx, x)); 5729566063dSJacob Faibussowitsch PetscCall(TaoComputeObjective(tao, xhat, &Jxhat_comp)); 573c4762a1bSJed Brown /* J(\hat(x)) \approx J(x) + g^T (xhat - x) = J(x) + h * g^T dx */ 574c4762a1bSJed Brown Jxhat_pred = Jx + h * gdotdx; 575c4762a1bSJed Brown /* Vector to dJdm scalar? Dot?*/ 576c4762a1bSJed Brown J = PetscAbsReal(Jxhat_comp - Jxhat_pred); 5779566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "J(xhat): %g, predicted: %g, diff %g\n", (double)Jxhat_comp, (double)Jxhat_pred, (double)J)); 578c4762a1bSJed Brown Js[i] = J; 579c4762a1bSJed Brown hs[i] = h; 580c4762a1bSJed Brown } 581c4762a1bSJed Brown for (j = 1; j < numValues; j++) { 582c4762a1bSJed Brown temp = PetscLogReal(Js[j] / Js[j - 1]) / PetscLogReal(hs[j] / hs[j - 1]); 58363a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Convergence rate step %" PetscInt_FMT ": %g\n", j - 1, (double)temp)); 584c4762a1bSJed Brown minrate = PetscMin(minrate, temp); 585c4762a1bSJed Brown } 586c4762a1bSJed Brown /* If O is not ~2, then the test is wrong */ 5879566063dSJacob Faibussowitsch PetscCall(PetscFree2(Js, hs)); 588c4762a1bSJed Brown *C = minrate; 5899566063dSJacob Faibussowitsch PetscCall(VecDestroy(&dx)); 5909566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xhat)); 5919566063dSJacob Faibussowitsch PetscCall(VecDestroy(&g)); 5923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 593c4762a1bSJed Brown } 594c4762a1bSJed Brown 595d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 596d71ae5a4SJacob Faibussowitsch { 597c4762a1bSJed Brown UserCtx ctx; 598c4762a1bSJed Brown Tao tao; 599c4762a1bSJed Brown Vec x; 600c4762a1bSJed Brown Mat H; 601c4762a1bSJed Brown 602327415f7SBarry Smith PetscFunctionBeginUser; 6039566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 6049566063dSJacob Faibussowitsch PetscCall(PetscNew(&ctx)); 6059566063dSJacob Faibussowitsch PetscCall(ConfigureContext(ctx)); 606a82e8c82SStefano Zampini /* Define two functions that could pass as objectives to TaoSetObjective(): one 607c4762a1bSJed Brown * for the misfit component, and one for the regularization component */ 608c4762a1bSJed Brown /* ObjectiveMisfit() and ObjectiveRegularization() */ 609c4762a1bSJed Brown 610c4762a1bSJed Brown /* Define a single function that calls both components adds them together: the complete objective, 611c4762a1bSJed Brown * in the absence of a Tao implementation that handles separability */ 612c4762a1bSJed Brown /* ObjectiveComplete() */ 6139566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao)); 6149566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao, TAONM)); 6159566063dSJacob Faibussowitsch PetscCall(TaoSetObjective(tao, ObjectiveComplete, (void *)ctx)); 6169566063dSJacob Faibussowitsch PetscCall(TaoSetGradient(tao, NULL, GradientComplete, (void *)ctx)); 6179566063dSJacob Faibussowitsch PetscCall(MatDuplicate(ctx->W, MAT_SHARE_NONZERO_PATTERN, &H)); 6189566063dSJacob Faibussowitsch PetscCall(TaoSetHessian(tao, H, H, HessianComplete, (void *)ctx)); 6199566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(ctx->F, NULL, &x)); 6209566063dSJacob Faibussowitsch PetscCall(VecSet(x, 0.)); 6219566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, x)); 6229566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao)); 6231baa6e33SBarry Smith if (ctx->use_admm) PetscCall(TaoSolveADMM(ctx, x)); 6241baa6e33SBarry Smith else PetscCall(TaoSolve(tao)); 625c4762a1bSJed Brown /* examine solution */ 6269566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(x, NULL, "-view_sol")); 627c4762a1bSJed Brown if (ctx->taylor) { 628c4762a1bSJed Brown PetscReal rate; 6299566063dSJacob Faibussowitsch PetscCall(TaylorTest(ctx, tao, x, &rate)); 630c4762a1bSJed Brown } 6319566063dSJacob Faibussowitsch PetscCall(MatDestroy(&H)); 6329566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao)); 6339566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 6349566063dSJacob Faibussowitsch PetscCall(DestroyContext(&ctx)); 6359566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 636b122ec5aSJacob Faibussowitsch return 0; 637c4762a1bSJed Brown } 638c4762a1bSJed Brown 639c4762a1bSJed Brown /*TEST 640c4762a1bSJed Brown 641c4762a1bSJed Brown build: 642c4762a1bSJed Brown requires: !complex 643c4762a1bSJed Brown 644c4762a1bSJed Brown test: 645c4762a1bSJed Brown suffix: 0 646c4762a1bSJed Brown args: 647c4762a1bSJed Brown 648c4762a1bSJed Brown test: 649c4762a1bSJed Brown suffix: l1_1 650c4762a1bSJed Brown args: -p 1 -tao_type lmvm -alpha 1. -epsilon 1.e-7 -m 64 -n 64 -view_sol -matrix_format 1 651c4762a1bSJed Brown 652c4762a1bSJed Brown test: 653c4762a1bSJed Brown suffix: hessian_1 654c5f5e425SStefano Zampini args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 1 -tao_type nls 655c4762a1bSJed Brown 656c4762a1bSJed Brown test: 657c4762a1bSJed Brown suffix: hessian_2 658c5f5e425SStefano Zampini args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 2 -tao_type nls 659c4762a1bSJed Brown 660c4762a1bSJed Brown test: 661c4762a1bSJed Brown suffix: nm_1 662c4762a1bSJed Brown args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 1 -tao_type nm -tao_max_it 50 663c4762a1bSJed Brown 664c4762a1bSJed Brown test: 665c4762a1bSJed Brown suffix: nm_2 666c4762a1bSJed Brown args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 2 -tao_type nm -tao_max_it 50 667c4762a1bSJed Brown 668c4762a1bSJed Brown test: 669c4762a1bSJed Brown suffix: lmvm_1 670c4762a1bSJed Brown args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 1 -tao_type lmvm -tao_max_it 40 671c4762a1bSJed Brown 672c4762a1bSJed Brown test: 673c4762a1bSJed Brown suffix: lmvm_2 674c4762a1bSJed Brown args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 2 -tao_type lmvm -tao_max_it 15 675c4762a1bSJed Brown 676c4762a1bSJed Brown test: 677c4762a1bSJed Brown suffix: soft_threshold_admm_1 678c4762a1bSJed Brown args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 1 -use_admm 679c4762a1bSJed Brown 680c4762a1bSJed Brown test: 681c4762a1bSJed Brown suffix: hessian_admm_1 682c4762a1bSJed Brown args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 1 -use_admm -reg_tao_type nls -misfit_tao_type nls 683c4762a1bSJed Brown 684c4762a1bSJed Brown test: 685c4762a1bSJed Brown suffix: hessian_admm_2 686c4762a1bSJed Brown args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 2 -use_admm -reg_tao_type nls -misfit_tao_type nls 687c4762a1bSJed Brown 688c4762a1bSJed Brown test: 689c4762a1bSJed Brown suffix: nm_admm_1 690c4762a1bSJed Brown args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 1 -use_admm -reg_tao_type nm -misfit_tao_type nm 691c4762a1bSJed Brown 692c4762a1bSJed Brown test: 693c4762a1bSJed Brown suffix: nm_admm_2 694c4762a1bSJed 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 695c4762a1bSJed Brown 696c4762a1bSJed Brown test: 697c4762a1bSJed Brown suffix: lmvm_admm_1 698c4762a1bSJed Brown args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 1 -use_admm -reg_tao_type lmvm -misfit_tao_type lmvm 699c4762a1bSJed Brown 700c4762a1bSJed Brown test: 701c4762a1bSJed Brown suffix: lmvm_admm_2 702c4762a1bSJed Brown args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 2 -use_admm -reg_tao_type lmvm -misfit_tao_type lmvm 703c4762a1bSJed Brown 704c4762a1bSJed Brown TEST*/ 705