Actual source code: parabolic.c
1: #include <petsc/private/taoimpl.h>
3: typedef struct {
4: PetscInt n; /* Number of variables */
5: PetscInt m; /* Number of constraints per time step */
6: PetscInt mx; /* grid points in each direction */
7: PetscInt nt; /* Number of time steps; as of now, must be divisible by 8 */
8: PetscInt ndata; /* Number of data points per sample */
9: PetscInt ns; /* Number of samples */
10: PetscInt *sample_times; /* Times of samples */
11: IS s_is;
12: IS d_is;
14: VecScatter state_scatter;
15: VecScatter design_scatter;
16: VecScatter *yi_scatter;
17: VecScatter *di_scatter;
19: Mat Js,Jd,JsBlockPrec,JsInv,JsBlock;
20: PetscBool jformed,dsg_formed;
22: PetscReal alpha; /* Regularization parameter */
23: PetscReal beta; /* Weight attributed to ||u||^2 in regularization functional */
24: PetscReal noise; /* Amount of noise to add to data */
25: PetscReal ht; /* Time step */
27: Mat Qblock,QblockT;
28: Mat L,LT;
29: Mat Div,Divwork;
30: Mat Grad;
31: Mat Av,Avwork,AvT;
32: Mat DSG;
33: Vec q;
34: Vec ur; /* reference */
36: Vec d;
37: Vec dwork;
38: Vec *di;
40: Vec y; /* state variables */
41: Vec ywork;
43: Vec ytrue;
44: Vec *yi,*yiwork;
46: Vec u; /* design variables */
47: Vec uwork;
49: Vec utrue;
50: Vec js_diag;
51: Vec c; /* constraint vector */
52: Vec cwork;
54: Vec lwork;
55: Vec S;
56: Vec Rwork,Swork,Twork;
57: Vec Av_u;
59: KSP solver;
60: PC prec;
62: PetscInt ksp_its;
63: PetscInt ksp_its_initial;
64: } AppCtx;
66: PetscErrorCode FormFunction(Tao, Vec, PetscReal*, void*);
67: PetscErrorCode FormGradient(Tao, Vec, Vec, void*);
68: PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal*, Vec, void*);
69: PetscErrorCode FormJacobianState(Tao, Vec, Mat, Mat, Mat, void*);
70: PetscErrorCode FormJacobianDesign(Tao, Vec, Mat, void*);
71: PetscErrorCode FormConstraints(Tao, Vec, Vec, void*);
72: PetscErrorCode FormHessian(Tao, Vec, Mat, Mat, void*);
73: PetscErrorCode Gather(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat);
74: PetscErrorCode Scatter(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat);
75: PetscErrorCode ParabolicInitialize(AppCtx *user);
76: PetscErrorCode ParabolicDestroy(AppCtx *user);
77: PetscErrorCode ParabolicMonitor(Tao, void*);
79: PetscErrorCode StateMatMult(Mat,Vec,Vec);
80: PetscErrorCode StateMatBlockMult(Mat,Vec,Vec);
81: PetscErrorCode StateMatMultTranspose(Mat,Vec,Vec);
82: PetscErrorCode StateMatGetDiagonal(Mat,Vec);
83: PetscErrorCode StateMatDuplicate(Mat,MatDuplicateOption,Mat*);
84: PetscErrorCode StateMatInvMult(Mat,Vec,Vec);
85: PetscErrorCode StateMatInvTransposeMult(Mat,Vec,Vec);
86: PetscErrorCode StateMatBlockPrecMult(PC,Vec,Vec);
88: PetscErrorCode DesignMatMult(Mat,Vec,Vec);
89: PetscErrorCode DesignMatMultTranspose(Mat,Vec,Vec);
91: PetscErrorCode Gather_i(Vec,Vec*,VecScatter*,PetscInt);
92: PetscErrorCode Scatter_i(Vec,Vec*,VecScatter*,PetscInt);
94: static char help[]="";
96: int main(int argc, char **argv)
97: {
98: PetscErrorCode ierr;
99: Vec x,x0;
100: Tao tao;
101: AppCtx user;
102: IS is_allstate,is_alldesign;
103: PetscInt lo,hi,hi2,lo2,ksp_old;
104: PetscInt ntests = 1;
105: PetscInt i;
106: #if defined(PETSC_USE_LOG)
107: PetscLogStage stages[1];
108: #endif
110: PetscInitialize(&argc, &argv, (char*)0,help);
111: user.mx = 8;
112: PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"parabolic example",NULL);
113: PetscOptionsInt("-mx","Number of grid points in each direction","",user.mx,&user.mx,NULL);
114: user.nt = 8;
115: PetscOptionsInt("-nt","Number of time steps","",user.nt,&user.nt,NULL);
116: user.ndata = 64;
117: PetscOptionsInt("-ndata","Numbers of data points per sample","",user.ndata,&user.ndata,NULL);
118: user.ns = 8;
119: PetscOptionsInt("-ns","Number of samples","",user.ns,&user.ns,NULL);
120: user.alpha = 1.0;
121: PetscOptionsReal("-alpha","Regularization parameter","",user.alpha,&user.alpha,NULL);
122: user.beta = 0.01;
123: PetscOptionsReal("-beta","Weight attributed to ||u||^2 in regularization functional","",user.beta,&user.beta,NULL);
124: user.noise = 0.01;
125: PetscOptionsReal("-noise","Amount of noise to add to data","",user.noise,&user.noise,NULL);
126: PetscOptionsInt("-ntests","Number of times to repeat TaoSolve","",ntests,&ntests,NULL);
127: PetscOptionsEnd();
129: user.m = user.mx*user.mx*user.mx; /* number of constraints per time step */
130: user.n = user.m*(user.nt+1); /* number of variables */
131: user.ht = (PetscReal)1/user.nt;
133: VecCreate(PETSC_COMM_WORLD,&user.u);
134: VecCreate(PETSC_COMM_WORLD,&user.y);
135: VecCreate(PETSC_COMM_WORLD,&user.c);
136: VecSetSizes(user.u,PETSC_DECIDE,user.n-user.m*user.nt);
137: VecSetSizes(user.y,PETSC_DECIDE,user.m*user.nt);
138: VecSetSizes(user.c,PETSC_DECIDE,user.m*user.nt);
139: VecSetFromOptions(user.u);
140: VecSetFromOptions(user.y);
141: VecSetFromOptions(user.c);
143: /* Create scatters for reduced spaces.
144: If the state vector y and design vector u are partitioned as
145: [y_1; y_2; ...; y_np] and [u_1; u_2; ...; u_np] (with np = # of processors),
146: then the solution vector x is organized as
147: [y_1; u_1; y_2; u_2; ...; y_np; u_np].
148: The index sets user.s_is and user.d_is correspond to the indices of the
149: state and design variables owned by the current processor.
150: */
151: VecCreate(PETSC_COMM_WORLD,&x);
153: VecGetOwnershipRange(user.y,&lo,&hi);
154: VecGetOwnershipRange(user.u,&lo2,&hi2);
156: ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_allstate);
157: ISCreateStride(PETSC_COMM_SELF,hi-lo,lo+lo2,1,&user.s_is);
159: ISCreateStride(PETSC_COMM_SELF,hi2-lo2,lo2,1,&is_alldesign);
160: ISCreateStride(PETSC_COMM_SELF,hi2-lo2,hi+lo2,1,&user.d_is);
162: VecSetSizes(x,hi-lo+hi2-lo2,user.n);
163: VecSetFromOptions(x);
165: VecScatterCreate(x,user.s_is,user.y,is_allstate,&user.state_scatter);
166: VecScatterCreate(x,user.d_is,user.u,is_alldesign,&user.design_scatter);
167: ISDestroy(&is_alldesign);
168: ISDestroy(&is_allstate);
170: /* Create TAO solver and set desired solution method */
171: TaoCreate(PETSC_COMM_WORLD,&tao);
172: TaoSetType(tao,TAOLCL);
174: /* Set up initial vectors and matrices */
175: ParabolicInitialize(&user);
177: Gather(x,user.y,user.state_scatter,user.u,user.design_scatter);
178: VecDuplicate(x,&x0);
179: VecCopy(x,x0);
181: /* Set solution vector with an initial guess */
182: TaoSetSolution(tao,x);
183: TaoSetObjective(tao, FormFunction, &user);
184: TaoSetGradient(tao, NULL, FormGradient, &user);
185: TaoSetConstraintsRoutine(tao, user.c, FormConstraints, &user);
187: TaoSetJacobianStateRoutine(tao, user.Js, user.JsBlockPrec, user.JsInv, FormJacobianState, &user);
188: TaoSetJacobianDesignRoutine(tao, user.Jd, FormJacobianDesign, &user);
190: TaoSetFromOptions(tao);
191: TaoSetStateDesignIS(tao,user.s_is,user.d_is);
193: /* SOLVE THE APPLICATION */
194: PetscLogStageRegister("Trials",&stages[0]);
195: PetscLogStagePush(stages[0]);
196: user.ksp_its_initial = user.ksp_its;
197: ksp_old = user.ksp_its;
198: for (i=0; i<ntests; i++) {
199: TaoSolve(tao);
200: PetscPrintf(PETSC_COMM_WORLD,"KSP Iterations = %D\n",user.ksp_its-ksp_old);
201: VecCopy(x0,x);
202: TaoSetSolution(tao,x);
203: }
204: PetscLogStagePop();
205: PetscBarrier((PetscObject)x);
206: PetscPrintf(PETSC_COMM_WORLD,"KSP iterations within initialization: ");
207: PetscPrintf(PETSC_COMM_WORLD,"%D\n",user.ksp_its_initial);
208: PetscPrintf(PETSC_COMM_WORLD,"Total KSP iterations over %D trial(s): ",ntests);
209: PetscPrintf(PETSC_COMM_WORLD,"%D\n",user.ksp_its);
210: PetscPrintf(PETSC_COMM_WORLD,"KSP iterations per trial: ");
211: PetscPrintf(PETSC_COMM_WORLD,"%D\n",(user.ksp_its-user.ksp_its_initial)/ntests);
213: TaoDestroy(&tao);
214: VecDestroy(&x);
215: VecDestroy(&x0);
216: ParabolicDestroy(&user);
218: PetscFinalize();
219: return 0;
220: }
221: /* ------------------------------------------------------------------- */
222: /*
223: dwork = Qy - d
224: lwork = L*(u-ur)
225: f = 1/2 * (dwork.dork + alpha*lwork.lwork)
226: */
227: PetscErrorCode FormFunction(Tao tao,Vec X,PetscReal *f,void *ptr)
228: {
229: PetscReal d1=0,d2=0;
230: PetscInt i,j;
231: AppCtx *user = (AppCtx*)ptr;
233: Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
234: Scatter_i(user->y,user->yi,user->yi_scatter,user->nt);
235: for (j=0; j<user->ns; j++) {
236: i = user->sample_times[j];
237: MatMult(user->Qblock,user->yi[i],user->di[j]);
238: }
239: Gather_i(user->dwork,user->di,user->di_scatter,user->ns);
240: VecAXPY(user->dwork,-1.0,user->d);
241: VecDot(user->dwork,user->dwork,&d1);
243: VecWAXPY(user->uwork,-1.0,user->ur,user->u);
244: MatMult(user->L,user->uwork,user->lwork);
245: VecDot(user->lwork,user->lwork,&d2);
247: *f = 0.5 * (d1 + user->alpha*d2);
248: return 0;
249: }
251: /* ------------------------------------------------------------------- */
252: /*
253: state: g_s = Q' *(Qy - d)
254: design: g_d = alpha*L'*L*(u-ur)
255: */
256: PetscErrorCode FormGradient(Tao tao,Vec X,Vec G,void *ptr)
257: {
258: PetscInt i,j;
259: AppCtx *user = (AppCtx*)ptr;
261: Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
262: Scatter_i(user->y,user->yi,user->yi_scatter,user->nt);
263: for (j=0; j<user->ns; j++) {
264: i = user->sample_times[j];
265: MatMult(user->Qblock,user->yi[i],user->di[j]);
266: }
267: Gather_i(user->dwork,user->di,user->di_scatter,user->ns);
268: VecAXPY(user->dwork,-1.0,user->d);
269: Scatter_i(user->dwork,user->di,user->di_scatter,user->ns);
270: VecSet(user->ywork,0.0);
271: Scatter_i(user->ywork,user->yiwork,user->yi_scatter,user->nt);
272: for (j=0; j<user->ns; j++) {
273: i = user->sample_times[j];
274: MatMult(user->QblockT,user->di[j],user->yiwork[i]);
275: }
276: Gather_i(user->ywork,user->yiwork,user->yi_scatter,user->nt);
278: VecWAXPY(user->uwork,-1.0,user->ur,user->u);
279: MatMult(user->L,user->uwork,user->lwork);
280: MatMult(user->LT,user->lwork,user->uwork);
281: VecScale(user->uwork, user->alpha);
282: Gather(G,user->ywork,user->state_scatter,user->uwork,user->design_scatter);
283: return 0;
284: }
286: PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, void *ptr)
287: {
288: PetscReal d1,d2;
289: PetscInt i,j;
290: AppCtx *user = (AppCtx*)ptr;
292: Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
293: Scatter_i(user->y,user->yi,user->yi_scatter,user->nt);
294: for (j=0; j<user->ns; j++) {
295: i = user->sample_times[j];
296: MatMult(user->Qblock,user->yi[i],user->di[j]);
297: }
298: Gather_i(user->dwork,user->di,user->di_scatter,user->ns);
299: VecAXPY(user->dwork,-1.0,user->d);
300: VecDot(user->dwork,user->dwork,&d1);
301: Scatter_i(user->dwork,user->di,user->di_scatter,user->ns);
302: VecSet(user->ywork,0.0);
303: Scatter_i(user->ywork,user->yiwork,user->yi_scatter,user->nt);
304: for (j=0; j<user->ns; j++) {
305: i = user->sample_times[j];
306: MatMult(user->QblockT,user->di[j],user->yiwork[i]);
307: }
308: Gather_i(user->ywork,user->yiwork,user->yi_scatter,user->nt);
310: VecWAXPY(user->uwork,-1.0,user->ur,user->u);
311: MatMult(user->L,user->uwork,user->lwork);
312: VecDot(user->lwork,user->lwork,&d2);
313: MatMult(user->LT,user->lwork,user->uwork);
314: VecScale(user->uwork, user->alpha);
315: *f = 0.5 * (d1 + user->alpha*d2);
317: Gather(G,user->ywork,user->state_scatter,user->uwork,user->design_scatter);
318: return 0;
319: }
321: /* ------------------------------------------------------------------- */
322: /* A
323: MatShell object
324: */
325: PetscErrorCode FormJacobianState(Tao tao, Vec X, Mat J, Mat JPre, Mat JInv, void *ptr)
326: {
327: AppCtx *user = (AppCtx*)ptr;
329: Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
330: VecSet(user->uwork,0);
331: VecAXPY(user->uwork,-1.0,user->u);
332: VecExp(user->uwork);
333: MatMult(user->Av,user->uwork,user->Av_u);
334: VecCopy(user->Av_u,user->Swork);
335: VecReciprocal(user->Swork);
336: MatCopy(user->Div,user->Divwork,SAME_NONZERO_PATTERN);
337: MatDiagonalScale(user->Divwork,NULL,user->Swork);
338: if (user->dsg_formed) {
339: MatProductNumeric(user->DSG);
340: } else {
341: MatMatMult(user->Divwork,user->Grad,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&user->DSG);
342: user->dsg_formed = PETSC_TRUE;
343: }
345: /* B = speye(nx^3) + ht*DSG; */
346: MatScale(user->DSG,user->ht);
347: MatShift(user->DSG,1.0);
348: return 0;
349: }
351: /* ------------------------------------------------------------------- */
352: /* B */
353: PetscErrorCode FormJacobianDesign(Tao tao, Vec X, Mat J, void *ptr)
354: {
355: AppCtx *user = (AppCtx*)ptr;
357: Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
358: return 0;
359: }
361: PetscErrorCode StateMatMult(Mat J_shell, Vec X, Vec Y)
362: {
363: PetscInt i;
364: AppCtx *user;
366: MatShellGetContext(J_shell,&user);
367: Scatter_i(X,user->yi,user->yi_scatter,user->nt);
368: MatMult(user->JsBlock,user->yi[0],user->yiwork[0]);
369: for (i=1; i<user->nt; i++) {
370: MatMult(user->JsBlock,user->yi[i],user->yiwork[i]);
371: VecAXPY(user->yiwork[i],-1.0,user->yi[i-1]);
372: }
373: Gather_i(Y,user->yiwork,user->yi_scatter,user->nt);
374: return 0;
375: }
377: PetscErrorCode StateMatMultTranspose(Mat J_shell, Vec X, Vec Y)
378: {
379: PetscInt i;
380: AppCtx *user;
382: MatShellGetContext(J_shell,&user);
383: Scatter_i(X,user->yi,user->yi_scatter,user->nt);
384: for (i=0; i<user->nt-1; i++) {
385: MatMult(user->JsBlock,user->yi[i],user->yiwork[i]);
386: VecAXPY(user->yiwork[i],-1.0,user->yi[i+1]);
387: }
388: i = user->nt-1;
389: MatMult(user->JsBlock,user->yi[i],user->yiwork[i]);
390: Gather_i(Y,user->yiwork,user->yi_scatter,user->nt);
391: return 0;
392: }
394: PetscErrorCode StateMatBlockMult(Mat J_shell, Vec X, Vec Y)
395: {
396: AppCtx *user;
398: MatShellGetContext(J_shell,&user);
399: MatMult(user->Grad,X,user->Swork);
400: VecPointwiseDivide(user->Swork,user->Swork,user->Av_u);
401: MatMult(user->Div,user->Swork,Y);
402: VecAYPX(Y,user->ht,X);
403: return 0;
404: }
406: PetscErrorCode DesignMatMult(Mat J_shell, Vec X, Vec Y)
407: {
408: PetscInt i;
409: AppCtx *user;
411: MatShellGetContext(J_shell,&user);
413: /* sdiag(1./v) */
414: VecSet(user->uwork,0);
415: VecAXPY(user->uwork,-1.0,user->u);
416: VecExp(user->uwork);
418: /* sdiag(1./((Av*(1./v)).^2)) */
419: MatMult(user->Av,user->uwork,user->Swork);
420: VecPointwiseMult(user->Swork,user->Swork,user->Swork);
421: VecReciprocal(user->Swork);
423: /* (Av * (sdiag(1./v) * b)) */
424: VecPointwiseMult(user->uwork,user->uwork,X);
425: MatMult(user->Av,user->uwork,user->Twork);
427: /* (sdiag(1./((Av*(1./v)).^2)) * (Av * (sdiag(1./v) * b))) */
428: VecPointwiseMult(user->Swork,user->Twork,user->Swork);
430: Scatter_i(user->y,user->yi,user->yi_scatter,user->nt);
431: for (i=0; i<user->nt; i++) {
432: /* (sdiag(Grad*y(:,i)) */
433: MatMult(user->Grad,user->yi[i],user->Twork);
435: /* ht * Div * (sdiag(Grad*y(:,i)) * (sdiag(1./((Av*(1./v)).^2)) * (Av * (sdiag(1./v) * b)))) */
436: VecPointwiseMult(user->Twork,user->Twork,user->Swork);
437: MatMult(user->Div,user->Twork,user->yiwork[i]);
438: VecScale(user->yiwork[i],user->ht);
439: }
440: Gather_i(Y,user->yiwork,user->yi_scatter,user->nt);
442: return 0;
443: }
445: PetscErrorCode DesignMatMultTranspose(Mat J_shell, Vec X, Vec Y)
446: {
447: PetscInt i;
448: AppCtx *user;
450: MatShellGetContext(J_shell,&user);
452: /* sdiag(1./((Av*(1./v)).^2)) */
453: VecSet(user->uwork,0);
454: VecAXPY(user->uwork,-1.0,user->u);
455: VecExp(user->uwork);
456: MatMult(user->Av,user->uwork,user->Rwork);
457: VecPointwiseMult(user->Rwork,user->Rwork,user->Rwork);
458: VecReciprocal(user->Rwork);
460: VecSet(Y,0.0);
461: Scatter_i(user->y,user->yi,user->yi_scatter,user->nt);
462: Scatter_i(X,user->yiwork,user->yi_scatter,user->nt);
463: for (i=0; i<user->nt; i++) {
464: /* (Div' * b(:,i)) */
465: MatMult(user->Grad,user->yiwork[i],user->Swork);
467: /* sdiag(Grad*y(:,i)) */
468: MatMult(user->Grad,user->yi[i],user->Twork);
470: /* (sdiag(Grad*y(:,i)) * (Div' * b(:,i))) */
471: VecPointwiseMult(user->Twork,user->Swork,user->Twork);
473: /* (sdiag(1./((Av*(1./v)).^2)) * (sdiag(Grad*y(:,i)) * (Div' * b(:,i)))) */
474: VecPointwiseMult(user->Twork,user->Rwork,user->Twork);
476: /* (Av' * (sdiag(1./((Av*(1./v)).^2)) * (sdiag(Grad*y(:,i)) * (Div' * b(:,i))))) */
477: MatMult(user->AvT,user->Twork,user->yiwork[i]);
479: /* sdiag(1./v) * (Av' * (sdiag(1./((Av*(1./v)).^2)) * (sdiag(Grad*y(:,i)) * (Div' * b(:,i))))) */
480: VecPointwiseMult(user->yiwork[i],user->uwork,user->yiwork[i]);
481: VecAXPY(Y,user->ht,user->yiwork[i]);
482: }
483: return 0;
484: }
486: PetscErrorCode StateMatBlockPrecMult(PC PC_shell, Vec X, Vec Y)
487: {
488: AppCtx *user;
490: PCShellGetContext(PC_shell,&user);
492: if (user->dsg_formed) {
493: MatSOR(user->DSG,X,1.0,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_SYMMETRIC_SWEEP),0.0,1,1,Y);
494: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"DSG not formed");
495: return 0;
496: }
498: PetscErrorCode StateMatInvMult(Mat J_shell, Vec X, Vec Y)
499: {
500: AppCtx *user;
501: PetscInt its,i;
503: MatShellGetContext(J_shell,&user);
505: if (Y == user->ytrue) {
506: /* First solve is done with true solution to set up problem */
507: KSPSetTolerances(user->solver,1e-8,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
508: } else {
509: KSPSetTolerances(user->solver,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
510: }
512: Scatter_i(X,user->yi,user->yi_scatter,user->nt);
513: KSPSolve(user->solver,user->yi[0],user->yiwork[0]);
514: KSPGetIterationNumber(user->solver,&its);
515: user->ksp_its = user->ksp_its + its;
517: for (i=1; i<user->nt; i++) {
518: VecAXPY(user->yi[i],1.0,user->yiwork[i-1]);
519: KSPSolve(user->solver,user->yi[i],user->yiwork[i]);
520: KSPGetIterationNumber(user->solver,&its);
521: user->ksp_its = user->ksp_its + its;
522: }
523: Gather_i(Y,user->yiwork,user->yi_scatter,user->nt);
524: return 0;
525: }
527: PetscErrorCode StateMatInvTransposeMult(Mat J_shell, Vec X, Vec Y)
528: {
529: AppCtx *user;
530: PetscInt its,i;
532: MatShellGetContext(J_shell,&user);
534: Scatter_i(X,user->yi,user->yi_scatter,user->nt);
536: i = user->nt - 1;
537: KSPSolve(user->solver,user->yi[i],user->yiwork[i]);
539: KSPGetIterationNumber(user->solver,&its);
540: user->ksp_its = user->ksp_its + its;
542: for (i=user->nt-2; i>=0; i--) {
543: VecAXPY(user->yi[i],1.0,user->yiwork[i+1]);
544: KSPSolve(user->solver,user->yi[i],user->yiwork[i]);
546: KSPGetIterationNumber(user->solver,&its);
547: user->ksp_its = user->ksp_its + its;
549: }
551: Gather_i(Y,user->yiwork,user->yi_scatter,user->nt);
552: return 0;
553: }
555: PetscErrorCode StateMatDuplicate(Mat J_shell, MatDuplicateOption opt, Mat *new_shell)
556: {
557: AppCtx *user;
559: MatShellGetContext(J_shell,&user);
561: MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m,user->m,user,new_shell);
562: MatShellSetOperation(*new_shell,MATOP_MULT,(void(*)(void))StateMatMult);
563: MatShellSetOperation(*new_shell,MATOP_DUPLICATE,(void(*)(void))StateMatDuplicate);
564: MatShellSetOperation(*new_shell,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatMultTranspose);
565: MatShellSetOperation(*new_shell,MATOP_GET_DIAGONAL,(void(*)(void))StateMatGetDiagonal);
566: return 0;
567: }
569: PetscErrorCode StateMatGetDiagonal(Mat J_shell, Vec X)
570: {
571: AppCtx *user;
573: MatShellGetContext(J_shell,&user);
574: VecCopy(user->js_diag,X);
575: return 0;
577: }
579: PetscErrorCode FormConstraints(Tao tao, Vec X, Vec C, void *ptr)
580: {
581: /* con = Ay - q, A = [B 0 0 ... 0;
582: -I B 0 ... 0;
583: 0 -I B ... 0;
584: ... ;
585: 0 ... -I B]
586: B = ht * Div * Sigma * Grad + eye */
587: PetscInt i;
588: AppCtx *user = (AppCtx*)ptr;
590: Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
591: Scatter_i(user->y,user->yi,user->yi_scatter,user->nt);
592: MatMult(user->JsBlock,user->yi[0],user->yiwork[0]);
593: for (i=1; i<user->nt; i++) {
594: MatMult(user->JsBlock,user->yi[i],user->yiwork[i]);
595: VecAXPY(user->yiwork[i],-1.0,user->yi[i-1]);
596: }
597: Gather_i(C,user->yiwork,user->yi_scatter,user->nt);
598: VecAXPY(C,-1.0,user->q);
599: return 0;
600: }
602: PetscErrorCode Scatter(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat)
603: {
604: VecScatterBegin(s_scat,x,state,INSERT_VALUES,SCATTER_FORWARD);
605: VecScatterEnd(s_scat,x,state,INSERT_VALUES,SCATTER_FORWARD);
606: VecScatterBegin(d_scat,x,design,INSERT_VALUES,SCATTER_FORWARD);
607: VecScatterEnd(d_scat,x,design,INSERT_VALUES,SCATTER_FORWARD);
608: return 0;
609: }
611: PetscErrorCode Scatter_i(Vec y, Vec *yi, VecScatter *scat, PetscInt nt)
612: {
613: PetscInt i;
615: for (i=0; i<nt; i++) {
616: VecScatterBegin(scat[i],y,yi[i],INSERT_VALUES,SCATTER_FORWARD);
617: VecScatterEnd(scat[i],y,yi[i],INSERT_VALUES,SCATTER_FORWARD);
618: }
619: return 0;
620: }
622: PetscErrorCode Gather(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat)
623: {
624: VecScatterBegin(s_scat,state,x,INSERT_VALUES,SCATTER_REVERSE);
625: VecScatterEnd(s_scat,state,x,INSERT_VALUES,SCATTER_REVERSE);
626: VecScatterBegin(d_scat,design,x,INSERT_VALUES,SCATTER_REVERSE);
627: VecScatterEnd(d_scat,design,x,INSERT_VALUES,SCATTER_REVERSE);
628: return 0;
629: }
631: PetscErrorCode Gather_i(Vec y, Vec *yi, VecScatter *scat, PetscInt nt)
632: {
633: PetscInt i;
635: for (i=0; i<nt; i++) {
636: VecScatterBegin(scat[i],yi[i],y,INSERT_VALUES,SCATTER_REVERSE);
637: VecScatterEnd(scat[i],yi[i],y,INSERT_VALUES,SCATTER_REVERSE);
638: }
639: return 0;
640: }
642: PetscErrorCode ParabolicInitialize(AppCtx *user)
643: {
644: PetscInt m,n,i,j,k,linear_index,istart,iend,iblock,lo,hi,lo2,hi2;
645: Vec XX,YY,ZZ,XXwork,YYwork,ZZwork,UTwork,yi,di,bc;
646: PetscReal *x, *y, *z;
647: PetscReal h,stime;
648: PetscScalar hinv,neg_hinv,half = 0.5,sqrt_beta;
649: PetscInt im,indx1,indx2,indy1,indy2,indz1,indz2,nx,ny,nz;
650: PetscReal xri,yri,zri,xim,yim,zim,dx1,dx2,dy1,dy2,dz1,dz2,Dx,Dy,Dz;
651: PetscScalar v,vx,vy,vz;
652: IS is_from_y,is_to_yi,is_from_d,is_to_di;
653: /* Data locations */
654: PetscScalar xr[64] = {0.4970, 0.8498, 0.7814, 0.6268, 0.7782, 0.6402, 0.3617, 0.3160,
655: 0.3610, 0.5298, 0.6987, 0.3331, 0.7962, 0.5596, 0.3866, 0.6774,
656: 0.5407, 0.4518, 0.6702, 0.6061, 0.7580, 0.8997, 0.5198, 0.8326,
657: 0.2138, 0.9198, 0.3000, 0.2833, 0.8288, 0.7076, 0.1820, 0.0728,
658: 0.8447, 0.2367, 0.3239, 0.6413, 0.3114, 0.4731, 0.1192, 0.9273,
659: 0.5724, 0.4331, 0.5136, 0.3547, 0.4413, 0.2602, 0.5698, 0.7278,
660: 0.5261, 0.6230, 0.2454, 0.3948, 0.7479, 0.6582, 0.4660, 0.5594,
661: 0.7574, 0.1143, 0.5900, 0.1065, 0.4260, 0.3294, 0.8276, 0.0756};
663: PetscScalar yr[64] = {0.7345, 0.9120, 0.9288, 0.7528, 0.4463, 0.4985, 0.2497, 0.6256,
664: 0.3425, 0.9026, 0.6983, 0.4230, 0.7140, 0.2970, 0.4474, 0.8792,
665: 0.6604, 0.2485, 0.7968, 0.6127, 0.1796, 0.2437, 0.5938, 0.6137,
666: 0.3867, 0.5658, 0.4575, 0.1009, 0.0863, 0.3361, 0.0738, 0.3985,
667: 0.6602, 0.1437, 0.0934, 0.5983, 0.5950, 0.0763, 0.0768, 0.2288,
668: 0.5761, 0.1129, 0.3841, 0.6150, 0.6904, 0.6686, 0.1361, 0.4601,
669: 0.4491, 0.3716, 0.1969, 0.6537, 0.6743, 0.6991, 0.4811, 0.5480,
670: 0.1684, 0.4569, 0.6889, 0.8437, 0.3015, 0.2854, 0.8199, 0.2658};
672: PetscScalar zr[64] = {0.7668, 0.8573, 0.2654, 0.2719, 0.1060, 0.1311, 0.6232, 0.2295,
673: 0.8009, 0.2147, 0.2119, 0.9325, 0.4473, 0.3600, 0.3374, 0.3819,
674: 0.4066, 0.5801, 0.1673, 0.0959, 0.4638, 0.8236, 0.8800, 0.2939,
675: 0.2028, 0.8262, 0.2706, 0.6276, 0.9085, 0.6443, 0.8241, 0.0712,
676: 0.1824, 0.7789, 0.4389, 0.8415, 0.7055, 0.6639, 0.3653, 0.2078,
677: 0.1987, 0.2297, 0.4321, 0.8115, 0.4915, 0.7764, 0.4657, 0.4627,
678: 0.4569, 0.4232, 0.8514, 0.0674, 0.3227, 0.1055, 0.6690, 0.6313,
679: 0.9226, 0.5461, 0.4126, 0.2364, 0.6096, 0.7042, 0.3914, 0.0711};
681: PetscMalloc1(user->mx,&x);
682: PetscMalloc1(user->mx,&y);
683: PetscMalloc1(user->mx,&z);
684: user->jformed = PETSC_FALSE;
685: user->dsg_formed = PETSC_FALSE;
687: n = user->mx * user->mx * user->mx;
688: m = 3 * user->mx * user->mx * (user->mx-1);
689: sqrt_beta = PetscSqrtScalar(user->beta);
691: user->ksp_its = 0;
692: user->ksp_its_initial = 0;
694: stime = (PetscReal)user->nt/user->ns;
695: PetscMalloc1(user->ns,&user->sample_times);
696: for (i=0; i<user->ns; i++) {
697: user->sample_times[i] = (PetscInt)(stime*((PetscReal)i+1.0)-0.5);
698: }
700: VecCreate(PETSC_COMM_WORLD,&XX);
701: VecCreate(PETSC_COMM_WORLD,&user->q);
702: VecSetSizes(XX,PETSC_DECIDE,n);
703: VecSetSizes(user->q,PETSC_DECIDE,n*user->nt);
704: VecSetFromOptions(XX);
705: VecSetFromOptions(user->q);
707: VecDuplicate(XX,&YY);
708: VecDuplicate(XX,&ZZ);
709: VecDuplicate(XX,&XXwork);
710: VecDuplicate(XX,&YYwork);
711: VecDuplicate(XX,&ZZwork);
712: VecDuplicate(XX,&UTwork);
713: VecDuplicate(XX,&user->utrue);
714: VecDuplicate(XX,&bc);
716: /* Generate 3D grid, and collect ns (1<=ns<=8) right-hand-side vectors into user->q */
717: h = 1.0/user->mx;
718: hinv = user->mx;
719: neg_hinv = -hinv;
721: VecGetOwnershipRange(XX,&istart,&iend);
722: for (linear_index=istart; linear_index<iend; linear_index++) {
723: i = linear_index % user->mx;
724: j = ((linear_index-i)/user->mx) % user->mx;
725: k = ((linear_index-i)/user->mx-j) / user->mx;
726: vx = h*(i+0.5);
727: vy = h*(j+0.5);
728: vz = h*(k+0.5);
729: VecSetValues(XX,1,&linear_index,&vx,INSERT_VALUES);
730: VecSetValues(YY,1,&linear_index,&vy,INSERT_VALUES);
731: VecSetValues(ZZ,1,&linear_index,&vz,INSERT_VALUES);
732: if ((vx<0.6) && (vx>0.4) && (vy<0.6) && (vy>0.4) && (vy<0.6) && (vz<0.6) && (vz>0.4)) {
733: v = 1000.0;
734: VecSetValues(bc,1,&linear_index,&v,INSERT_VALUES);
735: }
736: }
738: VecAssemblyBegin(XX);
739: VecAssemblyEnd(XX);
740: VecAssemblyBegin(YY);
741: VecAssemblyEnd(YY);
742: VecAssemblyBegin(ZZ);
743: VecAssemblyEnd(ZZ);
744: VecAssemblyBegin(bc);
745: VecAssemblyEnd(bc);
747: /* Compute true parameter function
748: ut = 0.5 + exp(-10*((x-0.5)^2+(y-0.5)^2+(z-0.5)^2)) */
749: VecCopy(XX,XXwork);
750: VecCopy(YY,YYwork);
751: VecCopy(ZZ,ZZwork);
753: VecShift(XXwork,-0.5);
754: VecShift(YYwork,-0.5);
755: VecShift(ZZwork,-0.5);
757: VecPointwiseMult(XXwork,XXwork,XXwork);
758: VecPointwiseMult(YYwork,YYwork,YYwork);
759: VecPointwiseMult(ZZwork,ZZwork,ZZwork);
761: VecCopy(XXwork,user->utrue);
762: VecAXPY(user->utrue,1.0,YYwork);
763: VecAXPY(user->utrue,1.0,ZZwork);
764: VecScale(user->utrue,-10.0);
765: VecExp(user->utrue);
766: VecShift(user->utrue,0.5);
768: VecDestroy(&XX);
769: VecDestroy(&YY);
770: VecDestroy(&ZZ);
771: VecDestroy(&XXwork);
772: VecDestroy(&YYwork);
773: VecDestroy(&ZZwork);
774: VecDestroy(&UTwork);
776: /* Initial guess and reference model */
777: VecDuplicate(user->utrue,&user->ur);
778: VecSet(user->ur,0.0);
780: /* Generate Grad matrix */
781: MatCreate(PETSC_COMM_WORLD,&user->Grad);
782: MatSetSizes(user->Grad,PETSC_DECIDE,PETSC_DECIDE,m,n);
783: MatSetFromOptions(user->Grad);
784: MatMPIAIJSetPreallocation(user->Grad,2,NULL,2,NULL);
785: MatSeqAIJSetPreallocation(user->Grad,2,NULL);
786: MatGetOwnershipRange(user->Grad,&istart,&iend);
788: for (i=istart; i<iend; i++) {
789: if (i<m/3) {
790: iblock = i / (user->mx-1);
791: j = iblock*user->mx + (i % (user->mx-1));
792: MatSetValues(user->Grad,1,&i,1,&j,&neg_hinv,INSERT_VALUES);
793: j = j+1;
794: MatSetValues(user->Grad,1,&i,1,&j,&hinv,INSERT_VALUES);
795: }
796: if (i>=m/3 && i<2*m/3) {
797: iblock = (i-m/3) / (user->mx*(user->mx-1));
798: j = iblock*user->mx*user->mx + ((i-m/3) % (user->mx*(user->mx-1)));
799: MatSetValues(user->Grad,1,&i,1,&j,&neg_hinv,INSERT_VALUES);
800: j = j + user->mx;
801: MatSetValues(user->Grad,1,&i,1,&j,&hinv,INSERT_VALUES);
802: }
803: if (i>=2*m/3) {
804: j = i-2*m/3;
805: MatSetValues(user->Grad,1,&i,1,&j,&neg_hinv,INSERT_VALUES);
806: j = j + user->mx*user->mx;
807: MatSetValues(user->Grad,1,&i,1,&j,&hinv,INSERT_VALUES);
808: }
809: }
811: MatAssemblyBegin(user->Grad,MAT_FINAL_ASSEMBLY);
812: MatAssemblyEnd(user->Grad,MAT_FINAL_ASSEMBLY);
814: /* Generate arithmetic averaging matrix Av */
815: MatCreate(PETSC_COMM_WORLD,&user->Av);
816: MatSetSizes(user->Av,PETSC_DECIDE,PETSC_DECIDE,m,n);
817: MatSetFromOptions(user->Av);
818: MatMPIAIJSetPreallocation(user->Av,2,NULL,2,NULL);
819: MatSeqAIJSetPreallocation(user->Av,2,NULL);
820: MatGetOwnershipRange(user->Av,&istart,&iend);
822: for (i=istart; i<iend; i++) {
823: if (i<m/3) {
824: iblock = i / (user->mx-1);
825: j = iblock*user->mx + (i % (user->mx-1));
826: MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES);
827: j = j+1;
828: MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES);
829: }
830: if (i>=m/3 && i<2*m/3) {
831: iblock = (i-m/3) / (user->mx*(user->mx-1));
832: j = iblock*user->mx*user->mx + ((i-m/3) % (user->mx*(user->mx-1)));
833: MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES);
834: j = j + user->mx;
835: MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES);
836: }
837: if (i>=2*m/3) {
838: j = i-2*m/3;
839: MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES);
840: j = j + user->mx*user->mx;
841: MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES);
842: }
843: }
845: MatAssemblyBegin(user->Av,MAT_FINAL_ASSEMBLY);
846: MatAssemblyEnd(user->Av,MAT_FINAL_ASSEMBLY);
848: /* Generate transpose of averaging matrix Av */
849: MatTranspose(user->Av,MAT_INITIAL_MATRIX,&user->AvT);
851: MatCreate(PETSC_COMM_WORLD,&user->L);
852: MatSetSizes(user->L,PETSC_DECIDE,PETSC_DECIDE,m+n,n);
853: MatSetFromOptions(user->L);
854: MatMPIAIJSetPreallocation(user->L,2,NULL,2,NULL);
855: MatSeqAIJSetPreallocation(user->L,2,NULL);
856: MatGetOwnershipRange(user->L,&istart,&iend);
858: for (i=istart; i<iend; i++) {
859: if (i<m/3) {
860: iblock = i / (user->mx-1);
861: j = iblock*user->mx + (i % (user->mx-1));
862: MatSetValues(user->L,1,&i,1,&j,&neg_hinv,INSERT_VALUES);
863: j = j+1;
864: MatSetValues(user->L,1,&i,1,&j,&hinv,INSERT_VALUES);
865: }
866: if (i>=m/3 && i<2*m/3) {
867: iblock = (i-m/3) / (user->mx*(user->mx-1));
868: j = iblock*user->mx*user->mx + ((i-m/3) % (user->mx*(user->mx-1)));
869: MatSetValues(user->L,1,&i,1,&j,&neg_hinv,INSERT_VALUES);
870: j = j + user->mx;
871: MatSetValues(user->L,1,&i,1,&j,&hinv,INSERT_VALUES);
872: }
873: if (i>=2*m/3 && i<m) {
874: j = i-2*m/3;
875: MatSetValues(user->L,1,&i,1,&j,&neg_hinv,INSERT_VALUES);
876: j = j + user->mx*user->mx;
877: MatSetValues(user->L,1,&i,1,&j,&hinv,INSERT_VALUES);
878: }
879: if (i>=m) {
880: j = i - m;
881: MatSetValues(user->L,1,&i,1,&j,&sqrt_beta,INSERT_VALUES);
882: }
883: }
885: MatAssemblyBegin(user->L,MAT_FINAL_ASSEMBLY);
886: MatAssemblyEnd(user->L,MAT_FINAL_ASSEMBLY);
888: MatScale(user->L,PetscPowScalar(h,1.5));
890: /* Generate Div matrix */
891: MatTranspose(user->Grad,MAT_INITIAL_MATRIX,&user->Div);
893: /* Build work vectors and matrices */
894: VecCreate(PETSC_COMM_WORLD,&user->S);
895: VecSetSizes(user->S, PETSC_DECIDE, user->mx*user->mx*(user->mx-1)*3);
896: VecSetFromOptions(user->S);
898: VecCreate(PETSC_COMM_WORLD,&user->lwork);
899: VecSetSizes(user->lwork,PETSC_DECIDE,m+user->mx*user->mx*user->mx);
900: VecSetFromOptions(user->lwork);
902: MatDuplicate(user->Div,MAT_SHARE_NONZERO_PATTERN,&user->Divwork);
903: MatDuplicate(user->Av,MAT_SHARE_NONZERO_PATTERN,&user->Avwork);
905: VecCreate(PETSC_COMM_WORLD,&user->d);
906: VecSetSizes(user->d,PETSC_DECIDE,user->ndata*user->nt);
907: VecSetFromOptions(user->d);
909: VecDuplicate(user->S,&user->Swork);
910: VecDuplicate(user->S,&user->Av_u);
911: VecDuplicate(user->S,&user->Twork);
912: VecDuplicate(user->S,&user->Rwork);
913: VecDuplicate(user->y,&user->ywork);
914: VecDuplicate(user->u,&user->uwork);
915: VecDuplicate(user->u,&user->js_diag);
916: VecDuplicate(user->c,&user->cwork);
917: VecDuplicate(user->d,&user->dwork);
919: /* Create matrix-free shell user->Js for computing A*x */
920: MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m*user->nt,user->m*user->nt,user,&user->Js);
921: MatShellSetOperation(user->Js,MATOP_MULT,(void(*)(void))StateMatMult);
922: MatShellSetOperation(user->Js,MATOP_DUPLICATE,(void(*)(void))StateMatDuplicate);
923: MatShellSetOperation(user->Js,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatMultTranspose);
924: MatShellSetOperation(user->Js,MATOP_GET_DIAGONAL,(void(*)(void))StateMatGetDiagonal);
926: /* Diagonal blocks of user->Js */
927: MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m,user->m,user,&user->JsBlock);
928: MatShellSetOperation(user->JsBlock,MATOP_MULT,(void(*)(void))StateMatBlockMult);
929: /* Blocks are symmetric */
930: MatShellSetOperation(user->JsBlock,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatBlockMult);
932: /* Create a matrix-free shell user->JsBlockPrec for computing (U+D)\D*(L+D)\x, where JsBlock = L+D+U,
933: D is diagonal, L is strictly lower triangular, and U is strictly upper triangular.
934: This is an SSOR preconditioner for user->JsBlock. */
935: MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m,user->m,user,&user->JsBlockPrec);
936: MatShellSetOperation(user->JsBlockPrec,MATOP_MULT,(void(*)(void))StateMatBlockPrecMult);
937: /* JsBlockPrec is symmetric */
938: MatShellSetOperation(user->JsBlockPrec,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatBlockPrecMult);
939: MatSetOption(user->JsBlockPrec,MAT_SYMMETRY_ETERNAL,PETSC_TRUE);
941: /* Create a matrix-free shell user->Jd for computing B*x */
942: MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m*user->nt,user->m,user,&user->Jd);
943: MatShellSetOperation(user->Jd,MATOP_MULT,(void(*)(void))DesignMatMult);
944: MatShellSetOperation(user->Jd,MATOP_MULT_TRANSPOSE,(void(*)(void))DesignMatMultTranspose);
946: /* User-defined routines for computing user->Js\x and user->Js^T\x*/
947: MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m*user->nt,user->m*user->nt,user,&user->JsInv);
948: MatShellSetOperation(user->JsInv,MATOP_MULT,(void(*)(void))StateMatInvMult);
949: MatShellSetOperation(user->JsInv,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatInvTransposeMult);
951: /* Solver options and tolerances */
952: KSPCreate(PETSC_COMM_WORLD,&user->solver);
953: KSPSetType(user->solver,KSPCG);
954: KSPSetOperators(user->solver,user->JsBlock,user->JsBlockPrec);
955: KSPSetInitialGuessNonzero(user->solver,PETSC_FALSE);
956: KSPSetTolerances(user->solver,1e-4,1e-20,1e3,500);
957: KSPSetFromOptions(user->solver);
958: KSPGetPC(user->solver,&user->prec);
959: PCSetType(user->prec,PCSHELL);
961: PCShellSetApply(user->prec,StateMatBlockPrecMult);
962: PCShellSetApplyTranspose(user->prec,StateMatBlockPrecMult);
963: PCShellSetContext(user->prec,user);
965: /* Create scatter from y to y_1,y_2,...,y_nt */
966: PetscMalloc1(user->nt*user->m,&user->yi_scatter);
967: VecCreate(PETSC_COMM_WORLD,&yi);
968: VecSetSizes(yi,PETSC_DECIDE,user->mx*user->mx*user->mx);
969: VecSetFromOptions(yi);
970: VecDuplicateVecs(yi,user->nt,&user->yi);
971: VecDuplicateVecs(yi,user->nt,&user->yiwork);
973: VecGetOwnershipRange(user->y,&lo2,&hi2);
974: istart = 0;
975: for (i=0; i<user->nt; i++) {
976: VecGetOwnershipRange(user->yi[i],&lo,&hi);
977: ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_to_yi);
978: ISCreateStride(PETSC_COMM_SELF,hi-lo,lo2+istart,1,&is_from_y);
979: VecScatterCreate(user->y,is_from_y,user->yi[i],is_to_yi,&user->yi_scatter[i]);
980: istart = istart + hi-lo;
981: ISDestroy(&is_to_yi);
982: ISDestroy(&is_from_y);
983: }
984: VecDestroy(&yi);
986: /* Create scatter from d to d_1,d_2,...,d_ns */
987: PetscMalloc1(user->ns*user->ndata,&user->di_scatter);
988: VecCreate(PETSC_COMM_WORLD,&di);
989: VecSetSizes(di,PETSC_DECIDE,user->ndata);
990: VecSetFromOptions(di);
991: VecDuplicateVecs(di,user->ns,&user->di);
992: VecGetOwnershipRange(user->d,&lo2,&hi2);
993: istart = 0;
994: for (i=0; i<user->ns; i++) {
995: VecGetOwnershipRange(user->di[i],&lo,&hi);
996: ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_to_di);
997: ISCreateStride(PETSC_COMM_SELF,hi-lo,lo2+istart,1,&is_from_d);
998: VecScatterCreate(user->d,is_from_d,user->di[i],is_to_di,&user->di_scatter[i]);
999: istart = istart + hi-lo;
1000: ISDestroy(&is_to_di);
1001: ISDestroy(&is_from_d);
1002: }
1003: VecDestroy(&di);
1005: /* Assemble RHS of forward problem */
1006: VecCopy(bc,user->yiwork[0]);
1007: for (i=1; i<user->nt; i++) {
1008: VecSet(user->yiwork[i],0.0);
1009: }
1010: Gather_i(user->q,user->yiwork,user->yi_scatter,user->nt);
1011: VecDestroy(&bc);
1013: /* Compute true state function yt given ut */
1014: VecCreate(PETSC_COMM_WORLD,&user->ytrue);
1015: VecSetSizes(user->ytrue,PETSC_DECIDE,n*user->nt);
1016: VecSetFromOptions(user->ytrue);
1018: /* First compute Av_u = Av*exp(-u) */
1019: VecSet(user->uwork,0);
1020: VecAXPY(user->uwork,-1.0,user->utrue); /* Note: user->utrue */
1021: VecExp(user->uwork);
1022: MatMult(user->Av,user->uwork,user->Av_u);
1024: /* Symbolic DSG = Div * Grad */
1025: MatProductCreate(user->Div,user->Grad,NULL,&user->DSG);
1026: MatProductSetType(user->DSG,MATPRODUCT_AB);
1027: MatProductSetAlgorithm(user->DSG,"default");
1028: MatProductSetFill(user->DSG,PETSC_DEFAULT);
1029: MatProductSetFromOptions(user->DSG);
1030: MatProductSymbolic(user->DSG);
1032: user->dsg_formed = PETSC_TRUE;
1034: /* Next form DSG = Div*Grad */
1035: MatCopy(user->Div,user->Divwork,SAME_NONZERO_PATTERN);
1036: MatDiagonalScale(user->Divwork,NULL,user->Av_u);
1037: if (user->dsg_formed) {
1038: MatProductNumeric(user->DSG);
1039: } else {
1040: MatMatMult(user->Div,user->Grad,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&user->DSG);
1041: user->dsg_formed = PETSC_TRUE;
1042: }
1043: /* B = speye(nx^3) + ht*DSG; */
1044: MatScale(user->DSG,user->ht);
1045: MatShift(user->DSG,1.0);
1047: /* Now solve for ytrue */
1048: StateMatInvMult(user->Js,user->q,user->ytrue);
1050: /* Initial guess y0 for state given u0 */
1052: /* First compute Av_u = Av*exp(-u) */
1053: VecSet(user->uwork,0);
1054: VecAXPY(user->uwork,-1.0,user->u); /* Note: user->u */
1055: VecExp(user->uwork);
1056: MatMult(user->Av,user->uwork,user->Av_u);
1058: /* Next form DSG = Div*S*Grad */
1059: MatCopy(user->Div,user->Divwork,SAME_NONZERO_PATTERN);
1060: MatDiagonalScale(user->Divwork,NULL,user->Av_u);
1061: if (user->dsg_formed) {
1062: MatProductNumeric(user->DSG);
1063: } else {
1064: MatMatMult(user->Div,user->Grad,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&user->DSG);
1066: user->dsg_formed = PETSC_TRUE;
1067: }
1068: /* B = speye(nx^3) + ht*DSG; */
1069: MatScale(user->DSG,user->ht);
1070: MatShift(user->DSG,1.0);
1072: /* Now solve for y */
1073: StateMatInvMult(user->Js,user->q,user->y);
1075: /* Construct projection matrix Q, a block diagonal matrix consisting of nt copies of Qblock along the diagonal */
1076: MatCreate(PETSC_COMM_WORLD,&user->Qblock);
1077: MatSetSizes(user->Qblock,PETSC_DECIDE,PETSC_DECIDE,user->ndata,n);
1078: MatSetFromOptions(user->Qblock);
1079: MatMPIAIJSetPreallocation(user->Qblock,8,NULL,8,NULL);
1080: MatSeqAIJSetPreallocation(user->Qblock,8,NULL);
1082: for (i=0; i<user->mx; i++) {
1083: x[i] = h*(i+0.5);
1084: y[i] = h*(i+0.5);
1085: z[i] = h*(i+0.5);
1086: }
1088: MatGetOwnershipRange(user->Qblock,&istart,&iend);
1089: nx = user->mx; ny = user->mx; nz = user->mx;
1090: for (i=istart; i<iend; i++) {
1091: xri = xr[i];
1092: im = 0;
1093: xim = x[im];
1094: while (xri>xim && im<nx) {
1095: im = im+1;
1096: xim = x[im];
1097: }
1098: indx1 = im-1;
1099: indx2 = im;
1100: dx1 = xri - x[indx1];
1101: dx2 = x[indx2] - xri;
1103: yri = yr[i];
1104: im = 0;
1105: yim = y[im];
1106: while (yri>yim && im<ny) {
1107: im = im+1;
1108: yim = y[im];
1109: }
1110: indy1 = im-1;
1111: indy2 = im;
1112: dy1 = yri - y[indy1];
1113: dy2 = y[indy2] - yri;
1115: zri = zr[i];
1116: im = 0;
1117: zim = z[im];
1118: while (zri>zim && im<nz) {
1119: im = im+1;
1120: zim = z[im];
1121: }
1122: indz1 = im-1;
1123: indz2 = im;
1124: dz1 = zri - z[indz1];
1125: dz2 = z[indz2] - zri;
1127: Dx = x[indx2] - x[indx1];
1128: Dy = y[indy2] - y[indy1];
1129: Dz = z[indz2] - z[indz1];
1131: j = indx1 + indy1*nx + indz1*nx*ny;
1132: v = (1-dx1/Dx)*(1-dy1/Dy)*(1-dz1/Dz);
1133: MatSetValues(user->Qblock,1,&i,1,&j,&v,INSERT_VALUES);
1135: j = indx1 + indy1*nx + indz2*nx*ny;
1136: v = (1-dx1/Dx)*(1-dy1/Dy)*(1-dz2/Dz);
1137: MatSetValues(user->Qblock,1,&i,1,&j,&v,INSERT_VALUES);
1139: j = indx1 + indy2*nx + indz1*nx*ny;
1140: v = (1-dx1/Dx)*(1-dy2/Dy)*(1-dz1/Dz);
1141: MatSetValues(user->Qblock,1,&i,1,&j,&v,INSERT_VALUES);
1143: j = indx1 + indy2*nx + indz2*nx*ny;
1144: v = (1-dx1/Dx)*(1-dy2/Dy)*(1-dz2/Dz);
1145: MatSetValues(user->Qblock,1,&i,1,&j,&v,INSERT_VALUES);
1147: j = indx2 + indy1*nx + indz1*nx*ny;
1148: v = (1-dx2/Dx)*(1-dy1/Dy)*(1-dz1/Dz);
1149: MatSetValues(user->Qblock,1,&i,1,&j,&v,INSERT_VALUES);
1151: j = indx2 + indy1*nx + indz2*nx*ny;
1152: v = (1-dx2/Dx)*(1-dy1/Dy)*(1-dz2/Dz);
1153: MatSetValues(user->Qblock,1,&i,1,&j,&v,INSERT_VALUES);
1155: j = indx2 + indy2*nx + indz1*nx*ny;
1156: v = (1-dx2/Dx)*(1-dy2/Dy)*(1-dz1/Dz);
1157: MatSetValues(user->Qblock,1,&i,1,&j,&v,INSERT_VALUES);
1159: j = indx2 + indy2*nx + indz2*nx*ny;
1160: v = (1-dx2/Dx)*(1-dy2/Dy)*(1-dz2/Dz);
1161: MatSetValues(user->Qblock,1,&i,1,&j,&v,INSERT_VALUES);
1162: }
1163: MatAssemblyBegin(user->Qblock,MAT_FINAL_ASSEMBLY);
1164: MatAssemblyEnd(user->Qblock,MAT_FINAL_ASSEMBLY);
1166: MatTranspose(user->Qblock,MAT_INITIAL_MATRIX,&user->QblockT);
1167: MatTranspose(user->L,MAT_INITIAL_MATRIX,&user->LT);
1169: /* Add noise to the measurement data */
1170: VecSet(user->ywork,1.0);
1171: VecAYPX(user->ywork,user->noise,user->ytrue);
1172: Scatter_i(user->ywork,user->yiwork,user->yi_scatter,user->nt);
1173: for (j=0; j<user->ns; j++) {
1174: i = user->sample_times[j];
1175: MatMult(user->Qblock,user->yiwork[i],user->di[j]);
1176: }
1177: Gather_i(user->d,user->di,user->di_scatter,user->ns);
1179: /* Now that initial conditions have been set, let the user pass tolerance options to the KSP solver */
1180: KSPSetFromOptions(user->solver);
1181: PetscFree(x);
1182: PetscFree(y);
1183: PetscFree(z);
1184: return 0;
1185: }
1187: PetscErrorCode ParabolicDestroy(AppCtx *user)
1188: {
1189: PetscInt i;
1191: MatDestroy(&user->Qblock);
1192: MatDestroy(&user->QblockT);
1193: MatDestroy(&user->Div);
1194: MatDestroy(&user->Divwork);
1195: MatDestroy(&user->Grad);
1196: MatDestroy(&user->Av);
1197: MatDestroy(&user->Avwork);
1198: MatDestroy(&user->AvT);
1199: MatDestroy(&user->DSG);
1200: MatDestroy(&user->L);
1201: MatDestroy(&user->LT);
1202: KSPDestroy(&user->solver);
1203: MatDestroy(&user->Js);
1204: MatDestroy(&user->Jd);
1205: MatDestroy(&user->JsInv);
1206: MatDestroy(&user->JsBlock);
1207: MatDestroy(&user->JsBlockPrec);
1208: VecDestroy(&user->u);
1209: VecDestroy(&user->uwork);
1210: VecDestroy(&user->utrue);
1211: VecDestroy(&user->y);
1212: VecDestroy(&user->ywork);
1213: VecDestroy(&user->ytrue);
1214: VecDestroyVecs(user->nt,&user->yi);
1215: VecDestroyVecs(user->nt,&user->yiwork);
1216: VecDestroyVecs(user->ns,&user->di);
1217: PetscFree(user->yi);
1218: PetscFree(user->yiwork);
1219: PetscFree(user->di);
1220: VecDestroy(&user->c);
1221: VecDestroy(&user->cwork);
1222: VecDestroy(&user->ur);
1223: VecDestroy(&user->q);
1224: VecDestroy(&user->d);
1225: VecDestroy(&user->dwork);
1226: VecDestroy(&user->lwork);
1227: VecDestroy(&user->S);
1228: VecDestroy(&user->Swork);
1229: VecDestroy(&user->Av_u);
1230: VecDestroy(&user->Twork);
1231: VecDestroy(&user->Rwork);
1232: VecDestroy(&user->js_diag);
1233: ISDestroy(&user->s_is);
1234: ISDestroy(&user->d_is);
1235: VecScatterDestroy(&user->state_scatter);
1236: VecScatterDestroy(&user->design_scatter);
1237: for (i=0; i<user->nt; i++) {
1238: VecScatterDestroy(&user->yi_scatter[i]);
1239: }
1240: for (i=0; i<user->ns; i++) {
1241: VecScatterDestroy(&user->di_scatter[i]);
1242: }
1243: PetscFree(user->yi_scatter);
1244: PetscFree(user->di_scatter);
1245: PetscFree(user->sample_times);
1246: return 0;
1247: }
1249: PetscErrorCode ParabolicMonitor(Tao tao, void *ptr)
1250: {
1251: Vec X;
1252: PetscReal unorm,ynorm;
1253: AppCtx *user = (AppCtx*)ptr;
1255: TaoGetSolution(tao,&X);
1256: Scatter(X,user->ywork,user->state_scatter,user->uwork,user->design_scatter);
1257: VecAXPY(user->ywork,-1.0,user->ytrue);
1258: VecAXPY(user->uwork,-1.0,user->utrue);
1259: VecNorm(user->uwork,NORM_2,&unorm);
1260: VecNorm(user->ywork,NORM_2,&ynorm);
1261: PetscPrintf(MPI_COMM_WORLD, "||u-ut||=%g ||y-yt||=%g\n",(double)unorm,(double)ynorm);
1262: return 0;
1263: }
1265: /*TEST
1267: build:
1268: requires: !complex
1270: test:
1271: args: -tao_cmonitor -tao_type lcl -ns 1 -tao_gatol 1.e-4 -ksp_max_it 30
1272: requires: !single
1274: TEST*/