Actual source code: elliptic.c

  1: #include <petsc/private/taoimpl.h>

  3: typedef struct {
  4:   PetscInt n; /* Number of total variables */
  5:   PetscInt m; /* Number of constraints */
  6:   PetscInt nstate;
  7:   PetscInt ndesign;
  8:   PetscInt mx;    /* grid points in each direction */
  9:   PetscInt ns;    /* Number of data samples (1<=ns<=8)
 10:                   Currently only ns=1 is supported */
 11:   PetscInt ndata; /* Number of data points per sample */
 12:   IS       s_is;
 13:   IS       d_is;

 15:   VecScatter  state_scatter;
 16:   VecScatter  design_scatter;
 17:   VecScatter *yi_scatter, *di_scatter;
 18:   Vec         suby, subq, subd;
 19:   Mat         Js, Jd, JsPrec, JsInv, JsBlock;

 21:   PetscReal  alpha; /* Regularization parameter */
 22:   PetscReal  beta;  /* Weight attributed to ||u||^2 in regularization functional */
 23:   PetscReal  noise; /* Amount of noise to add to data */
 24:   PetscReal *ones;
 25:   Mat        Q;
 26:   Mat        MQ;
 27:   Mat        L;

 29:   Mat Grad;
 30:   Mat Av, Avwork;
 31:   Mat Div, Divwork;
 32:   Mat DSG;
 33:   Mat Diag, Ones;

 35:   Vec q;
 36:   Vec ur; /* reference */

 38:   Vec d;
 39:   Vec dwork;

 41:   Vec x; /* super vec of y,u */

 43:   Vec y; /* state variables */
 44:   Vec ywork;

 46:   Vec ytrue;

 48:   Vec u; /* design variables */
 49:   Vec uwork;

 51:   Vec utrue;

 53:   Vec js_diag;

 55:   Vec c; /* constraint vector */
 56:   Vec cwork;

 58:   Vec lwork;
 59:   Vec S;
 60:   Vec Swork, Twork, Sdiag, Ywork;
 61:   Vec Av_u;

 63:   KSP solver;
 64:   PC  prec;

 66:   PetscReal     tola, tolb, tolc, told;
 67:   PetscInt      ksp_its;
 68:   PetscInt      ksp_its_initial;
 69:   PetscLogStage stages[10];
 70:   PetscBool     use_ptap;
 71:   PetscBool     use_lrc;
 72: } AppCtx;

 74: PetscErrorCode FormFunction(Tao, Vec, PetscReal *, void *);
 75: PetscErrorCode FormGradient(Tao, Vec, Vec, void *);
 76: PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal *, Vec, void *);
 77: PetscErrorCode FormJacobianState(Tao, Vec, Mat, Mat, Mat, void *);
 78: PetscErrorCode FormJacobianDesign(Tao, Vec, Mat, void *);
 79: PetscErrorCode FormConstraints(Tao, Vec, Vec, void *);
 80: PetscErrorCode FormHessian(Tao, Vec, Mat, Mat, void *);
 81: PetscErrorCode Gather(Vec, Vec, VecScatter, Vec, VecScatter);
 82: PetscErrorCode Scatter(Vec, Vec, VecScatter, Vec, VecScatter);
 83: PetscErrorCode EllipticInitialize(AppCtx *);
 84: PetscErrorCode EllipticDestroy(AppCtx *);
 85: PetscErrorCode EllipticMonitor(Tao, void *);

 87: PetscErrorCode StateBlockMatMult(Mat, Vec, Vec);
 88: PetscErrorCode StateMatMult(Mat, Vec, Vec);

 90: PetscErrorCode StateInvMatMult(Mat, Vec, Vec);
 91: PetscErrorCode DesignMatMult(Mat, Vec, Vec);
 92: PetscErrorCode DesignMatMultTranspose(Mat, Vec, Vec);

 94: PetscErrorCode QMatMult(Mat, Vec, Vec);
 95: PetscErrorCode QMatMultTranspose(Mat, Vec, Vec);

 97: static char help[] = "";

 99: int main(int argc, char **argv)
100: {
101:   Vec      x0;
102:   Tao      tao;
103:   AppCtx   user;
104:   PetscInt ntests = 1;

106:   PetscFunctionBeginUser;
107:   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
108:   user.mx = 8;
109:   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "elliptic example", NULL);
110:   PetscCall(PetscOptionsInt("-mx", "Number of grid points in each direction", "", user.mx, &user.mx, NULL));
111:   user.ns = 6;
112:   PetscCall(PetscOptionsInt("-ns", "Number of data samples (1<=ns<=8)", "", user.ns, &user.ns, NULL));
113:   user.ndata = 64;
114:   PetscCall(PetscOptionsInt("-ndata", "Numbers of data points per sample", "", user.ndata, &user.ndata, NULL));
115:   user.alpha = 0.1;
116:   PetscCall(PetscOptionsReal("-alpha", "Regularization parameter", "", user.alpha, &user.alpha, NULL));
117:   user.beta = 0.00001;
118:   PetscCall(PetscOptionsReal("-beta", "Weight attributed to ||u||^2 in regularization functional", "", user.beta, &user.beta, NULL));
119:   user.noise = 0.01;
120:   PetscCall(PetscOptionsReal("-noise", "Amount of noise to add to data", "", user.noise, &user.noise, NULL));

122:   user.use_ptap = PETSC_FALSE;
123:   PetscCall(PetscOptionsBool("-use_ptap", "Use ptap matrix for DSG", "", user.use_ptap, &user.use_ptap, NULL));
124:   user.use_lrc = PETSC_FALSE;
125:   PetscCall(PetscOptionsBool("-use_lrc", "Use lrc matrix for Js", "", user.use_lrc, &user.use_lrc, NULL));
126:   PetscCall(PetscOptionsInt("-ntests", "Number of times to repeat TaoSolve", "", ntests, &ntests, NULL));
127:   PetscOptionsEnd();

129:   user.m       = user.ns * user.mx * user.mx * user.mx; /* number of constraints */
130:   user.nstate  = user.m;
131:   user.ndesign = user.mx * user.mx * user.mx;
132:   user.n       = user.nstate + user.ndesign; /* number of variables */

134:   /* Create TAO solver and set desired solution method */
135:   PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao));
136:   PetscCall(TaoSetType(tao, TAOLCL));

138:   /* Set up initial vectors and matrices */
139:   PetscCall(EllipticInitialize(&user));

141:   PetscCall(Gather(user.x, user.y, user.state_scatter, user.u, user.design_scatter));
142:   PetscCall(VecDuplicate(user.x, &x0));
143:   PetscCall(VecCopy(user.x, x0));

145:   /* Set solution vector with an initial guess */
146:   PetscCall(TaoSetSolution(tao, user.x));
147:   PetscCall(TaoSetObjective(tao, FormFunction, &user));
148:   PetscCall(TaoSetGradient(tao, NULL, FormGradient, &user));
149:   PetscCall(TaoSetConstraintsRoutine(tao, user.c, FormConstraints, &user));

151:   PetscCall(TaoSetJacobianStateRoutine(tao, user.Js, NULL, user.JsInv, FormJacobianState, &user));
152:   PetscCall(TaoSetJacobianDesignRoutine(tao, user.Jd, FormJacobianDesign, &user));

154:   PetscCall(TaoSetStateDesignIS(tao, user.s_is, user.d_is));
155:   PetscCall(TaoSetFromOptions(tao));

157:   /* SOLVE THE APPLICATION */
158:   PetscCall(PetscLogStageRegister("Trials", &user.stages[1]));
159:   PetscCall(PetscLogStagePush(user.stages[1]));
160:   for (PetscInt i = 0; i < ntests; i++) {
161:     PetscCall(TaoSolve(tao));
162:     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "KSP Iterations = %" PetscInt_FMT "\n", user.ksp_its));
163:     PetscCall(VecCopy(x0, user.x));
164:   }
165:   PetscCall(PetscLogStagePop());
166:   PetscCall(PetscBarrier((PetscObject)user.x));
167:   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "KSP iterations within initialization: "));
168:   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%" PetscInt_FMT "\n", user.ksp_its_initial));

170:   PetscCall(TaoDestroy(&tao));
171:   PetscCall(VecDestroy(&x0));
172:   PetscCall(EllipticDestroy(&user));
173:   PetscCall(PetscFinalize());
174:   return 0;
175: }
176: /* ------------------------------------------------------------------- */
177: /*
178:    dwork = Qy - d
179:    lwork = L*(u-ur)
180:    f = 1/2 * (dwork.dwork + alpha*lwork.lwork)
181: */
182: PetscErrorCode FormFunction(Tao tao, Vec X, PetscReal *f, void *ptr)
183: {
184:   PetscReal d1 = 0, d2 = 0;
185:   AppCtx   *user = (AppCtx *)ptr;

187:   PetscFunctionBegin;
188:   PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter));
189:   PetscCall(MatMult(user->MQ, user->y, user->dwork));
190:   PetscCall(VecAXPY(user->dwork, -1.0, user->d));
191:   PetscCall(VecDot(user->dwork, user->dwork, &d1));
192:   PetscCall(VecWAXPY(user->uwork, -1.0, user->ur, user->u));
193:   PetscCall(MatMult(user->L, user->uwork, user->lwork));
194:   PetscCall(VecDot(user->lwork, user->lwork, &d2));
195:   *f = 0.5 * (d1 + user->alpha * d2);
196:   PetscFunctionReturn(PETSC_SUCCESS);
197: }

