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**); 36c4762a1bSJed Brown typedef struct {PetscReal kappa;PetscBool allencahn;PetscDrawViewPorts *ports;} UserCtx; 37c4762a1bSJed Brown 38c4762a1bSJed Brown int main(int argc,char **argv) 39c4762a1bSJed Brown { 40c4762a1bSJed Brown TS ts; /* time integrator */ 41c4762a1bSJed Brown Vec x,r; /* solution, residual vectors */ 42c4762a1bSJed Brown PetscInt steps,Mx; 43c4762a1bSJed Brown DM da; 44c4762a1bSJed Brown PetscReal dt; 45c4762a1bSJed Brown UserCtx ctx; 46c4762a1bSJed Brown PetscBool mymonitor; 47c4762a1bSJed Brown PetscViewer viewer; 48c4762a1bSJed Brown PetscBool flg; 49c4762a1bSJed Brown 50c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 51c4762a1bSJed Brown Initialize program 52c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 53*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help)); 54c4762a1bSJed Brown ctx.kappa = 1.0; 555f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-kappa",&ctx.kappa,NULL)); 56c4762a1bSJed Brown ctx.allencahn = PETSC_FALSE; 575f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-allen-cahn",&ctx.allencahn)); 585f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-mymonitor",&mymonitor)); 59c4762a1bSJed Brown 60c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 61c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 62c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 635f80ce2aSJacob Faibussowitsch CHKERRQ(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, 10,1,2,NULL,&da)); 645f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(da)); 655f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(da)); 665f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetFieldName(da,0,"Heat equation: u")); 675f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetInfo(da,0,&Mx,0,0,0,0,0,0,0,0,0,0,0)); 68c4762a1bSJed Brown dt = 1.0/(ctx.kappa*Mx*Mx); 69c4762a1bSJed Brown 70c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 71c4762a1bSJed Brown Extract global vectors from DMDA; then duplicate for remaining 72c4762a1bSJed Brown vectors that are the same types 73c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 745f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(da,&x)); 755f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&r)); 76c4762a1bSJed Brown 77c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 78c4762a1bSJed Brown Create timestepping solver context 79c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 805f80ce2aSJacob Faibussowitsch CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts)); 815f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetDM(ts,da)); 825f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetProblemType(ts,TS_NONLINEAR)); 835f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSFunction(ts,NULL,FormFunction,&ctx)); 84c4762a1bSJed Brown 85c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 86c4762a1bSJed Brown Customize nonlinear solver 87c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 885f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetType(ts,TSCN)); 89c4762a1bSJed Brown 90c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 91c4762a1bSJed Brown Set initial conditions 92c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 935f80ce2aSJacob Faibussowitsch CHKERRQ(FormInitialSolution(da,x)); 945f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTimeStep(ts,dt)); 955f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxTime(ts,.02)); 965f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_INTERPOLATE)); 975f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetSolution(ts,x)); 98c4762a1bSJed Brown 99c4762a1bSJed Brown if (mymonitor) { 100c4762a1bSJed Brown ctx.ports = NULL; 1015f80ce2aSJacob Faibussowitsch CHKERRQ(TSMonitorSet(ts,MyMonitor,&ctx,MyDestroy)); 102c4762a1bSJed Brown } 103c4762a1bSJed Brown 104c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 105c4762a1bSJed Brown Set runtime options 106c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1075f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(ts)); 108c4762a1bSJed Brown 109c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 110c4762a1bSJed Brown Solve nonlinear system 111c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1125f80ce2aSJacob Faibussowitsch CHKERRQ(TSSolve(ts,x)); 1135f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetStepNumber(ts,&steps)); 1145f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-square_initial",&flg)); 115c4762a1bSJed Brown if (flg) { 1165f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"InitialSolution.heat",FILE_MODE_WRITE,&viewer)); 1175f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(x,viewer)); 1185f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&viewer)); 119c4762a1bSJed Brown } 120c4762a1bSJed Brown 121c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 122c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 123c4762a1bSJed Brown are no longer needed. 124c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1255f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 1265f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&r)); 1275f80ce2aSJacob Faibussowitsch CHKERRQ(TSDestroy(&ts)); 1285f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&da)); 129c4762a1bSJed Brown 130*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 131*b122ec5aSJacob Faibussowitsch return 0; 132c4762a1bSJed Brown } 133c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 134c4762a1bSJed Brown /* 135c4762a1bSJed Brown FormFunction - Evaluates nonlinear function, F(x). 136c4762a1bSJed Brown 137c4762a1bSJed Brown Input Parameters: 138c4762a1bSJed Brown . ts - the TS context 139c4762a1bSJed Brown . X - input vector 140c4762a1bSJed Brown . ptr - optional user-defined context, as set by SNESSetFunction() 141c4762a1bSJed Brown 142c4762a1bSJed Brown Output Parameter: 143c4762a1bSJed Brown . F - function vector 144c4762a1bSJed Brown */ 145c4762a1bSJed Brown PetscErrorCode FormFunction(TS ts,PetscReal ftime,Vec X,Vec F,void *ptr) 146c4762a1bSJed Brown { 147c4762a1bSJed Brown DM da; 148c4762a1bSJed Brown PetscInt i,Mx,xs,xm; 149c4762a1bSJed Brown PetscReal hx,sx; 150c4762a1bSJed Brown PetscScalar *x,*f; 151c4762a1bSJed Brown Vec localX; 152c4762a1bSJed Brown UserCtx *ctx = (UserCtx*)ptr; 153c4762a1bSJed Brown 154c4762a1bSJed Brown PetscFunctionBegin; 1555f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts,&da)); 1565f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(da,&localX)); 1575f80ce2aSJacob Faibussowitsch CHKERRQ(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)); 158c4762a1bSJed Brown 159c4762a1bSJed Brown hx = 1.0/(PetscReal)Mx; sx = 1.0/(hx*hx); 160c4762a1bSJed Brown 161c4762a1bSJed Brown /* 162c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process 163c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 164c4762a1bSJed Brown By placing code between these two statements, computations can be 165c4762a1bSJed Brown done while messages are in transition. 166c4762a1bSJed Brown */ 1675f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX)); 1685f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX)); 169c4762a1bSJed Brown 170c4762a1bSJed Brown /* 171c4762a1bSJed Brown Get pointers to vector data 172c4762a1bSJed Brown */ 1735f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(da,localX,&x)); 1745f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da,F,&f)); 175c4762a1bSJed Brown 176c4762a1bSJed Brown /* 177c4762a1bSJed Brown Get local grid boundaries 178c4762a1bSJed Brown */ 1795f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL)); 180c4762a1bSJed Brown 181c4762a1bSJed Brown /* 182c4762a1bSJed Brown Compute function over the locally owned part of the grid 183c4762a1bSJed Brown */ 184c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 185c4762a1bSJed Brown f[i] = ctx->kappa*(x[i-1] + x[i+1] - 2.0*x[i])*sx; 186c4762a1bSJed Brown if (ctx->allencahn) f[i] += (x[i] - x[i]*x[i]*x[i]); 187c4762a1bSJed Brown } 188c4762a1bSJed Brown 189c4762a1bSJed Brown /* 190c4762a1bSJed Brown Restore vectors 191c4762a1bSJed Brown */ 1925f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(da,localX,&x)); 1935f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,F,&f)); 1945f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(da,&localX)); 195c4762a1bSJed Brown PetscFunctionReturn(0); 196c4762a1bSJed Brown } 197c4762a1bSJed Brown 198c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 199c4762a1bSJed Brown PetscErrorCode FormInitialSolution(DM da,Vec U) 200c4762a1bSJed Brown { 201c4762a1bSJed Brown PetscInt i,xs,xm,Mx,scale=1,N; 202c4762a1bSJed Brown PetscScalar *u; 203c4762a1bSJed Brown const PetscScalar *f; 204c4762a1bSJed Brown PetscReal hx,x,r; 205c4762a1bSJed Brown Vec finesolution; 206c4762a1bSJed Brown PetscViewer viewer; 207c4762a1bSJed Brown PetscBool flg; 208c4762a1bSJed Brown 209c4762a1bSJed Brown PetscFunctionBegin; 2105f80ce2aSJacob Faibussowitsch CHKERRQ(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)); 211c4762a1bSJed Brown 212c4762a1bSJed Brown hx = 1.0/(PetscReal)Mx; 213c4762a1bSJed Brown 214c4762a1bSJed Brown /* 215c4762a1bSJed Brown Get pointers to vector data 216c4762a1bSJed Brown */ 2175f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(da,U,&u)); 218c4762a1bSJed Brown 219c4762a1bSJed Brown /* 220c4762a1bSJed Brown Get local grid boundaries 221c4762a1bSJed Brown */ 2225f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL)); 223c4762a1bSJed Brown 224c4762a1bSJed Brown /* InitialSolution is obtained with 225c4762a1bSJed 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 226c4762a1bSJed Brown */ 2275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-square_initial",&flg)); 228c4762a1bSJed Brown if (!flg) { 2295f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"InitialSolution.heat",FILE_MODE_READ,&viewer)); 2305f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_WORLD,&finesolution)); 2315f80ce2aSJacob Faibussowitsch CHKERRQ(VecLoad(finesolution,viewer)); 2325f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&viewer)); 2335f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(finesolution,&N)); 234c4762a1bSJed Brown scale = N/Mx; 2355f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(finesolution,&f)); 236c4762a1bSJed Brown } 237c4762a1bSJed Brown 238c4762a1bSJed Brown /* 239c4762a1bSJed Brown Compute function over the locally owned part of the grid 240c4762a1bSJed Brown */ 241c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 242c4762a1bSJed Brown x = i*hx; 243c4762a1bSJed Brown r = PetscSqrtScalar((x-.5)*(x-.5)); 244c4762a1bSJed Brown if (r < .125) u[i] = 1.0; 245c4762a1bSJed Brown else u[i] = -.5; 246c4762a1bSJed Brown 247c4762a1bSJed Brown /* With the initial condition above the method is first order in space */ 248c4762a1bSJed Brown /* this is a smooth initial condition so the method becomes second order in space */ 249c4762a1bSJed Brown /*u[i] = PetscSinScalar(2*PETSC_PI*x); */ 250c4762a1bSJed Brown /* u[i] = f[scale*i];*/ 251c4762a1bSJed Brown if (!flg) u[i] = f[scale*i]; 252c4762a1bSJed Brown } 253c4762a1bSJed Brown if (!flg) { 2545f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(finesolution,&f)); 2555f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&finesolution)); 256c4762a1bSJed Brown } 257c4762a1bSJed Brown 258c4762a1bSJed Brown /* 259c4762a1bSJed Brown Restore vectors 260c4762a1bSJed Brown */ 2615f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,U,&u)); 262c4762a1bSJed Brown PetscFunctionReturn(0); 263c4762a1bSJed Brown } 264c4762a1bSJed Brown 265c4762a1bSJed Brown /* 266c4762a1bSJed Brown This routine is not parallel 267c4762a1bSJed Brown */ 268c4762a1bSJed Brown PetscErrorCode MyMonitor(TS ts,PetscInt step,PetscReal time,Vec U,void *ptr) 269c4762a1bSJed Brown { 270c4762a1bSJed Brown UserCtx *ctx = (UserCtx*)ptr; 271c4762a1bSJed Brown PetscDrawLG lg; 272c4762a1bSJed Brown PetscScalar *u; 273c4762a1bSJed Brown PetscInt Mx,i,xs,xm,cnt; 274c4762a1bSJed Brown PetscReal x,y,hx,pause,sx,len,max,xx[2],yy[2]; 275c4762a1bSJed Brown PetscDraw draw; 276c4762a1bSJed Brown Vec localU; 277c4762a1bSJed Brown DM da; 278c4762a1bSJed Brown int colors[] = {PETSC_DRAW_YELLOW,PETSC_DRAW_RED,PETSC_DRAW_BLUE}; 279c4762a1bSJed Brown const char*const legend[] = {"-kappa (\\grad u,\\grad u)","(1 - u^2)^2"}; 280c4762a1bSJed Brown PetscDrawAxis axis; 281c4762a1bSJed Brown PetscDrawViewPorts *ports; 282c4762a1bSJed Brown PetscReal vbounds[] = {-1.1,1.1}; 283c4762a1bSJed Brown 284c4762a1bSJed Brown PetscFunctionBegin; 2855f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),1,vbounds)); 2865f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDrawResize(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),1200,800)); 2875f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetDM(ts,&da)); 2885f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(da,&localU)); 2895f80ce2aSJacob Faibussowitsch CHKERRQ(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)); 2905f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL)); 291c4762a1bSJed Brown hx = 1.0/(PetscReal)Mx; sx = 1.0/(hx*hx); 2925f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU)); 2935f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU)); 2945f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(da,localU,&u)); 295c4762a1bSJed Brown 2965f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDrawGetDrawLG(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),1,&lg)); 2975f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawLGGetDraw(lg,&draw)); 2985f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawCheckResizedWindow(draw)); 299c4762a1bSJed Brown if (!ctx->ports) { 3005f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawViewPortsCreateRect(draw,1,3,&ctx->ports)); 301c4762a1bSJed Brown } 302c4762a1bSJed Brown ports = ctx->ports; 3035f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawLGGetAxis(lg,&axis)); 3045f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawLGReset(lg)); 305c4762a1bSJed Brown 306c4762a1bSJed Brown xx[0] = 0.0; xx[1] = 1.0; cnt = 2; 3075f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetRealArray(NULL,NULL,"-zoom",xx,&cnt,NULL)); 308c4762a1bSJed Brown xs = xx[0]/hx; xm = (xx[1] - xx[0])/hx; 309c4762a1bSJed Brown 310c4762a1bSJed Brown /* 311c4762a1bSJed Brown Plot the energies 312c4762a1bSJed Brown */ 3135f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawLGSetDimension(lg,1 + (ctx->allencahn ? 1 : 0))); 3145f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawLGSetColors(lg,colors+1)); 3155f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawViewPortsSet(ports,2)); 316c4762a1bSJed Brown x = hx*xs; 317c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 318c4762a1bSJed Brown xx[0] = xx[1] = x; 319c4762a1bSJed Brown yy[0] = PetscRealPart(.25*ctx->kappa*(u[i-1] - u[i+1])*(u[i-1] - u[i+1])*sx); 320c4762a1bSJed Brown if (ctx->allencahn) yy[1] = .25*PetscRealPart((1. - u[i]*u[i])*(1. - u[i]*u[i])); 3215f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawLGAddPoint(lg,xx,yy)); 322c4762a1bSJed Brown x += hx; 323c4762a1bSJed Brown } 3245f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawGetPause(draw,&pause)); 3255f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawSetPause(draw,0.0)); 3265f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawAxisSetLabels(axis,"Energy","","")); 3275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawLGSetLegend(lg,legend)); 3285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawLGDraw(lg)); 329c4762a1bSJed Brown 330c4762a1bSJed Brown /* 331c4762a1bSJed Brown Plot the forces 332c4762a1bSJed Brown */ 3335f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawViewPortsSet(ports,1)); 3345f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawLGReset(lg)); 335c4762a1bSJed Brown x = xs*hx; 336c4762a1bSJed Brown max = 0.; 337c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 338c4762a1bSJed Brown xx[0] = xx[1] = x; 339c4762a1bSJed Brown yy[0] = PetscRealPart(ctx->kappa*(u[i-1] + u[i+1] - 2.0*u[i])*sx); 340c4762a1bSJed Brown max = PetscMax(max,PetscAbs(yy[0])); 341c4762a1bSJed Brown if (ctx->allencahn) { 342c4762a1bSJed Brown yy[1] = PetscRealPart(u[i] - u[i]*u[i]*u[i]); 343c4762a1bSJed Brown max = PetscMax(max,PetscAbs(yy[1])); 344c4762a1bSJed Brown } 3455f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawLGAddPoint(lg,xx,yy)); 346c4762a1bSJed Brown x += hx; 347c4762a1bSJed Brown } 3485f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawAxisSetLabels(axis,"Right hand side","","")); 3495f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawLGSetLegend(lg,NULL)); 3505f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawLGDraw(lg)); 351c4762a1bSJed Brown 352c4762a1bSJed Brown /* 353c4762a1bSJed Brown Plot the solution 354c4762a1bSJed Brown */ 3555f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawLGSetDimension(lg,1)); 3565f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawViewPortsSet(ports,0)); 3575f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawLGReset(lg)); 358c4762a1bSJed Brown x = hx*xs; 3595f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawLGSetLimits(lg,x,x+(xm-1)*hx,-1.1,1.1)); 3605f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawLGSetColors(lg,colors)); 361c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 362c4762a1bSJed Brown xx[0] = x; 363c4762a1bSJed Brown yy[0] = PetscRealPart(u[i]); 3645f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawLGAddPoint(lg,xx,yy)); 365c4762a1bSJed Brown x += hx; 366c4762a1bSJed Brown } 3675f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawAxisSetLabels(axis,"Solution","","")); 3685f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawLGDraw(lg)); 369c4762a1bSJed Brown 370c4762a1bSJed Brown /* 371c4762a1bSJed Brown Print the forces as arrows on the solution 372c4762a1bSJed Brown */ 373c4762a1bSJed Brown x = hx*xs; 374c4762a1bSJed Brown cnt = xm/60; 375c4762a1bSJed Brown cnt = (!cnt) ? 1 : cnt; 376c4762a1bSJed Brown 377c4762a1bSJed Brown for (i=xs; i<xs+xm; i += cnt) { 378c4762a1bSJed Brown y = PetscRealPart(u[i]); 379c4762a1bSJed Brown len = .5*PetscRealPart(ctx->kappa*(u[i-1] + u[i+1] - 2.0*u[i])*sx)/max; 3805f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawArrow(draw,x,y,x,y+len,PETSC_DRAW_RED)); 381c4762a1bSJed Brown if (ctx->allencahn) { 382c4762a1bSJed Brown len = .5*PetscRealPart(u[i] - u[i]*u[i]*u[i])/max; 3835f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawArrow(draw,x,y,x,y+len,PETSC_DRAW_BLUE)); 384c4762a1bSJed Brown } 385c4762a1bSJed Brown x += cnt*hx; 386c4762a1bSJed Brown } 3875f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(da,localU,&x)); 3885f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(da,&localU)); 3895f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawStringSetSize(draw,.2,.2)); 3905f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawFlush(draw)); 3915f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawSetPause(draw,pause)); 3925f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawPause(draw)); 393c4762a1bSJed Brown PetscFunctionReturn(0); 394c4762a1bSJed Brown } 395c4762a1bSJed Brown 396c4762a1bSJed Brown PetscErrorCode MyDestroy(void **ptr) 397c4762a1bSJed Brown { 398c4762a1bSJed Brown UserCtx *ctx = *(UserCtx**)ptr; 399c4762a1bSJed Brown 400c4762a1bSJed Brown PetscFunctionBegin; 4015f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDrawViewPortsDestroy(ctx->ports)); 402c4762a1bSJed Brown PetscFunctionReturn(0); 403c4762a1bSJed Brown } 404c4762a1bSJed Brown 405c4762a1bSJed Brown /*TEST 406c4762a1bSJed Brown 407c4762a1bSJed Brown test: 408c4762a1bSJed Brown args: -ts_monitor -snes_monitor -pc_type lu -snes_converged_reason -ts_type cn -da_refine 2 -square_initial 409c4762a1bSJed Brown 410c4762a1bSJed Brown test: 411c4762a1bSJed Brown suffix: 2 412c4762a1bSJed 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 413c4762a1bSJed Brown requires: x 414c4762a1bSJed Brown 415c4762a1bSJed Brown TEST*/ 416