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