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; 5 6 Boundary conditions: 7 Zero dirichlet in y using ghosted values 8 Periodic in x 9 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)) 13 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. 16 17 */ 18 19 #include <petscdm.h> 20 #include <petscdmda.h> 21 #include <petscts.h> 22 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; 52 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); 58 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; 65 66 PetscCall(PetscInitialize(&argc,&argv,"petscopt_ex8", help)); 67 68 /* Get physics and time parameters */ 69 PetscCall(Parameter_settings(&user)); 70 /* Create a 2D DA with dof = 1 */ 71 PetscCall(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 PetscCall(DMSetFromOptions(user.da)); 73 PetscCall(DMSetUp(user.da)); 74 /* Set x and y coordinates */ 75 PetscCall(DMDASetUniformCoordinates(user.da,user.xmin,user.xmax,user.ymin,user.ymax,0,0)); 76 PetscCall(DMDASetCoordinateName(user.da,0,"X - the angle")); 77 PetscCall(DMDASetCoordinateName(user.da,1,"Y - the speed")); 78 79 /* Get global vector x from DM */ 80 PetscCall(DMCreateGlobalVector(user.da,&x)); 81 82 PetscCall(ini_bou(x,&user)); 83 PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"ini_x",FILE_MODE_WRITE,&viewer)); 84 PetscCall(VecView(x,viewer)); 85 PetscCall(PetscViewerDestroy(&viewer)); 86 87 PetscCall(TSCreate(PETSC_COMM_WORLD,&ts)); 88 PetscCall(TSSetDM(ts,user.da)); 89 PetscCall(TSSetProblemType(ts,TS_NONLINEAR)); 90 PetscCall(TSSetType(ts,TSARKIMEX)); 91 PetscCall(TSSetIFunction(ts,NULL,IFunction,&user)); 92 /* PetscCall(TSSetIJacobian(ts,NULL,NULL,IJacobian,&user)); */ 93 PetscCall(TSSetApplicationContext(ts,&user)); 94 PetscCall(TSSetTimeStep(ts,.005)); 95 PetscCall(TSSetFromOptions(ts)); 96 PetscCall(TSSetPostStep(ts,PostStep)); 97 PetscCall(TSSolve(ts,x)); 98 99 PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"fin_x",FILE_MODE_WRITE,&viewer)); 100 PetscCall(VecView(x,viewer)); 101 PetscCall(PetscViewerDestroy(&viewer)); 102 103 PetscCall(VecDestroy(&x)); 104 PetscCall(DMDestroy(&user.da)); 105 PetscCall(TSDestroy(&ts)); 106 PetscCall(PetscFinalize()); 107 return 0; 108 } 109 110 PetscErrorCode PostStep(TS ts) 111 { 112 Vec X; 113 AppCtx *user; 114 PetscReal t; 115 PetscScalar asum; 116 117 PetscFunctionBegin; 118 PetscCall(TSGetApplicationContext(ts,&user)); 119 PetscCall(TSGetTime(ts,&t)); 120 PetscCall(TSGetSolution(ts,&X)); 121 /* 122 if (t >= .2) { 123 PetscCall(TSGetSolution(ts,&X)); 124 PetscCall(VecView(X,PETSC_VIEWER_BINARY_WORLD)); 125 exit(0); 126 results in initial conditions after fault in binaryoutput 127 }*/ 128 129 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 */ 130 else user->Pmax = user->Pmax_s; 131 132 PetscCall(VecSum(X,&asum)); 133 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"sum(p) at t = %f = %f\n",(double)t,(double)(asum))); 134 PetscFunctionReturn(0); 135 } 136 137 PetscErrorCode ini_bou(Vec X,AppCtx* user) 138 { 139 DM cda; 140 DMDACoor2d **coors; 141 PetscScalar **p; 142 Vec gc; 143 PetscInt M,N,Ir,J; 144 PetscMPIInt rank; 145 146 PetscFunctionBeginUser; 147 PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 148 PetscCall(DMDAGetInfo(user->da,NULL,&M,&N,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL)); 149 user->dx = (user->xmax - user->xmin)/(M-1); user->dy = (user->ymax - user->ymin)/(N-1); 150 PetscCall(DMGetCoordinateDM(user->da,&cda)); 151 PetscCall(DMGetCoordinates(user->da,&gc)); 152 PetscCall(DMDAVecGetArrayRead(cda,gc,&coors)); 153 PetscCall(DMDAVecGetArray(user->da,X,&p)); 154 155 /* Point mass at (mux,muy) */ 156 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Original user->mux = %f, user->muy = %f\n",(double)PetscRealPart(user->mux),(double)PetscRealPart(user->muy))); 157 PetscCall(DMDAGetLogicalCoordinate(user->da,user->mux,user->muy,0.0,&Ir,&J,NULL,&user->mux,&user->muy,NULL)); 158 user->PM_min = user->Pmax*PetscSinScalar(user->mux); 159 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Corrected user->mux = %f, user->muy = %f user->PM_min = %f,user->dx = %f\n",(double)PetscRealPart(user->mux),(double)PetscRealPart(user->muy),(double)PetscRealPart(user->PM_min),(double)PetscRealPart(user->dx))); 160 if (Ir > -1 && J > -1) { 161 p[J][Ir] = 1.0; 162 } 163 164 PetscCall(DMDAVecRestoreArrayRead(cda,gc,&coors)); 165 PetscCall(DMDAVecRestoreArray(user->da,X,&p)); 166 PetscFunctionReturn(0); 167 } 168 169 /* First advection term */ 170 PetscErrorCode adv1(PetscScalar **p,PetscScalar y,PetscInt i,PetscInt j,PetscInt M,PetscScalar *p1,AppCtx *user) 171 { 172 PetscScalar f,fpos,fneg; 173 PetscFunctionBegin; 174 f = (y - user->ws); 175 fpos = PetscMax(f,0); 176 fneg = PetscMin(f,0); 177 if (user->st_width == 1) { 178 *p1 = fpos*(p[j][i] - p[j][i-1])/user->dx + fneg*(p[j][i+1] - p[j][i])/user->dx; 179 } else if (user->st_width == 2) { 180 *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); 181 } else if (user->st_width == 3) { 182 *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); 183 } 184 /* *p1 = f*(p[j][i+1] - p[j][i-1])/user->dx;*/ 185 PetscFunctionReturn(0); 186 } 187 188 /* Second advection term */ 189 PetscErrorCode adv2(PetscScalar **p,PetscScalar x,PetscScalar y,PetscInt i,PetscInt j,PetscInt N,PetscScalar *p2,AppCtx *user) 190 { 191 PetscScalar f,fpos,fneg; 192 PetscFunctionBegin; 193 f = (user->ws/(2*user->H))*(user->PM_min - user->Pmax*PetscSinScalar(x) - user->D*(y - user->ws)); 194 fpos = PetscMax(f,0); 195 fneg = PetscMin(f,0); 196 if (user->st_width == 1) { 197 *p2 = fpos*(p[j][i] - p[j-1][i])/user->dy + fneg*(p[j+1][i] - p[j][i])/user->dy; 198 } else if (user->st_width ==2) { 199 *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); 200 } else if (user->st_width == 3) { 201 *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); 202 } 203 204 /* *p2 = f*(p[j+1][i] - p[j-1][i])/user->dy;*/ 205 PetscFunctionReturn(0); 206 } 207 208 /* Diffusion term */ 209 PetscErrorCode diffuse(PetscScalar **p,PetscInt i,PetscInt j,PetscReal t,PetscScalar *p_diff,AppCtx * user) 210 { 211 PetscFunctionBeginUser; 212 if (user->st_width == 1) { 213 *p_diff = user->disper_coe*((p[j-1][i] - 2*p[j][i] + p[j+1][i])/(user->dy*user->dy)); 214 } else if (user->st_width == 2) { 215 *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)); 216 } else if (user->st_width == 3) { 217 *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)); 218 } 219 PetscFunctionReturn(0); 220 } 221 222 PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx) 223 { 224 AppCtx *user = (AppCtx*)ctx; 225 DM cda; 226 DMDACoor2d **coors; 227 PetscScalar **p,**f,**pdot; 228 PetscInt i,j; 229 PetscInt xs,ys,xm,ym,M,N; 230 Vec localX,gc,localXdot; 231 PetscScalar p_adv1 = 0.0,p_adv2 = 0.0,p_diff = 0; 232 PetscScalar diffuse1,gamma; 233 234 PetscFunctionBeginUser; 235 PetscCall(DMDAGetInfo(user->da,NULL,&M,&N,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL)); 236 PetscCall(DMGetCoordinateDM(user->da,&cda)); 237 PetscCall(DMDAGetCorners(cda,&xs,&ys,0,&xm,&ym,0)); 238 239 PetscCall(DMGetLocalVector(user->da,&localX)); 240 PetscCall(DMGetLocalVector(user->da,&localXdot)); 241 242 PetscCall(DMGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX)); 243 PetscCall(DMGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX)); 244 PetscCall(DMGlobalToLocalBegin(user->da,Xdot,INSERT_VALUES,localXdot)); 245 PetscCall(DMGlobalToLocalEnd(user->da,Xdot,INSERT_VALUES,localXdot)); 246 247 PetscCall(DMGetCoordinatesLocal(user->da,&gc)); 248 249 PetscCall(DMDAVecGetArrayRead(cda,gc,&coors)); 250 PetscCall(DMDAVecGetArrayRead(user->da,localX,&p)); 251 PetscCall(DMDAVecGetArrayRead(user->da,localXdot,&pdot)); 252 PetscCall(DMDAVecGetArray(user->da,F,&f)); 253 254 gamma = user->D*user->ws/(2*user->H); 255 diffuse1 = user->lambda*user->lambda*user->q/(user->lambda*gamma+1)*(1.0 - PetscExpScalar(-t*(gamma+1.0)/user->lambda)); 256 user->disper_coe = user->ws*user->ws/(4*user->H*user->H)*diffuse1; 257 258 for (i=xs; i < xs+xm; i++) { 259 for (j=ys; j < ys+ym; j++) { 260 PetscCall(adv1(p,coors[j][i].y,i,j,M,&p_adv1,user)); 261 PetscCall(adv2(p,coors[j][i].x,coors[j][i].y,i,j,N,&p_adv2,user)); 262 PetscCall(diffuse(p,i,j,t,&p_diff,user)); 263 f[j][i] = -p_adv1 - p_adv2 + p_diff - pdot[j][i]; 264 } 265 } 266 PetscCall(DMDAVecRestoreArrayRead(user->da,localX,&p)); 267 PetscCall(DMDAVecRestoreArrayRead(user->da,localX,&pdot)); 268 PetscCall(DMRestoreLocalVector(user->da,&localX)); 269 PetscCall(DMRestoreLocalVector(user->da,&localXdot)); 270 PetscCall(DMDAVecRestoreArray(user->da,F,&f)); 271 PetscCall(DMDAVecRestoreArrayRead(cda,gc,&coors)); 272 273 PetscFunctionReturn(0); 274 } 275 276 PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat J,Mat Jpre,void *ctx) 277 { 278 AppCtx *user=(AppCtx*)ctx; 279 DM cda; 280 DMDACoor2d **coors; 281 PetscInt i,j; 282 PetscInt xs,ys,xm,ym,M,N; 283 Vec gc; 284 PetscScalar val[5],xi,yi; 285 MatStencil row,col[5]; 286 PetscScalar c1,c3,c5,c1pos,c1neg,c3pos,c3neg; 287 288 PetscFunctionBeginUser; 289 PetscCall(DMDAGetInfo(user->da,NULL,&M,&N,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL)); 290 PetscCall(DMGetCoordinateDM(user->da,&cda)); 291 PetscCall(DMDAGetCorners(cda,&xs,&ys,0,&xm,&ym,0)); 292 293 PetscCall(DMGetCoordinatesLocal(user->da,&gc)); 294 PetscCall(DMDAVecGetArrayRead(cda,gc,&coors)); 295 for (i=xs; i < xs+xm; i++) { 296 for (j=ys; j < ys+ym; j++) { 297 PetscInt nc = 0; 298 xi = coors[j][i].x; yi = coors[j][i].y; 299 row.i = i; row.j = j; 300 c1 = (yi-user->ws)/user->dx; 301 c1pos = PetscMax(c1,0); 302 c1neg = PetscMin(c1,0); 303 c3 = (user->ws/(2.0*user->H))*(user->PM_min - user->Pmax*PetscSinScalar(xi) - user->D*(yi - user->ws))/user->dy; 304 c3pos = PetscMax(c3,0); 305 c3neg = PetscMin(c3,0); 306 c5 = (PetscPowScalar((user->lambda*user->ws)/(2*user->H),2)*user->q*(1.0-PetscExpScalar(-t/user->lambda)))/(user->dy*user->dy); 307 col[nc].i = i-1; col[nc].j = j; val[nc++] = c1pos; 308 col[nc].i = i+1; col[nc].j = j; val[nc++] = -c1neg; 309 col[nc].i = i; col[nc].j = j-1; val[nc++] = c3pos + c5; 310 col[nc].i = i; col[nc].j = j+1; val[nc++] = -c3neg + c5; 311 col[nc].i = i; col[nc].j = j; val[nc++] = -c1pos + c1neg -c3pos + c3neg -2*c5 -a; 312 PetscCall(MatSetValuesStencil(Jpre,1,&row,nc,col,val,INSERT_VALUES)); 313 } 314 } 315 PetscCall(DMDAVecRestoreArrayRead(cda,gc,&coors)); 316 317 PetscCall(MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY)); 318 PetscCall(MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY)); 319 if (J != Jpre) { 320 PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY)); 321 PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY)); 322 } 323 PetscFunctionReturn(0); 324 } 325 326 PetscErrorCode Parameter_settings(AppCtx *user) 327 { 328 PetscBool flg; 329 330 PetscFunctionBeginUser; 331 /* Set default parameters */ 332 user->ws = 1.0; 333 user->H = 5.0; 334 user->D = 0.0; 335 user->Pmax = user->Pmax_s = 2.1; 336 user->PM_min = 1.0; 337 user->lambda = 0.1; 338 user->q = 1.0; 339 user->mux = PetscAsinScalar(user->PM_min/user->Pmax); 340 user->sigmax = 0.1; 341 user->sigmay = 0.1; 342 user->rho = 0.0; 343 user->xmin = -PETSC_PI; 344 user->xmax = PETSC_PI; 345 user->bx = DM_BOUNDARY_PERIODIC; 346 user->by = DM_BOUNDARY_GHOSTED; 347 user->tf = user->tcl = -1; 348 user->ymin = -2.0; 349 user->ymax = 2.0; 350 user->st_width = 1; 351 352 PetscCall(PetscOptionsGetScalar(NULL,NULL,"-ws",&user->ws,&flg)); 353 PetscCall(PetscOptionsGetScalar(NULL,NULL,"-Inertia",&user->H,&flg)); 354 PetscCall(PetscOptionsGetScalar(NULL,NULL,"-D",&user->D,&flg)); 355 PetscCall(PetscOptionsGetScalar(NULL,NULL,"-Pmax",&user->Pmax,&flg)); 356 PetscCall(PetscOptionsGetScalar(NULL,NULL,"-PM_min",&user->PM_min,&flg)); 357 PetscCall(PetscOptionsGetScalar(NULL,NULL,"-lambda",&user->lambda,&flg)); 358 PetscCall(PetscOptionsGetScalar(NULL,NULL,"-q",&user->q,&flg)); 359 PetscCall(PetscOptionsGetScalar(NULL,NULL,"-mux",&user->mux,&flg)); 360 PetscCall(PetscOptionsGetScalar(NULL,NULL,"-muy",&user->muy,&flg)); 361 if (flg == 0) user->muy = user->ws; 362 PetscCall(PetscOptionsGetScalar(NULL,NULL,"-xmin",&user->xmin,&flg)); 363 PetscCall(PetscOptionsGetScalar(NULL,NULL,"-xmax",&user->xmax,&flg)); 364 PetscCall(PetscOptionsGetScalar(NULL,NULL,"-ymin",&user->ymin,&flg)); 365 PetscCall(PetscOptionsGetScalar(NULL,NULL,"-ymax",&user->ymax,&flg)); 366 PetscCall(PetscOptionsGetInt(NULL,NULL,"-stencil_width",&user->st_width,&flg)); 367 PetscCall(PetscOptionsGetEnum(NULL,NULL,"-bx",DMBoundaryTypes,(PetscEnum*)&user->bx,&flg)); 368 PetscCall(PetscOptionsGetEnum(NULL,NULL,"-by",DMBoundaryTypes,(PetscEnum*)&user->by,&flg)); 369 PetscCall(PetscOptionsGetReal(NULL,NULL,"-tf",&user->tf,&flg)); 370 PetscCall(PetscOptionsGetReal(NULL,NULL,"-tcl",&user->tcl,&flg)); 371 PetscFunctionReturn(0); 372 } 373 374 /*TEST 375 376 build: 377 requires: !complex x 378 379 test: 380 args: -ts_max_steps 1 381 localrunfiles: petscopt_ex8 382 timeoutfactor: 3 383 384 TEST*/ 385