1be1d678aSKris Buschelman 2b4319ba4SBarry Smith /* 3b4319ba4SBarry Smith Creates a matrix class for using the Neumann-Neumann type preconditioners. 4b4319ba4SBarry Smith This stores the matrices in globally unassembled form. Each processor 5b4319ba4SBarry Smith assembles only its local Neumann problem and the parallel matrix vector 6b4319ba4SBarry Smith product is handled "implicitly". 7b4319ba4SBarry Smith 8b4319ba4SBarry Smith We provide: 9b4319ba4SBarry Smith MatMult() 10b4319ba4SBarry Smith 11b4319ba4SBarry Smith Currently this allows for only one subdomain per processor. 12b4319ba4SBarry Smith 13b4319ba4SBarry Smith */ 14b4319ba4SBarry Smith 15c6db04a5SJed Brown #include <../src/mat/impls/is/matis.h> /*I "petscmat.h" I*/ 1628f4e0baSStefano Zampini 17*a72627d2SStefano Zampini static PetscErrorCode MatISComputeSF_Private(Mat); 18*a72627d2SStefano Zampini 1928f4e0baSStefano Zampini #undef __FUNCT__ 2028f4e0baSStefano Zampini #define __FUNCT__ "MatISComputeSF_Private" 21*a72627d2SStefano Zampini static PetscErrorCode MatISComputeSF_Private(Mat B) 2228f4e0baSStefano Zampini { 2328f4e0baSStefano Zampini Mat_IS *matis = (Mat_IS*)(B->data); 2428f4e0baSStefano Zampini const PetscInt *gidxs; 2528f4e0baSStefano Zampini PetscErrorCode ierr; 2628f4e0baSStefano Zampini 2728f4e0baSStefano Zampini PetscFunctionBegin; 2828f4e0baSStefano Zampini ierr = MatGetSize(matis->A,&matis->sf_nleaves,NULL);CHKERRQ(ierr); 2928f4e0baSStefano Zampini ierr = MatGetLocalSize(B,&matis->sf_nroots,NULL);CHKERRQ(ierr); 3028f4e0baSStefano Zampini ierr = PetscSFCreate(PetscObjectComm((PetscObject)B),&matis->sf);CHKERRQ(ierr); 3128f4e0baSStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(matis->mapping,&gidxs);CHKERRQ(ierr); 3228f4e0baSStefano Zampini /* PETSC_OWN_POINTER refers to ilocal which is NULL */ 3328f4e0baSStefano Zampini ierr = PetscSFSetGraphLayout(matis->sf,B->rmap,matis->sf_nleaves,NULL,PETSC_OWN_POINTER,gidxs);CHKERRQ(ierr); 3428f4e0baSStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(matis->mapping,&gidxs);CHKERRQ(ierr); 3528f4e0baSStefano Zampini ierr = PetscMalloc2(matis->sf_nroots,&matis->sf_rootdata,matis->sf_nleaves,&matis->sf_leafdata);CHKERRQ(ierr); 3628f4e0baSStefano Zampini PetscFunctionReturn(0); 3728f4e0baSStefano Zampini } 382e1947a5SStefano Zampini 392e1947a5SStefano Zampini #undef __FUNCT__ 402e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation" 41a88811baSStefano Zampini /* 42a88811baSStefano Zampini MatISSetPreallocation - Preallocates memory for a MATIS parallel matrix. 43a88811baSStefano Zampini 44a88811baSStefano Zampini Collective on MPI_Comm 45a88811baSStefano Zampini 46a88811baSStefano Zampini Input Parameters: 47a88811baSStefano Zampini + B - the matrix 48a88811baSStefano Zampini . d_nz - number of nonzeros per row in DIAGONAL portion of local submatrix 49a88811baSStefano Zampini (same value is used for all local rows) 50a88811baSStefano Zampini . d_nnz - array containing the number of nonzeros in the various rows of the 51a88811baSStefano Zampini DIAGONAL portion of the local submatrix (possibly different for each row) 52a88811baSStefano Zampini or NULL, if d_nz is used to specify the nonzero structure. 53a88811baSStefano Zampini The size of this array is equal to the number of local rows, i.e 'm'. 54a88811baSStefano Zampini For matrices that will be factored, you must leave room for (and set) 55a88811baSStefano Zampini the diagonal entry even if it is zero. 56a88811baSStefano Zampini . o_nz - number of nonzeros per row in the OFF-DIAGONAL portion of local 57a88811baSStefano Zampini submatrix (same value is used for all local rows). 58a88811baSStefano Zampini - o_nnz - array containing the number of nonzeros in the various rows of the 59a88811baSStefano Zampini OFF-DIAGONAL portion of the local submatrix (possibly different for 60a88811baSStefano Zampini each row) or NULL, if o_nz is used to specify the nonzero 61a88811baSStefano Zampini structure. The size of this array is equal to the number 62a88811baSStefano Zampini of local rows, i.e 'm'. 63a88811baSStefano Zampini 64a88811baSStefano Zampini If the *_nnz parameter is given then the *_nz parameter is ignored 65a88811baSStefano Zampini 66a88811baSStefano Zampini Level: intermediate 67a88811baSStefano Zampini 68a88811baSStefano Zampini Notes: This function has the same interface as the MPIAIJ preallocation routine in order to simplify the transition 69a88811baSStefano Zampini from the asssembled format to the unassembled one. It overestimates the preallocation of MATIS local 70a88811baSStefano Zampini matrices; for exact preallocation, the user should set the preallocation directly on local matrix objects. 71a88811baSStefano Zampini 72a88811baSStefano Zampini .keywords: matrix 73a88811baSStefano Zampini 74a88811baSStefano Zampini .seealso: MatCreate(), MatCreateIS(), MatMPIAIJSetPreallocation(), MatISGetLocalMat() 75a88811baSStefano Zampini @*/ 762e1947a5SStefano Zampini PetscErrorCode MatISSetPreallocation(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[]) 772e1947a5SStefano Zampini { 782e1947a5SStefano Zampini PetscErrorCode ierr; 792e1947a5SStefano Zampini 802e1947a5SStefano Zampini PetscFunctionBegin; 812e1947a5SStefano Zampini PetscValidHeaderSpecific(B,MAT_CLASSID,1); 822e1947a5SStefano Zampini PetscValidType(B,1); 832e1947a5SStefano Zampini ierr = PetscTryMethod(B,"MatISSetPreallocation_C",(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(B,d_nz,d_nnz,o_nz,o_nnz));CHKERRQ(ierr); 842e1947a5SStefano Zampini PetscFunctionReturn(0); 852e1947a5SStefano Zampini } 862e1947a5SStefano Zampini 872e1947a5SStefano Zampini #undef __FUNCT__ 882e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation_IS" 892e1947a5SStefano Zampini PetscErrorCode MatISSetPreallocation_IS(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[]) 902e1947a5SStefano Zampini { 912e1947a5SStefano Zampini Mat_IS *matis = (Mat_IS*)(B->data); 9228f4e0baSStefano Zampini PetscInt bs,i,nlocalcols; 932e1947a5SStefano Zampini PetscErrorCode ierr; 942e1947a5SStefano Zampini 952e1947a5SStefano Zampini PetscFunctionBegin; 962e1947a5SStefano Zampini if (!matis->A) { 972e1947a5SStefano Zampini SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"You should first call MatSetLocalToGlobalMapping"); 982e1947a5SStefano Zampini } 9928f4e0baSStefano Zampini if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */ 10028f4e0baSStefano Zampini ierr = MatISComputeSF_Private(B);CHKERRQ(ierr); 10128f4e0baSStefano Zampini } 1022e1947a5SStefano Zampini if (!d_nnz) { 10328f4e0baSStefano Zampini for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] = d_nz; 1042e1947a5SStefano Zampini } else { 10528f4e0baSStefano Zampini for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] = d_nnz[i]; 1062e1947a5SStefano Zampini } 1072e1947a5SStefano Zampini if (!o_nnz) { 10828f4e0baSStefano Zampini for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] += o_nz; 1092e1947a5SStefano Zampini } else { 11028f4e0baSStefano Zampini for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] += o_nnz[i]; 1112e1947a5SStefano Zampini } 11228f4e0baSStefano Zampini ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 11328f4e0baSStefano Zampini ierr = MatGetSize(matis->A,NULL,&nlocalcols);CHKERRQ(ierr); 11428f4e0baSStefano Zampini ierr = MatGetBlockSize(matis->A,&bs);CHKERRQ(ierr); 11528f4e0baSStefano Zampini ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 11628f4e0baSStefano Zampini for (i=0;i<matis->sf_nleaves;i++) { 11728f4e0baSStefano Zampini matis->sf_leafdata[i] = PetscMin(matis->sf_leafdata[i],nlocalcols); 1182e1947a5SStefano Zampini } 11928f4e0baSStefano Zampini ierr = MatSeqAIJSetPreallocation(matis->A,0,matis->sf_leafdata);CHKERRQ(ierr); 12028f4e0baSStefano Zampini for (i=0;i<matis->sf_nleaves/bs;i++) { 12128f4e0baSStefano Zampini matis->sf_leafdata[i] = matis->sf_leafdata[i*bs]/bs; 1222e1947a5SStefano Zampini } 12328f4e0baSStefano Zampini ierr = MatSeqBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr); 12428f4e0baSStefano Zampini for (i=0;i<matis->sf_nleaves/bs;i++) { 12528f4e0baSStefano Zampini matis->sf_leafdata[i] = matis->sf_leafdata[i]-i; 1262e1947a5SStefano Zampini } 12728f4e0baSStefano Zampini ierr = MatSeqSBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr); 1282e1947a5SStefano Zampini PetscFunctionReturn(0); 1292e1947a5SStefano Zampini } 130b4319ba4SBarry Smith 131b4319ba4SBarry Smith #undef __FUNCT__ 132b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ_IS" 133b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ_IS(Mat mat, MatReuse reuse, Mat *M) 134b7ce53b6SStefano Zampini { 135b7ce53b6SStefano Zampini Mat_IS *matis = (Mat_IS*)(mat->data); 136b7ce53b6SStefano Zampini /* info on mat */ 137b7ce53b6SStefano Zampini /* ISLocalToGlobalMapping rmapping,cmapping; */ 138b7ce53b6SStefano Zampini PetscInt bs,rows,cols; 139b7ce53b6SStefano Zampini PetscInt lrows,lcols; 140b7ce53b6SStefano Zampini PetscInt local_rows,local_cols; 141b7ce53b6SStefano Zampini PetscBool isdense,issbaij,issbaij_red; 1427c03b4e8SStefano Zampini PetscMPIInt nsubdomains; 143b7ce53b6SStefano Zampini /* values insertion */ 144b7ce53b6SStefano Zampini PetscScalar *array; 145b7ce53b6SStefano Zampini PetscInt *local_indices,*global_indices; 146b7ce53b6SStefano Zampini /* work */ 147b7ce53b6SStefano Zampini PetscInt i,j,index_row; 148b7ce53b6SStefano Zampini PetscErrorCode ierr; 149b7ce53b6SStefano Zampini 150b7ce53b6SStefano Zampini PetscFunctionBegin; 151b7ce53b6SStefano Zampini /* MISSING CHECKS 152b7ce53b6SStefano Zampini - rectangular case not covered (it is not allowed by MATIS) 153b7ce53b6SStefano Zampini */ 154b7ce53b6SStefano Zampini /* get info from mat */ 1557c03b4e8SStefano Zampini ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&nsubdomains);CHKERRQ(ierr); 1567c03b4e8SStefano Zampini if (nsubdomains == 1) { 1577c03b4e8SStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { 1587c03b4e8SStefano Zampini ierr = MatDuplicate(matis->A,MAT_COPY_VALUES,&(*M));CHKERRQ(ierr); 1597c03b4e8SStefano Zampini } else { 1607c03b4e8SStefano Zampini ierr = MatCopy(matis->A,*M,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 1617c03b4e8SStefano Zampini } 1627c03b4e8SStefano Zampini PetscFunctionReturn(0); 1637c03b4e8SStefano Zampini } 164b7ce53b6SStefano Zampini /* ierr = MatGetLocalToGlobalMapping(mat,&rmapping,&cmapping);CHKERRQ(ierr); */ 165b7ce53b6SStefano Zampini ierr = MatGetSize(mat,&rows,&cols);CHKERRQ(ierr); 166b7ce53b6SStefano Zampini ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 167b7ce53b6SStefano Zampini ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr); 168b7ce53b6SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr); 169b7ce53b6SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr); 170b7ce53b6SStefano Zampini 171b7ce53b6SStefano Zampini /* work */ 172b7ce53b6SStefano Zampini ierr = PetscMalloc1(local_rows,&local_indices);CHKERRQ(ierr); 173b7ce53b6SStefano Zampini for (i=0;i<local_rows;i++) local_indices[i]=i; 174b7ce53b6SStefano Zampini /* map indices of local mat to global values */ 175854ce69bSBarry Smith ierr = PetscMalloc1(PetscMax(local_cols,local_rows),&global_indices);CHKERRQ(ierr); 176b7ce53b6SStefano Zampini /* ierr = ISLocalToGlobalMappingApply(rmapping,local_rows,local_indices,global_indices);CHKERRQ(ierr); */ 177b7ce53b6SStefano Zampini ierr = ISLocalToGlobalMappingApply(matis->mapping,local_rows,local_indices,global_indices);CHKERRQ(ierr); 178b7ce53b6SStefano Zampini 179b7ce53b6SStefano Zampini if (issbaij) { 180b7ce53b6SStefano Zampini ierr = MatGetRowUpperTriangular(matis->A);CHKERRQ(ierr); 181b7ce53b6SStefano Zampini } 182b7ce53b6SStefano Zampini 183b7ce53b6SStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { 184b7ce53b6SStefano Zampini Mat new_mat; 185b7ce53b6SStefano Zampini MatType new_mat_type; 186*a72627d2SStefano Zampini PetscInt *my_dnz,*my_onz; 187b7ce53b6SStefano Zampini PetscInt *dnz,*onz,*mat_ranges,*row_ownership; 188b7ce53b6SStefano Zampini PetscInt index_col,owner; 189*a72627d2SStefano Zampini PetscBool maxreduce = PETSC_FALSE; 190b7ce53b6SStefano Zampini 191b7ce53b6SStefano Zampini /* determining new matrix type */ 192b7ce53b6SStefano Zampini ierr = MPI_Allreduce(&issbaij,&issbaij_red,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 193b7ce53b6SStefano Zampini if (issbaij_red) { 194b7ce53b6SStefano Zampini new_mat_type = MATSBAIJ; 195b7ce53b6SStefano Zampini } else { 196b7ce53b6SStefano Zampini if (bs>1) { 197b7ce53b6SStefano Zampini new_mat_type = MATBAIJ; 198b7ce53b6SStefano Zampini } else { 199b7ce53b6SStefano Zampini new_mat_type = MATAIJ; 200b7ce53b6SStefano Zampini } 201b7ce53b6SStefano Zampini } 202b7ce53b6SStefano Zampini 203b7ce53b6SStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)mat),&new_mat);CHKERRQ(ierr); 204b7ce53b6SStefano Zampini ierr = MatSetSizes(new_mat,PETSC_DECIDE,PETSC_DECIDE,rows,cols);CHKERRQ(ierr); 205b7ce53b6SStefano Zampini ierr = MatSetBlockSize(new_mat,bs);CHKERRQ(ierr); 206b7ce53b6SStefano Zampini ierr = MatSetType(new_mat,new_mat_type);CHKERRQ(ierr); 207b7ce53b6SStefano Zampini ierr = MatSetUp(new_mat);CHKERRQ(ierr); 208b7ce53b6SStefano Zampini ierr = MatGetLocalSize(new_mat,&lrows,&lcols);CHKERRQ(ierr); 209*a72627d2SStefano Zampini if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */ 210*a72627d2SStefano Zampini ierr = MatISComputeSF_Private(mat);CHKERRQ(ierr); 211*a72627d2SStefano Zampini } 212b7ce53b6SStefano Zampini 213b7ce53b6SStefano Zampini /* 214b7ce53b6SStefano Zampini preallocation 215b7ce53b6SStefano Zampini */ 216b7ce53b6SStefano Zampini 217b7ce53b6SStefano Zampini ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)new_mat),lrows,lcols,dnz,onz);CHKERRQ(ierr); 218b7ce53b6SStefano Zampini /* 219*a72627d2SStefano Zampini An SF reduce is needed to sum up properly on shared interface dofs. 220b7ce53b6SStefano Zampini Note that preallocation is not exact, since it overestimates nonzeros 221b7ce53b6SStefano Zampini */ 222b7ce53b6SStefano Zampini /* All processes need to compute entire row ownership */ 223b7ce53b6SStefano Zampini ierr = PetscMalloc1(rows,&row_ownership);CHKERRQ(ierr); 224b7ce53b6SStefano Zampini ierr = MatGetOwnershipRanges(new_mat,(const PetscInt**)&mat_ranges);CHKERRQ(ierr); 225b7ce53b6SStefano Zampini for (i=0;i<nsubdomains;i++) { 226b7ce53b6SStefano Zampini for (j=mat_ranges[i];j<mat_ranges[i+1];j++) { 227b7ce53b6SStefano Zampini row_ownership[j]=i; 228b7ce53b6SStefano Zampini } 229b7ce53b6SStefano Zampini } 230b7ce53b6SStefano Zampini 231b7ce53b6SStefano Zampini /* 232b7ce53b6SStefano Zampini my_dnz and my_onz contains exact contribution to preallocation from each local mat 233b7ce53b6SStefano Zampini then, they will be summed up properly. This way, preallocation is always sufficient 234b7ce53b6SStefano Zampini */ 235*a72627d2SStefano Zampini ierr = PetscCalloc2(local_rows,&my_dnz,local_rows,&my_onz);CHKERRQ(ierr); 236b7ce53b6SStefano Zampini /* preallocation as a MATAIJ */ 237b7ce53b6SStefano Zampini if (isdense) { /* special case for dense local matrices */ 238b7ce53b6SStefano Zampini for (i=0;i<local_rows;i++) { 239b7ce53b6SStefano Zampini index_row = global_indices[i]; 240b7ce53b6SStefano Zampini for (j=i;j<local_rows;j++) { 241b7ce53b6SStefano Zampini owner = row_ownership[index_row]; 242b7ce53b6SStefano Zampini index_col = global_indices[j]; 243b7ce53b6SStefano Zampini if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */ 244*a72627d2SStefano Zampini my_dnz[i] += 1; 245b7ce53b6SStefano Zampini } else { /* offdiag block */ 246*a72627d2SStefano Zampini my_onz[i] += 1; 247b7ce53b6SStefano Zampini } 248b7ce53b6SStefano Zampini /* same as before, interchanging rows and cols */ 249b7ce53b6SStefano Zampini if (i != j) { 250b7ce53b6SStefano Zampini owner = row_ownership[index_col]; 251b7ce53b6SStefano Zampini if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) { 252*a72627d2SStefano Zampini my_dnz[j] += 1; 253b7ce53b6SStefano Zampini } else { 254*a72627d2SStefano Zampini my_onz[j] += 1; 255b7ce53b6SStefano Zampini } 256b7ce53b6SStefano Zampini } 257b7ce53b6SStefano Zampini } 258b7ce53b6SStefano Zampini } 259b7ce53b6SStefano Zampini } else { 260b7ce53b6SStefano Zampini for (i=0;i<local_rows;i++) { 261b7ce53b6SStefano Zampini PetscInt ncols; 262b7ce53b6SStefano Zampini const PetscInt *cols; 263b7ce53b6SStefano Zampini index_row = global_indices[i]; 264b7ce53b6SStefano Zampini ierr = MatGetRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr); 265b7ce53b6SStefano Zampini for (j=0;j<ncols;j++) { 266b7ce53b6SStefano Zampini owner = row_ownership[index_row]; 267b7ce53b6SStefano Zampini index_col = global_indices[cols[j]]; 268b7ce53b6SStefano Zampini if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */ 269*a72627d2SStefano Zampini my_dnz[i] += 1; 270b7ce53b6SStefano Zampini } else { /* offdiag block */ 271*a72627d2SStefano Zampini my_onz[i] += 1; 272b7ce53b6SStefano Zampini } 273b7ce53b6SStefano Zampini /* same as before, interchanging rows and cols */ 274b7ce53b6SStefano Zampini if (issbaij) { 275b7ce53b6SStefano Zampini owner = row_ownership[index_col]; 276b7ce53b6SStefano Zampini if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) { 277*a72627d2SStefano Zampini my_dnz[j] += 1; 278b7ce53b6SStefano Zampini } else { 279*a72627d2SStefano Zampini my_onz[j] += 1; 280b7ce53b6SStefano Zampini } 281b7ce53b6SStefano Zampini } 282b7ce53b6SStefano Zampini } 283b7ce53b6SStefano Zampini ierr = MatRestoreRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr); 284b7ce53b6SStefano Zampini } 285b7ce53b6SStefano Zampini } 286b7ce53b6SStefano Zampini ierr = PetscFree(row_ownership);CHKERRQ(ierr); 287*a72627d2SStefano Zampini if (maxreduce) { 288*a72627d2SStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr); 289*a72627d2SStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr); 290*a72627d2SStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr); 291*a72627d2SStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr); 292*a72627d2SStefano Zampini } else { 293*a72627d2SStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr); 294*a72627d2SStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr); 295*a72627d2SStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr); 296*a72627d2SStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr); 297*a72627d2SStefano Zampini } 298*a72627d2SStefano Zampini ierr = PetscFree2(my_dnz,my_onz);CHKERRQ(ierr); 299b7ce53b6SStefano Zampini 300b7ce53b6SStefano Zampini /* Resize preallocation if overestimated */ 301b7ce53b6SStefano Zampini for (i=0;i<lrows;i++) { 302b7ce53b6SStefano Zampini dnz[i] = PetscMin(dnz[i],lcols); 303b7ce53b6SStefano Zampini onz[i] = PetscMin(onz[i],cols-lcols); 304b7ce53b6SStefano Zampini } 305b7ce53b6SStefano Zampini /* set preallocation */ 306b7ce53b6SStefano Zampini ierr = MatMPIAIJSetPreallocation(new_mat,0,dnz,0,onz);CHKERRQ(ierr); 307b7ce53b6SStefano Zampini for (i=0;i<lrows/bs;i++) { 308b7ce53b6SStefano Zampini dnz[i] = dnz[i*bs]/bs; 309b7ce53b6SStefano Zampini onz[i] = onz[i*bs]/bs; 310b7ce53b6SStefano Zampini } 311b7ce53b6SStefano Zampini ierr = MatMPIBAIJSetPreallocation(new_mat,bs,0,dnz,0,onz);CHKERRQ(ierr); 312b7ce53b6SStefano Zampini for (i=0;i<lrows/bs;i++) { 313b7ce53b6SStefano Zampini dnz[i] = dnz[i]-i; 314b7ce53b6SStefano Zampini } 315b7ce53b6SStefano Zampini ierr = MatMPISBAIJSetPreallocation(new_mat,bs,0,dnz,0,onz);CHKERRQ(ierr); 316b7ce53b6SStefano Zampini ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 317b7ce53b6SStefano Zampini *M = new_mat; 318b7ce53b6SStefano Zampini } else { 319b7ce53b6SStefano Zampini PetscInt mbs,mrows,mcols; 320b7ce53b6SStefano Zampini /* some checks */ 321b7ce53b6SStefano Zampini ierr = MatGetBlockSize(*M,&mbs);CHKERRQ(ierr); 322b7ce53b6SStefano Zampini ierr = MatGetSize(*M,&mrows,&mcols);CHKERRQ(ierr); 323b7ce53b6SStefano Zampini if (mrows != rows) { 324b7ce53b6SStefano Zampini SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of rows (%d != %d)",rows,mrows); 325b7ce53b6SStefano Zampini } 326b7ce53b6SStefano Zampini if (mrows != rows) { 327b7ce53b6SStefano Zampini SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of cols (%d != %d)",cols,mcols); 328b7ce53b6SStefano Zampini } 329b7ce53b6SStefano Zampini if (mbs != bs) { 330b7ce53b6SStefano Zampini SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong block size (%d != %d)",bs,mbs); 331b7ce53b6SStefano Zampini } 332b7ce53b6SStefano Zampini ierr = MatZeroEntries(*M);CHKERRQ(ierr); 333b7ce53b6SStefano Zampini } 334b7ce53b6SStefano Zampini /* set local to global mappings */ 335b7ce53b6SStefano Zampini /* ierr = MatSetLocalToGlobalMapping(*M,rmapping,cmapping);CHKERRQ(ierr); */ 336b7ce53b6SStefano Zampini /* Set values */ 337b7ce53b6SStefano Zampini if (isdense) { /* special case for dense local matrices */ 338b7ce53b6SStefano Zampini ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr); 339b7ce53b6SStefano Zampini ierr = MatDenseGetArray(matis->A,&array);CHKERRQ(ierr); 340b7ce53b6SStefano Zampini ierr = MatSetValues(*M,local_rows,global_indices,local_cols,global_indices,array,ADD_VALUES);CHKERRQ(ierr); 341b7ce53b6SStefano Zampini ierr = MatDenseRestoreArray(matis->A,&array);CHKERRQ(ierr); 342b7ce53b6SStefano Zampini ierr = PetscFree(local_indices);CHKERRQ(ierr); 343b7ce53b6SStefano Zampini ierr = PetscFree(global_indices);CHKERRQ(ierr); 344b7ce53b6SStefano Zampini } else { /* very basic values insertion for all other matrix types */ 345b7ce53b6SStefano Zampini ierr = PetscFree(local_indices);CHKERRQ(ierr); 346b7ce53b6SStefano Zampini for (i=0;i<local_rows;i++) { 347b7ce53b6SStefano Zampini ierr = MatGetRow(matis->A,i,&j,(const PetscInt**)&local_indices,(const PetscScalar**)&array);CHKERRQ(ierr); 348b7ce53b6SStefano Zampini /* ierr = MatSetValuesLocal(*M,1,&i,j,local_indices,array,ADD_VALUES);CHKERRQ(ierr); */ 349b7ce53b6SStefano Zampini ierr = ISLocalToGlobalMappingApply(matis->mapping,j,local_indices,global_indices);CHKERRQ(ierr); 350b7ce53b6SStefano Zampini ierr = ISLocalToGlobalMappingApply(matis->mapping,1,&i,&index_row);CHKERRQ(ierr); 351b7ce53b6SStefano Zampini ierr = MatSetValues(*M,1,&index_row,j,global_indices,array,ADD_VALUES);CHKERRQ(ierr); 352b7ce53b6SStefano Zampini ierr = MatRestoreRow(matis->A,i,&j,(const PetscInt**)&local_indices,(const PetscScalar**)&array);CHKERRQ(ierr); 353b7ce53b6SStefano Zampini } 354b7ce53b6SStefano Zampini ierr = PetscFree(global_indices);CHKERRQ(ierr); 355b7ce53b6SStefano Zampini } 356b7ce53b6SStefano Zampini ierr = MatAssemblyBegin(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 357b7ce53b6SStefano Zampini ierr = MatAssemblyEnd(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 358b7ce53b6SStefano Zampini if (isdense) { 359b7ce53b6SStefano Zampini ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr); 360b7ce53b6SStefano Zampini } 361b7ce53b6SStefano Zampini if (issbaij) { 362b7ce53b6SStefano Zampini ierr = MatRestoreRowUpperTriangular(matis->A);CHKERRQ(ierr); 363b7ce53b6SStefano Zampini } 364b7ce53b6SStefano Zampini PetscFunctionReturn(0); 365b7ce53b6SStefano Zampini } 366b7ce53b6SStefano Zampini 367b7ce53b6SStefano Zampini #undef __FUNCT__ 368b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ" 369b7ce53b6SStefano Zampini /*@ 370b7ce53b6SStefano Zampini MatISGetMPIXAIJ - Converts MATIS matrix into a parallel AIJ format 371b7ce53b6SStefano Zampini 372b7ce53b6SStefano Zampini Input Parameter: 373b7ce53b6SStefano Zampini . mat - the matrix (should be of type MATIS) 374b7ce53b6SStefano Zampini . reuse - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX 375b7ce53b6SStefano Zampini 376b7ce53b6SStefano Zampini Output Parameter: 377b7ce53b6SStefano Zampini . newmat - the matrix in AIJ format 378b7ce53b6SStefano Zampini 379b7ce53b6SStefano Zampini Level: developer 380b7ce53b6SStefano Zampini 381b7ce53b6SStefano Zampini Notes: 382b7ce53b6SStefano Zampini 383b7ce53b6SStefano Zampini .seealso: MATIS 384b7ce53b6SStefano Zampini @*/ 385b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ(Mat mat, MatReuse reuse, Mat *newmat) 386b7ce53b6SStefano Zampini { 387b7ce53b6SStefano Zampini PetscErrorCode ierr; 388b7ce53b6SStefano Zampini 389b7ce53b6SStefano Zampini PetscFunctionBegin; 390b7ce53b6SStefano Zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 391b7ce53b6SStefano Zampini PetscValidLogicalCollectiveEnum(mat,reuse,2); 392b7ce53b6SStefano Zampini PetscValidPointer(newmat,3); 393b7ce53b6SStefano Zampini if (reuse != MAT_INITIAL_MATRIX) { 394b7ce53b6SStefano Zampini PetscValidHeaderSpecific(*newmat,MAT_CLASSID,3); 395b7ce53b6SStefano Zampini PetscCheckSameComm(mat,1,*newmat,3); 396b7ce53b6SStefano Zampini } 397b7ce53b6SStefano Zampini ierr = PetscUseMethod(mat,"MatISGetMPIXAIJ_C",(Mat,MatReuse,Mat*),(mat,reuse,newmat));CHKERRQ(ierr); 398b7ce53b6SStefano Zampini PetscFunctionReturn(0); 399b7ce53b6SStefano Zampini } 400b7ce53b6SStefano Zampini 401b7ce53b6SStefano Zampini #undef __FUNCT__ 402ad6194a2SStefano Zampini #define __FUNCT__ "MatDuplicate_IS" 403ad6194a2SStefano Zampini PetscErrorCode MatDuplicate_IS(Mat mat,MatDuplicateOption op,Mat *newmat) 404ad6194a2SStefano Zampini { 405ad6194a2SStefano Zampini PetscErrorCode ierr; 406ad6194a2SStefano Zampini Mat_IS *matis = (Mat_IS*)(mat->data); 407ad6194a2SStefano Zampini PetscInt bs,m,n,M,N; 408ad6194a2SStefano Zampini Mat B,localmat; 409ad6194a2SStefano Zampini 410ad6194a2SStefano Zampini PetscFunctionBegin; 411ad6194a2SStefano Zampini ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 412ad6194a2SStefano Zampini ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 413ad6194a2SStefano Zampini ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 414ad6194a2SStefano Zampini ierr = MatCreateIS(PetscObjectComm((PetscObject)mat),bs,m,n,M,N,matis->mapping,&B);CHKERRQ(ierr); 415ad6194a2SStefano Zampini ierr = MatDuplicate(matis->A,op,&localmat);CHKERRQ(ierr); 416ad6194a2SStefano Zampini ierr = MatISSetLocalMat(B,localmat);CHKERRQ(ierr); 417b3317aa8SStefano Zampini ierr = MatDestroy(&localmat);CHKERRQ(ierr); 418ad6194a2SStefano Zampini ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 419ad6194a2SStefano Zampini ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 420ad6194a2SStefano Zampini *newmat = B; 421ad6194a2SStefano Zampini PetscFunctionReturn(0); 422ad6194a2SStefano Zampini } 423ad6194a2SStefano Zampini 424ad6194a2SStefano Zampini #undef __FUNCT__ 42569796d55SStefano Zampini #define __FUNCT__ "MatIsHermitian_IS" 42669796d55SStefano Zampini PetscErrorCode MatIsHermitian_IS(Mat A,PetscReal tol,PetscBool *flg) 42769796d55SStefano Zampini { 42869796d55SStefano Zampini PetscErrorCode ierr; 42969796d55SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 43069796d55SStefano Zampini PetscBool local_sym; 43169796d55SStefano Zampini 43269796d55SStefano Zampini PetscFunctionBegin; 43369796d55SStefano Zampini ierr = MatIsHermitian(matis->A,tol,&local_sym);CHKERRQ(ierr); 43469796d55SStefano Zampini ierr = MPI_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 43569796d55SStefano Zampini PetscFunctionReturn(0); 43669796d55SStefano Zampini } 43769796d55SStefano Zampini 43869796d55SStefano Zampini #undef __FUNCT__ 43969796d55SStefano Zampini #define __FUNCT__ "MatIsSymmetric_IS" 44069796d55SStefano Zampini PetscErrorCode MatIsSymmetric_IS(Mat A,PetscReal tol,PetscBool *flg) 44169796d55SStefano Zampini { 44269796d55SStefano Zampini PetscErrorCode ierr; 44369796d55SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 44469796d55SStefano Zampini PetscBool local_sym; 44569796d55SStefano Zampini 44669796d55SStefano Zampini PetscFunctionBegin; 44769796d55SStefano Zampini ierr = MatIsSymmetric(matis->A,tol,&local_sym);CHKERRQ(ierr); 44869796d55SStefano Zampini ierr = MPI_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 44969796d55SStefano Zampini PetscFunctionReturn(0); 45069796d55SStefano Zampini } 45169796d55SStefano Zampini 45269796d55SStefano Zampini #undef __FUNCT__ 453b4319ba4SBarry Smith #define __FUNCT__ "MatDestroy_IS" 454dfbe8321SBarry Smith PetscErrorCode MatDestroy_IS(Mat A) 455b4319ba4SBarry Smith { 456dfbe8321SBarry Smith PetscErrorCode ierr; 457b4319ba4SBarry Smith Mat_IS *b = (Mat_IS*)A->data; 458b4319ba4SBarry Smith 459b4319ba4SBarry Smith PetscFunctionBegin; 4606bf464f9SBarry Smith ierr = MatDestroy(&b->A);CHKERRQ(ierr); 4616bf464f9SBarry Smith ierr = VecScatterDestroy(&b->ctx);CHKERRQ(ierr); 4626bf464f9SBarry Smith ierr = VecDestroy(&b->x);CHKERRQ(ierr); 4636bf464f9SBarry Smith ierr = VecDestroy(&b->y);CHKERRQ(ierr); 4646bf464f9SBarry Smith ierr = ISLocalToGlobalMappingDestroy(&b->mapping);CHKERRQ(ierr); 46528f4e0baSStefano Zampini ierr = PetscSFDestroy(&b->sf);CHKERRQ(ierr); 46628f4e0baSStefano Zampini ierr = PetscFree2(b->sf_rootdata,b->sf_leafdata);CHKERRQ(ierr); 467bf0cc555SLisandro Dalcin ierr = PetscFree(A->data);CHKERRQ(ierr); 468dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr); 469bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",NULL);CHKERRQ(ierr); 470b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",NULL);CHKERRQ(ierr); 471b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",NULL);CHKERRQ(ierr); 4722e1947a5SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",NULL);CHKERRQ(ierr); 473b4319ba4SBarry Smith PetscFunctionReturn(0); 474b4319ba4SBarry Smith } 475b4319ba4SBarry Smith 476b4319ba4SBarry Smith #undef __FUNCT__ 477b4319ba4SBarry Smith #define __FUNCT__ "MatMult_IS" 478dfbe8321SBarry Smith PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y) 479b4319ba4SBarry Smith { 480dfbe8321SBarry Smith PetscErrorCode ierr; 481b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 482b4319ba4SBarry Smith PetscScalar zero = 0.0; 483b4319ba4SBarry Smith 484b4319ba4SBarry Smith PetscFunctionBegin; 485b4319ba4SBarry Smith /* scatter the global vector x into the local work vector */ 486ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 487ca9f406cSSatish Balay ierr = VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 488b4319ba4SBarry Smith 489b4319ba4SBarry Smith /* multiply the local matrix */ 490b4319ba4SBarry Smith ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr); 491b4319ba4SBarry Smith 492b4319ba4SBarry Smith /* scatter product back into global memory */ 4932dcb1b2aSMatthew Knepley ierr = VecSet(y,zero);CHKERRQ(ierr); 494ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 495ca9f406cSSatish Balay ierr = VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 496b4319ba4SBarry Smith PetscFunctionReturn(0); 497b4319ba4SBarry Smith } 498b4319ba4SBarry Smith 499b4319ba4SBarry Smith #undef __FUNCT__ 5002e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultAdd_IS" 501b7ce53b6SStefano Zampini PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 5022e74eeadSLisandro Dalcin { 503650997f4SStefano Zampini Vec temp_vec; 5042e74eeadSLisandro Dalcin PetscErrorCode ierr; 5052e74eeadSLisandro Dalcin 5062e74eeadSLisandro Dalcin PetscFunctionBegin; /* v3 = v2 + A * v1.*/ 507650997f4SStefano Zampini if (v3 != v2) { 508650997f4SStefano Zampini ierr = MatMult(A,v1,v3);CHKERRQ(ierr); 509650997f4SStefano Zampini ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 510650997f4SStefano Zampini } else { 511650997f4SStefano Zampini ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 512650997f4SStefano Zampini ierr = MatMult(A,v1,temp_vec);CHKERRQ(ierr); 513650997f4SStefano Zampini ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 514650997f4SStefano Zampini ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 515650997f4SStefano Zampini ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 516650997f4SStefano Zampini } 5172e74eeadSLisandro Dalcin PetscFunctionReturn(0); 5182e74eeadSLisandro Dalcin } 5192e74eeadSLisandro Dalcin 5202e74eeadSLisandro Dalcin #undef __FUNCT__ 5212e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTranspose_IS" 5222e74eeadSLisandro Dalcin PetscErrorCode MatMultTranspose_IS(Mat A,Vec x,Vec y) 5232e74eeadSLisandro Dalcin { 5242e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 5252e74eeadSLisandro Dalcin PetscErrorCode ierr; 5262e74eeadSLisandro Dalcin 5272e74eeadSLisandro Dalcin PetscFunctionBegin; /* y = A' * x */ 5282e74eeadSLisandro Dalcin /* scatter the global vector x into the local work vector */ 529ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 530ca9f406cSSatish Balay ierr = VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5312e74eeadSLisandro Dalcin 5322e74eeadSLisandro Dalcin /* multiply the local matrix */ 5332e74eeadSLisandro Dalcin ierr = MatMultTranspose(is->A,is->x,is->y);CHKERRQ(ierr); 5342e74eeadSLisandro Dalcin 5352e74eeadSLisandro Dalcin /* scatter product back into global vector */ 5362e74eeadSLisandro Dalcin ierr = VecSet(y,0);CHKERRQ(ierr); 537ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 538ca9f406cSSatish Balay ierr = VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 5392e74eeadSLisandro Dalcin PetscFunctionReturn(0); 5402e74eeadSLisandro Dalcin } 5412e74eeadSLisandro Dalcin 5422e74eeadSLisandro Dalcin #undef __FUNCT__ 5432e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTransposeAdd_IS" 5442e74eeadSLisandro Dalcin PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 5452e74eeadSLisandro Dalcin { 546650997f4SStefano Zampini Vec temp_vec; 5472e74eeadSLisandro Dalcin PetscErrorCode ierr; 5482e74eeadSLisandro Dalcin 5492e74eeadSLisandro Dalcin PetscFunctionBegin; /* v3 = v2 + A' * v1.*/ 550650997f4SStefano Zampini if (v3 != v2) { 551650997f4SStefano Zampini ierr = MatMultTranspose(A,v1,v3);CHKERRQ(ierr); 552650997f4SStefano Zampini ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 553650997f4SStefano Zampini } else { 554650997f4SStefano Zampini ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 555650997f4SStefano Zampini ierr = MatMultTranspose(A,v1,temp_vec);CHKERRQ(ierr); 556650997f4SStefano Zampini ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 557650997f4SStefano Zampini ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 558650997f4SStefano Zampini ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 559650997f4SStefano Zampini } 5602e74eeadSLisandro Dalcin PetscFunctionReturn(0); 5612e74eeadSLisandro Dalcin } 5622e74eeadSLisandro Dalcin 5632e74eeadSLisandro Dalcin #undef __FUNCT__ 564b4319ba4SBarry Smith #define __FUNCT__ "MatView_IS" 565dfbe8321SBarry Smith PetscErrorCode MatView_IS(Mat A,PetscViewer viewer) 566b4319ba4SBarry Smith { 567b4319ba4SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 568dfbe8321SBarry Smith PetscErrorCode ierr; 569b4319ba4SBarry Smith PetscViewer sviewer; 570b4319ba4SBarry Smith 571b4319ba4SBarry Smith PetscFunctionBegin; 572dd2fa690SBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr); 573b4319ba4SBarry Smith ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 574b4319ba4SBarry Smith ierr = MatView(a->A,sviewer);CHKERRQ(ierr); 575b4319ba4SBarry Smith ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 576b4319ba4SBarry Smith PetscFunctionReturn(0); 577b4319ba4SBarry Smith } 578b4319ba4SBarry Smith 579b4319ba4SBarry Smith #undef __FUNCT__ 580b4319ba4SBarry Smith #define __FUNCT__ "MatSetLocalToGlobalMapping_IS" 581784ac674SJed Brown PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping) 582b4319ba4SBarry Smith { 583dfbe8321SBarry Smith PetscErrorCode ierr; 5844e4c7dbeSStefano Zampini PetscInt n,bs; 585b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 586b4319ba4SBarry Smith IS from,to; 587b4319ba4SBarry Smith Vec global; 588b4319ba4SBarry Smith 589b4319ba4SBarry Smith PetscFunctionBegin; 590784ac674SJed Brown PetscCheckSameComm(A,1,rmapping,2); 591ce94432eSBarry Smith if (rmapping != cmapping) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"MATIS requires the row and column mappings to be identical"); 59270cf5478SStefano Zampini if (is->mapping) { /* Currenly destroys the objects that will be created by this routine. Is there anything else that should be checked? */ 59370cf5478SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&is->mapping);CHKERRQ(ierr); 59470cf5478SStefano Zampini ierr = VecDestroy(&is->x);CHKERRQ(ierr); 59570cf5478SStefano Zampini ierr = VecDestroy(&is->y);CHKERRQ(ierr); 59670cf5478SStefano Zampini ierr = VecScatterDestroy(&is->ctx);CHKERRQ(ierr); 5971c47cb0fSStefano Zampini ierr = MatDestroy(&is->A);CHKERRQ(ierr); 59828f4e0baSStefano Zampini ierr = PetscSFDestroy(&is->sf);CHKERRQ(ierr); 59928f4e0baSStefano Zampini ierr = PetscFree2(is->sf_rootdata,is->sf_leafdata);CHKERRQ(ierr); 60070cf5478SStefano Zampini } 601784ac674SJed Brown ierr = PetscObjectReference((PetscObject)rmapping);CHKERRQ(ierr); 6026bf464f9SBarry Smith ierr = ISLocalToGlobalMappingDestroy(&is->mapping);CHKERRQ(ierr); 603784ac674SJed Brown is->mapping = rmapping; 604fa7f1dd8SStefano Zampini /* 605fa7f1dd8SStefano Zampini ierr = PetscLayoutSetISLocalToGlobalMapping(A->rmap,rmapping);CHKERRQ(ierr); 606fa7f1dd8SStefano Zampini ierr = PetscLayoutSetISLocalToGlobalMapping(A->cmap,cmapping);CHKERRQ(ierr); 607fa7f1dd8SStefano Zampini */ 608b4319ba4SBarry Smith 609b4319ba4SBarry Smith /* Create the local matrix A */ 610784ac674SJed Brown ierr = ISLocalToGlobalMappingGetSize(rmapping,&n);CHKERRQ(ierr); 6112e1947a5SStefano Zampini ierr = ISLocalToGlobalMappingGetBlockSize(rmapping,&bs);CHKERRQ(ierr); 612f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr); 6132e1947a5SStefano Zampini if (bs > 1) { 6142e1947a5SStefano Zampini ierr = MatSetType(is->A,MATSEQBAIJ);CHKERRQ(ierr); 6152e1947a5SStefano Zampini } else { 6162e1947a5SStefano Zampini ierr = MatSetType(is->A,MATSEQAIJ);CHKERRQ(ierr); 6172e1947a5SStefano Zampini } 618f69a0ea3SMatthew Knepley ierr = MatSetSizes(is->A,n,n,n,n);CHKERRQ(ierr); 6194e4c7dbeSStefano Zampini ierr = MatSetBlockSize(is->A,bs);CHKERRQ(ierr); 620ff130e51SJed Brown ierr = MatSetOptionsPrefix(is->A,((PetscObject)A)->prefix);CHKERRQ(ierr); 621ff130e51SJed Brown ierr = MatAppendOptionsPrefix(is->A,"is_");CHKERRQ(ierr); 622b4319ba4SBarry Smith ierr = MatSetFromOptions(is->A);CHKERRQ(ierr); 623b4319ba4SBarry Smith 624b4319ba4SBarry Smith /* Create the local work vectors */ 6254e4c7dbeSStefano Zampini ierr = VecCreate(PETSC_COMM_SELF,&is->x);CHKERRQ(ierr); 6264e4c7dbeSStefano Zampini ierr = VecSetBlockSize(is->x,bs);CHKERRQ(ierr); 6274e4c7dbeSStefano Zampini ierr = VecSetSizes(is->x,n,n);CHKERRQ(ierr); 628ff130e51SJed Brown ierr = VecSetOptionsPrefix(is->x,((PetscObject)A)->prefix);CHKERRQ(ierr); 629ff130e51SJed Brown ierr = VecAppendOptionsPrefix(is->x,"is_");CHKERRQ(ierr); 6304e4c7dbeSStefano Zampini ierr = VecSetFromOptions(is->x);CHKERRQ(ierr); 631b4319ba4SBarry Smith ierr = VecDuplicate(is->x,&is->y);CHKERRQ(ierr); 632b4319ba4SBarry Smith 633b4319ba4SBarry Smith /* setup the global to local scatter */ 634b4319ba4SBarry Smith ierr = ISCreateStride(PETSC_COMM_SELF,n,0,1,&to);CHKERRQ(ierr); 635784ac674SJed Brown ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr); 6362a7a6963SBarry Smith ierr = MatCreateVecs(A,&global,NULL);CHKERRQ(ierr); 637b4319ba4SBarry Smith ierr = VecScatterCreate(global,from,is->x,to,&is->ctx);CHKERRQ(ierr); 6386bf464f9SBarry Smith ierr = VecDestroy(&global);CHKERRQ(ierr); 6396bf464f9SBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 6406bf464f9SBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 641b4319ba4SBarry Smith PetscFunctionReturn(0); 642b4319ba4SBarry Smith } 643b4319ba4SBarry Smith 6442e74eeadSLisandro Dalcin #undef __FUNCT__ 6452e74eeadSLisandro Dalcin #define __FUNCT__ "MatSetValues_IS" 6462e74eeadSLisandro Dalcin PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv) 6472e74eeadSLisandro Dalcin { 6482e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)mat->data; 6492e74eeadSLisandro Dalcin PetscInt rows_l[2048],cols_l[2048]; 6502e74eeadSLisandro Dalcin PetscErrorCode ierr; 6512e74eeadSLisandro Dalcin 6522e74eeadSLisandro Dalcin PetscFunctionBegin; 6532e74eeadSLisandro Dalcin #if defined(PETSC_USE_DEBUG) 654f23aa3ddSBarry Smith 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); 6552e74eeadSLisandro Dalcin #endif 6562e74eeadSLisandro Dalcin ierr = ISG2LMapApply(is->mapping,m,rows,rows_l);CHKERRQ(ierr); 6572e74eeadSLisandro Dalcin ierr = ISG2LMapApply(is->mapping,n,cols,cols_l);CHKERRQ(ierr); 6582e74eeadSLisandro Dalcin ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 6592e74eeadSLisandro Dalcin PetscFunctionReturn(0); 6602e74eeadSLisandro Dalcin } 6612e74eeadSLisandro Dalcin 6622e74eeadSLisandro Dalcin #undef ISG2LMapSetUp 6632e74eeadSLisandro Dalcin #undef ISG2LMapApply 664b4319ba4SBarry Smith 665b4319ba4SBarry Smith #undef __FUNCT__ 666b4319ba4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_IS" 66713f74950SBarry Smith PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 668b4319ba4SBarry Smith { 669dfbe8321SBarry Smith PetscErrorCode ierr; 670b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 671b4319ba4SBarry Smith 672b4319ba4SBarry Smith PetscFunctionBegin; 673b4319ba4SBarry Smith ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 674b4319ba4SBarry Smith PetscFunctionReturn(0); 675b4319ba4SBarry Smith } 676b4319ba4SBarry Smith 677b4319ba4SBarry Smith #undef __FUNCT__ 678f0006bf2SLisandro Dalcin #define __FUNCT__ "MatSetValuesBlockedLocal_IS" 679f0006bf2SLisandro Dalcin PetscErrorCode MatSetValuesBlockedLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 680f0006bf2SLisandro Dalcin { 681f0006bf2SLisandro Dalcin PetscErrorCode ierr; 682f0006bf2SLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 683f0006bf2SLisandro Dalcin 684f0006bf2SLisandro Dalcin PetscFunctionBegin; 685f0006bf2SLisandro Dalcin ierr = MatSetValuesBlocked(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 686f0006bf2SLisandro Dalcin PetscFunctionReturn(0); 687f0006bf2SLisandro Dalcin } 688f0006bf2SLisandro Dalcin 689f0006bf2SLisandro Dalcin #undef __FUNCT__ 6902e74eeadSLisandro Dalcin #define __FUNCT__ "MatZeroRows_IS" 6912b40b63fSBarry Smith PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 6922e74eeadSLisandro Dalcin { 6932e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 6940298fd71SBarry Smith PetscInt n_l = 0, *rows_l = NULL; 6952e74eeadSLisandro Dalcin PetscErrorCode ierr; 6962e74eeadSLisandro Dalcin 6972e74eeadSLisandro Dalcin PetscFunctionBegin; 698ce94432eSBarry Smith if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support"); 6992e74eeadSLisandro Dalcin if (n) { 700785e854fSJed Brown ierr = PetscMalloc1(n,&rows_l);CHKERRQ(ierr); 7012e74eeadSLisandro Dalcin ierr = ISGlobalToLocalMappingApply(is->mapping,IS_GTOLM_DROP,n,rows,&n_l,rows_l);CHKERRQ(ierr); 7022e74eeadSLisandro Dalcin } 7032b40b63fSBarry Smith ierr = MatZeroRowsLocal(A,n_l,rows_l,diag,x,b);CHKERRQ(ierr); 704c31cb41cSBarry Smith ierr = PetscFree(rows_l);CHKERRQ(ierr); 7052e74eeadSLisandro Dalcin PetscFunctionReturn(0); 7062e74eeadSLisandro Dalcin } 7072e74eeadSLisandro Dalcin 7082e74eeadSLisandro Dalcin #undef __FUNCT__ 709b4319ba4SBarry Smith #define __FUNCT__ "MatZeroRowsLocal_IS" 7102b40b63fSBarry Smith PetscErrorCode MatZeroRowsLocal_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 711b4319ba4SBarry Smith { 712b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 713dfbe8321SBarry Smith PetscErrorCode ierr; 714f4df32b1SMatthew Knepley PetscInt i; 715b4319ba4SBarry Smith PetscScalar *array; 716b4319ba4SBarry Smith 717b4319ba4SBarry Smith PetscFunctionBegin; 718ce94432eSBarry Smith if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support"); 719b4319ba4SBarry Smith { 720b4319ba4SBarry Smith /* 721b4319ba4SBarry Smith Set up is->x as a "counting vector". This is in order to MatMult_IS 722b4319ba4SBarry Smith work properly in the interface nodes. 723b4319ba4SBarry Smith */ 724b4319ba4SBarry Smith Vec counter; 725b4319ba4SBarry Smith PetscScalar one=1.0, zero=0.0; 7262a7a6963SBarry Smith ierr = MatCreateVecs(A,&counter,NULL);CHKERRQ(ierr); 7272dcb1b2aSMatthew Knepley ierr = VecSet(counter,zero);CHKERRQ(ierr); 7282dcb1b2aSMatthew Knepley ierr = VecSet(is->x,one);CHKERRQ(ierr); 729ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 730ca9f406cSSatish Balay ierr = VecScatterEnd (is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 731ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 732ca9f406cSSatish Balay ierr = VecScatterEnd (is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 7336bf464f9SBarry Smith ierr = VecDestroy(&counter);CHKERRQ(ierr); 734b4319ba4SBarry Smith } 735958c9bccSBarry Smith if (!n) { 736b4319ba4SBarry Smith is->pure_neumann = PETSC_TRUE; 737b4319ba4SBarry Smith } else { 738b4319ba4SBarry Smith is->pure_neumann = PETSC_FALSE; 7392205254eSKarl Rupp 740b4319ba4SBarry Smith ierr = VecGetArray(is->x,&array);CHKERRQ(ierr); 7412b40b63fSBarry Smith ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr); 742b4319ba4SBarry Smith for (i=0; i<n; i++) { 743f4df32b1SMatthew Knepley ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr); 744b4319ba4SBarry Smith } 745b4319ba4SBarry Smith ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 746b4319ba4SBarry Smith ierr = MatAssemblyEnd (is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 747b4319ba4SBarry Smith ierr = VecRestoreArray(is->x,&array);CHKERRQ(ierr); 748b4319ba4SBarry Smith } 749b4319ba4SBarry Smith PetscFunctionReturn(0); 750b4319ba4SBarry Smith } 751b4319ba4SBarry Smith 752b4319ba4SBarry Smith #undef __FUNCT__ 753b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyBegin_IS" 754dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type) 755b4319ba4SBarry Smith { 756b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 757dfbe8321SBarry Smith PetscErrorCode ierr; 758b4319ba4SBarry Smith 759b4319ba4SBarry Smith PetscFunctionBegin; 760b4319ba4SBarry Smith ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr); 761b4319ba4SBarry Smith PetscFunctionReturn(0); 762b4319ba4SBarry Smith } 763b4319ba4SBarry Smith 764b4319ba4SBarry Smith #undef __FUNCT__ 765b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_IS" 766dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type) 767b4319ba4SBarry Smith { 768b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 769dfbe8321SBarry Smith PetscErrorCode ierr; 770b4319ba4SBarry Smith 771b4319ba4SBarry Smith PetscFunctionBegin; 772b4319ba4SBarry Smith ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr); 773b4319ba4SBarry Smith PetscFunctionReturn(0); 774b4319ba4SBarry Smith } 775b4319ba4SBarry Smith 776b4319ba4SBarry Smith #undef __FUNCT__ 777b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat_IS" 7787087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local) 779b4319ba4SBarry Smith { 780b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)mat->data; 781b4319ba4SBarry Smith 782b4319ba4SBarry Smith PetscFunctionBegin; 783b4319ba4SBarry Smith *local = is->A; 784b4319ba4SBarry Smith PetscFunctionReturn(0); 785b4319ba4SBarry Smith } 786b4319ba4SBarry Smith 787b4319ba4SBarry Smith #undef __FUNCT__ 788b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat" 789b4319ba4SBarry Smith /*@ 790b4319ba4SBarry Smith MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix. 791b4319ba4SBarry Smith 792b4319ba4SBarry Smith Input Parameter: 793b4319ba4SBarry Smith . mat - the matrix 794b4319ba4SBarry Smith 795b4319ba4SBarry Smith Output Parameter: 796b4319ba4SBarry Smith . local - the local matrix usually MATSEQAIJ 797b4319ba4SBarry Smith 798b4319ba4SBarry Smith Level: advanced 799b4319ba4SBarry Smith 800b4319ba4SBarry Smith Notes: 801b4319ba4SBarry Smith This can be called if you have precomputed the nonzero structure of the 802b4319ba4SBarry Smith matrix and want to provide it to the inner matrix object to improve the performance 803b4319ba4SBarry Smith of the MatSetValues() operation. 804b4319ba4SBarry Smith 805b4319ba4SBarry Smith .seealso: MATIS 806b4319ba4SBarry Smith @*/ 8077087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local) 808b4319ba4SBarry Smith { 8094ac538c5SBarry Smith PetscErrorCode ierr; 810b4319ba4SBarry Smith 811b4319ba4SBarry Smith PetscFunctionBegin; 8120700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 813b4319ba4SBarry Smith PetscValidPointer(local,2); 8144ac538c5SBarry Smith ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat*),(mat,local));CHKERRQ(ierr); 815b4319ba4SBarry Smith PetscFunctionReturn(0); 816b4319ba4SBarry Smith } 817b4319ba4SBarry Smith 8183b03a366Sstefano_zampini #undef __FUNCT__ 8193b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat_IS" 8203b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat_IS(Mat mat,Mat local) 8213b03a366Sstefano_zampini { 8223b03a366Sstefano_zampini Mat_IS *is = (Mat_IS*)mat->data; 8233b03a366Sstefano_zampini PetscInt nrows,ncols,orows,ocols; 8243b03a366Sstefano_zampini PetscErrorCode ierr; 8253b03a366Sstefano_zampini 8263b03a366Sstefano_zampini PetscFunctionBegin; 8274e4c7dbeSStefano Zampini if (is->A) { 8283b03a366Sstefano_zampini ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr); 8293b03a366Sstefano_zampini ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr); 830f23aa3ddSBarry Smith 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); 8314e4c7dbeSStefano Zampini } 8323b03a366Sstefano_zampini ierr = PetscObjectReference((PetscObject)local);CHKERRQ(ierr); 8333b03a366Sstefano_zampini ierr = MatDestroy(&is->A);CHKERRQ(ierr); 8343b03a366Sstefano_zampini is->A = local; 8353b03a366Sstefano_zampini PetscFunctionReturn(0); 8363b03a366Sstefano_zampini } 8373b03a366Sstefano_zampini 8383b03a366Sstefano_zampini #undef __FUNCT__ 8393b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat" 8403b03a366Sstefano_zampini /*@ 8413b03a366Sstefano_zampini MatISSetLocalMat - Set the local matrix stored inside a MATIS matrix. 8423b03a366Sstefano_zampini 8433b03a366Sstefano_zampini Input Parameter: 8443b03a366Sstefano_zampini . mat - the matrix 8453b03a366Sstefano_zampini . local - the local matrix usually MATSEQAIJ 8463b03a366Sstefano_zampini 8473b03a366Sstefano_zampini Output Parameter: 8483b03a366Sstefano_zampini 8493b03a366Sstefano_zampini Level: advanced 8503b03a366Sstefano_zampini 8513b03a366Sstefano_zampini Notes: 8523b03a366Sstefano_zampini This can be called if you have precomputed the local matrix and 8533b03a366Sstefano_zampini want to provide it to the matrix object MATIS. 8543b03a366Sstefano_zampini 8553b03a366Sstefano_zampini .seealso: MATIS 8563b03a366Sstefano_zampini @*/ 8573b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat(Mat mat,Mat local) 8583b03a366Sstefano_zampini { 8593b03a366Sstefano_zampini PetscErrorCode ierr; 8603b03a366Sstefano_zampini 8613b03a366Sstefano_zampini PetscFunctionBegin; 8623b03a366Sstefano_zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 863b7ce53b6SStefano Zampini PetscValidHeaderSpecific(local,MAT_CLASSID,2); 8643b03a366Sstefano_zampini ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr); 8653b03a366Sstefano_zampini PetscFunctionReturn(0); 8663b03a366Sstefano_zampini } 8673b03a366Sstefano_zampini 8686726f965SBarry Smith #undef __FUNCT__ 8696726f965SBarry Smith #define __FUNCT__ "MatZeroEntries_IS" 8706726f965SBarry Smith PetscErrorCode MatZeroEntries_IS(Mat A) 8716726f965SBarry Smith { 8726726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 8736726f965SBarry Smith PetscErrorCode ierr; 8746726f965SBarry Smith 8756726f965SBarry Smith PetscFunctionBegin; 8766726f965SBarry Smith ierr = MatZeroEntries(a->A);CHKERRQ(ierr); 8776726f965SBarry Smith PetscFunctionReturn(0); 8786726f965SBarry Smith } 8796726f965SBarry Smith 8806726f965SBarry Smith #undef __FUNCT__ 8812e74eeadSLisandro Dalcin #define __FUNCT__ "MatScale_IS" 8822e74eeadSLisandro Dalcin PetscErrorCode MatScale_IS(Mat A,PetscScalar a) 8832e74eeadSLisandro Dalcin { 8842e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 8852e74eeadSLisandro Dalcin PetscErrorCode ierr; 8862e74eeadSLisandro Dalcin 8872e74eeadSLisandro Dalcin PetscFunctionBegin; 8882e74eeadSLisandro Dalcin ierr = MatScale(is->A,a);CHKERRQ(ierr); 8892e74eeadSLisandro Dalcin PetscFunctionReturn(0); 8902e74eeadSLisandro Dalcin } 8912e74eeadSLisandro Dalcin 8922e74eeadSLisandro Dalcin #undef __FUNCT__ 8932e74eeadSLisandro Dalcin #define __FUNCT__ "MatGetDiagonal_IS" 8942e74eeadSLisandro Dalcin PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v) 8952e74eeadSLisandro Dalcin { 8962e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 8972e74eeadSLisandro Dalcin PetscErrorCode ierr; 8982e74eeadSLisandro Dalcin 8992e74eeadSLisandro Dalcin PetscFunctionBegin; 9002e74eeadSLisandro Dalcin /* get diagonal of the local matrix */ 9012e74eeadSLisandro Dalcin ierr = MatGetDiagonal(is->A,is->x);CHKERRQ(ierr); 9022e74eeadSLisandro Dalcin 9032e74eeadSLisandro Dalcin /* scatter diagonal back into global vector */ 9042e74eeadSLisandro Dalcin ierr = VecSet(v,0);CHKERRQ(ierr); 905ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 906ca9f406cSSatish Balay ierr = VecScatterEnd (is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 9072e74eeadSLisandro Dalcin PetscFunctionReturn(0); 9082e74eeadSLisandro Dalcin } 9092e74eeadSLisandro Dalcin 9102e74eeadSLisandro Dalcin #undef __FUNCT__ 9116726f965SBarry Smith #define __FUNCT__ "MatSetOption_IS" 912ace3abfcSBarry Smith PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg) 9136726f965SBarry Smith { 9146726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 9156726f965SBarry Smith PetscErrorCode ierr; 9166726f965SBarry Smith 9176726f965SBarry Smith PetscFunctionBegin; 9184e0d8c25SBarry Smith ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr); 9196726f965SBarry Smith PetscFunctionReturn(0); 9206726f965SBarry Smith } 9216726f965SBarry Smith 922284134d9SBarry Smith #undef __FUNCT__ 923284134d9SBarry Smith #define __FUNCT__ "MatCreateIS" 924284134d9SBarry Smith /*@ 925284134d9SBarry Smith MatCreateIS - Creates a "process" unassmembled matrix, it is assembled on each 926284134d9SBarry Smith process but not across processes. 927284134d9SBarry Smith 928284134d9SBarry Smith Input Parameters: 929284134d9SBarry Smith + comm - MPI communicator that will share the matrix 9304e4c7dbeSStefano Zampini . bs - local and global block size of the matrix 931df3898eeSBarry Smith . m,n,M,N - local and/or global sizes of the left and right vector used in matrix vector products 932284134d9SBarry Smith - map - mapping that defines the global number for each local number 933284134d9SBarry Smith 934284134d9SBarry Smith Output Parameter: 935284134d9SBarry Smith . A - the resulting matrix 936284134d9SBarry Smith 9378e6c10adSSatish Balay Level: advanced 9388e6c10adSSatish Balay 939284134d9SBarry Smith Notes: See MATIS for more details 9408cda6cd7SBarry Smith m and n are NOT related to the size of the map, they are the size of the part of the vector owned 9418cda6cd7SBarry Smith by that process. m + nghosts (or n + nghosts) is the length of map since map maps all local points 9428cda6cd7SBarry Smith plus the ghost points to global indices. 943284134d9SBarry Smith 944284134d9SBarry Smith .seealso: MATIS, MatSetLocalToGlobalMapping() 945284134d9SBarry Smith @*/ 9464e4c7dbeSStefano Zampini PetscErrorCode MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping map,Mat *A) 947284134d9SBarry Smith { 948284134d9SBarry Smith PetscErrorCode ierr; 949284134d9SBarry Smith 950284134d9SBarry Smith PetscFunctionBegin; 951284134d9SBarry Smith ierr = MatCreate(comm,A);CHKERRQ(ierr); 952d3ee2243SStefano Zampini ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr); 953284134d9SBarry Smith ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 954284134d9SBarry Smith ierr = MatSetType(*A,MATIS);CHKERRQ(ierr); 9553b03a366Sstefano_zampini ierr = MatSetUp(*A);CHKERRQ(ierr); 956784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(*A,map,map);CHKERRQ(ierr); 957284134d9SBarry Smith PetscFunctionReturn(0); 958284134d9SBarry Smith } 959284134d9SBarry Smith 960b4319ba4SBarry Smith /*MC 9616a818285SBarry Smith MATIS - MATIS = "is" - A matrix type to be used for using the Neumann-Neumann type preconditioners, see PCBDDC. 962b4319ba4SBarry Smith This stores the matrices in globally unassembled form. Each processor 963b4319ba4SBarry Smith assembles only its local Neumann problem and the parallel matrix vector 964b4319ba4SBarry Smith product is handled "implicitly". 965b4319ba4SBarry Smith 966b4319ba4SBarry Smith Operations Provided: 9676726f965SBarry Smith + MatMult() 9682e74eeadSLisandro Dalcin . MatMultAdd() 9692e74eeadSLisandro Dalcin . MatMultTranspose() 9702e74eeadSLisandro Dalcin . MatMultTransposeAdd() 9716726f965SBarry Smith . MatZeroEntries() 9726726f965SBarry Smith . MatSetOption() 9732e74eeadSLisandro Dalcin . MatZeroRows() 9746726f965SBarry Smith . MatZeroRowsLocal() 9752e74eeadSLisandro Dalcin . MatSetValues() 9766726f965SBarry Smith . MatSetValuesLocal() 9772e74eeadSLisandro Dalcin . MatScale() 9782e74eeadSLisandro Dalcin . MatGetDiagonal() 9796726f965SBarry Smith - MatSetLocalToGlobalMapping() 980b4319ba4SBarry Smith 981b4319ba4SBarry Smith Options Database Keys: 982b4319ba4SBarry Smith . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions() 983b4319ba4SBarry Smith 984b4319ba4SBarry Smith Notes: Options prefix for the inner matrix are given by -is_mat_xxx 985b4319ba4SBarry Smith 986b4319ba4SBarry Smith You must call MatSetLocalToGlobalMapping() before using this matrix type. 987b4319ba4SBarry Smith 988b4319ba4SBarry Smith You can do matrix preallocation on the local matrix after you obtain it with 989b4319ba4SBarry Smith MatISGetLocalMat() 990b4319ba4SBarry Smith 991b4319ba4SBarry Smith Level: advanced 992b4319ba4SBarry Smith 9936a818285SBarry Smith .seealso: PC, MatISGetLocalMat(), MatSetLocalToGlobalMapping(), PCBDDC 994b4319ba4SBarry Smith 995b4319ba4SBarry Smith M*/ 996b4319ba4SBarry Smith 997b4319ba4SBarry Smith #undef __FUNCT__ 998b4319ba4SBarry Smith #define __FUNCT__ "MatCreate_IS" 9998cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_IS(Mat A) 1000b4319ba4SBarry Smith { 1001dfbe8321SBarry Smith PetscErrorCode ierr; 1002b4319ba4SBarry Smith Mat_IS *b; 1003b4319ba4SBarry Smith 1004b4319ba4SBarry Smith PetscFunctionBegin; 1005b00a9115SJed Brown ierr = PetscNewLog(A,&b);CHKERRQ(ierr); 1006b4319ba4SBarry Smith A->data = (void*)b; 1007b4319ba4SBarry Smith ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1008b4319ba4SBarry Smith 1009b4319ba4SBarry Smith A->ops->mult = MatMult_IS; 10102e74eeadSLisandro Dalcin A->ops->multadd = MatMultAdd_IS; 10112e74eeadSLisandro Dalcin A->ops->multtranspose = MatMultTranspose_IS; 10122e74eeadSLisandro Dalcin A->ops->multtransposeadd = MatMultTransposeAdd_IS; 1013b4319ba4SBarry Smith A->ops->destroy = MatDestroy_IS; 1014b4319ba4SBarry Smith A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS; 10152e74eeadSLisandro Dalcin A->ops->setvalues = MatSetValues_IS; 1016b4319ba4SBarry Smith A->ops->setvalueslocal = MatSetValuesLocal_IS; 1017f0006bf2SLisandro Dalcin A->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_IS; 10182e74eeadSLisandro Dalcin A->ops->zerorows = MatZeroRows_IS; 1019b4319ba4SBarry Smith A->ops->zerorowslocal = MatZeroRowsLocal_IS; 1020b4319ba4SBarry Smith A->ops->assemblybegin = MatAssemblyBegin_IS; 1021b4319ba4SBarry Smith A->ops->assemblyend = MatAssemblyEnd_IS; 1022b4319ba4SBarry Smith A->ops->view = MatView_IS; 10236726f965SBarry Smith A->ops->zeroentries = MatZeroEntries_IS; 10242e74eeadSLisandro Dalcin A->ops->scale = MatScale_IS; 10252e74eeadSLisandro Dalcin A->ops->getdiagonal = MatGetDiagonal_IS; 10266726f965SBarry Smith A->ops->setoption = MatSetOption_IS; 102769796d55SStefano Zampini A->ops->ishermitian = MatIsHermitian_IS; 102869796d55SStefano Zampini A->ops->issymmetric = MatIsSymmetric_IS; 1029ad6194a2SStefano Zampini A->ops->duplicate = MatDuplicate_IS; 1030b4319ba4SBarry Smith 103126283091SBarry Smith ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 103226283091SBarry Smith ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 1033b4319ba4SBarry Smith 1034b7ce53b6SStefano Zampini /* zeros pointer */ 1035b4319ba4SBarry Smith b->A = 0; 1036b4319ba4SBarry Smith b->ctx = 0; 1037b4319ba4SBarry Smith b->x = 0; 1038b4319ba4SBarry Smith b->y = 0; 1039b09366fdSStefano Zampini b->mapping = 0; 104028f4e0baSStefano Zampini b->sf = 0; 104128f4e0baSStefano Zampini b->sf_rootdata = 0; 104228f4e0baSStefano Zampini b->sf_leafdata = 0; 1043b7ce53b6SStefano Zampini 1044b7ce53b6SStefano Zampini /* special MATIS functions */ 1045bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",MatISGetLocalMat_IS);CHKERRQ(ierr); 1046bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",MatISSetLocalMat_IS);CHKERRQ(ierr); 1047b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",MatISGetMPIXAIJ_IS);CHKERRQ(ierr); 10482e1947a5SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",MatISSetPreallocation_IS);CHKERRQ(ierr); 104917667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr); 1050b4319ba4SBarry Smith PetscFunctionReturn(0); 1051b4319ba4SBarry Smith } 1052