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 19*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts)); 20*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetEquationType(ts,TS_EQ_ODE_IMPLICIT)); 21*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_WORLD,&f)); 22*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(f,1,PETSC_DECIDE)); 23*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(f)); 24*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetUp(f)); 25*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetIFunction(ts,f,IFunction,NULL)); 26*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&f)); 27c79dcfbdSBarry Smith 28*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 29*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A,1,1,PETSC_DECIDE,PETSC_DECIDE)); 30*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A)); 31*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(A)); 32*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 33*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 34c79dcfbdSBarry Smith /* ensure that the Jacobian matrix has diagonal entries since that is required by TS */ 35*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShift(A,(PetscReal)1)); 36*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShift(A,(PetscReal)-1)); 37*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetIJacobian(ts,A,A,IJacobian,NULL)); 38*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 39c79dcfbdSBarry Smith 40*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_WORLD,&x)); 41*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(x,1,PETSC_DECIDE)); 42*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(x)); 43*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetUp(x)); 44*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetSolution(ts,x)); 45*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 46*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(ts)); 47c79dcfbdSBarry Smith 48*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetStepNumber(ts,0)); 49*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTimeStep(ts,1)); 50*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTime(ts,0)); 51*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxTime(ts,PETSC_MAX_REAL)); 52*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxSteps(ts,3)); 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 61*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSDestroy(&ts)); 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 PetscFunctionBegin; 69*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(xdot,f)); 70*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScale(f,2.0)); 71*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecShift(f,-1.0)); 72c79dcfbdSBarry Smith PetscFunctionReturn(0); 73c79dcfbdSBarry Smith } 74c79dcfbdSBarry Smith 7591de4d3cSBarry Smith PetscErrorCode IJacobian(TS ts,PetscReal t,Vec x,Vec xdot,PetscReal shift,Mat A,Mat B,void *ctx) 76c79dcfbdSBarry Smith { 77c79dcfbdSBarry Smith PetscScalar j; 78c79dcfbdSBarry Smith 79c79dcfbdSBarry Smith PetscFunctionBegin; 80c79dcfbdSBarry Smith j = shift*2.0; 81*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValue(B,0,0,j,INSERT_VALUES)); 82*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 83*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 84c79dcfbdSBarry Smith PetscFunctionReturn(0); 85c79dcfbdSBarry Smith } 86c79dcfbdSBarry Smith 87c79dcfbdSBarry Smith /*TEST 88c79dcfbdSBarry Smith 89c79dcfbdSBarry Smith test: 90c79dcfbdSBarry Smith suffix: arkimex_explicit_stage 91dfd57a17SPierre Jolivet requires: defined(PETSC_USE_DEBUG) 92c79dcfbdSBarry Smith args: -ts_type arkimex -error_output_stdout 93184d557bSJunchao Zhang filter: egrep -v "(Petsc|on a| at |Configure)" 94c79dcfbdSBarry Smith 95c79dcfbdSBarry Smith test: 96c79dcfbdSBarry Smith suffix: arkimex_implicit_stage 97c79dcfbdSBarry Smith args: -ts_type arkimex -ts_arkimex_type l2 -ts_monitor_solution -ts_monitor 98c79dcfbdSBarry Smith 99c79dcfbdSBarry Smith TEST*/ 100