Actual source code: ex8.c

petsc-3.12.0 2019-09-29
Report Typos and Errors
  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: {
 62:   Vec            x;  /* Solution vector */
 63:   TS             ts;   /* Time-stepping context */
 64:   AppCtx         user; /* Application context */
 65:   PetscViewer    viewer;

 67:   PetscInitialize(&argc,&argv,"petscopt_ex8", help);if (ierr) return ierr;

 69:   /* Get physics and time parameters */
 70:   Parameter_settings(&user);
 71:   /* Create a 2D DA with dof = 1 */
 72:   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);
 73:   DMSetFromOptions(user.da);
 74:   DMSetUp(user.da);
 75:   /* Set x and y coordinates */
 76:   DMDASetUniformCoordinates(user.da,user.xmin,user.xmax,user.ymin,user.ymax,0,0);
 77:   DMDASetCoordinateName(user.da,0,"X - the angle");
 78:   DMDASetCoordinateName(user.da,1,"Y - the speed");

 80:   /* Get global vector x from DM  */
 81:   DMCreateGlobalVector(user.da,&x);

 83:   ini_bou(x,&user);
 84:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"ini_x",FILE_MODE_WRITE,&viewer);
 85:   VecView(x,viewer);
 86:   PetscViewerDestroy(&viewer);

 88:   TSCreate(PETSC_COMM_WORLD,&ts);
 89:   TSSetDM(ts,user.da);
 90:   TSSetProblemType(ts,TS_NONLINEAR);
 91:   TSSetType(ts,TSARKIMEX);
 92:   TSSetIFunction(ts,NULL,IFunction,&user);
 93:   /*  TSSetIJacobian(ts,NULL,NULL,IJacobian,&user); */
 94:   TSSetApplicationContext(ts,&user);
 95:   TSSetTimeStep(ts,.005);
 96:   TSSetFromOptions(ts);
 97:   TSSetPostStep(ts,PostStep);
 98:   TSSolve(ts,x);

100:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"fin_x",FILE_MODE_WRITE,&viewer);
101:   VecView(x,viewer);
102:   PetscViewerDestroy(&viewer);

104:   VecDestroy(&x);
105:   DMDestroy(&user.da);
106:   TSDestroy(&ts);
107:   PetscFinalize();
108:   return ierr;
109: }

111: PetscErrorCode PostStep(TS ts)
112: {
114:   Vec            X;
115:   AppCtx         *user;
116:   PetscReal      t;
117:   PetscScalar    asum;

120:   TSGetApplicationContext(ts,&user);
121:   TSGetTime(ts,&t);
122:   TSGetSolution(ts,&X);
123:   /*
124:   if (t >= .2) {
125:     TSGetSolution(ts,&X);
126:     VecView(X,PETSC_VIEWER_BINARY_WORLD);
127:     exit(0);
128:      results in initial conditions after fault in binaryoutput
129:   }*/

131:   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 */
132:   else user->Pmax = user->Pmax_s;

134:   VecSum(X,&asum);
135:   PetscPrintf(PETSC_COMM_WORLD,"sum(p) at t = %f = %f\n",(double)t,(double)(asum));
136:   return(0);
137: }

139: PetscErrorCode ini_bou(Vec X,AppCtx* user)
140: {
142:   DM             cda;
143:   DMDACoor2d     **coors;
144:   PetscScalar    **p;
145:   Vec            gc;
146:   PetscInt       M,N,Ir,J;
147:   PetscMPIInt    rank;

150:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
151:   DMDAGetInfo(user->da,NULL,&M,&N,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);
152:   user->dx = (user->xmax - user->xmin)/(M-1); user->dy = (user->ymax - user->ymin)/(N-1);
153:   DMGetCoordinateDM(user->da,&cda);
154:   DMGetCoordinates(user->da,&gc);
155:   DMDAVecGetArrayRead(cda,gc,&coors);
156:   DMDAVecGetArray(user->da,X,&p);

158:   /* Point mass at (mux,muy) */
159:   PetscPrintf(PETSC_COMM_WORLD,"Original user->mux = %f, user->muy = %f\n",user->mux,user->muy);
160:   DMDAGetLogicalCoordinate(user->da,user->mux,user->muy,0.0,&Ir,&J,NULL,&user->mux,&user->muy,NULL);
161:   user->PM_min = user->Pmax*PetscSinScalar(user->mux);
162:   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);
163:   if (Ir > -1 && J > -1) {
164:     p[J][Ir] = 1.0;
165:   }