199: /* ------------------------------------------------------------------- */
200: /*
201:     state: g_s = Q' *(Qy - d)
202:     design: g_d = alpha*L'*L*(u-ur)
203: */
204: PetscErrorCode FormGradient(Tao tao, Vec X, Vec G, void *ptr)
205: {
206:   AppCtx *user = (AppCtx *)ptr;

208:   PetscFunctionBegin;
209:   PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter));
210:   PetscCall(MatMult(user->MQ, user->y, user->dwork));
211:   PetscCall(VecAXPY(user->dwork, -1.0, user->d));
212:   PetscCall(MatMultTranspose(user->MQ, user->dwork, user->ywork));
213:   PetscCall(VecWAXPY(user->uwork, -1.0, user->ur, user->u));
214:   PetscCall(MatMult(user->L, user->uwork, user->lwork));
215:   PetscCall(MatMultTranspose(user->L, user->lwork, user->uwork));
216:   PetscCall(VecScale(user->uwork, user->alpha));
217:   PetscCall(Gather(G, user->ywork, user->state_scatter, user->uwork, user->design_scatter));
218:   PetscFunctionReturn(PETSC_SUCCESS);
219: }

221: PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, void *ptr)
222: {
223:   PetscReal d1, d2;
224:   AppCtx   *user = (AppCtx *)ptr;

226:   PetscFunctionBegin;
227:   PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter));
228:   PetscCall(MatMult(user->MQ, user->y, user->dwork));
229:   PetscCall(VecAXPY(user->dwork, -1.0, user->d));
230:   PetscCall(VecDot(user->dwork, user->dwork, &d1));
231:   PetscCall(MatMultTranspose(user->MQ, user->dwork, user->ywork));

233:   PetscCall(VecWAXPY(user->uwork, -1.0, user->ur, user->u));
234:   PetscCall(MatMult(user->L, user->uwork, user->lwork));
235:   PetscCall(VecDot(user->lwork, user->lwork, &d2));
236:   PetscCall(MatMultTranspose(user->L, user->lwork, user->uwork));
237:   PetscCall(VecScale(user->uwork, user->alpha));
238:   *f = 0.5 * (d1 + user->alpha * d2);
239:   PetscCall(Gather(G, user->ywork, user->state_scatter, user->uwork, user->design_scatter));
240:   PetscFunctionReturn(PETSC_SUCCESS);
241: }

243: /* ------------------------------------------------------------------- */
244: /* A
245: MatShell object
246: */
247: PetscErrorCode FormJacobianState(Tao tao, Vec X, Mat J, Mat JPre, Mat JInv, void *ptr)
248: {
249:   AppCtx *user = (AppCtx *)ptr;

251:   PetscFunctionBegin;
252:   PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter));
253:   /* DSG = Div * (1/Av_u) * Grad */
254:   PetscCall(VecSet(user->uwork, 0));
255:   PetscCall(VecAXPY(user->uwork, -1.0, user->u));
256:   PetscCall(VecExp(user->uwork));
257:   PetscCall(MatMult(user->Av, user->uwork, user->Av_u));
258:   PetscCall(VecCopy(user->Av_u, user->Swork));
259:   PetscCall(VecReciprocal(user->Swork));
260:   if (user->use_ptap) {
261:     PetscCall(MatDiagonalSet(user->Diag, user->Swork, INSERT_VALUES));
262:     PetscCall(MatPtAP(user->Diag, user->Grad, MAT_REUSE_MATRIX, 1.0, &user->DSG));
263:   } else {
264:     PetscCall(MatCopy(user->Div, user->Divwork, SAME_NONZERO_PATTERN));
265:     PetscCall(MatDiagonalScale(user->Divwork, NULL, user->Swork));
266:     PetscCall(MatProductNumeric(user->DSG));
267:   }
268:   PetscFunctionReturn(PETSC_SUCCESS);
269: }
270: /* ------------------------------------------------------------------- */
271: /* B */
272: PetscErrorCode FormJacobianDesign(Tao tao, Vec X, Mat J, void *ptr)
273: {
274:   AppCtx *user = (AppCtx *)ptr;

276:   PetscFunctionBegin;
277:   PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter));
278:   PetscFunctionReturn(PETSC_SUCCESS);
279: }

281: PetscErrorCode StateBlockMatMult(Mat J_shell, Vec X, Vec Y)
282: {
283:   PetscReal sum;
284:   AppCtx   *user;

286:   PetscFunctionBegin;
287:   PetscCall(MatShellGetContext(J_shell, &user));
288:   PetscCall(MatMult(user->DSG, X, Y));
289:   PetscCall(VecSum(X, &sum));
290:   sum /= user->ndesign;
291:   PetscCall(VecShift(Y, sum));
292:   PetscFunctionReturn(PETSC_SUCCESS);
293: }

295: PetscErrorCode StateMatMult(Mat J_shell, Vec X, Vec Y)
296: {
297:   PetscInt i;
298:   AppCtx  *user;

300:   PetscFunctionBegin;
301:   PetscCall(MatShellGetContext(J_shell, &user));
302:   if (user->ns == 1) {
303:     PetscCall(MatMult(user->JsBlock, X, Y));
304:   } else {
305:     for (i = 0; i < user->ns; i++) {
306:       PetscCall(Scatter(X, user->subq, user->yi_scatter[i], 0, 0));
307:       PetscCall(Scatter(Y, user->suby, user->yi_scatter[i], 0, 0));
308:       PetscCall(MatMult(user->JsBlock, user->subq, user->suby));
309:       PetscCall(Gather(Y, user->suby, user->yi_scatter[i], 0, 0));
310:     }
311:   }
312:   PetscFunctionReturn(PETSC_SUCCESS);
313: }

315: PetscErrorCode StateInvMatMult(Mat J_shell, Vec X, Vec Y)
316: {
317:   PetscInt its, i;
318:   AppCtx  *user;

320:   PetscFunctionBegin;
321:   PetscCall(MatShellGetContext(J_shell, &user));
322:   PetscCall(KSPSetOperators(user->solver, user->JsBlock, user->DSG));
323:   if (Y == user->ytrue) {
324:     /* First solve is done using true solution to set up problem */
325:     PetscCall(KSPSetTolerances(user->solver, 1e-8, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT));
326:   } else {
327:     PetscCall(KSPSetTolerances(user->solver, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT));
328:   }
329:   if (user->ns == 1) {
330:     PetscCall(KSPSolve(user->solver, X, Y));
331:     PetscCall(KSPGetIterationNumber(user->solver, &its));
332:     user->ksp_its += its;
333:   } else {
334:     for (i = 0; i < user->ns; i++) {
335:       PetscCall(Scatter(X, user->subq, user->yi_scatter[i], 0, 0));
336:       PetscCall(Scatter(Y, user->suby, user->yi_scatter[i], 0, 0));
337:       PetscCall(KSPSolve(user->solver, user->subq, user->suby));
338:       PetscCall(KSPGetIterationNumber(user->solver, &its));
339:       user->ksp_its += its;
340:       PetscCall(Gather(Y, user->suby, user->yi_scatter[i], 0, 0));
341:     }
342:   }
343:   PetscFunctionReturn(PETSC_SUCCESS);
344: }
345: PetscErrorCode QMatMult(Mat J_shell, Vec X, Vec Y)
346: {
347:   AppCtx  *user;
348:   PetscInt i;

350:   PetscFunctionBegin;
351:   PetscCall(MatShellGetContext(J_shell, &user));
352:   if (user->ns == 1) {
353:     PetscCall(MatMult(user->Q, X, Y));
354:   } else {
355:     for (i = 0; i < user->ns; i++) {
356:       PetscCall(Scatter(X, user->subq, user->yi_scatter[i], 0, 0));
357:       PetscCall(Scatter(Y, user->subd, user->di_scatter[i], 0, 0));
358:       PetscCall(MatMult(user->Q, user->subq, user->subd));
359:       PetscCall(Gather(Y, user->subd, user->di_scatter[i], 0, 0));
360:     }
361:   }
362:   PetscFunctionReturn(PETSC_SUCCESS);
363: }

365: PetscErrorCode QMatMultTranspose(Mat J_shell, Vec X, Vec Y)
366: {
367:   AppCtx  *user;
368:   PetscInt i;

370:   PetscFunctionBegin;
371:   PetscCall(MatShellGetContext(J_shell, &user));
372:   if (user->ns == 1) {
373:     PetscCall(MatMultTranspose(user->Q, X, Y));
374:   } else {
375:     for (i = 0; i < user->ns; i++) {
376:       PetscCall(Scatter(X, user->subd, user->di_scatter[i], 0, 0));
377:       PetscCall(Scatter(Y, user->suby, user->yi_scatter[i], 0, 0));
378:       PetscCall(MatMultTranspose(user->Q, user->subd, user->suby));
379:       PetscCall(Gather(Y, user->suby, user->yi_scatter[i], 0, 0));
380:     }
381:   }
382:   PetscFunctionReturn(PETSC_SUCCESS);
383: }

