1 #define PETSCKSP_DLL 2 3 /* 4 This file defines an additive Schwarz preconditioner for any Mat implementation. 5 6 Note that each processor may have any number of subdomains. But in order to 7 deal easily with the VecScatter(), we treat each processor as if it has the 8 same number of subdomains. 9 10 n - total number of true subdomains on all processors 11 n_local_true - actual number of subdomains on this processor 12 n_local = maximum over all processors of n_local_true 13 */ 14 #include "private/pcimpl.h" /*I "petscpc.h" I*/ 15 16 typedef struct { 17 PetscInt n, n_local, n_local_true; 18 PetscInt overlap; /* overlap requested by user */ 19 KSP *ksp; /* linear solvers for each block */ 20 VecScatter *restriction; /* mapping from global to subregion */ 21 VecScatter *prolongation; /* mapping from subregion to global */ 22 Vec *x,*y; /* work vectors */ 23 IS *is; /* index set that defines each overlapping subdomain */ 24 IS *is_local; /* index set that defines each non-overlapping subdomain, may be NULL */ 25 Mat *mat,*pmat; /* mat is not currently used */ 26 PCASMType type; /* use reduced interpolation, restriction or both */ 27 PetscTruth type_set; /* if user set this value (so won't change it for symmetric problems) */ 28 PetscTruth same_local_solves; /* flag indicating whether all local solvers are same */ 29 } PC_ASM; 30 31 #undef __FUNCT__ 32 #define __FUNCT__ "PCView_ASM" 33 static PetscErrorCode PCView_ASM(PC pc,PetscViewer viewer) 34 { 35 PC_ASM *osm = (PC_ASM*)pc->data; 36 PetscErrorCode ierr; 37 PetscMPIInt rank; 38 PetscInt i,bsz; 39 PetscTruth iascii,isstring; 40 PetscViewer sviewer; 41 42 43 PetscFunctionBegin; 44 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 45 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 46 if (iascii) { 47 if (osm->n > 0) { 48 ierr = PetscViewerASCIIPrintf(viewer," Additive Schwarz: total subdomain blocks = %D, amount of overlap = %D\n",osm->n,osm->overlap);CHKERRQ(ierr); 49 } else { 50 ierr = PetscViewerASCIIPrintf(viewer," Additive Schwarz: total subdomain blocks not yet set, amount of overlap = %D\n",osm->overlap);CHKERRQ(ierr); 51 } 52 ierr = PetscViewerASCIIPrintf(viewer," Additive Schwarz: restriction/interpolation type - %s\n",PCASMTypes[osm->type]);CHKERRQ(ierr); 53 ierr = MPI_Comm_rank(((PetscObject)pc)->comm,&rank);CHKERRQ(ierr); 54 if (osm->same_local_solves) { 55 if (osm->ksp) { 56 ierr = PetscViewerASCIIPrintf(viewer," Local solve is same for all blocks, in the following KSP and PC objects:\n");CHKERRQ(ierr); 57 ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 58 if (!rank) { 59 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 60 ierr = KSPView(osm->ksp[0],sviewer);CHKERRQ(ierr); 61 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 62 } 63 ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 64 } 65 } else { 66 ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] number of local blocks = %D\n",(int)rank,osm->n_local_true);CHKERRQ(ierr); 67 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 68 ierr = PetscViewerASCIIPrintf(viewer," Local solve info for each block is in the following KSP and PC objects:\n");CHKERRQ(ierr); 69 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 70 ierr = PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr); 71 for (i=0; i<osm->n_local; i++) { 72 ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 73 if (i < osm->n_local_true) { 74 ierr = ISGetLocalSize(osm->is[i],&bsz);CHKERRQ(ierr); 75 ierr = PetscViewerASCIISynchronizedPrintf(sviewer,"[%d] local block number %D, size = %D\n",(int)rank,i,bsz);CHKERRQ(ierr); 76 ierr = KSPView(osm->ksp[i],sviewer);CHKERRQ(ierr); 77 ierr = PetscViewerASCIISynchronizedPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr); 78 } 79 ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 80 } 81 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 82 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 83 } 84 } else if (isstring) { 85 ierr = PetscViewerStringSPrintf(viewer," blocks=%D, overlap=%D, type=%s",osm->n,osm->overlap,PCASMTypes[osm->type]);CHKERRQ(ierr); 86 ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 87 if (osm->ksp) {ierr = KSPView(osm->ksp[0],sviewer);CHKERRQ(ierr);} 88 ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 89 } else { 90 SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCASM",((PetscObject)viewer)->type_name); 91 } 92 PetscFunctionReturn(0); 93 } 94 95 #undef __FUNCT__ 96 #define __FUNCT__ "PCASMPrintSubdomains" 97 static PetscErrorCode PCASMPrintSubdomains(PC pc) 98 { 99 PC_ASM *osm = (PC_ASM*)pc->data; 100 const char *prefix; 101 char fname[PETSC_MAX_PATH_LEN+1]; 102 PetscViewer viewer; 103 PetscInt i,j,nidx; 104 const PetscInt *idx; 105 PetscErrorCode ierr; 106 107 PetscFunctionBegin; 108 ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 109 ierr = PetscOptionsGetString(prefix,"-pc_asm_print_subdomains",fname,PETSC_MAX_PATH_LEN,PETSC_NULL);CHKERRQ(ierr); 110 if (fname[0] == 0) { ierr = PetscStrcpy(fname,"stdout");CHKERRQ(ierr); }; 111 ierr = PetscViewerASCIIOpen(((PetscObject)pc)->comm,fname,&viewer);CHKERRQ(ierr); 112 for (i=0;i<osm->n_local_true;i++) { 113 ierr = ISGetLocalSize(osm->is[i],&nidx);CHKERRQ(ierr); 114 ierr = ISGetIndices(osm->is[i],&idx);CHKERRQ(ierr); 115 for (j=0; j<nidx; j++) { 116 ierr = PetscViewerASCIISynchronizedPrintf(viewer,"%D ",idx[j]);CHKERRQ(ierr); 117 } 118 ierr = ISRestoreIndices(osm->is[i],&idx);CHKERRQ(ierr); 119 ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr); 120 if (osm->is_local) { 121 ierr = ISGetLocalSize(osm->is_local[i],&nidx);CHKERRQ(ierr); 122 ierr = ISGetIndices(osm->is_local[i],&idx);CHKERRQ(ierr); 123 for (j=0; j<nidx; j++) { 124 ierr = PetscViewerASCIISynchronizedPrintf(viewer,"%D ",idx[j]);CHKERRQ(ierr); 125 } 126 ierr = ISRestoreIndices(osm->is_local[i],&idx);CHKERRQ(ierr); 127 ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr); 128 } 129 } 130 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 131 ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr); 132 PetscFunctionReturn(0); 133 } 134 135 #undef __FUNCT__ 136 #define __FUNCT__ "PCSetUp_ASM" 137 static PetscErrorCode PCSetUp_ASM(PC pc) 138 { 139 PC_ASM *osm = (PC_ASM*)pc->data; 140 PetscErrorCode ierr; 141 PetscTruth symset,flg; 142 PetscInt i,m; 143 PetscMPIInt size; 144 MatReuse scall = MAT_REUSE_MATRIX; 145 IS isl; 146 KSP ksp; 147 PC subpc; 148 const char *prefix,*pprefix; 149 Vec vec; 150 151 PetscFunctionBegin; 152 if (!pc->setupcalled) { 153 154 if (!osm->type_set) { 155 ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr); 156 if (symset && flg) { osm->type = PC_ASM_BASIC; } 157 } 158 159 if (osm->n == PETSC_DECIDE && osm->n_local_true < 1) { 160 /* no subdomains given, use one per processor */ 161 osm->n_local = osm->n_local_true = 1; 162 ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr); 163 osm->n = size; 164 } else if (osm->n == PETSC_DECIDE) { 165 /* determine global number of subdomains */ 166 PetscInt inwork[2],outwork[2]; 167 inwork[0] = inwork[1] = osm->n_local_true; 168 ierr = MPI_Allreduce(inwork,outwork,1,MPIU_2INT,PetscMaxSum_Op,((PetscObject)pc)->comm);CHKERRQ(ierr); 169 osm->n_local = outwork[0]; 170 osm->n = outwork[1]; 171 } 172 173 if (!osm->is){ /* create the index sets */ 174 ierr = PCASMCreateSubdomains(pc->pmat,osm->n_local_true,&osm->is);CHKERRQ(ierr); 175 } 176 ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 177 flg = PETSC_FALSE; 178 ierr = PetscOptionsGetTruth(prefix,"-pc_asm_print_subdomains",&flg,PETSC_NULL);CHKERRQ(ierr); 179 if (flg) { ierr = PCASMPrintSubdomains(pc);CHKERRQ(ierr); } 180 181 /* Extend the "overlapping" regions by a number of steps */ 182 ierr = MatIncreaseOverlap(pc->pmat,osm->n_local_true,osm->is,osm->overlap);CHKERRQ(ierr); 183 for (i=0; i<osm->n_local_true; i++) { 184 ierr = ISSort(osm->is[i]);CHKERRQ(ierr); 185 if (osm->is_local) { 186 ierr = ISSort(osm->is_local[i]);CHKERRQ(ierr); 187 } 188 } 189 190 /* Create the local work vectors and scatter contexts */ 191 ierr = MatGetVecs(pc->pmat,&vec,0);CHKERRQ(ierr); 192 ierr = PetscMalloc(osm->n_local*sizeof(VecScatter *),&osm->restriction);CHKERRQ(ierr); 193 ierr = PetscMalloc(osm->n_local*sizeof(VecScatter *),&osm->prolongation);CHKERRQ(ierr); 194 ierr = PetscMalloc(osm->n_local*sizeof(Vec *),&osm->x);CHKERRQ(ierr); 195 ierr = PetscMalloc(osm->n_local*sizeof(Vec *),&osm->y);CHKERRQ(ierr); 196 for (i=0; i<osm->n_local_true; i++) { 197 ierr = ISGetLocalSize(osm->is[i],&m);CHKERRQ(ierr); 198 ierr = VecCreateSeq(PETSC_COMM_SELF,m,&osm->x[i]);CHKERRQ(ierr); 199 ierr = ISCreateStride(PETSC_COMM_SELF,m,0,1,&isl);CHKERRQ(ierr); 200 ierr = VecScatterCreate(vec,osm->is[i],osm->x[i],isl,&osm->restriction[i]);CHKERRQ(ierr); 201 if (osm->is_local) { 202 ierr = ISGetLocalSize(osm->is_local[i],&m);CHKERRQ(ierr); 203 ierr = VecCreateSeq(PETSC_COMM_SELF,m,&osm->y[i]);CHKERRQ(ierr); 204 ierr = ISCreateStride(PETSC_COMM_SELF,m,0,1,&isl);CHKERRQ(ierr); 205 ierr = VecScatterCreate(vec,osm->is_local[i],osm->y[i],isl,&osm->prolongation[i]);CHKERRQ(ierr); 206 } else { 207 ierr = VecDuplicate(osm->x[i],&osm->y[i]);CHKERRQ(ierr); 208 osm->prolongation[i] = osm->restriction[i]; 209 ierr = PetscObjectReference((PetscObject) osm->restriction[i]);CHKERRQ(ierr); 210 } 211 ierr = ISDestroy(isl);CHKERRQ(ierr); 212 } 213 for (i=osm->n_local_true; i<osm->n_local; i++) { 214 ierr = VecCreateSeq(PETSC_COMM_SELF,0,&osm->x[i]);CHKERRQ(ierr); 215 ierr = VecDuplicate(osm->x[i],&osm->y[i]);CHKERRQ(ierr); 216 ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&isl);CHKERRQ(ierr); 217 ierr = VecScatterCreate(vec,isl,osm->x[i],isl,&osm->restriction[i]);CHKERRQ(ierr); 218 osm->prolongation[i] = osm->restriction[i]; 219 ierr = PetscObjectReference((PetscObject) osm->restriction[i]);CHKERRQ(ierr); 220 ierr = ISDestroy(isl);CHKERRQ(ierr); 221 } 222 ierr = VecDestroy(vec);CHKERRQ(ierr); 223 224 /* Create the local solvers */ 225 ierr = PetscMalloc(osm->n_local_true*sizeof(KSP *),&osm->ksp);CHKERRQ(ierr); 226 for (i=0; i<osm->n_local_true; i++) { 227 ierr = KSPCreate(PETSC_COMM_SELF,&ksp);CHKERRQ(ierr); 228 ierr = PetscLogObjectParent(pc,ksp);CHKERRQ(ierr); 229 ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);CHKERRQ(ierr); 230 ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr); 231 ierr = KSPGetPC(ksp,&subpc);CHKERRQ(ierr); 232 ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 233 ierr = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr); 234 ierr = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr); 235 osm->ksp[i] = ksp; 236 } 237 scall = MAT_INITIAL_MATRIX; 238 239 } else { 240 /* 241 Destroy the blocks from the previous iteration 242 */ 243 if (pc->flag == DIFFERENT_NONZERO_PATTERN) { 244 ierr = MatDestroyMatrices(osm->n_local_true,&osm->pmat);CHKERRQ(ierr); 245 scall = MAT_INITIAL_MATRIX; 246 } 247 } 248 249 /* 250 Extract out the submatrices 251 */ 252 ierr = MatGetSubMatrices(pc->pmat,osm->n_local_true,osm->is,osm->is,scall,&osm->pmat);CHKERRQ(ierr); 253 if (scall == MAT_INITIAL_MATRIX) { 254 ierr = PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);CHKERRQ(ierr); 255 for (i=0; i<osm->n_local_true; i++) { 256 ierr = PetscLogObjectParent(pc,osm->pmat[i]);CHKERRQ(ierr); 257 ierr = PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix);CHKERRQ(ierr); 258 } 259 } 260 261 /* Return control to the user so that the submatrices can be modified (e.g., to apply 262 different boundary conditions for the submatrices than for the global problem) */ 263 ierr = PCModifySubMatrices(pc,osm->n_local_true,osm->is,osm->is,osm->pmat,pc->modifysubmatricesP);CHKERRQ(ierr); 264 265 /* 266 Loop over subdomains putting them into local ksp 267 */ 268 for (i=0; i<osm->n_local_true; i++) { 269 ierr = KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i],pc->flag);CHKERRQ(ierr); 270 if (!pc->setupcalled) { 271 ierr = KSPSetFromOptions(osm->ksp[i]);CHKERRQ(ierr); 272 } 273 } 274 275 PetscFunctionReturn(0); 276 } 277 278 #undef __FUNCT__ 279 #define __FUNCT__ "PCSetUpOnBlocks_ASM" 280 static PetscErrorCode PCSetUpOnBlocks_ASM(PC pc) 281 { 282 PC_ASM *osm = (PC_ASM*)pc->data; 283 PetscErrorCode ierr; 284 PetscInt i; 285 286 PetscFunctionBegin; 287 for (i=0; i<osm->n_local_true; i++) { 288 ierr = KSPSetUp(osm->ksp[i]);CHKERRQ(ierr); 289 } 290 PetscFunctionReturn(0); 291 } 292 293 #undef __FUNCT__ 294 #define __FUNCT__ "PCApply_ASM" 295 static PetscErrorCode PCApply_ASM(PC pc,Vec x,Vec y) 296 { 297 PC_ASM *osm = (PC_ASM*)pc->data; 298 PetscErrorCode ierr; 299 PetscInt i,n_local = osm->n_local,n_local_true = osm->n_local_true; 300 ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 301 302 PetscFunctionBegin; 303 /* 304 Support for limiting the restriction or interpolation to only local 305 subdomain values (leaving the other values 0). 306 */ 307 if (!(osm->type & PC_ASM_RESTRICT)) { 308 forward = SCATTER_FORWARD_LOCAL; 309 /* have to zero the work RHS since scatter may leave some slots empty */ 310 for (i=0; i<n_local_true; i++) { 311 ierr = VecSet(osm->x[i],0.0);CHKERRQ(ierr); 312 } 313 } 314 if (!(osm->type & PC_ASM_INTERPOLATE)) { 315 reverse = SCATTER_REVERSE_LOCAL; 316 } 317 318 for (i=0; i<n_local; i++) { 319 ierr = VecScatterBegin(osm->restriction[i],x,osm->x[i],INSERT_VALUES,forward);CHKERRQ(ierr); 320 } 321 ierr = VecSet(y,0.0);CHKERRQ(ierr); 322 /* do the local solves */ 323 for (i=0; i<n_local_true; i++) { 324 ierr = VecScatterEnd(osm->restriction[i],x,osm->x[i],INSERT_VALUES,forward);CHKERRQ(ierr); 325 ierr = KSPSolve(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr); 326 ierr = VecScatterBegin(osm->prolongation[i],osm->y[i],y,ADD_VALUES,reverse);CHKERRQ(ierr); 327 } 328 /* handle the rest of the scatters that do not have local solves */ 329 for (i=n_local_true; i<n_local; i++) { 330 ierr = VecScatterEnd(osm->restriction[i],x,osm->x[i],INSERT_VALUES,forward);CHKERRQ(ierr); 331 ierr = VecScatterBegin(osm->prolongation[i],osm->y[i],y,ADD_VALUES,reverse);CHKERRQ(ierr); 332 } 333 for (i=0; i<n_local; i++) { 334 ierr = VecScatterEnd(osm->prolongation[i],osm->y[i],y,ADD_VALUES,reverse);CHKERRQ(ierr); 335 } 336 PetscFunctionReturn(0); 337 } 338 339 #undef __FUNCT__ 340 #define __FUNCT__ "PCApplyTranspose_ASM" 341 static PetscErrorCode PCApplyTranspose_ASM(PC pc,Vec x,Vec y) 342 { 343 PC_ASM *osm = (PC_ASM*)pc->data; 344 PetscErrorCode ierr; 345 PetscInt i,n_local = osm->n_local,n_local_true = osm->n_local_true; 346 ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 347 348 PetscFunctionBegin; 349 /* 350 Support for limiting the restriction or interpolation to only local 351 subdomain values (leaving the other values 0). 352 353 Note: these are reversed from the PCApply_ASM() because we are applying the 354 transpose of the three terms 355 */ 356 if (!(osm->type & PC_ASM_INTERPOLATE)) { 357 forward = SCATTER_FORWARD_LOCAL; 358 /* have to zero the work RHS since scatter may leave some slots empty */ 359 for (i=0; i<n_local_true; i++) { 360 ierr = VecSet(osm->x[i],0.0);CHKERRQ(ierr); 361 } 362 } 363 if (!(osm->type & PC_ASM_RESTRICT)) { 364 reverse = SCATTER_REVERSE_LOCAL; 365 } 366 367 for (i=0; i<n_local; i++) { 368 ierr = VecScatterBegin(osm->restriction[i],x,osm->x[i],INSERT_VALUES,forward);CHKERRQ(ierr); 369 } 370 ierr = VecSet(y,0.0);CHKERRQ(ierr); 371 /* do the local solves */ 372 for (i=0; i<n_local_true; i++) { 373 ierr = VecScatterEnd(osm->restriction[i],x,osm->x[i],INSERT_VALUES,forward);CHKERRQ(ierr); 374 ierr = KSPSolveTranspose(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr); 375 ierr = VecScatterBegin(osm->prolongation[i],osm->y[i],y,ADD_VALUES,reverse);CHKERRQ(ierr); 376 } 377 /* handle the rest of the scatters that do not have local solves */ 378 for (i=n_local_true; i<n_local; i++) { 379 ierr = VecScatterEnd(osm->restriction[i],x,osm->x[i],INSERT_VALUES,forward);CHKERRQ(ierr); 380 ierr = VecScatterBegin(osm->prolongation[i],osm->y[i],y,ADD_VALUES,reverse);CHKERRQ(ierr); 381 } 382 for (i=0; i<n_local; i++) { 383 ierr = VecScatterEnd(osm->prolongation[i],osm->y[i],y,ADD_VALUES,reverse);CHKERRQ(ierr); 384 } 385 PetscFunctionReturn(0); 386 } 387 388 #undef __FUNCT__ 389 #define __FUNCT__ "PCDestroy_ASM" 390 static PetscErrorCode PCDestroy_ASM(PC pc) 391 { 392 PC_ASM *osm = (PC_ASM*)pc->data; 393 PetscErrorCode ierr; 394 PetscInt i; 395 396 PetscFunctionBegin; 397 if (osm->ksp) { 398 for (i=0; i<osm->n_local_true; i++) { 399 ierr = KSPDestroy(osm->ksp[i]);CHKERRQ(ierr); 400 } 401 ierr = PetscFree(osm->ksp);CHKERRQ(ierr); 402 } 403 if (osm->pmat) { 404 if (osm->n_local_true > 0) { 405 ierr = MatDestroyMatrices(osm->n_local_true,&osm->pmat);CHKERRQ(ierr); 406 } 407 } 408 if (osm->restriction) { 409 for (i=0; i<osm->n_local; i++) { 410 ierr = VecScatterDestroy(osm->restriction[i]);CHKERRQ(ierr); 411 ierr = VecScatterDestroy(osm->prolongation[i]);CHKERRQ(ierr); 412 ierr = VecDestroy(osm->x[i]);CHKERRQ(ierr); 413 ierr = VecDestroy(osm->y[i]);CHKERRQ(ierr); 414 } 415 ierr = PetscFree(osm->restriction);CHKERRQ(ierr); 416 ierr = PetscFree(osm->prolongation);CHKERRQ(ierr); 417 ierr = PetscFree(osm->x);CHKERRQ(ierr); 418 ierr = PetscFree(osm->y);CHKERRQ(ierr); 419 } 420 if (osm->is) { 421 ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr); 422 } 423 ierr = PetscFree(osm);CHKERRQ(ierr); 424 PetscFunctionReturn(0); 425 } 426 427 #undef __FUNCT__ 428 #define __FUNCT__ "PCSetFromOptions_ASM" 429 static PetscErrorCode PCSetFromOptions_ASM(PC pc) 430 { 431 PC_ASM *osm = (PC_ASM*)pc->data; 432 PetscErrorCode ierr; 433 PetscInt blocks,ovl; 434 PetscTruth symset,flg; 435 PCASMType asmtype; 436 437 PetscFunctionBegin; 438 /* set the type to symmetric if matrix is symmetric */ 439 if (!osm->type_set && pc->pmat) { 440 ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr); 441 if (symset && flg) { osm->type = PC_ASM_BASIC; } 442 } 443 ierr = PetscOptionsHead("Additive Schwarz options");CHKERRQ(ierr); 444 ierr = PetscOptionsInt("-pc_asm_blocks","Number of subdomains","PCASMSetTotalSubdomains",osm->n,&blocks,&flg);CHKERRQ(ierr); 445 if (flg) {ierr = PCASMSetTotalSubdomains(pc,blocks,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); } 446 ierr = PetscOptionsInt("-pc_asm_overlap","Number of grid points overlap","PCASMSetOverlap",osm->overlap,&ovl,&flg);CHKERRQ(ierr); 447 if (flg) {ierr = PCASMSetOverlap(pc,ovl);CHKERRQ(ierr); } 448 flg = PETSC_FALSE; 449 ierr = PetscOptionsEnum("-pc_asm_type","Type of restriction/extension","PCASMSetType",PCASMTypes,(PetscEnum)osm->type,(PetscEnum*)&asmtype,&flg);CHKERRQ(ierr); 450 if (flg) {ierr = PCASMSetType(pc,asmtype);CHKERRQ(ierr); } 451 ierr = PetscOptionsTail();CHKERRQ(ierr); 452 PetscFunctionReturn(0); 453 } 454 455 /*------------------------------------------------------------------------------------*/ 456 457 EXTERN_C_BEGIN 458 #undef __FUNCT__ 459 #define __FUNCT__ "PCASMSetLocalSubdomains_ASM" 460 PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetLocalSubdomains_ASM(PC pc,PetscInt n,IS is[],IS is_local[]) 461 { 462 PC_ASM *osm = (PC_ASM*)pc->data; 463 PetscErrorCode ierr; 464 PetscInt i; 465 466 PetscFunctionBegin; 467 if (n < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more blocks, n = %D",n); 468 if (pc->setupcalled && (n != osm->n_local_true || is)) { 469 SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"PCASMSetLocalSubdomains() should be called before calling PCSetup()."); 470 } 471 if (!pc->setupcalled) { 472 if (is) { 473 for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is[i]);CHKERRQ(ierr);} 474 } 475 if (osm->is) { 476 ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr); 477 } 478 osm->n_local_true = n; 479 osm->is = 0; 480 osm->is_local = 0; 481 if (is) { 482 ierr = PetscMalloc(n*sizeof(IS *),&osm->is);CHKERRQ(ierr); 483 for (i=0; i<n; i++) { osm->is[i] = is[i]; } 484 } 485 if (is_local) { 486 ierr = PetscMalloc(n*sizeof(IS *),&osm->is_local);CHKERRQ(ierr); 487 for (i=0; i<n; i++) { osm->is_local[i] = is_local[i]; } 488 } 489 } 490 PetscFunctionReturn(0); 491 } 492 EXTERN_C_END 493 494 EXTERN_C_BEGIN 495 #undef __FUNCT__ 496 #define __FUNCT__ "PCASMSetTotalSubdomains_ASM" 497 PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetTotalSubdomains_ASM(PC pc,PetscInt N,IS *is) 498 { 499 PC_ASM *osm = (PC_ASM*)pc->data; 500 PetscErrorCode ierr; 501 PetscMPIInt rank,size; 502 PetscInt n; 503 504 PetscFunctionBegin; 505 if (N < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Number of total blocks must be > 0, N = %D",N); 506 if (is) SETERRQ(PETSC_ERR_SUP,"Use PCASMSetLocalSubdomains() to set specific index sets\n\they cannot be set globally yet."); 507 508 /* 509 Split the subdomains equally amoung all processors 510 */ 511 ierr = MPI_Comm_rank(((PetscObject)pc)->comm,&rank);CHKERRQ(ierr); 512 ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr); 513 n = N/size + ((N % size) > rank); 514 if (!n) SETERRQ3(PETSC_ERR_ARG_OUTOFRANGE,"Process %d must have at least one block: total processors %d total blocks %D",(int)rank,(int)size,N); 515 if (pc->setupcalled && n != osm->n_local_true) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"PCASMSetTotalSubdomains() should be called before PCSetup()."); 516 if (!pc->setupcalled) { 517 if (osm->is) { 518 ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr); 519 } 520 osm->n_local_true = n; 521 osm->is = 0; 522 osm->is_local = 0; 523 } 524 PetscFunctionReturn(0); 525 } 526 EXTERN_C_END 527 528 EXTERN_C_BEGIN 529 #undef __FUNCT__ 530 #define __FUNCT__ "PCASMSetOverlap_ASM" 531 PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetOverlap_ASM(PC pc,PetscInt ovl) 532 { 533 PC_ASM *osm = (PC_ASM*)pc->data; 534 535 PetscFunctionBegin; 536 if (ovl < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested"); 537 if (pc->setupcalled && ovl != osm->overlap) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"PCASMSetOverlap() should be called before PCSetup()."); 538 if (!pc->setupcalled) { 539 osm->overlap = ovl; 540 } 541 PetscFunctionReturn(0); 542 } 543 EXTERN_C_END 544 545 EXTERN_C_BEGIN 546 #undef __FUNCT__ 547 #define __FUNCT__ "PCASMSetType_ASM" 548 PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetType_ASM(PC pc,PCASMType type) 549 { 550 PC_ASM *osm = (PC_ASM*)pc->data; 551 552 PetscFunctionBegin; 553 osm->type = type; 554 osm->type_set = PETSC_TRUE; 555 PetscFunctionReturn(0); 556 } 557 EXTERN_C_END 558 559 EXTERN_C_BEGIN 560 #undef __FUNCT__ 561 #define __FUNCT__ "PCASMGetSubKSP_ASM" 562 PetscErrorCode PETSCKSP_DLLEXPORT PCASMGetSubKSP_ASM(PC pc,PetscInt *n_local,PetscInt *first_local,KSP **ksp) 563 { 564 PC_ASM *osm = (PC_ASM*)pc->data; 565 PetscErrorCode ierr; 566 567 PetscFunctionBegin; 568 if (osm->n_local_true < 1) { 569 SETERRQ(PETSC_ERR_ORDER,"Need to call PCSetUP() on PC (or KSPSetUp() on the outer KSP object) before calling here"); 570 } 571 572 if (n_local) { 573 *n_local = osm->n_local_true; 574 } 575 if (first_local) { 576 ierr = MPI_Scan(&osm->n_local_true,first_local,1,MPIU_INT,MPI_SUM,((PetscObject)pc)->comm);CHKERRQ(ierr); 577 *first_local -= osm->n_local_true; 578 } 579 if (ksp) { 580 /* Assume that local solves are now different; not necessarily 581 true though! This flag is used only for PCView_ASM() */ 582 *ksp = osm->ksp; 583 osm->same_local_solves = PETSC_FALSE; 584 } 585 PetscFunctionReturn(0); 586 } 587 EXTERN_C_END 588 589 590 #undef __FUNCT__ 591 #define __FUNCT__ "PCASMSetLocalSubdomains" 592 /*@C 593 PCASMSetLocalSubdomains - Sets the local subdomains (for this processor 594 only) for the additive Schwarz preconditioner. 595 596 Collective on PC 597 598 Input Parameters: 599 + pc - the preconditioner context 600 . n - the number of subdomains for this processor (default value = 1) 601 . is - the index sets that define the subdomains for this processor 602 (or PETSC_NULL for PETSc to determine subdomains) 603 - is_local - the index sets that define the local part of the subdomains for this processor 604 (or PETSC_NULL to use the default of 1 subdomain per process) 605 606 Notes: 607 The IS numbering is in the parallel, global numbering of the vector. 608 609 By default the ASM preconditioner uses 1 block per processor. 610 611 These index sets cannot be destroyed until after completion of the 612 linear solves for which the ASM preconditioner is being used. 613 614 Use PCASMSetTotalSubdomains() to set the subdomains for all processors. 615 616 Level: advanced 617 618 .keywords: PC, ASM, set, local, subdomains, additive Schwarz 619 620 .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(), 621 PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains() 622 @*/ 623 PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetLocalSubdomains(PC pc,PetscInt n,IS is[],IS is_local[]) 624 { 625 PetscErrorCode ierr,(*f)(PC,PetscInt,IS[],IS[]); 626 627 PetscFunctionBegin; 628 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 629 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCASMSetLocalSubdomains_C",(void (**)(void))&f);CHKERRQ(ierr); 630 if (f) { 631 ierr = (*f)(pc,n,is,is_local);CHKERRQ(ierr); 632 } 633 PetscFunctionReturn(0); 634 } 635 636 #undef __FUNCT__ 637 #define __FUNCT__ "PCASMSetTotalSubdomains" 638 /*@C 639 PCASMSetTotalSubdomains - Sets the subdomains for all processor for the 640 additive Schwarz preconditioner. Either all or no processors in the 641 PC communicator must call this routine, with the same index sets. 642 643 Collective on PC 644 645 Input Parameters: 646 + pc - the preconditioner context 647 . n - the number of subdomains for all processors 648 . is - the index sets that define the subdomains for all processor 649 (or PETSC_NULL for PETSc to determine subdomains) 650 - is_local - the index sets that define the local part of the subdomains for this processor 651 (or PETSC_NULL to use the default of 1 subdomain per process) 652 653 Options Database Key: 654 To set the total number of subdomain blocks rather than specify the 655 index sets, use the option 656 . -pc_asm_blocks <blks> - Sets total blocks 657 658 Notes: 659 Currently you cannot use this to set the actual subdomains with the argument is. 660 661 By default the ASM preconditioner uses 1 block per processor. 662 663 These index sets cannot be destroyed until after completion of the 664 linear solves for which the ASM preconditioner is being used. 665 666 Use PCASMSetLocalSubdomains() to set local subdomains. 667 668 Level: advanced 669 670 .keywords: PC, ASM, set, total, global, subdomains, additive Schwarz 671 672 .seealso: PCASMSetLocalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(), 673 PCASMCreateSubdomains2D() 674 @*/ 675 PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetTotalSubdomains(PC pc,PetscInt N,IS is[],IS is_local[]) 676 { 677 PetscErrorCode ierr,(*f)(PC,PetscInt,IS[],IS[]); 678 679 PetscFunctionBegin; 680 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 681 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCASMSetTotalSubdomains_C",(void (**)(void))&f);CHKERRQ(ierr); 682 if (f) { 683 ierr = (*f)(pc,N,is,is_local);CHKERRQ(ierr); 684 } 685 PetscFunctionReturn(0); 686 } 687 688 #undef __FUNCT__ 689 #define __FUNCT__ "PCASMSetOverlap" 690 /*@ 691 PCASMSetOverlap - Sets the overlap between a pair of subdomains for the 692 additive Schwarz preconditioner. Either all or no processors in the 693 PC communicator must call this routine. 694 695 Collective on PC 696 697 Input Parameters: 698 + pc - the preconditioner context 699 - ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1) 700 701 Options Database Key: 702 . -pc_asm_overlap <ovl> - Sets overlap 703 704 Notes: 705 By default the ASM preconditioner uses 1 block per processor. To use 706 multiple blocks per perocessor, see PCASMSetTotalSubdomains() and 707 PCASMSetLocalSubdomains() (and the option -pc_asm_blocks <blks>). 708 709 The overlap defaults to 1, so if one desires that no additional 710 overlap be computed beyond what may have been set with a call to 711 PCASMSetTotalSubdomains() or PCASMSetLocalSubdomains(), then ovl 712 must be set to be 0. In particular, if one does not explicitly set 713 the subdomains an application code, then all overlap would be computed 714 internally by PETSc, and using an overlap of 0 would result in an ASM 715 variant that is equivalent to the block Jacobi preconditioner. 716 717 Note that one can define initial index sets with any overlap via 718 PCASMSetTotalSubdomains() or PCASMSetLocalSubdomains(); the routine 719 PCASMSetOverlap() merely allows PETSc to extend that overlap further 720 if desired. 721 722 Level: intermediate 723 724 .keywords: PC, ASM, set, overlap 725 726 .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubKSP(), 727 PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains() 728 @*/ 729 PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetOverlap(PC pc,PetscInt ovl) 730 { 731 PetscErrorCode ierr,(*f)(PC,PetscInt); 732 733 PetscFunctionBegin; 734 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 735 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCASMSetOverlap_C",(void (**)(void))&f);CHKERRQ(ierr); 736 if (f) { 737 ierr = (*f)(pc,ovl);CHKERRQ(ierr); 738 } 739 PetscFunctionReturn(0); 740 } 741 742 #undef __FUNCT__ 743 #define __FUNCT__ "PCASMSetType" 744 /*@ 745 PCASMSetType - Sets the type of restriction and interpolation used 746 for local problems in the additive Schwarz method. 747 748 Collective on PC 749 750 Input Parameters: 751 + pc - the preconditioner context 752 - type - variant of ASM, one of 753 .vb 754 PC_ASM_BASIC - full interpolation and restriction 755 PC_ASM_RESTRICT - full restriction, local processor interpolation 756 PC_ASM_INTERPOLATE - full interpolation, local processor restriction 757 PC_ASM_NONE - local processor restriction and interpolation 758 .ve 759 760 Options Database Key: 761 . -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type 762 763 Level: intermediate 764 765 .keywords: PC, ASM, set, type 766 767 .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(), 768 PCASMCreateSubdomains2D() 769 @*/ 770 PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetType(PC pc,PCASMType type) 771 { 772 PetscErrorCode ierr,(*f)(PC,PCASMType); 773 774 PetscFunctionBegin; 775 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 776 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCASMSetType_C",(void (**)(void))&f);CHKERRQ(ierr); 777 if (f) { 778 ierr = (*f)(pc,type);CHKERRQ(ierr); 779 } 780 PetscFunctionReturn(0); 781 } 782 783 #undef __FUNCT__ 784 #define __FUNCT__ "PCASMGetSubKSP" 785 /*@C 786 PCASMGetSubKSP - Gets the local KSP contexts for all blocks on 787 this processor. 788 789 Collective on PC iff first_local is requested 790 791 Input Parameter: 792 . pc - the preconditioner context 793 794 Output Parameters: 795 + n_local - the number of blocks on this processor or PETSC_NULL 796 . first_local - the global number of the first block on this processor or PETSC_NULL, 797 all processors must request or all must pass PETSC_NULL 798 - ksp - the array of KSP contexts 799 800 Note: 801 After PCASMGetSubKSP() the array of KSPes is not to be freed 802 803 Currently for some matrix implementations only 1 block per processor 804 is supported. 805 806 You must call KSPSetUp() before calling PCASMGetSubKSP(). 807 808 Level: advanced 809 810 .keywords: PC, ASM, additive Schwarz, get, sub, KSP, context 811 812 .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap(), 813 PCASMCreateSubdomains2D(), 814 @*/ 815 PetscErrorCode PETSCKSP_DLLEXPORT PCASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[]) 816 { 817 PetscErrorCode ierr,(*f)(PC,PetscInt*,PetscInt*,KSP **); 818 819 PetscFunctionBegin; 820 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 821 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCASMGetSubKSP_C",(void (**)(void))&f);CHKERRQ(ierr); 822 if (f) { 823 ierr = (*f)(pc,n_local,first_local,ksp);CHKERRQ(ierr); 824 } else { 825 SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get subksp for this type of PC"); 826 } 827 828 PetscFunctionReturn(0); 829 } 830 831 /* -------------------------------------------------------------------------------------*/ 832 /*MC 833 PCASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with 834 its own KSP object. 835 836 Options Database Keys: 837 + -pc_asm_truelocal - Activates PCASMSetUseTrueLocal() 838 . -pc_asm_blocks <blks> - Sets total blocks 839 . -pc_asm_overlap <ovl> - Sets overlap 840 - -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type 841 842 IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you 843 will get a different convergence rate due to the default option of -pc_asm_type restrict. Use 844 -pc_asm_type basic to use the standard ASM. 845 846 Notes: Each processor can have one or more blocks, but a block cannot be shared by more 847 than one processor. Defaults to one block per processor. 848 849 To set options on the solvers for each block append -sub_ to all the KSP, and PC 850 options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly 851 852 To set the options on the solvers separate for each block call PCASMGetSubKSP() 853 and set the options directly on the resulting KSP object (you can access its PC 854 with KSPGetPC()) 855 856 857 Level: beginner 858 859 Concepts: additive Schwarz method 860 861 References: 862 An additive variant of the Schwarz alternating method for the case of many subregions 863 M Dryja, OB Widlund - Courant Institute, New York University Technical report 864 865 Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations, 866 Barry Smith, Petter Bjorstad, and William Gropp, Cambridge University Press, ISBN 0-521-49589-X. 867 868 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 869 PCBJACOBI, PCASMSetUseTrueLocal(), PCASMGetSubKSP(), PCASMSetLocalSubdomains(), 870 PCASMSetTotalSubdomains(), PCSetModifySubmatrices(), PCASMSetOverlap(), PCASMSetType() 871 872 M*/ 873 874 EXTERN_C_BEGIN 875 #undef __FUNCT__ 876 #define __FUNCT__ "PCCreate_ASM" 877 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ASM(PC pc) 878 { 879 PetscErrorCode ierr; 880 PC_ASM *osm; 881 882 PetscFunctionBegin; 883 ierr = PetscNewLog(pc,PC_ASM,&osm);CHKERRQ(ierr); 884 osm->n = PETSC_DECIDE; 885 osm->n_local = 0; 886 osm->n_local_true = 0; 887 osm->overlap = 1; 888 osm->ksp = 0; 889 osm->restriction = 0; 890 osm->prolongation = 0; 891 osm->x = 0; 892 osm->y = 0; 893 osm->is = 0; 894 osm->is_local = 0; 895 osm->mat = 0; 896 osm->pmat = 0; 897 osm->type = PC_ASM_RESTRICT; 898 osm->same_local_solves = PETSC_TRUE; 899 900 pc->data = (void*)osm; 901 pc->ops->apply = PCApply_ASM; 902 pc->ops->applytranspose = PCApplyTranspose_ASM; 903 pc->ops->setup = PCSetUp_ASM; 904 pc->ops->destroy = PCDestroy_ASM; 905 pc->ops->setfromoptions = PCSetFromOptions_ASM; 906 pc->ops->setuponblocks = PCSetUpOnBlocks_ASM; 907 pc->ops->view = PCView_ASM; 908 pc->ops->applyrichardson = 0; 909 910 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCASMSetLocalSubdomains_C","PCASMSetLocalSubdomains_ASM", 911 PCASMSetLocalSubdomains_ASM);CHKERRQ(ierr); 912 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCASMSetTotalSubdomains_C","PCASMSetTotalSubdomains_ASM", 913 PCASMSetTotalSubdomains_ASM);CHKERRQ(ierr); 914 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCASMSetOverlap_C","PCASMSetOverlap_ASM", 915 PCASMSetOverlap_ASM);CHKERRQ(ierr); 916 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCASMSetType_C","PCASMSetType_ASM", 917 PCASMSetType_ASM);CHKERRQ(ierr); 918 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCASMGetSubKSP_C","PCASMGetSubKSP_ASM", 919 PCASMGetSubKSP_ASM);CHKERRQ(ierr); 920 PetscFunctionReturn(0); 921 } 922 EXTERN_C_END 923 924 925 #undef __FUNCT__ 926 #define __FUNCT__ "PCASMCreateSubdomains" 927 /*@C 928 PCASMCreateSubdomains - Creates the index sets for the overlapping Schwarz 929 preconditioner for a any problem on a general grid. 930 931 Collective 932 933 Input Parameters: 934 + A - The global matrix operator 935 - n - the number of local blocks 936 937 Output Parameters: 938 . outis - the array of index sets defining the subdomains 939 940 Level: advanced 941 942 Note: this generates nonoverlapping subdomains; the PCASM will generate the overlap 943 from these if you use PCASMSetLocalSubdomains() 944 945 In the Fortran version you must provide the array outis[] already allocated of length n. 946 947 .keywords: PC, ASM, additive Schwarz, create, subdomains, unstructured grid 948 949 .seealso: PCASMSetLocalSubdomains(), PCASMDestroySubdomains() 950 @*/ 951 PetscErrorCode PETSCKSP_DLLEXPORT PCASMCreateSubdomains(Mat A, PetscInt n, IS* outis[]) 952 { 953 MatPartitioning mpart; 954 const char *prefix; 955 PetscErrorCode (*f)(Mat,PetscTruth*,MatReuse,Mat*); 956 PetscMPIInt size; 957 PetscInt i,j,rstart,rend,bs; 958 PetscTruth iscopy = PETSC_FALSE,isbaij = PETSC_FALSE,foundpart = PETSC_FALSE; 959 Mat Ad = PETSC_NULL, adj; 960 IS ispart,isnumb,*is; 961 PetscErrorCode ierr; 962 963 PetscFunctionBegin; 964 PetscValidHeaderSpecific(A,MAT_COOKIE,1); 965 PetscValidPointer(outis,3); 966 if (n < 1) {SETERRQ1(PETSC_ERR_ARG_WRONG,"number of local blocks must be > 0, n = %D",n);} 967 968 /* Get prefix, row distribution, and block size */ 969 ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr); 970 ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 971 ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 972 if (rstart/bs*bs != rstart || rend/bs*bs != rend) { 973 SETERRQ3(PETSC_ERR_ARG_WRONG,"bad row distribution [%D,%D) for matrix block size %D",rstart,rend,bs); 974 } 975 /* Get diagonal block from matrix if possible */ 976 ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr); 977 ierr = PetscObjectQueryFunction((PetscObject)A,"MatGetDiagonalBlock_C",(void (**)(void))&f);CHKERRQ(ierr); 978 if (f) { 979 ierr = (*f)(A,&iscopy,MAT_INITIAL_MATRIX,&Ad);CHKERRQ(ierr); 980 } else if (size == 1) { 981 iscopy = PETSC_FALSE; Ad = A; 982 } else { 983 iscopy = PETSC_FALSE; Ad = PETSC_NULL; 984 } 985 if (Ad) { 986 ierr = PetscTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);CHKERRQ(ierr); 987 if (!isbaij) {ierr = PetscTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);CHKERRQ(ierr);} 988 } 989 if (Ad && n > 1) { 990 PetscTruth match,done; 991 /* Try to setup a good matrix partitioning if available */ 992 ierr = MatPartitioningCreate(PETSC_COMM_SELF,&mpart);CHKERRQ(ierr); 993 ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr); 994 ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr); 995 ierr = PetscTypeCompare((PetscObject)mpart,MAT_PARTITIONING_CURRENT,&match);CHKERRQ(ierr); 996 if (!match) { 997 ierr = PetscTypeCompare((PetscObject)mpart,MAT_PARTITIONING_SQUARE,&match);CHKERRQ(ierr); 998 } 999 if (!match) { /* assume a "good" partitioner is available */ 1000 PetscInt na,*ia,*ja; 1001 ierr = MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr); 1002 if (done) { 1003 /* Build adjacency matrix by hand. Unfortunately a call to 1004 MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will 1005 remove the block-aij structure and we cannot expect 1006 MatPartitioning to split vertices as we need */ 1007 PetscInt i,j,*row,len,nnz,cnt,*iia=0,*jja=0; 1008 nnz = 0; 1009 for (i=0; i<na; i++) { /* count number of nonzeros */ 1010 len = ia[i+1] - ia[i]; 1011 row = ja + ia[i]; 1012 for (j=0; j<len; j++) { 1013 if (row[j] == i) { /* don't count diagonal */ 1014 len--; break; 1015 } 1016 } 1017 nnz += len; 1018 } 1019 ierr = PetscMalloc((na+1)*sizeof(PetscInt),&iia);CHKERRQ(ierr); 1020 ierr = PetscMalloc((nnz)*sizeof(PetscInt),&jja);CHKERRQ(ierr); 1021 nnz = 0; 1022 iia[0] = 0; 1023 for (i=0; i<na; i++) { /* fill adjacency */ 1024 cnt = 0; 1025 len = ia[i+1] - ia[i]; 1026 row = ja + ia[i]; 1027 for (j=0; j<len; j++) { 1028 if (row[j] != i) { /* if not diagonal */ 1029 jja[nnz+cnt++] = row[j]; 1030 } 1031 } 1032 nnz += cnt; 1033 iia[i+1] = nnz; 1034 } 1035 /* Partitioning of the adjacency matrix */ 1036 ierr = MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,PETSC_NULL,&adj);CHKERRQ(ierr); 1037 ierr = MatPartitioningSetAdjacency(mpart,adj);CHKERRQ(ierr); 1038 ierr = MatPartitioningSetNParts(mpart,n);CHKERRQ(ierr); 1039 ierr = MatPartitioningApply(mpart,&ispart);CHKERRQ(ierr); 1040 ierr = ISPartitioningToNumbering(ispart,&isnumb);CHKERRQ(ierr); 1041 ierr = MatDestroy(adj);CHKERRQ(ierr); 1042 foundpart = PETSC_TRUE; 1043 } 1044 ierr = MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr); 1045 } 1046 ierr = MatPartitioningDestroy(mpart);CHKERRQ(ierr); 1047 } 1048 if (iscopy) {ierr = MatDestroy(Ad);CHKERRQ(ierr);} 1049 1050 ierr = PetscMalloc(n*sizeof(IS),&is);CHKERRQ(ierr); 1051 *outis = is; 1052 1053 if (!foundpart) { 1054 1055 /* Partitioning by contiguous chunks of rows */ 1056 1057 PetscInt mbs = (rend-rstart)/bs; 1058 PetscInt start = rstart; 1059 for (i=0; i<n; i++) { 1060 PetscInt count = (mbs/n + ((mbs % n) > i)) * bs; 1061 ierr = ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);CHKERRQ(ierr); 1062 start += count; 1063 } 1064 1065 } else { 1066 1067 /* Partitioning by adjacency of diagonal block */ 1068 1069 const PetscInt *numbering; 1070 PetscInt *count,nidx,*indices,*newidx,start=0; 1071 /* Get node count in each partition */ 1072 ierr = PetscMalloc(n*sizeof(PetscInt),&count);CHKERRQ(ierr); 1073 ierr = ISPartitioningCount(ispart,n,count);CHKERRQ(ierr); 1074 if (isbaij && bs > 1) { /* adjust for the block-aij case */ 1075 for (i=0; i<n; i++) count[i] *= bs; 1076 } 1077 /* Build indices from node numbering */ 1078 ierr = ISGetLocalSize(isnumb,&nidx);CHKERRQ(ierr); 1079 ierr = PetscMalloc(nidx*sizeof(PetscInt),&indices);CHKERRQ(ierr); 1080 for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */ 1081 ierr = ISGetIndices(isnumb,&numbering);CHKERRQ(ierr); 1082 ierr = PetscSortIntWithPermutation(nidx,numbering,indices);CHKERRQ(ierr); 1083 ierr = ISRestoreIndices(isnumb,&numbering);CHKERRQ(ierr); 1084 if (isbaij && bs > 1) { /* adjust for the block-aij case */ 1085 ierr = PetscMalloc(nidx*bs*sizeof(PetscInt),&newidx);CHKERRQ(ierr); 1086 for (i=0; i<nidx; i++) 1087 for (j=0; j<bs; j++) 1088 newidx[i*bs+j] = indices[i]*bs + j; 1089 ierr = PetscFree(indices);CHKERRQ(ierr); 1090 nidx *= bs; 1091 indices = newidx; 1092 } 1093 /* Shift to get global indices */ 1094 for (i=0; i<nidx; i++) indices[i] += rstart; 1095 1096 /* Build the index sets for each block */ 1097 for (i=0; i<n; i++) { 1098 ierr = ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],&is[i]);CHKERRQ(ierr); 1099 ierr = ISSort(is[i]);CHKERRQ(ierr); 1100 start += count[i]; 1101 } 1102 1103 ierr = PetscFree(count); 1104 ierr = PetscFree(indices); 1105 ierr = ISDestroy(isnumb);CHKERRQ(ierr); 1106 ierr = ISDestroy(ispart);CHKERRQ(ierr); 1107 1108 } 1109 1110 PetscFunctionReturn(0); 1111 } 1112 1113 #undef __FUNCT__ 1114 #define __FUNCT__ "PCASMDestroySubdomains" 1115 /*@C 1116 PCASMDestroySubdomains - Destroys the index sets created with 1117 PCASMCreateSubdomains(). Should be called after setting subdomains 1118 with PCASMSetLocalSubdomains(). 1119 1120 Collective 1121 1122 Input Parameters: 1123 + n - the number of index sets 1124 . is - the array of index sets 1125 - is_local - the array of local index sets, can be PETSC_NULL 1126 1127 Level: advanced 1128 1129 .keywords: PC, ASM, additive Schwarz, create, subdomains, unstructured grid 1130 1131 .seealso: PCASMCreateSubdomains(), PCASMSetLocalSubdomains() 1132 @*/ 1133 PetscErrorCode PETSCKSP_DLLEXPORT PCASMDestroySubdomains(PetscInt n, IS is[], IS is_local[]) 1134 { 1135 PetscInt i; 1136 PetscErrorCode ierr; 1137 PetscFunctionBegin; 1138 if (n <= 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"n must be > 0: n = %D",n); 1139 PetscValidPointer(is,2); 1140 for (i=0; i<n; i++) { ierr = ISDestroy(is[i]);CHKERRQ(ierr); } 1141 ierr = PetscFree(is);CHKERRQ(ierr); 1142 if (is_local) { 1143 PetscValidPointer(is_local,3); 1144 for (i=0; i<n; i++) { ierr = ISDestroy(is_local[i]);CHKERRQ(ierr); } 1145 ierr = PetscFree(is_local);CHKERRQ(ierr); 1146 } 1147 PetscFunctionReturn(0); 1148 } 1149 1150 #undef __FUNCT__ 1151 #define __FUNCT__ "PCASMCreateSubdomains2D" 1152 /*@ 1153 PCASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz 1154 preconditioner for a two-dimensional problem on a regular grid. 1155 1156 Not Collective 1157 1158 Input Parameters: 1159 + m, n - the number of mesh points in the x and y directions 1160 . M, N - the number of subdomains in the x and y directions 1161 . dof - degrees of freedom per node 1162 - overlap - overlap in mesh lines 1163 1164 Output Parameters: 1165 + Nsub - the number of subdomains created 1166 - is - the array of index sets defining the subdomains 1167 1168 Note: 1169 Presently PCAMSCreateSubdomains2d() is valid only for sequential 1170 preconditioners. More general related routines are 1171 PCASMSetTotalSubdomains() and PCASMSetLocalSubdomains(). 1172 1173 Level: advanced 1174 1175 .keywords: PC, ASM, additive Schwarz, create, subdomains, 2D, regular grid 1176 1177 .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubKSP(), 1178 PCASMSetOverlap() 1179 @*/ 1180 PetscErrorCode PETSCKSP_DLLEXPORT PCASMCreateSubdomains2D(PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt dof,PetscInt overlap,PetscInt *Nsub,IS **is) 1181 { 1182 PetscInt i,j,height,width,ystart,xstart,yleft,yright,xleft,xright,loc_outter; 1183 PetscErrorCode ierr; 1184 PetscInt nidx,*idx,loc,ii,jj,count; 1185 1186 PetscFunctionBegin; 1187 if (dof != 1) SETERRQ(PETSC_ERR_SUP," "); 1188 1189 *Nsub = N*M; 1190 ierr = PetscMalloc((*Nsub)*sizeof(IS *),is);CHKERRQ(ierr); 1191 ystart = 0; 1192 loc_outter = 0; 1193 for (i=0; i<N; i++) { 1194 height = n/N + ((n % N) > i); /* height of subdomain */ 1195 if (height < 2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many N subdomains for mesh dimension n"); 1196 yleft = ystart - overlap; if (yleft < 0) yleft = 0; 1197 yright = ystart + height + overlap; if (yright > n) yright = n; 1198 xstart = 0; 1199 for (j=0; j<M; j++) { 1200 width = m/M + ((m % M) > j); /* width of subdomain */ 1201 if (width < 2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many M subdomains for mesh dimension m"); 1202 xleft = xstart - overlap; if (xleft < 0) xleft = 0; 1203 xright = xstart + width + overlap; if (xright > m) xright = m; 1204 nidx = (xright - xleft)*(yright - yleft); 1205 ierr = PetscMalloc(nidx*sizeof(PetscInt),&idx);CHKERRQ(ierr); 1206 loc = 0; 1207 for (ii=yleft; ii<yright; ii++) { 1208 count = m*ii + xleft; 1209 for (jj=xleft; jj<xright; jj++) { 1210 idx[loc++] = count++; 1211 } 1212 } 1213 ierr = ISCreateGeneral(PETSC_COMM_SELF,nidx,idx,(*is)+loc_outter++);CHKERRQ(ierr); 1214 ierr = PetscFree(idx);CHKERRQ(ierr); 1215 xstart += width; 1216 } 1217 ystart += height; 1218 } 1219 for (i=0; i<*Nsub; i++) { ierr = ISSort((*is)[i]);CHKERRQ(ierr); } 1220 PetscFunctionReturn(0); 1221 } 1222 1223 #undef __FUNCT__ 1224 #define __FUNCT__ "PCASMGetLocalSubdomains" 1225 /*@C 1226 PCASMGetLocalSubdomains - Gets the local subdomains (for this processor 1227 only) for the additive Schwarz preconditioner. 1228 1229 Collective on PC 1230 1231 Input Parameter: 1232 . pc - the preconditioner context 1233 1234 Output Parameters: 1235 + n - the number of subdomains for this processor (default value = 1) 1236 . is - the index sets that define the subdomains for this processor 1237 - is_local - the index sets that define the local part of the subdomains for this processor (can be PETSC_NULL) 1238 1239 1240 Notes: 1241 The IS numbering is in the parallel, global numbering of the vector. 1242 1243 Level: advanced 1244 1245 .keywords: PC, ASM, set, local, subdomains, additive Schwarz 1246 1247 .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(), 1248 PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubmatrices() 1249 @*/ 1250 PetscErrorCode PETSCKSP_DLLEXPORT PCASMGetLocalSubdomains(PC pc,PetscInt *n,IS *is[],IS *is_local[]) 1251 { 1252 PC_ASM *osm; 1253 PetscErrorCode ierr; 1254 PetscTruth match; 1255 1256 PetscFunctionBegin; 1257 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 1258 PetscValidIntPointer(n,2); 1259 if (is) PetscValidPointer(is,3); 1260 ierr = PetscTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr); 1261 if (!match) { 1262 if (n) *n = 0; 1263 if (is) *is = PETSC_NULL; 1264 } else { 1265 osm = (PC_ASM*)pc->data; 1266 if (n) *n = osm->n_local_true; 1267 if (is) *is = osm->is; 1268 if (is_local) *is_local = osm->is_local; 1269 } 1270 PetscFunctionReturn(0); 1271 } 1272 1273 #undef __FUNCT__ 1274 #define __FUNCT__ "PCASMGetLocalSubmatrices" 1275 /*@C 1276 PCASMGetLocalSubmatrices - Gets the local submatrices (for this processor 1277 only) for the additive Schwarz preconditioner. 1278 1279 Collective on PC 1280 1281 Input Parameter: 1282 . pc - the preconditioner context 1283 1284 Output Parameters: 1285 + n - the number of matrices for this processor (default value = 1) 1286 - mat - the matrices 1287 1288 1289 Level: advanced 1290 1291 .keywords: PC, ASM, set, local, subdomains, additive Schwarz, block Jacobi 1292 1293 .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(), 1294 PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains() 1295 @*/ 1296 PetscErrorCode PETSCKSP_DLLEXPORT PCASMGetLocalSubmatrices(PC pc,PetscInt *n,Mat *mat[]) 1297 { 1298 PC_ASM *osm; 1299 PetscErrorCode ierr; 1300 PetscTruth match; 1301 1302 PetscFunctionBegin; 1303 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 1304 PetscValidIntPointer(n,2); 1305 if (mat) PetscValidPointer(mat,3); 1306 if (!pc->setupcalled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUP() or PCSetUp()."); 1307 ierr = PetscTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr); 1308 if (!match) { 1309 if (n) *n = 0; 1310 if (mat) *mat = PETSC_NULL; 1311 } else { 1312 osm = (PC_ASM*)pc->data; 1313 if (n) *n = osm->n_local_true; 1314 if (mat) *mat = osm->pmat; 1315 } 1316 PetscFunctionReturn(0); 1317 } 1318