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: 7c4762a1bSJed Brown Zero dirichlet in y using ghosted values 8c4762a1bSJed Brown Periodic in x 9c4762a1bSJed Brown 10c4762a1bSJed 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. 11c4762a1bSJed Brown x_t = (y - ws) 12c4762a1bSJed Brown y_t = (ws/2H)*(Pm - Pmax*sin(x) - D*(w - ws)) 13c4762a1bSJed Brown 14c4762a1bSJed Brown In this example, we can see the effect of a fault, that zeroes the electrical power output 15c4762a1bSJed Brown Pmax*sin(x), on the PDF. The fault on/off times can be controlled by options -tf and -tcl respectively. 16c4762a1bSJed Brown 17c4762a1bSJed Brown */ 18c4762a1bSJed Brown 19c4762a1bSJed Brown #include <petscdm.h> 20c4762a1bSJed Brown #include <petscdmda.h> 21c4762a1bSJed Brown #include <petscts.h> 22c4762a1bSJed Brown 23c4762a1bSJed Brown /* 24c4762a1bSJed Brown User-defined data structures and routines 25c4762a1bSJed Brown */ 26c4762a1bSJed Brown typedef struct { 27c4762a1bSJed Brown PetscScalar ws; /* Synchronous speed */ 28c4762a1bSJed Brown PetscScalar H; /* Inertia constant */ 29c4762a1bSJed Brown PetscScalar D; /* Damping constant */ 30c4762a1bSJed Brown PetscScalar Pmax, Pmax_s; /* Maximum power output of generator */ 31c4762a1bSJed Brown PetscScalar PM_min; /* Mean mechanical power input */ 32c4762a1bSJed Brown PetscScalar lambda; /* correlation time */ 33c4762a1bSJed Brown PetscScalar q; /* noise strength */ 34c4762a1bSJed Brown PetscScalar mux; /* Initial average angle */ 35c4762a1bSJed Brown PetscScalar sigmax; /* Standard deviation of initial angle */ 36c4762a1bSJed Brown PetscScalar muy; /* Average speed */ 37c4762a1bSJed Brown PetscScalar sigmay; /* standard deviation of initial speed */ 38c4762a1bSJed Brown PetscScalar rho; /* Cross-correlation coefficient */ 39c4762a1bSJed Brown PetscScalar xmin; /* left boundary of angle */ 40c4762a1bSJed Brown PetscScalar xmax; /* right boundary of angle */ 41c4762a1bSJed Brown PetscScalar ymin; /* bottom boundary of speed */ 42c4762a1bSJed Brown PetscScalar ymax; /* top boundary of speed */ 43c4762a1bSJed Brown PetscScalar dx; /* x step size */ 44c4762a1bSJed Brown PetscScalar dy; /* y step size */ 45c4762a1bSJed Brown PetscScalar disper_coe; /* Dispersion coefficient */ 46c4762a1bSJed Brown DM da; 47c4762a1bSJed Brown PetscInt st_width; /* Stencil width */ 48c4762a1bSJed Brown DMBoundaryType bx; /* x boundary type */ 49c4762a1bSJed Brown DMBoundaryType by; /* y boundary type */ 50c4762a1bSJed Brown PetscReal tf, tcl; /* Fault incidence and clearing times */ 51c4762a1bSJed Brown } AppCtx; 52c4762a1bSJed Brown 53c4762a1bSJed Brown PetscErrorCode Parameter_settings(AppCtx *); 54c4762a1bSJed Brown PetscErrorCode ini_bou(Vec, AppCtx *); 55c4762a1bSJed Brown PetscErrorCode IFunction(TS, PetscReal, Vec, Vec, Vec, void *); 56c4762a1bSJed Brown PetscErrorCode IJacobian(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *); 57c4762a1bSJed Brown PetscErrorCode PostStep(TS); 58c4762a1bSJed Brown 59*9371c9d4SSatish Balay int main(int argc, char **argv) { 60c4762a1bSJed Brown Vec x; /* Solution vector */ 61c4762a1bSJed Brown TS ts; /* Time-stepping context */ 62c4762a1bSJed Brown AppCtx user; /* Application context */ 63c4762a1bSJed Brown PetscViewer viewer; 64c4762a1bSJed Brown 65327415f7SBarry Smith PetscFunctionBeginUser; 669566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, "petscopt_ex8", help)); 67c4762a1bSJed Brown 68c4762a1bSJed Brown /* Get physics and time parameters */ 699566063dSJacob Faibussowitsch PetscCall(Parameter_settings(&user)); 70c4762a1bSJed Brown /* Create a 2D DA with dof = 1 */ 719566063dSJacob 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)); 729566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(user.da)); 739566063dSJacob Faibussowitsch PetscCall(DMSetUp(user.da)); 74c4762a1bSJed Brown /* Set x and y coordinates */ 759566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(user.da, user.xmin, user.xmax, user.ymin, user.ymax, 0, 0)); 769566063dSJacob Faibussowitsch PetscCall(DMDASetCoordinateName(user.da, 0, "X - the angle")); 779566063dSJacob Faibussowitsch PetscCall(DMDASetCoordinateName(user.da, 1, "Y - the speed")); 78c4762a1bSJed Brown 79c4762a1bSJed Brown /* Get global vector x from DM */ 809566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(user.da, &x)); 81c4762a1bSJed Brown 829566063dSJacob Faibussowitsch PetscCall(ini_bou(x, &user)); 839566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "ini_x", FILE_MODE_WRITE, &viewer)); 849566063dSJacob Faibussowitsch PetscCall(VecView(x, viewer)); 859566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 86c4762a1bSJed Brown 879566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 889566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, user.da)); 899566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR)); 909566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSARKIMEX)); 919566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts, NULL, IFunction, &user)); 929566063dSJacob Faibussowitsch /* PetscCall(TSSetIJacobian(ts,NULL,NULL,IJacobian,&user)); */ 939566063dSJacob Faibussowitsch PetscCall(TSSetApplicationContext(ts, &user)); 949566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, .005)); 959566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 969566063dSJacob Faibussowitsch PetscCall(TSSetPostStep(ts, PostStep)); 979566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, x)); 98c4762a1bSJed Brown 999566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "fin_x", FILE_MODE_WRITE, &viewer)); 1009566063dSJacob Faibussowitsch PetscCall(VecView(x, viewer)); 1019566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 102c4762a1bSJed Brown 1039566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1049566063dSJacob Faibussowitsch PetscCall(DMDestroy(&user.da)); 1059566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 1069566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 107b122ec5aSJacob Faibussowitsch return 0; 108c4762a1bSJed Brown } 109c4762a1bSJed Brown 110*9371c9d4SSatish Balay PetscErrorCode PostStep(TS ts) { 111c4762a1bSJed Brown Vec X; 112c4762a1bSJed Brown AppCtx *user; 113c4762a1bSJed Brown PetscReal t; 114c4762a1bSJed Brown PetscScalar asum; 115c4762a1bSJed Brown 116c4762a1bSJed Brown PetscFunctionBegin; 1179566063dSJacob Faibussowitsch PetscCall(TSGetApplicationContext(ts, &user)); 1189566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &t)); 1199566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ts, &X)); 120c4762a1bSJed Brown /* 121c4762a1bSJed Brown if (t >= .2) { 1229566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ts,&X)); 1239566063dSJacob Faibussowitsch PetscCall(VecView(X,PETSC_VIEWER_BINARY_WORLD)); 124c4762a1bSJed Brown exit(0); 125c4762a1bSJed Brown results in initial conditions after fault in binaryoutput 126c4762a1bSJed Brown }*/ 127c4762a1bSJed Brown 128c4762a1bSJed Brown if ((t > user->tf) && (t < user->tcl)) user->Pmax = 0.0; /* A short-circuit that drives the electrical power output (Pmax*sin(delta)) to zero */ 129c4762a1bSJed Brown else user->Pmax = user->Pmax_s; 130c4762a1bSJed Brown 1319566063dSJacob Faibussowitsch PetscCall(VecSum(X, &asum)); 1329566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "sum(p) at t = %f = %f\n", (double)t, (double)(asum))); 133c4762a1bSJed Brown PetscFunctionReturn(0); 134c4762a1bSJed Brown } 135c4762a1bSJed Brown 136*9371c9d4SSatish Balay PetscErrorCode ini_bou(Vec X, AppCtx *user) { 137c4762a1bSJed Brown DM cda; 138c4762a1bSJed Brown DMDACoor2d **coors; 139c4762a1bSJed Brown PetscScalar **p; 140c4762a1bSJed Brown Vec gc; 141c4762a1bSJed Brown PetscInt M, N, Ir, J; 142c4762a1bSJed Brown PetscMPIInt rank; 143c4762a1bSJed Brown 144c4762a1bSJed Brown PetscFunctionBeginUser; 1459566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 1469566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(user->da, NULL, &M, &N, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 147*9371c9d4SSatish Balay user->dx = (user->xmax - user->xmin) / (M - 1); 148*9371c9d4SSatish Balay user->dy = (user->ymax - user->ymin) / (N - 1); 1499566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(user->da, &cda)); 1509566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(user->da, &gc)); 1519566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(cda, gc, &coors)); 1529566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(user->da, X, &p)); 153c4762a1bSJed Brown 154c4762a1bSJed Brown /* Point mass at (mux,muy) */ 15563a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Original user->mux = %f, user->muy = %f\n", (double)PetscRealPart(user->mux), (double)PetscRealPart(user->muy))); 1569566063dSJacob Faibussowitsch PetscCall(DMDAGetLogicalCoordinate(user->da, user->mux, user->muy, 0.0, &Ir, &J, NULL, &user->mux, &user->muy, NULL)); 157c4762a1bSJed Brown user->PM_min = user->Pmax * PetscSinScalar(user->mux); 158*9371c9d4SSatish Balay PetscCall( 159*9371c9d4SSatish Balay PetscPrintf(PETSC_COMM_WORLD, "Corrected user->mux = %f, user->muy = %f user->PM_min = %f,user->dx = %f\n", (double)PetscRealPart(user->mux), (double)PetscRealPart(user->muy), (double)PetscRealPart(user->PM_min), (double)PetscRealPart(user->dx))); 160*9371c9d4SSatish Balay if (Ir > -1 && J > -1) { p[J][Ir] = 1.0; } 161c4762a1bSJed Brown 1629566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(cda, gc, &coors)); 1639566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(user->da, X, &p)); 164c4762a1bSJed Brown PetscFunctionReturn(0); 165c4762a1bSJed Brown } 166c4762a1bSJed Brown 167c4762a1bSJed Brown /* First advection term */ 168*9371c9d4SSatish Balay PetscErrorCode adv1(PetscScalar **p, PetscScalar y, PetscInt i, PetscInt j, PetscInt M, PetscScalar *p1, AppCtx *user) { 169c4762a1bSJed Brown PetscScalar f, fpos, fneg; 170c4762a1bSJed Brown PetscFunctionBegin; 171c4762a1bSJed Brown f = (y - user->ws); 172c4762a1bSJed Brown fpos = PetscMax(f, 0); 173c4762a1bSJed Brown fneg = PetscMin(f, 0); 174c4762a1bSJed Brown if (user->st_width == 1) { 175c4762a1bSJed Brown *p1 = fpos * (p[j][i] - p[j][i - 1]) / user->dx + fneg * (p[j][i + 1] - p[j][i]) / user->dx; 176c4762a1bSJed Brown } else if (user->st_width == 2) { 177c4762a1bSJed 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); 178c4762a1bSJed Brown } else if (user->st_width == 3) { 179c4762a1bSJed 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); 180c4762a1bSJed Brown } 181c4762a1bSJed Brown /* *p1 = f*(p[j][i+1] - p[j][i-1])/user->dx;*/ 182c4762a1bSJed Brown PetscFunctionReturn(0); 183c4762a1bSJed Brown } 184c4762a1bSJed Brown 185c4762a1bSJed Brown /* Second advection term */ 186*9371c9d4SSatish Balay PetscErrorCode adv2(PetscScalar **p, PetscScalar x, PetscScalar y, PetscInt i, PetscInt j, PetscInt N, PetscScalar *p2, AppCtx *user) { 187c4762a1bSJed Brown PetscScalar f, fpos, fneg; 188c4762a1bSJed Brown PetscFunctionBegin; 189c4762a1bSJed Brown f = (user->ws / (2 * user->H)) * (user->PM_min - user->Pmax * PetscSinScalar(x) - user->D * (y - user->ws)); 190c4762a1bSJed Brown fpos = PetscMax(f, 0); 191c4762a1bSJed Brown fneg = PetscMin(f, 0); 192c4762a1bSJed Brown if (user->st_width == 1) { 193c4762a1bSJed Brown *p2 = fpos * (p[j][i] - p[j - 1][i]) / user->dy + fneg * (p[j + 1][i] - p[j][i]) / user->dy; 194c4762a1bSJed Brown } else if (user->st_width == 2) { 195c4762a1bSJed 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); 196c4762a1bSJed Brown } else if (user->st_width == 3) { 197c4762a1bSJed 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); 198c4762a1bSJed Brown } 199c4762a1bSJed Brown 200c4762a1bSJed Brown /* *p2 = f*(p[j+1][i] - p[j-1][i])/user->dy;*/ 201c4762a1bSJed Brown PetscFunctionReturn(0); 202c4762a1bSJed Brown } 203c4762a1bSJed Brown 204c4762a1bSJed Brown /* Diffusion term */ 205*9371c9d4SSatish Balay PetscErrorCode diffuse(PetscScalar **p, PetscInt i, PetscInt j, PetscReal t, PetscScalar *p_diff, AppCtx *user) { 206c4762a1bSJed Brown PetscFunctionBeginUser; 207c4762a1bSJed Brown if (user->st_width == 1) { 208c4762a1bSJed Brown *p_diff = user->disper_coe * ((p[j - 1][i] - 2 * p[j][i] + p[j + 1][i]) / (user->dy * user->dy)); 209c4762a1bSJed Brown } else if (user->st_width == 2) { 210c4762a1bSJed 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)); 211c4762a1bSJed Brown } else if (user->st_width == 3) { 212c4762a1bSJed 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)); 213c4762a1bSJed Brown } 214c4762a1bSJed Brown PetscFunctionReturn(0); 215c4762a1bSJed Brown } 216c4762a1bSJed Brown 217*9371c9d4SSatish Balay PetscErrorCode IFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx) { 218c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx; 219c4762a1bSJed Brown DM cda; 220c4762a1bSJed Brown DMDACoor2d **coors; 221c4762a1bSJed Brown PetscScalar **p, **f, **pdot; 222c4762a1bSJed Brown PetscInt i, j; 223c4762a1bSJed Brown PetscInt xs, ys, xm, ym, M, N; 224c4762a1bSJed Brown Vec localX, gc, localXdot; 225c4762a1bSJed Brown PetscScalar p_adv1 = 0.0, p_adv2 = 0.0, p_diff = 0; 226c4762a1bSJed Brown PetscScalar diffuse1, gamma; 227c4762a1bSJed Brown 228c4762a1bSJed Brown PetscFunctionBeginUser; 2299566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(user->da, NULL, &M, &N, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 2309566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(user->da, &cda)); 2319566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(cda, &xs, &ys, 0, &xm, &ym, 0)); 232c4762a1bSJed Brown 2339566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(user->da, &localX)); 2349566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(user->da, &localXdot)); 235c4762a1bSJed Brown 2369566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(user->da, X, INSERT_VALUES, localX)); 2379566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(user->da, X, INSERT_VALUES, localX)); 2389566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(user->da, Xdot, INSERT_VALUES, localXdot)); 2399566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(user->da, Xdot, INSERT_VALUES, localXdot)); 240c4762a1bSJed Brown 2419566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(user->da, &gc)); 242c4762a1bSJed Brown 2439566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(cda, gc, &coors)); 2449566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(user->da, localX, &p)); 2459566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(user->da, localXdot, &pdot)); 2469566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(user->da, F, &f)); 247c4762a1bSJed Brown 248c4762a1bSJed Brown gamma = user->D * user->ws / (2 * user->H); 249c4762a1bSJed Brown diffuse1 = user->lambda * user->lambda * user->q / (user->lambda * gamma + 1) * (1.0 - PetscExpScalar(-t * (gamma + 1.0) / user->lambda)); 250c4762a1bSJed Brown user->disper_coe = user->ws * user->ws / (4 * user->H * user->H) * diffuse1; 251c4762a1bSJed Brown 252c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 253c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 2549566063dSJacob Faibussowitsch PetscCall(adv1(p, coors[j][i].y, i, j, M, &p_adv1, user)); 2559566063dSJacob Faibussowitsch PetscCall(adv2(p, coors[j][i].x, coors[j][i].y, i, j, N, &p_adv2, user)); 2569566063dSJacob Faibussowitsch PetscCall(diffuse(p, i, j, t, &p_diff, user)); 257c4762a1bSJed Brown f[j][i] = -p_adv1 - p_adv2 + p_diff - pdot[j][i]; 258c4762a1bSJed Brown } 259c4762a1bSJed Brown } 2609566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(user->da, localX, &p)); 2619566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(user->da, localX, &pdot)); 2629566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(user->da, &localX)); 2639566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(user->da, &localXdot)); 2649566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(user->da, F, &f)); 2659566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(cda, gc, &coors)); 266c4762a1bSJed Brown 267c4762a1bSJed Brown PetscFunctionReturn(0); 268c4762a1bSJed Brown } 269c4762a1bSJed Brown 270*9371c9d4SSatish Balay PetscErrorCode IJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat J, Mat Jpre, void *ctx) { 271c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx; 272c4762a1bSJed Brown DM cda; 273c4762a1bSJed Brown DMDACoor2d **coors; 274c4762a1bSJed Brown PetscInt i, j; 275c4762a1bSJed Brown PetscInt xs, ys, xm, ym, M, N; 276c4762a1bSJed Brown Vec gc; 277c4762a1bSJed Brown PetscScalar val[5], xi, yi; 278c4762a1bSJed Brown MatStencil row, col[5]; 279c4762a1bSJed Brown PetscScalar c1, c3, c5, c1pos, c1neg, c3pos, c3neg; 280c4762a1bSJed Brown 281c4762a1bSJed Brown PetscFunctionBeginUser; 2829566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(user->da, NULL, &M, &N, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 2839566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(user->da, &cda)); 2849566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(cda, &xs, &ys, 0, &xm, &ym, 0)); 285c4762a1bSJed Brown 2869566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(user->da, &gc)); 2879566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(cda, gc, &coors)); 288c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 289c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 290c4762a1bSJed Brown PetscInt nc = 0; 291*9371c9d4SSatish Balay xi = coors[j][i].x; 292*9371c9d4SSatish Balay yi = coors[j][i].y; 293*9371c9d4SSatish Balay row.i = i; 294*9371c9d4SSatish Balay row.j = j; 295c4762a1bSJed Brown c1 = (yi - user->ws) / user->dx; 296c4762a1bSJed Brown c1pos = PetscMax(c1, 0); 297c4762a1bSJed Brown c1neg = PetscMin(c1, 0); 298c4762a1bSJed Brown c3 = (user->ws / (2.0 * user->H)) * (user->PM_min - user->Pmax * PetscSinScalar(xi) - user->D * (yi - user->ws)) / user->dy; 299c4762a1bSJed Brown c3pos = PetscMax(c3, 0); 300c4762a1bSJed Brown c3neg = PetscMin(c3, 0); 301c4762a1bSJed Brown c5 = (PetscPowScalar((user->lambda * user->ws) / (2 * user->H), 2) * user->q * (1.0 - PetscExpScalar(-t / user->lambda))) / (user->dy * user->dy); 302*9371c9d4SSatish Balay col[nc].i = i - 1; 303*9371c9d4SSatish Balay col[nc].j = j; 304*9371c9d4SSatish Balay val[nc++] = c1pos; 305*9371c9d4SSatish Balay col[nc].i = i + 1; 306*9371c9d4SSatish Balay col[nc].j = j; 307*9371c9d4SSatish Balay val[nc++] = -c1neg; 308*9371c9d4SSatish Balay col[nc].i = i; 309*9371c9d4SSatish Balay col[nc].j = j - 1; 310*9371c9d4SSatish Balay val[nc++] = c3pos + c5; 311*9371c9d4SSatish Balay col[nc].i = i; 312*9371c9d4SSatish Balay col[nc].j = j + 1; 313*9371c9d4SSatish Balay val[nc++] = -c3neg + c5; 314*9371c9d4SSatish Balay col[nc].i = i; 315*9371c9d4SSatish Balay col[nc].j = j; 316*9371c9d4SSatish Balay val[nc++] = -c1pos + c1neg - c3pos + c3neg - 2 * c5 - a; 3179566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(Jpre, 1, &row, nc, col, val, INSERT_VALUES)); 318c4762a1bSJed Brown } 319c4762a1bSJed Brown } 3209566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(cda, gc, &coors)); 321c4762a1bSJed Brown 3229566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Jpre, MAT_FINAL_ASSEMBLY)); 3239566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Jpre, MAT_FINAL_ASSEMBLY)); 324c4762a1bSJed Brown if (J != Jpre) { 3259566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 3269566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 327c4762a1bSJed Brown } 328c4762a1bSJed Brown PetscFunctionReturn(0); 329c4762a1bSJed Brown } 330c4762a1bSJed Brown 331*9371c9d4SSatish Balay PetscErrorCode Parameter_settings(AppCtx *user) { 332c4762a1bSJed Brown PetscBool flg; 333c4762a1bSJed Brown 334c4762a1bSJed Brown PetscFunctionBeginUser; 335c4762a1bSJed Brown /* Set default parameters */ 336c4762a1bSJed Brown user->ws = 1.0; 337c4762a1bSJed Brown user->H = 5.0; 338c4762a1bSJed Brown user->D = 0.0; 339c4762a1bSJed Brown user->Pmax = user->Pmax_s = 2.1; 340c4762a1bSJed Brown user->PM_min = 1.0; 341c4762a1bSJed Brown user->lambda = 0.1; 342c4762a1bSJed Brown user->q = 1.0; 343c4762a1bSJed Brown user->mux = PetscAsinScalar(user->PM_min / user->Pmax); 344c4762a1bSJed Brown user->sigmax = 0.1; 345c4762a1bSJed Brown user->sigmay = 0.1; 346c4762a1bSJed Brown user->rho = 0.0; 347c4762a1bSJed Brown user->xmin = -PETSC_PI; 348c4762a1bSJed Brown user->xmax = PETSC_PI; 349c4762a1bSJed Brown user->bx = DM_BOUNDARY_PERIODIC; 350c4762a1bSJed Brown user->by = DM_BOUNDARY_GHOSTED; 351c4762a1bSJed Brown user->tf = user->tcl = -1; 352c4762a1bSJed Brown user->ymin = -2.0; 353c4762a1bSJed Brown user->ymax = 2.0; 354c4762a1bSJed Brown user->st_width = 1; 355c4762a1bSJed Brown 3569566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-ws", &user->ws, &flg)); 3579566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-Inertia", &user->H, &flg)); 3589566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-D", &user->D, &flg)); 3599566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-Pmax", &user->Pmax, &flg)); 3609566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-PM_min", &user->PM_min, &flg)); 3619566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-lambda", &user->lambda, &flg)); 3629566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-q", &user->q, &flg)); 3639566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-mux", &user->mux, &flg)); 3649566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-muy", &user->muy, &flg)); 3655f80ce2aSJacob Faibussowitsch if (flg == 0) user->muy = user->ws; 3669566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-xmin", &user->xmin, &flg)); 3679566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-xmax", &user->xmax, &flg)); 3689566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-ymin", &user->ymin, &flg)); 3699566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-ymax", &user->ymax, &flg)); 3709566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-stencil_width", &user->st_width, &flg)); 3719566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetEnum(NULL, NULL, "-bx", DMBoundaryTypes, (PetscEnum *)&user->bx, &flg)); 3729566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetEnum(NULL, NULL, "-by", DMBoundaryTypes, (PetscEnum *)&user->by, &flg)); 3739566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-tf", &user->tf, &flg)); 3749566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-tcl", &user->tcl, &flg)); 375c4762a1bSJed Brown PetscFunctionReturn(0); 376c4762a1bSJed Brown } 377c4762a1bSJed Brown 378c4762a1bSJed Brown /*TEST 379c4762a1bSJed Brown 380c4762a1bSJed Brown build: 381c4762a1bSJed Brown requires: !complex x 382c4762a1bSJed Brown 383c4762a1bSJed Brown test: 384c4762a1bSJed Brown args: -ts_max_steps 1 385c4762a1bSJed Brown localrunfiles: petscopt_ex8 386c4762a1bSJed Brown timeoutfactor: 3 387c4762a1bSJed Brown 388c4762a1bSJed Brown TEST*/ 389