173f4d377SMatthew Knepley /*$Id: mpiadj.c,v 1.66 2001/08/07 03:02:59 balay Exp $*/ 2b97cf49bSBarry Smith 3b97cf49bSBarry Smith /* 4b97cf49bSBarry Smith Defines the basic matrix operations for the ADJ adjacency list matrix data-structure. 5b97cf49bSBarry Smith */ 63eda8832SBarry Smith #include "src/mat/impls/adj/mpi/mpiadj.h" 7325e03aeSBarry Smith #include "petscsys.h" 8b97cf49bSBarry Smith 94a2ae208SSatish Balay #undef __FUNCT__ 104a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIAdj_ASCII" 11b0a32e0cSBarry Smith int MatView_MPIAdj_ASCII(Mat A,PetscViewer viewer) 12b97cf49bSBarry Smith { 133eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 14fb9695e5SSatish Balay int ierr,i,j,m = A->m; 15fb9695e5SSatish Balay char *name; 16ce1411ecSBarry Smith PetscViewerFormat format; 17b97cf49bSBarry Smith 18433994e6SBarry Smith PetscFunctionBegin; 193a7fca6bSBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 20b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 21fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_INFO) { 223a40ed3dSBarry Smith PetscFunctionReturn(0); 23fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_MATLAB) { 24ffac6cdbSBarry Smith SETERRQ(PETSC_ERR_SUP,"Matlab format not supported"); 250752156aSBarry Smith } else { 26b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 27b97cf49bSBarry Smith for (i=0; i<m; i++) { 28b0a32e0cSBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"row %d:",i+a->rstart);CHKERRQ(ierr); 29b97cf49bSBarry Smith for (j=a->i[i]; j<a->i[i+1]; j++) { 30b0a32e0cSBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer," %d ",a->j[j]);CHKERRQ(ierr); 310752156aSBarry Smith } 32b0a32e0cSBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr); 33b97cf49bSBarry Smith } 34b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 35b97cf49bSBarry Smith } 36b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 373a40ed3dSBarry Smith PetscFunctionReturn(0); 38b97cf49bSBarry Smith } 39b97cf49bSBarry Smith 404a2ae208SSatish Balay #undef __FUNCT__ 414a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIAdj" 42b0a32e0cSBarry Smith int MatView_MPIAdj(Mat A,PetscViewer viewer) 43b97cf49bSBarry Smith { 44b97cf49bSBarry Smith int ierr; 456831982aSBarry Smith PetscTruth isascii; 46b97cf49bSBarry Smith 47433994e6SBarry Smith PetscFunctionBegin; 48b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 490f5bd95cSBarry Smith if (isascii) { 503eda8832SBarry Smith ierr = MatView_MPIAdj_ASCII(A,viewer);CHKERRQ(ierr); 515cd90555SBarry Smith } else { 5229bbc08cSBarry Smith SETERRQ1(1,"Viewer type %s not supported by MPIAdj",((PetscObject)viewer)->type_name); 53b97cf49bSBarry Smith } 543a40ed3dSBarry Smith PetscFunctionReturn(0); 55b97cf49bSBarry Smith } 56b97cf49bSBarry Smith 574a2ae208SSatish Balay #undef __FUNCT__ 584a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_MPIAdj" 593eda8832SBarry Smith int MatDestroy_MPIAdj(Mat mat) 60b97cf49bSBarry Smith { 613eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)mat->data; 6294d884c6SBarry Smith int ierr; 6394d884c6SBarry Smith 6494d884c6SBarry Smith PetscFunctionBegin; 65aa482453SBarry Smith #if defined(PETSC_USE_LOG) 66b0a32e0cSBarry Smith PetscLogObjectState((PetscObject)mat,"Rows=%d, Cols=%d, NZ=%d",mat->m,mat->n,a->nz); 67b97cf49bSBarry Smith #endif 68606d414cSSatish Balay if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);} 690400557aSBarry Smith if (a->freeaij) { 70606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 71606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 721198db28SBarry Smith if (a->values) {ierr = PetscFree(a->values);CHKERRQ(ierr);} 730400557aSBarry Smith } 740400557aSBarry Smith ierr = PetscFree(a->rowners);CHKERRQ(ierr); 751198db28SBarry Smith ierr = PetscFree(a);CHKERRQ(ierr); 763a40ed3dSBarry Smith PetscFunctionReturn(0); 77b97cf49bSBarry Smith } 78b97cf49bSBarry Smith 794a2ae208SSatish Balay #undef __FUNCT__ 804a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIAdj" 813eda8832SBarry Smith int MatSetOption_MPIAdj(Mat A,MatOption op) 82b97cf49bSBarry Smith { 833eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 84b97cf49bSBarry Smith 85433994e6SBarry Smith PetscFunctionBegin; 8612c028f9SKris Buschelman switch (op) { 87*77e54ba9SKris Buschelman case MAT_SYMMETRIC: 88*77e54ba9SKris Buschelman break; 8912c028f9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 90b97cf49bSBarry Smith a->symmetric = PETSC_TRUE; 9112c028f9SKris Buschelman break; 9212c028f9SKris Buschelman default: 93b0a32e0cSBarry Smith PetscLogInfo(A,"MatSetOption_MPIAdj:Option ignored\n"); 9412c028f9SKris Buschelman break; 95b97cf49bSBarry Smith } 963a40ed3dSBarry Smith PetscFunctionReturn(0); 97b97cf49bSBarry Smith } 98b97cf49bSBarry Smith 99b97cf49bSBarry Smith 100b97cf49bSBarry Smith /* 101b97cf49bSBarry Smith Adds diagonal pointers to sparse matrix structure. 102b97cf49bSBarry Smith */ 103b97cf49bSBarry Smith 1044a2ae208SSatish Balay #undef __FUNCT__ 1054a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_MPIAdj" 1063eda8832SBarry Smith int MatMarkDiagonal_MPIAdj(Mat A) 107b97cf49bSBarry Smith { 1083eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 10982502324SSatish Balay int i,j,*diag,m = A->m,ierr; 110b97cf49bSBarry Smith 111433994e6SBarry Smith PetscFunctionBegin; 112b0a32e0cSBarry Smith ierr = PetscMalloc((m+1)*sizeof(int),&diag);CHKERRQ(ierr); 113b0a32e0cSBarry Smith PetscLogObjectMemory(A,(m+1)*sizeof(int)); 114273d9f13SBarry Smith for (i=0; i<A->m; i++) { 115b97cf49bSBarry Smith for (j=a->i[i]; j<a->i[i+1]; j++) { 116b97cf49bSBarry Smith if (a->j[j] == i) { 117b97cf49bSBarry Smith diag[i] = j; 118b97cf49bSBarry Smith break; 119b97cf49bSBarry Smith } 120b97cf49bSBarry Smith } 121b97cf49bSBarry Smith } 122b97cf49bSBarry Smith a->diag = diag; 1233a40ed3dSBarry Smith PetscFunctionReturn(0); 124b97cf49bSBarry Smith } 125b97cf49bSBarry Smith 1264a2ae208SSatish Balay #undef __FUNCT__ 1274a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_MPIAdj" 12887828ca2SBarry Smith int MatGetRow_MPIAdj(Mat A,int row,int *nz,int **idx,PetscScalar **v) 129b97cf49bSBarry Smith { 1303eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 131b97cf49bSBarry Smith int *itmp; 132b97cf49bSBarry Smith 133433994e6SBarry Smith PetscFunctionBegin; 1340752156aSBarry Smith row -= a->rstart; 1350752156aSBarry Smith 136273d9f13SBarry Smith if (row < 0 || row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row out of range"); 137b97cf49bSBarry Smith 138b97cf49bSBarry Smith *nz = a->i[row+1] - a->i[row]; 139b97cf49bSBarry Smith if (v) *v = PETSC_NULL; 140b97cf49bSBarry Smith if (idx) { 141b97cf49bSBarry Smith itmp = a->j + a->i[row]; 142b97cf49bSBarry Smith if (*nz) { 143b97cf49bSBarry Smith *idx = itmp; 144b97cf49bSBarry Smith } 145b97cf49bSBarry Smith else *idx = 0; 146b97cf49bSBarry Smith } 1473a40ed3dSBarry Smith PetscFunctionReturn(0); 148b97cf49bSBarry Smith } 149b97cf49bSBarry Smith 1504a2ae208SSatish Balay #undef __FUNCT__ 1514a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_MPIAdj" 15287828ca2SBarry Smith int MatRestoreRow_MPIAdj(Mat A,int row,int *nz,int **idx,PetscScalar **v) 153b97cf49bSBarry Smith { 154433994e6SBarry Smith PetscFunctionBegin; 1553a40ed3dSBarry Smith PetscFunctionReturn(0); 156b97cf49bSBarry Smith } 157b97cf49bSBarry Smith 1584a2ae208SSatish Balay #undef __FUNCT__ 1594a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_MPIAdj" 1603eda8832SBarry Smith int MatGetBlockSize_MPIAdj(Mat A,int *bs) 161b97cf49bSBarry Smith { 162433994e6SBarry Smith PetscFunctionBegin; 163b97cf49bSBarry Smith *bs = 1; 1643a40ed3dSBarry Smith PetscFunctionReturn(0); 165b97cf49bSBarry Smith } 166b97cf49bSBarry Smith 167b97cf49bSBarry Smith 1684a2ae208SSatish Balay #undef __FUNCT__ 1694a2ae208SSatish Balay #define __FUNCT__ "MatEqual_MPIAdj" 1703eda8832SBarry Smith int MatEqual_MPIAdj(Mat A,Mat B,PetscTruth* flg) 171b97cf49bSBarry Smith { 1723eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data,*b = (Mat_MPIAdj *)B->data; 1730f5bd95cSBarry Smith int ierr; 1740f5bd95cSBarry Smith PetscTruth flag; 175b97cf49bSBarry Smith 176433994e6SBarry Smith PetscFunctionBegin; 177b97cf49bSBarry Smith /* If the matrix dimensions are not equal,or no of nonzeros */ 178273d9f13SBarry Smith if ((A->m != B->m) ||(a->nz != b->nz)) { 1790f5bd95cSBarry Smith flag = PETSC_FALSE; 180b97cf49bSBarry Smith } 181b97cf49bSBarry Smith 182b97cf49bSBarry Smith /* if the a->i are the same */ 183273d9f13SBarry Smith ierr = PetscMemcmp(a->i,b->i,(A->m+1)*sizeof(int),&flag);CHKERRQ(ierr); 184b97cf49bSBarry Smith 185b97cf49bSBarry Smith /* if a->j are the same */ 1860f5bd95cSBarry Smith ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(int),&flag);CHKERRQ(ierr); 187b97cf49bSBarry Smith 188ca161407SBarry Smith ierr = MPI_Allreduce(&flag,flg,1,MPI_INT,MPI_LAND,A->comm);CHKERRQ(ierr); 1893a40ed3dSBarry Smith PetscFunctionReturn(0); 190b97cf49bSBarry Smith } 191b97cf49bSBarry Smith 1924a2ae208SSatish Balay #undef __FUNCT__ 1934a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_MPIAdj" 1944e7234bfSBarry Smith int MatGetRowIJ_MPIAdj(Mat A,int oshift,PetscTruth symmetric,int *m,int *ia[],int *ja[],PetscTruth *done) 1956cd91f64SBarry Smith { 196d42735a1SBarry Smith int ierr,size,i; 1976cd91f64SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data; 1986cd91f64SBarry Smith 1996cd91f64SBarry Smith PetscFunctionBegin; 2006cd91f64SBarry Smith ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); 2016cd91f64SBarry Smith if (size > 1) {*done = PETSC_FALSE; PetscFunctionReturn(0);} 2026cd91f64SBarry Smith *m = A->m; 2036cd91f64SBarry Smith *ia = a->i; 2046cd91f64SBarry Smith *ja = a->j; 2056cd91f64SBarry Smith *done = PETSC_TRUE; 206d42735a1SBarry Smith if (oshift) { 207d42735a1SBarry Smith for (i=0; i<(*ia)[*m]; i++) { 208d42735a1SBarry Smith (*ja)[i]++; 209d42735a1SBarry Smith } 210d42735a1SBarry Smith for (i=0; i<=(*m); i++) (*ia)[i]++; 211d42735a1SBarry Smith } 212d42735a1SBarry Smith PetscFunctionReturn(0); 213d42735a1SBarry Smith } 214d42735a1SBarry Smith 2154a2ae208SSatish Balay #undef __FUNCT__ 2164a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_MPIAdj" 2174e7234bfSBarry Smith int MatRestoreRowIJ_MPIAdj(Mat A,int oshift,PetscTruth symmetric,int *m,int *ia[],int *ja[],PetscTruth *done) 218d42735a1SBarry Smith { 219d42735a1SBarry Smith int i; 220d42735a1SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data; 221d42735a1SBarry Smith 222d42735a1SBarry Smith PetscFunctionBegin; 223c5b3d447SBarry Smith if (ia && a->i != *ia) SETERRQ(1,"ia passed back is not one obtained with MatGetRowIJ()"); 224c5b3d447SBarry Smith if (ja && a->j != *ja) SETERRQ(1,"ja passed back is not one obtained with MatGetRowIJ()"); 225d42735a1SBarry Smith if (oshift) { 226d42735a1SBarry Smith for (i=0; i<=(*m); i++) (*ia)[i]--; 227d42735a1SBarry Smith for (i=0; i<(*ia)[*m]; i++) { 228d42735a1SBarry Smith (*ja)[i]--; 229d42735a1SBarry Smith } 230d42735a1SBarry Smith } 2316cd91f64SBarry Smith PetscFunctionReturn(0); 2326cd91f64SBarry Smith } 233b97cf49bSBarry Smith 234b97cf49bSBarry Smith /* -------------------------------------------------------------------*/ 2350b3b1487SBarry Smith static struct _MatOps MatOps_Values = {0, 2363eda8832SBarry Smith MatGetRow_MPIAdj, 2373eda8832SBarry Smith MatRestoreRow_MPIAdj, 238b97cf49bSBarry Smith 0, 23997304618SKris Buschelman /* 4*/ 0, 240b97cf49bSBarry Smith 0, 241b97cf49bSBarry Smith 0, 242b97cf49bSBarry Smith 0, 2430b3b1487SBarry Smith 0, 2440b3b1487SBarry Smith 0, 24597304618SKris Buschelman /*10*/ 0, 2460b3b1487SBarry Smith 0, 2470b3b1487SBarry Smith 0, 2480b3b1487SBarry Smith 0, 2490b3b1487SBarry Smith 0, 25097304618SKris Buschelman /*15*/ 0, 2513eda8832SBarry Smith MatEqual_MPIAdj, 2520b3b1487SBarry Smith 0, 2530b3b1487SBarry Smith 0, 2540b3b1487SBarry Smith 0, 25597304618SKris Buschelman /*20*/ 0, 2560b3b1487SBarry Smith 0, 2570b3b1487SBarry Smith 0, 2583eda8832SBarry Smith MatSetOption_MPIAdj, 2590b3b1487SBarry Smith 0, 26097304618SKris Buschelman /*25*/ 0, 2610b3b1487SBarry Smith 0, 2620b3b1487SBarry Smith 0, 2630b3b1487SBarry Smith 0, 2640b3b1487SBarry Smith 0, 26597304618SKris Buschelman /*30*/ 0, 2660b3b1487SBarry Smith 0, 267273d9f13SBarry Smith 0, 268273d9f13SBarry Smith 0, 2690b3b1487SBarry Smith 0, 27097304618SKris Buschelman /*35*/ 0, 2710b3b1487SBarry Smith 0, 2720b3b1487SBarry Smith 0, 2730b3b1487SBarry Smith 0, 2740b3b1487SBarry Smith 0, 27597304618SKris Buschelman /*40*/ 0, 2760b3b1487SBarry Smith 0, 2770b3b1487SBarry Smith 0, 2780b3b1487SBarry Smith 0, 2790b3b1487SBarry Smith 0, 28097304618SKris Buschelman /*45*/ 0, 2810b3b1487SBarry Smith 0, 2820b3b1487SBarry Smith 0, 2830b3b1487SBarry Smith 0, 2840b3b1487SBarry Smith 0, 28597304618SKris Buschelman /*50*/ MatGetBlockSize_MPIAdj, 2866cd91f64SBarry Smith MatGetRowIJ_MPIAdj, 287d42735a1SBarry Smith MatRestoreRowIJ_MPIAdj, 288b97cf49bSBarry Smith 0, 289b97cf49bSBarry Smith 0, 29097304618SKris Buschelman /*55*/ 0, 291b97cf49bSBarry Smith 0, 2920752156aSBarry Smith 0, 2930752156aSBarry Smith 0, 2940b3b1487SBarry Smith 0, 29597304618SKris Buschelman /*60*/ 0, 296b9b97703SBarry Smith MatDestroy_MPIAdj, 297b9b97703SBarry Smith MatView_MPIAdj, 29897304618SKris Buschelman MatGetPetscMaps_Petsc, 29997304618SKris Buschelman 0, 30097304618SKris Buschelman /*65*/ 0, 30197304618SKris Buschelman 0, 30297304618SKris Buschelman 0, 30397304618SKris Buschelman 0, 30497304618SKris Buschelman 0, 30597304618SKris Buschelman /*70*/ 0, 30697304618SKris Buschelman 0, 30797304618SKris Buschelman 0, 30897304618SKris Buschelman 0, 30997304618SKris Buschelman 0, 31097304618SKris Buschelman /*75*/ 0, 31197304618SKris Buschelman 0, 31297304618SKris Buschelman 0, 31397304618SKris Buschelman 0, 31497304618SKris Buschelman 0, 31597304618SKris Buschelman /*80*/ 0, 31697304618SKris Buschelman 0, 31797304618SKris Buschelman 0, 31897304618SKris Buschelman 0, 31997304618SKris Buschelman 0, 32097304618SKris Buschelman /*85*/ 0 32197304618SKris Buschelman }; 322b97cf49bSBarry Smith 323a23d5eceSKris Buschelman EXTERN_C_BEGIN 324a23d5eceSKris Buschelman #undef __FUNCT__ 325a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIAdjSetPreallocation_MPIAdj" 326a23d5eceSKris Buschelman int MatMPIAdjSetPreallocation_MPIAdj(Mat B,int *i,int *j,int *values) 327a23d5eceSKris Buschelman { 328a23d5eceSKris Buschelman Mat_MPIAdj *b = (Mat_MPIAdj *)B->data; 329a23d5eceSKris Buschelman int ierr; 330a23d5eceSKris Buschelman #if defined(PETSC_USE_BOPT_g) 331a23d5eceSKris Buschelman int ii; 332a23d5eceSKris Buschelman #endif 333a23d5eceSKris Buschelman 334a23d5eceSKris Buschelman PetscFunctionBegin; 335a23d5eceSKris Buschelman B->preallocated = PETSC_TRUE; 336a23d5eceSKris Buschelman #if defined(PETSC_USE_BOPT_g) 337a23d5eceSKris Buschelman if (i[0] != 0) SETERRQ1(1,"First i[] index must be zero, instead it is %d\n",i[0]); 338a23d5eceSKris Buschelman for (ii=1; ii<B->m; ii++) { 339a23d5eceSKris Buschelman if (i[ii] < 0 || i[ii] < i[ii-1]) { 340a23d5eceSKris Buschelman SETERRQ4(1,"i[%d]=%d index is out of range: i[%d]=%d",ii,i[ii],ii-1,i[ii-1]); 341a23d5eceSKris Buschelman } 342a23d5eceSKris Buschelman } 343a23d5eceSKris Buschelman for (ii=0; ii<i[B->m]; ii++) { 344a23d5eceSKris Buschelman if (j[ii] < 0 || j[ii] >= B->N) { 345a23d5eceSKris Buschelman SETERRQ2(1,"Column index %d out of range %d\n",ii,j[ii]); 346a23d5eceSKris Buschelman } 347a23d5eceSKris Buschelman } 348a23d5eceSKris Buschelman #endif 349a23d5eceSKris Buschelman 350a23d5eceSKris Buschelman b->j = j; 351a23d5eceSKris Buschelman b->i = i; 352a23d5eceSKris Buschelman b->values = values; 353a23d5eceSKris Buschelman 354a23d5eceSKris Buschelman b->nz = i[B->m]; 355a23d5eceSKris Buschelman b->diag = 0; 356a23d5eceSKris Buschelman b->symmetric = PETSC_FALSE; 357a23d5eceSKris Buschelman b->freeaij = PETSC_TRUE; 358a23d5eceSKris Buschelman 359a23d5eceSKris Buschelman ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 360a23d5eceSKris Buschelman ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 361a23d5eceSKris Buschelman PetscFunctionReturn(0); 362a23d5eceSKris Buschelman } 363a23d5eceSKris Buschelman EXTERN_C_END 364b97cf49bSBarry Smith 3650bad9183SKris Buschelman /*MC 366fafad747SKris Buschelman MATMPIADJ - MATMPIADJ = "mpiadj" - A matrix type to be used for distributed adjacency matrices, 3670bad9183SKris Buschelman intended for use constructing orderings and partitionings. 3680bad9183SKris Buschelman 3690bad9183SKris Buschelman Level: beginner 3700bad9183SKris Buschelman 3710bad9183SKris Buschelman .seealso: MatCreateMPIAdj 3720bad9183SKris Buschelman M*/ 3730bad9183SKris Buschelman 374273d9f13SBarry Smith EXTERN_C_BEGIN 3754a2ae208SSatish Balay #undef __FUNCT__ 3764a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIAdj" 377273d9f13SBarry Smith int MatCreate_MPIAdj(Mat B) 378273d9f13SBarry Smith { 379273d9f13SBarry Smith Mat_MPIAdj *b; 380273d9f13SBarry Smith int ii,ierr,size,rank; 381273d9f13SBarry Smith 382273d9f13SBarry Smith PetscFunctionBegin; 383273d9f13SBarry Smith ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); 384273d9f13SBarry Smith ierr = MPI_Comm_rank(B->comm,&rank);CHKERRQ(ierr); 385273d9f13SBarry Smith 386b0a32e0cSBarry Smith ierr = PetscNew(Mat_MPIAdj,&b);CHKERRQ(ierr); 387b0a32e0cSBarry Smith B->data = (void*)b; 388273d9f13SBarry Smith ierr = PetscMemzero(b,sizeof(Mat_MPIAdj));CHKERRQ(ierr); 389273d9f13SBarry Smith ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 390273d9f13SBarry Smith B->factor = 0; 391273d9f13SBarry Smith B->lupivotthreshold = 1.0; 392273d9f13SBarry Smith B->mapping = 0; 393273d9f13SBarry Smith B->assembled = PETSC_FALSE; 394273d9f13SBarry Smith 395273d9f13SBarry Smith ierr = MPI_Allreduce(&B->m,&B->M,1,MPI_INT,MPI_SUM,B->comm);CHKERRQ(ierr); 396273d9f13SBarry Smith B->n = B->N; 397273d9f13SBarry Smith 398273d9f13SBarry Smith /* the information in the maps duplicates the information computed below, eventually 399273d9f13SBarry Smith we should remove the duplicate information that is not contained in the maps */ 4008a124369SBarry Smith ierr = PetscMapCreateMPI(B->comm,B->m,B->M,&B->rmap);CHKERRQ(ierr); 401273d9f13SBarry Smith /* we don't know the "local columns" so just use the row information :-(*/ 4028a124369SBarry Smith ierr = PetscMapCreateMPI(B->comm,B->m,B->M,&B->cmap);CHKERRQ(ierr); 403273d9f13SBarry Smith 404b0a32e0cSBarry Smith ierr = PetscMalloc((size+1)*sizeof(int),&b->rowners);CHKERRQ(ierr); 405b0a32e0cSBarry Smith PetscLogObjectMemory(B,(size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIAdj)); 406273d9f13SBarry Smith ierr = MPI_Allgather(&B->m,1,MPI_INT,b->rowners+1,1,MPI_INT,B->comm);CHKERRQ(ierr); 407273d9f13SBarry Smith b->rowners[0] = 0; 408273d9f13SBarry Smith for (ii=2; ii<=size; ii++) { 409273d9f13SBarry Smith b->rowners[ii] += b->rowners[ii-1]; 410273d9f13SBarry Smith } 411273d9f13SBarry Smith b->rstart = b->rowners[rank]; 412273d9f13SBarry Smith b->rend = b->rowners[rank+1]; 413a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMPIAdjSetPreallocation_C", 414a23d5eceSKris Buschelman "MatMPIAdjSetPreallocation_MPIAdj", 415a23d5eceSKris Buschelman MatMPIAdjSetPreallocation_MPIAdj);CHKERRQ(ierr); 416273d9f13SBarry Smith PetscFunctionReturn(0); 417273d9f13SBarry Smith } 418273d9f13SBarry Smith EXTERN_C_END 419273d9f13SBarry Smith 4204a2ae208SSatish Balay #undef __FUNCT__ 4214a2ae208SSatish Balay #define __FUNCT__ "MatMPIAdjSetPreallocation" 422a23d5eceSKris Buschelman /*@C 423a23d5eceSKris Buschelman MatMPIAdjSetPreallocation - Sets the array used for storing the matrix elements 424a23d5eceSKris Buschelman 425a23d5eceSKris Buschelman Collective on MPI_Comm 426a23d5eceSKris Buschelman 427a23d5eceSKris Buschelman Input Parameters: 428a23d5eceSKris Buschelman + A - the matrix 429a23d5eceSKris Buschelman . i - the indices into j for the start of each row 430a23d5eceSKris Buschelman . j - the column indices for each row (sorted for each row). 431a23d5eceSKris Buschelman The indices in i and j start with zero (NOT with one). 432a23d5eceSKris Buschelman - values - [optional] edge weights 433a23d5eceSKris Buschelman 434a23d5eceSKris Buschelman Level: intermediate 435a23d5eceSKris Buschelman 436a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateMPIAdj(), MatSetValues() 437a23d5eceSKris Buschelman @*/ 438273d9f13SBarry Smith int MatMPIAdjSetPreallocation(Mat B,int *i,int *j,int *values) 439273d9f13SBarry Smith { 440a23d5eceSKris Buschelman int ierr,(*f)(Mat,int*,int*,int*); 441273d9f13SBarry Smith 442273d9f13SBarry Smith PetscFunctionBegin; 443a23d5eceSKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)B,"MatMPIAdjSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 444a23d5eceSKris Buschelman if (f) { 445a23d5eceSKris Buschelman ierr = (*f)(B,i,j,values);CHKERRQ(ierr); 446273d9f13SBarry Smith } 447273d9f13SBarry Smith PetscFunctionReturn(0); 448273d9f13SBarry Smith } 449273d9f13SBarry Smith 4504a2ae208SSatish Balay #undef __FUNCT__ 4514a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIAdj" 452b97cf49bSBarry Smith /*@C 4533eda8832SBarry Smith MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list. 454b97cf49bSBarry Smith The matrix does not have numerical values associated with it, but is 455b97cf49bSBarry Smith intended for ordering (to reduce bandwidth etc) and partitioning. 456b97cf49bSBarry Smith 457ef5ee4d1SLois Curfman McInnes Collective on MPI_Comm 458ef5ee4d1SLois Curfman McInnes 459b97cf49bSBarry Smith Input Parameters: 460c2e958feSBarry Smith + comm - MPI communicator 4610752156aSBarry Smith . m - number of local rows 462b97cf49bSBarry Smith . n - number of columns 463b97cf49bSBarry Smith . i - the indices into j for the start of each row 464329f5518SBarry Smith . j - the column indices for each row (sorted for each row). 465ef5ee4d1SLois Curfman McInnes The indices in i and j start with zero (NOT with one). 466329f5518SBarry Smith - values -[optional] edge weights 467b97cf49bSBarry Smith 468b97cf49bSBarry Smith Output Parameter: 469b97cf49bSBarry Smith . A - the matrix 470b97cf49bSBarry Smith 47115091d37SBarry Smith Level: intermediate 47215091d37SBarry Smith 4734bc6d8c0SBarry Smith Notes: This matrix object does not support most matrix operations, include 4744bc6d8c0SBarry Smith MatSetValues(). 475c2e958feSBarry Smith You must NOT free the ii, values and jj arrays yourself. PETSc will free them 476ededeb1aSvictorle when the matrix is destroyed; you must allocate them with PetscMalloc(). If you 4771198db28SBarry Smith call from Fortran you need not create the arrays with PetscMalloc(). 478327e1a73SBarry Smith Should not include the matrix diagonals. 479b97cf49bSBarry Smith 480e3f7e1d6Svictorle If you already have a matrix, you can create its adjacency matrix by a call 481ededeb1aSvictorle to MatConvert, specifying a type of MATMPIADJ. 482ededeb1aSvictorle 483ef5ee4d1SLois Curfman McInnes Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC 484b97cf49bSBarry Smith 485e3f7e1d6Svictorle .seealso: MatCreate(), MatConvert(), MatGetOrdering() 486b97cf49bSBarry Smith @*/ 4873eda8832SBarry Smith int MatCreateMPIAdj(MPI_Comm comm,int m,int n,int *i,int *j,int *values,Mat *A) 488b97cf49bSBarry Smith { 489273d9f13SBarry Smith int ierr; 490b97cf49bSBarry Smith 491433994e6SBarry Smith PetscFunctionBegin; 492273d9f13SBarry Smith ierr = MatCreate(comm,m,n,PETSC_DETERMINE,n,A);CHKERRQ(ierr); 493273d9f13SBarry Smith ierr = MatSetType(*A,MATMPIADJ);CHKERRQ(ierr); 494273d9f13SBarry Smith ierr = MatMPIAdjSetPreallocation(*A,i,j,values);CHKERRQ(ierr); 4953a40ed3dSBarry Smith PetscFunctionReturn(0); 496b97cf49bSBarry Smith } 497b97cf49bSBarry Smith 498273d9f13SBarry Smith EXTERN_C_BEGIN 4994a2ae208SSatish Balay #undef __FUNCT__ 5004a2ae208SSatish Balay #define __FUNCT__ "MatConvertTo_MPIAdj" 501676c34cdSKris Buschelman int MatConvertTo_MPIAdj(Mat A,MatType type,Mat *newmat) { 502676c34cdSKris Buschelman Mat B; 5034c49b128SBarry Smith int i,ierr,m,N,nzeros = 0,*ia,*ja,*rj,len,rstart,cnt,j,*a; 50487828ca2SBarry Smith PetscScalar *ra; 505c2e958feSBarry Smith MPI_Comm comm; 506c2e958feSBarry Smith 507c2e958feSBarry Smith PetscFunctionBegin; 508c2e958feSBarry Smith ierr = MatGetSize(A,PETSC_NULL,&N);CHKERRQ(ierr); 509c2e958feSBarry Smith ierr = MatGetLocalSize(A,&m,PETSC_NULL);CHKERRQ(ierr); 510c2e958feSBarry Smith ierr = MatGetOwnershipRange(A,&rstart,PETSC_NULL);CHKERRQ(ierr); 511c2e958feSBarry Smith 512c2e958feSBarry Smith /* count the number of nonzeros per row */ 513c2e958feSBarry Smith for (i=0; i<m; i++) { 514c2e958feSBarry Smith ierr = MatGetRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr); 515c2e958feSBarry Smith for (j=0; j<len; j++) { 516c2e958feSBarry Smith if (rj[j] == i+rstart) {len--; break;} /* don't count diagonal */ 517c2e958feSBarry Smith } 518c2e958feSBarry Smith ierr = MatRestoreRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr); 519c2e958feSBarry Smith nzeros += len; 520c2e958feSBarry Smith } 521c2e958feSBarry Smith 522c2e958feSBarry Smith /* malloc space for nonzeros */ 523b0a32e0cSBarry Smith ierr = PetscMalloc((nzeros+1)*sizeof(int),&a);CHKERRQ(ierr); 524b0a32e0cSBarry Smith ierr = PetscMalloc((N+1)*sizeof(int),&ia);CHKERRQ(ierr); 525b0a32e0cSBarry Smith ierr = PetscMalloc((nzeros+1)*sizeof(int),&ja);CHKERRQ(ierr); 526c2e958feSBarry Smith 527c2e958feSBarry Smith nzeros = 0; 528c2e958feSBarry Smith ia[0] = 0; 529c2e958feSBarry Smith for (i=0; i<m; i++) { 530c2e958feSBarry Smith ierr = MatGetRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr); 531c2e958feSBarry Smith cnt = 0; 532c2e958feSBarry Smith for (j=0; j<len; j++) { 533c2e958feSBarry Smith if (rj[j] != i+rstart) { /* if not diagonal */ 5344c49b128SBarry Smith a[nzeros+cnt] = (int) PetscAbsScalar(ra[j]); 535c2e958feSBarry Smith ja[nzeros+cnt++] = rj[j]; 536c2e958feSBarry Smith } 537c2e958feSBarry Smith } 538c2e958feSBarry Smith ierr = MatRestoreRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr); 539c2e958feSBarry Smith nzeros += cnt; 540c2e958feSBarry Smith ia[i+1] = nzeros; 541c2e958feSBarry Smith } 542c2e958feSBarry Smith 543c2e958feSBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 544676c34cdSKris Buschelman ierr = MatCreateMPIAdj(comm,m,N,ia,ja,a,&B);CHKERRQ(ierr); 545c2e958feSBarry Smith 546676c34cdSKris Buschelman /* Fake support for "inplace" convert. */ 547676c34cdSKris Buschelman if (*newmat == A) { 548676c34cdSKris Buschelman ierr = MatDestroy(A);CHKERRQ(ierr); 549676c34cdSKris Buschelman } 550676c34cdSKris Buschelman *newmat = B; 551c2e958feSBarry Smith PetscFunctionReturn(0); 552c2e958feSBarry Smith } 553273d9f13SBarry Smith EXTERN_C_END 554433994e6SBarry Smith 555433994e6SBarry Smith 556433994e6SBarry Smith 557433994e6SBarry Smith 558433994e6SBarry Smith 559