Actual source code: ex4.c

  1: static char help[] = "Simple example to test separable objective optimizers.\n";

  3: #include <petsc.h>
  4: #include <petsctao.h>
  5: #include <petscvec.h>
  6: #include <petscmath.h>

  8: #define NWORKLEFT  4
  9: #define NWORKRIGHT 12

 11: typedef struct _UserCtx {
 12:   PetscInt    m;       /* The row dimension of F */
 13:   PetscInt    n;       /* The column dimension of F */
 14:   PetscInt    matops;  /* Matrix format. 0 for stencil, 1 for random */
 15:   PetscInt    iter;    /* Number of iterations for ADMM */
 16:   PetscReal   hStart;  /* Starting point for Taylor test */
 17:   PetscReal   hFactor; /* Taylor test step factor */
 18:   PetscReal   hMin;    /* Taylor test end goal */
 19:   PetscReal   alpha;   /* regularization constant applied to || x ||_p */
 20:   PetscReal   eps;     /* small constant for approximating gradient of || x ||_1 */
 21:   PetscReal   mu;      /* the augmented Lagrangian term in ADMM */
 22:   PetscReal   abstol;
 23:   PetscReal   reltol;
 24:   Mat         F;                     /* matrix in least squares component $(1/2) * || F x - d ||_2^2$ */
 25:   Mat         W;                     /* Workspace matrix. ATA */
 26:   Mat         Hm;                    /* Hessian Misfit*/
 27:   Mat         Hr;                    /* Hessian Reg*/
 28:   Vec         d;                     /* RHS in least squares component $(1/2) * || F x - d ||_2^2$ */
 29:   Vec         workLeft[NWORKLEFT];   /* Workspace for temporary vec */
 30:   Vec         workRight[NWORKRIGHT]; /* Workspace for temporary vec */
 31:   NormType    p;
 32:   PetscRandom rctx;
 33:   PetscBool   taylor;   /* Flag to determine whether to run Taylor test or not */
 34:   PetscBool   use_admm; /* Flag to determine whether to run Taylor test or not */
 35: } *UserCtx;

 37: static PetscErrorCode CreateRHS(UserCtx ctx)
 38: {
 39:   PetscFunctionBegin;
 40:   /* build the rhs d in ctx */
 41:   PetscCall(VecCreate(PETSC_COMM_WORLD, &(ctx->d)));
 42:   PetscCall(VecSetSizes(ctx->d, PETSC_DECIDE, ctx->m));
 43:   PetscCall(VecSetFromOptions(ctx->d));
 44:   PetscCall(VecSetRandom(ctx->d, ctx->rctx));
 45:   PetscFunctionReturn(PETSC_SUCCESS);
 46: }

 48: static PetscErrorCode CreateMatrix(UserCtx ctx)
 49: {
 50:   PetscInt      Istart, Iend, i, j, Ii, gridN, I_n, I_s, I_e, I_w;
 51:   PetscLogStage stage;

 53:   PetscFunctionBegin;
 54:   /* build the matrix F in ctx */
 55:   PetscCall(MatCreate(PETSC_COMM_WORLD, &(ctx->F)));
 56:   PetscCall(MatSetSizes(ctx->F, PETSC_DECIDE, PETSC_DECIDE, ctx->m, ctx->n));
 57:   PetscCall(MatSetType(ctx->F, MATAIJ));                          /* TODO: Decide specific SetType other than dummy*/
 58:   PetscCall(MatMPIAIJSetPreallocation(ctx->F, 5, NULL, 5, NULL)); /*TODO: some number other than 5?*/
 59:   PetscCall(MatSeqAIJSetPreallocation(ctx->F, 5, NULL));
 60:   PetscCall(MatSetUp(ctx->F));
 61:   PetscCall(MatGetOwnershipRange(ctx->F, &Istart, &Iend));
 62:   PetscCall(PetscLogStageRegister("Assembly", &stage));
 63:   PetscCall(PetscLogStagePush(stage));

 65:   /* Set matrix elements in  2-D five point stencil format. */
 66:   if (!(ctx->matops)) {
 67:     PetscCheck(ctx->m == ctx->n, PETSC_COMM_WORLD, PETSC_ERR_ARG_SIZ, "Stencil matrix must be square");
 68:     gridN = (PetscInt)PetscSqrtReal((PetscReal)ctx->m);
 69:     PetscCheck(gridN * gridN == ctx->m, PETSC_COMM_WORLD, PETSC_ERR_ARG_SIZ, "Number of rows must be square");
 70:     for (Ii = Istart; Ii < Iend; Ii++) {
 71:       i   = Ii / gridN;
 72:       j   = Ii % gridN;
 73:       I_n = i * gridN + j + 1;
 74:       if (j + 1 >= gridN) I_n = -1;
 75:       I_s = i * gridN + j - 1;
 76:       if (j - 1 < 0) I_s = -1;
 77:       I_e = (i + 1) * gridN + j;
 78:       if (i + 1 >= gridN) I_e = -1;
 79:       I_w = (i - 1) * gridN + j;
 80:       if (i - 1 < 0) I_w = -1;
 81:       PetscCall(MatSetValue(ctx->F, Ii, Ii, 4., INSERT_VALUES));
 82:       PetscCall(MatSetValue(ctx->F, Ii, I_n, -1., INSERT_VALUES));
 83:       PetscCall(MatSetValue(ctx->F, Ii, I_s, -1., INSERT_VALUES));
 84:       PetscCall(MatSetValue(ctx->F, Ii, I_e, -1., INSERT_VALUES));
 85:       PetscCall(MatSetValue(ctx->F, Ii, I_w, -1., INSERT_VALUES));
 86:     }
 87:   } else PetscCall(MatSetRandom(ctx->F, ctx->rctx));
 88:   PetscCall(MatAssemblyBegin(ctx->F, MAT_FINAL_ASSEMBLY));
 89:   PetscCall(MatAssemblyEnd(ctx->F, MAT_FINAL_ASSEMBLY));
 90:   PetscCall(PetscLogStagePop());
 91:   /* Stencil matrix is symmetric. Setting symmetric flag for ICC/Cholesky preconditioner */
 92:   if (!(ctx->matops)) PetscCall(MatSetOption(ctx->F, MAT_SYMMETRIC, PETSC_TRUE));
 93:   PetscCall(MatTransposeMatMult(ctx->F, ctx->F, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &(ctx->W)));
 94:   /* Setup Hessian Workspace in same shape as W */
 95:   PetscCall(MatDuplicate(ctx->W, MAT_DO_NOT_COPY_VALUES, &(ctx->Hm)));
 96:   PetscCall(MatDuplicate(ctx->W, MAT_DO_NOT_COPY_VALUES, &(ctx->Hr)));
 97:   PetscFunctionReturn(PETSC_SUCCESS);
 98: }

