1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Newton's method for a two-variable system, sequential.\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown /*T 5c4762a1bSJed Brown Concepts: SNES^basic example 6c4762a1bSJed Brown T*/ 7c4762a1bSJed Brown 8c4762a1bSJed Brown 9c4762a1bSJed Brown 10c4762a1bSJed Brown /* 11c4762a1bSJed Brown Include "petscsnes.h" so that we can use SNES solvers. Note that this 12c4762a1bSJed Brown file automatically includes: 13c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 14c4762a1bSJed Brown petscmat.h - matrices 15c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods 16c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners 17c4762a1bSJed Brown petscksp.h - linear solvers 18c4762a1bSJed Brown */ 19c4762a1bSJed Brown /*F 20c4762a1bSJed Brown This examples solves either 21c4762a1bSJed Brown \begin{equation} 22c4762a1bSJed Brown F\genfrac{(}{)}{0pt}{}{x_0}{x_1} = \genfrac{(}{)}{0pt}{}{x^2_0 + x_0 x_1 - 3}{x_0 x_1 + x^2_1 - 6} 23c4762a1bSJed Brown \end{equation} 24c4762a1bSJed Brown or if the {\tt -hard} options is given 25c4762a1bSJed Brown \begin{equation} 26c4762a1bSJed Brown F\genfrac{(}{)}{0pt}{}{x_0}{x_1} = \genfrac{(}{)}{0pt}{}{\sin(3 x_0) + x_0}{x_1} 27c4762a1bSJed Brown \end{equation} 28c4762a1bSJed Brown F*/ 29c4762a1bSJed Brown #include <petscsnes.h> 30c4762a1bSJed Brown 31c4762a1bSJed Brown /* 32c4762a1bSJed Brown User-defined routines 33c4762a1bSJed Brown */ 34c4762a1bSJed Brown extern PetscErrorCode FormJacobian1(SNES,Vec,Mat,Mat,void*); 35c4762a1bSJed Brown extern PetscErrorCode FormFunction1(SNES,Vec,Vec,void*); 36c4762a1bSJed Brown extern PetscErrorCode FormJacobian2(SNES,Vec,Mat,Mat,void*); 37c4762a1bSJed Brown extern PetscErrorCode FormFunction2(SNES,Vec,Vec,void*); 38c4762a1bSJed Brown 39c4762a1bSJed Brown int main(int argc,char **argv) 40c4762a1bSJed Brown { 41c4762a1bSJed Brown SNES snes; /* nonlinear solver context */ 42c4762a1bSJed Brown KSP ksp; /* linear solver context */ 43c4762a1bSJed Brown PC pc; /* preconditioner context */ 44c4762a1bSJed Brown Vec x,r; /* solution, residual vectors */ 45c4762a1bSJed Brown Mat J; /* Jacobian matrix */ 46c4762a1bSJed Brown PetscErrorCode ierr; 47c4762a1bSJed Brown PetscMPIInt size; 48c4762a1bSJed Brown PetscScalar pfive = .5,*xx; 49c4762a1bSJed Brown PetscBool flg; 50c4762a1bSJed Brown 51c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 52*ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 53c4762a1bSJed Brown if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Example is only for sequential runs"); 54c4762a1bSJed Brown 55c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 56c4762a1bSJed Brown Create nonlinear solver context 57c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 58c4762a1bSJed Brown ierr = SNESCreate(PETSC_COMM_WORLD,&snes);CHKERRQ(ierr); 59c4762a1bSJed Brown 60c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 61c4762a1bSJed Brown Create matrix and vector data structures; set corresponding routines 62c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 63c4762a1bSJed Brown /* 64c4762a1bSJed Brown Create vectors for solution and nonlinear function 65c4762a1bSJed Brown */ 66c4762a1bSJed Brown ierr = VecCreate(PETSC_COMM_WORLD,&x);CHKERRQ(ierr); 67c4762a1bSJed Brown ierr = VecSetSizes(x,PETSC_DECIDE,2);CHKERRQ(ierr); 68c4762a1bSJed Brown ierr = VecSetFromOptions(x);CHKERRQ(ierr); 69c4762a1bSJed Brown ierr = VecDuplicate(x,&r);CHKERRQ(ierr); 70c4762a1bSJed Brown 71c4762a1bSJed Brown /* 72c4762a1bSJed Brown Create Jacobian matrix data structure 73c4762a1bSJed Brown */ 74c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&J);CHKERRQ(ierr); 75c4762a1bSJed Brown ierr = MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,2,2);CHKERRQ(ierr); 76c4762a1bSJed Brown ierr = MatSetFromOptions(J);CHKERRQ(ierr); 77c4762a1bSJed Brown ierr = MatSetUp(J);CHKERRQ(ierr); 78c4762a1bSJed Brown 79c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL,"-hard",&flg);CHKERRQ(ierr); 80c4762a1bSJed Brown if (!flg) { 81c4762a1bSJed Brown /* 82c4762a1bSJed Brown Set function evaluation routine and vector. 83c4762a1bSJed Brown */ 84c4762a1bSJed Brown ierr = SNESSetFunction(snes,r,FormFunction1,NULL);CHKERRQ(ierr); 85c4762a1bSJed Brown 86c4762a1bSJed Brown /* 87c4762a1bSJed Brown Set Jacobian matrix data structure and Jacobian evaluation routine 88c4762a1bSJed Brown */ 89c4762a1bSJed Brown ierr = SNESSetJacobian(snes,J,J,FormJacobian1,NULL);CHKERRQ(ierr); 90c4762a1bSJed Brown } else { 91c4762a1bSJed Brown ierr = SNESSetFunction(snes,r,FormFunction2,NULL);CHKERRQ(ierr); 92c4762a1bSJed Brown ierr = SNESSetJacobian(snes,J,J,FormJacobian2,NULL);CHKERRQ(ierr); 93c4762a1bSJed Brown } 94c4762a1bSJed Brown 95c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 96c4762a1bSJed Brown Customize nonlinear solver; set runtime options 97c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 98c4762a1bSJed Brown /* 99c4762a1bSJed Brown Set linear solver defaults for this problem. By extracting the 100c4762a1bSJed Brown KSP and PC contexts from the SNES context, we can then 101c4762a1bSJed Brown directly call any KSP and PC routines to set various options. 102c4762a1bSJed Brown */ 103c4762a1bSJed Brown ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 104c4762a1bSJed Brown ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 105c4762a1bSJed Brown ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 106c4762a1bSJed Brown ierr = KSPSetTolerances(ksp,1.e-4,PETSC_DEFAULT,PETSC_DEFAULT,20);CHKERRQ(ierr); 107c4762a1bSJed Brown 108c4762a1bSJed Brown /* 109c4762a1bSJed Brown Set SNES/KSP/KSP/PC runtime options, e.g., 110c4762a1bSJed Brown -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc> 111c4762a1bSJed Brown These options will override those specified above as long as 112c4762a1bSJed Brown SNESSetFromOptions() is called _after_ any other customization 113c4762a1bSJed Brown routines. 114c4762a1bSJed Brown */ 115c4762a1bSJed Brown ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 116c4762a1bSJed Brown 117c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 118c4762a1bSJed Brown Evaluate initial guess; then solve nonlinear system 119c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 120c4762a1bSJed Brown if (!flg) { 121c4762a1bSJed Brown ierr = VecSet(x,pfive);CHKERRQ(ierr); 122c4762a1bSJed Brown } else { 123c4762a1bSJed Brown ierr = VecGetArray(x,&xx);CHKERRQ(ierr); 124c4762a1bSJed Brown xx[0] = 2.0; xx[1] = 3.0; 125c4762a1bSJed Brown ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr); 126c4762a1bSJed Brown } 127c4762a1bSJed Brown /* 128c4762a1bSJed Brown Note: The user should initialize the vector, x, with the initial guess 129c4762a1bSJed Brown for the nonlinear solver prior to calling SNESSolve(). In particular, 130c4762a1bSJed Brown to employ an initial guess of zero, the user should explicitly set 131c4762a1bSJed Brown this vector to zero by calling VecSet(). 132c4762a1bSJed Brown */ 133c4762a1bSJed Brown 134c4762a1bSJed Brown ierr = SNESSolve(snes,NULL,x);CHKERRQ(ierr); 135c4762a1bSJed Brown if (flg) { 136c4762a1bSJed Brown Vec f; 137c4762a1bSJed Brown ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 138c4762a1bSJed Brown ierr = SNESGetFunction(snes,&f,0,0);CHKERRQ(ierr); 139c4762a1bSJed Brown ierr = VecView(r,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 140c4762a1bSJed Brown } 141c4762a1bSJed Brown 142c4762a1bSJed Brown 143c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 144c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 145c4762a1bSJed Brown are no longer needed. 146c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 147c4762a1bSJed Brown 148c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); ierr = VecDestroy(&r);CHKERRQ(ierr); 149c4762a1bSJed Brown ierr = MatDestroy(&J);CHKERRQ(ierr); ierr = SNESDestroy(&snes);CHKERRQ(ierr); 150c4762a1bSJed Brown ierr = PetscFinalize(); 151c4762a1bSJed Brown return ierr; 152c4762a1bSJed Brown } 153c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 154c4762a1bSJed Brown /* 155c4762a1bSJed Brown FormFunction1 - Evaluates nonlinear function, F(x). 156c4762a1bSJed Brown 157c4762a1bSJed Brown Input Parameters: 158c4762a1bSJed Brown . snes - the SNES context 159c4762a1bSJed Brown . x - input vector 160c4762a1bSJed Brown . ctx - optional user-defined context 161c4762a1bSJed Brown 162c4762a1bSJed Brown Output Parameter: 163c4762a1bSJed Brown . f - function vector 164c4762a1bSJed Brown */ 165c4762a1bSJed Brown PetscErrorCode FormFunction1(SNES snes,Vec x,Vec f,void *ctx) 166c4762a1bSJed Brown { 167c4762a1bSJed Brown PetscErrorCode ierr; 168c4762a1bSJed Brown const PetscScalar *xx; 169c4762a1bSJed Brown PetscScalar *ff; 170c4762a1bSJed Brown 171c4762a1bSJed Brown /* 172c4762a1bSJed Brown Get pointers to vector data. 173c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 174c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 175c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 176c4762a1bSJed Brown the array. 177c4762a1bSJed Brown */ 178c4762a1bSJed Brown ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr); 179c4762a1bSJed Brown ierr = VecGetArray(f,&ff);CHKERRQ(ierr); 180c4762a1bSJed Brown 181c4762a1bSJed Brown /* Compute function */ 182c4762a1bSJed Brown ff[0] = xx[0]*xx[0] + xx[0]*xx[1] - 3.0; 183c4762a1bSJed Brown ff[1] = xx[0]*xx[1] + xx[1]*xx[1] - 6.0; 184c4762a1bSJed Brown 185c4762a1bSJed Brown /* Restore vectors */ 186c4762a1bSJed Brown ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr); 187c4762a1bSJed Brown ierr = VecRestoreArray(f,&ff);CHKERRQ(ierr); 188c4762a1bSJed Brown return 0; 189c4762a1bSJed Brown } 190c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 191c4762a1bSJed Brown /* 192c4762a1bSJed Brown FormJacobian1 - Evaluates Jacobian matrix. 193c4762a1bSJed Brown 194c4762a1bSJed Brown Input Parameters: 195c4762a1bSJed Brown . snes - the SNES context 196c4762a1bSJed Brown . x - input vector 197c4762a1bSJed Brown . dummy - optional user-defined context (not used here) 198c4762a1bSJed Brown 199c4762a1bSJed Brown Output Parameters: 200c4762a1bSJed Brown . jac - Jacobian matrix 201c4762a1bSJed Brown . B - optionally different preconditioning matrix 202c4762a1bSJed Brown . flag - flag indicating matrix structure 203c4762a1bSJed Brown */ 204c4762a1bSJed Brown PetscErrorCode FormJacobian1(SNES snes,Vec x,Mat jac,Mat B,void *dummy) 205c4762a1bSJed Brown { 206c4762a1bSJed Brown const PetscScalar *xx; 207c4762a1bSJed Brown PetscScalar A[4]; 208c4762a1bSJed Brown PetscErrorCode ierr; 209c4762a1bSJed Brown PetscInt idx[2] = {0,1}; 210c4762a1bSJed Brown 211c4762a1bSJed Brown /* 212c4762a1bSJed Brown Get pointer to vector data 213c4762a1bSJed Brown */ 214c4762a1bSJed Brown ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr); 215c4762a1bSJed Brown 216c4762a1bSJed Brown /* 217c4762a1bSJed Brown Compute Jacobian entries and insert into matrix. 218c4762a1bSJed Brown - Since this is such a small problem, we set all entries for 219c4762a1bSJed Brown the matrix at once. 220c4762a1bSJed Brown */ 221c4762a1bSJed Brown A[0] = 2.0*xx[0] + xx[1]; A[1] = xx[0]; 222c4762a1bSJed Brown A[2] = xx[1]; A[3] = xx[0] + 2.0*xx[1]; 223c4762a1bSJed Brown ierr = MatSetValues(B,2,idx,2,idx,A,INSERT_VALUES);CHKERRQ(ierr); 224c4762a1bSJed Brown 225c4762a1bSJed Brown /* 226c4762a1bSJed Brown Restore vector 227c4762a1bSJed Brown */ 228c4762a1bSJed Brown ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr); 229c4762a1bSJed Brown 230c4762a1bSJed Brown /* 231c4762a1bSJed Brown Assemble matrix 232c4762a1bSJed Brown */ 233c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 234c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 235c4762a1bSJed Brown if (jac != B) { 236c4762a1bSJed Brown ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 237c4762a1bSJed Brown ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 238c4762a1bSJed Brown } 239c4762a1bSJed Brown return 0; 240c4762a1bSJed Brown } 241c4762a1bSJed Brown 242c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 243c4762a1bSJed Brown PetscErrorCode FormFunction2(SNES snes,Vec x,Vec f,void *dummy) 244c4762a1bSJed Brown { 245c4762a1bSJed Brown PetscErrorCode ierr; 246c4762a1bSJed Brown const PetscScalar *xx; 247c4762a1bSJed Brown PetscScalar *ff; 248c4762a1bSJed Brown 249c4762a1bSJed Brown /* 250c4762a1bSJed Brown Get pointers to vector data. 251c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 252c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 253c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 254c4762a1bSJed Brown the array. 255c4762a1bSJed Brown */ 256c4762a1bSJed Brown ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr); 257c4762a1bSJed Brown ierr = VecGetArray(f,&ff);CHKERRQ(ierr); 258c4762a1bSJed Brown 259c4762a1bSJed Brown /* 260c4762a1bSJed Brown Compute function 261c4762a1bSJed Brown */ 262c4762a1bSJed Brown ff[0] = PetscSinScalar(3.0*xx[0]) + xx[0]; 263c4762a1bSJed Brown ff[1] = xx[1]; 264c4762a1bSJed Brown 265c4762a1bSJed Brown /* 266c4762a1bSJed Brown Restore vectors 267c4762a1bSJed Brown */ 268c4762a1bSJed Brown ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr); 269c4762a1bSJed Brown ierr = VecRestoreArray(f,&ff);CHKERRQ(ierr); 270c4762a1bSJed Brown return 0; 271c4762a1bSJed Brown } 272c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 273c4762a1bSJed Brown PetscErrorCode FormJacobian2(SNES snes,Vec x,Mat jac,Mat B,void *dummy) 274c4762a1bSJed Brown { 275c4762a1bSJed Brown const PetscScalar *xx; 276c4762a1bSJed Brown PetscScalar A[4]; 277c4762a1bSJed Brown PetscErrorCode ierr; 278c4762a1bSJed Brown PetscInt idx[2] = {0,1}; 279c4762a1bSJed Brown 280c4762a1bSJed Brown /* 281c4762a1bSJed Brown Get pointer to vector data 282c4762a1bSJed Brown */ 283c4762a1bSJed Brown ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr); 284c4762a1bSJed Brown 285c4762a1bSJed Brown /* 286c4762a1bSJed Brown Compute Jacobian entries and insert into matrix. 287c4762a1bSJed Brown - Since this is such a small problem, we set all entries for 288c4762a1bSJed Brown the matrix at once. 289c4762a1bSJed Brown */ 290c4762a1bSJed Brown A[0] = 3.0*PetscCosScalar(3.0*xx[0]) + 1.0; A[1] = 0.0; 291c4762a1bSJed Brown A[2] = 0.0; A[3] = 1.0; 292c4762a1bSJed Brown ierr = MatSetValues(B,2,idx,2,idx,A,INSERT_VALUES);CHKERRQ(ierr); 293c4762a1bSJed Brown 294c4762a1bSJed Brown /* 295c4762a1bSJed Brown Restore vector 296c4762a1bSJed Brown */ 297c4762a1bSJed Brown ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr); 298c4762a1bSJed Brown 299c4762a1bSJed Brown /* 300c4762a1bSJed Brown Assemble matrix 301c4762a1bSJed Brown */ 302c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 303c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 304c4762a1bSJed Brown if (jac != B) { 305c4762a1bSJed Brown ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 306c4762a1bSJed Brown ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 307c4762a1bSJed Brown } 308c4762a1bSJed Brown return 0; 309c4762a1bSJed Brown } 310c4762a1bSJed Brown 311c4762a1bSJed Brown 312c4762a1bSJed Brown 313c4762a1bSJed Brown 314c4762a1bSJed Brown /*TEST 315c4762a1bSJed Brown 316c4762a1bSJed Brown test: 317c4762a1bSJed Brown args: -ksp_gmres_cgs_refinement_type refine_always -snes_monitor_short 318c4762a1bSJed Brown requires: !single 319c4762a1bSJed Brown 320c4762a1bSJed Brown test: 321c4762a1bSJed Brown suffix: 2 322c4762a1bSJed Brown requires: !single 323c4762a1bSJed Brown args: -snes_monitor_short 324c4762a1bSJed Brown output_file: output/ex1_1.out 325c4762a1bSJed Brown 326c4762a1bSJed Brown test: 327c4762a1bSJed Brown suffix: 3 328c4762a1bSJed Brown args: -ksp_view_solution ascii:ex1_2_sol.tmp:ascii_matlab -snes_monitor_short 329c4762a1bSJed Brown requires: !single 330c4762a1bSJed Brown output_file: output/ex1_1.out 331c4762a1bSJed Brown 332c4762a1bSJed Brown test: 333c4762a1bSJed Brown suffix: 4 334c4762a1bSJed Brown args: -ksp_view_solution ascii:ex1_2_sol.tmp::append -snes_monitor_short 335c4762a1bSJed Brown requires: !single 336c4762a1bSJed Brown output_file: output/ex1_1.out 337c4762a1bSJed Brown 338c4762a1bSJed Brown test: 339c4762a1bSJed Brown suffix: 5 340c4762a1bSJed Brown args: -ksp_view_solution ascii:ex1_2_sol.tmp:ascii_matlab:append -snes_monitor_short 341c4762a1bSJed Brown requires: !single 342c4762a1bSJed Brown output_file: output/ex1_1.out 343c4762a1bSJed Brown 344c4762a1bSJed Brown test: 345c4762a1bSJed Brown suffix: 6 346c4762a1bSJed Brown args: -ksp_view_solution ascii:ex1_2_sol.tmp:default:append -snes_monitor_short 347c4762a1bSJed Brown requires: !single 348c4762a1bSJed Brown output_file: output/ex1_1.out 349c4762a1bSJed Brown 350c4762a1bSJed Brown test: 351c4762a1bSJed Brown suffix: X 352c4762a1bSJed Brown args: -ksp_monitor_short -ksp_type gmres -ksp_gmres_krylov_monitor -snes_monitor_short -snes_rtol 1.e-4 353c4762a1bSJed Brown requires: !single x 354c4762a1bSJed Brown 355c4762a1bSJed Brown TEST*/ 356