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) { 8712c028f9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 88b97cf49bSBarry Smith a->symmetric = PETSC_TRUE; 8912c028f9SKris Buschelman break; 9012c028f9SKris Buschelman default: 91b0a32e0cSBarry Smith PetscLogInfo(A,"MatSetOption_MPIAdj:Option ignored\n"); 9212c028f9SKris Buschelman break; 93b97cf49bSBarry Smith } 943a40ed3dSBarry Smith PetscFunctionReturn(0); 95b97cf49bSBarry Smith } 96b97cf49bSBarry Smith 97b97cf49bSBarry Smith 98b97cf49bSBarry Smith /* 99b97cf49bSBarry Smith Adds diagonal pointers to sparse matrix structure. 100b97cf49bSBarry Smith */ 101b97cf49bSBarry Smith 1024a2ae208SSatish Balay #undef __FUNCT__ 1034a2ae208SSatish Balay #define __FUNCT__ "MatMarkDiagonal_MPIAdj" 1043eda8832SBarry Smith int MatMarkDiagonal_MPIAdj(Mat A) 105b97cf49bSBarry Smith { 1063eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 10782502324SSatish Balay int i,j,*diag,m = A->m,ierr; 108b97cf49bSBarry Smith 109433994e6SBarry Smith PetscFunctionBegin; 110b0a32e0cSBarry Smith ierr = PetscMalloc((m+1)*sizeof(int),&diag);CHKERRQ(ierr); 111b0a32e0cSBarry Smith PetscLogObjectMemory(A,(m+1)*sizeof(int)); 112273d9f13SBarry Smith for (i=0; i<A->m; i++) { 113b97cf49bSBarry Smith for (j=a->i[i]; j<a->i[i+1]; j++) { 114b97cf49bSBarry Smith if (a->j[j] == i) { 115b97cf49bSBarry Smith diag[i] = j; 116b97cf49bSBarry Smith break; 117b97cf49bSBarry Smith } 118b97cf49bSBarry Smith } 119b97cf49bSBarry Smith } 120b97cf49bSBarry Smith a->diag = diag; 1213a40ed3dSBarry Smith PetscFunctionReturn(0); 122b97cf49bSBarry Smith } 123b97cf49bSBarry Smith 1244a2ae208SSatish Balay #undef __FUNCT__ 1254a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_MPIAdj" 12687828ca2SBarry Smith int MatGetRow_MPIAdj(Mat A,int row,int *nz,int **idx,PetscScalar **v) 127b97cf49bSBarry Smith { 1283eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 129b97cf49bSBarry Smith int *itmp; 130b97cf49bSBarry Smith 131433994e6SBarry Smith PetscFunctionBegin; 1320752156aSBarry Smith row -= a->rstart; 1330752156aSBarry Smith 134273d9f13SBarry Smith if (row < 0 || row >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Row out of range"); 135b97cf49bSBarry Smith 136b97cf49bSBarry Smith *nz = a->i[row+1] - a->i[row]; 137b97cf49bSBarry Smith if (v) *v = PETSC_NULL; 138b97cf49bSBarry Smith if (idx) { 139b97cf49bSBarry Smith itmp = a->j + a->i[row]; 140b97cf49bSBarry Smith if (*nz) { 141b97cf49bSBarry Smith *idx = itmp; 142b97cf49bSBarry Smith } 143b97cf49bSBarry Smith else *idx = 0; 144b97cf49bSBarry Smith } 1453a40ed3dSBarry Smith PetscFunctionReturn(0); 146b97cf49bSBarry Smith } 147b97cf49bSBarry Smith 1484a2ae208SSatish Balay #undef __FUNCT__ 1494a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_MPIAdj" 15087828ca2SBarry Smith int MatRestoreRow_MPIAdj(Mat A,int row,int *nz,int **idx,PetscScalar **v) 151b97cf49bSBarry Smith { 152433994e6SBarry Smith PetscFunctionBegin; 1533a40ed3dSBarry Smith PetscFunctionReturn(0); 154b97cf49bSBarry Smith } 155b97cf49bSBarry Smith 1564a2ae208SSatish Balay #undef __FUNCT__ 1574a2ae208SSatish Balay #define __FUNCT__ "MatGetBlockSize_MPIAdj" 1583eda8832SBarry Smith int MatGetBlockSize_MPIAdj(Mat A,int *bs) 159b97cf49bSBarry Smith { 160433994e6SBarry Smith PetscFunctionBegin; 161b97cf49bSBarry Smith *bs = 1; 1623a40ed3dSBarry Smith PetscFunctionReturn(0); 163b97cf49bSBarry Smith } 164b97cf49bSBarry Smith 165b97cf49bSBarry Smith 1664a2ae208SSatish Balay #undef __FUNCT__ 1674a2ae208SSatish Balay #define __FUNCT__ "MatEqual_MPIAdj" 1683eda8832SBarry Smith int MatEqual_MPIAdj(Mat A,Mat B,PetscTruth* flg) 169b97cf49bSBarry Smith { 1703eda8832SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data,*b = (Mat_MPIAdj *)B->data; 1710f5bd95cSBarry Smith int ierr; 1720f5bd95cSBarry Smith PetscTruth flag; 173b97cf49bSBarry Smith 174433994e6SBarry Smith PetscFunctionBegin; 175273d9f13SBarry Smith ierr = PetscTypeCompare((PetscObject)B,MATMPIADJ,&flag);CHKERRQ(ierr); 176273d9f13SBarry Smith if (!flag) SETERRQ(PETSC_ERR_ARG_INCOMP,"Matrices must be same type"); 177b97cf49bSBarry Smith 178b97cf49bSBarry Smith /* If the matrix dimensions are not equal,or no of nonzeros */ 179273d9f13SBarry Smith if ((A->m != B->m) ||(a->nz != b->nz)) { 1800f5bd95cSBarry Smith flag = PETSC_FALSE; 181b97cf49bSBarry Smith } 182b97cf49bSBarry Smith 183b97cf49bSBarry Smith /* if the a->i are the same */ 184273d9f13SBarry Smith ierr = PetscMemcmp(a->i,b->i,(A->m+1)*sizeof(int),&flag);CHKERRQ(ierr); 185b97cf49bSBarry Smith 186b97cf49bSBarry Smith /* if a->j are the same */ 1870f5bd95cSBarry Smith ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(int),&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" 1956cd91f64SBarry Smith int MatGetRowIJ_MPIAdj(Mat A,int oshift,PetscTruth symmetric,int *m,int **ia,int **ja,PetscTruth *done) 1966cd91f64SBarry Smith { 197d42735a1SBarry Smith int ierr,size,i; 1986cd91f64SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data; 1996cd91f64SBarry Smith 2006cd91f64SBarry Smith PetscFunctionBegin; 2016cd91f64SBarry Smith ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); 2026cd91f64SBarry Smith if (size > 1) {*done = PETSC_FALSE; PetscFunctionReturn(0);} 2036cd91f64SBarry Smith *m = A->m; 2046cd91f64SBarry Smith *ia = a->i; 2056cd91f64SBarry Smith *ja = a->j; 2066cd91f64SBarry Smith *done = PETSC_TRUE; 207d42735a1SBarry Smith if (oshift) { 208d42735a1SBarry Smith for (i=0; i<(*ia)[*m]; i++) { 209d42735a1SBarry Smith (*ja)[i]++; 210d42735a1SBarry Smith } 211d42735a1SBarry Smith for (i=0; i<=(*m); i++) (*ia)[i]++; 212d42735a1SBarry Smith } 213d42735a1SBarry Smith PetscFunctionReturn(0); 214d42735a1SBarry Smith } 215d42735a1SBarry Smith 2164a2ae208SSatish Balay #undef __FUNCT__ 2174a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRowIJ_MPIAdj" 218d42735a1SBarry Smith int MatRestoreRowIJ_MPIAdj(Mat A,int oshift,PetscTruth symmetric,int *m,int **ia,int **ja,PetscTruth *done) 219d42735a1SBarry Smith { 220d42735a1SBarry Smith int i; 221d42735a1SBarry Smith Mat_MPIAdj *a = (Mat_MPIAdj *)A->data; 222d42735a1SBarry Smith 223d42735a1SBarry Smith PetscFunctionBegin; 224c5b3d447SBarry Smith if (ia && a->i != *ia) SETERRQ(1,"ia passed back is not one obtained with MatGetRowIJ()"); 225c5b3d447SBarry Smith if (ja && a->j != *ja) SETERRQ(1,"ja passed back is not one obtained with MatGetRowIJ()"); 226d42735a1SBarry Smith if (oshift) { 227d42735a1SBarry Smith for (i=0; i<=(*m); i++) (*ia)[i]--; 228d42735a1SBarry Smith for (i=0; i<(*ia)[*m]; i++) { 229d42735a1SBarry Smith (*ja)[i]--; 230d42735a1SBarry Smith } 231d42735a1SBarry Smith } 2326cd91f64SBarry Smith PetscFunctionReturn(0); 2336cd91f64SBarry Smith } 234b97cf49bSBarry Smith 235b97cf49bSBarry Smith /* -------------------------------------------------------------------*/ 2360b3b1487SBarry Smith static struct _MatOps MatOps_Values = {0, 2373eda8832SBarry Smith MatGetRow_MPIAdj, 2383eda8832SBarry Smith MatRestoreRow_MPIAdj, 239b97cf49bSBarry Smith 0, 240b97cf49bSBarry Smith 0, 241b97cf49bSBarry Smith 0, 242b97cf49bSBarry Smith 0, 2430b3b1487SBarry Smith 0, 2440b3b1487SBarry Smith 0, 2450b3b1487SBarry Smith 0, 2460b3b1487SBarry Smith 0, 2470b3b1487SBarry Smith 0, 2480b3b1487SBarry Smith 0, 2490b3b1487SBarry Smith 0, 2500b3b1487SBarry Smith 0, 2510b3b1487SBarry Smith 0, 2523eda8832SBarry Smith MatEqual_MPIAdj, 2530b3b1487SBarry Smith 0, 2540b3b1487SBarry Smith 0, 2550b3b1487SBarry Smith 0, 2560b3b1487SBarry Smith 0, 2570b3b1487SBarry Smith 0, 2580b3b1487SBarry Smith 0, 2593eda8832SBarry Smith MatSetOption_MPIAdj, 2600b3b1487SBarry Smith 0, 2610b3b1487SBarry Smith 0, 2620b3b1487SBarry Smith 0, 2630b3b1487SBarry Smith 0, 2640b3b1487SBarry Smith 0, 2650b3b1487SBarry Smith 0, 266273d9f13SBarry Smith 0, 267273d9f13SBarry Smith 0, 2680b3b1487SBarry Smith 0, 2690b3b1487SBarry Smith 0, 2700b3b1487SBarry Smith 0, 2710b3b1487SBarry Smith 0, 2720b3b1487SBarry Smith 0, 2730b3b1487SBarry Smith 0, 2740b3b1487SBarry Smith 0, 2750b3b1487SBarry Smith 0, 2760b3b1487SBarry Smith 0, 2770b3b1487SBarry Smith 0, 2780b3b1487SBarry Smith 0, 2790b3b1487SBarry Smith 0, 2800b3b1487SBarry Smith 0, 2810b3b1487SBarry Smith 0, 2820b3b1487SBarry Smith 0, 2830b3b1487SBarry Smith 0, 2840b3b1487SBarry Smith 0, 2850b3b1487SBarry Smith 0, 2863eda8832SBarry Smith MatGetBlockSize_MPIAdj, 2876cd91f64SBarry Smith MatGetRowIJ_MPIAdj, 288d42735a1SBarry Smith MatRestoreRowIJ_MPIAdj, 289b97cf49bSBarry Smith 0, 290b97cf49bSBarry Smith 0, 291b97cf49bSBarry Smith 0, 2920752156aSBarry Smith 0, 2930752156aSBarry Smith 0, 2940b3b1487SBarry Smith 0, 2950b3b1487SBarry Smith 0, 2960b3b1487SBarry Smith 0, 297b9b97703SBarry Smith MatDestroy_MPIAdj, 298b9b97703SBarry Smith MatView_MPIAdj, 2998a124369SBarry Smith MatGetPetscMaps_Petsc}; 300b97cf49bSBarry Smith 301*a23d5eceSKris Buschelman EXTERN_C_BEGIN 302*a23d5eceSKris Buschelman #undef __FUNCT__ 303*a23d5eceSKris Buschelman #define __FUNCT__ "MatMPIAdjSetPreallocation_MPIAdj" 304*a23d5eceSKris Buschelman int MatMPIAdjSetPreallocation_MPIAdj(Mat B,int *i,int *j,int *values) 305*a23d5eceSKris Buschelman { 306*a23d5eceSKris Buschelman Mat_MPIAdj *b = (Mat_MPIAdj *)B->data; 307*a23d5eceSKris Buschelman int ierr; 308*a23d5eceSKris Buschelman #if defined(PETSC_USE_BOPT_g) 309*a23d5eceSKris Buschelman int ii; 310*a23d5eceSKris Buschelman #endif 311*a23d5eceSKris Buschelman 312*a23d5eceSKris Buschelman PetscFunctionBegin; 313*a23d5eceSKris Buschelman B->preallocated = PETSC_TRUE; 314*a23d5eceSKris Buschelman #if defined(PETSC_USE_BOPT_g) 315*a23d5eceSKris Buschelman if (i[0] != 0) SETERRQ1(1,"First i[] index must be zero, instead it is %d\n",i[0]); 316*a23d5eceSKris Buschelman for (ii=1; ii<B->m; ii++) { 317*a23d5eceSKris Buschelman if (i[ii] < 0 || i[ii] < i[ii-1]) { 318*a23d5eceSKris Buschelman SETERRQ4(1,"i[%d]=%d index is out of range: i[%d]=%d",ii,i[ii],ii-1,i[ii-1]); 319*a23d5eceSKris Buschelman } 320*a23d5eceSKris Buschelman } 321*a23d5eceSKris Buschelman for (ii=0; ii<i[B->m]; ii++) { 322*a23d5eceSKris Buschelman if (j[ii] < 0 || j[ii] >= B->N) { 323*a23d5eceSKris Buschelman SETERRQ2(1,"Column index %d out of range %d\n",ii,j[ii]); 324*a23d5eceSKris Buschelman } 325*a23d5eceSKris Buschelman } 326*a23d5eceSKris Buschelman #endif 327*a23d5eceSKris Buschelman 328*a23d5eceSKris Buschelman b->j = j; 329*a23d5eceSKris Buschelman b->i = i; 330*a23d5eceSKris Buschelman b->values = values; 331*a23d5eceSKris Buschelman 332*a23d5eceSKris Buschelman b->nz = i[B->m]; 333*a23d5eceSKris Buschelman b->diag = 0; 334*a23d5eceSKris Buschelman b->symmetric = PETSC_FALSE; 335*a23d5eceSKris Buschelman b->freeaij = PETSC_TRUE; 336*a23d5eceSKris Buschelman 337*a23d5eceSKris Buschelman ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 338*a23d5eceSKris Buschelman ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 339*a23d5eceSKris Buschelman PetscFunctionReturn(0); 340*a23d5eceSKris Buschelman } 341*a23d5eceSKris Buschelman EXTERN_C_END 342b97cf49bSBarry Smith 343273d9f13SBarry Smith EXTERN_C_BEGIN 3444a2ae208SSatish Balay #undef __FUNCT__ 3454a2ae208SSatish Balay #define __FUNCT__ "MatCreate_MPIAdj" 346273d9f13SBarry Smith int MatCreate_MPIAdj(Mat B) 347273d9f13SBarry Smith { 348273d9f13SBarry Smith Mat_MPIAdj *b; 349273d9f13SBarry Smith int ii,ierr,size,rank; 350273d9f13SBarry Smith 351273d9f13SBarry Smith PetscFunctionBegin; 352273d9f13SBarry Smith ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); 353273d9f13SBarry Smith ierr = MPI_Comm_rank(B->comm,&rank);CHKERRQ(ierr); 354273d9f13SBarry Smith 355b0a32e0cSBarry Smith ierr = PetscNew(Mat_MPIAdj,&b);CHKERRQ(ierr); 356b0a32e0cSBarry Smith B->data = (void*)b; 357273d9f13SBarry Smith ierr = PetscMemzero(b,sizeof(Mat_MPIAdj));CHKERRQ(ierr); 358273d9f13SBarry Smith ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 359273d9f13SBarry Smith B->factor = 0; 360273d9f13SBarry Smith B->lupivotthreshold = 1.0; 361273d9f13SBarry Smith B->mapping = 0; 362273d9f13SBarry Smith B->assembled = PETSC_FALSE; 363273d9f13SBarry Smith 364273d9f13SBarry Smith ierr = MPI_Allreduce(&B->m,&B->M,1,MPI_INT,MPI_SUM,B->comm);CHKERRQ(ierr); 365273d9f13SBarry Smith B->n = B->N; 366273d9f13SBarry Smith 367273d9f13SBarry Smith /* the information in the maps duplicates the information computed below, eventually 368273d9f13SBarry Smith we should remove the duplicate information that is not contained in the maps */ 3698a124369SBarry Smith ierr = PetscMapCreateMPI(B->comm,B->m,B->M,&B->rmap);CHKERRQ(ierr); 370273d9f13SBarry Smith /* we don't know the "local columns" so just use the row information :-(*/ 3718a124369SBarry Smith ierr = PetscMapCreateMPI(B->comm,B->m,B->M,&B->cmap);CHKERRQ(ierr); 372273d9f13SBarry Smith 373b0a32e0cSBarry Smith ierr = PetscMalloc((size+1)*sizeof(int),&b->rowners);CHKERRQ(ierr); 374b0a32e0cSBarry Smith PetscLogObjectMemory(B,(size+2)*sizeof(int)+sizeof(struct _p_Mat)+sizeof(Mat_MPIAdj)); 375273d9f13SBarry Smith ierr = MPI_Allgather(&B->m,1,MPI_INT,b->rowners+1,1,MPI_INT,B->comm);CHKERRQ(ierr); 376273d9f13SBarry Smith b->rowners[0] = 0; 377273d9f13SBarry Smith for (ii=2; ii<=size; ii++) { 378273d9f13SBarry Smith b->rowners[ii] += b->rowners[ii-1]; 379273d9f13SBarry Smith } 380273d9f13SBarry Smith b->rstart = b->rowners[rank]; 381273d9f13SBarry Smith b->rend = b->rowners[rank+1]; 382*a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMPIAdjSetPreallocation_C", 383*a23d5eceSKris Buschelman "MatMPIAdjSetPreallocation_MPIAdj", 384*a23d5eceSKris Buschelman MatMPIAdjSetPreallocation_MPIAdj);CHKERRQ(ierr); 385273d9f13SBarry Smith PetscFunctionReturn(0); 386273d9f13SBarry Smith } 387273d9f13SBarry Smith EXTERN_C_END 388273d9f13SBarry Smith 3894a2ae208SSatish Balay #undef __FUNCT__ 3904a2ae208SSatish Balay #define __FUNCT__ "MatMPIAdjSetPreallocation" 391*a23d5eceSKris Buschelman /*@C 392*a23d5eceSKris Buschelman MatMPIAdjSetPreallocation - Sets the array used for storing the matrix elements 393*a23d5eceSKris Buschelman 394*a23d5eceSKris Buschelman Collective on MPI_Comm 395*a23d5eceSKris Buschelman 396*a23d5eceSKris Buschelman Input Parameters: 397*a23d5eceSKris Buschelman + A - the matrix 398*a23d5eceSKris Buschelman . i - the indices into j for the start of each row 399*a23d5eceSKris Buschelman . j - the column indices for each row (sorted for each row). 400*a23d5eceSKris Buschelman The indices in i and j start with zero (NOT with one). 401*a23d5eceSKris Buschelman - values - [optional] edge weights 402*a23d5eceSKris Buschelman 403*a23d5eceSKris Buschelman Level: intermediate 404*a23d5eceSKris Buschelman 405*a23d5eceSKris Buschelman .seealso: MatCreate(), MatCreateMPIAdj(), MatSetValues() 406*a23d5eceSKris Buschelman @*/ 407273d9f13SBarry Smith int MatMPIAdjSetPreallocation(Mat B,int *i,int *j,int *values) 408273d9f13SBarry Smith { 409*a23d5eceSKris Buschelman int ierr,(*f)(Mat,int*,int*,int*); 410273d9f13SBarry Smith 411273d9f13SBarry Smith PetscFunctionBegin; 412*a23d5eceSKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)B,"MatMPIAdjSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 413*a23d5eceSKris Buschelman if (f) { 414*a23d5eceSKris Buschelman ierr = (*f)(B,i,j,values);CHKERRQ(ierr); 415273d9f13SBarry Smith } 416273d9f13SBarry Smith PetscFunctionReturn(0); 417273d9f13SBarry Smith } 418273d9f13SBarry Smith 4194a2ae208SSatish Balay #undef __FUNCT__ 4204a2ae208SSatish Balay #define __FUNCT__ "MatCreateMPIAdj" 421b97cf49bSBarry Smith /*@C 4223eda8832SBarry Smith MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list. 423b97cf49bSBarry Smith The matrix does not have numerical values associated with it, but is 424b97cf49bSBarry Smith intended for ordering (to reduce bandwidth etc) and partitioning. 425b97cf49bSBarry Smith 426ef5ee4d1SLois Curfman McInnes Collective on MPI_Comm 427ef5ee4d1SLois Curfman McInnes 428b97cf49bSBarry Smith Input Parameters: 429c2e958feSBarry Smith + comm - MPI communicator 4300752156aSBarry Smith . m - number of local rows 431b97cf49bSBarry Smith . n - number of columns 432b97cf49bSBarry Smith . i - the indices into j for the start of each row 433329f5518SBarry Smith . j - the column indices for each row (sorted for each row). 434ef5ee4d1SLois Curfman McInnes The indices in i and j start with zero (NOT with one). 435329f5518SBarry Smith - values -[optional] edge weights 436b97cf49bSBarry Smith 437b97cf49bSBarry Smith Output Parameter: 438b97cf49bSBarry Smith . A - the matrix 439b97cf49bSBarry Smith 44015091d37SBarry Smith Level: intermediate 44115091d37SBarry Smith 4424bc6d8c0SBarry Smith Notes: This matrix object does not support most matrix operations, include 4434bc6d8c0SBarry Smith MatSetValues(). 444c2e958feSBarry Smith You must NOT free the ii, values and jj arrays yourself. PETSc will free them 445ededeb1aSvictorle when the matrix is destroyed; you must allocate them with PetscMalloc(). If you 4461198db28SBarry Smith call from Fortran you need not create the arrays with PetscMalloc(). 447327e1a73SBarry Smith Should not include the matrix diagonals. 448b97cf49bSBarry Smith 449e3f7e1d6Svictorle If you already have a matrix, you can create its adjacency matrix by a call 450ededeb1aSvictorle to MatConvert, specifying a type of MATMPIADJ. 451ededeb1aSvictorle 452ef5ee4d1SLois Curfman McInnes Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC 453b97cf49bSBarry Smith 454e3f7e1d6Svictorle .seealso: MatCreate(), MatConvert(), MatGetOrdering() 455b97cf49bSBarry Smith @*/ 4563eda8832SBarry Smith int MatCreateMPIAdj(MPI_Comm comm,int m,int n,int *i,int *j,int *values,Mat *A) 457b97cf49bSBarry Smith { 458273d9f13SBarry Smith int ierr; 459b97cf49bSBarry Smith 460433994e6SBarry Smith PetscFunctionBegin; 461273d9f13SBarry Smith ierr = MatCreate(comm,m,n,PETSC_DETERMINE,n,A);CHKERRQ(ierr); 462273d9f13SBarry Smith ierr = MatSetType(*A,MATMPIADJ);CHKERRQ(ierr); 463273d9f13SBarry Smith ierr = MatMPIAdjSetPreallocation(*A,i,j,values);CHKERRQ(ierr); 4643a40ed3dSBarry Smith PetscFunctionReturn(0); 465b97cf49bSBarry Smith } 466b97cf49bSBarry Smith 467273d9f13SBarry Smith EXTERN_C_BEGIN 4684a2ae208SSatish Balay #undef __FUNCT__ 4694a2ae208SSatish Balay #define __FUNCT__ "MatConvertTo_MPIAdj" 470273d9f13SBarry Smith int MatConvertTo_MPIAdj(Mat A,MatType type,Mat *B) 471c2e958feSBarry Smith { 4724c49b128SBarry Smith int i,ierr,m,N,nzeros = 0,*ia,*ja,*rj,len,rstart,cnt,j,*a; 47387828ca2SBarry Smith PetscScalar *ra; 474c2e958feSBarry Smith MPI_Comm comm; 475c2e958feSBarry Smith 476c2e958feSBarry Smith PetscFunctionBegin; 477c2e958feSBarry Smith ierr = MatGetSize(A,PETSC_NULL,&N);CHKERRQ(ierr); 478c2e958feSBarry Smith ierr = MatGetLocalSize(A,&m,PETSC_NULL);CHKERRQ(ierr); 479c2e958feSBarry Smith ierr = MatGetOwnershipRange(A,&rstart,PETSC_NULL);CHKERRQ(ierr); 480c2e958feSBarry Smith 481c2e958feSBarry Smith /* count the number of nonzeros per row */ 482c2e958feSBarry Smith for (i=0; i<m; i++) { 483c2e958feSBarry Smith ierr = MatGetRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr); 484c2e958feSBarry Smith for (j=0; j<len; j++) { 485c2e958feSBarry Smith if (rj[j] == i+rstart) {len--; break;} /* don't count diagonal */ 486c2e958feSBarry Smith } 487c2e958feSBarry Smith ierr = MatRestoreRow(A,i+rstart,&len,&rj,PETSC_NULL);CHKERRQ(ierr); 488c2e958feSBarry Smith nzeros += len; 489c2e958feSBarry Smith } 490c2e958feSBarry Smith 491c2e958feSBarry Smith /* malloc space for nonzeros */ 492b0a32e0cSBarry Smith ierr = PetscMalloc((nzeros+1)*sizeof(int),&a);CHKERRQ(ierr); 493b0a32e0cSBarry Smith ierr = PetscMalloc((N+1)*sizeof(int),&ia);CHKERRQ(ierr); 494b0a32e0cSBarry Smith ierr = PetscMalloc((nzeros+1)*sizeof(int),&ja);CHKERRQ(ierr); 495c2e958feSBarry Smith 496c2e958feSBarry Smith nzeros = 0; 497c2e958feSBarry Smith ia[0] = 0; 498c2e958feSBarry Smith for (i=0; i<m; i++) { 499c2e958feSBarry Smith ierr = MatGetRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr); 500c2e958feSBarry Smith cnt = 0; 501c2e958feSBarry Smith for (j=0; j<len; j++) { 502c2e958feSBarry Smith if (rj[j] != i+rstart) { /* if not diagonal */ 5034c49b128SBarry Smith a[nzeros+cnt] = (int) PetscAbsScalar(ra[j]); 504c2e958feSBarry Smith ja[nzeros+cnt++] = rj[j]; 505c2e958feSBarry Smith } 506c2e958feSBarry Smith } 507c2e958feSBarry Smith ierr = MatRestoreRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr); 508c2e958feSBarry Smith nzeros += cnt; 509c2e958feSBarry Smith ia[i+1] = nzeros; 510c2e958feSBarry Smith } 511c2e958feSBarry Smith 512c2e958feSBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 5133eda8832SBarry Smith ierr = MatCreateMPIAdj(comm,m,N,ia,ja,a,B);CHKERRQ(ierr); 514c2e958feSBarry Smith 515c2e958feSBarry Smith PetscFunctionReturn(0); 516c2e958feSBarry Smith } 517273d9f13SBarry Smith EXTERN_C_END 518433994e6SBarry Smith 519433994e6SBarry Smith 520433994e6SBarry Smith 521433994e6SBarry Smith 522433994e6SBarry Smith 523