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