1 #define PETSCMAT_DLL 2 3 #include "private/matimpl.h" /*I "petscmat.h" I*/ 4 5 #undef __FUNCT__ 6 #define __FUNCT__ "MatAXPY" 7 /*@ 8 MatAXPY - Computes Y = a*X + Y. 9 10 Logically 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, DIFFERENT_NONZERO_PATTERN 17 or SUBSET_NONZERO_PATTERN (nonzeros of X is a subset of Y's) 18 19 Notes: 20 Will only be efficient if one has the SAME_NONZERO_PATTERN or SUBSET_NONZERO_PATTERN 21 22 Level: intermediate 23 24 .keywords: matrix, add 25 26 .seealso: MatAYPX() 27 @*/ 28 PetscErrorCode MatAXPY(Mat Y,PetscScalar a,Mat X,MatStructure str) 29 { 30 PetscErrorCode ierr; 31 PetscInt m1,m2,n1,n2; 32 33 PetscFunctionBegin; 34 PetscValidHeaderSpecific(X,MAT_CLASSID,3); 35 PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 36 PetscValidLogicalCollectiveScalar(Y,a,2); 37 ierr = MatGetSize(X,&m1,&n1);CHKERRQ(ierr); 38 ierr = MatGetSize(Y,&m2,&n2);CHKERRQ(ierr); 39 if (m1 != m2 || n1 != n2) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Non conforming matrix add: %D %D %D %D",m1,m2,n1,n2); 40 41 ierr = PetscLogEventBegin(MAT_AXPY,Y,0,0,0);CHKERRQ(ierr); 42 if (Y->ops->axpy) { 43 ierr = (*Y->ops->axpy)(Y,a,X,str);CHKERRQ(ierr); 44 } else { 45 ierr = MatAXPY_Basic(Y,a,X,str);CHKERRQ(ierr); 46 } 47 ierr = PetscLogEventEnd(MAT_AXPY,Y,0,0,0);CHKERRQ(ierr); 48 PetscFunctionReturn(0); 49 } 50 51 #undef __FUNCT__ 52 #define __FUNCT__ "MatAXPY_Basic" 53 PetscErrorCode MatAXPY_Basic(Mat Y,PetscScalar a,Mat X,MatStructure str) 54 { 55 PetscInt i,start,end,j,ncols,m,n; 56 PetscErrorCode ierr; 57 const PetscInt *row; 58 PetscScalar *val; 59 const PetscScalar *vals; 60 61 PetscFunctionBegin; 62 ierr = MatGetSize(X,&m,&n);CHKERRQ(ierr); 63 ierr = MatGetOwnershipRange(X,&start,&end);CHKERRQ(ierr); 64 if (a == 1.0) { 65 for (i = start; i < end; i++) { 66 ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 67 ierr = MatSetValues(Y,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 68 ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 69 } 70 } else { 71 ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&val);CHKERRQ(ierr); 72 for (i=start; i<end; i++) { 73 ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 74 for (j=0; j<ncols; j++) { 75 val[j] = a*vals[j]; 76 } 77 ierr = MatSetValues(Y,1,&i,ncols,row,val,ADD_VALUES);CHKERRQ(ierr); 78 ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 79 } 80 ierr = PetscFree(val);CHKERRQ(ierr); 81 } 82 ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 83 ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 84 PetscFunctionReturn(0); 85 } 86 87 #undef __FUNCT__ 88 #define __FUNCT__ "MatAXPY_BasicWithPreallocation" 89 PetscErrorCode MatAXPY_BasicWithPreallocation(Mat B,Mat Y,PetscScalar a,Mat X,MatStructure str) 90 { 91 PetscInt i,start,end,j,ncols,m,n; 92 PetscErrorCode ierr; 93 const PetscInt *row; 94 PetscScalar *val; 95 const PetscScalar *vals; 96 97 PetscFunctionBegin; 98 ierr = MatGetSize(X,&m,&n);CHKERRQ(ierr); 99 ierr = MatGetOwnershipRange(X,&start,&end);CHKERRQ(ierr); 100 if (a == 1.0) { 101 for (i = start; i < end; i++) { 102 ierr = MatGetRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr); 103 ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 104 ierr = MatRestoreRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr); 105 106 ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 107 ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 108 ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 109 } 110 } else { 111 ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&val);CHKERRQ(ierr); 112 for (i=start; i<end; i++) { 113 ierr = MatGetRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr); 114 ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 115 ierr = MatRestoreRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr); 116 117 ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 118 for (j=0; j<ncols; j++) { 119 val[j] = a*vals[j]; 120 } 121 ierr = MatSetValues(B,1,&i,ncols,row,val,ADD_VALUES);CHKERRQ(ierr); 122 ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 123 } 124 ierr = PetscFree(val);CHKERRQ(ierr); 125 } 126 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 127 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 128 PetscFunctionReturn(0); 129 } 130 131 #undef __FUNCT__ 132 #define __FUNCT__ "MatShift" 133 /*@ 134 MatShift - Computes Y = Y + a I, where a is a PetscScalar and I is the identity matrix. 135 136 Neighbor-wise Collective on Mat 137 138 Input Parameters: 139 + Y - the matrices 140 - a - the PetscScalar 141 142 Level: intermediate 143 144 .keywords: matrix, add, shift 145 146 .seealso: MatDiagonalSet() 147 @*/ 148 PetscErrorCode MatShift(Mat Y,PetscScalar a) 149 { 150 PetscErrorCode ierr; 151 PetscInt i,start,end; 152 153 PetscFunctionBegin; 154 PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 155 if (!Y->assembled) SETERRQ(((PetscObject)Y)->comm,PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 156 if (Y->factortype) SETERRQ(((PetscObject)Y)->comm,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 157 ierr = MatPreallocated(Y);CHKERRQ(ierr); 158 159 if (Y->ops->shift) { 160 ierr = (*Y->ops->shift)(Y,a);CHKERRQ(ierr); 161 } else { 162 PetscScalar alpha = a; 163 ierr = MatGetOwnershipRange(Y,&start,&end);CHKERRQ(ierr); 164 for (i=start; i<end; i++) { 165 ierr = MatSetValues(Y,1,&i,1,&i,&alpha,ADD_VALUES);CHKERRQ(ierr); 166 } 167 ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 168 ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 169 } 170 PetscFunctionReturn(0); 171 } 172 173 #undef __FUNCT__ 174 #define __FUNCT__ "MatDiagonalSet_Default" 175 PetscErrorCode MatDiagonalSet_Default(Mat Y,Vec D,InsertMode is) 176 { 177 PetscErrorCode ierr; 178 PetscInt i,start,end,vstart,vend; 179 PetscScalar *v; 180 181 PetscFunctionBegin; 182 ierr = VecGetOwnershipRange(D,&vstart,&vend);CHKERRQ(ierr); 183 ierr = MatGetOwnershipRange(Y,&start,&end);CHKERRQ(ierr); 184 if (vstart != start || vend != end) { 185 SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Vector ownership range not compatible with matrix: %D %D vec %D %D mat",vstart,vend,start,end); 186 } 187 ierr = VecGetArray(D,&v);CHKERRQ(ierr); 188 for (i=start; i<end; i++) { 189 ierr = MatSetValues(Y,1,&i,1,&i,v+i-start,is);CHKERRQ(ierr); 190 } 191 ierr = VecRestoreArray(D,&v);CHKERRQ(ierr); 192 ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 193 ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 194 PetscFunctionReturn(0); 195 } 196 197 #undef __FUNCT__ 198 #define __FUNCT__ "MatDiagonalSet" 199 /*@ 200 MatDiagonalSet - Computes Y = Y + D, where D is a diagonal matrix 201 that is represented as a vector. Or Y[i,i] = D[i] if InsertMode is 202 INSERT_VALUES. 203 204 Input Parameters: 205 + Y - the input matrix 206 . D - the diagonal matrix, represented as a vector 207 - i - INSERT_VALUES or ADD_VALUES 208 209 Neighbor-wise Collective on Mat and Vec 210 211 Level: intermediate 212 213 .keywords: matrix, add, shift, diagonal 214 215 .seealso: MatShift() 216 @*/ 217 PetscErrorCode MatDiagonalSet(Mat Y,Vec D,InsertMode is) 218 { 219 PetscErrorCode ierr; 220 221 PetscFunctionBegin; 222 PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 223 PetscValidHeaderSpecific(D,VEC_CLASSID,2); 224 if (Y->ops->diagonalset) { 225 ierr = (*Y->ops->diagonalset)(Y,D,is);CHKERRQ(ierr); 226 } else { 227 ierr = MatDiagonalSet_Default(Y,D,is);CHKERRQ(ierr); 228 } 229 PetscFunctionReturn(0); 230 } 231 232 #undef __FUNCT__ 233 #define __FUNCT__ "MatAYPX" 234 /*@ 235 MatAYPX - Computes Y = a*Y + X. 236 237 Logically on Mat 238 239 Input Parameters: 240 + a - the PetscScalar multiplier 241 . Y - the first matrix 242 . X - the second matrix 243 - str - either SAME_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN or SUBSET_NONZERO_PATTERN 244 245 Notes: 246 Will only be efficient if one has the SAME_NONZERO_PATTERN or SUBSET_NONZERO_PATTERN 247 248 Level: intermediate 249 250 .keywords: matrix, add 251 252 .seealso: MatAXPY() 253 @*/ 254 PetscErrorCode MatAYPX(Mat Y,PetscScalar a,Mat X,MatStructure str) 255 { 256 PetscScalar one = 1.0; 257 PetscErrorCode ierr; 258 PetscInt mX,mY,nX,nY; 259 260 PetscFunctionBegin; 261 PetscValidHeaderSpecific(X,MAT_CLASSID,3); 262 PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 263 PetscValidLogicalCollectiveScalar(Y,a,2); 264 ierr = MatGetSize(X,&mX,&nX);CHKERRQ(ierr); 265 ierr = MatGetSize(X,&mY,&nY);CHKERRQ(ierr); 266 if (mX != mY || nX != nY) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Non conforming matrices: %D %D first %D %D second",mX,mY,nX,nY); 267 268 ierr = MatScale(Y,a);CHKERRQ(ierr); 269 ierr = MatAXPY(Y,one,X,str);CHKERRQ(ierr); 270 PetscFunctionReturn(0); 271 } 272 273 #undef __FUNCT__ 274 #define __FUNCT__ "MatComputeExplicitOperator" 275 /*@ 276 MatComputeExplicitOperator - Computes the explicit matrix 277 278 Collective on Mat 279 280 Input Parameter: 281 . inmat - the matrix 282 283 Output Parameter: 284 . mat - the explict preconditioned operator 285 286 Notes: 287 This computation is done by applying the operators to columns of the 288 identity matrix. 289 290 Currently, this routine uses a dense matrix format when 1 processor 291 is used and a sparse format otherwise. This routine is costly in general, 292 and is recommended for use only with relatively small systems. 293 294 Level: advanced 295 296 .keywords: Mat, compute, explicit, operator 297 298 @*/ 299 PetscErrorCode MatComputeExplicitOperator(Mat inmat,Mat *mat) 300 { 301 Vec in,out; 302 PetscErrorCode ierr; 303 PetscInt i,m,n,M,N,*rows,start,end; 304 MPI_Comm comm; 305 PetscScalar *array,zero = 0.0,one = 1.0; 306 PetscMPIInt size; 307 308 PetscFunctionBegin; 309 PetscValidHeaderSpecific(inmat,MAT_CLASSID,1); 310 PetscValidPointer(mat,2); 311 312 comm = ((PetscObject)inmat)->comm; 313 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 314 315 ierr = MatGetLocalSize(inmat,&m,&n);CHKERRQ(ierr); 316 ierr = MatGetSize(inmat,&M,&N);CHKERRQ(ierr); 317 ierr = MatGetVecs(inmat,&in,&out);CHKERRQ(ierr); 318 ierr = VecSetOption(in,VEC_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 319 ierr = VecGetOwnershipRange(out,&start,&end);CHKERRQ(ierr); 320 ierr = PetscMalloc(m*sizeof(PetscInt),&rows);CHKERRQ(ierr); 321 for (i=0; i<m; i++) {rows[i] = start + i;} 322 323 ierr = MatCreate(comm,mat);CHKERRQ(ierr); 324 ierr = MatSetSizes(*mat,m,n,M,N);CHKERRQ(ierr); 325 if (size == 1) { 326 ierr = MatSetType(*mat,MATSEQDENSE);CHKERRQ(ierr); 327 ierr = MatSeqDenseSetPreallocation(*mat,PETSC_NULL);CHKERRQ(ierr); 328 } else { 329 ierr = MatSetType(*mat,MATMPIAIJ);CHKERRQ(ierr); 330 ierr = MatMPIAIJSetPreallocation(*mat,n,PETSC_NULL,N-n,PETSC_NULL);CHKERRQ(ierr); 331 } 332 333 for (i=0; i<N; i++) { 334 335 ierr = VecSet(in,zero);CHKERRQ(ierr); 336 ierr = VecSetValues(in,1,&i,&one,INSERT_VALUES);CHKERRQ(ierr); 337 ierr = VecAssemblyBegin(in);CHKERRQ(ierr); 338 ierr = VecAssemblyEnd(in);CHKERRQ(ierr); 339 340 ierr = MatMult(inmat,in,out);CHKERRQ(ierr); 341 342 ierr = VecGetArray(out,&array);CHKERRQ(ierr); 343 ierr = MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);CHKERRQ(ierr); 344 ierr = VecRestoreArray(out,&array);CHKERRQ(ierr); 345 346 } 347 ierr = PetscFree(rows);CHKERRQ(ierr); 348 ierr = VecDestroy(out);CHKERRQ(ierr); 349 ierr = VecDestroy(in);CHKERRQ(ierr); 350 ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 351 ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 352 PetscFunctionReturn(0); 353 } 354 355 /* Get the map xtoy which is used by MatAXPY() in the case of SUBSET_NONZERO_PATTERN */ 356 #undef __FUNCT__ 357 #define __FUNCT__ "MatAXPYGetxtoy_Private" 358 PetscErrorCode MatAXPYGetxtoy_Private(PetscInt m,PetscInt *xi,PetscInt *xj,PetscInt *xgarray, PetscInt *yi,PetscInt *yj,PetscInt *ygarray, PetscInt **xtoy) 359 { 360 PetscErrorCode ierr; 361 PetscInt row,i,nz,xcol,ycol,jx,jy,*x2y; 362 363 PetscFunctionBegin; 364 ierr = PetscMalloc(xi[m]*sizeof(PetscInt),&x2y);CHKERRQ(ierr); 365 i = 0; 366 for (row=0; row<m; row++){ 367 nz = xi[1] - xi[0]; 368 jy = 0; 369 for (jx=0; jx<nz; jx++,jy++){ 370 if (xgarray && ygarray){ 371 xcol = xgarray[xj[*xi + jx]]; 372 ycol = ygarray[yj[*yi + jy]]; 373 } else { 374 xcol = xj[*xi + jx]; 375 ycol = yj[*yi + jy]; /* col index for y */ 376 } 377 while ( ycol < xcol ) { 378 jy++; 379 if (ygarray){ 380 ycol = ygarray[yj[*yi + jy]]; 381 } else { 382 ycol = yj[*yi + jy]; 383 } 384 } 385 if (xcol != ycol) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"X matrix entry (%D,%D) is not in Y matrix",row,ycol); 386 x2y[i++] = *yi + jy; 387 } 388 xi++; yi++; 389 } 390 *xtoy = x2y; 391 PetscFunctionReturn(0); 392 } 393