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 135147403d9SBarry Smith 1364b84eec9SHong Zhang Options database: 137147403d9SBarry Smith . -ts_mprk_type 2a22 - select this scheme 1384b84eec9SHong Zhang 1394b84eec9SHong Zhang Level: advanced 1404b84eec9SHong Zhang 1414b84eec9SHong Zhang .seealso: TSMPRK, TSMPRKType, TSMPRKSetType() 1424b84eec9SHong Zhang M*/ 1434b84eec9SHong Zhang /*MC 144f944a40eSHong Zhang TSMPRK2A23 - Second Order Multirate Partitioned Runge-Kutta scheme based on RK2A. 14519c2959aSHong Zhang 14619c2959aSHong Zhang This method has eight stages for slow and medium and fast parts. The refinement factor of the stepsize is 2. 14719c2959aSHong Zhang r = 2, np = 3 148147403d9SBarry Smith 14919c2959aSHong Zhang Options database: 150147403d9SBarry Smith . -ts_mprk_type 2a23 - select this scheme 15119c2959aSHong Zhang 15219c2959aSHong Zhang Level: advanced 15319c2959aSHong Zhang 15419c2959aSHong Zhang .seealso: TSMPRK, TSMPRKType, TSMPRKSetType() 15519c2959aSHong Zhang M*/ 15619c2959aSHong Zhang /*MC 157f944a40eSHong Zhang TSMPRK2A32 - Second Order Multirate Partitioned Runge-Kutta scheme based on RK2A. 15819c2959aSHong Zhang 15919c2959aSHong Zhang This method has four stages for slow and fast parts. The refinement factor of the stepsize is 3. 16019c2959aSHong Zhang r = 3, np = 2 161147403d9SBarry Smith 16219c2959aSHong Zhang Options database: 163147403d9SBarry Smith . -ts_mprk_type 2a32 - select this scheme 16419c2959aSHong Zhang 16519c2959aSHong Zhang Level: advanced 16619c2959aSHong Zhang 16719c2959aSHong Zhang .seealso: TSMPRK, TSMPRKType, TSMPRKSetType() 16819c2959aSHong Zhang M*/ 16919c2959aSHong Zhang /*MC 170f944a40eSHong Zhang TSMPRK2A33 - Second Order Multirate Partitioned Runge-Kutta scheme based on RK2A. 17119c2959aSHong Zhang 17219c2959aSHong Zhang This method has eight stages for slow and medium and fast parts. The refinement factor of the stepsize is 3. 17319c2959aSHong Zhang r = 3, np = 3 174147403d9SBarry Smith 17519c2959aSHong Zhang Options database: 176147403d9SBarry Smith . -ts_mprk_type 2a33- select this scheme 17719c2959aSHong Zhang 17819c2959aSHong Zhang Level: advanced 17919c2959aSHong Zhang 18019c2959aSHong Zhang .seealso: TSMPRK, TSMPRKType, TSMPRKSetType() 18119c2959aSHong Zhang M*/ 18219c2959aSHong Zhang /*MC 183f944a40eSHong Zhang TSMPRK3P2M - Third Order Multirate Partitioned Runge-Kutta scheme. 1844b84eec9SHong Zhang 1854b84eec9SHong Zhang This method has eight stages for both slow and fast parts. 1864b84eec9SHong Zhang 1874b84eec9SHong Zhang Options database: 188147403d9SBarry Smith . -ts_mprk_type pm3 - select this scheme 1894b84eec9SHong Zhang 1904b84eec9SHong Zhang Level: advanced 1914b84eec9SHong Zhang 1924b84eec9SHong Zhang .seealso: TSMPRK, TSMPRKType, TSMPRKSetType() 1934b84eec9SHong Zhang M*/ 1944b84eec9SHong Zhang /*MC 195f944a40eSHong Zhang TSMPRKP2 - Second Order Multirate Partitioned Runge-Kutta scheme. 1964b84eec9SHong Zhang 1974b84eec9SHong Zhang This method has five stages for both slow and fast parts. 1984b84eec9SHong Zhang 1994b84eec9SHong Zhang Options database: 200147403d9SBarry Smith . -ts_mprk_type p2 - select this scheme 2014b84eec9SHong Zhang 2024b84eec9SHong Zhang Level: advanced 2034b84eec9SHong Zhang 2044b84eec9SHong Zhang .seealso: TSMPRK, TSMPRKType, TSMPRKSetType() 2054b84eec9SHong Zhang M*/ 2064b84eec9SHong Zhang /*MC 207f944a40eSHong Zhang TSMPRKP3 - Third Order Multirate Partitioned Runge-Kutta scheme. 2084b84eec9SHong Zhang 2094b84eec9SHong Zhang This method has ten stages for both slow and fast parts. 2104b84eec9SHong Zhang 2114b84eec9SHong Zhang Options database: 212147403d9SBarry Smith . -ts_mprk_type p3 - select this scheme 2134b84eec9SHong Zhang 2144b84eec9SHong Zhang Level: advanced 2154b84eec9SHong Zhang 2164b84eec9SHong Zhang .seealso: TSMPRK, TSMPRKType, TSMPRKSetType() 2174b84eec9SHong Zhang M*/ 2184b84eec9SHong Zhang 2194b84eec9SHong Zhang /*@C 2204b84eec9SHong Zhang TSMPRKRegisterAll - Registers all of the Partirioned Runge-Kutta explicit methods in TSMPRK 2214b84eec9SHong Zhang 2224b84eec9SHong Zhang Not Collective, but should be called by all processes which will need the schemes to be registered 2234b84eec9SHong Zhang 2244b84eec9SHong Zhang Level: advanced 2254b84eec9SHong Zhang 2264b84eec9SHong Zhang .seealso: TSMPRKRegisterDestroy() 2274b84eec9SHong Zhang @*/ 2284b84eec9SHong Zhang PetscErrorCode TSMPRKRegisterAll(void) 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}; 2449566063dSJacob Faibussowitsch PetscCall(TSMPRKGenerateTableau2(2,2,&Abase[0][0],bbase,&Asb[0][0],bsb,&Af[0][0],bf)); 2459566063dSJacob Faibussowitsch PetscCall(TSMPRKRegister(TSMPRK2A22,2,2,2,1,&Asb[0][0],bsb,NULL,rsb,NULL,NULL,NULL,NULL,&Af[0][0],bf,NULL)); 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}; 2569566063dSJacob Faibussowitsch PetscCall(TSMPRKGenerateTableau3(2,2,&Abase[0][0],bbase,&Asb[0][0],bsb,&Amb[0][0],bmb,&Af[0][0],bf)); 2579566063dSJacob Faibussowitsch PetscCall(TSMPRKRegister(TSMPRK2A23,2,2,2,2,&Asb[0][0],bsb,NULL,rsb,&Amb[0][0],bmb,NULL,rmb,&Af[0][0],bf,NULL)); 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}; 2689566063dSJacob Faibussowitsch PetscCall(TSMPRKGenerateTableau2(3,2,&Abase[0][0],bbase,&Asb[0][0],bsb,&Af[0][0],bf)); 2699566063dSJacob Faibussowitsch PetscCall(TSMPRKRegister(TSMPRK2A32,2,2,3,1,&Asb[0][0],bsb,NULL,rsb,NULL,NULL,NULL,NULL,&Af[0][0],bf,NULL)); 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}; 2809566063dSJacob Faibussowitsch PetscCall(TSMPRKGenerateTableau3(3,2,&Abase[0][0],bbase,&Asb[0][0],bsb,&Amb[0][0],bmb,&Af[0][0],bf)); 2819566063dSJacob Faibussowitsch PetscCall(TSMPRKRegister(TSMPRK2A33,2,2,3,3,&Asb[0][0],bsb,NULL,rsb,&Amb[0][0],bmb,NULL,rmb,&Af[0][0],bf,NULL)); 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)}; 3339566063dSJacob Faibussowitsch PetscCall(TSMPRKRegister(TSMPRKPM3,3,8,&As[0][0],bs,NULL,&A[0][0],b,NULL)); 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}; 3529566063dSJacob Faibussowitsch PetscCall(TSMPRKRegister(TSMPRKP2,2,5,1,1,&Asb[0][0],bsb,NULL,rsb,NULL,NULL,NULL,NULL,&Af[0][0],bf,NULL)); 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}; 3819566063dSJacob Faibussowitsch PetscCall(TSMPRKRegister(TSMPRKP3,3,5,2,1,&Asb[0][0],bsb,NULL,rsb,NULL,NULL,NULL,NULL,&Af[0][0],bf,NULL)); 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 .seealso: TSMPRKRegister(), TSMPRKRegisterAll() 3954b84eec9SHong Zhang @*/ 3964b84eec9SHong Zhang PetscErrorCode TSMPRKRegisterDestroy(void) 3974b84eec9SHong Zhang { 3984b84eec9SHong Zhang MPRKTableauLink link; 3994b84eec9SHong Zhang 4004b84eec9SHong Zhang PetscFunctionBegin; 4014b84eec9SHong Zhang while ((link = MPRKTableauList)) { 4024b84eec9SHong Zhang MPRKTableau t = &link->tab; 4034b84eec9SHong Zhang MPRKTableauList = link->next; 4049566063dSJacob Faibussowitsch PetscCall(PetscFree3(t->Asb,t->bsb,t->csb)); 4059566063dSJacob Faibussowitsch PetscCall(PetscFree3(t->Amb,t->bmb,t->cmb)); 4069566063dSJacob Faibussowitsch PetscCall(PetscFree3(t->Af,t->bf,t->cf)); 4079566063dSJacob Faibussowitsch PetscCall(PetscFree(t->rsb)); 4089566063dSJacob Faibussowitsch PetscCall(PetscFree(t->rmb)); 4099566063dSJacob Faibussowitsch PetscCall(PetscFree(t->name)); 4109566063dSJacob Faibussowitsch PetscCall(PetscFree(link)); 4114b84eec9SHong Zhang } 4124b84eec9SHong Zhang TSMPRKRegisterAllCalled = PETSC_FALSE; 4134b84eec9SHong Zhang PetscFunctionReturn(0); 4144b84eec9SHong Zhang } 4154b84eec9SHong Zhang 4164b84eec9SHong Zhang /*@C 4174b84eec9SHong Zhang TSMPRKInitializePackage - This function initializes everything in the TSMPRK package. It is called 4184b84eec9SHong Zhang from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_MPRK() 4194b84eec9SHong Zhang when using static libraries. 4204b84eec9SHong Zhang 4214b84eec9SHong Zhang Level: developer 4224b84eec9SHong Zhang 4234b84eec9SHong Zhang .seealso: PetscInitialize() 4244b84eec9SHong Zhang @*/ 4254b84eec9SHong Zhang PetscErrorCode TSMPRKInitializePackage(void) 4264b84eec9SHong Zhang { 4274b84eec9SHong Zhang PetscFunctionBegin; 4284b84eec9SHong Zhang if (TSMPRKPackageInitialized) PetscFunctionReturn(0); 4294b84eec9SHong Zhang TSMPRKPackageInitialized = PETSC_TRUE; 4309566063dSJacob Faibussowitsch PetscCall(TSMPRKRegisterAll()); 4319566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalize(TSMPRKFinalizePackage)); 4324b84eec9SHong Zhang PetscFunctionReturn(0); 4334b84eec9SHong Zhang } 4344b84eec9SHong Zhang 4354b84eec9SHong Zhang /*@C 436f944a40eSHong Zhang TSMPRKFinalizePackage - This function destroys everything in the TSMPRK package. It is 4374b84eec9SHong Zhang called from PetscFinalize(). 4384b84eec9SHong Zhang 4394b84eec9SHong Zhang Level: developer 4404b84eec9SHong Zhang 4414b84eec9SHong Zhang .seealso: PetscFinalize() 4424b84eec9SHong Zhang @*/ 4434b84eec9SHong Zhang PetscErrorCode TSMPRKFinalizePackage(void) 4444b84eec9SHong Zhang { 4454b84eec9SHong Zhang PetscFunctionBegin; 4464b84eec9SHong Zhang TSMPRKPackageInitialized = PETSC_FALSE; 4479566063dSJacob Faibussowitsch PetscCall(TSMPRKRegisterDestroy()); 4484b84eec9SHong Zhang PetscFunctionReturn(0); 4494b84eec9SHong Zhang } 4504b84eec9SHong Zhang 4514b84eec9SHong Zhang /*@C 4524b84eec9SHong Zhang TSMPRKRegister - register a MPRK scheme by providing the entries in the Butcher tableau 4534b84eec9SHong Zhang 4544b84eec9SHong Zhang Not Collective, but the same schemes should be registered on all processes on which they will be used 4554b84eec9SHong Zhang 4564b84eec9SHong Zhang Input Parameters: 4574b84eec9SHong Zhang + name - identifier for method 4584b84eec9SHong Zhang . order - approximation order of method 45979bc8a77SHong Zhang . s - number of stages in the base methods 46079bc8a77SHong Zhang . ratio1 - stepsize ratio at 1st level (e.g. slow/medium) 46179bc8a77SHong Zhang . ratio2 - stepsize ratio at 2nd level (e.g. medium/fast) 4624b84eec9SHong Zhang . Af - stage coefficients for fast components(dimension s*s, row-major) 4634b84eec9SHong Zhang . bf - step completion table for fast components(dimension s) 4644b84eec9SHong Zhang . cf - abscissa for fast components(dimension s) 4654b84eec9SHong Zhang . As - stage coefficients for slow components(dimension s*s, row-major) 4664b84eec9SHong Zhang . bs - step completion table for slow components(dimension s) 4674b84eec9SHong Zhang - cs - abscissa for slow components(dimension s) 4684b84eec9SHong Zhang 4694b84eec9SHong Zhang Notes: 4704b84eec9SHong Zhang Several MPRK methods are provided, this function is only needed to create new methods. 4714b84eec9SHong Zhang 4724b84eec9SHong Zhang Level: advanced 4734b84eec9SHong Zhang 4744b84eec9SHong Zhang .seealso: TSMPRK 4754b84eec9SHong Zhang @*/ 47679bc8a77SHong Zhang PetscErrorCode TSMPRKRegister(TSMPRKType name,PetscInt order, 47779bc8a77SHong Zhang PetscInt sbase,PetscInt ratio1,PetscInt ratio2, 47819c2959aSHong Zhang const PetscReal Asb[],const PetscReal bsb[],const PetscReal csb[],const PetscInt rsb[], 47919c2959aSHong Zhang const PetscReal Amb[],const PetscReal bmb[],const PetscReal cmb[],const PetscInt rmb[], 4804b84eec9SHong Zhang const PetscReal Af[],const PetscReal bf[],const PetscReal cf[]) 4814b84eec9SHong Zhang { 4824b84eec9SHong Zhang MPRKTableauLink link; 4834b84eec9SHong Zhang MPRKTableau t; 48479bc8a77SHong Zhang PetscInt s,i,j; 4854b84eec9SHong Zhang 4864b84eec9SHong Zhang PetscFunctionBegin; 4874b84eec9SHong Zhang PetscValidCharPointer(name,1); 488064a246eSJacob Faibussowitsch PetscValidRealPointer(Asb,6); 48979bc8a77SHong Zhang if (bsb) PetscValidRealPointer(bsb,7); 49079bc8a77SHong Zhang if (csb) PetscValidRealPointer(csb,8); 491064a246eSJacob Faibussowitsch if (rsb) PetscValidIntPointer(rsb,9); 49279bc8a77SHong Zhang if (Amb) PetscValidRealPointer(Amb,10); 49379bc8a77SHong Zhang if (bmb) PetscValidRealPointer(bmb,11); 49479bc8a77SHong Zhang if (cmb) PetscValidRealPointer(cmb,12); 495064a246eSJacob Faibussowitsch if (rmb) PetscValidIntPointer(rmb,13); 49679bc8a77SHong Zhang PetscValidRealPointer(Af,14); 49779bc8a77SHong Zhang if (bf) PetscValidRealPointer(bf,15); 49879bc8a77SHong Zhang if (cf) PetscValidRealPointer(cf,16); 4994b84eec9SHong Zhang 5009566063dSJacob Faibussowitsch PetscCall(PetscNew(&link)); 5014b84eec9SHong Zhang t = &link->tab; 5024b84eec9SHong Zhang 5039566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(name,&t->name)); 50479bc8a77SHong Zhang s = sbase*ratio1*ratio2; /* this is the dimension of the matrices below */ 5054b84eec9SHong Zhang t->order = order; 50679bc8a77SHong Zhang t->sbase = sbase; 5074b84eec9SHong Zhang t->s = s; 50819c2959aSHong Zhang t->np = 2; 50979bc8a77SHong Zhang 5109566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(s*s,&t->Af,s,&t->bf,s,&t->cf)); 5119566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->Af,Af,s*s)); 5124b84eec9SHong Zhang if (bf) { 5139566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->bf,bf,s)); 51419c2959aSHong Zhang } else 5154b84eec9SHong Zhang for (i=0; i<s; i++) t->bf[i] = Af[(s-1)*s+i]; 5164b84eec9SHong Zhang if (cf) { 5179566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->cf,cf,s)); 51819c2959aSHong Zhang } else { 5194b84eec9SHong Zhang for (i=0; i<s; i++) 5204b84eec9SHong Zhang for (j=0,t->cf[i]=0; j<s; j++) 5214b84eec9SHong Zhang t->cf[i] += Af[i*s+j]; 5224b84eec9SHong Zhang } 52319c2959aSHong Zhang 52419c2959aSHong Zhang if (Amb) { 52519c2959aSHong Zhang t->np = 3; 5269566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(s*s,&t->Amb,s,&t->bmb,s,&t->cmb)); 5279566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->Amb,Amb,s*s)); 52819c2959aSHong Zhang if (bmb) { 5299566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->bmb,bmb,s)); 53019c2959aSHong Zhang } else { 53119c2959aSHong Zhang for (i=0; i<s; i++) t->bmb[i] = Amb[(s-1)*s+i]; 5324b84eec9SHong Zhang } 53319c2959aSHong Zhang if (cmb) { 5349566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->cmb,cmb,s)); 53519c2959aSHong Zhang } else { 5364b84eec9SHong Zhang for (i=0; i<s; i++) 53719c2959aSHong Zhang for (j=0,t->cmb[i]=0; j<s; j++) 53819c2959aSHong Zhang t->cmb[i] += Amb[i*s+j]; 53919c2959aSHong Zhang } 54019c2959aSHong Zhang if (rmb) { 5419566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(s,&t->rmb)); 5429566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->rmb,rmb,s)); 543071fcb05SBarry Smith } else { 5449566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(s,&t->rmb)); 54519c2959aSHong Zhang } 54619c2959aSHong Zhang } 54719c2959aSHong Zhang 5489566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(s*s,&t->Asb,s,&t->bsb,s,&t->csb)); 5499566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->Asb,Asb,s*s)); 55019c2959aSHong Zhang if (bsb) { 5519566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->bsb,bsb,s)); 55219c2959aSHong Zhang } else 55319c2959aSHong Zhang for (i=0; i<s; i++) t->bsb[i] = Asb[(s-1)*s+i]; 55419c2959aSHong Zhang if (csb) { 5559566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->csb,csb,s)); 55619c2959aSHong Zhang } else { 55719c2959aSHong Zhang for (i=0; i<s; i++) 55819c2959aSHong Zhang for (j=0,t->csb[i]=0; j<s; j++) 55919c2959aSHong Zhang t->csb[i] += Asb[i*s+j]; 56019c2959aSHong Zhang } 56119c2959aSHong Zhang if (rsb) { 5629566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(s,&t->rsb)); 5639566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->rsb,rsb,s)); 564071fcb05SBarry Smith } else { 5659566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(s,&t->rsb)); 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 5784b84eec9SHong Zhang PetscFunctionBegin; 5799566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetSubTS(ts,"slow",&mprk->subts_slow)); 5809566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetSubTS(ts,"fast",&mprk->subts_fast)); 5813c633725SBarry Smith PetscCheck(mprk->subts_slow && mprk->subts_fast,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"); 5824b84eec9SHong Zhang 58319c2959aSHong Zhang /* Only copy the DM */ 5849566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&dm)); 58519c2959aSHong Zhang 5869566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetSubTS(ts,"slowbuffer",&mprk->subts_slowbuffer)); 58719c2959aSHong Zhang if (!mprk->subts_slowbuffer) { 58819c2959aSHong Zhang mprk->subts_slowbuffer = mprk->subts_slow; 58919c2959aSHong Zhang mprk->subts_slow = NULL; 59019c2959aSHong Zhang } 59119c2959aSHong Zhang if (mprk->subts_slow) { 5929566063dSJacob Faibussowitsch PetscCall(DMClone(dm,&newdm)); 5939566063dSJacob Faibussowitsch PetscCall(TSGetDM(mprk->subts_slow,&subdm)); 5949566063dSJacob Faibussowitsch PetscCall(DMCopyDMTS(subdm,newdm)); 5959566063dSJacob Faibussowitsch PetscCall(DMCopyDMSNES(subdm,newdm)); 5969566063dSJacob Faibussowitsch PetscCall(TSSetDM(mprk->subts_slow,newdm)); 5979566063dSJacob Faibussowitsch PetscCall(DMDestroy(&newdm)); 59819c2959aSHong Zhang } 5999566063dSJacob Faibussowitsch PetscCall(DMClone(dm,&newdm)); 6009566063dSJacob Faibussowitsch PetscCall(TSGetDM(mprk->subts_slowbuffer,&subdm)); 6019566063dSJacob Faibussowitsch PetscCall(DMCopyDMTS(subdm,newdm)); 6029566063dSJacob Faibussowitsch PetscCall(DMCopyDMSNES(subdm,newdm)); 6039566063dSJacob Faibussowitsch PetscCall(TSSetDM(mprk->subts_slowbuffer,newdm)); 6049566063dSJacob Faibussowitsch PetscCall(DMDestroy(&newdm)); 60519c2959aSHong Zhang 6069566063dSJacob Faibussowitsch PetscCall(DMClone(dm,&newdm)); 6079566063dSJacob Faibussowitsch PetscCall(TSGetDM(mprk->subts_fast,&subdm)); 6089566063dSJacob Faibussowitsch PetscCall(DMCopyDMTS(subdm,newdm)); 6099566063dSJacob Faibussowitsch PetscCall(DMCopyDMSNES(subdm,newdm)); 6109566063dSJacob Faibussowitsch PetscCall(TSSetDM(mprk->subts_fast,newdm)); 6119566063dSJacob Faibussowitsch PetscCall(DMDestroy(&newdm)); 61219c2959aSHong Zhang 61319c2959aSHong Zhang if (tab->np == 3) { 6149566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetSubTS(ts,"medium",&mprk->subts_medium)); 6159566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetSubTS(ts,"mediumbuffer",&mprk->subts_mediumbuffer)); 61619c2959aSHong Zhang if (mprk->subts_medium && !mprk->subts_mediumbuffer) { 61719c2959aSHong Zhang mprk->subts_mediumbuffer = mprk->subts_medium; 61819c2959aSHong Zhang mprk->subts_medium = NULL; 61919c2959aSHong Zhang } 62019c2959aSHong Zhang if (mprk->subts_medium) { 6219566063dSJacob Faibussowitsch PetscCall(DMClone(dm,&newdm)); 6229566063dSJacob Faibussowitsch PetscCall(TSGetDM(mprk->subts_medium,&subdm)); 6239566063dSJacob Faibussowitsch PetscCall(DMCopyDMTS(subdm,newdm)); 6249566063dSJacob Faibussowitsch PetscCall(DMCopyDMSNES(subdm,newdm)); 6259566063dSJacob Faibussowitsch PetscCall(TSSetDM(mprk->subts_medium,newdm)); 6269566063dSJacob Faibussowitsch PetscCall(DMDestroy(&newdm)); 62719c2959aSHong Zhang } 6289566063dSJacob Faibussowitsch PetscCall(DMClone(dm,&newdm)); 6299566063dSJacob Faibussowitsch PetscCall(TSGetDM(mprk->subts_mediumbuffer,&subdm)); 6309566063dSJacob Faibussowitsch PetscCall(DMCopyDMTS(subdm,newdm)); 6319566063dSJacob Faibussowitsch PetscCall(DMCopyDMSNES(subdm,newdm)); 6329566063dSJacob Faibussowitsch PetscCall(TSSetDM(mprk->subts_mediumbuffer,newdm)); 6339566063dSJacob Faibussowitsch PetscCall(DMDestroy(&newdm)); 63419c2959aSHong Zhang } 6354b84eec9SHong Zhang PetscFunctionReturn(0); 6364b84eec9SHong Zhang } 6374b84eec9SHong Zhang 6384b84eec9SHong Zhang /* 6394b84eec9SHong Zhang This if for nonsplit RHS MPRK 6404b84eec9SHong Zhang The step completion formula is 6414b84eec9SHong Zhang 6424b84eec9SHong Zhang x1 = x0 + h b^T YdotRHS 6434b84eec9SHong Zhang 6444b84eec9SHong Zhang */ 6454b84eec9SHong Zhang static PetscErrorCode TSEvaluateStep_MPRK(TS ts,PetscInt order,Vec X,PetscBool *done) 6464b84eec9SHong Zhang { 6474b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 6484b84eec9SHong Zhang MPRKTableau tab = mprk->tableau; 6497c0df07dSHong Zhang PetscScalar *wf = mprk->work_fast; 6504b84eec9SHong Zhang PetscReal h = ts->time_step; 6514b84eec9SHong Zhang PetscInt s = tab->s,j; 6524b84eec9SHong Zhang 6534b84eec9SHong Zhang PetscFunctionBegin; 6544b84eec9SHong Zhang for (j=0; j<s; j++) wf[j] = h*tab->bf[j]; 6559566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol,X)); 6569566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X,s,wf,mprk->YdotRHS)); 6574b84eec9SHong Zhang PetscFunctionReturn(0); 6584b84eec9SHong Zhang } 6594b84eec9SHong Zhang 6604b84eec9SHong Zhang static PetscErrorCode TSStep_MPRK(TS ts) 6614b84eec9SHong Zhang { 6624b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 6637c0df07dSHong Zhang Vec *Y = mprk->Y,*YdotRHS = mprk->YdotRHS,*YdotRHS_fast = mprk->YdotRHS_fast,*YdotRHS_slow = mprk->YdotRHS_slow,*YdotRHS_slowbuffer = mprk->YdotRHS_slowbuffer; 6647c0df07dSHong Zhang Vec Yslow,Yslowbuffer,Yfast; 6654b84eec9SHong Zhang MPRKTableau tab = mprk->tableau; 6664b84eec9SHong Zhang const PetscInt s = tab->s; 66719c2959aSHong Zhang const PetscReal *Af = tab->Af,*cf = tab->cf,*Asb = tab->Asb,*csb = tab->csb; 668ebd5ed4eSHong Zhang PetscScalar *wf = mprk->work_fast,*wsb = mprk->work_slowbuffer; 669ebd5ed4eSHong Zhang PetscInt i,j; 6704b84eec9SHong Zhang PetscReal next_time_step = ts->time_step,t = ts->ptime,h = ts->time_step; 6714b84eec9SHong Zhang 6724b84eec9SHong Zhang PetscFunctionBegin; 6734b84eec9SHong Zhang for (i=0; i<s; i++) { 6744b84eec9SHong Zhang mprk->stage_time = t + h*cf[i]; 6759566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts,mprk->stage_time)); 6769566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol,Y[i])); 6774b84eec9SHong Zhang 6787c0df07dSHong Zhang /* slow buffer region */ 67919c2959aSHong Zhang for (j=0; j<i; j++) wsb[j] = h*Asb[i*s+j]; 6809849be05SHong Zhang for (j=0; j<i; j++) { 6819566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(YdotRHS[j],mprk->is_slowbuffer,&YdotRHS_slowbuffer[j])); 6827c0df07dSHong Zhang } 6839566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i],mprk->is_slowbuffer,&Yslowbuffer)); 6849566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Yslowbuffer,i,wsb,mprk->YdotRHS_slowbuffer)); 6859566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i],mprk->is_slowbuffer,&Yslowbuffer)); 6869849be05SHong Zhang for (j=0; j<i; j++) { 6879566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(YdotRHS[j],mprk->is_slowbuffer,&YdotRHS_slowbuffer[j])); 6887c0df07dSHong Zhang } 6897c0df07dSHong Zhang /* slow region */ 6907c0df07dSHong Zhang if (mprk->is_slow) { 6919849be05SHong Zhang for (j=0; j<i; j++) { 6929566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(YdotRHS[j],mprk->is_slow,&YdotRHS_slow[j])); 6937c0df07dSHong Zhang } 6949566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i],mprk->is_slow,&Yslow)); 6959566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Yslow,i,wsb,mprk->YdotRHS_slow)); 6969566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i],mprk->is_slow,&Yslow)); 697ebd5ed4eSHong Zhang for (j=0; j<i; j++) { 6989566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(YdotRHS[j],mprk->is_slow,&YdotRHS_slow[j])); 6997c0df07dSHong Zhang } 7007c0df07dSHong Zhang } 7014b84eec9SHong Zhang 7027c0df07dSHong Zhang /* fast region */ 7034b84eec9SHong Zhang for (j=0; j<i; j++) wf[j] = h*Af[i*s+j]; 7049849be05SHong Zhang for (j=0; j<i; j++) { 7059566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(YdotRHS[j],mprk->is_fast,&YdotRHS_fast[j])); 7067c0df07dSHong Zhang } 7079566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i],mprk->is_fast,&Yfast)); 7089566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Yfast,i,wf,mprk->YdotRHS_fast)); 7099566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i],mprk->is_fast,&Yfast)); 7109849be05SHong Zhang for (j=0; j<i; j++) { 7119566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(YdotRHS[j],mprk->is_fast,&YdotRHS_fast[j])); 7127c0df07dSHong Zhang } 7137c0df07dSHong Zhang if (tab->np == 3) { 7147c0df07dSHong Zhang Vec *YdotRHS_medium = mprk->YdotRHS_medium,*YdotRHS_mediumbuffer = mprk->YdotRHS_mediumbuffer; 7157c0df07dSHong Zhang Vec Ymedium,Ymediumbuffer; 716ebd5ed4eSHong Zhang PetscScalar *wmb = mprk->work_mediumbuffer; 717ebd5ed4eSHong Zhang 718ebd5ed4eSHong Zhang for (j=0; j<i; j++) wmb[j] = h*tab->Amb[i*s+j]; 7197c0df07dSHong Zhang /* medium region */ 7207c0df07dSHong Zhang if (mprk->is_medium) { 7217c0df07dSHong Zhang for (j=0; j<i; j++) { 7229566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(YdotRHS[j],mprk->is_medium,&YdotRHS_medium[j])); 7237c0df07dSHong Zhang } 7249566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i],mprk->is_medium,&Ymedium)); 7259566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Ymedium,i,wmb,mprk->YdotRHS_medium)); 7269566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i],mprk->is_medium,&Ymedium)); 727ebd5ed4eSHong Zhang for (j=0; j<i; j++) { 7289566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(YdotRHS[j],mprk->is_medium,&YdotRHS_medium[j])); 7297c0df07dSHong Zhang } 7307c0df07dSHong Zhang } 7317c0df07dSHong Zhang /* medium buffer region */ 7329849be05SHong Zhang for (j=0; j<i; j++) { 7339566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(YdotRHS[j],mprk->is_mediumbuffer,&YdotRHS_mediumbuffer[j])); 7347c0df07dSHong Zhang } 7359566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i],mprk->is_mediumbuffer,&Ymediumbuffer)); 7369566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Ymediumbuffer,i,wmb,mprk->YdotRHS_mediumbuffer)); 7379566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i],mprk->is_mediumbuffer,&Ymediumbuffer)); 7389849be05SHong Zhang for (j=0; j<i; j++) { 7399566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(YdotRHS[j],mprk->is_mediumbuffer,&YdotRHS_mediumbuffer[j])); 7407c0df07dSHong Zhang } 7417c0df07dSHong Zhang } 7429566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts,mprk->stage_time,i,Y)); 7434b84eec9SHong Zhang /* compute the stage RHS by fast and slow tableau respectively */ 7449566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(ts,t+h*csb[i],Y[i],YdotRHS[i])); 7454b84eec9SHong Zhang } 7469566063dSJacob Faibussowitsch PetscCall(TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL)); 7474b84eec9SHong Zhang ts->ptime += ts->time_step; 7484b84eec9SHong Zhang ts->time_step = next_time_step; 7494b84eec9SHong Zhang PetscFunctionReturn(0); 7504b84eec9SHong Zhang } 7514b84eec9SHong Zhang 7524b84eec9SHong Zhang /* 753f944a40eSHong Zhang This if for the case when split RHS is used 7544b84eec9SHong Zhang The step completion formula is 7554b84eec9SHong Zhang x1 = x0 + h b^T YdotRHS 7564b84eec9SHong Zhang */ 7574b84eec9SHong Zhang static PetscErrorCode TSEvaluateStep_MPRKSPLIT(TS ts,PetscInt order,Vec X,PetscBool *done) 7584b84eec9SHong Zhang { 7594b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 7604b84eec9SHong Zhang MPRKTableau tab = mprk->tableau; 761*6aad120cSJose E. Roman Vec Xslow,Xfast,Xslowbuffer; /* subvectors for slow and fast components in X respectively */ 76219c2959aSHong Zhang PetscScalar *wf = mprk->work_fast,*ws = mprk->work_slow,*wsb = mprk->work_slowbuffer; 7634b84eec9SHong Zhang PetscReal h = ts->time_step; 76479bc8a77SHong Zhang PetscInt s = tab->s,j,computedstages; 7654b84eec9SHong Zhang 7664b84eec9SHong Zhang PetscFunctionBegin; 7679566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol,X)); 76819c2959aSHong Zhang 76919c2959aSHong Zhang /* slow region */ 77019c2959aSHong Zhang if (mprk->is_slow) { 77179bc8a77SHong Zhang computedstages = 0; 77219c2959aSHong Zhang for (j=0; j<s; j++) { 7739849be05SHong Zhang if (tab->rsb[j]) ws[tab->rsb[j]-1] += h*tab->bsb[j]; 77479bc8a77SHong Zhang else ws[computedstages++] = h*tab->bsb[j]; 77519c2959aSHong Zhang } 7769566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(X,mprk->is_slow,&Xslow)); 7779566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Xslow,computedstages,ws,mprk->YdotRHS_slow)); 7789566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(X,mprk->is_slow,&Xslow)); 77919c2959aSHong Zhang } 78019c2959aSHong Zhang 7819d6e09e9SHong Zhang if (tab->np == 3 && mprk->is_medium) { 7829d6e09e9SHong Zhang computedstages = 0; 7839d6e09e9SHong Zhang for (j=0; j<s; j++) { 7849d6e09e9SHong Zhang if (tab->rmb[j]) wsb[computedstages-tab->sbase+(tab->rmb[j]-1)%tab->sbase] += h*tab->bsb[j]; 7859d6e09e9SHong Zhang else wsb[computedstages++] = h*tab->bsb[j]; 7869d6e09e9SHong Zhang } 7879566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(X,mprk->is_slowbuffer,&Xslowbuffer)); 7889566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Xslowbuffer,computedstages,wsb,mprk->YdotRHS_slowbuffer)); 7899566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(X,mprk->is_slowbuffer,&Xslowbuffer)); 7909d6e09e9SHong Zhang } else { 79119c2959aSHong Zhang /* slow buffer region */ 79219c2959aSHong Zhang for (j=0; j<s; j++) wsb[j] = h*tab->bsb[j]; 7939566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(X,mprk->is_slowbuffer,&Xslowbuffer)); 7949566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Xslowbuffer,s,wsb,mprk->YdotRHS_slowbuffer)); 7959566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(X,mprk->is_slowbuffer,&Xslowbuffer)); 7969d6e09e9SHong Zhang } 79719c2959aSHong Zhang if (tab->np == 3) { 79819c2959aSHong Zhang Vec Xmedium,Xmediumbuffer; 79919c2959aSHong Zhang PetscScalar *wm = mprk->work_medium,*wmb = mprk->work_mediumbuffer; 8009d6e09e9SHong Zhang /* medium region and slow buffer region */ 80119c2959aSHong Zhang if (mprk->is_medium) { 80279bc8a77SHong Zhang computedstages = 0; 80319c2959aSHong Zhang for (j=0; j<s; j++) { 80479bc8a77SHong Zhang if (tab->rmb[j]) wm[computedstages-tab->sbase+(tab->rmb[j]-1)%tab->sbase] += h*tab->bmb[j]; 80579bc8a77SHong Zhang else wm[computedstages++] = h*tab->bmb[j]; 80619c2959aSHong Zhang } 8079566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(X,mprk->is_medium,&Xmedium)); 8089566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Xmedium,computedstages,wm,mprk->YdotRHS_medium)); 8099566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(X,mprk->is_medium,&Xmedium)); 81019c2959aSHong Zhang } 81119c2959aSHong Zhang /* medium buffer region */ 81219c2959aSHong Zhang for (j=0; j<s; j++) wmb[j] = h*tab->bmb[j]; 8139566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(X,mprk->is_mediumbuffer,&Xmediumbuffer)); 8149566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Xmediumbuffer,s,wmb,mprk->YdotRHS_mediumbuffer)); 8159566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(X,mprk->is_mediumbuffer,&Xmediumbuffer)); 81619c2959aSHong Zhang } 81719c2959aSHong Zhang /* fast region */ 81819c2959aSHong Zhang for (j=0; j<s; j++) wf[j] = h*tab->bf[j]; 8199566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(X,mprk->is_fast,&Xfast)); 8209566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Xfast,s,wf,mprk->YdotRHS_fast)); 8219566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(X,mprk->is_fast,&Xfast)); 8224b84eec9SHong Zhang PetscFunctionReturn(0); 8234b84eec9SHong Zhang } 8244b84eec9SHong Zhang 8254b84eec9SHong Zhang static PetscErrorCode TSStep_MPRKSPLIT(TS ts) 8264b84eec9SHong Zhang { 8274b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 8284b84eec9SHong Zhang MPRKTableau tab = mprk->tableau; 82919c2959aSHong Zhang Vec *Y = mprk->Y,*YdotRHS_fast = mprk->YdotRHS_fast,*YdotRHS_slow = mprk->YdotRHS_slow,*YdotRHS_slowbuffer = mprk->YdotRHS_slowbuffer; 83019c2959aSHong Zhang Vec Yslow,Yslowbuffer,Yfast; /* subvectors for slow and fast components in Y[i] respectively */ 83119c2959aSHong Zhang PetscInt s = tab->s; 83219c2959aSHong Zhang const PetscReal *Af = tab->Af,*cf = tab->cf,*Asb = tab->Asb,*csb = tab->csb; 83319c2959aSHong Zhang PetscScalar *wf = mprk->work_fast,*ws = mprk->work_slow,*wsb = mprk->work_slowbuffer; 83479bc8a77SHong Zhang PetscInt i,j,computedstages; 8354b84eec9SHong Zhang PetscReal next_time_step = ts->time_step,t = ts->ptime,h = ts->time_step; 8364b84eec9SHong Zhang 8374b84eec9SHong Zhang PetscFunctionBegin; 8384b84eec9SHong Zhang for (i=0; i<s; i++) { 8394b84eec9SHong Zhang mprk->stage_time = t + h*cf[i]; 8409566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts,mprk->stage_time)); 8414b84eec9SHong Zhang /* calculate the stage value for fast and slow components respectively */ 8429566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol,Y[i])); 84319c2959aSHong Zhang for (j=0; j<i; j++) wsb[j] = h*Asb[i*s+j]; 84419c2959aSHong Zhang 84519c2959aSHong Zhang /* slow buffer region */ 8469d6e09e9SHong Zhang if (tab->np == 3 && mprk->is_medium) { 8479d6e09e9SHong Zhang if (tab->rmb[i]) { 8489566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i],mprk->is_slowbuffer,&Yslowbuffer)); 8499566063dSJacob Faibussowitsch PetscCall(VecISCopy(Y[tab->rmb[i]-1],mprk->is_slowbuffer,SCATTER_REVERSE,Yslowbuffer)); 8509566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i],mprk->is_slowbuffer,&Yslowbuffer)); 8519d6e09e9SHong Zhang } else { 8529d6e09e9SHong Zhang PetscScalar *wm = mprk->work_medium; 8539d6e09e9SHong Zhang computedstages = 0; 8549d6e09e9SHong Zhang for (j=0; j<i; j++) { 8559d6e09e9SHong Zhang if (tab->rmb[j]) wm[computedstages-tab->sbase+(tab->rmb[j]-1)%tab->sbase] += wsb[j]; 8569d6e09e9SHong Zhang else wm[computedstages++] = wsb[j]; 8579d6e09e9SHong Zhang } 8589566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i],mprk->is_slowbuffer,&Yslowbuffer)); 8599566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Yslowbuffer,computedstages,wm,YdotRHS_slowbuffer)); 8609566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i],mprk->is_slowbuffer,&Yslowbuffer)); 8619d6e09e9SHong Zhang } 8629d6e09e9SHong Zhang } else { 8639566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i],mprk->is_slowbuffer,&Yslowbuffer)); 8649566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Yslowbuffer,i,wsb,YdotRHS_slowbuffer)); 8659566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i],mprk->is_slowbuffer,&Yslowbuffer)); 8669d6e09e9SHong Zhang } 8679849be05SHong Zhang 86819c2959aSHong Zhang /* slow region */ 8699849be05SHong Zhang if (mprk->is_slow) { 8709849be05SHong Zhang if (tab->rsb[i]) { /* repeat previous stage */ 8719566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i],mprk->is_slow,&Yslow)); 8729566063dSJacob Faibussowitsch PetscCall(VecISCopy(Y[tab->rsb[i]-1],mprk->is_slow,SCATTER_REVERSE,Yslow)); 8739566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i],mprk->is_slow,&Yslow)); 8749849be05SHong Zhang } else { 87579bc8a77SHong Zhang computedstages = 0; 87619c2959aSHong Zhang for (j=0; j<i; j++) { 87779bc8a77SHong Zhang if (tab->rsb[j]) ws[tab->rsb[j]-1] += wsb[j]; 87879bc8a77SHong Zhang else ws[computedstages++] = wsb[j]; 87919c2959aSHong Zhang } 8809566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i],mprk->is_slow,&Yslow)); 8819566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Yslow,computedstages,ws,YdotRHS_slow)); 8829566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i],mprk->is_slow,&Yslow)); 88319c2959aSHong Zhang /* only depends on the slow buffer region */ 8849566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(mprk->subts_slow,t+h*csb[i],Y[i],YdotRHS_slow[computedstages])); 88519c2959aSHong Zhang } 8869849be05SHong Zhang } 88719c2959aSHong Zhang 88819c2959aSHong Zhang /* fast region */ 88919c2959aSHong Zhang for (j=0; j<i; j++) wf[j] = h*Af[i*s+j]; 8909566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i],mprk->is_fast,&Yfast)); 8919566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Yfast,i,wf,YdotRHS_fast)); 8929566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i],mprk->is_fast,&Yfast)); 89319c2959aSHong Zhang 89419c2959aSHong Zhang if (tab->np == 3) { 89519c2959aSHong Zhang Vec *YdotRHS_medium = mprk->YdotRHS_medium,*YdotRHS_mediumbuffer = mprk->YdotRHS_mediumbuffer; 89619c2959aSHong Zhang Vec Ymedium,Ymediumbuffer; 89719c2959aSHong Zhang const PetscReal *Amb = tab->Amb,*cmb = tab->cmb; 89819c2959aSHong Zhang PetscScalar *wm = mprk->work_medium,*wmb = mprk->work_mediumbuffer; 89919c2959aSHong Zhang 90019c2959aSHong Zhang for (j=0; j<i; j++) wmb[j] = h*Amb[i*s+j]; 90119c2959aSHong Zhang /* medium buffer region */ 9029566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i],mprk->is_mediumbuffer,&Ymediumbuffer)); 9039566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Ymediumbuffer,i,wmb,YdotRHS_mediumbuffer)); 9049566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i],mprk->is_mediumbuffer,&Ymediumbuffer)); 90519c2959aSHong Zhang /* medium region */ 90679bc8a77SHong Zhang if (mprk->is_medium) { 90779bc8a77SHong Zhang if (tab->rmb[i]) { /* repeat previous stage */ 9089566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i],mprk->is_medium,&Ymedium)); 9099566063dSJacob Faibussowitsch PetscCall(VecISCopy(Y[tab->rmb[i]-1],mprk->is_medium,SCATTER_REVERSE,Ymedium)); 9109566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i],mprk->is_medium,&Ymedium)); 91179bc8a77SHong Zhang } else { 91279bc8a77SHong Zhang computedstages = 0; 91379bc8a77SHong Zhang for (j=0; j<i; j++) { 91479bc8a77SHong Zhang if (tab->rmb[j]) wm[computedstages-tab->sbase+(tab->rmb[j]-1)%tab->sbase] += wmb[j]; 91579bc8a77SHong Zhang else wm[computedstages++] = wmb[j]; 91679bc8a77SHong Zhang 91779bc8a77SHong Zhang } 9189566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i],mprk->is_medium,&Ymedium)); 9199566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Ymedium,computedstages,wm,YdotRHS_medium)); 9209566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i],mprk->is_medium,&Ymedium)); 92179bc8a77SHong Zhang /* only depends on the medium buffer region */ 9229566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(mprk->subts_medium,t+h*cmb[i],Y[i],YdotRHS_medium[computedstages])); 9239d6e09e9SHong Zhang /* must be computed after all regions are updated in Y */ 9249566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(mprk->subts_slowbuffer,t+h*csb[i],Y[i],YdotRHS_slowbuffer[computedstages])); 92579bc8a77SHong Zhang } 92619c2959aSHong Zhang } 92719c2959aSHong Zhang /* must be computed after fast region and slow region are updated in Y */ 9289566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(mprk->subts_mediumbuffer,t+h*cmb[i],Y[i],YdotRHS_mediumbuffer[i])); 92919c2959aSHong Zhang } 9309d6e09e9SHong Zhang if (!(tab->np == 3 && mprk->is_medium)) { 9319566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(mprk->subts_slowbuffer,t+h*csb[i],Y[i],YdotRHS_slowbuffer[i])); 9329d6e09e9SHong Zhang } 9339566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(mprk->subts_fast,t+h*cf[i],Y[i],YdotRHS_fast[i])); 9344b84eec9SHong Zhang } 93579bc8a77SHong Zhang 9369566063dSJacob Faibussowitsch PetscCall(TSEvaluateStep(ts,tab->order,ts->vec_sol,NULL)); 9374b84eec9SHong Zhang ts->ptime += ts->time_step; 9384b84eec9SHong Zhang ts->time_step = next_time_step; 9394b84eec9SHong Zhang PetscFunctionReturn(0); 9404b84eec9SHong Zhang } 9414b84eec9SHong Zhang 9424b84eec9SHong Zhang static PetscErrorCode TSMPRKTableauReset(TS ts) 9434b84eec9SHong Zhang { 9444b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 9454b84eec9SHong Zhang MPRKTableau tab = mprk->tableau; 9464b84eec9SHong Zhang 9474b84eec9SHong Zhang PetscFunctionBegin; 9484b84eec9SHong Zhang if (!tab) PetscFunctionReturn(0); 9499566063dSJacob Faibussowitsch PetscCall(PetscFree(mprk->work_fast)); 9509566063dSJacob Faibussowitsch PetscCall(PetscFree(mprk->work_slow)); 9519566063dSJacob Faibussowitsch PetscCall(PetscFree(mprk->work_slowbuffer)); 9529566063dSJacob Faibussowitsch PetscCall(PetscFree(mprk->work_medium)); 9539566063dSJacob Faibussowitsch PetscCall(PetscFree(mprk->work_mediumbuffer)); 9549566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s,&mprk->Y)); 9550fe4d17eSHong Zhang if (ts->use_splitrhsfunction) { 9569566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s,&mprk->YdotRHS_fast)); 9579566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s,&mprk->YdotRHS_slow)); 9589566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s,&mprk->YdotRHS_slowbuffer)); 9599566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s,&mprk->YdotRHS_medium)); 9609566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s,&mprk->YdotRHS_mediumbuffer)); 9610fe4d17eSHong Zhang } else { 9629566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s,&mprk->YdotRHS)); 9637c0df07dSHong Zhang if (mprk->is_slow) { 9649566063dSJacob Faibussowitsch PetscCall(PetscFree(mprk->YdotRHS_slow)); 9654b84eec9SHong Zhang } 9669566063dSJacob Faibussowitsch PetscCall(PetscFree(mprk->YdotRHS_slowbuffer)); 9677c0df07dSHong Zhang if (tab->np == 3) { 9687c0df07dSHong Zhang if (mprk->is_medium) { 9699566063dSJacob Faibussowitsch PetscCall(PetscFree(mprk->YdotRHS_medium)); 9707c0df07dSHong Zhang } 9719566063dSJacob Faibussowitsch PetscCall(PetscFree(mprk->YdotRHS_mediumbuffer)); 9727c0df07dSHong Zhang } 9739566063dSJacob Faibussowitsch PetscCall(PetscFree(mprk->YdotRHS_fast)); 9747c0df07dSHong Zhang } 9754b84eec9SHong Zhang PetscFunctionReturn(0); 9764b84eec9SHong Zhang } 9774b84eec9SHong Zhang 9784b84eec9SHong Zhang static PetscErrorCode TSReset_MPRK(TS ts) 9794b84eec9SHong Zhang { 9804b84eec9SHong Zhang PetscFunctionBegin; 9819566063dSJacob Faibussowitsch PetscCall(TSMPRKTableauReset(ts)); 9824b84eec9SHong Zhang PetscFunctionReturn(0); 9834b84eec9SHong Zhang } 9844b84eec9SHong Zhang 9854b84eec9SHong Zhang static PetscErrorCode DMCoarsenHook_TSMPRK(DM fine,DM coarse,void *ctx) 9864b84eec9SHong Zhang { 9874b84eec9SHong Zhang PetscFunctionBegin; 9884b84eec9SHong Zhang PetscFunctionReturn(0); 9894b84eec9SHong Zhang } 9904b84eec9SHong Zhang 9914b84eec9SHong Zhang static PetscErrorCode DMRestrictHook_TSMPRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 9924b84eec9SHong Zhang { 9934b84eec9SHong Zhang PetscFunctionBegin; 9944b84eec9SHong Zhang PetscFunctionReturn(0); 9954b84eec9SHong Zhang } 9964b84eec9SHong Zhang 9974b84eec9SHong Zhang static PetscErrorCode DMSubDomainHook_TSMPRK(DM dm,DM subdm,void *ctx) 9984b84eec9SHong Zhang { 9994b84eec9SHong Zhang PetscFunctionBegin; 10004b84eec9SHong Zhang PetscFunctionReturn(0); 10014b84eec9SHong Zhang } 10024b84eec9SHong Zhang 10034b84eec9SHong Zhang static PetscErrorCode DMSubDomainRestrictHook_TSMPRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 10044b84eec9SHong Zhang { 10054b84eec9SHong Zhang PetscFunctionBegin; 10064b84eec9SHong Zhang PetscFunctionReturn(0); 10074b84eec9SHong Zhang } 10084b84eec9SHong Zhang 10094b84eec9SHong Zhang static PetscErrorCode TSMPRKTableauSetUp(TS ts) 10104b84eec9SHong Zhang { 10114b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 10124b84eec9SHong Zhang MPRKTableau tab = mprk->tableau; 101319c2959aSHong Zhang Vec YdotRHS_slow,YdotRHS_slowbuffer,YdotRHS_medium,YdotRHS_mediumbuffer,YdotRHS_fast; 10144b84eec9SHong Zhang 10154b84eec9SHong Zhang PetscFunctionBegin; 10169566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol,tab->s,&mprk->Y)); 10177c0df07dSHong Zhang if (mprk->is_slow) { 10189566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tab->s,&mprk->work_slow)); 10194b84eec9SHong Zhang } 10209566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tab->s,&mprk->work_slowbuffer)); 10217c0df07dSHong Zhang if (tab->np == 3) { 10227c0df07dSHong Zhang if (mprk->is_medium) { 10239566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tab->s,&mprk->work_medium)); 10247c0df07dSHong Zhang } 10259566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tab->s,&mprk->work_mediumbuffer)); 10267c0df07dSHong Zhang } 10279566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tab->s,&mprk->work_fast)); 10287c0df07dSHong Zhang 10290fe4d17eSHong Zhang if (ts->use_splitrhsfunction) { 103019c2959aSHong Zhang if (mprk->is_slow) { 10319566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(ts->vec_sol,mprk->is_slow,&YdotRHS_slow)); 10329566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(YdotRHS_slow,tab->s,&mprk->YdotRHS_slow)); 10339566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(ts->vec_sol,mprk->is_slow,&YdotRHS_slow)); 103419c2959aSHong Zhang } 10359566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(ts->vec_sol,mprk->is_slowbuffer,&YdotRHS_slowbuffer)); 10369566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(YdotRHS_slowbuffer,tab->s,&mprk->YdotRHS_slowbuffer)); 10379566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(ts->vec_sol,mprk->is_slowbuffer,&YdotRHS_slowbuffer)); 103819c2959aSHong Zhang if (tab->np == 3) { 103919c2959aSHong Zhang if (mprk->is_medium) { 10409566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(ts->vec_sol,mprk->is_medium,&YdotRHS_medium)); 10419566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(YdotRHS_medium,tab->s,&mprk->YdotRHS_medium)); 10429566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(ts->vec_sol,mprk->is_medium,&YdotRHS_medium)); 104319c2959aSHong Zhang } 10449566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(ts->vec_sol,mprk->is_mediumbuffer,&YdotRHS_mediumbuffer)); 10459566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(YdotRHS_mediumbuffer,tab->s,&mprk->YdotRHS_mediumbuffer)); 10469566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(ts->vec_sol,mprk->is_mediumbuffer,&YdotRHS_mediumbuffer)); 104719c2959aSHong Zhang } 10489566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(ts->vec_sol,mprk->is_fast,&YdotRHS_fast)); 10499566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(YdotRHS_fast,tab->s,&mprk->YdotRHS_fast)); 10509566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(ts->vec_sol,mprk->is_fast,&YdotRHS_fast)); 10510fe4d17eSHong Zhang } else { 10529566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol,tab->s,&mprk->YdotRHS)); 10530fe4d17eSHong Zhang if (mprk->is_slow) { 10549566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tab->s,&mprk->YdotRHS_slow)); 10550fe4d17eSHong Zhang } 10569566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tab->s,&mprk->YdotRHS_slowbuffer)); 10570fe4d17eSHong Zhang if (tab->np == 3) { 10580fe4d17eSHong Zhang if (mprk->is_medium) { 10599566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tab->s,&mprk->YdotRHS_medium)); 10600fe4d17eSHong Zhang } 10619566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tab->s,&mprk->YdotRHS_mediumbuffer)); 10620fe4d17eSHong Zhang } 10639566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tab->s,&mprk->YdotRHS_fast)); 10644b84eec9SHong Zhang } 10654b84eec9SHong Zhang PetscFunctionReturn(0); 10664b84eec9SHong Zhang } 10674b84eec9SHong Zhang 10684b84eec9SHong Zhang static PetscErrorCode TSSetUp_MPRK(TS ts) 10694b84eec9SHong Zhang { 10704b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 107119c2959aSHong Zhang MPRKTableau tab = mprk->tableau; 10724b84eec9SHong Zhang DM dm; 10734b84eec9SHong Zhang 10744b84eec9SHong Zhang PetscFunctionBegin; 10759566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetIS(ts,"slow",&mprk->is_slow)); 10769566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetIS(ts,"fast",&mprk->is_fast)); 10773c633725SBarry Smith PetscCheck(mprk->is_slow && mprk->is_fast,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); 107819c2959aSHong Zhang 107919c2959aSHong Zhang if (tab->np == 3) { 10809566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetIS(ts,"medium",&mprk->is_medium)); 10813c633725SBarry Smith PetscCheck(mprk->is_medium,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); 10829566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetIS(ts,"mediumbuffer",&mprk->is_mediumbuffer)); 108319c2959aSHong Zhang if (!mprk->is_mediumbuffer) { /* let medium buffer cover whole medium region */ 108419c2959aSHong Zhang mprk->is_mediumbuffer = mprk->is_medium; 108519c2959aSHong Zhang mprk->is_medium = NULL; 108619c2959aSHong Zhang } 108719c2959aSHong Zhang } 108819c2959aSHong Zhang 108919c2959aSHong Zhang /* If users do not provide buffer region settings, the solver will do them automatically, but with a performance penalty */ 10909566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetIS(ts,"slowbuffer",&mprk->is_slowbuffer)); 109119c2959aSHong Zhang if (!mprk->is_slowbuffer) { /* let slow buffer cover whole slow region */ 109219c2959aSHong Zhang mprk->is_slowbuffer = mprk->is_slow; 109319c2959aSHong Zhang mprk->is_slow = NULL; 109419c2959aSHong Zhang } 10959566063dSJacob Faibussowitsch PetscCall(TSCheckImplicitTerm(ts)); 10969566063dSJacob Faibussowitsch PetscCall(TSMPRKTableauSetUp(ts)); 10979566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&dm)); 10989566063dSJacob Faibussowitsch PetscCall(DMCoarsenHookAdd(dm,DMCoarsenHook_TSMPRK,DMRestrictHook_TSMPRK,ts)); 10999566063dSJacob Faibussowitsch PetscCall(DMSubDomainHookAdd(dm,DMSubDomainHook_TSMPRK,DMSubDomainRestrictHook_TSMPRK,ts)); 11000fe4d17eSHong Zhang if (ts->use_splitrhsfunction) { 11010fe4d17eSHong Zhang ts->ops->step = TSStep_MPRKSPLIT; 11020fe4d17eSHong Zhang ts->ops->evaluatestep = TSEvaluateStep_MPRKSPLIT; 11039566063dSJacob Faibussowitsch PetscCall(TSMPRKSetSplits(ts)); 11040fe4d17eSHong Zhang } else { 11050fe4d17eSHong Zhang ts->ops->step = TSStep_MPRK; 11060fe4d17eSHong Zhang ts->ops->evaluatestep = TSEvaluateStep_MPRK; 11070fe4d17eSHong Zhang } 11084b84eec9SHong Zhang PetscFunctionReturn(0); 11094b84eec9SHong Zhang } 11104b84eec9SHong Zhang 11114b84eec9SHong Zhang static PetscErrorCode TSSetFromOptions_MPRK(PetscOptionItems *PetscOptionsObject,TS ts) 11124b84eec9SHong Zhang { 11134b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 11144b84eec9SHong Zhang 11154b84eec9SHong Zhang PetscFunctionBegin; 11169566063dSJacob Faibussowitsch PetscCall(PetscOptionsHead(PetscOptionsObject,"PRK ODE solver options")); 11174b84eec9SHong Zhang { 11184b84eec9SHong Zhang MPRKTableauLink link; 11194b84eec9SHong Zhang PetscInt count,choice; 11204b84eec9SHong Zhang PetscBool flg; 11214b84eec9SHong Zhang const char **namelist; 11224b84eec9SHong Zhang for (link=MPRKTableauList,count=0; link; link=link->next,count++) ; 11239566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(count,(char***)&namelist)); 11244b84eec9SHong Zhang for (link=MPRKTableauList,count=0; link; link=link->next,count++) namelist[count] = link->tab.name; 11259566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-ts_mprk_type","Family of MPRK method","TSMPRKSetType",(const char*const*)namelist,count,mprk->tableau->name,&choice,&flg)); 11269566063dSJacob Faibussowitsch if (flg) PetscCall(TSMPRKSetType(ts,namelist[choice])); 11279566063dSJacob Faibussowitsch PetscCall(PetscFree(namelist)); 11284b84eec9SHong Zhang } 11299566063dSJacob Faibussowitsch PetscCall(PetscOptionsTail()); 11304b84eec9SHong Zhang PetscFunctionReturn(0); 11314b84eec9SHong Zhang } 11324b84eec9SHong Zhang 11334b84eec9SHong Zhang static PetscErrorCode TSView_MPRK(TS ts,PetscViewer viewer) 11344b84eec9SHong Zhang { 11354b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 11364b84eec9SHong Zhang PetscBool iascii; 11374b84eec9SHong Zhang 11384b84eec9SHong Zhang PetscFunctionBegin; 11399566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 11404b84eec9SHong Zhang if (iascii) { 11414b84eec9SHong Zhang MPRKTableau tab = mprk->tableau; 11424b84eec9SHong Zhang TSMPRKType mprktype; 11434b84eec9SHong Zhang char fbuf[512]; 11444b84eec9SHong Zhang char sbuf[512]; 114519c2959aSHong Zhang PetscInt i; 11469566063dSJacob Faibussowitsch PetscCall(TSMPRKGetType(ts,&mprktype)); 11479566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," MPRK type %s\n",mprktype)); 11489566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Order: %D\n",tab->order)); 114919c2959aSHong Zhang 11509566063dSJacob Faibussowitsch PetscCall(PetscFormatRealArray(fbuf,sizeof(fbuf),"% 8.6f",tab->s,tab->cf)); 11519566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Abscissa cf = %s\n",fbuf)); 11529566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Af = \n")); 115319c2959aSHong Zhang for (i=0; i<tab->s; i++) { 11549566063dSJacob Faibussowitsch PetscCall(PetscFormatRealArray(fbuf,sizeof(fbuf),"% 8.6f",tab->s,&tab->Af[i*tab->s])); 11559566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," %s\n",fbuf)); 115619c2959aSHong Zhang } 11579566063dSJacob Faibussowitsch PetscCall(PetscFormatRealArray(fbuf,sizeof(fbuf),"% 8.6f",tab->s,tab->bf)); 11589566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," bf = %s\n",fbuf)); 115919c2959aSHong Zhang 11609566063dSJacob Faibussowitsch PetscCall(PetscFormatRealArray(sbuf,sizeof(sbuf),"% 8.6f",tab->s,tab->csb)); 11619566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Abscissa csb = %s\n",sbuf)); 11629566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Asb = \n")); 116319c2959aSHong Zhang for (i=0; i<tab->s; i++) { 11649566063dSJacob Faibussowitsch PetscCall(PetscFormatRealArray(sbuf,sizeof(sbuf),"% 8.6f",tab->s,&tab->Asb[i*tab->s])); 11659566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," %s\n",sbuf)); 116619c2959aSHong Zhang } 11679566063dSJacob Faibussowitsch PetscCall(PetscFormatRealArray(sbuf,sizeof(sbuf),"% 8.6f",tab->s,tab->bsb)); 11689566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," bsb = %s\n",sbuf)); 116919c2959aSHong Zhang 117019c2959aSHong Zhang if (tab->np == 3) { 117119c2959aSHong Zhang char mbuf[512]; 11729566063dSJacob Faibussowitsch PetscCall(PetscFormatRealArray(mbuf,sizeof(mbuf),"% 8.6f",tab->s,tab->cmb)); 11739566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Abscissa cmb = %s\n",mbuf)); 11749566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Amb = \n")); 117519c2959aSHong Zhang for (i=0; i<tab->s; i++) { 11769566063dSJacob Faibussowitsch PetscCall(PetscFormatRealArray(mbuf,sizeof(mbuf),"% 8.6f",tab->s,&tab->Amb[i*tab->s])); 11779566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," %s\n",mbuf)); 117819c2959aSHong Zhang } 11799566063dSJacob Faibussowitsch PetscCall(PetscFormatRealArray(mbuf,sizeof(mbuf),"% 8.6f",tab->s,tab->bmb)); 11809566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," bmb = %s\n",mbuf)); 118119c2959aSHong Zhang } 11824b84eec9SHong Zhang } 11834b84eec9SHong Zhang PetscFunctionReturn(0); 11844b84eec9SHong Zhang } 11854b84eec9SHong Zhang 11864b84eec9SHong Zhang static PetscErrorCode TSLoad_MPRK(TS ts,PetscViewer viewer) 11874b84eec9SHong Zhang { 11884b84eec9SHong Zhang TSAdapt adapt; 11894b84eec9SHong Zhang 11904b84eec9SHong Zhang PetscFunctionBegin; 11919566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts,&adapt)); 11929566063dSJacob Faibussowitsch PetscCall(TSAdaptLoad(adapt,viewer)); 11934b84eec9SHong Zhang PetscFunctionReturn(0); 11944b84eec9SHong Zhang } 11954b84eec9SHong Zhang 11964b84eec9SHong Zhang /*@C 11974b84eec9SHong Zhang TSMPRKSetType - Set the type of MPRK scheme 11984b84eec9SHong Zhang 11990fe4d17eSHong Zhang Not collective 12004b84eec9SHong Zhang 1201d8d19677SJose E. Roman Input Parameters: 12024b84eec9SHong Zhang + ts - timestepping context 12034b84eec9SHong Zhang - mprktype - type of MPRK-scheme 12044b84eec9SHong Zhang 12054b84eec9SHong Zhang Options Database: 1206147403d9SBarry Smith . -ts_mprk_type - <pm2,p2,p3> - select the specific scheme 12074b84eec9SHong Zhang 12084b84eec9SHong Zhang Level: intermediate 12094b84eec9SHong Zhang 12104b84eec9SHong Zhang .seealso: TSMPRKGetType(), TSMPRK, TSMPRKType 12114b84eec9SHong Zhang @*/ 12124b84eec9SHong Zhang PetscErrorCode TSMPRKSetType(TS ts,TSMPRKType mprktype) 12134b84eec9SHong Zhang { 12144b84eec9SHong Zhang PetscFunctionBegin; 12154b84eec9SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 12164b84eec9SHong Zhang PetscValidCharPointer(mprktype,2); 1217cac4c232SBarry Smith PetscTryMethod(ts,"TSMPRKSetType_C",(TS,TSMPRKType),(ts,mprktype)); 12184b84eec9SHong Zhang PetscFunctionReturn(0); 12194b84eec9SHong Zhang } 12204b84eec9SHong Zhang 12214b84eec9SHong Zhang /*@C 12224b84eec9SHong Zhang TSMPRKGetType - Get the type of MPRK scheme 12234b84eec9SHong Zhang 12240fe4d17eSHong Zhang Not collective 12254b84eec9SHong Zhang 12264b84eec9SHong Zhang Input Parameter: 12274b84eec9SHong Zhang . ts - timestepping context 12284b84eec9SHong Zhang 12294b84eec9SHong Zhang Output Parameter: 12304b84eec9SHong Zhang . mprktype - type of MPRK-scheme 12314b84eec9SHong Zhang 12324b84eec9SHong Zhang Level: intermediate 12334b84eec9SHong Zhang 12344b84eec9SHong Zhang .seealso: TSMPRKGetType() 12354b84eec9SHong Zhang @*/ 12364b84eec9SHong Zhang PetscErrorCode TSMPRKGetType(TS ts,TSMPRKType *mprktype) 12374b84eec9SHong Zhang { 12384b84eec9SHong Zhang PetscFunctionBegin; 12394b84eec9SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1240cac4c232SBarry Smith PetscUseMethod(ts,"TSMPRKGetType_C",(TS,TSMPRKType*),(ts,mprktype)); 12414b84eec9SHong Zhang PetscFunctionReturn(0); 12424b84eec9SHong Zhang } 12434b84eec9SHong Zhang 12444b84eec9SHong Zhang static PetscErrorCode TSMPRKGetType_MPRK(TS ts,TSMPRKType *mprktype) 12454b84eec9SHong Zhang { 12464b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 12474b84eec9SHong Zhang 12484b84eec9SHong Zhang PetscFunctionBegin; 12494b84eec9SHong Zhang *mprktype = mprk->tableau->name; 12504b84eec9SHong Zhang PetscFunctionReturn(0); 12514b84eec9SHong Zhang } 12524b84eec9SHong Zhang 12534b84eec9SHong Zhang static PetscErrorCode TSMPRKSetType_MPRK(TS ts,TSMPRKType mprktype) 12544b84eec9SHong Zhang { 12554b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 12564b84eec9SHong Zhang PetscBool match; 12574b84eec9SHong Zhang MPRKTableauLink link; 12584b84eec9SHong Zhang 12594b84eec9SHong Zhang PetscFunctionBegin; 12604b84eec9SHong Zhang if (mprk->tableau) { 12619566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(mprk->tableau->name,mprktype,&match)); 12624b84eec9SHong Zhang if (match) PetscFunctionReturn(0); 12634b84eec9SHong Zhang } 12644b84eec9SHong Zhang for (link = MPRKTableauList; link; link=link->next) { 12659566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(link->tab.name,mprktype,&match)); 12664b84eec9SHong Zhang if (match) { 12679566063dSJacob Faibussowitsch if (ts->setupcalled) PetscCall(TSMPRKTableauReset(ts)); 12684b84eec9SHong Zhang mprk->tableau = &link->tab; 12699566063dSJacob Faibussowitsch if (ts->setupcalled) PetscCall(TSMPRKTableauSetUp(ts)); 12704b84eec9SHong Zhang PetscFunctionReturn(0); 12714b84eec9SHong Zhang } 12724b84eec9SHong Zhang } 127398921bdaSJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_UNKNOWN_TYPE,"Could not find '%s'",mprktype); 12744b84eec9SHong Zhang } 12754b84eec9SHong Zhang 12764b84eec9SHong Zhang static PetscErrorCode TSGetStages_MPRK(TS ts,PetscInt *ns,Vec **Y) 12774b84eec9SHong Zhang { 12784b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK*)ts->data; 12794b84eec9SHong Zhang 12804b84eec9SHong Zhang PetscFunctionBegin; 12814b84eec9SHong Zhang *ns = mprk->tableau->s; 12824b84eec9SHong Zhang if (Y) *Y = mprk->Y; 12834b84eec9SHong Zhang PetscFunctionReturn(0); 12844b84eec9SHong Zhang } 12854b84eec9SHong Zhang 12864b84eec9SHong Zhang static PetscErrorCode TSDestroy_MPRK(TS ts) 12874b84eec9SHong Zhang { 12884b84eec9SHong Zhang PetscFunctionBegin; 12899566063dSJacob Faibussowitsch PetscCall(TSReset_MPRK(ts)); 12904b84eec9SHong Zhang if (ts->dm) { 12919566063dSJacob Faibussowitsch PetscCall(DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSMPRK,DMRestrictHook_TSMPRK,ts)); 12929566063dSJacob Faibussowitsch PetscCall(DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSMPRK,DMSubDomainRestrictHook_TSMPRK,ts)); 12934b84eec9SHong Zhang } 12949566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data)); 12959566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSMPRKGetType_C",NULL)); 12969566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSMPRKSetType_C",NULL)); 12974b84eec9SHong Zhang PetscFunctionReturn(0); 12984b84eec9SHong Zhang } 12994b84eec9SHong Zhang 13004b84eec9SHong Zhang /*MC 1301f944a40eSHong Zhang TSMPRK - ODE solver using Multirate Partitioned Runge-Kutta schemes 13024b84eec9SHong Zhang 13034b84eec9SHong Zhang The user should provide the right hand side of the equation 13044b84eec9SHong Zhang using TSSetRHSFunction(). 13054b84eec9SHong Zhang 13064b84eec9SHong Zhang Notes: 1307f944a40eSHong Zhang The default is TSMPRKPM2, it can be changed with TSMPRKSetType() or -ts_mprk_type 13084b84eec9SHong Zhang 13094b84eec9SHong Zhang Level: beginner 13104b84eec9SHong Zhang 13114b84eec9SHong Zhang .seealso: TSCreate(), TS, TSSetType(), TSMPRKSetType(), TSMPRKGetType(), TSMPRKType, TSMPRKRegister(), TSMPRKSetMultirateType() 13124b84eec9SHong Zhang TSMPRKM2, TSMPRKM3, TSMPRKRFSMR3, TSMPRKRFSMR2 13134b84eec9SHong Zhang 13144b84eec9SHong Zhang M*/ 13154b84eec9SHong Zhang PETSC_EXTERN PetscErrorCode TSCreate_MPRK(TS ts) 13164b84eec9SHong Zhang { 13174b84eec9SHong Zhang TS_MPRK *mprk; 13184b84eec9SHong Zhang 13194b84eec9SHong Zhang PetscFunctionBegin; 13209566063dSJacob Faibussowitsch PetscCall(TSMPRKInitializePackage()); 13214b84eec9SHong Zhang 13224b84eec9SHong Zhang ts->ops->reset = TSReset_MPRK; 13234b84eec9SHong Zhang ts->ops->destroy = TSDestroy_MPRK; 13244b84eec9SHong Zhang ts->ops->view = TSView_MPRK; 13254b84eec9SHong Zhang ts->ops->load = TSLoad_MPRK; 13264b84eec9SHong Zhang ts->ops->setup = TSSetUp_MPRK; 13274b84eec9SHong Zhang ts->ops->step = TSStep_MPRK; 13284b84eec9SHong Zhang ts->ops->evaluatestep = TSEvaluateStep_MPRK; 13294b84eec9SHong Zhang ts->ops->setfromoptions = TSSetFromOptions_MPRK; 13304b84eec9SHong Zhang ts->ops->getstages = TSGetStages_MPRK; 13314b84eec9SHong Zhang 13329566063dSJacob Faibussowitsch PetscCall(PetscNewLog(ts,&mprk)); 13334b84eec9SHong Zhang ts->data = (void*)mprk; 13344b84eec9SHong Zhang 13359566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSMPRKGetType_C",TSMPRKGetType_MPRK)); 13369566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts,"TSMPRKSetType_C",TSMPRKSetType_MPRK)); 13374b84eec9SHong Zhang 13389566063dSJacob Faibussowitsch PetscCall(TSMPRKSetType(ts,TSMPRKDefault)); 13394b84eec9SHong Zhang PetscFunctionReturn(0); 13404b84eec9SHong Zhang } 1341