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