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 699371c9d4SSatish Balay static PetscErrorCode TSMPRKGenerateTableau2(PetscInt ratio, PetscInt s, const PetscReal Abase[], const PetscReal bbase[], PetscReal A1[], PetscReal b1[], PetscReal A2[], PetscReal b2[]) { 7019c2959aSHong Zhang PetscInt i, j, k, l; 7119c2959aSHong Zhang 7219c2959aSHong Zhang PetscFunctionBegin; 7319c2959aSHong Zhang for (k = 0; k < ratio; k++) { 7419c2959aSHong Zhang /* diagonal blocks */ 7519c2959aSHong Zhang for (i = 0; i < s; i++) 7619c2959aSHong Zhang for (j = 0; j < s; j++) { 7719c2959aSHong Zhang A1[(k * s + i) * ratio * s + k * s + j] = Abase[i * s + j]; 7819c2959aSHong Zhang A2[(k * s + i) * ratio * s + k * s + j] = Abase[i * s + j] / ratio; 7919c2959aSHong Zhang } 8019c2959aSHong Zhang /* off diagonal blocks */ 8119c2959aSHong Zhang for (l = 0; l < k; l++) 8219c2959aSHong Zhang for (i = 0; i < s; i++) 839371c9d4SSatish Balay for (j = 0; j < s; j++) A2[(k * s + i) * ratio * s + l * s + j] = bbase[j] / ratio; 8419c2959aSHong Zhang for (j = 0; j < s; j++) { 8519c2959aSHong Zhang b1[k * s + j] = bbase[j] / ratio; 8619c2959aSHong Zhang b2[k * s + j] = bbase[j] / ratio; 8719c2959aSHong Zhang } 8819c2959aSHong Zhang } 8919c2959aSHong Zhang PetscFunctionReturn(0); 9019c2959aSHong Zhang } 9119c2959aSHong Zhang 929371c9d4SSatish Balay 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[]) { 9319c2959aSHong Zhang PetscInt i, j, k, l, m, n; 9419c2959aSHong Zhang 9519c2959aSHong Zhang PetscFunctionBegin; 9619c2959aSHong Zhang for (k = 0; k < ratio; k++) { /* diagonal blocks of size ratio*s by ratio*s */ 9719c2959aSHong Zhang for (l = 0; l < ratio; l++) /* diagonal sub-blocks of size s by s */ 9819c2959aSHong Zhang for (i = 0; i < s; i++) 9919c2959aSHong Zhang for (j = 0; j < s; j++) { 10019c2959aSHong Zhang A1[((k * ratio + l) * s + i) * ratio * ratio * s + (k * ratio + l) * s + j] = Abase[i * s + j]; 10119c2959aSHong Zhang A2[((k * ratio + l) * s + i) * ratio * ratio * s + (k * ratio + l) * s + j] = Abase[i * s + j] / ratio; 10219c2959aSHong Zhang A3[((k * ratio + l) * s + i) * ratio * ratio * s + (k * ratio + l) * s + j] = Abase[i * s + j] / ratio / ratio; 10319c2959aSHong Zhang } 10419c2959aSHong Zhang for (l = 0; l < k; l++) /* off-diagonal blocks of size ratio*s by ratio*s */ 10519c2959aSHong Zhang for (m = 0; m < ratio; m++) 10619c2959aSHong Zhang for (n = 0; n < ratio; n++) 10719c2959aSHong Zhang for (i = 0; i < s; i++) 10819c2959aSHong Zhang for (j = 0; j < s; j++) { 10919c2959aSHong Zhang A2[((k * ratio + m) * s + i) * ratio * ratio * s + (l * ratio + n) * s + j] = bbase[j] / ratio / ratio; 11019c2959aSHong Zhang A3[((k * ratio + m) * s + i) * ratio * ratio * s + (l * ratio + n) * s + j] = bbase[j] / ratio / ratio; 11119c2959aSHong Zhang } 11219c2959aSHong Zhang for (m = 0; m < ratio; m++) 11319c2959aSHong Zhang for (n = 0; n < m; n++) /* off-diagonal sub-blocks of size s by s in the diagonal blocks */ 11419c2959aSHong Zhang for (i = 0; i < s; i++) 1159371c9d4SSatish Balay for (j = 0; j < s; j++) A3[((k * ratio + m) * s + i) * ratio * ratio * s + (k * ratio + n) * s + j] = bbase[j] / ratio / ratio; 11619c2959aSHong Zhang for (n = 0; n < ratio; n++) 11719c2959aSHong Zhang for (j = 0; j < s; j++) { 11819c2959aSHong Zhang b1[(k * ratio + n) * s + j] = bbase[j] / ratio / ratio; 11919c2959aSHong Zhang b2[(k * ratio + n) * s + j] = bbase[j] / ratio / ratio; 12019c2959aSHong Zhang b3[(k * ratio + n) * s + j] = bbase[j] / ratio / ratio; 12119c2959aSHong Zhang } 12219c2959aSHong Zhang } 12319c2959aSHong Zhang PetscFunctionReturn(0); 12419c2959aSHong Zhang } 12519c2959aSHong Zhang 1264b84eec9SHong Zhang /*MC 127f944a40eSHong Zhang TSMPRK2A22 - Second Order Multirate Partitioned Runge Kutta scheme based on RK2A. 1284b84eec9SHong Zhang 12919c2959aSHong Zhang This method has four stages for slow and fast parts. The refinement factor of the stepsize is 2. 13019c2959aSHong Zhang r = 2, np = 2 131147403d9SBarry Smith 1324b84eec9SHong Zhang Options database: 133147403d9SBarry Smith . -ts_mprk_type 2a22 - select this scheme 1344b84eec9SHong Zhang 1354b84eec9SHong Zhang Level: advanced 1364b84eec9SHong Zhang 137db781477SPatrick Sanan .seealso: `TSMPRK`, `TSMPRKType`, `TSMPRKSetType()` 1384b84eec9SHong Zhang M*/ 1394b84eec9SHong Zhang /*MC 140f944a40eSHong Zhang TSMPRK2A23 - Second Order Multirate Partitioned Runge-Kutta scheme based on RK2A. 14119c2959aSHong Zhang 14219c2959aSHong Zhang This method has eight stages for slow and medium and fast parts. The refinement factor of the stepsize is 2. 14319c2959aSHong Zhang r = 2, np = 3 144147403d9SBarry Smith 14519c2959aSHong Zhang Options database: 146147403d9SBarry Smith . -ts_mprk_type 2a23 - select this scheme 14719c2959aSHong Zhang 14819c2959aSHong Zhang Level: advanced 14919c2959aSHong Zhang 150db781477SPatrick Sanan .seealso: `TSMPRK`, `TSMPRKType`, `TSMPRKSetType()` 15119c2959aSHong Zhang M*/ 15219c2959aSHong Zhang /*MC 153f944a40eSHong Zhang TSMPRK2A32 - Second Order Multirate Partitioned Runge-Kutta scheme based on RK2A. 15419c2959aSHong Zhang 15519c2959aSHong Zhang This method has four stages for slow and fast parts. The refinement factor of the stepsize is 3. 15619c2959aSHong Zhang r = 3, np = 2 157147403d9SBarry Smith 15819c2959aSHong Zhang Options database: 159147403d9SBarry Smith . -ts_mprk_type 2a32 - select this scheme 16019c2959aSHong Zhang 16119c2959aSHong Zhang Level: advanced 16219c2959aSHong Zhang 163db781477SPatrick Sanan .seealso: `TSMPRK`, `TSMPRKType`, `TSMPRKSetType()` 16419c2959aSHong Zhang M*/ 16519c2959aSHong Zhang /*MC 166f944a40eSHong Zhang TSMPRK2A33 - Second Order Multirate Partitioned Runge-Kutta scheme based on RK2A. 16719c2959aSHong Zhang 16819c2959aSHong Zhang This method has eight stages for slow and medium and fast parts. The refinement factor of the stepsize is 3. 16919c2959aSHong Zhang r = 3, np = 3 170147403d9SBarry Smith 17119c2959aSHong Zhang Options database: 172147403d9SBarry Smith . -ts_mprk_type 2a33- select this scheme 17319c2959aSHong Zhang 17419c2959aSHong Zhang Level: advanced 17519c2959aSHong Zhang 176db781477SPatrick Sanan .seealso: `TSMPRK`, `TSMPRKType`, `TSMPRKSetType()` 17719c2959aSHong Zhang M*/ 17819c2959aSHong Zhang /*MC 179f944a40eSHong Zhang TSMPRK3P2M - Third Order Multirate Partitioned Runge-Kutta scheme. 1804b84eec9SHong Zhang 1814b84eec9SHong Zhang This method has eight stages for both slow and fast parts. 1824b84eec9SHong Zhang 1834b84eec9SHong Zhang Options database: 184147403d9SBarry Smith . -ts_mprk_type pm3 - select this scheme 1854b84eec9SHong Zhang 1864b84eec9SHong Zhang Level: advanced 1874b84eec9SHong Zhang 188db781477SPatrick Sanan .seealso: `TSMPRK`, `TSMPRKType`, `TSMPRKSetType()` 1894b84eec9SHong Zhang M*/ 1904b84eec9SHong Zhang /*MC 191f944a40eSHong Zhang TSMPRKP2 - Second Order Multirate Partitioned Runge-Kutta scheme. 1924b84eec9SHong Zhang 1934b84eec9SHong Zhang This method has five stages for both slow and fast parts. 1944b84eec9SHong Zhang 1954b84eec9SHong Zhang Options database: 196147403d9SBarry Smith . -ts_mprk_type p2 - select this scheme 1974b84eec9SHong Zhang 1984b84eec9SHong Zhang Level: advanced 1994b84eec9SHong Zhang 200db781477SPatrick Sanan .seealso: `TSMPRK`, `TSMPRKType`, `TSMPRKSetType()` 2014b84eec9SHong Zhang M*/ 2024b84eec9SHong Zhang /*MC 203f944a40eSHong Zhang TSMPRKP3 - Third Order Multirate Partitioned Runge-Kutta scheme. 2044b84eec9SHong Zhang 2054b84eec9SHong Zhang This method has ten stages for both slow and fast parts. 2064b84eec9SHong Zhang 2074b84eec9SHong Zhang Options database: 208147403d9SBarry Smith . -ts_mprk_type p3 - select this scheme 2094b84eec9SHong Zhang 2104b84eec9SHong Zhang Level: advanced 2114b84eec9SHong Zhang 212db781477SPatrick Sanan .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 222db781477SPatrick Sanan .seealso: `TSMPRKRegisterDestroy()` 2234b84eec9SHong Zhang @*/ 2249371c9d4SSatish Balay PetscErrorCode TSMPRKRegisterAll(void) { 2254b84eec9SHong Zhang PetscFunctionBegin; 2264b84eec9SHong Zhang if (TSMPRKRegisterAllCalled) PetscFunctionReturn(0); 2274b84eec9SHong Zhang TSMPRKRegisterAllCalled = PETSC_TRUE; 2284b84eec9SHong Zhang 2294b84eec9SHong Zhang #define RC PetscRealConstant 2304b84eec9SHong Zhang { 2319371c9d4SSatish Balay const PetscReal Abase[2][2] = 2329371c9d4SSatish Balay { 2339371c9d4SSatish Balay {0, 0}, 2349371c9d4SSatish Balay {RC(1.0), 0} 2359371c9d4SSatish Balay }, 23619c2959aSHong Zhang bbase[2] = {RC(0.5), RC(0.5)}; 2379371c9d4SSatish Balay PetscReal Asb[4][4] = {{0}}, Af[4][4] = {{0}}, bsb[4] = {0}, bf[4] = {0}; 2389371c9d4SSatish Balay PetscInt rsb[4] = {0, 0, 1, 2}; 2399566063dSJacob Faibussowitsch PetscCall(TSMPRKGenerateTableau2(2, 2, &Abase[0][0], bbase, &Asb[0][0], bsb, &Af[0][0], bf)); 2409566063dSJacob Faibussowitsch PetscCall(TSMPRKRegister(TSMPRK2A22, 2, 2, 2, 1, &Asb[0][0], bsb, NULL, rsb, NULL, NULL, NULL, NULL, &Af[0][0], bf, NULL)); 2414b84eec9SHong Zhang } 24219c2959aSHong Zhang { 2439371c9d4SSatish Balay const PetscReal Abase[2][2] = 2449371c9d4SSatish Balay { 2459371c9d4SSatish Balay {0, 0}, 2469371c9d4SSatish Balay {RC(1.0), 0} 2479371c9d4SSatish Balay }, 24819c2959aSHong Zhang bbase[2] = {RC(0.5), RC(0.5)}; 2499371c9d4SSatish Balay PetscReal Asb[8][8] = {{0}}, Amb[8][8] = {{0}}, Af[8][8] = {{0}}, bsb[8] = {0}, bmb[8] = {0}, bf[8] = {0}; 2509371c9d4SSatish Balay PetscInt rsb[8] = {0, 0, 1, 2, 1, 2, 1, 2}, rmb[8] = {0, 0, 1, 2, 0, 0, 5, 6}; 2519566063dSJacob Faibussowitsch PetscCall(TSMPRKGenerateTableau3(2, 2, &Abase[0][0], bbase, &Asb[0][0], bsb, &Amb[0][0], bmb, &Af[0][0], bf)); 2529566063dSJacob 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)); 25319c2959aSHong Zhang } 25419c2959aSHong Zhang { 2559371c9d4SSatish Balay const PetscReal Abase[2][2] = 2569371c9d4SSatish Balay { 2579371c9d4SSatish Balay {0, 0}, 2589371c9d4SSatish Balay {RC(1.0), 0} 2599371c9d4SSatish Balay }, 26019c2959aSHong Zhang bbase[2] = {RC(0.5), RC(0.5)}; 2619371c9d4SSatish Balay PetscReal Asb[6][6] = {{0}}, Af[6][6] = {{0}}, bsb[6] = {0}, bf[6] = {0}; 2629371c9d4SSatish Balay PetscInt rsb[6] = {0, 0, 1, 2, 1, 2}; 2639566063dSJacob Faibussowitsch PetscCall(TSMPRKGenerateTableau2(3, 2, &Abase[0][0], bbase, &Asb[0][0], bsb, &Af[0][0], bf)); 2649566063dSJacob Faibussowitsch PetscCall(TSMPRKRegister(TSMPRK2A32, 2, 2, 3, 1, &Asb[0][0], bsb, NULL, rsb, NULL, NULL, NULL, NULL, &Af[0][0], bf, NULL)); 26519c2959aSHong Zhang } 26619c2959aSHong Zhang { 2679371c9d4SSatish Balay const PetscReal Abase[2][2] = 2689371c9d4SSatish Balay { 2699371c9d4SSatish Balay {0, 0}, 2709371c9d4SSatish Balay {RC(1.0), 0} 2719371c9d4SSatish Balay }, 27219c2959aSHong Zhang bbase[2] = {RC(0.5), RC(0.5)}; 2739371c9d4SSatish Balay PetscReal Asb[18][18] = {{0}}, Amb[18][18] = {{0}}, Af[18][18] = {{0}}, bsb[18] = {0}, bmb[18] = {0}, bf[18] = {0}; 2749371c9d4SSatish Balay PetscInt 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}; 2759566063dSJacob Faibussowitsch PetscCall(TSMPRKGenerateTableau3(3, 2, &Abase[0][0], bbase, &Asb[0][0], bsb, &Amb[0][0], bmb, &Af[0][0], bf)); 2769566063dSJacob 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)); 27719c2959aSHong Zhang } 27819c2959aSHong Zhang /* 27919c2959aSHong Zhang PetscReal 28019c2959aSHong Zhang Asb[8][8] = {{Abase[0][0],Abase[0][1],0,0,0,0,0,0}, 28119c2959aSHong Zhang {Abase[1][0],Abase[1][1],0,0,0,0,0,0}, 28219c2959aSHong Zhang {0,0,Abase[0][0],Abase[0][1],0,0,0,0}, 28319c2959aSHong Zhang {0,0,Abase[1][0],Abase[1][1],0,0,0,0}, 28419c2959aSHong Zhang {0,0,0,0,Abase[0][0],Abase[0][1],0,0}, 28519c2959aSHong Zhang {0,0,0,0,Abase[1][0],Abase[1][1],0,0}, 28619c2959aSHong Zhang {0,0,0,0,0,0,Abase[0][0],Abase[0][1]}, 28719c2959aSHong Zhang {0,0,0,0,0,0,Abase[1][0],Abase[1][1]}}, 28819c2959aSHong Zhang Amb[8][8] = {{Abase[0][0]/m,Abase[0][1]/m,0,0,0,0,0,0}, 28919c2959aSHong Zhang {Abase[1][0]/m,Abase[1][1]/m,0,0,0,0,0,0}, 29019c2959aSHong Zhang {0,0,Abase[0][0]/m,Abase[0][1]/m,0,0,0,0}, 29119c2959aSHong Zhang {0,0,Abase[1][0]/m,Abase[1][1]/m,0,0,0,0}, 29219c2959aSHong Zhang {bbase[0]/m,bbase[1]/m,bbase[0]/m,bbase[1]/m,Abase[0][0]/m,Abase[0][1]/m,0,0}, 29319c2959aSHong Zhang {bbase[0]/m,bbase[1]/m,bbase[0]/m,bbase[1]/m,Abase[1][0]/m,Abase[1][1]/m,0,0}, 29419c2959aSHong Zhang {bbase[0]/m,bbase[1]/m,bbase[0]/m,bbase[1]/m,0,0,Abase[0][0]/m,Abase[0][1]/m}, 29519c2959aSHong Zhang {bbase[0]/m,bbase[1]/m,bbase[0]/m,bbase[1]/m,0,0,Abase[1][0]/m,Abase[1][1]/m}}, 29619c2959aSHong Zhang Af[8][8] = {{Abase[0][0]/m/m,Abase[0][1]/m/m,0,0,0,0,0,0}, 29719c2959aSHong Zhang {Abase[1][0]/m/m,Abase[1][1]/m/m,0,0,0,0,0,0}, 29819c2959aSHong Zhang {0,0,Abase[0][0]/m/m,Abase[0][1]/m/m,0,0,0,0}, 29919c2959aSHong Zhang {0,0,Abase[1][0]/m/m,Abase[1][1]/m/m,0,0,0,0}, 30019c2959aSHong 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}, 30119c2959aSHong 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}, 30219c2959aSHong 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}, 30319c2959aSHong 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}}, 30419c2959aSHong 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}, 30519c2959aSHong 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}, 30619c2959aSHong 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}, 30719c2959aSHong Zhang */ 3084b84eec9SHong Zhang /*{ 3094b84eec9SHong Zhang const PetscReal 3104b84eec9SHong Zhang As[8][8] = {{0,0,0,0,0,0,0,0}, 3114b84eec9SHong Zhang {RC(1.0)/RC(2.0),0,0,0,0,0,0,0}, 3124b84eec9SHong Zhang {RC(-1.0)/RC(6.0),RC(2.0)/RC(3.0),0,0,0,0,0,0}, 3134b84eec9SHong Zhang {RC(1.0)/RC(3.0),RC(-1.0)/RC(3.0),RC(1.0),0,0,0,0,0}, 3144b84eec9SHong Zhang {0,0,0,0,0,0,0,0}, 3154b84eec9SHong Zhang {0,0,0,0,RC(1.0)/RC(2.0),0,0,0}, 3164b84eec9SHong Zhang {0,0,0,0,RC(-1.0)/RC(6.0),RC(2.0)/RC(3.0),0,0}, 3174b84eec9SHong Zhang {0,0,0,0,RC(1.0)/RC(3.0),RC(-1.0)/RC(3.0),RC(1.0),0}}, 3184b84eec9SHong Zhang A[8][8] = {{0,0,0,0,0,0,0,0}, 3194b84eec9SHong Zhang {RC(1.0)/RC(4.0),0,0,0,0,0,0,0}, 3204b84eec9SHong Zhang {RC(-1.0)/RC(12.0),RC(1.0)/RC(3.0),0,0,0,0,0,0}, 3214b84eec9SHong Zhang {RC(1.0)/RC(6.0),RC(-1.0)/RC(6.0),RC(1.0)/RC(2.0),0,0,0,0,0}, 3224b84eec9SHong 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}, 3234b84eec9SHong 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}, 3244b84eec9SHong 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}, 3254b84eec9SHong Zhang {RC(1.0)/RC(12.0),RC(1.0)/RC(6.0),RC(1.0)/RC(6.0),RC(1.0)/RC(12.0),RC(1.0)/RC(6.0),RC(-1.0)/RC(6.0),RC(1.0)/RC(2.0),0}}, 3264b84eec9SHong 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)}, 3274b84eec9SHong 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)}; 3289566063dSJacob Faibussowitsch PetscCall(TSMPRKRegister(TSMPRKPM3,3,8,&As[0][0],bs,NULL,&A[0][0],b,NULL)); 3294b84eec9SHong Zhang }*/ 3304b84eec9SHong Zhang 3314b84eec9SHong Zhang { 3329371c9d4SSatish Balay const PetscReal Asb[5][5] = 3339371c9d4SSatish Balay { 3349371c9d4SSatish Balay {0, 0, 0, 0, 0}, 3354b84eec9SHong Zhang {RC(1.0) / RC(2.0), 0, 0, 0, 0}, 3364b84eec9SHong Zhang {RC(1.0) / RC(2.0), 0, 0, 0, 0}, 3374b84eec9SHong Zhang {RC(1.0), 0, 0, 0, 0}, 3389371c9d4SSatish Balay {RC(1.0), 0, 0, 0, 0} 3399371c9d4SSatish Balay }, 3409371c9d4SSatish Balay Af[5][5] = {{0, 0, 0, 0, 0}, {RC(1.0) / RC(2.0), 0, 0, 0, 0}, {RC(1.0) / RC(4.0), RC(1.0) / RC(4.0), 0, 0, 0}, {RC(1.0) / RC(4.0), RC(1.0) / RC(4.0), RC(1.0) / RC(2.0), 0, 0}, {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}}, bsb[5] = {RC(1.0) / RC(2.0), 0, 0, 0, RC(1.0) / RC(2.0)}, 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}; 3419371c9d4SSatish Balay const PetscInt rsb[5] = {0, 0, 2, 0, 4}; 3429566063dSJacob Faibussowitsch PetscCall(TSMPRKRegister(TSMPRKP2, 2, 5, 1, 1, &Asb[0][0], bsb, NULL, rsb, NULL, NULL, NULL, NULL, &Af[0][0], bf, NULL)); 3434b84eec9SHong Zhang } 3444b84eec9SHong Zhang 3454b84eec9SHong Zhang { 3469371c9d4SSatish Balay const PetscReal Asb[10][10] = 3479371c9d4SSatish Balay { 3489371c9d4SSatish Balay {0, 0, 0, 0, 0, 0, 0, 0, 0, 0}, 3494b84eec9SHong Zhang {RC(1.0) / RC(4.0), 0, 0, 0, 0, 0, 0, 0, 0, 0}, 3504b84eec9SHong Zhang {RC(1.0) / RC(4.0), 0, 0, 0, 0, 0, 0, 0, 0, 0}, 3514b84eec9SHong Zhang {RC(1.0) / RC(2.0), 0, 0, 0, 0, 0, 0, 0, 0, 0}, 3524b84eec9SHong Zhang {RC(1.0) / RC(2.0), 0, 0, 0, 0, 0, 0, 0, 0, 0}, 3534b84eec9SHong Zhang {RC(-1.0) / RC(6.0), 0, 0, 0, RC(2.0) / RC(3.0), 0, 0, 0, 0, 0}, 3544b84eec9SHong 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}, 3554b84eec9SHong 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}, 3564b84eec9SHong Zhang {RC(1.0) / RC(3.0), 0, 0, 0, RC(-1.0) / RC(3.0), RC(1.0), 0, 0, 0, 0}, 3579371c9d4SSatish Balay {RC(1.0) / RC(3.0), 0, 0, 0, RC(-1.0) / RC(3.0), RC(1.0), 0, 0, 0, 0} 3589371c9d4SSatish Balay }, 3599371c9d4SSatish Balay Af[10][10] = {{0, 0, 0, 0, 0, 0, 0, 0, 0, 0}, {RC(1.0) / RC(4.0), 0, 0, 0, 0, 0, 0, 0, 0, 0}, {RC(-1.0) / RC(12.0), RC(1.0) / RC(3.0), 0, 0, 0, 0, 0, 0, 0, 0}, {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}, {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}, {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}, {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}, {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}, {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}, {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}}, 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)}, 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}; 3609371c9d4SSatish Balay const PetscInt rsb[10] = {0, 0, 2, 0, 4, 0, 0, 7, 0, 9}; 3619566063dSJacob Faibussowitsch PetscCall(TSMPRKRegister(TSMPRKP3, 3, 5, 2, 1, &Asb[0][0], bsb, NULL, rsb, NULL, NULL, NULL, NULL, &Af[0][0], bf, NULL)); 3624b84eec9SHong Zhang } 3634b84eec9SHong Zhang #undef RC 3644b84eec9SHong Zhang PetscFunctionReturn(0); 3654b84eec9SHong Zhang } 3664b84eec9SHong Zhang 3674b84eec9SHong Zhang /*@C 3684b84eec9SHong Zhang TSMPRKRegisterDestroy - Frees the list of schemes that were registered by TSMPRKRegister(). 3694b84eec9SHong Zhang 3704b84eec9SHong Zhang Not Collective 3714b84eec9SHong Zhang 3724b84eec9SHong Zhang Level: advanced 3734b84eec9SHong Zhang 374db781477SPatrick Sanan .seealso: `TSMPRKRegister()`, `TSMPRKRegisterAll()` 3754b84eec9SHong Zhang @*/ 3769371c9d4SSatish Balay PetscErrorCode TSMPRKRegisterDestroy(void) { 3774b84eec9SHong Zhang MPRKTableauLink link; 3784b84eec9SHong Zhang 3794b84eec9SHong Zhang PetscFunctionBegin; 3804b84eec9SHong Zhang while ((link = MPRKTableauList)) { 3814b84eec9SHong Zhang MPRKTableau t = &link->tab; 3824b84eec9SHong Zhang MPRKTableauList = link->next; 3839566063dSJacob Faibussowitsch PetscCall(PetscFree3(t->Asb, t->bsb, t->csb)); 3849566063dSJacob Faibussowitsch PetscCall(PetscFree3(t->Amb, t->bmb, t->cmb)); 3859566063dSJacob Faibussowitsch PetscCall(PetscFree3(t->Af, t->bf, t->cf)); 3869566063dSJacob Faibussowitsch PetscCall(PetscFree(t->rsb)); 3879566063dSJacob Faibussowitsch PetscCall(PetscFree(t->rmb)); 3889566063dSJacob Faibussowitsch PetscCall(PetscFree(t->name)); 3899566063dSJacob Faibussowitsch PetscCall(PetscFree(link)); 3904b84eec9SHong Zhang } 3914b84eec9SHong Zhang TSMPRKRegisterAllCalled = PETSC_FALSE; 3924b84eec9SHong Zhang PetscFunctionReturn(0); 3934b84eec9SHong Zhang } 3944b84eec9SHong Zhang 3954b84eec9SHong Zhang /*@C 3964b84eec9SHong Zhang TSMPRKInitializePackage - This function initializes everything in the TSMPRK package. It is called 3974b84eec9SHong Zhang from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to TSCreate_MPRK() 3984b84eec9SHong Zhang when using static libraries. 3994b84eec9SHong Zhang 4004b84eec9SHong Zhang Level: developer 4014b84eec9SHong Zhang 402db781477SPatrick Sanan .seealso: `PetscInitialize()` 4034b84eec9SHong Zhang @*/ 4049371c9d4SSatish Balay PetscErrorCode TSMPRKInitializePackage(void) { 4054b84eec9SHong Zhang PetscFunctionBegin; 4064b84eec9SHong Zhang if (TSMPRKPackageInitialized) PetscFunctionReturn(0); 4074b84eec9SHong Zhang TSMPRKPackageInitialized = PETSC_TRUE; 4089566063dSJacob Faibussowitsch PetscCall(TSMPRKRegisterAll()); 4099566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalize(TSMPRKFinalizePackage)); 4104b84eec9SHong Zhang PetscFunctionReturn(0); 4114b84eec9SHong Zhang } 4124b84eec9SHong Zhang 4134b84eec9SHong Zhang /*@C 414f944a40eSHong Zhang TSMPRKFinalizePackage - This function destroys everything in the TSMPRK package. It is 4154b84eec9SHong Zhang called from PetscFinalize(). 4164b84eec9SHong Zhang 4174b84eec9SHong Zhang Level: developer 4184b84eec9SHong Zhang 419db781477SPatrick Sanan .seealso: `PetscFinalize()` 4204b84eec9SHong Zhang @*/ 4219371c9d4SSatish Balay PetscErrorCode TSMPRKFinalizePackage(void) { 4224b84eec9SHong Zhang PetscFunctionBegin; 4234b84eec9SHong Zhang TSMPRKPackageInitialized = PETSC_FALSE; 4249566063dSJacob Faibussowitsch PetscCall(TSMPRKRegisterDestroy()); 4254b84eec9SHong Zhang PetscFunctionReturn(0); 4264b84eec9SHong Zhang } 4274b84eec9SHong Zhang 4284b84eec9SHong Zhang /*@C 4294b84eec9SHong Zhang TSMPRKRegister - register a MPRK scheme by providing the entries in the Butcher tableau 4304b84eec9SHong Zhang 4314b84eec9SHong Zhang Not Collective, but the same schemes should be registered on all processes on which they will be used 4324b84eec9SHong Zhang 4334b84eec9SHong Zhang Input Parameters: 4344b84eec9SHong Zhang + name - identifier for method 4354b84eec9SHong Zhang . order - approximation order of method 43679bc8a77SHong Zhang . s - number of stages in the base methods 43779bc8a77SHong Zhang . ratio1 - stepsize ratio at 1st level (e.g. slow/medium) 43879bc8a77SHong Zhang . ratio2 - stepsize ratio at 2nd level (e.g. medium/fast) 4394b84eec9SHong Zhang . Af - stage coefficients for fast components(dimension s*s, row-major) 4404b84eec9SHong Zhang . bf - step completion table for fast components(dimension s) 4414b84eec9SHong Zhang . cf - abscissa for fast components(dimension s) 4424b84eec9SHong Zhang . As - stage coefficients for slow components(dimension s*s, row-major) 4434b84eec9SHong Zhang . bs - step completion table for slow components(dimension s) 4444b84eec9SHong Zhang - cs - abscissa for slow components(dimension s) 4454b84eec9SHong Zhang 4464b84eec9SHong Zhang Notes: 4474b84eec9SHong Zhang Several MPRK methods are provided, this function is only needed to create new methods. 4484b84eec9SHong Zhang 4494b84eec9SHong Zhang Level: advanced 4504b84eec9SHong Zhang 451db781477SPatrick Sanan .seealso: `TSMPRK` 4524b84eec9SHong Zhang @*/ 4539371c9d4SSatish Balay PetscErrorCode TSMPRKRegister(TSMPRKType name, PetscInt order, PetscInt sbase, PetscInt ratio1, PetscInt ratio2, const PetscReal Asb[], const PetscReal bsb[], const PetscReal csb[], const PetscInt rsb[], const PetscReal Amb[], const PetscReal bmb[], const PetscReal cmb[], const PetscInt rmb[], const PetscReal Af[], const PetscReal bf[], const PetscReal cf[]) { 4544b84eec9SHong Zhang MPRKTableauLink link; 4554b84eec9SHong Zhang MPRKTableau t; 45679bc8a77SHong Zhang PetscInt s, i, j; 4574b84eec9SHong Zhang 4584b84eec9SHong Zhang PetscFunctionBegin; 4594b84eec9SHong Zhang PetscValidCharPointer(name, 1); 460064a246eSJacob Faibussowitsch PetscValidRealPointer(Asb, 6); 46179bc8a77SHong Zhang if (bsb) PetscValidRealPointer(bsb, 7); 46279bc8a77SHong Zhang if (csb) PetscValidRealPointer(csb, 8); 463064a246eSJacob Faibussowitsch if (rsb) PetscValidIntPointer(rsb, 9); 46479bc8a77SHong Zhang if (Amb) PetscValidRealPointer(Amb, 10); 46579bc8a77SHong Zhang if (bmb) PetscValidRealPointer(bmb, 11); 46679bc8a77SHong Zhang if (cmb) PetscValidRealPointer(cmb, 12); 467064a246eSJacob Faibussowitsch if (rmb) PetscValidIntPointer(rmb, 13); 46879bc8a77SHong Zhang PetscValidRealPointer(Af, 14); 46979bc8a77SHong Zhang if (bf) PetscValidRealPointer(bf, 15); 47079bc8a77SHong Zhang if (cf) PetscValidRealPointer(cf, 16); 4714b84eec9SHong Zhang 4729566063dSJacob Faibussowitsch PetscCall(PetscNew(&link)); 4734b84eec9SHong Zhang t = &link->tab; 4744b84eec9SHong Zhang 4759566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(name, &t->name)); 47679bc8a77SHong Zhang s = sbase * ratio1 * ratio2; /* this is the dimension of the matrices below */ 4774b84eec9SHong Zhang t->order = order; 47879bc8a77SHong Zhang t->sbase = sbase; 4794b84eec9SHong Zhang t->s = s; 48019c2959aSHong Zhang t->np = 2; 48179bc8a77SHong Zhang 4829566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(s * s, &t->Af, s, &t->bf, s, &t->cf)); 4839566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->Af, Af, s * s)); 4844b84eec9SHong Zhang if (bf) { 4859566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->bf, bf, s)); 48619c2959aSHong Zhang } else 4874b84eec9SHong Zhang for (i = 0; i < s; i++) t->bf[i] = Af[(s - 1) * s + i]; 4884b84eec9SHong Zhang if (cf) { 4899566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->cf, cf, s)); 49019c2959aSHong Zhang } else { 4914b84eec9SHong Zhang for (i = 0; i < s; i++) 4929371c9d4SSatish Balay for (j = 0, t->cf[i] = 0; j < s; j++) t->cf[i] += Af[i * s + j]; 4934b84eec9SHong Zhang } 49419c2959aSHong Zhang 49519c2959aSHong Zhang if (Amb) { 49619c2959aSHong Zhang t->np = 3; 4979566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(s * s, &t->Amb, s, &t->bmb, s, &t->cmb)); 4989566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->Amb, Amb, s * s)); 49919c2959aSHong Zhang if (bmb) { 5009566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->bmb, bmb, s)); 50119c2959aSHong Zhang } else { 50219c2959aSHong Zhang for (i = 0; i < s; i++) t->bmb[i] = Amb[(s - 1) * s + i]; 5034b84eec9SHong Zhang } 50419c2959aSHong Zhang if (cmb) { 5059566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->cmb, cmb, s)); 50619c2959aSHong Zhang } else { 5074b84eec9SHong Zhang for (i = 0; i < s; i++) 5089371c9d4SSatish Balay for (j = 0, t->cmb[i] = 0; j < s; j++) t->cmb[i] += Amb[i * s + j]; 50919c2959aSHong Zhang } 51019c2959aSHong Zhang if (rmb) { 5119566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(s, &t->rmb)); 5129566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->rmb, rmb, s)); 513071fcb05SBarry Smith } else { 5149566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(s, &t->rmb)); 51519c2959aSHong Zhang } 51619c2959aSHong Zhang } 51719c2959aSHong Zhang 5189566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(s * s, &t->Asb, s, &t->bsb, s, &t->csb)); 5199566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->Asb, Asb, s * s)); 52019c2959aSHong Zhang if (bsb) { 5219566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->bsb, bsb, s)); 52219c2959aSHong Zhang } else 52319c2959aSHong Zhang for (i = 0; i < s; i++) t->bsb[i] = Asb[(s - 1) * s + i]; 52419c2959aSHong Zhang if (csb) { 5259566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->csb, csb, s)); 52619c2959aSHong Zhang } else { 52719c2959aSHong Zhang for (i = 0; i < s; i++) 5289371c9d4SSatish Balay for (j = 0, t->csb[i] = 0; j < s; j++) t->csb[i] += Asb[i * s + j]; 52919c2959aSHong Zhang } 53019c2959aSHong Zhang if (rsb) { 5319566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(s, &t->rsb)); 5329566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t->rsb, rsb, s)); 533071fcb05SBarry Smith } else { 5349566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(s, &t->rsb)); 5354b84eec9SHong Zhang } 5364b84eec9SHong Zhang link->next = MPRKTableauList; 5374b84eec9SHong Zhang MPRKTableauList = link; 5384b84eec9SHong Zhang PetscFunctionReturn(0); 5394b84eec9SHong Zhang } 5404b84eec9SHong Zhang 5419371c9d4SSatish Balay static PetscErrorCode TSMPRKSetSplits(TS ts) { 5424b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK *)ts->data; 54319c2959aSHong Zhang MPRKTableau tab = mprk->tableau; 5444b84eec9SHong Zhang DM dm, subdm, newdm; 5454b84eec9SHong Zhang 5464b84eec9SHong Zhang PetscFunctionBegin; 5479566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetSubTS(ts, "slow", &mprk->subts_slow)); 5489566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetSubTS(ts, "fast", &mprk->subts_fast)); 5493c633725SBarry 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"); 5504b84eec9SHong Zhang 55119c2959aSHong Zhang /* Only copy the DM */ 5529566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 55319c2959aSHong Zhang 5549566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetSubTS(ts, "slowbuffer", &mprk->subts_slowbuffer)); 55519c2959aSHong Zhang if (!mprk->subts_slowbuffer) { 55619c2959aSHong Zhang mprk->subts_slowbuffer = mprk->subts_slow; 55719c2959aSHong Zhang mprk->subts_slow = NULL; 55819c2959aSHong Zhang } 55919c2959aSHong Zhang if (mprk->subts_slow) { 5609566063dSJacob Faibussowitsch PetscCall(DMClone(dm, &newdm)); 5619566063dSJacob Faibussowitsch PetscCall(TSGetDM(mprk->subts_slow, &subdm)); 5629566063dSJacob Faibussowitsch PetscCall(DMCopyDMTS(subdm, newdm)); 5639566063dSJacob Faibussowitsch PetscCall(DMCopyDMSNES(subdm, newdm)); 5649566063dSJacob Faibussowitsch PetscCall(TSSetDM(mprk->subts_slow, newdm)); 5659566063dSJacob Faibussowitsch PetscCall(DMDestroy(&newdm)); 56619c2959aSHong Zhang } 5679566063dSJacob Faibussowitsch PetscCall(DMClone(dm, &newdm)); 5689566063dSJacob Faibussowitsch PetscCall(TSGetDM(mprk->subts_slowbuffer, &subdm)); 5699566063dSJacob Faibussowitsch PetscCall(DMCopyDMTS(subdm, newdm)); 5709566063dSJacob Faibussowitsch PetscCall(DMCopyDMSNES(subdm, newdm)); 5719566063dSJacob Faibussowitsch PetscCall(TSSetDM(mprk->subts_slowbuffer, newdm)); 5729566063dSJacob Faibussowitsch PetscCall(DMDestroy(&newdm)); 57319c2959aSHong Zhang 5749566063dSJacob Faibussowitsch PetscCall(DMClone(dm, &newdm)); 5759566063dSJacob Faibussowitsch PetscCall(TSGetDM(mprk->subts_fast, &subdm)); 5769566063dSJacob Faibussowitsch PetscCall(DMCopyDMTS(subdm, newdm)); 5779566063dSJacob Faibussowitsch PetscCall(DMCopyDMSNES(subdm, newdm)); 5789566063dSJacob Faibussowitsch PetscCall(TSSetDM(mprk->subts_fast, newdm)); 5799566063dSJacob Faibussowitsch PetscCall(DMDestroy(&newdm)); 58019c2959aSHong Zhang 58119c2959aSHong Zhang if (tab->np == 3) { 5829566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetSubTS(ts, "medium", &mprk->subts_medium)); 5839566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetSubTS(ts, "mediumbuffer", &mprk->subts_mediumbuffer)); 58419c2959aSHong Zhang if (mprk->subts_medium && !mprk->subts_mediumbuffer) { 58519c2959aSHong Zhang mprk->subts_mediumbuffer = mprk->subts_medium; 58619c2959aSHong Zhang mprk->subts_medium = NULL; 58719c2959aSHong Zhang } 58819c2959aSHong Zhang if (mprk->subts_medium) { 5899566063dSJacob Faibussowitsch PetscCall(DMClone(dm, &newdm)); 5909566063dSJacob Faibussowitsch PetscCall(TSGetDM(mprk->subts_medium, &subdm)); 5919566063dSJacob Faibussowitsch PetscCall(DMCopyDMTS(subdm, newdm)); 5929566063dSJacob Faibussowitsch PetscCall(DMCopyDMSNES(subdm, newdm)); 5939566063dSJacob Faibussowitsch PetscCall(TSSetDM(mprk->subts_medium, newdm)); 5949566063dSJacob Faibussowitsch PetscCall(DMDestroy(&newdm)); 59519c2959aSHong Zhang } 5969566063dSJacob Faibussowitsch PetscCall(DMClone(dm, &newdm)); 5979566063dSJacob Faibussowitsch PetscCall(TSGetDM(mprk->subts_mediumbuffer, &subdm)); 5989566063dSJacob Faibussowitsch PetscCall(DMCopyDMTS(subdm, newdm)); 5999566063dSJacob Faibussowitsch PetscCall(DMCopyDMSNES(subdm, newdm)); 6009566063dSJacob Faibussowitsch PetscCall(TSSetDM(mprk->subts_mediumbuffer, newdm)); 6019566063dSJacob Faibussowitsch PetscCall(DMDestroy(&newdm)); 60219c2959aSHong Zhang } 6034b84eec9SHong Zhang PetscFunctionReturn(0); 6044b84eec9SHong Zhang } 6054b84eec9SHong Zhang 6064b84eec9SHong Zhang /* 6074b84eec9SHong Zhang This if for nonsplit RHS MPRK 6084b84eec9SHong Zhang The step completion formula is 6094b84eec9SHong Zhang 6104b84eec9SHong Zhang x1 = x0 + h b^T YdotRHS 6114b84eec9SHong Zhang 6124b84eec9SHong Zhang */ 6139371c9d4SSatish Balay static PetscErrorCode TSEvaluateStep_MPRK(TS ts, PetscInt order, Vec X, PetscBool *done) { 6144b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK *)ts->data; 6154b84eec9SHong Zhang MPRKTableau tab = mprk->tableau; 6167c0df07dSHong Zhang PetscScalar *wf = mprk->work_fast; 6174b84eec9SHong Zhang PetscReal h = ts->time_step; 6184b84eec9SHong Zhang PetscInt s = tab->s, j; 6194b84eec9SHong Zhang 6204b84eec9SHong Zhang PetscFunctionBegin; 6214b84eec9SHong Zhang for (j = 0; j < s; j++) wf[j] = h * tab->bf[j]; 6229566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, X)); 6239566063dSJacob Faibussowitsch PetscCall(VecMAXPY(X, s, wf, mprk->YdotRHS)); 6244b84eec9SHong Zhang PetscFunctionReturn(0); 6254b84eec9SHong Zhang } 6264b84eec9SHong Zhang 6279371c9d4SSatish Balay static PetscErrorCode TSStep_MPRK(TS ts) { 6284b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK *)ts->data; 6297c0df07dSHong Zhang Vec *Y = mprk->Y, *YdotRHS = mprk->YdotRHS, *YdotRHS_fast = mprk->YdotRHS_fast, *YdotRHS_slow = mprk->YdotRHS_slow, *YdotRHS_slowbuffer = mprk->YdotRHS_slowbuffer; 6307c0df07dSHong Zhang Vec Yslow, Yslowbuffer, Yfast; 6314b84eec9SHong Zhang MPRKTableau tab = mprk->tableau; 6324b84eec9SHong Zhang const PetscInt s = tab->s; 63319c2959aSHong Zhang const PetscReal *Af = tab->Af, *cf = tab->cf, *Asb = tab->Asb, *csb = tab->csb; 634ebd5ed4eSHong Zhang PetscScalar *wf = mprk->work_fast, *wsb = mprk->work_slowbuffer; 635ebd5ed4eSHong Zhang PetscInt i, j; 6364b84eec9SHong Zhang PetscReal next_time_step = ts->time_step, t = ts->ptime, h = ts->time_step; 6374b84eec9SHong Zhang 6384b84eec9SHong Zhang PetscFunctionBegin; 6394b84eec9SHong Zhang for (i = 0; i < s; i++) { 6404b84eec9SHong Zhang mprk->stage_time = t + h * cf[i]; 6419566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, mprk->stage_time)); 6429566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, Y[i])); 6434b84eec9SHong Zhang 6447c0df07dSHong Zhang /* slow buffer region */ 64519c2959aSHong Zhang for (j = 0; j < i; j++) wsb[j] = h * Asb[i * s + j]; 646*48a46eb9SPierre Jolivet for (j = 0; j < i; j++) PetscCall(VecGetSubVector(YdotRHS[j], mprk->is_slowbuffer, &YdotRHS_slowbuffer[j])); 6479566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i], mprk->is_slowbuffer, &Yslowbuffer)); 6489566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Yslowbuffer, i, wsb, mprk->YdotRHS_slowbuffer)); 6499566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i], mprk->is_slowbuffer, &Yslowbuffer)); 650*48a46eb9SPierre Jolivet for (j = 0; j < i; j++) PetscCall(VecRestoreSubVector(YdotRHS[j], mprk->is_slowbuffer, &YdotRHS_slowbuffer[j])); 6517c0df07dSHong Zhang /* slow region */ 6527c0df07dSHong Zhang if (mprk->is_slow) { 653*48a46eb9SPierre Jolivet for (j = 0; j < i; j++) PetscCall(VecGetSubVector(YdotRHS[j], mprk->is_slow, &YdotRHS_slow[j])); 6549566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i], mprk->is_slow, &Yslow)); 6559566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Yslow, i, wsb, mprk->YdotRHS_slow)); 6569566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i], mprk->is_slow, &Yslow)); 657*48a46eb9SPierre Jolivet for (j = 0; j < i; j++) PetscCall(VecRestoreSubVector(YdotRHS[j], mprk->is_slow, &YdotRHS_slow[j])); 6587c0df07dSHong Zhang } 6594b84eec9SHong Zhang 6607c0df07dSHong Zhang /* fast region */ 6614b84eec9SHong Zhang for (j = 0; j < i; j++) wf[j] = h * Af[i * s + j]; 662*48a46eb9SPierre Jolivet for (j = 0; j < i; j++) PetscCall(VecGetSubVector(YdotRHS[j], mprk->is_fast, &YdotRHS_fast[j])); 6639566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i], mprk->is_fast, &Yfast)); 6649566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Yfast, i, wf, mprk->YdotRHS_fast)); 6659566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i], mprk->is_fast, &Yfast)); 666*48a46eb9SPierre Jolivet for (j = 0; j < i; j++) PetscCall(VecRestoreSubVector(YdotRHS[j], mprk->is_fast, &YdotRHS_fast[j])); 6677c0df07dSHong Zhang if (tab->np == 3) { 6687c0df07dSHong Zhang Vec *YdotRHS_medium = mprk->YdotRHS_medium, *YdotRHS_mediumbuffer = mprk->YdotRHS_mediumbuffer; 6697c0df07dSHong Zhang Vec Ymedium, Ymediumbuffer; 670ebd5ed4eSHong Zhang PetscScalar *wmb = mprk->work_mediumbuffer; 671ebd5ed4eSHong Zhang 672ebd5ed4eSHong Zhang for (j = 0; j < i; j++) wmb[j] = h * tab->Amb[i * s + j]; 6737c0df07dSHong Zhang /* medium region */ 6747c0df07dSHong Zhang if (mprk->is_medium) { 675*48a46eb9SPierre Jolivet for (j = 0; j < i; j++) PetscCall(VecGetSubVector(YdotRHS[j], mprk->is_medium, &YdotRHS_medium[j])); 6769566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i], mprk->is_medium, &Ymedium)); 6779566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Ymedium, i, wmb, mprk->YdotRHS_medium)); 6789566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i], mprk->is_medium, &Ymedium)); 679*48a46eb9SPierre Jolivet for (j = 0; j < i; j++) PetscCall(VecRestoreSubVector(YdotRHS[j], mprk->is_medium, &YdotRHS_medium[j])); 6807c0df07dSHong Zhang } 6817c0df07dSHong Zhang /* medium buffer region */ 682*48a46eb9SPierre Jolivet for (j = 0; j < i; j++) PetscCall(VecGetSubVector(YdotRHS[j], mprk->is_mediumbuffer, &YdotRHS_mediumbuffer[j])); 6839566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i], mprk->is_mediumbuffer, &Ymediumbuffer)); 6849566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Ymediumbuffer, i, wmb, mprk->YdotRHS_mediumbuffer)); 6859566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i], mprk->is_mediumbuffer, &Ymediumbuffer)); 686*48a46eb9SPierre Jolivet for (j = 0; j < i; j++) PetscCall(VecRestoreSubVector(YdotRHS[j], mprk->is_mediumbuffer, &YdotRHS_mediumbuffer[j])); 6877c0df07dSHong Zhang } 6889566063dSJacob Faibussowitsch PetscCall(TSPostStage(ts, mprk->stage_time, i, Y)); 6894b84eec9SHong Zhang /* compute the stage RHS by fast and slow tableau respectively */ 6909566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(ts, t + h * csb[i], Y[i], YdotRHS[i])); 6914b84eec9SHong Zhang } 6929566063dSJacob Faibussowitsch PetscCall(TSEvaluateStep(ts, tab->order, ts->vec_sol, NULL)); 6934b84eec9SHong Zhang ts->ptime += ts->time_step; 6944b84eec9SHong Zhang ts->time_step = next_time_step; 6954b84eec9SHong Zhang PetscFunctionReturn(0); 6964b84eec9SHong Zhang } 6974b84eec9SHong Zhang 6984b84eec9SHong Zhang /* 699f944a40eSHong Zhang This if for the case when split RHS is used 7004b84eec9SHong Zhang The step completion formula is 7014b84eec9SHong Zhang x1 = x0 + h b^T YdotRHS 7024b84eec9SHong Zhang */ 7039371c9d4SSatish Balay static PetscErrorCode TSEvaluateStep_MPRKSPLIT(TS ts, PetscInt order, Vec X, PetscBool *done) { 7044b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK *)ts->data; 7054b84eec9SHong Zhang MPRKTableau tab = mprk->tableau; 7066aad120cSJose E. Roman Vec Xslow, Xfast, Xslowbuffer; /* subvectors for slow and fast components in X respectively */ 70719c2959aSHong Zhang PetscScalar *wf = mprk->work_fast, *ws = mprk->work_slow, *wsb = mprk->work_slowbuffer; 7084b84eec9SHong Zhang PetscReal h = ts->time_step; 70979bc8a77SHong Zhang PetscInt s = tab->s, j, computedstages; 7104b84eec9SHong Zhang 7114b84eec9SHong Zhang PetscFunctionBegin; 7129566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, X)); 71319c2959aSHong Zhang 71419c2959aSHong Zhang /* slow region */ 71519c2959aSHong Zhang if (mprk->is_slow) { 71679bc8a77SHong Zhang computedstages = 0; 71719c2959aSHong Zhang for (j = 0; j < s; j++) { 7189849be05SHong Zhang if (tab->rsb[j]) ws[tab->rsb[j] - 1] += h * tab->bsb[j]; 71979bc8a77SHong Zhang else ws[computedstages++] = h * tab->bsb[j]; 72019c2959aSHong Zhang } 7219566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(X, mprk->is_slow, &Xslow)); 7229566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Xslow, computedstages, ws, mprk->YdotRHS_slow)); 7239566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(X, mprk->is_slow, &Xslow)); 72419c2959aSHong Zhang } 72519c2959aSHong Zhang 7269d6e09e9SHong Zhang if (tab->np == 3 && mprk->is_medium) { 7279d6e09e9SHong Zhang computedstages = 0; 7289d6e09e9SHong Zhang for (j = 0; j < s; j++) { 7299d6e09e9SHong Zhang if (tab->rmb[j]) wsb[computedstages - tab->sbase + (tab->rmb[j] - 1) % tab->sbase] += h * tab->bsb[j]; 7309d6e09e9SHong Zhang else wsb[computedstages++] = h * tab->bsb[j]; 7319d6e09e9SHong Zhang } 7329566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(X, mprk->is_slowbuffer, &Xslowbuffer)); 7339566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Xslowbuffer, computedstages, wsb, mprk->YdotRHS_slowbuffer)); 7349566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(X, mprk->is_slowbuffer, &Xslowbuffer)); 7359d6e09e9SHong Zhang } else { 73619c2959aSHong Zhang /* slow buffer region */ 73719c2959aSHong Zhang for (j = 0; j < s; j++) wsb[j] = h * tab->bsb[j]; 7389566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(X, mprk->is_slowbuffer, &Xslowbuffer)); 7399566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Xslowbuffer, s, wsb, mprk->YdotRHS_slowbuffer)); 7409566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(X, mprk->is_slowbuffer, &Xslowbuffer)); 7419d6e09e9SHong Zhang } 74219c2959aSHong Zhang if (tab->np == 3) { 74319c2959aSHong Zhang Vec Xmedium, Xmediumbuffer; 74419c2959aSHong Zhang PetscScalar *wm = mprk->work_medium, *wmb = mprk->work_mediumbuffer; 7459d6e09e9SHong Zhang /* medium region and slow buffer region */ 74619c2959aSHong Zhang if (mprk->is_medium) { 74779bc8a77SHong Zhang computedstages = 0; 74819c2959aSHong Zhang for (j = 0; j < s; j++) { 74979bc8a77SHong Zhang if (tab->rmb[j]) wm[computedstages - tab->sbase + (tab->rmb[j] - 1) % tab->sbase] += h * tab->bmb[j]; 75079bc8a77SHong Zhang else wm[computedstages++] = h * tab->bmb[j]; 75119c2959aSHong Zhang } 7529566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(X, mprk->is_medium, &Xmedium)); 7539566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Xmedium, computedstages, wm, mprk->YdotRHS_medium)); 7549566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(X, mprk->is_medium, &Xmedium)); 75519c2959aSHong Zhang } 75619c2959aSHong Zhang /* medium buffer region */ 75719c2959aSHong Zhang for (j = 0; j < s; j++) wmb[j] = h * tab->bmb[j]; 7589566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(X, mprk->is_mediumbuffer, &Xmediumbuffer)); 7599566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Xmediumbuffer, s, wmb, mprk->YdotRHS_mediumbuffer)); 7609566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(X, mprk->is_mediumbuffer, &Xmediumbuffer)); 76119c2959aSHong Zhang } 76219c2959aSHong Zhang /* fast region */ 76319c2959aSHong Zhang for (j = 0; j < s; j++) wf[j] = h * tab->bf[j]; 7649566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(X, mprk->is_fast, &Xfast)); 7659566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Xfast, s, wf, mprk->YdotRHS_fast)); 7669566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(X, mprk->is_fast, &Xfast)); 7674b84eec9SHong Zhang PetscFunctionReturn(0); 7684b84eec9SHong Zhang } 7694b84eec9SHong Zhang 7709371c9d4SSatish Balay static PetscErrorCode TSStep_MPRKSPLIT(TS ts) { 7714b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK *)ts->data; 7724b84eec9SHong Zhang MPRKTableau tab = mprk->tableau; 77319c2959aSHong Zhang Vec *Y = mprk->Y, *YdotRHS_fast = mprk->YdotRHS_fast, *YdotRHS_slow = mprk->YdotRHS_slow, *YdotRHS_slowbuffer = mprk->YdotRHS_slowbuffer; 77419c2959aSHong Zhang Vec Yslow, Yslowbuffer, Yfast; /* subvectors for slow and fast components in Y[i] respectively */ 77519c2959aSHong Zhang PetscInt s = tab->s; 77619c2959aSHong Zhang const PetscReal *Af = tab->Af, *cf = tab->cf, *Asb = tab->Asb, *csb = tab->csb; 77719c2959aSHong Zhang PetscScalar *wf = mprk->work_fast, *ws = mprk->work_slow, *wsb = mprk->work_slowbuffer; 77879bc8a77SHong Zhang PetscInt i, j, computedstages; 7794b84eec9SHong Zhang PetscReal next_time_step = ts->time_step, t = ts->ptime, h = ts->time_step; 7804b84eec9SHong Zhang 7814b84eec9SHong Zhang PetscFunctionBegin; 7824b84eec9SHong Zhang for (i = 0; i < s; i++) { 7834b84eec9SHong Zhang mprk->stage_time = t + h * cf[i]; 7849566063dSJacob Faibussowitsch PetscCall(TSPreStage(ts, mprk->stage_time)); 7854b84eec9SHong Zhang /* calculate the stage value for fast and slow components respectively */ 7869566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_sol, Y[i])); 78719c2959aSHong Zhang for (j = 0; j < i; j++) wsb[j] = h * Asb[i * s + j]; 78819c2959aSHong Zhang 78919c2959aSHong Zhang /* slow buffer region */ 7909d6e09e9SHong Zhang if (tab->np == 3 && mprk->is_medium) { 7919d6e09e9SHong Zhang if (tab->rmb[i]) { 7929566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i], mprk->is_slowbuffer, &Yslowbuffer)); 7939566063dSJacob Faibussowitsch PetscCall(VecISCopy(Y[tab->rmb[i] - 1], mprk->is_slowbuffer, SCATTER_REVERSE, Yslowbuffer)); 7949566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i], mprk->is_slowbuffer, &Yslowbuffer)); 7959d6e09e9SHong Zhang } else { 7969d6e09e9SHong Zhang PetscScalar *wm = mprk->work_medium; 7979d6e09e9SHong Zhang computedstages = 0; 7989d6e09e9SHong Zhang for (j = 0; j < i; j++) { 7999d6e09e9SHong Zhang if (tab->rmb[j]) wm[computedstages - tab->sbase + (tab->rmb[j] - 1) % tab->sbase] += wsb[j]; 8009d6e09e9SHong Zhang else wm[computedstages++] = wsb[j]; 8019d6e09e9SHong Zhang } 8029566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i], mprk->is_slowbuffer, &Yslowbuffer)); 8039566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Yslowbuffer, computedstages, wm, YdotRHS_slowbuffer)); 8049566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i], mprk->is_slowbuffer, &Yslowbuffer)); 8059d6e09e9SHong Zhang } 8069d6e09e9SHong Zhang } else { 8079566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i], mprk->is_slowbuffer, &Yslowbuffer)); 8089566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Yslowbuffer, i, wsb, YdotRHS_slowbuffer)); 8099566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i], mprk->is_slowbuffer, &Yslowbuffer)); 8109d6e09e9SHong Zhang } 8119849be05SHong Zhang 81219c2959aSHong Zhang /* slow region */ 8139849be05SHong Zhang if (mprk->is_slow) { 8149849be05SHong Zhang if (tab->rsb[i]) { /* repeat previous stage */ 8159566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i], mprk->is_slow, &Yslow)); 8169566063dSJacob Faibussowitsch PetscCall(VecISCopy(Y[tab->rsb[i] - 1], mprk->is_slow, SCATTER_REVERSE, Yslow)); 8179566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i], mprk->is_slow, &Yslow)); 8189849be05SHong Zhang } else { 81979bc8a77SHong Zhang computedstages = 0; 82019c2959aSHong Zhang for (j = 0; j < i; j++) { 82179bc8a77SHong Zhang if (tab->rsb[j]) ws[tab->rsb[j] - 1] += wsb[j]; 82279bc8a77SHong Zhang else ws[computedstages++] = wsb[j]; 82319c2959aSHong Zhang } 8249566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i], mprk->is_slow, &Yslow)); 8259566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Yslow, computedstages, ws, YdotRHS_slow)); 8269566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i], mprk->is_slow, &Yslow)); 82719c2959aSHong Zhang /* only depends on the slow buffer region */ 8289566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(mprk->subts_slow, t + h * csb[i], Y[i], YdotRHS_slow[computedstages])); 82919c2959aSHong Zhang } 8309849be05SHong Zhang } 83119c2959aSHong Zhang 83219c2959aSHong Zhang /* fast region */ 83319c2959aSHong Zhang for (j = 0; j < i; j++) wf[j] = h * Af[i * s + j]; 8349566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i], mprk->is_fast, &Yfast)); 8359566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Yfast, i, wf, YdotRHS_fast)); 8369566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i], mprk->is_fast, &Yfast)); 83719c2959aSHong Zhang 83819c2959aSHong Zhang if (tab->np == 3) { 83919c2959aSHong Zhang Vec *YdotRHS_medium = mprk->YdotRHS_medium, *YdotRHS_mediumbuffer = mprk->YdotRHS_mediumbuffer; 84019c2959aSHong Zhang Vec Ymedium, Ymediumbuffer; 84119c2959aSHong Zhang const PetscReal *Amb = tab->Amb, *cmb = tab->cmb; 84219c2959aSHong Zhang PetscScalar *wm = mprk->work_medium, *wmb = mprk->work_mediumbuffer; 84319c2959aSHong Zhang 84419c2959aSHong Zhang for (j = 0; j < i; j++) wmb[j] = h * Amb[i * s + j]; 84519c2959aSHong Zhang /* medium buffer region */ 8469566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i], mprk->is_mediumbuffer, &Ymediumbuffer)); 8479566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Ymediumbuffer, i, wmb, YdotRHS_mediumbuffer)); 8489566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i], mprk->is_mediumbuffer, &Ymediumbuffer)); 84919c2959aSHong Zhang /* medium region */ 85079bc8a77SHong Zhang if (mprk->is_medium) { 85179bc8a77SHong Zhang if (tab->rmb[i]) { /* repeat previous stage */ 8529566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i], mprk->is_medium, &Ymedium)); 8539566063dSJacob Faibussowitsch PetscCall(VecISCopy(Y[tab->rmb[i] - 1], mprk->is_medium, SCATTER_REVERSE, Ymedium)); 8549566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i], mprk->is_medium, &Ymedium)); 85579bc8a77SHong Zhang } else { 85679bc8a77SHong Zhang computedstages = 0; 85779bc8a77SHong Zhang for (j = 0; j < i; j++) { 85879bc8a77SHong Zhang if (tab->rmb[j]) wm[computedstages - tab->sbase + (tab->rmb[j] - 1) % tab->sbase] += wmb[j]; 85979bc8a77SHong Zhang else wm[computedstages++] = wmb[j]; 86079bc8a77SHong Zhang } 8619566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(Y[i], mprk->is_medium, &Ymedium)); 8629566063dSJacob Faibussowitsch PetscCall(VecMAXPY(Ymedium, computedstages, wm, YdotRHS_medium)); 8639566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(Y[i], mprk->is_medium, &Ymedium)); 86479bc8a77SHong Zhang /* only depends on the medium buffer region */ 8659566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(mprk->subts_medium, t + h * cmb[i], Y[i], YdotRHS_medium[computedstages])); 8669d6e09e9SHong Zhang /* must be computed after all regions are updated in Y */ 8679566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(mprk->subts_slowbuffer, t + h * csb[i], Y[i], YdotRHS_slowbuffer[computedstages])); 86879bc8a77SHong Zhang } 86919c2959aSHong Zhang } 87019c2959aSHong Zhang /* must be computed after fast region and slow region are updated in Y */ 8719566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(mprk->subts_mediumbuffer, t + h * cmb[i], Y[i], YdotRHS_mediumbuffer[i])); 87219c2959aSHong Zhang } 873*48a46eb9SPierre Jolivet if (!(tab->np == 3 && mprk->is_medium)) PetscCall(TSComputeRHSFunction(mprk->subts_slowbuffer, t + h * csb[i], Y[i], YdotRHS_slowbuffer[i])); 8749566063dSJacob Faibussowitsch PetscCall(TSComputeRHSFunction(mprk->subts_fast, t + h * cf[i], Y[i], YdotRHS_fast[i])); 8754b84eec9SHong Zhang } 87679bc8a77SHong Zhang 8779566063dSJacob Faibussowitsch PetscCall(TSEvaluateStep(ts, tab->order, ts->vec_sol, NULL)); 8784b84eec9SHong Zhang ts->ptime += ts->time_step; 8794b84eec9SHong Zhang ts->time_step = next_time_step; 8804b84eec9SHong Zhang PetscFunctionReturn(0); 8814b84eec9SHong Zhang } 8824b84eec9SHong Zhang 8839371c9d4SSatish Balay static PetscErrorCode TSMPRKTableauReset(TS ts) { 8844b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK *)ts->data; 8854b84eec9SHong Zhang MPRKTableau tab = mprk->tableau; 8864b84eec9SHong Zhang 8874b84eec9SHong Zhang PetscFunctionBegin; 8884b84eec9SHong Zhang if (!tab) PetscFunctionReturn(0); 8899566063dSJacob Faibussowitsch PetscCall(PetscFree(mprk->work_fast)); 8909566063dSJacob Faibussowitsch PetscCall(PetscFree(mprk->work_slow)); 8919566063dSJacob Faibussowitsch PetscCall(PetscFree(mprk->work_slowbuffer)); 8929566063dSJacob Faibussowitsch PetscCall(PetscFree(mprk->work_medium)); 8939566063dSJacob Faibussowitsch PetscCall(PetscFree(mprk->work_mediumbuffer)); 8949566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s, &mprk->Y)); 8950fe4d17eSHong Zhang if (ts->use_splitrhsfunction) { 8969566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s, &mprk->YdotRHS_fast)); 8979566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s, &mprk->YdotRHS_slow)); 8989566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s, &mprk->YdotRHS_slowbuffer)); 8999566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s, &mprk->YdotRHS_medium)); 9009566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s, &mprk->YdotRHS_mediumbuffer)); 9010fe4d17eSHong Zhang } else { 9029566063dSJacob Faibussowitsch PetscCall(VecDestroyVecs(tab->s, &mprk->YdotRHS)); 9031baa6e33SBarry Smith if (mprk->is_slow) PetscCall(PetscFree(mprk->YdotRHS_slow)); 9049566063dSJacob Faibussowitsch PetscCall(PetscFree(mprk->YdotRHS_slowbuffer)); 9057c0df07dSHong Zhang if (tab->np == 3) { 9061baa6e33SBarry Smith if (mprk->is_medium) PetscCall(PetscFree(mprk->YdotRHS_medium)); 9079566063dSJacob Faibussowitsch PetscCall(PetscFree(mprk->YdotRHS_mediumbuffer)); 9087c0df07dSHong Zhang } 9099566063dSJacob Faibussowitsch PetscCall(PetscFree(mprk->YdotRHS_fast)); 9107c0df07dSHong Zhang } 9114b84eec9SHong Zhang PetscFunctionReturn(0); 9124b84eec9SHong Zhang } 9134b84eec9SHong Zhang 9149371c9d4SSatish Balay static PetscErrorCode TSReset_MPRK(TS ts) { 9154b84eec9SHong Zhang PetscFunctionBegin; 9169566063dSJacob Faibussowitsch PetscCall(TSMPRKTableauReset(ts)); 9174b84eec9SHong Zhang PetscFunctionReturn(0); 9184b84eec9SHong Zhang } 9194b84eec9SHong Zhang 9209371c9d4SSatish Balay static PetscErrorCode DMCoarsenHook_TSMPRK(DM fine, DM coarse, void *ctx) { 9214b84eec9SHong Zhang PetscFunctionBegin; 9224b84eec9SHong Zhang PetscFunctionReturn(0); 9234b84eec9SHong Zhang } 9244b84eec9SHong Zhang 9259371c9d4SSatish Balay static PetscErrorCode DMRestrictHook_TSMPRK(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx) { 9264b84eec9SHong Zhang PetscFunctionBegin; 9274b84eec9SHong Zhang PetscFunctionReturn(0); 9284b84eec9SHong Zhang } 9294b84eec9SHong Zhang 9309371c9d4SSatish Balay static PetscErrorCode DMSubDomainHook_TSMPRK(DM dm, DM subdm, void *ctx) { 9314b84eec9SHong Zhang PetscFunctionBegin; 9324b84eec9SHong Zhang PetscFunctionReturn(0); 9334b84eec9SHong Zhang } 9344b84eec9SHong Zhang 9359371c9d4SSatish Balay static PetscErrorCode DMSubDomainRestrictHook_TSMPRK(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx) { 9364b84eec9SHong Zhang PetscFunctionBegin; 9374b84eec9SHong Zhang PetscFunctionReturn(0); 9384b84eec9SHong Zhang } 9394b84eec9SHong Zhang 9409371c9d4SSatish Balay static PetscErrorCode TSMPRKTableauSetUp(TS ts) { 9414b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK *)ts->data; 9424b84eec9SHong Zhang MPRKTableau tab = mprk->tableau; 94319c2959aSHong Zhang Vec YdotRHS_slow, YdotRHS_slowbuffer, YdotRHS_medium, YdotRHS_mediumbuffer, YdotRHS_fast; 9444b84eec9SHong Zhang 9454b84eec9SHong Zhang PetscFunctionBegin; 9469566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &mprk->Y)); 947*48a46eb9SPierre Jolivet if (mprk->is_slow) PetscCall(PetscMalloc1(tab->s, &mprk->work_slow)); 9489566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tab->s, &mprk->work_slowbuffer)); 9497c0df07dSHong Zhang if (tab->np == 3) { 950*48a46eb9SPierre Jolivet if (mprk->is_medium) PetscCall(PetscMalloc1(tab->s, &mprk->work_medium)); 9519566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tab->s, &mprk->work_mediumbuffer)); 9527c0df07dSHong Zhang } 9539566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tab->s, &mprk->work_fast)); 9547c0df07dSHong Zhang 9550fe4d17eSHong Zhang if (ts->use_splitrhsfunction) { 95619c2959aSHong Zhang if (mprk->is_slow) { 9579566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(ts->vec_sol, mprk->is_slow, &YdotRHS_slow)); 9589566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(YdotRHS_slow, tab->s, &mprk->YdotRHS_slow)); 9599566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(ts->vec_sol, mprk->is_slow, &YdotRHS_slow)); 96019c2959aSHong Zhang } 9619566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(ts->vec_sol, mprk->is_slowbuffer, &YdotRHS_slowbuffer)); 9629566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(YdotRHS_slowbuffer, tab->s, &mprk->YdotRHS_slowbuffer)); 9639566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(ts->vec_sol, mprk->is_slowbuffer, &YdotRHS_slowbuffer)); 96419c2959aSHong Zhang if (tab->np == 3) { 96519c2959aSHong Zhang if (mprk->is_medium) { 9669566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(ts->vec_sol, mprk->is_medium, &YdotRHS_medium)); 9679566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(YdotRHS_medium, tab->s, &mprk->YdotRHS_medium)); 9689566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(ts->vec_sol, mprk->is_medium, &YdotRHS_medium)); 96919c2959aSHong Zhang } 9709566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(ts->vec_sol, mprk->is_mediumbuffer, &YdotRHS_mediumbuffer)); 9719566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(YdotRHS_mediumbuffer, tab->s, &mprk->YdotRHS_mediumbuffer)); 9729566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(ts->vec_sol, mprk->is_mediumbuffer, &YdotRHS_mediumbuffer)); 97319c2959aSHong Zhang } 9749566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(ts->vec_sol, mprk->is_fast, &YdotRHS_fast)); 9759566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(YdotRHS_fast, tab->s, &mprk->YdotRHS_fast)); 9769566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(ts->vec_sol, mprk->is_fast, &YdotRHS_fast)); 9770fe4d17eSHong Zhang } else { 9789566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &mprk->YdotRHS)); 979*48a46eb9SPierre Jolivet if (mprk->is_slow) PetscCall(PetscMalloc1(tab->s, &mprk->YdotRHS_slow)); 9809566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tab->s, &mprk->YdotRHS_slowbuffer)); 9810fe4d17eSHong Zhang if (tab->np == 3) { 982*48a46eb9SPierre Jolivet if (mprk->is_medium) PetscCall(PetscMalloc1(tab->s, &mprk->YdotRHS_medium)); 9839566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tab->s, &mprk->YdotRHS_mediumbuffer)); 9840fe4d17eSHong Zhang } 9859566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(tab->s, &mprk->YdotRHS_fast)); 9864b84eec9SHong Zhang } 9874b84eec9SHong Zhang PetscFunctionReturn(0); 9884b84eec9SHong Zhang } 9894b84eec9SHong Zhang 9909371c9d4SSatish Balay static PetscErrorCode TSSetUp_MPRK(TS ts) { 9914b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK *)ts->data; 99219c2959aSHong Zhang MPRKTableau tab = mprk->tableau; 9934b84eec9SHong Zhang DM dm; 9944b84eec9SHong Zhang 9954b84eec9SHong Zhang PetscFunctionBegin; 9969566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetIS(ts, "slow", &mprk->is_slow)); 9979566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetIS(ts, "fast", &mprk->is_fast)); 9983c633725SBarry 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); 99919c2959aSHong Zhang 100019c2959aSHong Zhang if (tab->np == 3) { 10019566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetIS(ts, "medium", &mprk->is_medium)); 10023c633725SBarry 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); 10039566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetIS(ts, "mediumbuffer", &mprk->is_mediumbuffer)); 100419c2959aSHong Zhang if (!mprk->is_mediumbuffer) { /* let medium buffer cover whole medium region */ 100519c2959aSHong Zhang mprk->is_mediumbuffer = mprk->is_medium; 100619c2959aSHong Zhang mprk->is_medium = NULL; 100719c2959aSHong Zhang } 100819c2959aSHong Zhang } 100919c2959aSHong Zhang 101019c2959aSHong Zhang /* If users do not provide buffer region settings, the solver will do them automatically, but with a performance penalty */ 10119566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetIS(ts, "slowbuffer", &mprk->is_slowbuffer)); 101219c2959aSHong Zhang if (!mprk->is_slowbuffer) { /* let slow buffer cover whole slow region */ 101319c2959aSHong Zhang mprk->is_slowbuffer = mprk->is_slow; 101419c2959aSHong Zhang mprk->is_slow = NULL; 101519c2959aSHong Zhang } 10169566063dSJacob Faibussowitsch PetscCall(TSCheckImplicitTerm(ts)); 10179566063dSJacob Faibussowitsch PetscCall(TSMPRKTableauSetUp(ts)); 10189566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 10199566063dSJacob Faibussowitsch PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_TSMPRK, DMRestrictHook_TSMPRK, ts)); 10209566063dSJacob Faibussowitsch PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_TSMPRK, DMSubDomainRestrictHook_TSMPRK, ts)); 10210fe4d17eSHong Zhang if (ts->use_splitrhsfunction) { 10220fe4d17eSHong Zhang ts->ops->step = TSStep_MPRKSPLIT; 10230fe4d17eSHong Zhang ts->ops->evaluatestep = TSEvaluateStep_MPRKSPLIT; 10249566063dSJacob Faibussowitsch PetscCall(TSMPRKSetSplits(ts)); 10250fe4d17eSHong Zhang } else { 10260fe4d17eSHong Zhang ts->ops->step = TSStep_MPRK; 10270fe4d17eSHong Zhang ts->ops->evaluatestep = TSEvaluateStep_MPRK; 10280fe4d17eSHong Zhang } 10294b84eec9SHong Zhang PetscFunctionReturn(0); 10304b84eec9SHong Zhang } 10314b84eec9SHong Zhang 10329371c9d4SSatish Balay static PetscErrorCode TSSetFromOptions_MPRK(TS ts, PetscOptionItems *PetscOptionsObject) { 10334b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK *)ts->data; 10344b84eec9SHong Zhang 10354b84eec9SHong Zhang PetscFunctionBegin; 1036d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "PRK ODE solver options"); 10374b84eec9SHong Zhang { 10384b84eec9SHong Zhang MPRKTableauLink link; 10394b84eec9SHong Zhang PetscInt count, choice; 10404b84eec9SHong Zhang PetscBool flg; 10414b84eec9SHong Zhang const char **namelist; 10429371c9d4SSatish Balay for (link = MPRKTableauList, count = 0; link; link = link->next, count++) 10439371c9d4SSatish Balay ; 10449566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(count, (char ***)&namelist)); 10454b84eec9SHong Zhang for (link = MPRKTableauList, count = 0; link; link = link->next, count++) namelist[count] = link->tab.name; 10469566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-ts_mprk_type", "Family of MPRK method", "TSMPRKSetType", (const char *const *)namelist, count, mprk->tableau->name, &choice, &flg)); 10479566063dSJacob Faibussowitsch if (flg) PetscCall(TSMPRKSetType(ts, namelist[choice])); 10489566063dSJacob Faibussowitsch PetscCall(PetscFree(namelist)); 10494b84eec9SHong Zhang } 1050d0609cedSBarry Smith PetscOptionsHeadEnd(); 10514b84eec9SHong Zhang PetscFunctionReturn(0); 10524b84eec9SHong Zhang } 10534b84eec9SHong Zhang 10549371c9d4SSatish Balay static PetscErrorCode TSView_MPRK(TS ts, PetscViewer viewer) { 10554b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK *)ts->data; 10564b84eec9SHong Zhang PetscBool iascii; 10574b84eec9SHong Zhang 10584b84eec9SHong Zhang PetscFunctionBegin; 10599566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 10604b84eec9SHong Zhang if (iascii) { 10614b84eec9SHong Zhang MPRKTableau tab = mprk->tableau; 10624b84eec9SHong Zhang TSMPRKType mprktype; 10634b84eec9SHong Zhang char fbuf[512]; 10644b84eec9SHong Zhang char sbuf[512]; 106519c2959aSHong Zhang PetscInt i; 10669566063dSJacob Faibussowitsch PetscCall(TSMPRKGetType(ts, &mprktype)); 10679566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " MPRK type %s\n", mprktype)); 106863a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Order: %" PetscInt_FMT "\n", tab->order)); 106919c2959aSHong Zhang 10709566063dSJacob Faibussowitsch PetscCall(PetscFormatRealArray(fbuf, sizeof(fbuf), "% 8.6f", tab->s, tab->cf)); 10719566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Abscissa cf = %s\n", fbuf)); 10729566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Af = \n")); 107319c2959aSHong Zhang for (i = 0; i < tab->s; i++) { 10749566063dSJacob Faibussowitsch PetscCall(PetscFormatRealArray(fbuf, sizeof(fbuf), "% 8.6f", tab->s, &tab->Af[i * tab->s])); 10759566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " %s\n", fbuf)); 107619c2959aSHong Zhang } 10779566063dSJacob Faibussowitsch PetscCall(PetscFormatRealArray(fbuf, sizeof(fbuf), "% 8.6f", tab->s, tab->bf)); 10789566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " bf = %s\n", fbuf)); 107919c2959aSHong Zhang 10809566063dSJacob Faibussowitsch PetscCall(PetscFormatRealArray(sbuf, sizeof(sbuf), "% 8.6f", tab->s, tab->csb)); 10819566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Abscissa csb = %s\n", sbuf)); 10829566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Asb = \n")); 108319c2959aSHong Zhang for (i = 0; i < tab->s; i++) { 10849566063dSJacob Faibussowitsch PetscCall(PetscFormatRealArray(sbuf, sizeof(sbuf), "% 8.6f", tab->s, &tab->Asb[i * tab->s])); 10859566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " %s\n", sbuf)); 108619c2959aSHong Zhang } 10879566063dSJacob Faibussowitsch PetscCall(PetscFormatRealArray(sbuf, sizeof(sbuf), "% 8.6f", tab->s, tab->bsb)); 10889566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " bsb = %s\n", sbuf)); 108919c2959aSHong Zhang 109019c2959aSHong Zhang if (tab->np == 3) { 109119c2959aSHong Zhang char mbuf[512]; 10929566063dSJacob Faibussowitsch PetscCall(PetscFormatRealArray(mbuf, sizeof(mbuf), "% 8.6f", tab->s, tab->cmb)); 10939566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Abscissa cmb = %s\n", mbuf)); 10949566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Amb = \n")); 109519c2959aSHong Zhang for (i = 0; i < tab->s; i++) { 10969566063dSJacob Faibussowitsch PetscCall(PetscFormatRealArray(mbuf, sizeof(mbuf), "% 8.6f", tab->s, &tab->Amb[i * tab->s])); 10979566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " %s\n", mbuf)); 109819c2959aSHong Zhang } 10999566063dSJacob Faibussowitsch PetscCall(PetscFormatRealArray(mbuf, sizeof(mbuf), "% 8.6f", tab->s, tab->bmb)); 11009566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " bmb = %s\n", mbuf)); 110119c2959aSHong Zhang } 11024b84eec9SHong Zhang } 11034b84eec9SHong Zhang PetscFunctionReturn(0); 11044b84eec9SHong Zhang } 11054b84eec9SHong Zhang 11069371c9d4SSatish Balay static PetscErrorCode TSLoad_MPRK(TS ts, PetscViewer viewer) { 11074b84eec9SHong Zhang TSAdapt adapt; 11084b84eec9SHong Zhang 11094b84eec9SHong Zhang PetscFunctionBegin; 11109566063dSJacob Faibussowitsch PetscCall(TSGetAdapt(ts, &adapt)); 11119566063dSJacob Faibussowitsch PetscCall(TSAdaptLoad(adapt, viewer)); 11124b84eec9SHong Zhang PetscFunctionReturn(0); 11134b84eec9SHong Zhang } 11144b84eec9SHong Zhang 11154b84eec9SHong Zhang /*@C 11164b84eec9SHong Zhang TSMPRKSetType - Set the type of MPRK scheme 11174b84eec9SHong Zhang 11180fe4d17eSHong Zhang Not collective 11194b84eec9SHong Zhang 1120d8d19677SJose E. Roman Input Parameters: 11214b84eec9SHong Zhang + ts - timestepping context 11224b84eec9SHong Zhang - mprktype - type of MPRK-scheme 11234b84eec9SHong Zhang 11244b84eec9SHong Zhang Options Database: 1125147403d9SBarry Smith . -ts_mprk_type - <pm2,p2,p3> - select the specific scheme 11264b84eec9SHong Zhang 11274b84eec9SHong Zhang Level: intermediate 11284b84eec9SHong Zhang 1129db781477SPatrick Sanan .seealso: `TSMPRKGetType()`, `TSMPRK`, `TSMPRKType` 11304b84eec9SHong Zhang @*/ 11319371c9d4SSatish Balay PetscErrorCode TSMPRKSetType(TS ts, TSMPRKType mprktype) { 11324b84eec9SHong Zhang PetscFunctionBegin; 11334b84eec9SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 11344b84eec9SHong Zhang PetscValidCharPointer(mprktype, 2); 1135cac4c232SBarry Smith PetscTryMethod(ts, "TSMPRKSetType_C", (TS, TSMPRKType), (ts, mprktype)); 11364b84eec9SHong Zhang PetscFunctionReturn(0); 11374b84eec9SHong Zhang } 11384b84eec9SHong Zhang 11394b84eec9SHong Zhang /*@C 11404b84eec9SHong Zhang TSMPRKGetType - Get the type of MPRK scheme 11414b84eec9SHong Zhang 11420fe4d17eSHong Zhang Not collective 11434b84eec9SHong Zhang 11444b84eec9SHong Zhang Input Parameter: 11454b84eec9SHong Zhang . ts - timestepping context 11464b84eec9SHong Zhang 11474b84eec9SHong Zhang Output Parameter: 11484b84eec9SHong Zhang . mprktype - type of MPRK-scheme 11494b84eec9SHong Zhang 11504b84eec9SHong Zhang Level: intermediate 11514b84eec9SHong Zhang 1152db781477SPatrick Sanan .seealso: `TSMPRKGetType()` 11534b84eec9SHong Zhang @*/ 11549371c9d4SSatish Balay PetscErrorCode TSMPRKGetType(TS ts, TSMPRKType *mprktype) { 11554b84eec9SHong Zhang PetscFunctionBegin; 11564b84eec9SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1157cac4c232SBarry Smith PetscUseMethod(ts, "TSMPRKGetType_C", (TS, TSMPRKType *), (ts, mprktype)); 11584b84eec9SHong Zhang PetscFunctionReturn(0); 11594b84eec9SHong Zhang } 11604b84eec9SHong Zhang 11619371c9d4SSatish Balay static PetscErrorCode TSMPRKGetType_MPRK(TS ts, TSMPRKType *mprktype) { 11624b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK *)ts->data; 11634b84eec9SHong Zhang 11644b84eec9SHong Zhang PetscFunctionBegin; 11654b84eec9SHong Zhang *mprktype = mprk->tableau->name; 11664b84eec9SHong Zhang PetscFunctionReturn(0); 11674b84eec9SHong Zhang } 11684b84eec9SHong Zhang 11699371c9d4SSatish Balay static PetscErrorCode TSMPRKSetType_MPRK(TS ts, TSMPRKType mprktype) { 11704b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK *)ts->data; 11714b84eec9SHong Zhang PetscBool match; 11724b84eec9SHong Zhang MPRKTableauLink link; 11734b84eec9SHong Zhang 11744b84eec9SHong Zhang PetscFunctionBegin; 11754b84eec9SHong Zhang if (mprk->tableau) { 11769566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(mprk->tableau->name, mprktype, &match)); 11774b84eec9SHong Zhang if (match) PetscFunctionReturn(0); 11784b84eec9SHong Zhang } 11794b84eec9SHong Zhang for (link = MPRKTableauList; link; link = link->next) { 11809566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(link->tab.name, mprktype, &match)); 11814b84eec9SHong Zhang if (match) { 11829566063dSJacob Faibussowitsch if (ts->setupcalled) PetscCall(TSMPRKTableauReset(ts)); 11834b84eec9SHong Zhang mprk->tableau = &link->tab; 11849566063dSJacob Faibussowitsch if (ts->setupcalled) PetscCall(TSMPRKTableauSetUp(ts)); 11854b84eec9SHong Zhang PetscFunctionReturn(0); 11864b84eec9SHong Zhang } 11874b84eec9SHong Zhang } 118898921bdaSJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_UNKNOWN_TYPE, "Could not find '%s'", mprktype); 11894b84eec9SHong Zhang } 11904b84eec9SHong Zhang 11919371c9d4SSatish Balay static PetscErrorCode TSGetStages_MPRK(TS ts, PetscInt *ns, Vec **Y) { 11924b84eec9SHong Zhang TS_MPRK *mprk = (TS_MPRK *)ts->data; 11934b84eec9SHong Zhang 11944b84eec9SHong Zhang PetscFunctionBegin; 11954b84eec9SHong Zhang *ns = mprk->tableau->s; 11964b84eec9SHong Zhang if (Y) *Y = mprk->Y; 11974b84eec9SHong Zhang PetscFunctionReturn(0); 11984b84eec9SHong Zhang } 11994b84eec9SHong Zhang 12009371c9d4SSatish Balay static PetscErrorCode TSDestroy_MPRK(TS ts) { 12014b84eec9SHong Zhang PetscFunctionBegin; 12029566063dSJacob Faibussowitsch PetscCall(TSReset_MPRK(ts)); 12034b84eec9SHong Zhang if (ts->dm) { 12049566063dSJacob Faibussowitsch PetscCall(DMCoarsenHookRemove(ts->dm, DMCoarsenHook_TSMPRK, DMRestrictHook_TSMPRK, ts)); 12059566063dSJacob Faibussowitsch PetscCall(DMSubDomainHookRemove(ts->dm, DMSubDomainHook_TSMPRK, DMSubDomainRestrictHook_TSMPRK, ts)); 12064b84eec9SHong Zhang } 12079566063dSJacob Faibussowitsch PetscCall(PetscFree(ts->data)); 12089566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSMPRKGetType_C", NULL)); 12099566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSMPRKSetType_C", NULL)); 12104b84eec9SHong Zhang PetscFunctionReturn(0); 12114b84eec9SHong Zhang } 12124b84eec9SHong Zhang 12134b84eec9SHong Zhang /*MC 1214f944a40eSHong Zhang TSMPRK - ODE solver using Multirate Partitioned Runge-Kutta schemes 12154b84eec9SHong Zhang 12164b84eec9SHong Zhang The user should provide the right hand side of the equation 12174b84eec9SHong Zhang using TSSetRHSFunction(). 12184b84eec9SHong Zhang 12194b84eec9SHong Zhang Notes: 1220f944a40eSHong Zhang The default is TSMPRKPM2, it can be changed with TSMPRKSetType() or -ts_mprk_type 12214b84eec9SHong Zhang 12224b84eec9SHong Zhang Level: beginner 12234b84eec9SHong Zhang 1224db781477SPatrick Sanan .seealso: `TSCreate()`, `TS`, `TSSetType()`, `TSMPRKSetType()`, `TSMPRKGetType()`, `TSMPRKType`, `TSMPRKRegister()`, `TSMPRKSetMultirateType()` 1225db781477SPatrick Sanan `TSMPRKM2`, `TSMPRKM3`, `TSMPRKRFSMR3`, `TSMPRKRFSMR2` 12264b84eec9SHong Zhang 12274b84eec9SHong Zhang M*/ 12289371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode TSCreate_MPRK(TS ts) { 12294b84eec9SHong Zhang TS_MPRK *mprk; 12304b84eec9SHong Zhang 12314b84eec9SHong Zhang PetscFunctionBegin; 12329566063dSJacob Faibussowitsch PetscCall(TSMPRKInitializePackage()); 12334b84eec9SHong Zhang 12344b84eec9SHong Zhang ts->ops->reset = TSReset_MPRK; 12354b84eec9SHong Zhang ts->ops->destroy = TSDestroy_MPRK; 12364b84eec9SHong Zhang ts->ops->view = TSView_MPRK; 12374b84eec9SHong Zhang ts->ops->load = TSLoad_MPRK; 12384b84eec9SHong Zhang ts->ops->setup = TSSetUp_MPRK; 12394b84eec9SHong Zhang ts->ops->step = TSStep_MPRK; 12404b84eec9SHong Zhang ts->ops->evaluatestep = TSEvaluateStep_MPRK; 12414b84eec9SHong Zhang ts->ops->setfromoptions = TSSetFromOptions_MPRK; 12424b84eec9SHong Zhang ts->ops->getstages = TSGetStages_MPRK; 12434b84eec9SHong Zhang 12449566063dSJacob Faibussowitsch PetscCall(PetscNewLog(ts, &mprk)); 12454b84eec9SHong Zhang ts->data = (void *)mprk; 12464b84eec9SHong Zhang 12479566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSMPRKGetType_C", TSMPRKGetType_MPRK)); 12489566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSMPRKSetType_C", TSMPRKSetType_MPRK)); 12494b84eec9SHong Zhang 12509566063dSJacob Faibussowitsch PetscCall(TSMPRKSetType(ts, TSMPRKDefault)); 12514b84eec9SHong Zhang PetscFunctionReturn(0); 12524b84eec9SHong Zhang } 1253