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