Actual source code: parabolic.c
petsc-dev 2014-02-02
1: #include <petsc-private/taoimpl.h>
2: #include "../src/tao/pde_constrained/impls/lcl/lcl.h"
4: /*T
5: Concepts: TAO^Solving a system of nonlinear equations, nonlinear least squares
6: Routines: TaoCreate();
7: Routines: TaoSetType();
8: Routines: TaoSetInitialVector();
9: Routines: TaoSetObjectiveRoutine();
10: Routines: TaoSetGradientRoutine();
11: Routines: TaoSetConstraintsRoutine();
12: Routines: TaoSetJacobianStateRoutine();
13: Routines: TaoSetJacobianDesignRoutine();
14: Routines: TaoSetStateDesignIS();
15: Routines: TaoSetFromOptions();
16: Routines: TaoSetHistory(); TaoGetHistory();
17: Routines: TaoSolve();
18: Routines: TaoGetTerminationReason(); TaoDestroy();
19: Processors: n
20: T*/
22: typedef struct {
23: PetscInt n; /* Number of variables */
24: PetscInt m; /* Number of constraints per time step */
25: PetscInt mx; /* grid points in each direction */
26: PetscInt nt; /* Number of time steps; as of now, must be divisible by 8 */
27: PetscInt ndata; /* Number of data points per sample */
28: PetscInt ns; /* Number of samples */
29: PetscInt *sample_times; /* Times of samples */
30: IS s_is;
31: IS d_is;
33: VecScatter state_scatter;
34: VecScatter design_scatter;
35: VecScatter *yi_scatter;
36: VecScatter *di_scatter;
38: Mat Js,Jd,JsBlockPrec,JsInv,JsBlock;
39: PetscBool jformed,dsg_formed;
41: PetscReal alpha; /* Regularization parameter */
42: PetscReal beta; /* Weight attributed to ||u||^2 in regularization functional */
43: PetscReal noise; /* Amount of noise to add to data */
44: PetscReal ht; /* Time step */
46: Mat Qblock,QblockT;
47: Mat L,LT;
48: Mat Div,Divwork;
49: Mat Grad;
50: Mat Av,Avwork,AvT;
51: Mat DSG;
52: Vec q;
53: Vec ur; /* reference */
55: Vec d;
56: Vec dwork;
57: Vec *di;
59: Vec y; /* state variables */
60: Vec ywork;
62: Vec ytrue;
63: Vec *yi,*yiwork;
65: Vec u; /* design variables */
66: Vec uwork;
68: Vec utrue;
69: Vec js_diag;
70: Vec c; /* constraint vector */
71: Vec cwork;
73: Vec lwork;
74: Vec S;
75: Vec Rwork,Swork,Twork;
76: Vec Av_u;
78: KSP solver;
79: PC prec;
81: TAO_LCL *lcl;
82: PetscInt ksp_its;
83: PetscInt ksp_its_initial;
84: } AppCtx;
87: PetscErrorCode FormFunction(Tao, Vec, PetscReal*, void*);
88: PetscErrorCode FormGradient(Tao, Vec, Vec, void*);
89: PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal*, Vec, void*);
90: PetscErrorCode FormJacobianState(Tao, Vec, Mat*, Mat*, Mat*, MatStructure*,void*);
91: PetscErrorCode FormJacobianDesign(Tao, Vec, Mat*, void*);
92: PetscErrorCode FormConstraints(Tao, Vec, Vec, void*);
93: PetscErrorCode FormHessian(Tao, Vec, Mat*, Mat*, MatStructure*, void*);
94: PetscErrorCode Gather(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat);
95: PetscErrorCode Scatter(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat);
96: PetscErrorCode ParabolicInitialize(AppCtx *user);
97: PetscErrorCode ParabolicDestroy(AppCtx *user);
98: PetscErrorCode ParabolicMonitor(Tao, void*);
100: PetscErrorCode StateMatMult(Mat,Vec,Vec);
101: PetscErrorCode StateMatBlockMult(Mat,Vec,Vec);
102: PetscErrorCode StateMatMultTranspose(Mat,Vec,Vec);
103: PetscErrorCode StateMatGetDiagonal(Mat,Vec);
104: PetscErrorCode StateMatDuplicate(Mat,MatDuplicateOption,Mat*);
105: PetscErrorCode StateMatInvMult(Mat,Vec,Vec);
106: PetscErrorCode StateMatInvTransposeMult(Mat,Vec,Vec);
107: PetscErrorCode StateMatBlockPrecMult(PC,Vec,Vec);
109: PetscErrorCode DesignMatMult(Mat,Vec,Vec);
110: PetscErrorCode DesignMatMultTranspose(Mat,Vec,Vec);
112: PetscErrorCode Gather_i(Vec,Vec*,VecScatter*,PetscInt);
113: PetscErrorCode Scatter_i(Vec,Vec*,VecScatter*,PetscInt);
115: static char help[]="";
119: int main(int argc, char **argv)
120: {
121: PetscErrorCode ierr;
122: Vec x,x0;
123: Tao tao;
124: TaoTerminationReason reason;
125: AppCtx user;
126: IS is_allstate,is_alldesign;
127: PetscInt lo,hi,hi2,lo2,ksp_old;
128: PetscBool flag;
129: PetscInt ntests = 1;
130: PetscInt i;
131: int stages[1];
133: PetscInitialize(&argc, &argv, (char*)0,help);
135: user.mx = 8;
136: PetscOptionsInt("-mx","Number of grid points in each direction","",user.mx,&user.mx,&flag);
137: user.nt = 8;
138: PetscOptionsInt("-nt","Number of time steps","",user.nt,&user.nt,&flag);
139: user.ndata = 64;
140: PetscOptionsInt("-ndata","Numbers of data points per sample","",user.ndata,&user.ndata,&flag);
141: user.ns = 8;
142: PetscOptionsInt("-ns","Number of samples","",user.ns,&user.ns,&flag);
143: user.alpha = 1.0;
144: PetscOptionsReal("-alpha","Regularization parameter","",user.alpha,&user.alpha,&flag);
145: user.beta = 0.01;
146: PetscOptionsReal("-beta","Weight attributed to ||u||^2 in regularization functional","",user.beta,&user.beta,&flag);
147: user.noise = 0.01;
148: PetscOptionsReal("-noise","Amount of noise to add to data","",user.noise,&user.noise,&flag);
150: user.m = user.mx*user.mx*user.mx; /* number of constraints per time step */
151: user.n = user.m*(user.nt+1); /* number of variables */
152: user.ht = (PetscReal)1/user.nt;
154: VecCreate(PETSC_COMM_WORLD,&user.u);
155: VecCreate(PETSC_COMM_WORLD,&user.y);
156: VecCreate(PETSC_COMM_WORLD,&user.c);
157: VecSetSizes(user.u,PETSC_DECIDE,user.n-user.m*user.nt);
158: VecSetSizes(user.y,PETSC_DECIDE,user.m*user.nt);
159: VecSetSizes(user.c,PETSC_DECIDE,user.m*user.nt);
160: VecSetFromOptions(user.u);
161: VecSetFromOptions(user.y);
162: VecSetFromOptions(user.c);
164: /* Create scatters for reduced spaces.
165: If the state vector y and design vector u are partitioned as
166: [y_1; y_2; ...; y_np] and [u_1; u_2; ...; u_np] (with np = # of processors),
167: then the solution vector x is organized as
168: [y_1; u_1; y_2; u_2; ...; y_np; u_np].
169: The index sets user.s_is and user.d_is correspond to the indices of the
170: state and design variables owned by the current processor.
171: */
172: VecCreate(PETSC_COMM_WORLD,&x);
174: VecGetOwnershipRange(user.y,&lo,&hi);
175: VecGetOwnershipRange(user.u,&lo2,&hi2);
177: ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_allstate);
178: ISCreateStride(PETSC_COMM_SELF,hi-lo,lo+lo2,1,&user.s_is);
180: ISCreateStride(PETSC_COMM_SELF,hi2-lo2,lo2,1,&is_alldesign);
181: ISCreateStride(PETSC_COMM_SELF,hi2-lo2,hi+lo2,1,&user.d_is);
183: VecSetSizes(x,hi-lo+hi2-lo2,user.n);
184: VecSetFromOptions(x);
186: VecScatterCreate(x,user.s_is,user.y,is_allstate,&user.state_scatter);
187: VecScatterCreate(x,user.d_is,user.u,is_alldesign,&user.design_scatter);
188: ISDestroy(&is_alldesign);
189: ISDestroy(&is_allstate);
191: /* Create TAO solver and set desired solution method */
192: TaoCreate(PETSC_COMM_WORLD,&tao);
193: TaoSetType(tao,"tao_lcl");
194: user.lcl = (TAO_LCL*)(tao->data);
195: /* Set up initial vectors and matrices */
196: ParabolicInitialize(&user);
198: Gather(x,user.y,user.state_scatter,user.u,user.design_scatter);
199: VecDuplicate(x,&x0);
200: VecCopy(x,x0);
202: /* Set solution vector with an initial guess */
203: TaoSetInitialVector(tao,x);
204: TaoSetObjectiveRoutine(tao, FormFunction, (void *)&user);
205: TaoSetGradientRoutine(tao, FormGradient, (void *)&user);
206: TaoSetConstraintsRoutine(tao, user.c, FormConstraints, (void *)&user);
208: TaoSetJacobianStateRoutine(tao, user.Js, user.JsBlockPrec, user.JsInv, FormJacobianState, (void *)&user);
209: TaoSetJacobianDesignRoutine(tao, user.Jd, FormJacobianDesign, (void *)&user);
211: TaoSetFromOptions(tao);
212: TaoSetStateDesignIS(tao,user.s_is,user.d_is);
214: /* SOLVE THE APPLICATION */
215: PetscOptionsInt("-ntests","Number of times to repeat TaoSolve","",ntests,&ntests,&flag);
216: PetscLogStageRegister("Trials",&stages[0]);
217: PetscLogStagePush(stages[0]);
218: user.ksp_its_initial = user.ksp_its;
219: ksp_old = user.ksp_its;
220: for (i=0; i<ntests; i++){
221: TaoSolve(tao);
222: PetscPrintf(PETSC_COMM_WORLD,"KSP Iterations = %D\n",user.ksp_its-ksp_old);
223: VecCopy(x0,x);
224: TaoSetInitialVector(tao,x);
225: }
226: PetscLogStagePop();
227: PetscBarrier((PetscObject)x);
228: PetscPrintf(PETSC_COMM_WORLD,"KSP iterations within initialization: ");
229: PetscPrintf(PETSC_COMM_WORLD,"%D\n",user.ksp_its_initial);
230: PetscPrintf(PETSC_COMM_WORLD,"Total KSP iterations over %D trial(s): ",ntests);
231: PetscPrintf(PETSC_COMM_WORLD,"%D\n",user.ksp_its);
232: PetscPrintf(PETSC_COMM_WORLD,"KSP iterations per trial: ");
233: PetscPrintf(PETSC_COMM_WORLD,"%D\n",(user.ksp_its-user.ksp_its_initial)/ntests);
234: TaoGetTerminationReason(tao,&reason);
236: if (reason < 0) {
237: PetscPrintf(MPI_COMM_WORLD, "TAO failed to converge.\n");
238: } else {
239: PetscPrintf(MPI_COMM_WORLD, "Optimization terminated with status %D.\n", reason);
240: }
241: TaoDestroy(&tao);
242: VecDestroy(&x);
243: VecDestroy(&x0);
244: ParabolicDestroy(&user);
246: PetscFinalize();
247: return 0;
248: }
249: /* ------------------------------------------------------------------- */
252: /*
253: dwork = Qy - d
254: lwork = L*(u-ur)
255: f = 1/2 * (dwork.dork + alpha*lwork.lwork)
256: */
257: PetscErrorCode FormFunction(Tao tao,Vec X,PetscReal *f,void *ptr)
258: {
260: PetscReal d1=0,d2=0;
261: PetscInt i,j;
262: AppCtx *user = (AppCtx*)ptr;
265: Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
266: Scatter_i(user->y,user->yi,user->yi_scatter,user->nt);
267: for (j=0; j<user->ns; j++){
268: i = user->sample_times[j];
269: MatMult(user->Qblock,user->yi[i],user->di[j]);
270: }
271: Gather_i(user->dwork,user->di,user->di_scatter,user->ns);
272: VecAXPY(user->dwork,-1.0,user->d);
273: VecDot(user->dwork,user->dwork,&d1);
275: VecWAXPY(user->uwork,-1.0,user->ur,user->u);
276: MatMult(user->L,user->uwork,user->lwork);
277: VecDot(user->lwork,user->lwork,&d2);
279: *f = 0.5 * (d1 + user->alpha*d2);
280: return(0);
281: }
283: /* ------------------------------------------------------------------- */
286: /*
287: state: g_s = Q' *(Qy - d)
288: design: g_d = alpha*L'*L*(u-ur)
289: */
290: PetscErrorCode FormGradient(Tao tao,Vec X,Vec G,void *ptr)
291: {
293: PetscInt i,j;
294: AppCtx *user = (AppCtx*)ptr;
297: Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
298: Scatter_i(user->y,user->yi,user->yi_scatter,user->nt);
299: for (j=0; j<user->ns; j++){
300: i = user->sample_times[j];
301: MatMult(user->Qblock,user->yi[i],user->di[j]);
302: }
303: Gather_i(user->dwork,user->di,user->di_scatter,user->ns);
304: VecAXPY(user->dwork,-1.0,user->d);
305: Scatter_i(user->dwork,user->di,user->di_scatter,user->ns);
306: VecSet(user->ywork,0.0);
307: Scatter_i(user->ywork,user->yiwork,user->yi_scatter,user->nt);
308: for (j=0; j<user->ns; j++){
309: i = user->sample_times[j];
310: MatMult(user->QblockT,user->di[j],user->yiwork[i]);
311: }
312: Gather_i(user->ywork,user->yiwork,user->yi_scatter,user->nt);
314: VecWAXPY(user->uwork,-1.0,user->ur,user->u);
315: MatMult(user->L,user->uwork,user->lwork);
316: MatMult(user->LT,user->lwork,user->uwork);
317: VecScale(user->uwork, user->alpha);
318: Gather(G,user->ywork,user->state_scatter,user->uwork,user->design_scatter);
319: return(0);
320: }
324: PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, void *ptr)
325: {
327: PetscReal d1,d2;
328: PetscInt i,j;
329: AppCtx *user = (AppCtx*)ptr;
332: Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
333: Scatter_i(user->y,user->yi,user->yi_scatter,user->nt);
334: for (j=0; j<user->ns; j++){
335: i = user->sample_times[j];
336: MatMult(user->Qblock,user->yi[i],user->di[j]);
337: }
338: Gather_i(user->dwork,user->di,user->di_scatter,user->ns);
339: VecAXPY(user->dwork,-1.0,user->d);
340: VecDot(user->dwork,user->dwork,&d1);
341: Scatter_i(user->dwork,user->di,user->di_scatter,user->ns);
342: VecSet(user->ywork,0.0);
343: Scatter_i(user->ywork,user->yiwork,user->yi_scatter,user->nt);
344: for (j=0; j<user->ns; j++){
345: i = user->sample_times[j];
346: MatMult(user->QblockT,user->di[j],user->yiwork[i]);
347: }
348: Gather_i(user->ywork,user->yiwork,user->yi_scatter,user->nt);
350: VecWAXPY(user->uwork,-1.0,user->ur,user->u);
351: MatMult(user->L,user->uwork,user->lwork);
352: VecDot(user->lwork,user->lwork,&d2);
353: MatMult(user->LT,user->lwork,user->uwork);
354: VecScale(user->uwork, user->alpha);
355: *f = 0.5 * (d1 + user->alpha*d2);
357: Gather(G,user->ywork,user->state_scatter,user->uwork,user->design_scatter);
358: return(0);
359: }
361: /* ------------------------------------------------------------------- */
364: /* A
365: MatShell object
366: */
367: PetscErrorCode FormJacobianState(Tao tao, Vec X, Mat *J, Mat* JPre, Mat* JInv, MatStructure* flag, void *ptr)
368: {
370: AppCtx *user = (AppCtx*)ptr;
373: Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
374: VecSet(user->uwork,0);
375: VecAXPY(user->uwork,-1.0,user->u);
376: VecExp(user->uwork);
377: MatMult(user->Av,user->uwork,user->Av_u);
378: VecCopy(user->Av_u,user->Swork);
379: VecReciprocal(user->Swork);
380: MatCopy(user->Div,user->Divwork,SAME_NONZERO_PATTERN);
381: MatDiagonalScale(user->Divwork,NULL,user->Swork);
382: if (user->dsg_formed) {
383: MatMatMultNumeric(user->Divwork,user->Grad,user->DSG);
384: } else {
385: MatMatMultSymbolic(user->Divwork,user->Grad,PETSC_DECIDE,&user->DSG);
386: MatMatMultNumeric(user->Divwork,user->Grad,user->DSG);
387: user->dsg_formed = PETSC_TRUE;
388: }
390: /* B = speye(nx^3) + ht*DSG; */
391: MatScale(user->DSG,user->ht);
392: MatShift(user->DSG,1.0);
393: *JInv = user->JsInv;
394: return(0);
395: }
397: /* ------------------------------------------------------------------- */
400: /* B */
401: PetscErrorCode FormJacobianDesign(Tao tao, Vec X, Mat *J, void *ptr)
402: {
404: AppCtx *user = (AppCtx*)ptr;
407: Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
408: *J = user->Jd;
409: return(0);
410: }
414: PetscErrorCode StateMatMult(Mat J_shell, Vec X, Vec Y)
415: {
417: PetscInt i;
418: AppCtx *user;
421: MatShellGetContext(J_shell,(void**)&user);
422: Scatter_i(X,user->yi,user->yi_scatter,user->nt);
423: MatMult(user->JsBlock,user->yi[0],user->yiwork[0]);
424: for (i=1; i<user->nt; i++){
425: MatMult(user->JsBlock,user->yi[i],user->yiwork[i]);
426: VecAXPY(user->yiwork[i],-1.0,user->yi[i-1]);
427: }
428: Gather_i(Y,user->yiwork,user->yi_scatter,user->nt);
429: return(0);
430: }
434: PetscErrorCode StateMatMultTranspose(Mat J_shell, Vec X, Vec Y)
435: {
437: PetscInt i;
438: AppCtx *user;
441: MatShellGetContext(J_shell,(void**)&user);
442: Scatter_i(X,user->yi,user->yi_scatter,user->nt);
443: for (i=0; i<user->nt-1; i++){
444: MatMult(user->JsBlock,user->yi[i],user->yiwork[i]);
445: VecAXPY(user->yiwork[i],-1.0,user->yi[i+1]);
446: }
447: i = user->nt-1;
448: MatMult(user->JsBlock,user->yi[i],user->yiwork[i]);
449: Gather_i(Y,user->yiwork,user->yi_scatter,user->nt);
450: return(0);
451: }
455: PetscErrorCode StateMatBlockMult(Mat J_shell, Vec X, Vec Y)
456: {
458: AppCtx *user;
461: MatShellGetContext(J_shell,(void**)&user);
462: MatMult(user->Grad,X,user->Swork);
463: VecPointwiseDivide(user->Swork,user->Swork,user->Av_u);
464: MatMult(user->Div,user->Swork,Y);
465: VecAYPX(Y,user->ht,X);
466: return(0);
467: }
471: PetscErrorCode DesignMatMult(Mat J_shell, Vec X, Vec Y)
472: {
474: PetscInt i;
475: AppCtx *user;
478: MatShellGetContext(J_shell,(void**)&user);
480: /* sdiag(1./v) */
481: VecSet(user->uwork,0);
482: VecAXPY(user->uwork,-1.0,user->u);
483: VecExp(user->uwork);
485: /* sdiag(1./((Av*(1./v)).^2)) */
486: MatMult(user->Av,user->uwork,user->Swork);
487: VecPointwiseMult(user->Swork,user->Swork,user->Swork);
488: VecReciprocal(user->Swork);
490: /* (Av * (sdiag(1./v) * b)) */
491: VecPointwiseMult(user->uwork,user->uwork,X);
492: MatMult(user->Av,user->uwork,user->Twork);
494: /* (sdiag(1./((Av*(1./v)).^2)) * (Av * (sdiag(1./v) * b))) */
495: VecPointwiseMult(user->Swork,user->Twork,user->Swork);
497: Scatter_i(user->y,user->yi,user->yi_scatter,user->nt);
498: for (i=0; i<user->nt; i++){
499: /* (sdiag(Grad*y(:,i)) */
500: MatMult(user->Grad,user->yi[i],user->Twork);
502: /* ht * Div * (sdiag(Grad*y(:,i)) * (sdiag(1./((Av*(1./v)).^2)) * (Av * (sdiag(1./v) * b)))) */
503: VecPointwiseMult(user->Twork,user->Twork,user->Swork);
504: MatMult(user->Div,user->Twork,user->yiwork[i]);
505: VecScale(user->yiwork[i],user->ht);
506: }
507: Gather_i(Y,user->yiwork,user->yi_scatter,user->nt);
509: return(0);
510: }
514: PetscErrorCode DesignMatMultTranspose(Mat J_shell, Vec X, Vec Y)
515: {
517: PetscInt i;
518: AppCtx *user;
521: MatShellGetContext(J_shell,(void**)&user);
523: /* sdiag(1./((Av*(1./v)).^2)) */
524: VecSet(user->uwork,0);
525: VecAXPY(user->uwork,-1.0,user->u);
526: VecExp(user->uwork);
527: MatMult(user->Av,user->uwork,user->Rwork);
528: VecPointwiseMult(user->Rwork,user->Rwork,user->Rwork);
529: VecReciprocal(user->Rwork);
531: VecSet(Y,0.0);
532: Scatter_i(user->y,user->yi,user->yi_scatter,user->nt);
533: Scatter_i(X,user->yiwork,user->yi_scatter,user->nt);
534: for (i=0; i<user->nt; i++){
535: /* (Div' * b(:,i)) */
536: MatMult(user->Grad,user->yiwork[i],user->Swork);
538: /* sdiag(Grad*y(:,i)) */
539: MatMult(user->Grad,user->yi[i],user->Twork);
541: /* (sdiag(Grad*y(:,i)) * (Div' * b(:,i))) */
542: VecPointwiseMult(user->Twork,user->Swork,user->Twork);
544: /* (sdiag(1./((Av*(1./v)).^2)) * (sdiag(Grad*y(:,i)) * (Div' * b(:,i)))) */
545: VecPointwiseMult(user->Twork,user->Rwork,user->Twork);
547: /* (Av' * (sdiag(1./((Av*(1./v)).^2)) * (sdiag(Grad*y(:,i)) * (Div' * b(:,i))))) */
548: MatMult(user->AvT,user->Twork,user->yiwork[i]);
550: /* sdiag(1./v) * (Av' * (sdiag(1./((Av*(1./v)).^2)) * (sdiag(Grad*y(:,i)) * (Div' * b(:,i))))) */
551: VecPointwiseMult(user->yiwork[i],user->uwork,user->yiwork[i]);
552: VecAXPY(Y,user->ht,user->yiwork[i]);
553: }
554: return(0);
555: }
559: PetscErrorCode StateMatBlockPrecMult(PC PC_shell, Vec X, Vec Y)
560: {
562: AppCtx *user;
565: PCShellGetContext(PC_shell,(void**)&user);
567: if (user->dsg_formed) {
568: MatSOR(user->DSG,X,1.0,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_SYMMETRIC_SWEEP),0.0,1,1,Y);
569: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"DSG not formed");
570: return(0);
571: }
575: PetscErrorCode StateMatInvMult(Mat J_shell, Vec X, Vec Y)
576: {
578: AppCtx *user;
579: PetscReal tau;
580: PetscInt its,i;
583: MatShellGetContext(J_shell,(void**)&user);
585: if (Y == user->ytrue) {
586: KSPSetTolerances(user->solver,1e-8,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
587: } else if (user->lcl) {
588: tau = user->lcl->tau[user->lcl->solve_type];
589: KSPSetTolerances(user->solver,tau,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
590: }
592: Scatter_i(X,user->yi,user->yi_scatter,user->nt);
593: KSPSolve(user->solver,user->yi[0],user->yiwork[0]);
594: KSPGetIterationNumber(user->solver,&its);
595: user->ksp_its = user->ksp_its + its;
597: for (i=1; i<user->nt; i++){
598: VecAXPY(user->yi[i],1.0,user->yiwork[i-1]);
599: KSPSolve(user->solver,user->yi[i],user->yiwork[i]);
600: KSPGetIterationNumber(user->solver,&its);
601: user->ksp_its = user->ksp_its + its;
602: }
603: Gather_i(Y,user->yiwork,user->yi_scatter,user->nt);
604: return(0);
605: }
609: PetscErrorCode StateMatInvTransposeMult(Mat J_shell, Vec X, Vec Y)
610: {
612: PetscReal tau;
613: AppCtx *user;
614: PetscInt its,i;
617: MatShellGetContext(J_shell,(void**)&user);
618: if (user->lcl) {
619: tau = user->lcl->tau[user->lcl->solve_type];
620: KSPSetTolerances(user->solver,tau,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
621: }
622: Scatter_i(X,user->yi,user->yi_scatter,user->nt);
624: i = user->nt - 1;
625: KSPSolve(user->solver,user->yi[i],user->yiwork[i]);
627: KSPGetIterationNumber(user->solver,&its);
628: user->ksp_its = user->ksp_its + its;
630: for (i=user->nt-2; i>=0; i--){
631: VecAXPY(user->yi[i],1.0,user->yiwork[i+1]);
632: KSPSolve(user->solver,user->yi[i],user->yiwork[i]);
634: KSPGetIterationNumber(user->solver,&its);
635: user->ksp_its = user->ksp_its + its;
637: }
639: Gather_i(Y,user->yiwork,user->yi_scatter,user->nt);
640: return(0);
641: }
645: PetscErrorCode StateMatDuplicate(Mat J_shell, MatDuplicateOption opt, Mat *new_shell)
646: {
648: AppCtx *user;
651: MatShellGetContext(J_shell,(void**)&user);
653: MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m,user->m,user,new_shell);
654: MatShellSetOperation(*new_shell,MATOP_MULT,(void(*)(void))StateMatMult);
655: MatShellSetOperation(*new_shell,MATOP_DUPLICATE,(void(*)(void))StateMatDuplicate);
656: MatShellSetOperation(*new_shell,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatMultTranspose);
657: MatShellSetOperation(*new_shell,MATOP_GET_DIAGONAL,(void(*)(void))StateMatGetDiagonal);
658: return(0);
659: }
663: PetscErrorCode StateMatGetDiagonal(Mat J_shell, Vec X)
664: {
666: AppCtx *user;
669: MatShellGetContext(J_shell,(void**)&user);
670: VecCopy(user->js_diag,X);
671: return(0);
673: }
677: PetscErrorCode FormConstraints(Tao tao, Vec X, Vec C, void *ptr)
678: {
679: /* con = Ay - q, A = [B 0 0 ... 0;
680: -I B 0 ... 0;
681: 0 -I B ... 0;
682: ... ;
683: 0 ... -I B]
684: B = ht * Div * Sigma * Grad + eye */
686: PetscInt i;
687: AppCtx *user = (AppCtx*)ptr;
690: Scatter(X,user->y,user->state_scatter,user->u,user->design_scatter);
691: Scatter_i(user->y,user->yi,user->yi_scatter,user->nt);
692: MatMult(user->JsBlock,user->yi[0],user->yiwork[0]);
693: for (i=1; i<user->nt; i++){
694: MatMult(user->JsBlock,user->yi[i],user->yiwork[i]);
695: VecAXPY(user->yiwork[i],-1.0,user->yi[i-1]);
696: }
697: Gather_i(C,user->yiwork,user->yi_scatter,user->nt);
698: VecAXPY(C,-1.0,user->q);
699: return(0);
700: }
705: PetscErrorCode Scatter(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat)
706: {
710: VecScatterBegin(s_scat,x,state,INSERT_VALUES,SCATTER_FORWARD);
711: VecScatterEnd(s_scat,x,state,INSERT_VALUES,SCATTER_FORWARD);
712: VecScatterBegin(d_scat,x,design,INSERT_VALUES,SCATTER_FORWARD);
713: VecScatterEnd(d_scat,x,design,INSERT_VALUES,SCATTER_FORWARD);
714: return(0);
715: }
719: PetscErrorCode Scatter_i(Vec y, Vec *yi, VecScatter *scat, PetscInt nt)
720: {
722: PetscInt i;
725: for (i=0; i<nt; i++){
726: VecScatterBegin(scat[i],y,yi[i],INSERT_VALUES,SCATTER_FORWARD);
727: VecScatterEnd(scat[i],y,yi[i],INSERT_VALUES,SCATTER_FORWARD);
728: }
729: return(0);
730: }
735: PetscErrorCode Gather(Vec x, Vec state, VecScatter s_scat, Vec design, VecScatter d_scat)
736: {
740: VecScatterBegin(s_scat,state,x,INSERT_VALUES,SCATTER_REVERSE);
741: VecScatterEnd(s_scat,state,x,INSERT_VALUES,SCATTER_REVERSE);
742: VecScatterBegin(d_scat,design,x,INSERT_VALUES,SCATTER_REVERSE);
743: VecScatterEnd(d_scat,design,x,INSERT_VALUES,SCATTER_REVERSE);
744: return(0);
745: }
749: PetscErrorCode Gather_i(Vec y, Vec *yi, VecScatter *scat, PetscInt nt)
750: {
752: PetscInt i;
755: for (i=0; i<nt; i++){
756: VecScatterBegin(scat[i],yi[i],y,INSERT_VALUES,SCATTER_REVERSE);
757: VecScatterEnd(scat[i],yi[i],y,INSERT_VALUES,SCATTER_REVERSE);
758: }
759: return(0);
760: }
764: PetscErrorCode ParabolicInitialize(AppCtx *user)
765: {
767: PetscInt m,n,i,j,k,linear_index,istart,iend,iblock,lo,hi,lo2,hi2;
768: Vec XX,YY,ZZ,XXwork,YYwork,ZZwork,UTwork,yi,di,bc;
769: PetscReal *x, *y, *z;
770: PetscReal h,stime;
771: PetscScalar hinv,neg_hinv,half = 0.5,sqrt_beta;
772: PetscInt im,indx1,indx2,indy1,indy2,indz1,indz2,nx,ny,nz;
773: PetscReal xri,yri,zri,xim,yim,zim,dx1,dx2,dy1,dy2,dz1,dz2,Dx,Dy,Dz;
774: PetscScalar v,vx,vy,vz;
775: IS is_from_y,is_to_yi,is_from_d,is_to_di;
776: /* Data locations */
777: PetscScalar xr[64] = {0.4970, 0.8498, 0.7814, 0.6268, 0.7782, 0.6402, 0.3617, 0.3160,
778: 0.3610, 0.5298, 0.6987, 0.3331, 0.7962, 0.5596, 0.3866, 0.6774,
779: 0.5407, 0.4518, 0.6702, 0.6061, 0.7580, 0.8997, 0.5198, 0.8326,
780: 0.2138, 0.9198, 0.3000, 0.2833, 0.8288, 0.7076, 0.1820, 0.0728,
781: 0.8447, 0.2367, 0.3239, 0.6413, 0.3114, 0.4731, 0.1192, 0.9273,
782: 0.5724, 0.4331, 0.5136, 0.3547, 0.4413, 0.2602, 0.5698, 0.7278,
783: 0.5261, 0.6230, 0.2454, 0.3948, 0.7479, 0.6582, 0.4660, 0.5594,
784: 0.7574, 0.1143, 0.5900, 0.1065, 0.4260, 0.3294, 0.8276, 0.0756};
786: PetscScalar yr[64] = {0.7345, 0.9120, 0.9288, 0.7528, 0.4463, 0.4985, 0.2497, 0.6256,
787: 0.3425, 0.9026, 0.6983, 0.4230, 0.7140, 0.2970, 0.4474, 0.8792,
788: 0.6604, 0.2485, 0.7968, 0.6127, 0.1796, 0.2437, 0.5938, 0.6137,
789: 0.3867, 0.5658, 0.4575, 0.1009, 0.0863, 0.3361, 0.0738, 0.3985,
790: 0.6602, 0.1437, 0.0934, 0.5983, 0.5950, 0.0763, 0.0768, 0.2288,
791: 0.5761, 0.1129, 0.3841, 0.6150, 0.6904, 0.6686, 0.1361, 0.4601,
792: 0.4491, 0.3716, 0.1969, 0.6537, 0.6743, 0.6991, 0.4811, 0.5480,
793: 0.1684, 0.4569, 0.6889, 0.8437, 0.3015, 0.2854, 0.8199, 0.2658};
795: PetscScalar zr[64] = {0.7668, 0.8573, 0.2654, 0.2719, 0.1060, 0.1311, 0.6232, 0.2295,
796: 0.8009, 0.2147, 0.2119, 0.9325, 0.4473, 0.3600, 0.3374, 0.3819,
797: 0.4066, 0.5801, 0.1673, 0.0959, 0.4638, 0.8236, 0.8800, 0.2939,
798: 0.2028, 0.8262, 0.2706, 0.6276, 0.9085, 0.6443, 0.8241, 0.0712,
799: 0.1824, 0.7789, 0.4389, 0.8415, 0.7055, 0.6639, 0.3653, 0.2078,
800: 0.1987, 0.2297, 0.4321, 0.8115, 0.4915, 0.7764, 0.4657, 0.4627,
801: 0.4569, 0.4232, 0.8514, 0.0674, 0.3227, 0.1055, 0.6690, 0.6313,
802: 0.9226, 0.5461, 0.4126, 0.2364, 0.6096, 0.7042, 0.3914, 0.0711};
805: PetscMalloc1(user->mx,&x);
806: PetscMalloc1(user->mx,&y);
807: PetscMalloc1(user->mx,&z);
808: user->jformed = PETSC_FALSE;
809: user->dsg_formed = PETSC_FALSE;
811: n = user->mx * user->mx * user->mx;
812: m = 3 * user->mx * user->mx * (user->mx-1);
813: sqrt_beta = PetscSqrtScalar(user->beta);
815: user->ksp_its = 0;
816: user->ksp_its_initial = 0;
818: stime = (PetscReal)user->nt/user->ns;
819: PetscMalloc1(user->ns,&user->sample_times);
820: for (i=0; i<user->ns; i++){
821: user->sample_times[i] = (PetscInt)(stime*((PetscReal)i+1.0)-0.5);
822: }
824: VecCreate(PETSC_COMM_WORLD,&XX);
825: VecCreate(PETSC_COMM_WORLD,&user->q);
826: VecSetSizes(XX,PETSC_DECIDE,n);
827: VecSetSizes(user->q,PETSC_DECIDE,n*user->nt);
828: VecSetFromOptions(XX);
829: VecSetFromOptions(user->q);
831: VecDuplicate(XX,&YY);
832: VecDuplicate(XX,&ZZ);
833: VecDuplicate(XX,&XXwork);
834: VecDuplicate(XX,&YYwork);
835: VecDuplicate(XX,&ZZwork);
836: VecDuplicate(XX,&UTwork);
837: VecDuplicate(XX,&user->utrue);
838: VecDuplicate(XX,&bc);
840: /* Generate 3D grid, and collect ns (1<=ns<=8) right-hand-side vectors into user->q */
841: h = 1.0/user->mx;
842: hinv = user->mx;
843: neg_hinv = -hinv;
845: VecGetOwnershipRange(XX,&istart,&iend);
846: for (linear_index=istart; linear_index<iend; linear_index++){
847: i = linear_index % user->mx;
848: j = ((linear_index-i)/user->mx) % user->mx;
849: k = ((linear_index-i)/user->mx-j) / user->mx;
850: vx = h*(i+0.5);
851: vy = h*(j+0.5);
852: vz = h*(k+0.5);
853: VecSetValues(XX,1,&linear_index,&vx,INSERT_VALUES);
854: VecSetValues(YY,1,&linear_index,&vy,INSERT_VALUES);
855: VecSetValues(ZZ,1,&linear_index,&vz,INSERT_VALUES);
856: if ((vx<0.6) && (vx>0.4) && (vy<0.6) && (vy>0.4) && (vy<0.6) && (vz<0.6) && (vz>0.4)){
857: v = 1000.0;
858: VecSetValues(bc,1,&linear_index,&v,INSERT_VALUES);
859: }
860: }
862: VecAssemblyBegin(XX);
863: VecAssemblyEnd(XX);
864: VecAssemblyBegin(YY);
865: VecAssemblyEnd(YY);
866: VecAssemblyBegin(ZZ);
867: VecAssemblyEnd(ZZ);
868: VecAssemblyBegin(bc);
869: VecAssemblyEnd(bc);
871: /* Compute true parameter function
872: ut = 0.5 + exp(-10*((x-0.5)^2+(y-0.5)^2+(z-0.5)^2)) */
873: VecCopy(XX,XXwork);
874: VecCopy(YY,YYwork);
875: VecCopy(ZZ,ZZwork);
877: VecShift(XXwork,-0.5);
878: VecShift(YYwork,-0.5);
879: VecShift(ZZwork,-0.5);
881: VecPointwiseMult(XXwork,XXwork,XXwork);
882: VecPointwiseMult(YYwork,YYwork,YYwork);
883: VecPointwiseMult(ZZwork,ZZwork,ZZwork);
885: VecCopy(XXwork,user->utrue);
886: VecAXPY(user->utrue,1.0,YYwork);
887: VecAXPY(user->utrue,1.0,ZZwork);
888: VecScale(user->utrue,-10.0);
889: VecExp(user->utrue);
890: VecShift(user->utrue,0.5);
892: VecDestroy(&XX);
893: VecDestroy(&YY);
894: VecDestroy(&ZZ);
895: VecDestroy(&XXwork);
896: VecDestroy(&YYwork);
897: VecDestroy(&ZZwork);
898: VecDestroy(&UTwork);
900: /* Initial guess and reference model */
901: VecDuplicate(user->utrue,&user->ur);
902: VecSet(user->ur,0.0);
904: /* Generate Grad matrix */
905: MatCreate(PETSC_COMM_WORLD,&user->Grad);
906: MatSetSizes(user->Grad,PETSC_DECIDE,PETSC_DECIDE,m,n);
907: MatSetFromOptions(user->Grad);
908: MatMPIAIJSetPreallocation(user->Grad,2,NULL,2,NULL);
909: MatSeqAIJSetPreallocation(user->Grad,2,NULL);
910: MatGetOwnershipRange(user->Grad,&istart,&iend);
912: for (i=istart; i<iend; i++){
913: if (i<m/3){
914: iblock = i / (user->mx-1);
915: j = iblock*user->mx + (i % (user->mx-1));
916: MatSetValues(user->Grad,1,&i,1,&j,&neg_hinv,INSERT_VALUES);
917: j = j+1;
918: MatSetValues(user->Grad,1,&i,1,&j,&hinv,INSERT_VALUES);
919: }
920: if (i>=m/3 && i<2*m/3){
921: iblock = (i-m/3) / (user->mx*(user->mx-1));
922: j = iblock*user->mx*user->mx + ((i-m/3) % (user->mx*(user->mx-1)));
923: MatSetValues(user->Grad,1,&i,1,&j,&neg_hinv,INSERT_VALUES);
924: j = j + user->mx;
925: MatSetValues(user->Grad,1,&i,1,&j,&hinv,INSERT_VALUES);
926: }
927: if (i>=2*m/3){
928: j = i-2*m/3;
929: MatSetValues(user->Grad,1,&i,1,&j,&neg_hinv,INSERT_VALUES);
930: j = j + user->mx*user->mx;
931: MatSetValues(user->Grad,1,&i,1,&j,&hinv,INSERT_VALUES);
932: }
933: }
935: MatAssemblyBegin(user->Grad,MAT_FINAL_ASSEMBLY);
936: MatAssemblyEnd(user->Grad,MAT_FINAL_ASSEMBLY);
939: /* Generate arithmetic averaging matrix Av */
940: MatCreate(PETSC_COMM_WORLD,&user->Av);
941: MatSetSizes(user->Av,PETSC_DECIDE,PETSC_DECIDE,m,n);
942: MatSetFromOptions(user->Av);
943: MatMPIAIJSetPreallocation(user->Av,2,NULL,2,NULL);
944: MatSeqAIJSetPreallocation(user->Av,2,NULL);
945: MatGetOwnershipRange(user->Av,&istart,&iend);
947: for (i=istart; i<iend; i++){
948: if (i<m/3){
949: iblock = i / (user->mx-1);
950: j = iblock*user->mx + (i % (user->mx-1));
951: MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES);
952: j = j+1;
953: MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES);
954: }
955: if (i>=m/3 && i<2*m/3){
956: iblock = (i-m/3) / (user->mx*(user->mx-1));
957: j = iblock*user->mx*user->mx + ((i-m/3) % (user->mx*(user->mx-1)));
958: MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES);
959: j = j + user->mx;
960: MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES);
961: }
962: if (i>=2*m/3){
963: j = i-2*m/3;
964: MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES);
965: j = j + user->mx*user->mx;
966: MatSetValues(user->Av,1,&i,1,&j,&half,INSERT_VALUES);
967: }
968: }
970: MatAssemblyBegin(user->Av,MAT_FINAL_ASSEMBLY);
971: MatAssemblyEnd(user->Av,MAT_FINAL_ASSEMBLY);
973: /* Generate transpose of averaging matrix Av */
974: MatTranspose(user->Av,MAT_INITIAL_MATRIX,&user->AvT);
976: MatCreate(PETSC_COMM_WORLD,&user->L);
977: MatSetSizes(user->L,PETSC_DECIDE,PETSC_DECIDE,m+n,n);
978: MatSetFromOptions(user->L);
979: MatMPIAIJSetPreallocation(user->L,2,NULL,2,NULL);
980: MatSeqAIJSetPreallocation(user->L,2,NULL);
981: MatGetOwnershipRange(user->L,&istart,&iend);
983: for (i=istart; i<iend; i++){
984: if (i<m/3){
985: iblock = i / (user->mx-1);
986: j = iblock*user->mx + (i % (user->mx-1));
987: MatSetValues(user->L,1,&i,1,&j,&neg_hinv,INSERT_VALUES);
988: j = j+1;
989: MatSetValues(user->L,1,&i,1,&j,&hinv,INSERT_VALUES);
990: }
991: if (i>=m/3 && i<2*m/3){
992: iblock = (i-m/3) / (user->mx*(user->mx-1));
993: j = iblock*user->mx*user->mx + ((i-m/3) % (user->mx*(user->mx-1)));
994: MatSetValues(user->L,1,&i,1,&j,&neg_hinv,INSERT_VALUES);
995: j = j + user->mx;
996: MatSetValues(user->L,1,&i,1,&j,&hinv,INSERT_VALUES);
997: }
998: if (i>=2*m/3 && i<m){
999: j = i-2*m/3;
1000: MatSetValues(user->L,1,&i,1,&j,&neg_hinv,INSERT_VALUES);
1001: j = j + user->mx*user->mx;
1002: MatSetValues(user->L,1,&i,1,&j,&hinv,INSERT_VALUES);
1003: }
1004: if (i>=m){
1005: j = i - m;
1006: MatSetValues(user->L,1,&i,1,&j,&sqrt_beta,INSERT_VALUES);
1007: }
1008: }
1010: MatAssemblyBegin(user->L,MAT_FINAL_ASSEMBLY);
1011: MatAssemblyEnd(user->L,MAT_FINAL_ASSEMBLY);
1013: MatScale(user->L,PetscPowScalar(h,1.5));
1015: /* Generate Div matrix */
1016: MatTranspose(user->Grad,MAT_INITIAL_MATRIX,&user->Div);
1018: /* Build work vectors and matrices */
1019: VecCreate(PETSC_COMM_WORLD,&user->S);
1020: VecSetSizes(user->S, PETSC_DECIDE, user->mx*user->mx*(user->mx-1)*3);
1021: VecSetFromOptions(user->S);
1023: VecCreate(PETSC_COMM_WORLD,&user->lwork);
1024: VecSetSizes(user->lwork,PETSC_DECIDE,m+user->mx*user->mx*user->mx);
1025: VecSetFromOptions(user->lwork);
1027: MatDuplicate(user->Div,MAT_SHARE_NONZERO_PATTERN,&user->Divwork);
1028: MatDuplicate(user->Av,MAT_SHARE_NONZERO_PATTERN,&user->Avwork);
1030: VecCreate(PETSC_COMM_WORLD,&user->d);
1031: VecSetSizes(user->d,PETSC_DECIDE,user->ndata*user->nt);
1032: VecSetFromOptions(user->d);
1034: VecDuplicate(user->S,&user->Swork);
1035: VecDuplicate(user->S,&user->Av_u);
1036: VecDuplicate(user->S,&user->Twork);
1037: VecDuplicate(user->S,&user->Rwork);
1038: VecDuplicate(user->y,&user->ywork);
1039: VecDuplicate(user->u,&user->uwork);
1040: VecDuplicate(user->u,&user->js_diag);
1041: VecDuplicate(user->c,&user->cwork);
1042: VecDuplicate(user->d,&user->dwork);
1044: /* Create matrix-free shell user->Js for computing A*x */
1045: MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m*user->nt,user->m*user->nt,user,&user->Js);
1046: MatShellSetOperation(user->Js,MATOP_MULT,(void(*)(void))StateMatMult);
1047: MatShellSetOperation(user->Js,MATOP_DUPLICATE,(void(*)(void))StateMatDuplicate);
1048: MatShellSetOperation(user->Js,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatMultTranspose);
1049: MatShellSetOperation(user->Js,MATOP_GET_DIAGONAL,(void(*)(void))StateMatGetDiagonal);
1051: /* Diagonal blocks of user->Js */
1052: MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m,user->m,user,&user->JsBlock);
1053: MatShellSetOperation(user->JsBlock,MATOP_MULT,(void(*)(void))StateMatBlockMult);
1054: /* Blocks are symmetric */
1055: MatShellSetOperation(user->JsBlock,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatBlockMult);
1057: /* Create a matrix-free shell user->JsBlockPrec for computing (U+D)\D*(L+D)\x, where JsBlock = L+D+U,
1058: D is diagonal, L is strictly lower triangular, and U is strictly upper triangular.
1059: This is an SSOR preconditioner for user->JsBlock. */
1060: MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m,user->m,user,&user->JsBlockPrec);
1061: MatShellSetOperation(user->JsBlockPrec,MATOP_MULT,(void(*)(void))StateMatBlockPrecMult);
1062: /* JsBlockPrec is symmetric */
1063: MatShellSetOperation(user->JsBlockPrec,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatBlockPrecMult);
1064: MatSetOption(user->JsBlockPrec,MAT_SYMMETRY_ETERNAL,PETSC_TRUE);
1066: /* Create a matrix-free shell user->Jd for computing B*x */
1067: MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m*user->nt,user->m,user,&user->Jd);
1068: MatShellSetOperation(user->Jd,MATOP_MULT,(void(*)(void))DesignMatMult);
1069: MatShellSetOperation(user->Jd,MATOP_MULT_TRANSPOSE,(void(*)(void))DesignMatMultTranspose);
1071: /* User-defined routines for computing user->Js\x and user->Js^T\x*/
1072: MatCreateShell(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,user->m*user->nt,user->m*user->nt,user,&user->JsInv);
1073: MatShellSetOperation(user->JsInv,MATOP_MULT,(void(*)(void))StateMatInvMult);
1074: MatShellSetOperation(user->JsInv,MATOP_MULT_TRANSPOSE,(void(*)(void))StateMatInvTransposeMult);
1076: /* Solver options and tolerances */
1077: KSPCreate(PETSC_COMM_WORLD,&user->solver);
1078: KSPSetType(user->solver,KSPCG);
1079: KSPSetOperators(user->solver,user->JsBlock,user->JsBlockPrec,DIFFERENT_NONZERO_PATTERN);
1080: KSPSetInitialGuessNonzero(user->solver,PETSC_FALSE);
1081: KSPSetTolerances(user->solver,1e-4,1e-20,1e3,500);
1082: KSPGetPC(user->solver,&user->prec);
1083: PCSetType(user->prec,PCSHELL);
1085: PCShellSetApply(user->prec,StateMatBlockPrecMult);
1086: PCShellSetApplyTranspose(user->prec,StateMatBlockPrecMult);
1087: PCShellSetContext(user->prec,user);
1089: /* Create scatter from y to y_1,y_2,...,y_nt */
1090: PetscMalloc1(user->nt*user->m,&user->yi_scatter);
1091: VecCreate(PETSC_COMM_WORLD,&yi);
1092: VecSetSizes(yi,PETSC_DECIDE,user->mx*user->mx*user->mx);
1093: VecSetFromOptions(yi);
1094: VecDuplicateVecs(yi,user->nt,&user->yi);
1095: VecDuplicateVecs(yi,user->nt,&user->yiwork);
1097: VecGetOwnershipRange(user->y,&lo2,&hi2);
1098: istart = 0;
1099: for (i=0; i<user->nt; i++){
1100: VecGetOwnershipRange(user->yi[i],&lo,&hi);
1101: ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_to_yi);
1102: ISCreateStride(PETSC_COMM_SELF,hi-lo,lo2+istart,1,&is_from_y);
1103: VecScatterCreate(user->y,is_from_y,user->yi[i],is_to_yi,&user->yi_scatter[i]);
1104: istart = istart + hi-lo;
1105: ISDestroy(&is_to_yi);
1106: ISDestroy(&is_from_y);
1107: }
1108: VecDestroy(&yi);
1110: /* Create scatter from d to d_1,d_2,...,d_ns */
1111: PetscMalloc1(user->ns*user->ndata,&user->di_scatter);
1112: VecCreate(PETSC_COMM_WORLD,&di);
1113: VecSetSizes(di,PETSC_DECIDE,user->ndata);
1114: VecSetFromOptions(di);
1115: VecDuplicateVecs(di,user->ns,&user->di);
1116: VecGetOwnershipRange(user->d,&lo2,&hi2);
1117: istart = 0;
1118: for (i=0; i<user->ns; i++){
1119: VecGetOwnershipRange(user->di[i],&lo,&hi);
1120: ISCreateStride(PETSC_COMM_SELF,hi-lo,lo,1,&is_to_di);
1121: ISCreateStride(PETSC_COMM_SELF,hi-lo,lo2+istart,1,&is_from_d);
1122: VecScatterCreate(user->d,is_from_d,user->di[i],is_to_di,&user->di_scatter[i]);
1123: istart = istart + hi-lo;
1124: ISDestroy(&is_to_di);
1125: ISDestroy(&is_from_d);
1126: }
1127: VecDestroy(&di);
1129: /* Assemble RHS of forward problem */
1130: VecCopy(bc,user->yiwork[0]);
1131: for (i=1; i<user->nt; i++){
1132: VecSet(user->yiwork[i],0.0);
1133: }
1134: Gather_i(user->q,user->yiwork,user->yi_scatter,user->nt);
1135: VecDestroy(&bc);
1137: /* Compute true state function yt given ut */
1138: VecCreate(PETSC_COMM_WORLD,&user->ytrue);
1139: VecSetSizes(user->ytrue,PETSC_DECIDE,n*user->nt);
1140: VecSetFromOptions(user->ytrue);
1142: /* First compute Av_u = Av*exp(-u) */
1143: VecSet(user->uwork,0);
1144: VecAXPY(user->uwork,-1.0,user->utrue); /* Note: user->utrue */
1145: VecExp(user->uwork);
1146: MatMult(user->Av,user->uwork,user->Av_u);
1148: MatMatMultSymbolic(user->Div,user->Grad,PETSC_DEFAULT,&user->DSG);
1149: user->dsg_formed = PETSC_TRUE;
1151: /* Next form DSG = Div*S*Grad */
1152: MatCopy(user->Div,user->Divwork,SAME_NONZERO_PATTERN);
1153: MatDiagonalScale(user->Divwork,NULL,user->Av_u);
1154: if (user->dsg_formed) {
1155: MatMatMultNumeric(user->Divwork,user->Grad,user->DSG);
1156: } else {
1157: MatMatMultSymbolic(user->Div,user->Grad,PETSC_DEFAULT,&user->DSG);
1158: MatMatMultNumeric(user->Div,user->Grad,user->DSG);
1159: user->dsg_formed = PETSC_TRUE;
1160: }
1161: /* B = speye(nx^3) + ht*DSG; */
1162: MatScale(user->DSG,user->ht);
1163: MatShift(user->DSG,1.0);
1165: /* Now solve for ytrue */
1166: StateMatInvMult(user->Js,user->q,user->ytrue);
1168: /* Initial guess y0 for state given u0 */
1170: /* First compute Av_u = Av*exp(-u) */
1171: VecSet(user->uwork,0);
1172: VecAXPY(user->uwork,-1.0,user->u); /* Note: user->u */
1173: VecExp(user->uwork);
1174: MatMult(user->Av,user->uwork,user->Av_u);
1176: /* Next form DSG = Div*S*Grad */
1177: MatCopy(user->Div,user->Divwork,SAME_NONZERO_PATTERN);
1178: MatDiagonalScale(user->Divwork,NULL,user->Av_u);
1179: if (user->dsg_formed) {
1180: MatMatMultNumeric(user->Divwork,user->Grad,user->DSG);
1181: } else {
1182: MatMatMultSymbolic(user->Div,user->Grad,PETSC_DEFAULT,&user->DSG);
1183: MatMatMultNumeric(user->Div,user->Grad,user->DSG);
1184: user->dsg_formed = PETSC_TRUE;
1185: }
1186: /* B = speye(nx^3) + ht*DSG; */
1187: MatScale(user->DSG,user->ht);
1188: MatShift(user->DSG,1.0);
1190: /* Now solve for y */
1191: user->lcl->solve_type = LCL_FORWARD1;
1192: StateMatInvMult(user->Js,user->q,user->y);
1194: /* Construct projection matrix Q, a block diagonal matrix consisting of nt copies of Qblock along the diagonal */
1195: MatCreate(PETSC_COMM_WORLD,&user->Qblock);
1196: MatSetSizes(user->Qblock,PETSC_DECIDE,PETSC_DECIDE,user->ndata,n);
1197: MatSetFromOptions(user->Qblock);
1198: MatMPIAIJSetPreallocation(user->Qblock,8,NULL,8,NULL);
1199: MatSeqAIJSetPreallocation(user->Qblock,8,NULL);
1201: for (i=0; i<user->mx; i++){
1202: x[i] = h*(i+0.5);
1203: y[i] = h*(i+0.5);
1204: z[i] = h*(i+0.5);
1205: }
1207: MatGetOwnershipRange(user->Qblock,&istart,&iend);
1208: nx = user->mx; ny = user->mx; nz = user->mx;
1209: for (i=istart; i<iend; i++){
1210: xri = xr[i];
1211: im = 0;
1212: xim = x[im];
1213: while (xri>xim && im<nx){
1214: im = im+1;
1215: xim = x[im];
1216: }
1217: indx1 = im-1;
1218: indx2 = im;
1219: dx1 = xri - x[indx1];
1220: dx2 = x[indx2] - xri;
1222: yri = yr[i];
1223: im = 0;
1224: yim = y[im];
1225: while (yri>yim && im<ny){
1226: im = im+1;
1227: yim = y[im];
1228: }
1229: indy1 = im-1;
1230: indy2 = im;
1231: dy1 = yri - y[indy1];
1232: dy2 = y[indy2] - yri;
1234: zri = zr[i];
1235: im = 0;
1236: zim = z[im];
1237: while (zri>zim && im<nz){
1238: im = im+1;
1239: zim = z[im];
1240: }
1241: indz1 = im-1;
1242: indz2 = im;
1243: dz1 = zri - z[indz1];
1244: dz2 = z[indz2] - zri;
1246: Dx = x[indx2] - x[indx1];
1247: Dy = y[indy2] - y[indy1];
1248: Dz = z[indz2] - z[indz1];
1250: j = indx1 + indy1*nx + indz1*nx*ny;
1251: v = (1-dx1/Dx)*(1-dy1/Dy)*(1-dz1/Dz);
1252: MatSetValues(user->Qblock,1,&i,1,&j,&v,INSERT_VALUES);
1254: j = indx1 + indy1*nx + indz2*nx*ny;
1255: v = (1-dx1/Dx)*(1-dy1/Dy)*(1-dz2/Dz);
1256: MatSetValues(user->Qblock,1,&i,1,&j,&v,INSERT_VALUES);
1258: j = indx1 + indy2*nx + indz1*nx*ny;
1259: v = (1-dx1/Dx)*(1-dy2/Dy)*(1-dz1/Dz);
1260: MatSetValues(user->Qblock,1,&i,1,&j,&v,INSERT_VALUES);
1262: j = indx1 + indy2*nx + indz2*nx*ny;
1263: v = (1-dx1/Dx)*(1-dy2/Dy)*(1-dz2/Dz);
1264: MatSetValues(user->Qblock,1,&i,1,&j,&v,INSERT_VALUES);
1266: j = indx2 + indy1*nx + indz1*nx*ny;
1267: v = (1-dx2/Dx)*(1-dy1/Dy)*(1-dz1/Dz);
1268: MatSetValues(user->Qblock,1,&i,1,&j,&v,INSERT_VALUES);
1270: j = indx2 + indy1*nx + indz2*nx*ny;
1271: v = (1-dx2/Dx)*(1-dy1/Dy)*(1-dz2/Dz);
1272: MatSetValues(user->Qblock,1,&i,1,&j,&v,INSERT_VALUES);
1274: j = indx2 + indy2*nx + indz1*nx*ny;
1275: v = (1-dx2/Dx)*(1-dy2/Dy)*(1-dz1/Dz);
1276: MatSetValues(user->Qblock,1,&i,1,&j,&v,INSERT_VALUES);
1278: j = indx2 + indy2*nx + indz2*nx*ny;
1279: v = (1-dx2/Dx)*(1-dy2/Dy)*(1-dz2/Dz);
1280: MatSetValues(user->Qblock,1,&i,1,&j,&v,INSERT_VALUES);
1281: }
1282: MatAssemblyBegin(user->Qblock,MAT_FINAL_ASSEMBLY);
1283: MatAssemblyEnd(user->Qblock,MAT_FINAL_ASSEMBLY);
1285: MatTranspose(user->Qblock,MAT_INITIAL_MATRIX,&user->QblockT);
1286: MatTranspose(user->L,MAT_INITIAL_MATRIX,&user->LT);
1288: /* Add noise to the measurement data */
1289: VecSet(user->ywork,1.0);
1290: VecAYPX(user->ywork,user->noise,user->ytrue);
1291: Scatter_i(user->ywork,user->yiwork,user->yi_scatter,user->nt);
1292: for (j=0; j<user->ns; j++){
1293: i = user->sample_times[j];
1294: MatMult(user->Qblock,user->yiwork[i],user->di[j]);
1295: }
1296: Gather_i(user->d,user->di,user->di_scatter,user->ns);
1298: /* Now that initial conditions have been set, let the user pass tolerance options to the KSP solver */
1299: KSPSetFromOptions(user->solver);
1300: PetscFree(x);
1301: PetscFree(y);
1302: PetscFree(z);
1303: return(0);
1304: }
1308: PetscErrorCode ParabolicDestroy(AppCtx *user)
1309: {
1311: PetscInt i;
1314: MatDestroy(&user->Qblock);
1315: MatDestroy(&user->QblockT);
1316: MatDestroy(&user->Div);
1317: MatDestroy(&user->Divwork);
1318: MatDestroy(&user->Grad);
1319: MatDestroy(&user->Av);
1320: MatDestroy(&user->Avwork);
1321: MatDestroy(&user->AvT);
1322: MatDestroy(&user->DSG);
1323: MatDestroy(&user->L);
1324: MatDestroy(&user->LT);
1325: KSPDestroy(&user->solver);
1326: MatDestroy(&user->Js);
1327: MatDestroy(&user->Jd);
1328: MatDestroy(&user->JsInv);
1329: MatDestroy(&user->JsBlock);
1330: MatDestroy(&user->JsBlockPrec);
1331: VecDestroy(&user->u);
1332: VecDestroy(&user->uwork);
1333: VecDestroy(&user->utrue);
1334: VecDestroy(&user->y);
1335: VecDestroy(&user->ywork);
1336: VecDestroy(&user->ytrue);
1337: VecDestroyVecs(user->nt,&user->yi);
1338: VecDestroyVecs(user->nt,&user->yiwork);
1339: VecDestroyVecs(user->ns,&user->di);
1340: PetscFree(user->yi);
1341: PetscFree(user->yiwork);
1342: PetscFree(user->di);
1343: VecDestroy(&user->c);
1344: VecDestroy(&user->cwork);
1345: VecDestroy(&user->ur);
1346: VecDestroy(&user->q);
1347: VecDestroy(&user->d);
1348: VecDestroy(&user->dwork);
1349: VecDestroy(&user->lwork);
1350: VecDestroy(&user->S);
1351: VecDestroy(&user->Swork);
1352: VecDestroy(&user->Av_u);
1353: VecDestroy(&user->Twork);
1354: VecDestroy(&user->Rwork);
1355: VecDestroy(&user->js_diag);
1356: ISDestroy(&user->s_is);
1357: ISDestroy(&user->d_is);
1358: VecScatterDestroy(&user->state_scatter);
1359: VecScatterDestroy(&user->design_scatter);
1360: for (i=0; i<user->nt; i++){
1361: VecScatterDestroy(&user->yi_scatter[i]);
1362: }
1363: for (i=0; i<user->ns; i++){
1364: VecScatterDestroy(&user->di_scatter[i]);
1365: }
1366: PetscFree(user->yi_scatter);
1367: PetscFree(user->di_scatter);
1368: PetscFree(user->sample_times);
1369: return(0);
1370: }
1374: PetscErrorCode ParabolicMonitor(Tao tao, void *ptr)
1375: {
1377: Vec X;
1378: PetscReal unorm,ynorm;
1379: AppCtx *user = (AppCtx*)ptr;
1382: TaoGetSolutionVector(tao,&X);
1383: Scatter(X,user->ywork,user->state_scatter,user->uwork,user->design_scatter);
1384: VecAXPY(user->ywork,-1.0,user->ytrue);
1385: VecAXPY(user->uwork,-1.0,user->utrue);
1386: VecNorm(user->uwork,NORM_2,&unorm);
1387: VecNorm(user->ywork,NORM_2,&ynorm);
1388: PetscPrintf(MPI_COMM_WORLD, "||u-ut||=%g ||y-yt||=%g\n",(double)unorm,(double)ynorm);
1389: return(0);
1390: }