14b84eec9SHong Zhang /* 24b84eec9SHong Zhang Code for time stepping with the Partitioned Runge-Kutta method 34b84eec9SHong Zhang 44b84eec9SHong Zhang Notes: 54b84eec9SHong Zhang 1) The general system is written as 64b84eec9SHong Zhang Udot = F(t,U) for nonsplit RHS multi-rate RK, 74b84eec9SHong Zhang user should give the indexes for both slow and fast components; 84b84eec9SHong Zhang 2) The general system is written as 94b84eec9SHong Zhang Usdot = Fs(t,Us,Uf) 104b84eec9SHong Zhang Ufdot = Ff(t,Us,Uf) for split RHS multi-rate RK, 114b84eec9SHong Zhang user should partioned RHS by themselves and also provide the indexes for both slow and fast components. 1219c2959aSHong 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'. 1319c2959aSHong 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. 1419c2959aSHong Zhang 1519c2959aSHong Zhang Reference: 1619c2959aSHong Zhang Emil M. Constantinescu, Adrian Sandu, Multirate Timestepping Methods for Hyperbolic Conservation Laws, Journal of Scientific Computing 2007 1719c2959aSHong Zhang 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 */ 314b84eec9SHong Zhang PetscInt s; /* Number of stages */ 3219c2959aSHong Zhang PetscInt np; /* Number of partitions */ 334b84eec9SHong Zhang PetscReal *Af,*bf,*cf; /* Tableau for fast components */ 3419c2959aSHong Zhang PetscReal *Amb,*bmb,*cmb; /* Tableau for medium components */ 3519c2959aSHong Zhang PetscInt *rmb; /* Array of flags for repeated stages in medium method */ 3619c2959aSHong Zhang PetscReal *Asb,*bsb,*csb; /* Tableau for slow components */ 3719c2959aSHong Zhang PetscInt *rsb; /* Array of flags for repeated staged in slow method*/ 384b84eec9SHong Zhang }; 394b84eec9SHong Zhang typedef struct _MPRKTableauLink *MPRKTableauLink; 404b84eec9SHong Zhang struct _MPRKTableauLink { 414b84eec9SHong Zhang struct _MPRKTableau tab; 424b84eec9SHong Zhang MPRKTableauLink next; 434b84eec9SHong Zhang }; 444b84eec9SHong Zhang static MPRKTableauLink MPRKTableauList; 454b84eec9SHong Zhang 464b84eec9SHong Zhang typedef struct { 474b84eec9SHong Zhang MPRKTableau tableau; 484b84eec9SHong Zhang TSMPRKMultirateType mprkmtype; 494b84eec9SHong Zhang Vec *Y; /* States computed during the step */ 50*7c0df07dSHong 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 13119c2959aSHong Zhang TSMPRK2A22 - Second Order 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 14319c2959aSHong Zhang TSMPRK2A23 - Second Order 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 15519c2959aSHong Zhang TSMPRK2A32 - Second Order 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 16719c2959aSHong Zhang TSMPRK2A33 - Second Order 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 17919c2959aSHong Zhang TSMPRK3P2M - Third Order 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 1914b84eec9SHong Zhang TSMPRKP2 - Second Order 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 2034b84eec9SHong Zhang TSMPRKP3 - Third Order 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); 24519c2959aSHong Zhang ierr = TSMPRKRegister(TSMPRK2A22,2,4,&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); 25719c2959aSHong Zhang ierr = TSMPRKRegister(TSMPRK2A23,2,8,&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); 26919c2959aSHong Zhang ierr = TSMPRKRegister(TSMPRK2A32,2,6,&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); 28119c2959aSHong Zhang ierr = TSMPRKRegister(TSMPRK2A33,2,18,&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}; 35219c2959aSHong Zhang ierr = TSMPRKRegister(TSMPRKP2,2,5,&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}; 38119c2959aSHong Zhang ierr = TSMPRKRegister(TSMPRKP3,3,10,&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 4404b84eec9SHong Zhang TSRKFinalizePackage - 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 4664b84eec9SHong Zhang . s - number of stages, this is the dimension of the matrices below 4674b84eec9SHong Zhang . Af - stage coefficients for fast components(dimension s*s, row-major) 4684b84eec9SHong Zhang . bf - step completion table for fast components(dimension s) 4694b84eec9SHong Zhang . cf - abscissa for fast components(dimension s) 4704b84eec9SHong Zhang . As - stage coefficients for slow components(dimension s*s, row-major) 4714b84eec9SHong Zhang . bs - step completion table for slow components(dimension s) 4724b84eec9SHong Zhang - cs - abscissa for slow components(dimension s) 4734b84eec9SHong Zhang 4744b84eec9SHong Zhang Notes: 4754b84eec9SHong Zhang Several MPRK methods are provided, this function is only needed to create new methods. 4764b84eec9SHong Zhang 4774b84eec9SHong Zhang Level: advanced 4784b84eec9SHong Zhang 4794b84eec9SHong Zhang .keywords: TS, register 4804b84eec9SHong Zhang 4814b84eec9SHong Zhang .seealso: TSMPRK 4824b84eec9SHong Zhang @*/ 4834b84eec9SHong Zhang PetscErrorCode TSMPRKRegister(TSMPRKType name,PetscInt order,PetscInt s, 48419c2959aSHong Zhang const PetscReal Asb[],const PetscReal bsb[],const PetscReal csb[],const PetscInt rsb[], 48519c2959aSHong Zhang const PetscReal Amb[],const PetscReal bmb[],const PetscReal cmb[],const PetscInt rmb[], 4864b84eec9SHong Zhang const PetscReal Af[],const PetscReal bf[],const PetscReal cf[]) 4874b84eec9SHong Zhang { 4884b84eec9SHong Zhang MPRKTableauLink link; 4894b84eec9SHong Zhang MPRKTableau t; 4904b84eec9SHong Zhang PetscInt i,j; 4914b84eec9SHong Zhang PetscErrorCode ierr; 4924b84eec9SHong Zhang 4934b84eec9SHong Zhang PetscFunctionBegin; 4944b84eec9SHong Zhang PetscValidCharPointer(name,1); 49519c2959aSHong Zhang PetscValidRealPointer(Asb,4); 49619c2959aSHong Zhang if (bsb) PetscValidRealPointer(bsb,5); 49719c2959aSHong Zhang if (csb) PetscValidRealPointer(csb,6); 49819c2959aSHong Zhang if (rsb) PetscValidRealPointer(rsb,7); 49919c2959aSHong Zhang if (Amb) PetscValidRealPointer(Amb,8); 50019c2959aSHong Zhang if (bmb) PetscValidRealPointer(bmb,9); 50119c2959aSHong Zhang if (cmb) PetscValidRealPointer(cmb,10); 50219c2959aSHong Zhang if (rmb) PetscValidRealPointer(rmb,11); 50319c2959aSHong Zhang PetscValidRealPointer(Af,12); 5044b84eec9SHong Zhang if (bf) PetscValidRealPointer(bf,8); 5054b84eec9SHong Zhang if (cf) PetscValidRealPointer(cf,9); 5064b84eec9SHong Zhang 5074b84eec9SHong Zhang ierr = PetscNew(&link);CHKERRQ(ierr); 5084b84eec9SHong Zhang t = &link->tab; 5094b84eec9SHong Zhang 5104b84eec9SHong Zhang ierr = PetscStrallocpy(name,&t->name);CHKERRQ(ierr); 5114b84eec9SHong Zhang t->order = order; 5124b84eec9SHong Zhang t->s = s; 51319c2959aSHong Zhang t->np = 2; 5144b84eec9SHong Zhang ierr = PetscMalloc3(s*s,&t->Af,s,&t->bf,s,&t->cf);CHKERRQ(ierr); 5154b84eec9SHong Zhang ierr = PetscMemcpy(t->Af,Af,s*s*sizeof(Af[0]));CHKERRQ(ierr); 5164b84eec9SHong Zhang if (bf) { 5174b84eec9SHong Zhang ierr = PetscMemcpy(t->bf,bf,s*sizeof(bf[0]));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) { 5214b84eec9SHong Zhang ierr = PetscMemcpy(t->cf,cf,s*sizeof(cf[0]));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); 53119c2959aSHong Zhang ierr = PetscCalloc1(s,&t->rmb);CHKERRQ(ierr); 53219c2959aSHong Zhang ierr = PetscMemcpy(t->Amb,Amb,s*s*sizeof(Amb[0]));CHKERRQ(ierr); 53319c2959aSHong Zhang if (bmb) { 53419c2959aSHong Zhang ierr = PetscMemcpy(t->bmb,bmb,s*sizeof(bmb[0]));CHKERRQ(ierr); 53519c2959aSHong Zhang } else { 53619c2959aSHong Zhang for (i=0; i<s; i++) t->bmb[i] = Amb[(s-1)*s+i]; 5374b84eec9SHong Zhang } 53819c2959aSHong Zhang if (cmb) { 53919c2959aSHong Zhang ierr = PetscMemcpy(t->cmb,cmb,s*sizeof(cmb[0]));CHKERRQ(ierr); 54019c2959aSHong Zhang } else { 5414b84eec9SHong Zhang for (i=0; i<s; i++) 54219c2959aSHong Zhang for (j=0,t->cmb[i]=0; j<s; j++) 54319c2959aSHong Zhang t->cmb[i] += Amb[i*s+j]; 54419c2959aSHong Zhang } 54519c2959aSHong Zhang if (rmb) { 54619c2959aSHong Zhang ierr = PetscMemcpy(t->rmb,rmb,s*sizeof(rmb[0]));CHKERRQ(ierr); 54719c2959aSHong Zhang } 54819c2959aSHong Zhang } 54919c2959aSHong Zhang 55019c2959aSHong Zhang ierr = PetscMalloc3(s*s,&t->Asb,s,&t->bsb,s,&t->csb);CHKERRQ(ierr); 55119c2959aSHong Zhang ierr = PetscMemcpy(t->Asb,Asb,s*s*sizeof(Asb[0]));CHKERRQ(ierr); 55219c2959aSHong Zhang ierr = PetscCalloc1(s,&t->rsb);CHKERRQ(ierr); 55319c2959aSHong Zhang if (bsb) { 55419c2959aSHong Zhang ierr = PetscMemcpy(t->bsb,bsb,s*sizeof(bsb[0]));CHKERRQ(ierr); 55519c2959aSHong Zhang } else 55619c2959aSHong Zhang for (i=0; i<s; i++) t->bsb[i] = Asb[(s-1)*s+i]; 55719c2959aSHong Zhang if (csb) { 55819c2959aSHong Zhang ierr = PetscMemcpy(t->csb,csb,s*sizeof(csb[0]));CHKERRQ(ierr); 55919c2959aSHong Zhang } else { 56019c2959aSHong Zhang for (i=0; i<s; i++) 56119c2959aSHong Zhang for (j=0,t->csb[i]=0; j<s; j++) 56219c2959aSHong Zhang t->csb[i] += Asb[i*s+j]; 56319c2959aSHong Zhang } 56419c2959aSHong Zhang if (rsb) { 56519c2959aSHong Zhang ierr = PetscMemcpy(t->rsb,rsb,s*sizeof(rsb[0]));CHKERRQ(ierr); 5664b84eec9SHong Zhang } 5674b84eec9SHong Zhang link->next = MPRKTableauList; 5684b84eec9SHong Zhang MPRKTableauList = link; 5694b84eec9SHong Zhang PetscFunctionReturn(0); 5704b84eec9SHong Zhang } 5714b84eec9SHong Zhang 5724b84eec9SHong Zhang static PetscErrorCode TSMPRKSetSplits(TS ts) 5734b84eec9SHong Zhang { 5744b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 57519c2959aSHong Zhang MPRKTableau tab = mprk->tableau; 5764b84eec9SHong Zhang DM dm,subdm,newdm; 5774b84eec9SHong Zhang PetscErrorCode ierr; 5784b84eec9SHong Zhang 5794b84eec9SHong Zhang PetscFunctionBegin; 5804b84eec9SHong Zhang ierr = TSRHSSplitGetSubTS(ts,"slow",&mprk->subts_slow);CHKERRQ(ierr); 5814b84eec9SHong Zhang ierr = TSRHSSplitGetSubTS(ts,"fast",&mprk->subts_fast);CHKERRQ(ierr); 5824b84eec9SHong 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"); 5834b84eec9SHong Zhang 58419c2959aSHong Zhang /* Only copy the DM */ 5854b84eec9SHong Zhang ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 58619c2959aSHong Zhang 58719c2959aSHong Zhang ierr = TSRHSSplitGetSubTS(ts,"slowbuffer",&mprk->subts_slowbuffer);CHKERRQ(ierr); 58819c2959aSHong Zhang if (!mprk->subts_slowbuffer) { 58919c2959aSHong Zhang mprk->subts_slowbuffer = mprk->subts_slow; 59019c2959aSHong Zhang mprk->subts_slow = NULL; 59119c2959aSHong Zhang } 59219c2959aSHong Zhang if (mprk->subts_slow) { 5934b84eec9SHong Zhang ierr = DMClone(dm,&newdm);CHKERRQ(ierr); 5944b84eec9SHong Zhang ierr = TSGetDM(mprk->subts_slow,&subdm);CHKERRQ(ierr); 5954b84eec9SHong Zhang ierr = DMCopyDMTS(subdm,newdm);CHKERRQ(ierr); 5964b84eec9SHong Zhang ierr = DMCopyDMSNES(subdm,newdm);CHKERRQ(ierr); 5974b84eec9SHong Zhang ierr = TSSetDM(mprk->subts_slow,newdm);CHKERRQ(ierr); 5984b84eec9SHong Zhang ierr = DMDestroy(&newdm);CHKERRQ(ierr); 59919c2959aSHong Zhang } 60019c2959aSHong Zhang ierr = DMClone(dm,&newdm);CHKERRQ(ierr); 60119c2959aSHong Zhang ierr = TSGetDM(mprk->subts_slowbuffer,&subdm);CHKERRQ(ierr); 60219c2959aSHong Zhang ierr = DMCopyDMTS(subdm,newdm);CHKERRQ(ierr); 60319c2959aSHong Zhang ierr = DMCopyDMSNES(subdm,newdm);CHKERRQ(ierr); 60419c2959aSHong Zhang ierr = TSSetDM(mprk->subts_slowbuffer,newdm);CHKERRQ(ierr); 60519c2959aSHong Zhang ierr = DMDestroy(&newdm);CHKERRQ(ierr); 60619c2959aSHong Zhang 60719c2959aSHong Zhang ierr = DMClone(dm,&newdm);CHKERRQ(ierr); 60819c2959aSHong Zhang ierr = TSGetDM(mprk->subts_fast,&subdm);CHKERRQ(ierr); 60919c2959aSHong Zhang ierr = DMCopyDMTS(subdm,newdm);CHKERRQ(ierr); 61019c2959aSHong Zhang ierr = DMCopyDMSNES(subdm,newdm);CHKERRQ(ierr); 61119c2959aSHong Zhang ierr = TSSetDM(mprk->subts_fast,newdm);CHKERRQ(ierr); 61219c2959aSHong Zhang ierr = DMDestroy(&newdm);CHKERRQ(ierr); 61319c2959aSHong Zhang 61419c2959aSHong Zhang if (tab->np == 3) { 61519c2959aSHong Zhang ierr = TSRHSSplitGetSubTS(ts,"medium",&mprk->subts_medium);CHKERRQ(ierr); 61619c2959aSHong Zhang ierr = TSRHSSplitGetSubTS(ts,"mediumbuffer",&mprk->subts_mediumbuffer);CHKERRQ(ierr); 61719c2959aSHong Zhang if (mprk->subts_medium && !mprk->subts_mediumbuffer) { 61819c2959aSHong Zhang mprk->subts_mediumbuffer = mprk->subts_medium; 61919c2959aSHong Zhang mprk->subts_medium = NULL; 62019c2959aSHong Zhang } 62119c2959aSHong Zhang if (mprk->subts_medium) { 62219c2959aSHong Zhang ierr = DMClone(dm,&newdm);CHKERRQ(ierr); 62319c2959aSHong Zhang ierr = TSGetDM(mprk->subts_medium,&subdm);CHKERRQ(ierr); 62419c2959aSHong Zhang ierr = DMCopyDMTS(subdm,newdm);CHKERRQ(ierr); 62519c2959aSHong Zhang ierr = DMCopyDMSNES(subdm,newdm);CHKERRQ(ierr); 62619c2959aSHong Zhang ierr = TSSetDM(mprk->subts_medium,newdm);CHKERRQ(ierr); 62719c2959aSHong Zhang ierr = DMDestroy(&newdm);CHKERRQ(ierr); 62819c2959aSHong Zhang } 62919c2959aSHong Zhang ierr = DMClone(dm,&newdm);CHKERRQ(ierr); 63019c2959aSHong Zhang ierr = TSGetDM(mprk->subts_mediumbuffer,&subdm);CHKERRQ(ierr); 63119c2959aSHong Zhang ierr = DMCopyDMTS(subdm,newdm);CHKERRQ(ierr); 63219c2959aSHong Zhang ierr = DMCopyDMSNES(subdm,newdm);CHKERRQ(ierr); 63319c2959aSHong Zhang ierr = TSSetDM(mprk->subts_mediumbuffer,newdm);CHKERRQ(ierr); 63419c2959aSHong Zhang ierr = DMDestroy(&newdm);CHKERRQ(ierr); 63519c2959aSHong Zhang } 6364b84eec9SHong Zhang PetscFunctionReturn(0); 6374b84eec9SHong Zhang } 6384b84eec9SHong Zhang 6394b84eec9SHong Zhang /* 6404b84eec9SHong Zhang This if for nonsplit RHS MPRK 6414b84eec9SHong Zhang The step completion formula is 6424b84eec9SHong Zhang 6434b84eec9SHong Zhang x1 = x0 + h b^T YdotRHS 6444b84eec9SHong Zhang 6454b84eec9SHong Zhang */ 6464b84eec9SHong Zhang static PetscErrorCode TSEvaluateStep_MPRK(TS ts,PetscInt order,Vec X,PetscBool *done) 6474b84eec9SHong Zhang { 6484b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 6494b84eec9SHong Zhang MPRKTableau tab = mprk->tableau; 650*7c0df07dSHong Zhang PetscScalar *wf = mprk->work_fast; 6514b84eec9SHong Zhang PetscReal h = ts->time_step; 6524b84eec9SHong Zhang PetscInt s = tab->s,j; 6534b84eec9SHong Zhang PetscErrorCode ierr; 6544b84eec9SHong Zhang 6554b84eec9SHong Zhang PetscFunctionBegin; 6564b84eec9SHong Zhang for (j=0; j<s; j++) wf[j] = h*tab->bf[j]; 657*7c0df07dSHong Zhang ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 658*7c0df07dSHong Zhang ierr = VecMAXPY(X,s,wf,mprk->YdotRHS);CHKERRQ(ierr); 6594b84eec9SHong Zhang PetscFunctionReturn(0); 6604b84eec9SHong Zhang } 6614b84eec9SHong Zhang 6624b84eec9SHong Zhang static PetscErrorCode TSStep_MPRK(TS ts) 6634b84eec9SHong Zhang { 6644b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 665*7c0df07dSHong Zhang Vec *Y = mprk->Y,*YdotRHS = mprk->YdotRHS,*YdotRHS_fast = mprk->YdotRHS_fast,*YdotRHS_slow = mprk->YdotRHS_slow,*YdotRHS_slowbuffer = mprk->YdotRHS_slowbuffer; 666*7c0df07dSHong Zhang Vec Yslow,Yslowbuffer,Yfast; 6674b84eec9SHong Zhang MPRKTableau tab = mprk->tableau; 6684b84eec9SHong Zhang const PetscInt s = tab->s; 66919c2959aSHong Zhang const PetscReal *Af = tab->Af,*cf = tab->cf,*Asb = tab->Asb,*csb = tab->csb; 670*7c0df07dSHong Zhang PetscScalar *wf = mprk->work_fast,*ws = mprk->work_slow,*wsb = mprk->work_slowbuffer; 671*7c0df07dSHong Zhang PetscInt i,j,basestages; 6724b84eec9SHong Zhang PetscReal next_time_step = ts->time_step,t = ts->ptime,h = ts->time_step; 6734b84eec9SHong Zhang PetscErrorCode ierr; 6744b84eec9SHong Zhang 6754b84eec9SHong Zhang PetscFunctionBegin; 6764b84eec9SHong Zhang for (i=0; i<s; i++) { 6774b84eec9SHong Zhang mprk->stage_time = t + h*cf[i]; 6784b84eec9SHong Zhang ierr = TSPreStage(ts,mprk->stage_time);CHKERRQ(ierr); 679*7c0df07dSHong Zhang ierr = VecCopy(ts->vec_sol, Y[i]);CHKERRQ(ierr); 6804b84eec9SHong Zhang 681*7c0df07dSHong Zhang /* slow buffer region */ 68219c2959aSHong Zhang for (j=0; j<i; j++) wsb[j] = h*Asb[i*s+j]; 683*7c0df07dSHong Zhang for (j=0; j<s; j++) { 684*7c0df07dSHong Zhang ierr = VecGetSubVector(YdotRHS[j],mprk->is_slowbuffer,&YdotRHS_slowbuffer[j]);CHKERRQ(ierr); 685*7c0df07dSHong Zhang } 686*7c0df07dSHong Zhang ierr = VecGetSubVector(Y[i],mprk->is_slowbuffer,&Yslowbuffer);CHKERRQ(ierr); 687*7c0df07dSHong Zhang ierr = VecMAXPY(Yslowbuffer,i,wsb,mprk->YdotRHS_slowbuffer);CHKERRQ(ierr); 688*7c0df07dSHong Zhang ierr = VecRestoreSubVector(Y[i],mprk->is_slowbuffer,&Yslowbuffer);CHKERRQ(ierr); 689*7c0df07dSHong Zhang for (j=0; j<s; j++) { 690*7c0df07dSHong Zhang ierr = VecRestoreSubVector(YdotRHS[j],mprk->is_slowbuffer,&YdotRHS_slowbuffer[j]);CHKERRQ(ierr); 691*7c0df07dSHong Zhang } 692*7c0df07dSHong Zhang /* slow region */ 693*7c0df07dSHong Zhang if (mprk->is_slow) { 694*7c0df07dSHong Zhang basestages = 0; 695*7c0df07dSHong Zhang for (j=0; j<i; j++) 696*7c0df07dSHong Zhang { 697*7c0df07dSHong Zhang if (!tab->rsb[j]) ws[tab->rsb[j]-1] += h*tab->Asb[i*s+j]; 698*7c0df07dSHong Zhang else ws[basestages++] = h*tab->Asb[i*s+j]; 699*7c0df07dSHong Zhang } 700*7c0df07dSHong Zhang for (j=0; j<s; j++) { 701*7c0df07dSHong Zhang ierr = VecGetSubVector(YdotRHS[j],mprk->is_slow,&YdotRHS_slow[j]);CHKERRQ(ierr); 702*7c0df07dSHong Zhang } 703*7c0df07dSHong Zhang ierr = VecGetSubVector(Y[i],mprk->is_slow,&Yslow);CHKERRQ(ierr); 704*7c0df07dSHong Zhang ierr = VecMAXPY(Yslow,i,ws,mprk->YdotRHS_slow);CHKERRQ(ierr); 705*7c0df07dSHong Zhang ierr = VecRestoreSubVector(Y[i],mprk->is_slow,&Yslow);CHKERRQ(ierr); 706*7c0df07dSHong Zhang for (j=0; j<s; j++) { 707*7c0df07dSHong Zhang ierr = VecRestoreSubVector(YdotRHS[j],mprk->is_slow,&YdotRHS_slow[j]);CHKERRQ(ierr); 708*7c0df07dSHong Zhang } 709*7c0df07dSHong Zhang } 7104b84eec9SHong Zhang 711*7c0df07dSHong Zhang /* fast region */ 7124b84eec9SHong Zhang for (j=0; j<i; j++) wf[j] = h*Af[i*s+j]; 713*7c0df07dSHong Zhang for (j=0; j<s; j++) { 714*7c0df07dSHong Zhang ierr = VecGetSubVector(YdotRHS[j],mprk->is_fast,&YdotRHS_fast[j]);CHKERRQ(ierr); 715*7c0df07dSHong Zhang } 716*7c0df07dSHong Zhang ierr = VecGetSubVector(Y[i],mprk->is_fast,&Yfast);CHKERRQ(ierr); 717*7c0df07dSHong Zhang ierr = VecMAXPY(Yfast,i,wf,mprk->YdotRHS_fast);CHKERRQ(ierr); 718*7c0df07dSHong Zhang ierr = VecRestoreSubVector(Y[i],mprk->is_fast,&Yfast);CHKERRQ(ierr); 719*7c0df07dSHong Zhang for (j=0; j<s; j++) { 720*7c0df07dSHong Zhang ierr = VecRestoreSubVector(YdotRHS[j],mprk->is_fast,&YdotRHS_fast[j]);CHKERRQ(ierr); 721*7c0df07dSHong Zhang } 722*7c0df07dSHong Zhang if (tab->np == 3) { 723*7c0df07dSHong Zhang Vec *YdotRHS_medium = mprk->YdotRHS_medium,*YdotRHS_mediumbuffer = mprk->YdotRHS_mediumbuffer; 724*7c0df07dSHong Zhang Vec Ymedium,Ymediumbuffer; 725*7c0df07dSHong Zhang PetscScalar *wm = mprk->work_medium,*wmb = mprk->work_mediumbuffer; 726*7c0df07dSHong Zhang /* medium region */ 727*7c0df07dSHong Zhang if (mprk->is_medium) { 728*7c0df07dSHong Zhang basestages = 0; 729*7c0df07dSHong Zhang for (j=0; j<i; j++) { 730*7c0df07dSHong Zhang if (!tab->rmb[j]) wm[tab->rmb[j]-1] += h*tab->Amb[i*s+j]; 731*7c0df07dSHong Zhang else wm[basestages++] = h*tab->Amb[i*s+j]; 732*7c0df07dSHong Zhang } 733*7c0df07dSHong Zhang for (j=0; j<s; j++) { 734*7c0df07dSHong Zhang ierr = VecGetSubVector(YdotRHS[j],mprk->is_medium,&YdotRHS_medium[j]);CHKERRQ(ierr); 735*7c0df07dSHong Zhang } 736*7c0df07dSHong Zhang ierr = VecGetSubVector(Y[i],mprk->is_medium,&Ymedium);CHKERRQ(ierr); 737*7c0df07dSHong Zhang ierr = VecMAXPY(Ymedium,i,wm,mprk->YdotRHS_medium);CHKERRQ(ierr); 738*7c0df07dSHong Zhang ierr = VecRestoreSubVector(Y[i],mprk->is_medium,&Ymedium);CHKERRQ(ierr); 739*7c0df07dSHong Zhang for (j=0; j<s; j++) { 740*7c0df07dSHong Zhang ierr = VecRestoreSubVector(YdotRHS[j],mprk->is_medium,&YdotRHS_medium[j]);CHKERRQ(ierr); 741*7c0df07dSHong Zhang } 742*7c0df07dSHong Zhang } 743*7c0df07dSHong Zhang /* medium buffer region */ 744*7c0df07dSHong Zhang for (j=0; j<i; j++) wmb[j] = h*tab->Amb[i*s+j]; 745*7c0df07dSHong Zhang for (j=0; j<s; j++) { 746*7c0df07dSHong Zhang ierr = VecGetSubVector(YdotRHS[j],mprk->is_mediumbuffer,&YdotRHS_mediumbuffer[j]);CHKERRQ(ierr); 747*7c0df07dSHong Zhang } 748*7c0df07dSHong Zhang ierr = VecGetSubVector(Y[i],mprk->is_mediumbuffer,&Ymediumbuffer);CHKERRQ(ierr); 749*7c0df07dSHong Zhang ierr = VecMAXPY(Ymediumbuffer,i,wmb,mprk->YdotRHS_mediumbuffer);CHKERRQ(ierr); 750*7c0df07dSHong Zhang ierr = VecRestoreSubVector(Y[i],mprk->is_mediumbuffer,&Ymediumbuffer);CHKERRQ(ierr); 751*7c0df07dSHong Zhang for (j=0; j<s; j++) { 752*7c0df07dSHong Zhang ierr = VecRestoreSubVector(YdotRHS[j],mprk->is_mediumbuffer,&YdotRHS_mediumbuffer[j]);CHKERRQ(ierr); 753*7c0df07dSHong Zhang } 754*7c0df07dSHong Zhang } 7554b84eec9SHong Zhang ierr = TSPostStage(ts,mprk->stage_time,i,Y); CHKERRQ(ierr); 7564b84eec9SHong Zhang /* compute the stage RHS by fast and slow tableau respectively */ 757*7c0df07dSHong Zhang ierr = TSComputeRHSFunction(ts,t+h*csb[i],Y[i],YdotRHS[i]);CHKERRQ(ierr); 7584b84eec9SHong Zhang } 7594b84eec9SHong Zhang ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr); 7604b84eec9SHong Zhang ts->ptime += ts->time_step; 7614b84eec9SHong Zhang ts->time_step = next_time_step; 7624b84eec9SHong Zhang PetscFunctionReturn(0); 7634b84eec9SHong Zhang } 7644b84eec9SHong Zhang 7654b84eec9SHong Zhang /* 7664b84eec9SHong Zhang This if for partitioned RHS MPRK 7674b84eec9SHong Zhang The step completion formula is 7684b84eec9SHong Zhang 7694b84eec9SHong Zhang x1 = x0 + h b^T YdotRHS 7704b84eec9SHong Zhang 7714b84eec9SHong Zhang */ 7724b84eec9SHong Zhang static PetscErrorCode TSEvaluateStep_MPRKSPLIT(TS ts,PetscInt order,Vec X,PetscBool *done) 7734b84eec9SHong Zhang { 7744b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 7754b84eec9SHong Zhang MPRKTableau tab = mprk->tableau; 77619c2959aSHong Zhang Vec Xslow,Xfast,Xslowbuffer; /* subvectors for slow and fast componets in X respectively */ 77719c2959aSHong Zhang PetscScalar *wf = mprk->work_fast,*ws = mprk->work_slow,*wsb = mprk->work_slowbuffer; 7784b84eec9SHong Zhang PetscReal h = ts->time_step; 77919c2959aSHong Zhang PetscInt s = tab->s,j,basestages; 7804b84eec9SHong Zhang PetscErrorCode ierr; 7814b84eec9SHong Zhang 7824b84eec9SHong Zhang PetscFunctionBegin; 7834b84eec9SHong Zhang ierr = VecCopy(ts->vec_sol,X);CHKERRQ(ierr); 78419c2959aSHong Zhang 78519c2959aSHong Zhang /* slow region */ 78619c2959aSHong Zhang if (mprk->is_slow) { 78719c2959aSHong Zhang basestages = 0; 78819c2959aSHong Zhang for (j=0; j<s; j++) { 78919c2959aSHong Zhang if (!tab->rsb[j]) ws[tab->rsb[j]-1] += h*tab->bsb[j]; 79019c2959aSHong Zhang else ws[basestages++] = h*tab->bsb[j]; 79119c2959aSHong Zhang } 7924b84eec9SHong Zhang ierr = VecGetSubVector(X,mprk->is_slow,&Xslow);CHKERRQ(ierr); 79319c2959aSHong Zhang ierr = VecMAXPY(Xslow,basestages,ws,mprk->YdotRHS_slow);CHKERRQ(ierr); 79419c2959aSHong Zhang ierr = VecRestoreSubVector(X,mprk->is_slow,&Xslow);CHKERRQ(ierr); 79519c2959aSHong Zhang } 79619c2959aSHong Zhang 79719c2959aSHong Zhang /* slow buffer region */ 79819c2959aSHong Zhang for (j=0; j<s; j++) wsb[j] = h*tab->bsb[j]; 79919c2959aSHong Zhang ierr = VecGetSubVector(X,mprk->is_slowbuffer,&Xslowbuffer);CHKERRQ(ierr); 80019c2959aSHong Zhang ierr = VecMAXPY(Xslowbuffer,s,wsb,mprk->YdotRHS_slowbuffer);CHKERRQ(ierr); 80119c2959aSHong Zhang ierr = VecRestoreSubVector(X,mprk->is_slowbuffer,&Xslowbuffer);CHKERRQ(ierr); 80219c2959aSHong Zhang 80319c2959aSHong Zhang if (tab->np == 3) { 80419c2959aSHong Zhang Vec Xmedium,Xmediumbuffer; 80519c2959aSHong Zhang PetscScalar *wm = mprk->work_medium,*wmb = mprk->work_mediumbuffer; 80619c2959aSHong Zhang /* medium region */ 80719c2959aSHong Zhang if (mprk->is_medium) { 80819c2959aSHong Zhang basestages = 0; 80919c2959aSHong Zhang for (j=0; j<s; j++) { 81019c2959aSHong Zhang if (!tab->rmb[j]) wm[tab->rmb[j]-1] += h*tab->bmb[j]; 81119c2959aSHong Zhang else wm[basestages++] = h*tab->bmb[j]; 81219c2959aSHong Zhang } 81319c2959aSHong Zhang ierr = VecGetSubVector(X,mprk->is_medium,&Xmedium);CHKERRQ(ierr); 81419c2959aSHong Zhang ierr = VecMAXPY(Xmedium,basestages,wm,mprk->YdotRHS_medium);CHKERRQ(ierr); 81519c2959aSHong Zhang ierr = VecRestoreSubVector(X,mprk->is_medium,&Xmedium);CHKERRQ(ierr); 81619c2959aSHong Zhang } 81719c2959aSHong Zhang /* medium buffer region */ 81819c2959aSHong Zhang for (j=0; j<s; j++) wmb[j] = h*tab->bmb[j]; 81919c2959aSHong Zhang ierr = VecGetSubVector(X,mprk->is_mediumbuffer,&Xmediumbuffer);CHKERRQ(ierr); 82019c2959aSHong Zhang ierr = VecMAXPY(Xmediumbuffer,s,wmb,mprk->YdotRHS_mediumbuffer);CHKERRQ(ierr); 821*7c0df07dSHong Zhang 82219c2959aSHong Zhang ierr = VecRestoreSubVector(X,mprk->is_mediumbuffer,&Xmediumbuffer);CHKERRQ(ierr); 82319c2959aSHong Zhang } 82419c2959aSHong Zhang /* fast region */ 82519c2959aSHong Zhang for (j=0; j<s; j++) wf[j] = h*tab->bf[j]; 8264b84eec9SHong Zhang ierr = VecGetSubVector(X,mprk->is_fast,&Xfast);CHKERRQ(ierr); 8274b84eec9SHong Zhang ierr = VecMAXPY(Xfast,s,wf,mprk->YdotRHS_fast);CHKERRQ(ierr); 82819c2959aSHong Zhang ierr = VecRestoreSubVector(X,mprk->is_fast,&Xfast);CHKERRQ(ierr); 8294b84eec9SHong Zhang PetscFunctionReturn(0); 8304b84eec9SHong Zhang } 8314b84eec9SHong Zhang 8324b84eec9SHong Zhang static PetscErrorCode TSStep_MPRKSPLIT(TS ts) 8334b84eec9SHong Zhang { 8344b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 8354b84eec9SHong Zhang MPRKTableau tab = mprk->tableau; 83619c2959aSHong Zhang Vec *Y = mprk->Y,*YdotRHS_fast = mprk->YdotRHS_fast,*YdotRHS_slow = mprk->YdotRHS_slow,*YdotRHS_slowbuffer = mprk->YdotRHS_slowbuffer; 83719c2959aSHong Zhang Vec Yslow,Yslowbuffer,Yfast; /* subvectors for slow and fast components in Y[i] respectively */ 83819c2959aSHong Zhang PetscInt s = tab->s; 83919c2959aSHong Zhang const PetscReal *Af = tab->Af,*cf = tab->cf,*Asb = tab->Asb,*csb = tab->csb; 84019c2959aSHong Zhang PetscScalar *wf = mprk->work_fast,*ws = mprk->work_slow,*wsb = mprk->work_slowbuffer; 84119c2959aSHong Zhang PetscInt i,j,basestages; 8424b84eec9SHong Zhang PetscReal next_time_step = ts->time_step,t = ts->ptime,h = ts->time_step; 8434b84eec9SHong Zhang PetscErrorCode ierr; 8444b84eec9SHong Zhang 8454b84eec9SHong Zhang PetscFunctionBegin; 8464b84eec9SHong Zhang for (i=0; i<s; i++) { 8474b84eec9SHong Zhang mprk->stage_time = t + h*cf[i]; 8484b84eec9SHong Zhang ierr = TSPreStage(ts,mprk->stage_time);CHKERRQ(ierr); 8494b84eec9SHong Zhang /* calculate the stage value for fast and slow components respectively */ 8504b84eec9SHong Zhang ierr = VecCopy(ts->vec_sol,Y[i]);CHKERRQ(ierr); 85119c2959aSHong Zhang for (j=0; j<i; j++) wsb[j] = h*Asb[i*s+j]; 85219c2959aSHong Zhang 85319c2959aSHong Zhang /* slow buffer region */ 85419c2959aSHong Zhang ierr = VecGetSubVector(Y[i],mprk->is_slowbuffer,&Yslowbuffer);CHKERRQ(ierr); 85519c2959aSHong Zhang ierr = VecMAXPY(Yslowbuffer,i,wsb,YdotRHS_slowbuffer);CHKERRQ(ierr); 85619c2959aSHong Zhang ierr = VecRestoreSubVector(Y[i],mprk->is_slowbuffer,&Yslowbuffer);CHKERRQ(ierr); 85719c2959aSHong Zhang /* slow region */ 85819c2959aSHong Zhang if (!tab->rsb[i] && mprk->is_slow) { /* not a repeated stage */ 85919c2959aSHong Zhang basestages = 0; 86019c2959aSHong Zhang for (j=0; j<i; j++) { 86119c2959aSHong Zhang if (!tab->rsb[j]) ws[tab->rsb[j]-1] += h*Asb[i*s+j]; 86219c2959aSHong Zhang else ws[basestages++] = wsb[j]; 86319c2959aSHong Zhang } 8644b84eec9SHong Zhang ierr = VecGetSubVector(Y[i],mprk->is_slow,&Yslow);CHKERRQ(ierr); 8654b84eec9SHong Zhang ierr = VecMAXPY(Yslow,i,ws,YdotRHS_slow);CHKERRQ(ierr); 8664b84eec9SHong Zhang ierr = VecRestoreSubVector(Y[i],mprk->is_slow,&Yslow);CHKERRQ(ierr); 86719c2959aSHong Zhang /* only depends on the slow buffer region */ 86819c2959aSHong Zhang ierr = TSComputeRHSFunction(mprk->subts_slow,t+h*csb[i],Y[i],YdotRHS_slow[i]);CHKERRQ(ierr); 86919c2959aSHong Zhang } 87019c2959aSHong Zhang 87119c2959aSHong Zhang /* fast region */ 87219c2959aSHong Zhang for (j=0; j<i; j++) wf[j] = h*Af[i*s+j]; 87319c2959aSHong Zhang ierr = VecGetSubVector(Y[i],mprk->is_fast,&Yfast);CHKERRQ(ierr); 87419c2959aSHong Zhang ierr = VecMAXPY(Yfast,i,wf,YdotRHS_fast);CHKERRQ(ierr); 8754b84eec9SHong Zhang ierr = VecRestoreSubVector(Y[i],mprk->is_fast,&Yfast);CHKERRQ(ierr); 87619c2959aSHong Zhang 87719c2959aSHong Zhang if (tab->np == 3) { 87819c2959aSHong Zhang Vec *YdotRHS_medium = mprk->YdotRHS_medium,*YdotRHS_mediumbuffer = mprk->YdotRHS_mediumbuffer; 87919c2959aSHong Zhang Vec Ymedium,Ymediumbuffer; 88019c2959aSHong Zhang const PetscReal *Amb = tab->Amb,*cmb = tab->cmb; 88119c2959aSHong Zhang PetscScalar *wm = mprk->work_medium,*wmb = mprk->work_mediumbuffer; 88219c2959aSHong Zhang 88319c2959aSHong Zhang for (j=0; j<i; j++) wmb[j] = h*Amb[i*s+j]; 88419c2959aSHong Zhang /* medium buffer region */ 88519c2959aSHong Zhang ierr = VecGetSubVector(Y[i],mprk->is_mediumbuffer,&Ymediumbuffer);CHKERRQ(ierr); 88619c2959aSHong Zhang ierr = VecMAXPY(Ymediumbuffer,i,wmb,YdotRHS_mediumbuffer);CHKERRQ(ierr); 88719c2959aSHong Zhang ierr = VecRestoreSubVector(Y[i],mprk->is_mediumbuffer,&Ymediumbuffer);CHKERRQ(ierr); 88819c2959aSHong Zhang /* medium region */ 88919c2959aSHong Zhang if (!tab->rmb[i] && mprk->is_medium) { /* not a repeated stage */ 89019c2959aSHong Zhang basestages = 0; 89119c2959aSHong Zhang for (j=0; j<i; j++) { 89219c2959aSHong Zhang if (tab->rmb[j]) wm[tab->rmb[j]-1] += h*Amb[i*s+j]; 89319c2959aSHong Zhang else wm[basestages++] = wmb[j]; 89419c2959aSHong Zhang } 89519c2959aSHong Zhang ierr = VecGetSubVector(Y[basestages],mprk->is_medium,&Ymedium);CHKERRQ(ierr); 89619c2959aSHong Zhang ierr = VecMAXPY(Ymedium,i,wm,YdotRHS_medium);CHKERRQ(ierr); 89719c2959aSHong Zhang ierr = VecRestoreSubVector(Y[i],mprk->is_medium,&Ymedium);CHKERRQ(ierr); 89819c2959aSHong Zhang /* only depends on the medium buffer region and the slow buffer region */ 89919c2959aSHong Zhang ierr = TSComputeRHSFunction(mprk->subts_medium,t+h*cmb[i],Y[i],YdotRHS_medium[i]);CHKERRQ(ierr); 90019c2959aSHong Zhang } 90119c2959aSHong Zhang /* must be computed after fast region and slow region are updated in Y */ 90219c2959aSHong Zhang ierr = TSComputeRHSFunction(mprk->subts_mediumbuffer,t+h*cmb[i],Y[i],YdotRHS_mediumbuffer[i]);CHKERRQ(ierr); 90319c2959aSHong Zhang } 90419c2959aSHong Zhang /* must be computed after all regions are updated in Y */ 90519c2959aSHong Zhang ierr = TSComputeRHSFunction(mprk->subts_slowbuffer,t+h*csb[i],Y[i],YdotRHS_slowbuffer[i]);CHKERRQ(ierr); 906*7c0df07dSHong Zhang ierr = TSComputeRHSFunction(mprk->subts_fast,t+h*cf[i],Y[i],YdotRHS_fast[i]); 9074b84eec9SHong Zhang } 9084b84eec9SHong Zhang ierr = TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL);CHKERRQ(ierr); 9094b84eec9SHong Zhang ts->ptime += ts->time_step; 9104b84eec9SHong Zhang ts->time_step = next_time_step; 9114b84eec9SHong Zhang PetscFunctionReturn(0); 9124b84eec9SHong Zhang } 9134b84eec9SHong Zhang 9144b84eec9SHong Zhang static PetscErrorCode TSMPRKTableauReset(TS ts) 9154b84eec9SHong Zhang { 9164b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 9174b84eec9SHong Zhang MPRKTableau tab = mprk->tableau; 9184b84eec9SHong Zhang PetscErrorCode ierr; 9194b84eec9SHong Zhang 9204b84eec9SHong Zhang PetscFunctionBegin; 9214b84eec9SHong Zhang if (!tab) PetscFunctionReturn(0); 9224b84eec9SHong Zhang ierr = PetscFree(mprk->work_fast);CHKERRQ(ierr); 9234b84eec9SHong Zhang ierr = PetscFree(mprk->work_slow);CHKERRQ(ierr); 924*7c0df07dSHong Zhang ierr = PetscFree(mprk->work_slowbuffer);CHKERRQ(ierr); 925*7c0df07dSHong Zhang ierr = PetscFree(mprk->work_medium);CHKERRQ(ierr); 926*7c0df07dSHong Zhang ierr = PetscFree(mprk->work_mediumbuffer);CHKERRQ(ierr); 9274b84eec9SHong Zhang ierr = VecDestroyVecs(tab->s,&mprk->Y);CHKERRQ(ierr); 9284b84eec9SHong Zhang if (mprk->mprkmtype == TSMPRKNONSPLIT) { 929*7c0df07dSHong Zhang ierr = VecDestroyVecs(tab->s,&mprk->YdotRHS);CHKERRQ(ierr); 930*7c0df07dSHong Zhang if (mprk->is_slow) { 931*7c0df07dSHong Zhang ierr = PetscFree(mprk->YdotRHS_slow);CHKERRQ(ierr); 9324b84eec9SHong Zhang } 933*7c0df07dSHong Zhang ierr = PetscFree(mprk->YdotRHS_slowbuffer);CHKERRQ(ierr); 934*7c0df07dSHong Zhang if (tab->np == 3) { 935*7c0df07dSHong Zhang if (mprk->is_medium) { 936*7c0df07dSHong Zhang ierr = PetscFree(mprk->YdotRHS_medium);CHKERRQ(ierr); 937*7c0df07dSHong Zhang } 938*7c0df07dSHong Zhang ierr = PetscFree(mprk->YdotRHS_mediumbuffer);CHKERRQ(ierr); 939*7c0df07dSHong Zhang } 940*7c0df07dSHong Zhang ierr = PetscFree(mprk->YdotRHS_fast);CHKERRQ(ierr); 941*7c0df07dSHong Zhang 942*7c0df07dSHong Zhang } else { 9434b84eec9SHong Zhang ierr = VecDestroyVecs(tab->s,&mprk->YdotRHS_fast);CHKERRQ(ierr); 9444b84eec9SHong Zhang ierr = VecDestroyVecs(tab->s,&mprk->YdotRHS_slow);CHKERRQ(ierr); 945*7c0df07dSHong Zhang ierr = VecDestroyVecs(tab->s,&mprk->YdotRHS_slowbuffer);CHKERRQ(ierr); 946*7c0df07dSHong Zhang ierr = VecDestroyVecs(tab->s,&mprk->YdotRHS_medium);CHKERRQ(ierr); 947*7c0df07dSHong Zhang ierr = VecDestroyVecs(tab->s,&mprk->YdotRHS_mediumbuffer);CHKERRQ(ierr); 948*7c0df07dSHong Zhang } 9494b84eec9SHong Zhang PetscFunctionReturn(0); 9504b84eec9SHong Zhang } 9514b84eec9SHong Zhang 9524b84eec9SHong Zhang static PetscErrorCode TSReset_MPRK(TS ts) 9534b84eec9SHong Zhang { 9544b84eec9SHong Zhang PetscErrorCode ierr; 9554b84eec9SHong Zhang 9564b84eec9SHong Zhang PetscFunctionBegin; 9574b84eec9SHong Zhang ierr = TSMPRKTableauReset(ts);CHKERRQ(ierr); 9584b84eec9SHong Zhang PetscFunctionReturn(0); 9594b84eec9SHong Zhang } 9604b84eec9SHong Zhang 9614b84eec9SHong Zhang static PetscErrorCode DMCoarsenHook_TSMPRK(DM fine,DM coarse,void *ctx) 9624b84eec9SHong Zhang { 9634b84eec9SHong Zhang PetscFunctionBegin; 9644b84eec9SHong Zhang PetscFunctionReturn(0); 9654b84eec9SHong Zhang } 9664b84eec9SHong Zhang 9674b84eec9SHong Zhang static PetscErrorCode DMRestrictHook_TSMPRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 9684b84eec9SHong Zhang { 9694b84eec9SHong Zhang PetscFunctionBegin; 9704b84eec9SHong Zhang PetscFunctionReturn(0); 9714b84eec9SHong Zhang } 9724b84eec9SHong Zhang 9734b84eec9SHong Zhang static PetscErrorCode DMSubDomainHook_TSMPRK(DM dm,DM subdm,void *ctx) 9744b84eec9SHong Zhang { 9754b84eec9SHong Zhang PetscFunctionBegin; 9764b84eec9SHong Zhang PetscFunctionReturn(0); 9774b84eec9SHong Zhang } 9784b84eec9SHong Zhang 9794b84eec9SHong Zhang static PetscErrorCode DMSubDomainRestrictHook_TSMPRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 9804b84eec9SHong Zhang { 9814b84eec9SHong Zhang PetscFunctionBegin; 9824b84eec9SHong Zhang PetscFunctionReturn(0); 9834b84eec9SHong Zhang } 9844b84eec9SHong Zhang 9854b84eec9SHong Zhang static PetscErrorCode TSMPRKTableauSetUp(TS ts) 9864b84eec9SHong Zhang { 9874b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 9884b84eec9SHong Zhang MPRKTableau tab = mprk->tableau; 98919c2959aSHong Zhang Vec YdotRHS_slow,YdotRHS_slowbuffer,YdotRHS_medium,YdotRHS_mediumbuffer,YdotRHS_fast; 9904b84eec9SHong Zhang PetscErrorCode ierr; 9914b84eec9SHong Zhang 9924b84eec9SHong Zhang PetscFunctionBegin; 9934b84eec9SHong Zhang ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&mprk->Y);CHKERRQ(ierr); 994*7c0df07dSHong Zhang if (mprk->is_slow) { 99519c2959aSHong Zhang ierr = PetscMalloc1(tab->s,&mprk->work_slow);CHKERRQ(ierr); 9964b84eec9SHong Zhang } 997*7c0df07dSHong Zhang ierr = PetscMalloc1(tab->s,&mprk->work_slowbuffer);CHKERRQ(ierr); 998*7c0df07dSHong Zhang if (tab->np == 3) { 999*7c0df07dSHong Zhang if (mprk->is_medium) { 1000*7c0df07dSHong Zhang ierr = PetscMalloc1(tab->s,&mprk->work_medium);CHKERRQ(ierr); 1001*7c0df07dSHong Zhang } 1002*7c0df07dSHong Zhang ierr = PetscMalloc1(tab->s,&mprk->work_mediumbuffer);CHKERRQ(ierr); 1003*7c0df07dSHong Zhang } 1004*7c0df07dSHong Zhang ierr = PetscMalloc1(tab->s,&mprk->work_fast);CHKERRQ(ierr); 1005*7c0df07dSHong Zhang 1006*7c0df07dSHong Zhang if (mprk->mprkmtype == TSMPRKNONSPLIT) { 1007*7c0df07dSHong Zhang ierr = VecDuplicateVecs(ts->vec_sol,tab->s,&mprk->YdotRHS);CHKERRQ(ierr); 1008*7c0df07dSHong Zhang if (mprk->is_slow) { 1009*7c0df07dSHong Zhang ierr = PetscMalloc1(tab->s,&mprk->YdotRHS_slow);CHKERRQ(ierr); 1010*7c0df07dSHong Zhang } 1011*7c0df07dSHong Zhang ierr = PetscMalloc1(tab->s,&mprk->YdotRHS_slowbuffer);CHKERRQ(ierr); 1012*7c0df07dSHong Zhang if (tab->np == 3) { 1013*7c0df07dSHong Zhang if (mprk->is_medium) { 1014*7c0df07dSHong Zhang ierr = PetscMalloc1(tab->s,&mprk->YdotRHS_medium);CHKERRQ(ierr); 1015*7c0df07dSHong Zhang } 1016*7c0df07dSHong Zhang ierr = PetscMalloc1(tab->s,&mprk->YdotRHS_mediumbuffer);CHKERRQ(ierr); 1017*7c0df07dSHong Zhang } 1018*7c0df07dSHong Zhang ierr = PetscMalloc1(tab->s,&mprk->YdotRHS_fast);CHKERRQ(ierr); 1019*7c0df07dSHong Zhang } else { 102019c2959aSHong Zhang if (mprk->is_slow) { 10214b84eec9SHong Zhang ierr = VecGetSubVector(ts->vec_sol,mprk->is_slow,&YdotRHS_slow);CHKERRQ(ierr); 10224b84eec9SHong Zhang ierr = VecDuplicateVecs(YdotRHS_slow,tab->s,&mprk->YdotRHS_slow);CHKERRQ(ierr); 10234b84eec9SHong Zhang ierr = VecRestoreSubVector(ts->vec_sol,mprk->is_slow,&YdotRHS_slow);CHKERRQ(ierr); 102419c2959aSHong Zhang } 102519c2959aSHong Zhang ierr = VecGetSubVector(ts->vec_sol,mprk->is_slowbuffer,&YdotRHS_slowbuffer);CHKERRQ(ierr); 102619c2959aSHong Zhang ierr = VecDuplicateVecs(YdotRHS_slowbuffer,tab->s,&mprk->YdotRHS_slowbuffer);CHKERRQ(ierr); 102719c2959aSHong Zhang ierr = VecRestoreSubVector(ts->vec_sol,mprk->is_slowbuffer,&YdotRHS_slowbuffer);CHKERRQ(ierr); 102819c2959aSHong Zhang if (tab->np == 3) { 102919c2959aSHong Zhang if (mprk->is_medium) { 103019c2959aSHong Zhang ierr = VecGetSubVector(ts->vec_sol,mprk->is_medium,&YdotRHS_medium);CHKERRQ(ierr); 103119c2959aSHong Zhang ierr = VecDuplicateVecs(YdotRHS_medium,tab->s,&mprk->YdotRHS_medium);CHKERRQ(ierr); 103219c2959aSHong Zhang ierr = VecRestoreSubVector(ts->vec_sol,mprk->is_medium,&YdotRHS_medium);CHKERRQ(ierr); 103319c2959aSHong Zhang } 103419c2959aSHong Zhang ierr = VecGetSubVector(ts->vec_sol,mprk->is_mediumbuffer,&YdotRHS_mediumbuffer);CHKERRQ(ierr); 103519c2959aSHong Zhang ierr = VecDuplicateVecs(YdotRHS_mediumbuffer,tab->s,&mprk->YdotRHS_mediumbuffer);CHKERRQ(ierr); 103619c2959aSHong Zhang ierr = VecRestoreSubVector(ts->vec_sol,mprk->is_mediumbuffer,&YdotRHS_mediumbuffer);CHKERRQ(ierr); 103719c2959aSHong Zhang } 103819c2959aSHong Zhang ierr = VecGetSubVector(ts->vec_sol,mprk->is_fast,&YdotRHS_fast);CHKERRQ(ierr); 103919c2959aSHong Zhang ierr = VecDuplicateVecs(YdotRHS_fast,tab->s,&mprk->YdotRHS_fast);CHKERRQ(ierr); 10404b84eec9SHong Zhang ierr = VecRestoreSubVector(ts->vec_sol,mprk->is_fast,&YdotRHS_fast);CHKERRQ(ierr); 10414b84eec9SHong Zhang } 10424b84eec9SHong Zhang PetscFunctionReturn(0); 10434b84eec9SHong Zhang } 10444b84eec9SHong Zhang 10454b84eec9SHong Zhang static PetscErrorCode TSSetUp_MPRK(TS ts) 10464b84eec9SHong Zhang { 10474b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 104819c2959aSHong Zhang MPRKTableau tab = mprk->tableau; 10494b84eec9SHong Zhang DM dm; 10504b84eec9SHong Zhang PetscErrorCode ierr; 10514b84eec9SHong Zhang 10524b84eec9SHong Zhang PetscFunctionBegin; 10534b84eec9SHong Zhang ierr = TSRHSSplitGetIS(ts,"slow",&mprk->is_slow);CHKERRQ(ierr); 10544b84eec9SHong Zhang ierr = TSRHSSplitGetIS(ts,"fast",&mprk->is_fast);CHKERRQ(ierr); 105519c2959aSHong 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); 105619c2959aSHong Zhang 105719c2959aSHong Zhang if (tab->np == 3) { 105819c2959aSHong Zhang ierr = TSRHSSplitGetIS(ts,"medium",&mprk->is_medium);CHKERRQ(ierr); 105919c2959aSHong 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); 106019c2959aSHong Zhang ierr = TSRHSSplitGetIS(ts,"mediumbuffer",&mprk->is_mediumbuffer);CHKERRQ(ierr); 106119c2959aSHong Zhang if (!mprk->is_mediumbuffer) { /* let medium buffer cover whole medium region */ 106219c2959aSHong Zhang mprk->is_mediumbuffer = mprk->is_medium; 106319c2959aSHong Zhang mprk->is_medium = NULL; 106419c2959aSHong Zhang } 106519c2959aSHong Zhang } 106619c2959aSHong Zhang 106719c2959aSHong Zhang /* If users do not provide buffer region settings, the solver will do them automatically, but with a performance penalty */ 106819c2959aSHong Zhang ierr = TSRHSSplitGetIS(ts,"slowbuffer",&mprk->is_slowbuffer);CHKERRQ(ierr); 106919c2959aSHong Zhang if (!mprk->is_slowbuffer) { /* let slow buffer cover whole slow region */ 107019c2959aSHong Zhang mprk->is_slowbuffer = mprk->is_slow; 107119c2959aSHong Zhang mprk->is_slow = NULL; 107219c2959aSHong Zhang } 107319c2959aSHong Zhang /* 107419c2959aSHong Zhang if (!mprk->is_medium) { 107519c2959aSHong Zhang mprk->is_medium = mprk->is_fast; 107619c2959aSHong Zhang mprk->is_fast = NULL; 107719c2959aSHong Zhang } else { 107819c2959aSHong Zhang ierr = TSRHSSplitGetIS(ts,"mediumbuffer",&mprk->is_mediumbuffer);CHKERRQ(ierr); 107919c2959aSHong Zhang } 108019c2959aSHong Zhang */ 10814b84eec9SHong Zhang ierr = TSCheckImplicitTerm(ts);CHKERRQ(ierr); 10824b84eec9SHong Zhang ierr = TSMPRKTableauSetUp(ts);CHKERRQ(ierr); 10834b84eec9SHong Zhang ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 10844b84eec9SHong Zhang ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSMPRK,DMRestrictHook_TSMPRK,ts);CHKERRQ(ierr); 10854b84eec9SHong Zhang ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSMPRK,DMSubDomainRestrictHook_TSMPRK,ts);CHKERRQ(ierr); 10864b84eec9SHong Zhang ierr = PetscTryMethod(ts,"TSMPRKSetSplits_C",(TS),(ts));CHKERRQ(ierr); 10874b84eec9SHong Zhang PetscFunctionReturn(0); 10884b84eec9SHong Zhang } 10894b84eec9SHong Zhang 10904b84eec9SHong Zhang /* construct a database to chose nonsplit RHS mutirate mprk method or split RHS MPRK method */ 10914b84eec9SHong Zhang const char *const TSMPRKMultirateTypes[] = {"NONSPLIT","SPLIT","TSMPRKMultirateType","TSMPRK",0}; 10924b84eec9SHong Zhang 10934b84eec9SHong Zhang static PetscErrorCode TSSetFromOptions_MPRK(PetscOptionItems *PetscOptionsObject,TS ts) 10944b84eec9SHong Zhang { 10954b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 10964b84eec9SHong Zhang PetscErrorCode ierr; 10974b84eec9SHong Zhang 10984b84eec9SHong Zhang PetscFunctionBegin; 10994b84eec9SHong Zhang ierr = PetscOptionsHead(PetscOptionsObject,"PRK ODE solver options");CHKERRQ(ierr); 11004b84eec9SHong Zhang { 11014b84eec9SHong Zhang MPRKTableauLink link; 11024b84eec9SHong Zhang PetscInt count,choice; 11034b84eec9SHong Zhang PetscBool flg; 11044b84eec9SHong Zhang const char **namelist; 11054b84eec9SHong Zhang PetscInt mprkmtype = 0; 11064b84eec9SHong Zhang for (link=MPRKTableauList,count=0; link; link=link->next,count++) ; 11074b84eec9SHong Zhang ierr = PetscMalloc1(count,(char***)&namelist);CHKERRQ(ierr); 11084b84eec9SHong Zhang for (link=MPRKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 11094b84eec9SHong Zhang ierr = PetscOptionsEList("-ts_mprk_type","Family of MPRK method","TSMPRKSetType",(const char*const*)namelist,count,mprk->tableau->name,&choice,&flg);CHKERRQ(ierr); 11104b84eec9SHong Zhang if (flg) {ierr = TSMPRKSetType(ts,namelist[choice]);CHKERRQ(ierr);} 11114b84eec9SHong Zhang ierr = PetscFree(namelist);CHKERRQ(ierr); 11124b84eec9SHong Zhang ierr = PetscOptionsEList("-ts_mprk_multirate_type","Use Combined RHS Multirate or Partioned RHS Multirat MPRK method","TSMPRKSetMultirateType",TSMPRKMultirateTypes,2,TSMPRKMultirateTypes[0],&mprkmtype,&flg);CHKERRQ(ierr); 11134b84eec9SHong Zhang if (flg) {ierr = TSMPRKSetMultirateType(ts,mprkmtype);CHKERRQ(ierr);} 11144b84eec9SHong Zhang } 11154b84eec9SHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 11164b84eec9SHong Zhang PetscFunctionReturn(0); 11174b84eec9SHong Zhang } 11184b84eec9SHong Zhang 11194b84eec9SHong Zhang static PetscErrorCode TSView_MPRK(TS ts,PetscViewer viewer) 11204b84eec9SHong Zhang { 11214b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 11224b84eec9SHong Zhang PetscBool iascii; 11234b84eec9SHong Zhang PetscErrorCode ierr; 11244b84eec9SHong Zhang 11254b84eec9SHong Zhang PetscFunctionBegin; 11264b84eec9SHong Zhang ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 11274b84eec9SHong Zhang if (iascii) { 11284b84eec9SHong Zhang MPRKTableau tab = mprk->tableau; 11294b84eec9SHong Zhang TSMPRKType mprktype; 11304b84eec9SHong Zhang char fbuf[512]; 11314b84eec9SHong Zhang char sbuf[512]; 113219c2959aSHong Zhang PetscInt i; 11334b84eec9SHong Zhang ierr = TSMPRKGetType(ts,&mprktype);CHKERRQ(ierr); 11344b84eec9SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," MPRK type %s\n",mprktype);CHKERRQ(ierr); 11354b84eec9SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Order: %D\n",tab->order);CHKERRQ(ierr); 113619c2959aSHong Zhang 11374b84eec9SHong Zhang ierr = PetscFormatRealArray(fbuf,sizeof(fbuf),"% 8.6f",tab->s,tab->cf);CHKERRQ(ierr); 11384b84eec9SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Abscissa cf = %s\n",fbuf);CHKERRQ(ierr); 113919c2959aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Af = \n");CHKERRQ(ierr); 114019c2959aSHong Zhang for (i=0; i<tab->s; i++) { 114119c2959aSHong Zhang ierr = PetscFormatRealArray(fbuf,sizeof(fbuf),"% 8.6f",tab->s,&tab->Af[i*tab->s]);CHKERRQ(ierr); 114219c2959aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," %s\n",fbuf);CHKERRQ(ierr); 114319c2959aSHong Zhang } 114419c2959aSHong Zhang ierr = PetscFormatRealArray(fbuf,sizeof(fbuf),"% 8.6f",tab->s,tab->bf);CHKERRQ(ierr); 114519c2959aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," bf = %s\n",fbuf);CHKERRQ(ierr); 114619c2959aSHong Zhang 114719c2959aSHong Zhang ierr = PetscFormatRealArray(sbuf,sizeof(sbuf),"% 8.6f",tab->s,tab->csb);CHKERRQ(ierr); 114819c2959aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Abscissa csb = %s\n",sbuf);CHKERRQ(ierr); 114919c2959aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Asb = \n");CHKERRQ(ierr); 115019c2959aSHong Zhang for (i=0; i<tab->s; i++) { 115119c2959aSHong Zhang ierr = PetscFormatRealArray(sbuf,sizeof(sbuf),"% 8.6f",tab->s,&tab->Asb[i*tab->s]);CHKERRQ(ierr); 115219c2959aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," %s\n",sbuf);CHKERRQ(ierr); 115319c2959aSHong Zhang } 115419c2959aSHong Zhang ierr = PetscFormatRealArray(sbuf,sizeof(sbuf),"% 8.6f",tab->s,tab->bsb);CHKERRQ(ierr); 115519c2959aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," bsb = %s\n",sbuf);CHKERRQ(ierr); 115619c2959aSHong Zhang 115719c2959aSHong Zhang if (tab->np == 3) { 115819c2959aSHong Zhang char mbuf[512]; 115919c2959aSHong Zhang ierr = PetscFormatRealArray(mbuf,sizeof(mbuf),"% 8.6f",tab->s,tab->cmb);CHKERRQ(ierr); 116019c2959aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Abscissa cmb = %s\n",mbuf);CHKERRQ(ierr); 116119c2959aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Amb = \n");CHKERRQ(ierr); 116219c2959aSHong Zhang for (i=0; i<tab->s; i++) { 116319c2959aSHong Zhang ierr = PetscFormatRealArray(mbuf,sizeof(mbuf),"% 8.6f",tab->s,&tab->Amb[i*tab->s]);CHKERRQ(ierr); 116419c2959aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," %s\n",mbuf);CHKERRQ(ierr); 116519c2959aSHong Zhang } 116619c2959aSHong Zhang ierr = PetscFormatRealArray(mbuf,sizeof(mbuf),"% 8.6f",tab->s,tab->bmb);CHKERRQ(ierr); 116719c2959aSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," bmb = %s\n",mbuf);CHKERRQ(ierr); 116819c2959aSHong Zhang } 11694b84eec9SHong Zhang } 11704b84eec9SHong Zhang PetscFunctionReturn(0); 11714b84eec9SHong Zhang } 11724b84eec9SHong Zhang 11734b84eec9SHong Zhang static PetscErrorCode TSLoad_MPRK(TS ts,PetscViewer viewer) 11744b84eec9SHong Zhang { 11754b84eec9SHong Zhang PetscErrorCode ierr; 11764b84eec9SHong Zhang TSAdapt adapt; 11774b84eec9SHong Zhang 11784b84eec9SHong Zhang PetscFunctionBegin; 11794b84eec9SHong Zhang ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 11804b84eec9SHong Zhang ierr = TSAdaptLoad(adapt,viewer);CHKERRQ(ierr); 11814b84eec9SHong Zhang PetscFunctionReturn(0); 11824b84eec9SHong Zhang } 11834b84eec9SHong Zhang 11844b84eec9SHong Zhang /*@C 11854b84eec9SHong Zhang TSMPRKSetType - Set the type of MPRK scheme 11864b84eec9SHong Zhang 11874b84eec9SHong Zhang Logically collective 11884b84eec9SHong Zhang 11894b84eec9SHong Zhang Input Parameter: 11904b84eec9SHong Zhang + ts - timestepping context 11914b84eec9SHong Zhang - mprktype - type of MPRK-scheme 11924b84eec9SHong Zhang 11934b84eec9SHong Zhang Options Database: 11944b84eec9SHong Zhang . -ts_mprk_type - <pm2,p2,p3> 11954b84eec9SHong Zhang 11964b84eec9SHong Zhang Level: intermediate 11974b84eec9SHong Zhang 11984b84eec9SHong Zhang .seealso: TSMPRKGetType(), TSMPRK, TSMPRKType 11994b84eec9SHong Zhang @*/ 12004b84eec9SHong Zhang PetscErrorCode TSMPRKSetType(TS ts,TSMPRKType mprktype) 12014b84eec9SHong Zhang { 12024b84eec9SHong Zhang PetscErrorCode ierr; 12034b84eec9SHong Zhang 12044b84eec9SHong Zhang PetscFunctionBegin; 12054b84eec9SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 12064b84eec9SHong Zhang PetscValidCharPointer(mprktype,2); 12074b84eec9SHong Zhang ierr = PetscTryMethod(ts,"TSMPRKSetType_C",(TS,TSMPRKType),(ts,mprktype));CHKERRQ(ierr); 12084b84eec9SHong Zhang PetscFunctionReturn(0); 12094b84eec9SHong Zhang } 12104b84eec9SHong Zhang 12114b84eec9SHong Zhang /*@C 12124b84eec9SHong Zhang TSMPRKGetType - Get the type of MPRK scheme 12134b84eec9SHong Zhang 12144b84eec9SHong Zhang Logically collective 12154b84eec9SHong Zhang 12164b84eec9SHong Zhang Input Parameter: 12174b84eec9SHong Zhang . ts - timestepping context 12184b84eec9SHong Zhang 12194b84eec9SHong Zhang Output Parameter: 12204b84eec9SHong Zhang . mprktype - type of MPRK-scheme 12214b84eec9SHong Zhang 12224b84eec9SHong Zhang Level: intermediate 12234b84eec9SHong Zhang 12244b84eec9SHong Zhang .seealso: TSMPRKGetType() 12254b84eec9SHong Zhang @*/ 12264b84eec9SHong Zhang PetscErrorCode TSMPRKGetType(TS ts,TSMPRKType *mprktype) 12274b84eec9SHong Zhang { 12284b84eec9SHong Zhang PetscErrorCode ierr; 12294b84eec9SHong Zhang 12304b84eec9SHong Zhang PetscFunctionBegin; 12314b84eec9SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 12324b84eec9SHong Zhang ierr = PetscUseMethod(ts,"TSMPRKGetType_C",(TS,TSMPRKType*),(ts,mprktype));CHKERRQ(ierr); 12334b84eec9SHong Zhang PetscFunctionReturn(0); 12344b84eec9SHong Zhang } 12354b84eec9SHong Zhang 12364b84eec9SHong Zhang /*@C 12374b84eec9SHong Zhang TSMPRKSetMultirateType - Set the type of MPRK multirate scheme 12384b84eec9SHong Zhang 12394b84eec9SHong Zhang Logically collective 12404b84eec9SHong Zhang 12414b84eec9SHong Zhang Input Parameter: 12424b84eec9SHong Zhang + ts - timestepping context 12434b84eec9SHong Zhang - mprkmtype - type of the multirate configuration 12444b84eec9SHong Zhang 12454b84eec9SHong Zhang Options Database: 12464b84eec9SHong Zhang . -ts_mprk_multirate_type - <nonsplit,split> 12474b84eec9SHong Zhang 12484b84eec9SHong Zhang Level: intermediate 12494b84eec9SHong Zhang @*/ 12504b84eec9SHong Zhang PetscErrorCode TSMPRKSetMultirateType(TS ts, TSMPRKMultirateType mprkmtype) 12514b84eec9SHong Zhang { 12524b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 12534b84eec9SHong Zhang PetscErrorCode ierr; 12544b84eec9SHong Zhang 12554b84eec9SHong Zhang PetscFunctionBegin; 12564b84eec9SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 12574b84eec9SHong Zhang switch(mprkmtype){ 12584b84eec9SHong Zhang case TSMPRKNONSPLIT: 12594b84eec9SHong Zhang ts->ops->step = TSStep_MPRK; 12604b84eec9SHong Zhang ts->ops->evaluatestep = TSEvaluateStep_MPRK; 12614b84eec9SHong Zhang break; 12624b84eec9SHong Zhang case TSMPRKSPLIT: 12634b84eec9SHong Zhang ts->ops->step = TSStep_MPRKSPLIT; 12644b84eec9SHong Zhang ts->ops->evaluatestep = TSEvaluateStep_MPRKSPLIT; 12654b84eec9SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSMPRKSetSplits_C",TSMPRKSetSplits);CHKERRQ(ierr); 12664b84eec9SHong Zhang break; 12674b84eec9SHong Zhang default : 12684b84eec9SHong Zhang SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown type '%s'",mprkmtype); 12694b84eec9SHong Zhang } 12704b84eec9SHong Zhang mprk->mprkmtype = mprkmtype; 12714b84eec9SHong Zhang PetscFunctionReturn(0); 12724b84eec9SHong Zhang } 12734b84eec9SHong Zhang 12744b84eec9SHong Zhang static PetscErrorCode TSMPRKGetType_MPRK(TS ts,TSMPRKType *mprktype) 12754b84eec9SHong Zhang { 12764b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 12774b84eec9SHong Zhang 12784b84eec9SHong Zhang PetscFunctionBegin; 12794b84eec9SHong Zhang *mprktype = mprk->tableau->name; 12804b84eec9SHong Zhang PetscFunctionReturn(0); 12814b84eec9SHong Zhang } 12824b84eec9SHong Zhang 12834b84eec9SHong Zhang static PetscErrorCode TSMPRKSetType_MPRK(TS ts,TSMPRKType mprktype) 12844b84eec9SHong Zhang { 12854b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 12864b84eec9SHong Zhang PetscBool match; 12874b84eec9SHong Zhang MPRKTableauLink link; 12884b84eec9SHong Zhang PetscErrorCode ierr; 12894b84eec9SHong Zhang 12904b84eec9SHong Zhang PetscFunctionBegin; 12914b84eec9SHong Zhang if (mprk->tableau) { 12924b84eec9SHong Zhang ierr = PetscStrcmp(mprk->tableau->name,mprktype,&match);CHKERRQ(ierr); 12934b84eec9SHong Zhang if (match) PetscFunctionReturn(0); 12944b84eec9SHong Zhang } 12954b84eec9SHong Zhang for (link = MPRKTableauList; link; link=link->next) { 12964b84eec9SHong Zhang ierr = PetscStrcmp(link->tab.name,mprktype,&match);CHKERRQ(ierr); 12974b84eec9SHong Zhang if (match) { 12984b84eec9SHong Zhang if (ts->setupcalled) {ierr = TSMPRKTableauReset(ts);CHKERRQ(ierr);} 12994b84eec9SHong Zhang mprk->tableau = &link->tab; 13004b84eec9SHong Zhang if (ts->setupcalled) {ierr = TSMPRKTableauSetUp(ts);CHKERRQ(ierr);} 13014b84eec9SHong Zhang PetscFunctionReturn(0); 13024b84eec9SHong Zhang } 13034b84eec9SHong Zhang } 13044b84eec9SHong Zhang SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",mprktype); 13054b84eec9SHong Zhang PetscFunctionReturn(0); 13064b84eec9SHong Zhang } 13074b84eec9SHong Zhang 13084b84eec9SHong Zhang static PetscErrorCode TSGetStages_MPRK(TS ts,PetscInt *ns,Vec **Y) 13094b84eec9SHong Zhang { 13104b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 13114b84eec9SHong Zhang 13124b84eec9SHong Zhang PetscFunctionBegin; 13134b84eec9SHong Zhang *ns = mprk->tableau->s; 13144b84eec9SHong Zhang if (Y) *Y = mprk->Y; 13154b84eec9SHong Zhang PetscFunctionReturn(0); 13164b84eec9SHong Zhang } 13174b84eec9SHong Zhang 13184b84eec9SHong Zhang static PetscErrorCode TSDestroy_MPRK(TS ts) 13194b84eec9SHong Zhang { 13204b84eec9SHong Zhang PetscErrorCode ierr; 13214b84eec9SHong Zhang 13224b84eec9SHong Zhang PetscFunctionBegin; 13234b84eec9SHong Zhang ierr = TSReset_MPRK(ts);CHKERRQ(ierr); 13244b84eec9SHong Zhang if (ts->dm) { 13254b84eec9SHong Zhang ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSMPRK,DMRestrictHook_TSMPRK,ts);CHKERRQ(ierr); 13264b84eec9SHong Zhang ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSMPRK,DMSubDomainRestrictHook_TSMPRK,ts);CHKERRQ(ierr); 13274b84eec9SHong Zhang } 13284b84eec9SHong Zhang ierr = PetscFree(ts->data);CHKERRQ(ierr); 13294b84eec9SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSMPRKGetType_C",NULL);CHKERRQ(ierr); 13304b84eec9SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSMPRKSetType_C",NULL);CHKERRQ(ierr); 13314b84eec9SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSMPRKSetMultirateType_C",NULL);CHKERRQ(ierr); 13324b84eec9SHong Zhang PetscFunctionReturn(0); 13334b84eec9SHong Zhang } 13344b84eec9SHong Zhang 13354b84eec9SHong Zhang /*MC 13364b84eec9SHong Zhang TSMPRK - ODE solver using Partitioned Runge-Kutta schemes 13374b84eec9SHong Zhang 13384b84eec9SHong Zhang The user should provide the right hand side of the equation 13394b84eec9SHong Zhang using TSSetRHSFunction(). 13404b84eec9SHong Zhang 13414b84eec9SHong Zhang Notes: 13424b84eec9SHong Zhang The default is TSMPRKPM2, it can be changed with TSRKSetType() or -ts_mprk_type 13434b84eec9SHong Zhang 13444b84eec9SHong Zhang Level: beginner 13454b84eec9SHong Zhang 13464b84eec9SHong Zhang .seealso: TSCreate(), TS, TSSetType(), TSMPRKSetType(), TSMPRKGetType(), TSMPRKType, TSMPRKRegister(), TSMPRKSetMultirateType() 13474b84eec9SHong Zhang TSMPRKM2, TSMPRKM3, TSMPRKRFSMR3, TSMPRKRFSMR2 13484b84eec9SHong Zhang 13494b84eec9SHong Zhang M*/ 13504b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSCreate_MPRK(TS ts) 13514b84eec9SHong Zhang { 13524b84eec9SHong Zhang TS_MPRK *mprk; 13534b84eec9SHong Zhang PetscErrorCode ierr; 13544b84eec9SHong Zhang 13554b84eec9SHong Zhang PetscFunctionBegin; 13564b84eec9SHong Zhang ierr = TSMPRKInitializePackage();CHKERRQ(ierr); 13574b84eec9SHong Zhang 13584b84eec9SHong Zhang ts->ops->reset = TSReset_MPRK; 13594b84eec9SHong Zhang ts->ops->destroy = TSDestroy_MPRK; 13604b84eec9SHong Zhang ts->ops->view = TSView_MPRK; 13614b84eec9SHong Zhang ts->ops->load = TSLoad_MPRK; 13624b84eec9SHong Zhang ts->ops->setup = TSSetUp_MPRK; 13634b84eec9SHong Zhang ts->ops->step = TSStep_MPRK; 13644b84eec9SHong Zhang ts->ops->evaluatestep = TSEvaluateStep_MPRK; 13654b84eec9SHong Zhang ts->ops->setfromoptions = TSSetFromOptions_MPRK; 13664b84eec9SHong Zhang ts->ops->getstages = TSGetStages_MPRK; 13674b84eec9SHong Zhang 13684b84eec9SHong Zhang ierr = PetscNewLog(ts,&mprk);CHKERRQ(ierr); 13694b84eec9SHong Zhang ts->data = (void*)mprk; 13704b84eec9SHong Zhang 13714b84eec9SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSMPRKGetType_C",TSMPRKGetType_MPRK);CHKERRQ(ierr); 13724b84eec9SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSMPRKSetType_C",TSMPRKSetType_MPRK);CHKERRQ(ierr); 13734b84eec9SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSMPRKSetMultirateType_C",TSMPRKSetMultirateType);CHKERRQ(ierr); 13744b84eec9SHong Zhang 13754b84eec9SHong Zhang ierr = TSMPRKSetType(ts,TSMPRKDefault);CHKERRQ(ierr); 13764b84eec9SHong Zhang PetscFunctionReturn(0); 13774b84eec9SHong Zhang } 1378