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 */ 52c4762a1bSJed Brown PetscInt GetSize(const char *p) 53c4762a1bSJed Brown { 54c4762a1bSJed Brown PetscFunctionBegin; 55c4762a1bSJed Brown if ((!strcmp(p,"hull1972a1")) 56c4762a1bSJed Brown ||(!strcmp(p,"hull1972a2")) 57c4762a1bSJed Brown ||(!strcmp(p,"hull1972a3")) 58c4762a1bSJed Brown ||(!strcmp(p,"hull1972a4")) 59c4762a1bSJed Brown ||(!strcmp(p,"hull1972a5"))) PetscFunctionReturn(1); 60c4762a1bSJed Brown else if (!strcmp(p,"hull1972b1")) PetscFunctionReturn(2); 61c4762a1bSJed Brown else if ((!strcmp(p,"hull1972b2")) 62c4762a1bSJed Brown ||(!strcmp(p,"hull1972b3")) 63c4762a1bSJed Brown ||(!strcmp(p,"hull1972b4")) 64c4762a1bSJed Brown ||(!strcmp(p,"hull1972b5"))) PetscFunctionReturn(3); 65c4762a1bSJed Brown else if ((!strcmp(p,"kulik2013i"))) PetscFunctionReturn(4); 66c4762a1bSJed Brown else if ((!strcmp(p,"hull1972c1")) 67c4762a1bSJed Brown ||(!strcmp(p,"hull1972c2")) 68c4762a1bSJed Brown ||(!strcmp(p,"hull1972c3"))) PetscFunctionReturn(10); 69c4762a1bSJed Brown else if (!strcmp(p,"hull1972c4")) PetscFunctionReturn(51); 70c4762a1bSJed Brown else PetscFunctionReturn(-1); 71c4762a1bSJed Brown } 72c4762a1bSJed Brown 73c4762a1bSJed Brown /****************************************************************/ 74c4762a1bSJed Brown 75c4762a1bSJed Brown /* Problem specific functions */ 76c4762a1bSJed Brown 77c4762a1bSJed Brown /* Hull, 1972, Problem A1 */ 78c4762a1bSJed Brown 79c4762a1bSJed Brown PetscErrorCode RHSFunction_Hull1972A1(TS ts, PetscReal t, Vec Y, Vec F, void *s) 80c4762a1bSJed Brown { 81c4762a1bSJed Brown PetscScalar *f; 82c4762a1bSJed Brown const PetscScalar *y; 83c4762a1bSJed Brown 84c4762a1bSJed Brown PetscFunctionBegin; 859566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 869566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 87c4762a1bSJed Brown f[0] = -y[0]; 889566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 899566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 90c4762a1bSJed Brown PetscFunctionReturn(0); 91c4762a1bSJed Brown } 92c4762a1bSJed Brown 93c4762a1bSJed Brown PetscErrorCode RHSJacobian_Hull1972A1(TS ts, PetscReal t, Vec Y, Mat A, Mat B, void *s) 94c4762a1bSJed Brown { 95c4762a1bSJed Brown const PetscScalar *y; 96c4762a1bSJed Brown PetscInt row = 0,col = 0; 97c4762a1bSJed Brown PetscScalar value = -1.0; 98c4762a1bSJed Brown 99c4762a1bSJed Brown PetscFunctionBegin; 1009566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 1019566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&row,1,&col,&value,INSERT_VALUES)); 1029566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 1039566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY)); 1049566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 105c4762a1bSJed Brown PetscFunctionReturn(0); 106c4762a1bSJed Brown } 107c4762a1bSJed Brown 108c4762a1bSJed Brown PetscErrorCode IFunction_Hull1972A1(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) 109c4762a1bSJed Brown { 110c4762a1bSJed Brown const PetscScalar *y; 111c4762a1bSJed Brown PetscScalar *f; 112c4762a1bSJed Brown 113c4762a1bSJed Brown PetscFunctionBegin; 1149566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 1159566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 116c4762a1bSJed Brown f[0] = -y[0]; 1179566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 1189566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 119c4762a1bSJed Brown /* Left hand side = ydot - f(y) */ 1209566063dSJacob Faibussowitsch PetscCall(VecAYPX(F,-1.0,Ydot)); 121c4762a1bSJed Brown PetscFunctionReturn(0); 122c4762a1bSJed Brown } 123c4762a1bSJed Brown 124c4762a1bSJed Brown PetscErrorCode IJacobian_Hull1972A1(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) 125c4762a1bSJed Brown { 126c4762a1bSJed Brown const PetscScalar *y; 127c4762a1bSJed Brown PetscInt row = 0,col = 0; 128c4762a1bSJed Brown PetscScalar value = a - 1.0; 129c4762a1bSJed Brown 130c4762a1bSJed Brown PetscFunctionBegin; 1319566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 1329566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&row,1,&col,&value,INSERT_VALUES)); 1339566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 1349566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY)); 1359566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 136c4762a1bSJed Brown PetscFunctionReturn(0); 137c4762a1bSJed Brown } 138c4762a1bSJed Brown 139c4762a1bSJed Brown /* Hull, 1972, Problem A2 */ 140c4762a1bSJed Brown 141c4762a1bSJed Brown PetscErrorCode RHSFunction_Hull1972A2(TS ts, PetscReal t, Vec Y, Vec F, void *s) 142c4762a1bSJed Brown { 143c4762a1bSJed Brown const PetscScalar *y; 144c4762a1bSJed Brown PetscScalar *f; 145c4762a1bSJed Brown 146c4762a1bSJed Brown PetscFunctionBegin; 1479566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 1489566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 149c4762a1bSJed Brown f[0] = -0.5*y[0]*y[0]*y[0]; 1509566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 1519566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 152c4762a1bSJed Brown PetscFunctionReturn(0); 153c4762a1bSJed Brown } 154c4762a1bSJed Brown 155c4762a1bSJed Brown PetscErrorCode RHSJacobian_Hull1972A2(TS ts, PetscReal t, Vec Y, Mat A, Mat B, void *s) 156c4762a1bSJed Brown { 157c4762a1bSJed Brown const PetscScalar *y; 158c4762a1bSJed Brown PetscInt row = 0,col = 0; 159c4762a1bSJed Brown PetscScalar value; 160c4762a1bSJed Brown 161c4762a1bSJed Brown PetscFunctionBegin; 1629566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 163c4762a1bSJed Brown value = -0.5*3.0*y[0]*y[0]; 1649566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&row,1,&col,&value,INSERT_VALUES)); 1659566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 1669566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY)); 1679566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 168c4762a1bSJed Brown PetscFunctionReturn(0); 169c4762a1bSJed Brown } 170c4762a1bSJed Brown 171c4762a1bSJed Brown PetscErrorCode IFunction_Hull1972A2(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) 172c4762a1bSJed Brown { 173c4762a1bSJed Brown PetscScalar *f; 174c4762a1bSJed Brown const PetscScalar *y; 175c4762a1bSJed Brown 176c4762a1bSJed Brown PetscFunctionBegin; 1779566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 1789566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 179c4762a1bSJed Brown f[0] = -0.5*y[0]*y[0]*y[0]; 1809566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 1819566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 182c4762a1bSJed Brown /* Left hand side = ydot - f(y) */ 1839566063dSJacob Faibussowitsch PetscCall(VecAYPX(F,-1.0,Ydot)); 184c4762a1bSJed Brown PetscFunctionReturn(0); 185c4762a1bSJed Brown } 186c4762a1bSJed Brown 187c4762a1bSJed Brown PetscErrorCode IJacobian_Hull1972A2(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) 188c4762a1bSJed Brown { 189c4762a1bSJed Brown const PetscScalar *y; 190c4762a1bSJed Brown PetscInt row = 0,col = 0; 191c4762a1bSJed Brown PetscScalar value; 192c4762a1bSJed Brown 193c4762a1bSJed Brown PetscFunctionBegin; 1949566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 195c4762a1bSJed Brown value = a + 0.5*3.0*y[0]*y[0]; 1969566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&row,1,&col,&value,INSERT_VALUES)); 1979566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 1989566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY)); 1999566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 200c4762a1bSJed Brown PetscFunctionReturn(0); 201c4762a1bSJed Brown } 202c4762a1bSJed Brown 203c4762a1bSJed Brown /* Hull, 1972, Problem A3 */ 204c4762a1bSJed Brown 205c4762a1bSJed Brown PetscErrorCode RHSFunction_Hull1972A3(TS ts, PetscReal t, Vec Y, Vec F, void *s) 206c4762a1bSJed Brown { 207c4762a1bSJed Brown const PetscScalar *y; 208c4762a1bSJed Brown PetscScalar *f; 209c4762a1bSJed Brown 210c4762a1bSJed Brown PetscFunctionBegin; 2119566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 2129566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 213c4762a1bSJed Brown f[0] = y[0]*PetscCosReal(t); 2149566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 2159566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 216c4762a1bSJed Brown PetscFunctionReturn(0); 217c4762a1bSJed Brown } 218c4762a1bSJed Brown 219c4762a1bSJed Brown PetscErrorCode RHSJacobian_Hull1972A3(TS ts, PetscReal t, Vec Y, Mat A, Mat B, void *s) 220c4762a1bSJed Brown { 221c4762a1bSJed Brown const PetscScalar *y; 222c4762a1bSJed Brown PetscInt row = 0,col = 0; 223c4762a1bSJed Brown PetscScalar value = PetscCosReal(t); 224c4762a1bSJed Brown 225c4762a1bSJed Brown PetscFunctionBegin; 2269566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 2279566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&row,1,&col,&value,INSERT_VALUES)); 2289566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 2299566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY)); 2309566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 231c4762a1bSJed Brown PetscFunctionReturn(0); 232c4762a1bSJed Brown } 233c4762a1bSJed Brown 234c4762a1bSJed Brown PetscErrorCode IFunction_Hull1972A3(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) 235c4762a1bSJed Brown { 236c4762a1bSJed Brown PetscScalar *f; 237c4762a1bSJed Brown const PetscScalar *y; 238c4762a1bSJed Brown 239c4762a1bSJed Brown PetscFunctionBegin; 2409566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 2419566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 242c4762a1bSJed Brown f[0] = y[0]*PetscCosReal(t); 2439566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 2449566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 245c4762a1bSJed Brown /* Left hand side = ydot - f(y) */ 2469566063dSJacob Faibussowitsch PetscCall(VecAYPX(F,-1.0,Ydot)); 247c4762a1bSJed Brown PetscFunctionReturn(0); 248c4762a1bSJed Brown } 249c4762a1bSJed Brown 250c4762a1bSJed Brown PetscErrorCode IJacobian_Hull1972A3(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) 251c4762a1bSJed Brown { 252c4762a1bSJed Brown const PetscScalar *y; 253c4762a1bSJed Brown PetscInt row = 0,col = 0; 254c4762a1bSJed Brown PetscScalar value = a - PetscCosReal(t); 255c4762a1bSJed Brown 256c4762a1bSJed Brown PetscFunctionBegin; 2579566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 2589566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&row,1,&col,&value,INSERT_VALUES)); 2599566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 2609566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY)); 2619566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 262c4762a1bSJed Brown PetscFunctionReturn(0); 263c4762a1bSJed Brown } 264c4762a1bSJed Brown 265c4762a1bSJed Brown /* Hull, 1972, Problem A4 */ 266c4762a1bSJed Brown 267c4762a1bSJed Brown PetscErrorCode RHSFunction_Hull1972A4(TS ts, PetscReal t, Vec Y, Vec F, void *s) 268c4762a1bSJed Brown { 269c4762a1bSJed Brown PetscScalar *f; 270c4762a1bSJed Brown const PetscScalar *y; 271c4762a1bSJed Brown 272c4762a1bSJed Brown PetscFunctionBegin; 2739566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 2749566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 275c4762a1bSJed Brown f[0] = (0.25*y[0])*(1.0-0.05*y[0]); 2769566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 2779566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 278c4762a1bSJed Brown PetscFunctionReturn(0); 279c4762a1bSJed Brown } 280c4762a1bSJed Brown 281c4762a1bSJed Brown PetscErrorCode RHSJacobian_Hull1972A4(TS ts, PetscReal t, Vec Y, Mat A, Mat B, void *s) 282c4762a1bSJed Brown { 283c4762a1bSJed Brown const PetscScalar *y; 284c4762a1bSJed Brown PetscInt row = 0,col = 0; 285c4762a1bSJed Brown PetscScalar value; 286c4762a1bSJed Brown 287c4762a1bSJed Brown PetscFunctionBegin; 2889566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 289c4762a1bSJed Brown value = 0.25*(1.0-0.05*y[0]) - (0.25*y[0])*0.05; 2909566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&row,1,&col,&value,INSERT_VALUES)); 2919566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 2929566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY)); 2939566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 294c4762a1bSJed Brown PetscFunctionReturn(0); 295c4762a1bSJed Brown } 296c4762a1bSJed Brown 297c4762a1bSJed Brown PetscErrorCode IFunction_Hull1972A4(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) 298c4762a1bSJed Brown { 299c4762a1bSJed Brown PetscScalar *f; 300c4762a1bSJed Brown const PetscScalar *y; 301c4762a1bSJed Brown 302c4762a1bSJed Brown PetscFunctionBegin; 3039566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 3049566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 305c4762a1bSJed Brown f[0] = (0.25*y[0])*(1.0-0.05*y[0]); 3069566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 3079566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 308c4762a1bSJed Brown /* Left hand side = ydot - f(y) */ 3099566063dSJacob Faibussowitsch PetscCall(VecAYPX(F,-1.0,Ydot)); 310c4762a1bSJed Brown PetscFunctionReturn(0); 311c4762a1bSJed Brown } 312c4762a1bSJed Brown 313c4762a1bSJed Brown PetscErrorCode IJacobian_Hull1972A4(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) 314c4762a1bSJed Brown { 315c4762a1bSJed Brown const PetscScalar *y; 316c4762a1bSJed Brown PetscInt row = 0,col = 0; 317c4762a1bSJed Brown PetscScalar value; 318c4762a1bSJed Brown 319c4762a1bSJed Brown PetscFunctionBegin; 3209566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 321c4762a1bSJed Brown value = a - 0.25*(1.0-0.05*y[0]) + (0.25*y[0])*0.05; 3229566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&row,1,&col,&value,INSERT_VALUES)); 3239566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 3249566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY)); 3259566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 326c4762a1bSJed Brown PetscFunctionReturn(0); 327c4762a1bSJed Brown } 328c4762a1bSJed Brown 329c4762a1bSJed Brown /* Hull, 1972, Problem A5 */ 330c4762a1bSJed Brown 331c4762a1bSJed Brown PetscErrorCode RHSFunction_Hull1972A5(TS ts, PetscReal t, Vec Y, Vec F, void *s) 332c4762a1bSJed Brown { 333c4762a1bSJed Brown PetscScalar *f; 334c4762a1bSJed Brown const PetscScalar *y; 335c4762a1bSJed Brown 336c4762a1bSJed Brown PetscFunctionBegin; 3379566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 3389566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 339c4762a1bSJed Brown f[0] = (y[0]-t)/(y[0]+t); 3409566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 3419566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 342c4762a1bSJed Brown PetscFunctionReturn(0); 343c4762a1bSJed Brown } 344c4762a1bSJed Brown 345c4762a1bSJed Brown PetscErrorCode RHSJacobian_Hull1972A5(TS ts, PetscReal t, Vec Y, Mat A, Mat B, void *s) 346c4762a1bSJed Brown { 347c4762a1bSJed Brown const PetscScalar *y; 348c4762a1bSJed Brown PetscInt row = 0,col = 0; 349c4762a1bSJed Brown PetscScalar value; 350c4762a1bSJed Brown 351c4762a1bSJed Brown PetscFunctionBegin; 3529566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 353c4762a1bSJed Brown value = 2*t/((t+y[0])*(t+y[0])); 3549566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&row,1,&col,&value,INSERT_VALUES)); 3559566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 3569566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY)); 3579566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 358c4762a1bSJed Brown PetscFunctionReturn(0); 359c4762a1bSJed Brown } 360c4762a1bSJed Brown 361c4762a1bSJed Brown PetscErrorCode IFunction_Hull1972A5(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) 362c4762a1bSJed Brown { 363c4762a1bSJed Brown PetscScalar *f; 364c4762a1bSJed Brown const PetscScalar *y; 365c4762a1bSJed Brown 366c4762a1bSJed Brown PetscFunctionBegin; 3679566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 3689566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 369c4762a1bSJed Brown f[0] = (y[0]-t)/(y[0]+t); 3709566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 3719566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 372c4762a1bSJed Brown /* Left hand side = ydot - f(y) */ 3739566063dSJacob Faibussowitsch PetscCall(VecAYPX(F,-1.0,Ydot)); 374c4762a1bSJed Brown PetscFunctionReturn(0); 375c4762a1bSJed Brown } 376c4762a1bSJed Brown 377c4762a1bSJed Brown PetscErrorCode IJacobian_Hull1972A5(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) 378c4762a1bSJed Brown { 379c4762a1bSJed Brown const PetscScalar *y; 380c4762a1bSJed Brown PetscInt row = 0,col = 0; 381c4762a1bSJed Brown PetscScalar value; 382c4762a1bSJed Brown 383c4762a1bSJed Brown PetscFunctionBegin; 3849566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 385c4762a1bSJed Brown value = a - 2*t/((t+y[0])*(t+y[0])); 3869566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&row,1,&col,&value,INSERT_VALUES)); 3879566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 3889566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY)); 3899566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 390c4762a1bSJed Brown PetscFunctionReturn(0); 391c4762a1bSJed Brown } 392c4762a1bSJed Brown 393c4762a1bSJed Brown /* Hull, 1972, Problem B1 */ 394c4762a1bSJed Brown 395c4762a1bSJed Brown PetscErrorCode RHSFunction_Hull1972B1(TS ts, PetscReal t, Vec Y, Vec F, void *s) 396c4762a1bSJed Brown { 397c4762a1bSJed Brown PetscScalar *f; 398c4762a1bSJed Brown const PetscScalar *y; 399c4762a1bSJed Brown 400c4762a1bSJed Brown PetscFunctionBegin; 4019566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 4029566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 403c4762a1bSJed Brown f[0] = 2.0*(y[0] - y[0]*y[1]); 404c4762a1bSJed Brown f[1] = -(y[1]-y[0]*y[1]); 4059566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 4069566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 407c4762a1bSJed Brown PetscFunctionReturn(0); 408c4762a1bSJed Brown } 409c4762a1bSJed Brown 410c4762a1bSJed Brown PetscErrorCode IFunction_Hull1972B1(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) 411c4762a1bSJed Brown { 412c4762a1bSJed Brown PetscScalar *f; 413c4762a1bSJed Brown const PetscScalar *y; 414c4762a1bSJed Brown 415c4762a1bSJed Brown PetscFunctionBegin; 4169566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 4179566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 418c4762a1bSJed Brown f[0] = 2.0*(y[0] - y[0]*y[1]); 419c4762a1bSJed Brown f[1] = -(y[1]-y[0]*y[1]); 4209566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 4219566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 422c4762a1bSJed Brown /* Left hand side = ydot - f(y) */ 4239566063dSJacob Faibussowitsch PetscCall(VecAYPX(F,-1.0,Ydot)); 424c4762a1bSJed Brown PetscFunctionReturn(0); 425c4762a1bSJed Brown } 426c4762a1bSJed Brown 427c4762a1bSJed Brown PetscErrorCode IJacobian_Hull1972B1(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) 428c4762a1bSJed Brown { 429c4762a1bSJed Brown const PetscScalar *y; 430c4762a1bSJed Brown PetscInt row[2] = {0,1}; 431c4762a1bSJed Brown PetscScalar value[2][2]; 432c4762a1bSJed Brown 433c4762a1bSJed Brown PetscFunctionBegin; 4349566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 435c4762a1bSJed Brown value[0][0] = a - 2.0*(1.0-y[1]); value[0][1] = 2.0*y[0]; 436c4762a1bSJed Brown value[1][0] = -y[1]; value[1][1] = a + 1.0 - y[0]; 4379566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,2,&row[0],2,&row[0],&value[0][0],INSERT_VALUES)); 4389566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 4399566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY)); 4409566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 441c4762a1bSJed Brown PetscFunctionReturn(0); 442c4762a1bSJed Brown } 443c4762a1bSJed Brown 444c4762a1bSJed Brown /* Hull, 1972, Problem B2 */ 445c4762a1bSJed Brown 446c4762a1bSJed Brown PetscErrorCode RHSFunction_Hull1972B2(TS ts, PetscReal t, Vec Y, Vec F, void *s) 447c4762a1bSJed Brown { 448c4762a1bSJed Brown PetscScalar *f; 449c4762a1bSJed Brown const PetscScalar *y; 450c4762a1bSJed Brown 451c4762a1bSJed Brown PetscFunctionBegin; 4529566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 4539566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 454c4762a1bSJed Brown f[0] = -y[0] + y[1]; 455c4762a1bSJed Brown f[1] = y[0] - 2.0*y[1] + y[2]; 456c4762a1bSJed Brown f[2] = y[1] - y[2]; 4579566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 4589566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 459c4762a1bSJed Brown PetscFunctionReturn(0); 460c4762a1bSJed Brown } 461c4762a1bSJed Brown 462c4762a1bSJed Brown PetscErrorCode IFunction_Hull1972B2(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) 463c4762a1bSJed Brown { 464c4762a1bSJed Brown PetscScalar *f; 465c4762a1bSJed Brown const PetscScalar *y; 466c4762a1bSJed Brown 467c4762a1bSJed Brown PetscFunctionBegin; 4689566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 4699566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 470c4762a1bSJed Brown f[0] = -y[0] + y[1]; 471c4762a1bSJed Brown f[1] = y[0] - 2.0*y[1] + y[2]; 472c4762a1bSJed Brown f[2] = y[1] - y[2]; 4739566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 4749566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 475c4762a1bSJed Brown /* Left hand side = ydot - f(y) */ 4769566063dSJacob Faibussowitsch PetscCall(VecAYPX(F,-1.0,Ydot)); 477c4762a1bSJed Brown PetscFunctionReturn(0); 478c4762a1bSJed Brown } 479c4762a1bSJed Brown 480c4762a1bSJed Brown PetscErrorCode IJacobian_Hull1972B2(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) 481c4762a1bSJed Brown { 482c4762a1bSJed Brown const PetscScalar *y; 483c4762a1bSJed Brown PetscInt row[3] = {0,1,2}; 484c4762a1bSJed Brown PetscScalar value[3][3]; 485c4762a1bSJed Brown 486c4762a1bSJed Brown PetscFunctionBegin; 4879566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 488c4762a1bSJed Brown value[0][0] = a + 1.0; value[0][1] = -1.0; value[0][2] = 0; 489c4762a1bSJed Brown value[1][0] = -1.0; value[1][1] = a + 2.0; value[1][2] = -1.0; 490c4762a1bSJed Brown value[2][0] = 0; value[2][1] = -1.0; value[2][2] = a + 1.0; 4919566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,3,&row[0],3,&row[0],&value[0][0],INSERT_VALUES)); 4929566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 4939566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY)); 4949566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 495c4762a1bSJed Brown PetscFunctionReturn(0); 496c4762a1bSJed Brown } 497c4762a1bSJed Brown 498c4762a1bSJed Brown /* Hull, 1972, Problem B3 */ 499c4762a1bSJed Brown 500c4762a1bSJed Brown PetscErrorCode RHSFunction_Hull1972B3(TS ts, PetscReal t, Vec Y, Vec F, void *s) 501c4762a1bSJed Brown { 502c4762a1bSJed Brown PetscScalar *f; 503c4762a1bSJed Brown const PetscScalar *y; 504c4762a1bSJed Brown 505c4762a1bSJed Brown PetscFunctionBegin; 5069566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 5079566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 508c4762a1bSJed Brown f[0] = -y[0]; 509c4762a1bSJed Brown f[1] = y[0] - y[1]*y[1]; 510c4762a1bSJed Brown f[2] = y[1]*y[1]; 5119566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 5129566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 513c4762a1bSJed Brown PetscFunctionReturn(0); 514c4762a1bSJed Brown } 515c4762a1bSJed Brown 516c4762a1bSJed Brown PetscErrorCode IFunction_Hull1972B3(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) 517c4762a1bSJed Brown { 518c4762a1bSJed Brown PetscScalar *f; 519c4762a1bSJed Brown const PetscScalar *y; 520c4762a1bSJed Brown 521c4762a1bSJed Brown PetscFunctionBegin; 5229566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 5239566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 524c4762a1bSJed Brown f[0] = -y[0]; 525c4762a1bSJed Brown f[1] = y[0] - y[1]*y[1]; 526c4762a1bSJed Brown f[2] = y[1]*y[1]; 5279566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 5289566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 529c4762a1bSJed Brown /* Left hand side = ydot - f(y) */ 5309566063dSJacob Faibussowitsch PetscCall(VecAYPX(F,-1.0,Ydot)); 531c4762a1bSJed Brown PetscFunctionReturn(0); 532c4762a1bSJed Brown } 533c4762a1bSJed Brown 534c4762a1bSJed Brown PetscErrorCode IJacobian_Hull1972B3(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) 535c4762a1bSJed Brown { 536c4762a1bSJed Brown const PetscScalar *y; 537c4762a1bSJed Brown PetscInt row[3] = {0,1,2}; 538c4762a1bSJed Brown PetscScalar value[3][3]; 539c4762a1bSJed Brown 540c4762a1bSJed Brown PetscFunctionBegin; 5419566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 542c4762a1bSJed Brown value[0][0] = a + 1.0; value[0][1] = 0; value[0][2] = 0; 543c4762a1bSJed Brown value[1][0] = -1.0; value[1][1] = a + 2.0*y[1]; value[1][2] = 0; 544c4762a1bSJed Brown value[2][0] = 0; value[2][1] = -2.0*y[1]; value[2][2] = a; 5459566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,3,&row[0],3,&row[0],&value[0][0],INSERT_VALUES)); 5469566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 5479566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY)); 5489566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 549c4762a1bSJed Brown PetscFunctionReturn(0); 550c4762a1bSJed Brown } 551c4762a1bSJed Brown 552c4762a1bSJed Brown /* Hull, 1972, Problem B4 */ 553c4762a1bSJed Brown 554c4762a1bSJed Brown PetscErrorCode RHSFunction_Hull1972B4(TS ts, PetscReal t, Vec Y, Vec F, void *s) 555c4762a1bSJed Brown { 556c4762a1bSJed Brown PetscScalar *f; 557c4762a1bSJed Brown const PetscScalar *y; 558c4762a1bSJed Brown 559c4762a1bSJed Brown PetscFunctionBegin; 5609566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 5619566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 562c4762a1bSJed Brown f[0] = -y[1] - y[0]*y[2]/PetscSqrtScalar(y[0]*y[0]+y[1]*y[1]); 563c4762a1bSJed Brown f[1] = y[0] - y[1]*y[2]/PetscSqrtScalar(y[0]*y[0]+y[1]*y[1]); 564c4762a1bSJed Brown f[2] = y[0]/PetscSqrtScalar(y[0]*y[0]+y[1]*y[1]); 5659566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 5669566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 567c4762a1bSJed Brown PetscFunctionReturn(0); 568c4762a1bSJed Brown } 569c4762a1bSJed Brown 570c4762a1bSJed Brown PetscErrorCode IFunction_Hull1972B4(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) 571c4762a1bSJed Brown { 572c4762a1bSJed Brown PetscScalar *f; 573c4762a1bSJed Brown const PetscScalar *y; 574c4762a1bSJed Brown 575c4762a1bSJed Brown PetscFunctionBegin; 5769566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 5779566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 578c4762a1bSJed Brown f[0] = -y[1] - y[0]*y[2]/PetscSqrtScalar(y[0]*y[0]+y[1]*y[1]); 579c4762a1bSJed Brown f[1] = y[0] - y[1]*y[2]/PetscSqrtScalar(y[0]*y[0]+y[1]*y[1]); 580c4762a1bSJed Brown f[2] = y[0]/PetscSqrtScalar(y[0]*y[0]+y[1]*y[1]); 5819566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 5829566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 583c4762a1bSJed Brown /* Left hand side = ydot - f(y) */ 5849566063dSJacob Faibussowitsch PetscCall(VecAYPX(F,-1.0,Ydot)); 585c4762a1bSJed Brown PetscFunctionReturn(0); 586c4762a1bSJed Brown } 587c4762a1bSJed Brown 588c4762a1bSJed Brown PetscErrorCode IJacobian_Hull1972B4(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) 589c4762a1bSJed Brown { 590c4762a1bSJed Brown const PetscScalar *y; 591c4762a1bSJed Brown PetscInt row[3] = {0,1,2}; 592c4762a1bSJed Brown PetscScalar value[3][3],fac,fac2; 593c4762a1bSJed Brown 594c4762a1bSJed Brown PetscFunctionBegin; 5959566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 596c4762a1bSJed Brown fac = PetscPowScalar(y[0]*y[0]+y[1]*y[1],-1.5); 597c4762a1bSJed Brown fac2 = PetscPowScalar(y[0]*y[0]+y[1]*y[1],-0.5); 598c4762a1bSJed Brown value[0][0] = a + (y[1]*y[1]*y[2])*fac; 599c4762a1bSJed Brown value[0][1] = 1.0 - (y[0]*y[1]*y[2])*fac; 600c4762a1bSJed Brown value[0][2] = y[0]*fac2; 601c4762a1bSJed Brown value[1][0] = -1.0 - y[0]*y[1]*y[2]*fac; 602c4762a1bSJed Brown value[1][1] = a + y[0]*y[0]*y[2]*fac; 603c4762a1bSJed Brown value[1][2] = y[1]*fac2; 604c4762a1bSJed Brown value[2][0] = -y[1]*y[1]*fac; 605c4762a1bSJed Brown value[2][1] = y[0]*y[1]*fac; 606c4762a1bSJed Brown value[2][2] = a; 6079566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,3,&row[0],3,&row[0],&value[0][0],INSERT_VALUES)); 6089566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 6099566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY)); 6109566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 611c4762a1bSJed Brown PetscFunctionReturn(0); 612c4762a1bSJed Brown } 613c4762a1bSJed Brown 614c4762a1bSJed Brown /* Hull, 1972, Problem B5 */ 615c4762a1bSJed Brown 616c4762a1bSJed Brown PetscErrorCode RHSFunction_Hull1972B5(TS ts, PetscReal t, Vec Y, Vec F, void *s) 617c4762a1bSJed Brown { 618c4762a1bSJed Brown PetscScalar *f; 619c4762a1bSJed Brown const PetscScalar *y; 620c4762a1bSJed Brown 621c4762a1bSJed Brown PetscFunctionBegin; 6229566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 6239566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 624c4762a1bSJed Brown f[0] = y[1]*y[2]; 625c4762a1bSJed Brown f[1] = -y[0]*y[2]; 626c4762a1bSJed Brown f[2] = -0.51*y[0]*y[1]; 6279566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 6289566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 629c4762a1bSJed Brown PetscFunctionReturn(0); 630c4762a1bSJed Brown } 631c4762a1bSJed Brown 632c4762a1bSJed Brown PetscErrorCode IFunction_Hull1972B5(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) 633c4762a1bSJed Brown { 634c4762a1bSJed Brown PetscScalar *f; 635c4762a1bSJed Brown const PetscScalar *y; 636c4762a1bSJed Brown 637c4762a1bSJed Brown PetscFunctionBegin; 6389566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 6399566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 640c4762a1bSJed Brown f[0] = y[1]*y[2]; 641c4762a1bSJed Brown f[1] = -y[0]*y[2]; 642c4762a1bSJed Brown f[2] = -0.51*y[0]*y[1]; 6439566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 6449566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 645c4762a1bSJed Brown /* Left hand side = ydot - f(y) */ 6469566063dSJacob Faibussowitsch PetscCall(VecAYPX(F,-1.0,Ydot)); 647c4762a1bSJed Brown PetscFunctionReturn(0); 648c4762a1bSJed Brown } 649c4762a1bSJed Brown 650c4762a1bSJed Brown PetscErrorCode IJacobian_Hull1972B5(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) 651c4762a1bSJed Brown { 652c4762a1bSJed Brown const PetscScalar *y; 653c4762a1bSJed Brown PetscInt row[3] = {0,1,2}; 654c4762a1bSJed Brown PetscScalar value[3][3]; 655c4762a1bSJed Brown 656c4762a1bSJed Brown PetscFunctionBegin; 6579566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 658c4762a1bSJed Brown value[0][0] = a; value[0][1] = -y[2]; value[0][2] = -y[1]; 659c4762a1bSJed Brown value[1][0] = y[2]; value[1][1] = a; value[1][2] = y[0]; 660c4762a1bSJed Brown value[2][0] = 0.51*y[1]; value[2][1] = 0.51*y[0]; value[2][2] = a; 6619566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,3,&row[0],3,&row[0],&value[0][0],INSERT_VALUES)); 6629566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 6639566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY)); 6649566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 665c4762a1bSJed Brown PetscFunctionReturn(0); 666c4762a1bSJed Brown } 667c4762a1bSJed Brown 668c4762a1bSJed Brown /* Kulikov, 2013, Problem I */ 669c4762a1bSJed Brown 670c4762a1bSJed Brown PetscErrorCode RHSFunction_Kulikov2013I(TS ts, PetscReal t, Vec Y, Vec F, void *s) 671c4762a1bSJed Brown { 672c4762a1bSJed Brown PetscScalar *f; 673c4762a1bSJed Brown const PetscScalar *y; 674c4762a1bSJed Brown 675c4762a1bSJed Brown PetscFunctionBegin; 6769566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 6779566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 678c4762a1bSJed Brown f[0] = 2.*t*PetscPowScalar(y[1],1./5.)*y[3]; 679c4762a1bSJed Brown f[1] = 10.*t*y[3]*PetscExpScalar(5.0*(y[2]-1.)); 680c4762a1bSJed Brown f[2] = 2.*t*y[3]; 681c4762a1bSJed Brown f[3] = -2.*t*PetscLogScalar(y[0]); 6829566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 6839566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 684c4762a1bSJed Brown PetscFunctionReturn(0); 685c4762a1bSJed Brown } 686c4762a1bSJed Brown 687c4762a1bSJed Brown PetscErrorCode RHSJacobian_Kulikov2013I(TS ts, PetscReal t, Vec Y, Mat A, Mat B, void *s) 688c4762a1bSJed Brown { 689c4762a1bSJed Brown const PetscScalar *y; 690c4762a1bSJed Brown PetscInt row[4] = {0,1,2,3}; 691c4762a1bSJed Brown PetscScalar value[4][4]; 692c4762a1bSJed Brown PetscScalar m1,m2; 693c4762a1bSJed Brown PetscFunctionBegin; 6949566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 695c4762a1bSJed Brown m1=(2.*t*y[3])/(5.*PetscPowScalar(y[1],4./5.)); 696c4762a1bSJed Brown m2=2.*t*PetscPowScalar(y[1],1./5.); 697c4762a1bSJed Brown value[0][0] = 0. ; value[0][1] = m1; value[0][2] = 0.; value[0][3] = m2; 698c4762a1bSJed Brown m1=50.*t*y[3]*PetscExpScalar(5.0*(y[2]-1.)); 699c4762a1bSJed Brown m2=10.*t*PetscExpScalar(5.0*(y[2]-1.)); 700c4762a1bSJed Brown value[1][0] = 0.; value[1][1] = 0. ; value[1][2] = m1; value[1][3] = m2; 701c4762a1bSJed Brown value[2][0] = 0.; value[2][1] = 0.; value[2][2] = 0.; value[2][3] = 2*t; 702c4762a1bSJed Brown value[3][0] = -2.*t/y[0];value[3][1] = 0.; value[3][2] = 0.; value[3][3] = 0.; 7039566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,4,&row[0],4,&row[0],&value[0][0],INSERT_VALUES)); 7049566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 7059566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY)); 7069566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 707c4762a1bSJed Brown PetscFunctionReturn(0); 708c4762a1bSJed Brown } 709c4762a1bSJed Brown 710c4762a1bSJed Brown PetscErrorCode IFunction_Kulikov2013I(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) 711c4762a1bSJed Brown { 712c4762a1bSJed Brown PetscScalar *f; 713c4762a1bSJed Brown const PetscScalar *y; 714c4762a1bSJed Brown 715c4762a1bSJed Brown PetscFunctionBegin; 7169566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 7179566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 718c4762a1bSJed Brown f[0] = 2.*t*PetscPowScalar(y[1],1./5.)*y[3]; 719c4762a1bSJed Brown f[1] = 10.*t*y[3]*PetscExpScalar(5.0*(y[2]-1.)); 720c4762a1bSJed Brown f[2] = 2.*t*y[3]; 721c4762a1bSJed Brown f[3] = -2.*t*PetscLogScalar(y[0]); 7229566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 7239566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 724c4762a1bSJed Brown /* Left hand side = ydot - f(y) */ 7259566063dSJacob Faibussowitsch PetscCall(VecAYPX(F,-1.0,Ydot)); 726c4762a1bSJed Brown PetscFunctionReturn(0); 727c4762a1bSJed Brown } 728c4762a1bSJed Brown 729c4762a1bSJed Brown PetscErrorCode IJacobian_Kulikov2013I(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) 730c4762a1bSJed Brown { 731c4762a1bSJed Brown const PetscScalar *y; 732c4762a1bSJed Brown PetscInt row[4] = {0,1,2,3}; 733c4762a1bSJed Brown PetscScalar value[4][4]; 734c4762a1bSJed Brown PetscScalar m1,m2; 735c4762a1bSJed Brown 736c4762a1bSJed Brown PetscFunctionBegin; 7379566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 738c4762a1bSJed Brown m1=(2.*t*y[3])/(5.*PetscPowScalar(y[1],4./5.)); 739c4762a1bSJed Brown m2=2.*t*PetscPowScalar(y[1],1./5.); 740c4762a1bSJed Brown value[0][0] = a ; value[0][1] = m1; value[0][2] = 0.; value[0][3] = m2; 741c4762a1bSJed Brown m1=50.*t*y[3]*PetscExpScalar(5.0*(y[2]-1.)); 742c4762a1bSJed Brown m2=10.*t*PetscExpScalar(5.0*(y[2]-1.)); 743c4762a1bSJed Brown value[1][0] = 0.; value[1][1] = a ; value[1][2] = m1; value[1][3] = m2; 744c4762a1bSJed Brown value[2][0] = 0.; value[2][1] = 0.; value[2][2] = a; value[2][3] = 2*t; 745c4762a1bSJed Brown value[3][0] = -2.*t/y[0];value[3][1] = 0.; value[3][2] = 0.; value[3][3] = a; 7469566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,4,&row[0],4,&row[0],&value[0][0],INSERT_VALUES)); 7479566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 7489566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY)); 7499566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 750c4762a1bSJed Brown PetscFunctionReturn(0); 751c4762a1bSJed Brown } 752c4762a1bSJed Brown 753c4762a1bSJed Brown /* Hull, 1972, Problem C1 */ 754c4762a1bSJed Brown 755c4762a1bSJed Brown PetscErrorCode RHSFunction_Hull1972C1(TS ts, PetscReal t, Vec Y, Vec F, void *s) 756c4762a1bSJed Brown { 757c4762a1bSJed Brown PetscScalar *f; 758c4762a1bSJed Brown const PetscScalar *y; 759c4762a1bSJed Brown PetscInt N,i; 760c4762a1bSJed Brown 761c4762a1bSJed Brown PetscFunctionBegin; 7629566063dSJacob Faibussowitsch PetscCall(VecGetSize (Y,&N)); 7639566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 7649566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 765c4762a1bSJed Brown f[0] = -y[0]; 766c4762a1bSJed Brown for (i = 1; i < N-1; i++) { 767c4762a1bSJed Brown f[i] = y[i-1] - y[i]; 768c4762a1bSJed Brown } 769c4762a1bSJed Brown f[N-1] = y[N-2]; 7709566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 7719566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 772c4762a1bSJed Brown PetscFunctionReturn(0); 773c4762a1bSJed Brown } 774c4762a1bSJed Brown 775c4762a1bSJed Brown PetscErrorCode IFunction_Hull1972C1(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) 776c4762a1bSJed Brown { 777c4762a1bSJed Brown PetscScalar *f; 778c4762a1bSJed Brown const PetscScalar *y; 779c4762a1bSJed Brown PetscInt N,i; 780c4762a1bSJed Brown 781c4762a1bSJed Brown PetscFunctionBegin; 7829566063dSJacob Faibussowitsch PetscCall(VecGetSize (Y,&N)); 7839566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 7849566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 785c4762a1bSJed Brown f[0] = -y[0]; 786c4762a1bSJed Brown for (i = 1; i < N-1; i++) { 787c4762a1bSJed Brown f[i] = y[i-1] - y[i]; 788c4762a1bSJed Brown } 789c4762a1bSJed Brown f[N-1] = y[N-2]; 7909566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 7919566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 792c4762a1bSJed Brown /* Left hand side = ydot - f(y) */ 7939566063dSJacob Faibussowitsch PetscCall(VecAYPX(F,-1.0,Ydot)); 794c4762a1bSJed Brown PetscFunctionReturn(0); 795c4762a1bSJed Brown } 796c4762a1bSJed Brown 797c4762a1bSJed Brown PetscErrorCode IJacobian_Hull1972C1(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) 798c4762a1bSJed Brown { 799c4762a1bSJed Brown const PetscScalar *y; 800c4762a1bSJed Brown PetscInt N,i,col[2]; 801c4762a1bSJed Brown PetscScalar value[2]; 802c4762a1bSJed Brown 803c4762a1bSJed Brown PetscFunctionBegin; 8049566063dSJacob Faibussowitsch PetscCall(VecGetSize (Y,&N)); 8059566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 806c4762a1bSJed Brown i = 0; 807c4762a1bSJed Brown value[0] = a+1; col[0] = 0; 808c4762a1bSJed Brown value[1] = 0; col[1] = 1; 8099566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,2,col,value,INSERT_VALUES)); 810c4762a1bSJed Brown for (i = 0; i < N; i++) { 811c4762a1bSJed Brown value[0] = -1; col[0] = i-1; 812c4762a1bSJed Brown value[1] = a+1; col[1] = i; 8139566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,2,col,value,INSERT_VALUES)); 814c4762a1bSJed Brown } 815c4762a1bSJed Brown i = N-1; 816c4762a1bSJed Brown value[0] = -1; col[0] = N-2; 817c4762a1bSJed Brown value[1] = a; col[1] = N-1; 8189566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,2,col,value,INSERT_VALUES)); 8199566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 8209566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY)); 8219566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 822c4762a1bSJed Brown PetscFunctionReturn(0); 823c4762a1bSJed Brown } 824c4762a1bSJed Brown 825c4762a1bSJed Brown /* Hull, 1972, Problem C2 */ 826c4762a1bSJed Brown 827c4762a1bSJed Brown PetscErrorCode RHSFunction_Hull1972C2(TS ts, PetscReal t, Vec Y, Vec F, void *s) 828c4762a1bSJed Brown { 829c4762a1bSJed Brown const PetscScalar *y; 830c4762a1bSJed Brown PetscScalar *f; 831c4762a1bSJed Brown PetscInt N,i; 832c4762a1bSJed Brown 833c4762a1bSJed Brown PetscFunctionBegin; 8349566063dSJacob Faibussowitsch PetscCall(VecGetSize (Y,&N)); 8359566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 8369566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 837c4762a1bSJed Brown f[0] = -y[0]; 838c4762a1bSJed Brown for (i = 1; i < N-1; i++) { 839c4762a1bSJed Brown f[i] = (PetscReal)i*y[i-1] - (PetscReal)(i+1)*y[i]; 840c4762a1bSJed Brown } 841c4762a1bSJed Brown f[N-1] = (PetscReal)(N-1)*y[N-2]; 8429566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 8439566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 844c4762a1bSJed Brown PetscFunctionReturn(0); 845c4762a1bSJed Brown } 846c4762a1bSJed Brown 847c4762a1bSJed Brown PetscErrorCode IFunction_Hull1972C2(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) 848c4762a1bSJed Brown { 849c4762a1bSJed Brown PetscScalar *f; 850c4762a1bSJed Brown const PetscScalar *y; 851c4762a1bSJed Brown PetscInt N,i; 852c4762a1bSJed Brown 853c4762a1bSJed Brown PetscFunctionBegin; 8549566063dSJacob Faibussowitsch PetscCall(VecGetSize (Y,&N)); 8559566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 8569566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 857c4762a1bSJed Brown f[0] = -y[0]; 858c4762a1bSJed Brown for (i = 1; i < N-1; i++) { 859c4762a1bSJed Brown f[i] = (PetscReal)i*y[i-1] - (PetscReal)(i+1)*y[i]; 860c4762a1bSJed Brown } 861c4762a1bSJed Brown f[N-1] = (PetscReal)(N-1)*y[N-2]; 8629566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 8639566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 864c4762a1bSJed Brown /* Left hand side = ydot - f(y) */ 8659566063dSJacob Faibussowitsch PetscCall(VecAYPX(F,-1.0,Ydot)); 866c4762a1bSJed Brown PetscFunctionReturn(0); 867c4762a1bSJed Brown } 868c4762a1bSJed Brown 869c4762a1bSJed Brown PetscErrorCode IJacobian_Hull1972C2(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) 870c4762a1bSJed Brown { 871c4762a1bSJed Brown const PetscScalar *y; 872c4762a1bSJed Brown PetscInt N,i,col[2]; 873c4762a1bSJed Brown PetscScalar value[2]; 874c4762a1bSJed Brown 875c4762a1bSJed Brown PetscFunctionBegin; 8769566063dSJacob Faibussowitsch PetscCall(VecGetSize (Y,&N)); 8779566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 878c4762a1bSJed Brown i = 0; 879c4762a1bSJed Brown value[0] = a+1; col[0] = 0; 880c4762a1bSJed Brown value[1] = 0; col[1] = 1; 8819566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,2,col,value,INSERT_VALUES)); 882c4762a1bSJed Brown for (i = 0; i < N; i++) { 883c4762a1bSJed Brown value[0] = -(PetscReal) i; col[0] = i-1; 884c4762a1bSJed Brown value[1] = a+(PetscReal)(i+1); col[1] = i; 8859566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,2,col,value,INSERT_VALUES)); 886c4762a1bSJed Brown } 887c4762a1bSJed Brown i = N-1; 888c4762a1bSJed Brown value[0] = -(PetscReal) (N-1); col[0] = N-2; 889c4762a1bSJed Brown value[1] = a; col[1] = N-1; 8909566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,2,col,value,INSERT_VALUES)); 8919566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 8929566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY)); 8939566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 894c4762a1bSJed Brown PetscFunctionReturn(0); 895c4762a1bSJed Brown } 896c4762a1bSJed Brown 897c4762a1bSJed Brown /* Hull, 1972, Problem C3 and C4 */ 898c4762a1bSJed Brown 899c4762a1bSJed Brown PetscErrorCode RHSFunction_Hull1972C34(TS ts, PetscReal t, Vec Y, Vec F, void *s) 900c4762a1bSJed Brown { 901c4762a1bSJed Brown PetscScalar *f; 902c4762a1bSJed Brown const PetscScalar *y; 903c4762a1bSJed Brown PetscInt N,i; 904c4762a1bSJed Brown 905c4762a1bSJed Brown PetscFunctionBegin; 9069566063dSJacob Faibussowitsch PetscCall(VecGetSize (Y,&N)); 9079566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 9089566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 909c4762a1bSJed Brown f[0] = -2.0*y[0] + y[1]; 910c4762a1bSJed Brown for (i = 1; i < N-1; i++) { 911c4762a1bSJed Brown f[i] = y[i-1] - 2.0*y[i] + y[i+1]; 912c4762a1bSJed Brown } 913c4762a1bSJed Brown f[N-1] = y[N-2] - 2.0*y[N-1]; 9149566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 9159566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 916c4762a1bSJed Brown PetscFunctionReturn(0); 917c4762a1bSJed Brown } 918c4762a1bSJed Brown 919c4762a1bSJed Brown PetscErrorCode IFunction_Hull1972C34(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) 920c4762a1bSJed Brown { 921c4762a1bSJed Brown PetscScalar *f; 922c4762a1bSJed Brown const PetscScalar *y; 923c4762a1bSJed Brown PetscInt N,i; 924c4762a1bSJed Brown 925c4762a1bSJed Brown PetscFunctionBegin; 9269566063dSJacob Faibussowitsch PetscCall(VecGetSize (Y,&N)); 9279566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 9289566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 929c4762a1bSJed Brown f[0] = -2.0*y[0] + y[1]; 930c4762a1bSJed Brown for (i = 1; i < N-1; i++) { 931c4762a1bSJed Brown f[i] = y[i-1] - 2.0*y[i] + y[i+1]; 932c4762a1bSJed Brown } 933c4762a1bSJed Brown f[N-1] = y[N-2] - 2.0*y[N-1]; 9349566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 9359566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 936c4762a1bSJed Brown /* Left hand side = ydot - f(y) */ 9379566063dSJacob Faibussowitsch PetscCall(VecAYPX(F,-1.0,Ydot)); 938c4762a1bSJed Brown PetscFunctionReturn(0); 939c4762a1bSJed Brown } 940c4762a1bSJed Brown 941c4762a1bSJed Brown PetscErrorCode IJacobian_Hull1972C34(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) 942c4762a1bSJed Brown { 943c4762a1bSJed Brown const PetscScalar *y; 944c4762a1bSJed Brown PetscScalar value[3]; 945c4762a1bSJed Brown PetscInt N,i,col[3]; 946c4762a1bSJed Brown 947c4762a1bSJed Brown PetscFunctionBegin; 9489566063dSJacob Faibussowitsch PetscCall(VecGetSize (Y,&N)); 9499566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Y,&y)); 950c4762a1bSJed Brown for (i = 0; i < N; i++) { 951c4762a1bSJed Brown if (i == 0) { 952c4762a1bSJed Brown value[0] = a+2; col[0] = i; 953c4762a1bSJed Brown value[1] = -1; col[1] = i+1; 954c4762a1bSJed Brown value[2] = 0; col[2] = i+2; 955c4762a1bSJed Brown } else if (i == N-1) { 956c4762a1bSJed Brown value[0] = 0; col[0] = i-2; 957c4762a1bSJed Brown value[1] = -1; col[1] = i-1; 958c4762a1bSJed Brown value[2] = a+2; col[2] = i; 959c4762a1bSJed Brown } else { 960c4762a1bSJed Brown value[0] = -1; col[0] = i-1; 961c4762a1bSJed Brown value[1] = a+2; col[1] = i; 962c4762a1bSJed Brown value[2] = -1; col[2] = i+1; 963c4762a1bSJed Brown } 9649566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,3,col,value,INSERT_VALUES)); 965c4762a1bSJed Brown } 9669566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 9679566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY)); 9689566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Y,&y)); 969c4762a1bSJed Brown PetscFunctionReturn(0); 970c4762a1bSJed Brown } 971c4762a1bSJed Brown 972c4762a1bSJed Brown /***************************************************************************/ 973c4762a1bSJed Brown 974c4762a1bSJed Brown /* Sets the initial solution for the IVP and sets up the function pointers*/ 975c4762a1bSJed Brown PetscErrorCode Initialize(Vec Y, void* s) 976c4762a1bSJed Brown { 977c4762a1bSJed Brown char *p = (char*) s; 978c4762a1bSJed Brown PetscScalar *y; 979c4762a1bSJed Brown PetscReal t0; 980c4762a1bSJed Brown PetscInt N = GetSize((const char *)s); 981c4762a1bSJed Brown PetscBool flg; 982c4762a1bSJed Brown 983c4762a1bSJed Brown PetscFunctionBegin; 984c4762a1bSJed Brown VecZeroEntries(Y); 9859566063dSJacob Faibussowitsch PetscCall(VecGetArray(Y,&y)); 986c4762a1bSJed Brown if (!strcmp(p,"hull1972a1")) { 987c4762a1bSJed Brown y[0] = 1.0; 988c4762a1bSJed Brown RHSFunction = RHSFunction_Hull1972A1; 989c4762a1bSJed Brown RHSJacobian = RHSJacobian_Hull1972A1; 990c4762a1bSJed Brown IFunction = IFunction_Hull1972A1; 991c4762a1bSJed Brown IJacobian = IJacobian_Hull1972A1; 992c4762a1bSJed Brown } else if (!strcmp(p,"hull1972a2")) { 993c4762a1bSJed Brown y[0] = 1.0; 994c4762a1bSJed Brown RHSFunction = RHSFunction_Hull1972A2; 995c4762a1bSJed Brown RHSJacobian = RHSJacobian_Hull1972A2; 996c4762a1bSJed Brown IFunction = IFunction_Hull1972A2; 997c4762a1bSJed Brown IJacobian = IJacobian_Hull1972A2; 998c4762a1bSJed Brown } else if (!strcmp(p,"hull1972a3")) { 999c4762a1bSJed Brown y[0] = 1.0; 1000c4762a1bSJed Brown RHSFunction = RHSFunction_Hull1972A3; 1001c4762a1bSJed Brown RHSJacobian = RHSJacobian_Hull1972A3; 1002c4762a1bSJed Brown IFunction = IFunction_Hull1972A3; 1003c4762a1bSJed Brown IJacobian = IJacobian_Hull1972A3; 1004c4762a1bSJed Brown } else if (!strcmp(p,"hull1972a4")) { 1005c4762a1bSJed Brown y[0] = 1.0; 1006c4762a1bSJed Brown RHSFunction = RHSFunction_Hull1972A4; 1007c4762a1bSJed Brown RHSJacobian = RHSJacobian_Hull1972A4; 1008c4762a1bSJed Brown IFunction = IFunction_Hull1972A4; 1009c4762a1bSJed Brown IJacobian = IJacobian_Hull1972A4; 1010c4762a1bSJed Brown } else if (!strcmp(p,"hull1972a5")) { 1011c4762a1bSJed Brown y[0] = 4.0; 1012c4762a1bSJed Brown RHSFunction = RHSFunction_Hull1972A5; 1013c4762a1bSJed Brown RHSJacobian = RHSJacobian_Hull1972A5; 1014c4762a1bSJed Brown IFunction = IFunction_Hull1972A5; 1015c4762a1bSJed Brown IJacobian = IJacobian_Hull1972A5; 1016c4762a1bSJed Brown } else if (!strcmp(p,"hull1972b1")) { 1017c4762a1bSJed Brown y[0] = 1.0; 1018c4762a1bSJed Brown y[1] = 3.0; 1019c4762a1bSJed Brown RHSFunction = RHSFunction_Hull1972B1; 1020c4762a1bSJed Brown IFunction = IFunction_Hull1972B1; 1021c4762a1bSJed Brown IJacobian = IJacobian_Hull1972B1; 1022c4762a1bSJed Brown } else if (!strcmp(p,"hull1972b2")) { 1023c4762a1bSJed Brown y[0] = 2.0; 1024c4762a1bSJed Brown y[1] = 0.0; 1025c4762a1bSJed Brown y[2] = 1.0; 1026c4762a1bSJed Brown RHSFunction = RHSFunction_Hull1972B2; 1027c4762a1bSJed Brown IFunction = IFunction_Hull1972B2; 1028c4762a1bSJed Brown IJacobian = IJacobian_Hull1972B2; 1029c4762a1bSJed Brown } else if (!strcmp(p,"hull1972b3")) { 1030c4762a1bSJed Brown y[0] = 1.0; 1031c4762a1bSJed Brown y[1] = 0.0; 1032c4762a1bSJed Brown y[2] = 0.0; 1033c4762a1bSJed Brown RHSFunction = RHSFunction_Hull1972B3; 1034c4762a1bSJed Brown IFunction = IFunction_Hull1972B3; 1035c4762a1bSJed Brown IJacobian = IJacobian_Hull1972B3; 1036c4762a1bSJed Brown } else if (!strcmp(p,"hull1972b4")) { 1037c4762a1bSJed Brown y[0] = 3.0; 1038c4762a1bSJed Brown y[1] = 0.0; 1039c4762a1bSJed Brown y[2] = 0.0; 1040c4762a1bSJed Brown RHSFunction = RHSFunction_Hull1972B4; 1041c4762a1bSJed Brown IFunction = IFunction_Hull1972B4; 1042c4762a1bSJed Brown IJacobian = IJacobian_Hull1972B4; 1043c4762a1bSJed Brown } else if (!strcmp(p,"hull1972b5")) { 1044c4762a1bSJed Brown y[0] = 0.0; 1045c4762a1bSJed Brown y[1] = 1.0; 1046c4762a1bSJed Brown y[2] = 1.0; 1047c4762a1bSJed Brown RHSFunction = RHSFunction_Hull1972B5; 1048c4762a1bSJed Brown IFunction = IFunction_Hull1972B5; 1049c4762a1bSJed Brown IJacobian = IJacobian_Hull1972B5; 1050c4762a1bSJed Brown } else if (!strcmp(p,"kulik2013i")) { 1051c4762a1bSJed Brown t0=0.; 1052c4762a1bSJed Brown y[0] = PetscExpReal(PetscSinReal(t0*t0)); 1053c4762a1bSJed Brown y[1] = PetscExpReal(5.*PetscSinReal(t0*t0)); 1054c4762a1bSJed Brown y[2] = PetscSinReal(t0*t0)+1.0; 1055c4762a1bSJed Brown y[3] = PetscCosReal(t0*t0); 1056c4762a1bSJed Brown RHSFunction = RHSFunction_Kulikov2013I; 1057c4762a1bSJed Brown RHSJacobian = RHSJacobian_Kulikov2013I; 1058c4762a1bSJed Brown IFunction = IFunction_Kulikov2013I; 1059c4762a1bSJed Brown IJacobian = IJacobian_Kulikov2013I; 1060c4762a1bSJed Brown } else if (!strcmp(p,"hull1972c1")) { 1061c4762a1bSJed Brown y[0] = 1.0; 1062c4762a1bSJed Brown RHSFunction = RHSFunction_Hull1972C1; 1063c4762a1bSJed Brown IFunction = IFunction_Hull1972C1; 1064c4762a1bSJed Brown IJacobian = IJacobian_Hull1972C1; 1065c4762a1bSJed Brown } else if (!strcmp(p,"hull1972c2")) { 1066c4762a1bSJed Brown y[0] = 1.0; 1067c4762a1bSJed Brown RHSFunction = RHSFunction_Hull1972C2; 1068c4762a1bSJed Brown IFunction = IFunction_Hull1972C2; 1069c4762a1bSJed Brown IJacobian = IJacobian_Hull1972C2; 1070c4762a1bSJed Brown } else if ((!strcmp(p,"hull1972c3")) 1071c4762a1bSJed Brown ||(!strcmp(p,"hull1972c4"))) { 1072c4762a1bSJed Brown y[0] = 1.0; 1073c4762a1bSJed Brown RHSFunction = RHSFunction_Hull1972C34; 1074c4762a1bSJed Brown IFunction = IFunction_Hull1972C34; 1075c4762a1bSJed Brown IJacobian = IJacobian_Hull1972C34; 1076c4762a1bSJed Brown } 10779566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalarArray(NULL,NULL,"-yinit",y,&N,&flg)); 107863a3b9bcSJacob Faibussowitsch PetscCheck((N == GetSize((const char*)s)) || !flg,PETSC_COMM_WORLD,PETSC_ERR_ARG_SIZ,"Number of initial values %" PetscInt_FMT " does not match problem size %" PetscInt_FMT ".",N,GetSize((const char*)s)); 10799566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Y,&y)); 1080c4762a1bSJed Brown PetscFunctionReturn(0); 1081c4762a1bSJed Brown } 1082c4762a1bSJed Brown 1083c4762a1bSJed Brown /* Calculates the exact solution to problems that have one */ 1084c4762a1bSJed Brown PetscErrorCode ExactSolution(Vec Y, void* s, PetscReal t, PetscBool *flag) 1085c4762a1bSJed Brown { 1086c4762a1bSJed Brown char *p = (char*) s; 1087c4762a1bSJed Brown PetscScalar *y; 1088c4762a1bSJed Brown 1089c4762a1bSJed Brown PetscFunctionBegin; 1090c4762a1bSJed Brown if (!strcmp(p,"hull1972a1")) { 10919566063dSJacob Faibussowitsch PetscCall(VecGetArray(Y,&y)); 1092c4762a1bSJed Brown y[0] = PetscExpReal(-t); 1093c4762a1bSJed Brown *flag = PETSC_TRUE; 10949566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Y,&y)); 1095c4762a1bSJed Brown } else if (!strcmp(p,"hull1972a2")) { 10969566063dSJacob Faibussowitsch PetscCall(VecGetArray(Y,&y)); 1097c4762a1bSJed Brown y[0] = 1.0/PetscSqrtReal(t+1); 1098c4762a1bSJed Brown *flag = PETSC_TRUE; 10999566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Y,&y)); 1100c4762a1bSJed Brown } else if (!strcmp(p,"hull1972a3")) { 11019566063dSJacob Faibussowitsch PetscCall(VecGetArray(Y,&y)); 1102c4762a1bSJed Brown y[0] = PetscExpReal(PetscSinReal(t)); 1103c4762a1bSJed Brown *flag = PETSC_TRUE; 11049566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Y,&y)); 1105c4762a1bSJed Brown } else if (!strcmp(p,"hull1972a4")) { 11069566063dSJacob Faibussowitsch PetscCall(VecGetArray(Y,&y)); 1107c4762a1bSJed Brown y[0] = 20.0/(1+19.0*PetscExpReal(-t/4.0)); 1108c4762a1bSJed Brown *flag = PETSC_TRUE; 11099566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Y,&y)); 1110c4762a1bSJed Brown } else if (!strcmp(p,"kulik2013i")) { 11119566063dSJacob Faibussowitsch PetscCall(VecGetArray(Y,&y)); 1112c4762a1bSJed Brown y[0] = PetscExpReal(PetscSinReal(t*t)); 1113c4762a1bSJed Brown y[1] = PetscExpReal(5.*PetscSinReal(t*t)); 1114c4762a1bSJed Brown y[2] = PetscSinReal(t*t)+1.0; 1115c4762a1bSJed Brown y[3] = PetscCosReal(t*t); 1116c4762a1bSJed Brown *flag = PETSC_TRUE; 11179566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Y,&y)); 1118c4762a1bSJed Brown } else { 11199566063dSJacob Faibussowitsch PetscCall(VecSet(Y,0)); 1120c4762a1bSJed Brown *flag = PETSC_FALSE; 1121c4762a1bSJed Brown } 1122c4762a1bSJed Brown PetscFunctionReturn(0); 1123c4762a1bSJed Brown } 1124c4762a1bSJed Brown 1125c4762a1bSJed Brown /* Solves the specified ODE and computes the error if exact solution is available */ 1126c4762a1bSJed Brown PetscErrorCode SolveODE(char* ptype, PetscReal dt, PetscReal tfinal, PetscInt maxiter, PetscReal *error, PetscBool *exact_flag) 1127c4762a1bSJed Brown { 1128c4762a1bSJed Brown TS ts; /* time-integrator */ 1129c4762a1bSJed Brown Vec Y; /* Solution vector */ 1130c4762a1bSJed Brown Vec Yex; /* Exact solution */ 1131c4762a1bSJed Brown PetscInt N; /* Size of the system of equations */ 1132c4762a1bSJed Brown TSType time_scheme; /* Type of time-integration scheme */ 1133c4762a1bSJed Brown Mat Jac = NULL; /* Jacobian matrix */ 1134c4762a1bSJed Brown Vec Yerr; /* Auxiliary solution vector */ 1135c4762a1bSJed Brown PetscReal err_norm; /* Estimated error norm */ 1136c4762a1bSJed Brown PetscReal final_time; /* Actual final time from the integrator */ 1137c4762a1bSJed Brown 1138c4762a1bSJed Brown PetscFunctionBegin; 1139c4762a1bSJed Brown N = GetSize((const char *)&ptype[0]); 11403c633725SBarry Smith PetscCheck(N >= 0,PETSC_COMM_WORLD,PETSC_ERR_ARG_SIZ,"Illegal problem specification."); 11419566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&Y)); 11429566063dSJacob Faibussowitsch PetscCall(VecSetSizes(Y,N,PETSC_DECIDE)); 11439566063dSJacob Faibussowitsch PetscCall(VecSetUp(Y)); 11449566063dSJacob Faibussowitsch PetscCall(VecSet(Y,0)); 1145c4762a1bSJed Brown 1146c4762a1bSJed Brown /* Initialize the problem */ 11479566063dSJacob Faibussowitsch PetscCall(Initialize(Y,&ptype[0])); 1148c4762a1bSJed Brown 1149c4762a1bSJed Brown /* Create and initialize the time-integrator */ 11509566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD,&ts)); 1151c4762a1bSJed Brown /* Default time integration options */ 11529566063dSJacob Faibussowitsch PetscCall(TSSetType(ts,TSRK)); 11539566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(ts,maxiter)); 11549566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts,tfinal)); 11559566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts,dt)); 11569566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP)); 1157c4762a1bSJed Brown /* Read command line options for time integration */ 11589566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 1159c4762a1bSJed Brown /* Set solution vector */ 11609566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts,Y)); 1161c4762a1bSJed Brown /* Specify left/right-hand side functions */ 11629566063dSJacob Faibussowitsch PetscCall(TSGetType(ts,&time_scheme)); 1163c4762a1bSJed Brown 1164c4762a1bSJed Brown if ((!strcmp(time_scheme,TSEULER)) || (!strcmp(time_scheme,TSRK)) || (!strcmp(time_scheme,TSSSP) || (!strcmp(time_scheme,TSGLEE)))) { 1165c4762a1bSJed Brown /* Explicit time-integration -> specify right-hand side function ydot = f(y) */ 11669566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts,NULL,RHSFunction,&ptype[0])); 11679566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&Jac)); 11689566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Jac,PETSC_DECIDE,PETSC_DECIDE,N,N)); 11699566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(Jac)); 11709566063dSJacob Faibussowitsch PetscCall(MatSetUp(Jac)); 11719566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(ts,Jac,Jac,RHSJacobian,&ptype[0])); 1172c4762a1bSJed Brown } else if ((!strcmp(time_scheme,TSTHETA)) || (!strcmp(time_scheme,TSBEULER)) || (!strcmp(time_scheme,TSCN)) || (!strcmp(time_scheme,TSALPHA)) || (!strcmp(time_scheme,TSARKIMEX))) { 1173c4762a1bSJed Brown /* Implicit time-integration -> specify left-hand side function ydot-f(y) = 0 */ 1174c4762a1bSJed Brown /* and its Jacobian function */ 11759566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts,NULL,IFunction,&ptype[0])); 11769566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&Jac)); 11779566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Jac,PETSC_DECIDE,PETSC_DECIDE,N,N)); 11789566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(Jac)); 11799566063dSJacob Faibussowitsch PetscCall(MatSetUp(Jac)); 11809566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts,Jac,Jac,IJacobian,&ptype[0])); 1181c4762a1bSJed Brown } 1182c4762a1bSJed Brown 1183c4762a1bSJed Brown /* Solve */ 11849566063dSJacob Faibussowitsch PetscCall(TSSolve(ts,Y)); 11859566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts,&final_time)); 1186c4762a1bSJed Brown 1187c4762a1bSJed Brown /* Get the estimated error, if available */ 11889566063dSJacob Faibussowitsch PetscCall(VecDuplicate(Y,&Yerr)); 11899566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Yerr)); 11909566063dSJacob Faibussowitsch PetscCall(TSGetTimeError(ts,0,&Yerr)); 11919566063dSJacob Faibussowitsch PetscCall(VecNorm(Yerr,NORM_2,&err_norm)); 11929566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Yerr)); 119363a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Estimated Error = %e.\n",(double)err_norm)); 1194c4762a1bSJed Brown 1195c4762a1bSJed Brown /* Exact solution */ 11969566063dSJacob Faibussowitsch PetscCall(VecDuplicate(Y,&Yex)); 119763a3b9bcSJacob Faibussowitsch if (PetscAbsReal(final_time-tfinal)>2.*PETSC_MACHINE_EPSILON) { 11989566063dSJacob 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)); 1199c4762a1bSJed Brown } 12009566063dSJacob Faibussowitsch PetscCall(ExactSolution(Yex,&ptype[0],final_time,exact_flag)); 1201c4762a1bSJed Brown 1202c4762a1bSJed Brown /* Calculate Error */ 12039566063dSJacob Faibussowitsch PetscCall(VecAYPX(Yex,-1.0,Y)); 12049566063dSJacob Faibussowitsch PetscCall(VecNorm(Yex,NORM_2,error)); 1205c4762a1bSJed Brown *error = PetscSqrtReal(((*error)*(*error))/N); 1206c4762a1bSJed Brown 1207c4762a1bSJed Brown /* Clean up and finalize */ 12089566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Jac)); 12099566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 12109566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Yex)); 12119566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y)); 1212c4762a1bSJed Brown 1213c4762a1bSJed Brown PetscFunctionReturn(0); 1214c4762a1bSJed Brown } 1215c4762a1bSJed Brown 1216c4762a1bSJed Brown int main(int argc, char **argv) 1217c4762a1bSJed Brown { 1218c4762a1bSJed Brown char ptype[256] = "hull1972a1"; /* Problem specification */ 1219c4762a1bSJed Brown PetscInt n_refine = 1; /* Number of refinement levels for convergence analysis */ 1220c4762a1bSJed Brown PetscReal refine_fac = 2.0; /* Refinement factor for dt */ 1221c4762a1bSJed Brown PetscReal dt_initial = 0.01; /* Initial default value of dt */ 1222c4762a1bSJed Brown PetscReal dt; 1223c4762a1bSJed Brown PetscReal tfinal = 20.0; /* Final time for the time-integration */ 1224c4762a1bSJed Brown PetscInt maxiter = 100000; /* Maximum number of time-integration iterations */ 1225c4762a1bSJed Brown PetscReal *error; /* Array to store the errors for convergence analysis */ 1226c4762a1bSJed Brown PetscMPIInt size; /* No of processors */ 1227c4762a1bSJed Brown PetscBool flag; /* Flag denoting availability of exact solution */ 1228c4762a1bSJed Brown PetscInt r; 1229c4762a1bSJed Brown 1230c4762a1bSJed Brown /* Initialize program */ 1231*327415f7SBarry Smith PetscFunctionBeginUser; 12329566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 1233c4762a1bSJed Brown 1234c4762a1bSJed Brown /* Check if running with only 1 proc */ 12359566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 12363c633725SBarry Smith PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Only for sequential runs"); 1237c4762a1bSJed Brown 1238d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"ex31",NULL); 12399566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-problem","Problem specification","<hull1972a1>",ptype,ptype,sizeof(ptype),NULL)); 12409566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-refinement_levels","Number of refinement levels for convergence analysis","<1>",n_refine,&n_refine,NULL)); 12419566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-refinement_factor","Refinement factor for dt","<2.0>",refine_fac,&refine_fac,NULL)); 12429566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-dt","Time step size (for convergence analysis, initial time step)","<0.01>",dt_initial,&dt_initial,NULL)); 12439566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-final_time","Final time for the time-integration","<20.0>",tfinal,&tfinal,NULL)); 1244d0609cedSBarry Smith PetscOptionsEnd(); 1245c4762a1bSJed Brown 12469566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n_refine,&error)); 1247c4762a1bSJed Brown for (r = 0,dt = dt_initial; r < n_refine; r++) { 1248c4762a1bSJed Brown error[r] = 0; 1249c4762a1bSJed Brown if (r > 0) dt /= refine_fac; 1250c4762a1bSJed Brown 125163a3b9bcSJacob Faibussowitsch 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,GetSize(&ptype[0]))); 12529566063dSJacob Faibussowitsch PetscCall(SolveODE(&ptype[0],dt,tfinal,maxiter,&error[r],&flag)); 1253c4762a1bSJed Brown if (flag) { 1254c4762a1bSJed Brown /* If exact solution available for the specified ODE */ 1255c4762a1bSJed Brown if (r > 0) { 1256c4762a1bSJed Brown PetscReal conv_rate = (PetscLogReal(error[r]) - PetscLogReal(error[r-1])) / (-PetscLogReal(refine_fac)); 12579566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error = %E,\tConvergence rate = %f.\n",(double)error[r],(double)conv_rate)); 1258c4762a1bSJed Brown } else { 125963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error = %E.\n",(double)error[r])); 1260c4762a1bSJed Brown } 1261c4762a1bSJed Brown } 1262c4762a1bSJed Brown } 12639566063dSJacob Faibussowitsch PetscCall(PetscFree(error)); 12649566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 1265b122ec5aSJacob Faibussowitsch return 0; 1266c4762a1bSJed Brown } 1267c4762a1bSJed Brown 1268c4762a1bSJed Brown /*TEST 1269c4762a1bSJed Brown 1270c4762a1bSJed Brown test: 1271c4762a1bSJed Brown suffix: 2 1272c4762a1bSJed Brown args: -ts_type glee -final_time 5 -ts_adapt_type none 1273c4762a1bSJed Brown timeoutfactor: 3 1274c4762a1bSJed Brown requires: !single 1275c4762a1bSJed Brown 1276c4762a1bSJed Brown test: 1277c4762a1bSJed Brown suffix: 3 1278c4762a1bSJed 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 1279c4762a1bSJed Brown timeoutfactor: 3 1280c4762a1bSJed Brown requires: !single 1281c4762a1bSJed Brown 1282c4762a1bSJed Brown test: 1283c4762a1bSJed Brown suffix: 4 1284c4762a1bSJed 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 1285c4762a1bSJed Brown timeoutfactor: 3 1286c4762a1bSJed Brown requires: !single !__float128 1287c4762a1bSJed Brown 1288c4762a1bSJed Brown TEST*/ 1289