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; 12c4762a1bSJed Brown PetscErrorCode (*setfromoptions)(PetscOptionItems*,TSDAESimple); 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 22c4762a1bSJed Brown PetscErrorCode TSDAESimpleCreate(MPI_Comm comm,TSDAESimple *tsdae) 23c4762a1bSJed Brown { 24c4762a1bSJed Brown PetscErrorCode ierr; 25c4762a1bSJed Brown 26c4762a1bSJed Brown PetscFunctionBegin; 275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscNew(tsdae)); 28c4762a1bSJed Brown (*tsdae)->comm = comm; 29c4762a1bSJed Brown PetscFunctionReturn(0); 30c4762a1bSJed Brown } 31c4762a1bSJed Brown 32c4762a1bSJed Brown PetscErrorCode TSDAESimpleSetRHSFunction(TSDAESimple tsdae,Vec U,PetscErrorCode (*f)(PetscReal,Vec,Vec,Vec,void*),void *ctx) 33c4762a1bSJed Brown { 34c4762a1bSJed Brown PetscErrorCode ierr; 35c4762a1bSJed Brown 36c4762a1bSJed Brown PetscFunctionBegin; 37c4762a1bSJed Brown tsdae->f = f; 38c4762a1bSJed Brown tsdae->U = U; 395f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)U)); 40c4762a1bSJed Brown tsdae->fctx = ctx; 41c4762a1bSJed Brown PetscFunctionReturn(0); 42c4762a1bSJed Brown } 43c4762a1bSJed Brown 44c4762a1bSJed Brown PetscErrorCode TSDAESimpleSetIFunction(TSDAESimple tsdae,Vec V,PetscErrorCode (*F)(PetscReal,Vec,Vec,Vec,void*),void *ctx) 45c4762a1bSJed Brown { 46c4762a1bSJed Brown PetscErrorCode ierr; 47c4762a1bSJed Brown 48c4762a1bSJed Brown PetscFunctionBegin; 49c4762a1bSJed Brown tsdae->F = F; 50c4762a1bSJed Brown tsdae->V = V; 515f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)V)); 52c4762a1bSJed Brown tsdae->Fctx = ctx; 53c4762a1bSJed Brown PetscFunctionReturn(0); 54c4762a1bSJed Brown } 55c4762a1bSJed Brown 56c4762a1bSJed Brown PetscErrorCode TSDAESimpleDestroy(TSDAESimple *tsdae) 57c4762a1bSJed Brown { 58c4762a1bSJed Brown PetscErrorCode ierr; 59c4762a1bSJed Brown 60c4762a1bSJed Brown PetscFunctionBegin; 615f80ce2aSJacob Faibussowitsch CHKERRQ((*(*tsdae)->destroy)(*tsdae)); 625f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&(*tsdae)->U)); 635f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&(*tsdae)->V)); 645f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(*tsdae)); 65c4762a1bSJed Brown PetscFunctionReturn(0); 66c4762a1bSJed Brown } 67c4762a1bSJed Brown 68c4762a1bSJed Brown PetscErrorCode TSDAESimpleSolve(TSDAESimple tsdae,Vec Usolution) 69c4762a1bSJed Brown { 70c4762a1bSJed Brown PetscErrorCode ierr; 71c4762a1bSJed Brown 72c4762a1bSJed Brown PetscFunctionBegin; 735f80ce2aSJacob Faibussowitsch CHKERRQ((*tsdae->solve)(tsdae,Usolution)); 74c4762a1bSJed Brown PetscFunctionReturn(0); 75c4762a1bSJed Brown } 76c4762a1bSJed Brown 77c4762a1bSJed Brown PetscErrorCode TSDAESimpleSetFromOptions(TSDAESimple tsdae) 78c4762a1bSJed Brown { 79c4762a1bSJed Brown PetscErrorCode ierr; 80c4762a1bSJed Brown 81c4762a1bSJed Brown PetscFunctionBegin; 825f80ce2aSJacob Faibussowitsch CHKERRQ((*tsdae->setfromoptions)(PetscOptionsObject,tsdae)); 83c4762a1bSJed Brown PetscFunctionReturn(0); 84c4762a1bSJed Brown } 85c4762a1bSJed Brown 86c4762a1bSJed Brown /* ----------------------------------------------------------------------------*/ 87c4762a1bSJed Brown /* 88c4762a1bSJed Brown Integrates the system by integrating the reduced ODE system and solving the 89c4762a1bSJed Brown algebraic constraints at each stage with a separate SNES solve. 90c4762a1bSJed Brown */ 91c4762a1bSJed Brown 92c4762a1bSJed Brown typedef struct { 93c4762a1bSJed Brown PetscReal t; 94c4762a1bSJed Brown TS ts; 95c4762a1bSJed Brown SNES snes; 96c4762a1bSJed Brown Vec U; 97c4762a1bSJed Brown } TSDAESimple_Reduced; 98c4762a1bSJed Brown 99c4762a1bSJed Brown /* 100c4762a1bSJed Brown Defines the RHS function that is passed to the time-integrator. 101c4762a1bSJed Brown 102c4762a1bSJed Brown Solves F(U,V) for V and then computes f(U,V) 103c4762a1bSJed Brown 104c4762a1bSJed Brown */ 105c4762a1bSJed Brown PetscErrorCode TSDAESimple_Reduced_TSFunction(TS ts,PetscReal t,Vec U,Vec F,void *actx) 106c4762a1bSJed Brown { 107c4762a1bSJed Brown TSDAESimple tsdae = (TSDAESimple)actx; 108c4762a1bSJed Brown TSDAESimple_Reduced *red = (TSDAESimple_Reduced*)tsdae->data; 109c4762a1bSJed Brown PetscErrorCode ierr; 110c4762a1bSJed Brown 111c4762a1bSJed Brown PetscFunctionBeginUser; 112c4762a1bSJed Brown red->t = t; 113c4762a1bSJed Brown red->U = U; 1145f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSolve(red->snes,NULL,tsdae->V)); 1155f80ce2aSJacob Faibussowitsch CHKERRQ((*tsdae->f)(t,U,tsdae->V,F,tsdae->fctx)); 116c4762a1bSJed Brown PetscFunctionReturn(0); 117c4762a1bSJed Brown } 118c4762a1bSJed Brown 119c4762a1bSJed Brown /* 120c4762a1bSJed Brown Defines the nonlinear function that is passed to the nonlinear solver 121c4762a1bSJed Brown 122c4762a1bSJed Brown */ 123c4762a1bSJed Brown PetscErrorCode TSDAESimple_Reduced_SNESFunction(SNES snes,Vec V,Vec F,void *actx) 124c4762a1bSJed Brown { 125c4762a1bSJed Brown TSDAESimple tsdae = (TSDAESimple)actx; 126c4762a1bSJed Brown TSDAESimple_Reduced *red = (TSDAESimple_Reduced*)tsdae->data; 127c4762a1bSJed Brown PetscErrorCode ierr; 128c4762a1bSJed Brown 129c4762a1bSJed Brown PetscFunctionBeginUser; 1305f80ce2aSJacob Faibussowitsch CHKERRQ((*tsdae->F)(red->t,red->U,V,F,tsdae->Fctx)); 131c4762a1bSJed Brown PetscFunctionReturn(0); 132c4762a1bSJed Brown } 133c4762a1bSJed Brown 134c4762a1bSJed Brown PetscErrorCode TSDAESimpleSolve_Reduced(TSDAESimple tsdae,Vec U) 135c4762a1bSJed Brown { 136c4762a1bSJed Brown PetscErrorCode ierr; 137c4762a1bSJed Brown TSDAESimple_Reduced *red = (TSDAESimple_Reduced*)tsdae->data; 138c4762a1bSJed Brown 139c4762a1bSJed Brown PetscFunctionBegin; 1405f80ce2aSJacob Faibussowitsch CHKERRQ(TSSolve(red->ts,U)); 141c4762a1bSJed Brown PetscFunctionReturn(0); 142c4762a1bSJed Brown } 143c4762a1bSJed Brown 144c4762a1bSJed Brown PetscErrorCode TSDAESimpleSetFromOptions_Reduced(PetscOptionItems *PetscOptionsObject,TSDAESimple tsdae) 145c4762a1bSJed Brown { 146c4762a1bSJed Brown PetscErrorCode ierr; 147c4762a1bSJed Brown TSDAESimple_Reduced *red = (TSDAESimple_Reduced*)tsdae->data; 148c4762a1bSJed Brown 149c4762a1bSJed Brown PetscFunctionBegin; 1505f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(red->ts)); 1515f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetFromOptions(red->snes)); 152c4762a1bSJed Brown PetscFunctionReturn(0); 153c4762a1bSJed Brown } 154c4762a1bSJed Brown 155c4762a1bSJed Brown PetscErrorCode TSDAESimpleDestroy_Reduced(TSDAESimple tsdae) 156c4762a1bSJed Brown { 157c4762a1bSJed Brown PetscErrorCode ierr; 158c4762a1bSJed Brown TSDAESimple_Reduced *red = (TSDAESimple_Reduced*)tsdae->data; 159c4762a1bSJed Brown 160c4762a1bSJed Brown PetscFunctionBegin; 1615f80ce2aSJacob Faibussowitsch CHKERRQ(TSDestroy(&red->ts)); 1625f80ce2aSJacob Faibussowitsch CHKERRQ(SNESDestroy(&red->snes)); 1635f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(red)); 164c4762a1bSJed Brown PetscFunctionReturn(0); 165c4762a1bSJed Brown } 166c4762a1bSJed Brown 167c4762a1bSJed Brown PetscErrorCode TSDAESimpleSetUp_Reduced(TSDAESimple tsdae) 168c4762a1bSJed Brown { 169c4762a1bSJed Brown PetscErrorCode ierr; 170c4762a1bSJed Brown TSDAESimple_Reduced *red; 171c4762a1bSJed Brown Vec tsrhs; 172c4762a1bSJed Brown 173c4762a1bSJed Brown PetscFunctionBegin; 1745f80ce2aSJacob Faibussowitsch CHKERRQ(PetscNew(&red)); 175c4762a1bSJed Brown tsdae->data = red; 176c4762a1bSJed Brown 177c4762a1bSJed Brown tsdae->setfromoptions = TSDAESimpleSetFromOptions_Reduced; 178c4762a1bSJed Brown tsdae->solve = TSDAESimpleSolve_Reduced; 179c4762a1bSJed Brown tsdae->destroy = TSDAESimpleDestroy_Reduced; 180c4762a1bSJed Brown 1815f80ce2aSJacob Faibussowitsch CHKERRQ(TSCreate(tsdae->comm,&red->ts)); 1825f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetProblemType(red->ts,TS_NONLINEAR)); 1835f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetType(red->ts,TSEULER)); 1845f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetExactFinalTime(red->ts,TS_EXACTFINALTIME_STEPOVER)); 1855f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(tsdae->U,&tsrhs)); 1865f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSFunction(red->ts,tsrhs,TSDAESimple_Reduced_TSFunction,tsdae)); 1875f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&tsrhs)); 188c4762a1bSJed Brown 1895f80ce2aSJacob Faibussowitsch CHKERRQ(SNESCreate(tsdae->comm,&red->snes)); 1905f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetOptionsPrefix(red->snes,"tsdaesimple_")); 1915f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetFunction(red->snes,NULL,TSDAESimple_Reduced_SNESFunction,tsdae)); 1925f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetJacobian(red->snes,NULL,NULL,SNESComputeJacobianDefault,tsdae)); 193c4762a1bSJed Brown PetscFunctionReturn(0); 194c4762a1bSJed Brown } 195c4762a1bSJed Brown 196c4762a1bSJed Brown /* ----------------------------------------------------------------------------*/ 197c4762a1bSJed Brown 198c4762a1bSJed Brown /* 199c4762a1bSJed Brown Integrates the system by integrating directly the entire DAE system 200c4762a1bSJed Brown */ 201c4762a1bSJed Brown 202c4762a1bSJed Brown typedef struct { 203c4762a1bSJed Brown TS ts; 204c4762a1bSJed Brown Vec UV,UF,VF; 205c4762a1bSJed Brown VecScatter scatterU,scatterV; 206c4762a1bSJed Brown } TSDAESimple_Full; 207c4762a1bSJed Brown 208c4762a1bSJed Brown /* 209c4762a1bSJed Brown Defines the RHS function that is passed to the time-integrator. 210c4762a1bSJed Brown 211c4762a1bSJed Brown f(U,V) 212c4762a1bSJed Brown 0 213c4762a1bSJed Brown 214c4762a1bSJed Brown */ 215c4762a1bSJed Brown PetscErrorCode TSDAESimple_Full_TSRHSFunction(TS ts,PetscReal t,Vec UV,Vec F,void *actx) 216c4762a1bSJed Brown { 217c4762a1bSJed Brown TSDAESimple tsdae = (TSDAESimple)actx; 218c4762a1bSJed Brown TSDAESimple_Full *full = (TSDAESimple_Full*)tsdae->data; 219c4762a1bSJed Brown PetscErrorCode ierr; 220c4762a1bSJed Brown 221c4762a1bSJed Brown PetscFunctionBeginUser; 2225f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(F,0.0)); 2235f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(full->scatterU,UV,tsdae->U,INSERT_VALUES,SCATTER_REVERSE)); 2245f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(full->scatterU,UV,tsdae->U,INSERT_VALUES,SCATTER_REVERSE)); 2255f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(full->scatterV,UV,tsdae->V,INSERT_VALUES,SCATTER_REVERSE)); 2265f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(full->scatterV,UV,tsdae->V,INSERT_VALUES,SCATTER_REVERSE)); 2275f80ce2aSJacob Faibussowitsch CHKERRQ((*tsdae->f)(t,tsdae->U,tsdae->V,full->UF,tsdae->fctx)); 2285f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(full->scatterU,full->UF,F,INSERT_VALUES,SCATTER_FORWARD)); 2295f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(full->scatterU,full->UF,F,INSERT_VALUES,SCATTER_FORWARD)); 230c4762a1bSJed Brown PetscFunctionReturn(0); 231c4762a1bSJed Brown } 232c4762a1bSJed Brown 233c4762a1bSJed Brown /* 234c4762a1bSJed Brown Defines the nonlinear function that is passed to the nonlinear solver 235c4762a1bSJed Brown 236c4762a1bSJed Brown \dot{U} 237c4762a1bSJed Brown F(U,V) 238c4762a1bSJed Brown 239c4762a1bSJed Brown */ 240c4762a1bSJed Brown PetscErrorCode TSDAESimple_Full_TSIFunction(TS ts,PetscReal t,Vec UV,Vec UVdot,Vec F,void *actx) 241c4762a1bSJed Brown { 242c4762a1bSJed Brown TSDAESimple tsdae = (TSDAESimple)actx; 243c4762a1bSJed Brown TSDAESimple_Full *full = (TSDAESimple_Full*)tsdae->data; 244c4762a1bSJed Brown PetscErrorCode ierr; 245c4762a1bSJed Brown 246c4762a1bSJed Brown PetscFunctionBeginUser; 2475f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(UVdot,F)); 2485f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(full->scatterU,UV,tsdae->U,INSERT_VALUES,SCATTER_REVERSE)); 2495f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(full->scatterU,UV,tsdae->U,INSERT_VALUES,SCATTER_REVERSE)); 2505f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(full->scatterV,UV,tsdae->V,INSERT_VALUES,SCATTER_REVERSE)); 2515f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(full->scatterV,UV,tsdae->V,INSERT_VALUES,SCATTER_REVERSE)); 2525f80ce2aSJacob Faibussowitsch CHKERRQ((*tsdae->F)(t,tsdae->U,tsdae->V,full->VF,tsdae->Fctx)); 2535f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(full->scatterV,full->VF,F,INSERT_VALUES,SCATTER_FORWARD)); 2545f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(full->scatterV,full->VF,F,INSERT_VALUES,SCATTER_FORWARD)); 255c4762a1bSJed Brown PetscFunctionReturn(0); 256c4762a1bSJed Brown } 257c4762a1bSJed Brown 258c4762a1bSJed Brown PetscErrorCode TSDAESimpleSolve_Full(TSDAESimple tsdae,Vec U) 259c4762a1bSJed Brown { 260c4762a1bSJed Brown PetscErrorCode ierr; 261c4762a1bSJed Brown TSDAESimple_Full *full = (TSDAESimple_Full*)tsdae->data; 262c4762a1bSJed Brown 263c4762a1bSJed Brown PetscFunctionBegin; 2645f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(full->UV,1.0)); 2655f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(full->scatterU,U,full->UV,INSERT_VALUES,SCATTER_FORWARD)); 2665f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(full->scatterU,U,full->UV,INSERT_VALUES,SCATTER_FORWARD)); 2675f80ce2aSJacob Faibussowitsch CHKERRQ(TSSolve(full->ts,full->UV)); 2685f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(full->scatterU,full->UV,U,INSERT_VALUES,SCATTER_REVERSE)); 2695f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(full->scatterU,full->UV,U,INSERT_VALUES,SCATTER_REVERSE)); 270c4762a1bSJed Brown PetscFunctionReturn(0); 271c4762a1bSJed Brown } 272c4762a1bSJed Brown 273c4762a1bSJed Brown PetscErrorCode TSDAESimpleSetFromOptions_Full(PetscOptionItems *PetscOptionsObject,TSDAESimple tsdae) 274c4762a1bSJed Brown { 275c4762a1bSJed Brown PetscErrorCode ierr; 276c4762a1bSJed Brown TSDAESimple_Full *full = (TSDAESimple_Full*)tsdae->data; 277c4762a1bSJed Brown 278c4762a1bSJed Brown PetscFunctionBegin; 2795f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(full->ts)); 280c4762a1bSJed Brown PetscFunctionReturn(0); 281c4762a1bSJed Brown } 282c4762a1bSJed Brown 283c4762a1bSJed Brown PetscErrorCode TSDAESimpleDestroy_Full(TSDAESimple tsdae) 284c4762a1bSJed Brown { 285c4762a1bSJed Brown PetscErrorCode ierr; 286c4762a1bSJed Brown TSDAESimple_Full *full = (TSDAESimple_Full*)tsdae->data; 287c4762a1bSJed Brown 288c4762a1bSJed Brown PetscFunctionBegin; 2895f80ce2aSJacob Faibussowitsch CHKERRQ(TSDestroy(&full->ts)); 2905f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&full->UV)); 2915f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&full->UF)); 2925f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&full->VF)); 2935f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterDestroy(&full->scatterU)); 2945f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterDestroy(&full->scatterV)); 2955f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(full)); 296c4762a1bSJed Brown PetscFunctionReturn(0); 297c4762a1bSJed Brown } 298c4762a1bSJed Brown 299c4762a1bSJed Brown PetscErrorCode TSDAESimpleSetUp_Full(TSDAESimple tsdae) 300c4762a1bSJed Brown { 301c4762a1bSJed Brown PetscErrorCode ierr; 302c4762a1bSJed Brown TSDAESimple_Full *full; 303c4762a1bSJed Brown Vec tsrhs; 304c4762a1bSJed Brown PetscInt nU,nV,UVstart; 305c4762a1bSJed Brown IS is; 306c4762a1bSJed Brown 307c4762a1bSJed Brown PetscFunctionBegin; 3085f80ce2aSJacob Faibussowitsch CHKERRQ(PetscNew(&full)); 309c4762a1bSJed Brown tsdae->data = full; 310c4762a1bSJed Brown 311c4762a1bSJed Brown tsdae->setfromoptions = TSDAESimpleSetFromOptions_Full; 312c4762a1bSJed Brown tsdae->solve = TSDAESimpleSolve_Full; 313c4762a1bSJed Brown tsdae->destroy = TSDAESimpleDestroy_Full; 314c4762a1bSJed Brown 3155f80ce2aSJacob Faibussowitsch CHKERRQ(TSCreate(tsdae->comm,&full->ts)); 3165f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetProblemType(full->ts,TS_NONLINEAR)); 3175f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetType(full->ts,TSROSW)); 3185f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetExactFinalTime(full->ts,TS_EXACTFINALTIME_STEPOVER)); 3195f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(tsdae->U,&full->UF)); 3205f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(tsdae->V,&full->VF)); 321c4762a1bSJed Brown 3225f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(tsdae->U,&nU)); 3235f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(tsdae->V,&nV)); 3245f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateMPI(tsdae->comm,nU+nV,PETSC_DETERMINE,&tsrhs)); 3255f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(tsrhs,&full->UV)); 326c4762a1bSJed Brown 3275f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetOwnershipRange(tsrhs,&UVstart,NULL)); 3285f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateStride(tsdae->comm,nU,UVstart,1,&is)); 3295f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterCreate(tsdae->U,NULL,tsrhs,is,&full->scatterU)); 3305f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&is)); 3315f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateStride(tsdae->comm,nV,UVstart+nU,1,&is)); 3325f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterCreate(tsdae->V,NULL,tsrhs,is,&full->scatterV)); 3335f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&is)); 334c4762a1bSJed Brown 3355f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSFunction(full->ts,tsrhs,TSDAESimple_Full_TSRHSFunction,tsdae)); 3365f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetIFunction(full->ts,NULL,TSDAESimple_Full_TSIFunction,tsdae)); 3375f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&tsrhs)); 338c4762a1bSJed Brown PetscFunctionReturn(0); 339c4762a1bSJed Brown } 340c4762a1bSJed Brown 341c4762a1bSJed Brown /* ----------------------------------------------------------------------------*/ 342c4762a1bSJed Brown 343c4762a1bSJed Brown /* 344c4762a1bSJed Brown Simple example: f(U,V) = U + V 345c4762a1bSJed Brown 346c4762a1bSJed Brown */ 347c4762a1bSJed Brown PetscErrorCode f(PetscReal t,Vec U,Vec V,Vec F,void *ctx) 348c4762a1bSJed Brown { 349c4762a1bSJed Brown PetscErrorCode ierr; 350c4762a1bSJed Brown 351c4762a1bSJed Brown PetscFunctionBeginUser; 3525f80ce2aSJacob Faibussowitsch CHKERRQ(VecWAXPY(F,1.0,U,V)); 353c4762a1bSJed Brown PetscFunctionReturn(0); 354c4762a1bSJed Brown } 355c4762a1bSJed Brown 356c4762a1bSJed Brown /* 357c4762a1bSJed Brown Simple example: F(U,V) = U - V 358c4762a1bSJed Brown 359c4762a1bSJed Brown */ 360c4762a1bSJed Brown PetscErrorCode F(PetscReal t,Vec U,Vec V,Vec F,void *ctx) 361c4762a1bSJed Brown { 362c4762a1bSJed Brown PetscErrorCode ierr; 363c4762a1bSJed Brown 364c4762a1bSJed Brown PetscFunctionBeginUser; 3655f80ce2aSJacob Faibussowitsch CHKERRQ(VecWAXPY(F,-1.0,V,U)); 366c4762a1bSJed Brown PetscFunctionReturn(0); 367c4762a1bSJed Brown } 368c4762a1bSJed Brown 369c4762a1bSJed Brown int main(int argc,char **argv) 370c4762a1bSJed Brown { 371c4762a1bSJed Brown PetscErrorCode ierr; 372c4762a1bSJed Brown TSDAESimple tsdae; 373c4762a1bSJed Brown Vec U,V,Usolution; 374c4762a1bSJed Brown 375*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help)); 3765f80ce2aSJacob Faibussowitsch CHKERRQ(TSDAESimpleCreate(PETSC_COMM_WORLD,&tsdae)); 377c4762a1bSJed Brown 3785f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateMPI(PETSC_COMM_WORLD,1,PETSC_DETERMINE,&U)); 3795f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateMPI(PETSC_COMM_WORLD,1,PETSC_DETERMINE,&V)); 3805f80ce2aSJacob Faibussowitsch CHKERRQ(TSDAESimpleSetRHSFunction(tsdae,U,f,NULL)); 3815f80ce2aSJacob Faibussowitsch CHKERRQ(TSDAESimpleSetIFunction(tsdae,V,F,NULL)); 382c4762a1bSJed Brown 3835f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(U,&Usolution)); 3845f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(Usolution,1.0)); 385c4762a1bSJed Brown 3865f80ce2aSJacob Faibussowitsch /* CHKERRQ(TSDAESimpleSetUp_Full(tsdae)); */ 3875f80ce2aSJacob Faibussowitsch CHKERRQ(TSDAESimpleSetUp_Reduced(tsdae)); 388c4762a1bSJed Brown 3895f80ce2aSJacob Faibussowitsch CHKERRQ(TSDAESimpleSetFromOptions(tsdae)); 3905f80ce2aSJacob Faibussowitsch CHKERRQ(TSDAESimpleSolve(tsdae,Usolution)); 3915f80ce2aSJacob Faibussowitsch CHKERRQ(TSDAESimpleDestroy(&tsdae)); 392c4762a1bSJed Brown 3935f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&U)); 3945f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&Usolution)); 3955f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&V)); 396*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 397*b122ec5aSJacob Faibussowitsch return 0; 398c4762a1bSJed Brown } 399