385: PetscErrorCode DesignMatMult(Mat J_shell, Vec X, Vec Y)
386: {
387:   PetscInt i;
388:   AppCtx  *user;

390:   PetscFunctionBegin;
391:   PetscCall(MatShellGetContext(J_shell, &user));

393:   /* sdiag(1./v) */
394:   PetscCall(VecSet(user->uwork, 0));
395:   PetscCall(VecAXPY(user->uwork, -1.0, user->u));
396:   PetscCall(VecExp(user->uwork));

398:   /* sdiag(1./((Av*(1./v)).^2)) */
399:   PetscCall(MatMult(user->Av, user->uwork, user->Swork));
400:   PetscCall(VecPointwiseMult(user->Swork, user->Swork, user->Swork));
401:   PetscCall(VecReciprocal(user->Swork));

403:   /* (Av * (sdiag(1./v) * b)) */
404:   PetscCall(VecPointwiseMult(user->uwork, user->uwork, X));
405:   PetscCall(MatMult(user->Av, user->uwork, user->Twork));

407:   /* (sdiag(1./((Av*(1./v)).^2)) * (Av * (sdiag(1./v) * b))) */
408:   PetscCall(VecPointwiseMult(user->Swork, user->Twork, user->Swork));

410:   if (user->ns == 1) {
411:     /* (sdiag(Grad*y(:,i)) */
412:     PetscCall(MatMult(user->Grad, user->y, user->Twork));

414:     /* Div * (sdiag(Grad*y(:,i)) * (sdiag(1./((Av*(1./v)).^2)) * (Av * (sdiag(1./v) * b)))) */
415:     PetscCall(VecPointwiseMult(user->Swork, user->Twork, user->Swork));
416:     PetscCall(MatMultTranspose(user->Grad, user->Swork, Y));
417:   } else {
418:     for (i = 0; i < user->ns; i++) {
419:       PetscCall(Scatter(user->y, user->suby, user->yi_scatter[i], 0, 0));
420:       PetscCall(Scatter(Y, user->subq, user->yi_scatter[i], 0, 0));

422:       PetscCall(MatMult(user->Grad, user->suby, user->Twork));
423:       PetscCall(VecPointwiseMult(user->Twork, user->Twork, user->Swork));
424:       PetscCall(MatMultTranspose(user->Grad, user->Twork, user->subq));
425:       PetscCall(Gather(user->y, user->suby, user->yi_scatter[i], 0, 0));
426:       PetscCall(Gather(Y, user->subq, user->yi_scatter[i], 0, 0));
427:     }
428:   }
429:   PetscFunctionReturn(PETSC_SUCCESS);
430: }

432: PetscErrorCode DesignMatMultTranspose(Mat J_shell, Vec X, Vec Y)
433: {
434:   PetscInt i;
435:   AppCtx  *user;

437:   PetscFunctionBegin;
438:   PetscCall(MatShellGetContext(J_shell, &user));
439:   PetscCall(VecZeroEntries(Y));

441:   /* Sdiag = 1./((Av*(1./v)).^2) */
442:   PetscCall(VecSet(user->uwork, 0));
443:   PetscCall(VecAXPY(user->uwork, -1.0, user->u));
444:   PetscCall(VecExp(user->uwork));
445:   PetscCall(MatMult(user->Av, user->uwork, user->Swork));
446:   PetscCall(VecPointwiseMult(user->Sdiag, user->Swork, user->Swork));
447:   PetscCall(VecReciprocal(user->Sdiag));

449:   for (i = 0; i < user->ns; i++) {
450:     PetscCall(Scatter(X, user->subq, user->yi_scatter[i], 0, 0));
451:     PetscCall(Scatter(user->y, user->suby, user->yi_scatter[i], 0, 0));

453:     /* Swork = (Div' * b(:,i)) */
454:     PetscCall(MatMult(user->Grad, user->subq, user->Swork));

456:     /* Twork = Grad*y(:,i) */
457:     PetscCall(MatMult(user->Grad, user->suby, user->Twork));

459:     /* Twork = sdiag(Twork) * Swork */
460:     PetscCall(VecPointwiseMult(user->Twork, user->Swork, user->Twork));

462:     /* Swork = pointwisemult(Sdiag,Twork) */
463:     PetscCall(VecPointwiseMult(user->Swork, user->Twork, user->Sdiag));

465:     /* Ywork = Av' * Swork */
466:     PetscCall(MatMultTranspose(user->Av, user->Swork, user->Ywork));

468:     /* Ywork = pointwisemult(uwork,Ywork) */
469:     PetscCall(VecPointwiseMult(user->Ywork, user->uwork, user->Ywork));
470:     PetscCall(VecAXPY(Y, 1.0, user->Ywork));
471:     PetscCall(Gather(user->y, user->suby, user->yi_scatter[i], 0, 0));
472:   }
473:   PetscFunctionReturn(PETSC_SUCCESS);
474: }

476: PetscErrorCode FormConstraints(Tao tao, Vec X, Vec C, void *ptr)
477: {
478:   /* C=Ay - q      A = Div * Sigma * Grad + hx*hx*hx*ones(n,n) */
479:   PetscReal sum;
480:   PetscInt  i;
481:   AppCtx   *user = (AppCtx *)ptr;

483:   PetscFunctionBegin;
484:   PetscCall(Scatter(X, user->y, user->state_scatter, user->u, user->design_scatter));
485:   if (user->ns == 1) {
486:     PetscCall(MatMult(user->Grad, user->y, user->Swork));
487:     PetscCall(VecPointwiseDivide(user->Swork, user->Swork, user->Av_u));
488:     PetscCall(MatMultTranspose(user->Grad, user->Swork, C));
489:     PetscCall(VecSum(user->y, &sum));
490:     sum /= user->ndesign;
491:     PetscCall(VecShift(C, sum));
492:   } else {
493:     for (i = 0; i < user->ns; i++) {
494:       PetscCall(Scatter(user->y, user->suby, user->yi_scatter[i], 0, 0));
495:       PetscCall(Scatter(C, user->subq, user->yi_scatter[i], 0, 0));
496:       PetscCall(MatMult(user->Grad, user->suby, user->Swork));
497:       PetscCall(VecPointwiseDivide(user->Swork, user->Swork, user->Av_u));
498:       PetscCall(MatMultTranspose(user->Grad, user->Swork, user->subq));

500:       PetscCall(VecSum(user->suby, &sum));
501:       sum /= user->ndesign;
502:       PetscCall(VecShift(user->subq, sum));

504:       PetscCall(Gather(user->y, user->suby, user->yi_scatter[i], 0, 0));
505:       PetscCall(Gather(C, user->subq, user->yi_scatter[i], 0, 0));
506:     }
507:   }
508:   PetscCall(VecAXPY(C, -1.0, user->q));
509:   PetscFunctionReturn(PETSC_SUCCESS);
510: }

512: PetscErrorCode Scatter(Vec x, Vec sub1, VecScatter scat1, Vec sub2, VecScatter scat2)
513: {
514:   PetscFunctionBegin;
515:   PetscCall(VecScatterBegin(scat1, x, sub1, INSERT_VALUES, SCATTER_FORWARD));
516:   PetscCall(VecScatterEnd(scat1, x, sub1, INSERT_VALUES, SCATTER_FORWARD));
517:   if (sub2) {
518:     PetscCall(VecScatterBegin(scat2, x, sub2, INSERT_VALUES, SCATTER_FORWARD));
519:     PetscCall(VecScatterEnd(scat2, x, sub2, INSERT_VALUES, SCATTER_FORWARD));
520:   }
521:   PetscFunctionReturn(PETSC_SUCCESS);
522: }

524: PetscErrorCode Gather(Vec x, Vec sub1, VecScatter scat1, Vec sub2, VecScatter scat2)
525: {
526:   PetscFunctionBegin;
527:   PetscCall(VecScatterBegin(scat1, sub1, x, INSERT_VALUES, SCATTER_REVERSE));
528:   PetscCall(VecScatterEnd(scat1, sub1, x, INSERT_VALUES, SCATTER_REVERSE));
529:   if (sub2) {
530:     PetscCall(VecScatterBegin(scat2, sub2, x, INSERT_VALUES, SCATTER_REVERSE));
531:     PetscCall(VecScatterEnd(scat2, sub2, x, INSERT_VALUES, SCATTER_REVERSE));
532:   }
533:   PetscFunctionReturn(PETSC_SUCCESS);
534: }

