1c4762a1bSJed Brown static char help[] = "Time-dependent PDE in 2d for calculating joint PDF. \n"; 2c4762a1bSJed Brown /* 3c4762a1bSJed Brown p_t = -x_t*p_x -y_t*p_y + f(t)*p_yy 4c4762a1bSJed Brown xmin < x < xmax, ymin < y < ymax; 5c4762a1bSJed Brown 6c4762a1bSJed Brown Boundary conditions Neumman using mirror values 7c4762a1bSJed Brown 8c4762a1bSJed 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. 9c4762a1bSJed Brown x_t = (y - ws) y_t = (ws/2H)*(Pm - Pmax*sin(x)) 10c4762a1bSJed Brown 11c4762a1bSJed Brown */ 12c4762a1bSJed Brown 13c4762a1bSJed Brown #include <petscdm.h> 14c4762a1bSJed Brown #include <petscdmda.h> 15c4762a1bSJed Brown #include <petscts.h> 16c4762a1bSJed Brown 17c4762a1bSJed Brown /* 18c4762a1bSJed Brown User-defined data structures and routines 19c4762a1bSJed Brown */ 20c4762a1bSJed Brown typedef struct { 21c4762a1bSJed Brown PetscScalar ws; /* Synchronous speed */ 22c4762a1bSJed Brown PetscScalar H; /* Inertia constant */ 23c4762a1bSJed Brown PetscScalar D; /* Damping constant */ 24c4762a1bSJed Brown PetscScalar Pmax; /* Maximum power output of generator */ 25c4762a1bSJed Brown PetscScalar PM_min; /* Mean mechanical power input */ 26c4762a1bSJed Brown PetscScalar lambda; /* correlation time */ 27c4762a1bSJed Brown PetscScalar q; /* noise strength */ 28c4762a1bSJed Brown PetscScalar mux; /* Initial average angle */ 29c4762a1bSJed Brown PetscScalar sigmax; /* Standard deviation of initial angle */ 30c4762a1bSJed Brown PetscScalar muy; /* Average speed */ 31c4762a1bSJed Brown PetscScalar sigmay; /* standard deviation of initial speed */ 32c4762a1bSJed Brown PetscScalar rho; /* Cross-correlation coefficient */ 33c4762a1bSJed Brown PetscScalar xmin; /* left boundary of angle */ 34c4762a1bSJed Brown PetscScalar xmax; /* right boundary of angle */ 35c4762a1bSJed Brown PetscScalar ymin; /* bottom boundary of speed */ 36c4762a1bSJed Brown PetscScalar ymax; /* top boundary of speed */ 37c4762a1bSJed Brown PetscScalar dx; /* x step size */ 38c4762a1bSJed Brown PetscScalar dy; /* y step size */ 39c4762a1bSJed Brown PetscScalar disper_coe; /* Dispersion coefficient */ 40c4762a1bSJed Brown DM da; 41c4762a1bSJed Brown PetscInt st_width; /* Stencil width */ 42c4762a1bSJed Brown DMBoundaryType bx; /* x boundary type */ 43c4762a1bSJed Brown DMBoundaryType by; /* y boundary type */ 44c4762a1bSJed Brown PetscBool nonoiseinitial; 45c4762a1bSJed Brown } AppCtx; 46c4762a1bSJed Brown 47c4762a1bSJed Brown PetscErrorCode Parameter_settings(AppCtx *); 48c4762a1bSJed Brown PetscErrorCode ini_bou(Vec, AppCtx *); 49c4762a1bSJed Brown PetscErrorCode IFunction(TS, PetscReal, Vec, Vec, Vec, void *); 50c4762a1bSJed Brown PetscErrorCode IJacobian(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *); 51c4762a1bSJed Brown PetscErrorCode PostStep(TS); 52c4762a1bSJed Brown 53d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 54d71ae5a4SJacob Faibussowitsch { 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 105d71ae5a4SJacob Faibussowitsch PetscErrorCode PostStep(TS ts) 106d71ae5a4SJacob Faibussowitsch { 107c4762a1bSJed Brown Vec X, gc; 108c4762a1bSJed Brown AppCtx *user; 109c4762a1bSJed Brown PetscScalar sum = 0, asum; 110c4762a1bSJed Brown PetscReal t, **p; 111c4762a1bSJed Brown DMDACoor2d **coors; 112c4762a1bSJed Brown DM cda; 113c4762a1bSJed Brown PetscInt i, j, xs, ys, xm, ym; 114c4762a1bSJed Brown 115c4762a1bSJed Brown PetscFunctionBegin; 1169566063dSJacob Faibussowitsch PetscCall(TSGetApplicationContext(ts, &user)); 1179566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &t)); 1189566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ts, &X)); 119c4762a1bSJed Brown 1209566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(user->da, &cda)); 1219566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(cda, &xs, &ys, 0, &xm, &ym, 0)); 1229566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(user->da, &gc)); 1239566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(cda, gc, &coors)); 1249566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(user->da, X, &p)); 125c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 126c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 127c4762a1bSJed Brown if (coors[j][i].y < 5) sum += p[j][i]; 128c4762a1bSJed Brown } 129c4762a1bSJed Brown } 1309566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(cda, gc, &coors)); 1319566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(user->da, X, &p)); 132712fec58SPierre Jolivet PetscCall(MPIU_Allreduce(&sum, &asum, 1, MPIU_SCALAR, MPIU_SUM, PetscObjectComm((PetscObject)ts))); 1339566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "sum(p)*dw*dtheta at t = %f = %f\n", (double)t, (double)(asum))); 134c4762a1bSJed Brown if (sum < 1.0e-2) { 1359566063dSJacob Faibussowitsch PetscCall(TSSetConvergedReason(ts, TS_CONVERGED_USER)); 1369566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Exiting TS as the integral of PDF is almost zero\n")); 137c4762a1bSJed Brown } 1383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 139c4762a1bSJed Brown } 140c4762a1bSJed Brown 141d71ae5a4SJacob Faibussowitsch PetscErrorCode ini_bou(Vec X, AppCtx *user) 142d71ae5a4SJacob Faibussowitsch { 143c4762a1bSJed Brown DM cda; 144c4762a1bSJed Brown DMDACoor2d **coors; 145c4762a1bSJed Brown PetscScalar **p; 146c4762a1bSJed Brown Vec gc; 147c4762a1bSJed Brown PetscInt i, j; 148c4762a1bSJed Brown PetscInt xs, ys, xm, ym, M, N; 149c4762a1bSJed Brown PetscScalar xi, yi; 150c4762a1bSJed Brown PetscScalar sigmax = user->sigmax, sigmay = user->sigmay; 151c4762a1bSJed Brown PetscScalar rho = user->rho; 152c4762a1bSJed Brown PetscScalar muy = user->muy, mux; 153c4762a1bSJed Brown PetscMPIInt rank; 154c4762a1bSJed Brown PetscScalar sum; 155c4762a1bSJed Brown 156c4762a1bSJed Brown PetscFunctionBeginUser; 1579566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 1589566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(user->da, NULL, &M, &N, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 1599371c9d4SSatish Balay user->dx = (user->xmax - user->xmin) / (M - 1); 1609371c9d4SSatish Balay user->dy = (user->ymax - user->ymin) / (N - 1); 1619566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(user->da, &cda)); 1629566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(user->da, &gc)); 1639566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(cda, gc, &coors)); 1649566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(user->da, X, &p)); 1659566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(cda, &xs, &ys, 0, &xm, &ym, 0)); 166c4762a1bSJed Brown 167c4762a1bSJed Brown /* mux and muy need to be grid points in the x and y-direction otherwise the solution goes unstable 168c4762a1bSJed Brown muy is set by choosing the y domain, no. of grid points along y-direction so that muy is a grid point 169c4762a1bSJed Brown in the y-direction. We only modify mux here 170c4762a1bSJed Brown */ 171c4762a1bSJed Brown mux = user->mux = coors[0][M / 2 + 10].x; /* For -pi < x < pi, this should be some angle between 0 and pi/2 */ 172c4762a1bSJed Brown if (user->nonoiseinitial) { 173c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 174c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 1759371c9d4SSatish Balay xi = coors[j][i].x; 1769371c9d4SSatish Balay yi = coors[j][i].y; 177ad540459SPierre Jolivet if ((xi == mux) && (yi == muy)) p[j][i] = 1.0; 178c4762a1bSJed Brown } 179c4762a1bSJed Brown } 180c4762a1bSJed Brown } else { 181c4762a1bSJed Brown /* Change PM_min accordingly */ 182c4762a1bSJed Brown user->PM_min = user->Pmax * PetscSinScalar(mux); 183c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 184c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 1859371c9d4SSatish Balay xi = coors[j][i].x; 1869371c9d4SSatish Balay yi = coors[j][i].y; 187c4762a1bSJed 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))); 188c4762a1bSJed Brown } 189c4762a1bSJed Brown } 190c4762a1bSJed Brown } 1919566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(cda, gc, &coors)); 1929566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(user->da, X, &p)); 1939566063dSJacob Faibussowitsch PetscCall(VecSum(X, &sum)); 1949566063dSJacob Faibussowitsch PetscCall(VecScale(X, 1.0 / sum)); 1953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 196c4762a1bSJed Brown } 197c4762a1bSJed Brown 198c4762a1bSJed Brown /* First advection term */ 199d71ae5a4SJacob Faibussowitsch PetscErrorCode adv1(PetscScalar **p, PetscScalar y, PetscInt i, PetscInt j, PetscInt M, PetscScalar *p1, AppCtx *user) 200d71ae5a4SJacob Faibussowitsch { 201c4762a1bSJed Brown PetscScalar f, fpos, fneg; 202*4d86920dSPierre Jolivet 203c4762a1bSJed Brown PetscFunctionBegin; 204c4762a1bSJed Brown f = (y - user->ws); 205c4762a1bSJed Brown fpos = PetscMax(f, 0); 206c4762a1bSJed Brown fneg = PetscMin(f, 0); 207c4762a1bSJed Brown if (user->st_width == 1) { 208c4762a1bSJed Brown *p1 = fpos * (p[j][i] - p[j][i - 1]) / user->dx + fneg * (p[j][i + 1] - p[j][i]) / user->dx; 209c4762a1bSJed Brown } else if (user->st_width == 2) { 210c4762a1bSJed 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); 211c4762a1bSJed Brown } else if (user->st_width == 3) { 212c4762a1bSJed 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); 213c4762a1bSJed Brown } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No support for wider stencils"); 2143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 215c4762a1bSJed Brown } 216c4762a1bSJed Brown 217c4762a1bSJed Brown /* Second advection term */ 218d71ae5a4SJacob Faibussowitsch PetscErrorCode adv2(PetscScalar **p, PetscScalar x, PetscInt i, PetscInt j, PetscInt N, PetscScalar *p2, AppCtx *user) 219d71ae5a4SJacob Faibussowitsch { 220c4762a1bSJed Brown PetscScalar f, fpos, fneg; 221*4d86920dSPierre Jolivet 222c4762a1bSJed Brown PetscFunctionBegin; 223c4762a1bSJed Brown f = (user->ws / (2 * user->H)) * (user->PM_min - user->Pmax * PetscSinScalar(x)); 224c4762a1bSJed Brown fpos = PetscMax(f, 0); 225c4762a1bSJed Brown fneg = PetscMin(f, 0); 226c4762a1bSJed Brown if (user->st_width == 1) { 227c4762a1bSJed Brown *p2 = fpos * (p[j][i] - p[j - 1][i]) / user->dy + fneg * (p[j + 1][i] - p[j][i]) / user->dy; 228c4762a1bSJed Brown } else if (user->st_width == 2) { 229c4762a1bSJed 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); 230c4762a1bSJed Brown } else if (user->st_width == 3) { 231c4762a1bSJed 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); 232c4762a1bSJed Brown } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No support for wider stencils"); 2333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 234c4762a1bSJed Brown } 235c4762a1bSJed Brown 236c4762a1bSJed Brown /* Diffusion term */ 237d71ae5a4SJacob Faibussowitsch PetscErrorCode diffuse(PetscScalar **p, PetscInt i, PetscInt j, PetscReal t, PetscScalar *p_diff, AppCtx *user) 238d71ae5a4SJacob Faibussowitsch { 239c4762a1bSJed Brown PetscFunctionBeginUser; 240c4762a1bSJed Brown if (user->st_width == 1) { 241c4762a1bSJed Brown *p_diff = user->disper_coe * ((p[j - 1][i] - 2 * p[j][i] + p[j + 1][i]) / (user->dy * user->dy)); 242c4762a1bSJed Brown } else if (user->st_width == 2) { 243c4762a1bSJed 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)); 244c4762a1bSJed Brown } else if (user->st_width == 3) { 245c4762a1bSJed 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)); 246c4762a1bSJed Brown } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No support for wider stencils"); 2473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 248c4762a1bSJed Brown } 249c4762a1bSJed Brown 250d71ae5a4SJacob Faibussowitsch PetscErrorCode IFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx) 251d71ae5a4SJacob Faibussowitsch { 252c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx; 253c4762a1bSJed Brown DM cda; 254c4762a1bSJed Brown DMDACoor2d **coors; 255c4762a1bSJed Brown PetscScalar **p, **f, **pdot; 256c4762a1bSJed Brown PetscInt i, j; 257c4762a1bSJed Brown PetscInt xs, ys, xm, ym, M, N; 258c4762a1bSJed Brown Vec localX, gc, localXdot; 259bbadc9eeSBarry Smith PetscScalar p_adv1 = 0.0, p_adv2 = 0.0, p_diff; /* initialize to prevent incorrect compiler warnings */ 260c4762a1bSJed Brown 261c4762a1bSJed Brown PetscFunctionBeginUser; 2629566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(user->da, NULL, &M, &N, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 2639566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(user->da, &cda)); 2649566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(cda, &xs, &ys, 0, &xm, &ym, 0)); 265c4762a1bSJed Brown 2669566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(user->da, &localX)); 2679566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(user->da, &localXdot)); 268c4762a1bSJed Brown 2699566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(user->da, X, INSERT_VALUES, localX)); 2709566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(user->da, X, INSERT_VALUES, localX)); 2719566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(user->da, Xdot, INSERT_VALUES, localXdot)); 2729566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(user->da, Xdot, INSERT_VALUES, localXdot)); 273c4762a1bSJed Brown 2749566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(user->da, &gc)); 275c4762a1bSJed Brown 2769566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(cda, gc, &coors)); 2779566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(user->da, localX, &p)); 2789566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(user->da, localXdot, &pdot)); 2799566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(user->da, F, &f)); 280c4762a1bSJed Brown 281c4762a1bSJed Brown user->disper_coe = PetscPowScalar((user->lambda * user->ws) / (2 * user->H), 2) * user->q * (1.0 - PetscExpScalar(-t / user->lambda)); 282c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 283c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 2849566063dSJacob Faibussowitsch PetscCall(adv1(p, coors[j][i].y, i, j, M, &p_adv1, user)); 2859566063dSJacob Faibussowitsch PetscCall(adv2(p, coors[j][i].x, i, j, N, &p_adv2, user)); 2869566063dSJacob Faibussowitsch PetscCall(diffuse(p, i, j, t, &p_diff, user)); 287c4762a1bSJed Brown f[j][i] = -p_adv1 - p_adv2 + p_diff - pdot[j][i]; 288c4762a1bSJed Brown } 289c4762a1bSJed Brown } 2909566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(user->da, localX, &p)); 2919566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(user->da, localX, &pdot)); 2929566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(user->da, &localX)); 2939566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(user->da, &localXdot)); 2949566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(user->da, F, &f)); 2959566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(cda, gc, &coors)); 2963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 297c4762a1bSJed Brown } 298c4762a1bSJed Brown 299d71ae5a4SJacob Faibussowitsch PetscErrorCode IJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat J, Mat Jpre, void *ctx) 300d71ae5a4SJacob Faibussowitsch { 301c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx; 302c4762a1bSJed Brown DM cda; 303c4762a1bSJed Brown DMDACoor2d **coors; 304c4762a1bSJed Brown PetscInt i, j; 305c4762a1bSJed Brown PetscInt xs, ys, xm, ym, M, N; 306c4762a1bSJed Brown Vec gc; 307c4762a1bSJed Brown PetscScalar val[5], xi, yi; 308c4762a1bSJed Brown MatStencil row, col[5]; 309c4762a1bSJed Brown PetscScalar c1, c3, c5, c1pos, c1neg, c3pos, c3neg; 310c4762a1bSJed Brown 311c4762a1bSJed Brown PetscFunctionBeginUser; 3129566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(user->da, NULL, &M, &N, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 3139566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(user->da, &cda)); 3149566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(cda, &xs, &ys, 0, &xm, &ym, 0)); 315c4762a1bSJed Brown 3169566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(user->da, &gc)); 3179566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(cda, gc, &coors)); 318c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 319c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 320c4762a1bSJed Brown PetscInt nc = 0; 3219371c9d4SSatish Balay xi = coors[j][i].x; 3229371c9d4SSatish Balay yi = coors[j][i].y; 3239371c9d4SSatish Balay row.i = i; 3249371c9d4SSatish Balay row.j = j; 325c4762a1bSJed Brown c1 = (yi - user->ws) / user->dx; 326c4762a1bSJed Brown c1pos = PetscMax(c1, 0); 327c4762a1bSJed Brown c1neg = PetscMin(c1, 0); 328c4762a1bSJed Brown c3 = (user->ws / (2.0 * user->H)) * (user->PM_min - user->Pmax * PetscSinScalar(xi)) / user->dy; 329c4762a1bSJed Brown c3pos = PetscMax(c3, 0); 330c4762a1bSJed Brown c3neg = PetscMin(c3, 0); 331c4762a1bSJed Brown c5 = (PetscPowScalar((user->lambda * user->ws) / (2 * user->H), 2) * user->q * (1.0 - PetscExpScalar(-t / user->lambda))) / (user->dy * user->dy); 3329371c9d4SSatish Balay col[nc].i = i - 1; 3339371c9d4SSatish Balay col[nc].j = j; 3349371c9d4SSatish Balay val[nc++] = c1pos; 3359371c9d4SSatish Balay col[nc].i = i + 1; 3369371c9d4SSatish Balay col[nc].j = j; 3379371c9d4SSatish Balay val[nc++] = -c1neg; 3389371c9d4SSatish Balay col[nc].i = i; 3399371c9d4SSatish Balay col[nc].j = j - 1; 3409371c9d4SSatish Balay val[nc++] = c3pos + c5; 3419371c9d4SSatish Balay col[nc].i = i; 3429371c9d4SSatish Balay col[nc].j = j + 1; 3439371c9d4SSatish Balay val[nc++] = -c3neg + c5; 3449371c9d4SSatish Balay col[nc].i = i; 3459371c9d4SSatish Balay col[nc].j = j; 3469371c9d4SSatish Balay val[nc++] = -c1pos + c1neg - c3pos + c3neg - 2 * c5 - a; 3479566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(Jpre, 1, &row, nc, col, val, INSERT_VALUES)); 348c4762a1bSJed Brown } 349c4762a1bSJed Brown } 3509566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(cda, gc, &coors)); 351c4762a1bSJed Brown 3529566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Jpre, MAT_FINAL_ASSEMBLY)); 3539566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Jpre, MAT_FINAL_ASSEMBLY)); 354c4762a1bSJed Brown if (J != Jpre) { 3559566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 3569566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 357c4762a1bSJed Brown } 3583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 359c4762a1bSJed Brown } 360c4762a1bSJed Brown 361d71ae5a4SJacob Faibussowitsch PetscErrorCode Parameter_settings(AppCtx *user) 362d71ae5a4SJacob Faibussowitsch { 363c4762a1bSJed Brown PetscBool flg; 364c4762a1bSJed Brown 365c4762a1bSJed Brown PetscFunctionBeginUser; 366c4762a1bSJed Brown /* Set default parameters */ 367c4762a1bSJed Brown user->ws = 1.0; 368c4762a1bSJed Brown user->H = 5.0; 369c4762a1bSJed Brown user->Pmax = 2.1; 370c4762a1bSJed Brown user->PM_min = 1.0; 371c4762a1bSJed Brown user->lambda = 0.1; 372c4762a1bSJed Brown user->q = 1.0; 373c4762a1bSJed Brown user->mux = PetscAsinScalar(user->PM_min / user->Pmax); 374c4762a1bSJed Brown user->sigmax = 0.1; 375c4762a1bSJed Brown user->sigmay = 0.1; 376c4762a1bSJed Brown user->rho = 0.0; 377c4762a1bSJed Brown user->xmin = -PETSC_PI; 378c4762a1bSJed Brown user->xmax = PETSC_PI; 379c4762a1bSJed Brown user->bx = DM_BOUNDARY_PERIODIC; 380c4762a1bSJed Brown user->by = DM_BOUNDARY_MIRROR; 381c4762a1bSJed Brown user->nonoiseinitial = PETSC_FALSE; 382c4762a1bSJed Brown 383c4762a1bSJed Brown /* 384c4762a1bSJed Brown ymin of -3 seems to let the unstable solution move up and leave a zero in its wake 385c4762a1bSJed Brown with an ymin of -1 the wake is never exactly zero 386c4762a1bSJed Brown */ 387c4762a1bSJed Brown user->ymin = -3.0; 388c4762a1bSJed Brown user->ymax = 10.0; 389c4762a1bSJed Brown user->st_width = 1; 390c4762a1bSJed Brown 3919566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-ws", &user->ws, &flg)); 3929566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-Inertia", &user->H, &flg)); 3939566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-Pmax", &user->Pmax, &flg)); 3949566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-PM_min", &user->PM_min, &flg)); 3959566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-lambda", &user->lambda, &flg)); 3969566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-q", &user->q, &flg)); 3979566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-mux", &user->mux, &flg)); 3989566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-sigmax", &user->sigmax, &flg)); 3999566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-muy", &user->muy, &flg)); 400ad540459SPierre Jolivet if (flg == 0) user->muy = user->ws; 4019566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-sigmay", &user->sigmay, &flg)); 4029566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-rho", &user->rho, &flg)); 4039566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-xmin", &user->xmin, &flg)); 4049566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-xmax", &user->xmax, &flg)); 4059566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-ymin", &user->ymin, &flg)); 4069566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-ymax", &user->ymax, &flg)); 4079566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-stencil_width", &user->st_width, &flg)); 4089566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetEnum(NULL, NULL, "-bx", DMBoundaryTypes, (PetscEnum *)&user->bx, &flg)); 4099566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetEnum(NULL, NULL, "-by", DMBoundaryTypes, (PetscEnum *)&user->by, &flg)); 4109566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-nonoiseinitial", &user->nonoiseinitial, &flg)); 4113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 412c4762a1bSJed Brown } 413c4762a1bSJed Brown 414c4762a1bSJed Brown /*TEST 415c4762a1bSJed Brown 416c4762a1bSJed Brown build: 417c4762a1bSJed Brown requires: !complex !single 418c4762a1bSJed Brown 419c4762a1bSJed Brown test: 420c4762a1bSJed Brown args: -ts_max_steps 2 421c4762a1bSJed Brown localrunfiles: petscopt_ex7 422c4762a1bSJed Brown 423c4762a1bSJed Brown test: 424c4762a1bSJed Brown suffix: 2 425c4762a1bSJed Brown args: -ts_max_steps 2 -snes_mf_operator 426c4762a1bSJed Brown output_file: output/ex7_1.out 427c4762a1bSJed Brown localrunfiles: petscopt_ex7 428c4762a1bSJed Brown timeoutfactor: 2 429c4762a1bSJed Brown 430c4762a1bSJed Brown test: 431c4762a1bSJed Brown suffix: 3 432c4762a1bSJed Brown args: -ts_max_steps 2 -snes_mf -pc_type none 433c4762a1bSJed Brown output_file: output/ex7_1.out 434c4762a1bSJed Brown localrunfiles: petscopt_ex7 435c4762a1bSJed Brown timeoutfactor: 2 436c4762a1bSJed Brown 437c4762a1bSJed Brown TEST*/ 438