10925cdddSBarry Smith 2af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 30925cdddSBarry Smith 40925cdddSBarry Smith /*@ 582bf6240SBarry Smith MatGetColumnVector - Gets the values from a given column of a matrix. 60925cdddSBarry Smith 772257631SBarry Smith Not Collective 8fee21e36SBarry Smith 998a79cdbSBarry Smith Input Parameters: 105ed6d96aSBarry Smith + A - the matrix 115ed6d96aSBarry Smith . yy - the vector 12b7e53d48SBarry Smith - col - the column requested (in global numbering) 1398a79cdbSBarry Smith 1415091d37SBarry Smith Level: advanced 1515091d37SBarry Smith 169d006df2SBarry Smith Notes: 179a971ab9SStefano Zampini If a Mat type does not implement the operation, each processor for which this is called 189a971ab9SStefano Zampini gets the values for its rows using MatGetRow(). 199d006df2SBarry Smith 203d81755aSBarry Smith The vector must have the same parallel row layout as the matrix. 213d81755aSBarry Smith 2282bf6240SBarry Smith Contributed by: Denis Vanderstraeten 230925cdddSBarry Smith 24db781477SPatrick Sanan .seealso: `MatGetRow()`, `MatGetDiagonal()`, `MatMult()` 2515091d37SBarry Smith 260925cdddSBarry Smith @*/ 277087cfbeSBarry Smith PetscErrorCode MatGetColumnVector(Mat A,Vec yy,PetscInt col) 280925cdddSBarry Smith { 298aa348c1SBarry Smith PetscScalar *y; 30b3cc6726SBarry Smith const PetscScalar *v; 3138baddfdSBarry Smith PetscInt i,j,nz,N,Rs,Re,rs,re; 3238baddfdSBarry Smith const PetscInt *idx; 330925cdddSBarry Smith 340925cdddSBarry Smith PetscFunctionBegin; 350700a824SBarry Smith PetscValidHeaderSpecific(A,MAT_CLASSID,1); 360700a824SBarry Smith PetscValidHeaderSpecific(yy,VEC_CLASSID,2); 379a971ab9SStefano Zampini PetscValidLogicalCollectiveInt(A,col,3); 3808401ef6SPierre Jolivet PetscCheck(col >= 0,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Requested negative column: %" PetscInt_FMT,col); 399566063dSJacob Faibussowitsch PetscCall(MatGetSize(A,NULL,&N)); 4008401ef6SPierre Jolivet PetscCheck(col < N,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Requested column %" PetscInt_FMT " larger than number columns in matrix %" PetscInt_FMT,col,N); 419566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&Rs,&Re)); 429566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(yy,&rs,&re)); 43aed4548fSBarry Smith PetscCheck(Rs == rs && Re == re,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Matrix %" PetscInt_FMT " %" PetscInt_FMT " does not have same ownership range (size) as vector %" PetscInt_FMT " %" PetscInt_FMT,Rs,Re,rs,re); 4482bf6240SBarry Smith 45*dbbe0bcdSBarry Smith if (A->ops->getcolumnvector) PetscUseTypeMethod(A,getcolumnvector ,yy,col); 46*dbbe0bcdSBarry Smith else { 479566063dSJacob Faibussowitsch PetscCall(VecSet(yy,0.0)); 489566063dSJacob Faibussowitsch PetscCall(VecGetArray(yy,&y)); 499a971ab9SStefano Zampini /* TODO for general matrices */ 5082bf6240SBarry Smith for (i=Rs; i<Re; i++) { 519566063dSJacob Faibussowitsch PetscCall(MatGetRow(A,i,&nz,&idx,&v)); 5282bf6240SBarry Smith if (nz && idx[0] <= col) { 5382bf6240SBarry Smith /* 5482bf6240SBarry Smith Should use faster search here 5582bf6240SBarry Smith */ 5682bf6240SBarry Smith for (j=0; j<nz; j++) { 5782bf6240SBarry Smith if (idx[j] >= col) { 5882bf6240SBarry Smith if (idx[j] == col) y[i-rs] = v[j]; 5982bf6240SBarry Smith break; 600925cdddSBarry Smith } 610925cdddSBarry Smith } 620925cdddSBarry Smith } 639566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A,i,&nz,&idx,&v)); 640925cdddSBarry Smith } 659566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(yy,&y)); 668d0534beSBarry Smith } 670925cdddSBarry Smith PetscFunctionReturn(0); 680925cdddSBarry Smith } 69242f1d38SBarry Smith 7086819fdcSBarry Smith /*@ 710716a85fSBarry Smith MatGetColumnNorms - Gets the norms of each column of a sparse or dense matrix. 72242f1d38SBarry Smith 73d8d19677SJose E. Roman Input Parameters: 74242f1d38SBarry Smith + A - the matrix 75242f1d38SBarry Smith - type - NORM_2, NORM_1 or NORM_INFINITY 76242f1d38SBarry Smith 77242f1d38SBarry Smith Output Parameter: 78242f1d38SBarry Smith . norms - an array as large as the TOTAL number of columns in the matrix 79242f1d38SBarry Smith 80f6680f47SSatish Balay Level: intermediate 81f6680f47SSatish Balay 8295452b02SPatrick Sanan Notes: 8395452b02SPatrick Sanan Each process has ALL the column norms after the call. Because of the way this is computed each process gets all the values, 8486819fdcSBarry Smith if each process wants only some of the values it should extract the ones it wants from the array. 8586819fdcSBarry Smith 86db781477SPatrick Sanan .seealso: `NormType`, `MatNorm()` 8786819fdcSBarry Smith 8886819fdcSBarry Smith @*/ 896b9dab40SJed Brown PetscErrorCode MatGetColumnNorms(Mat A,NormType type,PetscReal norms[]) 90242f1d38SBarry Smith { 91857cbf51SRichard Tran Mills /* NOTE: MatGetColumnNorms() could simply be a macro that calls MatGetColumnReductions(). 92857cbf51SRichard Tran Mills * I've kept this as a function because it allows slightly more in the way of error checking, 93857cbf51SRichard Tran Mills * erroring out if MatGetColumnNorms() is not called with a valid NormType. */ 94242f1d38SBarry Smith 95242f1d38SBarry Smith PetscFunctionBegin; 96857cbf51SRichard Tran Mills if (type == NORM_2 || type == NORM_1 || type == NORM_FROBENIUS || type == NORM_INFINITY || type == NORM_1_AND_2) { 979566063dSJacob Faibussowitsch PetscCall(MatGetColumnReductions(A,type,norms)); 98857cbf51SRichard Tran Mills } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Unknown NormType"); 99857cbf51SRichard Tran Mills PetscFunctionReturn(0); 100a873a8cdSSam Reynolds } 101857cbf51SRichard Tran Mills 102857cbf51SRichard Tran Mills /*@ 103857cbf51SRichard Tran Mills MatGetColumnSumsRealPart - Gets the sums of the real part of each column of a sparse or dense matrix. 104857cbf51SRichard Tran Mills 105857cbf51SRichard Tran Mills Input Parameter: 106857cbf51SRichard Tran Mills . A - the matrix 107857cbf51SRichard Tran Mills 108857cbf51SRichard Tran Mills Output Parameter: 109857cbf51SRichard Tran Mills . sums - an array as large as the TOTAL number of columns in the matrix 110857cbf51SRichard Tran Mills 111857cbf51SRichard Tran Mills Level: intermediate 112857cbf51SRichard Tran Mills 113857cbf51SRichard Tran Mills Notes: 114857cbf51SRichard Tran Mills Each process has ALL the column sums after the call. Because of the way this is computed each process gets all the values, 115857cbf51SRichard Tran Mills if each process wants only some of the values it should extract the ones it wants from the array. 116857cbf51SRichard Tran Mills 117db781477SPatrick Sanan .seealso: `MatGetColumnSumsImaginaryPart()`, `VecSum()`, `MatGetColumnMeans()`, `MatGetColumnNorms()`, `MatGetColumnReductions()` 118857cbf51SRichard Tran Mills 119857cbf51SRichard Tran Mills @*/ 120857cbf51SRichard Tran Mills PetscErrorCode MatGetColumnSumsRealPart(Mat A,PetscReal sums[]) 121857cbf51SRichard Tran Mills { 122857cbf51SRichard Tran Mills PetscFunctionBegin; 1239566063dSJacob Faibussowitsch PetscCall(MatGetColumnReductions(A,REDUCTION_SUM_REALPART,sums)); 124857cbf51SRichard Tran Mills PetscFunctionReturn(0); 125857cbf51SRichard Tran Mills } 126857cbf51SRichard Tran Mills 127857cbf51SRichard Tran Mills /*@ 128857cbf51SRichard Tran Mills MatGetColumnSumsImaginaryPart - Gets the sums of the imaginary part of each column of a sparse or dense matrix. 129857cbf51SRichard Tran Mills 130857cbf51SRichard Tran Mills Input Parameter: 131857cbf51SRichard Tran Mills . A - the matrix 132857cbf51SRichard Tran Mills 133857cbf51SRichard Tran Mills Output Parameter: 134857cbf51SRichard Tran Mills . sums - an array as large as the TOTAL number of columns in the matrix 135857cbf51SRichard Tran Mills 136857cbf51SRichard Tran Mills Level: intermediate 137857cbf51SRichard Tran Mills 138857cbf51SRichard Tran Mills Notes: 139857cbf51SRichard Tran Mills Each process has ALL the column sums after the call. Because of the way this is computed each process gets all the values, 140857cbf51SRichard Tran Mills if each process wants only some of the values it should extract the ones it wants from the array. 141857cbf51SRichard Tran Mills 142db781477SPatrick Sanan .seealso: `MatGetColumnSumsRealPart()`, `VecSum()`, `MatGetColumnMeans()`, `MatGetColumnNorms()`, `MatGetColumnReductions()` 143857cbf51SRichard Tran Mills 144857cbf51SRichard Tran Mills @*/ 145857cbf51SRichard Tran Mills PetscErrorCode MatGetColumnSumsImaginaryPart(Mat A,PetscReal sums[]) 146857cbf51SRichard Tran Mills { 147857cbf51SRichard Tran Mills PetscFunctionBegin; 1489566063dSJacob Faibussowitsch PetscCall(MatGetColumnReductions(A,REDUCTION_SUM_IMAGINARYPART,sums)); 149857cbf51SRichard Tran Mills PetscFunctionReturn(0); 150857cbf51SRichard Tran Mills } 151857cbf51SRichard Tran Mills 152857cbf51SRichard Tran Mills /*@ 153857cbf51SRichard Tran Mills MatGetColumnSums - Gets the sums of each column of a sparse or dense matrix. 154857cbf51SRichard Tran Mills 155857cbf51SRichard Tran Mills Input Parameter: 156857cbf51SRichard Tran Mills . A - the matrix 157857cbf51SRichard Tran Mills 158857cbf51SRichard Tran Mills Output Parameter: 159857cbf51SRichard Tran Mills . sums - an array as large as the TOTAL number of columns in the matrix 160857cbf51SRichard Tran Mills 161857cbf51SRichard Tran Mills Level: intermediate 162857cbf51SRichard Tran Mills 163857cbf51SRichard Tran Mills Notes: 164857cbf51SRichard Tran Mills Each process has ALL the column sums after the call. Because of the way this is computed each process gets all the values, 165857cbf51SRichard Tran Mills if each process wants only some of the values it should extract the ones it wants from the array. 166857cbf51SRichard Tran Mills 167db781477SPatrick Sanan .seealso: `VecSum()`, `MatGetColumnMeans()`, `MatGetColumnNorms()`, `MatGetColumnReductions()` 168857cbf51SRichard Tran Mills 169857cbf51SRichard Tran Mills @*/ 170857cbf51SRichard Tran Mills PetscErrorCode MatGetColumnSums(Mat A,PetscScalar sums[]) 171857cbf51SRichard Tran Mills { 172857cbf51SRichard Tran Mills #if defined(PETSC_USE_COMPLEX) 173857cbf51SRichard Tran Mills PetscInt i,n; 174857cbf51SRichard Tran Mills PetscReal *work; 175857cbf51SRichard Tran Mills #endif 176857cbf51SRichard Tran Mills 177857cbf51SRichard Tran Mills PetscFunctionBegin; 178857cbf51SRichard Tran Mills 179857cbf51SRichard Tran Mills #if !defined(PETSC_USE_COMPLEX) 1809566063dSJacob Faibussowitsch PetscCall(MatGetColumnSumsRealPart(A,sums)); 181857cbf51SRichard Tran Mills #else 1829566063dSJacob Faibussowitsch PetscCall(MatGetSize(A,NULL,&n)); 1839566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(sums,n)); 1849566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(n,&work)); 1859566063dSJacob Faibussowitsch PetscCall(MatGetColumnSumsRealPart(A,work)); 186857cbf51SRichard Tran Mills for (i=0; i<n; i++) sums[i] = work[i]; 1879566063dSJacob Faibussowitsch PetscCall(MatGetColumnSumsImaginaryPart(A,work)); 188857cbf51SRichard Tran Mills for (i=0; i<n; i++) sums[i] += work[i]*PETSC_i; 1899566063dSJacob Faibussowitsch PetscCall(PetscFree(work)); 190857cbf51SRichard Tran Mills #endif 191857cbf51SRichard Tran Mills PetscFunctionReturn(0); 192857cbf51SRichard Tran Mills } 193857cbf51SRichard Tran Mills 194857cbf51SRichard Tran Mills /*@ 195857cbf51SRichard Tran Mills MatGetColumnMeansRealPart - Gets the arithmetic means of the real part of each column of a sparse or dense matrix. 196857cbf51SRichard Tran Mills 197857cbf51SRichard Tran Mills Input Parameter: 198857cbf51SRichard Tran Mills . A - the matrix 199857cbf51SRichard Tran Mills 200857cbf51SRichard Tran Mills Output Parameter: 201857cbf51SRichard Tran Mills . sums - an array as large as the TOTAL number of columns in the matrix 202857cbf51SRichard Tran Mills 203857cbf51SRichard Tran Mills Level: intermediate 204857cbf51SRichard Tran Mills 205857cbf51SRichard Tran Mills Notes: 206857cbf51SRichard Tran Mills Each process has ALL the column means after the call. Because of the way this is computed each process gets all the values, 207857cbf51SRichard Tran Mills if each process wants only some of the values it should extract the ones it wants from the array. 208857cbf51SRichard Tran Mills 209db781477SPatrick Sanan .seealso: `MatGetColumnMeansImaginaryPart()`, `VecSum()`, `MatGetColumnSums()`, `MatGetColumnNorms()`, `MatGetColumnReductions()` 210857cbf51SRichard Tran Mills 211857cbf51SRichard Tran Mills @*/ 212857cbf51SRichard Tran Mills PetscErrorCode MatGetColumnMeansRealPart(Mat A,PetscReal means[]) 213857cbf51SRichard Tran Mills { 214857cbf51SRichard Tran Mills PetscFunctionBegin; 2159566063dSJacob Faibussowitsch PetscCall(MatGetColumnReductions(A,REDUCTION_MEAN_REALPART,means)); 216857cbf51SRichard Tran Mills PetscFunctionReturn(0); 217857cbf51SRichard Tran Mills } 218857cbf51SRichard Tran Mills 219857cbf51SRichard Tran Mills /*@ 220857cbf51SRichard Tran Mills MatGetColumnMeansImaginaryPart - Gets the arithmetic means of the imaginary part of each column of a sparse or dense matrix. 221857cbf51SRichard Tran Mills 222857cbf51SRichard Tran Mills Input Parameter: 223857cbf51SRichard Tran Mills . A - the matrix 224857cbf51SRichard Tran Mills 225857cbf51SRichard Tran Mills Output Parameter: 226857cbf51SRichard Tran Mills . sums - an array as large as the TOTAL number of columns in the matrix 227857cbf51SRichard Tran Mills 228857cbf51SRichard Tran Mills Level: intermediate 229857cbf51SRichard Tran Mills 230857cbf51SRichard Tran Mills Notes: 231857cbf51SRichard Tran Mills Each process has ALL the column means after the call. Because of the way this is computed each process gets all the values, 232857cbf51SRichard Tran Mills if each process wants only some of the values it should extract the ones it wants from the array. 233857cbf51SRichard Tran Mills 234db781477SPatrick Sanan .seealso: `MatGetColumnMeansRealPart()`, `VecSum()`, `MatGetColumnSums()`, `MatGetColumnNorms()`, `MatGetColumnReductions()` 235857cbf51SRichard Tran Mills 236857cbf51SRichard Tran Mills @*/ 237857cbf51SRichard Tran Mills PetscErrorCode MatGetColumnMeansImaginaryPart(Mat A,PetscReal means[]) 238857cbf51SRichard Tran Mills { 239857cbf51SRichard Tran Mills PetscFunctionBegin; 2409566063dSJacob Faibussowitsch PetscCall(MatGetColumnReductions(A,REDUCTION_MEAN_IMAGINARYPART,means)); 241857cbf51SRichard Tran Mills PetscFunctionReturn(0); 242857cbf51SRichard Tran Mills } 243857cbf51SRichard Tran Mills 244857cbf51SRichard Tran Mills /*@ 245857cbf51SRichard Tran Mills MatGetColumnMeans - Gets the arithmetic means of each column of a sparse or dense matrix. 246857cbf51SRichard Tran Mills 247857cbf51SRichard Tran Mills Input Parameter: 248857cbf51SRichard Tran Mills . A - the matrix 249857cbf51SRichard Tran Mills 250857cbf51SRichard Tran Mills Output Parameter: 251857cbf51SRichard Tran Mills . means - an array as large as the TOTAL number of columns in the matrix 252857cbf51SRichard Tran Mills 253857cbf51SRichard Tran Mills Level: intermediate 254857cbf51SRichard Tran Mills 255857cbf51SRichard Tran Mills Notes: 256857cbf51SRichard Tran Mills Each process has ALL the column means after the call. Because of the way this is computed each process gets all the values, 257857cbf51SRichard Tran Mills if each process wants only some of the values it should extract the ones it wants from the array. 258857cbf51SRichard Tran Mills 259db781477SPatrick Sanan .seealso: `VecSum()`, `MatGetColumnSums()`, `MatGetColumnNorms()`, `MatGetColumnReductions()` 260857cbf51SRichard Tran Mills 261857cbf51SRichard Tran Mills @*/ 262857cbf51SRichard Tran Mills PetscErrorCode MatGetColumnMeans(Mat A,PetscScalar means[]) 263857cbf51SRichard Tran Mills { 264857cbf51SRichard Tran Mills #if defined(PETSC_USE_COMPLEX) 265857cbf51SRichard Tran Mills PetscInt i,n; 266857cbf51SRichard Tran Mills PetscReal *work; 267857cbf51SRichard Tran Mills #endif 268857cbf51SRichard Tran Mills 269857cbf51SRichard Tran Mills PetscFunctionBegin; 270857cbf51SRichard Tran Mills 271857cbf51SRichard Tran Mills #if !defined(PETSC_USE_COMPLEX) 2729566063dSJacob Faibussowitsch PetscCall(MatGetColumnMeansRealPart(A,means)); 273857cbf51SRichard Tran Mills #else 2749566063dSJacob Faibussowitsch PetscCall(MatGetSize(A,NULL,&n)); 2759566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(means,n)); 2769566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(n,&work)); 2779566063dSJacob Faibussowitsch PetscCall(MatGetColumnMeansRealPart(A,work)); 278857cbf51SRichard Tran Mills for (i=0; i<n; i++) means[i] = work[i]; 2799566063dSJacob Faibussowitsch PetscCall(MatGetColumnMeansImaginaryPart(A,work)); 280857cbf51SRichard Tran Mills for (i=0; i<n; i++) means[i] += work[i]*PETSC_i; 2819566063dSJacob Faibussowitsch PetscCall(PetscFree(work)); 282857cbf51SRichard Tran Mills #endif 283242f1d38SBarry Smith PetscFunctionReturn(0); 284242f1d38SBarry Smith } 285242f1d38SBarry Smith 286a873a8cdSSam Reynolds /*@ 287a873a8cdSSam Reynolds MatGetColumnReductions - Gets the reductions of each column of a sparse or dense matrix. 288a873a8cdSSam Reynolds 28902024617SSatish Balay Input Parameters: 290a873a8cdSSam Reynolds + A - the matrix 291857cbf51SRichard Tran Mills - type - A constant defined in NormType or ReductionType: NORM_2, NORM_1, NORM_INFINITY, REDUCTION_SUM_REALPART, 292857cbf51SRichard Tran Mills REDUCTION_SUM_IMAGINARYPART, REDUCTION_MEAN_REALPART, REDUCTION_MEAN_IMAGINARYPART 293a873a8cdSSam Reynolds 294a873a8cdSSam Reynolds Output Parameter: 295a873a8cdSSam Reynolds . reductions - an array as large as the TOTAL number of columns in the matrix 296a873a8cdSSam Reynolds 297857cbf51SRichard Tran Mills Level: developer 298a873a8cdSSam Reynolds 299a873a8cdSSam Reynolds Notes: 300a873a8cdSSam Reynolds Each process has ALL the column reductions after the call. Because of the way this is computed each process gets all the values, 301a873a8cdSSam Reynolds if each process wants only some of the values it should extract the ones it wants from the array. 302a873a8cdSSam Reynolds 303a873a8cdSSam Reynolds Developer Note: 304857cbf51SRichard Tran Mills This routine is primarily intended as a back-end. 305857cbf51SRichard Tran Mills MatGetColumnNorms(), MatGetColumnSums(), and MatGetColumnMeans() are implemented using this routine. 306a873a8cdSSam Reynolds 307db781477SPatrick Sanan .seealso: `ReductionType`, `NormType`, `MatGetColumnNorms()`, `MatGetColumnSums()`, `MatGetColumnMeans()` 308a873a8cdSSam Reynolds 309a873a8cdSSam Reynolds @*/ 310857cbf51SRichard Tran Mills PetscErrorCode MatGetColumnReductions(Mat A,PetscInt type,PetscReal reductions[]) 311a873a8cdSSam Reynolds { 312a873a8cdSSam Reynolds PetscFunctionBegin; 313a873a8cdSSam Reynolds PetscValidHeaderSpecific(A,MAT_CLASSID,1); 314*dbbe0bcdSBarry Smith PetscUseTypeMethod(A,getcolumnreductions ,type,reductions); 315a873a8cdSSam Reynolds PetscFunctionReturn(0); 316a873a8cdSSam Reynolds } 317