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