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