1c4762a1bSJed Brown static char help[] ="Solves the time dependent Allen-Cahn equation with IMEX methods"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown /* 4c4762a1bSJed Brown * allen_cahn.c 5c4762a1bSJed Brown * 6c4762a1bSJed Brown * Created on: Jun 8, 2012 7c4762a1bSJed Brown * Author: Hong Zhang 8c4762a1bSJed Brown */ 9c4762a1bSJed Brown 10c4762a1bSJed Brown #include <petscts.h> 11c4762a1bSJed Brown 12c4762a1bSJed Brown /* 13c4762a1bSJed Brown * application context 14c4762a1bSJed Brown */ 15c4762a1bSJed Brown typedef struct { 16c4762a1bSJed Brown PetscReal param; /* parameter */ 17c4762a1bSJed Brown PetscReal xleft,xright; /* range in x-direction */ 18c4762a1bSJed Brown PetscInt mx; /* Discretization in x-direction */ 19c4762a1bSJed Brown }AppCtx; 20c4762a1bSJed Brown 21c4762a1bSJed Brown static PetscErrorCode RHSFunction(TS,PetscReal,Vec,Vec,void*); 22c4762a1bSJed Brown static PetscErrorCode FormIFunction(TS,PetscReal,Vec,Vec,Vec,void*); 23c4762a1bSJed Brown static PetscErrorCode FormIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void *ctx); 24c4762a1bSJed Brown static PetscErrorCode FormInitialSolution(TS,Vec,void*); 25c4762a1bSJed Brown 26c4762a1bSJed Brown int main(int argc, char **argv) 27c4762a1bSJed Brown { 28c4762a1bSJed Brown TS ts; 29c4762a1bSJed Brown Vec x; /*solution vector*/ 30c4762a1bSJed Brown Mat A; /*Jacobian*/ 31c4762a1bSJed Brown PetscInt steps,mx; 32c4762a1bSJed Brown PetscReal ftime; 33c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 34c4762a1bSJed Brown 359566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,NULL,help)); 36c4762a1bSJed Brown 37c4762a1bSJed Brown /* Initialize user application context */ 38*d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Allen-Cahn equation",""); 39c4762a1bSJed Brown user.param = 9e-4; 40c4762a1bSJed Brown user.xleft = -1.; 41c4762a1bSJed Brown user.xright = 2.; 42c4762a1bSJed Brown user.mx = 400; 439566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-eps","parameter","",user.param,&user.param,NULL)); 44*d0609cedSBarry Smith PetscOptionsEnd(); 45c4762a1bSJed Brown 46c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 47c4762a1bSJed Brown Set runtime options 48c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 49c4762a1bSJed Brown /* 509566063dSJacob Faibussowitsch * PetscCall(PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL)); 51c4762a1bSJed Brown */ 52c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 53c4762a1bSJed Brown Create necessary matrix and vectors, solve same ODE on every process 54c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 559566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 569566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,user.mx,user.mx)); 579566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 589566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 599566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A,&x,NULL)); 60c4762a1bSJed Brown 61c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 62c4762a1bSJed Brown Create time stepping solver context 63c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 649566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD,&ts)); 659566063dSJacob Faibussowitsch PetscCall(TSSetType(ts,TSEIMEX)); 669566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts,NULL,RHSFunction,&user)); 679566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts,NULL,FormIFunction,&user)); 689566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts,A,A,FormIJacobian,&user)); 69c4762a1bSJed Brown ftime = 22; 709566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts,ftime)); 719566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); 72c4762a1bSJed Brown 73c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 74c4762a1bSJed Brown Set initial conditions 75c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 769566063dSJacob Faibussowitsch PetscCall(FormInitialSolution(ts,x,&user)); 779566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts,x)); 789566063dSJacob Faibussowitsch PetscCall(VecGetSize(x,&mx)); 79c4762a1bSJed Brown 80c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 81c4762a1bSJed Brown Set runtime options 82c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 839566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 84c4762a1bSJed Brown 85c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 86c4762a1bSJed Brown Solve nonlinear system 87c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 889566063dSJacob Faibussowitsch PetscCall(TSSolve(ts,x)); 899566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts,&ftime)); 909566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts,&steps)); 919566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"eps %g, steps %D, ftime %g\n",(double)user.param,steps,(double)ftime)); 929566063dSJacob Faibussowitsch /* PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD));*/ 93c4762a1bSJed Brown 94c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 95c4762a1bSJed Brown Free work space. 96c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 979566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 989566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 999566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 1009566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 101b122ec5aSJacob Faibussowitsch return 0; 102c4762a1bSJed Brown } 103c4762a1bSJed Brown 104c4762a1bSJed Brown static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr) 105c4762a1bSJed Brown { 106c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 107c4762a1bSJed Brown PetscScalar *f; 108c4762a1bSJed Brown const PetscScalar *x; 109c4762a1bSJed Brown PetscInt i,mx; 110c4762a1bSJed Brown PetscReal hx,eps; 111c4762a1bSJed Brown 112c4762a1bSJed Brown PetscFunctionBegin; 113c4762a1bSJed Brown mx = user->mx; 114c4762a1bSJed Brown eps = user->param; 115c4762a1bSJed Brown hx = (user->xright-user->xleft)/(mx-1); 1169566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&x)); 1179566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 118c4762a1bSJed Brown f[0] = 2.*eps*(x[1]-x[0])/(hx*hx); /*boundary*/ 119c4762a1bSJed Brown for (i=1;i<mx-1;i++) { 120c4762a1bSJed Brown f[i]= eps*(x[i+1]-2.*x[i]+x[i-1])/(hx*hx); 121c4762a1bSJed Brown } 122c4762a1bSJed Brown f[mx-1] = 2.*eps*(x[mx-2]- x[mx-1])/(hx*hx); /*boundary*/ 1239566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&x)); 1249566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 125c4762a1bSJed Brown PetscFunctionReturn(0); 126c4762a1bSJed Brown } 127c4762a1bSJed Brown 128c4762a1bSJed Brown static PetscErrorCode FormIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ptr) 129c4762a1bSJed Brown { 130c4762a1bSJed Brown AppCtx *user = (AppCtx*)ptr; 131c4762a1bSJed Brown PetscScalar *f; 132c4762a1bSJed Brown const PetscScalar *x,*xdot; 133c4762a1bSJed Brown PetscInt i,mx; 134c4762a1bSJed Brown 135c4762a1bSJed Brown PetscFunctionBegin; 136c4762a1bSJed Brown mx = user->mx; 1379566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&x)); 1389566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Xdot,&xdot)); 1399566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 140c4762a1bSJed Brown 141c4762a1bSJed Brown for (i=0;i<mx;i++) { 142c4762a1bSJed Brown f[i]= xdot[i] - x[i]*(1.-x[i]*x[i]); 143c4762a1bSJed Brown } 144c4762a1bSJed Brown 1459566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&x)); 1469566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Xdot,&xdot)); 1479566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 148c4762a1bSJed Brown PetscFunctionReturn(0); 149c4762a1bSJed Brown } 150c4762a1bSJed Brown 151c4762a1bSJed Brown static PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec U, Vec Udot, PetscReal a, Mat J,Mat Jpre,void *ctx) 152c4762a1bSJed Brown { 153c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx; 154c4762a1bSJed Brown PetscScalar v; 155c4762a1bSJed Brown const PetscScalar *x; 156c4762a1bSJed Brown PetscInt i,col; 157c4762a1bSJed Brown 158c4762a1bSJed Brown PetscFunctionBegin; 1599566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U,&x)); 160c4762a1bSJed Brown for (i=0; i < user->mx; i++) { 161c4762a1bSJed Brown v = a - 1. + 3.*x[i]*x[i]; 162c4762a1bSJed Brown col = i; 1639566063dSJacob Faibussowitsch PetscCall(MatSetValues(J,1,&i,1,&col,&v,INSERT_VALUES)); 164c4762a1bSJed Brown } 1659566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U,&x)); 166c4762a1bSJed Brown 1679566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY)); 1689566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY)); 169c4762a1bSJed Brown if (J != Jpre) { 1709566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY)); 1719566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY)); 172c4762a1bSJed Brown } 173c4762a1bSJed Brown /* MatView(J,PETSC_VIEWER_STDOUT_WORLD);*/ 174c4762a1bSJed Brown PetscFunctionReturn(0); 175c4762a1bSJed Brown } 176c4762a1bSJed Brown 177c4762a1bSJed Brown static PetscErrorCode FormInitialSolution(TS ts,Vec U,void *ctx) 178c4762a1bSJed Brown { 179c4762a1bSJed Brown AppCtx *user = (AppCtx*)ctx; 180c4762a1bSJed Brown PetscInt i; 181c4762a1bSJed Brown PetscScalar *x; 182c4762a1bSJed Brown PetscReal hx,x_map; 183c4762a1bSJed Brown 184c4762a1bSJed Brown PetscFunctionBegin; 185c4762a1bSJed Brown hx = (user->xright-user->xleft)/(PetscReal)(user->mx-1); 1869566063dSJacob Faibussowitsch PetscCall(VecGetArray(U,&x)); 187c4762a1bSJed Brown for (i=0;i<user->mx;i++) { 188c4762a1bSJed Brown x_map = user->xleft + i*hx; 189c4762a1bSJed Brown if (x_map >= 0.7065) { 190c4762a1bSJed Brown x[i] = PetscTanhReal((x_map-0.8)/(2.*PetscSqrtReal(user->param))); 191c4762a1bSJed Brown } else if (x_map >= 0.4865) { 192c4762a1bSJed Brown x[i] = PetscTanhReal((0.613-x_map)/(2.*PetscSqrtReal(user->param))); 193c4762a1bSJed Brown } else if (x_map >= 0.28) { 194c4762a1bSJed Brown x[i] = PetscTanhReal((x_map-0.36)/(2.*PetscSqrtReal(user->param))); 195c4762a1bSJed Brown } else if (x_map >= -0.7) { 196c4762a1bSJed Brown x[i] = PetscTanhReal((0.2-x_map)/(2.*PetscSqrtReal(user->param))); 197c4762a1bSJed Brown } else if (x_map >= -1) { 198c4762a1bSJed Brown x[i] = PetscTanhReal((x_map+0.9)/(2.*PetscSqrtReal(user->param))); 199c4762a1bSJed Brown } 200c4762a1bSJed Brown } 2019566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U,&x)); 202c4762a1bSJed Brown PetscFunctionReturn(0); 203c4762a1bSJed Brown } 204c4762a1bSJed Brown 205c4762a1bSJed Brown /*TEST 206c4762a1bSJed Brown 207c4762a1bSJed Brown test: 208c4762a1bSJed 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 209c4762a1bSJed Brown requires: x 210c4762a1bSJed Brown timeoutfactor: 3 211c4762a1bSJed Brown 212c4762a1bSJed Brown TEST*/ 213