167:   DMDAVecRestoreArrayRead(cda,gc,&coors);
168:   DMDAVecRestoreArray(user->da,X,&p);
169:   return(0);
170: }

172: /* First advection term */
173: PetscErrorCode adv1(PetscScalar **p,PetscScalar y,PetscInt i,PetscInt j,PetscInt M,PetscScalar *p1,AppCtx *user)
174: {
175:   PetscScalar f,fpos,fneg;
177:   f   =  (y - user->ws);
178:   fpos = PetscMax(f,0);
179:   fneg = PetscMin(f,0);
180:   if (user->st_width == 1) {
181:     *p1 = fpos*(p[j][i] - p[j][i-1])/user->dx + fneg*(p[j][i+1] - p[j][i])/user->dx;
182:   } else if (user->st_width == 2) {
183:     *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);
184:   } else if (user->st_width == 3) {
185:     *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);
186:   }
187:   /* *p1 = f*(p[j][i+1] - p[j][i-1])/user->dx;*/
188:   return(0);
189: }

191: /* Second advection term */
192: PetscErrorCode adv2(PetscScalar **p,PetscScalar x,PetscScalar y,PetscInt i,PetscInt j,PetscInt N,PetscScalar *p2,AppCtx *user)
193: {
194:   PetscScalar f,fpos,fneg;
196:   f   = (user->ws/(2*user->H))*(user->PM_min - user->Pmax*PetscSinScalar(x) - user->D*(y - user->ws));
197:   fpos = PetscMax(f,0);
198:   fneg = PetscMin(f,0);
199:   if (user->st_width == 1) {
200:     *p2 = fpos*(p[j][i] - p[j-1][i])/user->dy + fneg*(p[j+1][i] - p[j][i])/user->dy;
201:   } else if (user->st_width ==2) {
202:     *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);
203:   } else if (user->st_width == 3) {
204:     *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);
205:   }

207:   /* *p2 = f*(p[j+1][i] - p[j-1][i])/user->dy;*/
208:   return(0);
209: }

211: /* Diffusion term */
212: PetscErrorCode diffuse(PetscScalar **p,PetscInt i,PetscInt j,PetscReal t,PetscScalar *p_diff,AppCtx * user)
213: {
215:   if (user->st_width == 1) {
216:     *p_diff = user->disper_coe*((p[j-1][i] - 2*p[j][i] + p[j+1][i])/(user->dy*user->dy));
217:   } else if (user->st_width == 2) {
218:     *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));
219:   } else if (user->st_width == 3) {
220:     *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));
221:   }
222:   return(0);
223: }

225: PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
226: {
228:   AppCtx         *user=(AppCtx*)ctx;
229:   DM             cda;
230:   DMDACoor2d     **coors;
231:   PetscScalar    **p,**f,**pdot;
232:   PetscInt       i,j;
233:   PetscInt       xs,ys,xm,ym,M,N;
234:   Vec            localX,gc,localXdot;
235:   PetscScalar    p_adv1 = 0.0,p_adv2 = 0.0,p_diff = 0;
236:   PetscScalar    diffuse1,gamma;

239:   DMDAGetInfo(user->da,NULL,&M,&N,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);
240:   DMGetCoordinateDM(user->da,&cda);
241:   DMDAGetCorners(cda,&xs,&ys,0,&xm,&ym,0);

243:   DMGetLocalVector(user->da,&localX);
244:   DMGetLocalVector(user->da,&localXdot);

246:   DMGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);
247:   DMGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);
248:   DMGlobalToLocalBegin(user->da,Xdot,INSERT_VALUES,localXdot);
249:   DMGlobalToLocalEnd(user->da,Xdot,INSERT_VALUES,localXdot);