100: static PetscErrorCode SetupWorkspace(UserCtx ctx)
101: {
102:   PetscInt i;

104:   PetscFunctionBegin;
105:   PetscCall(MatCreateVecs(ctx->F, &ctx->workLeft[0], &ctx->workRight[0]));
106:   for (i = 1; i < NWORKLEFT; i++) PetscCall(VecDuplicate(ctx->workLeft[0], &(ctx->workLeft[i])));
107:   for (i = 1; i < NWORKRIGHT; i++) PetscCall(VecDuplicate(ctx->workRight[0], &(ctx->workRight[i])));
108:   PetscFunctionReturn(PETSC_SUCCESS);
109: }

111: static PetscErrorCode ConfigureContext(UserCtx ctx)
112: {
113:   PetscFunctionBegin;
114:   ctx->m        = 16;
115:   ctx->n        = 16;
116:   ctx->eps      = 1.e-3;
117:   ctx->abstol   = 1.e-4;
118:   ctx->reltol   = 1.e-2;
119:   ctx->hStart   = 1.;
120:   ctx->hMin     = 1.e-3;
121:   ctx->hFactor  = 0.5;
122:   ctx->alpha    = 1.;
123:   ctx->mu       = 1.0;
124:   ctx->matops   = 0;
125:   ctx->iter     = 10;
126:   ctx->p        = NORM_2;
127:   ctx->taylor   = PETSC_TRUE;
128:   ctx->use_admm = PETSC_FALSE;
129:   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Configure separable objection example", "ex4.c");
130:   PetscCall(PetscOptionsInt("-m", "The row dimension of matrix F", "ex4.c", ctx->m, &(ctx->m), NULL));
131:   PetscCall(PetscOptionsInt("-n", "The column dimension of matrix F", "ex4.c", ctx->n, &(ctx->n), NULL));
132:   PetscCall(PetscOptionsInt("-matrix_format", "Decide format of F matrix. 0 for stencil, 1 for random", "ex4.c", ctx->matops, &(ctx->matops), NULL));
133:   PetscCall(PetscOptionsInt("-iter", "Iteration number ADMM", "ex4.c", ctx->iter, &(ctx->iter), NULL));
134:   PetscCall(PetscOptionsReal("-alpha", "The regularization multiplier. 1 default", "ex4.c", ctx->alpha, &(ctx->alpha), NULL));
135:   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));
136:   PetscCall(PetscOptionsReal("-mu", "The augmented lagrangian multiplier in ADMM", "ex4.c", ctx->mu, &(ctx->mu), NULL));
137:   PetscCall(PetscOptionsReal("-hStart", "Taylor test starting point. 1 default.", "ex4.c", ctx->hStart, &(ctx->hStart), NULL));
138:   PetscCall(PetscOptionsReal("-hFactor", "Taylor test multiplier factor. 0.5 default", "ex4.c", ctx->hFactor, &(ctx->hFactor), NULL));
139:   PetscCall(PetscOptionsReal("-hMin", "Taylor test ending condition. 1.e-3 default", "ex4.c", ctx->hMin, &(ctx->hMin), NULL));
140:   PetscCall(PetscOptionsReal("-abstol", "Absolute stopping criterion for ADMM", "ex4.c", ctx->abstol, &(ctx->abstol), NULL));
141:   PetscCall(PetscOptionsReal("-reltol", "Relative stopping criterion for ADMM", "ex4.c", ctx->reltol, &(ctx->reltol), NULL));
142:   PetscCall(PetscOptionsBool("-taylor", "Flag for Taylor test. Default is true.", "ex4.c", ctx->taylor, &(ctx->taylor), NULL));
143:   PetscCall(PetscOptionsBool("-use_admm", "Use the ADMM solver in this example.", "ex4.c", ctx->use_admm, &(ctx->use_admm), NULL));
144:   PetscCall(PetscOptionsEnum("-p", "Norm type.", "ex4.c", NormTypes, (PetscEnum)ctx->p, (PetscEnum *)&(ctx->p), NULL));
145:   PetscOptionsEnd();
146:   /* Creating random ctx */
147:   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &(ctx->rctx)));
148:   PetscCall(PetscRandomSetFromOptions(ctx->rctx));
149:   PetscCall(CreateMatrix(ctx));
150:   PetscCall(CreateRHS(ctx));
151:   PetscCall(SetupWorkspace(ctx));
152:   PetscFunctionReturn(PETSC_SUCCESS);
153: }

