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