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