1 static char help[] = "Solves the trival ODE 2 du/dt = 1, u(0) = 0. \n\n"; 2 3 #include <petscts.h> 4 #include <petscpc.h> 5 6 PetscErrorCode IFunction(TS,PetscReal,Vec,Vec,Vec,void*); 7 PetscErrorCode IJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*); 8 9 int main(int argc,char **argv) 10 { 11 TS ts; 12 Vec x; 13 Vec f; 14 Mat A; 15 PetscErrorCode ierr; 16 17 ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 18 19 CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts)); 20 CHKERRQ(TSSetEquationType(ts,TS_EQ_ODE_IMPLICIT)); 21 CHKERRQ(VecCreate(PETSC_COMM_WORLD,&f)); 22 CHKERRQ(VecSetSizes(f,1,PETSC_DECIDE)); 23 CHKERRQ(VecSetFromOptions(f)); 24 CHKERRQ(VecSetUp(f)); 25 CHKERRQ(TSSetIFunction(ts,f,IFunction,NULL)); 26 CHKERRQ(VecDestroy(&f)); 27 28 CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 29 CHKERRQ(MatSetSizes(A,1,1,PETSC_DECIDE,PETSC_DECIDE)); 30 CHKERRQ(MatSetFromOptions(A)); 31 CHKERRQ(MatSetUp(A)); 32 CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 33 CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 34 /* ensure that the Jacobian matrix has diagonal entries since that is required by TS */ 35 CHKERRQ(MatShift(A,(PetscReal)1)); 36 CHKERRQ(MatShift(A,(PetscReal)-1)); 37 CHKERRQ(TSSetIJacobian(ts,A,A,IJacobian,NULL)); 38 CHKERRQ(MatDestroy(&A)); 39 40 CHKERRQ(VecCreate(PETSC_COMM_WORLD,&x)); 41 CHKERRQ(VecSetSizes(x,1,PETSC_DECIDE)); 42 CHKERRQ(VecSetFromOptions(x)); 43 CHKERRQ(VecSetUp(x)); 44 CHKERRQ(TSSetSolution(ts,x)); 45 CHKERRQ(VecDestroy(&x)); 46 CHKERRQ(TSSetFromOptions(ts)); 47 48 CHKERRQ(TSSetStepNumber(ts,0)); 49 CHKERRQ(TSSetTimeStep(ts,1)); 50 CHKERRQ(TSSetTime(ts,0)); 51 CHKERRQ(TSSetMaxTime(ts,PETSC_MAX_REAL)); 52 CHKERRQ(TSSetMaxSteps(ts,3)); 53 54 /* 55 When an ARKIMEX scheme with an explicit stage is used this will error with a message informing the user it is not possible to use 56 a non-trivial mass matrix with ARKIMEX schemes with explicit stages. 57 */ 58 ierr = TSSolve(ts,NULL); 59 if (ierr != PETSC_ERR_ARG_INCOMP) CHKERRQ(ierr); 60 61 CHKERRQ(TSDestroy(&ts)); 62 ierr = PetscFinalize(); 63 return ierr; 64 } 65 66 PetscErrorCode IFunction(TS ts,PetscReal t,Vec x,Vec xdot,Vec f,void *ctx) 67 { 68 PetscFunctionBegin; 69 CHKERRQ(VecCopy(xdot,f)); 70 CHKERRQ(VecScale(f,2.0)); 71 CHKERRQ(VecShift(f,-1.0)); 72 PetscFunctionReturn(0); 73 } 74 75 PetscErrorCode IJacobian(TS ts,PetscReal t,Vec x,Vec xdot,PetscReal shift,Mat A,Mat B,void *ctx) 76 { 77 PetscScalar j; 78 79 PetscFunctionBegin; 80 j = shift*2.0; 81 CHKERRQ(MatSetValue(B,0,0,j,INSERT_VALUES)); 82 CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 83 CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 84 PetscFunctionReturn(0); 85 } 86 87 /*TEST 88 89 test: 90 suffix: arkimex_explicit_stage 91 requires: defined(PETSC_USE_DEBUG) 92 args: -ts_type arkimex -error_output_stdout 93 filter: egrep -v "(Petsc|on a| at |Configure)" 94 95 test: 96 suffix: arkimex_implicit_stage 97 args: -ts_type arkimex -ts_arkimex_type l2 -ts_monitor_solution -ts_monitor 98 99 TEST*/ 100