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