xref: /petsc/src/ts/tests/ex26.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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 
17*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
18c79dcfbdSBarry Smith 
195f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts));
205f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetEquationType(ts,TS_EQ_ODE_IMPLICIT));
215f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&f));
225f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(f,1,PETSC_DECIDE));
235f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(f));
245f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetUp(f));
255f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetIFunction(ts,f,IFunction,NULL));
265f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&f));
27c79dcfbdSBarry Smith 
285f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
295f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A,1,1,PETSC_DECIDE,PETSC_DECIDE));
305f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(A));
315f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(A));
325f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
335f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
34c79dcfbdSBarry Smith   /* ensure that the Jacobian matrix has diagonal entries since that is required by TS */
355f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShift(A,(PetscReal)1));
365f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShift(A,(PetscReal)-1));
375f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetIJacobian(ts,A,A,IJacobian,NULL));
385f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
39c79dcfbdSBarry Smith 
405f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&x));
415f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(x,1,PETSC_DECIDE));
425f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(x));
435f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetUp(x));
445f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetSolution(ts,x));
455f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
465f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
47c79dcfbdSBarry Smith 
485f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetStepNumber(ts,0));
495f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(ts,1));
505f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTime(ts,0));
515f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxTime(ts,PETSC_MAX_REAL));
525f80ce2aSJacob 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 
615f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&ts));
62*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
63*b122ec5aSJacob Faibussowitsch   return 0;
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;
695f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCopy(xdot,f));
705f80ce2aSJacob Faibussowitsch   CHKERRQ(VecScale(f,2.0));
715f80ce2aSJacob 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;
815f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValue(B,0,0,j,INSERT_VALUES));
825f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
835f80ce2aSJacob 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