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*/ 16*28f4e0baSStefano Zampini 17*28f4e0baSStefano Zampini #undef __FUNCT__ 18*28f4e0baSStefano Zampini #define __FUNCT__ "MatISComputeSF_Private" 19*28f4e0baSStefano Zampini PetscErrorCode MatISComputeSF_Private(Mat B) 20*28f4e0baSStefano Zampini { 21*28f4e0baSStefano Zampini Mat_IS *matis = (Mat_IS*)(B->data); 22*28f4e0baSStefano Zampini const PetscInt *gidxs; 23*28f4e0baSStefano Zampini PetscErrorCode ierr; 24*28f4e0baSStefano Zampini 25*28f4e0baSStefano Zampini PetscFunctionBegin; 26*28f4e0baSStefano Zampini ierr = MatGetSize(matis->A,&matis->sf_nleaves,NULL);CHKERRQ(ierr); 27*28f4e0baSStefano Zampini ierr = MatGetLocalSize(B,&matis->sf_nroots,NULL);CHKERRQ(ierr); 28*28f4e0baSStefano Zampini ierr = PetscSFCreate(PetscObjectComm((PetscObject)B),&matis->sf);CHKERRQ(ierr); 29*28f4e0baSStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(matis->mapping,&gidxs);CHKERRQ(ierr); 30*28f4e0baSStefano Zampini /* PETSC_OWN_POINTER refers to ilocal which is NULL */ 31*28f4e0baSStefano Zampini ierr = PetscSFSetGraphLayout(matis->sf,B->rmap,matis->sf_nleaves,NULL,PETSC_OWN_POINTER,gidxs);CHKERRQ(ierr); 32*28f4e0baSStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(matis->mapping,&gidxs);CHKERRQ(ierr); 33*28f4e0baSStefano Zampini ierr = PetscMalloc2(matis->sf_nroots,&matis->sf_rootdata,matis->sf_nleaves,&matis->sf_leafdata);CHKERRQ(ierr); 34*28f4e0baSStefano Zampini PetscFunctionReturn(0); 35*28f4e0baSStefano Zampini } 362e1947a5SStefano Zampini 372e1947a5SStefano Zampini #undef __FUNCT__ 382e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation" 39a88811baSStefano Zampini /* 40a88811baSStefano Zampini MatISSetPreallocation - Preallocates memory for a MATIS parallel matrix. 41a88811baSStefano Zampini 42a88811baSStefano Zampini Collective on MPI_Comm 43a88811baSStefano Zampini 44a88811baSStefano Zampini Input Parameters: 45a88811baSStefano Zampini + B - the matrix 46a88811baSStefano Zampini . d_nz - number of nonzeros per row in DIAGONAL portion of local submatrix 47a88811baSStefano Zampini (same value is used for all local rows) 48a88811baSStefano Zampini . d_nnz - array containing the number of nonzeros in the various rows of the 49a88811baSStefano Zampini DIAGONAL portion of the local submatrix (possibly different for each row) 50a88811baSStefano Zampini or NULL, if d_nz is used to specify the nonzero structure. 51a88811baSStefano Zampini The size of this array is equal to the number of local rows, i.e 'm'. 52a88811baSStefano Zampini For matrices that will be factored, you must leave room for (and set) 53a88811baSStefano Zampini the diagonal entry even if it is zero. 54a88811baSStefano Zampini . o_nz - number of nonzeros per row in the OFF-DIAGONAL portion of local 55a88811baSStefano Zampini submatrix (same value is used for all local rows). 56a88811baSStefano Zampini - o_nnz - array containing the number of nonzeros in the various rows of the 57a88811baSStefano Zampini OFF-DIAGONAL portion of the local submatrix (possibly different for 58a88811baSStefano Zampini each row) or NULL, if o_nz is used to specify the nonzero 59a88811baSStefano Zampini structure. The size of this array is equal to the number 60a88811baSStefano Zampini of local rows, i.e 'm'. 61a88811baSStefano Zampini 62a88811baSStefano Zampini If the *_nnz parameter is given then the *_nz parameter is ignored 63a88811baSStefano Zampini 64a88811baSStefano Zampini Level: intermediate 65a88811baSStefano Zampini 66a88811baSStefano Zampini Notes: This function has the same interface as the MPIAIJ preallocation routine in order to simplify the transition 67a88811baSStefano Zampini from the asssembled format to the unassembled one. It overestimates the preallocation of MATIS local 68a88811baSStefano Zampini matrices; for exact preallocation, the user should set the preallocation directly on local matrix objects. 69a88811baSStefano Zampini 70a88811baSStefano Zampini .keywords: matrix 71a88811baSStefano Zampini 72a88811baSStefano Zampini .seealso: MatCreate(), MatCreateIS(), MatMPIAIJSetPreallocation(), MatISGetLocalMat() 73a88811baSStefano Zampini @*/ 742e1947a5SStefano Zampini PetscErrorCode MatISSetPreallocation(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[]) 752e1947a5SStefano Zampini { 762e1947a5SStefano Zampini PetscErrorCode ierr; 772e1947a5SStefano Zampini 782e1947a5SStefano Zampini PetscFunctionBegin; 792e1947a5SStefano Zampini PetscValidHeaderSpecific(B,MAT_CLASSID,1); 802e1947a5SStefano Zampini PetscValidType(B,1); 812e1947a5SStefano Zampini ierr = PetscTryMethod(B,"MatISSetPreallocation_C",(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(B,d_nz,d_nnz,o_nz,o_nnz));CHKERRQ(ierr); 822e1947a5SStefano Zampini PetscFunctionReturn(0); 832e1947a5SStefano Zampini } 842e1947a5SStefano Zampini 852e1947a5SStefano Zampini #undef __FUNCT__ 862e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation_IS" 872e1947a5SStefano Zampini PetscErrorCode MatISSetPreallocation_IS(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[]) 882e1947a5SStefano Zampini { 892e1947a5SStefano Zampini Mat_IS *matis = (Mat_IS*)(B->data); 90*28f4e0baSStefano Zampini PetscInt bs,i,nlocalcols; 912e1947a5SStefano Zampini PetscErrorCode ierr; 922e1947a5SStefano Zampini 932e1947a5SStefano Zampini PetscFunctionBegin; 942e1947a5SStefano Zampini if (!matis->A) { 952e1947a5SStefano Zampini SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"You should first call MatSetLocalToGlobalMapping"); 962e1947a5SStefano Zampini } 97*28f4e0baSStefano Zampini if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */ 98*28f4e0baSStefano Zampini ierr = MatISComputeSF_Private(B);CHKERRQ(ierr); 99*28f4e0baSStefano Zampini } 1002e1947a5SStefano Zampini if (!d_nnz) { 101*28f4e0baSStefano Zampini for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] = d_nz; 1022e1947a5SStefano Zampini } else { 103*28f4e0baSStefano Zampini for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] = d_nnz[i]; 1042e1947a5SStefano Zampini } 1052e1947a5SStefano Zampini if (!o_nnz) { 106*28f4e0baSStefano Zampini for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] += o_nz; 1072e1947a5SStefano Zampini } else { 108*28f4e0baSStefano Zampini for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] += o_nnz[i]; 1092e1947a5SStefano Zampini } 110*28f4e0baSStefano Zampini ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 111*28f4e0baSStefano Zampini ierr = MatGetSize(matis->A,NULL,&nlocalcols);CHKERRQ(ierr); 112*28f4e0baSStefano Zampini ierr = MatGetBlockSize(matis->A,&bs);CHKERRQ(ierr); 113*28f4e0baSStefano Zampini ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 114*28f4e0baSStefano Zampini for (i=0;i<matis->sf_nleaves;i++) { 115*28f4e0baSStefano Zampini matis->sf_leafdata[i] = PetscMin(matis->sf_leafdata[i],nlocalcols); 1162e1947a5SStefano Zampini } 117*28f4e0baSStefano Zampini ierr = MatSeqAIJSetPreallocation(matis->A,0,matis->sf_leafdata);CHKERRQ(ierr); 118*28f4e0baSStefano Zampini for (i=0;i<matis->sf_nleaves/bs;i++) { 119*28f4e0baSStefano Zampini matis->sf_leafdata[i] = matis->sf_leafdata[i*bs]/bs; 1202e1947a5SStefano Zampini } 121*28f4e0baSStefano Zampini ierr = MatSeqBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr); 122*28f4e0baSStefano Zampini for (i=0;i<matis->sf_nleaves/bs;i++) { 123*28f4e0baSStefano Zampini matis->sf_leafdata[i] = matis->sf_leafdata[i]-i; 1242e1947a5SStefano Zampini } 125*28f4e0baSStefano Zampini ierr = MatSeqSBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr); 1262e1947a5SStefano Zampini PetscFunctionReturn(0); 1272e1947a5SStefano Zampini } 128b4319ba4SBarry Smith 129b4319ba4SBarry Smith #undef __FUNCT__ 130b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ_IS" 131b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ_IS(Mat mat, MatReuse reuse, Mat *M) 132b7ce53b6SStefano Zampini { 133b7ce53b6SStefano Zampini Mat_IS *matis = (Mat_IS*)(mat->data); 134b7ce53b6SStefano Zampini /* info on mat */ 135b7ce53b6SStefano Zampini /* ISLocalToGlobalMapping rmapping,cmapping; */ 136b7ce53b6SStefano Zampini PetscInt bs,rows,cols; 137b7ce53b6SStefano Zampini PetscInt lrows,lcols; 138b7ce53b6SStefano Zampini PetscInt local_rows,local_cols; 139b7ce53b6SStefano Zampini PetscBool isdense,issbaij,issbaij_red; 1407c03b4e8SStefano Zampini PetscMPIInt nsubdomains; 141b7ce53b6SStefano Zampini /* values insertion */ 142b7ce53b6SStefano Zampini PetscScalar *array; 143b7ce53b6SStefano Zampini PetscInt *local_indices,*global_indices; 144b7ce53b6SStefano Zampini /* work */ 145b7ce53b6SStefano Zampini PetscInt i,j,index_row; 146b7ce53b6SStefano Zampini PetscErrorCode ierr; 147b7ce53b6SStefano Zampini 148b7ce53b6SStefano Zampini PetscFunctionBegin; 149b7ce53b6SStefano Zampini /* MISSING CHECKS 150b7ce53b6SStefano Zampini - rectangular case not covered (it is not allowed by MATIS) 151b7ce53b6SStefano Zampini */ 152b7ce53b6SStefano Zampini /* get info from mat */ 1537c03b4e8SStefano Zampini ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&nsubdomains);CHKERRQ(ierr); 1547c03b4e8SStefano Zampini if (nsubdomains == 1) { 1557c03b4e8SStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { 1567c03b4e8SStefano Zampini ierr = MatDuplicate(matis->A,MAT_COPY_VALUES,&(*M));CHKERRQ(ierr); 1577c03b4e8SStefano Zampini } else { 1587c03b4e8SStefano Zampini ierr = MatCopy(matis->A,*M,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 1597c03b4e8SStefano Zampini } 1607c03b4e8SStefano Zampini PetscFunctionReturn(0); 1617c03b4e8SStefano Zampini } 162b7ce53b6SStefano Zampini /* ierr = MatGetLocalToGlobalMapping(mat,&rmapping,&cmapping);CHKERRQ(ierr); */ 163b7ce53b6SStefano Zampini ierr = MatGetSize(mat,&rows,&cols);CHKERRQ(ierr); 164b7ce53b6SStefano Zampini ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 165b7ce53b6SStefano Zampini ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr); 166b7ce53b6SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr); 167b7ce53b6SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr); 168b7ce53b6SStefano Zampini 169b7ce53b6SStefano Zampini /* work */ 170b7ce53b6SStefano Zampini ierr = PetscMalloc1(local_rows,&local_indices);CHKERRQ(ierr); 171b7ce53b6SStefano Zampini for (i=0;i<local_rows;i++) local_indices[i]=i; 172b7ce53b6SStefano Zampini /* map indices of local mat to global values */ 173854ce69bSBarry Smith ierr = PetscMalloc1(PetscMax(local_cols,local_rows),&global_indices);CHKERRQ(ierr); 174b7ce53b6SStefano Zampini /* ierr = ISLocalToGlobalMappingApply(rmapping,local_rows,local_indices,global_indices);CHKERRQ(ierr); */ 175b7ce53b6SStefano Zampini ierr = ISLocalToGlobalMappingApply(matis->mapping,local_rows,local_indices,global_indices);CHKERRQ(ierr); 176b7ce53b6SStefano Zampini 177b7ce53b6SStefano Zampini if (issbaij) { 178b7ce53b6SStefano Zampini ierr = MatGetRowUpperTriangular(matis->A);CHKERRQ(ierr); 179b7ce53b6SStefano Zampini } 180b7ce53b6SStefano Zampini 181b7ce53b6SStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { 182b7ce53b6SStefano Zampini Mat new_mat; 183b7ce53b6SStefano Zampini MatType new_mat_type; 184b7ce53b6SStefano Zampini Vec vec_dnz,vec_onz; 185b7ce53b6SStefano Zampini PetscScalar *my_dnz,*my_onz; 186b7ce53b6SStefano Zampini PetscInt *dnz,*onz,*mat_ranges,*row_ownership; 187b7ce53b6SStefano Zampini PetscInt index_col,owner; 188b7ce53b6SStefano Zampini 189b7ce53b6SStefano Zampini /* determining new matrix type */ 190b7ce53b6SStefano Zampini ierr = MPI_Allreduce(&issbaij,&issbaij_red,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 191b7ce53b6SStefano Zampini if (issbaij_red) { 192b7ce53b6SStefano Zampini new_mat_type = MATSBAIJ; 193b7ce53b6SStefano Zampini } else { 194b7ce53b6SStefano Zampini if (bs>1) { 195b7ce53b6SStefano Zampini new_mat_type = MATBAIJ; 196b7ce53b6SStefano Zampini } else { 197b7ce53b6SStefano Zampini new_mat_type = MATAIJ; 198b7ce53b6SStefano Zampini } 199b7ce53b6SStefano Zampini } 200b7ce53b6SStefano Zampini 201b7ce53b6SStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)mat),&new_mat);CHKERRQ(ierr); 202b7ce53b6SStefano Zampini ierr = MatSetSizes(new_mat,PETSC_DECIDE,PETSC_DECIDE,rows,cols);CHKERRQ(ierr); 203b7ce53b6SStefano Zampini ierr = MatSetBlockSize(new_mat,bs);CHKERRQ(ierr); 204b7ce53b6SStefano Zampini ierr = MatSetType(new_mat,new_mat_type);CHKERRQ(ierr); 205b7ce53b6SStefano Zampini ierr = MatSetUp(new_mat);CHKERRQ(ierr); 206b7ce53b6SStefano Zampini ierr = MatGetLocalSize(new_mat,&lrows,&lcols);CHKERRQ(ierr); 207b7ce53b6SStefano Zampini 208b7ce53b6SStefano Zampini /* 209b7ce53b6SStefano Zampini preallocation 210b7ce53b6SStefano Zampini */ 211b7ce53b6SStefano Zampini 212b7ce53b6SStefano Zampini ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)new_mat),lrows,lcols,dnz,onz);CHKERRQ(ierr); 213b7ce53b6SStefano Zampini /* 214b7ce53b6SStefano Zampini Some vectors are needed to sum up properly on shared interface dofs. 215b7ce53b6SStefano Zampini Preallocation macros cannot do the job. 216b7ce53b6SStefano Zampini Note that preallocation is not exact, since it overestimates nonzeros 217b7ce53b6SStefano Zampini */ 2182a7a6963SBarry Smith ierr = MatCreateVecs(new_mat,NULL,&vec_dnz);CHKERRQ(ierr); 219b7ce53b6SStefano Zampini /* ierr = VecSetLocalToGlobalMapping(vec_dnz,rmapping);CHKERRQ(ierr); */ 220b7ce53b6SStefano Zampini ierr = VecSetLocalToGlobalMapping(vec_dnz,matis->mapping);CHKERRQ(ierr); 221b7ce53b6SStefano Zampini ierr = VecDuplicate(vec_dnz,&vec_onz);CHKERRQ(ierr); 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 */ 235b7ce53b6SStefano Zampini ierr = PetscMalloc1(local_rows,&my_dnz);CHKERRQ(ierr); 236b7ce53b6SStefano Zampini ierr = PetscMalloc1(local_rows,&my_onz);CHKERRQ(ierr); 237b7ce53b6SStefano Zampini ierr = PetscMemzero(my_dnz,local_rows*sizeof(*my_dnz));CHKERRQ(ierr); 238b7ce53b6SStefano Zampini ierr = PetscMemzero(my_onz,local_rows*sizeof(*my_onz));CHKERRQ(ierr); 239b7ce53b6SStefano Zampini /* preallocation as a MATAIJ */ 240b7ce53b6SStefano Zampini if (isdense) { /* special case for dense local matrices */ 241b7ce53b6SStefano Zampini for (i=0;i<local_rows;i++) { 242b7ce53b6SStefano Zampini index_row = global_indices[i]; 243b7ce53b6SStefano Zampini for (j=i;j<local_rows;j++) { 244b7ce53b6SStefano Zampini owner = row_ownership[index_row]; 245b7ce53b6SStefano Zampini index_col = global_indices[j]; 246b7ce53b6SStefano Zampini if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */ 247b7ce53b6SStefano Zampini my_dnz[i] += 1.0; 248b7ce53b6SStefano Zampini } else { /* offdiag block */ 249b7ce53b6SStefano Zampini my_onz[i] += 1.0; 250b7ce53b6SStefano Zampini } 251b7ce53b6SStefano Zampini /* same as before, interchanging rows and cols */ 252b7ce53b6SStefano Zampini if (i != j) { 253b7ce53b6SStefano Zampini owner = row_ownership[index_col]; 254b7ce53b6SStefano Zampini if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) { 255b7ce53b6SStefano Zampini my_dnz[j] += 1.0; 256b7ce53b6SStefano Zampini } else { 257b7ce53b6SStefano Zampini my_onz[j] += 1.0; 258b7ce53b6SStefano Zampini } 259b7ce53b6SStefano Zampini } 260b7ce53b6SStefano Zampini } 261b7ce53b6SStefano Zampini } 262b7ce53b6SStefano Zampini } else { 263b7ce53b6SStefano Zampini for (i=0;i<local_rows;i++) { 264b7ce53b6SStefano Zampini PetscInt ncols; 265b7ce53b6SStefano Zampini const PetscInt *cols; 266b7ce53b6SStefano Zampini index_row = global_indices[i]; 267b7ce53b6SStefano Zampini ierr = MatGetRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr); 268b7ce53b6SStefano Zampini for (j=0;j<ncols;j++) { 269b7ce53b6SStefano Zampini owner = row_ownership[index_row]; 270b7ce53b6SStefano Zampini index_col = global_indices[cols[j]]; 271b7ce53b6SStefano Zampini if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */ 272b7ce53b6SStefano Zampini my_dnz[i] += 1.0; 273b7ce53b6SStefano Zampini } else { /* offdiag block */ 274b7ce53b6SStefano Zampini my_onz[i] += 1.0; 275b7ce53b6SStefano Zampini } 276b7ce53b6SStefano Zampini /* same as before, interchanging rows and cols */ 277b7ce53b6SStefano Zampini if (issbaij) { 278b7ce53b6SStefano Zampini owner = row_ownership[index_col]; 279b7ce53b6SStefano Zampini if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) { 280b7ce53b6SStefano Zampini my_dnz[j] += 1.0; 281b7ce53b6SStefano Zampini } else { 282b7ce53b6SStefano Zampini my_onz[j] += 1.0; 283b7ce53b6SStefano Zampini } 284b7ce53b6SStefano Zampini } 285b7ce53b6SStefano Zampini } 286b7ce53b6SStefano Zampini ierr = MatRestoreRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr); 287b7ce53b6SStefano Zampini } 288b7ce53b6SStefano Zampini } 289b7ce53b6SStefano Zampini ierr = VecSet(vec_dnz,0.0);CHKERRQ(ierr); 290b7ce53b6SStefano Zampini ierr = VecSet(vec_onz,0.0);CHKERRQ(ierr); 291b7ce53b6SStefano Zampini if (local_rows) { /* multilevel guard */ 292b7ce53b6SStefano Zampini ierr = VecSetValuesLocal(vec_dnz,local_rows,local_indices,my_dnz,ADD_VALUES);CHKERRQ(ierr); 293b7ce53b6SStefano Zampini ierr = VecSetValuesLocal(vec_onz,local_rows,local_indices,my_onz,ADD_VALUES);CHKERRQ(ierr); 294b7ce53b6SStefano Zampini } 295b7ce53b6SStefano Zampini ierr = VecAssemblyBegin(vec_dnz);CHKERRQ(ierr); 296b7ce53b6SStefano Zampini ierr = VecAssemblyBegin(vec_onz);CHKERRQ(ierr); 297b7ce53b6SStefano Zampini ierr = VecAssemblyEnd(vec_dnz);CHKERRQ(ierr); 298b7ce53b6SStefano Zampini ierr = VecAssemblyEnd(vec_onz);CHKERRQ(ierr); 299b7ce53b6SStefano Zampini ierr = PetscFree(my_dnz);CHKERRQ(ierr); 300b7ce53b6SStefano Zampini ierr = PetscFree(my_onz);CHKERRQ(ierr); 301b7ce53b6SStefano Zampini ierr = PetscFree(row_ownership);CHKERRQ(ierr); 302b7ce53b6SStefano Zampini 303b7ce53b6SStefano Zampini /* set computed preallocation in dnz and onz */ 304b7ce53b6SStefano Zampini ierr = VecGetArray(vec_dnz,&array);CHKERRQ(ierr); 305b7ce53b6SStefano Zampini for (i=0; i<lrows; i++) dnz[i] = (PetscInt)PetscRealPart(array[i]); 306b7ce53b6SStefano Zampini ierr = VecRestoreArray(vec_dnz,&array);CHKERRQ(ierr); 307b7ce53b6SStefano Zampini ierr = VecGetArray(vec_onz,&array);CHKERRQ(ierr); 308b7ce53b6SStefano Zampini for (i=0;i<lrows;i++) onz[i] = (PetscInt)PetscRealPart(array[i]); 309b7ce53b6SStefano Zampini ierr = VecRestoreArray(vec_onz,&array);CHKERRQ(ierr); 310b7ce53b6SStefano Zampini ierr = VecDestroy(&vec_dnz);CHKERRQ(ierr); 311b7ce53b6SStefano Zampini ierr = VecDestroy(&vec_onz);CHKERRQ(ierr); 312b7ce53b6SStefano Zampini 313b7ce53b6SStefano Zampini /* Resize preallocation if overestimated */ 314b7ce53b6SStefano Zampini for (i=0;i<lrows;i++) { 315b7ce53b6SStefano Zampini dnz[i] = PetscMin(dnz[i],lcols); 316b7ce53b6SStefano Zampini onz[i] = PetscMin(onz[i],cols-lcols); 317b7ce53b6SStefano Zampini } 318b7ce53b6SStefano Zampini /* set preallocation */ 319b7ce53b6SStefano Zampini ierr = MatMPIAIJSetPreallocation(new_mat,0,dnz,0,onz);CHKERRQ(ierr); 320b7ce53b6SStefano Zampini for (i=0;i<lrows/bs;i++) { 321b7ce53b6SStefano Zampini dnz[i] = dnz[i*bs]/bs; 322b7ce53b6SStefano Zampini onz[i] = onz[i*bs]/bs; 323b7ce53b6SStefano Zampini } 324b7ce53b6SStefano Zampini ierr = MatMPIBAIJSetPreallocation(new_mat,bs,0,dnz,0,onz);CHKERRQ(ierr); 325b7ce53b6SStefano Zampini for (i=0;i<lrows/bs;i++) { 326b7ce53b6SStefano Zampini dnz[i] = dnz[i]-i; 327b7ce53b6SStefano Zampini } 328b7ce53b6SStefano Zampini ierr = MatMPISBAIJSetPreallocation(new_mat,bs,0,dnz,0,onz);CHKERRQ(ierr); 329b7ce53b6SStefano Zampini ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 330b7ce53b6SStefano Zampini *M = new_mat; 331b7ce53b6SStefano Zampini } else { 332b7ce53b6SStefano Zampini PetscInt mbs,mrows,mcols; 333b7ce53b6SStefano Zampini /* some checks */ 334b7ce53b6SStefano Zampini ierr = MatGetBlockSize(*M,&mbs);CHKERRQ(ierr); 335b7ce53b6SStefano Zampini ierr = MatGetSize(*M,&mrows,&mcols);CHKERRQ(ierr); 336b7ce53b6SStefano Zampini if (mrows != rows) { 337b7ce53b6SStefano Zampini SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of rows (%d != %d)",rows,mrows); 338b7ce53b6SStefano Zampini } 339b7ce53b6SStefano Zampini if (mrows != rows) { 340b7ce53b6SStefano Zampini SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of cols (%d != %d)",cols,mcols); 341b7ce53b6SStefano Zampini } 342b7ce53b6SStefano Zampini if (mbs != bs) { 343b7ce53b6SStefano Zampini SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong block size (%d != %d)",bs,mbs); 344b7ce53b6SStefano Zampini } 345b7ce53b6SStefano Zampini ierr = MatZeroEntries(*M);CHKERRQ(ierr); 346b7ce53b6SStefano Zampini } 347b7ce53b6SStefano Zampini /* set local to global mappings */ 348b7ce53b6SStefano Zampini /* ierr = MatSetLocalToGlobalMapping(*M,rmapping,cmapping);CHKERRQ(ierr); */ 349b7ce53b6SStefano Zampini /* Set values */ 350b7ce53b6SStefano Zampini if (isdense) { /* special case for dense local matrices */ 351b7ce53b6SStefano Zampini ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr); 352b7ce53b6SStefano Zampini ierr = MatDenseGetArray(matis->A,&array);CHKERRQ(ierr); 353b7ce53b6SStefano Zampini ierr = MatSetValues(*M,local_rows,global_indices,local_cols,global_indices,array,ADD_VALUES);CHKERRQ(ierr); 354b7ce53b6SStefano Zampini ierr = MatDenseRestoreArray(matis->A,&array);CHKERRQ(ierr); 355b7ce53b6SStefano Zampini ierr = PetscFree(local_indices);CHKERRQ(ierr); 356b7ce53b6SStefano Zampini ierr = PetscFree(global_indices);CHKERRQ(ierr); 357b7ce53b6SStefano Zampini } else { /* very basic values insertion for all other matrix types */ 358b7ce53b6SStefano Zampini ierr = PetscFree(local_indices);CHKERRQ(ierr); 359b7ce53b6SStefano Zampini for (i=0;i<local_rows;i++) { 360b7ce53b6SStefano Zampini ierr = MatGetRow(matis->A,i,&j,(const PetscInt**)&local_indices,(const PetscScalar**)&array);CHKERRQ(ierr); 361b7ce53b6SStefano Zampini /* ierr = MatSetValuesLocal(*M,1,&i,j,local_indices,array,ADD_VALUES);CHKERRQ(ierr); */ 362b7ce53b6SStefano Zampini ierr = ISLocalToGlobalMappingApply(matis->mapping,j,local_indices,global_indices);CHKERRQ(ierr); 363b7ce53b6SStefano Zampini ierr = ISLocalToGlobalMappingApply(matis->mapping,1,&i,&index_row);CHKERRQ(ierr); 364b7ce53b6SStefano Zampini ierr = MatSetValues(*M,1,&index_row,j,global_indices,array,ADD_VALUES);CHKERRQ(ierr); 365b7ce53b6SStefano Zampini ierr = MatRestoreRow(matis->A,i,&j,(const PetscInt**)&local_indices,(const PetscScalar**)&array);CHKERRQ(ierr); 366b7ce53b6SStefano Zampini } 367b7ce53b6SStefano Zampini ierr = PetscFree(global_indices);CHKERRQ(ierr); 368b7ce53b6SStefano Zampini } 369b7ce53b6SStefano Zampini ierr = MatAssemblyBegin(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 370b7ce53b6SStefano Zampini ierr = MatAssemblyEnd(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 371b7ce53b6SStefano Zampini if (isdense) { 372b7ce53b6SStefano Zampini ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr); 373b7ce53b6SStefano Zampini } 374b7ce53b6SStefano Zampini if (issbaij) { 375b7ce53b6SStefano Zampini ierr = MatRestoreRowUpperTriangular(matis->A);CHKERRQ(ierr); 376b7ce53b6SStefano Zampini } 377b7ce53b6SStefano Zampini PetscFunctionReturn(0); 378b7ce53b6SStefano Zampini } 379b7ce53b6SStefano Zampini 380b7ce53b6SStefano Zampini #undef __FUNCT__ 381b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ" 382b7ce53b6SStefano Zampini /*@ 383b7ce53b6SStefano Zampini MatISGetMPIXAIJ - Converts MATIS matrix into a parallel AIJ format 384b7ce53b6SStefano Zampini 385b7ce53b6SStefano Zampini Input Parameter: 386b7ce53b6SStefano Zampini . mat - the matrix (should be of type MATIS) 387b7ce53b6SStefano Zampini . reuse - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX 388b7ce53b6SStefano Zampini 389b7ce53b6SStefano Zampini Output Parameter: 390b7ce53b6SStefano Zampini . newmat - the matrix in AIJ format 391b7ce53b6SStefano Zampini 392b7ce53b6SStefano Zampini Level: developer 393b7ce53b6SStefano Zampini 394b7ce53b6SStefano Zampini Notes: 395b7ce53b6SStefano Zampini 396b7ce53b6SStefano Zampini .seealso: MATIS 397b7ce53b6SStefano Zampini @*/ 398b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ(Mat mat, MatReuse reuse, Mat *newmat) 399b7ce53b6SStefano Zampini { 400b7ce53b6SStefano Zampini PetscErrorCode ierr; 401b7ce53b6SStefano Zampini 402b7ce53b6SStefano Zampini PetscFunctionBegin; 403b7ce53b6SStefano Zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 404b7ce53b6SStefano Zampini PetscValidLogicalCollectiveEnum(mat,reuse,2); 405b7ce53b6SStefano Zampini PetscValidPointer(newmat,3); 406b7ce53b6SStefano Zampini if (reuse != MAT_INITIAL_MATRIX) { 407b7ce53b6SStefano Zampini PetscValidHeaderSpecific(*newmat,MAT_CLASSID,3); 408b7ce53b6SStefano Zampini PetscCheckSameComm(mat,1,*newmat,3); 409b7ce53b6SStefano Zampini } 410b7ce53b6SStefano Zampini ierr = PetscUseMethod(mat,"MatISGetMPIXAIJ_C",(Mat,MatReuse,Mat*),(mat,reuse,newmat));CHKERRQ(ierr); 411b7ce53b6SStefano Zampini PetscFunctionReturn(0); 412b7ce53b6SStefano Zampini } 413b7ce53b6SStefano Zampini 414b7ce53b6SStefano Zampini #undef __FUNCT__ 415ad6194a2SStefano Zampini #define __FUNCT__ "MatDuplicate_IS" 416ad6194a2SStefano Zampini PetscErrorCode MatDuplicate_IS(Mat mat,MatDuplicateOption op,Mat *newmat) 417ad6194a2SStefano Zampini { 418ad6194a2SStefano Zampini PetscErrorCode ierr; 419ad6194a2SStefano Zampini Mat_IS *matis = (Mat_IS*)(mat->data); 420ad6194a2SStefano Zampini PetscInt bs,m,n,M,N; 421ad6194a2SStefano Zampini Mat B,localmat; 422ad6194a2SStefano Zampini 423ad6194a2SStefano Zampini PetscFunctionBegin; 424ad6194a2SStefano Zampini ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 425ad6194a2SStefano Zampini ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 426ad6194a2SStefano Zampini ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 427ad6194a2SStefano Zampini ierr = MatCreateIS(PetscObjectComm((PetscObject)mat),bs,m,n,M,N,matis->mapping,&B);CHKERRQ(ierr); 428ad6194a2SStefano Zampini ierr = MatDuplicate(matis->A,op,&localmat);CHKERRQ(ierr); 429ad6194a2SStefano Zampini ierr = MatISSetLocalMat(B,localmat);CHKERRQ(ierr); 430b3317aa8SStefano Zampini ierr = MatDestroy(&localmat);CHKERRQ(ierr); 431ad6194a2SStefano Zampini ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 432ad6194a2SStefano Zampini ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 433ad6194a2SStefano Zampini *newmat = B; 434ad6194a2SStefano Zampini PetscFunctionReturn(0); 435ad6194a2SStefano Zampini } 436ad6194a2SStefano Zampini 437ad6194a2SStefano Zampini #undef __FUNCT__ 43869796d55SStefano Zampini #define __FUNCT__ "MatIsHermitian_IS" 43969796d55SStefano Zampini PetscErrorCode MatIsHermitian_IS(Mat A,PetscReal tol,PetscBool *flg) 44069796d55SStefano Zampini { 44169796d55SStefano Zampini PetscErrorCode ierr; 44269796d55SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 44369796d55SStefano Zampini PetscBool local_sym; 44469796d55SStefano Zampini 44569796d55SStefano Zampini PetscFunctionBegin; 44669796d55SStefano Zampini ierr = MatIsHermitian(matis->A,tol,&local_sym);CHKERRQ(ierr); 44769796d55SStefano Zampini ierr = MPI_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 44869796d55SStefano Zampini PetscFunctionReturn(0); 44969796d55SStefano Zampini } 45069796d55SStefano Zampini 45169796d55SStefano Zampini #undef __FUNCT__ 45269796d55SStefano Zampini #define __FUNCT__ "MatIsSymmetric_IS" 45369796d55SStefano Zampini PetscErrorCode MatIsSymmetric_IS(Mat A,PetscReal tol,PetscBool *flg) 45469796d55SStefano Zampini { 45569796d55SStefano Zampini PetscErrorCode ierr; 45669796d55SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 45769796d55SStefano Zampini PetscBool local_sym; 45869796d55SStefano Zampini 45969796d55SStefano Zampini PetscFunctionBegin; 46069796d55SStefano Zampini ierr = MatIsSymmetric(matis->A,tol,&local_sym);CHKERRQ(ierr); 46169796d55SStefano Zampini ierr = MPI_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 46269796d55SStefano Zampini PetscFunctionReturn(0); 46369796d55SStefano Zampini } 46469796d55SStefano Zampini 46569796d55SStefano Zampini #undef __FUNCT__ 466b4319ba4SBarry Smith #define __FUNCT__ "MatDestroy_IS" 467dfbe8321SBarry Smith PetscErrorCode MatDestroy_IS(Mat A) 468b4319ba4SBarry Smith { 469dfbe8321SBarry Smith PetscErrorCode ierr; 470b4319ba4SBarry Smith Mat_IS *b = (Mat_IS*)A->data; 471b4319ba4SBarry Smith 472b4319ba4SBarry Smith PetscFunctionBegin; 4736bf464f9SBarry Smith ierr = MatDestroy(&b->A);CHKERRQ(ierr); 4746bf464f9SBarry Smith ierr = VecScatterDestroy(&b->ctx);CHKERRQ(ierr); 4756bf464f9SBarry Smith ierr = VecDestroy(&b->x);CHKERRQ(ierr); 4766bf464f9SBarry Smith ierr = VecDestroy(&b->y);CHKERRQ(ierr); 4776bf464f9SBarry Smith ierr = ISLocalToGlobalMappingDestroy(&b->mapping);CHKERRQ(ierr); 478*28f4e0baSStefano Zampini ierr = PetscSFDestroy(&b->sf);CHKERRQ(ierr); 479*28f4e0baSStefano Zampini ierr = PetscFree2(b->sf_rootdata,b->sf_leafdata);CHKERRQ(ierr); 480bf0cc555SLisandro Dalcin ierr = PetscFree(A->data);CHKERRQ(ierr); 481dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr); 482bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",NULL);CHKERRQ(ierr); 483b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",NULL);CHKERRQ(ierr); 484b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",NULL);CHKERRQ(ierr); 4852e1947a5SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",NULL);CHKERRQ(ierr); 486b4319ba4SBarry Smith PetscFunctionReturn(0); 487b4319ba4SBarry Smith } 488b4319ba4SBarry Smith 489b4319ba4SBarry Smith #undef __FUNCT__ 490b4319ba4SBarry Smith #define __FUNCT__ "MatMult_IS" 491dfbe8321SBarry Smith PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y) 492b4319ba4SBarry Smith { 493dfbe8321SBarry Smith PetscErrorCode ierr; 494b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 495b4319ba4SBarry Smith PetscScalar zero = 0.0; 496b4319ba4SBarry Smith 497b4319ba4SBarry Smith PetscFunctionBegin; 498b4319ba4SBarry Smith /* scatter the global vector x into the local work vector */ 499ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 500ca9f406cSSatish Balay ierr = VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 501b4319ba4SBarry Smith 502b4319ba4SBarry Smith /* multiply the local matrix */ 503b4319ba4SBarry Smith ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr); 504b4319ba4SBarry Smith 505b4319ba4SBarry Smith /* scatter product back into global memory */ 5062dcb1b2aSMatthew Knepley ierr = VecSet(y,zero);CHKERRQ(ierr); 507ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 508ca9f406cSSatish Balay ierr = VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 509b4319ba4SBarry Smith PetscFunctionReturn(0); 510b4319ba4SBarry Smith } 511b4319ba4SBarry Smith 512b4319ba4SBarry Smith #undef __FUNCT__ 5132e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultAdd_IS" 514b7ce53b6SStefano Zampini PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 5152e74eeadSLisandro Dalcin { 516650997f4SStefano Zampini Vec temp_vec; 5172e74eeadSLisandro Dalcin PetscErrorCode ierr; 5182e74eeadSLisandro Dalcin 5192e74eeadSLisandro Dalcin PetscFunctionBegin; /* v3 = v2 + A * v1.*/ 520650997f4SStefano Zampini if (v3 != v2) { 521650997f4SStefano Zampini ierr = MatMult(A,v1,v3);CHKERRQ(ierr); 522650997f4SStefano Zampini ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 523650997f4SStefano Zampini } else { 524650997f4SStefano Zampini ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 525650997f4SStefano Zampini ierr = MatMult(A,v1,temp_vec);CHKERRQ(ierr); 526650997f4SStefano Zampini ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 527650997f4SStefano Zampini ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 528650997f4SStefano Zampini ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 529650997f4SStefano Zampini } 5302e74eeadSLisandro Dalcin PetscFunctionReturn(0); 5312e74eeadSLisandro Dalcin } 5322e74eeadSLisandro Dalcin 5332e74eeadSLisandro Dalcin #undef __FUNCT__ 5342e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTranspose_IS" 5352e74eeadSLisandro Dalcin PetscErrorCode MatMultTranspose_IS(Mat A,Vec x,Vec y) 5362e74eeadSLisandro Dalcin { 5372e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 5382e74eeadSLisandro Dalcin PetscErrorCode ierr; 5392e74eeadSLisandro Dalcin 5402e74eeadSLisandro Dalcin PetscFunctionBegin; /* y = A' * x */ 5412e74eeadSLisandro Dalcin /* scatter the global vector x into the local work vector */ 542ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 543ca9f406cSSatish Balay ierr = VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5442e74eeadSLisandro Dalcin 5452e74eeadSLisandro Dalcin /* multiply the local matrix */ 5462e74eeadSLisandro Dalcin ierr = MatMultTranspose(is->A,is->x,is->y);CHKERRQ(ierr); 5472e74eeadSLisandro Dalcin 5482e74eeadSLisandro Dalcin /* scatter product back into global vector */ 5492e74eeadSLisandro Dalcin ierr = VecSet(y,0);CHKERRQ(ierr); 550ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 551ca9f406cSSatish Balay ierr = VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 5522e74eeadSLisandro Dalcin PetscFunctionReturn(0); 5532e74eeadSLisandro Dalcin } 5542e74eeadSLisandro Dalcin 5552e74eeadSLisandro Dalcin #undef __FUNCT__ 5562e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTransposeAdd_IS" 5572e74eeadSLisandro Dalcin PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 5582e74eeadSLisandro Dalcin { 559650997f4SStefano Zampini Vec temp_vec; 5602e74eeadSLisandro Dalcin PetscErrorCode ierr; 5612e74eeadSLisandro Dalcin 5622e74eeadSLisandro Dalcin PetscFunctionBegin; /* v3 = v2 + A' * v1.*/ 563650997f4SStefano Zampini if (v3 != v2) { 564650997f4SStefano Zampini ierr = MatMultTranspose(A,v1,v3);CHKERRQ(ierr); 565650997f4SStefano Zampini ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 566650997f4SStefano Zampini } else { 567650997f4SStefano Zampini ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 568650997f4SStefano Zampini ierr = MatMultTranspose(A,v1,temp_vec);CHKERRQ(ierr); 569650997f4SStefano Zampini ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 570650997f4SStefano Zampini ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 571650997f4SStefano Zampini ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 572650997f4SStefano Zampini } 5732e74eeadSLisandro Dalcin PetscFunctionReturn(0); 5742e74eeadSLisandro Dalcin } 5752e74eeadSLisandro Dalcin 5762e74eeadSLisandro Dalcin #undef __FUNCT__ 577b4319ba4SBarry Smith #define __FUNCT__ "MatView_IS" 578dfbe8321SBarry Smith PetscErrorCode MatView_IS(Mat A,PetscViewer viewer) 579b4319ba4SBarry Smith { 580b4319ba4SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 581dfbe8321SBarry Smith PetscErrorCode ierr; 582b4319ba4SBarry Smith PetscViewer sviewer; 583b4319ba4SBarry Smith 584b4319ba4SBarry Smith PetscFunctionBegin; 585dd2fa690SBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr); 586b4319ba4SBarry Smith ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 587b4319ba4SBarry Smith ierr = MatView(a->A,sviewer);CHKERRQ(ierr); 588b4319ba4SBarry Smith ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 589b4319ba4SBarry Smith PetscFunctionReturn(0); 590b4319ba4SBarry Smith } 591b4319ba4SBarry Smith 592b4319ba4SBarry Smith #undef __FUNCT__ 593b4319ba4SBarry Smith #define __FUNCT__ "MatSetLocalToGlobalMapping_IS" 594784ac674SJed Brown PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping) 595b4319ba4SBarry Smith { 596dfbe8321SBarry Smith PetscErrorCode ierr; 5974e4c7dbeSStefano Zampini PetscInt n,bs; 598b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 599b4319ba4SBarry Smith IS from,to; 600b4319ba4SBarry Smith Vec global; 601b4319ba4SBarry Smith 602b4319ba4SBarry Smith PetscFunctionBegin; 603784ac674SJed Brown PetscCheckSameComm(A,1,rmapping,2); 604ce94432eSBarry Smith if (rmapping != cmapping) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"MATIS requires the row and column mappings to be identical"); 60570cf5478SStefano Zampini if (is->mapping) { /* Currenly destroys the objects that will be created by this routine. Is there anything else that should be checked? */ 60670cf5478SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&is->mapping);CHKERRQ(ierr); 60770cf5478SStefano Zampini ierr = VecDestroy(&is->x);CHKERRQ(ierr); 60870cf5478SStefano Zampini ierr = VecDestroy(&is->y);CHKERRQ(ierr); 60970cf5478SStefano Zampini ierr = VecScatterDestroy(&is->ctx);CHKERRQ(ierr); 6101c47cb0fSStefano Zampini ierr = MatDestroy(&is->A);CHKERRQ(ierr); 611*28f4e0baSStefano Zampini ierr = PetscSFDestroy(&is->sf);CHKERRQ(ierr); 612*28f4e0baSStefano Zampini ierr = PetscFree2(is->sf_rootdata,is->sf_leafdata);CHKERRQ(ierr); 61370cf5478SStefano Zampini } 614784ac674SJed Brown ierr = PetscObjectReference((PetscObject)rmapping);CHKERRQ(ierr); 6156bf464f9SBarry Smith ierr = ISLocalToGlobalMappingDestroy(&is->mapping);CHKERRQ(ierr); 616784ac674SJed Brown is->mapping = rmapping; 617fa7f1dd8SStefano Zampini /* 618fa7f1dd8SStefano Zampini ierr = PetscLayoutSetISLocalToGlobalMapping(A->rmap,rmapping);CHKERRQ(ierr); 619fa7f1dd8SStefano Zampini ierr = PetscLayoutSetISLocalToGlobalMapping(A->cmap,cmapping);CHKERRQ(ierr); 620fa7f1dd8SStefano Zampini */ 621b4319ba4SBarry Smith 622b4319ba4SBarry Smith /* Create the local matrix A */ 623784ac674SJed Brown ierr = ISLocalToGlobalMappingGetSize(rmapping,&n);CHKERRQ(ierr); 6242e1947a5SStefano Zampini ierr = ISLocalToGlobalMappingGetBlockSize(rmapping,&bs);CHKERRQ(ierr); 625f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr); 6262e1947a5SStefano Zampini if (bs > 1) { 6272e1947a5SStefano Zampini ierr = MatSetType(is->A,MATSEQBAIJ);CHKERRQ(ierr); 6282e1947a5SStefano Zampini } else { 6292e1947a5SStefano Zampini ierr = MatSetType(is->A,MATSEQAIJ);CHKERRQ(ierr); 6302e1947a5SStefano Zampini } 631f69a0ea3SMatthew Knepley ierr = MatSetSizes(is->A,n,n,n,n);CHKERRQ(ierr); 6324e4c7dbeSStefano Zampini ierr = MatSetBlockSize(is->A,bs);CHKERRQ(ierr); 633ff130e51SJed Brown ierr = MatSetOptionsPrefix(is->A,((PetscObject)A)->prefix);CHKERRQ(ierr); 634ff130e51SJed Brown ierr = MatAppendOptionsPrefix(is->A,"is_");CHKERRQ(ierr); 635b4319ba4SBarry Smith ierr = MatSetFromOptions(is->A);CHKERRQ(ierr); 636b4319ba4SBarry Smith 637b4319ba4SBarry Smith /* Create the local work vectors */ 6384e4c7dbeSStefano Zampini ierr = VecCreate(PETSC_COMM_SELF,&is->x);CHKERRQ(ierr); 6394e4c7dbeSStefano Zampini ierr = VecSetBlockSize(is->x,bs);CHKERRQ(ierr); 6404e4c7dbeSStefano Zampini ierr = VecSetSizes(is->x,n,n);CHKERRQ(ierr); 641ff130e51SJed Brown ierr = VecSetOptionsPrefix(is->x,((PetscObject)A)->prefix);CHKERRQ(ierr); 642ff130e51SJed Brown ierr = VecAppendOptionsPrefix(is->x,"is_");CHKERRQ(ierr); 6434e4c7dbeSStefano Zampini ierr = VecSetFromOptions(is->x);CHKERRQ(ierr); 644b4319ba4SBarry Smith ierr = VecDuplicate(is->x,&is->y);CHKERRQ(ierr); 645b4319ba4SBarry Smith 646b4319ba4SBarry Smith /* setup the global to local scatter */ 647b4319ba4SBarry Smith ierr = ISCreateStride(PETSC_COMM_SELF,n,0,1,&to);CHKERRQ(ierr); 648784ac674SJed Brown ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr); 6492a7a6963SBarry Smith ierr = MatCreateVecs(A,&global,NULL);CHKERRQ(ierr); 650b4319ba4SBarry Smith ierr = VecScatterCreate(global,from,is->x,to,&is->ctx);CHKERRQ(ierr); 6516bf464f9SBarry Smith ierr = VecDestroy(&global);CHKERRQ(ierr); 6526bf464f9SBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 6536bf464f9SBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 654b4319ba4SBarry Smith PetscFunctionReturn(0); 655b4319ba4SBarry Smith } 656b4319ba4SBarry Smith 6572e74eeadSLisandro Dalcin #undef __FUNCT__ 6582e74eeadSLisandro Dalcin #define __FUNCT__ "MatSetValues_IS" 6592e74eeadSLisandro Dalcin PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv) 6602e74eeadSLisandro Dalcin { 6612e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)mat->data; 6622e74eeadSLisandro Dalcin PetscInt rows_l[2048],cols_l[2048]; 6632e74eeadSLisandro Dalcin PetscErrorCode ierr; 6642e74eeadSLisandro Dalcin 6652e74eeadSLisandro Dalcin PetscFunctionBegin; 6662e74eeadSLisandro Dalcin #if defined(PETSC_USE_DEBUG) 667f23aa3ddSBarry 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); 6682e74eeadSLisandro Dalcin #endif 6692e74eeadSLisandro Dalcin ierr = ISG2LMapApply(is->mapping,m,rows,rows_l);CHKERRQ(ierr); 6702e74eeadSLisandro Dalcin ierr = ISG2LMapApply(is->mapping,n,cols,cols_l);CHKERRQ(ierr); 6712e74eeadSLisandro Dalcin ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 6722e74eeadSLisandro Dalcin PetscFunctionReturn(0); 6732e74eeadSLisandro Dalcin } 6742e74eeadSLisandro Dalcin 6752e74eeadSLisandro Dalcin #undef ISG2LMapSetUp 6762e74eeadSLisandro Dalcin #undef ISG2LMapApply 677b4319ba4SBarry Smith 678b4319ba4SBarry Smith #undef __FUNCT__ 679b4319ba4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_IS" 68013f74950SBarry Smith PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 681b4319ba4SBarry Smith { 682dfbe8321SBarry Smith PetscErrorCode ierr; 683b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 684b4319ba4SBarry Smith 685b4319ba4SBarry Smith PetscFunctionBegin; 686b4319ba4SBarry Smith ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 687b4319ba4SBarry Smith PetscFunctionReturn(0); 688b4319ba4SBarry Smith } 689b4319ba4SBarry Smith 690b4319ba4SBarry Smith #undef __FUNCT__ 691f0006bf2SLisandro Dalcin #define __FUNCT__ "MatSetValuesBlockedLocal_IS" 692f0006bf2SLisandro Dalcin PetscErrorCode MatSetValuesBlockedLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 693f0006bf2SLisandro Dalcin { 694f0006bf2SLisandro Dalcin PetscErrorCode ierr; 695f0006bf2SLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 696f0006bf2SLisandro Dalcin 697f0006bf2SLisandro Dalcin PetscFunctionBegin; 698f0006bf2SLisandro Dalcin ierr = MatSetValuesBlocked(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 699f0006bf2SLisandro Dalcin PetscFunctionReturn(0); 700f0006bf2SLisandro Dalcin } 701f0006bf2SLisandro Dalcin 702f0006bf2SLisandro Dalcin #undef __FUNCT__ 7032e74eeadSLisandro Dalcin #define __FUNCT__ "MatZeroRows_IS" 7042b40b63fSBarry Smith PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 7052e74eeadSLisandro Dalcin { 7062e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 7070298fd71SBarry Smith PetscInt n_l = 0, *rows_l = NULL; 7082e74eeadSLisandro Dalcin PetscErrorCode ierr; 7092e74eeadSLisandro Dalcin 7102e74eeadSLisandro Dalcin PetscFunctionBegin; 711ce94432eSBarry Smith if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support"); 7122e74eeadSLisandro Dalcin if (n) { 713785e854fSJed Brown ierr = PetscMalloc1(n,&rows_l);CHKERRQ(ierr); 7142e74eeadSLisandro Dalcin ierr = ISGlobalToLocalMappingApply(is->mapping,IS_GTOLM_DROP,n,rows,&n_l,rows_l);CHKERRQ(ierr); 7152e74eeadSLisandro Dalcin } 7162b40b63fSBarry Smith ierr = MatZeroRowsLocal(A,n_l,rows_l,diag,x,b);CHKERRQ(ierr); 717c31cb41cSBarry Smith ierr = PetscFree(rows_l);CHKERRQ(ierr); 7182e74eeadSLisandro Dalcin PetscFunctionReturn(0); 7192e74eeadSLisandro Dalcin } 7202e74eeadSLisandro Dalcin 7212e74eeadSLisandro Dalcin #undef __FUNCT__ 722b4319ba4SBarry Smith #define __FUNCT__ "MatZeroRowsLocal_IS" 7232b40b63fSBarry Smith PetscErrorCode MatZeroRowsLocal_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 724b4319ba4SBarry Smith { 725b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 726dfbe8321SBarry Smith PetscErrorCode ierr; 727f4df32b1SMatthew Knepley PetscInt i; 728b4319ba4SBarry Smith PetscScalar *array; 729b4319ba4SBarry Smith 730b4319ba4SBarry Smith PetscFunctionBegin; 731ce94432eSBarry Smith if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support"); 732b4319ba4SBarry Smith { 733b4319ba4SBarry Smith /* 734b4319ba4SBarry Smith Set up is->x as a "counting vector". This is in order to MatMult_IS 735b4319ba4SBarry Smith work properly in the interface nodes. 736b4319ba4SBarry Smith */ 737b4319ba4SBarry Smith Vec counter; 738b4319ba4SBarry Smith PetscScalar one=1.0, zero=0.0; 7392a7a6963SBarry Smith ierr = MatCreateVecs(A,&counter,NULL);CHKERRQ(ierr); 7402dcb1b2aSMatthew Knepley ierr = VecSet(counter,zero);CHKERRQ(ierr); 7412dcb1b2aSMatthew Knepley ierr = VecSet(is->x,one);CHKERRQ(ierr); 742ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 743ca9f406cSSatish Balay ierr = VecScatterEnd (is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 744ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 745ca9f406cSSatish Balay ierr = VecScatterEnd (is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 7466bf464f9SBarry Smith ierr = VecDestroy(&counter);CHKERRQ(ierr); 747b4319ba4SBarry Smith } 748958c9bccSBarry Smith if (!n) { 749b4319ba4SBarry Smith is->pure_neumann = PETSC_TRUE; 750b4319ba4SBarry Smith } else { 751b4319ba4SBarry Smith is->pure_neumann = PETSC_FALSE; 7522205254eSKarl Rupp 753b4319ba4SBarry Smith ierr = VecGetArray(is->x,&array);CHKERRQ(ierr); 7542b40b63fSBarry Smith ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr); 755b4319ba4SBarry Smith for (i=0; i<n; i++) { 756f4df32b1SMatthew Knepley ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr); 757b4319ba4SBarry Smith } 758b4319ba4SBarry Smith ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 759b4319ba4SBarry Smith ierr = MatAssemblyEnd (is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 760b4319ba4SBarry Smith ierr = VecRestoreArray(is->x,&array);CHKERRQ(ierr); 761b4319ba4SBarry Smith } 762b4319ba4SBarry Smith PetscFunctionReturn(0); 763b4319ba4SBarry Smith } 764b4319ba4SBarry Smith 765b4319ba4SBarry Smith #undef __FUNCT__ 766b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyBegin_IS" 767dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type) 768b4319ba4SBarry Smith { 769b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 770dfbe8321SBarry Smith PetscErrorCode ierr; 771b4319ba4SBarry Smith 772b4319ba4SBarry Smith PetscFunctionBegin; 773b4319ba4SBarry Smith ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr); 774b4319ba4SBarry Smith PetscFunctionReturn(0); 775b4319ba4SBarry Smith } 776b4319ba4SBarry Smith 777b4319ba4SBarry Smith #undef __FUNCT__ 778b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_IS" 779dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type) 780b4319ba4SBarry Smith { 781b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 782dfbe8321SBarry Smith PetscErrorCode ierr; 783b4319ba4SBarry Smith 784b4319ba4SBarry Smith PetscFunctionBegin; 785b4319ba4SBarry Smith ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr); 786b4319ba4SBarry Smith PetscFunctionReturn(0); 787b4319ba4SBarry Smith } 788b4319ba4SBarry Smith 789b4319ba4SBarry Smith #undef __FUNCT__ 790b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat_IS" 7917087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local) 792b4319ba4SBarry Smith { 793b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)mat->data; 794b4319ba4SBarry Smith 795b4319ba4SBarry Smith PetscFunctionBegin; 796b4319ba4SBarry Smith *local = is->A; 797b4319ba4SBarry Smith PetscFunctionReturn(0); 798b4319ba4SBarry Smith } 799b4319ba4SBarry Smith 800b4319ba4SBarry Smith #undef __FUNCT__ 801b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat" 802b4319ba4SBarry Smith /*@ 803b4319ba4SBarry Smith MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix. 804b4319ba4SBarry Smith 805b4319ba4SBarry Smith Input Parameter: 806b4319ba4SBarry Smith . mat - the matrix 807b4319ba4SBarry Smith 808b4319ba4SBarry Smith Output Parameter: 809b4319ba4SBarry Smith . local - the local matrix usually MATSEQAIJ 810b4319ba4SBarry Smith 811b4319ba4SBarry Smith Level: advanced 812b4319ba4SBarry Smith 813b4319ba4SBarry Smith Notes: 814b4319ba4SBarry Smith This can be called if you have precomputed the nonzero structure of the 815b4319ba4SBarry Smith matrix and want to provide it to the inner matrix object to improve the performance 816b4319ba4SBarry Smith of the MatSetValues() operation. 817b4319ba4SBarry Smith 818b4319ba4SBarry Smith .seealso: MATIS 819b4319ba4SBarry Smith @*/ 8207087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local) 821b4319ba4SBarry Smith { 8224ac538c5SBarry Smith PetscErrorCode ierr; 823b4319ba4SBarry Smith 824b4319ba4SBarry Smith PetscFunctionBegin; 8250700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 826b4319ba4SBarry Smith PetscValidPointer(local,2); 8274ac538c5SBarry Smith ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat*),(mat,local));CHKERRQ(ierr); 828b4319ba4SBarry Smith PetscFunctionReturn(0); 829b4319ba4SBarry Smith } 830b4319ba4SBarry Smith 8313b03a366Sstefano_zampini #undef __FUNCT__ 8323b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat_IS" 8333b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat_IS(Mat mat,Mat local) 8343b03a366Sstefano_zampini { 8353b03a366Sstefano_zampini Mat_IS *is = (Mat_IS*)mat->data; 8363b03a366Sstefano_zampini PetscInt nrows,ncols,orows,ocols; 8373b03a366Sstefano_zampini PetscErrorCode ierr; 8383b03a366Sstefano_zampini 8393b03a366Sstefano_zampini PetscFunctionBegin; 8404e4c7dbeSStefano Zampini if (is->A) { 8413b03a366Sstefano_zampini ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr); 8423b03a366Sstefano_zampini ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr); 843f23aa3ddSBarry 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); 8444e4c7dbeSStefano Zampini } 8453b03a366Sstefano_zampini ierr = PetscObjectReference((PetscObject)local);CHKERRQ(ierr); 8463b03a366Sstefano_zampini ierr = MatDestroy(&is->A);CHKERRQ(ierr); 8473b03a366Sstefano_zampini is->A = local; 8483b03a366Sstefano_zampini PetscFunctionReturn(0); 8493b03a366Sstefano_zampini } 8503b03a366Sstefano_zampini 8513b03a366Sstefano_zampini #undef __FUNCT__ 8523b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat" 8533b03a366Sstefano_zampini /*@ 8543b03a366Sstefano_zampini MatISSetLocalMat - Set the local matrix stored inside a MATIS matrix. 8553b03a366Sstefano_zampini 8563b03a366Sstefano_zampini Input Parameter: 8573b03a366Sstefano_zampini . mat - the matrix 8583b03a366Sstefano_zampini . local - the local matrix usually MATSEQAIJ 8593b03a366Sstefano_zampini 8603b03a366Sstefano_zampini Output Parameter: 8613b03a366Sstefano_zampini 8623b03a366Sstefano_zampini Level: advanced 8633b03a366Sstefano_zampini 8643b03a366Sstefano_zampini Notes: 8653b03a366Sstefano_zampini This can be called if you have precomputed the local matrix and 8663b03a366Sstefano_zampini want to provide it to the matrix object MATIS. 8673b03a366Sstefano_zampini 8683b03a366Sstefano_zampini .seealso: MATIS 8693b03a366Sstefano_zampini @*/ 8703b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat(Mat mat,Mat local) 8713b03a366Sstefano_zampini { 8723b03a366Sstefano_zampini PetscErrorCode ierr; 8733b03a366Sstefano_zampini 8743b03a366Sstefano_zampini PetscFunctionBegin; 8753b03a366Sstefano_zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 876b7ce53b6SStefano Zampini PetscValidHeaderSpecific(local,MAT_CLASSID,2); 8773b03a366Sstefano_zampini ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr); 8783b03a366Sstefano_zampini PetscFunctionReturn(0); 8793b03a366Sstefano_zampini } 8803b03a366Sstefano_zampini 8816726f965SBarry Smith #undef __FUNCT__ 8826726f965SBarry Smith #define __FUNCT__ "MatZeroEntries_IS" 8836726f965SBarry Smith PetscErrorCode MatZeroEntries_IS(Mat A) 8846726f965SBarry Smith { 8856726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 8866726f965SBarry Smith PetscErrorCode ierr; 8876726f965SBarry Smith 8886726f965SBarry Smith PetscFunctionBegin; 8896726f965SBarry Smith ierr = MatZeroEntries(a->A);CHKERRQ(ierr); 8906726f965SBarry Smith PetscFunctionReturn(0); 8916726f965SBarry Smith } 8926726f965SBarry Smith 8936726f965SBarry Smith #undef __FUNCT__ 8942e74eeadSLisandro Dalcin #define __FUNCT__ "MatScale_IS" 8952e74eeadSLisandro Dalcin PetscErrorCode MatScale_IS(Mat A,PetscScalar a) 8962e74eeadSLisandro Dalcin { 8972e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 8982e74eeadSLisandro Dalcin PetscErrorCode ierr; 8992e74eeadSLisandro Dalcin 9002e74eeadSLisandro Dalcin PetscFunctionBegin; 9012e74eeadSLisandro Dalcin ierr = MatScale(is->A,a);CHKERRQ(ierr); 9022e74eeadSLisandro Dalcin PetscFunctionReturn(0); 9032e74eeadSLisandro Dalcin } 9042e74eeadSLisandro Dalcin 9052e74eeadSLisandro Dalcin #undef __FUNCT__ 9062e74eeadSLisandro Dalcin #define __FUNCT__ "MatGetDiagonal_IS" 9072e74eeadSLisandro Dalcin PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v) 9082e74eeadSLisandro Dalcin { 9092e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 9102e74eeadSLisandro Dalcin PetscErrorCode ierr; 9112e74eeadSLisandro Dalcin 9122e74eeadSLisandro Dalcin PetscFunctionBegin; 9132e74eeadSLisandro Dalcin /* get diagonal of the local matrix */ 9142e74eeadSLisandro Dalcin ierr = MatGetDiagonal(is->A,is->x);CHKERRQ(ierr); 9152e74eeadSLisandro Dalcin 9162e74eeadSLisandro Dalcin /* scatter diagonal back into global vector */ 9172e74eeadSLisandro Dalcin ierr = VecSet(v,0);CHKERRQ(ierr); 918ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 919ca9f406cSSatish Balay ierr = VecScatterEnd (is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 9202e74eeadSLisandro Dalcin PetscFunctionReturn(0); 9212e74eeadSLisandro Dalcin } 9222e74eeadSLisandro Dalcin 9232e74eeadSLisandro Dalcin #undef __FUNCT__ 9246726f965SBarry Smith #define __FUNCT__ "MatSetOption_IS" 925ace3abfcSBarry Smith PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg) 9266726f965SBarry Smith { 9276726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 9286726f965SBarry Smith PetscErrorCode ierr; 9296726f965SBarry Smith 9306726f965SBarry Smith PetscFunctionBegin; 9314e0d8c25SBarry Smith ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr); 9326726f965SBarry Smith PetscFunctionReturn(0); 9336726f965SBarry Smith } 9346726f965SBarry Smith 935284134d9SBarry Smith #undef __FUNCT__ 936284134d9SBarry Smith #define __FUNCT__ "MatCreateIS" 937284134d9SBarry Smith /*@ 938284134d9SBarry Smith MatCreateIS - Creates a "process" unassmembled matrix, it is assembled on each 939284134d9SBarry Smith process but not across processes. 940284134d9SBarry Smith 941284134d9SBarry Smith Input Parameters: 942284134d9SBarry Smith + comm - MPI communicator that will share the matrix 9434e4c7dbeSStefano Zampini . bs - local and global block size of the matrix 944df3898eeSBarry Smith . m,n,M,N - local and/or global sizes of the left and right vector used in matrix vector products 945284134d9SBarry Smith - map - mapping that defines the global number for each local number 946284134d9SBarry Smith 947284134d9SBarry Smith Output Parameter: 948284134d9SBarry Smith . A - the resulting matrix 949284134d9SBarry Smith 9508e6c10adSSatish Balay Level: advanced 9518e6c10adSSatish Balay 952284134d9SBarry Smith Notes: See MATIS for more details 9538cda6cd7SBarry Smith m and n are NOT related to the size of the map, they are the size of the part of the vector owned 9548cda6cd7SBarry Smith by that process. m + nghosts (or n + nghosts) is the length of map since map maps all local points 9558cda6cd7SBarry Smith plus the ghost points to global indices. 956284134d9SBarry Smith 957284134d9SBarry Smith .seealso: MATIS, MatSetLocalToGlobalMapping() 958284134d9SBarry Smith @*/ 9594e4c7dbeSStefano Zampini PetscErrorCode MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping map,Mat *A) 960284134d9SBarry Smith { 961284134d9SBarry Smith PetscErrorCode ierr; 962284134d9SBarry Smith 963284134d9SBarry Smith PetscFunctionBegin; 964284134d9SBarry Smith ierr = MatCreate(comm,A);CHKERRQ(ierr); 965d3ee2243SStefano Zampini ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr); 966284134d9SBarry Smith ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 967284134d9SBarry Smith ierr = MatSetType(*A,MATIS);CHKERRQ(ierr); 9683b03a366Sstefano_zampini ierr = MatSetUp(*A);CHKERRQ(ierr); 969784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(*A,map,map);CHKERRQ(ierr); 970284134d9SBarry Smith PetscFunctionReturn(0); 971284134d9SBarry Smith } 972284134d9SBarry Smith 973b4319ba4SBarry Smith /*MC 9746a818285SBarry Smith MATIS - MATIS = "is" - A matrix type to be used for using the Neumann-Neumann type preconditioners, see PCBDDC. 975b4319ba4SBarry Smith This stores the matrices in globally unassembled form. Each processor 976b4319ba4SBarry Smith assembles only its local Neumann problem and the parallel matrix vector 977b4319ba4SBarry Smith product is handled "implicitly". 978b4319ba4SBarry Smith 979b4319ba4SBarry Smith Operations Provided: 9806726f965SBarry Smith + MatMult() 9812e74eeadSLisandro Dalcin . MatMultAdd() 9822e74eeadSLisandro Dalcin . MatMultTranspose() 9832e74eeadSLisandro Dalcin . MatMultTransposeAdd() 9846726f965SBarry Smith . MatZeroEntries() 9856726f965SBarry Smith . MatSetOption() 9862e74eeadSLisandro Dalcin . MatZeroRows() 9876726f965SBarry Smith . MatZeroRowsLocal() 9882e74eeadSLisandro Dalcin . MatSetValues() 9896726f965SBarry Smith . MatSetValuesLocal() 9902e74eeadSLisandro Dalcin . MatScale() 9912e74eeadSLisandro Dalcin . MatGetDiagonal() 9926726f965SBarry Smith - MatSetLocalToGlobalMapping() 993b4319ba4SBarry Smith 994b4319ba4SBarry Smith Options Database Keys: 995b4319ba4SBarry Smith . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions() 996b4319ba4SBarry Smith 997b4319ba4SBarry Smith Notes: Options prefix for the inner matrix are given by -is_mat_xxx 998b4319ba4SBarry Smith 999b4319ba4SBarry Smith You must call MatSetLocalToGlobalMapping() before using this matrix type. 1000b4319ba4SBarry Smith 1001b4319ba4SBarry Smith You can do matrix preallocation on the local matrix after you obtain it with 1002b4319ba4SBarry Smith MatISGetLocalMat() 1003b4319ba4SBarry Smith 1004b4319ba4SBarry Smith Level: advanced 1005b4319ba4SBarry Smith 10066a818285SBarry Smith .seealso: PC, MatISGetLocalMat(), MatSetLocalToGlobalMapping(), PCBDDC 1007b4319ba4SBarry Smith 1008b4319ba4SBarry Smith M*/ 1009b4319ba4SBarry Smith 1010b4319ba4SBarry Smith #undef __FUNCT__ 1011b4319ba4SBarry Smith #define __FUNCT__ "MatCreate_IS" 10128cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_IS(Mat A) 1013b4319ba4SBarry Smith { 1014dfbe8321SBarry Smith PetscErrorCode ierr; 1015b4319ba4SBarry Smith Mat_IS *b; 1016b4319ba4SBarry Smith 1017b4319ba4SBarry Smith PetscFunctionBegin; 1018b00a9115SJed Brown ierr = PetscNewLog(A,&b);CHKERRQ(ierr); 1019b4319ba4SBarry Smith A->data = (void*)b; 1020b4319ba4SBarry Smith ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1021b4319ba4SBarry Smith 1022b4319ba4SBarry Smith A->ops->mult = MatMult_IS; 10232e74eeadSLisandro Dalcin A->ops->multadd = MatMultAdd_IS; 10242e74eeadSLisandro Dalcin A->ops->multtranspose = MatMultTranspose_IS; 10252e74eeadSLisandro Dalcin A->ops->multtransposeadd = MatMultTransposeAdd_IS; 1026b4319ba4SBarry Smith A->ops->destroy = MatDestroy_IS; 1027b4319ba4SBarry Smith A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS; 10282e74eeadSLisandro Dalcin A->ops->setvalues = MatSetValues_IS; 1029b4319ba4SBarry Smith A->ops->setvalueslocal = MatSetValuesLocal_IS; 1030f0006bf2SLisandro Dalcin A->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_IS; 10312e74eeadSLisandro Dalcin A->ops->zerorows = MatZeroRows_IS; 1032b4319ba4SBarry Smith A->ops->zerorowslocal = MatZeroRowsLocal_IS; 1033b4319ba4SBarry Smith A->ops->assemblybegin = MatAssemblyBegin_IS; 1034b4319ba4SBarry Smith A->ops->assemblyend = MatAssemblyEnd_IS; 1035b4319ba4SBarry Smith A->ops->view = MatView_IS; 10366726f965SBarry Smith A->ops->zeroentries = MatZeroEntries_IS; 10372e74eeadSLisandro Dalcin A->ops->scale = MatScale_IS; 10382e74eeadSLisandro Dalcin A->ops->getdiagonal = MatGetDiagonal_IS; 10396726f965SBarry Smith A->ops->setoption = MatSetOption_IS; 104069796d55SStefano Zampini A->ops->ishermitian = MatIsHermitian_IS; 104169796d55SStefano Zampini A->ops->issymmetric = MatIsSymmetric_IS; 1042ad6194a2SStefano Zampini A->ops->duplicate = MatDuplicate_IS; 1043b4319ba4SBarry Smith 104426283091SBarry Smith ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 104526283091SBarry Smith ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 1046b4319ba4SBarry Smith 1047b7ce53b6SStefano Zampini /* zeros pointer */ 1048b4319ba4SBarry Smith b->A = 0; 1049b4319ba4SBarry Smith b->ctx = 0; 1050b4319ba4SBarry Smith b->x = 0; 1051b4319ba4SBarry Smith b->y = 0; 1052b09366fdSStefano Zampini b->mapping = 0; 1053*28f4e0baSStefano Zampini b->sf = 0; 1054*28f4e0baSStefano Zampini b->sf_rootdata = 0; 1055*28f4e0baSStefano Zampini b->sf_leafdata = 0; 1056b7ce53b6SStefano Zampini 1057b7ce53b6SStefano Zampini /* special MATIS functions */ 1058bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",MatISGetLocalMat_IS);CHKERRQ(ierr); 1059bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",MatISSetLocalMat_IS);CHKERRQ(ierr); 1060b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",MatISGetMPIXAIJ_IS);CHKERRQ(ierr); 10612e1947a5SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",MatISSetPreallocation_IS);CHKERRQ(ierr); 106217667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr); 1063b4319ba4SBarry Smith PetscFunctionReturn(0); 1064b4319ba4SBarry Smith } 1065