155: static PetscErrorCode DestroyContext(UserCtx *ctx)
156: {
157:   PetscInt i;

159:   PetscFunctionBegin;
160:   PetscCall(MatDestroy(&((*ctx)->F)));
161:   PetscCall(MatDestroy(&((*ctx)->W)));
162:   PetscCall(MatDestroy(&((*ctx)->Hm)));
163:   PetscCall(MatDestroy(&((*ctx)->Hr)));
164:   PetscCall(VecDestroy(&((*ctx)->d)));
165:   for (i = 0; i < NWORKLEFT; i++) PetscCall(VecDestroy(&((*ctx)->workLeft[i])));
166:   for (i = 0; i < NWORKRIGHT; i++) PetscCall(VecDestroy(&((*ctx)->workRight[i])));
167:   PetscCall(PetscRandomDestroy(&((*ctx)->rctx)));
168:   PetscCall(PetscFree(*ctx));
169:   PetscFunctionReturn(PETSC_SUCCESS);
170: }

172: /* compute (1/2) * ||F x - d||^2 */
173: static PetscErrorCode ObjectiveMisfit(Tao tao, Vec x, PetscReal *J, void *_ctx)
174: {
175:   UserCtx ctx = (UserCtx)_ctx;
176:   Vec     y;

178:   PetscFunctionBegin;
179:   y = ctx->workLeft[0];
180:   PetscCall(MatMult(ctx->F, x, y));
181:   PetscCall(VecAXPY(y, -1., ctx->d));
182:   PetscCall(VecDot(y, y, J));
183:   *J *= 0.5;
184:   PetscFunctionReturn(PETSC_SUCCESS);
185: }

187: /* compute V = FTFx - FTd */
188: static PetscErrorCode GradientMisfit(Tao tao, Vec x, Vec V, void *_ctx)
189: {
190:   UserCtx ctx = (UserCtx)_ctx;
191:   Vec     FTFx, FTd;

193:   PetscFunctionBegin;
194:   /* work1 is A^T Ax, work2 is Ab, W is A^T A*/
195:   FTFx = ctx->workRight[0];
196:   FTd  = ctx->workRight[1];
197:   PetscCall(MatMult(ctx->W, x, FTFx));
198:   PetscCall(MatMultTranspose(ctx->F, ctx->d, FTd));
199:   PetscCall(VecWAXPY(V, -1., FTd, FTFx));
200:   PetscFunctionReturn(PETSC_SUCCESS);
201: }

203: /* returns FTF */
204: static PetscErrorCode HessianMisfit(Tao tao, Vec x, Mat H, Mat Hpre, void *_ctx)
205: {
206:   UserCtx ctx = (UserCtx)_ctx;

208:   PetscFunctionBegin;
209:   if (H != ctx->W) PetscCall(MatCopy(ctx->W, H, DIFFERENT_NONZERO_PATTERN));
210:   if (Hpre != ctx->W) PetscCall(MatCopy(ctx->W, Hpre, DIFFERENT_NONZERO_PATTERN));
211:   PetscFunctionReturn(PETSC_SUCCESS);
212: }

