1 2 /* 3 Creates a matrix class for using the Neumann-Neumann type preconditioners. 4 This stores the matrices in globally unassembled form. Each processor 5 assembles only its local Neumann problem and the parallel matrix vector 6 product is handled "implicitly". 7 8 We provide: 9 MatMult() 10 11 Currently this allows for only one subdomain per processor. 12 13 */ 14 15 #include <../src/mat/impls/is/matis.h> /*I "petscmat.h" I*/ 16 #include <petscsf.h> 17 18 #undef __FUNCT__ 19 #define __FUNCT__ "MatISSetPreallocation" 20 /* 21 MatISSetPreallocation - Preallocates memory for a MATIS parallel matrix. 22 23 Collective on MPI_Comm 24 25 Input Parameters: 26 + B - the matrix 27 . d_nz - number of nonzeros per row in DIAGONAL portion of local submatrix 28 (same value is used for all local rows) 29 . d_nnz - array containing the number of nonzeros in the various rows of the 30 DIAGONAL portion of the local submatrix (possibly different for each row) 31 or NULL, if d_nz is used to specify the nonzero structure. 32 The size of this array is equal to the number of local rows, i.e 'm'. 33 For matrices that will be factored, you must leave room for (and set) 34 the diagonal entry even if it is zero. 35 . o_nz - number of nonzeros per row in the OFF-DIAGONAL portion of local 36 submatrix (same value is used for all local rows). 37 - o_nnz - array containing the number of nonzeros in the various rows of the 38 OFF-DIAGONAL portion of the local submatrix (possibly different for 39 each row) or NULL, if o_nz is used to specify the nonzero 40 structure. The size of this array is equal to the number 41 of local rows, i.e 'm'. 42 43 If the *_nnz parameter is given then the *_nz parameter is ignored 44 45 Level: intermediate 46 47 Notes: This function has the same interface as the MPIAIJ preallocation routine in order to simplify the transition 48 from the asssembled format to the unassembled one. It overestimates the preallocation of MATIS local 49 matrices; for exact preallocation, the user should set the preallocation directly on local matrix objects. 50 51 .keywords: matrix 52 53 .seealso: MatCreate(), MatCreateIS(), MatMPIAIJSetPreallocation(), MatISGetLocalMat() 54 @*/ 55 PetscErrorCode MatISSetPreallocation(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[]) 56 { 57 PetscErrorCode ierr; 58 59 PetscFunctionBegin; 60 PetscValidHeaderSpecific(B,MAT_CLASSID,1); 61 PetscValidType(B,1); 62 ierr = PetscTryMethod(B,"MatISSetPreallocation_C",(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(B,d_nz,d_nnz,o_nz,o_nnz));CHKERRQ(ierr); 63 PetscFunctionReturn(0); 64 } 65 66 #undef __FUNCT__ 67 #define __FUNCT__ "MatISSetPreallocation_IS" 68 PetscErrorCode MatISSetPreallocation_IS(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[]) 69 { 70 Mat_IS *matis = (Mat_IS*)(B->data); 71 PetscSF sf; 72 PetscInt bs,i,nroots,*rootdata,nleaves,*leafdata,nlocalcols; 73 const PetscInt *gidxs; 74 PetscErrorCode ierr; 75 76 PetscFunctionBegin; 77 if (!matis->A) { 78 SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"You should first call MatSetLocalToGlobalMapping"); 79 } 80 ierr = MatGetLocalSize(B,&nroots,NULL);CHKERRQ(ierr); 81 ierr = MatGetSize(matis->A,&nleaves,&nlocalcols);CHKERRQ(ierr); 82 ierr = MatGetBlockSize(matis->A,&bs);CHKERRQ(ierr); 83 ierr = PetscCalloc2(nroots,&rootdata,nleaves,&leafdata);CHKERRQ(ierr); 84 ierr = PetscSFCreate(PetscObjectComm((PetscObject)B),&sf);CHKERRQ(ierr); 85 ierr = PetscSFSetFromOptions(sf);CHKERRQ(ierr); 86 ierr = ISLocalToGlobalMappingGetIndices(matis->mapping,&gidxs);CHKERRQ(ierr); 87 ierr = PetscSFSetGraphLayout(sf,B->rmap,nleaves,NULL,PETSC_COPY_VALUES,gidxs);CHKERRQ(ierr); 88 ierr = ISLocalToGlobalMappingRestoreIndices(matis->mapping,&gidxs);CHKERRQ(ierr); 89 if (!d_nnz) { 90 for (i=0;i<nroots;i++) rootdata[i] += d_nz; 91 } else { 92 for (i=0;i<nroots;i++) rootdata[i] += d_nnz[i]; 93 } 94 if (!o_nnz) { 95 for (i=0;i<nroots;i++) rootdata[i] += o_nz; 96 } else { 97 for (i=0;i<nroots;i++) rootdata[i] += o_nnz[i]; 98 } 99 ierr = PetscSFBcastBegin(sf,MPIU_INT,rootdata,leafdata);CHKERRQ(ierr); 100 ierr = PetscSFBcastEnd(sf,MPIU_INT,rootdata,leafdata);CHKERRQ(ierr); 101 for (i=0;i<nleaves;i++) { 102 leafdata[i] = PetscMin(leafdata[i],nlocalcols); 103 } 104 ierr = MatSeqAIJSetPreallocation(matis->A,0,leafdata);CHKERRQ(ierr); 105 for (i=0;i<nleaves/bs;i++) { 106 leafdata[i] = leafdata[i*bs]/bs; 107 } 108 ierr = MatSeqBAIJSetPreallocation(matis->A,bs,0,leafdata);CHKERRQ(ierr); 109 for (i=0;i<nleaves/bs;i++) { 110 leafdata[i] = leafdata[i]-i; 111 } 112 ierr = MatSeqSBAIJSetPreallocation(matis->A,bs,0,leafdata);CHKERRQ(ierr); 113 ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 114 ierr = PetscFree2(rootdata,leafdata);CHKERRQ(ierr); 115 PetscFunctionReturn(0); 116 } 117 118 #undef __FUNCT__ 119 #define __FUNCT__ "MatISGetMPIXAIJ_IS" 120 PetscErrorCode MatISGetMPIXAIJ_IS(Mat mat, MatReuse reuse, Mat *M) 121 { 122 Mat_IS *matis = (Mat_IS*)(mat->data); 123 /* info on mat */ 124 /* ISLocalToGlobalMapping rmapping,cmapping; */ 125 PetscInt bs,rows,cols; 126 PetscInt lrows,lcols; 127 PetscInt local_rows,local_cols; 128 PetscBool isdense,issbaij,issbaij_red; 129 PetscMPIInt nsubdomains; 130 /* values insertion */ 131 PetscScalar *array; 132 PetscInt *local_indices,*global_indices; 133 /* work */ 134 PetscInt i,j,index_row; 135 PetscErrorCode ierr; 136 137 PetscFunctionBegin; 138 /* MISSING CHECKS 139 - rectangular case not covered (it is not allowed by MATIS) 140 */ 141 /* get info from mat */ 142 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&nsubdomains);CHKERRQ(ierr); 143 if (nsubdomains == 1) { 144 if (reuse == MAT_INITIAL_MATRIX) { 145 ierr = MatDuplicate(matis->A,MAT_COPY_VALUES,&(*M));CHKERRQ(ierr); 146 } else { 147 ierr = MatCopy(matis->A,*M,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 148 } 149 PetscFunctionReturn(0); 150 } 151 /* ierr = MatGetLocalToGlobalMapping(mat,&rmapping,&cmapping);CHKERRQ(ierr); */ 152 ierr = MatGetSize(mat,&rows,&cols);CHKERRQ(ierr); 153 ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 154 ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr); 155 ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr); 156 ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr); 157 158 /* work */ 159 ierr = PetscMalloc1(local_rows,&local_indices);CHKERRQ(ierr); 160 for (i=0;i<local_rows;i++) local_indices[i]=i; 161 /* map indices of local mat to global values */ 162 ierr = PetscMalloc1(PetscMax(local_cols,local_rows),&global_indices);CHKERRQ(ierr); 163 /* ierr = ISLocalToGlobalMappingApply(rmapping,local_rows,local_indices,global_indices);CHKERRQ(ierr); */ 164 ierr = ISLocalToGlobalMappingApply(matis->mapping,local_rows,local_indices,global_indices);CHKERRQ(ierr); 165 166 if (issbaij) { 167 ierr = MatGetRowUpperTriangular(matis->A);CHKERRQ(ierr); 168 } 169 170 if (reuse == MAT_INITIAL_MATRIX) { 171 Mat new_mat; 172 MatType new_mat_type; 173 Vec vec_dnz,vec_onz; 174 PetscScalar *my_dnz,*my_onz; 175 PetscInt *dnz,*onz,*mat_ranges,*row_ownership; 176 PetscInt index_col,owner; 177 178 /* determining new matrix type */ 179 ierr = MPI_Allreduce(&issbaij,&issbaij_red,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 180 if (issbaij_red) { 181 new_mat_type = MATSBAIJ; 182 } else { 183 if (bs>1) { 184 new_mat_type = MATBAIJ; 185 } else { 186 new_mat_type = MATAIJ; 187 } 188 } 189 190 ierr = MatCreate(PetscObjectComm((PetscObject)mat),&new_mat);CHKERRQ(ierr); 191 ierr = MatSetSizes(new_mat,PETSC_DECIDE,PETSC_DECIDE,rows,cols);CHKERRQ(ierr); 192 ierr = MatSetBlockSize(new_mat,bs);CHKERRQ(ierr); 193 ierr = MatSetType(new_mat,new_mat_type);CHKERRQ(ierr); 194 ierr = MatSetUp(new_mat);CHKERRQ(ierr); 195 ierr = MatGetLocalSize(new_mat,&lrows,&lcols);CHKERRQ(ierr); 196 197 /* 198 preallocation 199 */ 200 201 ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)new_mat),lrows,lcols,dnz,onz);CHKERRQ(ierr); 202 /* 203 Some vectors are needed to sum up properly on shared interface dofs. 204 Preallocation macros cannot do the job. 205 Note that preallocation is not exact, since it overestimates nonzeros 206 */ 207 ierr = MatCreateVecs(new_mat,NULL,&vec_dnz);CHKERRQ(ierr); 208 /* ierr = VecSetLocalToGlobalMapping(vec_dnz,rmapping);CHKERRQ(ierr); */ 209 ierr = VecSetLocalToGlobalMapping(vec_dnz,matis->mapping);CHKERRQ(ierr); 210 ierr = VecDuplicate(vec_dnz,&vec_onz);CHKERRQ(ierr); 211 /* All processes need to compute entire row ownership */ 212 ierr = PetscMalloc1(rows,&row_ownership);CHKERRQ(ierr); 213 ierr = MatGetOwnershipRanges(new_mat,(const PetscInt**)&mat_ranges);CHKERRQ(ierr); 214 for (i=0;i<nsubdomains;i++) { 215 for (j=mat_ranges[i];j<mat_ranges[i+1];j++) { 216 row_ownership[j]=i; 217 } 218 } 219 220 /* 221 my_dnz and my_onz contains exact contribution to preallocation from each local mat 222 then, they will be summed up properly. This way, preallocation is always sufficient 223 */ 224 ierr = PetscMalloc1(local_rows,&my_dnz);CHKERRQ(ierr); 225 ierr = PetscMalloc1(local_rows,&my_onz);CHKERRQ(ierr); 226 ierr = PetscMemzero(my_dnz,local_rows*sizeof(*my_dnz));CHKERRQ(ierr); 227 ierr = PetscMemzero(my_onz,local_rows*sizeof(*my_onz));CHKERRQ(ierr); 228 /* preallocation as a MATAIJ */ 229 if (isdense) { /* special case for dense local matrices */ 230 for (i=0;i<local_rows;i++) { 231 index_row = global_indices[i]; 232 for (j=i;j<local_rows;j++) { 233 owner = row_ownership[index_row]; 234 index_col = global_indices[j]; 235 if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */ 236 my_dnz[i] += 1.0; 237 } else { /* offdiag block */ 238 my_onz[i] += 1.0; 239 } 240 /* same as before, interchanging rows and cols */ 241 if (i != j) { 242 owner = row_ownership[index_col]; 243 if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) { 244 my_dnz[j] += 1.0; 245 } else { 246 my_onz[j] += 1.0; 247 } 248 } 249 } 250 } 251 } else { 252 for (i=0;i<local_rows;i++) { 253 PetscInt ncols; 254 const PetscInt *cols; 255 index_row = global_indices[i]; 256 ierr = MatGetRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr); 257 for (j=0;j<ncols;j++) { 258 owner = row_ownership[index_row]; 259 index_col = global_indices[cols[j]]; 260 if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */ 261 my_dnz[i] += 1.0; 262 } else { /* offdiag block */ 263 my_onz[i] += 1.0; 264 } 265 /* same as before, interchanging rows and cols */ 266 if (issbaij) { 267 owner = row_ownership[index_col]; 268 if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) { 269 my_dnz[j] += 1.0; 270 } else { 271 my_onz[j] += 1.0; 272 } 273 } 274 } 275 ierr = MatRestoreRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr); 276 } 277 } 278 ierr = VecSet(vec_dnz,0.0);CHKERRQ(ierr); 279 ierr = VecSet(vec_onz,0.0);CHKERRQ(ierr); 280 if (local_rows) { /* multilevel guard */ 281 ierr = VecSetValuesLocal(vec_dnz,local_rows,local_indices,my_dnz,ADD_VALUES);CHKERRQ(ierr); 282 ierr = VecSetValuesLocal(vec_onz,local_rows,local_indices,my_onz,ADD_VALUES);CHKERRQ(ierr); 283 } 284 ierr = VecAssemblyBegin(vec_dnz);CHKERRQ(ierr); 285 ierr = VecAssemblyBegin(vec_onz);CHKERRQ(ierr); 286 ierr = VecAssemblyEnd(vec_dnz);CHKERRQ(ierr); 287 ierr = VecAssemblyEnd(vec_onz);CHKERRQ(ierr); 288 ierr = PetscFree(my_dnz);CHKERRQ(ierr); 289 ierr = PetscFree(my_onz);CHKERRQ(ierr); 290 ierr = PetscFree(row_ownership);CHKERRQ(ierr); 291 292 /* set computed preallocation in dnz and onz */ 293 ierr = VecGetArray(vec_dnz,&array);CHKERRQ(ierr); 294 for (i=0; i<lrows; i++) dnz[i] = (PetscInt)PetscRealPart(array[i]); 295 ierr = VecRestoreArray(vec_dnz,&array);CHKERRQ(ierr); 296 ierr = VecGetArray(vec_onz,&array);CHKERRQ(ierr); 297 for (i=0;i<lrows;i++) onz[i] = (PetscInt)PetscRealPart(array[i]); 298 ierr = VecRestoreArray(vec_onz,&array);CHKERRQ(ierr); 299 ierr = VecDestroy(&vec_dnz);CHKERRQ(ierr); 300 ierr = VecDestroy(&vec_onz);CHKERRQ(ierr); 301 302 /* Resize preallocation if overestimated */ 303 for (i=0;i<lrows;i++) { 304 dnz[i] = PetscMin(dnz[i],lcols); 305 onz[i] = PetscMin(onz[i],cols-lcols); 306 } 307 /* set preallocation */ 308 ierr = MatMPIAIJSetPreallocation(new_mat,0,dnz,0,onz);CHKERRQ(ierr); 309 for (i=0;i<lrows/bs;i++) { 310 dnz[i] = dnz[i*bs]/bs; 311 onz[i] = onz[i*bs]/bs; 312 } 313 ierr = MatMPIBAIJSetPreallocation(new_mat,bs,0,dnz,0,onz);CHKERRQ(ierr); 314 for (i=0;i<lrows/bs;i++) { 315 dnz[i] = dnz[i]-i; 316 } 317 ierr = MatMPISBAIJSetPreallocation(new_mat,bs,0,dnz,0,onz);CHKERRQ(ierr); 318 ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 319 *M = new_mat; 320 } else { 321 PetscInt mbs,mrows,mcols; 322 /* some checks */ 323 ierr = MatGetBlockSize(*M,&mbs);CHKERRQ(ierr); 324 ierr = MatGetSize(*M,&mrows,&mcols);CHKERRQ(ierr); 325 if (mrows != rows) { 326 SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of rows (%d != %d)",rows,mrows); 327 } 328 if (mrows != rows) { 329 SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of cols (%d != %d)",cols,mcols); 330 } 331 if (mbs != bs) { 332 SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong block size (%d != %d)",bs,mbs); 333 } 334 ierr = MatZeroEntries(*M);CHKERRQ(ierr); 335 } 336 /* set local to global mappings */ 337 /* ierr = MatSetLocalToGlobalMapping(*M,rmapping,cmapping);CHKERRQ(ierr); */ 338 /* Set values */ 339 if (isdense) { /* special case for dense local matrices */ 340 ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr); 341 ierr = MatDenseGetArray(matis->A,&array);CHKERRQ(ierr); 342 ierr = MatSetValues(*M,local_rows,global_indices,local_cols,global_indices,array,ADD_VALUES);CHKERRQ(ierr); 343 ierr = MatDenseRestoreArray(matis->A,&array);CHKERRQ(ierr); 344 ierr = PetscFree(local_indices);CHKERRQ(ierr); 345 ierr = PetscFree(global_indices);CHKERRQ(ierr); 346 } else { /* very basic values insertion for all other matrix types */ 347 ierr = PetscFree(local_indices);CHKERRQ(ierr); 348 for (i=0;i<local_rows;i++) { 349 ierr = MatGetRow(matis->A,i,&j,(const PetscInt**)&local_indices,(const PetscScalar**)&array);CHKERRQ(ierr); 350 /* ierr = MatSetValuesLocal(*M,1,&i,j,local_indices,array,ADD_VALUES);CHKERRQ(ierr); */ 351 ierr = ISLocalToGlobalMappingApply(matis->mapping,j,local_indices,global_indices);CHKERRQ(ierr); 352 ierr = ISLocalToGlobalMappingApply(matis->mapping,1,&i,&index_row);CHKERRQ(ierr); 353 ierr = MatSetValues(*M,1,&index_row,j,global_indices,array,ADD_VALUES);CHKERRQ(ierr); 354 ierr = MatRestoreRow(matis->A,i,&j,(const PetscInt**)&local_indices,(const PetscScalar**)&array);CHKERRQ(ierr); 355 } 356 ierr = PetscFree(global_indices);CHKERRQ(ierr); 357 } 358 ierr = MatAssemblyBegin(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 359 ierr = MatAssemblyEnd(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 360 if (isdense) { 361 ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr); 362 } 363 if (issbaij) { 364 ierr = MatRestoreRowUpperTriangular(matis->A);CHKERRQ(ierr); 365 } 366 PetscFunctionReturn(0); 367 } 368 369 #undef __FUNCT__ 370 #define __FUNCT__ "MatISGetMPIXAIJ" 371 /*@ 372 MatISGetMPIXAIJ - Converts MATIS matrix into a parallel AIJ format 373 374 Input Parameter: 375 . mat - the matrix (should be of type MATIS) 376 . reuse - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX 377 378 Output Parameter: 379 . newmat - the matrix in AIJ format 380 381 Level: developer 382 383 Notes: 384 385 .seealso: MATIS 386 @*/ 387 PetscErrorCode MatISGetMPIXAIJ(Mat mat, MatReuse reuse, Mat *newmat) 388 { 389 PetscErrorCode ierr; 390 391 PetscFunctionBegin; 392 PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 393 PetscValidLogicalCollectiveEnum(mat,reuse,2); 394 PetscValidPointer(newmat,3); 395 if (reuse != MAT_INITIAL_MATRIX) { 396 PetscValidHeaderSpecific(*newmat,MAT_CLASSID,3); 397 PetscCheckSameComm(mat,1,*newmat,3); 398 } 399 ierr = PetscUseMethod(mat,"MatISGetMPIXAIJ_C",(Mat,MatReuse,Mat*),(mat,reuse,newmat));CHKERRQ(ierr); 400 PetscFunctionReturn(0); 401 } 402 403 #undef __FUNCT__ 404 #define __FUNCT__ "MatDuplicate_IS" 405 PetscErrorCode MatDuplicate_IS(Mat mat,MatDuplicateOption op,Mat *newmat) 406 { 407 PetscErrorCode ierr; 408 Mat_IS *matis = (Mat_IS*)(mat->data); 409 PetscInt bs,m,n,M,N; 410 Mat B,localmat; 411 412 PetscFunctionBegin; 413 ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 414 ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 415 ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 416 ierr = MatCreateIS(PetscObjectComm((PetscObject)mat),bs,m,n,M,N,matis->mapping,&B);CHKERRQ(ierr); 417 ierr = MatDuplicate(matis->A,op,&localmat);CHKERRQ(ierr); 418 ierr = MatISSetLocalMat(B,localmat);CHKERRQ(ierr); 419 ierr = MatDestroy(&localmat);CHKERRQ(ierr); 420 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 421 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 422 *newmat = B; 423 PetscFunctionReturn(0); 424 } 425 426 #undef __FUNCT__ 427 #define __FUNCT__ "MatIsHermitian_IS" 428 PetscErrorCode MatIsHermitian_IS(Mat A,PetscReal tol,PetscBool *flg) 429 { 430 PetscErrorCode ierr; 431 Mat_IS *matis = (Mat_IS*)A->data; 432 PetscBool local_sym; 433 434 PetscFunctionBegin; 435 ierr = MatIsHermitian(matis->A,tol,&local_sym);CHKERRQ(ierr); 436 ierr = MPI_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 437 PetscFunctionReturn(0); 438 } 439 440 #undef __FUNCT__ 441 #define __FUNCT__ "MatIsSymmetric_IS" 442 PetscErrorCode MatIsSymmetric_IS(Mat A,PetscReal tol,PetscBool *flg) 443 { 444 PetscErrorCode ierr; 445 Mat_IS *matis = (Mat_IS*)A->data; 446 PetscBool local_sym; 447 448 PetscFunctionBegin; 449 ierr = MatIsSymmetric(matis->A,tol,&local_sym);CHKERRQ(ierr); 450 ierr = MPI_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 451 PetscFunctionReturn(0); 452 } 453 454 #undef __FUNCT__ 455 #define __FUNCT__ "MatDestroy_IS" 456 PetscErrorCode MatDestroy_IS(Mat A) 457 { 458 PetscErrorCode ierr; 459 Mat_IS *b = (Mat_IS*)A->data; 460 461 PetscFunctionBegin; 462 ierr = MatDestroy(&b->A);CHKERRQ(ierr); 463 ierr = VecScatterDestroy(&b->ctx);CHKERRQ(ierr); 464 ierr = VecDestroy(&b->x);CHKERRQ(ierr); 465 ierr = VecDestroy(&b->y);CHKERRQ(ierr); 466 ierr = ISLocalToGlobalMappingDestroy(&b->mapping);CHKERRQ(ierr); 467 ierr = PetscFree(A->data);CHKERRQ(ierr); 468 ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr); 469 ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",NULL);CHKERRQ(ierr); 470 ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",NULL);CHKERRQ(ierr); 471 ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",NULL);CHKERRQ(ierr); 472 ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",NULL);CHKERRQ(ierr); 473 PetscFunctionReturn(0); 474 } 475 476 #undef __FUNCT__ 477 #define __FUNCT__ "MatMult_IS" 478 PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y) 479 { 480 PetscErrorCode ierr; 481 Mat_IS *is = (Mat_IS*)A->data; 482 PetscScalar zero = 0.0; 483 484 PetscFunctionBegin; 485 /* scatter the global vector x into the local work vector */ 486 ierr = VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 487 ierr = VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 488 489 /* multiply the local matrix */ 490 ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr); 491 492 /* scatter product back into global memory */ 493 ierr = VecSet(y,zero);CHKERRQ(ierr); 494 ierr = VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 495 ierr = VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 496 PetscFunctionReturn(0); 497 } 498 499 #undef __FUNCT__ 500 #define __FUNCT__ "MatMultAdd_IS" 501 PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 502 { 503 Vec temp_vec; 504 PetscErrorCode ierr; 505 506 PetscFunctionBegin; /* v3 = v2 + A * v1.*/ 507 if (v3 != v2) { 508 ierr = MatMult(A,v1,v3);CHKERRQ(ierr); 509 ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 510 } else { 511 ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 512 ierr = MatMult(A,v1,temp_vec);CHKERRQ(ierr); 513 ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 514 ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 515 ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 516 } 517 PetscFunctionReturn(0); 518 } 519 520 #undef __FUNCT__ 521 #define __FUNCT__ "MatMultTranspose_IS" 522 PetscErrorCode MatMultTranspose_IS(Mat A,Vec x,Vec y) 523 { 524 Mat_IS *is = (Mat_IS*)A->data; 525 PetscErrorCode ierr; 526 527 PetscFunctionBegin; /* y = A' * x */ 528 /* scatter the global vector x into the local work vector */ 529 ierr = VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 530 ierr = VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 531 532 /* multiply the local matrix */ 533 ierr = MatMultTranspose(is->A,is->x,is->y);CHKERRQ(ierr); 534 535 /* scatter product back into global vector */ 536 ierr = VecSet(y,0);CHKERRQ(ierr); 537 ierr = VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 538 ierr = VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 539 PetscFunctionReturn(0); 540 } 541 542 #undef __FUNCT__ 543 #define __FUNCT__ "MatMultTransposeAdd_IS" 544 PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 545 { 546 Vec temp_vec; 547 PetscErrorCode ierr; 548 549 PetscFunctionBegin; /* v3 = v2 + A' * v1.*/ 550 if (v3 != v2) { 551 ierr = MatMultTranspose(A,v1,v3);CHKERRQ(ierr); 552 ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 553 } else { 554 ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 555 ierr = MatMultTranspose(A,v1,temp_vec);CHKERRQ(ierr); 556 ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 557 ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 558 ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 559 } 560 PetscFunctionReturn(0); 561 } 562 563 #undef __FUNCT__ 564 #define __FUNCT__ "MatView_IS" 565 PetscErrorCode MatView_IS(Mat A,PetscViewer viewer) 566 { 567 Mat_IS *a = (Mat_IS*)A->data; 568 PetscErrorCode ierr; 569 PetscViewer sviewer; 570 571 PetscFunctionBegin; 572 ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr); 573 ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 574 ierr = MatView(a->A,sviewer);CHKERRQ(ierr); 575 ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 576 PetscFunctionReturn(0); 577 } 578 579 #undef __FUNCT__ 580 #define __FUNCT__ "MatSetLocalToGlobalMapping_IS" 581 PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping) 582 { 583 PetscErrorCode ierr; 584 PetscInt n,bs; 585 Mat_IS *is = (Mat_IS*)A->data; 586 IS from,to; 587 Vec global; 588 589 PetscFunctionBegin; 590 PetscCheckSameComm(A,1,rmapping,2); 591 if (rmapping != cmapping) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"MATIS requires the row and column mappings to be identical"); 592 if (is->mapping) { /* Currenly destroys the objects that will be created by this routine. Is there anything else that should be checked? */ 593 ierr = ISLocalToGlobalMappingDestroy(&is->mapping);CHKERRQ(ierr); 594 ierr = VecDestroy(&is->x);CHKERRQ(ierr); 595 ierr = VecDestroy(&is->y);CHKERRQ(ierr); 596 ierr = VecScatterDestroy(&is->ctx);CHKERRQ(ierr); 597 ierr = MatDestroy(&is->A);CHKERRQ(ierr); 598 } 599 ierr = PetscObjectReference((PetscObject)rmapping);CHKERRQ(ierr); 600 ierr = ISLocalToGlobalMappingDestroy(&is->mapping);CHKERRQ(ierr); 601 is->mapping = rmapping; 602 /* 603 ierr = PetscLayoutSetISLocalToGlobalMapping(A->rmap,rmapping);CHKERRQ(ierr); 604 ierr = PetscLayoutSetISLocalToGlobalMapping(A->cmap,cmapping);CHKERRQ(ierr); 605 */ 606 607 /* Create the local matrix A */ 608 ierr = ISLocalToGlobalMappingGetSize(rmapping,&n);CHKERRQ(ierr); 609 ierr = ISLocalToGlobalMappingGetBlockSize(rmapping,&bs);CHKERRQ(ierr); 610 ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr); 611 if (bs > 1) { 612 ierr = MatSetType(is->A,MATSEQBAIJ);CHKERRQ(ierr); 613 } else { 614 ierr = MatSetType(is->A,MATSEQAIJ);CHKERRQ(ierr); 615 } 616 ierr = MatSetSizes(is->A,n,n,n,n);CHKERRQ(ierr); 617 ierr = MatSetBlockSize(is->A,bs);CHKERRQ(ierr); 618 ierr = MatSetOptionsPrefix(is->A,((PetscObject)A)->prefix);CHKERRQ(ierr); 619 ierr = MatAppendOptionsPrefix(is->A,"is_");CHKERRQ(ierr); 620 ierr = MatSetFromOptions(is->A);CHKERRQ(ierr); 621 622 /* Create the local work vectors */ 623 ierr = VecCreate(PETSC_COMM_SELF,&is->x);CHKERRQ(ierr); 624 ierr = VecSetBlockSize(is->x,bs);CHKERRQ(ierr); 625 ierr = VecSetSizes(is->x,n,n);CHKERRQ(ierr); 626 ierr = VecSetOptionsPrefix(is->x,((PetscObject)A)->prefix);CHKERRQ(ierr); 627 ierr = VecAppendOptionsPrefix(is->x,"is_");CHKERRQ(ierr); 628 ierr = VecSetFromOptions(is->x);CHKERRQ(ierr); 629 ierr = VecDuplicate(is->x,&is->y);CHKERRQ(ierr); 630 631 /* setup the global to local scatter */ 632 ierr = ISCreateStride(PETSC_COMM_SELF,n,0,1,&to);CHKERRQ(ierr); 633 ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr); 634 ierr = MatCreateVecs(A,&global,NULL);CHKERRQ(ierr); 635 ierr = VecScatterCreate(global,from,is->x,to,&is->ctx);CHKERRQ(ierr); 636 ierr = VecDestroy(&global);CHKERRQ(ierr); 637 ierr = ISDestroy(&to);CHKERRQ(ierr); 638 ierr = ISDestroy(&from);CHKERRQ(ierr); 639 PetscFunctionReturn(0); 640 } 641 642 #undef __FUNCT__ 643 #define __FUNCT__ "MatSetValues_IS" 644 PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv) 645 { 646 Mat_IS *is = (Mat_IS*)mat->data; 647 PetscInt rows_l[2048],cols_l[2048]; 648 PetscErrorCode ierr; 649 650 PetscFunctionBegin; 651 #if defined(PETSC_USE_DEBUG) 652 if (m > 2048 || n > 2048) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of row/column indices must be <= 2048: they are %D %D",m,n); 653 #endif 654 ierr = ISG2LMapApply(is->mapping,m,rows,rows_l);CHKERRQ(ierr); 655 ierr = ISG2LMapApply(is->mapping,n,cols,cols_l);CHKERRQ(ierr); 656 ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 657 PetscFunctionReturn(0); 658 } 659 660 #undef ISG2LMapSetUp 661 #undef ISG2LMapApply 662 663 #undef __FUNCT__ 664 #define __FUNCT__ "MatSetValuesLocal_IS" 665 PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 666 { 667 PetscErrorCode ierr; 668 Mat_IS *is = (Mat_IS*)A->data; 669 670 PetscFunctionBegin; 671 ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 672 PetscFunctionReturn(0); 673 } 674 675 #undef __FUNCT__ 676 #define __FUNCT__ "MatSetValuesBlockedLocal_IS" 677 PetscErrorCode MatSetValuesBlockedLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 678 { 679 PetscErrorCode ierr; 680 Mat_IS *is = (Mat_IS*)A->data; 681 682 PetscFunctionBegin; 683 ierr = MatSetValuesBlocked(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 684 PetscFunctionReturn(0); 685 } 686 687 #undef __FUNCT__ 688 #define __FUNCT__ "MatZeroRows_IS" 689 PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 690 { 691 Mat_IS *is = (Mat_IS*)A->data; 692 PetscInt n_l = 0, *rows_l = NULL; 693 PetscErrorCode ierr; 694 695 PetscFunctionBegin; 696 if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support"); 697 if (n) { 698 ierr = PetscMalloc1(n,&rows_l);CHKERRQ(ierr); 699 ierr = ISGlobalToLocalMappingApply(is->mapping,IS_GTOLM_DROP,n,rows,&n_l,rows_l);CHKERRQ(ierr); 700 } 701 ierr = MatZeroRowsLocal(A,n_l,rows_l,diag,x,b);CHKERRQ(ierr); 702 ierr = PetscFree(rows_l);CHKERRQ(ierr); 703 PetscFunctionReturn(0); 704 } 705 706 #undef __FUNCT__ 707 #define __FUNCT__ "MatZeroRowsLocal_IS" 708 PetscErrorCode MatZeroRowsLocal_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 709 { 710 Mat_IS *is = (Mat_IS*)A->data; 711 PetscErrorCode ierr; 712 PetscInt i; 713 PetscScalar *array; 714 715 PetscFunctionBegin; 716 if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support"); 717 { 718 /* 719 Set up is->x as a "counting vector". This is in order to MatMult_IS 720 work properly in the interface nodes. 721 */ 722 Vec counter; 723 PetscScalar one=1.0, zero=0.0; 724 ierr = MatCreateVecs(A,&counter,NULL);CHKERRQ(ierr); 725 ierr = VecSet(counter,zero);CHKERRQ(ierr); 726 ierr = VecSet(is->x,one);CHKERRQ(ierr); 727 ierr = VecScatterBegin(is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 728 ierr = VecScatterEnd (is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 729 ierr = VecScatterBegin(is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 730 ierr = VecScatterEnd (is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 731 ierr = VecDestroy(&counter);CHKERRQ(ierr); 732 } 733 if (!n) { 734 is->pure_neumann = PETSC_TRUE; 735 } else { 736 is->pure_neumann = PETSC_FALSE; 737 738 ierr = VecGetArray(is->x,&array);CHKERRQ(ierr); 739 ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr); 740 for (i=0; i<n; i++) { 741 ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr); 742 } 743 ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 744 ierr = MatAssemblyEnd (is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 745 ierr = VecRestoreArray(is->x,&array);CHKERRQ(ierr); 746 } 747 PetscFunctionReturn(0); 748 } 749 750 #undef __FUNCT__ 751 #define __FUNCT__ "MatAssemblyBegin_IS" 752 PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type) 753 { 754 Mat_IS *is = (Mat_IS*)A->data; 755 PetscErrorCode ierr; 756 757 PetscFunctionBegin; 758 ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr); 759 PetscFunctionReturn(0); 760 } 761 762 #undef __FUNCT__ 763 #define __FUNCT__ "MatAssemblyEnd_IS" 764 PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type) 765 { 766 Mat_IS *is = (Mat_IS*)A->data; 767 PetscErrorCode ierr; 768 769 PetscFunctionBegin; 770 ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr); 771 PetscFunctionReturn(0); 772 } 773 774 #undef __FUNCT__ 775 #define __FUNCT__ "MatISGetLocalMat_IS" 776 PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local) 777 { 778 Mat_IS *is = (Mat_IS*)mat->data; 779 780 PetscFunctionBegin; 781 *local = is->A; 782 PetscFunctionReturn(0); 783 } 784 785 #undef __FUNCT__ 786 #define __FUNCT__ "MatISGetLocalMat" 787 /*@ 788 MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix. 789 790 Input Parameter: 791 . mat - the matrix 792 793 Output Parameter: 794 . local - the local matrix usually MATSEQAIJ 795 796 Level: advanced 797 798 Notes: 799 This can be called if you have precomputed the nonzero structure of the 800 matrix and want to provide it to the inner matrix object to improve the performance 801 of the MatSetValues() operation. 802 803 .seealso: MATIS 804 @*/ 805 PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local) 806 { 807 PetscErrorCode ierr; 808 809 PetscFunctionBegin; 810 PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 811 PetscValidPointer(local,2); 812 ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat*),(mat,local));CHKERRQ(ierr); 813 PetscFunctionReturn(0); 814 } 815 816 #undef __FUNCT__ 817 #define __FUNCT__ "MatISSetLocalMat_IS" 818 PetscErrorCode MatISSetLocalMat_IS(Mat mat,Mat local) 819 { 820 Mat_IS *is = (Mat_IS*)mat->data; 821 PetscInt nrows,ncols,orows,ocols; 822 PetscErrorCode ierr; 823 824 PetscFunctionBegin; 825 if (is->A) { 826 ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr); 827 ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr); 828 if (orows != nrows || ocols != ncols) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local MATIS matrix should be of size %dx%d (you passed a %dx%d matrix)\n",orows,ocols,nrows,ncols); 829 } 830 ierr = PetscObjectReference((PetscObject)local);CHKERRQ(ierr); 831 ierr = MatDestroy(&is->A);CHKERRQ(ierr); 832 is->A = local; 833 PetscFunctionReturn(0); 834 } 835 836 #undef __FUNCT__ 837 #define __FUNCT__ "MatISSetLocalMat" 838 /*@ 839 MatISSetLocalMat - Set the local matrix stored inside a MATIS matrix. 840 841 Input Parameter: 842 . mat - the matrix 843 . local - the local matrix usually MATSEQAIJ 844 845 Output Parameter: 846 847 Level: advanced 848 849 Notes: 850 This can be called if you have precomputed the local matrix and 851 want to provide it to the matrix object MATIS. 852 853 .seealso: MATIS 854 @*/ 855 PetscErrorCode MatISSetLocalMat(Mat mat,Mat local) 856 { 857 PetscErrorCode ierr; 858 859 PetscFunctionBegin; 860 PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 861 PetscValidHeaderSpecific(local,MAT_CLASSID,2); 862 ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr); 863 PetscFunctionReturn(0); 864 } 865 866 #undef __FUNCT__ 867 #define __FUNCT__ "MatZeroEntries_IS" 868 PetscErrorCode MatZeroEntries_IS(Mat A) 869 { 870 Mat_IS *a = (Mat_IS*)A->data; 871 PetscErrorCode ierr; 872 873 PetscFunctionBegin; 874 ierr = MatZeroEntries(a->A);CHKERRQ(ierr); 875 PetscFunctionReturn(0); 876 } 877 878 #undef __FUNCT__ 879 #define __FUNCT__ "MatScale_IS" 880 PetscErrorCode MatScale_IS(Mat A,PetscScalar a) 881 { 882 Mat_IS *is = (Mat_IS*)A->data; 883 PetscErrorCode ierr; 884 885 PetscFunctionBegin; 886 ierr = MatScale(is->A,a);CHKERRQ(ierr); 887 PetscFunctionReturn(0); 888 } 889 890 #undef __FUNCT__ 891 #define __FUNCT__ "MatGetDiagonal_IS" 892 PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v) 893 { 894 Mat_IS *is = (Mat_IS*)A->data; 895 PetscErrorCode ierr; 896 897 PetscFunctionBegin; 898 /* get diagonal of the local matrix */ 899 ierr = MatGetDiagonal(is->A,is->x);CHKERRQ(ierr); 900 901 /* scatter diagonal back into global vector */ 902 ierr = VecSet(v,0);CHKERRQ(ierr); 903 ierr = VecScatterBegin(is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 904 ierr = VecScatterEnd (is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 905 PetscFunctionReturn(0); 906 } 907 908 #undef __FUNCT__ 909 #define __FUNCT__ "MatSetOption_IS" 910 PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg) 911 { 912 Mat_IS *a = (Mat_IS*)A->data; 913 PetscErrorCode ierr; 914 915 PetscFunctionBegin; 916 ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr); 917 PetscFunctionReturn(0); 918 } 919 920 #undef __FUNCT__ 921 #define __FUNCT__ "MatCreateIS" 922 /*@ 923 MatCreateIS - Creates a "process" unassmembled matrix, it is assembled on each 924 process but not across processes. 925 926 Input Parameters: 927 + comm - MPI communicator that will share the matrix 928 . bs - local and global block size of the matrix 929 . m,n,M,N - local and/or global sizes of the left and right vector used in matrix vector products 930 - map - mapping that defines the global number for each local number 931 932 Output Parameter: 933 . A - the resulting matrix 934 935 Level: advanced 936 937 Notes: See MATIS for more details 938 m and n are NOT related to the size of the map, they are the size of the part of the vector owned 939 by that process. m + nghosts (or n + nghosts) is the length of map since map maps all local points 940 plus the ghost points to global indices. 941 942 .seealso: MATIS, MatSetLocalToGlobalMapping() 943 @*/ 944 PetscErrorCode MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping map,Mat *A) 945 { 946 PetscErrorCode ierr; 947 948 PetscFunctionBegin; 949 ierr = MatCreate(comm,A);CHKERRQ(ierr); 950 ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr); 951 ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 952 ierr = MatSetType(*A,MATIS);CHKERRQ(ierr); 953 ierr = MatSetUp(*A);CHKERRQ(ierr); 954 ierr = MatSetLocalToGlobalMapping(*A,map,map);CHKERRQ(ierr); 955 PetscFunctionReturn(0); 956 } 957 958 /*MC 959 MATIS - MATIS = "is" - A matrix type to be used for using the Neumann-Neumann type preconditioners, see PCBDDC. 960 This stores the matrices in globally unassembled form. Each processor 961 assembles only its local Neumann problem and the parallel matrix vector 962 product is handled "implicitly". 963 964 Operations Provided: 965 + MatMult() 966 . MatMultAdd() 967 . MatMultTranspose() 968 . MatMultTransposeAdd() 969 . MatZeroEntries() 970 . MatSetOption() 971 . MatZeroRows() 972 . MatZeroRowsLocal() 973 . MatSetValues() 974 . MatSetValuesLocal() 975 . MatScale() 976 . MatGetDiagonal() 977 - MatSetLocalToGlobalMapping() 978 979 Options Database Keys: 980 . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions() 981 982 Notes: Options prefix for the inner matrix are given by -is_mat_xxx 983 984 You must call MatSetLocalToGlobalMapping() before using this matrix type. 985 986 You can do matrix preallocation on the local matrix after you obtain it with 987 MatISGetLocalMat() 988 989 Level: advanced 990 991 .seealso: PC, MatISGetLocalMat(), MatSetLocalToGlobalMapping(), PCBDDC 992 993 M*/ 994 995 #undef __FUNCT__ 996 #define __FUNCT__ "MatCreate_IS" 997 PETSC_EXTERN PetscErrorCode MatCreate_IS(Mat A) 998 { 999 PetscErrorCode ierr; 1000 Mat_IS *b; 1001 1002 PetscFunctionBegin; 1003 ierr = PetscNewLog(A,&b);CHKERRQ(ierr); 1004 A->data = (void*)b; 1005 ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1006 1007 A->ops->mult = MatMult_IS; 1008 A->ops->multadd = MatMultAdd_IS; 1009 A->ops->multtranspose = MatMultTranspose_IS; 1010 A->ops->multtransposeadd = MatMultTransposeAdd_IS; 1011 A->ops->destroy = MatDestroy_IS; 1012 A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS; 1013 A->ops->setvalues = MatSetValues_IS; 1014 A->ops->setvalueslocal = MatSetValuesLocal_IS; 1015 A->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_IS; 1016 A->ops->zerorows = MatZeroRows_IS; 1017 A->ops->zerorowslocal = MatZeroRowsLocal_IS; 1018 A->ops->assemblybegin = MatAssemblyBegin_IS; 1019 A->ops->assemblyend = MatAssemblyEnd_IS; 1020 A->ops->view = MatView_IS; 1021 A->ops->zeroentries = MatZeroEntries_IS; 1022 A->ops->scale = MatScale_IS; 1023 A->ops->getdiagonal = MatGetDiagonal_IS; 1024 A->ops->setoption = MatSetOption_IS; 1025 A->ops->ishermitian = MatIsHermitian_IS; 1026 A->ops->issymmetric = MatIsSymmetric_IS; 1027 A->ops->duplicate = MatDuplicate_IS; 1028 1029 ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 1030 ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 1031 1032 /* zeros pointer */ 1033 b->A = 0; 1034 b->ctx = 0; 1035 b->x = 0; 1036 b->y = 0; 1037 b->mapping = 0; 1038 1039 /* special MATIS functions */ 1040 ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",MatISGetLocalMat_IS);CHKERRQ(ierr); 1041 ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",MatISSetLocalMat_IS);CHKERRQ(ierr); 1042 ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",MatISGetMPIXAIJ_IS);CHKERRQ(ierr); 1043 ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",MatISSetPreallocation_IS);CHKERRQ(ierr); 1044 ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr); 1045 PetscFunctionReturn(0); 1046 } 1047