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 x_t = (y - ws) y_t = (ws/2H)*(Pm - Pmax*sin(x)) 6 7 Boundary conditions: -bc_type 0 => Zero dirichlet boundary 8 -bc_type 1 => Steady state boundary condition 9 Steady state boundary condition found by setting p_t = 0 10 */ 11 12 #include <petscdm.h> 13 #include <petscdmda.h> 14 #include <petscts.h> 15 16 /* 17 User-defined data structures and routines 18 */ 19 typedef struct { 20 PetscScalar ws; /* Synchronous speed */ 21 PetscScalar H; /* Inertia constant */ 22 PetscScalar D; /* Damping constant */ 23 PetscScalar Pmax; /* Maximum power output of generator */ 24 PetscScalar PM_min; /* Mean mechanical power input */ 25 PetscScalar lambda; /* correlation time */ 26 PetscScalar q; /* noise strength */ 27 PetscScalar mux; /* Initial average angle */ 28 PetscScalar sigmax; /* Standard deviation of initial angle */ 29 PetscScalar muy; /* Average speed */ 30 PetscScalar sigmay; /* standard deviation of initial speed */ 31 PetscScalar rho; /* Cross-correlation coefficient */ 32 PetscScalar t0; /* Initial time */ 33 PetscScalar tmax; /* Final time */ 34 PetscScalar xmin; /* left boundary of angle */ 35 PetscScalar xmax; /* right boundary of angle */ 36 PetscScalar ymin; /* bottom boundary of speed */ 37 PetscScalar ymax; /* top boundary of speed */ 38 PetscScalar dx; /* x step size */ 39 PetscScalar dy; /* y step size */ 40 PetscInt bc; /* Boundary conditions */ 41 PetscScalar disper_coe; /* Dispersion coefficient */ 42 DM da; 43 } AppCtx; 44 45 PetscErrorCode Parameter_settings(AppCtx *); 46 PetscErrorCode ini_bou(Vec, AppCtx *); 47 PetscErrorCode IFunction(TS, PetscReal, Vec, Vec, Vec, void *); 48 PetscErrorCode IJacobian(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *); 49 PetscErrorCode PostStep(TS); 50 51 int main(int argc, char **argv) { 52 Vec x; /* Solution vector */ 53 TS ts; /* Time-stepping context */ 54 AppCtx user; /* Application context */ 55 Mat J; 56 PetscViewer viewer; 57 58 PetscFunctionBeginUser; 59 PetscCall(PetscInitialize(&argc, &argv, "petscopt_ex6", help)); 60 /* Get physics and time parameters */ 61 PetscCall(Parameter_settings(&user)); 62 /* Create a 2D DA with dof = 1 */ 63 PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, 4, 4, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, NULL, &user.da)); 64 PetscCall(DMSetFromOptions(user.da)); 65 PetscCall(DMSetUp(user.da)); 66 /* Set x and y coordinates */ 67 PetscCall(DMDASetUniformCoordinates(user.da, user.xmin, user.xmax, user.ymin, user.ymax, 0.0, 1.0)); 68 69 /* Get global vector x from DM */ 70 PetscCall(DMCreateGlobalVector(user.da, &x)); 71 72 PetscCall(ini_bou(x, &user)); 73 PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "ini_x", FILE_MODE_WRITE, &viewer)); 74 PetscCall(VecView(x, viewer)); 75 PetscCall(PetscViewerDestroy(&viewer)); 76 77 /* Get Jacobian matrix structure from the da */ 78 PetscCall(DMSetMatType(user.da, MATAIJ)); 79 PetscCall(DMCreateMatrix(user.da, &J)); 80 81 PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 82 PetscCall(TSSetProblemType(ts, TS_NONLINEAR)); 83 PetscCall(TSSetIFunction(ts, NULL, IFunction, &user)); 84 PetscCall(TSSetIJacobian(ts, J, J, IJacobian, &user)); 85 PetscCall(TSSetApplicationContext(ts, &user)); 86 PetscCall(TSSetMaxTime(ts, user.tmax)); 87 PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 88 PetscCall(TSSetTime(ts, user.t0)); 89 PetscCall(TSSetTimeStep(ts, .005)); 90 PetscCall(TSSetFromOptions(ts)); 91 PetscCall(TSSetPostStep(ts, PostStep)); 92 PetscCall(TSSolve(ts, x)); 93 94 PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "fin_x", FILE_MODE_WRITE, &viewer)); 95 PetscCall(VecView(x, viewer)); 96 PetscCall(PetscViewerDestroy(&viewer)); 97 98 PetscCall(VecDestroy(&x)); 99 PetscCall(MatDestroy(&J)); 100 PetscCall(DMDestroy(&user.da)); 101 PetscCall(TSDestroy(&ts)); 102 PetscCall(PetscFinalize()); 103 return 0; 104 } 105 106 PetscErrorCode PostStep(TS ts) { 107 Vec X; 108 AppCtx *user; 109 PetscScalar sum; 110 PetscReal t; 111 112 PetscFunctionBegin; 113 PetscCall(TSGetApplicationContext(ts, &user)); 114 PetscCall(TSGetTime(ts, &t)); 115 PetscCall(TSGetSolution(ts, &X)); 116 PetscCall(VecSum(X, &sum)); 117 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "sum(p)*dw*dtheta at t = %3.2f = %3.6f\n", (double)t, (double)(sum * user->dx * user->dy))); 118 PetscFunctionReturn(0); 119 } 120 121 PetscErrorCode ini_bou(Vec X, AppCtx *user) { 122 DM cda; 123 DMDACoor2d **coors; 124 PetscScalar **p; 125 Vec gc; 126 PetscInt i, j; 127 PetscInt xs, ys, xm, ym, M, N; 128 PetscScalar xi, yi; 129 PetscScalar sigmax = user->sigmax, sigmay = user->sigmay; 130 PetscScalar rho = user->rho; 131 PetscScalar mux = user->mux, muy = user->muy; 132 PetscMPIInt rank; 133 134 PetscFunctionBeginUser; 135 PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 136 PetscCall(DMDAGetInfo(user->da, NULL, &M, &N, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 137 user->dx = (user->xmax - user->xmin) / (M - 1); 138 user->dy = (user->ymax - user->ymin) / (N - 1); 139 PetscCall(DMGetCoordinateDM(user->da, &cda)); 140 PetscCall(DMGetCoordinates(user->da, &gc)); 141 PetscCall(DMDAVecGetArray(cda, gc, &coors)); 142 PetscCall(DMDAVecGetArray(user->da, X, &p)); 143 PetscCall(DMDAGetCorners(cda, &xs, &ys, 0, &xm, &ym, 0)); 144 for (i = xs; i < xs + xm; i++) { 145 for (j = ys; j < ys + ym; j++) { 146 xi = coors[j][i].x; 147 yi = coors[j][i].y; 148 if (i == 0 || j == 0 || i == M - 1 || j == N - 1) p[j][i] = 0.0; 149 else 150 p[j][i] = (0.5 / (PETSC_PI * sigmax * sigmay * PetscSqrtScalar(1.0 - rho * rho))) * PetscExpScalar(-0.5 / (1 - rho * rho) * (PetscPowScalar((xi - mux) / sigmax, 2) + PetscPowScalar((yi - muy) / sigmay, 2) - 2 * rho * (xi - mux) * (yi - muy) / (sigmax * sigmay))); 151 } 152 } 153 /* p[N/2+N%2][M/2+M%2] = 1/(user->dx*user->dy); */ 154 155 PetscCall(DMDAVecRestoreArray(cda, gc, &coors)); 156 PetscCall(DMDAVecRestoreArray(user->da, X, &p)); 157 PetscFunctionReturn(0); 158 } 159 160 /* First advection term */ 161 PetscErrorCode adv1(PetscScalar **p, PetscScalar y, PetscInt i, PetscInt j, PetscInt M, PetscScalar *p1, AppCtx *user) { 162 PetscScalar f; 163 /* PetscScalar v1,v2,v3,v4,v5,s1,s2,s3; */ 164 PetscFunctionBegin; 165 /* if (i > 2 && i < M-2) { 166 v1 = (y-user->ws)*(p[j][i-2] - p[j][i-3])/user->dx; 167 v2 = (y-user->ws)*(p[j][i-1] - p[j][i-2])/user->dx; 168 v3 = (y-user->ws)*(p[j][i] - p[j][i-1])/user->dx; 169 v4 = (y-user->ws)*(p[j][i+1] - p[j][i])/user->dx; 170 v5 = (y-user->ws)*(p[j][i+1] - p[j][i+2])/user->dx; 171 172 s1 = v1/3.0 - (7.0/6.0)*v2 + (11.0/6.0)*v3; 173 s2 =-v2/6.0 + (5.0/6.0)*v3 + (1.0/3.0)*v4; 174 s3 = v3/3.0 + (5.0/6.0)*v4 - (1.0/6.0)*v5; 175 176 *p1 = 0.1*s1 + 0.6*s2 + 0.3*s3; 177 } else *p1 = 0.0; */ 178 f = (y - user->ws); 179 *p1 = f * (p[j][i + 1] - p[j][i - 1]) / (2 * user->dx); 180 PetscFunctionReturn(0); 181 } 182 183 /* Second advection term */ 184 PetscErrorCode adv2(PetscScalar **p, PetscScalar x, PetscInt i, PetscInt j, PetscInt N, PetscScalar *p2, AppCtx *user) { 185 PetscScalar f; 186 /* PetscScalar v1,v2,v3,v4,v5,s1,s2,s3; */ 187 PetscFunctionBegin; 188 /* if (j > 2 && j < N-2) { 189 v1 = (user->ws/(2*user->H))*(user->PM_min - user->Pmax*sin(x))*(p[j-2][i] - p[j-3][i])/user->dy; 190 v2 = (user->ws/(2*user->H))*(user->PM_min - user->Pmax*sin(x))*(p[j-1][i] - p[j-2][i])/user->dy; 191 v3 = (user->ws/(2*user->H))*(user->PM_min - user->Pmax*sin(x))*(p[j][i] - p[j-1][i])/user->dy; 192 v4 = (user->ws/(2*user->H))*(user->PM_min - user->Pmax*sin(x))*(p[j+1][i] - p[j][i])/user->dy; 193 v5 = (user->ws/(2*user->H))*(user->PM_min - user->Pmax*sin(x))*(p[j+2][i] - p[j+1][i])/user->dy; 194 195 s1 = v1/3.0 - (7.0/6.0)*v2 + (11.0/6.0)*v3; 196 s2 =-v2/6.0 + (5.0/6.0)*v3 + (1.0/3.0)*v4; 197 s3 = v3/3.0 + (5.0/6.0)*v4 - (1.0/6.0)*v5; 198 199 *p2 = 0.1*s1 + 0.6*s2 + 0.3*s3; 200 } else *p2 = 0.0; */ 201 f = (user->ws / (2 * user->H)) * (user->PM_min - user->Pmax * PetscSinScalar(x)); 202 *p2 = f * (p[j + 1][i] - p[j - 1][i]) / (2 * user->dy); 203 PetscFunctionReturn(0); 204 } 205 206 /* Diffusion term */ 207 PetscErrorCode diffuse(PetscScalar **p, PetscInt i, PetscInt j, PetscReal t, PetscScalar *p_diff, AppCtx *user) { 208 PetscFunctionBeginUser; 209 210 *p_diff = user->disper_coe * ((p[j - 1][i] - 2 * p[j][i] + p[j + 1][i]) / (user->dy * user->dy)); 211 PetscFunctionReturn(0); 212 } 213 214 PetscErrorCode BoundaryConditions(PetscScalar **p, DMDACoor2d **coors, PetscInt i, PetscInt j, PetscInt M, PetscInt N, PetscScalar **f, AppCtx *user) { 215 PetscScalar fwc, fthetac; 216 PetscScalar w = coors[j][i].y, theta = coors[j][i].x; 217 218 PetscFunctionBeginUser; 219 if (user->bc == 0) { /* Natural boundary condition */ 220 f[j][i] = p[j][i]; 221 } else { /* Steady state boundary condition */ 222 fthetac = user->ws / (2 * user->H) * (user->PM_min - user->Pmax * PetscSinScalar(theta)); 223 fwc = (w * w / 2.0 - user->ws * w); 224 if (i == 0 && j == 0) { /* left bottom corner */ 225 f[j][i] = fwc * (p[j][i + 1] - p[j][i]) / user->dx + fthetac * p[j][i] - user->disper_coe * (p[j + 1][i] - p[j][i]) / user->dy; 226 } else if (i == 0 && j == N - 1) { /* right bottom corner */ 227 f[j][i] = fwc * (p[j][i + 1] - p[j][i]) / user->dx + fthetac * p[j][i] - user->disper_coe * (p[j][i] - p[j - 1][i]) / user->dy; 228 } else if (i == M - 1 && j == 0) { /* left top corner */ 229 f[j][i] = fwc * (p[j][i] - p[j][i - 1]) / user->dx + fthetac * p[j][i] - user->disper_coe * (p[j + 1][i] - p[j][i]) / user->dy; 230 } else if (i == M - 1 && j == N - 1) { /* right top corner */ 231 f[j][i] = fwc * (p[j][i] - p[j][i - 1]) / user->dx + fthetac * p[j][i] - user->disper_coe * (p[j][i] - p[j - 1][i]) / user->dy; 232 } else if (i == 0) { /* Bottom edge */ 233 f[j][i] = fwc * (p[j][i + 1] - p[j][i]) / (user->dx) + fthetac * p[j][i] - user->disper_coe * (p[j + 1][i] - p[j - 1][i]) / (2 * user->dy); 234 } else if (i == M - 1) { /* Top edge */ 235 f[j][i] = fwc * (p[j][i] - p[j][i - 1]) / (user->dx) + fthetac * p[j][i] - user->disper_coe * (p[j + 1][i] - p[j - 1][i]) / (2 * user->dy); 236 } else if (j == 0) { /* Left edge */ 237 f[j][i] = fwc * (p[j][i + 1] - p[j][i - 1]) / (2 * user->dx) + fthetac * p[j][i] - user->disper_coe * (p[j + 1][i] - p[j][i]) / (user->dy); 238 } else if (j == N - 1) { /* Right edge */ 239 f[j][i] = fwc * (p[j][i + 1] - p[j][i - 1]) / (2 * user->dx) + fthetac * p[j][i] - user->disper_coe * (p[j][i] - p[j - 1][i]) / (user->dy); 240 } 241 } 242 PetscFunctionReturn(0); 243 } 244 245 PetscErrorCode IFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx) { 246 AppCtx *user = (AppCtx *)ctx; 247 DM cda; 248 DMDACoor2d **coors; 249 PetscScalar **p, **f, **pdot; 250 PetscInt i, j; 251 PetscInt xs, ys, xm, ym, M, N; 252 Vec localX, gc, localXdot; 253 PetscScalar p_adv1, p_adv2, p_diff; 254 255 PetscFunctionBeginUser; 256 PetscCall(DMDAGetInfo(user->da, NULL, &M, &N, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 257 PetscCall(DMGetCoordinateDM(user->da, &cda)); 258 PetscCall(DMDAGetCorners(cda, &xs, &ys, 0, &xm, &ym, 0)); 259 260 PetscCall(DMGetLocalVector(user->da, &localX)); 261 PetscCall(DMGetLocalVector(user->da, &localXdot)); 262 263 PetscCall(DMGlobalToLocalBegin(user->da, X, INSERT_VALUES, localX)); 264 PetscCall(DMGlobalToLocalEnd(user->da, X, INSERT_VALUES, localX)); 265 PetscCall(DMGlobalToLocalBegin(user->da, Xdot, INSERT_VALUES, localXdot)); 266 PetscCall(DMGlobalToLocalEnd(user->da, Xdot, INSERT_VALUES, localXdot)); 267 268 PetscCall(DMGetCoordinatesLocal(user->da, &gc)); 269 270 PetscCall(DMDAVecGetArrayRead(cda, gc, &coors)); 271 PetscCall(DMDAVecGetArrayRead(user->da, localX, &p)); 272 PetscCall(DMDAVecGetArrayRead(user->da, localXdot, &pdot)); 273 PetscCall(DMDAVecGetArray(user->da, F, &f)); 274 275 user->disper_coe = PetscPowScalar((user->lambda * user->ws) / (2 * user->H), 2) * user->q * (1.0 - PetscExpScalar(-t / user->lambda)); 276 for (i = xs; i < xs + xm; i++) { 277 for (j = ys; j < ys + ym; j++) { 278 if (i == 0 || j == 0 || i == M - 1 || j == N - 1) { 279 PetscCall(BoundaryConditions(p, coors, i, j, M, N, f, user)); 280 } else { 281 PetscCall(adv1(p, coors[j][i].y, i, j, M, &p_adv1, user)); 282 PetscCall(adv2(p, coors[j][i].x, i, j, N, &p_adv2, user)); 283 PetscCall(diffuse(p, i, j, t, &p_diff, user)); 284 f[j][i] = -p_adv1 - p_adv2 + p_diff - pdot[j][i]; 285 } 286 } 287 } 288 PetscCall(DMDAVecRestoreArrayRead(user->da, localX, &p)); 289 PetscCall(DMDAVecRestoreArrayRead(user->da, localX, &pdot)); 290 PetscCall(DMRestoreLocalVector(user->da, &localX)); 291 PetscCall(DMRestoreLocalVector(user->da, &localXdot)); 292 PetscCall(DMDAVecRestoreArray(user->da, F, &f)); 293 PetscCall(DMDAVecRestoreArrayRead(cda, gc, &coors)); 294 295 PetscFunctionReturn(0); 296 } 297 298 PetscErrorCode IJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat J, Mat Jpre, void *ctx) { 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; 308 309 PetscFunctionBeginUser; 310 PetscCall(DMDAGetInfo(user->da, NULL, &M, &N, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 311 PetscCall(DMGetCoordinateDM(user->da, &cda)); 312 PetscCall(DMDAGetCorners(cda, &xs, &ys, 0, &xm, &ym, 0)); 313 314 PetscCall(DMGetCoordinatesLocal(user->da, &gc)); 315 PetscCall(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; 320 yi = coors[j][i].y; 321 row.i = i; 322 row.j = j; 323 if (i == 0 || j == 0 || i == M - 1 || j == N - 1) { 324 if (user->bc == 0) { 325 col[nc].i = i; 326 col[nc].j = j; 327 val[nc++] = 1.0; 328 } else { 329 PetscScalar fthetac, fwc; 330 fthetac = user->ws / (2 * user->H) * (user->PM_min - user->Pmax * PetscSinScalar(xi)); 331 fwc = (yi * yi / 2.0 - user->ws * yi); 332 if (i == 0 && j == 0) { 333 col[nc].i = i + 1; 334 col[nc].j = j; 335 val[nc++] = fwc / user->dx; 336 col[nc].i = i; 337 col[nc].j = j + 1; 338 val[nc++] = -user->disper_coe / user->dy; 339 col[nc].i = i; 340 col[nc].j = j; 341 val[nc++] = -fwc / user->dx + fthetac + user->disper_coe / user->dy; 342 } else if (i == 0 && j == N - 1) { 343 col[nc].i = i + 1; 344 col[nc].j = j; 345 val[nc++] = fwc / user->dx; 346 col[nc].i = i; 347 col[nc].j = j - 1; 348 val[nc++] = user->disper_coe / user->dy; 349 col[nc].i = i; 350 col[nc].j = j; 351 val[nc++] = -fwc / user->dx + fthetac - user->disper_coe / user->dy; 352 } else if (i == M - 1 && j == 0) { 353 col[nc].i = i - 1; 354 col[nc].j = j; 355 val[nc++] = -fwc / user->dx; 356 col[nc].i = i; 357 col[nc].j = j + 1; 358 val[nc++] = -user->disper_coe / user->dy; 359 col[nc].i = i; 360 col[nc].j = j; 361 val[nc++] = fwc / user->dx + fthetac + user->disper_coe / user->dy; 362 } else if (i == M - 1 && j == N - 1) { 363 col[nc].i = i - 1; 364 col[nc].j = j; 365 val[nc++] = -fwc / user->dx; 366 col[nc].i = i; 367 col[nc].j = j - 1; 368 val[nc++] = user->disper_coe / user->dy; 369 col[nc].i = i; 370 col[nc].j = j; 371 val[nc++] = fwc / user->dx + fthetac - user->disper_coe / user->dy; 372 } else if (i == 0) { 373 col[nc].i = i + 1; 374 col[nc].j = j; 375 val[nc++] = fwc / user->dx; 376 col[nc].i = i; 377 col[nc].j = j + 1; 378 val[nc++] = -user->disper_coe / (2 * user->dy); 379 col[nc].i = i; 380 col[nc].j = j - 1; 381 val[nc++] = user->disper_coe / (2 * user->dy); 382 col[nc].i = i; 383 col[nc].j = j; 384 val[nc++] = -fwc / user->dx + fthetac; 385 } else if (i == M - 1) { 386 col[nc].i = i - 1; 387 col[nc].j = j; 388 val[nc++] = -fwc / user->dx; 389 col[nc].i = i; 390 col[nc].j = j + 1; 391 val[nc++] = -user->disper_coe / (2 * user->dy); 392 col[nc].i = i; 393 col[nc].j = j - 1; 394 val[nc++] = user->disper_coe / (2 * user->dy); 395 col[nc].i = i; 396 col[nc].j = j; 397 val[nc++] = fwc / user->dx + fthetac; 398 } else if (j == 0) { 399 col[nc].i = i + 1; 400 col[nc].j = j; 401 val[nc++] = fwc / (2 * user->dx); 402 col[nc].i = i - 1; 403 col[nc].j = j; 404 val[nc++] = -fwc / (2 * user->dx); 405 col[nc].i = i; 406 col[nc].j = j + 1; 407 val[nc++] = -user->disper_coe / user->dy; 408 col[nc].i = i; 409 col[nc].j = j; 410 val[nc++] = user->disper_coe / user->dy + fthetac; 411 } else if (j == N - 1) { 412 col[nc].i = i + 1; 413 col[nc].j = j; 414 val[nc++] = fwc / (2 * user->dx); 415 col[nc].i = i - 1; 416 col[nc].j = j; 417 val[nc++] = -fwc / (2 * user->dx); 418 col[nc].i = i; 419 col[nc].j = j - 1; 420 val[nc++] = user->disper_coe / user->dy; 421 col[nc].i = i; 422 col[nc].j = j; 423 val[nc++] = -user->disper_coe / user->dy + fthetac; 424 } 425 } 426 } else { 427 c1 = (yi - user->ws) / (2 * user->dx); 428 c3 = (user->ws / (2.0 * user->H)) * (user->PM_min - user->Pmax * PetscSinScalar(xi)) / (2 * user->dy); 429 c5 = (PetscPowScalar((user->lambda * user->ws) / (2 * user->H), 2) * user->q * (1.0 - PetscExpScalar(-t / user->lambda))) / (user->dy * user->dy); 430 col[nc].i = i - 1; 431 col[nc].j = j; 432 val[nc++] = c1; 433 col[nc].i = i + 1; 434 col[nc].j = j; 435 val[nc++] = -c1; 436 col[nc].i = i; 437 col[nc].j = j - 1; 438 val[nc++] = c3 + c5; 439 col[nc].i = i; 440 col[nc].j = j + 1; 441 val[nc++] = -c3 + c5; 442 col[nc].i = i; 443 col[nc].j = j; 444 val[nc++] = -2 * c5 - a; 445 } 446 PetscCall(MatSetValuesStencil(Jpre, 1, &row, nc, col, val, INSERT_VALUES)); 447 } 448 } 449 PetscCall(DMDAVecRestoreArrayRead(cda, gc, &coors)); 450 451 PetscCall(MatAssemblyBegin(Jpre, MAT_FINAL_ASSEMBLY)); 452 PetscCall(MatAssemblyEnd(Jpre, MAT_FINAL_ASSEMBLY)); 453 if (J != Jpre) { 454 PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 455 PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 456 } 457 PetscFunctionReturn(0); 458 } 459 460 PetscErrorCode Parameter_settings(AppCtx *user) { 461 PetscBool flg; 462 463 PetscFunctionBeginUser; 464 465 /* Set default parameters */ 466 user->ws = 1.0; 467 user->H = 5.0; 468 user->Pmax = 2.1; 469 user->PM_min = 1.0; 470 user->lambda = 0.1; 471 user->q = 1.0; 472 user->mux = PetscAsinScalar(user->PM_min / user->Pmax); 473 user->sigmax = 0.1; 474 user->sigmay = 0.1; 475 user->rho = 0.0; 476 user->t0 = 0.0; 477 user->tmax = 2.0; 478 user->xmin = -1.0; 479 user->xmax = 10.0; 480 user->ymin = -1.0; 481 user->ymax = 10.0; 482 user->bc = 0; 483 484 PetscCall(PetscOptionsGetScalar(NULL, NULL, "-ws", &user->ws, &flg)); 485 PetscCall(PetscOptionsGetScalar(NULL, NULL, "-Inertia", &user->H, &flg)); 486 PetscCall(PetscOptionsGetScalar(NULL, NULL, "-Pmax", &user->Pmax, &flg)); 487 PetscCall(PetscOptionsGetScalar(NULL, NULL, "-PM_min", &user->PM_min, &flg)); 488 PetscCall(PetscOptionsGetScalar(NULL, NULL, "-lambda", &user->lambda, &flg)); 489 PetscCall(PetscOptionsGetScalar(NULL, NULL, "-q", &user->q, &flg)); 490 PetscCall(PetscOptionsGetScalar(NULL, NULL, "-mux", &user->mux, &flg)); 491 PetscCall(PetscOptionsGetScalar(NULL, NULL, "-sigmax", &user->sigmax, &flg)); 492 PetscCall(PetscOptionsGetScalar(NULL, NULL, "-muy", &user->muy, &flg)); 493 PetscCall(PetscOptionsGetScalar(NULL, NULL, "-sigmay", &user->sigmay, &flg)); 494 PetscCall(PetscOptionsGetScalar(NULL, NULL, "-rho", &user->rho, &flg)); 495 PetscCall(PetscOptionsGetScalar(NULL, NULL, "-t0", &user->t0, &flg)); 496 PetscCall(PetscOptionsGetScalar(NULL, NULL, "-tmax", &user->tmax, &flg)); 497 PetscCall(PetscOptionsGetScalar(NULL, NULL, "-xmin", &user->xmin, &flg)); 498 PetscCall(PetscOptionsGetScalar(NULL, NULL, "-xmax", &user->xmax, &flg)); 499 PetscCall(PetscOptionsGetScalar(NULL, NULL, "-ymin", &user->ymin, &flg)); 500 PetscCall(PetscOptionsGetScalar(NULL, NULL, "-ymax", &user->ymax, &flg)); 501 PetscCall(PetscOptionsGetInt(NULL, NULL, "-bc_type", &user->bc, &flg)); 502 user->muy = user->ws; 503 PetscFunctionReturn(0); 504 } 505 506 /*TEST 507 508 build: 509 requires: !complex 510 511 test: 512 args: -nox -ts_max_steps 2 513 localrunfiles: petscopt_ex6 514 timeoutfactor: 4 515 516 TEST*/ 517