1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*d64ed03dSBarry Smith static char vcid[] = "$Id: dense.c,v 1.132 1997/10/19 03:25:08 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; 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); 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); 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 } 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); 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 } 219e3372554SBarry Smith if (info) SETERRQ(1,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 ); 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); 2523a40ed3dSBarry 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 2643a40ed3dSBarry 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 ); 2753a40ed3dSBarry 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); 2823a40ed3dSBarry Smith } else { 28390f02eecSBarry Smith ierr = VecAXPY(&sone,zz,yy); CHKERRQ(ierr); 28490f02eecSBarry Smith } 28555659b69SBarry Smith PLogFlops(mat->n*mat->n - mat->n); 2863a40ed3dSBarry 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 2983a40ed3dSBarry Smith PetscFunctionBegin; 299289bc588SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 3003a40ed3dSBarry 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++ ) { 3073a40ed3dSBarry 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-- ) { 3243a40ed3dSBarry 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 } 3403a40ed3dSBarry 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; 3513a40ed3dSBarry Smith 3523a40ed3dSBarry 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); 3563a40ed3dSBarry 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; 3663a40ed3dSBarry Smith 3673a40ed3dSBarry 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); 3713a40ed3dSBarry 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; 3813a40ed3dSBarry Smith 3823a40ed3dSBarry 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); 3873a40ed3dSBarry 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; 3983a40ed3dSBarry Smith 3993a40ed3dSBarry 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); 4043a40ed3dSBarry 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 4163a40ed3dSBarry 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 } 4273a40ed3dSBarry 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); } 4363a40ed3dSBarry 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 4473a40ed3dSBarry Smith PetscFunctionBegin; 448289bc588SBarry Smith if (!mat->roworiented) { 449dbb450caSBarry Smith if (addv == INSERT_VALUES) { 450289bc588SBarry Smith for ( j=0; j<n; j++ ) { 4513a40ed3dSBarry 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++ ) { 4563a40ed3dSBarry 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 } 4633a40ed3dSBarry Smith } else { 464289bc588SBarry Smith for ( j=0; j<n; j++ ) { 4653a40ed3dSBarry 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++ ) { 4703a40ed3dSBarry 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 } 4783a40ed3dSBarry Smith } else { 479dbb450caSBarry Smith if (addv == INSERT_VALUES) { 480e8d4e0b9SBarry Smith for ( i=0; i<m; i++ ) { 4813a40ed3dSBarry 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++ ) { 4863a40ed3dSBarry 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 } 4933a40ed3dSBarry Smith } else { 494289bc588SBarry Smith for ( i=0; i<m; i++ ) { 4953a40ed3dSBarry 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++ ) { 5003a40ed3dSBarry 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 } 5093a40ed3dSBarry 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 5203a40ed3dSBarry 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 } 5273a40ed3dSBarry 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 5453a40ed3dSBarry 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); 563*d64ed03dSBarry 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 } 5913a40ed3dSBarry 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 6073a40ed3dSBarry 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) { 6123a40ed3dSBarry 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++ ) { 6193a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 620*d64ed03dSBarry Smith if (real(*v) != 0.0 && imag(*v) != 0.0) fprintf(fd," %d %g + %g i",j,real(*v),imag(*v)); 62180cd9d93SLois Curfman McInnes else if (real(*v)) fprintf(fd," %d %g ",j,real(*v)); 62280cd9d93SLois Curfman McInnes #else 62380cd9d93SLois Curfman McInnes if (*v) fprintf(fd," %d %g ",j,*v); 62480cd9d93SLois Curfman McInnes #endif 625*d64ed03dSBarry Smith v += a->m; 62680cd9d93SLois Curfman McInnes } 62780cd9d93SLois Curfman McInnes fprintf(fd,"\n"); 62880cd9d93SLois Curfman McInnes } 6293a40ed3dSBarry Smith } else { 6303a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 63147989497SBarry Smith int allreal = 1; 63247989497SBarry Smith /* determine if matrix has all real values */ 63347989497SBarry Smith v = a->v; 63447989497SBarry Smith for ( i=0; i<a->m*a->n; i++ ) { 63547989497SBarry Smith if (imag(v[i])) { allreal = 0; break ;} 63647989497SBarry Smith } 63747989497SBarry Smith #endif 638932b0c3eSLois Curfman McInnes for ( i=0; i<a->m; i++ ) { 639932b0c3eSLois Curfman McInnes v = a->v + i; 640932b0c3eSLois Curfman McInnes for ( j=0; j<a->n; j++ ) { 6413a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 64247989497SBarry Smith if (allreal) { 64347989497SBarry Smith fprintf(fd,"%6.4e ",real(*v)); v += a->m; 64447989497SBarry Smith } else { 645932b0c3eSLois Curfman McInnes fprintf(fd,"%6.4e + %6.4e i ",real(*v),imag(*v)); v += a->m; 64647989497SBarry Smith } 647289bc588SBarry Smith #else 648932b0c3eSLois Curfman McInnes fprintf(fd,"%6.4e ",*v); v += a->m; 649289bc588SBarry Smith #endif 650289bc588SBarry Smith } 6518e37dc09SBarry Smith fprintf(fd,"\n"); 652289bc588SBarry Smith } 653da3a660dSBarry Smith } 654932b0c3eSLois Curfman McInnes fflush(fd); 6553a40ed3dSBarry Smith PetscFunctionReturn(0); 656289bc588SBarry Smith } 657289bc588SBarry Smith 6585615d1e5SSatish Balay #undef __FUNC__ 6595615d1e5SSatish Balay #define __FUNC__ "MatView_SeqDense_Binary" 660932b0c3eSLois Curfman McInnes static int MatView_SeqDense_Binary(Mat A,Viewer viewer) 661932b0c3eSLois Curfman McInnes { 662932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense *) A->data; 663932b0c3eSLois Curfman McInnes int ict, j, n = a->n, m = a->m, i, fd, *col_lens, ierr, nz = m*n; 66490ace30eSBarry Smith int format; 66590ace30eSBarry Smith Scalar *v, *anonz,*vals; 666932b0c3eSLois Curfman McInnes 6673a40ed3dSBarry Smith PetscFunctionBegin; 66890ace30eSBarry Smith ierr = ViewerBinaryGetDescriptor(viewer,&fd); CHKERRQ(ierr); 66990ace30eSBarry Smith 67090ace30eSBarry Smith ierr = ViewerGetFormat(viewer,&format); CHKERRQ(ierr); 67102cad45dSBarry Smith if (format == VIEWER_FORMAT_BINARY_NATIVE) { 67290ace30eSBarry Smith /* store the matrix as a dense matrix */ 67390ace30eSBarry Smith col_lens = (int *) PetscMalloc( 4*sizeof(int) ); CHKPTRQ(col_lens); 67490ace30eSBarry Smith col_lens[0] = MAT_COOKIE; 67590ace30eSBarry Smith col_lens[1] = m; 67690ace30eSBarry Smith col_lens[2] = n; 67790ace30eSBarry Smith col_lens[3] = MATRIX_BINARY_FORMAT_DENSE; 6780752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,1); CHKERRQ(ierr); 67990ace30eSBarry Smith PetscFree(col_lens); 68090ace30eSBarry Smith 68190ace30eSBarry Smith /* write out matrix, by rows */ 68245d48df9SBarry Smith vals = (Scalar *) PetscMalloc((m*n+1)*sizeof(Scalar)); CHKPTRQ(vals); 68390ace30eSBarry Smith v = a->v; 68490ace30eSBarry Smith for ( i=0; i<m; i++ ) { 68590ace30eSBarry Smith for ( j=0; j<n; j++ ) { 68690ace30eSBarry Smith vals[i + j*m] = *v++; 68790ace30eSBarry Smith } 68890ace30eSBarry Smith } 6890752156aSBarry Smith ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,0); CHKERRQ(ierr); 69090ace30eSBarry Smith PetscFree(vals); 69190ace30eSBarry Smith } else { 6920452661fSBarry Smith col_lens = (int *) PetscMalloc( (4+nz)*sizeof(int) ); CHKPTRQ(col_lens); 693932b0c3eSLois Curfman McInnes col_lens[0] = MAT_COOKIE; 694932b0c3eSLois Curfman McInnes col_lens[1] = m; 695932b0c3eSLois Curfman McInnes col_lens[2] = n; 696932b0c3eSLois Curfman McInnes col_lens[3] = nz; 697932b0c3eSLois Curfman McInnes 698932b0c3eSLois Curfman McInnes /* store lengths of each row and write (including header) to file */ 699932b0c3eSLois Curfman McInnes for ( i=0; i<m; i++ ) col_lens[4+i] = n; 7000752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,1); CHKERRQ(ierr); 701932b0c3eSLois Curfman McInnes 702932b0c3eSLois Curfman McInnes /* Possibly should write in smaller increments, not whole matrix at once? */ 703932b0c3eSLois Curfman McInnes /* store column indices (zero start index) */ 704932b0c3eSLois Curfman McInnes ict = 0; 705932b0c3eSLois Curfman McInnes for ( i=0; i<m; i++ ) { 706932b0c3eSLois Curfman McInnes for ( j=0; j<n; j++ ) col_lens[ict++] = j; 707932b0c3eSLois Curfman McInnes } 7080752156aSBarry Smith ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,0); CHKERRQ(ierr); 7090452661fSBarry Smith PetscFree(col_lens); 710932b0c3eSLois Curfman McInnes 711932b0c3eSLois Curfman McInnes /* store nonzero values */ 71245d48df9SBarry Smith anonz = (Scalar *) PetscMalloc((nz+1)*sizeof(Scalar)); CHKPTRQ(anonz); 713932b0c3eSLois Curfman McInnes ict = 0; 714932b0c3eSLois Curfman McInnes for ( i=0; i<m; i++ ) { 715932b0c3eSLois Curfman McInnes v = a->v + i; 716932b0c3eSLois Curfman McInnes for ( j=0; j<n; j++ ) { 717932b0c3eSLois Curfman McInnes anonz[ict++] = *v; v += a->m; 718932b0c3eSLois Curfman McInnes } 719932b0c3eSLois Curfman McInnes } 7200752156aSBarry Smith ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,0); CHKERRQ(ierr); 7210452661fSBarry Smith PetscFree(anonz); 72290ace30eSBarry Smith } 7233a40ed3dSBarry Smith PetscFunctionReturn(0); 724932b0c3eSLois Curfman McInnes } 725932b0c3eSLois Curfman McInnes 7265615d1e5SSatish Balay #undef __FUNC__ 7275615d1e5SSatish Balay #define __FUNC__ "MatView_SeqDense" 7288f6be9afSLois Curfman McInnes int MatView_SeqDense(PetscObject obj,Viewer viewer) 729932b0c3eSLois Curfman McInnes { 730932b0c3eSLois Curfman McInnes Mat A = (Mat) obj; 731932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*) A->data; 732bcd2baecSBarry Smith ViewerType vtype; 733bcd2baecSBarry Smith int ierr; 734932b0c3eSLois Curfman McInnes 7353a40ed3dSBarry Smith PetscFunctionBegin; 736bcd2baecSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 737bcd2baecSBarry Smith 738bcd2baecSBarry Smith if (vtype == MATLAB_VIEWER) { 7393a40ed3dSBarry Smith ierr = ViewerMatlabPutArray_Private(viewer,a->m,a->n,a->v); CHKERRQ(ierr); 7403a40ed3dSBarry Smith } else if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) { 7413a40ed3dSBarry Smith ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr); 7423a40ed3dSBarry Smith } else if (vtype == BINARY_FILE_VIEWER) { 7433a40ed3dSBarry Smith ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr); 744932b0c3eSLois Curfman McInnes } 7453a40ed3dSBarry Smith PetscFunctionReturn(0); 746932b0c3eSLois Curfman McInnes } 747289bc588SBarry Smith 7485615d1e5SSatish Balay #undef __FUNC__ 7495615d1e5SSatish Balay #define __FUNC__ "MatDestroy_SeqDense" 7508f6be9afSLois Curfman McInnes int MatDestroy_SeqDense(PetscObject obj) 751289bc588SBarry Smith { 752289bc588SBarry Smith Mat mat = (Mat) obj; 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) 758a5a9c739SBarry Smith PLogObjectState(obj,"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 */ 782*d64ed03dSBarry Smith if (m != n) { /* malloc temp to hold transpose */ 783*d64ed03dSBarry Smith Scalar *w = (Scalar *) PetscMalloc((m+1)*(n+1)*sizeof(Scalar));CHKPTRQ(w); 784*d64ed03dSBarry Smith for ( j=0; j<m; j++ ) { 785*d64ed03dSBarry Smith for ( k=0; k<n; k++ ) { 786*d64ed03dSBarry Smith w[k + j*n] = v[j + k*m]; 787*d64ed03dSBarry Smith } 788*d64ed03dSBarry Smith } 789*d64ed03dSBarry Smith PetscMemcpy(v,w,m*n*sizeof(Scalar)); 790*d64ed03dSBarry Smith PetscFree(w); 791*d64ed03dSBarry 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 } 799*d64ed03dSBarry 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 || 9572b362799SSatish Balay op == MAT_IGNORE_OFF_PROC_ENTRIES) 95894a424c1SBarry Smith PLogInfo(A,"Info:MatSetOption_SeqDense:Option ignored\n"); 9593a40ed3dSBarry Smith else { 9603a40ed3dSBarry Smith SETERRQ(PETSC_ERR_SUP,0,"unknown option"); 9613a40ed3dSBarry Smith } 9623a40ed3dSBarry Smith PetscFunctionReturn(0); 963289bc588SBarry Smith } 964289bc588SBarry Smith 9655615d1e5SSatish Balay #undef __FUNC__ 9665615d1e5SSatish Balay #define __FUNC__ "MatZeroEntries_SeqDense" 9678f6be9afSLois Curfman McInnes int MatZeroEntries_SeqDense(Mat A) 9686f0a148fSBarry Smith { 969ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense *) A->data; 9703a40ed3dSBarry Smith 9713a40ed3dSBarry Smith PetscFunctionBegin; 972cddf8d76SBarry Smith PetscMemzero(l->v,l->m*l->n*sizeof(Scalar)); 9733a40ed3dSBarry Smith PetscFunctionReturn(0); 9746f0a148fSBarry Smith } 9756f0a148fSBarry Smith 9765615d1e5SSatish Balay #undef __FUNC__ 9775615d1e5SSatish Balay #define __FUNC__ "MatGetBlockSize_SeqDense" 9788f6be9afSLois Curfman McInnes int MatGetBlockSize_SeqDense(Mat A,int *bs) 9794e220ebcSLois Curfman McInnes { 9803a40ed3dSBarry Smith PetscFunctionBegin; 9814e220ebcSLois Curfman McInnes *bs = 1; 9823a40ed3dSBarry Smith PetscFunctionReturn(0); 9834e220ebcSLois Curfman McInnes } 9844e220ebcSLois Curfman McInnes 9855615d1e5SSatish Balay #undef __FUNC__ 9865615d1e5SSatish Balay #define __FUNC__ "MatZeroRows_SeqDense" 9878f6be9afSLois Curfman McInnes int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag) 9886f0a148fSBarry Smith { 989ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense *) A->data; 9906abc6512SBarry Smith int n = l->n, i, j,ierr,N, *rows; 9916f0a148fSBarry Smith Scalar *slot; 99255659b69SBarry Smith 9933a40ed3dSBarry Smith PetscFunctionBegin; 99477c4ece6SBarry Smith ierr = ISGetSize(is,&N); CHKERRQ(ierr); 99578b31e54SBarry Smith ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 9966f0a148fSBarry Smith for ( i=0; i<N; i++ ) { 9976f0a148fSBarry Smith slot = l->v + rows[i]; 9986f0a148fSBarry Smith for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;} 9996f0a148fSBarry Smith } 10006f0a148fSBarry Smith if (diag) { 10016f0a148fSBarry Smith for ( i=0; i<N; i++ ) { 10026f0a148fSBarry Smith slot = l->v + (n+1)*rows[i]; 1003260925b8SBarry Smith *slot = *diag; 10046f0a148fSBarry Smith } 10056f0a148fSBarry Smith } 1006260925b8SBarry Smith ISRestoreIndices(is,&rows); 10073a40ed3dSBarry Smith PetscFunctionReturn(0); 10086f0a148fSBarry Smith } 1009557bce09SLois Curfman McInnes 10105615d1e5SSatish Balay #undef __FUNC__ 10115615d1e5SSatish Balay #define __FUNC__ "MatGetSize_SeqDense" 10128f6be9afSLois Curfman McInnes int MatGetSize_SeqDense(Mat A,int *m,int *n) 1013557bce09SLois Curfman McInnes { 1014c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 10153a40ed3dSBarry Smith 10163a40ed3dSBarry Smith PetscFunctionBegin; 1017557bce09SLois Curfman McInnes *m = mat->m; *n = mat->n; 10183a40ed3dSBarry Smith PetscFunctionReturn(0); 1019557bce09SLois Curfman McInnes } 1020557bce09SLois Curfman McInnes 10215615d1e5SSatish Balay #undef __FUNC__ 10225615d1e5SSatish Balay #define __FUNC__ "MatGetOwnershipRange_SeqDense" 10238f6be9afSLois Curfman McInnes int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n) 1024d3e5ee88SLois Curfman McInnes { 1025c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 10263a40ed3dSBarry Smith 10273a40ed3dSBarry Smith PetscFunctionBegin; 1028d3e5ee88SLois Curfman McInnes *m = 0; *n = mat->m; 10293a40ed3dSBarry Smith PetscFunctionReturn(0); 1030d3e5ee88SLois Curfman McInnes } 1031d3e5ee88SLois Curfman McInnes 10325615d1e5SSatish Balay #undef __FUNC__ 10335615d1e5SSatish Balay #define __FUNC__ "MatGetArray_SeqDense" 10348f6be9afSLois Curfman McInnes int MatGetArray_SeqDense(Mat A,Scalar **array) 103564e87e97SBarry Smith { 1036c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 10373a40ed3dSBarry Smith 10383a40ed3dSBarry Smith PetscFunctionBegin; 103964e87e97SBarry Smith *array = mat->v; 10403a40ed3dSBarry Smith PetscFunctionReturn(0); 104164e87e97SBarry Smith } 10420754003eSLois Curfman McInnes 10435615d1e5SSatish Balay #undef __FUNC__ 10445615d1e5SSatish Balay #define __FUNC__ "MatRestoreArray_SeqDense" 10458f6be9afSLois Curfman McInnes int MatRestoreArray_SeqDense(Mat A,Scalar **array) 1046ff14e315SSatish Balay { 10473a40ed3dSBarry Smith PetscFunctionBegin; 10483a40ed3dSBarry Smith PetscFunctionReturn(0); 1049ff14e315SSatish Balay } 10500754003eSLois Curfman McInnes 10515615d1e5SSatish Balay #undef __FUNC__ 10525615d1e5SSatish Balay #define __FUNC__ "MatGetSubMatrix_SeqDense" 1053cddf8d76SBarry Smith static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall, 1054cddf8d76SBarry Smith 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" 11008f6be9afSLois Curfman McInnes int MatGetSubMatrices_SeqDense(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall, 1101905e6a2fSBarry Smith Mat **B) 1102905e6a2fSBarry Smith { 1103905e6a2fSBarry Smith int ierr,i; 1104905e6a2fSBarry Smith 11053a40ed3dSBarry Smith PetscFunctionBegin; 1106905e6a2fSBarry Smith if (scall == MAT_INITIAL_MATRIX) { 1107905e6a2fSBarry Smith *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B); 1108905e6a2fSBarry Smith } 1109905e6a2fSBarry Smith 1110905e6a2fSBarry Smith for ( i=0; i<n; i++ ) { 1111905e6a2fSBarry Smith ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],scall,&(*B)[i]);CHKERRQ(ierr); 1112905e6a2fSBarry Smith } 11133a40ed3dSBarry Smith PetscFunctionReturn(0); 1114905e6a2fSBarry Smith } 1115905e6a2fSBarry Smith 11165615d1e5SSatish Balay #undef __FUNC__ 11175615d1e5SSatish Balay #define __FUNC__ "MatCopy_SeqDense" 11188f6be9afSLois Curfman McInnes int MatCopy_SeqDense(Mat A, Mat B) 11194b0e389bSBarry Smith { 11204b0e389bSBarry Smith Mat_SeqDense *a = (Mat_SeqDense *) A->data, *b = (Mat_SeqDense *)B->data; 11213a40ed3dSBarry Smith int ierr; 11223a40ed3dSBarry Smith 11233a40ed3dSBarry Smith PetscFunctionBegin; 11243a40ed3dSBarry Smith if (B->type != MATSEQDENSE) { 11253a40ed3dSBarry Smith ierr = MatCopy_Basic(A,B);CHKERRQ(ierr); 11263a40ed3dSBarry Smith PetscFunctionReturn(0); 11273a40ed3dSBarry Smith } 1128e3372554SBarry Smith if (a->m != b->m || a->n != b->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"size(B) != size(A)"); 11294b0e389bSBarry Smith PetscMemcpy(b->v,a->v,a->m*a->n*sizeof(Scalar)); 11303a40ed3dSBarry Smith PetscFunctionReturn(0); 11314b0e389bSBarry Smith } 11324b0e389bSBarry Smith 1133289bc588SBarry Smith /* -------------------------------------------------------------------*/ 1134ae80bb75SLois Curfman McInnes static struct _MatOps MatOps = {MatSetValues_SeqDense, 1135905e6a2fSBarry Smith MatGetRow_SeqDense, 1136905e6a2fSBarry Smith MatRestoreRow_SeqDense, 1137905e6a2fSBarry Smith MatMult_SeqDense, 1138905e6a2fSBarry Smith MatMultAdd_SeqDense, 1139905e6a2fSBarry Smith MatMultTrans_SeqDense, 1140905e6a2fSBarry Smith MatMultTransAdd_SeqDense, 1141905e6a2fSBarry Smith MatSolve_SeqDense, 1142905e6a2fSBarry Smith MatSolveAdd_SeqDense, 1143905e6a2fSBarry Smith MatSolveTrans_SeqDense, 1144905e6a2fSBarry Smith MatSolveTransAdd_SeqDense, 1145905e6a2fSBarry Smith MatLUFactor_SeqDense, 1146905e6a2fSBarry Smith MatCholeskyFactor_SeqDense, 1147ec8511deSBarry Smith MatRelax_SeqDense, 1148ec8511deSBarry Smith MatTranspose_SeqDense, 1149905e6a2fSBarry Smith MatGetInfo_SeqDense, 1150905e6a2fSBarry Smith MatEqual_SeqDense, 1151905e6a2fSBarry Smith MatGetDiagonal_SeqDense, 1152905e6a2fSBarry Smith MatDiagonalScale_SeqDense, 1153905e6a2fSBarry Smith MatNorm_SeqDense, 1154905e6a2fSBarry Smith 0, 1155905e6a2fSBarry Smith 0, 1156905e6a2fSBarry Smith 0, 1157905e6a2fSBarry Smith MatSetOption_SeqDense, 1158905e6a2fSBarry Smith MatZeroEntries_SeqDense, 1159905e6a2fSBarry Smith MatZeroRows_SeqDense, 1160905e6a2fSBarry Smith MatLUFactorSymbolic_SeqDense, 1161905e6a2fSBarry Smith MatLUFactorNumeric_SeqDense, 1162905e6a2fSBarry Smith MatCholeskyFactorSymbolic_SeqDense, 1163905e6a2fSBarry Smith MatCholeskyFactorNumeric_SeqDense, 1164905e6a2fSBarry Smith MatGetSize_SeqDense, 1165905e6a2fSBarry Smith MatGetSize_SeqDense, 1166905e6a2fSBarry Smith MatGetOwnershipRange_SeqDense, 1167905e6a2fSBarry Smith 0, 1168905e6a2fSBarry Smith 0, 1169905e6a2fSBarry Smith MatGetArray_SeqDense, 1170905e6a2fSBarry Smith MatRestoreArray_SeqDense, 11713d1612f7SBarry Smith MatConvertSameType_SeqDense,0,0,0,0, 1172905e6a2fSBarry Smith MatAXPY_SeqDense,MatGetSubMatrices_SeqDense,0, 11734b0e389bSBarry Smith MatGetValues_SeqDense, 11744e220ebcSLois Curfman McInnes MatCopy_SeqDense,0,MatScale_SeqDense, 11754e220ebcSLois Curfman McInnes 0,0,0,MatGetBlockSize_SeqDense}; 117690ace30eSBarry Smith 11775615d1e5SSatish Balay #undef __FUNC__ 11785615d1e5SSatish Balay #define __FUNC__ "MatCreateSeqDense" 11794b828684SBarry Smith /*@C 1180fafbff53SBarry Smith MatCreateSeqDense - Creates a sequential dense matrix that 1181d65003e9SLois Curfman McInnes is stored in column major order (the usual Fortran 77 manner). Many 1182d65003e9SLois Curfman McInnes of the matrix operations use the BLAS and LAPACK routines. 1183289bc588SBarry Smith 118420563c6bSBarry Smith Input Parameters: 1185029af93fSBarry Smith . comm - MPI communicator, set to PETSC_COMM_SELF 11860c775827SLois Curfman McInnes . m - number of rows 118718f449edSLois Curfman McInnes . n - number of columns 1188b4fd4287SBarry Smith . data - optional location of matrix data. Set data=PETSC_NULL for PETSc 1189dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 119020563c6bSBarry Smith 119120563c6bSBarry Smith Output Parameter: 119244cd7ae7SLois Curfman McInnes . A - the matrix 119320563c6bSBarry Smith 119418f449edSLois Curfman McInnes Notes: 119518f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 119618f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 1197b4fd4287SBarry Smith set data=PETSC_NULL. 119818f449edSLois Curfman McInnes 1199dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS 1200d65003e9SLois Curfman McInnes 1201dbd7a890SLois Curfman McInnes .seealso: MatCreate(), MatSetValues() 120220563c6bSBarry Smith @*/ 120344cd7ae7SLois Curfman McInnes int MatCreateSeqDense(MPI_Comm comm,int m,int n,Scalar *data,Mat *A) 1204289bc588SBarry Smith { 120544cd7ae7SLois Curfman McInnes Mat B; 120644cd7ae7SLois Curfman McInnes Mat_SeqDense *b; 12073b2fbd54SBarry Smith int ierr,flg,size; 12083b2fbd54SBarry Smith 12093a40ed3dSBarry Smith PetscFunctionBegin; 12103b2fbd54SBarry Smith MPI_Comm_size(comm,&size); 1211e3372554SBarry Smith if (size > 1) SETERRQ(1,0,"Comm must be of size 1"); 121255659b69SBarry Smith 121344cd7ae7SLois Curfman McInnes *A = 0; 1214d4bb536fSBarry Smith PetscHeaderCreate(B,_p_Mat,MAT_COOKIE,MATSEQDENSE,comm,MatDestroy,MatView); 121544cd7ae7SLois Curfman McInnes PLogObjectCreate(B); 121644cd7ae7SLois Curfman McInnes b = (Mat_SeqDense *) PetscMalloc(sizeof(Mat_SeqDense)); CHKPTRQ(b); 121744cd7ae7SLois Curfman McInnes PetscMemzero(b,sizeof(Mat_SeqDense)); 121844cd7ae7SLois Curfman McInnes PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps)); 121944cd7ae7SLois Curfman McInnes B->destroy = MatDestroy_SeqDense; 122044cd7ae7SLois Curfman McInnes B->view = MatView_SeqDense; 122144cd7ae7SLois Curfman McInnes B->factor = 0; 122290f02eecSBarry Smith B->mapping = 0; 1223f09e8eb9SSatish Balay PLogObjectMemory(B,sizeof(struct _p_Mat)); 122444cd7ae7SLois Curfman McInnes B->data = (void *) b; 122518f449edSLois Curfman McInnes 122644cd7ae7SLois Curfman McInnes b->m = m; B->m = m; B->M = m; 122744cd7ae7SLois Curfman McInnes b->n = n; B->n = n; B->N = n; 122844cd7ae7SLois Curfman McInnes b->pivots = 0; 122944cd7ae7SLois Curfman McInnes b->roworiented = 1; 1230b4fd4287SBarry Smith if (data == PETSC_NULL) { 123144cd7ae7SLois Curfman McInnes b->v = (Scalar*) PetscMalloc((m*n+1)*sizeof(Scalar)); CHKPTRQ(b->v); 123244cd7ae7SLois Curfman McInnes PetscMemzero(b->v,m*n*sizeof(Scalar)); 123344cd7ae7SLois Curfman McInnes b->user_alloc = 0; 123444cd7ae7SLois Curfman McInnes PLogObjectMemory(B,n*m*sizeof(Scalar)); 123518f449edSLois Curfman McInnes } 12362dd2b441SLois Curfman McInnes else { /* user-allocated storage */ 123744cd7ae7SLois Curfman McInnes b->v = data; 123844cd7ae7SLois Curfman McInnes b->user_alloc = 1; 12392dd2b441SLois Curfman McInnes } 124025cdf11fSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr); 124144cd7ae7SLois Curfman McInnes if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr);} 12424e220ebcSLois Curfman McInnes 124344cd7ae7SLois Curfman McInnes *A = B; 12443a40ed3dSBarry Smith PetscFunctionReturn(0); 1245289bc588SBarry Smith } 1246289bc588SBarry Smith 12473b2fbd54SBarry Smith 12483b2fbd54SBarry Smith 12493b2fbd54SBarry Smith 12503b2fbd54SBarry Smith 12513b2fbd54SBarry Smith 1252