1b97cf49bSBarry Smith /* 2b97cf49bSBarry Smith Defines the basic matrix operations for the ADJ adjacency list matrix data-structure. 3b97cf49bSBarry Smith */ 43eda8832SBarry Smith #include "src/mat/impls/adj/mpi/mpiadj.h" 5325e03aeSBarry Smith #include "petscsys.h" 6b97cf49bSBarry Smith 74a2ae208SSatish Balay #undef __FUNCT__ 84a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIAdj_ASCII" 9dfbe8321SBarry Smith PetscErrorCode MatView_MPIAdj_ASCII(Mat A,PetscViewer viewer) 10b97cf49bSBarry Smith { 113eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 12dfbe8321SBarry Smith PetscErrorCode ierr; 13b24ad042SBarry Smith PetscInt i,j,m = A->m; 14fb9695e5SSatish Balay char *name; 15ce1411ecSBarry Smith PetscViewerFormat format; 16b97cf49bSBarry Smith 17433994e6SBarry Smith PetscFunctionBegin; 183a7fca6bSBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 19b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 20fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_INFO) { 213a40ed3dSBarry Smith PetscFunctionReturn(0); 22fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_MATLAB) { 23ffac6cdbSBarry Smith SETERRQ(PETSC_ERR_SUP,"Matlab format not supported"); 240752156aSBarry Smith } else { 25b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 26b97cf49bSBarry Smith for (i=0; i<m; i++) { 2777431f27SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"row %D:",i+a->rstart);CHKERRQ(ierr); 28b97cf49bSBarry Smith for (j=a->i[i]; j<a->i[i+1]; j++) { 2977431f27SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer," %D ",a->j[j]);CHKERRQ(ierr); 300752156aSBarry Smith } 31b0a32e0cSBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr); 32b97cf49bSBarry Smith } 33b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 34b97cf49bSBarry Smith } 35b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 363a40ed3dSBarry Smith PetscFunctionReturn(0); 37b97cf49bSBarry Smith } 38b97cf49bSBarry Smith 394a2ae208SSatish Balay #undef __FUNCT__ 404a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIAdj" 41dfbe8321SBarry Smith PetscErrorCode MatView_MPIAdj(Mat A,PetscViewer viewer) 42b97cf49bSBarry Smith { 43dfbe8321SBarry Smith PetscErrorCode ierr; 4432077d6dSBarry Smith PetscTruth iascii; 45b97cf49bSBarry Smith 46433994e6SBarry Smith PetscFunctionBegin; 4732077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 4832077d6dSBarry Smith if (iascii) { 493eda8832SBarry Smith ierr = MatView_MPIAdj_ASCII(A,viewer);CHKERRQ(ierr); 505cd90555SBarry Smith } else { 5179a5c55eSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by MPIAdj",((PetscObject)viewer)->type_name); 52b97cf49bSBarry Smith } 533a40ed3dSBarry Smith PetscFunctionReturn(0); 54b97cf49bSBarry Smith } 55b97cf49bSBarry Smith 564a2ae208SSatish Balay #undef __FUNCT__ 574a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_MPIAdj" 58dfbe8321SBarry Smith PetscErrorCode MatDestroy_MPIAdj(Mat mat) 59b97cf49bSBarry Smith { 603eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)mat->data; 61dfbe8321SBarry Smith PetscErrorCode ierr; 6294d884c6SBarry Smith 6394d884c6SBarry Smith PetscFunctionBegin; 64aa482453SBarry Smith #if defined(PETSC_USE_LOG) 6577431f27SBarry Smith PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D, NZ=%D",mat->m,mat->n,a->nz); 66b97cf49bSBarry Smith #endif 67606d414cSSatish Balay if (a->diag) {ierr = PetscFree(a->diag);CHKERRQ(ierr);} 680400557aSBarry Smith if (a->freeaij) { 69606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 70606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 711198db28SBarry Smith if (a->values) {ierr = PetscFree(a->values);CHKERRQ(ierr);} 720400557aSBarry Smith } 730400557aSBarry Smith ierr = PetscFree(a->rowners);CHKERRQ(ierr); 741198db28SBarry Smith ierr = PetscFree(a);CHKERRQ(ierr); 75901853e0SKris Buschelman 76901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr); 773a40ed3dSBarry Smith PetscFunctionReturn(0); 78b97cf49bSBarry Smith } 79b97cf49bSBarry Smith 804a2ae208SSatish Balay #undef __FUNCT__ 814a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIAdj" 82dfbe8321SBarry Smith PetscErrorCode MatSetOption_MPIAdj(Mat A,MatOption op) 83b97cf49bSBarry Smith { 843eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 85b97cf49bSBarry Smith 86433994e6SBarry Smith PetscFunctionBegin; 8712c028f9SKris Buschelman switch (op) { 8877e54ba9SKris Buschelman case MAT_SYMMETRIC: 8912c028f9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 909a4540c5SBarry Smith case MAT_HERMITIAN: 91b97cf49bSBarry Smith a->symmetric = PETSC_TRUE; 9212c028f9SKris Buschelman break; 939a4540c5SBarry Smith case MAT_NOT_SYMMETRIC: 949a4540c5SBarry Smith case MAT_NOT_STRUCTURALLY_SYMMETRIC: 959a4540c5SBarry Smith case MAT_NOT_HERMITIAN: 969a4540c5SBarry Smith a->symmetric = PETSC_FALSE; 979a4540c5SBarry Smith break; 989a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 999a4540c5SBarry Smith case MAT_NOT_SYMMETRY_ETERNAL: 1009a4540c5SBarry Smith break; 10112c028f9SKris Buschelman default: 102b0a32e0cSBarry Smith PetscLogInfo(A,"MatSetOption_MPIAdj:Option ignored\n"); 10312c028f9SKris Buschelman break; 104b97cf49bSBarry Smith } 1053a40ed3dSBarry Smith PetscFunctionReturn(0); 106b97cf49bSBarry Smith } 107b97cf49bSBarry Smith 108b97cf49bSBarry Smith 109b97cf49bSBarry Smith /* 110b97cf49bSBarry Smith Adds diagonal pointers to sparse matrix structure. 111b97cf49bSBarry Smith */ 112b97cf49bSBarry Smith 1134a2ae208SSatish Balay #undef __FUNCT__ 1144a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_MPIAdj" 115dfbe8321SBarry Smith PetscErrorCode MatMarkDiagonal_MPIAdj(Mat A) 116b97cf49bSBarry Smith { 1173eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 1186849ba73SBarry Smith PetscErrorCode ierr; 119b24ad042SBarry Smith PetscInt i,j,*diag,m = A->m; 120b97cf49bSBarry Smith 121433994e6SBarry Smith PetscFunctionBegin; 122b24ad042SBarry Smith ierr = PetscMalloc((m+1)*sizeof(PetscInt),&diag);CHKERRQ(ierr); 12352e6d16bSBarry Smith ierr = PetscLogObjectMemory(A,(m+1)*sizeof(PetscInt));CHKERRQ(ierr); 124273d9f13SBarry Smith for (i=0; i<A->m; i++) { 125b97cf49bSBarry Smith for (j=a->i[i]; j<a->i[i+1]; j++) { 126b97cf49bSBarry Smith if (a->j[j] == i) { 127b97cf49bSBarry Smith diag[i] = j; 128b97cf49bSBarry Smith break; 129b97cf49bSBarry Smith } 130b97cf49bSBarry Smith } 131b97cf49bSBarry Smith } 132b97cf49bSBarry Smith a->diag = diag; 1333a40ed3dSBarry Smith PetscFunctionReturn(0); 134b97cf49bSBarry Smith } 135b97cf49bSBarry Smith 1364a2ae208SSatish Balay #undef __FUNCT__ 1374a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_MPIAdj" 138b24ad042SBarry Smith PetscErrorCode MatGetRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 139b97cf49bSBarry Smith { 1403eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 141b24ad042SBarry Smith PetscInt *itmp; 142b97cf49bSBarry Smith 143433994e6SBarry Smith PetscFunctionBegin; 1440752156aSBarry Smith row -= a->rstart; 1450752156aSBarry Smith 146273d9f13SBarry Smith if (row < 0 || row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row out of range"); 147b97cf49bSBarry Smith 148b97cf49bSBarry Smith *nz = a->i[row+1] - a->i[row]; 149b97cf49bSBarry Smith if (v) *v = PETSC_NULL; 150b97cf49bSBarry Smith if (idx) { 151b97cf49bSBarry Smith itmp = a->j + a->i[row]; 152b97cf49bSBarry Smith if (*nz) { 153b97cf49bSBarry Smith *idx = itmp; 154b97cf49bSBarry Smith } 155b97cf49bSBarry Smith else *idx = 0; 156b97cf49bSBarry Smith } 1573a40ed3dSBarry Smith PetscFunctionReturn(0); 158b97cf49bSBarry Smith } 159b97cf49bSBarry Smith 1604a2ae208SSatish Balay #undef __FUNCT__ 1614a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_MPIAdj" 162b24ad042SBarry Smith PetscErrorCode MatRestoreRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 163b97cf49bSBarry Smith { 164433994e6SBarry Smith PetscFunctionBegin; 1653a40ed3dSBarry Smith PetscFunctionReturn(0); 166b97cf49bSBarry Smith } 167b97cf49bSBarry Smith 1684a2ae208SSatish Balay #undef __FUNCT__ 1694a2ae208SSatish Balay #define __FUNCT__ "MatEqual_MPIAdj" 170dfbe8321SBarry Smith PetscErrorCode MatEqual_MPIAdj(Mat A,Mat B,PetscTruth* flg) 171b97cf49bSBarry Smith { 1723eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data,*b = (Mat_MPIAdj *)B->data; 173dfbe8321SBarry Smith PetscErrorCode 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 */ 183b24ad042SBarry Smith ierr = PetscMemcmp(a->i,b->i,(A->m+1)*sizeof(PetscInt),&flag);CHKERRQ(ierr); 184b97cf49bSBarry Smith 185b97cf49bSBarry Smith /* if a->j are the same */ 186b24ad042SBarry Smith ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),&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" 194b24ad042SBarry Smith PetscErrorCode MatGetRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscTruth symmetric,PetscInt *m,PetscInt *ia[],PetscInt *ja[],PetscTruth *done) 1956cd91f64SBarry Smith { 196dfbe8321SBarry Smith PetscErrorCode ierr; 197b24ad042SBarry Smith PetscMPIInt size; 198b24ad042SBarry Smith PetscInt i; 1996cd91f64SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data; 2006cd91f64SBarry Smith 2016cd91f64SBarry Smith PetscFunctionBegin; 2026cd91f64SBarry Smith ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); 2036cd91f64SBarry Smith if (size > 1) {*done = PETSC_FALSE; PetscFunctionReturn(0);} 2046cd91f64SBarry Smith *m = A->m; 2056cd91f64SBarry Smith *ia = a->i; 2066cd91f64SBarry Smith *ja = a->j; 2076cd91f64SBarry Smith *done = PETSC_TRUE; 208d42735a1SBarry Smith if (oshift) { 209d42735a1SBarry Smith for (i=0; i<(*ia)[*m]; i++) { 210d42735a1SBarry Smith (*ja)[i]++; 211d42735a1SBarry Smith } 212d42735a1SBarry Smith for (i=0; i<=(*m); i++) (*ia)[i]++; 213d42735a1SBarry Smith } 214d42735a1SBarry Smith PetscFunctionReturn(0); 215d42735a1SBarry Smith } 216d42735a1SBarry Smith 2174a2ae208SSatish Balay #undef __FUNCT__ 2184a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_MPIAdj" 219b24ad042SBarry Smith PetscErrorCode MatRestoreRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscTruth symmetric,PetscInt *m,PetscInt *ia[],PetscInt *ja[],PetscTruth *done) 220d42735a1SBarry Smith { 221b24ad042SBarry Smith PetscInt i; 222d42735a1SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data; 223d42735a1SBarry Smith 224d42735a1SBarry Smith PetscFunctionBegin; 225e005ede5SBarry Smith if (ia && a->i != *ia) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"ia passed back is not one obtained with MatGetRowIJ()"); 226e005ede5SBarry Smith if (ja && a->j != *ja) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"ja passed back is not one obtained with MatGetRowIJ()"); 227d42735a1SBarry Smith if (oshift) { 228d42735a1SBarry Smith for (i=0; i<=(*m); i++) (*ia)[i]--; 229d42735a1SBarry Smith for (i=0; i<(*ia)[*m]; i++) { 230d42735a1SBarry Smith (*ja)[i]--; 231d42735a1SBarry Smith } 232d42735a1SBarry Smith } 2336cd91f64SBarry Smith PetscFunctionReturn(0); 2346cd91f64SBarry Smith } 235b97cf49bSBarry Smith 236b97cf49bSBarry Smith /* -------------------------------------------------------------------*/ 2370b3b1487SBarry Smith static struct _MatOps MatOps_Values = {0, 2383eda8832SBarry Smith MatGetRow_MPIAdj, 2393eda8832SBarry Smith MatRestoreRow_MPIAdj, 240b97cf49bSBarry Smith 0, 24197304618SKris Buschelman /* 4*/ 0, 242b97cf49bSBarry Smith 0, 243b97cf49bSBarry Smith 0, 244b97cf49bSBarry Smith 0, 2450b3b1487SBarry Smith 0, 2460b3b1487SBarry Smith 0, 24797304618SKris Buschelman /*10*/ 0, 2480b3b1487SBarry Smith 0, 2490b3b1487SBarry Smith 0, 2500b3b1487SBarry Smith 0, 2510b3b1487SBarry Smith 0, 25297304618SKris Buschelman /*15*/ 0, 2533eda8832SBarry Smith MatEqual_MPIAdj, 2540b3b1487SBarry Smith 0, 2550b3b1487SBarry Smith 0, 2560b3b1487SBarry Smith 0, 25797304618SKris Buschelman /*20*/ 0, 2580b3b1487SBarry Smith 0, 2590b3b1487SBarry Smith 0, 2603eda8832SBarry Smith MatSetOption_MPIAdj, 2610b3b1487SBarry Smith 0, 26297304618SKris Buschelman /*25*/ 0, 2630b3b1487SBarry Smith 0, 2640b3b1487SBarry Smith 0, 2650b3b1487SBarry Smith 0, 2660b3b1487SBarry Smith 0, 26797304618SKris Buschelman /*30*/ 0, 2680b3b1487SBarry Smith 0, 269273d9f13SBarry Smith 0, 270273d9f13SBarry Smith 0, 2710b3b1487SBarry Smith 0, 27297304618SKris Buschelman /*35*/ 0, 2730b3b1487SBarry Smith 0, 2740b3b1487SBarry Smith 0, 2750b3b1487SBarry Smith 0, 2760b3b1487SBarry Smith 0, 27797304618SKris Buschelman /*40*/ 0, 2780b3b1487SBarry Smith 0, 2790b3b1487SBarry Smith 0, 2800b3b1487SBarry Smith 0, 2810b3b1487SBarry Smith 0, 28297304618SKris Buschelman /*45*/ 0, 2830b3b1487SBarry Smith 0, 2840b3b1487SBarry Smith 0, 2850b3b1487SBarry Smith 0, 2860b3b1487SBarry Smith 0, 287521d7252SBarry Smith /*50*/ 0, 2886cd91f64SBarry Smith MatGetRowIJ_MPIAdj, 289d42735a1SBarry Smith MatRestoreRowIJ_MPIAdj, 290b97cf49bSBarry Smith 0, 291b97cf49bSBarry Smith 0, 29297304618SKris Buschelman /*55*/ 0, 293b97cf49bSBarry Smith 0, 2940752156aSBarry Smith 0, 2950752156aSBarry Smith 0, 2960b3b1487SBarry Smith 0, 29797304618SKris Buschelman /*60*/ 0, 298b9b97703SBarry Smith MatDestroy_MPIAdj, 299b9b97703SBarry Smith MatView_MPIAdj, 30097304618SKris Buschelman MatGetPetscMaps_Petsc, 30197304618SKris Buschelman 0, 30297304618SKris Buschelman /*65*/ 0, 30397304618SKris Buschelman 0, 30497304618SKris Buschelman 0, 30597304618SKris Buschelman 0, 30697304618SKris Buschelman 0, 30797304618SKris Buschelman /*70*/ 0, 30897304618SKris Buschelman 0, 30997304618SKris Buschelman 0, 31097304618SKris Buschelman 0, 31197304618SKris Buschelman 0, 31297304618SKris Buschelman /*75*/ 0, 31397304618SKris Buschelman 0, 31497304618SKris Buschelman 0, 31597304618SKris Buschelman 0, 31697304618SKris Buschelman 0, 31797304618SKris Buschelman /*80*/ 0, 31897304618SKris Buschelman 0, 31997304618SKris Buschelman 0, 32097304618SKris Buschelman 0, 321865e5f61SKris Buschelman 0, 322865e5f61SKris Buschelman /*85*/ 0, 323865e5f61SKris Buschelman 0, 324865e5f61SKris Buschelman 0, 325865e5f61SKris Buschelman 0, 326865e5f61SKris Buschelman 0, 327865e5f61SKris Buschelman /*90*/ 0, 328865e5f61SKris Buschelman 0, 329865e5f61SKris Buschelman 0, 330865e5f61SKris Buschelman 0, 331865e5f61SKris Buschelman 0, 332865e5f61SKris Buschelman /*95*/ 0, 333865e5f61SKris Buschelman 0, 334865e5f61SKris Buschelman 0, 335865e5f61SKris Buschelman 0}; 336b97cf49bSBarry Smith 337a23d5eceSKris Buschelman EXTERN_C_BEGIN 338a23d5eceSKris Buschelman #undef __FUNCT__ 339a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIAdjSetPreallocation_MPIAdj" 340b24ad042SBarry Smith PetscErrorCode MatMPIAdjSetPreallocation_MPIAdj(Mat B,PetscInt *i,PetscInt *j,PetscInt *values) 341a23d5eceSKris Buschelman { 342a23d5eceSKris Buschelman Mat_MPIAdj *b = (Mat_MPIAdj *)B->data; 343dfbe8321SBarry Smith PetscErrorCode ierr; 3442515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 345b24ad042SBarry Smith PetscInt ii; 346a23d5eceSKris Buschelman #endif 347a23d5eceSKris Buschelman 348a23d5eceSKris Buschelman PetscFunctionBegin; 349a23d5eceSKris Buschelman B->preallocated = PETSC_TRUE; 3502515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 35177431f27SBarry Smith if (i[0] != 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"First i[] index must be zero, instead it is %D\n",i[0]); 352a23d5eceSKris Buschelman for (ii=1; ii<B->m; ii++) { 353a23d5eceSKris Buschelman if (i[ii] < 0 || i[ii] < i[ii-1]) { 35477431f27SBarry Smith SETERRQ4(PETSC_ERR_ARG_OUTOFRANGE,"i[%D]=%D index is out of range: i[%D]=%D",ii,i[ii],ii-1,i[ii-1]); 355a23d5eceSKris Buschelman } 356a23d5eceSKris Buschelman } 357a23d5eceSKris Buschelman for (ii=0; ii<i[B->m]; ii++) { 358a23d5eceSKris Buschelman if (j[ii] < 0 || j[ii] >= B->N) { 35977431f27SBarry Smith SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column index %D out of range %D\n",ii,j[ii]); 360a23d5eceSKris Buschelman } 361a23d5eceSKris Buschelman } 362a23d5eceSKris Buschelman #endif 363a23d5eceSKris Buschelman 364a23d5eceSKris Buschelman b->j = j; 365a23d5eceSKris Buschelman b->i = i; 366a23d5eceSKris Buschelman b->values = values; 367a23d5eceSKris Buschelman 368a23d5eceSKris Buschelman b->nz = i[B->m]; 369a23d5eceSKris Buschelman b->diag = 0; 370a23d5eceSKris Buschelman b->symmetric = PETSC_FALSE; 371a23d5eceSKris Buschelman b->freeaij = PETSC_TRUE; 372a23d5eceSKris Buschelman 373a23d5eceSKris Buschelman ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 374a23d5eceSKris Buschelman ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 375a23d5eceSKris Buschelman PetscFunctionReturn(0); 376a23d5eceSKris Buschelman } 377a23d5eceSKris Buschelman EXTERN_C_END 378b97cf49bSBarry Smith 3790bad9183SKris Buschelman /*MC 380fafad747SKris Buschelman MATMPIADJ - MATMPIADJ = "mpiadj" - A matrix type to be used for distributed adjacency matrices, 3810bad9183SKris Buschelman intended for use constructing orderings and partitionings. 3820bad9183SKris Buschelman 3830bad9183SKris Buschelman Level: beginner 3840bad9183SKris Buschelman 3850bad9183SKris Buschelman .seealso: MatCreateMPIAdj 3860bad9183SKris Buschelman M*/ 3870bad9183SKris Buschelman 388273d9f13SBarry Smith EXTERN_C_BEGIN 3894a2ae208SSatish Balay #undef __FUNCT__ 3904a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIAdj" 391dfbe8321SBarry Smith PetscErrorCode MatCreate_MPIAdj(Mat B) 392273d9f13SBarry Smith { 393273d9f13SBarry Smith Mat_MPIAdj *b; 3946849ba73SBarry Smith PetscErrorCode ierr; 395b24ad042SBarry Smith PetscInt ii; 396b24ad042SBarry Smith PetscMPIInt size,rank; 397273d9f13SBarry Smith 398273d9f13SBarry Smith PetscFunctionBegin; 399273d9f13SBarry Smith ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); 400273d9f13SBarry Smith ierr = MPI_Comm_rank(B->comm,&rank);CHKERRQ(ierr); 401273d9f13SBarry Smith 402b0a32e0cSBarry Smith ierr = PetscNew(Mat_MPIAdj,&b);CHKERRQ(ierr); 403b0a32e0cSBarry Smith B->data = (void*)b; 404273d9f13SBarry Smith ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 405273d9f13SBarry Smith B->factor = 0; 406273d9f13SBarry Smith B->lupivotthreshold = 1.0; 407273d9f13SBarry Smith B->mapping = 0; 408273d9f13SBarry Smith B->assembled = PETSC_FALSE; 409273d9f13SBarry Smith 410ae77e121SBarry Smith ierr = PetscSplitOwnership(B->comm,&B->m,&B->M);CHKERRQ(ierr); 411ae77e121SBarry Smith B->n = B->N = PetscMax(B->N,B->n);CHKERRQ(ierr); 412273d9f13SBarry Smith 413273d9f13SBarry Smith /* the information in the maps duplicates the information computed below, eventually 414273d9f13SBarry Smith we should remove the duplicate information that is not contained in the maps */ 4158a124369SBarry Smith ierr = PetscMapCreateMPI(B->comm,B->m,B->M,&B->rmap);CHKERRQ(ierr); 416273d9f13SBarry Smith /* we don't know the "local columns" so just use the row information :-(*/ 4178a124369SBarry Smith ierr = PetscMapCreateMPI(B->comm,B->m,B->M,&B->cmap);CHKERRQ(ierr); 418273d9f13SBarry Smith 419b24ad042SBarry Smith ierr = PetscMalloc((size+1)*sizeof(PetscInt),&b->rowners);CHKERRQ(ierr); 42052e6d16bSBarry Smith ierr = PetscLogObjectMemory(B,(size+2)*sizeof(PetscInt)+sizeof(struct _p_Mat)+sizeof(Mat_MPIAdj));CHKERRQ(ierr); 421273d9f13SBarry Smith ierr = MPI_Allgather(&B->m,1,MPI_INT,b->rowners+1,1,MPI_INT,B->comm);CHKERRQ(ierr); 422273d9f13SBarry Smith b->rowners[0] = 0; 423273d9f13SBarry Smith for (ii=2; ii<=size; ii++) { 424273d9f13SBarry Smith b->rowners[ii] += b->rowners[ii-1]; 425273d9f13SBarry Smith } 426273d9f13SBarry Smith b->rstart = b->rowners[rank]; 427273d9f13SBarry Smith b->rend = b->rowners[rank+1]; 428a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMPIAdjSetPreallocation_C", 429a23d5eceSKris Buschelman "MatMPIAdjSetPreallocation_MPIAdj", 430a23d5eceSKris Buschelman MatMPIAdjSetPreallocation_MPIAdj);CHKERRQ(ierr); 431273d9f13SBarry Smith PetscFunctionReturn(0); 432273d9f13SBarry Smith } 433273d9f13SBarry Smith EXTERN_C_END 434273d9f13SBarry Smith 4354a2ae208SSatish Balay #undef __FUNCT__ 4364a2ae208SSatish Balay #define __FUNCT__ "MatMPIAdjSetPreallocation" 437a23d5eceSKris Buschelman /*@C 438a23d5eceSKris Buschelman MatMPIAdjSetPreallocation - Sets the array used for storing the matrix elements 439a23d5eceSKris Buschelman 440a23d5eceSKris Buschelman Collective on MPI_Comm 441a23d5eceSKris Buschelman 442a23d5eceSKris Buschelman Input Parameters: 443a23d5eceSKris Buschelman + A - the matrix 444a23d5eceSKris Buschelman . i - the indices into j for the start of each row 445a23d5eceSKris Buschelman . j - the column indices for each row (sorted for each row). 446a23d5eceSKris Buschelman The indices in i and j start with zero (NOT with one). 447a23d5eceSKris Buschelman - values - [optional] edge weights 448a23d5eceSKris Buschelman 449a23d5eceSKris Buschelman Level: intermediate 450a23d5eceSKris Buschelman 451a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateMPIAdj(), MatSetValues() 452a23d5eceSKris Buschelman @*/ 453b24ad042SBarry Smith PetscErrorCode MatMPIAdjSetPreallocation(Mat B,PetscInt *i,PetscInt *j,PetscInt *values) 454273d9f13SBarry Smith { 455b24ad042SBarry Smith PetscErrorCode ierr,(*f)(Mat,PetscInt*,PetscInt*,PetscInt*); 456273d9f13SBarry Smith 457273d9f13SBarry Smith PetscFunctionBegin; 458a23d5eceSKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)B,"MatMPIAdjSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 459a23d5eceSKris Buschelman if (f) { 460a23d5eceSKris Buschelman ierr = (*f)(B,i,j,values);CHKERRQ(ierr); 461273d9f13SBarry Smith } 462273d9f13SBarry Smith PetscFunctionReturn(0); 463273d9f13SBarry Smith } 464273d9f13SBarry Smith 4654a2ae208SSatish Balay #undef __FUNCT__ 4664a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIAdj" 467b97cf49bSBarry Smith /*@C 4683eda8832SBarry Smith MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list. 469b97cf49bSBarry Smith The matrix does not have numerical values associated with it, but is 470b97cf49bSBarry Smith intended for ordering (to reduce bandwidth etc) and partitioning. 471b97cf49bSBarry Smith 472ef5ee4d1SLois Curfman McInnes Collective on MPI_Comm 473ef5ee4d1SLois Curfman McInnes 474b97cf49bSBarry Smith Input Parameters: 475c2e958feSBarry Smith + comm - MPI communicator 4760752156aSBarry Smith . m - number of local rows 477b97cf49bSBarry Smith . n - number of columns 478b97cf49bSBarry Smith . i - the indices into j for the start of each row 479329f5518SBarry Smith . j - the column indices for each row (sorted for each row). 480ef5ee4d1SLois Curfman McInnes The indices in i and j start with zero (NOT with one). 481329f5518SBarry Smith - values -[optional] edge weights 482b97cf49bSBarry Smith 483b97cf49bSBarry Smith Output Parameter: 484b97cf49bSBarry Smith . A - the matrix 485b97cf49bSBarry Smith 48615091d37SBarry Smith Level: intermediate 48715091d37SBarry Smith 4884bc6d8c0SBarry Smith Notes: This matrix object does not support most matrix operations, include 4894bc6d8c0SBarry Smith MatSetValues(). 490c2e958feSBarry Smith You must NOT free the ii, values and jj arrays yourself. PETSc will free them 491ededeb1aSvictorle when the matrix is destroyed; you must allocate them with PetscMalloc(). If you 4921198db28SBarry Smith call from Fortran you need not create the arrays with PetscMalloc(). 493327e1a73SBarry Smith Should not include the matrix diagonals. 494b97cf49bSBarry Smith 495e3f7e1d6Svictorle If you already have a matrix, you can create its adjacency matrix by a call 496ededeb1aSvictorle to MatConvert, specifying a type of MATMPIADJ. 497ededeb1aSvictorle 498ef5ee4d1SLois Curfman McInnes Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC 499b97cf49bSBarry Smith 500e3f7e1d6Svictorle .seealso: MatCreate(), MatConvert(), MatGetOrdering() 501b97cf49bSBarry Smith @*/ 502b24ad042SBarry Smith PetscErrorCode MatCreateMPIAdj(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt *i,PetscInt *j,PetscInt *values,Mat *A) 503b97cf49bSBarry Smith { 504dfbe8321SBarry Smith PetscErrorCode ierr; 505b97cf49bSBarry Smith 506433994e6SBarry Smith PetscFunctionBegin; 507273d9f13SBarry Smith ierr = MatCreate(comm,m,n,PETSC_DETERMINE,n,A);CHKERRQ(ierr); 508273d9f13SBarry Smith ierr = MatSetType(*A,MATMPIADJ);CHKERRQ(ierr); 509273d9f13SBarry Smith ierr = MatMPIAdjSetPreallocation(*A,i,j,values);CHKERRQ(ierr); 5103a40ed3dSBarry Smith PetscFunctionReturn(0); 511b97cf49bSBarry Smith } 512b97cf49bSBarry Smith 513273d9f13SBarry Smith EXTERN_C_BEGIN 5144a2ae208SSatish Balay #undef __FUNCT__ 5154a2ae208SSatish Balay #define __FUNCT__ "MatConvertTo_MPIAdj" 516ceb03754SKris Buschelman PetscErrorCode MatConvertTo_MPIAdj(Mat A,MatType type,MatReuse reuse,Mat *newmat) 517b3cc6726SBarry Smith { 518676c34cdSKris Buschelman Mat B; 5196849ba73SBarry Smith PetscErrorCode ierr; 520b24ad042SBarry Smith PetscInt i,m,N,nzeros = 0,*ia,*ja,len,rstart,cnt,j,*a; 521b24ad042SBarry Smith const PetscInt *rj; 522b3cc6726SBarry Smith const PetscScalar *ra; 523c2e958feSBarry Smith MPI_Comm comm; 524c2e958feSBarry Smith 525c2e958feSBarry Smith PetscFunctionBegin; 526c2e958feSBarry Smith ierr = MatGetSize(A,PETSC_NULL,&N);CHKERRQ(ierr); 527c2e958feSBarry Smith ierr = MatGetLocalSize(A,&m,PETSC_NULL);CHKERRQ(ierr); 528c2e958feSBarry Smith ierr = MatGetOwnershipRange(A,&rstart,PETSC_NULL);CHKERRQ(ierr); 529c2e958feSBarry Smith 530c2e958feSBarry Smith /* count the number of nonzeros per row */ 531c2e958feSBarry Smith for (i=0; i<m; i++) { 532c2e958feSBarry Smith ierr = MatGetRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr); 533c2e958feSBarry Smith for (j=0; j<len; j++) { 534c2e958feSBarry Smith if (rj[j] == i+rstart) {len--; break;} /* don't count diagonal */ 535c2e958feSBarry Smith } 536c2e958feSBarry Smith ierr = MatRestoreRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr); 537c2e958feSBarry Smith nzeros += len; 538c2e958feSBarry Smith } 539c2e958feSBarry Smith 540c2e958feSBarry Smith /* malloc space for nonzeros */ 541b24ad042SBarry Smith ierr = PetscMalloc((nzeros+1)*sizeof(PetscInt),&a);CHKERRQ(ierr); 542b24ad042SBarry Smith ierr = PetscMalloc((N+1)*sizeof(PetscInt),&ia);CHKERRQ(ierr); 543b24ad042SBarry Smith ierr = PetscMalloc((nzeros+1)*sizeof(PetscInt),&ja);CHKERRQ(ierr); 544c2e958feSBarry Smith 545c2e958feSBarry Smith nzeros = 0; 546c2e958feSBarry Smith ia[0] = 0; 547c2e958feSBarry Smith for (i=0; i<m; i++) { 548c2e958feSBarry Smith ierr = MatGetRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr); 549c2e958feSBarry Smith cnt = 0; 550c2e958feSBarry Smith for (j=0; j<len; j++) { 551c2e958feSBarry Smith if (rj[j] != i+rstart) { /* if not diagonal */ 552b24ad042SBarry Smith a[nzeros+cnt] = (PetscInt) PetscAbsScalar(ra[j]); 553c2e958feSBarry Smith ja[nzeros+cnt++] = rj[j]; 554c2e958feSBarry Smith } 555c2e958feSBarry Smith } 556c2e958feSBarry Smith ierr = MatRestoreRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr); 557c2e958feSBarry Smith nzeros += cnt; 558c2e958feSBarry Smith ia[i+1] = nzeros; 559c2e958feSBarry Smith } 560c2e958feSBarry Smith 561c2e958feSBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 562f204ca49SKris Buschelman ierr = MatCreate(comm,m,N,PETSC_DETERMINE,N,&B);CHKERRQ(ierr); 563f204ca49SKris Buschelman ierr = MatSetType(B,type);CHKERRQ(ierr); 564f204ca49SKris Buschelman ierr = MatMPIAdjSetPreallocation(B,ia,ja,a);CHKERRQ(ierr); 565c2e958feSBarry Smith 566ceb03754SKris Buschelman if (reuse == MAT_REUSE_MATRIX) { 567*c3d102feSKris Buschelman ierr = MatHeaderCopy(A,B);CHKERRQ(ierr); 568*c3d102feSKris Buschelman } else { 569676c34cdSKris Buschelman *newmat = B; 570*c3d102feSKris Buschelman } 571c2e958feSBarry Smith PetscFunctionReturn(0); 572c2e958feSBarry Smith } 573273d9f13SBarry Smith EXTERN_C_END 574433994e6SBarry Smith 575433994e6SBarry Smith 576433994e6SBarry Smith 577433994e6SBarry Smith 578433994e6SBarry Smith 579