1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Solves heat equation in 1d.\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown /* 5c4762a1bSJed Brown Solves the equation 6c4762a1bSJed Brown 7c4762a1bSJed Brown u_t = kappa \Delta u 8c4762a1bSJed Brown Periodic boundary conditions 9c4762a1bSJed Brown 10c4762a1bSJed Brown Evolve the heat equation: 11c4762a1bSJed Brown --------------- 12c4762a1bSJed Brown ./heat -ts_monitor -snes_monitor -pc_type lu -draw_pause .1 -snes_converged_reason -ts_type cn -da_refine 5 -mymonitor 13c4762a1bSJed Brown 14c4762a1bSJed Brown Evolve the Allen-Cahn equation: 15c4762a1bSJed Brown --------------- 16c4762a1bSJed Brown ./heat -ts_monitor -snes_monitor -pc_type lu -draw_pause .1 -snes_converged_reason -ts_type cn -da_refine 5 -allen-cahn -kappa .001 -ts_max_time 5 -mymonitor 17c4762a1bSJed Brown 18c4762a1bSJed Brown Evolve the Allen-Cahn equation: zoom in on part of the domain 19c4762a1bSJed Brown --------------- 20c4762a1bSJed Brown ./heat -ts_monitor -snes_monitor -pc_type lu -snes_converged_reason -ts_type cn -da_refine 5 -allen-cahn -kappa .001 -ts_max_time 5 -zoom .25,.45 -mymonitor 21c4762a1bSJed Brown 22c4762a1bSJed Brown The option -square_initial indicates it should use a square wave initial condition otherwise it loads the file InitialSolution.heat as the initial solution. You should run with 23c4762a1bSJed Brown ./heat -square_initial -ts_monitor -snes_monitor -pc_type lu -snes_converged_reason -ts_type cn -da_refine 9 -ts_max_time 1.e-4 -ts_dt .125e-6 -snes_atol 1.e-25 -snes_rtol 1.e-25 -ts_max_steps 15 24c4762a1bSJed Brown to generate InitialSolution.heat 25c4762a1bSJed Brown 26c4762a1bSJed Brown */ 27c4762a1bSJed Brown #include <petscdm.h> 28c4762a1bSJed Brown #include <petscdmda.h> 29c4762a1bSJed Brown #include <petscts.h> 30c4762a1bSJed Brown #include <petscdraw.h> 31c4762a1bSJed Brown 32c4762a1bSJed Brown /* 33c4762a1bSJed Brown User-defined routines 34c4762a1bSJed Brown */ 35c4762a1bSJed Brown extern PetscErrorCode FormFunction(TS, PetscReal, Vec, Vec, void *), FormInitialSolution(DM, Vec), MyMonitor(TS, PetscInt, PetscReal, Vec, void *), MyDestroy(void **); 369371c9d4SSatish Balay typedef struct { 379371c9d4SSatish Balay PetscReal kappa; 389371c9d4SSatish Balay PetscBool allencahn; 399371c9d4SSatish Balay PetscDrawViewPorts *ports; 409371c9d4SSatish Balay } UserCtx; 41c4762a1bSJed Brown 429371c9d4SSatish Balay int main(int argc, char **argv) { 43c4762a1bSJed Brown TS ts; /* time integrator */ 44c4762a1bSJed Brown Vec x, r; /* solution, residual vectors */ 45c4762a1bSJed Brown PetscInt steps, Mx; 46c4762a1bSJed Brown DM da; 47c4762a1bSJed Brown PetscReal dt; 48c4762a1bSJed Brown UserCtx ctx; 49c4762a1bSJed Brown PetscBool mymonitor; 50c4762a1bSJed Brown PetscViewer viewer; 51c4762a1bSJed Brown PetscBool flg; 52c4762a1bSJed Brown 53c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 54c4762a1bSJed Brown Initialize program 55c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 56327415f7SBarry Smith PetscFunctionBeginUser; 579566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 58c4762a1bSJed Brown ctx.kappa = 1.0; 599566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-kappa", &ctx.kappa, NULL)); 60c4762a1bSJed Brown ctx.allencahn = PETSC_FALSE; 619566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-allen-cahn", &ctx.allencahn)); 629566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-mymonitor", &mymonitor)); 63c4762a1bSJed Brown 64c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 65c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 66c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 679566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, 10, 1, 2, NULL, &da)); 689566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 699566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 709566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da, 0, "Heat equation: u")); 719566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0)); 72c4762a1bSJed Brown dt = 1.0 / (ctx.kappa * Mx * Mx); 73c4762a1bSJed Brown 74c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 75c4762a1bSJed Brown Extract global vectors from DMDA; then duplicate for remaining 76c4762a1bSJed Brown vectors that are the same types 77c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 789566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &x)); 799566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &r)); 80c4762a1bSJed Brown 81c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 82c4762a1bSJed Brown Create timestepping solver context 83c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 849566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 859566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, da)); 869566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR)); 879566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, FormFunction, &ctx)); 88c4762a1bSJed Brown 89c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 90c4762a1bSJed Brown Customize nonlinear solver 91c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 929566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSCN)); 93c4762a1bSJed Brown 94c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 95c4762a1bSJed Brown Set initial conditions 96c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 979566063dSJacob Faibussowitsch PetscCall(FormInitialSolution(da, x)); 989566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, dt)); 999566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, .02)); 1009566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_INTERPOLATE)); 1019566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, x)); 102c4762a1bSJed Brown 103c4762a1bSJed Brown if (mymonitor) { 104c4762a1bSJed Brown ctx.ports = NULL; 1059566063dSJacob Faibussowitsch PetscCall(TSMonitorSet(ts, MyMonitor, &ctx, MyDestroy)); 106c4762a1bSJed Brown } 107c4762a1bSJed Brown 108c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 109c4762a1bSJed Brown Set runtime options 110c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1119566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 112c4762a1bSJed Brown 113c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 114c4762a1bSJed Brown Solve nonlinear system 115c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1169566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, x)); 1179566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &steps)); 1189566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-square_initial", &flg)); 119c4762a1bSJed Brown if (flg) { 1209566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "InitialSolution.heat", FILE_MODE_WRITE, &viewer)); 1219566063dSJacob Faibussowitsch PetscCall(VecView(x, viewer)); 1229566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 123c4762a1bSJed Brown } 124c4762a1bSJed Brown 125c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 126c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 127c4762a1bSJed Brown are no longer needed. 128c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1299566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1309566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 1319566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 1329566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 133c4762a1bSJed Brown 1349566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 135b122ec5aSJacob Faibussowitsch return 0; 136c4762a1bSJed Brown } 137c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 138c4762a1bSJed Brown /* 139c4762a1bSJed Brown FormFunction - Evaluates nonlinear function, F(x). 140c4762a1bSJed Brown 141c4762a1bSJed Brown Input Parameters: 142c4762a1bSJed Brown . ts - the TS context 143c4762a1bSJed Brown . X - input vector 144c4762a1bSJed Brown . ptr - optional user-defined context, as set by SNESSetFunction() 145c4762a1bSJed Brown 146c4762a1bSJed Brown Output Parameter: 147c4762a1bSJed Brown . F - function vector 148c4762a1bSJed Brown */ 1499371c9d4SSatish Balay PetscErrorCode FormFunction(TS ts, PetscReal ftime, Vec X, Vec F, void *ptr) { 150c4762a1bSJed Brown DM da; 151c4762a1bSJed Brown PetscInt i, Mx, xs, xm; 152c4762a1bSJed Brown PetscReal hx, sx; 153c4762a1bSJed Brown PetscScalar *x, *f; 154c4762a1bSJed Brown Vec localX; 155c4762a1bSJed Brown UserCtx *ctx = (UserCtx *)ptr; 156c4762a1bSJed Brown 157c4762a1bSJed Brown PetscFunctionBegin; 1589566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da)); 1599566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localX)); 1609566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &Mx, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE)); 161c4762a1bSJed Brown 1629371c9d4SSatish Balay hx = 1.0 / (PetscReal)Mx; 1639371c9d4SSatish Balay sx = 1.0 / (hx * hx); 164c4762a1bSJed Brown 165c4762a1bSJed Brown /* 166c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process 167c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 168c4762a1bSJed Brown By placing code between these two statements, computations can be 169c4762a1bSJed Brown done while messages are in transition. 170c4762a1bSJed Brown */ 1719566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX)); 1729566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX)); 173c4762a1bSJed Brown 174c4762a1bSJed Brown /* 175c4762a1bSJed Brown Get pointers to vector data 176c4762a1bSJed Brown */ 1779566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, localX, &x)); 1789566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, F, &f)); 179c4762a1bSJed Brown 180c4762a1bSJed Brown /* 181c4762a1bSJed Brown Get local grid boundaries 182c4762a1bSJed Brown */ 1839566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, NULL, NULL, &xm, NULL, NULL)); 184c4762a1bSJed Brown 185c4762a1bSJed Brown /* 186c4762a1bSJed Brown Compute function over the locally owned part of the grid 187c4762a1bSJed Brown */ 188c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 189c4762a1bSJed Brown f[i] = ctx->kappa * (x[i - 1] + x[i + 1] - 2.0 * x[i]) * sx; 190c4762a1bSJed Brown if (ctx->allencahn) f[i] += (x[i] - x[i] * x[i] * x[i]); 191c4762a1bSJed Brown } 192c4762a1bSJed Brown 193c4762a1bSJed Brown /* 194c4762a1bSJed Brown Restore vectors 195c4762a1bSJed Brown */ 1969566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, localX, &x)); 1979566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, F, &f)); 1989566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localX)); 199c4762a1bSJed Brown PetscFunctionReturn(0); 200c4762a1bSJed Brown } 201c4762a1bSJed Brown 202c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 2039371c9d4SSatish Balay PetscErrorCode FormInitialSolution(DM da, Vec U) { 204c4762a1bSJed Brown PetscInt i, xs, xm, Mx, scale = 1, N; 205c4762a1bSJed Brown PetscScalar *u; 206c4762a1bSJed Brown const PetscScalar *f; 207c4762a1bSJed Brown PetscReal hx, x, r; 208c4762a1bSJed Brown Vec finesolution; 209c4762a1bSJed Brown PetscViewer viewer; 210c4762a1bSJed Brown PetscBool flg; 211c4762a1bSJed Brown 212c4762a1bSJed Brown PetscFunctionBegin; 2139566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &Mx, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE)); 214c4762a1bSJed Brown 215c4762a1bSJed Brown hx = 1.0 / (PetscReal)Mx; 216c4762a1bSJed Brown 217c4762a1bSJed Brown /* 218c4762a1bSJed Brown Get pointers to vector data 219c4762a1bSJed Brown */ 2209566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, U, &u)); 221c4762a1bSJed Brown 222c4762a1bSJed Brown /* 223c4762a1bSJed Brown Get local grid boundaries 224c4762a1bSJed Brown */ 2259566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, NULL, NULL, &xm, NULL, NULL)); 226c4762a1bSJed Brown 227c4762a1bSJed Brown /* InitialSolution is obtained with 228c4762a1bSJed Brown ./heat -ts_monitor -snes_monitor -pc_type lu -snes_converged_reason -ts_type cn -da_refine 9 -ts_max_time 1.e-4 -ts_dt .125e-6 -snes_atol 1.e-25 -snes_rtol 1.e-25 -ts_max_steps 15 229c4762a1bSJed Brown */ 2309566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-square_initial", &flg)); 231c4762a1bSJed Brown if (!flg) { 2329566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "InitialSolution.heat", FILE_MODE_READ, &viewer)); 2339566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &finesolution)); 2349566063dSJacob Faibussowitsch PetscCall(VecLoad(finesolution, viewer)); 2359566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 2369566063dSJacob Faibussowitsch PetscCall(VecGetSize(finesolution, &N)); 237c4762a1bSJed Brown scale = N / Mx; 2389566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(finesolution, &f)); 239c4762a1bSJed Brown } 240c4762a1bSJed Brown 241c4762a1bSJed Brown /* 242c4762a1bSJed Brown Compute function over the locally owned part of the grid 243c4762a1bSJed Brown */ 244c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 245c4762a1bSJed Brown x = i * hx; 246c4762a1bSJed Brown r = PetscSqrtScalar((x - .5) * (x - .5)); 247c4762a1bSJed Brown if (r < .125) u[i] = 1.0; 248c4762a1bSJed Brown else u[i] = -.5; 249c4762a1bSJed Brown 250c4762a1bSJed Brown /* With the initial condition above the method is first order in space */ 251c4762a1bSJed Brown /* this is a smooth initial condition so the method becomes second order in space */ 252c4762a1bSJed Brown /*u[i] = PetscSinScalar(2*PETSC_PI*x); */ 253c4762a1bSJed Brown /* u[i] = f[scale*i];*/ 254c4762a1bSJed Brown if (!flg) u[i] = f[scale * i]; 255c4762a1bSJed Brown } 256c4762a1bSJed Brown if (!flg) { 2579566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(finesolution, &f)); 2589566063dSJacob Faibussowitsch PetscCall(VecDestroy(&finesolution)); 259c4762a1bSJed Brown } 260c4762a1bSJed Brown 261c4762a1bSJed Brown /* 262c4762a1bSJed Brown Restore vectors 263c4762a1bSJed Brown */ 2649566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, U, &u)); 265c4762a1bSJed Brown PetscFunctionReturn(0); 266c4762a1bSJed Brown } 267c4762a1bSJed Brown 268c4762a1bSJed Brown /* 269c4762a1bSJed Brown This routine is not parallel 270c4762a1bSJed Brown */ 2719371c9d4SSatish Balay PetscErrorCode MyMonitor(TS ts, PetscInt step, PetscReal time, Vec U, void *ptr) { 272c4762a1bSJed Brown UserCtx *ctx = (UserCtx *)ptr; 273c4762a1bSJed Brown PetscDrawLG lg; 274c4762a1bSJed Brown PetscScalar *u; 275c4762a1bSJed Brown PetscInt Mx, i, xs, xm, cnt; 276c4762a1bSJed Brown PetscReal x, y, hx, pause, sx, len, max, xx[2], yy[2]; 277c4762a1bSJed Brown PetscDraw draw; 278c4762a1bSJed Brown Vec localU; 279c4762a1bSJed Brown DM da; 280c4762a1bSJed Brown int colors[] = {PETSC_DRAW_YELLOW, PETSC_DRAW_RED, PETSC_DRAW_BLUE}; 281c4762a1bSJed Brown const char *const legend[] = {"-kappa (\\grad u,\\grad u)", "(1 - u^2)^2"}; 282c4762a1bSJed Brown PetscDrawAxis axis; 283c4762a1bSJed Brown PetscDrawViewPorts *ports; 284c4762a1bSJed Brown PetscReal vbounds[] = {-1.1, 1.1}; 285c4762a1bSJed Brown 286c4762a1bSJed Brown PetscFunctionBegin; 2879566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD), 1, vbounds)); 2889566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawResize(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD), 1200, 800)); 2899566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da)); 2909566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localU)); 2919566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &Mx, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE)); 2929566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, NULL, NULL, &xm, NULL, NULL)); 2939371c9d4SSatish Balay hx = 1.0 / (PetscReal)Mx; 2949371c9d4SSatish Balay sx = 1.0 / (hx * hx); 2959566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, U, INSERT_VALUES, localU)); 2969566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, U, INSERT_VALUES, localU)); 2979566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, localU, &u)); 298c4762a1bSJed Brown 2999566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDrawLG(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD), 1, &lg)); 3009566063dSJacob Faibussowitsch PetscCall(PetscDrawLGGetDraw(lg, &draw)); 3019566063dSJacob Faibussowitsch PetscCall(PetscDrawCheckResizedWindow(draw)); 302*48a46eb9SPierre Jolivet if (!ctx->ports) PetscCall(PetscDrawViewPortsCreateRect(draw, 1, 3, &ctx->ports)); 303c4762a1bSJed Brown ports = ctx->ports; 3049566063dSJacob Faibussowitsch PetscCall(PetscDrawLGGetAxis(lg, &axis)); 3059566063dSJacob Faibussowitsch PetscCall(PetscDrawLGReset(lg)); 306c4762a1bSJed Brown 3079371c9d4SSatish Balay xx[0] = 0.0; 3089371c9d4SSatish Balay xx[1] = 1.0; 3099371c9d4SSatish Balay cnt = 2; 3109566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetRealArray(NULL, NULL, "-zoom", xx, &cnt, NULL)); 3119371c9d4SSatish Balay xs = xx[0] / hx; 3129371c9d4SSatish Balay xm = (xx[1] - xx[0]) / hx; 313c4762a1bSJed Brown 314c4762a1bSJed Brown /* 315c4762a1bSJed Brown Plot the energies 316c4762a1bSJed Brown */ 3179566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetDimension(lg, 1 + (ctx->allencahn ? 1 : 0))); 3189566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetColors(lg, colors + 1)); 3199566063dSJacob Faibussowitsch PetscCall(PetscDrawViewPortsSet(ports, 2)); 320c4762a1bSJed Brown x = hx * xs; 321c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 322c4762a1bSJed Brown xx[0] = xx[1] = x; 323c4762a1bSJed Brown yy[0] = PetscRealPart(.25 * ctx->kappa * (u[i - 1] - u[i + 1]) * (u[i - 1] - u[i + 1]) * sx); 324c4762a1bSJed Brown if (ctx->allencahn) yy[1] = .25 * PetscRealPart((1. - u[i] * u[i]) * (1. - u[i] * u[i])); 3259566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddPoint(lg, xx, yy)); 326c4762a1bSJed Brown x += hx; 327c4762a1bSJed Brown } 3289566063dSJacob Faibussowitsch PetscCall(PetscDrawGetPause(draw, &pause)); 3299566063dSJacob Faibussowitsch PetscCall(PetscDrawSetPause(draw, 0.0)); 3309566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Energy", "", "")); 3319566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetLegend(lg, legend)); 3329566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDraw(lg)); 333c4762a1bSJed Brown 334c4762a1bSJed Brown /* 335c4762a1bSJed Brown Plot the forces 336c4762a1bSJed Brown */ 3379566063dSJacob Faibussowitsch PetscCall(PetscDrawViewPortsSet(ports, 1)); 3389566063dSJacob Faibussowitsch PetscCall(PetscDrawLGReset(lg)); 339c4762a1bSJed Brown x = xs * hx; 340c4762a1bSJed Brown max = 0.; 341c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 342c4762a1bSJed Brown xx[0] = xx[1] = x; 343c4762a1bSJed Brown yy[0] = PetscRealPart(ctx->kappa * (u[i - 1] + u[i + 1] - 2.0 * u[i]) * sx); 344c4762a1bSJed Brown max = PetscMax(max, PetscAbs(yy[0])); 345c4762a1bSJed Brown if (ctx->allencahn) { 346c4762a1bSJed Brown yy[1] = PetscRealPart(u[i] - u[i] * u[i] * u[i]); 347c4762a1bSJed Brown max = PetscMax(max, PetscAbs(yy[1])); 348c4762a1bSJed Brown } 3499566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddPoint(lg, xx, yy)); 350c4762a1bSJed Brown x += hx; 351c4762a1bSJed Brown } 3529566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Right hand side", "", "")); 3539566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetLegend(lg, NULL)); 3549566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDraw(lg)); 355c4762a1bSJed Brown 356c4762a1bSJed Brown /* 357c4762a1bSJed Brown Plot the solution 358c4762a1bSJed Brown */ 3599566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetDimension(lg, 1)); 3609566063dSJacob Faibussowitsch PetscCall(PetscDrawViewPortsSet(ports, 0)); 3619566063dSJacob Faibussowitsch PetscCall(PetscDrawLGReset(lg)); 362c4762a1bSJed Brown x = hx * xs; 3639566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetLimits(lg, x, x + (xm - 1) * hx, -1.1, 1.1)); 3649566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetColors(lg, colors)); 365c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 366c4762a1bSJed Brown xx[0] = x; 367c4762a1bSJed Brown yy[0] = PetscRealPart(u[i]); 3689566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddPoint(lg, xx, yy)); 369c4762a1bSJed Brown x += hx; 370c4762a1bSJed Brown } 3719566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Solution", "", "")); 3729566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDraw(lg)); 373c4762a1bSJed Brown 374c4762a1bSJed Brown /* 375c4762a1bSJed Brown Print the forces as arrows on the solution 376c4762a1bSJed Brown */ 377c4762a1bSJed Brown x = hx * xs; 378c4762a1bSJed Brown cnt = xm / 60; 379c4762a1bSJed Brown cnt = (!cnt) ? 1 : cnt; 380c4762a1bSJed Brown 381c4762a1bSJed Brown for (i = xs; i < xs + xm; i += cnt) { 382c4762a1bSJed Brown y = PetscRealPart(u[i]); 383c4762a1bSJed Brown len = .5 * PetscRealPart(ctx->kappa * (u[i - 1] + u[i + 1] - 2.0 * u[i]) * sx) / max; 3849566063dSJacob Faibussowitsch PetscCall(PetscDrawArrow(draw, x, y, x, y + len, PETSC_DRAW_RED)); 385c4762a1bSJed Brown if (ctx->allencahn) { 386c4762a1bSJed Brown len = .5 * PetscRealPart(u[i] - u[i] * u[i] * u[i]) / max; 3879566063dSJacob Faibussowitsch PetscCall(PetscDrawArrow(draw, x, y, x, y + len, PETSC_DRAW_BLUE)); 388c4762a1bSJed Brown } 389c4762a1bSJed Brown x += cnt * hx; 390c4762a1bSJed Brown } 3919566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, localU, &x)); 3929566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localU)); 3939566063dSJacob Faibussowitsch PetscCall(PetscDrawStringSetSize(draw, .2, .2)); 3949566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 3959566063dSJacob Faibussowitsch PetscCall(PetscDrawSetPause(draw, pause)); 3969566063dSJacob Faibussowitsch PetscCall(PetscDrawPause(draw)); 397c4762a1bSJed Brown PetscFunctionReturn(0); 398c4762a1bSJed Brown } 399c4762a1bSJed Brown 4009371c9d4SSatish Balay PetscErrorCode MyDestroy(void **ptr) { 401c4762a1bSJed Brown UserCtx *ctx = *(UserCtx **)ptr; 402c4762a1bSJed Brown 403c4762a1bSJed Brown PetscFunctionBegin; 4049566063dSJacob Faibussowitsch PetscCall(PetscDrawViewPortsDestroy(ctx->ports)); 405c4762a1bSJed Brown PetscFunctionReturn(0); 406c4762a1bSJed Brown } 407c4762a1bSJed Brown 408c4762a1bSJed Brown /*TEST 409c4762a1bSJed Brown 410c4762a1bSJed Brown test: 411c4762a1bSJed Brown args: -ts_monitor -snes_monitor -pc_type lu -snes_converged_reason -ts_type cn -da_refine 2 -square_initial 412c4762a1bSJed Brown 413c4762a1bSJed Brown test: 414c4762a1bSJed Brown suffix: 2 415c4762a1bSJed Brown args: -ts_monitor -snes_monitor -pc_type lu -snes_converged_reason -ts_type cn -da_refine 5 -mymonitor -square_initial -allen-cahn -kappa .001 416c4762a1bSJed Brown requires: x 417c4762a1bSJed Brown 418c4762a1bSJed Brown TEST*/ 419