536: PetscErrorCode EllipticInitialize(AppCtx *user)
537: {
538:   PetscInt        m, n, i, j, k, l, linear_index, is, js, ks, ls, istart, iend, iblock;
539:   Vec             XX, YY, ZZ, XXwork, YYwork, ZZwork, UTwork;
540:   PetscReal      *x, *y, *z;
541:   PetscReal       h, meanut;
542:   PetscScalar     hinv, neg_hinv, half = 0.5, sqrt_beta;
543:   PetscInt        im, indx1, indx2, indy1, indy2, indz1, indz2, nx, ny, nz;
544:   IS              is_alldesign, is_allstate;
545:   IS              is_from_d;
546:   IS              is_from_y;
547:   PetscInt        lo, hi, hi2, lo2, ysubnlocal, dsubnlocal;
548:   const PetscInt *ranges, *subranges;
549:   PetscMPIInt     size;
550:   PetscReal       xri, yri, zri, xim, yim, zim, dx1, dx2, dy1, dy2, dz1, dz2, Dx, Dy, Dz;
551:   PetscScalar     v, vx, vy, vz;
552:   PetscInt        offset, subindex, subvec, nrank, kk;

554:   PetscScalar xr[64] = {0.4970, 0.8498, 0.7814, 0.6268, 0.7782, 0.6402, 0.3617, 0.3160, 0.3610, 0.5298, 0.6987, 0.3331, 0.7962, 0.5596, 0.3866, 0.6774, 0.5407, 0.4518, 0.6702, 0.6061, 0.7580, 0.8997,
555:                         0.5198, 0.8326, 0.2138, 0.9198, 0.3000, 0.2833, 0.8288, 0.7076, 0.1820, 0.0728, 0.8447, 0.2367, 0.3239, 0.6413, 0.3114, 0.4731, 0.1192, 0.9273, 0.5724, 0.4331, 0.5136, 0.3547,
556:                         0.4413, 0.2602, 0.5698, 0.7278, 0.5261, 0.6230, 0.2454, 0.3948, 0.7479, 0.6582, 0.4660, 0.5594, 0.7574, 0.1143, 0.5900, 0.1065, 0.4260, 0.3294, 0.8276, 0.0756};

558:   PetscScalar yr[64] = {0.7345, 0.9120, 0.9288, 0.7528, 0.4463, 0.4985, 0.2497, 0.6256, 0.3425, 0.9026, 0.6983, 0.4230, 0.7140, 0.2970, 0.4474, 0.8792, 0.6604, 0.2485, 0.7968, 0.6127, 0.1796, 0.2437,
559:                         0.5938, 0.6137, 0.3867, 0.5658, 0.4575, 0.1009, 0.0863, 0.3361, 0.0738, 0.3985, 0.6602, 0.1437, 0.0934, 0.5983, 0.5950, 0.0763, 0.0768, 0.2288, 0.5761, 0.1129, 0.3841, 0.6150,
560:                         0.6904, 0.6686, 0.1361, 0.4601, 0.4491, 0.3716, 0.1969, 0.6537, 0.6743, 0.6991, 0.4811, 0.5480, 0.1684, 0.4569, 0.6889, 0.8437, 0.3015, 0.2854, 0.8199, 0.2658};

562:   PetscScalar zr[64] = {0.7668, 0.8573, 0.2654, 0.2719, 0.1060, 0.1311, 0.6232, 0.2295, 0.8009, 0.2147, 0.2119, 0.9325, 0.4473, 0.3600, 0.3374, 0.3819, 0.4066, 0.5801, 0.1673, 0.0959, 0.4638, 0.8236,
563:                         0.8800, 0.2939, 0.2028, 0.8262, 0.2706, 0.6276, 0.9085, 0.6443, 0.8241, 0.0712, 0.1824, 0.7789, 0.4389, 0.8415, 0.7055, 0.6639, 0.3653, 0.2078, 0.1987, 0.2297, 0.4321, 0.8115,
564:                         0.4915, 0.7764, 0.4657, 0.4627, 0.4569, 0.4232, 0.8514, 0.0674, 0.3227, 0.1055, 0.6690, 0.6313, 0.9226, 0.5461, 0.4126, 0.2364, 0.6096, 0.7042, 0.3914, 0.0711};

566:   PetscFunctionBegin;
567:   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
568:   PetscCall(PetscLogStageRegister("Elliptic Setup", &user->stages[0]));
569:   PetscCall(PetscLogStagePush(user->stages[0]));

571:   /* Create u,y,c,x */
572:   PetscCall(VecCreate(PETSC_COMM_WORLD, &user->u));
573:   PetscCall(VecCreate(PETSC_COMM_WORLD, &user->y));
574:   PetscCall(VecCreate(PETSC_COMM_WORLD, &user->c));
575:   PetscCall(VecSetSizes(user->u, PETSC_DECIDE, user->ndesign));
576:   PetscCall(VecSetFromOptions(user->u));
577:   PetscCall(VecGetLocalSize(user->u, &ysubnlocal));
578:   PetscCall(VecSetSizes(user->y, ysubnlocal * user->ns, user->nstate));
579:   PetscCall(VecSetSizes(user->c, ysubnlocal * user->ns, user->m));
580:   PetscCall(VecSetFromOptions(user->y));
581:   PetscCall(VecSetFromOptions(user->c));

583:   /*
584:      *******************************
585:      Create scatters for x <-> y,u
586:      *******************************

588:      If the state vector y and design vector u are partitioned as
589:      [y_1; y_2; ...; y_np] and [u_1; u_2; ...; u_np] (with np = # of processors),
590:      then the solution vector x is organized as
591:      [y_1; u_1; y_2; u_2; ...; y_np; u_np].
592:      The index sets user->s_is and user->d_is correspond to the indices of the
593:      state and design variables owned by the current processor.
594:   */
595:   PetscCall(VecCreate(PETSC_COMM_WORLD, &user->x));

597:   PetscCall(VecGetOwnershipRange(user->y, &lo, &hi));
598:   PetscCall(VecGetOwnershipRange(user->u, &lo2, &hi2));

600:   PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo, 1, &is_allstate));
601:   PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo + lo2, 1, &user->s_is));
602:   PetscCall(ISCreateStride(PETSC_COMM_SELF, hi2 - lo2, lo2, 1, &is_alldesign));
603:   PetscCall(ISCreateStride(PETSC_COMM_SELF, hi2 - lo2, hi + lo2, 1, &user->d_is));

605:   PetscCall(VecSetSizes(user->x, hi - lo + hi2 - lo2, user->n));
606:   PetscCall(VecSetFromOptions(user->x));

608:   PetscCall(VecScatterCreate(user->x, user->s_is, user->y, is_allstate, &user->state_scatter));
609:   PetscCall(VecScatterCreate(user->x, user->d_is, user->u, is_alldesign, &user->design_scatter));
610:   PetscCall(ISDestroy(&is_alldesign));
611:   PetscCall(ISDestroy(&is_allstate));
612:   /*
613:      *******************************
614:      Create scatter from y to y_1,y_2,...,y_ns
615:      *******************************
616:   */
617:   PetscCall(PetscMalloc1(user->ns, &user->yi_scatter));
618:   PetscCall(VecDuplicate(user->u, &user->suby));
619:   PetscCall(VecDuplicate(user->u, &user->subq));

621:   PetscCall(VecGetOwnershipRange(user->y, &lo2, &hi2));
622:   istart = 0;
623:   for (i = 0; i < user->ns; i++) {
624:     PetscCall(VecGetOwnershipRange(user->suby, &lo, &hi));
625:     PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo2 + istart, 1, &is_from_y));
626:     PetscCall(VecScatterCreate(user->y, is_from_y, user->suby, NULL, &user->yi_scatter[i]));
627:     istart = istart + hi - lo;
628:     PetscCall(ISDestroy(&is_from_y));
629:   }
630:   /*
631:      *******************************
632:      Create scatter from d to d_1,d_2,...,d_ns
633:      *******************************
634:   */
635:   PetscCall(VecCreate(PETSC_COMM_WORLD, &user->subd));
636:   PetscCall(VecSetSizes(user->subd, PETSC_DECIDE, user->ndata));
637:   PetscCall(VecSetFromOptions(user->subd));
638:   PetscCall(VecCreate(PETSC_COMM_WORLD, &user->d));
639:   PetscCall(VecGetLocalSize(user->subd, &dsubnlocal));
640:   PetscCall(VecSetSizes(user->d, dsubnlocal * user->ns, user->ndata * user->ns));
641:   PetscCall(VecSetFromOptions(user->d));
642:   PetscCall(PetscMalloc1(user->ns, &user->di_scatter));

644:   PetscCall(VecGetOwnershipRange(user->d, &lo2, &hi2));
645:   istart = 0;
646:   for (i = 0; i < user->ns; i++) {
647:     PetscCall(VecGetOwnershipRange(user->subd, &lo, &hi));
648:     PetscCall(ISCreateStride(PETSC_COMM_SELF, hi - lo, lo2 + istart, 1, &is_from_d));
649:     PetscCall(VecScatterCreate(user->d, is_from_d, user->subd, NULL, &user->di_scatter[i]));
650:     istart = istart + hi - lo;
651:     PetscCall(ISDestroy(&is_from_d));
652:   }

654:   PetscCall(PetscMalloc1(user->mx, &x));
655:   PetscCall(PetscMalloc1(user->mx, &y));
656:   PetscCall(PetscMalloc1(user->mx, &z));

658:   user->ksp_its         = 0;
659:   user->ksp_its_initial = 0;

661:   n         = user->mx * user->mx * user->mx;
662:   m         = 3 * user->mx * user->mx * (user->mx - 1);
663:   sqrt_beta = PetscSqrtScalar(user->beta);

