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