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 } 77242f1d38SBarry Smith 78242f1d38SBarry Smith #include "../src/mat/impls/aij/seq/aij.h" 79242f1d38SBarry Smith 80242f1d38SBarry Smith #undef __FUNCT__ 81242f1d38SBarry Smith #define __FUNCT__ "MatGetColumnNorms_SeqAIJ" 82242f1d38SBarry Smith PetscErrorCode MatGetColumnNorms_SeqAIJ(Mat A,NormType type,PetscReal *norms) 83242f1d38SBarry Smith { 84242f1d38SBarry Smith PetscErrorCode ierr; 85242f1d38SBarry Smith PetscInt i,m,n; 86242f1d38SBarry Smith Mat_SeqAIJ *aij = (Mat_SeqAIJ*)A->data; 87242f1d38SBarry Smith 88242f1d38SBarry Smith PetscFunctionBegin; 89242f1d38SBarry Smith ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr); 90242f1d38SBarry Smith ierr = PetscMemzero(norms,n*sizeof(PetscReal));CHKERRQ(ierr); 91242f1d38SBarry Smith if (type == NORM_2) { 92242f1d38SBarry Smith for (i=0; i<aij->i[m]; i++) { 93242f1d38SBarry Smith norms[aij->j[i]] += PetscAbsScalar(aij->a[i]*aij->a[i]); 94242f1d38SBarry Smith } 95242f1d38SBarry Smith } else if (type == NORM_1) { 96242f1d38SBarry Smith for (i=0; i<aij->i[m]; i++) { 97242f1d38SBarry Smith norms[aij->j[i]] += PetscAbsScalar(aij->a[i]); 98242f1d38SBarry Smith } 99242f1d38SBarry Smith } else if (type == NORM_INFINITY) { 100242f1d38SBarry Smith for (i=0; i<aij->i[m]; i++) { 101242f1d38SBarry Smith norms[aij->j[i]] = PetscMax(PetscAbsScalar(aij->a[i]),norms[aij->j[i]]); 102242f1d38SBarry Smith } 103242f1d38SBarry Smith } else SETERRQ(PETSC_ERR_ARG_WRONG,"Unknown NormType"); 104242f1d38SBarry Smith 105242f1d38SBarry Smith if (type == NORM_2) { 106242f1d38SBarry Smith for (i=0; i<n; i++) norms[i] = sqrt(norms[i]); 107242f1d38SBarry Smith } 108242f1d38SBarry Smith PetscFunctionReturn(0); 109242f1d38SBarry Smith } 110242f1d38SBarry Smith 111242f1d38SBarry Smith #include "../src/mat/impls/aij/mpi/mpiaij.h" 112242f1d38SBarry Smith 113242f1d38SBarry Smith #undef __FUNCT__ 114242f1d38SBarry Smith #define __FUNCT__ "MatGetColumnNorms_MPIAIJ" 115242f1d38SBarry Smith PetscErrorCode MatGetColumnNorms_MPIAIJ(Mat A,NormType type,PetscReal *norms) 116242f1d38SBarry Smith { 117242f1d38SBarry Smith PetscErrorCode ierr; 118242f1d38SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 119242f1d38SBarry Smith PetscInt i,n,*garray = aij->garray; 120242f1d38SBarry Smith Mat_SeqAIJ *a_aij = (Mat_SeqAIJ*) aij->A->data; 121242f1d38SBarry Smith Mat_SeqAIJ *b_aij = (Mat_SeqAIJ*) aij->B->data; 122242f1d38SBarry Smith PetscReal *work; 123242f1d38SBarry Smith 124242f1d38SBarry Smith PetscFunctionBegin; 125242f1d38SBarry Smith ierr = MatGetSize(A,PETSC_NULL,&n);CHKERRQ(ierr); 126242f1d38SBarry Smith ierr = PetscMalloc(n*sizeof(PetscReal),&work);CHKERRQ(ierr); 127242f1d38SBarry Smith ierr = PetscMemzero(work,n*sizeof(PetscReal));CHKERRQ(ierr); 128242f1d38SBarry Smith if (type == NORM_2) { 129242f1d38SBarry Smith for (i=0; i<a_aij->i[aij->A->rmap->n]; i++) { 130242f1d38SBarry Smith work[A->rmap->rstart + a_aij->j[i]] += PetscAbsScalar(a_aij->a[i]*a_aij->a[i]); 131242f1d38SBarry Smith } 132242f1d38SBarry Smith for (i=0; i<b_aij->i[aij->B->rmap->n]; i++) { 133242f1d38SBarry Smith work[garray[b_aij->j[i]]] += PetscAbsScalar(b_aij->a[i]*b_aij->a[i]); 134242f1d38SBarry Smith } 135242f1d38SBarry Smith } else if (type == NORM_1) { 136242f1d38SBarry Smith for (i=0; i<a_aij->i[aij->A->rmap->n]; i++) { 137242f1d38SBarry Smith work[A->rmap->rstart + a_aij->j[i]] += PetscAbsScalar(a_aij->a[i]); 138242f1d38SBarry Smith } 139242f1d38SBarry Smith for (i=0; i<b_aij->i[aij->B->rmap->n]; i++) { 140242f1d38SBarry Smith work[garray[b_aij->j[i]]] += PetscAbsScalar(b_aij->a[i]); 141242f1d38SBarry Smith } 142242f1d38SBarry Smith } else if (type == NORM_INFINITY) { 143242f1d38SBarry Smith for (i=0; i<a_aij->i[aij->A->rmap->n]; i++) { 144242f1d38SBarry Smith work[A->rmap->rstart + a_aij->j[i]] = PetscMax(PetscAbsScalar(a_aij->a[i]), work[A->rmap->rstart + a_aij->j[i]]); 145242f1d38SBarry Smith } 146242f1d38SBarry Smith for (i=0; i<b_aij->i[aij->B->rmap->n]; i++) { 147242f1d38SBarry Smith work[garray[b_aij->j[i]]] = PetscMax(PetscAbsScalar(b_aij->a[i]),work[garray[b_aij->j[i]]]); 148242f1d38SBarry Smith } 149242f1d38SBarry Smith 150242f1d38SBarry Smith } else SETERRQ(PETSC_ERR_ARG_WRONG,"Unknown NormType"); 151242f1d38SBarry Smith if (type == NORM_INFINITY) { 152242f1d38SBarry Smith ierr = MPI_Allreduce(work,norms,n,MPIU_REAL,MPI_MAX,A->hdr.comm);CHKERRQ(ierr); 153242f1d38SBarry Smith } else { 154242f1d38SBarry Smith ierr = MPI_Allreduce(work,norms,n,MPIU_REAL,MPI_SUM,A->hdr.comm);CHKERRQ(ierr); 155242f1d38SBarry Smith } 156242f1d38SBarry Smith ierr = PetscFree(work);CHKERRQ(ierr); 157242f1d38SBarry Smith if (type == NORM_2) { 158242f1d38SBarry Smith for (i=0; i<n; i++) norms[i] = sqrt(norms[i]); 159242f1d38SBarry Smith } 160242f1d38SBarry Smith PetscFunctionReturn(0); 161242f1d38SBarry Smith } 162242f1d38SBarry Smith 163242f1d38SBarry Smith #undef __FUNCT__ 164242f1d38SBarry Smith #define __FUNCT__ "MatGetColumnNorms_SeqDense" 165242f1d38SBarry Smith PetscErrorCode MatGetColumnNorms_SeqDense(Mat A,NormType type,PetscReal *norms) 166242f1d38SBarry Smith { 167242f1d38SBarry Smith PetscErrorCode ierr; 168242f1d38SBarry Smith PetscInt i,j,m,n; 169242f1d38SBarry Smith PetscScalar *a; 170242f1d38SBarry Smith 171242f1d38SBarry Smith PetscFunctionBegin; 172242f1d38SBarry Smith ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr); 173242f1d38SBarry Smith ierr = PetscMemzero(norms,n*sizeof(PetscReal));CHKERRQ(ierr); 174242f1d38SBarry Smith ierr = MatGetArray(A,&a);CHKERRQ(ierr); 175242f1d38SBarry Smith if (type == NORM_2) { 176242f1d38SBarry Smith for (i=0; i<n; i++ ){ 177242f1d38SBarry Smith for (j=0; j<m; j++) { 178242f1d38SBarry Smith norms[i] += PetscAbsScalar(a[j]*a[j]); 179242f1d38SBarry Smith } 180242f1d38SBarry Smith a += m; 181242f1d38SBarry Smith } 182242f1d38SBarry Smith } else if (type == NORM_1) { 183242f1d38SBarry Smith for (i=0; i<n; i++ ){ 184242f1d38SBarry Smith for (j=0; j<m; j++) { 185242f1d38SBarry Smith norms[i] += PetscAbsScalar(a[j]); 186242f1d38SBarry Smith } 187242f1d38SBarry Smith a += m; 188242f1d38SBarry Smith } 189242f1d38SBarry Smith } else if (type == NORM_INFINITY) { 190242f1d38SBarry Smith for (i=0; i<n; i++ ){ 191242f1d38SBarry Smith for (j=0; j<m; j++) { 192242f1d38SBarry Smith norms[i] = PetscMax(PetscAbsScalar(a[j]),norms[i]); 193242f1d38SBarry Smith } 194242f1d38SBarry Smith a += m; 195242f1d38SBarry Smith } 196242f1d38SBarry Smith } else SETERRQ(PETSC_ERR_ARG_WRONG,"Unknown NormType"); 197242f1d38SBarry Smith if (type == NORM_2) { 198242f1d38SBarry Smith for (i=0; i<n; i++) norms[i] = sqrt(norms[i]); 199242f1d38SBarry Smith } 200242f1d38SBarry Smith PetscFunctionReturn(0); 201242f1d38SBarry Smith } 202242f1d38SBarry Smith 203242f1d38SBarry Smith #include "../src/mat/impls/dense/mpi/mpidense.h" 204242f1d38SBarry Smith #undef __FUNCT__ 205242f1d38SBarry Smith #define __FUNCT__ "MatGetColumnNorms_MPIDense" 206242f1d38SBarry Smith PetscErrorCode MatGetColumnNorms_MPIDense(Mat A,NormType type,PetscReal *norms) 207242f1d38SBarry Smith { 208242f1d38SBarry Smith PetscErrorCode ierr; 209242f1d38SBarry Smith PetscInt i,n; 210242f1d38SBarry Smith Mat_MPIDense *a = (Mat_MPIDense*) A->data; 211242f1d38SBarry Smith PetscReal *work; 212242f1d38SBarry Smith 213242f1d38SBarry Smith PetscFunctionBegin; 214242f1d38SBarry Smith ierr = MatGetSize(A,PETSC_NULL,&n);CHKERRQ(ierr); 215242f1d38SBarry Smith ierr = PetscMalloc(n*sizeof(PetscReal),&work);CHKERRQ(ierr); 216242f1d38SBarry Smith ierr = MatGetColumnNorms_SeqDense(a->A,type,work);CHKERRQ(ierr); 217242f1d38SBarry Smith if (type == NORM_2) { 218242f1d38SBarry Smith for (i=0; i<n; i++) work[i] *= work[i]; 219242f1d38SBarry Smith } 220242f1d38SBarry Smith if (type == NORM_INFINITY) { 221242f1d38SBarry Smith ierr = MPI_Allreduce(work,norms,n,MPIU_REAL,MPI_MAX,A->hdr.comm);CHKERRQ(ierr); 222242f1d38SBarry Smith } else { 223242f1d38SBarry Smith ierr = MPI_Allreduce(work,norms,n,MPIU_REAL,MPI_SUM,A->hdr.comm);CHKERRQ(ierr); 224242f1d38SBarry Smith } 225242f1d38SBarry Smith ierr = PetscFree(work);CHKERRQ(ierr); 226242f1d38SBarry Smith if (type == NORM_2) { 227242f1d38SBarry Smith for (i=0; i<n; i++) norms[i] = sqrt(norms[i]); 228242f1d38SBarry Smith } 229242f1d38SBarry Smith PetscFunctionReturn(0); 230242f1d38SBarry Smith } 231242f1d38SBarry Smith 232242f1d38SBarry Smith #undef __FUNCT__ 233242f1d38SBarry Smith #define __FUNCT__ "MatGetColumnNorms" 234*86819fdcSBarry Smith /*@ 235242f1d38SBarry Smith MatGetColumnNorms - Gets the 2 norms of each column of a sparse or dense matrix. 236242f1d38SBarry Smith 237242f1d38SBarry Smith Input Parameter: 238242f1d38SBarry Smith + A - the matrix 239242f1d38SBarry Smith - type - NORM_2, NORM_1 or NORM_INFINITY 240242f1d38SBarry Smith 241242f1d38SBarry Smith Output Parameter: 242242f1d38SBarry Smith . norms - an array as large as the TOTAL number of columns in the matrix 243242f1d38SBarry Smith 244*86819fdcSBarry Smith Notes: Each process has ALL the column norms after the call. Because of the way this is computed each process gets all the values, 245*86819fdcSBarry Smith if each process wants only some of the values it should extract the ones it wants from the array. 246*86819fdcSBarry Smith 247*86819fdcSBarry Smith .seealso: MatGetColumns() 248*86819fdcSBarry Smith 249*86819fdcSBarry Smith @*/ 250242f1d38SBarry Smith PetscErrorCode MatGetColumnNorms(Mat A,NormType type,PetscReal *norms) 251242f1d38SBarry Smith { 252242f1d38SBarry Smith PetscErrorCode ierr; 253242f1d38SBarry Smith PetscTruth flg; 254242f1d38SBarry Smith 255242f1d38SBarry Smith PetscFunctionBegin; 256242f1d38SBarry Smith ierr = PetscTypeCompare((PetscObject)A,MATSEQAIJ,&flg);CHKERRQ(ierr); 257242f1d38SBarry Smith if (flg) { 258242f1d38SBarry Smith ierr = MatGetColumnNorms_SeqAIJ(A,type,norms);CHKERRQ(ierr); 259242f1d38SBarry Smith } else { 260242f1d38SBarry Smith ierr = PetscTypeCompare((PetscObject)A,MATSEQDENSE,&flg);CHKERRQ(ierr); 261242f1d38SBarry Smith if (flg) { 262242f1d38SBarry Smith ierr = MatGetColumnNorms_SeqDense(A,type,norms);CHKERRQ(ierr); 263242f1d38SBarry Smith } else { 264242f1d38SBarry Smith ierr = PetscTypeCompare((PetscObject)A,MATMPIDENSE,&flg);CHKERRQ(ierr); 265242f1d38SBarry Smith if (flg) { 266242f1d38SBarry Smith ierr = MatGetColumnNorms_MPIDense(A,type,norms);CHKERRQ(ierr); 267242f1d38SBarry Smith } else { 268242f1d38SBarry Smith ierr = PetscTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr); 269242f1d38SBarry Smith if (flg) { 270242f1d38SBarry Smith ierr = MatGetColumnNorms_MPIAIJ(A,type,norms);CHKERRQ(ierr); 271242f1d38SBarry Smith } else SETERRQ(PETSC_ERR_SUP,"Not coded for this matrix type"); 272242f1d38SBarry Smith } 273242f1d38SBarry Smith } 274242f1d38SBarry Smith } 275242f1d38SBarry Smith PetscFunctionReturn(0); 276242f1d38SBarry Smith } 277242f1d38SBarry Smith 278