214: /* computes augment Lagrangian objective (with scaled dual):
215:  * 0.5 * ||F x - d||^2  + 0.5 * mu ||x - z + u||^2 */
216: static PetscErrorCode ObjectiveMisfitADMM(Tao tao, Vec x, PetscReal *J, void *_ctx)
217: {
218:   UserCtx   ctx = (UserCtx)_ctx;
219:   PetscReal mu, workNorm, misfit;
220:   Vec       z, u, temp;

222:   PetscFunctionBegin;
223:   mu   = ctx->mu;
224:   z    = ctx->workRight[5];
225:   u    = ctx->workRight[6];
226:   temp = ctx->workRight[10];
227:   /* misfit = f(x) */
228:   PetscCall(ObjectiveMisfit(tao, x, &misfit, _ctx));
229:   PetscCall(VecCopy(x, temp));
230:   /* temp = x - z + u */
231:   PetscCall(VecAXPBYPCZ(temp, -1., 1., 1., z, u));
232:   /* workNorm = ||x - z + u||^2 */
233:   PetscCall(VecDot(temp, temp, &workNorm));
234:   /* augment Lagrangian objective (with scaled dual): f(x) + 0.5 * mu ||x -z + u||^2 */
235:   *J = misfit + 0.5 * mu * workNorm;
236:   PetscFunctionReturn(PETSC_SUCCESS);
237: }

239: /* computes FTFx - FTd  mu*(x - z + u) */
240: static PetscErrorCode GradientMisfitADMM(Tao tao, Vec x, Vec V, void *_ctx)
241: {
242:   UserCtx   ctx = (UserCtx)_ctx;
243:   PetscReal mu;
244:   Vec       z, u, temp;

246:   PetscFunctionBegin;
247:   mu   = ctx->mu;
248:   z    = ctx->workRight[5];
249:   u    = ctx->workRight[6];
250:   temp = ctx->workRight[10];
251:   PetscCall(GradientMisfit(tao, x, V, _ctx));
252:   PetscCall(VecCopy(x, temp));
253:   /* temp = x - z + u */
254:   PetscCall(VecAXPBYPCZ(temp, -1., 1., 1., z, u));
255:   /* V =  FTFx - FTd  mu*(x - z + u) */
256:   PetscCall(VecAXPY(V, mu, temp));
257:   PetscFunctionReturn(PETSC_SUCCESS);
258: }

260: /* returns FTF + diag(mu) */
261: static PetscErrorCode HessianMisfitADMM(Tao tao, Vec x, Mat H, Mat Hpre, void *_ctx)
262: {
263:   UserCtx ctx = (UserCtx)_ctx;

265:   PetscFunctionBegin;
266:   PetscCall(MatCopy(ctx->W, H, DIFFERENT_NONZERO_PATTERN));
267:   PetscCall(MatShift(H, ctx->mu));
268:   if (Hpre != H) PetscCall(MatCopy(H, Hpre, DIFFERENT_NONZERO_PATTERN));
269:   PetscFunctionReturn(PETSC_SUCCESS);
270: }

272: /* computes || x ||_p (mult by 0.5 in case of NORM_2) */
273: static PetscErrorCode ObjectiveRegularization(Tao tao, Vec x, PetscReal *J, void *_ctx)
274: {
275:   UserCtx   ctx = (UserCtx)_ctx;
276:   PetscReal norm;

278:   PetscFunctionBegin;
279:   *J = 0;
280:   PetscCall(VecNorm(x, ctx->p, &norm));
281:   if (ctx->p == NORM_2) norm = 0.5 * norm * norm;
282:   *J = ctx->alpha * norm;
283:   PetscFunctionReturn(PETSC_SUCCESS);
284: }

286: /* NORM_2 Case: return x
287:  * NORM_1 Case: x/(|x| + eps)
288:  * Else: TODO */
289: static PetscErrorCode GradientRegularization(Tao tao, Vec x, Vec V, void *_ctx)
290: {
291:   UserCtx   ctx = (UserCtx)_ctx;
292:   PetscReal eps = ctx->eps;

294:   PetscFunctionBegin;
295:   if (ctx->p == NORM_2) {
296:     PetscCall(VecCopy(x, V));
297:   } else if (ctx->p == NORM_1) {
298:     PetscCall(VecCopy(x, ctx->workRight[1]));
299:     PetscCall(VecAbs(ctx->workRight[1]));
300:     PetscCall(VecShift(ctx->workRight[1], eps));
301:     PetscCall(VecPointwiseDivide(V, x, ctx->workRight[1]));
302:   } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Example only works for NORM_1 and NORM_2");
303:   PetscFunctionReturn(PETSC_SUCCESS);
304: }

