1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*3a40ed3dSBarry Smith static char vcid[] = "$Id: dense.c,v 1.131 1997/09/26 02:18:55 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" 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; 18*3a40ed3dSBarry Smith 19*3a40ed3dSBarry Smith PetscFunctionBegin; 201987afe7SBarry Smith BLaxpy_( &N, alpha, x->v, &one, y->v, &one ); 211987afe7SBarry Smith PLogFlops(2*N-1); 22*3a40ed3dSBarry 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; 32*3a40ed3dSBarry Smith 33*3a40ed3dSBarry 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 51*3a40ed3dSBarry 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 61*3a40ed3dSBarry Smith PetscFunctionBegin; 6280cd9d93SLois Curfman McInnes nz = a->m*a->n; 6380cd9d93SLois Curfman McInnes BLscal_( &nz, alpha, a->v, &one ); 6480cd9d93SLois Curfman McInnes PLogFlops(nz); 65*3a40ed3dSBarry 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 78*3a40ed3dSBarry 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); 84e3372554SBarry Smith if (info<0) SETERRQ(1,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); 88*3a40ed3dSBarry 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 99*3a40ed3dSBarry 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; 107*3a40ed3dSBarry 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 { 114*3a40ed3dSBarry Smith int ierr; 115*3a40ed3dSBarry Smith 116*3a40ed3dSBarry Smith PetscFunctionBegin; 117*3a40ed3dSBarry Smith ierr = MatConvertSameType_SeqDense(A,fact,PETSC_FALSE);CHKERRQ(ierr); 118*3a40ed3dSBarry 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; 126*3a40ed3dSBarry Smith int ierr; 127*3a40ed3dSBarry Smith 128*3a40ed3dSBarry 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; 132*3a40ed3dSBarry Smith ierr = MatLUFactor(*fact,0,0,1.0);CHKERRQ(ierr); 133*3a40ed3dSBarry 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 { 140*3a40ed3dSBarry Smith int ierr; 141*3a40ed3dSBarry Smith 142*3a40ed3dSBarry Smith PetscFunctionBegin; 143*3a40ed3dSBarry Smith ierr = MatConvert(A,MATSAME,fact);CHKERRQ(ierr); 144*3a40ed3dSBarry 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 { 151*3a40ed3dSBarry Smith int ierr; 152*3a40ed3dSBarry Smith 153*3a40ed3dSBarry Smith PetscFunctionBegin; 154*3a40ed3dSBarry Smith ierr = MatCholeskyFactor(*fact,0,1.0);CHKERRQ(ierr); 155*3a40ed3dSBarry 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 165*3a40ed3dSBarry 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); 175*3a40ed3dSBarry 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 186*3a40ed3dSBarry 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 } 195e3372554SBarry Smith else SETERRQ(1,0,"Matrix must be factored to solve"); 196e3372554SBarry Smith if (info) SETERRQ(1,0,"MBad solve"); 19755659b69SBarry Smith PLogFlops(mat->n*mat->n - mat->n); 198*3a40ed3dSBarry 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 209*3a40ed3dSBarry 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 } 219e3372554SBarry Smith if (info) SETERRQ(1,0,"Bad solve"); 22055659b69SBarry Smith PLogFlops(mat->n*mat->n - mat->n); 221*3a40ed3dSBarry 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 233*3a40ed3dSBarry 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 ); 244da3a660dSBarry Smith } 245da3a660dSBarry Smith else { 24648d91487SBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 247da3a660dSBarry Smith } 248e3372554SBarry Smith if (info) SETERRQ(1,0,"Bad solve"); 249da3a660dSBarry Smith if (tmp) {VecAXPY(&sone,tmp,yy); VecDestroy(tmp);} 250da3a660dSBarry Smith else VecAXPY(&sone,zz,yy); 25155659b69SBarry Smith PLogFlops(mat->n*mat->n - mat->n); 252*3a40ed3dSBarry Smith PetscFunctionReturn(0); 253da3a660dSBarry Smith } 25467e560aaSBarry Smith 2555615d1e5SSatish Balay #undef __FUNC__ 2565615d1e5SSatish Balay #define __FUNC__ "MatSolveTransAdd_SeqDense" 2578f6be9afSLois Curfman McInnes int MatSolveTransAdd_SeqDense(Mat A,Vec xx,Vec zz, Vec yy) 258da3a660dSBarry Smith { 259c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 2606abc6512SBarry Smith int one = 1, info,ierr; 2616abc6512SBarry Smith Scalar *x, *y, sone = 1.0; 262da3a660dSBarry Smith Vec tmp; 26367e560aaSBarry Smith 264*3a40ed3dSBarry Smith PetscFunctionBegin; 265da3a660dSBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 266da3a660dSBarry Smith if (yy == zz) { 26778b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp); CHKERRQ(ierr); 268c0bbcb79SLois Curfman McInnes PLogObjectParent(A,tmp); 26978b31e54SBarry Smith ierr = VecCopy(yy,tmp); CHKERRQ(ierr); 270da3a660dSBarry Smith } 271416022c9SBarry Smith PetscMemcpy(y,x,mat->m*sizeof(Scalar)); 272752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 273da3a660dSBarry Smith if (mat->pivots) { 27448d91487SBarry Smith LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 275*3a40ed3dSBarry Smith } else { 27648d91487SBarry Smith LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 277da3a660dSBarry Smith } 278e3372554SBarry Smith if (info) SETERRQ(1,0,"Bad solve"); 27990f02eecSBarry Smith if (tmp) { 28090f02eecSBarry Smith ierr = VecAXPY(&sone,tmp,yy); CHKERRQ(ierr); 28190f02eecSBarry Smith ierr = VecDestroy(tmp); CHKERRQ(ierr); 282*3a40ed3dSBarry Smith } else { 28390f02eecSBarry Smith ierr = VecAXPY(&sone,zz,yy); CHKERRQ(ierr); 28490f02eecSBarry Smith } 28555659b69SBarry Smith PLogFlops(mat->n*mat->n - mat->n); 286*3a40ed3dSBarry Smith PetscFunctionReturn(0); 287da3a660dSBarry Smith } 288289bc588SBarry Smith /* ------------------------------------------------------------------*/ 2895615d1e5SSatish Balay #undef __FUNC__ 2905615d1e5SSatish Balay #define __FUNC__ "MatRelax_SeqDense" 2918f6be9afSLois Curfman McInnes int MatRelax_SeqDense(Mat A,Vec bb,double omega,MatSORType flag, 29220e1a0d4SLois Curfman McInnes double shift,int its,Vec xx) 293289bc588SBarry Smith { 294c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 295289bc588SBarry Smith Scalar *x, *b, *v = mat->v, zero = 0.0, xt; 2966abc6512SBarry Smith int o = 1,ierr, m = mat->m, i; 297289bc588SBarry Smith 298*3a40ed3dSBarry Smith PetscFunctionBegin; 299289bc588SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 300*3a40ed3dSBarry Smith /* this is a hack fix, should have another version without the second BLdot */ 301bbb6d6a8SBarry Smith ierr = VecSet(&zero,xx); CHKERRQ(ierr); 302289bc588SBarry Smith } 303289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(bb,&b); 304289bc588SBarry Smith while (its--) { 305289bc588SBarry Smith if (flag & SOR_FORWARD_SWEEP){ 306289bc588SBarry Smith for ( i=0; i<m; i++ ) { 307*3a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 308f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 309f1747703SBarry Smith not happy about returning a double complex */ 310f1747703SBarry Smith int _i; 311f1747703SBarry Smith Scalar sum = b[i]; 312f1747703SBarry Smith for ( _i=0; _i<m; _i++ ) { 313f1747703SBarry Smith sum -= conj(v[i+_i*m])*x[_i]; 314f1747703SBarry Smith } 315f1747703SBarry Smith xt = sum; 316f1747703SBarry Smith #else 317289bc588SBarry Smith xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 318f1747703SBarry Smith #endif 31955a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 320289bc588SBarry Smith } 321289bc588SBarry Smith } 322289bc588SBarry Smith if (flag & SOR_BACKWARD_SWEEP) { 323289bc588SBarry Smith for ( i=m-1; i>=0; i-- ) { 324*3a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 325f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 326f1747703SBarry Smith not happy about returning a double complex */ 327f1747703SBarry Smith int _i; 328f1747703SBarry Smith Scalar sum = b[i]; 329f1747703SBarry Smith for ( _i=0; _i<m; _i++ ) { 330f1747703SBarry Smith sum -= conj(v[i+_i*m])*x[_i]; 331f1747703SBarry Smith } 332f1747703SBarry Smith xt = sum; 333f1747703SBarry Smith #else 334289bc588SBarry Smith xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 335f1747703SBarry Smith #endif 33655a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 337289bc588SBarry Smith } 338289bc588SBarry Smith } 339289bc588SBarry Smith } 340*3a40ed3dSBarry Smith PetscFunctionReturn(0); 341289bc588SBarry Smith } 342289bc588SBarry Smith 343289bc588SBarry Smith /* -----------------------------------------------------------------*/ 3445615d1e5SSatish Balay #undef __FUNC__ 3455615d1e5SSatish Balay #define __FUNC__ "MatMultTrans_SeqDense" 34644cd7ae7SLois Curfman McInnes int MatMultTrans_SeqDense(Mat A,Vec xx,Vec yy) 347289bc588SBarry Smith { 348c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 349289bc588SBarry Smith Scalar *v = mat->v, *x, *y; 350289bc588SBarry Smith int _One=1;Scalar _DOne=1.0, _DZero=0.0; 351*3a40ed3dSBarry Smith 352*3a40ed3dSBarry Smith PetscFunctionBegin; 353289bc588SBarry Smith VecGetArray(xx,&x), VecGetArray(yy,&y); 35448d91487SBarry Smith LAgemv_("T",&(mat->m),&(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DZero,y,&_One); 35555659b69SBarry Smith PLogFlops(2*mat->m*mat->n - mat->n); 356*3a40ed3dSBarry Smith PetscFunctionReturn(0); 357289bc588SBarry Smith } 3586ee01492SSatish Balay 3595615d1e5SSatish Balay #undef __FUNC__ 3605615d1e5SSatish Balay #define __FUNC__ "MatMult_SeqDense" 36144cd7ae7SLois Curfman McInnes int MatMult_SeqDense(Mat A,Vec xx,Vec yy) 362289bc588SBarry Smith { 363c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 364289bc588SBarry Smith Scalar *v = mat->v, *x, *y; 365289bc588SBarry Smith int _One=1;Scalar _DOne=1.0, _DZero=0.0; 366*3a40ed3dSBarry Smith 367*3a40ed3dSBarry Smith PetscFunctionBegin; 368289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); 36948d91487SBarry Smith LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DZero,y,&_One); 37055659b69SBarry Smith PLogFlops(2*mat->m*mat->n - mat->m); 371*3a40ed3dSBarry Smith PetscFunctionReturn(0); 372289bc588SBarry Smith } 3736ee01492SSatish Balay 3745615d1e5SSatish Balay #undef __FUNC__ 3755615d1e5SSatish Balay #define __FUNC__ "MatMultAdd_SeqDense" 37644cd7ae7SLois Curfman McInnes int MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 377289bc588SBarry Smith { 378c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 379289bc588SBarry Smith Scalar *v = mat->v, *x, *y, *z; 3806abc6512SBarry Smith int _One=1; Scalar _DOne=1.0; 381*3a40ed3dSBarry Smith 382*3a40ed3dSBarry Smith PetscFunctionBegin; 383289bc588SBarry Smith VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z); 384416022c9SBarry Smith if (zz != yy) PetscMemcpy(y,z,mat->m*sizeof(Scalar)); 38548d91487SBarry Smith LAgemv_( "N", &(mat->m), &(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DOne,y,&_One); 38655659b69SBarry Smith PLogFlops(2*mat->m*mat->n); 387*3a40ed3dSBarry Smith PetscFunctionReturn(0); 388289bc588SBarry Smith } 3896ee01492SSatish Balay 3905615d1e5SSatish Balay #undef __FUNC__ 3915615d1e5SSatish Balay #define __FUNC__ "MatMultTransAdd_SeqDense" 39244cd7ae7SLois Curfman McInnes int MatMultTransAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 393289bc588SBarry Smith { 394c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 395289bc588SBarry Smith Scalar *v = mat->v, *x, *y, *z; 396289bc588SBarry Smith int _One=1; 3976abc6512SBarry Smith Scalar _DOne=1.0; 398*3a40ed3dSBarry Smith 399*3a40ed3dSBarry Smith PetscFunctionBegin; 40044cd7ae7SLois Curfman McInnes VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z); 401717eeb74SLois Curfman McInnes if (zz != yy) PetscMemcpy(y,z,mat->n*sizeof(Scalar)); 40248d91487SBarry Smith LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DOne,y,&_One); 40355659b69SBarry Smith PLogFlops(2*mat->m*mat->n); 404*3a40ed3dSBarry Smith PetscFunctionReturn(0); 405289bc588SBarry Smith } 406289bc588SBarry Smith 407289bc588SBarry Smith /* -----------------------------------------------------------------*/ 4085615d1e5SSatish Balay #undef __FUNC__ 4095615d1e5SSatish Balay #define __FUNC__ "MatGetRow_SeqDense" 4108f6be9afSLois Curfman McInnes int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals) 411289bc588SBarry Smith { 412c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 413289bc588SBarry Smith Scalar *v; 414289bc588SBarry Smith int i; 41567e560aaSBarry Smith 416*3a40ed3dSBarry Smith PetscFunctionBegin; 417289bc588SBarry Smith *ncols = mat->n; 418289bc588SBarry Smith if (cols) { 41945d48df9SBarry Smith *cols = (int *) PetscMalloc((mat->n+1)*sizeof(int)); CHKPTRQ(*cols); 42073c3ce57SBarry Smith for ( i=0; i<mat->n; i++ ) (*cols)[i] = i; 421289bc588SBarry Smith } 422289bc588SBarry Smith if (vals) { 42345d48df9SBarry Smith *vals = (Scalar *) PetscMalloc((mat->n+1)*sizeof(Scalar)); CHKPTRQ(*vals); 424289bc588SBarry Smith v = mat->v + row; 42573c3ce57SBarry Smith for ( i=0; i<mat->n; i++ ) {(*vals)[i] = *v; v += mat->m;} 426289bc588SBarry Smith } 427*3a40ed3dSBarry Smith PetscFunctionReturn(0); 428289bc588SBarry Smith } 4296ee01492SSatish Balay 4305615d1e5SSatish Balay #undef __FUNC__ 4315615d1e5SSatish Balay #define __FUNC__ "MatRestoreRow_SeqDense" 4328f6be9afSLois Curfman McInnes int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals) 433289bc588SBarry Smith { 4340452661fSBarry Smith if (cols) { PetscFree(*cols); } 4350452661fSBarry Smith if (vals) { PetscFree(*vals); } 436*3a40ed3dSBarry Smith PetscFunctionReturn(0); 437289bc588SBarry Smith } 438289bc588SBarry Smith /* ----------------------------------------------------------------*/ 4395615d1e5SSatish Balay #undef __FUNC__ 4405615d1e5SSatish Balay #define __FUNC__ "MatSetValues_SeqDense" 4418f6be9afSLois Curfman McInnes int MatSetValues_SeqDense(Mat A,int m,int *indexm,int n, 442e8d4e0b9SBarry Smith int *indexn,Scalar *v,InsertMode addv) 443289bc588SBarry Smith { 444c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 445289bc588SBarry Smith int i,j; 446d6dfbf8fSBarry Smith 447*3a40ed3dSBarry Smith PetscFunctionBegin; 448289bc588SBarry Smith if (!mat->roworiented) { 449dbb450caSBarry Smith if (addv == INSERT_VALUES) { 450289bc588SBarry Smith for ( j=0; j<n; j++ ) { 451*3a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g) 45258804f6dSBarry Smith if (indexn[j] < 0) SETERRQ(1,0,"Negative column"); 45358804f6dSBarry Smith if (indexn[j] >= A->n) SETERRQ(1,0,"Column too large"); 45458804f6dSBarry Smith #endif 455289bc588SBarry Smith for ( i=0; i<m; i++ ) { 456*3a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g) 45758804f6dSBarry Smith if (indexm[i] < 0) SETERRQ(1,0,"Negative row"); 45858804f6dSBarry Smith if (indexm[i] >= A->m) SETERRQ(1,0,"Row too large"); 45958804f6dSBarry Smith #endif 460289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 461289bc588SBarry Smith } 462289bc588SBarry Smith } 463*3a40ed3dSBarry Smith } else { 464289bc588SBarry Smith for ( j=0; j<n; j++ ) { 465*3a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g) 46658804f6dSBarry Smith if (indexn[j] < 0) SETERRQ(1,0,"Negative column"); 46758804f6dSBarry Smith if (indexn[j] >= A->n) SETERRQ(1,0,"Column too large"); 46858804f6dSBarry Smith #endif 469289bc588SBarry Smith for ( i=0; i<m; i++ ) { 470*3a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g) 47158804f6dSBarry Smith if (indexm[i] < 0) SETERRQ(1,0,"Negative row"); 47258804f6dSBarry Smith if (indexm[i] >= A->m) SETERRQ(1,0,"Row too large"); 47358804f6dSBarry Smith #endif 474289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 475289bc588SBarry Smith } 476289bc588SBarry Smith } 477289bc588SBarry Smith } 478*3a40ed3dSBarry Smith } else { 479dbb450caSBarry Smith if (addv == INSERT_VALUES) { 480e8d4e0b9SBarry Smith for ( i=0; i<m; i++ ) { 481*3a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g) 48258804f6dSBarry Smith if (indexm[i] < 0) SETERRQ(1,0,"Negative row"); 48358804f6dSBarry Smith if (indexm[i] >= A->m) SETERRQ(1,0,"Row too large"); 48458804f6dSBarry Smith #endif 485e8d4e0b9SBarry Smith for ( j=0; j<n; j++ ) { 486*3a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g) 48758804f6dSBarry Smith if (indexn[j] < 0) SETERRQ(1,0,"Negative column"); 48858804f6dSBarry Smith if (indexn[j] >= A->n) SETERRQ(1,0,"Column too large"); 48958804f6dSBarry Smith #endif 490e8d4e0b9SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 491e8d4e0b9SBarry Smith } 492e8d4e0b9SBarry Smith } 493*3a40ed3dSBarry Smith } else { 494289bc588SBarry Smith for ( i=0; i<m; i++ ) { 495*3a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g) 49658804f6dSBarry Smith if (indexm[i] < 0) SETERRQ(1,0,"Negative row"); 49758804f6dSBarry Smith if (indexm[i] >= A->m) SETERRQ(1,0,"Row too large"); 49858804f6dSBarry Smith #endif 499289bc588SBarry Smith for ( j=0; j<n; j++ ) { 500*3a40ed3dSBarry Smith #if defined(USE_PETSC_BOPT_g) 50158804f6dSBarry Smith if (indexn[j] < 0) SETERRQ(1,0,"Negative column"); 50258804f6dSBarry Smith if (indexn[j] >= A->n) SETERRQ(1,0,"Column too large"); 50358804f6dSBarry Smith #endif 504289bc588SBarry Smith mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 505289bc588SBarry Smith } 506289bc588SBarry Smith } 507289bc588SBarry Smith } 508e8d4e0b9SBarry Smith } 509*3a40ed3dSBarry Smith PetscFunctionReturn(0); 510289bc588SBarry Smith } 511e8d4e0b9SBarry Smith 5125615d1e5SSatish Balay #undef __FUNC__ 5135615d1e5SSatish Balay #define __FUNC__ "MatGetValues_SeqDense" 5148f6be9afSLois Curfman McInnes int MatGetValues_SeqDense(Mat A,int m,int *indexm,int n,int *indexn,Scalar *v) 515ae80bb75SLois Curfman McInnes { 516ae80bb75SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 517ae80bb75SLois Curfman McInnes int i, j; 518ae80bb75SLois Curfman McInnes Scalar *vpt = v; 519ae80bb75SLois Curfman McInnes 520*3a40ed3dSBarry Smith PetscFunctionBegin; 521ae80bb75SLois Curfman McInnes /* row-oriented output */ 522ae80bb75SLois Curfman McInnes for ( i=0; i<m; i++ ) { 523ae80bb75SLois Curfman McInnes for ( j=0; j<n; j++ ) { 524ae80bb75SLois Curfman McInnes *vpt++ = mat->v[indexn[j]*mat->m + indexm[i]]; 525ae80bb75SLois Curfman McInnes } 526ae80bb75SLois Curfman McInnes } 527*3a40ed3dSBarry Smith PetscFunctionReturn(0); 528ae80bb75SLois Curfman McInnes } 529ae80bb75SLois Curfman McInnes 530289bc588SBarry Smith /* -----------------------------------------------------------------*/ 531289bc588SBarry Smith 53277c4ece6SBarry Smith #include "sys.h" 53388e32edaSLois Curfman McInnes 5345615d1e5SSatish Balay #undef __FUNC__ 5355615d1e5SSatish Balay #define __FUNC__ "MatLoad_SeqDense" 53619bcc07fSBarry Smith int MatLoad_SeqDense(Viewer viewer,MatType type,Mat *A) 53788e32edaSLois Curfman McInnes { 53888e32edaSLois Curfman McInnes Mat_SeqDense *a; 53988e32edaSLois Curfman McInnes Mat B; 54088e32edaSLois Curfman McInnes int *scols, i, j, nz, ierr, fd, header[4], size; 54188e32edaSLois Curfman McInnes int *rowlengths = 0, M, N, *cols; 54288e32edaSLois Curfman McInnes Scalar *vals, *svals, *v; 54319bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 54488e32edaSLois Curfman McInnes 545*3a40ed3dSBarry Smith PetscFunctionBegin; 54688e32edaSLois Curfman McInnes MPI_Comm_size(comm,&size); 547e3372554SBarry Smith if (size > 1) SETERRQ(1,0,"view must have one processor"); 54890ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 5490752156aSBarry Smith ierr = PetscBinaryRead(fd,header,4,PETSC_INT); CHKERRQ(ierr); 550e3372554SBarry Smith if (header[0] != MAT_COOKIE) SETERRQ(1,0,"Not matrix object"); 55188e32edaSLois Curfman McInnes M = header[1]; N = header[2]; nz = header[3]; 55288e32edaSLois Curfman McInnes 55390ace30eSBarry Smith if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */ 55490ace30eSBarry Smith ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr); 55590ace30eSBarry Smith B = *A; 55690ace30eSBarry Smith a = (Mat_SeqDense *) B->data; 55790ace30eSBarry Smith 55890ace30eSBarry Smith /* read in nonzero values */ 5590752156aSBarry Smith ierr = PetscBinaryRead(fd,a->v,M*N,PETSC_SCALAR); CHKERRQ(ierr); 56090ace30eSBarry Smith 5616d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 5626d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 56390ace30eSBarry Smith ierr = MatTranspose(B,PETSC_NULL); CHKERRQ(ierr); 56490ace30eSBarry Smith } else { 56588e32edaSLois Curfman McInnes /* read row lengths */ 56645d48df9SBarry Smith rowlengths = (int*) PetscMalloc( (M+1)*sizeof(int) ); CHKPTRQ(rowlengths); 5670752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT); CHKERRQ(ierr); 56888e32edaSLois Curfman McInnes 56988e32edaSLois Curfman McInnes /* create our matrix */ 570b4fd4287SBarry Smith ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A); CHKERRQ(ierr); 57188e32edaSLois Curfman McInnes B = *A; 57288e32edaSLois Curfman McInnes a = (Mat_SeqDense *) B->data; 57388e32edaSLois Curfman McInnes v = a->v; 57488e32edaSLois Curfman McInnes 57588e32edaSLois Curfman McInnes /* read column indices and nonzeros */ 57645d48df9SBarry Smith cols = scols = (int *) PetscMalloc( (nz+1)*sizeof(int) ); CHKPTRQ(cols); 5770752156aSBarry Smith ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT); CHKERRQ(ierr); 57845d48df9SBarry Smith vals = svals = (Scalar *) PetscMalloc( (nz+1)*sizeof(Scalar) ); CHKPTRQ(vals); 5790752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR); CHKERRQ(ierr); 58088e32edaSLois Curfman McInnes 58188e32edaSLois Curfman McInnes /* insert into matrix */ 58288e32edaSLois Curfman McInnes for ( i=0; i<M; i++ ) { 58388e32edaSLois Curfman McInnes for ( j=0; j<rowlengths[i]; j++ ) v[i+M*scols[j]] = svals[j]; 58488e32edaSLois Curfman McInnes svals += rowlengths[i]; scols += rowlengths[i]; 58588e32edaSLois Curfman McInnes } 5860452661fSBarry Smith PetscFree(vals); PetscFree(cols); PetscFree(rowlengths); 58788e32edaSLois Curfman McInnes 5886d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 5896d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 59090ace30eSBarry Smith } 591*3a40ed3dSBarry Smith PetscFunctionReturn(0); 59288e32edaSLois Curfman McInnes } 59388e32edaSLois Curfman McInnes 594932b0c3eSLois Curfman McInnes #include "pinclude/pviewer.h" 59577c4ece6SBarry Smith #include "sys.h" 596932b0c3eSLois Curfman McInnes 5975615d1e5SSatish Balay #undef __FUNC__ 5985615d1e5SSatish Balay #define __FUNC__ "MatView_SeqDense_ASCII" 599932b0c3eSLois Curfman McInnes static int MatView_SeqDense_ASCII(Mat A,Viewer viewer) 600289bc588SBarry Smith { 601932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense *) A->data; 602932b0c3eSLois Curfman McInnes int ierr, i, j, format; 603d636dbe3SBarry Smith FILE *fd; 604932b0c3eSLois Curfman McInnes char *outputname; 605932b0c3eSLois Curfman McInnes Scalar *v; 606932b0c3eSLois Curfman McInnes 607*3a40ed3dSBarry Smith PetscFunctionBegin; 60890ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 609932b0c3eSLois Curfman McInnes ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr); 61090ace30eSBarry Smith ierr = ViewerGetFormat(viewer,&format); 611639f9d9dSBarry Smith if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) { 612*3a40ed3dSBarry Smith PetscFunctionReturn(0); /* do nothing for now */ 613f72585f2SLois Curfman McInnes } 61402cad45dSBarry Smith else if (format == VIEWER_FORMAT_ASCII_COMMON) { 61580cd9d93SLois Curfman McInnes for ( i=0; i<a->m; i++ ) { 61644cd7ae7SLois Curfman McInnes v = a->v + i; 61780cd9d93SLois Curfman McInnes fprintf(fd,"row %d:",i); 61880cd9d93SLois Curfman McInnes for ( j=0; j<a->n; j++ ) { 619*3a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 62080cd9d93SLois Curfman McInnes if (real(*v) != 0.0 && imag(*v) != 0.0) 62180cd9d93SLois Curfman McInnes fprintf(fd," %d %g + %g i",j,real(*v),imag(*v)); 62280cd9d93SLois Curfman McInnes else if (real(*v)) fprintf(fd," %d %g ",j,real(*v)); 62380cd9d93SLois Curfman McInnes v += a->m; 62480cd9d93SLois Curfman McInnes #else 62580cd9d93SLois Curfman McInnes if (*v) fprintf(fd," %d %g ",j,*v); 62680cd9d93SLois Curfman McInnes v += a->m; 62780cd9d93SLois Curfman McInnes #endif 62880cd9d93SLois Curfman McInnes } 62980cd9d93SLois Curfman McInnes fprintf(fd,"\n"); 63080cd9d93SLois Curfman McInnes } 631*3a40ed3dSBarry Smith } else { 632*3a40ed3dSBarry 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++ ) { 643*3a40ed3dSBarry 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); 657*3a40ed3dSBarry 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 669*3a40ed3dSBarry 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 } 725*3a40ed3dSBarry Smith PetscFunctionReturn(0); 726932b0c3eSLois Curfman McInnes } 727932b0c3eSLois Curfman McInnes 7285615d1e5SSatish Balay #undef __FUNC__ 7295615d1e5SSatish Balay #define __FUNC__ "MatView_SeqDense" 7308f6be9afSLois Curfman McInnes int MatView_SeqDense(PetscObject obj,Viewer viewer) 731932b0c3eSLois Curfman McInnes { 732932b0c3eSLois Curfman McInnes Mat A = (Mat) obj; 733932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*) A->data; 734bcd2baecSBarry Smith ViewerType vtype; 735bcd2baecSBarry Smith int ierr; 736932b0c3eSLois Curfman McInnes 737*3a40ed3dSBarry Smith PetscFunctionBegin; 738bcd2baecSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 739bcd2baecSBarry Smith 740bcd2baecSBarry Smith if (vtype == MATLAB_VIEWER) { 741*3a40ed3dSBarry Smith ierr = ViewerMatlabPutArray_Private(viewer,a->m,a->n,a->v); CHKERRQ(ierr); 742*3a40ed3dSBarry Smith } else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) { 743*3a40ed3dSBarry Smith ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr); 744*3a40ed3dSBarry Smith } else if (vtype == BINARY_FILE_VIEWER) { 745*3a40ed3dSBarry Smith ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr); 746932b0c3eSLois Curfman McInnes } 747*3a40ed3dSBarry Smith PetscFunctionReturn(0); 748932b0c3eSLois Curfman McInnes } 749289bc588SBarry Smith 7505615d1e5SSatish Balay #undef __FUNC__ 7515615d1e5SSatish Balay #define __FUNC__ "MatDestroy_SeqDense" 7528f6be9afSLois Curfman McInnes int MatDestroy_SeqDense(PetscObject obj) 753289bc588SBarry Smith { 754289bc588SBarry Smith Mat mat = (Mat) obj; 755ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense *) mat->data; 75690f02eecSBarry Smith int ierr; 75790f02eecSBarry Smith 758*3a40ed3dSBarry Smith PetscFunctionBegin; 759*3a40ed3dSBarry Smith #if defined(USE_PETSC_LOG) 760a5a9c739SBarry Smith PLogObjectState(obj,"Rows %d Cols %d",l->m,l->n); 761a5a9c739SBarry Smith #endif 7620452661fSBarry Smith if (l->pivots) PetscFree(l->pivots); 7633439631bSBarry Smith if (!l->user_alloc) PetscFree(l->v); 7640452661fSBarry Smith PetscFree(l); 76590f02eecSBarry Smith if (mat->mapping) { 76690f02eecSBarry Smith ierr = ISLocalToGlobalMappingDestroy(mat->mapping); CHKERRQ(ierr); 76790f02eecSBarry Smith } 768a5a9c739SBarry Smith PLogObjectDestroy(mat); 7690452661fSBarry Smith PetscHeaderDestroy(mat); 770*3a40ed3dSBarry Smith PetscFunctionReturn(0); 771289bc588SBarry Smith } 772289bc588SBarry Smith 7735615d1e5SSatish Balay #undef __FUNC__ 7745615d1e5SSatish Balay #define __FUNC__ "MatTranspose_SeqDense" 7758f6be9afSLois Curfman McInnes int MatTranspose_SeqDense(Mat A,Mat *matout) 776289bc588SBarry Smith { 777c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 778d3e5ee88SLois Curfman McInnes int k, j, m, n; 779d3e5ee88SLois Curfman McInnes Scalar *v, tmp; 78048b35521SBarry Smith 781*3a40ed3dSBarry Smith PetscFunctionBegin; 782d3e5ee88SLois Curfman McInnes v = mat->v; m = mat->m; n = mat->n; 7833638b69dSLois Curfman McInnes if (matout == PETSC_NULL) { /* in place transpose */ 784e3372554SBarry Smith if (m != n) SETERRQ(PETSC_ERR_SUP,0,"Supports square matrix only in-place"); 785d3e5ee88SLois Curfman McInnes for ( j=0; j<m; j++ ) { 786289bc588SBarry Smith for ( k=0; k<j; k++ ) { 787d3e5ee88SLois Curfman McInnes tmp = v[j + k*n]; 788d3e5ee88SLois Curfman McInnes v[j + k*n] = v[k + j*n]; 789d3e5ee88SLois Curfman McInnes v[k + j*n] = tmp; 790289bc588SBarry Smith } 791289bc588SBarry Smith } 792*3a40ed3dSBarry Smith } else { /* out-of-place transpose */ 793d3e5ee88SLois Curfman McInnes int ierr; 794d3e5ee88SLois Curfman McInnes Mat tmat; 795ec8511deSBarry Smith Mat_SeqDense *tmatd; 796d3e5ee88SLois Curfman McInnes Scalar *v2; 797b4fd4287SBarry Smith ierr = MatCreateSeqDense(A->comm,mat->n,mat->m,PETSC_NULL,&tmat); CHKERRQ(ierr); 798ec8511deSBarry Smith tmatd = (Mat_SeqDense *) tmat->data; 7990de55854SLois Curfman McInnes v = mat->v; v2 = tmatd->v; 800d3e5ee88SLois Curfman McInnes for ( j=0; j<n; j++ ) { 8010de55854SLois Curfman McInnes for ( k=0; k<m; k++ ) v2[j + k*n] = v[k + j*m]; 802d3e5ee88SLois Curfman McInnes } 8036d4a8577SBarry Smith ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 8046d4a8577SBarry Smith ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 805d3e5ee88SLois Curfman McInnes *matout = tmat; 80648b35521SBarry Smith } 807*3a40ed3dSBarry Smith PetscFunctionReturn(0); 808289bc588SBarry Smith } 809289bc588SBarry Smith 8105615d1e5SSatish Balay #undef __FUNC__ 8115615d1e5SSatish Balay #define __FUNC__ "MatEqual_SeqDense" 8128f6be9afSLois Curfman McInnes int MatEqual_SeqDense(Mat A1,Mat A2, PetscTruth *flg) 813289bc588SBarry Smith { 814c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat1 = (Mat_SeqDense *) A1->data; 815c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat2 = (Mat_SeqDense *) A2->data; 816289bc588SBarry Smith int i; 817289bc588SBarry Smith Scalar *v1 = mat1->v, *v2 = mat2->v; 8189ea5d5aeSSatish Balay 819*3a40ed3dSBarry Smith PetscFunctionBegin; 820*3a40ed3dSBarry Smith if (A2->type != MATSEQDENSE) SETERRQ(PETSC_ERR_SUP,0,"Matrices should be of type MATSEQDENSE"); 821*3a40ed3dSBarry Smith if (mat1->m != mat2->m) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 822*3a40ed3dSBarry Smith if (mat1->n != mat2->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 823289bc588SBarry Smith for ( i=0; i<mat1->m*mat1->n; i++ ) { 824*3a40ed3dSBarry Smith if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 825289bc588SBarry Smith v1++; v2++; 826289bc588SBarry Smith } 82777c4ece6SBarry Smith *flg = PETSC_TRUE; 828*3a40ed3dSBarry Smith PetscFunctionReturn(0); 829289bc588SBarry Smith } 830289bc588SBarry Smith 8315615d1e5SSatish Balay #undef __FUNC__ 8325615d1e5SSatish Balay #define __FUNC__ "MatGetDiagonal_SeqDense" 8338f6be9afSLois Curfman McInnes int MatGetDiagonal_SeqDense(Mat A,Vec v) 834289bc588SBarry Smith { 835c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 83644cd7ae7SLois Curfman McInnes int i, n, len; 83744cd7ae7SLois Curfman McInnes Scalar *x, zero = 0.0; 83844cd7ae7SLois Curfman McInnes 839*3a40ed3dSBarry Smith PetscFunctionBegin; 84044cd7ae7SLois Curfman McInnes VecSet(&zero,v); 841289bc588SBarry Smith VecGetArray(v,&x); VecGetSize(v,&n); 84244cd7ae7SLois Curfman McInnes len = PetscMin(mat->m,mat->n); 843e3372554SBarry Smith if (n != mat->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Nonconforming mat and vec"); 84444cd7ae7SLois Curfman McInnes for ( i=0; i<len; i++ ) { 845289bc588SBarry Smith x[i] = mat->v[i*mat->m + i]; 846289bc588SBarry Smith } 847*3a40ed3dSBarry Smith PetscFunctionReturn(0); 848289bc588SBarry Smith } 849289bc588SBarry Smith 8505615d1e5SSatish Balay #undef __FUNC__ 8515615d1e5SSatish Balay #define __FUNC__ "MatDiagonalScale_SeqDense" 8528f6be9afSLois Curfman McInnes int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr) 853289bc588SBarry Smith { 854c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 855da3a660dSBarry Smith Scalar *l,*r,x,*v; 856da3a660dSBarry Smith int i,j,m = mat->m, n = mat->n; 85755659b69SBarry Smith 858*3a40ed3dSBarry Smith PetscFunctionBegin; 85928988994SBarry Smith if (ll) { 860da3a660dSBarry Smith VecGetArray(ll,&l); VecGetSize(ll,&m); 861e3372554SBarry Smith if (m != mat->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Left scaling vec wrong size"); 862da3a660dSBarry Smith for ( i=0; i<m; i++ ) { 863da3a660dSBarry Smith x = l[i]; 864da3a660dSBarry Smith v = mat->v + i; 865da3a660dSBarry Smith for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;} 866da3a660dSBarry Smith } 86744cd7ae7SLois Curfman McInnes PLogFlops(n*m); 868da3a660dSBarry Smith } 86928988994SBarry Smith if (rr) { 870da3a660dSBarry Smith VecGetArray(rr,&r); VecGetSize(rr,&n); 871e3372554SBarry Smith if (n != mat->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Right scaling vec wrong size"); 872da3a660dSBarry Smith for ( i=0; i<n; i++ ) { 873da3a660dSBarry Smith x = r[i]; 874da3a660dSBarry Smith v = mat->v + i*m; 875da3a660dSBarry Smith for ( j=0; j<m; j++ ) { (*v++) *= x;} 876da3a660dSBarry Smith } 87744cd7ae7SLois Curfman McInnes PLogFlops(n*m); 878da3a660dSBarry Smith } 879*3a40ed3dSBarry Smith PetscFunctionReturn(0); 880289bc588SBarry Smith } 881289bc588SBarry Smith 8825615d1e5SSatish Balay #undef __FUNC__ 8835615d1e5SSatish Balay #define __FUNC__ "MatNorm_SeqDense" 8848f6be9afSLois Curfman McInnes int MatNorm_SeqDense(Mat A,NormType type,double *norm) 885289bc588SBarry Smith { 886c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 887289bc588SBarry Smith Scalar *v = mat->v; 888289bc588SBarry Smith double sum = 0.0; 889289bc588SBarry Smith int i, j; 89055659b69SBarry Smith 891*3a40ed3dSBarry Smith PetscFunctionBegin; 892289bc588SBarry Smith if (type == NORM_FROBENIUS) { 893289bc588SBarry Smith for (i=0; i<mat->n*mat->m; i++ ) { 894*3a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 895289bc588SBarry Smith sum += real(conj(*v)*(*v)); v++; 896289bc588SBarry Smith #else 897289bc588SBarry Smith sum += (*v)*(*v); v++; 898289bc588SBarry Smith #endif 899289bc588SBarry Smith } 900289bc588SBarry Smith *norm = sqrt(sum); 90155659b69SBarry Smith PLogFlops(2*mat->n*mat->m); 902*3a40ed3dSBarry Smith } else if (type == NORM_1) { 903289bc588SBarry Smith *norm = 0.0; 904289bc588SBarry Smith for ( j=0; j<mat->n; j++ ) { 905289bc588SBarry Smith sum = 0.0; 906289bc588SBarry Smith for ( i=0; i<mat->m; i++ ) { 90733a8263dSBarry Smith sum += PetscAbsScalar(*v); v++; 908289bc588SBarry Smith } 909289bc588SBarry Smith if (sum > *norm) *norm = sum; 910289bc588SBarry Smith } 91155659b69SBarry Smith PLogFlops(mat->n*mat->m); 912*3a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { 913289bc588SBarry Smith *norm = 0.0; 914289bc588SBarry Smith for ( j=0; j<mat->m; j++ ) { 915289bc588SBarry Smith v = mat->v + j; 916289bc588SBarry Smith sum = 0.0; 917289bc588SBarry Smith for ( i=0; i<mat->n; i++ ) { 918cddf8d76SBarry Smith sum += PetscAbsScalar(*v); v += mat->m; 919289bc588SBarry Smith } 920289bc588SBarry Smith if (sum > *norm) *norm = sum; 921289bc588SBarry Smith } 92255659b69SBarry Smith PLogFlops(mat->n*mat->m); 923*3a40ed3dSBarry Smith } else { 924e3372554SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"No two norm"); 925289bc588SBarry Smith } 926*3a40ed3dSBarry Smith PetscFunctionReturn(0); 927289bc588SBarry Smith } 928289bc588SBarry Smith 9295615d1e5SSatish Balay #undef __FUNC__ 9305615d1e5SSatish Balay #define __FUNC__ "MatSetOption_SeqDense" 9318f6be9afSLois Curfman McInnes int MatSetOption_SeqDense(Mat A,MatOption op) 932289bc588SBarry Smith { 933c0bbcb79SLois Curfman McInnes Mat_SeqDense *aij = (Mat_SeqDense *) A->data; 93467e560aaSBarry Smith 935*3a40ed3dSBarry Smith PetscFunctionBegin; 9366d4a8577SBarry Smith if (op == MAT_ROW_ORIENTED) aij->roworiented = 1; 9376d4a8577SBarry Smith else if (op == MAT_COLUMN_ORIENTED) aij->roworiented = 0; 9386d4a8577SBarry Smith else if (op == MAT_ROWS_SORTED || 939219d9a1aSLois Curfman McInnes op == MAT_ROWS_UNSORTED || 9406d4a8577SBarry Smith op == MAT_COLUMNS_SORTED || 941219d9a1aSLois Curfman McInnes op == MAT_COLUMNS_UNSORTED || 9426d4a8577SBarry Smith op == MAT_SYMMETRIC || 9436d4a8577SBarry Smith op == MAT_STRUCTURALLY_SYMMETRIC || 9446d4a8577SBarry Smith op == MAT_NO_NEW_NONZERO_LOCATIONS || 9456d4a8577SBarry Smith op == MAT_YES_NEW_NONZERO_LOCATIONS || 946c2653b3dSLois Curfman McInnes op == MAT_NEW_NONZERO_LOCATION_ERROR || 9476d4a8577SBarry Smith op == MAT_NO_NEW_DIAGONALS || 94890f02eecSBarry Smith op == MAT_YES_NEW_DIAGONALS || 9492b362799SSatish Balay op == MAT_IGNORE_OFF_PROC_ENTRIES) 95094a424c1SBarry Smith PLogInfo(A,"Info:MatSetOption_SeqDense:Option ignored\n"); 951*3a40ed3dSBarry Smith else { 952*3a40ed3dSBarry Smith SETERRQ(PETSC_ERR_SUP,0,"unknown option"); 953*3a40ed3dSBarry Smith } 954*3a40ed3dSBarry Smith PetscFunctionReturn(0); 955289bc588SBarry Smith } 956289bc588SBarry Smith 9575615d1e5SSatish Balay #undef __FUNC__ 9585615d1e5SSatish Balay #define __FUNC__ "MatZeroEntries_SeqDense" 9598f6be9afSLois Curfman McInnes int MatZeroEntries_SeqDense(Mat A) 9606f0a148fSBarry Smith { 961ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense *) A->data; 962*3a40ed3dSBarry Smith 963*3a40ed3dSBarry Smith PetscFunctionBegin; 964cddf8d76SBarry Smith PetscMemzero(l->v,l->m*l->n*sizeof(Scalar)); 965*3a40ed3dSBarry Smith PetscFunctionReturn(0); 9666f0a148fSBarry Smith } 9676f0a148fSBarry Smith 9685615d1e5SSatish Balay #undef __FUNC__ 9695615d1e5SSatish Balay #define __FUNC__ "MatGetBlockSize_SeqDense" 9708f6be9afSLois Curfman McInnes int MatGetBlockSize_SeqDense(Mat A,int *bs) 9714e220ebcSLois Curfman McInnes { 972*3a40ed3dSBarry Smith PetscFunctionBegin; 9734e220ebcSLois Curfman McInnes *bs = 1; 974*3a40ed3dSBarry Smith PetscFunctionReturn(0); 9754e220ebcSLois Curfman McInnes } 9764e220ebcSLois Curfman McInnes 9775615d1e5SSatish Balay #undef __FUNC__ 9785615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqDense" 9798f6be9afSLois Curfman McInnes int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag) 9806f0a148fSBarry Smith { 981ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense *) A->data; 9826abc6512SBarry Smith int n = l->n, i, j,ierr,N, *rows; 9836f0a148fSBarry Smith Scalar *slot; 98455659b69SBarry Smith 985*3a40ed3dSBarry Smith PetscFunctionBegin; 98677c4ece6SBarry Smith ierr = ISGetSize(is,&N); CHKERRQ(ierr); 98778b31e54SBarry Smith ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 9886f0a148fSBarry Smith for ( i=0; i<N; i++ ) { 9896f0a148fSBarry Smith slot = l->v + rows[i]; 9906f0a148fSBarry Smith for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;} 9916f0a148fSBarry Smith } 9926f0a148fSBarry Smith if (diag) { 9936f0a148fSBarry Smith for ( i=0; i<N; i++ ) { 9946f0a148fSBarry Smith slot = l->v + (n+1)*rows[i]; 995260925b8SBarry Smith *slot = *diag; 9966f0a148fSBarry Smith } 9976f0a148fSBarry Smith } 998260925b8SBarry Smith ISRestoreIndices(is,&rows); 999*3a40ed3dSBarry Smith PetscFunctionReturn(0); 10006f0a148fSBarry Smith } 1001557bce09SLois Curfman McInnes 10025615d1e5SSatish Balay #undef __FUNC__ 10035615d1e5SSatish Balay #define __FUNC__ "MatGetSize_SeqDense" 10048f6be9afSLois Curfman McInnes int MatGetSize_SeqDense(Mat A,int *m,int *n) 1005557bce09SLois Curfman McInnes { 1006c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 1007*3a40ed3dSBarry Smith 1008*3a40ed3dSBarry Smith PetscFunctionBegin; 1009557bce09SLois Curfman McInnes *m = mat->m; *n = mat->n; 1010*3a40ed3dSBarry Smith PetscFunctionReturn(0); 1011557bce09SLois Curfman McInnes } 1012557bce09SLois Curfman McInnes 10135615d1e5SSatish Balay #undef __FUNC__ 10145615d1e5SSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqDense" 10158f6be9afSLois Curfman McInnes int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n) 1016d3e5ee88SLois Curfman McInnes { 1017c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 1018*3a40ed3dSBarry Smith 1019*3a40ed3dSBarry Smith PetscFunctionBegin; 1020d3e5ee88SLois Curfman McInnes *m = 0; *n = mat->m; 1021*3a40ed3dSBarry Smith PetscFunctionReturn(0); 1022d3e5ee88SLois Curfman McInnes } 1023d3e5ee88SLois Curfman McInnes 10245615d1e5SSatish Balay #undef __FUNC__ 10255615d1e5SSatish Balay #define __FUNC__ "MatGetArray_SeqDense" 10268f6be9afSLois Curfman McInnes int MatGetArray_SeqDense(Mat A,Scalar **array) 102764e87e97SBarry Smith { 1028c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 1029*3a40ed3dSBarry Smith 1030*3a40ed3dSBarry Smith PetscFunctionBegin; 103164e87e97SBarry Smith *array = mat->v; 1032*3a40ed3dSBarry Smith PetscFunctionReturn(0); 103364e87e97SBarry Smith } 10340754003eSLois Curfman McInnes 10355615d1e5SSatish Balay #undef __FUNC__ 10365615d1e5SSatish Balay #define __FUNC__ "MatRestoreArray_SeqDense" 10378f6be9afSLois Curfman McInnes int MatRestoreArray_SeqDense(Mat A,Scalar **array) 1038ff14e315SSatish Balay { 1039*3a40ed3dSBarry Smith PetscFunctionBegin; 1040*3a40ed3dSBarry Smith PetscFunctionReturn(0); 1041ff14e315SSatish Balay } 10420754003eSLois Curfman McInnes 10435615d1e5SSatish Balay #undef __FUNC__ 10445615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrix_SeqDense" 1045cddf8d76SBarry Smith static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall, 1046cddf8d76SBarry Smith Mat *submat) 10470754003eSLois Curfman McInnes { 1048c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 10490754003eSLois Curfman McInnes int nznew, *smap, i, j, ierr, oldcols = mat->n; 1050160018dcSBarry Smith int *irow, *icol, nrows, ncols, *cwork; 10510754003eSLois Curfman McInnes Scalar *vwork, *val; 10520754003eSLois Curfman McInnes Mat newmat; 10530754003eSLois Curfman McInnes 1054*3a40ed3dSBarry Smith PetscFunctionBegin; 105578b31e54SBarry Smith ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr); 105678b31e54SBarry Smith ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr); 105778b31e54SBarry Smith ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr); 105878b31e54SBarry Smith ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr); 10590754003eSLois Curfman McInnes 10600452661fSBarry Smith smap = (int *) PetscMalloc(oldcols*sizeof(int)); CHKPTRQ(smap); 10610452661fSBarry Smith cwork = (int *) PetscMalloc(ncols*sizeof(int)); CHKPTRQ(cwork); 10620452661fSBarry Smith vwork = (Scalar *) PetscMalloc(ncols*sizeof(Scalar)); CHKPTRQ(vwork); 1063cddf8d76SBarry Smith PetscMemzero((char*)smap,oldcols*sizeof(int)); 10640754003eSLois Curfman McInnes for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1; 10650754003eSLois Curfman McInnes 10660754003eSLois Curfman McInnes /* Create and fill new matrix */ 1067b4fd4287SBarry Smith ierr = MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat); CHKERRQ(ierr); 10680754003eSLois Curfman McInnes for (i=0; i<nrows; i++) { 10690754003eSLois Curfman McInnes nznew = 0; 10700754003eSLois Curfman McInnes val = mat->v + irow[i]; 10710754003eSLois Curfman McInnes for (j=0; j<oldcols; j++) { 10720754003eSLois Curfman McInnes if (smap[j]) { 10730754003eSLois Curfman McInnes cwork[nznew] = smap[j] - 1; 10740754003eSLois Curfman McInnes vwork[nznew++] = val[j * mat->m]; 10750754003eSLois Curfman McInnes } 10760754003eSLois Curfman McInnes } 1077*3a40ed3dSBarry Smith ierr = MatSetValues(newmat,1,&i,nznew,cwork,vwork,INSERT_VALUES);CHKERRQ(ierr); 10780754003eSLois Curfman McInnes } 10796d4a8577SBarry Smith ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 10806d4a8577SBarry Smith ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY); CHKERRQ(ierr); 10810754003eSLois Curfman McInnes 10820754003eSLois Curfman McInnes /* Free work space */ 10830452661fSBarry Smith PetscFree(smap); PetscFree(cwork); PetscFree(vwork); 108478b31e54SBarry Smith ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr); 108578b31e54SBarry Smith ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr); 10860754003eSLois Curfman McInnes *submat = newmat; 1087*3a40ed3dSBarry Smith PetscFunctionReturn(0); 10880754003eSLois Curfman McInnes } 10890754003eSLois Curfman McInnes 10905615d1e5SSatish Balay #undef __FUNC__ 10915615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrices_SeqDense" 10928f6be9afSLois Curfman McInnes int MatGetSubMatrices_SeqDense(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall, 1093905e6a2fSBarry Smith Mat **B) 1094905e6a2fSBarry Smith { 1095905e6a2fSBarry Smith int ierr,i; 1096905e6a2fSBarry Smith 1097*3a40ed3dSBarry Smith PetscFunctionBegin; 1098905e6a2fSBarry Smith if (scall == MAT_INITIAL_MATRIX) { 1099905e6a2fSBarry Smith *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B); 1100905e6a2fSBarry Smith } 1101905e6a2fSBarry Smith 1102905e6a2fSBarry Smith for ( i=0; i<n; i++ ) { 1103905e6a2fSBarry Smith ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],scall,&(*B)[i]);CHKERRQ(ierr); 1104905e6a2fSBarry Smith } 1105*3a40ed3dSBarry Smith PetscFunctionReturn(0); 1106905e6a2fSBarry Smith } 1107905e6a2fSBarry Smith 11085615d1e5SSatish Balay #undef __FUNC__ 11095615d1e5SSatish Balay #define __FUNC__ "MatCopy_SeqDense" 11108f6be9afSLois Curfman McInnes int MatCopy_SeqDense(Mat A, Mat B) 11114b0e389bSBarry Smith { 11124b0e389bSBarry Smith Mat_SeqDense *a = (Mat_SeqDense *) A->data, *b = (Mat_SeqDense *)B->data; 1113*3a40ed3dSBarry Smith int ierr; 1114*3a40ed3dSBarry Smith 1115*3a40ed3dSBarry Smith PetscFunctionBegin; 1116*3a40ed3dSBarry Smith if (B->type != MATSEQDENSE) { 1117*3a40ed3dSBarry Smith ierr = MatCopy_Basic(A,B);CHKERRQ(ierr); 1118*3a40ed3dSBarry Smith PetscFunctionReturn(0); 1119*3a40ed3dSBarry Smith } 1120e3372554SBarry Smith if (a->m != b->m || a->n != b->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"size(B) != size(A)"); 11214b0e389bSBarry Smith PetscMemcpy(b->v,a->v,a->m*a->n*sizeof(Scalar)); 1122*3a40ed3dSBarry Smith PetscFunctionReturn(0); 11234b0e389bSBarry Smith } 11244b0e389bSBarry Smith 1125289bc588SBarry Smith /* -------------------------------------------------------------------*/ 1126ae80bb75SLois Curfman McInnes static struct _MatOps MatOps = {MatSetValues_SeqDense, 1127905e6a2fSBarry Smith MatGetRow_SeqDense, 1128905e6a2fSBarry Smith MatRestoreRow_SeqDense, 1129905e6a2fSBarry Smith MatMult_SeqDense, 1130905e6a2fSBarry Smith MatMultAdd_SeqDense, 1131905e6a2fSBarry Smith MatMultTrans_SeqDense, 1132905e6a2fSBarry Smith MatMultTransAdd_SeqDense, 1133905e6a2fSBarry Smith MatSolve_SeqDense, 1134905e6a2fSBarry Smith MatSolveAdd_SeqDense, 1135905e6a2fSBarry Smith MatSolveTrans_SeqDense, 1136905e6a2fSBarry Smith MatSolveTransAdd_SeqDense, 1137905e6a2fSBarry Smith MatLUFactor_SeqDense, 1138905e6a2fSBarry Smith MatCholeskyFactor_SeqDense, 1139ec8511deSBarry Smith MatRelax_SeqDense, 1140ec8511deSBarry Smith MatTranspose_SeqDense, 1141905e6a2fSBarry Smith MatGetInfo_SeqDense, 1142905e6a2fSBarry Smith MatEqual_SeqDense, 1143905e6a2fSBarry Smith MatGetDiagonal_SeqDense, 1144905e6a2fSBarry Smith MatDiagonalScale_SeqDense, 1145905e6a2fSBarry Smith MatNorm_SeqDense, 1146905e6a2fSBarry Smith 0, 1147905e6a2fSBarry Smith 0, 1148905e6a2fSBarry Smith 0, 1149905e6a2fSBarry Smith MatSetOption_SeqDense, 1150905e6a2fSBarry Smith MatZeroEntries_SeqDense, 1151905e6a2fSBarry Smith MatZeroRows_SeqDense, 1152905e6a2fSBarry Smith MatLUFactorSymbolic_SeqDense, 1153905e6a2fSBarry Smith MatLUFactorNumeric_SeqDense, 1154905e6a2fSBarry Smith MatCholeskyFactorSymbolic_SeqDense, 1155905e6a2fSBarry Smith MatCholeskyFactorNumeric_SeqDense, 1156905e6a2fSBarry Smith MatGetSize_SeqDense, 1157905e6a2fSBarry Smith MatGetSize_SeqDense, 1158905e6a2fSBarry Smith MatGetOwnershipRange_SeqDense, 1159905e6a2fSBarry Smith 0, 1160905e6a2fSBarry Smith 0, 1161905e6a2fSBarry Smith MatGetArray_SeqDense, 1162905e6a2fSBarry Smith MatRestoreArray_SeqDense, 11633d1612f7SBarry Smith MatConvertSameType_SeqDense,0,0,0,0, 1164905e6a2fSBarry Smith MatAXPY_SeqDense,MatGetSubMatrices_SeqDense,0, 11654b0e389bSBarry Smith MatGetValues_SeqDense, 11664e220ebcSLois Curfman McInnes MatCopy_SeqDense,0,MatScale_SeqDense, 11674e220ebcSLois Curfman McInnes 0,0,0,MatGetBlockSize_SeqDense}; 116890ace30eSBarry Smith 11695615d1e5SSatish Balay #undef __FUNC__ 11705615d1e5SSatish Balay #define __FUNC__ "MatCreateSeqDense" 11714b828684SBarry Smith /*@C 1172fafbff53SBarry Smith MatCreateSeqDense - Creates a sequential dense matrix that 1173d65003e9SLois Curfman McInnes is stored in column major order (the usual Fortran 77 manner). Many 1174d65003e9SLois Curfman McInnes of the matrix operations use the BLAS and LAPACK routines. 1175289bc588SBarry Smith 117620563c6bSBarry Smith Input Parameters: 1177029af93fSBarry Smith . comm - MPI communicator, set to PETSC_COMM_SELF 11780c775827SLois Curfman McInnes . m - number of rows 117918f449edSLois Curfman McInnes . n - number of columns 1180b4fd4287SBarry Smith . data - optional location of matrix data. Set data=PETSC_NULL for PETSc 1181dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 118220563c6bSBarry Smith 118320563c6bSBarry Smith Output Parameter: 118444cd7ae7SLois Curfman McInnes . A - the matrix 118520563c6bSBarry Smith 118618f449edSLois Curfman McInnes Notes: 118718f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 118818f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 1189b4fd4287SBarry Smith set data=PETSC_NULL. 119018f449edSLois Curfman McInnes 1191dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS 1192d65003e9SLois Curfman McInnes 1193dbd7a890SLois Curfman McInnes .seealso: MatCreate(), MatSetValues() 119420563c6bSBarry Smith @*/ 119544cd7ae7SLois Curfman McInnes int MatCreateSeqDense(MPI_Comm comm,int m,int n,Scalar *data,Mat *A) 1196289bc588SBarry Smith { 119744cd7ae7SLois Curfman McInnes Mat B; 119844cd7ae7SLois Curfman McInnes Mat_SeqDense *b; 11993b2fbd54SBarry Smith int ierr,flg,size; 12003b2fbd54SBarry Smith 1201*3a40ed3dSBarry Smith PetscFunctionBegin; 12023b2fbd54SBarry Smith MPI_Comm_size(comm,&size); 1203e3372554SBarry Smith if (size > 1) SETERRQ(1,0,"Comm must be of size 1"); 120455659b69SBarry Smith 120544cd7ae7SLois Curfman McInnes *A = 0; 1206d4bb536fSBarry Smith PetscHeaderCreate(B,_p_Mat,MAT_COOKIE,MATSEQDENSE,comm,MatDestroy,MatView); 120744cd7ae7SLois Curfman McInnes PLogObjectCreate(B); 120844cd7ae7SLois Curfman McInnes b = (Mat_SeqDense *) PetscMalloc(sizeof(Mat_SeqDense)); CHKPTRQ(b); 120944cd7ae7SLois Curfman McInnes PetscMemzero(b,sizeof(Mat_SeqDense)); 121044cd7ae7SLois Curfman McInnes PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps)); 121144cd7ae7SLois Curfman McInnes B->destroy = MatDestroy_SeqDense; 121244cd7ae7SLois Curfman McInnes B->view = MatView_SeqDense; 121344cd7ae7SLois Curfman McInnes B->factor = 0; 121490f02eecSBarry Smith B->mapping = 0; 1215f09e8eb9SSatish Balay PLogObjectMemory(B,sizeof(struct _p_Mat)); 121644cd7ae7SLois Curfman McInnes B->data = (void *) b; 121718f449edSLois Curfman McInnes 121844cd7ae7SLois Curfman McInnes b->m = m; B->m = m; B->M = m; 121944cd7ae7SLois Curfman McInnes b->n = n; B->n = n; B->N = n; 122044cd7ae7SLois Curfman McInnes b->pivots = 0; 122144cd7ae7SLois Curfman McInnes b->roworiented = 1; 1222b4fd4287SBarry Smith if (data == PETSC_NULL) { 122344cd7ae7SLois Curfman McInnes b->v = (Scalar*) PetscMalloc((m*n+1)*sizeof(Scalar)); CHKPTRQ(b->v); 122444cd7ae7SLois Curfman McInnes PetscMemzero(b->v,m*n*sizeof(Scalar)); 122544cd7ae7SLois Curfman McInnes b->user_alloc = 0; 122644cd7ae7SLois Curfman McInnes PLogObjectMemory(B,n*m*sizeof(Scalar)); 122718f449edSLois Curfman McInnes } 12282dd2b441SLois Curfman McInnes else { /* user-allocated storage */ 122944cd7ae7SLois Curfman McInnes b->v = data; 123044cd7ae7SLois Curfman McInnes b->user_alloc = 1; 12312dd2b441SLois Curfman McInnes } 123225cdf11fSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr); 123344cd7ae7SLois Curfman McInnes if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr);} 12344e220ebcSLois Curfman McInnes 123544cd7ae7SLois Curfman McInnes *A = B; 1236*3a40ed3dSBarry Smith PetscFunctionReturn(0); 1237289bc588SBarry Smith } 1238289bc588SBarry Smith 12393b2fbd54SBarry Smith 12403b2fbd54SBarry Smith 12413b2fbd54SBarry Smith 12423b2fbd54SBarry Smith 12433b2fbd54SBarry Smith 1244