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