1 #ifdef PETSC_RCS_HEADER 2 static char vcid[] = "$Id: dense.c,v 1.173 1999/09/02 14:53:18 bsmith Exp bsmith $"; 3 #endif 4 /* 5 Defines the basic matrix operations for sequential dense. 6 */ 7 8 #include "src/mat/impls/dense/seq/dense.h" 9 #include "pinclude/blaslapack.h" 10 11 #undef __FUNC__ 12 #define __FUNC__ "MatAXPY_SeqDense" 13 int MatAXPY_SeqDense(Scalar *alpha,Mat X,Mat Y) 14 { 15 Mat_SeqDense *x = (Mat_SeqDense*) X->data,*y = (Mat_SeqDense*) Y->data; 16 int N = x->m*x->n, one = 1; 17 18 PetscFunctionBegin; 19 BLaxpy_( &N, alpha, x->v, &one, y->v, &one ); 20 PLogFlops(2*N-1); 21 PetscFunctionReturn(0); 22 } 23 24 #undef __FUNC__ 25 #define __FUNC__ "MatGetInfo_SeqDense" 26 int MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info) 27 { 28 Mat_SeqDense *a = (Mat_SeqDense *) A->data; 29 int i,N = a->m*a->n,count = 0; 30 Scalar *v = a->v; 31 32 PetscFunctionBegin; 33 for ( i=0; i<N; i++ ) {if (*v != 0.0) count++; v++;} 34 35 info->rows_global = (double)a->m; 36 info->columns_global = (double)a->n; 37 info->rows_local = (double)a->m; 38 info->columns_local = (double)a->n; 39 info->block_size = 1.0; 40 info->nz_allocated = (double)N; 41 info->nz_used = (double)count; 42 info->nz_unneeded = (double)(N-count); 43 info->assemblies = (double)A->num_ass; 44 info->mallocs = 0; 45 info->memory = A->mem; 46 info->fill_ratio_given = 0; 47 info->fill_ratio_needed = 0; 48 info->factor_mallocs = 0; 49 50 PetscFunctionReturn(0); 51 } 52 53 #undef __FUNC__ 54 #define __FUNC__ "MatScale_SeqDense" 55 int MatScale_SeqDense(Scalar *alpha,Mat inA) 56 { 57 Mat_SeqDense *a = (Mat_SeqDense *) inA->data; 58 int one = 1, nz; 59 60 PetscFunctionBegin; 61 nz = a->m*a->n; 62 BLscal_( &nz, alpha, a->v, &one ); 63 PLogFlops(nz); 64 PetscFunctionReturn(0); 65 } 66 67 /* ---------------------------------------------------------------*/ 68 /* COMMENT: I have chosen to hide column permutation in the pivots, 69 rather than put it in the Mat->col slot.*/ 70 #undef __FUNC__ 71 #define __FUNC__ "MatLUFactor_SeqDense" 72 int MatLUFactor_SeqDense(Mat A,IS row,IS col,double f) 73 { 74 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 75 int info; 76 77 PetscFunctionBegin; 78 if (!mat->pivots) { 79 mat->pivots = (int *) PetscMalloc((mat->m+1)*sizeof(int));CHKPTRQ(mat->pivots); 80 PLogObjectMemory(A,mat->m*sizeof(int)); 81 } 82 if (!mat->m || !mat->n) PetscFunctionReturn(0); 83 LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info); 84 if (info<0) SETERRQ(PETSC_ERR_LIB,0,"Bad argument to LU factorization"); 85 if (info>0) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,0,"Bad LU factorization"); 86 A->factor = FACTOR_LU; 87 PLogFlops((2*mat->n*mat->n*mat->n)/3); 88 PetscFunctionReturn(0); 89 } 90 91 #undef __FUNC__ 92 #define __FUNC__ "MatDuplicate_SeqDense" 93 int MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat) 94 { 95 Mat_SeqDense *mat = (Mat_SeqDense *) A->data, *l; 96 int ierr; 97 Mat newi; 98 99 PetscFunctionBegin; 100 ierr = MatCreateSeqDense(A->comm,mat->m,mat->n,PETSC_NULL,&newi);CHKERRQ(ierr); 101 l = (Mat_SeqDense *) newi->data; 102 if (cpvalues == MAT_COPY_VALUES) { 103 ierr = PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar));CHKERRQ(ierr); 104 } 105 newi->assembled = PETSC_TRUE; 106 *newmat = newi; 107 PetscFunctionReturn(0); 108 } 109 110 #undef __FUNC__ 111 #define __FUNC__ "MatLUFactorSymbolic_SeqDense" 112 int MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,double f,Mat *fact) 113 { 114 int ierr; 115 116 PetscFunctionBegin; 117 ierr = MatDuplicate_SeqDense(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr); 118 PetscFunctionReturn(0); 119 } 120 121 #undef __FUNC__ 122 #define __FUNC__ "MatLUFactorNumeric_SeqDense" 123 int MatLUFactorNumeric_SeqDense(Mat A,Mat *fact) 124 { 125 Mat_SeqDense *mat = (Mat_SeqDense*) A->data, *l = (Mat_SeqDense*) (*fact)->data; 126 int ierr; 127 128 PetscFunctionBegin; 129 /* copy the numerical values */ 130 ierr = PetscMemcpy(l->v,mat->v,mat->m*mat->n*sizeof(Scalar));CHKERRQ(ierr); 131 (*fact)->factor = 0; 132 ierr = MatLUFactor(*fact,0,0,1.0);CHKERRQ(ierr); 133 PetscFunctionReturn(0); 134 } 135 136 #undef __FUNC__ 137 #define __FUNC__ "MatCholeskyFactorSymbolic_SeqDense" 138 int MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,double f,Mat *fact) 139 { 140 int ierr; 141 142 PetscFunctionBegin; 143 ierr = MatConvert(A,MATSAME,fact);CHKERRQ(ierr); 144 PetscFunctionReturn(0); 145 } 146 147 #undef __FUNC__ 148 #define __FUNC__ "MatCholeskyFactorNumeric_SeqDense" 149 int MatCholeskyFactorNumeric_SeqDense(Mat A,Mat *fact) 150 { 151 int ierr; 152 153 PetscFunctionBegin; 154 ierr = MatCholeskyFactor(*fact,0,1.0);CHKERRQ(ierr); 155 PetscFunctionReturn(0); 156 } 157 158 #undef __FUNC__ 159 #define __FUNC__ "MatCholeskyFactor_SeqDense" 160 int MatCholeskyFactor_SeqDense(Mat A,IS perm,double f) 161 { 162 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 163 int info,ierr; 164 165 PetscFunctionBegin; 166 if (mat->pivots) { 167 ierr = PetscFree(mat->pivots);CHKERRQ(ierr); 168 PLogObjectMemory(A,-mat->m*sizeof(int)); 169 mat->pivots = 0; 170 } 171 if (!mat->m || !mat->n) PetscFunctionReturn(0); 172 LApotrf_("L",&mat->n,mat->v,&mat->m,&info); 173 if (info) SETERRQ(PETSC_ERR_MAT_CH_ZRPVT,0,"Bad factorization"); 174 A->factor = FACTOR_CHOLESKY; 175 PLogFlops((mat->n*mat->n*mat->n)/3); 176 PetscFunctionReturn(0); 177 } 178 179 #undef __FUNC__ 180 #define __FUNC__ "MatSolve_SeqDense" 181 int MatSolve_SeqDense(Mat A,Vec xx,Vec yy) 182 { 183 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 184 int one = 1, info, ierr; 185 Scalar *x, *y; 186 187 PetscFunctionBegin; 188 if (!mat->m || !mat->n) PetscFunctionReturn(0); 189 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 190 ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 191 ierr = PetscMemcpy(y,x,mat->m*sizeof(Scalar));CHKERRQ(ierr); 192 if (A->factor == FACTOR_LU) { 193 LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 194 } else if (A->factor == FACTOR_CHOLESKY){ 195 LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 196 } 197 else SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Matrix must be factored to solve"); 198 if (info) SETERRQ(PETSC_ERR_LIB,0,"MBad solve"); 199 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 200 ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 201 PLogFlops(mat->n*mat->n - mat->n); 202 PetscFunctionReturn(0); 203 } 204 205 #undef __FUNC__ 206 #define __FUNC__ "MatSolveTrans_SeqDense" 207 int MatSolveTrans_SeqDense(Mat A,Vec xx,Vec yy) 208 { 209 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 210 int ierr,one = 1, info; 211 Scalar *x, *y; 212 213 PetscFunctionBegin; 214 if (!mat->m || !mat->n) PetscFunctionReturn(0); 215 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 216 ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 217 ierr = PetscMemcpy(y,x,mat->m*sizeof(Scalar));CHKERRQ(ierr); 218 /* assume if pivots exist then use LU; else Cholesky */ 219 if (mat->pivots) { 220 LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 221 } else { 222 LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 223 } 224 if (info) SETERRQ(PETSC_ERR_LIB,0,"Bad solve"); 225 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 226 ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 227 PLogFlops(mat->n*mat->n - mat->n); 228 PetscFunctionReturn(0); 229 } 230 231 #undef __FUNC__ 232 #define __FUNC__ "MatSolveAdd_SeqDense" 233 int MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 234 { 235 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 236 int ierr,one = 1, info; 237 Scalar *x, *y, sone = 1.0; 238 Vec tmp = 0; 239 240 PetscFunctionBegin; 241 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 242 ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 243 if (!mat->m || !mat->n) PetscFunctionReturn(0); 244 if (yy == zz) { 245 ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr); 246 PLogObjectParent(A,tmp); 247 ierr = VecCopy(yy,tmp);CHKERRQ(ierr); 248 } 249 ierr = PetscMemcpy(y,x,mat->m*sizeof(Scalar));CHKERRQ(ierr); 250 /* assume if pivots exist then use LU; else Cholesky */ 251 if (mat->pivots) { 252 LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 253 } else { 254 LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 255 } 256 if (info) SETERRQ(PETSC_ERR_LIB,0,"Bad solve"); 257 if (tmp) {ierr = VecAXPY(&sone,tmp,yy);CHKERRQ(ierr); ierr = VecDestroy(tmp);CHKERRQ(ierr);} 258 else {ierr = VecAXPY(&sone,zz,yy);CHKERRQ(ierr);} 259 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 260 ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 261 PLogFlops(mat->n*mat->n - mat->n); 262 PetscFunctionReturn(0); 263 } 264 265 #undef __FUNC__ 266 #define __FUNC__ "MatSolveTransAdd_SeqDense" 267 int MatSolveTransAdd_SeqDense(Mat A,Vec xx,Vec zz, Vec yy) 268 { 269 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 270 int one = 1, info,ierr; 271 Scalar *x, *y, sone = 1.0; 272 Vec tmp; 273 274 PetscFunctionBegin; 275 if (!mat->m || !mat->n) PetscFunctionReturn(0); 276 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 277 ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 278 if (yy == zz) { 279 ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr); 280 PLogObjectParent(A,tmp); 281 ierr = VecCopy(yy,tmp);CHKERRQ(ierr); 282 } 283 ierr = PetscMemcpy(y,x,mat->m*sizeof(Scalar));CHKERRQ(ierr); 284 /* assume if pivots exist then use LU; else Cholesky */ 285 if (mat->pivots) { 286 LAgetrs_( "T", &mat->m, &one, mat->v, &mat->m, mat->pivots,y, &mat->m, &info ); 287 } else { 288 LApotrs_( "L", &mat->m, &one, mat->v, &mat->m,y, &mat->m, &info ); 289 } 290 if (info) SETERRQ(PETSC_ERR_LIB,0,"Bad solve"); 291 if (tmp) { 292 ierr = VecAXPY(&sone,tmp,yy);CHKERRQ(ierr); 293 ierr = VecDestroy(tmp);CHKERRQ(ierr); 294 } else { 295 ierr = VecAXPY(&sone,zz,yy);CHKERRQ(ierr); 296 } 297 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 298 ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 299 PLogFlops(mat->n*mat->n - mat->n); 300 PetscFunctionReturn(0); 301 } 302 /* ------------------------------------------------------------------*/ 303 #undef __FUNC__ 304 #define __FUNC__ "MatRelax_SeqDense" 305 int MatRelax_SeqDense(Mat A,Vec bb,double omega,MatSORType flag, 306 double shift,int its,Vec xx) 307 { 308 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 309 Scalar *x, *b, *v = mat->v, zero = 0.0, xt; 310 int ierr, m = mat->m, i; 311 #if !defined(PETSC_USE_COMPLEX) 312 int o = 1; 313 #endif 314 315 PetscFunctionBegin; 316 if (flag & SOR_ZERO_INITIAL_GUESS) { 317 /* this is a hack fix, should have another version without the second BLdot */ 318 ierr = VecSet(&zero,xx);CHKERRQ(ierr); 319 } 320 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 321 ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 322 while (its--) { 323 if (flag & SOR_FORWARD_SWEEP){ 324 for ( i=0; i<m; i++ ) { 325 #if defined(PETSC_USE_COMPLEX) 326 /* cannot use BLAS dot for complex because compiler/linker is 327 not happy about returning a double complex */ 328 int _i; 329 Scalar sum = b[i]; 330 for ( _i=0; _i<m; _i++ ) { 331 sum -= PetscConj(v[i+_i*m])*x[_i]; 332 } 333 xt = sum; 334 #else 335 xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 336 #endif 337 x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 338 } 339 } 340 if (flag & SOR_BACKWARD_SWEEP) { 341 for ( i=m-1; i>=0; i-- ) { 342 #if defined(PETSC_USE_COMPLEX) 343 /* cannot use BLAS dot for complex because compiler/linker is 344 not happy about returning a double complex */ 345 int _i; 346 Scalar sum = b[i]; 347 for ( _i=0; _i<m; _i++ ) { 348 sum -= PetscConj(v[i+_i*m])*x[_i]; 349 } 350 xt = sum; 351 #else 352 xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 353 #endif 354 x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 355 } 356 } 357 } 358 ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 359 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 360 PetscFunctionReturn(0); 361 } 362 363 /* -----------------------------------------------------------------*/ 364 #undef __FUNC__ 365 #define __FUNC__ "MatMultTrans_SeqDense" 366 int MatMultTrans_SeqDense(Mat A,Vec xx,Vec yy) 367 { 368 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 369 Scalar *v = mat->v, *x, *y; 370 int ierr,_One=1;Scalar _DOne=1.0, _DZero=0.0; 371 372 PetscFunctionBegin; 373 if (!mat->m || !mat->n) PetscFunctionReturn(0); 374 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 375 ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 376 LAgemv_("T",&(mat->m),&(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DZero,y,&_One); 377 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 378 ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 379 PLogFlops(2*mat->m*mat->n - mat->n); 380 PetscFunctionReturn(0); 381 } 382 383 #undef __FUNC__ 384 #define __FUNC__ "MatMult_SeqDense" 385 int MatMult_SeqDense(Mat A,Vec xx,Vec yy) 386 { 387 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 388 Scalar *v = mat->v, *x, *y; 389 int ierr,_One=1;Scalar _DOne=1.0, _DZero=0.0; 390 391 PetscFunctionBegin; 392 if (!mat->m || !mat->n) PetscFunctionReturn(0); 393 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 394 ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 395 LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DZero,y,&_One); 396 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 397 ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 398 PLogFlops(2*mat->m*mat->n - mat->m); 399 PetscFunctionReturn(0); 400 } 401 402 #undef __FUNC__ 403 #define __FUNC__ "MatMultAdd_SeqDense" 404 int MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 405 { 406 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 407 Scalar *v = mat->v, *x, *y; 408 int ierr,_One=1; Scalar _DOne=1.0; 409 410 PetscFunctionBegin; 411 if (!mat->m || !mat->n) PetscFunctionReturn(0); 412 if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 413 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 414 ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 415 LAgemv_( "N", &(mat->m), &(mat->n),&_DOne,v,&(mat->m),x,&_One,&_DOne,y,&_One); 416 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 417 ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 418 PLogFlops(2*mat->m*mat->n); 419 PetscFunctionReturn(0); 420 } 421 422 #undef __FUNC__ 423 #define __FUNC__ "MatMultTransAdd_SeqDense" 424 int MatMultTransAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 425 { 426 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 427 Scalar *v = mat->v, *x, *y; 428 int ierr,_One=1; 429 Scalar _DOne=1.0; 430 431 PetscFunctionBegin; 432 if (!mat->m || !mat->n) PetscFunctionReturn(0); 433 if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 434 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 435 ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 436 LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m),x,&_One,&_DOne,y,&_One); 437 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 438 ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 439 PLogFlops(2*mat->m*mat->n); 440 PetscFunctionReturn(0); 441 } 442 443 /* -----------------------------------------------------------------*/ 444 #undef __FUNC__ 445 #define __FUNC__ "MatGetRow_SeqDense" 446 int MatGetRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals) 447 { 448 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 449 Scalar *v; 450 int i; 451 452 PetscFunctionBegin; 453 *ncols = mat->n; 454 if (cols) { 455 *cols = (int *) PetscMalloc((mat->n+1)*sizeof(int));CHKPTRQ(*cols); 456 for ( i=0; i<mat->n; i++ ) (*cols)[i] = i; 457 } 458 if (vals) { 459 *vals = (Scalar *) PetscMalloc((mat->n+1)*sizeof(Scalar));CHKPTRQ(*vals); 460 v = mat->v + row; 461 for ( i=0; i<mat->n; i++ ) {(*vals)[i] = *v; v += mat->m;} 462 } 463 PetscFunctionReturn(0); 464 } 465 466 #undef __FUNC__ 467 #define __FUNC__ "MatRestoreRow_SeqDense" 468 int MatRestoreRow_SeqDense(Mat A,int row,int *ncols,int **cols,Scalar **vals) 469 { 470 int ierr; 471 PetscFunctionBegin; 472 if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);} 473 if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); } 474 PetscFunctionReturn(0); 475 } 476 /* ----------------------------------------------------------------*/ 477 #undef __FUNC__ 478 #define __FUNC__ "MatSetValues_SeqDense" 479 int MatSetValues_SeqDense(Mat A,int m,int *indexm,int n, 480 int *indexn,Scalar *v,InsertMode addv) 481 { 482 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 483 int i,j; 484 485 PetscFunctionBegin; 486 if (!mat->roworiented) { 487 if (addv == INSERT_VALUES) { 488 for ( j=0; j<n; j++ ) { 489 if (indexn[j] < 0) {v += m; continue;} 490 #if defined(PETSC_USE_BOPT_g) 491 if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 492 #endif 493 for ( i=0; i<m; i++ ) { 494 if (indexm[i] < 0) {v++; continue;} 495 #if defined(PETSC_USE_BOPT_g) 496 if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 497 #endif 498 mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 499 } 500 } 501 } else { 502 for ( j=0; j<n; j++ ) { 503 if (indexn[j] < 0) {v += m; continue;} 504 #if defined(PETSC_USE_BOPT_g) 505 if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 506 #endif 507 for ( i=0; i<m; i++ ) { 508 if (indexm[i] < 0) {v++; continue;} 509 #if defined(PETSC_USE_BOPT_g) 510 if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 511 #endif 512 mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 513 } 514 } 515 } 516 } else { 517 if (addv == INSERT_VALUES) { 518 for ( i=0; i<m; i++ ) { 519 if (indexm[i] < 0) { v += n; continue;} 520 #if defined(PETSC_USE_BOPT_g) 521 if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 522 #endif 523 for ( j=0; j<n; j++ ) { 524 if (indexn[j] < 0) { v++; continue;} 525 #if defined(PETSC_USE_BOPT_g) 526 if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 527 #endif 528 mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 529 } 530 } 531 } else { 532 for ( i=0; i<m; i++ ) { 533 if (indexm[i] < 0) { v += n; continue;} 534 #if defined(PETSC_USE_BOPT_g) 535 if (indexm[i] >= A->m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Row too large"); 536 #endif 537 for ( j=0; j<n; j++ ) { 538 if (indexn[j] < 0) { v++; continue;} 539 #if defined(PETSC_USE_BOPT_g) 540 if (indexn[j] >= A->n) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Column too large"); 541 #endif 542 mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 543 } 544 } 545 } 546 } 547 PetscFunctionReturn(0); 548 } 549 550 #undef __FUNC__ 551 #define __FUNC__ "MatGetValues_SeqDense" 552 int MatGetValues_SeqDense(Mat A,int m,int *indexm,int n,int *indexn,Scalar *v) 553 { 554 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 555 int i, j; 556 Scalar *vpt = v; 557 558 PetscFunctionBegin; 559 /* row-oriented output */ 560 for ( i=0; i<m; i++ ) { 561 for ( j=0; j<n; j++ ) { 562 *vpt++ = mat->v[indexn[j]*mat->m + indexm[i]]; 563 } 564 } 565 PetscFunctionReturn(0); 566 } 567 568 /* -----------------------------------------------------------------*/ 569 570 #include "sys.h" 571 572 #undef __FUNC__ 573 #define __FUNC__ "MatLoad_SeqDense" 574 int MatLoad_SeqDense(Viewer viewer,MatType type,Mat *A) 575 { 576 Mat_SeqDense *a; 577 Mat B; 578 int *scols, i, j, nz, ierr, fd, header[4], size; 579 int *rowlengths = 0, M, N, *cols; 580 Scalar *vals, *svals, *v,*w; 581 MPI_Comm comm = ((PetscObject)viewer)->comm; 582 583 PetscFunctionBegin; 584 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 585 if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"view must have one processor"); 586 ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 587 ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 588 if (header[0] != MAT_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,0,"Not matrix object"); 589 M = header[1]; N = header[2]; nz = header[3]; 590 591 if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */ 592 ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A);CHKERRQ(ierr); 593 B = *A; 594 a = (Mat_SeqDense *) B->data; 595 v = a->v; 596 /* Allocate some temp space to read in the values and then flip them 597 from row major to column major */ 598 w = (Scalar *) PetscMalloc((M*N+1)*sizeof(Scalar));CHKPTRQ(w); 599 /* read in nonzero values */ 600 ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr); 601 /* now flip the values and store them in the matrix*/ 602 for ( j=0; j<N; j++ ) { 603 for ( i=0; i<M; i++ ) { 604 *v++ =w[i*N+j]; 605 } 606 } 607 ierr = PetscFree(w);CHKERRQ(ierr); 608 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 609 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 610 } else { 611 /* read row lengths */ 612 rowlengths = (int*) PetscMalloc( (M+1)*sizeof(int) );CHKPTRQ(rowlengths); 613 ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 614 615 /* create our matrix */ 616 ierr = MatCreateSeqDense(comm,M,N,PETSC_NULL,A);CHKERRQ(ierr); 617 B = *A; 618 a = (Mat_SeqDense *) B->data; 619 v = a->v; 620 621 /* read column indices and nonzeros */ 622 cols = scols = (int *) PetscMalloc( (nz+1)*sizeof(int) );CHKPTRQ(cols); 623 ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr); 624 vals = svals = (Scalar *) PetscMalloc( (nz+1)*sizeof(Scalar) );CHKPTRQ(vals); 625 ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 626 627 /* insert into matrix */ 628 for ( i=0; i<M; i++ ) { 629 for ( j=0; j<rowlengths[i]; j++ ) v[i+M*scols[j]] = svals[j]; 630 svals += rowlengths[i]; scols += rowlengths[i]; 631 } 632 ierr = PetscFree(vals);CHKERRQ(ierr); 633 ierr = PetscFree(cols);CHKERRQ(ierr); 634 ierr = PetscFree(rowlengths);CHKERRQ(ierr); 635 636 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 637 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 638 } 639 PetscFunctionReturn(0); 640 } 641 642 #include "sys.h" 643 644 #undef __FUNC__ 645 #define __FUNC__ "MatView_SeqDense_ASCII" 646 static int MatView_SeqDense_ASCII(Mat A,Viewer viewer) 647 { 648 Mat_SeqDense *a = (Mat_SeqDense *) A->data; 649 int ierr, i, j, format; 650 FILE *fd; 651 char *outputname; 652 Scalar *v; 653 654 PetscFunctionBegin; 655 ierr = ViewerASCIIGetPointer(viewer,&fd);CHKERRQ(ierr); 656 ierr = ViewerGetOutputname(viewer,&outputname);CHKERRQ(ierr); 657 ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr); 658 if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) { 659 PetscFunctionReturn(0); /* do nothing for now */ 660 } 661 else if (format == VIEWER_FORMAT_ASCII_COMMON) { 662 for ( i=0; i<a->m; i++ ) { 663 v = a->v + i; 664 fprintf(fd,"row %d:",i); 665 for ( j=0; j<a->n; j++ ) { 666 #if defined(PETSC_USE_COMPLEX) 667 if (PetscReal(*v) != 0.0 && PetscImaginary(*v) != 0.0) fprintf(fd," %d %g + %g i",j,PetscReal(*v),PetscImaginary(*v)); 668 else if (PetscReal(*v)) fprintf(fd," %d %g ",j,PetscReal(*v)); 669 #else 670 if (*v) fprintf(fd," %d %g ",j,*v); 671 #endif 672 v += a->m; 673 } 674 fprintf(fd,"\n"); 675 } 676 } else { 677 #if defined(PETSC_USE_COMPLEX) 678 int allreal = 1; 679 /* determine if matrix has all real values */ 680 v = a->v; 681 for ( i=0; i<a->m*a->n; i++ ) { 682 if (PetscImaginary(v[i])) { allreal = 0; break ;} 683 } 684 #endif 685 for ( i=0; i<a->m; i++ ) { 686 v = a->v + i; 687 for ( j=0; j<a->n; j++ ) { 688 #if defined(PETSC_USE_COMPLEX) 689 if (allreal) { 690 fprintf(fd,"%6.4e ",PetscReal(*v)); v += a->m; 691 } else { 692 fprintf(fd,"%6.4e + %6.4e i ",PetscReal(*v),PetscImaginary(*v)); v += a->m; 693 } 694 #else 695 fprintf(fd,"%6.4e ",*v); v += a->m; 696 #endif 697 } 698 fprintf(fd,"\n"); 699 } 700 } 701 fflush(fd); 702 PetscFunctionReturn(0); 703 } 704 705 #undef __FUNC__ 706 #define __FUNC__ "MatView_SeqDense_Binary" 707 static int MatView_SeqDense_Binary(Mat A,Viewer viewer) 708 { 709 Mat_SeqDense *a = (Mat_SeqDense *) A->data; 710 int ict, j, n = a->n, m = a->m, i, fd, *col_lens, ierr, nz = m*n; 711 int format; 712 Scalar *v, *anonz,*vals; 713 714 PetscFunctionBegin; 715 ierr = ViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 716 717 ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr); 718 if (format == VIEWER_FORMAT_BINARY_NATIVE) { 719 /* store the matrix as a dense matrix */ 720 col_lens = (int *) PetscMalloc( 4*sizeof(int) );CHKPTRQ(col_lens); 721 col_lens[0] = MAT_COOKIE; 722 col_lens[1] = m; 723 col_lens[2] = n; 724 col_lens[3] = MATRIX_BINARY_FORMAT_DENSE; 725 ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,1);CHKERRQ(ierr); 726 ierr = PetscFree(col_lens);CHKERRQ(ierr); 727 728 /* write out matrix, by rows */ 729 vals = (Scalar *) PetscMalloc((m*n+1)*sizeof(Scalar));CHKPTRQ(vals); 730 v = a->v; 731 for ( i=0; i<m; i++ ) { 732 for ( j=0; j<n; j++ ) { 733 vals[i + j*m] = *v++; 734 } 735 } 736 ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,0);CHKERRQ(ierr); 737 ierr = PetscFree(vals);CHKERRQ(ierr); 738 } else { 739 col_lens = (int *) PetscMalloc( (4+nz)*sizeof(int) );CHKPTRQ(col_lens); 740 col_lens[0] = MAT_COOKIE; 741 col_lens[1] = m; 742 col_lens[2] = n; 743 col_lens[3] = nz; 744 745 /* store lengths of each row and write (including header) to file */ 746 for ( i=0; i<m; i++ ) col_lens[4+i] = n; 747 ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,1);CHKERRQ(ierr); 748 749 /* Possibly should write in smaller increments, not whole matrix at once? */ 750 /* store column indices (zero start index) */ 751 ict = 0; 752 for ( i=0; i<m; i++ ) { 753 for ( j=0; j<n; j++ ) col_lens[ict++] = j; 754 } 755 ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,0);CHKERRQ(ierr); 756 ierr = PetscFree(col_lens);CHKERRQ(ierr); 757 758 /* store nonzero values */ 759 anonz = (Scalar *) PetscMalloc((nz+1)*sizeof(Scalar));CHKPTRQ(anonz); 760 ict = 0; 761 for ( i=0; i<m; i++ ) { 762 v = a->v + i; 763 for ( j=0; j<n; j++ ) { 764 anonz[ict++] = *v; v += a->m; 765 } 766 } 767 ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,0);CHKERRQ(ierr); 768 ierr = PetscFree(anonz);CHKERRQ(ierr); 769 } 770 PetscFunctionReturn(0); 771 } 772 773 #undef __FUNC__ 774 #define __FUNC__ "MatView_SeqDense" 775 int MatView_SeqDense(Mat A,Viewer viewer) 776 { 777 Mat_SeqDense *a = (Mat_SeqDense*) A->data; 778 int ierr; 779 780 PetscFunctionBegin; 781 if (PetscTypeCompare(viewer,SOCKET_VIEWER)) { 782 ierr = ViewerSocketPutScalar_Private(viewer,a->m,a->n,a->v);CHKERRQ(ierr); 783 } else if (PetscTypeCompare(viewer,ASCII_VIEWER)) { 784 ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr); 785 } else if (PetscTypeCompare(viewer,BINARY_VIEWER)) { 786 ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr); 787 } else { 788 SETERRQ(1,1,"Viewer type not supported by PETSc object"); 789 } 790 PetscFunctionReturn(0); 791 } 792 793 #undef __FUNC__ 794 #define __FUNC__ "MatDestroy_SeqDense" 795 int MatDestroy_SeqDense(Mat mat) 796 { 797 Mat_SeqDense *l = (Mat_SeqDense *) mat->data; 798 int ierr; 799 800 PetscFunctionBegin; 801 802 if (mat->mapping) { 803 ierr = ISLocalToGlobalMappingDestroy(mat->mapping);CHKERRQ(ierr); 804 } 805 if (mat->bmapping) { 806 ierr = ISLocalToGlobalMappingDestroy(mat->bmapping);CHKERRQ(ierr); 807 } 808 #if defined(PETSC_USE_LOG) 809 PLogObjectState((PetscObject)mat,"Rows %d Cols %d",l->m,l->n); 810 #endif 811 if (l->pivots) {ierr = PetscFree(l->pivots);CHKERRQ(ierr);} 812 if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);} 813 ierr = PetscFree(l);CHKERRQ(ierr); 814 if (mat->rmap) { 815 ierr = MapDestroy(mat->rmap);CHKERRQ(ierr); 816 } 817 if (mat->cmap) { 818 ierr = MapDestroy(mat->cmap);CHKERRQ(ierr); 819 } 820 PLogObjectDestroy(mat); 821 PetscHeaderDestroy(mat); 822 PetscFunctionReturn(0); 823 } 824 825 #undef __FUNC__ 826 #define __FUNC__ "MatTranspose_SeqDense" 827 int MatTranspose_SeqDense(Mat A,Mat *matout) 828 { 829 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 830 int k, j, m, n, ierr; 831 Scalar *v, tmp; 832 833 PetscFunctionBegin; 834 v = mat->v; m = mat->m; n = mat->n; 835 if (matout == PETSC_NULL) { /* in place transpose */ 836 if (m != n) { /* malloc temp to hold transpose */ 837 Scalar *w = (Scalar *) PetscMalloc((m+1)*(n+1)*sizeof(Scalar));CHKPTRQ(w); 838 for ( j=0; j<m; j++ ) { 839 for ( k=0; k<n; k++ ) { 840 w[k + j*n] = v[j + k*m]; 841 } 842 } 843 ierr = PetscMemcpy(v,w,m*n*sizeof(Scalar));CHKERRQ(ierr); 844 ierr = PetscFree(w);CHKERRQ(ierr); 845 } else { 846 for ( j=0; j<m; j++ ) { 847 for ( k=0; k<j; k++ ) { 848 tmp = v[j + k*n]; 849 v[j + k*n] = v[k + j*n]; 850 v[k + j*n] = tmp; 851 } 852 } 853 } 854 } else { /* out-of-place transpose */ 855 Mat tmat; 856 Mat_SeqDense *tmatd; 857 Scalar *v2; 858 ierr = MatCreateSeqDense(A->comm,mat->n,mat->m,PETSC_NULL,&tmat);CHKERRQ(ierr); 859 tmatd = (Mat_SeqDense *) tmat->data; 860 v = mat->v; v2 = tmatd->v; 861 for ( j=0; j<n; j++ ) { 862 for ( k=0; k<m; k++ ) v2[j + k*n] = v[k + j*m]; 863 } 864 ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 865 ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 866 *matout = tmat; 867 } 868 PetscFunctionReturn(0); 869 } 870 871 #undef __FUNC__ 872 #define __FUNC__ "MatEqual_SeqDense" 873 int MatEqual_SeqDense(Mat A1,Mat A2, PetscTruth *flg) 874 { 875 Mat_SeqDense *mat1 = (Mat_SeqDense *) A1->data; 876 Mat_SeqDense *mat2 = (Mat_SeqDense *) A2->data; 877 int i; 878 Scalar *v1 = mat1->v, *v2 = mat2->v; 879 880 PetscFunctionBegin; 881 if (A2->type != MATSEQDENSE) SETERRQ(PETSC_ERR_SUP,0,"Matrices should be of type MATSEQDENSE"); 882 if (mat1->m != mat2->m) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 883 if (mat1->n != mat2->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 884 for ( i=0; i<mat1->m*mat1->n; i++ ) { 885 if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 886 v1++; v2++; 887 } 888 *flg = PETSC_TRUE; 889 PetscFunctionReturn(0); 890 } 891 892 #undef __FUNC__ 893 #define __FUNC__ "MatGetDiagonal_SeqDense" 894 int MatGetDiagonal_SeqDense(Mat A,Vec v) 895 { 896 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 897 int ierr,i, n, len; 898 Scalar *x, zero = 0.0; 899 900 PetscFunctionBegin; 901 ierr = VecSet(&zero,v);CHKERRQ(ierr); 902 ierr = VecGetSize(v,&n);CHKERRQ(ierr); 903 ierr = VecGetArray(v,&x);CHKERRQ(ierr); 904 len = PetscMin(mat->m,mat->n); 905 if (n != mat->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Nonconforming mat and vec"); 906 for ( i=0; i<len; i++ ) { 907 x[i] = mat->v[i*mat->m + i]; 908 } 909 ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 910 PetscFunctionReturn(0); 911 } 912 913 #undef __FUNC__ 914 #define __FUNC__ "MatDiagonalScale_SeqDense" 915 int MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr) 916 { 917 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 918 Scalar *l,*r,x,*v; 919 int ierr,i,j,m = mat->m, n = mat->n; 920 921 PetscFunctionBegin; 922 if (ll) { 923 ierr = VecGetSize(ll,&m);CHKERRQ(ierr); 924 ierr = VecGetArray(ll,&l);CHKERRQ(ierr); 925 if (m != mat->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Left scaling vec wrong size"); 926 for ( i=0; i<m; i++ ) { 927 x = l[i]; 928 v = mat->v + i; 929 for ( j=0; j<n; j++ ) { (*v) *= x; v+= m;} 930 } 931 ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr); 932 PLogFlops(n*m); 933 } 934 if (rr) { 935 ierr = VecGetSize(rr,&n);CHKERRQ(ierr); 936 ierr = VecGetArray(rr,&r);CHKERRQ(ierr); 937 if (n != mat->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Right scaling vec wrong size"); 938 for ( i=0; i<n; i++ ) { 939 x = r[i]; 940 v = mat->v + i*m; 941 for ( j=0; j<m; j++ ) { (*v++) *= x;} 942 } 943 ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr); 944 PLogFlops(n*m); 945 } 946 PetscFunctionReturn(0); 947 } 948 949 #undef __FUNC__ 950 #define __FUNC__ "MatNorm_SeqDense" 951 int MatNorm_SeqDense(Mat A,NormType type,double *norm) 952 { 953 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 954 Scalar *v = mat->v; 955 double sum = 0.0; 956 int i, j; 957 958 PetscFunctionBegin; 959 if (type == NORM_FROBENIUS) { 960 for (i=0; i<mat->n*mat->m; i++ ) { 961 #if defined(PETSC_USE_COMPLEX) 962 sum += PetscReal(PetscConj(*v)*(*v)); v++; 963 #else 964 sum += (*v)*(*v); v++; 965 #endif 966 } 967 *norm = sqrt(sum); 968 PLogFlops(2*mat->n*mat->m); 969 } else if (type == NORM_1) { 970 *norm = 0.0; 971 for ( j=0; j<mat->n; j++ ) { 972 sum = 0.0; 973 for ( i=0; i<mat->m; i++ ) { 974 sum += PetscAbsScalar(*v); v++; 975 } 976 if (sum > *norm) *norm = sum; 977 } 978 PLogFlops(mat->n*mat->m); 979 } else if (type == NORM_INFINITY) { 980 *norm = 0.0; 981 for ( j=0; j<mat->m; j++ ) { 982 v = mat->v + j; 983 sum = 0.0; 984 for ( i=0; i<mat->n; i++ ) { 985 sum += PetscAbsScalar(*v); v += mat->m; 986 } 987 if (sum > *norm) *norm = sum; 988 } 989 PLogFlops(mat->n*mat->m); 990 } else { 991 SETERRQ(PETSC_ERR_SUP,0,"No two norm"); 992 } 993 PetscFunctionReturn(0); 994 } 995 996 #undef __FUNC__ 997 #define __FUNC__ "MatSetOption_SeqDense" 998 int MatSetOption_SeqDense(Mat A,MatOption op) 999 { 1000 Mat_SeqDense *aij = (Mat_SeqDense *) A->data; 1001 1002 PetscFunctionBegin; 1003 if (op == MAT_ROW_ORIENTED) aij->roworiented = 1; 1004 else if (op == MAT_COLUMN_ORIENTED) aij->roworiented = 0; 1005 else if (op == MAT_ROWS_SORTED || 1006 op == MAT_ROWS_UNSORTED || 1007 op == MAT_COLUMNS_SORTED || 1008 op == MAT_COLUMNS_UNSORTED || 1009 op == MAT_SYMMETRIC || 1010 op == MAT_STRUCTURALLY_SYMMETRIC || 1011 op == MAT_NO_NEW_NONZERO_LOCATIONS || 1012 op == MAT_YES_NEW_NONZERO_LOCATIONS || 1013 op == MAT_NEW_NONZERO_LOCATION_ERR || 1014 op == MAT_NO_NEW_DIAGONALS || 1015 op == MAT_YES_NEW_DIAGONALS || 1016 op == MAT_IGNORE_OFF_PROC_ENTRIES || 1017 op == MAT_USE_HASH_TABLE) 1018 PLogInfo(A,"MatSetOption_SeqDense:Option ignored\n"); 1019 else { 1020 SETERRQ(PETSC_ERR_SUP,0,"unknown option"); 1021 } 1022 PetscFunctionReturn(0); 1023 } 1024 1025 #undef __FUNC__ 1026 #define __FUNC__ "MatZeroEntries_SeqDense" 1027 int MatZeroEntries_SeqDense(Mat A) 1028 { 1029 Mat_SeqDense *l = (Mat_SeqDense *) A->data; 1030 int ierr; 1031 1032 PetscFunctionBegin; 1033 ierr = PetscMemzero(l->v,l->m*l->n*sizeof(Scalar));CHKERRQ(ierr); 1034 PetscFunctionReturn(0); 1035 } 1036 1037 #undef __FUNC__ 1038 #define __FUNC__ "MatGetBlockSize_SeqDense" 1039 int MatGetBlockSize_SeqDense(Mat A,int *bs) 1040 { 1041 PetscFunctionBegin; 1042 *bs = 1; 1043 PetscFunctionReturn(0); 1044 } 1045 1046 #undef __FUNC__ 1047 #define __FUNC__ "MatZeroRows_SeqDense" 1048 int MatZeroRows_SeqDense(Mat A,IS is,Scalar *diag) 1049 { 1050 Mat_SeqDense *l = (Mat_SeqDense *) A->data; 1051 int n = l->n, i, j,ierr,N, *rows; 1052 Scalar *slot; 1053 1054 PetscFunctionBegin; 1055 ierr = ISGetSize(is,&N);CHKERRQ(ierr); 1056 ierr = ISGetIndices(is,&rows);CHKERRQ(ierr); 1057 for ( i=0; i<N; i++ ) { 1058 slot = l->v + rows[i]; 1059 for ( j=0; j<n; j++ ) { *slot = 0.0; slot += n;} 1060 } 1061 if (diag) { 1062 for ( i=0; i<N; i++ ) { 1063 slot = l->v + (n+1)*rows[i]; 1064 *slot = *diag; 1065 } 1066 } 1067 ISRestoreIndices(is,&rows); 1068 PetscFunctionReturn(0); 1069 } 1070 1071 #undef __FUNC__ 1072 #define __FUNC__ "MatGetSize_SeqDense" 1073 int MatGetSize_SeqDense(Mat A,int *m,int *n) 1074 { 1075 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 1076 1077 PetscFunctionBegin; 1078 *m = mat->m; *n = mat->n; 1079 PetscFunctionReturn(0); 1080 } 1081 1082 #undef __FUNC__ 1083 #define __FUNC__ "MatGetOwnershipRange_SeqDense" 1084 int MatGetOwnershipRange_SeqDense(Mat A,int *m,int *n) 1085 { 1086 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 1087 1088 PetscFunctionBegin; 1089 *m = 0; *n = mat->m; 1090 PetscFunctionReturn(0); 1091 } 1092 1093 #undef __FUNC__ 1094 #define __FUNC__ "MatGetArray_SeqDense" 1095 int MatGetArray_SeqDense(Mat A,Scalar **array) 1096 { 1097 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 1098 1099 PetscFunctionBegin; 1100 *array = mat->v; 1101 PetscFunctionReturn(0); 1102 } 1103 1104 #undef __FUNC__ 1105 #define __FUNC__ "MatRestoreArray_SeqDense" 1106 int MatRestoreArray_SeqDense(Mat A,Scalar **array) 1107 { 1108 PetscFunctionBegin; 1109 *array = 0; /* user cannot accidently use the array later */ 1110 PetscFunctionReturn(0); 1111 } 1112 1113 #undef __FUNC__ 1114 #define __FUNC__ "MatGetSubMatrix_SeqDense" 1115 static int MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B) 1116 { 1117 Mat_SeqDense *mat = (Mat_SeqDense *) A->data; 1118 int i, j, ierr, m = mat->m, *irow, *icol, nrows, ncols; 1119 Scalar *av, *bv, *v = mat->v; 1120 Mat newmat; 1121 1122 PetscFunctionBegin; 1123 ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 1124 ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 1125 ierr = ISGetSize(isrow,&nrows);CHKERRQ(ierr); 1126 ierr = ISGetSize(iscol,&ncols);CHKERRQ(ierr); 1127 1128 /* Check submatrixcall */ 1129 if (scall == MAT_REUSE_MATRIX) { 1130 int n_cols,n_rows; 1131 ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr); 1132 if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Reused submatrix wrong size"); 1133 newmat = *B; 1134 } else { 1135 /* Create and fill new matrix */ 1136 ierr = MatCreateSeqDense(A->comm,nrows,ncols,PETSC_NULL,&newmat);CHKERRQ(ierr); 1137 } 1138 1139 /* Now extract the data pointers and do the copy, column at a time */ 1140 bv = ((Mat_SeqDense*)newmat->data)->v; 1141 1142 for ( i=0; i<ncols; i++ ) { 1143 av = v + m*icol[i]; 1144 for (j=0; j<nrows; j++ ) { 1145 *bv++ = av[irow[j]]; 1146 } 1147 } 1148 1149 /* Assemble the matrices so that the correct flags are set */ 1150 ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1151 ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1152 1153 /* Free work space */ 1154 ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 1155 ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 1156 *B = newmat; 1157 PetscFunctionReturn(0); 1158 } 1159 1160 #undef __FUNC__ 1161 #define __FUNC__ "MatGetSubMatrices_SeqDense" 1162 int MatGetSubMatrices_SeqDense(Mat A,int n, IS *irow,IS *icol,MatReuse scall,Mat **B) 1163 { 1164 int ierr,i; 1165 1166 PetscFunctionBegin; 1167 if (scall == MAT_INITIAL_MATRIX) { 1168 *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) );CHKPTRQ(*B); 1169 } 1170 1171 for ( i=0; i<n; i++ ) { 1172 ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr); 1173 } 1174 PetscFunctionReturn(0); 1175 } 1176 1177 #undef __FUNC__ 1178 #define __FUNC__ "MatCopy_SeqDense" 1179 int MatCopy_SeqDense(Mat A, Mat B,MatStructure str) 1180 { 1181 Mat_SeqDense *a = (Mat_SeqDense *) A->data, *b = (Mat_SeqDense *)B->data; 1182 int ierr; 1183 1184 PetscFunctionBegin; 1185 if (B->type != MATSEQDENSE) { 1186 ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr); 1187 PetscFunctionReturn(0); 1188 } 1189 if (a->m != b->m || a->n != b->n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"size(B) != size(A)"); 1190 ierr = PetscMemcpy(b->v,a->v,a->m*a->n*sizeof(Scalar));CHKERRQ(ierr); 1191 PetscFunctionReturn(0); 1192 } 1193 1194 /* -------------------------------------------------------------------*/ 1195 static struct _MatOps MatOps_Values = {MatSetValues_SeqDense, 1196 MatGetRow_SeqDense, 1197 MatRestoreRow_SeqDense, 1198 MatMult_SeqDense, 1199 MatMultAdd_SeqDense, 1200 MatMultTrans_SeqDense, 1201 MatMultTransAdd_SeqDense, 1202 MatSolve_SeqDense, 1203 MatSolveAdd_SeqDense, 1204 MatSolveTrans_SeqDense, 1205 MatSolveTransAdd_SeqDense, 1206 MatLUFactor_SeqDense, 1207 MatCholeskyFactor_SeqDense, 1208 MatRelax_SeqDense, 1209 MatTranspose_SeqDense, 1210 MatGetInfo_SeqDense, 1211 MatEqual_SeqDense, 1212 MatGetDiagonal_SeqDense, 1213 MatDiagonalScale_SeqDense, 1214 MatNorm_SeqDense, 1215 0, 1216 0, 1217 0, 1218 MatSetOption_SeqDense, 1219 MatZeroEntries_SeqDense, 1220 MatZeroRows_SeqDense, 1221 MatLUFactorSymbolic_SeqDense, 1222 MatLUFactorNumeric_SeqDense, 1223 MatCholeskyFactorSymbolic_SeqDense, 1224 MatCholeskyFactorNumeric_SeqDense, 1225 MatGetSize_SeqDense, 1226 MatGetSize_SeqDense, 1227 MatGetOwnershipRange_SeqDense, 1228 0, 1229 0, 1230 MatGetArray_SeqDense, 1231 MatRestoreArray_SeqDense, 1232 MatDuplicate_SeqDense, 1233 0, 1234 0, 1235 0, 1236 0, 1237 MatAXPY_SeqDense, 1238 MatGetSubMatrices_SeqDense, 1239 0, 1240 MatGetValues_SeqDense, 1241 MatCopy_SeqDense, 1242 0, 1243 MatScale_SeqDense, 1244 0, 1245 0, 1246 0, 1247 MatGetBlockSize_SeqDense, 1248 0, 1249 0, 1250 0, 1251 0, 1252 0, 1253 0, 1254 0, 1255 0, 1256 0, 1257 0, 1258 0, 1259 0, 1260 MatGetMaps_Petsc}; 1261 1262 #undef __FUNC__ 1263 #define __FUNC__ "MatCreateSeqDense" 1264 /*@C 1265 MatCreateSeqDense - Creates a sequential dense matrix that 1266 is stored in column major order (the usual Fortran 77 manner). Many 1267 of the matrix operations use the BLAS and LAPACK routines. 1268 1269 Collective on MPI_Comm 1270 1271 Input Parameters: 1272 + comm - MPI communicator, set to PETSC_COMM_SELF 1273 . m - number of rows 1274 . n - number of columns 1275 - data - optional location of matrix data. Set data=PETSC_NULL for PETSc 1276 to control all matrix memory allocation. 1277 1278 Output Parameter: 1279 . A - the matrix 1280 1281 Notes: 1282 The data input variable is intended primarily for Fortran programmers 1283 who wish to allocate their own matrix memory space. Most users should 1284 set data=PETSC_NULL. 1285 1286 Level: intermediate 1287 1288 .keywords: dense, matrix, LAPACK, BLAS 1289 1290 .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues() 1291 @*/ 1292 int MatCreateSeqDense(MPI_Comm comm,int m,int n,Scalar *data,Mat *A) 1293 { 1294 Mat B; 1295 Mat_SeqDense *b; 1296 int ierr,flg,size; 1297 1298 PetscFunctionBegin; 1299 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1300 if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Comm must be of size 1"); 1301 1302 *A = 0; 1303 PetscHeaderCreate(B,_p_Mat,struct _MatOps,MAT_COOKIE,MATSEQDENSE,"Mat",comm,MatDestroy,MatView); 1304 PLogObjectCreate(B); 1305 b = (Mat_SeqDense *) PetscMalloc(sizeof(Mat_SeqDense));CHKPTRQ(b); 1306 ierr = PetscMemzero(b,sizeof(Mat_SeqDense));CHKERRQ(ierr); 1307 ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 1308 B->ops->destroy = MatDestroy_SeqDense; 1309 B->ops->view = MatView_SeqDense; 1310 B->factor = 0; 1311 B->mapping = 0; 1312 PLogObjectMemory(B,sizeof(struct _p_Mat)); 1313 B->data = (void *) b; 1314 1315 b->m = m; B->m = m; B->M = m; 1316 b->n = n; B->n = n; B->N = n; 1317 1318 ierr = MapCreateMPI(comm,m,m,&B->rmap);CHKERRQ(ierr); 1319 ierr = MapCreateMPI(comm,n,n,&B->cmap);CHKERRQ(ierr); 1320 1321 b->pivots = 0; 1322 b->roworiented = 1; 1323 if (data == PETSC_NULL) { 1324 b->v = (Scalar*) PetscMalloc((m*n+1)*sizeof(Scalar));CHKPTRQ(b->v); 1325 ierr = PetscMemzero(b->v,m*n*sizeof(Scalar));CHKERRQ(ierr); 1326 b->user_alloc = 0; 1327 PLogObjectMemory(B,n*m*sizeof(Scalar)); 1328 } else { /* user-allocated storage */ 1329 b->v = data; 1330 b->user_alloc = 1; 1331 } 1332 ierr = OptionsHasName(PETSC_NULL,"-help",&flg);CHKERRQ(ierr); 1333 if (flg) {ierr = MatPrintHelp(B);CHKERRQ(ierr);} 1334 1335 *A = B; 1336 PetscFunctionReturn(0); 1337 } 1338 1339 1340 1341 1342 1343 1344