1 /*$Id: axpy.c,v 1.54 2001/08/06 21:16:10 bsmith Exp bsmith $*/ 2 3 #include "src/mat/matimpl.h" /*I "petscmat.h" I*/ 4 5 #undef __FUNCT__ 6 #define __FUNCT__ "MatAXPY" 7 /*@ 8 MatAXPY - Computes Y = a*X + Y. 9 10 Collective on Mat 11 12 Input Parameters: 13 + a - the scalar multiplier 14 . X - the first matrix 15 . Y - the second matrix 16 - str - either SAME_NONZERO_PATTERN or DIFFERENT_NONZERO_PATTERN 17 18 Contributed by: Matthew Knepley 19 20 Notes: 21 Will only be efficient if one has the SAME_NONZERO_PATTERN 22 23 Level: intermediate 24 25 .keywords: matrix, add 26 27 .seealso: MatAYPX() 28 @*/ 29 int MatAXPY(PetscScalar *a,Mat X,Mat Y,MatStructure str) 30 { 31 int m1,m2,n1,n2,i,*row,start,end,j,ncols,ierr; 32 PetscScalar *val,*vals; 33 34 PetscFunctionBegin; 35 PetscValidHeaderSpecific(X,MAT_COOKIE); 36 PetscValidHeaderSpecific(Y,MAT_COOKIE); 37 PetscValidScalarPointer(a); 38 39 ierr = MatGetSize(X,&m1,&n1);CHKERRQ(ierr); 40 ierr = MatGetSize(Y,&m2,&n2);CHKERRQ(ierr); 41 if (m1 != m2 || n1 != n2) SETERRQ4(PETSC_ERR_ARG_SIZ,"Non conforming matrix add: %d %d %d %d",m1,m2,n1,n2); 42 43 if (X->ops->axpy) { 44 ierr = (*X->ops->axpy)(a,X,Y,str);CHKERRQ(ierr); 45 } else { 46 ierr = MatAXPY_Basic(a,X,Y,str);CHKERRQ(ierr); 47 } 48 PetscFunctionReturn(0); 49 } 50 51 52 #undef __FUNCT__ 53 #define __FUNCT__ "MatAXPY_Basic" 54 int MatAXPY_Basic(PetscScalar *a,Mat X,Mat Y,MatStructure str) 55 { 56 int i,*row,start,end,j,ncols,ierr,m,n; 57 PetscScalar *val,*vals; 58 59 PetscFunctionBegin; 60 ierr = MatGetSize(X,&m,&n);CHKERRQ(ierr); 61 ierr = MatGetOwnershipRange(X,&start,&end);CHKERRQ(ierr); 62 if (*a == 1.0) { 63 for (i = start; i < end; i++) { 64 ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 65 ierr = MatSetValues(Y,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 66 ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 67 } 68 } else { 69 ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 70 for (i=start; i<end; i++) { 71 ierr = MatGetRow(X,i,&ncols,&row,&val);CHKERRQ(ierr); 72 for (j=0; j<ncols; j++) { 73 vals[j] = (*a)*val[j]; 74 } 75 ierr = MatSetValues(Y,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 76 ierr = MatRestoreRow(X,i,&ncols,&row,&val);CHKERRQ(ierr); 77 } 78 ierr = PetscFree(vals);CHKERRQ(ierr); 79 } 80 ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 81 ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 82 PetscFunctionReturn(0); 83 } 84 85 #undef __FUNCT__ 86 #define __FUNCT__ "MatShift" 87 /*@ 88 MatShift - Computes Y = Y + a I, where a is a PetscScalar and I is the identity matrix. 89 90 Collective on Mat 91 92 Input Parameters: 93 + Y - the matrices 94 - a - the PetscScalar 95 96 Level: intermediate 97 98 .keywords: matrix, add, shift 99 100 .seealso: MatDiagonalSet() 101 @*/ 102 int MatShift(PetscScalar *a,Mat Y) 103 { 104 int i,start,end,ierr; 105 106 PetscFunctionBegin; 107 PetscValidHeaderSpecific(Y,MAT_COOKIE); 108 PetscValidScalarPointer(a); 109 if (Y->ops->shift) { 110 ierr = (*Y->ops->shift)(a,Y);CHKERRQ(ierr); 111 } else { 112 ierr = MatGetOwnershipRange(Y,&start,&end);CHKERRQ(ierr); 113 for (i=start; i<end; i++) { 114 ierr = MatSetValues(Y,1,&i,1,&i,a,ADD_VALUES);CHKERRQ(ierr); 115 } 116 ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 117 ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 118 } 119 PetscFunctionReturn(0); 120 } 121 122 #undef __FUNCT__ 123 #define __FUNCT__ "MatDiagonalSet" 124 /*@ 125 MatDiagonalSet - Computes Y = Y + D, where D is a diagonal matrix 126 that is represented as a vector. Or Y[i,i] = D[i] if InsertMode is 127 INSERT_VALUES. 128 129 Input Parameters: 130 + Y - the input matrix 131 . D - the diagonal matrix, represented as a vector 132 - i - INSERT_VALUES or ADD_VALUES 133 134 Collective on Mat and Vec 135 136 Level: intermediate 137 138 .keywords: matrix, add, shift, diagonal 139 140 .seealso: MatShift() 141 @*/ 142 int MatDiagonalSet(Mat Y,Vec D,InsertMode is) 143 { 144 int i,start,end,ierr; 145 146 PetscFunctionBegin; 147 PetscValidHeaderSpecific(Y,MAT_COOKIE); 148 PetscValidHeaderSpecific(D,VEC_COOKIE); 149 if (Y->ops->diagonalset) { 150 ierr = (*Y->ops->diagonalset)(Y,D,is);CHKERRQ(ierr); 151 } else { 152 int vstart,vend; 153 PetscScalar *v; 154 ierr = VecGetOwnershipRange(D,&vstart,&vend);CHKERRQ(ierr); 155 ierr = MatGetOwnershipRange(Y,&start,&end);CHKERRQ(ierr); 156 if (vstart != start || vend != end) { 157 SETERRQ4(PETSC_ERR_ARG_SIZ,"Vector ownership range not compatible with matrix: %d %d vec %d %d mat",vstart,vend,start,end); 158 } 159 ierr = VecGetArray(D,&v);CHKERRQ(ierr); 160 for (i=start; i<end; i++) { 161 ierr = MatSetValues(Y,1,&i,1,&i,v+i-start,is);CHKERRQ(ierr); 162 } 163 ierr = VecRestoreArray(D,&v);CHKERRQ(ierr); 164 ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 165 ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 166 } 167 PetscFunctionReturn(0); 168 } 169 170 #undef __FUNCT__ 171 #define __FUNCT__ "MatAYPX" 172 /*@ 173 MatAYPX - Computes Y = X + a*Y. 174 175 Collective on Mat 176 177 Input Parameters: 178 + X,Y - the matrices 179 - a - the PetscScalar multiplier 180 181 Contributed by: Matthew Knepley 182 183 Notes: 184 This routine currently uses the MatAXPY() implementation. 185 186 This is slow, if you need it fast send email to petsc-maint@mcs.anl.gov 187 188 Level: intermediate 189 190 .keywords: matrix, add 191 192 .seealso: MatAXPY() 193 @*/ 194 int MatAYPX(PetscScalar *a,Mat X,Mat Y) 195 { 196 PetscScalar one = 1.0; 197 int mX,mY,nX,nY,ierr; 198 199 PetscFunctionBegin; 200 PetscValidHeaderSpecific(X,MAT_COOKIE); 201 PetscValidHeaderSpecific(Y,MAT_COOKIE); 202 PetscValidScalarPointer(a); 203 204 ierr = MatGetSize(X,&mX,&nX);CHKERRQ(ierr); 205 ierr = MatGetSize(X,&mY,&nY);CHKERRQ(ierr); 206 if (mX != mY || nX != nY) SETERRQ4(PETSC_ERR_ARG_SIZ,"Non conforming matrices: %d %d first %d %d second",mX,mY,nX,nY); 207 208 ierr = MatScale(a,Y);CHKERRQ(ierr); 209 ierr = MatAXPY(&one,X,Y,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 210 PetscFunctionReturn(0); 211 } 212 213 #undef __FUNCT__ 214 #define __FUNCT__ "MatComputeExplicitOperator" 215 /*@ 216 MatComputeExplicitOperator - Computes the explicit matrix 217 218 Collective on Mat 219 220 Input Parameter: 221 . inmat - the matrix 222 223 Output Parameter: 224 . mat - the explict preconditioned operator 225 226 Notes: 227 This computation is done by applying the operators to columns of the 228 identity matrix. 229 230 Currently, this routine uses a dense matrix format when 1 processor 231 is used and a sparse format otherwise. This routine is costly in general, 232 and is recommended for use only with relatively small systems. 233 234 Level: advanced 235 236 .keywords: Mat, compute, explicit, operator 237 238 @*/ 239 int MatComputeExplicitOperator(Mat inmat,Mat *mat) 240 { 241 Vec in,out; 242 int ierr,i,M,m,size,*rows,start,end; 243 MPI_Comm comm; 244 PetscScalar *array,zero = 0.0,one = 1.0; 245 246 PetscFunctionBegin; 247 PetscValidHeaderSpecific(inmat,MAT_COOKIE); 248 PetscValidPointer(mat); 249 250 comm = inmat->comm; 251 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 252 253 ierr = MatGetLocalSize(inmat,&m,0);CHKERRQ(ierr); 254 ierr = MatGetSize(inmat,&M,0);CHKERRQ(ierr); 255 ierr = VecCreateMPI(comm,m,M,&in);CHKERRQ(ierr); 256 ierr = VecDuplicate(in,&out);CHKERRQ(ierr); 257 ierr = VecGetOwnershipRange(in,&start,&end);CHKERRQ(ierr); 258 ierr = PetscMalloc((m+1)*sizeof(int),&rows);CHKERRQ(ierr); 259 for (i=0; i<m; i++) {rows[i] = start + i;} 260 261 if (size == 1) { 262 ierr = MatCreateSeqDense(comm,M,M,PETSC_NULL,mat);CHKERRQ(ierr); 263 } else { 264 ierr = MatCreateMPIAIJ(comm,m,m,M,M,0,0,0,0,mat);CHKERRQ(ierr); 265 } 266 267 for (i=0; i<M; i++) { 268 269 ierr = VecSet(&zero,in);CHKERRQ(ierr); 270 ierr = VecSetValues(in,1,&i,&one,INSERT_VALUES);CHKERRQ(ierr); 271 ierr = VecAssemblyBegin(in);CHKERRQ(ierr); 272 ierr = VecAssemblyEnd(in);CHKERRQ(ierr); 273 274 ierr = MatMult(inmat,in,out);CHKERRQ(ierr); 275 276 ierr = VecGetArray(out,&array);CHKERRQ(ierr); 277 ierr = MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);CHKERRQ(ierr); 278 ierr = VecRestoreArray(out,&array);CHKERRQ(ierr); 279 280 } 281 ierr = PetscFree(rows);CHKERRQ(ierr); 282 ierr = VecDestroy(out);CHKERRQ(ierr); 283 ierr = VecDestroy(in);CHKERRQ(ierr); 284 ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 285 ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 286 PetscFunctionReturn(0); 287 } 288 289