665:   PetscCall(VecCreate(PETSC_COMM_WORLD, &XX));
666:   PetscCall(VecCreate(PETSC_COMM_WORLD, &user->q));
667:   PetscCall(VecSetSizes(XX, ysubnlocal, n));
668:   PetscCall(VecSetSizes(user->q, ysubnlocal * user->ns, user->m));
669:   PetscCall(VecSetFromOptions(XX));
670:   PetscCall(VecSetFromOptions(user->q));

672:   PetscCall(VecDuplicate(XX, &YY));
673:   PetscCall(VecDuplicate(XX, &ZZ));
674:   PetscCall(VecDuplicate(XX, &XXwork));
675:   PetscCall(VecDuplicate(XX, &YYwork));
676:   PetscCall(VecDuplicate(XX, &ZZwork));
677:   PetscCall(VecDuplicate(XX, &UTwork));
678:   PetscCall(VecDuplicate(XX, &user->utrue));

680:   /* map for striding q */
681:   PetscCall(VecGetOwnershipRanges(user->q, &ranges));
682:   PetscCall(VecGetOwnershipRanges(user->u, &subranges));

684:   PetscCall(VecGetOwnershipRange(user->q, &lo2, &hi2));
685:   PetscCall(VecGetOwnershipRange(user->u, &lo, &hi));
686:   /* Generate 3D grid, and collect ns (1<=ns<=8) right-hand-side vectors into user->q */
687:   h        = 1.0 / user->mx;
688:   hinv     = user->mx;
689:   neg_hinv = -hinv;

691:   PetscCall(VecGetOwnershipRange(XX, &istart, &iend));
692:   for (linear_index = istart; linear_index < iend; linear_index++) {
693:     i  = linear_index % user->mx;
694:     j  = ((linear_index - i) / user->mx) % user->mx;
695:     k  = ((linear_index - i) / user->mx - j) / user->mx;
696:     vx = h * (i + 0.5);
697:     vy = h * (j + 0.5);
698:     vz = h * (k + 0.5);
699:     PetscCall(VecSetValues(XX, 1, &linear_index, &vx, INSERT_VALUES));
700:     PetscCall(VecSetValues(YY, 1, &linear_index, &vy, INSERT_VALUES));
701:     PetscCall(VecSetValues(ZZ, 1, &linear_index, &vz, INSERT_VALUES));
702:     for (is = 0; is < 2; is++) {
703:       for (js = 0; js < 2; js++) {
704:         for (ks = 0; ks < 2; ks++) {
705:           ls = is * 4 + js * 2 + ks;
706:           if (ls < user->ns) {
707:             l = ls * n + linear_index;
708:             /* remap */
709:             subindex = l % n;
710:             subvec   = l / n;
711:             nrank    = 0;
712:             while (subindex >= subranges[nrank + 1]) nrank++;
713:             offset = subindex - subranges[nrank];
714:             istart = 0;
715:             for (kk = 0; kk < nrank; kk++) istart += user->ns * (subranges[kk + 1] - subranges[kk]);
716:             istart += (subranges[nrank + 1] - subranges[nrank]) * subvec;
717:             l = istart + offset;
718:             v = 100 * PetscSinScalar(2 * PETSC_PI * (vx + 0.25 * is)) * PetscSinScalar(2 * PETSC_PI * (vy + 0.25 * js)) * PetscSinScalar(2 * PETSC_PI * (vz + 0.25 * ks));
719:             PetscCall(VecSetValues(user->q, 1, &l, &v, INSERT_VALUES));
720:           }
721:         }
722:       }
723:     }
724:   }

726:   PetscCall(VecAssemblyBegin(XX));
727:   PetscCall(VecAssemblyEnd(XX));
728:   PetscCall(VecAssemblyBegin(YY));
729:   PetscCall(VecAssemblyEnd(YY));
730:   PetscCall(VecAssemblyBegin(ZZ));
731:   PetscCall(VecAssemblyEnd(ZZ));
732:   PetscCall(VecAssemblyBegin(user->q));
733:   PetscCall(VecAssemblyEnd(user->q));

735:   /* Compute true parameter function
736:      ut = exp(-((x-0.25)^2+(y-0.25)^2+(z-0.25)^2)/0.05) - exp((x-0.75)^2-(y-0.75)^2-(z-0.75))^2/0.05) */
737:   PetscCall(VecCopy(XX, XXwork));
738:   PetscCall(VecCopy(YY, YYwork));
739:   PetscCall(VecCopy(ZZ, ZZwork));

741:   PetscCall(VecShift(XXwork, -0.25));
742:   PetscCall(VecShift(YYwork, -0.25));
743:   PetscCall(VecShift(ZZwork, -0.25));

745:   PetscCall(VecPointwiseMult(XXwork, XXwork, XXwork));
746:   PetscCall(VecPointwiseMult(YYwork, YYwork, YYwork));
747:   PetscCall(VecPointwiseMult(ZZwork, ZZwork, ZZwork));

749:   PetscCall(VecCopy(XXwork, UTwork));
750:   PetscCall(VecAXPY(UTwork, 1.0, YYwork));
751:   PetscCall(VecAXPY(UTwork, 1.0, ZZwork));
752:   PetscCall(VecScale(UTwork, -20.0));
753:   PetscCall(VecExp(UTwork));
754:   PetscCall(VecCopy(UTwork, user->utrue));

756:   PetscCall(VecCopy(XX, XXwork));
757:   PetscCall(VecCopy(YY, YYwork));
758:   PetscCall(VecCopy(ZZ, ZZwork));

760:   PetscCall(VecShift(XXwork, -0.75));
761:   PetscCall(VecShift(YYwork, -0.75));
762:   PetscCall(VecShift(ZZwork, -0.75));

764:   PetscCall(VecPointwiseMult(XXwork, XXwork, XXwork));
765:   PetscCall(VecPointwiseMult(YYwork, YYwork, YYwork));
766:   PetscCall(VecPointwiseMult(ZZwork, ZZwork, ZZwork));

768:   PetscCall(VecCopy(XXwork, UTwork));
769:   PetscCall(VecAXPY(UTwork, 1.0, YYwork));
770:   PetscCall(VecAXPY(UTwork, 1.0, ZZwork));
771:   PetscCall(VecScale(UTwork, -20.0));
772:   PetscCall(VecExp(UTwork));

774:   PetscCall(VecAXPY(user->utrue, -1.0, UTwork));

776:   PetscCall(VecDestroy(&XX));
777:   PetscCall(VecDestroy(&YY));
778:   PetscCall(VecDestroy(&ZZ));
779:   PetscCall(VecDestroy(&XXwork));
780:   PetscCall(VecDestroy(&YYwork));
781:   PetscCall(VecDestroy(&ZZwork));
782:   PetscCall(VecDestroy(&UTwork));

784:   /* Initial guess and reference model */
785:   PetscCall(VecDuplicate(user->utrue, &user->ur));
786:   PetscCall(VecSum(user->utrue, &meanut));
787:   meanut = meanut / n;
788:   PetscCall(VecSet(user->ur, meanut));
789:   PetscCall(VecCopy(user->ur, user->u));

791:   /* Generate Grad matrix */
792:   PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Grad));
793:   PetscCall(MatSetSizes(user->Grad, PETSC_DECIDE, ysubnlocal, m, n));
794:   PetscCall(MatSetFromOptions(user->Grad));
795:   PetscCall(MatMPIAIJSetPreallocation(user->Grad, 2, NULL, 2, NULL));
796:   PetscCall(MatSeqAIJSetPreallocation(user->Grad, 2, NULL));
797:   PetscCall(MatGetOwnershipRange(user->Grad, &istart, &iend));

799:   for (i = istart; i < iend; i++) {
800:     if (i < m / 3) {
801:       iblock = i / (user->mx - 1);
802:       j      = iblock * user->mx + (i % (user->mx - 1));
803:       PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &neg_hinv, INSERT_VALUES));
804:       j = j + 1;
805:       PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &hinv, INSERT_VALUES));
806:     }
807:     if (i >= m / 3 && i < 2 * m / 3) {
808:       iblock = (i - m / 3) / (user->mx * (user->mx - 1));
809:       j      = iblock * user->mx * user->mx + ((i - m / 3) % (user->mx * (user->mx - 1)));
810:       PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &neg_hinv, INSERT_VALUES));
811:       j = j + user->mx;
812:       PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &hinv, INSERT_VALUES));
813:     }
814:     if (i >= 2 * m / 3) {
815:       j = i - 2 * m / 3;
816:       PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &neg_hinv, INSERT_VALUES));
817:       j = j + user->mx * user->mx;
818:       PetscCall(MatSetValues(user->Grad, 1, &i, 1, &j, &hinv, INSERT_VALUES));
819:     }
820:   }

822:   PetscCall(MatAssemblyBegin(user->Grad, MAT_FINAL_ASSEMBLY));
823:   PetscCall(MatAssemblyEnd(user->Grad, MAT_FINAL_ASSEMBLY));

