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 /* 5c4762a1bSJed Brown Include "petscsnes.h" so that we can use SNES solvers. Note that this 6c4762a1bSJed Brown file automatically includes: 7c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 8c4762a1bSJed Brown petscsys.h - system routines petscmat.h - matrices 9c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods 10c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners 11c4762a1bSJed Brown petscksp.h - linear solvers 12c4762a1bSJed Brown */ 13c4762a1bSJed Brown #include <petscsnes.h> 14c4762a1bSJed Brown 15c4762a1bSJed Brown /* 16c4762a1bSJed Brown This example is block version of the test found at 17c4762a1bSJed Brown ${PETSC_DIR}/src/snes/tutorials/ex1.c 18c4762a1bSJed Brown In this test we replace the Jacobian systems 19c4762a1bSJed Brown [J]{x} = {F} 20c4762a1bSJed Brown where 21c4762a1bSJed Brown 22c4762a1bSJed Brown [J] = (j_00, j_01), {x} = (x_0, x_1)^T, {F} = (f_0, f_1)^T 23c4762a1bSJed Brown (j_10, j_11) 24c4762a1bSJed Brown where [J] \in \mathbb^{2 \times 2}, {x},{F} \in \mathbb^{2 \times 1}, 25c4762a1bSJed Brown 26c4762a1bSJed Brown with a block system in which each block is of length 1. 27c4762a1bSJed Brown i.e. The block system is thus 28c4762a1bSJed Brown 29c4762a1bSJed Brown [J] = ([j00], [j01]), {x} = ({x0}, {x1})^T, {F} = ({f0}, {f1})^T 30c4762a1bSJed Brown ([j10], [j11]) 31c4762a1bSJed Brown where 32c4762a1bSJed Brown [j00], [j01], [j10], [j11] \in \mathbb^{1 \times 1} 33c4762a1bSJed Brown {x0}, {x1}, {f0}, {f1} \in \mathbb^{1 \times 1} 34c4762a1bSJed Brown 35c4762a1bSJed Brown In practice we would not bother defing blocks of size one, and would instead assemble the 36c4762a1bSJed Brown full system. This is just a simple test to illustrate how to manipulate the blocks and 37c4762a1bSJed Brown to confirm the implementation is correct. 38c4762a1bSJed Brown */ 39c4762a1bSJed Brown 40c4762a1bSJed Brown /* 41c4762a1bSJed Brown User-defined routines 42c4762a1bSJed Brown */ 43c4762a1bSJed Brown static PetscErrorCode FormJacobian1(SNES,Vec,Mat,Mat,void*); 44c4762a1bSJed Brown static PetscErrorCode FormFunction1(SNES,Vec,Vec,void*); 45c4762a1bSJed Brown static PetscErrorCode FormJacobian2(SNES,Vec,Mat,Mat,void*); 46c4762a1bSJed Brown static PetscErrorCode FormFunction2(SNES,Vec,Vec,void*); 47c4762a1bSJed Brown static PetscErrorCode FormJacobian1_block(SNES,Vec,Mat,Mat,void*); 48c4762a1bSJed Brown static PetscErrorCode FormFunction1_block(SNES,Vec,Vec,void*); 49c4762a1bSJed Brown static PetscErrorCode FormJacobian2_block(SNES,Vec,Mat,Mat,void*); 50c4762a1bSJed Brown static PetscErrorCode FormFunction2_block(SNES,Vec,Vec,void*); 51c4762a1bSJed Brown 52c4762a1bSJed Brown static PetscErrorCode assembled_system(void) 53c4762a1bSJed Brown { 54c4762a1bSJed Brown SNES snes; /* nonlinear solver context */ 55c4762a1bSJed Brown KSP ksp; /* linear solver context */ 56c4762a1bSJed Brown PC pc; /* preconditioner context */ 57c4762a1bSJed Brown Vec x,r; /* solution, residual vectors */ 58c4762a1bSJed Brown Mat J; /* Jacobian matrix */ 59c4762a1bSJed Brown PetscInt its; 60c4762a1bSJed Brown PetscScalar pfive = .5,*xx; 61c4762a1bSJed Brown PetscBool flg; 62c4762a1bSJed Brown 63c4762a1bSJed Brown PetscFunctionBeginUser; 649566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n\n========================= Assembled system =========================\n\n")); 65c4762a1bSJed Brown 66c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 67c4762a1bSJed Brown Create nonlinear solver context 68c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 69c4762a1bSJed Brown 709566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD,&snes)); 71c4762a1bSJed Brown 72c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 73c4762a1bSJed Brown Create matrix and vector data structures; set corresponding routines 74c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 75c4762a1bSJed Brown 76c4762a1bSJed Brown /* 77c4762a1bSJed Brown Create vectors for solution and nonlinear function 78c4762a1bSJed Brown */ 799566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,2,&x)); 809566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&r)); 81c4762a1bSJed Brown 82c4762a1bSJed Brown /* 83c4762a1bSJed Brown Create Jacobian matrix data structure 84c4762a1bSJed Brown */ 859566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF,&J)); 869566063dSJacob Faibussowitsch PetscCall(MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,2,2)); 879566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(J)); 889566063dSJacob Faibussowitsch PetscCall(MatSetUp(J)); 89c4762a1bSJed Brown 909566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-hard",&flg)); 91c4762a1bSJed Brown if (!flg) { 92c4762a1bSJed Brown /* 93c4762a1bSJed Brown Set function evaluation routine and vector. 94c4762a1bSJed Brown */ 959566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes,r,FormFunction1,NULL)); 96c4762a1bSJed Brown 97c4762a1bSJed Brown /* 98c4762a1bSJed Brown Set Jacobian matrix data structure and Jacobian evaluation routine 99c4762a1bSJed Brown */ 1009566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes,J,J,FormJacobian1,NULL)); 101c4762a1bSJed Brown } else { 1029566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes,r,FormFunction2,NULL)); 1039566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes,J,J,FormJacobian2,NULL)); 104c4762a1bSJed Brown } 105c4762a1bSJed Brown 106c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 107c4762a1bSJed Brown Customize nonlinear solver; set runtime options 108c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 109c4762a1bSJed Brown 110c4762a1bSJed Brown /* 111c4762a1bSJed Brown Set linear solver defaults for this problem. By extracting the 112c4762a1bSJed Brown KSP, KSP, and PC contexts from the SNES context, we can then 113c4762a1bSJed Brown directly call any KSP, KSP, and PC routines to set various options. 114c4762a1bSJed Brown */ 1159566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes,&ksp)); 1169566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp,&pc)); 1179566063dSJacob Faibussowitsch PetscCall(PCSetType(pc,PCNONE)); 1189566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(ksp,1.e-4,PETSC_DEFAULT,PETSC_DEFAULT,20)); 119c4762a1bSJed Brown 120c4762a1bSJed Brown /* 121c4762a1bSJed Brown Set SNES/KSP/KSP/PC runtime options, e.g., 122c4762a1bSJed Brown -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc> 123c4762a1bSJed Brown These options will override those specified above as long as 124c4762a1bSJed Brown SNESSetFromOptions() is called _after_ any other customization 125c4762a1bSJed Brown routines. 126c4762a1bSJed Brown */ 1279566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 128c4762a1bSJed Brown 129c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 130c4762a1bSJed Brown Evaluate initial guess; then solve nonlinear system 131c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 132c4762a1bSJed Brown if (!flg) { 1339566063dSJacob Faibussowitsch PetscCall(VecSet(x,pfive)); 134c4762a1bSJed Brown } else { 1359566063dSJacob Faibussowitsch PetscCall(VecGetArray(x,&xx)); 136c4762a1bSJed Brown xx[0] = 2.0; xx[1] = 3.0; 1379566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x,&xx)); 138c4762a1bSJed Brown } 139c4762a1bSJed Brown /* 140c4762a1bSJed Brown Note: The user should initialize the vector, x, with the initial guess 141c4762a1bSJed Brown for the nonlinear solver prior to calling SNESSolve(). In particular, 142c4762a1bSJed Brown to employ an initial guess of zero, the user should explicitly set 143c4762a1bSJed Brown this vector to zero by calling VecSet(). 144c4762a1bSJed Brown */ 145c4762a1bSJed Brown 1469566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes,NULL,x)); 1479566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes,&its)); 148c4762a1bSJed Brown if (flg) { 149c4762a1bSJed Brown Vec f; 1509566063dSJacob Faibussowitsch PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); 1519566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(snes,&f,0,0)); 1529566063dSJacob Faibussowitsch PetscCall(VecView(r,PETSC_VIEWER_STDOUT_WORLD)); 153c4762a1bSJed Brown } 15463a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"number of SNES iterations = %" PetscInt_FMT "\n\n",its)); 155c4762a1bSJed Brown 156c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 157c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 158c4762a1bSJed Brown are no longer needed. 159c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 160c4762a1bSJed Brown 1619566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); PetscCall(VecDestroy(&r)); 1629566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); PetscCall(SNESDestroy(&snes)); 163c4762a1bSJed Brown PetscFunctionReturn(0); 164c4762a1bSJed Brown } 165c4762a1bSJed Brown 166c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 167c4762a1bSJed Brown /* 168c4762a1bSJed Brown FormFunction1 - Evaluates nonlinear function, F(x). 169c4762a1bSJed Brown 170c4762a1bSJed Brown Input Parameters: 171c4762a1bSJed Brown . snes - the SNES context 172c4762a1bSJed Brown . x - input vector 173c4762a1bSJed Brown . dummy - optional user-defined context (not used here) 174c4762a1bSJed Brown 175c4762a1bSJed Brown Output Parameter: 176c4762a1bSJed Brown . f - function vector 177c4762a1bSJed Brown */ 178c4762a1bSJed Brown static PetscErrorCode FormFunction1(SNES snes,Vec x,Vec f,void *dummy) 179c4762a1bSJed Brown { 180c4762a1bSJed Brown const PetscScalar *xx; 181c4762a1bSJed Brown PetscScalar *ff; 182c4762a1bSJed Brown 183c4762a1bSJed Brown PetscFunctionBeginUser; 184c4762a1bSJed Brown /* 185c4762a1bSJed Brown Get pointers to vector data. 186c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 187c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 188c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 189c4762a1bSJed Brown the array. 190c4762a1bSJed Brown */ 1919566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x,&xx)); 1929566063dSJacob Faibussowitsch PetscCall(VecGetArray(f,&ff)); 193c4762a1bSJed Brown 194c4762a1bSJed Brown /* 195c4762a1bSJed Brown Compute function 196c4762a1bSJed Brown */ 197c4762a1bSJed Brown ff[0] = xx[0]*xx[0] + xx[0]*xx[1] - 3.0; 198c4762a1bSJed Brown ff[1] = xx[0]*xx[1] + xx[1]*xx[1] - 6.0; 199c4762a1bSJed Brown 200c4762a1bSJed Brown /* 201c4762a1bSJed Brown Restore vectors 202c4762a1bSJed Brown */ 2039566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x,&xx)); 2049566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(f,&ff)); 205c4762a1bSJed Brown PetscFunctionReturn(0); 206c4762a1bSJed Brown } 207c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 208c4762a1bSJed Brown /* 209c4762a1bSJed Brown FormJacobian1 - Evaluates Jacobian matrix. 210c4762a1bSJed Brown 211c4762a1bSJed Brown Input Parameters: 212c4762a1bSJed Brown . snes - the SNES context 213c4762a1bSJed Brown . x - input vector 214c4762a1bSJed Brown . dummy - optional user-defined context (not used here) 215c4762a1bSJed Brown 216c4762a1bSJed Brown Output Parameters: 217c4762a1bSJed Brown . jac - Jacobian matrix 218c4762a1bSJed Brown . B - optionally different preconditioning matrix 219c4762a1bSJed Brown . flag - flag indicating matrix structure 220c4762a1bSJed Brown */ 221c4762a1bSJed Brown static PetscErrorCode FormJacobian1(SNES snes,Vec x,Mat jac,Mat B,void *dummy) 222c4762a1bSJed Brown { 223c4762a1bSJed Brown const PetscScalar *xx; 224c4762a1bSJed Brown PetscScalar A[4]; 225c4762a1bSJed Brown PetscInt idx[2] = {0,1}; 226c4762a1bSJed Brown 227c4762a1bSJed Brown PetscFunctionBeginUser; 228c4762a1bSJed Brown /* 229c4762a1bSJed Brown Get pointer to vector data 230c4762a1bSJed Brown */ 2319566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x,&xx)); 232c4762a1bSJed Brown 233c4762a1bSJed Brown /* 234c4762a1bSJed Brown Compute Jacobian entries and insert into matrix. 235c4762a1bSJed Brown - Since this is such a small problem, we set all entries for 236c4762a1bSJed Brown the matrix at once. 237c4762a1bSJed Brown */ 238c4762a1bSJed Brown A[0] = 2.0*xx[0] + xx[1]; A[1] = xx[0]; 239c4762a1bSJed Brown A[2] = xx[1]; A[3] = xx[0] + 2.0*xx[1]; 2409566063dSJacob Faibussowitsch PetscCall(MatSetValues(jac,2,idx,2,idx,A,INSERT_VALUES)); 241c4762a1bSJed Brown 242c4762a1bSJed Brown /* 243c4762a1bSJed Brown Restore vector 244c4762a1bSJed Brown */ 2459566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x,&xx)); 246c4762a1bSJed Brown 247c4762a1bSJed Brown /* 248c4762a1bSJed Brown Assemble matrix 249c4762a1bSJed Brown */ 2509566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY)); 2519566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY)); 252c4762a1bSJed Brown PetscFunctionReturn(0); 253c4762a1bSJed Brown } 254c4762a1bSJed Brown 255c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 256c4762a1bSJed Brown static PetscErrorCode FormFunction2(SNES snes,Vec x,Vec f,void *dummy) 257c4762a1bSJed Brown { 258c4762a1bSJed Brown const PetscScalar *xx; 259c4762a1bSJed Brown PetscScalar *ff; 260c4762a1bSJed Brown 261c4762a1bSJed Brown PetscFunctionBeginUser; 262c4762a1bSJed Brown /* 263c4762a1bSJed Brown Get pointers to vector data. 264c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 265c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 266c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 267c4762a1bSJed Brown the array. 268c4762a1bSJed Brown */ 2699566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x,&xx)); 2709566063dSJacob Faibussowitsch PetscCall(VecGetArray(f,&ff)); 271c4762a1bSJed Brown 272c4762a1bSJed Brown /* 273c4762a1bSJed Brown Compute function 274c4762a1bSJed Brown */ 275c4762a1bSJed Brown ff[0] = PetscSinScalar(3.0*xx[0]) + xx[0]; 276c4762a1bSJed Brown ff[1] = xx[1]; 277c4762a1bSJed Brown 278c4762a1bSJed Brown /* 279c4762a1bSJed Brown Restore vectors 280c4762a1bSJed Brown */ 2819566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x,&xx)); 2829566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(f,&ff)); 283c4762a1bSJed Brown PetscFunctionReturn(0); 284c4762a1bSJed Brown } 285c4762a1bSJed Brown 286c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 287c4762a1bSJed Brown static PetscErrorCode FormJacobian2(SNES snes,Vec x,Mat jac,Mat B,void *dummy) 288c4762a1bSJed Brown { 289c4762a1bSJed Brown const PetscScalar *xx; 290c4762a1bSJed Brown PetscScalar A[4]; 291c4762a1bSJed Brown PetscInt idx[2] = {0,1}; 292c4762a1bSJed Brown 293c4762a1bSJed Brown PetscFunctionBeginUser; 294c4762a1bSJed Brown /* 295c4762a1bSJed Brown Get pointer to vector data 296c4762a1bSJed Brown */ 2979566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x,&xx)); 298c4762a1bSJed Brown 299c4762a1bSJed Brown /* 300c4762a1bSJed Brown Compute Jacobian entries and insert into matrix. 301c4762a1bSJed Brown - Since this is such a small problem, we set all entries for 302c4762a1bSJed Brown the matrix at once. 303c4762a1bSJed Brown */ 304c4762a1bSJed Brown A[0] = 3.0*PetscCosScalar(3.0*xx[0]) + 1.0; A[1] = 0.0; 305c4762a1bSJed Brown A[2] = 0.0; A[3] = 1.0; 3069566063dSJacob Faibussowitsch PetscCall(MatSetValues(jac,2,idx,2,idx,A,INSERT_VALUES)); 307c4762a1bSJed Brown 308c4762a1bSJed Brown /* 309c4762a1bSJed Brown Restore vector 310c4762a1bSJed Brown */ 3119566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x,&xx)); 312c4762a1bSJed Brown 313c4762a1bSJed Brown /* 314c4762a1bSJed Brown Assemble matrix 315c4762a1bSJed Brown */ 3169566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY)); 3179566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY)); 318c4762a1bSJed Brown PetscFunctionReturn(0); 319c4762a1bSJed Brown } 320c4762a1bSJed Brown 321b5675b0fSBarry Smith static PetscErrorCode block_system(void) 322c4762a1bSJed Brown { 323c4762a1bSJed Brown SNES snes; /* nonlinear solver context */ 324c4762a1bSJed Brown KSP ksp; /* linear solver context */ 325c4762a1bSJed Brown PC pc; /* preconditioner context */ 326c4762a1bSJed Brown Vec x,r; /* solution, residual vectors */ 327c4762a1bSJed Brown Mat J; /* Jacobian matrix */ 328c4762a1bSJed Brown PetscInt its; 329c4762a1bSJed Brown PetscScalar pfive = .5; 330c4762a1bSJed Brown PetscBool flg; 331c4762a1bSJed Brown 332c4762a1bSJed Brown Mat j11, j12, j21, j22; 333c4762a1bSJed Brown Vec x1, x2, r1, r2; 334c4762a1bSJed Brown Vec bv; 335c4762a1bSJed Brown Vec bx[2]; 336c4762a1bSJed Brown Mat bA[2][2]; 337c4762a1bSJed Brown 338c4762a1bSJed Brown PetscFunctionBeginUser; 3399566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n\n========================= Block system =========================\n\n")); 340c4762a1bSJed Brown 3419566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD,&snes)); 342c4762a1bSJed Brown 343c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 344c4762a1bSJed Brown Create matrix and vector data structures; set corresponding routines 345c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 346c4762a1bSJed Brown 347c4762a1bSJed Brown /* 348c4762a1bSJed Brown Create sub vectors for solution and nonlinear function 349c4762a1bSJed Brown */ 3509566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,1,&x1)); 3519566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x1,&r1)); 352c4762a1bSJed Brown 3539566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,1,&x2)); 3549566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x2,&r2)); 355c4762a1bSJed Brown 356c4762a1bSJed Brown /* 357c4762a1bSJed Brown Create the block vectors 358c4762a1bSJed Brown */ 359c4762a1bSJed Brown bx[0] = x1; 360c4762a1bSJed Brown bx[1] = x2; 3619566063dSJacob Faibussowitsch PetscCall(VecCreateNest(PETSC_COMM_WORLD,2,NULL,bx,&x)); 3629566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(x)); 3639566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(x)); 3649566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x1)); 3659566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x2)); 366c4762a1bSJed Brown 367c4762a1bSJed Brown bx[0] = r1; 368c4762a1bSJed Brown bx[1] = r2; 3699566063dSJacob Faibussowitsch PetscCall(VecCreateNest(PETSC_COMM_WORLD,2,NULL,bx,&r)); 3709566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r1)); 3719566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r2)); 3729566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(r)); 3739566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(r)); 374c4762a1bSJed Brown 375c4762a1bSJed Brown /* 376c4762a1bSJed Brown Create sub Jacobian matrix data structure 377c4762a1bSJed Brown */ 3789566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &j11)); 3799566063dSJacob Faibussowitsch PetscCall(MatSetSizes(j11, 1, 1, 1, 1)); 3809566063dSJacob Faibussowitsch PetscCall(MatSetType(j11, MATSEQAIJ)); 3819566063dSJacob Faibussowitsch PetscCall(MatSetUp(j11)); 382c4762a1bSJed Brown 3839566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &j12)); 3849566063dSJacob Faibussowitsch PetscCall(MatSetSizes(j12, 1, 1, 1, 1)); 3859566063dSJacob Faibussowitsch PetscCall(MatSetType(j12, MATSEQAIJ)); 3869566063dSJacob Faibussowitsch PetscCall(MatSetUp(j12)); 387c4762a1bSJed Brown 3889566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &j21)); 3899566063dSJacob Faibussowitsch PetscCall(MatSetSizes(j21, 1, 1, 1, 1)); 3909566063dSJacob Faibussowitsch PetscCall(MatSetType(j21, MATSEQAIJ)); 3919566063dSJacob Faibussowitsch PetscCall(MatSetUp(j21)); 392c4762a1bSJed Brown 3939566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &j22)); 3949566063dSJacob Faibussowitsch PetscCall(MatSetSizes(j22, PETSC_DECIDE, PETSC_DECIDE, 1, 1)); 3959566063dSJacob Faibussowitsch PetscCall(MatSetType(j22, MATSEQAIJ)); 3969566063dSJacob Faibussowitsch PetscCall(MatSetUp(j22)); 397c4762a1bSJed Brown /* 398c4762a1bSJed Brown Create block Jacobian matrix data structure 399c4762a1bSJed Brown */ 400c4762a1bSJed Brown bA[0][0] = j11; 401c4762a1bSJed Brown bA[0][1] = j12; 402c4762a1bSJed Brown bA[1][0] = j21; 403c4762a1bSJed Brown bA[1][1] = j22; 404c4762a1bSJed Brown 4059566063dSJacob Faibussowitsch PetscCall(MatCreateNest(PETSC_COMM_WORLD,2,NULL,2,NULL,&bA[0][0],&J)); 4069566063dSJacob Faibussowitsch PetscCall(MatSetUp(J)); 4079566063dSJacob Faibussowitsch PetscCall(MatNestSetVecType(J,VECNEST)); 4089566063dSJacob Faibussowitsch PetscCall(MatDestroy(&j11)); 4099566063dSJacob Faibussowitsch PetscCall(MatDestroy(&j12)); 4109566063dSJacob Faibussowitsch PetscCall(MatDestroy(&j21)); 4119566063dSJacob Faibussowitsch PetscCall(MatDestroy(&j22)); 412c4762a1bSJed Brown 4139566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY)); 4149566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY)); 415c4762a1bSJed Brown 4169566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-hard",&flg)); 417c4762a1bSJed Brown if (!flg) { 418c4762a1bSJed Brown /* 419c4762a1bSJed Brown Set function evaluation routine and vector. 420c4762a1bSJed Brown */ 4219566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes,r,FormFunction1_block,NULL)); 422c4762a1bSJed Brown 423c4762a1bSJed Brown /* 424c4762a1bSJed Brown Set Jacobian matrix data structure and Jacobian evaluation routine 425c4762a1bSJed Brown */ 4269566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes,J,J,FormJacobian1_block,NULL)); 427c4762a1bSJed Brown } else { 4289566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes,r,FormFunction2_block,NULL)); 4299566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes,J,J,FormJacobian2_block,NULL)); 430c4762a1bSJed Brown } 431c4762a1bSJed Brown 432c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 433c4762a1bSJed Brown Customize nonlinear solver; set runtime options 434c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 435c4762a1bSJed Brown 436c4762a1bSJed Brown /* 437c4762a1bSJed Brown Set linear solver defaults for this problem. By extracting the 438c4762a1bSJed Brown KSP, KSP, and PC contexts from the SNES context, we can then 439c4762a1bSJed Brown directly call any KSP, KSP, and PC routines to set various options. 440c4762a1bSJed Brown */ 4419566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes,&ksp)); 4429566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp,&pc)); 4439566063dSJacob Faibussowitsch PetscCall(PCSetType(pc,PCNONE)); 4449566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(ksp,1.e-4,PETSC_DEFAULT,PETSC_DEFAULT,20)); 445c4762a1bSJed Brown 446c4762a1bSJed Brown /* 447c4762a1bSJed Brown Set SNES/KSP/KSP/PC runtime options, e.g., 448c4762a1bSJed Brown -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc> 449c4762a1bSJed Brown These options will override those specified above as long as 450c4762a1bSJed Brown SNESSetFromOptions() is called _after_ any other customization 451c4762a1bSJed Brown routines. 452c4762a1bSJed Brown */ 4539566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 454c4762a1bSJed Brown 455c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 456c4762a1bSJed Brown Evaluate initial guess; then solve nonlinear system 457c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 458c4762a1bSJed Brown if (!flg) { 4599566063dSJacob Faibussowitsch PetscCall(VecSet(x,pfive)); 460c4762a1bSJed Brown } else { 461c4762a1bSJed Brown Vec *vecs; 4629566063dSJacob Faibussowitsch PetscCall(VecNestGetSubVecs(x, NULL, &vecs)); 463c4762a1bSJed Brown bv = vecs[0]; 4649566063dSJacob Faibussowitsch /* PetscCall(VecBlockGetSubVec(x, 0, &bv)); */ 4659566063dSJacob Faibussowitsch PetscCall(VecSetValue(bv, 0, 2.0, INSERT_VALUES)); /* xx[0] = 2.0; */ 4669566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(bv)); 4679566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(bv)); 468c4762a1bSJed Brown 4699566063dSJacob Faibussowitsch /* PetscCall(VecBlockGetSubVec(x, 1, &bv)); */ 470c4762a1bSJed Brown bv = vecs[1]; 4719566063dSJacob Faibussowitsch PetscCall(VecSetValue(bv, 0, 3.0, INSERT_VALUES)); /* xx[1] = 3.0; */ 4729566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(bv)); 4739566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(bv)); 474c4762a1bSJed Brown } 475c4762a1bSJed Brown /* 476c4762a1bSJed Brown Note: The user should initialize the vector, x, with the initial guess 477c4762a1bSJed Brown for the nonlinear solver prior to calling SNESSolve(). In particular, 478c4762a1bSJed Brown to employ an initial guess of zero, the user should explicitly set 479c4762a1bSJed Brown this vector to zero by calling VecSet(). 480c4762a1bSJed Brown */ 4819566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes,NULL,x)); 4829566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes,&its)); 483c4762a1bSJed Brown if (flg) { 484c4762a1bSJed Brown Vec f; 4859566063dSJacob Faibussowitsch PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); 4869566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(snes,&f,0,0)); 4879566063dSJacob Faibussowitsch PetscCall(VecView(r,PETSC_VIEWER_STDOUT_WORLD)); 488c4762a1bSJed Brown } 489c4762a1bSJed Brown 49063a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"number of SNES iterations = %" PetscInt_FMT "\n\n",its)); 491c4762a1bSJed Brown 492c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 493c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 494c4762a1bSJed Brown are no longer needed. 495c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 4969566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); PetscCall(VecDestroy(&r)); 4979566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); PetscCall(SNESDestroy(&snes)); 498c4762a1bSJed Brown PetscFunctionReturn(0); 499c4762a1bSJed Brown } 500c4762a1bSJed Brown 501c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 502c4762a1bSJed Brown static PetscErrorCode FormFunction1_block(SNES snes,Vec x,Vec f,void *dummy) 503c4762a1bSJed Brown { 504c4762a1bSJed Brown Vec *xx, *ff, x1,x2, f1,f2; 505c4762a1bSJed Brown PetscScalar ff_0, ff_1; 506c4762a1bSJed Brown PetscScalar xx_0, xx_1; 507c4762a1bSJed Brown PetscInt index,nb; 508c4762a1bSJed Brown 509c4762a1bSJed Brown PetscFunctionBeginUser; 510c4762a1bSJed Brown /* get blocks for function */ 5119566063dSJacob Faibussowitsch PetscCall(VecNestGetSubVecs(f, &nb, &ff)); 512c4762a1bSJed Brown f1 = ff[0]; f2 = ff[1]; 513c4762a1bSJed Brown 514c4762a1bSJed Brown /* get blocks for solution */ 5159566063dSJacob Faibussowitsch PetscCall(VecNestGetSubVecs(x, &nb, &xx)); 516c4762a1bSJed Brown x1 = xx[0]; x2 = xx[1]; 517c4762a1bSJed Brown 518c4762a1bSJed Brown /* get solution values */ 519c4762a1bSJed Brown index = 0; 5209566063dSJacob Faibussowitsch PetscCall(VecGetValues(x1,1, &index, &xx_0)); 5219566063dSJacob Faibussowitsch PetscCall(VecGetValues(x2,1, &index, &xx_1)); 522c4762a1bSJed Brown 523c4762a1bSJed Brown /* Compute function */ 524c4762a1bSJed Brown ff_0 = xx_0*xx_0 + xx_0*xx_1 - 3.0; 525c4762a1bSJed Brown ff_1 = xx_0*xx_1 + xx_1*xx_1 - 6.0; 526c4762a1bSJed Brown 527c4762a1bSJed Brown /* set function values */ 5289566063dSJacob Faibussowitsch PetscCall(VecSetValue(f1, index, ff_0, INSERT_VALUES)); 529c4762a1bSJed Brown 5309566063dSJacob Faibussowitsch PetscCall(VecSetValue(f2, index, ff_1, INSERT_VALUES)); 531c4762a1bSJed Brown 5329566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(f)); 5339566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(f)); 534c4762a1bSJed Brown PetscFunctionReturn(0); 535c4762a1bSJed Brown } 536c4762a1bSJed Brown 537c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 538c4762a1bSJed Brown static PetscErrorCode FormJacobian1_block(SNES snes,Vec x,Mat jac,Mat B,void *dummy) 539c4762a1bSJed Brown { 540c4762a1bSJed Brown Vec *xx, x1,x2; 541c4762a1bSJed Brown PetscScalar xx_0, xx_1; 542c4762a1bSJed Brown PetscInt index,nb; 543c4762a1bSJed Brown PetscScalar A_00, A_01, A_10, A_11; 544c4762a1bSJed Brown Mat j11, j12, j21, j22; 545c4762a1bSJed Brown Mat **mats; 546c4762a1bSJed Brown 547c4762a1bSJed Brown PetscFunctionBeginUser; 548c4762a1bSJed Brown /* get blocks for solution */ 5499566063dSJacob Faibussowitsch PetscCall(VecNestGetSubVecs(x, &nb, &xx)); 550c4762a1bSJed Brown x1 = xx[0]; x2 = xx[1]; 551c4762a1bSJed Brown 552c4762a1bSJed Brown /* get solution values */ 553c4762a1bSJed Brown index = 0; 5549566063dSJacob Faibussowitsch PetscCall(VecGetValues(x1,1, &index, &xx_0)); 5559566063dSJacob Faibussowitsch PetscCall(VecGetValues(x2,1, &index, &xx_1)); 556c4762a1bSJed Brown 557c4762a1bSJed Brown /* get block matrices */ 5589566063dSJacob Faibussowitsch PetscCall(MatNestGetSubMats(jac,NULL,NULL,&mats)); 559c4762a1bSJed Brown j11 = mats[0][0]; 560c4762a1bSJed Brown j12 = mats[0][1]; 561c4762a1bSJed Brown j21 = mats[1][0]; 562c4762a1bSJed Brown j22 = mats[1][1]; 563c4762a1bSJed Brown 564c4762a1bSJed Brown /* compute jacobian entries */ 565c4762a1bSJed Brown A_00 = 2.0*xx_0 + xx_1; 566c4762a1bSJed Brown A_01 = xx_0; 567c4762a1bSJed Brown A_10 = xx_1; 568c4762a1bSJed Brown A_11 = xx_0 + 2.0*xx_1; 569c4762a1bSJed Brown 570c4762a1bSJed Brown /* set jacobian values */ 5719566063dSJacob Faibussowitsch PetscCall(MatSetValue(j11, 0,0, A_00, INSERT_VALUES)); 5729566063dSJacob Faibussowitsch PetscCall(MatSetValue(j12, 0,0, A_01, INSERT_VALUES)); 5739566063dSJacob Faibussowitsch PetscCall(MatSetValue(j21, 0,0, A_10, INSERT_VALUES)); 5749566063dSJacob Faibussowitsch PetscCall(MatSetValue(j22, 0,0, A_11, INSERT_VALUES)); 575c4762a1bSJed Brown 576c4762a1bSJed Brown /* Assemble sub matrix */ 5779566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY)); 5789566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY)); 579c4762a1bSJed Brown PetscFunctionReturn(0); 580c4762a1bSJed Brown } 581c4762a1bSJed Brown 582c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 583c4762a1bSJed Brown static PetscErrorCode FormFunction2_block(SNES snes,Vec x,Vec f,void *dummy) 584c4762a1bSJed Brown { 585c4762a1bSJed Brown PetscScalar *ff; 586c4762a1bSJed Brown const PetscScalar *xx; 587c4762a1bSJed Brown 588c4762a1bSJed Brown PetscFunctionBeginUser; 589c4762a1bSJed Brown /* 590c4762a1bSJed Brown Get pointers to vector data. 591c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 592c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 593c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 594c4762a1bSJed Brown the array. 595c4762a1bSJed Brown */ 5969566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x,&xx)); 5979566063dSJacob Faibussowitsch PetscCall(VecGetArray(f,&ff)); 598c4762a1bSJed Brown 599c4762a1bSJed Brown /* 600c4762a1bSJed Brown Compute function 601c4762a1bSJed Brown */ 602c4762a1bSJed Brown ff[0] = PetscSinScalar(3.0*xx[0]) + xx[0]; 603c4762a1bSJed Brown ff[1] = xx[1]; 604c4762a1bSJed Brown 605c4762a1bSJed Brown /* 606c4762a1bSJed Brown Restore vectors 607c4762a1bSJed Brown */ 6089566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x,&xx)); 6099566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(f,&ff)); 610c4762a1bSJed Brown PetscFunctionReturn(0); 611c4762a1bSJed Brown } 612c4762a1bSJed Brown 613c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 614c4762a1bSJed Brown static PetscErrorCode FormJacobian2_block(SNES snes,Vec x,Mat jac,Mat B,void *dummy) 615c4762a1bSJed Brown { 616c4762a1bSJed Brown const PetscScalar *xx; 617c4762a1bSJed Brown PetscScalar A[4]; 618c4762a1bSJed Brown PetscInt idx[2] = {0,1}; 619c4762a1bSJed Brown 620c4762a1bSJed Brown PetscFunctionBeginUser; 621c4762a1bSJed Brown /* 622c4762a1bSJed Brown Get pointer to vector data 623c4762a1bSJed Brown */ 6249566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x,&xx)); 625c4762a1bSJed Brown 626c4762a1bSJed Brown /* 627c4762a1bSJed Brown Compute Jacobian entries and insert into matrix. 628c4762a1bSJed Brown - Since this is such a small problem, we set all entries for 629c4762a1bSJed Brown the matrix at once. 630c4762a1bSJed Brown */ 631c4762a1bSJed Brown A[0] = 3.0*PetscCosScalar(3.0*xx[0]) + 1.0; A[1] = 0.0; 632c4762a1bSJed Brown A[2] = 0.0; A[3] = 1.0; 6339566063dSJacob Faibussowitsch PetscCall(MatSetValues(jac,2,idx,2,idx,A,INSERT_VALUES)); 634c4762a1bSJed Brown 635c4762a1bSJed Brown /* 636c4762a1bSJed Brown Restore vector 637c4762a1bSJed Brown */ 6389566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x,&xx)); 639c4762a1bSJed Brown 640c4762a1bSJed Brown /* 641c4762a1bSJed Brown Assemble matrix 642c4762a1bSJed Brown */ 6439566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY)); 6449566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY)); 645c4762a1bSJed Brown PetscFunctionReturn(0); 646c4762a1bSJed Brown } 647c4762a1bSJed Brown 648c4762a1bSJed Brown int main(int argc,char **argv) 649c4762a1bSJed Brown { 650c4762a1bSJed Brown PetscMPIInt size; 651c4762a1bSJed Brown 652*327415f7SBarry Smith PetscFunctionBeginUser; 6539566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 6549566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 655be096a46SBarry Smith PetscCheck(size == 1,PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 6569566063dSJacob Faibussowitsch PetscCall(assembled_system()); 6579566063dSJacob Faibussowitsch PetscCall(block_system()); 6589566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 659b122ec5aSJacob Faibussowitsch return 0; 660c4762a1bSJed Brown } 661c4762a1bSJed Brown 662c4762a1bSJed Brown /*TEST 663c4762a1bSJed Brown 664c4762a1bSJed Brown test: 665c4762a1bSJed Brown args: -snes_monitor_short 666c4762a1bSJed Brown requires: !single 667c4762a1bSJed Brown 668c4762a1bSJed Brown TEST*/ 669