14b84eec9SHong Zhang /* 2*f944a40eSHong Zhang Code for time stepping with the Multirate Partitioned Runge-Kutta method 34b84eec9SHong Zhang 44b84eec9SHong Zhang Notes: 54b84eec9SHong Zhang 1) The general system is written as 6*f944a40eSHong Zhang Udot = F(t,U) 7*f944a40eSHong Zhang if one does not split the RHS function, but gives the indexes for both slow and fast components; 84b84eec9SHong Zhang 2) The general system is written as 94b84eec9SHong Zhang Usdot = Fs(t,Us,Uf) 10*f944a40eSHong Zhang Ufdot = Ff(t,Us,Uf) 11*f944a40eSHong Zhang for component-wise partitioned system, 12*f944a40eSHong Zhang users should split the RHS function themselves and also provide the indexes for both slow and fast components. 1319c2959aSHong Zhang 3) To correct The confusing terminology in the paper, we use 'slow method', 'slow buffer method' and 'fast method' to denote the methods applied to 'slow region', 'slow buffer region' and 'fast region' respectively. The 'slow method' in the original paper actually means the 'slow buffer method'. 1419c2959aSHong Zhang 4) Why does the buffer region have to be inside the slow region? The buffer region is treated with a slow method essentially. Applying the slow method to a region with a fast characteristic time scale is apparently not a good choice. 1519c2959aSHong Zhang 1619c2959aSHong Zhang Reference: 1719c2959aSHong Zhang Emil M. Constantinescu, Adrian Sandu, Multirate Timestepping Methods for Hyperbolic Conservation Laws, Journal of Scientific Computing 2007 184b84eec9SHong Zhang */ 194b84eec9SHong Zhang 204b84eec9SHong Zhang #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 214b84eec9SHong Zhang #include <petscdm.h> 224b84eec9SHong Zhang 2319c2959aSHong Zhang static TSMPRKType TSMPRKDefault = TSMPRK2A22; 244b84eec9SHong Zhang static PetscBool TSMPRKRegisterAllCalled; 254b84eec9SHong Zhang static PetscBool TSMPRKPackageInitialized; 264b84eec9SHong Zhang 274b84eec9SHong Zhang typedef struct _MPRKTableau *MPRKTableau; 284b84eec9SHong Zhang struct _MPRKTableau { 294b84eec9SHong Zhang char *name; 304b84eec9SHong Zhang PetscInt order; /* Classical approximation order of the method i */ 3179bc8a77SHong Zhang PetscInt sbase; /* Number of stages in the base method*/ 324b84eec9SHong Zhang PetscInt s; /* Number of stages */ 3319c2959aSHong Zhang PetscInt np; /* Number of partitions */ 344b84eec9SHong Zhang PetscReal *Af,*bf,*cf; /* Tableau for fast components */ 3519c2959aSHong Zhang PetscReal *Amb,*bmb,*cmb; /* Tableau for medium components */ 3619c2959aSHong Zhang PetscInt *rmb; /* Array of flags for repeated stages in medium method */ 3719c2959aSHong Zhang PetscReal *Asb,*bsb,*csb; /* Tableau for slow components */ 3819c2959aSHong Zhang PetscInt *rsb; /* Array of flags for repeated staged in slow method*/ 394b84eec9SHong Zhang }; 404b84eec9SHong Zhang typedef struct _MPRKTableauLink *MPRKTableauLink; 414b84eec9SHong Zhang struct _MPRKTableauLink { 424b84eec9SHong Zhang struct _MPRKTableau tab; 434b84eec9SHong Zhang MPRKTableauLink next; 444b84eec9SHong Zhang }; 454b84eec9SHong Zhang static MPRKTableauLink MPRKTableauList; 464b84eec9SHong Zhang 474b84eec9SHong Zhang typedef struct { 484b84eec9SHong Zhang MPRKTableau tableau; 494b84eec9SHong Zhang Vec *Y; /* States computed during the step */ 507c0df07dSHong Zhang Vec *YdotRHS; 514b84eec9SHong Zhang Vec *YdotRHS_slow; /* Function evaluations by slow tableau for slow components */ 5219c2959aSHong Zhang Vec *YdotRHS_slowbuffer; /* Function evaluations by slow tableau for slow components */ 5319c2959aSHong Zhang Vec *YdotRHS_medium; /* Function evaluations by slow tableau for slow components */ 5419c2959aSHong Zhang Vec *YdotRHS_mediumbuffer; /* Function evaluations by slow tableau for slow components */ 5519c2959aSHong Zhang Vec *YdotRHS_fast; /* Function evaluations by fast tableau for fast components */ 564b84eec9SHong Zhang PetscScalar *work_slow; /* Scalar work_slow by slow tableau */ 5719c2959aSHong Zhang PetscScalar *work_slowbuffer; /* Scalar work_slow by slow tableau */ 5819c2959aSHong Zhang PetscScalar *work_medium; /* Scalar work_slow by medium tableau */ 5919c2959aSHong Zhang PetscScalar *work_mediumbuffer; /* Scalar work_slow by medium tableau */ 6019c2959aSHong Zhang PetscScalar *work_fast; /* Scalar work_fast by fast tableau */ 614b84eec9SHong Zhang PetscReal stage_time; 624b84eec9SHong Zhang TSStepStatus status; 634b84eec9SHong Zhang PetscReal ptime; 644b84eec9SHong Zhang PetscReal time_step; 6519c2959aSHong Zhang IS is_slow,is_slowbuffer,is_medium,is_mediumbuffer,is_fast; 6619c2959aSHong Zhang TS subts_slow,subts_slowbuffer,subts_medium,subts_mediumbuffer,subts_fast; 674b84eec9SHong Zhang } TS_MPRK; 684b84eec9SHong Zhang 6919c2959aSHong Zhang static PetscErrorCode TSMPRKGenerateTableau2(PetscInt ratio,PetscInt s,const PetscReal Abase[],const PetscReal bbase[],PetscReal A1[],PetscReal b1[],PetscReal A2[],PetscReal b2[]) 7019c2959aSHong Zhang { 7119c2959aSHong Zhang PetscInt i,j,k,l; 7219c2959aSHong Zhang 7319c2959aSHong Zhang PetscFunctionBegin; 7419c2959aSHong Zhang for (k=0; k<ratio; k++) { 7519c2959aSHong Zhang /* diagonal blocks */ 7619c2959aSHong Zhang for (i=0; i<s; i++) 7719c2959aSHong Zhang for (j=0; j<s; j++) { 7819c2959aSHong Zhang A1[(k*s+i)*ratio*s+k*s+j] = Abase[i*s+j]; 7919c2959aSHong Zhang A2[(k*s+i)*ratio*s+k*s+j] = Abase[i*s+j]/ratio; 8019c2959aSHong Zhang } 8119c2959aSHong Zhang /* off diagonal blocks */ 8219c2959aSHong Zhang for (l=0; l<k; l++) 8319c2959aSHong Zhang for (i=0; i<s; i++) 8419c2959aSHong Zhang for (j=0; j<s; j++) 8519c2959aSHong Zhang A2[(k*s+i)*ratio*s+l*s+j] = bbase[j]/ratio; 8619c2959aSHong Zhang for (j=0; j<s; j++) { 8719c2959aSHong Zhang b1[k*s+j] = bbase[j]/ratio; 8819c2959aSHong Zhang b2[k*s+j] = bbase[j]/ratio; 8919c2959aSHong Zhang } 9019c2959aSHong Zhang } 9119c2959aSHong Zhang PetscFunctionReturn(0); 9219c2959aSHong Zhang } 9319c2959aSHong Zhang 9419c2959aSHong Zhang static PetscErrorCode TSMPRKGenerateTableau3(PetscInt ratio,PetscInt s,const PetscReal Abase[],const PetscReal bbase[],PetscReal A1[],PetscReal b1[],PetscReal A2[],PetscReal b2[],PetscReal A3[],PetscReal b3[]) 9519c2959aSHong Zhang { 9619c2959aSHong Zhang PetscInt i,j,k,l,m,n; 9719c2959aSHong Zhang 9819c2959aSHong Zhang PetscFunctionBegin; 9919c2959aSHong Zhang for (k=0; k<ratio; k++) { /* diagonal blocks of size ratio*s by ratio*s */ 10019c2959aSHong Zhang for (l=0; l<ratio; l++) /* diagonal sub-blocks of size s by s */ 10119c2959aSHong Zhang for (i=0; i<s; i++) 10219c2959aSHong Zhang for (j=0; j<s; j++) { 10319c2959aSHong Zhang A1[((k*ratio+l)*s+i)*ratio*ratio*s+(k*ratio+l)*s+j] = Abase[i*s+j]; 10419c2959aSHong Zhang A2[((k*ratio+l)*s+i)*ratio*ratio*s+(k*ratio+l)*s+j] = Abase[i*s+j]/ratio; 10519c2959aSHong Zhang A3[((k*ratio+l)*s+i)*ratio*ratio*s+(k*ratio+l)*s+j] = Abase[i*s+j]/ratio/ratio; 10619c2959aSHong Zhang } 10719c2959aSHong Zhang for (l=0; l<k; l++) /* off-diagonal blocks of size ratio*s by ratio*s */ 10819c2959aSHong Zhang for (m=0; m<ratio; m++) 10919c2959aSHong Zhang for (n=0; n<ratio; n++) 11019c2959aSHong Zhang for (i=0; i<s; i++) 11119c2959aSHong Zhang for (j=0; j<s; j++) { 11219c2959aSHong Zhang A2[((k*ratio+m)*s+i)*ratio*ratio*s+(l*ratio+n)*s+j] = bbase[j]/ratio/ratio; 11319c2959aSHong Zhang A3[((k*ratio+m)*s+i)*ratio*ratio*s+(l*ratio+n)*s+j] = bbase[j]/ratio/ratio; 11419c2959aSHong Zhang } 11519c2959aSHong Zhang for (m=0; m<ratio; m++) 11619c2959aSHong Zhang for (n=0; n<m; n++) /* off-diagonal sub-blocks of size s by s in the diagonal blocks */ 11719c2959aSHong Zhang for (i=0; i<s; i++) 11819c2959aSHong Zhang for (j=0; j<s; j++) 11919c2959aSHong Zhang A3[((k*ratio+m)*s+i)*ratio*ratio*s+(k*ratio+n)*s+j] = bbase[j]/ratio/ratio; 12019c2959aSHong Zhang for (n=0; n<ratio; n++) 12119c2959aSHong Zhang for (j=0; j<s; j++) { 12219c2959aSHong Zhang b1[(k*ratio+n)*s+j] = bbase[j]/ratio/ratio; 12319c2959aSHong Zhang b2[(k*ratio+n)*s+j] = bbase[j]/ratio/ratio; 12419c2959aSHong Zhang b3[(k*ratio+n)*s+j] = bbase[j]/ratio/ratio; 12519c2959aSHong Zhang } 12619c2959aSHong Zhang } 12719c2959aSHong Zhang PetscFunctionReturn(0); 12819c2959aSHong Zhang } 12919c2959aSHong Zhang 1304b84eec9SHong Zhang /*MC 131*f944a40eSHong Zhang TSMPRK2A22 - Second Order Multirate Partitioned Runge Kutta scheme based on RK2A. 1324b84eec9SHong Zhang 13319c2959aSHong Zhang This method has four stages for slow and fast parts. The refinement factor of the stepsize is 2. 13419c2959aSHong Zhang r = 2, np = 2 1354b84eec9SHong Zhang Options database: 13619c2959aSHong Zhang . -ts_mprk_type 2a22 1374b84eec9SHong Zhang 1384b84eec9SHong Zhang Level: advanced 1394b84eec9SHong Zhang 1404b84eec9SHong Zhang .seealso: TSMPRK, TSMPRKType, TSMPRKSetType() 1414b84eec9SHong Zhang M*/ 1424b84eec9SHong Zhang /*MC 143*f944a40eSHong Zhang TSMPRK2A23 - Second Order Multirate Partitioned Runge-Kutta scheme based on RK2A. 14419c2959aSHong Zhang 14519c2959aSHong Zhang This method has eight stages for slow and medium and fast parts. The refinement factor of the stepsize is 2. 14619c2959aSHong Zhang r = 2, np = 3 14719c2959aSHong Zhang Options database: 14819c2959aSHong Zhang . -ts_mprk_type 2a23 14919c2959aSHong Zhang 15019c2959aSHong Zhang Level: advanced 15119c2959aSHong Zhang 15219c2959aSHong Zhang .seealso: TSMPRK, TSMPRKType, TSMPRKSetType() 15319c2959aSHong Zhang M*/ 15419c2959aSHong Zhang /*MC 155*f944a40eSHong Zhang TSMPRK2A32 - Second Order Multirate Partitioned Runge-Kutta scheme based on RK2A. 15619c2959aSHong Zhang 15719c2959aSHong Zhang This method has four stages for slow and fast parts. The refinement factor of the stepsize is 3. 15819c2959aSHong Zhang r = 3, np = 2 15919c2959aSHong Zhang Options database: 16019c2959aSHong Zhang . -ts_mprk_type 2a32 16119c2959aSHong Zhang 16219c2959aSHong Zhang Level: advanced 16319c2959aSHong Zhang 16419c2959aSHong Zhang .seealso: TSMPRK, TSMPRKType, TSMPRKSetType() 16519c2959aSHong Zhang M*/ 16619c2959aSHong Zhang /*MC 167*f944a40eSHong Zhang TSMPRK2A33 - Second Order Multirate Partitioned Runge-Kutta scheme based on RK2A. 16819c2959aSHong Zhang 16919c2959aSHong Zhang This method has eight stages for slow and medium and fast parts. The refinement factor of the stepsize is 3. 17019c2959aSHong Zhang r = 3, np = 3 17119c2959aSHong Zhang Options database: 17219c2959aSHong Zhang . -ts_mprk_type 2a33 17319c2959aSHong Zhang 17419c2959aSHong Zhang Level: advanced 17519c2959aSHong Zhang 17619c2959aSHong Zhang .seealso: TSMPRK, TSMPRKType, TSMPRKSetType() 17719c2959aSHong Zhang M*/ 17819c2959aSHong Zhang /*MC 179*f944a40eSHong Zhang TSMPRK3P2M - Third Order Multirate Partitioned Runge-Kutta scheme. 1804b84eec9SHong Zhang 1814b84eec9SHong Zhang This method has eight stages for both slow and fast parts. 1824b84eec9SHong Zhang 1834b84eec9SHong Zhang Options database: 1844b84eec9SHong Zhang . -ts_mprk_type pm3 (put here temporarily) 1854b84eec9SHong Zhang 1864b84eec9SHong Zhang Level: advanced 1874b84eec9SHong Zhang 1884b84eec9SHong Zhang .seealso: TSMPRK, TSMPRKType, TSMPRKSetType() 1894b84eec9SHong Zhang M*/ 1904b84eec9SHong Zhang /*MC 191*f944a40eSHong Zhang TSMPRKP2 - Second Order Multirate Partitioned Runge-Kutta scheme. 1924b84eec9SHong Zhang 1934b84eec9SHong Zhang This method has five stages for both slow and fast parts. 1944b84eec9SHong Zhang 1954b84eec9SHong Zhang Options database: 1964b84eec9SHong Zhang . -ts_mprk_type p2 1974b84eec9SHong Zhang 1984b84eec9SHong Zhang Level: advanced 1994b84eec9SHong Zhang 2004b84eec9SHong Zhang .seealso: TSMPRK, TSMPRKType, TSMPRKSetType() 2014b84eec9SHong Zhang M*/ 2024b84eec9SHong Zhang /*MC 203*f944a40eSHong Zhang TSMPRKP3 - Third Order Multirate Partitioned Runge-Kutta scheme. 2044b84eec9SHong Zhang 2054b84eec9SHong Zhang This method has ten stages for both slow and fast parts. 2064b84eec9SHong Zhang 2074b84eec9SHong Zhang Options database: 2084b84eec9SHong Zhang . -ts_mprk_type p3 2094b84eec9SHong Zhang 2104b84eec9SHong Zhang Level: advanced 2114b84eec9SHong Zhang 2124b84eec9SHong Zhang .seealso: TSMPRK, TSMPRKType, TSMPRKSetType() 2134b84eec9SHong Zhang M*/ 2144b84eec9SHong Zhang 2154b84eec9SHong Zhang /*@C 2164b84eec9SHong Zhang TSMPRKRegisterAll - Registers all of the Partirioned Runge-Kutta explicit methods in TSMPRK 2174b84eec9SHong Zhang 2184b84eec9SHong Zhang Not Collective, but should be called by all processes which will need the schemes to be registered 2194b84eec9SHong Zhang 2204b84eec9SHong Zhang Level: advanced 2214b84eec9SHong Zhang 2224b84eec9SHong Zhang .keywords: TS, TSMPRK, register, all 2234b84eec9SHong Zhang 2244b84eec9SHong Zhang .seealso: TSMPRKRegisterDestroy() 2254b84eec9SHong Zhang @*/ 2264b84eec9SHong Zhang PetscErrorCode TSMPRKRegisterAll(void) 2274b84eec9SHong Zhang { 2284b84eec9SHong Zhang PetscErrorCode ierr; 2294b84eec9SHong Zhang 2304b84eec9SHong Zhang PetscFunctionBegin; 2314b84eec9SHong Zhang if (TSMPRKRegisterAllCalled) PetscFunctionReturn(0); 2324b84eec9SHong Zhang TSMPRKRegisterAllCalled = PETSC_TRUE; 2334b84eec9SHong Zhang 2344b84eec9SHong Zhang #define RC PetscRealConstant 2354b84eec9SHong Zhang { 2364b84eec9SHong Zhang const PetscReal 23719c2959aSHong Zhang Abase[2][2] = {{0,0}, 23819c2959aSHong Zhang {RC(1.0),0}}, 23919c2959aSHong Zhang bbase[2] = {RC(0.5),RC(0.5)}; 24019c2959aSHong Zhang PetscReal 24119c2959aSHong Zhang Asb[4][4] = {{0}},Af[4][4] = {{0}},bsb[4] = {0},bf[4] = {0}; 24219c2959aSHong Zhang PetscInt 24319c2959aSHong Zhang rsb[4] = {0,0,1,2}; 24419c2959aSHong Zhang ierr = TSMPRKGenerateTableau2(2,2,&Abase[0][0],bbase,&Asb[0][0],bsb,&Af[0][0],bf);CHKERRQ(ierr); 24579bc8a77SHong Zhang ierr = TSMPRKRegister(TSMPRK2A22,2,2,2,1,&Asb[0][0],bsb,NULL,rsb,NULL,NULL,NULL,NULL,&Af[0][0],bf,NULL);CHKERRQ(ierr); 2464b84eec9SHong Zhang } 24719c2959aSHong Zhang { 24819c2959aSHong Zhang const PetscReal 24919c2959aSHong Zhang Abase[2][2] = {{0,0}, 25019c2959aSHong Zhang {RC(1.0),0}}, 25119c2959aSHong Zhang bbase[2] = {RC(0.5),RC(0.5)}; 25219c2959aSHong Zhang PetscReal 25319c2959aSHong Zhang Asb[8][8] = {{0}},Amb[8][8] = {{0}},Af[8][8] = {{0}},bsb[8] ={0},bmb[8] = {0},bf[8] = {0}; 25419c2959aSHong Zhang PetscInt 25519c2959aSHong Zhang rsb[8] = {0,0,1,2,1,2,1,2},rmb[8] = {0,0,1,2,0,0,5,6}; 25619c2959aSHong Zhang ierr = TSMPRKGenerateTableau3(2,2,&Abase[0][0],bbase,&Asb[0][0],bsb,&Amb[0][0],bmb,&Af[0][0],bf);CHKERRQ(ierr); 25779bc8a77SHong Zhang ierr = TSMPRKRegister(TSMPRK2A23,2,2,2,2,&Asb[0][0],bsb,NULL,rsb,&Amb[0][0],bmb,NULL,rmb,&Af[0][0],bf,NULL);CHKERRQ(ierr); 25819c2959aSHong Zhang } 25919c2959aSHong Zhang { 26019c2959aSHong Zhang const PetscReal 26119c2959aSHong Zhang Abase[2][2] = {{0,0}, 26219c2959aSHong Zhang {RC(1.0),0}}, 26319c2959aSHong Zhang bbase[2] = {RC(0.5),RC(0.5)}; 26419c2959aSHong Zhang PetscReal 26519c2959aSHong Zhang Asb[6][6] = {{0}},Af[6][6] = {{0}},bsb[6] = {0},bf[6] = {0}; 26619c2959aSHong Zhang PetscInt 26719c2959aSHong Zhang rsb[6] = {0,0,1,2,1,2}; 26819c2959aSHong Zhang ierr = TSMPRKGenerateTableau2(3,2,&Abase[0][0],bbase,&Asb[0][0],bsb,&Af[0][0],bf);CHKERRQ(ierr); 26979bc8a77SHong Zhang ierr = TSMPRKRegister(TSMPRK2A32,2,2,3,1,&Asb[0][0],bsb,NULL,rsb,NULL,NULL,NULL,NULL,&Af[0][0],bf,NULL);CHKERRQ(ierr); 27019c2959aSHong Zhang } 27119c2959aSHong Zhang { 27219c2959aSHong Zhang const PetscReal 27319c2959aSHong Zhang Abase[2][2] = {{0,0}, 27419c2959aSHong Zhang {RC(1.0),0}}, 27519c2959aSHong Zhang bbase[2] = {RC(0.5),RC(0.5)}; 27619c2959aSHong Zhang PetscReal 27719c2959aSHong Zhang Asb[18][18] = {{0}},Amb[18][18] = {{0}},Af[18][18] = {{0}},bsb[18] ={0},bmb[18] = {0},bf[18] = {0}; 27819c2959aSHong Zhang PetscInt 27919c2959aSHong Zhang rsb[18] = {0,0,1,2,1,2,1,2,1,2,1,2,1,2,1,2,1,2},rmb[18] = {0,0,1,2,1,2,0,0,7,8,7,8,0,0,13,14,13,14}; 28019c2959aSHong Zhang ierr = TSMPRKGenerateTableau3(3,2,&Abase[0][0],bbase,&Asb[0][0],bsb,&Amb[0][0],bmb,&Af[0][0],bf);CHKERRQ(ierr); 28179bc8a77SHong Zhang ierr = TSMPRKRegister(TSMPRK2A33,2,2,3,3,&Asb[0][0],bsb,NULL,rsb,&Amb[0][0],bmb,NULL,rmb,&Af[0][0],bf,NULL);CHKERRQ(ierr); 28219c2959aSHong Zhang } 28319c2959aSHong Zhang /* 28419c2959aSHong Zhang PetscReal 28519c2959aSHong Zhang Asb[8][8] = {{Abase[0][0],Abase[0][1],0,0,0,0,0,0}, 28619c2959aSHong Zhang {Abase[1][0],Abase[1][1],0,0,0,0,0,0}, 28719c2959aSHong Zhang {0,0,Abase[0][0],Abase[0][1],0,0,0,0}, 28819c2959aSHong Zhang {0,0,Abase[1][0],Abase[1][1],0,0,0,0}, 28919c2959aSHong Zhang {0,0,0,0,Abase[0][0],Abase[0][1],0,0}, 29019c2959aSHong Zhang {0,0,0,0,Abase[1][0],Abase[1][1],0,0}, 29119c2959aSHong Zhang {0,0,0,0,0,0,Abase[0][0],Abase[0][1]}, 29219c2959aSHong Zhang {0,0,0,0,0,0,Abase[1][0],Abase[1][1]}}, 29319c2959aSHong Zhang Amb[8][8] = {{Abase[0][0]/m,Abase[0][1]/m,0,0,0,0,0,0}, 29419c2959aSHong Zhang {Abase[1][0]/m,Abase[1][1]/m,0,0,0,0,0,0}, 29519c2959aSHong Zhang {0,0,Abase[0][0]/m,Abase[0][1]/m,0,0,0,0}, 29619c2959aSHong Zhang {0,0,Abase[1][0]/m,Abase[1][1]/m,0,0,0,0}, 29719c2959aSHong Zhang {bbase[0]/m,bbase[1]/m,bbase[0]/m,bbase[1]/m,Abase[0][0]/m,Abase[0][1]/m,0,0}, 29819c2959aSHong Zhang {bbase[0]/m,bbase[1]/m,bbase[0]/m,bbase[1]/m,Abase[1][0]/m,Abase[1][1]/m,0,0}, 29919c2959aSHong Zhang {bbase[0]/m,bbase[1]/m,bbase[0]/m,bbase[1]/m,0,0,Abase[0][0]/m,Abase[0][1]/m}, 30019c2959aSHong Zhang {bbase[0]/m,bbase[1]/m,bbase[0]/m,bbase[1]/m,0,0,Abase[1][0]/m,Abase[1][1]/m}}, 30119c2959aSHong Zhang Af[8][8] = {{Abase[0][0]/m/m,Abase[0][1]/m/m,0,0,0,0,0,0}, 30219c2959aSHong Zhang {Abase[1][0]/m/m,Abase[1][1]/m/m,0,0,0,0,0,0}, 30319c2959aSHong Zhang {0,0,Abase[0][0]/m/m,Abase[0][1]/m/m,0,0,0,0}, 30419c2959aSHong Zhang {0,0,Abase[1][0]/m/m,Abase[1][1]/m/m,0,0,0,0}, 30519c2959aSHong Zhang {bbase[0]/m/m,bbase[1]/m/m,bbase[0]/m/m,bbase[1]/m/m,Abase[0][0]/m/m,Abase[0][1]/m/m,0,0}, 30619c2959aSHong Zhang {bbase[0]/m/m,bbase[1]/m/m,bbase[0]/m/m,bbase[1]/m/m,Abase[1][0]/m/m,Abase[1][1]/m/m,0,0}, 30719c2959aSHong Zhang {bbase[0]/m/m,bbase[1]/m/m,bbase[0]/m/m,bbase[1]/m/m,bbase[0]/m/m,bbase[1]/m/m,Abase[0][0]/m,Abase[0][1]/m}, 30819c2959aSHong Zhang {bbase[0]/m/m,bbase[1]/m/m,bbase[0]/m/m,bbase[1]/m/m,bbase[0]/m/m,bbase[1]/m/m,Abase[1][0]/m,Abase[1][1]/m}}, 30919c2959aSHong Zhang bsb[8] = {bbase[0]/m/m,bbase[1]/m/m,bbase[0]/m/m,bbase[1]/m/m,bbase[0]/m/m,bbase[1]/m/m,bbase[0]/m/m,bbase[1]/m/m}, 31019c2959aSHong Zhang bmb[8] = {bbase[0]/m/m,bbase[1]/m/m,bbase[0]/m/m,bbase[1]/m/m,bbase[0]/m/m,bbase[1]/m/m,bbase[0]/m,bbase[1]/m/m}, 31119c2959aSHong Zhang bf[8] = {bbase[0]/m/m,bbase[1]/m/m,bbase[0]/m/m,bbase[1]/m/m,bbase[0]/m/m,bbase[1]/m,bbase[0]/m/m,bbase[1]/m/m}, 31219c2959aSHong Zhang */ 3134b84eec9SHong Zhang /*{ 3144b84eec9SHong Zhang const PetscReal 3154b84eec9SHong Zhang As[8][8] = {{0,0,0,0,0,0,0,0}, 3164b84eec9SHong Zhang {RC(1.0)/RC(2.0),0,0,0,0,0,0,0}, 3174b84eec9SHong Zhang {RC(-1.0)/RC(6.0),RC(2.0)/RC(3.0),0,0,0,0,0,0}, 3184b84eec9SHong Zhang {RC(1.0)/RC(3.0),RC(-1.0)/RC(3.0),RC(1.0),0,0,0,0,0}, 3194b84eec9SHong Zhang {0,0,0,0,0,0,0,0}, 3204b84eec9SHong Zhang {0,0,0,0,RC(1.0)/RC(2.0),0,0,0}, 3214b84eec9SHong Zhang {0,0,0,0,RC(-1.0)/RC(6.0),RC(2.0)/RC(3.0),0,0}, 3224b84eec9SHong Zhang {0,0,0,0,RC(1.0)/RC(3.0),RC(-1.0)/RC(3.0),RC(1.0),0}}, 3234b84eec9SHong Zhang A[8][8] = {{0,0,0,0,0,0,0,0}, 3244b84eec9SHong Zhang {RC(1.0)/RC(4.0),0,0,0,0,0,0,0}, 3254b84eec9SHong Zhang {RC(-1.0)/RC(12.0),RC(1.0)/RC(3.0),0,0,0,0,0,0}, 3264b84eec9SHong Zhang {RC(1.0)/RC(6.0),RC(-1.0)/RC(6.0),RC(1.0)/RC(2.0),0,0,0,0,0}, 3274b84eec9SHong Zhang {RC(1.0)/RC(12.0),RC(1.0)/RC(6.0),RC(1.0)/RC(6.0),RC(1.0)/RC(12.0),0,0,0,0}, 3284b84eec9SHong Zhang {RC(1.0)/RC(12.0),RC(1.0)/RC(6.0),RC(1.0)/RC(6.0),RC(1.0)/RC(12.0),RC(1.0)/RC(4.0),0,0,0}, 3294b84eec9SHong Zhang {RC(1.0)/RC(12.0),RC(1.0)/RC(6.0),RC(1.0)/RC(6.0),RC(1.0)/RC(12.0),RC(-1.0)/RC(12.0),RC(1.0)/RC(3.0),0,0}, 3304b84eec9SHong Zhang {RC(1.0)/RC(12.0),RC(1.0)/RC(6.0),RC(1.0)/RC(6.0),RC(1.0)/RC(12.0),RC(1.0)/RC(6.0),RC(-1.0)/RC(6.0),RC(1.0)/RC(2.0),0}}, 3314b84eec9SHong Zhang bs[8] = {RC(1.0)/RC(12.0),RC(1.0)/RC(6.0),RC(1.0)/RC(6.0),RC(1.0)/RC(12.0),RC(1.0)/RC(12.0),RC(1.0)/RC(6.0),RC(1.0)/RC(6.0),RC(1.0)/RC(12.0)}, 3324b84eec9SHong Zhang b[8] = {RC(1.0)/RC(12.0),RC(1.0)/RC(6.0),RC(1.0)/RC(6.0),RC(1.0)/RC(12.0),RC(1.0)/RC(12.0),RC(1.0)/RC(6.0),RC(1.0)/RC(6.0),RC(1.0)/RC(12.0)}; 3334b84eec9SHong Zhang ierr = TSMPRKRegister(TSMPRKPM3,3,8,&As[0][0],bs,NULL,&A[0][0],b,NULL);CHKERRQ(ierr); 3344b84eec9SHong Zhang }*/ 3354b84eec9SHong Zhang 3364b84eec9SHong Zhang { 3374b84eec9SHong Zhang const PetscReal 33819c2959aSHong Zhang Asb[5][5] = {{0,0,0,0,0}, 3394b84eec9SHong Zhang {RC(1.0)/RC(2.0),0,0,0,0}, 3404b84eec9SHong Zhang {RC(1.0)/RC(2.0),0,0,0,0}, 3414b84eec9SHong Zhang {RC(1.0),0,0,0,0}, 3424b84eec9SHong Zhang {RC(1.0),0,0,0,0}}, 34319c2959aSHong Zhang Af[5][5] = {{0,0,0,0,0}, 3444b84eec9SHong Zhang {RC(1.0)/RC(2.0),0,0,0,0}, 3454b84eec9SHong Zhang {RC(1.0)/RC(4.0),RC(1.0)/RC(4.0),0,0,0}, 3464b84eec9SHong Zhang {RC(1.0)/RC(4.0),RC(1.0)/RC(4.0),RC(1.0)/RC(2.0),0,0}, 3474b84eec9SHong Zhang {RC(1.0)/RC(4.0),RC(1.0)/RC(4.0),RC(1.0)/RC(4.0),RC(1.0)/RC(4.0),0}}, 34819c2959aSHong Zhang bsb[5] = {RC(1.0)/RC(2.0),0,0,0,RC(1.0)/RC(2.0)}, 34919c2959aSHong Zhang bf[5] = {RC(1.0)/RC(4.0),RC(1.0)/RC(4.0),RC(1.0)/RC(4.0),RC(1.0)/RC(4.0),0}; 35019c2959aSHong Zhang const PetscInt 35119c2959aSHong Zhang rsb[5] = {0,0,2,0,4}; 35279bc8a77SHong Zhang ierr = TSMPRKRegister(TSMPRKP2,2,5,1,1,&Asb[0][0],bsb,NULL,rsb,NULL,NULL,NULL,NULL,&Af[0][0],bf,NULL);CHKERRQ(ierr); 3534b84eec9SHong Zhang } 3544b84eec9SHong Zhang 3554b84eec9SHong Zhang { 3564b84eec9SHong Zhang const PetscReal 35719c2959aSHong Zhang Asb[10][10] = {{0,0,0,0,0,0,0,0,0,0}, 3584b84eec9SHong Zhang {RC(1.0)/RC(4.0),0,0,0,0,0,0,0,0,0}, 3594b84eec9SHong Zhang {RC(1.0)/RC(4.0),0,0,0,0,0,0,0,0,0}, 3604b84eec9SHong Zhang {RC(1.0)/RC(2.0),0,0,0,0,0,0,0,0,0}, 3614b84eec9SHong Zhang {RC(1.0)/RC(2.0),0,0,0,0,0,0,0,0,0}, 3624b84eec9SHong Zhang {RC(-1.0)/RC(6.0),0,0,0,RC(2.0)/RC(3.0),0,0,0,0,0}, 3634b84eec9SHong Zhang {RC(1.0)/RC(12.0),0,0,0,RC(1.0)/RC(6.0),RC(1.0)/RC(2.0),0,0,0,0}, 3644b84eec9SHong Zhang {RC(1.0)/RC(12.0),0,0,0,RC(1.0)/RC(6.0),RC(1.0)/RC(2.0),0,0,0,0}, 3654b84eec9SHong Zhang {RC(1.0)/RC(3.0),0,0,0,RC(-1.0)/RC(3.0),RC(1.0),0,0,0,0}, 3664b84eec9SHong Zhang {RC(1.0)/RC(3.0),0,0,0,RC(-1.0)/RC(3.0),RC(1.0),0,0,0,0}}, 36719c2959aSHong Zhang Af[10][10] = {{0,0,0,0,0,0,0,0,0,0}, 3684b84eec9SHong Zhang {RC(1.0)/RC(4.0),0,0,0,0,0,0,0,0,0}, 3694b84eec9SHong Zhang {RC(-1.0)/RC(12.0),RC(1.0)/RC(3.0),0,0,0,0,0,0,0,0}, 3704b84eec9SHong Zhang {RC(1.0)/RC(6.0),RC(-1.0)/RC(6.0),RC(1.0)/RC(2.0),0,0,0,0,0,0,0}, 3714b84eec9SHong Zhang {RC(1.0)/RC(12.0),RC(1.0)/RC(6.0),RC(1.0)/RC(6.0),RC(1.0)/RC(12.0),0,0,0,0,0,0}, 3724b84eec9SHong Zhang {RC(1.0)/RC(12.0),RC(1.0)/RC(6.0),RC(1.0)/RC(6.0),RC(1.0)/RC(12.0),0,0,0,0,0,0}, 3734b84eec9SHong Zhang {RC(1.0)/RC(12.0),RC(1.0)/RC(6.0),RC(1.0)/RC(6.0),RC(1.0)/RC(12.0),0,RC(1.0)/RC(4.0),0,0,0,0}, 3744b84eec9SHong Zhang {RC(1.0)/RC(12.0),RC(1.0)/RC(6.0),RC(1.0)/RC(6.0),RC(1.0)/RC(12.0),0,RC(-1.0)/RC(12.0),RC(1.0)/RC(3.0),0,0,0}, 3754b84eec9SHong Zhang {RC(1.0)/RC(12.0),RC(1.0)/RC(6.0),RC(1.0)/RC(6.0),RC(1.0)/RC(12.0),0,RC(1.0)/RC(6.0),RC(-1.0)/RC(6.0),RC(1.0)/RC(2.0),0,0}, 3764b84eec9SHong Zhang {RC(1.0)/RC(12.0),RC(1.0)/RC(6.0),RC(1.0)/RC(6.0),RC(1.0)/RC(12.0),0,RC(1.0)/RC(12.0),RC(1.0)/RC(6.0),RC(1.0)/RC(6.0),RC(1.0)/RC(12.0),0}}, 37719c2959aSHong Zhang bsb[10] = {RC(1.0)/RC(6.0),0,0,0,RC(1.0)/RC(3.0),RC(1.0)/RC(3.0),0,0,0,RC(1.0)/RC(6.0)}, 37819c2959aSHong Zhang bf[10] = {RC(1.0)/RC(12.0),RC(1.0)/RC(6.0),RC(1.0)/RC(6.0),RC(1.0)/RC(12.0),0,RC(1.0)/RC(12.0),RC(1.0)/RC(6.0),RC(1.0)/RC(6.0),RC(1.0)/RC(12.0),0}; 37919c2959aSHong Zhang const PetscInt 38019c2959aSHong Zhang rsb[10] = {0,0,2,0,4,0,0,7,0,9}; 38179bc8a77SHong Zhang ierr = TSMPRKRegister(TSMPRKP3,3,5,2,1,&Asb[0][0],bsb,NULL,rsb,NULL,NULL,NULL,NULL,&Af[0][0],bf,NULL);CHKERRQ(ierr); 3824b84eec9SHong Zhang } 3834b84eec9SHong Zhang #undef RC 3844b84eec9SHong Zhang PetscFunctionReturn(0); 3854b84eec9SHong Zhang } 3864b84eec9SHong Zhang 3874b84eec9SHong Zhang /*@C 3884b84eec9SHong Zhang TSMPRKRegisterDestroy - Frees the list of schemes that were registered by TSMPRKRegister(). 3894b84eec9SHong Zhang 3904b84eec9SHong Zhang Not Collective 3914b84eec9SHong Zhang 3924b84eec9SHong Zhang Level: advanced 3934b84eec9SHong Zhang 3944b84eec9SHong Zhang .keywords: TSMPRK, register, destroy 3954b84eec9SHong Zhang .seealso: TSMPRKRegister(), TSMPRKRegisterAll() 3964b84eec9SHong Zhang @*/ 3974b84eec9SHong Zhang PetscErrorCode TSMPRKRegisterDestroy(void) 3984b84eec9SHong Zhang { 3994b84eec9SHong Zhang PetscErrorCode ierr; 4004b84eec9SHong Zhang MPRKTableauLink link; 4014b84eec9SHong Zhang 4024b84eec9SHong Zhang PetscFunctionBegin; 4034b84eec9SHong Zhang while ((link = MPRKTableauList)) { 4044b84eec9SHong Zhang MPRKTableau t = &link->tab; 4054b84eec9SHong Zhang MPRKTableauList = link->next; 40619c2959aSHong Zhang ierr = PetscFree3(t->Asb,t->bsb,t->csb);CHKERRQ(ierr); 40719c2959aSHong Zhang ierr = PetscFree3(t->Amb,t->bmb,t->cmb);CHKERRQ(ierr); 4084b84eec9SHong Zhang ierr = PetscFree3(t->Af,t->bf,t->cf);CHKERRQ(ierr); 40919c2959aSHong Zhang ierr = PetscFree2(t->rsb,t->rmb);CHKERRQ(ierr); 4104b84eec9SHong Zhang ierr = PetscFree (t->name);CHKERRQ(ierr); 4114b84eec9SHong Zhang ierr = PetscFree (link);CHKERRQ(ierr); 4124b84eec9SHong Zhang } 4134b84eec9SHong Zhang TSMPRKRegisterAllCalled = PETSC_FALSE; 4144b84eec9SHong Zhang PetscFunctionReturn(0); 4154b84eec9SHong Zhang } 4164b84eec9SHong Zhang 4174b84eec9SHong Zhang /*@C 4184b84eec9SHong Zhang TSMPRKInitializePackage - This function initializes everything in the TSMPRK package. It is called 4194b84eec9SHong Zhang from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_MPRK() 4204b84eec9SHong Zhang when using static libraries. 4214b84eec9SHong Zhang 4224b84eec9SHong Zhang Level: developer 4234b84eec9SHong Zhang 4244b84eec9SHong Zhang .keywords: TS, TSMPRK, initialize, package 4254b84eec9SHong Zhang .seealso: PetscInitialize() 4264b84eec9SHong Zhang @*/ 4274b84eec9SHong Zhang PetscErrorCode TSMPRKInitializePackage(void) 4284b84eec9SHong Zhang { 4294b84eec9SHong Zhang PetscErrorCode ierr; 4304b84eec9SHong Zhang 4314b84eec9SHong Zhang PetscFunctionBegin; 4324b84eec9SHong Zhang if (TSMPRKPackageInitialized) PetscFunctionReturn(0); 4334b84eec9SHong Zhang TSMPRKPackageInitialized = PETSC_TRUE; 4344b84eec9SHong Zhang ierr = TSMPRKRegisterAll();CHKERRQ(ierr); 4354b84eec9SHong Zhang ierr = PetscRegisterFinalize(TSMPRKFinalizePackage);CHKERRQ(ierr); 4364b84eec9SHong Zhang PetscFunctionReturn(0); 4374b84eec9SHong Zhang } 4384b84eec9SHong Zhang 4394b84eec9SHong Zhang /*@C 440*f944a40eSHong Zhang TSMPRKFinalizePackage - This function destroys everything in the TSMPRK package. It is 4414b84eec9SHong Zhang called from PetscFinalize(). 4424b84eec9SHong Zhang 4434b84eec9SHong Zhang Level: developer 4444b84eec9SHong Zhang 4454b84eec9SHong Zhang .keywords: Petsc, destroy, package 4464b84eec9SHong Zhang .seealso: PetscFinalize() 4474b84eec9SHong Zhang @*/ 4484b84eec9SHong Zhang PetscErrorCode TSMPRKFinalizePackage(void) 4494b84eec9SHong Zhang { 4504b84eec9SHong Zhang PetscErrorCode ierr; 4514b84eec9SHong Zhang 4524b84eec9SHong Zhang PetscFunctionBegin; 4534b84eec9SHong Zhang TSMPRKPackageInitialized = PETSC_FALSE; 4544b84eec9SHong Zhang ierr = TSMPRKRegisterDestroy();CHKERRQ(ierr); 4554b84eec9SHong Zhang PetscFunctionReturn(0); 4564b84eec9SHong Zhang } 4574b84eec9SHong Zhang 4584b84eec9SHong Zhang /*@C 4594b84eec9SHong Zhang TSMPRKRegister - register a MPRK scheme by providing the entries in the Butcher tableau 4604b84eec9SHong Zhang 4614b84eec9SHong Zhang Not Collective, but the same schemes should be registered on all processes on which they will be used 4624b84eec9SHong Zhang 4634b84eec9SHong Zhang Input Parameters: 4644b84eec9SHong Zhang + name - identifier for method 4654b84eec9SHong Zhang . order - approximation order of method 46679bc8a77SHong Zhang . s - number of stages in the base methods 46779bc8a77SHong Zhang . ratio1 - stepsize ratio at 1st level (e.g. slow/medium) 46879bc8a77SHong Zhang . ratio2 - stepsize ratio at 2nd level (e.g. medium/fast) 4694b84eec9SHong Zhang . Af - stage coefficients for fast components(dimension s*s, row-major) 4704b84eec9SHong Zhang . bf - step completion table for fast components(dimension s) 4714b84eec9SHong Zhang . cf - abscissa for fast components(dimension s) 4724b84eec9SHong Zhang . As - stage coefficients for slow components(dimension s*s, row-major) 4734b84eec9SHong Zhang . bs - step completion table for slow components(dimension s) 4744b84eec9SHong Zhang - cs - abscissa for slow components(dimension s) 4754b84eec9SHong Zhang 4764b84eec9SHong Zhang Notes: 4774b84eec9SHong Zhang Several MPRK methods are provided, this function is only needed to create new methods. 4784b84eec9SHong Zhang 4794b84eec9SHong Zhang Level: advanced 4804b84eec9SHong Zhang 4814b84eec9SHong Zhang .keywords: TS, register 4824b84eec9SHong Zhang 4834b84eec9SHong Zhang .seealso: TSMPRK 4844b84eec9SHong Zhang @*/ 48579bc8a77SHong Zhang PetscErrorCode TSMPRKRegister(TSMPRKType name,PetscInt order, 48679bc8a77SHong Zhang PetscInt sbase,PetscInt ratio1,PetscInt ratio2, 48719c2959aSHong Zhang const PetscReal Asb[],const PetscReal bsb[],const PetscReal csb[],const PetscInt rsb[], 48819c2959aSHong Zhang const PetscReal Amb[],const PetscReal bmb[],const PetscReal cmb[],const PetscInt rmb[], 4894b84eec9SHong Zhang const PetscReal Af[],const PetscReal bf[],const PetscReal cf[]) 4904b84eec9SHong Zhang { 4914b84eec9SHong Zhang MPRKTableauLink link; 4924b84eec9SHong Zhang MPRKTableau t; 49379bc8a77SHong Zhang PetscInt s,i,j; 4944b84eec9SHong Zhang PetscErrorCode ierr; 4954b84eec9SHong Zhang 4964b84eec9SHong Zhang PetscFunctionBegin; 4974b84eec9SHong Zhang PetscValidCharPointer(name,1); 49819c2959aSHong Zhang PetscValidRealPointer(Asb,4); 49979bc8a77SHong Zhang if (bsb) PetscValidRealPointer(bsb,7); 50079bc8a77SHong Zhang if (csb) PetscValidRealPointer(csb,8); 50179bc8a77SHong Zhang if (rsb) PetscValidRealPointer(rsb,9); 50279bc8a77SHong Zhang if (Amb) PetscValidRealPointer(Amb,10); 50379bc8a77SHong Zhang if (bmb) PetscValidRealPointer(bmb,11); 50479bc8a77SHong Zhang if (cmb) PetscValidRealPointer(cmb,12); 50579bc8a77SHong Zhang if (rmb) PetscValidRealPointer(rmb,13); 50679bc8a77SHong Zhang PetscValidRealPointer(Af,14); 50779bc8a77SHong Zhang if (bf) PetscValidRealPointer(bf,15); 50879bc8a77SHong Zhang if (cf) PetscValidRealPointer(cf,16); 5094b84eec9SHong Zhang 5104b84eec9SHong Zhang ierr = PetscNew(&link);CHKERRQ(ierr); 5114b84eec9SHong Zhang t = &link->tab; 5124b84eec9SHong Zhang 5134b84eec9SHong Zhang ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr); 51479bc8a77SHong Zhang s = sbase*ratio1*ratio2; /* this is the dimension of the matrices below */ 5154b84eec9SHong Zhang t->order = order; 51679bc8a77SHong Zhang t->sbase = sbase; 5174b84eec9SHong Zhang t->s = s; 51819c2959aSHong Zhang t->np = 2; 51979bc8a77SHong Zhang 5204b84eec9SHong Zhang ierr = PetscMalloc3(s*s,&t->Af,s,&t->bf,s,&t->cf);CHKERRQ(ierr); 5214b84eec9SHong Zhang ierr = PetscMemcpy(t->Af,Af,s*s*sizeof(Af[0]));CHKERRQ(ierr); 5224b84eec9SHong Zhang if (bf) { 5234b84eec9SHong Zhang ierr = PetscMemcpy(t->bf,bf,s*sizeof(bf[0]));CHKERRQ(ierr); 52419c2959aSHong Zhang } else 5254b84eec9SHong Zhang for (i=0; i<s; i++) t->bf[i] = Af[(s-1)*s+i]; 5264b84eec9SHong Zhang if (cf) { 5274b84eec9SHong Zhang ierr = PetscMemcpy(t->cf,cf,s*sizeof(cf[0]));CHKERRQ(ierr); 52819c2959aSHong Zhang } else { 5294b84eec9SHong Zhang for (i=0; i<s; i++) 5304b84eec9SHong Zhang for (j=0,t->cf[i]=0; j<s; j++) 5314b84eec9SHong Zhang t->cf[i] += Af[i*s+j]; 5324b84eec9SHong Zhang } 53319c2959aSHong Zhang 53419c2959aSHong Zhang if (Amb) { 53519c2959aSHong Zhang t->np = 3; 53619c2959aSHong Zhang ierr = PetscMalloc3(s*s,&t->Amb,s,&t->bmb,s,&t->cmb);CHKERRQ(ierr); 53719c2959aSHong Zhang ierr = PetscCalloc1(s,&t->rmb);CHKERRQ(ierr); 53819c2959aSHong Zhang ierr = PetscMemcpy(t->Amb,Amb,s*s*sizeof(Amb[0]));CHKERRQ(ierr); 53919c2959aSHong Zhang if (bmb) { 54019c2959aSHong Zhang ierr = PetscMemcpy(t->bmb,bmb,s*sizeof(bmb[0]));CHKERRQ(ierr); 54119c2959aSHong Zhang } else { 54219c2959aSHong Zhang for (i=0; i<s; i++) t->bmb[i] = Amb[(s-1)*s+i]; 5434b84eec9SHong Zhang } 54419c2959aSHong Zhang if (cmb) { 54519c2959aSHong Zhang ierr = PetscMemcpy(t->cmb,cmb,s*sizeof(cmb[0]));CHKERRQ(ierr); 54619c2959aSHong Zhang } else { 5474b84eec9SHong Zhang for (i=0; i<s; i++) 54819c2959aSHong Zhang for (j=0,t->cmb[i]=0; j<s; j++) 54919c2959aSHong Zhang t->cmb[i] += Amb[i*s+j]; 55019c2959aSHong Zhang } 55119c2959aSHong Zhang if (rmb) { 55219c2959aSHong Zhang ierr = PetscMemcpy(t->rmb,rmb,s*sizeof(rmb[0]));CHKERRQ(ierr); 55319c2959aSHong Zhang } 55419c2959aSHong Zhang } 55519c2959aSHong Zhang 55619c2959aSHong Zhang ierr = PetscMalloc3(s*s,&t->Asb,s,&t->bsb,s,&t->csb);CHKERRQ(ierr); 55719c2959aSHong Zhang ierr = PetscMemcpy(t->Asb,Asb,s*s*sizeof(Asb[0]));CHKERRQ(ierr); 55819c2959aSHong Zhang ierr = PetscCalloc1(s,&t->rsb);CHKERRQ(ierr); 55919c2959aSHong Zhang if (bsb) { 56019c2959aSHong Zhang ierr = PetscMemcpy(t->bsb,bsb,s*sizeof(bsb[0]));CHKERRQ(ierr); 56119c2959aSHong Zhang } else 56219c2959aSHong Zhang for (i=0; i<s; i++) t->bsb[i] = Asb[(s-1)*s+i]; 56319c2959aSHong Zhang if (csb) { 56419c2959aSHong Zhang ierr = PetscMemcpy(t->csb,csb,s*sizeof(csb[0]));CHKERRQ(ierr); 56519c2959aSHong Zhang } else { 56619c2959aSHong Zhang for (i=0; i<s; i++) 56719c2959aSHong Zhang for (j=0,t->csb[i]=0; j<s; j++) 56819c2959aSHong Zhang t->csb[i] += Asb[i*s+j]; 56919c2959aSHong Zhang } 57019c2959aSHong Zhang if (rsb) { 57119c2959aSHong Zhang ierr = PetscMemcpy(t->rsb,rsb,s*sizeof(rsb[0]));CHKERRQ(ierr); 5724b84eec9SHong Zhang } 5734b84eec9SHong Zhang link->next = MPRKTableauList; 5744b84eec9SHong Zhang MPRKTableauList = link; 5754b84eec9SHong Zhang PetscFunctionReturn(0); 5764b84eec9SHong Zhang } 5774b84eec9SHong Zhang 5784b84eec9SHong Zhang static PetscErrorCode TSMPRKSetSplits(TS ts) 5794b84eec9SHong Zhang { 5804b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 58119c2959aSHong Zhang MPRKTableau tab = mprk->tableau; 5824b84eec9SHong Zhang DM dm,subdm,newdm; 5834b84eec9SHong Zhang PetscErrorCode ierr; 5844b84eec9SHong Zhang 5854b84eec9SHong Zhang PetscFunctionBegin; 5864b84eec9SHong Zhang ierr = TSRHSSplitGetSubTS(ts,"slow",&mprk->subts_slow);CHKERRQ(ierr); 5874b84eec9SHong Zhang ierr = TSRHSSplitGetSubTS(ts,"fast",&mprk->subts_fast);CHKERRQ(ierr); 5884b84eec9SHong Zhang if (!mprk->subts_slow || !mprk->subts_fast) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must set up the RHSFunctions for 'slow' and 'fast' components using TSRHSSplitSetRHSFunction() or calling TSSetRHSFunction() for each sub-TS"); 5894b84eec9SHong Zhang 59019c2959aSHong Zhang /* Only copy the DM */ 5914b84eec9SHong Zhang ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 59219c2959aSHong Zhang 59319c2959aSHong Zhang ierr = TSRHSSplitGetSubTS(ts,"slowbuffer",&mprk->subts_slowbuffer);CHKERRQ(ierr); 59419c2959aSHong Zhang if (!mprk->subts_slowbuffer) { 59519c2959aSHong Zhang mprk->subts_slowbuffer = mprk->subts_slow; 59619c2959aSHong Zhang mprk->subts_slow = NULL; 59719c2959aSHong Zhang } 59819c2959aSHong Zhang if (mprk->subts_slow) { 5994b84eec9SHong Zhang ierr = DMClone(dm,&newdm);CHKERRQ(ierr); 6004b84eec9SHong Zhang ierr = TSGetDM(mprk->subts_slow,&subdm);CHKERRQ(ierr); 6014b84eec9SHong Zhang ierr = DMCopyDMTS(subdm,newdm);CHKERRQ(ierr); 6024b84eec9SHong Zhang ierr = DMCopyDMSNES(subdm,newdm);CHKERRQ(ierr); 6034b84eec9SHong Zhang ierr = TSSetDM(mprk->subts_slow,newdm);CHKERRQ(ierr); 6044b84eec9SHong Zhang ierr = DMDestroy(&newdm);CHKERRQ(ierr); 60519c2959aSHong Zhang } 60619c2959aSHong Zhang ierr = DMClone(dm,&newdm);CHKERRQ(ierr); 60719c2959aSHong Zhang ierr = TSGetDM(mprk->subts_slowbuffer,&subdm);CHKERRQ(ierr); 60819c2959aSHong Zhang ierr = DMCopyDMTS(subdm,newdm);CHKERRQ(ierr); 60919c2959aSHong Zhang ierr = DMCopyDMSNES(subdm,newdm);CHKERRQ(ierr); 61019c2959aSHong Zhang ierr = TSSetDM(mprk->subts_slowbuffer,newdm);CHKERRQ(ierr); 61119c2959aSHong Zhang ierr = DMDestroy(&newdm);CHKERRQ(ierr); 61219c2959aSHong Zhang 61319c2959aSHong Zhang ierr = DMClone(dm,&newdm);CHKERRQ(ierr); 61419c2959aSHong Zhang ierr = TSGetDM(mprk->subts_fast,&subdm);CHKERRQ(ierr); 61519c2959aSHong Zhang ierr = DMCopyDMTS(subdm,newdm);CHKERRQ(ierr); 61619c2959aSHong Zhang ierr = DMCopyDMSNES(subdm,newdm);CHKERRQ(ierr); 61719c2959aSHong Zhang ierr = TSSetDM(mprk->subts_fast,newdm);CHKERRQ(ierr); 61819c2959aSHong Zhang ierr = DMDestroy(&newdm);CHKERRQ(ierr); 61919c2959aSHong Zhang 62019c2959aSHong Zhang if (tab->np == 3) { 62119c2959aSHong Zhang ierr = TSRHSSplitGetSubTS(ts,"medium",&mprk->subts_medium);CHKERRQ(ierr); 62219c2959aSHong Zhang ierr = TSRHSSplitGetSubTS(ts,"mediumbuffer",&mprk->subts_mediumbuffer);CHKERRQ(ierr); 62319c2959aSHong Zhang if (mprk->subts_medium && !mprk->subts_mediumbuffer) { 62419c2959aSHong Zhang mprk->subts_mediumbuffer = mprk->subts_medium; 62519c2959aSHong Zhang mprk->subts_medium = NULL; 62619c2959aSHong Zhang } 62719c2959aSHong Zhang if (mprk->subts_medium) { 62819c2959aSHong Zhang ierr = DMClone(dm,&newdm);CHKERRQ(ierr); 62919c2959aSHong Zhang ierr = TSGetDM(mprk->subts_medium,&subdm);CHKERRQ(ierr); 63019c2959aSHong Zhang ierr = DMCopyDMTS(subdm,newdm);CHKERRQ(ierr); 63119c2959aSHong Zhang ierr = DMCopyDMSNES(subdm,newdm);CHKERRQ(ierr); 63219c2959aSHong Zhang ierr = TSSetDM(mprk->subts_medium,newdm);CHKERRQ(ierr); 63319c2959aSHong Zhang ierr = DMDestroy(&newdm);CHKERRQ(ierr); 63419c2959aSHong Zhang } 63519c2959aSHong Zhang ierr = DMClone(dm,&newdm);CHKERRQ(ierr); 63619c2959aSHong Zhang ierr = TSGetDM(mprk->subts_mediumbuffer,&subdm);CHKERRQ(ierr); 63719c2959aSHong Zhang ierr = DMCopyDMTS(subdm,newdm);CHKERRQ(ierr); 63819c2959aSHong Zhang ierr = DMCopyDMSNES(subdm,newdm);CHKERRQ(ierr); 63919c2959aSHong Zhang ierr = TSSetDM(mprk->subts_mediumbuffer,newdm);CHKERRQ(ierr); 64019c2959aSHong Zhang ierr = DMDestroy(&newdm);CHKERRQ(ierr); 64119c2959aSHong Zhang } 6424b84eec9SHong Zhang PetscFunctionReturn(0); 6434b84eec9SHong Zhang } 6444b84eec9SHong Zhang 6454b84eec9SHong Zhang /* 6464b84eec9SHong Zhang This if for nonsplit RHS MPRK 6474b84eec9SHong Zhang The step completion formula is 6484b84eec9SHong Zhang 6494b84eec9SHong Zhang x1 = x0 + h b^T YdotRHS 6504b84eec9SHong Zhang 6514b84eec9SHong Zhang */ 6524b84eec9SHong Zhang static PetscErrorCode TSEvaluateStep_MPRK(TS ts,PetscInt order,Vec X,PetscBool *done) 6534b84eec9SHong Zhang { 6544b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 6554b84eec9SHong Zhang MPRKTableau tab = mprk->tableau; 6567c0df07dSHong Zhang PetscScalar *wf = mprk->work_fast; 6574b84eec9SHong Zhang PetscReal h = ts->time_step; 6584b84eec9SHong Zhang PetscInt s = tab->s,j; 6594b84eec9SHong Zhang PetscErrorCode ierr; 6604b84eec9SHong Zhang 6614b84eec9SHong Zhang PetscFunctionBegin; 6624b84eec9SHong Zhang for (j=0; j<s; j++) wf[j] = h*tab->bf[j]; 6637c0df07dSHong Zhang ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 6647c0df07dSHong Zhang ierr = VecMAXPY(X,s,wf,mprk->YdotRHS);CHKERRQ(ierr); 6654b84eec9SHong Zhang PetscFunctionReturn(0); 6664b84eec9SHong Zhang } 6674b84eec9SHong Zhang 6684b84eec9SHong Zhang static PetscErrorCode TSStep_MPRK(TS ts) 6694b84eec9SHong Zhang { 6704b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 6717c0df07dSHong Zhang Vec *Y = mprk->Y,*YdotRHS = mprk->YdotRHS,*YdotRHS_fast = mprk->YdotRHS_fast,*YdotRHS_slow = mprk->YdotRHS_slow,*YdotRHS_slowbuffer = mprk->YdotRHS_slowbuffer; 6727c0df07dSHong Zhang Vec Yslow,Yslowbuffer,Yfast; 6734b84eec9SHong Zhang MPRKTableau tab = mprk->tableau; 6744b84eec9SHong Zhang const PetscInt s = tab->s; 67519c2959aSHong Zhang const PetscReal *Af = tab->Af,*cf = tab->cf,*Asb = tab->Asb,*csb = tab->csb; 676ebd5ed4eSHong Zhang PetscScalar *wf = mprk->work_fast,*wsb = mprk->work_slowbuffer; 677ebd5ed4eSHong Zhang PetscInt i,j; 6784b84eec9SHong Zhang PetscReal next_time_step = ts->time_step,t = ts->ptime,h = ts->time_step; 6794b84eec9SHong Zhang PetscErrorCode ierr; 6804b84eec9SHong Zhang 6814b84eec9SHong Zhang PetscFunctionBegin; 6824b84eec9SHong Zhang for (i=0; i<s; i++) { 6834b84eec9SHong Zhang mprk->stage_time = t + h*cf[i]; 6844b84eec9SHong Zhang ierr = TSPreStage(ts,mprk->stage_time);CHKERRQ(ierr); 6857c0df07dSHong Zhang ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr); 6864b84eec9SHong Zhang 6877c0df07dSHong Zhang /* slow buffer region */ 68819c2959aSHong Zhang for (j=0; j<i; j++) wsb[j] = h*Asb[i*s+j]; 6899849be05SHong Zhang for (j=0; j<i; j++) { 6907c0df07dSHong Zhang ierr = VecGetSubVector(YdotRHS[j],mprk->is_slowbuffer,&YdotRHS_slowbuffer[j]);CHKERRQ(ierr); 6917c0df07dSHong Zhang } 6927c0df07dSHong Zhang ierr = VecGetSubVector(Y[i],mprk->is_slowbuffer,&Yslowbuffer);CHKERRQ(ierr); 6937c0df07dSHong Zhang ierr = VecMAXPY(Yslowbuffer,i,wsb,mprk->YdotRHS_slowbuffer);CHKERRQ(ierr); 6947c0df07dSHong Zhang ierr = VecRestoreSubVector(Y[i],mprk->is_slowbuffer,&Yslowbuffer);CHKERRQ(ierr); 6959849be05SHong Zhang for (j=0; j<i; j++) { 6967c0df07dSHong Zhang ierr = VecRestoreSubVector(YdotRHS[j],mprk->is_slowbuffer,&YdotRHS_slowbuffer[j]);CHKERRQ(ierr); 6977c0df07dSHong Zhang } 6987c0df07dSHong Zhang /* slow region */ 6997c0df07dSHong Zhang if (mprk->is_slow) { 7009849be05SHong Zhang for (j=0; j<i; j++) { 7017c0df07dSHong Zhang ierr = VecGetSubVector(YdotRHS[j],mprk->is_slow,&YdotRHS_slow[j]);CHKERRQ(ierr); 7027c0df07dSHong Zhang } 7037c0df07dSHong Zhang ierr = VecGetSubVector(Y[i],mprk->is_slow,&Yslow);CHKERRQ(ierr); 704ebd5ed4eSHong Zhang ierr = VecMAXPY(Yslow,i,wsb,mprk->YdotRHS_slow);CHKERRQ(ierr); 7057c0df07dSHong Zhang ierr = VecRestoreSubVector(Y[i],mprk->is_slow,&Yslow);CHKERRQ(ierr); 706ebd5ed4eSHong Zhang for (j=0; j<i; j++) { 7077c0df07dSHong Zhang ierr = VecRestoreSubVector(YdotRHS[j],mprk->is_slow,&YdotRHS_slow[j]);CHKERRQ(ierr); 7087c0df07dSHong Zhang } 7097c0df07dSHong Zhang } 7104b84eec9SHong Zhang 7117c0df07dSHong Zhang /* fast region */ 7124b84eec9SHong Zhang for (j=0; j<i; j++) wf[j] = h*Af[i*s+j]; 7139849be05SHong Zhang for (j=0; j<i; j++) { 7147c0df07dSHong Zhang ierr = VecGetSubVector(YdotRHS[j],mprk->is_fast,&YdotRHS_fast[j]);CHKERRQ(ierr); 7157c0df07dSHong Zhang } 7167c0df07dSHong Zhang ierr = VecGetSubVector(Y[i],mprk->is_fast,&Yfast);CHKERRQ(ierr); 7177c0df07dSHong Zhang ierr = VecMAXPY(Yfast,i,wf,mprk->YdotRHS_fast);CHKERRQ(ierr); 7187c0df07dSHong Zhang ierr = VecRestoreSubVector(Y[i],mprk->is_fast,&Yfast);CHKERRQ(ierr); 7199849be05SHong Zhang for (j=0; j<i; j++) { 7207c0df07dSHong Zhang ierr = VecRestoreSubVector(YdotRHS[j],mprk->is_fast,&YdotRHS_fast[j]);CHKERRQ(ierr); 7217c0df07dSHong Zhang } 7227c0df07dSHong Zhang if (tab->np == 3) { 7237c0df07dSHong Zhang Vec *YdotRHS_medium = mprk->YdotRHS_medium,*YdotRHS_mediumbuffer = mprk->YdotRHS_mediumbuffer; 7247c0df07dSHong Zhang Vec Ymedium,Ymediumbuffer; 725ebd5ed4eSHong Zhang PetscScalar *wmb = mprk->work_mediumbuffer; 726ebd5ed4eSHong Zhang 727ebd5ed4eSHong Zhang for (j=0; j<i; j++) wmb[j] = h*tab->Amb[i*s+j]; 7287c0df07dSHong Zhang /* medium region */ 7297c0df07dSHong Zhang if (mprk->is_medium) { 7307c0df07dSHong Zhang for (j=0; j<i; j++) { 7317c0df07dSHong Zhang ierr = VecGetSubVector(YdotRHS[j],mprk->is_medium,&YdotRHS_medium[j]);CHKERRQ(ierr); 7327c0df07dSHong Zhang } 7337c0df07dSHong Zhang ierr = VecGetSubVector(Y[i],mprk->is_medium,&Ymedium);CHKERRQ(ierr); 734ebd5ed4eSHong Zhang ierr = VecMAXPY(Ymedium,i,wmb,mprk->YdotRHS_medium);CHKERRQ(ierr); 7357c0df07dSHong Zhang ierr = VecRestoreSubVector(Y[i],mprk->is_medium,&Ymedium);CHKERRQ(ierr); 736ebd5ed4eSHong Zhang for (j=0; j<i; j++) { 7377c0df07dSHong Zhang ierr = VecRestoreSubVector(YdotRHS[j],mprk->is_medium,&YdotRHS_medium[j]);CHKERRQ(ierr); 7387c0df07dSHong Zhang } 7397c0df07dSHong Zhang } 7407c0df07dSHong Zhang /* medium buffer region */ 7419849be05SHong Zhang for (j=0; j<i; j++) { 7427c0df07dSHong Zhang ierr = VecGetSubVector(YdotRHS[j],mprk->is_mediumbuffer,&YdotRHS_mediumbuffer[j]);CHKERRQ(ierr); 7437c0df07dSHong Zhang } 7447c0df07dSHong Zhang ierr = VecGetSubVector(Y[i],mprk->is_mediumbuffer,&Ymediumbuffer);CHKERRQ(ierr); 7457c0df07dSHong Zhang ierr = VecMAXPY(Ymediumbuffer,i,wmb,mprk->YdotRHS_mediumbuffer);CHKERRQ(ierr); 7467c0df07dSHong Zhang ierr = VecRestoreSubVector(Y[i],mprk->is_mediumbuffer,&Ymediumbuffer);CHKERRQ(ierr); 7479849be05SHong Zhang for (j=0; j<i; j++) { 7487c0df07dSHong Zhang ierr = VecRestoreSubVector(YdotRHS[j],mprk->is_mediumbuffer,&YdotRHS_mediumbuffer[j]);CHKERRQ(ierr); 7497c0df07dSHong Zhang } 7507c0df07dSHong Zhang } 7514b84eec9SHong Zhang ierr = TSPostStage(ts,mprk->stage_time,i,Y); CHKERRQ(ierr); 7524b84eec9SHong Zhang /* compute the stage RHS by fast and slow tableau respectively */ 7537c0df07dSHong Zhang ierr = TSComputeRHSFunction(ts,t+h*csb[i],Y[i],YdotRHS[i]);CHKERRQ(ierr); 7544b84eec9SHong Zhang } 7554b84eec9SHong Zhang ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr); 7564b84eec9SHong Zhang ts->ptime += ts->time_step; 7574b84eec9SHong Zhang ts->time_step = next_time_step; 7584b84eec9SHong Zhang PetscFunctionReturn(0); 7594b84eec9SHong Zhang } 7604b84eec9SHong Zhang 7614b84eec9SHong Zhang /* 762*f944a40eSHong Zhang This if for the case when split RHS is used 7634b84eec9SHong Zhang The step completion formula is 7644b84eec9SHong Zhang x1 = x0 + h b^T YdotRHS 7654b84eec9SHong Zhang */ 7664b84eec9SHong Zhang static PetscErrorCode TSEvaluateStep_MPRKSPLIT(TS ts,PetscInt order,Vec X,PetscBool *done) 7674b84eec9SHong Zhang { 7684b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 7694b84eec9SHong Zhang MPRKTableau tab = mprk->tableau; 77019c2959aSHong Zhang Vec Xslow,Xfast,Xslowbuffer; /* subvectors for slow and fast componets in X respectively */ 77119c2959aSHong Zhang PetscScalar *wf = mprk->work_fast,*ws = mprk->work_slow,*wsb = mprk->work_slowbuffer; 7724b84eec9SHong Zhang PetscReal h = ts->time_step; 77379bc8a77SHong Zhang PetscInt s = tab->s,j,computedstages; 7744b84eec9SHong Zhang PetscErrorCode ierr; 7754b84eec9SHong Zhang 7764b84eec9SHong Zhang PetscFunctionBegin; 7774b84eec9SHong Zhang ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 77819c2959aSHong Zhang 77919c2959aSHong Zhang /* slow region */ 78019c2959aSHong Zhang if (mprk->is_slow) { 78179bc8a77SHong Zhang computedstages = 0; 78219c2959aSHong Zhang for (j=0; j<s; j++) { 7839849be05SHong Zhang if (tab->rsb[j]) ws[tab->rsb[j]-1] += h*tab->bsb[j]; 78479bc8a77SHong Zhang else ws[computedstages++] = h*tab->bsb[j]; 78519c2959aSHong Zhang } 7864b84eec9SHong Zhang ierr = VecGetSubVector(X,mprk->is_slow,&Xslow);CHKERRQ(ierr); 78779bc8a77SHong Zhang ierr = VecMAXPY(Xslow,computedstages,ws,mprk->YdotRHS_slow);CHKERRQ(ierr); 78819c2959aSHong Zhang ierr = VecRestoreSubVector(X,mprk->is_slow,&Xslow);CHKERRQ(ierr); 78919c2959aSHong Zhang } 79019c2959aSHong Zhang 7919d6e09e9SHong Zhang if (tab->np == 3 && mprk->is_medium) { 7929d6e09e9SHong Zhang computedstages = 0; 7939d6e09e9SHong Zhang for (j=0; j<s; j++) { 7949d6e09e9SHong Zhang if (tab->rmb[j]) wsb[computedstages-tab->sbase+(tab->rmb[j]-1)%tab->sbase] += h*tab->bsb[j]; 7959d6e09e9SHong Zhang else wsb[computedstages++] = h*tab->bsb[j]; 7969d6e09e9SHong Zhang } 7979d6e09e9SHong Zhang ierr = VecGetSubVector(X,mprk->is_slowbuffer,&Xslowbuffer);CHKERRQ(ierr); 7989d6e09e9SHong Zhang ierr = VecMAXPY(Xslowbuffer,computedstages,wsb,mprk->YdotRHS_slowbuffer);CHKERRQ(ierr); 7999d6e09e9SHong Zhang ierr = VecRestoreSubVector(X,mprk->is_slowbuffer,&Xslowbuffer);CHKERRQ(ierr); 8009d6e09e9SHong Zhang } else { 80119c2959aSHong Zhang /* slow buffer region */ 80219c2959aSHong Zhang for (j=0; j<s; j++) wsb[j] = h*tab->bsb[j]; 80319c2959aSHong Zhang ierr = VecGetSubVector(X,mprk->is_slowbuffer,&Xslowbuffer);CHKERRQ(ierr); 80419c2959aSHong Zhang ierr = VecMAXPY(Xslowbuffer,s,wsb,mprk->YdotRHS_slowbuffer);CHKERRQ(ierr); 80519c2959aSHong Zhang ierr = VecRestoreSubVector(X,mprk->is_slowbuffer,&Xslowbuffer);CHKERRQ(ierr); 8069d6e09e9SHong Zhang } 80719c2959aSHong Zhang if (tab->np == 3) { 80819c2959aSHong Zhang Vec Xmedium,Xmediumbuffer; 80919c2959aSHong Zhang PetscScalar *wm = mprk->work_medium,*wmb = mprk->work_mediumbuffer; 8109d6e09e9SHong Zhang /* medium region and slow buffer region */ 81119c2959aSHong Zhang if (mprk->is_medium) { 81279bc8a77SHong Zhang computedstages = 0; 81319c2959aSHong Zhang for (j=0; j<s; j++) { 81479bc8a77SHong Zhang if (tab->rmb[j]) wm[computedstages-tab->sbase+(tab->rmb[j]-1)%tab->sbase] += h*tab->bmb[j]; 81579bc8a77SHong Zhang else wm[computedstages++] = h*tab->bmb[j]; 81619c2959aSHong Zhang } 81719c2959aSHong Zhang ierr = VecGetSubVector(X,mprk->is_medium,&Xmedium);CHKERRQ(ierr); 81879bc8a77SHong Zhang ierr = VecMAXPY(Xmedium,computedstages,wm,mprk->YdotRHS_medium);CHKERRQ(ierr); 81919c2959aSHong Zhang ierr = VecRestoreSubVector(X,mprk->is_medium,&Xmedium);CHKERRQ(ierr); 82019c2959aSHong Zhang } 82119c2959aSHong Zhang /* medium buffer region */ 82219c2959aSHong Zhang for (j=0; j<s; j++) wmb[j] = h*tab->bmb[j]; 82319c2959aSHong Zhang ierr = VecGetSubVector(X,mprk->is_mediumbuffer,&Xmediumbuffer);CHKERRQ(ierr); 82419c2959aSHong Zhang ierr = VecMAXPY(Xmediumbuffer,s,wmb,mprk->YdotRHS_mediumbuffer);CHKERRQ(ierr); 82519c2959aSHong Zhang ierr = VecRestoreSubVector(X,mprk->is_mediumbuffer,&Xmediumbuffer);CHKERRQ(ierr); 82619c2959aSHong Zhang } 82719c2959aSHong Zhang /* fast region */ 82819c2959aSHong Zhang for (j=0; j<s; j++) wf[j] = h*tab->bf[j]; 8294b84eec9SHong Zhang ierr = VecGetSubVector(X,mprk->is_fast,&Xfast);CHKERRQ(ierr); 8304b84eec9SHong Zhang ierr = VecMAXPY(Xfast,s,wf,mprk->YdotRHS_fast);CHKERRQ(ierr); 83119c2959aSHong Zhang ierr = VecRestoreSubVector(X,mprk->is_fast,&Xfast);CHKERRQ(ierr); 8324b84eec9SHong Zhang PetscFunctionReturn(0); 8334b84eec9SHong Zhang } 8344b84eec9SHong Zhang 8354b84eec9SHong Zhang static PetscErrorCode TSStep_MPRKSPLIT(TS ts) 8364b84eec9SHong Zhang { 8374b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 8384b84eec9SHong Zhang MPRKTableau tab = mprk->tableau; 83919c2959aSHong Zhang Vec *Y = mprk->Y,*YdotRHS_fast = mprk->YdotRHS_fast,*YdotRHS_slow = mprk->YdotRHS_slow,*YdotRHS_slowbuffer = mprk->YdotRHS_slowbuffer; 84019c2959aSHong Zhang Vec Yslow,Yslowbuffer,Yfast; /* subvectors for slow and fast components in Y[i] respectively */ 84119c2959aSHong Zhang PetscInt s = tab->s; 84219c2959aSHong Zhang const PetscReal *Af = tab->Af,*cf = tab->cf,*Asb = tab->Asb,*csb = tab->csb; 84319c2959aSHong Zhang PetscScalar *wf = mprk->work_fast,*ws = mprk->work_slow,*wsb = mprk->work_slowbuffer; 84479bc8a77SHong Zhang PetscInt i,j,computedstages; 8454b84eec9SHong Zhang PetscReal next_time_step = ts->time_step,t = ts->ptime,h = ts->time_step; 8464b84eec9SHong Zhang PetscErrorCode ierr; 8474b84eec9SHong Zhang 8484b84eec9SHong Zhang PetscFunctionBegin; 8494b84eec9SHong Zhang for (i=0; i<s; i++) { 8504b84eec9SHong Zhang mprk->stage_time = t + h*cf[i]; 8514b84eec9SHong Zhang ierr = TSPreStage(ts,mprk->stage_time);CHKERRQ(ierr); 8524b84eec9SHong Zhang /* calculate the stage value for fast and slow components respectively */ 8534b84eec9SHong Zhang ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr); 85419c2959aSHong Zhang for (j=0; j<i; j++) wsb[j] = h*Asb[i*s+j]; 85519c2959aSHong Zhang 85619c2959aSHong Zhang /* slow buffer region */ 8579d6e09e9SHong Zhang if (tab->np == 3 && mprk->is_medium) { 8589d6e09e9SHong Zhang if (tab->rmb[i]) { 8599d6e09e9SHong Zhang ierr = VecGetSubVector(Y[i],mprk->is_slowbuffer,&Yslowbuffer);CHKERRQ(ierr); 8609d6e09e9SHong Zhang ierr = VecISCopy(Y[tab->rmb[i]-1],mprk->is_slowbuffer,SCATTER_REVERSE,Yslowbuffer);CHKERRQ(ierr); 8619d6e09e9SHong Zhang ierr = VecRestoreSubVector(Y[i],mprk->is_slowbuffer,&Yslowbuffer);CHKERRQ(ierr); 8629d6e09e9SHong Zhang } else { 8639d6e09e9SHong Zhang PetscScalar *wm = mprk->work_medium; 8649d6e09e9SHong Zhang computedstages = 0; 8659d6e09e9SHong Zhang for (j=0; j<i; j++) { 8669d6e09e9SHong Zhang if (tab->rmb[j]) wm[computedstages-tab->sbase+(tab->rmb[j]-1)%tab->sbase] += wsb[j]; 8679d6e09e9SHong Zhang else wm[computedstages++] = wsb[j]; 8689d6e09e9SHong Zhang } 8699d6e09e9SHong Zhang ierr = VecGetSubVector(Y[i],mprk->is_slowbuffer,&Yslowbuffer);CHKERRQ(ierr); 8709d6e09e9SHong Zhang ierr = VecMAXPY(Yslowbuffer,computedstages,wm,YdotRHS_slowbuffer);CHKERRQ(ierr); 8719d6e09e9SHong Zhang ierr = VecRestoreSubVector(Y[i],mprk->is_slowbuffer,&Yslowbuffer);CHKERRQ(ierr); 8729d6e09e9SHong Zhang } 8739d6e09e9SHong Zhang } else { 87419c2959aSHong Zhang ierr = VecGetSubVector(Y[i],mprk->is_slowbuffer,&Yslowbuffer);CHKERRQ(ierr); 87519c2959aSHong Zhang ierr = VecMAXPY(Yslowbuffer,i,wsb,YdotRHS_slowbuffer);CHKERRQ(ierr); 87619c2959aSHong Zhang ierr = VecRestoreSubVector(Y[i],mprk->is_slowbuffer,&Yslowbuffer);CHKERRQ(ierr); 8779d6e09e9SHong Zhang } 8789849be05SHong Zhang 87919c2959aSHong Zhang /* slow region */ 8809849be05SHong Zhang if (mprk->is_slow) { 8819849be05SHong Zhang if (tab->rsb[i]) { /* repeat previous stage */ 8829849be05SHong Zhang ierr = VecGetSubVector(Y[i],mprk->is_slow,&Yslow);CHKERRQ(ierr); 8839849be05SHong Zhang ierr = VecISCopy(Y[tab->rsb[i]-1],mprk->is_slow,SCATTER_REVERSE,Yslow);CHKERRQ(ierr); 8849849be05SHong Zhang ierr = VecRestoreSubVector(Y[i],mprk->is_slow,&Yslow);CHKERRQ(ierr); 8859849be05SHong Zhang } else { 88679bc8a77SHong Zhang computedstages = 0; 88719c2959aSHong Zhang for (j=0; j<i; j++) { 88879bc8a77SHong Zhang if (tab->rsb[j]) ws[tab->rsb[j]-1] += wsb[j]; 88979bc8a77SHong Zhang else ws[computedstages++] = wsb[j]; 89019c2959aSHong Zhang } 8914b84eec9SHong Zhang ierr = VecGetSubVector(Y[i],mprk->is_slow,&Yslow);CHKERRQ(ierr); 89279bc8a77SHong Zhang ierr = VecMAXPY(Yslow,computedstages,ws,YdotRHS_slow);CHKERRQ(ierr); 8934b84eec9SHong Zhang ierr = VecRestoreSubVector(Y[i],mprk->is_slow,&Yslow);CHKERRQ(ierr); 89419c2959aSHong Zhang /* only depends on the slow buffer region */ 89579bc8a77SHong Zhang ierr = TSComputeRHSFunction(mprk->subts_slow,t+h*csb[i],Y[i],YdotRHS_slow[computedstages]);CHKERRQ(ierr); 89619c2959aSHong Zhang } 8979849be05SHong Zhang } 89819c2959aSHong Zhang 89919c2959aSHong Zhang /* fast region */ 90019c2959aSHong Zhang for (j=0; j<i; j++) wf[j] = h*Af[i*s+j]; 90119c2959aSHong Zhang ierr = VecGetSubVector(Y[i],mprk->is_fast,&Yfast);CHKERRQ(ierr); 90219c2959aSHong Zhang ierr = VecMAXPY(Yfast,i,wf,YdotRHS_fast);CHKERRQ(ierr); 9034b84eec9SHong Zhang ierr = VecRestoreSubVector(Y[i],mprk->is_fast,&Yfast);CHKERRQ(ierr); 90419c2959aSHong Zhang 90519c2959aSHong Zhang if (tab->np == 3) { 90619c2959aSHong Zhang Vec *YdotRHS_medium = mprk->YdotRHS_medium,*YdotRHS_mediumbuffer = mprk->YdotRHS_mediumbuffer; 90719c2959aSHong Zhang Vec Ymedium,Ymediumbuffer; 90819c2959aSHong Zhang const PetscReal *Amb = tab->Amb,*cmb = tab->cmb; 90919c2959aSHong Zhang PetscScalar *wm = mprk->work_medium,*wmb = mprk->work_mediumbuffer; 91019c2959aSHong Zhang 91119c2959aSHong Zhang for (j=0; j<i; j++) wmb[j] = h*Amb[i*s+j]; 91219c2959aSHong Zhang /* medium buffer region */ 91319c2959aSHong Zhang ierr = VecGetSubVector(Y[i],mprk->is_mediumbuffer,&Ymediumbuffer);CHKERRQ(ierr); 91419c2959aSHong Zhang ierr = VecMAXPY(Ymediumbuffer,i,wmb,YdotRHS_mediumbuffer);CHKERRQ(ierr); 91519c2959aSHong Zhang ierr = VecRestoreSubVector(Y[i],mprk->is_mediumbuffer,&Ymediumbuffer);CHKERRQ(ierr); 91619c2959aSHong Zhang /* medium region */ 91779bc8a77SHong Zhang if (mprk->is_medium) { 91879bc8a77SHong Zhang if (tab->rmb[i]) { /* repeat previous stage */ 91979bc8a77SHong Zhang ierr = VecGetSubVector(Y[i],mprk->is_medium,&Ymedium);CHKERRQ(ierr); 92079bc8a77SHong Zhang ierr = VecISCopy(Y[tab->rmb[i]-1],mprk->is_medium,SCATTER_REVERSE,Ymedium);CHKERRQ(ierr); 92119c2959aSHong Zhang ierr = VecRestoreSubVector(Y[i],mprk->is_medium,&Ymedium);CHKERRQ(ierr); 92279bc8a77SHong Zhang } else { 92379bc8a77SHong Zhang computedstages = 0; 92479bc8a77SHong Zhang for (j=0; j<i; j++) { 92579bc8a77SHong Zhang if (tab->rmb[j]) wm[computedstages-tab->sbase+(tab->rmb[j]-1)%tab->sbase] += wmb[j]; 92679bc8a77SHong Zhang else wm[computedstages++] = wmb[j]; 92779bc8a77SHong Zhang 92879bc8a77SHong Zhang } 92979bc8a77SHong Zhang ierr = VecGetSubVector(Y[i],mprk->is_medium,&Ymedium);CHKERRQ(ierr); 93079bc8a77SHong Zhang ierr = VecMAXPY(Ymedium,computedstages,wm,YdotRHS_medium);CHKERRQ(ierr); 93179bc8a77SHong Zhang ierr = VecRestoreSubVector(Y[i],mprk->is_medium,&Ymedium);CHKERRQ(ierr); 93279bc8a77SHong Zhang /* only depends on the medium buffer region */ 93379bc8a77SHong Zhang ierr = TSComputeRHSFunction(mprk->subts_medium,t+h*cmb[i],Y[i],YdotRHS_medium[computedstages]);CHKERRQ(ierr); 9349d6e09e9SHong Zhang /* must be computed after all regions are updated in Y */ 9359d6e09e9SHong Zhang ierr = TSComputeRHSFunction(mprk->subts_slowbuffer,t+h*csb[i],Y[i],YdotRHS_slowbuffer[computedstages]);CHKERRQ(ierr); 93679bc8a77SHong Zhang } 93719c2959aSHong Zhang } 93819c2959aSHong Zhang /* must be computed after fast region and slow region are updated in Y */ 93919c2959aSHong Zhang ierr = TSComputeRHSFunction(mprk->subts_mediumbuffer,t+h*cmb[i],Y[i],YdotRHS_mediumbuffer[i]);CHKERRQ(ierr); 94019c2959aSHong Zhang } 9419d6e09e9SHong Zhang if (!(tab->np == 3 && mprk->is_medium)) { 94219c2959aSHong Zhang ierr = TSComputeRHSFunction(mprk->subts_slowbuffer,t+h*csb[i],Y[i],YdotRHS_slowbuffer[i]);CHKERRQ(ierr); 9439d6e09e9SHong Zhang } 9447c0df07dSHong Zhang ierr = TSComputeRHSFunction(mprk->subts_fast,t+h*cf[i],Y[i],YdotRHS_fast[i]); 9454b84eec9SHong Zhang } 94679bc8a77SHong Zhang 9474b84eec9SHong Zhang ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr); 9484b84eec9SHong Zhang ts->ptime += ts->time_step; 9494b84eec9SHong Zhang ts->time_step = next_time_step; 9504b84eec9SHong Zhang PetscFunctionReturn(0); 9514b84eec9SHong Zhang } 9524b84eec9SHong Zhang 9534b84eec9SHong Zhang static PetscErrorCode TSMPRKTableauReset(TS ts) 9544b84eec9SHong Zhang { 9554b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 9564b84eec9SHong Zhang MPRKTableau tab = mprk->tableau; 9574b84eec9SHong Zhang PetscErrorCode ierr; 9584b84eec9SHong Zhang 9594b84eec9SHong Zhang PetscFunctionBegin; 9604b84eec9SHong Zhang if (!tab) PetscFunctionReturn(0); 9614b84eec9SHong Zhang ierr = PetscFree(mprk->work_fast);CHKERRQ(ierr); 9624b84eec9SHong Zhang ierr = PetscFree(mprk->work_slow);CHKERRQ(ierr); 9637c0df07dSHong Zhang ierr = PetscFree(mprk->work_slowbuffer);CHKERRQ(ierr); 9647c0df07dSHong Zhang ierr = PetscFree(mprk->work_medium);CHKERRQ(ierr); 9657c0df07dSHong Zhang ierr = PetscFree(mprk->work_mediumbuffer);CHKERRQ(ierr); 9664b84eec9SHong Zhang ierr = VecDestroyVecs(tab->s,&mprk->Y);CHKERRQ(ierr); 9670fe4d17eSHong Zhang if (ts->use_splitrhsfunction) { 9680fe4d17eSHong Zhang ierr = VecDestroyVecs(tab->s,&mprk->YdotRHS_fast);CHKERRQ(ierr); 9690fe4d17eSHong Zhang ierr = VecDestroyVecs(tab->s,&mprk->YdotRHS_slow);CHKERRQ(ierr); 9700fe4d17eSHong Zhang ierr = VecDestroyVecs(tab->s,&mprk->YdotRHS_slowbuffer);CHKERRQ(ierr); 9710fe4d17eSHong Zhang ierr = VecDestroyVecs(tab->s,&mprk->YdotRHS_medium);CHKERRQ(ierr); 9720fe4d17eSHong Zhang ierr = VecDestroyVecs(tab->s,&mprk->YdotRHS_mediumbuffer);CHKERRQ(ierr); 9730fe4d17eSHong Zhang } else { 9747c0df07dSHong Zhang ierr = VecDestroyVecs(tab->s,&mprk->YdotRHS);CHKERRQ(ierr); 9757c0df07dSHong Zhang if (mprk->is_slow) { 9767c0df07dSHong Zhang ierr = PetscFree(mprk->YdotRHS_slow);CHKERRQ(ierr); 9774b84eec9SHong Zhang } 9787c0df07dSHong Zhang ierr = PetscFree(mprk->YdotRHS_slowbuffer);CHKERRQ(ierr); 9797c0df07dSHong Zhang if (tab->np == 3) { 9807c0df07dSHong Zhang if (mprk->is_medium) { 9817c0df07dSHong Zhang ierr = PetscFree(mprk->YdotRHS_medium);CHKERRQ(ierr); 9827c0df07dSHong Zhang } 9837c0df07dSHong Zhang ierr = PetscFree(mprk->YdotRHS_mediumbuffer);CHKERRQ(ierr); 9847c0df07dSHong Zhang } 9857c0df07dSHong Zhang ierr = PetscFree(mprk->YdotRHS_fast);CHKERRQ(ierr); 9867c0df07dSHong Zhang } 9874b84eec9SHong Zhang PetscFunctionReturn(0); 9884b84eec9SHong Zhang } 9894b84eec9SHong Zhang 9904b84eec9SHong Zhang static PetscErrorCode TSReset_MPRK(TS ts) 9914b84eec9SHong Zhang { 9924b84eec9SHong Zhang PetscErrorCode ierr; 9934b84eec9SHong Zhang 9944b84eec9SHong Zhang PetscFunctionBegin; 9954b84eec9SHong Zhang ierr = TSMPRKTableauReset(ts);CHKERRQ(ierr); 9964b84eec9SHong Zhang PetscFunctionReturn(0); 9974b84eec9SHong Zhang } 9984b84eec9SHong Zhang 9994b84eec9SHong Zhang static PetscErrorCode DMCoarsenHook_TSMPRK(DM fine,DM coarse,void *ctx) 10004b84eec9SHong Zhang { 10014b84eec9SHong Zhang PetscFunctionBegin; 10024b84eec9SHong Zhang PetscFunctionReturn(0); 10034b84eec9SHong Zhang } 10044b84eec9SHong Zhang 10054b84eec9SHong Zhang static PetscErrorCode DMRestrictHook_TSMPRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 10064b84eec9SHong Zhang { 10074b84eec9SHong Zhang PetscFunctionBegin; 10084b84eec9SHong Zhang PetscFunctionReturn(0); 10094b84eec9SHong Zhang } 10104b84eec9SHong Zhang 10114b84eec9SHong Zhang static PetscErrorCode DMSubDomainHook_TSMPRK(DM dm,DM subdm,void *ctx) 10124b84eec9SHong Zhang { 10134b84eec9SHong Zhang PetscFunctionBegin; 10144b84eec9SHong Zhang PetscFunctionReturn(0); 10154b84eec9SHong Zhang } 10164b84eec9SHong Zhang 10174b84eec9SHong Zhang static PetscErrorCode DMSubDomainRestrictHook_TSMPRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 10184b84eec9SHong Zhang { 10194b84eec9SHong Zhang PetscFunctionBegin; 10204b84eec9SHong Zhang PetscFunctionReturn(0); 10214b84eec9SHong Zhang } 10224b84eec9SHong Zhang 10234b84eec9SHong Zhang static PetscErrorCode TSMPRKTableauSetUp(TS ts) 10244b84eec9SHong Zhang { 10254b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 10264b84eec9SHong Zhang MPRKTableau tab = mprk->tableau; 102719c2959aSHong Zhang Vec YdotRHS_slow,YdotRHS_slowbuffer,YdotRHS_medium,YdotRHS_mediumbuffer,YdotRHS_fast; 10284b84eec9SHong Zhang PetscErrorCode ierr; 10294b84eec9SHong Zhang 10304b84eec9SHong Zhang PetscFunctionBegin; 10314b84eec9SHong Zhang ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&mprk->Y);CHKERRQ(ierr); 10327c0df07dSHong Zhang if (mprk->is_slow) { 103319c2959aSHong Zhang ierr = PetscMalloc1(tab->s,&mprk->work_slow);CHKERRQ(ierr); 10344b84eec9SHong Zhang } 10357c0df07dSHong Zhang ierr = PetscMalloc1(tab->s,&mprk->work_slowbuffer);CHKERRQ(ierr); 10367c0df07dSHong Zhang if (tab->np == 3) { 10377c0df07dSHong Zhang if (mprk->is_medium) { 10387c0df07dSHong Zhang ierr = PetscMalloc1(tab->s,&mprk->work_medium);CHKERRQ(ierr); 10397c0df07dSHong Zhang } 10407c0df07dSHong Zhang ierr = PetscMalloc1(tab->s,&mprk->work_mediumbuffer);CHKERRQ(ierr); 10417c0df07dSHong Zhang } 10427c0df07dSHong Zhang ierr = PetscMalloc1(tab->s,&mprk->work_fast);CHKERRQ(ierr); 10437c0df07dSHong Zhang 10440fe4d17eSHong Zhang if (ts->use_splitrhsfunction) { 104519c2959aSHong Zhang if (mprk->is_slow) { 10464b84eec9SHong Zhang ierr = VecGetSubVector(ts->vec_sol,mprk->is_slow,&YdotRHS_slow);CHKERRQ(ierr); 10474b84eec9SHong Zhang ierr = VecDuplicateVecs(YdotRHS_slow,tab->s,&mprk->YdotRHS_slow);CHKERRQ(ierr); 10484b84eec9SHong Zhang ierr = VecRestoreSubVector(ts->vec_sol,mprk->is_slow,&YdotRHS_slow);CHKERRQ(ierr); 104919c2959aSHong Zhang } 105019c2959aSHong Zhang ierr = VecGetSubVector(ts->vec_sol,mprk->is_slowbuffer,&YdotRHS_slowbuffer);CHKERRQ(ierr); 105119c2959aSHong Zhang ierr = VecDuplicateVecs(YdotRHS_slowbuffer,tab->s,&mprk->YdotRHS_slowbuffer);CHKERRQ(ierr); 105219c2959aSHong Zhang ierr = VecRestoreSubVector(ts->vec_sol,mprk->is_slowbuffer,&YdotRHS_slowbuffer);CHKERRQ(ierr); 105319c2959aSHong Zhang if (tab->np == 3) { 105419c2959aSHong Zhang if (mprk->is_medium) { 105519c2959aSHong Zhang ierr = VecGetSubVector(ts->vec_sol,mprk->is_medium,&YdotRHS_medium);CHKERRQ(ierr); 105619c2959aSHong Zhang ierr = VecDuplicateVecs(YdotRHS_medium,tab->s,&mprk->YdotRHS_medium);CHKERRQ(ierr); 105719c2959aSHong Zhang ierr = VecRestoreSubVector(ts->vec_sol,mprk->is_medium,&YdotRHS_medium);CHKERRQ(ierr); 105819c2959aSHong Zhang } 105919c2959aSHong Zhang ierr = VecGetSubVector(ts->vec_sol,mprk->is_mediumbuffer,&YdotRHS_mediumbuffer);CHKERRQ(ierr); 106019c2959aSHong Zhang ierr = VecDuplicateVecs(YdotRHS_mediumbuffer,tab->s,&mprk->YdotRHS_mediumbuffer);CHKERRQ(ierr); 106119c2959aSHong Zhang ierr = VecRestoreSubVector(ts->vec_sol,mprk->is_mediumbuffer,&YdotRHS_mediumbuffer);CHKERRQ(ierr); 106219c2959aSHong Zhang } 106319c2959aSHong Zhang ierr = VecGetSubVector(ts->vec_sol,mprk->is_fast,&YdotRHS_fast);CHKERRQ(ierr); 106419c2959aSHong Zhang ierr = VecDuplicateVecs(YdotRHS_fast,tab->s,&mprk->YdotRHS_fast);CHKERRQ(ierr); 10654b84eec9SHong Zhang ierr = VecRestoreSubVector(ts->vec_sol,mprk->is_fast,&YdotRHS_fast);CHKERRQ(ierr); 10660fe4d17eSHong Zhang } else { 10670fe4d17eSHong Zhang ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&mprk->YdotRHS);CHKERRQ(ierr); 10680fe4d17eSHong Zhang if (mprk->is_slow) { 10690fe4d17eSHong Zhang ierr = PetscMalloc1(tab->s,&mprk->YdotRHS_slow);CHKERRQ(ierr); 10700fe4d17eSHong Zhang } 10710fe4d17eSHong Zhang ierr = PetscMalloc1(tab->s,&mprk->YdotRHS_slowbuffer);CHKERRQ(ierr); 10720fe4d17eSHong Zhang if (tab->np == 3) { 10730fe4d17eSHong Zhang if (mprk->is_medium) { 10740fe4d17eSHong Zhang ierr = PetscMalloc1(tab->s,&mprk->YdotRHS_medium);CHKERRQ(ierr); 10750fe4d17eSHong Zhang } 10760fe4d17eSHong Zhang ierr = PetscMalloc1(tab->s,&mprk->YdotRHS_mediumbuffer);CHKERRQ(ierr); 10770fe4d17eSHong Zhang } 10780fe4d17eSHong Zhang ierr = PetscMalloc1(tab->s,&mprk->YdotRHS_fast);CHKERRQ(ierr); 10794b84eec9SHong Zhang } 10804b84eec9SHong Zhang PetscFunctionReturn(0); 10814b84eec9SHong Zhang } 10824b84eec9SHong Zhang 10834b84eec9SHong Zhang static PetscErrorCode TSSetUp_MPRK(TS ts) 10844b84eec9SHong Zhang { 10854b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 108619c2959aSHong Zhang MPRKTableau tab = mprk->tableau; 10874b84eec9SHong Zhang DM dm; 10884b84eec9SHong Zhang PetscErrorCode ierr; 10894b84eec9SHong Zhang 10904b84eec9SHong Zhang PetscFunctionBegin; 10914b84eec9SHong Zhang ierr = TSRHSSplitGetIS(ts,"slow",&mprk->is_slow);CHKERRQ(ierr); 10924b84eec9SHong Zhang ierr = TSRHSSplitGetIS(ts,"fast",&mprk->is_fast);CHKERRQ(ierr); 109319c2959aSHong Zhang if (!mprk->is_slow || !mprk->is_fast) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must set up RHSSplits with TSRHSSplitSetIS() using split names 'slow' and 'fast' respectively in order to use the method '%s'",tab->name); 109419c2959aSHong Zhang 109519c2959aSHong Zhang if (tab->np == 3) { 109619c2959aSHong Zhang ierr = TSRHSSplitGetIS(ts,"medium",&mprk->is_medium);CHKERRQ(ierr); 109719c2959aSHong Zhang if (!mprk->is_medium) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must set up RHSSplits with TSRHSSplitSetIS() using split names 'slow' and 'medium' and 'fast' respectively in order to use the method '%s'",tab->name); 109819c2959aSHong Zhang ierr = TSRHSSplitGetIS(ts,"mediumbuffer",&mprk->is_mediumbuffer);CHKERRQ(ierr); 109919c2959aSHong Zhang if (!mprk->is_mediumbuffer) { /* let medium buffer cover whole medium region */ 110019c2959aSHong Zhang mprk->is_mediumbuffer = mprk->is_medium; 110119c2959aSHong Zhang mprk->is_medium = NULL; 110219c2959aSHong Zhang } 110319c2959aSHong Zhang } 110419c2959aSHong Zhang 110519c2959aSHong Zhang /* If users do not provide buffer region settings, the solver will do them automatically, but with a performance penalty */ 110619c2959aSHong Zhang ierr = TSRHSSplitGetIS(ts,"slowbuffer",&mprk->is_slowbuffer);CHKERRQ(ierr); 110719c2959aSHong Zhang if (!mprk->is_slowbuffer) { /* let slow buffer cover whole slow region */ 110819c2959aSHong Zhang mprk->is_slowbuffer = mprk->is_slow; 110919c2959aSHong Zhang mprk->is_slow = NULL; 111019c2959aSHong Zhang } 11114b84eec9SHong Zhang ierr = TSCheckImplicitTerm(ts);CHKERRQ(ierr); 11124b84eec9SHong Zhang ierr = TSMPRKTableauSetUp(ts);CHKERRQ(ierr); 11134b84eec9SHong Zhang ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 11144b84eec9SHong Zhang ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSMPRK,DMRestrictHook_TSMPRK,ts);CHKERRQ(ierr); 11154b84eec9SHong Zhang ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSMPRK,DMSubDomainRestrictHook_TSMPRK,ts);CHKERRQ(ierr); 11160fe4d17eSHong Zhang if (ts->use_splitrhsfunction) { 11170fe4d17eSHong Zhang ts->ops->step = TSStep_MPRKSPLIT; 11180fe4d17eSHong Zhang ts->ops->evaluatestep = TSEvaluateStep_MPRKSPLIT; 11190fe4d17eSHong Zhang ierr = TSMPRKSetSplits(ts);CHKERRQ(ierr); 11200fe4d17eSHong Zhang } else { 11210fe4d17eSHong Zhang ts->ops->step = TSStep_MPRK; 11220fe4d17eSHong Zhang ts->ops->evaluatestep = TSEvaluateStep_MPRK; 11230fe4d17eSHong Zhang } 11244b84eec9SHong Zhang PetscFunctionReturn(0); 11254b84eec9SHong Zhang } 11264b84eec9SHong Zhang 11274b84eec9SHong Zhang static PetscErrorCode TSSetFromOptions_MPRK(PetscOptionItems *PetscOptionsObject,TS ts) 11284b84eec9SHong Zhang { 11294b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 11304b84eec9SHong Zhang PetscErrorCode ierr; 11314b84eec9SHong Zhang 11324b84eec9SHong Zhang PetscFunctionBegin; 11334b84eec9SHong Zhang ierr = PetscOptionsHead(PetscOptionsObject,"PRK ODE solver options");CHKERRQ(ierr); 11344b84eec9SHong Zhang { 11354b84eec9SHong Zhang MPRKTableauLink link; 11364b84eec9SHong Zhang PetscInt count,choice; 11374b84eec9SHong Zhang PetscBool flg; 11384b84eec9SHong Zhang const char **namelist; 11394b84eec9SHong Zhang for (link=MPRKTableauList,count=0; link; link=link->next,count++) ; 11404b84eec9SHong Zhang ierr = PetscMalloc1(count,(char***)&namelist);CHKERRQ(ierr); 11414b84eec9SHong Zhang for (link=MPRKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 11424b84eec9SHong Zhang ierr = PetscOptionsEList("-ts_mprk_type","Family of MPRK method","TSMPRKSetType",(const char*const*)namelist,count,mprk->tableau->name,&choice,&flg);CHKERRQ(ierr); 11434b84eec9SHong Zhang if (flg) {ierr = TSMPRKSetType(ts,namelist[choice]);CHKERRQ(ierr);} 11444b84eec9SHong Zhang ierr = PetscFree(namelist);CHKERRQ(ierr); 11454b84eec9SHong Zhang } 11464b84eec9SHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 11474b84eec9SHong Zhang PetscFunctionReturn(0); 11484b84eec9SHong Zhang } 11494b84eec9SHong Zhang 11504b84eec9SHong Zhang static PetscErrorCode TSView_MPRK(TS ts,PetscViewer viewer) 11514b84eec9SHong Zhang { 11524b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 11534b84eec9SHong Zhang PetscBool iascii; 11544b84eec9SHong Zhang PetscErrorCode ierr; 11554b84eec9SHong Zhang 11564b84eec9SHong Zhang PetscFunctionBegin; 11574b84eec9SHong Zhang ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 11584b84eec9SHong Zhang if (iascii) { 11594b84eec9SHong Zhang MPRKTableau tab = mprk->tableau; 11604b84eec9SHong Zhang TSMPRKType mprktype; 11614b84eec9SHong Zhang char fbuf[512]; 11624b84eec9SHong Zhang char sbuf[512]; 116319c2959aSHong Zhang PetscInt i; 11644b84eec9SHong Zhang ierr = TSMPRKGetType(ts,&mprktype);CHKERRQ(ierr); 11654b84eec9SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," MPRK type %s\n",mprktype);CHKERRQ(ierr); 11664b84eec9SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Order: %D\n",tab->order);CHKERRQ(ierr); 116719c2959aSHong Zhang 11684b84eec9SHong Zhang ierr = PetscFormatRealArray(fbuf,sizeof(fbuf),"% 8.6f",tab->s,tab->cf);CHKERRQ(ierr); 11694b84eec9SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Abscissa cf = %s\n",fbuf);CHKERRQ(ierr); 117019c2959aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Af = \n");CHKERRQ(ierr); 117119c2959aSHong Zhang for (i=0; i<tab->s; i++) { 117219c2959aSHong Zhang ierr = PetscFormatRealArray(fbuf,sizeof(fbuf),"% 8.6f",tab->s,&tab->Af[i*tab->s]);CHKERRQ(ierr); 117319c2959aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," %s\n",fbuf);CHKERRQ(ierr); 117419c2959aSHong Zhang } 117519c2959aSHong Zhang ierr = PetscFormatRealArray(fbuf,sizeof(fbuf),"% 8.6f",tab->s,tab->bf);CHKERRQ(ierr); 117619c2959aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," bf = %s\n",fbuf);CHKERRQ(ierr); 117719c2959aSHong Zhang 117819c2959aSHong Zhang ierr = PetscFormatRealArray(sbuf,sizeof(sbuf),"% 8.6f",tab->s,tab->csb);CHKERRQ(ierr); 117919c2959aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Abscissa csb = %s\n",sbuf);CHKERRQ(ierr); 118019c2959aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Asb = \n");CHKERRQ(ierr); 118119c2959aSHong Zhang for (i=0; i<tab->s; i++) { 118219c2959aSHong Zhang ierr = PetscFormatRealArray(sbuf,sizeof(sbuf),"% 8.6f",tab->s,&tab->Asb[i*tab->s]);CHKERRQ(ierr); 118319c2959aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," %s\n",sbuf);CHKERRQ(ierr); 118419c2959aSHong Zhang } 118519c2959aSHong Zhang ierr = PetscFormatRealArray(sbuf,sizeof(sbuf),"% 8.6f",tab->s,tab->bsb);CHKERRQ(ierr); 118619c2959aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," bsb = %s\n",sbuf);CHKERRQ(ierr); 118719c2959aSHong Zhang 118819c2959aSHong Zhang if (tab->np == 3) { 118919c2959aSHong Zhang char mbuf[512]; 119019c2959aSHong Zhang ierr = PetscFormatRealArray(mbuf,sizeof(mbuf),"% 8.6f",tab->s,tab->cmb);CHKERRQ(ierr); 119119c2959aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Abscissa cmb = %s\n",mbuf);CHKERRQ(ierr); 119219c2959aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Amb = \n");CHKERRQ(ierr); 119319c2959aSHong Zhang for (i=0; i<tab->s; i++) { 119419c2959aSHong Zhang ierr = PetscFormatRealArray(mbuf,sizeof(mbuf),"% 8.6f",tab->s,&tab->Amb[i*tab->s]);CHKERRQ(ierr); 119519c2959aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," %s\n",mbuf);CHKERRQ(ierr); 119619c2959aSHong Zhang } 119719c2959aSHong Zhang ierr = PetscFormatRealArray(mbuf,sizeof(mbuf),"% 8.6f",tab->s,tab->bmb);CHKERRQ(ierr); 119819c2959aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," bmb = %s\n",mbuf);CHKERRQ(ierr); 119919c2959aSHong Zhang } 12004b84eec9SHong Zhang } 12014b84eec9SHong Zhang PetscFunctionReturn(0); 12024b84eec9SHong Zhang } 12034b84eec9SHong Zhang 12044b84eec9SHong Zhang static PetscErrorCode TSLoad_MPRK(TS ts,PetscViewer viewer) 12054b84eec9SHong Zhang { 12064b84eec9SHong Zhang PetscErrorCode ierr; 12074b84eec9SHong Zhang TSAdapt adapt; 12084b84eec9SHong Zhang 12094b84eec9SHong Zhang PetscFunctionBegin; 12104b84eec9SHong Zhang ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 12114b84eec9SHong Zhang ierr = TSAdaptLoad(adapt,viewer);CHKERRQ(ierr); 12124b84eec9SHong Zhang PetscFunctionReturn(0); 12134b84eec9SHong Zhang } 12144b84eec9SHong Zhang 12154b84eec9SHong Zhang /*@C 12164b84eec9SHong Zhang TSMPRKSetType - Set the type of MPRK scheme 12174b84eec9SHong Zhang 12180fe4d17eSHong Zhang Not collective 12194b84eec9SHong Zhang 12204b84eec9SHong Zhang Input Parameter: 12214b84eec9SHong Zhang + ts - timestepping context 12224b84eec9SHong Zhang - mprktype - type of MPRK-scheme 12234b84eec9SHong Zhang 12244b84eec9SHong Zhang Options Database: 12254b84eec9SHong Zhang . -ts_mprk_type - <pm2,p2,p3> 12264b84eec9SHong Zhang 12274b84eec9SHong Zhang Level: intermediate 12284b84eec9SHong Zhang 12294b84eec9SHong Zhang .seealso: TSMPRKGetType(), TSMPRK, TSMPRKType 12304b84eec9SHong Zhang @*/ 12314b84eec9SHong Zhang PetscErrorCode TSMPRKSetType(TS ts,TSMPRKType mprktype) 12324b84eec9SHong Zhang { 12334b84eec9SHong Zhang PetscErrorCode ierr; 12344b84eec9SHong Zhang 12354b84eec9SHong Zhang PetscFunctionBegin; 12364b84eec9SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 12374b84eec9SHong Zhang PetscValidCharPointer(mprktype,2); 12384b84eec9SHong Zhang ierr = PetscTryMethod(ts,"TSMPRKSetType_C",(TS,TSMPRKType),(ts,mprktype));CHKERRQ(ierr); 12394b84eec9SHong Zhang PetscFunctionReturn(0); 12404b84eec9SHong Zhang } 12414b84eec9SHong Zhang 12424b84eec9SHong Zhang /*@C 12434b84eec9SHong Zhang TSMPRKGetType - Get the type of MPRK scheme 12444b84eec9SHong Zhang 12450fe4d17eSHong Zhang Not collective 12464b84eec9SHong Zhang 12474b84eec9SHong Zhang Input Parameter: 12484b84eec9SHong Zhang . ts - timestepping context 12494b84eec9SHong Zhang 12504b84eec9SHong Zhang Output Parameter: 12514b84eec9SHong Zhang . mprktype - type of MPRK-scheme 12524b84eec9SHong Zhang 12534b84eec9SHong Zhang Level: intermediate 12544b84eec9SHong Zhang 12554b84eec9SHong Zhang .seealso: TSMPRKGetType() 12564b84eec9SHong Zhang @*/ 12574b84eec9SHong Zhang PetscErrorCode TSMPRKGetType(TS ts,TSMPRKType *mprktype) 12584b84eec9SHong Zhang { 12594b84eec9SHong Zhang PetscErrorCode ierr; 12604b84eec9SHong Zhang 12614b84eec9SHong Zhang PetscFunctionBegin; 12624b84eec9SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 12634b84eec9SHong Zhang ierr = PetscUseMethod(ts,"TSMPRKGetType_C",(TS,TSMPRKType*),(ts,mprktype));CHKERRQ(ierr); 12644b84eec9SHong Zhang PetscFunctionReturn(0); 12654b84eec9SHong Zhang } 12664b84eec9SHong Zhang 12674b84eec9SHong Zhang static PetscErrorCode TSMPRKGetType_MPRK(TS ts,TSMPRKType *mprktype) 12684b84eec9SHong Zhang { 12694b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 12704b84eec9SHong Zhang 12714b84eec9SHong Zhang PetscFunctionBegin; 12724b84eec9SHong Zhang *mprktype = mprk->tableau->name; 12734b84eec9SHong Zhang PetscFunctionReturn(0); 12744b84eec9SHong Zhang } 12754b84eec9SHong Zhang 12764b84eec9SHong Zhang static PetscErrorCode TSMPRKSetType_MPRK(TS ts,TSMPRKType mprktype) 12774b84eec9SHong Zhang { 12784b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 12794b84eec9SHong Zhang PetscBool match; 12804b84eec9SHong Zhang MPRKTableauLink link; 12814b84eec9SHong Zhang PetscErrorCode ierr; 12824b84eec9SHong Zhang 12834b84eec9SHong Zhang PetscFunctionBegin; 12844b84eec9SHong Zhang if (mprk->tableau) { 12854b84eec9SHong Zhang ierr = PetscStrcmp(mprk->tableau->name,mprktype,&match);CHKERRQ(ierr); 12864b84eec9SHong Zhang if (match) PetscFunctionReturn(0); 12874b84eec9SHong Zhang } 12884b84eec9SHong Zhang for (link = MPRKTableauList; link; link=link->next) { 12894b84eec9SHong Zhang ierr = PetscStrcmp(link->tab.name,mprktype,&match);CHKERRQ(ierr); 12904b84eec9SHong Zhang if (match) { 12914b84eec9SHong Zhang if (ts->setupcalled) {ierr = TSMPRKTableauReset(ts);CHKERRQ(ierr);} 12924b84eec9SHong Zhang mprk->tableau = &link->tab; 12934b84eec9SHong Zhang if (ts->setupcalled) {ierr = TSMPRKTableauSetUp(ts);CHKERRQ(ierr);} 12944b84eec9SHong Zhang PetscFunctionReturn(0); 12954b84eec9SHong Zhang } 12964b84eec9SHong Zhang } 12974b84eec9SHong Zhang SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",mprktype); 12984b84eec9SHong Zhang PetscFunctionReturn(0); 12994b84eec9SHong Zhang } 13004b84eec9SHong Zhang 13014b84eec9SHong Zhang static PetscErrorCode TSGetStages_MPRK(TS ts,PetscInt *ns,Vec **Y) 13024b84eec9SHong Zhang { 13034b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 13044b84eec9SHong Zhang 13054b84eec9SHong Zhang PetscFunctionBegin; 13064b84eec9SHong Zhang *ns = mprk->tableau->s; 13074b84eec9SHong Zhang if (Y) *Y = mprk->Y; 13084b84eec9SHong Zhang PetscFunctionReturn(0); 13094b84eec9SHong Zhang } 13104b84eec9SHong Zhang 13114b84eec9SHong Zhang static PetscErrorCode TSDestroy_MPRK(TS ts) 13124b84eec9SHong Zhang { 13134b84eec9SHong Zhang PetscErrorCode ierr; 13144b84eec9SHong Zhang 13154b84eec9SHong Zhang PetscFunctionBegin; 13164b84eec9SHong Zhang ierr = TSReset_MPRK(ts);CHKERRQ(ierr); 13174b84eec9SHong Zhang if (ts->dm) { 13184b84eec9SHong Zhang ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSMPRK,DMRestrictHook_TSMPRK,ts);CHKERRQ(ierr); 13194b84eec9SHong Zhang ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSMPRK,DMSubDomainRestrictHook_TSMPRK,ts);CHKERRQ(ierr); 13204b84eec9SHong Zhang } 13214b84eec9SHong Zhang ierr = PetscFree(ts->data);CHKERRQ(ierr); 13224b84eec9SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSMPRKGetType_C",NULL);CHKERRQ(ierr); 13234b84eec9SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSMPRKSetType_C",NULL);CHKERRQ(ierr); 13244b84eec9SHong Zhang PetscFunctionReturn(0); 13254b84eec9SHong Zhang } 13264b84eec9SHong Zhang 13274b84eec9SHong Zhang /*MC 1328*f944a40eSHong Zhang TSMPRK - ODE solver using Multirate Partitioned Runge-Kutta schemes 13294b84eec9SHong Zhang 13304b84eec9SHong Zhang The user should provide the right hand side of the equation 13314b84eec9SHong Zhang using TSSetRHSFunction(). 13324b84eec9SHong Zhang 13334b84eec9SHong Zhang Notes: 1334*f944a40eSHong Zhang The default is TSMPRKPM2, it can be changed with TSMPRKSetType() or -ts_mprk_type 13354b84eec9SHong Zhang 13364b84eec9SHong Zhang Level: beginner 13374b84eec9SHong Zhang 13384b84eec9SHong Zhang .seealso: TSCreate(), TS, TSSetType(), TSMPRKSetType(), TSMPRKGetType(), TSMPRKType, TSMPRKRegister(), TSMPRKSetMultirateType() 13394b84eec9SHong Zhang TSMPRKM2, TSMPRKM3, TSMPRKRFSMR3, TSMPRKRFSMR2 13404b84eec9SHong Zhang 13414b84eec9SHong Zhang M*/ 13424b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSCreate_MPRK(TS ts) 13434b84eec9SHong Zhang { 13444b84eec9SHong Zhang TS_MPRK *mprk; 13454b84eec9SHong Zhang PetscErrorCode ierr; 13464b84eec9SHong Zhang 13474b84eec9SHong Zhang PetscFunctionBegin; 13484b84eec9SHong Zhang ierr = TSMPRKInitializePackage();CHKERRQ(ierr); 13494b84eec9SHong Zhang 13504b84eec9SHong Zhang ts->ops->reset = TSReset_MPRK; 13514b84eec9SHong Zhang ts->ops->destroy = TSDestroy_MPRK; 13524b84eec9SHong Zhang ts->ops->view = TSView_MPRK; 13534b84eec9SHong Zhang ts->ops->load = TSLoad_MPRK; 13544b84eec9SHong Zhang ts->ops->setup = TSSetUp_MPRK; 13554b84eec9SHong Zhang ts->ops->step = TSStep_MPRK; 13564b84eec9SHong Zhang ts->ops->evaluatestep = TSEvaluateStep_MPRK; 13574b84eec9SHong Zhang ts->ops->setfromoptions = TSSetFromOptions_MPRK; 13584b84eec9SHong Zhang ts->ops->getstages = TSGetStages_MPRK; 13594b84eec9SHong Zhang 13604b84eec9SHong Zhang ierr = PetscNewLog(ts,&mprk);CHKERRQ(ierr); 13614b84eec9SHong Zhang ts->data = (void*)mprk; 13624b84eec9SHong Zhang 13634b84eec9SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSMPRKGetType_C",TSMPRKGetType_MPRK);CHKERRQ(ierr); 13644b84eec9SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSMPRKSetType_C",TSMPRKSetType_MPRK);CHKERRQ(ierr); 13654b84eec9SHong Zhang 13664b84eec9SHong Zhang ierr = TSMPRKSetType(ts,TSMPRKDefault);CHKERRQ(ierr); 13674b84eec9SHong Zhang PetscFunctionReturn(0); 13684b84eec9SHong Zhang } 1369