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 249a971ab9SStefano Zampini .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; 31dfbe8321SBarry Smith PetscErrorCode ierr; 3238baddfdSBarry Smith PetscInt i,j,nz,N,Rs,Re,rs,re; 3338baddfdSBarry Smith const PetscInt *idx; 340925cdddSBarry Smith 350925cdddSBarry Smith PetscFunctionBegin; 360700a824SBarry Smith PetscValidHeaderSpecific(A,MAT_CLASSID,1); 370700a824SBarry Smith PetscValidHeaderSpecific(yy,VEC_CLASSID,2); 389a971ab9SStefano Zampini PetscValidLogicalCollectiveInt(A,col,3); 39*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(col < 0,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Requested negative column: %" PetscInt_FMT,col); 400298fd71SBarry Smith ierr = MatGetSize(A,NULL,&N);CHKERRQ(ierr); 41*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(col >= N,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Requested column %" PetscInt_FMT " larger than number columns in matrix %" PetscInt_FMT,col,N); 4282bf6240SBarry Smith ierr = MatGetOwnershipRange(A,&Rs,&Re);CHKERRQ(ierr); 4382bf6240SBarry Smith ierr = VecGetOwnershipRange(yy,&rs,&re);CHKERRQ(ierr); 44*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(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); 4582bf6240SBarry Smith 468d0534beSBarry Smith if (A->ops->getcolumnvector) { 478d0534beSBarry Smith ierr = (*A->ops->getcolumnvector)(A,yy,col);CHKERRQ(ierr); 488d0534beSBarry Smith } else { 498aa348c1SBarry Smith ierr = VecSet(yy,0.0);CHKERRQ(ierr); 5082bf6240SBarry Smith ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 519a971ab9SStefano Zampini /* TODO for general matrices */ 5282bf6240SBarry Smith for (i=Rs; i<Re; i++) { 5382bf6240SBarry Smith ierr = MatGetRow(A,i,&nz,&idx,&v);CHKERRQ(ierr); 5482bf6240SBarry Smith if (nz && idx[0] <= col) { 5582bf6240SBarry Smith /* 5682bf6240SBarry Smith Should use faster search here 5782bf6240SBarry Smith */ 5882bf6240SBarry Smith for (j=0; j<nz; j++) { 5982bf6240SBarry Smith if (idx[j] >= col) { 6082bf6240SBarry Smith if (idx[j] == col) y[i-rs] = v[j]; 6182bf6240SBarry Smith break; 620925cdddSBarry Smith } 630925cdddSBarry Smith } 640925cdddSBarry Smith } 6582bf6240SBarry Smith ierr = MatRestoreRow(A,i,&nz,&idx,&v);CHKERRQ(ierr); 660925cdddSBarry Smith } 6782bf6240SBarry Smith ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 688d0534beSBarry Smith } 690925cdddSBarry Smith PetscFunctionReturn(0); 700925cdddSBarry Smith } 71242f1d38SBarry Smith 7286819fdcSBarry Smith /*@ 730716a85fSBarry Smith MatGetColumnNorms - Gets the norms of each column of a sparse or dense matrix. 74242f1d38SBarry Smith 75d8d19677SJose E. Roman Input Parameters: 76242f1d38SBarry Smith + A - the matrix 77242f1d38SBarry Smith - type - NORM_2, NORM_1 or NORM_INFINITY 78242f1d38SBarry Smith 79242f1d38SBarry Smith Output Parameter: 80242f1d38SBarry Smith . norms - an array as large as the TOTAL number of columns in the matrix 81242f1d38SBarry Smith 82f6680f47SSatish Balay Level: intermediate 83f6680f47SSatish Balay 8495452b02SPatrick Sanan Notes: 8595452b02SPatrick Sanan Each process has ALL the column norms after the call. Because of the way this is computed each process gets all the values, 8686819fdcSBarry Smith if each process wants only some of the values it should extract the ones it wants from the array. 8786819fdcSBarry Smith 8840b1048aSBarry Smith .seealso: NormType, MatNorm() 8986819fdcSBarry Smith 9086819fdcSBarry Smith @*/ 916b9dab40SJed Brown PetscErrorCode MatGetColumnNorms(Mat A,NormType type,PetscReal norms[]) 92242f1d38SBarry Smith { 93857cbf51SRichard Tran Mills /* NOTE: MatGetColumnNorms() could simply be a macro that calls MatGetColumnReductions(). 94857cbf51SRichard Tran Mills * I've kept this as a function because it allows slightly more in the way of error checking, 95857cbf51SRichard Tran Mills * erroring out if MatGetColumnNorms() is not called with a valid NormType. */ 96242f1d38SBarry Smith PetscErrorCode ierr; 97242f1d38SBarry Smith 98242f1d38SBarry Smith PetscFunctionBegin; 99857cbf51SRichard Tran Mills if (type == NORM_2 || type == NORM_1 || type == NORM_FROBENIUS || type == NORM_INFINITY || type == NORM_1_AND_2) { 100857cbf51SRichard Tran Mills ierr = MatGetColumnReductions(A,type,norms);CHKERRQ(ierr); 101857cbf51SRichard Tran Mills } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Unknown NormType"); 102857cbf51SRichard Tran Mills PetscFunctionReturn(0); 103a873a8cdSSam Reynolds } 104857cbf51SRichard Tran Mills 105857cbf51SRichard Tran Mills /*@ 106857cbf51SRichard Tran Mills MatGetColumnSumsRealPart - Gets the sums of the real part of each column of a sparse or dense matrix. 107857cbf51SRichard Tran Mills 108857cbf51SRichard Tran Mills Input Parameter: 109857cbf51SRichard Tran Mills . A - the matrix 110857cbf51SRichard Tran Mills 111857cbf51SRichard Tran Mills Output Parameter: 112857cbf51SRichard Tran Mills . sums - an array as large as the TOTAL number of columns in the matrix 113857cbf51SRichard Tran Mills 114857cbf51SRichard Tran Mills Level: intermediate 115857cbf51SRichard Tran Mills 116857cbf51SRichard Tran Mills Notes: 117857cbf51SRichard 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, 118857cbf51SRichard Tran Mills if each process wants only some of the values it should extract the ones it wants from the array. 119857cbf51SRichard Tran Mills 120857cbf51SRichard Tran Mills .seealso: MatGetColumnSumsImaginaryPart(), VecSum(), MatGetColumnMeans(), MatGetColumnNorms(), MatGetColumnReductions() 121857cbf51SRichard Tran Mills 122857cbf51SRichard Tran Mills @*/ 123857cbf51SRichard Tran Mills PetscErrorCode MatGetColumnSumsRealPart(Mat A,PetscReal sums[]) 124857cbf51SRichard Tran Mills { 125857cbf51SRichard Tran Mills PetscErrorCode ierr; 126857cbf51SRichard Tran Mills 127857cbf51SRichard Tran Mills PetscFunctionBegin; 128857cbf51SRichard Tran Mills ierr = MatGetColumnReductions(A,REDUCTION_SUM_REALPART,sums);CHKERRQ(ierr); 129857cbf51SRichard Tran Mills PetscFunctionReturn(0); 130857cbf51SRichard Tran Mills } 131857cbf51SRichard Tran Mills 132857cbf51SRichard Tran Mills /*@ 133857cbf51SRichard Tran Mills MatGetColumnSumsImaginaryPart - Gets the sums of the imaginary part of each column of a sparse or dense matrix. 134857cbf51SRichard Tran Mills 135857cbf51SRichard Tran Mills Input Parameter: 136857cbf51SRichard Tran Mills . A - the matrix 137857cbf51SRichard Tran Mills 138857cbf51SRichard Tran Mills Output Parameter: 139857cbf51SRichard Tran Mills . sums - an array as large as the TOTAL number of columns in the matrix 140857cbf51SRichard Tran Mills 141857cbf51SRichard Tran Mills Level: intermediate 142857cbf51SRichard Tran Mills 143857cbf51SRichard Tran Mills Notes: 144857cbf51SRichard 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, 145857cbf51SRichard Tran Mills if each process wants only some of the values it should extract the ones it wants from the array. 146857cbf51SRichard Tran Mills 147857cbf51SRichard Tran Mills .seealso: MatGetColumnSumsRealPart(), VecSum(), MatGetColumnMeans(), MatGetColumnNorms(), MatGetColumnReductions() 148857cbf51SRichard Tran Mills 149857cbf51SRichard Tran Mills @*/ 150857cbf51SRichard Tran Mills PetscErrorCode MatGetColumnSumsImaginaryPart(Mat A,PetscReal sums[]) 151857cbf51SRichard Tran Mills { 152857cbf51SRichard Tran Mills PetscErrorCode ierr; 153857cbf51SRichard Tran Mills 154857cbf51SRichard Tran Mills PetscFunctionBegin; 155857cbf51SRichard Tran Mills ierr = MatGetColumnReductions(A,REDUCTION_SUM_IMAGINARYPART,sums);CHKERRQ(ierr); 156857cbf51SRichard Tran Mills PetscFunctionReturn(0); 157857cbf51SRichard Tran Mills } 158857cbf51SRichard Tran Mills 159857cbf51SRichard Tran Mills /*@ 160857cbf51SRichard Tran Mills MatGetColumnSums - Gets the sums of each column of a sparse or dense matrix. 161857cbf51SRichard Tran Mills 162857cbf51SRichard Tran Mills Input Parameter: 163857cbf51SRichard Tran Mills . A - the matrix 164857cbf51SRichard Tran Mills 165857cbf51SRichard Tran Mills Output Parameter: 166857cbf51SRichard Tran Mills . sums - an array as large as the TOTAL number of columns in the matrix 167857cbf51SRichard Tran Mills 168857cbf51SRichard Tran Mills Level: intermediate 169857cbf51SRichard Tran Mills 170857cbf51SRichard Tran Mills Notes: 171857cbf51SRichard 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, 172857cbf51SRichard Tran Mills if each process wants only some of the values it should extract the ones it wants from the array. 173857cbf51SRichard Tran Mills 174857cbf51SRichard Tran Mills .seealso: VecSum(), MatGetColumnMeans(), MatGetColumnNorms(), MatGetColumnReductions() 175857cbf51SRichard Tran Mills 176857cbf51SRichard Tran Mills @*/ 177857cbf51SRichard Tran Mills PetscErrorCode MatGetColumnSums(Mat A,PetscScalar sums[]) 178857cbf51SRichard Tran Mills { 179857cbf51SRichard Tran Mills PetscErrorCode ierr; 180857cbf51SRichard Tran Mills #if defined(PETSC_USE_COMPLEX) 181857cbf51SRichard Tran Mills PetscInt i,n; 182857cbf51SRichard Tran Mills PetscReal *work; 183857cbf51SRichard Tran Mills #endif 184857cbf51SRichard Tran Mills 185857cbf51SRichard Tran Mills PetscFunctionBegin; 186857cbf51SRichard Tran Mills 187857cbf51SRichard Tran Mills #if !defined(PETSC_USE_COMPLEX) 188857cbf51SRichard Tran Mills ierr = MatGetColumnSumsRealPart(A,sums);CHKERRQ(ierr); 189857cbf51SRichard Tran Mills #else 190857cbf51SRichard Tran Mills ierr = MatGetSize(A,NULL,&n);CHKERRQ(ierr); 191857cbf51SRichard Tran Mills ierr = PetscArrayzero(sums,n);CHKERRQ(ierr); 192857cbf51SRichard Tran Mills ierr = PetscCalloc1(n,&work);CHKERRQ(ierr); 193857cbf51SRichard Tran Mills ierr = MatGetColumnSumsRealPart(A,work);CHKERRQ(ierr); 194857cbf51SRichard Tran Mills for (i=0; i<n; i++) sums[i] = work[i]; 195857cbf51SRichard Tran Mills ierr = MatGetColumnSumsImaginaryPart(A,work);CHKERRQ(ierr); 196857cbf51SRichard Tran Mills for (i=0; i<n; i++) sums[i] += work[i]*PETSC_i; 197857cbf51SRichard Tran Mills ierr = PetscFree(work);CHKERRQ(ierr); 198857cbf51SRichard Tran Mills #endif 199857cbf51SRichard Tran Mills PetscFunctionReturn(0); 200857cbf51SRichard Tran Mills } 201857cbf51SRichard Tran Mills 202857cbf51SRichard Tran Mills /*@ 203857cbf51SRichard Tran Mills MatGetColumnMeansRealPart - Gets the arithmetic means of the real part of each column of a sparse or dense matrix. 204857cbf51SRichard Tran Mills 205857cbf51SRichard Tran Mills Input Parameter: 206857cbf51SRichard Tran Mills . A - the matrix 207857cbf51SRichard Tran Mills 208857cbf51SRichard Tran Mills Output Parameter: 209857cbf51SRichard Tran Mills . sums - an array as large as the TOTAL number of columns in the matrix 210857cbf51SRichard Tran Mills 211857cbf51SRichard Tran Mills Level: intermediate 212857cbf51SRichard Tran Mills 213857cbf51SRichard Tran Mills Notes: 214857cbf51SRichard 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, 215857cbf51SRichard Tran Mills if each process wants only some of the values it should extract the ones it wants from the array. 216857cbf51SRichard Tran Mills 217857cbf51SRichard Tran Mills .seealso: MatGetColumnMeansImaginaryPart(), VecSum(), MatGetColumnSums(), MatGetColumnNorms(), MatGetColumnReductions() 218857cbf51SRichard Tran Mills 219857cbf51SRichard Tran Mills @*/ 220857cbf51SRichard Tran Mills PetscErrorCode MatGetColumnMeansRealPart(Mat A,PetscReal means[]) 221857cbf51SRichard Tran Mills { 222857cbf51SRichard Tran Mills PetscErrorCode ierr; 223857cbf51SRichard Tran Mills 224857cbf51SRichard Tran Mills PetscFunctionBegin; 225857cbf51SRichard Tran Mills ierr = MatGetColumnReductions(A,REDUCTION_MEAN_REALPART,means);CHKERRQ(ierr); 226857cbf51SRichard Tran Mills PetscFunctionReturn(0); 227857cbf51SRichard Tran Mills } 228857cbf51SRichard Tran Mills 229857cbf51SRichard Tran Mills /*@ 230857cbf51SRichard Tran Mills MatGetColumnMeansImaginaryPart - Gets the arithmetic means of the imaginary part of each column of a sparse or dense matrix. 231857cbf51SRichard Tran Mills 232857cbf51SRichard Tran Mills Input Parameter: 233857cbf51SRichard Tran Mills . A - the matrix 234857cbf51SRichard Tran Mills 235857cbf51SRichard Tran Mills Output Parameter: 236857cbf51SRichard Tran Mills . sums - an array as large as the TOTAL number of columns in the matrix 237857cbf51SRichard Tran Mills 238857cbf51SRichard Tran Mills Level: intermediate 239857cbf51SRichard Tran Mills 240857cbf51SRichard Tran Mills Notes: 241857cbf51SRichard 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, 242857cbf51SRichard Tran Mills if each process wants only some of the values it should extract the ones it wants from the array. 243857cbf51SRichard Tran Mills 244857cbf51SRichard Tran Mills .seealso: MatGetColumnMeansRealPart(), VecSum(), MatGetColumnSums(), MatGetColumnNorms(), MatGetColumnReductions() 245857cbf51SRichard Tran Mills 246857cbf51SRichard Tran Mills @*/ 247857cbf51SRichard Tran Mills PetscErrorCode MatGetColumnMeansImaginaryPart(Mat A,PetscReal means[]) 248857cbf51SRichard Tran Mills { 249857cbf51SRichard Tran Mills PetscErrorCode ierr; 250857cbf51SRichard Tran Mills 251857cbf51SRichard Tran Mills PetscFunctionBegin; 252857cbf51SRichard Tran Mills ierr = MatGetColumnReductions(A,REDUCTION_MEAN_IMAGINARYPART,means);CHKERRQ(ierr); 253857cbf51SRichard Tran Mills PetscFunctionReturn(0); 254857cbf51SRichard Tran Mills } 255857cbf51SRichard Tran Mills 256857cbf51SRichard Tran Mills /*@ 257857cbf51SRichard Tran Mills MatGetColumnMeans - Gets the arithmetic means of each column of a sparse or dense matrix. 258857cbf51SRichard Tran Mills 259857cbf51SRichard Tran Mills Input Parameter: 260857cbf51SRichard Tran Mills . A - the matrix 261857cbf51SRichard Tran Mills 262857cbf51SRichard Tran Mills Output Parameter: 263857cbf51SRichard Tran Mills . means - an array as large as the TOTAL number of columns in the matrix 264857cbf51SRichard Tran Mills 265857cbf51SRichard Tran Mills Level: intermediate 266857cbf51SRichard Tran Mills 267857cbf51SRichard Tran Mills Notes: 268857cbf51SRichard 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, 269857cbf51SRichard Tran Mills if each process wants only some of the values it should extract the ones it wants from the array. 270857cbf51SRichard Tran Mills 271857cbf51SRichard Tran Mills .seealso: VecSum(), MatGetColumnSums(), MatGetColumnNorms(), MatGetColumnReductions() 272857cbf51SRichard Tran Mills 273857cbf51SRichard Tran Mills @*/ 274857cbf51SRichard Tran Mills PetscErrorCode MatGetColumnMeans(Mat A,PetscScalar means[]) 275857cbf51SRichard Tran Mills { 276857cbf51SRichard Tran Mills PetscErrorCode ierr; 277857cbf51SRichard Tran Mills #if defined(PETSC_USE_COMPLEX) 278857cbf51SRichard Tran Mills PetscInt i,n; 279857cbf51SRichard Tran Mills PetscReal *work; 280857cbf51SRichard Tran Mills #endif 281857cbf51SRichard Tran Mills 282857cbf51SRichard Tran Mills PetscFunctionBegin; 283857cbf51SRichard Tran Mills 284857cbf51SRichard Tran Mills #if !defined(PETSC_USE_COMPLEX) 285857cbf51SRichard Tran Mills ierr = MatGetColumnMeansRealPart(A,means);CHKERRQ(ierr); 286857cbf51SRichard Tran Mills #else 287857cbf51SRichard Tran Mills ierr = MatGetSize(A,NULL,&n);CHKERRQ(ierr); 288857cbf51SRichard Tran Mills ierr = PetscArrayzero(means,n);CHKERRQ(ierr); 289857cbf51SRichard Tran Mills ierr = PetscCalloc1(n,&work);CHKERRQ(ierr); 290857cbf51SRichard Tran Mills ierr = MatGetColumnMeansRealPart(A,work);CHKERRQ(ierr); 291857cbf51SRichard Tran Mills for (i=0; i<n; i++) means[i] = work[i]; 292857cbf51SRichard Tran Mills ierr = MatGetColumnMeansImaginaryPart(A,work);CHKERRQ(ierr); 293857cbf51SRichard Tran Mills for (i=0; i<n; i++) means[i] += work[i]*PETSC_i; 294857cbf51SRichard Tran Mills ierr = PetscFree(work);CHKERRQ(ierr); 295857cbf51SRichard Tran Mills #endif 296242f1d38SBarry Smith PetscFunctionReturn(0); 297242f1d38SBarry Smith } 298242f1d38SBarry Smith 299a873a8cdSSam Reynolds /*@ 300a873a8cdSSam Reynolds MatGetColumnReductions - Gets the reductions of each column of a sparse or dense matrix. 301a873a8cdSSam Reynolds 30202024617SSatish Balay Input Parameters: 303a873a8cdSSam Reynolds + A - the matrix 304857cbf51SRichard Tran Mills - type - A constant defined in NormType or ReductionType: NORM_2, NORM_1, NORM_INFINITY, REDUCTION_SUM_REALPART, 305857cbf51SRichard Tran Mills REDUCTION_SUM_IMAGINARYPART, REDUCTION_MEAN_REALPART, REDUCTION_MEAN_IMAGINARYPART 306a873a8cdSSam Reynolds 307a873a8cdSSam Reynolds Output Parameter: 308a873a8cdSSam Reynolds . reductions - an array as large as the TOTAL number of columns in the matrix 309a873a8cdSSam Reynolds 310857cbf51SRichard Tran Mills Level: developer 311a873a8cdSSam Reynolds 312a873a8cdSSam Reynolds Notes: 313a873a8cdSSam Reynolds Each process has ALL the column reductions after the call. Because of the way this is computed each process gets all the values, 314a873a8cdSSam Reynolds if each process wants only some of the values it should extract the ones it wants from the array. 315a873a8cdSSam Reynolds 316a873a8cdSSam Reynolds Developer Note: 317857cbf51SRichard Tran Mills This routine is primarily intended as a back-end. 318857cbf51SRichard Tran Mills MatGetColumnNorms(), MatGetColumnSums(), and MatGetColumnMeans() are implemented using this routine. 319a873a8cdSSam Reynolds 320857cbf51SRichard Tran Mills .seealso: ReductionType, NormType, MatGetColumnNorms(), MatGetColumnSums(), MatGetColumnMeans() 321a873a8cdSSam Reynolds 322a873a8cdSSam Reynolds @*/ 323857cbf51SRichard Tran Mills PetscErrorCode MatGetColumnReductions(Mat A,PetscInt type,PetscReal reductions[]) 324a873a8cdSSam Reynolds { 325a873a8cdSSam Reynolds PetscErrorCode ierr; 326a873a8cdSSam Reynolds 327a873a8cdSSam Reynolds PetscFunctionBegin; 328a873a8cdSSam Reynolds PetscValidHeaderSpecific(A,MAT_CLASSID,1); 329a873a8cdSSam Reynolds if (A->ops->getcolumnreductions) { 330a873a8cdSSam Reynolds ierr = (*A->ops->getcolumnreductions)(A,type,reductions);CHKERRQ(ierr); 331a873a8cdSSam Reynolds } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Not coded for this matrix type"); 332a873a8cdSSam Reynolds PetscFunctionReturn(0); 333a873a8cdSSam Reynolds } 334