Actual source code: ex8.c
petsc-3.9.2 2018-05-20
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: static const char *const BoundaryTypes[] = {"NONE","GHOSTED","MIRROR","PERIODIC","DMBoundaryType","DM_BOUNDARY_",0};
25: /*
26: User-defined data structures and routines
27: */
28: typedef struct {
29: PetscScalar ws; /* Synchronous speed */
30: PetscScalar H; /* Inertia constant */
31: PetscScalar D; /* Damping constant */
32: PetscScalar Pmax,Pmax_s; /* Maximum power output of generator */
33: PetscScalar PM_min; /* Mean mechanical power input */
34: PetscScalar lambda; /* correlation time */
35: PetscScalar q; /* noise strength */
36: PetscScalar mux; /* Initial average angle */
37: PetscScalar sigmax; /* Standard deviation of initial angle */
38: PetscScalar muy; /* Average speed */
39: PetscScalar sigmay; /* standard deviation of initial speed */
40: PetscScalar rho; /* Cross-correlation coefficient */
41: PetscScalar xmin; /* left boundary of angle */
42: PetscScalar xmax; /* right boundary of angle */
43: PetscScalar ymin; /* bottom boundary of speed */
44: PetscScalar ymax; /* top boundary of speed */
45: PetscScalar dx; /* x step size */
46: PetscScalar dy; /* y step size */
47: PetscScalar disper_coe; /* Dispersion coefficient */
48: DM da;
49: PetscInt st_width; /* Stencil width */
50: DMBoundaryType bx; /* x boundary type */
51: DMBoundaryType by; /* y boundary type */
52: PetscReal tf,tcl; /* Fault incidence and clearing times */
53: } AppCtx;
55: PetscErrorCode Parameter_settings(AppCtx*);
56: PetscErrorCode ini_bou(Vec,AppCtx*);
57: PetscErrorCode IFunction(TS,PetscReal,Vec,Vec,Vec,void*);
58: PetscErrorCode IJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
59: PetscErrorCode PostStep(TS);
61: int main(int argc, char **argv)
62: {
64: Vec x; /* Solution vector */
65: TS ts; /* Time-stepping context */
66: AppCtx user; /* Application context */
67: PetscViewer viewer;
69: PetscInitialize(&argc,&argv,"petscopt_ex8", help);if (ierr) return ierr;
71: /* Get physics and time parameters */
72: Parameter_settings(&user);
73: /* Create a 2D DA with dof = 1 */
74: 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);
75: DMSetFromOptions(user.da);
76: DMSetUp(user.da);
77: /* Set x and y coordinates */
78: DMDASetUniformCoordinates(user.da,user.xmin,user.xmax,user.ymin,user.ymax,0,0);
79: DMDASetCoordinateName(user.da,0,"X - the angle");
80: DMDASetCoordinateName(user.da,1,"Y - the speed");
82: /* Get global vector x from DM */
83: DMCreateGlobalVector(user.da,&x);
85: ini_bou(x,&user);
86: PetscViewerBinaryOpen(PETSC_COMM_WORLD,"ini_x",FILE_MODE_WRITE,&viewer);
87: VecView(x,viewer);
88: PetscViewerDestroy(&viewer);
90: TSCreate(PETSC_COMM_WORLD,&ts);
91: TSSetDM(ts,user.da);
92: TSSetProblemType(ts,TS_NONLINEAR);
93: TSSetType(ts,TSARKIMEX);
94: TSSetIFunction(ts,NULL,IFunction,&user);
95: /* TSSetIJacobian(ts,NULL,NULL,IJacobian,&user); */
96: TSSetApplicationContext(ts,&user);
97: TSSetTimeStep(ts,.005);
98: TSSetFromOptions(ts);
99: TSSetPostStep(ts,PostStep);
100: TSSolve(ts,x);
102: PetscViewerBinaryOpen(PETSC_COMM_WORLD,"fin_x",FILE_MODE_WRITE,&viewer);
103: VecView(x,viewer);
104: PetscViewerDestroy(&viewer);
106: VecDestroy(&x);
107: DMDestroy(&user.da);
108: TSDestroy(&ts);
109: PetscFinalize();
110: return ierr;
111: }
113: PetscErrorCode PostStep(TS ts)
114: {
116: Vec X;
117: AppCtx *user;
118: PetscReal t;
119: PetscScalar asum;
122: TSGetApplicationContext(ts,&user);
123: TSGetTime(ts,&t);
124: TSGetSolution(ts,&X);
125: /*
126: if (t >= .2) {
127: TSGetSolution(ts,&X);
128: VecView(X,PETSC_VIEWER_BINARY_WORLD);
129: exit(0);
130: results in initial conditions after fault in binaryoutput
131: }*/
133: 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 */
134: else user->Pmax = user->Pmax_s;
136: VecSum(X,&asum);
137: PetscPrintf(PETSC_COMM_WORLD,"sum(p) at t = %f = %f\n",(double)t,(double)(asum));
138: return(0);
139: }
141: PetscErrorCode ini_bou(Vec X,AppCtx* user)
142: {
144: DM cda;
145: DMDACoor2d **coors;
146: PetscScalar **p;
147: Vec gc;
148: PetscInt M,N,Ir,J;
149: PetscMPIInt rank;
152: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
153: DMDAGetInfo(user->da,NULL,&M,&N,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);
154: user->dx = (user->xmax - user->xmin)/(M-1); user->dy = (user->ymax - user->ymin)/(N-1);
155: DMGetCoordinateDM(user->da,&cda);
156: DMGetCoordinates(user->da,&gc);
157: DMDAVecGetArrayRead(cda,gc,&coors);
158: DMDAVecGetArray(user->da,X,&p);
160: /* Point mass at (mux,muy) */
161: PetscPrintf(PETSC_COMM_WORLD,"Original user->mux = %f, user->muy = %f\n",user->mux,user->muy);
162: DMDAGetLogicalCoordinate(user->da,user->mux,user->muy,0.0,&Ir,&J,NULL,&user->mux,&user->muy,NULL);
163: user->PM_min = user->Pmax*PetscSinScalar(user->mux);
164: 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);
165: if (Ir > -1 && J > -1) {
166: p[J][Ir] = 1.0;
167: }
169: DMDAVecRestoreArrayRead(cda,gc,&coors);
170: DMDAVecRestoreArray(user->da,X,&p);
171: return(0);
172: }
174: /* First advection term */
175: PetscErrorCode adv1(PetscScalar **p,PetscScalar y,PetscInt i,PetscInt j,PetscInt M,PetscScalar *p1,AppCtx *user)
176: {
177: PetscScalar f,fpos,fneg;
179: f = (y - user->ws);
180: fpos = PetscMax(f,0);
181: fneg = PetscMin(f,0);
182: if (user->st_width == 1) {
183: *p1 = fpos*(p[j][i] - p[j][i-1])/user->dx + fneg*(p[j][i+1] - p[j][i])/user->dx;
184: } else if (user->st_width == 2) {
185: *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);
186: } else if (user->st_width == 3) {
187: *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);
188: }
189: /* *p1 = f*(p[j][i+1] - p[j][i-1])/user->dx;*/
190: return(0);
191: }
193: /* Second advection term */
194: PetscErrorCode adv2(PetscScalar **p,PetscScalar x,PetscScalar y,PetscInt i,PetscInt j,PetscInt N,PetscScalar *p2,AppCtx *user)
195: {
196: PetscScalar f,fpos,fneg;
198: f = (user->ws/(2*user->H))*(user->PM_min - user->Pmax*PetscSinScalar(x) - user->D*(y - user->ws));
199: fpos = PetscMax(f,0);
200: fneg = PetscMin(f,0);
201: if (user->st_width == 1) {
202: *p2 = fpos*(p[j][i] - p[j-1][i])/user->dy + fneg*(p[j+1][i] - p[j][i])/user->dy;
203: } else if (user->st_width ==2) {
204: *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);
205: } else if (user->st_width == 3) {
206: *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);
207: }
209: /* *p2 = f*(p[j+1][i] - p[j-1][i])/user->dy;*/
210: return(0);
211: }
213: /* Diffusion term */
214: PetscErrorCode diffuse(PetscScalar **p,PetscInt i,PetscInt j,PetscReal t,PetscScalar *p_diff,AppCtx * user)
215: {
217: if (user->st_width == 1) {
218: *p_diff = user->disper_coe*((p[j-1][i] - 2*p[j][i] + p[j+1][i])/(user->dy*user->dy));
219: } else if (user->st_width == 2) {
220: *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));
221: } else if (user->st_width == 3) {
222: *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));
223: }
224: return(0);
225: }
227: PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
228: {
230: AppCtx *user=(AppCtx*)ctx;
231: DM cda;
232: DMDACoor2d **coors;
233: PetscScalar **p,**f,**pdot;
234: PetscInt i,j;
235: PetscInt xs,ys,xm,ym,M,N;
236: Vec localX,gc,localXdot;
237: PetscScalar p_adv1 = 0.0,p_adv2 = 0.0,p_diff = 0;
238: PetscScalar diffuse1,gamma;
241: DMDAGetInfo(user->da,NULL,&M,&N,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);
242: DMGetCoordinateDM(user->da,&cda);
243: DMDAGetCorners(cda,&xs,&ys,0,&xm,&ym,0);
245: DMGetLocalVector(user->da,&localX);
246: DMGetLocalVector(user->da,&localXdot);
248: DMGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);
249: DMGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);
250: DMGlobalToLocalBegin(user->da,Xdot,INSERT_VALUES,localXdot);
251: DMGlobalToLocalEnd(user->da,Xdot,INSERT_VALUES,localXdot);
253: DMGetCoordinatesLocal(user->da,&gc);
255: DMDAVecGetArrayRead(cda,gc,&coors);
256: DMDAVecGetArrayRead(user->da,localX,&p);
257: DMDAVecGetArrayRead(user->da,localXdot,&pdot);
258: DMDAVecGetArray(user->da,F,&f);
260: gamma = user->D*user->ws/(2*user->H);
261: diffuse1 = user->lambda*user->lambda*user->q/(user->lambda*gamma+1)*(1.0 - PetscExpScalar(-t*(gamma+1.0)/user->lambda));
262: user->disper_coe = user->ws*user->ws/(4*user->H*user->H)*diffuse1;
264: for (i=xs; i < xs+xm; i++) {
265: for (j=ys; j < ys+ym; j++) {
266: adv1(p,coors[j][i].y,i,j,M,&p_adv1,user);
267: adv2(p,coors[j][i].x,coors[j][i].y,i,j,N,&p_adv2,user);
268: diffuse(p,i,j,t,&p_diff,user);
269: f[j][i] = -p_adv1 - p_adv2 + p_diff - pdot[j][i];
270: }
271: }
272: DMDAVecRestoreArrayRead(user->da,localX,&p);
273: DMDAVecRestoreArrayRead(user->da,localX,&pdot);
274: DMRestoreLocalVector(user->da,&localX);
275: DMRestoreLocalVector(user->da,&localXdot);
276: DMDAVecRestoreArray(user->da,F,&f);
277: DMDAVecRestoreArrayRead(cda,gc,&coors);
279: return(0);
280: }
282: PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat J,Mat Jpre,void *ctx)
283: {
285: AppCtx *user=(AppCtx*)ctx;
286: DM cda;
287: DMDACoor2d **coors;
288: PetscInt i,j;
289: PetscInt xs,ys,xm,ym,M,N;
290: Vec gc;
291: PetscScalar val[5],xi,yi;
292: MatStencil row,col[5];
293: PetscScalar c1,c3,c5,c1pos,c1neg,c3pos,c3neg;
296: DMDAGetInfo(user->da,NULL,&M,&N,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);
297: DMGetCoordinateDM(user->da,&cda);
298: DMDAGetCorners(cda,&xs,&ys,0,&xm,&ym,0);
300: DMGetCoordinatesLocal(user->da,&gc);
301: DMDAVecGetArrayRead(cda,gc,&coors);
302: for (i=xs; i < xs+xm; i++) {
303: for (j=ys; j < ys+ym; j++) {
304: PetscInt nc = 0;
305: xi = coors[j][i].x; yi = coors[j][i].y;
306: row.i = i; row.j = j;
307: c1 = (yi-user->ws)/user->dx;
308: c1pos = PetscMax(c1,0);
309: c1neg = PetscMin(c1,0);
310: c3 = (user->ws/(2.0*user->H))*(user->PM_min - user->Pmax*PetscSinScalar(xi) - user->D*(yi - user->ws))/user->dy;
311: c3pos = PetscMax(c3,0);
312: c3neg = PetscMin(c3,0);
313: c5 = (PetscPowScalar((user->lambda*user->ws)/(2*user->H),2)*user->q*(1.0-PetscExpScalar(-t/user->lambda)))/(user->dy*user->dy);
314: col[nc].i = i-1; col[nc].j = j; val[nc++] = c1pos;
315: col[nc].i = i+1; col[nc].j = j; val[nc++] = -c1neg;
316: col[nc].i = i; col[nc].j = j-1; val[nc++] = c3pos + c5;
317: col[nc].i = i; col[nc].j = j+1; val[nc++] = -c3neg + c5;
318: col[nc].i = i; col[nc].j = j; val[nc++] = -c1pos + c1neg -c3pos + c3neg -2*c5 -a;
319: MatSetValuesStencil(Jpre,1,&row,nc,col,val,INSERT_VALUES);
320: }
321: }
322: DMDAVecRestoreArrayRead(cda,gc,&coors);
324: MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY);
325: MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY);
326: if (J != Jpre) {
327: MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
328: MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
329: }
330: return(0);
331: }
333: PetscErrorCode Parameter_settings(AppCtx *user)
334: {
336: PetscBool flg;
340: /* Set default parameters */
341: user->ws = 1.0;
342: user->H = 5.0;
343: user->D = 0.0;
344: user->Pmax = user->Pmax_s = 2.1;
345: user->PM_min = 1.0;
346: user->lambda = 0.1;
347: user->q = 1.0;
348: user->mux = PetscAsinScalar(user->PM_min/user->Pmax);
349: user->sigmax = 0.1;
350: user->sigmay = 0.1;
351: user->rho = 0.0;
352: user->xmin = -PETSC_PI;
353: user->xmax = PETSC_PI;
354: user->bx = DM_BOUNDARY_PERIODIC;
355: user->by = DM_BOUNDARY_GHOSTED;
356: user->tf = user->tcl = -1;
357: user->ymin = -2.0;
358: user->ymax = 2.0;
359: user->st_width = 1;
361: PetscOptionsGetScalar(NULL,NULL,"-ws",&user->ws,&flg);
362: PetscOptionsGetScalar(NULL,NULL,"-Inertia",&user->H,&flg);
363: PetscOptionsGetScalar(NULL,NULL,"-D",&user->D,&flg);
364: PetscOptionsGetScalar(NULL,NULL,"-Pmax",&user->Pmax,&flg);
365: PetscOptionsGetScalar(NULL,NULL,"-PM_min",&user->PM_min,&flg);
366: PetscOptionsGetScalar(NULL,NULL,"-lambda",&user->lambda,&flg);
367: PetscOptionsGetScalar(NULL,NULL,"-q",&user->q,&flg);
368: PetscOptionsGetScalar(NULL,NULL,"-mux",&user->mux,&flg);
369: PetscOptionsGetScalar(NULL,NULL,"-muy",&user->muy,&flg);
370: if (flg == 0) {
371: user->muy = user->ws;
372: }
373: PetscOptionsGetScalar(NULL,NULL,"-xmin",&user->xmin,&flg);
374: PetscOptionsGetScalar(NULL,NULL,"-xmax",&user->xmax,&flg);
375: PetscOptionsGetScalar(NULL,NULL,"-ymin",&user->ymin,&flg);
376: PetscOptionsGetScalar(NULL,NULL,"-ymax",&user->ymax,&flg);
377: PetscOptionsGetInt(NULL,NULL,"-stencil_width",&user->st_width,&flg);
378: PetscOptionsGetEnum(NULL,NULL,"-bx",BoundaryTypes,(PetscEnum*)&user->bx,&flg);
379: PetscOptionsGetEnum(NULL,NULL,"-by",BoundaryTypes,(PetscEnum*)&user->by,&flg);
380: PetscOptionsGetReal(NULL,NULL,"-tf",&user->tf,&flg);
381: PetscOptionsGetReal(NULL,NULL,"-tcl",&user->tcl,&flg);
382: return(0);
383: }
385: /*TEST
387: build:
388: requires: !complex x
390: test:
391: args: -ts_max_steps 1
392: localrunfiles: petscopt_ex8
393: timeoutfactor: 3
395: TEST*/