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 x_t = (y - ws) y_t = (ws/2H)*(Pm - Pmax*sin(x)) 6c4762a1bSJed Brown 7c4762a1bSJed Brown Boundary conditions: -bc_type 0 => Zero dirichlet boundary 8c4762a1bSJed Brown -bc_type 1 => Steady state boundary condition 9c4762a1bSJed Brown Steady state boundary condition found by setting p_t = 0 10c4762a1bSJed Brown */ 11c4762a1bSJed Brown 12c4762a1bSJed Brown #include <petscdm.h> 13c4762a1bSJed Brown #include <petscdmda.h> 14c4762a1bSJed Brown #include <petscts.h> 15c4762a1bSJed Brown 16c4762a1bSJed Brown /* 17c4762a1bSJed Brown User-defined data structures and routines 18c4762a1bSJed Brown */ 19c4762a1bSJed Brown typedef struct { 20c4762a1bSJed Brown PetscScalar ws; /* Synchronous speed */ 21c4762a1bSJed Brown PetscScalar H; /* Inertia constant */ 22c4762a1bSJed Brown PetscScalar D; /* Damping constant */ 23c4762a1bSJed Brown PetscScalar Pmax; /* Maximum power output of generator */ 24c4762a1bSJed Brown PetscScalar PM_min; /* Mean mechanical power input */ 25c4762a1bSJed Brown PetscScalar lambda; /* correlation time */ 26c4762a1bSJed Brown PetscScalar q; /* noise strength */ 27c4762a1bSJed Brown PetscScalar mux; /* Initial average angle */ 28c4762a1bSJed Brown PetscScalar sigmax; /* Standard deviation of initial angle */ 29c4762a1bSJed Brown PetscScalar muy; /* Average speed */ 30c4762a1bSJed Brown PetscScalar sigmay; /* standard deviation of initial speed */ 31c4762a1bSJed Brown PetscScalar rho; /* Cross-correlation coefficient */ 32c4762a1bSJed Brown PetscScalar t0; /* Initial time */ 33c4762a1bSJed Brown PetscScalar tmax; /* Final time */ 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 PetscInt bc; /* Boundary conditions */ 41c4762a1bSJed Brown PetscScalar disper_coe; /* Dispersion coefficient */ 42c4762a1bSJed Brown DM da; 43c4762a1bSJed Brown } AppCtx; 44c4762a1bSJed Brown 45c4762a1bSJed Brown PetscErrorCode Parameter_settings(AppCtx *); 46c4762a1bSJed Brown PetscErrorCode ini_bou(Vec, AppCtx *); 47c4762a1bSJed Brown PetscErrorCode IFunction(TS, PetscReal, Vec, Vec, Vec, void *); 48c4762a1bSJed Brown PetscErrorCode IJacobian(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *); 49c4762a1bSJed Brown PetscErrorCode PostStep(TS); 50c4762a1bSJed Brown 51d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 52d71ae5a4SJacob Faibussowitsch { 53c4762a1bSJed Brown Vec x; /* Solution vector */ 54c4762a1bSJed Brown TS ts; /* Time-stepping context */ 55c4762a1bSJed Brown AppCtx user; /* Application context */ 56c4762a1bSJed Brown Mat J; 57c4762a1bSJed Brown PetscViewer viewer; 58c4762a1bSJed Brown 59327415f7SBarry Smith PetscFunctionBeginUser; 609566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, "petscopt_ex6", help)); 61c4762a1bSJed Brown /* Get physics and time parameters */ 629566063dSJacob Faibussowitsch PetscCall(Parameter_settings(&user)); 63c4762a1bSJed Brown /* Create a 2D DA with dof = 1 */ 649566063dSJacob Faibussowitsch 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)); 659566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(user.da)); 669566063dSJacob Faibussowitsch PetscCall(DMSetUp(user.da)); 67c4762a1bSJed Brown /* Set x and y coordinates */ 689566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(user.da, user.xmin, user.xmax, user.ymin, user.ymax, 0.0, 1.0)); 69c4762a1bSJed Brown 70c4762a1bSJed Brown /* Get global vector x from DM */ 719566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(user.da, &x)); 72c4762a1bSJed Brown 739566063dSJacob Faibussowitsch PetscCall(ini_bou(x, &user)); 749566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "ini_x", FILE_MODE_WRITE, &viewer)); 759566063dSJacob Faibussowitsch PetscCall(VecView(x, viewer)); 769566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 77c4762a1bSJed Brown 78c4762a1bSJed Brown /* Get Jacobian matrix structure from the da */ 799566063dSJacob Faibussowitsch PetscCall(DMSetMatType(user.da, MATAIJ)); 809566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(user.da, &J)); 81c4762a1bSJed Brown 829566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 839566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR)); 849566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts, NULL, IFunction, &user)); 859566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts, J, J, IJacobian, &user)); 869566063dSJacob Faibussowitsch PetscCall(TSSetApplicationContext(ts, &user)); 879566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, user.tmax)); 889566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 899566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts, user.t0)); 909566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, .005)); 919566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 929566063dSJacob Faibussowitsch PetscCall(TSSetPostStep(ts, PostStep)); 939566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, x)); 94c4762a1bSJed Brown 959566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "fin_x", FILE_MODE_WRITE, &viewer)); 969566063dSJacob Faibussowitsch PetscCall(VecView(x, viewer)); 979566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 98c4762a1bSJed Brown 999566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1009566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 1019566063dSJacob Faibussowitsch PetscCall(DMDestroy(&user.da)); 1029566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 1039566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 104b122ec5aSJacob Faibussowitsch return 0; 105c4762a1bSJed Brown } 106c4762a1bSJed Brown 107d71ae5a4SJacob Faibussowitsch PetscErrorCode PostStep(TS ts) 108d71ae5a4SJacob Faibussowitsch { 109c4762a1bSJed Brown Vec X; 110c4762a1bSJed Brown AppCtx *user; 111c4762a1bSJed Brown PetscScalar sum; 112c4762a1bSJed Brown PetscReal t; 113c4762a1bSJed Brown 114c4762a1bSJed Brown PetscFunctionBegin; 1159566063dSJacob Faibussowitsch PetscCall(TSGetApplicationContext(ts, &user)); 1169566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &t)); 1179566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ts, &X)); 1189566063dSJacob Faibussowitsch PetscCall(VecSum(X, &sum)); 11963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "sum(p)*dw*dtheta at t = %3.2f = %3.6f\n", (double)t, (double)(sum * user->dx * user->dy))); 1203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 121c4762a1bSJed Brown } 122c4762a1bSJed Brown 123d71ae5a4SJacob Faibussowitsch PetscErrorCode ini_bou(Vec X, AppCtx *user) 124d71ae5a4SJacob Faibussowitsch { 125c4762a1bSJed Brown DM cda; 126c4762a1bSJed Brown DMDACoor2d **coors; 127c4762a1bSJed Brown PetscScalar **p; 128c4762a1bSJed Brown Vec gc; 129c4762a1bSJed Brown PetscInt i, j; 130c4762a1bSJed Brown PetscInt xs, ys, xm, ym, M, N; 131c4762a1bSJed Brown PetscScalar xi, yi; 132c4762a1bSJed Brown PetscScalar sigmax = user->sigmax, sigmay = user->sigmay; 133c4762a1bSJed Brown PetscScalar rho = user->rho; 134c4762a1bSJed Brown PetscScalar mux = user->mux, muy = user->muy; 135c4762a1bSJed Brown PetscMPIInt rank; 136c4762a1bSJed Brown 137c4762a1bSJed Brown PetscFunctionBeginUser; 1389566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 1399566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(user->da, NULL, &M, &N, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 1409371c9d4SSatish Balay user->dx = (user->xmax - user->xmin) / (M - 1); 1419371c9d4SSatish Balay user->dy = (user->ymax - user->ymin) / (N - 1); 1429566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(user->da, &cda)); 1439566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(user->da, &gc)); 1449566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(cda, gc, &coors)); 1459566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(user->da, X, &p)); 1469566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(cda, &xs, &ys, 0, &xm, &ym, 0)); 147c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 148c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 1499371c9d4SSatish Balay xi = coors[j][i].x; 1509371c9d4SSatish Balay yi = coors[j][i].y; 151c4762a1bSJed Brown if (i == 0 || j == 0 || i == M - 1 || j == N - 1) p[j][i] = 0.0; 1529371c9d4SSatish Balay else 1539371c9d4SSatish Balay 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))); 154c4762a1bSJed Brown } 155c4762a1bSJed Brown } 156c4762a1bSJed Brown /* p[N/2+N%2][M/2+M%2] = 1/(user->dx*user->dy); */ 157c4762a1bSJed Brown 1589566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(cda, gc, &coors)); 1599566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(user->da, X, &p)); 1603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 161c4762a1bSJed Brown } 162c4762a1bSJed Brown 163c4762a1bSJed Brown /* First advection term */ 164d71ae5a4SJacob Faibussowitsch PetscErrorCode adv1(PetscScalar **p, PetscScalar y, PetscInt i, PetscInt j, PetscInt M, PetscScalar *p1, AppCtx *user) 165d71ae5a4SJacob Faibussowitsch { 166*4d86920dSPierre Jolivet PetscScalar f; /* v1,v2,v3,v4,v5,s1,s2,s3; */ 167*4d86920dSPierre Jolivet 168c4762a1bSJed Brown PetscFunctionBegin; 169c4762a1bSJed Brown /* if (i > 2 && i < M-2) { 170c4762a1bSJed Brown v1 = (y-user->ws)*(p[j][i-2] - p[j][i-3])/user->dx; 171c4762a1bSJed Brown v2 = (y-user->ws)*(p[j][i-1] - p[j][i-2])/user->dx; 172c4762a1bSJed Brown v3 = (y-user->ws)*(p[j][i] - p[j][i-1])/user->dx; 173c4762a1bSJed Brown v4 = (y-user->ws)*(p[j][i+1] - p[j][i])/user->dx; 174c4762a1bSJed Brown v5 = (y-user->ws)*(p[j][i+1] - p[j][i+2])/user->dx; 175c4762a1bSJed Brown 176c4762a1bSJed Brown s1 = v1/3.0 - (7.0/6.0)*v2 + (11.0/6.0)*v3; 177c4762a1bSJed Brown s2 =-v2/6.0 + (5.0/6.0)*v3 + (1.0/3.0)*v4; 178c4762a1bSJed Brown s3 = v3/3.0 + (5.0/6.0)*v4 - (1.0/6.0)*v5; 179c4762a1bSJed Brown 180c4762a1bSJed Brown *p1 = 0.1*s1 + 0.6*s2 + 0.3*s3; 181c4762a1bSJed Brown } else *p1 = 0.0; */ 182c4762a1bSJed Brown f = (y - user->ws); 183c4762a1bSJed Brown *p1 = f * (p[j][i + 1] - p[j][i - 1]) / (2 * user->dx); 1843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 185c4762a1bSJed Brown } 186c4762a1bSJed Brown 187c4762a1bSJed Brown /* Second advection term */ 188d71ae5a4SJacob Faibussowitsch PetscErrorCode adv2(PetscScalar **p, PetscScalar x, PetscInt i, PetscInt j, PetscInt N, PetscScalar *p2, AppCtx *user) 189d71ae5a4SJacob Faibussowitsch { 190*4d86920dSPierre Jolivet PetscScalar f; /* v1,v2,v3,v4,v5,s1,s2,s3; */ 191*4d86920dSPierre Jolivet 192c4762a1bSJed Brown PetscFunctionBegin; 193c4762a1bSJed Brown /* if (j > 2 && j < N-2) { 194c4762a1bSJed Brown v1 = (user->ws/(2*user->H))*(user->PM_min - user->Pmax*sin(x))*(p[j-2][i] - p[j-3][i])/user->dy; 195c4762a1bSJed Brown v2 = (user->ws/(2*user->H))*(user->PM_min - user->Pmax*sin(x))*(p[j-1][i] - p[j-2][i])/user->dy; 196c4762a1bSJed Brown v3 = (user->ws/(2*user->H))*(user->PM_min - user->Pmax*sin(x))*(p[j][i] - p[j-1][i])/user->dy; 197c4762a1bSJed Brown v4 = (user->ws/(2*user->H))*(user->PM_min - user->Pmax*sin(x))*(p[j+1][i] - p[j][i])/user->dy; 198c4762a1bSJed Brown v5 = (user->ws/(2*user->H))*(user->PM_min - user->Pmax*sin(x))*(p[j+2][i] - p[j+1][i])/user->dy; 199c4762a1bSJed Brown 200c4762a1bSJed Brown s1 = v1/3.0 - (7.0/6.0)*v2 + (11.0/6.0)*v3; 201c4762a1bSJed Brown s2 =-v2/6.0 + (5.0/6.0)*v3 + (1.0/3.0)*v4; 202c4762a1bSJed Brown s3 = v3/3.0 + (5.0/6.0)*v4 - (1.0/6.0)*v5; 203c4762a1bSJed Brown 204c4762a1bSJed Brown *p2 = 0.1*s1 + 0.6*s2 + 0.3*s3; 205c4762a1bSJed Brown } else *p2 = 0.0; */ 206c4762a1bSJed Brown f = (user->ws / (2 * user->H)) * (user->PM_min - user->Pmax * PetscSinScalar(x)); 207c4762a1bSJed Brown *p2 = f * (p[j + 1][i] - p[j - 1][i]) / (2 * user->dy); 2083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 209c4762a1bSJed Brown } 210c4762a1bSJed Brown 211c4762a1bSJed Brown /* Diffusion term */ 212d71ae5a4SJacob Faibussowitsch PetscErrorCode diffuse(PetscScalar **p, PetscInt i, PetscInt j, PetscReal t, PetscScalar *p_diff, AppCtx *user) 213d71ae5a4SJacob Faibussowitsch { 214c4762a1bSJed Brown PetscFunctionBeginUser; 215c4762a1bSJed Brown *p_diff = user->disper_coe * ((p[j - 1][i] - 2 * p[j][i] + p[j + 1][i]) / (user->dy * user->dy)); 2163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 217c4762a1bSJed Brown } 218c4762a1bSJed Brown 219d71ae5a4SJacob Faibussowitsch PetscErrorCode BoundaryConditions(PetscScalar **p, DMDACoor2d **coors, PetscInt i, PetscInt j, PetscInt M, PetscInt N, PetscScalar **f, AppCtx *user) 220d71ae5a4SJacob Faibussowitsch { 221c4762a1bSJed Brown PetscScalar fwc, fthetac; 222c4762a1bSJed Brown PetscScalar w = coors[j][i].y, theta = coors[j][i].x; 223c4762a1bSJed Brown 224c4762a1bSJed Brown PetscFunctionBeginUser; 225c4762a1bSJed Brown if (user->bc == 0) { /* Natural boundary condition */ 226c4762a1bSJed Brown f[j][i] = p[j][i]; 227c4762a1bSJed Brown } else { /* Steady state boundary condition */ 228c4762a1bSJed Brown fthetac = user->ws / (2 * user->H) * (user->PM_min - user->Pmax * PetscSinScalar(theta)); 229c4762a1bSJed Brown fwc = (w * w / 2.0 - user->ws * w); 230c4762a1bSJed Brown if (i == 0 && j == 0) { /* left bottom corner */ 231c4762a1bSJed Brown 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; 232c4762a1bSJed Brown } else if (i == 0 && j == N - 1) { /* right bottom corner */ 233c4762a1bSJed Brown 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; 234c4762a1bSJed Brown } else if (i == M - 1 && j == 0) { /* left top corner */ 235c4762a1bSJed Brown 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; 236c4762a1bSJed Brown } else if (i == M - 1 && j == N - 1) { /* right top corner */ 237c4762a1bSJed Brown 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; 238c4762a1bSJed Brown } else if (i == 0) { /* Bottom edge */ 239c4762a1bSJed Brown 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); 240c4762a1bSJed Brown } else if (i == M - 1) { /* Top edge */ 241c4762a1bSJed Brown 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); 242c4762a1bSJed Brown } else if (j == 0) { /* Left edge */ 243c4762a1bSJed Brown 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); 244c4762a1bSJed Brown } else if (j == N - 1) { /* Right edge */ 245c4762a1bSJed Brown 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); 246c4762a1bSJed Brown } 247c4762a1bSJed Brown } 2483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 249c4762a1bSJed Brown } 250c4762a1bSJed Brown 251d71ae5a4SJacob Faibussowitsch PetscErrorCode IFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx) 252d71ae5a4SJacob Faibussowitsch { 253c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx; 254c4762a1bSJed Brown DM cda; 255c4762a1bSJed Brown DMDACoor2d **coors; 256c4762a1bSJed Brown PetscScalar **p, **f, **pdot; 257c4762a1bSJed Brown PetscInt i, j; 258c4762a1bSJed Brown PetscInt xs, ys, xm, ym, M, N; 259c4762a1bSJed Brown Vec localX, gc, localXdot; 260c4762a1bSJed Brown PetscScalar p_adv1, p_adv2, p_diff; 261c4762a1bSJed Brown 262c4762a1bSJed Brown PetscFunctionBeginUser; 2639566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(user->da, NULL, &M, &N, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 2649566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(user->da, &cda)); 2659566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(cda, &xs, &ys, 0, &xm, &ym, 0)); 266c4762a1bSJed Brown 2679566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(user->da, &localX)); 2689566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(user->da, &localXdot)); 269c4762a1bSJed Brown 2709566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(user->da, X, INSERT_VALUES, localX)); 2719566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(user->da, X, INSERT_VALUES, localX)); 2729566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(user->da, Xdot, INSERT_VALUES, localXdot)); 2739566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(user->da, Xdot, INSERT_VALUES, localXdot)); 274c4762a1bSJed Brown 2759566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(user->da, &gc)); 276c4762a1bSJed Brown 2779566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(cda, gc, &coors)); 2789566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(user->da, localX, &p)); 2799566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(user->da, localXdot, &pdot)); 2809566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(user->da, F, &f)); 281c4762a1bSJed Brown 282c4762a1bSJed Brown user->disper_coe = PetscPowScalar((user->lambda * user->ws) / (2 * user->H), 2) * user->q * (1.0 - PetscExpScalar(-t / user->lambda)); 283c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 284c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 285c4762a1bSJed Brown if (i == 0 || j == 0 || i == M - 1 || j == N - 1) { 2869566063dSJacob Faibussowitsch PetscCall(BoundaryConditions(p, coors, i, j, M, N, f, user)); 287c4762a1bSJed Brown } else { 2889566063dSJacob Faibussowitsch PetscCall(adv1(p, coors[j][i].y, i, j, M, &p_adv1, user)); 2899566063dSJacob Faibussowitsch PetscCall(adv2(p, coors[j][i].x, i, j, N, &p_adv2, user)); 2909566063dSJacob Faibussowitsch PetscCall(diffuse(p, i, j, t, &p_diff, user)); 291c4762a1bSJed Brown f[j][i] = -p_adv1 - p_adv2 + p_diff - pdot[j][i]; 292c4762a1bSJed Brown } 293c4762a1bSJed Brown } 294c4762a1bSJed Brown } 2959566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(user->da, localX, &p)); 2969566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(user->da, localX, &pdot)); 2979566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(user->da, &localX)); 2989566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(user->da, &localXdot)); 2999566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(user->da, F, &f)); 3009566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(cda, gc, &coors)); 3013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 302c4762a1bSJed Brown } 303c4762a1bSJed Brown 304d71ae5a4SJacob Faibussowitsch PetscErrorCode IJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat J, Mat Jpre, void *ctx) 305d71ae5a4SJacob Faibussowitsch { 306c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx; 307c4762a1bSJed Brown DM cda; 308c4762a1bSJed Brown DMDACoor2d **coors; 309c4762a1bSJed Brown PetscInt i, j; 310c4762a1bSJed Brown PetscInt xs, ys, xm, ym, M, N; 311c4762a1bSJed Brown Vec gc; 312c4762a1bSJed Brown PetscScalar val[5], xi, yi; 313c4762a1bSJed Brown MatStencil row, col[5]; 314c4762a1bSJed Brown PetscScalar c1, c3, c5; 315c4762a1bSJed Brown 316c4762a1bSJed Brown PetscFunctionBeginUser; 3179566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(user->da, NULL, &M, &N, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 3189566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(user->da, &cda)); 3199566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(cda, &xs, &ys, 0, &xm, &ym, 0)); 320c4762a1bSJed Brown 3219566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(user->da, &gc)); 3229566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(cda, gc, &coors)); 323c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 324c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 325c4762a1bSJed Brown PetscInt nc = 0; 3269371c9d4SSatish Balay xi = coors[j][i].x; 3279371c9d4SSatish Balay yi = coors[j][i].y; 3289371c9d4SSatish Balay row.i = i; 3299371c9d4SSatish Balay row.j = j; 330c4762a1bSJed Brown if (i == 0 || j == 0 || i == M - 1 || j == N - 1) { 331c4762a1bSJed Brown if (user->bc == 0) { 3329371c9d4SSatish Balay col[nc].i = i; 3339371c9d4SSatish Balay col[nc].j = j; 3349371c9d4SSatish Balay val[nc++] = 1.0; 335c4762a1bSJed Brown } else { 336c4762a1bSJed Brown PetscScalar fthetac, fwc; 337c4762a1bSJed Brown fthetac = user->ws / (2 * user->H) * (user->PM_min - user->Pmax * PetscSinScalar(xi)); 338c4762a1bSJed Brown fwc = (yi * yi / 2.0 - user->ws * yi); 339c4762a1bSJed Brown if (i == 0 && j == 0) { 3409371c9d4SSatish Balay col[nc].i = i + 1; 3419371c9d4SSatish Balay col[nc].j = j; 3429371c9d4SSatish Balay val[nc++] = fwc / user->dx; 3439371c9d4SSatish Balay col[nc].i = i; 3449371c9d4SSatish Balay col[nc].j = j + 1; 3459371c9d4SSatish Balay val[nc++] = -user->disper_coe / user->dy; 3469371c9d4SSatish Balay col[nc].i = i; 3479371c9d4SSatish Balay col[nc].j = j; 3489371c9d4SSatish Balay val[nc++] = -fwc / user->dx + fthetac + user->disper_coe / user->dy; 349c4762a1bSJed Brown } else if (i == 0 && j == N - 1) { 3509371c9d4SSatish Balay col[nc].i = i + 1; 3519371c9d4SSatish Balay col[nc].j = j; 3529371c9d4SSatish Balay val[nc++] = fwc / user->dx; 3539371c9d4SSatish Balay col[nc].i = i; 3549371c9d4SSatish Balay col[nc].j = j - 1; 3559371c9d4SSatish Balay val[nc++] = user->disper_coe / user->dy; 3569371c9d4SSatish Balay col[nc].i = i; 3579371c9d4SSatish Balay col[nc].j = j; 3589371c9d4SSatish Balay val[nc++] = -fwc / user->dx + fthetac - user->disper_coe / user->dy; 359c4762a1bSJed Brown } else if (i == M - 1 && j == 0) { 3609371c9d4SSatish Balay col[nc].i = i - 1; 3619371c9d4SSatish Balay col[nc].j = j; 3629371c9d4SSatish Balay val[nc++] = -fwc / user->dx; 3639371c9d4SSatish Balay col[nc].i = i; 3649371c9d4SSatish Balay col[nc].j = j + 1; 3659371c9d4SSatish Balay val[nc++] = -user->disper_coe / user->dy; 3669371c9d4SSatish Balay col[nc].i = i; 3679371c9d4SSatish Balay col[nc].j = j; 3689371c9d4SSatish Balay val[nc++] = fwc / user->dx + fthetac + user->disper_coe / user->dy; 369c4762a1bSJed Brown } else if (i == M - 1 && j == N - 1) { 3709371c9d4SSatish Balay col[nc].i = i - 1; 3719371c9d4SSatish Balay col[nc].j = j; 3729371c9d4SSatish Balay val[nc++] = -fwc / user->dx; 3739371c9d4SSatish Balay col[nc].i = i; 3749371c9d4SSatish Balay col[nc].j = j - 1; 3759371c9d4SSatish Balay val[nc++] = user->disper_coe / user->dy; 3769371c9d4SSatish Balay col[nc].i = i; 3779371c9d4SSatish Balay col[nc].j = j; 3789371c9d4SSatish Balay val[nc++] = fwc / user->dx + fthetac - user->disper_coe / user->dy; 379c4762a1bSJed Brown } else if (i == 0) { 3809371c9d4SSatish Balay col[nc].i = i + 1; 3819371c9d4SSatish Balay col[nc].j = j; 3829371c9d4SSatish Balay val[nc++] = fwc / user->dx; 3839371c9d4SSatish Balay col[nc].i = i; 3849371c9d4SSatish Balay col[nc].j = j + 1; 3859371c9d4SSatish Balay val[nc++] = -user->disper_coe / (2 * user->dy); 3869371c9d4SSatish Balay col[nc].i = i; 3879371c9d4SSatish Balay col[nc].j = j - 1; 3889371c9d4SSatish Balay val[nc++] = user->disper_coe / (2 * user->dy); 3899371c9d4SSatish Balay col[nc].i = i; 3909371c9d4SSatish Balay col[nc].j = j; 3919371c9d4SSatish Balay val[nc++] = -fwc / user->dx + fthetac; 392c4762a1bSJed Brown } else if (i == M - 1) { 3939371c9d4SSatish Balay col[nc].i = i - 1; 3949371c9d4SSatish Balay col[nc].j = j; 3959371c9d4SSatish Balay val[nc++] = -fwc / user->dx; 3969371c9d4SSatish Balay col[nc].i = i; 3979371c9d4SSatish Balay col[nc].j = j + 1; 3989371c9d4SSatish Balay val[nc++] = -user->disper_coe / (2 * user->dy); 3999371c9d4SSatish Balay col[nc].i = i; 4009371c9d4SSatish Balay col[nc].j = j - 1; 4019371c9d4SSatish Balay val[nc++] = user->disper_coe / (2 * user->dy); 4029371c9d4SSatish Balay col[nc].i = i; 4039371c9d4SSatish Balay col[nc].j = j; 4049371c9d4SSatish Balay val[nc++] = fwc / user->dx + fthetac; 405c4762a1bSJed Brown } else if (j == 0) { 4069371c9d4SSatish Balay col[nc].i = i + 1; 4079371c9d4SSatish Balay col[nc].j = j; 4089371c9d4SSatish Balay val[nc++] = fwc / (2 * user->dx); 4099371c9d4SSatish Balay col[nc].i = i - 1; 4109371c9d4SSatish Balay col[nc].j = j; 4119371c9d4SSatish Balay val[nc++] = -fwc / (2 * user->dx); 4129371c9d4SSatish Balay col[nc].i = i; 4139371c9d4SSatish Balay col[nc].j = j + 1; 4149371c9d4SSatish Balay val[nc++] = -user->disper_coe / user->dy; 4159371c9d4SSatish Balay col[nc].i = i; 4169371c9d4SSatish Balay col[nc].j = j; 4179371c9d4SSatish Balay val[nc++] = user->disper_coe / user->dy + fthetac; 418c4762a1bSJed Brown } else if (j == N - 1) { 4199371c9d4SSatish Balay col[nc].i = i + 1; 4209371c9d4SSatish Balay col[nc].j = j; 4219371c9d4SSatish Balay val[nc++] = fwc / (2 * user->dx); 4229371c9d4SSatish Balay col[nc].i = i - 1; 4239371c9d4SSatish Balay col[nc].j = j; 4249371c9d4SSatish Balay val[nc++] = -fwc / (2 * user->dx); 4259371c9d4SSatish Balay col[nc].i = i; 4269371c9d4SSatish Balay col[nc].j = j - 1; 4279371c9d4SSatish Balay val[nc++] = user->disper_coe / user->dy; 4289371c9d4SSatish Balay col[nc].i = i; 4299371c9d4SSatish Balay col[nc].j = j; 4309371c9d4SSatish Balay val[nc++] = -user->disper_coe / user->dy + fthetac; 431c4762a1bSJed Brown } 432c4762a1bSJed Brown } 433c4762a1bSJed Brown } else { 434c4762a1bSJed Brown c1 = (yi - user->ws) / (2 * user->dx); 435c4762a1bSJed Brown c3 = (user->ws / (2.0 * user->H)) * (user->PM_min - user->Pmax * PetscSinScalar(xi)) / (2 * user->dy); 436c4762a1bSJed Brown c5 = (PetscPowScalar((user->lambda * user->ws) / (2 * user->H), 2) * user->q * (1.0 - PetscExpScalar(-t / user->lambda))) / (user->dy * user->dy); 4379371c9d4SSatish Balay col[nc].i = i - 1; 4389371c9d4SSatish Balay col[nc].j = j; 4399371c9d4SSatish Balay val[nc++] = c1; 4409371c9d4SSatish Balay col[nc].i = i + 1; 4419371c9d4SSatish Balay col[nc].j = j; 4429371c9d4SSatish Balay val[nc++] = -c1; 4439371c9d4SSatish Balay col[nc].i = i; 4449371c9d4SSatish Balay col[nc].j = j - 1; 4459371c9d4SSatish Balay val[nc++] = c3 + c5; 4469371c9d4SSatish Balay col[nc].i = i; 4479371c9d4SSatish Balay col[nc].j = j + 1; 4489371c9d4SSatish Balay val[nc++] = -c3 + c5; 4499371c9d4SSatish Balay col[nc].i = i; 4509371c9d4SSatish Balay col[nc].j = j; 4519371c9d4SSatish Balay val[nc++] = -2 * c5 - a; 452c4762a1bSJed Brown } 4539566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(Jpre, 1, &row, nc, col, val, INSERT_VALUES)); 454c4762a1bSJed Brown } 455c4762a1bSJed Brown } 4569566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(cda, gc, &coors)); 457c4762a1bSJed Brown 4589566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Jpre, MAT_FINAL_ASSEMBLY)); 4599566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Jpre, MAT_FINAL_ASSEMBLY)); 460c4762a1bSJed Brown if (J != Jpre) { 4619566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 4629566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 463c4762a1bSJed Brown } 4643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 465c4762a1bSJed Brown } 466c4762a1bSJed Brown 467d71ae5a4SJacob Faibussowitsch PetscErrorCode Parameter_settings(AppCtx *user) 468d71ae5a4SJacob Faibussowitsch { 469c4762a1bSJed Brown PetscBool flg; 470c4762a1bSJed Brown 471c4762a1bSJed Brown PetscFunctionBeginUser; 472c4762a1bSJed Brown /* Set default parameters */ 473c4762a1bSJed Brown user->ws = 1.0; 4749371c9d4SSatish Balay user->H = 5.0; 4759371c9d4SSatish Balay user->Pmax = 2.1; 4769371c9d4SSatish Balay user->PM_min = 1.0; 4779371c9d4SSatish Balay user->lambda = 0.1; 4789371c9d4SSatish Balay user->q = 1.0; 4799371c9d4SSatish Balay user->mux = PetscAsinScalar(user->PM_min / user->Pmax); 480c4762a1bSJed Brown user->sigmax = 0.1; 4819371c9d4SSatish Balay user->sigmay = 0.1; 4829371c9d4SSatish Balay user->rho = 0.0; 4839371c9d4SSatish Balay user->t0 = 0.0; 4849371c9d4SSatish Balay user->tmax = 2.0; 4859371c9d4SSatish Balay user->xmin = -1.0; 4869371c9d4SSatish Balay user->xmax = 10.0; 4879371c9d4SSatish Balay user->ymin = -1.0; 4889371c9d4SSatish Balay user->ymax = 10.0; 489c4762a1bSJed Brown user->bc = 0; 490c4762a1bSJed Brown 4919566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-ws", &user->ws, &flg)); 4929566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-Inertia", &user->H, &flg)); 4939566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-Pmax", &user->Pmax, &flg)); 4949566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-PM_min", &user->PM_min, &flg)); 4959566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-lambda", &user->lambda, &flg)); 4969566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-q", &user->q, &flg)); 4979566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-mux", &user->mux, &flg)); 4989566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-sigmax", &user->sigmax, &flg)); 4999566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-muy", &user->muy, &flg)); 5009566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-sigmay", &user->sigmay, &flg)); 5019566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-rho", &user->rho, &flg)); 5029566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-t0", &user->t0, &flg)); 5039566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-tmax", &user->tmax, &flg)); 5049566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-xmin", &user->xmin, &flg)); 5059566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-xmax", &user->xmax, &flg)); 5069566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-ymin", &user->ymin, &flg)); 5079566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-ymax", &user->ymax, &flg)); 5089566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-bc_type", &user->bc, &flg)); 509c4762a1bSJed Brown user->muy = user->ws; 5103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 511c4762a1bSJed Brown } 512c4762a1bSJed Brown 513c4762a1bSJed Brown /*TEST 514c4762a1bSJed Brown 515c4762a1bSJed Brown build: 516c4762a1bSJed Brown requires: !complex 517c4762a1bSJed Brown 518c4762a1bSJed Brown test: 519c4762a1bSJed Brown args: -nox -ts_max_steps 2 520c4762a1bSJed Brown localrunfiles: petscopt_ex6 521c4762a1bSJed Brown timeoutfactor: 4 522c4762a1bSJed Brown 523c4762a1bSJed Brown TEST*/ 524