14b84eec9SHong Zhang /* 2f944a40eSHong 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 6f944a40eSHong Zhang Udot = F(t,U) 7f944a40eSHong 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) 10f944a40eSHong Zhang Ufdot = Ff(t,Us,Uf) 11f944a40eSHong Zhang for component-wise partitioned system, 12f944a40eSHong 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 131f944a40eSHong 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 143f944a40eSHong 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 155f944a40eSHong 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 167f944a40eSHong 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 179f944a40eSHong 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 191f944a40eSHong 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 203f944a40eSHong 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 .seealso: TSMPRKRegisterDestroy() 2234b84eec9SHong Zhang @*/ 2244b84eec9SHong Zhang PetscErrorCode TSMPRKRegisterAll(void) 2254b84eec9SHong Zhang { 2264b84eec9SHong Zhang PetscErrorCode ierr; 2274b84eec9SHong Zhang 2284b84eec9SHong Zhang PetscFunctionBegin; 2294b84eec9SHong Zhang if (TSMPRKRegisterAllCalled) PetscFunctionReturn(0); 2304b84eec9SHong Zhang TSMPRKRegisterAllCalled = PETSC_TRUE; 2314b84eec9SHong Zhang 2324b84eec9SHong Zhang #define RC PetscRealConstant 2334b84eec9SHong Zhang { 2344b84eec9SHong Zhang const PetscReal 23519c2959aSHong Zhang Abase[2][2] = {{0,0}, 23619c2959aSHong Zhang {RC(1.0),0}}, 23719c2959aSHong Zhang bbase[2] = {RC(0.5),RC(0.5)}; 23819c2959aSHong Zhang PetscReal 23919c2959aSHong Zhang Asb[4][4] = {{0}},Af[4][4] = {{0}},bsb[4] = {0},bf[4] = {0}; 24019c2959aSHong Zhang PetscInt 24119c2959aSHong Zhang rsb[4] = {0,0,1,2}; 24219c2959aSHong Zhang ierr = TSMPRKGenerateTableau2(2,2,&Abase[0][0],bbase,&Asb[0][0],bsb,&Af[0][0],bf);CHKERRQ(ierr); 24379bc8a77SHong 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); 2444b84eec9SHong Zhang } 24519c2959aSHong Zhang { 24619c2959aSHong Zhang const PetscReal 24719c2959aSHong Zhang Abase[2][2] = {{0,0}, 24819c2959aSHong Zhang {RC(1.0),0}}, 24919c2959aSHong Zhang bbase[2] = {RC(0.5),RC(0.5)}; 25019c2959aSHong Zhang PetscReal 25119c2959aSHong Zhang Asb[8][8] = {{0}},Amb[8][8] = {{0}},Af[8][8] = {{0}},bsb[8] ={0},bmb[8] = {0},bf[8] = {0}; 25219c2959aSHong Zhang PetscInt 25319c2959aSHong Zhang rsb[8] = {0,0,1,2,1,2,1,2},rmb[8] = {0,0,1,2,0,0,5,6}; 25419c2959aSHong Zhang ierr = TSMPRKGenerateTableau3(2,2,&Abase[0][0],bbase,&Asb[0][0],bsb,&Amb[0][0],bmb,&Af[0][0],bf);CHKERRQ(ierr); 25579bc8a77SHong 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); 25619c2959aSHong Zhang } 25719c2959aSHong Zhang { 25819c2959aSHong Zhang const PetscReal 25919c2959aSHong Zhang Abase[2][2] = {{0,0}, 26019c2959aSHong Zhang {RC(1.0),0}}, 26119c2959aSHong Zhang bbase[2] = {RC(0.5),RC(0.5)}; 26219c2959aSHong Zhang PetscReal 26319c2959aSHong Zhang Asb[6][6] = {{0}},Af[6][6] = {{0}},bsb[6] = {0},bf[6] = {0}; 26419c2959aSHong Zhang PetscInt 26519c2959aSHong Zhang rsb[6] = {0,0,1,2,1,2}; 26619c2959aSHong Zhang ierr = TSMPRKGenerateTableau2(3,2,&Abase[0][0],bbase,&Asb[0][0],bsb,&Af[0][0],bf);CHKERRQ(ierr); 26779bc8a77SHong 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); 26819c2959aSHong Zhang } 26919c2959aSHong Zhang { 27019c2959aSHong Zhang const PetscReal 27119c2959aSHong Zhang Abase[2][2] = {{0,0}, 27219c2959aSHong Zhang {RC(1.0),0}}, 27319c2959aSHong Zhang bbase[2] = {RC(0.5),RC(0.5)}; 27419c2959aSHong Zhang PetscReal 27519c2959aSHong Zhang Asb[18][18] = {{0}},Amb[18][18] = {{0}},Af[18][18] = {{0}},bsb[18] ={0},bmb[18] = {0},bf[18] = {0}; 27619c2959aSHong Zhang PetscInt 27719c2959aSHong 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}; 27819c2959aSHong Zhang ierr = TSMPRKGenerateTableau3(3,2,&Abase[0][0],bbase,&Asb[0][0],bsb,&Amb[0][0],bmb,&Af[0][0],bf);CHKERRQ(ierr); 27979bc8a77SHong 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); 28019c2959aSHong Zhang } 28119c2959aSHong Zhang /* 28219c2959aSHong Zhang PetscReal 28319c2959aSHong Zhang Asb[8][8] = {{Abase[0][0],Abase[0][1],0,0,0,0,0,0}, 28419c2959aSHong Zhang {Abase[1][0],Abase[1][1],0,0,0,0,0,0}, 28519c2959aSHong Zhang {0,0,Abase[0][0],Abase[0][1],0,0,0,0}, 28619c2959aSHong Zhang {0,0,Abase[1][0],Abase[1][1],0,0,0,0}, 28719c2959aSHong Zhang {0,0,0,0,Abase[0][0],Abase[0][1],0,0}, 28819c2959aSHong Zhang {0,0,0,0,Abase[1][0],Abase[1][1],0,0}, 28919c2959aSHong Zhang {0,0,0,0,0,0,Abase[0][0],Abase[0][1]}, 29019c2959aSHong Zhang {0,0,0,0,0,0,Abase[1][0],Abase[1][1]}}, 29119c2959aSHong Zhang Amb[8][8] = {{Abase[0][0]/m,Abase[0][1]/m,0,0,0,0,0,0}, 29219c2959aSHong Zhang {Abase[1][0]/m,Abase[1][1]/m,0,0,0,0,0,0}, 29319c2959aSHong Zhang {0,0,Abase[0][0]/m,Abase[0][1]/m,0,0,0,0}, 29419c2959aSHong Zhang {0,0,Abase[1][0]/m,Abase[1][1]/m,0,0,0,0}, 29519c2959aSHong Zhang {bbase[0]/m,bbase[1]/m,bbase[0]/m,bbase[1]/m,Abase[0][0]/m,Abase[0][1]/m,0,0}, 29619c2959aSHong Zhang {bbase[0]/m,bbase[1]/m,bbase[0]/m,bbase[1]/m,Abase[1][0]/m,Abase[1][1]/m,0,0}, 29719c2959aSHong Zhang {bbase[0]/m,bbase[1]/m,bbase[0]/m,bbase[1]/m,0,0,Abase[0][0]/m,Abase[0][1]/m}, 29819c2959aSHong Zhang {bbase[0]/m,bbase[1]/m,bbase[0]/m,bbase[1]/m,0,0,Abase[1][0]/m,Abase[1][1]/m}}, 29919c2959aSHong Zhang Af[8][8] = {{Abase[0][0]/m/m,Abase[0][1]/m/m,0,0,0,0,0,0}, 30019c2959aSHong Zhang {Abase[1][0]/m/m,Abase[1][1]/m/m,0,0,0,0,0,0}, 30119c2959aSHong Zhang {0,0,Abase[0][0]/m/m,Abase[0][1]/m/m,0,0,0,0}, 30219c2959aSHong Zhang {0,0,Abase[1][0]/m/m,Abase[1][1]/m/m,0,0,0,0}, 30319c2959aSHong 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}, 30419c2959aSHong 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}, 30519c2959aSHong 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}, 30619c2959aSHong 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}}, 30719c2959aSHong 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}, 30819c2959aSHong 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}, 30919c2959aSHong 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}, 31019c2959aSHong Zhang */ 3114b84eec9SHong Zhang /*{ 3124b84eec9SHong Zhang const PetscReal 3134b84eec9SHong Zhang As[8][8] = {{0,0,0,0,0,0,0,0}, 3144b84eec9SHong Zhang {RC(1.0)/RC(2.0),0,0,0,0,0,0,0}, 3154b84eec9SHong Zhang {RC(-1.0)/RC(6.0),RC(2.0)/RC(3.0),0,0,0,0,0,0}, 3164b84eec9SHong Zhang {RC(1.0)/RC(3.0),RC(-1.0)/RC(3.0),RC(1.0),0,0,0,0,0}, 3174b84eec9SHong Zhang {0,0,0,0,0,0,0,0}, 3184b84eec9SHong Zhang {0,0,0,0,RC(1.0)/RC(2.0),0,0,0}, 3194b84eec9SHong Zhang {0,0,0,0,RC(-1.0)/RC(6.0),RC(2.0)/RC(3.0),0,0}, 3204b84eec9SHong Zhang {0,0,0,0,RC(1.0)/RC(3.0),RC(-1.0)/RC(3.0),RC(1.0),0}}, 3214b84eec9SHong Zhang A[8][8] = {{0,0,0,0,0,0,0,0}, 3224b84eec9SHong Zhang {RC(1.0)/RC(4.0),0,0,0,0,0,0,0}, 3234b84eec9SHong Zhang {RC(-1.0)/RC(12.0),RC(1.0)/RC(3.0),0,0,0,0,0,0}, 3244b84eec9SHong Zhang {RC(1.0)/RC(6.0),RC(-1.0)/RC(6.0),RC(1.0)/RC(2.0),0,0,0,0,0}, 3254b84eec9SHong 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}, 3264b84eec9SHong 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}, 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),RC(-1.0)/RC(12.0),RC(1.0)/RC(3.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(6.0),RC(-1.0)/RC(6.0),RC(1.0)/RC(2.0),0}}, 3294b84eec9SHong 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)}, 3304b84eec9SHong 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)}; 3314b84eec9SHong Zhang ierr = TSMPRKRegister(TSMPRKPM3,3,8,&As[0][0],bs,NULL,&A[0][0],b,NULL);CHKERRQ(ierr); 3324b84eec9SHong Zhang }*/ 3334b84eec9SHong Zhang 3344b84eec9SHong Zhang { 3354b84eec9SHong Zhang const PetscReal 33619c2959aSHong Zhang Asb[5][5] = {{0,0,0,0,0}, 3374b84eec9SHong Zhang {RC(1.0)/RC(2.0),0,0,0,0}, 3384b84eec9SHong Zhang {RC(1.0)/RC(2.0),0,0,0,0}, 3394b84eec9SHong Zhang {RC(1.0),0,0,0,0}, 3404b84eec9SHong Zhang {RC(1.0),0,0,0,0}}, 34119c2959aSHong Zhang Af[5][5] = {{0,0,0,0,0}, 3424b84eec9SHong Zhang {RC(1.0)/RC(2.0),0,0,0,0}, 3434b84eec9SHong Zhang {RC(1.0)/RC(4.0),RC(1.0)/RC(4.0),0,0,0}, 3444b84eec9SHong Zhang {RC(1.0)/RC(4.0),RC(1.0)/RC(4.0),RC(1.0)/RC(2.0),0,0}, 3454b84eec9SHong 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}}, 34619c2959aSHong Zhang bsb[5] = {RC(1.0)/RC(2.0),0,0,0,RC(1.0)/RC(2.0)}, 34719c2959aSHong 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}; 34819c2959aSHong Zhang const PetscInt 34919c2959aSHong Zhang rsb[5] = {0,0,2,0,4}; 35079bc8a77SHong 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); 3514b84eec9SHong Zhang } 3524b84eec9SHong Zhang 3534b84eec9SHong Zhang { 3544b84eec9SHong Zhang const PetscReal 35519c2959aSHong Zhang Asb[10][10] = {{0,0,0,0,0,0,0,0,0,0}, 3564b84eec9SHong Zhang {RC(1.0)/RC(4.0),0,0,0,0,0,0,0,0,0}, 3574b84eec9SHong Zhang {RC(1.0)/RC(4.0),0,0,0,0,0,0,0,0,0}, 3584b84eec9SHong Zhang {RC(1.0)/RC(2.0),0,0,0,0,0,0,0,0,0}, 3594b84eec9SHong Zhang {RC(1.0)/RC(2.0),0,0,0,0,0,0,0,0,0}, 3604b84eec9SHong Zhang {RC(-1.0)/RC(6.0),0,0,0,RC(2.0)/RC(3.0),0,0,0,0,0}, 3614b84eec9SHong 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}, 3624b84eec9SHong 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}, 3634b84eec9SHong Zhang {RC(1.0)/RC(3.0),0,0,0,RC(-1.0)/RC(3.0),RC(1.0),0,0,0,0}, 3644b84eec9SHong Zhang {RC(1.0)/RC(3.0),0,0,0,RC(-1.0)/RC(3.0),RC(1.0),0,0,0,0}}, 36519c2959aSHong Zhang Af[10][10] = {{0,0,0,0,0,0,0,0,0,0}, 3664b84eec9SHong Zhang {RC(1.0)/RC(4.0),0,0,0,0,0,0,0,0,0}, 3674b84eec9SHong Zhang {RC(-1.0)/RC(12.0),RC(1.0)/RC(3.0),0,0,0,0,0,0,0,0}, 3684b84eec9SHong 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}, 3694b84eec9SHong 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}, 3704b84eec9SHong 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}, 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,RC(1.0)/RC(4.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,RC(-1.0)/RC(12.0),RC(1.0)/RC(3.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(6.0),RC(-1.0)/RC(6.0),RC(1.0)/RC(2.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(6.0),RC(1.0)/RC(6.0),RC(1.0)/RC(12.0),0}}, 37519c2959aSHong 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)}, 37619c2959aSHong 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}; 37719c2959aSHong Zhang const PetscInt 37819c2959aSHong Zhang rsb[10] = {0,0,2,0,4,0,0,7,0,9}; 37979bc8a77SHong 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); 3804b84eec9SHong Zhang } 3814b84eec9SHong Zhang #undef RC 3824b84eec9SHong Zhang PetscFunctionReturn(0); 3834b84eec9SHong Zhang } 3844b84eec9SHong Zhang 3854b84eec9SHong Zhang /*@C 3864b84eec9SHong Zhang TSMPRKRegisterDestroy - Frees the list of schemes that were registered by TSMPRKRegister(). 3874b84eec9SHong Zhang 3884b84eec9SHong Zhang Not Collective 3894b84eec9SHong Zhang 3904b84eec9SHong Zhang Level: advanced 3914b84eec9SHong Zhang 3924b84eec9SHong Zhang .seealso: TSMPRKRegister(), TSMPRKRegisterAll() 3934b84eec9SHong Zhang @*/ 3944b84eec9SHong Zhang PetscErrorCode TSMPRKRegisterDestroy(void) 3954b84eec9SHong Zhang { 3964b84eec9SHong Zhang PetscErrorCode ierr; 3974b84eec9SHong Zhang MPRKTableauLink link; 3984b84eec9SHong Zhang 3994b84eec9SHong Zhang PetscFunctionBegin; 4004b84eec9SHong Zhang while ((link = MPRKTableauList)) { 4014b84eec9SHong Zhang MPRKTableau t = &link->tab; 4024b84eec9SHong Zhang MPRKTableauList = link->next; 40319c2959aSHong Zhang ierr = PetscFree3(t->Asb,t->bsb,t->csb);CHKERRQ(ierr); 40419c2959aSHong Zhang ierr = PetscFree3(t->Amb,t->bmb,t->cmb);CHKERRQ(ierr); 4054b84eec9SHong Zhang ierr = PetscFree3(t->Af,t->bf,t->cf);CHKERRQ(ierr); 40641eb3033SHong Zhang ierr = PetscFree(t->rsb);CHKERRQ(ierr); 40741eb3033SHong Zhang ierr = PetscFree(t->rmb);CHKERRQ(ierr); 4084b84eec9SHong Zhang ierr = PetscFree(t->name);CHKERRQ(ierr); 4094b84eec9SHong Zhang ierr = PetscFree(link);CHKERRQ(ierr); 4104b84eec9SHong Zhang } 4114b84eec9SHong Zhang TSMPRKRegisterAllCalled = PETSC_FALSE; 4124b84eec9SHong Zhang PetscFunctionReturn(0); 4134b84eec9SHong Zhang } 4144b84eec9SHong Zhang 4154b84eec9SHong Zhang /*@C 4164b84eec9SHong Zhang TSMPRKInitializePackage - This function initializes everything in the TSMPRK package. It is called 4174b84eec9SHong Zhang from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_MPRK() 4184b84eec9SHong Zhang when using static libraries. 4194b84eec9SHong Zhang 4204b84eec9SHong Zhang Level: developer 4214b84eec9SHong Zhang 4224b84eec9SHong Zhang .seealso: PetscInitialize() 4234b84eec9SHong Zhang @*/ 4244b84eec9SHong Zhang PetscErrorCode TSMPRKInitializePackage(void) 4254b84eec9SHong Zhang { 4264b84eec9SHong Zhang PetscErrorCode ierr; 4274b84eec9SHong Zhang 4284b84eec9SHong Zhang PetscFunctionBegin; 4294b84eec9SHong Zhang if (TSMPRKPackageInitialized) PetscFunctionReturn(0); 4304b84eec9SHong Zhang TSMPRKPackageInitialized = PETSC_TRUE; 4314b84eec9SHong Zhang ierr = TSMPRKRegisterAll();CHKERRQ(ierr); 4324b84eec9SHong Zhang ierr = PetscRegisterFinalize(TSMPRKFinalizePackage);CHKERRQ(ierr); 4334b84eec9SHong Zhang PetscFunctionReturn(0); 4344b84eec9SHong Zhang } 4354b84eec9SHong Zhang 4364b84eec9SHong Zhang /*@C 437f944a40eSHong Zhang TSMPRKFinalizePackage - This function destroys everything in the TSMPRK package. It is 4384b84eec9SHong Zhang called from PetscFinalize(). 4394b84eec9SHong Zhang 4404b84eec9SHong Zhang Level: developer 4414b84eec9SHong Zhang 4424b84eec9SHong Zhang .seealso: PetscFinalize() 4434b84eec9SHong Zhang @*/ 4444b84eec9SHong Zhang PetscErrorCode TSMPRKFinalizePackage(void) 4454b84eec9SHong Zhang { 4464b84eec9SHong Zhang PetscErrorCode ierr; 4474b84eec9SHong Zhang 4484b84eec9SHong Zhang PetscFunctionBegin; 4494b84eec9SHong Zhang TSMPRKPackageInitialized = PETSC_FALSE; 4504b84eec9SHong Zhang ierr = TSMPRKRegisterDestroy();CHKERRQ(ierr); 4514b84eec9SHong Zhang PetscFunctionReturn(0); 4524b84eec9SHong Zhang } 4534b84eec9SHong Zhang 4544b84eec9SHong Zhang /*@C 4554b84eec9SHong Zhang TSMPRKRegister - register a MPRK scheme by providing the entries in the Butcher tableau 4564b84eec9SHong Zhang 4574b84eec9SHong Zhang Not Collective, but the same schemes should be registered on all processes on which they will be used 4584b84eec9SHong Zhang 4594b84eec9SHong Zhang Input Parameters: 4604b84eec9SHong Zhang + name - identifier for method 4614b84eec9SHong Zhang . order - approximation order of method 46279bc8a77SHong Zhang . s - number of stages in the base methods 46379bc8a77SHong Zhang . ratio1 - stepsize ratio at 1st level (e.g. slow/medium) 46479bc8a77SHong Zhang . ratio2 - stepsize ratio at 2nd level (e.g. medium/fast) 4654b84eec9SHong Zhang . Af - stage coefficients for fast components(dimension s*s, row-major) 4664b84eec9SHong Zhang . bf - step completion table for fast components(dimension s) 4674b84eec9SHong Zhang . cf - abscissa for fast components(dimension s) 4684b84eec9SHong Zhang . As - stage coefficients for slow components(dimension s*s, row-major) 4694b84eec9SHong Zhang . bs - step completion table for slow components(dimension s) 4704b84eec9SHong Zhang - cs - abscissa for slow components(dimension s) 4714b84eec9SHong Zhang 4724b84eec9SHong Zhang Notes: 4734b84eec9SHong Zhang Several MPRK methods are provided, this function is only needed to create new methods. 4744b84eec9SHong Zhang 4754b84eec9SHong Zhang Level: advanced 4764b84eec9SHong Zhang 4774b84eec9SHong Zhang .seealso: TSMPRK 4784b84eec9SHong Zhang @*/ 47979bc8a77SHong Zhang PetscErrorCode TSMPRKRegister(TSMPRKType name,PetscInt order, 48079bc8a77SHong Zhang PetscInt sbase,PetscInt ratio1,PetscInt ratio2, 48119c2959aSHong Zhang const PetscReal Asb[],const PetscReal bsb[],const PetscReal csb[],const PetscInt rsb[], 48219c2959aSHong Zhang const PetscReal Amb[],const PetscReal bmb[],const PetscReal cmb[],const PetscInt rmb[], 4834b84eec9SHong Zhang const PetscReal Af[],const PetscReal bf[],const PetscReal cf[]) 4844b84eec9SHong Zhang { 4854b84eec9SHong Zhang MPRKTableauLink link; 4864b84eec9SHong Zhang MPRKTableau t; 48779bc8a77SHong Zhang PetscInt s,i,j; 4884b84eec9SHong Zhang PetscErrorCode ierr; 4894b84eec9SHong Zhang 4904b84eec9SHong Zhang PetscFunctionBegin; 4914b84eec9SHong Zhang PetscValidCharPointer(name,1); 492064a246eSJacob Faibussowitsch PetscValidRealPointer(Asb,6); 49379bc8a77SHong Zhang if (bsb) PetscValidRealPointer(bsb,7); 49479bc8a77SHong Zhang if (csb) PetscValidRealPointer(csb,8); 495064a246eSJacob Faibussowitsch if (rsb) PetscValidIntPointer(rsb,9); 49679bc8a77SHong Zhang if (Amb) PetscValidRealPointer(Amb,10); 49779bc8a77SHong Zhang if (bmb) PetscValidRealPointer(bmb,11); 49879bc8a77SHong Zhang if (cmb) PetscValidRealPointer(cmb,12); 499064a246eSJacob Faibussowitsch if (rmb) PetscValidIntPointer(rmb,13); 50079bc8a77SHong Zhang PetscValidRealPointer(Af,14); 50179bc8a77SHong Zhang if (bf) PetscValidRealPointer(bf,15); 50279bc8a77SHong Zhang if (cf) PetscValidRealPointer(cf,16); 5034b84eec9SHong Zhang 5044b84eec9SHong Zhang ierr = PetscNew(&link);CHKERRQ(ierr); 5054b84eec9SHong Zhang t = &link->tab; 5064b84eec9SHong Zhang 5074b84eec9SHong Zhang ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr); 50879bc8a77SHong Zhang s = sbase*ratio1*ratio2; /* this is the dimension of the matrices below */ 5094b84eec9SHong Zhang t->order = order; 51079bc8a77SHong Zhang t->sbase = sbase; 5114b84eec9SHong Zhang t->s = s; 51219c2959aSHong Zhang t->np = 2; 51379bc8a77SHong Zhang 5144b84eec9SHong Zhang ierr = PetscMalloc3(s*s,&t->Af,s,&t->bf,s,&t->cf);CHKERRQ(ierr); 515580bdb30SBarry Smith ierr = PetscArraycpy(t->Af,Af,s*s);CHKERRQ(ierr); 5164b84eec9SHong Zhang if (bf) { 517580bdb30SBarry Smith ierr = PetscArraycpy(t->bf,bf,s);CHKERRQ(ierr); 51819c2959aSHong Zhang } else 5194b84eec9SHong Zhang for (i=0; i<s; i++) t->bf[i] = Af[(s-1)*s+i]; 5204b84eec9SHong Zhang if (cf) { 521580bdb30SBarry Smith ierr = PetscArraycpy(t->cf,cf,s);CHKERRQ(ierr); 52219c2959aSHong Zhang } else { 5234b84eec9SHong Zhang for (i=0; i<s; i++) 5244b84eec9SHong Zhang for (j=0,t->cf[i]=0; j<s; j++) 5254b84eec9SHong Zhang t->cf[i] += Af[i*s+j]; 5264b84eec9SHong Zhang } 52719c2959aSHong Zhang 52819c2959aSHong Zhang if (Amb) { 52919c2959aSHong Zhang t->np = 3; 53019c2959aSHong Zhang ierr = PetscMalloc3(s*s,&t->Amb,s,&t->bmb,s,&t->cmb);CHKERRQ(ierr); 531580bdb30SBarry Smith ierr = PetscArraycpy(t->Amb,Amb,s*s);CHKERRQ(ierr); 53219c2959aSHong Zhang if (bmb) { 533580bdb30SBarry Smith ierr = PetscArraycpy(t->bmb,bmb,s);CHKERRQ(ierr); 53419c2959aSHong Zhang } else { 53519c2959aSHong Zhang for (i=0; i<s; i++) t->bmb[i] = Amb[(s-1)*s+i]; 5364b84eec9SHong Zhang } 53719c2959aSHong Zhang if (cmb) { 538580bdb30SBarry Smith ierr = PetscArraycpy(t->cmb,cmb,s);CHKERRQ(ierr); 53919c2959aSHong Zhang } else { 5404b84eec9SHong Zhang for (i=0; i<s; i++) 54119c2959aSHong Zhang for (j=0,t->cmb[i]=0; j<s; j++) 54219c2959aSHong Zhang t->cmb[i] += Amb[i*s+j]; 54319c2959aSHong Zhang } 54419c2959aSHong Zhang if (rmb) { 545071fcb05SBarry Smith ierr = PetscMalloc1(s,&t->rmb);CHKERRQ(ierr); 546580bdb30SBarry Smith ierr = PetscArraycpy(t->rmb,rmb,s);CHKERRQ(ierr); 547071fcb05SBarry Smith } else { 548071fcb05SBarry Smith ierr = PetscCalloc1(s,&t->rmb);CHKERRQ(ierr); 54919c2959aSHong Zhang } 55019c2959aSHong Zhang } 55119c2959aSHong Zhang 55219c2959aSHong Zhang ierr = PetscMalloc3(s*s,&t->Asb,s,&t->bsb,s,&t->csb);CHKERRQ(ierr); 553580bdb30SBarry Smith ierr = PetscArraycpy(t->Asb,Asb,s*s);CHKERRQ(ierr); 55419c2959aSHong Zhang if (bsb) { 555580bdb30SBarry Smith ierr = PetscArraycpy(t->bsb,bsb,s);CHKERRQ(ierr); 55619c2959aSHong Zhang } else 55719c2959aSHong Zhang for (i=0; i<s; i++) t->bsb[i] = Asb[(s-1)*s+i]; 55819c2959aSHong Zhang if (csb) { 559580bdb30SBarry Smith ierr = PetscArraycpy(t->csb,csb,s);CHKERRQ(ierr); 56019c2959aSHong Zhang } else { 56119c2959aSHong Zhang for (i=0; i<s; i++) 56219c2959aSHong Zhang for (j=0,t->csb[i]=0; j<s; j++) 56319c2959aSHong Zhang t->csb[i] += Asb[i*s+j]; 56419c2959aSHong Zhang } 56519c2959aSHong Zhang if (rsb) { 566071fcb05SBarry Smith ierr = PetscMalloc1(s,&t->rsb);CHKERRQ(ierr); 567580bdb30SBarry Smith ierr = PetscArraycpy(t->rsb,rsb,s);CHKERRQ(ierr); 568071fcb05SBarry Smith } else { 569071fcb05SBarry Smith ierr = PetscCalloc1(s,&t->rsb);CHKERRQ(ierr); 5704b84eec9SHong Zhang } 5714b84eec9SHong Zhang link->next = MPRKTableauList; 5724b84eec9SHong Zhang MPRKTableauList = link; 5734b84eec9SHong Zhang PetscFunctionReturn(0); 5744b84eec9SHong Zhang } 5754b84eec9SHong Zhang 5764b84eec9SHong Zhang static PetscErrorCode TSMPRKSetSplits(TS ts) 5774b84eec9SHong Zhang { 5784b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 57919c2959aSHong Zhang MPRKTableau tab = mprk->tableau; 5804b84eec9SHong Zhang DM dm,subdm,newdm; 5814b84eec9SHong Zhang PetscErrorCode ierr; 5824b84eec9SHong Zhang 5834b84eec9SHong Zhang PetscFunctionBegin; 5844b84eec9SHong Zhang ierr = TSRHSSplitGetSubTS(ts,"slow",&mprk->subts_slow);CHKERRQ(ierr); 5854b84eec9SHong Zhang ierr = TSRHSSplitGetSubTS(ts,"fast",&mprk->subts_fast);CHKERRQ(ierr); 5864b84eec9SHong 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"); 5874b84eec9SHong Zhang 58819c2959aSHong Zhang /* Only copy the DM */ 5894b84eec9SHong Zhang ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 59019c2959aSHong Zhang 59119c2959aSHong Zhang ierr = TSRHSSplitGetSubTS(ts,"slowbuffer",&mprk->subts_slowbuffer);CHKERRQ(ierr); 59219c2959aSHong Zhang if (!mprk->subts_slowbuffer) { 59319c2959aSHong Zhang mprk->subts_slowbuffer = mprk->subts_slow; 59419c2959aSHong Zhang mprk->subts_slow = NULL; 59519c2959aSHong Zhang } 59619c2959aSHong Zhang if (mprk->subts_slow) { 5974b84eec9SHong Zhang ierr = DMClone(dm,&newdm);CHKERRQ(ierr); 5984b84eec9SHong Zhang ierr = TSGetDM(mprk->subts_slow,&subdm);CHKERRQ(ierr); 5994b84eec9SHong Zhang ierr = DMCopyDMTS(subdm,newdm);CHKERRQ(ierr); 6004b84eec9SHong Zhang ierr = DMCopyDMSNES(subdm,newdm);CHKERRQ(ierr); 6014b84eec9SHong Zhang ierr = TSSetDM(mprk->subts_slow,newdm);CHKERRQ(ierr); 6024b84eec9SHong Zhang ierr = DMDestroy(&newdm);CHKERRQ(ierr); 60319c2959aSHong Zhang } 60419c2959aSHong Zhang ierr = DMClone(dm,&newdm);CHKERRQ(ierr); 60519c2959aSHong Zhang ierr = TSGetDM(mprk->subts_slowbuffer,&subdm);CHKERRQ(ierr); 60619c2959aSHong Zhang ierr = DMCopyDMTS(subdm,newdm);CHKERRQ(ierr); 60719c2959aSHong Zhang ierr = DMCopyDMSNES(subdm,newdm);CHKERRQ(ierr); 60819c2959aSHong Zhang ierr = TSSetDM(mprk->subts_slowbuffer,newdm);CHKERRQ(ierr); 60919c2959aSHong Zhang ierr = DMDestroy(&newdm);CHKERRQ(ierr); 61019c2959aSHong Zhang 61119c2959aSHong Zhang ierr = DMClone(dm,&newdm);CHKERRQ(ierr); 61219c2959aSHong Zhang ierr = TSGetDM(mprk->subts_fast,&subdm);CHKERRQ(ierr); 61319c2959aSHong Zhang ierr = DMCopyDMTS(subdm,newdm);CHKERRQ(ierr); 61419c2959aSHong Zhang ierr = DMCopyDMSNES(subdm,newdm);CHKERRQ(ierr); 61519c2959aSHong Zhang ierr = TSSetDM(mprk->subts_fast,newdm);CHKERRQ(ierr); 61619c2959aSHong Zhang ierr = DMDestroy(&newdm);CHKERRQ(ierr); 61719c2959aSHong Zhang 61819c2959aSHong Zhang if (tab->np == 3) { 61919c2959aSHong Zhang ierr = TSRHSSplitGetSubTS(ts,"medium",&mprk->subts_medium);CHKERRQ(ierr); 62019c2959aSHong Zhang ierr = TSRHSSplitGetSubTS(ts,"mediumbuffer",&mprk->subts_mediumbuffer);CHKERRQ(ierr); 62119c2959aSHong Zhang if (mprk->subts_medium && !mprk->subts_mediumbuffer) { 62219c2959aSHong Zhang mprk->subts_mediumbuffer = mprk->subts_medium; 62319c2959aSHong Zhang mprk->subts_medium = NULL; 62419c2959aSHong Zhang } 62519c2959aSHong Zhang if (mprk->subts_medium) { 62619c2959aSHong Zhang ierr = DMClone(dm,&newdm);CHKERRQ(ierr); 62719c2959aSHong Zhang ierr = TSGetDM(mprk->subts_medium,&subdm);CHKERRQ(ierr); 62819c2959aSHong Zhang ierr = DMCopyDMTS(subdm,newdm);CHKERRQ(ierr); 62919c2959aSHong Zhang ierr = DMCopyDMSNES(subdm,newdm);CHKERRQ(ierr); 63019c2959aSHong Zhang ierr = TSSetDM(mprk->subts_medium,newdm);CHKERRQ(ierr); 63119c2959aSHong Zhang ierr = DMDestroy(&newdm);CHKERRQ(ierr); 63219c2959aSHong Zhang } 63319c2959aSHong Zhang ierr = DMClone(dm,&newdm);CHKERRQ(ierr); 63419c2959aSHong Zhang ierr = TSGetDM(mprk->subts_mediumbuffer,&subdm);CHKERRQ(ierr); 63519c2959aSHong Zhang ierr = DMCopyDMTS(subdm,newdm);CHKERRQ(ierr); 63619c2959aSHong Zhang ierr = DMCopyDMSNES(subdm,newdm);CHKERRQ(ierr); 63719c2959aSHong Zhang ierr = TSSetDM(mprk->subts_mediumbuffer,newdm);CHKERRQ(ierr); 63819c2959aSHong Zhang ierr = DMDestroy(&newdm);CHKERRQ(ierr); 63919c2959aSHong Zhang } 6404b84eec9SHong Zhang PetscFunctionReturn(0); 6414b84eec9SHong Zhang } 6424b84eec9SHong Zhang 6434b84eec9SHong Zhang /* 6444b84eec9SHong Zhang This if for nonsplit RHS MPRK 6454b84eec9SHong Zhang The step completion formula is 6464b84eec9SHong Zhang 6474b84eec9SHong Zhang x1 = x0 + h b^T YdotRHS 6484b84eec9SHong Zhang 6494b84eec9SHong Zhang */ 6504b84eec9SHong Zhang static PetscErrorCode TSEvaluateStep_MPRK(TS ts,PetscInt order,Vec X,PetscBool *done) 6514b84eec9SHong Zhang { 6524b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 6534b84eec9SHong Zhang MPRKTableau tab = mprk->tableau; 6547c0df07dSHong Zhang PetscScalar *wf = mprk->work_fast; 6554b84eec9SHong Zhang PetscReal h = ts->time_step; 6564b84eec9SHong Zhang PetscInt s = tab->s,j; 6574b84eec9SHong Zhang PetscErrorCode ierr; 6584b84eec9SHong Zhang 6594b84eec9SHong Zhang PetscFunctionBegin; 6604b84eec9SHong Zhang for (j=0; j<s; j++) wf[j] = h*tab->bf[j]; 6617c0df07dSHong Zhang ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 6627c0df07dSHong Zhang ierr = VecMAXPY(X,s,wf,mprk->YdotRHS);CHKERRQ(ierr); 6634b84eec9SHong Zhang PetscFunctionReturn(0); 6644b84eec9SHong Zhang } 6654b84eec9SHong Zhang 6664b84eec9SHong Zhang static PetscErrorCode TSStep_MPRK(TS ts) 6674b84eec9SHong Zhang { 6684b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 6697c0df07dSHong Zhang Vec *Y = mprk->Y,*YdotRHS = mprk->YdotRHS,*YdotRHS_fast = mprk->YdotRHS_fast,*YdotRHS_slow = mprk->YdotRHS_slow,*YdotRHS_slowbuffer = mprk->YdotRHS_slowbuffer; 6707c0df07dSHong Zhang Vec Yslow,Yslowbuffer,Yfast; 6714b84eec9SHong Zhang MPRKTableau tab = mprk->tableau; 6724b84eec9SHong Zhang const PetscInt s = tab->s; 67319c2959aSHong Zhang const PetscReal *Af = tab->Af,*cf = tab->cf,*Asb = tab->Asb,*csb = tab->csb; 674ebd5ed4eSHong Zhang PetscScalar *wf = mprk->work_fast,*wsb = mprk->work_slowbuffer; 675ebd5ed4eSHong Zhang PetscInt i,j; 6764b84eec9SHong Zhang PetscReal next_time_step = ts->time_step,t = ts->ptime,h = ts->time_step; 6774b84eec9SHong Zhang PetscErrorCode ierr; 6784b84eec9SHong Zhang 6794b84eec9SHong Zhang PetscFunctionBegin; 6804b84eec9SHong Zhang for (i=0; i<s; i++) { 6814b84eec9SHong Zhang mprk->stage_time = t + h*cf[i]; 6824b84eec9SHong Zhang ierr = TSPreStage(ts,mprk->stage_time);CHKERRQ(ierr); 6837c0df07dSHong Zhang ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr); 6844b84eec9SHong Zhang 6857c0df07dSHong Zhang /* slow buffer region */ 68619c2959aSHong Zhang for (j=0; j<i; j++) wsb[j] = h*Asb[i*s+j]; 6879849be05SHong Zhang for (j=0; j<i; j++) { 6887c0df07dSHong Zhang ierr = VecGetSubVector(YdotRHS[j],mprk->is_slowbuffer,&YdotRHS_slowbuffer[j]);CHKERRQ(ierr); 6897c0df07dSHong Zhang } 6907c0df07dSHong Zhang ierr = VecGetSubVector(Y[i],mprk->is_slowbuffer,&Yslowbuffer);CHKERRQ(ierr); 6917c0df07dSHong Zhang ierr = VecMAXPY(Yslowbuffer,i,wsb,mprk->YdotRHS_slowbuffer);CHKERRQ(ierr); 6927c0df07dSHong Zhang ierr = VecRestoreSubVector(Y[i],mprk->is_slowbuffer,&Yslowbuffer);CHKERRQ(ierr); 6939849be05SHong Zhang for (j=0; j<i; j++) { 6947c0df07dSHong Zhang ierr = VecRestoreSubVector(YdotRHS[j],mprk->is_slowbuffer,&YdotRHS_slowbuffer[j]);CHKERRQ(ierr); 6957c0df07dSHong Zhang } 6967c0df07dSHong Zhang /* slow region */ 6977c0df07dSHong Zhang if (mprk->is_slow) { 6989849be05SHong Zhang for (j=0; j<i; j++) { 6997c0df07dSHong Zhang ierr = VecGetSubVector(YdotRHS[j],mprk->is_slow,&YdotRHS_slow[j]);CHKERRQ(ierr); 7007c0df07dSHong Zhang } 7017c0df07dSHong Zhang ierr = VecGetSubVector(Y[i],mprk->is_slow,&Yslow);CHKERRQ(ierr); 702ebd5ed4eSHong Zhang ierr = VecMAXPY(Yslow,i,wsb,mprk->YdotRHS_slow);CHKERRQ(ierr); 7037c0df07dSHong Zhang ierr = VecRestoreSubVector(Y[i],mprk->is_slow,&Yslow);CHKERRQ(ierr); 704ebd5ed4eSHong Zhang for (j=0; j<i; j++) { 7057c0df07dSHong Zhang ierr = VecRestoreSubVector(YdotRHS[j],mprk->is_slow,&YdotRHS_slow[j]);CHKERRQ(ierr); 7067c0df07dSHong Zhang } 7077c0df07dSHong Zhang } 7084b84eec9SHong Zhang 7097c0df07dSHong Zhang /* fast region */ 7104b84eec9SHong Zhang for (j=0; j<i; j++) wf[j] = h*Af[i*s+j]; 7119849be05SHong Zhang for (j=0; j<i; j++) { 7127c0df07dSHong Zhang ierr = VecGetSubVector(YdotRHS[j],mprk->is_fast,&YdotRHS_fast[j]);CHKERRQ(ierr); 7137c0df07dSHong Zhang } 7147c0df07dSHong Zhang ierr = VecGetSubVector(Y[i],mprk->is_fast,&Yfast);CHKERRQ(ierr); 7157c0df07dSHong Zhang ierr = VecMAXPY(Yfast,i,wf,mprk->YdotRHS_fast);CHKERRQ(ierr); 7167c0df07dSHong Zhang ierr = VecRestoreSubVector(Y[i],mprk->is_fast,&Yfast);CHKERRQ(ierr); 7179849be05SHong Zhang for (j=0; j<i; j++) { 7187c0df07dSHong Zhang ierr = VecRestoreSubVector(YdotRHS[j],mprk->is_fast,&YdotRHS_fast[j]);CHKERRQ(ierr); 7197c0df07dSHong Zhang } 7207c0df07dSHong Zhang if (tab->np == 3) { 7217c0df07dSHong Zhang Vec *YdotRHS_medium = mprk->YdotRHS_medium,*YdotRHS_mediumbuffer = mprk->YdotRHS_mediumbuffer; 7227c0df07dSHong Zhang Vec Ymedium,Ymediumbuffer; 723ebd5ed4eSHong Zhang PetscScalar *wmb = mprk->work_mediumbuffer; 724ebd5ed4eSHong Zhang 725ebd5ed4eSHong Zhang for (j=0; j<i; j++) wmb[j] = h*tab->Amb[i*s+j]; 7267c0df07dSHong Zhang /* medium region */ 7277c0df07dSHong Zhang if (mprk->is_medium) { 7287c0df07dSHong Zhang for (j=0; j<i; j++) { 7297c0df07dSHong Zhang ierr = VecGetSubVector(YdotRHS[j],mprk->is_medium,&YdotRHS_medium[j]);CHKERRQ(ierr); 7307c0df07dSHong Zhang } 7317c0df07dSHong Zhang ierr = VecGetSubVector(Y[i],mprk->is_medium,&Ymedium);CHKERRQ(ierr); 732ebd5ed4eSHong Zhang ierr = VecMAXPY(Ymedium,i,wmb,mprk->YdotRHS_medium);CHKERRQ(ierr); 7337c0df07dSHong Zhang ierr = VecRestoreSubVector(Y[i],mprk->is_medium,&Ymedium);CHKERRQ(ierr); 734ebd5ed4eSHong Zhang for (j=0; j<i; j++) { 7357c0df07dSHong Zhang ierr = VecRestoreSubVector(YdotRHS[j],mprk->is_medium,&YdotRHS_medium[j]);CHKERRQ(ierr); 7367c0df07dSHong Zhang } 7377c0df07dSHong Zhang } 7387c0df07dSHong Zhang /* medium buffer region */ 7399849be05SHong Zhang for (j=0; j<i; j++) { 7407c0df07dSHong Zhang ierr = VecGetSubVector(YdotRHS[j],mprk->is_mediumbuffer,&YdotRHS_mediumbuffer[j]);CHKERRQ(ierr); 7417c0df07dSHong Zhang } 7427c0df07dSHong Zhang ierr = VecGetSubVector(Y[i],mprk->is_mediumbuffer,&Ymediumbuffer);CHKERRQ(ierr); 7437c0df07dSHong Zhang ierr = VecMAXPY(Ymediumbuffer,i,wmb,mprk->YdotRHS_mediumbuffer);CHKERRQ(ierr); 7447c0df07dSHong Zhang ierr = VecRestoreSubVector(Y[i],mprk->is_mediumbuffer,&Ymediumbuffer);CHKERRQ(ierr); 7459849be05SHong Zhang for (j=0; j<i; j++) { 7467c0df07dSHong Zhang ierr = VecRestoreSubVector(YdotRHS[j],mprk->is_mediumbuffer,&YdotRHS_mediumbuffer[j]);CHKERRQ(ierr); 7477c0df07dSHong Zhang } 7487c0df07dSHong Zhang } 7494b84eec9SHong Zhang ierr = TSPostStage(ts,mprk->stage_time,i,Y);CHKERRQ(ierr); 7504b84eec9SHong Zhang /* compute the stage RHS by fast and slow tableau respectively */ 7517c0df07dSHong Zhang ierr = TSComputeRHSFunction(ts,t+h*csb[i],Y[i],YdotRHS[i]);CHKERRQ(ierr); 7524b84eec9SHong Zhang } 7534b84eec9SHong Zhang ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr); 7544b84eec9SHong Zhang ts->ptime += ts->time_step; 7554b84eec9SHong Zhang ts->time_step = next_time_step; 7564b84eec9SHong Zhang PetscFunctionReturn(0); 7574b84eec9SHong Zhang } 7584b84eec9SHong Zhang 7594b84eec9SHong Zhang /* 760f944a40eSHong Zhang This if for the case when split RHS is used 7614b84eec9SHong Zhang The step completion formula is 7624b84eec9SHong Zhang x1 = x0 + h b^T YdotRHS 7634b84eec9SHong Zhang */ 7644b84eec9SHong Zhang static PetscErrorCode TSEvaluateStep_MPRKSPLIT(TS ts,PetscInt order,Vec X,PetscBool *done) 7654b84eec9SHong Zhang { 7664b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 7674b84eec9SHong Zhang MPRKTableau tab = mprk->tableau; 76819c2959aSHong Zhang Vec Xslow,Xfast,Xslowbuffer; /* subvectors for slow and fast componets in X respectively */ 76919c2959aSHong Zhang PetscScalar *wf = mprk->work_fast,*ws = mprk->work_slow,*wsb = mprk->work_slowbuffer; 7704b84eec9SHong Zhang PetscReal h = ts->time_step; 77179bc8a77SHong Zhang PetscInt s = tab->s,j,computedstages; 7724b84eec9SHong Zhang PetscErrorCode ierr; 7734b84eec9SHong Zhang 7744b84eec9SHong Zhang PetscFunctionBegin; 7754b84eec9SHong Zhang ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 77619c2959aSHong Zhang 77719c2959aSHong Zhang /* slow region */ 77819c2959aSHong Zhang if (mprk->is_slow) { 77979bc8a77SHong Zhang computedstages = 0; 78019c2959aSHong Zhang for (j=0; j<s; j++) { 7819849be05SHong Zhang if (tab->rsb[j]) ws[tab->rsb[j]-1] += h*tab->bsb[j]; 78279bc8a77SHong Zhang else ws[computedstages++] = h*tab->bsb[j]; 78319c2959aSHong Zhang } 7844b84eec9SHong Zhang ierr = VecGetSubVector(X,mprk->is_slow,&Xslow);CHKERRQ(ierr); 78579bc8a77SHong Zhang ierr = VecMAXPY(Xslow,computedstages,ws,mprk->YdotRHS_slow);CHKERRQ(ierr); 78619c2959aSHong Zhang ierr = VecRestoreSubVector(X,mprk->is_slow,&Xslow);CHKERRQ(ierr); 78719c2959aSHong Zhang } 78819c2959aSHong Zhang 7899d6e09e9SHong Zhang if (tab->np == 3 && mprk->is_medium) { 7909d6e09e9SHong Zhang computedstages = 0; 7919d6e09e9SHong Zhang for (j=0; j<s; j++) { 7929d6e09e9SHong Zhang if (tab->rmb[j]) wsb[computedstages-tab->sbase+(tab->rmb[j]-1)%tab->sbase] += h*tab->bsb[j]; 7939d6e09e9SHong Zhang else wsb[computedstages++] = h*tab->bsb[j]; 7949d6e09e9SHong Zhang } 7959d6e09e9SHong Zhang ierr = VecGetSubVector(X,mprk->is_slowbuffer,&Xslowbuffer);CHKERRQ(ierr); 7969d6e09e9SHong Zhang ierr = VecMAXPY(Xslowbuffer,computedstages,wsb,mprk->YdotRHS_slowbuffer);CHKERRQ(ierr); 7979d6e09e9SHong Zhang ierr = VecRestoreSubVector(X,mprk->is_slowbuffer,&Xslowbuffer);CHKERRQ(ierr); 7989d6e09e9SHong Zhang } else { 79919c2959aSHong Zhang /* slow buffer region */ 80019c2959aSHong Zhang for (j=0; j<s; j++) wsb[j] = h*tab->bsb[j]; 80119c2959aSHong Zhang ierr = VecGetSubVector(X,mprk->is_slowbuffer,&Xslowbuffer);CHKERRQ(ierr); 80219c2959aSHong Zhang ierr = VecMAXPY(Xslowbuffer,s,wsb,mprk->YdotRHS_slowbuffer);CHKERRQ(ierr); 80319c2959aSHong Zhang ierr = VecRestoreSubVector(X,mprk->is_slowbuffer,&Xslowbuffer);CHKERRQ(ierr); 8049d6e09e9SHong Zhang } 80519c2959aSHong Zhang if (tab->np == 3) { 80619c2959aSHong Zhang Vec Xmedium,Xmediumbuffer; 80719c2959aSHong Zhang PetscScalar *wm = mprk->work_medium,*wmb = mprk->work_mediumbuffer; 8089d6e09e9SHong Zhang /* medium region and slow buffer region */ 80919c2959aSHong Zhang if (mprk->is_medium) { 81079bc8a77SHong Zhang computedstages = 0; 81119c2959aSHong Zhang for (j=0; j<s; j++) { 81279bc8a77SHong Zhang if (tab->rmb[j]) wm[computedstages-tab->sbase+(tab->rmb[j]-1)%tab->sbase] += h*tab->bmb[j]; 81379bc8a77SHong Zhang else wm[computedstages++] = h*tab->bmb[j]; 81419c2959aSHong Zhang } 81519c2959aSHong Zhang ierr = VecGetSubVector(X,mprk->is_medium,&Xmedium);CHKERRQ(ierr); 81679bc8a77SHong Zhang ierr = VecMAXPY(Xmedium,computedstages,wm,mprk->YdotRHS_medium);CHKERRQ(ierr); 81719c2959aSHong Zhang ierr = VecRestoreSubVector(X,mprk->is_medium,&Xmedium);CHKERRQ(ierr); 81819c2959aSHong Zhang } 81919c2959aSHong Zhang /* medium buffer region */ 82019c2959aSHong Zhang for (j=0; j<s; j++) wmb[j] = h*tab->bmb[j]; 82119c2959aSHong Zhang ierr = VecGetSubVector(X,mprk->is_mediumbuffer,&Xmediumbuffer);CHKERRQ(ierr); 82219c2959aSHong Zhang ierr = VecMAXPY(Xmediumbuffer,s,wmb,mprk->YdotRHS_mediumbuffer);CHKERRQ(ierr); 82319c2959aSHong Zhang ierr = VecRestoreSubVector(X,mprk->is_mediumbuffer,&Xmediumbuffer);CHKERRQ(ierr); 82419c2959aSHong Zhang } 82519c2959aSHong Zhang /* fast region */ 82619c2959aSHong Zhang for (j=0; j<s; j++) wf[j] = h*tab->bf[j]; 8274b84eec9SHong Zhang ierr = VecGetSubVector(X,mprk->is_fast,&Xfast);CHKERRQ(ierr); 8284b84eec9SHong Zhang ierr = VecMAXPY(Xfast,s,wf,mprk->YdotRHS_fast);CHKERRQ(ierr); 82919c2959aSHong Zhang ierr = VecRestoreSubVector(X,mprk->is_fast,&Xfast);CHKERRQ(ierr); 8304b84eec9SHong Zhang PetscFunctionReturn(0); 8314b84eec9SHong Zhang } 8324b84eec9SHong Zhang 8334b84eec9SHong Zhang static PetscErrorCode TSStep_MPRKSPLIT(TS ts) 8344b84eec9SHong Zhang { 8354b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 8364b84eec9SHong Zhang MPRKTableau tab = mprk->tableau; 83719c2959aSHong Zhang Vec *Y = mprk->Y,*YdotRHS_fast = mprk->YdotRHS_fast,*YdotRHS_slow = mprk->YdotRHS_slow,*YdotRHS_slowbuffer = mprk->YdotRHS_slowbuffer; 83819c2959aSHong Zhang Vec Yslow,Yslowbuffer,Yfast; /* subvectors for slow and fast components in Y[i] respectively */ 83919c2959aSHong Zhang PetscInt s = tab->s; 84019c2959aSHong Zhang const PetscReal *Af = tab->Af,*cf = tab->cf,*Asb = tab->Asb,*csb = tab->csb; 84119c2959aSHong Zhang PetscScalar *wf = mprk->work_fast,*ws = mprk->work_slow,*wsb = mprk->work_slowbuffer; 84279bc8a77SHong Zhang PetscInt i,j,computedstages; 8434b84eec9SHong Zhang PetscReal next_time_step = ts->time_step,t = ts->ptime,h = ts->time_step; 8444b84eec9SHong Zhang PetscErrorCode ierr; 8454b84eec9SHong Zhang 8464b84eec9SHong Zhang PetscFunctionBegin; 8474b84eec9SHong Zhang for (i=0; i<s; i++) { 8484b84eec9SHong Zhang mprk->stage_time = t + h*cf[i]; 8494b84eec9SHong Zhang ierr = TSPreStage(ts,mprk->stage_time);CHKERRQ(ierr); 8504b84eec9SHong Zhang /* calculate the stage value for fast and slow components respectively */ 8514b84eec9SHong Zhang ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr); 85219c2959aSHong Zhang for (j=0; j<i; j++) wsb[j] = h*Asb[i*s+j]; 85319c2959aSHong Zhang 85419c2959aSHong Zhang /* slow buffer region */ 8559d6e09e9SHong Zhang if (tab->np == 3 && mprk->is_medium) { 8569d6e09e9SHong Zhang if (tab->rmb[i]) { 8579d6e09e9SHong Zhang ierr = VecGetSubVector(Y[i],mprk->is_slowbuffer,&Yslowbuffer);CHKERRQ(ierr); 8589d6e09e9SHong Zhang ierr = VecISCopy(Y[tab->rmb[i]-1],mprk->is_slowbuffer,SCATTER_REVERSE,Yslowbuffer);CHKERRQ(ierr); 8599d6e09e9SHong Zhang ierr = VecRestoreSubVector(Y[i],mprk->is_slowbuffer,&Yslowbuffer);CHKERRQ(ierr); 8609d6e09e9SHong Zhang } else { 8619d6e09e9SHong Zhang PetscScalar *wm = mprk->work_medium; 8629d6e09e9SHong Zhang computedstages = 0; 8639d6e09e9SHong Zhang for (j=0; j<i; j++) { 8649d6e09e9SHong Zhang if (tab->rmb[j]) wm[computedstages-tab->sbase+(tab->rmb[j]-1)%tab->sbase] += wsb[j]; 8659d6e09e9SHong Zhang else wm[computedstages++] = wsb[j]; 8669d6e09e9SHong Zhang } 8679d6e09e9SHong Zhang ierr = VecGetSubVector(Y[i],mprk->is_slowbuffer,&Yslowbuffer);CHKERRQ(ierr); 8689d6e09e9SHong Zhang ierr = VecMAXPY(Yslowbuffer,computedstages,wm,YdotRHS_slowbuffer);CHKERRQ(ierr); 8699d6e09e9SHong Zhang ierr = VecRestoreSubVector(Y[i],mprk->is_slowbuffer,&Yslowbuffer);CHKERRQ(ierr); 8709d6e09e9SHong Zhang } 8719d6e09e9SHong Zhang } else { 87219c2959aSHong Zhang ierr = VecGetSubVector(Y[i],mprk->is_slowbuffer,&Yslowbuffer);CHKERRQ(ierr); 87319c2959aSHong Zhang ierr = VecMAXPY(Yslowbuffer,i,wsb,YdotRHS_slowbuffer);CHKERRQ(ierr); 87419c2959aSHong Zhang ierr = VecRestoreSubVector(Y[i],mprk->is_slowbuffer,&Yslowbuffer);CHKERRQ(ierr); 8759d6e09e9SHong Zhang } 8769849be05SHong Zhang 87719c2959aSHong Zhang /* slow region */ 8789849be05SHong Zhang if (mprk->is_slow) { 8799849be05SHong Zhang if (tab->rsb[i]) { /* repeat previous stage */ 8809849be05SHong Zhang ierr = VecGetSubVector(Y[i],mprk->is_slow,&Yslow);CHKERRQ(ierr); 8819849be05SHong Zhang ierr = VecISCopy(Y[tab->rsb[i]-1],mprk->is_slow,SCATTER_REVERSE,Yslow);CHKERRQ(ierr); 8829849be05SHong Zhang ierr = VecRestoreSubVector(Y[i],mprk->is_slow,&Yslow);CHKERRQ(ierr); 8839849be05SHong Zhang } else { 88479bc8a77SHong Zhang computedstages = 0; 88519c2959aSHong Zhang for (j=0; j<i; j++) { 88679bc8a77SHong Zhang if (tab->rsb[j]) ws[tab->rsb[j]-1] += wsb[j]; 88779bc8a77SHong Zhang else ws[computedstages++] = wsb[j]; 88819c2959aSHong Zhang } 8894b84eec9SHong Zhang ierr = VecGetSubVector(Y[i],mprk->is_slow,&Yslow);CHKERRQ(ierr); 89079bc8a77SHong Zhang ierr = VecMAXPY(Yslow,computedstages,ws,YdotRHS_slow);CHKERRQ(ierr); 8914b84eec9SHong Zhang ierr = VecRestoreSubVector(Y[i],mprk->is_slow,&Yslow);CHKERRQ(ierr); 89219c2959aSHong Zhang /* only depends on the slow buffer region */ 89379bc8a77SHong Zhang ierr = TSComputeRHSFunction(mprk->subts_slow,t+h*csb[i],Y[i],YdotRHS_slow[computedstages]);CHKERRQ(ierr); 89419c2959aSHong Zhang } 8959849be05SHong Zhang } 89619c2959aSHong Zhang 89719c2959aSHong Zhang /* fast region */ 89819c2959aSHong Zhang for (j=0; j<i; j++) wf[j] = h*Af[i*s+j]; 89919c2959aSHong Zhang ierr = VecGetSubVector(Y[i],mprk->is_fast,&Yfast);CHKERRQ(ierr); 90019c2959aSHong Zhang ierr = VecMAXPY(Yfast,i,wf,YdotRHS_fast);CHKERRQ(ierr); 9014b84eec9SHong Zhang ierr = VecRestoreSubVector(Y[i],mprk->is_fast,&Yfast);CHKERRQ(ierr); 90219c2959aSHong Zhang 90319c2959aSHong Zhang if (tab->np == 3) { 90419c2959aSHong Zhang Vec *YdotRHS_medium = mprk->YdotRHS_medium,*YdotRHS_mediumbuffer = mprk->YdotRHS_mediumbuffer; 90519c2959aSHong Zhang Vec Ymedium,Ymediumbuffer; 90619c2959aSHong Zhang const PetscReal *Amb = tab->Amb,*cmb = tab->cmb; 90719c2959aSHong Zhang PetscScalar *wm = mprk->work_medium,*wmb = mprk->work_mediumbuffer; 90819c2959aSHong Zhang 90919c2959aSHong Zhang for (j=0; j<i; j++) wmb[j] = h*Amb[i*s+j]; 91019c2959aSHong Zhang /* medium buffer region */ 91119c2959aSHong Zhang ierr = VecGetSubVector(Y[i],mprk->is_mediumbuffer,&Ymediumbuffer);CHKERRQ(ierr); 91219c2959aSHong Zhang ierr = VecMAXPY(Ymediumbuffer,i,wmb,YdotRHS_mediumbuffer);CHKERRQ(ierr); 91319c2959aSHong Zhang ierr = VecRestoreSubVector(Y[i],mprk->is_mediumbuffer,&Ymediumbuffer);CHKERRQ(ierr); 91419c2959aSHong Zhang /* medium region */ 91579bc8a77SHong Zhang if (mprk->is_medium) { 91679bc8a77SHong Zhang if (tab->rmb[i]) { /* repeat previous stage */ 91779bc8a77SHong Zhang ierr = VecGetSubVector(Y[i],mprk->is_medium,&Ymedium);CHKERRQ(ierr); 91879bc8a77SHong Zhang ierr = VecISCopy(Y[tab->rmb[i]-1],mprk->is_medium,SCATTER_REVERSE,Ymedium);CHKERRQ(ierr); 91919c2959aSHong Zhang ierr = VecRestoreSubVector(Y[i],mprk->is_medium,&Ymedium);CHKERRQ(ierr); 92079bc8a77SHong Zhang } else { 92179bc8a77SHong Zhang computedstages = 0; 92279bc8a77SHong Zhang for (j=0; j<i; j++) { 92379bc8a77SHong Zhang if (tab->rmb[j]) wm[computedstages-tab->sbase+(tab->rmb[j]-1)%tab->sbase] += wmb[j]; 92479bc8a77SHong Zhang else wm[computedstages++] = wmb[j]; 92579bc8a77SHong Zhang 92679bc8a77SHong Zhang } 92779bc8a77SHong Zhang ierr = VecGetSubVector(Y[i],mprk->is_medium,&Ymedium);CHKERRQ(ierr); 92879bc8a77SHong Zhang ierr = VecMAXPY(Ymedium,computedstages,wm,YdotRHS_medium);CHKERRQ(ierr); 92979bc8a77SHong Zhang ierr = VecRestoreSubVector(Y[i],mprk->is_medium,&Ymedium);CHKERRQ(ierr); 93079bc8a77SHong Zhang /* only depends on the medium buffer region */ 93179bc8a77SHong Zhang ierr = TSComputeRHSFunction(mprk->subts_medium,t+h*cmb[i],Y[i],YdotRHS_medium[computedstages]);CHKERRQ(ierr); 9329d6e09e9SHong Zhang /* must be computed after all regions are updated in Y */ 9339d6e09e9SHong Zhang ierr = TSComputeRHSFunction(mprk->subts_slowbuffer,t+h*csb[i],Y[i],YdotRHS_slowbuffer[computedstages]);CHKERRQ(ierr); 93479bc8a77SHong Zhang } 93519c2959aSHong Zhang } 93619c2959aSHong Zhang /* must be computed after fast region and slow region are updated in Y */ 93719c2959aSHong Zhang ierr = TSComputeRHSFunction(mprk->subts_mediumbuffer,t+h*cmb[i],Y[i],YdotRHS_mediumbuffer[i]);CHKERRQ(ierr); 93819c2959aSHong Zhang } 9399d6e09e9SHong Zhang if (!(tab->np == 3 && mprk->is_medium)) { 94019c2959aSHong Zhang ierr = TSComputeRHSFunction(mprk->subts_slowbuffer,t+h*csb[i],Y[i],YdotRHS_slowbuffer[i]);CHKERRQ(ierr); 9419d6e09e9SHong Zhang } 942bf0cca7dSHong Zhang ierr = TSComputeRHSFunction(mprk->subts_fast,t+h*cf[i],Y[i],YdotRHS_fast[i]);CHKERRQ(ierr); 9434b84eec9SHong Zhang } 94479bc8a77SHong Zhang 9454b84eec9SHong Zhang ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr); 9464b84eec9SHong Zhang ts->ptime += ts->time_step; 9474b84eec9SHong Zhang ts->time_step = next_time_step; 9484b84eec9SHong Zhang PetscFunctionReturn(0); 9494b84eec9SHong Zhang } 9504b84eec9SHong Zhang 9514b84eec9SHong Zhang static PetscErrorCode TSMPRKTableauReset(TS ts) 9524b84eec9SHong Zhang { 9534b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 9544b84eec9SHong Zhang MPRKTableau tab = mprk->tableau; 9554b84eec9SHong Zhang PetscErrorCode ierr; 9564b84eec9SHong Zhang 9574b84eec9SHong Zhang PetscFunctionBegin; 9584b84eec9SHong Zhang if (!tab) PetscFunctionReturn(0); 9594b84eec9SHong Zhang ierr = PetscFree(mprk->work_fast);CHKERRQ(ierr); 9604b84eec9SHong Zhang ierr = PetscFree(mprk->work_slow);CHKERRQ(ierr); 9617c0df07dSHong Zhang ierr = PetscFree(mprk->work_slowbuffer);CHKERRQ(ierr); 9627c0df07dSHong Zhang ierr = PetscFree(mprk->work_medium);CHKERRQ(ierr); 9637c0df07dSHong Zhang ierr = PetscFree(mprk->work_mediumbuffer);CHKERRQ(ierr); 9644b84eec9SHong Zhang ierr = VecDestroyVecs(tab->s,&mprk->Y);CHKERRQ(ierr); 9650fe4d17eSHong Zhang if (ts->use_splitrhsfunction) { 9660fe4d17eSHong Zhang ierr = VecDestroyVecs(tab->s,&mprk->YdotRHS_fast);CHKERRQ(ierr); 9670fe4d17eSHong Zhang ierr = VecDestroyVecs(tab->s,&mprk->YdotRHS_slow);CHKERRQ(ierr); 9680fe4d17eSHong Zhang ierr = VecDestroyVecs(tab->s,&mprk->YdotRHS_slowbuffer);CHKERRQ(ierr); 9690fe4d17eSHong Zhang ierr = VecDestroyVecs(tab->s,&mprk->YdotRHS_medium);CHKERRQ(ierr); 9700fe4d17eSHong Zhang ierr = VecDestroyVecs(tab->s,&mprk->YdotRHS_mediumbuffer);CHKERRQ(ierr); 9710fe4d17eSHong Zhang } else { 9727c0df07dSHong Zhang ierr = VecDestroyVecs(tab->s,&mprk->YdotRHS);CHKERRQ(ierr); 9737c0df07dSHong Zhang if (mprk->is_slow) { 9747c0df07dSHong Zhang ierr = PetscFree(mprk->YdotRHS_slow);CHKERRQ(ierr); 9754b84eec9SHong Zhang } 9767c0df07dSHong Zhang ierr = PetscFree(mprk->YdotRHS_slowbuffer);CHKERRQ(ierr); 9777c0df07dSHong Zhang if (tab->np == 3) { 9787c0df07dSHong Zhang if (mprk->is_medium) { 9797c0df07dSHong Zhang ierr = PetscFree(mprk->YdotRHS_medium);CHKERRQ(ierr); 9807c0df07dSHong Zhang } 9817c0df07dSHong Zhang ierr = PetscFree(mprk->YdotRHS_mediumbuffer);CHKERRQ(ierr); 9827c0df07dSHong Zhang } 9837c0df07dSHong Zhang ierr = PetscFree(mprk->YdotRHS_fast);CHKERRQ(ierr); 9847c0df07dSHong Zhang } 9854b84eec9SHong Zhang PetscFunctionReturn(0); 9864b84eec9SHong Zhang } 9874b84eec9SHong Zhang 9884b84eec9SHong Zhang static PetscErrorCode TSReset_MPRK(TS ts) 9894b84eec9SHong Zhang { 9904b84eec9SHong Zhang PetscErrorCode ierr; 9914b84eec9SHong Zhang 9924b84eec9SHong Zhang PetscFunctionBegin; 9934b84eec9SHong Zhang ierr = TSMPRKTableauReset(ts);CHKERRQ(ierr); 9944b84eec9SHong Zhang PetscFunctionReturn(0); 9954b84eec9SHong Zhang } 9964b84eec9SHong Zhang 9974b84eec9SHong Zhang static PetscErrorCode DMCoarsenHook_TSMPRK(DM fine,DM coarse,void *ctx) 9984b84eec9SHong Zhang { 9994b84eec9SHong Zhang PetscFunctionBegin; 10004b84eec9SHong Zhang PetscFunctionReturn(0); 10014b84eec9SHong Zhang } 10024b84eec9SHong Zhang 10034b84eec9SHong Zhang static PetscErrorCode DMRestrictHook_TSMPRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 10044b84eec9SHong Zhang { 10054b84eec9SHong Zhang PetscFunctionBegin; 10064b84eec9SHong Zhang PetscFunctionReturn(0); 10074b84eec9SHong Zhang } 10084b84eec9SHong Zhang 10094b84eec9SHong Zhang static PetscErrorCode DMSubDomainHook_TSMPRK(DM dm,DM subdm,void *ctx) 10104b84eec9SHong Zhang { 10114b84eec9SHong Zhang PetscFunctionBegin; 10124b84eec9SHong Zhang PetscFunctionReturn(0); 10134b84eec9SHong Zhang } 10144b84eec9SHong Zhang 10154b84eec9SHong Zhang static PetscErrorCode DMSubDomainRestrictHook_TSMPRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 10164b84eec9SHong Zhang { 10174b84eec9SHong Zhang PetscFunctionBegin; 10184b84eec9SHong Zhang PetscFunctionReturn(0); 10194b84eec9SHong Zhang } 10204b84eec9SHong Zhang 10214b84eec9SHong Zhang static PetscErrorCode TSMPRKTableauSetUp(TS ts) 10224b84eec9SHong Zhang { 10234b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 10244b84eec9SHong Zhang MPRKTableau tab = mprk->tableau; 102519c2959aSHong Zhang Vec YdotRHS_slow,YdotRHS_slowbuffer,YdotRHS_medium,YdotRHS_mediumbuffer,YdotRHS_fast; 10264b84eec9SHong Zhang PetscErrorCode ierr; 10274b84eec9SHong Zhang 10284b84eec9SHong Zhang PetscFunctionBegin; 10294b84eec9SHong Zhang ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&mprk->Y);CHKERRQ(ierr); 10307c0df07dSHong Zhang if (mprk->is_slow) { 103119c2959aSHong Zhang ierr = PetscMalloc1(tab->s,&mprk->work_slow);CHKERRQ(ierr); 10324b84eec9SHong Zhang } 10337c0df07dSHong Zhang ierr = PetscMalloc1(tab->s,&mprk->work_slowbuffer);CHKERRQ(ierr); 10347c0df07dSHong Zhang if (tab->np == 3) { 10357c0df07dSHong Zhang if (mprk->is_medium) { 10367c0df07dSHong Zhang ierr = PetscMalloc1(tab->s,&mprk->work_medium);CHKERRQ(ierr); 10377c0df07dSHong Zhang } 10387c0df07dSHong Zhang ierr = PetscMalloc1(tab->s,&mprk->work_mediumbuffer);CHKERRQ(ierr); 10397c0df07dSHong Zhang } 10407c0df07dSHong Zhang ierr = PetscMalloc1(tab->s,&mprk->work_fast);CHKERRQ(ierr); 10417c0df07dSHong Zhang 10420fe4d17eSHong Zhang if (ts->use_splitrhsfunction) { 104319c2959aSHong Zhang if (mprk->is_slow) { 10444b84eec9SHong Zhang ierr = VecGetSubVector(ts->vec_sol,mprk->is_slow,&YdotRHS_slow);CHKERRQ(ierr); 10454b84eec9SHong Zhang ierr = VecDuplicateVecs(YdotRHS_slow,tab->s,&mprk->YdotRHS_slow);CHKERRQ(ierr); 10464b84eec9SHong Zhang ierr = VecRestoreSubVector(ts->vec_sol,mprk->is_slow,&YdotRHS_slow);CHKERRQ(ierr); 104719c2959aSHong Zhang } 104819c2959aSHong Zhang ierr = VecGetSubVector(ts->vec_sol,mprk->is_slowbuffer,&YdotRHS_slowbuffer);CHKERRQ(ierr); 104919c2959aSHong Zhang ierr = VecDuplicateVecs(YdotRHS_slowbuffer,tab->s,&mprk->YdotRHS_slowbuffer);CHKERRQ(ierr); 105019c2959aSHong Zhang ierr = VecRestoreSubVector(ts->vec_sol,mprk->is_slowbuffer,&YdotRHS_slowbuffer);CHKERRQ(ierr); 105119c2959aSHong Zhang if (tab->np == 3) { 105219c2959aSHong Zhang if (mprk->is_medium) { 105319c2959aSHong Zhang ierr = VecGetSubVector(ts->vec_sol,mprk->is_medium,&YdotRHS_medium);CHKERRQ(ierr); 105419c2959aSHong Zhang ierr = VecDuplicateVecs(YdotRHS_medium,tab->s,&mprk->YdotRHS_medium);CHKERRQ(ierr); 105519c2959aSHong Zhang ierr = VecRestoreSubVector(ts->vec_sol,mprk->is_medium,&YdotRHS_medium);CHKERRQ(ierr); 105619c2959aSHong Zhang } 105719c2959aSHong Zhang ierr = VecGetSubVector(ts->vec_sol,mprk->is_mediumbuffer,&YdotRHS_mediumbuffer);CHKERRQ(ierr); 105819c2959aSHong Zhang ierr = VecDuplicateVecs(YdotRHS_mediumbuffer,tab->s,&mprk->YdotRHS_mediumbuffer);CHKERRQ(ierr); 105919c2959aSHong Zhang ierr = VecRestoreSubVector(ts->vec_sol,mprk->is_mediumbuffer,&YdotRHS_mediumbuffer);CHKERRQ(ierr); 106019c2959aSHong Zhang } 106119c2959aSHong Zhang ierr = VecGetSubVector(ts->vec_sol,mprk->is_fast,&YdotRHS_fast);CHKERRQ(ierr); 106219c2959aSHong Zhang ierr = VecDuplicateVecs(YdotRHS_fast,tab->s,&mprk->YdotRHS_fast);CHKERRQ(ierr); 10634b84eec9SHong Zhang ierr = VecRestoreSubVector(ts->vec_sol,mprk->is_fast,&YdotRHS_fast);CHKERRQ(ierr); 10640fe4d17eSHong Zhang } else { 10650fe4d17eSHong Zhang ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&mprk->YdotRHS);CHKERRQ(ierr); 10660fe4d17eSHong Zhang if (mprk->is_slow) { 10670fe4d17eSHong Zhang ierr = PetscMalloc1(tab->s,&mprk->YdotRHS_slow);CHKERRQ(ierr); 10680fe4d17eSHong Zhang } 10690fe4d17eSHong Zhang ierr = PetscMalloc1(tab->s,&mprk->YdotRHS_slowbuffer);CHKERRQ(ierr); 10700fe4d17eSHong Zhang if (tab->np == 3) { 10710fe4d17eSHong Zhang if (mprk->is_medium) { 10720fe4d17eSHong Zhang ierr = PetscMalloc1(tab->s,&mprk->YdotRHS_medium);CHKERRQ(ierr); 10730fe4d17eSHong Zhang } 10740fe4d17eSHong Zhang ierr = PetscMalloc1(tab->s,&mprk->YdotRHS_mediumbuffer);CHKERRQ(ierr); 10750fe4d17eSHong Zhang } 10760fe4d17eSHong Zhang ierr = PetscMalloc1(tab->s,&mprk->YdotRHS_fast);CHKERRQ(ierr); 10774b84eec9SHong Zhang } 10784b84eec9SHong Zhang PetscFunctionReturn(0); 10794b84eec9SHong Zhang } 10804b84eec9SHong Zhang 10814b84eec9SHong Zhang static PetscErrorCode TSSetUp_MPRK(TS ts) 10824b84eec9SHong Zhang { 10834b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 108419c2959aSHong Zhang MPRKTableau tab = mprk->tableau; 10854b84eec9SHong Zhang DM dm; 10864b84eec9SHong Zhang PetscErrorCode ierr; 10874b84eec9SHong Zhang 10884b84eec9SHong Zhang PetscFunctionBegin; 10894b84eec9SHong Zhang ierr = TSRHSSplitGetIS(ts,"slow",&mprk->is_slow);CHKERRQ(ierr); 10904b84eec9SHong Zhang ierr = TSRHSSplitGetIS(ts,"fast",&mprk->is_fast);CHKERRQ(ierr); 109119c2959aSHong 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); 109219c2959aSHong Zhang 109319c2959aSHong Zhang if (tab->np == 3) { 109419c2959aSHong Zhang ierr = TSRHSSplitGetIS(ts,"medium",&mprk->is_medium);CHKERRQ(ierr); 109519c2959aSHong 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); 109619c2959aSHong Zhang ierr = TSRHSSplitGetIS(ts,"mediumbuffer",&mprk->is_mediumbuffer);CHKERRQ(ierr); 109719c2959aSHong Zhang if (!mprk->is_mediumbuffer) { /* let medium buffer cover whole medium region */ 109819c2959aSHong Zhang mprk->is_mediumbuffer = mprk->is_medium; 109919c2959aSHong Zhang mprk->is_medium = NULL; 110019c2959aSHong Zhang } 110119c2959aSHong Zhang } 110219c2959aSHong Zhang 110319c2959aSHong Zhang /* If users do not provide buffer region settings, the solver will do them automatically, but with a performance penalty */ 110419c2959aSHong Zhang ierr = TSRHSSplitGetIS(ts,"slowbuffer",&mprk->is_slowbuffer);CHKERRQ(ierr); 110519c2959aSHong Zhang if (!mprk->is_slowbuffer) { /* let slow buffer cover whole slow region */ 110619c2959aSHong Zhang mprk->is_slowbuffer = mprk->is_slow; 110719c2959aSHong Zhang mprk->is_slow = NULL; 110819c2959aSHong Zhang } 11094b84eec9SHong Zhang ierr = TSCheckImplicitTerm(ts);CHKERRQ(ierr); 11104b84eec9SHong Zhang ierr = TSMPRKTableauSetUp(ts);CHKERRQ(ierr); 11114b84eec9SHong Zhang ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 11124b84eec9SHong Zhang ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSMPRK,DMRestrictHook_TSMPRK,ts);CHKERRQ(ierr); 11134b84eec9SHong Zhang ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSMPRK,DMSubDomainRestrictHook_TSMPRK,ts);CHKERRQ(ierr); 11140fe4d17eSHong Zhang if (ts->use_splitrhsfunction) { 11150fe4d17eSHong Zhang ts->ops->step = TSStep_MPRKSPLIT; 11160fe4d17eSHong Zhang ts->ops->evaluatestep = TSEvaluateStep_MPRKSPLIT; 11170fe4d17eSHong Zhang ierr = TSMPRKSetSplits(ts);CHKERRQ(ierr); 11180fe4d17eSHong Zhang } else { 11190fe4d17eSHong Zhang ts->ops->step = TSStep_MPRK; 11200fe4d17eSHong Zhang ts->ops->evaluatestep = TSEvaluateStep_MPRK; 11210fe4d17eSHong Zhang } 11224b84eec9SHong Zhang PetscFunctionReturn(0); 11234b84eec9SHong Zhang } 11244b84eec9SHong Zhang 11254b84eec9SHong Zhang static PetscErrorCode TSSetFromOptions_MPRK(PetscOptionItems *PetscOptionsObject,TS ts) 11264b84eec9SHong Zhang { 11274b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 11284b84eec9SHong Zhang PetscErrorCode ierr; 11294b84eec9SHong Zhang 11304b84eec9SHong Zhang PetscFunctionBegin; 11314b84eec9SHong Zhang ierr = PetscOptionsHead(PetscOptionsObject,"PRK ODE solver options");CHKERRQ(ierr); 11324b84eec9SHong Zhang { 11334b84eec9SHong Zhang MPRKTableauLink link; 11344b84eec9SHong Zhang PetscInt count,choice; 11354b84eec9SHong Zhang PetscBool flg; 11364b84eec9SHong Zhang const char **namelist; 11374b84eec9SHong Zhang for (link=MPRKTableauList,count=0; link; link=link->next,count++) ; 11384b84eec9SHong Zhang ierr = PetscMalloc1(count,(char***)&namelist);CHKERRQ(ierr); 11394b84eec9SHong Zhang for (link=MPRKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 11404b84eec9SHong Zhang ierr = PetscOptionsEList("-ts_mprk_type","Family of MPRK method","TSMPRKSetType",(const char*const*)namelist,count,mprk->tableau->name,&choice,&flg);CHKERRQ(ierr); 11414b84eec9SHong Zhang if (flg) {ierr = TSMPRKSetType(ts,namelist[choice]);CHKERRQ(ierr);} 11424b84eec9SHong Zhang ierr = PetscFree(namelist);CHKERRQ(ierr); 11434b84eec9SHong Zhang } 11444b84eec9SHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 11454b84eec9SHong Zhang PetscFunctionReturn(0); 11464b84eec9SHong Zhang } 11474b84eec9SHong Zhang 11484b84eec9SHong Zhang static PetscErrorCode TSView_MPRK(TS ts,PetscViewer viewer) 11494b84eec9SHong Zhang { 11504b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 11514b84eec9SHong Zhang PetscBool iascii; 11524b84eec9SHong Zhang PetscErrorCode ierr; 11534b84eec9SHong Zhang 11544b84eec9SHong Zhang PetscFunctionBegin; 11554b84eec9SHong Zhang ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 11564b84eec9SHong Zhang if (iascii) { 11574b84eec9SHong Zhang MPRKTableau tab = mprk->tableau; 11584b84eec9SHong Zhang TSMPRKType mprktype; 11594b84eec9SHong Zhang char fbuf[512]; 11604b84eec9SHong Zhang char sbuf[512]; 116119c2959aSHong Zhang PetscInt i; 11624b84eec9SHong Zhang ierr = TSMPRKGetType(ts,&mprktype);CHKERRQ(ierr); 11634b84eec9SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," MPRK type %s\n",mprktype);CHKERRQ(ierr); 11644b84eec9SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Order: %D\n",tab->order);CHKERRQ(ierr); 116519c2959aSHong Zhang 11664b84eec9SHong Zhang ierr = PetscFormatRealArray(fbuf,sizeof(fbuf),"% 8.6f",tab->s,tab->cf);CHKERRQ(ierr); 11674b84eec9SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Abscissa cf = %s\n",fbuf);CHKERRQ(ierr); 116819c2959aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Af = \n");CHKERRQ(ierr); 116919c2959aSHong Zhang for (i=0; i<tab->s; i++) { 117019c2959aSHong Zhang ierr = PetscFormatRealArray(fbuf,sizeof(fbuf),"% 8.6f",tab->s,&tab->Af[i*tab->s]);CHKERRQ(ierr); 117119c2959aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," %s\n",fbuf);CHKERRQ(ierr); 117219c2959aSHong Zhang } 117319c2959aSHong Zhang ierr = PetscFormatRealArray(fbuf,sizeof(fbuf),"% 8.6f",tab->s,tab->bf);CHKERRQ(ierr); 117419c2959aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," bf = %s\n",fbuf);CHKERRQ(ierr); 117519c2959aSHong Zhang 117619c2959aSHong Zhang ierr = PetscFormatRealArray(sbuf,sizeof(sbuf),"% 8.6f",tab->s,tab->csb);CHKERRQ(ierr); 117719c2959aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Abscissa csb = %s\n",sbuf);CHKERRQ(ierr); 117819c2959aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Asb = \n");CHKERRQ(ierr); 117919c2959aSHong Zhang for (i=0; i<tab->s; i++) { 118019c2959aSHong Zhang ierr = PetscFormatRealArray(sbuf,sizeof(sbuf),"% 8.6f",tab->s,&tab->Asb[i*tab->s]);CHKERRQ(ierr); 118119c2959aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," %s\n",sbuf);CHKERRQ(ierr); 118219c2959aSHong Zhang } 118319c2959aSHong Zhang ierr = PetscFormatRealArray(sbuf,sizeof(sbuf),"% 8.6f",tab->s,tab->bsb);CHKERRQ(ierr); 118419c2959aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," bsb = %s\n",sbuf);CHKERRQ(ierr); 118519c2959aSHong Zhang 118619c2959aSHong Zhang if (tab->np == 3) { 118719c2959aSHong Zhang char mbuf[512]; 118819c2959aSHong Zhang ierr = PetscFormatRealArray(mbuf,sizeof(mbuf),"% 8.6f",tab->s,tab->cmb);CHKERRQ(ierr); 118919c2959aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Abscissa cmb = %s\n",mbuf);CHKERRQ(ierr); 119019c2959aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Amb = \n");CHKERRQ(ierr); 119119c2959aSHong Zhang for (i=0; i<tab->s; i++) { 119219c2959aSHong Zhang ierr = PetscFormatRealArray(mbuf,sizeof(mbuf),"% 8.6f",tab->s,&tab->Amb[i*tab->s]);CHKERRQ(ierr); 119319c2959aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," %s\n",mbuf);CHKERRQ(ierr); 119419c2959aSHong Zhang } 119519c2959aSHong Zhang ierr = PetscFormatRealArray(mbuf,sizeof(mbuf),"% 8.6f",tab->s,tab->bmb);CHKERRQ(ierr); 119619c2959aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," bmb = %s\n",mbuf);CHKERRQ(ierr); 119719c2959aSHong Zhang } 11984b84eec9SHong Zhang } 11994b84eec9SHong Zhang PetscFunctionReturn(0); 12004b84eec9SHong Zhang } 12014b84eec9SHong Zhang 12024b84eec9SHong Zhang static PetscErrorCode TSLoad_MPRK(TS ts,PetscViewer viewer) 12034b84eec9SHong Zhang { 12044b84eec9SHong Zhang PetscErrorCode ierr; 12054b84eec9SHong Zhang TSAdapt adapt; 12064b84eec9SHong Zhang 12074b84eec9SHong Zhang PetscFunctionBegin; 12084b84eec9SHong Zhang ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 12094b84eec9SHong Zhang ierr = TSAdaptLoad(adapt,viewer);CHKERRQ(ierr); 12104b84eec9SHong Zhang PetscFunctionReturn(0); 12114b84eec9SHong Zhang } 12124b84eec9SHong Zhang 12134b84eec9SHong Zhang /*@C 12144b84eec9SHong Zhang TSMPRKSetType - Set the type of MPRK scheme 12154b84eec9SHong Zhang 12160fe4d17eSHong Zhang Not collective 12174b84eec9SHong Zhang 1218*d8d19677SJose E. Roman Input Parameters: 12194b84eec9SHong Zhang + ts - timestepping context 12204b84eec9SHong Zhang - mprktype - type of MPRK-scheme 12214b84eec9SHong Zhang 12224b84eec9SHong Zhang Options Database: 12234b84eec9SHong Zhang . -ts_mprk_type - <pm2,p2,p3> 12244b84eec9SHong Zhang 12254b84eec9SHong Zhang Level: intermediate 12264b84eec9SHong Zhang 12274b84eec9SHong Zhang .seealso: TSMPRKGetType(), TSMPRK, TSMPRKType 12284b84eec9SHong Zhang @*/ 12294b84eec9SHong Zhang PetscErrorCode TSMPRKSetType(TS ts,TSMPRKType mprktype) 12304b84eec9SHong Zhang { 12314b84eec9SHong Zhang PetscErrorCode ierr; 12324b84eec9SHong Zhang 12334b84eec9SHong Zhang PetscFunctionBegin; 12344b84eec9SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 12354b84eec9SHong Zhang PetscValidCharPointer(mprktype,2); 12364b84eec9SHong Zhang ierr = PetscTryMethod(ts,"TSMPRKSetType_C",(TS,TSMPRKType),(ts,mprktype));CHKERRQ(ierr); 12374b84eec9SHong Zhang PetscFunctionReturn(0); 12384b84eec9SHong Zhang } 12394b84eec9SHong Zhang 12404b84eec9SHong Zhang /*@C 12414b84eec9SHong Zhang TSMPRKGetType - Get the type of MPRK scheme 12424b84eec9SHong Zhang 12430fe4d17eSHong Zhang Not collective 12444b84eec9SHong Zhang 12454b84eec9SHong Zhang Input Parameter: 12464b84eec9SHong Zhang . ts - timestepping context 12474b84eec9SHong Zhang 12484b84eec9SHong Zhang Output Parameter: 12494b84eec9SHong Zhang . mprktype - type of MPRK-scheme 12504b84eec9SHong Zhang 12514b84eec9SHong Zhang Level: intermediate 12524b84eec9SHong Zhang 12534b84eec9SHong Zhang .seealso: TSMPRKGetType() 12544b84eec9SHong Zhang @*/ 12554b84eec9SHong Zhang PetscErrorCode TSMPRKGetType(TS ts,TSMPRKType *mprktype) 12564b84eec9SHong Zhang { 12574b84eec9SHong Zhang PetscErrorCode ierr; 12584b84eec9SHong Zhang 12594b84eec9SHong Zhang PetscFunctionBegin; 12604b84eec9SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 12614b84eec9SHong Zhang ierr = PetscUseMethod(ts,"TSMPRKGetType_C",(TS,TSMPRKType*),(ts,mprktype));CHKERRQ(ierr); 12624b84eec9SHong Zhang PetscFunctionReturn(0); 12634b84eec9SHong Zhang } 12644b84eec9SHong Zhang 12654b84eec9SHong Zhang static PetscErrorCode TSMPRKGetType_MPRK(TS ts,TSMPRKType *mprktype) 12664b84eec9SHong Zhang { 12674b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 12684b84eec9SHong Zhang 12694b84eec9SHong Zhang PetscFunctionBegin; 12704b84eec9SHong Zhang *mprktype = mprk->tableau->name; 12714b84eec9SHong Zhang PetscFunctionReturn(0); 12724b84eec9SHong Zhang } 12734b84eec9SHong Zhang 12744b84eec9SHong Zhang static PetscErrorCode TSMPRKSetType_MPRK(TS ts,TSMPRKType mprktype) 12754b84eec9SHong Zhang { 12764b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 12774b84eec9SHong Zhang PetscBool match; 12784b84eec9SHong Zhang MPRKTableauLink link; 12794b84eec9SHong Zhang PetscErrorCode ierr; 12804b84eec9SHong Zhang 12814b84eec9SHong Zhang PetscFunctionBegin; 12824b84eec9SHong Zhang if (mprk->tableau) { 12834b84eec9SHong Zhang ierr = PetscStrcmp(mprk->tableau->name,mprktype,&match);CHKERRQ(ierr); 12844b84eec9SHong Zhang if (match) PetscFunctionReturn(0); 12854b84eec9SHong Zhang } 12864b84eec9SHong Zhang for (link = MPRKTableauList; link; link=link->next) { 12874b84eec9SHong Zhang ierr = PetscStrcmp(link->tab.name,mprktype,&match);CHKERRQ(ierr); 12884b84eec9SHong Zhang if (match) { 12894b84eec9SHong Zhang if (ts->setupcalled) {ierr = TSMPRKTableauReset(ts);CHKERRQ(ierr);} 12904b84eec9SHong Zhang mprk->tableau = &link->tab; 12914b84eec9SHong Zhang if (ts->setupcalled) {ierr = TSMPRKTableauSetUp(ts);CHKERRQ(ierr);} 12924b84eec9SHong Zhang PetscFunctionReturn(0); 12934b84eec9SHong Zhang } 12944b84eec9SHong Zhang } 12954b84eec9SHong Zhang SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",mprktype); 12964b84eec9SHong Zhang } 12974b84eec9SHong Zhang 12984b84eec9SHong Zhang static PetscErrorCode TSGetStages_MPRK(TS ts,PetscInt *ns,Vec **Y) 12994b84eec9SHong Zhang { 13004b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 13014b84eec9SHong Zhang 13024b84eec9SHong Zhang PetscFunctionBegin; 13034b84eec9SHong Zhang *ns = mprk->tableau->s; 13044b84eec9SHong Zhang if (Y) *Y = mprk->Y; 13054b84eec9SHong Zhang PetscFunctionReturn(0); 13064b84eec9SHong Zhang } 13074b84eec9SHong Zhang 13084b84eec9SHong Zhang static PetscErrorCode TSDestroy_MPRK(TS ts) 13094b84eec9SHong Zhang { 13104b84eec9SHong Zhang PetscErrorCode ierr; 13114b84eec9SHong Zhang 13124b84eec9SHong Zhang PetscFunctionBegin; 13134b84eec9SHong Zhang ierr = TSReset_MPRK(ts);CHKERRQ(ierr); 13144b84eec9SHong Zhang if (ts->dm) { 13154b84eec9SHong Zhang ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSMPRK,DMRestrictHook_TSMPRK,ts);CHKERRQ(ierr); 13164b84eec9SHong Zhang ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSMPRK,DMSubDomainRestrictHook_TSMPRK,ts);CHKERRQ(ierr); 13174b84eec9SHong Zhang } 13184b84eec9SHong Zhang ierr = PetscFree(ts->data);CHKERRQ(ierr); 13194b84eec9SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSMPRKGetType_C",NULL);CHKERRQ(ierr); 13204b84eec9SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSMPRKSetType_C",NULL);CHKERRQ(ierr); 13214b84eec9SHong Zhang PetscFunctionReturn(0); 13224b84eec9SHong Zhang } 13234b84eec9SHong Zhang 13244b84eec9SHong Zhang /*MC 1325f944a40eSHong Zhang TSMPRK - ODE solver using Multirate Partitioned Runge-Kutta schemes 13264b84eec9SHong Zhang 13274b84eec9SHong Zhang The user should provide the right hand side of the equation 13284b84eec9SHong Zhang using TSSetRHSFunction(). 13294b84eec9SHong Zhang 13304b84eec9SHong Zhang Notes: 1331f944a40eSHong Zhang The default is TSMPRKPM2, it can be changed with TSMPRKSetType() or -ts_mprk_type 13324b84eec9SHong Zhang 13334b84eec9SHong Zhang Level: beginner 13344b84eec9SHong Zhang 13354b84eec9SHong Zhang .seealso: TSCreate(), TS, TSSetType(), TSMPRKSetType(), TSMPRKGetType(), TSMPRKType, TSMPRKRegister(), TSMPRKSetMultirateType() 13364b84eec9SHong Zhang TSMPRKM2, TSMPRKM3, TSMPRKRFSMR3, TSMPRKRFSMR2 13374b84eec9SHong Zhang 13384b84eec9SHong Zhang M*/ 13394b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSCreate_MPRK(TS ts) 13404b84eec9SHong Zhang { 13414b84eec9SHong Zhang TS_MPRK *mprk; 13424b84eec9SHong Zhang PetscErrorCode ierr; 13434b84eec9SHong Zhang 13444b84eec9SHong Zhang PetscFunctionBegin; 13454b84eec9SHong Zhang ierr = TSMPRKInitializePackage();CHKERRQ(ierr); 13464b84eec9SHong Zhang 13474b84eec9SHong Zhang ts->ops->reset = TSReset_MPRK; 13484b84eec9SHong Zhang ts->ops->destroy = TSDestroy_MPRK; 13494b84eec9SHong Zhang ts->ops->view = TSView_MPRK; 13504b84eec9SHong Zhang ts->ops->load = TSLoad_MPRK; 13514b84eec9SHong Zhang ts->ops->setup = TSSetUp_MPRK; 13524b84eec9SHong Zhang ts->ops->step = TSStep_MPRK; 13534b84eec9SHong Zhang ts->ops->evaluatestep = TSEvaluateStep_MPRK; 13544b84eec9SHong Zhang ts->ops->setfromoptions = TSSetFromOptions_MPRK; 13554b84eec9SHong Zhang ts->ops->getstages = TSGetStages_MPRK; 13564b84eec9SHong Zhang 13574b84eec9SHong Zhang ierr = PetscNewLog(ts,&mprk);CHKERRQ(ierr); 13584b84eec9SHong Zhang ts->data = (void*)mprk; 13594b84eec9SHong Zhang 13604b84eec9SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSMPRKGetType_C",TSMPRKGetType_MPRK);CHKERRQ(ierr); 13614b84eec9SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSMPRKSetType_C",TSMPRKSetType_MPRK);CHKERRQ(ierr); 13624b84eec9SHong Zhang 13634b84eec9SHong Zhang ierr = TSMPRKSetType(ts,TSMPRKDefault);CHKERRQ(ierr); 13644b84eec9SHong Zhang PetscFunctionReturn(0); 13654b84eec9SHong Zhang } 1366