1 2 /* 3 Standard Fortran style matrices 4 */ 5 #include "ptscimpl.h" 6 #include "plapack.h" 7 #include "matimpl.h" 8 #include "math.h" 9 #include "vec/vecimpl.h" 10 11 typedef struct { 12 Scalar *v; 13 int roworiented; 14 int m,n,pad; 15 int *pivots; /* pivots in LU factorization */ 16 } MatiSD; 17 18 19 static int MatiSDnz(Mat matin,int *nz) 20 { 21 MatiSD *mat = (MatiSD *) matin->data; 22 int i,N = mat->m*mat->n,count = 0; 23 Scalar *v = mat->v; 24 for ( i=0; i<N; i++ ) {if (*v != 0.0) count++; v++;} 25 *nz = count; return 0; 26 } 27 static int MatiSDmemory(Mat matin,int *mem) 28 { 29 MatiSD *mat = (MatiSD *) matin->data; 30 *mem = mat->m*mat->n*sizeof(Scalar); return 0; 31 } 32 33 /* ---------------------------------------------------------------*/ 34 /* COMMENT: I have chosen to hide column permutation in the pivots, 35 rather than put it in the Mat->col slot.*/ 36 static int MatiSDlufactor(Mat matin,IS row,IS col) 37 { 38 MatiSD *mat = (MatiSD *) matin->data; 39 int ierr, one = 1, info; 40 if (!mat->pivots) { 41 mat->pivots = (int *) MALLOC( mat->m*sizeof(int) ); 42 CHKPTR(mat->pivots); 43 } 44 LAgetrf_(&mat->m,&mat->n,mat->v,&mat->m,mat->pivots,&info); 45 if (info) SETERR(1,"Bad LU factorization"); 46 matin->factor = FACTOR_LU; 47 return 0; 48 } 49 static int MatiSDlufactorsymbolic(Mat matin,IS row,IS col,Mat *fact) 50 { 51 int ierr; 52 if (ierr = MatCopy(matin,fact)) SETERR(ierr,0); 53 return 0; 54 } 55 static int MatiSDlufactornumeric(Mat matin,Mat fact) 56 { 57 return MatLUFactor(fact,0,0); 58 } 59 static int MatiSDchfactorsymbolic(Mat matin,IS row,Mat *fact) 60 { 61 int ierr; 62 if (ierr = MatCopy(matin,fact)) SETERR(ierr,0); 63 return 0; 64 } 65 static int MatiSDchfactornumeric(Mat matin,Mat fact) 66 { 67 return MatCholeskyFactor(fact,0); 68 } 69 static int MatiSDchfactor(Mat matin,IS perm) 70 { 71 MatiSD *mat = (MatiSD *) matin->data; 72 int ierr, one = 1, info; 73 if (mat->pivots) {FREE(mat->pivots); mat->pivots = 0;} 74 LApotrf_("L",&mat->n,mat->v,&mat->m,&info); 75 if (info) SETERR(1,"Bad Cholesky factorization"); 76 matin->factor = FACTOR_CHOLESKY; 77 return 0; 78 } 79 80 static int MatiSDsolve(Mat matin,Vec xx,Vec yy) 81 { 82 MatiSD *mat = (MatiSD *) matin->data; 83 int i,j, one = 1, info; 84 Scalar *v = mat->v, *x, *y; 85 VecGetArray(xx,&x); VecGetArray(yy,&y); 86 MEMCPY(y,x,mat->m*sizeof(Scalar)); 87 /* assume if pivots exist then LU else Cholesky */ 88 if (mat->pivots) { 89 LAgetrs_( "N", &mat->m, &one, mat->v, &mat->m, mat->pivots, 90 y, &mat->m, &info ); 91 } 92 else { 93 LApotrs_( "L", &mat->m, &one, mat->v, &mat->m, 94 y, &mat->m, &info ); 95 } 96 if (info) SETERR(1,"Bad solve"); 97 return 0; 98 } 99 static int MatiSDsolveadd(Mat matin,Vec xx,Vec yy,Vec zz) 100 { 101 MatiSD *mat = (MatiSD *) matin->data; 102 int i; 103 Scalar *y, *z; 104 VecGetArray(yy,&y); VecGetArray(zz,&z); 105 MatiSDsolve(matin,xx,zz); 106 for ( i=0; i<mat->m; i++ ) z[i] += y[i]; 107 return 0; 108 } 109 static int MatiSDsolvetrans(Mat matin,Vec xx,Vec yy) 110 {return 0;} 111 112 static int MatiSDsolvetransadd(Mat matin,Vec xx,Vec yy,Vec zz) 113 {return 0;} 114 115 /* ------------------------------------------------------------------*/ 116 static int MatiSDrelax(Mat matin,Vec bb,double omega,int flag,IS is, 117 int its,Vec xx) 118 { 119 MatiSD *mat = (MatiSD *) matin->data; 120 Scalar *x, *b, *v = mat->v, zero = 0.0, xt; 121 int o = 1, tmp,n = mat->n,ierr, m = mat->m, i, j; 122 123 if (is) SETERR(1,"No support for ordering in relaxation"); 124 if (flag & SOR_ZERO_INITIAL_GUESS) { 125 /* this is a hack fix, should have another version without 126 the second BLdot */ 127 if (ierr = VecSet(&zero,xx)) SETERR(ierr,0); 128 } 129 VecGetArray(xx,&x); VecGetArray(bb,&b); 130 while (its--) { 131 if (flag & SOR_FORWARD_SWEEP){ 132 for ( i=0; i<m; i++ ) { 133 xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 134 x[i] = (1. - omega)*x[i] + omega*(xt/v[i + i*m] + x[i]); 135 } 136 } 137 if (flag & SOR_BACKWARD_SWEEP) { 138 for ( i=m-1; i>=0; i-- ) { 139 xt = b[i]-BLdot_(&m,v+i,&m,x,&o); 140 x[i] = (1. - omega)*x[i] + omega*(xt/v[i + i*m] + x[i]); 141 } 142 } 143 } 144 return 0; 145 } 146 147 /* -----------------------------------------------------------------*/ 148 static int MatiSDmulttrans(Mat matin,Vec xx,Vec yy) 149 { 150 MatiSD *mat = (MatiSD *) matin->data; 151 Scalar *v = mat->v, *x, *y; 152 int _One=1;Scalar _DOne=1.0, _DZero=0.0; 153 VecGetArray(xx,&x), VecGetArray(yy,&y); 154 LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 155 x, &_One, &_DZero, y, &_One ); 156 return 0; 157 } 158 static int MatiSDmult(Mat matin,Vec xx,Vec yy) 159 { 160 MatiSD *mat = (MatiSD *) matin->data; 161 Scalar *v = mat->v, *x, *y; 162 int _One=1;Scalar _DOne=1.0, _DZero=0.0; 163 VecGetArray(xx,&x); VecGetArray(yy,&y); 164 LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 165 x, &_One, &_DZero, y, &_One ); 166 return 0; 167 } 168 static int MatiSDmultadd(Mat matin,Vec xx,Vec zz,Vec yy) 169 { 170 MatiSD *mat = (MatiSD *) matin->data; 171 Scalar *v = mat->v, *x, *y, *z; 172 int _One=1; Scalar _DOne=1.0, _DZero=0.0; 173 VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z); 174 if (zz != yy) MEMCPY(y,z,mat->m*sizeof(Scalar)); 175 LAgemv_( "N", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 176 x, &_One, &_DOne, y, &_One ); 177 return 0; 178 } 179 static int MatiSDmulttransadd(Mat matin,Vec xx,Vec zz,Vec yy) 180 { 181 MatiSD *mat = (MatiSD *) matin->data; 182 Scalar *v = mat->v, *x, *y, *z; 183 int _One=1; 184 Scalar _DOne=1.0, _DZero=0.0; 185 VecGetArray(xx,&x); VecGetArray(yy,&y); 186 VecGetArray(zz,&z); 187 if (zz != yy) MEMCPY(y,z,mat->m*sizeof(Scalar)); 188 LAgemv_( "T", &(mat->m), &(mat->n), &_DOne, v, &(mat->m), 189 x, &_One, &_DOne, y, &_One ); 190 return 0; 191 } 192 193 /* -----------------------------------------------------------------*/ 194 static int MatiSDgetrow(Mat matin,int row,int *ncols,int **cols, 195 Scalar **vals) 196 { 197 MatiSD *mat = (MatiSD *) matin->data; 198 Scalar *v; 199 int i; 200 *ncols = mat->n; 201 if (cols) { 202 *cols = (int *) MALLOC(mat->n*sizeof(int)); CHKPTR(*cols); 203 for ( i=0; i<mat->n; i++ ) *cols[i] = i; 204 } 205 if (vals) { 206 *vals = (Scalar *) MALLOC(mat->n*sizeof(Scalar)); CHKPTR(*vals); 207 v = mat->v + row; 208 for ( i=0; i<mat->n; i++ ) {*vals[i] = *v; v += mat->m;} 209 } 210 return 0; 211 } 212 static int MatiSDrestorerow(Mat matin,int row,int *ncols,int **cols, 213 Scalar **vals) 214 { 215 MatiSD *mat = (MatiSD *) matin->data; 216 if (cols) { FREE(*cols); } 217 if (vals) { FREE(*vals); } 218 return 0; 219 } 220 /* ----------------------------------------------------------------*/ 221 static int MatiSDinsert(Mat matin,Scalar *v,int m,int *indexm,int n, 222 int *indexn) 223 { 224 MatiSD *mat = (MatiSD *) matin->data; 225 int i,j; 226 if (!mat->roworiented) { 227 for ( j=0; j<n; j++ ) { 228 for ( i=0; i<m; i++ ) { 229 mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 230 } 231 } 232 } 233 else { 234 for ( i=0; i<m; i++ ) { 235 for ( j=0; j<n; j++ ) { 236 mat->v[indexn[j]*mat->m + indexm[i]] = *v++; 237 } 238 } 239 } 240 return 0; 241 } 242 static int MatiSDinsertadd(Mat matin,Scalar *v,int m,int *indexm, 243 int n,int *indexn) 244 { 245 MatiSD *mat = (MatiSD *) matin->data; 246 int i,j; 247 if (!mat->roworiented) { 248 for ( j=0; j<n; j++ ) { 249 for ( i=0; i<m; i++ ) { 250 mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 251 } 252 } 253 } 254 else { 255 for ( i=0; i<m; i++ ) { 256 for ( j=0; j<n; j++ ) { 257 mat->v[indexn[j]*mat->m + indexm[i]] += *v++; 258 } 259 } 260 } 261 return 0; 262 } 263 /* -----------------------------------------------------------------*/ 264 static int MatiSDcopy(Mat matin,Mat *newmat) 265 { 266 MatiSD *mat = (MatiSD *) matin->data; 267 int ierr; 268 Mat newi; 269 MatiSD *l; 270 if (ierr = MatCreateSequentialDense(mat->m,mat->n,&newi)) SETERR(ierr,0); 271 l = (MatiSD *) newi->data; 272 MEMCPY(l->v,mat->v,mat->m*mat->n*sizeof(Scalar)); 273 *newmat = newi; 274 return 0; 275 } 276 277 int MatiSDview(PetscObject obj,Viewer ptr) 278 { 279 Mat matin = (Mat) obj; 280 MatiSD *mat = (MatiSD *) matin->data; 281 Scalar *v; 282 int i,j; 283 for ( i=0; i<mat->m; i++ ) { 284 v = mat->v + i; 285 for ( j=0; j<mat->n; j++ ) { 286 #if defined(PETSC_COMPLEX) 287 printf("%6.4e + %6.4e i ",real(*v),imag(*v)); v += mat->m; 288 #else 289 printf("%6.4e ",*v); v += mat->m; 290 #endif 291 } 292 printf("\n"); 293 } 294 return 0; 295 } 296 297 298 static int MatiSDdestroy(PetscObject obj) 299 { 300 Mat mat = (Mat) obj; 301 MatiSD *l = (MatiSD *) mat->data; 302 if (l->pivots) FREE(l->pivots); 303 FREE(l); 304 FREE(mat); 305 return 0; 306 } 307 308 static int MatiSDtrans(Mat matin) 309 { 310 MatiSD *mat = (MatiSD *) matin->data; 311 int k,j; 312 Scalar *v = mat->v, tmp; 313 if (mat->m != mat->n) { 314 SETERR(1,"Cannot transpose rectangular dense matrix"); 315 } 316 for ( j=0; j<mat->m; j++ ) { 317 for ( k=0; k<j; k++ ) { 318 tmp = v[j + k*mat->n]; 319 v[j + k*mat->n] = v[k + j*mat->n]; 320 v[k + j*mat->n] = tmp; 321 } 322 } 323 return 0; 324 } 325 326 static int MatiSDequal(Mat matin1,Mat matin2) 327 { 328 MatiSD *mat1 = (MatiSD *) matin1->data; 329 MatiSD *mat2 = (MatiSD *) matin2->data; 330 int i; 331 Scalar *v1 = mat1->v, *v2 = mat2->v; 332 if (mat1->m != mat2->m) return 0; 333 if (mat1->n != mat2->n) return 0; 334 for ( i=0; i<mat1->m*mat1->n; i++ ) { 335 if (*v1 != *v2) return 0; 336 v1++; v2++; 337 } 338 return 1; 339 } 340 341 static int MatiSDgetdiag(Mat matin,Vec v) 342 { 343 MatiSD *mat = (MatiSD *) matin->data; 344 int i,j, n; 345 Scalar *x, zero = 0.0; 346 CHKTYPE(v,SEQVECTOR); 347 VecGetArray(v,&x); VecGetSize(v,&n); 348 if (n != mat->m) SETERR(1,"Nonconforming matrix and vector"); 349 for ( i=0; i<mat->m; i++ ) { 350 x[i] = mat->v[i*mat->m + i]; 351 } 352 return 0; 353 } 354 355 static int MatiSDscale(Mat matin,Vec l,Vec r) 356 { 357 return 0; 358 } 359 360 static int MatiSDnorm(Mat matin,int type,double *norm) 361 { 362 MatiSD *mat = (MatiSD *) matin->data; 363 Scalar *v = mat->v; 364 double sum = 0.0; 365 int i, j; 366 if (type == NORM_FROBENIUS) { 367 for (i=0; i<mat->n*mat->m; i++ ) { 368 #if defined(PETSC_COMPLEX) 369 sum += real(conj(*v)*(*v)); v++; 370 #else 371 sum += (*v)*(*v); v++; 372 #endif 373 } 374 *norm = sqrt(sum); 375 } 376 else if (type == NORM_1) { 377 *norm = 0.0; 378 for ( j=0; j<mat->n; j++ ) { 379 sum = 0.0; 380 for ( i=0; i<mat->m; i++ ) { 381 #if defined(PETSC_COMPLEX) 382 sum += abs(*v++); 383 #else 384 sum += fabs(*v++); 385 #endif 386 } 387 if (sum > *norm) *norm = sum; 388 } 389 } 390 else if (type == NORM_INFINITY) { 391 *norm = 0.0; 392 for ( j=0; j<mat->m; j++ ) { 393 v = mat->v + j; 394 sum = 0.0; 395 for ( i=0; i<mat->n; i++ ) { 396 #if defined(PETSC_COMPLEX) 397 sum += abs(*v); v += mat->m; 398 #else 399 sum += fabs(*v); v += mat->m; 400 #endif 401 } 402 if (sum > *norm) *norm = sum; 403 } 404 } 405 else { 406 SETERR(1,"No support for the two norm yet"); 407 } 408 return 0; 409 } 410 411 static int MatiDenseinsopt(Mat aijin,int op) 412 { 413 MatiSD *aij = (MatiSD *) aijin->data; 414 if (op == ROW_ORIENTED) aij->roworiented = 1; 415 else if (op == COLUMN_ORIENTED) aij->roworiented = 0; 416 /* doesn't care about sorted rows or columns */ 417 return 0; 418 } 419 420 /* -------------------------------------------------------------------*/ 421 static struct _MatOps MatOps = {MatiSDinsert,MatiSDinsert, 422 MatiSDgetrow, MatiSDrestorerow, 423 MatiSDmult, MatiSDmultadd, MatiSDmulttrans, MatiSDmulttransadd, 424 MatiSDsolve,MatiSDsolveadd,MatiSDsolvetrans,MatiSDsolvetransadd, 425 MatiSDlufactor,MatiSDchfactor, 426 MatiSDrelax, 427 MatiSDtrans, 428 MatiSDnz,MatiSDmemory,MatiSDequal, 429 MatiSDcopy, 430 MatiSDgetdiag,MatiSDscale,MatiSDnorm, 431 0,0, 432 0, MatiDenseinsopt,0,0,0, 433 MatiSDlufactorsymbolic,MatiSDlufactornumeric, 434 MatiSDchfactorsymbolic,MatiSDchfactornumeric 435 }; 436 437 int MatCreateSequentialDense(int m,int n,Mat *newmat) 438 { 439 int size = sizeof(MatiSD) + m*n*sizeof(Scalar); 440 Mat mat; 441 MatiSD *l; 442 *newmat = 0; 443 CREATEHEADER(mat,_Mat); 444 l = (MatiSD *) MALLOC(size); CHKPTR(l); 445 mat->cookie = MAT_COOKIE; 446 mat->ops = &MatOps; 447 mat->destroy = MatiSDdestroy; 448 mat->view = MatiSDview; 449 mat->data = (void *) l; 450 mat->type = MATDENSESEQ; 451 mat->factor = 0; 452 mat->col = 0; 453 mat->row = 0; 454 l->m = m; 455 l->n = n; 456 l->v = (Scalar *) (l + 1); 457 l->pivots = 0; 458 l->roworiented = 1; 459 MEMSET(l->v,0,m*n*sizeof(Scalar)); 460 *newmat = mat; 461 return 0; 462 } 463 464 int MatiSDCreate(Mat matin,Mat *newmat) 465 { 466 MatiSD *m = (MatiSD *) matin->data; 467 return MatCreateSequentialDense(m->m,m->n,newmat); 468 } 469