306: /* NORM_2 Case: returns diag(mu)
307:  * NORM_1 Case: diag(mu* 1/sqrt(x_i^2 + eps) * (1 - x_i^2/ABS(x_i^2+eps)))  */
308: static PetscErrorCode HessianRegularization(Tao tao, Vec x, Mat H, Mat Hpre, void *_ctx)
309: {
310:   UserCtx   ctx = (UserCtx)_ctx;
311:   PetscReal eps = ctx->eps;
312:   Vec       copy1, copy2, copy3;

314:   PetscFunctionBegin;
315:   if (ctx->p == NORM_2) {
316:     /* Identity matrix scaled by mu */
317:     PetscCall(MatZeroEntries(H));
318:     PetscCall(MatShift(H, ctx->mu));
319:     if (Hpre != H) {
320:       PetscCall(MatZeroEntries(Hpre));
321:       PetscCall(MatShift(Hpre, ctx->mu));
322:     }
323:   } else if (ctx->p == NORM_1) {
324:     /* 1/sqrt(x_i^2 + eps) * (1 - x_i^2/ABS(x_i^2+eps)) */
325:     copy1 = ctx->workRight[1];
326:     copy2 = ctx->workRight[2];
327:     copy3 = ctx->workRight[3];
328:     /* copy1 : 1/sqrt(x_i^2 + eps) */
329:     PetscCall(VecCopy(x, copy1));
330:     PetscCall(VecPow(copy1, 2));
331:     PetscCall(VecShift(copy1, eps));
332:     PetscCall(VecSqrtAbs(copy1));
333:     PetscCall(VecReciprocal(copy1));
334:     /* copy2:  x_i^2.*/
335:     PetscCall(VecCopy(x, copy2));
336:     PetscCall(VecPow(copy2, 2));
337:     /* copy3: abs(x_i^2 + eps) */
338:     PetscCall(VecCopy(x, copy3));
339:     PetscCall(VecPow(copy3, 2));
340:     PetscCall(VecShift(copy3, eps));
341:     PetscCall(VecAbs(copy3));
342:     /* copy2: 1 - x_i^2/abs(x_i^2 + eps) */
343:     PetscCall(VecPointwiseDivide(copy2, copy2, copy3));
344:     PetscCall(VecScale(copy2, -1.));
345:     PetscCall(VecShift(copy2, 1.));
346:     PetscCall(VecAXPY(copy1, 1., copy2));
347:     PetscCall(VecScale(copy1, ctx->mu));
348:     PetscCall(MatZeroEntries(H));
349:     PetscCall(MatDiagonalSet(H, copy1, INSERT_VALUES));
350:     if (Hpre != H) {
351:       PetscCall(MatZeroEntries(Hpre));
352:       PetscCall(MatDiagonalSet(Hpre, copy1, INSERT_VALUES));
353:     }
354:   } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Example only works for NORM_1 and NORM_2");
355:   PetscFunctionReturn(PETSC_SUCCESS);
356: }

358: /* NORM_2 Case: 0.5 || x ||_2 + 0.5 * mu * ||x + u - z||^2
359:  * Else : || x ||_2 + 0.5 * mu * ||x + u - z||^2 */
360: static PetscErrorCode ObjectiveRegularizationADMM(Tao tao, Vec z, PetscReal *J, void *_ctx)
361: {
362:   UserCtx   ctx = (UserCtx)_ctx;
363:   PetscReal mu, workNorm, reg;
364:   Vec       x, u, temp;

366:   PetscFunctionBegin;
367:   mu   = ctx->mu;
368:   x    = ctx->workRight[4];
369:   u    = ctx->workRight[6];
370:   temp = ctx->workRight[10];
371:   PetscCall(ObjectiveRegularization(tao, z, &reg, _ctx));
372:   PetscCall(VecCopy(z, temp));
373:   /* temp = x + u -z */
374:   PetscCall(VecAXPBYPCZ(temp, 1., 1., -1., x, u));
375:   /* workNorm = ||x + u - z ||^2 */
376:   PetscCall(VecDot(temp, temp, &workNorm));
377:   *J = reg + 0.5 * mu * workNorm;
378:   PetscFunctionReturn(PETSC_SUCCESS);
379: }

381: /* NORM_2 Case: x - mu*(x + u - z)
382:  * NORM_1 Case: x/(|x| + eps) - mu*(x + u - z)
383:  * Else: TODO */
384: static PetscErrorCode GradientRegularizationADMM(Tao tao, Vec z, Vec V, void *_ctx)
385: {
386:   UserCtx   ctx = (UserCtx)_ctx;
387:   PetscReal mu;
388:   Vec       x, u, temp;

390:   PetscFunctionBegin;
391:   mu   = ctx->mu;
392:   x    = ctx->workRight[4];
393:   u    = ctx->workRight[6];
394:   temp = ctx->workRight[10];
395:   PetscCall(GradientRegularization(tao, z, V, _ctx));
396:   PetscCall(VecCopy(z, temp));
397:   /* temp = x + u -z */
398:   PetscCall(VecAXPBYPCZ(temp, 1., 1., -1., x, u));
399:   PetscCall(VecAXPY(V, -mu, temp));
400:   PetscFunctionReturn(PETSC_SUCCESS);
401: }

403: /* NORM_2 Case: returns diag(mu)
404:  * NORM_1 Case: FTF + diag(mu) */
405: static PetscErrorCode HessianRegularizationADMM(Tao tao, Vec x, Mat H, Mat Hpre, void *_ctx)
406: {
407:   UserCtx ctx = (UserCtx)_ctx;

409:   PetscFunctionBegin;
410:   if (ctx->p == NORM_2) {
411:     /* Identity matrix scaled by mu */
412:     PetscCall(MatZeroEntries(H));
413:     PetscCall(MatShift(H, ctx->mu));
414:     if (Hpre != H) {
415:       PetscCall(MatZeroEntries(Hpre));
416:       PetscCall(MatShift(Hpre, ctx->mu));
417:     }
418:   } else if (ctx->p == NORM_1) {
419:     PetscCall(HessianMisfit(tao, x, H, Hpre, (void *)ctx));
420:     PetscCall(MatShift(H, ctx->mu));
421:     if (Hpre != H) PetscCall(MatShift(Hpre, ctx->mu));
422:   } else SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_OUTOFRANGE, "Example only works for NORM_1 and NORM_2");
423:   PetscFunctionReturn(PETSC_SUCCESS);
424: }

