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 Concepts: TS 6c4762a1bSJed Brown Useful command line parameters: 7c4762a1bSJed Brown -problem <hull1972a1>: choose which problem to solve (see references 8c4762a1bSJed Brown for complete listing of problems). 9c4762a1bSJed Brown -ts_type <euler>: specify time-integrator 10c4762a1bSJed Brown -ts_adapt_type <basic>: specify time-step adapting (none,basic,advanced) 11c4762a1bSJed Brown -refinement_levels <1>: number of refinement levels for convergence analysis 12c4762a1bSJed Brown -refinement_factor <2.0>: factor to refine time step size by for convergence analysis 13c4762a1bSJed Brown -dt <0.01>: specify time step (initial time step for convergence analysis) 14c4762a1bSJed Brown 15c4762a1bSJed Brown */ 16c4762a1bSJed Brown 17c4762a1bSJed Brown /* 18c4762a1bSJed Brown List of cases and their names in the code:- 19c4762a1bSJed Brown From Hull, T.E., Enright, W.H., Fellen, B.M., and Sedgwick, A.E., 20c4762a1bSJed Brown "Comparing Numerical Methods for Ordinary Differential 21c4762a1bSJed Brown Equations", SIAM J. Numer. Anal., 9(4), 1972, pp. 603 - 635 22c4762a1bSJed Brown A1 -> "hull1972a1" (exact solution available) 23c4762a1bSJed Brown A2 -> "hull1972a2" (exact solution available) 24c4762a1bSJed Brown A3 -> "hull1972a3" (exact solution available) 25c4762a1bSJed Brown A4 -> "hull1972a4" (exact solution available) 26c4762a1bSJed Brown A5 -> "hull1972a5" 27c4762a1bSJed Brown B1 -> "hull1972b1" 28c4762a1bSJed Brown B2 -> "hull1972b2" 29c4762a1bSJed Brown B3 -> "hull1972b3" 30c4762a1bSJed Brown B4 -> "hull1972b4" 31c4762a1bSJed Brown B5 -> "hull1972b5" 32c4762a1bSJed Brown C1 -> "hull1972c1" 33c4762a1bSJed Brown C2 -> "hull1972c2" 34c4762a1bSJed Brown C3 -> "hull1972c3" 35c4762a1bSJed Brown C4 -> "hull1972c4" 36c4762a1bSJed Brown 37c4762a1bSJed Brown From Constantinescu, E. "Estimating Global Errors in Time Stepping" ArXiv e-prints, 38c4762a1bSJed Brown https://arxiv.org/abs/1503.05166, 2016 39c4762a1bSJed Brown 40c4762a1bSJed Brown Kulikov2013I -> "kulik2013i" 41c4762a1bSJed Brown 42c4762a1bSJed Brown */ 43c4762a1bSJed Brown 44c4762a1bSJed Brown #include <petscts.h> 45c4762a1bSJed Brown 46c4762a1bSJed Brown /* Function declarations */ 47c4762a1bSJed Brown PetscErrorCode (*RHSFunction) (TS,PetscReal,Vec,Vec,void*); 48c4762a1bSJed Brown PetscErrorCode (*RHSJacobian) (TS,PetscReal,Vec,Mat,Mat,void*); 49c4762a1bSJed Brown PetscErrorCode (*IFunction) (TS,PetscReal,Vec,Vec,Vec,void*); 50c4762a1bSJed Brown PetscErrorCode (*IJacobian) (TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*); 51c4762a1bSJed Brown 52c4762a1bSJed Brown /* Returns the size of the system of equations depending on problem specification */ 53c4762a1bSJed Brown PetscInt GetSize(const char *p) 54c4762a1bSJed Brown { 55c4762a1bSJed Brown PetscFunctionBegin; 56c4762a1bSJed Brown if ((!strcmp(p,"hull1972a1")) 57c4762a1bSJed Brown ||(!strcmp(p,"hull1972a2")) 58c4762a1bSJed Brown ||(!strcmp(p,"hull1972a3")) 59c4762a1bSJed Brown ||(!strcmp(p,"hull1972a4")) 60c4762a1bSJed Brown ||(!strcmp(p,"hull1972a5"))) PetscFunctionReturn(1); 61c4762a1bSJed Brown else if (!strcmp(p,"hull1972b1")) PetscFunctionReturn(2); 62c4762a1bSJed Brown else if ((!strcmp(p,"hull1972b2")) 63c4762a1bSJed Brown ||(!strcmp(p,"hull1972b3")) 64c4762a1bSJed Brown ||(!strcmp(p,"hull1972b4")) 65c4762a1bSJed Brown ||(!strcmp(p,"hull1972b5"))) PetscFunctionReturn(3); 66c4762a1bSJed Brown else if ((!strcmp(p,"kulik2013i"))) PetscFunctionReturn(4); 67c4762a1bSJed Brown else if ((!strcmp(p,"hull1972c1")) 68c4762a1bSJed Brown ||(!strcmp(p,"hull1972c2")) 69c4762a1bSJed Brown ||(!strcmp(p,"hull1972c3"))) PetscFunctionReturn(10); 70c4762a1bSJed Brown else if (!strcmp(p,"hull1972c4")) PetscFunctionReturn(51); 71c4762a1bSJed Brown else PetscFunctionReturn(-1); 72c4762a1bSJed Brown } 73c4762a1bSJed Brown 74c4762a1bSJed Brown /****************************************************************/ 75c4762a1bSJed Brown 76c4762a1bSJed Brown /* Problem specific functions */ 77c4762a1bSJed Brown 78c4762a1bSJed Brown /* Hull, 1972, Problem A1 */ 79c4762a1bSJed Brown 80c4762a1bSJed Brown PetscErrorCode RHSFunction_Hull1972A1(TS ts, PetscReal t, Vec Y, Vec F, void *s) 81c4762a1bSJed Brown { 82c4762a1bSJed Brown PetscErrorCode ierr; 83c4762a1bSJed Brown PetscScalar *f; 84c4762a1bSJed Brown const PetscScalar *y; 85c4762a1bSJed Brown 86c4762a1bSJed Brown PetscFunctionBegin; 87c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 88c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 89c4762a1bSJed Brown f[0] = -y[0]; 90c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 91c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 92c4762a1bSJed Brown PetscFunctionReturn(0); 93c4762a1bSJed Brown } 94c4762a1bSJed Brown 95c4762a1bSJed Brown PetscErrorCode RHSJacobian_Hull1972A1(TS ts, PetscReal t, Vec Y, Mat A, Mat B, void *s) 96c4762a1bSJed Brown { 97c4762a1bSJed Brown PetscErrorCode ierr; 98c4762a1bSJed Brown const PetscScalar *y; 99c4762a1bSJed Brown PetscInt row = 0,col = 0; 100c4762a1bSJed Brown PetscScalar value = -1.0; 101c4762a1bSJed Brown 102c4762a1bSJed Brown PetscFunctionBegin; 103c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 104c4762a1bSJed Brown ierr = MatSetValues(A,1,&row,1,&col,&value,INSERT_VALUES);CHKERRQ(ierr); 105c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 106c4762a1bSJed Brown ierr = MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 107c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 108c4762a1bSJed Brown PetscFunctionReturn(0); 109c4762a1bSJed Brown } 110c4762a1bSJed Brown 111c4762a1bSJed Brown PetscErrorCode IFunction_Hull1972A1(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) 112c4762a1bSJed Brown { 113c4762a1bSJed Brown PetscErrorCode ierr; 114c4762a1bSJed Brown const PetscScalar *y; 115c4762a1bSJed Brown PetscScalar *f; 116c4762a1bSJed Brown 117c4762a1bSJed Brown PetscFunctionBegin; 118c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 119c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 120c4762a1bSJed Brown f[0] = -y[0]; 121c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 122c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 123c4762a1bSJed Brown /* Left hand side = ydot - f(y) */ 124c4762a1bSJed Brown ierr = VecAYPX(F,-1.0,Ydot);CHKERRQ(ierr); 125c4762a1bSJed Brown PetscFunctionReturn(0); 126c4762a1bSJed Brown } 127c4762a1bSJed Brown 128c4762a1bSJed Brown PetscErrorCode IJacobian_Hull1972A1(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) 129c4762a1bSJed Brown { 130c4762a1bSJed Brown PetscErrorCode ierr; 131c4762a1bSJed Brown const PetscScalar *y; 132c4762a1bSJed Brown PetscInt row = 0,col = 0; 133c4762a1bSJed Brown PetscScalar value = a - 1.0; 134c4762a1bSJed Brown 135c4762a1bSJed Brown PetscFunctionBegin; 136c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 137c4762a1bSJed Brown ierr = MatSetValues(A,1,&row,1,&col,&value,INSERT_VALUES);CHKERRQ(ierr); 138c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 139c4762a1bSJed Brown ierr = MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 140c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 141c4762a1bSJed Brown PetscFunctionReturn(0); 142c4762a1bSJed Brown } 143c4762a1bSJed Brown 144c4762a1bSJed Brown /* Hull, 1972, Problem A2 */ 145c4762a1bSJed Brown 146c4762a1bSJed Brown PetscErrorCode RHSFunction_Hull1972A2(TS ts, PetscReal t, Vec Y, Vec F, void *s) 147c4762a1bSJed Brown { 148c4762a1bSJed Brown PetscErrorCode ierr; 149c4762a1bSJed Brown const PetscScalar *y; 150c4762a1bSJed Brown PetscScalar *f; 151c4762a1bSJed Brown 152c4762a1bSJed Brown PetscFunctionBegin; 153c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 154c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 155c4762a1bSJed Brown f[0] = -0.5*y[0]*y[0]*y[0]; 156c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 157c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 158c4762a1bSJed Brown PetscFunctionReturn(0); 159c4762a1bSJed Brown } 160c4762a1bSJed Brown 161c4762a1bSJed Brown PetscErrorCode RHSJacobian_Hull1972A2(TS ts, PetscReal t, Vec Y, Mat A, Mat B, void *s) 162c4762a1bSJed Brown { 163c4762a1bSJed Brown PetscErrorCode ierr; 164c4762a1bSJed Brown const PetscScalar *y; 165c4762a1bSJed Brown PetscInt row = 0,col = 0; 166c4762a1bSJed Brown PetscScalar value; 167c4762a1bSJed Brown 168c4762a1bSJed Brown PetscFunctionBegin; 169c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 170c4762a1bSJed Brown value = -0.5*3.0*y[0]*y[0]; 171c4762a1bSJed Brown ierr = MatSetValues(A,1,&row,1,&col,&value,INSERT_VALUES);CHKERRQ(ierr); 172c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 173c4762a1bSJed Brown ierr = MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 174c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 175c4762a1bSJed Brown PetscFunctionReturn(0); 176c4762a1bSJed Brown } 177c4762a1bSJed Brown 178c4762a1bSJed Brown PetscErrorCode IFunction_Hull1972A2(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) 179c4762a1bSJed Brown { 180c4762a1bSJed Brown PetscErrorCode ierr; 181c4762a1bSJed Brown PetscScalar *f; 182c4762a1bSJed Brown const PetscScalar *y; 183c4762a1bSJed Brown 184c4762a1bSJed Brown PetscFunctionBegin; 185c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 186c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 187c4762a1bSJed Brown f[0] = -0.5*y[0]*y[0]*y[0]; 188c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 189c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 190c4762a1bSJed Brown /* Left hand side = ydot - f(y) */ 191c4762a1bSJed Brown ierr = VecAYPX(F,-1.0,Ydot);CHKERRQ(ierr); 192c4762a1bSJed Brown PetscFunctionReturn(0); 193c4762a1bSJed Brown } 194c4762a1bSJed Brown 195c4762a1bSJed Brown PetscErrorCode IJacobian_Hull1972A2(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) 196c4762a1bSJed Brown { 197c4762a1bSJed Brown PetscErrorCode ierr; 198c4762a1bSJed Brown const PetscScalar *y; 199c4762a1bSJed Brown PetscInt row = 0,col = 0; 200c4762a1bSJed Brown PetscScalar value; 201c4762a1bSJed Brown 202c4762a1bSJed Brown PetscFunctionBegin; 203c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 204c4762a1bSJed Brown value = a + 0.5*3.0*y[0]*y[0]; 205c4762a1bSJed Brown ierr = MatSetValues(A,1,&row,1,&col,&value,INSERT_VALUES);CHKERRQ(ierr); 206c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 207c4762a1bSJed Brown ierr = MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 208c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 209c4762a1bSJed Brown PetscFunctionReturn(0); 210c4762a1bSJed Brown } 211c4762a1bSJed Brown 212c4762a1bSJed Brown /* Hull, 1972, Problem A3 */ 213c4762a1bSJed Brown 214c4762a1bSJed Brown PetscErrorCode RHSFunction_Hull1972A3(TS ts, PetscReal t, Vec Y, Vec F, void *s) 215c4762a1bSJed Brown { 216c4762a1bSJed Brown PetscErrorCode ierr; 217c4762a1bSJed Brown const PetscScalar *y; 218c4762a1bSJed Brown PetscScalar *f; 219c4762a1bSJed Brown 220c4762a1bSJed Brown PetscFunctionBegin; 221c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 222c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 223c4762a1bSJed Brown f[0] = y[0]*PetscCosReal(t); 224c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 225c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 226c4762a1bSJed Brown PetscFunctionReturn(0); 227c4762a1bSJed Brown } 228c4762a1bSJed Brown 229c4762a1bSJed Brown PetscErrorCode RHSJacobian_Hull1972A3(TS ts, PetscReal t, Vec Y, Mat A, Mat B, void *s) 230c4762a1bSJed Brown { 231c4762a1bSJed Brown PetscErrorCode ierr; 232c4762a1bSJed Brown const PetscScalar *y; 233c4762a1bSJed Brown PetscInt row = 0,col = 0; 234c4762a1bSJed Brown PetscScalar value = PetscCosReal(t); 235c4762a1bSJed Brown 236c4762a1bSJed Brown PetscFunctionBegin; 237c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 238c4762a1bSJed Brown ierr = MatSetValues(A,1,&row,1,&col,&value,INSERT_VALUES);CHKERRQ(ierr); 239c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 240c4762a1bSJed Brown ierr = MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 241c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 242c4762a1bSJed Brown PetscFunctionReturn(0); 243c4762a1bSJed Brown } 244c4762a1bSJed Brown 245c4762a1bSJed Brown PetscErrorCode IFunction_Hull1972A3(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) 246c4762a1bSJed Brown { 247c4762a1bSJed Brown PetscErrorCode ierr; 248c4762a1bSJed Brown PetscScalar *f; 249c4762a1bSJed Brown const PetscScalar *y; 250c4762a1bSJed Brown 251c4762a1bSJed Brown PetscFunctionBegin; 252c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 253c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 254c4762a1bSJed Brown f[0] = y[0]*PetscCosReal(t); 255c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 256c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 257c4762a1bSJed Brown /* Left hand side = ydot - f(y) */ 258c4762a1bSJed Brown ierr = VecAYPX(F,-1.0,Ydot);CHKERRQ(ierr); 259c4762a1bSJed Brown PetscFunctionReturn(0); 260c4762a1bSJed Brown } 261c4762a1bSJed Brown 262c4762a1bSJed Brown PetscErrorCode IJacobian_Hull1972A3(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) 263c4762a1bSJed Brown { 264c4762a1bSJed Brown PetscErrorCode ierr; 265c4762a1bSJed Brown const PetscScalar *y; 266c4762a1bSJed Brown PetscInt row = 0,col = 0; 267c4762a1bSJed Brown PetscScalar value = a - PetscCosReal(t); 268c4762a1bSJed Brown 269c4762a1bSJed Brown PetscFunctionBegin; 270c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 271c4762a1bSJed Brown ierr = MatSetValues(A,1,&row,1,&col,&value,INSERT_VALUES);CHKERRQ(ierr); 272c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 273c4762a1bSJed Brown ierr = MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 274c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 275c4762a1bSJed Brown PetscFunctionReturn(0); 276c4762a1bSJed Brown } 277c4762a1bSJed Brown 278c4762a1bSJed Brown /* Hull, 1972, Problem A4 */ 279c4762a1bSJed Brown 280c4762a1bSJed Brown PetscErrorCode RHSFunction_Hull1972A4(TS ts, PetscReal t, Vec Y, Vec F, void *s) 281c4762a1bSJed Brown { 282c4762a1bSJed Brown PetscErrorCode ierr; 283c4762a1bSJed Brown PetscScalar *f; 284c4762a1bSJed Brown const PetscScalar *y; 285c4762a1bSJed Brown 286c4762a1bSJed Brown PetscFunctionBegin; 287c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 288c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 289c4762a1bSJed Brown f[0] = (0.25*y[0])*(1.0-0.05*y[0]); 290c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 291c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 292c4762a1bSJed Brown PetscFunctionReturn(0); 293c4762a1bSJed Brown } 294c4762a1bSJed Brown 295c4762a1bSJed Brown PetscErrorCode RHSJacobian_Hull1972A4(TS ts, PetscReal t, Vec Y, Mat A, Mat B, void *s) 296c4762a1bSJed Brown { 297c4762a1bSJed Brown PetscErrorCode ierr; 298c4762a1bSJed Brown const PetscScalar *y; 299c4762a1bSJed Brown PetscInt row = 0,col = 0; 300c4762a1bSJed Brown PetscScalar value; 301c4762a1bSJed Brown 302c4762a1bSJed Brown PetscFunctionBegin; 303c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 304c4762a1bSJed Brown value = 0.25*(1.0-0.05*y[0]) - (0.25*y[0])*0.05; 305c4762a1bSJed Brown ierr = MatSetValues(A,1,&row,1,&col,&value,INSERT_VALUES);CHKERRQ(ierr); 306c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 307c4762a1bSJed Brown ierr = MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 308c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 309c4762a1bSJed Brown PetscFunctionReturn(0); 310c4762a1bSJed Brown } 311c4762a1bSJed Brown 312c4762a1bSJed Brown PetscErrorCode IFunction_Hull1972A4(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) 313c4762a1bSJed Brown { 314c4762a1bSJed Brown PetscErrorCode ierr; 315c4762a1bSJed Brown PetscScalar *f; 316c4762a1bSJed Brown const PetscScalar *y; 317c4762a1bSJed Brown 318c4762a1bSJed Brown PetscFunctionBegin; 319c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 320c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 321c4762a1bSJed Brown f[0] = (0.25*y[0])*(1.0-0.05*y[0]); 322c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 323c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 324c4762a1bSJed Brown /* Left hand side = ydot - f(y) */ 325c4762a1bSJed Brown ierr = VecAYPX(F,-1.0,Ydot);CHKERRQ(ierr); 326c4762a1bSJed Brown PetscFunctionReturn(0); 327c4762a1bSJed Brown } 328c4762a1bSJed Brown 329c4762a1bSJed Brown PetscErrorCode IJacobian_Hull1972A4(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) 330c4762a1bSJed Brown { 331c4762a1bSJed Brown PetscErrorCode ierr; 332c4762a1bSJed Brown const PetscScalar *y; 333c4762a1bSJed Brown PetscInt row = 0,col = 0; 334c4762a1bSJed Brown PetscScalar value; 335c4762a1bSJed Brown 336c4762a1bSJed Brown PetscFunctionBegin; 337c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 338c4762a1bSJed Brown value = a - 0.25*(1.0-0.05*y[0]) + (0.25*y[0])*0.05; 339c4762a1bSJed Brown ierr = MatSetValues(A,1,&row,1,&col,&value,INSERT_VALUES);CHKERRQ(ierr); 340c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 341c4762a1bSJed Brown ierr = MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 342c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 343c4762a1bSJed Brown PetscFunctionReturn(0); 344c4762a1bSJed Brown } 345c4762a1bSJed Brown 346c4762a1bSJed Brown /* Hull, 1972, Problem A5 */ 347c4762a1bSJed Brown 348c4762a1bSJed Brown PetscErrorCode RHSFunction_Hull1972A5(TS ts, PetscReal t, Vec Y, Vec F, void *s) 349c4762a1bSJed Brown { 350c4762a1bSJed Brown PetscErrorCode ierr; 351c4762a1bSJed Brown PetscScalar *f; 352c4762a1bSJed Brown const PetscScalar *y; 353c4762a1bSJed Brown 354c4762a1bSJed Brown PetscFunctionBegin; 355c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 356c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 357c4762a1bSJed Brown f[0] = (y[0]-t)/(y[0]+t); 358c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 359c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 360c4762a1bSJed Brown PetscFunctionReturn(0); 361c4762a1bSJed Brown } 362c4762a1bSJed Brown 363c4762a1bSJed Brown PetscErrorCode RHSJacobian_Hull1972A5(TS ts, PetscReal t, Vec Y, Mat A, Mat B, void *s) 364c4762a1bSJed Brown { 365c4762a1bSJed Brown PetscErrorCode ierr; 366c4762a1bSJed Brown const PetscScalar *y; 367c4762a1bSJed Brown PetscInt row = 0,col = 0; 368c4762a1bSJed Brown PetscScalar value; 369c4762a1bSJed Brown 370c4762a1bSJed Brown PetscFunctionBegin; 371c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 372c4762a1bSJed Brown value = 2*t/((t+y[0])*(t+y[0])); 373c4762a1bSJed Brown ierr = MatSetValues(A,1,&row,1,&col,&value,INSERT_VALUES);CHKERRQ(ierr); 374c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 375c4762a1bSJed Brown ierr = MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 376c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 377c4762a1bSJed Brown PetscFunctionReturn(0); 378c4762a1bSJed Brown } 379c4762a1bSJed Brown 380c4762a1bSJed Brown PetscErrorCode IFunction_Hull1972A5(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) 381c4762a1bSJed Brown { 382c4762a1bSJed Brown PetscErrorCode ierr; 383c4762a1bSJed Brown PetscScalar *f; 384c4762a1bSJed Brown const PetscScalar *y; 385c4762a1bSJed Brown 386c4762a1bSJed Brown PetscFunctionBegin; 387c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 388c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 389c4762a1bSJed Brown f[0] = (y[0]-t)/(y[0]+t); 390c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 391c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 392c4762a1bSJed Brown /* Left hand side = ydot - f(y) */ 393c4762a1bSJed Brown ierr = VecAYPX(F,-1.0,Ydot);CHKERRQ(ierr); 394c4762a1bSJed Brown PetscFunctionReturn(0); 395c4762a1bSJed Brown } 396c4762a1bSJed Brown 397c4762a1bSJed Brown PetscErrorCode IJacobian_Hull1972A5(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) 398c4762a1bSJed Brown { 399c4762a1bSJed Brown PetscErrorCode ierr; 400c4762a1bSJed Brown const PetscScalar *y; 401c4762a1bSJed Brown PetscInt row = 0,col = 0; 402c4762a1bSJed Brown PetscScalar value; 403c4762a1bSJed Brown 404c4762a1bSJed Brown PetscFunctionBegin; 405c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 406c4762a1bSJed Brown value = a - 2*t/((t+y[0])*(t+y[0])); 407c4762a1bSJed Brown ierr = MatSetValues(A,1,&row,1,&col,&value,INSERT_VALUES);CHKERRQ(ierr); 408c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 409c4762a1bSJed Brown ierr = MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 410c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 411c4762a1bSJed Brown PetscFunctionReturn(0); 412c4762a1bSJed Brown } 413c4762a1bSJed Brown 414c4762a1bSJed Brown /* Hull, 1972, Problem B1 */ 415c4762a1bSJed Brown 416c4762a1bSJed Brown PetscErrorCode RHSFunction_Hull1972B1(TS ts, PetscReal t, Vec Y, Vec F, void *s) 417c4762a1bSJed Brown { 418c4762a1bSJed Brown PetscErrorCode ierr; 419c4762a1bSJed Brown PetscScalar *f; 420c4762a1bSJed Brown const PetscScalar *y; 421c4762a1bSJed Brown 422c4762a1bSJed Brown PetscFunctionBegin; 423c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 424c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 425c4762a1bSJed Brown f[0] = 2.0*(y[0] - y[0]*y[1]); 426c4762a1bSJed Brown f[1] = -(y[1]-y[0]*y[1]); 427c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 428c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 429c4762a1bSJed Brown PetscFunctionReturn(0); 430c4762a1bSJed Brown } 431c4762a1bSJed Brown 432c4762a1bSJed Brown PetscErrorCode IFunction_Hull1972B1(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) 433c4762a1bSJed Brown { 434c4762a1bSJed Brown PetscErrorCode ierr; 435c4762a1bSJed Brown PetscScalar *f; 436c4762a1bSJed Brown const PetscScalar *y; 437c4762a1bSJed Brown 438c4762a1bSJed Brown PetscFunctionBegin; 439c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 440c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 441c4762a1bSJed Brown f[0] = 2.0*(y[0] - y[0]*y[1]); 442c4762a1bSJed Brown f[1] = -(y[1]-y[0]*y[1]); 443c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 444c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 445c4762a1bSJed Brown /* Left hand side = ydot - f(y) */ 446c4762a1bSJed Brown ierr = VecAYPX(F,-1.0,Ydot);CHKERRQ(ierr); 447c4762a1bSJed Brown PetscFunctionReturn(0); 448c4762a1bSJed Brown } 449c4762a1bSJed Brown 450c4762a1bSJed Brown PetscErrorCode IJacobian_Hull1972B1(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) 451c4762a1bSJed Brown { 452c4762a1bSJed Brown PetscErrorCode ierr; 453c4762a1bSJed Brown const PetscScalar *y; 454c4762a1bSJed Brown PetscInt row[2] = {0,1}; 455c4762a1bSJed Brown PetscScalar value[2][2]; 456c4762a1bSJed Brown 457c4762a1bSJed Brown PetscFunctionBegin; 458c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 459c4762a1bSJed Brown value[0][0] = a - 2.0*(1.0-y[1]); value[0][1] = 2.0*y[0]; 460c4762a1bSJed Brown value[1][0] = -y[1]; value[1][1] = a + 1.0 - y[0]; 461c4762a1bSJed Brown ierr = MatSetValues(A,2,&row[0],2,&row[0],&value[0][0],INSERT_VALUES);CHKERRQ(ierr); 462c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 463c4762a1bSJed Brown ierr = MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 464c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 465c4762a1bSJed Brown PetscFunctionReturn(0); 466c4762a1bSJed Brown } 467c4762a1bSJed Brown 468c4762a1bSJed Brown /* Hull, 1972, Problem B2 */ 469c4762a1bSJed Brown 470c4762a1bSJed Brown PetscErrorCode RHSFunction_Hull1972B2(TS ts, PetscReal t, Vec Y, Vec F, void *s) 471c4762a1bSJed Brown { 472c4762a1bSJed Brown PetscErrorCode ierr; 473c4762a1bSJed Brown PetscScalar *f; 474c4762a1bSJed Brown const PetscScalar *y; 475c4762a1bSJed Brown 476c4762a1bSJed Brown PetscFunctionBegin; 477c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 478c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 479c4762a1bSJed Brown f[0] = -y[0] + y[1]; 480c4762a1bSJed Brown f[1] = y[0] - 2.0*y[1] + y[2]; 481c4762a1bSJed Brown f[2] = y[1] - y[2]; 482c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 483c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 484c4762a1bSJed Brown PetscFunctionReturn(0); 485c4762a1bSJed Brown } 486c4762a1bSJed Brown 487c4762a1bSJed Brown PetscErrorCode IFunction_Hull1972B2(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) 488c4762a1bSJed Brown { 489c4762a1bSJed Brown PetscErrorCode ierr; 490c4762a1bSJed Brown PetscScalar *f; 491c4762a1bSJed Brown const PetscScalar *y; 492c4762a1bSJed Brown 493c4762a1bSJed Brown PetscFunctionBegin; 494c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 495c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 496c4762a1bSJed Brown f[0] = -y[0] + y[1]; 497c4762a1bSJed Brown f[1] = y[0] - 2.0*y[1] + y[2]; 498c4762a1bSJed Brown f[2] = y[1] - y[2]; 499c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 500c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 501c4762a1bSJed Brown /* Left hand side = ydot - f(y) */ 502c4762a1bSJed Brown ierr = VecAYPX(F,-1.0,Ydot);CHKERRQ(ierr); 503c4762a1bSJed Brown PetscFunctionReturn(0); 504c4762a1bSJed Brown } 505c4762a1bSJed Brown 506c4762a1bSJed Brown PetscErrorCode IJacobian_Hull1972B2(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) 507c4762a1bSJed Brown { 508c4762a1bSJed Brown PetscErrorCode ierr; 509c4762a1bSJed Brown const PetscScalar *y; 510c4762a1bSJed Brown PetscInt row[3] = {0,1,2}; 511c4762a1bSJed Brown PetscScalar value[3][3]; 512c4762a1bSJed Brown 513c4762a1bSJed Brown PetscFunctionBegin; 514c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 515c4762a1bSJed Brown value[0][0] = a + 1.0; value[0][1] = -1.0; value[0][2] = 0; 516c4762a1bSJed Brown value[1][0] = -1.0; value[1][1] = a + 2.0; value[1][2] = -1.0; 517c4762a1bSJed Brown value[2][0] = 0; value[2][1] = -1.0; value[2][2] = a + 1.0; 518c4762a1bSJed Brown ierr = MatSetValues(A,3,&row[0],3,&row[0],&value[0][0],INSERT_VALUES);CHKERRQ(ierr); 519c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 520c4762a1bSJed Brown ierr = MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 521c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 522c4762a1bSJed Brown PetscFunctionReturn(0); 523c4762a1bSJed Brown } 524c4762a1bSJed Brown 525c4762a1bSJed Brown /* Hull, 1972, Problem B3 */ 526c4762a1bSJed Brown 527c4762a1bSJed Brown PetscErrorCode RHSFunction_Hull1972B3(TS ts, PetscReal t, Vec Y, Vec F, void *s) 528c4762a1bSJed Brown { 529c4762a1bSJed Brown PetscErrorCode ierr; 530c4762a1bSJed Brown PetscScalar *f; 531c4762a1bSJed Brown const PetscScalar *y; 532c4762a1bSJed Brown 533c4762a1bSJed Brown PetscFunctionBegin; 534c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 535c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 536c4762a1bSJed Brown f[0] = -y[0]; 537c4762a1bSJed Brown f[1] = y[0] - y[1]*y[1]; 538c4762a1bSJed Brown f[2] = y[1]*y[1]; 539c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 540c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 541c4762a1bSJed Brown PetscFunctionReturn(0); 542c4762a1bSJed Brown } 543c4762a1bSJed Brown 544c4762a1bSJed Brown PetscErrorCode IFunction_Hull1972B3(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) 545c4762a1bSJed Brown { 546c4762a1bSJed Brown PetscErrorCode ierr; 547c4762a1bSJed Brown PetscScalar *f; 548c4762a1bSJed Brown const PetscScalar *y; 549c4762a1bSJed Brown 550c4762a1bSJed Brown PetscFunctionBegin; 551c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 552c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 553c4762a1bSJed Brown f[0] = -y[0]; 554c4762a1bSJed Brown f[1] = y[0] - y[1]*y[1]; 555c4762a1bSJed Brown f[2] = y[1]*y[1]; 556c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 557c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 558c4762a1bSJed Brown /* Left hand side = ydot - f(y) */ 559c4762a1bSJed Brown ierr = VecAYPX(F,-1.0,Ydot);CHKERRQ(ierr); 560c4762a1bSJed Brown PetscFunctionReturn(0); 561c4762a1bSJed Brown } 562c4762a1bSJed Brown 563c4762a1bSJed Brown PetscErrorCode IJacobian_Hull1972B3(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) 564c4762a1bSJed Brown { 565c4762a1bSJed Brown PetscErrorCode ierr; 566c4762a1bSJed Brown const PetscScalar *y; 567c4762a1bSJed Brown PetscInt row[3] = {0,1,2}; 568c4762a1bSJed Brown PetscScalar value[3][3]; 569c4762a1bSJed Brown 570c4762a1bSJed Brown PetscFunctionBegin; 571c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 572c4762a1bSJed Brown value[0][0] = a + 1.0; value[0][1] = 0; value[0][2] = 0; 573c4762a1bSJed Brown value[1][0] = -1.0; value[1][1] = a + 2.0*y[1]; value[1][2] = 0; 574c4762a1bSJed Brown value[2][0] = 0; value[2][1] = -2.0*y[1]; value[2][2] = a; 575c4762a1bSJed Brown ierr = MatSetValues(A,3,&row[0],3,&row[0],&value[0][0],INSERT_VALUES);CHKERRQ(ierr); 576c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 577c4762a1bSJed Brown ierr = MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 578c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 579c4762a1bSJed Brown PetscFunctionReturn(0); 580c4762a1bSJed Brown } 581c4762a1bSJed Brown 582c4762a1bSJed Brown /* Hull, 1972, Problem B4 */ 583c4762a1bSJed Brown 584c4762a1bSJed Brown PetscErrorCode RHSFunction_Hull1972B4(TS ts, PetscReal t, Vec Y, Vec F, void *s) 585c4762a1bSJed Brown { 586c4762a1bSJed Brown PetscErrorCode ierr; 587c4762a1bSJed Brown PetscScalar *f; 588c4762a1bSJed Brown const PetscScalar *y; 589c4762a1bSJed Brown 590c4762a1bSJed Brown PetscFunctionBegin; 591c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 592c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 593c4762a1bSJed Brown f[0] = -y[1] - y[0]*y[2]/PetscSqrtScalar(y[0]*y[0]+y[1]*y[1]); 594c4762a1bSJed Brown f[1] = y[0] - y[1]*y[2]/PetscSqrtScalar(y[0]*y[0]+y[1]*y[1]); 595c4762a1bSJed Brown f[2] = y[0]/PetscSqrtScalar(y[0]*y[0]+y[1]*y[1]); 596c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 597c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 598c4762a1bSJed Brown PetscFunctionReturn(0); 599c4762a1bSJed Brown } 600c4762a1bSJed Brown 601c4762a1bSJed Brown PetscErrorCode IFunction_Hull1972B4(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) 602c4762a1bSJed Brown { 603c4762a1bSJed Brown PetscErrorCode ierr; 604c4762a1bSJed Brown PetscScalar *f; 605c4762a1bSJed Brown const PetscScalar *y; 606c4762a1bSJed Brown 607c4762a1bSJed Brown PetscFunctionBegin; 608c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 609c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 610c4762a1bSJed Brown f[0] = -y[1] - y[0]*y[2]/PetscSqrtScalar(y[0]*y[0]+y[1]*y[1]); 611c4762a1bSJed Brown f[1] = y[0] - y[1]*y[2]/PetscSqrtScalar(y[0]*y[0]+y[1]*y[1]); 612c4762a1bSJed Brown f[2] = y[0]/PetscSqrtScalar(y[0]*y[0]+y[1]*y[1]); 613c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 614c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 615c4762a1bSJed Brown /* Left hand side = ydot - f(y) */ 616c4762a1bSJed Brown ierr = VecAYPX(F,-1.0,Ydot);CHKERRQ(ierr); 617c4762a1bSJed Brown PetscFunctionReturn(0); 618c4762a1bSJed Brown } 619c4762a1bSJed Brown 620c4762a1bSJed Brown PetscErrorCode IJacobian_Hull1972B4(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) 621c4762a1bSJed Brown { 622c4762a1bSJed Brown PetscErrorCode ierr; 623c4762a1bSJed Brown const PetscScalar *y; 624c4762a1bSJed Brown PetscInt row[3] = {0,1,2}; 625c4762a1bSJed Brown PetscScalar value[3][3],fac,fac2; 626c4762a1bSJed Brown 627c4762a1bSJed Brown PetscFunctionBegin; 628c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 629c4762a1bSJed Brown fac = PetscPowScalar(y[0]*y[0]+y[1]*y[1],-1.5); 630c4762a1bSJed Brown fac2 = PetscPowScalar(y[0]*y[0]+y[1]*y[1],-0.5); 631c4762a1bSJed Brown value[0][0] = a + (y[1]*y[1]*y[2])*fac; 632c4762a1bSJed Brown value[0][1] = 1.0 - (y[0]*y[1]*y[2])*fac; 633c4762a1bSJed Brown value[0][2] = y[0]*fac2; 634c4762a1bSJed Brown value[1][0] = -1.0 - y[0]*y[1]*y[2]*fac; 635c4762a1bSJed Brown value[1][1] = a + y[0]*y[0]*y[2]*fac; 636c4762a1bSJed Brown value[1][2] = y[1]*fac2; 637c4762a1bSJed Brown value[2][0] = -y[1]*y[1]*fac; 638c4762a1bSJed Brown value[2][1] = y[0]*y[1]*fac; 639c4762a1bSJed Brown value[2][2] = a; 640c4762a1bSJed Brown ierr = MatSetValues(A,3,&row[0],3,&row[0],&value[0][0],INSERT_VALUES);CHKERRQ(ierr); 641c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 642c4762a1bSJed Brown ierr = MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 643c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 644c4762a1bSJed Brown PetscFunctionReturn(0); 645c4762a1bSJed Brown } 646c4762a1bSJed Brown 647c4762a1bSJed Brown /* Hull, 1972, Problem B5 */ 648c4762a1bSJed Brown 649c4762a1bSJed Brown PetscErrorCode RHSFunction_Hull1972B5(TS ts, PetscReal t, Vec Y, Vec F, void *s) 650c4762a1bSJed Brown { 651c4762a1bSJed Brown PetscErrorCode ierr; 652c4762a1bSJed Brown PetscScalar *f; 653c4762a1bSJed Brown const PetscScalar *y; 654c4762a1bSJed Brown 655c4762a1bSJed Brown PetscFunctionBegin; 656c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 657c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 658c4762a1bSJed Brown f[0] = y[1]*y[2]; 659c4762a1bSJed Brown f[1] = -y[0]*y[2]; 660c4762a1bSJed Brown f[2] = -0.51*y[0]*y[1]; 661c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 662c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 663c4762a1bSJed Brown PetscFunctionReturn(0); 664c4762a1bSJed Brown } 665c4762a1bSJed Brown 666c4762a1bSJed Brown PetscErrorCode IFunction_Hull1972B5(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) 667c4762a1bSJed Brown { 668c4762a1bSJed Brown PetscErrorCode ierr; 669c4762a1bSJed Brown PetscScalar *f; 670c4762a1bSJed Brown const PetscScalar *y; 671c4762a1bSJed Brown 672c4762a1bSJed Brown PetscFunctionBegin; 673c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 674c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 675c4762a1bSJed Brown f[0] = y[1]*y[2]; 676c4762a1bSJed Brown f[1] = -y[0]*y[2]; 677c4762a1bSJed Brown f[2] = -0.51*y[0]*y[1]; 678c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 679c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 680c4762a1bSJed Brown /* Left hand side = ydot - f(y) */ 681c4762a1bSJed Brown ierr = VecAYPX(F,-1.0,Ydot);CHKERRQ(ierr); 682c4762a1bSJed Brown PetscFunctionReturn(0); 683c4762a1bSJed Brown } 684c4762a1bSJed Brown 685c4762a1bSJed Brown PetscErrorCode IJacobian_Hull1972B5(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) 686c4762a1bSJed Brown { 687c4762a1bSJed Brown PetscErrorCode ierr; 688c4762a1bSJed Brown const PetscScalar *y; 689c4762a1bSJed Brown PetscInt row[3] = {0,1,2}; 690c4762a1bSJed Brown PetscScalar value[3][3]; 691c4762a1bSJed Brown 692c4762a1bSJed Brown PetscFunctionBegin; 693c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 694c4762a1bSJed Brown value[0][0] = a; value[0][1] = -y[2]; value[0][2] = -y[1]; 695c4762a1bSJed Brown value[1][0] = y[2]; value[1][1] = a; value[1][2] = y[0]; 696c4762a1bSJed Brown value[2][0] = 0.51*y[1]; value[2][1] = 0.51*y[0]; value[2][2] = a; 697c4762a1bSJed Brown ierr = MatSetValues(A,3,&row[0],3,&row[0],&value[0][0],INSERT_VALUES);CHKERRQ(ierr); 698c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 699c4762a1bSJed Brown ierr = MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 700c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 701c4762a1bSJed Brown PetscFunctionReturn(0); 702c4762a1bSJed Brown } 703c4762a1bSJed Brown 704c4762a1bSJed Brown /* Kulikov, 2013, Problem I */ 705c4762a1bSJed Brown 706c4762a1bSJed Brown PetscErrorCode RHSFunction_Kulikov2013I(TS ts, PetscReal t, Vec Y, Vec F, void *s) 707c4762a1bSJed Brown { 708c4762a1bSJed Brown PetscErrorCode ierr; 709c4762a1bSJed Brown PetscScalar *f; 710c4762a1bSJed Brown const PetscScalar *y; 711c4762a1bSJed Brown 712c4762a1bSJed Brown PetscFunctionBegin; 713c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 714c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 715c4762a1bSJed Brown f[0] = 2.*t*PetscPowScalar(y[1],1./5.)*y[3]; 716c4762a1bSJed Brown f[1] = 10.*t*y[3]*PetscExpScalar(5.0*(y[2]-1.)); 717c4762a1bSJed Brown f[2] = 2.*t*y[3]; 718c4762a1bSJed Brown f[3] = -2.*t*PetscLogScalar(y[0]); 719c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 720c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 721c4762a1bSJed Brown PetscFunctionReturn(0); 722c4762a1bSJed Brown } 723c4762a1bSJed Brown 724c4762a1bSJed Brown PetscErrorCode RHSJacobian_Kulikov2013I(TS ts, PetscReal t, Vec Y, Mat A, Mat B, void *s) 725c4762a1bSJed Brown { 726c4762a1bSJed Brown PetscErrorCode ierr; 727c4762a1bSJed Brown const PetscScalar *y; 728c4762a1bSJed Brown PetscInt row[4] = {0,1,2,3}; 729c4762a1bSJed Brown PetscScalar value[4][4]; 730c4762a1bSJed Brown PetscScalar m1,m2; 731c4762a1bSJed Brown PetscFunctionBegin; 732c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 733c4762a1bSJed Brown m1=(2.*t*y[3])/(5.*PetscPowScalar(y[1],4./5.)); 734c4762a1bSJed Brown m2=2.*t*PetscPowScalar(y[1],1./5.); 735c4762a1bSJed Brown value[0][0] = 0. ; value[0][1] = m1; value[0][2] = 0.; value[0][3] = m2; 736c4762a1bSJed Brown m1=50.*t*y[3]*PetscExpScalar(5.0*(y[2]-1.)); 737c4762a1bSJed Brown m2=10.*t*PetscExpScalar(5.0*(y[2]-1.)); 738c4762a1bSJed Brown value[1][0] = 0.; value[1][1] = 0. ; value[1][2] = m1; value[1][3] = m2; 739c4762a1bSJed Brown value[2][0] = 0.; value[2][1] = 0.; value[2][2] = 0.; value[2][3] = 2*t; 740c4762a1bSJed Brown value[3][0] = -2.*t/y[0];value[3][1] = 0.; value[3][2] = 0.; value[3][3] = 0.; 741c4762a1bSJed Brown ierr = MatSetValues(A,4,&row[0],4,&row[0],&value[0][0],INSERT_VALUES);CHKERRQ(ierr); 742c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 743c4762a1bSJed Brown ierr = MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 744c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 745c4762a1bSJed Brown PetscFunctionReturn(0); 746c4762a1bSJed Brown } 747c4762a1bSJed Brown 748c4762a1bSJed Brown PetscErrorCode IFunction_Kulikov2013I(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) 749c4762a1bSJed Brown { 750c4762a1bSJed Brown PetscErrorCode ierr; 751c4762a1bSJed Brown PetscScalar *f; 752c4762a1bSJed Brown const PetscScalar *y; 753c4762a1bSJed Brown 754c4762a1bSJed Brown PetscFunctionBegin; 755c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 756c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 757c4762a1bSJed Brown f[0] = 2.*t*PetscPowScalar(y[1],1./5.)*y[3]; 758c4762a1bSJed Brown f[1] = 10.*t*y[3]*PetscExpScalar(5.0*(y[2]-1.)); 759c4762a1bSJed Brown f[2] = 2.*t*y[3]; 760c4762a1bSJed Brown f[3] = -2.*t*PetscLogScalar(y[0]); 761c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 762c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 763c4762a1bSJed Brown /* Left hand side = ydot - f(y) */ 764c4762a1bSJed Brown ierr = VecAYPX(F,-1.0,Ydot);CHKERRQ(ierr); 765c4762a1bSJed Brown PetscFunctionReturn(0); 766c4762a1bSJed Brown } 767c4762a1bSJed Brown 768c4762a1bSJed Brown PetscErrorCode IJacobian_Kulikov2013I(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) 769c4762a1bSJed Brown { 770c4762a1bSJed Brown PetscErrorCode ierr; 771c4762a1bSJed Brown const PetscScalar *y; 772c4762a1bSJed Brown PetscInt row[4] = {0,1,2,3}; 773c4762a1bSJed Brown PetscScalar value[4][4]; 774c4762a1bSJed Brown PetscScalar m1,m2; 775c4762a1bSJed Brown 776c4762a1bSJed Brown PetscFunctionBegin; 777c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 778c4762a1bSJed Brown m1=(2.*t*y[3])/(5.*PetscPowScalar(y[1],4./5.)); 779c4762a1bSJed Brown m2=2.*t*PetscPowScalar(y[1],1./5.); 780c4762a1bSJed Brown value[0][0] = a ; value[0][1] = m1; value[0][2] = 0.; value[0][3] = m2; 781c4762a1bSJed Brown m1=50.*t*y[3]*PetscExpScalar(5.0*(y[2]-1.)); 782c4762a1bSJed Brown m2=10.*t*PetscExpScalar(5.0*(y[2]-1.)); 783c4762a1bSJed Brown value[1][0] = 0.; value[1][1] = a ; value[1][2] = m1; value[1][3] = m2; 784c4762a1bSJed Brown value[2][0] = 0.; value[2][1] = 0.; value[2][2] = a; value[2][3] = 2*t; 785c4762a1bSJed Brown value[3][0] = -2.*t/y[0];value[3][1] = 0.; value[3][2] = 0.; value[3][3] = a; 786c4762a1bSJed Brown ierr = MatSetValues(A,4,&row[0],4,&row[0],&value[0][0],INSERT_VALUES);CHKERRQ(ierr); 787c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 788c4762a1bSJed Brown ierr = MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 789c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 790c4762a1bSJed Brown PetscFunctionReturn(0); 791c4762a1bSJed Brown } 792c4762a1bSJed Brown 793c4762a1bSJed Brown /* Hull, 1972, Problem C1 */ 794c4762a1bSJed Brown 795c4762a1bSJed Brown PetscErrorCode RHSFunction_Hull1972C1(TS ts, PetscReal t, Vec Y, Vec F, void *s) 796c4762a1bSJed Brown { 797c4762a1bSJed Brown PetscErrorCode ierr; 798c4762a1bSJed Brown PetscScalar *f; 799c4762a1bSJed Brown const PetscScalar *y; 800c4762a1bSJed Brown PetscInt N,i; 801c4762a1bSJed Brown 802c4762a1bSJed Brown PetscFunctionBegin; 803c4762a1bSJed Brown ierr = VecGetSize (Y,&N);CHKERRQ(ierr); 804c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 805c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 806c4762a1bSJed Brown f[0] = -y[0]; 807c4762a1bSJed Brown for (i = 1; i < N-1; i++) { 808c4762a1bSJed Brown f[i] = y[i-1] - y[i]; 809c4762a1bSJed Brown } 810c4762a1bSJed Brown f[N-1] = y[N-2]; 811c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 812c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 813c4762a1bSJed Brown PetscFunctionReturn(0); 814c4762a1bSJed Brown } 815c4762a1bSJed Brown 816c4762a1bSJed Brown PetscErrorCode IFunction_Hull1972C1(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) 817c4762a1bSJed Brown { 818c4762a1bSJed Brown PetscErrorCode ierr; 819c4762a1bSJed Brown PetscScalar *f; 820c4762a1bSJed Brown const PetscScalar *y; 821c4762a1bSJed Brown PetscInt N,i; 822c4762a1bSJed Brown 823c4762a1bSJed Brown PetscFunctionBegin; 824c4762a1bSJed Brown ierr = VecGetSize (Y,&N);CHKERRQ(ierr); 825c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 826c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 827c4762a1bSJed Brown f[0] = -y[0]; 828c4762a1bSJed Brown for (i = 1; i < N-1; i++) { 829c4762a1bSJed Brown f[i] = y[i-1] - y[i]; 830c4762a1bSJed Brown } 831c4762a1bSJed Brown f[N-1] = y[N-2]; 832c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 833c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 834c4762a1bSJed Brown /* Left hand side = ydot - f(y) */ 835c4762a1bSJed Brown ierr = VecAYPX(F,-1.0,Ydot);CHKERRQ(ierr); 836c4762a1bSJed Brown PetscFunctionReturn(0); 837c4762a1bSJed Brown } 838c4762a1bSJed Brown 839c4762a1bSJed Brown PetscErrorCode IJacobian_Hull1972C1(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) 840c4762a1bSJed Brown { 841c4762a1bSJed Brown PetscErrorCode ierr; 842c4762a1bSJed Brown const PetscScalar *y; 843c4762a1bSJed Brown PetscInt N,i,col[2]; 844c4762a1bSJed Brown PetscScalar value[2]; 845c4762a1bSJed Brown 846c4762a1bSJed Brown PetscFunctionBegin; 847c4762a1bSJed Brown ierr = VecGetSize (Y,&N);CHKERRQ(ierr); 848c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 849c4762a1bSJed Brown i = 0; 850c4762a1bSJed Brown value[0] = a+1; col[0] = 0; 851c4762a1bSJed Brown value[1] = 0; col[1] = 1; 852c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr); 853c4762a1bSJed Brown for (i = 0; i < N; i++) { 854c4762a1bSJed Brown value[0] = -1; col[0] = i-1; 855c4762a1bSJed Brown value[1] = a+1; col[1] = i; 856c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr); 857c4762a1bSJed Brown } 858c4762a1bSJed Brown i = N-1; 859c4762a1bSJed Brown value[0] = -1; col[0] = N-2; 860c4762a1bSJed Brown value[1] = a; col[1] = N-1; 861c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr); 862c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 863c4762a1bSJed Brown ierr = MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 864c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 865c4762a1bSJed Brown PetscFunctionReturn(0); 866c4762a1bSJed Brown } 867c4762a1bSJed Brown 868c4762a1bSJed Brown /* Hull, 1972, Problem C2 */ 869c4762a1bSJed Brown 870c4762a1bSJed Brown PetscErrorCode RHSFunction_Hull1972C2(TS ts, PetscReal t, Vec Y, Vec F, void *s) 871c4762a1bSJed Brown { 872c4762a1bSJed Brown PetscErrorCode ierr; 873c4762a1bSJed Brown const PetscScalar *y; 874c4762a1bSJed Brown PetscScalar *f; 875c4762a1bSJed Brown PetscInt N,i; 876c4762a1bSJed Brown 877c4762a1bSJed Brown PetscFunctionBegin; 878c4762a1bSJed Brown ierr = VecGetSize (Y,&N);CHKERRQ(ierr); 879c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 880c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 881c4762a1bSJed Brown f[0] = -y[0]; 882c4762a1bSJed Brown for (i = 1; i < N-1; i++) { 883c4762a1bSJed Brown f[i] = (PetscReal)i*y[i-1] - (PetscReal)(i+1)*y[i]; 884c4762a1bSJed Brown } 885c4762a1bSJed Brown f[N-1] = (PetscReal)(N-1)*y[N-2]; 886c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 887c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 888c4762a1bSJed Brown PetscFunctionReturn(0); 889c4762a1bSJed Brown } 890c4762a1bSJed Brown 891c4762a1bSJed Brown PetscErrorCode IFunction_Hull1972C2(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) 892c4762a1bSJed Brown { 893c4762a1bSJed Brown PetscErrorCode ierr; 894c4762a1bSJed Brown PetscScalar *f; 895c4762a1bSJed Brown const PetscScalar *y; 896c4762a1bSJed Brown PetscInt N,i; 897c4762a1bSJed Brown 898c4762a1bSJed Brown PetscFunctionBegin; 899c4762a1bSJed Brown ierr = VecGetSize (Y,&N);CHKERRQ(ierr); 900c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 901c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 902c4762a1bSJed Brown f[0] = -y[0]; 903c4762a1bSJed Brown for (i = 1; i < N-1; i++) { 904c4762a1bSJed Brown f[i] = (PetscReal)i*y[i-1] - (PetscReal)(i+1)*y[i]; 905c4762a1bSJed Brown } 906c4762a1bSJed Brown f[N-1] = (PetscReal)(N-1)*y[N-2]; 907c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 908c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 909c4762a1bSJed Brown /* Left hand side = ydot - f(y) */ 910c4762a1bSJed Brown ierr = VecAYPX(F,-1.0,Ydot);CHKERRQ(ierr); 911c4762a1bSJed Brown PetscFunctionReturn(0); 912c4762a1bSJed Brown } 913c4762a1bSJed Brown 914c4762a1bSJed Brown PetscErrorCode IJacobian_Hull1972C2(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) 915c4762a1bSJed Brown { 916c4762a1bSJed Brown PetscErrorCode ierr; 917c4762a1bSJed Brown const PetscScalar *y; 918c4762a1bSJed Brown PetscInt N,i,col[2]; 919c4762a1bSJed Brown PetscScalar value[2]; 920c4762a1bSJed Brown 921c4762a1bSJed Brown PetscFunctionBegin; 922c4762a1bSJed Brown ierr = VecGetSize (Y,&N);CHKERRQ(ierr); 923c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 924c4762a1bSJed Brown i = 0; 925c4762a1bSJed Brown value[0] = a+1; col[0] = 0; 926c4762a1bSJed Brown value[1] = 0; col[1] = 1; 927c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr); 928c4762a1bSJed Brown for (i = 0; i < N; i++) { 929c4762a1bSJed Brown value[0] = -(PetscReal) i; col[0] = i-1; 930c4762a1bSJed Brown value[1] = a+(PetscReal)(i+1); col[1] = i; 931c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr); 932c4762a1bSJed Brown } 933c4762a1bSJed Brown i = N-1; 934c4762a1bSJed Brown value[0] = -(PetscReal) (N-1); col[0] = N-2; 935c4762a1bSJed Brown value[1] = a; col[1] = N-1; 936c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr); 937c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 938c4762a1bSJed Brown ierr = MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 939c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 940c4762a1bSJed Brown PetscFunctionReturn(0); 941c4762a1bSJed Brown } 942c4762a1bSJed Brown 943c4762a1bSJed Brown /* Hull, 1972, Problem C3 and C4 */ 944c4762a1bSJed Brown 945c4762a1bSJed Brown PetscErrorCode RHSFunction_Hull1972C34(TS ts, PetscReal t, Vec Y, Vec F, void *s) 946c4762a1bSJed Brown { 947c4762a1bSJed Brown PetscErrorCode ierr; 948c4762a1bSJed Brown PetscScalar *f; 949c4762a1bSJed Brown const PetscScalar *y; 950c4762a1bSJed Brown PetscInt N,i; 951c4762a1bSJed Brown 952c4762a1bSJed Brown PetscFunctionBegin; 953c4762a1bSJed Brown ierr = VecGetSize (Y,&N);CHKERRQ(ierr); 954c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 955c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 956c4762a1bSJed Brown f[0] = -2.0*y[0] + y[1]; 957c4762a1bSJed Brown for (i = 1; i < N-1; i++) { 958c4762a1bSJed Brown f[i] = y[i-1] - 2.0*y[i] + y[i+1]; 959c4762a1bSJed Brown } 960c4762a1bSJed Brown f[N-1] = y[N-2] - 2.0*y[N-1]; 961c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 962c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 963c4762a1bSJed Brown PetscFunctionReturn(0); 964c4762a1bSJed Brown } 965c4762a1bSJed Brown 966c4762a1bSJed Brown PetscErrorCode IFunction_Hull1972C34(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) 967c4762a1bSJed Brown { 968c4762a1bSJed Brown PetscErrorCode ierr; 969c4762a1bSJed Brown PetscScalar *f; 970c4762a1bSJed Brown const PetscScalar *y; 971c4762a1bSJed Brown PetscInt N,i; 972c4762a1bSJed Brown 973c4762a1bSJed Brown PetscFunctionBegin; 974c4762a1bSJed Brown ierr = VecGetSize (Y,&N);CHKERRQ(ierr); 975c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 976c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 977c4762a1bSJed Brown f[0] = -2.0*y[0] + y[1]; 978c4762a1bSJed Brown for (i = 1; i < N-1; i++) { 979c4762a1bSJed Brown f[i] = y[i-1] - 2.0*y[i] + y[i+1]; 980c4762a1bSJed Brown } 981c4762a1bSJed Brown f[N-1] = y[N-2] - 2.0*y[N-1]; 982c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 983c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 984c4762a1bSJed Brown /* Left hand side = ydot - f(y) */ 985c4762a1bSJed Brown ierr = VecAYPX(F,-1.0,Ydot);CHKERRQ(ierr); 986c4762a1bSJed Brown PetscFunctionReturn(0); 987c4762a1bSJed Brown } 988c4762a1bSJed Brown 989c4762a1bSJed Brown PetscErrorCode IJacobian_Hull1972C34(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) 990c4762a1bSJed Brown { 991c4762a1bSJed Brown PetscErrorCode ierr; 992c4762a1bSJed Brown const PetscScalar *y; 993c4762a1bSJed Brown PetscScalar value[3]; 994c4762a1bSJed Brown PetscInt N,i,col[3]; 995c4762a1bSJed Brown 996c4762a1bSJed Brown PetscFunctionBegin; 997c4762a1bSJed Brown ierr = VecGetSize (Y,&N);CHKERRQ(ierr); 998c4762a1bSJed Brown ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr); 999c4762a1bSJed Brown for (i = 0; i < N; i++) { 1000c4762a1bSJed Brown if (i == 0) { 1001c4762a1bSJed Brown value[0] = a+2; col[0] = i; 1002c4762a1bSJed Brown value[1] = -1; col[1] = i+1; 1003c4762a1bSJed Brown value[2] = 0; col[2] = i+2; 1004c4762a1bSJed Brown } else if (i == N-1) { 1005c4762a1bSJed Brown value[0] = 0; col[0] = i-2; 1006c4762a1bSJed Brown value[1] = -1; col[1] = i-1; 1007c4762a1bSJed Brown value[2] = a+2; col[2] = i; 1008c4762a1bSJed Brown } else { 1009c4762a1bSJed Brown value[0] = -1; col[0] = i-1; 1010c4762a1bSJed Brown value[1] = a+2; col[1] = i; 1011c4762a1bSJed Brown value[2] = -1; col[2] = i+1; 1012c4762a1bSJed Brown } 1013c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);CHKERRQ(ierr); 1014c4762a1bSJed Brown } 1015c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1016c4762a1bSJed Brown ierr = MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1017c4762a1bSJed Brown ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr); 1018c4762a1bSJed Brown PetscFunctionReturn(0); 1019c4762a1bSJed Brown } 1020c4762a1bSJed Brown 1021c4762a1bSJed Brown /***************************************************************************/ 1022c4762a1bSJed Brown 1023c4762a1bSJed Brown /* Sets the initial solution for the IVP and sets up the function pointers*/ 1024c4762a1bSJed Brown PetscErrorCode Initialize(Vec Y, void* s) 1025c4762a1bSJed Brown { 1026c4762a1bSJed Brown PetscErrorCode ierr; 1027c4762a1bSJed Brown char *p = (char*) s; 1028c4762a1bSJed Brown PetscScalar *y; 1029c4762a1bSJed Brown PetscReal t0; 1030c4762a1bSJed Brown PetscInt N = GetSize((const char *)s); 1031c4762a1bSJed Brown PetscBool flg; 1032c4762a1bSJed Brown 1033c4762a1bSJed Brown PetscFunctionBegin; 1034c4762a1bSJed Brown VecZeroEntries(Y); 1035c4762a1bSJed Brown ierr = VecGetArray(Y,&y);CHKERRQ(ierr); 1036c4762a1bSJed Brown if (!strcmp(p,"hull1972a1")) { 1037c4762a1bSJed Brown y[0] = 1.0; 1038c4762a1bSJed Brown RHSFunction = RHSFunction_Hull1972A1; 1039c4762a1bSJed Brown RHSJacobian = RHSJacobian_Hull1972A1; 1040c4762a1bSJed Brown IFunction = IFunction_Hull1972A1; 1041c4762a1bSJed Brown IJacobian = IJacobian_Hull1972A1; 1042c4762a1bSJed Brown } else if (!strcmp(p,"hull1972a2")) { 1043c4762a1bSJed Brown y[0] = 1.0; 1044c4762a1bSJed Brown RHSFunction = RHSFunction_Hull1972A2; 1045c4762a1bSJed Brown RHSJacobian = RHSJacobian_Hull1972A2; 1046c4762a1bSJed Brown IFunction = IFunction_Hull1972A2; 1047c4762a1bSJed Brown IJacobian = IJacobian_Hull1972A2; 1048c4762a1bSJed Brown } else if (!strcmp(p,"hull1972a3")) { 1049c4762a1bSJed Brown y[0] = 1.0; 1050c4762a1bSJed Brown RHSFunction = RHSFunction_Hull1972A3; 1051c4762a1bSJed Brown RHSJacobian = RHSJacobian_Hull1972A3; 1052c4762a1bSJed Brown IFunction = IFunction_Hull1972A3; 1053c4762a1bSJed Brown IJacobian = IJacobian_Hull1972A3; 1054c4762a1bSJed Brown } else if (!strcmp(p,"hull1972a4")) { 1055c4762a1bSJed Brown y[0] = 1.0; 1056c4762a1bSJed Brown RHSFunction = RHSFunction_Hull1972A4; 1057c4762a1bSJed Brown RHSJacobian = RHSJacobian_Hull1972A4; 1058c4762a1bSJed Brown IFunction = IFunction_Hull1972A4; 1059c4762a1bSJed Brown IJacobian = IJacobian_Hull1972A4; 1060c4762a1bSJed Brown } else if (!strcmp(p,"hull1972a5")) { 1061c4762a1bSJed Brown y[0] = 4.0; 1062c4762a1bSJed Brown RHSFunction = RHSFunction_Hull1972A5; 1063c4762a1bSJed Brown RHSJacobian = RHSJacobian_Hull1972A5; 1064c4762a1bSJed Brown IFunction = IFunction_Hull1972A5; 1065c4762a1bSJed Brown IJacobian = IJacobian_Hull1972A5; 1066c4762a1bSJed Brown } else if (!strcmp(p,"hull1972b1")) { 1067c4762a1bSJed Brown y[0] = 1.0; 1068c4762a1bSJed Brown y[1] = 3.0; 1069c4762a1bSJed Brown RHSFunction = RHSFunction_Hull1972B1; 1070c4762a1bSJed Brown IFunction = IFunction_Hull1972B1; 1071c4762a1bSJed Brown IJacobian = IJacobian_Hull1972B1; 1072c4762a1bSJed Brown } else if (!strcmp(p,"hull1972b2")) { 1073c4762a1bSJed Brown y[0] = 2.0; 1074c4762a1bSJed Brown y[1] = 0.0; 1075c4762a1bSJed Brown y[2] = 1.0; 1076c4762a1bSJed Brown RHSFunction = RHSFunction_Hull1972B2; 1077c4762a1bSJed Brown IFunction = IFunction_Hull1972B2; 1078c4762a1bSJed Brown IJacobian = IJacobian_Hull1972B2; 1079c4762a1bSJed Brown } else if (!strcmp(p,"hull1972b3")) { 1080c4762a1bSJed Brown y[0] = 1.0; 1081c4762a1bSJed Brown y[1] = 0.0; 1082c4762a1bSJed Brown y[2] = 0.0; 1083c4762a1bSJed Brown RHSFunction = RHSFunction_Hull1972B3; 1084c4762a1bSJed Brown IFunction = IFunction_Hull1972B3; 1085c4762a1bSJed Brown IJacobian = IJacobian_Hull1972B3; 1086c4762a1bSJed Brown } else if (!strcmp(p,"hull1972b4")) { 1087c4762a1bSJed Brown y[0] = 3.0; 1088c4762a1bSJed Brown y[1] = 0.0; 1089c4762a1bSJed Brown y[2] = 0.0; 1090c4762a1bSJed Brown RHSFunction = RHSFunction_Hull1972B4; 1091c4762a1bSJed Brown IFunction = IFunction_Hull1972B4; 1092c4762a1bSJed Brown IJacobian = IJacobian_Hull1972B4; 1093c4762a1bSJed Brown } else if (!strcmp(p,"hull1972b5")) { 1094c4762a1bSJed Brown y[0] = 0.0; 1095c4762a1bSJed Brown y[1] = 1.0; 1096c4762a1bSJed Brown y[2] = 1.0; 1097c4762a1bSJed Brown RHSFunction = RHSFunction_Hull1972B5; 1098c4762a1bSJed Brown IFunction = IFunction_Hull1972B5; 1099c4762a1bSJed Brown IJacobian = IJacobian_Hull1972B5; 1100c4762a1bSJed Brown } else if (!strcmp(p,"kulik2013i")) { 1101c4762a1bSJed Brown t0=0.; 1102c4762a1bSJed Brown y[0] = PetscExpReal(PetscSinReal(t0*t0)); 1103c4762a1bSJed Brown y[1] = PetscExpReal(5.*PetscSinReal(t0*t0)); 1104c4762a1bSJed Brown y[2] = PetscSinReal(t0*t0)+1.0; 1105c4762a1bSJed Brown y[3] = PetscCosReal(t0*t0); 1106c4762a1bSJed Brown RHSFunction = RHSFunction_Kulikov2013I; 1107c4762a1bSJed Brown RHSJacobian = RHSJacobian_Kulikov2013I; 1108c4762a1bSJed Brown IFunction = IFunction_Kulikov2013I; 1109c4762a1bSJed Brown IJacobian = IJacobian_Kulikov2013I; 1110c4762a1bSJed Brown } else if (!strcmp(p,"hull1972c1")) { 1111c4762a1bSJed Brown y[0] = 1.0; 1112c4762a1bSJed Brown RHSFunction = RHSFunction_Hull1972C1; 1113c4762a1bSJed Brown IFunction = IFunction_Hull1972C1; 1114c4762a1bSJed Brown IJacobian = IJacobian_Hull1972C1; 1115c4762a1bSJed Brown } else if (!strcmp(p,"hull1972c2")) { 1116c4762a1bSJed Brown y[0] = 1.0; 1117c4762a1bSJed Brown RHSFunction = RHSFunction_Hull1972C2; 1118c4762a1bSJed Brown IFunction = IFunction_Hull1972C2; 1119c4762a1bSJed Brown IJacobian = IJacobian_Hull1972C2; 1120c4762a1bSJed Brown } else if ((!strcmp(p,"hull1972c3")) 1121c4762a1bSJed Brown ||(!strcmp(p,"hull1972c4"))) { 1122c4762a1bSJed Brown y[0] = 1.0; 1123c4762a1bSJed Brown RHSFunction = RHSFunction_Hull1972C34; 1124c4762a1bSJed Brown IFunction = IFunction_Hull1972C34; 1125c4762a1bSJed Brown IJacobian = IJacobian_Hull1972C34; 1126c4762a1bSJed Brown } 1127c4762a1bSJed Brown ierr = PetscOptionsGetScalarArray(NULL,NULL,"-yinit",y,&N,&flg);CHKERRQ(ierr); 1128*3c633725SBarry Smith PetscCheck((N == GetSize((const char*)s)) || !flg,PETSC_COMM_WORLD,PETSC_ERR_ARG_SIZ,"Number of initial values %D does not match problem size %D.",N,GetSize((const char*)s)); 1129c4762a1bSJed Brown ierr = VecRestoreArray(Y,&y);CHKERRQ(ierr); 1130c4762a1bSJed Brown PetscFunctionReturn(0); 1131c4762a1bSJed Brown } 1132c4762a1bSJed Brown 1133c4762a1bSJed Brown /* Calculates the exact solution to problems that have one */ 1134c4762a1bSJed Brown PetscErrorCode ExactSolution(Vec Y, void* s, PetscReal t, PetscBool *flag) 1135c4762a1bSJed Brown { 1136c4762a1bSJed Brown PetscErrorCode ierr; 1137c4762a1bSJed Brown char *p = (char*) s; 1138c4762a1bSJed Brown PetscScalar *y; 1139c4762a1bSJed Brown 1140c4762a1bSJed Brown PetscFunctionBegin; 1141c4762a1bSJed Brown if (!strcmp(p,"hull1972a1")) { 1142c4762a1bSJed Brown ierr = VecGetArray(Y,&y);CHKERRQ(ierr); 1143c4762a1bSJed Brown y[0] = PetscExpReal(-t); 1144c4762a1bSJed Brown *flag = PETSC_TRUE; 1145c4762a1bSJed Brown ierr = VecRestoreArray(Y,&y);CHKERRQ(ierr); 1146c4762a1bSJed Brown } else if (!strcmp(p,"hull1972a2")) { 1147c4762a1bSJed Brown ierr = VecGetArray(Y,&y);CHKERRQ(ierr); 1148c4762a1bSJed Brown y[0] = 1.0/PetscSqrtReal(t+1); 1149c4762a1bSJed Brown *flag = PETSC_TRUE; 1150c4762a1bSJed Brown ierr = VecRestoreArray(Y,&y);CHKERRQ(ierr); 1151c4762a1bSJed Brown } else if (!strcmp(p,"hull1972a3")) { 1152c4762a1bSJed Brown ierr = VecGetArray(Y,&y);CHKERRQ(ierr); 1153c4762a1bSJed Brown y[0] = PetscExpReal(PetscSinReal(t)); 1154c4762a1bSJed Brown *flag = PETSC_TRUE; 1155c4762a1bSJed Brown ierr = VecRestoreArray(Y,&y);CHKERRQ(ierr); 1156c4762a1bSJed Brown } else if (!strcmp(p,"hull1972a4")) { 1157c4762a1bSJed Brown ierr = VecGetArray(Y,&y);CHKERRQ(ierr); 1158c4762a1bSJed Brown y[0] = 20.0/(1+19.0*PetscExpReal(-t/4.0)); 1159c4762a1bSJed Brown *flag = PETSC_TRUE; 1160c4762a1bSJed Brown ierr = VecRestoreArray(Y,&y);CHKERRQ(ierr); 1161c4762a1bSJed Brown } else if (!strcmp(p,"kulik2013i")) { 1162c4762a1bSJed Brown ierr = VecGetArray(Y,&y);CHKERRQ(ierr); 1163c4762a1bSJed Brown y[0] = PetscExpReal(PetscSinReal(t*t)); 1164c4762a1bSJed Brown y[1] = PetscExpReal(5.*PetscSinReal(t*t)); 1165c4762a1bSJed Brown y[2] = PetscSinReal(t*t)+1.0; 1166c4762a1bSJed Brown y[3] = PetscCosReal(t*t); 1167c4762a1bSJed Brown *flag = PETSC_TRUE; 1168c4762a1bSJed Brown ierr = VecRestoreArray(Y,&y);CHKERRQ(ierr); 1169c4762a1bSJed Brown } else { 1170c4762a1bSJed Brown ierr = VecSet(Y,0);CHKERRQ(ierr); 1171c4762a1bSJed Brown *flag = PETSC_FALSE; 1172c4762a1bSJed Brown } 1173c4762a1bSJed Brown PetscFunctionReturn(0); 1174c4762a1bSJed Brown } 1175c4762a1bSJed Brown 1176c4762a1bSJed Brown /* Solves the specified ODE and computes the error if exact solution is available */ 1177c4762a1bSJed Brown PetscErrorCode SolveODE(char* ptype, PetscReal dt, PetscReal tfinal, PetscInt maxiter, PetscReal *error, PetscBool *exact_flag) 1178c4762a1bSJed Brown { 1179c4762a1bSJed Brown PetscErrorCode ierr; /* Error code */ 1180c4762a1bSJed Brown TS ts; /* time-integrator */ 1181c4762a1bSJed Brown Vec Y; /* Solution vector */ 1182c4762a1bSJed Brown Vec Yex; /* Exact solution */ 1183c4762a1bSJed Brown PetscInt N; /* Size of the system of equations */ 1184c4762a1bSJed Brown TSType time_scheme; /* Type of time-integration scheme */ 1185c4762a1bSJed Brown Mat Jac = NULL; /* Jacobian matrix */ 1186c4762a1bSJed Brown Vec Yerr; /* Auxiliary solution vector */ 1187c4762a1bSJed Brown PetscReal err_norm; /* Estimated error norm */ 1188c4762a1bSJed Brown PetscReal final_time; /* Actual final time from the integrator */ 1189c4762a1bSJed Brown 1190c4762a1bSJed Brown PetscFunctionBegin; 1191c4762a1bSJed Brown N = GetSize((const char *)&ptype[0]); 1192*3c633725SBarry Smith PetscCheck(N >= 0,PETSC_COMM_WORLD,PETSC_ERR_ARG_SIZ,"Illegal problem specification."); 1193c4762a1bSJed Brown ierr = VecCreate(PETSC_COMM_WORLD,&Y);CHKERRQ(ierr); 1194c4762a1bSJed Brown ierr = VecSetSizes(Y,N,PETSC_DECIDE);CHKERRQ(ierr); 1195c4762a1bSJed Brown ierr = VecSetUp(Y);CHKERRQ(ierr); 1196c4762a1bSJed Brown ierr = VecSet(Y,0);CHKERRQ(ierr); 1197c4762a1bSJed Brown 1198c4762a1bSJed Brown /* Initialize the problem */ 1199c4762a1bSJed Brown ierr = Initialize(Y,&ptype[0]);CHKERRQ(ierr); 1200c4762a1bSJed Brown 1201c4762a1bSJed Brown /* Create and initialize the time-integrator */ 1202c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 1203c4762a1bSJed Brown /* Default time integration options */ 1204c4762a1bSJed Brown ierr = TSSetType(ts,TSRK);CHKERRQ(ierr); 1205c4762a1bSJed Brown ierr = TSSetMaxSteps(ts,maxiter);CHKERRQ(ierr); 1206c4762a1bSJed Brown ierr = TSSetMaxTime(ts,tfinal);CHKERRQ(ierr); 1207c4762a1bSJed Brown ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr); 1208c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr); 1209c4762a1bSJed Brown /* Read command line options for time integration */ 1210c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 1211c4762a1bSJed Brown /* Set solution vector */ 1212c4762a1bSJed Brown ierr = TSSetSolution(ts,Y);CHKERRQ(ierr); 1213c4762a1bSJed Brown /* Specify left/right-hand side functions */ 1214c4762a1bSJed Brown ierr = TSGetType(ts,&time_scheme);CHKERRQ(ierr); 1215c4762a1bSJed Brown 1216c4762a1bSJed Brown if ((!strcmp(time_scheme,TSEULER)) || (!strcmp(time_scheme,TSRK)) || (!strcmp(time_scheme,TSSSP) || (!strcmp(time_scheme,TSGLEE)))) { 1217c4762a1bSJed Brown /* Explicit time-integration -> specify right-hand side function ydot = f(y) */ 1218c4762a1bSJed Brown ierr = TSSetRHSFunction(ts,NULL,RHSFunction,&ptype[0]);CHKERRQ(ierr); 1219c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&Jac);CHKERRQ(ierr); 1220c4762a1bSJed Brown ierr = MatSetSizes(Jac,PETSC_DECIDE,PETSC_DECIDE,N,N);CHKERRQ(ierr); 1221c4762a1bSJed Brown ierr = MatSetFromOptions(Jac);CHKERRQ(ierr); 1222c4762a1bSJed Brown ierr = MatSetUp(Jac);CHKERRQ(ierr); 1223c4762a1bSJed Brown ierr = TSSetRHSJacobian(ts,Jac,Jac,RHSJacobian,&ptype[0]);CHKERRQ(ierr); 1224c4762a1bSJed Brown } else if ((!strcmp(time_scheme,TSTHETA)) || (!strcmp(time_scheme,TSBEULER)) || (!strcmp(time_scheme,TSCN)) || (!strcmp(time_scheme,TSALPHA)) || (!strcmp(time_scheme,TSARKIMEX))) { 1225c4762a1bSJed Brown /* Implicit time-integration -> specify left-hand side function ydot-f(y) = 0 */ 1226c4762a1bSJed Brown /* and its Jacobian function */ 1227c4762a1bSJed Brown ierr = TSSetIFunction(ts,NULL,IFunction,&ptype[0]);CHKERRQ(ierr); 1228c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&Jac);CHKERRQ(ierr); 1229c4762a1bSJed Brown ierr = MatSetSizes(Jac,PETSC_DECIDE,PETSC_DECIDE,N,N);CHKERRQ(ierr); 1230c4762a1bSJed Brown ierr = MatSetFromOptions(Jac);CHKERRQ(ierr); 1231c4762a1bSJed Brown ierr = MatSetUp(Jac);CHKERRQ(ierr); 1232c4762a1bSJed Brown ierr = TSSetIJacobian(ts,Jac,Jac,IJacobian,&ptype[0]);CHKERRQ(ierr); 1233c4762a1bSJed Brown } 1234c4762a1bSJed Brown 1235c4762a1bSJed Brown /* Solve */ 1236c4762a1bSJed Brown ierr = TSSolve(ts,Y);CHKERRQ(ierr); 1237c4762a1bSJed Brown ierr = TSGetTime(ts,&final_time);CHKERRQ(ierr); 1238c4762a1bSJed Brown 1239c4762a1bSJed Brown /* Get the estimated error, if available */ 1240c4762a1bSJed Brown ierr = VecDuplicate(Y,&Yerr);CHKERRQ(ierr); 1241c4762a1bSJed Brown ierr = VecZeroEntries(Yerr);CHKERRQ(ierr); 1242c4762a1bSJed Brown ierr = TSGetTimeError(ts,0,&Yerr);CHKERRQ(ierr); 1243c4762a1bSJed Brown ierr = VecNorm(Yerr,NORM_2,&err_norm);CHKERRQ(ierr); 1244c4762a1bSJed Brown ierr = VecDestroy(&Yerr);CHKERRQ(ierr); 1245c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Estimated Error = %E.\n",err_norm);CHKERRQ(ierr); 1246c4762a1bSJed Brown 1247c4762a1bSJed Brown /* Exact solution */ 1248c4762a1bSJed Brown ierr = VecDuplicate(Y,&Yex);CHKERRQ(ierr); 1249c4762a1bSJed Brown if (PetscAbsScalar(final_time-tfinal)>2.*PETSC_MACHINE_EPSILON) { 1250c4762a1bSJed Brown ierr = 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);CHKERRQ(ierr); 1251c4762a1bSJed Brown } 1252c4762a1bSJed Brown ierr = ExactSolution(Yex,&ptype[0],final_time,exact_flag);CHKERRQ(ierr); 1253c4762a1bSJed Brown 1254c4762a1bSJed Brown /* Calculate Error */ 1255c4762a1bSJed Brown ierr = VecAYPX(Yex,-1.0,Y);CHKERRQ(ierr); 1256c4762a1bSJed Brown ierr = VecNorm(Yex,NORM_2,error);CHKERRQ(ierr); 1257c4762a1bSJed Brown *error = PetscSqrtReal(((*error)*(*error))/N); 1258c4762a1bSJed Brown 1259c4762a1bSJed Brown /* Clean up and finalize */ 1260c4762a1bSJed Brown ierr = MatDestroy(&Jac);CHKERRQ(ierr); 1261c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 1262c4762a1bSJed Brown ierr = VecDestroy(&Yex);CHKERRQ(ierr); 1263c4762a1bSJed Brown ierr = VecDestroy(&Y);CHKERRQ(ierr); 1264c4762a1bSJed Brown 1265c4762a1bSJed Brown PetscFunctionReturn(0); 1266c4762a1bSJed Brown } 1267c4762a1bSJed Brown 1268c4762a1bSJed Brown int main(int argc, char **argv) 1269c4762a1bSJed Brown { 1270c4762a1bSJed Brown PetscErrorCode ierr; /* Error code */ 1271c4762a1bSJed Brown char ptype[256] = "hull1972a1"; /* Problem specification */ 1272c4762a1bSJed Brown PetscInt n_refine = 1; /* Number of refinement levels for convergence analysis */ 1273c4762a1bSJed Brown PetscReal refine_fac = 2.0; /* Refinement factor for dt */ 1274c4762a1bSJed Brown PetscReal dt_initial = 0.01; /* Initial default value of dt */ 1275c4762a1bSJed Brown PetscReal dt; 1276c4762a1bSJed Brown PetscReal tfinal = 20.0; /* Final time for the time-integration */ 1277c4762a1bSJed Brown PetscInt maxiter = 100000; /* Maximum number of time-integration iterations */ 1278c4762a1bSJed Brown PetscReal *error; /* Array to store the errors for convergence analysis */ 1279c4762a1bSJed Brown PetscMPIInt size; /* No of processors */ 1280c4762a1bSJed Brown PetscBool flag; /* Flag denoting availability of exact solution */ 1281c4762a1bSJed Brown PetscInt r; 1282c4762a1bSJed Brown 1283c4762a1bSJed Brown /* Initialize program */ 1284c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 1285c4762a1bSJed Brown 1286c4762a1bSJed Brown /* Check if running with only 1 proc */ 1287ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 1288*3c633725SBarry Smith PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Only for sequential runs"); 1289c4762a1bSJed Brown 129076280437SVaclav Hapla ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"ex31",NULL);CHKERRQ(ierr); 1291c4762a1bSJed Brown ierr = PetscOptionsString("-problem","Problem specification","<hull1972a1>",ptype,ptype,sizeof(ptype),NULL);CHKERRQ(ierr); 1292c4762a1bSJed Brown ierr = PetscOptionsInt("-refinement_levels","Number of refinement levels for convergence analysis","<1>",n_refine,&n_refine,NULL);CHKERRQ(ierr); 1293c4762a1bSJed Brown ierr = PetscOptionsReal("-refinement_factor","Refinement factor for dt","<2.0>",refine_fac,&refine_fac,NULL);CHKERRQ(ierr); 1294c4762a1bSJed Brown ierr = PetscOptionsReal("-dt","Time step size (for convergence analysis, initial time step)","<0.01>",dt_initial,&dt_initial,NULL);CHKERRQ(ierr); 1295c4762a1bSJed Brown ierr = PetscOptionsReal("-final_time","Final time for the time-integration","<20.0>",tfinal,&tfinal,NULL);CHKERRQ(ierr); 1296c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 1297c4762a1bSJed Brown 1298c4762a1bSJed Brown ierr = PetscMalloc1(n_refine,&error);CHKERRQ(ierr); 1299c4762a1bSJed Brown for (r = 0,dt = dt_initial; r < n_refine; r++) { 1300c4762a1bSJed Brown error[r] = 0; 1301c4762a1bSJed Brown if (r > 0) dt /= refine_fac; 1302c4762a1bSJed Brown 1303c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Solving ODE \"%s\" with dt %f, final time %f and system size %D.\n",ptype,(double)dt,(double)tfinal,GetSize(&ptype[0]));CHKERRQ(ierr); 1304c4762a1bSJed Brown ierr = SolveODE(&ptype[0],dt,tfinal,maxiter,&error[r],&flag);CHKERRQ(ierr); 1305c4762a1bSJed Brown if (flag) { 1306c4762a1bSJed Brown /* If exact solution available for the specified ODE */ 1307c4762a1bSJed Brown if (r > 0) { 1308c4762a1bSJed Brown PetscReal conv_rate = (PetscLogReal(error[r]) - PetscLogReal(error[r-1])) / (-PetscLogReal(refine_fac)); 1309c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Error = %E,\tConvergence rate = %f.\n",(double)error[r],(double)conv_rate);CHKERRQ(ierr); 1310c4762a1bSJed Brown } else { 1311c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Error = %E.\n",error[r]);CHKERRQ(ierr); 1312c4762a1bSJed Brown } 1313c4762a1bSJed Brown } 1314c4762a1bSJed Brown } 1315c4762a1bSJed Brown ierr = PetscFree(error);CHKERRQ(ierr); 1316c4762a1bSJed Brown ierr = PetscFinalize(); 1317c4762a1bSJed Brown return ierr; 1318c4762a1bSJed Brown } 1319c4762a1bSJed Brown 1320c4762a1bSJed Brown /*TEST 1321c4762a1bSJed Brown 1322c4762a1bSJed Brown test: 1323c4762a1bSJed Brown suffix: 2 1324c4762a1bSJed Brown args: -ts_type glee -final_time 5 -ts_adapt_type none 1325c4762a1bSJed Brown timeoutfactor: 3 1326c4762a1bSJed Brown requires: !single 1327c4762a1bSJed Brown 1328c4762a1bSJed Brown test: 1329c4762a1bSJed Brown suffix: 3 1330c4762a1bSJed 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 1331c4762a1bSJed Brown timeoutfactor: 3 1332c4762a1bSJed Brown requires: !single 1333c4762a1bSJed Brown 1334c4762a1bSJed Brown test: 1335c4762a1bSJed Brown suffix: 4 1336c4762a1bSJed 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 1337c4762a1bSJed Brown timeoutfactor: 3 1338c4762a1bSJed Brown requires: !single !__float128 1339c4762a1bSJed Brown 1340c4762a1bSJed Brown TEST*/ 1341