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