426: /* NORM_2 Case : (1/2) * ||F x - d||^2 + 0.5 * || x ||_p
427: *  NORM_1 Case : (1/2) * ||F x - d||^2 + || x ||_p */
428: static PetscErrorCode ObjectiveComplete(Tao tao, Vec x, PetscReal *J, void *ctx)
429: {
430:   PetscReal Jm, Jr;

432:   PetscFunctionBegin;
433:   PetscCall(ObjectiveMisfit(tao, x, &Jm, ctx));
434:   PetscCall(ObjectiveRegularization(tao, x, &Jr, ctx));
435:   *J = Jm + Jr;
436:   PetscFunctionReturn(PETSC_SUCCESS);
437: }

439: /* NORM_2 Case: FTFx - FTd + x
440:  * NORM_1 Case: FTFx - FTd + x/(|x| + eps) */
441: static PetscErrorCode GradientComplete(Tao tao, Vec x, Vec V, void *ctx)
442: {
443:   UserCtx cntx = (UserCtx)ctx;

445:   PetscFunctionBegin;
446:   PetscCall(GradientMisfit(tao, x, cntx->workRight[2], ctx));
447:   PetscCall(GradientRegularization(tao, x, cntx->workRight[3], ctx));
448:   PetscCall(VecWAXPY(V, 1, cntx->workRight[2], cntx->workRight[3]));
449:   PetscFunctionReturn(PETSC_SUCCESS);
450: }

452: /* NORM_2 Case: diag(mu) + FTF
453:  * NORM_1 Case: diag(mu* 1/sqrt(x_i^2 + eps) * (1 - x_i^2/ABS(x_i^2+eps))) + FTF  */
454: static PetscErrorCode HessianComplete(Tao tao, Vec x, Mat H, Mat Hpre, void *ctx)
455: {
456:   Mat tempH;

458:   PetscFunctionBegin;
459:   PetscCall(MatDuplicate(H, MAT_SHARE_NONZERO_PATTERN, &tempH));
460:   PetscCall(HessianMisfit(tao, x, H, H, ctx));
461:   PetscCall(HessianRegularization(tao, x, tempH, tempH, ctx));
462:   PetscCall(MatAXPY(H, 1., tempH, DIFFERENT_NONZERO_PATTERN));
463:   if (Hpre != H) PetscCall(MatCopy(H, Hpre, DIFFERENT_NONZERO_PATTERN));
464:   PetscCall(MatDestroy(&tempH));
465:   PetscFunctionReturn(PETSC_SUCCESS);
466: }

468: static PetscErrorCode TaoSolveADMM(UserCtx ctx, Vec x)
469: {
470:   PetscInt  i;
471:   PetscReal u_norm, r_norm, s_norm, primal, dual, x_norm, z_norm;
472:   Tao       tao1, tao2;
473:   Vec       xk, z, u, diff, zold, zdiff, temp;
474:   PetscReal mu;

476:   PetscFunctionBegin;
477:   xk    = ctx->workRight[4];
478:   z     = ctx->workRight[5];
479:   u     = ctx->workRight[6];
480:   diff  = ctx->workRight[7];
481:   zold  = ctx->workRight[8];
482:   zdiff = ctx->workRight[9];
483:   temp  = ctx->workRight[11];
484:   mu    = ctx->mu;
485:   PetscCall(VecSet(u, 0.));
486:   PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao1));
487:   PetscCall(TaoSetType(tao1, TAONLS));
488:   PetscCall(TaoSetObjective(tao1, ObjectiveMisfitADMM, (void *)ctx));
489:   PetscCall(TaoSetGradient(tao1, NULL, GradientMisfitADMM, (void *)ctx));
490:   PetscCall(TaoSetHessian(tao1, ctx->Hm, ctx->Hm, HessianMisfitADMM, (void *)ctx));
491:   PetscCall(VecSet(xk, 0.));
492:   PetscCall(TaoSetSolution(tao1, xk));
493:   PetscCall(TaoSetOptionsPrefix(tao1, "misfit_"));
494:   PetscCall(TaoSetFromOptions(tao1));
495:   PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao2));
496:   if (ctx->p == NORM_2) {
497:     PetscCall(TaoSetType(tao2, TAONLS));
498:     PetscCall(TaoSetObjective(tao2, ObjectiveRegularizationADMM, (void *)ctx));
499:     PetscCall(TaoSetGradient(tao2, NULL, GradientRegularizationADMM, (void *)ctx));
500:     PetscCall(TaoSetHessian(tao2, ctx->Hr, ctx->Hr, HessianRegularizationADMM, (void *)ctx));
501:   }
502:   PetscCall(VecSet(z, 0.));
503:   PetscCall(TaoSetSolution(tao2, z));
504:   PetscCall(TaoSetOptionsPrefix(tao2, "reg_"));
505:   PetscCall(TaoSetFromOptions(tao2));