251:   DMGetCoordinatesLocal(user->da,&gc);

253:   DMDAVecGetArrayRead(cda,gc,&coors);
254:   DMDAVecGetArrayRead(user->da,localX,&p);
255:   DMDAVecGetArrayRead(user->da,localXdot,&pdot);
256:   DMDAVecGetArray(user->da,F,&f);

258:   gamma = user->D*user->ws/(2*user->H);
259:   diffuse1 = user->lambda*user->lambda*user->q/(user->lambda*gamma+1)*(1.0 - PetscExpScalar(-t*(gamma+1.0)/user->lambda));
260:   user->disper_coe = user->ws*user->ws/(4*user->H*user->H)*diffuse1;

262:   for (i=xs; i < xs+xm; i++) {
263:     for (j=ys; j < ys+ym; j++) {
264:       adv1(p,coors[j][i].y,i,j,M,&p_adv1,user);
265:       adv2(p,coors[j][i].x,coors[j][i].y,i,j,N,&p_adv2,user);
266:       diffuse(p,i,j,t,&p_diff,user);
267:       f[j][i] = -p_adv1 - p_adv2  + p_diff - pdot[j][i];
268:     }
269:   }
270:   DMDAVecRestoreArrayRead(user->da,localX,&p);
271:   DMDAVecRestoreArrayRead(user->da,localX,&pdot);
272:   DMRestoreLocalVector(user->da,&localX);
273:   DMRestoreLocalVector(user->da,&localXdot);
274:   DMDAVecRestoreArray(user->da,F,&f);
275:   DMDAVecRestoreArrayRead(cda,gc,&coors);

277:   return(0);
278: }

280: PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat J,Mat Jpre,void *ctx)
281: {
283:   AppCtx         *user=(AppCtx*)ctx;
284:   DM             cda;
285:   DMDACoor2d     **coors;
286:   PetscInt       i,j;
287:   PetscInt       xs,ys,xm,ym,M,N;
288:   Vec            gc;
289:   PetscScalar    val[5],xi,yi;
290:   MatStencil     row,col[5];
291:   PetscScalar    c1,c3,c5,c1pos,c1neg,c3pos,c3neg;

294:   DMDAGetInfo(user->da,NULL,&M,&N,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);
295:   DMGetCoordinateDM(user->da,&cda);
296:   DMDAGetCorners(cda,&xs,&ys,0,&xm,&ym,0);

298:   DMGetCoordinatesLocal(user->da,&gc);
299:   DMDAVecGetArrayRead(cda,gc,&coors);
300:   for (i=xs; i < xs+xm; i++) {
301:     for (j=ys; j < ys+ym; j++) {
302:       PetscInt nc = 0;
303:       xi = coors[j][i].x; yi = coors[j][i].y;
304:       row.i = i; row.j = j;
305:       c1        = (yi-user->ws)/user->dx;
306:       c1pos    = PetscMax(c1,0);
307:       c1neg    = PetscMin(c1,0);
308:       c3        = (user->ws/(2.0*user->H))*(user->PM_min - user->Pmax*PetscSinScalar(xi) - user->D*(yi - user->ws))/user->dy;
309:       c3pos    = PetscMax(c3,0);
310:       c3neg    = PetscMin(c3,0);
311:       c5        = (PetscPowScalar((user->lambda*user->ws)/(2*user->H),2)*user->q*(1.0-PetscExpScalar(-t/user->lambda)))/(user->dy*user->dy);
312:       col[nc].i = i-1; col[nc].j = j;   val[nc++] = c1pos;
313:       col[nc].i = i+1; col[nc].j = j;   val[nc++] = -c1neg;
314:       col[nc].i = i;   col[nc].j = j-1; val[nc++] = c3pos + c5;
315:       col[nc].i = i;   col[nc].j = j+1; val[nc++] = -c3neg + c5;
316:       col[nc].i = i;   col[nc].j = j;   val[nc++] = -c1pos + c1neg -c3pos + c3neg -2*c5 -a;
317:       MatSetValuesStencil(Jpre,1,&row,nc,col,val,INSERT_VALUES);
318:     }
319:   }
320:   DMDAVecRestoreArrayRead(cda,gc,&coors);

322:    MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY);
323:   MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY);
324:   if (J != Jpre) {
325:     MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
326:     MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
327:   }
328:   return(0);
329: }

