1c4762a1bSJed Brown static char help[] = "Simple wrapper object to solve DAE of the form:\n\ 2c4762a1bSJed Brown \\dot{U} = f(U,V)\n\ 3c4762a1bSJed Brown F(U,V) = 0\n\n"; 4c4762a1bSJed Brown 5c4762a1bSJed Brown #include <petscts.h> 6c4762a1bSJed Brown 7c4762a1bSJed Brown /* ----------------------------------------------------------------------------*/ 8c4762a1bSJed Brown 9c4762a1bSJed Brown typedef struct _p_TSDAESimple *TSDAESimple; 10c4762a1bSJed Brown struct _p_TSDAESimple { 11c4762a1bSJed Brown MPI_Comm comm; 12ce78bad3SBarry Smith PetscErrorCode (*setfromoptions)(TSDAESimple, PetscOptionItems); 13c4762a1bSJed Brown PetscErrorCode (*solve)(TSDAESimple, Vec); 14c4762a1bSJed Brown PetscErrorCode (*destroy)(TSDAESimple); 15c4762a1bSJed Brown Vec U, V; 16c4762a1bSJed Brown PetscErrorCode (*f)(PetscReal, Vec, Vec, Vec, void *); 17c4762a1bSJed Brown PetscErrorCode (*F)(PetscReal, Vec, Vec, Vec, void *); 18c4762a1bSJed Brown void *fctx, *Fctx; 19c4762a1bSJed Brown void *data; 20c4762a1bSJed Brown }; 21c4762a1bSJed Brown 22d71ae5a4SJacob Faibussowitsch PetscErrorCode TSDAESimpleCreate(MPI_Comm comm, TSDAESimple *tsdae) 23d71ae5a4SJacob Faibussowitsch { 247510d9b0SBarry Smith PetscFunctionBeginUser; 259566063dSJacob Faibussowitsch PetscCall(PetscNew(tsdae)); 26c4762a1bSJed Brown (*tsdae)->comm = comm; 273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 28c4762a1bSJed Brown } 29c4762a1bSJed Brown 30*2a8381b2SBarry Smith PetscErrorCode TSDAESimpleSetRHSFunction(TSDAESimple tsdae, Vec U, PetscErrorCode (*f)(PetscReal, Vec, Vec, Vec, void *), PetscCtx ctx) 31d71ae5a4SJacob Faibussowitsch { 327510d9b0SBarry Smith PetscFunctionBeginUser; 33c4762a1bSJed Brown tsdae->f = f; 34c4762a1bSJed Brown tsdae->U = U; 359566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)U)); 36c4762a1bSJed Brown tsdae->fctx = ctx; 373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 38c4762a1bSJed Brown } 39c4762a1bSJed Brown 40*2a8381b2SBarry Smith PetscErrorCode TSDAESimpleSetIFunction(TSDAESimple tsdae, Vec V, PetscErrorCode (*F)(PetscReal, Vec, Vec, Vec, void *), PetscCtx ctx) 41d71ae5a4SJacob Faibussowitsch { 427510d9b0SBarry Smith PetscFunctionBeginUser; 43c4762a1bSJed Brown tsdae->F = F; 44c4762a1bSJed Brown tsdae->V = V; 459566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)V)); 46c4762a1bSJed Brown tsdae->Fctx = ctx; 473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 48c4762a1bSJed Brown } 49c4762a1bSJed Brown 50d71ae5a4SJacob Faibussowitsch PetscErrorCode TSDAESimpleDestroy(TSDAESimple *tsdae) 51d71ae5a4SJacob Faibussowitsch { 527510d9b0SBarry Smith PetscFunctionBeginUser; 539566063dSJacob Faibussowitsch PetscCall((*(*tsdae)->destroy)(*tsdae)); 549566063dSJacob Faibussowitsch PetscCall(VecDestroy(&(*tsdae)->U)); 559566063dSJacob Faibussowitsch PetscCall(VecDestroy(&(*tsdae)->V)); 569566063dSJacob Faibussowitsch PetscCall(PetscFree(*tsdae)); 573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 58c4762a1bSJed Brown } 59c4762a1bSJed Brown 60d71ae5a4SJacob Faibussowitsch PetscErrorCode TSDAESimpleSolve(TSDAESimple tsdae, Vec Usolution) 61d71ae5a4SJacob Faibussowitsch { 627510d9b0SBarry Smith PetscFunctionBeginUser; 639566063dSJacob Faibussowitsch PetscCall((*tsdae->solve)(tsdae, Usolution)); 643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 65c4762a1bSJed Brown } 66c4762a1bSJed Brown 67d71ae5a4SJacob Faibussowitsch PetscErrorCode TSDAESimpleSetFromOptions(TSDAESimple tsdae) 68d71ae5a4SJacob Faibussowitsch { 697510d9b0SBarry Smith PetscFunctionBeginUser; 709566063dSJacob Faibussowitsch PetscCall((*tsdae->setfromoptions)(PetscOptionsObject, tsdae)); 713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 72c4762a1bSJed Brown } 73c4762a1bSJed Brown 74c4762a1bSJed Brown /* ----------------------------------------------------------------------------*/ 75c4762a1bSJed Brown /* 76c4762a1bSJed Brown Integrates the system by integrating the reduced ODE system and solving the 77c4762a1bSJed Brown algebraic constraints at each stage with a separate SNES solve. 78c4762a1bSJed Brown */ 79c4762a1bSJed Brown 80c4762a1bSJed Brown typedef struct { 81c4762a1bSJed Brown PetscReal t; 82c4762a1bSJed Brown TS ts; 83c4762a1bSJed Brown SNES snes; 84c4762a1bSJed Brown Vec U; 85c4762a1bSJed Brown } TSDAESimple_Reduced; 86c4762a1bSJed Brown 87c4762a1bSJed Brown /* 88c4762a1bSJed Brown Defines the RHS function that is passed to the time-integrator. 89c4762a1bSJed Brown 90c4762a1bSJed Brown Solves F(U,V) for V and then computes f(U,V) 91c4762a1bSJed Brown 92c4762a1bSJed Brown */ 93d71ae5a4SJacob Faibussowitsch PetscErrorCode TSDAESimple_Reduced_TSFunction(TS ts, PetscReal t, Vec U, Vec F, void *actx) 94d71ae5a4SJacob Faibussowitsch { 95c4762a1bSJed Brown TSDAESimple tsdae = (TSDAESimple)actx; 96c4762a1bSJed Brown TSDAESimple_Reduced *red = (TSDAESimple_Reduced *)tsdae->data; 97c4762a1bSJed Brown 98c4762a1bSJed Brown PetscFunctionBeginUser; 99c4762a1bSJed Brown red->t = t; 100c4762a1bSJed Brown red->U = U; 1019566063dSJacob Faibussowitsch PetscCall(SNESSolve(red->snes, NULL, tsdae->V)); 1029566063dSJacob Faibussowitsch PetscCall((*tsdae->f)(t, U, tsdae->V, F, tsdae->fctx)); 1033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 104c4762a1bSJed Brown } 105c4762a1bSJed Brown 106c4762a1bSJed Brown /* 107c4762a1bSJed Brown Defines the nonlinear function that is passed to the nonlinear solver 108c4762a1bSJed Brown 109c4762a1bSJed Brown */ 110d71ae5a4SJacob Faibussowitsch PetscErrorCode TSDAESimple_Reduced_SNESFunction(SNES snes, Vec V, Vec F, void *actx) 111d71ae5a4SJacob Faibussowitsch { 112c4762a1bSJed Brown TSDAESimple tsdae = (TSDAESimple)actx; 113c4762a1bSJed Brown TSDAESimple_Reduced *red = (TSDAESimple_Reduced *)tsdae->data; 114c4762a1bSJed Brown 115c4762a1bSJed Brown PetscFunctionBeginUser; 1169566063dSJacob Faibussowitsch PetscCall((*tsdae->F)(red->t, red->U, V, F, tsdae->Fctx)); 1173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 118c4762a1bSJed Brown } 119c4762a1bSJed Brown 120d71ae5a4SJacob Faibussowitsch PetscErrorCode TSDAESimpleSolve_Reduced(TSDAESimple tsdae, Vec U) 121d71ae5a4SJacob Faibussowitsch { 122c4762a1bSJed Brown TSDAESimple_Reduced *red = (TSDAESimple_Reduced *)tsdae->data; 123c4762a1bSJed Brown 1247510d9b0SBarry Smith PetscFunctionBeginUser; 1259566063dSJacob Faibussowitsch PetscCall(TSSolve(red->ts, U)); 1263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 127c4762a1bSJed Brown } 128c4762a1bSJed Brown 129ce78bad3SBarry Smith PetscErrorCode TSDAESimpleSetFromOptions_Reduced(TSDAESimple tsdae, PetscOptionItems PetscOptionsObject) 130d71ae5a4SJacob Faibussowitsch { 131c4762a1bSJed Brown TSDAESimple_Reduced *red = (TSDAESimple_Reduced *)tsdae->data; 132c4762a1bSJed Brown 1337510d9b0SBarry Smith PetscFunctionBeginUser; 1349566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(red->ts)); 1359566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(red->snes)); 1363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 137c4762a1bSJed Brown } 138c4762a1bSJed Brown 139d71ae5a4SJacob Faibussowitsch PetscErrorCode TSDAESimpleDestroy_Reduced(TSDAESimple tsdae) 140d71ae5a4SJacob Faibussowitsch { 141c4762a1bSJed Brown TSDAESimple_Reduced *red = (TSDAESimple_Reduced *)tsdae->data; 142c4762a1bSJed Brown 1437510d9b0SBarry Smith PetscFunctionBeginUser; 1449566063dSJacob Faibussowitsch PetscCall(TSDestroy(&red->ts)); 1459566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&red->snes)); 1469566063dSJacob Faibussowitsch PetscCall(PetscFree(red)); 1473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 148c4762a1bSJed Brown } 149c4762a1bSJed Brown 150d71ae5a4SJacob Faibussowitsch PetscErrorCode TSDAESimpleSetUp_Reduced(TSDAESimple tsdae) 151d71ae5a4SJacob Faibussowitsch { 152c4762a1bSJed Brown TSDAESimple_Reduced *red; 153c4762a1bSJed Brown Vec tsrhs; 154c4762a1bSJed Brown 1557510d9b0SBarry Smith PetscFunctionBeginUser; 1569566063dSJacob Faibussowitsch PetscCall(PetscNew(&red)); 157c4762a1bSJed Brown tsdae->data = red; 158c4762a1bSJed Brown 159c4762a1bSJed Brown tsdae->setfromoptions = TSDAESimpleSetFromOptions_Reduced; 160c4762a1bSJed Brown tsdae->solve = TSDAESimpleSolve_Reduced; 161c4762a1bSJed Brown tsdae->destroy = TSDAESimpleDestroy_Reduced; 162c4762a1bSJed Brown 1639566063dSJacob Faibussowitsch PetscCall(TSCreate(tsdae->comm, &red->ts)); 1649566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(red->ts, TS_NONLINEAR)); 1659566063dSJacob Faibussowitsch PetscCall(TSSetType(red->ts, TSEULER)); 1669566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(red->ts, TS_EXACTFINALTIME_STEPOVER)); 1679566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tsdae->U, &tsrhs)); 1689566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(red->ts, tsrhs, TSDAESimple_Reduced_TSFunction, tsdae)); 1699566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tsrhs)); 170c4762a1bSJed Brown 1719566063dSJacob Faibussowitsch PetscCall(SNESCreate(tsdae->comm, &red->snes)); 1729566063dSJacob Faibussowitsch PetscCall(SNESSetOptionsPrefix(red->snes, "tsdaesimple_")); 1739566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(red->snes, NULL, TSDAESimple_Reduced_SNESFunction, tsdae)); 1749566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(red->snes, NULL, NULL, SNESComputeJacobianDefault, tsdae)); 1753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 176c4762a1bSJed Brown } 177c4762a1bSJed Brown 178c4762a1bSJed Brown /* ----------------------------------------------------------------------------*/ 179c4762a1bSJed Brown 180c4762a1bSJed Brown /* 181c4762a1bSJed Brown Integrates the system by integrating directly the entire DAE system 182c4762a1bSJed Brown */ 183c4762a1bSJed Brown 184c4762a1bSJed Brown typedef struct { 185c4762a1bSJed Brown TS ts; 186c4762a1bSJed Brown Vec UV, UF, VF; 187c4762a1bSJed Brown VecScatter scatterU, scatterV; 188c4762a1bSJed Brown } TSDAESimple_Full; 189c4762a1bSJed Brown 190c4762a1bSJed Brown /* 191c4762a1bSJed Brown Defines the RHS function that is passed to the time-integrator. 192c4762a1bSJed Brown 193c4762a1bSJed Brown f(U,V) 194c4762a1bSJed Brown 0 195c4762a1bSJed Brown 196c4762a1bSJed Brown */ 197d71ae5a4SJacob Faibussowitsch PetscErrorCode TSDAESimple_Full_TSRHSFunction(TS ts, PetscReal t, Vec UV, Vec F, void *actx) 198d71ae5a4SJacob Faibussowitsch { 199c4762a1bSJed Brown TSDAESimple tsdae = (TSDAESimple)actx; 200c4762a1bSJed Brown TSDAESimple_Full *full = (TSDAESimple_Full *)tsdae->data; 201c4762a1bSJed Brown 202c4762a1bSJed Brown PetscFunctionBeginUser; 2039566063dSJacob Faibussowitsch PetscCall(VecSet(F, 0.0)); 2049566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(full->scatterU, UV, tsdae->U, INSERT_VALUES, SCATTER_REVERSE)); 2059566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(full->scatterU, UV, tsdae->U, INSERT_VALUES, SCATTER_REVERSE)); 2069566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(full->scatterV, UV, tsdae->V, INSERT_VALUES, SCATTER_REVERSE)); 2079566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(full->scatterV, UV, tsdae->V, INSERT_VALUES, SCATTER_REVERSE)); 2089566063dSJacob Faibussowitsch PetscCall((*tsdae->f)(t, tsdae->U, tsdae->V, full->UF, tsdae->fctx)); 2099566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(full->scatterU, full->UF, F, INSERT_VALUES, SCATTER_FORWARD)); 2109566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(full->scatterU, full->UF, F, INSERT_VALUES, SCATTER_FORWARD)); 2113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 212c4762a1bSJed Brown } 213c4762a1bSJed Brown 214c4762a1bSJed Brown /* 215c4762a1bSJed Brown Defines the nonlinear function that is passed to the nonlinear solver 216c4762a1bSJed Brown 217c4762a1bSJed Brown \dot{U} 218c4762a1bSJed Brown F(U,V) 219c4762a1bSJed Brown 220c4762a1bSJed Brown */ 221d71ae5a4SJacob Faibussowitsch PetscErrorCode TSDAESimple_Full_TSIFunction(TS ts, PetscReal t, Vec UV, Vec UVdot, Vec F, void *actx) 222d71ae5a4SJacob Faibussowitsch { 223c4762a1bSJed Brown TSDAESimple tsdae = (TSDAESimple)actx; 224c4762a1bSJed Brown TSDAESimple_Full *full = (TSDAESimple_Full *)tsdae->data; 225c4762a1bSJed Brown 226c4762a1bSJed Brown PetscFunctionBeginUser; 2279566063dSJacob Faibussowitsch PetscCall(VecCopy(UVdot, F)); 2289566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(full->scatterU, UV, tsdae->U, INSERT_VALUES, SCATTER_REVERSE)); 2299566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(full->scatterU, UV, tsdae->U, INSERT_VALUES, SCATTER_REVERSE)); 2309566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(full->scatterV, UV, tsdae->V, INSERT_VALUES, SCATTER_REVERSE)); 2319566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(full->scatterV, UV, tsdae->V, INSERT_VALUES, SCATTER_REVERSE)); 2329566063dSJacob Faibussowitsch PetscCall((*tsdae->F)(t, tsdae->U, tsdae->V, full->VF, tsdae->Fctx)); 2339566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(full->scatterV, full->VF, F, INSERT_VALUES, SCATTER_FORWARD)); 2349566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(full->scatterV, full->VF, F, INSERT_VALUES, SCATTER_FORWARD)); 2353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 236c4762a1bSJed Brown } 237c4762a1bSJed Brown 238d71ae5a4SJacob Faibussowitsch PetscErrorCode TSDAESimpleSolve_Full(TSDAESimple tsdae, Vec U) 239d71ae5a4SJacob Faibussowitsch { 240c4762a1bSJed Brown TSDAESimple_Full *full = (TSDAESimple_Full *)tsdae->data; 241c4762a1bSJed Brown 2427510d9b0SBarry Smith PetscFunctionBeginUser; 2439566063dSJacob Faibussowitsch PetscCall(VecSet(full->UV, 1.0)); 2449566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(full->scatterU, U, full->UV, INSERT_VALUES, SCATTER_FORWARD)); 2459566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(full->scatterU, U, full->UV, INSERT_VALUES, SCATTER_FORWARD)); 2469566063dSJacob Faibussowitsch PetscCall(TSSolve(full->ts, full->UV)); 2479566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(full->scatterU, full->UV, U, INSERT_VALUES, SCATTER_REVERSE)); 2489566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(full->scatterU, full->UV, U, INSERT_VALUES, SCATTER_REVERSE)); 2493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 250c4762a1bSJed Brown } 251c4762a1bSJed Brown 252ce78bad3SBarry Smith PetscErrorCode TSDAESimpleSetFromOptions_Full(TSDAESimple tsdae, PetscOptionItems PetscOptionsObject) 253d71ae5a4SJacob Faibussowitsch { 254c4762a1bSJed Brown TSDAESimple_Full *full = (TSDAESimple_Full *)tsdae->data; 255c4762a1bSJed Brown 2567510d9b0SBarry Smith PetscFunctionBeginUser; 2579566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(full->ts)); 2583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 259c4762a1bSJed Brown } 260c4762a1bSJed Brown 261d71ae5a4SJacob Faibussowitsch PetscErrorCode TSDAESimpleDestroy_Full(TSDAESimple tsdae) 262d71ae5a4SJacob Faibussowitsch { 263c4762a1bSJed Brown TSDAESimple_Full *full = (TSDAESimple_Full *)tsdae->data; 264c4762a1bSJed Brown 2657510d9b0SBarry Smith PetscFunctionBeginUser; 2669566063dSJacob Faibussowitsch PetscCall(TSDestroy(&full->ts)); 2679566063dSJacob Faibussowitsch PetscCall(VecDestroy(&full->UV)); 2689566063dSJacob Faibussowitsch PetscCall(VecDestroy(&full->UF)); 2699566063dSJacob Faibussowitsch PetscCall(VecDestroy(&full->VF)); 2709566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&full->scatterU)); 2719566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&full->scatterV)); 2729566063dSJacob Faibussowitsch PetscCall(PetscFree(full)); 2733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 274c4762a1bSJed Brown } 275c4762a1bSJed Brown 276d71ae5a4SJacob Faibussowitsch PetscErrorCode TSDAESimpleSetUp_Full(TSDAESimple tsdae) 277d71ae5a4SJacob Faibussowitsch { 278c4762a1bSJed Brown TSDAESimple_Full *full; 279c4762a1bSJed Brown Vec tsrhs; 280c4762a1bSJed Brown PetscInt nU, nV, UVstart; 281c4762a1bSJed Brown IS is; 282c4762a1bSJed Brown 2837510d9b0SBarry Smith PetscFunctionBeginUser; 2849566063dSJacob Faibussowitsch PetscCall(PetscNew(&full)); 285c4762a1bSJed Brown tsdae->data = full; 286c4762a1bSJed Brown 287c4762a1bSJed Brown tsdae->setfromoptions = TSDAESimpleSetFromOptions_Full; 288c4762a1bSJed Brown tsdae->solve = TSDAESimpleSolve_Full; 289c4762a1bSJed Brown tsdae->destroy = TSDAESimpleDestroy_Full; 290c4762a1bSJed Brown 2919566063dSJacob Faibussowitsch PetscCall(TSCreate(tsdae->comm, &full->ts)); 2929566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(full->ts, TS_NONLINEAR)); 2939566063dSJacob Faibussowitsch PetscCall(TSSetType(full->ts, TSROSW)); 2949566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(full->ts, TS_EXACTFINALTIME_STEPOVER)); 2959566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tsdae->U, &full->UF)); 2969566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tsdae->V, &full->VF)); 297c4762a1bSJed Brown 2989566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(tsdae->U, &nU)); 2999566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(tsdae->V, &nV)); 30077433607SBarry Smith PetscCall(VecCreateFromOptions(tsdae->comm, NULL, nU + nV, PETSC_DETERMINE, &tsrhs)); 3019566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tsrhs, &full->UV)); 302c4762a1bSJed Brown 3039566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(tsrhs, &UVstart, NULL)); 3049566063dSJacob Faibussowitsch PetscCall(ISCreateStride(tsdae->comm, nU, UVstart, 1, &is)); 3059566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(tsdae->U, NULL, tsrhs, is, &full->scatterU)); 3069566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 3079566063dSJacob Faibussowitsch PetscCall(ISCreateStride(tsdae->comm, nV, UVstart + nU, 1, &is)); 3089566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(tsdae->V, NULL, tsrhs, is, &full->scatterV)); 3099566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 310c4762a1bSJed Brown 3119566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(full->ts, tsrhs, TSDAESimple_Full_TSRHSFunction, tsdae)); 3129566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(full->ts, NULL, TSDAESimple_Full_TSIFunction, tsdae)); 3139566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tsrhs)); 3143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 315c4762a1bSJed Brown } 316c4762a1bSJed Brown 317c4762a1bSJed Brown /* ----------------------------------------------------------------------------*/ 318c4762a1bSJed Brown 319c4762a1bSJed Brown /* 320c4762a1bSJed Brown Simple example: f(U,V) = U + V 321c4762a1bSJed Brown 322c4762a1bSJed Brown */ 323*2a8381b2SBarry Smith PetscErrorCode f(PetscReal t, Vec U, Vec V, Vec F, PetscCtx ctx) 324d71ae5a4SJacob Faibussowitsch { 325c4762a1bSJed Brown PetscFunctionBeginUser; 3269566063dSJacob Faibussowitsch PetscCall(VecWAXPY(F, 1.0, U, V)); 3273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 328c4762a1bSJed Brown } 329c4762a1bSJed Brown 330c4762a1bSJed Brown /* 331c4762a1bSJed Brown Simple example: F(U,V) = U - V 332c4762a1bSJed Brown 333c4762a1bSJed Brown */ 334*2a8381b2SBarry Smith PetscErrorCode F(PetscReal t, Vec U, Vec V, Vec F, PetscCtx ctx) 335d71ae5a4SJacob Faibussowitsch { 336c4762a1bSJed Brown PetscFunctionBeginUser; 3379566063dSJacob Faibussowitsch PetscCall(VecWAXPY(F, -1.0, V, U)); 3383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 339c4762a1bSJed Brown } 340c4762a1bSJed Brown 341d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 342d71ae5a4SJacob Faibussowitsch { 343c4762a1bSJed Brown TSDAESimple tsdae; 344c4762a1bSJed Brown Vec U, V, Usolution; 345c4762a1bSJed Brown 346327415f7SBarry Smith PetscFunctionBeginUser; 347c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 3489566063dSJacob Faibussowitsch PetscCall(TSDAESimpleCreate(PETSC_COMM_WORLD, &tsdae)); 349c4762a1bSJed Brown 35077433607SBarry Smith PetscCall(VecCreateFromOptions(PETSC_COMM_WORLD, NULL, 1, 1, PETSC_DETERMINE, &U)); 35177433607SBarry Smith PetscCall(VecCreateFromOptions(PETSC_COMM_WORLD, NULL, 1, 1, PETSC_DETERMINE, &V)); 3529566063dSJacob Faibussowitsch PetscCall(TSDAESimpleSetRHSFunction(tsdae, U, f, NULL)); 3539566063dSJacob Faibussowitsch PetscCall(TSDAESimpleSetIFunction(tsdae, V, F, NULL)); 354c4762a1bSJed Brown 3559566063dSJacob Faibussowitsch PetscCall(VecDuplicate(U, &Usolution)); 3569566063dSJacob Faibussowitsch PetscCall(VecSet(Usolution, 1.0)); 357c4762a1bSJed Brown 3589566063dSJacob Faibussowitsch /* PetscCall(TSDAESimpleSetUp_Full(tsdae)); */ 3599566063dSJacob Faibussowitsch PetscCall(TSDAESimpleSetUp_Reduced(tsdae)); 360c4762a1bSJed Brown 3619566063dSJacob Faibussowitsch PetscCall(TSDAESimpleSetFromOptions(tsdae)); 3629566063dSJacob Faibussowitsch PetscCall(TSDAESimpleSolve(tsdae, Usolution)); 3639566063dSJacob Faibussowitsch PetscCall(TSDAESimpleDestroy(&tsdae)); 364c4762a1bSJed Brown 3659566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U)); 3669566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Usolution)); 3679566063dSJacob Faibussowitsch PetscCall(VecDestroy(&V)); 3689566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 369b122ec5aSJacob Faibussowitsch return 0; 370c4762a1bSJed Brown } 371