1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*6831982aSBarry Smith static char vcid[] = "$Id: dense.c,v 1.175 1999/10/13 20:37:16 bsmith Exp bsmith $"; 3cb512458SBarry Smith #endif 467e560aaSBarry Smith /* 567e560aaSBarry Smith Defines the basic matrix operations for sequential dense. 667e560aaSBarry Smith */ 7289bc588SBarry Smith 870f55243SBarry Smith #include "src/mat/impls/dense/seq/dense.h" 9eadb2fb4SBarry Smith #include "pinclude/blaslapack.h" 10289bc588SBarry Smith 115615d1e5SSatish Balay #undef __FUNC__ 125615d1e5SSatish Balay #define __FUNC__ "MatAXPY_SeqDense" 131987afe7SBarry Smith int MatAXPY_SeqDense(Scalar *alpha,Mat X,Mat Y) 141987afe7SBarry Smith { 151987afe7SBarry Smith Mat_SeqDense *x = (Mat_SeqDense*) X->data,*y = (Mat_SeqDense*) Y->data; 161987afe7SBarry Smith int N = x->m*x->n, one = 1; 173a40ed3dSBarry Smith 183a40ed3dSBarry Smith PetscFunctionBegin; 191987afe7SBarry Smith BLaxpy_( &N, alpha, x->v, &one, y->v, &one ); 201987afe7SBarry Smith PLogFlops(2*N-1); 213a40ed3dSBarry Smith PetscFunctionReturn(0); 221987afe7SBarry Smith } 231987afe7SBarry Smith 245615d1e5SSatish Balay #undef __FUNC__ 255615d1e5SSatish Balay #define __FUNC__ "MatGetInfo_SeqDense" 268f6be9afSLois Curfman McInnes int MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info) 27289bc588SBarry Smith { 284e220ebcSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense *) A->data; 294e220ebcSLois Curfman McInnes int i,N = a->m*a->n,count = 0; 304e220ebcSLois Curfman McInnes Scalar *v = a->v; 313a40ed3dSBarry Smith 323a40ed3dSBarry Smith PetscFunctionBegin; 33289bc588SBarry Smith for ( i=0; i<N; i++ ) {if (*v != 0.0) count++; v++;} 344e220ebcSLois Curfman McInnes 354e220ebcSLois Curfman McInnes info->rows_global = (double)a->m; 364e220ebcSLois Curfman McInnes info->columns_global = (double)a->n; 374e220ebcSLois Curfman McInnes info->rows_local = (double)a->m; 384e220ebcSLois Curfman McInnes info->columns_local = (double)a->n; 394e220ebcSLois Curfman McInnes info->block_size = 1.0; 404e220ebcSLois Curfman McInnes info->nz_allocated = (double)N; 414e220ebcSLois Curfman McInnes info->nz_used = (double)count; 424e220ebcSLois Curfman McInnes info->nz_unneeded = (double)(N-count); 434e220ebcSLois Curfman McInnes info->assemblies = (double)A->num_ass; 444e220ebcSLois Curfman McInnes info->mallocs = 0; 454e220ebcSLois Curfman McInnes info->memory = A->mem; 464e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; 474e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 484e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 494e220ebcSLois Curfman McInnes 503a40ed3dSBarry Smith PetscFunctionReturn(0); 51289bc588SBarry Smith } 52289bc588SBarry Smith 535615d1e5SSatish Balay #undef __FUNC__ 545615d1e5SSatish Balay #define __FUNC__ "MatScale_SeqDense" 558f6be9afSLois Curfman McInnes int MatScale_SeqDense(Scalar *alpha,Mat inA) 5680cd9d93SLois Curfman McInnes { 5780cd9d93SLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense *) inA->data; 5880cd9d93SLois Curfman McInnes int one = 1, nz; 5980cd9d93SLois Curfman McInnes 603a40ed3dSBarry Smith PetscFunctionBegin; 6180cd9d93SLois Curfman McInnes nz = a->m*a->n; 6280cd9d93SLois Curfman McInnes BLscal_( &nz, alpha, a->v, &one ); 6380cd9d93SLois Curfman McInnes PLogFlops(nz); 643a40ed3dSBarry Smith PetscFunctionReturn(0); 6580cd9d93SLois Curfman McInnes } 6680cd9d93SLois Curfman McInnes 67289bc588SBarry Smith /* ---------------------------------------------------------------*/ 68*6831982aSBarry Smith /* COMMENT: I have chosen to hide row permutation in the pivots, 69*6831982aSBarry Smith rather than put it in the Mat->row slot.*/ 705615d1e5SSatish Balay #undef __FUNC__ 715615d1e5SSatish Balay #define __FUNC__ "MatLUFactor_SeqDense" 728f6be9afSLois Curfman McInnes int MatLUFactor_SeqDense(Mat A,IS row,IS col,double f) 73289bc588SBarry Smith { 74c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 756abc6512SBarry Smith int info; 7655659b69SBarry Smith 773a40ed3dSBarry Smith PetscFunctionBegin; 78289bc588SBarry Smith if (!mat->pivots) { 7945d48df9SBarry Smith mat->pivots = (int *) PetscMalloc((mat->m+1)*sizeof(int));CHKPTRQ(mat->pivots); 80c0bbcb79SLois Curfman McInnes PLogObjectMemory(A,mat->m*sizeof(int)); 81289bc588SBarry Smith } 827a97a34bSBarry Smith if (!mat->m || !mat->n) PetscFunctionReturn(0); 83289bc588SBarry Smith LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info); 84a8c6a408SBarry Smith if (info<0) SETERRQ(PETSC_ERR_LIB,0,"Bad argument to LU factorization"); 85e3372554SBarry Smith if (info>0) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,0,"Bad LU factorization"); 86c0bbcb79SLois Curfman McInnes A->factor = FACTOR_LU; 8755659b69SBarry Smith PLogFlops((2*mat->n*mat->n*mat->n)/3); 883a40ed3dSBarry Smith PetscFunctionReturn(0); 89289bc588SBarry Smith } 906ee01492SSatish Balay 915615d1e5SSatish Balay #undef __FUNC__ 925609ef8eSBarry Smith #define __FUNC__ "MatDuplicate_SeqDense" 935609ef8eSBarry Smith int MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat) 9402cad45dSBarry Smith { 9502cad45dSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense *) A->data, *l; 9602cad45dSBarry Smith int ierr; 9702cad45dSBarry Smith Mat newi; 9802cad45dSBarry Smith 993a40ed3dSBarry Smith PetscFunctionBegin; 10002cad45dSBarry Smith ierr = MatCreateSeqDense(A->comm,mat->m,mat->n,PETSC_NULL,&newi);CHKERRQ(ierr); 10102cad45dSBarry Smith l = (Mat_SeqDense *) newi->data; 1025609ef8eSBarry Smith if (cpvalues == MAT_COPY_VALUES) { 103549d3d68SSatish Balay ierr = PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar));CHKERRQ(ierr); 10402cad45dSBarry Smith } 10502cad45dSBarry Smith newi->assembled = PETSC_TRUE; 10602cad45dSBarry Smith *newmat = newi; 1073a40ed3dSBarry Smith PetscFunctionReturn(0); 10802cad45dSBarry Smith } 10902cad45dSBarry Smith 1105615d1e5SSatish Balay #undef __FUNC__ 1115615d1e5SSatish Balay #define __FUNC__ "MatLUFactorSymbolic_SeqDense" 1128f6be9afSLois Curfman McInnes int MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,double f,Mat *fact) 113289bc588SBarry Smith { 1143a40ed3dSBarry Smith int ierr; 1153a40ed3dSBarry Smith 1163a40ed3dSBarry Smith PetscFunctionBegin; 1175609ef8eSBarry Smith ierr = MatDuplicate_SeqDense(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr); 1183a40ed3dSBarry Smith PetscFunctionReturn(0); 119289bc588SBarry Smith } 1206ee01492SSatish Balay 1215615d1e5SSatish Balay #undef __FUNC__ 1225615d1e5SSatish Balay #define __FUNC__ "MatLUFactorNumeric_SeqDense" 1238f6be9afSLois Curfman McInnes int MatLUFactorNumeric_SeqDense(Mat A,Mat *fact) 124289bc588SBarry Smith { 12502cad45dSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense*) A->data, *l = (Mat_SeqDense*) (*fact)->data; 1263a40ed3dSBarry Smith int ierr; 1273a40ed3dSBarry Smith 1283a40ed3dSBarry Smith PetscFunctionBegin; 12902cad45dSBarry Smith /* copy the numerical values */ 130549d3d68SSatish Balay ierr = PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar));CHKERRQ(ierr); 13102cad45dSBarry Smith (*fact)->factor = 0; 1323a40ed3dSBarry Smith ierr = MatLUFactor(*fact,0,0,1.0);CHKERRQ(ierr); 1333a40ed3dSBarry Smith PetscFunctionReturn(0); 134289bc588SBarry Smith } 1356ee01492SSatish Balay 1365615d1e5SSatish Balay #undef __FUNC__ 1375615d1e5SSatish Balay #define __FUNC__ "MatCholeskyFactorSymbolic_SeqDense" 1388f6be9afSLois Curfman McInnes int MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,double f,Mat *fact) 139289bc588SBarry Smith { 1403a40ed3dSBarry Smith int ierr; 1413a40ed3dSBarry Smith 1423a40ed3dSBarry Smith PetscFunctionBegin; 1433a40ed3dSBarry Smith ierr = MatConvert(A,MATSAME,fact);CHKERRQ(ierr); 1443a40ed3dSBarry Smith PetscFunctionReturn(0); 145289bc588SBarry Smith } 1466ee01492SSatish Balay 1475615d1e5SSatish Balay #undef __FUNC__ 1485615d1e5SSatish Balay #define __FUNC__ "MatCholeskyFactorNumeric_SeqDense" 1498f6be9afSLois Curfman McInnes int MatCholeskyFactorNumeric_SeqDense(Mat A,Mat *fact) 150289bc588SBarry Smith { 1513a40ed3dSBarry Smith int ierr; 1523a40ed3dSBarry Smith 1533a40ed3dSBarry Smith PetscFunctionBegin; 1543a40ed3dSBarry Smith ierr = MatCholeskyFactor(*fact,0,1.0);CHKERRQ(ierr); 1553a40ed3dSBarry Smith PetscFunctionReturn(0); 156289bc588SBarry Smith } 1576ee01492SSatish Balay 1585615d1e5SSatish Balay #undef __FUNC__ 1595615d1e5SSatish Balay #define __FUNC__ "MatCholeskyFactor_SeqDense" 1608f6be9afSLois Curfman McInnes int MatCholeskyFactor_SeqDense(Mat A,IS perm,double f) 161289bc588SBarry Smith { 162c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 163606d414cSSatish Balay int info,ierr; 16455659b69SBarry Smith 1653a40ed3dSBarry Smith PetscFunctionBegin; 166752f5567SLois Curfman McInnes if (mat->pivots) { 167606d414cSSatish Balay ierr = PetscFree(mat->pivots);CHKERRQ(ierr); 168c0bbcb79SLois Curfman McInnes PLogObjectMemory(A,-mat->m*sizeof(int)); 169752f5567SLois Curfman McInnes mat->pivots = 0; 170752f5567SLois Curfman McInnes } 1717a97a34bSBarry Smith if (!mat->m || !mat->n) PetscFunctionReturn(0); 172289bc588SBarry Smith LApotrf_("L",&mat->n,mat->v,&mat->m,&info); 173e3372554SBarry Smith if (info) SETERRQ(PETSC_ERR_MAT_CH_ZRPVT,0,"Bad factorization"); 174c0bbcb79SLois Curfman McInnes A->factor = FACTOR_CHOLESKY; 17555659b69SBarry Smith PLogFlops((mat->n*mat->n*mat->n)/3); 1763a40ed3dSBarry Smith PetscFunctionReturn(0); 177289bc588SBarry Smith } 178289bc588SBarry Smith 1795615d1e5SSatish Balay #undef __FUNC__ 1805615d1e5SSatish Balay #define __FUNC__ "MatSolve_SeqDense" 1818f6be9afSLois Curfman McInnes int MatSolve_SeqDense(Mat A,Vec xx,Vec yy) 182289bc588SBarry Smith { 183c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 184a57ad546SLois Curfman McInnes int one = 1, info, ierr; 1856abc6512SBarry Smith Scalar *x, *y; 18667e560aaSBarry Smith 1873a40ed3dSBarry Smith PetscFunctionBegin; 1887a97a34bSBarry Smith if (!mat->m || !mat->n) PetscFunctionReturn(0); 1897a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 1907a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 191549d3d68SSatish Balay ierr = PetscMemcpy(y,x,mat->m*sizeof(Scalar));CHKERRQ(ierr); 192c0bbcb79SLois Curfman McInnes if (A->factor == FACTOR_LU) { 19348d91487SBarry Smith LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 1947a97a34bSBarry Smith } else if (A->factor == FACTOR_CHOLESKY){ 19548d91487SBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 196289bc588SBarry Smith } 197a8c6a408SBarry Smith else SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Matrix must be factored to solve"); 198a8c6a408SBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,0,"MBad solve"); 1997a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2007a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 20155659b69SBarry Smith PLogFlops(mat->n*mat->n - mat->n); 2023a40ed3dSBarry Smith PetscFunctionReturn(0); 203289bc588SBarry Smith } 2046ee01492SSatish Balay 2055615d1e5SSatish Balay #undef __FUNC__ 2065615d1e5SSatish Balay #define __FUNC__ "MatSolveTrans_SeqDense" 2078f6be9afSLois Curfman McInnes int MatSolveTrans_SeqDense(Mat A,Vec xx,Vec yy) 208da3a660dSBarry Smith { 209c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 2107a97a34bSBarry Smith int ierr,one = 1, info; 2116abc6512SBarry Smith Scalar *x, *y; 21267e560aaSBarry Smith 2133a40ed3dSBarry Smith PetscFunctionBegin; 2147a97a34bSBarry Smith if (!mat->m || !mat->n) PetscFunctionReturn(0); 2157a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2167a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 217549d3d68SSatish Balay ierr = PetscMemcpy(y,x,mat->m*sizeof(Scalar));CHKERRQ(ierr); 218752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 219da3a660dSBarry Smith if (mat->pivots) { 22048d91487SBarry Smith LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 2217a97a34bSBarry Smith } else { 22248d91487SBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 223da3a660dSBarry Smith } 224a8c6a408SBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,0,"Bad solve"); 2257a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2267a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 22755659b69SBarry Smith PLogFlops(mat->n*mat->n - mat->n); 2283a40ed3dSBarry Smith PetscFunctionReturn(0); 229da3a660dSBarry Smith } 2306ee01492SSatish Balay 2315615d1e5SSatish Balay #undef __FUNC__ 2325615d1e5SSatish Balay #define __FUNC__ "MatSolveAdd_SeqDense" 2338f6be9afSLois Curfman McInnes int MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 234da3a660dSBarry Smith { 235c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 2367a97a34bSBarry Smith int ierr,one = 1, info; 2376abc6512SBarry Smith Scalar *x, *y, sone = 1.0; 238da3a660dSBarry Smith Vec tmp = 0; 23967e560aaSBarry Smith 2403a40ed3dSBarry Smith PetscFunctionBegin; 2417a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2427a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 2437a97a34bSBarry Smith if (!mat->m || !mat->n) PetscFunctionReturn(0); 244da3a660dSBarry Smith if (yy == zz) { 24578b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr); 246c0bbcb79SLois Curfman McInnes PLogObjectParent(A,tmp); 24778b31e54SBarry Smith ierr = VecCopy(yy,tmp);CHKERRQ(ierr); 248da3a660dSBarry Smith } 249549d3d68SSatish Balay ierr = PetscMemcpy(y,x,mat->m*sizeof(Scalar));CHKERRQ(ierr); 250752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 251da3a660dSBarry Smith if (mat->pivots) { 25248d91487SBarry Smith LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 253a8c6a408SBarry Smith } else { 25448d91487SBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 255da3a660dSBarry Smith } 256a8c6a408SBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,0,"Bad solve"); 257a8c6a408SBarry Smith if (tmp) {ierr = VecAXPY(&sone,tmp,yy);CHKERRQ(ierr); ierr = VecDestroy(tmp);CHKERRQ(ierr);} 258a8c6a408SBarry Smith else {ierr = VecAXPY(&sone,zz,yy);CHKERRQ(ierr);} 2597a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2607a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 26155659b69SBarry Smith PLogFlops(mat->n*mat->n - mat->n); 2623a40ed3dSBarry Smith PetscFunctionReturn(0); 263da3a660dSBarry Smith } 26467e560aaSBarry Smith 2655615d1e5SSatish Balay #undef __FUNC__ 2665615d1e5SSatish Balay #define __FUNC__ "MatSolveTransAdd_SeqDense" 2678f6be9afSLois Curfman McInnes int MatSolveTransAdd_SeqDense(Mat A,Vec xx,Vec zz, Vec yy) 268da3a660dSBarry Smith { 269c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 2706abc6512SBarry Smith int one = 1, info,ierr; 2716abc6512SBarry Smith Scalar *x, *y, sone = 1.0; 272da3a660dSBarry Smith Vec tmp; 27367e560aaSBarry Smith 2743a40ed3dSBarry Smith PetscFunctionBegin; 2757a97a34bSBarry Smith if (!mat->m || !mat->n) PetscFunctionReturn(0); 2767a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2777a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 278da3a660dSBarry Smith if (yy == zz) { 27978b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr); 280c0bbcb79SLois Curfman McInnes PLogObjectParent(A,tmp); 28178b31e54SBarry Smith ierr = VecCopy(yy,tmp);CHKERRQ(ierr); 282da3a660dSBarry Smith } 283549d3d68SSatish Balay ierr = PetscMemcpy(y,x,mat->m*sizeof(Scalar));CHKERRQ(ierr); 284752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 285da3a660dSBarry Smith if (mat->pivots) { 28648d91487SBarry Smith LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 2873a40ed3dSBarry Smith } else { 28848d91487SBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 289da3a660dSBarry Smith } 290a8c6a408SBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,0,"Bad solve"); 29190f02eecSBarry Smith if (tmp) { 29290f02eecSBarry Smith ierr = VecAXPY(&sone,tmp,yy);CHKERRQ(ierr); 29390f02eecSBarry Smith ierr = VecDestroy(tmp);CHKERRQ(ierr); 2943a40ed3dSBarry Smith } else { 29590f02eecSBarry Smith ierr = VecAXPY(&sone,zz,yy);CHKERRQ(ierr); 29690f02eecSBarry Smith } 2977a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2987a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 29955659b69SBarry Smith PLogFlops(mat->n*mat->n - mat->n); 3003a40ed3dSBarry Smith PetscFunctionReturn(0); 301da3a660dSBarry Smith } 302289bc588SBarry Smith /* ------------------------------------------------------------------*/ 3035615d1e5SSatish Balay #undef __FUNC__ 3045615d1e5SSatish Balay #define __FUNC__ "MatRelax_SeqDense" 3058f6be9afSLois Curfman McInnes int MatRelax_SeqDense(Mat A,Vec bb,double omega,MatSORType flag, 30620e1a0d4SLois Curfman McInnes double shift,int its,Vec xx) 307289bc588SBarry Smith { 308c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 309289bc588SBarry Smith Scalar *x, *b, *v = mat->v, zero = 0.0, xt; 310bc1b551cSSatish Balay int ierr, m = mat->m, i; 311aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX) 312bc1b551cSSatish Balay int o = 1; 313bc1b551cSSatish Balay #endif 314289bc588SBarry Smith 3153a40ed3dSBarry Smith PetscFunctionBegin; 316289bc588SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 3173a40ed3dSBarry Smith /* this is a hack fix, should have another version without the second BLdot */ 318bbb6d6a8SBarry Smith ierr = VecSet(&zero,xx);CHKERRQ(ierr); 319289bc588SBarry Smith } 3207a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 3217a97a34bSBarry Smith ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 322289bc588SBarry Smith while (its--) { 323289bc588SBarry Smith if (flag & SOR_FORWARD_SWEEP){ 324289bc588SBarry Smith for ( i=0; i<m; i++ ) { 325aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 326f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 327f1747703SBarry Smith not happy about returning a double complex */ 328f1747703SBarry Smith int _i; 329f1747703SBarry Smith Scalar sum = b[i]; 330f1747703SBarry Smith for ( _i=0; _i<m; _i++ ) { 3313f6de6efSSatish Balay sum -= PetscConj(v[i+_i*m])*x[_i]; 332f1747703SBarry Smith } 333f1747703SBarry Smith xt = sum; 334f1747703SBarry Smith #else 335289bc588SBarry Smith xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 336f1747703SBarry Smith #endif 33755a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 338289bc588SBarry Smith } 339289bc588SBarry Smith } 340289bc588SBarry Smith if (flag & SOR_BACKWARD_SWEEP) { 341289bc588SBarry Smith for ( i=m-1; i>=0; i-- ) { 342aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 343f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 344f1747703SBarry Smith not happy about returning a double complex */ 345f1747703SBarry Smith int _i; 346f1747703SBarry Smith Scalar sum = b[i]; 347f1747703SBarry Smith for ( _i=0; _i<m; _i++ ) { 3483f6de6efSSatish Balay sum -= PetscConj(v[i+_i*m])*x[_i]; 349f1747703SBarry Smith } 350f1747703SBarry Smith xt = sum; 351f1747703SBarry Smith #else 352289bc588SBarry Smith xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 353f1747703SBarry Smith #endif 35455a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 355289bc588SBarry Smith } 356289bc588SBarry Smith } 357289bc588SBarry Smith } 358600d6d0bSBarry Smith ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 3597a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 3603a40ed3dSBarry Smith PetscFunctionReturn(0); 361289bc588SBarry Smith } 362289bc588SBarry Smith 363289bc588SBarry Smith /* -----------------------------------------------------------------*/ 3645615d1e5SSatish Balay #undef __FUNC__ 3655615d1e5SSatish Balay #define __FUNC__ "MatMultTrans_SeqDense" 36644cd7ae7SLois Curfman McInnes int MatMultTrans_SeqDense(Mat A,Vec xx,Vec yy) 367289bc588SBarry Smith { 368c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 369289bc588SBarry Smith Scalar *v = mat->v, *x, *y; 3707a97a34bSBarry Smith int ierr,_One=1;Scalar _DOne=1.0, _DZero=0.0; 3713a40ed3dSBarry Smith 3723a40ed3dSBarry Smith PetscFunctionBegin; 3737a97a34bSBarry Smith if (!mat->m || !mat->n) PetscFunctionReturn(0); 3747a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 3757a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 37648d91487SBarry Smith LAgemv_("T",&(mat->m),&(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DZero,y,&_One); 3777a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 3787a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 37955659b69SBarry Smith PLogFlops(2*mat->m*mat->n - mat->n); 3803a40ed3dSBarry Smith PetscFunctionReturn(0); 381289bc588SBarry Smith } 3826ee01492SSatish Balay 3835615d1e5SSatish Balay #undef __FUNC__ 3845615d1e5SSatish Balay #define __FUNC__ "MatMult_SeqDense" 38544cd7ae7SLois Curfman McInnes int MatMult_SeqDense(Mat A,Vec xx,Vec yy) 386289bc588SBarry Smith { 387c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 388289bc588SBarry Smith Scalar *v = mat->v, *x, *y; 3897a97a34bSBarry Smith int ierr,_One=1;Scalar _DOne=1.0, _DZero=0.0; 3903a40ed3dSBarry Smith 3913a40ed3dSBarry Smith PetscFunctionBegin; 3927a97a34bSBarry Smith if (!mat->m || !mat->n) PetscFunctionReturn(0); 3937a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 3947a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 39548d91487SBarry Smith LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DZero,y,&_One); 3967a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 3977a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 39855659b69SBarry Smith PLogFlops(2*mat->m*mat->n - mat->m); 3993a40ed3dSBarry Smith PetscFunctionReturn(0); 400289bc588SBarry Smith } 4016ee01492SSatish Balay 4025615d1e5SSatish Balay #undef __FUNC__ 4035615d1e5SSatish Balay #define __FUNC__ "MatMultAdd_SeqDense" 40444cd7ae7SLois Curfman McInnes int MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 405289bc588SBarry Smith { 406c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 407600d6d0bSBarry Smith Scalar *v = mat->v, *x, *y; 4087a97a34bSBarry Smith int ierr,_One=1; Scalar _DOne=1.0; 4093a40ed3dSBarry Smith 4103a40ed3dSBarry Smith PetscFunctionBegin; 4117a97a34bSBarry Smith if (!mat->m || !mat->n) PetscFunctionReturn(0); 412600d6d0bSBarry Smith if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 4137a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 4147a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 41548d91487SBarry Smith LAgemv_( "N", &(mat->m), &(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DOne,y,&_One); 4167a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 4177a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 41855659b69SBarry Smith PLogFlops(2*mat->m*mat->n); 4193a40ed3dSBarry Smith PetscFunctionReturn(0); 420289bc588SBarry Smith } 4216ee01492SSatish Balay 4225615d1e5SSatish Balay #undef __FUNC__ 4235615d1e5SSatish Balay #define __FUNC__ "MatMultTransAdd_SeqDense" 42444cd7ae7SLois Curfman McInnes int MatMultTransAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 425289bc588SBarry Smith { 426c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 427600d6d0bSBarry Smith Scalar *v = mat->v, *x, *y; 4287a97a34bSBarry Smith int ierr,_One=1; 4296abc6512SBarry Smith Scalar _DOne=1.0; 4303a40ed3dSBarry Smith 4313a40ed3dSBarry Smith PetscFunctionBegin; 4327a97a34bSBarry Smith if (!mat->m || !mat->n) PetscFunctionReturn(0); 433600d6d0bSBarry Smith if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 4347a97a34bSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 4357a97a34bSBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 43648d91487SBarry Smith LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DOne,y,&_One); 4377a97a34bSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 4387a97a34bSBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 43955659b69SBarry Smith PLogFlops(2*mat->m*mat->n); 4403a40ed3dSBarry Smith PetscFunctionReturn(0); 441289bc588SBarry Smith } 442289bc588SBarry Smith 443289bc588SBarry Smith /* -----------------------------------------------------------------*/ 4445615d1e5SSatish Balay #undef __FUNC__ 4455615d1e5SSatish Balay #define __FUNC__ "MatGetRow_SeqDense" 4468f6be9afSLois Curfman McInnes int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals) 447289bc588SBarry Smith { 448c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 449289bc588SBarry Smith Scalar *v; 450289bc588SBarry Smith int i; 45167e560aaSBarry Smith 4523a40ed3dSBarry Smith PetscFunctionBegin; 453289bc588SBarry Smith *ncols = mat->n; 454289bc588SBarry Smith if (cols) { 45545d48df9SBarry Smith *cols = (int *) PetscMalloc((mat->n+1)*sizeof(int));CHKPTRQ(*cols); 45673c3ce57SBarry Smith for ( i=0; i<mat->n; i++ ) (*cols)[i] = i; 457289bc588SBarry Smith } 458289bc588SBarry Smith if (vals) { 45945d48df9SBarry Smith *vals = (Scalar *) PetscMalloc((mat->n+1)*sizeof(Scalar));CHKPTRQ(*vals); 460289bc588SBarry Smith v = mat->v + row; 46173c3ce57SBarry Smith for ( i=0; i<mat->n; i++ ) {(*vals)[i] = *v; v += mat->m;} 462289bc588SBarry Smith } 4633a40ed3dSBarry Smith PetscFunctionReturn(0); 464289bc588SBarry Smith } 4656ee01492SSatish Balay 4665615d1e5SSatish Balay #undef __FUNC__ 4675615d1e5SSatish Balay #define __FUNC__ "MatRestoreRow_SeqDense" 4688f6be9afSLois Curfman McInnes int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals) 469289bc588SBarry Smith { 470606d414cSSatish Balay int ierr; 471606d414cSSatish Balay PetscFunctionBegin; 472606d414cSSatish Balay if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);} 473606d414cSSatish Balay if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); } 4743a40ed3dSBarry Smith PetscFunctionReturn(0); 475289bc588SBarry Smith } 476289bc588SBarry Smith /* ----------------------------------------------------------------*/ 4775615d1e5SSatish Balay #undef __FUNC__ 4785615d1e5SSatish Balay #define __FUNC__ "MatSetValues_SeqDense" 4798f6be9afSLois Curfman McInnes int MatSetValues_SeqDense(Mat A,int m,int *indexm,int n, 480e8d4e0b9SBarry Smith int *indexn,Scalar *v,InsertMode addv) 481289bc588SBarry Smith { 482c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 483289bc588SBarry Smith int i,j; 484d6dfbf8fSBarry Smith 4853a40ed3dSBarry Smith PetscFunctionBegin; 486289bc588SBarry Smith if (!mat->roworiented) { 487dbb450caSBarry Smith if (addv == INSERT_VALUES) { 488289bc588SBarry Smith for ( j=0; j<n; j++ ) { 4895ef9f2a5SBarry Smith if (indexn[j] < 0) {v += m; continue;} 490aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 491a8c6a408SBarry Smith if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 49258804f6dSBarry Smith #endif 493289bc588SBarry Smith for ( i=0; i<m; i++ ) { 4945ef9f2a5SBarry Smith if (indexm[i] < 0) {v++; continue;} 495aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 496a8c6a408SBarry Smith if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 49758804f6dSBarry Smith #endif 498289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 499289bc588SBarry Smith } 500289bc588SBarry Smith } 5013a40ed3dSBarry Smith } else { 502289bc588SBarry Smith for ( j=0; j<n; j++ ) { 5035ef9f2a5SBarry Smith if (indexn[j] < 0) {v += m; continue;} 504aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 505a8c6a408SBarry Smith if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 50658804f6dSBarry Smith #endif 507289bc588SBarry Smith for ( i=0; i<m; i++ ) { 5085ef9f2a5SBarry Smith if (indexm[i] < 0) {v++; continue;} 509aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 510a8c6a408SBarry Smith if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 51158804f6dSBarry Smith #endif 512289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 513289bc588SBarry Smith } 514289bc588SBarry Smith } 515289bc588SBarry Smith } 5163a40ed3dSBarry Smith } else { 517dbb450caSBarry Smith if (addv == INSERT_VALUES) { 518e8d4e0b9SBarry Smith for ( i=0; i<m; i++ ) { 5195ef9f2a5SBarry Smith if (indexm[i] < 0) { v += n; continue;} 520aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 521a8c6a408SBarry Smith if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 52258804f6dSBarry Smith #endif 523e8d4e0b9SBarry Smith for ( j=0; j<n; j++ ) { 5245ef9f2a5SBarry Smith if (indexn[j] < 0) { v++; continue;} 525aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 526a8c6a408SBarry Smith if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 52758804f6dSBarry Smith #endif 528e8d4e0b9SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 529e8d4e0b9SBarry Smith } 530e8d4e0b9SBarry Smith } 5313a40ed3dSBarry Smith } else { 532289bc588SBarry Smith for ( i=0; i<m; i++ ) { 5335ef9f2a5SBarry Smith if (indexm[i] < 0) { v += n; continue;} 534aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 535a8c6a408SBarry Smith if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 53658804f6dSBarry Smith #endif 537289bc588SBarry Smith for ( j=0; j<n; j++ ) { 5385ef9f2a5SBarry Smith if (indexn[j] < 0) { v++; continue;} 539aa482453SBarry Smith #if defined(PETSC_USE_BOPT_g) 540a8c6a408SBarry Smith if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 54158804f6dSBarry Smith #endif 542289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 543289bc588SBarry Smith } 544289bc588SBarry Smith } 545289bc588SBarry Smith } 546e8d4e0b9SBarry Smith } 5473a40ed3dSBarry Smith PetscFunctionReturn(0); 548289bc588SBarry Smith } 549e8d4e0b9SBarry Smith 5505615d1e5SSatish Balay #undef __FUNC__ 5515615d1e5SSatish Balay #define __FUNC__ "MatGetValues_SeqDense" 5528f6be9afSLois Curfman McInnes int MatGetValues_SeqDense(Mat A,int m,int *indexm,int n,int *indexn,Scalar *v) 553ae80bb75SLois Curfman McInnes { 554ae80bb75SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 555ae80bb75SLois Curfman McInnes int i, j; 556ae80bb75SLois Curfman McInnes Scalar *vpt = v; 557ae80bb75SLois Curfman McInnes 5583a40ed3dSBarry Smith PetscFunctionBegin; 559ae80bb75SLois Curfman McInnes /* row-oriented output */ 560ae80bb75SLois Curfman McInnes for ( i=0; i<m; i++ ) { 561ae80bb75SLois Curfman McInnes for ( j=0; j<n; j++ ) { 562ae80bb75SLois Curfman McInnes *vpt++ = mat->v[indexn[j]*mat->m + indexm[i]]; 563ae80bb75SLois Curfman McInnes } 564ae80bb75SLois Curfman McInnes } 5653a40ed3dSBarry Smith PetscFunctionReturn(0); 566ae80bb75SLois Curfman McInnes } 567ae80bb75SLois Curfman McInnes 568289bc588SBarry Smith /* -----------------------------------------------------------------*/ 569289bc588SBarry Smith 57077c4ece6SBarry Smith #include "sys.h" 57188e32edaSLois Curfman McInnes 5725615d1e5SSatish Balay #undef __FUNC__ 5735615d1e5SSatish Balay #define __FUNC__ "MatLoad_SeqDense" 57419bcc07fSBarry Smith int MatLoad_SeqDense(Viewer viewer,MatType type,Mat *A) 57588e32edaSLois Curfman McInnes { 57688e32edaSLois Curfman McInnes Mat_SeqDense *a; 57788e32edaSLois Curfman McInnes Mat B; 57888e32edaSLois Curfman McInnes int *scols, i, j, nz, ierr, fd, header[4], size; 57988e32edaSLois Curfman McInnes int *rowlengths = 0, M, N, *cols; 5804a41a4d3SSatish Balay Scalar *vals, *svals, *v,*w; 58119bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 58288e32edaSLois Curfman McInnes 5833a40ed3dSBarry Smith PetscFunctionBegin; 584d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 585a8c6a408SBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"view must have one processor"); 58690ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 5870752156aSBarry Smith ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 588a8c6a408SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"Not matrix object"); 58988e32edaSLois Curfman McInnes M = header[1]; N = header[2]; nz = header[3]; 59088e32edaSLois Curfman McInnes 59190ace30eSBarry Smith if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */ 59290ace30eSBarry Smith ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A);CHKERRQ(ierr); 59390ace30eSBarry Smith B = *A; 59490ace30eSBarry Smith a = (Mat_SeqDense *) B->data; 5954a41a4d3SSatish Balay v = a->v; 5964a41a4d3SSatish Balay /* Allocate some temp space to read in the values and then flip them 5974a41a4d3SSatish Balay from row major to column major */ 5984a41a4d3SSatish Balay w = (Scalar *) PetscMalloc((M*N+1)*sizeof(Scalar));CHKPTRQ(w); 59990ace30eSBarry Smith /* read in nonzero values */ 6004a41a4d3SSatish Balay ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr); 6014a41a4d3SSatish Balay /* now flip the values and store them in the matrix*/ 6024a41a4d3SSatish Balay for ( j=0; j<N; j++ ) { 6034a41a4d3SSatish Balay for ( i=0; i<M; i++ ) { 6044a41a4d3SSatish Balay *v++ =w[i*N+j]; 6054a41a4d3SSatish Balay } 6064a41a4d3SSatish Balay } 607606d414cSSatish Balay ierr = PetscFree(w);CHKERRQ(ierr); 6086d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 6096d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 61090ace30eSBarry Smith } else { 61188e32edaSLois Curfman McInnes /* read row lengths */ 61245d48df9SBarry Smith rowlengths = (int*) PetscMalloc( (M+1)*sizeof(int) );CHKPTRQ(rowlengths); 6130752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 61488e32edaSLois Curfman McInnes 61588e32edaSLois Curfman McInnes /* create our matrix */ 616b4fd4287SBarry Smith ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A);CHKERRQ(ierr); 61788e32edaSLois Curfman McInnes B = *A; 61888e32edaSLois Curfman McInnes a = (Mat_SeqDense *) B->data; 61988e32edaSLois Curfman McInnes v = a->v; 62088e32edaSLois Curfman McInnes 62188e32edaSLois Curfman McInnes /* read column indices and nonzeros */ 62245d48df9SBarry Smith cols = scols = (int *) PetscMalloc( (nz+1)*sizeof(int) );CHKPTRQ(cols); 6230752156aSBarry Smith ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr); 62445d48df9SBarry Smith vals = svals = (Scalar *) PetscMalloc( (nz+1)*sizeof(Scalar) );CHKPTRQ(vals); 6250752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 62688e32edaSLois Curfman McInnes 62788e32edaSLois Curfman McInnes /* insert into matrix */ 62888e32edaSLois Curfman McInnes for ( i=0; i<M; i++ ) { 62988e32edaSLois Curfman McInnes for ( j=0; j<rowlengths[i]; j++ ) v[i+M*scols[j]] = svals[j]; 63088e32edaSLois Curfman McInnes svals += rowlengths[i]; scols += rowlengths[i]; 63188e32edaSLois Curfman McInnes } 632606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 633606d414cSSatish Balay ierr = PetscFree(cols);CHKERRQ(ierr); 634606d414cSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 63588e32edaSLois Curfman McInnes 6366d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 6376d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 63890ace30eSBarry Smith } 6393a40ed3dSBarry Smith PetscFunctionReturn(0); 64088e32edaSLois Curfman McInnes } 64188e32edaSLois Curfman McInnes 64277c4ece6SBarry Smith #include "sys.h" 643932b0c3eSLois Curfman McInnes 6445615d1e5SSatish Balay #undef __FUNC__ 6455615d1e5SSatish Balay #define __FUNC__ "MatView_SeqDense_ASCII" 646932b0c3eSLois Curfman McInnes static int MatView_SeqDense_ASCII(Mat A,Viewer viewer) 647289bc588SBarry Smith { 648932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense *) A->data; 649932b0c3eSLois Curfman McInnes int ierr, i, j, format; 650932b0c3eSLois Curfman McInnes char *outputname; 651932b0c3eSLois Curfman McInnes Scalar *v; 652932b0c3eSLois Curfman McInnes 6533a40ed3dSBarry Smith PetscFunctionBegin; 65477ed5343SBarry Smith ierr = ViewerGetOutputname(viewer,&outputname);CHKERRQ(ierr); 655888f2ed8SSatish Balay ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr); 656639f9d9dSBarry Smith if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) { 6573a40ed3dSBarry Smith PetscFunctionReturn(0); /* do nothing for now */ 658*6831982aSBarry Smith } else if (format == VIEWER_FORMAT_ASCII_COMMON) { 659*6831982aSBarry Smith ierr = ViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 66080cd9d93SLois Curfman McInnes for ( i=0; i<a->m; i++ ) { 66144cd7ae7SLois Curfman McInnes v = a->v + i; 662*6831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer,"row %d:",i);CHKERRQ(ierr); 66380cd9d93SLois Curfman McInnes for ( j=0; j<a->n; j++ ) { 664aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 665*6831982aSBarry Smith if (PetscReal(*v) != 0.0 && PetscImaginary(*v) != 0.0) { 666*6831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g + %g i",j,PetscReal(*v),PetscImaginary(*v));CHKERRQ(ierr); 667*6831982aSBarry Smith } else if (PetscReal(*v)) { 668*6831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g ",j,PetscReal(*v));CHKERRQ(ierr); 669*6831982aSBarry Smith } 67080cd9d93SLois Curfman McInnes #else 671*6831982aSBarry Smith if (*v) { 672*6831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer," %d %g ",j,*v);CHKERRQ(ierr); 673*6831982aSBarry Smith } 67480cd9d93SLois Curfman McInnes #endif 675d64ed03dSBarry Smith v += a->m; 67680cd9d93SLois Curfman McInnes } 677*6831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 67880cd9d93SLois Curfman McInnes } 679*6831982aSBarry Smith ierr = ViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 6803a40ed3dSBarry Smith } else { 681*6831982aSBarry Smith ierr = ViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 682aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 68347989497SBarry Smith int allreal = 1; 68447989497SBarry Smith /* determine if matrix has all real values */ 68547989497SBarry Smith v = a->v; 68647989497SBarry Smith for ( i=0; i<a->m*a->n; i++ ) { 687e20fef11SSatish Balay if (PetscImaginary(v[i])) { allreal = 0; break ;} 68847989497SBarry Smith } 68947989497SBarry Smith #endif 690932b0c3eSLois Curfman McInnes for ( i=0; i<a->m; i++ ) { 691932b0c3eSLois Curfman McInnes v = a->v + i; 692932b0c3eSLois Curfman McInnes for ( j=0; j<a->n; j++ ) { 693aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 69447989497SBarry Smith if (allreal) { 695*6831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer,"%6.4e ",PetscReal(*v));CHKERRQ(ierr); v += a->m; 69647989497SBarry Smith } else { 697*6831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer,"%6.4e + %6.4e i ",PetscReal(*v),PetscImaginary(*v));CHKERRQ(ierr); v += a->m; 69847989497SBarry Smith } 699289bc588SBarry Smith #else 700*6831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer,"%6.4e ",*v);CHKERRQ(ierr); v += a->m; 701289bc588SBarry Smith #endif 702289bc588SBarry Smith } 703*6831982aSBarry Smith ierr = ViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 704289bc588SBarry Smith } 705*6831982aSBarry Smith ierr = ViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 706da3a660dSBarry Smith } 707*6831982aSBarry Smith ierr = ViewerFlush(viewer);CHKERRQ(ierr); 7083a40ed3dSBarry Smith PetscFunctionReturn(0); 709289bc588SBarry Smith } 710289bc588SBarry Smith 7115615d1e5SSatish Balay #undef __FUNC__ 7125615d1e5SSatish Balay #define __FUNC__ "MatView_SeqDense_Binary" 713932b0c3eSLois Curfman McInnes static int MatView_SeqDense_Binary(Mat A,Viewer viewer) 714932b0c3eSLois Curfman McInnes { 715932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense *) A->data; 716932b0c3eSLois Curfman McInnes int ict, j, n = a->n, m = a->m, i, fd, *col_lens, ierr, nz = m*n; 71790ace30eSBarry Smith int format; 71890ace30eSBarry Smith Scalar *v, *anonz,*vals; 719932b0c3eSLois Curfman McInnes 7203a40ed3dSBarry Smith PetscFunctionBegin; 72190ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 72290ace30eSBarry Smith 72390ace30eSBarry Smith ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr); 72402cad45dSBarry Smith if (format == VIEWER_FORMAT_BINARY_NATIVE) { 72590ace30eSBarry Smith /* store the matrix as a dense matrix */ 72690ace30eSBarry Smith col_lens = (int *) PetscMalloc( 4*sizeof(int) );CHKPTRQ(col_lens); 72790ace30eSBarry Smith col_lens[0] = MAT_COOKIE; 72890ace30eSBarry Smith col_lens[1] = m; 72990ace30eSBarry Smith col_lens[2] = n; 73090ace30eSBarry Smith col_lens[3] = MATRIX_BINARY_FORMAT_DENSE; 7310752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,1);CHKERRQ(ierr); 732606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 73390ace30eSBarry Smith 73490ace30eSBarry Smith /* write out matrix, by rows */ 73545d48df9SBarry Smith vals = (Scalar *) PetscMalloc((m*n+1)*sizeof(Scalar));CHKPTRQ(vals); 73690ace30eSBarry Smith v = a->v; 73790ace30eSBarry Smith for ( i=0; i<m; i++ ) { 73890ace30eSBarry Smith for ( j=0; j<n; j++ ) { 73990ace30eSBarry Smith vals[i + j*m] = *v++; 74090ace30eSBarry Smith } 74190ace30eSBarry Smith } 7420752156aSBarry Smith ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,0);CHKERRQ(ierr); 743606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 74490ace30eSBarry Smith } else { 7450452661fSBarry Smith col_lens = (int *) PetscMalloc( (4+nz)*sizeof(int) );CHKPTRQ(col_lens); 746932b0c3eSLois Curfman McInnes col_lens[0] = MAT_COOKIE; 747932b0c3eSLois Curfman McInnes col_lens[1] = m; 748932b0c3eSLois Curfman McInnes col_lens[2] = n; 749932b0c3eSLois Curfman McInnes col_lens[3] = nz; 750932b0c3eSLois Curfman McInnes 751932b0c3eSLois Curfman McInnes /* store lengths of each row and write (including header) to file */ 752932b0c3eSLois Curfman McInnes for ( i=0; i<m; i++ ) col_lens[4+i] = n; 7530752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,1);CHKERRQ(ierr); 754932b0c3eSLois Curfman McInnes 755932b0c3eSLois Curfman McInnes /* Possibly should write in smaller increments, not whole matrix at once? */ 756932b0c3eSLois Curfman McInnes /* store column indices (zero start index) */ 757932b0c3eSLois Curfman McInnes ict = 0; 758932b0c3eSLois Curfman McInnes for ( i=0; i<m; i++ ) { 759932b0c3eSLois Curfman McInnes for ( j=0; j<n; j++ ) col_lens[ict++] = j; 760932b0c3eSLois Curfman McInnes } 7610752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,0);CHKERRQ(ierr); 762606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 763932b0c3eSLois Curfman McInnes 764932b0c3eSLois Curfman McInnes /* store nonzero values */ 76545d48df9SBarry Smith anonz = (Scalar *) PetscMalloc((nz+1)*sizeof(Scalar));CHKPTRQ(anonz); 766932b0c3eSLois Curfman McInnes ict = 0; 767932b0c3eSLois Curfman McInnes for ( i=0; i<m; i++ ) { 768932b0c3eSLois Curfman McInnes v = a->v + i; 769932b0c3eSLois Curfman McInnes for ( j=0; j<n; j++ ) { 770932b0c3eSLois Curfman McInnes anonz[ict++] = *v; v += a->m; 771932b0c3eSLois Curfman McInnes } 772932b0c3eSLois Curfman McInnes } 7730752156aSBarry Smith ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,0);CHKERRQ(ierr); 774606d414cSSatish Balay ierr = PetscFree(anonz);CHKERRQ(ierr); 77590ace30eSBarry Smith } 7763a40ed3dSBarry Smith PetscFunctionReturn(0); 777932b0c3eSLois Curfman McInnes } 778932b0c3eSLois Curfman McInnes 7795615d1e5SSatish Balay #undef __FUNC__ 7805615d1e5SSatish Balay #define __FUNC__ "MatView_SeqDense" 781c6e7dd08SBarry Smith int MatView_SeqDense(Mat A,Viewer viewer) 782932b0c3eSLois Curfman McInnes { 783932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*) A->data; 784bcd2baecSBarry Smith int ierr; 785*6831982aSBarry Smith PetscTruth issocket,isascii,isbinary; 786932b0c3eSLois Curfman McInnes 7873a40ed3dSBarry Smith PetscFunctionBegin; 788*6831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,SOCKET_VIEWER,&issocket);CHKERRQ(ierr); 789*6831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr); 790*6831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,BINARY_VIEWER,&isbinary);CHKERRQ(ierr); 7910f5bd95cSBarry Smith 7920f5bd95cSBarry Smith if (issocket) { 7937b2a1423SBarry Smith ierr = ViewerSocketPutScalar_Private(viewer,a->m,a->n,a->v);CHKERRQ(ierr); 7940f5bd95cSBarry Smith } else if (isascii) { 7953a40ed3dSBarry Smith ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr); 7960f5bd95cSBarry Smith } else if (isbinary) { 7973a40ed3dSBarry Smith ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr); 7985cd90555SBarry Smith } else { 7990f5bd95cSBarry Smith SETERRQ1(1,1,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name); 800932b0c3eSLois Curfman McInnes } 8013a40ed3dSBarry Smith PetscFunctionReturn(0); 802932b0c3eSLois Curfman McInnes } 803289bc588SBarry Smith 8045615d1e5SSatish Balay #undef __FUNC__ 8055615d1e5SSatish Balay #define __FUNC__ "MatDestroy_SeqDense" 806c6e7dd08SBarry Smith int MatDestroy_SeqDense(Mat mat) 807289bc588SBarry Smith { 808ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense *) mat->data; 80990f02eecSBarry Smith int ierr; 81090f02eecSBarry Smith 8113a40ed3dSBarry Smith PetscFunctionBegin; 81294d884c6SBarry Smith 81394d884c6SBarry Smith if (mat->mapping) { 81494d884c6SBarry Smith ierr = ISLocalToGlobalMappingDestroy(mat->mapping);CHKERRQ(ierr); 81594d884c6SBarry Smith } 81694d884c6SBarry Smith if (mat->bmapping) { 81794d884c6SBarry Smith ierr = ISLocalToGlobalMappingDestroy(mat->bmapping);CHKERRQ(ierr); 81894d884c6SBarry Smith } 819aa482453SBarry Smith #if defined(PETSC_USE_LOG) 820c6e7dd08SBarry Smith PLogObjectState((PetscObject)mat,"Rows %d Cols %d",l->m,l->n); 821a5a9c739SBarry Smith #endif 822606d414cSSatish Balay if (l->pivots) {ierr = PetscFree(l->pivots);CHKERRQ(ierr);} 823606d414cSSatish Balay if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);} 824606d414cSSatish Balay ierr = PetscFree(l);CHKERRQ(ierr); 825a5ae1ecdSBarry Smith if (mat->rmap) { 826a5ae1ecdSBarry Smith ierr = MapDestroy(mat->rmap);CHKERRQ(ierr); 827a5ae1ecdSBarry Smith } 828a5ae1ecdSBarry Smith if (mat->cmap) { 829a5ae1ecdSBarry Smith ierr = MapDestroy(mat->cmap);CHKERRQ(ierr); 83090f02eecSBarry Smith } 831a5a9c739SBarry Smith PLogObjectDestroy(mat); 8320452661fSBarry Smith PetscHeaderDestroy(mat); 8333a40ed3dSBarry Smith PetscFunctionReturn(0); 834289bc588SBarry Smith } 835289bc588SBarry Smith 8365615d1e5SSatish Balay #undef __FUNC__ 8375615d1e5SSatish Balay #define __FUNC__ "MatTranspose_SeqDense" 8388f6be9afSLois Curfman McInnes int MatTranspose_SeqDense(Mat A,Mat *matout) 839289bc588SBarry Smith { 840c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 841549d3d68SSatish Balay int k, j, m, n, ierr; 842d3e5ee88SLois Curfman McInnes Scalar *v, tmp; 84348b35521SBarry Smith 8443a40ed3dSBarry Smith PetscFunctionBegin; 845d3e5ee88SLois Curfman McInnes v = mat->v; m = mat->m; n = mat->n; 8463638b69dSLois Curfman McInnes if (matout == PETSC_NULL) { /* in place transpose */ 847d64ed03dSBarry Smith if (m != n) { /* malloc temp to hold transpose */ 848d64ed03dSBarry Smith Scalar *w = (Scalar *) PetscMalloc((m+1)*(n+1)*sizeof(Scalar));CHKPTRQ(w); 849d64ed03dSBarry Smith for ( j=0; j<m; j++ ) { 850d64ed03dSBarry Smith for ( k=0; k<n; k++ ) { 851d64ed03dSBarry Smith w[k + j*n] = v[j + k*m]; 852d64ed03dSBarry Smith } 853d64ed03dSBarry Smith } 854549d3d68SSatish Balay ierr = PetscMemcpy(v,w,m*n*sizeof(Scalar));CHKERRQ(ierr); 855606d414cSSatish Balay ierr = PetscFree(w);CHKERRQ(ierr); 856d64ed03dSBarry Smith } else { 857d3e5ee88SLois Curfman McInnes for ( j=0; j<m; j++ ) { 858289bc588SBarry Smith for ( k=0; k<j; k++ ) { 859d3e5ee88SLois Curfman McInnes tmp = v[j + k*n]; 860d3e5ee88SLois Curfman McInnes v[j + k*n] = v[k + j*n]; 861d3e5ee88SLois Curfman McInnes v[k + j*n] = tmp; 862289bc588SBarry Smith } 863289bc588SBarry Smith } 864d64ed03dSBarry Smith } 8653a40ed3dSBarry Smith } else { /* out-of-place transpose */ 866d3e5ee88SLois Curfman McInnes Mat tmat; 867ec8511deSBarry Smith Mat_SeqDense *tmatd; 868d3e5ee88SLois Curfman McInnes Scalar *v2; 869b4fd4287SBarry Smith ierr = MatCreateSeqDense(A->comm,mat->n,mat->m,PETSC_NULL,&tmat);CHKERRQ(ierr); 870ec8511deSBarry Smith tmatd = (Mat_SeqDense *) tmat->data; 8710de55854SLois Curfman McInnes v = mat->v; v2 = tmatd->v; 872d3e5ee88SLois Curfman McInnes for ( j=0; j<n; j++ ) { 8730de55854SLois Curfman McInnes for ( k=0; k<m; k++ ) v2[j + k*n] = v[k + j*m]; 874d3e5ee88SLois Curfman McInnes } 8756d4a8577SBarry Smith ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 8766d4a8577SBarry Smith ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 877d3e5ee88SLois Curfman McInnes *matout = tmat; 87848b35521SBarry Smith } 8793a40ed3dSBarry Smith PetscFunctionReturn(0); 880289bc588SBarry Smith } 881289bc588SBarry Smith 8825615d1e5SSatish Balay #undef __FUNC__ 8835615d1e5SSatish Balay #define __FUNC__ "MatEqual_SeqDense" 8848f6be9afSLois Curfman McInnes int MatEqual_SeqDense(Mat A1,Mat A2, PetscTruth *flg) 885289bc588SBarry Smith { 886c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat1 = (Mat_SeqDense *) A1->data; 887c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat2 = (Mat_SeqDense *) A2->data; 888289bc588SBarry Smith int i; 889289bc588SBarry Smith Scalar *v1 = mat1->v, *v2 = mat2->v; 8909ea5d5aeSSatish Balay 8913a40ed3dSBarry Smith PetscFunctionBegin; 8923a40ed3dSBarry Smith if (A2->type != MATSEQDENSE) SETERRQ(PETSC_ERR_SUP,0,"Matrices should be of type MATSEQDENSE"); 8933a40ed3dSBarry Smith if (mat1->m != mat2->m) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 8943a40ed3dSBarry Smith if (mat1->n != mat2->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 895289bc588SBarry Smith for ( i=0; i<mat1->m*mat1->n; i++ ) { 8963a40ed3dSBarry Smith if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 897289bc588SBarry Smith v1++; v2++; 898289bc588SBarry Smith } 89977c4ece6SBarry Smith *flg = PETSC_TRUE; 9003a40ed3dSBarry Smith PetscFunctionReturn(0); 901289bc588SBarry Smith } 902289bc588SBarry Smith 9035615d1e5SSatish Balay #undef __FUNC__ 9045615d1e5SSatish Balay #define __FUNC__ "MatGetDiagonal_SeqDense" 9058f6be9afSLois Curfman McInnes int MatGetDiagonal_SeqDense(Mat A,Vec v) 906289bc588SBarry Smith { 907c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 9087a97a34bSBarry Smith int ierr,i, n, len; 90944cd7ae7SLois Curfman McInnes Scalar *x, zero = 0.0; 91044cd7ae7SLois Curfman McInnes 9113a40ed3dSBarry Smith PetscFunctionBegin; 9127a97a34bSBarry Smith ierr = VecSet(&zero,v);CHKERRQ(ierr); 9137a97a34bSBarry Smith ierr = VecGetSize(v,&n);CHKERRQ(ierr); 914600d6d0bSBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 91544cd7ae7SLois Curfman McInnes len = PetscMin(mat->m,mat->n); 916e3372554SBarry Smith if (n != mat->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Nonconforming mat and vec"); 91744cd7ae7SLois Curfman McInnes for ( i=0; i<len; i++ ) { 918289bc588SBarry Smith x[i] = mat->v[i*mat->m + i]; 919289bc588SBarry Smith } 9207a97a34bSBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 9213a40ed3dSBarry Smith PetscFunctionReturn(0); 922289bc588SBarry Smith } 923289bc588SBarry Smith 9245615d1e5SSatish Balay #undef __FUNC__ 9255615d1e5SSatish Balay #define __FUNC__ "MatDiagonalScale_SeqDense" 9268f6be9afSLois Curfman McInnes int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr) 927289bc588SBarry Smith { 928c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 929da3a660dSBarry Smith Scalar *l,*r,x,*v; 9307a97a34bSBarry Smith int ierr,i,j,m = mat->m, n = mat->n; 93155659b69SBarry Smith 9323a40ed3dSBarry Smith PetscFunctionBegin; 93328988994SBarry Smith if (ll) { 9347a97a34bSBarry Smith ierr = VecGetSize(ll,&m);CHKERRQ(ierr); 935600d6d0bSBarry Smith ierr = VecGetArray(ll,&l);CHKERRQ(ierr); 936e3372554SBarry Smith if (m != mat->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Left scaling vec wrong size"); 937da3a660dSBarry Smith for ( i=0; i<m; i++ ) { 938da3a660dSBarry Smith x = l[i]; 939da3a660dSBarry Smith v = mat->v + i; 940da3a660dSBarry Smith for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;} 941da3a660dSBarry Smith } 9427a97a34bSBarry Smith ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr); 94344cd7ae7SLois Curfman McInnes PLogFlops(n*m); 944da3a660dSBarry Smith } 94528988994SBarry Smith if (rr) { 9467a97a34bSBarry Smith ierr = VecGetSize(rr,&n);CHKERRQ(ierr); 947600d6d0bSBarry Smith ierr = VecGetArray(rr,&r);CHKERRQ(ierr); 948e3372554SBarry Smith if (n != mat->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Right scaling vec wrong size"); 949da3a660dSBarry Smith for ( i=0; i<n; i++ ) { 950da3a660dSBarry Smith x = r[i]; 951da3a660dSBarry Smith v = mat->v + i*m; 952da3a660dSBarry Smith for ( j=0; j<m; j++ ) { (*v++) *= x;} 953da3a660dSBarry Smith } 9547a97a34bSBarry Smith ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr); 95544cd7ae7SLois Curfman McInnes PLogFlops(n*m); 956da3a660dSBarry Smith } 9573a40ed3dSBarry Smith PetscFunctionReturn(0); 958289bc588SBarry Smith } 959289bc588SBarry Smith 9605615d1e5SSatish Balay #undef __FUNC__ 9615615d1e5SSatish Balay #define __FUNC__ "MatNorm_SeqDense" 9628f6be9afSLois Curfman McInnes int MatNorm_SeqDense(Mat A,NormType type,double *norm) 963289bc588SBarry Smith { 964c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 965289bc588SBarry Smith Scalar *v = mat->v; 966289bc588SBarry Smith double sum = 0.0; 967289bc588SBarry Smith int i, j; 96855659b69SBarry Smith 9693a40ed3dSBarry Smith PetscFunctionBegin; 970289bc588SBarry Smith if (type == NORM_FROBENIUS) { 971289bc588SBarry Smith for (i=0; i<mat->n*mat->m; i++ ) { 972aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 9733f6de6efSSatish Balay sum += PetscReal(PetscConj(*v)*(*v)); v++; 974289bc588SBarry Smith #else 975289bc588SBarry Smith sum += (*v)*(*v); v++; 976289bc588SBarry Smith #endif 977289bc588SBarry Smith } 978289bc588SBarry Smith *norm = sqrt(sum); 97955659b69SBarry Smith PLogFlops(2*mat->n*mat->m); 9803a40ed3dSBarry Smith } else if (type == NORM_1) { 981289bc588SBarry Smith *norm = 0.0; 982289bc588SBarry Smith for ( j=0; j<mat->n; j++ ) { 983289bc588SBarry Smith sum = 0.0; 984289bc588SBarry Smith for ( i=0; i<mat->m; i++ ) { 98533a8263dSBarry Smith sum += PetscAbsScalar(*v); v++; 986289bc588SBarry Smith } 987289bc588SBarry Smith if (sum > *norm) *norm = sum; 988289bc588SBarry Smith } 98955659b69SBarry Smith PLogFlops(mat->n*mat->m); 9903a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { 991289bc588SBarry Smith *norm = 0.0; 992289bc588SBarry Smith for ( j=0; j<mat->m; j++ ) { 993289bc588SBarry Smith v = mat->v + j; 994289bc588SBarry Smith sum = 0.0; 995289bc588SBarry Smith for ( i=0; i<mat->n; i++ ) { 996cddf8d76SBarry Smith sum += PetscAbsScalar(*v); v += mat->m; 997289bc588SBarry Smith } 998289bc588SBarry Smith if (sum > *norm) *norm = sum; 999289bc588SBarry Smith } 100055659b69SBarry Smith PLogFlops(mat->n*mat->m); 10013a40ed3dSBarry Smith } else { 1002e3372554SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"No two norm"); 1003289bc588SBarry Smith } 10043a40ed3dSBarry Smith PetscFunctionReturn(0); 1005289bc588SBarry Smith } 1006289bc588SBarry Smith 10075615d1e5SSatish Balay #undef __FUNC__ 10085615d1e5SSatish Balay #define __FUNC__ "MatSetOption_SeqDense" 10098f6be9afSLois Curfman McInnes int MatSetOption_SeqDense(Mat A,MatOption op) 1010289bc588SBarry Smith { 1011c0bbcb79SLois Curfman McInnes Mat_SeqDense *aij = (Mat_SeqDense *) A->data; 101267e560aaSBarry Smith 10133a40ed3dSBarry Smith PetscFunctionBegin; 10146d4a8577SBarry Smith if (op == MAT_ROW_ORIENTED) aij->roworiented = 1; 10156d4a8577SBarry Smith else if (op == MAT_COLUMN_ORIENTED) aij->roworiented = 0; 10166d4a8577SBarry Smith else if (op == MAT_ROWS_SORTED || 1017219d9a1aSLois Curfman McInnes op == MAT_ROWS_UNSORTED || 10186d4a8577SBarry Smith op == MAT_COLUMNS_SORTED || 1019219d9a1aSLois Curfman McInnes op == MAT_COLUMNS_UNSORTED || 10206d4a8577SBarry Smith op == MAT_SYMMETRIC || 10216d4a8577SBarry Smith op == MAT_STRUCTURALLY_SYMMETRIC || 10226d4a8577SBarry Smith op == MAT_NO_NEW_NONZERO_LOCATIONS || 10236d4a8577SBarry Smith op == MAT_YES_NEW_NONZERO_LOCATIONS || 10244787f768SSatish Balay op == MAT_NEW_NONZERO_LOCATION_ERR || 10256d4a8577SBarry Smith op == MAT_NO_NEW_DIAGONALS || 102690f02eecSBarry Smith op == MAT_YES_NEW_DIAGONALS || 1027b51ba29fSSatish Balay op == MAT_IGNORE_OFF_PROC_ENTRIES || 1028b51ba29fSSatish Balay op == MAT_USE_HASH_TABLE) 1029981c4779SBarry Smith PLogInfo(A,"MatSetOption_SeqDense:Option ignored\n"); 10303a40ed3dSBarry Smith else { 10313a40ed3dSBarry Smith SETERRQ(PETSC_ERR_SUP,0,"unknown option"); 10323a40ed3dSBarry Smith } 10333a40ed3dSBarry Smith PetscFunctionReturn(0); 1034289bc588SBarry Smith } 1035289bc588SBarry Smith 10365615d1e5SSatish Balay #undef __FUNC__ 10375615d1e5SSatish Balay #define __FUNC__ "MatZeroEntries_SeqDense" 10388f6be9afSLois Curfman McInnes int MatZeroEntries_SeqDense(Mat A) 10396f0a148fSBarry Smith { 1040ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense *) A->data; 1041549d3d68SSatish Balay int ierr; 10423a40ed3dSBarry Smith 10433a40ed3dSBarry Smith PetscFunctionBegin; 1044549d3d68SSatish Balay ierr = PetscMemzero(l->v,l->m*l->n*sizeof(Scalar));CHKERRQ(ierr); 10453a40ed3dSBarry Smith PetscFunctionReturn(0); 10466f0a148fSBarry Smith } 10476f0a148fSBarry Smith 10485615d1e5SSatish Balay #undef __FUNC__ 10495615d1e5SSatish Balay #define __FUNC__ "MatGetBlockSize_SeqDense" 10508f6be9afSLois Curfman McInnes int MatGetBlockSize_SeqDense(Mat A,int *bs) 10514e220ebcSLois Curfman McInnes { 10523a40ed3dSBarry Smith PetscFunctionBegin; 10534e220ebcSLois Curfman McInnes *bs = 1; 10543a40ed3dSBarry Smith PetscFunctionReturn(0); 10554e220ebcSLois Curfman McInnes } 10564e220ebcSLois Curfman McInnes 10575615d1e5SSatish Balay #undef __FUNC__ 10585615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqDense" 10598f6be9afSLois Curfman McInnes int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag) 10606f0a148fSBarry Smith { 1061ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense *) A->data; 10626abc6512SBarry Smith int n = l->n, i, j,ierr,N, *rows; 10636f0a148fSBarry Smith Scalar *slot; 106455659b69SBarry Smith 10653a40ed3dSBarry Smith PetscFunctionBegin; 106677c4ece6SBarry Smith ierr = ISGetSize(is,&N);CHKERRQ(ierr); 106778b31e54SBarry Smith ierr = ISGetIndices(is,&rows);CHKERRQ(ierr); 10686f0a148fSBarry Smith for ( i=0; i<N; i++ ) { 10696f0a148fSBarry Smith slot = l->v + rows[i]; 10706f0a148fSBarry Smith for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;} 10716f0a148fSBarry Smith } 10726f0a148fSBarry Smith if (diag) { 10736f0a148fSBarry Smith for ( i=0; i<N; i++ ) { 10746f0a148fSBarry Smith slot = l->v + (n+1)*rows[i]; 1075260925b8SBarry Smith *slot = *diag; 10766f0a148fSBarry Smith } 10776f0a148fSBarry Smith } 1078260925b8SBarry Smith ISRestoreIndices(is,&rows); 10793a40ed3dSBarry Smith PetscFunctionReturn(0); 10806f0a148fSBarry Smith } 1081557bce09SLois Curfman McInnes 10825615d1e5SSatish Balay #undef __FUNC__ 10835615d1e5SSatish Balay #define __FUNC__ "MatGetSize_SeqDense" 10848f6be9afSLois Curfman McInnes int MatGetSize_SeqDense(Mat A,int *m,int *n) 1085557bce09SLois Curfman McInnes { 1086c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 10873a40ed3dSBarry Smith 10883a40ed3dSBarry Smith PetscFunctionBegin; 1089557bce09SLois Curfman McInnes *m = mat->m; *n = mat->n; 10903a40ed3dSBarry Smith PetscFunctionReturn(0); 1091557bce09SLois Curfman McInnes } 1092557bce09SLois Curfman McInnes 10935615d1e5SSatish Balay #undef __FUNC__ 10945615d1e5SSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqDense" 10958f6be9afSLois Curfman McInnes int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n) 1096d3e5ee88SLois Curfman McInnes { 1097c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 10983a40ed3dSBarry Smith 10993a40ed3dSBarry Smith PetscFunctionBegin; 1100d3e5ee88SLois Curfman McInnes *m = 0; *n = mat->m; 11013a40ed3dSBarry Smith PetscFunctionReturn(0); 1102d3e5ee88SLois Curfman McInnes } 1103d3e5ee88SLois Curfman McInnes 11045615d1e5SSatish Balay #undef __FUNC__ 11055615d1e5SSatish Balay #define __FUNC__ "MatGetArray_SeqDense" 11068f6be9afSLois Curfman McInnes int MatGetArray_SeqDense(Mat A,Scalar **array) 110764e87e97SBarry Smith { 1108c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 11093a40ed3dSBarry Smith 11103a40ed3dSBarry Smith PetscFunctionBegin; 111164e87e97SBarry Smith *array = mat->v; 11123a40ed3dSBarry Smith PetscFunctionReturn(0); 111364e87e97SBarry Smith } 11140754003eSLois Curfman McInnes 11155615d1e5SSatish Balay #undef __FUNC__ 11165615d1e5SSatish Balay #define __FUNC__ "MatRestoreArray_SeqDense" 11178f6be9afSLois Curfman McInnes int MatRestoreArray_SeqDense(Mat A,Scalar **array) 1118ff14e315SSatish Balay { 11193a40ed3dSBarry Smith PetscFunctionBegin; 112009b544d4SBarry Smith *array = 0; /* user cannot accidently use the array later */ 11213a40ed3dSBarry Smith PetscFunctionReturn(0); 1122ff14e315SSatish Balay } 11230754003eSLois Curfman McInnes 11245615d1e5SSatish Balay #undef __FUNC__ 11255615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrix_SeqDense" 11267b2a1423SBarry Smith static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B) 11270754003eSLois Curfman McInnes { 1128c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 1129182d2002SSatish Balay int i, j, ierr, m = mat->m, *irow, *icol, nrows, ncols; 1130182d2002SSatish Balay Scalar *av, *bv, *v = mat->v; 11310754003eSLois Curfman McInnes Mat newmat; 11320754003eSLois Curfman McInnes 11333a40ed3dSBarry Smith PetscFunctionBegin; 113478b31e54SBarry Smith ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 113578b31e54SBarry Smith ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 113678b31e54SBarry Smith ierr = ISGetSize(isrow,&nrows);CHKERRQ(ierr); 113778b31e54SBarry Smith ierr = ISGetSize(iscol,&ncols);CHKERRQ(ierr); 11380754003eSLois Curfman McInnes 1139182d2002SSatish Balay /* Check submatrixcall */ 1140182d2002SSatish Balay if (scall == MAT_REUSE_MATRIX) { 1141182d2002SSatish Balay int n_cols,n_rows; 1142182d2002SSatish Balay ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr); 1143182d2002SSatish Balay if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Reused submatrix wrong size"); 1144182d2002SSatish Balay newmat = *B; 1145182d2002SSatish Balay } else { 11460754003eSLois Curfman McInnes /* Create and fill new matrix */ 1147b4fd4287SBarry Smith ierr = MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat);CHKERRQ(ierr); 1148182d2002SSatish Balay } 1149182d2002SSatish Balay 1150182d2002SSatish Balay /* Now extract the data pointers and do the copy, column at a time */ 1151182d2002SSatish Balay bv = ((Mat_SeqDense*)newmat->data)->v; 1152182d2002SSatish Balay 1153182d2002SSatish Balay for ( i=0; i<ncols; i++ ) { 1154182d2002SSatish Balay av = v + m*icol[i]; 1155182d2002SSatish Balay for (j=0; j<nrows; j++ ) { 1156182d2002SSatish Balay *bv++ = av[irow[j]]; 11570754003eSLois Curfman McInnes } 11580754003eSLois Curfman McInnes } 1159182d2002SSatish Balay 1160182d2002SSatish Balay /* Assemble the matrices so that the correct flags are set */ 11616d4a8577SBarry Smith ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 11626d4a8577SBarry Smith ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 11630754003eSLois Curfman McInnes 11640754003eSLois Curfman McInnes /* Free work space */ 116578b31e54SBarry Smith ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 116678b31e54SBarry Smith ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 1167182d2002SSatish Balay *B = newmat; 11683a40ed3dSBarry Smith PetscFunctionReturn(0); 11690754003eSLois Curfman McInnes } 11700754003eSLois Curfman McInnes 11715615d1e5SSatish Balay #undef __FUNC__ 11725615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrices_SeqDense" 11737b2a1423SBarry Smith int MatGetSubMatrices_SeqDense(Mat A,int n, IS *irow,IS *icol,MatReuse scall,Mat **B) 1174905e6a2fSBarry Smith { 1175905e6a2fSBarry Smith int ierr,i; 1176905e6a2fSBarry Smith 11773a40ed3dSBarry Smith PetscFunctionBegin; 1178905e6a2fSBarry Smith if (scall == MAT_INITIAL_MATRIX) { 1179905e6a2fSBarry Smith *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) );CHKPTRQ(*B); 1180905e6a2fSBarry Smith } 1181905e6a2fSBarry Smith 1182905e6a2fSBarry Smith for ( i=0; i<n; i++ ) { 11836a6a5d1dSBarry Smith ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr); 1184905e6a2fSBarry Smith } 11853a40ed3dSBarry Smith PetscFunctionReturn(0); 1186905e6a2fSBarry Smith } 1187905e6a2fSBarry Smith 11885615d1e5SSatish Balay #undef __FUNC__ 11895615d1e5SSatish Balay #define __FUNC__ "MatCopy_SeqDense" 1190cb5b572fSBarry Smith int MatCopy_SeqDense(Mat A, Mat B,MatStructure str) 11914b0e389bSBarry Smith { 11924b0e389bSBarry Smith Mat_SeqDense *a = (Mat_SeqDense *) A->data, *b = (Mat_SeqDense *)B->data; 11933a40ed3dSBarry Smith int ierr; 11943a40ed3dSBarry Smith 11953a40ed3dSBarry Smith PetscFunctionBegin; 11963a40ed3dSBarry Smith if (B->type != MATSEQDENSE) { 1197cb5b572fSBarry Smith ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr); 11983a40ed3dSBarry Smith PetscFunctionReturn(0); 11993a40ed3dSBarry Smith } 1200e3372554SBarry Smith if (a->m != b->m || a->n != b->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"size(B) != size(A)"); 1201549d3d68SSatish Balay ierr = PetscMemcpy(b->v,a->v,a->m*a->n*sizeof(Scalar));CHKERRQ(ierr); 12023a40ed3dSBarry Smith PetscFunctionReturn(0); 12034b0e389bSBarry Smith } 12044b0e389bSBarry Smith 1205289bc588SBarry Smith /* -------------------------------------------------------------------*/ 1206a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqDense, 1207905e6a2fSBarry Smith MatGetRow_SeqDense, 1208905e6a2fSBarry Smith MatRestoreRow_SeqDense, 1209905e6a2fSBarry Smith MatMult_SeqDense, 1210905e6a2fSBarry Smith MatMultAdd_SeqDense, 1211905e6a2fSBarry Smith MatMultTrans_SeqDense, 1212905e6a2fSBarry Smith MatMultTransAdd_SeqDense, 1213905e6a2fSBarry Smith MatSolve_SeqDense, 1214905e6a2fSBarry Smith MatSolveAdd_SeqDense, 1215905e6a2fSBarry Smith MatSolveTrans_SeqDense, 1216905e6a2fSBarry Smith MatSolveTransAdd_SeqDense, 1217905e6a2fSBarry Smith MatLUFactor_SeqDense, 1218905e6a2fSBarry Smith MatCholeskyFactor_SeqDense, 1219ec8511deSBarry Smith MatRelax_SeqDense, 1220ec8511deSBarry Smith MatTranspose_SeqDense, 1221905e6a2fSBarry Smith MatGetInfo_SeqDense, 1222905e6a2fSBarry Smith MatEqual_SeqDense, 1223905e6a2fSBarry Smith MatGetDiagonal_SeqDense, 1224905e6a2fSBarry Smith MatDiagonalScale_SeqDense, 1225905e6a2fSBarry Smith MatNorm_SeqDense, 1226905e6a2fSBarry Smith 0, 1227905e6a2fSBarry Smith 0, 1228905e6a2fSBarry Smith 0, 1229905e6a2fSBarry Smith MatSetOption_SeqDense, 1230905e6a2fSBarry Smith MatZeroEntries_SeqDense, 1231905e6a2fSBarry Smith MatZeroRows_SeqDense, 1232905e6a2fSBarry Smith MatLUFactorSymbolic_SeqDense, 1233905e6a2fSBarry Smith MatLUFactorNumeric_SeqDense, 1234905e6a2fSBarry Smith MatCholeskyFactorSymbolic_SeqDense, 1235905e6a2fSBarry Smith MatCholeskyFactorNumeric_SeqDense, 1236905e6a2fSBarry Smith MatGetSize_SeqDense, 1237905e6a2fSBarry Smith MatGetSize_SeqDense, 1238905e6a2fSBarry Smith MatGetOwnershipRange_SeqDense, 1239905e6a2fSBarry Smith 0, 1240905e6a2fSBarry Smith 0, 1241905e6a2fSBarry Smith MatGetArray_SeqDense, 1242905e6a2fSBarry Smith MatRestoreArray_SeqDense, 12435609ef8eSBarry Smith MatDuplicate_SeqDense, 1244a5ae1ecdSBarry Smith 0, 1245a5ae1ecdSBarry Smith 0, 1246a5ae1ecdSBarry Smith 0, 1247a5ae1ecdSBarry Smith 0, 1248a5ae1ecdSBarry Smith MatAXPY_SeqDense, 1249a5ae1ecdSBarry Smith MatGetSubMatrices_SeqDense, 1250a5ae1ecdSBarry Smith 0, 12514b0e389bSBarry Smith MatGetValues_SeqDense, 1252a5ae1ecdSBarry Smith MatCopy_SeqDense, 1253a5ae1ecdSBarry Smith 0, 1254a5ae1ecdSBarry Smith MatScale_SeqDense, 1255a5ae1ecdSBarry Smith 0, 1256a5ae1ecdSBarry Smith 0, 1257a5ae1ecdSBarry Smith 0, 1258a5ae1ecdSBarry Smith MatGetBlockSize_SeqDense, 1259a5ae1ecdSBarry Smith 0, 1260a5ae1ecdSBarry Smith 0, 1261a5ae1ecdSBarry Smith 0, 1262a5ae1ecdSBarry Smith 0, 1263a5ae1ecdSBarry Smith 0, 1264a5ae1ecdSBarry Smith 0, 1265a5ae1ecdSBarry Smith 0, 1266a5ae1ecdSBarry Smith 0, 1267a5ae1ecdSBarry Smith 0, 1268a5ae1ecdSBarry Smith 0, 1269a5ae1ecdSBarry Smith 0, 1270a5ae1ecdSBarry Smith 0, 1271a5ae1ecdSBarry Smith MatGetMaps_Petsc}; 127290ace30eSBarry Smith 12735615d1e5SSatish Balay #undef __FUNC__ 12745615d1e5SSatish Balay #define __FUNC__ "MatCreateSeqDense" 12754b828684SBarry Smith /*@C 1276fafbff53SBarry Smith MatCreateSeqDense - Creates a sequential dense matrix that 1277d65003e9SLois Curfman McInnes is stored in column major order (the usual Fortran 77 manner). Many 1278d65003e9SLois Curfman McInnes of the matrix operations use the BLAS and LAPACK routines. 1279289bc588SBarry Smith 1280db81eaa0SLois Curfman McInnes Collective on MPI_Comm 1281db81eaa0SLois Curfman McInnes 128220563c6bSBarry Smith Input Parameters: 1283db81eaa0SLois Curfman McInnes + comm - MPI communicator, set to PETSC_COMM_SELF 12840c775827SLois Curfman McInnes . m - number of rows 128518f449edSLois Curfman McInnes . n - number of columns 1286db81eaa0SLois Curfman McInnes - data - optional location of matrix data. Set data=PETSC_NULL for PETSc 1287dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 128820563c6bSBarry Smith 128920563c6bSBarry Smith Output Parameter: 129044cd7ae7SLois Curfman McInnes . A - the matrix 129120563c6bSBarry Smith 1292b259b22eSLois Curfman McInnes Notes: 129318f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 129418f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 1295b4fd4287SBarry Smith set data=PETSC_NULL. 129618f449edSLois Curfman McInnes 1297027ccd11SLois Curfman McInnes Level: intermediate 1298027ccd11SLois Curfman McInnes 1299dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS 1300d65003e9SLois Curfman McInnes 1301db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues() 130220563c6bSBarry Smith @*/ 130344cd7ae7SLois Curfman McInnes int MatCreateSeqDense(MPI_Comm comm,int m,int n,Scalar *data,Mat *A) 1304289bc588SBarry Smith { 130544cd7ae7SLois Curfman McInnes Mat B; 130644cd7ae7SLois Curfman McInnes Mat_SeqDense *b; 13073b2fbd54SBarry Smith int ierr,flg,size; 13083b2fbd54SBarry Smith 13093a40ed3dSBarry Smith PetscFunctionBegin; 1310d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1311a8c6a408SBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Comm must be of size 1"); 131255659b69SBarry Smith 131344cd7ae7SLois Curfman McInnes *A = 0; 13143f1db9ecSBarry Smith PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQDENSE,"Mat",comm,MatDestroy,MatView); 131544cd7ae7SLois Curfman McInnes PLogObjectCreate(B); 131644cd7ae7SLois Curfman McInnes b = (Mat_SeqDense *) PetscMalloc(sizeof(Mat_SeqDense));CHKPTRQ(b); 1317549d3d68SSatish Balay ierr = PetscMemzero(b,sizeof(Mat_SeqDense));CHKERRQ(ierr); 1318549d3d68SSatish Balay ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1319c6e7dd08SBarry Smith B->ops->destroy = MatDestroy_SeqDense; 1320c6e7dd08SBarry Smith B->ops->view = MatView_SeqDense; 132144cd7ae7SLois Curfman McInnes B->factor = 0; 132290f02eecSBarry Smith B->mapping = 0; 1323f09e8eb9SSatish Balay PLogObjectMemory(B,sizeof(struct _p_Mat)); 132444cd7ae7SLois Curfman McInnes B->data = (void *) b; 132518f449edSLois Curfman McInnes 132644cd7ae7SLois Curfman McInnes b->m = m; B->m = m; B->M = m; 132744cd7ae7SLois Curfman McInnes b->n = n; B->n = n; B->N = n; 1328a5ae1ecdSBarry Smith 1329488ecbafSBarry Smith ierr = MapCreateMPI(comm,m,m,&B->rmap);CHKERRQ(ierr); 1330488ecbafSBarry Smith ierr = MapCreateMPI(comm,n,n,&B->cmap);CHKERRQ(ierr); 1331a5ae1ecdSBarry Smith 133244cd7ae7SLois Curfman McInnes b->pivots = 0; 133344cd7ae7SLois Curfman McInnes b->roworiented = 1; 1334b4fd4287SBarry Smith if (data == PETSC_NULL) { 133544cd7ae7SLois Curfman McInnes b->v = (Scalar*) PetscMalloc((m*n+1)*sizeof(Scalar));CHKPTRQ(b->v); 1336549d3d68SSatish Balay ierr = PetscMemzero(b->v,m*n*sizeof(Scalar));CHKERRQ(ierr); 133744cd7ae7SLois Curfman McInnes b->user_alloc = 0; 133844cd7ae7SLois Curfman McInnes PLogObjectMemory(B,n*m*sizeof(Scalar)); 1339273fdc76SBarry Smith } else { /* user-allocated storage */ 134044cd7ae7SLois Curfman McInnes b->v = data; 134144cd7ae7SLois Curfman McInnes b->user_alloc = 1; 13422dd2b441SLois Curfman McInnes } 134325cdf11fSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help",&flg);CHKERRQ(ierr); 134444cd7ae7SLois Curfman McInnes if (flg) {ierr = MatPrintHelp(B);CHKERRQ(ierr);} 13454e220ebcSLois Curfman McInnes 134644cd7ae7SLois Curfman McInnes *A = B; 13473a40ed3dSBarry Smith PetscFunctionReturn(0); 1348289bc588SBarry Smith } 1349289bc588SBarry Smith 13503b2fbd54SBarry Smith 13513b2fbd54SBarry Smith 13523b2fbd54SBarry Smith 13533b2fbd54SBarry Smith 13543b2fbd54SBarry Smith 1355