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 51*9371c9d4SSatish Balay int main(int argc, char **argv) { 52c4762a1bSJed Brown Vec x; /* Solution vector */ 53c4762a1bSJed Brown TS ts; /* Time-stepping context */ 54c4762a1bSJed Brown AppCtx user; /* Application context */ 55c4762a1bSJed Brown Mat J; 56c4762a1bSJed Brown PetscViewer viewer; 57c4762a1bSJed Brown 58327415f7SBarry Smith PetscFunctionBeginUser; 599566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, "petscopt_ex6", help)); 60c4762a1bSJed Brown /* Get physics and time parameters */ 619566063dSJacob Faibussowitsch PetscCall(Parameter_settings(&user)); 62c4762a1bSJed Brown /* Create a 2D DA with dof = 1 */ 639566063dSJacob 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)); 649566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(user.da)); 659566063dSJacob Faibussowitsch PetscCall(DMSetUp(user.da)); 66c4762a1bSJed Brown /* Set x and y coordinates */ 679566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(user.da, user.xmin, user.xmax, user.ymin, user.ymax, 0.0, 1.0)); 68c4762a1bSJed Brown 69c4762a1bSJed Brown /* Get global vector x from DM */ 709566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(user.da, &x)); 71c4762a1bSJed Brown 729566063dSJacob Faibussowitsch PetscCall(ini_bou(x, &user)); 739566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "ini_x", FILE_MODE_WRITE, &viewer)); 749566063dSJacob Faibussowitsch PetscCall(VecView(x, viewer)); 759566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 76c4762a1bSJed Brown 77c4762a1bSJed Brown /* Get Jacobian matrix structure from the da */ 789566063dSJacob Faibussowitsch PetscCall(DMSetMatType(user.da, MATAIJ)); 799566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(user.da, &J)); 80c4762a1bSJed Brown 819566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 829566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR)); 839566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts, NULL, IFunction, &user)); 849566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts, J, J, IJacobian, &user)); 859566063dSJacob Faibussowitsch PetscCall(TSSetApplicationContext(ts, &user)); 869566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, user.tmax)); 879566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 889566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts, user.t0)); 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(MatDestroy(&J)); 1009566063dSJacob Faibussowitsch PetscCall(DMDestroy(&user.da)); 1019566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 1029566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 103b122ec5aSJacob Faibussowitsch return 0; 104c4762a1bSJed Brown } 105c4762a1bSJed Brown 106*9371c9d4SSatish Balay PetscErrorCode PostStep(TS ts) { 107c4762a1bSJed Brown Vec X; 108c4762a1bSJed Brown AppCtx *user; 109c4762a1bSJed Brown PetscScalar sum; 110c4762a1bSJed Brown PetscReal t; 111c4762a1bSJed Brown 112c4762a1bSJed Brown PetscFunctionBegin; 1139566063dSJacob Faibussowitsch PetscCall(TSGetApplicationContext(ts, &user)); 1149566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &t)); 1159566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ts, &X)); 1169566063dSJacob Faibussowitsch PetscCall(VecSum(X, &sum)); 11763a3b9bcSJacob 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))); 118c4762a1bSJed Brown PetscFunctionReturn(0); 119c4762a1bSJed Brown } 120c4762a1bSJed Brown 121*9371c9d4SSatish Balay PetscErrorCode ini_bou(Vec X, AppCtx *user) { 122c4762a1bSJed Brown DM cda; 123c4762a1bSJed Brown DMDACoor2d **coors; 124c4762a1bSJed Brown PetscScalar **p; 125c4762a1bSJed Brown Vec gc; 126c4762a1bSJed Brown PetscInt i, j; 127c4762a1bSJed Brown PetscInt xs, ys, xm, ym, M, N; 128c4762a1bSJed Brown PetscScalar xi, yi; 129c4762a1bSJed Brown PetscScalar sigmax = user->sigmax, sigmay = user->sigmay; 130c4762a1bSJed Brown PetscScalar rho = user->rho; 131c4762a1bSJed Brown PetscScalar mux = user->mux, muy = user->muy; 132c4762a1bSJed Brown PetscMPIInt rank; 133c4762a1bSJed Brown 134c4762a1bSJed Brown PetscFunctionBeginUser; 1359566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 1369566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(user->da, NULL, &M, &N, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 137*9371c9d4SSatish Balay user->dx = (user->xmax - user->xmin) / (M - 1); 138*9371c9d4SSatish Balay user->dy = (user->ymax - user->ymin) / (N - 1); 1399566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(user->da, &cda)); 1409566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(user->da, &gc)); 1419566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(cda, gc, &coors)); 1429566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(user->da, X, &p)); 1439566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(cda, &xs, &ys, 0, &xm, &ym, 0)); 144c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 145c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 146*9371c9d4SSatish Balay xi = coors[j][i].x; 147*9371c9d4SSatish Balay yi = coors[j][i].y; 148c4762a1bSJed Brown if (i == 0 || j == 0 || i == M - 1 || j == N - 1) p[j][i] = 0.0; 149*9371c9d4SSatish Balay else 150*9371c9d4SSatish 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))); 151c4762a1bSJed Brown } 152c4762a1bSJed Brown } 153c4762a1bSJed Brown /* p[N/2+N%2][M/2+M%2] = 1/(user->dx*user->dy); */ 154c4762a1bSJed Brown 1559566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(cda, gc, &coors)); 1569566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(user->da, X, &p)); 157c4762a1bSJed Brown PetscFunctionReturn(0); 158c4762a1bSJed Brown } 159c4762a1bSJed Brown 160c4762a1bSJed Brown /* First advection term */ 161*9371c9d4SSatish Balay PetscErrorCode adv1(PetscScalar **p, PetscScalar y, PetscInt i, PetscInt j, PetscInt M, PetscScalar *p1, AppCtx *user) { 162c4762a1bSJed Brown PetscScalar f; 163c4762a1bSJed Brown /* PetscScalar v1,v2,v3,v4,v5,s1,s2,s3; */ 164c4762a1bSJed Brown PetscFunctionBegin; 165c4762a1bSJed Brown /* if (i > 2 && i < M-2) { 166c4762a1bSJed Brown v1 = (y-user->ws)*(p[j][i-2] - p[j][i-3])/user->dx; 167c4762a1bSJed Brown v2 = (y-user->ws)*(p[j][i-1] - p[j][i-2])/user->dx; 168c4762a1bSJed Brown v3 = (y-user->ws)*(p[j][i] - p[j][i-1])/user->dx; 169c4762a1bSJed Brown v4 = (y-user->ws)*(p[j][i+1] - p[j][i])/user->dx; 170c4762a1bSJed Brown v5 = (y-user->ws)*(p[j][i+1] - p[j][i+2])/user->dx; 171c4762a1bSJed Brown 172c4762a1bSJed Brown s1 = v1/3.0 - (7.0/6.0)*v2 + (11.0/6.0)*v3; 173c4762a1bSJed Brown s2 =-v2/6.0 + (5.0/6.0)*v3 + (1.0/3.0)*v4; 174c4762a1bSJed Brown s3 = v3/3.0 + (5.0/6.0)*v4 - (1.0/6.0)*v5; 175c4762a1bSJed Brown 176c4762a1bSJed Brown *p1 = 0.1*s1 + 0.6*s2 + 0.3*s3; 177c4762a1bSJed Brown } else *p1 = 0.0; */ 178c4762a1bSJed Brown f = (y - user->ws); 179c4762a1bSJed Brown *p1 = f * (p[j][i + 1] - p[j][i - 1]) / (2 * user->dx); 180c4762a1bSJed Brown PetscFunctionReturn(0); 181c4762a1bSJed Brown } 182c4762a1bSJed Brown 183c4762a1bSJed Brown /* Second advection term */ 184*9371c9d4SSatish Balay PetscErrorCode adv2(PetscScalar **p, PetscScalar x, PetscInt i, PetscInt j, PetscInt N, PetscScalar *p2, AppCtx *user) { 185c4762a1bSJed Brown PetscScalar f; 186c4762a1bSJed Brown /* PetscScalar v1,v2,v3,v4,v5,s1,s2,s3; */ 187c4762a1bSJed Brown PetscFunctionBegin; 188c4762a1bSJed Brown /* if (j > 2 && j < N-2) { 189c4762a1bSJed Brown v1 = (user->ws/(2*user->H))*(user->PM_min - user->Pmax*sin(x))*(p[j-2][i] - p[j-3][i])/user->dy; 190c4762a1bSJed Brown v2 = (user->ws/(2*user->H))*(user->PM_min - user->Pmax*sin(x))*(p[j-1][i] - p[j-2][i])/user->dy; 191c4762a1bSJed Brown v3 = (user->ws/(2*user->H))*(user->PM_min - user->Pmax*sin(x))*(p[j][i] - p[j-1][i])/user->dy; 192c4762a1bSJed Brown v4 = (user->ws/(2*user->H))*(user->PM_min - user->Pmax*sin(x))*(p[j+1][i] - p[j][i])/user->dy; 193c4762a1bSJed Brown v5 = (user->ws/(2*user->H))*(user->PM_min - user->Pmax*sin(x))*(p[j+2][i] - p[j+1][i])/user->dy; 194c4762a1bSJed Brown 195c4762a1bSJed Brown s1 = v1/3.0 - (7.0/6.0)*v2 + (11.0/6.0)*v3; 196c4762a1bSJed Brown s2 =-v2/6.0 + (5.0/6.0)*v3 + (1.0/3.0)*v4; 197c4762a1bSJed Brown s3 = v3/3.0 + (5.0/6.0)*v4 - (1.0/6.0)*v5; 198c4762a1bSJed Brown 199c4762a1bSJed Brown *p2 = 0.1*s1 + 0.6*s2 + 0.3*s3; 200c4762a1bSJed Brown } else *p2 = 0.0; */ 201c4762a1bSJed Brown f = (user->ws / (2 * user->H)) * (user->PM_min - user->Pmax * PetscSinScalar(x)); 202c4762a1bSJed Brown *p2 = f * (p[j + 1][i] - p[j - 1][i]) / (2 * user->dy); 203c4762a1bSJed Brown PetscFunctionReturn(0); 204c4762a1bSJed Brown } 205c4762a1bSJed Brown 206c4762a1bSJed Brown /* Diffusion term */ 207*9371c9d4SSatish Balay PetscErrorCode diffuse(PetscScalar **p, PetscInt i, PetscInt j, PetscReal t, PetscScalar *p_diff, AppCtx *user) { 208c4762a1bSJed Brown PetscFunctionBeginUser; 209c4762a1bSJed Brown 210c4762a1bSJed Brown *p_diff = user->disper_coe * ((p[j - 1][i] - 2 * p[j][i] + p[j + 1][i]) / (user->dy * user->dy)); 211c4762a1bSJed Brown PetscFunctionReturn(0); 212c4762a1bSJed Brown } 213c4762a1bSJed Brown 214*9371c9d4SSatish Balay PetscErrorCode BoundaryConditions(PetscScalar **p, DMDACoor2d **coors, PetscInt i, PetscInt j, PetscInt M, PetscInt N, PetscScalar **f, AppCtx *user) { 215c4762a1bSJed Brown PetscScalar fwc, fthetac; 216c4762a1bSJed Brown PetscScalar w = coors[j][i].y, theta = coors[j][i].x; 217c4762a1bSJed Brown 218c4762a1bSJed Brown PetscFunctionBeginUser; 219c4762a1bSJed Brown if (user->bc == 0) { /* Natural boundary condition */ 220c4762a1bSJed Brown f[j][i] = p[j][i]; 221c4762a1bSJed Brown } else { /* Steady state boundary condition */ 222c4762a1bSJed Brown fthetac = user->ws / (2 * user->H) * (user->PM_min - user->Pmax * PetscSinScalar(theta)); 223c4762a1bSJed Brown fwc = (w * w / 2.0 - user->ws * w); 224c4762a1bSJed Brown if (i == 0 && j == 0) { /* left bottom corner */ 225c4762a1bSJed 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; 226c4762a1bSJed Brown } else if (i == 0 && j == N - 1) { /* right bottom corner */ 227c4762a1bSJed 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; 228c4762a1bSJed Brown } else if (i == M - 1 && j == 0) { /* left top corner */ 229c4762a1bSJed 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; 230c4762a1bSJed Brown } else if (i == M - 1 && j == N - 1) { /* right top corner */ 231c4762a1bSJed 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; 232c4762a1bSJed Brown } else if (i == 0) { /* Bottom edge */ 233c4762a1bSJed 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); 234c4762a1bSJed Brown } else if (i == M - 1) { /* Top edge */ 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 - 1][i]) / (2 * user->dy); 236c4762a1bSJed Brown } else if (j == 0) { /* Left edge */ 237c4762a1bSJed 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); 238c4762a1bSJed Brown } else if (j == N - 1) { /* Right edge */ 239c4762a1bSJed 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); 240c4762a1bSJed Brown } 241c4762a1bSJed Brown } 242c4762a1bSJed Brown PetscFunctionReturn(0); 243c4762a1bSJed Brown } 244c4762a1bSJed Brown 245*9371c9d4SSatish Balay PetscErrorCode IFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx) { 246c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx; 247c4762a1bSJed Brown DM cda; 248c4762a1bSJed Brown DMDACoor2d **coors; 249c4762a1bSJed Brown PetscScalar **p, **f, **pdot; 250c4762a1bSJed Brown PetscInt i, j; 251c4762a1bSJed Brown PetscInt xs, ys, xm, ym, M, N; 252c4762a1bSJed Brown Vec localX, gc, localXdot; 253c4762a1bSJed Brown PetscScalar p_adv1, p_adv2, p_diff; 254c4762a1bSJed Brown 255c4762a1bSJed Brown PetscFunctionBeginUser; 2569566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(user->da, NULL, &M, &N, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 2579566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(user->da, &cda)); 2589566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(cda, &xs, &ys, 0, &xm, &ym, 0)); 259c4762a1bSJed Brown 2609566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(user->da, &localX)); 2619566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(user->da, &localXdot)); 262c4762a1bSJed Brown 2639566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(user->da, X, INSERT_VALUES, localX)); 2649566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(user->da, X, INSERT_VALUES, localX)); 2659566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(user->da, Xdot, INSERT_VALUES, localXdot)); 2669566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(user->da, Xdot, INSERT_VALUES, localXdot)); 267c4762a1bSJed Brown 2689566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(user->da, &gc)); 269c4762a1bSJed Brown 2709566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(cda, gc, &coors)); 2719566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(user->da, localX, &p)); 2729566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(user->da, localXdot, &pdot)); 2739566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(user->da, F, &f)); 274c4762a1bSJed Brown 275c4762a1bSJed Brown user->disper_coe = PetscPowScalar((user->lambda * user->ws) / (2 * user->H), 2) * user->q * (1.0 - PetscExpScalar(-t / user->lambda)); 276c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 277c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 278c4762a1bSJed Brown if (i == 0 || j == 0 || i == M - 1 || j == N - 1) { 2799566063dSJacob Faibussowitsch PetscCall(BoundaryConditions(p, coors, i, j, M, N, f, user)); 280c4762a1bSJed Brown } else { 2819566063dSJacob Faibussowitsch PetscCall(adv1(p, coors[j][i].y, i, j, M, &p_adv1, user)); 2829566063dSJacob Faibussowitsch PetscCall(adv2(p, coors[j][i].x, i, j, N, &p_adv2, user)); 2839566063dSJacob Faibussowitsch PetscCall(diffuse(p, i, j, t, &p_diff, user)); 284c4762a1bSJed Brown f[j][i] = -p_adv1 - p_adv2 + p_diff - pdot[j][i]; 285c4762a1bSJed Brown } 286c4762a1bSJed Brown } 287c4762a1bSJed Brown } 2889566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(user->da, localX, &p)); 2899566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(user->da, localX, &pdot)); 2909566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(user->da, &localX)); 2919566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(user->da, &localXdot)); 2929566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(user->da, F, &f)); 2939566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(cda, gc, &coors)); 294c4762a1bSJed Brown 295c4762a1bSJed Brown PetscFunctionReturn(0); 296c4762a1bSJed Brown } 297c4762a1bSJed Brown 298*9371c9d4SSatish Balay PetscErrorCode IJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat J, Mat Jpre, void *ctx) { 299c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx; 300c4762a1bSJed Brown DM cda; 301c4762a1bSJed Brown DMDACoor2d **coors; 302c4762a1bSJed Brown PetscInt i, j; 303c4762a1bSJed Brown PetscInt xs, ys, xm, ym, M, N; 304c4762a1bSJed Brown Vec gc; 305c4762a1bSJed Brown PetscScalar val[5], xi, yi; 306c4762a1bSJed Brown MatStencil row, col[5]; 307c4762a1bSJed Brown PetscScalar c1, c3, c5; 308c4762a1bSJed Brown 309c4762a1bSJed Brown PetscFunctionBeginUser; 3109566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(user->da, NULL, &M, &N, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 3119566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(user->da, &cda)); 3129566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(cda, &xs, &ys, 0, &xm, &ym, 0)); 313c4762a1bSJed Brown 3149566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(user->da, &gc)); 3159566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(cda, gc, &coors)); 316c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 317c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 318c4762a1bSJed Brown PetscInt nc = 0; 319*9371c9d4SSatish Balay xi = coors[j][i].x; 320*9371c9d4SSatish Balay yi = coors[j][i].y; 321*9371c9d4SSatish Balay row.i = i; 322*9371c9d4SSatish Balay row.j = j; 323c4762a1bSJed Brown if (i == 0 || j == 0 || i == M - 1 || j == N - 1) { 324c4762a1bSJed Brown if (user->bc == 0) { 325*9371c9d4SSatish Balay col[nc].i = i; 326*9371c9d4SSatish Balay col[nc].j = j; 327*9371c9d4SSatish Balay val[nc++] = 1.0; 328c4762a1bSJed Brown } else { 329c4762a1bSJed Brown PetscScalar fthetac, fwc; 330c4762a1bSJed Brown fthetac = user->ws / (2 * user->H) * (user->PM_min - user->Pmax * PetscSinScalar(xi)); 331c4762a1bSJed Brown fwc = (yi * yi / 2.0 - user->ws * yi); 332c4762a1bSJed Brown if (i == 0 && j == 0) { 333*9371c9d4SSatish Balay col[nc].i = i + 1; 334*9371c9d4SSatish Balay col[nc].j = j; 335*9371c9d4SSatish Balay val[nc++] = fwc / user->dx; 336*9371c9d4SSatish Balay col[nc].i = i; 337*9371c9d4SSatish Balay col[nc].j = j + 1; 338*9371c9d4SSatish Balay val[nc++] = -user->disper_coe / user->dy; 339*9371c9d4SSatish Balay col[nc].i = i; 340*9371c9d4SSatish Balay col[nc].j = j; 341*9371c9d4SSatish Balay val[nc++] = -fwc / user->dx + fthetac + user->disper_coe / user->dy; 342c4762a1bSJed Brown } else if (i == 0 && j == N - 1) { 343*9371c9d4SSatish Balay col[nc].i = i + 1; 344*9371c9d4SSatish Balay col[nc].j = j; 345*9371c9d4SSatish Balay val[nc++] = fwc / user->dx; 346*9371c9d4SSatish Balay col[nc].i = i; 347*9371c9d4SSatish Balay col[nc].j = j - 1; 348*9371c9d4SSatish Balay val[nc++] = user->disper_coe / user->dy; 349*9371c9d4SSatish Balay col[nc].i = i; 350*9371c9d4SSatish Balay col[nc].j = j; 351*9371c9d4SSatish Balay val[nc++] = -fwc / user->dx + fthetac - user->disper_coe / user->dy; 352c4762a1bSJed Brown } else if (i == M - 1 && j == 0) { 353*9371c9d4SSatish Balay col[nc].i = i - 1; 354*9371c9d4SSatish Balay col[nc].j = j; 355*9371c9d4SSatish Balay val[nc++] = -fwc / user->dx; 356*9371c9d4SSatish Balay col[nc].i = i; 357*9371c9d4SSatish Balay col[nc].j = j + 1; 358*9371c9d4SSatish Balay val[nc++] = -user->disper_coe / user->dy; 359*9371c9d4SSatish Balay col[nc].i = i; 360*9371c9d4SSatish Balay col[nc].j = j; 361*9371c9d4SSatish Balay val[nc++] = fwc / user->dx + fthetac + user->disper_coe / user->dy; 362c4762a1bSJed Brown } else if (i == M - 1 && j == N - 1) { 363*9371c9d4SSatish Balay col[nc].i = i - 1; 364*9371c9d4SSatish Balay col[nc].j = j; 365*9371c9d4SSatish Balay val[nc++] = -fwc / user->dx; 366*9371c9d4SSatish Balay col[nc].i = i; 367*9371c9d4SSatish Balay col[nc].j = j - 1; 368*9371c9d4SSatish Balay val[nc++] = user->disper_coe / user->dy; 369*9371c9d4SSatish Balay col[nc].i = i; 370*9371c9d4SSatish Balay col[nc].j = j; 371*9371c9d4SSatish Balay val[nc++] = fwc / user->dx + fthetac - user->disper_coe / user->dy; 372c4762a1bSJed Brown } else if (i == 0) { 373*9371c9d4SSatish Balay col[nc].i = i + 1; 374*9371c9d4SSatish Balay col[nc].j = j; 375*9371c9d4SSatish Balay val[nc++] = fwc / user->dx; 376*9371c9d4SSatish Balay col[nc].i = i; 377*9371c9d4SSatish Balay col[nc].j = j + 1; 378*9371c9d4SSatish Balay val[nc++] = -user->disper_coe / (2 * user->dy); 379*9371c9d4SSatish Balay col[nc].i = i; 380*9371c9d4SSatish Balay col[nc].j = j - 1; 381*9371c9d4SSatish Balay val[nc++] = user->disper_coe / (2 * user->dy); 382*9371c9d4SSatish Balay col[nc].i = i; 383*9371c9d4SSatish Balay col[nc].j = j; 384*9371c9d4SSatish Balay val[nc++] = -fwc / user->dx + fthetac; 385c4762a1bSJed Brown } else if (i == M - 1) { 386*9371c9d4SSatish Balay col[nc].i = i - 1; 387*9371c9d4SSatish Balay col[nc].j = j; 388*9371c9d4SSatish Balay val[nc++] = -fwc / user->dx; 389*9371c9d4SSatish Balay col[nc].i = i; 390*9371c9d4SSatish Balay col[nc].j = j + 1; 391*9371c9d4SSatish Balay val[nc++] = -user->disper_coe / (2 * user->dy); 392*9371c9d4SSatish Balay col[nc].i = i; 393*9371c9d4SSatish Balay col[nc].j = j - 1; 394*9371c9d4SSatish Balay val[nc++] = user->disper_coe / (2 * user->dy); 395*9371c9d4SSatish Balay col[nc].i = i; 396*9371c9d4SSatish Balay col[nc].j = j; 397*9371c9d4SSatish Balay val[nc++] = fwc / user->dx + fthetac; 398c4762a1bSJed Brown } else if (j == 0) { 399*9371c9d4SSatish Balay col[nc].i = i + 1; 400*9371c9d4SSatish Balay col[nc].j = j; 401*9371c9d4SSatish Balay val[nc++] = fwc / (2 * user->dx); 402*9371c9d4SSatish Balay col[nc].i = i - 1; 403*9371c9d4SSatish Balay col[nc].j = j; 404*9371c9d4SSatish Balay val[nc++] = -fwc / (2 * user->dx); 405*9371c9d4SSatish Balay col[nc].i = i; 406*9371c9d4SSatish Balay col[nc].j = j + 1; 407*9371c9d4SSatish Balay val[nc++] = -user->disper_coe / user->dy; 408*9371c9d4SSatish Balay col[nc].i = i; 409*9371c9d4SSatish Balay col[nc].j = j; 410*9371c9d4SSatish Balay val[nc++] = user->disper_coe / user->dy + fthetac; 411c4762a1bSJed Brown } else if (j == N - 1) { 412*9371c9d4SSatish Balay col[nc].i = i + 1; 413*9371c9d4SSatish Balay col[nc].j = j; 414*9371c9d4SSatish Balay val[nc++] = fwc / (2 * user->dx); 415*9371c9d4SSatish Balay col[nc].i = i - 1; 416*9371c9d4SSatish Balay col[nc].j = j; 417*9371c9d4SSatish Balay val[nc++] = -fwc / (2 * user->dx); 418*9371c9d4SSatish Balay col[nc].i = i; 419*9371c9d4SSatish Balay col[nc].j = j - 1; 420*9371c9d4SSatish Balay val[nc++] = user->disper_coe / user->dy; 421*9371c9d4SSatish Balay col[nc].i = i; 422*9371c9d4SSatish Balay col[nc].j = j; 423*9371c9d4SSatish Balay val[nc++] = -user->disper_coe / user->dy + fthetac; 424c4762a1bSJed Brown } 425c4762a1bSJed Brown } 426c4762a1bSJed Brown } else { 427c4762a1bSJed Brown c1 = (yi - user->ws) / (2 * user->dx); 428c4762a1bSJed Brown c3 = (user->ws / (2.0 * user->H)) * (user->PM_min - user->Pmax * PetscSinScalar(xi)) / (2 * user->dy); 429c4762a1bSJed Brown c5 = (PetscPowScalar((user->lambda * user->ws) / (2 * user->H), 2) * user->q * (1.0 - PetscExpScalar(-t / user->lambda))) / (user->dy * user->dy); 430*9371c9d4SSatish Balay col[nc].i = i - 1; 431*9371c9d4SSatish Balay col[nc].j = j; 432*9371c9d4SSatish Balay val[nc++] = c1; 433*9371c9d4SSatish Balay col[nc].i = i + 1; 434*9371c9d4SSatish Balay col[nc].j = j; 435*9371c9d4SSatish Balay val[nc++] = -c1; 436*9371c9d4SSatish Balay col[nc].i = i; 437*9371c9d4SSatish Balay col[nc].j = j - 1; 438*9371c9d4SSatish Balay val[nc++] = c3 + c5; 439*9371c9d4SSatish Balay col[nc].i = i; 440*9371c9d4SSatish Balay col[nc].j = j + 1; 441*9371c9d4SSatish Balay val[nc++] = -c3 + c5; 442*9371c9d4SSatish Balay col[nc].i = i; 443*9371c9d4SSatish Balay col[nc].j = j; 444*9371c9d4SSatish Balay val[nc++] = -2 * c5 - a; 445c4762a1bSJed Brown } 4469566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(Jpre, 1, &row, nc, col, val, INSERT_VALUES)); 447c4762a1bSJed Brown } 448c4762a1bSJed Brown } 4499566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(cda, gc, &coors)); 450c4762a1bSJed Brown 4519566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Jpre, MAT_FINAL_ASSEMBLY)); 4529566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Jpre, MAT_FINAL_ASSEMBLY)); 453c4762a1bSJed Brown if (J != Jpre) { 4549566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 4559566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 456c4762a1bSJed Brown } 457c4762a1bSJed Brown PetscFunctionReturn(0); 458c4762a1bSJed Brown } 459c4762a1bSJed Brown 460*9371c9d4SSatish Balay PetscErrorCode Parameter_settings(AppCtx *user) { 461c4762a1bSJed Brown PetscBool flg; 462c4762a1bSJed Brown 463c4762a1bSJed Brown PetscFunctionBeginUser; 464c4762a1bSJed Brown 465c4762a1bSJed Brown /* Set default parameters */ 466c4762a1bSJed Brown user->ws = 1.0; 467*9371c9d4SSatish Balay user->H = 5.0; 468*9371c9d4SSatish Balay user->Pmax = 2.1; 469*9371c9d4SSatish Balay user->PM_min = 1.0; 470*9371c9d4SSatish Balay user->lambda = 0.1; 471*9371c9d4SSatish Balay user->q = 1.0; 472*9371c9d4SSatish Balay user->mux = PetscAsinScalar(user->PM_min / user->Pmax); 473c4762a1bSJed Brown user->sigmax = 0.1; 474*9371c9d4SSatish Balay user->sigmay = 0.1; 475*9371c9d4SSatish Balay user->rho = 0.0; 476*9371c9d4SSatish Balay user->t0 = 0.0; 477*9371c9d4SSatish Balay user->tmax = 2.0; 478*9371c9d4SSatish Balay user->xmin = -1.0; 479*9371c9d4SSatish Balay user->xmax = 10.0; 480*9371c9d4SSatish Balay user->ymin = -1.0; 481*9371c9d4SSatish Balay user->ymax = 10.0; 482c4762a1bSJed Brown user->bc = 0; 483c4762a1bSJed Brown 4849566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-ws", &user->ws, &flg)); 4859566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-Inertia", &user->H, &flg)); 4869566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-Pmax", &user->Pmax, &flg)); 4879566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-PM_min", &user->PM_min, &flg)); 4889566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-lambda", &user->lambda, &flg)); 4899566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-q", &user->q, &flg)); 4909566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-mux", &user->mux, &flg)); 4919566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-sigmax", &user->sigmax, &flg)); 4929566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-muy", &user->muy, &flg)); 4939566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-sigmay", &user->sigmay, &flg)); 4949566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-rho", &user->rho, &flg)); 4959566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-t0", &user->t0, &flg)); 4969566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-tmax", &user->tmax, &flg)); 4979566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-xmin", &user->xmin, &flg)); 4989566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-xmax", &user->xmax, &flg)); 4999566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-ymin", &user->ymin, &flg)); 5009566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-ymax", &user->ymax, &flg)); 5019566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-bc_type", &user->bc, &flg)); 502c4762a1bSJed Brown user->muy = user->ws; 503c4762a1bSJed Brown PetscFunctionReturn(0); 504c4762a1bSJed Brown } 505c4762a1bSJed Brown 506c4762a1bSJed Brown /*TEST 507c4762a1bSJed Brown 508c4762a1bSJed Brown build: 509c4762a1bSJed Brown requires: !complex 510c4762a1bSJed Brown 511c4762a1bSJed Brown test: 512c4762a1bSJed Brown args: -nox -ts_max_steps 2 513c4762a1bSJed Brown localrunfiles: petscopt_ex6 514c4762a1bSJed Brown timeoutfactor: 4 515c4762a1bSJed Brown 516c4762a1bSJed Brown TEST*/ 517