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 Vec x; /* Solution vector */ 61 TS ts; /* Time-stepping context */ 62 AppCtx user; /* Application context */ 63 PetscViewer viewer; 64 65 PetscFunctionBeginUser; 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 Vec X; 112 AppCtx *user; 113 PetscReal t; 114 PetscScalar asum; 115 116 PetscFunctionBegin; 117 PetscCall(TSGetApplicationContext(ts, &user)); 118 PetscCall(TSGetTime(ts, &t)); 119 PetscCall(TSGetSolution(ts, &X)); 120 /* 121 if (t >= .2) { 122 PetscCall(TSGetSolution(ts,&X)); 123 PetscCall(VecView(X,PETSC_VIEWER_BINARY_WORLD)); 124 exit(0); 125 results in initial conditions after fault in binaryoutput 126 }*/ 127 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; 130 131 PetscCall(VecSum(X, &asum)); 132 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "sum(p) at t = %f = %f\n", (double)t, (double)(asum))); 133 PetscFunctionReturn(0); 134 } 135 136 PetscErrorCode ini_bou(Vec X, AppCtx *user) { 137 DM cda; 138 DMDACoor2d **coors; 139 PetscScalar **p; 140 Vec gc; 141 PetscInt M, N, Ir, J; 142 PetscMPIInt rank; 143 144 PetscFunctionBeginUser; 145 PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 146 PetscCall(DMDAGetInfo(user->da, NULL, &M, &N, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 147 user->dx = (user->xmax - user->xmin) / (M - 1); 148 user->dy = (user->ymax - user->ymin) / (N - 1); 149 PetscCall(DMGetCoordinateDM(user->da, &cda)); 150 PetscCall(DMGetCoordinates(user->da, &gc)); 151 PetscCall(DMDAVecGetArrayRead(cda, gc, &coors)); 152 PetscCall(DMDAVecGetArray(user->da, X, &p)); 153 154 /* Point mass at (mux,muy) */ 155 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Original user->mux = %f, user->muy = %f\n", (double)PetscRealPart(user->mux), (double)PetscRealPart(user->muy))); 156 PetscCall(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 PetscCall( 159 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) { p[J][Ir] = 1.0; } 161 162 PetscCall(DMDAVecRestoreArrayRead(cda, gc, &coors)); 163 PetscCall(DMDAVecRestoreArray(user->da, X, &p)); 164 PetscFunctionReturn(0); 165 } 166 167 /* First advection term */ 168 PetscErrorCode adv1(PetscScalar **p, PetscScalar y, PetscInt i, PetscInt j, PetscInt M, PetscScalar *p1, AppCtx *user) { 169 PetscScalar f, fpos, fneg; 170 PetscFunctionBegin; 171 f = (y - user->ws); 172 fpos = PetscMax(f, 0); 173 fneg = PetscMin(f, 0); 174 if (user->st_width == 1) { 175 *p1 = fpos * (p[j][i] - p[j][i - 1]) / user->dx + fneg * (p[j][i + 1] - p[j][i]) / user->dx; 176 } else if (user->st_width == 2) { 177 *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); 178 } else if (user->st_width == 3) { 179 *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); 180 } 181 /* *p1 = f*(p[j][i+1] - p[j][i-1])/user->dx;*/ 182 PetscFunctionReturn(0); 183 } 184 185 /* Second advection term */ 186 PetscErrorCode adv2(PetscScalar **p, PetscScalar x, PetscScalar y, PetscInt i, PetscInt j, PetscInt N, PetscScalar *p2, AppCtx *user) { 187 PetscScalar f, fpos, fneg; 188 PetscFunctionBegin; 189 f = (user->ws / (2 * user->H)) * (user->PM_min - user->Pmax * PetscSinScalar(x) - user->D * (y - user->ws)); 190 fpos = PetscMax(f, 0); 191 fneg = PetscMin(f, 0); 192 if (user->st_width == 1) { 193 *p2 = fpos * (p[j][i] - p[j - 1][i]) / user->dy + fneg * (p[j + 1][i] - p[j][i]) / user->dy; 194 } else if (user->st_width == 2) { 195 *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); 196 } else if (user->st_width == 3) { 197 *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); 198 } 199 200 /* *p2 = f*(p[j+1][i] - p[j-1][i])/user->dy;*/ 201 PetscFunctionReturn(0); 202 } 203 204 /* Diffusion term */ 205 PetscErrorCode diffuse(PetscScalar **p, PetscInt i, PetscInt j, PetscReal t, PetscScalar *p_diff, AppCtx *user) { 206 PetscFunctionBeginUser; 207 if (user->st_width == 1) { 208 *p_diff = user->disper_coe * ((p[j - 1][i] - 2 * p[j][i] + p[j + 1][i]) / (user->dy * user->dy)); 209 } else if (user->st_width == 2) { 210 *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)); 211 } else if (user->st_width == 3) { 212 *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)); 213 } 214 PetscFunctionReturn(0); 215 } 216 217 PetscErrorCode IFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx) { 218 AppCtx *user = (AppCtx *)ctx; 219 DM cda; 220 DMDACoor2d **coors; 221 PetscScalar **p, **f, **pdot; 222 PetscInt i, j; 223 PetscInt xs, ys, xm, ym, M, N; 224 Vec localX, gc, localXdot; 225 PetscScalar p_adv1 = 0.0, p_adv2 = 0.0, p_diff = 0; 226 PetscScalar diffuse1, gamma; 227 228 PetscFunctionBeginUser; 229 PetscCall(DMDAGetInfo(user->da, NULL, &M, &N, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 230 PetscCall(DMGetCoordinateDM(user->da, &cda)); 231 PetscCall(DMDAGetCorners(cda, &xs, &ys, 0, &xm, &ym, 0)); 232 233 PetscCall(DMGetLocalVector(user->da, &localX)); 234 PetscCall(DMGetLocalVector(user->da, &localXdot)); 235 236 PetscCall(DMGlobalToLocalBegin(user->da, X, INSERT_VALUES, localX)); 237 PetscCall(DMGlobalToLocalEnd(user->da, X, INSERT_VALUES, localX)); 238 PetscCall(DMGlobalToLocalBegin(user->da, Xdot, INSERT_VALUES, localXdot)); 239 PetscCall(DMGlobalToLocalEnd(user->da, Xdot, INSERT_VALUES, localXdot)); 240 241 PetscCall(DMGetCoordinatesLocal(user->da, &gc)); 242 243 PetscCall(DMDAVecGetArrayRead(cda, gc, &coors)); 244 PetscCall(DMDAVecGetArrayRead(user->da, localX, &p)); 245 PetscCall(DMDAVecGetArrayRead(user->da, localXdot, &pdot)); 246 PetscCall(DMDAVecGetArray(user->da, F, &f)); 247 248 gamma = user->D * user->ws / (2 * user->H); 249 diffuse1 = user->lambda * user->lambda * user->q / (user->lambda * gamma + 1) * (1.0 - PetscExpScalar(-t * (gamma + 1.0) / user->lambda)); 250 user->disper_coe = user->ws * user->ws / (4 * user->H * user->H) * diffuse1; 251 252 for (i = xs; i < xs + xm; i++) { 253 for (j = ys; j < ys + ym; j++) { 254 PetscCall(adv1(p, coors[j][i].y, i, j, M, &p_adv1, user)); 255 PetscCall(adv2(p, coors[j][i].x, coors[j][i].y, i, j, N, &p_adv2, user)); 256 PetscCall(diffuse(p, i, j, t, &p_diff, user)); 257 f[j][i] = -p_adv1 - p_adv2 + p_diff - pdot[j][i]; 258 } 259 } 260 PetscCall(DMDAVecRestoreArrayRead(user->da, localX, &p)); 261 PetscCall(DMDAVecRestoreArrayRead(user->da, localX, &pdot)); 262 PetscCall(DMRestoreLocalVector(user->da, &localX)); 263 PetscCall(DMRestoreLocalVector(user->da, &localXdot)); 264 PetscCall(DMDAVecRestoreArray(user->da, F, &f)); 265 PetscCall(DMDAVecRestoreArrayRead(cda, gc, &coors)); 266 267 PetscFunctionReturn(0); 268 } 269 270 PetscErrorCode IJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat J, Mat Jpre, void *ctx) { 271 AppCtx *user = (AppCtx *)ctx; 272 DM cda; 273 DMDACoor2d **coors; 274 PetscInt i, j; 275 PetscInt xs, ys, xm, ym, M, N; 276 Vec gc; 277 PetscScalar val[5], xi, yi; 278 MatStencil row, col[5]; 279 PetscScalar c1, c3, c5, c1pos, c1neg, c3pos, c3neg; 280 281 PetscFunctionBeginUser; 282 PetscCall(DMDAGetInfo(user->da, NULL, &M, &N, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 283 PetscCall(DMGetCoordinateDM(user->da, &cda)); 284 PetscCall(DMDAGetCorners(cda, &xs, &ys, 0, &xm, &ym, 0)); 285 286 PetscCall(DMGetCoordinatesLocal(user->da, &gc)); 287 PetscCall(DMDAVecGetArrayRead(cda, gc, &coors)); 288 for (i = xs; i < xs + xm; i++) { 289 for (j = ys; j < ys + ym; j++) { 290 PetscInt nc = 0; 291 xi = coors[j][i].x; 292 yi = coors[j][i].y; 293 row.i = i; 294 row.j = j; 295 c1 = (yi - user->ws) / user->dx; 296 c1pos = PetscMax(c1, 0); 297 c1neg = PetscMin(c1, 0); 298 c3 = (user->ws / (2.0 * user->H)) * (user->PM_min - user->Pmax * PetscSinScalar(xi) - user->D * (yi - user->ws)) / user->dy; 299 c3pos = PetscMax(c3, 0); 300 c3neg = PetscMin(c3, 0); 301 c5 = (PetscPowScalar((user->lambda * user->ws) / (2 * user->H), 2) * user->q * (1.0 - PetscExpScalar(-t / user->lambda))) / (user->dy * user->dy); 302 col[nc].i = i - 1; 303 col[nc].j = j; 304 val[nc++] = c1pos; 305 col[nc].i = i + 1; 306 col[nc].j = j; 307 val[nc++] = -c1neg; 308 col[nc].i = i; 309 col[nc].j = j - 1; 310 val[nc++] = c3pos + c5; 311 col[nc].i = i; 312 col[nc].j = j + 1; 313 val[nc++] = -c3neg + c5; 314 col[nc].i = i; 315 col[nc].j = j; 316 val[nc++] = -c1pos + c1neg - c3pos + c3neg - 2 * c5 - a; 317 PetscCall(MatSetValuesStencil(Jpre, 1, &row, nc, col, val, INSERT_VALUES)); 318 } 319 } 320 PetscCall(DMDAVecRestoreArrayRead(cda, gc, &coors)); 321 322 PetscCall(MatAssemblyBegin(Jpre, MAT_FINAL_ASSEMBLY)); 323 PetscCall(MatAssemblyEnd(Jpre, MAT_FINAL_ASSEMBLY)); 324 if (J != Jpre) { 325 PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 326 PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 327 } 328 PetscFunctionReturn(0); 329 } 330 331 PetscErrorCode Parameter_settings(AppCtx *user) { 332 PetscBool flg; 333 334 PetscFunctionBeginUser; 335 /* Set default parameters */ 336 user->ws = 1.0; 337 user->H = 5.0; 338 user->D = 0.0; 339 user->Pmax = user->Pmax_s = 2.1; 340 user->PM_min = 1.0; 341 user->lambda = 0.1; 342 user->q = 1.0; 343 user->mux = PetscAsinScalar(user->PM_min / user->Pmax); 344 user->sigmax = 0.1; 345 user->sigmay = 0.1; 346 user->rho = 0.0; 347 user->xmin = -PETSC_PI; 348 user->xmax = PETSC_PI; 349 user->bx = DM_BOUNDARY_PERIODIC; 350 user->by = DM_BOUNDARY_GHOSTED; 351 user->tf = user->tcl = -1; 352 user->ymin = -2.0; 353 user->ymax = 2.0; 354 user->st_width = 1; 355 356 PetscCall(PetscOptionsGetScalar(NULL, NULL, "-ws", &user->ws, &flg)); 357 PetscCall(PetscOptionsGetScalar(NULL, NULL, "-Inertia", &user->H, &flg)); 358 PetscCall(PetscOptionsGetScalar(NULL, NULL, "-D", &user->D, &flg)); 359 PetscCall(PetscOptionsGetScalar(NULL, NULL, "-Pmax", &user->Pmax, &flg)); 360 PetscCall(PetscOptionsGetScalar(NULL, NULL, "-PM_min", &user->PM_min, &flg)); 361 PetscCall(PetscOptionsGetScalar(NULL, NULL, "-lambda", &user->lambda, &flg)); 362 PetscCall(PetscOptionsGetScalar(NULL, NULL, "-q", &user->q, &flg)); 363 PetscCall(PetscOptionsGetScalar(NULL, NULL, "-mux", &user->mux, &flg)); 364 PetscCall(PetscOptionsGetScalar(NULL, NULL, "-muy", &user->muy, &flg)); 365 if (flg == 0) user->muy = user->ws; 366 PetscCall(PetscOptionsGetScalar(NULL, NULL, "-xmin", &user->xmin, &flg)); 367 PetscCall(PetscOptionsGetScalar(NULL, NULL, "-xmax", &user->xmax, &flg)); 368 PetscCall(PetscOptionsGetScalar(NULL, NULL, "-ymin", &user->ymin, &flg)); 369 PetscCall(PetscOptionsGetScalar(NULL, NULL, "-ymax", &user->ymax, &flg)); 370 PetscCall(PetscOptionsGetInt(NULL, NULL, "-stencil_width", &user->st_width, &flg)); 371 PetscCall(PetscOptionsGetEnum(NULL, NULL, "-bx", DMBoundaryTypes, (PetscEnum *)&user->bx, &flg)); 372 PetscCall(PetscOptionsGetEnum(NULL, NULL, "-by", DMBoundaryTypes, (PetscEnum *)&user->by, &flg)); 373 PetscCall(PetscOptionsGetReal(NULL, NULL, "-tf", &user->tf, &flg)); 374 PetscCall(PetscOptionsGetReal(NULL, NULL, "-tcl", &user->tcl, &flg)); 375 PetscFunctionReturn(0); 376 } 377 378 /*TEST 379 380 build: 381 requires: !complex x 382 383 test: 384 args: -ts_max_steps 1 385 localrunfiles: petscopt_ex8 386 timeoutfactor: 3 387 388 TEST*/ 389