/* Formatted test for TS routines. Solves U_t=F(t,u) Where: [2*u1+u2 F(t,u)= [u1+2*u2+u3 [ u2+2*u3 We can compare the solutions from euler, beuler and SUNDIALS to see what is the difference. */ static char help[] = "Solves a linear ODE. \n\n"; #include #include extern PetscErrorCode RHSFunction(TS,PetscReal,Vec,Vec,void*); extern PetscErrorCode RHSJacobian(TS,PetscReal,Vec,Mat,Mat,void*); extern PetscErrorCode Monitor(TS,PetscInt,PetscReal,Vec,void*); extern PetscErrorCode Initial(Vec,void*); extern PetscErrorCode MyMatMult(Mat,Vec,Vec); extern PetscReal solx(PetscReal); extern PetscReal soly(PetscReal); extern PetscReal solz(PetscReal); int main(int argc,char **argv) { PetscErrorCode ierr; PetscInt time_steps = 100,steps; Vec global; PetscReal dt,ftime; TS ts; Mat A = 0,S; ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-time",&time_steps,NULL)); /* set initial conditions */ CHKERRQ(VecCreate(PETSC_COMM_WORLD,&global)); CHKERRQ(VecSetSizes(global,PETSC_DECIDE,3)); CHKERRQ(VecSetFromOptions(global)); CHKERRQ(Initial(global,NULL)); /* make timestep context */ CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts)); CHKERRQ(TSSetProblemType(ts,TS_NONLINEAR)); CHKERRQ(TSMonitorSet(ts,Monitor,NULL,NULL)); dt = 0.001; /* The user provides the RHS and Jacobian */ CHKERRQ(TSSetRHSFunction(ts,NULL,RHSFunction,NULL)); CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,3,3)); CHKERRQ(MatSetFromOptions(A)); CHKERRQ(MatSetUp(A)); CHKERRQ(RHSJacobian(ts,0.0,global,A,A,NULL)); CHKERRQ(TSSetRHSJacobian(ts,A,A,RHSJacobian,NULL)); CHKERRQ(MatCreateShell(PETSC_COMM_WORLD,3,3,3,3,NULL,&S)); CHKERRQ(MatShellSetOperation(S,MATOP_MULT,(void (*)(void))MyMatMult)); CHKERRQ(TSSetRHSJacobian(ts,S,A,RHSJacobian,NULL)); CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP)); CHKERRQ(TSSetFromOptions(ts)); CHKERRQ(TSSetTimeStep(ts,dt)); CHKERRQ(TSSetMaxSteps(ts,time_steps)); CHKERRQ(TSSetMaxTime(ts,1)); CHKERRQ(TSSetSolution(ts,global)); CHKERRQ(TSSolve(ts,global)); CHKERRQ(TSGetSolveTime(ts,&ftime)); CHKERRQ(TSGetStepNumber(ts,&steps)); /* free the memories */ CHKERRQ(TSDestroy(&ts)); CHKERRQ(VecDestroy(&global)); CHKERRQ(MatDestroy(&A)); CHKERRQ(MatDestroy(&S)); ierr = PetscFinalize(); return ierr; } PetscErrorCode MyMatMult(Mat S,Vec x,Vec y) { const PetscScalar *inptr; PetscScalar *outptr; PetscFunctionBegin; CHKERRQ(VecGetArrayRead(x,&inptr)); CHKERRQ(VecGetArrayWrite(y,&outptr)); outptr[0] = 2.0*inptr[0]+inptr[1]; outptr[1] = inptr[0]+2.0*inptr[1]+inptr[2]; outptr[2] = inptr[1]+2.0*inptr[2]; CHKERRQ(VecRestoreArrayRead(x,&inptr)); CHKERRQ(VecRestoreArrayWrite(y,&outptr)); PetscFunctionReturn(0); } /* this test problem has initial values (1,1,1). */ PetscErrorCode Initial(Vec global,void *ctx) { PetscScalar *localptr; PetscInt i,mybase,myend,locsize; /* determine starting point of each processor */ CHKERRQ(VecGetOwnershipRange(global,&mybase,&myend)); CHKERRQ(VecGetLocalSize(global,&locsize)); /* Initialize the array */ CHKERRQ(VecGetArrayWrite(global,&localptr)); for (i=0; i