1*3d177a5cSEmil Constantinescu 2*3d177a5cSEmil Constantinescu #include <../src/ts/impls/implicit/glle/glle.h> /*I "petscts.h" I*/ 3*3d177a5cSEmil Constantinescu #include <petscdm.h> 4*3d177a5cSEmil Constantinescu #include <petscblaslapack.h> 5*3d177a5cSEmil Constantinescu 6*3d177a5cSEmil Constantinescu static const char *TSGLLEErrorDirections[] = {"FORWARD","BACKWARD","TSGLLEErrorDirection","TSGLLEERROR_",0}; 7*3d177a5cSEmil Constantinescu static PetscFunctionList TSGLLEList; 8*3d177a5cSEmil Constantinescu static PetscFunctionList TSGLLEAcceptList; 9*3d177a5cSEmil Constantinescu static PetscBool TSGLLEPackageInitialized; 10*3d177a5cSEmil Constantinescu static PetscBool TSGLLERegisterAllCalled; 11*3d177a5cSEmil Constantinescu 12*3d177a5cSEmil Constantinescu /* This function is pure */ 13*3d177a5cSEmil Constantinescu static PetscScalar Factorial(PetscInt n) 14*3d177a5cSEmil Constantinescu { 15*3d177a5cSEmil Constantinescu PetscInt i; 16*3d177a5cSEmil Constantinescu if (n < 12) { /* Can compute with 32-bit integers */ 17*3d177a5cSEmil Constantinescu PetscInt f = 1; 18*3d177a5cSEmil Constantinescu for (i=2; i<=n; i++) f *= i; 19*3d177a5cSEmil Constantinescu return (PetscScalar)f; 20*3d177a5cSEmil Constantinescu } else { 21*3d177a5cSEmil Constantinescu PetscScalar f = 1.; 22*3d177a5cSEmil Constantinescu for (i=2; i<=n; i++) f *= (PetscScalar)i; 23*3d177a5cSEmil Constantinescu return f; 24*3d177a5cSEmil Constantinescu } 25*3d177a5cSEmil Constantinescu } 26*3d177a5cSEmil Constantinescu 27*3d177a5cSEmil Constantinescu /* This function is pure */ 28*3d177a5cSEmil Constantinescu static PetscScalar CPowF(PetscScalar c,PetscInt p) 29*3d177a5cSEmil Constantinescu { 30*3d177a5cSEmil Constantinescu return PetscPowRealInt(PetscRealPart(c),p)/Factorial(p); 31*3d177a5cSEmil Constantinescu } 32*3d177a5cSEmil Constantinescu 33*3d177a5cSEmil Constantinescu #undef __FUNCT__ 34*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEGetVecs" 35*3d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEGetVecs(TS ts,DM dm,Vec *Z,Vec *Ydotstage) 36*3d177a5cSEmil Constantinescu { 37*3d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 38*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 39*3d177a5cSEmil Constantinescu 40*3d177a5cSEmil Constantinescu PetscFunctionBegin; 41*3d177a5cSEmil Constantinescu if (Z) { 42*3d177a5cSEmil Constantinescu if (dm && dm != ts->dm) { 43*3d177a5cSEmil Constantinescu ierr = DMGetNamedGlobalVector(dm,"TSGLLE_Z",Z);CHKERRQ(ierr); 44*3d177a5cSEmil Constantinescu } else *Z = gl->Z; 45*3d177a5cSEmil Constantinescu } 46*3d177a5cSEmil Constantinescu if (Ydotstage) { 47*3d177a5cSEmil Constantinescu if (dm && dm != ts->dm) { 48*3d177a5cSEmil Constantinescu ierr = DMGetNamedGlobalVector(dm,"TSGLLE_Ydot",Ydotstage);CHKERRQ(ierr); 49*3d177a5cSEmil Constantinescu } else *Ydotstage = gl->Ydot[gl->stage]; 50*3d177a5cSEmil Constantinescu } 51*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 52*3d177a5cSEmil Constantinescu } 53*3d177a5cSEmil Constantinescu 54*3d177a5cSEmil Constantinescu 55*3d177a5cSEmil Constantinescu #undef __FUNCT__ 56*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLERestoreVecs" 57*3d177a5cSEmil Constantinescu static PetscErrorCode TSGLLERestoreVecs(TS ts,DM dm,Vec *Z,Vec *Ydotstage) 58*3d177a5cSEmil Constantinescu { 59*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 60*3d177a5cSEmil Constantinescu 61*3d177a5cSEmil Constantinescu PetscFunctionBegin; 62*3d177a5cSEmil Constantinescu if (Z) { 63*3d177a5cSEmil Constantinescu if (dm && dm != ts->dm) { 64*3d177a5cSEmil Constantinescu ierr = DMRestoreNamedGlobalVector(dm,"TSGLLE_Z",Z);CHKERRQ(ierr); 65*3d177a5cSEmil Constantinescu } 66*3d177a5cSEmil Constantinescu } 67*3d177a5cSEmil Constantinescu if (Ydotstage) { 68*3d177a5cSEmil Constantinescu 69*3d177a5cSEmil Constantinescu if (dm && dm != ts->dm) { 70*3d177a5cSEmil Constantinescu ierr = DMRestoreNamedGlobalVector(dm,"TSGLLE_Ydot",Ydotstage);CHKERRQ(ierr); 71*3d177a5cSEmil Constantinescu } 72*3d177a5cSEmil Constantinescu } 73*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 74*3d177a5cSEmil Constantinescu } 75*3d177a5cSEmil Constantinescu 76*3d177a5cSEmil Constantinescu #undef __FUNCT__ 77*3d177a5cSEmil Constantinescu #define __FUNCT__ "DMCoarsenHook_TSGLLE" 78*3d177a5cSEmil Constantinescu static PetscErrorCode DMCoarsenHook_TSGLLE(DM fine,DM coarse,void *ctx) 79*3d177a5cSEmil Constantinescu { 80*3d177a5cSEmil Constantinescu PetscFunctionBegin; 81*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 82*3d177a5cSEmil Constantinescu } 83*3d177a5cSEmil Constantinescu 84*3d177a5cSEmil Constantinescu #undef __FUNCT__ 85*3d177a5cSEmil Constantinescu #define __FUNCT__ "DMRestrictHook_TSGLLE" 86*3d177a5cSEmil Constantinescu static PetscErrorCode DMRestrictHook_TSGLLE(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 87*3d177a5cSEmil Constantinescu { 88*3d177a5cSEmil Constantinescu TS ts = (TS)ctx; 89*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 90*3d177a5cSEmil Constantinescu Vec Ydot,Ydot_c; 91*3d177a5cSEmil Constantinescu 92*3d177a5cSEmil Constantinescu PetscFunctionBegin; 93*3d177a5cSEmil Constantinescu ierr = TSGLLEGetVecs(ts,fine,NULL,&Ydot);CHKERRQ(ierr); 94*3d177a5cSEmil Constantinescu ierr = TSGLLEGetVecs(ts,coarse,NULL,&Ydot_c);CHKERRQ(ierr); 95*3d177a5cSEmil Constantinescu ierr = MatRestrict(restrct,Ydot,Ydot_c);CHKERRQ(ierr); 96*3d177a5cSEmil Constantinescu ierr = VecPointwiseMult(Ydot_c,rscale,Ydot_c);CHKERRQ(ierr); 97*3d177a5cSEmil Constantinescu ierr = TSGLLERestoreVecs(ts,fine,NULL,&Ydot);CHKERRQ(ierr); 98*3d177a5cSEmil Constantinescu ierr = TSGLLERestoreVecs(ts,coarse,NULL,&Ydot_c);CHKERRQ(ierr); 99*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 100*3d177a5cSEmil Constantinescu } 101*3d177a5cSEmil Constantinescu 102*3d177a5cSEmil Constantinescu #undef __FUNCT__ 103*3d177a5cSEmil Constantinescu #define __FUNCT__ "DMSubDomainHook_TSGLLE" 104*3d177a5cSEmil Constantinescu static PetscErrorCode DMSubDomainHook_TSGLLE(DM dm,DM subdm,void *ctx) 105*3d177a5cSEmil Constantinescu { 106*3d177a5cSEmil Constantinescu PetscFunctionBegin; 107*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 108*3d177a5cSEmil Constantinescu } 109*3d177a5cSEmil Constantinescu 110*3d177a5cSEmil Constantinescu #undef __FUNCT__ 111*3d177a5cSEmil Constantinescu #define __FUNCT__ "DMSubDomainRestrictHook_TSGLLE" 112*3d177a5cSEmil Constantinescu static PetscErrorCode DMSubDomainRestrictHook_TSGLLE(DM dm,VecScatter gscat, VecScatter lscat,DM subdm,void *ctx) 113*3d177a5cSEmil Constantinescu { 114*3d177a5cSEmil Constantinescu TS ts = (TS)ctx; 115*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 116*3d177a5cSEmil Constantinescu Vec Ydot,Ydot_s; 117*3d177a5cSEmil Constantinescu 118*3d177a5cSEmil Constantinescu PetscFunctionBegin; 119*3d177a5cSEmil Constantinescu ierr = TSGLLEGetVecs(ts,dm,NULL,&Ydot);CHKERRQ(ierr); 120*3d177a5cSEmil Constantinescu ierr = TSGLLEGetVecs(ts,subdm,NULL,&Ydot_s);CHKERRQ(ierr); 121*3d177a5cSEmil Constantinescu 122*3d177a5cSEmil Constantinescu ierr = VecScatterBegin(gscat,Ydot,Ydot_s,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 123*3d177a5cSEmil Constantinescu ierr = VecScatterEnd(gscat,Ydot,Ydot_s,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 124*3d177a5cSEmil Constantinescu 125*3d177a5cSEmil Constantinescu ierr = TSGLLERestoreVecs(ts,dm,NULL,&Ydot);CHKERRQ(ierr); 126*3d177a5cSEmil Constantinescu ierr = TSGLLERestoreVecs(ts,subdm,NULL,&Ydot_s);CHKERRQ(ierr); 127*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 128*3d177a5cSEmil Constantinescu } 129*3d177a5cSEmil Constantinescu 130*3d177a5cSEmil Constantinescu #undef __FUNCT__ 131*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLESchemeCreate" 132*3d177a5cSEmil Constantinescu static PetscErrorCode TSGLLESchemeCreate(PetscInt p,PetscInt q,PetscInt r,PetscInt s,const PetscScalar *c, 133*3d177a5cSEmil Constantinescu const PetscScalar *a,const PetscScalar *b,const PetscScalar *u,const PetscScalar *v,TSGLLEScheme *inscheme) 134*3d177a5cSEmil Constantinescu { 135*3d177a5cSEmil Constantinescu TSGLLEScheme scheme; 136*3d177a5cSEmil Constantinescu PetscInt j; 137*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 138*3d177a5cSEmil Constantinescu 139*3d177a5cSEmil Constantinescu PetscFunctionBegin; 140*3d177a5cSEmil Constantinescu if (p < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Scheme order must be positive"); 141*3d177a5cSEmil Constantinescu if (r < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"At least one item must be carried between steps"); 142*3d177a5cSEmil Constantinescu if (s < 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"At least one stage is required"); 143*3d177a5cSEmil Constantinescu PetscValidPointer(inscheme,4); 144*3d177a5cSEmil Constantinescu *inscheme = 0; 145*3d177a5cSEmil Constantinescu ierr = PetscNew(&scheme);CHKERRQ(ierr); 146*3d177a5cSEmil Constantinescu scheme->p = p; 147*3d177a5cSEmil Constantinescu scheme->q = q; 148*3d177a5cSEmil Constantinescu scheme->r = r; 149*3d177a5cSEmil Constantinescu scheme->s = s; 150*3d177a5cSEmil Constantinescu 151*3d177a5cSEmil Constantinescu ierr = PetscMalloc5(s,&scheme->c,s*s,&scheme->a,r*s,&scheme->b,r*s,&scheme->u,r*r,&scheme->v);CHKERRQ(ierr); 152*3d177a5cSEmil Constantinescu ierr = PetscMemcpy(scheme->c,c,s*sizeof(PetscScalar));CHKERRQ(ierr); 153*3d177a5cSEmil Constantinescu for (j=0; j<s*s; j++) scheme->a[j] = (PetscAbsScalar(a[j]) < 1e-12) ? 0 : a[j]; 154*3d177a5cSEmil Constantinescu for (j=0; j<r*s; j++) scheme->b[j] = (PetscAbsScalar(b[j]) < 1e-12) ? 0 : b[j]; 155*3d177a5cSEmil Constantinescu for (j=0; j<s*r; j++) scheme->u[j] = (PetscAbsScalar(u[j]) < 1e-12) ? 0 : u[j]; 156*3d177a5cSEmil Constantinescu for (j=0; j<r*r; j++) scheme->v[j] = (PetscAbsScalar(v[j]) < 1e-12) ? 0 : v[j]; 157*3d177a5cSEmil Constantinescu 158*3d177a5cSEmil Constantinescu ierr = PetscMalloc6(r,&scheme->alpha,r,&scheme->beta,r,&scheme->gamma,3*s,&scheme->phi,3*r,&scheme->psi,r,&scheme->stage_error);CHKERRQ(ierr); 159*3d177a5cSEmil Constantinescu { 160*3d177a5cSEmil Constantinescu PetscInt i,j,k,ss=s+2; 161*3d177a5cSEmil Constantinescu PetscBLASInt m,n,one=1,*ipiv,lwork=4*((s+3)*3+3),info,ldb; 162*3d177a5cSEmil Constantinescu PetscReal rcond,*sing,*workreal; 163*3d177a5cSEmil Constantinescu PetscScalar *ImV,*H,*bmat,*workscalar,*c=scheme->c,*a=scheme->a,*b=scheme->b,*u=scheme->u,*v=scheme->v; 164*3d177a5cSEmil Constantinescu #if !defined(PETSC_MISSING_LAPACK_GELSS) 165*3d177a5cSEmil Constantinescu PetscBLASInt rank; 166*3d177a5cSEmil Constantinescu #endif 167*3d177a5cSEmil Constantinescu ierr = PetscMalloc7(PetscSqr(r),&ImV,3*s,&H,3*ss,&bmat,lwork,&workscalar,5*(3+r),&workreal,r+s,&sing,r+s,&ipiv);CHKERRQ(ierr); 168*3d177a5cSEmil Constantinescu 169*3d177a5cSEmil Constantinescu /* column-major input */ 170*3d177a5cSEmil Constantinescu for (i=0; i<r-1; i++) { 171*3d177a5cSEmil Constantinescu for (j=0; j<r-1; j++) ImV[i+j*r] = 1.0*(i==j) - v[(i+1)*r+j+1]; 172*3d177a5cSEmil Constantinescu } 173*3d177a5cSEmil Constantinescu /* Build right hand side for alpha (tp - glm.B(2:end,:)*(glm.c.^(p)./factorial(p))) */ 174*3d177a5cSEmil Constantinescu for (i=1; i<r; i++) { 175*3d177a5cSEmil Constantinescu scheme->alpha[i] = 1./Factorial(p+1-i); 176*3d177a5cSEmil Constantinescu for (j=0; j<s; j++) scheme->alpha[i] -= b[i*s+j]*CPowF(c[j],p); 177*3d177a5cSEmil Constantinescu } 178*3d177a5cSEmil Constantinescu ierr = PetscBLASIntCast(r-1,&m);CHKERRQ(ierr); 179*3d177a5cSEmil Constantinescu ierr = PetscBLASIntCast(r,&n);CHKERRQ(ierr); 180*3d177a5cSEmil Constantinescu PetscStackCallBLAS("LAPACKgesv",LAPACKgesv_(&m,&one,ImV,&n,ipiv,scheme->alpha+1,&n,&info)); 181*3d177a5cSEmil Constantinescu if (info < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to GESV"); 182*3d177a5cSEmil Constantinescu if (info > 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization"); 183*3d177a5cSEmil Constantinescu 184*3d177a5cSEmil Constantinescu /* Build right hand side for beta (tp1 - glm.B(2:end,:)*(glm.c.^(p+1)./factorial(p+1)) - e.alpha) */ 185*3d177a5cSEmil Constantinescu for (i=1; i<r; i++) { 186*3d177a5cSEmil Constantinescu scheme->beta[i] = 1./Factorial(p+2-i) - scheme->alpha[i]; 187*3d177a5cSEmil Constantinescu for (j=0; j<s; j++) scheme->beta[i] -= b[i*s+j]*CPowF(c[j],p+1); 188*3d177a5cSEmil Constantinescu } 189*3d177a5cSEmil Constantinescu PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("No transpose",&m,&one,ImV,&n,ipiv,scheme->beta+1,&n,&info)); 190*3d177a5cSEmil Constantinescu if (info < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to GETRS"); 191*3d177a5cSEmil Constantinescu if (info > 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Should not happen"); 192*3d177a5cSEmil Constantinescu 193*3d177a5cSEmil Constantinescu /* Build stage_error vector 194*3d177a5cSEmil Constantinescu xi = glm.c.^(p+1)/factorial(p+1) - glm.A*glm.c.^p/factorial(p) + glm.U(:,2:end)*e.alpha; 195*3d177a5cSEmil Constantinescu */ 196*3d177a5cSEmil Constantinescu for (i=0; i<s; i++) { 197*3d177a5cSEmil Constantinescu scheme->stage_error[i] = CPowF(c[i],p+1); 198*3d177a5cSEmil Constantinescu for (j=0; j<s; j++) scheme->stage_error[i] -= a[i*s+j]*CPowF(c[j],p); 199*3d177a5cSEmil Constantinescu for (j=1; j<r; j++) scheme->stage_error[i] += u[i*r+j]*scheme->alpha[j]; 200*3d177a5cSEmil Constantinescu } 201*3d177a5cSEmil Constantinescu 202*3d177a5cSEmil Constantinescu /* alpha[0] (epsilon in B,J,W 2007) 203*3d177a5cSEmil Constantinescu epsilon = 1/factorial(p+1) - B(1,:)*c.^p/factorial(p) + V(1,2:end)*e.alpha; 204*3d177a5cSEmil Constantinescu */ 205*3d177a5cSEmil Constantinescu scheme->alpha[0] = 1./Factorial(p+1); 206*3d177a5cSEmil Constantinescu for (j=0; j<s; j++) scheme->alpha[0] -= b[0*s+j]*CPowF(c[j],p); 207*3d177a5cSEmil Constantinescu for (j=1; j<r; j++) scheme->alpha[0] += v[0*r+j]*scheme->alpha[j]; 208*3d177a5cSEmil Constantinescu 209*3d177a5cSEmil Constantinescu /* right hand side for gamma (glm.B(2:end,:)*e.xi - e.epsilon*eye(s-1,1)) */ 210*3d177a5cSEmil Constantinescu for (i=1; i<r; i++) { 211*3d177a5cSEmil Constantinescu scheme->gamma[i] = (i==1 ? -1. : 0)*scheme->alpha[0]; 212*3d177a5cSEmil Constantinescu for (j=0; j<s; j++) scheme->gamma[i] += b[i*s+j]*scheme->stage_error[j]; 213*3d177a5cSEmil Constantinescu } 214*3d177a5cSEmil Constantinescu PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("No transpose",&m,&one,ImV,&n,ipiv,scheme->gamma+1,&n,&info)); 215*3d177a5cSEmil Constantinescu if (info < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to GETRS"); 216*3d177a5cSEmil Constantinescu if (info > 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Should not happen"); 217*3d177a5cSEmil Constantinescu 218*3d177a5cSEmil Constantinescu /* beta[0] (rho in B,J,W 2007) 219*3d177a5cSEmil Constantinescu e.rho = 1/factorial(p+2) - glm.B(1,:)*glm.c.^(p+1)/factorial(p+1) ... 220*3d177a5cSEmil Constantinescu + glm.V(1,2:end)*e.beta;% - e.epsilon; 221*3d177a5cSEmil Constantinescu % Note: The paper (B,J,W 2007) includes the last term in their definition 222*3d177a5cSEmil Constantinescu * */ 223*3d177a5cSEmil Constantinescu scheme->beta[0] = 1./Factorial(p+2); 224*3d177a5cSEmil Constantinescu for (j=0; j<s; j++) scheme->beta[0] -= b[0*s+j]*CPowF(c[j],p+1); 225*3d177a5cSEmil Constantinescu for (j=1; j<r; j++) scheme->beta[0] += v[0*r+j]*scheme->beta[j]; 226*3d177a5cSEmil Constantinescu 227*3d177a5cSEmil Constantinescu /* gamma[0] (sigma in B,J,W 2007) 228*3d177a5cSEmil Constantinescu * e.sigma = glm.B(1,:)*e.xi + glm.V(1,2:end)*e.gamma; 229*3d177a5cSEmil Constantinescu * */ 230*3d177a5cSEmil Constantinescu scheme->gamma[0] = 0.0; 231*3d177a5cSEmil Constantinescu for (j=0; j<s; j++) scheme->gamma[0] += b[0*s+j]*scheme->stage_error[j]; 232*3d177a5cSEmil Constantinescu for (j=1; j<r; j++) scheme->gamma[0] += v[0*s+j]*scheme->gamma[j]; 233*3d177a5cSEmil Constantinescu 234*3d177a5cSEmil Constantinescu /* Assemble H 235*3d177a5cSEmil Constantinescu * % Determine the error estimators phi 236*3d177a5cSEmil Constantinescu H = [[cpow(glm.c,p) + C*e.alpha] [cpow(glm.c,p+1) + C*e.beta] ... 237*3d177a5cSEmil Constantinescu [e.xi - C*(e.gamma + 0*e.epsilon*eye(s-1,1))]]'; 238*3d177a5cSEmil Constantinescu % Paper has formula above without the 0, but that term must be left 239*3d177a5cSEmil Constantinescu % out to satisfy the conditions they propose and to make the 240*3d177a5cSEmil Constantinescu % example schemes work 241*3d177a5cSEmil Constantinescu e.H = H; 242*3d177a5cSEmil Constantinescu e.phi = (H \ [1 0 0;1 1 0;0 0 -1])'; 243*3d177a5cSEmil Constantinescu e.psi = -e.phi*C; 244*3d177a5cSEmil Constantinescu * */ 245*3d177a5cSEmil Constantinescu for (j=0; j<s; j++) { 246*3d177a5cSEmil Constantinescu H[0+j*3] = CPowF(c[j],p); 247*3d177a5cSEmil Constantinescu H[1+j*3] = CPowF(c[j],p+1); 248*3d177a5cSEmil Constantinescu H[2+j*3] = scheme->stage_error[j]; 249*3d177a5cSEmil Constantinescu for (k=1; k<r; k++) { 250*3d177a5cSEmil Constantinescu H[0+j*3] += CPowF(c[j],k-1)*scheme->alpha[k]; 251*3d177a5cSEmil Constantinescu H[1+j*3] += CPowF(c[j],k-1)*scheme->beta[k]; 252*3d177a5cSEmil Constantinescu H[2+j*3] -= CPowF(c[j],k-1)*scheme->gamma[k]; 253*3d177a5cSEmil Constantinescu } 254*3d177a5cSEmil Constantinescu } 255*3d177a5cSEmil Constantinescu bmat[0+0*ss] = 1.; bmat[0+1*ss] = 0.; bmat[0+2*ss] = 0.; 256*3d177a5cSEmil Constantinescu bmat[1+0*ss] = 1.; bmat[1+1*ss] = 1.; bmat[1+2*ss] = 0.; 257*3d177a5cSEmil Constantinescu bmat[2+0*ss] = 0.; bmat[2+1*ss] = 0.; bmat[2+2*ss] = -1.; 258*3d177a5cSEmil Constantinescu m = 3; 259*3d177a5cSEmil Constantinescu ierr = PetscBLASIntCast(s,&n);CHKERRQ(ierr); 260*3d177a5cSEmil Constantinescu ierr = PetscBLASIntCast(ss,&ldb);CHKERRQ(ierr); 261*3d177a5cSEmil Constantinescu rcond = 1e-12; 262*3d177a5cSEmil Constantinescu #if defined(PETSC_MISSING_LAPACK_GELSS) 263*3d177a5cSEmil Constantinescu /* ESSL does not have this routine */ 264*3d177a5cSEmil Constantinescu SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GELSS - Lapack routine is unavailable\nNot able to run GL time stepping."); 265*3d177a5cSEmil Constantinescu #else 266*3d177a5cSEmil Constantinescu #if defined(PETSC_USE_COMPLEX) 267*3d177a5cSEmil Constantinescu /* ZGELSS( M, N, NRHS, A, LDA, B, LDB, S, RCOND, RANK, WORK, LWORK, RWORK, INFO) */ 268*3d177a5cSEmil Constantinescu PetscStackCallBLAS("LAPACKgelss",LAPACKgelss_(&m,&n,&m,H,&m,bmat,&ldb,sing,&rcond,&rank,workscalar,&lwork,workreal,&info)); 269*3d177a5cSEmil Constantinescu #else 270*3d177a5cSEmil Constantinescu /* DGELSS( M, N, NRHS, A, LDA, B, LDB, S, RCOND, RANK, WORK, LWORK, INFO) */ 271*3d177a5cSEmil Constantinescu PetscStackCallBLAS("LAPACKgelss",LAPACKgelss_(&m,&n,&m,H,&m,bmat,&ldb,sing,&rcond,&rank,workscalar,&lwork,&info)); 272*3d177a5cSEmil Constantinescu #endif 273*3d177a5cSEmil Constantinescu if (info < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to GELSS"); 274*3d177a5cSEmil Constantinescu if (info > 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"SVD failed to converge"); 275*3d177a5cSEmil Constantinescu #endif 276*3d177a5cSEmil Constantinescu 277*3d177a5cSEmil Constantinescu for (j=0; j<3; j++) { 278*3d177a5cSEmil Constantinescu for (k=0; k<s; k++) scheme->phi[k+j*s] = bmat[k+j*ss]; 279*3d177a5cSEmil Constantinescu } 280*3d177a5cSEmil Constantinescu 281*3d177a5cSEmil Constantinescu /* the other part of the error estimator, psi in B,J,W 2007 */ 282*3d177a5cSEmil Constantinescu scheme->psi[0*r+0] = 0.; 283*3d177a5cSEmil Constantinescu scheme->psi[1*r+0] = 0.; 284*3d177a5cSEmil Constantinescu scheme->psi[2*r+0] = 0.; 285*3d177a5cSEmil Constantinescu for (j=1; j<r; j++) { 286*3d177a5cSEmil Constantinescu scheme->psi[0*r+j] = 0.; 287*3d177a5cSEmil Constantinescu scheme->psi[1*r+j] = 0.; 288*3d177a5cSEmil Constantinescu scheme->psi[2*r+j] = 0.; 289*3d177a5cSEmil Constantinescu for (k=0; k<s; k++) { 290*3d177a5cSEmil Constantinescu scheme->psi[0*r+j] -= CPowF(c[k],j-1)*scheme->phi[0*s+k]; 291*3d177a5cSEmil Constantinescu scheme->psi[1*r+j] -= CPowF(c[k],j-1)*scheme->phi[1*s+k]; 292*3d177a5cSEmil Constantinescu scheme->psi[2*r+j] -= CPowF(c[k],j-1)*scheme->phi[2*s+k]; 293*3d177a5cSEmil Constantinescu } 294*3d177a5cSEmil Constantinescu } 295*3d177a5cSEmil Constantinescu ierr = PetscFree7(ImV,H,bmat,workscalar,workreal,sing,ipiv);CHKERRQ(ierr); 296*3d177a5cSEmil Constantinescu } 297*3d177a5cSEmil Constantinescu /* Check which properties are satisfied */ 298*3d177a5cSEmil Constantinescu scheme->stiffly_accurate = PETSC_TRUE; 299*3d177a5cSEmil Constantinescu if (scheme->c[s-1] != 1.) scheme->stiffly_accurate = PETSC_FALSE; 300*3d177a5cSEmil Constantinescu for (j=0; j<s; j++) if (a[(s-1)*s+j] != b[j]) scheme->stiffly_accurate = PETSC_FALSE; 301*3d177a5cSEmil Constantinescu for (j=0; j<r; j++) if (u[(s-1)*r+j] != v[j]) scheme->stiffly_accurate = PETSC_FALSE; 302*3d177a5cSEmil Constantinescu scheme->fsal = scheme->stiffly_accurate; /* FSAL is stronger */ 303*3d177a5cSEmil Constantinescu for (j=0; j<s-1; j++) if (r>1 && b[1*s+j] != 0.) scheme->fsal = PETSC_FALSE; 304*3d177a5cSEmil Constantinescu if (b[1*s+r-1] != 1.) scheme->fsal = PETSC_FALSE; 305*3d177a5cSEmil Constantinescu for (j=0; j<r; j++) if (r>1 && v[1*r+j] != 0.) scheme->fsal = PETSC_FALSE; 306*3d177a5cSEmil Constantinescu 307*3d177a5cSEmil Constantinescu *inscheme = scheme; 308*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 309*3d177a5cSEmil Constantinescu } 310*3d177a5cSEmil Constantinescu 311*3d177a5cSEmil Constantinescu #undef __FUNCT__ 312*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLESchemeDestroy" 313*3d177a5cSEmil Constantinescu static PetscErrorCode TSGLLESchemeDestroy(TSGLLEScheme sc) 314*3d177a5cSEmil Constantinescu { 315*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 316*3d177a5cSEmil Constantinescu 317*3d177a5cSEmil Constantinescu PetscFunctionBegin; 318*3d177a5cSEmil Constantinescu ierr = PetscFree5(sc->c,sc->a,sc->b,sc->u,sc->v);CHKERRQ(ierr); 319*3d177a5cSEmil Constantinescu ierr = PetscFree6(sc->alpha,sc->beta,sc->gamma,sc->phi,sc->psi,sc->stage_error);CHKERRQ(ierr); 320*3d177a5cSEmil Constantinescu ierr = PetscFree(sc);CHKERRQ(ierr); 321*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 322*3d177a5cSEmil Constantinescu } 323*3d177a5cSEmil Constantinescu 324*3d177a5cSEmil Constantinescu #undef __FUNCT__ 325*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEDestroy_Default" 326*3d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEDestroy_Default(TS_GLLE *gl) 327*3d177a5cSEmil Constantinescu { 328*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 329*3d177a5cSEmil Constantinescu PetscInt i; 330*3d177a5cSEmil Constantinescu 331*3d177a5cSEmil Constantinescu PetscFunctionBegin; 332*3d177a5cSEmil Constantinescu for (i=0; i<gl->nschemes; i++) { 333*3d177a5cSEmil Constantinescu if (gl->schemes[i]) {ierr = TSGLLESchemeDestroy(gl->schemes[i]);CHKERRQ(ierr);} 334*3d177a5cSEmil Constantinescu } 335*3d177a5cSEmil Constantinescu ierr = PetscFree(gl->schemes);CHKERRQ(ierr); 336*3d177a5cSEmil Constantinescu gl->nschemes = 0; 337*3d177a5cSEmil Constantinescu ierr = PetscMemzero(gl->type_name,sizeof(gl->type_name));CHKERRQ(ierr); 338*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 339*3d177a5cSEmil Constantinescu } 340*3d177a5cSEmil Constantinescu 341*3d177a5cSEmil Constantinescu #undef __FUNCT__ 342*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEViewTable_Private" 343*3d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEViewTable_Private(PetscViewer viewer,PetscInt m,PetscInt n,const PetscScalar a[],const char name[]) 344*3d177a5cSEmil Constantinescu { 345*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 346*3d177a5cSEmil Constantinescu PetscBool iascii; 347*3d177a5cSEmil Constantinescu PetscInt i,j; 348*3d177a5cSEmil Constantinescu 349*3d177a5cSEmil Constantinescu PetscFunctionBegin; 350*3d177a5cSEmil Constantinescu ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 351*3d177a5cSEmil Constantinescu if (iascii) { 352*3d177a5cSEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer,"%30s = [",name);CHKERRQ(ierr); 353*3d177a5cSEmil Constantinescu for (i=0; i<m; i++) { 354*3d177a5cSEmil Constantinescu if (i) {ierr = PetscViewerASCIIPrintf(viewer,"%30s [","");CHKERRQ(ierr);} 355*3d177a5cSEmil Constantinescu ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 356*3d177a5cSEmil Constantinescu for (j=0; j<n; j++) { 357*3d177a5cSEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer," %12.8g",PetscRealPart(a[i*n+j]));CHKERRQ(ierr); 358*3d177a5cSEmil Constantinescu } 359*3d177a5cSEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer,"]\n");CHKERRQ(ierr); 360*3d177a5cSEmil Constantinescu ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 361*3d177a5cSEmil Constantinescu } 362*3d177a5cSEmil Constantinescu } 363*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 364*3d177a5cSEmil Constantinescu } 365*3d177a5cSEmil Constantinescu 366*3d177a5cSEmil Constantinescu 367*3d177a5cSEmil Constantinescu #undef __FUNCT__ 368*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLESchemeView" 369*3d177a5cSEmil Constantinescu static PetscErrorCode TSGLLESchemeView(TSGLLEScheme sc,PetscBool view_details,PetscViewer viewer) 370*3d177a5cSEmil Constantinescu { 371*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 372*3d177a5cSEmil Constantinescu PetscBool iascii; 373*3d177a5cSEmil Constantinescu 374*3d177a5cSEmil Constantinescu PetscFunctionBegin; 375*3d177a5cSEmil Constantinescu ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 376*3d177a5cSEmil Constantinescu if (iascii) { 377*3d177a5cSEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer,"GL scheme p,q,r,s = %d,%d,%d,%d\n",sc->p,sc->q,sc->r,sc->s);CHKERRQ(ierr); 378*3d177a5cSEmil Constantinescu ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 379*3d177a5cSEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer,"Stiffly accurate: %s, FSAL: %s\n",sc->stiffly_accurate ? "yes" : "no",sc->fsal ? "yes" : "no");CHKERRQ(ierr); 380*3d177a5cSEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer,"Leading error constants: %10.3e %10.3e %10.3e\n", 381*3d177a5cSEmil Constantinescu PetscRealPart(sc->alpha[0]),PetscRealPart(sc->beta[0]),PetscRealPart(sc->gamma[0]));CHKERRQ(ierr); 382*3d177a5cSEmil Constantinescu ierr = TSGLLEViewTable_Private(viewer,1,sc->s,sc->c,"Abscissas c");CHKERRQ(ierr); 383*3d177a5cSEmil Constantinescu if (view_details) { 384*3d177a5cSEmil Constantinescu ierr = TSGLLEViewTable_Private(viewer,sc->s,sc->s,sc->a,"A");CHKERRQ(ierr); 385*3d177a5cSEmil Constantinescu ierr = TSGLLEViewTable_Private(viewer,sc->r,sc->s,sc->b,"B");CHKERRQ(ierr); 386*3d177a5cSEmil Constantinescu ierr = TSGLLEViewTable_Private(viewer,sc->s,sc->r,sc->u,"U");CHKERRQ(ierr); 387*3d177a5cSEmil Constantinescu ierr = TSGLLEViewTable_Private(viewer,sc->r,sc->r,sc->v,"V");CHKERRQ(ierr); 388*3d177a5cSEmil Constantinescu 389*3d177a5cSEmil Constantinescu ierr = TSGLLEViewTable_Private(viewer,3,sc->s,sc->phi,"Error estimate phi");CHKERRQ(ierr); 390*3d177a5cSEmil Constantinescu ierr = TSGLLEViewTable_Private(viewer,3,sc->r,sc->psi,"Error estimate psi");CHKERRQ(ierr); 391*3d177a5cSEmil Constantinescu ierr = TSGLLEViewTable_Private(viewer,1,sc->r,sc->alpha,"Modify alpha");CHKERRQ(ierr); 392*3d177a5cSEmil Constantinescu ierr = TSGLLEViewTable_Private(viewer,1,sc->r,sc->beta,"Modify beta");CHKERRQ(ierr); 393*3d177a5cSEmil Constantinescu ierr = TSGLLEViewTable_Private(viewer,1,sc->r,sc->gamma,"Modify gamma");CHKERRQ(ierr); 394*3d177a5cSEmil Constantinescu ierr = TSGLLEViewTable_Private(viewer,1,sc->s,sc->stage_error,"Stage error xi");CHKERRQ(ierr); 395*3d177a5cSEmil Constantinescu } 396*3d177a5cSEmil Constantinescu ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 397*3d177a5cSEmil Constantinescu } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported",((PetscObject)viewer)->type_name); 398*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 399*3d177a5cSEmil Constantinescu } 400*3d177a5cSEmil Constantinescu 401*3d177a5cSEmil Constantinescu #undef __FUNCT__ 402*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEEstimateHigherMoments_Default" 403*3d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEEstimateHigherMoments_Default(TSGLLEScheme sc,PetscReal h,Vec Ydot[],Vec Xold[],Vec hm[]) 404*3d177a5cSEmil Constantinescu { 405*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 406*3d177a5cSEmil Constantinescu PetscInt i; 407*3d177a5cSEmil Constantinescu 408*3d177a5cSEmil Constantinescu PetscFunctionBegin; 409*3d177a5cSEmil Constantinescu if (sc->r > 64 || sc->s > 64) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Ridiculous number of stages or items passed between stages"); 410*3d177a5cSEmil Constantinescu /* build error vectors*/ 411*3d177a5cSEmil Constantinescu for (i=0; i<3; i++) { 412*3d177a5cSEmil Constantinescu PetscScalar phih[64]; 413*3d177a5cSEmil Constantinescu PetscInt j; 414*3d177a5cSEmil Constantinescu for (j=0; j<sc->s; j++) phih[j] = sc->phi[i*sc->s+j]*h; 415*3d177a5cSEmil Constantinescu ierr = VecZeroEntries(hm[i]);CHKERRQ(ierr); 416*3d177a5cSEmil Constantinescu ierr = VecMAXPY(hm[i],sc->s,phih,Ydot);CHKERRQ(ierr); 417*3d177a5cSEmil Constantinescu ierr = VecMAXPY(hm[i],sc->r,&sc->psi[i*sc->r],Xold);CHKERRQ(ierr); 418*3d177a5cSEmil Constantinescu } 419*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 420*3d177a5cSEmil Constantinescu } 421*3d177a5cSEmil Constantinescu 422*3d177a5cSEmil Constantinescu #undef __FUNCT__ 423*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLECompleteStep_Rescale" 424*3d177a5cSEmil Constantinescu static PetscErrorCode TSGLLECompleteStep_Rescale(TSGLLEScheme sc,PetscReal h,TSGLLEScheme next_sc,PetscReal next_h,Vec Ydot[],Vec Xold[],Vec X[]) 425*3d177a5cSEmil Constantinescu { 426*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 427*3d177a5cSEmil Constantinescu PetscScalar brow[32],vrow[32]; 428*3d177a5cSEmil Constantinescu PetscInt i,j,r,s; 429*3d177a5cSEmil Constantinescu 430*3d177a5cSEmil Constantinescu PetscFunctionBegin; 431*3d177a5cSEmil Constantinescu /* Build the new solution from (X,Ydot) */ 432*3d177a5cSEmil Constantinescu r = sc->r; 433*3d177a5cSEmil Constantinescu s = sc->s; 434*3d177a5cSEmil Constantinescu for (i=0; i<r; i++) { 435*3d177a5cSEmil Constantinescu ierr = VecZeroEntries(X[i]);CHKERRQ(ierr); 436*3d177a5cSEmil Constantinescu for (j=0; j<s; j++) brow[j] = h*sc->b[i*s+j]; 437*3d177a5cSEmil Constantinescu ierr = VecMAXPY(X[i],s,brow,Ydot);CHKERRQ(ierr); 438*3d177a5cSEmil Constantinescu for (j=0; j<r; j++) vrow[j] = sc->v[i*r+j]; 439*3d177a5cSEmil Constantinescu ierr = VecMAXPY(X[i],r,vrow,Xold);CHKERRQ(ierr); 440*3d177a5cSEmil Constantinescu } 441*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 442*3d177a5cSEmil Constantinescu } 443*3d177a5cSEmil Constantinescu 444*3d177a5cSEmil Constantinescu #undef __FUNCT__ 445*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLECompleteStep_RescaleAndModify" 446*3d177a5cSEmil Constantinescu static PetscErrorCode TSGLLECompleteStep_RescaleAndModify(TSGLLEScheme sc,PetscReal h,TSGLLEScheme next_sc,PetscReal next_h,Vec Ydot[],Vec Xold[],Vec X[]) 447*3d177a5cSEmil Constantinescu { 448*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 449*3d177a5cSEmil Constantinescu PetscScalar brow[32],vrow[32]; 450*3d177a5cSEmil Constantinescu PetscReal ratio; 451*3d177a5cSEmil Constantinescu PetscInt i,j,p,r,s; 452*3d177a5cSEmil Constantinescu 453*3d177a5cSEmil Constantinescu PetscFunctionBegin; 454*3d177a5cSEmil Constantinescu /* Build the new solution from (X,Ydot) */ 455*3d177a5cSEmil Constantinescu p = sc->p; 456*3d177a5cSEmil Constantinescu r = sc->r; 457*3d177a5cSEmil Constantinescu s = sc->s; 458*3d177a5cSEmil Constantinescu ratio = next_h/h; 459*3d177a5cSEmil Constantinescu for (i=0; i<r; i++) { 460*3d177a5cSEmil Constantinescu ierr = VecZeroEntries(X[i]);CHKERRQ(ierr); 461*3d177a5cSEmil Constantinescu for (j=0; j<s; j++) { 462*3d177a5cSEmil Constantinescu brow[j] = h*(PetscPowRealInt(ratio,i)*sc->b[i*s+j] 463*3d177a5cSEmil Constantinescu + (PetscPowRealInt(ratio,i) - PetscPowRealInt(ratio,p+1))*(+ sc->alpha[i]*sc->phi[0*s+j]) 464*3d177a5cSEmil Constantinescu + (PetscPowRealInt(ratio,i) - PetscPowRealInt(ratio,p+2))*(+ sc->beta [i]*sc->phi[1*s+j] 465*3d177a5cSEmil Constantinescu + sc->gamma[i]*sc->phi[2*s+j])); 466*3d177a5cSEmil Constantinescu } 467*3d177a5cSEmil Constantinescu ierr = VecMAXPY(X[i],s,brow,Ydot);CHKERRQ(ierr); 468*3d177a5cSEmil Constantinescu for (j=0; j<r; j++) { 469*3d177a5cSEmil Constantinescu vrow[j] = (PetscPowRealInt(ratio,i)*sc->v[i*r+j] 470*3d177a5cSEmil Constantinescu + (PetscPowRealInt(ratio,i) - PetscPowRealInt(ratio,p+1))*(+ sc->alpha[i]*sc->psi[0*r+j]) 471*3d177a5cSEmil Constantinescu + (PetscPowRealInt(ratio,i) - PetscPowRealInt(ratio,p+2))*(+ sc->beta [i]*sc->psi[1*r+j] 472*3d177a5cSEmil Constantinescu + sc->gamma[i]*sc->psi[2*r+j])); 473*3d177a5cSEmil Constantinescu } 474*3d177a5cSEmil Constantinescu ierr = VecMAXPY(X[i],r,vrow,Xold);CHKERRQ(ierr); 475*3d177a5cSEmil Constantinescu } 476*3d177a5cSEmil Constantinescu if (r < next_sc->r) { 477*3d177a5cSEmil Constantinescu if (r+1 != next_sc->r) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot accommodate jump in r greater than 1"); 478*3d177a5cSEmil Constantinescu ierr = VecZeroEntries(X[r]);CHKERRQ(ierr); 479*3d177a5cSEmil Constantinescu for (j=0; j<s; j++) brow[j] = h*PetscPowRealInt(ratio,p+1)*sc->phi[0*s+j]; 480*3d177a5cSEmil Constantinescu ierr = VecMAXPY(X[r],s,brow,Ydot);CHKERRQ(ierr); 481*3d177a5cSEmil Constantinescu for (j=0; j<r; j++) vrow[j] = PetscPowRealInt(ratio,p+1)*sc->psi[0*r+j]; 482*3d177a5cSEmil Constantinescu ierr = VecMAXPY(X[r],r,vrow,Xold);CHKERRQ(ierr); 483*3d177a5cSEmil Constantinescu } 484*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 485*3d177a5cSEmil Constantinescu } 486*3d177a5cSEmil Constantinescu 487*3d177a5cSEmil Constantinescu #undef __FUNCT__ 488*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLECreate_IRKS" 489*3d177a5cSEmil Constantinescu static PetscErrorCode TSGLLECreate_IRKS(TS ts) 490*3d177a5cSEmil Constantinescu { 491*3d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 492*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 493*3d177a5cSEmil Constantinescu 494*3d177a5cSEmil Constantinescu PetscFunctionBegin; 495*3d177a5cSEmil Constantinescu gl->Destroy = TSGLLEDestroy_Default; 496*3d177a5cSEmil Constantinescu gl->EstimateHigherMoments = TSGLLEEstimateHigherMoments_Default; 497*3d177a5cSEmil Constantinescu gl->CompleteStep = TSGLLECompleteStep_RescaleAndModify; 498*3d177a5cSEmil Constantinescu ierr = PetscMalloc1(10,&gl->schemes);CHKERRQ(ierr); 499*3d177a5cSEmil Constantinescu gl->nschemes = 0; 500*3d177a5cSEmil Constantinescu 501*3d177a5cSEmil Constantinescu { 502*3d177a5cSEmil Constantinescu /* p=1,q=1, r=s=2, A- and L-stable with error estimates of order 2 and 3 503*3d177a5cSEmil Constantinescu * Listed in Butcher & Podhaisky 2006. On error estimation in general linear methods for stiff ODE. 504*3d177a5cSEmil Constantinescu * irks(0.3,0,[.3,1],[1],1) 505*3d177a5cSEmil Constantinescu * Note: can be made to have classical order (not stage order) 2 by replacing 0.3 with 1-sqrt(1/2) 506*3d177a5cSEmil Constantinescu * but doing so would sacrifice the error estimator. 507*3d177a5cSEmil Constantinescu */ 508*3d177a5cSEmil Constantinescu const PetscScalar c[2] = {3./10., 1.}; 509*3d177a5cSEmil Constantinescu const PetscScalar a[2][2] = {{3./10., 0}, {7./10., 3./10.}}; 510*3d177a5cSEmil Constantinescu const PetscScalar b[2][2] = {{7./10., 3./10.}, {0,1}}; 511*3d177a5cSEmil Constantinescu const PetscScalar u[2][2] = {{1,0},{1,0}}; 512*3d177a5cSEmil Constantinescu const PetscScalar v[2][2] = {{1,0},{0,0}}; 513*3d177a5cSEmil Constantinescu ierr = TSGLLESchemeCreate(1,1,2,2,c,*a,*b,*u,*v,&gl->schemes[gl->nschemes++]);CHKERRQ(ierr); 514*3d177a5cSEmil Constantinescu } 515*3d177a5cSEmil Constantinescu 516*3d177a5cSEmil Constantinescu { 517*3d177a5cSEmil Constantinescu /* p=q=2, r=s=3: irks(4/9,0,[1:3]/3,[0.33852],1) */ 518*3d177a5cSEmil Constantinescu /* http://www.math.auckland.ac.nz/~hpod/atlas/i2a.html */ 519*3d177a5cSEmil Constantinescu const PetscScalar c[3] = {1./3., 2./3., 1} 520*3d177a5cSEmil Constantinescu ,a[3][3] = {{4./9. ,0 , 0}, 521*3d177a5cSEmil Constantinescu {1.03750643704090e+00 , 4./9., 0}, 522*3d177a5cSEmil Constantinescu {7.67024779410304e-01 , -3.81140216918943e-01, 4./9.}} 523*3d177a5cSEmil Constantinescu ,b[3][3] = {{0.767024779410304, -0.381140216918943, 4./9.}, 524*3d177a5cSEmil Constantinescu {0.000000000000000, 0.000000000000000, 1.000000000000000}, 525*3d177a5cSEmil Constantinescu {-2.075048385225385, 0.621728385225383, 1.277197204924873}} 526*3d177a5cSEmil Constantinescu ,u[3][3] = {{1.0000000000000000, -0.1111111111111109, -0.0925925925925922}, 527*3d177a5cSEmil Constantinescu {1.0000000000000000, -0.8152842148186744, -0.4199095530877056}, 528*3d177a5cSEmil Constantinescu {1.0000000000000000, 0.1696709930641948, 0.0539741070314165}} 529*3d177a5cSEmil Constantinescu ,v[3][3] = {{1.0000000000000000, 0.1696709930641948, 0.0539741070314165}, 530*3d177a5cSEmil Constantinescu {0.000000000000000, 0.000000000000000, 0.000000000000000}, 531*3d177a5cSEmil Constantinescu {0.000000000000000, 0.176122795075129, 0.000000000000000}}; 532*3d177a5cSEmil Constantinescu ierr = TSGLLESchemeCreate(2,2,3,3,c,*a,*b,*u,*v,&gl->schemes[gl->nschemes++]);CHKERRQ(ierr); 533*3d177a5cSEmil Constantinescu } 534*3d177a5cSEmil Constantinescu { 535*3d177a5cSEmil Constantinescu /* p=q=3, r=s=4: irks(9/40,0,[1:4]/4,[0.3312 1.0050],[0.49541 1;1 0]) */ 536*3d177a5cSEmil Constantinescu const PetscScalar c[4] = {0.25,0.5,0.75,1.0} 537*3d177a5cSEmil Constantinescu ,a[4][4] = {{9./40. , 0, 0, 0}, 538*3d177a5cSEmil Constantinescu {2.11286958887701e-01 , 9./40. , 0, 0}, 539*3d177a5cSEmil Constantinescu {9.46338294287584e-01 , -3.42942861246094e-01, 9./40. , 0}, 540*3d177a5cSEmil Constantinescu {0.521490453970721 , -0.662474225622980, 0.490476425459734, 9./40. }} 541*3d177a5cSEmil Constantinescu ,b[4][4] = {{0.521490453970721 , -0.662474225622980, 0.490476425459734, 9./40. }, 542*3d177a5cSEmil Constantinescu {0.000000000000000 , 0.000000000000000, 0.000000000000000, 1.000000000000000}, 543*3d177a5cSEmil Constantinescu {-0.084677029310348 , 1.390757514776085, -1.568157386206001, 2.023192696767826}, 544*3d177a5cSEmil Constantinescu {0.465383797936408 , 1.478273530625148, -1.930836081010182, 1.644872111193354}} 545*3d177a5cSEmil Constantinescu ,u[4][4] = {{1.00000000000000000 , 0.02500000000001035, -0.02499999999999053, -0.00442708333332865}, 546*3d177a5cSEmil Constantinescu {1.00000000000000000 , 0.06371304111232945, -0.04032173972189845, -0.01389438413189452}, 547*3d177a5cSEmil Constantinescu {1.00000000000000000 , -0.07839543304147778, 0.04738685705116663, 0.02032603595928376}, 548*3d177a5cSEmil Constantinescu {1.00000000000000000 , 0.42550734619251651, 0.10800718022400080, -0.01726712647760034}} 549*3d177a5cSEmil Constantinescu ,v[4][4] = {{1.00000000000000000 , 0.42550734619251651, 0.10800718022400080, -0.01726712647760034}, 550*3d177a5cSEmil Constantinescu {0.000000000000000 , 0.000000000000000, 0.000000000000000, 0.000000000000000}, 551*3d177a5cSEmil Constantinescu {0.000000000000000 , -1.761115796027561, -0.521284157173780, 0.258249384305463}, 552*3d177a5cSEmil Constantinescu {0.000000000000000 , -1.657693358744728, -1.052227765232394, 0.521284157173780}}; 553*3d177a5cSEmil Constantinescu ierr = TSGLLESchemeCreate(3,3,4,4,c,*a,*b,*u,*v,&gl->schemes[gl->nschemes++]);CHKERRQ(ierr); 554*3d177a5cSEmil Constantinescu } 555*3d177a5cSEmil Constantinescu { 556*3d177a5cSEmil Constantinescu /* p=q=4, r=s=5: 557*3d177a5cSEmil Constantinescu irks(3/11,0,[1:5]/5, [0.1715 -0.1238 0.6617],... 558*3d177a5cSEmil Constantinescu [ -0.0812 0.4079 1.0000 559*3d177a5cSEmil Constantinescu 1.0000 0 0 560*3d177a5cSEmil Constantinescu 0.8270 1.0000 0]) 561*3d177a5cSEmil Constantinescu */ 562*3d177a5cSEmil Constantinescu const PetscScalar c[5] = {0.2,0.4,0.6,0.8,1.0} 563*3d177a5cSEmil Constantinescu ,a[5][5] = {{2.72727272727352e-01 , 0.00000000000000e+00, 0.00000000000000e+00 , 0.00000000000000e+00 , 0.00000000000000e+00}, 564*3d177a5cSEmil Constantinescu {-1.03980153733431e-01, 2.72727272727405e-01, 0.00000000000000e+00, 0.00000000000000e+00 , 0.00000000000000e+00}, 565*3d177a5cSEmil Constantinescu {-1.58615400341492e+00, 7.44168951881122e-01, 2.72727272727309e-01, 0.00000000000000e+00 , 0.00000000000000e+00}, 566*3d177a5cSEmil Constantinescu {-8.73658042865628e-01, 5.37884671894595e-01, -1.63298538799523e-01, 2.72727272726996e-01 , 0.00000000000000e+00}, 567*3d177a5cSEmil Constantinescu {2.95489397443992e-01 , -1.18481693910097e+00 , -6.68029812659953e-01 , 1.00716687860943e+00 , 2.72727272727288e-01}} 568*3d177a5cSEmil Constantinescu ,b[5][5] = {{2.95489397443992e-01 , -1.18481693910097e+00 , -6.68029812659953e-01 , 1.00716687860943e+00 , 2.72727272727288e-01}, 569*3d177a5cSEmil Constantinescu {0.00000000000000e+00 , 1.11022302462516e-16 , -2.22044604925031e-16 , 0.00000000000000e+00 , 1.00000000000000e+00}, 570*3d177a5cSEmil Constantinescu {-4.05882503986005e+00, -4.00924006567769e+00, -1.38930610972481e+00, 4.45223930308488e+00 , 6.32331093108427e-01}, 571*3d177a5cSEmil Constantinescu {8.35690179937017e+00 , -2.26640927349732e+00 , 6.86647884973826e+00 , -5.22595158025740e+00 , 4.50893068837431e+00}, 572*3d177a5cSEmil Constantinescu {1.27656267027479e+01 , 2.80882153840821e+00 , 8.91173096522890e+00 , -1.07936444078906e+01 , 4.82534148988854e+00}} 573*3d177a5cSEmil Constantinescu ,u[5][5] = {{1.00000000000000e+00 , -7.27272727273551e-02 , -3.45454545454419e-02 , -4.12121212119565e-03 ,-2.96969696964014e-04}, 574*3d177a5cSEmil Constantinescu {1.00000000000000e+00 , 2.31252881006154e-01 , -8.29487834416481e-03 , -9.07191207681020e-03 ,-1.70378403743473e-03}, 575*3d177a5cSEmil Constantinescu {1.00000000000000e+00 , 1.16925777880663e+00 , 3.59268562942635e-02 , -4.09013451730615e-02 ,-1.02411119670164e-02}, 576*3d177a5cSEmil Constantinescu {1.00000000000000e+00 , 1.02634463704356e+00 , 1.59375044913405e-01 , 1.89673015035370e-03 ,-4.89987231897569e-03}, 577*3d177a5cSEmil Constantinescu {1.00000000000000e+00 , 1.27746320298021e+00 , 2.37186008132728e-01 , -8.28694373940065e-02 ,-5.34396510196430e-02}} 578*3d177a5cSEmil Constantinescu ,v[5][5] = {{1.00000000000000e+00 , 1.27746320298021e+00 , 2.37186008132728e-01 , -8.28694373940065e-02 ,-5.34396510196430e-02}, 579*3d177a5cSEmil Constantinescu {0.00000000000000e+00 , -1.77635683940025e-15 , -1.99840144432528e-15 , -9.99200722162641e-16 ,-3.33066907387547e-16}, 580*3d177a5cSEmil Constantinescu {0.00000000000000e+00 , 4.37280081906924e+00 , 5.49221645016377e-02 , -8.88913177394943e-02 , 1.12879077989154e-01}, 581*3d177a5cSEmil Constantinescu {0.00000000000000e+00 , -1.22399504837280e+01 , -5.21287338448645e+00 , -8.03952325565291e-01 , 4.60298678047147e-01}, 582*3d177a5cSEmil Constantinescu {0.00000000000000e+00 , -1.85178762883829e+01 , -5.21411849862624e+00 , -1.04283436528809e+00 , 7.49030161063651e-01}}; 583*3d177a5cSEmil Constantinescu ierr = TSGLLESchemeCreate(4,4,5,5,c,*a,*b,*u,*v,&gl->schemes[gl->nschemes++]);CHKERRQ(ierr); 584*3d177a5cSEmil Constantinescu } 585*3d177a5cSEmil Constantinescu { 586*3d177a5cSEmil Constantinescu /* p=q=5, r=s=6; 587*3d177a5cSEmil Constantinescu irks(1/3,0,[1:6]/6,... 588*3d177a5cSEmil Constantinescu [-0.0489 0.4228 -0.8814 0.9021],... 589*3d177a5cSEmil Constantinescu [-0.3474 -0.6617 0.6294 0.2129 590*3d177a5cSEmil Constantinescu 0.0044 -0.4256 -0.1427 -0.8936 591*3d177a5cSEmil Constantinescu -0.8267 0.4821 0.1371 -0.2557 592*3d177a5cSEmil Constantinescu -0.4426 -0.3855 -0.7514 0.3014]) 593*3d177a5cSEmil Constantinescu */ 594*3d177a5cSEmil Constantinescu const PetscScalar c[6] = {1./6, 2./6, 3./6, 4./6, 5./6, 1.} 595*3d177a5cSEmil Constantinescu ,a[6][6] = {{ 3.33333333333940e-01, 0 , 0 , 0 , 0 , 0 }, 596*3d177a5cSEmil Constantinescu { -8.64423857333350e-02, 3.33333333332888e-01, 0 , 0 , 0 , 0 }, 597*3d177a5cSEmil Constantinescu { -2.16850174258252e+00, -2.23619072028839e+00, 3.33333333335204e-01, 0 , 0 , 0 }, 598*3d177a5cSEmil Constantinescu { -4.73160970138997e+00, -3.89265344629268e+00, -2.76318716520933e-01, 3.33333333335759e-01, 0 , 0 }, 599*3d177a5cSEmil Constantinescu { -6.75187540297338e+00, -7.90756533769377e+00, 7.90245051802259e-01, -4.48352364517632e-01, 3.33333333328483e-01, 0 }, 600*3d177a5cSEmil Constantinescu { -4.26488287921548e+00, -1.19320395589302e+01, 3.38924509887755e+00, -2.23969848002481e+00, 6.62807710124007e-01, 3.33333333335440e-01}} 601*3d177a5cSEmil Constantinescu ,b[6][6] = {{ -4.26488287921548e+00, -1.19320395589302e+01, 3.38924509887755e+00, -2.23969848002481e+00, 6.62807710124007e-01, 3.33333333335440e-01}, 602*3d177a5cSEmil Constantinescu { -8.88178419700125e-16, 4.44089209850063e-16, -1.54737334057131e-15, -8.88178419700125e-16, 0.00000000000000e+00, 1.00000000000001e+00}, 603*3d177a5cSEmil Constantinescu { -2.87780425770651e+01, -1.13520448264971e+01, 2.62002318943161e+01, 2.56943874812797e+01, -3.06702268304488e+01, 6.68067773510103e+00}, 604*3d177a5cSEmil Constantinescu { 5.47971245256474e+01, 6.80366875868284e+01, -6.50952588861999e+01, -8.28643975339097e+01, 8.17416943896414e+01, -1.17819043489036e+01}, 605*3d177a5cSEmil Constantinescu { -2.33332114788869e+02, 6.12942539462634e+01, -4.91850135865944e+01, 1.82716844135480e+02, -1.29788173979395e+02, 3.09968095651099e+01}, 606*3d177a5cSEmil Constantinescu { -1.72049132343751e+02, 8.60194713593999e+00, 7.98154219170200e-01, 1.50371386053218e+02, -1.18515423962066e+02, 2.50898277784663e+01}} 607*3d177a5cSEmil Constantinescu ,u[6][6] = {{ 1.00000000000000e+00, -1.66666666666870e-01, -4.16666666664335e-02, -3.85802469124815e-03, -2.25051440302250e-04, -9.64506172339142e-06}, 608*3d177a5cSEmil Constantinescu { 1.00000000000000e+00, 8.64423857327162e-02, -4.11484912671353e-02, -1.11450903217645e-02, -1.47651050487126e-03, -1.34395070766826e-04}, 609*3d177a5cSEmil Constantinescu { 1.00000000000000e+00, 4.57135912953434e+00, 1.06514719719137e+00, 1.33517564218007e-01, 1.11365952968659e-02, 6.12382756769504e-04}, 610*3d177a5cSEmil Constantinescu { 1.00000000000000e+00, 9.23391519753404e+00, 2.22431212392095e+00, 2.91823807741891e-01, 2.52058456411084e-02, 1.22800542949647e-03}, 611*3d177a5cSEmil Constantinescu { 1.00000000000000e+00, 1.48175480533865e+01, 3.73439117461835e+00, 5.14648336541804e-01, 4.76430038853402e-02, 2.56798515502156e-03}, 612*3d177a5cSEmil Constantinescu { 1.00000000000000e+00, 1.50512347758335e+01, 4.10099701165164e+00, 5.66039141003603e-01, 3.91213893800891e-02, -2.99136269067853e-03}} 613*3d177a5cSEmil Constantinescu ,v[6][6] = {{ 1.00000000000000e+00, 1.50512347758335e+01, 4.10099701165164e+00, 5.66039141003603e-01, 3.91213893800891e-02, -2.99136269067853e-03}, 614*3d177a5cSEmil Constantinescu { 0.00000000000000e+00, -4.88498130835069e-15, -6.43929354282591e-15, -3.55271367880050e-15, -1.22124532708767e-15, -3.12250225675825e-16}, 615*3d177a5cSEmil Constantinescu { 0.00000000000000e+00, 1.22250171233141e+01, -1.77150760606169e+00, 3.54516769879390e-01, 6.22298845883398e-01, 2.31647447450276e-01}, 616*3d177a5cSEmil Constantinescu { 0.00000000000000e+00, -4.48339457331040e+01, -3.57363126641880e-01, 5.18750173123425e-01, 6.55727990241799e-02, 1.63175368287079e-01}, 617*3d177a5cSEmil Constantinescu { 0.00000000000000e+00, 1.37297394708005e+02, -1.60145272991317e+00, -5.05319555199441e+00, 1.55328940390990e-01, 9.16629423682464e-01}, 618*3d177a5cSEmil Constantinescu { 0.00000000000000e+00, 1.05703241119022e+02, -1.16610260983038e+00, -2.99767252773859e+00, -1.13472315553890e-01, 1.09742849254729e+00}}; 619*3d177a5cSEmil Constantinescu ierr = TSGLLESchemeCreate(5,5,6,6,c,*a,*b,*u,*v,&gl->schemes[gl->nschemes++]);CHKERRQ(ierr); 620*3d177a5cSEmil Constantinescu } 621*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 622*3d177a5cSEmil Constantinescu } 623*3d177a5cSEmil Constantinescu 624*3d177a5cSEmil Constantinescu #undef __FUNCT__ 625*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLESetType" 626*3d177a5cSEmil Constantinescu /*@C 627*3d177a5cSEmil Constantinescu TSGLLESetType - sets the class of general linear method to use for time-stepping 628*3d177a5cSEmil Constantinescu 629*3d177a5cSEmil Constantinescu Collective on TS 630*3d177a5cSEmil Constantinescu 631*3d177a5cSEmil Constantinescu Input Parameters: 632*3d177a5cSEmil Constantinescu + ts - the TS context 633*3d177a5cSEmil Constantinescu - type - a method 634*3d177a5cSEmil Constantinescu 635*3d177a5cSEmil Constantinescu Options Database Key: 636*3d177a5cSEmil Constantinescu . -ts_gl_type <type> - sets the method, use -help for a list of available method (e.g. irks) 637*3d177a5cSEmil Constantinescu 638*3d177a5cSEmil Constantinescu Notes: 639*3d177a5cSEmil Constantinescu See "petsc/include/petscts.h" for available methods (for instance) 640*3d177a5cSEmil Constantinescu . TSGLLE_IRKS - Diagonally implicit methods with inherent Runge-Kutta stability (for stiff problems) 641*3d177a5cSEmil Constantinescu 642*3d177a5cSEmil Constantinescu Normally, it is best to use the TSSetFromOptions() command and 643*3d177a5cSEmil Constantinescu then set the TSGLLE type from the options database rather than by using 644*3d177a5cSEmil Constantinescu this routine. Using the options database provides the user with 645*3d177a5cSEmil Constantinescu maximum flexibility in evaluating the many different solvers. 646*3d177a5cSEmil Constantinescu The TSGLLESetType() routine is provided for those situations where it 647*3d177a5cSEmil Constantinescu is necessary to set the timestepping solver independently of the 648*3d177a5cSEmil Constantinescu command line or options database. This might be the case, for example, 649*3d177a5cSEmil Constantinescu when the choice of solver changes during the execution of the 650*3d177a5cSEmil Constantinescu program, and the user's application is taking responsibility for 651*3d177a5cSEmil Constantinescu choosing the appropriate method. 652*3d177a5cSEmil Constantinescu 653*3d177a5cSEmil Constantinescu Level: intermediate 654*3d177a5cSEmil Constantinescu 655*3d177a5cSEmil Constantinescu .keywords: TS, TSGLLE, set, type 656*3d177a5cSEmil Constantinescu @*/ 657*3d177a5cSEmil Constantinescu PetscErrorCode TSGLLESetType(TS ts,TSGLLEType type) 658*3d177a5cSEmil Constantinescu { 659*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 660*3d177a5cSEmil Constantinescu 661*3d177a5cSEmil Constantinescu PetscFunctionBegin; 662*3d177a5cSEmil Constantinescu PetscValidHeaderSpecific(ts,TS_CLASSID,1); 663*3d177a5cSEmil Constantinescu PetscValidCharPointer(type,2); 664*3d177a5cSEmil Constantinescu ierr = PetscTryMethod(ts,"TSGLLESetType_C",(TS,TSGLLEType),(ts,type));CHKERRQ(ierr); 665*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 666*3d177a5cSEmil Constantinescu } 667*3d177a5cSEmil Constantinescu 668*3d177a5cSEmil Constantinescu #undef __FUNCT__ 669*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLESetAcceptType" 670*3d177a5cSEmil Constantinescu /*@C 671*3d177a5cSEmil Constantinescu TSGLLESetAcceptType - sets the acceptance test 672*3d177a5cSEmil Constantinescu 673*3d177a5cSEmil Constantinescu Time integrators that need to control error must have the option to reject a time step based on local error 674*3d177a5cSEmil Constantinescu estimates. This function allows different schemes to be set. 675*3d177a5cSEmil Constantinescu 676*3d177a5cSEmil Constantinescu Logically Collective on TS 677*3d177a5cSEmil Constantinescu 678*3d177a5cSEmil Constantinescu Input Parameters: 679*3d177a5cSEmil Constantinescu + ts - the TS context 680*3d177a5cSEmil Constantinescu - type - the type 681*3d177a5cSEmil Constantinescu 682*3d177a5cSEmil Constantinescu Options Database Key: 683*3d177a5cSEmil Constantinescu . -ts_gl_accept_type <type> - sets the method used to determine whether to accept or reject a step 684*3d177a5cSEmil Constantinescu 685*3d177a5cSEmil Constantinescu Level: intermediate 686*3d177a5cSEmil Constantinescu 687*3d177a5cSEmil Constantinescu .seealso: TS, TSGLLE, TSGLLEAcceptRegister(), TSGLLEAdapt, set type 688*3d177a5cSEmil Constantinescu @*/ 689*3d177a5cSEmil Constantinescu PetscErrorCode TSGLLESetAcceptType(TS ts,TSGLLEAcceptType type) 690*3d177a5cSEmil Constantinescu { 691*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 692*3d177a5cSEmil Constantinescu 693*3d177a5cSEmil Constantinescu PetscFunctionBegin; 694*3d177a5cSEmil Constantinescu PetscValidHeaderSpecific(ts,TS_CLASSID,1); 695*3d177a5cSEmil Constantinescu PetscValidCharPointer(type,2); 696*3d177a5cSEmil Constantinescu ierr = PetscTryMethod(ts,"TSGLLESetAcceptType_C",(TS,TSGLLEAcceptType),(ts,type));CHKERRQ(ierr); 697*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 698*3d177a5cSEmil Constantinescu } 699*3d177a5cSEmil Constantinescu 700*3d177a5cSEmil Constantinescu #undef __FUNCT__ 701*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEGetAdapt" 702*3d177a5cSEmil Constantinescu /*@C 703*3d177a5cSEmil Constantinescu TSGLLEGetAdapt - gets the TSGLLEAdapt object from the TS 704*3d177a5cSEmil Constantinescu 705*3d177a5cSEmil Constantinescu Not Collective 706*3d177a5cSEmil Constantinescu 707*3d177a5cSEmil Constantinescu Input Parameter: 708*3d177a5cSEmil Constantinescu . ts - the TS context 709*3d177a5cSEmil Constantinescu 710*3d177a5cSEmil Constantinescu Output Parameter: 711*3d177a5cSEmil Constantinescu . adapt - the TSGLLEAdapt context 712*3d177a5cSEmil Constantinescu 713*3d177a5cSEmil Constantinescu Notes: 714*3d177a5cSEmil Constantinescu This allows the user set options on the TSGLLEAdapt object. Usually it is better to do this using the options 715*3d177a5cSEmil Constantinescu database, so this function is rarely needed. 716*3d177a5cSEmil Constantinescu 717*3d177a5cSEmil Constantinescu Level: advanced 718*3d177a5cSEmil Constantinescu 719*3d177a5cSEmil Constantinescu .seealso: TSGLLEAdapt, TSGLLEAdaptRegister() 720*3d177a5cSEmil Constantinescu @*/ 721*3d177a5cSEmil Constantinescu PetscErrorCode TSGLLEGetAdapt(TS ts,TSGLLEAdapt *adapt) 722*3d177a5cSEmil Constantinescu { 723*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 724*3d177a5cSEmil Constantinescu 725*3d177a5cSEmil Constantinescu PetscFunctionBegin; 726*3d177a5cSEmil Constantinescu PetscValidHeaderSpecific(ts,TS_CLASSID,1); 727*3d177a5cSEmil Constantinescu PetscValidPointer(adapt,2); 728*3d177a5cSEmil Constantinescu ierr = PetscUseMethod(ts,"TSGLLEGetAdapt_C",(TS,TSGLLEAdapt*),(ts,adapt));CHKERRQ(ierr); 729*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 730*3d177a5cSEmil Constantinescu } 731*3d177a5cSEmil Constantinescu 732*3d177a5cSEmil Constantinescu #undef __FUNCT__ 733*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEAccept_Always" 734*3d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEAccept_Always(TS ts,PetscReal tleft,PetscReal h,const PetscReal enorms[],PetscBool *accept) 735*3d177a5cSEmil Constantinescu { 736*3d177a5cSEmil Constantinescu PetscFunctionBegin; 737*3d177a5cSEmil Constantinescu *accept = PETSC_TRUE; 738*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 739*3d177a5cSEmil Constantinescu } 740*3d177a5cSEmil Constantinescu 741*3d177a5cSEmil Constantinescu #undef __FUNCT__ 742*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEUpdateWRMS" 743*3d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEUpdateWRMS(TS ts) 744*3d177a5cSEmil Constantinescu { 745*3d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 746*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 747*3d177a5cSEmil Constantinescu PetscScalar *x,*w; 748*3d177a5cSEmil Constantinescu PetscInt n,i; 749*3d177a5cSEmil Constantinescu 750*3d177a5cSEmil Constantinescu PetscFunctionBegin; 751*3d177a5cSEmil Constantinescu ierr = VecGetArray(gl->X[0],&x);CHKERRQ(ierr); 752*3d177a5cSEmil Constantinescu ierr = VecGetArray(gl->W,&w);CHKERRQ(ierr); 753*3d177a5cSEmil Constantinescu ierr = VecGetLocalSize(gl->W,&n);CHKERRQ(ierr); 754*3d177a5cSEmil Constantinescu for (i=0; i<n; i++) w[i] = 1./(gl->wrms_atol + gl->wrms_rtol*PetscAbsScalar(x[i])); 755*3d177a5cSEmil Constantinescu ierr = VecRestoreArray(gl->X[0],&x);CHKERRQ(ierr); 756*3d177a5cSEmil Constantinescu ierr = VecRestoreArray(gl->W,&w);CHKERRQ(ierr); 757*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 758*3d177a5cSEmil Constantinescu } 759*3d177a5cSEmil Constantinescu 760*3d177a5cSEmil Constantinescu #undef __FUNCT__ 761*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEVecNormWRMS" 762*3d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEVecNormWRMS(TS ts,Vec X,PetscReal *nrm) 763*3d177a5cSEmil Constantinescu { 764*3d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 765*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 766*3d177a5cSEmil Constantinescu PetscScalar *x,*w; 767*3d177a5cSEmil Constantinescu PetscReal sum = 0.0,gsum; 768*3d177a5cSEmil Constantinescu PetscInt n,N,i; 769*3d177a5cSEmil Constantinescu 770*3d177a5cSEmil Constantinescu PetscFunctionBegin; 771*3d177a5cSEmil Constantinescu ierr = VecGetArray(X,&x);CHKERRQ(ierr); 772*3d177a5cSEmil Constantinescu ierr = VecGetArray(gl->W,&w);CHKERRQ(ierr); 773*3d177a5cSEmil Constantinescu ierr = VecGetLocalSize(gl->W,&n);CHKERRQ(ierr); 774*3d177a5cSEmil Constantinescu for (i=0; i<n; i++) sum += PetscAbsScalar(PetscSqr(x[i]*w[i])); 775*3d177a5cSEmil Constantinescu ierr = VecRestoreArray(X,&x);CHKERRQ(ierr); 776*3d177a5cSEmil Constantinescu ierr = VecRestoreArray(gl->W,&w);CHKERRQ(ierr); 777*3d177a5cSEmil Constantinescu ierr = MPIU_Allreduce(&sum,&gsum,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)ts));CHKERRQ(ierr); 778*3d177a5cSEmil Constantinescu ierr = VecGetSize(gl->W,&N);CHKERRQ(ierr); 779*3d177a5cSEmil Constantinescu *nrm = PetscSqrtReal(gsum/(1.*N)); 780*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 781*3d177a5cSEmil Constantinescu } 782*3d177a5cSEmil Constantinescu 783*3d177a5cSEmil Constantinescu #undef __FUNCT__ 784*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLESetType_GLLE" 785*3d177a5cSEmil Constantinescu static PetscErrorCode TSGLLESetType_GLLE(TS ts,TSGLLEType type) 786*3d177a5cSEmil Constantinescu { 787*3d177a5cSEmil Constantinescu PetscErrorCode ierr,(*r)(TS); 788*3d177a5cSEmil Constantinescu PetscBool same; 789*3d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 790*3d177a5cSEmil Constantinescu 791*3d177a5cSEmil Constantinescu PetscFunctionBegin; 792*3d177a5cSEmil Constantinescu if (gl->type_name[0]) { 793*3d177a5cSEmil Constantinescu ierr = PetscStrcmp(gl->type_name,type,&same);CHKERRQ(ierr); 794*3d177a5cSEmil Constantinescu if (same) PetscFunctionReturn(0); 795*3d177a5cSEmil Constantinescu ierr = (*gl->Destroy)(gl);CHKERRQ(ierr); 796*3d177a5cSEmil Constantinescu } 797*3d177a5cSEmil Constantinescu 798*3d177a5cSEmil Constantinescu ierr = PetscFunctionListFind(TSGLLEList,type,&r);CHKERRQ(ierr); 799*3d177a5cSEmil Constantinescu if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSGLLE type \"%s\" given",type); 800*3d177a5cSEmil Constantinescu ierr = (*r)(ts);CHKERRQ(ierr); 801*3d177a5cSEmil Constantinescu ierr = PetscStrcpy(gl->type_name,type);CHKERRQ(ierr); 802*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 803*3d177a5cSEmil Constantinescu } 804*3d177a5cSEmil Constantinescu 805*3d177a5cSEmil Constantinescu #undef __FUNCT__ 806*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLESetAcceptType_GLLE" 807*3d177a5cSEmil Constantinescu static PetscErrorCode TSGLLESetAcceptType_GLLE(TS ts,TSGLLEAcceptType type) 808*3d177a5cSEmil Constantinescu { 809*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 810*3d177a5cSEmil Constantinescu TSGLLEAcceptFunction r; 811*3d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 812*3d177a5cSEmil Constantinescu 813*3d177a5cSEmil Constantinescu PetscFunctionBegin; 814*3d177a5cSEmil Constantinescu ierr = PetscFunctionListFind(TSGLLEAcceptList,type,&r);CHKERRQ(ierr); 815*3d177a5cSEmil Constantinescu if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSGLLEAccept type \"%s\" given",type); 816*3d177a5cSEmil Constantinescu gl->Accept = r; 817*3d177a5cSEmil Constantinescu ierr = PetscStrncpy(gl->accept_name,type,sizeof(gl->accept_name));CHKERRQ(ierr); 818*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 819*3d177a5cSEmil Constantinescu } 820*3d177a5cSEmil Constantinescu 821*3d177a5cSEmil Constantinescu #undef __FUNCT__ 822*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEGetAdapt_GLLE" 823*3d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEGetAdapt_GLLE(TS ts,TSGLLEAdapt *adapt) 824*3d177a5cSEmil Constantinescu { 825*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 826*3d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 827*3d177a5cSEmil Constantinescu 828*3d177a5cSEmil Constantinescu PetscFunctionBegin; 829*3d177a5cSEmil Constantinescu if (!gl->adapt) { 830*3d177a5cSEmil Constantinescu ierr = TSGLLEAdaptCreate(PetscObjectComm((PetscObject)ts),&gl->adapt);CHKERRQ(ierr); 831*3d177a5cSEmil Constantinescu ierr = PetscObjectIncrementTabLevel((PetscObject)gl->adapt,(PetscObject)ts,1);CHKERRQ(ierr); 832*3d177a5cSEmil Constantinescu ierr = PetscLogObjectParent((PetscObject)ts,(PetscObject)gl->adapt);CHKERRQ(ierr); 833*3d177a5cSEmil Constantinescu } 834*3d177a5cSEmil Constantinescu *adapt = gl->adapt; 835*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 836*3d177a5cSEmil Constantinescu } 837*3d177a5cSEmil Constantinescu 838*3d177a5cSEmil Constantinescu #undef __FUNCT__ 839*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEChooseNextScheme" 840*3d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEChooseNextScheme(TS ts,PetscReal h,const PetscReal hmnorm[],PetscInt *next_scheme,PetscReal *next_h,PetscBool *finish) 841*3d177a5cSEmil Constantinescu { 842*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 843*3d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 844*3d177a5cSEmil Constantinescu PetscInt i,n,cur_p,cur,next_sc,candidates[64],orders[64]; 845*3d177a5cSEmil Constantinescu PetscReal errors[64],costs[64],tleft; 846*3d177a5cSEmil Constantinescu 847*3d177a5cSEmil Constantinescu PetscFunctionBegin; 848*3d177a5cSEmil Constantinescu cur = -1; 849*3d177a5cSEmil Constantinescu cur_p = gl->schemes[gl->current_scheme]->p; 850*3d177a5cSEmil Constantinescu tleft = ts->max_time - (ts->ptime + ts->time_step); 851*3d177a5cSEmil Constantinescu for (i=0,n=0; i<gl->nschemes; i++) { 852*3d177a5cSEmil Constantinescu TSGLLEScheme sc = gl->schemes[i]; 853*3d177a5cSEmil Constantinescu if (sc->p < gl->min_order || gl->max_order < sc->p) continue; 854*3d177a5cSEmil Constantinescu if (sc->p == cur_p - 1) errors[n] = PetscAbsScalar(sc->alpha[0])*hmnorm[0]; 855*3d177a5cSEmil Constantinescu else if (sc->p == cur_p) errors[n] = PetscAbsScalar(sc->alpha[0])*hmnorm[1]; 856*3d177a5cSEmil Constantinescu else if (sc->p == cur_p+1) errors[n] = PetscAbsScalar(sc->alpha[0])*(hmnorm[2]+hmnorm[3]); 857*3d177a5cSEmil Constantinescu else continue; 858*3d177a5cSEmil Constantinescu candidates[n] = i; 859*3d177a5cSEmil Constantinescu orders[n] = PetscMin(sc->p,sc->q); /* order of global truncation error */ 860*3d177a5cSEmil Constantinescu costs[n] = sc->s; /* estimate the cost as the number of stages */ 861*3d177a5cSEmil Constantinescu if (i == gl->current_scheme) cur = n; 862*3d177a5cSEmil Constantinescu n++; 863*3d177a5cSEmil Constantinescu } 864*3d177a5cSEmil Constantinescu if (cur < 0 || gl->nschemes <= cur) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Current scheme not found in scheme list"); 865*3d177a5cSEmil Constantinescu ierr = TSGLLEAdaptChoose(gl->adapt,n,orders,errors,costs,cur,h,tleft,&next_sc,next_h,finish);CHKERRQ(ierr); 866*3d177a5cSEmil Constantinescu *next_scheme = candidates[next_sc]; 867*3d177a5cSEmil Constantinescu ierr = PetscInfo7(ts,"Adapt chose scheme %d (%d,%d,%d,%d) with step size %6.2e, finish=%d\n",*next_scheme,gl->schemes[*next_scheme]->p,gl->schemes[*next_scheme]->q,gl->schemes[*next_scheme]->r,gl->schemes[*next_scheme]->s,*next_h,*finish);CHKERRQ(ierr); 868*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 869*3d177a5cSEmil Constantinescu } 870*3d177a5cSEmil Constantinescu 871*3d177a5cSEmil Constantinescu #undef __FUNCT__ 872*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEGetMaxSizes" 873*3d177a5cSEmil Constantinescu static PetscErrorCode TSGLLEGetMaxSizes(TS ts,PetscInt *max_r,PetscInt *max_s) 874*3d177a5cSEmil Constantinescu { 875*3d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 876*3d177a5cSEmil Constantinescu 877*3d177a5cSEmil Constantinescu PetscFunctionBegin; 878*3d177a5cSEmil Constantinescu *max_r = gl->schemes[gl->nschemes-1]->r; 879*3d177a5cSEmil Constantinescu *max_s = gl->schemes[gl->nschemes-1]->s; 880*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 881*3d177a5cSEmil Constantinescu } 882*3d177a5cSEmil Constantinescu 883*3d177a5cSEmil Constantinescu #undef __FUNCT__ 884*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSSolve_GLLE" 885*3d177a5cSEmil Constantinescu static PetscErrorCode TSSolve_GLLE(TS ts) 886*3d177a5cSEmil Constantinescu { 887*3d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 888*3d177a5cSEmil Constantinescu PetscInt i,k,its,lits,max_r,max_s; 889*3d177a5cSEmil Constantinescu PetscBool final_step,finish; 890*3d177a5cSEmil Constantinescu SNESConvergedReason snesreason; 891*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 892*3d177a5cSEmil Constantinescu 893*3d177a5cSEmil Constantinescu PetscFunctionBegin; 894*3d177a5cSEmil Constantinescu ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 895*3d177a5cSEmil Constantinescu 896*3d177a5cSEmil Constantinescu ierr = TSGLLEGetMaxSizes(ts,&max_r,&max_s);CHKERRQ(ierr); 897*3d177a5cSEmil Constantinescu ierr = VecCopy(ts->vec_sol,gl->X[0]);CHKERRQ(ierr); 898*3d177a5cSEmil Constantinescu for (i=1; i<max_r; i++) { 899*3d177a5cSEmil Constantinescu ierr = VecZeroEntries(gl->X[i]);CHKERRQ(ierr); 900*3d177a5cSEmil Constantinescu } 901*3d177a5cSEmil Constantinescu ierr = TSGLLEUpdateWRMS(ts);CHKERRQ(ierr); 902*3d177a5cSEmil Constantinescu 903*3d177a5cSEmil Constantinescu if (0) { 904*3d177a5cSEmil Constantinescu /* Find consistent initial data for DAE */ 905*3d177a5cSEmil Constantinescu gl->stage_time = ts->ptime + ts->time_step; 906*3d177a5cSEmil Constantinescu gl->scoeff = 1.; 907*3d177a5cSEmil Constantinescu gl->stage = 0; 908*3d177a5cSEmil Constantinescu 909*3d177a5cSEmil Constantinescu ierr = VecCopy(ts->vec_sol,gl->Z);CHKERRQ(ierr); 910*3d177a5cSEmil Constantinescu ierr = VecCopy(ts->vec_sol,gl->Y);CHKERRQ(ierr); 911*3d177a5cSEmil Constantinescu ierr = SNESSolve(ts->snes,NULL,gl->Y);CHKERRQ(ierr); 912*3d177a5cSEmil Constantinescu ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr); 913*3d177a5cSEmil Constantinescu ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 914*3d177a5cSEmil Constantinescu ierr = SNESGetConvergedReason(ts->snes,&snesreason);CHKERRQ(ierr); 915*3d177a5cSEmil Constantinescu 916*3d177a5cSEmil Constantinescu ts->snes_its += its; ts->ksp_its += lits; 917*3d177a5cSEmil Constantinescu if (snesreason < 0 && ts->max_snes_failures > 0 && ++ts->num_snes_failures >= ts->max_snes_failures) { 918*3d177a5cSEmil Constantinescu ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 919*3d177a5cSEmil Constantinescu ierr = PetscInfo2(ts,"Step=%D, nonlinear solve solve failures %D greater than current TS allowed, stopping solve\n",ts->steps,ts->num_snes_failures);CHKERRQ(ierr); 920*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 921*3d177a5cSEmil Constantinescu } 922*3d177a5cSEmil Constantinescu } 923*3d177a5cSEmil Constantinescu 924*3d177a5cSEmil Constantinescu if (gl->current_scheme < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"A starting scheme has not been provided"); 925*3d177a5cSEmil Constantinescu 926*3d177a5cSEmil Constantinescu for (k=0,final_step=PETSC_FALSE,finish=PETSC_FALSE; k<ts->max_steps && !finish; k++) { 927*3d177a5cSEmil Constantinescu PetscInt j,r,s,next_scheme = 0,rejections; 928*3d177a5cSEmil Constantinescu PetscReal h,hmnorm[4],enorm[3],next_h; 929*3d177a5cSEmil Constantinescu PetscBool accept; 930*3d177a5cSEmil Constantinescu const PetscScalar *c,*a,*u; 931*3d177a5cSEmil Constantinescu Vec *X,*Ydot,Y; 932*3d177a5cSEmil Constantinescu TSGLLEScheme scheme = gl->schemes[gl->current_scheme]; 933*3d177a5cSEmil Constantinescu 934*3d177a5cSEmil Constantinescu r = scheme->r; s = scheme->s; 935*3d177a5cSEmil Constantinescu c = scheme->c; 936*3d177a5cSEmil Constantinescu a = scheme->a; u = scheme->u; 937*3d177a5cSEmil Constantinescu h = ts->time_step; 938*3d177a5cSEmil Constantinescu X = gl->X; Ydot = gl->Ydot; Y = gl->Y; 939*3d177a5cSEmil Constantinescu 940*3d177a5cSEmil Constantinescu if (ts->ptime > ts->max_time) break; 941*3d177a5cSEmil Constantinescu 942*3d177a5cSEmil Constantinescu /* 943*3d177a5cSEmil Constantinescu We only call PreStep at the start of each STEP, not each STAGE. This is because it is 944*3d177a5cSEmil Constantinescu possible to fail (have to restart a step) after multiple stages. 945*3d177a5cSEmil Constantinescu */ 946*3d177a5cSEmil Constantinescu ierr = TSPreStep(ts);CHKERRQ(ierr); 947*3d177a5cSEmil Constantinescu 948*3d177a5cSEmil Constantinescu rejections = 0; 949*3d177a5cSEmil Constantinescu while (1) { 950*3d177a5cSEmil Constantinescu for (i=0; i<s; i++) { 951*3d177a5cSEmil Constantinescu PetscScalar shift; 952*3d177a5cSEmil Constantinescu gl->scoeff = 1./PetscRealPart(a[i*s+i]); 953*3d177a5cSEmil Constantinescu shift = gl->scoeff/ts->time_step; 954*3d177a5cSEmil Constantinescu gl->stage = i; 955*3d177a5cSEmil Constantinescu gl->stage_time = ts->ptime + PetscRealPart(c[i])*h; 956*3d177a5cSEmil Constantinescu 957*3d177a5cSEmil Constantinescu /* 958*3d177a5cSEmil Constantinescu * Stage equation: Y = h A Y' + U X 959*3d177a5cSEmil Constantinescu * We assume that A is lower-triangular so that we can solve the stages (Y,Y') sequentially 960*3d177a5cSEmil Constantinescu * Build the affine vector z_i = -[1/(h a_ii)](h sum_j a_ij y'_j + sum_j u_ij x_j) 961*3d177a5cSEmil Constantinescu * Then y'_i = z + 1/(h a_ii) y_i 962*3d177a5cSEmil Constantinescu */ 963*3d177a5cSEmil Constantinescu ierr = VecZeroEntries(gl->Z);CHKERRQ(ierr); 964*3d177a5cSEmil Constantinescu for (j=0; j<r; j++) { 965*3d177a5cSEmil Constantinescu ierr = VecAXPY(gl->Z,-shift*u[i*r+j],X[j]);CHKERRQ(ierr); 966*3d177a5cSEmil Constantinescu } 967*3d177a5cSEmil Constantinescu for (j=0; j<i; j++) { 968*3d177a5cSEmil Constantinescu ierr = VecAXPY(gl->Z,-shift*h*a[i*s+j],Ydot[j]);CHKERRQ(ierr); 969*3d177a5cSEmil Constantinescu } 970*3d177a5cSEmil Constantinescu /* Note: Z is used within function evaluation, Ydot = Z + shift*Y */ 971*3d177a5cSEmil Constantinescu 972*3d177a5cSEmil Constantinescu /* Compute an estimate of Y to start Newton iteration */ 973*3d177a5cSEmil Constantinescu if (gl->extrapolate) { 974*3d177a5cSEmil Constantinescu if (i==0) { 975*3d177a5cSEmil Constantinescu /* Linear extrapolation on the first stage */ 976*3d177a5cSEmil Constantinescu ierr = VecWAXPY(Y,c[i]*h,X[1],X[0]);CHKERRQ(ierr); 977*3d177a5cSEmil Constantinescu } else { 978*3d177a5cSEmil Constantinescu /* Linear extrapolation from the last stage */ 979*3d177a5cSEmil Constantinescu ierr = VecAXPY(Y,(c[i]-c[i-1])*h,Ydot[i-1]);CHKERRQ(ierr); 980*3d177a5cSEmil Constantinescu } 981*3d177a5cSEmil Constantinescu } else if (i==0) { /* Directly use solution from the last step, otherwise reuse the last stage (do nothing) */ 982*3d177a5cSEmil Constantinescu ierr = VecCopy(X[0],Y);CHKERRQ(ierr); 983*3d177a5cSEmil Constantinescu } 984*3d177a5cSEmil Constantinescu 985*3d177a5cSEmil Constantinescu /* Solve this stage (Ydot[i] is computed during function evaluation) */ 986*3d177a5cSEmil Constantinescu ierr = SNESSolve(ts->snes,NULL,Y);CHKERRQ(ierr); 987*3d177a5cSEmil Constantinescu ierr = SNESGetIterationNumber(ts->snes,&its);CHKERRQ(ierr); 988*3d177a5cSEmil Constantinescu ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 989*3d177a5cSEmil Constantinescu ierr = SNESGetConvergedReason(ts->snes,&snesreason);CHKERRQ(ierr); 990*3d177a5cSEmil Constantinescu ts->snes_its += its; ts->ksp_its += lits; 991*3d177a5cSEmil Constantinescu if (snesreason < 0 && ts->max_snes_failures > 0 && ++ts->num_snes_failures >= ts->max_snes_failures) { 992*3d177a5cSEmil Constantinescu ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 993*3d177a5cSEmil Constantinescu ierr = PetscInfo2(ts,"Step=%D, nonlinear solve solve failures %D greater than current TS allowed, stopping solve\n",ts->steps,ts->num_snes_failures);CHKERRQ(ierr); 994*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 995*3d177a5cSEmil Constantinescu } 996*3d177a5cSEmil Constantinescu } 997*3d177a5cSEmil Constantinescu 998*3d177a5cSEmil Constantinescu gl->stage_time = ts->ptime + ts->time_step; 999*3d177a5cSEmil Constantinescu 1000*3d177a5cSEmil Constantinescu ierr = (*gl->EstimateHigherMoments)(scheme,h,Ydot,gl->X,gl->himom);CHKERRQ(ierr); 1001*3d177a5cSEmil Constantinescu /* hmnorm[i] = h^{p+i}x^{(p+i)} with i=0,1,2; hmnorm[3] = h^{p+2}(dx'/dx) x^{(p+1)} */ 1002*3d177a5cSEmil Constantinescu for (i=0; i<3; i++) { 1003*3d177a5cSEmil Constantinescu ierr = TSGLLEVecNormWRMS(ts,gl->himom[i],&hmnorm[i+1]);CHKERRQ(ierr); 1004*3d177a5cSEmil Constantinescu } 1005*3d177a5cSEmil Constantinescu enorm[0] = PetscRealPart(scheme->alpha[0])*hmnorm[1]; 1006*3d177a5cSEmil Constantinescu enorm[1] = PetscRealPart(scheme->beta[0]) *hmnorm[2]; 1007*3d177a5cSEmil Constantinescu enorm[2] = PetscRealPart(scheme->gamma[0])*hmnorm[3]; 1008*3d177a5cSEmil Constantinescu ierr = (*gl->Accept)(ts,ts->max_time-gl->stage_time,h,enorm,&accept);CHKERRQ(ierr); 1009*3d177a5cSEmil Constantinescu if (accept) goto accepted; 1010*3d177a5cSEmil Constantinescu rejections++; 1011*3d177a5cSEmil Constantinescu ierr = PetscInfo3(ts,"Step %D (t=%g) not accepted, rejections=%D\n",k,gl->stage_time,rejections);CHKERRQ(ierr); 1012*3d177a5cSEmil Constantinescu if (rejections > gl->max_step_rejections) break; 1013*3d177a5cSEmil Constantinescu /* 1014*3d177a5cSEmil Constantinescu There are lots of reasons why a step might be rejected, including solvers not converging and other factors that 1015*3d177a5cSEmil Constantinescu TSGLLEChooseNextScheme does not support. Additionally, the error estimates may be very screwed up, so I'm not 1016*3d177a5cSEmil Constantinescu convinced that it's safe to just compute a new error estimate using the same interface as the current adaptor 1017*3d177a5cSEmil Constantinescu (the adaptor interface probably has to change). Here we make an arbitrary and naive choice. This assumes that 1018*3d177a5cSEmil Constantinescu steps were written in Nordsieck form. The "correct" method would be to re-complete the previous time step with 1019*3d177a5cSEmil Constantinescu the correct "next" step size. It is unclear to me whether the present ad-hoc method of rescaling X is stable. 1020*3d177a5cSEmil Constantinescu */ 1021*3d177a5cSEmil Constantinescu h *= 0.5; 1022*3d177a5cSEmil Constantinescu for (i=1; i<scheme->r; i++) { 1023*3d177a5cSEmil Constantinescu ierr = VecScale(X[i],PetscPowRealInt(0.5,i));CHKERRQ(ierr); 1024*3d177a5cSEmil Constantinescu } 1025*3d177a5cSEmil Constantinescu } 1026*3d177a5cSEmil Constantinescu SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_CONV_FAILED,"Time step %D (t=%g) not accepted after %D failures\n",k,gl->stage_time,rejections); 1027*3d177a5cSEmil Constantinescu 1028*3d177a5cSEmil Constantinescu accepted: 1029*3d177a5cSEmil Constantinescu /* This term is not error, but it *would* be the leading term for a lower order method */ 1030*3d177a5cSEmil Constantinescu ierr = TSGLLEVecNormWRMS(ts,gl->X[scheme->r-1],&hmnorm[0]);CHKERRQ(ierr); 1031*3d177a5cSEmil Constantinescu /* Correct scaling so that these are equivalent to norms of the Nordsieck vectors */ 1032*3d177a5cSEmil Constantinescu 1033*3d177a5cSEmil Constantinescu ierr = PetscInfo4(ts,"Last moment norm %10.2e, estimated error norms %10.2e %10.2e %10.2e\n",hmnorm[0],enorm[0],enorm[1],enorm[2]);CHKERRQ(ierr); 1034*3d177a5cSEmil Constantinescu if (!final_step) { 1035*3d177a5cSEmil Constantinescu ierr = TSGLLEChooseNextScheme(ts,h,hmnorm,&next_scheme,&next_h,&final_step);CHKERRQ(ierr); 1036*3d177a5cSEmil Constantinescu } else { 1037*3d177a5cSEmil Constantinescu /* Dummy values to complete the current step in a consistent manner */ 1038*3d177a5cSEmil Constantinescu next_scheme = gl->current_scheme; 1039*3d177a5cSEmil Constantinescu next_h = h; 1040*3d177a5cSEmil Constantinescu finish = PETSC_TRUE; 1041*3d177a5cSEmil Constantinescu } 1042*3d177a5cSEmil Constantinescu 1043*3d177a5cSEmil Constantinescu X = gl->Xold; 1044*3d177a5cSEmil Constantinescu gl->Xold = gl->X; 1045*3d177a5cSEmil Constantinescu gl->X = X; 1046*3d177a5cSEmil Constantinescu ierr = (*gl->CompleteStep)(scheme,h,gl->schemes[next_scheme],next_h,Ydot,gl->Xold,gl->X);CHKERRQ(ierr); 1047*3d177a5cSEmil Constantinescu 1048*3d177a5cSEmil Constantinescu ierr = TSGLLEUpdateWRMS(ts);CHKERRQ(ierr); 1049*3d177a5cSEmil Constantinescu 1050*3d177a5cSEmil Constantinescu /* Post the solution for the user, we could avoid this copy with a small bit of cleverness */ 1051*3d177a5cSEmil Constantinescu ierr = VecCopy(gl->X[0],ts->vec_sol);CHKERRQ(ierr); 1052*3d177a5cSEmil Constantinescu ts->ptime += h; 1053*3d177a5cSEmil Constantinescu ts->steps++; 1054*3d177a5cSEmil Constantinescu ts->total_steps++; 1055*3d177a5cSEmil Constantinescu 1056*3d177a5cSEmil Constantinescu ierr = TSPostStep(ts);CHKERRQ(ierr); 1057*3d177a5cSEmil Constantinescu ierr = TSMonitor(ts,ts->steps,ts->ptime,ts->vec_sol);CHKERRQ(ierr); 1058*3d177a5cSEmil Constantinescu 1059*3d177a5cSEmil Constantinescu gl->current_scheme = next_scheme; 1060*3d177a5cSEmil Constantinescu ts->time_step = next_h; 1061*3d177a5cSEmil Constantinescu } 1062*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 1063*3d177a5cSEmil Constantinescu } 1064*3d177a5cSEmil Constantinescu 1065*3d177a5cSEmil Constantinescu /*------------------------------------------------------------*/ 1066*3d177a5cSEmil Constantinescu 1067*3d177a5cSEmil Constantinescu #undef __FUNCT__ 1068*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSReset_GLLE" 1069*3d177a5cSEmil Constantinescu static PetscErrorCode TSReset_GLLE(TS ts) 1070*3d177a5cSEmil Constantinescu { 1071*3d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 1072*3d177a5cSEmil Constantinescu PetscInt max_r,max_s; 1073*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 1074*3d177a5cSEmil Constantinescu 1075*3d177a5cSEmil Constantinescu PetscFunctionBegin; 1076*3d177a5cSEmil Constantinescu if (gl->setupcalled) { 1077*3d177a5cSEmil Constantinescu ierr = TSGLLEGetMaxSizes(ts,&max_r,&max_s);CHKERRQ(ierr); 1078*3d177a5cSEmil Constantinescu ierr = VecDestroyVecs(max_r,&gl->Xold);CHKERRQ(ierr); 1079*3d177a5cSEmil Constantinescu ierr = VecDestroyVecs(max_r,&gl->X);CHKERRQ(ierr); 1080*3d177a5cSEmil Constantinescu ierr = VecDestroyVecs(max_s,&gl->Ydot);CHKERRQ(ierr); 1081*3d177a5cSEmil Constantinescu ierr = VecDestroyVecs(3,&gl->himom);CHKERRQ(ierr); 1082*3d177a5cSEmil Constantinescu ierr = VecDestroy(&gl->W);CHKERRQ(ierr); 1083*3d177a5cSEmil Constantinescu ierr = VecDestroy(&gl->Y);CHKERRQ(ierr); 1084*3d177a5cSEmil Constantinescu ierr = VecDestroy(&gl->Z);CHKERRQ(ierr); 1085*3d177a5cSEmil Constantinescu } 1086*3d177a5cSEmil Constantinescu gl->setupcalled = PETSC_FALSE; 1087*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 1088*3d177a5cSEmil Constantinescu } 1089*3d177a5cSEmil Constantinescu 1090*3d177a5cSEmil Constantinescu #undef __FUNCT__ 1091*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSDestroy_GLLE" 1092*3d177a5cSEmil Constantinescu static PetscErrorCode TSDestroy_GLLE(TS ts) 1093*3d177a5cSEmil Constantinescu { 1094*3d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 1095*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 1096*3d177a5cSEmil Constantinescu 1097*3d177a5cSEmil Constantinescu PetscFunctionBegin; 1098*3d177a5cSEmil Constantinescu ierr = TSReset_GLLE(ts);CHKERRQ(ierr); 1099*3d177a5cSEmil Constantinescu if (gl->adapt) {ierr = TSGLLEAdaptDestroy(&gl->adapt);CHKERRQ(ierr);} 1100*3d177a5cSEmil Constantinescu if (gl->Destroy) {ierr = (*gl->Destroy)(gl);CHKERRQ(ierr);} 1101*3d177a5cSEmil Constantinescu ierr = PetscFree(ts->data);CHKERRQ(ierr); 1102*3d177a5cSEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSGLLESetType_C",NULL);CHKERRQ(ierr); 1103*3d177a5cSEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSGLLESetAcceptType_C",NULL);CHKERRQ(ierr); 1104*3d177a5cSEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSGLLEGetAdapt_C",NULL);CHKERRQ(ierr); 1105*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 1106*3d177a5cSEmil Constantinescu } 1107*3d177a5cSEmil Constantinescu 1108*3d177a5cSEmil Constantinescu /* 1109*3d177a5cSEmil Constantinescu This defines the nonlinear equation that is to be solved with SNES 1110*3d177a5cSEmil Constantinescu g(x) = f(t,x,z+shift*x) = 0 1111*3d177a5cSEmil Constantinescu */ 1112*3d177a5cSEmil Constantinescu #undef __FUNCT__ 1113*3d177a5cSEmil Constantinescu #define __FUNCT__ "SNESTSFormFunction_GLLE" 1114*3d177a5cSEmil Constantinescu static PetscErrorCode SNESTSFormFunction_GLLE(SNES snes,Vec x,Vec f,TS ts) 1115*3d177a5cSEmil Constantinescu { 1116*3d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 1117*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 1118*3d177a5cSEmil Constantinescu Vec Z,Ydot; 1119*3d177a5cSEmil Constantinescu DM dm,dmsave; 1120*3d177a5cSEmil Constantinescu 1121*3d177a5cSEmil Constantinescu PetscFunctionBegin; 1122*3d177a5cSEmil Constantinescu ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 1123*3d177a5cSEmil Constantinescu ierr = TSGLLEGetVecs(ts,dm,&Z,&Ydot);CHKERRQ(ierr); 1124*3d177a5cSEmil Constantinescu ierr = VecWAXPY(Ydot,gl->scoeff/ts->time_step,x,Z);CHKERRQ(ierr); 1125*3d177a5cSEmil Constantinescu dmsave = ts->dm; 1126*3d177a5cSEmil Constantinescu ts->dm = dm; 1127*3d177a5cSEmil Constantinescu ierr = TSComputeIFunction(ts,gl->stage_time,x,Ydot,f,PETSC_FALSE);CHKERRQ(ierr); 1128*3d177a5cSEmil Constantinescu ts->dm = dmsave; 1129*3d177a5cSEmil Constantinescu ierr = TSGLLERestoreVecs(ts,dm,&Z,&Ydot);CHKERRQ(ierr); 1130*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 1131*3d177a5cSEmil Constantinescu } 1132*3d177a5cSEmil Constantinescu 1133*3d177a5cSEmil Constantinescu #undef __FUNCT__ 1134*3d177a5cSEmil Constantinescu #define __FUNCT__ "SNESTSFormJacobian_GLLE" 1135*3d177a5cSEmil Constantinescu static PetscErrorCode SNESTSFormJacobian_GLLE(SNES snes,Vec x,Mat A,Mat B,TS ts) 1136*3d177a5cSEmil Constantinescu { 1137*3d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 1138*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 1139*3d177a5cSEmil Constantinescu Vec Z,Ydot; 1140*3d177a5cSEmil Constantinescu DM dm,dmsave; 1141*3d177a5cSEmil Constantinescu 1142*3d177a5cSEmil Constantinescu PetscFunctionBegin; 1143*3d177a5cSEmil Constantinescu ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 1144*3d177a5cSEmil Constantinescu ierr = TSGLLEGetVecs(ts,dm,&Z,&Ydot);CHKERRQ(ierr); 1145*3d177a5cSEmil Constantinescu dmsave = ts->dm; 1146*3d177a5cSEmil Constantinescu ts->dm = dm; 1147*3d177a5cSEmil Constantinescu /* gl->Xdot will have already been computed in SNESTSFormFunction_GLLE */ 1148*3d177a5cSEmil Constantinescu ierr = TSComputeIJacobian(ts,gl->stage_time,x,gl->Ydot[gl->stage],gl->scoeff/ts->time_step,A,B,PETSC_FALSE);CHKERRQ(ierr); 1149*3d177a5cSEmil Constantinescu ts->dm = dmsave; 1150*3d177a5cSEmil Constantinescu ierr = TSGLLERestoreVecs(ts,dm,&Z,&Ydot);CHKERRQ(ierr); 1151*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 1152*3d177a5cSEmil Constantinescu } 1153*3d177a5cSEmil Constantinescu 1154*3d177a5cSEmil Constantinescu 1155*3d177a5cSEmil Constantinescu #undef __FUNCT__ 1156*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSSetUp_GLLE" 1157*3d177a5cSEmil Constantinescu static PetscErrorCode TSSetUp_GLLE(TS ts) 1158*3d177a5cSEmil Constantinescu { 1159*3d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 1160*3d177a5cSEmil Constantinescu PetscInt max_r,max_s; 1161*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 1162*3d177a5cSEmil Constantinescu DM dm; 1163*3d177a5cSEmil Constantinescu 1164*3d177a5cSEmil Constantinescu PetscFunctionBegin; 1165*3d177a5cSEmil Constantinescu gl->setupcalled = PETSC_TRUE; 1166*3d177a5cSEmil Constantinescu ierr = TSGLLEGetMaxSizes(ts,&max_r,&max_s);CHKERRQ(ierr); 1167*3d177a5cSEmil Constantinescu ierr = VecDuplicateVecs(ts->vec_sol,max_r,&gl->X);CHKERRQ(ierr); 1168*3d177a5cSEmil Constantinescu ierr = VecDuplicateVecs(ts->vec_sol,max_r,&gl->Xold);CHKERRQ(ierr); 1169*3d177a5cSEmil Constantinescu ierr = VecDuplicateVecs(ts->vec_sol,max_s,&gl->Ydot);CHKERRQ(ierr); 1170*3d177a5cSEmil Constantinescu ierr = VecDuplicateVecs(ts->vec_sol,3,&gl->himom);CHKERRQ(ierr); 1171*3d177a5cSEmil Constantinescu ierr = VecDuplicate(ts->vec_sol,&gl->W);CHKERRQ(ierr); 1172*3d177a5cSEmil Constantinescu ierr = VecDuplicate(ts->vec_sol,&gl->Y);CHKERRQ(ierr); 1173*3d177a5cSEmil Constantinescu ierr = VecDuplicate(ts->vec_sol,&gl->Z);CHKERRQ(ierr); 1174*3d177a5cSEmil Constantinescu 1175*3d177a5cSEmil Constantinescu /* Default acceptance tests and adaptivity */ 1176*3d177a5cSEmil Constantinescu if (!gl->Accept) {ierr = TSGLLESetAcceptType(ts,TSGLLEACCEPT_ALWAYS);CHKERRQ(ierr);} 1177*3d177a5cSEmil Constantinescu if (!gl->adapt) {ierr = TSGLLEGetAdapt(ts,&gl->adapt);CHKERRQ(ierr);} 1178*3d177a5cSEmil Constantinescu 1179*3d177a5cSEmil Constantinescu if (gl->current_scheme < 0) { 1180*3d177a5cSEmil Constantinescu PetscInt i; 1181*3d177a5cSEmil Constantinescu for (i=0;; i++) { 1182*3d177a5cSEmil Constantinescu if (gl->schemes[i]->p == gl->start_order) break; 1183*3d177a5cSEmil Constantinescu if (i+1 == gl->nschemes) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"No schemes available with requested start order %d",i); 1184*3d177a5cSEmil Constantinescu } 1185*3d177a5cSEmil Constantinescu gl->current_scheme = i; 1186*3d177a5cSEmil Constantinescu } 1187*3d177a5cSEmil Constantinescu ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 1188*3d177a5cSEmil Constantinescu if (dm) { 1189*3d177a5cSEmil Constantinescu ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSGLLE,DMRestrictHook_TSGLLE,ts);CHKERRQ(ierr); 1190*3d177a5cSEmil Constantinescu ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSGLLE,DMSubDomainRestrictHook_TSGLLE,ts);CHKERRQ(ierr); 1191*3d177a5cSEmil Constantinescu } 1192*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 1193*3d177a5cSEmil Constantinescu } 1194*3d177a5cSEmil Constantinescu /*------------------------------------------------------------*/ 1195*3d177a5cSEmil Constantinescu 1196*3d177a5cSEmil Constantinescu #undef __FUNCT__ 1197*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSSetFromOptions_GLLE" 1198*3d177a5cSEmil Constantinescu static PetscErrorCode TSSetFromOptions_GLLE(PetscOptionItems *PetscOptionsObject,TS ts) 1199*3d177a5cSEmil Constantinescu { 1200*3d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 1201*3d177a5cSEmil Constantinescu char tname[256] = TSGLLE_IRKS,completef[256] = "rescale-and-modify"; 1202*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 1203*3d177a5cSEmil Constantinescu 1204*3d177a5cSEmil Constantinescu PetscFunctionBegin; 1205*3d177a5cSEmil Constantinescu ierr = PetscOptionsHead(PetscOptionsObject,"General Linear ODE solver options");CHKERRQ(ierr); 1206*3d177a5cSEmil Constantinescu { 1207*3d177a5cSEmil Constantinescu PetscBool flg; 1208*3d177a5cSEmil Constantinescu ierr = PetscOptionsFList("-ts_gl_type","Type of GL method","TSGLLESetType",TSGLLEList,gl->type_name[0] ? gl->type_name : tname,tname,sizeof(tname),&flg);CHKERRQ(ierr); 1209*3d177a5cSEmil Constantinescu if (flg || !gl->type_name[0]) { 1210*3d177a5cSEmil Constantinescu ierr = TSGLLESetType(ts,tname);CHKERRQ(ierr); 1211*3d177a5cSEmil Constantinescu } 1212*3d177a5cSEmil Constantinescu ierr = PetscOptionsInt("-ts_gl_max_step_rejections","Maximum number of times to attempt a step","None",gl->max_step_rejections,&gl->max_step_rejections,NULL);CHKERRQ(ierr); 1213*3d177a5cSEmil Constantinescu ierr = PetscOptionsInt("-ts_gl_max_order","Maximum order to try","TSGLLESetMaxOrder",gl->max_order,&gl->max_order,NULL);CHKERRQ(ierr); 1214*3d177a5cSEmil Constantinescu ierr = PetscOptionsInt("-ts_gl_min_order","Minimum order to try","TSGLLESetMinOrder",gl->min_order,&gl->min_order,NULL);CHKERRQ(ierr); 1215*3d177a5cSEmil Constantinescu ierr = PetscOptionsInt("-ts_gl_start_order","Initial order to try","TSGLLESetMinOrder",gl->start_order,&gl->start_order,NULL);CHKERRQ(ierr); 1216*3d177a5cSEmil Constantinescu ierr = PetscOptionsEnum("-ts_gl_error_direction","Which direction to look when estimating error","TSGLLESetErrorDirection",TSGLLEErrorDirections,(PetscEnum)gl->error_direction,(PetscEnum*)&gl->error_direction,NULL);CHKERRQ(ierr); 1217*3d177a5cSEmil Constantinescu ierr = PetscOptionsBool("-ts_gl_extrapolate","Extrapolate stage solution from previous solution (sometimes unstable)","TSGLLESetExtrapolate",gl->extrapolate,&gl->extrapolate,NULL);CHKERRQ(ierr); 1218*3d177a5cSEmil Constantinescu ierr = PetscOptionsReal("-ts_gl_atol","Absolute tolerance","TSGLLESetTolerances",gl->wrms_atol,&gl->wrms_atol,NULL);CHKERRQ(ierr); 1219*3d177a5cSEmil Constantinescu ierr = PetscOptionsReal("-ts_gl_rtol","Relative tolerance","TSGLLESetTolerances",gl->wrms_rtol,&gl->wrms_rtol,NULL);CHKERRQ(ierr); 1220*3d177a5cSEmil Constantinescu ierr = PetscOptionsString("-ts_gl_complete","Method to use for completing the step","none",completef,completef,sizeof(completef),&flg);CHKERRQ(ierr); 1221*3d177a5cSEmil Constantinescu if (flg) { 1222*3d177a5cSEmil Constantinescu PetscBool match1,match2; 1223*3d177a5cSEmil Constantinescu ierr = PetscStrcmp(completef,"rescale",&match1);CHKERRQ(ierr); 1224*3d177a5cSEmil Constantinescu ierr = PetscStrcmp(completef,"rescale-and-modify",&match2);CHKERRQ(ierr); 1225*3d177a5cSEmil Constantinescu if (match1) gl->CompleteStep = TSGLLECompleteStep_Rescale; 1226*3d177a5cSEmil Constantinescu else if (match2) gl->CompleteStep = TSGLLECompleteStep_RescaleAndModify; 1227*3d177a5cSEmil Constantinescu else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"%s",completef); 1228*3d177a5cSEmil Constantinescu } 1229*3d177a5cSEmil Constantinescu { 1230*3d177a5cSEmil Constantinescu char type[256] = TSGLLEACCEPT_ALWAYS; 1231*3d177a5cSEmil Constantinescu ierr = PetscOptionsFList("-ts_gl_accept_type","Method to use for determining whether to accept a step","TSGLLESetAcceptType",TSGLLEAcceptList,gl->accept_name[0] ? gl->accept_name : type,type,sizeof(type),&flg);CHKERRQ(ierr); 1232*3d177a5cSEmil Constantinescu if (flg || !gl->accept_name[0]) { 1233*3d177a5cSEmil Constantinescu ierr = TSGLLESetAcceptType(ts,type);CHKERRQ(ierr); 1234*3d177a5cSEmil Constantinescu } 1235*3d177a5cSEmil Constantinescu } 1236*3d177a5cSEmil Constantinescu { 1237*3d177a5cSEmil Constantinescu TSGLLEAdapt adapt; 1238*3d177a5cSEmil Constantinescu ierr = TSGLLEGetAdapt(ts,&adapt);CHKERRQ(ierr); 1239*3d177a5cSEmil Constantinescu ierr = TSGLLEAdaptSetFromOptions(PetscOptionsObject,adapt);CHKERRQ(ierr); 1240*3d177a5cSEmil Constantinescu } 1241*3d177a5cSEmil Constantinescu } 1242*3d177a5cSEmil Constantinescu ierr = PetscOptionsTail();CHKERRQ(ierr); 1243*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 1244*3d177a5cSEmil Constantinescu } 1245*3d177a5cSEmil Constantinescu 1246*3d177a5cSEmil Constantinescu #undef __FUNCT__ 1247*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSView_GLLE" 1248*3d177a5cSEmil Constantinescu static PetscErrorCode TSView_GLLE(TS ts,PetscViewer viewer) 1249*3d177a5cSEmil Constantinescu { 1250*3d177a5cSEmil Constantinescu TS_GLLE *gl = (TS_GLLE*)ts->data; 1251*3d177a5cSEmil Constantinescu PetscInt i; 1252*3d177a5cSEmil Constantinescu PetscBool iascii,details; 1253*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 1254*3d177a5cSEmil Constantinescu 1255*3d177a5cSEmil Constantinescu PetscFunctionBegin; 1256*3d177a5cSEmil Constantinescu ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1257*3d177a5cSEmil Constantinescu if (iascii) { 1258*3d177a5cSEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer," min order %D, max order %D, current order %D\n",gl->min_order,gl->max_order,gl->schemes[gl->current_scheme]->p);CHKERRQ(ierr); 1259*3d177a5cSEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer," Error estimation: %s\n",TSGLLEErrorDirections[gl->error_direction]);CHKERRQ(ierr); 1260*3d177a5cSEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer," Extrapolation: %s\n",gl->extrapolate ? "yes" : "no");CHKERRQ(ierr); 1261*3d177a5cSEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer," Acceptance test: %s\n",gl->accept_name[0] ? gl->accept_name : "(not yet set)");CHKERRQ(ierr); 1262*3d177a5cSEmil Constantinescu ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1263*3d177a5cSEmil Constantinescu ierr = TSGLLEAdaptView(gl->adapt,viewer);CHKERRQ(ierr); 1264*3d177a5cSEmil Constantinescu ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1265*3d177a5cSEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer," type: %s\n",gl->type_name[0] ? gl->type_name : "(not yet set)");CHKERRQ(ierr); 1266*3d177a5cSEmil Constantinescu ierr = PetscViewerASCIIPrintf(viewer,"Schemes within family (%d):\n",gl->nschemes);CHKERRQ(ierr); 1267*3d177a5cSEmil Constantinescu details = PETSC_FALSE; 1268*3d177a5cSEmil Constantinescu ierr = PetscOptionsGetBool(((PetscObject)ts)->options,((PetscObject)ts)->prefix,"-ts_gl_view_detailed",&details,NULL);CHKERRQ(ierr); 1269*3d177a5cSEmil Constantinescu ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1270*3d177a5cSEmil Constantinescu for (i=0; i<gl->nschemes; i++) { 1271*3d177a5cSEmil Constantinescu ierr = TSGLLESchemeView(gl->schemes[i],details,viewer);CHKERRQ(ierr); 1272*3d177a5cSEmil Constantinescu } 1273*3d177a5cSEmil Constantinescu if (gl->View) { 1274*3d177a5cSEmil Constantinescu ierr = (*gl->View)(gl,viewer);CHKERRQ(ierr); 1275*3d177a5cSEmil Constantinescu } 1276*3d177a5cSEmil Constantinescu ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 1277*3d177a5cSEmil Constantinescu } 1278*3d177a5cSEmil Constantinescu ierr = SNESView(ts->snes,viewer);CHKERRQ(ierr); 1279*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 1280*3d177a5cSEmil Constantinescu } 1281*3d177a5cSEmil Constantinescu 1282*3d177a5cSEmil Constantinescu #undef __FUNCT__ 1283*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLERegister" 1284*3d177a5cSEmil Constantinescu /*@C 1285*3d177a5cSEmil Constantinescu TSGLLERegister - adds a TSGLLE implementation 1286*3d177a5cSEmil Constantinescu 1287*3d177a5cSEmil Constantinescu Not Collective 1288*3d177a5cSEmil Constantinescu 1289*3d177a5cSEmil Constantinescu Input Parameters: 1290*3d177a5cSEmil Constantinescu + name_scheme - name of user-defined general linear scheme 1291*3d177a5cSEmil Constantinescu - routine_create - routine to create method context 1292*3d177a5cSEmil Constantinescu 1293*3d177a5cSEmil Constantinescu Notes: 1294*3d177a5cSEmil Constantinescu TSGLLERegister() may be called multiple times to add several user-defined families. 1295*3d177a5cSEmil Constantinescu 1296*3d177a5cSEmil Constantinescu Sample usage: 1297*3d177a5cSEmil Constantinescu .vb 1298*3d177a5cSEmil Constantinescu TSGLLERegister("my_scheme",MySchemeCreate); 1299*3d177a5cSEmil Constantinescu .ve 1300*3d177a5cSEmil Constantinescu 1301*3d177a5cSEmil Constantinescu Then, your scheme can be chosen with the procedural interface via 1302*3d177a5cSEmil Constantinescu $ TSGLLESetType(ts,"my_scheme") 1303*3d177a5cSEmil Constantinescu or at runtime via the option 1304*3d177a5cSEmil Constantinescu $ -ts_gl_type my_scheme 1305*3d177a5cSEmil Constantinescu 1306*3d177a5cSEmil Constantinescu Level: advanced 1307*3d177a5cSEmil Constantinescu 1308*3d177a5cSEmil Constantinescu .keywords: TSGLLE, register 1309*3d177a5cSEmil Constantinescu 1310*3d177a5cSEmil Constantinescu .seealso: TSGLLERegisterAll() 1311*3d177a5cSEmil Constantinescu @*/ 1312*3d177a5cSEmil Constantinescu PetscErrorCode TSGLLERegister(const char sname[],PetscErrorCode (*function)(TS)) 1313*3d177a5cSEmil Constantinescu { 1314*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 1315*3d177a5cSEmil Constantinescu 1316*3d177a5cSEmil Constantinescu PetscFunctionBegin; 1317*3d177a5cSEmil Constantinescu ierr = PetscFunctionListAdd(&TSGLLEList,sname,function);CHKERRQ(ierr); 1318*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 1319*3d177a5cSEmil Constantinescu } 1320*3d177a5cSEmil Constantinescu 1321*3d177a5cSEmil Constantinescu #undef __FUNCT__ 1322*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEAcceptRegister" 1323*3d177a5cSEmil Constantinescu /*@C 1324*3d177a5cSEmil Constantinescu TSGLLEAcceptRegister - adds a TSGLLE acceptance scheme 1325*3d177a5cSEmil Constantinescu 1326*3d177a5cSEmil Constantinescu Not Collective 1327*3d177a5cSEmil Constantinescu 1328*3d177a5cSEmil Constantinescu Input Parameters: 1329*3d177a5cSEmil Constantinescu + name_scheme - name of user-defined acceptance scheme 1330*3d177a5cSEmil Constantinescu - routine_create - routine to create method context 1331*3d177a5cSEmil Constantinescu 1332*3d177a5cSEmil Constantinescu Notes: 1333*3d177a5cSEmil Constantinescu TSGLLEAcceptRegister() may be called multiple times to add several user-defined families. 1334*3d177a5cSEmil Constantinescu 1335*3d177a5cSEmil Constantinescu Sample usage: 1336*3d177a5cSEmil Constantinescu .vb 1337*3d177a5cSEmil Constantinescu TSGLLEAcceptRegister("my_scheme",MySchemeCreate); 1338*3d177a5cSEmil Constantinescu .ve 1339*3d177a5cSEmil Constantinescu 1340*3d177a5cSEmil Constantinescu Then, your scheme can be chosen with the procedural interface via 1341*3d177a5cSEmil Constantinescu $ TSGLLESetAcceptType(ts,"my_scheme") 1342*3d177a5cSEmil Constantinescu or at runtime via the option 1343*3d177a5cSEmil Constantinescu $ -ts_gl_accept_type my_scheme 1344*3d177a5cSEmil Constantinescu 1345*3d177a5cSEmil Constantinescu Level: advanced 1346*3d177a5cSEmil Constantinescu 1347*3d177a5cSEmil Constantinescu .keywords: TSGLLE, TSGLLEAcceptType, register 1348*3d177a5cSEmil Constantinescu 1349*3d177a5cSEmil Constantinescu .seealso: TSGLLERegisterAll() 1350*3d177a5cSEmil Constantinescu @*/ 1351*3d177a5cSEmil Constantinescu PetscErrorCode TSGLLEAcceptRegister(const char sname[],TSGLLEAcceptFunction function) 1352*3d177a5cSEmil Constantinescu { 1353*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 1354*3d177a5cSEmil Constantinescu 1355*3d177a5cSEmil Constantinescu PetscFunctionBegin; 1356*3d177a5cSEmil Constantinescu ierr = PetscFunctionListAdd(&TSGLLEAcceptList,sname,function);CHKERRQ(ierr); 1357*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 1358*3d177a5cSEmil Constantinescu } 1359*3d177a5cSEmil Constantinescu 1360*3d177a5cSEmil Constantinescu #undef __FUNCT__ 1361*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLERegisterAll" 1362*3d177a5cSEmil Constantinescu /*@C 1363*3d177a5cSEmil Constantinescu TSGLLERegisterAll - Registers all of the general linear methods in TSGLLE 1364*3d177a5cSEmil Constantinescu 1365*3d177a5cSEmil Constantinescu Not Collective 1366*3d177a5cSEmil Constantinescu 1367*3d177a5cSEmil Constantinescu Level: advanced 1368*3d177a5cSEmil Constantinescu 1369*3d177a5cSEmil Constantinescu .keywords: TS, TSGLLE, register, all 1370*3d177a5cSEmil Constantinescu 1371*3d177a5cSEmil Constantinescu .seealso: TSGLLERegisterDestroy() 1372*3d177a5cSEmil Constantinescu @*/ 1373*3d177a5cSEmil Constantinescu PetscErrorCode TSGLLERegisterAll(void) 1374*3d177a5cSEmil Constantinescu { 1375*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 1376*3d177a5cSEmil Constantinescu 1377*3d177a5cSEmil Constantinescu PetscFunctionBegin; 1378*3d177a5cSEmil Constantinescu if (TSGLLERegisterAllCalled) PetscFunctionReturn(0); 1379*3d177a5cSEmil Constantinescu TSGLLERegisterAllCalled = PETSC_TRUE; 1380*3d177a5cSEmil Constantinescu 1381*3d177a5cSEmil Constantinescu ierr = TSGLLERegister(TSGLLE_IRKS, TSGLLECreate_IRKS);CHKERRQ(ierr); 1382*3d177a5cSEmil Constantinescu ierr = TSGLLEAcceptRegister(TSGLLEACCEPT_ALWAYS,TSGLLEAccept_Always);CHKERRQ(ierr); 1383*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 1384*3d177a5cSEmil Constantinescu } 1385*3d177a5cSEmil Constantinescu 1386*3d177a5cSEmil Constantinescu #undef __FUNCT__ 1387*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEInitializePackage" 1388*3d177a5cSEmil Constantinescu /*@C 1389*3d177a5cSEmil Constantinescu TSGLLEInitializePackage - This function initializes everything in the TSGLLE package. It is called 1390*3d177a5cSEmil Constantinescu from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_GLLE() 1391*3d177a5cSEmil Constantinescu when using static libraries. 1392*3d177a5cSEmil Constantinescu 1393*3d177a5cSEmil Constantinescu Level: developer 1394*3d177a5cSEmil Constantinescu 1395*3d177a5cSEmil Constantinescu .keywords: TS, TSGLLE, initialize, package 1396*3d177a5cSEmil Constantinescu .seealso: PetscInitialize() 1397*3d177a5cSEmil Constantinescu @*/ 1398*3d177a5cSEmil Constantinescu PetscErrorCode TSGLLEInitializePackage(void) 1399*3d177a5cSEmil Constantinescu { 1400*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 1401*3d177a5cSEmil Constantinescu 1402*3d177a5cSEmil Constantinescu PetscFunctionBegin; 1403*3d177a5cSEmil Constantinescu if (TSGLLEPackageInitialized) PetscFunctionReturn(0); 1404*3d177a5cSEmil Constantinescu TSGLLEPackageInitialized = PETSC_TRUE; 1405*3d177a5cSEmil Constantinescu ierr = TSGLLERegisterAll();CHKERRQ(ierr); 1406*3d177a5cSEmil Constantinescu ierr = PetscRegisterFinalize(TSGLLEFinalizePackage);CHKERRQ(ierr); 1407*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 1408*3d177a5cSEmil Constantinescu } 1409*3d177a5cSEmil Constantinescu 1410*3d177a5cSEmil Constantinescu #undef __FUNCT__ 1411*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSGLLEFinalizePackage" 1412*3d177a5cSEmil Constantinescu /*@C 1413*3d177a5cSEmil Constantinescu TSGLLEFinalizePackage - This function destroys everything in the TSGLLE package. It is 1414*3d177a5cSEmil Constantinescu called from PetscFinalize(). 1415*3d177a5cSEmil Constantinescu 1416*3d177a5cSEmil Constantinescu Level: developer 1417*3d177a5cSEmil Constantinescu 1418*3d177a5cSEmil Constantinescu .keywords: Petsc, destroy, package 1419*3d177a5cSEmil Constantinescu .seealso: PetscFinalize() 1420*3d177a5cSEmil Constantinescu @*/ 1421*3d177a5cSEmil Constantinescu PetscErrorCode TSGLLEFinalizePackage(void) 1422*3d177a5cSEmil Constantinescu { 1423*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 1424*3d177a5cSEmil Constantinescu 1425*3d177a5cSEmil Constantinescu PetscFunctionBegin; 1426*3d177a5cSEmil Constantinescu ierr = PetscFunctionListDestroy(&TSGLLEList);CHKERRQ(ierr); 1427*3d177a5cSEmil Constantinescu ierr = PetscFunctionListDestroy(&TSGLLEAcceptList);CHKERRQ(ierr); 1428*3d177a5cSEmil Constantinescu TSGLLEPackageInitialized = PETSC_FALSE; 1429*3d177a5cSEmil Constantinescu TSGLLERegisterAllCalled = PETSC_FALSE; 1430*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 1431*3d177a5cSEmil Constantinescu } 1432*3d177a5cSEmil Constantinescu 1433*3d177a5cSEmil Constantinescu /* ------------------------------------------------------------ */ 1434*3d177a5cSEmil Constantinescu /*MC 1435*3d177a5cSEmil Constantinescu TSGLLE - DAE solver using implicit General Linear methods 1436*3d177a5cSEmil Constantinescu 1437*3d177a5cSEmil Constantinescu These methods contain Runge-Kutta and multistep schemes as special cases. These special cases have some fundamental 1438*3d177a5cSEmil Constantinescu limitations. For example, diagonally implicit Runge-Kutta cannot have stage order greater than 1 which limits their 1439*3d177a5cSEmil Constantinescu applicability to very stiff systems. Meanwhile, multistep methods cannot be A-stable for order greater than 2 and BDF 1440*3d177a5cSEmil Constantinescu are not 0-stable for order greater than 6. GL methods can be A- and L-stable with arbitrarily high stage order and 1441*3d177a5cSEmil Constantinescu reliable error estimates for both 1 and 2 orders higher to facilitate adaptive step sizes and adaptive order schemes. 1442*3d177a5cSEmil Constantinescu All this is possible while preserving a singly diagonally implicit structure. 1443*3d177a5cSEmil Constantinescu 1444*3d177a5cSEmil Constantinescu Options database keys: 1445*3d177a5cSEmil Constantinescu + -ts_gl_type <type> - the class of general linear method (irks) 1446*3d177a5cSEmil Constantinescu . -ts_gl_rtol <tol> - relative error 1447*3d177a5cSEmil Constantinescu . -ts_gl_atol <tol> - absolute error 1448*3d177a5cSEmil Constantinescu . -ts_gl_min_order <p> - minimum order method to consider (default=1) 1449*3d177a5cSEmil Constantinescu . -ts_gl_max_order <p> - maximum order method to consider (default=3) 1450*3d177a5cSEmil Constantinescu . -ts_gl_start_order <p> - order of starting method (default=1) 1451*3d177a5cSEmil Constantinescu . -ts_gl_complete <method> - method to use for completing the step (rescale-and-modify or rescale) 1452*3d177a5cSEmil Constantinescu - -ts_adapt_type <method> - adaptive controller to use (none step both) 1453*3d177a5cSEmil Constantinescu 1454*3d177a5cSEmil Constantinescu Notes: 1455*3d177a5cSEmil Constantinescu This integrator can be applied to DAE. 1456*3d177a5cSEmil Constantinescu 1457*3d177a5cSEmil Constantinescu Diagonally implicit general linear (DIGL) methods are a generalization of diagonally implicit Runge-Kutta (DIRK). 1458*3d177a5cSEmil Constantinescu They are represented by the tableau 1459*3d177a5cSEmil Constantinescu 1460*3d177a5cSEmil Constantinescu .vb 1461*3d177a5cSEmil Constantinescu A | U 1462*3d177a5cSEmil Constantinescu ------- 1463*3d177a5cSEmil Constantinescu B | V 1464*3d177a5cSEmil Constantinescu .ve 1465*3d177a5cSEmil Constantinescu 1466*3d177a5cSEmil Constantinescu combined with a vector c of abscissa. "Diagonally implicit" means that A is lower triangular. 1467*3d177a5cSEmil Constantinescu A step of the general method reads 1468*3d177a5cSEmil Constantinescu 1469*3d177a5cSEmil Constantinescu .vb 1470*3d177a5cSEmil Constantinescu [ Y ] = [A U] [ Y' ] 1471*3d177a5cSEmil Constantinescu [X^k] = [B V] [X^{k-1}] 1472*3d177a5cSEmil Constantinescu .ve 1473*3d177a5cSEmil Constantinescu 1474*3d177a5cSEmil Constantinescu where Y is the multivector of stage values, Y' is the multivector of stage derivatives, X^k is the Nordsieck vector of 1475*3d177a5cSEmil Constantinescu the solution at step k. The Nordsieck vector consists of the first r moments of the solution, given by 1476*3d177a5cSEmil Constantinescu 1477*3d177a5cSEmil Constantinescu .vb 1478*3d177a5cSEmil Constantinescu X = [x_0,x_1,...,x_{r-1}] = [x, h x', h^2 x'', ..., h^{r-1} x^{(r-1)} ] 1479*3d177a5cSEmil Constantinescu .ve 1480*3d177a5cSEmil Constantinescu 1481*3d177a5cSEmil Constantinescu If A is lower triangular, we can solve the stages (Y,Y') sequentially 1482*3d177a5cSEmil Constantinescu 1483*3d177a5cSEmil Constantinescu .vb 1484*3d177a5cSEmil Constantinescu y_i = h sum_{j=0}^{s-1} (a_ij y'_j) + sum_{j=0}^{r-1} u_ij x_j, i=0,...,{s-1} 1485*3d177a5cSEmil Constantinescu .ve 1486*3d177a5cSEmil Constantinescu 1487*3d177a5cSEmil Constantinescu and then construct the pieces to carry to the next step 1488*3d177a5cSEmil Constantinescu 1489*3d177a5cSEmil Constantinescu .vb 1490*3d177a5cSEmil Constantinescu xx_i = h sum_{j=0}^{s-1} b_ij y'_j + sum_{j=0}^{r-1} v_ij x_j, i=0,...,{r-1} 1491*3d177a5cSEmil Constantinescu .ve 1492*3d177a5cSEmil Constantinescu 1493*3d177a5cSEmil Constantinescu Note that when the equations are cast in implicit form, we are using the stage equation to define y'_i 1494*3d177a5cSEmil Constantinescu in terms of y_i and known stuff (y_j for j<i and x_j for all j). 1495*3d177a5cSEmil Constantinescu 1496*3d177a5cSEmil Constantinescu 1497*3d177a5cSEmil Constantinescu Error estimation 1498*3d177a5cSEmil Constantinescu 1499*3d177a5cSEmil Constantinescu At present, the most attractive GL methods for stiff problems are singly diagonally implicit schemes which posses 1500*3d177a5cSEmil Constantinescu Inherent Runge-Kutta Stability (IRKS). These methods have r=s, the number of items passed between steps is equal to 1501*3d177a5cSEmil Constantinescu the number of stages. The order and stage-order are one less than the number of stages. We use the error estimates 1502*3d177a5cSEmil Constantinescu in the 2007 paper which provide the following estimates 1503*3d177a5cSEmil Constantinescu 1504*3d177a5cSEmil Constantinescu .vb 1505*3d177a5cSEmil Constantinescu h^{p+1} X^{(p+1)} = phi_0^T Y' + [0 psi_0^T] Xold 1506*3d177a5cSEmil Constantinescu h^{p+2} X^{(p+2)} = phi_1^T Y' + [0 psi_1^T] Xold 1507*3d177a5cSEmil Constantinescu h^{p+2} (dx'/dx) X^{(p+1)} = phi_2^T Y' + [0 psi_2^T] Xold 1508*3d177a5cSEmil Constantinescu .ve 1509*3d177a5cSEmil Constantinescu 1510*3d177a5cSEmil Constantinescu These estimates are accurate to O(h^{p+3}). 1511*3d177a5cSEmil Constantinescu 1512*3d177a5cSEmil Constantinescu Changing the step size 1513*3d177a5cSEmil Constantinescu 1514*3d177a5cSEmil Constantinescu We use the generalized "rescale and modify" scheme, see equation (4.5) of the 2007 paper. 1515*3d177a5cSEmil Constantinescu 1516*3d177a5cSEmil Constantinescu Level: beginner 1517*3d177a5cSEmil Constantinescu 1518*3d177a5cSEmil Constantinescu References: 1519*3d177a5cSEmil Constantinescu + 1. - John Butcher and Z. Jackieweicz and W. Wright, On error propagation in general linear methods for 1520*3d177a5cSEmil Constantinescu ordinary differential equations, Journal of Complexity, Vol 23, 2007. 1521*3d177a5cSEmil Constantinescu - 2. - John Butcher, Numerical methods for ordinary differential equations, second edition, Wiley, 2009. 1522*3d177a5cSEmil Constantinescu 1523*3d177a5cSEmil Constantinescu .seealso: TSCreate(), TS, TSSetType() 1524*3d177a5cSEmil Constantinescu 1525*3d177a5cSEmil Constantinescu M*/ 1526*3d177a5cSEmil Constantinescu #undef __FUNCT__ 1527*3d177a5cSEmil Constantinescu #define __FUNCT__ "TSCreate_GLLE" 1528*3d177a5cSEmil Constantinescu PETSC_EXTERN PetscErrorCode TSCreate_GLLE(TS ts) 1529*3d177a5cSEmil Constantinescu { 1530*3d177a5cSEmil Constantinescu TS_GLLE *gl; 1531*3d177a5cSEmil Constantinescu PetscErrorCode ierr; 1532*3d177a5cSEmil Constantinescu 1533*3d177a5cSEmil Constantinescu PetscFunctionBegin; 1534*3d177a5cSEmil Constantinescu ierr = TSGLLEInitializePackage();CHKERRQ(ierr); 1535*3d177a5cSEmil Constantinescu 1536*3d177a5cSEmil Constantinescu ierr = PetscNewLog(ts,&gl);CHKERRQ(ierr); 1537*3d177a5cSEmil Constantinescu ts->data = (void*)gl; 1538*3d177a5cSEmil Constantinescu 1539*3d177a5cSEmil Constantinescu ts->ops->reset = TSReset_GLLE; 1540*3d177a5cSEmil Constantinescu ts->ops->destroy = TSDestroy_GLLE; 1541*3d177a5cSEmil Constantinescu ts->ops->view = TSView_GLLE; 1542*3d177a5cSEmil Constantinescu ts->ops->setup = TSSetUp_GLLE; 1543*3d177a5cSEmil Constantinescu ts->ops->solve = TSSolve_GLLE; 1544*3d177a5cSEmil Constantinescu ts->ops->setfromoptions = TSSetFromOptions_GLLE; 1545*3d177a5cSEmil Constantinescu ts->ops->snesfunction = SNESTSFormFunction_GLLE; 1546*3d177a5cSEmil Constantinescu ts->ops->snesjacobian = SNESTSFormJacobian_GLLE; 1547*3d177a5cSEmil Constantinescu 1548*3d177a5cSEmil Constantinescu gl->max_step_rejections = 1; 1549*3d177a5cSEmil Constantinescu gl->min_order = 1; 1550*3d177a5cSEmil Constantinescu gl->max_order = 3; 1551*3d177a5cSEmil Constantinescu gl->start_order = 1; 1552*3d177a5cSEmil Constantinescu gl->current_scheme = -1; 1553*3d177a5cSEmil Constantinescu gl->extrapolate = PETSC_FALSE; 1554*3d177a5cSEmil Constantinescu 1555*3d177a5cSEmil Constantinescu gl->wrms_atol = 1e-8; 1556*3d177a5cSEmil Constantinescu gl->wrms_rtol = 1e-5; 1557*3d177a5cSEmil Constantinescu 1558*3d177a5cSEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSGLLESetType_C", &TSGLLESetType_GLLE);CHKERRQ(ierr); 1559*3d177a5cSEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSGLLESetAcceptType_C",&TSGLLESetAcceptType_GLLE);CHKERRQ(ierr); 1560*3d177a5cSEmil Constantinescu ierr = PetscObjectComposeFunction((PetscObject)ts,"TSGLLEGetAdapt_C", &TSGLLEGetAdapt_GLLE);CHKERRQ(ierr); 1561*3d177a5cSEmil Constantinescu PetscFunctionReturn(0); 1562*3d177a5cSEmil Constantinescu } 1563