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*/ 162e1947a5SStefano Zampini #include <petscsf.h> 172e1947a5SStefano Zampini 182e1947a5SStefano Zampini #undef __FUNCT__ 192e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation" 20*a88811baSStefano Zampini /* 21*a88811baSStefano Zampini MatISSetPreallocation - Preallocates memory for a MATIS parallel matrix. 22*a88811baSStefano Zampini 23*a88811baSStefano Zampini Collective on MPI_Comm 24*a88811baSStefano Zampini 25*a88811baSStefano Zampini Input Parameters: 26*a88811baSStefano Zampini + B - the matrix 27*a88811baSStefano Zampini . d_nz - number of nonzeros per row in DIAGONAL portion of local submatrix 28*a88811baSStefano Zampini (same value is used for all local rows) 29*a88811baSStefano Zampini . d_nnz - array containing the number of nonzeros in the various rows of the 30*a88811baSStefano Zampini DIAGONAL portion of the local submatrix (possibly different for each row) 31*a88811baSStefano Zampini or NULL, if d_nz is used to specify the nonzero structure. 32*a88811baSStefano Zampini The size of this array is equal to the number of local rows, i.e 'm'. 33*a88811baSStefano Zampini For matrices that will be factored, you must leave room for (and set) 34*a88811baSStefano Zampini the diagonal entry even if it is zero. 35*a88811baSStefano Zampini . o_nz - number of nonzeros per row in the OFF-DIAGONAL portion of local 36*a88811baSStefano Zampini submatrix (same value is used for all local rows). 37*a88811baSStefano Zampini - o_nnz - array containing the number of nonzeros in the various rows of the 38*a88811baSStefano Zampini OFF-DIAGONAL portion of the local submatrix (possibly different for 39*a88811baSStefano Zampini each row) or NULL, if o_nz is used to specify the nonzero 40*a88811baSStefano Zampini structure. The size of this array is equal to the number 41*a88811baSStefano Zampini of local rows, i.e 'm'. 42*a88811baSStefano Zampini 43*a88811baSStefano Zampini If the *_nnz parameter is given then the *_nz parameter is ignored 44*a88811baSStefano Zampini 45*a88811baSStefano Zampini Level: intermediate 46*a88811baSStefano Zampini 47*a88811baSStefano Zampini Notes: This function has the same interface as the MPIAIJ preallocation routine in order to simplify the transition 48*a88811baSStefano Zampini from the asssembled format to the unassembled one. It overestimates the preallocation of MATIS local 49*a88811baSStefano Zampini matrices; for exact preallocation, the user should set the preallocation directly on local matrix objects. 50*a88811baSStefano Zampini 51*a88811baSStefano Zampini .keywords: matrix 52*a88811baSStefano Zampini 53*a88811baSStefano Zampini .seealso: MatCreate(), MatCreateIS(), MatMPIAIJSetPreallocation(), MatISGetLocalMat() 54*a88811baSStefano Zampini @*/ 552e1947a5SStefano Zampini PetscErrorCode MatISSetPreallocation(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[]) 562e1947a5SStefano Zampini { 572e1947a5SStefano Zampini PetscErrorCode ierr; 582e1947a5SStefano Zampini 592e1947a5SStefano Zampini PetscFunctionBegin; 602e1947a5SStefano Zampini PetscValidHeaderSpecific(B,MAT_CLASSID,1); 612e1947a5SStefano Zampini PetscValidType(B,1); 622e1947a5SStefano Zampini ierr = PetscTryMethod(B,"MatISSetPreallocation_C",(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(B,d_nz,d_nnz,o_nz,o_nnz));CHKERRQ(ierr); 632e1947a5SStefano Zampini PetscFunctionReturn(0); 642e1947a5SStefano Zampini } 652e1947a5SStefano Zampini 662e1947a5SStefano Zampini #undef __FUNCT__ 672e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation_IS" 682e1947a5SStefano Zampini PetscErrorCode MatISSetPreallocation_IS(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[]) 692e1947a5SStefano Zampini { 702e1947a5SStefano Zampini Mat_IS *matis = (Mat_IS*)(B->data); 712e1947a5SStefano Zampini PetscSF sf; 722e1947a5SStefano Zampini PetscInt bs,i,nroots,*rootdata,nleaves,*leafdata,nlocalcols; 732e1947a5SStefano Zampini const PetscInt *gidxs; 742e1947a5SStefano Zampini PetscErrorCode ierr; 752e1947a5SStefano Zampini 762e1947a5SStefano Zampini PetscFunctionBegin; 772e1947a5SStefano Zampini if (!matis->A) { 782e1947a5SStefano Zampini SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"You should first call MatSetLocalToGlobalMapping"); 792e1947a5SStefano Zampini } 802e1947a5SStefano Zampini ierr = MatGetLocalSize(B,&nroots,NULL);CHKERRQ(ierr); 812e1947a5SStefano Zampini ierr = MatGetSize(matis->A,&nleaves,&nlocalcols);CHKERRQ(ierr); 822e1947a5SStefano Zampini ierr = MatGetBlockSize(matis->A,&bs);CHKERRQ(ierr); 832e1947a5SStefano Zampini ierr = PetscCalloc2(nroots,&rootdata,nleaves,&leafdata);CHKERRQ(ierr); 842e1947a5SStefano Zampini ierr = PetscSFCreate(PetscObjectComm((PetscObject)B),&sf);CHKERRQ(ierr); 852e1947a5SStefano Zampini ierr = PetscSFSetFromOptions(sf);CHKERRQ(ierr); 862e1947a5SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(matis->mapping,&gidxs);CHKERRQ(ierr); 872e1947a5SStefano Zampini ierr = PetscSFSetGraphLayout(sf,B->rmap,nleaves,NULL,PETSC_COPY_VALUES,gidxs);CHKERRQ(ierr); 882e1947a5SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(matis->mapping,&gidxs);CHKERRQ(ierr); 892e1947a5SStefano Zampini if (!d_nnz) { 902e1947a5SStefano Zampini for (i=0;i<nroots;i++) rootdata[i] += d_nz; 912e1947a5SStefano Zampini } else { 922e1947a5SStefano Zampini for (i=0;i<nroots;i++) rootdata[i] += d_nnz[i]; 932e1947a5SStefano Zampini } 942e1947a5SStefano Zampini if (!o_nnz) { 952e1947a5SStefano Zampini for (i=0;i<nroots;i++) rootdata[i] += o_nz; 962e1947a5SStefano Zampini } else { 972e1947a5SStefano Zampini for (i=0;i<nroots;i++) rootdata[i] += o_nnz[i]; 982e1947a5SStefano Zampini } 992e1947a5SStefano Zampini ierr = PetscSFBcastBegin(sf,MPIU_INT,rootdata,leafdata);CHKERRQ(ierr); 1002e1947a5SStefano Zampini ierr = PetscSFBcastEnd(sf,MPIU_INT,rootdata,leafdata);CHKERRQ(ierr); 1012e1947a5SStefano Zampini for (i=0;i<nleaves;i++) { 1022e1947a5SStefano Zampini leafdata[i] = PetscMin(leafdata[i],nlocalcols); 1032e1947a5SStefano Zampini } 1042e1947a5SStefano Zampini ierr = MatSeqAIJSetPreallocation(matis->A,0,leafdata);CHKERRQ(ierr); 1052e1947a5SStefano Zampini for (i=0;i<nleaves/bs;i++) { 1062e1947a5SStefano Zampini leafdata[i] = leafdata[i*bs]/bs; 1072e1947a5SStefano Zampini } 1082e1947a5SStefano Zampini ierr = MatSeqBAIJSetPreallocation(matis->A,bs,0,leafdata);CHKERRQ(ierr); 1092e1947a5SStefano Zampini for (i=0;i<nleaves/bs;i++) { 1102e1947a5SStefano Zampini leafdata[i] = leafdata[i]-i; 1112e1947a5SStefano Zampini } 1122e1947a5SStefano Zampini ierr = MatSeqSBAIJSetPreallocation(matis->A,bs,0,leafdata);CHKERRQ(ierr); 1132e1947a5SStefano Zampini ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 1142e1947a5SStefano Zampini ierr = PetscFree2(rootdata,leafdata);CHKERRQ(ierr); 1152e1947a5SStefano Zampini PetscFunctionReturn(0); 1162e1947a5SStefano Zampini } 117b4319ba4SBarry Smith 118b4319ba4SBarry Smith #undef __FUNCT__ 119b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ_IS" 120b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ_IS(Mat mat, MatReuse reuse, Mat *M) 121b7ce53b6SStefano Zampini { 122b7ce53b6SStefano Zampini Mat_IS *matis = (Mat_IS*)(mat->data); 123b7ce53b6SStefano Zampini /* info on mat */ 124b7ce53b6SStefano Zampini /* ISLocalToGlobalMapping rmapping,cmapping; */ 125b7ce53b6SStefano Zampini PetscInt bs,rows,cols; 126b7ce53b6SStefano Zampini PetscInt lrows,lcols; 127b7ce53b6SStefano Zampini PetscInt local_rows,local_cols; 128b7ce53b6SStefano Zampini PetscBool isdense,issbaij,issbaij_red; 129b7ce53b6SStefano Zampini /* values insertion */ 130b7ce53b6SStefano Zampini PetscScalar *array; 131b7ce53b6SStefano Zampini PetscInt *local_indices,*global_indices; 132b7ce53b6SStefano Zampini /* work */ 133b7ce53b6SStefano Zampini PetscInt i,j,index_row; 134b7ce53b6SStefano Zampini PetscErrorCode ierr; 135b7ce53b6SStefano Zampini 136b7ce53b6SStefano Zampini PetscFunctionBegin; 137b7ce53b6SStefano Zampini /* MISSING CHECKS 138b7ce53b6SStefano Zampini - rectangular case not covered (it is not allowed by MATIS) 139b7ce53b6SStefano Zampini */ 140b7ce53b6SStefano Zampini /* get info from mat */ 141b7ce53b6SStefano Zampini /* ierr = MatGetLocalToGlobalMapping(mat,&rmapping,&cmapping);CHKERRQ(ierr); */ 142b7ce53b6SStefano Zampini ierr = MatGetSize(mat,&rows,&cols);CHKERRQ(ierr); 143b7ce53b6SStefano Zampini ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 144b7ce53b6SStefano Zampini ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr); 145b7ce53b6SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr); 146b7ce53b6SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr); 147b7ce53b6SStefano Zampini 148b7ce53b6SStefano Zampini /* work */ 149b7ce53b6SStefano Zampini ierr = PetscMalloc1(local_rows,&local_indices);CHKERRQ(ierr); 150b7ce53b6SStefano Zampini for (i=0;i<local_rows;i++) local_indices[i]=i; 151b7ce53b6SStefano Zampini /* map indices of local mat to global values */ 152b7ce53b6SStefano Zampini ierr = PetscMalloc(PetscMax(local_cols,local_rows)*sizeof(*global_indices),&global_indices);CHKERRQ(ierr); 153b7ce53b6SStefano Zampini /* ierr = ISLocalToGlobalMappingApply(rmapping,local_rows,local_indices,global_indices);CHKERRQ(ierr); */ 154b7ce53b6SStefano Zampini ierr = ISLocalToGlobalMappingApply(matis->mapping,local_rows,local_indices,global_indices);CHKERRQ(ierr); 155b7ce53b6SStefano Zampini 156b7ce53b6SStefano Zampini if (issbaij) { 157b7ce53b6SStefano Zampini ierr = MatGetRowUpperTriangular(matis->A);CHKERRQ(ierr); 158b7ce53b6SStefano Zampini } 159b7ce53b6SStefano Zampini 160b7ce53b6SStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { 161b7ce53b6SStefano Zampini Mat new_mat; 162b7ce53b6SStefano Zampini MatType new_mat_type; 163b7ce53b6SStefano Zampini Vec vec_dnz,vec_onz; 164b7ce53b6SStefano Zampini PetscScalar *my_dnz,*my_onz; 165b7ce53b6SStefano Zampini PetscInt *dnz,*onz,*mat_ranges,*row_ownership; 166b7ce53b6SStefano Zampini PetscInt index_col,owner; 167b7ce53b6SStefano Zampini PetscMPIInt nsubdomains; 168b7ce53b6SStefano Zampini 169b7ce53b6SStefano Zampini /* determining new matrix type */ 170b7ce53b6SStefano Zampini ierr = MPI_Allreduce(&issbaij,&issbaij_red,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 171b7ce53b6SStefano Zampini if (issbaij_red) { 172b7ce53b6SStefano Zampini new_mat_type = MATSBAIJ; 173b7ce53b6SStefano Zampini } else { 174b7ce53b6SStefano Zampini if (bs>1) { 175b7ce53b6SStefano Zampini new_mat_type = MATBAIJ; 176b7ce53b6SStefano Zampini } else { 177b7ce53b6SStefano Zampini new_mat_type = MATAIJ; 178b7ce53b6SStefano Zampini } 179b7ce53b6SStefano Zampini } 180b7ce53b6SStefano Zampini 181b7ce53b6SStefano Zampini ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&nsubdomains);CHKERRQ(ierr); 182b7ce53b6SStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)mat),&new_mat);CHKERRQ(ierr); 183b7ce53b6SStefano Zampini ierr = MatSetSizes(new_mat,PETSC_DECIDE,PETSC_DECIDE,rows,cols);CHKERRQ(ierr); 184b7ce53b6SStefano Zampini ierr = MatSetBlockSize(new_mat,bs);CHKERRQ(ierr); 185b7ce53b6SStefano Zampini ierr = MatSetType(new_mat,new_mat_type);CHKERRQ(ierr); 186b7ce53b6SStefano Zampini ierr = MatSetUp(new_mat);CHKERRQ(ierr); 187b7ce53b6SStefano Zampini ierr = MatGetLocalSize(new_mat,&lrows,&lcols);CHKERRQ(ierr); 188b7ce53b6SStefano Zampini 189b7ce53b6SStefano Zampini /* 190b7ce53b6SStefano Zampini preallocation 191b7ce53b6SStefano Zampini */ 192b7ce53b6SStefano Zampini 193b7ce53b6SStefano Zampini ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)new_mat),lrows,lcols,dnz,onz);CHKERRQ(ierr); 194b7ce53b6SStefano Zampini /* 195b7ce53b6SStefano Zampini Some vectors are needed to sum up properly on shared interface dofs. 196b7ce53b6SStefano Zampini Preallocation macros cannot do the job. 197b7ce53b6SStefano Zampini Note that preallocation is not exact, since it overestimates nonzeros 198b7ce53b6SStefano Zampini */ 199b7ce53b6SStefano Zampini ierr = MatGetVecs(new_mat,NULL,&vec_dnz);CHKERRQ(ierr); 200b7ce53b6SStefano Zampini /* ierr = VecSetLocalToGlobalMapping(vec_dnz,rmapping);CHKERRQ(ierr); */ 201b7ce53b6SStefano Zampini ierr = VecSetLocalToGlobalMapping(vec_dnz,matis->mapping);CHKERRQ(ierr); 202b7ce53b6SStefano Zampini ierr = VecDuplicate(vec_dnz,&vec_onz);CHKERRQ(ierr); 203b7ce53b6SStefano Zampini /* All processes need to compute entire row ownership */ 204b7ce53b6SStefano Zampini ierr = PetscMalloc1(rows,&row_ownership);CHKERRQ(ierr); 205b7ce53b6SStefano Zampini ierr = MatGetOwnershipRanges(new_mat,(const PetscInt**)&mat_ranges);CHKERRQ(ierr); 206b7ce53b6SStefano Zampini for (i=0;i<nsubdomains;i++) { 207b7ce53b6SStefano Zampini for (j=mat_ranges[i];j<mat_ranges[i+1];j++) { 208b7ce53b6SStefano Zampini row_ownership[j]=i; 209b7ce53b6SStefano Zampini } 210b7ce53b6SStefano Zampini } 211b7ce53b6SStefano Zampini 212b7ce53b6SStefano Zampini /* 213b7ce53b6SStefano Zampini my_dnz and my_onz contains exact contribution to preallocation from each local mat 214b7ce53b6SStefano Zampini then, they will be summed up properly. This way, preallocation is always sufficient 215b7ce53b6SStefano Zampini */ 216b7ce53b6SStefano Zampini ierr = PetscMalloc1(local_rows,&my_dnz);CHKERRQ(ierr); 217b7ce53b6SStefano Zampini ierr = PetscMalloc1(local_rows,&my_onz);CHKERRQ(ierr); 218b7ce53b6SStefano Zampini ierr = PetscMemzero(my_dnz,local_rows*sizeof(*my_dnz));CHKERRQ(ierr); 219b7ce53b6SStefano Zampini ierr = PetscMemzero(my_onz,local_rows*sizeof(*my_onz));CHKERRQ(ierr); 220b7ce53b6SStefano Zampini /* preallocation as a MATAIJ */ 221b7ce53b6SStefano Zampini if (isdense) { /* special case for dense local matrices */ 222b7ce53b6SStefano Zampini for (i=0;i<local_rows;i++) { 223b7ce53b6SStefano Zampini index_row = global_indices[i]; 224b7ce53b6SStefano Zampini for (j=i;j<local_rows;j++) { 225b7ce53b6SStefano Zampini owner = row_ownership[index_row]; 226b7ce53b6SStefano Zampini index_col = global_indices[j]; 227b7ce53b6SStefano Zampini if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */ 228b7ce53b6SStefano Zampini my_dnz[i] += 1.0; 229b7ce53b6SStefano Zampini } else { /* offdiag block */ 230b7ce53b6SStefano Zampini my_onz[i] += 1.0; 231b7ce53b6SStefano Zampini } 232b7ce53b6SStefano Zampini /* same as before, interchanging rows and cols */ 233b7ce53b6SStefano Zampini if (i != j) { 234b7ce53b6SStefano Zampini owner = row_ownership[index_col]; 235b7ce53b6SStefano Zampini if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) { 236b7ce53b6SStefano Zampini my_dnz[j] += 1.0; 237b7ce53b6SStefano Zampini } else { 238b7ce53b6SStefano Zampini my_onz[j] += 1.0; 239b7ce53b6SStefano Zampini } 240b7ce53b6SStefano Zampini } 241b7ce53b6SStefano Zampini } 242b7ce53b6SStefano Zampini } 243b7ce53b6SStefano Zampini } else { 244b7ce53b6SStefano Zampini for (i=0;i<local_rows;i++) { 245b7ce53b6SStefano Zampini PetscInt ncols; 246b7ce53b6SStefano Zampini const PetscInt *cols; 247b7ce53b6SStefano Zampini index_row = global_indices[i]; 248b7ce53b6SStefano Zampini ierr = MatGetRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr); 249b7ce53b6SStefano Zampini for (j=0;j<ncols;j++) { 250b7ce53b6SStefano Zampini owner = row_ownership[index_row]; 251b7ce53b6SStefano Zampini index_col = global_indices[cols[j]]; 252b7ce53b6SStefano Zampini if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */ 253b7ce53b6SStefano Zampini my_dnz[i] += 1.0; 254b7ce53b6SStefano Zampini } else { /* offdiag block */ 255b7ce53b6SStefano Zampini my_onz[i] += 1.0; 256b7ce53b6SStefano Zampini } 257b7ce53b6SStefano Zampini /* same as before, interchanging rows and cols */ 258b7ce53b6SStefano Zampini if (issbaij) { 259b7ce53b6SStefano Zampini owner = row_ownership[index_col]; 260b7ce53b6SStefano Zampini if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) { 261b7ce53b6SStefano Zampini my_dnz[j] += 1.0; 262b7ce53b6SStefano Zampini } else { 263b7ce53b6SStefano Zampini my_onz[j] += 1.0; 264b7ce53b6SStefano Zampini } 265b7ce53b6SStefano Zampini } 266b7ce53b6SStefano Zampini } 267b7ce53b6SStefano Zampini ierr = MatRestoreRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr); 268b7ce53b6SStefano Zampini } 269b7ce53b6SStefano Zampini } 270b7ce53b6SStefano Zampini ierr = VecSet(vec_dnz,0.0);CHKERRQ(ierr); 271b7ce53b6SStefano Zampini ierr = VecSet(vec_onz,0.0);CHKERRQ(ierr); 272b7ce53b6SStefano Zampini if (local_rows) { /* multilevel guard */ 273b7ce53b6SStefano Zampini ierr = VecSetValuesLocal(vec_dnz,local_rows,local_indices,my_dnz,ADD_VALUES);CHKERRQ(ierr); 274b7ce53b6SStefano Zampini ierr = VecSetValuesLocal(vec_onz,local_rows,local_indices,my_onz,ADD_VALUES);CHKERRQ(ierr); 275b7ce53b6SStefano Zampini } 276b7ce53b6SStefano Zampini ierr = VecAssemblyBegin(vec_dnz);CHKERRQ(ierr); 277b7ce53b6SStefano Zampini ierr = VecAssemblyBegin(vec_onz);CHKERRQ(ierr); 278b7ce53b6SStefano Zampini ierr = VecAssemblyEnd(vec_dnz);CHKERRQ(ierr); 279b7ce53b6SStefano Zampini ierr = VecAssemblyEnd(vec_onz);CHKERRQ(ierr); 280b7ce53b6SStefano Zampini ierr = PetscFree(my_dnz);CHKERRQ(ierr); 281b7ce53b6SStefano Zampini ierr = PetscFree(my_onz);CHKERRQ(ierr); 282b7ce53b6SStefano Zampini ierr = PetscFree(row_ownership);CHKERRQ(ierr); 283b7ce53b6SStefano Zampini 284b7ce53b6SStefano Zampini /* set computed preallocation in dnz and onz */ 285b7ce53b6SStefano Zampini ierr = VecGetArray(vec_dnz,&array);CHKERRQ(ierr); 286b7ce53b6SStefano Zampini for (i=0; i<lrows; i++) dnz[i] = (PetscInt)PetscRealPart(array[i]); 287b7ce53b6SStefano Zampini ierr = VecRestoreArray(vec_dnz,&array);CHKERRQ(ierr); 288b7ce53b6SStefano Zampini ierr = VecGetArray(vec_onz,&array);CHKERRQ(ierr); 289b7ce53b6SStefano Zampini for (i=0;i<lrows;i++) onz[i] = (PetscInt)PetscRealPart(array[i]); 290b7ce53b6SStefano Zampini ierr = VecRestoreArray(vec_onz,&array);CHKERRQ(ierr); 291b7ce53b6SStefano Zampini ierr = VecDestroy(&vec_dnz);CHKERRQ(ierr); 292b7ce53b6SStefano Zampini ierr = VecDestroy(&vec_onz);CHKERRQ(ierr); 293b7ce53b6SStefano Zampini 294b7ce53b6SStefano Zampini /* Resize preallocation if overestimated */ 295b7ce53b6SStefano Zampini for (i=0;i<lrows;i++) { 296b7ce53b6SStefano Zampini dnz[i] = PetscMin(dnz[i],lcols); 297b7ce53b6SStefano Zampini onz[i] = PetscMin(onz[i],cols-lcols); 298b7ce53b6SStefano Zampini } 299b7ce53b6SStefano Zampini /* set preallocation */ 300b7ce53b6SStefano Zampini ierr = MatMPIAIJSetPreallocation(new_mat,0,dnz,0,onz);CHKERRQ(ierr); 301b7ce53b6SStefano Zampini for (i=0;i<lrows/bs;i++) { 302b7ce53b6SStefano Zampini dnz[i] = dnz[i*bs]/bs; 303b7ce53b6SStefano Zampini onz[i] = onz[i*bs]/bs; 304b7ce53b6SStefano Zampini } 305b7ce53b6SStefano Zampini ierr = MatMPIBAIJSetPreallocation(new_mat,bs,0,dnz,0,onz);CHKERRQ(ierr); 306b7ce53b6SStefano Zampini for (i=0;i<lrows/bs;i++) { 307b7ce53b6SStefano Zampini dnz[i] = dnz[i]-i; 308b7ce53b6SStefano Zampini } 309b7ce53b6SStefano Zampini ierr = MatMPISBAIJSetPreallocation(new_mat,bs,0,dnz,0,onz);CHKERRQ(ierr); 310b7ce53b6SStefano Zampini ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 311b7ce53b6SStefano Zampini *M = new_mat; 312b7ce53b6SStefano Zampini } else { 313b7ce53b6SStefano Zampini PetscInt mbs,mrows,mcols; 314b7ce53b6SStefano Zampini /* some checks */ 315b7ce53b6SStefano Zampini ierr = MatGetBlockSize(*M,&mbs);CHKERRQ(ierr); 316b7ce53b6SStefano Zampini ierr = MatGetSize(*M,&mrows,&mcols);CHKERRQ(ierr); 317b7ce53b6SStefano Zampini if (mrows != rows) { 318b7ce53b6SStefano Zampini SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of rows (%d != %d)",rows,mrows); 319b7ce53b6SStefano Zampini } 320b7ce53b6SStefano Zampini if (mrows != rows) { 321b7ce53b6SStefano Zampini SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of cols (%d != %d)",cols,mcols); 322b7ce53b6SStefano Zampini } 323b7ce53b6SStefano Zampini if (mbs != bs) { 324b7ce53b6SStefano Zampini SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong block size (%d != %d)",bs,mbs); 325b7ce53b6SStefano Zampini } 326b7ce53b6SStefano Zampini ierr = MatZeroEntries(*M);CHKERRQ(ierr); 327b7ce53b6SStefano Zampini } 328b7ce53b6SStefano Zampini /* set local to global mappings */ 329b7ce53b6SStefano Zampini /* ierr = MatSetLocalToGlobalMapping(*M,rmapping,cmapping);CHKERRQ(ierr); */ 330b7ce53b6SStefano Zampini /* Set values */ 331b7ce53b6SStefano Zampini if (isdense) { /* special case for dense local matrices */ 332b7ce53b6SStefano Zampini ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr); 333b7ce53b6SStefano Zampini ierr = MatDenseGetArray(matis->A,&array);CHKERRQ(ierr); 334b7ce53b6SStefano Zampini ierr = MatSetValues(*M,local_rows,global_indices,local_cols,global_indices,array,ADD_VALUES);CHKERRQ(ierr); 335b7ce53b6SStefano Zampini ierr = MatDenseRestoreArray(matis->A,&array);CHKERRQ(ierr); 336b7ce53b6SStefano Zampini ierr = PetscFree(local_indices);CHKERRQ(ierr); 337b7ce53b6SStefano Zampini ierr = PetscFree(global_indices);CHKERRQ(ierr); 338b7ce53b6SStefano Zampini } else { /* very basic values insertion for all other matrix types */ 339b7ce53b6SStefano Zampini ierr = PetscFree(local_indices);CHKERRQ(ierr); 340b7ce53b6SStefano Zampini for (i=0;i<local_rows;i++) { 341b7ce53b6SStefano Zampini ierr = MatGetRow(matis->A,i,&j,(const PetscInt**)&local_indices,(const PetscScalar**)&array);CHKERRQ(ierr); 342b7ce53b6SStefano Zampini /* ierr = MatSetValuesLocal(*M,1,&i,j,local_indices,array,ADD_VALUES);CHKERRQ(ierr); */ 343b7ce53b6SStefano Zampini ierr = ISLocalToGlobalMappingApply(matis->mapping,j,local_indices,global_indices);CHKERRQ(ierr); 344b7ce53b6SStefano Zampini ierr = ISLocalToGlobalMappingApply(matis->mapping,1,&i,&index_row);CHKERRQ(ierr); 345b7ce53b6SStefano Zampini ierr = MatSetValues(*M,1,&index_row,j,global_indices,array,ADD_VALUES);CHKERRQ(ierr); 346b7ce53b6SStefano Zampini ierr = MatRestoreRow(matis->A,i,&j,(const PetscInt**)&local_indices,(const PetscScalar**)&array);CHKERRQ(ierr); 347b7ce53b6SStefano Zampini } 348b7ce53b6SStefano Zampini ierr = PetscFree(global_indices);CHKERRQ(ierr); 349b7ce53b6SStefano Zampini } 350b7ce53b6SStefano Zampini ierr = MatAssemblyBegin(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 351b7ce53b6SStefano Zampini ierr = MatAssemblyEnd(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 352b7ce53b6SStefano Zampini if (isdense) { 353b7ce53b6SStefano Zampini ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr); 354b7ce53b6SStefano Zampini } 355b7ce53b6SStefano Zampini if (issbaij) { 356b7ce53b6SStefano Zampini ierr = MatRestoreRowUpperTriangular(matis->A);CHKERRQ(ierr); 357b7ce53b6SStefano Zampini } 358b7ce53b6SStefano Zampini PetscFunctionReturn(0); 359b7ce53b6SStefano Zampini } 360b7ce53b6SStefano Zampini 361b7ce53b6SStefano Zampini #undef __FUNCT__ 362b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ" 363b7ce53b6SStefano Zampini /*@ 364b7ce53b6SStefano Zampini MatISGetMPIXAIJ - Converts MATIS matrix into a parallel AIJ format 365b7ce53b6SStefano Zampini 366b7ce53b6SStefano Zampini Input Parameter: 367b7ce53b6SStefano Zampini . mat - the matrix (should be of type MATIS) 368b7ce53b6SStefano Zampini . reuse - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX 369b7ce53b6SStefano Zampini 370b7ce53b6SStefano Zampini Output Parameter: 371b7ce53b6SStefano Zampini . newmat - the matrix in AIJ format 372b7ce53b6SStefano Zampini 373b7ce53b6SStefano Zampini Level: developer 374b7ce53b6SStefano Zampini 375b7ce53b6SStefano Zampini Notes: 376b7ce53b6SStefano Zampini 377b7ce53b6SStefano Zampini .seealso: MATIS 378b7ce53b6SStefano Zampini @*/ 379b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ(Mat mat, MatReuse reuse, Mat *newmat) 380b7ce53b6SStefano Zampini { 381b7ce53b6SStefano Zampini PetscErrorCode ierr; 382b7ce53b6SStefano Zampini 383b7ce53b6SStefano Zampini PetscFunctionBegin; 384b7ce53b6SStefano Zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 385b7ce53b6SStefano Zampini PetscValidLogicalCollectiveEnum(mat,reuse,2); 386b7ce53b6SStefano Zampini PetscValidPointer(newmat,3); 387b7ce53b6SStefano Zampini if (reuse != MAT_INITIAL_MATRIX) { 388b7ce53b6SStefano Zampini PetscValidHeaderSpecific(*newmat,MAT_CLASSID,3); 389b7ce53b6SStefano Zampini PetscCheckSameComm(mat,1,*newmat,3); 390b7ce53b6SStefano Zampini } 391b7ce53b6SStefano Zampini ierr = PetscUseMethod(mat,"MatISGetMPIXAIJ_C",(Mat,MatReuse,Mat*),(mat,reuse,newmat));CHKERRQ(ierr); 392b7ce53b6SStefano Zampini PetscFunctionReturn(0); 393b7ce53b6SStefano Zampini } 394b7ce53b6SStefano Zampini 395b7ce53b6SStefano Zampini #undef __FUNCT__ 396ad6194a2SStefano Zampini #define __FUNCT__ "MatDuplicate_IS" 397ad6194a2SStefano Zampini PetscErrorCode MatDuplicate_IS(Mat mat,MatDuplicateOption op,Mat *newmat) 398ad6194a2SStefano Zampini { 399ad6194a2SStefano Zampini PetscErrorCode ierr; 400ad6194a2SStefano Zampini Mat_IS *matis = (Mat_IS*)(mat->data); 401ad6194a2SStefano Zampini PetscInt bs,m,n,M,N; 402ad6194a2SStefano Zampini Mat B,localmat; 403ad6194a2SStefano Zampini 404ad6194a2SStefano Zampini PetscFunctionBegin; 405ad6194a2SStefano Zampini ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 406ad6194a2SStefano Zampini ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 407ad6194a2SStefano Zampini ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 408ad6194a2SStefano Zampini ierr = MatCreateIS(PetscObjectComm((PetscObject)mat),bs,m,n,M,N,matis->mapping,&B);CHKERRQ(ierr); 409ad6194a2SStefano Zampini ierr = MatDuplicate(matis->A,op,&localmat);CHKERRQ(ierr); 410ad6194a2SStefano Zampini ierr = MatISSetLocalMat(B,localmat);CHKERRQ(ierr); 411ad6194a2SStefano Zampini ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 412ad6194a2SStefano Zampini ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 413ad6194a2SStefano Zampini *newmat = B; 414ad6194a2SStefano Zampini PetscFunctionReturn(0); 415ad6194a2SStefano Zampini } 416ad6194a2SStefano Zampini 417ad6194a2SStefano Zampini #undef __FUNCT__ 41869796d55SStefano Zampini #define __FUNCT__ "MatIsHermitian_IS" 41969796d55SStefano Zampini PetscErrorCode MatIsHermitian_IS(Mat A,PetscReal tol,PetscBool *flg) 42069796d55SStefano Zampini { 42169796d55SStefano Zampini PetscErrorCode ierr; 42269796d55SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 42369796d55SStefano Zampini PetscBool local_sym; 42469796d55SStefano Zampini 42569796d55SStefano Zampini PetscFunctionBegin; 42669796d55SStefano Zampini ierr = MatIsHermitian(matis->A,tol,&local_sym);CHKERRQ(ierr); 42769796d55SStefano Zampini ierr = MPI_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 42869796d55SStefano Zampini PetscFunctionReturn(0); 42969796d55SStefano Zampini } 43069796d55SStefano Zampini 43169796d55SStefano Zampini #undef __FUNCT__ 43269796d55SStefano Zampini #define __FUNCT__ "MatIsSymmetric_IS" 43369796d55SStefano Zampini PetscErrorCode MatIsSymmetric_IS(Mat A,PetscReal tol,PetscBool *flg) 43469796d55SStefano Zampini { 43569796d55SStefano Zampini PetscErrorCode ierr; 43669796d55SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 43769796d55SStefano Zampini PetscBool local_sym; 43869796d55SStefano Zampini 43969796d55SStefano Zampini PetscFunctionBegin; 44069796d55SStefano Zampini ierr = MatIsSymmetric(matis->A,tol,&local_sym);CHKERRQ(ierr); 44169796d55SStefano Zampini ierr = MPI_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 44269796d55SStefano Zampini PetscFunctionReturn(0); 44369796d55SStefano Zampini } 44469796d55SStefano Zampini 44569796d55SStefano Zampini #undef __FUNCT__ 446b4319ba4SBarry Smith #define __FUNCT__ "MatDestroy_IS" 447dfbe8321SBarry Smith PetscErrorCode MatDestroy_IS(Mat A) 448b4319ba4SBarry Smith { 449dfbe8321SBarry Smith PetscErrorCode ierr; 450b4319ba4SBarry Smith Mat_IS *b = (Mat_IS*)A->data; 451b4319ba4SBarry Smith 452b4319ba4SBarry Smith PetscFunctionBegin; 4536bf464f9SBarry Smith ierr = MatDestroy(&b->A);CHKERRQ(ierr); 4546bf464f9SBarry Smith ierr = VecScatterDestroy(&b->ctx);CHKERRQ(ierr); 4556bf464f9SBarry Smith ierr = VecDestroy(&b->x);CHKERRQ(ierr); 4566bf464f9SBarry Smith ierr = VecDestroy(&b->y);CHKERRQ(ierr); 4576bf464f9SBarry Smith ierr = ISLocalToGlobalMappingDestroy(&b->mapping);CHKERRQ(ierr); 458bf0cc555SLisandro Dalcin ierr = PetscFree(A->data);CHKERRQ(ierr); 459dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr); 460bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",NULL);CHKERRQ(ierr); 461b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",NULL);CHKERRQ(ierr); 462b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",NULL);CHKERRQ(ierr); 4632e1947a5SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",NULL);CHKERRQ(ierr); 464b4319ba4SBarry Smith PetscFunctionReturn(0); 465b4319ba4SBarry Smith } 466b4319ba4SBarry Smith 467b4319ba4SBarry Smith #undef __FUNCT__ 468b4319ba4SBarry Smith #define __FUNCT__ "MatMult_IS" 469dfbe8321SBarry Smith PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y) 470b4319ba4SBarry Smith { 471dfbe8321SBarry Smith PetscErrorCode ierr; 472b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 473b4319ba4SBarry Smith PetscScalar zero = 0.0; 474b4319ba4SBarry Smith 475b4319ba4SBarry Smith PetscFunctionBegin; 476b4319ba4SBarry Smith /* scatter the global vector x into the local work vector */ 477ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 478ca9f406cSSatish Balay ierr = VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 479b4319ba4SBarry Smith 480b4319ba4SBarry Smith /* multiply the local matrix */ 481b4319ba4SBarry Smith ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr); 482b4319ba4SBarry Smith 483b4319ba4SBarry Smith /* scatter product back into global memory */ 4842dcb1b2aSMatthew Knepley ierr = VecSet(y,zero);CHKERRQ(ierr); 485ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 486ca9f406cSSatish Balay ierr = VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 487b4319ba4SBarry Smith PetscFunctionReturn(0); 488b4319ba4SBarry Smith } 489b4319ba4SBarry Smith 490b4319ba4SBarry Smith #undef __FUNCT__ 4912e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultAdd_IS" 492b7ce53b6SStefano Zampini PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 4932e74eeadSLisandro Dalcin { 494650997f4SStefano Zampini Vec temp_vec; 4952e74eeadSLisandro Dalcin PetscErrorCode ierr; 4962e74eeadSLisandro Dalcin 4972e74eeadSLisandro Dalcin PetscFunctionBegin; /* v3 = v2 + A * v1.*/ 498650997f4SStefano Zampini if (v3 != v2) { 499650997f4SStefano Zampini ierr = MatMult(A,v1,v3);CHKERRQ(ierr); 500650997f4SStefano Zampini ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 501650997f4SStefano Zampini } else { 502650997f4SStefano Zampini ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 503650997f4SStefano Zampini ierr = MatMult(A,v1,temp_vec);CHKERRQ(ierr); 504650997f4SStefano Zampini ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 505650997f4SStefano Zampini ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 506650997f4SStefano Zampini ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 507650997f4SStefano Zampini } 5082e74eeadSLisandro Dalcin PetscFunctionReturn(0); 5092e74eeadSLisandro Dalcin } 5102e74eeadSLisandro Dalcin 5112e74eeadSLisandro Dalcin #undef __FUNCT__ 5122e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTranspose_IS" 5132e74eeadSLisandro Dalcin PetscErrorCode MatMultTranspose_IS(Mat A,Vec x,Vec y) 5142e74eeadSLisandro Dalcin { 5152e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 5162e74eeadSLisandro Dalcin PetscErrorCode ierr; 5172e74eeadSLisandro Dalcin 5182e74eeadSLisandro Dalcin PetscFunctionBegin; /* y = A' * x */ 5192e74eeadSLisandro Dalcin /* scatter the global vector x into the local work vector */ 520ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 521ca9f406cSSatish Balay ierr = VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5222e74eeadSLisandro Dalcin 5232e74eeadSLisandro Dalcin /* multiply the local matrix */ 5242e74eeadSLisandro Dalcin ierr = MatMultTranspose(is->A,is->x,is->y);CHKERRQ(ierr); 5252e74eeadSLisandro Dalcin 5262e74eeadSLisandro Dalcin /* scatter product back into global vector */ 5272e74eeadSLisandro Dalcin ierr = VecSet(y,0);CHKERRQ(ierr); 528ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 529ca9f406cSSatish Balay ierr = VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 5302e74eeadSLisandro Dalcin PetscFunctionReturn(0); 5312e74eeadSLisandro Dalcin } 5322e74eeadSLisandro Dalcin 5332e74eeadSLisandro Dalcin #undef __FUNCT__ 5342e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTransposeAdd_IS" 5352e74eeadSLisandro Dalcin PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 5362e74eeadSLisandro Dalcin { 537650997f4SStefano Zampini Vec temp_vec; 5382e74eeadSLisandro Dalcin PetscErrorCode ierr; 5392e74eeadSLisandro Dalcin 5402e74eeadSLisandro Dalcin PetscFunctionBegin; /* v3 = v2 + A' * v1.*/ 541650997f4SStefano Zampini if (v3 != v2) { 542650997f4SStefano Zampini ierr = MatMultTranspose(A,v1,v3);CHKERRQ(ierr); 543650997f4SStefano Zampini ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 544650997f4SStefano Zampini } else { 545650997f4SStefano Zampini ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 546650997f4SStefano Zampini ierr = MatMultTranspose(A,v1,temp_vec);CHKERRQ(ierr); 547650997f4SStefano Zampini ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 548650997f4SStefano Zampini ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 549650997f4SStefano Zampini ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 550650997f4SStefano Zampini } 5512e74eeadSLisandro Dalcin PetscFunctionReturn(0); 5522e74eeadSLisandro Dalcin } 5532e74eeadSLisandro Dalcin 5542e74eeadSLisandro Dalcin #undef __FUNCT__ 555b4319ba4SBarry Smith #define __FUNCT__ "MatView_IS" 556dfbe8321SBarry Smith PetscErrorCode MatView_IS(Mat A,PetscViewer viewer) 557b4319ba4SBarry Smith { 558b4319ba4SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 559dfbe8321SBarry Smith PetscErrorCode ierr; 560b4319ba4SBarry Smith PetscViewer sviewer; 561b4319ba4SBarry Smith 562b4319ba4SBarry Smith PetscFunctionBegin; 563dd2fa690SBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr); 564b4319ba4SBarry Smith ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 565b4319ba4SBarry Smith ierr = MatView(a->A,sviewer);CHKERRQ(ierr); 566b4319ba4SBarry Smith ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 567b4319ba4SBarry Smith PetscFunctionReturn(0); 568b4319ba4SBarry Smith } 569b4319ba4SBarry Smith 570b4319ba4SBarry Smith #undef __FUNCT__ 571b4319ba4SBarry Smith #define __FUNCT__ "MatSetLocalToGlobalMapping_IS" 572784ac674SJed Brown PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping) 573b4319ba4SBarry Smith { 574dfbe8321SBarry Smith PetscErrorCode ierr; 5754e4c7dbeSStefano Zampini PetscInt n,bs; 576b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 577b4319ba4SBarry Smith IS from,to; 578b4319ba4SBarry Smith Vec global; 579b4319ba4SBarry Smith 580b4319ba4SBarry Smith PetscFunctionBegin; 581784ac674SJed Brown PetscCheckSameComm(A,1,rmapping,2); 582ce94432eSBarry Smith if (rmapping != cmapping) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"MATIS requires the row and column mappings to be identical"); 58370cf5478SStefano Zampini if (is->mapping) { /* Currenly destroys the objects that will be created by this routine. Is there anything else that should be checked? */ 58470cf5478SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&is->mapping);CHKERRQ(ierr); 58570cf5478SStefano Zampini ierr = VecDestroy(&is->x);CHKERRQ(ierr); 58670cf5478SStefano Zampini ierr = VecDestroy(&is->y);CHKERRQ(ierr); 58770cf5478SStefano Zampini ierr = VecScatterDestroy(&is->ctx);CHKERRQ(ierr); 5881c47cb0fSStefano Zampini ierr = MatDestroy(&is->A);CHKERRQ(ierr); 58970cf5478SStefano Zampini } 590784ac674SJed Brown ierr = PetscObjectReference((PetscObject)rmapping);CHKERRQ(ierr); 5916bf464f9SBarry Smith ierr = ISLocalToGlobalMappingDestroy(&is->mapping);CHKERRQ(ierr); 592784ac674SJed Brown is->mapping = rmapping; 593fa7f1dd8SStefano Zampini /* 594fa7f1dd8SStefano Zampini ierr = PetscLayoutSetISLocalToGlobalMapping(A->rmap,rmapping);CHKERRQ(ierr); 595fa7f1dd8SStefano Zampini ierr = PetscLayoutSetISLocalToGlobalMapping(A->cmap,cmapping);CHKERRQ(ierr); 596fa7f1dd8SStefano Zampini */ 597b4319ba4SBarry Smith 598b4319ba4SBarry Smith /* Create the local matrix A */ 599784ac674SJed Brown ierr = ISLocalToGlobalMappingGetSize(rmapping,&n);CHKERRQ(ierr); 6002e1947a5SStefano Zampini ierr = ISLocalToGlobalMappingGetBlockSize(rmapping,&bs);CHKERRQ(ierr); 601f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr); 6022e1947a5SStefano Zampini if (bs > 1) { 6032e1947a5SStefano Zampini ierr = MatSetType(is->A,MATSEQBAIJ);CHKERRQ(ierr); 6042e1947a5SStefano Zampini } else { 6052e1947a5SStefano Zampini ierr = MatSetType(is->A,MATSEQAIJ);CHKERRQ(ierr); 6062e1947a5SStefano Zampini } 607f69a0ea3SMatthew Knepley ierr = MatSetSizes(is->A,n,n,n,n);CHKERRQ(ierr); 6084e4c7dbeSStefano Zampini ierr = MatSetBlockSize(is->A,bs);CHKERRQ(ierr); 609ff130e51SJed Brown ierr = MatSetOptionsPrefix(is->A,((PetscObject)A)->prefix);CHKERRQ(ierr); 610ff130e51SJed Brown ierr = MatAppendOptionsPrefix(is->A,"is_");CHKERRQ(ierr); 611b4319ba4SBarry Smith ierr = MatSetFromOptions(is->A);CHKERRQ(ierr); 612b4319ba4SBarry Smith 613b4319ba4SBarry Smith /* Create the local work vectors */ 6144e4c7dbeSStefano Zampini ierr = VecCreate(PETSC_COMM_SELF,&is->x);CHKERRQ(ierr); 6154e4c7dbeSStefano Zampini ierr = VecSetBlockSize(is->x,bs);CHKERRQ(ierr); 6164e4c7dbeSStefano Zampini ierr = VecSetSizes(is->x,n,n);CHKERRQ(ierr); 617ff130e51SJed Brown ierr = VecSetOptionsPrefix(is->x,((PetscObject)A)->prefix);CHKERRQ(ierr); 618ff130e51SJed Brown ierr = VecAppendOptionsPrefix(is->x,"is_");CHKERRQ(ierr); 6194e4c7dbeSStefano Zampini ierr = VecSetFromOptions(is->x);CHKERRQ(ierr); 620b4319ba4SBarry Smith ierr = VecDuplicate(is->x,&is->y);CHKERRQ(ierr); 621b4319ba4SBarry Smith 622b4319ba4SBarry Smith /* setup the global to local scatter */ 623b4319ba4SBarry Smith ierr = ISCreateStride(PETSC_COMM_SELF,n,0,1,&to);CHKERRQ(ierr); 624784ac674SJed Brown ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr); 6250298fd71SBarry Smith ierr = MatGetVecs(A,&global,NULL);CHKERRQ(ierr); 626b4319ba4SBarry Smith ierr = VecScatterCreate(global,from,is->x,to,&is->ctx);CHKERRQ(ierr); 6276bf464f9SBarry Smith ierr = VecDestroy(&global);CHKERRQ(ierr); 6286bf464f9SBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 6296bf464f9SBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 630b4319ba4SBarry Smith PetscFunctionReturn(0); 631b4319ba4SBarry Smith } 632b4319ba4SBarry Smith 6332e74eeadSLisandro Dalcin #undef __FUNCT__ 6342e74eeadSLisandro Dalcin #define __FUNCT__ "MatSetValues_IS" 6352e74eeadSLisandro Dalcin PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv) 6362e74eeadSLisandro Dalcin { 6372e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)mat->data; 6382e74eeadSLisandro Dalcin PetscInt rows_l[2048],cols_l[2048]; 6392e74eeadSLisandro Dalcin PetscErrorCode ierr; 6402e74eeadSLisandro Dalcin 6412e74eeadSLisandro Dalcin PetscFunctionBegin; 6422e74eeadSLisandro Dalcin #if defined(PETSC_USE_DEBUG) 643f23aa3ddSBarry 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); 6442e74eeadSLisandro Dalcin #endif 6452e74eeadSLisandro Dalcin ierr = ISG2LMapApply(is->mapping,m,rows,rows_l);CHKERRQ(ierr); 6462e74eeadSLisandro Dalcin ierr = ISG2LMapApply(is->mapping,n,cols,cols_l);CHKERRQ(ierr); 6472e74eeadSLisandro Dalcin ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 6482e74eeadSLisandro Dalcin PetscFunctionReturn(0); 6492e74eeadSLisandro Dalcin } 6502e74eeadSLisandro Dalcin 6512e74eeadSLisandro Dalcin #undef ISG2LMapSetUp 6522e74eeadSLisandro Dalcin #undef ISG2LMapApply 653b4319ba4SBarry Smith 654b4319ba4SBarry Smith #undef __FUNCT__ 655b4319ba4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_IS" 65613f74950SBarry Smith PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 657b4319ba4SBarry Smith { 658dfbe8321SBarry Smith PetscErrorCode ierr; 659b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 660b4319ba4SBarry Smith 661b4319ba4SBarry Smith PetscFunctionBegin; 662b4319ba4SBarry Smith ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 663b4319ba4SBarry Smith PetscFunctionReturn(0); 664b4319ba4SBarry Smith } 665b4319ba4SBarry Smith 666b4319ba4SBarry Smith #undef __FUNCT__ 667f0006bf2SLisandro Dalcin #define __FUNCT__ "MatSetValuesBlockedLocal_IS" 668f0006bf2SLisandro Dalcin PetscErrorCode MatSetValuesBlockedLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 669f0006bf2SLisandro Dalcin { 670f0006bf2SLisandro Dalcin PetscErrorCode ierr; 671f0006bf2SLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 672f0006bf2SLisandro Dalcin 673f0006bf2SLisandro Dalcin PetscFunctionBegin; 674f0006bf2SLisandro Dalcin ierr = MatSetValuesBlocked(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 675f0006bf2SLisandro Dalcin PetscFunctionReturn(0); 676f0006bf2SLisandro Dalcin } 677f0006bf2SLisandro Dalcin 678f0006bf2SLisandro Dalcin #undef __FUNCT__ 6792e74eeadSLisandro Dalcin #define __FUNCT__ "MatZeroRows_IS" 6802b40b63fSBarry Smith PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 6812e74eeadSLisandro Dalcin { 6822e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 6830298fd71SBarry Smith PetscInt n_l = 0, *rows_l = NULL; 6842e74eeadSLisandro Dalcin PetscErrorCode ierr; 6852e74eeadSLisandro Dalcin 6862e74eeadSLisandro Dalcin PetscFunctionBegin; 687ce94432eSBarry Smith if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support"); 6882e74eeadSLisandro Dalcin if (n) { 689785e854fSJed Brown ierr = PetscMalloc1(n,&rows_l);CHKERRQ(ierr); 6902e74eeadSLisandro Dalcin ierr = ISGlobalToLocalMappingApply(is->mapping,IS_GTOLM_DROP,n,rows,&n_l,rows_l);CHKERRQ(ierr); 6912e74eeadSLisandro Dalcin } 6922b40b63fSBarry Smith ierr = MatZeroRowsLocal(A,n_l,rows_l,diag,x,b);CHKERRQ(ierr); 693c31cb41cSBarry Smith ierr = PetscFree(rows_l);CHKERRQ(ierr); 6942e74eeadSLisandro Dalcin PetscFunctionReturn(0); 6952e74eeadSLisandro Dalcin } 6962e74eeadSLisandro Dalcin 6972e74eeadSLisandro Dalcin #undef __FUNCT__ 698b4319ba4SBarry Smith #define __FUNCT__ "MatZeroRowsLocal_IS" 6992b40b63fSBarry Smith PetscErrorCode MatZeroRowsLocal_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 700b4319ba4SBarry Smith { 701b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 702dfbe8321SBarry Smith PetscErrorCode ierr; 703f4df32b1SMatthew Knepley PetscInt i; 704b4319ba4SBarry Smith PetscScalar *array; 705b4319ba4SBarry Smith 706b4319ba4SBarry Smith PetscFunctionBegin; 707ce94432eSBarry Smith if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support"); 708b4319ba4SBarry Smith { 709b4319ba4SBarry Smith /* 710b4319ba4SBarry Smith Set up is->x as a "counting vector". This is in order to MatMult_IS 711b4319ba4SBarry Smith work properly in the interface nodes. 712b4319ba4SBarry Smith */ 713b4319ba4SBarry Smith Vec counter; 714b4319ba4SBarry Smith PetscScalar one=1.0, zero=0.0; 7150298fd71SBarry Smith ierr = MatGetVecs(A,&counter,NULL);CHKERRQ(ierr); 7162dcb1b2aSMatthew Knepley ierr = VecSet(counter,zero);CHKERRQ(ierr); 7172dcb1b2aSMatthew Knepley ierr = VecSet(is->x,one);CHKERRQ(ierr); 718ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 719ca9f406cSSatish Balay ierr = VecScatterEnd (is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 720ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 721ca9f406cSSatish Balay ierr = VecScatterEnd (is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 7226bf464f9SBarry Smith ierr = VecDestroy(&counter);CHKERRQ(ierr); 723b4319ba4SBarry Smith } 724958c9bccSBarry Smith if (!n) { 725b4319ba4SBarry Smith is->pure_neumann = PETSC_TRUE; 726b4319ba4SBarry Smith } else { 727b4319ba4SBarry Smith is->pure_neumann = PETSC_FALSE; 7282205254eSKarl Rupp 729b4319ba4SBarry Smith ierr = VecGetArray(is->x,&array);CHKERRQ(ierr); 7302b40b63fSBarry Smith ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr); 731b4319ba4SBarry Smith for (i=0; i<n; i++) { 732f4df32b1SMatthew Knepley ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr); 733b4319ba4SBarry Smith } 734b4319ba4SBarry Smith ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 735b4319ba4SBarry Smith ierr = MatAssemblyEnd (is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 736b4319ba4SBarry Smith ierr = VecRestoreArray(is->x,&array);CHKERRQ(ierr); 737b4319ba4SBarry Smith } 738b4319ba4SBarry Smith PetscFunctionReturn(0); 739b4319ba4SBarry Smith } 740b4319ba4SBarry Smith 741b4319ba4SBarry Smith #undef __FUNCT__ 742b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyBegin_IS" 743dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type) 744b4319ba4SBarry Smith { 745b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 746dfbe8321SBarry Smith PetscErrorCode ierr; 747b4319ba4SBarry Smith 748b4319ba4SBarry Smith PetscFunctionBegin; 749b4319ba4SBarry Smith ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr); 750b4319ba4SBarry Smith PetscFunctionReturn(0); 751b4319ba4SBarry Smith } 752b4319ba4SBarry Smith 753b4319ba4SBarry Smith #undef __FUNCT__ 754b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_IS" 755dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type) 756b4319ba4SBarry Smith { 757b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 758dfbe8321SBarry Smith PetscErrorCode ierr; 759b4319ba4SBarry Smith 760b4319ba4SBarry Smith PetscFunctionBegin; 761b4319ba4SBarry Smith ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr); 762b4319ba4SBarry Smith PetscFunctionReturn(0); 763b4319ba4SBarry Smith } 764b4319ba4SBarry Smith 765b4319ba4SBarry Smith #undef __FUNCT__ 766b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat_IS" 7677087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local) 768b4319ba4SBarry Smith { 769b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)mat->data; 770b4319ba4SBarry Smith 771b4319ba4SBarry Smith PetscFunctionBegin; 772b4319ba4SBarry Smith *local = is->A; 773b4319ba4SBarry Smith PetscFunctionReturn(0); 774b4319ba4SBarry Smith } 775b4319ba4SBarry Smith 776b4319ba4SBarry Smith #undef __FUNCT__ 777b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat" 778b4319ba4SBarry Smith /*@ 779b4319ba4SBarry Smith MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix. 780b4319ba4SBarry Smith 781b4319ba4SBarry Smith Input Parameter: 782b4319ba4SBarry Smith . mat - the matrix 783b4319ba4SBarry Smith 784b4319ba4SBarry Smith Output Parameter: 785b4319ba4SBarry Smith . local - the local matrix usually MATSEQAIJ 786b4319ba4SBarry Smith 787b4319ba4SBarry Smith Level: advanced 788b4319ba4SBarry Smith 789b4319ba4SBarry Smith Notes: 790b4319ba4SBarry Smith This can be called if you have precomputed the nonzero structure of the 791b4319ba4SBarry Smith matrix and want to provide it to the inner matrix object to improve the performance 792b4319ba4SBarry Smith of the MatSetValues() operation. 793b4319ba4SBarry Smith 794b4319ba4SBarry Smith .seealso: MATIS 795b4319ba4SBarry Smith @*/ 7967087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local) 797b4319ba4SBarry Smith { 7984ac538c5SBarry Smith PetscErrorCode ierr; 799b4319ba4SBarry Smith 800b4319ba4SBarry Smith PetscFunctionBegin; 8010700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 802b4319ba4SBarry Smith PetscValidPointer(local,2); 8034ac538c5SBarry Smith ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat*),(mat,local));CHKERRQ(ierr); 804b4319ba4SBarry Smith PetscFunctionReturn(0); 805b4319ba4SBarry Smith } 806b4319ba4SBarry Smith 8073b03a366Sstefano_zampini #undef __FUNCT__ 8083b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat_IS" 8093b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat_IS(Mat mat,Mat local) 8103b03a366Sstefano_zampini { 8113b03a366Sstefano_zampini Mat_IS *is = (Mat_IS*)mat->data; 8123b03a366Sstefano_zampini PetscInt nrows,ncols,orows,ocols; 8133b03a366Sstefano_zampini PetscErrorCode ierr; 8143b03a366Sstefano_zampini 8153b03a366Sstefano_zampini PetscFunctionBegin; 8164e4c7dbeSStefano Zampini if (is->A) { 8173b03a366Sstefano_zampini ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr); 8183b03a366Sstefano_zampini ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr); 819f23aa3ddSBarry 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); 8204e4c7dbeSStefano Zampini } 8213b03a366Sstefano_zampini ierr = PetscObjectReference((PetscObject)local);CHKERRQ(ierr); 8223b03a366Sstefano_zampini ierr = MatDestroy(&is->A);CHKERRQ(ierr); 8233b03a366Sstefano_zampini is->A = local; 8243b03a366Sstefano_zampini PetscFunctionReturn(0); 8253b03a366Sstefano_zampini } 8263b03a366Sstefano_zampini 8273b03a366Sstefano_zampini #undef __FUNCT__ 8283b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat" 8293b03a366Sstefano_zampini /*@ 8303b03a366Sstefano_zampini MatISSetLocalMat - Set the local matrix stored inside a MATIS matrix. 8313b03a366Sstefano_zampini 8323b03a366Sstefano_zampini Input Parameter: 8333b03a366Sstefano_zampini . mat - the matrix 8343b03a366Sstefano_zampini . local - the local matrix usually MATSEQAIJ 8353b03a366Sstefano_zampini 8363b03a366Sstefano_zampini Output Parameter: 8373b03a366Sstefano_zampini 8383b03a366Sstefano_zampini Level: advanced 8393b03a366Sstefano_zampini 8403b03a366Sstefano_zampini Notes: 8413b03a366Sstefano_zampini This can be called if you have precomputed the local matrix and 8423b03a366Sstefano_zampini want to provide it to the matrix object MATIS. 8433b03a366Sstefano_zampini 8443b03a366Sstefano_zampini .seealso: MATIS 8453b03a366Sstefano_zampini @*/ 8463b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat(Mat mat,Mat local) 8473b03a366Sstefano_zampini { 8483b03a366Sstefano_zampini PetscErrorCode ierr; 8493b03a366Sstefano_zampini 8503b03a366Sstefano_zampini PetscFunctionBegin; 8513b03a366Sstefano_zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 852b7ce53b6SStefano Zampini PetscValidHeaderSpecific(local,MAT_CLASSID,2); 8533b03a366Sstefano_zampini ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr); 8543b03a366Sstefano_zampini PetscFunctionReturn(0); 8553b03a366Sstefano_zampini } 8563b03a366Sstefano_zampini 8576726f965SBarry Smith #undef __FUNCT__ 8586726f965SBarry Smith #define __FUNCT__ "MatZeroEntries_IS" 8596726f965SBarry Smith PetscErrorCode MatZeroEntries_IS(Mat A) 8606726f965SBarry Smith { 8616726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 8626726f965SBarry Smith PetscErrorCode ierr; 8636726f965SBarry Smith 8646726f965SBarry Smith PetscFunctionBegin; 8656726f965SBarry Smith ierr = MatZeroEntries(a->A);CHKERRQ(ierr); 8666726f965SBarry Smith PetscFunctionReturn(0); 8676726f965SBarry Smith } 8686726f965SBarry Smith 8696726f965SBarry Smith #undef __FUNCT__ 8702e74eeadSLisandro Dalcin #define __FUNCT__ "MatScale_IS" 8712e74eeadSLisandro Dalcin PetscErrorCode MatScale_IS(Mat A,PetscScalar a) 8722e74eeadSLisandro Dalcin { 8732e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 8742e74eeadSLisandro Dalcin PetscErrorCode ierr; 8752e74eeadSLisandro Dalcin 8762e74eeadSLisandro Dalcin PetscFunctionBegin; 8772e74eeadSLisandro Dalcin ierr = MatScale(is->A,a);CHKERRQ(ierr); 8782e74eeadSLisandro Dalcin PetscFunctionReturn(0); 8792e74eeadSLisandro Dalcin } 8802e74eeadSLisandro Dalcin 8812e74eeadSLisandro Dalcin #undef __FUNCT__ 8822e74eeadSLisandro Dalcin #define __FUNCT__ "MatGetDiagonal_IS" 8832e74eeadSLisandro Dalcin PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v) 8842e74eeadSLisandro Dalcin { 8852e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 8862e74eeadSLisandro Dalcin PetscErrorCode ierr; 8872e74eeadSLisandro Dalcin 8882e74eeadSLisandro Dalcin PetscFunctionBegin; 8892e74eeadSLisandro Dalcin /* get diagonal of the local matrix */ 8902e74eeadSLisandro Dalcin ierr = MatGetDiagonal(is->A,is->x);CHKERRQ(ierr); 8912e74eeadSLisandro Dalcin 8922e74eeadSLisandro Dalcin /* scatter diagonal back into global vector */ 8932e74eeadSLisandro Dalcin ierr = VecSet(v,0);CHKERRQ(ierr); 894ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 895ca9f406cSSatish Balay ierr = VecScatterEnd (is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8962e74eeadSLisandro Dalcin PetscFunctionReturn(0); 8972e74eeadSLisandro Dalcin } 8982e74eeadSLisandro Dalcin 8992e74eeadSLisandro Dalcin #undef __FUNCT__ 9006726f965SBarry Smith #define __FUNCT__ "MatSetOption_IS" 901ace3abfcSBarry Smith PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg) 9026726f965SBarry Smith { 9036726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 9046726f965SBarry Smith PetscErrorCode ierr; 9056726f965SBarry Smith 9066726f965SBarry Smith PetscFunctionBegin; 9074e0d8c25SBarry Smith ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr); 9086726f965SBarry Smith PetscFunctionReturn(0); 9096726f965SBarry Smith } 9106726f965SBarry Smith 911284134d9SBarry Smith #undef __FUNCT__ 912284134d9SBarry Smith #define __FUNCT__ "MatCreateIS" 913284134d9SBarry Smith /*@ 914284134d9SBarry Smith MatCreateIS - Creates a "process" unassmembled matrix, it is assembled on each 915284134d9SBarry Smith process but not across processes. 916284134d9SBarry Smith 917284134d9SBarry Smith Input Parameters: 918284134d9SBarry Smith + comm - MPI communicator that will share the matrix 9194e4c7dbeSStefano Zampini . bs - local and global block size of the matrix 920df3898eeSBarry Smith . m,n,M,N - local and/or global sizes of the left and right vector used in matrix vector products 921284134d9SBarry Smith - map - mapping that defines the global number for each local number 922284134d9SBarry Smith 923284134d9SBarry Smith Output Parameter: 924284134d9SBarry Smith . A - the resulting matrix 925284134d9SBarry Smith 9268e6c10adSSatish Balay Level: advanced 9278e6c10adSSatish Balay 928284134d9SBarry Smith Notes: See MATIS for more details 9298cda6cd7SBarry Smith m and n are NOT related to the size of the map, they are the size of the part of the vector owned 9308cda6cd7SBarry Smith by that process. m + nghosts (or n + nghosts) is the length of map since map maps all local points 9318cda6cd7SBarry Smith plus the ghost points to global indices. 932284134d9SBarry Smith 933284134d9SBarry Smith .seealso: MATIS, MatSetLocalToGlobalMapping() 934284134d9SBarry Smith @*/ 9354e4c7dbeSStefano Zampini PetscErrorCode MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping map,Mat *A) 936284134d9SBarry Smith { 937284134d9SBarry Smith PetscErrorCode ierr; 938284134d9SBarry Smith 939284134d9SBarry Smith PetscFunctionBegin; 940284134d9SBarry Smith ierr = MatCreate(comm,A);CHKERRQ(ierr); 941d3ee2243SStefano Zampini ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr); 942284134d9SBarry Smith ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 943284134d9SBarry Smith ierr = MatSetType(*A,MATIS);CHKERRQ(ierr); 9443b03a366Sstefano_zampini ierr = MatSetUp(*A);CHKERRQ(ierr); 945784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(*A,map,map);CHKERRQ(ierr); 946284134d9SBarry Smith PetscFunctionReturn(0); 947284134d9SBarry Smith } 948284134d9SBarry Smith 949b4319ba4SBarry Smith /*MC 9506a818285SBarry Smith MATIS - MATIS = "is" - A matrix type to be used for using the Neumann-Neumann type preconditioners, see PCBDDC. 951b4319ba4SBarry Smith This stores the matrices in globally unassembled form. Each processor 952b4319ba4SBarry Smith assembles only its local Neumann problem and the parallel matrix vector 953b4319ba4SBarry Smith product is handled "implicitly". 954b4319ba4SBarry Smith 955b4319ba4SBarry Smith Operations Provided: 9566726f965SBarry Smith + MatMult() 9572e74eeadSLisandro Dalcin . MatMultAdd() 9582e74eeadSLisandro Dalcin . MatMultTranspose() 9592e74eeadSLisandro Dalcin . MatMultTransposeAdd() 9606726f965SBarry Smith . MatZeroEntries() 9616726f965SBarry Smith . MatSetOption() 9622e74eeadSLisandro Dalcin . MatZeroRows() 9636726f965SBarry Smith . MatZeroRowsLocal() 9642e74eeadSLisandro Dalcin . MatSetValues() 9656726f965SBarry Smith . MatSetValuesLocal() 9662e74eeadSLisandro Dalcin . MatScale() 9672e74eeadSLisandro Dalcin . MatGetDiagonal() 9686726f965SBarry Smith - MatSetLocalToGlobalMapping() 969b4319ba4SBarry Smith 970b4319ba4SBarry Smith Options Database Keys: 971b4319ba4SBarry Smith . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions() 972b4319ba4SBarry Smith 973b4319ba4SBarry Smith Notes: Options prefix for the inner matrix are given by -is_mat_xxx 974b4319ba4SBarry Smith 975b4319ba4SBarry Smith You must call MatSetLocalToGlobalMapping() before using this matrix type. 976b4319ba4SBarry Smith 977b4319ba4SBarry Smith You can do matrix preallocation on the local matrix after you obtain it with 978b4319ba4SBarry Smith MatISGetLocalMat() 979b4319ba4SBarry Smith 980b4319ba4SBarry Smith Level: advanced 981b4319ba4SBarry Smith 9826a818285SBarry Smith .seealso: PC, MatISGetLocalMat(), MatSetLocalToGlobalMapping(), PCBDDC 983b4319ba4SBarry Smith 984b4319ba4SBarry Smith M*/ 985b4319ba4SBarry Smith 986b4319ba4SBarry Smith #undef __FUNCT__ 987b4319ba4SBarry Smith #define __FUNCT__ "MatCreate_IS" 9888cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_IS(Mat A) 989b4319ba4SBarry Smith { 990dfbe8321SBarry Smith PetscErrorCode ierr; 991b4319ba4SBarry Smith Mat_IS *b; 992b4319ba4SBarry Smith 993b4319ba4SBarry Smith PetscFunctionBegin; 994b00a9115SJed Brown ierr = PetscNewLog(A,&b);CHKERRQ(ierr); 995b4319ba4SBarry Smith A->data = (void*)b; 996b4319ba4SBarry Smith ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 997b4319ba4SBarry Smith 998b4319ba4SBarry Smith A->ops->mult = MatMult_IS; 9992e74eeadSLisandro Dalcin A->ops->multadd = MatMultAdd_IS; 10002e74eeadSLisandro Dalcin A->ops->multtranspose = MatMultTranspose_IS; 10012e74eeadSLisandro Dalcin A->ops->multtransposeadd = MatMultTransposeAdd_IS; 1002b4319ba4SBarry Smith A->ops->destroy = MatDestroy_IS; 1003b4319ba4SBarry Smith A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS; 10042e74eeadSLisandro Dalcin A->ops->setvalues = MatSetValues_IS; 1005b4319ba4SBarry Smith A->ops->setvalueslocal = MatSetValuesLocal_IS; 1006f0006bf2SLisandro Dalcin A->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_IS; 10072e74eeadSLisandro Dalcin A->ops->zerorows = MatZeroRows_IS; 1008b4319ba4SBarry Smith A->ops->zerorowslocal = MatZeroRowsLocal_IS; 1009b4319ba4SBarry Smith A->ops->assemblybegin = MatAssemblyBegin_IS; 1010b4319ba4SBarry Smith A->ops->assemblyend = MatAssemblyEnd_IS; 1011b4319ba4SBarry Smith A->ops->view = MatView_IS; 10126726f965SBarry Smith A->ops->zeroentries = MatZeroEntries_IS; 10132e74eeadSLisandro Dalcin A->ops->scale = MatScale_IS; 10142e74eeadSLisandro Dalcin A->ops->getdiagonal = MatGetDiagonal_IS; 10156726f965SBarry Smith A->ops->setoption = MatSetOption_IS; 101669796d55SStefano Zampini A->ops->ishermitian = MatIsHermitian_IS; 101769796d55SStefano Zampini A->ops->issymmetric = MatIsSymmetric_IS; 1018ad6194a2SStefano Zampini A->ops->duplicate = MatDuplicate_IS; 1019b4319ba4SBarry Smith 102026283091SBarry Smith ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 102126283091SBarry Smith ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 1022b4319ba4SBarry Smith 1023b7ce53b6SStefano Zampini /* zeros pointer */ 1024b4319ba4SBarry Smith b->A = 0; 1025b4319ba4SBarry Smith b->ctx = 0; 1026b4319ba4SBarry Smith b->x = 0; 1027b4319ba4SBarry Smith b->y = 0; 1028b09366fdSStefano Zampini b->mapping = 0; 1029b7ce53b6SStefano Zampini 1030b7ce53b6SStefano Zampini /* special MATIS functions */ 1031bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",MatISGetLocalMat_IS);CHKERRQ(ierr); 1032bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",MatISSetLocalMat_IS);CHKERRQ(ierr); 1033b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",MatISGetMPIXAIJ_IS);CHKERRQ(ierr); 10342e1947a5SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",MatISSetPreallocation_IS);CHKERRQ(ierr); 103517667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr); 1036b4319ba4SBarry Smith PetscFunctionReturn(0); 1037b4319ba4SBarry Smith } 1038