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