507:   for (i = 0; i < ctx->iter; i++) {
508:     PetscCall(VecCopy(z, zold));
509:     PetscCall(TaoSolve(tao1)); /* Updates xk */
510:     if (ctx->p == NORM_1) {
511:       PetscCall(VecWAXPY(temp, 1., xk, u));
512:       PetscCall(TaoSoftThreshold(temp, -ctx->alpha / mu, ctx->alpha / mu, z));
513:     } else {
514:       PetscCall(TaoSolve(tao2)); /* Update zk */
515:     }
516:     /* u = u + xk -z */
517:     PetscCall(VecAXPBYPCZ(u, 1., -1., 1., xk, z));
518:     /* r_norm : norm(x-z) */
519:     PetscCall(VecWAXPY(diff, -1., z, xk));
520:     PetscCall(VecNorm(diff, NORM_2, &r_norm));
521:     /* s_norm : norm(-mu(z-zold)) */
522:     PetscCall(VecWAXPY(zdiff, -1., zold, z));
523:     PetscCall(VecNorm(zdiff, NORM_2, &s_norm));
524:     s_norm = s_norm * mu;
525:     /* primal : sqrt(n)*ABSTOL + RELTOL*max(norm(x), norm(-z))*/
526:     PetscCall(VecNorm(xk, NORM_2, &x_norm));
527:     PetscCall(VecNorm(z, NORM_2, &z_norm));
528:     primal = PetscSqrtReal(ctx->n) * ctx->abstol + ctx->reltol * PetscMax(x_norm, z_norm);
529:     /* Duality : sqrt(n)*ABSTOL + RELTOL*norm(mu*u)*/
530:     PetscCall(VecNorm(u, NORM_2, &u_norm));
531:     dual = PetscSqrtReal(ctx->n) * ctx->abstol + ctx->reltol * u_norm * mu;
532:     PetscCall(PetscPrintf(PetscObjectComm((PetscObject)tao1), "Iter %" PetscInt_FMT " : ||x-z||: %g, mu*||z-zold||: %g\n", i, (double)r_norm, (double)s_norm));
533:     if (r_norm < primal && s_norm < dual) break;
534:   }
535:   PetscCall(VecCopy(xk, x));
536:   PetscCall(TaoDestroy(&tao1));
537:   PetscCall(TaoDestroy(&tao2));
538:   PetscFunctionReturn(PETSC_SUCCESS);
539: }

541: /* Second order Taylor remainder convergence test */
542: static PetscErrorCode TaylorTest(UserCtx ctx, Tao tao, Vec x, PetscReal *C)
543: {
544:   PetscReal  h, J, temp;
545:   PetscInt   i, j;
546:   PetscInt   numValues;
547:   PetscReal  Jx, Jxhat_comp, Jxhat_pred;
548:   PetscReal *Js, *hs;
549:   PetscReal  gdotdx;
550:   PetscReal  minrate = PETSC_MAX_REAL;
551:   MPI_Comm   comm    = PetscObjectComm((PetscObject)x);
552:   Vec        g, dx, xhat;

554:   PetscFunctionBegin;
555:   PetscCall(VecDuplicate(x, &g));
556:   PetscCall(VecDuplicate(x, &xhat));
557:   /* choose a perturbation direction */
558:   PetscCall(VecDuplicate(x, &dx));
559:   PetscCall(VecSetRandom(dx, ctx->rctx));
560:   /* evaluate objective at x: J(x) */
561:   PetscCall(TaoComputeObjective(tao, x, &Jx));
562:   /* evaluate gradient at x, save in vector g */
563:   PetscCall(TaoComputeGradient(tao, x, g));
564:   PetscCall(VecDot(g, dx, &gdotdx));

566:   for (numValues = 0, h = ctx->hStart; h >= ctx->hMin; h *= ctx->hFactor) numValues++;
567:   PetscCall(PetscCalloc2(numValues, &Js, numValues, &hs));
568:   for (i = 0, h = ctx->hStart; h >= ctx->hMin; h *= ctx->hFactor, i++) {
569:     PetscCall(VecWAXPY(xhat, h, dx, x));
570:     PetscCall(TaoComputeObjective(tao, xhat, &Jxhat_comp));
571:     /* J(\hat(x)) \approx J(x) + g^T (xhat - x) = J(x) + h * g^T dx */
572:     Jxhat_pred = Jx + h * gdotdx;
573:     /* Vector to dJdm scalar? Dot?*/
574:     J = PetscAbsReal(Jxhat_comp - Jxhat_pred);
575:     PetscCall(PetscPrintf(comm, "J(xhat): %g, predicted: %g, diff %g\n", (double)Jxhat_comp, (double)Jxhat_pred, (double)J));
576:     Js[i] = J;
577:     hs[i] = h;
578:   }
579:   for (j = 1; j < numValues; j++) {
580:     temp = PetscLogReal(Js[j] / Js[j - 1]) / PetscLogReal(hs[j] / hs[j - 1]);
581:     PetscCall(PetscPrintf(comm, "Convergence rate step %" PetscInt_FMT ": %g\n", j - 1, (double)temp));
582:     minrate = PetscMin(minrate, temp);
583:   }
584:   /* If O is not ~2, then the test is wrong */
585:   PetscCall(PetscFree2(Js, hs));
586:   *C = minrate;
587:   PetscCall(VecDestroy(&dx));
588:   PetscCall(VecDestroy(&xhat));
589:   PetscCall(VecDestroy(&g));
590:   PetscFunctionReturn(PETSC_SUCCESS);
591: }

