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