1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*c6e7dd08SBarry Smith static char vcid[] = "$Id: dense.c,v 1.139 1998/03/16 18:40:12 balay 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" 9b16a3bb1SBarry Smith #include "pinclude/plapack.h" 10b16a3bb1SBarry Smith #include "pinclude/pviewer.h" 11289bc588SBarry Smith 125615d1e5SSatish Balay #undef __FUNC__ 135615d1e5SSatish Balay #define __FUNC__ "MatAXPY_SeqDense" 141987afe7SBarry Smith int MatAXPY_SeqDense(Scalar *alpha,Mat X,Mat Y) 151987afe7SBarry Smith { 161987afe7SBarry Smith Mat_SeqDense *x = (Mat_SeqDense*) X->data,*y = (Mat_SeqDense*) Y->data; 171987afe7SBarry Smith int N = x->m*x->n, one = 1; 183a40ed3dSBarry Smith 193a40ed3dSBarry Smith PetscFunctionBegin; 201987afe7SBarry Smith BLaxpy_( &N, alpha, x->v, &one, y->v, &one ); 211987afe7SBarry Smith PLogFlops(2*N-1); 223a40ed3dSBarry Smith PetscFunctionReturn(0); 231987afe7SBarry Smith } 241987afe7SBarry Smith 255615d1e5SSatish Balay #undef __FUNC__ 265615d1e5SSatish Balay #define __FUNC__ "MatGetInfo_SeqDense" 278f6be9afSLois Curfman McInnes int MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info) 28289bc588SBarry Smith { 294e220ebcSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense *) A->data; 304e220ebcSLois Curfman McInnes int i,N = a->m*a->n,count = 0; 314e220ebcSLois Curfman McInnes Scalar *v = a->v; 323a40ed3dSBarry Smith 333a40ed3dSBarry Smith PetscFunctionBegin; 34289bc588SBarry Smith for ( i=0; i<N; i++ ) {if (*v != 0.0) count++; v++;} 354e220ebcSLois Curfman McInnes 364e220ebcSLois Curfman McInnes info->rows_global = (double)a->m; 374e220ebcSLois Curfman McInnes info->columns_global = (double)a->n; 384e220ebcSLois Curfman McInnes info->rows_local = (double)a->m; 394e220ebcSLois Curfman McInnes info->columns_local = (double)a->n; 404e220ebcSLois Curfman McInnes info->block_size = 1.0; 414e220ebcSLois Curfman McInnes info->nz_allocated = (double)N; 424e220ebcSLois Curfman McInnes info->nz_used = (double)count; 434e220ebcSLois Curfman McInnes info->nz_unneeded = (double)(N-count); 444e220ebcSLois Curfman McInnes info->assemblies = (double)A->num_ass; 454e220ebcSLois Curfman McInnes info->mallocs = 0; 464e220ebcSLois Curfman McInnes info->memory = A->mem; 474e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; 484e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 494e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 504e220ebcSLois Curfman McInnes 513a40ed3dSBarry Smith PetscFunctionReturn(0); 52289bc588SBarry Smith } 53289bc588SBarry Smith 545615d1e5SSatish Balay #undef __FUNC__ 555615d1e5SSatish Balay #define __FUNC__ "MatScale_SeqDense" 568f6be9afSLois Curfman McInnes int MatScale_SeqDense(Scalar *alpha,Mat inA) 5780cd9d93SLois Curfman McInnes { 5880cd9d93SLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense *) inA->data; 5980cd9d93SLois Curfman McInnes int one = 1, nz; 6080cd9d93SLois Curfman McInnes 613a40ed3dSBarry Smith PetscFunctionBegin; 6280cd9d93SLois Curfman McInnes nz = a->m*a->n; 6380cd9d93SLois Curfman McInnes BLscal_( &nz, alpha, a->v, &one ); 6480cd9d93SLois Curfman McInnes PLogFlops(nz); 653a40ed3dSBarry Smith PetscFunctionReturn(0); 6680cd9d93SLois Curfman McInnes } 6780cd9d93SLois Curfman McInnes 68289bc588SBarry Smith /* ---------------------------------------------------------------*/ 69289bc588SBarry Smith /* COMMENT: I have chosen to hide column permutation in the pivots, 70289bc588SBarry Smith rather than put it in the Mat->col slot.*/ 715615d1e5SSatish Balay #undef __FUNC__ 725615d1e5SSatish Balay #define __FUNC__ "MatLUFactor_SeqDense" 738f6be9afSLois Curfman McInnes int MatLUFactor_SeqDense(Mat A,IS row,IS col,double f) 74289bc588SBarry Smith { 75c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 766abc6512SBarry Smith int info; 7755659b69SBarry Smith 783a40ed3dSBarry Smith PetscFunctionBegin; 79289bc588SBarry Smith if (!mat->pivots) { 8045d48df9SBarry Smith mat->pivots = (int *) PetscMalloc((mat->m+1)*sizeof(int));CHKPTRQ(mat->pivots); 81c0bbcb79SLois Curfman McInnes PLogObjectMemory(A,mat->m*sizeof(int)); 82289bc588SBarry Smith } 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__ 925615d1e5SSatish Balay #define __FUNC__ "MatConvertSameType_SeqDense" 938f6be9afSLois Curfman McInnes int MatConvertSameType_SeqDense(Mat A,Mat *newmat,int cpvalues) 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; 10202cad45dSBarry Smith if (cpvalues == COPY_VALUES) { 10302cad45dSBarry Smith PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar)); 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; 1173a40ed3dSBarry Smith ierr = MatConvertSameType_SeqDense(A,fact,PETSC_FALSE);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 */ 13002cad45dSBarry Smith PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar)); 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; 1636abc6512SBarry Smith int info; 16455659b69SBarry Smith 1653a40ed3dSBarry Smith PetscFunctionBegin; 166752f5567SLois Curfman McInnes if (mat->pivots) { 1670452661fSBarry Smith PetscFree(mat->pivots); 168c0bbcb79SLois Curfman McInnes PLogObjectMemory(A,-mat->m*sizeof(int)); 169752f5567SLois Curfman McInnes mat->pivots = 0; 170752f5567SLois Curfman McInnes } 171289bc588SBarry Smith LApotrf_("L",&mat->n,mat->v,&mat->m,&info); 172e3372554SBarry Smith if (info) SETERRQ(PETSC_ERR_MAT_CH_ZRPVT,0,"Bad factorization"); 173c0bbcb79SLois Curfman McInnes A->factor = FACTOR_CHOLESKY; 17455659b69SBarry Smith PLogFlops((mat->n*mat->n*mat->n)/3); 1753a40ed3dSBarry Smith PetscFunctionReturn(0); 176289bc588SBarry Smith } 177289bc588SBarry Smith 1785615d1e5SSatish Balay #undef __FUNC__ 1795615d1e5SSatish Balay #define __FUNC__ "MatSolve_SeqDense" 1808f6be9afSLois Curfman McInnes int MatSolve_SeqDense(Mat A,Vec xx,Vec yy) 181289bc588SBarry Smith { 182c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 183a57ad546SLois Curfman McInnes int one = 1, info, ierr; 1846abc6512SBarry Smith Scalar *x, *y; 18567e560aaSBarry Smith 1863a40ed3dSBarry Smith PetscFunctionBegin; 187a57ad546SLois Curfman McInnes ierr = VecGetArray(xx,&x); CHKERRQ(ierr); ierr = VecGetArray(yy,&y); CHKERRQ(ierr); 188416022c9SBarry Smith PetscMemcpy(y,x,mat->m*sizeof(Scalar)); 189c0bbcb79SLois Curfman McInnes if (A->factor == FACTOR_LU) { 19048d91487SBarry Smith LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 191289bc588SBarry Smith } 192c0bbcb79SLois Curfman McInnes else if (A->factor == FACTOR_CHOLESKY){ 19348d91487SBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 194289bc588SBarry Smith } 195a8c6a408SBarry Smith else SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Matrix must be factored to solve"); 196a8c6a408SBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,0,"MBad solve"); 19755659b69SBarry Smith PLogFlops(mat->n*mat->n - mat->n); 1983a40ed3dSBarry Smith PetscFunctionReturn(0); 199289bc588SBarry Smith } 2006ee01492SSatish Balay 2015615d1e5SSatish Balay #undef __FUNC__ 2025615d1e5SSatish Balay #define __FUNC__ "MatSolveTrans_SeqDense" 2038f6be9afSLois Curfman McInnes int MatSolveTrans_SeqDense(Mat A,Vec xx,Vec yy) 204da3a660dSBarry Smith { 205c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 2066abc6512SBarry Smith int one = 1, info; 2076abc6512SBarry Smith Scalar *x, *y; 20867e560aaSBarry Smith 2093a40ed3dSBarry Smith PetscFunctionBegin; 210da3a660dSBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 211416022c9SBarry Smith PetscMemcpy(y,x,mat->m*sizeof(Scalar)); 212752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 213da3a660dSBarry Smith if (mat->pivots) { 21448d91487SBarry Smith LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 215da3a660dSBarry Smith } 216da3a660dSBarry Smith else { 21748d91487SBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 218da3a660dSBarry Smith } 219a8c6a408SBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,0,"Bad solve"); 22055659b69SBarry Smith PLogFlops(mat->n*mat->n - mat->n); 2213a40ed3dSBarry Smith PetscFunctionReturn(0); 222da3a660dSBarry Smith } 2236ee01492SSatish Balay 2245615d1e5SSatish Balay #undef __FUNC__ 2255615d1e5SSatish Balay #define __FUNC__ "MatSolveAdd_SeqDense" 2268f6be9afSLois Curfman McInnes int MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 227da3a660dSBarry Smith { 228c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 2296abc6512SBarry Smith int one = 1, info,ierr; 2306abc6512SBarry Smith Scalar *x, *y, sone = 1.0; 231da3a660dSBarry Smith Vec tmp = 0; 23267e560aaSBarry Smith 2333a40ed3dSBarry Smith PetscFunctionBegin; 234da3a660dSBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 235da3a660dSBarry Smith if (yy == zz) { 23678b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr); 237c0bbcb79SLois Curfman McInnes PLogObjectParent(A,tmp); 23878b31e54SBarry Smith ierr = VecCopy(yy,tmp); CHKERRQ(ierr); 239da3a660dSBarry Smith } 240416022c9SBarry Smith PetscMemcpy(y,x,mat->m*sizeof(Scalar)); 241752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 242da3a660dSBarry Smith if (mat->pivots) { 24348d91487SBarry Smith LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 244a8c6a408SBarry Smith } else { 24548d91487SBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 246da3a660dSBarry Smith } 247a8c6a408SBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,0,"Bad solve"); 248a8c6a408SBarry Smith if (tmp) {ierr = VecAXPY(&sone,tmp,yy); CHKERRQ(ierr); ierr = VecDestroy(tmp);CHKERRQ(ierr);} 249a8c6a408SBarry Smith else {ierr = VecAXPY(&sone,zz,yy); CHKERRQ(ierr);} 25055659b69SBarry Smith PLogFlops(mat->n*mat->n - mat->n); 2513a40ed3dSBarry Smith PetscFunctionReturn(0); 252da3a660dSBarry Smith } 25367e560aaSBarry Smith 2545615d1e5SSatish Balay #undef __FUNC__ 2555615d1e5SSatish Balay #define __FUNC__ "MatSolveTransAdd_SeqDense" 2568f6be9afSLois Curfman McInnes int MatSolveTransAdd_SeqDense(Mat A,Vec xx,Vec zz, Vec yy) 257da3a660dSBarry Smith { 258c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 2596abc6512SBarry Smith int one = 1, info,ierr; 2606abc6512SBarry Smith Scalar *x, *y, sone = 1.0; 261da3a660dSBarry Smith Vec tmp; 26267e560aaSBarry Smith 2633a40ed3dSBarry Smith PetscFunctionBegin; 264da3a660dSBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 265da3a660dSBarry Smith if (yy == zz) { 26678b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr); 267c0bbcb79SLois Curfman McInnes PLogObjectParent(A,tmp); 26878b31e54SBarry Smith ierr = VecCopy(yy,tmp); CHKERRQ(ierr); 269da3a660dSBarry Smith } 270416022c9SBarry Smith PetscMemcpy(y,x,mat->m*sizeof(Scalar)); 271752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 272da3a660dSBarry Smith if (mat->pivots) { 27348d91487SBarry Smith LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 2743a40ed3dSBarry Smith } else { 27548d91487SBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 276da3a660dSBarry Smith } 277a8c6a408SBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,0,"Bad solve"); 27890f02eecSBarry Smith if (tmp) { 27990f02eecSBarry Smith ierr = VecAXPY(&sone,tmp,yy); CHKERRQ(ierr); 28090f02eecSBarry Smith ierr = VecDestroy(tmp); CHKERRQ(ierr); 2813a40ed3dSBarry Smith } else { 28290f02eecSBarry Smith ierr = VecAXPY(&sone,zz,yy); CHKERRQ(ierr); 28390f02eecSBarry Smith } 28455659b69SBarry Smith PLogFlops(mat->n*mat->n - mat->n); 2853a40ed3dSBarry Smith PetscFunctionReturn(0); 286da3a660dSBarry Smith } 287289bc588SBarry Smith /* ------------------------------------------------------------------*/ 2885615d1e5SSatish Balay #undef __FUNC__ 2895615d1e5SSatish Balay #define __FUNC__ "MatRelax_SeqDense" 2908f6be9afSLois Curfman McInnes int MatRelax_SeqDense(Mat A,Vec bb,double omega,MatSORType flag, 29120e1a0d4SLois Curfman McInnes double shift,int its,Vec xx) 292289bc588SBarry Smith { 293c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 294289bc588SBarry Smith Scalar *x, *b, *v = mat->v, zero = 0.0, xt; 295bc1b551cSSatish Balay int ierr, m = mat->m, i; 296bc1b551cSSatish Balay #if !defined(USE_PETSC_COMPLEX) 297bc1b551cSSatish Balay int o = 1; 298bc1b551cSSatish Balay #endif 299289bc588SBarry Smith 3003a40ed3dSBarry Smith PetscFunctionBegin; 301289bc588SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 3023a40ed3dSBarry Smith /* this is a hack fix, should have another version without the second BLdot */ 303bbb6d6a8SBarry Smith ierr = VecSet(&zero,xx); CHKERRQ(ierr); 304289bc588SBarry Smith } 305289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(bb,&b); 306289bc588SBarry Smith while (its--) { 307289bc588SBarry Smith if (flag & SOR_FORWARD_SWEEP){ 308289bc588SBarry Smith for ( i=0; i<m; i++ ) { 3093a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 310f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 311f1747703SBarry Smith not happy about returning a double complex */ 312f1747703SBarry Smith int _i; 313f1747703SBarry Smith Scalar sum = b[i]; 314f1747703SBarry Smith for ( _i=0; _i<m; _i++ ) { 315f1747703SBarry Smith sum -= conj(v[i+_i*m])*x[_i]; 316f1747703SBarry Smith } 317f1747703SBarry Smith xt = sum; 318f1747703SBarry Smith #else 319289bc588SBarry Smith xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 320f1747703SBarry Smith #endif 32155a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 322289bc588SBarry Smith } 323289bc588SBarry Smith } 324289bc588SBarry Smith if (flag & SOR_BACKWARD_SWEEP) { 325289bc588SBarry Smith for ( i=m-1; i>=0; i-- ) { 3263a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 327f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 328f1747703SBarry Smith not happy about returning a double complex */ 329f1747703SBarry Smith int _i; 330f1747703SBarry Smith Scalar sum = b[i]; 331f1747703SBarry Smith for ( _i=0; _i<m; _i++ ) { 332f1747703SBarry Smith sum -= conj(v[i+_i*m])*x[_i]; 333f1747703SBarry Smith } 334f1747703SBarry Smith xt = sum; 335f1747703SBarry Smith #else 336289bc588SBarry Smith xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 337f1747703SBarry Smith #endif 33855a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 339289bc588SBarry Smith } 340289bc588SBarry Smith } 341289bc588SBarry Smith } 3423a40ed3dSBarry Smith PetscFunctionReturn(0); 343289bc588SBarry Smith } 344289bc588SBarry Smith 345289bc588SBarry Smith /* -----------------------------------------------------------------*/ 3465615d1e5SSatish Balay #undef __FUNC__ 3475615d1e5SSatish Balay #define __FUNC__ "MatMultTrans_SeqDense" 34844cd7ae7SLois Curfman McInnes int MatMultTrans_SeqDense(Mat A,Vec xx,Vec yy) 349289bc588SBarry Smith { 350c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 351289bc588SBarry Smith Scalar *v = mat->v, *x, *y; 352289bc588SBarry Smith int _One=1;Scalar _DOne=1.0, _DZero=0.0; 3533a40ed3dSBarry Smith 3543a40ed3dSBarry Smith PetscFunctionBegin; 355289bc588SBarry Smith VecGetArray(xx,&x), VecGetArray(yy,&y); 35648d91487SBarry Smith LAgemv_("T",&(mat->m),&(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DZero,y,&_One); 35755659b69SBarry Smith PLogFlops(2*mat->m*mat->n - mat->n); 3583a40ed3dSBarry Smith PetscFunctionReturn(0); 359289bc588SBarry Smith } 3606ee01492SSatish Balay 3615615d1e5SSatish Balay #undef __FUNC__ 3625615d1e5SSatish Balay #define __FUNC__ "MatMult_SeqDense" 36344cd7ae7SLois Curfman McInnes int MatMult_SeqDense(Mat A,Vec xx,Vec yy) 364289bc588SBarry Smith { 365c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 366289bc588SBarry Smith Scalar *v = mat->v, *x, *y; 367289bc588SBarry Smith int _One=1;Scalar _DOne=1.0, _DZero=0.0; 3683a40ed3dSBarry Smith 3693a40ed3dSBarry Smith PetscFunctionBegin; 370289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 37148d91487SBarry Smith LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DZero,y,&_One); 37255659b69SBarry Smith PLogFlops(2*mat->m*mat->n - mat->m); 3733a40ed3dSBarry Smith PetscFunctionReturn(0); 374289bc588SBarry Smith } 3756ee01492SSatish Balay 3765615d1e5SSatish Balay #undef __FUNC__ 3775615d1e5SSatish Balay #define __FUNC__ "MatMultAdd_SeqDense" 37844cd7ae7SLois Curfman McInnes int MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 379289bc588SBarry Smith { 380c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 381289bc588SBarry Smith Scalar *v = mat->v, *x, *y, *z; 3826abc6512SBarry Smith int _One=1; Scalar _DOne=1.0; 3833a40ed3dSBarry Smith 3843a40ed3dSBarry Smith PetscFunctionBegin; 385289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z); 386416022c9SBarry Smith if (zz != yy) PetscMemcpy(y,z,mat->m*sizeof(Scalar)); 38748d91487SBarry Smith LAgemv_( "N", &(mat->m), &(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DOne,y,&_One); 38855659b69SBarry Smith PLogFlops(2*mat->m*mat->n); 3893a40ed3dSBarry Smith PetscFunctionReturn(0); 390289bc588SBarry Smith } 3916ee01492SSatish Balay 3925615d1e5SSatish Balay #undef __FUNC__ 3935615d1e5SSatish Balay #define __FUNC__ "MatMultTransAdd_SeqDense" 39444cd7ae7SLois Curfman McInnes int MatMultTransAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 395289bc588SBarry Smith { 396c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 397289bc588SBarry Smith Scalar *v = mat->v, *x, *y, *z; 398289bc588SBarry Smith int _One=1; 3996abc6512SBarry Smith Scalar _DOne=1.0; 4003a40ed3dSBarry Smith 4013a40ed3dSBarry Smith PetscFunctionBegin; 40244cd7ae7SLois Curfman McInnes VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z); 403717eeb74SLois Curfman McInnes if (zz != yy) PetscMemcpy(y,z,mat->n*sizeof(Scalar)); 40448d91487SBarry Smith LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DOne,y,&_One); 40555659b69SBarry Smith PLogFlops(2*mat->m*mat->n); 4063a40ed3dSBarry Smith PetscFunctionReturn(0); 407289bc588SBarry Smith } 408289bc588SBarry Smith 409289bc588SBarry Smith /* -----------------------------------------------------------------*/ 4105615d1e5SSatish Balay #undef __FUNC__ 4115615d1e5SSatish Balay #define __FUNC__ "MatGetRow_SeqDense" 4128f6be9afSLois Curfman McInnes int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals) 413289bc588SBarry Smith { 414c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 415289bc588SBarry Smith Scalar *v; 416289bc588SBarry Smith int i; 41767e560aaSBarry Smith 4183a40ed3dSBarry Smith PetscFunctionBegin; 419289bc588SBarry Smith *ncols = mat->n; 420289bc588SBarry Smith if (cols) { 42145d48df9SBarry Smith *cols = (int *) PetscMalloc((mat->n+1)*sizeof(int)); CHKPTRQ(*cols); 42273c3ce57SBarry Smith for ( i=0; i<mat->n; i++ ) (*cols)[i] = i; 423289bc588SBarry Smith } 424289bc588SBarry Smith if (vals) { 42545d48df9SBarry Smith *vals = (Scalar *) PetscMalloc((mat->n+1)*sizeof(Scalar)); CHKPTRQ(*vals); 426289bc588SBarry Smith v = mat->v + row; 42773c3ce57SBarry Smith for ( i=0; i<mat->n; i++ ) {(*vals)[i] = *v; v += mat->m;} 428289bc588SBarry Smith } 4293a40ed3dSBarry Smith PetscFunctionReturn(0); 430289bc588SBarry Smith } 4316ee01492SSatish Balay 4325615d1e5SSatish Balay #undef __FUNC__ 4335615d1e5SSatish Balay #define __FUNC__ "MatRestoreRow_SeqDense" 4348f6be9afSLois Curfman McInnes int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals) 435289bc588SBarry Smith { 4360452661fSBarry Smith if (cols) { PetscFree(*cols); } 4370452661fSBarry Smith if (vals) { PetscFree(*vals); } 4383a40ed3dSBarry Smith PetscFunctionReturn(0); 439289bc588SBarry Smith } 440289bc588SBarry Smith /* ----------------------------------------------------------------*/ 4415615d1e5SSatish Balay #undef __FUNC__ 4425615d1e5SSatish Balay #define __FUNC__ "MatSetValues_SeqDense" 4438f6be9afSLois Curfman McInnes int MatSetValues_SeqDense(Mat A,int m,int *indexm,int n, 444e8d4e0b9SBarry Smith int *indexn,Scalar *v,InsertMode addv) 445289bc588SBarry Smith { 446c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 447289bc588SBarry Smith int i,j; 448d6dfbf8fSBarry Smith 4493a40ed3dSBarry Smith PetscFunctionBegin; 450289bc588SBarry Smith if (!mat->roworiented) { 451dbb450caSBarry Smith if (addv == INSERT_VALUES) { 452289bc588SBarry Smith for ( j=0; j<n; j++ ) { 4533a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g) 454a8c6a408SBarry Smith if (indexn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column"); 455a8c6a408SBarry Smith if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 45658804f6dSBarry Smith #endif 457289bc588SBarry Smith for ( i=0; i<m; i++ ) { 4583a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g) 459a8c6a408SBarry Smith if (indexm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row"); 460a8c6a408SBarry Smith if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 46158804f6dSBarry Smith #endif 462289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 463289bc588SBarry Smith } 464289bc588SBarry Smith } 4653a40ed3dSBarry Smith } else { 466289bc588SBarry Smith for ( j=0; j<n; j++ ) { 4673a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g) 468a8c6a408SBarry Smith if (indexn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column"); 469a8c6a408SBarry Smith if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 47058804f6dSBarry Smith #endif 471289bc588SBarry Smith for ( i=0; i<m; i++ ) { 4723a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g) 473a8c6a408SBarry Smith if (indexm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row"); 474a8c6a408SBarry Smith if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 47558804f6dSBarry Smith #endif 476289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 477289bc588SBarry Smith } 478289bc588SBarry Smith } 479289bc588SBarry Smith } 4803a40ed3dSBarry Smith } else { 481dbb450caSBarry Smith if (addv == INSERT_VALUES) { 482e8d4e0b9SBarry Smith for ( i=0; i<m; i++ ) { 4833a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g) 484a8c6a408SBarry Smith if (indexm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row"); 485a8c6a408SBarry Smith if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 48658804f6dSBarry Smith #endif 487e8d4e0b9SBarry Smith for ( j=0; j<n; j++ ) { 4883a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g) 489a8c6a408SBarry Smith if (indexn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column"); 490a8c6a408SBarry Smith if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 49158804f6dSBarry Smith #endif 492e8d4e0b9SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 493e8d4e0b9SBarry Smith } 494e8d4e0b9SBarry Smith } 4953a40ed3dSBarry Smith } else { 496289bc588SBarry Smith for ( i=0; i<m; i++ ) { 4973a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g) 498a8c6a408SBarry Smith if (indexm[i] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative row"); 499a8c6a408SBarry Smith if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 50058804f6dSBarry Smith #endif 501289bc588SBarry Smith for ( j=0; j<n; j++ ) { 5023a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g) 503a8c6a408SBarry Smith if (indexn[j] < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative column"); 504a8c6a408SBarry Smith if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 50558804f6dSBarry Smith #endif 506289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 507289bc588SBarry Smith } 508289bc588SBarry Smith } 509289bc588SBarry Smith } 510e8d4e0b9SBarry Smith } 5113a40ed3dSBarry Smith PetscFunctionReturn(0); 512289bc588SBarry Smith } 513e8d4e0b9SBarry Smith 5145615d1e5SSatish Balay #undef __FUNC__ 5155615d1e5SSatish Balay #define __FUNC__ "MatGetValues_SeqDense" 5168f6be9afSLois Curfman McInnes int MatGetValues_SeqDense(Mat A,int m,int *indexm,int n,int *indexn,Scalar *v) 517ae80bb75SLois Curfman McInnes { 518ae80bb75SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 519ae80bb75SLois Curfman McInnes int i, j; 520ae80bb75SLois Curfman McInnes Scalar *vpt = v; 521ae80bb75SLois Curfman McInnes 5223a40ed3dSBarry Smith PetscFunctionBegin; 523ae80bb75SLois Curfman McInnes /* row-oriented output */ 524ae80bb75SLois Curfman McInnes for ( i=0; i<m; i++ ) { 525ae80bb75SLois Curfman McInnes for ( j=0; j<n; j++ ) { 526ae80bb75SLois Curfman McInnes *vpt++ = mat->v[indexn[j]*mat->m + indexm[i]]; 527ae80bb75SLois Curfman McInnes } 528ae80bb75SLois Curfman McInnes } 5293a40ed3dSBarry Smith PetscFunctionReturn(0); 530ae80bb75SLois Curfman McInnes } 531ae80bb75SLois Curfman McInnes 532289bc588SBarry Smith /* -----------------------------------------------------------------*/ 533289bc588SBarry Smith 53477c4ece6SBarry Smith #include "sys.h" 53588e32edaSLois Curfman McInnes 5365615d1e5SSatish Balay #undef __FUNC__ 5375615d1e5SSatish Balay #define __FUNC__ "MatLoad_SeqDense" 53819bcc07fSBarry Smith int MatLoad_SeqDense(Viewer viewer,MatType type,Mat *A) 53988e32edaSLois Curfman McInnes { 54088e32edaSLois Curfman McInnes Mat_SeqDense *a; 54188e32edaSLois Curfman McInnes Mat B; 54288e32edaSLois Curfman McInnes int *scols, i, j, nz, ierr, fd, header[4], size; 54388e32edaSLois Curfman McInnes int *rowlengths = 0, M, N, *cols; 54488e32edaSLois Curfman McInnes Scalar *vals, *svals, *v; 54519bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 54688e32edaSLois Curfman McInnes 5473a40ed3dSBarry Smith PetscFunctionBegin; 54888e32edaSLois Curfman McInnes MPI_Comm_size(comm,&size); 549a8c6a408SBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"view must have one processor"); 55090ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 5510752156aSBarry Smith ierr = PetscBinaryRead(fd,header,4,PETSC_INT); CHKERRQ(ierr); 552a8c6a408SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"Not matrix object"); 55388e32edaSLois Curfman McInnes M = header[1]; N = header[2]; nz = header[3]; 55488e32edaSLois Curfman McInnes 55590ace30eSBarry Smith if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */ 55690ace30eSBarry Smith ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr); 55790ace30eSBarry Smith B = *A; 55890ace30eSBarry Smith a = (Mat_SeqDense *) B->data; 55990ace30eSBarry Smith 56090ace30eSBarry Smith /* read in nonzero values */ 5610752156aSBarry Smith ierr = PetscBinaryRead(fd,a->v,M*N,PETSC_SCALAR); CHKERRQ(ierr); 56290ace30eSBarry Smith 5636d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 5646d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 565d64ed03dSBarry Smith /* ierr = MatTranspose(B,PETSC_NULL); CHKERRQ(ierr); */ 56690ace30eSBarry Smith } else { 56788e32edaSLois Curfman McInnes /* read row lengths */ 56845d48df9SBarry Smith rowlengths = (int*) PetscMalloc( (M+1)*sizeof(int) ); CHKPTRQ(rowlengths); 5690752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT); CHKERRQ(ierr); 57088e32edaSLois Curfman McInnes 57188e32edaSLois Curfman McInnes /* create our matrix */ 572b4fd4287SBarry Smith ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr); 57388e32edaSLois Curfman McInnes B = *A; 57488e32edaSLois Curfman McInnes a = (Mat_SeqDense *) B->data; 57588e32edaSLois Curfman McInnes v = a->v; 57688e32edaSLois Curfman McInnes 57788e32edaSLois Curfman McInnes /* read column indices and nonzeros */ 57845d48df9SBarry Smith cols = scols = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(cols); 5790752156aSBarry Smith ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT); CHKERRQ(ierr); 58045d48df9SBarry Smith vals = svals = (Scalar *) PetscMalloc( (nz+1)*sizeof(Scalar) ); CHKPTRQ(vals); 5810752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR); CHKERRQ(ierr); 58288e32edaSLois Curfman McInnes 58388e32edaSLois Curfman McInnes /* insert into matrix */ 58488e32edaSLois Curfman McInnes for ( i=0; i<M; i++ ) { 58588e32edaSLois Curfman McInnes for ( j=0; j<rowlengths[i]; j++ ) v[i+M*scols[j]] = svals[j]; 58688e32edaSLois Curfman McInnes svals += rowlengths[i]; scols += rowlengths[i]; 58788e32edaSLois Curfman McInnes } 5880452661fSBarry Smith PetscFree(vals); PetscFree(cols); PetscFree(rowlengths); 58988e32edaSLois Curfman McInnes 5906d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 5916d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 59290ace30eSBarry Smith } 5933a40ed3dSBarry Smith PetscFunctionReturn(0); 59488e32edaSLois Curfman McInnes } 59588e32edaSLois Curfman McInnes 596932b0c3eSLois Curfman McInnes #include "pinclude/pviewer.h" 59777c4ece6SBarry Smith #include "sys.h" 598932b0c3eSLois Curfman McInnes 5995615d1e5SSatish Balay #undef __FUNC__ 6005615d1e5SSatish Balay #define __FUNC__ "MatView_SeqDense_ASCII" 601932b0c3eSLois Curfman McInnes static int MatView_SeqDense_ASCII(Mat A,Viewer viewer) 602289bc588SBarry Smith { 603932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense *) A->data; 604932b0c3eSLois Curfman McInnes int ierr, i, j, format; 605d636dbe3SBarry Smith FILE *fd; 606932b0c3eSLois Curfman McInnes char *outputname; 607932b0c3eSLois Curfman McInnes Scalar *v; 608932b0c3eSLois Curfman McInnes 6093a40ed3dSBarry Smith PetscFunctionBegin; 61090ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 611932b0c3eSLois Curfman McInnes ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr); 61290ace30eSBarry Smith ierr = ViewerGetFormat(viewer,&format); 613639f9d9dSBarry Smith if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) { 6143a40ed3dSBarry Smith PetscFunctionReturn(0); /* do nothing for now */ 615f72585f2SLois Curfman McInnes } 61602cad45dSBarry Smith else if (format == VIEWER_FORMAT_ASCII_COMMON) { 61780cd9d93SLois Curfman McInnes for ( i=0; i<a->m; i++ ) { 61844cd7ae7SLois Curfman McInnes v = a->v + i; 61980cd9d93SLois Curfman McInnes fprintf(fd,"row %d:",i); 62080cd9d93SLois Curfman McInnes for ( j=0; j<a->n; j++ ) { 6213a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 622d64ed03dSBarry Smith if (real(*v) != 0.0 && imag(*v) != 0.0) fprintf(fd," %d %g + %g i",j,real(*v),imag(*v)); 62380cd9d93SLois Curfman McInnes else if (real(*v)) fprintf(fd," %d %g ",j,real(*v)); 62480cd9d93SLois Curfman McInnes #else 62580cd9d93SLois Curfman McInnes if (*v) fprintf(fd," %d %g ",j,*v); 62680cd9d93SLois Curfman McInnes #endif 627d64ed03dSBarry Smith v += a->m; 62880cd9d93SLois Curfman McInnes } 62980cd9d93SLois Curfman McInnes fprintf(fd,"\n"); 63080cd9d93SLois Curfman McInnes } 6313a40ed3dSBarry Smith } else { 6323a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 63347989497SBarry Smith int allreal = 1; 63447989497SBarry Smith /* determine if matrix has all real values */ 63547989497SBarry Smith v = a->v; 63647989497SBarry Smith for ( i=0; i<a->m*a->n; i++ ) { 63747989497SBarry Smith if (imag(v[i])) { allreal = 0; break ;} 63847989497SBarry Smith } 63947989497SBarry Smith #endif 640932b0c3eSLois Curfman McInnes for ( i=0; i<a->m; i++ ) { 641932b0c3eSLois Curfman McInnes v = a->v + i; 642932b0c3eSLois Curfman McInnes for ( j=0; j<a->n; j++ ) { 6433a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 64447989497SBarry Smith if (allreal) { 64547989497SBarry Smith fprintf(fd,"%6.4e ",real(*v)); v += a->m; 64647989497SBarry Smith } else { 647932b0c3eSLois Curfman McInnes fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += a->m; 64847989497SBarry Smith } 649289bc588SBarry Smith #else 650932b0c3eSLois Curfman McInnes fprintf(fd,"%6.4e ",*v); v += a->m; 651289bc588SBarry Smith #endif 652289bc588SBarry Smith } 6538e37dc09SBarry Smith fprintf(fd,"\n"); 654289bc588SBarry Smith } 655da3a660dSBarry Smith } 656932b0c3eSLois Curfman McInnes fflush(fd); 6573a40ed3dSBarry Smith PetscFunctionReturn(0); 658289bc588SBarry Smith } 659289bc588SBarry Smith 6605615d1e5SSatish Balay #undef __FUNC__ 6615615d1e5SSatish Balay #define __FUNC__ "MatView_SeqDense_Binary" 662932b0c3eSLois Curfman McInnes static int MatView_SeqDense_Binary(Mat A,Viewer viewer) 663932b0c3eSLois Curfman McInnes { 664932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense *) A->data; 665932b0c3eSLois Curfman McInnes int ict, j, n = a->n, m = a->m, i, fd, *col_lens, ierr, nz = m*n; 66690ace30eSBarry Smith int format; 66790ace30eSBarry Smith Scalar *v, *anonz,*vals; 668932b0c3eSLois Curfman McInnes 6693a40ed3dSBarry Smith PetscFunctionBegin; 67090ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 67190ace30eSBarry Smith 67290ace30eSBarry Smith ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr); 67302cad45dSBarry Smith if (format == VIEWER_FORMAT_BINARY_NATIVE) { 67490ace30eSBarry Smith /* store the matrix as a dense matrix */ 67590ace30eSBarry Smith col_lens = (int *) PetscMalloc( 4*sizeof(int) ); CHKPTRQ(col_lens); 67690ace30eSBarry Smith col_lens[0] = MAT_COOKIE; 67790ace30eSBarry Smith col_lens[1] = m; 67890ace30eSBarry Smith col_lens[2] = n; 67990ace30eSBarry Smith col_lens[3] = MATRIX_BINARY_FORMAT_DENSE; 6800752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,1); CHKERRQ(ierr); 68190ace30eSBarry Smith PetscFree(col_lens); 68290ace30eSBarry Smith 68390ace30eSBarry Smith /* write out matrix, by rows */ 68445d48df9SBarry Smith vals = (Scalar *) PetscMalloc((m*n+1)*sizeof(Scalar)); CHKPTRQ(vals); 68590ace30eSBarry Smith v = a->v; 68690ace30eSBarry Smith for ( i=0; i<m; i++ ) { 68790ace30eSBarry Smith for ( j=0; j<n; j++ ) { 68890ace30eSBarry Smith vals[i + j*m] = *v++; 68990ace30eSBarry Smith } 69090ace30eSBarry Smith } 6910752156aSBarry Smith ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,0); CHKERRQ(ierr); 69290ace30eSBarry Smith PetscFree(vals); 69390ace30eSBarry Smith } else { 6940452661fSBarry Smith col_lens = (int *) PetscMalloc( (4+nz)*sizeof(int) ); CHKPTRQ(col_lens); 695932b0c3eSLois Curfman McInnes col_lens[0] = MAT_COOKIE; 696932b0c3eSLois Curfman McInnes col_lens[1] = m; 697932b0c3eSLois Curfman McInnes col_lens[2] = n; 698932b0c3eSLois Curfman McInnes col_lens[3] = nz; 699932b0c3eSLois Curfman McInnes 700932b0c3eSLois Curfman McInnes /* store lengths of each row and write (including header) to file */ 701932b0c3eSLois Curfman McInnes for ( i=0; i<m; i++ ) col_lens[4+i] = n; 7020752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,1); CHKERRQ(ierr); 703932b0c3eSLois Curfman McInnes 704932b0c3eSLois Curfman McInnes /* Possibly should write in smaller increments, not whole matrix at once? */ 705932b0c3eSLois Curfman McInnes /* store column indices (zero start index) */ 706932b0c3eSLois Curfman McInnes ict = 0; 707932b0c3eSLois Curfman McInnes for ( i=0; i<m; i++ ) { 708932b0c3eSLois Curfman McInnes for ( j=0; j<n; j++ ) col_lens[ict++] = j; 709932b0c3eSLois Curfman McInnes } 7100752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,0); CHKERRQ(ierr); 7110452661fSBarry Smith PetscFree(col_lens); 712932b0c3eSLois Curfman McInnes 713932b0c3eSLois Curfman McInnes /* store nonzero values */ 71445d48df9SBarry Smith anonz = (Scalar *) PetscMalloc((nz+1)*sizeof(Scalar)); CHKPTRQ(anonz); 715932b0c3eSLois Curfman McInnes ict = 0; 716932b0c3eSLois Curfman McInnes for ( i=0; i<m; i++ ) { 717932b0c3eSLois Curfman McInnes v = a->v + i; 718932b0c3eSLois Curfman McInnes for ( j=0; j<n; j++ ) { 719932b0c3eSLois Curfman McInnes anonz[ict++] = *v; v += a->m; 720932b0c3eSLois Curfman McInnes } 721932b0c3eSLois Curfman McInnes } 7220752156aSBarry Smith ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,0); CHKERRQ(ierr); 7230452661fSBarry Smith PetscFree(anonz); 72490ace30eSBarry Smith } 7253a40ed3dSBarry Smith PetscFunctionReturn(0); 726932b0c3eSLois Curfman McInnes } 727932b0c3eSLois Curfman McInnes 7285615d1e5SSatish Balay #undef __FUNC__ 7295615d1e5SSatish Balay #define __FUNC__ "MatView_SeqDense" 730*c6e7dd08SBarry Smith int MatView_SeqDense(Mat A,Viewer viewer) 731932b0c3eSLois Curfman McInnes { 732932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*) A->data; 733bcd2baecSBarry Smith ViewerType vtype; 734bcd2baecSBarry Smith int ierr; 735932b0c3eSLois Curfman McInnes 7363a40ed3dSBarry Smith PetscFunctionBegin; 737bcd2baecSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 738bcd2baecSBarry Smith 739bcd2baecSBarry Smith if (vtype == MATLAB_VIEWER) { 7403a40ed3dSBarry Smith ierr = ViewerMatlabPutArray_Private(viewer,a->m,a->n,a->v); CHKERRQ(ierr); 7413a40ed3dSBarry Smith } else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) { 7423a40ed3dSBarry Smith ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr); 7433a40ed3dSBarry Smith } else if (vtype == BINARY_FILE_VIEWER) { 7443a40ed3dSBarry Smith ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr); 745932b0c3eSLois Curfman McInnes } 7463a40ed3dSBarry Smith PetscFunctionReturn(0); 747932b0c3eSLois Curfman McInnes } 748289bc588SBarry Smith 7495615d1e5SSatish Balay #undef __FUNC__ 7505615d1e5SSatish Balay #define __FUNC__ "MatDestroy_SeqDense" 751*c6e7dd08SBarry Smith int MatDestroy_SeqDense(Mat mat) 752289bc588SBarry Smith { 753ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense *) mat->data; 75490f02eecSBarry Smith int ierr; 75590f02eecSBarry Smith 7563a40ed3dSBarry Smith PetscFunctionBegin; 7573a40ed3dSBarry Smith #if defined(USE_PETSC_LOG) 758*c6e7dd08SBarry Smith PLogObjectState((PetscObject)mat,"Rows %d Cols %d",l->m,l->n); 759a5a9c739SBarry Smith #endif 7600452661fSBarry Smith if (l->pivots) PetscFree(l->pivots); 7613439631bSBarry Smith if (!l->user_alloc) PetscFree(l->v); 7620452661fSBarry Smith PetscFree(l); 76390f02eecSBarry Smith if (mat->mapping) { 76490f02eecSBarry Smith ierr = ISLocalToGlobalMappingDestroy(mat->mapping); CHKERRQ(ierr); 76590f02eecSBarry Smith } 766a5a9c739SBarry Smith PLogObjectDestroy(mat); 7670452661fSBarry Smith PetscHeaderDestroy(mat); 7683a40ed3dSBarry Smith PetscFunctionReturn(0); 769289bc588SBarry Smith } 770289bc588SBarry Smith 7715615d1e5SSatish Balay #undef __FUNC__ 7725615d1e5SSatish Balay #define __FUNC__ "MatTranspose_SeqDense" 7738f6be9afSLois Curfman McInnes int MatTranspose_SeqDense(Mat A,Mat *matout) 774289bc588SBarry Smith { 775c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 776d3e5ee88SLois Curfman McInnes int k, j, m, n; 777d3e5ee88SLois Curfman McInnes Scalar *v, tmp; 77848b35521SBarry Smith 7793a40ed3dSBarry Smith PetscFunctionBegin; 780d3e5ee88SLois Curfman McInnes v = mat->v; m = mat->m; n = mat->n; 7813638b69dSLois Curfman McInnes if (matout == PETSC_NULL) { /* in place transpose */ 782d64ed03dSBarry Smith if (m != n) { /* malloc temp to hold transpose */ 783d64ed03dSBarry Smith Scalar *w = (Scalar *) PetscMalloc((m+1)*(n+1)*sizeof(Scalar));CHKPTRQ(w); 784d64ed03dSBarry Smith for ( j=0; j<m; j++ ) { 785d64ed03dSBarry Smith for ( k=0; k<n; k++ ) { 786d64ed03dSBarry Smith w[k + j*n] = v[j + k*m]; 787d64ed03dSBarry Smith } 788d64ed03dSBarry Smith } 789d64ed03dSBarry Smith PetscMemcpy(v,w,m*n*sizeof(Scalar)); 790d64ed03dSBarry Smith PetscFree(w); 791d64ed03dSBarry Smith } else { 792d3e5ee88SLois Curfman McInnes for ( j=0; j<m; j++ ) { 793289bc588SBarry Smith for ( k=0; k<j; k++ ) { 794d3e5ee88SLois Curfman McInnes tmp = v[j + k*n]; 795d3e5ee88SLois Curfman McInnes v[j + k*n] = v[k + j*n]; 796d3e5ee88SLois Curfman McInnes v[k + j*n] = tmp; 797289bc588SBarry Smith } 798289bc588SBarry Smith } 799d64ed03dSBarry Smith } 8003a40ed3dSBarry Smith } else { /* out-of-place transpose */ 801d3e5ee88SLois Curfman McInnes int ierr; 802d3e5ee88SLois Curfman McInnes Mat tmat; 803ec8511deSBarry Smith Mat_SeqDense *tmatd; 804d3e5ee88SLois Curfman McInnes Scalar *v2; 805b4fd4287SBarry Smith ierr = MatCreateSeqDense(A->comm,mat->n,mat->m,PETSC_NULL,&tmat); CHKERRQ(ierr); 806ec8511deSBarry Smith tmatd = (Mat_SeqDense *) tmat->data; 8070de55854SLois Curfman McInnes v = mat->v; v2 = tmatd->v; 808d3e5ee88SLois Curfman McInnes for ( j=0; j<n; j++ ) { 8090de55854SLois Curfman McInnes for ( k=0; k<m; k++ ) v2[j + k*n] = v[k + j*m]; 810d3e5ee88SLois Curfman McInnes } 8116d4a8577SBarry Smith ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 8126d4a8577SBarry Smith ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 813d3e5ee88SLois Curfman McInnes *matout = tmat; 81448b35521SBarry Smith } 8153a40ed3dSBarry Smith PetscFunctionReturn(0); 816289bc588SBarry Smith } 817289bc588SBarry Smith 8185615d1e5SSatish Balay #undef __FUNC__ 8195615d1e5SSatish Balay #define __FUNC__ "MatEqual_SeqDense" 8208f6be9afSLois Curfman McInnes int MatEqual_SeqDense(Mat A1,Mat A2, PetscTruth *flg) 821289bc588SBarry Smith { 822c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat1 = (Mat_SeqDense *) A1->data; 823c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat2 = (Mat_SeqDense *) A2->data; 824289bc588SBarry Smith int i; 825289bc588SBarry Smith Scalar *v1 = mat1->v, *v2 = mat2->v; 8269ea5d5aeSSatish Balay 8273a40ed3dSBarry Smith PetscFunctionBegin; 8283a40ed3dSBarry Smith if (A2->type != MATSEQDENSE) SETERRQ(PETSC_ERR_SUP,0,"Matrices should be of type MATSEQDENSE"); 8293a40ed3dSBarry Smith if (mat1->m != mat2->m) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 8303a40ed3dSBarry Smith if (mat1->n != mat2->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 831289bc588SBarry Smith for ( i=0; i<mat1->m*mat1->n; i++ ) { 8323a40ed3dSBarry Smith if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 833289bc588SBarry Smith v1++; v2++; 834289bc588SBarry Smith } 83577c4ece6SBarry Smith *flg = PETSC_TRUE; 8363a40ed3dSBarry Smith PetscFunctionReturn(0); 837289bc588SBarry Smith } 838289bc588SBarry Smith 8395615d1e5SSatish Balay #undef __FUNC__ 8405615d1e5SSatish Balay #define __FUNC__ "MatGetDiagonal_SeqDense" 8418f6be9afSLois Curfman McInnes int MatGetDiagonal_SeqDense(Mat A,Vec v) 842289bc588SBarry Smith { 843c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 84444cd7ae7SLois Curfman McInnes int i, n, len; 84544cd7ae7SLois Curfman McInnes Scalar *x, zero = 0.0; 84644cd7ae7SLois Curfman McInnes 8473a40ed3dSBarry Smith PetscFunctionBegin; 84844cd7ae7SLois Curfman McInnes VecSet(&zero,v); 849289bc588SBarry Smith VecGetArray(v,&x); VecGetSize(v,&n); 85044cd7ae7SLois Curfman McInnes len = PetscMin(mat->m,mat->n); 851e3372554SBarry Smith if (n != mat->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Nonconforming mat and vec"); 85244cd7ae7SLois Curfman McInnes for ( i=0; i<len; i++ ) { 853289bc588SBarry Smith x[i] = mat->v[i*mat->m + i]; 854289bc588SBarry Smith } 8553a40ed3dSBarry Smith PetscFunctionReturn(0); 856289bc588SBarry Smith } 857289bc588SBarry Smith 8585615d1e5SSatish Balay #undef __FUNC__ 8595615d1e5SSatish Balay #define __FUNC__ "MatDiagonalScale_SeqDense" 8608f6be9afSLois Curfman McInnes int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr) 861289bc588SBarry Smith { 862c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 863da3a660dSBarry Smith Scalar *l,*r,x,*v; 864da3a660dSBarry Smith int i,j,m = mat->m, n = mat->n; 86555659b69SBarry Smith 8663a40ed3dSBarry Smith PetscFunctionBegin; 86728988994SBarry Smith if (ll) { 868da3a660dSBarry Smith VecGetArray(ll,&l); VecGetSize(ll,&m); 869e3372554SBarry Smith if (m != mat->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Left scaling vec wrong size"); 870da3a660dSBarry Smith for ( i=0; i<m; i++ ) { 871da3a660dSBarry Smith x = l[i]; 872da3a660dSBarry Smith v = mat->v + i; 873da3a660dSBarry Smith for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;} 874da3a660dSBarry Smith } 87544cd7ae7SLois Curfman McInnes PLogFlops(n*m); 876da3a660dSBarry Smith } 87728988994SBarry Smith if (rr) { 878da3a660dSBarry Smith VecGetArray(rr,&r); VecGetSize(rr,&n); 879e3372554SBarry Smith if (n != mat->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Right scaling vec wrong size"); 880da3a660dSBarry Smith for ( i=0; i<n; i++ ) { 881da3a660dSBarry Smith x = r[i]; 882da3a660dSBarry Smith v = mat->v + i*m; 883da3a660dSBarry Smith for ( j=0; j<m; j++ ) { (*v++) *= x;} 884da3a660dSBarry Smith } 88544cd7ae7SLois Curfman McInnes PLogFlops(n*m); 886da3a660dSBarry Smith } 8873a40ed3dSBarry Smith PetscFunctionReturn(0); 888289bc588SBarry Smith } 889289bc588SBarry Smith 8905615d1e5SSatish Balay #undef __FUNC__ 8915615d1e5SSatish Balay #define __FUNC__ "MatNorm_SeqDense" 8928f6be9afSLois Curfman McInnes int MatNorm_SeqDense(Mat A,NormType type,double *norm) 893289bc588SBarry Smith { 894c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 895289bc588SBarry Smith Scalar *v = mat->v; 896289bc588SBarry Smith double sum = 0.0; 897289bc588SBarry Smith int i, j; 89855659b69SBarry Smith 8993a40ed3dSBarry Smith PetscFunctionBegin; 900289bc588SBarry Smith if (type == NORM_FROBENIUS) { 901289bc588SBarry Smith for (i=0; i<mat->n*mat->m; i++ ) { 9023a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 903289bc588SBarry Smith sum += real(conj(*v)*(*v)); v++; 904289bc588SBarry Smith #else 905289bc588SBarry Smith sum += (*v)*(*v); v++; 906289bc588SBarry Smith #endif 907289bc588SBarry Smith } 908289bc588SBarry Smith *norm = sqrt(sum); 90955659b69SBarry Smith PLogFlops(2*mat->n*mat->m); 9103a40ed3dSBarry Smith } else if (type == NORM_1) { 911289bc588SBarry Smith *norm = 0.0; 912289bc588SBarry Smith for ( j=0; j<mat->n; j++ ) { 913289bc588SBarry Smith sum = 0.0; 914289bc588SBarry Smith for ( i=0; i<mat->m; i++ ) { 91533a8263dSBarry Smith sum += PetscAbsScalar(*v); v++; 916289bc588SBarry Smith } 917289bc588SBarry Smith if (sum > *norm) *norm = sum; 918289bc588SBarry Smith } 91955659b69SBarry Smith PLogFlops(mat->n*mat->m); 9203a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { 921289bc588SBarry Smith *norm = 0.0; 922289bc588SBarry Smith for ( j=0; j<mat->m; j++ ) { 923289bc588SBarry Smith v = mat->v + j; 924289bc588SBarry Smith sum = 0.0; 925289bc588SBarry Smith for ( i=0; i<mat->n; i++ ) { 926cddf8d76SBarry Smith sum += PetscAbsScalar(*v); v += mat->m; 927289bc588SBarry Smith } 928289bc588SBarry Smith if (sum > *norm) *norm = sum; 929289bc588SBarry Smith } 93055659b69SBarry Smith PLogFlops(mat->n*mat->m); 9313a40ed3dSBarry Smith } else { 932e3372554SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"No two norm"); 933289bc588SBarry Smith } 9343a40ed3dSBarry Smith PetscFunctionReturn(0); 935289bc588SBarry Smith } 936289bc588SBarry Smith 9375615d1e5SSatish Balay #undef __FUNC__ 9385615d1e5SSatish Balay #define __FUNC__ "MatSetOption_SeqDense" 9398f6be9afSLois Curfman McInnes int MatSetOption_SeqDense(Mat A,MatOption op) 940289bc588SBarry Smith { 941c0bbcb79SLois Curfman McInnes Mat_SeqDense *aij = (Mat_SeqDense *) A->data; 94267e560aaSBarry Smith 9433a40ed3dSBarry Smith PetscFunctionBegin; 9446d4a8577SBarry Smith if (op == MAT_ROW_ORIENTED) aij->roworiented = 1; 9456d4a8577SBarry Smith else if (op == MAT_COLUMN_ORIENTED) aij->roworiented = 0; 9466d4a8577SBarry Smith else if (op == MAT_ROWS_SORTED || 947219d9a1aSLois Curfman McInnes op == MAT_ROWS_UNSORTED || 9486d4a8577SBarry Smith op == MAT_COLUMNS_SORTED || 949219d9a1aSLois Curfman McInnes op == MAT_COLUMNS_UNSORTED || 9506d4a8577SBarry Smith op == MAT_SYMMETRIC || 9516d4a8577SBarry Smith op == MAT_STRUCTURALLY_SYMMETRIC || 9526d4a8577SBarry Smith op == MAT_NO_NEW_NONZERO_LOCATIONS || 9536d4a8577SBarry Smith op == MAT_YES_NEW_NONZERO_LOCATIONS || 954c2653b3dSLois Curfman McInnes op == MAT_NEW_NONZERO_LOCATION_ERROR || 9556d4a8577SBarry Smith op == MAT_NO_NEW_DIAGONALS || 95690f02eecSBarry Smith op == MAT_YES_NEW_DIAGONALS || 957b51ba29fSSatish Balay op == MAT_IGNORE_OFF_PROC_ENTRIES || 958b51ba29fSSatish Balay op == MAT_USE_HASH_TABLE) 959981c4779SBarry Smith PLogInfo(A,"MatSetOption_SeqDense:Option ignored\n"); 9603a40ed3dSBarry Smith else { 9613a40ed3dSBarry Smith SETERRQ(PETSC_ERR_SUP,0,"unknown option"); 9623a40ed3dSBarry Smith } 9633a40ed3dSBarry Smith PetscFunctionReturn(0); 964289bc588SBarry Smith } 965289bc588SBarry Smith 9665615d1e5SSatish Balay #undef __FUNC__ 9675615d1e5SSatish Balay #define __FUNC__ "MatZeroEntries_SeqDense" 9688f6be9afSLois Curfman McInnes int MatZeroEntries_SeqDense(Mat A) 9696f0a148fSBarry Smith { 970ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense *) A->data; 9713a40ed3dSBarry Smith 9723a40ed3dSBarry Smith PetscFunctionBegin; 973cddf8d76SBarry Smith PetscMemzero(l->v,l->m*l->n*sizeof(Scalar)); 9743a40ed3dSBarry Smith PetscFunctionReturn(0); 9756f0a148fSBarry Smith } 9766f0a148fSBarry Smith 9775615d1e5SSatish Balay #undef __FUNC__ 9785615d1e5SSatish Balay #define __FUNC__ "MatGetBlockSize_SeqDense" 9798f6be9afSLois Curfman McInnes int MatGetBlockSize_SeqDense(Mat A,int *bs) 9804e220ebcSLois Curfman McInnes { 9813a40ed3dSBarry Smith PetscFunctionBegin; 9824e220ebcSLois Curfman McInnes *bs = 1; 9833a40ed3dSBarry Smith PetscFunctionReturn(0); 9844e220ebcSLois Curfman McInnes } 9854e220ebcSLois Curfman McInnes 9865615d1e5SSatish Balay #undef __FUNC__ 9875615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqDense" 9888f6be9afSLois Curfman McInnes int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag) 9896f0a148fSBarry Smith { 990ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense *) A->data; 9916abc6512SBarry Smith int n = l->n, i, j,ierr,N, *rows; 9926f0a148fSBarry Smith Scalar *slot; 99355659b69SBarry Smith 9943a40ed3dSBarry Smith PetscFunctionBegin; 99577c4ece6SBarry Smith ierr = ISGetSize(is,&N); CHKERRQ(ierr); 99678b31e54SBarry Smith ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 9976f0a148fSBarry Smith for ( i=0; i<N; i++ ) { 9986f0a148fSBarry Smith slot = l->v + rows[i]; 9996f0a148fSBarry Smith for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;} 10006f0a148fSBarry Smith } 10016f0a148fSBarry Smith if (diag) { 10026f0a148fSBarry Smith for ( i=0; i<N; i++ ) { 10036f0a148fSBarry Smith slot = l->v + (n+1)*rows[i]; 1004260925b8SBarry Smith *slot = *diag; 10056f0a148fSBarry Smith } 10066f0a148fSBarry Smith } 1007260925b8SBarry Smith ISRestoreIndices(is,&rows); 10083a40ed3dSBarry Smith PetscFunctionReturn(0); 10096f0a148fSBarry Smith } 1010557bce09SLois Curfman McInnes 10115615d1e5SSatish Balay #undef __FUNC__ 10125615d1e5SSatish Balay #define __FUNC__ "MatGetSize_SeqDense" 10138f6be9afSLois Curfman McInnes int MatGetSize_SeqDense(Mat A,int *m,int *n) 1014557bce09SLois Curfman McInnes { 1015c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 10163a40ed3dSBarry Smith 10173a40ed3dSBarry Smith PetscFunctionBegin; 1018557bce09SLois Curfman McInnes *m = mat->m; *n = mat->n; 10193a40ed3dSBarry Smith PetscFunctionReturn(0); 1020557bce09SLois Curfman McInnes } 1021557bce09SLois Curfman McInnes 10225615d1e5SSatish Balay #undef __FUNC__ 10235615d1e5SSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqDense" 10248f6be9afSLois Curfman McInnes int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n) 1025d3e5ee88SLois Curfman McInnes { 1026c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 10273a40ed3dSBarry Smith 10283a40ed3dSBarry Smith PetscFunctionBegin; 1029d3e5ee88SLois Curfman McInnes *m = 0; *n = mat->m; 10303a40ed3dSBarry Smith PetscFunctionReturn(0); 1031d3e5ee88SLois Curfman McInnes } 1032d3e5ee88SLois Curfman McInnes 10335615d1e5SSatish Balay #undef __FUNC__ 10345615d1e5SSatish Balay #define __FUNC__ "MatGetArray_SeqDense" 10358f6be9afSLois Curfman McInnes int MatGetArray_SeqDense(Mat A,Scalar **array) 103664e87e97SBarry Smith { 1037c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 10383a40ed3dSBarry Smith 10393a40ed3dSBarry Smith PetscFunctionBegin; 104064e87e97SBarry Smith *array = mat->v; 10413a40ed3dSBarry Smith PetscFunctionReturn(0); 104264e87e97SBarry Smith } 10430754003eSLois Curfman McInnes 10445615d1e5SSatish Balay #undef __FUNC__ 10455615d1e5SSatish Balay #define __FUNC__ "MatRestoreArray_SeqDense" 10468f6be9afSLois Curfman McInnes int MatRestoreArray_SeqDense(Mat A,Scalar **array) 1047ff14e315SSatish Balay { 10483a40ed3dSBarry Smith PetscFunctionBegin; 10493a40ed3dSBarry Smith PetscFunctionReturn(0); 1050ff14e315SSatish Balay } 10510754003eSLois Curfman McInnes 10525615d1e5SSatish Balay #undef __FUNC__ 10535615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrix_SeqDense" 10546a6a5d1dSBarry Smith static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,int cs,MatGetSubMatrixCall scall,Mat *submat) 10550754003eSLois Curfman McInnes { 1056c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 10570754003eSLois Curfman McInnes int nznew, *smap, i, j, ierr, oldcols = mat->n; 1058160018dcSBarry Smith int *irow, *icol, nrows, ncols, *cwork; 10590754003eSLois Curfman McInnes Scalar *vwork, *val; 10600754003eSLois Curfman McInnes Mat newmat; 10610754003eSLois Curfman McInnes 10623a40ed3dSBarry Smith PetscFunctionBegin; 106378b31e54SBarry Smith ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr); 106478b31e54SBarry Smith ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr); 106578b31e54SBarry Smith ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr); 106678b31e54SBarry Smith ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr); 10670754003eSLois Curfman McInnes 10680452661fSBarry Smith smap = (int *) PetscMalloc(oldcols*sizeof(int)); CHKPTRQ(smap); 10690452661fSBarry Smith cwork = (int *) PetscMalloc(ncols*sizeof(int)); CHKPTRQ(cwork); 10700452661fSBarry Smith vwork = (Scalar *) PetscMalloc(ncols*sizeof(Scalar)); CHKPTRQ(vwork); 1071cddf8d76SBarry Smith PetscMemzero((char*)smap,oldcols*sizeof(int)); 10720754003eSLois Curfman McInnes for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1; 10730754003eSLois Curfman McInnes 10740754003eSLois Curfman McInnes /* Create and fill new matrix */ 1075b4fd4287SBarry Smith ierr = MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat); CHKERRQ(ierr); 10760754003eSLois Curfman McInnes for (i=0; i<nrows; i++) { 10770754003eSLois Curfman McInnes nznew = 0; 10780754003eSLois Curfman McInnes val = mat->v + irow[i]; 10790754003eSLois Curfman McInnes for (j=0; j<oldcols; j++) { 10800754003eSLois Curfman McInnes if (smap[j]) { 10810754003eSLois Curfman McInnes cwork[nznew] = smap[j] - 1; 10820754003eSLois Curfman McInnes vwork[nznew++] = val[j * mat->m]; 10830754003eSLois Curfman McInnes } 10840754003eSLois Curfman McInnes } 10853a40ed3dSBarry Smith ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERT_VALUES);CHKERRQ(ierr); 10860754003eSLois Curfman McInnes } 10876d4a8577SBarry Smith ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 10886d4a8577SBarry Smith ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 10890754003eSLois Curfman McInnes 10900754003eSLois Curfman McInnes /* Free work space */ 10910452661fSBarry Smith PetscFree(smap); PetscFree(cwork); PetscFree(vwork); 109278b31e54SBarry Smith ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr); 109378b31e54SBarry Smith ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr); 10940754003eSLois Curfman McInnes *submat = newmat; 10953a40ed3dSBarry Smith PetscFunctionReturn(0); 10960754003eSLois Curfman McInnes } 10970754003eSLois Curfman McInnes 10985615d1e5SSatish Balay #undef __FUNC__ 10995615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrices_SeqDense" 11006a6a5d1dSBarry Smith int MatGetSubMatrices_SeqDense(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall,Mat **B) 1101905e6a2fSBarry Smith { 1102905e6a2fSBarry Smith int ierr,i; 1103905e6a2fSBarry Smith 11043a40ed3dSBarry Smith PetscFunctionBegin; 1105905e6a2fSBarry Smith if (scall == MAT_INITIAL_MATRIX) { 1106905e6a2fSBarry Smith *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B); 1107905e6a2fSBarry Smith } 1108905e6a2fSBarry Smith 1109905e6a2fSBarry Smith for ( i=0; i<n; i++ ) { 11106a6a5d1dSBarry Smith ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr); 1111905e6a2fSBarry Smith } 11123a40ed3dSBarry Smith PetscFunctionReturn(0); 1113905e6a2fSBarry Smith } 1114905e6a2fSBarry Smith 11155615d1e5SSatish Balay #undef __FUNC__ 11165615d1e5SSatish Balay #define __FUNC__ "MatCopy_SeqDense" 11178f6be9afSLois Curfman McInnes int MatCopy_SeqDense(Mat A, Mat B) 11184b0e389bSBarry Smith { 11194b0e389bSBarry Smith Mat_SeqDense *a = (Mat_SeqDense *) A->data, *b = (Mat_SeqDense *)B->data; 11203a40ed3dSBarry Smith int ierr; 11213a40ed3dSBarry Smith 11223a40ed3dSBarry Smith PetscFunctionBegin; 11233a40ed3dSBarry Smith if (B->type != MATSEQDENSE) { 11243a40ed3dSBarry Smith ierr = MatCopy_Basic(A,B);CHKERRQ(ierr); 11253a40ed3dSBarry Smith PetscFunctionReturn(0); 11263a40ed3dSBarry Smith } 1127e3372554SBarry Smith if (a->m != b->m || a->n != b->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"size(B) != size(A)"); 11284b0e389bSBarry Smith PetscMemcpy(b->v,a->v,a->m*a->n*sizeof(Scalar)); 11293a40ed3dSBarry Smith PetscFunctionReturn(0); 11304b0e389bSBarry Smith } 11314b0e389bSBarry Smith 1132289bc588SBarry Smith /* -------------------------------------------------------------------*/ 1133ae80bb75SLois Curfman McInnes static struct _MatOps MatOps = {MatSetValues_SeqDense, 1134905e6a2fSBarry Smith MatGetRow_SeqDense, 1135905e6a2fSBarry Smith MatRestoreRow_SeqDense, 1136905e6a2fSBarry Smith MatMult_SeqDense, 1137905e6a2fSBarry Smith MatMultAdd_SeqDense, 1138905e6a2fSBarry Smith MatMultTrans_SeqDense, 1139905e6a2fSBarry Smith MatMultTransAdd_SeqDense, 1140905e6a2fSBarry Smith MatSolve_SeqDense, 1141905e6a2fSBarry Smith MatSolveAdd_SeqDense, 1142905e6a2fSBarry Smith MatSolveTrans_SeqDense, 1143905e6a2fSBarry Smith MatSolveTransAdd_SeqDense, 1144905e6a2fSBarry Smith MatLUFactor_SeqDense, 1145905e6a2fSBarry Smith MatCholeskyFactor_SeqDense, 1146ec8511deSBarry Smith MatRelax_SeqDense, 1147ec8511deSBarry Smith MatTranspose_SeqDense, 1148905e6a2fSBarry Smith MatGetInfo_SeqDense, 1149905e6a2fSBarry Smith MatEqual_SeqDense, 1150905e6a2fSBarry Smith MatGetDiagonal_SeqDense, 1151905e6a2fSBarry Smith MatDiagonalScale_SeqDense, 1152905e6a2fSBarry Smith MatNorm_SeqDense, 1153905e6a2fSBarry Smith 0, 1154905e6a2fSBarry Smith 0, 1155905e6a2fSBarry Smith 0, 1156905e6a2fSBarry Smith MatSetOption_SeqDense, 1157905e6a2fSBarry Smith MatZeroEntries_SeqDense, 1158905e6a2fSBarry Smith MatZeroRows_SeqDense, 1159905e6a2fSBarry Smith MatLUFactorSymbolic_SeqDense, 1160905e6a2fSBarry Smith MatLUFactorNumeric_SeqDense, 1161905e6a2fSBarry Smith MatCholeskyFactorSymbolic_SeqDense, 1162905e6a2fSBarry Smith MatCholeskyFactorNumeric_SeqDense, 1163905e6a2fSBarry Smith MatGetSize_SeqDense, 1164905e6a2fSBarry Smith MatGetSize_SeqDense, 1165905e6a2fSBarry Smith MatGetOwnershipRange_SeqDense, 1166905e6a2fSBarry Smith 0, 1167905e6a2fSBarry Smith 0, 1168905e6a2fSBarry Smith MatGetArray_SeqDense, 1169905e6a2fSBarry Smith MatRestoreArray_SeqDense, 11703d1612f7SBarry Smith MatConvertSameType_SeqDense,0,0,0,0, 1171905e6a2fSBarry Smith MatAXPY_SeqDense,MatGetSubMatrices_SeqDense,0, 11724b0e389bSBarry Smith MatGetValues_SeqDense, 11734e220ebcSLois Curfman McInnes MatCopy_SeqDense,0,MatScale_SeqDense, 11744e220ebcSLois Curfman McInnes 0,0,0,MatGetBlockSize_SeqDense}; 117590ace30eSBarry Smith 11765615d1e5SSatish Balay #undef __FUNC__ 11775615d1e5SSatish Balay #define __FUNC__ "MatCreateSeqDense" 11784b828684SBarry Smith /*@C 1179fafbff53SBarry Smith MatCreateSeqDense - Creates a sequential dense matrix that 1180d65003e9SLois Curfman McInnes is stored in column major order (the usual Fortran 77 manner). Many 1181d65003e9SLois Curfman McInnes of the matrix operations use the BLAS and LAPACK routines. 1182289bc588SBarry Smith 118320563c6bSBarry Smith Input Parameters: 1184029af93fSBarry Smith . comm - MPI communicator, set to PETSC_COMM_SELF 11850c775827SLois Curfman McInnes . m - number of rows 118618f449edSLois Curfman McInnes . n - number of columns 1187b4fd4287SBarry Smith . data - optional location of matrix data. Set data=PETSC_NULL for PETSc 1188dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 118920563c6bSBarry Smith 119020563c6bSBarry Smith Output Parameter: 119144cd7ae7SLois Curfman McInnes . A - the matrix 119220563c6bSBarry Smith 119318f449edSLois Curfman McInnes Notes: 119418f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 119518f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 1196b4fd4287SBarry Smith set data=PETSC_NULL. 119718f449edSLois Curfman McInnes 1198dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS 1199d65003e9SLois Curfman McInnes 1200dbd7a890SLois Curfman McInnes .seealso: MatCreate(), MatSetValues() 120120563c6bSBarry Smith @*/ 120244cd7ae7SLois Curfman McInnes int MatCreateSeqDense(MPI_Comm comm,int m,int n,Scalar *data,Mat *A) 1203289bc588SBarry Smith { 120444cd7ae7SLois Curfman McInnes Mat B; 120544cd7ae7SLois Curfman McInnes Mat_SeqDense *b; 12063b2fbd54SBarry Smith int ierr,flg,size; 12073b2fbd54SBarry Smith 12083a40ed3dSBarry Smith PetscFunctionBegin; 12093b2fbd54SBarry Smith MPI_Comm_size(comm,&size); 1210a8c6a408SBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Comm must be of size 1"); 121155659b69SBarry Smith 121244cd7ae7SLois Curfman McInnes *A = 0; 1213f830108cSBarry Smith PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQDENSE,comm,MatDestroy,MatView); 121444cd7ae7SLois Curfman McInnes PLogObjectCreate(B); 121544cd7ae7SLois Curfman McInnes b = (Mat_SeqDense *) PetscMalloc(sizeof(Mat_SeqDense)); CHKPTRQ(b); 121644cd7ae7SLois Curfman McInnes PetscMemzero(b,sizeof(Mat_SeqDense)); 1217f830108cSBarry Smith PetscMemcpy(B->ops,&MatOps,sizeof(struct _MatOps)); 1218*c6e7dd08SBarry Smith B->ops->destroy = MatDestroy_SeqDense; 1219*c6e7dd08SBarry Smith B->ops->view = MatView_SeqDense; 122044cd7ae7SLois Curfman McInnes B->factor = 0; 122190f02eecSBarry Smith B->mapping = 0; 1222f09e8eb9SSatish Balay PLogObjectMemory(B,sizeof(struct _p_Mat)); 122344cd7ae7SLois Curfman McInnes B->data = (void *) b; 122418f449edSLois Curfman McInnes 122544cd7ae7SLois Curfman McInnes b->m = m; B->m = m; B->M = m; 122644cd7ae7SLois Curfman McInnes b->n = n; B->n = n; B->N = n; 122744cd7ae7SLois Curfman McInnes b->pivots = 0; 122844cd7ae7SLois Curfman McInnes b->roworiented = 1; 1229b4fd4287SBarry Smith if (data == PETSC_NULL) { 123044cd7ae7SLois Curfman McInnes b->v = (Scalar*) PetscMalloc((m*n+1)*sizeof(Scalar)); CHKPTRQ(b->v); 123144cd7ae7SLois Curfman McInnes PetscMemzero(b->v,m*n*sizeof(Scalar)); 123244cd7ae7SLois Curfman McInnes b->user_alloc = 0; 123344cd7ae7SLois Curfman McInnes PLogObjectMemory(B,n*m*sizeof(Scalar)); 123418f449edSLois Curfman McInnes } 12352dd2b441SLois Curfman McInnes else { /* user-allocated storage */ 123644cd7ae7SLois Curfman McInnes b->v = data; 123744cd7ae7SLois Curfman McInnes b->user_alloc = 1; 12382dd2b441SLois Curfman McInnes } 123925cdf11fSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr); 124044cd7ae7SLois Curfman McInnes if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr);} 12414e220ebcSLois Curfman McInnes 124244cd7ae7SLois Curfman McInnes *A = B; 12433a40ed3dSBarry Smith PetscFunctionReturn(0); 1244289bc588SBarry Smith } 1245289bc588SBarry Smith 12463b2fbd54SBarry Smith 12473b2fbd54SBarry Smith 12483b2fbd54SBarry Smith 12493b2fbd54SBarry Smith 12503b2fbd54SBarry Smith 1251