1*8a381b04SJed Brown /* 2*8a381b04SJed Brown Code for timestepping with additive Runge-Kutta IMEX method 3*8a381b04SJed Brown 4*8a381b04SJed Brown Notes: 5*8a381b04SJed Brown The general system is written as 6*8a381b04SJed Brown 7*8a381b04SJed Brown F(t,X,Xdot) = G(t,X) 8*8a381b04SJed Brown 9*8a381b04SJed Brown where F represents the stiff part of the physics and G represents the non-stiff part. 10*8a381b04SJed Brown 11*8a381b04SJed Brown */ 12*8a381b04SJed Brown #include <private/tsimpl.h> /*I "petscts.h" I*/ 13*8a381b04SJed Brown 14*8a381b04SJed Brown static const TSARKIMEXType TSARKIMEXDefault = TSARKIMEX2D; 15*8a381b04SJed Brown static PetscBool TSARKIMEXRegisterAllCalled; 16*8a381b04SJed Brown static PetscBool TSARKIMEXPackageInitialized; 17*8a381b04SJed Brown 18*8a381b04SJed Brown typedef struct _ARKTableau *ARKTableau; 19*8a381b04SJed Brown struct _ARKTableau { 20*8a381b04SJed Brown char *name; 21*8a381b04SJed Brown PetscInt order; 22*8a381b04SJed Brown PetscInt s; 23*8a381b04SJed Brown PetscReal *At,*bt,*ct; 24*8a381b04SJed Brown PetscReal *A,*b,*c; /* Non-stiff tableau */ 25*8a381b04SJed Brown }; 26*8a381b04SJed Brown typedef struct _ARKTableauLink *ARKTableauLink; 27*8a381b04SJed Brown struct _ARKTableauLink { 28*8a381b04SJed Brown struct _ARKTableau tab; 29*8a381b04SJed Brown ARKTableauLink next; 30*8a381b04SJed Brown }; 31*8a381b04SJed Brown static ARKTableauLink ARKTableauList; 32*8a381b04SJed Brown 33*8a381b04SJed Brown typedef struct { 34*8a381b04SJed Brown ARKTableau tableau; 35*8a381b04SJed Brown Vec *Y; /* States computed during the step */ 36*8a381b04SJed Brown Vec *YdotI; /* Time derivatives for the stiff part */ 37*8a381b04SJed Brown Vec *YdotRHS; /* Function evaluations for the non-stiff part */ 38*8a381b04SJed Brown Vec Ydot; /* Work vector holding Ydot during residual evaluation */ 39*8a381b04SJed Brown Vec Work; /* Generic work vector */ 40*8a381b04SJed Brown Vec Z; /* Ydot = shift(Y-Z) */ 41*8a381b04SJed Brown PetscScalar *work; /* Scalar work */ 42*8a381b04SJed Brown PetscReal shift; 43*8a381b04SJed Brown PetscReal stage_time; 44*8a381b04SJed Brown } TS_ARKIMEX; 45*8a381b04SJed Brown 46*8a381b04SJed Brown #undef __FUNCT__ 47*8a381b04SJed Brown #define __FUNCT__ "TSARKIMEXRegisterAll" 48*8a381b04SJed Brown /*@C 49*8a381b04SJed Brown TSARKIMEXRegisterAll - Registers all of the additive Runge-Kutta implicit-explicit methods in TSARKIMEX 50*8a381b04SJed Brown 51*8a381b04SJed Brown Not Collective 52*8a381b04SJed Brown 53*8a381b04SJed Brown Level: advanced 54*8a381b04SJed Brown 55*8a381b04SJed Brown .keywords: TS, TSARKIMEX, register, all 56*8a381b04SJed Brown 57*8a381b04SJed Brown .seealso: TSARKIMEXRegisterDestroy() 58*8a381b04SJed Brown @*/ 59*8a381b04SJed Brown PetscErrorCode TSARKIMEXRegisterAll(void) 60*8a381b04SJed Brown { 61*8a381b04SJed Brown PetscErrorCode ierr; 62*8a381b04SJed Brown 63*8a381b04SJed Brown PetscFunctionBegin; 64*8a381b04SJed Brown if (TSARKIMEXRegisterAllCalled) PetscFunctionReturn(0); 65*8a381b04SJed Brown TSARKIMEXRegisterAllCalled = PETSC_TRUE; 66*8a381b04SJed Brown 67*8a381b04SJed Brown { 68*8a381b04SJed Brown const PetscReal 69*8a381b04SJed Brown A[3][3] = {{0,0,0}, 70*8a381b04SJed Brown {0.41421356237309504880,0,0}, 71*8a381b04SJed Brown {0.75,0.25,0}}, 72*8a381b04SJed Brown At[3][3] = {{0,0,0}, 73*8a381b04SJed Brown {0.12132034355964257320,0.29289321881345247560,0}, 74*8a381b04SJed Brown {0.20710678118654752440,0.50000000000000000000,0.29289321881345247560}}; 75*8a381b04SJed Brown ierr = TSARKIMEXRegister(TSARKIMEX2D,2,3,&At[0][0],PETSC_NULL,PETSC_NULL,&A[0][0],PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 76*8a381b04SJed Brown } 77*8a381b04SJed Brown PetscFunctionReturn(0); 78*8a381b04SJed Brown } 79*8a381b04SJed Brown 80*8a381b04SJed Brown #undef __FUNCT__ 81*8a381b04SJed Brown #define __FUNCT__ "TSARKIMEXRegisterDestroy" 82*8a381b04SJed Brown /*@C 83*8a381b04SJed Brown TSARKIMEXRegisterDestroy - Frees the list of schemes that were registered by TSARKIMEXRegister(). 84*8a381b04SJed Brown 85*8a381b04SJed Brown Not Collective 86*8a381b04SJed Brown 87*8a381b04SJed Brown Level: advanced 88*8a381b04SJed Brown 89*8a381b04SJed Brown .keywords: TSARKIMEX, register, destroy 90*8a381b04SJed Brown .seealso: TSARKIMEXRegister(), TSARKIMEXRegisterAll(), TSARKIMEXRegisterDynamic() 91*8a381b04SJed Brown @*/ 92*8a381b04SJed Brown PetscErrorCode TSARKIMEXRegisterDestroy(void) 93*8a381b04SJed Brown { 94*8a381b04SJed Brown PetscErrorCode ierr; 95*8a381b04SJed Brown ARKTableauLink link; 96*8a381b04SJed Brown 97*8a381b04SJed Brown PetscFunctionBegin; 98*8a381b04SJed Brown while ((link = ARKTableauList)) { 99*8a381b04SJed Brown ARKTableau t = &link->tab; 100*8a381b04SJed Brown ARKTableauList = link->next; 101*8a381b04SJed Brown ierr = PetscFree6(t->At,t->bt,t->ct,t->A,t->b,t->c);CHKERRQ(ierr); 102*8a381b04SJed Brown ierr = PetscFree(t->name);CHKERRQ(ierr); 103*8a381b04SJed Brown ierr = PetscFree(link);CHKERRQ(ierr); 104*8a381b04SJed Brown } 105*8a381b04SJed Brown TSARKIMEXRegisterAllCalled = PETSC_FALSE; 106*8a381b04SJed Brown PetscFunctionReturn(0); 107*8a381b04SJed Brown } 108*8a381b04SJed Brown 109*8a381b04SJed Brown #undef __FUNCT__ 110*8a381b04SJed Brown #define __FUNCT__ "TSARKIMEXInitializePackage" 111*8a381b04SJed Brown /*@C 112*8a381b04SJed Brown TSARKIMEXInitializePackage - This function initializes everything in the TSARKIMEX package. It is called 113*8a381b04SJed Brown from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_ARKIMEX() 114*8a381b04SJed Brown when using static libraries. 115*8a381b04SJed Brown 116*8a381b04SJed Brown Input Parameter: 117*8a381b04SJed Brown path - The dynamic library path, or PETSC_NULL 118*8a381b04SJed Brown 119*8a381b04SJed Brown Level: developer 120*8a381b04SJed Brown 121*8a381b04SJed Brown .keywords: TS, TSARKIMEX, initialize, package 122*8a381b04SJed Brown .seealso: PetscInitialize() 123*8a381b04SJed Brown @*/ 124*8a381b04SJed Brown PetscErrorCode TSARKIMEXInitializePackage(const char path[]) 125*8a381b04SJed Brown { 126*8a381b04SJed Brown PetscErrorCode ierr; 127*8a381b04SJed Brown 128*8a381b04SJed Brown PetscFunctionBegin; 129*8a381b04SJed Brown if (TSARKIMEXPackageInitialized) PetscFunctionReturn(0); 130*8a381b04SJed Brown TSARKIMEXPackageInitialized = PETSC_TRUE; 131*8a381b04SJed Brown ierr = TSARKIMEXRegisterAll();CHKERRQ(ierr); 132*8a381b04SJed Brown ierr = PetscRegisterFinalize(TSARKIMEXFinalizePackage);CHKERRQ(ierr); 133*8a381b04SJed Brown PetscFunctionReturn(0); 134*8a381b04SJed Brown } 135*8a381b04SJed Brown 136*8a381b04SJed Brown #undef __FUNCT__ 137*8a381b04SJed Brown #define __FUNCT__ "TSARKIMEXFinalizePackage" 138*8a381b04SJed Brown /*@C 139*8a381b04SJed Brown TSARKIMEXFinalizePackage - This function destroys everything in the TSARKIMEX package. It is 140*8a381b04SJed Brown called from PetscFinalize(). 141*8a381b04SJed Brown 142*8a381b04SJed Brown Level: developer 143*8a381b04SJed Brown 144*8a381b04SJed Brown .keywords: Petsc, destroy, package 145*8a381b04SJed Brown .seealso: PetscFinalize() 146*8a381b04SJed Brown @*/ 147*8a381b04SJed Brown PetscErrorCode TSARKIMEXFinalizePackage(void) 148*8a381b04SJed Brown { 149*8a381b04SJed Brown PetscErrorCode ierr; 150*8a381b04SJed Brown 151*8a381b04SJed Brown PetscFunctionBegin; 152*8a381b04SJed Brown TSARKIMEXPackageInitialized = PETSC_FALSE; 153*8a381b04SJed Brown ierr = TSARKIMEXRegisterDestroy();CHKERRQ(ierr); 154*8a381b04SJed Brown PetscFunctionReturn(0); 155*8a381b04SJed Brown } 156*8a381b04SJed Brown 157*8a381b04SJed Brown #undef __FUNCT__ 158*8a381b04SJed Brown #define __FUNCT__ "TSARKIMEXRegister" 159*8a381b04SJed Brown PetscErrorCode TSARKIMEXRegister(const TSARKIMEXType name,PetscInt order,PetscInt s, 160*8a381b04SJed Brown const PetscReal At[],const PetscReal bt[],const PetscReal ct[], 161*8a381b04SJed Brown const PetscReal A[],const PetscReal b[],const PetscReal c[]) 162*8a381b04SJed Brown { 163*8a381b04SJed Brown PetscErrorCode ierr; 164*8a381b04SJed Brown ARKTableauLink link; 165*8a381b04SJed Brown ARKTableau t; 166*8a381b04SJed Brown PetscInt i,j; 167*8a381b04SJed Brown 168*8a381b04SJed Brown PetscFunctionBegin; 169*8a381b04SJed Brown ierr = PetscMalloc(sizeof(*link),&link);CHKERRQ(ierr); 170*8a381b04SJed Brown t = &link->tab; 171*8a381b04SJed Brown ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr); 172*8a381b04SJed Brown t->order = order; 173*8a381b04SJed Brown t->s = s; 174*8a381b04SJed Brown ierr = PetscMalloc6(s*s,PetscReal,&t->At,s,PetscReal,&t->bt,s,PetscReal,&t->ct,s*s,PetscReal,&t->A,s,PetscReal,&t->b,s,PetscReal,&t->c);CHKERRQ(ierr); 175*8a381b04SJed Brown ierr = PetscMemcpy(t->At,At,s*s*sizeof(At[0]));CHKERRQ(ierr); 176*8a381b04SJed Brown ierr = PetscMemcpy(t->A,A,s*s*sizeof(A[0]));CHKERRQ(ierr); 177*8a381b04SJed Brown if (bt) {ierr = PetscMemcpy(t->bt,bt,s*sizeof(bt[0]));CHKERRQ(ierr);} 178*8a381b04SJed Brown else for (i=0; i<s; i++) t->bt[i] = At[(s-1)*s+i]; 179*8a381b04SJed Brown if (b) {ierr = PetscMemcpy(t->b,b,s*sizeof(b[0]));CHKERRQ(ierr);} 180*8a381b04SJed Brown else for (i=0; i<s; i++) t->b[i] = At[(s-1)*s+i]; 181*8a381b04SJed Brown if (ct) {ierr = PetscMemcpy(t->ct,ct,s*sizeof(ct[0]));CHKERRQ(ierr);} 182*8a381b04SJed Brown else for (i=0; i<s; i++) for (j=0,t->ct[i]=0; j<s; j++) t->ct[i] += At[i*s+j]; 183*8a381b04SJed Brown if (c) {ierr = PetscMemcpy(t->c,c,s*sizeof(c[0]));CHKERRQ(ierr);} 184*8a381b04SJed Brown else for (i=0; i<s; i++) for (j=0,t->c[i]=0; j<s; j++) t->c[i] += A[i*s+j]; 185*8a381b04SJed Brown link->next = ARKTableauList; 186*8a381b04SJed Brown ARKTableauList = link; 187*8a381b04SJed Brown PetscFunctionReturn(0); 188*8a381b04SJed Brown } 189*8a381b04SJed Brown 190*8a381b04SJed Brown #undef __FUNCT__ 191*8a381b04SJed Brown #define __FUNCT__ "TSStep_ARKIMEX" 192*8a381b04SJed Brown static PetscErrorCode TSStep_ARKIMEX(TS ts) 193*8a381b04SJed Brown { 194*8a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 195*8a381b04SJed Brown ARKTableau tab = ark->tableau; 196*8a381b04SJed Brown const PetscInt s = tab->s; 197*8a381b04SJed Brown const PetscReal *At = tab->At,*A = tab->A,*bt = tab->bt,*b = tab->b,*ct = tab->ct,*c = tab->c; 198*8a381b04SJed Brown PetscReal *w = ark->work; 199*8a381b04SJed Brown Vec *Y = ark->Y,*YdotI = ark->YdotI,*YdotRHS = ark->YdotRHS,Ydot = ark->Ydot,W = ark->Work,Z = ark->Z; 200*8a381b04SJed Brown SNES snes; 201*8a381b04SJed Brown PetscInt i,j,its,lits; 202*8a381b04SJed Brown PetscReal h,t; 203*8a381b04SJed Brown PetscErrorCode ierr; 204*8a381b04SJed Brown 205*8a381b04SJed Brown PetscFunctionBegin; 206*8a381b04SJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 207*8a381b04SJed Brown h = ts->time_step = ts->next_time_step; 208*8a381b04SJed Brown t = ts->ptime; 209*8a381b04SJed Brown 210*8a381b04SJed Brown for (i=0; i<s; i++) { 211*8a381b04SJed Brown if (At[i*s+i] == 0) { /* This stage is explicit */ 212*8a381b04SJed Brown ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr); 213*8a381b04SJed Brown for (j=0; j<i; j++) w[j] = -h*At[i*s+j]; 214*8a381b04SJed Brown ierr = VecMAXPY(Y[i],i,w,YdotI);CHKERRQ(ierr); 215*8a381b04SJed Brown for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 216*8a381b04SJed Brown ierr = VecMAXPY(Y[i],i,w,YdotRHS);CHKERRQ(ierr); 217*8a381b04SJed Brown } else { 218*8a381b04SJed Brown ark->stage_time = t + h*ct[i]; 219*8a381b04SJed Brown ark->shift = 1./(h*At[i*s+i]); 220*8a381b04SJed Brown /* Affine part */ 221*8a381b04SJed Brown ierr = VecZeroEntries(W);CHKERRQ(ierr); 222*8a381b04SJed Brown for (j=0; j<i; j++) w[j] = h*A[i*s+j]; 223*8a381b04SJed Brown ierr = VecMAXPY(W,i,w,YdotRHS);CHKERRQ(ierr); 224*8a381b04SJed Brown /* Ydot = shift*(Y-Z) */ 225*8a381b04SJed Brown ierr = VecCopy(ts->vec_sol,Z);CHKERRQ(ierr); 226*8a381b04SJed Brown for (j=0; j<i; j++) w[j] = h*At[i*s+j]; 227*8a381b04SJed Brown ierr = VecMAXPY(Z,i,w,YdotRHS);CHKERRQ(ierr); 228*8a381b04SJed Brown /* Initial guess taken from last stage */ 229*8a381b04SJed Brown ierr = VecCopy(i>0?Y[i-1]:ts->vec_sol,Y[i]);CHKERRQ(ierr); 230*8a381b04SJed Brown ierr = SNESSolve(snes,W,Y[i]);CHKERRQ(ierr); 231*8a381b04SJed Brown ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr); 232*8a381b04SJed Brown ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr); 233*8a381b04SJed Brown ts->nonlinear_its += its; ts->linear_its += lits; 234*8a381b04SJed Brown } 235*8a381b04SJed Brown ierr = VecZeroEntries(Ydot);CHKERRQ(ierr); 236*8a381b04SJed Brown ierr = TSComputeIFunction(ts,t+h*ct[i],Y[i],Ydot,YdotI[i],PETSC_TRUE);CHKERRQ(ierr); 237*8a381b04SJed Brown ierr = TSComputeRHSFunction(ts,t+h*c[i],Y[i],YdotRHS[i]);CHKERRQ(ierr); 238*8a381b04SJed Brown } 239*8a381b04SJed Brown for (j=0; j<s; j++) w[j] = h*bt[j]; 240*8a381b04SJed Brown ierr = VecMAXPY(ts->vec_sol,s,w,YdotI);CHKERRQ(ierr); 241*8a381b04SJed Brown for (j=0; j<s; j++) w[j] = h*b[j]; 242*8a381b04SJed Brown ierr = VecMAXPY(ts->vec_sol,s,w,YdotRHS);CHKERRQ(ierr); 243*8a381b04SJed Brown 244*8a381b04SJed Brown ts->ptime += ts->time_step; 245*8a381b04SJed Brown ts->next_time_step = ts->time_step; 246*8a381b04SJed Brown ts->steps++; 247*8a381b04SJed Brown PetscFunctionReturn(0); 248*8a381b04SJed Brown } 249*8a381b04SJed Brown 250*8a381b04SJed Brown /*------------------------------------------------------------*/ 251*8a381b04SJed Brown #undef __FUNCT__ 252*8a381b04SJed Brown #define __FUNCT__ "TSReset_ARKIMEX" 253*8a381b04SJed Brown static PetscErrorCode TSReset_ARKIMEX(TS ts) 254*8a381b04SJed Brown { 255*8a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 256*8a381b04SJed Brown PetscInt s; 257*8a381b04SJed Brown PetscErrorCode ierr; 258*8a381b04SJed Brown 259*8a381b04SJed Brown PetscFunctionBegin; 260*8a381b04SJed Brown if (!ark->tableau) PetscFunctionReturn(0); 261*8a381b04SJed Brown s = ark->tableau->s; 262*8a381b04SJed Brown ierr = VecDestroyVecs(s,&ark->Y);CHKERRQ(ierr); 263*8a381b04SJed Brown ierr = VecDestroyVecs(s,&ark->YdotI);CHKERRQ(ierr); 264*8a381b04SJed Brown ierr = VecDestroyVecs(s,&ark->YdotRHS);CHKERRQ(ierr); 265*8a381b04SJed Brown ierr = VecDestroy(&ark->Ydot);CHKERRQ(ierr); 266*8a381b04SJed Brown ierr = VecDestroy(&ark->Work);CHKERRQ(ierr); 267*8a381b04SJed Brown ierr = VecDestroy(&ark->Z);CHKERRQ(ierr); 268*8a381b04SJed Brown ierr = PetscFree(ark->work);CHKERRQ(ierr); 269*8a381b04SJed Brown PetscFunctionReturn(0); 270*8a381b04SJed Brown } 271*8a381b04SJed Brown 272*8a381b04SJed Brown #undef __FUNCT__ 273*8a381b04SJed Brown #define __FUNCT__ "TSDestroy_ARKIMEX" 274*8a381b04SJed Brown static PetscErrorCode TSDestroy_ARKIMEX(TS ts) 275*8a381b04SJed Brown { 276*8a381b04SJed Brown PetscErrorCode ierr; 277*8a381b04SJed Brown 278*8a381b04SJed Brown PetscFunctionBegin; 279*8a381b04SJed Brown ierr = TSReset_ARKIMEX(ts);CHKERRQ(ierr); 280*8a381b04SJed Brown ierr = PetscFree(ts->data);CHKERRQ(ierr); 281*8a381b04SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXGetType_C","",PETSC_NULL);CHKERRQ(ierr); 282*8a381b04SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXSetType_C","",PETSC_NULL);CHKERRQ(ierr); 283*8a381b04SJed Brown PetscFunctionReturn(0); 284*8a381b04SJed Brown } 285*8a381b04SJed Brown 286*8a381b04SJed Brown /* 287*8a381b04SJed Brown This defines the nonlinear equation that is to be solved with SNES 288*8a381b04SJed Brown G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 289*8a381b04SJed Brown */ 290*8a381b04SJed Brown #undef __FUNCT__ 291*8a381b04SJed Brown #define __FUNCT__ "SNESTSFormFunction_ARKIMEX" 292*8a381b04SJed Brown static PetscErrorCode SNESTSFormFunction_ARKIMEX(SNES snes,Vec X,Vec F,TS ts) 293*8a381b04SJed Brown { 294*8a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 295*8a381b04SJed Brown PetscErrorCode ierr; 296*8a381b04SJed Brown 297*8a381b04SJed Brown PetscFunctionBegin; 298*8a381b04SJed Brown ierr = VecAXPBYPCZ(ark->Ydot,-ark->shift,ark->shift,0,ark->Z,X);CHKERRQ(ierr); /* Ydot = shift*(X-Z) */ 299*8a381b04SJed Brown ierr = TSComputeIFunction(ts,ark->stage_time,X,ark->Ydot,X,PETSC_TRUE);CHKERRQ(ierr); 300*8a381b04SJed Brown PetscFunctionReturn(0); 301*8a381b04SJed Brown } 302*8a381b04SJed Brown 303*8a381b04SJed Brown #undef __FUNCT__ 304*8a381b04SJed Brown #define __FUNCT__ "SNESTSFormJacobian_ARKIMEX" 305*8a381b04SJed Brown static PetscErrorCode SNESTSFormJacobian_ARKIMEX(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *str,TS ts) 306*8a381b04SJed Brown { 307*8a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 308*8a381b04SJed Brown PetscErrorCode ierr; 309*8a381b04SJed Brown 310*8a381b04SJed Brown PetscFunctionBegin; 311*8a381b04SJed Brown /* ark->Ydot has already been computed in SNESTSFormFunction_ARKIMEX (SNES guarantees this) */ 312*8a381b04SJed Brown ierr = TSComputeIJacobian(ts,ark->stage_time,X,ark->Ydot,ark->shift,A,B,str,PETSC_TRUE);CHKERRQ(ierr); 313*8a381b04SJed Brown PetscFunctionReturn(0); 314*8a381b04SJed Brown } 315*8a381b04SJed Brown 316*8a381b04SJed Brown #undef __FUNCT__ 317*8a381b04SJed Brown #define __FUNCT__ "TSSetUp_ARKIMEX" 318*8a381b04SJed Brown static PetscErrorCode TSSetUp_ARKIMEX(TS ts) 319*8a381b04SJed Brown { 320*8a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 321*8a381b04SJed Brown ARKTableau tab = ark->tableau; 322*8a381b04SJed Brown PetscInt s = tab->s; 323*8a381b04SJed Brown PetscErrorCode ierr; 324*8a381b04SJed Brown 325*8a381b04SJed Brown PetscFunctionBegin; 326*8a381b04SJed Brown if (!ark->tableau) { 327*8a381b04SJed Brown ierr = TSSetType(ts,TSARKIMEXDefault);CHKERRQ(ierr); 328*8a381b04SJed Brown } 329*8a381b04SJed Brown ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->Y);CHKERRQ(ierr); 330*8a381b04SJed Brown ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->YdotI);CHKERRQ(ierr); 331*8a381b04SJed Brown ierr = VecDuplicateVecs(ts->vec_sol,s,&ark->YdotRHS);CHKERRQ(ierr); 332*8a381b04SJed Brown ierr = VecDuplicate(ts->vec_sol,&ark->Ydot);CHKERRQ(ierr); 333*8a381b04SJed Brown ierr = VecDuplicate(ts->vec_sol,&ark->Work);CHKERRQ(ierr); 334*8a381b04SJed Brown ierr = VecDuplicate(ts->vec_sol,&ark->Z);CHKERRQ(ierr); 335*8a381b04SJed Brown ierr = PetscMalloc(s*sizeof(ark->work[0]),&ark->work);CHKERRQ(ierr); 336*8a381b04SJed Brown PetscFunctionReturn(0); 337*8a381b04SJed Brown } 338*8a381b04SJed Brown /*------------------------------------------------------------*/ 339*8a381b04SJed Brown 340*8a381b04SJed Brown #undef __FUNCT__ 341*8a381b04SJed Brown #define __FUNCT__ "TSSetFromOptions_ARKIMEX" 342*8a381b04SJed Brown static PetscErrorCode TSSetFromOptions_ARKIMEX(TS ts) 343*8a381b04SJed Brown { 344*8a381b04SJed Brown PetscErrorCode ierr; 345*8a381b04SJed Brown char arktype[256]; 346*8a381b04SJed Brown 347*8a381b04SJed Brown PetscFunctionBegin; 348*8a381b04SJed Brown ierr = PetscOptionsHead("ARKIMEX ODE solver options");CHKERRQ(ierr); 349*8a381b04SJed Brown { 350*8a381b04SJed Brown ARKTableauLink link; 351*8a381b04SJed Brown PetscInt count,choice; 352*8a381b04SJed Brown PetscBool flg; 353*8a381b04SJed Brown const char **namelist; 354*8a381b04SJed Brown ierr = PetscStrncpy(arktype,TSARKIMEXDefault,sizeof arktype);CHKERRQ(ierr); 355*8a381b04SJed Brown for (link=ARKTableauList,count=0; link; link=link->next,count++) ; 356*8a381b04SJed Brown ierr = PetscMalloc(count*sizeof(char*),&namelist);CHKERRQ(ierr); 357*8a381b04SJed Brown for (link=ARKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 358*8a381b04SJed Brown ierr = PetscOptionsEList("-ts_arkimex_type","Family of ARK IMEX method","TSARKIMEXSetType",(const char*const*)namelist,count,arktype,&choice,&flg);CHKERRQ(ierr); 359*8a381b04SJed Brown ierr = TSARKIMEXSetType(ts,flg ? namelist[choice] : arktype);CHKERRQ(ierr); 360*8a381b04SJed Brown ierr = PetscFree(namelist);CHKERRQ(ierr); 361*8a381b04SJed Brown } 362*8a381b04SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 363*8a381b04SJed Brown PetscFunctionReturn(0); 364*8a381b04SJed Brown } 365*8a381b04SJed Brown 366*8a381b04SJed Brown #undef __FUNCT__ 367*8a381b04SJed Brown #define __FUNCT__ "PetscFormatRealArray" 368*8a381b04SJed Brown static PetscErrorCode PetscFormatRealArray(char buf[],size_t len,const char *fmt,PetscInt n,const PetscReal x[]) 369*8a381b04SJed Brown { 370*8a381b04SJed Brown int i,left,count; 371*8a381b04SJed Brown char *p; 372*8a381b04SJed Brown 373*8a381b04SJed Brown PetscFunctionBegin; 374*8a381b04SJed Brown for (i=0,p=buf,left=(int)len; i<n; i++) { 375*8a381b04SJed Brown count = snprintf(p,left,fmt,x[i]); 376*8a381b04SJed Brown if (count < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"snprintf()"); 377*8a381b04SJed Brown if (count >= left) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Insufficient space in buffer"); 378*8a381b04SJed Brown left -= count; 379*8a381b04SJed Brown p += count; 380*8a381b04SJed Brown *p++ = ' '; 381*8a381b04SJed Brown } 382*8a381b04SJed Brown p[i ? 0 : -1] = 0; 383*8a381b04SJed Brown PetscFunctionReturn(0); 384*8a381b04SJed Brown } 385*8a381b04SJed Brown 386*8a381b04SJed Brown #undef __FUNCT__ 387*8a381b04SJed Brown #define __FUNCT__ "TSView_ARKIMEX" 388*8a381b04SJed Brown static PetscErrorCode TSView_ARKIMEX(TS ts,PetscViewer viewer) 389*8a381b04SJed Brown { 390*8a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 391*8a381b04SJed Brown ARKTableau tab = ark->tableau; 392*8a381b04SJed Brown PetscBool iascii; 393*8a381b04SJed Brown PetscErrorCode ierr; 394*8a381b04SJed Brown 395*8a381b04SJed Brown PetscFunctionBegin; 396*8a381b04SJed Brown ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 397*8a381b04SJed Brown if (iascii) { 398*8a381b04SJed Brown const TSARKIMEXType arktype; 399*8a381b04SJed Brown char buf[512]; 400*8a381b04SJed Brown ierr = TSARKIMEXGetType(ts,&arktype);CHKERRQ(ierr); 401*8a381b04SJed Brown ierr = PetscViewerASCIIPrintf(viewer," ARK IMEX %s\n",arktype);CHKERRQ(ierr); 402*8a381b04SJed Brown ierr = PetscFormatRealArray(buf,sizeof buf,"%8.6f",tab->s,tab->ct);CHKERRQ(ierr); 403*8a381b04SJed Brown ierr = PetscViewerASCIIPrintf(viewer," Stiff abscissa at ct = %s\n",buf);CHKERRQ(ierr); 404*8a381b04SJed Brown } else { 405*8a381b04SJed Brown SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for TS_ARKIMEX",((PetscObject)viewer)->type_name); 406*8a381b04SJed Brown } 407*8a381b04SJed Brown PetscFunctionReturn(0); 408*8a381b04SJed Brown } 409*8a381b04SJed Brown 410*8a381b04SJed Brown #undef __FUNCT__ 411*8a381b04SJed Brown #define __FUNCT__ "TSARKIMEXSetType" 412*8a381b04SJed Brown /*@C 413*8a381b04SJed Brown TSARKIMEXSetType - Set the type of ARK IMEX scheme 414*8a381b04SJed Brown 415*8a381b04SJed Brown Logically collective 416*8a381b04SJed Brown 417*8a381b04SJed Brown Input Parameter: 418*8a381b04SJed Brown + ts - timestepping context 419*8a381b04SJed Brown - arktype - type of ARK-IMEX scheme 420*8a381b04SJed Brown 421*8a381b04SJed Brown Level: intermediate 422*8a381b04SJed Brown 423*8a381b04SJed Brown .seealso: TSARKIMEXGetType() 424*8a381b04SJed Brown @*/ 425*8a381b04SJed Brown PetscErrorCode TSARKIMEXSetType(TS ts,const TSARKIMEXType arktype) 426*8a381b04SJed Brown { 427*8a381b04SJed Brown PetscErrorCode ierr; 428*8a381b04SJed Brown 429*8a381b04SJed Brown PetscFunctionBegin; 430*8a381b04SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 431*8a381b04SJed Brown ierr = PetscTryMethod(ts,"TSARKIMEXSetType_C",(TS,const TSARKIMEXType),(ts,arktype));CHKERRQ(ierr); 432*8a381b04SJed Brown PetscFunctionReturn(0); 433*8a381b04SJed Brown } 434*8a381b04SJed Brown 435*8a381b04SJed Brown #undef __FUNCT__ 436*8a381b04SJed Brown #define __FUNCT__ "TSARKIMEXGetType" 437*8a381b04SJed Brown /*@C 438*8a381b04SJed Brown TSARKIMEXGetType - Get the type of ARK IMEX scheme 439*8a381b04SJed Brown 440*8a381b04SJed Brown Logically collective 441*8a381b04SJed Brown 442*8a381b04SJed Brown Input Parameter: 443*8a381b04SJed Brown . ts - timestepping context 444*8a381b04SJed Brown 445*8a381b04SJed Brown Output Parameter: 446*8a381b04SJed Brown . arktype - type of ARK-IMEX scheme 447*8a381b04SJed Brown 448*8a381b04SJed Brown Level: intermediate 449*8a381b04SJed Brown 450*8a381b04SJed Brown .seealso: TSARKIMEXGetType() 451*8a381b04SJed Brown @*/ 452*8a381b04SJed Brown PetscErrorCode TSARKIMEXGetType(TS ts,const TSARKIMEXType *arktype) 453*8a381b04SJed Brown { 454*8a381b04SJed Brown PetscErrorCode ierr; 455*8a381b04SJed Brown 456*8a381b04SJed Brown PetscFunctionBegin; 457*8a381b04SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,1); 458*8a381b04SJed Brown ierr = PetscUseMethod(ts,"TSARKIMEXGetType_C",(TS,const TSARKIMEXType*),(ts,arktype));CHKERRQ(ierr); 459*8a381b04SJed Brown PetscFunctionReturn(0); 460*8a381b04SJed Brown } 461*8a381b04SJed Brown 462*8a381b04SJed Brown EXTERN_C_BEGIN 463*8a381b04SJed Brown #undef __FUNCT__ 464*8a381b04SJed Brown #define __FUNCT__ "TSARKIMEXGetType_ARKIMEX" 465*8a381b04SJed Brown PetscErrorCode TSARKIMEXGetType_ARKIMEX(TS ts,const TSARKIMEXType *arktype) 466*8a381b04SJed Brown { 467*8a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 468*8a381b04SJed Brown PetscErrorCode ierr; 469*8a381b04SJed Brown 470*8a381b04SJed Brown PetscFunctionBegin; 471*8a381b04SJed Brown if (!ark->tableau) {ierr = TSARKIMEXSetType(ts,TSARKIMEXDefault);CHKERRQ(ierr);} 472*8a381b04SJed Brown *arktype = ark->tableau->name; 473*8a381b04SJed Brown PetscFunctionReturn(0); 474*8a381b04SJed Brown } 475*8a381b04SJed Brown #undef __FUNCT__ 476*8a381b04SJed Brown #define __FUNCT__ "TSARKIMEXSetType_ARKIMEX" 477*8a381b04SJed Brown PetscErrorCode TSARKIMEXSetType_ARKIMEX(TS ts,const TSARKIMEXType arktype) 478*8a381b04SJed Brown { 479*8a381b04SJed Brown TS_ARKIMEX *ark = (TS_ARKIMEX*)ts->data; 480*8a381b04SJed Brown PetscErrorCode ierr; 481*8a381b04SJed Brown PetscBool match; 482*8a381b04SJed Brown ARKTableauLink link; 483*8a381b04SJed Brown 484*8a381b04SJed Brown PetscFunctionBegin; 485*8a381b04SJed Brown if (ark->tableau) { 486*8a381b04SJed Brown ierr = PetscStrcmp(ark->tableau->name,arktype,&match);CHKERRQ(ierr); 487*8a381b04SJed Brown if (match) PetscFunctionReturn(0); 488*8a381b04SJed Brown } 489*8a381b04SJed Brown for (link = ARKTableauList; link; link=link->next) { 490*8a381b04SJed Brown ierr = PetscStrcmp(link->tab.name,arktype,&match);CHKERRQ(ierr); 491*8a381b04SJed Brown if (match) { 492*8a381b04SJed Brown ierr = TSReset_ARKIMEX(ts);CHKERRQ(ierr); 493*8a381b04SJed Brown ark->tableau = &link->tab; 494*8a381b04SJed Brown PetscFunctionReturn(0); 495*8a381b04SJed Brown } 496*8a381b04SJed Brown } 497*8a381b04SJed Brown SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",arktype); 498*8a381b04SJed Brown PetscFunctionReturn(0); 499*8a381b04SJed Brown } 500*8a381b04SJed Brown EXTERN_C_END 501*8a381b04SJed Brown 502*8a381b04SJed Brown /* ------------------------------------------------------------ */ 503*8a381b04SJed Brown /*MC 504*8a381b04SJed Brown TSARKIMEX - ODE solver using Additive Runge-Kutta IMEX schemes 505*8a381b04SJed Brown 506*8a381b04SJed Brown Level: beginner 507*8a381b04SJed Brown 508*8a381b04SJed Brown .seealso: TSCreate(), TS, TSSetType() 509*8a381b04SJed Brown 510*8a381b04SJed Brown M*/ 511*8a381b04SJed Brown EXTERN_C_BEGIN 512*8a381b04SJed Brown #undef __FUNCT__ 513*8a381b04SJed Brown #define __FUNCT__ "TSCreate_ARKIMEX" 514*8a381b04SJed Brown PetscErrorCode TSCreate_ARKIMEX(TS ts) 515*8a381b04SJed Brown { 516*8a381b04SJed Brown TS_ARKIMEX *th; 517*8a381b04SJed Brown PetscErrorCode ierr; 518*8a381b04SJed Brown 519*8a381b04SJed Brown PetscFunctionBegin; 520*8a381b04SJed Brown #if !defined(PETSC_USE_DYNAMIC_LIBRARIES) 521*8a381b04SJed Brown ierr = TSARKIMEXInitializePackage(PETSC_NULL);CHKERRQ(ierr); 522*8a381b04SJed Brown #endif 523*8a381b04SJed Brown 524*8a381b04SJed Brown ts->ops->reset = TSReset_ARKIMEX; 525*8a381b04SJed Brown ts->ops->destroy = TSDestroy_ARKIMEX; 526*8a381b04SJed Brown ts->ops->view = TSView_ARKIMEX; 527*8a381b04SJed Brown ts->ops->setup = TSSetUp_ARKIMEX; 528*8a381b04SJed Brown ts->ops->step = TSStep_ARKIMEX; 529*8a381b04SJed Brown ts->ops->setfromoptions = TSSetFromOptions_ARKIMEX; 530*8a381b04SJed Brown ts->ops->snesfunction = SNESTSFormFunction_ARKIMEX; 531*8a381b04SJed Brown ts->ops->snesjacobian = SNESTSFormJacobian_ARKIMEX; 532*8a381b04SJed Brown 533*8a381b04SJed Brown ierr = PetscNewLog(ts,TS_ARKIMEX,&th);CHKERRQ(ierr); 534*8a381b04SJed Brown ts->data = (void*)th; 535*8a381b04SJed Brown 536*8a381b04SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXGetType_C","TSARKIMEXGetType_ARKIMEX",TSARKIMEXGetType_ARKIMEX);CHKERRQ(ierr); 537*8a381b04SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)ts,"TSARKIMEXSetType_C","TSARKIMEXSetType_ARKIMEX",TSARKIMEXSetType_ARKIMEX);CHKERRQ(ierr); 538*8a381b04SJed Brown PetscFunctionReturn(0); 539*8a381b04SJed Brown } 540*8a381b04SJed Brown EXTERN_C_END 541