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 .keywords: matrix, add, shift 151 152 .seealso: MatDiagonalSet() 153 @*/ 154 PetscErrorCode MatShift(Mat Y,PetscScalar a) 155 { 156 PetscErrorCode ierr; 157 158 PetscFunctionBegin; 159 PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 160 if (!Y->assembled) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 161 if (Y->factortype) SETERRQ(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 162 MatCheckPreallocated(Y,1); 163 164 ierr = (*Y->ops->shift)(Y,a);CHKERRQ(ierr); 165 166 #if defined(PETSC_HAVE_CUSP) 167 if (Y->valid_GPU_matrix != PETSC_CUSP_UNALLOCATED) { 168 Y->valid_GPU_matrix = PETSC_CUSP_CPU; 169 } 170 #endif 171 #if defined(PETSC_HAVE_VIENNACL) 172 if (Y->valid_GPU_matrix != PETSC_VIENNACL_UNALLOCATED) { 173 Y->valid_GPU_matrix = PETSC_VIENNACL_CPU; 174 } 175 #endif 176 PetscFunctionReturn(0); 177 } 178 179 #undef __FUNCT__ 180 #define __FUNCT__ "MatDiagonalSet_Default" 181 PetscErrorCode MatDiagonalSet_Default(Mat Y,Vec D,InsertMode is) 182 { 183 PetscErrorCode ierr; 184 PetscInt i,start,end; 185 PetscScalar *v; 186 187 PetscFunctionBegin; 188 ierr = MatGetOwnershipRange(Y,&start,&end);CHKERRQ(ierr); 189 ierr = VecGetArray(D,&v);CHKERRQ(ierr); 190 for (i=start; i<end; i++) { 191 ierr = MatSetValues(Y,1,&i,1,&i,v+i-start,is);CHKERRQ(ierr); 192 } 193 ierr = VecRestoreArray(D,&v);CHKERRQ(ierr); 194 ierr = MatAssemblyBegin(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 195 ierr = MatAssemblyEnd(Y,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 196 PetscFunctionReturn(0); 197 } 198 199 #undef __FUNCT__ 200 #define __FUNCT__ "MatDiagonalSet" 201 /*@ 202 MatDiagonalSet - Computes Y = Y + D, where D is a diagonal matrix 203 that is represented as a vector. Or Y[i,i] = D[i] if InsertMode is 204 INSERT_VALUES. 205 206 Input Parameters: 207 + Y - the input matrix 208 . D - the diagonal matrix, represented as a vector 209 - i - INSERT_VALUES or ADD_VALUES 210 211 Neighbor-wise Collective on Mat and Vec 212 213 Level: intermediate 214 215 .keywords: matrix, add, shift, diagonal 216 217 .seealso: MatShift() 218 @*/ 219 PetscErrorCode MatDiagonalSet(Mat Y,Vec D,InsertMode is) 220 { 221 PetscErrorCode ierr; 222 PetscInt matlocal,veclocal; 223 224 PetscFunctionBegin; 225 PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 226 PetscValidHeaderSpecific(D,VEC_CLASSID,2); 227 ierr = MatGetLocalSize(Y,&matlocal,NULL);CHKERRQ(ierr); 228 ierr = VecGetLocalSize(D,&veclocal);CHKERRQ(ierr); 229 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); 230 if (Y->ops->diagonalset) { 231 ierr = (*Y->ops->diagonalset)(Y,D,is);CHKERRQ(ierr); 232 } else { 233 ierr = MatDiagonalSet_Default(Y,D,is);CHKERRQ(ierr); 234 } 235 PetscFunctionReturn(0); 236 } 237 238 #undef __FUNCT__ 239 #define __FUNCT__ "MatAYPX" 240 /*@ 241 MatAYPX - Computes Y = a*Y + X. 242 243 Logically on Mat 244 245 Input Parameters: 246 + a - the PetscScalar multiplier 247 . Y - the first matrix 248 . X - the second matrix 249 - str - either SAME_NONZERO_PATTERN, DIFFERENT_NONZERO_PATTERN or SUBSET_NONZERO_PATTERN 250 251 Level: intermediate 252 253 .keywords: matrix, add 254 255 .seealso: MatAXPY() 256 @*/ 257 PetscErrorCode MatAYPX(Mat Y,PetscScalar a,Mat X,MatStructure str) 258 { 259 PetscScalar one = 1.0; 260 PetscErrorCode ierr; 261 PetscInt mX,mY,nX,nY; 262 263 PetscFunctionBegin; 264 PetscValidHeaderSpecific(X,MAT_CLASSID,3); 265 PetscValidHeaderSpecific(Y,MAT_CLASSID,1); 266 PetscValidLogicalCollectiveScalar(Y,a,2); 267 ierr = MatGetSize(X,&mX,&nX);CHKERRQ(ierr); 268 ierr = MatGetSize(X,&mY,&nY);CHKERRQ(ierr); 269 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); 270 271 ierr = MatScale(Y,a);CHKERRQ(ierr); 272 ierr = MatAXPY(Y,one,X,str);CHKERRQ(ierr); 273 PetscFunctionReturn(0); 274 } 275 276 #undef __FUNCT__ 277 #define __FUNCT__ "MatComputeExplicitOperator" 278 /*@ 279 MatComputeExplicitOperator - Computes the explicit matrix 280 281 Collective on Mat 282 283 Input Parameter: 284 . inmat - the matrix 285 286 Output Parameter: 287 . mat - the explict preconditioned operator 288 289 Notes: 290 This computation is done by applying the operators to columns of the 291 identity matrix. 292 293 Currently, this routine uses a dense matrix format when 1 processor 294 is used and a sparse format otherwise. This routine is costly in general, 295 and is recommended for use only with relatively small systems. 296 297 Level: advanced 298 299 .keywords: Mat, compute, explicit, operator 300 @*/ 301 PetscErrorCode MatComputeExplicitOperator(Mat inmat,Mat *mat) 302 { 303 Vec in,out; 304 PetscErrorCode ierr; 305 PetscInt i,m,n,M,N,*rows,start,end; 306 MPI_Comm comm; 307 PetscScalar *array,zero = 0.0,one = 1.0; 308 PetscMPIInt size; 309 310 PetscFunctionBegin; 311 PetscValidHeaderSpecific(inmat,MAT_CLASSID,1); 312 PetscValidPointer(mat,2); 313 314 ierr = PetscObjectGetComm((PetscObject)inmat,&comm);CHKERRQ(ierr); 315 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 316 317 ierr = MatGetLocalSize(inmat,&m,&n);CHKERRQ(ierr); 318 ierr = MatGetSize(inmat,&M,&N);CHKERRQ(ierr); 319 ierr = MatCreateVecs(inmat,&in,&out);CHKERRQ(ierr); 320 ierr = VecSetOption(in,VEC_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 321 ierr = VecGetOwnershipRange(out,&start,&end);CHKERRQ(ierr); 322 ierr = PetscMalloc1(m,&rows);CHKERRQ(ierr); 323 for (i=0; i<m; i++) rows[i] = start + i; 324 325 ierr = MatCreate(comm,mat);CHKERRQ(ierr); 326 ierr = MatSetSizes(*mat,m,n,M,N);CHKERRQ(ierr); 327 if (size == 1) { 328 ierr = MatSetType(*mat,MATSEQDENSE);CHKERRQ(ierr); 329 ierr = MatSeqDenseSetPreallocation(*mat,NULL);CHKERRQ(ierr); 330 } else { 331 ierr = MatSetType(*mat,MATMPIAIJ);CHKERRQ(ierr); 332 ierr = MatMPIAIJSetPreallocation(*mat,n,NULL,N-n,NULL);CHKERRQ(ierr); 333 } 334 335 for (i=0; i<N; i++) { 336 337 ierr = VecSet(in,zero);CHKERRQ(ierr); 338 ierr = VecSetValues(in,1,&i,&one,INSERT_VALUES);CHKERRQ(ierr); 339 ierr = VecAssemblyBegin(in);CHKERRQ(ierr); 340 ierr = VecAssemblyEnd(in);CHKERRQ(ierr); 341 342 ierr = MatMult(inmat,in,out);CHKERRQ(ierr); 343 344 ierr = VecGetArray(out,&array);CHKERRQ(ierr); 345 ierr = MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);CHKERRQ(ierr); 346 ierr = VecRestoreArray(out,&array);CHKERRQ(ierr); 347 348 } 349 ierr = PetscFree(rows);CHKERRQ(ierr); 350 ierr = VecDestroy(&out);CHKERRQ(ierr); 351 ierr = VecDestroy(&in);CHKERRQ(ierr); 352 ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 353 ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 354 PetscFunctionReturn(0); 355 } 356 357 #undef __FUNCT__ 358 #define __FUNCT__ "MatChop" 359 /*@ 360 MatChop - Set all values in the matrix less than the tolerance to zero 361 362 Input Parameters: 363 + A - The matrix 364 - tol - The zero tolerance 365 366 Output Parameters: 367 . A - The chopped matrix 368 369 Level: intermediate 370 371 .seealso: MatCreate(), MatZeroEntries() 372 @*/ 373 PetscErrorCode MatChop(Mat A, PetscReal tol) 374 { 375 PetscScalar *newVals; 376 PetscInt *newCols; 377 PetscInt rStart, rEnd, numRows, maxRows, r, colMax = 0; 378 PetscErrorCode ierr; 379 380 PetscFunctionBegin; 381 ierr = MatGetOwnershipRange(A, &rStart, &rEnd);CHKERRQ(ierr); 382 for (r = rStart; r < rEnd; ++r) { 383 PetscInt ncols; 384 385 ierr = MatGetRow(A, r, &ncols, NULL, NULL);CHKERRQ(ierr); 386 colMax = PetscMax(colMax, ncols);CHKERRQ(ierr); 387 ierr = MatRestoreRow(A, r, &ncols, NULL, NULL);CHKERRQ(ierr); 388 } 389 numRows = rEnd - rStart; 390 ierr = MPI_Allreduce(&numRows, &maxRows, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 391 ierr = PetscMalloc2(colMax,&newCols,colMax,&newVals);CHKERRQ(ierr); 392 for (r = rStart; r < rStart+maxRows; ++r) { 393 const PetscScalar *vals; 394 const PetscInt *cols; 395 PetscInt ncols, newcols, c; 396 397 if (r < rEnd) { 398 ierr = MatGetRow(A, r, &ncols, &cols, &vals);CHKERRQ(ierr); 399 for (c = 0; c < ncols; ++c) { 400 newCols[c] = cols[c]; 401 newVals[c] = PetscAbsScalar(vals[c]) < tol ? 0.0 : vals[c]; 402 } 403 newcols = ncols; 404 ierr = MatRestoreRow(A, r, &ncols, &cols, &vals);CHKERRQ(ierr); 405 ierr = MatSetValues(A, 1, &r, newcols, newCols, newVals, INSERT_VALUES);CHKERRQ(ierr); 406 } 407 ierr = MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 408 ierr = MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 409 } 410 ierr = PetscFree2(newCols,newVals);CHKERRQ(ierr); 411 PetscFunctionReturn(0); 412 } 413