1c4762a1bSJed Brown static const char help[] = "Newton's method to solve a two-variable system, sequentially.\n" 2c4762a1bSJed Brown "The same problem is solved twice - i) fully assembled system + ii) block system\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown /*T 5c4762a1bSJed Brown Concepts: SNES^basic uniprocessor example, block objects 6c4762a1bSJed Brown Processors: 1 7c4762a1bSJed Brown T*/ 8c4762a1bSJed Brown 9c4762a1bSJed Brown /* 10c4762a1bSJed Brown Include "petscsnes.h" so that we can use SNES solvers. Note that this 11c4762a1bSJed Brown file automatically includes: 12c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 13c4762a1bSJed Brown petscsys.h - system routines petscmat.h - matrices 14c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods 15c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners 16c4762a1bSJed Brown petscksp.h - linear solvers 17c4762a1bSJed Brown */ 18c4762a1bSJed Brown #include <petscsnes.h> 19c4762a1bSJed Brown 20c4762a1bSJed Brown /* 21c4762a1bSJed Brown This example is block version of the test found at 22c4762a1bSJed Brown ${PETSC_DIR}/src/snes/tutorials/ex1.c 23c4762a1bSJed Brown In this test we replace the Jacobian systems 24c4762a1bSJed Brown [J]{x} = {F} 25c4762a1bSJed Brown where 26c4762a1bSJed Brown 27c4762a1bSJed Brown [J] = (j_00, j_01), {x} = (x_0, x_1)^T, {F} = (f_0, f_1)^T 28c4762a1bSJed Brown (j_10, j_11) 29c4762a1bSJed Brown where [J] \in \mathbb^{2 \times 2}, {x},{F} \in \mathbb^{2 \times 1}, 30c4762a1bSJed Brown 31c4762a1bSJed Brown with a block system in which each block is of length 1. 32c4762a1bSJed Brown i.e. The block system is thus 33c4762a1bSJed Brown 34c4762a1bSJed Brown [J] = ([j00], [j01]), {x} = ({x0}, {x1})^T, {F} = ({f0}, {f1})^T 35c4762a1bSJed Brown ([j10], [j11]) 36c4762a1bSJed Brown where 37c4762a1bSJed Brown [j00], [j01], [j10], [j11] \in \mathbb^{1 \times 1} 38c4762a1bSJed Brown {x0}, {x1}, {f0}, {f1} \in \mathbb^{1 \times 1} 39c4762a1bSJed Brown 40c4762a1bSJed Brown In practice we would not bother defing blocks of size one, and would instead assemble the 41c4762a1bSJed Brown full system. This is just a simple test to illustrate how to manipulate the blocks and 42c4762a1bSJed Brown to confirm the implementation is correct. 43c4762a1bSJed Brown */ 44c4762a1bSJed Brown 45c4762a1bSJed Brown /* 46c4762a1bSJed Brown User-defined routines 47c4762a1bSJed Brown */ 48c4762a1bSJed Brown static PetscErrorCode FormJacobian1(SNES,Vec,Mat,Mat,void*); 49c4762a1bSJed Brown static PetscErrorCode FormFunction1(SNES,Vec,Vec,void*); 50c4762a1bSJed Brown static PetscErrorCode FormJacobian2(SNES,Vec,Mat,Mat,void*); 51c4762a1bSJed Brown static PetscErrorCode FormFunction2(SNES,Vec,Vec,void*); 52c4762a1bSJed Brown static PetscErrorCode FormJacobian1_block(SNES,Vec,Mat,Mat,void*); 53c4762a1bSJed Brown static PetscErrorCode FormFunction1_block(SNES,Vec,Vec,void*); 54c4762a1bSJed Brown static PetscErrorCode FormJacobian2_block(SNES,Vec,Mat,Mat,void*); 55c4762a1bSJed Brown static PetscErrorCode FormFunction2_block(SNES,Vec,Vec,void*); 56c4762a1bSJed Brown 57c4762a1bSJed Brown static PetscErrorCode assembled_system(void) 58c4762a1bSJed Brown { 59c4762a1bSJed Brown SNES snes; /* nonlinear solver context */ 60c4762a1bSJed Brown KSP ksp; /* linear solver context */ 61c4762a1bSJed Brown PC pc; /* preconditioner context */ 62c4762a1bSJed Brown Vec x,r; /* solution, residual vectors */ 63c4762a1bSJed Brown Mat J; /* Jacobian matrix */ 64c4762a1bSJed Brown PetscInt its; 65c4762a1bSJed Brown PetscScalar pfive = .5,*xx; 66c4762a1bSJed Brown PetscBool flg; 67c4762a1bSJed Brown 68c4762a1bSJed Brown PetscFunctionBeginUser; 695f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "\n\n========================= Assembled system =========================\n\n")); 70c4762a1bSJed Brown 71c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 72c4762a1bSJed Brown Create nonlinear solver context 73c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 74c4762a1bSJed Brown 755f80ce2aSJacob Faibussowitsch CHKERRQ(SNESCreate(PETSC_COMM_WORLD,&snes)); 76c4762a1bSJed Brown 77c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 78c4762a1bSJed Brown Create matrix and vector data structures; set corresponding routines 79c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 80c4762a1bSJed Brown 81c4762a1bSJed Brown /* 82c4762a1bSJed Brown Create vectors for solution and nonlinear function 83c4762a1bSJed Brown */ 845f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,2,&x)); 855f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&r)); 86c4762a1bSJed Brown 87c4762a1bSJed Brown /* 88c4762a1bSJed Brown Create Jacobian matrix data structure 89c4762a1bSJed Brown */ 905f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_SELF,&J)); 915f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,2,2)); 925f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(J)); 935f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(J)); 94c4762a1bSJed Brown 955f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-hard",&flg)); 96c4762a1bSJed Brown if (!flg) { 97c4762a1bSJed Brown /* 98c4762a1bSJed Brown Set function evaluation routine and vector. 99c4762a1bSJed Brown */ 1005f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetFunction(snes,r,FormFunction1,NULL)); 101c4762a1bSJed Brown 102c4762a1bSJed Brown /* 103c4762a1bSJed Brown Set Jacobian matrix data structure and Jacobian evaluation routine 104c4762a1bSJed Brown */ 1055f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetJacobian(snes,J,J,FormJacobian1,NULL)); 106c4762a1bSJed Brown } else { 1075f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetFunction(snes,r,FormFunction2,NULL)); 1085f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetJacobian(snes,J,J,FormJacobian2,NULL)); 109c4762a1bSJed Brown } 110c4762a1bSJed Brown 111c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 112c4762a1bSJed Brown Customize nonlinear solver; set runtime options 113c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 114c4762a1bSJed Brown 115c4762a1bSJed Brown /* 116c4762a1bSJed Brown Set linear solver defaults for this problem. By extracting the 117c4762a1bSJed Brown KSP, KSP, and PC contexts from the SNES context, we can then 118c4762a1bSJed Brown directly call any KSP, KSP, and PC routines to set various options. 119c4762a1bSJed Brown */ 1205f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetKSP(snes,&ksp)); 1215f80ce2aSJacob Faibussowitsch CHKERRQ(KSPGetPC(ksp,&pc)); 1225f80ce2aSJacob Faibussowitsch CHKERRQ(PCSetType(pc,PCNONE)); 1235f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetTolerances(ksp,1.e-4,PETSC_DEFAULT,PETSC_DEFAULT,20)); 124c4762a1bSJed Brown 125c4762a1bSJed Brown /* 126c4762a1bSJed Brown Set SNES/KSP/KSP/PC runtime options, e.g., 127c4762a1bSJed Brown -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc> 128c4762a1bSJed Brown These options will override those specified above as long as 129c4762a1bSJed Brown SNESSetFromOptions() is called _after_ any other customization 130c4762a1bSJed Brown routines. 131c4762a1bSJed Brown */ 1325f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetFromOptions(snes)); 133c4762a1bSJed Brown 134c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 135c4762a1bSJed Brown Evaluate initial guess; then solve nonlinear system 136c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 137c4762a1bSJed Brown if (!flg) { 1385f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(x,pfive)); 139c4762a1bSJed Brown } else { 1405f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(x,&xx)); 141c4762a1bSJed Brown xx[0] = 2.0; xx[1] = 3.0; 1425f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(x,&xx)); 143c4762a1bSJed Brown } 144c4762a1bSJed Brown /* 145c4762a1bSJed Brown Note: The user should initialize the vector, x, with the initial guess 146c4762a1bSJed Brown for the nonlinear solver prior to calling SNESSolve(). In particular, 147c4762a1bSJed Brown to employ an initial guess of zero, the user should explicitly set 148c4762a1bSJed Brown this vector to zero by calling VecSet(). 149c4762a1bSJed Brown */ 150c4762a1bSJed Brown 1515f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSolve(snes,NULL,x)); 1525f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetIterationNumber(snes,&its)); 153c4762a1bSJed Brown if (flg) { 154c4762a1bSJed Brown Vec f; 1555f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); 1565f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetFunction(snes,&f,0,0)); 1575f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(r,PETSC_VIEWER_STDOUT_WORLD)); 158c4762a1bSJed Brown } 1595f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"number of SNES iterations = %D\n\n",its)); 160c4762a1bSJed Brown 161c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 162c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 163c4762a1bSJed Brown are no longer needed. 164c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 165c4762a1bSJed Brown 1665f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); CHKERRQ(VecDestroy(&r)); 1675f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&J)); CHKERRQ(SNESDestroy(&snes)); 168c4762a1bSJed Brown PetscFunctionReturn(0); 169c4762a1bSJed Brown } 170c4762a1bSJed Brown 171c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 172c4762a1bSJed Brown /* 173c4762a1bSJed Brown FormFunction1 - Evaluates nonlinear function, F(x). 174c4762a1bSJed Brown 175c4762a1bSJed Brown Input Parameters: 176c4762a1bSJed Brown . snes - the SNES context 177c4762a1bSJed Brown . x - input vector 178c4762a1bSJed Brown . dummy - optional user-defined context (not used here) 179c4762a1bSJed Brown 180c4762a1bSJed Brown Output Parameter: 181c4762a1bSJed Brown . f - function vector 182c4762a1bSJed Brown */ 183c4762a1bSJed Brown static PetscErrorCode FormFunction1(SNES snes,Vec x,Vec f,void *dummy) 184c4762a1bSJed Brown { 185c4762a1bSJed Brown const PetscScalar *xx; 186c4762a1bSJed Brown PetscScalar *ff; 187c4762a1bSJed Brown 188c4762a1bSJed Brown PetscFunctionBeginUser; 189c4762a1bSJed Brown /* 190c4762a1bSJed Brown Get pointers to vector data. 191c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 192c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 193c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 194c4762a1bSJed Brown the array. 195c4762a1bSJed Brown */ 1965f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(x,&xx)); 1975f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(f,&ff)); 198c4762a1bSJed Brown 199c4762a1bSJed Brown /* 200c4762a1bSJed Brown Compute function 201c4762a1bSJed Brown */ 202c4762a1bSJed Brown ff[0] = xx[0]*xx[0] + xx[0]*xx[1] - 3.0; 203c4762a1bSJed Brown ff[1] = xx[0]*xx[1] + xx[1]*xx[1] - 6.0; 204c4762a1bSJed Brown 205c4762a1bSJed Brown /* 206c4762a1bSJed Brown Restore vectors 207c4762a1bSJed Brown */ 2085f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(x,&xx)); 2095f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(f,&ff)); 210c4762a1bSJed Brown PetscFunctionReturn(0); 211c4762a1bSJed Brown } 212c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 213c4762a1bSJed Brown /* 214c4762a1bSJed Brown FormJacobian1 - Evaluates Jacobian matrix. 215c4762a1bSJed Brown 216c4762a1bSJed Brown Input Parameters: 217c4762a1bSJed Brown . snes - the SNES context 218c4762a1bSJed Brown . x - input vector 219c4762a1bSJed Brown . dummy - optional user-defined context (not used here) 220c4762a1bSJed Brown 221c4762a1bSJed Brown Output Parameters: 222c4762a1bSJed Brown . jac - Jacobian matrix 223c4762a1bSJed Brown . B - optionally different preconditioning matrix 224c4762a1bSJed Brown . flag - flag indicating matrix structure 225c4762a1bSJed Brown */ 226c4762a1bSJed Brown static PetscErrorCode FormJacobian1(SNES snes,Vec x,Mat jac,Mat B,void *dummy) 227c4762a1bSJed Brown { 228c4762a1bSJed Brown const PetscScalar *xx; 229c4762a1bSJed Brown PetscScalar A[4]; 230c4762a1bSJed Brown PetscInt idx[2] = {0,1}; 231c4762a1bSJed Brown 232c4762a1bSJed Brown PetscFunctionBeginUser; 233c4762a1bSJed Brown /* 234c4762a1bSJed Brown Get pointer to vector data 235c4762a1bSJed Brown */ 2365f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(x,&xx)); 237c4762a1bSJed Brown 238c4762a1bSJed Brown /* 239c4762a1bSJed Brown Compute Jacobian entries and insert into matrix. 240c4762a1bSJed Brown - Since this is such a small problem, we set all entries for 241c4762a1bSJed Brown the matrix at once. 242c4762a1bSJed Brown */ 243c4762a1bSJed Brown A[0] = 2.0*xx[0] + xx[1]; A[1] = xx[0]; 244c4762a1bSJed Brown A[2] = xx[1]; A[3] = xx[0] + 2.0*xx[1]; 2455f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(jac,2,idx,2,idx,A,INSERT_VALUES)); 246c4762a1bSJed Brown 247c4762a1bSJed Brown /* 248c4762a1bSJed Brown Restore vector 249c4762a1bSJed Brown */ 2505f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(x,&xx)); 251c4762a1bSJed Brown 252c4762a1bSJed Brown /* 253c4762a1bSJed Brown Assemble matrix 254c4762a1bSJed Brown */ 2555f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY)); 2565f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY)); 257c4762a1bSJed Brown PetscFunctionReturn(0); 258c4762a1bSJed Brown } 259c4762a1bSJed Brown 260c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 261c4762a1bSJed Brown static PetscErrorCode FormFunction2(SNES snes,Vec x,Vec f,void *dummy) 262c4762a1bSJed Brown { 263c4762a1bSJed Brown const PetscScalar *xx; 264c4762a1bSJed Brown PetscScalar *ff; 265c4762a1bSJed Brown 266c4762a1bSJed Brown PetscFunctionBeginUser; 267c4762a1bSJed Brown /* 268c4762a1bSJed Brown Get pointers to vector data. 269c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 270c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 271c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 272c4762a1bSJed Brown the array. 273c4762a1bSJed Brown */ 2745f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(x,&xx)); 2755f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(f,&ff)); 276c4762a1bSJed Brown 277c4762a1bSJed Brown /* 278c4762a1bSJed Brown Compute function 279c4762a1bSJed Brown */ 280c4762a1bSJed Brown ff[0] = PetscSinScalar(3.0*xx[0]) + xx[0]; 281c4762a1bSJed Brown ff[1] = xx[1]; 282c4762a1bSJed Brown 283c4762a1bSJed Brown /* 284c4762a1bSJed Brown Restore vectors 285c4762a1bSJed Brown */ 2865f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(x,&xx)); 2875f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(f,&ff)); 288c4762a1bSJed Brown PetscFunctionReturn(0); 289c4762a1bSJed Brown } 290c4762a1bSJed Brown 291c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 292c4762a1bSJed Brown static PetscErrorCode FormJacobian2(SNES snes,Vec x,Mat jac,Mat B,void *dummy) 293c4762a1bSJed Brown { 294c4762a1bSJed Brown const PetscScalar *xx; 295c4762a1bSJed Brown PetscScalar A[4]; 296c4762a1bSJed Brown PetscInt idx[2] = {0,1}; 297c4762a1bSJed Brown 298c4762a1bSJed Brown PetscFunctionBeginUser; 299c4762a1bSJed Brown /* 300c4762a1bSJed Brown Get pointer to vector data 301c4762a1bSJed Brown */ 3025f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(x,&xx)); 303c4762a1bSJed Brown 304c4762a1bSJed Brown /* 305c4762a1bSJed Brown Compute Jacobian entries and insert into matrix. 306c4762a1bSJed Brown - Since this is such a small problem, we set all entries for 307c4762a1bSJed Brown the matrix at once. 308c4762a1bSJed Brown */ 309c4762a1bSJed Brown A[0] = 3.0*PetscCosScalar(3.0*xx[0]) + 1.0; A[1] = 0.0; 310c4762a1bSJed Brown A[2] = 0.0; A[3] = 1.0; 3115f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(jac,2,idx,2,idx,A,INSERT_VALUES)); 312c4762a1bSJed Brown 313c4762a1bSJed Brown /* 314c4762a1bSJed Brown Restore vector 315c4762a1bSJed Brown */ 3165f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(x,&xx)); 317c4762a1bSJed Brown 318c4762a1bSJed Brown /* 319c4762a1bSJed Brown Assemble matrix 320c4762a1bSJed Brown */ 3215f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY)); 3225f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY)); 323c4762a1bSJed Brown PetscFunctionReturn(0); 324c4762a1bSJed Brown } 325c4762a1bSJed Brown 326b5675b0fSBarry Smith static PetscErrorCode block_system(void) 327c4762a1bSJed Brown { 328c4762a1bSJed Brown SNES snes; /* nonlinear solver context */ 329c4762a1bSJed Brown KSP ksp; /* linear solver context */ 330c4762a1bSJed Brown PC pc; /* preconditioner context */ 331c4762a1bSJed Brown Vec x,r; /* solution, residual vectors */ 332c4762a1bSJed Brown Mat J; /* Jacobian matrix */ 333c4762a1bSJed Brown PetscInt its; 334c4762a1bSJed Brown PetscScalar pfive = .5; 335c4762a1bSJed Brown PetscBool flg; 336c4762a1bSJed Brown 337c4762a1bSJed Brown Mat j11, j12, j21, j22; 338c4762a1bSJed Brown Vec x1, x2, r1, r2; 339c4762a1bSJed Brown Vec bv; 340c4762a1bSJed Brown Vec bx[2]; 341c4762a1bSJed Brown Mat bA[2][2]; 342c4762a1bSJed Brown 343c4762a1bSJed Brown PetscFunctionBeginUser; 3445f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "\n\n========================= Block system =========================\n\n")); 345c4762a1bSJed Brown 3465f80ce2aSJacob Faibussowitsch CHKERRQ(SNESCreate(PETSC_COMM_WORLD,&snes)); 347c4762a1bSJed Brown 348c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 349c4762a1bSJed Brown Create matrix and vector data structures; set corresponding routines 350c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 351c4762a1bSJed Brown 352c4762a1bSJed Brown /* 353c4762a1bSJed Brown Create sub vectors for solution and nonlinear function 354c4762a1bSJed Brown */ 3555f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,1,&x1)); 3565f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x1,&r1)); 357c4762a1bSJed Brown 3585f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,1,&x2)); 3595f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x2,&r2)); 360c4762a1bSJed Brown 361c4762a1bSJed Brown /* 362c4762a1bSJed Brown Create the block vectors 363c4762a1bSJed Brown */ 364c4762a1bSJed Brown bx[0] = x1; 365c4762a1bSJed Brown bx[1] = x2; 3665f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateNest(PETSC_COMM_WORLD,2,NULL,bx,&x)); 3675f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyBegin(x)); 3685f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyEnd(x)); 3695f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x1)); 3705f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x2)); 371c4762a1bSJed Brown 372c4762a1bSJed Brown bx[0] = r1; 373c4762a1bSJed Brown bx[1] = r2; 3745f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateNest(PETSC_COMM_WORLD,2,NULL,bx,&r)); 3755f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&r1)); 3765f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&r2)); 3775f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyBegin(r)); 3785f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyEnd(r)); 379c4762a1bSJed Brown 380c4762a1bSJed Brown /* 381c4762a1bSJed Brown Create sub Jacobian matrix data structure 382c4762a1bSJed Brown */ 3835f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD, &j11)); 3845f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(j11, 1, 1, 1, 1)); 3855f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(j11, MATSEQAIJ)); 3865f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(j11)); 387c4762a1bSJed Brown 3885f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD, &j12)); 3895f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(j12, 1, 1, 1, 1)); 3905f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(j12, MATSEQAIJ)); 3915f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(j12)); 392c4762a1bSJed Brown 3935f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD, &j21)); 3945f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(j21, 1, 1, 1, 1)); 3955f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(j21, MATSEQAIJ)); 3965f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(j21)); 397c4762a1bSJed Brown 3985f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD, &j22)); 3995f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(j22, PETSC_DECIDE, PETSC_DECIDE, 1, 1)); 4005f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(j22, MATSEQAIJ)); 4015f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(j22)); 402c4762a1bSJed Brown /* 403c4762a1bSJed Brown Create block Jacobian matrix data structure 404c4762a1bSJed Brown */ 405c4762a1bSJed Brown bA[0][0] = j11; 406c4762a1bSJed Brown bA[0][1] = j12; 407c4762a1bSJed Brown bA[1][0] = j21; 408c4762a1bSJed Brown bA[1][1] = j22; 409c4762a1bSJed Brown 4105f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateNest(PETSC_COMM_WORLD,2,NULL,2,NULL,&bA[0][0],&J)); 4115f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(J)); 4125f80ce2aSJacob Faibussowitsch CHKERRQ(MatNestSetVecType(J,VECNEST)); 4135f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&j11)); 4145f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&j12)); 4155f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&j21)); 4165f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&j22)); 417c4762a1bSJed Brown 4185f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY)); 4195f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY)); 420c4762a1bSJed Brown 4215f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-hard",&flg)); 422c4762a1bSJed Brown if (!flg) { 423c4762a1bSJed Brown /* 424c4762a1bSJed Brown Set function evaluation routine and vector. 425c4762a1bSJed Brown */ 4265f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetFunction(snes,r,FormFunction1_block,NULL)); 427c4762a1bSJed Brown 428c4762a1bSJed Brown /* 429c4762a1bSJed Brown Set Jacobian matrix data structure and Jacobian evaluation routine 430c4762a1bSJed Brown */ 4315f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetJacobian(snes,J,J,FormJacobian1_block,NULL)); 432c4762a1bSJed Brown } else { 4335f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetFunction(snes,r,FormFunction2_block,NULL)); 4345f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetJacobian(snes,J,J,FormJacobian2_block,NULL)); 435c4762a1bSJed Brown } 436c4762a1bSJed Brown 437c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 438c4762a1bSJed Brown Customize nonlinear solver; set runtime options 439c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 440c4762a1bSJed Brown 441c4762a1bSJed Brown /* 442c4762a1bSJed Brown Set linear solver defaults for this problem. By extracting the 443c4762a1bSJed Brown KSP, KSP, and PC contexts from the SNES context, we can then 444c4762a1bSJed Brown directly call any KSP, KSP, and PC routines to set various options. 445c4762a1bSJed Brown */ 4465f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetKSP(snes,&ksp)); 4475f80ce2aSJacob Faibussowitsch CHKERRQ(KSPGetPC(ksp,&pc)); 4485f80ce2aSJacob Faibussowitsch CHKERRQ(PCSetType(pc,PCNONE)); 4495f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetTolerances(ksp,1.e-4,PETSC_DEFAULT,PETSC_DEFAULT,20)); 450c4762a1bSJed Brown 451c4762a1bSJed Brown /* 452c4762a1bSJed Brown Set SNES/KSP/KSP/PC runtime options, e.g., 453c4762a1bSJed Brown -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc> 454c4762a1bSJed Brown These options will override those specified above as long as 455c4762a1bSJed Brown SNESSetFromOptions() is called _after_ any other customization 456c4762a1bSJed Brown routines. 457c4762a1bSJed Brown */ 4585f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetFromOptions(snes)); 459c4762a1bSJed Brown 460c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 461c4762a1bSJed Brown Evaluate initial guess; then solve nonlinear system 462c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 463c4762a1bSJed Brown if (!flg) { 4645f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(x,pfive)); 465c4762a1bSJed Brown } else { 466c4762a1bSJed Brown Vec *vecs; 4675f80ce2aSJacob Faibussowitsch CHKERRQ(VecNestGetSubVecs(x, NULL, &vecs)); 468c4762a1bSJed Brown bv = vecs[0]; 4695f80ce2aSJacob Faibussowitsch /* CHKERRQ(VecBlockGetSubVec(x, 0, &bv)); */ 4705f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetValue(bv, 0, 2.0, INSERT_VALUES)); /* xx[0] = 2.0; */ 4715f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyBegin(bv)); 4725f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyEnd(bv)); 473c4762a1bSJed Brown 4745f80ce2aSJacob Faibussowitsch /* CHKERRQ(VecBlockGetSubVec(x, 1, &bv)); */ 475c4762a1bSJed Brown bv = vecs[1]; 4765f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetValue(bv, 0, 3.0, INSERT_VALUES)); /* xx[1] = 3.0; */ 4775f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyBegin(bv)); 4785f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyEnd(bv)); 479c4762a1bSJed Brown } 480c4762a1bSJed Brown /* 481c4762a1bSJed Brown Note: The user should initialize the vector, x, with the initial guess 482c4762a1bSJed Brown for the nonlinear solver prior to calling SNESSolve(). In particular, 483c4762a1bSJed Brown to employ an initial guess of zero, the user should explicitly set 484c4762a1bSJed Brown this vector to zero by calling VecSet(). 485c4762a1bSJed Brown */ 4865f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSolve(snes,NULL,x)); 4875f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetIterationNumber(snes,&its)); 488c4762a1bSJed Brown if (flg) { 489c4762a1bSJed Brown Vec f; 4905f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); 4915f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetFunction(snes,&f,0,0)); 4925f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(r,PETSC_VIEWER_STDOUT_WORLD)); 493c4762a1bSJed Brown } 494c4762a1bSJed Brown 4955f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"number of SNES iterations = %D\n\n",its)); 496c4762a1bSJed Brown 497c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 498c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 499c4762a1bSJed Brown are no longer needed. 500c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 5015f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); CHKERRQ(VecDestroy(&r)); 5025f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&J)); CHKERRQ(SNESDestroy(&snes)); 503c4762a1bSJed Brown PetscFunctionReturn(0); 504c4762a1bSJed Brown } 505c4762a1bSJed Brown 506c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 507c4762a1bSJed Brown static PetscErrorCode FormFunction1_block(SNES snes,Vec x,Vec f,void *dummy) 508c4762a1bSJed Brown { 509c4762a1bSJed Brown Vec *xx, *ff, x1,x2, f1,f2; 510c4762a1bSJed Brown PetscScalar ff_0, ff_1; 511c4762a1bSJed Brown PetscScalar xx_0, xx_1; 512c4762a1bSJed Brown PetscInt index,nb; 513c4762a1bSJed Brown 514c4762a1bSJed Brown PetscFunctionBeginUser; 515c4762a1bSJed Brown /* get blocks for function */ 5165f80ce2aSJacob Faibussowitsch CHKERRQ(VecNestGetSubVecs(f, &nb, &ff)); 517c4762a1bSJed Brown f1 = ff[0]; f2 = ff[1]; 518c4762a1bSJed Brown 519c4762a1bSJed Brown /* get blocks for solution */ 5205f80ce2aSJacob Faibussowitsch CHKERRQ(VecNestGetSubVecs(x, &nb, &xx)); 521c4762a1bSJed Brown x1 = xx[0]; x2 = xx[1]; 522c4762a1bSJed Brown 523c4762a1bSJed Brown /* get solution values */ 524c4762a1bSJed Brown index = 0; 5255f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetValues(x1,1, &index, &xx_0)); 5265f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetValues(x2,1, &index, &xx_1)); 527c4762a1bSJed Brown 528c4762a1bSJed Brown /* Compute function */ 529c4762a1bSJed Brown ff_0 = xx_0*xx_0 + xx_0*xx_1 - 3.0; 530c4762a1bSJed Brown ff_1 = xx_0*xx_1 + xx_1*xx_1 - 6.0; 531c4762a1bSJed Brown 532c4762a1bSJed Brown /* set function values */ 5335f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetValue(f1, index, ff_0, INSERT_VALUES)); 534c4762a1bSJed Brown 5355f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetValue(f2, index, ff_1, INSERT_VALUES)); 536c4762a1bSJed Brown 5375f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyBegin(f)); 5385f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyEnd(f)); 539c4762a1bSJed Brown PetscFunctionReturn(0); 540c4762a1bSJed Brown } 541c4762a1bSJed Brown 542c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 543c4762a1bSJed Brown static PetscErrorCode FormJacobian1_block(SNES snes,Vec x,Mat jac,Mat B,void *dummy) 544c4762a1bSJed Brown { 545c4762a1bSJed Brown Vec *xx, x1,x2; 546c4762a1bSJed Brown PetscScalar xx_0, xx_1; 547c4762a1bSJed Brown PetscInt index,nb; 548c4762a1bSJed Brown PetscScalar A_00, A_01, A_10, A_11; 549c4762a1bSJed Brown Mat j11, j12, j21, j22; 550c4762a1bSJed Brown Mat **mats; 551c4762a1bSJed Brown 552c4762a1bSJed Brown PetscFunctionBeginUser; 553c4762a1bSJed Brown /* get blocks for solution */ 5545f80ce2aSJacob Faibussowitsch CHKERRQ(VecNestGetSubVecs(x, &nb, &xx)); 555c4762a1bSJed Brown x1 = xx[0]; x2 = xx[1]; 556c4762a1bSJed Brown 557c4762a1bSJed Brown /* get solution values */ 558c4762a1bSJed Brown index = 0; 5595f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetValues(x1,1, &index, &xx_0)); 5605f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetValues(x2,1, &index, &xx_1)); 561c4762a1bSJed Brown 562c4762a1bSJed Brown /* get block matrices */ 5635f80ce2aSJacob Faibussowitsch CHKERRQ(MatNestGetSubMats(jac,NULL,NULL,&mats)); 564c4762a1bSJed Brown j11 = mats[0][0]; 565c4762a1bSJed Brown j12 = mats[0][1]; 566c4762a1bSJed Brown j21 = mats[1][0]; 567c4762a1bSJed Brown j22 = mats[1][1]; 568c4762a1bSJed Brown 569c4762a1bSJed Brown /* compute jacobian entries */ 570c4762a1bSJed Brown A_00 = 2.0*xx_0 + xx_1; 571c4762a1bSJed Brown A_01 = xx_0; 572c4762a1bSJed Brown A_10 = xx_1; 573c4762a1bSJed Brown A_11 = xx_0 + 2.0*xx_1; 574c4762a1bSJed Brown 575c4762a1bSJed Brown /* set jacobian values */ 5765f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValue(j11, 0,0, A_00, INSERT_VALUES)); 5775f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValue(j12, 0,0, A_01, INSERT_VALUES)); 5785f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValue(j21, 0,0, A_10, INSERT_VALUES)); 5795f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValue(j22, 0,0, A_11, INSERT_VALUES)); 580c4762a1bSJed Brown 581c4762a1bSJed Brown /* Assemble sub matrix */ 5825f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY)); 5835f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY)); 584c4762a1bSJed Brown PetscFunctionReturn(0); 585c4762a1bSJed Brown } 586c4762a1bSJed Brown 587c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 588c4762a1bSJed Brown static PetscErrorCode FormFunction2_block(SNES snes,Vec x,Vec f,void *dummy) 589c4762a1bSJed Brown { 590c4762a1bSJed Brown PetscScalar *ff; 591c4762a1bSJed Brown const PetscScalar *xx; 592c4762a1bSJed Brown 593c4762a1bSJed Brown PetscFunctionBeginUser; 594c4762a1bSJed Brown /* 595c4762a1bSJed Brown Get pointers to vector data. 596c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 597c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 598c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 599c4762a1bSJed Brown the array. 600c4762a1bSJed Brown */ 6015f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(x,&xx)); 6025f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(f,&ff)); 603c4762a1bSJed Brown 604c4762a1bSJed Brown /* 605c4762a1bSJed Brown Compute function 606c4762a1bSJed Brown */ 607c4762a1bSJed Brown ff[0] = PetscSinScalar(3.0*xx[0]) + xx[0]; 608c4762a1bSJed Brown ff[1] = xx[1]; 609c4762a1bSJed Brown 610c4762a1bSJed Brown /* 611c4762a1bSJed Brown Restore vectors 612c4762a1bSJed Brown */ 6135f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(x,&xx)); 6145f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(f,&ff)); 615c4762a1bSJed Brown PetscFunctionReturn(0); 616c4762a1bSJed Brown } 617c4762a1bSJed Brown 618c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 619c4762a1bSJed Brown static PetscErrorCode FormJacobian2_block(SNES snes,Vec x,Mat jac,Mat B,void *dummy) 620c4762a1bSJed Brown { 621c4762a1bSJed Brown const PetscScalar *xx; 622c4762a1bSJed Brown PetscScalar A[4]; 623c4762a1bSJed Brown PetscInt idx[2] = {0,1}; 624c4762a1bSJed Brown 625c4762a1bSJed Brown PetscFunctionBeginUser; 626c4762a1bSJed Brown /* 627c4762a1bSJed Brown Get pointer to vector data 628c4762a1bSJed Brown */ 6295f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(x,&xx)); 630c4762a1bSJed Brown 631c4762a1bSJed Brown /* 632c4762a1bSJed Brown Compute Jacobian entries and insert into matrix. 633c4762a1bSJed Brown - Since this is such a small problem, we set all entries for 634c4762a1bSJed Brown the matrix at once. 635c4762a1bSJed Brown */ 636c4762a1bSJed Brown A[0] = 3.0*PetscCosScalar(3.0*xx[0]) + 1.0; A[1] = 0.0; 637c4762a1bSJed Brown A[2] = 0.0; A[3] = 1.0; 6385f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(jac,2,idx,2,idx,A,INSERT_VALUES)); 639c4762a1bSJed Brown 640c4762a1bSJed Brown /* 641c4762a1bSJed Brown Restore vector 642c4762a1bSJed Brown */ 6435f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(x,&xx)); 644c4762a1bSJed Brown 645c4762a1bSJed Brown /* 646c4762a1bSJed Brown Assemble matrix 647c4762a1bSJed Brown */ 6485f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY)); 6495f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY)); 650c4762a1bSJed Brown PetscFunctionReturn(0); 651c4762a1bSJed Brown } 652c4762a1bSJed Brown 653c4762a1bSJed Brown int main(int argc,char **argv) 654c4762a1bSJed Brown { 655c4762a1bSJed Brown PetscMPIInt size; 656c4762a1bSJed Brown 657*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help)); 6585f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 6592c71b3e2SJacob Faibussowitsch PetscCheckFalse(size != 1,PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 6605f80ce2aSJacob Faibussowitsch CHKERRQ(assembled_system()); 6615f80ce2aSJacob Faibussowitsch CHKERRQ(block_system()); 662*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 663*b122ec5aSJacob Faibussowitsch return 0; 664c4762a1bSJed Brown } 665c4762a1bSJed Brown 666c4762a1bSJed Brown /*TEST 667c4762a1bSJed Brown 668c4762a1bSJed Brown test: 669c4762a1bSJed Brown args: -snes_monitor_short 670c4762a1bSJed Brown requires: !single 671c4762a1bSJed Brown 672c4762a1bSJed Brown TEST*/ 673