xref: /petsc/src/mat/utils/getcolv.c (revision c6db04a5321582041def2b1e244c75985478b3ef)
10925cdddSBarry Smith 
2*c6db04a5SJed Brown #include <private/matimpl.h>  /*I   "petscmat.h"  I*/
30925cdddSBarry Smith 
44a2ae208SSatish Balay #undef __FUNCT__
54a2ae208SSatish Balay #define __FUNCT__ "MatGetColumnVector"
60925cdddSBarry Smith /*@
782bf6240SBarry Smith    MatGetColumnVector - Gets the values from a given column of a matrix.
80925cdddSBarry Smith 
972257631SBarry Smith    Not Collective
10fee21e36SBarry Smith 
1198a79cdbSBarry Smith    Input Parameters:
125ed6d96aSBarry Smith +  A - the matrix
135ed6d96aSBarry Smith .  yy - the vector
145ed6d96aSBarry Smith -  c - the column requested (in global numbering)
1598a79cdbSBarry Smith 
1615091d37SBarry Smith    Level: advanced
1715091d37SBarry Smith 
189d006df2SBarry Smith    Notes:
199d006df2SBarry Smith    Each processor for which this is called gets the values for its rows.
209d006df2SBarry Smith 
219d006df2SBarry Smith    Since PETSc matrices are usually stored in compressed row format, this routine
229d006df2SBarry Smith    will generally be slow.
239d006df2SBarry Smith 
243d81755aSBarry Smith    The vector must have the same parallel row layout as the matrix.
253d81755aSBarry Smith 
2682bf6240SBarry Smith    Contributed by: Denis Vanderstraeten
270925cdddSBarry Smith 
2882bf6240SBarry Smith .keywords: matrix, column, get
290925cdddSBarry Smith 
3015091d37SBarry Smith .seealso: MatGetRow(), MatGetDiagonal()
3115091d37SBarry Smith 
320925cdddSBarry Smith @*/
337087cfbeSBarry Smith PetscErrorCode  MatGetColumnVector(Mat A,Vec yy,PetscInt col)
340925cdddSBarry Smith {
358aa348c1SBarry Smith   PetscScalar        *y;
36b3cc6726SBarry Smith   const PetscScalar  *v;
37dfbe8321SBarry Smith   PetscErrorCode     ierr;
3838baddfdSBarry Smith   PetscInt           i,j,nz,N,Rs,Re,rs,re;
3938baddfdSBarry Smith   const PetscInt     *idx;
400925cdddSBarry Smith 
410925cdddSBarry Smith   PetscFunctionBegin;
420700a824SBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
430700a824SBarry Smith   PetscValidHeaderSpecific(yy,VEC_CLASSID,2);
44e32f2f54SBarry Smith   if (col < 0)  SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Requested negative column: %D",col);
45011b8408SBarry Smith   ierr = MatGetSize(A,PETSC_NULL,&N);CHKERRQ(ierr);
46e32f2f54SBarry Smith   if (col >= N)  SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Requested column %D larger than number columns in matrix %D",col,N);
4782bf6240SBarry Smith   ierr = MatGetOwnershipRange(A,&Rs,&Re);CHKERRQ(ierr);
4882bf6240SBarry Smith   ierr = VecGetOwnershipRange(yy,&rs,&re);CHKERRQ(ierr);
49e32f2f54SBarry Smith   if (Rs != rs || Re != re) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Matrix %D %D does not have same ownership range (size) as vector %D %D",Rs,Re,rs,re);
5082bf6240SBarry Smith 
518d0534beSBarry Smith   if (A->ops->getcolumnvector) {
528d0534beSBarry Smith     ierr = (*A->ops->getcolumnvector)(A,yy,col);CHKERRQ(ierr);
538d0534beSBarry Smith   } else {
548aa348c1SBarry Smith     ierr = VecSet(yy,0.0);CHKERRQ(ierr);
5582bf6240SBarry Smith     ierr = VecGetArray(yy,&y);CHKERRQ(ierr);
5682bf6240SBarry Smith 
5782bf6240SBarry Smith     for (i=Rs; i<Re; i++) {
5882bf6240SBarry Smith       ierr = MatGetRow(A,i,&nz,&idx,&v);CHKERRQ(ierr);
5982bf6240SBarry Smith       if (nz && idx[0] <= col) {
6082bf6240SBarry Smith 	/*
6182bf6240SBarry Smith           Should use faster search here
6282bf6240SBarry Smith 	*/
6382bf6240SBarry Smith 	for (j=0; j<nz; j++) {
6482bf6240SBarry Smith 	  if (idx[j] >= col) {
6582bf6240SBarry Smith 	    if (idx[j] == col) y[i-rs] = v[j];
6682bf6240SBarry Smith 	    break;
670925cdddSBarry Smith 	  }
680925cdddSBarry Smith 	}
690925cdddSBarry Smith       }
7082bf6240SBarry Smith       ierr = MatRestoreRow(A,i,&nz,&idx,&v);CHKERRQ(ierr);
710925cdddSBarry Smith     }
7282bf6240SBarry Smith     ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr);
738d0534beSBarry Smith   }
740925cdddSBarry Smith   PetscFunctionReturn(0);
750925cdddSBarry Smith }
76242f1d38SBarry Smith 
77*c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h>
78242f1d38SBarry Smith 
79242f1d38SBarry Smith #undef __FUNCT__
80242f1d38SBarry Smith #define __FUNCT__ "MatGetColumnNorms_SeqAIJ"
81242f1d38SBarry Smith PetscErrorCode MatGetColumnNorms_SeqAIJ(Mat A,NormType type,PetscReal *norms)
82242f1d38SBarry Smith {
83242f1d38SBarry Smith   PetscErrorCode ierr;
84242f1d38SBarry Smith   PetscInt       i,m,n;
85242f1d38SBarry Smith   Mat_SeqAIJ     *aij = (Mat_SeqAIJ*)A->data;
86242f1d38SBarry Smith 
87242f1d38SBarry Smith   PetscFunctionBegin;
88242f1d38SBarry Smith   ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr);
89242f1d38SBarry Smith   ierr = PetscMemzero(norms,n*sizeof(PetscReal));CHKERRQ(ierr);
90242f1d38SBarry Smith   if (type == NORM_2) {
91242f1d38SBarry Smith     for (i=0; i<aij->i[m]; i++) {
92242f1d38SBarry Smith       norms[aij->j[i]] += PetscAbsScalar(aij->a[i]*aij->a[i]);
93242f1d38SBarry Smith     }
94242f1d38SBarry Smith   } else if (type == NORM_1) {
95242f1d38SBarry Smith     for (i=0; i<aij->i[m]; i++) {
96242f1d38SBarry Smith       norms[aij->j[i]] += PetscAbsScalar(aij->a[i]);
97242f1d38SBarry Smith     }
98242f1d38SBarry Smith   } else if (type == NORM_INFINITY) {
99242f1d38SBarry Smith     for (i=0; i<aij->i[m]; i++) {
100242f1d38SBarry Smith       norms[aij->j[i]] = PetscMax(PetscAbsScalar(aij->a[i]),norms[aij->j[i]]);
101242f1d38SBarry Smith     }
102e32f2f54SBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Unknown NormType");
103242f1d38SBarry Smith 
104242f1d38SBarry Smith   if (type == NORM_2) {
105242f1d38SBarry Smith     for (i=0; i<n; i++) norms[i] = sqrt(norms[i]);
106242f1d38SBarry Smith   }
107242f1d38SBarry Smith   PetscFunctionReturn(0);
108242f1d38SBarry Smith }
109242f1d38SBarry Smith 
110*c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h>
111242f1d38SBarry Smith 
112242f1d38SBarry Smith #undef __FUNCT__
113242f1d38SBarry Smith #define __FUNCT__ "MatGetColumnNorms_MPIAIJ"
114242f1d38SBarry Smith PetscErrorCode MatGetColumnNorms_MPIAIJ(Mat A,NormType type,PetscReal *norms)
115242f1d38SBarry Smith {
116242f1d38SBarry Smith   PetscErrorCode ierr;
117242f1d38SBarry Smith   Mat_MPIAIJ     *aij = (Mat_MPIAIJ*)A->data;
118242f1d38SBarry Smith   PetscInt       i,n,*garray = aij->garray;
119242f1d38SBarry Smith   Mat_SeqAIJ     *a_aij = (Mat_SeqAIJ*) aij->A->data;
120242f1d38SBarry Smith   Mat_SeqAIJ     *b_aij = (Mat_SeqAIJ*) aij->B->data;
121242f1d38SBarry Smith   PetscReal      *work;
122242f1d38SBarry Smith 
123242f1d38SBarry Smith   PetscFunctionBegin;
124242f1d38SBarry Smith   ierr = MatGetSize(A,PETSC_NULL,&n);CHKERRQ(ierr);
125242f1d38SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscReal),&work);CHKERRQ(ierr);
126242f1d38SBarry Smith   ierr = PetscMemzero(work,n*sizeof(PetscReal));CHKERRQ(ierr);
127242f1d38SBarry Smith   if (type == NORM_2) {
128242f1d38SBarry Smith     for (i=0; i<a_aij->i[aij->A->rmap->n]; i++) {
129c51dce65SBarry Smith       work[A->cmap->rstart + a_aij->j[i]] += PetscAbsScalar(a_aij->a[i]*a_aij->a[i]);
130242f1d38SBarry Smith     }
131242f1d38SBarry Smith     for (i=0; i<b_aij->i[aij->B->rmap->n]; i++) {
132242f1d38SBarry Smith       work[garray[b_aij->j[i]]] += PetscAbsScalar(b_aij->a[i]*b_aij->a[i]);
133242f1d38SBarry Smith     }
134242f1d38SBarry Smith   } else if (type == NORM_1) {
135242f1d38SBarry Smith     for (i=0; i<a_aij->i[aij->A->rmap->n]; i++) {
136c51dce65SBarry Smith       work[A->cmap->rstart + a_aij->j[i]] += PetscAbsScalar(a_aij->a[i]);
137242f1d38SBarry Smith     }
138242f1d38SBarry Smith     for (i=0; i<b_aij->i[aij->B->rmap->n]; i++) {
139242f1d38SBarry Smith       work[garray[b_aij->j[i]]] += PetscAbsScalar(b_aij->a[i]);
140242f1d38SBarry Smith     }
141242f1d38SBarry Smith   } else if (type == NORM_INFINITY) {
142242f1d38SBarry Smith     for (i=0; i<a_aij->i[aij->A->rmap->n]; i++) {
143c51dce65SBarry Smith       work[A->cmap->rstart + a_aij->j[i]] = PetscMax(PetscAbsScalar(a_aij->a[i]), work[A->cmap->rstart + a_aij->j[i]]);
144242f1d38SBarry Smith     }
145242f1d38SBarry Smith     for (i=0; i<b_aij->i[aij->B->rmap->n]; i++) {
146242f1d38SBarry Smith       work[garray[b_aij->j[i]]] = PetscMax(PetscAbsScalar(b_aij->a[i]),work[garray[b_aij->j[i]]]);
147242f1d38SBarry Smith     }
148242f1d38SBarry Smith 
14917186662SBarry Smith   } else SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Unknown NormType");
150242f1d38SBarry Smith   if (type == NORM_INFINITY) {
151242f1d38SBarry Smith     ierr = MPI_Allreduce(work,norms,n,MPIU_REAL,MPI_MAX,A->hdr.comm);CHKERRQ(ierr);
152242f1d38SBarry Smith   } else {
153242f1d38SBarry Smith     ierr = MPI_Allreduce(work,norms,n,MPIU_REAL,MPI_SUM,A->hdr.comm);CHKERRQ(ierr);
154242f1d38SBarry Smith   }
155242f1d38SBarry Smith   ierr = PetscFree(work);CHKERRQ(ierr);
156242f1d38SBarry Smith   if (type == NORM_2) {
157242f1d38SBarry Smith     for (i=0; i<n; i++) norms[i] = sqrt(norms[i]);
158242f1d38SBarry Smith   }
159242f1d38SBarry Smith   PetscFunctionReturn(0);
160242f1d38SBarry Smith }
161242f1d38SBarry Smith 
162242f1d38SBarry Smith #undef __FUNCT__
163242f1d38SBarry Smith #define __FUNCT__ "MatGetColumnNorms_SeqDense"
164242f1d38SBarry Smith PetscErrorCode MatGetColumnNorms_SeqDense(Mat A,NormType type,PetscReal *norms)
165242f1d38SBarry Smith {
166242f1d38SBarry Smith   PetscErrorCode ierr;
167242f1d38SBarry Smith   PetscInt       i,j,m,n;
168242f1d38SBarry Smith   PetscScalar    *a;
169242f1d38SBarry Smith 
170242f1d38SBarry Smith   PetscFunctionBegin;
171242f1d38SBarry Smith   ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr);
172242f1d38SBarry Smith   ierr = PetscMemzero(norms,n*sizeof(PetscReal));CHKERRQ(ierr);
173242f1d38SBarry Smith   ierr = MatGetArray(A,&a);CHKERRQ(ierr);
174242f1d38SBarry Smith   if (type == NORM_2) {
175242f1d38SBarry Smith     for (i=0; i<n; i++ ){
176242f1d38SBarry Smith       for (j=0; j<m; j++) {
177242f1d38SBarry Smith 	norms[i] += PetscAbsScalar(a[j]*a[j]);
178242f1d38SBarry Smith       }
179242f1d38SBarry Smith       a += m;
180242f1d38SBarry Smith     }
181242f1d38SBarry Smith   } else if (type == NORM_1) {
182242f1d38SBarry Smith     for (i=0; i<n; i++ ){
183242f1d38SBarry Smith       for (j=0; j<m; j++) {
184242f1d38SBarry Smith 	norms[i] += PetscAbsScalar(a[j]);
185242f1d38SBarry Smith       }
186242f1d38SBarry Smith       a += m;
187242f1d38SBarry Smith     }
188242f1d38SBarry Smith   } else if (type == NORM_INFINITY) {
189242f1d38SBarry Smith     for (i=0; i<n; i++ ){
190242f1d38SBarry Smith       for (j=0; j<m; j++) {
191242f1d38SBarry Smith 	norms[i] = PetscMax(PetscAbsScalar(a[j]),norms[i]);
192242f1d38SBarry Smith       }
193242f1d38SBarry Smith       a += m;
194242f1d38SBarry Smith     }
19517186662SBarry Smith   } else SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Unknown NormType");
196242f1d38SBarry Smith   if (type == NORM_2) {
197242f1d38SBarry Smith     for (i=0; i<n; i++) norms[i] = sqrt(norms[i]);
198242f1d38SBarry Smith   }
199242f1d38SBarry Smith   PetscFunctionReturn(0);
200242f1d38SBarry Smith }
201242f1d38SBarry Smith 
202*c6db04a5SJed Brown #include <../src/mat/impls/dense/mpi/mpidense.h>
203242f1d38SBarry Smith #undef __FUNCT__
204242f1d38SBarry Smith #define __FUNCT__ "MatGetColumnNorms_MPIDense"
205242f1d38SBarry Smith PetscErrorCode MatGetColumnNorms_MPIDense(Mat A,NormType type,PetscReal *norms)
206242f1d38SBarry Smith {
207242f1d38SBarry Smith   PetscErrorCode ierr;
208242f1d38SBarry Smith   PetscInt       i,n;
209242f1d38SBarry Smith   Mat_MPIDense   *a = (Mat_MPIDense*) A->data;
210242f1d38SBarry Smith   PetscReal      *work;
211242f1d38SBarry Smith 
212242f1d38SBarry Smith   PetscFunctionBegin;
213242f1d38SBarry Smith   ierr = MatGetSize(A,PETSC_NULL,&n);CHKERRQ(ierr);
214242f1d38SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscReal),&work);CHKERRQ(ierr);
215242f1d38SBarry Smith   ierr = MatGetColumnNorms_SeqDense(a->A,type,work);CHKERRQ(ierr);
216242f1d38SBarry Smith   if (type == NORM_2) {
217242f1d38SBarry Smith     for (i=0; i<n; i++) work[i] *= work[i];
218242f1d38SBarry Smith   }
219242f1d38SBarry Smith   if (type == NORM_INFINITY) {
220242f1d38SBarry Smith     ierr = MPI_Allreduce(work,norms,n,MPIU_REAL,MPI_MAX,A->hdr.comm);CHKERRQ(ierr);
221242f1d38SBarry Smith   } else {
222242f1d38SBarry Smith     ierr = MPI_Allreduce(work,norms,n,MPIU_REAL,MPI_SUM,A->hdr.comm);CHKERRQ(ierr);
223242f1d38SBarry Smith   }
224242f1d38SBarry Smith   ierr = PetscFree(work);CHKERRQ(ierr);
225242f1d38SBarry Smith   if (type == NORM_2) {
226242f1d38SBarry Smith     for (i=0; i<n; i++) norms[i] = sqrt(norms[i]);
227242f1d38SBarry Smith   }
228242f1d38SBarry Smith   PetscFunctionReturn(0);
229242f1d38SBarry Smith }
230242f1d38SBarry Smith 
231242f1d38SBarry Smith #undef __FUNCT__
232242f1d38SBarry Smith #define __FUNCT__ "MatGetColumnNorms"
23386819fdcSBarry Smith /*@
234242f1d38SBarry Smith     MatGetColumnNorms - Gets the 2 norms of each column of a sparse or dense matrix.
235242f1d38SBarry Smith 
236242f1d38SBarry Smith   Input Parameter:
237242f1d38SBarry Smith +  A - the matrix
238242f1d38SBarry Smith -  type - NORM_2, NORM_1 or NORM_INFINITY
239242f1d38SBarry Smith 
240242f1d38SBarry Smith   Output Parameter:
241242f1d38SBarry Smith .  norms - an array as large as the TOTAL number of columns in the matrix
242242f1d38SBarry Smith 
243f6680f47SSatish Balay    Level: intermediate
244f6680f47SSatish Balay 
24586819fdcSBarry 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,
24686819fdcSBarry Smith     if each process wants only some of the values it should extract the ones it wants from the array.
24786819fdcSBarry Smith 
24886819fdcSBarry Smith .seealso: MatGetColumns()
24986819fdcSBarry Smith 
25086819fdcSBarry Smith @*/
251242f1d38SBarry Smith PetscErrorCode MatGetColumnNorms(Mat A,NormType type,PetscReal *norms)
252242f1d38SBarry Smith {
253242f1d38SBarry Smith   PetscErrorCode ierr;
254ace3abfcSBarry Smith   PetscBool      flg;
255242f1d38SBarry Smith 
256242f1d38SBarry Smith   PetscFunctionBegin;
257242f1d38SBarry Smith   ierr = PetscTypeCompare((PetscObject)A,MATSEQAIJ,&flg);CHKERRQ(ierr);
258242f1d38SBarry Smith   if (flg) {
259242f1d38SBarry Smith     ierr = MatGetColumnNorms_SeqAIJ(A,type,norms);CHKERRQ(ierr);
260242f1d38SBarry Smith   } else {
261242f1d38SBarry Smith     ierr = PetscTypeCompare((PetscObject)A,MATSEQDENSE,&flg);CHKERRQ(ierr);
262242f1d38SBarry Smith     if (flg) {
263242f1d38SBarry Smith       ierr = MatGetColumnNorms_SeqDense(A,type,norms);CHKERRQ(ierr);
264242f1d38SBarry Smith     } else {
265242f1d38SBarry Smith       ierr = PetscTypeCompare((PetscObject)A,MATMPIDENSE,&flg);CHKERRQ(ierr);
266242f1d38SBarry Smith       if (flg) {
267242f1d38SBarry Smith         ierr = MatGetColumnNorms_MPIDense(A,type,norms);CHKERRQ(ierr);
268242f1d38SBarry Smith       } else {
269242f1d38SBarry Smith         ierr = PetscTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr);
270242f1d38SBarry Smith         if (flg) {
271242f1d38SBarry Smith           ierr = MatGetColumnNorms_MPIAIJ(A,type,norms);CHKERRQ(ierr);
27217186662SBarry Smith         } else SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Not coded for this matrix type");
273242f1d38SBarry Smith       }
274242f1d38SBarry Smith     }
275242f1d38SBarry Smith   }
276242f1d38SBarry Smith   PetscFunctionReturn(0);
277242f1d38SBarry Smith }
278242f1d38SBarry Smith 
279