825:   /* Generate arithmetic averaging matrix Av */
826:   PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Av));
827:   PetscCall(MatSetSizes(user->Av, PETSC_DECIDE, ysubnlocal, m, n));
828:   PetscCall(MatSetFromOptions(user->Av));
829:   PetscCall(MatMPIAIJSetPreallocation(user->Av, 2, NULL, 2, NULL));
830:   PetscCall(MatSeqAIJSetPreallocation(user->Av, 2, NULL));
831:   PetscCall(MatGetOwnershipRange(user->Av, &istart, &iend));

833:   for (i = istart; i < iend; i++) {
834:     if (i < m / 3) {
835:       iblock = i / (user->mx - 1);
836:       j      = iblock * user->mx + (i % (user->mx - 1));
837:       PetscCall(MatSetValues(user->Av, 1, &i, 1, &j, &half, INSERT_VALUES));
838:       j = j + 1;
839:       PetscCall(MatSetValues(user->Av, 1, &i, 1, &j, &half, INSERT_VALUES));
840:     }
841:     if (i >= m / 3 && i < 2 * m / 3) {
842:       iblock = (i - m / 3) / (user->mx * (user->mx - 1));
843:       j      = iblock * user->mx * user->mx + ((i - m / 3) % (user->mx * (user->mx - 1)));
844:       PetscCall(MatSetValues(user->Av, 1, &i, 1, &j, &half, INSERT_VALUES));
845:       j = j + user->mx;
846:       PetscCall(MatSetValues(user->Av, 1, &i, 1, &j, &half, INSERT_VALUES));
847:     }
848:     if (i >= 2 * m / 3) {
849:       j = i - 2 * m / 3;
850:       PetscCall(MatSetValues(user->Av, 1, &i, 1, &j, &half, INSERT_VALUES));
851:       j = j + user->mx * user->mx;
852:       PetscCall(MatSetValues(user->Av, 1, &i, 1, &j, &half, INSERT_VALUES));
853:     }
854:   }

856:   PetscCall(MatAssemblyBegin(user->Av, MAT_FINAL_ASSEMBLY));
857:   PetscCall(MatAssemblyEnd(user->Av, MAT_FINAL_ASSEMBLY));

859:   PetscCall(MatCreate(PETSC_COMM_WORLD, &user->L));
860:   PetscCall(MatSetSizes(user->L, PETSC_DECIDE, ysubnlocal, m + n, n));
861:   PetscCall(MatSetFromOptions(user->L));
862:   PetscCall(MatMPIAIJSetPreallocation(user->L, 2, NULL, 2, NULL));
863:   PetscCall(MatSeqAIJSetPreallocation(user->L, 2, NULL));
864:   PetscCall(MatGetOwnershipRange(user->L, &istart, &iend));

866:   for (i = istart; i < iend; i++) {
867:     if (i < m / 3) {
868:       iblock = i / (user->mx - 1);
869:       j      = iblock * user->mx + (i % (user->mx - 1));
870:       PetscCall(MatSetValues(user->L, 1, &i, 1, &j, &neg_hinv, INSERT_VALUES));
871:       j = j + 1;
872:       PetscCall(MatSetValues(user->L, 1, &i, 1, &j, &hinv, INSERT_VALUES));
873:     }
874:     if (i >= m / 3 && i < 2 * m / 3) {
875:       iblock = (i - m / 3) / (user->mx * (user->mx - 1));
876:       j      = iblock * user->mx * user->mx + ((i - m / 3) % (user->mx * (user->mx - 1)));
877:       PetscCall(MatSetValues(user->L, 1, &i, 1, &j, &neg_hinv, INSERT_VALUES));
878:       j = j + user->mx;
879:       PetscCall(MatSetValues(user->L, 1, &i, 1, &j, &hinv, INSERT_VALUES));
880:     }
881:     if (i >= 2 * m / 3 && i < m) {
882:       j = i - 2 * m / 3;
883:       PetscCall(MatSetValues(user->L, 1, &i, 1, &j, &neg_hinv, INSERT_VALUES));
884:       j = j + user->mx * user->mx;
885:       PetscCall(MatSetValues(user->L, 1, &i, 1, &j, &hinv, INSERT_VALUES));
886:     }
887:     if (i >= m) {
888:       j = i - m;
889:       PetscCall(MatSetValues(user->L, 1, &i, 1, &j, &sqrt_beta, INSERT_VALUES));
890:     }
891:   }
892:   PetscCall(MatAssemblyBegin(user->L, MAT_FINAL_ASSEMBLY));
893:   PetscCall(MatAssemblyEnd(user->L, MAT_FINAL_ASSEMBLY));
894:   PetscCall(MatScale(user->L, PetscPowScalar(h, 1.5)));

896:   /* Generate Div matrix */
897:   if (!user->use_ptap) {
898:     /* Generate Div matrix */
899:     PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Div));
900:     PetscCall(MatSetSizes(user->Div, ysubnlocal, PETSC_DECIDE, n, m));
901:     PetscCall(MatSetFromOptions(user->Div));
902:     PetscCall(MatMPIAIJSetPreallocation(user->Div, 4, NULL, 4, NULL));
903:     PetscCall(MatSeqAIJSetPreallocation(user->Div, 6, NULL));
904:     PetscCall(MatGetOwnershipRange(user->Grad, &istart, &iend));

906:     for (i = istart; i < iend; i++) {
907:       if (i < m / 3) {
908:         iblock = i / (user->mx - 1);
909:         j      = iblock * user->mx + (i % (user->mx - 1));
910:         PetscCall(MatSetValues(user->Div, 1, &j, 1, &i, &neg_hinv, INSERT_VALUES));
911:         j = j + 1;
912:         PetscCall(MatSetValues(user->Div, 1, &j, 1, &i, &hinv, INSERT_VALUES));
913:       }
914:       if (i >= m / 3 && i < 2 * m / 3) {
915:         iblock = (i - m / 3) / (user->mx * (user->mx - 1));
916:         j      = iblock * user->mx * user->mx + ((i - m / 3) % (user->mx * (user->mx - 1)));
917:         PetscCall(MatSetValues(user->Div, 1, &j, 1, &i, &neg_hinv, INSERT_VALUES));
918:         j = j + user->mx;
919:         PetscCall(MatSetValues(user->Div, 1, &j, 1, &i, &hinv, INSERT_VALUES));
920:       }
921:       if (i >= 2 * m / 3) {
922:         j = i - 2 * m / 3;
923:         PetscCall(MatSetValues(user->Div, 1, &j, 1, &i, &neg_hinv, INSERT_VALUES));
924:         j = j + user->mx * user->mx;
925:         PetscCall(MatSetValues(user->Div, 1, &j, 1, &i, &hinv, INSERT_VALUES));
926:       }
927:     }

929:     PetscCall(MatAssemblyBegin(user->Div, MAT_FINAL_ASSEMBLY));
930:     PetscCall(MatAssemblyEnd(user->Div, MAT_FINAL_ASSEMBLY));
931:     PetscCall(MatDuplicate(user->Div, MAT_SHARE_NONZERO_PATTERN, &user->Divwork));
932:   } else {
933:     PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Diag));
934:     PetscCall(MatSetSizes(user->Diag, PETSC_DECIDE, PETSC_DECIDE, m, m));
935:     PetscCall(MatSetFromOptions(user->Diag));
936:     PetscCall(MatMPIAIJSetPreallocation(user->Diag, 1, NULL, 0, NULL));
937:     PetscCall(MatSeqAIJSetPreallocation(user->Diag, 1, NULL));
938:   }

940:   /* Build work vectors and matrices */
941:   PetscCall(VecCreate(PETSC_COMM_WORLD, &user->S));
942:   PetscCall(VecSetSizes(user->S, PETSC_DECIDE, m));
943:   PetscCall(VecSetFromOptions(user->S));

945:   PetscCall(VecCreate(PETSC_COMM_WORLD, &user->lwork));
946:   PetscCall(VecSetSizes(user->lwork, PETSC_DECIDE, m + user->mx * user->mx * user->mx));
947:   PetscCall(VecSetFromOptions(user->lwork));

949:   PetscCall(MatDuplicate(user->Av, MAT_SHARE_NONZERO_PATTERN, &user->Avwork));

951:   PetscCall(VecDuplicate(user->S, &user->Swork));
952:   PetscCall(VecDuplicate(user->S, &user->Sdiag));
953:   PetscCall(VecDuplicate(user->S, &user->Av_u));
954:   PetscCall(VecDuplicate(user->S, &user->Twork));
955:   PetscCall(VecDuplicate(user->y, &user->ywork));
956:   PetscCall(VecDuplicate(user->u, &user->Ywork));
957:   PetscCall(VecDuplicate(user->u, &user->uwork));
958:   PetscCall(VecDuplicate(user->u, &user->js_diag));
959:   PetscCall(VecDuplicate(user->c, &user->cwork));
960:   PetscCall(VecDuplicate(user->d, &user->dwork));

962:   /* Create a matrix-free shell user->Jd for computing B*x */
963:   PetscCall(MatCreateShell(PETSC_COMM_WORLD, ysubnlocal * user->ns, ysubnlocal, user->nstate, user->ndesign, user, &user->Jd));
964:   PetscCall(MatShellSetOperation(user->Jd, MATOP_MULT, (PetscErrorCodeFn *)DesignMatMult));
965:   PetscCall(MatShellSetOperation(user->Jd, MATOP_MULT_TRANSPOSE, (PetscErrorCodeFn *)DesignMatMultTranspose));

