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