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 52*9371c9d4SSatish Balay static PetscErrorCode assembled_system(void) { 53c4762a1bSJed Brown SNES snes; /* nonlinear solver context */ 54c4762a1bSJed Brown KSP ksp; /* linear solver context */ 55c4762a1bSJed Brown PC pc; /* preconditioner context */ 56c4762a1bSJed Brown Vec x, r; /* solution, residual vectors */ 57c4762a1bSJed Brown Mat J; /* Jacobian matrix */ 58c4762a1bSJed Brown PetscInt its; 59c4762a1bSJed Brown PetscScalar pfive = .5, *xx; 60c4762a1bSJed Brown PetscBool flg; 61c4762a1bSJed Brown 62c4762a1bSJed Brown PetscFunctionBeginUser; 639566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n\n========================= Assembled system =========================\n\n")); 64c4762a1bSJed Brown 65c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 66c4762a1bSJed Brown Create nonlinear solver context 67c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 68c4762a1bSJed Brown 699566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 70c4762a1bSJed Brown 71c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 72c4762a1bSJed Brown Create matrix and vector data structures; set corresponding routines 73c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 74c4762a1bSJed Brown 75c4762a1bSJed Brown /* 76c4762a1bSJed Brown Create vectors for solution and nonlinear function 77c4762a1bSJed Brown */ 789566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, 2, &x)); 799566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &r)); 80c4762a1bSJed Brown 81c4762a1bSJed Brown /* 82c4762a1bSJed Brown Create Jacobian matrix data structure 83c4762a1bSJed Brown */ 849566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &J)); 859566063dSJacob Faibussowitsch PetscCall(MatSetSizes(J, PETSC_DECIDE, PETSC_DECIDE, 2, 2)); 869566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(J)); 879566063dSJacob Faibussowitsch PetscCall(MatSetUp(J)); 88c4762a1bSJed Brown 899566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-hard", &flg)); 90c4762a1bSJed Brown if (!flg) { 91c4762a1bSJed Brown /* 92c4762a1bSJed Brown Set function evaluation routine and vector. 93c4762a1bSJed Brown */ 949566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes, r, FormFunction1, NULL)); 95c4762a1bSJed Brown 96c4762a1bSJed Brown /* 97c4762a1bSJed Brown Set Jacobian matrix data structure and Jacobian evaluation routine 98c4762a1bSJed Brown */ 999566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, J, J, FormJacobian1, NULL)); 100c4762a1bSJed Brown } else { 1019566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes, r, FormFunction2, NULL)); 1029566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, J, J, FormJacobian2, NULL)); 103c4762a1bSJed Brown } 104c4762a1bSJed Brown 105c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 106c4762a1bSJed Brown Customize nonlinear solver; set runtime options 107c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 108c4762a1bSJed Brown 109c4762a1bSJed Brown /* 110c4762a1bSJed Brown Set linear solver defaults for this problem. By extracting the 111c4762a1bSJed Brown KSP, KSP, and PC contexts from the SNES context, we can then 112c4762a1bSJed Brown directly call any KSP, KSP, and PC routines to set various options. 113c4762a1bSJed Brown */ 1149566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes, &ksp)); 1159566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &pc)); 1169566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCNONE)); 1179566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(ksp, 1.e-4, PETSC_DEFAULT, PETSC_DEFAULT, 20)); 118c4762a1bSJed Brown 119c4762a1bSJed Brown /* 120c4762a1bSJed Brown Set SNES/KSP/KSP/PC runtime options, e.g., 121c4762a1bSJed Brown -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc> 122c4762a1bSJed Brown These options will override those specified above as long as 123c4762a1bSJed Brown SNESSetFromOptions() is called _after_ any other customization 124c4762a1bSJed Brown routines. 125c4762a1bSJed Brown */ 1269566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 127c4762a1bSJed Brown 128c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 129c4762a1bSJed Brown Evaluate initial guess; then solve nonlinear system 130c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 131c4762a1bSJed Brown if (!flg) { 1329566063dSJacob Faibussowitsch PetscCall(VecSet(x, pfive)); 133c4762a1bSJed Brown } else { 1349566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &xx)); 135*9371c9d4SSatish Balay xx[0] = 2.0; 136*9371c9d4SSatish Balay 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 161*9371c9d4SSatish Balay PetscCall(VecDestroy(&x)); 162*9371c9d4SSatish Balay PetscCall(VecDestroy(&r)); 163*9371c9d4SSatish Balay PetscCall(MatDestroy(&J)); 164*9371c9d4SSatish Balay PetscCall(SNESDestroy(&snes)); 165c4762a1bSJed Brown PetscFunctionReturn(0); 166c4762a1bSJed Brown } 167c4762a1bSJed Brown 168c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 169c4762a1bSJed Brown /* 170c4762a1bSJed Brown FormFunction1 - Evaluates nonlinear function, F(x). 171c4762a1bSJed Brown 172c4762a1bSJed Brown Input Parameters: 173c4762a1bSJed Brown . snes - the SNES context 174c4762a1bSJed Brown . x - input vector 175c4762a1bSJed Brown . dummy - optional user-defined context (not used here) 176c4762a1bSJed Brown 177c4762a1bSJed Brown Output Parameter: 178c4762a1bSJed Brown . f - function vector 179c4762a1bSJed Brown */ 180*9371c9d4SSatish Balay static PetscErrorCode FormFunction1(SNES snes, Vec x, Vec f, void *dummy) { 181c4762a1bSJed Brown const PetscScalar *xx; 182c4762a1bSJed Brown PetscScalar *ff; 183c4762a1bSJed Brown 184c4762a1bSJed Brown PetscFunctionBeginUser; 185c4762a1bSJed Brown /* 186c4762a1bSJed Brown Get pointers to vector data. 187c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 188c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 189c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 190c4762a1bSJed Brown the array. 191c4762a1bSJed Brown */ 1929566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xx)); 1939566063dSJacob Faibussowitsch PetscCall(VecGetArray(f, &ff)); 194c4762a1bSJed Brown 195c4762a1bSJed Brown /* 196c4762a1bSJed Brown Compute function 197c4762a1bSJed Brown */ 198c4762a1bSJed Brown ff[0] = xx[0] * xx[0] + xx[0] * xx[1] - 3.0; 199c4762a1bSJed Brown ff[1] = xx[0] * xx[1] + xx[1] * xx[1] - 6.0; 200c4762a1bSJed Brown 201c4762a1bSJed Brown /* 202c4762a1bSJed Brown Restore vectors 203c4762a1bSJed Brown */ 2049566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xx)); 2059566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(f, &ff)); 206c4762a1bSJed Brown PetscFunctionReturn(0); 207c4762a1bSJed Brown } 208c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 209c4762a1bSJed Brown /* 210c4762a1bSJed Brown FormJacobian1 - Evaluates Jacobian matrix. 211c4762a1bSJed Brown 212c4762a1bSJed Brown Input Parameters: 213c4762a1bSJed Brown . snes - the SNES context 214c4762a1bSJed Brown . x - input vector 215c4762a1bSJed Brown . dummy - optional user-defined context (not used here) 216c4762a1bSJed Brown 217c4762a1bSJed Brown Output Parameters: 218c4762a1bSJed Brown . jac - Jacobian matrix 219c4762a1bSJed Brown . B - optionally different preconditioning matrix 220c4762a1bSJed Brown . flag - flag indicating matrix structure 221c4762a1bSJed Brown */ 222*9371c9d4SSatish Balay static PetscErrorCode FormJacobian1(SNES snes, Vec x, Mat jac, Mat B, void *dummy) { 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 */ 238*9371c9d4SSatish Balay A[0] = 2.0 * xx[0] + xx[1]; 239*9371c9d4SSatish Balay A[1] = xx[0]; 240*9371c9d4SSatish Balay A[2] = xx[1]; 241*9371c9d4SSatish Balay A[3] = xx[0] + 2.0 * xx[1]; 2429566063dSJacob Faibussowitsch PetscCall(MatSetValues(jac, 2, idx, 2, idx, A, INSERT_VALUES)); 243c4762a1bSJed Brown 244c4762a1bSJed Brown /* 245c4762a1bSJed Brown Restore vector 246c4762a1bSJed Brown */ 2479566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xx)); 248c4762a1bSJed Brown 249c4762a1bSJed Brown /* 250c4762a1bSJed Brown Assemble matrix 251c4762a1bSJed Brown */ 2529566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY)); 2539566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY)); 254c4762a1bSJed Brown PetscFunctionReturn(0); 255c4762a1bSJed Brown } 256c4762a1bSJed Brown 257c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 258*9371c9d4SSatish Balay static PetscErrorCode FormFunction2(SNES snes, Vec x, Vec f, void *dummy) { 259c4762a1bSJed Brown const PetscScalar *xx; 260c4762a1bSJed Brown PetscScalar *ff; 261c4762a1bSJed Brown 262c4762a1bSJed Brown PetscFunctionBeginUser; 263c4762a1bSJed Brown /* 264c4762a1bSJed Brown Get pointers to vector data. 265c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 266c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 267c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 268c4762a1bSJed Brown the array. 269c4762a1bSJed Brown */ 2709566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xx)); 2719566063dSJacob Faibussowitsch PetscCall(VecGetArray(f, &ff)); 272c4762a1bSJed Brown 273c4762a1bSJed Brown /* 274c4762a1bSJed Brown Compute function 275c4762a1bSJed Brown */ 276c4762a1bSJed Brown ff[0] = PetscSinScalar(3.0 * xx[0]) + xx[0]; 277c4762a1bSJed Brown ff[1] = xx[1]; 278c4762a1bSJed Brown 279c4762a1bSJed Brown /* 280c4762a1bSJed Brown Restore vectors 281c4762a1bSJed Brown */ 2829566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xx)); 2839566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(f, &ff)); 284c4762a1bSJed Brown PetscFunctionReturn(0); 285c4762a1bSJed Brown } 286c4762a1bSJed Brown 287c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 288*9371c9d4SSatish Balay static PetscErrorCode FormJacobian2(SNES snes, Vec x, Mat jac, Mat B, void *dummy) { 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 */ 304*9371c9d4SSatish Balay A[0] = 3.0 * PetscCosScalar(3.0 * xx[0]) + 1.0; 305*9371c9d4SSatish Balay A[1] = 0.0; 306*9371c9d4SSatish Balay A[2] = 0.0; 307*9371c9d4SSatish Balay A[3] = 1.0; 3089566063dSJacob Faibussowitsch PetscCall(MatSetValues(jac, 2, idx, 2, idx, A, INSERT_VALUES)); 309c4762a1bSJed Brown 310c4762a1bSJed Brown /* 311c4762a1bSJed Brown Restore vector 312c4762a1bSJed Brown */ 3139566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xx)); 314c4762a1bSJed Brown 315c4762a1bSJed Brown /* 316c4762a1bSJed Brown Assemble matrix 317c4762a1bSJed Brown */ 3189566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY)); 3199566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY)); 320c4762a1bSJed Brown PetscFunctionReturn(0); 321c4762a1bSJed Brown } 322c4762a1bSJed Brown 323*9371c9d4SSatish Balay static PetscErrorCode block_system(void) { 324c4762a1bSJed Brown SNES snes; /* nonlinear solver context */ 325c4762a1bSJed Brown KSP ksp; /* linear solver context */ 326c4762a1bSJed Brown PC pc; /* preconditioner context */ 327c4762a1bSJed Brown Vec x, r; /* solution, residual vectors */ 328c4762a1bSJed Brown Mat J; /* Jacobian matrix */ 329c4762a1bSJed Brown PetscInt its; 330c4762a1bSJed Brown PetscScalar pfive = .5; 331c4762a1bSJed Brown PetscBool flg; 332c4762a1bSJed Brown 333c4762a1bSJed Brown Mat j11, j12, j21, j22; 334c4762a1bSJed Brown Vec x1, x2, r1, r2; 335c4762a1bSJed Brown Vec bv; 336c4762a1bSJed Brown Vec bx[2]; 337c4762a1bSJed Brown Mat bA[2][2]; 338c4762a1bSJed Brown 339c4762a1bSJed Brown PetscFunctionBeginUser; 3409566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n\n========================= Block system =========================\n\n")); 341c4762a1bSJed Brown 3429566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 343c4762a1bSJed Brown 344c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 345c4762a1bSJed Brown Create matrix and vector data structures; set corresponding routines 346c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 347c4762a1bSJed Brown 348c4762a1bSJed Brown /* 349c4762a1bSJed Brown Create sub vectors for solution and nonlinear function 350c4762a1bSJed Brown */ 3519566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, 1, &x1)); 3529566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x1, &r1)); 353c4762a1bSJed Brown 3549566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, 1, &x2)); 3559566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x2, &r2)); 356c4762a1bSJed Brown 357c4762a1bSJed Brown /* 358c4762a1bSJed Brown Create the block vectors 359c4762a1bSJed Brown */ 360c4762a1bSJed Brown bx[0] = x1; 361c4762a1bSJed Brown bx[1] = x2; 3629566063dSJacob Faibussowitsch PetscCall(VecCreateNest(PETSC_COMM_WORLD, 2, NULL, bx, &x)); 3639566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(x)); 3649566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(x)); 3659566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x1)); 3669566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x2)); 367c4762a1bSJed Brown 368c4762a1bSJed Brown bx[0] = r1; 369c4762a1bSJed Brown bx[1] = r2; 3709566063dSJacob Faibussowitsch PetscCall(VecCreateNest(PETSC_COMM_WORLD, 2, NULL, bx, &r)); 3719566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r1)); 3729566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r2)); 3739566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(r)); 3749566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(r)); 375c4762a1bSJed Brown 376c4762a1bSJed Brown /* 377c4762a1bSJed Brown Create sub Jacobian matrix data structure 378c4762a1bSJed Brown */ 3799566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &j11)); 3809566063dSJacob Faibussowitsch PetscCall(MatSetSizes(j11, 1, 1, 1, 1)); 3819566063dSJacob Faibussowitsch PetscCall(MatSetType(j11, MATSEQAIJ)); 3829566063dSJacob Faibussowitsch PetscCall(MatSetUp(j11)); 383c4762a1bSJed Brown 3849566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &j12)); 3859566063dSJacob Faibussowitsch PetscCall(MatSetSizes(j12, 1, 1, 1, 1)); 3869566063dSJacob Faibussowitsch PetscCall(MatSetType(j12, MATSEQAIJ)); 3879566063dSJacob Faibussowitsch PetscCall(MatSetUp(j12)); 388c4762a1bSJed Brown 3899566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &j21)); 3909566063dSJacob Faibussowitsch PetscCall(MatSetSizes(j21, 1, 1, 1, 1)); 3919566063dSJacob Faibussowitsch PetscCall(MatSetType(j21, MATSEQAIJ)); 3929566063dSJacob Faibussowitsch PetscCall(MatSetUp(j21)); 393c4762a1bSJed Brown 3949566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &j22)); 3959566063dSJacob Faibussowitsch PetscCall(MatSetSizes(j22, PETSC_DECIDE, PETSC_DECIDE, 1, 1)); 3969566063dSJacob Faibussowitsch PetscCall(MatSetType(j22, MATSEQAIJ)); 3979566063dSJacob Faibussowitsch PetscCall(MatSetUp(j22)); 398c4762a1bSJed Brown /* 399c4762a1bSJed Brown Create block Jacobian matrix data structure 400c4762a1bSJed Brown */ 401c4762a1bSJed Brown bA[0][0] = j11; 402c4762a1bSJed Brown bA[0][1] = j12; 403c4762a1bSJed Brown bA[1][0] = j21; 404c4762a1bSJed Brown bA[1][1] = j22; 405c4762a1bSJed Brown 4069566063dSJacob Faibussowitsch PetscCall(MatCreateNest(PETSC_COMM_WORLD, 2, NULL, 2, NULL, &bA[0][0], &J)); 4079566063dSJacob Faibussowitsch PetscCall(MatSetUp(J)); 4089566063dSJacob Faibussowitsch PetscCall(MatNestSetVecType(J, VECNEST)); 4099566063dSJacob Faibussowitsch PetscCall(MatDestroy(&j11)); 4109566063dSJacob Faibussowitsch PetscCall(MatDestroy(&j12)); 4119566063dSJacob Faibussowitsch PetscCall(MatDestroy(&j21)); 4129566063dSJacob Faibussowitsch PetscCall(MatDestroy(&j22)); 413c4762a1bSJed Brown 4149566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 4159566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 416c4762a1bSJed Brown 4179566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-hard", &flg)); 418c4762a1bSJed Brown if (!flg) { 419c4762a1bSJed Brown /* 420c4762a1bSJed Brown Set function evaluation routine and vector. 421c4762a1bSJed Brown */ 4229566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes, r, FormFunction1_block, NULL)); 423c4762a1bSJed Brown 424c4762a1bSJed Brown /* 425c4762a1bSJed Brown Set Jacobian matrix data structure and Jacobian evaluation routine 426c4762a1bSJed Brown */ 4279566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, J, J, FormJacobian1_block, NULL)); 428c4762a1bSJed Brown } else { 4299566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes, r, FormFunction2_block, NULL)); 4309566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, J, J, FormJacobian2_block, NULL)); 431c4762a1bSJed Brown } 432c4762a1bSJed Brown 433c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 434c4762a1bSJed Brown Customize nonlinear solver; set runtime options 435c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 436c4762a1bSJed Brown 437c4762a1bSJed Brown /* 438c4762a1bSJed Brown Set linear solver defaults for this problem. By extracting the 439c4762a1bSJed Brown KSP, KSP, and PC contexts from the SNES context, we can then 440c4762a1bSJed Brown directly call any KSP, KSP, and PC routines to set various options. 441c4762a1bSJed Brown */ 4429566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes, &ksp)); 4439566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &pc)); 4449566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCNONE)); 4459566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(ksp, 1.e-4, PETSC_DEFAULT, PETSC_DEFAULT, 20)); 446c4762a1bSJed Brown 447c4762a1bSJed Brown /* 448c4762a1bSJed Brown Set SNES/KSP/KSP/PC runtime options, e.g., 449c4762a1bSJed Brown -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc> 450c4762a1bSJed Brown These options will override those specified above as long as 451c4762a1bSJed Brown SNESSetFromOptions() is called _after_ any other customization 452c4762a1bSJed Brown routines. 453c4762a1bSJed Brown */ 4549566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 455c4762a1bSJed Brown 456c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 457c4762a1bSJed Brown Evaluate initial guess; then solve nonlinear system 458c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 459c4762a1bSJed Brown if (!flg) { 4609566063dSJacob Faibussowitsch PetscCall(VecSet(x, pfive)); 461c4762a1bSJed Brown } else { 462c4762a1bSJed Brown Vec *vecs; 4639566063dSJacob Faibussowitsch PetscCall(VecNestGetSubVecs(x, NULL, &vecs)); 464c4762a1bSJed Brown bv = vecs[0]; 4659566063dSJacob Faibussowitsch /* PetscCall(VecBlockGetSubVec(x, 0, &bv)); */ 4669566063dSJacob Faibussowitsch PetscCall(VecSetValue(bv, 0, 2.0, INSERT_VALUES)); /* xx[0] = 2.0; */ 4679566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(bv)); 4689566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(bv)); 469c4762a1bSJed Brown 4709566063dSJacob Faibussowitsch /* PetscCall(VecBlockGetSubVec(x, 1, &bv)); */ 471c4762a1bSJed Brown bv = vecs[1]; 4729566063dSJacob Faibussowitsch PetscCall(VecSetValue(bv, 0, 3.0, INSERT_VALUES)); /* xx[1] = 3.0; */ 4739566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(bv)); 4749566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(bv)); 475c4762a1bSJed Brown } 476c4762a1bSJed Brown /* 477c4762a1bSJed Brown Note: The user should initialize the vector, x, with the initial guess 478c4762a1bSJed Brown for the nonlinear solver prior to calling SNESSolve(). In particular, 479c4762a1bSJed Brown to employ an initial guess of zero, the user should explicitly set 480c4762a1bSJed Brown this vector to zero by calling VecSet(). 481c4762a1bSJed Brown */ 4829566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, x)); 4839566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes, &its)); 484c4762a1bSJed Brown if (flg) { 485c4762a1bSJed Brown Vec f; 4869566063dSJacob Faibussowitsch PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD)); 4879566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(snes, &f, 0, 0)); 4889566063dSJacob Faibussowitsch PetscCall(VecView(r, PETSC_VIEWER_STDOUT_WORLD)); 489c4762a1bSJed Brown } 490c4762a1bSJed Brown 49163a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "number of SNES iterations = %" PetscInt_FMT "\n\n", its)); 492c4762a1bSJed Brown 493c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 494c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 495c4762a1bSJed Brown are no longer needed. 496c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 497*9371c9d4SSatish Balay PetscCall(VecDestroy(&x)); 498*9371c9d4SSatish Balay PetscCall(VecDestroy(&r)); 499*9371c9d4SSatish Balay PetscCall(MatDestroy(&J)); 500*9371c9d4SSatish Balay PetscCall(SNESDestroy(&snes)); 501c4762a1bSJed Brown PetscFunctionReturn(0); 502c4762a1bSJed Brown } 503c4762a1bSJed Brown 504c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 505*9371c9d4SSatish Balay static PetscErrorCode FormFunction1_block(SNES snes, Vec x, Vec f, void *dummy) { 506c4762a1bSJed Brown Vec *xx, *ff, x1, x2, f1, f2; 507c4762a1bSJed Brown PetscScalar ff_0, ff_1; 508c4762a1bSJed Brown PetscScalar xx_0, xx_1; 509c4762a1bSJed Brown PetscInt index, nb; 510c4762a1bSJed Brown 511c4762a1bSJed Brown PetscFunctionBeginUser; 512c4762a1bSJed Brown /* get blocks for function */ 5139566063dSJacob Faibussowitsch PetscCall(VecNestGetSubVecs(f, &nb, &ff)); 514*9371c9d4SSatish Balay f1 = ff[0]; 515*9371c9d4SSatish Balay f2 = ff[1]; 516c4762a1bSJed Brown 517c4762a1bSJed Brown /* get blocks for solution */ 5189566063dSJacob Faibussowitsch PetscCall(VecNestGetSubVecs(x, &nb, &xx)); 519*9371c9d4SSatish Balay x1 = xx[0]; 520*9371c9d4SSatish Balay x2 = xx[1]; 521c4762a1bSJed Brown 522c4762a1bSJed Brown /* get solution values */ 523c4762a1bSJed Brown index = 0; 5249566063dSJacob Faibussowitsch PetscCall(VecGetValues(x1, 1, &index, &xx_0)); 5259566063dSJacob Faibussowitsch PetscCall(VecGetValues(x2, 1, &index, &xx_1)); 526c4762a1bSJed Brown 527c4762a1bSJed Brown /* Compute function */ 528c4762a1bSJed Brown ff_0 = xx_0 * xx_0 + xx_0 * xx_1 - 3.0; 529c4762a1bSJed Brown ff_1 = xx_0 * xx_1 + xx_1 * xx_1 - 6.0; 530c4762a1bSJed Brown 531c4762a1bSJed Brown /* set function values */ 5329566063dSJacob Faibussowitsch PetscCall(VecSetValue(f1, index, ff_0, INSERT_VALUES)); 533c4762a1bSJed Brown 5349566063dSJacob Faibussowitsch PetscCall(VecSetValue(f2, index, ff_1, INSERT_VALUES)); 535c4762a1bSJed Brown 5369566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(f)); 5379566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(f)); 538c4762a1bSJed Brown PetscFunctionReturn(0); 539c4762a1bSJed Brown } 540c4762a1bSJed Brown 541c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 542*9371c9d4SSatish Balay static PetscErrorCode FormJacobian1_block(SNES snes, Vec x, Mat jac, Mat B, void *dummy) { 543c4762a1bSJed Brown Vec *xx, x1, x2; 544c4762a1bSJed Brown PetscScalar xx_0, xx_1; 545c4762a1bSJed Brown PetscInt index, nb; 546c4762a1bSJed Brown PetscScalar A_00, A_01, A_10, A_11; 547c4762a1bSJed Brown Mat j11, j12, j21, j22; 548c4762a1bSJed Brown Mat **mats; 549c4762a1bSJed Brown 550c4762a1bSJed Brown PetscFunctionBeginUser; 551c4762a1bSJed Brown /* get blocks for solution */ 5529566063dSJacob Faibussowitsch PetscCall(VecNestGetSubVecs(x, &nb, &xx)); 553*9371c9d4SSatish Balay x1 = xx[0]; 554*9371c9d4SSatish Balay x2 = xx[1]; 555c4762a1bSJed Brown 556c4762a1bSJed Brown /* get solution values */ 557c4762a1bSJed Brown index = 0; 5589566063dSJacob Faibussowitsch PetscCall(VecGetValues(x1, 1, &index, &xx_0)); 5599566063dSJacob Faibussowitsch PetscCall(VecGetValues(x2, 1, &index, &xx_1)); 560c4762a1bSJed Brown 561c4762a1bSJed Brown /* get block matrices */ 5629566063dSJacob Faibussowitsch PetscCall(MatNestGetSubMats(jac, NULL, NULL, &mats)); 563c4762a1bSJed Brown j11 = mats[0][0]; 564c4762a1bSJed Brown j12 = mats[0][1]; 565c4762a1bSJed Brown j21 = mats[1][0]; 566c4762a1bSJed Brown j22 = mats[1][1]; 567c4762a1bSJed Brown 568c4762a1bSJed Brown /* compute jacobian entries */ 569c4762a1bSJed Brown A_00 = 2.0 * xx_0 + xx_1; 570c4762a1bSJed Brown A_01 = xx_0; 571c4762a1bSJed Brown A_10 = xx_1; 572c4762a1bSJed Brown A_11 = xx_0 + 2.0 * xx_1; 573c4762a1bSJed Brown 574c4762a1bSJed Brown /* set jacobian values */ 5759566063dSJacob Faibussowitsch PetscCall(MatSetValue(j11, 0, 0, A_00, INSERT_VALUES)); 5769566063dSJacob Faibussowitsch PetscCall(MatSetValue(j12, 0, 0, A_01, INSERT_VALUES)); 5779566063dSJacob Faibussowitsch PetscCall(MatSetValue(j21, 0, 0, A_10, INSERT_VALUES)); 5789566063dSJacob Faibussowitsch PetscCall(MatSetValue(j22, 0, 0, A_11, INSERT_VALUES)); 579c4762a1bSJed Brown 580c4762a1bSJed Brown /* Assemble sub matrix */ 5819566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY)); 5829566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY)); 583c4762a1bSJed Brown PetscFunctionReturn(0); 584c4762a1bSJed Brown } 585c4762a1bSJed Brown 586c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 587*9371c9d4SSatish Balay static PetscErrorCode FormFunction2_block(SNES snes, Vec x, Vec f, void *dummy) { 588c4762a1bSJed Brown PetscScalar *ff; 589c4762a1bSJed Brown const PetscScalar *xx; 590c4762a1bSJed Brown 591c4762a1bSJed Brown PetscFunctionBeginUser; 592c4762a1bSJed Brown /* 593c4762a1bSJed Brown Get pointers to vector data. 594c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 595c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 596c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 597c4762a1bSJed Brown the array. 598c4762a1bSJed Brown */ 5999566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xx)); 6009566063dSJacob Faibussowitsch PetscCall(VecGetArray(f, &ff)); 601c4762a1bSJed Brown 602c4762a1bSJed Brown /* 603c4762a1bSJed Brown Compute function 604c4762a1bSJed Brown */ 605c4762a1bSJed Brown ff[0] = PetscSinScalar(3.0 * xx[0]) + xx[0]; 606c4762a1bSJed Brown ff[1] = xx[1]; 607c4762a1bSJed Brown 608c4762a1bSJed Brown /* 609c4762a1bSJed Brown Restore vectors 610c4762a1bSJed Brown */ 6119566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xx)); 6129566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(f, &ff)); 613c4762a1bSJed Brown PetscFunctionReturn(0); 614c4762a1bSJed Brown } 615c4762a1bSJed Brown 616c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 617*9371c9d4SSatish Balay static PetscErrorCode FormJacobian2_block(SNES snes, Vec x, Mat jac, Mat B, void *dummy) { 618c4762a1bSJed Brown const PetscScalar *xx; 619c4762a1bSJed Brown PetscScalar A[4]; 620c4762a1bSJed Brown PetscInt idx[2] = {0, 1}; 621c4762a1bSJed Brown 622c4762a1bSJed Brown PetscFunctionBeginUser; 623c4762a1bSJed Brown /* 624c4762a1bSJed Brown Get pointer to vector data 625c4762a1bSJed Brown */ 6269566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x, &xx)); 627c4762a1bSJed Brown 628c4762a1bSJed Brown /* 629c4762a1bSJed Brown Compute Jacobian entries and insert into matrix. 630c4762a1bSJed Brown - Since this is such a small problem, we set all entries for 631c4762a1bSJed Brown the matrix at once. 632c4762a1bSJed Brown */ 633*9371c9d4SSatish Balay A[0] = 3.0 * PetscCosScalar(3.0 * xx[0]) + 1.0; 634*9371c9d4SSatish Balay A[1] = 0.0; 635*9371c9d4SSatish Balay A[2] = 0.0; 636*9371c9d4SSatish Balay A[3] = 1.0; 6379566063dSJacob Faibussowitsch PetscCall(MatSetValues(jac, 2, idx, 2, idx, A, INSERT_VALUES)); 638c4762a1bSJed Brown 639c4762a1bSJed Brown /* 640c4762a1bSJed Brown Restore vector 641c4762a1bSJed Brown */ 6429566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x, &xx)); 643c4762a1bSJed Brown 644c4762a1bSJed Brown /* 645c4762a1bSJed Brown Assemble matrix 646c4762a1bSJed Brown */ 6479566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY)); 6489566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY)); 649c4762a1bSJed Brown PetscFunctionReturn(0); 650c4762a1bSJed Brown } 651c4762a1bSJed Brown 652*9371c9d4SSatish Balay int main(int argc, char **argv) { 653c4762a1bSJed Brown PetscMPIInt size; 654c4762a1bSJed Brown 655327415f7SBarry Smith PetscFunctionBeginUser; 6569566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 6579566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 658be096a46SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!"); 6599566063dSJacob Faibussowitsch PetscCall(assembled_system()); 6609566063dSJacob Faibussowitsch PetscCall(block_system()); 6619566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 662b122ec5aSJacob Faibussowitsch return 0; 663c4762a1bSJed Brown } 664c4762a1bSJed Brown 665c4762a1bSJed Brown /*TEST 666c4762a1bSJed Brown 667c4762a1bSJed Brown test: 668c4762a1bSJed Brown args: -snes_monitor_short 669c4762a1bSJed Brown requires: !single 670c4762a1bSJed Brown 671c4762a1bSJed Brown TEST*/ 672