331: PetscErrorCode Parameter_settings(AppCtx *user)
332: {
334:   PetscBool      flg;


338:   /* Set default parameters */
339:   user->ws     = 1.0;
340:   user->H      = 5.0;
341:   user->D      = 0.0;
342:   user->Pmax = user->Pmax_s  = 2.1;
343:   user->PM_min = 1.0;
344:   user->lambda = 0.1;
345:   user->q      = 1.0;
346:   user->mux    = PetscAsinScalar(user->PM_min/user->Pmax);
347:   user->sigmax = 0.1;
348:   user->sigmay = 0.1;
349:   user->rho    = 0.0;
350:   user->xmin   = -PETSC_PI;
351:   user->xmax   = PETSC_PI;
352:   user->bx     = DM_BOUNDARY_PERIODIC;
353:   user->by     = DM_BOUNDARY_GHOSTED;
354:   user->tf = user->tcl = -1;
355:   user->ymin   = -2.0;
356:   user->ymax   = 2.0;
357:   user->st_width = 1;

359:   PetscOptionsGetScalar(NULL,NULL,"-ws",&user->ws,&flg);
360:   PetscOptionsGetScalar(NULL,NULL,"-Inertia",&user->H,&flg);
361:   PetscOptionsGetScalar(NULL,NULL,"-D",&user->D,&flg);
362:   PetscOptionsGetScalar(NULL,NULL,"-Pmax",&user->Pmax,&flg);
363:   PetscOptionsGetScalar(NULL,NULL,"-PM_min",&user->PM_min,&flg);
364:   PetscOptionsGetScalar(NULL,NULL,"-lambda",&user->lambda,&flg);
365:   PetscOptionsGetScalar(NULL,NULL,"-q",&user->q,&flg);
366:   PetscOptionsGetScalar(NULL,NULL,"-mux",&user->mux,&flg);
367:   PetscOptionsGetScalar(NULL,NULL,"-muy",&user->muy,&flg);
368:   if (flg == 0) {
369:     user->muy = user->ws;
370:   }
371:   PetscOptionsGetScalar(NULL,NULL,"-xmin",&user->xmin,&flg);
372:   PetscOptionsGetScalar(NULL,NULL,"-xmax",&user->xmax,&flg);
373:   PetscOptionsGetScalar(NULL,NULL,"-ymin",&user->ymin,&flg);
374:   PetscOptionsGetScalar(NULL,NULL,"-ymax",&user->ymax,&flg);
375:   PetscOptionsGetInt(NULL,NULL,"-stencil_width",&user->st_width,&flg);
376:   PetscOptionsGetEnum(NULL,NULL,"-bx",DMBoundaryTypes,(PetscEnum*)&user->bx,&flg);
377:   PetscOptionsGetEnum(NULL,NULL,"-by",DMBoundaryTypes,(PetscEnum*)&user->by,&flg);
378:   PetscOptionsGetReal(NULL,NULL,"-tf",&user->tf,&flg);
379:   PetscOptionsGetReal(NULL,NULL,"-tcl",&user->tcl,&flg);
380:   return(0);
381: }

383: /*TEST

385:    build:
386:       requires: !complex x

388:    test:
389:       args: -ts_max_steps 1
390:       localrunfiles: petscopt_ex8
391:       timeoutfactor: 3

393: TEST*/