1c4762a1bSJed Brown static char help[] = "Small ODE to test TS accuracy.\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown /* 4c4762a1bSJed Brown The ODE 5c4762a1bSJed Brown u1_t = cos(t), 6c4762a1bSJed Brown u2_t = sin(u2) 7c4762a1bSJed Brown with analytical solution 8c4762a1bSJed Brown u1(t) = sin(t), 9c4762a1bSJed Brown u2(t) = 2 * atan(exp(t) * tan(0.5)) 10c4762a1bSJed Brown is used to test the accuracy of TS schemes. 11c4762a1bSJed Brown */ 12c4762a1bSJed Brown 13c4762a1bSJed Brown #include <petscts.h> 14c4762a1bSJed Brown 15c4762a1bSJed Brown /* 16c4762a1bSJed Brown Defines the ODE passed to the ODE solver in explicit form: U_t = F(U) 17c4762a1bSJed Brown */ 18c4762a1bSJed Brown static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec F, void *s) 19c4762a1bSJed Brown { 20c4762a1bSJed Brown PetscScalar *f; 21c4762a1bSJed Brown const PetscScalar *u; 22c4762a1bSJed Brown 23c4762a1bSJed Brown PetscFunctionBegin; 245f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(U,&u)); 255f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(F,&f)); 26c4762a1bSJed Brown 27c4762a1bSJed Brown f[0] = PetscCosReal(t); 28c4762a1bSJed Brown f[1] = PetscSinReal(u[1]); 29c4762a1bSJed Brown 305f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(U,&u)); 315f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(F,&f)); 32c4762a1bSJed Brown PetscFunctionReturn(0); 33c4762a1bSJed Brown } 34c4762a1bSJed Brown 35c4762a1bSJed Brown /* 36c4762a1bSJed Brown Defines the exact solution. 37c4762a1bSJed Brown */ 38c4762a1bSJed Brown static PetscErrorCode ExactSolution(PetscReal t, Vec U) 39c4762a1bSJed Brown { 40c4762a1bSJed Brown PetscScalar *u; 41c4762a1bSJed Brown 42c4762a1bSJed Brown PetscFunctionBegin; 435f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(U,&u)); 44c4762a1bSJed Brown 45c4762a1bSJed Brown u[0] = PetscSinReal(t); 46c4762a1bSJed Brown u[1] = 2 * PetscAtanReal(PetscExpReal(t) * PetscTanReal(0.5)); 47c4762a1bSJed Brown 485f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(U,&u)); 49c4762a1bSJed Brown PetscFunctionReturn(0); 50c4762a1bSJed Brown } 51c4762a1bSJed Brown 52c4762a1bSJed Brown int main(int argc,char **argv) 53c4762a1bSJed Brown { 54c4762a1bSJed Brown TS ts; /* ODE integrator */ 55c4762a1bSJed Brown Vec U; /* numerical solution will be stored here */ 56c4762a1bSJed Brown Vec Uex; /* analytical (exact) solution will be stored here */ 57c4762a1bSJed Brown PetscMPIInt size; 58c4762a1bSJed Brown PetscInt n = 2; 59c4762a1bSJed Brown PetscScalar *u; 60c4762a1bSJed Brown PetscReal t, final_time = 1.0, dt = 0.25; 61c4762a1bSJed Brown PetscReal error; 62c4762a1bSJed Brown TSAdapt adapt; 63c4762a1bSJed Brown 64c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 65c4762a1bSJed Brown Initialize program 66c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 67*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help)); 685f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 693c633725SBarry Smith PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Only for sequential runs"); 70c4762a1bSJed Brown 71c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 72c4762a1bSJed Brown Create timestepping solver context 73c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 745f80ce2aSJacob Faibussowitsch CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts)); 755f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetType(ts,TSROSW)); 76c4762a1bSJed Brown 77c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 78c4762a1bSJed Brown Set ODE routines 79c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 805f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetProblemType(ts,TS_NONLINEAR)); 815f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSFunction(ts,NULL,RHSFunction,NULL)); 82c4762a1bSJed Brown 83c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 84c4762a1bSJed Brown Set initial conditions 85c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 865f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_WORLD,&U)); 875f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(U,n,PETSC_DETERMINE)); 885f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetUp(U)); 895f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(U,&u)); 90c4762a1bSJed Brown u[0] = 0.0; 91c4762a1bSJed Brown u[1] = 1.0; 925f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(U,&u)); 935f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetSolution(ts,U)); 94c4762a1bSJed Brown 95c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 96c4762a1bSJed Brown Set solver options 97c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 985f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetSaveTrajectory(ts)); 995f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxTime(ts,final_time)); 1005f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); 1015f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTimeStep(ts,dt)); 102a5b23f4aSJose E. Roman /* The adaptive time step controller is forced to take constant time steps. */ 1035f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetAdapt(ts,&adapt)); 1045f80ce2aSJacob Faibussowitsch CHKERRQ(TSAdaptSetType(adapt,TSADAPTNONE)); 105c4762a1bSJed Brown 1065f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(ts)); 107c4762a1bSJed Brown 108c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 109c4762a1bSJed Brown Run timestepping solver and compute error 110c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1115f80ce2aSJacob Faibussowitsch CHKERRQ(TSSolve(ts,U)); 1125f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetTime(ts,&t)); 113c4762a1bSJed Brown 114c4762a1bSJed Brown if (PetscAbsReal(t-final_time)>100*PETSC_MACHINE_EPSILON) { 1155f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Note: There is a difference of %g between the prescribed final time %g and the actual final time.\n",(double)(final_time-t),(double)final_time)); 116c4762a1bSJed Brown } 1175f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(U,&Uex)); 1185f80ce2aSJacob Faibussowitsch CHKERRQ(ExactSolution(t,Uex)); 119c4762a1bSJed Brown 1205f80ce2aSJacob Faibussowitsch CHKERRQ(VecAYPX(Uex,-1.0,U)); 1215f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(Uex,NORM_2,&error)); 1225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Error at final time: %.2E\n",(double)error)); 123c4762a1bSJed Brown 124c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 125c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they are no longer needed. 126c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1275f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&U)); 1285f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&Uex)); 1295f80ce2aSJacob Faibussowitsch CHKERRQ(TSDestroy(&ts)); 130c4762a1bSJed Brown 131*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 132*b122ec5aSJacob Faibussowitsch return 0; 133c4762a1bSJed Brown } 134c4762a1bSJed Brown 135c4762a1bSJed Brown /*TEST 136c4762a1bSJed Brown 137c4762a1bSJed Brown test: 138c4762a1bSJed Brown suffix: 3bs 139c4762a1bSJed Brown args: -ts_type rk -ts_rk_type 3bs 140c4762a1bSJed Brown requires: !single 141c4762a1bSJed Brown 142c4762a1bSJed Brown test: 143c4762a1bSJed Brown suffix: 5bs 144c4762a1bSJed Brown args: -ts_type rk -ts_rk_type 5bs 145c4762a1bSJed Brown requires: !single 146c4762a1bSJed Brown 147c4762a1bSJed Brown test: 148c4762a1bSJed Brown suffix: 5dp 149c4762a1bSJed Brown args: -ts_type rk -ts_rk_type 5dp 150c4762a1bSJed Brown requires: !single 151c4762a1bSJed Brown 152c4762a1bSJed Brown test: 153c4762a1bSJed Brown suffix: 6vr 154c4762a1bSJed Brown args: -ts_type rk -ts_rk_type 6vr 155c4762a1bSJed Brown requires: !single 156c4762a1bSJed Brown 157c4762a1bSJed Brown test: 158c4762a1bSJed Brown suffix: 7vr 159c4762a1bSJed Brown args: -ts_type rk -ts_rk_type 7vr 160c4762a1bSJed Brown requires: !single 161c4762a1bSJed Brown 162c4762a1bSJed Brown test: 163c4762a1bSJed Brown suffix: 8vr 164c4762a1bSJed Brown args: -ts_type rk -ts_rk_type 8vr 165c4762a1bSJed Brown requires: !single 166c4762a1bSJed Brown 167c4762a1bSJed Brown TEST*/ 168