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++) { 18712dfadf8SStefano Zampini PetscInt owner = row_ownership[global_indices_r[i]]; 18812dfadf8SStefano 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 */ 343*65066ba5SStefano Zampini PetscInt i,*dummy; 344ecf5a873SStefano Zampini 345*65066ba5SStefano Zampini ierr = PetscMalloc1(PetscMax(local_rows,local_cols),&dummy);CHKERRQ(ierr); 346*65066ba5SStefano Zampini for (i=0;i<PetscMax(local_rows,local_cols);i++) dummy[i] = i; 347b7ce53b6SStefano Zampini ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr); 348d9a9e74cSStefano Zampini ierr = MatDenseGetArray(local_mat,&array);CHKERRQ(ierr); 349*65066ba5SStefano Zampini ierr = MatSetValuesLocal(*M,local_rows,dummy,local_cols,dummy,array,ADD_VALUES);CHKERRQ(ierr); 350d9a9e74cSStefano Zampini ierr = MatDenseRestoreArray(local_mat,&array);CHKERRQ(ierr); 351*65066ba5SStefano Zampini ierr = PetscFree(dummy);CHKERRQ(ierr); 352686e3a49SStefano Zampini } else if (isseqaij) { 353ecf5a873SStefano Zampini PetscInt i,nvtxs,*xadj,*adjncy; 354686e3a49SStefano Zampini PetscBool done; 355686e3a49SStefano Zampini 356d9a9e74cSStefano Zampini ierr = MatGetRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr); 3576c4ed002SBarry Smith if (!done) SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__); 358d9a9e74cSStefano Zampini ierr = MatSeqAIJGetArray(local_mat,&array);CHKERRQ(ierr); 359686e3a49SStefano Zampini for (i=0;i<nvtxs;i++) { 360ecf5a873SStefano Zampini ierr = MatSetValuesLocal(*M,1,&i,xadj[i+1]-xadj[i],adjncy+xadj[i],array+xadj[i],ADD_VALUES);CHKERRQ(ierr); 361686e3a49SStefano Zampini } 362d9a9e74cSStefano Zampini ierr = MatRestoreRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr); 3636c4ed002SBarry Smith if (!done) SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__); 364d9a9e74cSStefano Zampini ierr = MatSeqAIJRestoreArray(local_mat,&array);CHKERRQ(ierr); 365686e3a49SStefano Zampini } else { /* very basic values insertion for all other matrix types */ 366ecf5a873SStefano Zampini PetscInt i; 367c0962df8SStefano Zampini 368686e3a49SStefano Zampini for (i=0;i<local_rows;i++) { 369686e3a49SStefano Zampini PetscInt j; 370ecf5a873SStefano Zampini const PetscInt *local_indices_cols; 371686e3a49SStefano Zampini 372ecf5a873SStefano Zampini ierr = MatGetRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr); 373ecf5a873SStefano Zampini ierr = MatSetValuesLocal(*M,1,&i,j,local_indices_cols,array,ADD_VALUES);CHKERRQ(ierr); 374ecf5a873SStefano Zampini ierr = MatRestoreRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr); 375686e3a49SStefano Zampini } 376b7ce53b6SStefano Zampini } 377b7ce53b6SStefano Zampini ierr = MatAssemblyBegin(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 378d9a9e74cSStefano Zampini ierr = MatDestroy(&local_mat);CHKERRQ(ierr); 379b7ce53b6SStefano Zampini ierr = MatAssemblyEnd(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 380b7ce53b6SStefano Zampini if (isdense) { 381b7ce53b6SStefano Zampini ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr); 382b7ce53b6SStefano Zampini } 383b7ce53b6SStefano Zampini PetscFunctionReturn(0); 384b7ce53b6SStefano Zampini } 385b7ce53b6SStefano Zampini 386b7ce53b6SStefano Zampini #undef __FUNCT__ 387b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ" 388b7ce53b6SStefano Zampini /*@ 389b7ce53b6SStefano Zampini MatISGetMPIXAIJ - Converts MATIS matrix into a parallel AIJ format 390b7ce53b6SStefano Zampini 391b7ce53b6SStefano Zampini Input Parameter: 392b7ce53b6SStefano Zampini . mat - the matrix (should be of type MATIS) 393b7ce53b6SStefano Zampini . reuse - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX 394b7ce53b6SStefano Zampini 395b7ce53b6SStefano Zampini Output Parameter: 396b7ce53b6SStefano Zampini . newmat - the matrix in AIJ format 397b7ce53b6SStefano Zampini 398b7ce53b6SStefano Zampini Level: developer 399b7ce53b6SStefano Zampini 400eb82efa4SStefano Zampini Notes: mat and *newmat cannot be the same object when MAT_REUSE_MATRIX is requested. 401b7ce53b6SStefano Zampini 402b7ce53b6SStefano Zampini .seealso: MATIS 403b7ce53b6SStefano Zampini @*/ 404b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ(Mat mat, MatReuse reuse, Mat *newmat) 405b7ce53b6SStefano Zampini { 406b7ce53b6SStefano Zampini PetscErrorCode ierr; 407b7ce53b6SStefano Zampini 408b7ce53b6SStefano Zampini PetscFunctionBegin; 409b7ce53b6SStefano Zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 410b7ce53b6SStefano Zampini PetscValidLogicalCollectiveEnum(mat,reuse,2); 411b7ce53b6SStefano Zampini PetscValidPointer(newmat,3); 412b7ce53b6SStefano Zampini if (reuse != MAT_INITIAL_MATRIX) { 413b7ce53b6SStefano Zampini PetscValidHeaderSpecific(*newmat,MAT_CLASSID,3); 414b7ce53b6SStefano Zampini PetscCheckSameComm(mat,1,*newmat,3); 4156c4ed002SBarry Smith if (mat == *newmat) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse the same matrix"); 416b7ce53b6SStefano Zampini } 417b7ce53b6SStefano Zampini ierr = PetscUseMethod(mat,"MatISGetMPIXAIJ_C",(Mat,MatReuse,Mat*),(mat,reuse,newmat));CHKERRQ(ierr); 418b7ce53b6SStefano Zampini PetscFunctionReturn(0); 419b7ce53b6SStefano Zampini } 420b7ce53b6SStefano Zampini 421b7ce53b6SStefano Zampini #undef __FUNCT__ 422ad6194a2SStefano Zampini #define __FUNCT__ "MatDuplicate_IS" 423ad6194a2SStefano Zampini PetscErrorCode MatDuplicate_IS(Mat mat,MatDuplicateOption op,Mat *newmat) 424ad6194a2SStefano Zampini { 425ad6194a2SStefano Zampini PetscErrorCode ierr; 426ad6194a2SStefano Zampini Mat_IS *matis = (Mat_IS*)(mat->data); 427ad6194a2SStefano Zampini PetscInt bs,m,n,M,N; 428ad6194a2SStefano Zampini Mat B,localmat; 429ad6194a2SStefano Zampini 430ad6194a2SStefano Zampini PetscFunctionBegin; 431ad6194a2SStefano Zampini ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 432ad6194a2SStefano Zampini ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 433ad6194a2SStefano Zampini ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 434e176bc59SStefano Zampini ierr = MatCreateIS(PetscObjectComm((PetscObject)mat),bs,m,n,M,N,mat->rmap->mapping,mat->cmap->mapping,&B);CHKERRQ(ierr); 435ad6194a2SStefano Zampini ierr = MatDuplicate(matis->A,op,&localmat);CHKERRQ(ierr); 436ad6194a2SStefano Zampini ierr = MatISSetLocalMat(B,localmat);CHKERRQ(ierr); 437b3317aa8SStefano Zampini ierr = MatDestroy(&localmat);CHKERRQ(ierr); 438ad6194a2SStefano Zampini ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 439ad6194a2SStefano Zampini ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 440ad6194a2SStefano Zampini *newmat = B; 441ad6194a2SStefano Zampini PetscFunctionReturn(0); 442ad6194a2SStefano Zampini } 443ad6194a2SStefano Zampini 444ad6194a2SStefano Zampini #undef __FUNCT__ 44569796d55SStefano Zampini #define __FUNCT__ "MatIsHermitian_IS" 44669796d55SStefano Zampini PetscErrorCode MatIsHermitian_IS(Mat A,PetscReal tol,PetscBool *flg) 44769796d55SStefano Zampini { 44869796d55SStefano Zampini PetscErrorCode ierr; 44969796d55SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 45069796d55SStefano Zampini PetscBool local_sym; 45169796d55SStefano Zampini 45269796d55SStefano Zampini PetscFunctionBegin; 45369796d55SStefano Zampini ierr = MatIsHermitian(matis->A,tol,&local_sym);CHKERRQ(ierr); 454b2566f29SBarry Smith ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 45569796d55SStefano Zampini PetscFunctionReturn(0); 45669796d55SStefano Zampini } 45769796d55SStefano Zampini 45869796d55SStefano Zampini #undef __FUNCT__ 45969796d55SStefano Zampini #define __FUNCT__ "MatIsSymmetric_IS" 46069796d55SStefano Zampini PetscErrorCode MatIsSymmetric_IS(Mat A,PetscReal tol,PetscBool *flg) 46169796d55SStefano Zampini { 46269796d55SStefano Zampini PetscErrorCode ierr; 46369796d55SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 46469796d55SStefano Zampini PetscBool local_sym; 46569796d55SStefano Zampini 46669796d55SStefano Zampini PetscFunctionBegin; 46769796d55SStefano Zampini ierr = MatIsSymmetric(matis->A,tol,&local_sym);CHKERRQ(ierr); 468b2566f29SBarry Smith ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 46969796d55SStefano Zampini PetscFunctionReturn(0); 47069796d55SStefano Zampini } 47169796d55SStefano Zampini 47269796d55SStefano Zampini #undef __FUNCT__ 473b4319ba4SBarry Smith #define __FUNCT__ "MatDestroy_IS" 474dfbe8321SBarry Smith PetscErrorCode MatDestroy_IS(Mat A) 475b4319ba4SBarry Smith { 476dfbe8321SBarry Smith PetscErrorCode ierr; 477b4319ba4SBarry Smith Mat_IS *b = (Mat_IS*)A->data; 478b4319ba4SBarry Smith 479b4319ba4SBarry Smith PetscFunctionBegin; 4806bf464f9SBarry Smith ierr = MatDestroy(&b->A);CHKERRQ(ierr); 481e176bc59SStefano Zampini ierr = VecScatterDestroy(&b->cctx);CHKERRQ(ierr); 482e176bc59SStefano Zampini ierr = VecScatterDestroy(&b->rctx);CHKERRQ(ierr); 4836bf464f9SBarry Smith ierr = VecDestroy(&b->x);CHKERRQ(ierr); 4846bf464f9SBarry Smith ierr = VecDestroy(&b->y);CHKERRQ(ierr); 48528f4e0baSStefano Zampini ierr = PetscSFDestroy(&b->sf);CHKERRQ(ierr); 48628f4e0baSStefano Zampini ierr = PetscFree2(b->sf_rootdata,b->sf_leafdata);CHKERRQ(ierr); 487bf0cc555SLisandro Dalcin ierr = PetscFree(A->data);CHKERRQ(ierr); 488dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr); 489bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",NULL);CHKERRQ(ierr); 490b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",NULL);CHKERRQ(ierr); 491b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",NULL);CHKERRQ(ierr); 4922e1947a5SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",NULL);CHKERRQ(ierr); 493b4319ba4SBarry Smith PetscFunctionReturn(0); 494b4319ba4SBarry Smith } 495b4319ba4SBarry Smith 496b4319ba4SBarry Smith #undef __FUNCT__ 497b4319ba4SBarry Smith #define __FUNCT__ "MatMult_IS" 498dfbe8321SBarry Smith PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y) 499b4319ba4SBarry Smith { 500dfbe8321SBarry Smith PetscErrorCode ierr; 501b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 502b4319ba4SBarry Smith PetscScalar zero = 0.0; 503b4319ba4SBarry Smith 504b4319ba4SBarry Smith PetscFunctionBegin; 505b4319ba4SBarry Smith /* scatter the global vector x into the local work vector */ 506e176bc59SStefano Zampini ierr = VecScatterBegin(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 507e176bc59SStefano Zampini ierr = VecScatterEnd(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 508b4319ba4SBarry Smith 509b4319ba4SBarry Smith /* multiply the local matrix */ 510b4319ba4SBarry Smith ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr); 511b4319ba4SBarry Smith 512b4319ba4SBarry Smith /* scatter product back into global memory */ 5132dcb1b2aSMatthew Knepley ierr = VecSet(y,zero);CHKERRQ(ierr); 514e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 515e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 516b4319ba4SBarry Smith PetscFunctionReturn(0); 517b4319ba4SBarry Smith } 518b4319ba4SBarry Smith 519b4319ba4SBarry Smith #undef __FUNCT__ 5202e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultAdd_IS" 521b7ce53b6SStefano Zampini PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 5222e74eeadSLisandro Dalcin { 523650997f4SStefano Zampini Vec temp_vec; 5242e74eeadSLisandro Dalcin PetscErrorCode ierr; 5252e74eeadSLisandro Dalcin 5262e74eeadSLisandro Dalcin PetscFunctionBegin; /* v3 = v2 + A * v1.*/ 527650997f4SStefano Zampini if (v3 != v2) { 528650997f4SStefano Zampini ierr = MatMult(A,v1,v3);CHKERRQ(ierr); 529650997f4SStefano Zampini ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 530650997f4SStefano Zampini } else { 531650997f4SStefano Zampini ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 532650997f4SStefano Zampini ierr = MatMult(A,v1,temp_vec);CHKERRQ(ierr); 533650997f4SStefano Zampini ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 534650997f4SStefano Zampini ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 535650997f4SStefano Zampini ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 536650997f4SStefano Zampini } 5372e74eeadSLisandro Dalcin PetscFunctionReturn(0); 5382e74eeadSLisandro Dalcin } 5392e74eeadSLisandro Dalcin 5402e74eeadSLisandro Dalcin #undef __FUNCT__ 5412e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTranspose_IS" 542e176bc59SStefano Zampini PetscErrorCode MatMultTranspose_IS(Mat A,Vec y,Vec x) 5432e74eeadSLisandro Dalcin { 5442e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 5452e74eeadSLisandro Dalcin PetscErrorCode ierr; 5462e74eeadSLisandro Dalcin 547e176bc59SStefano Zampini PetscFunctionBegin; 5482e74eeadSLisandro Dalcin /* scatter the global vector x into the local work vector */ 549e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 550e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5512e74eeadSLisandro Dalcin 5522e74eeadSLisandro Dalcin /* multiply the local matrix */ 553e176bc59SStefano Zampini ierr = MatMultTranspose(is->A,is->y,is->x);CHKERRQ(ierr); 5542e74eeadSLisandro Dalcin 5552e74eeadSLisandro Dalcin /* scatter product back into global vector */ 556e176bc59SStefano Zampini ierr = VecSet(x,0);CHKERRQ(ierr); 557e176bc59SStefano Zampini ierr = VecScatterBegin(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 558e176bc59SStefano Zampini ierr = VecScatterEnd(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 5592e74eeadSLisandro Dalcin PetscFunctionReturn(0); 5602e74eeadSLisandro Dalcin } 5612e74eeadSLisandro Dalcin 5622e74eeadSLisandro Dalcin #undef __FUNCT__ 5632e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTransposeAdd_IS" 5642e74eeadSLisandro Dalcin PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 5652e74eeadSLisandro Dalcin { 566650997f4SStefano Zampini Vec temp_vec; 5672e74eeadSLisandro Dalcin PetscErrorCode ierr; 5682e74eeadSLisandro Dalcin 5692e74eeadSLisandro Dalcin PetscFunctionBegin; /* v3 = v2 + A' * v1.*/ 570650997f4SStefano Zampini if (v3 != v2) { 571650997f4SStefano Zampini ierr = MatMultTranspose(A,v1,v3);CHKERRQ(ierr); 572650997f4SStefano Zampini ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 573650997f4SStefano Zampini } else { 574650997f4SStefano Zampini ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 575650997f4SStefano Zampini ierr = MatMultTranspose(A,v1,temp_vec);CHKERRQ(ierr); 576650997f4SStefano Zampini ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 577650997f4SStefano Zampini ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 578650997f4SStefano Zampini ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 579650997f4SStefano Zampini } 5802e74eeadSLisandro Dalcin PetscFunctionReturn(0); 5812e74eeadSLisandro Dalcin } 5822e74eeadSLisandro Dalcin 5832e74eeadSLisandro Dalcin #undef __FUNCT__ 584b4319ba4SBarry Smith #define __FUNCT__ "MatView_IS" 585dfbe8321SBarry Smith PetscErrorCode MatView_IS(Mat A,PetscViewer viewer) 586b4319ba4SBarry Smith { 587b4319ba4SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 588dfbe8321SBarry Smith PetscErrorCode ierr; 589b4319ba4SBarry Smith PetscViewer sviewer; 590b4319ba4SBarry Smith 591b4319ba4SBarry Smith PetscFunctionBegin; 5923f08860eSBarry Smith ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 593b4319ba4SBarry Smith ierr = MatView(a->A,sviewer);CHKERRQ(ierr); 5943f08860eSBarry Smith ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 595b4319ba4SBarry Smith PetscFunctionReturn(0); 596b4319ba4SBarry Smith } 597b4319ba4SBarry Smith 598b4319ba4SBarry Smith #undef __FUNCT__ 599b4319ba4SBarry Smith #define __FUNCT__ "MatSetLocalToGlobalMapping_IS" 600784ac674SJed Brown PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping) 601b4319ba4SBarry Smith { 602dfbe8321SBarry Smith PetscErrorCode ierr; 603e176bc59SStefano Zampini PetscInt nr,rbs,nc,cbs; 604b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 605b4319ba4SBarry Smith IS from,to; 606e176bc59SStefano Zampini Vec cglobal,rglobal; 607b4319ba4SBarry Smith 608b4319ba4SBarry Smith PetscFunctionBegin; 609784ac674SJed Brown PetscCheckSameComm(A,1,rmapping,2); 610e176bc59SStefano Zampini PetscCheckSameComm(A,1,cmapping,3); 6113bbff08aSStefano Zampini /* Destroy any previous data */ 61270cf5478SStefano Zampini ierr = VecDestroy(&is->x);CHKERRQ(ierr); 61370cf5478SStefano Zampini ierr = VecDestroy(&is->y);CHKERRQ(ierr); 614e176bc59SStefano Zampini ierr = VecScatterDestroy(&is->rctx);CHKERRQ(ierr); 615e176bc59SStefano Zampini ierr = VecScatterDestroy(&is->cctx);CHKERRQ(ierr); 6161c47cb0fSStefano Zampini ierr = MatDestroy(&is->A);CHKERRQ(ierr); 61728f4e0baSStefano Zampini ierr = PetscSFDestroy(&is->sf);CHKERRQ(ierr); 61828f4e0baSStefano Zampini ierr = PetscFree2(is->sf_rootdata,is->sf_leafdata);CHKERRQ(ierr); 6193bbff08aSStefano Zampini 6203bbff08aSStefano Zampini /* Setup Layout and set local to global maps */ 621fc27028aSStefano Zampini ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 622fc27028aSStefano Zampini ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 623fa7f1dd8SStefano Zampini ierr = PetscLayoutSetISLocalToGlobalMapping(A->rmap,rmapping);CHKERRQ(ierr); 624fa7f1dd8SStefano Zampini ierr = PetscLayoutSetISLocalToGlobalMapping(A->cmap,cmapping);CHKERRQ(ierr); 625b4319ba4SBarry Smith 626b4319ba4SBarry Smith /* Create the local matrix A */ 627e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(rmapping,&nr);CHKERRQ(ierr); 628e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetBlockSize(rmapping,&rbs);CHKERRQ(ierr); 629e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(cmapping,&nc);CHKERRQ(ierr); 630e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetBlockSize(cmapping,&cbs);CHKERRQ(ierr); 631f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr); 632e176bc59SStefano Zampini ierr = MatSetType(is->A,MATAIJ);CHKERRQ(ierr); 633e176bc59SStefano Zampini ierr = MatSetSizes(is->A,nr,nc,nr,nc);CHKERRQ(ierr); 634e176bc59SStefano Zampini ierr = MatSetBlockSizes(is->A,rbs,cbs);CHKERRQ(ierr); 635ff130e51SJed Brown ierr = MatSetOptionsPrefix(is->A,((PetscObject)A)->prefix);CHKERRQ(ierr); 636ff130e51SJed Brown ierr = MatAppendOptionsPrefix(is->A,"is_");CHKERRQ(ierr); 637b4319ba4SBarry Smith ierr = MatSetFromOptions(is->A);CHKERRQ(ierr); 638b4319ba4SBarry Smith 639b4319ba4SBarry Smith /* Create the local work vectors */ 6403bbff08aSStefano Zampini ierr = MatCreateVecs(is->A,&is->x,&is->y);CHKERRQ(ierr); 641b4319ba4SBarry Smith 642e176bc59SStefano Zampini /* setup the global to local scatters */ 643e176bc59SStefano Zampini ierr = MatCreateVecs(A,&cglobal,&rglobal);CHKERRQ(ierr); 644e176bc59SStefano Zampini ierr = ISCreateStride(PETSC_COMM_SELF,nr,0,1,&to);CHKERRQ(ierr); 645784ac674SJed Brown ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr); 646e176bc59SStefano Zampini ierr = VecScatterCreate(rglobal,from,is->y,to,&is->rctx);CHKERRQ(ierr); 647e176bc59SStefano Zampini if (rmapping != cmapping) { 648e176bc59SStefano Zampini ierr = ISDestroy(&to);CHKERRQ(ierr); 649e176bc59SStefano Zampini ierr = ISDestroy(&from);CHKERRQ(ierr); 650e176bc59SStefano Zampini ierr = ISCreateStride(PETSC_COMM_SELF,nc,0,1,&to);CHKERRQ(ierr); 651e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingApplyIS(cmapping,to,&from);CHKERRQ(ierr); 652e176bc59SStefano Zampini ierr = VecScatterCreate(cglobal,from,is->x,to,&is->cctx);CHKERRQ(ierr); 653e176bc59SStefano Zampini } else { 654e176bc59SStefano Zampini ierr = PetscObjectReference((PetscObject)is->rctx);CHKERRQ(ierr); 655e176bc59SStefano Zampini is->cctx = is->rctx; 656e176bc59SStefano Zampini } 657e176bc59SStefano Zampini ierr = VecDestroy(&rglobal);CHKERRQ(ierr); 658e176bc59SStefano Zampini ierr = VecDestroy(&cglobal);CHKERRQ(ierr); 6596bf464f9SBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 6606bf464f9SBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 661b4319ba4SBarry Smith PetscFunctionReturn(0); 662b4319ba4SBarry Smith } 663b4319ba4SBarry Smith 6642e74eeadSLisandro Dalcin #undef __FUNCT__ 6652e74eeadSLisandro Dalcin #define __FUNCT__ "MatSetValues_IS" 6662e74eeadSLisandro Dalcin PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv) 6672e74eeadSLisandro Dalcin { 6682e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)mat->data; 6692e74eeadSLisandro Dalcin PetscInt rows_l[2048],cols_l[2048]; 6702e74eeadSLisandro Dalcin PetscErrorCode ierr; 6712e74eeadSLisandro Dalcin 6722e74eeadSLisandro Dalcin PetscFunctionBegin; 6732e74eeadSLisandro Dalcin #if defined(PETSC_USE_DEBUG) 674f23aa3ddSBarry 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); 6752e74eeadSLisandro Dalcin #endif 6763bbff08aSStefano Zampini ierr = ISG2LMapApply(mat->rmap->mapping,m,rows,rows_l);CHKERRQ(ierr); 6773bbff08aSStefano Zampini ierr = ISG2LMapApply(mat->cmap->mapping,n,cols,cols_l);CHKERRQ(ierr); 6782e74eeadSLisandro Dalcin ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 6792e74eeadSLisandro Dalcin PetscFunctionReturn(0); 6802e74eeadSLisandro Dalcin } 6812e74eeadSLisandro Dalcin 682b4319ba4SBarry Smith #undef __FUNCT__ 683b4319ba4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_IS" 68413f74950SBarry Smith PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 685b4319ba4SBarry Smith { 686dfbe8321SBarry Smith PetscErrorCode ierr; 687b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 688b4319ba4SBarry Smith 689b4319ba4SBarry Smith PetscFunctionBegin; 690b4319ba4SBarry Smith ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 691b4319ba4SBarry Smith PetscFunctionReturn(0); 692b4319ba4SBarry Smith } 693b4319ba4SBarry Smith 694b4319ba4SBarry Smith #undef __FUNCT__ 695f0006bf2SLisandro Dalcin #define __FUNCT__ "MatSetValuesBlockedLocal_IS" 696f0006bf2SLisandro Dalcin PetscErrorCode MatSetValuesBlockedLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 697f0006bf2SLisandro Dalcin { 698f0006bf2SLisandro Dalcin PetscErrorCode ierr; 699f0006bf2SLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 700f0006bf2SLisandro Dalcin 701f0006bf2SLisandro Dalcin PetscFunctionBegin; 702f0006bf2SLisandro Dalcin ierr = MatSetValuesBlocked(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 703f0006bf2SLisandro Dalcin PetscFunctionReturn(0); 704f0006bf2SLisandro Dalcin } 705f0006bf2SLisandro Dalcin 706f0006bf2SLisandro Dalcin #undef __FUNCT__ 7072e74eeadSLisandro Dalcin #define __FUNCT__ "MatZeroRows_IS" 7082b40b63fSBarry Smith PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 7092e74eeadSLisandro Dalcin { 7100298fd71SBarry Smith PetscInt n_l = 0, *rows_l = NULL; 7112e74eeadSLisandro Dalcin PetscErrorCode ierr; 7122e74eeadSLisandro Dalcin 7132e74eeadSLisandro Dalcin PetscFunctionBegin; 714ce94432eSBarry Smith if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support"); 7152e74eeadSLisandro Dalcin if (n) { 716785e854fSJed Brown ierr = PetscMalloc1(n,&rows_l);CHKERRQ(ierr); 7173bbff08aSStefano Zampini ierr = ISGlobalToLocalMappingApply(A->rmap->mapping,IS_GTOLM_DROP,n,rows,&n_l,rows_l);CHKERRQ(ierr); 7182e74eeadSLisandro Dalcin } 7192b40b63fSBarry Smith ierr = MatZeroRowsLocal(A,n_l,rows_l,diag,x,b);CHKERRQ(ierr); 720c31cb41cSBarry Smith ierr = PetscFree(rows_l);CHKERRQ(ierr); 7212e74eeadSLisandro Dalcin PetscFunctionReturn(0); 7222e74eeadSLisandro Dalcin } 7232e74eeadSLisandro Dalcin 7242e74eeadSLisandro Dalcin #undef __FUNCT__ 725b4319ba4SBarry Smith #define __FUNCT__ "MatZeroRowsLocal_IS" 7262b40b63fSBarry Smith PetscErrorCode MatZeroRowsLocal_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 727b4319ba4SBarry Smith { 728b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 729dfbe8321SBarry Smith PetscErrorCode ierr; 730f4df32b1SMatthew Knepley PetscInt i; 731b4319ba4SBarry Smith PetscScalar *array; 732b4319ba4SBarry Smith 733b4319ba4SBarry Smith PetscFunctionBegin; 734ce94432eSBarry Smith if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support"); 735b4319ba4SBarry Smith { 736b4319ba4SBarry Smith /* 737b4319ba4SBarry Smith Set up is->x as a "counting vector". This is in order to MatMult_IS 738b4319ba4SBarry Smith work properly in the interface nodes. 739b4319ba4SBarry Smith */ 740b4319ba4SBarry Smith Vec counter; 741e176bc59SStefano Zampini ierr = MatCreateVecs(A,NULL,&counter);CHKERRQ(ierr); 742e176bc59SStefano Zampini ierr = VecSet(counter,0.);CHKERRQ(ierr); 743e176bc59SStefano Zampini ierr = VecSet(is->y,1.);CHKERRQ(ierr); 744e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,is->y,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 745e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,is->y,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 746e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,counter,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 747e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,counter,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 7486bf464f9SBarry Smith ierr = VecDestroy(&counter);CHKERRQ(ierr); 749b4319ba4SBarry Smith } 750958c9bccSBarry Smith if (!n) { 751b4319ba4SBarry Smith is->pure_neumann = PETSC_TRUE; 752b4319ba4SBarry Smith } else { 753b4319ba4SBarry Smith is->pure_neumann = PETSC_FALSE; 7542205254eSKarl Rupp 755e176bc59SStefano Zampini ierr = VecGetArray(is->y,&array);CHKERRQ(ierr); 7562b40b63fSBarry Smith ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr); 757b4319ba4SBarry Smith for (i=0; i<n; i++) { 758f4df32b1SMatthew Knepley ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr); 759b4319ba4SBarry Smith } 760b4319ba4SBarry Smith ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 761b4319ba4SBarry Smith ierr = MatAssemblyEnd(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 762e176bc59SStefano Zampini ierr = VecRestoreArray(is->y,&array);CHKERRQ(ierr); 763b4319ba4SBarry Smith } 764b4319ba4SBarry Smith PetscFunctionReturn(0); 765b4319ba4SBarry Smith } 766b4319ba4SBarry Smith 767b4319ba4SBarry Smith #undef __FUNCT__ 768b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyBegin_IS" 769dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type) 770b4319ba4SBarry Smith { 771b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 772dfbe8321SBarry Smith PetscErrorCode ierr; 773b4319ba4SBarry Smith 774b4319ba4SBarry Smith PetscFunctionBegin; 775b4319ba4SBarry Smith ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr); 776b4319ba4SBarry Smith PetscFunctionReturn(0); 777b4319ba4SBarry Smith } 778b4319ba4SBarry Smith 779b4319ba4SBarry Smith #undef __FUNCT__ 780b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_IS" 781dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type) 782b4319ba4SBarry Smith { 783b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 784dfbe8321SBarry Smith PetscErrorCode ierr; 785b4319ba4SBarry Smith 786b4319ba4SBarry Smith PetscFunctionBegin; 787b4319ba4SBarry Smith ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr); 788b4319ba4SBarry Smith PetscFunctionReturn(0); 789b4319ba4SBarry Smith } 790b4319ba4SBarry Smith 791b4319ba4SBarry Smith #undef __FUNCT__ 792b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat_IS" 7937087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local) 794b4319ba4SBarry Smith { 795b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)mat->data; 796b4319ba4SBarry Smith 797b4319ba4SBarry Smith PetscFunctionBegin; 798b4319ba4SBarry Smith *local = is->A; 799b4319ba4SBarry Smith PetscFunctionReturn(0); 800b4319ba4SBarry Smith } 801b4319ba4SBarry Smith 802b4319ba4SBarry Smith #undef __FUNCT__ 803b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat" 804b4319ba4SBarry Smith /*@ 805b4319ba4SBarry Smith MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix. 806b4319ba4SBarry Smith 807b4319ba4SBarry Smith Input Parameter: 808b4319ba4SBarry Smith . mat - the matrix 809b4319ba4SBarry Smith 810b4319ba4SBarry Smith Output Parameter: 811eb82efa4SStefano Zampini . local - the local matrix 812b4319ba4SBarry Smith 813b4319ba4SBarry Smith Level: advanced 814b4319ba4SBarry Smith 815b4319ba4SBarry Smith Notes: 816b4319ba4SBarry Smith This can be called if you have precomputed the nonzero structure of the 817b4319ba4SBarry Smith matrix and want to provide it to the inner matrix object to improve the performance 818b4319ba4SBarry Smith of the MatSetValues() operation. 819b4319ba4SBarry Smith 82096a6f129SJed Brown This function does not increase the reference count for the local Mat. Do not destroy it and do not attempt to use 82196a6f129SJed Brown your reference after destroying the parent mat. 82296a6f129SJed Brown 823b4319ba4SBarry Smith .seealso: MATIS 824b4319ba4SBarry Smith @*/ 8257087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local) 826b4319ba4SBarry Smith { 8274ac538c5SBarry Smith PetscErrorCode ierr; 828b4319ba4SBarry Smith 829b4319ba4SBarry Smith PetscFunctionBegin; 8300700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 831b4319ba4SBarry Smith PetscValidPointer(local,2); 8324ac538c5SBarry Smith ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat*),(mat,local));CHKERRQ(ierr); 833b4319ba4SBarry Smith PetscFunctionReturn(0); 834b4319ba4SBarry Smith } 835b4319ba4SBarry Smith 8363b03a366Sstefano_zampini #undef __FUNCT__ 8373b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat_IS" 8383b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat_IS(Mat mat,Mat local) 8393b03a366Sstefano_zampini { 8403b03a366Sstefano_zampini Mat_IS *is = (Mat_IS*)mat->data; 8413b03a366Sstefano_zampini PetscInt nrows,ncols,orows,ocols; 8423b03a366Sstefano_zampini PetscErrorCode ierr; 8433b03a366Sstefano_zampini 8443b03a366Sstefano_zampini PetscFunctionBegin; 8454e4c7dbeSStefano Zampini if (is->A) { 8463b03a366Sstefano_zampini ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr); 8473b03a366Sstefano_zampini ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr); 848f23aa3ddSBarry 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); 8494e4c7dbeSStefano Zampini } 8503b03a366Sstefano_zampini ierr = PetscObjectReference((PetscObject)local);CHKERRQ(ierr); 8513b03a366Sstefano_zampini ierr = MatDestroy(&is->A);CHKERRQ(ierr); 8523b03a366Sstefano_zampini is->A = local; 8533b03a366Sstefano_zampini PetscFunctionReturn(0); 8543b03a366Sstefano_zampini } 8553b03a366Sstefano_zampini 8563b03a366Sstefano_zampini #undef __FUNCT__ 8573b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat" 8583b03a366Sstefano_zampini /*@ 859eb82efa4SStefano Zampini MatISSetLocalMat - Replace the local matrix stored inside a MATIS object. 8603b03a366Sstefano_zampini 8613b03a366Sstefano_zampini Input Parameter: 8623b03a366Sstefano_zampini . mat - the matrix 863eb82efa4SStefano Zampini . local - the local matrix 8643b03a366Sstefano_zampini 8653b03a366Sstefano_zampini Output Parameter: 8663b03a366Sstefano_zampini 8673b03a366Sstefano_zampini Level: advanced 8683b03a366Sstefano_zampini 8693b03a366Sstefano_zampini Notes: 8703b03a366Sstefano_zampini This can be called if you have precomputed the local matrix and 8713b03a366Sstefano_zampini want to provide it to the matrix object MATIS. 8723b03a366Sstefano_zampini 8733b03a366Sstefano_zampini .seealso: MATIS 8743b03a366Sstefano_zampini @*/ 8753b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat(Mat mat,Mat local) 8763b03a366Sstefano_zampini { 8773b03a366Sstefano_zampini PetscErrorCode ierr; 8783b03a366Sstefano_zampini 8793b03a366Sstefano_zampini PetscFunctionBegin; 8803b03a366Sstefano_zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 881b7ce53b6SStefano Zampini PetscValidHeaderSpecific(local,MAT_CLASSID,2); 8823b03a366Sstefano_zampini ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr); 8833b03a366Sstefano_zampini PetscFunctionReturn(0); 8843b03a366Sstefano_zampini } 8853b03a366Sstefano_zampini 8866726f965SBarry Smith #undef __FUNCT__ 8876726f965SBarry Smith #define __FUNCT__ "MatZeroEntries_IS" 8886726f965SBarry Smith PetscErrorCode MatZeroEntries_IS(Mat A) 8896726f965SBarry Smith { 8906726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 8916726f965SBarry Smith PetscErrorCode ierr; 8926726f965SBarry Smith 8936726f965SBarry Smith PetscFunctionBegin; 8946726f965SBarry Smith ierr = MatZeroEntries(a->A);CHKERRQ(ierr); 8956726f965SBarry Smith PetscFunctionReturn(0); 8966726f965SBarry Smith } 8976726f965SBarry Smith 8986726f965SBarry Smith #undef __FUNCT__ 8992e74eeadSLisandro Dalcin #define __FUNCT__ "MatScale_IS" 9002e74eeadSLisandro Dalcin PetscErrorCode MatScale_IS(Mat A,PetscScalar a) 9012e74eeadSLisandro Dalcin { 9022e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 9032e74eeadSLisandro Dalcin PetscErrorCode ierr; 9042e74eeadSLisandro Dalcin 9052e74eeadSLisandro Dalcin PetscFunctionBegin; 9062e74eeadSLisandro Dalcin ierr = MatScale(is->A,a);CHKERRQ(ierr); 9072e74eeadSLisandro Dalcin PetscFunctionReturn(0); 9082e74eeadSLisandro Dalcin } 9092e74eeadSLisandro Dalcin 9102e74eeadSLisandro Dalcin #undef __FUNCT__ 9112e74eeadSLisandro Dalcin #define __FUNCT__ "MatGetDiagonal_IS" 9122e74eeadSLisandro Dalcin PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v) 9132e74eeadSLisandro Dalcin { 9142e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 9152e74eeadSLisandro Dalcin PetscErrorCode ierr; 9162e74eeadSLisandro Dalcin 9172e74eeadSLisandro Dalcin PetscFunctionBegin; 9182e74eeadSLisandro Dalcin /* get diagonal of the local matrix */ 919e176bc59SStefano Zampini ierr = MatGetDiagonal(is->A,is->y);CHKERRQ(ierr); 9202e74eeadSLisandro Dalcin 9212e74eeadSLisandro Dalcin /* scatter diagonal back into global vector */ 9222e74eeadSLisandro Dalcin ierr = VecSet(v,0);CHKERRQ(ierr); 923e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 924e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 9252e74eeadSLisandro Dalcin PetscFunctionReturn(0); 9262e74eeadSLisandro Dalcin } 9272e74eeadSLisandro Dalcin 9282e74eeadSLisandro Dalcin #undef __FUNCT__ 9296726f965SBarry Smith #define __FUNCT__ "MatSetOption_IS" 930ace3abfcSBarry Smith PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg) 9316726f965SBarry Smith { 9326726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 9336726f965SBarry Smith PetscErrorCode ierr; 9346726f965SBarry Smith 9356726f965SBarry Smith PetscFunctionBegin; 9364e0d8c25SBarry Smith ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr); 9376726f965SBarry Smith PetscFunctionReturn(0); 9386726f965SBarry Smith } 9396726f965SBarry Smith 940284134d9SBarry Smith #undef __FUNCT__ 941284134d9SBarry Smith #define __FUNCT__ "MatCreateIS" 942284134d9SBarry Smith /*@ 943284134d9SBarry Smith MatCreateIS - Creates a "process" unassmembled matrix, it is assembled on each 944284134d9SBarry Smith process but not across processes. 945284134d9SBarry Smith 946284134d9SBarry Smith Input Parameters: 947284134d9SBarry Smith + comm - MPI communicator that will share the matrix 948e176bc59SStefano Zampini . bs - block size of the matrix 949df3898eeSBarry Smith . m,n,M,N - local and/or global sizes of the left and right vector used in matrix vector products 950e176bc59SStefano Zampini . rmap - local to global map for rows 951e176bc59SStefano Zampini - cmap - local to global map for cols 952284134d9SBarry Smith 953284134d9SBarry Smith Output Parameter: 954284134d9SBarry Smith . A - the resulting matrix 955284134d9SBarry Smith 9568e6c10adSSatish Balay Level: advanced 9578e6c10adSSatish Balay 958284134d9SBarry Smith Notes: See MATIS for more details 9598cda6cd7SBarry Smith m and n are NOT related to the size of the map, they are the size of the part of the vector owned 960e176bc59SStefano Zampini by that process. The sizes of rmap and cmap define the size of the local matrices. 961e176bc59SStefano Zampini If either rmap or cmap are NULL, than the matrix is assumed to be square 962284134d9SBarry Smith 963284134d9SBarry Smith .seealso: MATIS, MatSetLocalToGlobalMapping() 964284134d9SBarry Smith @*/ 965e176bc59SStefano Zampini PetscErrorCode MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping rmap,ISLocalToGlobalMapping cmap,Mat *A) 966284134d9SBarry Smith { 967284134d9SBarry Smith PetscErrorCode ierr; 968284134d9SBarry Smith 969284134d9SBarry Smith PetscFunctionBegin; 970e176bc59SStefano Zampini if (!rmap && !cmap) SETERRQ(comm,PETSC_ERR_USER,"You need to provide at least one of the mapping"); 971284134d9SBarry Smith ierr = MatCreate(comm,A);CHKERRQ(ierr); 972284134d9SBarry Smith ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 9734f619741Sstefano_zampini ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr); 974284134d9SBarry Smith ierr = MatSetType(*A,MATIS);CHKERRQ(ierr); 9753b03a366Sstefano_zampini ierr = MatSetUp(*A);CHKERRQ(ierr); 976e176bc59SStefano Zampini if (rmap && cmap) { 977e176bc59SStefano Zampini ierr = MatSetLocalToGlobalMapping(*A,rmap,cmap);CHKERRQ(ierr); 978e176bc59SStefano Zampini } else if (!rmap) { 979e176bc59SStefano Zampini ierr = MatSetLocalToGlobalMapping(*A,cmap,cmap);CHKERRQ(ierr); 980e176bc59SStefano Zampini } else { 981e176bc59SStefano Zampini ierr = MatSetLocalToGlobalMapping(*A,rmap,rmap);CHKERRQ(ierr); 982e176bc59SStefano Zampini } 983284134d9SBarry Smith PetscFunctionReturn(0); 984284134d9SBarry Smith } 985284134d9SBarry Smith 986b4319ba4SBarry Smith /*MC 987eb82efa4SStefano Zampini MATIS - MATIS = "is" - A matrix type to be used for using the non-overlapping domain decomposition type preconditioners (e.g. PCBDDC). 988b4319ba4SBarry Smith This stores the matrices in globally unassembled form. Each processor 989b4319ba4SBarry Smith assembles only its local Neumann problem and the parallel matrix vector 990b4319ba4SBarry Smith product is handled "implicitly". 991b4319ba4SBarry Smith 992b4319ba4SBarry Smith Operations Provided: 9936726f965SBarry Smith + MatMult() 9942e74eeadSLisandro Dalcin . MatMultAdd() 9952e74eeadSLisandro Dalcin . MatMultTranspose() 9962e74eeadSLisandro Dalcin . MatMultTransposeAdd() 9976726f965SBarry Smith . MatZeroEntries() 9986726f965SBarry Smith . MatSetOption() 9992e74eeadSLisandro Dalcin . MatZeroRows() 10006726f965SBarry Smith . MatZeroRowsLocal() 10012e74eeadSLisandro Dalcin . MatSetValues() 10026726f965SBarry Smith . MatSetValuesLocal() 10032e74eeadSLisandro Dalcin . MatScale() 10042e74eeadSLisandro Dalcin . MatGetDiagonal() 10056726f965SBarry Smith - MatSetLocalToGlobalMapping() 1006b4319ba4SBarry Smith 1007b4319ba4SBarry Smith Options Database Keys: 1008b4319ba4SBarry Smith . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions() 1009b4319ba4SBarry Smith 1010b4319ba4SBarry Smith Notes: Options prefix for the inner matrix are given by -is_mat_xxx 1011b4319ba4SBarry Smith 1012b4319ba4SBarry Smith You must call MatSetLocalToGlobalMapping() before using this matrix type. 1013b4319ba4SBarry Smith 1014b4319ba4SBarry Smith You can do matrix preallocation on the local matrix after you obtain it with 1015eb82efa4SStefano Zampini MatISGetLocalMat(); otherwise, you could use MatISSetPreallocation() 1016b4319ba4SBarry Smith 1017b4319ba4SBarry Smith Level: advanced 1018b4319ba4SBarry Smith 1019eb82efa4SStefano Zampini .seealso: Mat, MatISGetLocalMat(), MatSetLocalToGlobalMapping(), MatISSetPreallocation(), PCBDDC 1020b4319ba4SBarry Smith 1021b4319ba4SBarry Smith M*/ 1022b4319ba4SBarry Smith 1023b4319ba4SBarry Smith #undef __FUNCT__ 1024b4319ba4SBarry Smith #define __FUNCT__ "MatCreate_IS" 10258cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_IS(Mat A) 1026b4319ba4SBarry Smith { 1027dfbe8321SBarry Smith PetscErrorCode ierr; 1028b4319ba4SBarry Smith Mat_IS *b; 1029b4319ba4SBarry Smith 1030b4319ba4SBarry Smith PetscFunctionBegin; 1031b00a9115SJed Brown ierr = PetscNewLog(A,&b);CHKERRQ(ierr); 1032b4319ba4SBarry Smith A->data = (void*)b; 1033b4319ba4SBarry Smith 1034e176bc59SStefano Zampini /* matrix ops */ 1035e176bc59SStefano Zampini ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1036b4319ba4SBarry Smith A->ops->mult = MatMult_IS; 10372e74eeadSLisandro Dalcin A->ops->multadd = MatMultAdd_IS; 10382e74eeadSLisandro Dalcin A->ops->multtranspose = MatMultTranspose_IS; 10392e74eeadSLisandro Dalcin A->ops->multtransposeadd = MatMultTransposeAdd_IS; 1040b4319ba4SBarry Smith A->ops->destroy = MatDestroy_IS; 1041b4319ba4SBarry Smith A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS; 10422e74eeadSLisandro Dalcin A->ops->setvalues = MatSetValues_IS; 1043b4319ba4SBarry Smith A->ops->setvalueslocal = MatSetValuesLocal_IS; 1044f0006bf2SLisandro Dalcin A->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_IS; 10452e74eeadSLisandro Dalcin A->ops->zerorows = MatZeroRows_IS; 1046b4319ba4SBarry Smith A->ops->zerorowslocal = MatZeroRowsLocal_IS; 1047b4319ba4SBarry Smith A->ops->assemblybegin = MatAssemblyBegin_IS; 1048b4319ba4SBarry Smith A->ops->assemblyend = MatAssemblyEnd_IS; 1049b4319ba4SBarry Smith A->ops->view = MatView_IS; 10506726f965SBarry Smith A->ops->zeroentries = MatZeroEntries_IS; 10512e74eeadSLisandro Dalcin A->ops->scale = MatScale_IS; 10522e74eeadSLisandro Dalcin A->ops->getdiagonal = MatGetDiagonal_IS; 10536726f965SBarry Smith A->ops->setoption = MatSetOption_IS; 105469796d55SStefano Zampini A->ops->ishermitian = MatIsHermitian_IS; 105569796d55SStefano Zampini A->ops->issymmetric = MatIsSymmetric_IS; 1056ad6194a2SStefano Zampini A->ops->duplicate = MatDuplicate_IS; 1057b4319ba4SBarry Smith 1058b7ce53b6SStefano Zampini /* special MATIS functions */ 1059bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",MatISGetLocalMat_IS);CHKERRQ(ierr); 1060bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",MatISSetLocalMat_IS);CHKERRQ(ierr); 1061b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",MatISGetMPIXAIJ_IS);CHKERRQ(ierr); 10622e1947a5SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",MatISSetPreallocation_IS);CHKERRQ(ierr); 106317667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr); 1064b4319ba4SBarry Smith PetscFunctionReturn(0); 1065b4319ba4SBarry Smith } 1066