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*2e1947a5SStefano Zampini #include <petscsf.h> 17*2e1947a5SStefano Zampini 18*2e1947a5SStefano Zampini #undef __FUNCT__ 19*2e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation" 20*2e1947a5SStefano Zampini PetscErrorCode MatISSetPreallocation(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[]) 21*2e1947a5SStefano Zampini { 22*2e1947a5SStefano Zampini PetscErrorCode ierr; 23*2e1947a5SStefano Zampini 24*2e1947a5SStefano Zampini PetscFunctionBegin; 25*2e1947a5SStefano Zampini PetscValidHeaderSpecific(B,MAT_CLASSID,1); 26*2e1947a5SStefano Zampini PetscValidType(B,1); 27*2e1947a5SStefano Zampini ierr = PetscTryMethod(B,"MatISSetPreallocation_C",(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(B,d_nz,d_nnz,o_nz,o_nnz));CHKERRQ(ierr); 28*2e1947a5SStefano Zampini PetscFunctionReturn(0); 29*2e1947a5SStefano Zampini } 30*2e1947a5SStefano Zampini 31*2e1947a5SStefano Zampini #undef __FUNCT__ 32*2e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation_IS" 33*2e1947a5SStefano Zampini PetscErrorCode MatISSetPreallocation_IS(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[]) 34*2e1947a5SStefano Zampini { 35*2e1947a5SStefano Zampini Mat_IS *matis = (Mat_IS*)(B->data); 36*2e1947a5SStefano Zampini PetscSF sf; 37*2e1947a5SStefano Zampini PetscInt bs,i,nroots,*rootdata,nleaves,*leafdata,nlocalcols; 38*2e1947a5SStefano Zampini const PetscInt *gidxs; 39*2e1947a5SStefano Zampini PetscErrorCode ierr; 40*2e1947a5SStefano Zampini 41*2e1947a5SStefano Zampini PetscFunctionBegin; 42*2e1947a5SStefano Zampini if (!matis->A) { 43*2e1947a5SStefano Zampini SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"You should first call MatSetLocalToGlobalMapping"); 44*2e1947a5SStefano Zampini } 45*2e1947a5SStefano Zampini ierr = MatGetLocalSize(B,&nroots,NULL);CHKERRQ(ierr); 46*2e1947a5SStefano Zampini ierr = MatGetSize(matis->A,&nleaves,&nlocalcols);CHKERRQ(ierr); 47*2e1947a5SStefano Zampini ierr = MatGetBlockSize(matis->A,&bs);CHKERRQ(ierr); 48*2e1947a5SStefano Zampini ierr = PetscCalloc2(nroots,&rootdata,nleaves,&leafdata);CHKERRQ(ierr); 49*2e1947a5SStefano Zampini ierr = PetscSFCreate(PetscObjectComm((PetscObject)B),&sf);CHKERRQ(ierr); 50*2e1947a5SStefano Zampini ierr = PetscSFSetFromOptions(sf);CHKERRQ(ierr); 51*2e1947a5SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(matis->mapping,&gidxs);CHKERRQ(ierr); 52*2e1947a5SStefano Zampini ierr = PetscSFSetGraphLayout(sf,B->rmap,nleaves,NULL,PETSC_COPY_VALUES,gidxs);CHKERRQ(ierr); 53*2e1947a5SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(matis->mapping,&gidxs);CHKERRQ(ierr); 54*2e1947a5SStefano Zampini if (!d_nnz) { 55*2e1947a5SStefano Zampini for (i=0;i<nroots;i++) rootdata[i] += d_nz; 56*2e1947a5SStefano Zampini } else { 57*2e1947a5SStefano Zampini for (i=0;i<nroots;i++) rootdata[i] += d_nnz[i]; 58*2e1947a5SStefano Zampini } 59*2e1947a5SStefano Zampini if (!o_nnz) { 60*2e1947a5SStefano Zampini for (i=0;i<nroots;i++) rootdata[i] += o_nz; 61*2e1947a5SStefano Zampini } else { 62*2e1947a5SStefano Zampini for (i=0;i<nroots;i++) rootdata[i] += o_nnz[i]; 63*2e1947a5SStefano Zampini } 64*2e1947a5SStefano Zampini ierr = PetscSFBcastBegin(sf,MPIU_INT,rootdata,leafdata);CHKERRQ(ierr); 65*2e1947a5SStefano Zampini ierr = PetscSFBcastEnd(sf,MPIU_INT,rootdata,leafdata);CHKERRQ(ierr); 66*2e1947a5SStefano Zampini for (i=0;i<nleaves;i++) { 67*2e1947a5SStefano Zampini leafdata[i] = PetscMin(leafdata[i],nlocalcols); 68*2e1947a5SStefano Zampini } 69*2e1947a5SStefano Zampini ierr = MatSeqAIJSetPreallocation(matis->A,0,leafdata);CHKERRQ(ierr); 70*2e1947a5SStefano Zampini for (i=0;i<nleaves/bs;i++) { 71*2e1947a5SStefano Zampini leafdata[i] = leafdata[i*bs]/bs; 72*2e1947a5SStefano Zampini } 73*2e1947a5SStefano Zampini ierr = MatSeqBAIJSetPreallocation(matis->A,bs,0,leafdata);CHKERRQ(ierr); 74*2e1947a5SStefano Zampini for (i=0;i<nleaves/bs;i++) { 75*2e1947a5SStefano Zampini leafdata[i] = leafdata[i]-i; 76*2e1947a5SStefano Zampini } 77*2e1947a5SStefano Zampini ierr = MatSeqSBAIJSetPreallocation(matis->A,bs,0,leafdata);CHKERRQ(ierr); 78*2e1947a5SStefano Zampini ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 79*2e1947a5SStefano Zampini ierr = PetscFree2(rootdata,leafdata);CHKERRQ(ierr); 80*2e1947a5SStefano Zampini PetscFunctionReturn(0); 81*2e1947a5SStefano Zampini } 82b4319ba4SBarry Smith 83b4319ba4SBarry Smith #undef __FUNCT__ 84b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ_IS" 85b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ_IS(Mat mat, MatReuse reuse, Mat *M) 86b7ce53b6SStefano Zampini { 87b7ce53b6SStefano Zampini Mat_IS *matis = (Mat_IS*)(mat->data); 88b7ce53b6SStefano Zampini /* info on mat */ 89b7ce53b6SStefano Zampini /* ISLocalToGlobalMapping rmapping,cmapping; */ 90b7ce53b6SStefano Zampini PetscInt bs,rows,cols; 91b7ce53b6SStefano Zampini PetscInt lrows,lcols; 92b7ce53b6SStefano Zampini PetscInt local_rows,local_cols; 93b7ce53b6SStefano Zampini PetscBool isdense,issbaij,issbaij_red; 94b7ce53b6SStefano Zampini /* values insertion */ 95b7ce53b6SStefano Zampini PetscScalar *array; 96b7ce53b6SStefano Zampini PetscInt *local_indices,*global_indices; 97b7ce53b6SStefano Zampini /* work */ 98b7ce53b6SStefano Zampini PetscInt i,j,index_row; 99b7ce53b6SStefano Zampini PetscErrorCode ierr; 100b7ce53b6SStefano Zampini 101b7ce53b6SStefano Zampini PetscFunctionBegin; 102b7ce53b6SStefano Zampini /* MISSING CHECKS 103b7ce53b6SStefano Zampini - rectangular case not covered (it is not allowed by MATIS) 104b7ce53b6SStefano Zampini */ 105b7ce53b6SStefano Zampini /* get info from mat */ 106b7ce53b6SStefano Zampini /* ierr = MatGetLocalToGlobalMapping(mat,&rmapping,&cmapping);CHKERRQ(ierr); */ 107b7ce53b6SStefano Zampini ierr = MatGetSize(mat,&rows,&cols);CHKERRQ(ierr); 108b7ce53b6SStefano Zampini ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 109b7ce53b6SStefano Zampini ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr); 110b7ce53b6SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr); 111b7ce53b6SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr); 112b7ce53b6SStefano Zampini 113b7ce53b6SStefano Zampini /* work */ 114b7ce53b6SStefano Zampini ierr = PetscMalloc1(local_rows,&local_indices);CHKERRQ(ierr); 115b7ce53b6SStefano Zampini for (i=0;i<local_rows;i++) local_indices[i]=i; 116b7ce53b6SStefano Zampini /* map indices of local mat to global values */ 117b7ce53b6SStefano Zampini ierr = PetscMalloc(PetscMax(local_cols,local_rows)*sizeof(*global_indices),&global_indices);CHKERRQ(ierr); 118b7ce53b6SStefano Zampini /* ierr = ISLocalToGlobalMappingApply(rmapping,local_rows,local_indices,global_indices);CHKERRQ(ierr); */ 119b7ce53b6SStefano Zampini ierr = ISLocalToGlobalMappingApply(matis->mapping,local_rows,local_indices,global_indices);CHKERRQ(ierr); 120b7ce53b6SStefano Zampini 121b7ce53b6SStefano Zampini if (issbaij) { 122b7ce53b6SStefano Zampini ierr = MatGetRowUpperTriangular(matis->A);CHKERRQ(ierr); 123b7ce53b6SStefano Zampini } 124b7ce53b6SStefano Zampini 125b7ce53b6SStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { 126b7ce53b6SStefano Zampini Mat new_mat; 127b7ce53b6SStefano Zampini MatType new_mat_type; 128b7ce53b6SStefano Zampini Vec vec_dnz,vec_onz; 129b7ce53b6SStefano Zampini PetscScalar *my_dnz,*my_onz; 130b7ce53b6SStefano Zampini PetscInt *dnz,*onz,*mat_ranges,*row_ownership; 131b7ce53b6SStefano Zampini PetscInt index_col,owner; 132b7ce53b6SStefano Zampini PetscMPIInt nsubdomains; 133b7ce53b6SStefano Zampini 134b7ce53b6SStefano Zampini /* determining new matrix type */ 135b7ce53b6SStefano Zampini ierr = MPI_Allreduce(&issbaij,&issbaij_red,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 136b7ce53b6SStefano Zampini if (issbaij_red) { 137b7ce53b6SStefano Zampini new_mat_type = MATSBAIJ; 138b7ce53b6SStefano Zampini } else { 139b7ce53b6SStefano Zampini if (bs>1) { 140b7ce53b6SStefano Zampini new_mat_type = MATBAIJ; 141b7ce53b6SStefano Zampini } else { 142b7ce53b6SStefano Zampini new_mat_type = MATAIJ; 143b7ce53b6SStefano Zampini } 144b7ce53b6SStefano Zampini } 145b7ce53b6SStefano Zampini 146b7ce53b6SStefano Zampini ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&nsubdomains);CHKERRQ(ierr); 147b7ce53b6SStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)mat),&new_mat);CHKERRQ(ierr); 148b7ce53b6SStefano Zampini ierr = MatSetSizes(new_mat,PETSC_DECIDE,PETSC_DECIDE,rows,cols);CHKERRQ(ierr); 149b7ce53b6SStefano Zampini ierr = MatSetBlockSize(new_mat,bs);CHKERRQ(ierr); 150b7ce53b6SStefano Zampini ierr = MatSetType(new_mat,new_mat_type);CHKERRQ(ierr); 151b7ce53b6SStefano Zampini ierr = MatSetUp(new_mat);CHKERRQ(ierr); 152b7ce53b6SStefano Zampini ierr = MatGetLocalSize(new_mat,&lrows,&lcols);CHKERRQ(ierr); 153b7ce53b6SStefano Zampini 154b7ce53b6SStefano Zampini /* 155b7ce53b6SStefano Zampini preallocation 156b7ce53b6SStefano Zampini */ 157b7ce53b6SStefano Zampini 158b7ce53b6SStefano Zampini ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)new_mat),lrows,lcols,dnz,onz);CHKERRQ(ierr); 159b7ce53b6SStefano Zampini /* 160b7ce53b6SStefano Zampini Some vectors are needed to sum up properly on shared interface dofs. 161b7ce53b6SStefano Zampini Preallocation macros cannot do the job. 162b7ce53b6SStefano Zampini Note that preallocation is not exact, since it overestimates nonzeros 163b7ce53b6SStefano Zampini */ 164b7ce53b6SStefano Zampini ierr = MatGetVecs(new_mat,NULL,&vec_dnz);CHKERRQ(ierr); 165b7ce53b6SStefano Zampini /* ierr = VecSetLocalToGlobalMapping(vec_dnz,rmapping);CHKERRQ(ierr); */ 166b7ce53b6SStefano Zampini ierr = VecSetLocalToGlobalMapping(vec_dnz,matis->mapping);CHKERRQ(ierr); 167b7ce53b6SStefano Zampini ierr = VecDuplicate(vec_dnz,&vec_onz);CHKERRQ(ierr); 168b7ce53b6SStefano Zampini /* All processes need to compute entire row ownership */ 169b7ce53b6SStefano Zampini ierr = PetscMalloc1(rows,&row_ownership);CHKERRQ(ierr); 170b7ce53b6SStefano Zampini ierr = MatGetOwnershipRanges(new_mat,(const PetscInt**)&mat_ranges);CHKERRQ(ierr); 171b7ce53b6SStefano Zampini for (i=0;i<nsubdomains;i++) { 172b7ce53b6SStefano Zampini for (j=mat_ranges[i];j<mat_ranges[i+1];j++) { 173b7ce53b6SStefano Zampini row_ownership[j]=i; 174b7ce53b6SStefano Zampini } 175b7ce53b6SStefano Zampini } 176b7ce53b6SStefano Zampini 177b7ce53b6SStefano Zampini /* 178b7ce53b6SStefano Zampini my_dnz and my_onz contains exact contribution to preallocation from each local mat 179b7ce53b6SStefano Zampini then, they will be summed up properly. This way, preallocation is always sufficient 180b7ce53b6SStefano Zampini */ 181b7ce53b6SStefano Zampini ierr = PetscMalloc1(local_rows,&my_dnz);CHKERRQ(ierr); 182b7ce53b6SStefano Zampini ierr = PetscMalloc1(local_rows,&my_onz);CHKERRQ(ierr); 183b7ce53b6SStefano Zampini ierr = PetscMemzero(my_dnz,local_rows*sizeof(*my_dnz));CHKERRQ(ierr); 184b7ce53b6SStefano Zampini ierr = PetscMemzero(my_onz,local_rows*sizeof(*my_onz));CHKERRQ(ierr); 185b7ce53b6SStefano Zampini /* preallocation as a MATAIJ */ 186b7ce53b6SStefano Zampini if (isdense) { /* special case for dense local matrices */ 187b7ce53b6SStefano Zampini for (i=0;i<local_rows;i++) { 188b7ce53b6SStefano Zampini index_row = global_indices[i]; 189b7ce53b6SStefano Zampini for (j=i;j<local_rows;j++) { 190b7ce53b6SStefano Zampini owner = row_ownership[index_row]; 191b7ce53b6SStefano Zampini index_col = global_indices[j]; 192b7ce53b6SStefano Zampini if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */ 193b7ce53b6SStefano Zampini my_dnz[i] += 1.0; 194b7ce53b6SStefano Zampini } else { /* offdiag block */ 195b7ce53b6SStefano Zampini my_onz[i] += 1.0; 196b7ce53b6SStefano Zampini } 197b7ce53b6SStefano Zampini /* same as before, interchanging rows and cols */ 198b7ce53b6SStefano Zampini if (i != j) { 199b7ce53b6SStefano Zampini owner = row_ownership[index_col]; 200b7ce53b6SStefano Zampini if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) { 201b7ce53b6SStefano Zampini my_dnz[j] += 1.0; 202b7ce53b6SStefano Zampini } else { 203b7ce53b6SStefano Zampini my_onz[j] += 1.0; 204b7ce53b6SStefano Zampini } 205b7ce53b6SStefano Zampini } 206b7ce53b6SStefano Zampini } 207b7ce53b6SStefano Zampini } 208b7ce53b6SStefano Zampini } else { 209b7ce53b6SStefano Zampini for (i=0;i<local_rows;i++) { 210b7ce53b6SStefano Zampini PetscInt ncols; 211b7ce53b6SStefano Zampini const PetscInt *cols; 212b7ce53b6SStefano Zampini index_row = global_indices[i]; 213b7ce53b6SStefano Zampini ierr = MatGetRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr); 214b7ce53b6SStefano Zampini for (j=0;j<ncols;j++) { 215b7ce53b6SStefano Zampini owner = row_ownership[index_row]; 216b7ce53b6SStefano Zampini index_col = global_indices[cols[j]]; 217b7ce53b6SStefano Zampini if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */ 218b7ce53b6SStefano Zampini my_dnz[i] += 1.0; 219b7ce53b6SStefano Zampini } else { /* offdiag block */ 220b7ce53b6SStefano Zampini my_onz[i] += 1.0; 221b7ce53b6SStefano Zampini } 222b7ce53b6SStefano Zampini /* same as before, interchanging rows and cols */ 223b7ce53b6SStefano Zampini if (issbaij) { 224b7ce53b6SStefano Zampini owner = row_ownership[index_col]; 225b7ce53b6SStefano Zampini if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) { 226b7ce53b6SStefano Zampini my_dnz[j] += 1.0; 227b7ce53b6SStefano Zampini } else { 228b7ce53b6SStefano Zampini my_onz[j] += 1.0; 229b7ce53b6SStefano Zampini } 230b7ce53b6SStefano Zampini } 231b7ce53b6SStefano Zampini } 232b7ce53b6SStefano Zampini ierr = MatRestoreRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr); 233b7ce53b6SStefano Zampini } 234b7ce53b6SStefano Zampini } 235b7ce53b6SStefano Zampini ierr = VecSet(vec_dnz,0.0);CHKERRQ(ierr); 236b7ce53b6SStefano Zampini ierr = VecSet(vec_onz,0.0);CHKERRQ(ierr); 237b7ce53b6SStefano Zampini if (local_rows) { /* multilevel guard */ 238b7ce53b6SStefano Zampini ierr = VecSetValuesLocal(vec_dnz,local_rows,local_indices,my_dnz,ADD_VALUES);CHKERRQ(ierr); 239b7ce53b6SStefano Zampini ierr = VecSetValuesLocal(vec_onz,local_rows,local_indices,my_onz,ADD_VALUES);CHKERRQ(ierr); 240b7ce53b6SStefano Zampini } 241b7ce53b6SStefano Zampini ierr = VecAssemblyBegin(vec_dnz);CHKERRQ(ierr); 242b7ce53b6SStefano Zampini ierr = VecAssemblyBegin(vec_onz);CHKERRQ(ierr); 243b7ce53b6SStefano Zampini ierr = VecAssemblyEnd(vec_dnz);CHKERRQ(ierr); 244b7ce53b6SStefano Zampini ierr = VecAssemblyEnd(vec_onz);CHKERRQ(ierr); 245b7ce53b6SStefano Zampini ierr = PetscFree(my_dnz);CHKERRQ(ierr); 246b7ce53b6SStefano Zampini ierr = PetscFree(my_onz);CHKERRQ(ierr); 247b7ce53b6SStefano Zampini ierr = PetscFree(row_ownership);CHKERRQ(ierr); 248b7ce53b6SStefano Zampini 249b7ce53b6SStefano Zampini /* set computed preallocation in dnz and onz */ 250b7ce53b6SStefano Zampini ierr = VecGetArray(vec_dnz,&array);CHKERRQ(ierr); 251b7ce53b6SStefano Zampini for (i=0; i<lrows; i++) dnz[i] = (PetscInt)PetscRealPart(array[i]); 252b7ce53b6SStefano Zampini ierr = VecRestoreArray(vec_dnz,&array);CHKERRQ(ierr); 253b7ce53b6SStefano Zampini ierr = VecGetArray(vec_onz,&array);CHKERRQ(ierr); 254b7ce53b6SStefano Zampini for (i=0;i<lrows;i++) onz[i] = (PetscInt)PetscRealPart(array[i]); 255b7ce53b6SStefano Zampini ierr = VecRestoreArray(vec_onz,&array);CHKERRQ(ierr); 256b7ce53b6SStefano Zampini ierr = VecDestroy(&vec_dnz);CHKERRQ(ierr); 257b7ce53b6SStefano Zampini ierr = VecDestroy(&vec_onz);CHKERRQ(ierr); 258b7ce53b6SStefano Zampini 259b7ce53b6SStefano Zampini /* Resize preallocation if overestimated */ 260b7ce53b6SStefano Zampini for (i=0;i<lrows;i++) { 261b7ce53b6SStefano Zampini dnz[i] = PetscMin(dnz[i],lcols); 262b7ce53b6SStefano Zampini onz[i] = PetscMin(onz[i],cols-lcols); 263b7ce53b6SStefano Zampini } 264b7ce53b6SStefano Zampini /* set preallocation */ 265b7ce53b6SStefano Zampini ierr = MatMPIAIJSetPreallocation(new_mat,0,dnz,0,onz);CHKERRQ(ierr); 266b7ce53b6SStefano Zampini for (i=0;i<lrows/bs;i++) { 267b7ce53b6SStefano Zampini dnz[i] = dnz[i*bs]/bs; 268b7ce53b6SStefano Zampini onz[i] = onz[i*bs]/bs; 269b7ce53b6SStefano Zampini } 270b7ce53b6SStefano Zampini ierr = MatMPIBAIJSetPreallocation(new_mat,bs,0,dnz,0,onz);CHKERRQ(ierr); 271b7ce53b6SStefano Zampini for (i=0;i<lrows/bs;i++) { 272b7ce53b6SStefano Zampini dnz[i] = dnz[i]-i; 273b7ce53b6SStefano Zampini } 274b7ce53b6SStefano Zampini ierr = MatMPISBAIJSetPreallocation(new_mat,bs,0,dnz,0,onz);CHKERRQ(ierr); 275b7ce53b6SStefano Zampini ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 276b7ce53b6SStefano Zampini *M = new_mat; 277b7ce53b6SStefano Zampini } else { 278b7ce53b6SStefano Zampini PetscInt mbs,mrows,mcols; 279b7ce53b6SStefano Zampini /* some checks */ 280b7ce53b6SStefano Zampini ierr = MatGetBlockSize(*M,&mbs);CHKERRQ(ierr); 281b7ce53b6SStefano Zampini ierr = MatGetSize(*M,&mrows,&mcols);CHKERRQ(ierr); 282b7ce53b6SStefano Zampini if (mrows != rows) { 283b7ce53b6SStefano Zampini SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of rows (%d != %d)",rows,mrows); 284b7ce53b6SStefano Zampini } 285b7ce53b6SStefano Zampini if (mrows != rows) { 286b7ce53b6SStefano Zampini SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of cols (%d != %d)",cols,mcols); 287b7ce53b6SStefano Zampini } 288b7ce53b6SStefano Zampini if (mbs != bs) { 289b7ce53b6SStefano Zampini SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong block size (%d != %d)",bs,mbs); 290b7ce53b6SStefano Zampini } 291b7ce53b6SStefano Zampini ierr = MatZeroEntries(*M);CHKERRQ(ierr); 292b7ce53b6SStefano Zampini } 293b7ce53b6SStefano Zampini /* set local to global mappings */ 294b7ce53b6SStefano Zampini /* ierr = MatSetLocalToGlobalMapping(*M,rmapping,cmapping);CHKERRQ(ierr); */ 295b7ce53b6SStefano Zampini /* Set values */ 296b7ce53b6SStefano Zampini if (isdense) { /* special case for dense local matrices */ 297b7ce53b6SStefano Zampini ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr); 298b7ce53b6SStefano Zampini ierr = MatDenseGetArray(matis->A,&array);CHKERRQ(ierr); 299b7ce53b6SStefano Zampini ierr = MatSetValues(*M,local_rows,global_indices,local_cols,global_indices,array,ADD_VALUES);CHKERRQ(ierr); 300b7ce53b6SStefano Zampini ierr = MatDenseRestoreArray(matis->A,&array);CHKERRQ(ierr); 301b7ce53b6SStefano Zampini ierr = PetscFree(local_indices);CHKERRQ(ierr); 302b7ce53b6SStefano Zampini ierr = PetscFree(global_indices);CHKERRQ(ierr); 303b7ce53b6SStefano Zampini } else { /* very basic values insertion for all other matrix types */ 304b7ce53b6SStefano Zampini ierr = PetscFree(local_indices);CHKERRQ(ierr); 305b7ce53b6SStefano Zampini for (i=0;i<local_rows;i++) { 306b7ce53b6SStefano Zampini ierr = MatGetRow(matis->A,i,&j,(const PetscInt**)&local_indices,(const PetscScalar**)&array);CHKERRQ(ierr); 307b7ce53b6SStefano Zampini /* ierr = MatSetValuesLocal(*M,1,&i,j,local_indices,array,ADD_VALUES);CHKERRQ(ierr); */ 308b7ce53b6SStefano Zampini ierr = ISLocalToGlobalMappingApply(matis->mapping,j,local_indices,global_indices);CHKERRQ(ierr); 309b7ce53b6SStefano Zampini ierr = ISLocalToGlobalMappingApply(matis->mapping,1,&i,&index_row);CHKERRQ(ierr); 310b7ce53b6SStefano Zampini ierr = MatSetValues(*M,1,&index_row,j,global_indices,array,ADD_VALUES);CHKERRQ(ierr); 311b7ce53b6SStefano Zampini ierr = MatRestoreRow(matis->A,i,&j,(const PetscInt**)&local_indices,(const PetscScalar**)&array);CHKERRQ(ierr); 312b7ce53b6SStefano Zampini } 313b7ce53b6SStefano Zampini ierr = PetscFree(global_indices);CHKERRQ(ierr); 314b7ce53b6SStefano Zampini } 315b7ce53b6SStefano Zampini ierr = MatAssemblyBegin(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 316b7ce53b6SStefano Zampini ierr = MatAssemblyEnd(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 317b7ce53b6SStefano Zampini if (isdense) { 318b7ce53b6SStefano Zampini ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr); 319b7ce53b6SStefano Zampini } 320b7ce53b6SStefano Zampini if (issbaij) { 321b7ce53b6SStefano Zampini ierr = MatRestoreRowUpperTriangular(matis->A);CHKERRQ(ierr); 322b7ce53b6SStefano Zampini } 323b7ce53b6SStefano Zampini PetscFunctionReturn(0); 324b7ce53b6SStefano Zampini } 325b7ce53b6SStefano Zampini 326b7ce53b6SStefano Zampini #undef __FUNCT__ 327b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ" 328b7ce53b6SStefano Zampini /*@ 329b7ce53b6SStefano Zampini MatISGetMPIXAIJ - Converts MATIS matrix into a parallel AIJ format 330b7ce53b6SStefano Zampini 331b7ce53b6SStefano Zampini Input Parameter: 332b7ce53b6SStefano Zampini . mat - the matrix (should be of type MATIS) 333b7ce53b6SStefano Zampini . reuse - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX 334b7ce53b6SStefano Zampini 335b7ce53b6SStefano Zampini Output Parameter: 336b7ce53b6SStefano Zampini . newmat - the matrix in AIJ format 337b7ce53b6SStefano Zampini 338b7ce53b6SStefano Zampini Level: developer 339b7ce53b6SStefano Zampini 340b7ce53b6SStefano Zampini Notes: 341b7ce53b6SStefano Zampini 342b7ce53b6SStefano Zampini .seealso: MATIS 343b7ce53b6SStefano Zampini @*/ 344b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ(Mat mat, MatReuse reuse, Mat *newmat) 345b7ce53b6SStefano Zampini { 346b7ce53b6SStefano Zampini PetscErrorCode ierr; 347b7ce53b6SStefano Zampini 348b7ce53b6SStefano Zampini PetscFunctionBegin; 349b7ce53b6SStefano Zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 350b7ce53b6SStefano Zampini PetscValidLogicalCollectiveEnum(mat,reuse,2); 351b7ce53b6SStefano Zampini PetscValidPointer(newmat,3); 352b7ce53b6SStefano Zampini if (reuse != MAT_INITIAL_MATRIX) { 353b7ce53b6SStefano Zampini PetscValidHeaderSpecific(*newmat,MAT_CLASSID,3); 354b7ce53b6SStefano Zampini PetscCheckSameComm(mat,1,*newmat,3); 355b7ce53b6SStefano Zampini } 356b7ce53b6SStefano Zampini ierr = PetscUseMethod(mat,"MatISGetMPIXAIJ_C",(Mat,MatReuse,Mat*),(mat,reuse,newmat));CHKERRQ(ierr); 357b7ce53b6SStefano Zampini PetscFunctionReturn(0); 358b7ce53b6SStefano Zampini } 359b7ce53b6SStefano Zampini 360b7ce53b6SStefano Zampini #undef __FUNCT__ 361ad6194a2SStefano Zampini #define __FUNCT__ "MatDuplicate_IS" 362ad6194a2SStefano Zampini PetscErrorCode MatDuplicate_IS(Mat mat,MatDuplicateOption op,Mat *newmat) 363ad6194a2SStefano Zampini { 364ad6194a2SStefano Zampini PetscErrorCode ierr; 365ad6194a2SStefano Zampini Mat_IS *matis = (Mat_IS*)(mat->data); 366ad6194a2SStefano Zampini PetscInt bs,m,n,M,N; 367ad6194a2SStefano Zampini Mat B,localmat; 368ad6194a2SStefano Zampini 369ad6194a2SStefano Zampini PetscFunctionBegin; 370ad6194a2SStefano Zampini ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 371ad6194a2SStefano Zampini ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 372ad6194a2SStefano Zampini ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 373ad6194a2SStefano Zampini ierr = MatCreateIS(PetscObjectComm((PetscObject)mat),bs,m,n,M,N,matis->mapping,&B);CHKERRQ(ierr); 374ad6194a2SStefano Zampini ierr = MatDuplicate(matis->A,op,&localmat);CHKERRQ(ierr); 375ad6194a2SStefano Zampini ierr = MatISSetLocalMat(B,localmat);CHKERRQ(ierr); 376ad6194a2SStefano Zampini ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 377ad6194a2SStefano Zampini ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 378ad6194a2SStefano Zampini *newmat = B; 379ad6194a2SStefano Zampini PetscFunctionReturn(0); 380ad6194a2SStefano Zampini } 381ad6194a2SStefano Zampini 382ad6194a2SStefano Zampini #undef __FUNCT__ 38369796d55SStefano Zampini #define __FUNCT__ "MatIsHermitian_IS" 38469796d55SStefano Zampini PetscErrorCode MatIsHermitian_IS(Mat A,PetscReal tol,PetscBool *flg) 38569796d55SStefano Zampini { 38669796d55SStefano Zampini PetscErrorCode ierr; 38769796d55SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 38869796d55SStefano Zampini PetscBool local_sym; 38969796d55SStefano Zampini 39069796d55SStefano Zampini PetscFunctionBegin; 39169796d55SStefano Zampini ierr = MatIsHermitian(matis->A,tol,&local_sym);CHKERRQ(ierr); 39269796d55SStefano Zampini ierr = MPI_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 39369796d55SStefano Zampini PetscFunctionReturn(0); 39469796d55SStefano Zampini } 39569796d55SStefano Zampini 39669796d55SStefano Zampini #undef __FUNCT__ 39769796d55SStefano Zampini #define __FUNCT__ "MatIsSymmetric_IS" 39869796d55SStefano Zampini PetscErrorCode MatIsSymmetric_IS(Mat A,PetscReal tol,PetscBool *flg) 39969796d55SStefano Zampini { 40069796d55SStefano Zampini PetscErrorCode ierr; 40169796d55SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 40269796d55SStefano Zampini PetscBool local_sym; 40369796d55SStefano Zampini 40469796d55SStefano Zampini PetscFunctionBegin; 40569796d55SStefano Zampini ierr = MatIsSymmetric(matis->A,tol,&local_sym);CHKERRQ(ierr); 40669796d55SStefano Zampini ierr = MPI_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 40769796d55SStefano Zampini PetscFunctionReturn(0); 40869796d55SStefano Zampini } 40969796d55SStefano Zampini 41069796d55SStefano Zampini #undef __FUNCT__ 411b4319ba4SBarry Smith #define __FUNCT__ "MatDestroy_IS" 412dfbe8321SBarry Smith PetscErrorCode MatDestroy_IS(Mat A) 413b4319ba4SBarry Smith { 414dfbe8321SBarry Smith PetscErrorCode ierr; 415b4319ba4SBarry Smith Mat_IS *b = (Mat_IS*)A->data; 416b4319ba4SBarry Smith 417b4319ba4SBarry Smith PetscFunctionBegin; 4186bf464f9SBarry Smith ierr = MatDestroy(&b->A);CHKERRQ(ierr); 4196bf464f9SBarry Smith ierr = VecScatterDestroy(&b->ctx);CHKERRQ(ierr); 4206bf464f9SBarry Smith ierr = VecDestroy(&b->x);CHKERRQ(ierr); 4216bf464f9SBarry Smith ierr = VecDestroy(&b->y);CHKERRQ(ierr); 4226bf464f9SBarry Smith ierr = ISLocalToGlobalMappingDestroy(&b->mapping);CHKERRQ(ierr); 423bf0cc555SLisandro Dalcin ierr = PetscFree(A->data);CHKERRQ(ierr); 424dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr); 425bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",NULL);CHKERRQ(ierr); 426b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",NULL);CHKERRQ(ierr); 427b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",NULL);CHKERRQ(ierr); 428*2e1947a5SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",NULL);CHKERRQ(ierr); 429b4319ba4SBarry Smith PetscFunctionReturn(0); 430b4319ba4SBarry Smith } 431b4319ba4SBarry Smith 432b4319ba4SBarry Smith #undef __FUNCT__ 433b4319ba4SBarry Smith #define __FUNCT__ "MatMult_IS" 434dfbe8321SBarry Smith PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y) 435b4319ba4SBarry Smith { 436dfbe8321SBarry Smith PetscErrorCode ierr; 437b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 438b4319ba4SBarry Smith PetscScalar zero = 0.0; 439b4319ba4SBarry Smith 440b4319ba4SBarry Smith PetscFunctionBegin; 441b4319ba4SBarry Smith /* scatter the global vector x into the local work vector */ 442ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 443ca9f406cSSatish Balay ierr = VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 444b4319ba4SBarry Smith 445b4319ba4SBarry Smith /* multiply the local matrix */ 446b4319ba4SBarry Smith ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr); 447b4319ba4SBarry Smith 448b4319ba4SBarry Smith /* scatter product back into global memory */ 4492dcb1b2aSMatthew Knepley ierr = VecSet(y,zero);CHKERRQ(ierr); 450ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 451ca9f406cSSatish Balay ierr = VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 452b4319ba4SBarry Smith PetscFunctionReturn(0); 453b4319ba4SBarry Smith } 454b4319ba4SBarry Smith 455b4319ba4SBarry Smith #undef __FUNCT__ 4562e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultAdd_IS" 457b7ce53b6SStefano Zampini PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 4582e74eeadSLisandro Dalcin { 459650997f4SStefano Zampini Vec temp_vec; 4602e74eeadSLisandro Dalcin PetscErrorCode ierr; 4612e74eeadSLisandro Dalcin 4622e74eeadSLisandro Dalcin PetscFunctionBegin; /* v3 = v2 + A * v1.*/ 463650997f4SStefano Zampini if (v3 != v2) { 464650997f4SStefano Zampini ierr = MatMult(A,v1,v3);CHKERRQ(ierr); 465650997f4SStefano Zampini ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 466650997f4SStefano Zampini } else { 467650997f4SStefano Zampini ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 468650997f4SStefano Zampini ierr = MatMult(A,v1,temp_vec);CHKERRQ(ierr); 469650997f4SStefano Zampini ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 470650997f4SStefano Zampini ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 471650997f4SStefano Zampini ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 472650997f4SStefano Zampini } 4732e74eeadSLisandro Dalcin PetscFunctionReturn(0); 4742e74eeadSLisandro Dalcin } 4752e74eeadSLisandro Dalcin 4762e74eeadSLisandro Dalcin #undef __FUNCT__ 4772e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTranspose_IS" 4782e74eeadSLisandro Dalcin PetscErrorCode MatMultTranspose_IS(Mat A,Vec x,Vec y) 4792e74eeadSLisandro Dalcin { 4802e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 4812e74eeadSLisandro Dalcin PetscErrorCode ierr; 4822e74eeadSLisandro Dalcin 4832e74eeadSLisandro Dalcin PetscFunctionBegin; /* y = A' * x */ 4842e74eeadSLisandro Dalcin /* scatter the global vector x into the local work vector */ 485ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 486ca9f406cSSatish Balay ierr = VecScatterEnd(is->ctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 4872e74eeadSLisandro Dalcin 4882e74eeadSLisandro Dalcin /* multiply the local matrix */ 4892e74eeadSLisandro Dalcin ierr = MatMultTranspose(is->A,is->x,is->y);CHKERRQ(ierr); 4902e74eeadSLisandro Dalcin 4912e74eeadSLisandro Dalcin /* scatter product back into global vector */ 4922e74eeadSLisandro Dalcin ierr = VecSet(y,0);CHKERRQ(ierr); 493ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 494ca9f406cSSatish Balay ierr = VecScatterEnd(is->ctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 4952e74eeadSLisandro Dalcin PetscFunctionReturn(0); 4962e74eeadSLisandro Dalcin } 4972e74eeadSLisandro Dalcin 4982e74eeadSLisandro Dalcin #undef __FUNCT__ 4992e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTransposeAdd_IS" 5002e74eeadSLisandro Dalcin PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 5012e74eeadSLisandro Dalcin { 502650997f4SStefano Zampini Vec temp_vec; 5032e74eeadSLisandro Dalcin PetscErrorCode ierr; 5042e74eeadSLisandro Dalcin 5052e74eeadSLisandro Dalcin PetscFunctionBegin; /* v3 = v2 + A' * v1.*/ 506650997f4SStefano Zampini if (v3 != v2) { 507650997f4SStefano Zampini ierr = MatMultTranspose(A,v1,v3);CHKERRQ(ierr); 508650997f4SStefano Zampini ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 509650997f4SStefano Zampini } else { 510650997f4SStefano Zampini ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 511650997f4SStefano Zampini ierr = MatMultTranspose(A,v1,temp_vec);CHKERRQ(ierr); 512650997f4SStefano Zampini ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 513650997f4SStefano Zampini ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 514650997f4SStefano Zampini ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 515650997f4SStefano Zampini } 5162e74eeadSLisandro Dalcin PetscFunctionReturn(0); 5172e74eeadSLisandro Dalcin } 5182e74eeadSLisandro Dalcin 5192e74eeadSLisandro Dalcin #undef __FUNCT__ 520b4319ba4SBarry Smith #define __FUNCT__ "MatView_IS" 521dfbe8321SBarry Smith PetscErrorCode MatView_IS(Mat A,PetscViewer viewer) 522b4319ba4SBarry Smith { 523b4319ba4SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 524dfbe8321SBarry Smith PetscErrorCode ierr; 525b4319ba4SBarry Smith PetscViewer sviewer; 526b4319ba4SBarry Smith 527b4319ba4SBarry Smith PetscFunctionBegin; 528dd2fa690SBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr); 529b4319ba4SBarry Smith ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 530b4319ba4SBarry Smith ierr = MatView(a->A,sviewer);CHKERRQ(ierr); 531b4319ba4SBarry Smith ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 532b4319ba4SBarry Smith PetscFunctionReturn(0); 533b4319ba4SBarry Smith } 534b4319ba4SBarry Smith 535b4319ba4SBarry Smith #undef __FUNCT__ 536b4319ba4SBarry Smith #define __FUNCT__ "MatSetLocalToGlobalMapping_IS" 537784ac674SJed Brown PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping) 538b4319ba4SBarry Smith { 539dfbe8321SBarry Smith PetscErrorCode ierr; 5404e4c7dbeSStefano Zampini PetscInt n,bs; 541b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 542b4319ba4SBarry Smith IS from,to; 543b4319ba4SBarry Smith Vec global; 544b4319ba4SBarry Smith 545b4319ba4SBarry Smith PetscFunctionBegin; 546784ac674SJed Brown PetscCheckSameComm(A,1,rmapping,2); 547ce94432eSBarry Smith if (rmapping != cmapping) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"MATIS requires the row and column mappings to be identical"); 54870cf5478SStefano Zampini if (is->mapping) { /* Currenly destroys the objects that will be created by this routine. Is there anything else that should be checked? */ 54970cf5478SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&is->mapping);CHKERRQ(ierr); 55070cf5478SStefano Zampini ierr = VecDestroy(&is->x);CHKERRQ(ierr); 55170cf5478SStefano Zampini ierr = VecDestroy(&is->y);CHKERRQ(ierr); 55270cf5478SStefano Zampini ierr = VecScatterDestroy(&is->ctx);CHKERRQ(ierr); 5531c47cb0fSStefano Zampini ierr = MatDestroy(&is->A);CHKERRQ(ierr); 55470cf5478SStefano Zampini } 555784ac674SJed Brown ierr = PetscObjectReference((PetscObject)rmapping);CHKERRQ(ierr); 5566bf464f9SBarry Smith ierr = ISLocalToGlobalMappingDestroy(&is->mapping);CHKERRQ(ierr); 557784ac674SJed Brown is->mapping = rmapping; 558fa7f1dd8SStefano Zampini /* 559fa7f1dd8SStefano Zampini ierr = PetscLayoutSetISLocalToGlobalMapping(A->rmap,rmapping);CHKERRQ(ierr); 560fa7f1dd8SStefano Zampini ierr = PetscLayoutSetISLocalToGlobalMapping(A->cmap,cmapping);CHKERRQ(ierr); 561fa7f1dd8SStefano Zampini */ 562b4319ba4SBarry Smith 563b4319ba4SBarry Smith /* Create the local matrix A */ 564784ac674SJed Brown ierr = ISLocalToGlobalMappingGetSize(rmapping,&n);CHKERRQ(ierr); 565*2e1947a5SStefano Zampini ierr = ISLocalToGlobalMappingGetBlockSize(rmapping,&bs);CHKERRQ(ierr); 566f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr); 567*2e1947a5SStefano Zampini if (bs > 1) { 568*2e1947a5SStefano Zampini ierr = MatSetType(is->A,MATSEQBAIJ);CHKERRQ(ierr); 569*2e1947a5SStefano Zampini } else { 570*2e1947a5SStefano Zampini ierr = MatSetType(is->A,MATSEQAIJ);CHKERRQ(ierr); 571*2e1947a5SStefano Zampini } 572f69a0ea3SMatthew Knepley ierr = MatSetSizes(is->A,n,n,n,n);CHKERRQ(ierr); 5734e4c7dbeSStefano Zampini ierr = MatSetBlockSize(is->A,bs);CHKERRQ(ierr); 574ff130e51SJed Brown ierr = MatSetOptionsPrefix(is->A,((PetscObject)A)->prefix);CHKERRQ(ierr); 575ff130e51SJed Brown ierr = MatAppendOptionsPrefix(is->A,"is_");CHKERRQ(ierr); 576b4319ba4SBarry Smith ierr = MatSetFromOptions(is->A);CHKERRQ(ierr); 577b4319ba4SBarry Smith 578b4319ba4SBarry Smith /* Create the local work vectors */ 5794e4c7dbeSStefano Zampini ierr = VecCreate(PETSC_COMM_SELF,&is->x);CHKERRQ(ierr); 5804e4c7dbeSStefano Zampini ierr = VecSetBlockSize(is->x,bs);CHKERRQ(ierr); 5814e4c7dbeSStefano Zampini ierr = VecSetSizes(is->x,n,n);CHKERRQ(ierr); 582ff130e51SJed Brown ierr = VecSetOptionsPrefix(is->x,((PetscObject)A)->prefix);CHKERRQ(ierr); 583ff130e51SJed Brown ierr = VecAppendOptionsPrefix(is->x,"is_");CHKERRQ(ierr); 5844e4c7dbeSStefano Zampini ierr = VecSetFromOptions(is->x);CHKERRQ(ierr); 585b4319ba4SBarry Smith ierr = VecDuplicate(is->x,&is->y);CHKERRQ(ierr); 586b4319ba4SBarry Smith 587b4319ba4SBarry Smith /* setup the global to local scatter */ 588b4319ba4SBarry Smith ierr = ISCreateStride(PETSC_COMM_SELF,n,0,1,&to);CHKERRQ(ierr); 589784ac674SJed Brown ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr); 5900298fd71SBarry Smith ierr = MatGetVecs(A,&global,NULL);CHKERRQ(ierr); 591b4319ba4SBarry Smith ierr = VecScatterCreate(global,from,is->x,to,&is->ctx);CHKERRQ(ierr); 5926bf464f9SBarry Smith ierr = VecDestroy(&global);CHKERRQ(ierr); 5936bf464f9SBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 5946bf464f9SBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 595b4319ba4SBarry Smith PetscFunctionReturn(0); 596b4319ba4SBarry Smith } 597b4319ba4SBarry Smith 5982e74eeadSLisandro Dalcin #undef __FUNCT__ 5992e74eeadSLisandro Dalcin #define __FUNCT__ "MatSetValues_IS" 6002e74eeadSLisandro Dalcin PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv) 6012e74eeadSLisandro Dalcin { 6022e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)mat->data; 6032e74eeadSLisandro Dalcin PetscInt rows_l[2048],cols_l[2048]; 6042e74eeadSLisandro Dalcin PetscErrorCode ierr; 6052e74eeadSLisandro Dalcin 6062e74eeadSLisandro Dalcin PetscFunctionBegin; 6072e74eeadSLisandro Dalcin #if defined(PETSC_USE_DEBUG) 608f23aa3ddSBarry 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); 6092e74eeadSLisandro Dalcin #endif 6102e74eeadSLisandro Dalcin ierr = ISG2LMapApply(is->mapping,m,rows,rows_l);CHKERRQ(ierr); 6112e74eeadSLisandro Dalcin ierr = ISG2LMapApply(is->mapping,n,cols,cols_l);CHKERRQ(ierr); 6122e74eeadSLisandro Dalcin ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 6132e74eeadSLisandro Dalcin PetscFunctionReturn(0); 6142e74eeadSLisandro Dalcin } 6152e74eeadSLisandro Dalcin 6162e74eeadSLisandro Dalcin #undef ISG2LMapSetUp 6172e74eeadSLisandro Dalcin #undef ISG2LMapApply 618b4319ba4SBarry Smith 619b4319ba4SBarry Smith #undef __FUNCT__ 620b4319ba4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_IS" 62113f74950SBarry Smith PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 622b4319ba4SBarry Smith { 623dfbe8321SBarry Smith PetscErrorCode ierr; 624b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 625b4319ba4SBarry Smith 626b4319ba4SBarry Smith PetscFunctionBegin; 627b4319ba4SBarry Smith ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 628b4319ba4SBarry Smith PetscFunctionReturn(0); 629b4319ba4SBarry Smith } 630b4319ba4SBarry Smith 631b4319ba4SBarry Smith #undef __FUNCT__ 632f0006bf2SLisandro Dalcin #define __FUNCT__ "MatSetValuesBlockedLocal_IS" 633f0006bf2SLisandro Dalcin PetscErrorCode MatSetValuesBlockedLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 634f0006bf2SLisandro Dalcin { 635f0006bf2SLisandro Dalcin PetscErrorCode ierr; 636f0006bf2SLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 637f0006bf2SLisandro Dalcin 638f0006bf2SLisandro Dalcin PetscFunctionBegin; 639f0006bf2SLisandro Dalcin ierr = MatSetValuesBlocked(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 640f0006bf2SLisandro Dalcin PetscFunctionReturn(0); 641f0006bf2SLisandro Dalcin } 642f0006bf2SLisandro Dalcin 643f0006bf2SLisandro Dalcin #undef __FUNCT__ 6442e74eeadSLisandro Dalcin #define __FUNCT__ "MatZeroRows_IS" 6452b40b63fSBarry Smith PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 6462e74eeadSLisandro Dalcin { 6472e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 6480298fd71SBarry Smith PetscInt n_l = 0, *rows_l = NULL; 6492e74eeadSLisandro Dalcin PetscErrorCode ierr; 6502e74eeadSLisandro Dalcin 6512e74eeadSLisandro Dalcin PetscFunctionBegin; 652ce94432eSBarry Smith if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support"); 6532e74eeadSLisandro Dalcin if (n) { 654785e854fSJed Brown ierr = PetscMalloc1(n,&rows_l);CHKERRQ(ierr); 6552e74eeadSLisandro Dalcin ierr = ISGlobalToLocalMappingApply(is->mapping,IS_GTOLM_DROP,n,rows,&n_l,rows_l);CHKERRQ(ierr); 6562e74eeadSLisandro Dalcin } 6572b40b63fSBarry Smith ierr = MatZeroRowsLocal(A,n_l,rows_l,diag,x,b);CHKERRQ(ierr); 658c31cb41cSBarry Smith ierr = PetscFree(rows_l);CHKERRQ(ierr); 6592e74eeadSLisandro Dalcin PetscFunctionReturn(0); 6602e74eeadSLisandro Dalcin } 6612e74eeadSLisandro Dalcin 6622e74eeadSLisandro Dalcin #undef __FUNCT__ 663b4319ba4SBarry Smith #define __FUNCT__ "MatZeroRowsLocal_IS" 6642b40b63fSBarry Smith PetscErrorCode MatZeroRowsLocal_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 665b4319ba4SBarry Smith { 666b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 667dfbe8321SBarry Smith PetscErrorCode ierr; 668f4df32b1SMatthew Knepley PetscInt i; 669b4319ba4SBarry Smith PetscScalar *array; 670b4319ba4SBarry Smith 671b4319ba4SBarry Smith PetscFunctionBegin; 672ce94432eSBarry Smith if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support"); 673b4319ba4SBarry Smith { 674b4319ba4SBarry Smith /* 675b4319ba4SBarry Smith Set up is->x as a "counting vector". This is in order to MatMult_IS 676b4319ba4SBarry Smith work properly in the interface nodes. 677b4319ba4SBarry Smith */ 678b4319ba4SBarry Smith Vec counter; 679b4319ba4SBarry Smith PetscScalar one=1.0, zero=0.0; 6800298fd71SBarry Smith ierr = MatGetVecs(A,&counter,NULL);CHKERRQ(ierr); 6812dcb1b2aSMatthew Knepley ierr = VecSet(counter,zero);CHKERRQ(ierr); 6822dcb1b2aSMatthew Knepley ierr = VecSet(is->x,one);CHKERRQ(ierr); 683ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 684ca9f406cSSatish Balay ierr = VecScatterEnd (is->ctx,is->x,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 685ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 686ca9f406cSSatish Balay ierr = VecScatterEnd (is->ctx,counter,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 6876bf464f9SBarry Smith ierr = VecDestroy(&counter);CHKERRQ(ierr); 688b4319ba4SBarry Smith } 689958c9bccSBarry Smith if (!n) { 690b4319ba4SBarry Smith is->pure_neumann = PETSC_TRUE; 691b4319ba4SBarry Smith } else { 692b4319ba4SBarry Smith is->pure_neumann = PETSC_FALSE; 6932205254eSKarl Rupp 694b4319ba4SBarry Smith ierr = VecGetArray(is->x,&array);CHKERRQ(ierr); 6952b40b63fSBarry Smith ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr); 696b4319ba4SBarry Smith for (i=0; i<n; i++) { 697f4df32b1SMatthew Knepley ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr); 698b4319ba4SBarry Smith } 699b4319ba4SBarry Smith ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 700b4319ba4SBarry Smith ierr = MatAssemblyEnd (is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 701b4319ba4SBarry Smith ierr = VecRestoreArray(is->x,&array);CHKERRQ(ierr); 702b4319ba4SBarry Smith } 703b4319ba4SBarry Smith PetscFunctionReturn(0); 704b4319ba4SBarry Smith } 705b4319ba4SBarry Smith 706b4319ba4SBarry Smith #undef __FUNCT__ 707b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyBegin_IS" 708dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type) 709b4319ba4SBarry Smith { 710b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 711dfbe8321SBarry Smith PetscErrorCode ierr; 712b4319ba4SBarry Smith 713b4319ba4SBarry Smith PetscFunctionBegin; 714b4319ba4SBarry Smith ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr); 715b4319ba4SBarry Smith PetscFunctionReturn(0); 716b4319ba4SBarry Smith } 717b4319ba4SBarry Smith 718b4319ba4SBarry Smith #undef __FUNCT__ 719b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_IS" 720dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type) 721b4319ba4SBarry Smith { 722b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 723dfbe8321SBarry Smith PetscErrorCode ierr; 724b4319ba4SBarry Smith 725b4319ba4SBarry Smith PetscFunctionBegin; 726b4319ba4SBarry Smith ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr); 727b4319ba4SBarry Smith PetscFunctionReturn(0); 728b4319ba4SBarry Smith } 729b4319ba4SBarry Smith 730b4319ba4SBarry Smith #undef __FUNCT__ 731b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat_IS" 7327087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local) 733b4319ba4SBarry Smith { 734b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)mat->data; 735b4319ba4SBarry Smith 736b4319ba4SBarry Smith PetscFunctionBegin; 737b4319ba4SBarry Smith *local = is->A; 738b4319ba4SBarry Smith PetscFunctionReturn(0); 739b4319ba4SBarry Smith } 740b4319ba4SBarry Smith 741b4319ba4SBarry Smith #undef __FUNCT__ 742b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat" 743b4319ba4SBarry Smith /*@ 744b4319ba4SBarry Smith MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix. 745b4319ba4SBarry Smith 746b4319ba4SBarry Smith Input Parameter: 747b4319ba4SBarry Smith . mat - the matrix 748b4319ba4SBarry Smith 749b4319ba4SBarry Smith Output Parameter: 750b4319ba4SBarry Smith . local - the local matrix usually MATSEQAIJ 751b4319ba4SBarry Smith 752b4319ba4SBarry Smith Level: advanced 753b4319ba4SBarry Smith 754b4319ba4SBarry Smith Notes: 755b4319ba4SBarry Smith This can be called if you have precomputed the nonzero structure of the 756b4319ba4SBarry Smith matrix and want to provide it to the inner matrix object to improve the performance 757b4319ba4SBarry Smith of the MatSetValues() operation. 758b4319ba4SBarry Smith 759b4319ba4SBarry Smith .seealso: MATIS 760b4319ba4SBarry Smith @*/ 7617087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local) 762b4319ba4SBarry Smith { 7634ac538c5SBarry Smith PetscErrorCode ierr; 764b4319ba4SBarry Smith 765b4319ba4SBarry Smith PetscFunctionBegin; 7660700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 767b4319ba4SBarry Smith PetscValidPointer(local,2); 7684ac538c5SBarry Smith ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat*),(mat,local));CHKERRQ(ierr); 769b4319ba4SBarry Smith PetscFunctionReturn(0); 770b4319ba4SBarry Smith } 771b4319ba4SBarry Smith 7723b03a366Sstefano_zampini #undef __FUNCT__ 7733b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat_IS" 7743b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat_IS(Mat mat,Mat local) 7753b03a366Sstefano_zampini { 7763b03a366Sstefano_zampini Mat_IS *is = (Mat_IS*)mat->data; 7773b03a366Sstefano_zampini PetscInt nrows,ncols,orows,ocols; 7783b03a366Sstefano_zampini PetscErrorCode ierr; 7793b03a366Sstefano_zampini 7803b03a366Sstefano_zampini PetscFunctionBegin; 7814e4c7dbeSStefano Zampini if (is->A) { 7823b03a366Sstefano_zampini ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr); 7833b03a366Sstefano_zampini ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr); 784f23aa3ddSBarry 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); 7854e4c7dbeSStefano Zampini } 7863b03a366Sstefano_zampini ierr = PetscObjectReference((PetscObject)local);CHKERRQ(ierr); 7873b03a366Sstefano_zampini ierr = MatDestroy(&is->A);CHKERRQ(ierr); 7883b03a366Sstefano_zampini is->A = local; 7893b03a366Sstefano_zampini PetscFunctionReturn(0); 7903b03a366Sstefano_zampini } 7913b03a366Sstefano_zampini 7923b03a366Sstefano_zampini #undef __FUNCT__ 7933b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat" 7943b03a366Sstefano_zampini /*@ 7953b03a366Sstefano_zampini MatISSetLocalMat - Set the local matrix stored inside a MATIS matrix. 7963b03a366Sstefano_zampini 7973b03a366Sstefano_zampini Input Parameter: 7983b03a366Sstefano_zampini . mat - the matrix 7993b03a366Sstefano_zampini . local - the local matrix usually MATSEQAIJ 8003b03a366Sstefano_zampini 8013b03a366Sstefano_zampini Output Parameter: 8023b03a366Sstefano_zampini 8033b03a366Sstefano_zampini Level: advanced 8043b03a366Sstefano_zampini 8053b03a366Sstefano_zampini Notes: 8063b03a366Sstefano_zampini This can be called if you have precomputed the local matrix and 8073b03a366Sstefano_zampini want to provide it to the matrix object MATIS. 8083b03a366Sstefano_zampini 8093b03a366Sstefano_zampini .seealso: MATIS 8103b03a366Sstefano_zampini @*/ 8113b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat(Mat mat,Mat local) 8123b03a366Sstefano_zampini { 8133b03a366Sstefano_zampini PetscErrorCode ierr; 8143b03a366Sstefano_zampini 8153b03a366Sstefano_zampini PetscFunctionBegin; 8163b03a366Sstefano_zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 817b7ce53b6SStefano Zampini PetscValidHeaderSpecific(local,MAT_CLASSID,2); 8183b03a366Sstefano_zampini ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr); 8193b03a366Sstefano_zampini PetscFunctionReturn(0); 8203b03a366Sstefano_zampini } 8213b03a366Sstefano_zampini 8226726f965SBarry Smith #undef __FUNCT__ 8236726f965SBarry Smith #define __FUNCT__ "MatZeroEntries_IS" 8246726f965SBarry Smith PetscErrorCode MatZeroEntries_IS(Mat A) 8256726f965SBarry Smith { 8266726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 8276726f965SBarry Smith PetscErrorCode ierr; 8286726f965SBarry Smith 8296726f965SBarry Smith PetscFunctionBegin; 8306726f965SBarry Smith ierr = MatZeroEntries(a->A);CHKERRQ(ierr); 8316726f965SBarry Smith PetscFunctionReturn(0); 8326726f965SBarry Smith } 8336726f965SBarry Smith 8346726f965SBarry Smith #undef __FUNCT__ 8352e74eeadSLisandro Dalcin #define __FUNCT__ "MatScale_IS" 8362e74eeadSLisandro Dalcin PetscErrorCode MatScale_IS(Mat A,PetscScalar a) 8372e74eeadSLisandro Dalcin { 8382e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 8392e74eeadSLisandro Dalcin PetscErrorCode ierr; 8402e74eeadSLisandro Dalcin 8412e74eeadSLisandro Dalcin PetscFunctionBegin; 8422e74eeadSLisandro Dalcin ierr = MatScale(is->A,a);CHKERRQ(ierr); 8432e74eeadSLisandro Dalcin PetscFunctionReturn(0); 8442e74eeadSLisandro Dalcin } 8452e74eeadSLisandro Dalcin 8462e74eeadSLisandro Dalcin #undef __FUNCT__ 8472e74eeadSLisandro Dalcin #define __FUNCT__ "MatGetDiagonal_IS" 8482e74eeadSLisandro Dalcin PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v) 8492e74eeadSLisandro Dalcin { 8502e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 8512e74eeadSLisandro Dalcin PetscErrorCode ierr; 8522e74eeadSLisandro Dalcin 8532e74eeadSLisandro Dalcin PetscFunctionBegin; 8542e74eeadSLisandro Dalcin /* get diagonal of the local matrix */ 8552e74eeadSLisandro Dalcin ierr = MatGetDiagonal(is->A,is->x);CHKERRQ(ierr); 8562e74eeadSLisandro Dalcin 8572e74eeadSLisandro Dalcin /* scatter diagonal back into global vector */ 8582e74eeadSLisandro Dalcin ierr = VecSet(v,0);CHKERRQ(ierr); 859ca9f406cSSatish Balay ierr = VecScatterBegin(is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 860ca9f406cSSatish Balay ierr = VecScatterEnd (is->ctx,is->x,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 8612e74eeadSLisandro Dalcin PetscFunctionReturn(0); 8622e74eeadSLisandro Dalcin } 8632e74eeadSLisandro Dalcin 8642e74eeadSLisandro Dalcin #undef __FUNCT__ 8656726f965SBarry Smith #define __FUNCT__ "MatSetOption_IS" 866ace3abfcSBarry Smith PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg) 8676726f965SBarry Smith { 8686726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 8696726f965SBarry Smith PetscErrorCode ierr; 8706726f965SBarry Smith 8716726f965SBarry Smith PetscFunctionBegin; 8724e0d8c25SBarry Smith ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr); 8736726f965SBarry Smith PetscFunctionReturn(0); 8746726f965SBarry Smith } 8756726f965SBarry Smith 876284134d9SBarry Smith #undef __FUNCT__ 877284134d9SBarry Smith #define __FUNCT__ "MatCreateIS" 878284134d9SBarry Smith /*@ 879284134d9SBarry Smith MatCreateIS - Creates a "process" unassmembled matrix, it is assembled on each 880284134d9SBarry Smith process but not across processes. 881284134d9SBarry Smith 882284134d9SBarry Smith Input Parameters: 883284134d9SBarry Smith + comm - MPI communicator that will share the matrix 8844e4c7dbeSStefano Zampini . bs - local and global block size of the matrix 885df3898eeSBarry Smith . m,n,M,N - local and/or global sizes of the left and right vector used in matrix vector products 886284134d9SBarry Smith - map - mapping that defines the global number for each local number 887284134d9SBarry Smith 888284134d9SBarry Smith Output Parameter: 889284134d9SBarry Smith . A - the resulting matrix 890284134d9SBarry Smith 8918e6c10adSSatish Balay Level: advanced 8928e6c10adSSatish Balay 893284134d9SBarry Smith Notes: See MATIS for more details 8948cda6cd7SBarry Smith m and n are NOT related to the size of the map, they are the size of the part of the vector owned 8958cda6cd7SBarry Smith by that process. m + nghosts (or n + nghosts) is the length of map since map maps all local points 8968cda6cd7SBarry Smith plus the ghost points to global indices. 897284134d9SBarry Smith 898284134d9SBarry Smith .seealso: MATIS, MatSetLocalToGlobalMapping() 899284134d9SBarry Smith @*/ 9004e4c7dbeSStefano Zampini PetscErrorCode MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping map,Mat *A) 901284134d9SBarry Smith { 902284134d9SBarry Smith PetscErrorCode ierr; 903284134d9SBarry Smith 904284134d9SBarry Smith PetscFunctionBegin; 905284134d9SBarry Smith ierr = MatCreate(comm,A);CHKERRQ(ierr); 906d3ee2243SStefano Zampini ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr); 907284134d9SBarry Smith ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 908284134d9SBarry Smith ierr = MatSetType(*A,MATIS);CHKERRQ(ierr); 9093b03a366Sstefano_zampini ierr = MatSetUp(*A);CHKERRQ(ierr); 910784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(*A,map,map);CHKERRQ(ierr); 911284134d9SBarry Smith PetscFunctionReturn(0); 912284134d9SBarry Smith } 913284134d9SBarry Smith 914b4319ba4SBarry Smith /*MC 9156a818285SBarry Smith MATIS - MATIS = "is" - A matrix type to be used for using the Neumann-Neumann type preconditioners, see PCBDDC. 916b4319ba4SBarry Smith This stores the matrices in globally unassembled form. Each processor 917b4319ba4SBarry Smith assembles only its local Neumann problem and the parallel matrix vector 918b4319ba4SBarry Smith product is handled "implicitly". 919b4319ba4SBarry Smith 920b4319ba4SBarry Smith Operations Provided: 9216726f965SBarry Smith + MatMult() 9222e74eeadSLisandro Dalcin . MatMultAdd() 9232e74eeadSLisandro Dalcin . MatMultTranspose() 9242e74eeadSLisandro Dalcin . MatMultTransposeAdd() 9256726f965SBarry Smith . MatZeroEntries() 9266726f965SBarry Smith . MatSetOption() 9272e74eeadSLisandro Dalcin . MatZeroRows() 9286726f965SBarry Smith . MatZeroRowsLocal() 9292e74eeadSLisandro Dalcin . MatSetValues() 9306726f965SBarry Smith . MatSetValuesLocal() 9312e74eeadSLisandro Dalcin . MatScale() 9322e74eeadSLisandro Dalcin . MatGetDiagonal() 9336726f965SBarry Smith - MatSetLocalToGlobalMapping() 934b4319ba4SBarry Smith 935b4319ba4SBarry Smith Options Database Keys: 936b4319ba4SBarry Smith . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions() 937b4319ba4SBarry Smith 938b4319ba4SBarry Smith Notes: Options prefix for the inner matrix are given by -is_mat_xxx 939b4319ba4SBarry Smith 940b4319ba4SBarry Smith You must call MatSetLocalToGlobalMapping() before using this matrix type. 941b4319ba4SBarry Smith 942b4319ba4SBarry Smith You can do matrix preallocation on the local matrix after you obtain it with 943b4319ba4SBarry Smith MatISGetLocalMat() 944b4319ba4SBarry Smith 945b4319ba4SBarry Smith Level: advanced 946b4319ba4SBarry Smith 9476a818285SBarry Smith .seealso: PC, MatISGetLocalMat(), MatSetLocalToGlobalMapping(), PCBDDC 948b4319ba4SBarry Smith 949b4319ba4SBarry Smith M*/ 950b4319ba4SBarry Smith 951b4319ba4SBarry Smith #undef __FUNCT__ 952b4319ba4SBarry Smith #define __FUNCT__ "MatCreate_IS" 9538cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_IS(Mat A) 954b4319ba4SBarry Smith { 955dfbe8321SBarry Smith PetscErrorCode ierr; 956b4319ba4SBarry Smith Mat_IS *b; 957b4319ba4SBarry Smith 958b4319ba4SBarry Smith PetscFunctionBegin; 959b00a9115SJed Brown ierr = PetscNewLog(A,&b);CHKERRQ(ierr); 960b4319ba4SBarry Smith A->data = (void*)b; 961b4319ba4SBarry Smith ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 962b4319ba4SBarry Smith 963b4319ba4SBarry Smith A->ops->mult = MatMult_IS; 9642e74eeadSLisandro Dalcin A->ops->multadd = MatMultAdd_IS; 9652e74eeadSLisandro Dalcin A->ops->multtranspose = MatMultTranspose_IS; 9662e74eeadSLisandro Dalcin A->ops->multtransposeadd = MatMultTransposeAdd_IS; 967b4319ba4SBarry Smith A->ops->destroy = MatDestroy_IS; 968b4319ba4SBarry Smith A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS; 9692e74eeadSLisandro Dalcin A->ops->setvalues = MatSetValues_IS; 970b4319ba4SBarry Smith A->ops->setvalueslocal = MatSetValuesLocal_IS; 971f0006bf2SLisandro Dalcin A->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_IS; 9722e74eeadSLisandro Dalcin A->ops->zerorows = MatZeroRows_IS; 973b4319ba4SBarry Smith A->ops->zerorowslocal = MatZeroRowsLocal_IS; 974b4319ba4SBarry Smith A->ops->assemblybegin = MatAssemblyBegin_IS; 975b4319ba4SBarry Smith A->ops->assemblyend = MatAssemblyEnd_IS; 976b4319ba4SBarry Smith A->ops->view = MatView_IS; 9776726f965SBarry Smith A->ops->zeroentries = MatZeroEntries_IS; 9782e74eeadSLisandro Dalcin A->ops->scale = MatScale_IS; 9792e74eeadSLisandro Dalcin A->ops->getdiagonal = MatGetDiagonal_IS; 9806726f965SBarry Smith A->ops->setoption = MatSetOption_IS; 98169796d55SStefano Zampini A->ops->ishermitian = MatIsHermitian_IS; 98269796d55SStefano Zampini A->ops->issymmetric = MatIsSymmetric_IS; 983ad6194a2SStefano Zampini A->ops->duplicate = MatDuplicate_IS; 984b4319ba4SBarry Smith 98526283091SBarry Smith ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 98626283091SBarry Smith ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 987b4319ba4SBarry Smith 988b7ce53b6SStefano Zampini /* zeros pointer */ 989b4319ba4SBarry Smith b->A = 0; 990b4319ba4SBarry Smith b->ctx = 0; 991b4319ba4SBarry Smith b->x = 0; 992b4319ba4SBarry Smith b->y = 0; 993b09366fdSStefano Zampini b->mapping = 0; 994b7ce53b6SStefano Zampini 995b7ce53b6SStefano Zampini /* special MATIS functions */ 996bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",MatISGetLocalMat_IS);CHKERRQ(ierr); 997bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",MatISSetLocalMat_IS);CHKERRQ(ierr); 998b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",MatISGetMPIXAIJ_IS);CHKERRQ(ierr); 999*2e1947a5SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",MatISSetPreallocation_IS);CHKERRQ(ierr); 100017667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr); 1001b4319ba4SBarry Smith PetscFunctionReturn(0); 1002b4319ba4SBarry Smith } 1003