1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Time-dependent PDE in 2d for calculating joint PDF. \n"; 3c4762a1bSJed Brown /* 4c4762a1bSJed Brown p_t = -x_t*p_x -y_t*p_y + f(t)*p_yy 5c4762a1bSJed Brown xmin < x < xmax, ymin < y < ymax; 6c4762a1bSJed Brown 7c4762a1bSJed Brown Boundary conditions Neumman using mirror values 8c4762a1bSJed Brown 9c4762a1bSJed Brown 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. 10c4762a1bSJed Brown x_t = (y - ws) y_t = (ws/2H)*(Pm - Pmax*sin(x)) 11c4762a1bSJed Brown 12c4762a1bSJed Brown */ 13c4762a1bSJed Brown 14c4762a1bSJed Brown #include <petscdm.h> 15c4762a1bSJed Brown #include <petscdmda.h> 16c4762a1bSJed Brown #include <petscts.h> 17c4762a1bSJed Brown 18c4762a1bSJed Brown /* 19c4762a1bSJed Brown User-defined data structures and routines 20c4762a1bSJed Brown */ 21c4762a1bSJed Brown typedef struct { 22c4762a1bSJed Brown PetscScalar ws; /* Synchronous speed */ 23c4762a1bSJed Brown PetscScalar H; /* Inertia constant */ 24c4762a1bSJed Brown PetscScalar D; /* Damping constant */ 25c4762a1bSJed Brown PetscScalar Pmax; /* Maximum power output of generator */ 26c4762a1bSJed Brown PetscScalar PM_min; /* Mean mechanical power input */ 27c4762a1bSJed Brown PetscScalar lambda; /* correlation time */ 28c4762a1bSJed Brown PetscScalar q; /* noise strength */ 29c4762a1bSJed Brown PetscScalar mux; /* Initial average angle */ 30c4762a1bSJed Brown PetscScalar sigmax; /* Standard deviation of initial angle */ 31c4762a1bSJed Brown PetscScalar muy; /* Average speed */ 32c4762a1bSJed Brown PetscScalar sigmay; /* standard deviation of initial speed */ 33c4762a1bSJed Brown PetscScalar rho; /* Cross-correlation coefficient */ 34c4762a1bSJed Brown PetscScalar xmin; /* left boundary of angle */ 35c4762a1bSJed Brown PetscScalar xmax; /* right boundary of angle */ 36c4762a1bSJed Brown PetscScalar ymin; /* bottom boundary of speed */ 37c4762a1bSJed Brown PetscScalar ymax; /* top boundary of speed */ 38c4762a1bSJed Brown PetscScalar dx; /* x step size */ 39c4762a1bSJed Brown PetscScalar dy; /* y step size */ 40c4762a1bSJed Brown PetscScalar disper_coe; /* Dispersion coefficient */ 41c4762a1bSJed Brown DM da; 42c4762a1bSJed Brown PetscInt st_width; /* Stencil width */ 43c4762a1bSJed Brown DMBoundaryType bx; /* x boundary type */ 44c4762a1bSJed Brown DMBoundaryType by; /* y boundary type */ 45c4762a1bSJed Brown PetscBool nonoiseinitial; 46c4762a1bSJed Brown } AppCtx; 47c4762a1bSJed Brown 48c4762a1bSJed Brown PetscErrorCode Parameter_settings(AppCtx *); 49c4762a1bSJed Brown PetscErrorCode ini_bou(Vec, AppCtx *); 50c4762a1bSJed Brown PetscErrorCode IFunction(TS, PetscReal, Vec, Vec, Vec, void *); 51c4762a1bSJed Brown PetscErrorCode IJacobian(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *); 52c4762a1bSJed Brown PetscErrorCode PostStep(TS); 53c4762a1bSJed Brown 549371c9d4SSatish Balay int main(int argc, char **argv) { 55c4762a1bSJed Brown Vec x; /* Solution vector */ 56c4762a1bSJed Brown TS ts; /* Time-stepping context */ 57c4762a1bSJed Brown AppCtx user; /* Application context */ 58c4762a1bSJed Brown PetscViewer viewer; 59c4762a1bSJed Brown 60327415f7SBarry Smith PetscFunctionBeginUser; 619566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, "petscopt_ex7", help)); 62c4762a1bSJed Brown 63c4762a1bSJed Brown /* Get physics and time parameters */ 649566063dSJacob Faibussowitsch PetscCall(Parameter_settings(&user)); 65c4762a1bSJed Brown /* Create a 2D DA with dof = 1 */ 669566063dSJacob Faibussowitsch 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)); 679566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(user.da)); 689566063dSJacob Faibussowitsch PetscCall(DMSetUp(user.da)); 69c4762a1bSJed Brown /* Set x and y coordinates */ 709566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(user.da, user.xmin, user.xmax, user.ymin, user.ymax, 0, 0)); 719566063dSJacob Faibussowitsch PetscCall(DMDASetCoordinateName(user.da, 0, "X - the angle")); 729566063dSJacob Faibussowitsch PetscCall(DMDASetCoordinateName(user.da, 1, "Y - the speed")); 73c4762a1bSJed Brown 74c4762a1bSJed Brown /* Get global vector x from DM */ 759566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(user.da, &x)); 76c4762a1bSJed Brown 779566063dSJacob Faibussowitsch PetscCall(ini_bou(x, &user)); 789566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "ini_x", FILE_MODE_WRITE, &viewer)); 799566063dSJacob Faibussowitsch PetscCall(VecView(x, viewer)); 809566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 81c4762a1bSJed Brown 829566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 839566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, user.da)); 849566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR)); 859566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSARKIMEX)); 869566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts, NULL, IFunction, &user)); 879566063dSJacob Faibussowitsch /* PetscCall(TSSetIJacobian(ts,NULL,NULL,IJacobian,&user)); */ 889566063dSJacob Faibussowitsch PetscCall(TSSetApplicationContext(ts, &user)); 899566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, .005)); 909566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 919566063dSJacob Faibussowitsch PetscCall(TSSetPostStep(ts, PostStep)); 929566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, x)); 93c4762a1bSJed Brown 949566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "fin_x", FILE_MODE_WRITE, &viewer)); 959566063dSJacob Faibussowitsch PetscCall(VecView(x, viewer)); 969566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 97c4762a1bSJed Brown 989566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 999566063dSJacob Faibussowitsch PetscCall(DMDestroy(&user.da)); 1009566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 1019566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 102b122ec5aSJacob Faibussowitsch return 0; 103c4762a1bSJed Brown } 104c4762a1bSJed Brown 1059371c9d4SSatish Balay PetscErrorCode PostStep(TS ts) { 106c4762a1bSJed Brown Vec X, gc; 107c4762a1bSJed Brown AppCtx *user; 108c4762a1bSJed Brown PetscScalar sum = 0, asum; 109c4762a1bSJed Brown PetscReal t, **p; 110c4762a1bSJed Brown DMDACoor2d **coors; 111c4762a1bSJed Brown DM cda; 112c4762a1bSJed Brown PetscInt i, j, xs, ys, xm, ym; 113c4762a1bSJed Brown 114c4762a1bSJed Brown PetscFunctionBegin; 1159566063dSJacob Faibussowitsch PetscCall(TSGetApplicationContext(ts, &user)); 1169566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &t)); 1179566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ts, &X)); 118c4762a1bSJed Brown 1199566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(user->da, &cda)); 1209566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(cda, &xs, &ys, 0, &xm, &ym, 0)); 1219566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(user->da, &gc)); 1229566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(cda, gc, &coors)); 1239566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(user->da, X, &p)); 124c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 125c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 126c4762a1bSJed Brown if (coors[j][i].y < 5) sum += p[j][i]; 127c4762a1bSJed Brown } 128c4762a1bSJed Brown } 1299566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(cda, gc, &coors)); 1309566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(user->da, X, &p)); 1319566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(&sum, &asum, 1, MPIU_SCALAR, MPIU_SUM, PetscObjectComm((PetscObject)ts))); 1329566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "sum(p)*dw*dtheta at t = %f = %f\n", (double)t, (double)(asum))); 133c4762a1bSJed Brown if (sum < 1.0e-2) { 1349566063dSJacob Faibussowitsch PetscCall(TSSetConvergedReason(ts, TS_CONVERGED_USER)); 1359566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Exiting TS as the integral of PDF is almost zero\n")); 136c4762a1bSJed Brown } 137c4762a1bSJed Brown PetscFunctionReturn(0); 138c4762a1bSJed Brown } 139c4762a1bSJed Brown 1409371c9d4SSatish Balay PetscErrorCode ini_bou(Vec X, AppCtx *user) { 141c4762a1bSJed Brown DM cda; 142c4762a1bSJed Brown DMDACoor2d **coors; 143c4762a1bSJed Brown PetscScalar **p; 144c4762a1bSJed Brown Vec gc; 145c4762a1bSJed Brown PetscInt i, j; 146c4762a1bSJed Brown PetscInt xs, ys, xm, ym, M, N; 147c4762a1bSJed Brown PetscScalar xi, yi; 148c4762a1bSJed Brown PetscScalar sigmax = user->sigmax, sigmay = user->sigmay; 149c4762a1bSJed Brown PetscScalar rho = user->rho; 150c4762a1bSJed Brown PetscScalar muy = user->muy, mux; 151c4762a1bSJed Brown PetscMPIInt rank; 152c4762a1bSJed Brown PetscScalar sum; 153c4762a1bSJed Brown 154c4762a1bSJed Brown PetscFunctionBeginUser; 1559566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 1569566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(user->da, NULL, &M, &N, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 1579371c9d4SSatish Balay user->dx = (user->xmax - user->xmin) / (M - 1); 1589371c9d4SSatish Balay user->dy = (user->ymax - user->ymin) / (N - 1); 1599566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(user->da, &cda)); 1609566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(user->da, &gc)); 1619566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(cda, gc, &coors)); 1629566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(user->da, X, &p)); 1639566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(cda, &xs, &ys, 0, &xm, &ym, 0)); 164c4762a1bSJed Brown 165c4762a1bSJed Brown /* mux and muy need to be grid points in the x and y-direction otherwise the solution goes unstable 166c4762a1bSJed Brown muy is set by choosing the y domain, no. of grid points along y-direction so that muy is a grid point 167c4762a1bSJed Brown in the y-direction. We only modify mux here 168c4762a1bSJed Brown */ 169c4762a1bSJed Brown mux = user->mux = coors[0][M / 2 + 10].x; /* For -pi < x < pi, this should be some angle between 0 and pi/2 */ 170c4762a1bSJed Brown if (user->nonoiseinitial) { 171c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 172c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 1739371c9d4SSatish Balay xi = coors[j][i].x; 1749371c9d4SSatish Balay yi = coors[j][i].y; 175*ad540459SPierre Jolivet if ((xi == mux) && (yi == muy)) p[j][i] = 1.0; 176c4762a1bSJed Brown } 177c4762a1bSJed Brown } 178c4762a1bSJed Brown } else { 179c4762a1bSJed Brown /* Change PM_min accordingly */ 180c4762a1bSJed Brown user->PM_min = user->Pmax * PetscSinScalar(mux); 181c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 182c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 1839371c9d4SSatish Balay xi = coors[j][i].x; 1849371c9d4SSatish Balay yi = coors[j][i].y; 185c4762a1bSJed Brown 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))); 186c4762a1bSJed Brown } 187c4762a1bSJed Brown } 188c4762a1bSJed Brown } 1899566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(cda, gc, &coors)); 1909566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(user->da, X, &p)); 1919566063dSJacob Faibussowitsch PetscCall(VecSum(X, &sum)); 1929566063dSJacob Faibussowitsch PetscCall(VecScale(X, 1.0 / sum)); 193c4762a1bSJed Brown PetscFunctionReturn(0); 194c4762a1bSJed Brown } 195c4762a1bSJed Brown 196c4762a1bSJed Brown /* First advection term */ 1979371c9d4SSatish Balay PetscErrorCode adv1(PetscScalar **p, PetscScalar y, PetscInt i, PetscInt j, PetscInt M, PetscScalar *p1, AppCtx *user) { 198c4762a1bSJed Brown PetscScalar f, fpos, fneg; 199c4762a1bSJed Brown PetscFunctionBegin; 200c4762a1bSJed Brown f = (y - user->ws); 201c4762a1bSJed Brown fpos = PetscMax(f, 0); 202c4762a1bSJed Brown fneg = PetscMin(f, 0); 203c4762a1bSJed Brown if (user->st_width == 1) { 204c4762a1bSJed Brown *p1 = fpos * (p[j][i] - p[j][i - 1]) / user->dx + fneg * (p[j][i + 1] - p[j][i]) / user->dx; 205c4762a1bSJed Brown } else if (user->st_width == 2) { 206c4762a1bSJed Brown *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); 207c4762a1bSJed Brown } else if (user->st_width == 3) { 208c4762a1bSJed Brown *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); 209c4762a1bSJed Brown } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No support for wider stencils"); 210c4762a1bSJed Brown PetscFunctionReturn(0); 211c4762a1bSJed Brown } 212c4762a1bSJed Brown 213c4762a1bSJed Brown /* Second advection term */ 2149371c9d4SSatish Balay PetscErrorCode adv2(PetscScalar **p, PetscScalar x, PetscInt i, PetscInt j, PetscInt N, PetscScalar *p2, AppCtx *user) { 215c4762a1bSJed Brown PetscScalar f, fpos, fneg; 216c4762a1bSJed Brown PetscFunctionBegin; 217c4762a1bSJed Brown f = (user->ws / (2 * user->H)) * (user->PM_min - user->Pmax * PetscSinScalar(x)); 218c4762a1bSJed Brown fpos = PetscMax(f, 0); 219c4762a1bSJed Brown fneg = PetscMin(f, 0); 220c4762a1bSJed Brown if (user->st_width == 1) { 221c4762a1bSJed Brown *p2 = fpos * (p[j][i] - p[j - 1][i]) / user->dy + fneg * (p[j + 1][i] - p[j][i]) / user->dy; 222c4762a1bSJed Brown } else if (user->st_width == 2) { 223c4762a1bSJed Brown *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); 224c4762a1bSJed Brown } else if (user->st_width == 3) { 225c4762a1bSJed Brown *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); 226c4762a1bSJed Brown } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No support for wider stencils"); 227c4762a1bSJed Brown PetscFunctionReturn(0); 228c4762a1bSJed Brown } 229c4762a1bSJed Brown 230c4762a1bSJed Brown /* Diffusion term */ 2319371c9d4SSatish Balay PetscErrorCode diffuse(PetscScalar **p, PetscInt i, PetscInt j, PetscReal t, PetscScalar *p_diff, AppCtx *user) { 232c4762a1bSJed Brown PetscFunctionBeginUser; 233c4762a1bSJed Brown if (user->st_width == 1) { 234c4762a1bSJed Brown *p_diff = user->disper_coe * ((p[j - 1][i] - 2 * p[j][i] + p[j + 1][i]) / (user->dy * user->dy)); 235c4762a1bSJed Brown } else if (user->st_width == 2) { 236c4762a1bSJed Brown *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)); 237c4762a1bSJed Brown } else if (user->st_width == 3) { 238c4762a1bSJed Brown *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)); 239c4762a1bSJed Brown } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No support for wider stencils"); 240c4762a1bSJed Brown PetscFunctionReturn(0); 241c4762a1bSJed Brown } 242c4762a1bSJed Brown 2439371c9d4SSatish Balay PetscErrorCode IFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx) { 244c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx; 245c4762a1bSJed Brown DM cda; 246c4762a1bSJed Brown DMDACoor2d **coors; 247c4762a1bSJed Brown PetscScalar **p, **f, **pdot; 248c4762a1bSJed Brown PetscInt i, j; 249c4762a1bSJed Brown PetscInt xs, ys, xm, ym, M, N; 250c4762a1bSJed Brown Vec localX, gc, localXdot; 251bbadc9eeSBarry Smith PetscScalar p_adv1 = 0.0, p_adv2 = 0.0, p_diff; /* initialize to prevent incorrect compiler warnings */ 252c4762a1bSJed Brown 253c4762a1bSJed Brown PetscFunctionBeginUser; 2549566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(user->da, NULL, &M, &N, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 2559566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(user->da, &cda)); 2569566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(cda, &xs, &ys, 0, &xm, &ym, 0)); 257c4762a1bSJed Brown 2589566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(user->da, &localX)); 2599566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(user->da, &localXdot)); 260c4762a1bSJed Brown 2619566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(user->da, X, INSERT_VALUES, localX)); 2629566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(user->da, X, INSERT_VALUES, localX)); 2639566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(user->da, Xdot, INSERT_VALUES, localXdot)); 2649566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(user->da, Xdot, INSERT_VALUES, localXdot)); 265c4762a1bSJed Brown 2669566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(user->da, &gc)); 267c4762a1bSJed Brown 2689566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(cda, gc, &coors)); 2699566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(user->da, localX, &p)); 2709566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(user->da, localXdot, &pdot)); 2719566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(user->da, F, &f)); 272c4762a1bSJed Brown 273c4762a1bSJed Brown user->disper_coe = PetscPowScalar((user->lambda * user->ws) / (2 * user->H), 2) * user->q * (1.0 - PetscExpScalar(-t / user->lambda)); 274c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 275c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 2769566063dSJacob Faibussowitsch PetscCall(adv1(p, coors[j][i].y, i, j, M, &p_adv1, user)); 2779566063dSJacob Faibussowitsch PetscCall(adv2(p, coors[j][i].x, i, j, N, &p_adv2, user)); 2789566063dSJacob Faibussowitsch PetscCall(diffuse(p, i, j, t, &p_diff, user)); 279c4762a1bSJed Brown f[j][i] = -p_adv1 - p_adv2 + p_diff - pdot[j][i]; 280c4762a1bSJed Brown } 281c4762a1bSJed Brown } 2829566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(user->da, localX, &p)); 2839566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(user->da, localX, &pdot)); 2849566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(user->da, &localX)); 2859566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(user->da, &localXdot)); 2869566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(user->da, F, &f)); 2879566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(cda, gc, &coors)); 288c4762a1bSJed Brown 289c4762a1bSJed Brown PetscFunctionReturn(0); 290c4762a1bSJed Brown } 291c4762a1bSJed Brown 2929371c9d4SSatish Balay PetscErrorCode IJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat J, Mat Jpre, void *ctx) { 293c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx; 294c4762a1bSJed Brown DM cda; 295c4762a1bSJed Brown DMDACoor2d **coors; 296c4762a1bSJed Brown PetscInt i, j; 297c4762a1bSJed Brown PetscInt xs, ys, xm, ym, M, N; 298c4762a1bSJed Brown Vec gc; 299c4762a1bSJed Brown PetscScalar val[5], xi, yi; 300c4762a1bSJed Brown MatStencil row, col[5]; 301c4762a1bSJed Brown PetscScalar c1, c3, c5, c1pos, c1neg, c3pos, c3neg; 302c4762a1bSJed Brown 303c4762a1bSJed Brown PetscFunctionBeginUser; 3049566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(user->da, NULL, &M, &N, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 3059566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(user->da, &cda)); 3069566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(cda, &xs, &ys, 0, &xm, &ym, 0)); 307c4762a1bSJed Brown 3089566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(user->da, &gc)); 3099566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(cda, gc, &coors)); 310c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 311c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 312c4762a1bSJed Brown PetscInt nc = 0; 3139371c9d4SSatish Balay xi = coors[j][i].x; 3149371c9d4SSatish Balay yi = coors[j][i].y; 3159371c9d4SSatish Balay row.i = i; 3169371c9d4SSatish Balay row.j = j; 317c4762a1bSJed Brown c1 = (yi - user->ws) / user->dx; 318c4762a1bSJed Brown c1pos = PetscMax(c1, 0); 319c4762a1bSJed Brown c1neg = PetscMin(c1, 0); 320c4762a1bSJed Brown c3 = (user->ws / (2.0 * user->H)) * (user->PM_min - user->Pmax * PetscSinScalar(xi)) / user->dy; 321c4762a1bSJed Brown c3pos = PetscMax(c3, 0); 322c4762a1bSJed Brown c3neg = PetscMin(c3, 0); 323c4762a1bSJed Brown c5 = (PetscPowScalar((user->lambda * user->ws) / (2 * user->H), 2) * user->q * (1.0 - PetscExpScalar(-t / user->lambda))) / (user->dy * user->dy); 3249371c9d4SSatish Balay col[nc].i = i - 1; 3259371c9d4SSatish Balay col[nc].j = j; 3269371c9d4SSatish Balay val[nc++] = c1pos; 3279371c9d4SSatish Balay col[nc].i = i + 1; 3289371c9d4SSatish Balay col[nc].j = j; 3299371c9d4SSatish Balay val[nc++] = -c1neg; 3309371c9d4SSatish Balay col[nc].i = i; 3319371c9d4SSatish Balay col[nc].j = j - 1; 3329371c9d4SSatish Balay val[nc++] = c3pos + c5; 3339371c9d4SSatish Balay col[nc].i = i; 3349371c9d4SSatish Balay col[nc].j = j + 1; 3359371c9d4SSatish Balay val[nc++] = -c3neg + c5; 3369371c9d4SSatish Balay col[nc].i = i; 3379371c9d4SSatish Balay col[nc].j = j; 3389371c9d4SSatish Balay val[nc++] = -c1pos + c1neg - c3pos + c3neg - 2 * c5 - a; 3399566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(Jpre, 1, &row, nc, col, val, INSERT_VALUES)); 340c4762a1bSJed Brown } 341c4762a1bSJed Brown } 3429566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(cda, gc, &coors)); 343c4762a1bSJed Brown 3449566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Jpre, MAT_FINAL_ASSEMBLY)); 3459566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Jpre, MAT_FINAL_ASSEMBLY)); 346c4762a1bSJed Brown if (J != Jpre) { 3479566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 3489566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 349c4762a1bSJed Brown } 350c4762a1bSJed Brown PetscFunctionReturn(0); 351c4762a1bSJed Brown } 352c4762a1bSJed Brown 3539371c9d4SSatish Balay PetscErrorCode Parameter_settings(AppCtx *user) { 354c4762a1bSJed Brown PetscBool flg; 355c4762a1bSJed Brown 356c4762a1bSJed Brown PetscFunctionBeginUser; 357c4762a1bSJed Brown 358c4762a1bSJed Brown /* Set default parameters */ 359c4762a1bSJed Brown user->ws = 1.0; 360c4762a1bSJed Brown user->H = 5.0; 361c4762a1bSJed Brown user->Pmax = 2.1; 362c4762a1bSJed Brown user->PM_min = 1.0; 363c4762a1bSJed Brown user->lambda = 0.1; 364c4762a1bSJed Brown user->q = 1.0; 365c4762a1bSJed Brown user->mux = PetscAsinScalar(user->PM_min / user->Pmax); 366c4762a1bSJed Brown user->sigmax = 0.1; 367c4762a1bSJed Brown user->sigmay = 0.1; 368c4762a1bSJed Brown user->rho = 0.0; 369c4762a1bSJed Brown user->xmin = -PETSC_PI; 370c4762a1bSJed Brown user->xmax = PETSC_PI; 371c4762a1bSJed Brown user->bx = DM_BOUNDARY_PERIODIC; 372c4762a1bSJed Brown user->by = DM_BOUNDARY_MIRROR; 373c4762a1bSJed Brown user->nonoiseinitial = PETSC_FALSE; 374c4762a1bSJed Brown 375c4762a1bSJed Brown /* 376c4762a1bSJed Brown ymin of -3 seems to let the unstable solution move up and leave a zero in its wake 377c4762a1bSJed Brown with an ymin of -1 the wake is never exactly zero 378c4762a1bSJed Brown */ 379c4762a1bSJed Brown user->ymin = -3.0; 380c4762a1bSJed Brown user->ymax = 10.0; 381c4762a1bSJed Brown user->st_width = 1; 382c4762a1bSJed Brown 3839566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-ws", &user->ws, &flg)); 3849566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-Inertia", &user->H, &flg)); 3859566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-Pmax", &user->Pmax, &flg)); 3869566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-PM_min", &user->PM_min, &flg)); 3879566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-lambda", &user->lambda, &flg)); 3889566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-q", &user->q, &flg)); 3899566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-mux", &user->mux, &flg)); 3909566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-sigmax", &user->sigmax, &flg)); 3919566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-muy", &user->muy, &flg)); 392*ad540459SPierre Jolivet if (flg == 0) user->muy = user->ws; 3939566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-sigmay", &user->sigmay, &flg)); 3949566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-rho", &user->rho, &flg)); 3959566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-xmin", &user->xmin, &flg)); 3969566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-xmax", &user->xmax, &flg)); 3979566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-ymin", &user->ymin, &flg)); 3989566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-ymax", &user->ymax, &flg)); 3999566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-stencil_width", &user->st_width, &flg)); 4009566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetEnum(NULL, NULL, "-bx", DMBoundaryTypes, (PetscEnum *)&user->bx, &flg)); 4019566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetEnum(NULL, NULL, "-by", DMBoundaryTypes, (PetscEnum *)&user->by, &flg)); 4029566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-nonoiseinitial", &user->nonoiseinitial, &flg)); 403c4762a1bSJed Brown 404c4762a1bSJed Brown PetscFunctionReturn(0); 405c4762a1bSJed Brown } 406c4762a1bSJed Brown 407c4762a1bSJed Brown /*TEST 408c4762a1bSJed Brown 409c4762a1bSJed Brown build: 410c4762a1bSJed Brown requires: !complex !single 411c4762a1bSJed Brown 412c4762a1bSJed Brown test: 413c4762a1bSJed Brown args: -ts_max_steps 2 414c4762a1bSJed Brown localrunfiles: petscopt_ex7 415c4762a1bSJed Brown 416c4762a1bSJed Brown test: 417c4762a1bSJed Brown suffix: 2 418c4762a1bSJed Brown args: -ts_max_steps 2 -snes_mf_operator 419c4762a1bSJed Brown output_file: output/ex7_1.out 420c4762a1bSJed Brown localrunfiles: petscopt_ex7 421c4762a1bSJed Brown timeoutfactor: 2 422c4762a1bSJed Brown 423c4762a1bSJed Brown test: 424c4762a1bSJed Brown suffix: 3 425c4762a1bSJed Brown args: -ts_max_steps 2 -snes_mf -pc_type none 426c4762a1bSJed Brown output_file: output/ex7_1.out 427c4762a1bSJed Brown localrunfiles: petscopt_ex7 428c4762a1bSJed Brown timeoutfactor: 2 429c4762a1bSJed Brown 430c4762a1bSJed Brown TEST*/ 431