593: int main(int argc, char **argv)
594: {
595:   UserCtx ctx;
596:   Tao     tao;
597:   Vec     x;
598:   Mat     H;

600:   PetscFunctionBeginUser;
601:   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
602:   PetscCall(PetscNew(&ctx));
603:   PetscCall(ConfigureContext(ctx));
604:   /* Define two functions that could pass as objectives to TaoSetObjective(): one
605:    * for the misfit component, and one for the regularization component */
606:   /* ObjectiveMisfit() and ObjectiveRegularization() */

608:   /* Define a single function that calls both components adds them together: the complete objective,
609:    * in the absence of a Tao implementation that handles separability */
610:   /* ObjectiveComplete() */
611:   PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao));
612:   PetscCall(TaoSetType(tao, TAONM));
613:   PetscCall(TaoSetObjective(tao, ObjectiveComplete, (void *)ctx));
614:   PetscCall(TaoSetGradient(tao, NULL, GradientComplete, (void *)ctx));
615:   PetscCall(MatDuplicate(ctx->W, MAT_SHARE_NONZERO_PATTERN, &H));
616:   PetscCall(TaoSetHessian(tao, H, H, HessianComplete, (void *)ctx));
617:   PetscCall(MatCreateVecs(ctx->F, NULL, &x));
618:   PetscCall(VecSet(x, 0.));
619:   PetscCall(TaoSetSolution(tao, x));
620:   PetscCall(TaoSetFromOptions(tao));
621:   if (ctx->use_admm) PetscCall(TaoSolveADMM(ctx, x));
622:   else PetscCall(TaoSolve(tao));
623:   /* examine solution */
624:   PetscCall(VecViewFromOptions(x, NULL, "-view_sol"));
625:   if (ctx->taylor) {
626:     PetscReal rate;
627:     PetscCall(TaylorTest(ctx, tao, x, &rate));
628:   }
629:   PetscCall(MatDestroy(&H));
630:   PetscCall(TaoDestroy(&tao));
631:   PetscCall(VecDestroy(&x));
632:   PetscCall(DestroyContext(&ctx));
633:   PetscCall(PetscFinalize());
634:   return 0;
635: }

637: /*TEST

639:   build:
640:     requires: !complex

642:   test:
643:     suffix: 0
644:     args:

646:   test:
647:     suffix: l1_1
648:     args: -p 1 -tao_type lmvm -alpha 1. -epsilon 1.e-7 -m 64 -n 64 -view_sol -matrix_format 1

650:   test:
651:     suffix: hessian_1
652:     args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 1 -tao_type nls

654:   test:
655:     suffix: hessian_2
656:     args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 2 -tao_type nls

658:   test:
659:     suffix: nm_1
660:     args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 1 -tao_type nm -tao_max_it 50

662:   test:
663:     suffix: nm_2
664:     args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 2 -tao_type nm -tao_max_it 50

666:   test:
667:     suffix: lmvm_1
668:     args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 1 -tao_type lmvm -tao_max_it 40

670:   test:
671:     suffix: lmvm_2
672:     args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 2 -tao_type lmvm -tao_max_it 15

674:   test:
675:     suffix: soft_threshold_admm_1
676:     args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 1 -use_admm

678:   test:
679:     suffix: hessian_admm_1
680:     args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 1 -use_admm -reg_tao_type nls -misfit_tao_type nls

682:   test:
683:     suffix: hessian_admm_2
684:     args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 2 -use_admm -reg_tao_type nls -misfit_tao_type nls

686:   test:
687:     suffix: nm_admm_1
688:     args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 1 -use_admm -reg_tao_type nm -misfit_tao_type nm

690:   test:
691:     suffix: nm_admm_2
692:     args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 2 -use_admm -reg_tao_type nm -misfit_tao_type nm -iter 7

694:   test:
695:     suffix: lmvm_admm_1
696:     args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 1 -use_admm -reg_tao_type lmvm -misfit_tao_type lmvm

698:   test:
699:     suffix: lmvm_admm_2
700:     args: -matrix_format 1 -m 100 -n 100 -tao_monitor -p 2 -use_admm -reg_tao_type lmvm -misfit_tao_type lmvm

702: TEST*/