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 *); 8d27334e2SStefano Zampini PetscErrorCode RHSFunction(TS, PetscReal, Vec, Vec, void *); 9c79dcfbdSBarry Smith 10d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 11d71ae5a4SJacob Faibussowitsch { 12c79dcfbdSBarry Smith TS ts; 13c79dcfbdSBarry Smith Vec x; 14c79dcfbdSBarry Smith Mat A; 15d27334e2SStefano Zampini PetscBool flg = PETSC_FALSE, usingimex; 16c79dcfbdSBarry Smith 17327415f7SBarry Smith PetscFunctionBeginUser; 18c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 19c79dcfbdSBarry Smith 209566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 21d27334e2SStefano Zampini PetscCall(PetscOptionsGetBool(NULL, NULL, "-set_implicit", &flg, NULL)); 22d27334e2SStefano Zampini if (flg) PetscCall(TSSetEquationType(ts, TS_EQ_ODE_IMPLICIT)); 23d27334e2SStefano Zampini PetscCall(TSSetIFunction(ts, NULL, IFunction, NULL)); 24d27334e2SStefano Zampini PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, &usingimex)); 25c79dcfbdSBarry Smith 269566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 279566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, 1, 1, PETSC_DECIDE, PETSC_DECIDE)); 289566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 299566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 309566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 319566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 329566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts, A, A, IJacobian, NULL)); 339566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 34c79dcfbdSBarry Smith 359566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &x)); 369566063dSJacob Faibussowitsch PetscCall(VecSetSizes(x, 1, PETSC_DECIDE)); 379566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(x)); 389566063dSJacob Faibussowitsch PetscCall(VecSetUp(x)); 399566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, x)); 409566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 419566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 42c79dcfbdSBarry Smith 43d27334e2SStefano Zampini /* Need to know if we are using an IMEX scheme to decide on the form 44d27334e2SStefano Zampini of the RHS function */ 45d27334e2SStefano Zampini PetscCall(PetscObjectTypeCompare((PetscObject)ts, TSARKIMEX, &usingimex)); 46d27334e2SStefano Zampini if (usingimex) { 47d27334e2SStefano Zampini PetscCall(TSARKIMEXGetFullyImplicit(ts, &flg)); 48d27334e2SStefano Zampini if (flg) usingimex = PETSC_FALSE; 49d27334e2SStefano Zampini } 509566063dSJacob Faibussowitsch PetscCall(TSSetStepNumber(ts, 0)); 519566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, 1)); 529566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts, 0)); 539566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, PETSC_MAX_REAL)); 549566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(ts, 3)); 55c79dcfbdSBarry Smith 56cab2b192SBarry Smith PetscCall(TSSolve(ts, NULL)); 57c79dcfbdSBarry Smith 589566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 599566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 60b122ec5aSJacob Faibussowitsch return 0; 61c79dcfbdSBarry Smith } 62c79dcfbdSBarry Smith 63*2a8381b2SBarry Smith PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec x, Vec f, PetscCtx ctx) 64d27334e2SStefano Zampini { 65d27334e2SStefano Zampini PetscBool usingimex = *(PetscBool *)ctx; 66d27334e2SStefano Zampini 67d27334e2SStefano Zampini PetscFunctionBeginUser; 68d27334e2SStefano Zampini PetscCall(VecSet(f, usingimex ? 0.5 : 1)); 69d27334e2SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 70d27334e2SStefano Zampini } 71d27334e2SStefano Zampini 72*2a8381b2SBarry Smith PetscErrorCode IFunction(TS ts, PetscReal t, Vec x, Vec xdot, Vec f, PetscCtx ctx) 73d71ae5a4SJacob Faibussowitsch { 747510d9b0SBarry Smith PetscFunctionBeginUser; 759566063dSJacob Faibussowitsch PetscCall(VecCopy(xdot, f)); 769566063dSJacob Faibussowitsch PetscCall(VecScale(f, 2.0)); 773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 78c79dcfbdSBarry Smith } 79c79dcfbdSBarry Smith 80*2a8381b2SBarry Smith PetscErrorCode IJacobian(TS ts, PetscReal t, Vec x, Vec xdot, PetscReal shift, Mat A, Mat B, PetscCtx ctx) 81d71ae5a4SJacob Faibussowitsch { 82c79dcfbdSBarry Smith PetscScalar j; 83c79dcfbdSBarry Smith 847510d9b0SBarry Smith PetscFunctionBeginUser; 85c79dcfbdSBarry Smith j = shift * 2.0; 869566063dSJacob Faibussowitsch PetscCall(MatSetValue(B, 0, 0, j, INSERT_VALUES)); 879566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 889566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 90c79dcfbdSBarry Smith } 91c79dcfbdSBarry Smith 92c79dcfbdSBarry Smith /*TEST 93c79dcfbdSBarry Smith 94c79dcfbdSBarry Smith test: 95c79dcfbdSBarry Smith suffix: arkimex_explicit_stage 960ef292d3SStefano Zampini requires: !defined(PETSCTEST_VALGRIND) defined(PETSC_USE_DEBUG) !defined(PETSC_HAVE_SANITIZER) 97d27334e2SStefano Zampini args: -ts_type arkimex -petsc_ci_portable_error_output -error_output_stdout -set_implicit 98e9a33e21SBarry Smith filter: grep -E -v "(memory block|leaked context|not freed before MPI_Finalize|Could be the program crashed)" 99c79dcfbdSBarry Smith 100c79dcfbdSBarry Smith test: 101c79dcfbdSBarry Smith suffix: arkimex_implicit_stage 102d27334e2SStefano Zampini args: -ts_type arkimex -ts_arkimex_type {{3 l2}} -ts_monitor_solution -ts_monitor 103c79dcfbdSBarry Smith 104c79dcfbdSBarry Smith TEST*/ 105