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