967:   /* Compute true state function ytrue given utrue */
968:   PetscCall(VecDuplicate(user->y, &user->ytrue));

970:   /* First compute Av_u = Av*exp(-u) */
971:   PetscCall(VecSet(user->uwork, 0));
972:   PetscCall(VecAXPY(user->uwork, -1.0, user->utrue)); /* Note: user->utrue */
973:   PetscCall(VecExp(user->uwork));
974:   PetscCall(MatMult(user->Av, user->uwork, user->Av_u));

976:   /* Next form DSG = Div*S*Grad */
977:   PetscCall(VecCopy(user->Av_u, user->Swork));
978:   PetscCall(VecReciprocal(user->Swork));
979:   if (user->use_ptap) {
980:     PetscCall(MatDiagonalSet(user->Diag, user->Swork, INSERT_VALUES));
981:     PetscCall(MatPtAP(user->Diag, user->Grad, MAT_INITIAL_MATRIX, 1.0, &user->DSG));
982:   } else {
983:     PetscCall(MatCopy(user->Div, user->Divwork, SAME_NONZERO_PATTERN));
984:     PetscCall(MatDiagonalScale(user->Divwork, NULL, user->Swork));

986:     PetscCall(MatMatMult(user->Divwork, user->Grad, MAT_INITIAL_MATRIX, 1.0, &user->DSG));
987:   }

989:   PetscCall(MatSetOption(user->DSG, MAT_SYMMETRIC, PETSC_TRUE));
990:   PetscCall(MatSetOption(user->DSG, MAT_SYMMETRY_ETERNAL, PETSC_TRUE));

992:   if (user->use_lrc == PETSC_TRUE) {
993:     v = PetscSqrtReal(1.0 / user->ndesign);
994:     PetscCall(PetscMalloc1(user->ndesign, &user->ones));

996:     for (i = 0; i < user->ndesign; i++) user->ones[i] = v;
997:     PetscCall(MatCreateDense(PETSC_COMM_WORLD, ysubnlocal, PETSC_DECIDE, user->ndesign, 1, user->ones, &user->Ones));
998:     PetscCall(MatCreateLRC(user->DSG, user->Ones, NULL, user->Ones, &user->JsBlock));
999:     PetscCall(MatSetUp(user->JsBlock));
1000:   } else {
1001:     /* Create matrix-free shell user->Js for computing (A + h^3*e*e^T)*x */
1002:     PetscCall(MatCreateShell(PETSC_COMM_WORLD, ysubnlocal, ysubnlocal, user->ndesign, user->ndesign, user, &user->JsBlock));
1003:     PetscCall(MatShellSetOperation(user->JsBlock, MATOP_MULT, (PetscErrorCodeFn *)StateBlockMatMult));
1004:     PetscCall(MatShellSetOperation(user->JsBlock, MATOP_MULT_TRANSPOSE, (PetscErrorCodeFn *)StateBlockMatMult));
1005:   }
1006:   PetscCall(MatSetOption(user->JsBlock, MAT_SYMMETRIC, PETSC_TRUE));
1007:   PetscCall(MatSetOption(user->JsBlock, MAT_SYMMETRY_ETERNAL, PETSC_TRUE));
1008:   PetscCall(MatCreateShell(PETSC_COMM_WORLD, ysubnlocal * user->ns, ysubnlocal * user->ns, user->nstate, user->nstate, user, &user->Js));
1009:   PetscCall(MatShellSetOperation(user->Js, MATOP_MULT, (PetscErrorCodeFn *)StateMatMult));
1010:   PetscCall(MatShellSetOperation(user->Js, MATOP_MULT_TRANSPOSE, (PetscErrorCodeFn *)StateMatMult));
1011:   PetscCall(MatSetOption(user->Js, MAT_SYMMETRIC, PETSC_TRUE));
1012:   PetscCall(MatSetOption(user->Js, MAT_SYMMETRY_ETERNAL, PETSC_TRUE));

1014:   PetscCall(MatCreateShell(PETSC_COMM_WORLD, ysubnlocal * user->ns, ysubnlocal * user->ns, user->nstate, user->nstate, user, &user->JsInv));
1015:   PetscCall(MatShellSetOperation(user->JsInv, MATOP_MULT, (PetscErrorCodeFn *)StateInvMatMult));
1016:   PetscCall(MatShellSetOperation(user->JsInv, MATOP_MULT_TRANSPOSE, (PetscErrorCodeFn *)StateInvMatMult));
1017:   PetscCall(MatSetOption(user->JsInv, MAT_SYMMETRIC, PETSC_TRUE));
1018:   PetscCall(MatSetOption(user->JsInv, MAT_SYMMETRY_ETERNAL, PETSC_TRUE));

1020:   PetscCall(MatSetOption(user->DSG, MAT_SYMMETRIC, PETSC_TRUE));
1021:   PetscCall(MatSetOption(user->DSG, MAT_SYMMETRY_ETERNAL, PETSC_TRUE));
1022:   /* Now solve for ytrue */
1023:   PetscCall(KSPCreate(PETSC_COMM_WORLD, &user->solver));
1024:   PetscCall(KSPSetFromOptions(user->solver));

1026:   PetscCall(KSPSetOperators(user->solver, user->JsBlock, user->DSG));

1028:   PetscCall(MatMult(user->JsInv, user->q, user->ytrue));
1029:   /* First compute Av_u = Av*exp(-u) */
1030:   PetscCall(VecSet(user->uwork, 0));
1031:   PetscCall(VecAXPY(user->uwork, -1.0, user->u)); /* Note: user->u */
1032:   PetscCall(VecExp(user->uwork));
1033:   PetscCall(MatMult(user->Av, user->uwork, user->Av_u));

1035:   /* Next update DSG = Div*S*Grad  with user->u */
1036:   PetscCall(VecCopy(user->Av_u, user->Swork));
1037:   PetscCall(VecReciprocal(user->Swork));
1038:   if (user->use_ptap) {
1039:     PetscCall(MatDiagonalSet(user->Diag, user->Swork, INSERT_VALUES));
1040:     PetscCall(MatPtAP(user->Diag, user->Grad, MAT_REUSE_MATRIX, 1.0, &user->DSG));
1041:   } else {
1042:     PetscCall(MatCopy(user->Div, user->Divwork, SAME_NONZERO_PATTERN));
1043:     PetscCall(MatDiagonalScale(user->Divwork, NULL, user->Av_u));
1044:     PetscCall(MatProductNumeric(user->DSG));
1045:   }

1047:   /* Now solve for y */

1049:   PetscCall(MatMult(user->JsInv, user->q, user->y));

1051:   user->ksp_its_initial = user->ksp_its;
1052:   user->ksp_its         = 0;
1053:   /* Construct projection matrix Q (blocks) */
1054:   PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Q));
1055:   PetscCall(MatSetSizes(user->Q, dsubnlocal, ysubnlocal, user->ndata, user->ndesign));
1056:   PetscCall(MatSetFromOptions(user->Q));
1057:   PetscCall(MatMPIAIJSetPreallocation(user->Q, 8, NULL, 8, NULL));
1058:   PetscCall(MatSeqAIJSetPreallocation(user->Q, 8, NULL));

1060:   for (i = 0; i < user->mx; i++) {
1061:     x[i] = h * (i + 0.5);
1062:     y[i] = h * (i + 0.5);
1063:     z[i] = h * (i + 0.5);
1064:   }
1065:   PetscCall(MatGetOwnershipRange(user->Q, &istart, &iend));

