Actual source code: ex8.c

petsc-3.7.4 2016-10-02
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: 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);

 63: int main(int argc, char **argv)
 64: {
 66:   Vec            x;  /* Solution vector */
 67:   TS             ts;   /* Time-stepping context */
 68:   AppCtx         user; /* Application context */
 69:   PetscViewer    viewer;

 71:   PetscInitialize(&argc,&argv,"petscopt_ex8", help);

 73:   /* Get physics and time parameters */
 74:   Parameter_settings(&user);
 75:   /* Create a 2D DA with dof = 1 */
 76:   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);
 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:   TSSetInitialTimeStep(ts,0.0,.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 0;
111: }

115: PetscErrorCode PostStep(TS ts)
116: {
118:   Vec            X;
119:   AppCtx         *user;
120:   PetscReal      t;
121:   PetscScalar    asum;

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

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

138:   VecSum(X,&asum);
139:   PetscPrintf(PETSC_COMM_WORLD,"sum(p) at t = %f = %f\n",(double)t,(double)(asum));
140:   return(0);
141: }

145: PetscErrorCode ini_bou(Vec X,AppCtx* user)
146: {
148:   DM             cda;
149:   DMDACoor2d     **coors;
150:   PetscScalar    **p;
151:   Vec            gc;
152:   PetscInt       M,N,I,J;
153:   PetscMPIInt    rank;

156:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
157:   DMDAGetInfo(user->da,NULL,&M,&N,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);
158:   user->dx = (user->xmax - user->xmin)/(M-1); user->dy = (user->ymax - user->ymin)/(N-1);
159:   DMGetCoordinateDM(user->da,&cda);
160:   DMGetCoordinates(user->da,&gc);
161:   DMDAVecGetArrayRead(cda,gc,&coors);
162:   DMDAVecGetArray(user->da,X,&p);

164:   /* Point mass at (mux,muy) */
165:   PetscPrintf(PETSC_COMM_WORLD,"Original user->mux = %f, user->muy = %f\n",user->mux,user->muy);
166:   DMDAGetLogicalCoordinate(user->da,user->mux,user->muy,0.0,&I,&J,NULL,&user->mux,&user->muy,NULL);
167:   user->PM_min = user->Pmax*PetscSinScalar(user->mux);
168:   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);
169:   if (I > -1 && J > -1) {
170:     p[J][I] = 1.0;
171:   }

173:   DMDAVecRestoreArrayRead(cda,gc,&coors);
174:   DMDAVecRestoreArray(user->da,X,&p);
175:   return(0);
176: }

178: /* First advection term */
181: PetscErrorCode adv1(PetscScalar **p,PetscScalar y,PetscInt i,PetscInt j,PetscInt M,PetscScalar *p1,AppCtx *user)
182: {
183:   PetscScalar f,fpos,fneg;
185:   f   =  (y - user->ws);
186:   fpos = PetscMax(f,0);
187:   fneg = PetscMin(f,0);
188:   if (user->st_width == 1) {
189:     *p1 = fpos*(p[j][i] - p[j][i-1])/user->dx + fneg*(p[j][i+1] - p[j][i])/user->dx;
190:   } else if (user->st_width == 2) {
191:     *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);
192:   } else if (user->st_width == 3) {
193:     *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);
194:   }
195:   /* *p1 = f*(p[j][i+1] - p[j][i-1])/user->dx;*/
196:   return(0);
197: }

199: /* Second advection term */
202: PetscErrorCode adv2(PetscScalar **p,PetscScalar x,PetscScalar y,PetscInt i,PetscInt j,PetscInt N,PetscScalar *p2,AppCtx *user)
203: {
204:   PetscScalar f,fpos,fneg;
206:   f   = (user->ws/(2*user->H))*(user->PM_min - user->Pmax*PetscSinScalar(x) - user->D*(y - user->ws));
207:   fpos = PetscMax(f,0);
208:   fneg = PetscMin(f,0);
209:   if (user->st_width == 1) {
210:     *p2 = fpos*(p[j][i] - p[j-1][i])/user->dy + fneg*(p[j+1][i] - p[j][i])/user->dy;
211:   } else if (user->st_width ==2) {
212:     *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);
213:   } else if (user->st_width == 3) {
214:     *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);
215:   }

217:   /* *p2 = f*(p[j+1][i] - p[j-1][i])/user->dy;*/
218:   return(0);
219: }

