1be1d678aSKris Buschelman 2b4319ba4SBarry Smith /* 3b4319ba4SBarry Smith Creates a matrix class for using the Neumann-Neumann type preconditioners. 4b4319ba4SBarry Smith This stores the matrices in globally unassembled form. Each processor 5b4319ba4SBarry Smith assembles only its local Neumann problem and the parallel matrix vector 6b4319ba4SBarry Smith product is handled "implicitly". 7b4319ba4SBarry Smith 8b4319ba4SBarry Smith We provide: 9b4319ba4SBarry Smith MatMult() 10b4319ba4SBarry Smith 11b4319ba4SBarry Smith Currently this allows for only one subdomain per processor. 12b4319ba4SBarry Smith 13b4319ba4SBarry Smith */ 14b4319ba4SBarry Smith 15c6db04a5SJed Brown #include <../src/mat/impls/is/matis.h> /*I "petscmat.h" I*/ 1628f4e0baSStefano Zampini 17a72627d2SStefano Zampini static PetscErrorCode MatISComputeSF_Private(Mat); 18a72627d2SStefano Zampini 1928f4e0baSStefano Zampini #undef __FUNCT__ 2028f4e0baSStefano Zampini #define __FUNCT__ "MatISComputeSF_Private" 21a72627d2SStefano Zampini static PetscErrorCode MatISComputeSF_Private(Mat B) 2228f4e0baSStefano Zampini { 2328f4e0baSStefano Zampini Mat_IS *matis = (Mat_IS*)(B->data); 2428f4e0baSStefano Zampini const PetscInt *gidxs; 2528f4e0baSStefano Zampini PetscErrorCode ierr; 2628f4e0baSStefano Zampini 2728f4e0baSStefano Zampini PetscFunctionBegin; 2828f4e0baSStefano Zampini ierr = MatGetSize(matis->A,&matis->sf_nleaves,NULL);CHKERRQ(ierr); 2928f4e0baSStefano Zampini ierr = MatGetLocalSize(B,&matis->sf_nroots,NULL);CHKERRQ(ierr); 3028f4e0baSStefano Zampini ierr = PetscSFCreate(PetscObjectComm((PetscObject)B),&matis->sf);CHKERRQ(ierr); 313bbff08aSStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(B->rmap->mapping,&gidxs);CHKERRQ(ierr); 3228f4e0baSStefano Zampini /* PETSC_OWN_POINTER refers to ilocal which is NULL */ 3328f4e0baSStefano Zampini ierr = PetscSFSetGraphLayout(matis->sf,B->rmap,matis->sf_nleaves,NULL,PETSC_OWN_POINTER,gidxs);CHKERRQ(ierr); 343bbff08aSStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(B->rmap->mapping,&gidxs);CHKERRQ(ierr); 3528f4e0baSStefano Zampini ierr = PetscMalloc2(matis->sf_nroots,&matis->sf_rootdata,matis->sf_nleaves,&matis->sf_leafdata);CHKERRQ(ierr); 3628f4e0baSStefano Zampini PetscFunctionReturn(0); 3728f4e0baSStefano Zampini } 382e1947a5SStefano Zampini 392e1947a5SStefano Zampini #undef __FUNCT__ 402e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation" 41eb82efa4SStefano Zampini /*@ 42a88811baSStefano Zampini MatISSetPreallocation - Preallocates memory for a MATIS parallel matrix. 43a88811baSStefano Zampini 44a88811baSStefano Zampini Collective on MPI_Comm 45a88811baSStefano Zampini 46a88811baSStefano Zampini Input Parameters: 47a88811baSStefano Zampini + B - the matrix 48a88811baSStefano Zampini . d_nz - number of nonzeros per row in DIAGONAL portion of local submatrix 49a88811baSStefano Zampini (same value is used for all local rows) 50a88811baSStefano Zampini . d_nnz - array containing the number of nonzeros in the various rows of the 51a88811baSStefano Zampini DIAGONAL portion of the local submatrix (possibly different for each row) 52a88811baSStefano Zampini or NULL, if d_nz is used to specify the nonzero structure. 53a88811baSStefano Zampini The size of this array is equal to the number of local rows, i.e 'm'. 54a88811baSStefano Zampini For matrices that will be factored, you must leave room for (and set) 55a88811baSStefano Zampini the diagonal entry even if it is zero. 56a88811baSStefano Zampini . o_nz - number of nonzeros per row in the OFF-DIAGONAL portion of local 57a88811baSStefano Zampini submatrix (same value is used for all local rows). 58a88811baSStefano Zampini - o_nnz - array containing the number of nonzeros in the various rows of the 59a88811baSStefano Zampini OFF-DIAGONAL portion of the local submatrix (possibly different for 60a88811baSStefano Zampini each row) or NULL, if o_nz is used to specify the nonzero 61a88811baSStefano Zampini structure. The size of this array is equal to the number 62a88811baSStefano Zampini of local rows, i.e 'm'. 63a88811baSStefano Zampini 64a88811baSStefano Zampini If the *_nnz parameter is given then the *_nz parameter is ignored 65a88811baSStefano Zampini 66a88811baSStefano Zampini Level: intermediate 67a88811baSStefano Zampini 68a88811baSStefano Zampini Notes: This function has the same interface as the MPIAIJ preallocation routine in order to simplify the transition 69a88811baSStefano Zampini from the asssembled format to the unassembled one. It overestimates the preallocation of MATIS local 70a88811baSStefano Zampini matrices; for exact preallocation, the user should set the preallocation directly on local matrix objects. 71a88811baSStefano Zampini 72a88811baSStefano Zampini .keywords: matrix 73a88811baSStefano Zampini 74a88811baSStefano Zampini .seealso: MatCreate(), MatCreateIS(), MatMPIAIJSetPreallocation(), MatISGetLocalMat() 75a88811baSStefano Zampini @*/ 762e1947a5SStefano Zampini PetscErrorCode MatISSetPreallocation(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[]) 772e1947a5SStefano Zampini { 782e1947a5SStefano Zampini PetscErrorCode ierr; 792e1947a5SStefano Zampini 802e1947a5SStefano Zampini PetscFunctionBegin; 812e1947a5SStefano Zampini PetscValidHeaderSpecific(B,MAT_CLASSID,1); 822e1947a5SStefano Zampini PetscValidType(B,1); 832e1947a5SStefano Zampini ierr = PetscTryMethod(B,"MatISSetPreallocation_C",(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(B,d_nz,d_nnz,o_nz,o_nnz));CHKERRQ(ierr); 842e1947a5SStefano Zampini PetscFunctionReturn(0); 852e1947a5SStefano Zampini } 862e1947a5SStefano Zampini 872e1947a5SStefano Zampini #undef __FUNCT__ 882e1947a5SStefano Zampini #define __FUNCT__ "MatISSetPreallocation_IS" 892e1947a5SStefano Zampini PetscErrorCode MatISSetPreallocation_IS(Mat B,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[]) 902e1947a5SStefano Zampini { 912e1947a5SStefano Zampini Mat_IS *matis = (Mat_IS*)(B->data); 9228f4e0baSStefano Zampini PetscInt bs,i,nlocalcols; 932e1947a5SStefano Zampini PetscErrorCode ierr; 942e1947a5SStefano Zampini 952e1947a5SStefano Zampini PetscFunctionBegin; 966c4ed002SBarry Smith if (!matis->A) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"You should first call MatSetLocalToGlobalMapping"); 9728f4e0baSStefano Zampini if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */ 9828f4e0baSStefano Zampini ierr = MatISComputeSF_Private(B);CHKERRQ(ierr); 9928f4e0baSStefano Zampini } 1002e1947a5SStefano Zampini if (!d_nnz) { 10128f4e0baSStefano Zampini for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] = d_nz; 1022e1947a5SStefano Zampini } else { 10328f4e0baSStefano Zampini for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] = d_nnz[i]; 1042e1947a5SStefano Zampini } 1052e1947a5SStefano Zampini if (!o_nnz) { 10628f4e0baSStefano Zampini for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] += o_nz; 1072e1947a5SStefano Zampini } else { 10828f4e0baSStefano Zampini for (i=0;i<matis->sf_nroots;i++) matis->sf_rootdata[i] += o_nnz[i]; 1092e1947a5SStefano Zampini } 11028f4e0baSStefano Zampini ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 11128f4e0baSStefano Zampini ierr = MatGetSize(matis->A,NULL,&nlocalcols);CHKERRQ(ierr); 11228f4e0baSStefano Zampini ierr = MatGetBlockSize(matis->A,&bs);CHKERRQ(ierr); 11328f4e0baSStefano Zampini ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 11428f4e0baSStefano Zampini for (i=0;i<matis->sf_nleaves;i++) { 11528f4e0baSStefano Zampini matis->sf_leafdata[i] = PetscMin(matis->sf_leafdata[i],nlocalcols); 1162e1947a5SStefano Zampini } 11728f4e0baSStefano Zampini ierr = MatSeqAIJSetPreallocation(matis->A,0,matis->sf_leafdata);CHKERRQ(ierr); 11828f4e0baSStefano Zampini for (i=0;i<matis->sf_nleaves/bs;i++) { 11928f4e0baSStefano Zampini matis->sf_leafdata[i] = matis->sf_leafdata[i*bs]/bs; 1202e1947a5SStefano Zampini } 12128f4e0baSStefano Zampini ierr = MatSeqBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr); 12228f4e0baSStefano Zampini for (i=0;i<matis->sf_nleaves/bs;i++) { 12328f4e0baSStefano Zampini matis->sf_leafdata[i] = matis->sf_leafdata[i]-i; 1242e1947a5SStefano Zampini } 12528f4e0baSStefano Zampini ierr = MatSeqSBAIJSetPreallocation(matis->A,bs,0,matis->sf_leafdata);CHKERRQ(ierr); 1262e1947a5SStefano Zampini PetscFunctionReturn(0); 1272e1947a5SStefano Zampini } 128b4319ba4SBarry Smith 129b4319ba4SBarry Smith #undef __FUNCT__ 1303927de2eSStefano Zampini #define __FUNCT__ "MatISSetMPIXAIJPreallocation_Private" 1313927de2eSStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetMPIXAIJPreallocation_Private(Mat A, Mat B, PetscBool maxreduce) 1323927de2eSStefano Zampini { 1333927de2eSStefano Zampini Mat_IS *matis = (Mat_IS*)(A->data); 1343927de2eSStefano Zampini PetscInt *my_dnz,*my_onz,*dnz,*onz,*mat_ranges,*row_ownership; 135ecf5a873SStefano Zampini const PetscInt *global_indices_r,*global_indices_c; 1363927de2eSStefano Zampini PetscInt i,j,bs,rows,cols; 1373927de2eSStefano Zampini PetscInt lrows,lcols; 1383927de2eSStefano Zampini PetscInt local_rows,local_cols; 1393927de2eSStefano Zampini PetscMPIInt nsubdomains; 1403927de2eSStefano Zampini PetscBool isdense,issbaij; 1413927de2eSStefano Zampini PetscErrorCode ierr; 1423927de2eSStefano Zampini 1433927de2eSStefano Zampini PetscFunctionBegin; 1443927de2eSStefano Zampini ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&nsubdomains);CHKERRQ(ierr); 1453927de2eSStefano Zampini ierr = MatGetSize(A,&rows,&cols);CHKERRQ(ierr); 1463927de2eSStefano Zampini ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 1473927de2eSStefano Zampini ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr); 1483927de2eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr); 1493927de2eSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr); 150ecf5a873SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr); 151ecf5a873SStefano Zampini if (A->rmap->mapping != A->cmap->mapping) { 1520ea065fbSStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(A->cmap->mapping,&global_indices_c);CHKERRQ(ierr); 153ecf5a873SStefano Zampini } else { 154ecf5a873SStefano Zampini global_indices_c = global_indices_r; 155ecf5a873SStefano Zampini } 156ecf5a873SStefano Zampini 1573927de2eSStefano Zampini if (issbaij) { 1583927de2eSStefano Zampini ierr = MatGetRowUpperTriangular(matis->A);CHKERRQ(ierr); 1593927de2eSStefano Zampini } 1603927de2eSStefano Zampini /* 161ecf5a873SStefano Zampini An SF reduce is needed to sum up properly on shared rows. 1623927de2eSStefano Zampini Note that generally preallocation is not exact, since it overestimates nonzeros 1633927de2eSStefano Zampini */ 1643927de2eSStefano Zampini if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */ 1653927de2eSStefano Zampini ierr = MatISComputeSF_Private(A);CHKERRQ(ierr); 1663927de2eSStefano Zampini } 1673927de2eSStefano Zampini ierr = MatGetLocalSize(A,&lrows,&lcols);CHKERRQ(ierr); 1683927de2eSStefano Zampini ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)A),lrows,lcols,dnz,onz);CHKERRQ(ierr); 1693927de2eSStefano Zampini /* All processes need to compute entire row ownership */ 1703927de2eSStefano Zampini ierr = PetscMalloc1(rows,&row_ownership);CHKERRQ(ierr); 1713927de2eSStefano Zampini ierr = MatGetOwnershipRanges(A,(const PetscInt**)&mat_ranges);CHKERRQ(ierr); 1723927de2eSStefano Zampini for (i=0;i<nsubdomains;i++) { 1733927de2eSStefano Zampini for (j=mat_ranges[i];j<mat_ranges[i+1];j++) { 1743927de2eSStefano Zampini row_ownership[j] = i; 1753927de2eSStefano Zampini } 1763927de2eSStefano Zampini } 1770ea065fbSStefano Zampini ierr = MatGetOwnershipRangesColumn(A,(const PetscInt**)&mat_ranges);CHKERRQ(ierr); 1783927de2eSStefano Zampini 1793927de2eSStefano Zampini /* 1803927de2eSStefano Zampini my_dnz and my_onz contains exact contribution to preallocation from each local mat 1813927de2eSStefano Zampini then, they will be summed up properly. This way, preallocation is always sufficient 1823927de2eSStefano Zampini */ 1833927de2eSStefano Zampini ierr = PetscCalloc2(local_rows,&my_dnz,local_rows,&my_onz);CHKERRQ(ierr); 1843927de2eSStefano Zampini /* preallocation as a MATAIJ */ 1853927de2eSStefano Zampini if (isdense) { /* special case for dense local matrices */ 1863927de2eSStefano Zampini for (i=0;i<local_rows;i++) { 187*12dfadf8SStefano Zampini PetscInt owner = row_ownership[global_indices_r[i]]; 188*12dfadf8SStefano Zampini for (j=0;j<local_cols;j++) { 189ecf5a873SStefano Zampini PetscInt index_col = global_indices_c[j]; 1903927de2eSStefano Zampini if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */ 1913927de2eSStefano Zampini my_dnz[i] += 1; 1923927de2eSStefano Zampini } else { /* offdiag block */ 1933927de2eSStefano Zampini my_onz[i] += 1; 1943927de2eSStefano Zampini } 1953927de2eSStefano Zampini } 1963927de2eSStefano Zampini } 197ecf5a873SStefano Zampini } else { /* TODO: this could be optimized using MatGetRowIJ */ 1983927de2eSStefano Zampini for (i=0;i<local_rows;i++) { 1993927de2eSStefano Zampini const PetscInt *cols; 200ecf5a873SStefano Zampini PetscInt ncols,index_row = global_indices_r[i]; 2013927de2eSStefano Zampini ierr = MatGetRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr); 2023927de2eSStefano Zampini for (j=0;j<ncols;j++) { 2033927de2eSStefano Zampini PetscInt owner = row_ownership[index_row]; 204ecf5a873SStefano Zampini PetscInt index_col = global_indices_c[cols[j]]; 2053927de2eSStefano Zampini if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */ 2063927de2eSStefano Zampini my_dnz[i] += 1; 2073927de2eSStefano Zampini } else { /* offdiag block */ 2083927de2eSStefano Zampini my_onz[i] += 1; 2093927de2eSStefano Zampini } 2103927de2eSStefano Zampini /* same as before, interchanging rows and cols */ 211d9a9e74cSStefano Zampini if (issbaij && index_col != index_row) { 2123927de2eSStefano Zampini owner = row_ownership[index_col]; 2133927de2eSStefano Zampini if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) { 214d9a9e74cSStefano Zampini my_dnz[cols[j]] += 1; 2153927de2eSStefano Zampini } else { 216d9a9e74cSStefano Zampini my_onz[cols[j]] += 1; 2173927de2eSStefano Zampini } 2183927de2eSStefano Zampini } 2193927de2eSStefano Zampini } 2203927de2eSStefano Zampini ierr = MatRestoreRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr); 2213927de2eSStefano Zampini } 2223927de2eSStefano Zampini } 223ecf5a873SStefano Zampini if (global_indices_c != global_indices_r) { 2240ea065fbSStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(A->cmap->mapping,&global_indices_c);CHKERRQ(ierr); 225ecf5a873SStefano Zampini } 2264f619741Sstefano_zampini ierr = ISLocalToGlobalMappingRestoreIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr); 2273927de2eSStefano Zampini ierr = PetscFree(row_ownership);CHKERRQ(ierr); 228ecf5a873SStefano Zampini 229ecf5a873SStefano Zampini /* Reduce my_dnz and my_onz */ 2303927de2eSStefano Zampini if (maxreduce) { 2313927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr); 2323927de2eSStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr); 2333927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr); 2343927de2eSStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr); 2353927de2eSStefano Zampini } else { 2363927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr); 2373927de2eSStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr); 2383927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr); 2393927de2eSStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr); 2403927de2eSStefano Zampini } 2413927de2eSStefano Zampini ierr = PetscFree2(my_dnz,my_onz);CHKERRQ(ierr); 2423927de2eSStefano Zampini 2433927de2eSStefano Zampini /* Resize preallocation if overestimated */ 2443927de2eSStefano Zampini for (i=0;i<lrows;i++) { 2453927de2eSStefano Zampini dnz[i] = PetscMin(dnz[i],lcols); 2463927de2eSStefano Zampini onz[i] = PetscMin(onz[i],cols-lcols); 2473927de2eSStefano Zampini } 2483927de2eSStefano Zampini /* set preallocation */ 2493927de2eSStefano Zampini ierr = MatMPIAIJSetPreallocation(B,0,dnz,0,onz);CHKERRQ(ierr); 2503927de2eSStefano Zampini for (i=0;i<lrows/bs;i++) { 2513927de2eSStefano Zampini dnz[i] = dnz[i*bs]/bs; 2523927de2eSStefano Zampini onz[i] = onz[i*bs]/bs; 2533927de2eSStefano Zampini } 2543927de2eSStefano Zampini ierr = MatMPIBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr); 2553927de2eSStefano Zampini ierr = MatMPISBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr); 2563927de2eSStefano Zampini ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 2573927de2eSStefano Zampini if (issbaij) { 2583927de2eSStefano Zampini ierr = MatRestoreRowUpperTriangular(matis->A);CHKERRQ(ierr); 2593927de2eSStefano Zampini } 2603927de2eSStefano Zampini PetscFunctionReturn(0); 2613927de2eSStefano Zampini } 2623927de2eSStefano Zampini 2633927de2eSStefano Zampini #undef __FUNCT__ 264b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ_IS" 265b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ_IS(Mat mat, MatReuse reuse, Mat *M) 266b7ce53b6SStefano Zampini { 267b7ce53b6SStefano Zampini Mat_IS *matis = (Mat_IS*)(mat->data); 268d9a9e74cSStefano Zampini Mat local_mat; 269b7ce53b6SStefano Zampini /* info on mat */ 2703cfa4ea4SStefano Zampini PetscInt bs,rows,cols,lrows,lcols; 271b7ce53b6SStefano Zampini PetscInt local_rows,local_cols; 272686e3a49SStefano Zampini PetscBool isdense,issbaij,isseqaij; 2737c03b4e8SStefano Zampini PetscMPIInt nsubdomains; 274b7ce53b6SStefano Zampini /* values insertion */ 275b7ce53b6SStefano Zampini PetscScalar *array; 276b7ce53b6SStefano Zampini /* work */ 277b7ce53b6SStefano Zampini PetscErrorCode ierr; 278b7ce53b6SStefano Zampini 279b7ce53b6SStefano Zampini PetscFunctionBegin; 280b7ce53b6SStefano Zampini /* get info from mat */ 2817c03b4e8SStefano Zampini ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&nsubdomains);CHKERRQ(ierr); 2827c03b4e8SStefano Zampini if (nsubdomains == 1) { 2837c03b4e8SStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { 2847c03b4e8SStefano Zampini ierr = MatDuplicate(matis->A,MAT_COPY_VALUES,&(*M));CHKERRQ(ierr); 2857c03b4e8SStefano Zampini } else { 2867c03b4e8SStefano Zampini ierr = MatCopy(matis->A,*M,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 2877c03b4e8SStefano Zampini } 2887c03b4e8SStefano Zampini PetscFunctionReturn(0); 2897c03b4e8SStefano Zampini } 290b7ce53b6SStefano Zampini ierr = MatGetSize(mat,&rows,&cols);CHKERRQ(ierr); 291b7ce53b6SStefano Zampini ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 2923cfa4ea4SStefano Zampini ierr = MatGetLocalSize(mat,&lrows,&lcols);CHKERRQ(ierr); 293b7ce53b6SStefano Zampini ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr); 294b7ce53b6SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr); 295686e3a49SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr); 296b7ce53b6SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr); 297b7ce53b6SStefano Zampini 298b7ce53b6SStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { 299b7ce53b6SStefano Zampini MatType new_mat_type; 3003927de2eSStefano Zampini PetscBool issbaij_red; 301b7ce53b6SStefano Zampini 302b7ce53b6SStefano Zampini /* determining new matrix type */ 303b2566f29SBarry Smith ierr = MPIU_Allreduce(&issbaij,&issbaij_red,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 304b7ce53b6SStefano Zampini if (issbaij_red) { 305b7ce53b6SStefano Zampini new_mat_type = MATSBAIJ; 306b7ce53b6SStefano Zampini } else { 307b7ce53b6SStefano Zampini if (bs>1) { 308b7ce53b6SStefano Zampini new_mat_type = MATBAIJ; 309b7ce53b6SStefano Zampini } else { 310b7ce53b6SStefano Zampini new_mat_type = MATAIJ; 311b7ce53b6SStefano Zampini } 312b7ce53b6SStefano Zampini } 313b7ce53b6SStefano Zampini 3143927de2eSStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)mat),M);CHKERRQ(ierr); 3153cfa4ea4SStefano Zampini ierr = MatSetSizes(*M,lrows,lcols,rows,cols);CHKERRQ(ierr); 3163927de2eSStefano Zampini ierr = MatSetBlockSize(*M,bs);CHKERRQ(ierr); 3173927de2eSStefano Zampini ierr = MatSetType(*M,new_mat_type);CHKERRQ(ierr); 3183927de2eSStefano Zampini ierr = MatISSetMPIXAIJPreallocation_Private(mat,*M,PETSC_FALSE);CHKERRQ(ierr); 319b7ce53b6SStefano Zampini } else { 3203cfa4ea4SStefano Zampini PetscInt mbs,mrows,mcols,mlrows,mlcols; 321b7ce53b6SStefano Zampini /* some checks */ 322b7ce53b6SStefano Zampini ierr = MatGetBlockSize(*M,&mbs);CHKERRQ(ierr); 323b7ce53b6SStefano Zampini ierr = MatGetSize(*M,&mrows,&mcols);CHKERRQ(ierr); 3243cfa4ea4SStefano Zampini ierr = MatGetLocalSize(*M,&mlrows,&mlcols);CHKERRQ(ierr); 3256c4ed002SBarry Smith if (mrows != rows) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of rows (%d != %d)",rows,mrows); 3266c4ed002SBarry Smith if (mcols != cols) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of cols (%d != %d)",cols,mcols); 3276c4ed002SBarry Smith if (mlrows != lrows) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local rows (%d != %d)",lrows,mlrows); 3286c4ed002SBarry Smith if (mlcols != lcols) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local cols (%d != %d)",lcols,mlcols); 3296c4ed002SBarry Smith if (mbs != bs) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong block size (%d != %d)",bs,mbs); 330b7ce53b6SStefano Zampini ierr = MatZeroEntries(*M);CHKERRQ(ierr); 331b7ce53b6SStefano Zampini } 332d9a9e74cSStefano Zampini 333d9a9e74cSStefano Zampini if (issbaij) { 334d9a9e74cSStefano Zampini ierr = MatConvert(matis->A,MATSEQBAIJ,MAT_INITIAL_MATRIX,&local_mat);CHKERRQ(ierr); 335d9a9e74cSStefano Zampini } else { 336d9a9e74cSStefano Zampini ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr); 337d9a9e74cSStefano Zampini local_mat = matis->A; 338d9a9e74cSStefano Zampini } 339686e3a49SStefano Zampini 340b7ce53b6SStefano Zampini /* Set values */ 341ecf5a873SStefano Zampini ierr = MatSetLocalToGlobalMapping(*M,mat->rmap->mapping,mat->cmap->mapping);CHKERRQ(ierr); 342b7ce53b6SStefano Zampini if (isdense) { /* special case for dense local matrices */ 343ecf5a873SStefano Zampini PetscInt i,*dummy_rows,*dummy_cols; 344ecf5a873SStefano Zampini 345ecf5a873SStefano Zampini if (local_rows != local_cols) { 346ecf5a873SStefano Zampini ierr = PetscMalloc2(local_rows,&dummy_rows,local_cols,&dummy_cols);CHKERRQ(ierr); 347ecf5a873SStefano Zampini for (i=0;i<local_rows;i++) dummy_rows[i] = i; 348ecf5a873SStefano Zampini for (i=0;i<local_cols;i++) dummy_cols[i] = i; 349ecf5a873SStefano Zampini } else { 350ecf5a873SStefano Zampini ierr = PetscMalloc1(local_rows,&dummy_rows);CHKERRQ(ierr); 351ecf5a873SStefano Zampini for (i=0;i<local_rows;i++) dummy_rows[i] = i; 352ecf5a873SStefano Zampini dummy_cols = dummy_rows; 353ecf5a873SStefano Zampini } 354b7ce53b6SStefano Zampini ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr); 355d9a9e74cSStefano Zampini ierr = MatDenseGetArray(local_mat,&array);CHKERRQ(ierr); 356ecf5a873SStefano Zampini ierr = MatSetValuesLocal(*M,local_rows,dummy_rows,local_cols,dummy_cols,array,ADD_VALUES);CHKERRQ(ierr); 357d9a9e74cSStefano Zampini ierr = MatDenseRestoreArray(local_mat,&array);CHKERRQ(ierr); 358ecf5a873SStefano Zampini if (dummy_rows != dummy_cols) { 359ecf5a873SStefano Zampini ierr = PetscFree2(dummy_rows,dummy_cols);CHKERRQ(ierr); 360ecf5a873SStefano Zampini } else { 361ecf5a873SStefano Zampini ierr = PetscFree(dummy_rows);CHKERRQ(ierr); 362ecf5a873SStefano Zampini } 363686e3a49SStefano Zampini } else if (isseqaij) { 364ecf5a873SStefano Zampini PetscInt i,nvtxs,*xadj,*adjncy; 365686e3a49SStefano Zampini PetscBool done; 366686e3a49SStefano Zampini 367d9a9e74cSStefano Zampini ierr = MatGetRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr); 3686c4ed002SBarry Smith if (!done) SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__); 369d9a9e74cSStefano Zampini ierr = MatSeqAIJGetArray(local_mat,&array);CHKERRQ(ierr); 370686e3a49SStefano Zampini for (i=0;i<nvtxs;i++) { 371ecf5a873SStefano Zampini ierr = MatSetValuesLocal(*M,1,&i,xadj[i+1]-xadj[i],adjncy+xadj[i],array+xadj[i],ADD_VALUES);CHKERRQ(ierr); 372686e3a49SStefano Zampini } 373d9a9e74cSStefano Zampini ierr = MatRestoreRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr); 3746c4ed002SBarry Smith if (!done) SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__); 375d9a9e74cSStefano Zampini ierr = MatSeqAIJRestoreArray(local_mat,&array);CHKERRQ(ierr); 376686e3a49SStefano Zampini } else { /* very basic values insertion for all other matrix types */ 377ecf5a873SStefano Zampini PetscInt i; 378c0962df8SStefano Zampini 379686e3a49SStefano Zampini for (i=0;i<local_rows;i++) { 380686e3a49SStefano Zampini PetscInt j; 381ecf5a873SStefano Zampini const PetscInt *local_indices_cols; 382686e3a49SStefano Zampini 383ecf5a873SStefano Zampini ierr = MatGetRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr); 384ecf5a873SStefano Zampini ierr = MatSetValuesLocal(*M,1,&i,j,local_indices_cols,array,ADD_VALUES);CHKERRQ(ierr); 385ecf5a873SStefano Zampini ierr = MatRestoreRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr); 386686e3a49SStefano Zampini } 387b7ce53b6SStefano Zampini } 388b7ce53b6SStefano Zampini ierr = MatAssemblyBegin(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 389d9a9e74cSStefano Zampini ierr = MatDestroy(&local_mat);CHKERRQ(ierr); 390b7ce53b6SStefano Zampini ierr = MatAssemblyEnd(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 391b7ce53b6SStefano Zampini if (isdense) { 392b7ce53b6SStefano Zampini ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr); 393b7ce53b6SStefano Zampini } 394b7ce53b6SStefano Zampini PetscFunctionReturn(0); 395b7ce53b6SStefano Zampini } 396b7ce53b6SStefano Zampini 397b7ce53b6SStefano Zampini #undef __FUNCT__ 398b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ" 399b7ce53b6SStefano Zampini /*@ 400b7ce53b6SStefano Zampini MatISGetMPIXAIJ - Converts MATIS matrix into a parallel AIJ format 401b7ce53b6SStefano Zampini 402b7ce53b6SStefano Zampini Input Parameter: 403b7ce53b6SStefano Zampini . mat - the matrix (should be of type MATIS) 404b7ce53b6SStefano Zampini . reuse - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX 405b7ce53b6SStefano Zampini 406b7ce53b6SStefano Zampini Output Parameter: 407b7ce53b6SStefano Zampini . newmat - the matrix in AIJ format 408b7ce53b6SStefano Zampini 409b7ce53b6SStefano Zampini Level: developer 410b7ce53b6SStefano Zampini 411eb82efa4SStefano Zampini Notes: mat and *newmat cannot be the same object when MAT_REUSE_MATRIX is requested. 412b7ce53b6SStefano Zampini 413b7ce53b6SStefano Zampini .seealso: MATIS 414b7ce53b6SStefano Zampini @*/ 415b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ(Mat mat, MatReuse reuse, Mat *newmat) 416b7ce53b6SStefano Zampini { 417b7ce53b6SStefano Zampini PetscErrorCode ierr; 418b7ce53b6SStefano Zampini 419b7ce53b6SStefano Zampini PetscFunctionBegin; 420b7ce53b6SStefano Zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 421b7ce53b6SStefano Zampini PetscValidLogicalCollectiveEnum(mat,reuse,2); 422b7ce53b6SStefano Zampini PetscValidPointer(newmat,3); 423b7ce53b6SStefano Zampini if (reuse != MAT_INITIAL_MATRIX) { 424b7ce53b6SStefano Zampini PetscValidHeaderSpecific(*newmat,MAT_CLASSID,3); 425b7ce53b6SStefano Zampini PetscCheckSameComm(mat,1,*newmat,3); 4266c4ed002SBarry Smith if (mat == *newmat) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse the same matrix"); 427b7ce53b6SStefano Zampini } 428b7ce53b6SStefano Zampini ierr = PetscUseMethod(mat,"MatISGetMPIXAIJ_C",(Mat,MatReuse,Mat*),(mat,reuse,newmat));CHKERRQ(ierr); 429b7ce53b6SStefano Zampini PetscFunctionReturn(0); 430b7ce53b6SStefano Zampini } 431b7ce53b6SStefano Zampini 432b7ce53b6SStefano Zampini #undef __FUNCT__ 433ad6194a2SStefano Zampini #define __FUNCT__ "MatDuplicate_IS" 434ad6194a2SStefano Zampini PetscErrorCode MatDuplicate_IS(Mat mat,MatDuplicateOption op,Mat *newmat) 435ad6194a2SStefano Zampini { 436ad6194a2SStefano Zampini PetscErrorCode ierr; 437ad6194a2SStefano Zampini Mat_IS *matis = (Mat_IS*)(mat->data); 438ad6194a2SStefano Zampini PetscInt bs,m,n,M,N; 439ad6194a2SStefano Zampini Mat B,localmat; 440ad6194a2SStefano Zampini 441ad6194a2SStefano Zampini PetscFunctionBegin; 442ad6194a2SStefano Zampini ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 443ad6194a2SStefano Zampini ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 444ad6194a2SStefano Zampini ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 445e176bc59SStefano Zampini ierr = MatCreateIS(PetscObjectComm((PetscObject)mat),bs,m,n,M,N,mat->rmap->mapping,mat->cmap->mapping,&B);CHKERRQ(ierr); 446ad6194a2SStefano Zampini ierr = MatDuplicate(matis->A,op,&localmat);CHKERRQ(ierr); 447ad6194a2SStefano Zampini ierr = MatISSetLocalMat(B,localmat);CHKERRQ(ierr); 448b3317aa8SStefano Zampini ierr = MatDestroy(&localmat);CHKERRQ(ierr); 449ad6194a2SStefano Zampini ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 450ad6194a2SStefano Zampini ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 451ad6194a2SStefano Zampini *newmat = B; 452ad6194a2SStefano Zampini PetscFunctionReturn(0); 453ad6194a2SStefano Zampini } 454ad6194a2SStefano Zampini 455ad6194a2SStefano Zampini #undef __FUNCT__ 45669796d55SStefano Zampini #define __FUNCT__ "MatIsHermitian_IS" 45769796d55SStefano Zampini PetscErrorCode MatIsHermitian_IS(Mat A,PetscReal tol,PetscBool *flg) 45869796d55SStefano Zampini { 45969796d55SStefano Zampini PetscErrorCode ierr; 46069796d55SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 46169796d55SStefano Zampini PetscBool local_sym; 46269796d55SStefano Zampini 46369796d55SStefano Zampini PetscFunctionBegin; 46469796d55SStefano Zampini ierr = MatIsHermitian(matis->A,tol,&local_sym);CHKERRQ(ierr); 465b2566f29SBarry Smith ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 46669796d55SStefano Zampini PetscFunctionReturn(0); 46769796d55SStefano Zampini } 46869796d55SStefano Zampini 46969796d55SStefano Zampini #undef __FUNCT__ 47069796d55SStefano Zampini #define __FUNCT__ "MatIsSymmetric_IS" 47169796d55SStefano Zampini PetscErrorCode MatIsSymmetric_IS(Mat A,PetscReal tol,PetscBool *flg) 47269796d55SStefano Zampini { 47369796d55SStefano Zampini PetscErrorCode ierr; 47469796d55SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 47569796d55SStefano Zampini PetscBool local_sym; 47669796d55SStefano Zampini 47769796d55SStefano Zampini PetscFunctionBegin; 47869796d55SStefano Zampini ierr = MatIsSymmetric(matis->A,tol,&local_sym);CHKERRQ(ierr); 479b2566f29SBarry Smith ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 48069796d55SStefano Zampini PetscFunctionReturn(0); 48169796d55SStefano Zampini } 48269796d55SStefano Zampini 48369796d55SStefano Zampini #undef __FUNCT__ 484b4319ba4SBarry Smith #define __FUNCT__ "MatDestroy_IS" 485dfbe8321SBarry Smith PetscErrorCode MatDestroy_IS(Mat A) 486b4319ba4SBarry Smith { 487dfbe8321SBarry Smith PetscErrorCode ierr; 488b4319ba4SBarry Smith Mat_IS *b = (Mat_IS*)A->data; 489b4319ba4SBarry Smith 490b4319ba4SBarry Smith PetscFunctionBegin; 4916bf464f9SBarry Smith ierr = MatDestroy(&b->A);CHKERRQ(ierr); 492e176bc59SStefano Zampini ierr = VecScatterDestroy(&b->cctx);CHKERRQ(ierr); 493e176bc59SStefano Zampini ierr = VecScatterDestroy(&b->rctx);CHKERRQ(ierr); 4946bf464f9SBarry Smith ierr = VecDestroy(&b->x);CHKERRQ(ierr); 4956bf464f9SBarry Smith ierr = VecDestroy(&b->y);CHKERRQ(ierr); 49628f4e0baSStefano Zampini ierr = PetscSFDestroy(&b->sf);CHKERRQ(ierr); 49728f4e0baSStefano Zampini ierr = PetscFree2(b->sf_rootdata,b->sf_leafdata);CHKERRQ(ierr); 498bf0cc555SLisandro Dalcin ierr = PetscFree(A->data);CHKERRQ(ierr); 499dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr); 500bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",NULL);CHKERRQ(ierr); 501b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",NULL);CHKERRQ(ierr); 502b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",NULL);CHKERRQ(ierr); 5032e1947a5SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",NULL);CHKERRQ(ierr); 504b4319ba4SBarry Smith PetscFunctionReturn(0); 505b4319ba4SBarry Smith } 506b4319ba4SBarry Smith 507b4319ba4SBarry Smith #undef __FUNCT__ 508b4319ba4SBarry Smith #define __FUNCT__ "MatMult_IS" 509dfbe8321SBarry Smith PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y) 510b4319ba4SBarry Smith { 511dfbe8321SBarry Smith PetscErrorCode ierr; 512b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 513b4319ba4SBarry Smith PetscScalar zero = 0.0; 514b4319ba4SBarry Smith 515b4319ba4SBarry Smith PetscFunctionBegin; 516b4319ba4SBarry Smith /* scatter the global vector x into the local work vector */ 517e176bc59SStefano Zampini ierr = VecScatterBegin(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 518e176bc59SStefano Zampini ierr = VecScatterEnd(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 519b4319ba4SBarry Smith 520b4319ba4SBarry Smith /* multiply the local matrix */ 521b4319ba4SBarry Smith ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr); 522b4319ba4SBarry Smith 523b4319ba4SBarry Smith /* scatter product back into global memory */ 5242dcb1b2aSMatthew Knepley ierr = VecSet(y,zero);CHKERRQ(ierr); 525e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 526e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 527b4319ba4SBarry Smith PetscFunctionReturn(0); 528b4319ba4SBarry Smith } 529b4319ba4SBarry Smith 530b4319ba4SBarry Smith #undef __FUNCT__ 5312e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultAdd_IS" 532b7ce53b6SStefano Zampini PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 5332e74eeadSLisandro Dalcin { 534650997f4SStefano Zampini Vec temp_vec; 5352e74eeadSLisandro Dalcin PetscErrorCode ierr; 5362e74eeadSLisandro Dalcin 5372e74eeadSLisandro Dalcin PetscFunctionBegin; /* v3 = v2 + A * v1.*/ 538650997f4SStefano Zampini if (v3 != v2) { 539650997f4SStefano Zampini ierr = MatMult(A,v1,v3);CHKERRQ(ierr); 540650997f4SStefano Zampini ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 541650997f4SStefano Zampini } else { 542650997f4SStefano Zampini ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 543650997f4SStefano Zampini ierr = MatMult(A,v1,temp_vec);CHKERRQ(ierr); 544650997f4SStefano Zampini ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 545650997f4SStefano Zampini ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 546650997f4SStefano Zampini ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 547650997f4SStefano Zampini } 5482e74eeadSLisandro Dalcin PetscFunctionReturn(0); 5492e74eeadSLisandro Dalcin } 5502e74eeadSLisandro Dalcin 5512e74eeadSLisandro Dalcin #undef __FUNCT__ 5522e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTranspose_IS" 553e176bc59SStefano Zampini PetscErrorCode MatMultTranspose_IS(Mat A,Vec y,Vec x) 5542e74eeadSLisandro Dalcin { 5552e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 5562e74eeadSLisandro Dalcin PetscErrorCode ierr; 5572e74eeadSLisandro Dalcin 558e176bc59SStefano Zampini PetscFunctionBegin; 5592e74eeadSLisandro Dalcin /* scatter the global vector x into the local work vector */ 560e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 561e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5622e74eeadSLisandro Dalcin 5632e74eeadSLisandro Dalcin /* multiply the local matrix */ 564e176bc59SStefano Zampini ierr = MatMultTranspose(is->A,is->y,is->x);CHKERRQ(ierr); 5652e74eeadSLisandro Dalcin 5662e74eeadSLisandro Dalcin /* scatter product back into global vector */ 567e176bc59SStefano Zampini ierr = VecSet(x,0);CHKERRQ(ierr); 568e176bc59SStefano Zampini ierr = VecScatterBegin(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 569e176bc59SStefano Zampini ierr = VecScatterEnd(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 5702e74eeadSLisandro Dalcin PetscFunctionReturn(0); 5712e74eeadSLisandro Dalcin } 5722e74eeadSLisandro Dalcin 5732e74eeadSLisandro Dalcin #undef __FUNCT__ 5742e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTransposeAdd_IS" 5752e74eeadSLisandro Dalcin PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 5762e74eeadSLisandro Dalcin { 577650997f4SStefano Zampini Vec temp_vec; 5782e74eeadSLisandro Dalcin PetscErrorCode ierr; 5792e74eeadSLisandro Dalcin 5802e74eeadSLisandro Dalcin PetscFunctionBegin; /* v3 = v2 + A' * v1.*/ 581650997f4SStefano Zampini if (v3 != v2) { 582650997f4SStefano Zampini ierr = MatMultTranspose(A,v1,v3);CHKERRQ(ierr); 583650997f4SStefano Zampini ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 584650997f4SStefano Zampini } else { 585650997f4SStefano Zampini ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 586650997f4SStefano Zampini ierr = MatMultTranspose(A,v1,temp_vec);CHKERRQ(ierr); 587650997f4SStefano Zampini ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 588650997f4SStefano Zampini ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 589650997f4SStefano Zampini ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 590650997f4SStefano Zampini } 5912e74eeadSLisandro Dalcin PetscFunctionReturn(0); 5922e74eeadSLisandro Dalcin } 5932e74eeadSLisandro Dalcin 5942e74eeadSLisandro Dalcin #undef __FUNCT__ 595b4319ba4SBarry Smith #define __FUNCT__ "MatView_IS" 596dfbe8321SBarry Smith PetscErrorCode MatView_IS(Mat A,PetscViewer viewer) 597b4319ba4SBarry Smith { 598b4319ba4SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 599dfbe8321SBarry Smith PetscErrorCode ierr; 600b4319ba4SBarry Smith PetscViewer sviewer; 601b4319ba4SBarry Smith 602b4319ba4SBarry Smith PetscFunctionBegin; 6033f08860eSBarry Smith ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 604b4319ba4SBarry Smith ierr = MatView(a->A,sviewer);CHKERRQ(ierr); 6053f08860eSBarry Smith ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 606b4319ba4SBarry Smith PetscFunctionReturn(0); 607b4319ba4SBarry Smith } 608b4319ba4SBarry Smith 609b4319ba4SBarry Smith #undef __FUNCT__ 610b4319ba4SBarry Smith #define __FUNCT__ "MatSetLocalToGlobalMapping_IS" 611784ac674SJed Brown PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping) 612b4319ba4SBarry Smith { 613dfbe8321SBarry Smith PetscErrorCode ierr; 614e176bc59SStefano Zampini PetscInt nr,rbs,nc,cbs; 615b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 616b4319ba4SBarry Smith IS from,to; 617e176bc59SStefano Zampini Vec cglobal,rglobal; 618b4319ba4SBarry Smith 619b4319ba4SBarry Smith PetscFunctionBegin; 620784ac674SJed Brown PetscCheckSameComm(A,1,rmapping,2); 621e176bc59SStefano Zampini PetscCheckSameComm(A,1,cmapping,3); 6223bbff08aSStefano Zampini /* Destroy any previous data */ 62370cf5478SStefano Zampini ierr = VecDestroy(&is->x);CHKERRQ(ierr); 62470cf5478SStefano Zampini ierr = VecDestroy(&is->y);CHKERRQ(ierr); 625e176bc59SStefano Zampini ierr = VecScatterDestroy(&is->rctx);CHKERRQ(ierr); 626e176bc59SStefano Zampini ierr = VecScatterDestroy(&is->cctx);CHKERRQ(ierr); 6271c47cb0fSStefano Zampini ierr = MatDestroy(&is->A);CHKERRQ(ierr); 62828f4e0baSStefano Zampini ierr = PetscSFDestroy(&is->sf);CHKERRQ(ierr); 62928f4e0baSStefano Zampini ierr = PetscFree2(is->sf_rootdata,is->sf_leafdata);CHKERRQ(ierr); 6303bbff08aSStefano Zampini 6313bbff08aSStefano Zampini /* Setup Layout and set local to global maps */ 632fc27028aSStefano Zampini ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 633fc27028aSStefano Zampini ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 634fa7f1dd8SStefano Zampini ierr = PetscLayoutSetISLocalToGlobalMapping(A->rmap,rmapping);CHKERRQ(ierr); 635fa7f1dd8SStefano Zampini ierr = PetscLayoutSetISLocalToGlobalMapping(A->cmap,cmapping);CHKERRQ(ierr); 636b4319ba4SBarry Smith 637b4319ba4SBarry Smith /* Create the local matrix A */ 638e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(rmapping,&nr);CHKERRQ(ierr); 639e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetBlockSize(rmapping,&rbs);CHKERRQ(ierr); 640e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(cmapping,&nc);CHKERRQ(ierr); 641e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetBlockSize(cmapping,&cbs);CHKERRQ(ierr); 642f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr); 643e176bc59SStefano Zampini ierr = MatSetType(is->A,MATAIJ);CHKERRQ(ierr); 644e176bc59SStefano Zampini ierr = MatSetSizes(is->A,nr,nc,nr,nc);CHKERRQ(ierr); 645e176bc59SStefano Zampini ierr = MatSetBlockSizes(is->A,rbs,cbs);CHKERRQ(ierr); 646ff130e51SJed Brown ierr = MatSetOptionsPrefix(is->A,((PetscObject)A)->prefix);CHKERRQ(ierr); 647ff130e51SJed Brown ierr = MatAppendOptionsPrefix(is->A,"is_");CHKERRQ(ierr); 648b4319ba4SBarry Smith ierr = MatSetFromOptions(is->A);CHKERRQ(ierr); 649b4319ba4SBarry Smith 650b4319ba4SBarry Smith /* Create the local work vectors */ 6513bbff08aSStefano Zampini ierr = MatCreateVecs(is->A,&is->x,&is->y);CHKERRQ(ierr); 652b4319ba4SBarry Smith 653e176bc59SStefano Zampini /* setup the global to local scatters */ 654e176bc59SStefano Zampini ierr = MatCreateVecs(A,&cglobal,&rglobal);CHKERRQ(ierr); 655e176bc59SStefano Zampini ierr = ISCreateStride(PETSC_COMM_SELF,nr,0,1,&to);CHKERRQ(ierr); 656784ac674SJed Brown ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr); 657e176bc59SStefano Zampini ierr = VecScatterCreate(rglobal,from,is->y,to,&is->rctx);CHKERRQ(ierr); 658e176bc59SStefano Zampini if (rmapping != cmapping) { 659e176bc59SStefano Zampini ierr = ISDestroy(&to);CHKERRQ(ierr); 660e176bc59SStefano Zampini ierr = ISDestroy(&from);CHKERRQ(ierr); 661e176bc59SStefano Zampini ierr = ISCreateStride(PETSC_COMM_SELF,nc,0,1,&to);CHKERRQ(ierr); 662e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingApplyIS(cmapping,to,&from);CHKERRQ(ierr); 663e176bc59SStefano Zampini ierr = VecScatterCreate(cglobal,from,is->x,to,&is->cctx);CHKERRQ(ierr); 664e176bc59SStefano Zampini } else { 665e176bc59SStefano Zampini ierr = PetscObjectReference((PetscObject)is->rctx);CHKERRQ(ierr); 666e176bc59SStefano Zampini is->cctx = is->rctx; 667e176bc59SStefano Zampini } 668e176bc59SStefano Zampini ierr = VecDestroy(&rglobal);CHKERRQ(ierr); 669e176bc59SStefano Zampini ierr = VecDestroy(&cglobal);CHKERRQ(ierr); 6706bf464f9SBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 6716bf464f9SBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 672b4319ba4SBarry Smith PetscFunctionReturn(0); 673b4319ba4SBarry Smith } 674b4319ba4SBarry Smith 6752e74eeadSLisandro Dalcin #undef __FUNCT__ 6762e74eeadSLisandro Dalcin #define __FUNCT__ "MatSetValues_IS" 6772e74eeadSLisandro Dalcin PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv) 6782e74eeadSLisandro Dalcin { 6792e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)mat->data; 6802e74eeadSLisandro Dalcin PetscInt rows_l[2048],cols_l[2048]; 6812e74eeadSLisandro Dalcin PetscErrorCode ierr; 6822e74eeadSLisandro Dalcin 6832e74eeadSLisandro Dalcin PetscFunctionBegin; 6842e74eeadSLisandro Dalcin #if defined(PETSC_USE_DEBUG) 685f23aa3ddSBarry 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); 6862e74eeadSLisandro Dalcin #endif 6873bbff08aSStefano Zampini ierr = ISG2LMapApply(mat->rmap->mapping,m,rows,rows_l);CHKERRQ(ierr); 6883bbff08aSStefano Zampini ierr = ISG2LMapApply(mat->cmap->mapping,n,cols,cols_l);CHKERRQ(ierr); 6892e74eeadSLisandro Dalcin ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 6902e74eeadSLisandro Dalcin PetscFunctionReturn(0); 6912e74eeadSLisandro Dalcin } 6922e74eeadSLisandro Dalcin 693b4319ba4SBarry Smith #undef __FUNCT__ 694b4319ba4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_IS" 69513f74950SBarry Smith PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 696b4319ba4SBarry Smith { 697dfbe8321SBarry Smith PetscErrorCode ierr; 698b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 699b4319ba4SBarry Smith 700b4319ba4SBarry Smith PetscFunctionBegin; 701b4319ba4SBarry Smith ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 702b4319ba4SBarry Smith PetscFunctionReturn(0); 703b4319ba4SBarry Smith } 704b4319ba4SBarry Smith 705b4319ba4SBarry Smith #undef __FUNCT__ 706f0006bf2SLisandro Dalcin #define __FUNCT__ "MatSetValuesBlockedLocal_IS" 707f0006bf2SLisandro Dalcin PetscErrorCode MatSetValuesBlockedLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 708f0006bf2SLisandro Dalcin { 709f0006bf2SLisandro Dalcin PetscErrorCode ierr; 710f0006bf2SLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 711f0006bf2SLisandro Dalcin 712f0006bf2SLisandro Dalcin PetscFunctionBegin; 713f0006bf2SLisandro Dalcin ierr = MatSetValuesBlocked(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 714f0006bf2SLisandro Dalcin PetscFunctionReturn(0); 715f0006bf2SLisandro Dalcin } 716f0006bf2SLisandro Dalcin 717f0006bf2SLisandro Dalcin #undef __FUNCT__ 7182e74eeadSLisandro Dalcin #define __FUNCT__ "MatZeroRows_IS" 7192b40b63fSBarry Smith PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 7202e74eeadSLisandro Dalcin { 7210298fd71SBarry Smith PetscInt n_l = 0, *rows_l = NULL; 7222e74eeadSLisandro Dalcin PetscErrorCode ierr; 7232e74eeadSLisandro Dalcin 7242e74eeadSLisandro Dalcin PetscFunctionBegin; 725ce94432eSBarry Smith if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support"); 7262e74eeadSLisandro Dalcin if (n) { 727785e854fSJed Brown ierr = PetscMalloc1(n,&rows_l);CHKERRQ(ierr); 7283bbff08aSStefano Zampini ierr = ISGlobalToLocalMappingApply(A->rmap->mapping,IS_GTOLM_DROP,n,rows,&n_l,rows_l);CHKERRQ(ierr); 7292e74eeadSLisandro Dalcin } 7302b40b63fSBarry Smith ierr = MatZeroRowsLocal(A,n_l,rows_l,diag,x,b);CHKERRQ(ierr); 731c31cb41cSBarry Smith ierr = PetscFree(rows_l);CHKERRQ(ierr); 7322e74eeadSLisandro Dalcin PetscFunctionReturn(0); 7332e74eeadSLisandro Dalcin } 7342e74eeadSLisandro Dalcin 7352e74eeadSLisandro Dalcin #undef __FUNCT__ 736b4319ba4SBarry Smith #define __FUNCT__ "MatZeroRowsLocal_IS" 7372b40b63fSBarry Smith PetscErrorCode MatZeroRowsLocal_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 738b4319ba4SBarry Smith { 739b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 740dfbe8321SBarry Smith PetscErrorCode ierr; 741f4df32b1SMatthew Knepley PetscInt i; 742b4319ba4SBarry Smith PetscScalar *array; 743b4319ba4SBarry Smith 744b4319ba4SBarry Smith PetscFunctionBegin; 745ce94432eSBarry Smith if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support"); 746b4319ba4SBarry Smith { 747b4319ba4SBarry Smith /* 748b4319ba4SBarry Smith Set up is->x as a "counting vector". This is in order to MatMult_IS 749b4319ba4SBarry Smith work properly in the interface nodes. 750b4319ba4SBarry Smith */ 751b4319ba4SBarry Smith Vec counter; 752e176bc59SStefano Zampini ierr = MatCreateVecs(A,NULL,&counter);CHKERRQ(ierr); 753e176bc59SStefano Zampini ierr = VecSet(counter,0.);CHKERRQ(ierr); 754e176bc59SStefano Zampini ierr = VecSet(is->y,1.);CHKERRQ(ierr); 755e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,is->y,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 756e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,is->y,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 757e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,counter,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 758e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,counter,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 7596bf464f9SBarry Smith ierr = VecDestroy(&counter);CHKERRQ(ierr); 760b4319ba4SBarry Smith } 761958c9bccSBarry Smith if (!n) { 762b4319ba4SBarry Smith is->pure_neumann = PETSC_TRUE; 763b4319ba4SBarry Smith } else { 764b4319ba4SBarry Smith is->pure_neumann = PETSC_FALSE; 7652205254eSKarl Rupp 766e176bc59SStefano Zampini ierr = VecGetArray(is->y,&array);CHKERRQ(ierr); 7672b40b63fSBarry Smith ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr); 768b4319ba4SBarry Smith for (i=0; i<n; i++) { 769f4df32b1SMatthew Knepley ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr); 770b4319ba4SBarry Smith } 771b4319ba4SBarry Smith ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 772b4319ba4SBarry Smith ierr = MatAssemblyEnd(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 773e176bc59SStefano Zampini ierr = VecRestoreArray(is->y,&array);CHKERRQ(ierr); 774b4319ba4SBarry Smith } 775b4319ba4SBarry Smith PetscFunctionReturn(0); 776b4319ba4SBarry Smith } 777b4319ba4SBarry Smith 778b4319ba4SBarry Smith #undef __FUNCT__ 779b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyBegin_IS" 780dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type) 781b4319ba4SBarry Smith { 782b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 783dfbe8321SBarry Smith PetscErrorCode ierr; 784b4319ba4SBarry Smith 785b4319ba4SBarry Smith PetscFunctionBegin; 786b4319ba4SBarry Smith ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr); 787b4319ba4SBarry Smith PetscFunctionReturn(0); 788b4319ba4SBarry Smith } 789b4319ba4SBarry Smith 790b4319ba4SBarry Smith #undef __FUNCT__ 791b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_IS" 792dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type) 793b4319ba4SBarry Smith { 794b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 795dfbe8321SBarry Smith PetscErrorCode ierr; 796b4319ba4SBarry Smith 797b4319ba4SBarry Smith PetscFunctionBegin; 798b4319ba4SBarry Smith ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr); 799b4319ba4SBarry Smith PetscFunctionReturn(0); 800b4319ba4SBarry Smith } 801b4319ba4SBarry Smith 802b4319ba4SBarry Smith #undef __FUNCT__ 803b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat_IS" 8047087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local) 805b4319ba4SBarry Smith { 806b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)mat->data; 807b4319ba4SBarry Smith 808b4319ba4SBarry Smith PetscFunctionBegin; 809b4319ba4SBarry Smith *local = is->A; 810b4319ba4SBarry Smith PetscFunctionReturn(0); 811b4319ba4SBarry Smith } 812b4319ba4SBarry Smith 813b4319ba4SBarry Smith #undef __FUNCT__ 814b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat" 815b4319ba4SBarry Smith /*@ 816b4319ba4SBarry Smith MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix. 817b4319ba4SBarry Smith 818b4319ba4SBarry Smith Input Parameter: 819b4319ba4SBarry Smith . mat - the matrix 820b4319ba4SBarry Smith 821b4319ba4SBarry Smith Output Parameter: 822eb82efa4SStefano Zampini . local - the local matrix 823b4319ba4SBarry Smith 824b4319ba4SBarry Smith Level: advanced 825b4319ba4SBarry Smith 826b4319ba4SBarry Smith Notes: 827b4319ba4SBarry Smith This can be called if you have precomputed the nonzero structure of the 828b4319ba4SBarry Smith matrix and want to provide it to the inner matrix object to improve the performance 829b4319ba4SBarry Smith of the MatSetValues() operation. 830b4319ba4SBarry Smith 83196a6f129SJed Brown This function does not increase the reference count for the local Mat. Do not destroy it and do not attempt to use 83296a6f129SJed Brown your reference after destroying the parent mat. 83396a6f129SJed Brown 834b4319ba4SBarry Smith .seealso: MATIS 835b4319ba4SBarry Smith @*/ 8367087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local) 837b4319ba4SBarry Smith { 8384ac538c5SBarry Smith PetscErrorCode ierr; 839b4319ba4SBarry Smith 840b4319ba4SBarry Smith PetscFunctionBegin; 8410700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 842b4319ba4SBarry Smith PetscValidPointer(local,2); 8434ac538c5SBarry Smith ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat*),(mat,local));CHKERRQ(ierr); 844b4319ba4SBarry Smith PetscFunctionReturn(0); 845b4319ba4SBarry Smith } 846b4319ba4SBarry Smith 8473b03a366Sstefano_zampini #undef __FUNCT__ 8483b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat_IS" 8493b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat_IS(Mat mat,Mat local) 8503b03a366Sstefano_zampini { 8513b03a366Sstefano_zampini Mat_IS *is = (Mat_IS*)mat->data; 8523b03a366Sstefano_zampini PetscInt nrows,ncols,orows,ocols; 8533b03a366Sstefano_zampini PetscErrorCode ierr; 8543b03a366Sstefano_zampini 8553b03a366Sstefano_zampini PetscFunctionBegin; 8564e4c7dbeSStefano Zampini if (is->A) { 8573b03a366Sstefano_zampini ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr); 8583b03a366Sstefano_zampini ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr); 859f23aa3ddSBarry 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); 8604e4c7dbeSStefano Zampini } 8613b03a366Sstefano_zampini ierr = PetscObjectReference((PetscObject)local);CHKERRQ(ierr); 8623b03a366Sstefano_zampini ierr = MatDestroy(&is->A);CHKERRQ(ierr); 8633b03a366Sstefano_zampini is->A = local; 8643b03a366Sstefano_zampini PetscFunctionReturn(0); 8653b03a366Sstefano_zampini } 8663b03a366Sstefano_zampini 8673b03a366Sstefano_zampini #undef __FUNCT__ 8683b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat" 8693b03a366Sstefano_zampini /*@ 870eb82efa4SStefano Zampini MatISSetLocalMat - Replace the local matrix stored inside a MATIS object. 8713b03a366Sstefano_zampini 8723b03a366Sstefano_zampini Input Parameter: 8733b03a366Sstefano_zampini . mat - the matrix 874eb82efa4SStefano Zampini . local - the local matrix 8753b03a366Sstefano_zampini 8763b03a366Sstefano_zampini Output Parameter: 8773b03a366Sstefano_zampini 8783b03a366Sstefano_zampini Level: advanced 8793b03a366Sstefano_zampini 8803b03a366Sstefano_zampini Notes: 8813b03a366Sstefano_zampini This can be called if you have precomputed the local matrix and 8823b03a366Sstefano_zampini want to provide it to the matrix object MATIS. 8833b03a366Sstefano_zampini 8843b03a366Sstefano_zampini .seealso: MATIS 8853b03a366Sstefano_zampini @*/ 8863b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat(Mat mat,Mat local) 8873b03a366Sstefano_zampini { 8883b03a366Sstefano_zampini PetscErrorCode ierr; 8893b03a366Sstefano_zampini 8903b03a366Sstefano_zampini PetscFunctionBegin; 8913b03a366Sstefano_zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 892b7ce53b6SStefano Zampini PetscValidHeaderSpecific(local,MAT_CLASSID,2); 8933b03a366Sstefano_zampini ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr); 8943b03a366Sstefano_zampini PetscFunctionReturn(0); 8953b03a366Sstefano_zampini } 8963b03a366Sstefano_zampini 8976726f965SBarry Smith #undef __FUNCT__ 8986726f965SBarry Smith #define __FUNCT__ "MatZeroEntries_IS" 8996726f965SBarry Smith PetscErrorCode MatZeroEntries_IS(Mat A) 9006726f965SBarry Smith { 9016726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 9026726f965SBarry Smith PetscErrorCode ierr; 9036726f965SBarry Smith 9046726f965SBarry Smith PetscFunctionBegin; 9056726f965SBarry Smith ierr = MatZeroEntries(a->A);CHKERRQ(ierr); 9066726f965SBarry Smith PetscFunctionReturn(0); 9076726f965SBarry Smith } 9086726f965SBarry Smith 9096726f965SBarry Smith #undef __FUNCT__ 9102e74eeadSLisandro Dalcin #define __FUNCT__ "MatScale_IS" 9112e74eeadSLisandro Dalcin PetscErrorCode MatScale_IS(Mat A,PetscScalar a) 9122e74eeadSLisandro Dalcin { 9132e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 9142e74eeadSLisandro Dalcin PetscErrorCode ierr; 9152e74eeadSLisandro Dalcin 9162e74eeadSLisandro Dalcin PetscFunctionBegin; 9172e74eeadSLisandro Dalcin ierr = MatScale(is->A,a);CHKERRQ(ierr); 9182e74eeadSLisandro Dalcin PetscFunctionReturn(0); 9192e74eeadSLisandro Dalcin } 9202e74eeadSLisandro Dalcin 9212e74eeadSLisandro Dalcin #undef __FUNCT__ 9222e74eeadSLisandro Dalcin #define __FUNCT__ "MatGetDiagonal_IS" 9232e74eeadSLisandro Dalcin PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v) 9242e74eeadSLisandro Dalcin { 9252e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 9262e74eeadSLisandro Dalcin PetscErrorCode ierr; 9272e74eeadSLisandro Dalcin 9282e74eeadSLisandro Dalcin PetscFunctionBegin; 9292e74eeadSLisandro Dalcin /* get diagonal of the local matrix */ 930e176bc59SStefano Zampini ierr = MatGetDiagonal(is->A,is->y);CHKERRQ(ierr); 9312e74eeadSLisandro Dalcin 9322e74eeadSLisandro Dalcin /* scatter diagonal back into global vector */ 9332e74eeadSLisandro Dalcin ierr = VecSet(v,0);CHKERRQ(ierr); 934e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 935e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 9362e74eeadSLisandro Dalcin PetscFunctionReturn(0); 9372e74eeadSLisandro Dalcin } 9382e74eeadSLisandro Dalcin 9392e74eeadSLisandro Dalcin #undef __FUNCT__ 9406726f965SBarry Smith #define __FUNCT__ "MatSetOption_IS" 941ace3abfcSBarry Smith PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg) 9426726f965SBarry Smith { 9436726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 9446726f965SBarry Smith PetscErrorCode ierr; 9456726f965SBarry Smith 9466726f965SBarry Smith PetscFunctionBegin; 9474e0d8c25SBarry Smith ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr); 9486726f965SBarry Smith PetscFunctionReturn(0); 9496726f965SBarry Smith } 9506726f965SBarry Smith 951284134d9SBarry Smith #undef __FUNCT__ 952284134d9SBarry Smith #define __FUNCT__ "MatCreateIS" 953284134d9SBarry Smith /*@ 954284134d9SBarry Smith MatCreateIS - Creates a "process" unassmembled matrix, it is assembled on each 955284134d9SBarry Smith process but not across processes. 956284134d9SBarry Smith 957284134d9SBarry Smith Input Parameters: 958284134d9SBarry Smith + comm - MPI communicator that will share the matrix 959e176bc59SStefano Zampini . bs - block size of the matrix 960df3898eeSBarry Smith . m,n,M,N - local and/or global sizes of the left and right vector used in matrix vector products 961e176bc59SStefano Zampini . rmap - local to global map for rows 962e176bc59SStefano Zampini - cmap - local to global map for cols 963284134d9SBarry Smith 964284134d9SBarry Smith Output Parameter: 965284134d9SBarry Smith . A - the resulting matrix 966284134d9SBarry Smith 9678e6c10adSSatish Balay Level: advanced 9688e6c10adSSatish Balay 969284134d9SBarry Smith Notes: See MATIS for more details 9708cda6cd7SBarry Smith m and n are NOT related to the size of the map, they are the size of the part of the vector owned 971e176bc59SStefano Zampini by that process. The sizes of rmap and cmap define the size of the local matrices. 972e176bc59SStefano Zampini If either rmap or cmap are NULL, than the matrix is assumed to be square 973284134d9SBarry Smith 974284134d9SBarry Smith .seealso: MATIS, MatSetLocalToGlobalMapping() 975284134d9SBarry Smith @*/ 976e176bc59SStefano Zampini PetscErrorCode MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping rmap,ISLocalToGlobalMapping cmap,Mat *A) 977284134d9SBarry Smith { 978284134d9SBarry Smith PetscErrorCode ierr; 979284134d9SBarry Smith 980284134d9SBarry Smith PetscFunctionBegin; 981e176bc59SStefano Zampini if (!rmap && !cmap) SETERRQ(comm,PETSC_ERR_USER,"You need to provide at least one of the mapping"); 982284134d9SBarry Smith ierr = MatCreate(comm,A);CHKERRQ(ierr); 983284134d9SBarry Smith ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 9844f619741Sstefano_zampini ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr); 985284134d9SBarry Smith ierr = MatSetType(*A,MATIS);CHKERRQ(ierr); 9863b03a366Sstefano_zampini ierr = MatSetUp(*A);CHKERRQ(ierr); 987e176bc59SStefano Zampini if (rmap && cmap) { 988e176bc59SStefano Zampini ierr = MatSetLocalToGlobalMapping(*A,rmap,cmap);CHKERRQ(ierr); 989e176bc59SStefano Zampini } else if (!rmap) { 990e176bc59SStefano Zampini ierr = MatSetLocalToGlobalMapping(*A,cmap,cmap);CHKERRQ(ierr); 991e176bc59SStefano Zampini } else { 992e176bc59SStefano Zampini ierr = MatSetLocalToGlobalMapping(*A,rmap,rmap);CHKERRQ(ierr); 993e176bc59SStefano Zampini } 994284134d9SBarry Smith PetscFunctionReturn(0); 995284134d9SBarry Smith } 996284134d9SBarry Smith 997b4319ba4SBarry Smith /*MC 998eb82efa4SStefano Zampini MATIS - MATIS = "is" - A matrix type to be used for using the non-overlapping domain decomposition type preconditioners (e.g. PCBDDC). 999b4319ba4SBarry Smith This stores the matrices in globally unassembled form. Each processor 1000b4319ba4SBarry Smith assembles only its local Neumann problem and the parallel matrix vector 1001b4319ba4SBarry Smith product is handled "implicitly". 1002b4319ba4SBarry Smith 1003b4319ba4SBarry Smith Operations Provided: 10046726f965SBarry Smith + MatMult() 10052e74eeadSLisandro Dalcin . MatMultAdd() 10062e74eeadSLisandro Dalcin . MatMultTranspose() 10072e74eeadSLisandro Dalcin . MatMultTransposeAdd() 10086726f965SBarry Smith . MatZeroEntries() 10096726f965SBarry Smith . MatSetOption() 10102e74eeadSLisandro Dalcin . MatZeroRows() 10116726f965SBarry Smith . MatZeroRowsLocal() 10122e74eeadSLisandro Dalcin . MatSetValues() 10136726f965SBarry Smith . MatSetValuesLocal() 10142e74eeadSLisandro Dalcin . MatScale() 10152e74eeadSLisandro Dalcin . MatGetDiagonal() 10166726f965SBarry Smith - MatSetLocalToGlobalMapping() 1017b4319ba4SBarry Smith 1018b4319ba4SBarry Smith Options Database Keys: 1019b4319ba4SBarry Smith . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions() 1020b4319ba4SBarry Smith 1021b4319ba4SBarry Smith Notes: Options prefix for the inner matrix are given by -is_mat_xxx 1022b4319ba4SBarry Smith 1023b4319ba4SBarry Smith You must call MatSetLocalToGlobalMapping() before using this matrix type. 1024b4319ba4SBarry Smith 1025b4319ba4SBarry Smith You can do matrix preallocation on the local matrix after you obtain it with 1026eb82efa4SStefano Zampini MatISGetLocalMat(); otherwise, you could use MatISSetPreallocation() 1027b4319ba4SBarry Smith 1028b4319ba4SBarry Smith Level: advanced 1029b4319ba4SBarry Smith 1030eb82efa4SStefano Zampini .seealso: Mat, MatISGetLocalMat(), MatSetLocalToGlobalMapping(), MatISSetPreallocation(), PCBDDC 1031b4319ba4SBarry Smith 1032b4319ba4SBarry Smith M*/ 1033b4319ba4SBarry Smith 1034b4319ba4SBarry Smith #undef __FUNCT__ 1035b4319ba4SBarry Smith #define __FUNCT__ "MatCreate_IS" 10368cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_IS(Mat A) 1037b4319ba4SBarry Smith { 1038dfbe8321SBarry Smith PetscErrorCode ierr; 1039b4319ba4SBarry Smith Mat_IS *b; 1040b4319ba4SBarry Smith 1041b4319ba4SBarry Smith PetscFunctionBegin; 1042b00a9115SJed Brown ierr = PetscNewLog(A,&b);CHKERRQ(ierr); 1043b4319ba4SBarry Smith A->data = (void*)b; 1044b4319ba4SBarry Smith 1045e176bc59SStefano Zampini /* matrix ops */ 1046e176bc59SStefano Zampini ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1047b4319ba4SBarry Smith A->ops->mult = MatMult_IS; 10482e74eeadSLisandro Dalcin A->ops->multadd = MatMultAdd_IS; 10492e74eeadSLisandro Dalcin A->ops->multtranspose = MatMultTranspose_IS; 10502e74eeadSLisandro Dalcin A->ops->multtransposeadd = MatMultTransposeAdd_IS; 1051b4319ba4SBarry Smith A->ops->destroy = MatDestroy_IS; 1052b4319ba4SBarry Smith A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS; 10532e74eeadSLisandro Dalcin A->ops->setvalues = MatSetValues_IS; 1054b4319ba4SBarry Smith A->ops->setvalueslocal = MatSetValuesLocal_IS; 1055f0006bf2SLisandro Dalcin A->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_IS; 10562e74eeadSLisandro Dalcin A->ops->zerorows = MatZeroRows_IS; 1057b4319ba4SBarry Smith A->ops->zerorowslocal = MatZeroRowsLocal_IS; 1058b4319ba4SBarry Smith A->ops->assemblybegin = MatAssemblyBegin_IS; 1059b4319ba4SBarry Smith A->ops->assemblyend = MatAssemblyEnd_IS; 1060b4319ba4SBarry Smith A->ops->view = MatView_IS; 10616726f965SBarry Smith A->ops->zeroentries = MatZeroEntries_IS; 10622e74eeadSLisandro Dalcin A->ops->scale = MatScale_IS; 10632e74eeadSLisandro Dalcin A->ops->getdiagonal = MatGetDiagonal_IS; 10646726f965SBarry Smith A->ops->setoption = MatSetOption_IS; 106569796d55SStefano Zampini A->ops->ishermitian = MatIsHermitian_IS; 106669796d55SStefano Zampini A->ops->issymmetric = MatIsSymmetric_IS; 1067ad6194a2SStefano Zampini A->ops->duplicate = MatDuplicate_IS; 1068b4319ba4SBarry Smith 1069b7ce53b6SStefano Zampini /* special MATIS functions */ 1070bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",MatISGetLocalMat_IS);CHKERRQ(ierr); 1071bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",MatISSetLocalMat_IS);CHKERRQ(ierr); 1072b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",MatISGetMPIXAIJ_IS);CHKERRQ(ierr); 10732e1947a5SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",MatISSetPreallocation_IS);CHKERRQ(ierr); 107417667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr); 1075b4319ba4SBarry Smith PetscFunctionReturn(0); 1076b4319ba4SBarry Smith } 1077