1*c4762a1bSJed Brown static char help[] = "Simple wrapper object to solve DAE of the form:\n\ 2*c4762a1bSJed Brown \\dot{U} = f(U,V)\n\ 3*c4762a1bSJed Brown F(U,V) = 0\n\n"; 4*c4762a1bSJed Brown 5*c4762a1bSJed Brown #include <petscts.h> 6*c4762a1bSJed Brown 7*c4762a1bSJed Brown /* ----------------------------------------------------------------------------*/ 8*c4762a1bSJed Brown 9*c4762a1bSJed Brown typedef struct _p_TSDAESimple *TSDAESimple; 10*c4762a1bSJed Brown struct _p_TSDAESimple { 11*c4762a1bSJed Brown MPI_Comm comm; 12*c4762a1bSJed Brown PetscErrorCode (*setfromoptions)(PetscOptionItems*,TSDAESimple); 13*c4762a1bSJed Brown PetscErrorCode (*solve)(TSDAESimple,Vec); 14*c4762a1bSJed Brown PetscErrorCode (*destroy)(TSDAESimple); 15*c4762a1bSJed Brown Vec U,V; 16*c4762a1bSJed Brown PetscErrorCode (*f)(PetscReal,Vec,Vec,Vec,void*); 17*c4762a1bSJed Brown PetscErrorCode (*F)(PetscReal,Vec,Vec,Vec,void*); 18*c4762a1bSJed Brown void *fctx,*Fctx; 19*c4762a1bSJed Brown void *data; 20*c4762a1bSJed Brown }; 21*c4762a1bSJed Brown 22*c4762a1bSJed Brown PetscErrorCode TSDAESimpleCreate(MPI_Comm comm,TSDAESimple *tsdae) 23*c4762a1bSJed Brown { 24*c4762a1bSJed Brown PetscErrorCode ierr; 25*c4762a1bSJed Brown 26*c4762a1bSJed Brown PetscFunctionBegin; 27*c4762a1bSJed Brown ierr = PetscNew(tsdae);CHKERRQ(ierr); 28*c4762a1bSJed Brown (*tsdae)->comm = comm; 29*c4762a1bSJed Brown PetscFunctionReturn(0); 30*c4762a1bSJed Brown } 31*c4762a1bSJed Brown 32*c4762a1bSJed Brown PetscErrorCode TSDAESimpleSetRHSFunction(TSDAESimple tsdae,Vec U,PetscErrorCode (*f)(PetscReal,Vec,Vec,Vec,void*),void *ctx) 33*c4762a1bSJed Brown { 34*c4762a1bSJed Brown PetscErrorCode ierr; 35*c4762a1bSJed Brown 36*c4762a1bSJed Brown PetscFunctionBegin; 37*c4762a1bSJed Brown tsdae->f = f; 38*c4762a1bSJed Brown tsdae->U = U; 39*c4762a1bSJed Brown ierr = PetscObjectReference((PetscObject)U);CHKERRQ(ierr); 40*c4762a1bSJed Brown tsdae->fctx = ctx; 41*c4762a1bSJed Brown PetscFunctionReturn(0); 42*c4762a1bSJed Brown } 43*c4762a1bSJed Brown 44*c4762a1bSJed Brown PetscErrorCode TSDAESimpleSetIFunction(TSDAESimple tsdae,Vec V,PetscErrorCode (*F)(PetscReal,Vec,Vec,Vec,void*),void *ctx) 45*c4762a1bSJed Brown { 46*c4762a1bSJed Brown PetscErrorCode ierr; 47*c4762a1bSJed Brown 48*c4762a1bSJed Brown PetscFunctionBegin; 49*c4762a1bSJed Brown tsdae->F = F; 50*c4762a1bSJed Brown tsdae->V = V; 51*c4762a1bSJed Brown ierr = PetscObjectReference((PetscObject)V);CHKERRQ(ierr); 52*c4762a1bSJed Brown tsdae->Fctx = ctx; 53*c4762a1bSJed Brown PetscFunctionReturn(0); 54*c4762a1bSJed Brown } 55*c4762a1bSJed Brown 56*c4762a1bSJed Brown PetscErrorCode TSDAESimpleDestroy(TSDAESimple *tsdae) 57*c4762a1bSJed Brown { 58*c4762a1bSJed Brown PetscErrorCode ierr; 59*c4762a1bSJed Brown 60*c4762a1bSJed Brown PetscFunctionBegin; 61*c4762a1bSJed Brown ierr = (*(*tsdae)->destroy)(*tsdae);CHKERRQ(ierr); 62*c4762a1bSJed Brown ierr = VecDestroy(&(*tsdae)->U);CHKERRQ(ierr); 63*c4762a1bSJed Brown ierr = VecDestroy(&(*tsdae)->V);CHKERRQ(ierr); 64*c4762a1bSJed Brown ierr = PetscFree(*tsdae);CHKERRQ(ierr); 65*c4762a1bSJed Brown PetscFunctionReturn(0); 66*c4762a1bSJed Brown } 67*c4762a1bSJed Brown 68*c4762a1bSJed Brown PetscErrorCode TSDAESimpleSolve(TSDAESimple tsdae,Vec Usolution) 69*c4762a1bSJed Brown { 70*c4762a1bSJed Brown PetscErrorCode ierr; 71*c4762a1bSJed Brown 72*c4762a1bSJed Brown PetscFunctionBegin; 73*c4762a1bSJed Brown ierr = (*tsdae->solve)(tsdae,Usolution);CHKERRQ(ierr); 74*c4762a1bSJed Brown PetscFunctionReturn(0); 75*c4762a1bSJed Brown } 76*c4762a1bSJed Brown 77*c4762a1bSJed Brown PetscErrorCode TSDAESimpleSetFromOptions(TSDAESimple tsdae) 78*c4762a1bSJed Brown { 79*c4762a1bSJed Brown PetscErrorCode ierr; 80*c4762a1bSJed Brown 81*c4762a1bSJed Brown PetscFunctionBegin; 82*c4762a1bSJed Brown ierr = (*tsdae->setfromoptions)(PetscOptionsObject,tsdae);CHKERRQ(ierr); 83*c4762a1bSJed Brown PetscFunctionReturn(0); 84*c4762a1bSJed Brown } 85*c4762a1bSJed Brown 86*c4762a1bSJed Brown /* ----------------------------------------------------------------------------*/ 87*c4762a1bSJed Brown /* 88*c4762a1bSJed Brown Integrates the system by integrating the reduced ODE system and solving the 89*c4762a1bSJed Brown algebraic constraints at each stage with a separate SNES solve. 90*c4762a1bSJed Brown */ 91*c4762a1bSJed Brown 92*c4762a1bSJed Brown typedef struct { 93*c4762a1bSJed Brown PetscReal t; 94*c4762a1bSJed Brown TS ts; 95*c4762a1bSJed Brown SNES snes; 96*c4762a1bSJed Brown Vec U; 97*c4762a1bSJed Brown } TSDAESimple_Reduced; 98*c4762a1bSJed Brown 99*c4762a1bSJed Brown /* 100*c4762a1bSJed Brown Defines the RHS function that is passed to the time-integrator. 101*c4762a1bSJed Brown 102*c4762a1bSJed Brown Solves F(U,V) for V and then computes f(U,V) 103*c4762a1bSJed Brown 104*c4762a1bSJed Brown */ 105*c4762a1bSJed Brown PetscErrorCode TSDAESimple_Reduced_TSFunction(TS ts,PetscReal t,Vec U,Vec F,void *actx) 106*c4762a1bSJed Brown { 107*c4762a1bSJed Brown TSDAESimple tsdae = (TSDAESimple)actx; 108*c4762a1bSJed Brown TSDAESimple_Reduced *red = (TSDAESimple_Reduced*)tsdae->data; 109*c4762a1bSJed Brown PetscErrorCode ierr; 110*c4762a1bSJed Brown 111*c4762a1bSJed Brown PetscFunctionBeginUser; 112*c4762a1bSJed Brown red->t = t; 113*c4762a1bSJed Brown red->U = U; 114*c4762a1bSJed Brown ierr = SNESSolve(red->snes,NULL,tsdae->V);CHKERRQ(ierr); 115*c4762a1bSJed Brown ierr = (*tsdae->f)(t,U,tsdae->V,F,tsdae->fctx);CHKERRQ(ierr); 116*c4762a1bSJed Brown PetscFunctionReturn(0); 117*c4762a1bSJed Brown } 118*c4762a1bSJed Brown 119*c4762a1bSJed Brown /* 120*c4762a1bSJed Brown Defines the nonlinear function that is passed to the nonlinear solver 121*c4762a1bSJed Brown 122*c4762a1bSJed Brown */ 123*c4762a1bSJed Brown PetscErrorCode TSDAESimple_Reduced_SNESFunction(SNES snes,Vec V,Vec F,void *actx) 124*c4762a1bSJed Brown { 125*c4762a1bSJed Brown TSDAESimple tsdae = (TSDAESimple)actx; 126*c4762a1bSJed Brown TSDAESimple_Reduced *red = (TSDAESimple_Reduced*)tsdae->data; 127*c4762a1bSJed Brown PetscErrorCode ierr; 128*c4762a1bSJed Brown 129*c4762a1bSJed Brown PetscFunctionBeginUser; 130*c4762a1bSJed Brown ierr = (*tsdae->F)(red->t,red->U,V,F,tsdae->Fctx);CHKERRQ(ierr); 131*c4762a1bSJed Brown PetscFunctionReturn(0); 132*c4762a1bSJed Brown } 133*c4762a1bSJed Brown 134*c4762a1bSJed Brown 135*c4762a1bSJed Brown PetscErrorCode TSDAESimpleSolve_Reduced(TSDAESimple tsdae,Vec U) 136*c4762a1bSJed Brown { 137*c4762a1bSJed Brown PetscErrorCode ierr; 138*c4762a1bSJed Brown TSDAESimple_Reduced *red = (TSDAESimple_Reduced*)tsdae->data; 139*c4762a1bSJed Brown 140*c4762a1bSJed Brown PetscFunctionBegin; 141*c4762a1bSJed Brown ierr = TSSolve(red->ts,U);CHKERRQ(ierr); 142*c4762a1bSJed Brown PetscFunctionReturn(0); 143*c4762a1bSJed Brown } 144*c4762a1bSJed Brown 145*c4762a1bSJed Brown PetscErrorCode TSDAESimpleSetFromOptions_Reduced(PetscOptionItems *PetscOptionsObject,TSDAESimple tsdae) 146*c4762a1bSJed Brown { 147*c4762a1bSJed Brown PetscErrorCode ierr; 148*c4762a1bSJed Brown TSDAESimple_Reduced *red = (TSDAESimple_Reduced*)tsdae->data; 149*c4762a1bSJed Brown 150*c4762a1bSJed Brown PetscFunctionBegin; 151*c4762a1bSJed Brown ierr = TSSetFromOptions(red->ts);CHKERRQ(ierr); 152*c4762a1bSJed Brown ierr = SNESSetFromOptions(red->snes);CHKERRQ(ierr); 153*c4762a1bSJed Brown PetscFunctionReturn(0); 154*c4762a1bSJed Brown } 155*c4762a1bSJed Brown 156*c4762a1bSJed Brown PetscErrorCode TSDAESimpleDestroy_Reduced(TSDAESimple tsdae) 157*c4762a1bSJed Brown { 158*c4762a1bSJed Brown PetscErrorCode ierr; 159*c4762a1bSJed Brown TSDAESimple_Reduced *red = (TSDAESimple_Reduced*)tsdae->data; 160*c4762a1bSJed Brown 161*c4762a1bSJed Brown PetscFunctionBegin; 162*c4762a1bSJed Brown ierr = TSDestroy(&red->ts);CHKERRQ(ierr); 163*c4762a1bSJed Brown ierr = SNESDestroy(&red->snes);CHKERRQ(ierr); 164*c4762a1bSJed Brown ierr = PetscFree(red);CHKERRQ(ierr); 165*c4762a1bSJed Brown PetscFunctionReturn(0); 166*c4762a1bSJed Brown } 167*c4762a1bSJed Brown 168*c4762a1bSJed Brown PetscErrorCode TSDAESimpleSetUp_Reduced(TSDAESimple tsdae) 169*c4762a1bSJed Brown { 170*c4762a1bSJed Brown PetscErrorCode ierr; 171*c4762a1bSJed Brown TSDAESimple_Reduced *red; 172*c4762a1bSJed Brown Vec tsrhs; 173*c4762a1bSJed Brown 174*c4762a1bSJed Brown PetscFunctionBegin; 175*c4762a1bSJed Brown ierr = PetscNew(&red);CHKERRQ(ierr); 176*c4762a1bSJed Brown tsdae->data = red; 177*c4762a1bSJed Brown 178*c4762a1bSJed Brown tsdae->setfromoptions = TSDAESimpleSetFromOptions_Reduced; 179*c4762a1bSJed Brown tsdae->solve = TSDAESimpleSolve_Reduced; 180*c4762a1bSJed Brown tsdae->destroy = TSDAESimpleDestroy_Reduced; 181*c4762a1bSJed Brown 182*c4762a1bSJed Brown ierr = TSCreate(tsdae->comm,&red->ts);CHKERRQ(ierr); 183*c4762a1bSJed Brown ierr = TSSetProblemType(red->ts,TS_NONLINEAR);CHKERRQ(ierr); 184*c4762a1bSJed Brown ierr = TSSetType(red->ts,TSEULER);CHKERRQ(ierr); 185*c4762a1bSJed Brown ierr = TSSetExactFinalTime(red->ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); 186*c4762a1bSJed Brown ierr = VecDuplicate(tsdae->U,&tsrhs);CHKERRQ(ierr); 187*c4762a1bSJed Brown ierr = TSSetRHSFunction(red->ts,tsrhs,TSDAESimple_Reduced_TSFunction,tsdae);CHKERRQ(ierr); 188*c4762a1bSJed Brown ierr = VecDestroy(&tsrhs);CHKERRQ(ierr); 189*c4762a1bSJed Brown 190*c4762a1bSJed Brown ierr = SNESCreate(tsdae->comm,&red->snes);CHKERRQ(ierr); 191*c4762a1bSJed Brown ierr = SNESSetOptionsPrefix(red->snes,"tsdaesimple_");CHKERRQ(ierr); 192*c4762a1bSJed Brown ierr = SNESSetFunction(red->snes,NULL,TSDAESimple_Reduced_SNESFunction,tsdae);CHKERRQ(ierr); 193*c4762a1bSJed Brown ierr = SNESSetJacobian(red->snes,NULL,NULL,SNESComputeJacobianDefault,tsdae);CHKERRQ(ierr); 194*c4762a1bSJed Brown PetscFunctionReturn(0); 195*c4762a1bSJed Brown } 196*c4762a1bSJed Brown 197*c4762a1bSJed Brown 198*c4762a1bSJed Brown /* ----------------------------------------------------------------------------*/ 199*c4762a1bSJed Brown 200*c4762a1bSJed Brown /* 201*c4762a1bSJed Brown Integrates the system by integrating directly the entire DAE system 202*c4762a1bSJed Brown */ 203*c4762a1bSJed Brown 204*c4762a1bSJed Brown typedef struct { 205*c4762a1bSJed Brown TS ts; 206*c4762a1bSJed Brown Vec UV,UF,VF; 207*c4762a1bSJed Brown VecScatter scatterU,scatterV; 208*c4762a1bSJed Brown } TSDAESimple_Full; 209*c4762a1bSJed Brown 210*c4762a1bSJed Brown /* 211*c4762a1bSJed Brown Defines the RHS function that is passed to the time-integrator. 212*c4762a1bSJed Brown 213*c4762a1bSJed Brown f(U,V) 214*c4762a1bSJed Brown 0 215*c4762a1bSJed Brown 216*c4762a1bSJed Brown */ 217*c4762a1bSJed Brown PetscErrorCode TSDAESimple_Full_TSRHSFunction(TS ts,PetscReal t,Vec UV,Vec F,void *actx) 218*c4762a1bSJed Brown { 219*c4762a1bSJed Brown TSDAESimple tsdae = (TSDAESimple)actx; 220*c4762a1bSJed Brown TSDAESimple_Full *full = (TSDAESimple_Full*)tsdae->data; 221*c4762a1bSJed Brown PetscErrorCode ierr; 222*c4762a1bSJed Brown 223*c4762a1bSJed Brown PetscFunctionBeginUser; 224*c4762a1bSJed Brown ierr = VecSet(F,0.0);CHKERRQ(ierr); 225*c4762a1bSJed Brown ierr = VecScatterBegin(full->scatterU,UV,tsdae->U,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 226*c4762a1bSJed Brown ierr = VecScatterEnd(full->scatterU,UV,tsdae->U,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 227*c4762a1bSJed Brown ierr = VecScatterBegin(full->scatterV,UV,tsdae->V,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 228*c4762a1bSJed Brown ierr = VecScatterEnd(full->scatterV,UV,tsdae->V,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 229*c4762a1bSJed Brown ierr = (*tsdae->f)(t,tsdae->U,tsdae->V,full->UF,tsdae->fctx);CHKERRQ(ierr); 230*c4762a1bSJed Brown ierr = VecScatterBegin(full->scatterU,full->UF,F,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 231*c4762a1bSJed Brown ierr = VecScatterEnd(full->scatterU,full->UF,F,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 232*c4762a1bSJed Brown PetscFunctionReturn(0); 233*c4762a1bSJed Brown } 234*c4762a1bSJed Brown 235*c4762a1bSJed Brown /* 236*c4762a1bSJed Brown Defines the nonlinear function that is passed to the nonlinear solver 237*c4762a1bSJed Brown 238*c4762a1bSJed Brown \dot{U} 239*c4762a1bSJed Brown F(U,V) 240*c4762a1bSJed Brown 241*c4762a1bSJed Brown */ 242*c4762a1bSJed Brown PetscErrorCode TSDAESimple_Full_TSIFunction(TS ts,PetscReal t,Vec UV,Vec UVdot,Vec F,void *actx) 243*c4762a1bSJed Brown { 244*c4762a1bSJed Brown TSDAESimple tsdae = (TSDAESimple)actx; 245*c4762a1bSJed Brown TSDAESimple_Full *full = (TSDAESimple_Full*)tsdae->data; 246*c4762a1bSJed Brown PetscErrorCode ierr; 247*c4762a1bSJed Brown 248*c4762a1bSJed Brown PetscFunctionBeginUser; 249*c4762a1bSJed Brown ierr = VecCopy(UVdot,F);CHKERRQ(ierr); 250*c4762a1bSJed Brown ierr = VecScatterBegin(full->scatterU,UV,tsdae->U,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 251*c4762a1bSJed Brown ierr = VecScatterEnd(full->scatterU,UV,tsdae->U,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 252*c4762a1bSJed Brown ierr = VecScatterBegin(full->scatterV,UV,tsdae->V,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 253*c4762a1bSJed Brown ierr = VecScatterEnd(full->scatterV,UV,tsdae->V,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 254*c4762a1bSJed Brown ierr = (*tsdae->F)(t,tsdae->U,tsdae->V,full->VF,tsdae->Fctx);CHKERRQ(ierr); 255*c4762a1bSJed Brown ierr = VecScatterBegin(full->scatterV,full->VF,F,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 256*c4762a1bSJed Brown ierr = VecScatterEnd(full->scatterV,full->VF,F,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 257*c4762a1bSJed Brown PetscFunctionReturn(0); 258*c4762a1bSJed Brown } 259*c4762a1bSJed Brown 260*c4762a1bSJed Brown 261*c4762a1bSJed Brown PetscErrorCode TSDAESimpleSolve_Full(TSDAESimple tsdae,Vec U) 262*c4762a1bSJed Brown { 263*c4762a1bSJed Brown PetscErrorCode ierr; 264*c4762a1bSJed Brown TSDAESimple_Full *full = (TSDAESimple_Full*)tsdae->data; 265*c4762a1bSJed Brown 266*c4762a1bSJed Brown PetscFunctionBegin; 267*c4762a1bSJed Brown ierr = VecSet(full->UV,1.0);CHKERRQ(ierr); 268*c4762a1bSJed Brown ierr = VecScatterBegin(full->scatterU,U,full->UV,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 269*c4762a1bSJed Brown ierr = VecScatterEnd(full->scatterU,U,full->UV,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 270*c4762a1bSJed Brown ierr = TSSolve(full->ts,full->UV);CHKERRQ(ierr); 271*c4762a1bSJed Brown ierr = VecScatterBegin(full->scatterU,full->UV,U,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 272*c4762a1bSJed Brown ierr = VecScatterEnd(full->scatterU,full->UV,U,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 273*c4762a1bSJed Brown PetscFunctionReturn(0); 274*c4762a1bSJed Brown } 275*c4762a1bSJed Brown 276*c4762a1bSJed Brown PetscErrorCode TSDAESimpleSetFromOptions_Full(PetscOptionItems *PetscOptionsObject,TSDAESimple tsdae) 277*c4762a1bSJed Brown { 278*c4762a1bSJed Brown PetscErrorCode ierr; 279*c4762a1bSJed Brown TSDAESimple_Full *full = (TSDAESimple_Full*)tsdae->data; 280*c4762a1bSJed Brown 281*c4762a1bSJed Brown PetscFunctionBegin; 282*c4762a1bSJed Brown ierr = TSSetFromOptions(full->ts);CHKERRQ(ierr); 283*c4762a1bSJed Brown PetscFunctionReturn(0); 284*c4762a1bSJed Brown } 285*c4762a1bSJed Brown 286*c4762a1bSJed Brown PetscErrorCode TSDAESimpleDestroy_Full(TSDAESimple tsdae) 287*c4762a1bSJed Brown { 288*c4762a1bSJed Brown PetscErrorCode ierr; 289*c4762a1bSJed Brown TSDAESimple_Full *full = (TSDAESimple_Full*)tsdae->data; 290*c4762a1bSJed Brown 291*c4762a1bSJed Brown PetscFunctionBegin; 292*c4762a1bSJed Brown ierr = TSDestroy(&full->ts);CHKERRQ(ierr); 293*c4762a1bSJed Brown ierr = VecDestroy(&full->UV);CHKERRQ(ierr); 294*c4762a1bSJed Brown ierr = VecDestroy(&full->UF);CHKERRQ(ierr); 295*c4762a1bSJed Brown ierr = VecDestroy(&full->VF);CHKERRQ(ierr); 296*c4762a1bSJed Brown ierr = VecScatterDestroy(&full->scatterU);CHKERRQ(ierr); 297*c4762a1bSJed Brown ierr = VecScatterDestroy(&full->scatterV);CHKERRQ(ierr); 298*c4762a1bSJed Brown ierr = PetscFree(full);CHKERRQ(ierr); 299*c4762a1bSJed Brown PetscFunctionReturn(0); 300*c4762a1bSJed Brown } 301*c4762a1bSJed Brown 302*c4762a1bSJed Brown PetscErrorCode TSDAESimpleSetUp_Full(TSDAESimple tsdae) 303*c4762a1bSJed Brown { 304*c4762a1bSJed Brown PetscErrorCode ierr; 305*c4762a1bSJed Brown TSDAESimple_Full *full; 306*c4762a1bSJed Brown Vec tsrhs; 307*c4762a1bSJed Brown PetscInt nU,nV,UVstart; 308*c4762a1bSJed Brown IS is; 309*c4762a1bSJed Brown 310*c4762a1bSJed Brown PetscFunctionBegin; 311*c4762a1bSJed Brown ierr = PetscNew(&full);CHKERRQ(ierr); 312*c4762a1bSJed Brown tsdae->data = full; 313*c4762a1bSJed Brown 314*c4762a1bSJed Brown tsdae->setfromoptions = TSDAESimpleSetFromOptions_Full; 315*c4762a1bSJed Brown tsdae->solve = TSDAESimpleSolve_Full; 316*c4762a1bSJed Brown tsdae->destroy = TSDAESimpleDestroy_Full; 317*c4762a1bSJed Brown 318*c4762a1bSJed Brown ierr = TSCreate(tsdae->comm,&full->ts);CHKERRQ(ierr); 319*c4762a1bSJed Brown ierr = TSSetProblemType(full->ts,TS_NONLINEAR);CHKERRQ(ierr); 320*c4762a1bSJed Brown ierr = TSSetType(full->ts,TSROSW);CHKERRQ(ierr); 321*c4762a1bSJed Brown ierr = TSSetExactFinalTime(full->ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); 322*c4762a1bSJed Brown ierr = VecDuplicate(tsdae->U,&full->UF);CHKERRQ(ierr); 323*c4762a1bSJed Brown ierr = VecDuplicate(tsdae->V,&full->VF);CHKERRQ(ierr); 324*c4762a1bSJed Brown 325*c4762a1bSJed Brown ierr = VecGetLocalSize(tsdae->U,&nU);CHKERRQ(ierr); 326*c4762a1bSJed Brown ierr = VecGetLocalSize(tsdae->V,&nV);CHKERRQ(ierr); 327*c4762a1bSJed Brown ierr = VecCreateMPI(tsdae->comm,nU+nV,PETSC_DETERMINE,&tsrhs);CHKERRQ(ierr); 328*c4762a1bSJed Brown ierr = VecDuplicate(tsrhs,&full->UV);CHKERRQ(ierr); 329*c4762a1bSJed Brown 330*c4762a1bSJed Brown ierr = VecGetOwnershipRange(tsrhs,&UVstart,NULL);CHKERRQ(ierr); 331*c4762a1bSJed Brown ierr = ISCreateStride(tsdae->comm,nU,UVstart,1,&is);CHKERRQ(ierr); 332*c4762a1bSJed Brown ierr = VecScatterCreate(tsdae->U,NULL,tsrhs,is,&full->scatterU);CHKERRQ(ierr); 333*c4762a1bSJed Brown ierr = ISDestroy(&is);CHKERRQ(ierr); 334*c4762a1bSJed Brown ierr = ISCreateStride(tsdae->comm,nV,UVstart+nU,1,&is);CHKERRQ(ierr); 335*c4762a1bSJed Brown ierr = VecScatterCreate(tsdae->V,NULL,tsrhs,is,&full->scatterV);CHKERRQ(ierr); 336*c4762a1bSJed Brown ierr = ISDestroy(&is);CHKERRQ(ierr); 337*c4762a1bSJed Brown 338*c4762a1bSJed Brown ierr = TSSetRHSFunction(full->ts,tsrhs,TSDAESimple_Full_TSRHSFunction,tsdae);CHKERRQ(ierr); 339*c4762a1bSJed Brown ierr = TSSetIFunction(full->ts,NULL,TSDAESimple_Full_TSIFunction,tsdae);CHKERRQ(ierr); 340*c4762a1bSJed Brown ierr = VecDestroy(&tsrhs);CHKERRQ(ierr); 341*c4762a1bSJed Brown PetscFunctionReturn(0); 342*c4762a1bSJed Brown } 343*c4762a1bSJed Brown 344*c4762a1bSJed Brown 345*c4762a1bSJed Brown /* ----------------------------------------------------------------------------*/ 346*c4762a1bSJed Brown 347*c4762a1bSJed Brown 348*c4762a1bSJed Brown /* 349*c4762a1bSJed Brown Simple example: f(U,V) = U + V 350*c4762a1bSJed Brown 351*c4762a1bSJed Brown */ 352*c4762a1bSJed Brown PetscErrorCode f(PetscReal t,Vec U,Vec V,Vec F,void *ctx) 353*c4762a1bSJed Brown { 354*c4762a1bSJed Brown PetscErrorCode ierr; 355*c4762a1bSJed Brown 356*c4762a1bSJed Brown PetscFunctionBeginUser; 357*c4762a1bSJed Brown ierr = VecWAXPY(F,1.0,U,V);CHKERRQ(ierr); 358*c4762a1bSJed Brown PetscFunctionReturn(0); 359*c4762a1bSJed Brown } 360*c4762a1bSJed Brown 361*c4762a1bSJed Brown /* 362*c4762a1bSJed Brown Simple example: F(U,V) = U - V 363*c4762a1bSJed Brown 364*c4762a1bSJed Brown */ 365*c4762a1bSJed Brown PetscErrorCode F(PetscReal t,Vec U,Vec V,Vec F,void *ctx) 366*c4762a1bSJed Brown { 367*c4762a1bSJed Brown PetscErrorCode ierr; 368*c4762a1bSJed Brown 369*c4762a1bSJed Brown PetscFunctionBeginUser; 370*c4762a1bSJed Brown ierr = VecWAXPY(F,-1.0,V,U);CHKERRQ(ierr); 371*c4762a1bSJed Brown PetscFunctionReturn(0); 372*c4762a1bSJed Brown } 373*c4762a1bSJed Brown 374*c4762a1bSJed Brown int main(int argc,char **argv) 375*c4762a1bSJed Brown { 376*c4762a1bSJed Brown PetscErrorCode ierr; 377*c4762a1bSJed Brown TSDAESimple tsdae; 378*c4762a1bSJed Brown Vec U,V,Usolution; 379*c4762a1bSJed Brown 380*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 381*c4762a1bSJed Brown ierr = TSDAESimpleCreate(PETSC_COMM_WORLD,&tsdae);CHKERRQ(ierr); 382*c4762a1bSJed Brown 383*c4762a1bSJed Brown ierr = VecCreateMPI(PETSC_COMM_WORLD,1,PETSC_DETERMINE,&U);CHKERRQ(ierr); 384*c4762a1bSJed Brown ierr = VecCreateMPI(PETSC_COMM_WORLD,1,PETSC_DETERMINE,&V);CHKERRQ(ierr); 385*c4762a1bSJed Brown ierr = TSDAESimpleSetRHSFunction(tsdae,U,f,NULL);CHKERRQ(ierr); 386*c4762a1bSJed Brown ierr = TSDAESimpleSetIFunction(tsdae,V,F,NULL);CHKERRQ(ierr); 387*c4762a1bSJed Brown 388*c4762a1bSJed Brown ierr = VecDuplicate(U,&Usolution);CHKERRQ(ierr); 389*c4762a1bSJed Brown ierr = VecSet(Usolution,1.0);CHKERRQ(ierr); 390*c4762a1bSJed Brown 391*c4762a1bSJed Brown /* ierr = TSDAESimpleSetUp_Full(tsdae);CHKERRQ(ierr); */ 392*c4762a1bSJed Brown ierr = TSDAESimpleSetUp_Reduced(tsdae);CHKERRQ(ierr); 393*c4762a1bSJed Brown 394*c4762a1bSJed Brown ierr = TSDAESimpleSetFromOptions(tsdae);CHKERRQ(ierr); 395*c4762a1bSJed Brown ierr = TSDAESimpleSolve(tsdae,Usolution);CHKERRQ(ierr); 396*c4762a1bSJed Brown ierr = TSDAESimpleDestroy(&tsdae);CHKERRQ(ierr); 397*c4762a1bSJed Brown 398*c4762a1bSJed Brown ierr = VecDestroy(&U);CHKERRQ(ierr); 399*c4762a1bSJed Brown ierr = VecDestroy(&Usolution);CHKERRQ(ierr); 400*c4762a1bSJed Brown ierr = VecDestroy(&V);CHKERRQ(ierr); 401*c4762a1bSJed Brown ierr = PetscFinalize(); 402*c4762a1bSJed Brown return ierr; 403*c4762a1bSJed Brown } 404*c4762a1bSJed Brown 405*c4762a1bSJed Brown 406*c4762a1bSJed Brown 407