1 2 #include "src/ksp/pc/impls/is/nn/nn.h" 3 4 /* -------------------------------------------------------------------------- */ 5 /* 6 PCSetUp_NN - Prepares for the use of the NN preconditioner 7 by setting data structures and options. 8 9 Input Parameter: 10 . pc - the preconditioner context 11 12 Application Interface Routine: PCSetUp() 13 14 Notes: 15 The interface routine PCSetUp() is not usually called directly by 16 the user, but instead is called by PCApply() if necessary. 17 */ 18 #undef __FUNCT__ 19 #define __FUNCT__ "PCSetUp_NN" 20 static PetscErrorCode PCSetUp_NN(PC pc) 21 { 22 PetscErrorCode ierr; 23 24 PetscFunctionBegin; 25 if (!pc->setupcalled) { 26 /* Set up all the "iterative substructuring" common block */ 27 ierr = PCISSetUp(pc);CHKERRQ(ierr); 28 /* Create the coarse matrix. */ 29 ierr = PCNNCreateCoarseMatrix(pc);CHKERRQ(ierr); 30 } 31 PetscFunctionReturn(0); 32 } 33 34 /* -------------------------------------------------------------------------- */ 35 /* 36 PCApply_NN - Applies the NN preconditioner to a vector. 37 38 Input Parameters: 39 . pc - the preconditioner context 40 . r - input vector (global) 41 42 Output Parameter: 43 . z - output vector (global) 44 45 Application Interface Routine: PCApply() 46 */ 47 #undef __FUNCT__ 48 #define __FUNCT__ "PCApply_NN" 49 static PetscErrorCode PCApply_NN(PC pc,Vec r,Vec z) 50 { 51 PC_IS *pcis = (PC_IS*)(pc->data); 52 PetscErrorCode ierr; 53 PetscScalar m_one = -1.0; 54 Vec w = pcis->vec1_global; 55 56 PetscFunctionBegin; 57 58 /* 59 Dirichlet solvers. 60 Solving $ B_I^{(i)}r_I^{(i)} $ at each processor. 61 Storing the local results at vec2_D 62 */ 63 ierr = VecScatterBegin(r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD,pcis->global_to_D);CHKERRQ(ierr); 64 ierr = VecScatterEnd (r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD,pcis->global_to_D);CHKERRQ(ierr); 65 ierr = KSPSolve(pcis->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 66 67 /* 68 Computing $ r_B - \sum_j \tilde R_j^T A_{BI}^{(j)} (B_I^{(j)}r_I^{(j)}) $ . 69 Storing the result in the interface portion of the global vector w. 70 */ 71 ierr = MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);CHKERRQ(ierr); 72 ierr = VecScale(&m_one,pcis->vec1_B);CHKERRQ(ierr); 73 ierr = VecCopy(r,w);CHKERRQ(ierr); 74 ierr = VecScatterBegin(pcis->vec1_B,w,ADD_VALUES,SCATTER_REVERSE,pcis->global_to_B);CHKERRQ(ierr); 75 ierr = VecScatterEnd (pcis->vec1_B,w,ADD_VALUES,SCATTER_REVERSE,pcis->global_to_B);CHKERRQ(ierr); 76 77 /* 78 Apply the interface preconditioner 79 */ 80 ierr = PCNNApplyInterfacePreconditioner(pc,w,z,pcis->work_N,pcis->vec1_B,pcis->vec2_B,pcis->vec3_B,pcis->vec1_D, 81 pcis->vec3_D,pcis->vec1_N,pcis->vec2_N);CHKERRQ(ierr); 82 83 /* 84 Computing $ t_I^{(i)} = A_{IB}^{(i)} \tilde R_i z_B $ 85 The result is stored in vec1_D. 86 */ 87 ierr = VecScatterBegin(z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD,pcis->global_to_B);CHKERRQ(ierr); 88 ierr = VecScatterEnd (z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD,pcis->global_to_B);CHKERRQ(ierr); 89 ierr = MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);CHKERRQ(ierr); 90 91 /* 92 Dirichlet solvers. 93 Computing $ B_I^{(i)}t_I^{(i)} $ and sticking into the global vector the blocks 94 $ B_I^{(i)}r_I^{(i)} - B_I^{(i)}t_I^{(i)} $. 95 */ 96 ierr = VecScatterBegin(pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE,pcis->global_to_D);CHKERRQ(ierr); 97 ierr = VecScatterEnd (pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE,pcis->global_to_D);CHKERRQ(ierr); 98 ierr = KSPSolve(pcis->ksp_D,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 99 ierr = VecScale(&m_one,pcis->vec2_D);CHKERRQ(ierr); 100 ierr = VecScatterBegin(pcis->vec2_D,z,ADD_VALUES,SCATTER_REVERSE,pcis->global_to_D);CHKERRQ(ierr); 101 ierr = VecScatterEnd (pcis->vec2_D,z,ADD_VALUES,SCATTER_REVERSE,pcis->global_to_D);CHKERRQ(ierr); 102 103 PetscFunctionReturn(0); 104 } 105 106 /* -------------------------------------------------------------------------- */ 107 /* 108 PCDestroy_NN - Destroys the private context for the NN preconditioner 109 that was created with PCCreate_NN(). 110 111 Input Parameter: 112 . pc - the preconditioner context 113 114 Application Interface Routine: PCDestroy() 115 */ 116 #undef __FUNCT__ 117 #define __FUNCT__ "PCDestroy_NN" 118 static PetscErrorCode PCDestroy_NN(PC pc) 119 { 120 PC_NN *pcnn = (PC_NN*)pc->data; 121 PetscErrorCode ierr; 122 123 PetscFunctionBegin; 124 125 ierr = PCISDestroy(pc);CHKERRQ(ierr); 126 127 if (pcnn->coarse_mat) {ierr = MatDestroy(pcnn->coarse_mat);CHKERRQ(ierr);} 128 if (pcnn->coarse_x) {ierr = VecDestroy(pcnn->coarse_x);CHKERRQ(ierr);} 129 if (pcnn->coarse_b) {ierr = VecDestroy(pcnn->coarse_b);CHKERRQ(ierr);} 130 if (pcnn->ksp_coarse) {ierr = KSPDestroy(pcnn->ksp_coarse);CHKERRQ(ierr);} 131 if (pcnn->DZ_IN) { 132 if (pcnn->DZ_IN[0]) {ierr = PetscFree(pcnn->DZ_IN[0]);CHKERRQ(ierr);} 133 ierr = PetscFree(pcnn->DZ_IN);CHKERRQ(ierr); 134 } 135 136 /* 137 Free the private data structure that was hanging off the PC 138 */ 139 ierr = PetscFree(pcnn);CHKERRQ(ierr); 140 PetscFunctionReturn(0); 141 } 142 143 /* -------------------------------------------------------------------------- */ 144 /*MC 145 PCNN - Balancing Neumann-Neumann for scalar elliptic PDEs. 146 147 Options Database Keys: 148 + -pc_nn_turn_off_first_balancing - do not balance the residual before solving the local Neumann problems 149 (this skips the first coarse grid solve in the preconditioner) 150 . -pc_nn_turn_off_second_balancing - do not balance the solution solving the local Neumann problems 151 (this skips the second coarse grid solve in the preconditioner) 152 . -pc_is_damp_fixed <fact> - 153 . -pc_is_remove_nullspace_fixed - 154 . -pc_is_set_damping_factor_floating <fact> - 155 . -pc_is_not_damp_floating - 156 + -pc_is_not_remove_nullspace_floating - 157 158 Level: intermediate 159 160 Notes: The matrix used with this preconditioner must be of type MATIS 161 162 Unlike more 'conventional' Neumann-Neumann preconditioners this iterates over ALL the 163 degrees of freedom, NOT just those on the interface (this allows the use of approximate solvers 164 on the subdomains; though in our experience using approximate solvers is slower.). 165 166 Options for the coarse grid preconditioner can be set with -nn_coarse_pc_xxx 167 Options for the Dirichlet subproblem preconditioner can be set with -is_localD_pc_xxx 168 Options for the Neumann subproblem preconditioner can be set with -is_localN_pc_xxx 169 170 Contributed by Paulo Goldfeld 171 172 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, MatIS 173 M*/ 174 EXTERN_C_BEGIN 175 #undef __FUNCT__ 176 #define __FUNCT__ "PCCreate_NN" 177 PetscErrorCode PCCreate_NN(PC pc) 178 { 179 PetscErrorCode ierr; 180 PC_NN *pcnn; 181 182 PetscFunctionBegin; 183 184 /* 185 Creates the private data structure for this preconditioner and 186 attach it to the PC object. 187 */ 188 ierr = PetscNew(PC_NN,&pcnn);CHKERRQ(ierr); 189 pc->data = (void*)pcnn; 190 191 /* 192 Logs the memory usage; this is not needed but allows PETSc to 193 monitor how much memory is being used for various purposes. 194 */ 195 PetscLogObjectMemory(pc,sizeof(PC_NN)+sizeof(PC_IS)); /* Is this the right thing to do? */ 196 197 ierr = PCISCreate(pc);CHKERRQ(ierr); 198 pcnn->coarse_mat = 0; 199 pcnn->coarse_x = 0; 200 pcnn->coarse_b = 0; 201 pcnn->ksp_coarse = 0; 202 pcnn->DZ_IN = 0; 203 204 /* 205 Set the pointers for the functions that are provided above. 206 Now when the user-level routines (such as PCApply(), PCDestroy(), etc.) 207 are called, they will automatically call these functions. Note we 208 choose not to provide a couple of these functions since they are 209 not needed. 210 */ 211 pc->ops->apply = PCApply_NN; 212 pc->ops->applytranspose = 0; 213 pc->ops->setup = PCSetUp_NN; 214 pc->ops->destroy = PCDestroy_NN; 215 pc->ops->view = 0; 216 pc->ops->applyrichardson = 0; 217 pc->ops->applysymmetricleft = 0; 218 pc->ops->applysymmetricright = 0; 219 220 PetscFunctionReturn(0); 221 } 222 EXTERN_C_END 223 224 225 /* -------------------------------------------------------------------------- */ 226 /* 227 PCNNCreateCoarseMatrix - 228 */ 229 #undef __FUNCT__ 230 #define __FUNCT__ "PCNNCreateCoarseMatrix" 231 PetscErrorCode PCNNCreateCoarseMatrix (PC pc) 232 { 233 MPI_Request *send_request, *recv_request; 234 PetscErrorCode ierr; 235 int i, j, k; 236 237 PetscScalar* mat; /* Sub-matrix with this subdomain's contribution to the coarse matrix */ 238 PetscScalar** DZ_OUT; /* proc[k].DZ_OUT[i][] = bit of vector to be sent from processor k to processor i */ 239 240 /* aliasing some names */ 241 PC_IS* pcis = (PC_IS*)(pc->data); 242 PC_NN* pcnn = (PC_NN*)pc->data; 243 int n_neigh = pcis->n_neigh; 244 int* neigh = pcis->neigh; 245 int* n_shared = pcis->n_shared; 246 int** shared = pcis->shared; 247 PetscScalar** DZ_IN; /* Must be initialized after memory allocation. */ 248 249 PetscFunctionBegin; 250 251 /* Allocate memory for mat (the +1 is to handle the case n_neigh equal to zero) */ 252 ierr = PetscMalloc((n_neigh*n_neigh+1)*sizeof(PetscScalar),&mat);CHKERRQ(ierr); 253 254 /* Allocate memory for DZ */ 255 /* Notice that DZ_OUT[0] is allocated some space that is never used. */ 256 /* This is just in order to DZ_OUT and DZ_IN to have exactly the same form. */ 257 { 258 int size_of_Z = 0; 259 ierr = PetscMalloc ((n_neigh+1)*sizeof(PetscScalar*),&pcnn->DZ_IN);CHKERRQ(ierr); 260 DZ_IN = pcnn->DZ_IN; 261 ierr = PetscMalloc ((n_neigh+1)*sizeof(PetscScalar*),&DZ_OUT);CHKERRQ(ierr); 262 for (i=0; i<n_neigh; i++) { 263 size_of_Z += n_shared[i]; 264 } 265 ierr = PetscMalloc ((size_of_Z+1)*sizeof(PetscScalar),&DZ_IN[0]);CHKERRQ(ierr); 266 ierr = PetscMalloc ((size_of_Z+1)*sizeof(PetscScalar),&DZ_OUT[0]);CHKERRQ(ierr); 267 } 268 for (i=1; i<n_neigh; i++) { 269 DZ_IN[i] = DZ_IN [i-1] + n_shared[i-1]; 270 DZ_OUT[i] = DZ_OUT[i-1] + n_shared[i-1]; 271 } 272 273 /* Set the values of DZ_OUT, in order to send this info to the neighbours */ 274 /* First, set the auxiliary array pcis->work_N. */ 275 ierr = PCISScatterArrayNToVecB(pcis->work_N,pcis->D,INSERT_VALUES,SCATTER_REVERSE,pc);CHKERRQ(ierr); 276 for (i=1; i<n_neigh; i++){ 277 for (j=0; j<n_shared[i]; j++) { 278 DZ_OUT[i][j] = pcis->work_N[shared[i][j]]; 279 } 280 } 281 282 /* Non-blocking send/receive the common-interface chunks of scaled nullspaces */ 283 /* Notice that send_request[] and recv_request[] could have one less element. */ 284 /* We make them longer to have request[i] corresponding to neigh[i]. */ 285 { 286 int tag; 287 ierr = PetscObjectGetNewTag((PetscObject)pc,&tag);CHKERRQ(ierr); 288 ierr = PetscMalloc((2*(n_neigh)+1)*sizeof(MPI_Request),&send_request);CHKERRQ(ierr); 289 recv_request = send_request + (n_neigh); 290 for (i=1; i<n_neigh; i++) { 291 ierr = MPI_Isend((void*)(DZ_OUT[i]),n_shared[i],MPIU_SCALAR,neigh[i],tag,pc->comm,&(send_request[i]));CHKERRQ(ierr); 292 ierr = MPI_Irecv((void*)(DZ_IN [i]),n_shared[i],MPIU_SCALAR,neigh[i],tag,pc->comm,&(recv_request[i]));CHKERRQ(ierr); 293 } 294 } 295 296 /* Set DZ_IN[0][] (recall that neigh[0]==rank, always) */ 297 for(j=0; j<n_shared[0]; j++) { 298 DZ_IN[0][j] = pcis->work_N[shared[0][j]]; 299 } 300 301 /* Start computing with local D*Z while communication goes on. */ 302 /* Apply Schur complement. The result is "stored" in vec (more */ 303 /* precisely, vec points to the result, stored in pc_nn->vec1_B) */ 304 /* and also scattered to pcnn->work_N. */ 305 ierr = PCNNApplySchurToChunk(pc,n_shared[0],shared[0],DZ_IN[0],pcis->work_N,pcis->vec1_B, 306 pcis->vec2_B,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 307 308 /* Compute the first column, while completing the receiving. */ 309 for (i=0; i<n_neigh; i++) { 310 MPI_Status stat; 311 int ind=0; 312 if (i>0) { ierr = MPI_Waitany(n_neigh-1,recv_request+1,&ind,&stat);CHKERRQ(ierr); ind++;} 313 mat[ind*n_neigh+0] = 0.0; 314 for (k=0; k<n_shared[ind]; k++) { 315 mat[ind*n_neigh+0] += DZ_IN[ind][k] * pcis->work_N[shared[ind][k]]; 316 } 317 } 318 319 /* Compute the remaining of the columns */ 320 for (j=1; j<n_neigh; j++) { 321 ierr = PCNNApplySchurToChunk(pc,n_shared[j],shared[j],DZ_IN[j],pcis->work_N,pcis->vec1_B, 322 pcis->vec2_B,pcis->vec1_D,pcis->vec2_D);CHKERRQ(ierr); 323 for (i=0; i<n_neigh; i++) { 324 mat[i*n_neigh+j] = 0.0; 325 for (k=0; k<n_shared[i]; k++) { 326 mat[i*n_neigh+j] += DZ_IN[i][k] * pcis->work_N[shared[i][k]]; 327 } 328 } 329 } 330 331 /* Complete the sending. */ 332 if (n_neigh>1) { 333 MPI_Status *stat; 334 ierr = PetscMalloc((n_neigh-1)*sizeof(MPI_Status),&stat);CHKERRQ(ierr); 335 ierr = MPI_Waitall(n_neigh-1,&(send_request[1]),stat);CHKERRQ(ierr); 336 ierr = PetscFree(stat);CHKERRQ(ierr); 337 } 338 339 /* Free the memory for the MPI requests */ 340 ierr = PetscFree(send_request);CHKERRQ(ierr); 341 342 /* Free the memory for DZ_OUT */ 343 if (DZ_OUT) { 344 if (DZ_OUT[0]) { ierr = PetscFree(DZ_OUT[0]);CHKERRQ(ierr); } 345 ierr = PetscFree(DZ_OUT);CHKERRQ(ierr); 346 } 347 348 { 349 int size,n_neigh_m1; 350 ierr = MPI_Comm_size(pc->comm,&size);CHKERRQ(ierr); 351 n_neigh_m1 = (n_neigh) ? n_neigh-1 : 0; 352 /* Create the global coarse vectors (rhs and solution). */ 353 ierr = VecCreateMPI(pc->comm,1,size,&(pcnn->coarse_b));CHKERRQ(ierr); 354 ierr = VecDuplicate(pcnn->coarse_b,&(pcnn->coarse_x));CHKERRQ(ierr); 355 /* Create and set the global coarse AIJ matrix. */ 356 ierr = MatCreate(pc->comm,1,1,size,size,&(pcnn->coarse_mat));CHKERRQ(ierr); 357 ierr = MatSetType(pcnn->coarse_mat,MATAIJ);CHKERRQ(ierr); 358 ierr = MatSeqAIJSetPreallocation(pcnn->coarse_mat,1,PETSC_NULL);CHKERRQ(ierr); 359 ierr = MatMPIAIJSetPreallocation(pcnn->coarse_mat,1,PETSC_NULL,1,PETSC_NULL);CHKERRQ(ierr); 360 ierr = MatSetValues(pcnn->coarse_mat,n_neigh,neigh,n_neigh,neigh,mat,ADD_VALUES);CHKERRQ(ierr); 361 ierr = MatAssemblyBegin(pcnn->coarse_mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 362 ierr = MatAssemblyEnd (pcnn->coarse_mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 363 } 364 365 { 366 int rank; 367 PetscScalar one = 1.0; 368 IS is; 369 ierr = MPI_Comm_rank(pc->comm,&rank);CHKERRQ(ierr); 370 /* "Zero out" rows of not-purely-Neumann subdomains */ 371 if (pcis->pure_neumann) { /* does NOT zero the row; create an empty index set. The reason is that MatZeroRows() is collective. */ 372 ierr = ISCreateStride(pc->comm,0,0,0,&is);CHKERRQ(ierr); 373 } else { /* here it DOES zero the row, since it's not a floating subdomain. */ 374 ierr = ISCreateStride(pc->comm,1,rank,0,&is);CHKERRQ(ierr); 375 } 376 ierr = MatZeroRows(pcnn->coarse_mat,is,&one);CHKERRQ(ierr); 377 ierr = ISDestroy(is);CHKERRQ(ierr); 378 } 379 380 /* Create the coarse linear solver context */ 381 { 382 PC pc_ctx, inner_pc; 383 ierr = KSPCreate(pc->comm,&pcnn->ksp_coarse);CHKERRQ(ierr); 384 ierr = KSPSetOperators(pcnn->ksp_coarse,pcnn->coarse_mat,pcnn->coarse_mat,SAME_PRECONDITIONER);CHKERRQ(ierr); 385 ierr = KSPGetPC(pcnn->ksp_coarse,&pc_ctx);CHKERRQ(ierr); 386 ierr = PCSetType(pc_ctx,PCREDUNDANT);CHKERRQ(ierr); 387 ierr = KSPSetType(pcnn->ksp_coarse,KSPPREONLY);CHKERRQ(ierr); 388 ierr = PCRedundantGetPC(pc_ctx,&inner_pc);CHKERRQ(ierr); 389 ierr = PCSetType(inner_pc,PCLU);CHKERRQ(ierr); 390 ierr = KSPSetOptionsPrefix(pcnn->ksp_coarse,"nn_coarse_");CHKERRQ(ierr); 391 ierr = KSPSetFromOptions(pcnn->ksp_coarse);CHKERRQ(ierr); 392 /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */ 393 ierr = KSPSetUp(pcnn->ksp_coarse);CHKERRQ(ierr); 394 } 395 396 /* Free the memory for mat */ 397 ierr = PetscFree(mat);CHKERRQ(ierr); 398 399 /* for DEBUGGING, save the coarse matrix to a file. */ 400 { 401 PetscTruth flg; 402 ierr = PetscOptionsHasName(PETSC_NULL,"-pc_nn_save_coarse_matrix",&flg);CHKERRQ(ierr); 403 if (flg) { 404 PetscViewer viewer; 405 ierr = PetscViewerASCIIOpen(PETSC_COMM_WORLD,"coarse.m",&viewer);CHKERRQ(ierr); 406 ierr = PetscViewerSetFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr); 407 ierr = MatView(pcnn->coarse_mat,viewer);CHKERRQ(ierr); 408 ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr); 409 } 410 } 411 412 /* Set the variable pcnn->factor_coarse_rhs. */ 413 pcnn->factor_coarse_rhs = (pcis->pure_neumann) ? 1.0 : 0.0; 414 415 /* See historical note 02, at the bottom of this file. */ 416 417 PetscFunctionReturn(0); 418 } 419 420 /* -------------------------------------------------------------------------- */ 421 /* 422 PCNNApplySchurToChunk - 423 424 Input parameters: 425 . pcnn 426 . n - size of chunk 427 . idx - indices of chunk 428 . chunk - values 429 430 Output parameters: 431 . array_N - result of Schur complement applied to chunk, scattered to big array 432 . vec1_B - result of Schur complement applied to chunk 433 . vec2_B - garbage (used as work space) 434 . vec1_D - garbage (used as work space) 435 . vec2_D - garbage (used as work space) 436 437 */ 438 #undef __FUNCT__ 439 #define __FUNCT__ "PCNNApplySchurToChunk" 440 PetscErrorCode PCNNApplySchurToChunk(PC pc, int n, int* idx, PetscScalar *chunk, PetscScalar* array_N, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D) 441 { 442 PetscErrorCode ierr; 443 int i; 444 PC_IS *pcis = (PC_IS*)(pc->data); 445 446 PetscFunctionBegin; 447 448 ierr = PetscMemzero((void*)array_N, pcis->n*sizeof(PetscScalar));CHKERRQ(ierr); 449 for (i=0; i<n; i++) { array_N[idx[i]] = chunk[i]; } 450 ierr = PCISScatterArrayNToVecB(array_N,vec2_B,INSERT_VALUES,SCATTER_FORWARD,pc);CHKERRQ(ierr); 451 ierr = PCISApplySchur(pc,vec2_B,vec1_B,(Vec)0,vec1_D,vec2_D);CHKERRQ(ierr); 452 ierr = PCISScatterArrayNToVecB(array_N,vec1_B,INSERT_VALUES,SCATTER_REVERSE,pc);CHKERRQ(ierr); 453 454 PetscFunctionReturn(0); 455 } 456 457 /* -------------------------------------------------------------------------- */ 458 /* 459 PCNNApplyInterfacePreconditioner - Apply the interface preconditioner, i.e., 460 the preconditioner for the Schur complement. 461 462 Input parameter: 463 . r - global vector of interior and interface nodes. The values on the interior nodes are NOT used. 464 465 Output parameters: 466 . z - global vector of interior and interface nodes. The values on the interface are the result of 467 the application of the interface preconditioner to the interface part of r. The values on the 468 interior nodes are garbage. 469 . work_N - array of local nodes (interior and interface, including ghosts); returns garbage (used as work space) 470 . vec1_B - vector of local interface nodes (including ghosts); returns garbage (used as work space) 471 . vec2_B - vector of local interface nodes (including ghosts); returns garbage (used as work space) 472 . vec3_B - vector of local interface nodes (including ghosts); returns garbage (used as work space) 473 . vec1_D - vector of local interior nodes; returns garbage (used as work space) 474 . vec2_D - vector of local interior nodes; returns garbage (used as work space) 475 . vec1_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space) 476 . vec2_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space) 477 478 */ 479 #undef __FUNCT__ 480 #define __FUNCT__ "PCNNApplyInterfacePreconditioner" 481 PetscErrorCode PCNNApplyInterfacePreconditioner (PC pc, Vec r, Vec z, PetscScalar* work_N, Vec vec1_B, Vec vec2_B, Vec vec3_B, Vec vec1_D, 482 Vec vec2_D, Vec vec1_N, Vec vec2_N) 483 { 484 PetscErrorCode ierr; 485 PC_IS* pcis = (PC_IS*)(pc->data); 486 487 PetscFunctionBegin; 488 489 /* 490 First balancing step. 491 */ 492 { 493 PetscTruth flg; 494 ierr = PetscOptionsHasName(PETSC_NULL,"-pc_nn_turn_off_first_balancing",&flg);CHKERRQ(ierr); 495 if (!flg) { 496 ierr = PCNNBalancing(pc,r,(Vec)0,z,vec1_B,vec2_B,(Vec)0,vec1_D,vec2_D,work_N);CHKERRQ(ierr); 497 } else { 498 ierr = VecCopy(r,z);CHKERRQ(ierr); 499 } 500 } 501 502 /* 503 Extract the local interface part of z and scale it by D 504 */ 505 ierr = VecScatterBegin(z,vec1_B,INSERT_VALUES,SCATTER_FORWARD,pcis->global_to_B);CHKERRQ(ierr); 506 ierr = VecScatterEnd (z,vec1_B,INSERT_VALUES,SCATTER_FORWARD,pcis->global_to_B);CHKERRQ(ierr); 507 ierr = VecPointwiseMult(pcis->D,vec1_B,vec2_B);CHKERRQ(ierr); 508 509 /* Neumann Solver */ 510 ierr = PCISApplyInvSchur(pc,vec2_B,vec1_B,vec1_N,vec2_N);CHKERRQ(ierr); 511 512 /* 513 Second balancing step. 514 */ 515 { 516 PetscTruth flg; 517 ierr = PetscOptionsHasName(PETSC_NULL,"-pc_turn_off_second_balancing",&flg);CHKERRQ(ierr); 518 if (!flg) { 519 ierr = PCNNBalancing(pc,r,vec1_B,z,vec2_B,vec3_B,(Vec)0,vec1_D,vec2_D,work_N);CHKERRQ(ierr); 520 } else { 521 PetscScalar zero = 0.0; 522 ierr = VecPointwiseMult(pcis->D,vec1_B,vec2_B);CHKERRQ(ierr); 523 ierr = VecSet(&zero,z);CHKERRQ(ierr); 524 ierr = VecScatterBegin(vec2_B,z,ADD_VALUES,SCATTER_REVERSE,pcis->global_to_B);CHKERRQ(ierr); 525 ierr = VecScatterEnd (vec2_B,z,ADD_VALUES,SCATTER_REVERSE,pcis->global_to_B);CHKERRQ(ierr); 526 } 527 } 528 529 PetscFunctionReturn(0); 530 } 531 532 /* -------------------------------------------------------------------------- */ 533 /* 534 PCNNBalancing - Computes z, as given in equations (15) and (16) (if the 535 input argument u is provided), or s, as given in equations 536 (12) and (13), if the input argument u is a null vector. 537 Notice that the input argument u plays the role of u_i in 538 equation (14). The equation numbers refer to [Man93]. 539 540 Input Parameters: 541 . pcnn - NN preconditioner context. 542 . r - MPI vector of all nodes (interior and interface). It's preserved. 543 . u - (Optional) sequential vector of local interface nodes. It's preserved UNLESS vec3_B is null. 544 545 Output Parameters: 546 . z - MPI vector of interior and interface nodes. Returns s or z (see description above). 547 . vec1_B - Sequential vector of local interface nodes. Workspace. 548 . vec2_B - Sequential vector of local interface nodes. Workspace. 549 . vec3_B - (Optional) sequential vector of local interface nodes. Workspace. 550 . vec1_D - Sequential vector of local interior nodes. Workspace. 551 . vec2_D - Sequential vector of local interior nodes. Workspace. 552 . work_N - Array of all local nodes (interior and interface). Workspace. 553 554 */ 555 #undef __FUNCT__ 556 #define __FUNCT__ "PCNNBalancing" 557 PetscErrorCode PCNNBalancing (PC pc, Vec r, Vec u, Vec z, Vec vec1_B, Vec vec2_B, Vec vec3_B, 558 Vec vec1_D, Vec vec2_D, PetscScalar *work_N) 559 { 560 PetscErrorCode ierr; 561 int k; 562 PetscScalar zero = 0.0; 563 PetscScalar m_one = -1.0; 564 PetscScalar value; 565 PetscScalar* lambda; 566 PC_NN* pcnn = (PC_NN*)(pc->data); 567 PC_IS* pcis = (PC_IS*)(pc->data); 568 569 PetscFunctionBegin; 570 ierr = PetscLogEventBegin(PC_ApplyCoarse,0,0,0,0);CHKERRQ(ierr); 571 572 if (u) { 573 if (!vec3_B) { vec3_B = u; } 574 ierr = VecPointwiseMult(pcis->D,u,vec1_B);CHKERRQ(ierr); 575 ierr = VecSet(&zero,z);CHKERRQ(ierr); 576 ierr = VecScatterBegin(vec1_B,z,ADD_VALUES,SCATTER_REVERSE,pcis->global_to_B);CHKERRQ(ierr); 577 ierr = VecScatterEnd (vec1_B,z,ADD_VALUES,SCATTER_REVERSE,pcis->global_to_B);CHKERRQ(ierr); 578 ierr = VecScatterBegin(z,vec2_B,INSERT_VALUES,SCATTER_FORWARD,pcis->global_to_B);CHKERRQ(ierr); 579 ierr = VecScatterEnd (z,vec2_B,INSERT_VALUES,SCATTER_FORWARD,pcis->global_to_B);CHKERRQ(ierr); 580 ierr = PCISApplySchur(pc,vec2_B,vec3_B,(Vec)0,vec1_D,vec2_D);CHKERRQ(ierr); 581 ierr = VecScale(&m_one,vec3_B);CHKERRQ(ierr); 582 ierr = VecCopy(r,z);CHKERRQ(ierr); 583 ierr = VecScatterBegin(vec3_B,z,ADD_VALUES,SCATTER_REVERSE,pcis->global_to_B);CHKERRQ(ierr); 584 ierr = VecScatterEnd (vec3_B,z,ADD_VALUES,SCATTER_REVERSE,pcis->global_to_B);CHKERRQ(ierr); 585 } else { 586 ierr = VecCopy(r,z);CHKERRQ(ierr); 587 } 588 ierr = VecScatterBegin(z,vec2_B,INSERT_VALUES,SCATTER_FORWARD,pcis->global_to_B);CHKERRQ(ierr); 589 ierr = VecScatterEnd (z,vec2_B,INSERT_VALUES,SCATTER_FORWARD,pcis->global_to_B);CHKERRQ(ierr); 590 ierr = PCISScatterArrayNToVecB(work_N,vec2_B,INSERT_VALUES,SCATTER_REVERSE,pc);CHKERRQ(ierr); 591 for (k=0, value=0.0; k<pcis->n_shared[0]; k++) { value += pcnn->DZ_IN[0][k] * work_N[pcis->shared[0][k]]; } 592 value *= pcnn->factor_coarse_rhs; /* This factor is set in CreateCoarseMatrix(). */ 593 { 594 int rank; 595 ierr = MPI_Comm_rank(pc->comm,&rank);CHKERRQ(ierr); 596 ierr = VecSetValue(pcnn->coarse_b,rank,value,INSERT_VALUES);CHKERRQ(ierr); 597 /* 598 Since we are only inserting local values (one value actually) we don't need to do the 599 reduction that tells us there is no data that needs to be moved. Hence we comment out these 600 ierr = VecAssemblyBegin(pcnn->coarse_b);CHKERRQ(ierr); 601 ierr = VecAssemblyEnd (pcnn->coarse_b);CHKERRQ(ierr); 602 */ 603 } 604 ierr = KSPSolve(pcnn->ksp_coarse,pcnn->coarse_b,pcnn->coarse_x);CHKERRQ(ierr); 605 if (!u) { ierr = VecScale(&m_one,pcnn->coarse_x);CHKERRQ(ierr); } 606 ierr = VecGetArray(pcnn->coarse_x,&lambda);CHKERRQ(ierr); 607 for (k=0; k<pcis->n_shared[0]; k++) { work_N[pcis->shared[0][k]] = *lambda * pcnn->DZ_IN[0][k]; } 608 ierr = VecRestoreArray(pcnn->coarse_x,&lambda);CHKERRQ(ierr); 609 ierr = PCISScatterArrayNToVecB(work_N,vec2_B,INSERT_VALUES,SCATTER_FORWARD,pc);CHKERRQ(ierr); 610 ierr = VecSet(&zero,z);CHKERRQ(ierr); 611 ierr = VecScatterBegin(vec2_B,z,ADD_VALUES,SCATTER_REVERSE,pcis->global_to_B);CHKERRQ(ierr); 612 ierr = VecScatterEnd (vec2_B,z,ADD_VALUES,SCATTER_REVERSE,pcis->global_to_B);CHKERRQ(ierr); 613 if (!u) { 614 ierr = VecScatterBegin(z,vec2_B,INSERT_VALUES,SCATTER_FORWARD,pcis->global_to_B);CHKERRQ(ierr); 615 ierr = VecScatterEnd (z,vec2_B,INSERT_VALUES,SCATTER_FORWARD,pcis->global_to_B);CHKERRQ(ierr); 616 ierr = PCISApplySchur(pc,vec2_B,vec1_B,(Vec)0,vec1_D,vec2_D);CHKERRQ(ierr); 617 ierr = VecCopy(r,z);CHKERRQ(ierr); 618 } 619 ierr = VecScatterBegin(vec1_B,z,ADD_VALUES,SCATTER_REVERSE,pcis->global_to_B);CHKERRQ(ierr); 620 ierr = VecScatterEnd (vec1_B,z,ADD_VALUES,SCATTER_REVERSE,pcis->global_to_B);CHKERRQ(ierr); 621 ierr = PetscLogEventEnd(PC_ApplyCoarse,0,0,0,0);CHKERRQ(ierr); 622 623 PetscFunctionReturn(0); 624 } 625 626 #undef __FUNCT__ 627 628 629 630 /* ------- E N D O F T H E C O D E ------- */ 631 /* */ 632 /* From now on, "footnotes" (or "historical notes"). */ 633 /* */ 634 /* ------------------------------------------------- */ 635 636 637 638 /* -------------------------------------------------------------------------- 639 Historical note 01 640 -------------------------------------------------------------------------- */ 641 /* 642 We considered the possibility of an alternative D_i that would still 643 provide a partition of unity (i.e., $ \sum_i N_i D_i N_i^T = I $). 644 The basic principle was still the pseudo-inverse of the counting 645 function; the difference was that we would not count subdomains 646 that do not contribute to the coarse space (i.e., not pure-Neumann 647 subdomains). 648 649 This turned out to be a bad idea: we would solve trivial Neumann 650 problems in the not pure-Neumann subdomains, since we would be scaling 651 the balanced residual by zero. 652 */ 653 654 655 656 657 /* -------------------------------------------------------------------------- 658 Historical note 02 659 -------------------------------------------------------------------------- */ 660 /* 661 We tried an alternative coarse problem, that would eliminate exactly a 662 constant error. Turned out not to improve the overall convergence. 663 */ 664 665 666