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