16f79c3a4SBarry Smith 2e090d566SSatish Balay #include "src/mat/matimpl.h" /*I "petscmat.h" I*/ 36f79c3a4SBarry Smith 44a2ae208SSatish Balay #undef __FUNCT__ 54a2ae208SSatish Balay #define __FUNCT__ "MatAXPY" 606be10caSBarry Smith /*@ 721c89e3eSBarry Smith MatAXPY - Computes Y = a*X + Y. 86f79c3a4SBarry Smith 9fee21e36SBarry Smith Collective on Mat 10fee21e36SBarry Smith 1198a79cdbSBarry Smith Input Parameters: 12607cd303SBarry Smith + a - the scalar multiplier 13607cd303SBarry Smith . X - the first matrix 14607cd303SBarry Smith . Y - the second matrix 151a5ee19eSHong Zhang - str - either SAME_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN or SUBSET_NONZERO_PATTERN 1698a79cdbSBarry Smith 17e182c471SBarry Smith Contributed by: Matthew Knepley 18d4bb536fSBarry Smith 192860a424SLois Curfman McInnes Notes: 201a5ee19eSHong Zhang Will only be efficient if one has the SAME_NONZERO_PATTERN or SUBSET_NONZERO_PATTERN 212860a424SLois Curfman McInnes 222860a424SLois Curfman McInnes Level: intermediate 232860a424SLois Curfman McInnes 249cf4f1e8SLois Curfman McInnes .keywords: matrix, add 25d4bb536fSBarry Smith 262860a424SLois Curfman McInnes .seealso: MatAYPX() 2706be10caSBarry Smith @*/ 28268466fbSBarry Smith int MatAXPY(const PetscScalar *a,Mat X,Mat Y,MatStructure str) 296f79c3a4SBarry Smith { 309f12381dSSatish Balay int m1,m2,n1,n2,ierr; 316f79c3a4SBarry Smith 323a40ed3dSBarry Smith PetscFunctionBegin; 334482741eSBarry Smith PetscValidScalarPointer(a,1); 344482741eSBarry Smith PetscValidHeaderSpecific(X,MAT_COOKIE,2); 354482741eSBarry Smith PetscValidHeaderSpecific(Y,MAT_COOKIE,3); 3690f02eecSBarry Smith 37273d9f13SBarry Smith ierr = MatGetSize(X,&m1,&n1);CHKERRQ(ierr); 38273d9f13SBarry Smith ierr = MatGetSize(Y,&m2,&n2);CHKERRQ(ierr); 3929bbc08cSBarry Smith if (m1 != m2 || n1 != n2) SETERRQ4(PETSC_ERR_ARG_SIZ,"Non conforming matrix add: %d %d %d %d",m1,m2,n1,n2); 401987afe7SBarry Smith 41f830108cSBarry Smith if (X->ops->axpy) { 42607cd303SBarry Smith ierr = (*X->ops->axpy)(a,X,Y,str);CHKERRQ(ierr); 43d4bb536fSBarry Smith } else { 44607cd303SBarry Smith ierr = MatAXPY_Basic(a,X,Y,str);CHKERRQ(ierr); 45607cd303SBarry Smith } 46607cd303SBarry Smith PetscFunctionReturn(0); 47607cd303SBarry Smith } 48607cd303SBarry Smith 49607cd303SBarry Smith 50607cd303SBarry Smith #undef __FUNCT__ 51607cd303SBarry Smith #define __FUNCT__ "MatAXPY_Basic" 52268466fbSBarry Smith int MatAXPY_Basic(const PetscScalar *a,Mat X,Mat Y,MatStructure str) 53607cd303SBarry Smith { 54*b3cc6726SBarry Smith int i,start,end,j,ncols,ierr,m,n; 55*b3cc6726SBarry Smith const int *row; 56*b3cc6726SBarry Smith PetscScalar *val; 57*b3cc6726SBarry Smith const PetscScalar *vals; 58607cd303SBarry Smith 59607cd303SBarry Smith PetscFunctionBegin; 608dadbd76SSatish Balay ierr = MatGetSize(X,&m,&n);CHKERRQ(ierr); 6190f02eecSBarry Smith ierr = MatGetOwnershipRange(X,&start,&end);CHKERRQ(ierr); 62d4bb536fSBarry Smith if (*a == 1.0) { 63d4bb536fSBarry Smith for (i = start; i < end; i++) { 64d4bb536fSBarry Smith ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 65d4bb536fSBarry Smith ierr = MatSetValues(Y,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 66d4bb536fSBarry Smith ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 67d4bb536fSBarry Smith } 68d4bb536fSBarry Smith } else { 69*b3cc6726SBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&val);CHKERRQ(ierr); 7006be10caSBarry Smith for (i=start; i<end; i++) { 71*b3cc6726SBarry Smith ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 7206be10caSBarry Smith for (j=0; j<ncols; j++) { 73*b3cc6726SBarry Smith val[j] = (*a)*vals[j]; 746f79c3a4SBarry Smith } 75*b3cc6726SBarry Smith ierr = MatSetValues(Y,1,&i,ncols,row,val,ADD_VALUES);CHKERRQ(ierr); 76*b3cc6726SBarry Smith ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 776f79c3a4SBarry Smith } 78*b3cc6726SBarry Smith ierr = PetscFree(val);CHKERRQ(ierr); 79d4bb536fSBarry Smith } 806d4a8577SBarry Smith ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 816d4a8577SBarry Smith ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 823a40ed3dSBarry Smith PetscFunctionReturn(0); 836f79c3a4SBarry Smith } 84052efed2SBarry Smith 854a2ae208SSatish Balay #undef __FUNCT__ 864a2ae208SSatish Balay #define __FUNCT__ "MatShift" 87052efed2SBarry Smith /*@ 8887828ca2SBarry Smith MatShift - Computes Y = Y + a I, where a is a PetscScalar and I is the identity matrix. 89052efed2SBarry Smith 90fee21e36SBarry Smith Collective on Mat 91fee21e36SBarry Smith 9298a79cdbSBarry Smith Input Parameters: 9398a79cdbSBarry Smith + Y - the matrices 9487828ca2SBarry Smith - a - the PetscScalar 9598a79cdbSBarry Smith 962860a424SLois Curfman McInnes Level: intermediate 972860a424SLois Curfman McInnes 98052efed2SBarry Smith .keywords: matrix, add, shift 996b9ee512SLois Curfman McInnes 100f56f2b3fSBarry Smith .seealso: MatDiagonalSet() 101052efed2SBarry Smith @*/ 102268466fbSBarry Smith int MatShift(const PetscScalar *a,Mat Y) 103052efed2SBarry Smith { 104052efed2SBarry Smith int i,start,end,ierr; 105052efed2SBarry Smith 1063a40ed3dSBarry Smith PetscFunctionBegin; 1074482741eSBarry Smith PetscValidScalarPointer(a,1); 1084482741eSBarry Smith PetscValidHeaderSpecific(Y,MAT_COOKIE,2); 109f830108cSBarry Smith if (Y->ops->shift) { 110f830108cSBarry Smith ierr = (*Y->ops->shift)(a,Y);CHKERRQ(ierr); 11162d58ce1SBarry Smith } else { 112d4bb536fSBarry Smith ierr = MatGetOwnershipRange(Y,&start,&end);CHKERRQ(ierr); 113052efed2SBarry Smith for (i=start; i<end; i++) { 114052efed2SBarry Smith ierr = MatSetValues(Y,1,&i,1,&i,a,ADD_VALUES);CHKERRQ(ierr); 115052efed2SBarry Smith } 1166d4a8577SBarry Smith ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1176d4a8577SBarry Smith ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 118052efed2SBarry Smith } 1193a40ed3dSBarry Smith PetscFunctionReturn(0); 120052efed2SBarry Smith } 1216d84be18SBarry Smith 1224a2ae208SSatish Balay #undef __FUNCT__ 1234a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalSet" 1246d84be18SBarry Smith /*@ 125f56f2b3fSBarry Smith MatDiagonalSet - Computes Y = Y + D, where D is a diagonal matrix 126f56f2b3fSBarry Smith that is represented as a vector. Or Y[i,i] = D[i] if InsertMode is 127f56f2b3fSBarry Smith INSERT_VALUES. 1286d84be18SBarry Smith 1296d84be18SBarry Smith Input Parameters: 13098a79cdbSBarry Smith + Y - the input matrix 131f56f2b3fSBarry Smith . D - the diagonal matrix, represented as a vector 132f56f2b3fSBarry Smith - i - INSERT_VALUES or ADD_VALUES 1336d84be18SBarry Smith 134fee21e36SBarry Smith Collective on Mat and Vec 135fee21e36SBarry Smith 1362860a424SLois Curfman McInnes Level: intermediate 1372860a424SLois Curfman McInnes 1386b9ee512SLois Curfman McInnes .keywords: matrix, add, shift, diagonal 1396b9ee512SLois Curfman McInnes 1406b9ee512SLois Curfman McInnes .seealso: MatShift() 1416d84be18SBarry Smith @*/ 142f56f2b3fSBarry Smith int MatDiagonalSet(Mat Y,Vec D,InsertMode is) 1436d84be18SBarry Smith { 1446d84be18SBarry Smith int i,start,end,ierr; 1456d84be18SBarry Smith 1463a40ed3dSBarry Smith PetscFunctionBegin; 1474482741eSBarry Smith PetscValidHeaderSpecific(Y,MAT_COOKIE,1); 1484482741eSBarry Smith PetscValidHeaderSpecific(D,VEC_COOKIE,2); 149f56f2b3fSBarry Smith if (Y->ops->diagonalset) { 150f56f2b3fSBarry Smith ierr = (*Y->ops->diagonalset)(Y,D,is);CHKERRQ(ierr); 15194d884c6SBarry Smith } else { 1526d84be18SBarry Smith int vstart,vend; 15387828ca2SBarry Smith PetscScalar *v; 1546d84be18SBarry Smith ierr = VecGetOwnershipRange(D,&vstart,&vend);CHKERRQ(ierr); 1556d84be18SBarry Smith ierr = MatGetOwnershipRange(Y,&start,&end);CHKERRQ(ierr); 156d4bb536fSBarry Smith if (vstart != start || vend != end) { 15729bbc08cSBarry Smith SETERRQ4(PETSC_ERR_ARG_SIZ,"Vector ownership range not compatible with matrix: %d %d vec %d %d mat",vstart,vend,start,end); 158d4bb536fSBarry Smith } 1596d84be18SBarry Smith ierr = VecGetArray(D,&v);CHKERRQ(ierr); 1606d84be18SBarry Smith for (i=start; i<end; i++) { 161f56f2b3fSBarry Smith ierr = MatSetValues(Y,1,&i,1,&i,v+i-start,is);CHKERRQ(ierr); 1626d84be18SBarry Smith } 1632e8a6d31SBarry Smith ierr = VecRestoreArray(D,&v);CHKERRQ(ierr); 1646d4a8577SBarry Smith ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1656d4a8577SBarry Smith ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1666d84be18SBarry Smith } 1673a40ed3dSBarry Smith PetscFunctionReturn(0); 1686d84be18SBarry Smith } 169d4bb536fSBarry Smith 1704a2ae208SSatish Balay #undef __FUNCT__ 1714a2ae208SSatish Balay #define __FUNCT__ "MatAYPX" 172d4bb536fSBarry Smith /*@ 173d4bb536fSBarry Smith MatAYPX - Computes Y = X + a*Y. 174d4bb536fSBarry Smith 175fee21e36SBarry Smith Collective on Mat 176fee21e36SBarry Smith 17798a79cdbSBarry Smith Input Parameters: 17898a79cdbSBarry Smith + X,Y - the matrices 17987828ca2SBarry Smith - a - the PetscScalar multiplier 18098a79cdbSBarry Smith 181e182c471SBarry Smith Contributed by: Matthew Knepley 182d4bb536fSBarry Smith 1832860a424SLois Curfman McInnes Notes: 1842860a424SLois Curfman McInnes This routine currently uses the MatAXPY() implementation. 1852860a424SLois Curfman McInnes 186607cd303SBarry Smith This is slow, if you need it fast send email to petsc-maint@mcs.anl.gov 187607cd303SBarry Smith 1882860a424SLois Curfman McInnes Level: intermediate 1892860a424SLois Curfman McInnes 190d4bb536fSBarry Smith .keywords: matrix, add 191d4bb536fSBarry Smith 1922860a424SLois Curfman McInnes .seealso: MatAXPY() 193d4bb536fSBarry Smith @*/ 194268466fbSBarry Smith int MatAYPX(const PetscScalar *a,Mat X,Mat Y) 195d4bb536fSBarry Smith { 19687828ca2SBarry Smith PetscScalar one = 1.0; 197d4bb536fSBarry Smith int mX,mY,nX,nY,ierr; 198d4bb536fSBarry Smith 1993a40ed3dSBarry Smith PetscFunctionBegin; 2004482741eSBarry Smith PetscValidScalarPointer(a,1); 2014482741eSBarry Smith PetscValidHeaderSpecific(X,MAT_COOKIE,2); 2024482741eSBarry Smith PetscValidHeaderSpecific(Y,MAT_COOKIE,3); 203d4bb536fSBarry Smith 204329f5518SBarry Smith ierr = MatGetSize(X,&mX,&nX);CHKERRQ(ierr); 205329f5518SBarry Smith ierr = MatGetSize(X,&mY,&nY);CHKERRQ(ierr); 20629bbc08cSBarry Smith if (mX != mY || nX != nY) SETERRQ4(PETSC_ERR_ARG_SIZ,"Non conforming matrices: %d %d first %d %d second",mX,mY,nX,nY); 207d4bb536fSBarry Smith 208d4bb536fSBarry Smith ierr = MatScale(a,Y);CHKERRQ(ierr); 209607cd303SBarry Smith ierr = MatAXPY(&one,X,Y,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 2103a40ed3dSBarry Smith PetscFunctionReturn(0); 211d4bb536fSBarry Smith } 212b0a32e0cSBarry Smith 2134a2ae208SSatish Balay #undef __FUNCT__ 2144a2ae208SSatish Balay #define __FUNCT__ "MatComputeExplicitOperator" 215b0a32e0cSBarry Smith /*@ 216b0a32e0cSBarry Smith MatComputeExplicitOperator - Computes the explicit matrix 217b0a32e0cSBarry Smith 218b0a32e0cSBarry Smith Collective on Mat 219b0a32e0cSBarry Smith 220b0a32e0cSBarry Smith Input Parameter: 221b0a32e0cSBarry Smith . inmat - the matrix 222b0a32e0cSBarry Smith 223b0a32e0cSBarry Smith Output Parameter: 224b0a32e0cSBarry Smith . mat - the explict preconditioned operator 225b0a32e0cSBarry Smith 226b0a32e0cSBarry Smith Notes: 227b0a32e0cSBarry Smith This computation is done by applying the operators to columns of the 228b0a32e0cSBarry Smith identity matrix. 229b0a32e0cSBarry Smith 230b0a32e0cSBarry Smith Currently, this routine uses a dense matrix format when 1 processor 231b0a32e0cSBarry Smith is used and a sparse format otherwise. This routine is costly in general, 232b0a32e0cSBarry Smith and is recommended for use only with relatively small systems. 233b0a32e0cSBarry Smith 234b0a32e0cSBarry Smith Level: advanced 235b0a32e0cSBarry Smith 236b0a32e0cSBarry Smith .keywords: Mat, compute, explicit, operator 237b0a32e0cSBarry Smith 238b0a32e0cSBarry Smith @*/ 239b0a32e0cSBarry Smith int MatComputeExplicitOperator(Mat inmat,Mat *mat) 240b0a32e0cSBarry Smith { 241b0a32e0cSBarry Smith Vec in,out; 242b0a32e0cSBarry Smith int ierr,i,M,m,size,*rows,start,end; 243b0a32e0cSBarry Smith MPI_Comm comm; 24487828ca2SBarry Smith PetscScalar *array,zero = 0.0,one = 1.0; 245b0a32e0cSBarry Smith 246b0a32e0cSBarry Smith PetscFunctionBegin; 2474482741eSBarry Smith PetscValidHeaderSpecific(inmat,MAT_COOKIE,1); 2484482741eSBarry Smith PetscValidPointer(mat,2); 249b0a32e0cSBarry Smith 250b0a32e0cSBarry Smith comm = inmat->comm; 251b0a32e0cSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 252b0a32e0cSBarry Smith 253b22afee1SSatish Balay ierr = MatGetLocalSize(inmat,&m,0);CHKERRQ(ierr); 254b22afee1SSatish Balay ierr = MatGetSize(inmat,&M,0);CHKERRQ(ierr); 255b0a32e0cSBarry Smith ierr = VecCreateMPI(comm,m,M,&in);CHKERRQ(ierr); 256b0a32e0cSBarry Smith ierr = VecDuplicate(in,&out);CHKERRQ(ierr); 257b0a32e0cSBarry Smith ierr = VecGetOwnershipRange(in,&start,&end);CHKERRQ(ierr); 258b0a32e0cSBarry Smith ierr = PetscMalloc((m+1)*sizeof(int),&rows);CHKERRQ(ierr); 259b0a32e0cSBarry Smith for (i=0; i<m; i++) {rows[i] = start + i;} 260b0a32e0cSBarry Smith 261be5d1d56SKris Buschelman ierr = MatCreate(comm,m,m,M,M,mat);CHKERRQ(ierr); 262b0a32e0cSBarry Smith if (size == 1) { 263be5d1d56SKris Buschelman ierr = MatSetType(*mat,MATSEQDENSE);CHKERRQ(ierr); 264be5d1d56SKris Buschelman ierr = MatSeqDenseSetPreallocation(*mat,PETSC_NULL);CHKERRQ(ierr); 265b0a32e0cSBarry Smith } else { 266be5d1d56SKris Buschelman ierr = MatSetType(*mat,MATMPIAIJ);CHKERRQ(ierr); 267be5d1d56SKris Buschelman ierr = MatMPIAIJSetPreallocation(*mat,0,PETSC_NULL,0,PETSC_NULL);CHKERRQ(ierr); 268b0a32e0cSBarry Smith } 269b0a32e0cSBarry Smith 270b0a32e0cSBarry Smith for (i=0; i<M; i++) { 271b0a32e0cSBarry Smith 272b0a32e0cSBarry Smith ierr = VecSet(&zero,in);CHKERRQ(ierr); 273b0a32e0cSBarry Smith ierr = VecSetValues(in,1,&i,&one,INSERT_VALUES);CHKERRQ(ierr); 274b0a32e0cSBarry Smith ierr = VecAssemblyBegin(in);CHKERRQ(ierr); 275b0a32e0cSBarry Smith ierr = VecAssemblyEnd(in);CHKERRQ(ierr); 276b0a32e0cSBarry Smith 277b0a32e0cSBarry Smith ierr = MatMult(inmat,in,out);CHKERRQ(ierr); 278b0a32e0cSBarry Smith 279b0a32e0cSBarry Smith ierr = VecGetArray(out,&array);CHKERRQ(ierr); 280b0a32e0cSBarry Smith ierr = MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);CHKERRQ(ierr); 281b0a32e0cSBarry Smith ierr = VecRestoreArray(out,&array);CHKERRQ(ierr); 282b0a32e0cSBarry Smith 283b0a32e0cSBarry Smith } 284b0a32e0cSBarry Smith ierr = PetscFree(rows);CHKERRQ(ierr); 285b0a32e0cSBarry Smith ierr = VecDestroy(out);CHKERRQ(ierr); 286b0a32e0cSBarry Smith ierr = VecDestroy(in);CHKERRQ(ierr); 287b0a32e0cSBarry Smith ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 288b0a32e0cSBarry Smith ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 289b0a32e0cSBarry Smith PetscFunctionReturn(0); 290b0a32e0cSBarry Smith } 291b0a32e0cSBarry Smith 29224f910e3SHong Zhang /* Get the map xtoy which is used by MatAXPY() in the case of SUBSET_NONZERO_PATTERN */ 29324f910e3SHong Zhang #undef __FUNCT__ 29424f910e3SHong Zhang #define __FUNCT__ "MatAXPYGetxtoy_Private" 29524f910e3SHong Zhang int MatAXPYGetxtoy_Private(int m,int *xi,int *xj,int *xgarray, int *yi,int *yj,int *ygarray, int **xtoy) 29624f910e3SHong Zhang { 29724f910e3SHong Zhang int ierr,row,i,nz,xcol,ycol,jx,jy,*x2y; 29824f910e3SHong Zhang 29924f910e3SHong Zhang PetscFunctionBegin; 30024f910e3SHong Zhang ierr = PetscMalloc(xi[m]*sizeof(int),&x2y);CHKERRQ(ierr); 30124f910e3SHong Zhang i = 0; 30224f910e3SHong Zhang for (row=0; row<m; row++){ 30324f910e3SHong Zhang nz = xi[1] - xi[0]; 30424f910e3SHong Zhang jy = 0; 30524f910e3SHong Zhang for (jx=0; jx<nz; jx++,jy++){ 30624f910e3SHong Zhang if (xgarray && ygarray){ 30724f910e3SHong Zhang xcol = xgarray[xj[*xi + jx]]; 30824f910e3SHong Zhang ycol = ygarray[yj[*yi + jy]]; 30924f910e3SHong Zhang } else { 31024f910e3SHong Zhang xcol = xj[*xi + jx]; 31124f910e3SHong Zhang ycol = yj[*yi + jy]; /* col index for y */ 31224f910e3SHong Zhang } 31324f910e3SHong Zhang while ( ycol < xcol ) { 31424f910e3SHong Zhang jy++; 31524f910e3SHong Zhang if (ygarray){ 31624f910e3SHong Zhang ycol = ygarray[yj[*yi + jy]]; 31724f910e3SHong Zhang } else { 31824f910e3SHong Zhang ycol = yj[*yi + jy]; 31924f910e3SHong Zhang } 32024f910e3SHong Zhang } 32124f910e3SHong Zhang if (xcol != ycol) SETERRQ2(PETSC_ERR_ARG_WRONG,"X matrix entry (%d,%d) is not in Y matrix",row,ycol); 32224f910e3SHong Zhang x2y[i++] = *yi + jy; 32324f910e3SHong Zhang } 32424f910e3SHong Zhang xi++; yi++; 32524f910e3SHong Zhang } 32624f910e3SHong Zhang *xtoy = x2y; 32724f910e3SHong Zhang PetscFunctionReturn(0); 32824f910e3SHong Zhang } 329