1*c4762a1bSJed Brown static char help[] ="Solves the time dependent Allen-Cahn equation with IMEX methods"; 2*c4762a1bSJed Brown 3*c4762a1bSJed Brown /* 4*c4762a1bSJed Brown * allen_cahn.c 5*c4762a1bSJed Brown * 6*c4762a1bSJed Brown * Created on: Jun 8, 2012 7*c4762a1bSJed Brown * Author: Hong Zhang 8*c4762a1bSJed Brown */ 9*c4762a1bSJed Brown 10*c4762a1bSJed Brown #include <petscts.h> 11*c4762a1bSJed Brown 12*c4762a1bSJed Brown /* 13*c4762a1bSJed Brown * application context 14*c4762a1bSJed Brown */ 15*c4762a1bSJed Brown typedef struct { 16*c4762a1bSJed Brown PetscReal param; /* parameter */ 17*c4762a1bSJed Brown PetscReal xleft,xright; /* range in x-direction */ 18*c4762a1bSJed Brown PetscInt mx; /* Discretization in x-direction */ 19*c4762a1bSJed Brown }AppCtx; 20*c4762a1bSJed Brown 21*c4762a1bSJed Brown 22*c4762a1bSJed Brown static PetscErrorCode RHSFunction(TS,PetscReal,Vec,Vec,void*); 23*c4762a1bSJed Brown static PetscErrorCode FormIFunction(TS,PetscReal,Vec,Vec,Vec,void*); 24*c4762a1bSJed Brown static PetscErrorCode FormIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void *ctx); 25*c4762a1bSJed Brown static PetscErrorCode FormInitialSolution(TS,Vec,void*); 26*c4762a1bSJed Brown 27*c4762a1bSJed Brown 28*c4762a1bSJed Brown int main(int argc, char **argv) 29*c4762a1bSJed Brown { 30*c4762a1bSJed Brown TS ts; 31*c4762a1bSJed Brown Vec x; /*solution vector*/ 32*c4762a1bSJed Brown Mat A; /*Jacobian*/ 33*c4762a1bSJed Brown PetscInt steps,mx; 34*c4762a1bSJed Brown PetscErrorCode ierr; 35*c4762a1bSJed Brown PetscReal ftime; 36*c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 37*c4762a1bSJed Brown 38*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; 39*c4762a1bSJed Brown 40*c4762a1bSJed Brown /* Initialize user application context */ 41*c4762a1bSJed Brown ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Allen-Cahn equation",""); 42*c4762a1bSJed Brown user.param = 9e-4; 43*c4762a1bSJed Brown user.xleft = -1.; 44*c4762a1bSJed Brown user.xright = 2.; 45*c4762a1bSJed Brown user.mx = 400; 46*c4762a1bSJed Brown ierr = PetscOptionsReal("-eps","parameter","",user.param,&user.param,NULL);CHKERRQ(ierr); 47*c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 48*c4762a1bSJed Brown 49*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 50*c4762a1bSJed Brown Set runtime options 51*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 52*c4762a1bSJed Brown /* 53*c4762a1bSJed Brown * ierr = PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL);CHKERRQ(ierr); 54*c4762a1bSJed Brown */ 55*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 56*c4762a1bSJed Brown Create necessary matrix and vectors, solve same ODE on every process 57*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 58*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 59*c4762a1bSJed Brown ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,user.mx,user.mx);CHKERRQ(ierr); 60*c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 61*c4762a1bSJed Brown ierr = MatSetUp(A);CHKERRQ(ierr); 62*c4762a1bSJed Brown ierr = MatCreateVecs(A,&x,NULL);CHKERRQ(ierr); 63*c4762a1bSJed Brown 64*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 65*c4762a1bSJed Brown Create time stepping solver context 66*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 67*c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 68*c4762a1bSJed Brown ierr = TSSetType(ts,TSEIMEX);CHKERRQ(ierr); 69*c4762a1bSJed Brown ierr = TSSetRHSFunction(ts,NULL,RHSFunction,&user);CHKERRQ(ierr); 70*c4762a1bSJed Brown ierr = TSSetIFunction(ts,NULL,FormIFunction,&user);CHKERRQ(ierr); 71*c4762a1bSJed Brown ierr = TSSetIJacobian(ts,A,A,FormIJacobian,&user);CHKERRQ(ierr); 72*c4762a1bSJed Brown ftime = 22; 73*c4762a1bSJed Brown ierr = TSSetMaxTime(ts,ftime);CHKERRQ(ierr); 74*c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); 75*c4762a1bSJed Brown 76*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 77*c4762a1bSJed Brown Set initial conditions 78*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 79*c4762a1bSJed Brown ierr = FormInitialSolution(ts,x,&user);CHKERRQ(ierr); 80*c4762a1bSJed Brown ierr = TSSetSolution(ts,x);CHKERRQ(ierr); 81*c4762a1bSJed Brown ierr = VecGetSize(x,&mx);CHKERRQ(ierr); 82*c4762a1bSJed Brown 83*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 84*c4762a1bSJed Brown Set runtime options 85*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 86*c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 87*c4762a1bSJed Brown 88*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 89*c4762a1bSJed Brown Solve nonlinear system 90*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 91*c4762a1bSJed Brown ierr = TSSolve(ts,x);CHKERRQ(ierr); 92*c4762a1bSJed Brown ierr = TSGetTime(ts,&ftime);CHKERRQ(ierr); 93*c4762a1bSJed Brown ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr); 94*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"eps %g, steps %D, ftime %g\n",(double)user.param,steps,(double)ftime);CHKERRQ(ierr); 95*c4762a1bSJed Brown /* ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);*/ 96*c4762a1bSJed Brown 97*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 98*c4762a1bSJed Brown Free work space. 99*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 100*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 101*c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 102*c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 103*c4762a1bSJed Brown ierr = PetscFinalize(); 104*c4762a1bSJed Brown return ierr; 105*c4762a1bSJed Brown } 106*c4762a1bSJed Brown 107*c4762a1bSJed Brown 108*c4762a1bSJed Brown static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr) 109*c4762a1bSJed Brown { 110*c4762a1bSJed Brown PetscErrorCode ierr; 111*c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 112*c4762a1bSJed Brown PetscScalar *f; 113*c4762a1bSJed Brown const PetscScalar *x; 114*c4762a1bSJed Brown PetscInt i,mx; 115*c4762a1bSJed Brown PetscReal hx,eps; 116*c4762a1bSJed Brown 117*c4762a1bSJed Brown PetscFunctionBegin; 118*c4762a1bSJed Brown mx = user->mx; 119*c4762a1bSJed Brown eps = user->param; 120*c4762a1bSJed Brown hx = (user->xright-user->xleft)/(mx-1); 121*c4762a1bSJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 122*c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 123*c4762a1bSJed Brown f[0] = 2.*eps*(x[1]-x[0])/(hx*hx); /*boundary*/ 124*c4762a1bSJed Brown for(i=1;i<mx-1;i++) { 125*c4762a1bSJed Brown f[i]= eps*(x[i+1]-2.*x[i]+x[i-1])/(hx*hx); 126*c4762a1bSJed Brown } 127*c4762a1bSJed Brown f[mx-1] = 2.*eps*(x[mx-2]- x[mx-1])/(hx*hx); /*boundary*/ 128*c4762a1bSJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 129*c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 130*c4762a1bSJed Brown PetscFunctionReturn(0); 131*c4762a1bSJed Brown } 132*c4762a1bSJed Brown 133*c4762a1bSJed Brown 134*c4762a1bSJed Brown static PetscErrorCode FormIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ptr) 135*c4762a1bSJed Brown { 136*c4762a1bSJed Brown PetscErrorCode ierr; 137*c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 138*c4762a1bSJed Brown PetscScalar *f; 139*c4762a1bSJed Brown const PetscScalar *x,*xdot; 140*c4762a1bSJed Brown PetscInt i,mx; 141*c4762a1bSJed Brown 142*c4762a1bSJed Brown PetscFunctionBegin; 143*c4762a1bSJed Brown mx = user->mx; 144*c4762a1bSJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 145*c4762a1bSJed Brown ierr = VecGetArrayRead(Xdot,&xdot);CHKERRQ(ierr); 146*c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 147*c4762a1bSJed Brown 148*c4762a1bSJed Brown for(i=0;i<mx;i++) { 149*c4762a1bSJed Brown f[i]= xdot[i] - x[i]*(1.-x[i]*x[i]); 150*c4762a1bSJed Brown } 151*c4762a1bSJed Brown 152*c4762a1bSJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 153*c4762a1bSJed Brown ierr = VecRestoreArrayRead(Xdot,&xdot);CHKERRQ(ierr); 154*c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 155*c4762a1bSJed Brown PetscFunctionReturn(0); 156*c4762a1bSJed Brown } 157*c4762a1bSJed Brown 158*c4762a1bSJed Brown 159*c4762a1bSJed Brown static PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec U, Vec Udot, PetscReal a, Mat J,Mat Jpre,void *ctx) 160*c4762a1bSJed Brown { 161*c4762a1bSJed Brown PetscErrorCode ierr; 162*c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx; 163*c4762a1bSJed Brown PetscScalar v; 164*c4762a1bSJed Brown const PetscScalar *x; 165*c4762a1bSJed Brown PetscInt i,col; 166*c4762a1bSJed Brown 167*c4762a1bSJed Brown PetscFunctionBegin; 168*c4762a1bSJed Brown ierr = VecGetArrayRead(U,&x);CHKERRQ(ierr); 169*c4762a1bSJed Brown for(i=0; i < user->mx; i++) { 170*c4762a1bSJed Brown v = a - 1. + 3.*x[i]*x[i]; 171*c4762a1bSJed Brown col = i; 172*c4762a1bSJed Brown ierr = MatSetValues(J,1,&i,1,&col,&v,INSERT_VALUES);CHKERRQ(ierr); 173*c4762a1bSJed Brown } 174*c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&x);CHKERRQ(ierr); 175*c4762a1bSJed Brown 176*c4762a1bSJed Brown ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 177*c4762a1bSJed Brown ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 178*c4762a1bSJed Brown if (J != Jpre) { 179*c4762a1bSJed Brown ierr = MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 180*c4762a1bSJed Brown ierr = MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 181*c4762a1bSJed Brown } 182*c4762a1bSJed Brown /* MatView(J,PETSC_VIEWER_STDOUT_WORLD);*/ 183*c4762a1bSJed Brown PetscFunctionReturn(0); 184*c4762a1bSJed Brown } 185*c4762a1bSJed Brown 186*c4762a1bSJed Brown 187*c4762a1bSJed Brown static PetscErrorCode FormInitialSolution(TS ts,Vec U,void *ctx) 188*c4762a1bSJed Brown { 189*c4762a1bSJed Brown AppCtx *user = (AppCtx*)ctx; 190*c4762a1bSJed Brown PetscInt i; 191*c4762a1bSJed Brown PetscScalar *x; 192*c4762a1bSJed Brown PetscReal hx,x_map; 193*c4762a1bSJed Brown PetscErrorCode ierr; 194*c4762a1bSJed Brown 195*c4762a1bSJed Brown PetscFunctionBegin; 196*c4762a1bSJed Brown hx = (user->xright-user->xleft)/(PetscReal)(user->mx-1); 197*c4762a1bSJed Brown ierr = VecGetArray(U,&x);CHKERRQ(ierr); 198*c4762a1bSJed Brown for(i=0;i<user->mx;i++) { 199*c4762a1bSJed Brown x_map = user->xleft + i*hx; 200*c4762a1bSJed Brown if(x_map >= 0.7065) { 201*c4762a1bSJed Brown x[i] = PetscTanhReal((x_map-0.8)/(2.*PetscSqrtReal(user->param))); 202*c4762a1bSJed Brown } else if(x_map >= 0.4865) { 203*c4762a1bSJed Brown x[i] = PetscTanhReal((0.613-x_map)/(2.*PetscSqrtReal(user->param))); 204*c4762a1bSJed Brown } else if(x_map >= 0.28) { 205*c4762a1bSJed Brown x[i] = PetscTanhReal((x_map-0.36)/(2.*PetscSqrtReal(user->param))); 206*c4762a1bSJed Brown } else if(x_map >= -0.7) { 207*c4762a1bSJed Brown x[i] = PetscTanhReal((0.2-x_map)/(2.*PetscSqrtReal(user->param))); 208*c4762a1bSJed Brown } else if(x_map >= -1) { 209*c4762a1bSJed Brown x[i] = PetscTanhReal((x_map+0.9)/(2.*PetscSqrtReal(user->param))); 210*c4762a1bSJed Brown } 211*c4762a1bSJed Brown } 212*c4762a1bSJed Brown ierr = VecRestoreArray(U,&x);CHKERRQ(ierr); 213*c4762a1bSJed Brown PetscFunctionReturn(0); 214*c4762a1bSJed Brown } 215*c4762a1bSJed Brown 216*c4762a1bSJed Brown /*TEST 217*c4762a1bSJed Brown 218*c4762a1bSJed Brown test: 219*c4762a1bSJed Brown args: -ts_rtol 1e-04 -ts_dt 0.025 -pc_type lu -ksp_error_if_not_converged TRUE -ts_type eimex -ts_adapt_type none -ts_eimex_order_adapt -ts_eimex_max_rows 7 -ts_monitor_draw_solution 220*c4762a1bSJed Brown requires: x 221*c4762a1bSJed Brown timeoutfactor: 3 222*c4762a1bSJed Brown 223*c4762a1bSJed Brown TEST*/ 224