Actual source code: ex8.c
1: static char help[] = "Time-dependent PDE in 2d for calculating joint PDF. \n";
2: /*
3: p_t = -x_t*p_x -y_t*p_y + f(t)*p_yy
4: xmin < x < xmax, ymin < y < ymax;
6: Boundary conditions:
7: Zero dirichlet in y using ghosted values
8: Periodic in x
10: Note that x_t and y_t in the above are given functions of x and y; they are not derivatives of x and y.
11: x_t = (y - ws)
12: y_t = (ws/2H)*(Pm - Pmax*sin(x) - D*(w - ws))
14: In this example, we can see the effect of a fault, that zeroes the electrical power output
15: Pmax*sin(x), on the PDF. The fault on/off times can be controlled by options -tf and -tcl respectively.
17: */
19: #include <petscdm.h>
20: #include <petscdmda.h>
21: #include <petscts.h>
23: /*
24: User-defined data structures and routines
25: */
26: typedef struct {
27: PetscScalar ws; /* Synchronous speed */
28: PetscScalar H; /* Inertia constant */
29: PetscScalar D; /* Damping constant */
30: PetscScalar Pmax,Pmax_s; /* Maximum power output of generator */
31: PetscScalar PM_min; /* Mean mechanical power input */
32: PetscScalar lambda; /* correlation time */
33: PetscScalar q; /* noise strength */
34: PetscScalar mux; /* Initial average angle */
35: PetscScalar sigmax; /* Standard deviation of initial angle */
36: PetscScalar muy; /* Average speed */
37: PetscScalar sigmay; /* standard deviation of initial speed */
38: PetscScalar rho; /* Cross-correlation coefficient */
39: PetscScalar xmin; /* left boundary of angle */
40: PetscScalar xmax; /* right boundary of angle */
41: PetscScalar ymin; /* bottom boundary of speed */
42: PetscScalar ymax; /* top boundary of speed */
43: PetscScalar dx; /* x step size */
44: PetscScalar dy; /* y step size */
45: PetscScalar disper_coe; /* Dispersion coefficient */
46: DM da;
47: PetscInt st_width; /* Stencil width */
48: DMBoundaryType bx; /* x boundary type */
49: DMBoundaryType by; /* y boundary type */
50: PetscReal tf,tcl; /* Fault incidence and clearing times */
51: } AppCtx;
53: PetscErrorCode Parameter_settings(AppCtx*);
54: PetscErrorCode ini_bou(Vec,AppCtx*);
55: PetscErrorCode IFunction(TS,PetscReal,Vec,Vec,Vec,void*);
56: PetscErrorCode IJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
57: PetscErrorCode PostStep(TS);
59: int main(int argc, char **argv)
60: {
61: Vec x; /* Solution vector */
62: TS ts; /* Time-stepping context */
63: AppCtx user; /* Application context */
64: PetscViewer viewer;
66: PetscInitialize(&argc,&argv,"petscopt_ex8", help);
68: /* Get physics and time parameters */
69: Parameter_settings(&user);
70: /* Create a 2D DA with dof = 1 */
71: DMDACreate2d(PETSC_COMM_WORLD,user.bx,user.by,DMDA_STENCIL_STAR,4,4,PETSC_DECIDE,PETSC_DECIDE,1,user.st_width,NULL,NULL,&user.da);
72: DMSetFromOptions(user.da);
73: DMSetUp(user.da);
74: /* Set x and y coordinates */
75: DMDASetUniformCoordinates(user.da,user.xmin,user.xmax,user.ymin,user.ymax,0,0);
76: DMDASetCoordinateName(user.da,0,"X - the angle");
77: DMDASetCoordinateName(user.da,1,"Y - the speed");
79: /* Get global vector x from DM */
80: DMCreateGlobalVector(user.da,&x);
82: ini_bou(x,&user);
83: PetscViewerBinaryOpen(PETSC_COMM_WORLD,"ini_x",FILE_MODE_WRITE,&viewer);
84: VecView(x,viewer);
85: PetscViewerDestroy(&viewer);
87: TSCreate(PETSC_COMM_WORLD,&ts);
88: TSSetDM(ts,user.da);
89: TSSetProblemType(ts,TS_NONLINEAR);
90: TSSetType(ts,TSARKIMEX);
91: TSSetIFunction(ts,NULL,IFunction,&user);
92: /* TSSetIJacobian(ts,NULL,NULL,IJacobian,&user); */
93: TSSetApplicationContext(ts,&user);
94: TSSetTimeStep(ts,.005);
95: TSSetFromOptions(ts);
96: TSSetPostStep(ts,PostStep);
97: TSSolve(ts,x);
99: PetscViewerBinaryOpen(PETSC_COMM_WORLD,"fin_x",FILE_MODE_WRITE,&viewer);
100: VecView(x,viewer);
101: PetscViewerDestroy(&viewer);
103: VecDestroy(&x);
104: DMDestroy(&user.da);
105: TSDestroy(&ts);
106: PetscFinalize();
107: return 0;
108: }
110: PetscErrorCode PostStep(TS ts)
111: {
112: Vec X;
113: AppCtx *user;
114: PetscReal t;
115: PetscScalar asum;
117: TSGetApplicationContext(ts,&user);
118: TSGetTime(ts,&t);
119: TSGetSolution(ts,&X);
120: /*
121: if (t >= .2) {
122: TSGetSolution(ts,&X);
123: VecView(X,PETSC_VIEWER_BINARY_WORLD);
124: exit(0);
125: results in initial conditions after fault in binaryoutput
126: }*/
128: if ((t > user->tf) && (t < user->tcl)) user->Pmax = 0.0; /* A short-circuit that drives the electrical power output (Pmax*sin(delta)) to zero */
129: else user->Pmax = user->Pmax_s;
131: VecSum(X,&asum);
132: PetscPrintf(PETSC_COMM_WORLD,"sum(p) at t = %f = %f\n",(double)t,(double)(asum));
133: return 0;
134: }
136: PetscErrorCode ini_bou(Vec X,AppCtx* user)
137: {
138: DM cda;
139: DMDACoor2d **coors;
140: PetscScalar **p;
141: Vec gc;
142: PetscInt M,N,Ir,J;
143: PetscMPIInt rank;
146: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
147: DMDAGetInfo(user->da,NULL,&M,&N,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);
148: user->dx = (user->xmax - user->xmin)/(M-1); user->dy = (user->ymax - user->ymin)/(N-1);
149: DMGetCoordinateDM(user->da,&cda);
150: DMGetCoordinates(user->da,&gc);
151: DMDAVecGetArrayRead(cda,gc,&coors);
152: DMDAVecGetArray(user->da,X,&p);
154: /* Point mass at (mux,muy) */
155: PetscPrintf(PETSC_COMM_WORLD,"Original user->mux = %f, user->muy = %f\n",user->mux,user->muy);
156: DMDAGetLogicalCoordinate(user->da,user->mux,user->muy,0.0,&Ir,&J,NULL,&user->mux,&user->muy,NULL);
157: user->PM_min = user->Pmax*PetscSinScalar(user->mux);
158: PetscPrintf(PETSC_COMM_WORLD,"Corrected user->mux = %f, user->muy = %f user->PM_min = %f,user->dx = %f\n",user->mux,user->muy,user->PM_min,user->dx);
159: if (Ir > -1 && J > -1) {
160: p[J][Ir] = 1.0;
161: }
163: DMDAVecRestoreArrayRead(cda,gc,&coors);
164: DMDAVecRestoreArray(user->da,X,&p);
165: return 0;
166: }
168: /* First advection term */
169: PetscErrorCode adv1(PetscScalar **p,PetscScalar y,PetscInt i,PetscInt j,PetscInt M,PetscScalar *p1,AppCtx *user)
170: {
171: PetscScalar f,fpos,fneg;
172: f = (y - user->ws);
173: fpos = PetscMax(f,0);
174: fneg = PetscMin(f,0);
175: if (user->st_width == 1) {
176: *p1 = fpos*(p[j][i] - p[j][i-1])/user->dx + fneg*(p[j][i+1] - p[j][i])/user->dx;
177: } else if (user->st_width == 2) {
178: *p1 = fpos*(3*p[j][i] - 4*p[j][i-1] + p[j][i-2])/(2*user->dx) + fneg*(-p[j][i+2] + 4*p[j][i+1] - 3*p[j][i])/(2*user->dx);
179: } else if (user->st_width == 3) {
180: *p1 = fpos*(2*p[j][i+1] + 3*p[j][i] - 6*p[j][i-1] + p[j][i-2])/(6*user->dx) + fneg*(-p[j][i+2] + 6*p[j][i+1] - 3*p[j][i] - 2*p[j][i-1])/(6*user->dx);
181: }
182: /* *p1 = f*(p[j][i+1] - p[j][i-1])/user->dx;*/
183: return 0;
184: }
186: /* Second advection term */
187: PetscErrorCode adv2(PetscScalar **p,PetscScalar x,PetscScalar y,PetscInt i,PetscInt j,PetscInt N,PetscScalar *p2,AppCtx *user)
188: {
189: PetscScalar f,fpos,fneg;
190: f = (user->ws/(2*user->H))*(user->PM_min - user->Pmax*PetscSinScalar(x) - user->D*(y - user->ws));
191: fpos = PetscMax(f,0);
192: fneg = PetscMin(f,0);
193: if (user->st_width == 1) {
194: *p2 = fpos*(p[j][i] - p[j-1][i])/user->dy + fneg*(p[j+1][i] - p[j][i])/user->dy;
195: } else if (user->st_width ==2) {
196: *p2 = fpos*(3*p[j][i] - 4*p[j-1][i] + p[j-2][i])/(2*user->dy) + fneg*(-p[j+2][i] + 4*p[j+1][i] - 3*p[j][i])/(2*user->dy);
197: } else if (user->st_width == 3) {
198: *p2 = fpos*(2*p[j+1][i] + 3*p[j][i] - 6*p[j-1][i] + p[j-2][i])/(6*user->dy) + fneg*(-p[j+2][i] + 6*p[j+1][i] - 3*p[j][i] - 2*p[j-1][i])/(6*user->dy);
199: }
201: /* *p2 = f*(p[j+1][i] - p[j-1][i])/user->dy;*/
202: return 0;
203: }
205: /* Diffusion term */
206: PetscErrorCode diffuse(PetscScalar **p,PetscInt i,PetscInt j,PetscReal t,PetscScalar *p_diff,AppCtx * user)
207: {
209: if (user->st_width == 1) {
210: *p_diff = user->disper_coe*((p[j-1][i] - 2*p[j][i] + p[j+1][i])/(user->dy*user->dy));
211: } else if (user->st_width == 2) {
212: *p_diff = user->disper_coe*((-p[j-2][i] + 16*p[j-1][i] - 30*p[j][i] + 16*p[j+1][i] - p[j+2][i])/(12.0*user->dy*user->dy));
213: } else if (user->st_width == 3) {
214: *p_diff = user->disper_coe*((2*p[j-3][i] - 27*p[j-2][i] + 270*p[j-1][i] - 490*p[j][i] + 270*p[j+1][i] - 27*p[j+2][i] + 2*p[j+3][i])/(180.0*user->dy*user->dy));
215: }
216: return 0;
217: }
219: PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
220: {
221: AppCtx *user = (AppCtx*)ctx;
222: DM cda;
223: DMDACoor2d **coors;
224: PetscScalar **p,**f,**pdot;
225: PetscInt i,j;
226: PetscInt xs,ys,xm,ym,M,N;
227: Vec localX,gc,localXdot;
228: PetscScalar p_adv1 = 0.0,p_adv2 = 0.0,p_diff = 0;
229: PetscScalar diffuse1,gamma;
232: DMDAGetInfo(user->da,NULL,&M,&N,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);
233: DMGetCoordinateDM(user->da,&cda);
234: DMDAGetCorners(cda,&xs,&ys,0,&xm,&ym,0);
236: DMGetLocalVector(user->da,&localX);
237: DMGetLocalVector(user->da,&localXdot);
239: DMGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);
240: DMGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);
241: DMGlobalToLocalBegin(user->da,Xdot,INSERT_VALUES,localXdot);
242: DMGlobalToLocalEnd(user->da,Xdot,INSERT_VALUES,localXdot);
244: DMGetCoordinatesLocal(user->da,&gc);
246: DMDAVecGetArrayRead(cda,gc,&coors);
247: DMDAVecGetArrayRead(user->da,localX,&p);
248: DMDAVecGetArrayRead(user->da,localXdot,&pdot);
249: DMDAVecGetArray(user->da,F,&f);
251: gamma = user->D*user->ws/(2*user->H);
252: diffuse1 = user->lambda*user->lambda*user->q/(user->lambda*gamma+1)*(1.0 - PetscExpScalar(-t*(gamma+1.0)/user->lambda));
253: user->disper_coe = user->ws*user->ws/(4*user->H*user->H)*diffuse1;
255: for (i=xs; i < xs+xm; i++) {
256: for (j=ys; j < ys+ym; j++) {
257: adv1(p,coors[j][i].y,i,j,M,&p_adv1,user);
258: adv2(p,coors[j][i].x,coors[j][i].y,i,j,N,&p_adv2,user);
259: diffuse(p,i,j,t,&p_diff,user);
260: f[j][i] = -p_adv1 - p_adv2 + p_diff - pdot[j][i];
261: }
262: }
263: DMDAVecRestoreArrayRead(user->da,localX,&p);
264: DMDAVecRestoreArrayRead(user->da,localX,&pdot);
265: DMRestoreLocalVector(user->da,&localX);
266: DMRestoreLocalVector(user->da,&localXdot);
267: DMDAVecRestoreArray(user->da,F,&f);
268: DMDAVecRestoreArrayRead(cda,gc,&coors);
270: return 0;
271: }
273: PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat J,Mat Jpre,void *ctx)
274: {
275: AppCtx *user=(AppCtx*)ctx;
276: DM cda;
277: DMDACoor2d **coors;
278: PetscInt i,j;
279: PetscInt xs,ys,xm,ym,M,N;
280: Vec gc;
281: PetscScalar val[5],xi,yi;
282: MatStencil row,col[5];
283: PetscScalar c1,c3,c5,c1pos,c1neg,c3pos,c3neg;
286: DMDAGetInfo(user->da,NULL,&M,&N,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);
287: DMGetCoordinateDM(user->da,&cda);
288: DMDAGetCorners(cda,&xs,&ys,0,&xm,&ym,0);
290: DMGetCoordinatesLocal(user->da,&gc);
291: DMDAVecGetArrayRead(cda,gc,&coors);
292: for (i=xs; i < xs+xm; i++) {
293: for (j=ys; j < ys+ym; j++) {
294: PetscInt nc = 0;
295: xi = coors[j][i].x; yi = coors[j][i].y;
296: row.i = i; row.j = j;
297: c1 = (yi-user->ws)/user->dx;
298: c1pos = PetscMax(c1,0);
299: c1neg = PetscMin(c1,0);
300: c3 = (user->ws/(2.0*user->H))*(user->PM_min - user->Pmax*PetscSinScalar(xi) - user->D*(yi - user->ws))/user->dy;
301: c3pos = PetscMax(c3,0);
302: c3neg = PetscMin(c3,0);
303: c5 = (PetscPowScalar((user->lambda*user->ws)/(2*user->H),2)*user->q*(1.0-PetscExpScalar(-t/user->lambda)))/(user->dy*user->dy);
304: col[nc].i = i-1; col[nc].j = j; val[nc++] = c1pos;
305: col[nc].i = i+1; col[nc].j = j; val[nc++] = -c1neg;
306: col[nc].i = i; col[nc].j = j-1; val[nc++] = c3pos + c5;
307: col[nc].i = i; col[nc].j = j+1; val[nc++] = -c3neg + c5;
308: col[nc].i = i; col[nc].j = j; val[nc++] = -c1pos + c1neg -c3pos + c3neg -2*c5 -a;
309: MatSetValuesStencil(Jpre,1,&row,nc,col,val,INSERT_VALUES);
310: }
311: }
312: DMDAVecRestoreArrayRead(cda,gc,&coors);
314: MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY);
315: MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY);
316: if (J != Jpre) {
317: MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
318: MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
319: }
320: return 0;
321: }
323: PetscErrorCode Parameter_settings(AppCtx *user)
324: {
325: PetscBool flg;
328: /* Set default parameters */
329: user->ws = 1.0;
330: user->H = 5.0;
331: user->D = 0.0;
332: user->Pmax = user->Pmax_s = 2.1;
333: user->PM_min = 1.0;
334: user->lambda = 0.1;
335: user->q = 1.0;
336: user->mux = PetscAsinScalar(user->PM_min/user->Pmax);
337: user->sigmax = 0.1;
338: user->sigmay = 0.1;
339: user->rho = 0.0;
340: user->xmin = -PETSC_PI;
341: user->xmax = PETSC_PI;
342: user->bx = DM_BOUNDARY_PERIODIC;
343: user->by = DM_BOUNDARY_GHOSTED;
344: user->tf = user->tcl = -1;
345: user->ymin = -2.0;
346: user->ymax = 2.0;
347: user->st_width = 1;
349: PetscOptionsGetScalar(NULL,NULL,"-ws",&user->ws,&flg);
350: PetscOptionsGetScalar(NULL,NULL,"-Inertia",&user->H,&flg);
351: PetscOptionsGetScalar(NULL,NULL,"-D",&user->D,&flg);
352: PetscOptionsGetScalar(NULL,NULL,"-Pmax",&user->Pmax,&flg);
353: PetscOptionsGetScalar(NULL,NULL,"-PM_min",&user->PM_min,&flg);
354: PetscOptionsGetScalar(NULL,NULL,"-lambda",&user->lambda,&flg);
355: PetscOptionsGetScalar(NULL,NULL,"-q",&user->q,&flg);
356: PetscOptionsGetScalar(NULL,NULL,"-mux",&user->mux,&flg);
357: PetscOptionsGetScalar(NULL,NULL,"-muy",&user->muy,&flg);
358: if (flg == 0) user->muy = user->ws;
359: PetscOptionsGetScalar(NULL,NULL,"-xmin",&user->xmin,&flg);
360: PetscOptionsGetScalar(NULL,NULL,"-xmax",&user->xmax,&flg);
361: PetscOptionsGetScalar(NULL,NULL,"-ymin",&user->ymin,&flg);
362: PetscOptionsGetScalar(NULL,NULL,"-ymax",&user->ymax,&flg);
363: PetscOptionsGetInt(NULL,NULL,"-stencil_width",&user->st_width,&flg);
364: PetscOptionsGetEnum(NULL,NULL,"-bx",DMBoundaryTypes,(PetscEnum*)&user->bx,&flg);
365: PetscOptionsGetEnum(NULL,NULL,"-by",DMBoundaryTypes,(PetscEnum*)&user->by,&flg);
366: PetscOptionsGetReal(NULL,NULL,"-tf",&user->tf,&flg);
367: PetscOptionsGetReal(NULL,NULL,"-tcl",&user->tcl,&flg);
368: return 0;
369: }
371: /*TEST
373: build:
374: requires: !complex x
376: test:
377: args: -ts_max_steps 1
378: localrunfiles: petscopt_ex8
379: timeoutfactor: 3
381: TEST*/