1be1d678aSKris Buschelman #define PETSCMAT_DLL 2be1d678aSKris Buschelman 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" 11dfbe8321SBarry Smith PetscErrorCode MatView_MPIAdj_ASCII(Mat A,PetscViewer viewer) 12b97cf49bSBarry Smith { 133eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 14dfbe8321SBarry Smith PetscErrorCode ierr; 15899cda47SBarry Smith PetscInt i,j,m = A->rmap.n; 162dcb1b2aSMatthew Knepley const char *name; 17ce1411ecSBarry Smith PetscViewerFormat format; 18b97cf49bSBarry Smith 19433994e6SBarry Smith PetscFunctionBegin; 203a7fca6bSBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 21b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 22fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_INFO) { 233a40ed3dSBarry Smith PetscFunctionReturn(0); 24fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_MATLAB) { 25ffac6cdbSBarry Smith SETERRQ(PETSC_ERR_SUP,"Matlab format not supported"); 260752156aSBarry Smith } else { 27b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 28b97cf49bSBarry Smith for (i=0; i<m; i++) { 29899cda47SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"row %D:",i+A->rmap.rstart);CHKERRQ(ierr); 30b97cf49bSBarry Smith for (j=a->i[i]; j<a->i[i+1]; j++) { 3177431f27SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer," %D ",a->j[j]);CHKERRQ(ierr); 320752156aSBarry Smith } 33b0a32e0cSBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr); 34b97cf49bSBarry Smith } 35b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 36b97cf49bSBarry Smith } 37b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 383a40ed3dSBarry Smith PetscFunctionReturn(0); 39b97cf49bSBarry Smith } 40b97cf49bSBarry Smith 414a2ae208SSatish Balay #undef __FUNCT__ 424a2ae208SSatish Balay #define __FUNCT__ "MatView_MPIAdj" 43dfbe8321SBarry Smith PetscErrorCode MatView_MPIAdj(Mat A,PetscViewer viewer) 44b97cf49bSBarry Smith { 45dfbe8321SBarry Smith PetscErrorCode ierr; 4632077d6dSBarry Smith PetscTruth iascii; 47b97cf49bSBarry Smith 48433994e6SBarry Smith PetscFunctionBegin; 4932077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 5032077d6dSBarry Smith if (iascii) { 513eda8832SBarry Smith ierr = MatView_MPIAdj_ASCII(A,viewer);CHKERRQ(ierr); 525cd90555SBarry Smith } else { 5379a5c55eSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by MPIAdj",((PetscObject)viewer)->type_name); 54b97cf49bSBarry Smith } 553a40ed3dSBarry Smith PetscFunctionReturn(0); 56b97cf49bSBarry Smith } 57b97cf49bSBarry Smith 584a2ae208SSatish Balay #undef __FUNCT__ 594a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_MPIAdj" 60dfbe8321SBarry Smith PetscErrorCode MatDestroy_MPIAdj(Mat mat) 61b97cf49bSBarry Smith { 623eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)mat->data; 63dfbe8321SBarry Smith PetscErrorCode ierr; 6494d884c6SBarry Smith 6594d884c6SBarry Smith PetscFunctionBegin; 66aa482453SBarry Smith #if defined(PETSC_USE_LOG) 67899cda47SBarry Smith PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D, NZ=%D",mat->rmap.n,mat->cmap.n,a->nz); 68b97cf49bSBarry Smith #endif 6905b42c5fSBarry Smith ierr = PetscFree(a->diag);CHKERRQ(ierr); 700400557aSBarry Smith if (a->freeaij) { 71606d414cSSatish Balay ierr = PetscFree(a->i);CHKERRQ(ierr); 72606d414cSSatish Balay ierr = PetscFree(a->j);CHKERRQ(ierr); 7305b42c5fSBarry Smith ierr = PetscFree(a->values);CHKERRQ(ierr); 740400557aSBarry Smith } 751198db28SBarry Smith ierr = PetscFree(a);CHKERRQ(ierr); 76901853e0SKris Buschelman 77901853e0SKris Buschelman ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr); 783a40ed3dSBarry Smith PetscFunctionReturn(0); 79b97cf49bSBarry Smith } 80b97cf49bSBarry Smith 814a2ae208SSatish Balay #undef __FUNCT__ 824a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_MPIAdj" 83dfbe8321SBarry Smith PetscErrorCode MatSetOption_MPIAdj(Mat A,MatOption op) 84b97cf49bSBarry Smith { 853eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 8663ba0a88SBarry Smith PetscErrorCode ierr; 87b97cf49bSBarry Smith 88433994e6SBarry Smith PetscFunctionBegin; 8912c028f9SKris Buschelman switch (op) { 9077e54ba9SKris Buschelman case MAT_SYMMETRIC: 9112c028f9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 929a4540c5SBarry Smith case MAT_HERMITIAN: 93b97cf49bSBarry Smith a->symmetric = PETSC_TRUE; 9412c028f9SKris Buschelman break; 959a4540c5SBarry Smith case MAT_NOT_SYMMETRIC: 969a4540c5SBarry Smith case MAT_NOT_STRUCTURALLY_SYMMETRIC: 979a4540c5SBarry Smith case MAT_NOT_HERMITIAN: 989a4540c5SBarry Smith a->symmetric = PETSC_FALSE; 999a4540c5SBarry Smith break; 1009a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 1019a4540c5SBarry Smith case MAT_NOT_SYMMETRY_ETERNAL: 1029a4540c5SBarry Smith break; 10312c028f9SKris Buschelman default: 104ae15b995SBarry Smith ierr = PetscInfo(A,"Option ignored\n");CHKERRQ(ierr); 10512c028f9SKris Buschelman break; 106b97cf49bSBarry Smith } 1073a40ed3dSBarry Smith PetscFunctionReturn(0); 108b97cf49bSBarry Smith } 109b97cf49bSBarry Smith 110b97cf49bSBarry Smith 111b97cf49bSBarry Smith /* 112b97cf49bSBarry Smith Adds diagonal pointers to sparse matrix structure. 113b97cf49bSBarry Smith */ 114b97cf49bSBarry Smith 1154a2ae208SSatish Balay #undef __FUNCT__ 1164a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_MPIAdj" 117dfbe8321SBarry Smith PetscErrorCode MatMarkDiagonal_MPIAdj(Mat A) 118b97cf49bSBarry Smith { 1193eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 1206849ba73SBarry Smith PetscErrorCode ierr; 121*09f38230SBarry Smith PetscInt i,j,m = A->rmap.n; 122b97cf49bSBarry Smith 123433994e6SBarry Smith PetscFunctionBegin; 124*09f38230SBarry Smith ierr = PetscMalloc(m*sizeof(PetscInt),&a->diag);CHKERRQ(ierr); 125*09f38230SBarry Smith ierr = PetscLogObjectMemory(A,m*sizeof(PetscInt));CHKERRQ(ierr); 126899cda47SBarry Smith for (i=0; i<A->rmap.n; i++) { 127b97cf49bSBarry Smith for (j=a->i[i]; j<a->i[i+1]; j++) { 128b97cf49bSBarry Smith if (a->j[j] == i) { 129*09f38230SBarry Smith a->diag[i] = j; 130b97cf49bSBarry Smith break; 131b97cf49bSBarry Smith } 132b97cf49bSBarry Smith } 133b97cf49bSBarry Smith } 1343a40ed3dSBarry Smith PetscFunctionReturn(0); 135b97cf49bSBarry Smith } 136b97cf49bSBarry Smith 1374a2ae208SSatish Balay #undef __FUNCT__ 1384a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_MPIAdj" 139b24ad042SBarry Smith PetscErrorCode MatGetRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 140b97cf49bSBarry Smith { 1413eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 142b24ad042SBarry Smith PetscInt *itmp; 143b97cf49bSBarry Smith 144433994e6SBarry Smith PetscFunctionBegin; 145899cda47SBarry Smith row -= A->rmap.rstart; 1460752156aSBarry Smith 147899cda47SBarry Smith if (row < 0 || row >= A->rmap.n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row out of range"); 148b97cf49bSBarry Smith 149b97cf49bSBarry Smith *nz = a->i[row+1] - a->i[row]; 150b97cf49bSBarry Smith if (v) *v = PETSC_NULL; 151b97cf49bSBarry Smith if (idx) { 152b97cf49bSBarry Smith itmp = a->j + a->i[row]; 153b97cf49bSBarry Smith if (*nz) { 154b97cf49bSBarry Smith *idx = itmp; 155b97cf49bSBarry Smith } 156b97cf49bSBarry Smith else *idx = 0; 157b97cf49bSBarry Smith } 1583a40ed3dSBarry Smith PetscFunctionReturn(0); 159b97cf49bSBarry Smith } 160b97cf49bSBarry Smith 1614a2ae208SSatish Balay #undef __FUNCT__ 1624a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_MPIAdj" 163b24ad042SBarry Smith PetscErrorCode MatRestoreRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 164b97cf49bSBarry Smith { 165433994e6SBarry Smith PetscFunctionBegin; 1663a40ed3dSBarry Smith PetscFunctionReturn(0); 167b97cf49bSBarry Smith } 168b97cf49bSBarry Smith 1694a2ae208SSatish Balay #undef __FUNCT__ 1704a2ae208SSatish Balay #define __FUNCT__ "MatEqual_MPIAdj" 171dfbe8321SBarry Smith PetscErrorCode MatEqual_MPIAdj(Mat A,Mat B,PetscTruth* flg) 172b97cf49bSBarry Smith { 1733eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data,*b = (Mat_MPIAdj *)B->data; 174dfbe8321SBarry Smith PetscErrorCode ierr; 1750f5bd95cSBarry Smith PetscTruth flag; 176b97cf49bSBarry Smith 177433994e6SBarry Smith PetscFunctionBegin; 178b97cf49bSBarry Smith /* If the matrix dimensions are not equal,or no of nonzeros */ 179899cda47SBarry Smith if ((A->rmap.n != B->rmap.n) ||(a->nz != b->nz)) { 1800f5bd95cSBarry Smith flag = PETSC_FALSE; 181b97cf49bSBarry Smith } 182b97cf49bSBarry Smith 183b97cf49bSBarry Smith /* if the a->i are the same */ 184899cda47SBarry Smith ierr = PetscMemcmp(a->i,b->i,(A->rmap.n+1)*sizeof(PetscInt),&flag);CHKERRQ(ierr); 185b97cf49bSBarry Smith 186b97cf49bSBarry Smith /* if a->j are the same */ 187b24ad042SBarry Smith ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),&flag);CHKERRQ(ierr); 188b97cf49bSBarry Smith 189ca161407SBarry Smith ierr = MPI_Allreduce(&flag,flg,1,MPI_INT,MPI_LAND,A->comm);CHKERRQ(ierr); 1903a40ed3dSBarry Smith PetscFunctionReturn(0); 191b97cf49bSBarry Smith } 192b97cf49bSBarry Smith 1934a2ae208SSatish Balay #undef __FUNCT__ 1944a2ae208SSatish Balay #define __FUNCT__ "MatGetRowIJ_MPIAdj" 195b24ad042SBarry Smith PetscErrorCode MatGetRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscTruth symmetric,PetscInt *m,PetscInt *ia[],PetscInt *ja[],PetscTruth *done) 1966cd91f64SBarry Smith { 197dfbe8321SBarry Smith PetscErrorCode ierr; 198b24ad042SBarry Smith PetscMPIInt size; 199b24ad042SBarry Smith PetscInt i; 2006cd91f64SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data; 2016cd91f64SBarry Smith 2026cd91f64SBarry Smith PetscFunctionBegin; 2036cd91f64SBarry Smith ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); 2046cd91f64SBarry Smith if (size > 1) {*done = PETSC_FALSE; PetscFunctionReturn(0);} 205899cda47SBarry Smith *m = A->rmap.n; 2066cd91f64SBarry Smith *ia = a->i; 2076cd91f64SBarry Smith *ja = a->j; 2086cd91f64SBarry Smith *done = PETSC_TRUE; 209d42735a1SBarry Smith if (oshift) { 210d42735a1SBarry Smith for (i=0; i<(*ia)[*m]; i++) { 211d42735a1SBarry Smith (*ja)[i]++; 212d42735a1SBarry Smith } 213d42735a1SBarry Smith for (i=0; i<=(*m); i++) (*ia)[i]++; 214d42735a1SBarry Smith } 215d42735a1SBarry Smith PetscFunctionReturn(0); 216d42735a1SBarry Smith } 217d42735a1SBarry Smith 2184a2ae208SSatish Balay #undef __FUNCT__ 2194a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_MPIAdj" 220b24ad042SBarry Smith PetscErrorCode MatRestoreRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscTruth symmetric,PetscInt *m,PetscInt *ia[],PetscInt *ja[],PetscTruth *done) 221d42735a1SBarry Smith { 222b24ad042SBarry Smith PetscInt i; 223d42735a1SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data; 224d42735a1SBarry Smith 225d42735a1SBarry Smith PetscFunctionBegin; 226e005ede5SBarry Smith if (ia && a->i != *ia) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"ia passed back is not one obtained with MatGetRowIJ()"); 227e005ede5SBarry Smith if (ja && a->j != *ja) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"ja passed back is not one obtained with MatGetRowIJ()"); 228d42735a1SBarry Smith if (oshift) { 229d42735a1SBarry Smith for (i=0; i<=(*m); i++) (*ia)[i]--; 230d42735a1SBarry Smith for (i=0; i<(*ia)[*m]; i++) { 231d42735a1SBarry Smith (*ja)[i]--; 232d42735a1SBarry Smith } 233d42735a1SBarry Smith } 2346cd91f64SBarry Smith PetscFunctionReturn(0); 2356cd91f64SBarry Smith } 236b97cf49bSBarry Smith 237b97cf49bSBarry Smith /* -------------------------------------------------------------------*/ 2380b3b1487SBarry Smith static struct _MatOps MatOps_Values = {0, 2393eda8832SBarry Smith MatGetRow_MPIAdj, 2403eda8832SBarry Smith MatRestoreRow_MPIAdj, 241b97cf49bSBarry Smith 0, 24297304618SKris Buschelman /* 4*/ 0, 243b97cf49bSBarry Smith 0, 244b97cf49bSBarry Smith 0, 245b97cf49bSBarry Smith 0, 2460b3b1487SBarry Smith 0, 2470b3b1487SBarry Smith 0, 24897304618SKris Buschelman /*10*/ 0, 2490b3b1487SBarry Smith 0, 2500b3b1487SBarry Smith 0, 2510b3b1487SBarry Smith 0, 2520b3b1487SBarry Smith 0, 25397304618SKris Buschelman /*15*/ 0, 2543eda8832SBarry Smith MatEqual_MPIAdj, 2550b3b1487SBarry Smith 0, 2560b3b1487SBarry Smith 0, 2570b3b1487SBarry Smith 0, 25897304618SKris Buschelman /*20*/ 0, 2590b3b1487SBarry Smith 0, 2600b3b1487SBarry Smith 0, 2613eda8832SBarry Smith MatSetOption_MPIAdj, 2620b3b1487SBarry Smith 0, 26397304618SKris Buschelman /*25*/ 0, 2640b3b1487SBarry Smith 0, 2650b3b1487SBarry Smith 0, 2660b3b1487SBarry Smith 0, 2670b3b1487SBarry Smith 0, 26897304618SKris Buschelman /*30*/ 0, 2690b3b1487SBarry Smith 0, 270273d9f13SBarry Smith 0, 271273d9f13SBarry Smith 0, 2720b3b1487SBarry Smith 0, 27397304618SKris Buschelman /*35*/ 0, 2740b3b1487SBarry Smith 0, 2750b3b1487SBarry Smith 0, 2760b3b1487SBarry Smith 0, 2770b3b1487SBarry Smith 0, 27897304618SKris Buschelman /*40*/ 0, 2790b3b1487SBarry Smith 0, 2800b3b1487SBarry Smith 0, 2810b3b1487SBarry Smith 0, 2820b3b1487SBarry Smith 0, 28397304618SKris Buschelman /*45*/ 0, 2840b3b1487SBarry Smith 0, 2850b3b1487SBarry Smith 0, 2860b3b1487SBarry Smith 0, 2870b3b1487SBarry Smith 0, 288521d7252SBarry Smith /*50*/ 0, 2896cd91f64SBarry Smith MatGetRowIJ_MPIAdj, 290d42735a1SBarry Smith MatRestoreRowIJ_MPIAdj, 291b97cf49bSBarry Smith 0, 292b97cf49bSBarry Smith 0, 29397304618SKris Buschelman /*55*/ 0, 294b97cf49bSBarry Smith 0, 2950752156aSBarry Smith 0, 2960752156aSBarry Smith 0, 2970b3b1487SBarry Smith 0, 29897304618SKris Buschelman /*60*/ 0, 299b9b97703SBarry Smith MatDestroy_MPIAdj, 300b9b97703SBarry Smith MatView_MPIAdj, 301357abbc8SBarry Smith 0, 30297304618SKris Buschelman 0, 30397304618SKris Buschelman /*65*/ 0, 30497304618SKris Buschelman 0, 30597304618SKris Buschelman 0, 30697304618SKris Buschelman 0, 30797304618SKris Buschelman 0, 30897304618SKris Buschelman /*70*/ 0, 30997304618SKris Buschelman 0, 31097304618SKris Buschelman 0, 31197304618SKris Buschelman 0, 31297304618SKris Buschelman 0, 31397304618SKris Buschelman /*75*/ 0, 31497304618SKris Buschelman 0, 31597304618SKris Buschelman 0, 31697304618SKris Buschelman 0, 31797304618SKris Buschelman 0, 31897304618SKris Buschelman /*80*/ 0, 31997304618SKris Buschelman 0, 32097304618SKris Buschelman 0, 32197304618SKris Buschelman 0, 322865e5f61SKris Buschelman 0, 323865e5f61SKris Buschelman /*85*/ 0, 324865e5f61SKris Buschelman 0, 325865e5f61SKris Buschelman 0, 326865e5f61SKris Buschelman 0, 327865e5f61SKris Buschelman 0, 328865e5f61SKris Buschelman /*90*/ 0, 329865e5f61SKris Buschelman 0, 330865e5f61SKris Buschelman 0, 331865e5f61SKris Buschelman 0, 332865e5f61SKris Buschelman 0, 333865e5f61SKris Buschelman /*95*/ 0, 334865e5f61SKris Buschelman 0, 335865e5f61SKris Buschelman 0, 336865e5f61SKris Buschelman 0}; 337b97cf49bSBarry Smith 338a23d5eceSKris Buschelman EXTERN_C_BEGIN 339a23d5eceSKris Buschelman #undef __FUNCT__ 340a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIAdjSetPreallocation_MPIAdj" 341be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAdjSetPreallocation_MPIAdj(Mat B,PetscInt *i,PetscInt *j,PetscInt *values) 342a23d5eceSKris Buschelman { 343a23d5eceSKris Buschelman Mat_MPIAdj *b = (Mat_MPIAdj *)B->data; 344dfbe8321SBarry Smith PetscErrorCode ierr; 3452515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 346b24ad042SBarry Smith PetscInt ii; 347a23d5eceSKris Buschelman #endif 348a23d5eceSKris Buschelman 349a23d5eceSKris Buschelman PetscFunctionBegin; 350a23d5eceSKris Buschelman B->preallocated = PETSC_TRUE; 3512515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 35277431f27SBarry Smith if (i[0] != 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"First i[] index must be zero, instead it is %D\n",i[0]); 353899cda47SBarry Smith for (ii=1; ii<B->rmap.n; ii++) { 354a23d5eceSKris Buschelman if (i[ii] < 0 || i[ii] < i[ii-1]) { 35577431f27SBarry 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]); 356a23d5eceSKris Buschelman } 357a23d5eceSKris Buschelman } 358899cda47SBarry Smith for (ii=0; ii<i[B->rmap.n]; ii++) { 359899cda47SBarry Smith if (j[ii] < 0 || j[ii] >= B->cmap.N) { 36077431f27SBarry Smith SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column index %D out of range %D\n",ii,j[ii]); 361a23d5eceSKris Buschelman } 362a23d5eceSKris Buschelman } 363a23d5eceSKris Buschelman #endif 364a23d5eceSKris Buschelman 365a23d5eceSKris Buschelman b->j = j; 366a23d5eceSKris Buschelman b->i = i; 367a23d5eceSKris Buschelman b->values = values; 368a23d5eceSKris Buschelman 369899cda47SBarry Smith b->nz = i[B->rmap.n]; 370a23d5eceSKris Buschelman b->diag = 0; 371a23d5eceSKris Buschelman b->symmetric = PETSC_FALSE; 372a23d5eceSKris Buschelman b->freeaij = PETSC_TRUE; 373a23d5eceSKris Buschelman 374a23d5eceSKris Buschelman ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 375a23d5eceSKris Buschelman ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 376a23d5eceSKris Buschelman PetscFunctionReturn(0); 377a23d5eceSKris Buschelman } 378a23d5eceSKris Buschelman EXTERN_C_END 379b97cf49bSBarry Smith 3800bad9183SKris Buschelman /*MC 381fafad747SKris Buschelman MATMPIADJ - MATMPIADJ = "mpiadj" - A matrix type to be used for distributed adjacency matrices, 3820bad9183SKris Buschelman intended for use constructing orderings and partitionings. 3830bad9183SKris Buschelman 3840bad9183SKris Buschelman Level: beginner 3850bad9183SKris Buschelman 3860bad9183SKris Buschelman .seealso: MatCreateMPIAdj 3870bad9183SKris Buschelman M*/ 3880bad9183SKris Buschelman 389273d9f13SBarry Smith EXTERN_C_BEGIN 3904a2ae208SSatish Balay #undef __FUNCT__ 3914a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIAdj" 392be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_MPIAdj(Mat B) 393273d9f13SBarry Smith { 394273d9f13SBarry Smith Mat_MPIAdj *b; 3956849ba73SBarry Smith PetscErrorCode ierr; 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->mapping = 0; 407273d9f13SBarry Smith B->assembled = PETSC_FALSE; 408273d9f13SBarry Smith 409899cda47SBarry Smith ierr = PetscMapInitialize(B->comm,&B->rmap);CHKERRQ(ierr); 410e2ee2a47SBarry Smith if (B->cmap.n < 0) B->cmap.n = B->cmap.N; 411e2ee2a47SBarry Smith if (B->cmap.N < 0) B->cmap.N = B->cmap.n; 412273d9f13SBarry Smith 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 @*/ 438be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAdjSetPreallocation(Mat B,PetscInt *i,PetscInt *j,PetscInt *values) 439273d9f13SBarry Smith { 440b24ad042SBarry Smith PetscErrorCode ierr,(*f)(Mat,PetscInt*,PetscInt*,PetscInt*); 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 @*/ 487be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAdj(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt *i,PetscInt *j,PetscInt *values,Mat *A) 488b97cf49bSBarry Smith { 489dfbe8321SBarry Smith PetscErrorCode ierr; 490b97cf49bSBarry Smith 491433994e6SBarry Smith PetscFunctionBegin; 492f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 493e2ee2a47SBarry Smith ierr = MatSetSizes(*A,m,n,PETSC_DETERMINE,n);CHKERRQ(ierr); 494273d9f13SBarry Smith ierr = MatSetType(*A,MATMPIADJ);CHKERRQ(ierr); 495273d9f13SBarry Smith ierr = MatMPIAdjSetPreallocation(*A,i,j,values);CHKERRQ(ierr); 4963a40ed3dSBarry Smith PetscFunctionReturn(0); 497b97cf49bSBarry Smith } 498b97cf49bSBarry Smith 499273d9f13SBarry Smith EXTERN_C_BEGIN 5004a2ae208SSatish Balay #undef __FUNCT__ 5014a2ae208SSatish Balay #define __FUNCT__ "MatConvertTo_MPIAdj" 502be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatConvertTo_MPIAdj(Mat A,MatType type,MatReuse reuse,Mat *newmat) 503b3cc6726SBarry Smith { 504676c34cdSKris Buschelman Mat B; 5056849ba73SBarry Smith PetscErrorCode ierr; 506b24ad042SBarry Smith PetscInt i,m,N,nzeros = 0,*ia,*ja,len,rstart,cnt,j,*a; 507b24ad042SBarry Smith const PetscInt *rj; 508b3cc6726SBarry Smith const PetscScalar *ra; 509c2e958feSBarry Smith MPI_Comm comm; 510c2e958feSBarry Smith 511c2e958feSBarry Smith PetscFunctionBegin; 512c2e958feSBarry Smith ierr = MatGetSize(A,PETSC_NULL,&N);CHKERRQ(ierr); 513c2e958feSBarry Smith ierr = MatGetLocalSize(A,&m,PETSC_NULL);CHKERRQ(ierr); 514c2e958feSBarry Smith ierr = MatGetOwnershipRange(A,&rstart,PETSC_NULL);CHKERRQ(ierr); 515c2e958feSBarry Smith 516c2e958feSBarry Smith /* count the number of nonzeros per row */ 517c2e958feSBarry Smith for (i=0; i<m; i++) { 518c2e958feSBarry Smith ierr = MatGetRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr); 519c2e958feSBarry Smith for (j=0; j<len; j++) { 520c2e958feSBarry Smith if (rj[j] == i+rstart) {len--; break;} /* don't count diagonal */ 521c2e958feSBarry Smith } 522c2e958feSBarry Smith ierr = MatRestoreRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr); 523c2e958feSBarry Smith nzeros += len; 524c2e958feSBarry Smith } 525c2e958feSBarry Smith 526c2e958feSBarry Smith /* malloc space for nonzeros */ 527b24ad042SBarry Smith ierr = PetscMalloc((nzeros+1)*sizeof(PetscInt),&a);CHKERRQ(ierr); 528b24ad042SBarry Smith ierr = PetscMalloc((N+1)*sizeof(PetscInt),&ia);CHKERRQ(ierr); 529b24ad042SBarry Smith ierr = PetscMalloc((nzeros+1)*sizeof(PetscInt),&ja);CHKERRQ(ierr); 530c2e958feSBarry Smith 531c2e958feSBarry Smith nzeros = 0; 532c2e958feSBarry Smith ia[0] = 0; 533c2e958feSBarry Smith for (i=0; i<m; i++) { 534c2e958feSBarry Smith ierr = MatGetRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr); 535c2e958feSBarry Smith cnt = 0; 536c2e958feSBarry Smith for (j=0; j<len; j++) { 537c2e958feSBarry Smith if (rj[j] != i+rstart) { /* if not diagonal */ 538b24ad042SBarry Smith a[nzeros+cnt] = (PetscInt) PetscAbsScalar(ra[j]); 539c2e958feSBarry Smith ja[nzeros+cnt++] = rj[j]; 540c2e958feSBarry Smith } 541c2e958feSBarry Smith } 542c2e958feSBarry Smith ierr = MatRestoreRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr); 543c2e958feSBarry Smith nzeros += cnt; 544c2e958feSBarry Smith ia[i+1] = nzeros; 545c2e958feSBarry Smith } 546c2e958feSBarry Smith 547c2e958feSBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 548f69a0ea3SMatthew Knepley ierr = MatCreate(comm,&B);CHKERRQ(ierr); 549f69a0ea3SMatthew Knepley ierr = MatSetSizes(B,m,N,PETSC_DETERMINE,N);CHKERRQ(ierr); 550f204ca49SKris Buschelman ierr = MatSetType(B,type);CHKERRQ(ierr); 551f204ca49SKris Buschelman ierr = MatMPIAdjSetPreallocation(B,ia,ja,a);CHKERRQ(ierr); 552c2e958feSBarry Smith 553ceb03754SKris Buschelman if (reuse == MAT_REUSE_MATRIX) { 5548ab5b326SKris Buschelman ierr = MatHeaderReplace(A,B);CHKERRQ(ierr); 555c3d102feSKris Buschelman } else { 556676c34cdSKris Buschelman *newmat = B; 557c3d102feSKris Buschelman } 558c2e958feSBarry Smith PetscFunctionReturn(0); 559c2e958feSBarry Smith } 560273d9f13SBarry Smith EXTERN_C_END 561433994e6SBarry Smith 562433994e6SBarry Smith 563433994e6SBarry Smith 564433994e6SBarry Smith 565433994e6SBarry Smith 566