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