221: /* Diffusion term */
224: PetscErrorCode diffuse(PetscScalar **p,PetscInt i,PetscInt j,PetscReal t,PetscScalar *p_diff,AppCtx * user)
225: {
227:   if (user->st_width == 1) {
228:     *p_diff = user->disper_coe*((p[j-1][i] - 2*p[j][i] + p[j+1][i])/(user->dy*user->dy));
229:   } else if (user->st_width == 2) {
230:     *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));
231:   } else if (user->st_width == 3) {
232:     *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));
233:   }
234:   return(0);
235: }

239: PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
240: {
242:   AppCtx         *user=(AppCtx*)ctx;
243:   DM             cda;
244:   DMDACoor2d     **coors;
245:   PetscScalar    **p,**f,**pdot;
246:   PetscInt       i,j;
247:   PetscInt       xs,ys,xm,ym,M,N;
248:   Vec            localX,gc,localXdot;
249:   PetscScalar    p_adv1,p_adv2,p_diff;

252:   DMDAGetInfo(user->da,NULL,&M,&N,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);
253:   DMGetCoordinateDM(user->da,&cda);
254:   DMDAGetCorners(cda,&xs,&ys,0,&xm,&ym,0);

256:   DMGetLocalVector(user->da,&localX);
257:   DMGetLocalVector(user->da,&localXdot);

259:   DMGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);
260:   DMGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);
261:   DMGlobalToLocalBegin(user->da,Xdot,INSERT_VALUES,localXdot);
262:   DMGlobalToLocalEnd(user->da,Xdot,INSERT_VALUES,localXdot);

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

266:   DMDAVecGetArrayRead(cda,gc,&coors);
267:   DMDAVecGetArrayRead(user->da,localX,&p);
268:   DMDAVecGetArrayRead(user->da,localXdot,&pdot);
269:   DMDAVecGetArray(user->da,F,&f);

271:   PetscScalar diffuse1,gamma;
272:   gamma = user->D*user->ws/(2*user->H);
273:   diffuse1 = user->lambda*user->lambda*user->q/(user->lambda*gamma+1)*(1.0 - PetscExpScalar(-t*(gamma+1.0)/user->lambda));
274:   user->disper_coe = user->ws*user->ws/(4*user->H*user->H)*diffuse1;

276:   for (i=xs; i < xs+xm; i++) {
277:     for (j=ys; j < ys+ym; j++) {
278:       adv1(p,coors[j][i].y,i,j,M,&p_adv1,user);
279:       adv2(p,coors[j][i].x,coors[j][i].y,i,j,N,&p_adv2,user);
280:       diffuse(p,i,j,t,&p_diff,user);
281:       f[j][i] = -p_adv1 - p_adv2  + p_diff - pdot[j][i];
282:     }
283:   }
284:   DMDAVecRestoreArrayRead(user->da,localX,&p);
285:   DMDAVecRestoreArrayRead(user->da,localX,&pdot);
286:   DMRestoreLocalVector(user->da,&localX);
287:   DMRestoreLocalVector(user->da,&localXdot);
288:   DMDAVecRestoreArray(user->da,F,&f);
289:   DMDAVecRestoreArrayRead(cda,gc,&coors);

291:   return(0);
292: }

