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 743c212e90SHong Zhang .seealso: MatCreate(), MatCreateIS(), MatMPIAIJSetPreallocation(), MatISGetLocalMat(), MATIS 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) { 152ecf5a873SStefano Zampini ierr = ISLocalToGlobalMappingGetIndices(A->rmap->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 } 1773927de2eSStefano Zampini 1783927de2eSStefano Zampini /* 1793927de2eSStefano Zampini my_dnz and my_onz contains exact contribution to preallocation from each local mat 1803927de2eSStefano Zampini then, they will be summed up properly. This way, preallocation is always sufficient 1813927de2eSStefano Zampini */ 1823927de2eSStefano Zampini ierr = PetscCalloc2(local_rows,&my_dnz,local_rows,&my_onz);CHKERRQ(ierr); 1833927de2eSStefano Zampini /* preallocation as a MATAIJ */ 1843927de2eSStefano Zampini if (isdense) { /* special case for dense local matrices */ 1853927de2eSStefano Zampini for (i=0;i<local_rows;i++) { 186ecf5a873SStefano Zampini PetscInt index_row = global_indices_r[i]; 1873927de2eSStefano Zampini for (j=i;j<local_rows;j++) { 1883927de2eSStefano Zampini PetscInt owner = row_ownership[index_row]; 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 /* same as before, interchanging rows and cols */ 1963927de2eSStefano Zampini if (i != j) { 1973927de2eSStefano Zampini owner = row_ownership[index_col]; 1983927de2eSStefano Zampini if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) { 1993927de2eSStefano Zampini my_dnz[j] += 1; 2003927de2eSStefano Zampini } else { 2013927de2eSStefano Zampini my_onz[j] += 1; 2023927de2eSStefano Zampini } 2033927de2eSStefano Zampini } 2043927de2eSStefano Zampini } 2053927de2eSStefano Zampini } 206ecf5a873SStefano Zampini } else { /* TODO: this could be optimized using MatGetRowIJ */ 2073927de2eSStefano Zampini for (i=0;i<local_rows;i++) { 2083927de2eSStefano Zampini const PetscInt *cols; 209ecf5a873SStefano Zampini PetscInt ncols,index_row = global_indices_r[i]; 2103927de2eSStefano Zampini ierr = MatGetRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr); 2113927de2eSStefano Zampini for (j=0;j<ncols;j++) { 2123927de2eSStefano Zampini PetscInt owner = row_ownership[index_row]; 213ecf5a873SStefano Zampini PetscInt index_col = global_indices_c[cols[j]]; 2143927de2eSStefano Zampini if (index_col > mat_ranges[owner]-1 && index_col < mat_ranges[owner+1] ) { /* diag block */ 2153927de2eSStefano Zampini my_dnz[i] += 1; 2163927de2eSStefano Zampini } else { /* offdiag block */ 2173927de2eSStefano Zampini my_onz[i] += 1; 2183927de2eSStefano Zampini } 2193927de2eSStefano Zampini /* same as before, interchanging rows and cols */ 220d9a9e74cSStefano Zampini if (issbaij && index_col != index_row) { 2213927de2eSStefano Zampini owner = row_ownership[index_col]; 2223927de2eSStefano Zampini if (index_row > mat_ranges[owner]-1 && index_row < mat_ranges[owner+1] ) { 223d9a9e74cSStefano Zampini my_dnz[cols[j]] += 1; 2243927de2eSStefano Zampini } else { 225d9a9e74cSStefano Zampini my_onz[cols[j]] += 1; 2263927de2eSStefano Zampini } 2273927de2eSStefano Zampini } 2283927de2eSStefano Zampini } 2293927de2eSStefano Zampini ierr = MatRestoreRow(matis->A,i,&ncols,&cols,NULL);CHKERRQ(ierr); 2303927de2eSStefano Zampini } 2313927de2eSStefano Zampini } 232ecf5a873SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(A->rmap->mapping,&global_indices_r);CHKERRQ(ierr); 233ecf5a873SStefano Zampini if (global_indices_c != global_indices_r) { 234ecf5a873SStefano Zampini ierr = ISLocalToGlobalMappingRestoreIndices(A->rmap->mapping,&global_indices_c);CHKERRQ(ierr); 235ecf5a873SStefano Zampini } 2363927de2eSStefano Zampini ierr = PetscFree(row_ownership);CHKERRQ(ierr); 237ecf5a873SStefano Zampini 238ecf5a873SStefano Zampini /* Reduce my_dnz and my_onz */ 2393927de2eSStefano Zampini if (maxreduce) { 2403927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr); 2413927de2eSStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_MAX);CHKERRQ(ierr); 2423927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr); 2433927de2eSStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_MAX);CHKERRQ(ierr); 2443927de2eSStefano Zampini } else { 2453927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr); 2463927de2eSStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_dnz,dnz,MPI_SUM);CHKERRQ(ierr); 2473927de2eSStefano Zampini ierr = PetscSFReduceBegin(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr); 2483927de2eSStefano Zampini ierr = PetscSFReduceEnd(matis->sf,MPIU_INT,my_onz,onz,MPI_SUM);CHKERRQ(ierr); 2493927de2eSStefano Zampini } 2503927de2eSStefano Zampini ierr = PetscFree2(my_dnz,my_onz);CHKERRQ(ierr); 2513927de2eSStefano Zampini 2523927de2eSStefano Zampini /* Resize preallocation if overestimated */ 2533927de2eSStefano Zampini for (i=0;i<lrows;i++) { 2543927de2eSStefano Zampini dnz[i] = PetscMin(dnz[i],lcols); 2553927de2eSStefano Zampini onz[i] = PetscMin(onz[i],cols-lcols); 2563927de2eSStefano Zampini } 2573927de2eSStefano Zampini /* set preallocation */ 2583927de2eSStefano Zampini ierr = MatMPIAIJSetPreallocation(B,0,dnz,0,onz);CHKERRQ(ierr); 2593927de2eSStefano Zampini for (i=0;i<lrows/bs;i++) { 2603927de2eSStefano Zampini dnz[i] = dnz[i*bs]/bs; 2613927de2eSStefano Zampini onz[i] = onz[i*bs]/bs; 2623927de2eSStefano Zampini } 2633927de2eSStefano Zampini ierr = MatMPIBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr); 2643927de2eSStefano Zampini ierr = MatMPISBAIJSetPreallocation(B,bs,0,dnz,0,onz);CHKERRQ(ierr); 2653927de2eSStefano Zampini ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 2663927de2eSStefano Zampini if (issbaij) { 2673927de2eSStefano Zampini ierr = MatRestoreRowUpperTriangular(matis->A);CHKERRQ(ierr); 2683927de2eSStefano Zampini } 2693927de2eSStefano Zampini PetscFunctionReturn(0); 2703927de2eSStefano Zampini } 2713927de2eSStefano Zampini 2723927de2eSStefano Zampini #undef __FUNCT__ 273b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ_IS" 274b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ_IS(Mat mat, MatReuse reuse, Mat *M) 275b7ce53b6SStefano Zampini { 276b7ce53b6SStefano Zampini Mat_IS *matis = (Mat_IS*)(mat->data); 277d9a9e74cSStefano Zampini Mat local_mat; 278b7ce53b6SStefano Zampini /* info on mat */ 2793cfa4ea4SStefano Zampini PetscInt bs,rows,cols,lrows,lcols; 280b7ce53b6SStefano Zampini PetscInt local_rows,local_cols; 281686e3a49SStefano Zampini PetscBool isdense,issbaij,isseqaij; 2827c03b4e8SStefano Zampini PetscMPIInt nsubdomains; 283b7ce53b6SStefano Zampini /* values insertion */ 284b7ce53b6SStefano Zampini PetscScalar *array; 285b7ce53b6SStefano Zampini /* work */ 286b7ce53b6SStefano Zampini PetscErrorCode ierr; 287b7ce53b6SStefano Zampini 288b7ce53b6SStefano Zampini PetscFunctionBegin; 289b7ce53b6SStefano Zampini /* get info from mat */ 2907c03b4e8SStefano Zampini ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&nsubdomains);CHKERRQ(ierr); 2917c03b4e8SStefano Zampini if (nsubdomains == 1) { 2927c03b4e8SStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { 2937c03b4e8SStefano Zampini ierr = MatDuplicate(matis->A,MAT_COPY_VALUES,&(*M));CHKERRQ(ierr); 2947c03b4e8SStefano Zampini } else { 2957c03b4e8SStefano Zampini ierr = MatCopy(matis->A,*M,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 2967c03b4e8SStefano Zampini } 2977c03b4e8SStefano Zampini PetscFunctionReturn(0); 2987c03b4e8SStefano Zampini } 299b7ce53b6SStefano Zampini ierr = MatGetSize(mat,&rows,&cols);CHKERRQ(ierr); 300b7ce53b6SStefano Zampini ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 3013cfa4ea4SStefano Zampini ierr = MatGetLocalSize(mat,&lrows,&lcols);CHKERRQ(ierr); 302b7ce53b6SStefano Zampini ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr); 303b7ce53b6SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr); 304686e3a49SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr); 305b7ce53b6SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr); 306b7ce53b6SStefano Zampini 307b7ce53b6SStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { 308b7ce53b6SStefano Zampini MatType new_mat_type; 3093927de2eSStefano Zampini PetscBool issbaij_red; 310b7ce53b6SStefano Zampini 311b7ce53b6SStefano Zampini /* determining new matrix type */ 312b2566f29SBarry Smith ierr = MPIU_Allreduce(&issbaij,&issbaij_red,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 313b7ce53b6SStefano Zampini if (issbaij_red) { 314b7ce53b6SStefano Zampini new_mat_type = MATSBAIJ; 315b7ce53b6SStefano Zampini } else { 316b7ce53b6SStefano Zampini if (bs>1) { 317b7ce53b6SStefano Zampini new_mat_type = MATBAIJ; 318b7ce53b6SStefano Zampini } else { 319b7ce53b6SStefano Zampini new_mat_type = MATAIJ; 320b7ce53b6SStefano Zampini } 321b7ce53b6SStefano Zampini } 322b7ce53b6SStefano Zampini 3233927de2eSStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)mat),M);CHKERRQ(ierr); 3243cfa4ea4SStefano Zampini ierr = MatSetSizes(*M,lrows,lcols,rows,cols);CHKERRQ(ierr); 3253927de2eSStefano Zampini ierr = MatSetBlockSize(*M,bs);CHKERRQ(ierr); 3263927de2eSStefano Zampini ierr = MatSetType(*M,new_mat_type);CHKERRQ(ierr); 3273927de2eSStefano Zampini ierr = MatISSetMPIXAIJPreallocation_Private(mat,*M,PETSC_FALSE);CHKERRQ(ierr); 328b7ce53b6SStefano Zampini } else { 3293cfa4ea4SStefano Zampini PetscInt mbs,mrows,mcols,mlrows,mlcols; 330b7ce53b6SStefano Zampini /* some checks */ 331b7ce53b6SStefano Zampini ierr = MatGetBlockSize(*M,&mbs);CHKERRQ(ierr); 332b7ce53b6SStefano Zampini ierr = MatGetSize(*M,&mrows,&mcols);CHKERRQ(ierr); 3333cfa4ea4SStefano Zampini ierr = MatGetLocalSize(*M,&mlrows,&mlcols);CHKERRQ(ierr); 3346c4ed002SBarry Smith if (mrows != rows) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of rows (%d != %d)",rows,mrows); 3356c4ed002SBarry Smith if (mcols != cols) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of cols (%d != %d)",cols,mcols); 3366c4ed002SBarry Smith if (mlrows != lrows) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local rows (%d != %d)",lrows,mlrows); 3376c4ed002SBarry Smith if (mlcols != lcols) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local cols (%d != %d)",lcols,mlcols); 3386c4ed002SBarry Smith if (mbs != bs) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong block size (%d != %d)",bs,mbs); 339b7ce53b6SStefano Zampini ierr = MatZeroEntries(*M);CHKERRQ(ierr); 340b7ce53b6SStefano Zampini } 341d9a9e74cSStefano Zampini 342d9a9e74cSStefano Zampini if (issbaij) { 343d9a9e74cSStefano Zampini ierr = MatConvert(matis->A,MATSEQBAIJ,MAT_INITIAL_MATRIX,&local_mat);CHKERRQ(ierr); 344d9a9e74cSStefano Zampini } else { 345d9a9e74cSStefano Zampini ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr); 346d9a9e74cSStefano Zampini local_mat = matis->A; 347d9a9e74cSStefano Zampini } 348686e3a49SStefano Zampini 349b7ce53b6SStefano Zampini /* Set values */ 350ecf5a873SStefano Zampini ierr = MatSetLocalToGlobalMapping(*M,mat->rmap->mapping,mat->cmap->mapping);CHKERRQ(ierr); 351b7ce53b6SStefano Zampini if (isdense) { /* special case for dense local matrices */ 352ecf5a873SStefano Zampini PetscInt i,*dummy_rows,*dummy_cols; 353ecf5a873SStefano Zampini 354ecf5a873SStefano Zampini if (local_rows != local_cols) { 355ecf5a873SStefano Zampini ierr = PetscMalloc2(local_rows,&dummy_rows,local_cols,&dummy_cols);CHKERRQ(ierr); 356ecf5a873SStefano Zampini for (i=0;i<local_rows;i++) dummy_rows[i] = i; 357ecf5a873SStefano Zampini for (i=0;i<local_cols;i++) dummy_cols[i] = i; 358ecf5a873SStefano Zampini } else { 359ecf5a873SStefano Zampini ierr = PetscMalloc1(local_rows,&dummy_rows);CHKERRQ(ierr); 360ecf5a873SStefano Zampini for (i=0;i<local_rows;i++) dummy_rows[i] = i; 361ecf5a873SStefano Zampini dummy_cols = dummy_rows; 362ecf5a873SStefano Zampini } 363b7ce53b6SStefano Zampini ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr); 364d9a9e74cSStefano Zampini ierr = MatDenseGetArray(local_mat,&array);CHKERRQ(ierr); 365ecf5a873SStefano Zampini ierr = MatSetValuesLocal(*M,local_rows,dummy_rows,local_cols,dummy_cols,array,ADD_VALUES);CHKERRQ(ierr); 366d9a9e74cSStefano Zampini ierr = MatDenseRestoreArray(local_mat,&array);CHKERRQ(ierr); 367ecf5a873SStefano Zampini if (dummy_rows != dummy_cols) { 368ecf5a873SStefano Zampini ierr = PetscFree2(dummy_rows,dummy_cols);CHKERRQ(ierr); 369ecf5a873SStefano Zampini } else { 370ecf5a873SStefano Zampini ierr = PetscFree(dummy_rows);CHKERRQ(ierr); 371ecf5a873SStefano Zampini } 372686e3a49SStefano Zampini } else if (isseqaij) { 373ecf5a873SStefano Zampini PetscInt i,nvtxs,*xadj,*adjncy; 374686e3a49SStefano Zampini PetscBool done; 375686e3a49SStefano Zampini 376d9a9e74cSStefano Zampini ierr = MatGetRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr); 3776c4ed002SBarry Smith if (!done) SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__); 378d9a9e74cSStefano Zampini ierr = MatSeqAIJGetArray(local_mat,&array);CHKERRQ(ierr); 379686e3a49SStefano Zampini for (i=0;i<nvtxs;i++) { 380ecf5a873SStefano Zampini ierr = MatSetValuesLocal(*M,1,&i,xadj[i+1]-xadj[i],adjncy+xadj[i],array+xadj[i],ADD_VALUES);CHKERRQ(ierr); 381686e3a49SStefano Zampini } 382d9a9e74cSStefano Zampini ierr = MatRestoreRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr); 3836c4ed002SBarry Smith if (!done) SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__); 384d9a9e74cSStefano Zampini ierr = MatSeqAIJRestoreArray(local_mat,&array);CHKERRQ(ierr); 385686e3a49SStefano Zampini } else { /* very basic values insertion for all other matrix types */ 386ecf5a873SStefano Zampini PetscInt i; 387c0962df8SStefano Zampini 388686e3a49SStefano Zampini for (i=0;i<local_rows;i++) { 389686e3a49SStefano Zampini PetscInt j; 390ecf5a873SStefano Zampini const PetscInt *local_indices_cols; 391686e3a49SStefano Zampini 392ecf5a873SStefano Zampini ierr = MatGetRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr); 393ecf5a873SStefano Zampini ierr = MatSetValuesLocal(*M,1,&i,j,local_indices_cols,array,ADD_VALUES);CHKERRQ(ierr); 394ecf5a873SStefano Zampini ierr = MatRestoreRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr); 395686e3a49SStefano Zampini } 396b7ce53b6SStefano Zampini } 397b7ce53b6SStefano Zampini ierr = MatAssemblyBegin(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 398d9a9e74cSStefano Zampini ierr = MatDestroy(&local_mat);CHKERRQ(ierr); 399b7ce53b6SStefano Zampini ierr = MatAssemblyEnd(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 400b7ce53b6SStefano Zampini if (isdense) { 401b7ce53b6SStefano Zampini ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr); 402b7ce53b6SStefano Zampini } 403b7ce53b6SStefano Zampini PetscFunctionReturn(0); 404b7ce53b6SStefano Zampini } 405b7ce53b6SStefano Zampini 406b7ce53b6SStefano Zampini #undef __FUNCT__ 407b7ce53b6SStefano Zampini #define __FUNCT__ "MatISGetMPIXAIJ" 408b7ce53b6SStefano Zampini /*@ 409b7ce53b6SStefano Zampini MatISGetMPIXAIJ - Converts MATIS matrix into a parallel AIJ format 410b7ce53b6SStefano Zampini 411b7ce53b6SStefano Zampini Input Parameter: 412b7ce53b6SStefano Zampini . mat - the matrix (should be of type MATIS) 413b7ce53b6SStefano Zampini . reuse - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX 414b7ce53b6SStefano Zampini 415b7ce53b6SStefano Zampini Output Parameter: 416b7ce53b6SStefano Zampini . newmat - the matrix in AIJ format 417b7ce53b6SStefano Zampini 418b7ce53b6SStefano Zampini Level: developer 419b7ce53b6SStefano Zampini 420eb82efa4SStefano Zampini Notes: mat and *newmat cannot be the same object when MAT_REUSE_MATRIX is requested. 421b7ce53b6SStefano Zampini 422b7ce53b6SStefano Zampini .seealso: MATIS 423b7ce53b6SStefano Zampini @*/ 424b7ce53b6SStefano Zampini PetscErrorCode MatISGetMPIXAIJ(Mat mat, MatReuse reuse, Mat *newmat) 425b7ce53b6SStefano Zampini { 426b7ce53b6SStefano Zampini PetscErrorCode ierr; 427b7ce53b6SStefano Zampini 428b7ce53b6SStefano Zampini PetscFunctionBegin; 429b7ce53b6SStefano Zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 430b7ce53b6SStefano Zampini PetscValidLogicalCollectiveEnum(mat,reuse,2); 431b7ce53b6SStefano Zampini PetscValidPointer(newmat,3); 432b7ce53b6SStefano Zampini if (reuse != MAT_INITIAL_MATRIX) { 433b7ce53b6SStefano Zampini PetscValidHeaderSpecific(*newmat,MAT_CLASSID,3); 434b7ce53b6SStefano Zampini PetscCheckSameComm(mat,1,*newmat,3); 4356c4ed002SBarry Smith if (mat == *newmat) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse the same matrix"); 436b7ce53b6SStefano Zampini } 437b7ce53b6SStefano Zampini ierr = PetscUseMethod(mat,"MatISGetMPIXAIJ_C",(Mat,MatReuse,Mat*),(mat,reuse,newmat));CHKERRQ(ierr); 438b7ce53b6SStefano Zampini PetscFunctionReturn(0); 439b7ce53b6SStefano Zampini } 440b7ce53b6SStefano Zampini 441b7ce53b6SStefano Zampini #undef __FUNCT__ 442ad6194a2SStefano Zampini #define __FUNCT__ "MatDuplicate_IS" 443ad6194a2SStefano Zampini PetscErrorCode MatDuplicate_IS(Mat mat,MatDuplicateOption op,Mat *newmat) 444ad6194a2SStefano Zampini { 445ad6194a2SStefano Zampini PetscErrorCode ierr; 446ad6194a2SStefano Zampini Mat_IS *matis = (Mat_IS*)(mat->data); 447ad6194a2SStefano Zampini PetscInt bs,m,n,M,N; 448ad6194a2SStefano Zampini Mat B,localmat; 449ad6194a2SStefano Zampini 450ad6194a2SStefano Zampini PetscFunctionBegin; 451ad6194a2SStefano Zampini ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 452ad6194a2SStefano Zampini ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr); 453ad6194a2SStefano Zampini ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr); 454e176bc59SStefano Zampini ierr = MatCreateIS(PetscObjectComm((PetscObject)mat),bs,m,n,M,N,mat->rmap->mapping,mat->cmap->mapping,&B);CHKERRQ(ierr); 455ad6194a2SStefano Zampini ierr = MatDuplicate(matis->A,op,&localmat);CHKERRQ(ierr); 456ad6194a2SStefano Zampini ierr = MatISSetLocalMat(B,localmat);CHKERRQ(ierr); 457b3317aa8SStefano Zampini ierr = MatDestroy(&localmat);CHKERRQ(ierr); 458ad6194a2SStefano Zampini ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 459ad6194a2SStefano Zampini ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 460ad6194a2SStefano Zampini *newmat = B; 461ad6194a2SStefano Zampini PetscFunctionReturn(0); 462ad6194a2SStefano Zampini } 463ad6194a2SStefano Zampini 464ad6194a2SStefano Zampini #undef __FUNCT__ 46569796d55SStefano Zampini #define __FUNCT__ "MatIsHermitian_IS" 46669796d55SStefano Zampini PetscErrorCode MatIsHermitian_IS(Mat A,PetscReal tol,PetscBool *flg) 46769796d55SStefano Zampini { 46869796d55SStefano Zampini PetscErrorCode ierr; 46969796d55SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 47069796d55SStefano Zampini PetscBool local_sym; 47169796d55SStefano Zampini 47269796d55SStefano Zampini PetscFunctionBegin; 47369796d55SStefano Zampini ierr = MatIsHermitian(matis->A,tol,&local_sym);CHKERRQ(ierr); 474b2566f29SBarry Smith ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 47569796d55SStefano Zampini PetscFunctionReturn(0); 47669796d55SStefano Zampini } 47769796d55SStefano Zampini 47869796d55SStefano Zampini #undef __FUNCT__ 47969796d55SStefano Zampini #define __FUNCT__ "MatIsSymmetric_IS" 48069796d55SStefano Zampini PetscErrorCode MatIsSymmetric_IS(Mat A,PetscReal tol,PetscBool *flg) 48169796d55SStefano Zampini { 48269796d55SStefano Zampini PetscErrorCode ierr; 48369796d55SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 48469796d55SStefano Zampini PetscBool local_sym; 48569796d55SStefano Zampini 48669796d55SStefano Zampini PetscFunctionBegin; 48769796d55SStefano Zampini ierr = MatIsSymmetric(matis->A,tol,&local_sym);CHKERRQ(ierr); 488b2566f29SBarry Smith ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 48969796d55SStefano Zampini PetscFunctionReturn(0); 49069796d55SStefano Zampini } 49169796d55SStefano Zampini 49269796d55SStefano Zampini #undef __FUNCT__ 493b4319ba4SBarry Smith #define __FUNCT__ "MatDestroy_IS" 494dfbe8321SBarry Smith PetscErrorCode MatDestroy_IS(Mat A) 495b4319ba4SBarry Smith { 496dfbe8321SBarry Smith PetscErrorCode ierr; 497b4319ba4SBarry Smith Mat_IS *b = (Mat_IS*)A->data; 498b4319ba4SBarry Smith 499b4319ba4SBarry Smith PetscFunctionBegin; 5006bf464f9SBarry Smith ierr = MatDestroy(&b->A);CHKERRQ(ierr); 501e176bc59SStefano Zampini ierr = VecScatterDestroy(&b->cctx);CHKERRQ(ierr); 502e176bc59SStefano Zampini ierr = VecScatterDestroy(&b->rctx);CHKERRQ(ierr); 5036bf464f9SBarry Smith ierr = VecDestroy(&b->x);CHKERRQ(ierr); 5046bf464f9SBarry Smith ierr = VecDestroy(&b->y);CHKERRQ(ierr); 50528f4e0baSStefano Zampini ierr = PetscSFDestroy(&b->sf);CHKERRQ(ierr); 50628f4e0baSStefano Zampini ierr = PetscFree2(b->sf_rootdata,b->sf_leafdata);CHKERRQ(ierr); 507bf0cc555SLisandro Dalcin ierr = PetscFree(A->data);CHKERRQ(ierr); 508dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr); 509bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",NULL);CHKERRQ(ierr); 510b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",NULL);CHKERRQ(ierr); 511b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",NULL);CHKERRQ(ierr); 5122e1947a5SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",NULL);CHKERRQ(ierr); 513b4319ba4SBarry Smith PetscFunctionReturn(0); 514b4319ba4SBarry Smith } 515b4319ba4SBarry Smith 516b4319ba4SBarry Smith #undef __FUNCT__ 517b4319ba4SBarry Smith #define __FUNCT__ "MatMult_IS" 518dfbe8321SBarry Smith PetscErrorCode MatMult_IS(Mat A,Vec x,Vec y) 519b4319ba4SBarry Smith { 520dfbe8321SBarry Smith PetscErrorCode ierr; 521b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 522b4319ba4SBarry Smith PetscScalar zero = 0.0; 523b4319ba4SBarry Smith 524b4319ba4SBarry Smith PetscFunctionBegin; 525b4319ba4SBarry Smith /* scatter the global vector x into the local work vector */ 526e176bc59SStefano Zampini ierr = VecScatterBegin(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 527e176bc59SStefano Zampini ierr = VecScatterEnd(is->cctx,x,is->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 528b4319ba4SBarry Smith 529b4319ba4SBarry Smith /* multiply the local matrix */ 530b4319ba4SBarry Smith ierr = MatMult(is->A,is->x,is->y);CHKERRQ(ierr); 531b4319ba4SBarry Smith 532b4319ba4SBarry Smith /* scatter product back into global memory */ 5332dcb1b2aSMatthew Knepley ierr = VecSet(y,zero);CHKERRQ(ierr); 534e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 535e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,is->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 536b4319ba4SBarry Smith PetscFunctionReturn(0); 537b4319ba4SBarry Smith } 538b4319ba4SBarry Smith 539b4319ba4SBarry Smith #undef __FUNCT__ 5402e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultAdd_IS" 541b7ce53b6SStefano Zampini PetscErrorCode MatMultAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 5422e74eeadSLisandro Dalcin { 543650997f4SStefano Zampini Vec temp_vec; 5442e74eeadSLisandro Dalcin PetscErrorCode ierr; 5452e74eeadSLisandro Dalcin 5462e74eeadSLisandro Dalcin PetscFunctionBegin; /* v3 = v2 + A * v1.*/ 547650997f4SStefano Zampini if (v3 != v2) { 548650997f4SStefano Zampini ierr = MatMult(A,v1,v3);CHKERRQ(ierr); 549650997f4SStefano Zampini ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 550650997f4SStefano Zampini } else { 551650997f4SStefano Zampini ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 552650997f4SStefano Zampini ierr = MatMult(A,v1,temp_vec);CHKERRQ(ierr); 553650997f4SStefano Zampini ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 554650997f4SStefano Zampini ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 555650997f4SStefano Zampini ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 556650997f4SStefano Zampini } 5572e74eeadSLisandro Dalcin PetscFunctionReturn(0); 5582e74eeadSLisandro Dalcin } 5592e74eeadSLisandro Dalcin 5602e74eeadSLisandro Dalcin #undef __FUNCT__ 5612e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTranspose_IS" 562e176bc59SStefano Zampini PetscErrorCode MatMultTranspose_IS(Mat A,Vec y,Vec x) 5632e74eeadSLisandro Dalcin { 5642e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 5652e74eeadSLisandro Dalcin PetscErrorCode ierr; 5662e74eeadSLisandro Dalcin 567e176bc59SStefano Zampini PetscFunctionBegin; 5682e74eeadSLisandro Dalcin /* scatter the global vector x into the local work vector */ 569e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 570e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,y,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5712e74eeadSLisandro Dalcin 5722e74eeadSLisandro Dalcin /* multiply the local matrix */ 573e176bc59SStefano Zampini ierr = MatMultTranspose(is->A,is->y,is->x);CHKERRQ(ierr); 5742e74eeadSLisandro Dalcin 5752e74eeadSLisandro Dalcin /* scatter product back into global vector */ 576e176bc59SStefano Zampini ierr = VecSet(x,0);CHKERRQ(ierr); 577e176bc59SStefano Zampini ierr = VecScatterBegin(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 578e176bc59SStefano Zampini ierr = VecScatterEnd(is->cctx,is->x,x,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 5792e74eeadSLisandro Dalcin PetscFunctionReturn(0); 5802e74eeadSLisandro Dalcin } 5812e74eeadSLisandro Dalcin 5822e74eeadSLisandro Dalcin #undef __FUNCT__ 5832e74eeadSLisandro Dalcin #define __FUNCT__ "MatMultTransposeAdd_IS" 5842e74eeadSLisandro Dalcin PetscErrorCode MatMultTransposeAdd_IS(Mat A,Vec v1,Vec v2,Vec v3) 5852e74eeadSLisandro Dalcin { 586650997f4SStefano Zampini Vec temp_vec; 5872e74eeadSLisandro Dalcin PetscErrorCode ierr; 5882e74eeadSLisandro Dalcin 5892e74eeadSLisandro Dalcin PetscFunctionBegin; /* v3 = v2 + A' * v1.*/ 590650997f4SStefano Zampini if (v3 != v2) { 591650997f4SStefano Zampini ierr = MatMultTranspose(A,v1,v3);CHKERRQ(ierr); 592650997f4SStefano Zampini ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 593650997f4SStefano Zampini } else { 594650997f4SStefano Zampini ierr = VecDuplicate(v2,&temp_vec);CHKERRQ(ierr); 595650997f4SStefano Zampini ierr = MatMultTranspose(A,v1,temp_vec);CHKERRQ(ierr); 596650997f4SStefano Zampini ierr = VecAXPY(temp_vec,1.0,v2);CHKERRQ(ierr); 597650997f4SStefano Zampini ierr = VecCopy(temp_vec,v3);CHKERRQ(ierr); 598650997f4SStefano Zampini ierr = VecDestroy(&temp_vec);CHKERRQ(ierr); 599650997f4SStefano Zampini } 6002e74eeadSLisandro Dalcin PetscFunctionReturn(0); 6012e74eeadSLisandro Dalcin } 6022e74eeadSLisandro Dalcin 6032e74eeadSLisandro Dalcin #undef __FUNCT__ 604b4319ba4SBarry Smith #define __FUNCT__ "MatView_IS" 605dfbe8321SBarry Smith PetscErrorCode MatView_IS(Mat A,PetscViewer viewer) 606b4319ba4SBarry Smith { 607b4319ba4SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 608dfbe8321SBarry Smith PetscErrorCode ierr; 609b4319ba4SBarry Smith PetscViewer sviewer; 610b4319ba4SBarry Smith 611b4319ba4SBarry Smith PetscFunctionBegin; 6123f08860eSBarry Smith ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 613b4319ba4SBarry Smith ierr = MatView(a->A,sviewer);CHKERRQ(ierr); 6143f08860eSBarry Smith ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 615*6e520ac8SStefano Zampini ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 616b4319ba4SBarry Smith PetscFunctionReturn(0); 617b4319ba4SBarry Smith } 618b4319ba4SBarry Smith 619b4319ba4SBarry Smith #undef __FUNCT__ 620b4319ba4SBarry Smith #define __FUNCT__ "MatSetLocalToGlobalMapping_IS" 621784ac674SJed Brown PetscErrorCode MatSetLocalToGlobalMapping_IS(Mat A,ISLocalToGlobalMapping rmapping,ISLocalToGlobalMapping cmapping) 622b4319ba4SBarry Smith { 623dfbe8321SBarry Smith PetscErrorCode ierr; 624e176bc59SStefano Zampini PetscInt nr,rbs,nc,cbs; 625b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 626b4319ba4SBarry Smith IS from,to; 627e176bc59SStefano Zampini Vec cglobal,rglobal; 628b4319ba4SBarry Smith 629b4319ba4SBarry Smith PetscFunctionBegin; 630784ac674SJed Brown PetscCheckSameComm(A,1,rmapping,2); 631e176bc59SStefano Zampini PetscCheckSameComm(A,1,cmapping,3); 6323bbff08aSStefano Zampini /* Destroy any previous data */ 63370cf5478SStefano Zampini ierr = VecDestroy(&is->x);CHKERRQ(ierr); 63470cf5478SStefano Zampini ierr = VecDestroy(&is->y);CHKERRQ(ierr); 635e176bc59SStefano Zampini ierr = VecScatterDestroy(&is->rctx);CHKERRQ(ierr); 636e176bc59SStefano Zampini ierr = VecScatterDestroy(&is->cctx);CHKERRQ(ierr); 6371c47cb0fSStefano Zampini ierr = MatDestroy(&is->A);CHKERRQ(ierr); 63828f4e0baSStefano Zampini ierr = PetscSFDestroy(&is->sf);CHKERRQ(ierr); 63928f4e0baSStefano Zampini ierr = PetscFree2(is->sf_rootdata,is->sf_leafdata);CHKERRQ(ierr); 6403bbff08aSStefano Zampini 6413bbff08aSStefano Zampini /* Setup Layout and set local to global maps */ 642fc27028aSStefano Zampini ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 643fc27028aSStefano Zampini ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 644fa7f1dd8SStefano Zampini ierr = PetscLayoutSetISLocalToGlobalMapping(A->rmap,rmapping);CHKERRQ(ierr); 645fa7f1dd8SStefano Zampini ierr = PetscLayoutSetISLocalToGlobalMapping(A->cmap,cmapping);CHKERRQ(ierr); 646b4319ba4SBarry Smith 647b4319ba4SBarry Smith /* Create the local matrix A */ 648e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(rmapping,&nr);CHKERRQ(ierr); 649e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetBlockSize(rmapping,&rbs);CHKERRQ(ierr); 650e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetSize(cmapping,&nc);CHKERRQ(ierr); 651e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingGetBlockSize(cmapping,&cbs);CHKERRQ(ierr); 652f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,&is->A);CHKERRQ(ierr); 653e176bc59SStefano Zampini ierr = MatSetType(is->A,MATAIJ);CHKERRQ(ierr); 654e176bc59SStefano Zampini ierr = MatSetSizes(is->A,nr,nc,nr,nc);CHKERRQ(ierr); 655e176bc59SStefano Zampini ierr = MatSetBlockSizes(is->A,rbs,cbs);CHKERRQ(ierr); 656ff130e51SJed Brown ierr = MatSetOptionsPrefix(is->A,((PetscObject)A)->prefix);CHKERRQ(ierr); 657ff130e51SJed Brown ierr = MatAppendOptionsPrefix(is->A,"is_");CHKERRQ(ierr); 658b4319ba4SBarry Smith ierr = MatSetFromOptions(is->A);CHKERRQ(ierr); 659b4319ba4SBarry Smith 660b4319ba4SBarry Smith /* Create the local work vectors */ 6613bbff08aSStefano Zampini ierr = MatCreateVecs(is->A,&is->x,&is->y);CHKERRQ(ierr); 662b4319ba4SBarry Smith 663e176bc59SStefano Zampini /* setup the global to local scatters */ 664e176bc59SStefano Zampini ierr = MatCreateVecs(A,&cglobal,&rglobal);CHKERRQ(ierr); 665e176bc59SStefano Zampini ierr = ISCreateStride(PETSC_COMM_SELF,nr,0,1,&to);CHKERRQ(ierr); 666784ac674SJed Brown ierr = ISLocalToGlobalMappingApplyIS(rmapping,to,&from);CHKERRQ(ierr); 667e176bc59SStefano Zampini ierr = VecScatterCreate(rglobal,from,is->y,to,&is->rctx);CHKERRQ(ierr); 668e176bc59SStefano Zampini if (rmapping != cmapping) { 669e176bc59SStefano Zampini ierr = ISDestroy(&to);CHKERRQ(ierr); 670e176bc59SStefano Zampini ierr = ISDestroy(&from);CHKERRQ(ierr); 671e176bc59SStefano Zampini ierr = ISCreateStride(PETSC_COMM_SELF,nc,0,1,&to);CHKERRQ(ierr); 672e176bc59SStefano Zampini ierr = ISLocalToGlobalMappingApplyIS(cmapping,to,&from);CHKERRQ(ierr); 673e176bc59SStefano Zampini ierr = VecScatterCreate(cglobal,from,is->x,to,&is->cctx);CHKERRQ(ierr); 674e176bc59SStefano Zampini } else { 675e176bc59SStefano Zampini ierr = PetscObjectReference((PetscObject)is->rctx);CHKERRQ(ierr); 676e176bc59SStefano Zampini is->cctx = is->rctx; 677e176bc59SStefano Zampini } 678e176bc59SStefano Zampini ierr = VecDestroy(&rglobal);CHKERRQ(ierr); 679e176bc59SStefano Zampini ierr = VecDestroy(&cglobal);CHKERRQ(ierr); 6806bf464f9SBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 6816bf464f9SBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 68248ff6bf3SStefano Zampini ierr = MatSetUp(A);CHKERRQ(ierr); 683b4319ba4SBarry Smith PetscFunctionReturn(0); 684b4319ba4SBarry Smith } 685b4319ba4SBarry Smith 6862e74eeadSLisandro Dalcin #undef __FUNCT__ 6872e74eeadSLisandro Dalcin #define __FUNCT__ "MatSetValues_IS" 6882e74eeadSLisandro Dalcin PetscErrorCode MatSetValues_IS(Mat mat, PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols, const PetscScalar *values, InsertMode addv) 6892e74eeadSLisandro Dalcin { 6902e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)mat->data; 6912e74eeadSLisandro Dalcin PetscInt rows_l[2048],cols_l[2048]; 6922e74eeadSLisandro Dalcin PetscErrorCode ierr; 6932e74eeadSLisandro Dalcin 6942e74eeadSLisandro Dalcin PetscFunctionBegin; 6952e74eeadSLisandro Dalcin #if defined(PETSC_USE_DEBUG) 696f23aa3ddSBarry 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); 6972e74eeadSLisandro Dalcin #endif 6983bbff08aSStefano Zampini ierr = ISG2LMapApply(mat->rmap->mapping,m,rows,rows_l);CHKERRQ(ierr); 6993bbff08aSStefano Zampini ierr = ISG2LMapApply(mat->cmap->mapping,n,cols,cols_l);CHKERRQ(ierr); 7002e74eeadSLisandro Dalcin ierr = MatSetValues(is->A,m,rows_l,n,cols_l,values,addv);CHKERRQ(ierr); 7012e74eeadSLisandro Dalcin PetscFunctionReturn(0); 7022e74eeadSLisandro Dalcin } 7032e74eeadSLisandro Dalcin 704b4319ba4SBarry Smith #undef __FUNCT__ 705b4319ba4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_IS" 70613f74950SBarry Smith PetscErrorCode MatSetValuesLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 707b4319ba4SBarry Smith { 708dfbe8321SBarry Smith PetscErrorCode ierr; 709b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 710b4319ba4SBarry Smith 711b4319ba4SBarry Smith PetscFunctionBegin; 712b4319ba4SBarry Smith ierr = MatSetValues(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 713b4319ba4SBarry Smith PetscFunctionReturn(0); 714b4319ba4SBarry Smith } 715b4319ba4SBarry Smith 716b4319ba4SBarry Smith #undef __FUNCT__ 717f0006bf2SLisandro Dalcin #define __FUNCT__ "MatSetValuesBlockedLocal_IS" 718f0006bf2SLisandro Dalcin PetscErrorCode MatSetValuesBlockedLocal_IS(Mat A,PetscInt m,const PetscInt *rows, PetscInt n,const PetscInt *cols,const PetscScalar *values,InsertMode addv) 719f0006bf2SLisandro Dalcin { 720f0006bf2SLisandro Dalcin PetscErrorCode ierr; 721f0006bf2SLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 722f0006bf2SLisandro Dalcin 723f0006bf2SLisandro Dalcin PetscFunctionBegin; 724f0006bf2SLisandro Dalcin ierr = MatSetValuesBlocked(is->A,m,rows,n,cols,values,addv);CHKERRQ(ierr); 725f0006bf2SLisandro Dalcin PetscFunctionReturn(0); 726f0006bf2SLisandro Dalcin } 727f0006bf2SLisandro Dalcin 728f0006bf2SLisandro Dalcin #undef __FUNCT__ 7292e74eeadSLisandro Dalcin #define __FUNCT__ "MatZeroRows_IS" 7302b40b63fSBarry Smith PetscErrorCode MatZeroRows_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 7312e74eeadSLisandro Dalcin { 732*6e520ac8SStefano Zampini Mat_IS *matis = (Mat_IS*)A->data; 733*6e520ac8SStefano Zampini PetscInt nr,nl,len,i; 734*6e520ac8SStefano Zampini PetscInt *lrows; 7352e74eeadSLisandro Dalcin PetscErrorCode ierr; 7362e74eeadSLisandro Dalcin 7372e74eeadSLisandro Dalcin PetscFunctionBegin; 738*6e520ac8SStefano Zampini /* get locally owned rows */ 739*6e520ac8SStefano Zampini ierr = MatZeroRowsMapLocal_Private(A,n,rows,&len,&lrows);CHKERRQ(ierr); 740*6e520ac8SStefano Zampini /* fix right hand side if needed */ 741*6e520ac8SStefano Zampini if (x && b) { 742*6e520ac8SStefano Zampini const PetscScalar *xx; 743*6e520ac8SStefano Zampini PetscScalar *bb; 744*6e520ac8SStefano Zampini 745*6e520ac8SStefano Zampini ierr = VecGetArrayRead(x, &xx);CHKERRQ(ierr); 746*6e520ac8SStefano Zampini ierr = VecGetArray(b, &bb);CHKERRQ(ierr); 747*6e520ac8SStefano Zampini for (i=0;i<len;++i) bb[lrows[i]] = diag*xx[lrows[i]]; 748*6e520ac8SStefano Zampini ierr = VecRestoreArrayRead(x, &xx);CHKERRQ(ierr); 749*6e520ac8SStefano Zampini ierr = VecRestoreArray(b, &bb);CHKERRQ(ierr); 7502e74eeadSLisandro Dalcin } 751*6e520ac8SStefano Zampini /* get rows associated to the local matrices */ 752*6e520ac8SStefano Zampini if (!matis->sf) { /* setup SF if not yet created and allocate rootdata and leafdata */ 753*6e520ac8SStefano Zampini ierr = MatISComputeSF_Private(A);CHKERRQ(ierr); 754*6e520ac8SStefano Zampini } 755*6e520ac8SStefano Zampini ierr = MatGetSize(matis->A,&nl,NULL);CHKERRQ(ierr); 756*6e520ac8SStefano Zampini ierr = PetscMemzero(matis->sf_leafdata,nl*sizeof(PetscInt));CHKERRQ(ierr); 757*6e520ac8SStefano Zampini ierr = PetscMemzero(matis->sf_rootdata,A->rmap->n*sizeof(PetscInt));CHKERRQ(ierr); 758*6e520ac8SStefano Zampini for (i=0;i<len;i++) matis->sf_rootdata[lrows[i]] = 1; 759*6e520ac8SStefano Zampini ierr = PetscFree(lrows);CHKERRQ(ierr); 760*6e520ac8SStefano Zampini ierr = PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 761*6e520ac8SStefano Zampini ierr = PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);CHKERRQ(ierr); 762*6e520ac8SStefano Zampini ierr = PetscMalloc1(nl,&lrows);CHKERRQ(ierr); 763*6e520ac8SStefano Zampini for (i=0,nr=0;i<nl;i++) if (matis->sf_leafdata[i]) lrows[nr++] = i; 764*6e520ac8SStefano Zampini ierr = MatZeroRowsLocal(A,nr,lrows,diag,NULL,NULL);CHKERRQ(ierr); 765*6e520ac8SStefano Zampini ierr = PetscFree(lrows);CHKERRQ(ierr); 7662e74eeadSLisandro Dalcin PetscFunctionReturn(0); 7672e74eeadSLisandro Dalcin } 7682e74eeadSLisandro Dalcin 7692e74eeadSLisandro Dalcin #undef __FUNCT__ 770b4319ba4SBarry Smith #define __FUNCT__ "MatZeroRowsLocal_IS" 7712b40b63fSBarry Smith PetscErrorCode MatZeroRowsLocal_IS(Mat A,PetscInt n,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 772b4319ba4SBarry Smith { 773b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 774dfbe8321SBarry Smith PetscErrorCode ierr; 775f4df32b1SMatthew Knepley PetscInt i; 776b4319ba4SBarry Smith 777b4319ba4SBarry Smith PetscFunctionBegin; 778ce94432eSBarry Smith if (x && b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support"); 779*6e520ac8SStefano Zampini if (diag != 0.) { 780b4319ba4SBarry Smith /* 781b4319ba4SBarry Smith Set up is->x as a "counting vector". This is in order to MatMult_IS 782b4319ba4SBarry Smith work properly in the interface nodes. 783b4319ba4SBarry Smith */ 784b4319ba4SBarry Smith Vec counter; 785e176bc59SStefano Zampini ierr = MatCreateVecs(A,NULL,&counter);CHKERRQ(ierr); 786e176bc59SStefano Zampini ierr = VecSet(counter,0.);CHKERRQ(ierr); 787e176bc59SStefano Zampini ierr = VecSet(is->y,1.);CHKERRQ(ierr); 788e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,is->y,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 789e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,is->y,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 790e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,counter,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 791e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,counter,is->y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 7926bf464f9SBarry Smith ierr = VecDestroy(&counter);CHKERRQ(ierr); 793b4319ba4SBarry Smith } 794958c9bccSBarry Smith if (!n) { 795b4319ba4SBarry Smith is->pure_neumann = PETSC_TRUE; 796b4319ba4SBarry Smith } else { 797b4319ba4SBarry Smith is->pure_neumann = PETSC_FALSE; 7982205254eSKarl Rupp 7992b40b63fSBarry Smith ierr = MatZeroRows(is->A,n,rows,diag,0,0);CHKERRQ(ierr); 800*6e520ac8SStefano Zampini if (diag != 0.) { 801*6e520ac8SStefano Zampini const PetscScalar *array; 802*6e520ac8SStefano Zampini ierr = VecGetArrayRead(is->y,&array);CHKERRQ(ierr); 803b4319ba4SBarry Smith for (i=0; i<n; i++) { 804f4df32b1SMatthew Knepley ierr = MatSetValue(is->A,rows[i],rows[i],diag/(array[rows[i]]),INSERT_VALUES);CHKERRQ(ierr); 805b4319ba4SBarry Smith } 806*6e520ac8SStefano Zampini ierr = VecRestoreArrayRead(is->y,&array);CHKERRQ(ierr); 807*6e520ac8SStefano Zampini } 808b4319ba4SBarry Smith ierr = MatAssemblyBegin(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 809b4319ba4SBarry Smith ierr = MatAssemblyEnd(is->A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 810b4319ba4SBarry Smith } 811b4319ba4SBarry Smith PetscFunctionReturn(0); 812b4319ba4SBarry Smith } 813b4319ba4SBarry Smith 814b4319ba4SBarry Smith #undef __FUNCT__ 815b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyBegin_IS" 816dfbe8321SBarry Smith PetscErrorCode MatAssemblyBegin_IS(Mat A,MatAssemblyType type) 817b4319ba4SBarry Smith { 818b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 819dfbe8321SBarry Smith PetscErrorCode ierr; 820b4319ba4SBarry Smith 821b4319ba4SBarry Smith PetscFunctionBegin; 822b4319ba4SBarry Smith ierr = MatAssemblyBegin(is->A,type);CHKERRQ(ierr); 823b4319ba4SBarry Smith PetscFunctionReturn(0); 824b4319ba4SBarry Smith } 825b4319ba4SBarry Smith 826b4319ba4SBarry Smith #undef __FUNCT__ 827b4319ba4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_IS" 828dfbe8321SBarry Smith PetscErrorCode MatAssemblyEnd_IS(Mat A,MatAssemblyType type) 829b4319ba4SBarry Smith { 830b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)A->data; 831dfbe8321SBarry Smith PetscErrorCode ierr; 832b4319ba4SBarry Smith 833b4319ba4SBarry Smith PetscFunctionBegin; 834b4319ba4SBarry Smith ierr = MatAssemblyEnd(is->A,type);CHKERRQ(ierr); 835b4319ba4SBarry Smith PetscFunctionReturn(0); 836b4319ba4SBarry Smith } 837b4319ba4SBarry Smith 838b4319ba4SBarry Smith #undef __FUNCT__ 839b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat_IS" 8407087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat_IS(Mat mat,Mat *local) 841b4319ba4SBarry Smith { 842b4319ba4SBarry Smith Mat_IS *is = (Mat_IS*)mat->data; 843b4319ba4SBarry Smith 844b4319ba4SBarry Smith PetscFunctionBegin; 845b4319ba4SBarry Smith *local = is->A; 846b4319ba4SBarry Smith PetscFunctionReturn(0); 847b4319ba4SBarry Smith } 848b4319ba4SBarry Smith 849b4319ba4SBarry Smith #undef __FUNCT__ 850b4319ba4SBarry Smith #define __FUNCT__ "MatISGetLocalMat" 851b4319ba4SBarry Smith /*@ 852b4319ba4SBarry Smith MatISGetLocalMat - Gets the local matrix stored inside a MATIS matrix. 853b4319ba4SBarry Smith 854b4319ba4SBarry Smith Input Parameter: 855b4319ba4SBarry Smith . mat - the matrix 856b4319ba4SBarry Smith 857b4319ba4SBarry Smith Output Parameter: 858eb82efa4SStefano Zampini . local - the local matrix 859b4319ba4SBarry Smith 860b4319ba4SBarry Smith Level: advanced 861b4319ba4SBarry Smith 862b4319ba4SBarry Smith Notes: 863b4319ba4SBarry Smith This can be called if you have precomputed the nonzero structure of the 864b4319ba4SBarry Smith matrix and want to provide it to the inner matrix object to improve the performance 865b4319ba4SBarry Smith of the MatSetValues() operation. 866b4319ba4SBarry Smith 867b4319ba4SBarry Smith .seealso: MATIS 868b4319ba4SBarry Smith @*/ 8697087cfbeSBarry Smith PetscErrorCode MatISGetLocalMat(Mat mat,Mat *local) 870b4319ba4SBarry Smith { 8714ac538c5SBarry Smith PetscErrorCode ierr; 872b4319ba4SBarry Smith 873b4319ba4SBarry Smith PetscFunctionBegin; 8740700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 875b4319ba4SBarry Smith PetscValidPointer(local,2); 8764ac538c5SBarry Smith ierr = PetscUseMethod(mat,"MatISGetLocalMat_C",(Mat,Mat*),(mat,local));CHKERRQ(ierr); 877b4319ba4SBarry Smith PetscFunctionReturn(0); 878b4319ba4SBarry Smith } 879b4319ba4SBarry Smith 8803b03a366Sstefano_zampini #undef __FUNCT__ 8813b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat_IS" 8823b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat_IS(Mat mat,Mat local) 8833b03a366Sstefano_zampini { 8843b03a366Sstefano_zampini Mat_IS *is = (Mat_IS*)mat->data; 8853b03a366Sstefano_zampini PetscInt nrows,ncols,orows,ocols; 8863b03a366Sstefano_zampini PetscErrorCode ierr; 8873b03a366Sstefano_zampini 8883b03a366Sstefano_zampini PetscFunctionBegin; 8894e4c7dbeSStefano Zampini if (is->A) { 8903b03a366Sstefano_zampini ierr = MatGetSize(is->A,&orows,&ocols);CHKERRQ(ierr); 8913b03a366Sstefano_zampini ierr = MatGetSize(local,&nrows,&ncols);CHKERRQ(ierr); 892f23aa3ddSBarry 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); 8934e4c7dbeSStefano Zampini } 8943b03a366Sstefano_zampini ierr = PetscObjectReference((PetscObject)local);CHKERRQ(ierr); 8953b03a366Sstefano_zampini ierr = MatDestroy(&is->A);CHKERRQ(ierr); 8963b03a366Sstefano_zampini is->A = local; 8973b03a366Sstefano_zampini PetscFunctionReturn(0); 8983b03a366Sstefano_zampini } 8993b03a366Sstefano_zampini 9003b03a366Sstefano_zampini #undef __FUNCT__ 9013b03a366Sstefano_zampini #define __FUNCT__ "MatISSetLocalMat" 9023b03a366Sstefano_zampini /*@ 903eb82efa4SStefano Zampini MatISSetLocalMat - Replace the local matrix stored inside a MATIS object. 9043b03a366Sstefano_zampini 9053b03a366Sstefano_zampini Input Parameter: 9063b03a366Sstefano_zampini . mat - the matrix 907eb82efa4SStefano Zampini . local - the local matrix 9083b03a366Sstefano_zampini 9093b03a366Sstefano_zampini Output Parameter: 9103b03a366Sstefano_zampini 9113b03a366Sstefano_zampini Level: advanced 9123b03a366Sstefano_zampini 9133b03a366Sstefano_zampini Notes: 9143b03a366Sstefano_zampini This can be called if you have precomputed the local matrix and 9153b03a366Sstefano_zampini want to provide it to the matrix object MATIS. 9163b03a366Sstefano_zampini 9173b03a366Sstefano_zampini .seealso: MATIS 9183b03a366Sstefano_zampini @*/ 9193b03a366Sstefano_zampini PetscErrorCode MatISSetLocalMat(Mat mat,Mat local) 9203b03a366Sstefano_zampini { 9213b03a366Sstefano_zampini PetscErrorCode ierr; 9223b03a366Sstefano_zampini 9233b03a366Sstefano_zampini PetscFunctionBegin; 9243b03a366Sstefano_zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 925b7ce53b6SStefano Zampini PetscValidHeaderSpecific(local,MAT_CLASSID,2); 9263b03a366Sstefano_zampini ierr = PetscUseMethod(mat,"MatISSetLocalMat_C",(Mat,Mat),(mat,local));CHKERRQ(ierr); 9273b03a366Sstefano_zampini PetscFunctionReturn(0); 9283b03a366Sstefano_zampini } 9293b03a366Sstefano_zampini 9306726f965SBarry Smith #undef __FUNCT__ 9316726f965SBarry Smith #define __FUNCT__ "MatZeroEntries_IS" 9326726f965SBarry Smith PetscErrorCode MatZeroEntries_IS(Mat A) 9336726f965SBarry Smith { 9346726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 9356726f965SBarry Smith PetscErrorCode ierr; 9366726f965SBarry Smith 9376726f965SBarry Smith PetscFunctionBegin; 9386726f965SBarry Smith ierr = MatZeroEntries(a->A);CHKERRQ(ierr); 9396726f965SBarry Smith PetscFunctionReturn(0); 9406726f965SBarry Smith } 9416726f965SBarry Smith 9426726f965SBarry Smith #undef __FUNCT__ 9432e74eeadSLisandro Dalcin #define __FUNCT__ "MatScale_IS" 9442e74eeadSLisandro Dalcin PetscErrorCode MatScale_IS(Mat A,PetscScalar a) 9452e74eeadSLisandro Dalcin { 9462e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 9472e74eeadSLisandro Dalcin PetscErrorCode ierr; 9482e74eeadSLisandro Dalcin 9492e74eeadSLisandro Dalcin PetscFunctionBegin; 9502e74eeadSLisandro Dalcin ierr = MatScale(is->A,a);CHKERRQ(ierr); 9512e74eeadSLisandro Dalcin PetscFunctionReturn(0); 9522e74eeadSLisandro Dalcin } 9532e74eeadSLisandro Dalcin 9542e74eeadSLisandro Dalcin #undef __FUNCT__ 9552e74eeadSLisandro Dalcin #define __FUNCT__ "MatGetDiagonal_IS" 9562e74eeadSLisandro Dalcin PetscErrorCode MatGetDiagonal_IS(Mat A, Vec v) 9572e74eeadSLisandro Dalcin { 9582e74eeadSLisandro Dalcin Mat_IS *is = (Mat_IS*)A->data; 9592e74eeadSLisandro Dalcin PetscErrorCode ierr; 9602e74eeadSLisandro Dalcin 9612e74eeadSLisandro Dalcin PetscFunctionBegin; 9622e74eeadSLisandro Dalcin /* get diagonal of the local matrix */ 963e176bc59SStefano Zampini ierr = MatGetDiagonal(is->A,is->y);CHKERRQ(ierr); 9642e74eeadSLisandro Dalcin 9652e74eeadSLisandro Dalcin /* scatter diagonal back into global vector */ 9662e74eeadSLisandro Dalcin ierr = VecSet(v,0);CHKERRQ(ierr); 967e176bc59SStefano Zampini ierr = VecScatterBegin(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 968e176bc59SStefano Zampini ierr = VecScatterEnd(is->rctx,is->y,v,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 9692e74eeadSLisandro Dalcin PetscFunctionReturn(0); 9702e74eeadSLisandro Dalcin } 9712e74eeadSLisandro Dalcin 9722e74eeadSLisandro Dalcin #undef __FUNCT__ 9736726f965SBarry Smith #define __FUNCT__ "MatSetOption_IS" 974ace3abfcSBarry Smith PetscErrorCode MatSetOption_IS(Mat A,MatOption op,PetscBool flg) 9756726f965SBarry Smith { 9766726f965SBarry Smith Mat_IS *a = (Mat_IS*)A->data; 9776726f965SBarry Smith PetscErrorCode ierr; 9786726f965SBarry Smith 9796726f965SBarry Smith PetscFunctionBegin; 9804e0d8c25SBarry Smith ierr = MatSetOption(a->A,op,flg);CHKERRQ(ierr); 9816726f965SBarry Smith PetscFunctionReturn(0); 9826726f965SBarry Smith } 9836726f965SBarry Smith 984284134d9SBarry Smith #undef __FUNCT__ 985284134d9SBarry Smith #define __FUNCT__ "MatCreateIS" 986284134d9SBarry Smith /*@ 9873c212e90SHong Zhang MatCreateIS - Creates a "process" unassembled matrix, assembled on each 988284134d9SBarry Smith process but not across processes. 989284134d9SBarry Smith 990284134d9SBarry Smith Input Parameters: 991284134d9SBarry Smith + comm - MPI communicator that will share the matrix 992e176bc59SStefano Zampini . bs - block size of the matrix 993df3898eeSBarry Smith . m,n,M,N - local and/or global sizes of the left and right vector used in matrix vector products 994e176bc59SStefano Zampini . rmap - local to global map for rows 995e176bc59SStefano Zampini - cmap - local to global map for cols 996284134d9SBarry Smith 997284134d9SBarry Smith Output Parameter: 998284134d9SBarry Smith . A - the resulting matrix 999284134d9SBarry Smith 10008e6c10adSSatish Balay Level: advanced 10018e6c10adSSatish Balay 10023c212e90SHong Zhang Notes: See MATIS for more details. 10033c212e90SHong Zhang m and n are NOT related to the size of the map; they are the size of the part of the vector owned 1004e176bc59SStefano Zampini by that process. The sizes of rmap and cmap define the size of the local matrices. 10053c212e90SHong Zhang If either rmap or cmap are NULL, then the matrix is assumed to be square. 1006284134d9SBarry Smith 1007284134d9SBarry Smith .seealso: MATIS, MatSetLocalToGlobalMapping() 1008284134d9SBarry Smith @*/ 1009e176bc59SStefano Zampini PetscErrorCode MatCreateIS(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,ISLocalToGlobalMapping rmap,ISLocalToGlobalMapping cmap,Mat *A) 1010284134d9SBarry Smith { 1011284134d9SBarry Smith PetscErrorCode ierr; 1012284134d9SBarry Smith 1013284134d9SBarry Smith PetscFunctionBegin; 1014e176bc59SStefano Zampini if (!rmap && !cmap) SETERRQ(comm,PETSC_ERR_USER,"You need to provide at least one of the mapping"); 1015284134d9SBarry Smith ierr = MatCreate(comm,A);CHKERRQ(ierr); 1016d3ee2243SStefano Zampini ierr = MatSetBlockSize(*A,bs);CHKERRQ(ierr); 1017284134d9SBarry Smith ierr = MatSetSizes(*A,m,n,M,N);CHKERRQ(ierr); 1018284134d9SBarry Smith ierr = MatSetType(*A,MATIS);CHKERRQ(ierr); 1019e176bc59SStefano Zampini if (rmap && cmap) { 1020e176bc59SStefano Zampini ierr = MatSetLocalToGlobalMapping(*A,rmap,cmap);CHKERRQ(ierr); 1021e176bc59SStefano Zampini } else if (!rmap) { 1022e176bc59SStefano Zampini ierr = MatSetLocalToGlobalMapping(*A,cmap,cmap);CHKERRQ(ierr); 1023e176bc59SStefano Zampini } else { 1024e176bc59SStefano Zampini ierr = MatSetLocalToGlobalMapping(*A,rmap,rmap);CHKERRQ(ierr); 1025e176bc59SStefano Zampini } 1026284134d9SBarry Smith PetscFunctionReturn(0); 1027284134d9SBarry Smith } 1028284134d9SBarry Smith 1029b4319ba4SBarry Smith /*MC 1030eb82efa4SStefano Zampini MATIS - MATIS = "is" - A matrix type to be used for using the non-overlapping domain decomposition type preconditioners (e.g. PCBDDC). 1031b4319ba4SBarry Smith This stores the matrices in globally unassembled form. Each processor 1032b4319ba4SBarry Smith assembles only its local Neumann problem and the parallel matrix vector 1033b4319ba4SBarry Smith product is handled "implicitly". 1034b4319ba4SBarry Smith 1035b4319ba4SBarry Smith Operations Provided: 10366726f965SBarry Smith + MatMult() 10372e74eeadSLisandro Dalcin . MatMultAdd() 10382e74eeadSLisandro Dalcin . MatMultTranspose() 10392e74eeadSLisandro Dalcin . MatMultTransposeAdd() 10406726f965SBarry Smith . MatZeroEntries() 10416726f965SBarry Smith . MatSetOption() 10422e74eeadSLisandro Dalcin . MatZeroRows() 10436726f965SBarry Smith . MatZeroRowsLocal() 10442e74eeadSLisandro Dalcin . MatSetValues() 10456726f965SBarry Smith . MatSetValuesLocal() 10462e74eeadSLisandro Dalcin . MatScale() 10472e74eeadSLisandro Dalcin . MatGetDiagonal() 10486726f965SBarry Smith - MatSetLocalToGlobalMapping() 1049b4319ba4SBarry Smith 1050b4319ba4SBarry Smith Options Database Keys: 1051b4319ba4SBarry Smith . -mat_type is - sets the matrix type to "is" during a call to MatSetFromOptions() 1052b4319ba4SBarry Smith 1053b4319ba4SBarry Smith Notes: Options prefix for the inner matrix are given by -is_mat_xxx 1054b4319ba4SBarry Smith 1055b4319ba4SBarry Smith You must call MatSetLocalToGlobalMapping() before using this matrix type. 1056b4319ba4SBarry Smith 1057b4319ba4SBarry Smith You can do matrix preallocation on the local matrix after you obtain it with 1058eb82efa4SStefano Zampini MatISGetLocalMat(); otherwise, you could use MatISSetPreallocation() 1059b4319ba4SBarry Smith 1060b4319ba4SBarry Smith Level: advanced 1061b4319ba4SBarry Smith 10623c212e90SHong Zhang .seealso: Mat, MatISGetLocalMat(), MatSetLocalToGlobalMapping(), MatISSetPreallocation(), PCBDDC, MatCreateIS() 1063b4319ba4SBarry Smith 1064b4319ba4SBarry Smith M*/ 1065b4319ba4SBarry Smith 1066b4319ba4SBarry Smith #undef __FUNCT__ 1067b4319ba4SBarry Smith #define __FUNCT__ "MatCreate_IS" 10688cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_IS(Mat A) 1069b4319ba4SBarry Smith { 1070dfbe8321SBarry Smith PetscErrorCode ierr; 1071b4319ba4SBarry Smith Mat_IS *b; 1072b4319ba4SBarry Smith 1073b4319ba4SBarry Smith PetscFunctionBegin; 1074b00a9115SJed Brown ierr = PetscNewLog(A,&b);CHKERRQ(ierr); 1075b4319ba4SBarry Smith A->data = (void*)b; 1076b4319ba4SBarry Smith 1077e176bc59SStefano Zampini /* matrix ops */ 1078e176bc59SStefano Zampini ierr = PetscMemzero(A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 1079b4319ba4SBarry Smith A->ops->mult = MatMult_IS; 10802e74eeadSLisandro Dalcin A->ops->multadd = MatMultAdd_IS; 10812e74eeadSLisandro Dalcin A->ops->multtranspose = MatMultTranspose_IS; 10822e74eeadSLisandro Dalcin A->ops->multtransposeadd = MatMultTransposeAdd_IS; 1083b4319ba4SBarry Smith A->ops->destroy = MatDestroy_IS; 1084b4319ba4SBarry Smith A->ops->setlocaltoglobalmapping = MatSetLocalToGlobalMapping_IS; 10852e74eeadSLisandro Dalcin A->ops->setvalues = MatSetValues_IS; 1086b4319ba4SBarry Smith A->ops->setvalueslocal = MatSetValuesLocal_IS; 1087f0006bf2SLisandro Dalcin A->ops->setvaluesblockedlocal = MatSetValuesBlockedLocal_IS; 10882e74eeadSLisandro Dalcin A->ops->zerorows = MatZeroRows_IS; 1089b4319ba4SBarry Smith A->ops->zerorowslocal = MatZeroRowsLocal_IS; 1090b4319ba4SBarry Smith A->ops->assemblybegin = MatAssemblyBegin_IS; 1091b4319ba4SBarry Smith A->ops->assemblyend = MatAssemblyEnd_IS; 1092b4319ba4SBarry Smith A->ops->view = MatView_IS; 10936726f965SBarry Smith A->ops->zeroentries = MatZeroEntries_IS; 10942e74eeadSLisandro Dalcin A->ops->scale = MatScale_IS; 10952e74eeadSLisandro Dalcin A->ops->getdiagonal = MatGetDiagonal_IS; 10966726f965SBarry Smith A->ops->setoption = MatSetOption_IS; 109769796d55SStefano Zampini A->ops->ishermitian = MatIsHermitian_IS; 109869796d55SStefano Zampini A->ops->issymmetric = MatIsSymmetric_IS; 1099ad6194a2SStefano Zampini A->ops->duplicate = MatDuplicate_IS; 1100b4319ba4SBarry Smith 1101b7ce53b6SStefano Zampini /* special MATIS functions */ 1102bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetLocalMat_C",MatISGetLocalMat_IS);CHKERRQ(ierr); 1103bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetLocalMat_C",MatISSetLocalMat_IS);CHKERRQ(ierr); 1104b7ce53b6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISGetMPIXAIJ_C",MatISGetMPIXAIJ_IS);CHKERRQ(ierr); 11052e1947a5SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatISSetPreallocation_C",MatISSetPreallocation_IS);CHKERRQ(ierr); 110617667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATIS);CHKERRQ(ierr); 1107b4319ba4SBarry Smith PetscFunctionReturn(0); 1108b4319ba4SBarry Smith } 1109