16aad120cSJose E. Roman static char help[] = "Solves the trivial 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 99371c9d4SSatish Balay int main(int argc, char **argv) { 10c79dcfbdSBarry Smith TS ts; 11c79dcfbdSBarry Smith Vec x; 12c79dcfbdSBarry Smith Vec f; 13c79dcfbdSBarry Smith Mat A; 14c79dcfbdSBarry Smith 15327415f7SBarry Smith PetscFunctionBeginUser; 169566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 17c79dcfbdSBarry Smith 189566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 199566063dSJacob Faibussowitsch PetscCall(TSSetEquationType(ts, TS_EQ_ODE_IMPLICIT)); 209566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &f)); 219566063dSJacob Faibussowitsch PetscCall(VecSetSizes(f, 1, PETSC_DECIDE)); 229566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(f)); 239566063dSJacob Faibussowitsch PetscCall(VecSetUp(f)); 249566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts, f, IFunction, NULL)); 259566063dSJacob Faibussowitsch PetscCall(VecDestroy(&f)); 26c79dcfbdSBarry Smith 279566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 289566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, 1, 1, PETSC_DECIDE, PETSC_DECIDE)); 299566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 309566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 319566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 329566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 33c79dcfbdSBarry Smith /* ensure that the Jacobian matrix has diagonal entries since that is required by TS */ 349566063dSJacob Faibussowitsch PetscCall(MatShift(A, (PetscReal)1)); 359566063dSJacob Faibussowitsch PetscCall(MatShift(A, (PetscReal)-1)); 369566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts, A, A, IJacobian, NULL)); 379566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 38c79dcfbdSBarry Smith 399566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &x)); 409566063dSJacob Faibussowitsch PetscCall(VecSetSizes(x, 1, PETSC_DECIDE)); 419566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(x)); 429566063dSJacob Faibussowitsch PetscCall(VecSetUp(x)); 439566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, x)); 449566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 459566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 46c79dcfbdSBarry Smith 479566063dSJacob Faibussowitsch PetscCall(TSSetStepNumber(ts, 0)); 489566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, 1)); 499566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts, 0)); 509566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, PETSC_MAX_REAL)); 519566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(ts, 3)); 52c79dcfbdSBarry Smith 53c79dcfbdSBarry Smith /* 54c79dcfbdSBarry 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 55c79dcfbdSBarry Smith a non-trivial mass matrix with ARKIMEX schemes with explicit stages. 56c79dcfbdSBarry Smith */ 57cab2b192SBarry Smith PetscCall(TSSolve(ts, NULL)); 58c79dcfbdSBarry Smith 599566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 609566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 61b122ec5aSJacob Faibussowitsch return 0; 62c79dcfbdSBarry Smith } 63c79dcfbdSBarry Smith 649371c9d4SSatish Balay PetscErrorCode IFunction(TS ts, PetscReal t, Vec x, Vec xdot, Vec f, void *ctx) { 657510d9b0SBarry Smith PetscFunctionBeginUser; 669566063dSJacob Faibussowitsch PetscCall(VecCopy(xdot, f)); 679566063dSJacob Faibussowitsch PetscCall(VecScale(f, 2.0)); 689566063dSJacob Faibussowitsch PetscCall(VecShift(f, -1.0)); 69c79dcfbdSBarry Smith PetscFunctionReturn(0); 70c79dcfbdSBarry Smith } 71c79dcfbdSBarry Smith 729371c9d4SSatish Balay PetscErrorCode IJacobian(TS ts, PetscReal t, Vec x, Vec xdot, PetscReal shift, Mat A, Mat B, void *ctx) { 73c79dcfbdSBarry Smith PetscScalar j; 74c79dcfbdSBarry Smith 757510d9b0SBarry Smith PetscFunctionBeginUser; 76c79dcfbdSBarry Smith j = shift * 2.0; 779566063dSJacob Faibussowitsch PetscCall(MatSetValue(B, 0, 0, j, INSERT_VALUES)); 789566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 799566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 80c79dcfbdSBarry Smith PetscFunctionReturn(0); 81c79dcfbdSBarry Smith } 82c79dcfbdSBarry Smith 83c79dcfbdSBarry Smith /*TEST 84c79dcfbdSBarry Smith 85c79dcfbdSBarry Smith test: 86c79dcfbdSBarry Smith suffix: arkimex_explicit_stage 87cab2b192SBarry Smith requires: !defined(PETSCTEST_VALGRIND) defined(PETSC_USE_DEBUG) 88cab2b192SBarry Smith args: -ts_type arkimex -petsc_ci_portable_error_output -error_output_stdout 89*36f0be53SJacob Faibussowitsch filter: egrep -v "(options_left|memory block|leaked context|not freed before MPI_Finalize|Could be the program crashed)" 90c79dcfbdSBarry Smith 91c79dcfbdSBarry Smith test: 92c79dcfbdSBarry Smith suffix: arkimex_implicit_stage 93c79dcfbdSBarry Smith args: -ts_type arkimex -ts_arkimex_type l2 -ts_monitor_solution -ts_monitor 94c79dcfbdSBarry Smith 95c79dcfbdSBarry Smith TEST*/ 96