1c4762a1bSJed Brown static char help[] = "Solves the ordinary differential equations (IVPs) using explicit and implicit time-integration methods.\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown /* 4c4762a1bSJed Brown 5c4762a1bSJed Brown Useful command line parameters: 6c4762a1bSJed Brown -problem <hull1972a1>: choose which problem to solve (see references 7c4762a1bSJed Brown for complete listing of problems). 8c4762a1bSJed Brown -ts_type <euler>: specify time-integrator 9c4762a1bSJed Brown -ts_adapt_type <basic>: specify time-step adapting (none,basic,advanced) 10c4762a1bSJed Brown -refinement_levels <1>: number of refinement levels for convergence analysis 11c4762a1bSJed Brown -refinement_factor <2.0>: factor to refine time step size by for convergence analysis 12c4762a1bSJed Brown -dt <0.01>: specify time step (initial time step for convergence analysis) 13c4762a1bSJed Brown 14c4762a1bSJed Brown */ 15c4762a1bSJed Brown 16c4762a1bSJed Brown /* 17c4762a1bSJed Brown List of cases and their names in the code:- 18c4762a1bSJed Brown From Hull, T.E., Enright, W.H., Fellen, B.M., and Sedgwick, A.E., 19c4762a1bSJed Brown "Comparing Numerical Methods for Ordinary Differential 20c4762a1bSJed Brown Equations", SIAM J. Numer. Anal., 9(4), 1972, pp. 603 - 635 21c4762a1bSJed Brown A1 -> "hull1972a1" (exact solution available) 22c4762a1bSJed Brown A2 -> "hull1972a2" (exact solution available) 23c4762a1bSJed Brown A3 -> "hull1972a3" (exact solution available) 24c4762a1bSJed Brown A4 -> "hull1972a4" (exact solution available) 25c4762a1bSJed Brown A5 -> "hull1972a5" 26c4762a1bSJed Brown B1 -> "hull1972b1" 27c4762a1bSJed Brown B2 -> "hull1972b2" 28c4762a1bSJed Brown B3 -> "hull1972b3" 29c4762a1bSJed Brown B4 -> "hull1972b4" 30c4762a1bSJed Brown B5 -> "hull1972b5" 31c4762a1bSJed Brown C1 -> "hull1972c1" 32c4762a1bSJed Brown C2 -> "hull1972c2" 33c4762a1bSJed Brown C3 -> "hull1972c3" 34c4762a1bSJed Brown C4 -> "hull1972c4" 35c4762a1bSJed Brown 36c4762a1bSJed Brown From Constantinescu, E. "Estimating Global Errors in Time Stepping" ArXiv e-prints, 37c4762a1bSJed Brown https://arxiv.org/abs/1503.05166, 2016 38c4762a1bSJed Brown 39c4762a1bSJed Brown Kulikov2013I -> "kulik2013i" 40c4762a1bSJed Brown 41c4762a1bSJed Brown */ 42c4762a1bSJed Brown 43c4762a1bSJed Brown #include <petscts.h> 44c4762a1bSJed Brown 45c4762a1bSJed Brown /* Function declarations */ 46c4762a1bSJed Brown PetscErrorCode (*RHSFunction)(TS, PetscReal, Vec, Vec, void *); 47c4762a1bSJed Brown PetscErrorCode (*RHSJacobian)(TS, PetscReal, Vec, Mat, Mat, void *); 48c4762a1bSJed Brown PetscErrorCode (*IFunction)(TS, PetscReal, Vec, Vec, Vec, void *); 49c4762a1bSJed Brown PetscErrorCode (*IJacobian)(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *); 50c4762a1bSJed Brown 51c4762a1bSJed Brown /* Returns the size of the system of equations depending on problem specification */ 52*11cc89d2SBarry Smith PetscErrorCode GetSize(const char *p, PetscInt *sz) { 537510d9b0SBarry Smith PetscFunctionBeginUser; 54*11cc89d2SBarry Smith 55*11cc89d2SBarry Smith if (!strcmp(p, "hull1972a1") || !strcmp(p, "hull1972a2") || !strcmp(p, "hull1972a3") || !strcmp(p, "hull1972a4") || !strcmp(p, "hull1972a5")) *sz = 1; 56*11cc89d2SBarry Smith else if (!strcmp(p, "hull1972b1")) *sz = 2; 57*11cc89d2SBarry Smith else if (!strcmp(p, "hull1972b2") || !strcmp(p, "hull1972b3") || !strcmp(p, "hull1972b4") || !strcmp(p, "hull1972b5")) *sz = 3; 58*11cc89d2SBarry Smith else if (!strcmp(p, "kulik2013i")) *sz = 4; 59*11cc89d2SBarry Smith else if (!strcmp(p, "hull1972c1") || !strcmp(p, "hull1972c2") || !strcmp(p, "hull1972c3")) *sz = 10; 60*11cc89d2SBarry Smith else if (!strcmp(p, "hull1972c4")) *sz = 52; 61*11cc89d2SBarry Smith else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Incorrect problem name %s", p); 62*11cc89d2SBarry Smith PetscFunctionReturn(0); 63c4762a1bSJed Brown } 64c4762a1bSJed Brown 65c4762a1bSJed Brown /****************************************************************/ 66c4762a1bSJed Brown 67c4762a1bSJed Brown /* Problem specific functions */ 68c4762a1bSJed Brown 69c4762a1bSJed Brown /* Hull, 1972, Problem A1 */ 70c4762a1bSJed Brown 719371c9d4SSatish Balay PetscErrorCode RHSFunction_Hull1972A1(TS ts, PetscReal t, Vec Y, Vec F, void *s) { 72c4762a1bSJed Brown PetscScalar *f; 73c4762a1bSJed Brown const PetscScalar *y; 74c4762a1bSJed Brown 757510d9b0SBarry Smith PetscFunctionBeginUser; 769566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 779566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 78c4762a1bSJed Brown f[0] = -y[0]; 799566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 809566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 81c4762a1bSJed Brown PetscFunctionReturn(0); 82c4762a1bSJed Brown } 83c4762a1bSJed Brown 849371c9d4SSatish Balay PetscErrorCode RHSJacobian_Hull1972A1(TS ts, PetscReal t, Vec Y, Mat A, Mat B, void *s) { 85c4762a1bSJed Brown const PetscScalar *y; 86c4762a1bSJed Brown PetscInt row = 0, col = 0; 87c4762a1bSJed Brown PetscScalar value = -1.0; 88c4762a1bSJed Brown 897510d9b0SBarry Smith PetscFunctionBeginUser; 909566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 919566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 1, &col, &value, INSERT_VALUES)); 929566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 939566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 949566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 95c4762a1bSJed Brown PetscFunctionReturn(0); 96c4762a1bSJed Brown } 97c4762a1bSJed Brown 989371c9d4SSatish Balay PetscErrorCode IFunction_Hull1972A1(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) { 99c4762a1bSJed Brown const PetscScalar *y; 100c4762a1bSJed Brown PetscScalar *f; 101c4762a1bSJed Brown 1027510d9b0SBarry Smith PetscFunctionBeginUser; 1039566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 1049566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 105c4762a1bSJed Brown f[0] = -y[0]; 1069566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 1079566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 108c4762a1bSJed Brown /* Left hand side = ydot - f(y) */ 1099566063dSJacob Faibussowitsch PetscCall(VecAYPX(F, -1.0, Ydot)); 110c4762a1bSJed Brown PetscFunctionReturn(0); 111c4762a1bSJed Brown } 112c4762a1bSJed Brown 1139371c9d4SSatish Balay PetscErrorCode IJacobian_Hull1972A1(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) { 114c4762a1bSJed Brown const PetscScalar *y; 115c4762a1bSJed Brown PetscInt row = 0, col = 0; 116c4762a1bSJed Brown PetscScalar value = a - 1.0; 117c4762a1bSJed Brown 1187510d9b0SBarry Smith PetscFunctionBeginUser; 1199566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 1209566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 1, &col, &value, INSERT_VALUES)); 1219566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1229566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 1239566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 124c4762a1bSJed Brown PetscFunctionReturn(0); 125c4762a1bSJed Brown } 126c4762a1bSJed Brown 127c4762a1bSJed Brown /* Hull, 1972, Problem A2 */ 128c4762a1bSJed Brown 1299371c9d4SSatish Balay PetscErrorCode RHSFunction_Hull1972A2(TS ts, PetscReal t, Vec Y, Vec F, void *s) { 130c4762a1bSJed Brown const PetscScalar *y; 131c4762a1bSJed Brown PetscScalar *f; 132c4762a1bSJed Brown 1337510d9b0SBarry Smith PetscFunctionBeginUser; 1349566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 1359566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 136c4762a1bSJed Brown f[0] = -0.5 * y[0] * y[0] * y[0]; 1379566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 1389566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 139c4762a1bSJed Brown PetscFunctionReturn(0); 140c4762a1bSJed Brown } 141c4762a1bSJed Brown 1429371c9d4SSatish Balay PetscErrorCode RHSJacobian_Hull1972A2(TS ts, PetscReal t, Vec Y, Mat A, Mat B, void *s) { 143c4762a1bSJed Brown const PetscScalar *y; 144c4762a1bSJed Brown PetscInt row = 0, col = 0; 145c4762a1bSJed Brown PetscScalar value; 146c4762a1bSJed Brown 1477510d9b0SBarry Smith PetscFunctionBeginUser; 1489566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 149c4762a1bSJed Brown value = -0.5 * 3.0 * y[0] * y[0]; 1509566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 1, &col, &value, INSERT_VALUES)); 1519566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1529566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 1539566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 154c4762a1bSJed Brown PetscFunctionReturn(0); 155c4762a1bSJed Brown } 156c4762a1bSJed Brown 1579371c9d4SSatish Balay PetscErrorCode IFunction_Hull1972A2(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) { 158c4762a1bSJed Brown PetscScalar *f; 159c4762a1bSJed Brown const PetscScalar *y; 160c4762a1bSJed Brown 1617510d9b0SBarry Smith PetscFunctionBeginUser; 1629566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 1639566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 164c4762a1bSJed Brown f[0] = -0.5 * y[0] * y[0] * y[0]; 1659566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 1669566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 167c4762a1bSJed Brown /* Left hand side = ydot - f(y) */ 1689566063dSJacob Faibussowitsch PetscCall(VecAYPX(F, -1.0, Ydot)); 169c4762a1bSJed Brown PetscFunctionReturn(0); 170c4762a1bSJed Brown } 171c4762a1bSJed Brown 1729371c9d4SSatish Balay PetscErrorCode IJacobian_Hull1972A2(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) { 173c4762a1bSJed Brown const PetscScalar *y; 174c4762a1bSJed Brown PetscInt row = 0, col = 0; 175c4762a1bSJed Brown PetscScalar value; 176c4762a1bSJed Brown 1777510d9b0SBarry Smith PetscFunctionBeginUser; 1789566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 179c4762a1bSJed Brown value = a + 0.5 * 3.0 * y[0] * y[0]; 1809566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 1, &col, &value, INSERT_VALUES)); 1819566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1829566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 1839566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 184c4762a1bSJed Brown PetscFunctionReturn(0); 185c4762a1bSJed Brown } 186c4762a1bSJed Brown 187c4762a1bSJed Brown /* Hull, 1972, Problem A3 */ 188c4762a1bSJed Brown 1899371c9d4SSatish Balay PetscErrorCode RHSFunction_Hull1972A3(TS ts, PetscReal t, Vec Y, Vec F, void *s) { 190c4762a1bSJed Brown const PetscScalar *y; 191c4762a1bSJed Brown PetscScalar *f; 192c4762a1bSJed Brown 1937510d9b0SBarry Smith PetscFunctionBeginUser; 1949566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 1959566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 196c4762a1bSJed Brown f[0] = y[0] * PetscCosReal(t); 1979566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 1989566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 199c4762a1bSJed Brown PetscFunctionReturn(0); 200c4762a1bSJed Brown } 201c4762a1bSJed Brown 2029371c9d4SSatish Balay PetscErrorCode RHSJacobian_Hull1972A3(TS ts, PetscReal t, Vec Y, Mat A, Mat B, void *s) { 203c4762a1bSJed Brown const PetscScalar *y; 204c4762a1bSJed Brown PetscInt row = 0, col = 0; 205c4762a1bSJed Brown PetscScalar value = PetscCosReal(t); 206c4762a1bSJed Brown 2077510d9b0SBarry Smith PetscFunctionBeginUser; 2089566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 2099566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 1, &col, &value, INSERT_VALUES)); 2109566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 2119566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 2129566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 213c4762a1bSJed Brown PetscFunctionReturn(0); 214c4762a1bSJed Brown } 215c4762a1bSJed Brown 2169371c9d4SSatish Balay PetscErrorCode IFunction_Hull1972A3(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) { 217c4762a1bSJed Brown PetscScalar *f; 218c4762a1bSJed Brown const PetscScalar *y; 219c4762a1bSJed Brown 2207510d9b0SBarry Smith PetscFunctionBeginUser; 2219566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 2229566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 223c4762a1bSJed Brown f[0] = y[0] * PetscCosReal(t); 2249566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 2259566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 226c4762a1bSJed Brown /* Left hand side = ydot - f(y) */ 2279566063dSJacob Faibussowitsch PetscCall(VecAYPX(F, -1.0, Ydot)); 228c4762a1bSJed Brown PetscFunctionReturn(0); 229c4762a1bSJed Brown } 230c4762a1bSJed Brown 2319371c9d4SSatish Balay PetscErrorCode IJacobian_Hull1972A3(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) { 232c4762a1bSJed Brown const PetscScalar *y; 233c4762a1bSJed Brown PetscInt row = 0, col = 0; 234c4762a1bSJed Brown PetscScalar value = a - PetscCosReal(t); 235c4762a1bSJed Brown 2367510d9b0SBarry Smith PetscFunctionBeginUser; 2379566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 2389566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 1, &col, &value, INSERT_VALUES)); 2399566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 2409566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 2419566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 242c4762a1bSJed Brown PetscFunctionReturn(0); 243c4762a1bSJed Brown } 244c4762a1bSJed Brown 245c4762a1bSJed Brown /* Hull, 1972, Problem A4 */ 246c4762a1bSJed Brown 2479371c9d4SSatish Balay PetscErrorCode RHSFunction_Hull1972A4(TS ts, PetscReal t, Vec Y, Vec F, void *s) { 248c4762a1bSJed Brown PetscScalar *f; 249c4762a1bSJed Brown const PetscScalar *y; 250c4762a1bSJed Brown 2517510d9b0SBarry Smith PetscFunctionBeginUser; 2529566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 2539566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 254c4762a1bSJed Brown f[0] = (0.25 * y[0]) * (1.0 - 0.05 * y[0]); 2559566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 2569566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 257c4762a1bSJed Brown PetscFunctionReturn(0); 258c4762a1bSJed Brown } 259c4762a1bSJed Brown 2609371c9d4SSatish Balay PetscErrorCode RHSJacobian_Hull1972A4(TS ts, PetscReal t, Vec Y, Mat A, Mat B, void *s) { 261c4762a1bSJed Brown const PetscScalar *y; 262c4762a1bSJed Brown PetscInt row = 0, col = 0; 263c4762a1bSJed Brown PetscScalar value; 264c4762a1bSJed Brown 2657510d9b0SBarry Smith PetscFunctionBeginUser; 2669566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 267c4762a1bSJed Brown value = 0.25 * (1.0 - 0.05 * y[0]) - (0.25 * y[0]) * 0.05; 2689566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 1, &col, &value, INSERT_VALUES)); 2699566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 2709566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 2719566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 272c4762a1bSJed Brown PetscFunctionReturn(0); 273c4762a1bSJed Brown } 274c4762a1bSJed Brown 2759371c9d4SSatish Balay PetscErrorCode IFunction_Hull1972A4(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) { 276c4762a1bSJed Brown PetscScalar *f; 277c4762a1bSJed Brown const PetscScalar *y; 278c4762a1bSJed Brown 2797510d9b0SBarry Smith PetscFunctionBeginUser; 2809566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 2819566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 282c4762a1bSJed Brown f[0] = (0.25 * y[0]) * (1.0 - 0.05 * y[0]); 2839566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 2849566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 285c4762a1bSJed Brown /* Left hand side = ydot - f(y) */ 2869566063dSJacob Faibussowitsch PetscCall(VecAYPX(F, -1.0, Ydot)); 287c4762a1bSJed Brown PetscFunctionReturn(0); 288c4762a1bSJed Brown } 289c4762a1bSJed Brown 2909371c9d4SSatish Balay PetscErrorCode IJacobian_Hull1972A4(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) { 291c4762a1bSJed Brown const PetscScalar *y; 292c4762a1bSJed Brown PetscInt row = 0, col = 0; 293c4762a1bSJed Brown PetscScalar value; 294c4762a1bSJed Brown 2957510d9b0SBarry Smith PetscFunctionBeginUser; 2969566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 297c4762a1bSJed Brown value = a - 0.25 * (1.0 - 0.05 * y[0]) + (0.25 * y[0]) * 0.05; 2989566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 1, &col, &value, INSERT_VALUES)); 2999566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 3009566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 3019566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 302c4762a1bSJed Brown PetscFunctionReturn(0); 303c4762a1bSJed Brown } 304c4762a1bSJed Brown 305c4762a1bSJed Brown /* Hull, 1972, Problem A5 */ 306c4762a1bSJed Brown 3079371c9d4SSatish Balay PetscErrorCode RHSFunction_Hull1972A5(TS ts, PetscReal t, Vec Y, Vec F, void *s) { 308c4762a1bSJed Brown PetscScalar *f; 309c4762a1bSJed Brown const PetscScalar *y; 310c4762a1bSJed Brown 3117510d9b0SBarry Smith PetscFunctionBeginUser; 3129566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 3139566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 314c4762a1bSJed Brown f[0] = (y[0] - t) / (y[0] + t); 3159566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 3169566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 317c4762a1bSJed Brown PetscFunctionReturn(0); 318c4762a1bSJed Brown } 319c4762a1bSJed Brown 3209371c9d4SSatish Balay PetscErrorCode RHSJacobian_Hull1972A5(TS ts, PetscReal t, Vec Y, Mat A, Mat B, void *s) { 321c4762a1bSJed Brown const PetscScalar *y; 322c4762a1bSJed Brown PetscInt row = 0, col = 0; 323c4762a1bSJed Brown PetscScalar value; 324c4762a1bSJed Brown 3257510d9b0SBarry Smith PetscFunctionBeginUser; 3269566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 327c4762a1bSJed Brown value = 2 * t / ((t + y[0]) * (t + y[0])); 3289566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 1, &col, &value, INSERT_VALUES)); 3299566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 3309566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 3319566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 332c4762a1bSJed Brown PetscFunctionReturn(0); 333c4762a1bSJed Brown } 334c4762a1bSJed Brown 3359371c9d4SSatish Balay PetscErrorCode IFunction_Hull1972A5(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) { 336c4762a1bSJed Brown PetscScalar *f; 337c4762a1bSJed Brown const PetscScalar *y; 338c4762a1bSJed Brown 3397510d9b0SBarry Smith PetscFunctionBeginUser; 3409566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 3419566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 342c4762a1bSJed Brown f[0] = (y[0] - t) / (y[0] + t); 3439566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 3449566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 345c4762a1bSJed Brown /* Left hand side = ydot - f(y) */ 3469566063dSJacob Faibussowitsch PetscCall(VecAYPX(F, -1.0, Ydot)); 347c4762a1bSJed Brown PetscFunctionReturn(0); 348c4762a1bSJed Brown } 349c4762a1bSJed Brown 3509371c9d4SSatish Balay PetscErrorCode IJacobian_Hull1972A5(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) { 351c4762a1bSJed Brown const PetscScalar *y; 352c4762a1bSJed Brown PetscInt row = 0, col = 0; 353c4762a1bSJed Brown PetscScalar value; 354c4762a1bSJed Brown 3557510d9b0SBarry Smith PetscFunctionBeginUser; 3569566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 357c4762a1bSJed Brown value = a - 2 * t / ((t + y[0]) * (t + y[0])); 3589566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 1, &col, &value, INSERT_VALUES)); 3599566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 3609566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 3619566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 362c4762a1bSJed Brown PetscFunctionReturn(0); 363c4762a1bSJed Brown } 364c4762a1bSJed Brown 365c4762a1bSJed Brown /* Hull, 1972, Problem B1 */ 366c4762a1bSJed Brown 3679371c9d4SSatish Balay PetscErrorCode RHSFunction_Hull1972B1(TS ts, PetscReal t, Vec Y, Vec F, void *s) { 368c4762a1bSJed Brown PetscScalar *f; 369c4762a1bSJed Brown const PetscScalar *y; 370c4762a1bSJed Brown 3717510d9b0SBarry Smith PetscFunctionBeginUser; 3729566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 3739566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 374c4762a1bSJed Brown f[0] = 2.0 * (y[0] - y[0] * y[1]); 375c4762a1bSJed Brown f[1] = -(y[1] - y[0] * y[1]); 3769566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 3779566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 378c4762a1bSJed Brown PetscFunctionReturn(0); 379c4762a1bSJed Brown } 380c4762a1bSJed Brown 3819371c9d4SSatish Balay PetscErrorCode IFunction_Hull1972B1(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) { 382c4762a1bSJed Brown PetscScalar *f; 383c4762a1bSJed Brown const PetscScalar *y; 384c4762a1bSJed Brown 3857510d9b0SBarry Smith PetscFunctionBeginUser; 3869566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 3879566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 388c4762a1bSJed Brown f[0] = 2.0 * (y[0] - y[0] * y[1]); 389c4762a1bSJed Brown f[1] = -(y[1] - y[0] * y[1]); 3909566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 3919566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 392c4762a1bSJed Brown /* Left hand side = ydot - f(y) */ 3939566063dSJacob Faibussowitsch PetscCall(VecAYPX(F, -1.0, Ydot)); 394c4762a1bSJed Brown PetscFunctionReturn(0); 395c4762a1bSJed Brown } 396c4762a1bSJed Brown 3979371c9d4SSatish Balay PetscErrorCode IJacobian_Hull1972B1(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) { 398c4762a1bSJed Brown const PetscScalar *y; 399c4762a1bSJed Brown PetscInt row[2] = {0, 1}; 400c4762a1bSJed Brown PetscScalar value[2][2]; 401c4762a1bSJed Brown 4027510d9b0SBarry Smith PetscFunctionBeginUser; 4039566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 4049371c9d4SSatish Balay value[0][0] = a - 2.0 * (1.0 - y[1]); 4059371c9d4SSatish Balay value[0][1] = 2.0 * y[0]; 4069371c9d4SSatish Balay value[1][0] = -y[1]; 4079371c9d4SSatish Balay value[1][1] = a + 1.0 - y[0]; 4089566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 2, &row[0], 2, &row[0], &value[0][0], INSERT_VALUES)); 4099566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 4109566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 4119566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 412c4762a1bSJed Brown PetscFunctionReturn(0); 413c4762a1bSJed Brown } 414c4762a1bSJed Brown 415c4762a1bSJed Brown /* Hull, 1972, Problem B2 */ 416c4762a1bSJed Brown 4179371c9d4SSatish Balay PetscErrorCode RHSFunction_Hull1972B2(TS ts, PetscReal t, Vec Y, Vec F, void *s) { 418c4762a1bSJed Brown PetscScalar *f; 419c4762a1bSJed Brown const PetscScalar *y; 420c4762a1bSJed Brown 4217510d9b0SBarry Smith PetscFunctionBeginUser; 4229566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 4239566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 424c4762a1bSJed Brown f[0] = -y[0] + y[1]; 425c4762a1bSJed Brown f[1] = y[0] - 2.0 * y[1] + y[2]; 426c4762a1bSJed Brown f[2] = y[1] - y[2]; 4279566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 4289566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 429c4762a1bSJed Brown PetscFunctionReturn(0); 430c4762a1bSJed Brown } 431c4762a1bSJed Brown 4329371c9d4SSatish Balay PetscErrorCode IFunction_Hull1972B2(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) { 433c4762a1bSJed Brown PetscScalar *f; 434c4762a1bSJed Brown const PetscScalar *y; 435c4762a1bSJed Brown 4367510d9b0SBarry Smith PetscFunctionBeginUser; 4379566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 4389566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 439c4762a1bSJed Brown f[0] = -y[0] + y[1]; 440c4762a1bSJed Brown f[1] = y[0] - 2.0 * y[1] + y[2]; 441c4762a1bSJed Brown f[2] = y[1] - y[2]; 4429566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 4439566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 444c4762a1bSJed Brown /* Left hand side = ydot - f(y) */ 4459566063dSJacob Faibussowitsch PetscCall(VecAYPX(F, -1.0, Ydot)); 446c4762a1bSJed Brown PetscFunctionReturn(0); 447c4762a1bSJed Brown } 448c4762a1bSJed Brown 4499371c9d4SSatish Balay PetscErrorCode IJacobian_Hull1972B2(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) { 450c4762a1bSJed Brown const PetscScalar *y; 451c4762a1bSJed Brown PetscInt row[3] = {0, 1, 2}; 452c4762a1bSJed Brown PetscScalar value[3][3]; 453c4762a1bSJed Brown 4547510d9b0SBarry Smith PetscFunctionBeginUser; 4559566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 4569371c9d4SSatish Balay value[0][0] = a + 1.0; 4579371c9d4SSatish Balay value[0][1] = -1.0; 4589371c9d4SSatish Balay value[0][2] = 0; 4599371c9d4SSatish Balay value[1][0] = -1.0; 4609371c9d4SSatish Balay value[1][1] = a + 2.0; 4619371c9d4SSatish Balay value[1][2] = -1.0; 4629371c9d4SSatish Balay value[2][0] = 0; 4639371c9d4SSatish Balay value[2][1] = -1.0; 4649371c9d4SSatish Balay value[2][2] = a + 1.0; 4659566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 3, &row[0], 3, &row[0], &value[0][0], INSERT_VALUES)); 4669566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 4679566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 4689566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 469c4762a1bSJed Brown PetscFunctionReturn(0); 470c4762a1bSJed Brown } 471c4762a1bSJed Brown 472c4762a1bSJed Brown /* Hull, 1972, Problem B3 */ 473c4762a1bSJed Brown 4749371c9d4SSatish Balay PetscErrorCode RHSFunction_Hull1972B3(TS ts, PetscReal t, Vec Y, Vec F, void *s) { 475c4762a1bSJed Brown PetscScalar *f; 476c4762a1bSJed Brown const PetscScalar *y; 477c4762a1bSJed Brown 4787510d9b0SBarry Smith PetscFunctionBeginUser; 4799566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 4809566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 481c4762a1bSJed Brown f[0] = -y[0]; 482c4762a1bSJed Brown f[1] = y[0] - y[1] * y[1]; 483c4762a1bSJed Brown f[2] = y[1] * y[1]; 4849566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 4859566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 486c4762a1bSJed Brown PetscFunctionReturn(0); 487c4762a1bSJed Brown } 488c4762a1bSJed Brown 4899371c9d4SSatish Balay PetscErrorCode IFunction_Hull1972B3(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) { 490c4762a1bSJed Brown PetscScalar *f; 491c4762a1bSJed Brown const PetscScalar *y; 492c4762a1bSJed Brown 4937510d9b0SBarry Smith PetscFunctionBeginUser; 4949566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 4959566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 496c4762a1bSJed Brown f[0] = -y[0]; 497c4762a1bSJed Brown f[1] = y[0] - y[1] * y[1]; 498c4762a1bSJed Brown f[2] = y[1] * y[1]; 4999566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 5009566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 501c4762a1bSJed Brown /* Left hand side = ydot - f(y) */ 5029566063dSJacob Faibussowitsch PetscCall(VecAYPX(F, -1.0, Ydot)); 503c4762a1bSJed Brown PetscFunctionReturn(0); 504c4762a1bSJed Brown } 505c4762a1bSJed Brown 5069371c9d4SSatish Balay PetscErrorCode IJacobian_Hull1972B3(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) { 507c4762a1bSJed Brown const PetscScalar *y; 508c4762a1bSJed Brown PetscInt row[3] = {0, 1, 2}; 509c4762a1bSJed Brown PetscScalar value[3][3]; 510c4762a1bSJed Brown 5117510d9b0SBarry Smith PetscFunctionBeginUser; 5129566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 5139371c9d4SSatish Balay value[0][0] = a + 1.0; 5149371c9d4SSatish Balay value[0][1] = 0; 5159371c9d4SSatish Balay value[0][2] = 0; 5169371c9d4SSatish Balay value[1][0] = -1.0; 5179371c9d4SSatish Balay value[1][1] = a + 2.0 * y[1]; 5189371c9d4SSatish Balay value[1][2] = 0; 5199371c9d4SSatish Balay value[2][0] = 0; 5209371c9d4SSatish Balay value[2][1] = -2.0 * y[1]; 5219371c9d4SSatish Balay value[2][2] = a; 5229566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 3, &row[0], 3, &row[0], &value[0][0], INSERT_VALUES)); 5239566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 5249566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 5259566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 526c4762a1bSJed Brown PetscFunctionReturn(0); 527c4762a1bSJed Brown } 528c4762a1bSJed Brown 529c4762a1bSJed Brown /* Hull, 1972, Problem B4 */ 530c4762a1bSJed Brown 5319371c9d4SSatish Balay PetscErrorCode RHSFunction_Hull1972B4(TS ts, PetscReal t, Vec Y, Vec F, void *s) { 532c4762a1bSJed Brown PetscScalar *f; 533c4762a1bSJed Brown const PetscScalar *y; 534c4762a1bSJed Brown 5357510d9b0SBarry Smith PetscFunctionBeginUser; 5369566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 5379566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 538c4762a1bSJed Brown f[0] = -y[1] - y[0] * y[2] / PetscSqrtScalar(y[0] * y[0] + y[1] * y[1]); 539c4762a1bSJed Brown f[1] = y[0] - y[1] * y[2] / PetscSqrtScalar(y[0] * y[0] + y[1] * y[1]); 540c4762a1bSJed Brown f[2] = y[0] / PetscSqrtScalar(y[0] * y[0] + y[1] * y[1]); 5419566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 5429566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 543c4762a1bSJed Brown PetscFunctionReturn(0); 544c4762a1bSJed Brown } 545c4762a1bSJed Brown 5469371c9d4SSatish Balay PetscErrorCode IFunction_Hull1972B4(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) { 547c4762a1bSJed Brown PetscScalar *f; 548c4762a1bSJed Brown const PetscScalar *y; 549c4762a1bSJed Brown 5507510d9b0SBarry Smith PetscFunctionBeginUser; 5519566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 5529566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 553c4762a1bSJed Brown f[0] = -y[1] - y[0] * y[2] / PetscSqrtScalar(y[0] * y[0] + y[1] * y[1]); 554c4762a1bSJed Brown f[1] = y[0] - y[1] * y[2] / PetscSqrtScalar(y[0] * y[0] + y[1] * y[1]); 555c4762a1bSJed Brown f[2] = y[0] / PetscSqrtScalar(y[0] * y[0] + y[1] * y[1]); 5569566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 5579566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 558c4762a1bSJed Brown /* Left hand side = ydot - f(y) */ 5599566063dSJacob Faibussowitsch PetscCall(VecAYPX(F, -1.0, Ydot)); 560c4762a1bSJed Brown PetscFunctionReturn(0); 561c4762a1bSJed Brown } 562c4762a1bSJed Brown 5639371c9d4SSatish Balay PetscErrorCode IJacobian_Hull1972B4(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) { 564c4762a1bSJed Brown const PetscScalar *y; 565c4762a1bSJed Brown PetscInt row[3] = {0, 1, 2}; 566c4762a1bSJed Brown PetscScalar value[3][3], fac, fac2; 567c4762a1bSJed Brown 5687510d9b0SBarry Smith PetscFunctionBeginUser; 5699566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 570c4762a1bSJed Brown fac = PetscPowScalar(y[0] * y[0] + y[1] * y[1], -1.5); 571c4762a1bSJed Brown fac2 = PetscPowScalar(y[0] * y[0] + y[1] * y[1], -0.5); 572c4762a1bSJed Brown value[0][0] = a + (y[1] * y[1] * y[2]) * fac; 573c4762a1bSJed Brown value[0][1] = 1.0 - (y[0] * y[1] * y[2]) * fac; 574c4762a1bSJed Brown value[0][2] = y[0] * fac2; 575c4762a1bSJed Brown value[1][0] = -1.0 - y[0] * y[1] * y[2] * fac; 576c4762a1bSJed Brown value[1][1] = a + y[0] * y[0] * y[2] * fac; 577c4762a1bSJed Brown value[1][2] = y[1] * fac2; 578c4762a1bSJed Brown value[2][0] = -y[1] * y[1] * fac; 579c4762a1bSJed Brown value[2][1] = y[0] * y[1] * fac; 580c4762a1bSJed Brown value[2][2] = a; 5819566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 3, &row[0], 3, &row[0], &value[0][0], INSERT_VALUES)); 5829566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 5839566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 5849566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 585c4762a1bSJed Brown PetscFunctionReturn(0); 586c4762a1bSJed Brown } 587c4762a1bSJed Brown 588c4762a1bSJed Brown /* Hull, 1972, Problem B5 */ 589c4762a1bSJed Brown 5909371c9d4SSatish Balay PetscErrorCode RHSFunction_Hull1972B5(TS ts, PetscReal t, Vec Y, Vec F, void *s) { 591c4762a1bSJed Brown PetscScalar *f; 592c4762a1bSJed Brown const PetscScalar *y; 593c4762a1bSJed Brown 5947510d9b0SBarry Smith PetscFunctionBeginUser; 5959566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 5969566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 597c4762a1bSJed Brown f[0] = y[1] * y[2]; 598c4762a1bSJed Brown f[1] = -y[0] * y[2]; 599c4762a1bSJed Brown f[2] = -0.51 * y[0] * y[1]; 6009566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 6019566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 602c4762a1bSJed Brown PetscFunctionReturn(0); 603c4762a1bSJed Brown } 604c4762a1bSJed Brown 6059371c9d4SSatish Balay PetscErrorCode IFunction_Hull1972B5(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) { 606c4762a1bSJed Brown PetscScalar *f; 607c4762a1bSJed Brown const PetscScalar *y; 608c4762a1bSJed Brown 6097510d9b0SBarry Smith PetscFunctionBeginUser; 6109566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 6119566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 612c4762a1bSJed Brown f[0] = y[1] * y[2]; 613c4762a1bSJed Brown f[1] = -y[0] * y[2]; 614c4762a1bSJed Brown f[2] = -0.51 * y[0] * y[1]; 6159566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 6169566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 617c4762a1bSJed Brown /* Left hand side = ydot - f(y) */ 6189566063dSJacob Faibussowitsch PetscCall(VecAYPX(F, -1.0, Ydot)); 619c4762a1bSJed Brown PetscFunctionReturn(0); 620c4762a1bSJed Brown } 621c4762a1bSJed Brown 6229371c9d4SSatish Balay PetscErrorCode IJacobian_Hull1972B5(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) { 623c4762a1bSJed Brown const PetscScalar *y; 624c4762a1bSJed Brown PetscInt row[3] = {0, 1, 2}; 625c4762a1bSJed Brown PetscScalar value[3][3]; 626c4762a1bSJed Brown 6277510d9b0SBarry Smith PetscFunctionBeginUser; 6289566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 6299371c9d4SSatish Balay value[0][0] = a; 6309371c9d4SSatish Balay value[0][1] = -y[2]; 6319371c9d4SSatish Balay value[0][2] = -y[1]; 6329371c9d4SSatish Balay value[1][0] = y[2]; 6339371c9d4SSatish Balay value[1][1] = a; 6349371c9d4SSatish Balay value[1][2] = y[0]; 6359371c9d4SSatish Balay value[2][0] = 0.51 * y[1]; 6369371c9d4SSatish Balay value[2][1] = 0.51 * y[0]; 6379371c9d4SSatish Balay value[2][2] = a; 6389566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 3, &row[0], 3, &row[0], &value[0][0], INSERT_VALUES)); 6399566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 6409566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 6419566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 642c4762a1bSJed Brown PetscFunctionReturn(0); 643c4762a1bSJed Brown } 644c4762a1bSJed Brown 645c4762a1bSJed Brown /* Kulikov, 2013, Problem I */ 646c4762a1bSJed Brown 6479371c9d4SSatish Balay PetscErrorCode RHSFunction_Kulikov2013I(TS ts, PetscReal t, Vec Y, Vec F, void *s) { 648c4762a1bSJed Brown PetscScalar *f; 649c4762a1bSJed Brown const PetscScalar *y; 650c4762a1bSJed Brown 6517510d9b0SBarry Smith PetscFunctionBeginUser; 6529566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 6539566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 654c4762a1bSJed Brown f[0] = 2. * t * PetscPowScalar(y[1], 1. / 5.) * y[3]; 655c4762a1bSJed Brown f[1] = 10. * t * y[3] * PetscExpScalar(5.0 * (y[2] - 1.)); 656c4762a1bSJed Brown f[2] = 2. * t * y[3]; 657c4762a1bSJed Brown f[3] = -2. * t * PetscLogScalar(y[0]); 6589566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 6599566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 660c4762a1bSJed Brown PetscFunctionReturn(0); 661c4762a1bSJed Brown } 662c4762a1bSJed Brown 6639371c9d4SSatish Balay PetscErrorCode RHSJacobian_Kulikov2013I(TS ts, PetscReal t, Vec Y, Mat A, Mat B, void *s) { 664c4762a1bSJed Brown const PetscScalar *y; 665c4762a1bSJed Brown PetscInt row[4] = {0, 1, 2, 3}; 666c4762a1bSJed Brown PetscScalar value[4][4]; 667c4762a1bSJed Brown PetscScalar m1, m2; 6687510d9b0SBarry Smith PetscFunctionBeginUser; 6699566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 670c4762a1bSJed Brown m1 = (2. * t * y[3]) / (5. * PetscPowScalar(y[1], 4. / 5.)); 671c4762a1bSJed Brown m2 = 2. * t * PetscPowScalar(y[1], 1. / 5.); 6729371c9d4SSatish Balay value[0][0] = 0.; 6739371c9d4SSatish Balay value[0][1] = m1; 6749371c9d4SSatish Balay value[0][2] = 0.; 6759371c9d4SSatish Balay value[0][3] = m2; 676c4762a1bSJed Brown m1 = 50. * t * y[3] * PetscExpScalar(5.0 * (y[2] - 1.)); 677c4762a1bSJed Brown m2 = 10. * t * PetscExpScalar(5.0 * (y[2] - 1.)); 6789371c9d4SSatish Balay value[1][0] = 0.; 6799371c9d4SSatish Balay value[1][1] = 0.; 6809371c9d4SSatish Balay value[1][2] = m1; 6819371c9d4SSatish Balay value[1][3] = m2; 6829371c9d4SSatish Balay value[2][0] = 0.; 6839371c9d4SSatish Balay value[2][1] = 0.; 6849371c9d4SSatish Balay value[2][2] = 0.; 6859371c9d4SSatish Balay value[2][3] = 2 * t; 6869371c9d4SSatish Balay value[3][0] = -2. * t / y[0]; 6879371c9d4SSatish Balay value[3][1] = 0.; 6889371c9d4SSatish Balay value[3][2] = 0.; 6899371c9d4SSatish Balay value[3][3] = 0.; 6909566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 4, &row[0], 4, &row[0], &value[0][0], INSERT_VALUES)); 6919566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 6929566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 6939566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 694c4762a1bSJed Brown PetscFunctionReturn(0); 695c4762a1bSJed Brown } 696c4762a1bSJed Brown 6979371c9d4SSatish Balay PetscErrorCode IFunction_Kulikov2013I(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) { 698c4762a1bSJed Brown PetscScalar *f; 699c4762a1bSJed Brown const PetscScalar *y; 700c4762a1bSJed Brown 7017510d9b0SBarry Smith PetscFunctionBeginUser; 7029566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 7039566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 704c4762a1bSJed Brown f[0] = 2. * t * PetscPowScalar(y[1], 1. / 5.) * y[3]; 705c4762a1bSJed Brown f[1] = 10. * t * y[3] * PetscExpScalar(5.0 * (y[2] - 1.)); 706c4762a1bSJed Brown f[2] = 2. * t * y[3]; 707c4762a1bSJed Brown f[3] = -2. * t * PetscLogScalar(y[0]); 7089566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 7099566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 710c4762a1bSJed Brown /* Left hand side = ydot - f(y) */ 7119566063dSJacob Faibussowitsch PetscCall(VecAYPX(F, -1.0, Ydot)); 712c4762a1bSJed Brown PetscFunctionReturn(0); 713c4762a1bSJed Brown } 714c4762a1bSJed Brown 7159371c9d4SSatish Balay PetscErrorCode IJacobian_Kulikov2013I(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) { 716c4762a1bSJed Brown const PetscScalar *y; 717c4762a1bSJed Brown PetscInt row[4] = {0, 1, 2, 3}; 718c4762a1bSJed Brown PetscScalar value[4][4]; 719c4762a1bSJed Brown PetscScalar m1, m2; 720c4762a1bSJed Brown 7217510d9b0SBarry Smith PetscFunctionBeginUser; 7229566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 723c4762a1bSJed Brown m1 = (2. * t * y[3]) / (5. * PetscPowScalar(y[1], 4. / 5.)); 724c4762a1bSJed Brown m2 = 2. * t * PetscPowScalar(y[1], 1. / 5.); 7259371c9d4SSatish Balay value[0][0] = a; 7269371c9d4SSatish Balay value[0][1] = m1; 7279371c9d4SSatish Balay value[0][2] = 0.; 7289371c9d4SSatish Balay value[0][3] = m2; 729c4762a1bSJed Brown m1 = 50. * t * y[3] * PetscExpScalar(5.0 * (y[2] - 1.)); 730c4762a1bSJed Brown m2 = 10. * t * PetscExpScalar(5.0 * (y[2] - 1.)); 7319371c9d4SSatish Balay value[1][0] = 0.; 7329371c9d4SSatish Balay value[1][1] = a; 7339371c9d4SSatish Balay value[1][2] = m1; 7349371c9d4SSatish Balay value[1][3] = m2; 7359371c9d4SSatish Balay value[2][0] = 0.; 7369371c9d4SSatish Balay value[2][1] = 0.; 7379371c9d4SSatish Balay value[2][2] = a; 7389371c9d4SSatish Balay value[2][3] = 2 * t; 7399371c9d4SSatish Balay value[3][0] = -2. * t / y[0]; 7409371c9d4SSatish Balay value[3][1] = 0.; 7419371c9d4SSatish Balay value[3][2] = 0.; 7429371c9d4SSatish Balay value[3][3] = a; 7439566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 4, &row[0], 4, &row[0], &value[0][0], INSERT_VALUES)); 7449566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 7459566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 7469566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 747c4762a1bSJed Brown PetscFunctionReturn(0); 748c4762a1bSJed Brown } 749c4762a1bSJed Brown 750c4762a1bSJed Brown /* Hull, 1972, Problem C1 */ 751c4762a1bSJed Brown 7529371c9d4SSatish Balay PetscErrorCode RHSFunction_Hull1972C1(TS ts, PetscReal t, Vec Y, Vec F, void *s) { 753c4762a1bSJed Brown PetscScalar *f; 754c4762a1bSJed Brown const PetscScalar *y; 755c4762a1bSJed Brown PetscInt N, i; 756c4762a1bSJed Brown 7577510d9b0SBarry Smith PetscFunctionBeginUser; 7589566063dSJacob Faibussowitsch PetscCall(VecGetSize(Y, &N)); 7599566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 7609566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 761c4762a1bSJed Brown f[0] = -y[0]; 762ad540459SPierre Jolivet for (i = 1; i < N - 1; i++) f[i] = y[i - 1] - y[i]; 763c4762a1bSJed Brown f[N - 1] = y[N - 2]; 7649566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 7659566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 766c4762a1bSJed Brown PetscFunctionReturn(0); 767c4762a1bSJed Brown } 768c4762a1bSJed Brown 7699371c9d4SSatish Balay PetscErrorCode IFunction_Hull1972C1(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) { 770c4762a1bSJed Brown PetscScalar *f; 771c4762a1bSJed Brown const PetscScalar *y; 772c4762a1bSJed Brown PetscInt N, i; 773c4762a1bSJed Brown 7747510d9b0SBarry Smith PetscFunctionBeginUser; 7759566063dSJacob Faibussowitsch PetscCall(VecGetSize(Y, &N)); 7769566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 7779566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 778c4762a1bSJed Brown f[0] = -y[0]; 779ad540459SPierre Jolivet for (i = 1; i < N - 1; i++) f[i] = y[i - 1] - y[i]; 780c4762a1bSJed Brown f[N - 1] = y[N - 2]; 7819566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 7829566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 783c4762a1bSJed Brown /* Left hand side = ydot - f(y) */ 7849566063dSJacob Faibussowitsch PetscCall(VecAYPX(F, -1.0, Ydot)); 785c4762a1bSJed Brown PetscFunctionReturn(0); 786c4762a1bSJed Brown } 787c4762a1bSJed Brown 7889371c9d4SSatish Balay PetscErrorCode IJacobian_Hull1972C1(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) { 789c4762a1bSJed Brown const PetscScalar *y; 790c4762a1bSJed Brown PetscInt N, i, col[2]; 791c4762a1bSJed Brown PetscScalar value[2]; 792c4762a1bSJed Brown 7937510d9b0SBarry Smith PetscFunctionBeginUser; 7949566063dSJacob Faibussowitsch PetscCall(VecGetSize(Y, &N)); 7959566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 796c4762a1bSJed Brown i = 0; 7979371c9d4SSatish Balay value[0] = a + 1; 7989371c9d4SSatish Balay col[0] = 0; 7999371c9d4SSatish Balay value[1] = 0; 8009371c9d4SSatish Balay col[1] = 1; 8019566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 2, col, value, INSERT_VALUES)); 802c4762a1bSJed Brown for (i = 0; i < N; i++) { 8039371c9d4SSatish Balay value[0] = -1; 8049371c9d4SSatish Balay col[0] = i - 1; 8059371c9d4SSatish Balay value[1] = a + 1; 8069371c9d4SSatish Balay col[1] = i; 8079566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 2, col, value, INSERT_VALUES)); 808c4762a1bSJed Brown } 809c4762a1bSJed Brown i = N - 1; 8109371c9d4SSatish Balay value[0] = -1; 8119371c9d4SSatish Balay col[0] = N - 2; 8129371c9d4SSatish Balay value[1] = a; 8139371c9d4SSatish Balay col[1] = N - 1; 8149566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 2, col, value, INSERT_VALUES)); 8159566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 8169566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 8179566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 818c4762a1bSJed Brown PetscFunctionReturn(0); 819c4762a1bSJed Brown } 820c4762a1bSJed Brown 821c4762a1bSJed Brown /* Hull, 1972, Problem C2 */ 822c4762a1bSJed Brown 8239371c9d4SSatish Balay PetscErrorCode RHSFunction_Hull1972C2(TS ts, PetscReal t, Vec Y, Vec F, void *s) { 824c4762a1bSJed Brown const PetscScalar *y; 825c4762a1bSJed Brown PetscScalar *f; 826c4762a1bSJed Brown PetscInt N, i; 827c4762a1bSJed Brown 8287510d9b0SBarry Smith PetscFunctionBeginUser; 8299566063dSJacob Faibussowitsch PetscCall(VecGetSize(Y, &N)); 8309566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 8319566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 832c4762a1bSJed Brown f[0] = -y[0]; 833ad540459SPierre Jolivet for (i = 1; i < N - 1; i++) f[i] = (PetscReal)i * y[i - 1] - (PetscReal)(i + 1) * y[i]; 834c4762a1bSJed Brown f[N - 1] = (PetscReal)(N - 1) * y[N - 2]; 8359566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 8369566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 837c4762a1bSJed Brown PetscFunctionReturn(0); 838c4762a1bSJed Brown } 839c4762a1bSJed Brown 8409371c9d4SSatish Balay PetscErrorCode IFunction_Hull1972C2(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) { 841c4762a1bSJed Brown PetscScalar *f; 842c4762a1bSJed Brown const PetscScalar *y; 843c4762a1bSJed Brown PetscInt N, i; 844c4762a1bSJed Brown 8457510d9b0SBarry Smith PetscFunctionBeginUser; 8469566063dSJacob Faibussowitsch PetscCall(VecGetSize(Y, &N)); 8479566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 8489566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 849c4762a1bSJed Brown f[0] = -y[0]; 850ad540459SPierre Jolivet for (i = 1; i < N - 1; i++) f[i] = (PetscReal)i * y[i - 1] - (PetscReal)(i + 1) * y[i]; 851c4762a1bSJed Brown f[N - 1] = (PetscReal)(N - 1) * y[N - 2]; 8529566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 8539566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 854c4762a1bSJed Brown /* Left hand side = ydot - f(y) */ 8559566063dSJacob Faibussowitsch PetscCall(VecAYPX(F, -1.0, Ydot)); 856c4762a1bSJed Brown PetscFunctionReturn(0); 857c4762a1bSJed Brown } 858c4762a1bSJed Brown 8599371c9d4SSatish Balay PetscErrorCode IJacobian_Hull1972C2(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) { 860c4762a1bSJed Brown const PetscScalar *y; 861c4762a1bSJed Brown PetscInt N, i, col[2]; 862c4762a1bSJed Brown PetscScalar value[2]; 863c4762a1bSJed Brown 8647510d9b0SBarry Smith PetscFunctionBeginUser; 8659566063dSJacob Faibussowitsch PetscCall(VecGetSize(Y, &N)); 8669566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 867c4762a1bSJed Brown i = 0; 8689371c9d4SSatish Balay value[0] = a + 1; 8699371c9d4SSatish Balay col[0] = 0; 8709371c9d4SSatish Balay value[1] = 0; 8719371c9d4SSatish Balay col[1] = 1; 8729566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 2, col, value, INSERT_VALUES)); 873c4762a1bSJed Brown for (i = 0; i < N; i++) { 8749371c9d4SSatish Balay value[0] = -(PetscReal)i; 8759371c9d4SSatish Balay col[0] = i - 1; 8769371c9d4SSatish Balay value[1] = a + (PetscReal)(i + 1); 8779371c9d4SSatish Balay col[1] = i; 8789566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 2, col, value, INSERT_VALUES)); 879c4762a1bSJed Brown } 880c4762a1bSJed Brown i = N - 1; 8819371c9d4SSatish Balay value[0] = -(PetscReal)(N - 1); 8829371c9d4SSatish Balay col[0] = N - 2; 8839371c9d4SSatish Balay value[1] = a; 8849371c9d4SSatish Balay col[1] = N - 1; 8859566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 2, col, value, INSERT_VALUES)); 8869566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 8879566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 8889566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 889c4762a1bSJed Brown PetscFunctionReturn(0); 890c4762a1bSJed Brown } 891c4762a1bSJed Brown 892c4762a1bSJed Brown /* Hull, 1972, Problem C3 and C4 */ 893c4762a1bSJed Brown 8949371c9d4SSatish Balay PetscErrorCode RHSFunction_Hull1972C34(TS ts, PetscReal t, Vec Y, Vec F, void *s) { 895c4762a1bSJed Brown PetscScalar *f; 896c4762a1bSJed Brown const PetscScalar *y; 897c4762a1bSJed Brown PetscInt N, i; 898c4762a1bSJed Brown 8997510d9b0SBarry Smith PetscFunctionBeginUser; 9009566063dSJacob Faibussowitsch PetscCall(VecGetSize(Y, &N)); 9019566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 9029566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 903c4762a1bSJed Brown f[0] = -2.0 * y[0] + y[1]; 904ad540459SPierre Jolivet for (i = 1; i < N - 1; i++) f[i] = y[i - 1] - 2.0 * y[i] + y[i + 1]; 905c4762a1bSJed Brown f[N - 1] = y[N - 2] - 2.0 * y[N - 1]; 9069566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 9079566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 908c4762a1bSJed Brown PetscFunctionReturn(0); 909c4762a1bSJed Brown } 910c4762a1bSJed Brown 9119371c9d4SSatish Balay PetscErrorCode IFunction_Hull1972C34(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) { 912c4762a1bSJed Brown PetscScalar *f; 913c4762a1bSJed Brown const PetscScalar *y; 914c4762a1bSJed Brown PetscInt N, i; 915c4762a1bSJed Brown 9167510d9b0SBarry Smith PetscFunctionBeginUser; 9179566063dSJacob Faibussowitsch PetscCall(VecGetSize(Y, &N)); 9189566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 9199566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 920c4762a1bSJed Brown f[0] = -2.0 * y[0] + y[1]; 921ad540459SPierre Jolivet for (i = 1; i < N - 1; i++) f[i] = y[i - 1] - 2.0 * y[i] + y[i + 1]; 922c4762a1bSJed Brown f[N - 1] = y[N - 2] - 2.0 * y[N - 1]; 9239566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 9249566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 925c4762a1bSJed Brown /* Left hand side = ydot - f(y) */ 9269566063dSJacob Faibussowitsch PetscCall(VecAYPX(F, -1.0, Ydot)); 927c4762a1bSJed Brown PetscFunctionReturn(0); 928c4762a1bSJed Brown } 929c4762a1bSJed Brown 9309371c9d4SSatish Balay PetscErrorCode IJacobian_Hull1972C34(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) { 931c4762a1bSJed Brown const PetscScalar *y; 932c4762a1bSJed Brown PetscScalar value[3]; 933c4762a1bSJed Brown PetscInt N, i, col[3]; 934c4762a1bSJed Brown 9357510d9b0SBarry Smith PetscFunctionBeginUser; 9369566063dSJacob Faibussowitsch PetscCall(VecGetSize(Y, &N)); 9379566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y, &y)); 938c4762a1bSJed Brown for (i = 0; i < N; i++) { 939c4762a1bSJed Brown if (i == 0) { 9409371c9d4SSatish Balay value[0] = a + 2; 9419371c9d4SSatish Balay col[0] = i; 9429371c9d4SSatish Balay value[1] = -1; 9439371c9d4SSatish Balay col[1] = i + 1; 9449371c9d4SSatish Balay value[2] = 0; 9459371c9d4SSatish Balay col[2] = i + 2; 946c4762a1bSJed Brown } else if (i == N - 1) { 9479371c9d4SSatish Balay value[0] = 0; 9489371c9d4SSatish Balay col[0] = i - 2; 9499371c9d4SSatish Balay value[1] = -1; 9509371c9d4SSatish Balay col[1] = i - 1; 9519371c9d4SSatish Balay value[2] = a + 2; 9529371c9d4SSatish Balay col[2] = i; 953c4762a1bSJed Brown } else { 9549371c9d4SSatish Balay value[0] = -1; 9559371c9d4SSatish Balay col[0] = i - 1; 9569371c9d4SSatish Balay value[1] = a + 2; 9579371c9d4SSatish Balay col[1] = i; 9589371c9d4SSatish Balay value[2] = -1; 9599371c9d4SSatish Balay col[2] = i + 1; 960c4762a1bSJed Brown } 9619566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 3, col, value, INSERT_VALUES)); 962c4762a1bSJed Brown } 9639566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 9649566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 9659566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y, &y)); 966c4762a1bSJed Brown PetscFunctionReturn(0); 967c4762a1bSJed Brown } 968c4762a1bSJed Brown 969c4762a1bSJed Brown /***************************************************************************/ 970c4762a1bSJed Brown 971c4762a1bSJed Brown /* Sets the initial solution for the IVP and sets up the function pointers*/ 9729371c9d4SSatish Balay PetscErrorCode Initialize(Vec Y, void *s) { 973c4762a1bSJed Brown char *p = (char *)s; 974c4762a1bSJed Brown PetscScalar *y; 975c4762a1bSJed Brown PetscReal t0; 976*11cc89d2SBarry Smith PetscInt N; 977c4762a1bSJed Brown PetscBool flg; 978c4762a1bSJed Brown 9797510d9b0SBarry Smith PetscFunctionBeginUser; 980*11cc89d2SBarry Smith PetscCall(GetSize((const char *)s, &N)); 981c4762a1bSJed Brown VecZeroEntries(Y); 9829566063dSJacob Faibussowitsch PetscCall(VecGetArray(Y, &y)); 983c4762a1bSJed Brown if (!strcmp(p, "hull1972a1")) { 984c4762a1bSJed Brown y[0] = 1.0; 985c4762a1bSJed Brown RHSFunction = RHSFunction_Hull1972A1; 986c4762a1bSJed Brown RHSJacobian = RHSJacobian_Hull1972A1; 987c4762a1bSJed Brown IFunction = IFunction_Hull1972A1; 988c4762a1bSJed Brown IJacobian = IJacobian_Hull1972A1; 989c4762a1bSJed Brown } else if (!strcmp(p, "hull1972a2")) { 990c4762a1bSJed Brown y[0] = 1.0; 991c4762a1bSJed Brown RHSFunction = RHSFunction_Hull1972A2; 992c4762a1bSJed Brown RHSJacobian = RHSJacobian_Hull1972A2; 993c4762a1bSJed Brown IFunction = IFunction_Hull1972A2; 994c4762a1bSJed Brown IJacobian = IJacobian_Hull1972A2; 995c4762a1bSJed Brown } else if (!strcmp(p, "hull1972a3")) { 996c4762a1bSJed Brown y[0] = 1.0; 997c4762a1bSJed Brown RHSFunction = RHSFunction_Hull1972A3; 998c4762a1bSJed Brown RHSJacobian = RHSJacobian_Hull1972A3; 999c4762a1bSJed Brown IFunction = IFunction_Hull1972A3; 1000c4762a1bSJed Brown IJacobian = IJacobian_Hull1972A3; 1001c4762a1bSJed Brown } else if (!strcmp(p, "hull1972a4")) { 1002c4762a1bSJed Brown y[0] = 1.0; 1003c4762a1bSJed Brown RHSFunction = RHSFunction_Hull1972A4; 1004c4762a1bSJed Brown RHSJacobian = RHSJacobian_Hull1972A4; 1005c4762a1bSJed Brown IFunction = IFunction_Hull1972A4; 1006c4762a1bSJed Brown IJacobian = IJacobian_Hull1972A4; 1007c4762a1bSJed Brown } else if (!strcmp(p, "hull1972a5")) { 1008c4762a1bSJed Brown y[0] = 4.0; 1009c4762a1bSJed Brown RHSFunction = RHSFunction_Hull1972A5; 1010c4762a1bSJed Brown RHSJacobian = RHSJacobian_Hull1972A5; 1011c4762a1bSJed Brown IFunction = IFunction_Hull1972A5; 1012c4762a1bSJed Brown IJacobian = IJacobian_Hull1972A5; 1013c4762a1bSJed Brown } else if (!strcmp(p, "hull1972b1")) { 1014c4762a1bSJed Brown y[0] = 1.0; 1015c4762a1bSJed Brown y[1] = 3.0; 1016c4762a1bSJed Brown RHSFunction = RHSFunction_Hull1972B1; 1017c4762a1bSJed Brown IFunction = IFunction_Hull1972B1; 1018c4762a1bSJed Brown IJacobian = IJacobian_Hull1972B1; 1019c4762a1bSJed Brown } else if (!strcmp(p, "hull1972b2")) { 1020c4762a1bSJed Brown y[0] = 2.0; 1021c4762a1bSJed Brown y[1] = 0.0; 1022c4762a1bSJed Brown y[2] = 1.0; 1023c4762a1bSJed Brown RHSFunction = RHSFunction_Hull1972B2; 1024c4762a1bSJed Brown IFunction = IFunction_Hull1972B2; 1025c4762a1bSJed Brown IJacobian = IJacobian_Hull1972B2; 1026c4762a1bSJed Brown } else if (!strcmp(p, "hull1972b3")) { 1027c4762a1bSJed Brown y[0] = 1.0; 1028c4762a1bSJed Brown y[1] = 0.0; 1029c4762a1bSJed Brown y[2] = 0.0; 1030c4762a1bSJed Brown RHSFunction = RHSFunction_Hull1972B3; 1031c4762a1bSJed Brown IFunction = IFunction_Hull1972B3; 1032c4762a1bSJed Brown IJacobian = IJacobian_Hull1972B3; 1033c4762a1bSJed Brown } else if (!strcmp(p, "hull1972b4")) { 1034c4762a1bSJed Brown y[0] = 3.0; 1035c4762a1bSJed Brown y[1] = 0.0; 1036c4762a1bSJed Brown y[2] = 0.0; 1037c4762a1bSJed Brown RHSFunction = RHSFunction_Hull1972B4; 1038c4762a1bSJed Brown IFunction = IFunction_Hull1972B4; 1039c4762a1bSJed Brown IJacobian = IJacobian_Hull1972B4; 1040c4762a1bSJed Brown } else if (!strcmp(p, "hull1972b5")) { 1041c4762a1bSJed Brown y[0] = 0.0; 1042c4762a1bSJed Brown y[1] = 1.0; 1043c4762a1bSJed Brown y[2] = 1.0; 1044c4762a1bSJed Brown RHSFunction = RHSFunction_Hull1972B5; 1045c4762a1bSJed Brown IFunction = IFunction_Hull1972B5; 1046c4762a1bSJed Brown IJacobian = IJacobian_Hull1972B5; 1047c4762a1bSJed Brown } else if (!strcmp(p, "kulik2013i")) { 1048c4762a1bSJed Brown t0 = 0.; 1049c4762a1bSJed Brown y[0] = PetscExpReal(PetscSinReal(t0 * t0)); 1050c4762a1bSJed Brown y[1] = PetscExpReal(5. * PetscSinReal(t0 * t0)); 1051c4762a1bSJed Brown y[2] = PetscSinReal(t0 * t0) + 1.0; 1052c4762a1bSJed Brown y[3] = PetscCosReal(t0 * t0); 1053c4762a1bSJed Brown RHSFunction = RHSFunction_Kulikov2013I; 1054c4762a1bSJed Brown RHSJacobian = RHSJacobian_Kulikov2013I; 1055c4762a1bSJed Brown IFunction = IFunction_Kulikov2013I; 1056c4762a1bSJed Brown IJacobian = IJacobian_Kulikov2013I; 1057c4762a1bSJed Brown } else if (!strcmp(p, "hull1972c1")) { 1058c4762a1bSJed Brown y[0] = 1.0; 1059c4762a1bSJed Brown RHSFunction = RHSFunction_Hull1972C1; 1060c4762a1bSJed Brown IFunction = IFunction_Hull1972C1; 1061c4762a1bSJed Brown IJacobian = IJacobian_Hull1972C1; 1062c4762a1bSJed Brown } else if (!strcmp(p, "hull1972c2")) { 1063c4762a1bSJed Brown y[0] = 1.0; 1064c4762a1bSJed Brown RHSFunction = RHSFunction_Hull1972C2; 1065c4762a1bSJed Brown IFunction = IFunction_Hull1972C2; 1066c4762a1bSJed Brown IJacobian = IJacobian_Hull1972C2; 1067*11cc89d2SBarry Smith } else if (!strcmp(p, "hull1972c3") || !strcmp(p, "hull1972c4")) { 1068c4762a1bSJed Brown y[0] = 1.0; 1069c4762a1bSJed Brown RHSFunction = RHSFunction_Hull1972C34; 1070c4762a1bSJed Brown IFunction = IFunction_Hull1972C34; 1071c4762a1bSJed Brown IJacobian = IJacobian_Hull1972C34; 1072c4762a1bSJed Brown } 10739566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalarArray(NULL, NULL, "-yinit", y, &N, &flg)); 10749566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Y, &y)); 1075c4762a1bSJed Brown PetscFunctionReturn(0); 1076c4762a1bSJed Brown } 1077c4762a1bSJed Brown 1078c4762a1bSJed Brown /* Calculates the exact solution to problems that have one */ 10799371c9d4SSatish Balay PetscErrorCode ExactSolution(Vec Y, void *s, PetscReal t, PetscBool *flag) { 1080c4762a1bSJed Brown char *p = (char *)s; 1081c4762a1bSJed Brown PetscScalar *y; 1082c4762a1bSJed Brown 10837510d9b0SBarry Smith PetscFunctionBeginUser; 1084c4762a1bSJed Brown if (!strcmp(p, "hull1972a1")) { 10859566063dSJacob Faibussowitsch PetscCall(VecGetArray(Y, &y)); 1086c4762a1bSJed Brown y[0] = PetscExpReal(-t); 1087c4762a1bSJed Brown *flag = PETSC_TRUE; 10889566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Y, &y)); 1089c4762a1bSJed Brown } else if (!strcmp(p, "hull1972a2")) { 10909566063dSJacob Faibussowitsch PetscCall(VecGetArray(Y, &y)); 1091c4762a1bSJed Brown y[0] = 1.0 / PetscSqrtReal(t + 1); 1092c4762a1bSJed Brown *flag = PETSC_TRUE; 10939566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Y, &y)); 1094c4762a1bSJed Brown } else if (!strcmp(p, "hull1972a3")) { 10959566063dSJacob Faibussowitsch PetscCall(VecGetArray(Y, &y)); 1096c4762a1bSJed Brown y[0] = PetscExpReal(PetscSinReal(t)); 1097c4762a1bSJed Brown *flag = PETSC_TRUE; 10989566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Y, &y)); 1099c4762a1bSJed Brown } else if (!strcmp(p, "hull1972a4")) { 11009566063dSJacob Faibussowitsch PetscCall(VecGetArray(Y, &y)); 1101c4762a1bSJed Brown y[0] = 20.0 / (1 + 19.0 * PetscExpReal(-t / 4.0)); 1102c4762a1bSJed Brown *flag = PETSC_TRUE; 11039566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Y, &y)); 1104c4762a1bSJed Brown } else if (!strcmp(p, "kulik2013i")) { 11059566063dSJacob Faibussowitsch PetscCall(VecGetArray(Y, &y)); 1106c4762a1bSJed Brown y[0] = PetscExpReal(PetscSinReal(t * t)); 1107c4762a1bSJed Brown y[1] = PetscExpReal(5. * PetscSinReal(t * t)); 1108c4762a1bSJed Brown y[2] = PetscSinReal(t * t) + 1.0; 1109c4762a1bSJed Brown y[3] = PetscCosReal(t * t); 1110c4762a1bSJed Brown *flag = PETSC_TRUE; 11119566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Y, &y)); 1112c4762a1bSJed Brown } else { 11139566063dSJacob Faibussowitsch PetscCall(VecSet(Y, 0)); 1114c4762a1bSJed Brown *flag = PETSC_FALSE; 1115c4762a1bSJed Brown } 1116c4762a1bSJed Brown PetscFunctionReturn(0); 1117c4762a1bSJed Brown } 1118c4762a1bSJed Brown 1119c4762a1bSJed Brown /* Solves the specified ODE and computes the error if exact solution is available */ 11209371c9d4SSatish Balay PetscErrorCode SolveODE(char *ptype, PetscReal dt, PetscReal tfinal, PetscInt maxiter, PetscReal *error, PetscBool *exact_flag) { 1121c4762a1bSJed Brown TS ts; /* time-integrator */ 1122c4762a1bSJed Brown Vec Y; /* Solution vector */ 1123c4762a1bSJed Brown Vec Yex; /* Exact solution */ 1124c4762a1bSJed Brown PetscInt N; /* Size of the system of equations */ 1125c4762a1bSJed Brown TSType time_scheme; /* Type of time-integration scheme */ 1126c4762a1bSJed Brown Mat Jac = NULL; /* Jacobian matrix */ 1127c4762a1bSJed Brown Vec Yerr; /* Auxiliary solution vector */ 1128c4762a1bSJed Brown PetscReal err_norm; /* Estimated error norm */ 1129c4762a1bSJed Brown PetscReal final_time; /* Actual final time from the integrator */ 1130c4762a1bSJed Brown 11317510d9b0SBarry Smith PetscFunctionBeginUser; 1132*11cc89d2SBarry Smith PetscCall(GetSize((const char *)&ptype[0], &N)); 11333c633725SBarry Smith PetscCheck(N >= 0, PETSC_COMM_WORLD, PETSC_ERR_ARG_SIZ, "Illegal problem specification."); 11349566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &Y)); 11359566063dSJacob Faibussowitsch PetscCall(VecSetSizes(Y, N, PETSC_DECIDE)); 11369566063dSJacob Faibussowitsch PetscCall(VecSetUp(Y)); 11379566063dSJacob Faibussowitsch PetscCall(VecSet(Y, 0)); 1138c4762a1bSJed Brown 1139c4762a1bSJed Brown /* Initialize the problem */ 11409566063dSJacob Faibussowitsch PetscCall(Initialize(Y, &ptype[0])); 1141c4762a1bSJed Brown 1142c4762a1bSJed Brown /* Create and initialize the time-integrator */ 11439566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 1144c4762a1bSJed Brown /* Default time integration options */ 11459566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSRK)); 11469566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(ts, maxiter)); 11479566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, tfinal)); 11489566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, dt)); 11499566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 1150c4762a1bSJed Brown /* Read command line options for time integration */ 11519566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 1152c4762a1bSJed Brown /* Set solution vector */ 11539566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, Y)); 1154c4762a1bSJed Brown /* Specify left/right-hand side functions */ 11559566063dSJacob Faibussowitsch PetscCall(TSGetType(ts, &time_scheme)); 1156c4762a1bSJed Brown 1157c4762a1bSJed Brown if ((!strcmp(time_scheme, TSEULER)) || (!strcmp(time_scheme, TSRK)) || (!strcmp(time_scheme, TSSSP) || (!strcmp(time_scheme, TSGLEE)))) { 1158c4762a1bSJed Brown /* Explicit time-integration -> specify right-hand side function ydot = f(y) */ 11599566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, &ptype[0])); 11609566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &Jac)); 11619566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Jac, PETSC_DECIDE, PETSC_DECIDE, N, N)); 11629566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(Jac)); 11639566063dSJacob Faibussowitsch PetscCall(MatSetUp(Jac)); 11649566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(ts, Jac, Jac, RHSJacobian, &ptype[0])); 1165c4762a1bSJed Brown } else if ((!strcmp(time_scheme, TSTHETA)) || (!strcmp(time_scheme, TSBEULER)) || (!strcmp(time_scheme, TSCN)) || (!strcmp(time_scheme, TSALPHA)) || (!strcmp(time_scheme, TSARKIMEX))) { 1166c4762a1bSJed Brown /* Implicit time-integration -> specify left-hand side function ydot-f(y) = 0 */ 1167c4762a1bSJed Brown /* and its Jacobian function */ 11689566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts, NULL, IFunction, &ptype[0])); 11699566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &Jac)); 11709566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Jac, PETSC_DECIDE, PETSC_DECIDE, N, N)); 11719566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(Jac)); 11729566063dSJacob Faibussowitsch PetscCall(MatSetUp(Jac)); 11739566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts, Jac, Jac, IJacobian, &ptype[0])); 1174c4762a1bSJed Brown } 1175c4762a1bSJed Brown 1176c4762a1bSJed Brown /* Solve */ 11779566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, Y)); 11789566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &final_time)); 1179c4762a1bSJed Brown 1180c4762a1bSJed Brown /* Get the estimated error, if available */ 11819566063dSJacob Faibussowitsch PetscCall(VecDuplicate(Y, &Yerr)); 11829566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Yerr)); 11839566063dSJacob Faibussowitsch PetscCall(TSGetTimeError(ts, 0, &Yerr)); 11849566063dSJacob Faibussowitsch PetscCall(VecNorm(Yerr, NORM_2, &err_norm)); 11859566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Yerr)); 118663a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Estimated Error = %e.\n", (double)err_norm)); 1187c4762a1bSJed Brown 1188c4762a1bSJed Brown /* Exact solution */ 11899566063dSJacob Faibussowitsch PetscCall(VecDuplicate(Y, &Yex)); 119063a3b9bcSJacob Faibussowitsch if (PetscAbsReal(final_time - tfinal) > 2. * PETSC_MACHINE_EPSILON) { 11919566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Note: There is a difference between the prescribed final time %g and the actual final time, %g.\n", (double)tfinal, (double)final_time)); 1192c4762a1bSJed Brown } 11939566063dSJacob Faibussowitsch PetscCall(ExactSolution(Yex, &ptype[0], final_time, exact_flag)); 1194c4762a1bSJed Brown 1195c4762a1bSJed Brown /* Calculate Error */ 11969566063dSJacob Faibussowitsch PetscCall(VecAYPX(Yex, -1.0, Y)); 11979566063dSJacob Faibussowitsch PetscCall(VecNorm(Yex, NORM_2, error)); 1198c4762a1bSJed Brown *error = PetscSqrtReal(((*error) * (*error)) / N); 1199c4762a1bSJed Brown 1200c4762a1bSJed Brown /* Clean up and finalize */ 12019566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Jac)); 12029566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 12039566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Yex)); 12049566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y)); 1205c4762a1bSJed Brown 1206c4762a1bSJed Brown PetscFunctionReturn(0); 1207c4762a1bSJed Brown } 1208c4762a1bSJed Brown 12099371c9d4SSatish Balay int main(int argc, char **argv) { 1210c4762a1bSJed Brown char ptype[256] = "hull1972a1"; /* Problem specification */ 1211c4762a1bSJed Brown PetscInt n_refine = 1; /* Number of refinement levels for convergence analysis */ 1212c4762a1bSJed Brown PetscReal refine_fac = 2.0; /* Refinement factor for dt */ 1213c4762a1bSJed Brown PetscReal dt_initial = 0.01; /* Initial default value of dt */ 1214c4762a1bSJed Brown PetscReal dt; 1215c4762a1bSJed Brown PetscReal tfinal = 20.0; /* Final time for the time-integration */ 1216c4762a1bSJed Brown PetscInt maxiter = 100000; /* Maximum number of time-integration iterations */ 1217c4762a1bSJed Brown PetscReal *error; /* Array to store the errors for convergence analysis */ 1218c4762a1bSJed Brown PetscMPIInt size; /* No of processors */ 1219c4762a1bSJed Brown PetscBool flag; /* Flag denoting availability of exact solution */ 1220*11cc89d2SBarry Smith PetscInt r, N; 1221c4762a1bSJed Brown 1222c4762a1bSJed Brown /* Initialize program */ 1223327415f7SBarry Smith PetscFunctionBeginUser; 1224*11cc89d2SBarry Smith PetscCall(GetSize(&ptype[0], &N)); 12259566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 1226c4762a1bSJed Brown 1227c4762a1bSJed Brown /* Check if running with only 1 proc */ 12289566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 12293c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs"); 1230c4762a1bSJed Brown 1231d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "ex31", NULL); 12329566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-problem", "Problem specification", "<hull1972a1>", ptype, ptype, sizeof(ptype), NULL)); 12339566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-refinement_levels", "Number of refinement levels for convergence analysis", "<1>", n_refine, &n_refine, NULL)); 12349566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-refinement_factor", "Refinement factor for dt", "<2.0>", refine_fac, &refine_fac, NULL)); 12359566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-dt", "Time step size (for convergence analysis, initial time step)", "<0.01>", dt_initial, &dt_initial, NULL)); 12369566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-final_time", "Final time for the time-integration", "<20.0>", tfinal, &tfinal, NULL)); 1237d0609cedSBarry Smith PetscOptionsEnd(); 1238c4762a1bSJed Brown 12399566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n_refine, &error)); 1240c4762a1bSJed Brown for (r = 0, dt = dt_initial; r < n_refine; r++) { 1241c4762a1bSJed Brown error[r] = 0; 1242c4762a1bSJed Brown if (r > 0) dt /= refine_fac; 1243c4762a1bSJed Brown 1244*11cc89d2SBarry Smith PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Solving ODE \"%s\" with dt %f, final time %f and system size %" PetscInt_FMT ".\n", ptype, (double)dt, (double)tfinal, N)); 12459566063dSJacob Faibussowitsch PetscCall(SolveODE(&ptype[0], dt, tfinal, maxiter, &error[r], &flag)); 1246c4762a1bSJed Brown if (flag) { 1247c4762a1bSJed Brown /* If exact solution available for the specified ODE */ 1248c4762a1bSJed Brown if (r > 0) { 1249c4762a1bSJed Brown PetscReal conv_rate = (PetscLogReal(error[r]) - PetscLogReal(error[r - 1])) / (-PetscLogReal(refine_fac)); 12509566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error = %E,\tConvergence rate = %f.\n", (double)error[r], (double)conv_rate)); 1251c4762a1bSJed Brown } else { 125263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error = %E.\n", (double)error[r])); 1253c4762a1bSJed Brown } 1254c4762a1bSJed Brown } 1255c4762a1bSJed Brown } 12569566063dSJacob Faibussowitsch PetscCall(PetscFree(error)); 12579566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 1258b122ec5aSJacob Faibussowitsch return 0; 1259c4762a1bSJed Brown } 1260c4762a1bSJed Brown 1261c4762a1bSJed Brown /*TEST 1262c4762a1bSJed Brown 1263c4762a1bSJed Brown test: 1264c4762a1bSJed Brown suffix: 2 1265c4762a1bSJed Brown args: -ts_type glee -final_time 5 -ts_adapt_type none 1266c4762a1bSJed Brown timeoutfactor: 3 1267c4762a1bSJed Brown requires: !single 1268c4762a1bSJed Brown 1269c4762a1bSJed Brown test: 1270c4762a1bSJed Brown suffix: 3 1271c4762a1bSJed Brown args: -ts_type glee -final_time 5 -ts_adapt_type glee -ts_adapt_monitor -ts_max_steps 50 -problem hull1972a3 -ts_adapt_glee_use_local 1 1272c4762a1bSJed Brown timeoutfactor: 3 1273c4762a1bSJed Brown requires: !single 1274c4762a1bSJed Brown 1275c4762a1bSJed Brown test: 1276c4762a1bSJed Brown suffix: 4 1277c4762a1bSJed Brown args: -ts_type glee -final_time 5 -ts_adapt_type glee -ts_adapt_monitor -ts_max_steps 50 -problem hull1972a3 -ts_max_reject 100 -ts_adapt_glee_use_local 0 1278c4762a1bSJed Brown timeoutfactor: 3 1279c4762a1bSJed Brown requires: !single !__float128 1280c4762a1bSJed Brown 1281c4762a1bSJed Brown TEST*/ 1282