1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Solves biharmonic equation in 1d.\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown /* 5c4762a1bSJed Brown Solves the equation 6c4762a1bSJed Brown 7c4762a1bSJed Brown u_t = - kappa \Delta \Delta u 8c4762a1bSJed Brown Periodic boundary conditions 9c4762a1bSJed Brown 10c4762a1bSJed Brown Evolve the biharmonic heat equation: 11c4762a1bSJed Brown --------------- 12c4762a1bSJed Brown ./biharmonic -ts_monitor -snes_monitor -pc_type lu -draw_pause .1 -snes_converged_reason -draw_pause -2 -ts_type cn -da_refine 5 -mymonitor 13c4762a1bSJed Brown 14c4762a1bSJed Brown Evolve with the restriction that -1 <= u <= 1; i.e. as a variational inequality 15c4762a1bSJed Brown --------------- 16c4762a1bSJed Brown ./biharmonic -ts_monitor -snes_monitor -pc_type lu -draw_pause .1 -snes_converged_reason -draw_pause -2 -ts_type cn -da_refine 5 -mymonitor 17c4762a1bSJed Brown 18c4762a1bSJed Brown u_t = kappa \Delta \Delta u + 6.*u*(u_x)^2 + (3*u^2 - 12) \Delta u 19c4762a1bSJed Brown -1 <= u <= 1 20c4762a1bSJed Brown Periodic boundary conditions 21c4762a1bSJed Brown 22c4762a1bSJed Brown Evolve the Cahn-Hillard equations: double well Initial hump shrinks then grows 23c4762a1bSJed Brown --------------- 24c4762a1bSJed Brown ./biharmonic -ts_monitor -snes_monitor -pc_type lu -draw_pause .1 -snes_converged_reason -draw_pause -2 -ts_type cn -da_refine 6 -kappa .00001 -ts_dt 5.96046e-06 -cahn-hillard -ts_monitor_draw_solution --mymonitor 25c4762a1bSJed Brown 26c4762a1bSJed Brown Initial hump neither shrinks nor grows when degenerate (otherwise similar solution) 27c4762a1bSJed Brown 28c4762a1bSJed Brown ./biharmonic -ts_monitor -snes_monitor -pc_type lu -draw_pause .1 -snes_converged_reason -draw_pause -2 -ts_type cn -da_refine 6 -kappa .00001 -ts_dt 5.96046e-06 -cahn-hillard -degenerate -ts_monitor_draw_solution --mymonitor 29c4762a1bSJed Brown 30c4762a1bSJed Brown ./biharmonic -ts_monitor -snes_monitor -pc_type lu -draw_pause .1 -snes_converged_reason -draw_pause -2 -ts_type cn -da_refine 6 -kappa .00001 -ts_dt 5.96046e-06 -cahn-hillard -snes_vi_ignore_function_sign -ts_monitor_draw_solution --mymonitor 31c4762a1bSJed Brown 32c4762a1bSJed Brown Evolve the Cahn-Hillard equations: double obstacle 33c4762a1bSJed Brown --------------- 34c4762a1bSJed Brown ./biharmonic -ts_monitor -snes_monitor -pc_type lu -draw_pause .1 -snes_converged_reason -draw_pause -2 -ts_type cn -da_refine 5 -kappa .00001 -ts_dt 5.96046e-06 -cahn-hillard -energy 2 -snes_linesearch_monitor -ts_monitor_draw_solution --mymonitor 35c4762a1bSJed Brown 36c4762a1bSJed Brown Evolve the Cahn-Hillard equations: logarithmic + double well (never shrinks and then grows) 37c4762a1bSJed Brown --------------- 38c4762a1bSJed Brown ./biharmonic -ts_monitor -snes_monitor -pc_type lu --snes_converged_reason -draw_pause -2 -ts_type cn -da_refine 5 -kappa .0001 -ts_dt 5.96046e-06 -cahn-hillard -energy 3 -snes_linesearch_monitor -theta .00000001 -ts_monitor_draw_solution --ts_max_time 1. -mymonitor 39c4762a1bSJed Brown 40c4762a1bSJed Brown ./biharmonic -ts_monitor -snes_monitor -pc_type lu --snes_converged_reason -draw_pause -2 -ts_type cn -da_refine 5 -kappa .0001 -ts_dt 5.96046e-06 -cahn-hillard -energy 3 -snes_linesearch_monitor -theta .00000001 -ts_monitor_draw_solution --ts_max_time 1. -degenerate -mymonitor 41c4762a1bSJed Brown 42c4762a1bSJed Brown Evolve the Cahn-Hillard equations: logarithmic + double obstacle (never shrinks, never grows) 43c4762a1bSJed Brown --------------- 44c4762a1bSJed Brown ./biharmonic -ts_monitor -snes_monitor -pc_type lu --snes_converged_reason -draw_pause -2 -ts_type cn -da_refine 5 -kappa .00001 -ts_dt 5.96046e-06 -cahn-hillard -energy 4 -snes_linesearch_monitor -theta .00000001 -ts_monitor_draw_solution --mymonitor 45c4762a1bSJed Brown 46c4762a1bSJed Brown */ 47c4762a1bSJed Brown #include <petscdm.h> 48c4762a1bSJed Brown #include <petscdmda.h> 49c4762a1bSJed Brown #include <petscts.h> 50c4762a1bSJed Brown #include <petscdraw.h> 51c4762a1bSJed Brown 52c4762a1bSJed Brown extern PetscErrorCode FormFunction(TS,PetscReal,Vec,Vec,void*),FormInitialSolution(DM,Vec),MyMonitor(TS,PetscInt,PetscReal,Vec,void*),MyDestroy(void**),FormJacobian(TS,PetscReal,Vec,Mat,Mat,void*); 53c4762a1bSJed Brown typedef struct {PetscBool cahnhillard;PetscBool degenerate;PetscReal kappa;PetscInt energy;PetscReal tol;PetscReal theta,theta_c;PetscInt truncation;PetscBool netforce; PetscDrawViewPorts *ports;} UserCtx; 54c4762a1bSJed Brown 55c4762a1bSJed Brown int main(int argc,char **argv) 56c4762a1bSJed Brown { 57c4762a1bSJed Brown TS ts; /* nonlinear solver */ 58c4762a1bSJed Brown Vec x,r; /* solution, residual vectors */ 59c4762a1bSJed Brown Mat J; /* Jacobian matrix */ 60c4762a1bSJed Brown PetscInt steps,Mx; 61c4762a1bSJed Brown DM da; 62c4762a1bSJed Brown PetscReal dt; 63c4762a1bSJed Brown PetscBool mymonitor; 64c4762a1bSJed Brown UserCtx ctx; 65c4762a1bSJed Brown 66c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 67c4762a1bSJed Brown Initialize program 68c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 69*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 70c4762a1bSJed Brown ctx.kappa = 1.0; 71*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-kappa",&ctx.kappa,NULL)); 72c4762a1bSJed Brown ctx.degenerate = PETSC_FALSE; 73*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-degenerate",&ctx.degenerate,NULL)); 74c4762a1bSJed Brown ctx.cahnhillard = PETSC_FALSE; 75*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-cahn-hillard",&ctx.cahnhillard,NULL)); 76c4762a1bSJed Brown ctx.netforce = PETSC_FALSE; 77*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-netforce",&ctx.netforce,NULL)); 78c4762a1bSJed Brown ctx.energy = 1; 79*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-energy",&ctx.energy,NULL)); 80c4762a1bSJed Brown ctx.tol = 1.0e-8; 81*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-tol",&ctx.tol,NULL)); 82c4762a1bSJed Brown ctx.theta = .001; 83c4762a1bSJed Brown ctx.theta_c = 1.0; 84*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-theta",&ctx.theta,NULL)); 85*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-theta_c",&ctx.theta_c,NULL)); 86c4762a1bSJed Brown ctx.truncation = 1; 87*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-truncation",&ctx.truncation,NULL)); 88*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-mymonitor",&mymonitor)); 89c4762a1bSJed Brown 90c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 91c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 92c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 93*9566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, 10,1,2,NULL,&da)); 94*9566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 95*9566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 96*9566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da,0,"Biharmonic heat equation: u")); 97*9566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,0,&Mx,0,0,0,0,0,0,0,0,0,0,0)); 98c4762a1bSJed Brown dt = 1.0/(10.*ctx.kappa*Mx*Mx*Mx*Mx); 99c4762a1bSJed Brown 100c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 101c4762a1bSJed Brown Extract global vectors from DMDA; then duplicate for remaining 102c4762a1bSJed Brown vectors that are the same types 103c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 104*9566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da,&x)); 105*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&r)); 106c4762a1bSJed Brown 107c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 108c4762a1bSJed Brown Create timestepping solver context 109c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 110*9566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD,&ts)); 111*9566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts,da)); 112*9566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts,TS_NONLINEAR)); 113*9566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts,NULL,FormFunction,&ctx)); 114*9566063dSJacob Faibussowitsch PetscCall(DMSetMatType(da,MATAIJ)); 115*9566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(da,&J)); 116*9566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(ts,J,J,FormJacobian,&ctx)); 117*9566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts,.02)); 118*9566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_INTERPOLATE)); 119c4762a1bSJed Brown 120c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 121c4762a1bSJed Brown Create matrix data structure; set Jacobian evaluation routine 122c4762a1bSJed Brown 123c4762a1bSJed Brown Set Jacobian matrix data structure and default Jacobian evaluation 124c4762a1bSJed Brown routine. User can override with: 125c4762a1bSJed Brown -snes_mf : matrix-free Newton-Krylov method with no preconditioning 126c4762a1bSJed Brown (unless user explicitly sets preconditioner) 127c4762a1bSJed Brown -snes_mf_operator : form preconditioning matrix as set by the user, 128c4762a1bSJed Brown but use matrix-free approx for Jacobian-vector 129c4762a1bSJed Brown products within Newton-Krylov method 130c4762a1bSJed Brown 131c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 132c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 133c4762a1bSJed Brown Customize nonlinear solver 134c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 135*9566063dSJacob Faibussowitsch PetscCall(TSSetType(ts,TSCN)); 136c4762a1bSJed Brown 137c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 138c4762a1bSJed Brown Set initial conditions 139c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 140*9566063dSJacob Faibussowitsch PetscCall(FormInitialSolution(da,x)); 141*9566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts,dt)); 142*9566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts,x)); 143c4762a1bSJed Brown 144c4762a1bSJed Brown if (mymonitor) { 145c4762a1bSJed Brown ctx.ports = NULL; 146*9566063dSJacob Faibussowitsch PetscCall(TSMonitorSet(ts,MyMonitor,&ctx,MyDestroy)); 147c4762a1bSJed Brown } 148c4762a1bSJed Brown 149c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 150c4762a1bSJed Brown Set runtime options 151c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 152*9566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 153c4762a1bSJed Brown 154c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 155c4762a1bSJed Brown Solve nonlinear system 156c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 157*9566063dSJacob Faibussowitsch PetscCall(TSSolve(ts,x)); 158*9566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts,&steps)); 159*9566063dSJacob Faibussowitsch PetscCall(VecView(x,PETSC_VIEWER_BINARY_WORLD)); 160c4762a1bSJed Brown 161c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 162c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 163c4762a1bSJed Brown are no longer needed. 164c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 165*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 166*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 167*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 168*9566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 169*9566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 170c4762a1bSJed Brown 171*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 172b122ec5aSJacob Faibussowitsch return 0; 173c4762a1bSJed Brown } 174c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 175c4762a1bSJed Brown /* 176c4762a1bSJed Brown FormFunction - Evaluates nonlinear function, F(x). 177c4762a1bSJed Brown 178c4762a1bSJed Brown Input Parameters: 179c4762a1bSJed Brown . ts - the TS context 180c4762a1bSJed Brown . X - input vector 181c4762a1bSJed Brown . ptr - optional user-defined context, as set by SNESSetFunction() 182c4762a1bSJed Brown 183c4762a1bSJed Brown Output Parameter: 184c4762a1bSJed Brown . F - function vector 185c4762a1bSJed Brown */ 186c4762a1bSJed Brown PetscErrorCode FormFunction(TS ts,PetscReal ftime,Vec X,Vec F,void *ptr) 187c4762a1bSJed Brown { 188c4762a1bSJed Brown DM da; 189c4762a1bSJed Brown PetscInt i,Mx,xs,xm; 190c4762a1bSJed Brown PetscReal hx,sx; 191c4762a1bSJed Brown PetscScalar *x,*f,c,r,l; 192c4762a1bSJed Brown Vec localX; 193c4762a1bSJed Brown UserCtx *ctx = (UserCtx*)ptr; 194c4762a1bSJed Brown PetscReal tol = ctx->tol, theta=ctx->theta,theta_c=ctx->theta_c,a,b; /* a and b are used in the cubic truncation of the log function */ 195c4762a1bSJed Brown 196c4762a1bSJed Brown PetscFunctionBegin; 197*9566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&da)); 198*9566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da,&localX)); 199*9566063dSJacob 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)); 200c4762a1bSJed Brown 201c4762a1bSJed Brown hx = 1.0/(PetscReal)Mx; sx = 1.0/(hx*hx); 202c4762a1bSJed Brown 203c4762a1bSJed Brown /* 204c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process 205c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 206c4762a1bSJed Brown By placing code between these two statements, computations can be 207c4762a1bSJed Brown done while messages are in transition. 208c4762a1bSJed Brown */ 209*9566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX)); 210*9566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX)); 211c4762a1bSJed Brown 212c4762a1bSJed Brown /* 213c4762a1bSJed Brown Get pointers to vector data 214c4762a1bSJed Brown */ 215*9566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da,localX,&x)); 216*9566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,F,&f)); 217c4762a1bSJed Brown 218c4762a1bSJed Brown /* 219c4762a1bSJed Brown Get local grid boundaries 220c4762a1bSJed Brown */ 221*9566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL)); 222c4762a1bSJed Brown 223c4762a1bSJed Brown /* 224c4762a1bSJed Brown Compute function over the locally owned part of the grid 225c4762a1bSJed Brown */ 226c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 227c4762a1bSJed Brown if (ctx->degenerate) { 228c4762a1bSJed Brown c = (1. - x[i]*x[i])*(x[i-1] + x[i+1] - 2.0*x[i])*sx; 229c4762a1bSJed Brown r = (1. - x[i+1]*x[i+1])*(x[i] + x[i+2] - 2.0*x[i+1])*sx; 230c4762a1bSJed Brown l = (1. - x[i-1]*x[i-1])*(x[i-2] + x[i] - 2.0*x[i-1])*sx; 231c4762a1bSJed Brown } else { 232c4762a1bSJed Brown c = (x[i-1] + x[i+1] - 2.0*x[i])*sx; 233c4762a1bSJed Brown r = (x[i] + x[i+2] - 2.0*x[i+1])*sx; 234c4762a1bSJed Brown l = (x[i-2] + x[i] - 2.0*x[i-1])*sx; 235c4762a1bSJed Brown } 236c4762a1bSJed Brown f[i] = -ctx->kappa*(l + r - 2.0*c)*sx; 237c4762a1bSJed Brown if (ctx->cahnhillard) { 238c4762a1bSJed Brown switch (ctx->energy) { 239c4762a1bSJed Brown case 1: /* double well */ 240c4762a1bSJed Brown f[i] += 6.*.25*x[i]*(x[i+1] - x[i-1])*(x[i+1] - x[i-1])*sx + (3.*x[i]*x[i] - 1.)*(x[i-1] + x[i+1] - 2.0*x[i])*sx; 241c4762a1bSJed Brown break; 242c4762a1bSJed Brown case 2: /* double obstacle */ 243c4762a1bSJed Brown f[i] += -(x[i-1] + x[i+1] - 2.0*x[i])*sx; 244c4762a1bSJed Brown break; 245c4762a1bSJed Brown case 3: /* logarithmic + double well */ 246c4762a1bSJed Brown f[i] += 6.*.25*x[i]*(x[i+1] - x[i-1])*(x[i+1] - x[i-1])*sx + (3.*x[i]*x[i] - 1.)*(x[i-1] + x[i+1] - 2.0*x[i])*sx; 247c4762a1bSJed Brown if (ctx->truncation==2) { /* log function with approximated with a quadratic polynomial outside -1.0+2*tol, 1.0-2*tol */ 248c4762a1bSJed Brown if (PetscRealPart(x[i]) < -1.0 + 2.0*tol) f[i] += (.25*theta/(tol-tol*tol))*(x[i-1] + x[i+1] - 2.0*x[i])*sx; 249c4762a1bSJed Brown else if (PetscRealPart(x[i]) > 1.0 - 2.0*tol) f[i] += (.25*theta/(tol-tol*tol))*(x[i-1] + x[i+1] - 2.0*x[i])*sx; 250c4762a1bSJed Brown else f[i] += 2.0*theta*x[i]/((1.0-x[i]*x[i])*(1.0-x[i]*x[i]))*.25*(x[i+1] - x[i-1])*(x[i+1] - x[i-1])*sx + (theta/(1.0-x[i]*x[i]))*(x[i-1] + x[i+1] - 2.0*x[i])*sx; 251c4762a1bSJed Brown } else { /* log function is approximated with a cubic polynomial outside -1.0+2*tol, 1.0-2*tol */ 252c4762a1bSJed Brown a = 2.0*theta*(1.0-2.0*tol)/(16.0*tol*tol*(1.0-tol)*(1.0-tol)); 253c4762a1bSJed Brown b = theta/(4.0*tol*(1.0-tol)) - a*(1.0-2.0*tol); 254c4762a1bSJed Brown if (PetscRealPart(x[i]) < -1.0 + 2.0*tol) f[i] += -1.0*a*.25*(x[i+1] - x[i-1])*(x[i+1] - x[i-1])*sx + (-1.0*a*x[i] + b)*(x[i-1] + x[i+1] - 2.0*x[i])*sx; 255c4762a1bSJed Brown else if (PetscRealPart(x[i]) > 1.0 - 2.0*tol) f[i] += 1.0*a*.25*(x[i+1] - x[i-1])*(x[i+1] - x[i-1])*sx + ( a*x[i] + b)*(x[i-1] + x[i+1] - 2.0*x[i])*sx; 256c4762a1bSJed Brown else f[i] += 2.0*theta*x[i]/((1.0-x[i]*x[i])*(1.0-x[i]*x[i]))*.25*(x[i+1] - x[i-1])*(x[i+1] - x[i-1])*sx + (theta/(1.0-x[i]*x[i]))*(x[i-1] + x[i+1] - 2.0*x[i])*sx; 257c4762a1bSJed Brown } 258c4762a1bSJed Brown break; 259c4762a1bSJed Brown case 4: /* logarithmic + double obstacle */ 260c4762a1bSJed Brown f[i] += -theta_c*(x[i-1] + x[i+1] - 2.0*x[i])*sx; 261c4762a1bSJed Brown if (ctx->truncation==2) { /* quadratic */ 262c4762a1bSJed Brown if (PetscRealPart(x[i]) < -1.0 + 2.0*tol) f[i] += (.25*theta/(tol-tol*tol))*(x[i-1] + x[i+1] - 2.0*x[i])*sx; 263c4762a1bSJed Brown else if (PetscRealPart(x[i]) > 1.0 - 2.0*tol) f[i] += (.25*theta/(tol-tol*tol))*(x[i-1] + x[i+1] - 2.0*x[i])*sx; 264c4762a1bSJed Brown else f[i] += 2.0*theta*x[i]/((1.0-x[i]*x[i])*(1.0-x[i]*x[i]))*.25*(x[i+1] - x[i-1])*(x[i+1] - x[i-1])*sx + (theta/(1.0-x[i]*x[i]))*(x[i-1] + x[i+1] - 2.0*x[i])*sx; 265c4762a1bSJed Brown } else { /* cubic */ 266c4762a1bSJed Brown a = 2.0*theta*(1.0-2.0*tol)/(16.0*tol*tol*(1.0-tol)*(1.0-tol)); 267c4762a1bSJed Brown b = theta/(4.0*tol*(1.0-tol)) - a*(1.0-2.0*tol); 268c4762a1bSJed Brown if (PetscRealPart(x[i]) < -1.0 + 2.0*tol) f[i] += -1.0*a*.25*(x[i+1] - x[i-1])*(x[i+1] - x[i-1])*sx + (-1.0*a*x[i] + b)*(x[i-1] + x[i+1] - 2.0*x[i])*sx; 269c4762a1bSJed Brown else if (PetscRealPart(x[i]) > 1.0 - 2.0*tol) f[i] += 1.0*a*.25*(x[i+1] - x[i-1])*(x[i+1] - x[i-1])*sx + ( a*x[i] + b)*(x[i-1] + x[i+1] - 2.0*x[i])*sx; 270c4762a1bSJed Brown else f[i] += 2.0*theta*x[i]/((1.0-x[i]*x[i])*(1.0-x[i]*x[i]))*.25*(x[i+1] - x[i-1])*(x[i+1] - x[i-1])*sx + (theta/(1.0-x[i]*x[i]))*(x[i-1] + x[i+1] - 2.0*x[i])*sx; 271c4762a1bSJed Brown } 272c4762a1bSJed Brown break; 273c4762a1bSJed Brown } 274c4762a1bSJed Brown } 275c4762a1bSJed Brown 276c4762a1bSJed Brown } 277c4762a1bSJed Brown 278c4762a1bSJed Brown /* 279c4762a1bSJed Brown Restore vectors 280c4762a1bSJed Brown */ 281*9566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da,localX,&x)); 282*9566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,F,&f)); 283*9566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da,&localX)); 284c4762a1bSJed Brown PetscFunctionReturn(0); 285c4762a1bSJed Brown } 286c4762a1bSJed Brown 287c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 288c4762a1bSJed Brown /* 289c4762a1bSJed Brown FormJacobian - Evaluates nonlinear function's Jacobian 290c4762a1bSJed Brown 291c4762a1bSJed Brown */ 292c4762a1bSJed Brown PetscErrorCode FormJacobian(TS ts,PetscReal ftime,Vec X,Mat A,Mat B,void *ptr) 293c4762a1bSJed Brown { 294c4762a1bSJed Brown DM da; 295c4762a1bSJed Brown PetscInt i,Mx,xs,xm; 296c4762a1bSJed Brown MatStencil row,cols[5]; 297c4762a1bSJed Brown PetscReal hx,sx; 298c4762a1bSJed Brown PetscScalar *x,vals[5]; 299c4762a1bSJed Brown Vec localX; 300c4762a1bSJed Brown UserCtx *ctx = (UserCtx*)ptr; 301c4762a1bSJed Brown 302c4762a1bSJed Brown PetscFunctionBegin; 303*9566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&da)); 304*9566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da,&localX)); 305*9566063dSJacob 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)); 306c4762a1bSJed Brown 307c4762a1bSJed Brown hx = 1.0/(PetscReal)Mx; sx = 1.0/(hx*hx); 308c4762a1bSJed Brown 309c4762a1bSJed Brown /* 310c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process 311c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 312c4762a1bSJed Brown By placing code between these two statements, computations can be 313c4762a1bSJed Brown done while messages are in transition. 314c4762a1bSJed Brown */ 315*9566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX)); 316*9566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX)); 317c4762a1bSJed Brown 318c4762a1bSJed Brown /* 319c4762a1bSJed Brown Get pointers to vector data 320c4762a1bSJed Brown */ 321*9566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da,localX,&x)); 322c4762a1bSJed Brown 323c4762a1bSJed Brown /* 324c4762a1bSJed Brown Get local grid boundaries 325c4762a1bSJed Brown */ 326*9566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL)); 327c4762a1bSJed Brown 328c4762a1bSJed Brown /* 329c4762a1bSJed Brown Compute function over the locally owned part of the grid 330c4762a1bSJed Brown */ 331c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 332c4762a1bSJed Brown row.i = i; 333c4762a1bSJed Brown if (ctx->degenerate) { 334c4762a1bSJed Brown /*PetscScalar c,r,l; 335c4762a1bSJed Brown c = (1. - x[i]*x[i])*(x[i-1] + x[i+1] - 2.0*x[i])*sx; 336c4762a1bSJed Brown r = (1. - x[i+1]*x[i+1])*(x[i] + x[i+2] - 2.0*x[i+1])*sx; 337c4762a1bSJed Brown l = (1. - x[i-1]*x[i-1])*(x[i-2] + x[i] - 2.0*x[i-1])*sx; */ 338c4762a1bSJed Brown } else { 339c4762a1bSJed Brown cols[0].i = i - 2; vals[0] = -ctx->kappa*sx*sx; 340c4762a1bSJed Brown cols[1].i = i - 1; vals[1] = 4.0*ctx->kappa*sx*sx; 341c4762a1bSJed Brown cols[2].i = i ; vals[2] = -6.0*ctx->kappa*sx*sx; 342c4762a1bSJed Brown cols[3].i = i + 1; vals[3] = 4.0*ctx->kappa*sx*sx; 343c4762a1bSJed Brown cols[4].i = i + 2; vals[4] = -ctx->kappa*sx*sx; 344c4762a1bSJed Brown } 345*9566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(B,1,&row,5,cols,vals,INSERT_VALUES)); 346c4762a1bSJed Brown 347c4762a1bSJed Brown if (ctx->cahnhillard) { 348c4762a1bSJed Brown switch (ctx->energy) { 349c4762a1bSJed Brown case 1: /* double well */ 350c4762a1bSJed Brown /* f[i] += 6.*.25*x[i]*(x[i+1] - x[i-1])*(x[i+1] - x[i-1])*sx + (3.*x[i]*x[i] - 1.)*(x[i-1] + x[i+1] - 2.0*x[i])*sx; */ 351c4762a1bSJed Brown break; 352c4762a1bSJed Brown case 2: /* double obstacle */ 353c4762a1bSJed Brown /* f[i] += -(x[i-1] + x[i+1] - 2.0*x[i])*sx; */ 354c4762a1bSJed Brown break; 355c4762a1bSJed Brown case 3: /* logarithmic + double well */ 356c4762a1bSJed Brown break; 357c4762a1bSJed Brown case 4: /* logarithmic + double obstacle */ 358c4762a1bSJed Brown break; 359c4762a1bSJed Brown } 360c4762a1bSJed Brown } 361c4762a1bSJed Brown 362c4762a1bSJed Brown } 363c4762a1bSJed Brown 364c4762a1bSJed Brown /* 365c4762a1bSJed Brown Restore vectors 366c4762a1bSJed Brown */ 367*9566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da,localX,&x)); 368*9566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da,&localX)); 369*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 370*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 371c4762a1bSJed Brown if (A != B) { 372*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 373*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 374c4762a1bSJed Brown } 375c4762a1bSJed Brown PetscFunctionReturn(0); 376c4762a1bSJed Brown } 377c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 378c4762a1bSJed Brown PetscErrorCode FormInitialSolution(DM da,Vec U) 379c4762a1bSJed Brown { 380c4762a1bSJed Brown PetscInt i,xs,xm,Mx,N,scale; 381c4762a1bSJed Brown PetscScalar *u; 382c4762a1bSJed Brown PetscReal r,hx,x; 383c4762a1bSJed Brown const PetscScalar *f; 384c4762a1bSJed Brown Vec finesolution; 385c4762a1bSJed Brown PetscViewer viewer; 386c4762a1bSJed Brown 387c4762a1bSJed Brown PetscFunctionBegin; 388*9566063dSJacob 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)); 389c4762a1bSJed Brown 390c4762a1bSJed Brown hx = 1.0/(PetscReal)Mx; 391c4762a1bSJed Brown 392c4762a1bSJed Brown /* 393c4762a1bSJed Brown Get pointers to vector data 394c4762a1bSJed Brown */ 395*9566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da,U,&u)); 396c4762a1bSJed Brown 397c4762a1bSJed Brown /* 398c4762a1bSJed Brown Get local grid boundaries 399c4762a1bSJed Brown */ 400*9566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL)); 401c4762a1bSJed Brown 402c4762a1bSJed Brown /* 403c4762a1bSJed Brown Seee heat.c for how to generate InitialSolution.heat 404c4762a1bSJed Brown */ 405*9566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"InitialSolution.heat",FILE_MODE_READ,&viewer)); 406*9566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&finesolution)); 407*9566063dSJacob Faibussowitsch PetscCall(VecLoad(finesolution,viewer)); 408*9566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 409*9566063dSJacob Faibussowitsch PetscCall(VecGetSize(finesolution,&N)); 410c4762a1bSJed Brown scale = N/Mx; 411*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(finesolution,&f)); 412c4762a1bSJed Brown 413c4762a1bSJed Brown /* 414c4762a1bSJed Brown Compute function over the locally owned part of the grid 415c4762a1bSJed Brown */ 416c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 417c4762a1bSJed Brown x = i*hx; 418c4762a1bSJed Brown r = PetscSqrtReal((x-.5)*(x-.5)); 419c4762a1bSJed Brown if (r < .125) u[i] = 1.0; 420c4762a1bSJed Brown else u[i] = -.5; 421c4762a1bSJed Brown 422c4762a1bSJed Brown /* With the initial condition above the method is first order in space */ 423c4762a1bSJed Brown /* this is a smooth initial condition so the method becomes second order in space */ 424c4762a1bSJed Brown /*u[i] = PetscSinScalar(2*PETSC_PI*x); */ 425c4762a1bSJed Brown u[i] = f[scale*i]; 426c4762a1bSJed Brown } 427*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(finesolution,&f)); 428*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&finesolution)); 429c4762a1bSJed Brown 430c4762a1bSJed Brown /* 431c4762a1bSJed Brown Restore vectors 432c4762a1bSJed Brown */ 433*9566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da,U,&u)); 434c4762a1bSJed Brown PetscFunctionReturn(0); 435c4762a1bSJed Brown } 436c4762a1bSJed Brown 437c4762a1bSJed Brown /* 438c4762a1bSJed Brown This routine is not parallel 439c4762a1bSJed Brown */ 440c4762a1bSJed Brown PetscErrorCode MyMonitor(TS ts,PetscInt step,PetscReal time,Vec U,void *ptr) 441c4762a1bSJed Brown { 442c4762a1bSJed Brown UserCtx *ctx = (UserCtx*)ptr; 443c4762a1bSJed Brown PetscDrawLG lg; 444c4762a1bSJed Brown PetscScalar *u,l,r,c; 445c4762a1bSJed Brown PetscInt Mx,i,xs,xm,cnt; 446c4762a1bSJed Brown PetscReal x,y,hx,pause,sx,len,max,xx[4],yy[4],xx_netforce,yy_netforce,yup,ydown,y2,len2; 447c4762a1bSJed Brown PetscDraw draw; 448c4762a1bSJed Brown Vec localU; 449c4762a1bSJed Brown DM da; 450c4762a1bSJed Brown int colors[] = {PETSC_DRAW_YELLOW,PETSC_DRAW_RED,PETSC_DRAW_BLUE,PETSC_DRAW_PLUM,PETSC_DRAW_BLACK}; 451c4762a1bSJed Brown /* 452c4762a1bSJed Brown const char *const legend[3][3] = {{"-kappa (\\grad u,\\grad u)","(1 - u^2)^2"},{"-kappa (\\grad u,\\grad u)","(1 - u^2)"},{"-kappa (\\grad u,\\grad u)","logarithmic"}}; 453c4762a1bSJed Brown */ 454c4762a1bSJed Brown PetscDrawAxis axis; 455c4762a1bSJed Brown PetscDrawViewPorts *ports; 456c4762a1bSJed Brown PetscReal tol = ctx->tol, theta=ctx->theta,theta_c=ctx->theta_c,a,b; /* a and b are used in the cubic truncation of the log function */ 457c4762a1bSJed Brown PetscReal vbounds[] = {-1.1,1.1}; 458c4762a1bSJed Brown 459c4762a1bSJed Brown PetscFunctionBegin; 460*9566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawSetBounds(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),1,vbounds)); 461*9566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawResize(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),800,600)); 462*9566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&da)); 463*9566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da,&localU)); 464*9566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE, 465b122ec5aSJacob Faibussowitsch PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE)); 466*9566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL)); 467c4762a1bSJed Brown hx = 1.0/(PetscReal)Mx; sx = 1.0/(hx*hx); 468*9566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU)); 469*9566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU)); 470*9566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da,localU,&u)); 471c4762a1bSJed Brown 472*9566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDrawLG(PETSC_VIEWER_DRAW_(PETSC_COMM_WORLD),1,&lg)); 473*9566063dSJacob Faibussowitsch PetscCall(PetscDrawLGGetDraw(lg,&draw)); 474*9566063dSJacob Faibussowitsch PetscCall(PetscDrawCheckResizedWindow(draw)); 475c4762a1bSJed Brown if (!ctx->ports) { 476*9566063dSJacob Faibussowitsch PetscCall(PetscDrawViewPortsCreateRect(draw,1,3,&ctx->ports)); 477c4762a1bSJed Brown } 478c4762a1bSJed Brown ports = ctx->ports; 479*9566063dSJacob Faibussowitsch PetscCall(PetscDrawLGGetAxis(lg,&axis)); 480*9566063dSJacob Faibussowitsch PetscCall(PetscDrawLGReset(lg)); 481c4762a1bSJed Brown 482c4762a1bSJed Brown xx[0] = 0.0; xx[1] = 1.0; cnt = 2; 483*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetRealArray(NULL,NULL,"-zoom",xx,&cnt,NULL)); 484c4762a1bSJed Brown xs = xx[0]/hx; xm = (xx[1] - xx[0])/hx; 485c4762a1bSJed Brown 486c4762a1bSJed Brown /* 487c4762a1bSJed Brown Plot the energies 488c4762a1bSJed Brown */ 489*9566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetDimension(lg,1 + (ctx->cahnhillard ? 1 : 0) + (ctx->energy == 3))); 490*9566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetColors(lg,colors+1)); 491*9566063dSJacob Faibussowitsch PetscCall(PetscDrawViewPortsSet(ports,2)); 492c4762a1bSJed Brown x = hx*xs; 493c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 494c4762a1bSJed Brown xx[0] = xx[1] = xx[2] = x; 495c4762a1bSJed Brown if (ctx->degenerate) yy[0] = PetscRealPart(.25*(1. - u[i]*u[i])*ctx->kappa*(u[i-1] - u[i+1])*(u[i-1] - u[i+1])*sx); 496c4762a1bSJed Brown else yy[0] = PetscRealPart(.25*ctx->kappa*(u[i-1] - u[i+1])*(u[i-1] - u[i+1])*sx); 497c4762a1bSJed Brown 498c4762a1bSJed Brown if (ctx->cahnhillard) { 499c4762a1bSJed Brown switch (ctx->energy) { 500c4762a1bSJed Brown case 1: /* double well */ 501c4762a1bSJed Brown yy[1] = .25*PetscRealPart((1. - u[i]*u[i])*(1. - u[i]*u[i])); 502c4762a1bSJed Brown break; 503c4762a1bSJed Brown case 2: /* double obstacle */ 504c4762a1bSJed Brown yy[1] = .5*PetscRealPart(1. - u[i]*u[i]); 505c4762a1bSJed Brown break; 506c4762a1bSJed Brown case 3: /* logarithm + double well */ 507c4762a1bSJed Brown yy[1] = .25*PetscRealPart((1. - u[i]*u[i])*(1. - u[i]*u[i])); 508c4762a1bSJed Brown if (PetscRealPart(u[i]) < -1.0 + 2.0*tol) yy[2] = .5*theta*(2.0*tol*PetscLogReal(tol) + PetscRealPart(1.0-u[i])*PetscLogReal(PetscRealPart(1.-u[i])/2.0)); 509c4762a1bSJed Brown else if (PetscRealPart(u[i]) > 1.0 - 2.0*tol) yy[2] = .5*theta*(PetscRealPart(1.0+u[i])*PetscLogReal(PetscRealPart(1.0+u[i])/2.0) + 2.0*tol*PetscLogReal(tol)); 510c4762a1bSJed Brown else yy[2] = .5*theta*(PetscRealPart(1.0+u[i])*PetscLogReal(PetscRealPart(1.0+u[i])/2.0) + PetscRealPart(1.0-u[i])*PetscLogReal(PetscRealPart(1.0-u[i])/2.0)); 511c4762a1bSJed Brown break; 512c4762a1bSJed Brown case 4: /* logarithm + double obstacle */ 513c4762a1bSJed Brown yy[1] = .5*theta_c*PetscRealPart(1.0-u[i]*u[i]); 514c4762a1bSJed Brown if (PetscRealPart(u[i]) < -1.0 + 2.0*tol) yy[2] = .5*theta*(2.0*tol*PetscLogReal(tol) + PetscRealPart(1.0-u[i])*PetscLogReal(PetscRealPart(1.-u[i])/2.0)); 515c4762a1bSJed Brown else if (PetscRealPart(u[i]) > 1.0 - 2.0*tol) yy[2] = .5*theta*(PetscRealPart(1.0+u[i])*PetscLogReal(PetscRealPart(1.0+u[i])/2.0) + 2.0*tol*PetscLogReal(tol)); 516c4762a1bSJed Brown else yy[2] = .5*theta*(PetscRealPart(1.0+u[i])*PetscLogReal(PetscRealPart(1.0+u[i])/2.0) + PetscRealPart(1.0-u[i])*PetscLogReal(PetscRealPart(1.0-u[i])/2.0)); 517c4762a1bSJed Brown break; 518c4762a1bSJed Brown default: 519c4762a1bSJed Brown SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"It will always be one of the values"); 520c4762a1bSJed Brown } 521c4762a1bSJed Brown } 522*9566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddPoint(lg,xx,yy)); 523c4762a1bSJed Brown x += hx; 524c4762a1bSJed Brown } 525*9566063dSJacob Faibussowitsch PetscCall(PetscDrawGetPause(draw,&pause)); 526*9566063dSJacob Faibussowitsch PetscCall(PetscDrawSetPause(draw,0.0)); 527*9566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis,"Energy","","")); 528*9566063dSJacob Faibussowitsch /* PetscCall(PetscDrawLGSetLegend(lg,legend[ctx->energy-1])); */ 529*9566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDraw(lg)); 530c4762a1bSJed Brown 531c4762a1bSJed Brown /* 532c4762a1bSJed Brown Plot the forces 533c4762a1bSJed Brown */ 534*9566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetDimension(lg,0 + (ctx->cahnhillard ? 2 : 0) + (ctx->energy == 3))); 535*9566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetColors(lg,colors+1)); 536*9566063dSJacob Faibussowitsch PetscCall(PetscDrawViewPortsSet(ports,1)); 537*9566063dSJacob Faibussowitsch PetscCall(PetscDrawLGReset(lg)); 538c4762a1bSJed Brown x = xs*hx; 539c4762a1bSJed Brown max = 0.; 540c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 541c4762a1bSJed Brown xx[0] = xx[1] = xx[2] = xx[3] = x; 542c4762a1bSJed Brown xx_netforce = x; 543c4762a1bSJed Brown if (ctx->degenerate) { 544c4762a1bSJed Brown c = (1. - u[i]*u[i])*(u[i-1] + u[i+1] - 2.0*u[i])*sx; 545c4762a1bSJed Brown r = (1. - u[i+1]*u[i+1])*(u[i] + u[i+2] - 2.0*u[i+1])*sx; 546c4762a1bSJed Brown l = (1. - u[i-1]*u[i-1])*(u[i-2] + u[i] - 2.0*u[i-1])*sx; 547c4762a1bSJed Brown } else { 548c4762a1bSJed Brown c = (u[i-1] + u[i+1] - 2.0*u[i])*sx; 549c4762a1bSJed Brown r = (u[i] + u[i+2] - 2.0*u[i+1])*sx; 550c4762a1bSJed Brown l = (u[i-2] + u[i] - 2.0*u[i-1])*sx; 551c4762a1bSJed Brown } 552c4762a1bSJed Brown yy[0] = PetscRealPart(-ctx->kappa*(l + r - 2.0*c)*sx); 553c4762a1bSJed Brown yy_netforce = yy[0]; 554c4762a1bSJed Brown max = PetscMax(max,PetscAbs(yy[0])); 555c4762a1bSJed Brown if (ctx->cahnhillard) { 556c4762a1bSJed Brown switch (ctx->energy) { 557c4762a1bSJed Brown case 1: /* double well */ 558c4762a1bSJed Brown yy[1] = PetscRealPart(6.*.25*u[i]*(u[i+1] - u[i-1])*(u[i+1] - u[i-1])*sx + (3.*u[i]*u[i] - 1.)*(u[i-1] + u[i+1] - 2.0*u[i])*sx); 559c4762a1bSJed Brown break; 560c4762a1bSJed Brown case 2: /* double obstacle */ 561c4762a1bSJed Brown yy[1] = -PetscRealPart(u[i-1] + u[i+1] - 2.0*u[i])*sx; 562c4762a1bSJed Brown break; 563c4762a1bSJed Brown case 3: /* logarithmic + double well */ 564c4762a1bSJed Brown yy[1] = PetscRealPart(6.*.25*u[i]*(u[i+1] - u[i-1])*(u[i+1] - u[i-1])*sx + (3.*u[i]*u[i] - 1.)*(u[i-1] + u[i+1] - 2.0*u[i])*sx); 565c4762a1bSJed Brown if (ctx->truncation==2) { /* quadratic */ 566c4762a1bSJed Brown if (PetscRealPart(u[i]) < -1.0 + 2.0*tol) yy[2] = (.25*theta/(tol-tol*tol))*PetscRealPart(u[i-1] + u[i+1] - 2.0*u[i])*sx; 567c4762a1bSJed Brown else if (PetscRealPart(u[i]) > 1.0 - 2.0*tol) yy[2] = (.25*theta/(tol-tol*tol))*PetscRealPart(u[i-1] + u[i+1] - 2.0*u[i])*sx; 568c4762a1bSJed Brown else yy[2] = PetscRealPart(2.0*theta*u[i]/((1.0-u[i]*u[i])*(1.0-u[i]*u[i]))*.25*(u[i+1] - u[i-1])*(u[i+1] - u[i-1])*sx + (theta/(1.0-u[i]*u[i]))*(u[i-1] + u[i+1] - 2.0*u[i])*sx); 569c4762a1bSJed Brown } else { /* cubic */ 570c4762a1bSJed Brown a = 2.0*theta*(1.0-2.0*tol)/(16.0*tol*tol*(1.0-tol)*(1.0-tol)); 571c4762a1bSJed Brown b = theta/(4.0*tol*(1.0-tol)) - a*(1.0-2.0*tol); 572c4762a1bSJed Brown if (PetscRealPart(u[i]) < -1.0 + 2.0*tol) yy[2] = PetscRealPart(-1.0*a*.25*(u[i+1] - u[i-1])*(u[i+1] - u[i-1])*sx + (-1.0*a*u[i] + b)*(u[i-1] + u[i+1] - 2.0*u[i])*sx); 573c4762a1bSJed Brown else if (PetscRealPart(u[i]) > 1.0 - 2.0*tol) yy[2] = PetscRealPart(1.0*a*.25*(u[i+1] - u[i-1])*(u[i+1] - u[i-1])*sx + ( a*u[i] + b)*(u[i-1] + u[i+1] - 2.0*u[i])*sx); 574c4762a1bSJed Brown else yy[2] = PetscRealPart(2.0*theta*u[i]/((1.0-u[i]*u[i])*(1.0-u[i]*u[i]))*.25*(u[i+1] - u[i-1])*(u[i+1] - u[i-1])*sx + (theta/(1.0-u[i]*u[i]))*(u[i-1] + u[i+1] - 2.0*u[i])*sx); 575c4762a1bSJed Brown } 576c4762a1bSJed Brown break; 577c4762a1bSJed Brown case 4: /* logarithmic + double obstacle */ 578c4762a1bSJed Brown yy[1] = theta_c*PetscRealPart(-(u[i-1] + u[i+1] - 2.0*u[i]))*sx; 579c4762a1bSJed Brown if (ctx->truncation==2) { 580c4762a1bSJed Brown if (PetscRealPart(u[i]) < -1.0 + 2.0*tol) yy[2] = (.25*theta/(tol-tol*tol))*PetscRealPart(u[i-1] + u[i+1] - 2.0*u[i])*sx; 581c4762a1bSJed Brown else if (PetscRealPart(u[i]) > 1.0 - 2.0*tol) yy[2] = (.25*theta/(tol-tol*tol))*PetscRealPart(u[i-1] + u[i+1] - 2.0*u[i])*sx; 582c4762a1bSJed Brown else yy[2] = PetscRealPart(2.0*theta*u[i]/((1.0-u[i]*u[i])*(1.0-u[i]*u[i]))*.25*(u[i+1] - u[i-1])*(u[i+1] - u[i-1])*sx + (theta/(1.0-u[i]*u[i]))*(u[i-1] + u[i+1] - 2.0*u[i])*sx); 583c4762a1bSJed Brown } 584c4762a1bSJed Brown else { 585c4762a1bSJed Brown a = 2.0*theta*(1.0-2.0*tol)/(16.0*tol*tol*(1.0-tol)*(1.0-tol)); 586c4762a1bSJed Brown b = theta/(4.0*tol*(1.0-tol)) - a*(1.0-2.0*tol); 587c4762a1bSJed Brown if (PetscRealPart(u[i]) < -1.0 + 2.0*tol) yy[2] = PetscRealPart(-1.0*a*.25*(u[i+1] - u[i-1])*(u[i+1] - u[i-1])*sx + (-1.0*a*u[i] + b)*(u[i-1] + u[i+1] - 2.0*u[i])*sx); 588c4762a1bSJed Brown else if (PetscRealPart(u[i]) > 1.0 - 2.0*tol) yy[2] = PetscRealPart(1.0*a*.25*(u[i+1] - u[i-1])*(u[i+1] - u[i-1])*sx + ( a*u[i] + b)*(u[i-1] + u[i+1] - 2.0*u[i])*sx); 589c4762a1bSJed Brown else yy[2] = PetscRealPart(2.0*theta*u[i]/((1.0-u[i]*u[i])*(1.0-u[i]*u[i]))*.25*(u[i+1] - u[i-1])*(u[i+1] - u[i-1])*sx + (theta/(1.0-u[i]*u[i]))*(u[i-1] + u[i+1] - 2.0*u[i])*sx); 590c4762a1bSJed Brown } 591c4762a1bSJed Brown break; 592c4762a1bSJed Brown default: 593c4762a1bSJed Brown SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"It will always be one of the values"); 594c4762a1bSJed Brown } 595c4762a1bSJed Brown if (ctx->energy < 3) { 596c4762a1bSJed Brown max = PetscMax(max,PetscAbs(yy[1])); 597c4762a1bSJed Brown yy[2] = yy[0]+yy[1]; 598c4762a1bSJed Brown yy_netforce = yy[2]; 599c4762a1bSJed Brown } else { 600c4762a1bSJed Brown max = PetscMax(max,PetscAbs(yy[1]+yy[2])); 601c4762a1bSJed Brown yy[3] = yy[0]+yy[1]+yy[2]; 602c4762a1bSJed Brown yy_netforce = yy[3]; 603c4762a1bSJed Brown } 604c4762a1bSJed Brown } 605c4762a1bSJed Brown if (ctx->netforce) { 606*9566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddPoint(lg,&xx_netforce,&yy_netforce)); 607c4762a1bSJed Brown } else { 608*9566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddPoint(lg,xx,yy)); 609c4762a1bSJed Brown } 610c4762a1bSJed Brown x += hx; 611c4762a1bSJed Brown /*if (max > 7200150000.0) */ 612c4762a1bSJed Brown /* printf("max very big when i = %d\n",i); */ 613c4762a1bSJed Brown } 614*9566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis,"Right hand side","","")); 615*9566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetLegend(lg,NULL)); 616*9566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDraw(lg)); 617c4762a1bSJed Brown 618c4762a1bSJed Brown /* 619c4762a1bSJed Brown Plot the solution 620c4762a1bSJed Brown */ 621*9566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetDimension(lg,1)); 622*9566063dSJacob Faibussowitsch PetscCall(PetscDrawViewPortsSet(ports,0)); 623*9566063dSJacob Faibussowitsch PetscCall(PetscDrawLGReset(lg)); 624c4762a1bSJed Brown x = hx*xs; 625*9566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetLimits(lg,x,x+(xm-1)*hx,-1.1,1.1)); 626*9566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetColors(lg,colors)); 627c4762a1bSJed Brown for (i=xs; i<xs+xm; i++) { 628c4762a1bSJed Brown xx[0] = x; 629c4762a1bSJed Brown yy[0] = PetscRealPart(u[i]); 630*9566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddPoint(lg,xx,yy)); 631c4762a1bSJed Brown x += hx; 632c4762a1bSJed Brown } 633*9566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis,"Solution","","")); 634*9566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDraw(lg)); 635c4762a1bSJed Brown 636c4762a1bSJed Brown /* 637c4762a1bSJed Brown Print the forces as arrows on the solution 638c4762a1bSJed Brown */ 639c4762a1bSJed Brown x = hx*xs; 640c4762a1bSJed Brown cnt = xm/60; 641c4762a1bSJed Brown cnt = (!cnt) ? 1 : cnt; 642c4762a1bSJed Brown 643c4762a1bSJed Brown for (i=xs; i<xs+xm; i += cnt) { 644c4762a1bSJed Brown y = yup = ydown = PetscRealPart(u[i]); 645c4762a1bSJed Brown c = (u[i-1] + u[i+1] - 2.0*u[i])*sx; 646c4762a1bSJed Brown r = (u[i] + u[i+2] - 2.0*u[i+1])*sx; 647c4762a1bSJed Brown l = (u[i-2] + u[i] - 2.0*u[i-1])*sx; 648c4762a1bSJed Brown len = -.5*PetscRealPart(ctx->kappa*(l + r - 2.0*c)*sx)/max; 649*9566063dSJacob Faibussowitsch PetscCall(PetscDrawArrow(draw,x,y,x,y+len,PETSC_DRAW_RED)); 650c4762a1bSJed Brown if (ctx->cahnhillard) { 651c4762a1bSJed Brown if (len < 0.) ydown += len; 652c4762a1bSJed Brown else yup += len; 653c4762a1bSJed Brown 654c4762a1bSJed Brown switch (ctx->energy) { 655c4762a1bSJed Brown case 1: /* double well */ 656c4762a1bSJed Brown len = .5*PetscRealPart(6.*.25*u[i]*(u[i+1] - u[i-1])*(u[i+1] - u[i-1])*sx + (3.*u[i]*u[i] - 1.)*(u[i-1] + u[i+1] - 2.0*u[i])*sx)/max; 657c4762a1bSJed Brown break; 658c4762a1bSJed Brown case 2: /* double obstacle */ 659c4762a1bSJed Brown len = -.5*PetscRealPart(u[i-1] + u[i+1] - 2.0*u[i])*sx/max; 660c4762a1bSJed Brown break; 661c4762a1bSJed Brown case 3: /* logarithmic + double well */ 662c4762a1bSJed Brown len = .5*PetscRealPart(6.*.25*u[i]*(u[i+1] - u[i-1])*(u[i+1] - u[i-1])*sx + (3.*u[i]*u[i] - 1.)*(u[i-1] + u[i+1] - 2.0*u[i])*sx)/max; 663c4762a1bSJed Brown if (len < 0.) ydown += len; 664c4762a1bSJed Brown else yup += len; 665c4762a1bSJed Brown 666c4762a1bSJed Brown if (ctx->truncation==2) { /* quadratic */ 667c4762a1bSJed Brown if (PetscRealPart(u[i]) < -1.0 + 2.0*tol) len2 = .5*(.25*theta/(tol-tol*tol))*PetscRealPart(u[i-1] + u[i+1] - 2.0*u[i])*sx/max; 668c4762a1bSJed Brown else if (PetscRealPart(u[i]) > 1.0 - 2.0*tol) len2 = .5*(.25*theta/(tol-tol*tol))*PetscRealPart(u[i-1] + u[i+1] - 2.0*u[i])*sx/max; 669c4762a1bSJed Brown else len2 = PetscRealPart(.5*(2.0*theta*u[i]/((1.0-u[i]*u[i])*(1.0-u[i]*u[i]))*.25*(u[i+1] - u[i-1])*(u[i+1] - u[i-1])*sx + (theta/(1.0-u[i]*u[i]))*(u[i-1] + u[i+1] - 2.0*u[i])*sx)/max); 670c4762a1bSJed Brown } else { /* cubic */ 671c4762a1bSJed Brown a = 2.0*theta*(1.0-2.0*tol)/(16.0*tol*tol*(1.0-tol)*(1.0-tol)); 672c4762a1bSJed Brown b = theta/(4.0*tol*(1.0-tol)) - a*(1.0-2.0*tol); 673c4762a1bSJed Brown if (PetscRealPart(u[i]) < -1.0 + 2.0*tol) len2 = PetscRealPart(.5*(-1.0*a*.25*(u[i+1] - u[i-1])*(u[i+1] - u[i-1])*sx + (-1.0*a*u[i] + b)*(u[i-1] + u[i+1] - 2.0*u[i])*sx)/max); 674c4762a1bSJed Brown else if (PetscRealPart(u[i]) > 1.0 - 2.0*tol) len2 = PetscRealPart(.5*(a*.25*(u[i+1] - u[i-1])*(u[i+1] - u[i-1])*sx + ( a*u[i] + b)*(u[i-1] + u[i+1] - 2.0*u[i])*sx)/max); 675c4762a1bSJed Brown else len2 = PetscRealPart(.5*(2.0*theta*u[i]/((1.0-u[i]*u[i])*(1.0-u[i]*u[i]))*.25*(u[i+1] - u[i-1])*(u[i+1] - u[i-1])*sx + (theta/(1.0-u[i]*u[i]))*(u[i-1] + u[i+1] - 2.0*u[i])*sx)/max); 676c4762a1bSJed Brown } 677c4762a1bSJed Brown y2 = len < 0 ? ydown : yup; 678*9566063dSJacob Faibussowitsch PetscCall(PetscDrawArrow(draw,x,y2,x,y2+len2,PETSC_DRAW_PLUM)); 679c4762a1bSJed Brown break; 680c4762a1bSJed Brown case 4: /* logarithmic + double obstacle */ 681c4762a1bSJed Brown len = -.5*theta_c*PetscRealPart(-(u[i-1] + u[i+1] - 2.0*u[i])*sx/max); 682c4762a1bSJed Brown if (len < 0.) ydown += len; 683c4762a1bSJed Brown else yup += len; 684c4762a1bSJed Brown 685c4762a1bSJed Brown if (ctx->truncation==2) { /* quadratic */ 686c4762a1bSJed Brown if (PetscRealPart(u[i]) < -1.0 + 2.0*tol) len2 = .5*(.25*theta/(tol-tol*tol))*PetscRealPart(u[i-1] + u[i+1] - 2.0*u[i])*sx/max; 687c4762a1bSJed Brown else if (PetscRealPart(u[i]) > 1.0 - 2.0*tol) len2 = .5*(.25*theta/(tol-tol*tol))*PetscRealPart(u[i-1] + u[i+1] - 2.0*u[i])*sx/max; 688c4762a1bSJed Brown else len2 = PetscRealPart(.5*(2.0*theta*u[i]/((1.0-u[i]*u[i])*(1.0-u[i]*u[i]))*.25*(u[i+1] - u[i-1])*(u[i+1] - u[i-1])*sx + (theta/(1.0-u[i]*u[i]))*(u[i-1] + u[i+1] - 2.0*u[i])*sx)/max); 689c4762a1bSJed Brown } else { /* cubic */ 690c4762a1bSJed Brown a = 2.0*theta*(1.0-2.0*tol)/(16.0*tol*tol*(1.0-tol)*(1.0-tol)); 691c4762a1bSJed Brown b = theta/(4.0*tol*(1.0-tol)) - a*(1.0-2.0*tol); 692c4762a1bSJed Brown if (PetscRealPart(u[i]) < -1.0 + 2.0*tol) len2 = .5*PetscRealPart(-1.0*a*.25*(u[i+1] - u[i-1])*(u[i+1] - u[i-1])*sx + (-1.0*a*u[i] + b)*(u[i-1] + u[i+1] - 2.0*u[i])*sx)/max; 693c4762a1bSJed Brown else if (PetscRealPart(u[i]) > 1.0 - 2.0*tol) len2 = .5*PetscRealPart(a*.25*(u[i+1] - u[i-1])*(u[i+1] - u[i-1])*sx + ( a*u[i] + b)*(u[i-1] + u[i+1] - 2.0*u[i])*sx)/max; 694c4762a1bSJed Brown else len2 = .5*PetscRealPart(2.0*theta*u[i]/((1.0-u[i]*u[i])*(1.0-u[i]*u[i]))*.25*(u[i+1] - u[i-1])*(u[i+1] - u[i-1])*sx + (theta/(1.0-u[i]*u[i]))*(u[i-1] + u[i+1] - 2.0*u[i])*sx)/max; 695c4762a1bSJed Brown } 696c4762a1bSJed Brown y2 = len < 0 ? ydown : yup; 697*9566063dSJacob Faibussowitsch PetscCall(PetscDrawArrow(draw,x,y2,x,y2+len2,PETSC_DRAW_PLUM)); 698c4762a1bSJed Brown break; 699c4762a1bSJed Brown } 700*9566063dSJacob Faibussowitsch PetscCall(PetscDrawArrow(draw,x,y,x,y+len,PETSC_DRAW_BLUE)); 701c4762a1bSJed Brown } 702c4762a1bSJed Brown x += cnt*hx; 703c4762a1bSJed Brown } 704*9566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da,localU,&x)); 705*9566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da,&localU)); 706*9566063dSJacob Faibussowitsch PetscCall(PetscDrawStringSetSize(draw,.2,.2)); 707*9566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 708*9566063dSJacob Faibussowitsch PetscCall(PetscDrawSetPause(draw,pause)); 709*9566063dSJacob Faibussowitsch PetscCall(PetscDrawPause(draw)); 710c4762a1bSJed Brown PetscFunctionReturn(0); 711c4762a1bSJed Brown } 712c4762a1bSJed Brown 713c4762a1bSJed Brown PetscErrorCode MyDestroy(void **ptr) 714c4762a1bSJed Brown { 715c4762a1bSJed Brown UserCtx *ctx = *(UserCtx**)ptr; 716c4762a1bSJed Brown 717c4762a1bSJed Brown PetscFunctionBegin; 718*9566063dSJacob Faibussowitsch PetscCall(PetscDrawViewPortsDestroy(ctx->ports)); 719c4762a1bSJed Brown PetscFunctionReturn(0); 720c4762a1bSJed Brown } 721c4762a1bSJed Brown 722c4762a1bSJed Brown /*TEST 723c4762a1bSJed Brown 724c4762a1bSJed Brown test: 725c4762a1bSJed Brown TODO: currently requires initial condition file generated by heat 726c4762a1bSJed Brown 727c4762a1bSJed Brown TEST*/ 728