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 11c4762a1bSJed Brown /* 12c4762a1bSJed Brown Include "petscsnes.h" so that we can use SNES solvers. Note that this 13c4762a1bSJed Brown file automatically includes: 14c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 15c4762a1bSJed Brown petscsys.h - system routines petscmat.h - matrices 16c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods 17c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners 18c4762a1bSJed Brown petscksp.h - linear solvers 19c4762a1bSJed Brown */ 20c4762a1bSJed Brown #include <petscsnes.h> 21c4762a1bSJed Brown 22c4762a1bSJed Brown /* 23c4762a1bSJed Brown This example is block version of the test found at 24c4762a1bSJed Brown ${PETSC_DIR}/src/snes/tutorials/ex1.c 25c4762a1bSJed Brown In this test we replace the Jacobian systems 26c4762a1bSJed Brown [J]{x} = {F} 27c4762a1bSJed Brown where 28c4762a1bSJed Brown 29c4762a1bSJed Brown [J] = (j_00, j_01), {x} = (x_0, x_1)^T, {F} = (f_0, f_1)^T 30c4762a1bSJed Brown (j_10, j_11) 31c4762a1bSJed Brown where [J] \in \mathbb^{2 \times 2}, {x},{F} \in \mathbb^{2 \times 1}, 32c4762a1bSJed Brown 33c4762a1bSJed Brown with a block system in which each block is of length 1. 34c4762a1bSJed Brown i.e. The block system is thus 35c4762a1bSJed Brown 36c4762a1bSJed Brown [J] = ([j00], [j01]), {x} = ({x0}, {x1})^T, {F} = ({f0}, {f1})^T 37c4762a1bSJed Brown ([j10], [j11]) 38c4762a1bSJed Brown where 39c4762a1bSJed Brown [j00], [j01], [j10], [j11] \in \mathbb^{1 \times 1} 40c4762a1bSJed Brown {x0}, {x1}, {f0}, {f1} \in \mathbb^{1 \times 1} 41c4762a1bSJed Brown 42c4762a1bSJed Brown In practice we would not bother defing blocks of size one, and would instead assemble the 43c4762a1bSJed Brown full system. This is just a simple test to illustrate how to manipulate the blocks and 44c4762a1bSJed Brown to confirm the implementation is correct. 45c4762a1bSJed Brown */ 46c4762a1bSJed Brown 47c4762a1bSJed Brown /* 48c4762a1bSJed Brown User-defined routines 49c4762a1bSJed Brown */ 50c4762a1bSJed Brown static PetscErrorCode FormJacobian1(SNES,Vec,Mat,Mat,void*); 51c4762a1bSJed Brown static PetscErrorCode FormFunction1(SNES,Vec,Vec,void*); 52c4762a1bSJed Brown static PetscErrorCode FormJacobian2(SNES,Vec,Mat,Mat,void*); 53c4762a1bSJed Brown static PetscErrorCode FormFunction2(SNES,Vec,Vec,void*); 54c4762a1bSJed Brown static PetscErrorCode FormJacobian1_block(SNES,Vec,Mat,Mat,void*); 55c4762a1bSJed Brown static PetscErrorCode FormFunction1_block(SNES,Vec,Vec,void*); 56c4762a1bSJed Brown static PetscErrorCode FormJacobian2_block(SNES,Vec,Mat,Mat,void*); 57c4762a1bSJed Brown static PetscErrorCode FormFunction2_block(SNES,Vec,Vec,void*); 58c4762a1bSJed Brown 59c4762a1bSJed Brown 60c4762a1bSJed Brown static PetscErrorCode assembled_system(void) 61c4762a1bSJed Brown { 62c4762a1bSJed Brown SNES snes; /* nonlinear solver context */ 63c4762a1bSJed Brown KSP ksp; /* linear solver context */ 64c4762a1bSJed Brown PC pc; /* preconditioner context */ 65c4762a1bSJed Brown Vec x,r; /* solution, residual vectors */ 66c4762a1bSJed Brown Mat J; /* Jacobian matrix */ 67c4762a1bSJed Brown PetscErrorCode ierr; 68c4762a1bSJed Brown PetscInt its; 69c4762a1bSJed Brown PetscScalar pfive = .5,*xx; 70c4762a1bSJed Brown PetscBool flg; 71c4762a1bSJed Brown 72c4762a1bSJed Brown PetscFunctionBeginUser; 73*b5675b0fSBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD, "\n\n========================= Assembled system =========================\n\n");CHKERRQ(ierr); 74c4762a1bSJed Brown 75c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 76c4762a1bSJed Brown Create nonlinear solver context 77c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 78c4762a1bSJed Brown 79c4762a1bSJed Brown ierr = SNESCreate(PETSC_COMM_WORLD,&snes);CHKERRQ(ierr); 80c4762a1bSJed Brown 81c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 82c4762a1bSJed Brown Create matrix and vector data structures; set corresponding routines 83c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 84c4762a1bSJed Brown 85c4762a1bSJed Brown /* 86c4762a1bSJed Brown Create vectors for solution and nonlinear function 87c4762a1bSJed Brown */ 88c4762a1bSJed Brown ierr = VecCreateSeq(PETSC_COMM_SELF,2,&x);CHKERRQ(ierr); 89c4762a1bSJed Brown ierr = VecDuplicate(x,&r);CHKERRQ(ierr); 90c4762a1bSJed Brown 91c4762a1bSJed Brown /* 92c4762a1bSJed Brown Create Jacobian matrix data structure 93c4762a1bSJed Brown */ 94c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_SELF,&J);CHKERRQ(ierr); 95c4762a1bSJed Brown ierr = MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,2,2);CHKERRQ(ierr); 96c4762a1bSJed Brown ierr = MatSetFromOptions(J);CHKERRQ(ierr); 97c4762a1bSJed Brown ierr = MatSetUp(J);CHKERRQ(ierr); 98c4762a1bSJed Brown 99c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL,"-hard",&flg);CHKERRQ(ierr); 100c4762a1bSJed Brown if (!flg) { 101c4762a1bSJed Brown /* 102c4762a1bSJed Brown Set function evaluation routine and vector. 103c4762a1bSJed Brown */ 104c4762a1bSJed Brown ierr = SNESSetFunction(snes,r,FormFunction1,NULL);CHKERRQ(ierr); 105c4762a1bSJed Brown 106c4762a1bSJed Brown /* 107c4762a1bSJed Brown Set Jacobian matrix data structure and Jacobian evaluation routine 108c4762a1bSJed Brown */ 109c4762a1bSJed Brown ierr = SNESSetJacobian(snes,J,J,FormJacobian1,NULL);CHKERRQ(ierr); 110c4762a1bSJed Brown } else { 111c4762a1bSJed Brown ierr = SNESSetFunction(snes,r,FormFunction2,NULL);CHKERRQ(ierr); 112c4762a1bSJed Brown ierr = SNESSetJacobian(snes,J,J,FormJacobian2,NULL);CHKERRQ(ierr); 113c4762a1bSJed Brown } 114c4762a1bSJed Brown 115c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 116c4762a1bSJed Brown Customize nonlinear solver; set runtime options 117c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 118c4762a1bSJed Brown 119c4762a1bSJed Brown /* 120c4762a1bSJed Brown Set linear solver defaults for this problem. By extracting the 121c4762a1bSJed Brown KSP, KSP, and PC contexts from the SNES context, we can then 122c4762a1bSJed Brown directly call any KSP, KSP, and PC routines to set various options. 123c4762a1bSJed Brown */ 124c4762a1bSJed Brown ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 125c4762a1bSJed Brown ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 126c4762a1bSJed Brown ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 127c4762a1bSJed Brown ierr = KSPSetTolerances(ksp,1.e-4,PETSC_DEFAULT,PETSC_DEFAULT,20);CHKERRQ(ierr); 128c4762a1bSJed Brown 129c4762a1bSJed Brown /* 130c4762a1bSJed Brown Set SNES/KSP/KSP/PC runtime options, e.g., 131c4762a1bSJed Brown -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc> 132c4762a1bSJed Brown These options will override those specified above as long as 133c4762a1bSJed Brown SNESSetFromOptions() is called _after_ any other customization 134c4762a1bSJed Brown routines. 135c4762a1bSJed Brown */ 136c4762a1bSJed Brown ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 137c4762a1bSJed Brown 138c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 139c4762a1bSJed Brown Evaluate initial guess; then solve nonlinear system 140c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 141c4762a1bSJed Brown if (!flg) { 142c4762a1bSJed Brown ierr = VecSet(x,pfive);CHKERRQ(ierr); 143c4762a1bSJed Brown } else { 144c4762a1bSJed Brown ierr = VecGetArray(x,&xx);CHKERRQ(ierr); 145c4762a1bSJed Brown xx[0] = 2.0; xx[1] = 3.0; 146c4762a1bSJed Brown ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr); 147c4762a1bSJed Brown } 148c4762a1bSJed Brown /* 149c4762a1bSJed Brown Note: The user should initialize the vector, x, with the initial guess 150c4762a1bSJed Brown for the nonlinear solver prior to calling SNESSolve(). In particular, 151c4762a1bSJed Brown to employ an initial guess of zero, the user should explicitly set 152c4762a1bSJed Brown this vector to zero by calling VecSet(). 153c4762a1bSJed Brown */ 154c4762a1bSJed Brown 155c4762a1bSJed Brown ierr = SNESSolve(snes,NULL,x);CHKERRQ(ierr); 156c4762a1bSJed Brown ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr); 157c4762a1bSJed Brown if (flg) { 158c4762a1bSJed Brown Vec f; 159c4762a1bSJed Brown ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 160c4762a1bSJed Brown ierr = SNESGetFunction(snes,&f,0,0);CHKERRQ(ierr); 161c4762a1bSJed Brown ierr = VecView(r,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 162c4762a1bSJed Brown } 163c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"number of SNES iterations = %D\n\n",its);CHKERRQ(ierr); 164c4762a1bSJed Brown 165c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 166c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 167c4762a1bSJed Brown are no longer needed. 168c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 169c4762a1bSJed Brown 170c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); ierr = VecDestroy(&r);CHKERRQ(ierr); 171c4762a1bSJed Brown ierr = MatDestroy(&J);CHKERRQ(ierr); ierr = SNESDestroy(&snes);CHKERRQ(ierr); 172c4762a1bSJed Brown PetscFunctionReturn(0); 173c4762a1bSJed Brown } 174c4762a1bSJed Brown 175c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 176c4762a1bSJed Brown /* 177c4762a1bSJed Brown FormFunction1 - Evaluates nonlinear function, F(x). 178c4762a1bSJed Brown 179c4762a1bSJed Brown Input Parameters: 180c4762a1bSJed Brown . snes - the SNES context 181c4762a1bSJed Brown . x - input vector 182c4762a1bSJed Brown . dummy - optional user-defined context (not used here) 183c4762a1bSJed Brown 184c4762a1bSJed Brown Output Parameter: 185c4762a1bSJed Brown . f - function vector 186c4762a1bSJed Brown */ 187c4762a1bSJed Brown static PetscErrorCode FormFunction1(SNES snes,Vec x,Vec f,void *dummy) 188c4762a1bSJed Brown { 189c4762a1bSJed Brown PetscErrorCode ierr; 190c4762a1bSJed Brown const PetscScalar *xx; 191c4762a1bSJed Brown PetscScalar *ff; 192c4762a1bSJed Brown 193c4762a1bSJed Brown PetscFunctionBeginUser; 194c4762a1bSJed Brown /* 195c4762a1bSJed Brown Get pointers to vector data. 196c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 197c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 198c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 199c4762a1bSJed Brown the array. 200c4762a1bSJed Brown */ 201c4762a1bSJed Brown ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr); 202c4762a1bSJed Brown ierr = VecGetArray(f,&ff);CHKERRQ(ierr); 203c4762a1bSJed Brown 204c4762a1bSJed Brown /* 205c4762a1bSJed Brown Compute function 206c4762a1bSJed Brown */ 207c4762a1bSJed Brown ff[0] = xx[0]*xx[0] + xx[0]*xx[1] - 3.0; 208c4762a1bSJed Brown ff[1] = xx[0]*xx[1] + xx[1]*xx[1] - 6.0; 209c4762a1bSJed Brown 210c4762a1bSJed Brown 211c4762a1bSJed Brown /* 212c4762a1bSJed Brown Restore vectors 213c4762a1bSJed Brown */ 214c4762a1bSJed Brown ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr); 215c4762a1bSJed Brown ierr = VecRestoreArray(f,&ff);CHKERRQ(ierr); 216c4762a1bSJed Brown PetscFunctionReturn(0); 217c4762a1bSJed Brown } 218c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 219c4762a1bSJed Brown /* 220c4762a1bSJed Brown FormJacobian1 - Evaluates Jacobian matrix. 221c4762a1bSJed Brown 222c4762a1bSJed Brown Input Parameters: 223c4762a1bSJed Brown . snes - the SNES context 224c4762a1bSJed Brown . x - input vector 225c4762a1bSJed Brown . dummy - optional user-defined context (not used here) 226c4762a1bSJed Brown 227c4762a1bSJed Brown Output Parameters: 228c4762a1bSJed Brown . jac - Jacobian matrix 229c4762a1bSJed Brown . B - optionally different preconditioning matrix 230c4762a1bSJed Brown . flag - flag indicating matrix structure 231c4762a1bSJed Brown */ 232c4762a1bSJed Brown static PetscErrorCode FormJacobian1(SNES snes,Vec x,Mat jac,Mat B,void *dummy) 233c4762a1bSJed Brown { 234c4762a1bSJed Brown const PetscScalar *xx; 235c4762a1bSJed Brown PetscScalar A[4]; 236c4762a1bSJed Brown PetscErrorCode ierr; 237c4762a1bSJed Brown PetscInt idx[2] = {0,1}; 238c4762a1bSJed Brown 239c4762a1bSJed Brown PetscFunctionBeginUser; 240c4762a1bSJed Brown /* 241c4762a1bSJed Brown Get pointer to vector data 242c4762a1bSJed Brown */ 243c4762a1bSJed Brown ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr); 244c4762a1bSJed Brown 245c4762a1bSJed Brown /* 246c4762a1bSJed Brown Compute Jacobian entries and insert into matrix. 247c4762a1bSJed Brown - Since this is such a small problem, we set all entries for 248c4762a1bSJed Brown the matrix at once. 249c4762a1bSJed Brown */ 250c4762a1bSJed Brown A[0] = 2.0*xx[0] + xx[1]; A[1] = xx[0]; 251c4762a1bSJed Brown A[2] = xx[1]; A[3] = xx[0] + 2.0*xx[1]; 252c4762a1bSJed Brown ierr = MatSetValues(jac,2,idx,2,idx,A,INSERT_VALUES);CHKERRQ(ierr); 253c4762a1bSJed Brown 254c4762a1bSJed Brown /* 255c4762a1bSJed Brown Restore vector 256c4762a1bSJed Brown */ 257c4762a1bSJed Brown ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr); 258c4762a1bSJed Brown 259c4762a1bSJed Brown /* 260c4762a1bSJed Brown Assemble matrix 261c4762a1bSJed Brown */ 262c4762a1bSJed Brown ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 263c4762a1bSJed Brown ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 264c4762a1bSJed Brown PetscFunctionReturn(0); 265c4762a1bSJed Brown } 266c4762a1bSJed Brown 267c4762a1bSJed Brown 268c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 269c4762a1bSJed Brown static PetscErrorCode FormFunction2(SNES snes,Vec x,Vec f,void *dummy) 270c4762a1bSJed Brown { 271c4762a1bSJed Brown PetscErrorCode ierr; 272c4762a1bSJed Brown const PetscScalar *xx; 273c4762a1bSJed Brown PetscScalar *ff; 274c4762a1bSJed Brown 275c4762a1bSJed Brown PetscFunctionBeginUser; 276c4762a1bSJed Brown /* 277c4762a1bSJed Brown Get pointers to vector data. 278c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 279c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 280c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 281c4762a1bSJed Brown the array. 282c4762a1bSJed Brown */ 283c4762a1bSJed Brown ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr); 284c4762a1bSJed Brown ierr = VecGetArray(f,&ff);CHKERRQ(ierr); 285c4762a1bSJed Brown 286c4762a1bSJed Brown /* 287c4762a1bSJed Brown Compute function 288c4762a1bSJed Brown */ 289c4762a1bSJed Brown ff[0] = PetscSinScalar(3.0*xx[0]) + xx[0]; 290c4762a1bSJed Brown ff[1] = xx[1]; 291c4762a1bSJed Brown 292c4762a1bSJed Brown /* 293c4762a1bSJed Brown Restore vectors 294c4762a1bSJed Brown */ 295c4762a1bSJed Brown ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr); 296c4762a1bSJed Brown ierr = VecRestoreArray(f,&ff);CHKERRQ(ierr); 297c4762a1bSJed Brown PetscFunctionReturn(0); 298c4762a1bSJed Brown } 299c4762a1bSJed Brown 300c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 301c4762a1bSJed Brown static PetscErrorCode FormJacobian2(SNES snes,Vec x,Mat jac,Mat B,void *dummy) 302c4762a1bSJed Brown { 303c4762a1bSJed Brown const PetscScalar *xx; 304c4762a1bSJed Brown PetscScalar A[4]; 305c4762a1bSJed Brown PetscErrorCode ierr; 306c4762a1bSJed Brown PetscInt idx[2] = {0,1}; 307c4762a1bSJed Brown 308c4762a1bSJed Brown PetscFunctionBeginUser; 309c4762a1bSJed Brown /* 310c4762a1bSJed Brown Get pointer to vector data 311c4762a1bSJed Brown */ 312c4762a1bSJed Brown ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr); 313c4762a1bSJed Brown 314c4762a1bSJed Brown /* 315c4762a1bSJed Brown Compute Jacobian entries and insert into matrix. 316c4762a1bSJed Brown - Since this is such a small problem, we set all entries for 317c4762a1bSJed Brown the matrix at once. 318c4762a1bSJed Brown */ 319c4762a1bSJed Brown A[0] = 3.0*PetscCosScalar(3.0*xx[0]) + 1.0; A[1] = 0.0; 320c4762a1bSJed Brown A[2] = 0.0; A[3] = 1.0; 321c4762a1bSJed Brown ierr = MatSetValues(jac,2,idx,2,idx,A,INSERT_VALUES);CHKERRQ(ierr); 322c4762a1bSJed Brown 323c4762a1bSJed Brown /* 324c4762a1bSJed Brown Restore vector 325c4762a1bSJed Brown */ 326c4762a1bSJed Brown ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr); 327c4762a1bSJed Brown 328c4762a1bSJed Brown /* 329c4762a1bSJed Brown Assemble matrix 330c4762a1bSJed Brown */ 331c4762a1bSJed Brown ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 332c4762a1bSJed Brown ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 333c4762a1bSJed Brown PetscFunctionReturn(0); 334c4762a1bSJed Brown } 335c4762a1bSJed Brown 336*b5675b0fSBarry Smith static PetscErrorCode block_system(void) 337c4762a1bSJed Brown { 338c4762a1bSJed Brown SNES snes; /* nonlinear solver context */ 339c4762a1bSJed Brown KSP ksp; /* linear solver context */ 340c4762a1bSJed Brown PC pc; /* preconditioner context */ 341c4762a1bSJed Brown Vec x,r; /* solution, residual vectors */ 342c4762a1bSJed Brown Mat J; /* Jacobian matrix */ 343c4762a1bSJed Brown PetscErrorCode ierr; 344c4762a1bSJed Brown PetscInt its; 345c4762a1bSJed Brown PetscScalar pfive = .5; 346c4762a1bSJed Brown PetscBool flg; 347c4762a1bSJed Brown 348c4762a1bSJed Brown Mat j11, j12, j21, j22; 349c4762a1bSJed Brown Vec x1, x2, r1, r2; 350c4762a1bSJed Brown Vec bv; 351c4762a1bSJed Brown Vec bx[2]; 352c4762a1bSJed Brown Mat bA[2][2]; 353c4762a1bSJed Brown 354c4762a1bSJed Brown PetscFunctionBeginUser; 355*b5675b0fSBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD, "\n\n========================= Block system =========================\n\n");CHKERRQ(ierr); 356c4762a1bSJed Brown 357c4762a1bSJed Brown ierr = SNESCreate(PETSC_COMM_WORLD,&snes);CHKERRQ(ierr); 358c4762a1bSJed Brown 359c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 360c4762a1bSJed Brown Create matrix and vector data structures; set corresponding routines 361c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 362c4762a1bSJed Brown 363c4762a1bSJed Brown /* 364c4762a1bSJed Brown Create sub vectors for solution and nonlinear function 365c4762a1bSJed Brown */ 366c4762a1bSJed Brown ierr = VecCreateSeq(PETSC_COMM_SELF,1,&x1);CHKERRQ(ierr); 367c4762a1bSJed Brown ierr = VecDuplicate(x1,&r1);CHKERRQ(ierr); 368c4762a1bSJed Brown 369c4762a1bSJed Brown ierr = VecCreateSeq(PETSC_COMM_SELF,1,&x2);CHKERRQ(ierr); 370c4762a1bSJed Brown ierr = VecDuplicate(x2,&r2);CHKERRQ(ierr); 371c4762a1bSJed Brown 372c4762a1bSJed Brown /* 373c4762a1bSJed Brown Create the block vectors 374c4762a1bSJed Brown */ 375c4762a1bSJed Brown bx[0] = x1; 376c4762a1bSJed Brown bx[1] = x2; 377c4762a1bSJed Brown ierr = VecCreateNest(PETSC_COMM_WORLD,2,NULL,bx,&x);CHKERRQ(ierr); 378c4762a1bSJed Brown ierr = VecAssemblyBegin(x);CHKERRQ(ierr); 379c4762a1bSJed Brown ierr = VecAssemblyEnd(x);CHKERRQ(ierr); 380c4762a1bSJed Brown ierr = VecDestroy(&x1);CHKERRQ(ierr); 381c4762a1bSJed Brown ierr = VecDestroy(&x2);CHKERRQ(ierr); 382c4762a1bSJed Brown 383c4762a1bSJed Brown bx[0] = r1; 384c4762a1bSJed Brown bx[1] = r2; 385c4762a1bSJed Brown ierr = VecCreateNest(PETSC_COMM_WORLD,2,NULL,bx,&r);CHKERRQ(ierr); 386c4762a1bSJed Brown ierr = VecDestroy(&r1);CHKERRQ(ierr); 387c4762a1bSJed Brown ierr = VecDestroy(&r2);CHKERRQ(ierr); 388c4762a1bSJed Brown ierr = VecAssemblyBegin(r);CHKERRQ(ierr); 389c4762a1bSJed Brown ierr = VecAssemblyEnd(r);CHKERRQ(ierr); 390c4762a1bSJed Brown 391c4762a1bSJed Brown /* 392c4762a1bSJed Brown Create sub Jacobian matrix data structure 393c4762a1bSJed Brown */ 394c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD, &j11);CHKERRQ(ierr); 395c4762a1bSJed Brown ierr = MatSetSizes(j11, 1, 1, 1, 1);CHKERRQ(ierr); 396c4762a1bSJed Brown ierr = MatSetType(j11, MATSEQAIJ);CHKERRQ(ierr); 397c4762a1bSJed Brown ierr = MatSetUp(j11);CHKERRQ(ierr); 398c4762a1bSJed Brown 399c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD, &j12);CHKERRQ(ierr); 400c4762a1bSJed Brown ierr = MatSetSizes(j12, 1, 1, 1, 1);CHKERRQ(ierr); 401c4762a1bSJed Brown ierr = MatSetType(j12, MATSEQAIJ);CHKERRQ(ierr); 402c4762a1bSJed Brown ierr = MatSetUp(j12);CHKERRQ(ierr); 403c4762a1bSJed Brown 404c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD, &j21);CHKERRQ(ierr); 405c4762a1bSJed Brown ierr = MatSetSizes(j21, 1, 1, 1, 1);CHKERRQ(ierr); 406c4762a1bSJed Brown ierr = MatSetType(j21, MATSEQAIJ);CHKERRQ(ierr); 407c4762a1bSJed Brown ierr = MatSetUp(j21);CHKERRQ(ierr); 408c4762a1bSJed Brown 409c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD, &j22);CHKERRQ(ierr); 410c4762a1bSJed Brown ierr = MatSetSizes(j22, PETSC_DECIDE, PETSC_DECIDE, 1, 1);CHKERRQ(ierr); 411c4762a1bSJed Brown ierr = MatSetType(j22, MATSEQAIJ);CHKERRQ(ierr); 412c4762a1bSJed Brown ierr = MatSetUp(j22);CHKERRQ(ierr); 413c4762a1bSJed Brown /* 414c4762a1bSJed Brown Create block Jacobian matrix data structure 415c4762a1bSJed Brown */ 416c4762a1bSJed Brown bA[0][0] = j11; 417c4762a1bSJed Brown bA[0][1] = j12; 418c4762a1bSJed Brown bA[1][0] = j21; 419c4762a1bSJed Brown bA[1][1] = j22; 420c4762a1bSJed Brown 421c4762a1bSJed Brown ierr = MatCreateNest(PETSC_COMM_WORLD,2,NULL,2,NULL,&bA[0][0],&J);CHKERRQ(ierr); 422c4762a1bSJed Brown ierr = MatSetUp(J);CHKERRQ(ierr); 423c4762a1bSJed Brown ierr = MatNestSetVecType(J,VECNEST);CHKERRQ(ierr); 424c4762a1bSJed Brown ierr = MatDestroy(&j11);CHKERRQ(ierr); 425c4762a1bSJed Brown ierr = MatDestroy(&j12);CHKERRQ(ierr); 426c4762a1bSJed Brown ierr = MatDestroy(&j21);CHKERRQ(ierr); 427c4762a1bSJed Brown ierr = MatDestroy(&j22);CHKERRQ(ierr); 428c4762a1bSJed Brown 429c4762a1bSJed Brown ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 430c4762a1bSJed Brown ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 431c4762a1bSJed Brown 432c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL,"-hard",&flg);CHKERRQ(ierr); 433c4762a1bSJed Brown if (!flg) { 434c4762a1bSJed Brown /* 435c4762a1bSJed Brown Set function evaluation routine and vector. 436c4762a1bSJed Brown */ 437c4762a1bSJed Brown ierr = SNESSetFunction(snes,r,FormFunction1_block,NULL);CHKERRQ(ierr); 438c4762a1bSJed Brown 439c4762a1bSJed Brown /* 440c4762a1bSJed Brown Set Jacobian matrix data structure and Jacobian evaluation routine 441c4762a1bSJed Brown */ 442c4762a1bSJed Brown ierr = SNESSetJacobian(snes,J,J,FormJacobian1_block,NULL);CHKERRQ(ierr); 443c4762a1bSJed Brown } else { 444c4762a1bSJed Brown ierr = SNESSetFunction(snes,r,FormFunction2_block,NULL);CHKERRQ(ierr); 445c4762a1bSJed Brown ierr = SNESSetJacobian(snes,J,J,FormJacobian2_block,NULL);CHKERRQ(ierr); 446c4762a1bSJed Brown } 447c4762a1bSJed Brown 448c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 449c4762a1bSJed Brown Customize nonlinear solver; set runtime options 450c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 451c4762a1bSJed Brown 452c4762a1bSJed Brown /* 453c4762a1bSJed Brown Set linear solver defaults for this problem. By extracting the 454c4762a1bSJed Brown KSP, KSP, and PC contexts from the SNES context, we can then 455c4762a1bSJed Brown directly call any KSP, KSP, and PC routines to set various options. 456c4762a1bSJed Brown */ 457c4762a1bSJed Brown ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 458c4762a1bSJed Brown ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 459c4762a1bSJed Brown ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 460c4762a1bSJed Brown ierr = KSPSetTolerances(ksp,1.e-4,PETSC_DEFAULT,PETSC_DEFAULT,20);CHKERRQ(ierr); 461c4762a1bSJed Brown 462c4762a1bSJed Brown /* 463c4762a1bSJed Brown Set SNES/KSP/KSP/PC runtime options, e.g., 464c4762a1bSJed Brown -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc> 465c4762a1bSJed Brown These options will override those specified above as long as 466c4762a1bSJed Brown SNESSetFromOptions() is called _after_ any other customization 467c4762a1bSJed Brown routines. 468c4762a1bSJed Brown */ 469c4762a1bSJed Brown ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 470c4762a1bSJed Brown 471c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 472c4762a1bSJed Brown Evaluate initial guess; then solve nonlinear system 473c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 474c4762a1bSJed Brown if (!flg) { 475c4762a1bSJed Brown ierr = VecSet(x,pfive);CHKERRQ(ierr); 476c4762a1bSJed Brown } else { 477c4762a1bSJed Brown Vec *vecs; 478c4762a1bSJed Brown ierr = VecNestGetSubVecs(x, NULL, &vecs);CHKERRQ(ierr); 479c4762a1bSJed Brown bv = vecs[0]; 480c4762a1bSJed Brown /* ierr = VecBlockGetSubVec(x, 0, &bv);CHKERRQ(ierr); */ 481c4762a1bSJed Brown ierr = VecSetValue(bv, 0, 2.0, INSERT_VALUES);CHKERRQ(ierr); /* xx[0] = 2.0; */ 482c4762a1bSJed Brown ierr = VecAssemblyBegin(bv);CHKERRQ(ierr); 483c4762a1bSJed Brown ierr = VecAssemblyEnd(bv);CHKERRQ(ierr); 484c4762a1bSJed Brown 485c4762a1bSJed Brown /* ierr = VecBlockGetSubVec(x, 1, &bv);CHKERRQ(ierr); */ 486c4762a1bSJed Brown bv = vecs[1]; 487c4762a1bSJed Brown ierr = VecSetValue(bv, 0, 3.0, INSERT_VALUES);CHKERRQ(ierr); /* xx[1] = 3.0; */ 488c4762a1bSJed Brown ierr = VecAssemblyBegin(bv);CHKERRQ(ierr); 489c4762a1bSJed Brown ierr = VecAssemblyEnd(bv);CHKERRQ(ierr); 490c4762a1bSJed Brown } 491c4762a1bSJed Brown /* 492c4762a1bSJed Brown Note: The user should initialize the vector, x, with the initial guess 493c4762a1bSJed Brown for the nonlinear solver prior to calling SNESSolve(). In particular, 494c4762a1bSJed Brown to employ an initial guess of zero, the user should explicitly set 495c4762a1bSJed Brown this vector to zero by calling VecSet(). 496c4762a1bSJed Brown */ 497c4762a1bSJed Brown ierr = SNESSolve(snes,NULL,x);CHKERRQ(ierr); 498c4762a1bSJed Brown ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr); 499c4762a1bSJed Brown if (flg) { 500c4762a1bSJed Brown Vec f; 501c4762a1bSJed Brown ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 502c4762a1bSJed Brown ierr = SNESGetFunction(snes,&f,0,0);CHKERRQ(ierr); 503c4762a1bSJed Brown ierr = VecView(r,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 504c4762a1bSJed Brown } 505c4762a1bSJed Brown 506c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"number of SNES iterations = %D\n\n",its);CHKERRQ(ierr); 507c4762a1bSJed Brown 508c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 509c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 510c4762a1bSJed Brown are no longer needed. 511c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 512c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); ierr = VecDestroy(&r);CHKERRQ(ierr); 513c4762a1bSJed Brown ierr = MatDestroy(&J);CHKERRQ(ierr); ierr = SNESDestroy(&snes);CHKERRQ(ierr); 514c4762a1bSJed Brown PetscFunctionReturn(0); 515c4762a1bSJed Brown } 516c4762a1bSJed Brown 517c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 518c4762a1bSJed Brown static PetscErrorCode FormFunction1_block(SNES snes,Vec x,Vec f,void *dummy) 519c4762a1bSJed Brown { 520c4762a1bSJed Brown Vec *xx, *ff, x1,x2, f1,f2; 521c4762a1bSJed Brown PetscScalar ff_0, ff_1; 522c4762a1bSJed Brown PetscScalar xx_0, xx_1; 523c4762a1bSJed Brown PetscInt index,nb; 524c4762a1bSJed Brown PetscErrorCode ierr; 525c4762a1bSJed Brown 526c4762a1bSJed Brown PetscFunctionBeginUser; 527c4762a1bSJed Brown /* get blocks for function */ 528c4762a1bSJed Brown ierr = VecNestGetSubVecs(f, &nb, &ff);CHKERRQ(ierr); 529c4762a1bSJed Brown f1 = ff[0]; f2 = ff[1]; 530c4762a1bSJed Brown 531c4762a1bSJed Brown /* get blocks for solution */ 532c4762a1bSJed Brown ierr = VecNestGetSubVecs(x, &nb, &xx);CHKERRQ(ierr); 533c4762a1bSJed Brown x1 = xx[0]; x2 = xx[1]; 534c4762a1bSJed Brown 535c4762a1bSJed Brown /* get solution values */ 536c4762a1bSJed Brown index = 0; 537c4762a1bSJed Brown ierr = VecGetValues(x1,1, &index, &xx_0);CHKERRQ(ierr); 538c4762a1bSJed Brown ierr = VecGetValues(x2,1, &index, &xx_1);CHKERRQ(ierr); 539c4762a1bSJed Brown 540c4762a1bSJed Brown /* Compute function */ 541c4762a1bSJed Brown ff_0 = xx_0*xx_0 + xx_0*xx_1 - 3.0; 542c4762a1bSJed Brown ff_1 = xx_0*xx_1 + xx_1*xx_1 - 6.0; 543c4762a1bSJed Brown 544c4762a1bSJed Brown /* set function values */ 545c4762a1bSJed Brown ierr = VecSetValue(f1, index, ff_0, INSERT_VALUES);CHKERRQ(ierr); 546c4762a1bSJed Brown 547c4762a1bSJed Brown ierr = VecSetValue(f2, index, ff_1, INSERT_VALUES);CHKERRQ(ierr); 548c4762a1bSJed Brown 549c4762a1bSJed Brown ierr = VecAssemblyBegin(f);CHKERRQ(ierr); 550c4762a1bSJed Brown ierr = VecAssemblyEnd(f);CHKERRQ(ierr); 551c4762a1bSJed Brown PetscFunctionReturn(0); 552c4762a1bSJed Brown } 553c4762a1bSJed Brown 554c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 555c4762a1bSJed Brown static PetscErrorCode FormJacobian1_block(SNES snes,Vec x,Mat jac,Mat B,void *dummy) 556c4762a1bSJed Brown { 557c4762a1bSJed Brown Vec *xx, x1,x2; 558c4762a1bSJed Brown PetscScalar xx_0, xx_1; 559c4762a1bSJed Brown PetscInt index,nb; 560c4762a1bSJed Brown PetscScalar A_00, A_01, A_10, A_11; 561c4762a1bSJed Brown Mat j11, j12, j21, j22; 562c4762a1bSJed Brown Mat **mats; 563c4762a1bSJed Brown PetscErrorCode ierr; 564c4762a1bSJed Brown 565c4762a1bSJed Brown PetscFunctionBeginUser; 566c4762a1bSJed Brown /* get blocks for solution */ 567c4762a1bSJed Brown ierr = VecNestGetSubVecs(x, &nb, &xx);CHKERRQ(ierr); 568c4762a1bSJed Brown x1 = xx[0]; x2 = xx[1]; 569c4762a1bSJed Brown 570c4762a1bSJed Brown /* get solution values */ 571c4762a1bSJed Brown index = 0; 572c4762a1bSJed Brown ierr = VecGetValues(x1,1, &index, &xx_0);CHKERRQ(ierr); 573c4762a1bSJed Brown ierr = VecGetValues(x2,1, &index, &xx_1);CHKERRQ(ierr); 574c4762a1bSJed Brown 575c4762a1bSJed Brown /* get block matrices */ 576c4762a1bSJed Brown ierr = MatNestGetSubMats(jac,NULL,NULL,&mats);CHKERRQ(ierr); 577c4762a1bSJed Brown j11 = mats[0][0]; 578c4762a1bSJed Brown j12 = mats[0][1]; 579c4762a1bSJed Brown j21 = mats[1][0]; 580c4762a1bSJed Brown j22 = mats[1][1]; 581c4762a1bSJed Brown 582c4762a1bSJed Brown /* compute jacobian entries */ 583c4762a1bSJed Brown A_00 = 2.0*xx_0 + xx_1; 584c4762a1bSJed Brown A_01 = xx_0; 585c4762a1bSJed Brown A_10 = xx_1; 586c4762a1bSJed Brown A_11 = xx_0 + 2.0*xx_1; 587c4762a1bSJed Brown 588c4762a1bSJed Brown /* set jacobian values */ 589c4762a1bSJed Brown ierr = MatSetValue(j11, 0,0, A_00, INSERT_VALUES);CHKERRQ(ierr); 590c4762a1bSJed Brown ierr = MatSetValue(j12, 0,0, A_01, INSERT_VALUES);CHKERRQ(ierr); 591c4762a1bSJed Brown ierr = MatSetValue(j21, 0,0, A_10, INSERT_VALUES);CHKERRQ(ierr); 592c4762a1bSJed Brown ierr = MatSetValue(j22, 0,0, A_11, INSERT_VALUES);CHKERRQ(ierr); 593c4762a1bSJed Brown 594c4762a1bSJed Brown /* Assemble sub matrix */ 595c4762a1bSJed Brown ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 596c4762a1bSJed Brown ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 597c4762a1bSJed Brown PetscFunctionReturn(0); 598c4762a1bSJed Brown } 599c4762a1bSJed Brown 600c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 601c4762a1bSJed Brown static PetscErrorCode FormFunction2_block(SNES snes,Vec x,Vec f,void *dummy) 602c4762a1bSJed Brown { 603c4762a1bSJed Brown PetscErrorCode ierr; 604c4762a1bSJed Brown PetscScalar *ff; 605c4762a1bSJed Brown const PetscScalar *xx; 606c4762a1bSJed Brown 607c4762a1bSJed Brown PetscFunctionBeginUser; 608c4762a1bSJed Brown /* 609c4762a1bSJed Brown Get pointers to vector data. 610c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 611c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 612c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 613c4762a1bSJed Brown the array. 614c4762a1bSJed Brown */ 615c4762a1bSJed Brown ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr); 616c4762a1bSJed Brown ierr = VecGetArray(f,&ff);CHKERRQ(ierr); 617c4762a1bSJed Brown 618c4762a1bSJed Brown /* 619c4762a1bSJed Brown Compute function 620c4762a1bSJed Brown */ 621c4762a1bSJed Brown ff[0] = PetscSinScalar(3.0*xx[0]) + xx[0]; 622c4762a1bSJed Brown ff[1] = xx[1]; 623c4762a1bSJed Brown 624c4762a1bSJed Brown /* 625c4762a1bSJed Brown Restore vectors 626c4762a1bSJed Brown */ 627c4762a1bSJed Brown ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr); 628c4762a1bSJed Brown ierr = VecRestoreArray(f,&ff);CHKERRQ(ierr); 629c4762a1bSJed Brown PetscFunctionReturn(0); 630c4762a1bSJed Brown } 631c4762a1bSJed Brown 632c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 633c4762a1bSJed Brown static PetscErrorCode FormJacobian2_block(SNES snes,Vec x,Mat jac,Mat B,void *dummy) 634c4762a1bSJed Brown { 635c4762a1bSJed Brown const PetscScalar *xx; 636c4762a1bSJed Brown PetscScalar A[4]; 637c4762a1bSJed Brown PetscErrorCode ierr; 638c4762a1bSJed Brown PetscInt idx[2] = {0,1}; 639c4762a1bSJed Brown 640c4762a1bSJed Brown PetscFunctionBeginUser; 641c4762a1bSJed Brown /* 642c4762a1bSJed Brown Get pointer to vector data 643c4762a1bSJed Brown */ 644c4762a1bSJed Brown ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr); 645c4762a1bSJed Brown 646c4762a1bSJed Brown /* 647c4762a1bSJed Brown Compute Jacobian entries and insert into matrix. 648c4762a1bSJed Brown - Since this is such a small problem, we set all entries for 649c4762a1bSJed Brown the matrix at once. 650c4762a1bSJed Brown */ 651c4762a1bSJed Brown A[0] = 3.0*PetscCosScalar(3.0*xx[0]) + 1.0; A[1] = 0.0; 652c4762a1bSJed Brown A[2] = 0.0; A[3] = 1.0; 653c4762a1bSJed Brown ierr = MatSetValues(jac,2,idx,2,idx,A,INSERT_VALUES);CHKERRQ(ierr); 654c4762a1bSJed Brown 655c4762a1bSJed Brown /* 656c4762a1bSJed Brown Restore vector 657c4762a1bSJed Brown */ 658c4762a1bSJed Brown ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr); 659c4762a1bSJed Brown 660c4762a1bSJed Brown /* 661c4762a1bSJed Brown Assemble matrix 662c4762a1bSJed Brown */ 663c4762a1bSJed Brown ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 664c4762a1bSJed Brown ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 665c4762a1bSJed Brown PetscFunctionReturn(0); 666c4762a1bSJed Brown } 667c4762a1bSJed Brown 668c4762a1bSJed Brown int main(int argc,char **argv) 669c4762a1bSJed Brown { 670c4762a1bSJed Brown PetscMPIInt size; 671c4762a1bSJed Brown PetscErrorCode ierr; 672c4762a1bSJed Brown 673c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 674c4762a1bSJed Brown 675c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 676c4762a1bSJed Brown if (size != 1) SETERRQ(PETSC_COMM_WORLD, 1,"This is a uniprocessor example only!"); 677c4762a1bSJed Brown 678c4762a1bSJed Brown ierr = assembled_system();CHKERRQ(ierr); 679c4762a1bSJed Brown 680c4762a1bSJed Brown ierr = block_system();CHKERRQ(ierr); 681c4762a1bSJed Brown 682c4762a1bSJed Brown ierr = PetscFinalize(); 683c4762a1bSJed Brown return ierr; 684c4762a1bSJed Brown } 685c4762a1bSJed Brown 686c4762a1bSJed Brown 687c4762a1bSJed Brown /*TEST 688c4762a1bSJed Brown 689c4762a1bSJed Brown test: 690c4762a1bSJed Brown args: -snes_monitor_short 691c4762a1bSJed Brown requires: !single 692c4762a1bSJed Brown 693c4762a1bSJed Brown TEST*/ 694