296: PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat J,Mat Jpre,void *ctx)
297: {
299:   AppCtx         *user=(AppCtx*)ctx;
300:   DM             cda;
301:   DMDACoor2d     **coors;
302:   PetscInt       i,j;
303:   PetscInt       xs,ys,xm,ym,M,N;
304:   Vec            gc;
305:   PetscScalar    val[5],xi,yi;
306:   MatStencil     row,col[5];
307:   PetscScalar    c1,c3,c5,c1pos,c1neg,c3pos,c3neg;

310:   DMDAGetInfo(user->da,NULL,&M,&N,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);
311:   DMGetCoordinateDM(user->da,&cda);
312:   DMDAGetCorners(cda,&xs,&ys,0,&xm,&ym,0);

314:   DMGetCoordinatesLocal(user->da,&gc);
315:   DMDAVecGetArrayRead(cda,gc,&coors);
316:   for (i=xs; i < xs+xm; i++) {
317:     for (j=ys; j < ys+ym; j++) {
318:       PetscInt nc = 0;
319:       xi = coors[j][i].x; yi = coors[j][i].y;
320:       row.i = i; row.j = j;
321:       c1        = (yi-user->ws)/user->dx;
322:       c1pos    = PetscMax(c1,0);
323:       c1neg    = PetscMin(c1,0);
324:       c3        = (user->ws/(2.0*user->H))*(user->PM_min - user->Pmax*PetscSinScalar(xi) - user->D*(yi - user->ws))/user->dy;
325:       c3pos    = PetscMax(c3,0);
326:       c3neg    = PetscMin(c3,0);
327:       c5        = (PetscPowScalar((user->lambda*user->ws)/(2*user->H),2)*user->q*(1.0-PetscExpScalar(-t/user->lambda)))/(user->dy*user->dy);
328:       col[nc].i = i-1; col[nc].j = j;   val[nc++] = c1pos;
329:       col[nc].i = i+1; col[nc].j = j;   val[nc++] = -c1neg;
330:       col[nc].i = i;   col[nc].j = j-1; val[nc++] = c3pos + c5;
331:       col[nc].i = i;   col[nc].j = j+1; val[nc++] = -c3neg + c5;
332:       col[nc].i = i;   col[nc].j = j;   val[nc++] = -c1pos + c1neg -c3pos + c3neg -2*c5 -a;
333:       MatSetValuesStencil(Jpre,1,&row,nc,col,val,INSERT_VALUES);
334:     }
335:   }
336:   DMDAVecRestoreArrayRead(cda,gc,&coors);

338:    MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY);
339:   MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY);
340:   if (J != Jpre) {
341:     MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
342:     MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
343:   }
344:   return(0);
345: }

349: PetscErrorCode Parameter_settings(AppCtx *user)
350: {
352:   PetscBool      flg;


356:   /* Set default parameters */
357:   user->ws     = 1.0;
358:   user->H      = 5.0;
359:   user->D      = 0.0;
360:   user->Pmax = user->Pmax_s  = 2.1;
361:   user->PM_min = 1.0;
362:   user->lambda = 0.1;
363:   user->q      = 1.0;
364:   user->mux    = PetscAsinScalar(user->PM_min/user->Pmax);
365:   user->sigmax = 0.1;
366:   user->sigmay = 0.1;
367:   user->rho    = 0.0;
368:   user->xmin   = -PETSC_PI;
369:   user->xmax   = PETSC_PI;
370:   user->bx     = DM_BOUNDARY_PERIODIC;
371:   user->by     = DM_BOUNDARY_GHOSTED;
372:   user->tf = user->tcl = -1;
373:   user->ymin   = -2.0;
374:   user->ymax   = 2.0;
375:   user->st_width = 1;

377:   PetscOptionsGetScalar(NULL,NULL,"-ws",&user->ws,&flg);
378:   PetscOptionsGetScalar(NULL,NULL,"-Inertia",&user->H,&flg);
379:   PetscOptionsGetScalar(NULL,NULL,"-D",&user->D,&flg);
380:   PetscOptionsGetScalar(NULL,NULL,"-Pmax",&user->Pmax,&flg);
381:   PetscOptionsGetScalar(NULL,NULL,"-PM_min",&user->PM_min,&flg);
382:   PetscOptionsGetScalar(NULL,NULL,"-lambda",&user->lambda,&flg);
383:   PetscOptionsGetScalar(NULL,NULL,"-q",&user->q,&flg);
384:   PetscOptionsGetScalar(NULL,NULL,"-mux",&user->mux,&flg);
385:   PetscOptionsGetScalar(NULL,NULL,"-muy",&user->muy,&flg);
386:   if (flg == 0) {
387:     user->muy = user->ws;
388:   }
389:   PetscOptionsGetScalar(NULL,NULL,"-xmin",&user->xmin,&flg);
390:   PetscOptionsGetScalar(NULL,NULL,"-xmax",&user->xmax,&flg);
391:   PetscOptionsGetScalar(NULL,NULL,"-ymin",&user->ymin,&flg);
392:   PetscOptionsGetScalar(NULL,NULL,"-ymax",&user->ymax,&flg);
393:   PetscOptionsGetInt(NULL,NULL,"-stencil_width",&user->st_width,&flg);
394:   PetscOptionsGetEnum("","-bx",BoundaryTypes,(PetscEnum*)&user->bx,&flg);
395:   PetscOptionsGetEnum("","-by",BoundaryTypes,(PetscEnum*)&user->by,&flg);
396:   PetscOptionsGetReal(NULL,"-tf",&user->tf,&flg);
397:   PetscOptionsGetReal(NULL,"-tcl",&user->tcl,&flg);
398:   return(0);
399: }