1be1d678aSKris Buschelman #define PETSCMAT_DLL 20925cdddSBarry Smith 37c4f633dSBarry Smith #include "private/matimpl.h" /*I "petscmat.h" I*/ 40925cdddSBarry Smith 54a2ae208SSatish Balay #undef __FUNCT__ 64a2ae208SSatish Balay #define __FUNCT__ "MatGetColumnVector" 70925cdddSBarry Smith /*@ 882bf6240SBarry Smith MatGetColumnVector - Gets the values from a given column of a matrix. 90925cdddSBarry Smith 1072257631SBarry Smith Not Collective 11fee21e36SBarry Smith 1298a79cdbSBarry Smith Input Parameters: 135ed6d96aSBarry Smith + A - the matrix 145ed6d96aSBarry Smith . yy - the vector 155ed6d96aSBarry Smith - c - the column requested (in global numbering) 1698a79cdbSBarry Smith 1715091d37SBarry Smith Level: advanced 1815091d37SBarry Smith 199d006df2SBarry Smith Notes: 209d006df2SBarry Smith Each processor for which this is called gets the values for its rows. 219d006df2SBarry Smith 229d006df2SBarry Smith Since PETSc matrices are usually stored in compressed row format, this routine 239d006df2SBarry Smith will generally be slow. 249d006df2SBarry Smith 253d81755aSBarry Smith The vector must have the same parallel row layout as the matrix. 263d81755aSBarry Smith 2782bf6240SBarry Smith Contributed by: Denis Vanderstraeten 280925cdddSBarry Smith 2982bf6240SBarry Smith .keywords: matrix, column, get 300925cdddSBarry Smith 3115091d37SBarry Smith .seealso: MatGetRow(), MatGetDiagonal() 3215091d37SBarry Smith 330925cdddSBarry Smith @*/ 34be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnVector(Mat A,Vec yy,PetscInt col) 350925cdddSBarry Smith { 368aa348c1SBarry Smith PetscScalar *y; 37b3cc6726SBarry Smith const PetscScalar *v; 38dfbe8321SBarry Smith PetscErrorCode ierr; 3938baddfdSBarry Smith PetscInt i,j,nz,N,Rs,Re,rs,re; 4038baddfdSBarry Smith const PetscInt *idx; 410925cdddSBarry Smith 420925cdddSBarry Smith PetscFunctionBegin; 434482741eSBarry Smith PetscValidHeaderSpecific(A,MAT_COOKIE,1); 444482741eSBarry Smith PetscValidHeaderSpecific(yy,VEC_COOKIE,2); 4577431f27SBarry Smith if (col < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Requested negative column: %D",col); 46011b8408SBarry Smith ierr = MatGetSize(A,PETSC_NULL,&N);CHKERRQ(ierr); 4777431f27SBarry Smith if (col >= N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Requested column %D larger than number columns in matrix %D",col,N); 4882bf6240SBarry Smith ierr = MatGetOwnershipRange(A,&Rs,&Re);CHKERRQ(ierr); 4982bf6240SBarry Smith ierr = VecGetOwnershipRange(yy,&rs,&re);CHKERRQ(ierr); 5077431f27SBarry Smith if (Rs != rs || Re != re) SETERRQ4(PETSC_ERR_ARG_INCOMP,"Matrix %D %D does not have same ownership range (size) as vector %D %D",Rs,Re,rs,re); 5182bf6240SBarry Smith 528d0534beSBarry Smith if (A->ops->getcolumnvector) { 538d0534beSBarry Smith ierr = (*A->ops->getcolumnvector)(A,yy,col);CHKERRQ(ierr); 548d0534beSBarry Smith } else { 558aa348c1SBarry Smith ierr = VecSet(yy,0.0);CHKERRQ(ierr); 5682bf6240SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 5782bf6240SBarry Smith 5882bf6240SBarry Smith for (i=Rs; i<Re; i++) { 5982bf6240SBarry Smith ierr = MatGetRow(A,i,&nz,&idx,&v);CHKERRQ(ierr); 6082bf6240SBarry Smith if (nz && idx[0] <= col) { 6182bf6240SBarry Smith /* 6282bf6240SBarry Smith Should use faster search here 6382bf6240SBarry Smith */ 6482bf6240SBarry Smith for (j=0; j<nz; j++) { 6582bf6240SBarry Smith if (idx[j] >= col) { 6682bf6240SBarry Smith if (idx[j] == col) y[i-rs] = v[j]; 6782bf6240SBarry Smith break; 680925cdddSBarry Smith } 690925cdddSBarry Smith } 700925cdddSBarry Smith } 7182bf6240SBarry Smith ierr = MatRestoreRow(A,i,&nz,&idx,&v);CHKERRQ(ierr); 720925cdddSBarry Smith } 7382bf6240SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 748d0534beSBarry Smith } 750925cdddSBarry Smith PetscFunctionReturn(0); 760925cdddSBarry Smith } 77*242f1d38SBarry Smith 78*242f1d38SBarry Smith #include "../src/mat/impls/aij/seq/aij.h" 79*242f1d38SBarry Smith 80*242f1d38SBarry Smith #undef __FUNCT__ 81*242f1d38SBarry Smith #define __FUNCT__ "MatGetColumnNorms_SeqAIJ" 82*242f1d38SBarry Smith PetscErrorCode MatGetColumnNorms_SeqAIJ(Mat A,NormType type,PetscReal *norms) 83*242f1d38SBarry Smith { 84*242f1d38SBarry Smith PetscErrorCode ierr; 85*242f1d38SBarry Smith PetscInt i,m,n; 86*242f1d38SBarry Smith Mat_SeqAIJ *aij = (Mat_SeqAIJ*)A->data; 87*242f1d38SBarry Smith 88*242f1d38SBarry Smith PetscFunctionBegin; 89*242f1d38SBarry Smith ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr); 90*242f1d38SBarry Smith ierr = PetscMemzero(norms,n*sizeof(PetscReal));CHKERRQ(ierr); 91*242f1d38SBarry Smith if (type == NORM_2) { 92*242f1d38SBarry Smith for (i=0; i<aij->i[m]; i++) { 93*242f1d38SBarry Smith norms[aij->j[i]] += PetscAbsScalar(aij->a[i]*aij->a[i]); 94*242f1d38SBarry Smith } 95*242f1d38SBarry Smith } else if (type == NORM_1) { 96*242f1d38SBarry Smith for (i=0; i<aij->i[m]; i++) { 97*242f1d38SBarry Smith norms[aij->j[i]] += PetscAbsScalar(aij->a[i]); 98*242f1d38SBarry Smith } 99*242f1d38SBarry Smith } else if (type == NORM_INFINITY) { 100*242f1d38SBarry Smith for (i=0; i<aij->i[m]; i++) { 101*242f1d38SBarry Smith norms[aij->j[i]] = PetscMax(PetscAbsScalar(aij->a[i]),norms[aij->j[i]]); 102*242f1d38SBarry Smith } 103*242f1d38SBarry Smith } else SETERRQ(PETSC_ERR_ARG_WRONG,"Unknown NormType"); 104*242f1d38SBarry Smith 105*242f1d38SBarry Smith if (type == NORM_2) { 106*242f1d38SBarry Smith for (i=0; i<n; i++) norms[i] = sqrt(norms[i]); 107*242f1d38SBarry Smith } 108*242f1d38SBarry Smith PetscFunctionReturn(0); 109*242f1d38SBarry Smith } 110*242f1d38SBarry Smith 111*242f1d38SBarry Smith #include "../src/mat/impls/aij/mpi/mpiaij.h" 112*242f1d38SBarry Smith 113*242f1d38SBarry Smith #undef __FUNCT__ 114*242f1d38SBarry Smith #define __FUNCT__ "MatGetColumnNorms_MPIAIJ" 115*242f1d38SBarry Smith PetscErrorCode MatGetColumnNorms_MPIAIJ(Mat A,NormType type,PetscReal *norms) 116*242f1d38SBarry Smith { 117*242f1d38SBarry Smith PetscErrorCode ierr; 118*242f1d38SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 119*242f1d38SBarry Smith PetscInt i,n,*garray = aij->garray; 120*242f1d38SBarry Smith Mat_SeqAIJ *a_aij = (Mat_SeqAIJ*) aij->A->data; 121*242f1d38SBarry Smith Mat_SeqAIJ *b_aij = (Mat_SeqAIJ*) aij->B->data; 122*242f1d38SBarry Smith PetscReal *work; 123*242f1d38SBarry Smith 124*242f1d38SBarry Smith PetscFunctionBegin; 125*242f1d38SBarry Smith ierr = MatGetSize(A,PETSC_NULL,&n);CHKERRQ(ierr); 126*242f1d38SBarry Smith ierr = PetscMalloc(n*sizeof(PetscReal),&work);CHKERRQ(ierr); 127*242f1d38SBarry Smith ierr = PetscMemzero(work,n*sizeof(PetscReal));CHKERRQ(ierr); 128*242f1d38SBarry Smith if (type == NORM_2) { 129*242f1d38SBarry Smith for (i=0; i<a_aij->i[aij->A->rmap->n]; i++) { 130*242f1d38SBarry Smith work[A->rmap->rstart + a_aij->j[i]] += PetscAbsScalar(a_aij->a[i]*a_aij->a[i]); 131*242f1d38SBarry Smith } 132*242f1d38SBarry Smith for (i=0; i<b_aij->i[aij->B->rmap->n]; i++) { 133*242f1d38SBarry Smith work[garray[b_aij->j[i]]] += PetscAbsScalar(b_aij->a[i]*b_aij->a[i]); 134*242f1d38SBarry Smith } 135*242f1d38SBarry Smith } else if (type == NORM_1) { 136*242f1d38SBarry Smith for (i=0; i<a_aij->i[aij->A->rmap->n]; i++) { 137*242f1d38SBarry Smith work[A->rmap->rstart + a_aij->j[i]] += PetscAbsScalar(a_aij->a[i]); 138*242f1d38SBarry Smith } 139*242f1d38SBarry Smith for (i=0; i<b_aij->i[aij->B->rmap->n]; i++) { 140*242f1d38SBarry Smith work[garray[b_aij->j[i]]] += PetscAbsScalar(b_aij->a[i]); 141*242f1d38SBarry Smith } 142*242f1d38SBarry Smith } else if (type == NORM_INFINITY) { 143*242f1d38SBarry Smith for (i=0; i<a_aij->i[aij->A->rmap->n]; i++) { 144*242f1d38SBarry Smith work[A->rmap->rstart + a_aij->j[i]] = PetscMax(PetscAbsScalar(a_aij->a[i]), work[A->rmap->rstart + a_aij->j[i]]); 145*242f1d38SBarry Smith } 146*242f1d38SBarry Smith for (i=0; i<b_aij->i[aij->B->rmap->n]; i++) { 147*242f1d38SBarry Smith work[garray[b_aij->j[i]]] = PetscMax(PetscAbsScalar(b_aij->a[i]),work[garray[b_aij->j[i]]]); 148*242f1d38SBarry Smith } 149*242f1d38SBarry Smith 150*242f1d38SBarry Smith } else SETERRQ(PETSC_ERR_ARG_WRONG,"Unknown NormType"); 151*242f1d38SBarry Smith if (type == NORM_INFINITY) { 152*242f1d38SBarry Smith ierr = MPI_Allreduce(work,norms,n,MPIU_REAL,MPI_MAX,A->hdr.comm);CHKERRQ(ierr); 153*242f1d38SBarry Smith } else { 154*242f1d38SBarry Smith ierr = MPI_Allreduce(work,norms,n,MPIU_REAL,MPI_SUM,A->hdr.comm);CHKERRQ(ierr); 155*242f1d38SBarry Smith } 156*242f1d38SBarry Smith ierr = PetscFree(work);CHKERRQ(ierr); 157*242f1d38SBarry Smith if (type == NORM_2) { 158*242f1d38SBarry Smith for (i=0; i<n; i++) norms[i] = sqrt(norms[i]); 159*242f1d38SBarry Smith } 160*242f1d38SBarry Smith PetscFunctionReturn(0); 161*242f1d38SBarry Smith } 162*242f1d38SBarry Smith 163*242f1d38SBarry Smith #undef __FUNCT__ 164*242f1d38SBarry Smith #define __FUNCT__ "MatGetColumnNorms_SeqDense" 165*242f1d38SBarry Smith PetscErrorCode MatGetColumnNorms_SeqDense(Mat A,NormType type,PetscReal *norms) 166*242f1d38SBarry Smith { 167*242f1d38SBarry Smith PetscErrorCode ierr; 168*242f1d38SBarry Smith PetscInt i,j,m,n; 169*242f1d38SBarry Smith PetscScalar *a; 170*242f1d38SBarry Smith 171*242f1d38SBarry Smith PetscFunctionBegin; 172*242f1d38SBarry Smith ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr); 173*242f1d38SBarry Smith ierr = PetscMemzero(norms,n*sizeof(PetscReal));CHKERRQ(ierr); 174*242f1d38SBarry Smith ierr = MatGetArray(A,&a);CHKERRQ(ierr); 175*242f1d38SBarry Smith if (type == NORM_2) { 176*242f1d38SBarry Smith for (i=0; i<n; i++ ){ 177*242f1d38SBarry Smith for (j=0; j<m; j++) { 178*242f1d38SBarry Smith norms[i] += PetscAbsScalar(a[j]*a[j]); 179*242f1d38SBarry Smith } 180*242f1d38SBarry Smith a += m; 181*242f1d38SBarry Smith } 182*242f1d38SBarry Smith } else if (type == NORM_1) { 183*242f1d38SBarry Smith for (i=0; i<n; i++ ){ 184*242f1d38SBarry Smith for (j=0; j<m; j++) { 185*242f1d38SBarry Smith norms[i] += PetscAbsScalar(a[j]); 186*242f1d38SBarry Smith } 187*242f1d38SBarry Smith a += m; 188*242f1d38SBarry Smith } 189*242f1d38SBarry Smith } else if (type == NORM_INFINITY) { 190*242f1d38SBarry Smith for (i=0; i<n; i++ ){ 191*242f1d38SBarry Smith for (j=0; j<m; j++) { 192*242f1d38SBarry Smith norms[i] = PetscMax(PetscAbsScalar(a[j]),norms[i]); 193*242f1d38SBarry Smith } 194*242f1d38SBarry Smith a += m; 195*242f1d38SBarry Smith } 196*242f1d38SBarry Smith } else SETERRQ(PETSC_ERR_ARG_WRONG,"Unknown NormType"); 197*242f1d38SBarry Smith if (type == NORM_2) { 198*242f1d38SBarry Smith for (i=0; i<n; i++) norms[i] = sqrt(norms[i]); 199*242f1d38SBarry Smith } 200*242f1d38SBarry Smith PetscFunctionReturn(0); 201*242f1d38SBarry Smith } 202*242f1d38SBarry Smith 203*242f1d38SBarry Smith #include "../src/mat/impls/dense/mpi/mpidense.h" 204*242f1d38SBarry Smith #undef __FUNCT__ 205*242f1d38SBarry Smith #define __FUNCT__ "MatGetColumnNorms_MPIDense" 206*242f1d38SBarry Smith PetscErrorCode MatGetColumnNorms_MPIDense(Mat A,NormType type,PetscReal *norms) 207*242f1d38SBarry Smith { 208*242f1d38SBarry Smith PetscErrorCode ierr; 209*242f1d38SBarry Smith PetscInt i,n; 210*242f1d38SBarry Smith Mat_MPIDense *a = (Mat_MPIDense*) A->data; 211*242f1d38SBarry Smith PetscReal *work; 212*242f1d38SBarry Smith 213*242f1d38SBarry Smith PetscFunctionBegin; 214*242f1d38SBarry Smith ierr = MatGetSize(A,PETSC_NULL,&n);CHKERRQ(ierr); 215*242f1d38SBarry Smith ierr = PetscMalloc(n*sizeof(PetscReal),&work);CHKERRQ(ierr); 216*242f1d38SBarry Smith ierr = MatGetColumnNorms_SeqDense(a->A,type,work);CHKERRQ(ierr); 217*242f1d38SBarry Smith if (type == NORM_2) { 218*242f1d38SBarry Smith for (i=0; i<n; i++) work[i] *= work[i]; 219*242f1d38SBarry Smith } 220*242f1d38SBarry Smith if (type == NORM_INFINITY) { 221*242f1d38SBarry Smith ierr = MPI_Allreduce(work,norms,n,MPIU_REAL,MPI_MAX,A->hdr.comm);CHKERRQ(ierr); 222*242f1d38SBarry Smith } else { 223*242f1d38SBarry Smith ierr = MPI_Allreduce(work,norms,n,MPIU_REAL,MPI_SUM,A->hdr.comm);CHKERRQ(ierr); 224*242f1d38SBarry Smith } 225*242f1d38SBarry Smith ierr = PetscFree(work);CHKERRQ(ierr); 226*242f1d38SBarry Smith if (type == NORM_2) { 227*242f1d38SBarry Smith for (i=0; i<n; i++) norms[i] = sqrt(norms[i]); 228*242f1d38SBarry Smith } 229*242f1d38SBarry Smith PetscFunctionReturn(0); 230*242f1d38SBarry Smith } 231*242f1d38SBarry Smith 232*242f1d38SBarry Smith #undef __FUNCT__ 233*242f1d38SBarry Smith #define __FUNCT__ "MatGetColumnNorms" 234*242f1d38SBarry Smith /* 235*242f1d38SBarry Smith MatGetColumnNorms - Gets the 2 norms of each column of a sparse or dense matrix. 236*242f1d38SBarry Smith 237*242f1d38SBarry Smith Input Parameter: 238*242f1d38SBarry Smith + A - the matrix 239*242f1d38SBarry Smith - type - NORM_2, NORM_1 or NORM_INFINITY 240*242f1d38SBarry Smith 241*242f1d38SBarry Smith Output Parameter: 242*242f1d38SBarry Smith . norms - an array as large as the TOTAL number of columns in the matrix 243*242f1d38SBarry Smith 244*242f1d38SBarry Smith Notes: Each process has ALL the column norms after the call. 245*242f1d38SBarry Smith */ 246*242f1d38SBarry Smith PetscErrorCode MatGetColumnNorms(Mat A,NormType type,PetscReal *norms) 247*242f1d38SBarry Smith { 248*242f1d38SBarry Smith PetscErrorCode ierr; 249*242f1d38SBarry Smith PetscTruth flg; 250*242f1d38SBarry Smith 251*242f1d38SBarry Smith PetscFunctionBegin; 252*242f1d38SBarry Smith ierr = PetscTypeCompare((PetscObject)A,MATSEQAIJ,&flg);CHKERRQ(ierr); 253*242f1d38SBarry Smith if (flg) { 254*242f1d38SBarry Smith ierr = MatGetColumnNorms_SeqAIJ(A,type,norms);CHKERRQ(ierr); 255*242f1d38SBarry Smith } else { 256*242f1d38SBarry Smith ierr = PetscTypeCompare((PetscObject)A,MATSEQDENSE,&flg);CHKERRQ(ierr); 257*242f1d38SBarry Smith if (flg) { 258*242f1d38SBarry Smith ierr = MatGetColumnNorms_SeqDense(A,type,norms);CHKERRQ(ierr); 259*242f1d38SBarry Smith } else { 260*242f1d38SBarry Smith ierr = PetscTypeCompare((PetscObject)A,MATMPIDENSE,&flg);CHKERRQ(ierr); 261*242f1d38SBarry Smith if (flg) { 262*242f1d38SBarry Smith ierr = MatGetColumnNorms_MPIDense(A,type,norms);CHKERRQ(ierr); 263*242f1d38SBarry Smith } else { 264*242f1d38SBarry Smith ierr = PetscTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr); 265*242f1d38SBarry Smith if (flg) { 266*242f1d38SBarry Smith ierr = MatGetColumnNorms_MPIAIJ(A,type,norms);CHKERRQ(ierr); 267*242f1d38SBarry Smith } else SETERRQ(PETSC_ERR_SUP,"Not coded for this matrix type"); 268*242f1d38SBarry Smith } 269*242f1d38SBarry Smith } 270*242f1d38SBarry Smith } 271*242f1d38SBarry Smith PetscFunctionReturn(0); 272*242f1d38SBarry Smith } 273*242f1d38SBarry Smith 274