1067:   nx = user->mx;
1068:   ny = user->mx;
1069:   nz = user->mx;
1070:   for (i = istart; i < iend; i++) {
1071:     xri = xr[i];
1072:     im  = 0;
1073:     xim = x[im];
1074:     while (xri > xim && im < nx) {
1075:       im  = im + 1;
1076:       xim = x[im];
1077:     }
1078:     indx1 = im - 1;
1079:     indx2 = im;
1080:     dx1   = xri - x[indx1];
1081:     dx2   = x[indx2] - xri;

1083:     yri = yr[i];
1084:     im  = 0;
1085:     yim = y[im];
1086:     while (yri > yim && im < ny) {
1087:       im  = im + 1;
1088:       yim = y[im];
1089:     }
1090:     indy1 = im - 1;
1091:     indy2 = im;
1092:     dy1   = yri - y[indy1];
1093:     dy2   = y[indy2] - yri;

1095:     zri = zr[i];
1096:     im  = 0;
1097:     zim = z[im];
1098:     while (zri > zim && im < nz) {
1099:       im  = im + 1;
1100:       zim = z[im];
1101:     }
1102:     indz1 = im - 1;
1103:     indz2 = im;
1104:     dz1   = zri - z[indz1];
1105:     dz2   = z[indz2] - zri;

1107:     Dx = x[indx2] - x[indx1];
1108:     Dy = y[indy2] - y[indy1];
1109:     Dz = z[indz2] - z[indz1];

1111:     j = indx1 + indy1 * nx + indz1 * nx * ny;
1112:     v = (1 - dx1 / Dx) * (1 - dy1 / Dy) * (1 - dz1 / Dz);
1113:     PetscCall(MatSetValues(user->Q, 1, &i, 1, &j, &v, INSERT_VALUES));

1115:     j = indx1 + indy1 * nx + indz2 * nx * ny;
1116:     v = (1 - dx1 / Dx) * (1 - dy1 / Dy) * (1 - dz2 / Dz);
1117:     PetscCall(MatSetValues(user->Q, 1, &i, 1, &j, &v, INSERT_VALUES));

1119:     j = indx1 + indy2 * nx + indz1 * nx * ny;
1120:     v = (1 - dx1 / Dx) * (1 - dy2 / Dy) * (1 - dz1 / Dz);
1121:     PetscCall(MatSetValues(user->Q, 1, &i, 1, &j, &v, INSERT_VALUES));

1123:     j = indx1 + indy2 * nx + indz2 * nx * ny;
1124:     v = (1 - dx1 / Dx) * (1 - dy2 / Dy) * (1 - dz2 / Dz);
1125:     PetscCall(MatSetValues(user->Q, 1, &i, 1, &j, &v, INSERT_VALUES));

1127:     j = indx2 + indy1 * nx + indz1 * nx * ny;
1128:     v = (1 - dx2 / Dx) * (1 - dy1 / Dy) * (1 - dz1 / Dz);
1129:     PetscCall(MatSetValues(user->Q, 1, &i, 1, &j, &v, INSERT_VALUES));

1131:     j = indx2 + indy1 * nx + indz2 * nx * ny;
1132:     v = (1 - dx2 / Dx) * (1 - dy1 / Dy) * (1 - dz2 / Dz);
1133:     PetscCall(MatSetValues(user->Q, 1, &i, 1, &j, &v, INSERT_VALUES));

1135:     j = indx2 + indy2 * nx + indz1 * nx * ny;
1136:     v = (1 - dx2 / Dx) * (1 - dy2 / Dy) * (1 - dz1 / Dz);
1137:     PetscCall(MatSetValues(user->Q, 1, &i, 1, &j, &v, INSERT_VALUES));

1139:     j = indx2 + indy2 * nx + indz2 * nx * ny;
1140:     v = (1 - dx2 / Dx) * (1 - dy2 / Dy) * (1 - dz2 / Dz);
1141:     PetscCall(MatSetValues(user->Q, 1, &i, 1, &j, &v, INSERT_VALUES));
1142:   }

1144:   PetscCall(MatAssemblyBegin(user->Q, MAT_FINAL_ASSEMBLY));
1145:   PetscCall(MatAssemblyEnd(user->Q, MAT_FINAL_ASSEMBLY));
1146:   /* Create MQ (composed of blocks of Q */
1147:   PetscCall(MatCreateShell(PETSC_COMM_WORLD, dsubnlocal * user->ns, PETSC_DECIDE, user->ndata * user->ns, user->nstate, user, &user->MQ));
1148:   PetscCall(MatShellSetOperation(user->MQ, MATOP_MULT, (PetscErrorCodeFn *)QMatMult));
1149:   PetscCall(MatShellSetOperation(user->MQ, MATOP_MULT_TRANSPOSE, (PetscErrorCodeFn *)QMatMultTranspose));

1151:   /* Add noise to the measurement data */
1152:   PetscCall(VecSet(user->ywork, 1.0));
1153:   PetscCall(VecAYPX(user->ywork, user->noise, user->ytrue));
1154:   PetscCall(MatMult(user->MQ, user->ywork, user->d));

1156:   /* Now that initial conditions have been set, let the user pass tolerance options to the KSP solver */
1157:   PetscCall(PetscFree(x));
1158:   PetscCall(PetscFree(y));
1159:   PetscCall(PetscFree(z));
1160:   PetscCall(PetscLogStagePop());
1161:   PetscFunctionReturn(PETSC_SUCCESS);
1162: }

1164: PetscErrorCode EllipticDestroy(AppCtx *user)
1165: {
1166:   PetscFunctionBegin;
1167:   PetscCall(MatDestroy(&user->DSG));
1168:   PetscCall(KSPDestroy(&user->solver));
1169:   PetscCall(MatDestroy(&user->Q));
1170:   PetscCall(MatDestroy(&user->MQ));
1171:   if (!user->use_ptap) {
1172:     PetscCall(MatDestroy(&user->Div));
1173:     PetscCall(MatDestroy(&user->Divwork));
1174:   } else {
1175:     PetscCall(MatDestroy(&user->Diag));
1176:   }
1177:   if (user->use_lrc) PetscCall(MatDestroy(&user->Ones));

1179:   PetscCall(MatDestroy(&user->Grad));
1180:   PetscCall(MatDestroy(&user->Av));
1181:   PetscCall(MatDestroy(&user->Avwork));
1182:   PetscCall(MatDestroy(&user->L));
1183:   PetscCall(MatDestroy(&user->Js));
1184:   PetscCall(MatDestroy(&user->Jd));
1185:   PetscCall(MatDestroy(&user->JsBlock));
1186:   PetscCall(MatDestroy(&user->JsInv));

1188:   PetscCall(VecDestroy(&user->x));
1189:   PetscCall(VecDestroy(&user->u));
1190:   PetscCall(VecDestroy(&user->uwork));
1191:   PetscCall(VecDestroy(&user->utrue));
1192:   PetscCall(VecDestroy(&user->y));
1193:   PetscCall(VecDestroy(&user->ywork));
1194:   PetscCall(VecDestroy(&user->ytrue));
1195:   PetscCall(VecDestroy(&user->c));
1196:   PetscCall(VecDestroy(&user->cwork));
1197:   PetscCall(VecDestroy(&user->ur));
1198:   PetscCall(VecDestroy(&user->q));
1199:   PetscCall(VecDestroy(&user->d));
1200:   PetscCall(VecDestroy(&user->dwork));
1201:   PetscCall(VecDestroy(&user->lwork));
1202:   PetscCall(VecDestroy(&user->S));
1203:   PetscCall(VecDestroy(&user->Swork));
1204:   PetscCall(VecDestroy(&user->Sdiag));
1205:   PetscCall(VecDestroy(&user->Ywork));
1206:   PetscCall(VecDestroy(&user->Twork));
1207:   PetscCall(VecDestroy(&user->Av_u));
1208:   PetscCall(VecDestroy(&user->js_diag));
1209:   PetscCall(ISDestroy(&user->s_is));
1210:   PetscCall(ISDestroy(&user->d_is));
1211:   PetscCall(VecDestroy(&user->suby));
1212:   PetscCall(VecDestroy(&user->subd));
1213:   PetscCall(VecDestroy(&user->subq));
1214:   PetscCall(VecScatterDestroy(&user->state_scatter));
1215:   PetscCall(VecScatterDestroy(&user->design_scatter));
1216:   for (PetscInt i = 0; i < user->ns; i++) {
1217:     PetscCall(VecScatterDestroy(&user->yi_scatter[i]));
1218:     PetscCall(VecScatterDestroy(&user->di_scatter[i]));
1219:   }
1220:   PetscCall(PetscFree(user->yi_scatter));
1221:   PetscCall(PetscFree(user->di_scatter));
1222:   if (user->use_lrc) {
1223:     PetscCall(PetscFree(user->ones));
1224:     PetscCall(MatDestroy(&user->Ones));
1225:   }
1226:   PetscFunctionReturn(PETSC_SUCCESS);
1227: }

1229: PetscErrorCode EllipticMonitor(Tao tao, void *ptr)
1230: {
1231:   Vec       X;
1232:   PetscReal unorm, ynorm;
1233:   AppCtx   *user = (AppCtx *)ptr;

1235:   PetscFunctionBegin;
1236:   PetscCall(TaoGetSolution(tao, &X));
1237:   PetscCall(Scatter(X, user->ywork, user->state_scatter, user->uwork, user->design_scatter));
1238:   PetscCall(VecAXPY(user->ywork, -1.0, user->ytrue));
1239:   PetscCall(VecAXPY(user->uwork, -1.0, user->utrue));
1240:   PetscCall(VecNorm(user->uwork, NORM_2, &unorm));
1241:   PetscCall(VecNorm(user->ywork, NORM_2, &ynorm));
1242:   PetscCall(PetscPrintf(MPI_COMM_WORLD, "||u-ut||=%g ||y-yt||=%g\n", (double)unorm, (double)ynorm));
1243:   PetscFunctionReturn(PETSC_SUCCESS);
1244: }

1246: /*TEST

1248:    build:
1249:       requires: !complex

1251:    test:
1252:       args: -tao_monitor_constraint_norm -ns 1 -tao_type lcl -tao_gatol 1.e-3 -tao_max_it 11
1253:       requires: !single

1255:    test:
1256:       suffix: 2
1257:       args: -tao_monitor_constraint_norm -tao_type lcl -tao_max_it 11 -use_ptap -use_lrc -ns 1 -tao_gatol 1.e-3
1258:       requires: !single

1260: TEST*/