1 2 #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 3 4 #undef __FUNCT__ 5 #define __FUNCT__ "MatAXPY" 6 /*@ 7 MatAXPY - Computes Y = a*X + Y. 8 9 Logically Collective on Mat 10 11 Input Parameters: 12 + a - the scalar multiplier 13 . X - the first matrix 14 . Y - the second matrix 15 - str - either SAME_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN 16 or SUBSET_NONZERO_PATTERN (nonzeros of X is a subset of Y's) 17 18 Level: intermediate 19 20 .keywords: matrix, add 21 22 .seealso: MatAYPX() 23 @*/ 24 PetscErrorCode MatAXPY(Mat Y,PetscScalar a,Mat X,MatStructure str) 25 { 26 PetscErrorCode ierr; 27 PetscInt m1,m2,n1,n2; 28 29 PetscFunctionBegin; 30 PetscValidHeaderSpecific(X,MAT_CLASSID,3); 31 PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 32 PetscValidLogicalCollectiveScalar(Y,a,2); 33 ierr = MatGetSize(X,&m1,&n1);CHKERRQ(ierr); 34 ierr = MatGetSize(Y,&m2,&n2);CHKERRQ(ierr); 35 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); 36 37 ierr = PetscLogEventBegin(MAT_AXPY,Y,0,0,0);CHKERRQ(ierr); 38 if (Y->ops->axpy) { 39 ierr = (*Y->ops->axpy)(Y,a,X,str);CHKERRQ(ierr); 40 } else { 41 ierr = MatAXPY_Basic(Y,a,X,str);CHKERRQ(ierr); 42 } 43 ierr = PetscLogEventEnd(MAT_AXPY,Y,0,0,0);CHKERRQ(ierr); 44 #if defined(PETSC_HAVE_CUSP) 45 if (Y->valid_GPU_matrix != PETSC_CUSP_UNALLOCATED) { 46 Y->valid_GPU_matrix = PETSC_CUSP_CPU; 47 } 48 #endif 49 #if defined(PETSC_HAVE_VIENNACL) 50 if (Y->valid_GPU_matrix != PETSC_VIENNACL_UNALLOCATED) { 51 Y->valid_GPU_matrix = PETSC_VIENNACL_CPU; 52 } 53 #endif 54 PetscFunctionReturn(0); 55 } 56 57 #undef __FUNCT__ 58 #define __FUNCT__ "MatAXPY_Basic" 59 PetscErrorCode MatAXPY_Basic(Mat Y,PetscScalar a,Mat X,MatStructure str) 60 { 61 PetscInt i,start,end,j,ncols,m,n; 62 PetscErrorCode ierr; 63 const PetscInt *row; 64 PetscScalar *val; 65 const PetscScalar *vals; 66 67 PetscFunctionBegin; 68 ierr = MatGetSize(X,&m,&n);CHKERRQ(ierr); 69 ierr = MatGetOwnershipRange(X,&start,&end);CHKERRQ(ierr); 70 if (a == 1.0) { 71 for (i = start; i < end; i++) { 72 ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 73 ierr = MatSetValues(Y,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 74 ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 75 } 76 } else { 77 ierr = PetscMalloc1(n+1,&val);CHKERRQ(ierr); 78 for (i=start; i<end; i++) { 79 ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 80 for (j=0; j<ncols; j++) { 81 val[j] = a*vals[j]; 82 } 83 ierr = MatSetValues(Y,1,&i,ncols,row,val,ADD_VALUES);CHKERRQ(ierr); 84 ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 85 } 86 ierr = PetscFree(val);CHKERRQ(ierr); 87 } 88 ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 89 ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 90 PetscFunctionReturn(0); 91 } 92 93 #undef __FUNCT__ 94 #define __FUNCT__ "MatAXPY_BasicWithPreallocation" 95 PetscErrorCode MatAXPY_BasicWithPreallocation(Mat B,Mat Y,PetscScalar a,Mat X,MatStructure str) 96 { 97 PetscInt i,start,end,j,ncols,m,n; 98 PetscErrorCode ierr; 99 const PetscInt *row; 100 PetscScalar *val; 101 const PetscScalar *vals; 102 103 PetscFunctionBegin; 104 ierr = MatGetSize(X,&m,&n);CHKERRQ(ierr); 105 ierr = MatGetOwnershipRange(X,&start,&end);CHKERRQ(ierr); 106 if (a == 1.0) { 107 for (i = start; i < end; i++) { 108 ierr = MatGetRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr); 109 ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 110 ierr = MatRestoreRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr); 111 112 ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 113 ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 114 ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 115 } 116 } else { 117 ierr = PetscMalloc1(n+1,&val);CHKERRQ(ierr); 118 for (i=start; i<end; i++) { 119 ierr = MatGetRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr); 120 ierr = MatSetValues(B,1,&i,ncols,row,vals,ADD_VALUES);CHKERRQ(ierr); 121 ierr = MatRestoreRow(Y,i,&ncols,&row,&vals);CHKERRQ(ierr); 122 123 ierr = MatGetRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 124 for (j=0; j<ncols; j++) { 125 val[j] = a*vals[j]; 126 } 127 ierr = MatSetValues(B,1,&i,ncols,row,val,ADD_VALUES);CHKERRQ(ierr); 128 ierr = MatRestoreRow(X,i,&ncols,&row,&vals);CHKERRQ(ierr); 129 } 130 ierr = PetscFree(val);CHKERRQ(ierr); 131 } 132 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 133 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 134 PetscFunctionReturn(0); 135 } 136 137 #undef __FUNCT__ 138 #define __FUNCT__ "MatShift" 139 /*@ 140 MatShift - Computes Y = Y + a I, where a is a PetscScalar and I is the identity matrix. 141 142 Neighbor-wise Collective on Mat 143 144 Input Parameters: 145 + Y - the matrices 146 - a - the PetscScalar 147 148 Level: intermediate 149 150 Notes: If the matrix Y is missing some diagonal entries this routine can be very slow. To make it fast one should initially 151 fill the matrix so that all diagonal entries have a value (with a value of zero for those locations that would not have an 152 entry). 153 154 Developers Note: If the local "diagonal part" of the matrix Y has no entries then the local diagonal part is 155 preallocated with 1 nonzero per row for the to be added values. This allows for fast shifting of an empty matrix. 156 157 .keywords: matrix, add, shift 158 159 .seealso: MatDiagonalSet() 160 @*/ 161 PetscErrorCode MatShift(Mat Y,PetscScalar a) 162 { 163 PetscErrorCode ierr; 164 165 PetscFunctionBegin; 166 PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 167 if (!Y->assembled) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 168 if (Y->factortype) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 169 MatCheckPreallocated(Y,1); 170 171 ierr = (*Y->ops->shift)(Y,a);CHKERRQ(ierr); 172 173 #if defined(PETSC_HAVE_CUSP) 174 if (Y->valid_GPU_matrix != PETSC_CUSP_UNALLOCATED) { 175 Y->valid_GPU_matrix = PETSC_CUSP_CPU; 176 } 177 #endif 178 #if defined(PETSC_HAVE_VIENNACL) 179 if (Y->valid_GPU_matrix != PETSC_VIENNACL_UNALLOCATED) { 180 Y->valid_GPU_matrix = PETSC_VIENNACL_CPU; 181 } 182 #endif 183 PetscFunctionReturn(0); 184 } 185 186 #undef __FUNCT__ 187 #define __FUNCT__ "MatDiagonalSet_Default" 188 PetscErrorCode MatDiagonalSet_Default(Mat Y,Vec D,InsertMode is) 189 { 190 PetscErrorCode ierr; 191 PetscInt i,start,end; 192 PetscScalar *v; 193 194 PetscFunctionBegin; 195 ierr = MatGetOwnershipRange(Y,&start,&end);CHKERRQ(ierr); 196 ierr = VecGetArray(D,&v);CHKERRQ(ierr); 197 for (i=start; i<end; i++) { 198 ierr = MatSetValues(Y,1,&i,1,&i,v+i-start,is);CHKERRQ(ierr); 199 } 200 ierr = VecRestoreArray(D,&v);CHKERRQ(ierr); 201 ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 202 ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 203 PetscFunctionReturn(0); 204 } 205 206 #undef __FUNCT__ 207 #define __FUNCT__ "MatDiagonalSet" 208 /*@ 209 MatDiagonalSet - Computes Y = Y + D, where D is a diagonal matrix 210 that is represented as a vector. Or Y[i,i] = D[i] if InsertMode is 211 INSERT_VALUES. 212 213 Input Parameters: 214 + Y - the input matrix 215 . D - the diagonal matrix, represented as a vector 216 - i - INSERT_VALUES or ADD_VALUES 217 218 Neighbor-wise Collective on Mat and Vec 219 220 Notes: If the matrix Y is missing some diagonal entries this routine can be very slow. To make it fast one should initially 221 fill the matrix so that all diagonal entries have a value (with a value of zero for those locations that would not have an 222 entry). 223 224 Level: intermediate 225 226 .keywords: matrix, add, shift, diagonal 227 228 .seealso: MatShift() 229 @*/ 230 PetscErrorCode MatDiagonalSet(Mat Y,Vec D,InsertMode is) 231 { 232 PetscErrorCode ierr; 233 PetscInt matlocal,veclocal; 234 235 PetscFunctionBegin; 236 PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 237 PetscValidHeaderSpecific(D,VEC_CLASSID,2); 238 ierr = MatGetLocalSize(Y,&matlocal,NULL);CHKERRQ(ierr); 239 ierr = VecGetLocalSize(D,&veclocal);CHKERRQ(ierr); 240 if (matlocal != veclocal) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number local rows of matrix %D does not match that of vector for diagonal %D",matlocal,veclocal); 241 if (Y->ops->diagonalset) { 242 ierr = (*Y->ops->diagonalset)(Y,D,is);CHKERRQ(ierr); 243 } else { 244 ierr = MatDiagonalSet_Default(Y,D,is);CHKERRQ(ierr); 245 } 246 PetscFunctionReturn(0); 247 } 248 249 #undef __FUNCT__ 250 #define __FUNCT__ "MatAYPX" 251 /*@ 252 MatAYPX - Computes Y = a*Y + X. 253 254 Logically on Mat 255 256 Input Parameters: 257 + a - the PetscScalar multiplier 258 . Y - the first matrix 259 . X - the second matrix 260 - str - either SAME_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN or SUBSET_NONZERO_PATTERN 261 262 Level: intermediate 263 264 .keywords: matrix, add 265 266 .seealso: MatAXPY() 267 @*/ 268 PetscErrorCode MatAYPX(Mat Y,PetscScalar a,Mat X,MatStructure str) 269 { 270 PetscScalar one = 1.0; 271 PetscErrorCode ierr; 272 PetscInt mX,mY,nX,nY; 273 274 PetscFunctionBegin; 275 PetscValidHeaderSpecific(X,MAT_CLASSID,3); 276 PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 277 PetscValidLogicalCollectiveScalar(Y,a,2); 278 ierr = MatGetSize(X,&mX,&nX);CHKERRQ(ierr); 279 ierr = MatGetSize(X,&mY,&nY);CHKERRQ(ierr); 280 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); 281 282 ierr = MatScale(Y,a);CHKERRQ(ierr); 283 ierr = MatAXPY(Y,one,X,str);CHKERRQ(ierr); 284 PetscFunctionReturn(0); 285 } 286 287 #undef __FUNCT__ 288 #define __FUNCT__ "MatComputeExplicitOperator" 289 /*@ 290 MatComputeExplicitOperator - Computes the explicit matrix 291 292 Collective on Mat 293 294 Input Parameter: 295 . inmat - the matrix 296 297 Output Parameter: 298 . mat - the explict preconditioned operator 299 300 Notes: 301 This computation is done by applying the operators to columns of the 302 identity matrix. 303 304 Currently, this routine uses a dense matrix format when 1 processor 305 is used and a sparse format otherwise. This routine is costly in general, 306 and is recommended for use only with relatively small systems. 307 308 Level: advanced 309 310 .keywords: Mat, compute, explicit, operator 311 @*/ 312 PetscErrorCode MatComputeExplicitOperator(Mat inmat,Mat *mat) 313 { 314 Vec in,out; 315 PetscErrorCode ierr; 316 PetscInt i,m,n,M,N,*rows,start,end; 317 MPI_Comm comm; 318 PetscScalar *array,zero = 0.0,one = 1.0; 319 PetscMPIInt size; 320 321 PetscFunctionBegin; 322 PetscValidHeaderSpecific(inmat,MAT_CLASSID,1); 323 PetscValidPointer(mat,2); 324 325 ierr = PetscObjectGetComm((PetscObject)inmat,&comm);CHKERRQ(ierr); 326 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 327 328 ierr = MatGetLocalSize(inmat,&m,&n);CHKERRQ(ierr); 329 ierr = MatGetSize(inmat,&M,&N);CHKERRQ(ierr); 330 ierr = MatCreateVecs(inmat,&in,&out);CHKERRQ(ierr); 331 ierr = VecSetOption(in,VEC_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 332 ierr = VecGetOwnershipRange(out,&start,&end);CHKERRQ(ierr); 333 ierr = PetscMalloc1(m,&rows);CHKERRQ(ierr); 334 for (i=0; i<m; i++) rows[i] = start + i; 335 336 ierr = MatCreate(comm,mat);CHKERRQ(ierr); 337 ierr = MatSetSizes(*mat,m,n,M,N);CHKERRQ(ierr); 338 if (size == 1) { 339 ierr = MatSetType(*mat,MATSEQDENSE);CHKERRQ(ierr); 340 ierr = MatSeqDenseSetPreallocation(*mat,NULL);CHKERRQ(ierr); 341 } else { 342 ierr = MatSetType(*mat,MATMPIAIJ);CHKERRQ(ierr); 343 ierr = MatMPIAIJSetPreallocation(*mat,n,NULL,N-n,NULL);CHKERRQ(ierr); 344 } 345 346 for (i=0; i<N; i++) { 347 348 ierr = VecSet(in,zero);CHKERRQ(ierr); 349 ierr = VecSetValues(in,1,&i,&one,INSERT_VALUES);CHKERRQ(ierr); 350 ierr = VecAssemblyBegin(in);CHKERRQ(ierr); 351 ierr = VecAssemblyEnd(in);CHKERRQ(ierr); 352 353 ierr = MatMult(inmat,in,out);CHKERRQ(ierr); 354 355 ierr = VecGetArray(out,&array);CHKERRQ(ierr); 356 ierr = MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);CHKERRQ(ierr); 357 ierr = VecRestoreArray(out,&array);CHKERRQ(ierr); 358 359 } 360 ierr = PetscFree(rows);CHKERRQ(ierr); 361 ierr = VecDestroy(&out);CHKERRQ(ierr); 362 ierr = VecDestroy(&in);CHKERRQ(ierr); 363 ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 364 ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 365 PetscFunctionReturn(0); 366 } 367 368 #undef __FUNCT__ 369 #define __FUNCT__ "MatChop" 370 /*@ 371 MatChop - Set all values in the matrix less than the tolerance to zero 372 373 Input Parameters: 374 + A - The matrix 375 - tol - The zero tolerance 376 377 Output Parameters: 378 . A - The chopped matrix 379 380 Level: intermediate 381 382 .seealso: MatCreate(), MatZeroEntries() 383 @*/ 384 PetscErrorCode MatChop(Mat A, PetscReal tol) 385 { 386 PetscScalar *newVals; 387 PetscInt *newCols; 388 PetscInt rStart, rEnd, numRows, maxRows, r, colMax = 0; 389 PetscErrorCode ierr; 390 391 PetscFunctionBegin; 392 ierr = MatGetOwnershipRange(A, &rStart, &rEnd);CHKERRQ(ierr); 393 for (r = rStart; r < rEnd; ++r) { 394 PetscInt ncols; 395 396 ierr = MatGetRow(A, r, &ncols, NULL, NULL);CHKERRQ(ierr); 397 colMax = PetscMax(colMax, ncols);CHKERRQ(ierr); 398 ierr = MatRestoreRow(A, r, &ncols, NULL, NULL);CHKERRQ(ierr); 399 } 400 numRows = rEnd - rStart; 401 ierr = MPIU_Allreduce(&numRows, &maxRows, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 402 ierr = PetscMalloc2(colMax,&newCols,colMax,&newVals);CHKERRQ(ierr); 403 for (r = rStart; r < rStart+maxRows; ++r) { 404 const PetscScalar *vals; 405 const PetscInt *cols; 406 PetscInt ncols, newcols, c; 407 408 if (r < rEnd) { 409 ierr = MatGetRow(A, r, &ncols, &cols, &vals);CHKERRQ(ierr); 410 for (c = 0; c < ncols; ++c) { 411 newCols[c] = cols[c]; 412 newVals[c] = PetscAbsScalar(vals[c]) < tol ? 0.0 : vals[c]; 413 } 414 newcols = ncols; 415 ierr = MatRestoreRow(A, r, &ncols, &cols, &vals);CHKERRQ(ierr); 416 ierr = MatSetValues(A, 1, &r, newcols, newCols, newVals, INSERT_VALUES);CHKERRQ(ierr); 417 } 418 ierr = MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 419 ierr = MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 420 } 421 ierr = PetscFree2(newCols,newVals);CHKERRQ(ierr); 422 PetscFunctionReturn(0); 423 } 424