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 ierr = PetscFree(domain_names[d]); CHKERRQ(ierr); 227 ierr = ISDestroy(&inner_domain_is[d]); CHKERRQ(ierr); 228 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 if (osm->is) { 562 ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr); 563 osm->is = 0; 564 osm->is_local = 0; 565 } 566 PetscFunctionReturn(0); 567 } 568 569 #undef __FUNCT__ 570 #define __FUNCT__ "PCDestroy_ASM" 571 static PetscErrorCode PCDestroy_ASM(PC pc) 572 { 573 PC_ASM *osm = (PC_ASM*)pc->data; 574 PetscErrorCode ierr; 575 PetscInt i; 576 577 PetscFunctionBegin; 578 ierr = PCReset_ASM(pc);CHKERRQ(ierr); 579 if (osm->ksp) { 580 for (i=0; i<osm->n_local_true; i++) { 581 ierr = KSPDestroy(&osm->ksp[i]);CHKERRQ(ierr); 582 } 583 ierr = PetscFree(osm->ksp);CHKERRQ(ierr); 584 } 585 ierr = PetscFree(pc->data);CHKERRQ(ierr); 586 PetscFunctionReturn(0); 587 } 588 589 #undef __FUNCT__ 590 #define __FUNCT__ "PCSetFromOptions_ASM" 591 static PetscErrorCode PCSetFromOptions_ASM(PC pc) 592 { 593 PC_ASM *osm = (PC_ASM*)pc->data; 594 PetscErrorCode ierr; 595 PetscInt blocks,ovl; 596 PetscBool symset,flg; 597 PCASMType asmtype; 598 599 PetscFunctionBegin; 600 /* set the type to symmetric if matrix is symmetric */ 601 if (!osm->type_set && pc->pmat) { 602 ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr); 603 if (symset && flg) { osm->type = PC_ASM_BASIC; } 604 } 605 ierr = PetscOptionsHead("Additive Schwarz options");CHKERRQ(ierr); 606 ierr = PetscOptionsInt("-pc_asm_blocks","Number of subdomains","PCASMSetTotalSubdomains",osm->n,&blocks,&flg);CHKERRQ(ierr); 607 if (flg) {ierr = PCASMSetTotalSubdomains(pc,blocks,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); } 608 ierr = PetscOptionsInt("-pc_asm_overlap","Number of grid points overlap","PCASMSetOverlap",osm->overlap,&ovl,&flg);CHKERRQ(ierr); 609 if (flg) {ierr = PCASMSetOverlap(pc,ovl);CHKERRQ(ierr); } 610 flg = PETSC_FALSE; 611 ierr = PetscOptionsEnum("-pc_asm_type","Type of restriction/extension","PCASMSetType",PCASMTypes,(PetscEnum)osm->type,(PetscEnum*)&asmtype,&flg);CHKERRQ(ierr); 612 if (flg) {ierr = PCASMSetType(pc,asmtype);CHKERRQ(ierr); } 613 ierr = PetscOptionsTail();CHKERRQ(ierr); 614 PetscFunctionReturn(0); 615 } 616 617 /*------------------------------------------------------------------------------------*/ 618 619 EXTERN_C_BEGIN 620 #undef __FUNCT__ 621 #define __FUNCT__ "PCASMSetLocalSubdomains_ASM" 622 PetscErrorCode PCASMSetLocalSubdomains_ASM(PC pc,PetscInt n,IS is[],IS is_local[]) 623 { 624 PC_ASM *osm = (PC_ASM*)pc->data; 625 PetscErrorCode ierr; 626 PetscInt i; 627 628 PetscFunctionBegin; 629 if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more blocks, n = %D",n); 630 if (pc->setupcalled && (n != osm->n_local_true || is)) 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 if (osm->is) { 640 ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr); 641 } 642 osm->n_local_true = n; 643 osm->is = 0; 644 osm->is_local = 0; 645 if (is) { 646 ierr = PetscMalloc(n*sizeof(IS),&osm->is);CHKERRQ(ierr); 647 for (i=0; i<n; i++) { osm->is[i] = is[i]; } 648 /* Flag indicating that the user has set overlapping subdomains so PCASM should not increase their size. */ 649 osm->overlap = -1; 650 } 651 if (is_local) { 652 ierr = PetscMalloc(n*sizeof(IS),&osm->is_local);CHKERRQ(ierr); 653 for (i=0; i<n; i++) { osm->is_local[i] = is_local[i]; } 654 } 655 } 656 PetscFunctionReturn(0); 657 } 658 EXTERN_C_END 659 660 EXTERN_C_BEGIN 661 #undef __FUNCT__ 662 #define __FUNCT__ "PCASMSetTotalSubdomains_ASM" 663 PetscErrorCode PCASMSetTotalSubdomains_ASM(PC pc,PetscInt N,IS *is,IS *is_local) 664 { 665 PC_ASM *osm = (PC_ASM*)pc->data; 666 PetscErrorCode ierr; 667 PetscMPIInt rank,size; 668 PetscInt n; 669 670 PetscFunctionBegin; 671 if (N < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of total blocks must be > 0, N = %D",N); 672 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."); 673 674 /* 675 Split the subdomains equally among all processors 676 */ 677 ierr = MPI_Comm_rank(((PetscObject)pc)->comm,&rank);CHKERRQ(ierr); 678 ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr); 679 n = N/size + ((N % size) > rank); 680 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); 681 if (pc->setupcalled && n != osm->n_local_true) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PCASMSetTotalSubdomains() should be called before PCSetUp()."); 682 if (!pc->setupcalled) { 683 if (osm->is) { 684 ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr); 685 } 686 osm->n_local_true = n; 687 osm->is = 0; 688 osm->is_local = 0; 689 } 690 PetscFunctionReturn(0); 691 } 692 EXTERN_C_END 693 694 EXTERN_C_BEGIN 695 #undef __FUNCT__ 696 #define __FUNCT__ "PCASMSetOverlap_ASM" 697 PetscErrorCode PCASMSetOverlap_ASM(PC pc,PetscInt ovl) 698 { 699 PC_ASM *osm = (PC_ASM*)pc->data; 700 701 PetscFunctionBegin; 702 if (ovl < 0) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested"); 703 if (pc->setupcalled && ovl != osm->overlap) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"PCASMSetOverlap() should be called before PCSetUp()."); 704 if (!pc->setupcalled) { 705 osm->overlap = ovl; 706 } 707 PetscFunctionReturn(0); 708 } 709 EXTERN_C_END 710 711 EXTERN_C_BEGIN 712 #undef __FUNCT__ 713 #define __FUNCT__ "PCASMSetType_ASM" 714 PetscErrorCode PCASMSetType_ASM(PC pc,PCASMType type) 715 { 716 PC_ASM *osm = (PC_ASM*)pc->data; 717 718 PetscFunctionBegin; 719 osm->type = type; 720 osm->type_set = PETSC_TRUE; 721 PetscFunctionReturn(0); 722 } 723 EXTERN_C_END 724 725 EXTERN_C_BEGIN 726 #undef __FUNCT__ 727 #define __FUNCT__ "PCASMSetSortIndices_ASM" 728 PetscErrorCode PCASMSetSortIndices_ASM(PC pc,PetscBool doSort) 729 { 730 PC_ASM *osm = (PC_ASM*)pc->data; 731 732 PetscFunctionBegin; 733 osm->sort_indices = doSort; 734 PetscFunctionReturn(0); 735 } 736 EXTERN_C_END 737 738 EXTERN_C_BEGIN 739 #undef __FUNCT__ 740 #define __FUNCT__ "PCASMGetSubKSP_ASM" 741 PetscErrorCode PCASMGetSubKSP_ASM(PC pc,PetscInt *n_local,PetscInt *first_local,KSP **ksp) 742 { 743 PC_ASM *osm = (PC_ASM*)pc->data; 744 PetscErrorCode ierr; 745 746 PetscFunctionBegin; 747 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"); 748 749 if (n_local) { 750 *n_local = osm->n_local_true; 751 } 752 if (first_local) { 753 ierr = MPI_Scan(&osm->n_local_true,first_local,1,MPIU_INT,MPI_SUM,((PetscObject)pc)->comm);CHKERRQ(ierr); 754 *first_local -= osm->n_local_true; 755 } 756 if (ksp) { 757 /* Assume that local solves are now different; not necessarily 758 true though! This flag is used only for PCView_ASM() */ 759 *ksp = osm->ksp; 760 osm->same_local_solves = PETSC_FALSE; 761 } 762 PetscFunctionReturn(0); 763 } 764 EXTERN_C_END 765 766 767 #undef __FUNCT__ 768 #define __FUNCT__ "PCASMSetLocalSubdomains" 769 /*@C 770 PCASMSetLocalSubdomains - Sets the local subdomains (for this processor only) for the additive Schwarz preconditioner. 771 772 Collective on PC 773 774 Input Parameters: 775 + pc - the preconditioner context 776 . n - the number of subdomains for this processor (default value = 1) 777 . is - the index set that defines the subdomains for this processor 778 (or PETSC_NULL for PETSc to determine subdomains) 779 - is_local - the index sets that define the local part of the subdomains for this processor 780 (or PETSC_NULL to use the default of 1 subdomain per process) 781 782 Notes: 783 The IS numbering is in the parallel, global numbering of the vector for both is and is_local 784 785 By default the ASM preconditioner uses 1 block per processor. 786 787 Use PCASMSetTotalSubdomains() to set the subdomains for all processors. 788 789 Level: advanced 790 791 .keywords: PC, ASM, set, local, subdomains, additive Schwarz 792 793 .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(), 794 PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains() 795 @*/ 796 PetscErrorCode PCASMSetLocalSubdomains(PC pc,PetscInt n,IS is[],IS is_local[]) 797 { 798 PetscErrorCode ierr; 799 800 PetscFunctionBegin; 801 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 802 ierr = PetscTryMethod(pc,"PCASMSetLocalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,is,is_local));CHKERRQ(ierr); 803 PetscFunctionReturn(0); 804 } 805 806 #undef __FUNCT__ 807 #define __FUNCT__ "PCASMSetTotalSubdomains" 808 /*@C 809 PCASMSetTotalSubdomains - Sets the subdomains for all processors for the 810 additive Schwarz preconditioner. Either all or no processors in the 811 PC communicator must call this routine, with the same index sets. 812 813 Collective on PC 814 815 Input Parameters: 816 + pc - the preconditioner context 817 . N - the number of subdomains for all processors 818 . is - the index sets that define the subdomains for all processors 819 (or PETSC_NULL to ask PETSc to compe up with subdomains) 820 - is_local - the index sets that define the local part of the subdomains for this processor 821 (or PETSC_NULL to use the default of 1 subdomain per process) 822 823 Options Database Key: 824 To set the total number of subdomain blocks rather than specify the 825 index sets, use the option 826 . -pc_asm_blocks <blks> - Sets total blocks 827 828 Notes: 829 Currently you cannot use this to set the actual subdomains with the argument is. 830 831 By default the ASM preconditioner uses 1 block per processor. 832 833 These index sets cannot be destroyed until after completion of the 834 linear solves for which the ASM preconditioner is being used. 835 836 Use PCASMSetLocalSubdomains() to set local subdomains. 837 838 The IS numbering is in the parallel, global numbering of the vector for both is and is_local 839 840 Level: advanced 841 842 .keywords: PC, ASM, set, total, global, subdomains, additive Schwarz 843 844 .seealso: PCASMSetLocalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(), 845 PCASMCreateSubdomains2D() 846 @*/ 847 PetscErrorCode PCASMSetTotalSubdomains(PC pc,PetscInt N,IS is[],IS is_local[]) 848 { 849 PetscErrorCode ierr; 850 851 PetscFunctionBegin; 852 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 853 ierr = PetscTryMethod(pc,"PCASMSetTotalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,N,is,is_local));CHKERRQ(ierr); 854 PetscFunctionReturn(0); 855 } 856 857 #undef __FUNCT__ 858 #define __FUNCT__ "PCASMSetOverlap" 859 /*@ 860 PCASMSetOverlap - Sets the overlap between a pair of subdomains for the 861 additive Schwarz preconditioner. Either all or no processors in the 862 PC communicator must call this routine. 863 864 Logically Collective on PC 865 866 Input Parameters: 867 + pc - the preconditioner context 868 - ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1) 869 870 Options Database Key: 871 . -pc_asm_overlap <ovl> - Sets overlap 872 873 Notes: 874 By default the ASM preconditioner uses 1 block per processor. To use 875 multiple blocks per perocessor, see PCASMSetTotalSubdomains() and 876 PCASMSetLocalSubdomains() (and the option -pc_asm_blocks <blks>). 877 878 The overlap defaults to 1, so if one desires that no additional 879 overlap be computed beyond what may have been set with a call to 880 PCASMSetTotalSubdomains() or PCASMSetLocalSubdomains(), then ovl 881 must be set to be 0. In particular, if one does not explicitly set 882 the subdomains an application code, then all overlap would be computed 883 internally by PETSc, and using an overlap of 0 would result in an ASM 884 variant that is equivalent to the block Jacobi preconditioner. 885 886 Note that one can define initial index sets with any overlap via 887 PCASMSetTotalSubdomains() or PCASMSetLocalSubdomains(); the routine 888 PCASMSetOverlap() merely allows PETSc to extend that overlap further 889 if desired. 890 891 Level: intermediate 892 893 .keywords: PC, ASM, set, overlap 894 895 .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubKSP(), 896 PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains() 897 @*/ 898 PetscErrorCode PCASMSetOverlap(PC pc,PetscInt ovl) 899 { 900 PetscErrorCode ierr; 901 902 PetscFunctionBegin; 903 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 904 PetscValidLogicalCollectiveInt(pc,ovl,2); 905 ierr = PetscTryMethod(pc,"PCASMSetOverlap_C",(PC,PetscInt),(pc,ovl));CHKERRQ(ierr); 906 PetscFunctionReturn(0); 907 } 908 909 #undef __FUNCT__ 910 #define __FUNCT__ "PCASMSetType" 911 /*@ 912 PCASMSetType - Sets the type of restriction and interpolation used 913 for local problems in the additive Schwarz method. 914 915 Logically Collective on PC 916 917 Input Parameters: 918 + pc - the preconditioner context 919 - type - variant of ASM, one of 920 .vb 921 PC_ASM_BASIC - full interpolation and restriction 922 PC_ASM_RESTRICT - full restriction, local processor interpolation 923 PC_ASM_INTERPOLATE - full interpolation, local processor restriction 924 PC_ASM_NONE - local processor restriction and interpolation 925 .ve 926 927 Options Database Key: 928 . -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type 929 930 Level: intermediate 931 932 .keywords: PC, ASM, set, type 933 934 .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(), 935 PCASMCreateSubdomains2D() 936 @*/ 937 PetscErrorCode PCASMSetType(PC pc,PCASMType type) 938 { 939 PetscErrorCode ierr; 940 941 PetscFunctionBegin; 942 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 943 PetscValidLogicalCollectiveEnum(pc,type,2); 944 ierr = PetscTryMethod(pc,"PCASMSetType_C",(PC,PCASMType),(pc,type));CHKERRQ(ierr); 945 PetscFunctionReturn(0); 946 } 947 948 #undef __FUNCT__ 949 #define __FUNCT__ "PCASMSetSortIndices" 950 /*@ 951 PCASMSetSortIndices - Determines whether subdomain indices are sorted. 952 953 Logically Collective on PC 954 955 Input Parameters: 956 + pc - the preconditioner context 957 - doSort - sort the subdomain indices 958 959 Level: intermediate 960 961 .keywords: PC, ASM, set, type 962 963 .seealso: PCASMSetLocalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(), 964 PCASMCreateSubdomains2D() 965 @*/ 966 PetscErrorCode PCASMSetSortIndices(PC pc,PetscBool doSort) 967 { 968 PetscErrorCode ierr; 969 970 PetscFunctionBegin; 971 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 972 PetscValidLogicalCollectiveBool(pc,doSort,2); 973 ierr = PetscTryMethod(pc,"PCASMSetSortIndices_C",(PC,PetscBool),(pc,doSort));CHKERRQ(ierr); 974 PetscFunctionReturn(0); 975 } 976 977 #undef __FUNCT__ 978 #define __FUNCT__ "PCASMGetSubKSP" 979 /*@C 980 PCASMGetSubKSP - Gets the local KSP contexts for all blocks on 981 this processor. 982 983 Collective on PC iff first_local is requested 984 985 Input Parameter: 986 . pc - the preconditioner context 987 988 Output Parameters: 989 + n_local - the number of blocks on this processor or PETSC_NULL 990 . first_local - the global number of the first block on this processor or PETSC_NULL, 991 all processors must request or all must pass PETSC_NULL 992 - ksp - the array of KSP contexts 993 994 Note: 995 After PCASMGetSubKSP() the array of KSPes is not to be freed 996 997 Currently for some matrix implementations only 1 block per processor 998 is supported. 999 1000 You must call KSPSetUp() before calling PCASMGetSubKSP(). 1001 1002 Level: advanced 1003 1004 .keywords: PC, ASM, additive Schwarz, get, sub, KSP, context 1005 1006 .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap(), 1007 PCASMCreateSubdomains2D(), 1008 @*/ 1009 PetscErrorCode PCASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[]) 1010 { 1011 PetscErrorCode ierr; 1012 1013 PetscFunctionBegin; 1014 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1015 ierr = PetscUseMethod(pc,"PCASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));CHKERRQ(ierr); 1016 PetscFunctionReturn(0); 1017 } 1018 1019 /* -------------------------------------------------------------------------------------*/ 1020 /*MC 1021 PCASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with 1022 its own KSP object. 1023 1024 Options Database Keys: 1025 + -pc_asm_truelocal - Activates PCASMSetUseTrueLocal() 1026 . -pc_asm_blocks <blks> - Sets total blocks 1027 . -pc_asm_overlap <ovl> - Sets overlap 1028 - -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type 1029 1030 IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you 1031 will get a different convergence rate due to the default option of -pc_asm_type restrict. Use 1032 -pc_asm_type basic to use the standard ASM. 1033 1034 Notes: Each processor can have one or more blocks, but a block cannot be shared by more 1035 than one processor. Defaults to one block per processor. 1036 1037 To set options on the solvers for each block append -sub_ to all the KSP, and PC 1038 options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly 1039 1040 To set the options on the solvers separate for each block call PCASMGetSubKSP() 1041 and set the options directly on the resulting KSP object (you can access its PC 1042 with KSPGetPC()) 1043 1044 1045 Level: beginner 1046 1047 Concepts: additive Schwarz method 1048 1049 References: 1050 An additive variant of the Schwarz alternating method for the case of many subregions 1051 M Dryja, OB Widlund - Courant Institute, New York University Technical report 1052 1053 Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations, 1054 Barry Smith, Petter Bjorstad, and William Gropp, Cambridge University Press, ISBN 0-521-49589-X. 1055 1056 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 1057 PCBJACOBI, PCASMSetUseTrueLocal(), PCASMGetSubKSP(), PCASMSetLocalSubdomains(), 1058 PCASMSetTotalSubdomains(), PCSetModifySubmatrices(), PCASMSetOverlap(), PCASMSetType() 1059 1060 M*/ 1061 1062 EXTERN_C_BEGIN 1063 #undef __FUNCT__ 1064 #define __FUNCT__ "PCCreate_ASM" 1065 PetscErrorCode PCCreate_ASM(PC pc) 1066 { 1067 PetscErrorCode ierr; 1068 PC_ASM *osm; 1069 1070 PetscFunctionBegin; 1071 ierr = PetscNewLog(pc,PC_ASM,&osm);CHKERRQ(ierr); 1072 osm->n = PETSC_DECIDE; 1073 osm->n_local = 0; 1074 osm->n_local_true = PETSC_DECIDE; 1075 osm->overlap = 1; 1076 osm->ksp = 0; 1077 osm->restriction = 0; 1078 osm->localization = 0; 1079 osm->prolongation = 0; 1080 osm->x = 0; 1081 osm->y = 0; 1082 osm->y_local = 0; 1083 osm->is = 0; 1084 osm->is_local = 0; 1085 osm->mat = 0; 1086 osm->pmat = 0; 1087 osm->type = PC_ASM_RESTRICT; 1088 osm->same_local_solves = PETSC_TRUE; 1089 osm->sort_indices = PETSC_TRUE; 1090 1091 pc->data = (void*)osm; 1092 pc->ops->apply = PCApply_ASM; 1093 pc->ops->applytranspose = PCApplyTranspose_ASM; 1094 pc->ops->setup = PCSetUp_ASM; 1095 pc->ops->reset = PCReset_ASM; 1096 pc->ops->destroy = PCDestroy_ASM; 1097 pc->ops->setfromoptions = PCSetFromOptions_ASM; 1098 pc->ops->setuponblocks = PCSetUpOnBlocks_ASM; 1099 pc->ops->view = PCView_ASM; 1100 pc->ops->applyrichardson = 0; 1101 1102 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCASMSetLocalSubdomains_C","PCASMSetLocalSubdomains_ASM", 1103 PCASMSetLocalSubdomains_ASM);CHKERRQ(ierr); 1104 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCASMSetTotalSubdomains_C","PCASMSetTotalSubdomains_ASM", 1105 PCASMSetTotalSubdomains_ASM);CHKERRQ(ierr); 1106 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCASMSetOverlap_C","PCASMSetOverlap_ASM", 1107 PCASMSetOverlap_ASM);CHKERRQ(ierr); 1108 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCASMSetType_C","PCASMSetType_ASM", 1109 PCASMSetType_ASM);CHKERRQ(ierr); 1110 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCASMSetSortIndices_C","PCASMSetSortIndices_ASM", 1111 PCASMSetSortIndices_ASM);CHKERRQ(ierr); 1112 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCASMGetSubKSP_C","PCASMGetSubKSP_ASM", 1113 PCASMGetSubKSP_ASM);CHKERRQ(ierr); 1114 PetscFunctionReturn(0); 1115 } 1116 EXTERN_C_END 1117 1118 1119 #undef __FUNCT__ 1120 #define __FUNCT__ "PCASMCreateSubdomains" 1121 /*@C 1122 PCASMCreateSubdomains - Creates the index sets for the overlapping Schwarz 1123 preconditioner for a any problem on a general grid. 1124 1125 Collective 1126 1127 Input Parameters: 1128 + A - The global matrix operator 1129 - n - the number of local blocks 1130 1131 Output Parameters: 1132 . outis - the array of index sets defining the subdomains 1133 1134 Level: advanced 1135 1136 Note: this generates nonoverlapping subdomains; the PCASM will generate the overlap 1137 from these if you use PCASMSetLocalSubdomains() 1138 1139 In the Fortran version you must provide the array outis[] already allocated of length n. 1140 1141 .keywords: PC, ASM, additive Schwarz, create, subdomains, unstructured grid 1142 1143 .seealso: PCASMSetLocalSubdomains(), PCASMDestroySubdomains() 1144 @*/ 1145 PetscErrorCode PCASMCreateSubdomains(Mat A, PetscInt n, IS* outis[]) 1146 { 1147 MatPartitioning mpart; 1148 const char *prefix; 1149 PetscErrorCode (*f)(Mat,Mat*); 1150 PetscMPIInt size; 1151 PetscInt i,j,rstart,rend,bs; 1152 PetscBool isbaij = PETSC_FALSE,foundpart = PETSC_FALSE; 1153 Mat Ad = PETSC_NULL, adj; 1154 IS ispart,isnumb,*is; 1155 PetscErrorCode ierr; 1156 1157 PetscFunctionBegin; 1158 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 1159 PetscValidPointer(outis,3); 1160 if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local blocks must be > 0, n = %D",n); 1161 1162 /* Get prefix, row distribution, and block size */ 1163 ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr); 1164 ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 1165 ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 1166 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); 1167 1168 /* Get diagonal block from matrix if possible */ 1169 ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr); 1170 ierr = PetscObjectQueryFunction((PetscObject)A,"MatGetDiagonalBlock_C",(void (**)(void))&f);CHKERRQ(ierr); 1171 if (f) { 1172 ierr = MatGetDiagonalBlock(A,&Ad);CHKERRQ(ierr); 1173 } else if (size == 1) { 1174 Ad = A; 1175 } 1176 if (Ad) { 1177 ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);CHKERRQ(ierr); 1178 if (!isbaij) {ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);CHKERRQ(ierr);} 1179 } 1180 if (Ad && n > 1) { 1181 PetscBool match,done; 1182 /* Try to setup a good matrix partitioning if available */ 1183 ierr = MatPartitioningCreate(PETSC_COMM_SELF,&mpart);CHKERRQ(ierr); 1184 ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr); 1185 ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr); 1186 ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match);CHKERRQ(ierr); 1187 if (!match) { 1188 ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match);CHKERRQ(ierr); 1189 } 1190 if (!match) { /* assume a "good" partitioner is available */ 1191 PetscInt na,*ia,*ja; 1192 ierr = MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr); 1193 if (done) { 1194 /* Build adjacency matrix by hand. Unfortunately a call to 1195 MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will 1196 remove the block-aij structure and we cannot expect 1197 MatPartitioning to split vertices as we need */ 1198 PetscInt i,j,*row,len,nnz,cnt,*iia=0,*jja=0; 1199 nnz = 0; 1200 for (i=0; i<na; i++) { /* count number of nonzeros */ 1201 len = ia[i+1] - ia[i]; 1202 row = ja + ia[i]; 1203 for (j=0; j<len; j++) { 1204 if (row[j] == i) { /* don't count diagonal */ 1205 len--; break; 1206 } 1207 } 1208 nnz += len; 1209 } 1210 ierr = PetscMalloc((na+1)*sizeof(PetscInt),&iia);CHKERRQ(ierr); 1211 ierr = PetscMalloc((nnz)*sizeof(PetscInt),&jja);CHKERRQ(ierr); 1212 nnz = 0; 1213 iia[0] = 0; 1214 for (i=0; i<na; i++) { /* fill adjacency */ 1215 cnt = 0; 1216 len = ia[i+1] - ia[i]; 1217 row = ja + ia[i]; 1218 for (j=0; j<len; j++) { 1219 if (row[j] != i) { /* if not diagonal */ 1220 jja[nnz+cnt++] = row[j]; 1221 } 1222 } 1223 nnz += cnt; 1224 iia[i+1] = nnz; 1225 } 1226 /* Partitioning of the adjacency matrix */ 1227 ierr = MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,PETSC_NULL,&adj);CHKERRQ(ierr); 1228 ierr = MatPartitioningSetAdjacency(mpart,adj);CHKERRQ(ierr); 1229 ierr = MatPartitioningSetNParts(mpart,n);CHKERRQ(ierr); 1230 ierr = MatPartitioningApply(mpart,&ispart);CHKERRQ(ierr); 1231 ierr = ISPartitioningToNumbering(ispart,&isnumb);CHKERRQ(ierr); 1232 ierr = MatDestroy(&adj);CHKERRQ(ierr); 1233 foundpart = PETSC_TRUE; 1234 } 1235 ierr = MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr); 1236 } 1237 ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr); 1238 } 1239 1240 ierr = PetscMalloc(n*sizeof(IS),&is);CHKERRQ(ierr); 1241 *outis = is; 1242 1243 if (!foundpart) { 1244 1245 /* Partitioning by contiguous chunks of rows */ 1246 1247 PetscInt mbs = (rend-rstart)/bs; 1248 PetscInt start = rstart; 1249 for (i=0; i<n; i++) { 1250 PetscInt count = (mbs/n + ((mbs % n) > i)) * bs; 1251 ierr = ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);CHKERRQ(ierr); 1252 start += count; 1253 } 1254 1255 } else { 1256 1257 /* Partitioning by adjacency of diagonal block */ 1258 1259 const PetscInt *numbering; 1260 PetscInt *count,nidx,*indices,*newidx,start=0; 1261 /* Get node count in each partition */ 1262 ierr = PetscMalloc(n*sizeof(PetscInt),&count);CHKERRQ(ierr); 1263 ierr = ISPartitioningCount(ispart,n,count);CHKERRQ(ierr); 1264 if (isbaij && bs > 1) { /* adjust for the block-aij case */ 1265 for (i=0; i<n; i++) count[i] *= bs; 1266 } 1267 /* Build indices from node numbering */ 1268 ierr = ISGetLocalSize(isnumb,&nidx);CHKERRQ(ierr); 1269 ierr = PetscMalloc(nidx*sizeof(PetscInt),&indices);CHKERRQ(ierr); 1270 for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */ 1271 ierr = ISGetIndices(isnumb,&numbering);CHKERRQ(ierr); 1272 ierr = PetscSortIntWithPermutation(nidx,numbering,indices);CHKERRQ(ierr); 1273 ierr = ISRestoreIndices(isnumb,&numbering);CHKERRQ(ierr); 1274 if (isbaij && bs > 1) { /* adjust for the block-aij case */ 1275 ierr = PetscMalloc(nidx*bs*sizeof(PetscInt),&newidx);CHKERRQ(ierr); 1276 for (i=0; i<nidx; i++) 1277 for (j=0; j<bs; j++) 1278 newidx[i*bs+j] = indices[i]*bs + j; 1279 ierr = PetscFree(indices);CHKERRQ(ierr); 1280 nidx *= bs; 1281 indices = newidx; 1282 } 1283 /* Shift to get global indices */ 1284 for (i=0; i<nidx; i++) indices[i] += rstart; 1285 1286 /* Build the index sets for each block */ 1287 for (i=0; i<n; i++) { 1288 ierr = ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]);CHKERRQ(ierr); 1289 ierr = ISSort(is[i]);CHKERRQ(ierr); 1290 start += count[i]; 1291 } 1292 1293 ierr = PetscFree(count); 1294 ierr = PetscFree(indices); 1295 ierr = ISDestroy(&isnumb);CHKERRQ(ierr); 1296 ierr = ISDestroy(&ispart);CHKERRQ(ierr); 1297 1298 } 1299 1300 PetscFunctionReturn(0); 1301 } 1302 1303 #undef __FUNCT__ 1304 #define __FUNCT__ "PCASMDestroySubdomains" 1305 /*@C 1306 PCASMDestroySubdomains - Destroys the index sets created with 1307 PCASMCreateSubdomains(). Should be called after setting subdomains 1308 with PCASMSetLocalSubdomains(). 1309 1310 Collective 1311 1312 Input Parameters: 1313 + n - the number of index sets 1314 . is - the array of index sets 1315 - is_local - the array of local index sets, can be PETSC_NULL 1316 1317 Level: advanced 1318 1319 .keywords: PC, ASM, additive Schwarz, create, subdomains, unstructured grid 1320 1321 .seealso: PCASMCreateSubdomains(), PCASMSetLocalSubdomains() 1322 @*/ 1323 PetscErrorCode PCASMDestroySubdomains(PetscInt n, IS is[], IS is_local[]) 1324 { 1325 PetscInt i; 1326 PetscErrorCode ierr; 1327 PetscFunctionBegin; 1328 if (n <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"n must be > 0: n = %D",n); 1329 PetscValidPointer(is,2); 1330 for (i=0; i<n; i++) { ierr = ISDestroy(&is[i]);CHKERRQ(ierr); } 1331 ierr = PetscFree(is);CHKERRQ(ierr); 1332 if (is_local) { 1333 PetscValidPointer(is_local,3); 1334 for (i=0; i<n; i++) { ierr = ISDestroy(&is_local[i]);CHKERRQ(ierr); } 1335 ierr = PetscFree(is_local);CHKERRQ(ierr); 1336 } 1337 PetscFunctionReturn(0); 1338 } 1339 1340 #undef __FUNCT__ 1341 #define __FUNCT__ "PCASMCreateSubdomains2D" 1342 /*@ 1343 PCASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz 1344 preconditioner for a two-dimensional problem on a regular grid. 1345 1346 Not Collective 1347 1348 Input Parameters: 1349 + m, n - the number of mesh points in the x and y directions 1350 . M, N - the number of subdomains in the x and y directions 1351 . dof - degrees of freedom per node 1352 - overlap - overlap in mesh lines 1353 1354 Output Parameters: 1355 + Nsub - the number of subdomains created 1356 . is - array of index sets defining overlapping (if overlap > 0) subdomains 1357 - is_local - array of index sets defining non-overlapping subdomains 1358 1359 Note: 1360 Presently PCAMSCreateSubdomains2d() is valid only for sequential 1361 preconditioners. More general related routines are 1362 PCASMSetTotalSubdomains() and PCASMSetLocalSubdomains(). 1363 1364 Level: advanced 1365 1366 .keywords: PC, ASM, additive Schwarz, create, subdomains, 2D, regular grid 1367 1368 .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubKSP(), 1369 PCASMSetOverlap() 1370 @*/ 1371 PetscErrorCode PCASMCreateSubdomains2D(PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt dof,PetscInt overlap,PetscInt *Nsub,IS **is,IS **is_local) 1372 { 1373 PetscInt i,j,height,width,ystart,xstart,yleft,yright,xleft,xright,loc_outer; 1374 PetscErrorCode ierr; 1375 PetscInt nidx,*idx,loc,ii,jj,count; 1376 1377 PetscFunctionBegin; 1378 if (dof != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP," "); 1379 1380 *Nsub = N*M; 1381 ierr = PetscMalloc((*Nsub)*sizeof(IS*),is);CHKERRQ(ierr); 1382 ierr = PetscMalloc((*Nsub)*sizeof(IS*),is_local);CHKERRQ(ierr); 1383 ystart = 0; 1384 loc_outer = 0; 1385 for (i=0; i<N; i++) { 1386 height = n/N + ((n % N) > i); /* height of subdomain */ 1387 if (height < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many N subdomains for mesh dimension n"); 1388 yleft = ystart - overlap; if (yleft < 0) yleft = 0; 1389 yright = ystart + height + overlap; if (yright > n) yright = n; 1390 xstart = 0; 1391 for (j=0; j<M; j++) { 1392 width = m/M + ((m % M) > j); /* width of subdomain */ 1393 if (width < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many M subdomains for mesh dimension m"); 1394 xleft = xstart - overlap; if (xleft < 0) xleft = 0; 1395 xright = xstart + width + overlap; if (xright > m) xright = m; 1396 nidx = (xright - xleft)*(yright - yleft); 1397 ierr = PetscMalloc(nidx*sizeof(PetscInt),&idx);CHKERRQ(ierr); 1398 loc = 0; 1399 for (ii=yleft; ii<yright; ii++) { 1400 count = m*ii + xleft; 1401 for (jj=xleft; jj<xright; jj++) { 1402 idx[loc++] = count++; 1403 } 1404 } 1405 ierr = ISCreateGeneral(PETSC_COMM_SELF,nidx,idx,PETSC_COPY_VALUES,(*is)+loc_outer);CHKERRQ(ierr); 1406 if (overlap == 0) { 1407 ierr = PetscObjectReference((PetscObject)(*is)[loc_outer]);CHKERRQ(ierr); 1408 (*is_local)[loc_outer] = (*is)[loc_outer]; 1409 } else { 1410 for (loc=0,ii=ystart; ii<ystart+height; ii++) { 1411 for (jj=xstart; jj<xstart+width; jj++) { 1412 idx[loc++] = m*ii + jj; 1413 } 1414 } 1415 ierr = ISCreateGeneral(PETSC_COMM_SELF,loc,idx,PETSC_COPY_VALUES,*is_local+loc_outer);CHKERRQ(ierr); 1416 } 1417 ierr = PetscFree(idx);CHKERRQ(ierr); 1418 xstart += width; 1419 loc_outer++; 1420 } 1421 ystart += height; 1422 } 1423 for (i=0; i<*Nsub; i++) { ierr = ISSort((*is)[i]);CHKERRQ(ierr); } 1424 PetscFunctionReturn(0); 1425 } 1426 1427 #undef __FUNCT__ 1428 #define __FUNCT__ "PCASMGetLocalSubdomains" 1429 /*@C 1430 PCASMGetLocalSubdomains - Gets the local subdomains (for this processor 1431 only) for the additive Schwarz preconditioner. 1432 1433 Not Collective 1434 1435 Input Parameter: 1436 . pc - the preconditioner context 1437 1438 Output Parameters: 1439 + n - the number of subdomains for this processor (default value = 1) 1440 . is - the index sets that define the subdomains for this processor 1441 - is_local - the index sets that define the local part of the subdomains for this processor (can be PETSC_NULL) 1442 1443 1444 Notes: 1445 The IS numbering is in the parallel, global numbering of the vector. 1446 1447 Level: advanced 1448 1449 .keywords: PC, ASM, set, local, subdomains, additive Schwarz 1450 1451 .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(), 1452 PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubmatrices() 1453 @*/ 1454 PetscErrorCode PCASMGetLocalSubdomains(PC pc,PetscInt *n,IS *is[],IS *is_local[]) 1455 { 1456 PC_ASM *osm; 1457 PetscErrorCode ierr; 1458 PetscBool match; 1459 1460 PetscFunctionBegin; 1461 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1462 PetscValidIntPointer(n,2); 1463 if (is) PetscValidPointer(is,3); 1464 ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr); 1465 if (!match) { 1466 if (n) *n = 0; 1467 if (is) *is = PETSC_NULL; 1468 } else { 1469 osm = (PC_ASM*)pc->data; 1470 if (n) *n = osm->n_local_true; 1471 if (is) *is = osm->is; 1472 if (is_local) *is_local = osm->is_local; 1473 } 1474 PetscFunctionReturn(0); 1475 } 1476 1477 #undef __FUNCT__ 1478 #define __FUNCT__ "PCASMGetLocalSubmatrices" 1479 /*@C 1480 PCASMGetLocalSubmatrices - Gets the local submatrices (for this processor 1481 only) for the additive Schwarz preconditioner. 1482 1483 Not Collective 1484 1485 Input Parameter: 1486 . pc - the preconditioner context 1487 1488 Output Parameters: 1489 + n - the number of matrices for this processor (default value = 1) 1490 - mat - the matrices 1491 1492 1493 Level: advanced 1494 1495 Notes: Call after PCSetUp() (or KSPSetUp()) but before PCApply() (or KSPApply()) and before PCSetUpOnBlocks()) 1496 1497 Usually one would use PCSetModifySubmatrices() to change the submatrices in building the preconditioner. 1498 1499 .keywords: PC, ASM, set, local, subdomains, additive Schwarz, block Jacobi 1500 1501 .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(), 1502 PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains(), PCSetModifySubmatrices() 1503 @*/ 1504 PetscErrorCode PCASMGetLocalSubmatrices(PC pc,PetscInt *n,Mat *mat[]) 1505 { 1506 PC_ASM *osm; 1507 PetscErrorCode ierr; 1508 PetscBool match; 1509 1510 PetscFunctionBegin; 1511 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1512 PetscValidIntPointer(n,2); 1513 if (mat) PetscValidPointer(mat,3); 1514 if (!pc->setupcalled) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUP() or PCSetUp()."); 1515 ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr); 1516 if (!match) { 1517 if (n) *n = 0; 1518 if (mat) *mat = PETSC_NULL; 1519 } else { 1520 osm = (PC_ASM*)pc->data; 1521 if (n) *n = osm->n_local_true; 1522 if (mat) *mat = osm->pmat; 1523 } 1524 PetscFunctionReturn(0); 1525 } 1526