1 2 #ifndef lint 3 static char vcid[] = "$Id: aij.c,v 1.152 1996/03/01 03:43:59 bsmith Exp bsmith $"; 4 #endif 5 6 /* 7 Defines the basic matrix operations for the AIJ (compressed row) 8 matrix storage format. 9 */ 10 #include "aij.h" 11 #include "vec/vecimpl.h" 12 #include "inline/spops.h" 13 #include "petsc.h" 14 #include "inline/bitarray.h" 15 16 extern int MatToSymmetricIJ_SeqAIJ(int,int*,int*,int,int**,int**); 17 18 static int MatGetReordering_SeqAIJ(Mat A,MatOrdering type,IS *rperm, IS *cperm) 19 { 20 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 21 int ierr, *ia, *ja,n,*idx,i; 22 /*Viewer V1, V2;*/ 23 24 /* 25 this is tacky: In the future when we have written special factorization 26 and solve routines for the identity permutation we should use a 27 stride index set instead of the general one. 28 */ 29 if (type == ORDER_NATURAL) { 30 n = a->n; 31 idx = (int *) PetscMalloc( n*sizeof(int) ); CHKPTRQ(idx); 32 for ( i=0; i<n; i++ ) idx[i] = i; 33 ierr = ISCreateSeq(MPI_COMM_SELF,n,idx,rperm); CHKERRQ(ierr); 34 ierr = ISCreateSeq(MPI_COMM_SELF,n,idx,cperm); CHKERRQ(ierr); 35 PetscFree(idx); 36 ISSetPermutation(*rperm); 37 ISSetPermutation(*cperm); 38 ISSetIdentity(*rperm); 39 ISSetIdentity(*cperm); 40 return 0; 41 } 42 43 ierr = MatToSymmetricIJ_SeqAIJ(a->n,a->i,a->j,a->indexshift, &ia, &ja);CHKERRQ(ierr); 44 ierr = MatGetReordering_IJ(a->n,ia,ja,type,rperm,cperm); CHKERRQ(ierr); 45 46 PetscFree(ia); PetscFree(ja); 47 return 0; 48 } 49 50 #define CHUNKSIZE 15 51 52 /* This version has row oriented v */ 53 static int MatSetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v,InsertMode is) 54 { 55 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 56 int *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax, N, sorted = a->sorted; 57 int *imax = a->imax, *ai = a->i, *ailen = a->ilen,roworiented = a->roworiented; 58 int *aj = a->j, nonew = a->nonew,shift = a->indexshift; 59 Scalar *ap,value, *aa = a->a; 60 61 for ( k=0; k<m; k++ ) { /* loop over added rows */ 62 row = im[k]; 63 if (row < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative row"); 64 if (row >= a->m) SETERRQ(1,"MatSetValues_SeqAIJ:Row too large"); 65 rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 66 rmax = imax[row]; nrow = ailen[row]; 67 low = 0; 68 for ( l=0; l<n; l++ ) { /* loop over added columns */ 69 if (in[l] < 0) SETERRQ(1,"MatSetValues_SeqAIJ:Negative column"); 70 if (in[l] >= a->n) SETERRQ(1,"MatSetValues_SeqAIJ:Column too large"); 71 col = in[l] - shift; 72 if (roworiented) { 73 value = *v++; 74 } 75 else { 76 value = v[k + l*m]; 77 } 78 if (!sorted) low = 0; high = nrow; 79 while (high-low > 5) { 80 t = (low+high)/2; 81 if (rp[t] > col) high = t; 82 else low = t; 83 } 84 for ( i=low; i<high; i++ ) { 85 if (rp[i] > col) break; 86 if (rp[i] == col) { 87 if (is == ADD_VALUES) ap[i] += value; 88 else ap[i] = value; 89 goto noinsert; 90 } 91 } 92 if (nonew) goto noinsert; 93 if (nrow >= rmax) { 94 /* there is no extra room in row, therefore enlarge */ 95 int new_nz = ai[a->m] + CHUNKSIZE,len,*new_i,*new_j; 96 Scalar *new_a; 97 98 /* malloc new storage space */ 99 len = new_nz*(sizeof(int)+sizeof(Scalar))+(a->m+1)*sizeof(int); 100 new_a = (Scalar *) PetscMalloc( len ); CHKPTRQ(new_a); 101 new_j = (int *) (new_a + new_nz); 102 new_i = new_j + new_nz; 103 104 /* copy over old data into new slots */ 105 for ( ii=0; ii<row+1; ii++ ) {new_i[ii] = ai[ii];} 106 for ( ii=row+1; ii<a->m+1; ii++ ) {new_i[ii] = ai[ii]+CHUNKSIZE;} 107 PetscMemcpy(new_j,aj,(ai[row]+nrow+shift)*sizeof(int)); 108 len = (new_nz - CHUNKSIZE - ai[row] - nrow - shift); 109 PetscMemcpy(new_j+ai[row]+shift+nrow+CHUNKSIZE,aj+ai[row]+shift+nrow, 110 len*sizeof(int)); 111 PetscMemcpy(new_a,aa,(ai[row]+nrow+shift)*sizeof(Scalar)); 112 PetscMemcpy(new_a+ai[row]+shift+nrow+CHUNKSIZE,aa+ai[row]+shift+nrow, 113 len*sizeof(Scalar)); 114 /* free up old matrix storage */ 115 PetscFree(a->a); 116 if (!a->singlemalloc) {PetscFree(a->i);PetscFree(a->j);} 117 aa = a->a = new_a; ai = a->i = new_i; aj = a->j = new_j; 118 a->singlemalloc = 1; 119 120 rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 121 rmax = imax[row] = imax[row] + CHUNKSIZE; 122 PLogObjectMemory(A,CHUNKSIZE*(sizeof(int) + sizeof(Scalar))); 123 a->maxnz += CHUNKSIZE; 124 a->reallocs++; 125 } 126 N = nrow++ - 1; a->nz++; 127 /* shift up all the later entries in this row */ 128 for ( ii=N; ii>=i; ii-- ) { 129 rp[ii+1] = rp[ii]; 130 ap[ii+1] = ap[ii]; 131 } 132 rp[i] = col; 133 ap[i] = value; 134 noinsert:; 135 low = i + 1; 136 } 137 ailen[row] = nrow; 138 } 139 return 0; 140 } 141 142 static int MatGetValues_SeqAIJ(Mat A,int m,int *im,int n,int *in,Scalar *v) 143 { 144 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 145 int *rp, k, low, high, t, row, nrow, i, col, l, *aj = a->j; 146 int *ai = a->i, *ailen = a->ilen, shift = a->indexshift; 147 Scalar *ap, *aa = a->a, zero = 0.0; 148 149 for ( k=0; k<m; k++ ) { /* loop over rows */ 150 row = im[k]; 151 if (row < 0) SETERRQ(1,"MatGetValues_SeqAIJ:Negative row"); 152 if (row >= a->m) SETERRQ(1,"MatGetValues_SeqAIJ:Row too large"); 153 rp = aj + ai[row] + shift; ap = aa + ai[row] + shift; 154 nrow = ailen[row]; 155 for ( l=0; l<n; l++ ) { /* loop over columns */ 156 if (in[l] < 0) SETERRQ(1,"MatGetValues_SeqAIJ:Negative column"); 157 if (in[l] >= a->n) SETERRQ(1,"MatGetValues_SeqAIJ:Column too large"); 158 col = in[l] - shift; 159 high = nrow; low = 0; /* assume unsorted */ 160 while (high-low > 5) { 161 t = (low+high)/2; 162 if (rp[t] > col) high = t; 163 else low = t; 164 } 165 for ( i=low; i<high; i++ ) { 166 if (rp[i] > col) break; 167 if (rp[i] == col) { 168 *v++ = ap[i]; 169 goto finished; 170 } 171 } 172 *v++ = zero; 173 finished:; 174 } 175 } 176 return 0; 177 } 178 179 #include "draw.h" 180 #include "pinclude/pviewer.h" 181 #include "sysio.h" 182 183 static int MatView_SeqAIJ_Binary(Mat A,Viewer viewer) 184 { 185 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 186 int i, fd, *col_lens, ierr; 187 188 ierr = ViewerFileGetDescriptor_Private(viewer,&fd); CHKERRQ(ierr); 189 col_lens = (int *) PetscMalloc( (4+a->m)*sizeof(int) ); CHKPTRQ(col_lens); 190 col_lens[0] = MAT_COOKIE; 191 col_lens[1] = a->m; 192 col_lens[2] = a->n; 193 col_lens[3] = a->nz; 194 195 /* store lengths of each row and write (including header) to file */ 196 for ( i=0; i<a->m; i++ ) { 197 col_lens[4+i] = a->i[i+1] - a->i[i]; 198 } 199 ierr = SYWrite(fd,col_lens,4+a->m,SYINT,1); CHKERRQ(ierr); 200 PetscFree(col_lens); 201 202 /* store column indices (zero start index) */ 203 if (a->indexshift) { 204 for ( i=0; i<a->nz; i++ ) a->j[i]--; 205 } 206 ierr = SYWrite(fd,a->j,a->nz,SYINT,0); CHKERRQ(ierr); 207 if (a->indexshift) { 208 for ( i=0; i<a->nz; i++ ) a->j[i]++; 209 } 210 211 /* store nonzero values */ 212 ierr = SYWrite(fd,a->a,a->nz,SYSCALAR,0); CHKERRQ(ierr); 213 return 0; 214 } 215 216 static int MatView_SeqAIJ_ASCII(Mat A,Viewer viewer) 217 { 218 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 219 int ierr, i,j, m = a->m, shift = a->indexshift, format, flg; 220 FILE *fd; 221 char *outputname; 222 223 ierr = ViewerFileGetPointer(viewer,&fd); CHKERRQ(ierr); 224 ierr = ViewerFileGetOutputname_Private(viewer,&outputname); CHKERRQ(ierr); 225 ierr = ViewerFileGetFormat_Private(viewer,&format); 226 if (format == FILE_FORMAT_INFO) { 227 return 0; 228 } 229 else if (format == FILE_FORMAT_INFO_DETAILED) { 230 ierr = OptionsHasName(PETSC_NULL,"-mat_aij_no_inode",&flg); CHKERRQ(ierr); 231 if (flg) fprintf(fd," not using I-node routines\n"); 232 else fprintf(fd," using I-node routines: found %d nodes, limit used is %d\n", 233 a->inode.node_count,a->inode.limit); 234 } 235 else if (format == FILE_FORMAT_MATLAB) { 236 int nz, nzalloc, mem; 237 MatGetInfo(A,MAT_LOCAL,&nz,&nzalloc,&mem); 238 fprintf(fd,"%% Size = %d %d \n",m,a->n); 239 fprintf(fd,"%% Nonzeros = %d \n",nz); 240 fprintf(fd,"zzz = zeros(%d,3);\n",nz); 241 fprintf(fd,"zzz = [\n"); 242 243 for (i=0; i<m; i++) { 244 for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 245 #if defined(PETSC_COMPLEX) 246 fprintf(fd,"%d %d %18.16e %18.16e \n",i+1,a->j[j]+!shift,real(a->a[j]), 247 imag(a->a[j])); 248 #else 249 fprintf(fd,"%d %d %18.16e\n", i+1, a->j[j]+!shift, a->a[j]); 250 #endif 251 } 252 } 253 fprintf(fd,"];\n %s = spconvert(zzz);\n",outputname); 254 } 255 else { 256 for ( i=0; i<m; i++ ) { 257 fprintf(fd,"row %d:",i); 258 for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 259 #if defined(PETSC_COMPLEX) 260 if (imag(a->a[j]) != 0.0) { 261 fprintf(fd," %d %g + %g i",a->j[j]+shift,real(a->a[j]),imag(a->a[j])); 262 } 263 else { 264 fprintf(fd," %d %g ",a->j[j]+shift,real(a->a[j])); 265 } 266 #else 267 fprintf(fd," %d %g ",a->j[j]+shift,a->a[j]); 268 #endif 269 } 270 fprintf(fd,"\n"); 271 } 272 } 273 fflush(fd); 274 return 0; 275 } 276 277 static int MatView_SeqAIJ_Draw(Mat A,Viewer viewer) 278 { 279 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 280 int ierr, i,j, m = a->m, shift = a->indexshift,pause,color; 281 double xl,yl,xr,yr,w,h,xc,yc,scale = 1.0,x_l,x_r,y_l,y_r; 282 Draw draw = (Draw) viewer; 283 DrawButton button; 284 285 xr = a->n; yr = a->m; h = yr/10.0; w = xr/10.0; 286 xr += w; yr += h; xl = -w; yl = -h; 287 ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 288 /* loop over matrix elements drawing boxes */ 289 color = DRAW_BLUE; 290 for ( i=0; i<m; i++ ) { 291 y_l = m - i - 1.0; y_r = y_l + 1.0; 292 for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 293 x_l = a->j[j] + shift; x_r = x_l + 1.0; 294 #if defined(PETSC_COMPLEX) 295 if (real(a->a[j]) >= 0.) continue; 296 #else 297 if (a->a[j] >= 0.) continue; 298 #endif 299 DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 300 } 301 } 302 color = DRAW_CYAN; 303 for ( i=0; i<m; i++ ) { 304 y_l = m - i - 1.0; y_r = y_l + 1.0; 305 for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 306 x_l = a->j[j] + shift; x_r = x_l + 1.0; 307 if (a->a[j] != 0.) continue; 308 DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 309 } 310 } 311 color = DRAW_RED; 312 for ( i=0; i<m; i++ ) { 313 y_l = m - i - 1.0; y_r = y_l + 1.0; 314 for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 315 x_l = a->j[j] + shift; x_r = x_l + 1.0; 316 #if defined(PETSC_COMPLEX) 317 if (real(a->a[j]) <= 0.) continue; 318 #else 319 if (a->a[j] <= 0.) continue; 320 #endif 321 DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 322 } 323 } 324 DrawFlush(draw); 325 DrawGetPause(draw,&pause); 326 if (pause >= 0) { PetscSleep(pause); return 0;} 327 328 /* allow the matrix to zoom or shrink */ 329 ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0); 330 while (button != BUTTON_RIGHT) { 331 DrawClear(draw); 332 if (button == BUTTON_LEFT) scale = .5; 333 else if (button == BUTTON_CENTER) scale = 2.; 334 xl = scale*(xl + w - xc) + xc - w*scale; 335 xr = scale*(xr - w - xc) + xc + w*scale; 336 yl = scale*(yl + h - yc) + yc - h*scale; 337 yr = scale*(yr - h - yc) + yc + h*scale; 338 w *= scale; h *= scale; 339 ierr = DrawSetCoordinates(draw,xl,yl,xr,yr); CHKERRQ(ierr); 340 color = DRAW_BLUE; 341 for ( i=0; i<m; i++ ) { 342 y_l = m - i - 1.0; y_r = y_l + 1.0; 343 for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 344 x_l = a->j[j] + shift; x_r = x_l + 1.0; 345 #if defined(PETSC_COMPLEX) 346 if (real(a->a[j]) >= 0.) continue; 347 #else 348 if (a->a[j] >= 0.) continue; 349 #endif 350 DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 351 } 352 } 353 color = DRAW_CYAN; 354 for ( i=0; i<m; i++ ) { 355 y_l = m - i - 1.0; y_r = y_l + 1.0; 356 for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 357 x_l = a->j[j] + shift; x_r = x_l + 1.0; 358 if (a->a[j] != 0.) continue; 359 DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 360 } 361 } 362 color = DRAW_RED; 363 for ( i=0; i<m; i++ ) { 364 y_l = m - i - 1.0; y_r = y_l + 1.0; 365 for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 366 x_l = a->j[j] + shift; x_r = x_l + 1.0; 367 #if defined(PETSC_COMPLEX) 368 if (real(a->a[j]) <= 0.) continue; 369 #else 370 if (a->a[j] <= 0.) continue; 371 #endif 372 DrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color); 373 } 374 } 375 ierr = DrawGetMouseButton(draw,&button,&xc,&yc,0,0); 376 } 377 return 0; 378 } 379 380 static int MatView_SeqAIJ(PetscObject obj,Viewer viewer) 381 { 382 Mat A = (Mat) obj; 383 Mat_SeqAIJ *a = (Mat_SeqAIJ*) A->data; 384 PetscObject vobj = (PetscObject) viewer; 385 386 if (!viewer) { 387 viewer = STDOUT_VIEWER_SELF; vobj = (PetscObject) viewer; 388 } 389 if (vobj->cookie == VIEWER_COOKIE) { 390 if (vobj->type == MATLAB_VIEWER) { 391 return ViewerMatlabPutSparse_Private(viewer,a->m,a->n,a->nz,a->a,a->i,a->j); 392 } 393 else if (vobj->type == ASCII_FILE_VIEWER || vobj->type == ASCII_FILES_VIEWER){ 394 return MatView_SeqAIJ_ASCII(A,viewer); 395 } 396 else if (vobj->type == BINARY_FILE_VIEWER) { 397 return MatView_SeqAIJ_Binary(A,viewer); 398 } 399 } 400 else if (vobj->cookie == DRAW_COOKIE) { 401 if (vobj->type == NULLWINDOW) return 0; 402 else return MatView_SeqAIJ_Draw(A,viewer); 403 } 404 return 0; 405 } 406 extern int Mat_AIJ_CheckInode(Mat); 407 static int MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode) 408 { 409 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 410 int fshift = 0,i,j,*ai = a->i, *aj = a->j, *imax = a->imax,ierr; 411 int m = a->m, *ip, N, *ailen = a->ilen,shift = a->indexshift; 412 Scalar *aa = a->a, *ap; 413 414 if (mode == FLUSH_ASSEMBLY) return 0; 415 416 for ( i=1; i<m; i++ ) { 417 /* move each row back by the amount of empty slots (fshift) before it*/ 418 fshift += imax[i-1] - ailen[i-1]; 419 if (fshift) { 420 ip = aj + ai[i] + shift; ap = aa + ai[i] + shift; 421 N = ailen[i]; 422 for ( j=0; j<N; j++ ) { 423 ip[j-fshift] = ip[j]; 424 ap[j-fshift] = ap[j]; 425 } 426 } 427 ai[i] = ai[i-1] + ailen[i-1]; 428 } 429 if (m) { 430 fshift += imax[m-1] - ailen[m-1]; 431 ai[m] = ai[m-1] + ailen[m-1]; 432 } 433 /* reset ilen and imax for each row */ 434 for ( i=0; i<m; i++ ) { 435 ailen[i] = imax[i] = ai[i+1] - ai[i]; 436 } 437 a->nz = ai[m] + shift; 438 439 /* diagonals may have moved, so kill the diagonal pointers */ 440 if (fshift && a->diag) { 441 PetscFree(a->diag); 442 PLogObjectMemory(A,-(m+1)*sizeof(int)); 443 a->diag = 0; 444 } 445 PLogInfo((PetscObject)A,"MatAssemblyEnd_SeqAIJ:Unneeded storage space %d used %d rows %d\n", 446 fshift,a->nz,m); 447 PLogInfo((PetscObject)A,"MatAssemblyEnd_SeqAIJ:Number of mallocs during MatSetValues %d\n", 448 a->reallocs); 449 /* check out for identical nodes. If found, use inode functions */ 450 ierr = Mat_AIJ_CheckInode(A); CHKERRQ(ierr); 451 return 0; 452 } 453 454 static int MatZeroEntries_SeqAIJ(Mat A) 455 { 456 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 457 PetscMemzero(a->a,(a->i[a->m]+a->indexshift)*sizeof(Scalar)); 458 return 0; 459 } 460 461 int MatDestroy_SeqAIJ(PetscObject obj) 462 { 463 Mat A = (Mat) obj; 464 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 465 466 #if defined(PETSC_LOG) 467 PLogObjectState(obj,"Rows=%d, Cols=%d, NZ=%d",a->m,a->n,a->nz); 468 #endif 469 PetscFree(a->a); 470 if (!a->singlemalloc) { PetscFree(a->i); PetscFree(a->j);} 471 if (a->diag) PetscFree(a->diag); 472 if (a->ilen) PetscFree(a->ilen); 473 if (a->imax) PetscFree(a->imax); 474 if (a->solve_work) PetscFree(a->solve_work); 475 if (a->inode.size) PetscFree(a->inode.size); 476 PetscFree(a); 477 PLogObjectDestroy(A); 478 PetscHeaderDestroy(A); 479 return 0; 480 } 481 482 static int MatCompress_SeqAIJ(Mat A) 483 { 484 return 0; 485 } 486 487 static int MatSetOption_SeqAIJ(Mat A,MatOption op) 488 { 489 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 490 if (op == ROW_ORIENTED) a->roworiented = 1; 491 else if (op == COLUMN_ORIENTED) a->roworiented = 0; 492 else if (op == COLUMNS_SORTED) a->sorted = 1; 493 else if (op == NO_NEW_NONZERO_LOCATIONS) a->nonew = 1; 494 else if (op == YES_NEW_NONZERO_LOCATIONS) a->nonew = 0; 495 else if (op == ROWS_SORTED || 496 op == SYMMETRIC_MATRIX || 497 op == STRUCTURALLY_SYMMETRIC_MATRIX || 498 op == YES_NEW_DIAGONALS) 499 PLogInfo((PetscObject)A,"Info:MatSetOption_SeqAIJ:Option ignored\n"); 500 else if (op == NO_NEW_DIAGONALS) 501 {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:NO_NEW_DIAGONALS");} 502 else if (op == INODE_LIMIT_1) a->inode.limit = 1; 503 else if (op == INODE_LIMIT_2) a->inode.limit = 2; 504 else if (op == INODE_LIMIT_3) a->inode.limit = 3; 505 else if (op == INODE_LIMIT_4) a->inode.limit = 4; 506 else if (op == INODE_LIMIT_5) a->inode.limit = 5; 507 else 508 {SETERRQ(PETSC_ERR_SUP,"MatSetOption_SeqAIJ:unknown option");} 509 return 0; 510 } 511 512 static int MatGetDiagonal_SeqAIJ(Mat A,Vec v) 513 { 514 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 515 int i,j, n,shift = a->indexshift; 516 Scalar *x, zero = 0.0; 517 518 VecSet(&zero,v); 519 VecGetArray(v,&x); VecGetLocalSize(v,&n); 520 if (n != a->m) SETERRQ(1,"MatGetDiagonal_SeqAIJ:Nonconforming matrix and vector"); 521 for ( i=0; i<a->m; i++ ) { 522 for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 523 if (a->j[j]+shift == i) { 524 x[i] = a->a[j]; 525 break; 526 } 527 } 528 } 529 return 0; 530 } 531 532 /* -------------------------------------------------------*/ 533 /* Should check that shapes of vectors and matrices match */ 534 /* -------------------------------------------------------*/ 535 static int MatMultTrans_SeqAIJ(Mat A,Vec xx,Vec yy) 536 { 537 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 538 Scalar *x, *y, *v, alpha; 539 int m = a->m, n, i, *idx, shift = a->indexshift; 540 541 VecGetArray(xx,&x); VecGetArray(yy,&y); 542 PetscMemzero(y,a->n*sizeof(Scalar)); 543 y = y + shift; /* shift for Fortran start by 1 indexing */ 544 for ( i=0; i<m; i++ ) { 545 idx = a->j + a->i[i] + shift; 546 v = a->a + a->i[i] + shift; 547 n = a->i[i+1] - a->i[i]; 548 alpha = x[i]; 549 while (n-->0) {y[*idx++] += alpha * *v++;} 550 } 551 PLogFlops(2*a->nz - a->n); 552 return 0; 553 } 554 555 static int MatMultTransAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy) 556 { 557 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 558 Scalar *x, *y, *v, alpha; 559 int m = a->m, n, i, *idx,shift = a->indexshift; 560 561 VecGetArray(xx,&x); VecGetArray(yy,&y); 562 if (zz != yy) VecCopy(zz,yy); 563 y = y + shift; /* shift for Fortran start by 1 indexing */ 564 for ( i=0; i<m; i++ ) { 565 idx = a->j + a->i[i] + shift; 566 v = a->a + a->i[i] + shift; 567 n = a->i[i+1] - a->i[i]; 568 alpha = x[i]; 569 while (n-->0) {y[*idx++] += alpha * *v++;} 570 } 571 return 0; 572 } 573 574 static int MatMult_SeqAIJ(Mat A,Vec xx,Vec yy) 575 { 576 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 577 Scalar *x, *y, *v, sum; 578 int m = a->m, n, i, *idx, shift = a->indexshift,*ii; 579 580 VecGetArray(xx,&x); VecGetArray(yy,&y); 581 x = x + shift; /* shift for Fortran start by 1 indexing */ 582 idx = a->j; 583 v = a->a; 584 ii = a->i; 585 for ( i=0; i<m; i++ ) { 586 n = ii[1] - ii[0]; ii++; 587 sum = 0.0; 588 /* SPARSEDENSEDOT(sum,x,v,idx,n); */ 589 /* for ( j=n-1; j>-1; j--) sum += v[j]*x[idx[j]]; */ 590 while (n--) sum += *v++ * x[*idx++]; 591 y[i] = sum; 592 } 593 PLogFlops(2*a->nz - m); 594 return 0; 595 } 596 597 static int MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz) 598 { 599 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 600 Scalar *x, *y, *z, *v, sum; 601 int m = a->m, n, i, *idx, shift = a->indexshift,*ii; 602 603 VecGetArray(xx,&x); VecGetArray(yy,&y); VecGetArray(zz,&z); 604 x = x + shift; /* shift for Fortran start by 1 indexing */ 605 idx = a->j; 606 v = a->a; 607 ii = a->i; 608 for ( i=0; i<m; i++ ) { 609 n = ii[1] - ii[0]; ii++; 610 sum = y[i]; 611 /* SPARSEDENSEDOT(sum,x,v,idx,n); */ 612 while (n--) sum += *v++ * x[*idx++]; 613 z[i] = sum; 614 } 615 PLogFlops(2*a->nz); 616 return 0; 617 } 618 619 /* 620 Adds diagonal pointers to sparse matrix structure. 621 */ 622 623 int MatMarkDiag_SeqAIJ(Mat A) 624 { 625 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 626 int i,j, *diag, m = a->m,shift = a->indexshift; 627 628 diag = (int *) PetscMalloc( (m+1)*sizeof(int)); CHKPTRQ(diag); 629 PLogObjectMemory(A,(m+1)*sizeof(int)); 630 for ( i=0; i<a->m; i++ ) { 631 for ( j=a->i[i]+shift; j<a->i[i+1]+shift; j++ ) { 632 if (a->j[j]+shift == i) { 633 diag[i] = j - shift; 634 break; 635 } 636 } 637 } 638 a->diag = diag; 639 return 0; 640 } 641 642 static int MatRelax_SeqAIJ(Mat A,Vec bb,double omega,MatSORType flag, 643 double fshift,int its,Vec xx) 644 { 645 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 646 Scalar *x, *b, *bs, d, *xs, sum, *v = a->a,*t,scale,*ts, *xb; 647 int ierr, *idx, *diag,n = a->n, m = a->m, i, shift = a->indexshift; 648 649 VecGetArray(xx,&x); VecGetArray(bb,&b); 650 if (!a->diag) {if ((ierr = MatMarkDiag_SeqAIJ(A))) return ierr;} 651 diag = a->diag; 652 xs = x + shift; /* shifted by one for index start of a or a->j*/ 653 if (flag == SOR_APPLY_UPPER) { 654 /* apply ( U + D/omega) to the vector */ 655 bs = b + shift; 656 for ( i=0; i<m; i++ ) { 657 d = fshift + a->a[diag[i] + shift]; 658 n = a->i[i+1] - diag[i] - 1; 659 idx = a->j + diag[i] + (!shift); 660 v = a->a + diag[i] + (!shift); 661 sum = b[i]*d/omega; 662 SPARSEDENSEDOT(sum,bs,v,idx,n); 663 x[i] = sum; 664 } 665 return 0; 666 } 667 if (flag == SOR_APPLY_LOWER) { 668 SETERRQ(1,"MatRelax_SeqAIJ:SOR_APPLY_LOWER is not done"); 669 } 670 else if (flag & SOR_EISENSTAT) { 671 /* Let A = L + U + D; where L is lower trianglar, 672 U is upper triangular, E is diagonal; This routine applies 673 674 (L + E)^{-1} A (U + E)^{-1} 675 676 to a vector efficiently using Eisenstat's trick. This is for 677 the case of SSOR preconditioner, so E is D/omega where omega 678 is the relaxation factor. 679 */ 680 t = (Scalar *) PetscMalloc( m*sizeof(Scalar) ); CHKPTRQ(t); 681 scale = (2.0/omega) - 1.0; 682 683 /* x = (E + U)^{-1} b */ 684 for ( i=m-1; i>=0; i-- ) { 685 d = fshift + a->a[diag[i] + shift]; 686 n = a->i[i+1] - diag[i] - 1; 687 idx = a->j + diag[i] + (!shift); 688 v = a->a + diag[i] + (!shift); 689 sum = b[i]; 690 SPARSEDENSEMDOT(sum,xs,v,idx,n); 691 x[i] = omega*(sum/d); 692 } 693 694 /* t = b - (2*E - D)x */ 695 v = a->a; 696 for ( i=0; i<m; i++ ) { t[i] = b[i] - scale*(v[*diag++ + shift])*x[i]; } 697 698 /* t = (E + L)^{-1}t */ 699 ts = t + shift; /* shifted by one for index start of a or a->j*/ 700 diag = a->diag; 701 for ( i=0; i<m; i++ ) { 702 d = fshift + a->a[diag[i]+shift]; 703 n = diag[i] - a->i[i]; 704 idx = a->j + a->i[i] + shift; 705 v = a->a + a->i[i] + shift; 706 sum = t[i]; 707 SPARSEDENSEMDOT(sum,ts,v,idx,n); 708 t[i] = omega*(sum/d); 709 } 710 711 /* x = x + t */ 712 for ( i=0; i<m; i++ ) { x[i] += t[i]; } 713 PetscFree(t); 714 return 0; 715 } 716 if (flag & SOR_ZERO_INITIAL_GUESS) { 717 if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 718 for ( i=0; i<m; i++ ) { 719 d = fshift + a->a[diag[i]+shift]; 720 n = diag[i] - a->i[i]; 721 idx = a->j + a->i[i] + shift; 722 v = a->a + a->i[i] + shift; 723 sum = b[i]; 724 SPARSEDENSEMDOT(sum,xs,v,idx,n); 725 x[i] = omega*(sum/d); 726 } 727 xb = x; 728 } 729 else xb = b; 730 if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) && 731 (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) { 732 for ( i=0; i<m; i++ ) { 733 x[i] *= a->a[diag[i]+shift]; 734 } 735 } 736 if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 737 for ( i=m-1; i>=0; i-- ) { 738 d = fshift + a->a[diag[i] + shift]; 739 n = a->i[i+1] - diag[i] - 1; 740 idx = a->j + diag[i] + (!shift); 741 v = a->a + diag[i] + (!shift); 742 sum = xb[i]; 743 SPARSEDENSEMDOT(sum,xs,v,idx,n); 744 x[i] = omega*(sum/d); 745 } 746 } 747 its--; 748 } 749 while (its--) { 750 if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 751 for ( i=0; i<m; i++ ) { 752 d = fshift + a->a[diag[i]+shift]; 753 n = a->i[i+1] - a->i[i]; 754 idx = a->j + a->i[i] + shift; 755 v = a->a + a->i[i] + shift; 756 sum = b[i]; 757 SPARSEDENSEMDOT(sum,xs,v,idx,n); 758 x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 759 } 760 } 761 if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 762 for ( i=m-1; i>=0; i-- ) { 763 d = fshift + a->a[diag[i] + shift]; 764 n = a->i[i+1] - a->i[i]; 765 idx = a->j + a->i[i] + shift; 766 v = a->a + a->i[i] + shift; 767 sum = b[i]; 768 SPARSEDENSEMDOT(sum,xs,v,idx,n); 769 x[i] = (1. - omega)*x[i] + omega*(sum/d + x[i]); 770 } 771 } 772 } 773 return 0; 774 } 775 776 static int MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,int *nz,int *nzalloc,int *mem) 777 { 778 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 779 *nz = a->nz; 780 *nzalloc = a->maxnz; 781 *mem = (int)A->mem; 782 return 0; 783 } 784 785 extern int MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,double,Mat*); 786 extern int MatLUFactorNumeric_SeqAIJ(Mat,Mat*); 787 extern int MatLUFactor_SeqAIJ(Mat,IS,IS,double); 788 extern int MatSolve_SeqAIJ(Mat,Vec,Vec); 789 extern int MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec); 790 extern int MatSolveTrans_SeqAIJ(Mat,Vec,Vec); 791 extern int MatSolveTransAdd_SeqAIJ(Mat,Vec,Vec,Vec); 792 793 static int MatZeroRows_SeqAIJ(Mat A,IS is,Scalar *diag) 794 { 795 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 796 int i,ierr,N, *rows,m = a->m - 1,shift = a->indexshift; 797 798 ierr = ISGetLocalSize(is,&N); CHKERRQ(ierr); 799 ierr = ISGetIndices(is,&rows); CHKERRQ(ierr); 800 if (diag) { 801 for ( i=0; i<N; i++ ) { 802 if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range"); 803 if (a->ilen[rows[i]] > 0) { /* in case row was completely empty */ 804 a->ilen[rows[i]] = 1; 805 a->a[a->i[rows[i]]+shift] = *diag; 806 a->j[a->i[rows[i]]+shift] = rows[i]+shift; 807 } 808 else { 809 ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],diag,INSERT_VALUES); 810 CHKERRQ(ierr); 811 } 812 } 813 } 814 else { 815 for ( i=0; i<N; i++ ) { 816 if (rows[i] < 0 || rows[i] > m) SETERRQ(1,"MatZeroRows_SeqAIJ:row out of range"); 817 a->ilen[rows[i]] = 0; 818 } 819 } 820 ISRestoreIndices(is,&rows); 821 ierr = MatAssemblyBegin(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 822 ierr = MatAssemblyEnd(A,FINAL_ASSEMBLY); CHKERRQ(ierr); 823 return 0; 824 } 825 826 static int MatGetSize_SeqAIJ(Mat A,int *m,int *n) 827 { 828 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 829 *m = a->m; *n = a->n; 830 return 0; 831 } 832 833 static int MatGetOwnershipRange_SeqAIJ(Mat A,int *m,int *n) 834 { 835 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 836 *m = 0; *n = a->m; 837 return 0; 838 } 839 int MatGetRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 840 { 841 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 842 int *itmp,i,shift = a->indexshift; 843 844 if (row < 0 || row >= a->m) SETERRQ(1,"MatGetRow_SeqAIJ:Row out of range"); 845 846 *nz = a->i[row+1] - a->i[row]; 847 if (v) *v = a->a + a->i[row] + shift; 848 if (idx) { 849 itmp = a->j + a->i[row] + shift; 850 if (*nz && shift) { 851 *idx = (int *) PetscMalloc( (*nz)*sizeof(int) ); CHKPTRQ(*idx); 852 for ( i=0; i<(*nz); i++ ) {(*idx)[i] = itmp[i] + shift;} 853 } else if (*nz) { 854 *idx = itmp; 855 } 856 else *idx = 0; 857 } 858 return 0; 859 } 860 861 int MatRestoreRow_SeqAIJ(Mat A,int row,int *nz,int **idx,Scalar **v) 862 { 863 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 864 if (idx) {if (*idx && a->indexshift) PetscFree(*idx);} 865 return 0; 866 } 867 868 static int MatNorm_SeqAIJ(Mat A,NormType type,double *norm) 869 { 870 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 871 Scalar *v = a->a; 872 double sum = 0.0; 873 int i, j,shift = a->indexshift; 874 875 if (type == NORM_FROBENIUS) { 876 for (i=0; i<a->nz; i++ ) { 877 #if defined(PETSC_COMPLEX) 878 sum += real(conj(*v)*(*v)); v++; 879 #else 880 sum += (*v)*(*v); v++; 881 #endif 882 } 883 *norm = sqrt(sum); 884 } 885 else if (type == NORM_1) { 886 double *tmp; 887 int *jj = a->j; 888 tmp = (double *) PetscMalloc( a->n*sizeof(double) ); CHKPTRQ(tmp); 889 PetscMemzero(tmp,a->n*sizeof(double)); 890 *norm = 0.0; 891 for ( j=0; j<a->nz; j++ ) { 892 tmp[*jj++ + shift] += PetscAbsScalar(*v); v++; 893 } 894 for ( j=0; j<a->n; j++ ) { 895 if (tmp[j] > *norm) *norm = tmp[j]; 896 } 897 PetscFree(tmp); 898 } 899 else if (type == NORM_INFINITY) { 900 *norm = 0.0; 901 for ( j=0; j<a->m; j++ ) { 902 v = a->a + a->i[j] + shift; 903 sum = 0.0; 904 for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) { 905 sum += PetscAbsScalar(*v); v++; 906 } 907 if (sum > *norm) *norm = sum; 908 } 909 } 910 else { 911 SETERRQ(1,"MatNorm_SeqAIJ:No support for two norm yet"); 912 } 913 return 0; 914 } 915 916 static int MatTranspose_SeqAIJ(Mat A,Mat *B) 917 { 918 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 919 Mat C; 920 int i, ierr, *aj = a->j, *ai = a->i, m = a->m, len, *col; 921 int shift = a->indexshift; 922 Scalar *array = a->a; 923 924 if (B == PETSC_NULL && m != a->n) 925 SETERRQ(1,"MatTranspose_MPIAIJ:Square matrix only for in-place"); 926 col = (int *) PetscMalloc((1+a->n)*sizeof(int)); CHKPTRQ(col); 927 PetscMemzero(col,(1+a->n)*sizeof(int)); 928 if (shift) { 929 for ( i=0; i<ai[m]-1; i++ ) aj[i] -= 1; 930 } 931 for ( i=0; i<ai[m]+shift; i++ ) col[aj[i]] += 1; 932 ierr = MatCreateSeqAIJ(A->comm,a->n,m,0,col,&C); CHKERRQ(ierr); 933 PetscFree(col); 934 for ( i=0; i<m; i++ ) { 935 len = ai[i+1]-ai[i]; 936 ierr = MatSetValues(C,len,aj,1,&i,array,INSERT_VALUES); CHKERRQ(ierr); 937 array += len; aj += len; 938 } 939 if (shift) { 940 for ( i=0; i<ai[m]-1; i++ ) aj[i] += 1; 941 } 942 943 ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 944 ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 945 946 if (B != PETSC_NULL) { 947 *B = C; 948 } else { 949 /* This isn't really an in-place transpose */ 950 PetscFree(a->a); 951 if (!a->singlemalloc) {PetscFree(a->i); PetscFree(a->j);} 952 if (a->diag) PetscFree(a->diag); 953 if (a->ilen) PetscFree(a->ilen); 954 if (a->imax) PetscFree(a->imax); 955 if (a->solve_work) PetscFree(a->solve_work); 956 PetscFree(a); 957 PetscMemcpy(A,C,sizeof(struct _Mat)); 958 PetscHeaderDestroy(C); 959 } 960 return 0; 961 } 962 963 static int MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr) 964 { 965 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 966 Scalar *l,*r,x,*v; 967 int i,j,m = a->m, n = a->n, M, nz = a->nz, *jj,shift = a->indexshift; 968 969 if (ll) { 970 VecGetArray(ll,&l); VecGetSize(ll,&m); 971 if (m != a->m) SETERRQ(1,"MatDiagonalScale_SeqAIJ:Left scaling vector wrong length"); 972 v = a->a; 973 for ( i=0; i<m; i++ ) { 974 x = l[i]; 975 M = a->i[i+1] - a->i[i]; 976 for ( j=0; j<M; j++ ) { (*v++) *= x;} 977 } 978 } 979 if (rr) { 980 VecGetArray(rr,&r); VecGetSize(rr,&n); 981 if (n != a->n) SETERRQ(1,"MatDiagonalScale_SeqAIJ:Right scaling vector wrong length"); 982 v = a->a; jj = a->j; 983 for ( i=0; i<nz; i++ ) { 984 (*v++) *= r[*jj++ + shift]; 985 } 986 } 987 return 0; 988 } 989 990 static int MatGetSubMatrix_SeqAIJ(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,Mat *B) 991 { 992 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data,*c; 993 int nznew, *smap, i, k, kstart, kend, ierr, oldcols = a->n,*lens; 994 register int sum,lensi; 995 int *irow, *icol, nrows, ncols, *cwork, shift = a->indexshift,*ssmap; 996 int first,step,*starts,*j_new,*i_new,*aj = a->j, *ai = a->i,ii,*ailen = a->ilen; 997 Scalar *vwork,*a_new; 998 Mat C; 999 1000 ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr); 1001 ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr); 1002 ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr); 1003 1004 if (ISStrideGetInfo(iscol,&first,&step) && step == 1) { /* no need to sort */ 1005 /* special case of contiguous rows */ 1006 lens = (int *) PetscMalloc((2*ncols+1)*sizeof(int)); CHKPTRQ(lens); 1007 starts = lens + ncols; 1008 /* loop over new rows determining lens and starting points */ 1009 for (i=0; i<nrows; i++) { 1010 kstart = ai[irow[i]]+shift; 1011 kend = kstart + ailen[irow[i]]; 1012 for ( k=kstart; k<kend; k++ ) { 1013 if (aj[k]+shift >= first) { 1014 starts[i] = k; 1015 break; 1016 } 1017 } 1018 sum = 0; 1019 while (k < kend) { 1020 if (aj[k++]+shift >= first+ncols) break; 1021 sum++; 1022 } 1023 lens[i] = sum; 1024 } 1025 /* create submatrix */ 1026 if (scall == MAT_REUSE_MATRIX) { 1027 int n_cols,n_rows; 1028 ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr); 1029 if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:"); 1030 ierr = MatZeroEntries(*B); CHKERRQ(ierr); 1031 C = *B; 1032 } 1033 else { 1034 ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr); 1035 } 1036 c = (Mat_SeqAIJ*) C->data; 1037 1038 /* loop over rows inserting into submatrix */ 1039 a_new = c->a; 1040 j_new = c->j; 1041 i_new = c->i; 1042 i_new[0] = -shift; 1043 for (i=0; i<nrows; i++) { 1044 ii = starts[i]; 1045 lensi = lens[i]; 1046 for ( k=0; k<lensi; k++ ) { 1047 *j_new++ = aj[ii+k] - first; 1048 } 1049 PetscMemcpy(a_new,a->a + starts[i],lensi*sizeof(Scalar)); 1050 a_new += lensi; 1051 i_new[i+1] = i_new[i] + lensi; 1052 c->ilen[i] = lensi; 1053 } 1054 PetscFree(lens); 1055 } 1056 else { 1057 ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr); 1058 ierr = SYIsort(ncols,icol); CHKERRQ(ierr); 1059 smap = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap); 1060 ssmap = smap + shift; 1061 cwork = (int *) PetscMalloc((1+nrows+ncols)*sizeof(int)); CHKPTRQ(cwork); 1062 lens = cwork + ncols; 1063 vwork = (Scalar *) PetscMalloc((1+ncols)*sizeof(Scalar)); CHKPTRQ(vwork); 1064 PetscMemzero(smap,oldcols*sizeof(int)); 1065 for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1; 1066 /* determine lens of each row */ 1067 for (i=0; i<nrows; i++) { 1068 kstart = ai[irow[i]]+shift; 1069 kend = kstart + a->ilen[irow[i]]; 1070 lens[i] = 0; 1071 for ( k=kstart; k<kend; k++ ) { 1072 if (ssmap[aj[k]]) { 1073 lens[i]++; 1074 } 1075 } 1076 } 1077 /* Create and fill new matrix */ 1078 if (scall == MAT_REUSE_MATRIX) { 1079 int n_cols,n_rows; 1080 ierr = MatGetSize(*B,&n_rows,&n_cols); CHKERRQ(ierr); 1081 if (n_rows != nrows || n_cols != ncols) SETERRQ(1,"MatGetSubMatrix_SeqAIJ:"); 1082 ierr = MatZeroEntries(*B); CHKERRQ(ierr); 1083 C = *B; 1084 } 1085 else { 1086 ierr = MatCreateSeqAIJ(A->comm,nrows,ncols,0,lens,&C);CHKERRQ(ierr); 1087 } 1088 for (i=0; i<nrows; i++) { 1089 nznew = 0; 1090 kstart = ai[irow[i]]+shift; 1091 kend = kstart + a->ilen[irow[i]]; 1092 for ( k=kstart; k<kend; k++ ) { 1093 if (ssmap[a->j[k]]) { 1094 cwork[nznew] = ssmap[a->j[k]] - 1; 1095 vwork[nznew++] = a->a[k]; 1096 } 1097 } 1098 ierr = MatSetValues(C,1,&i,nznew,cwork,vwork,INSERT_VALUES); CHKERRQ(ierr); 1099 } 1100 /* Free work space */ 1101 ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr); 1102 PetscFree(smap); PetscFree(cwork); PetscFree(vwork); 1103 } 1104 ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 1105 ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 1106 1107 ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr); 1108 *B = C; 1109 return 0; 1110 } 1111 1112 /* 1113 note: This can only work for identity for row and col. It would 1114 be good to check this and otherwise generate an error. 1115 */ 1116 static int MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,double efill,int fill) 1117 { 1118 Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data; 1119 int ierr; 1120 Mat outA; 1121 1122 if (fill != 0) SETERRQ(1,"MatILUFactor_SeqAIJ:Only fill=0 supported"); 1123 1124 outA = inA; 1125 inA->factor = FACTOR_LU; 1126 a->row = row; 1127 a->col = col; 1128 1129 a->solve_work = (Scalar *) PetscMalloc( (a->m+1)*sizeof(Scalar)); CHKPTRQ(a->solve_work); 1130 1131 if (!a->diag) { 1132 ierr = MatMarkDiag_SeqAIJ(inA); CHKERRQ(ierr); 1133 } 1134 ierr = MatLUFactorNumeric_SeqAIJ(inA,&outA); CHKERRQ(ierr); 1135 return 0; 1136 } 1137 1138 #include "pinclude/plapack.h" 1139 static int MatScale_SeqAIJ(Scalar *alpha,Mat inA) 1140 { 1141 Mat_SeqAIJ *a = (Mat_SeqAIJ *) inA->data; 1142 int one = 1; 1143 BLscal_( &a->nz, alpha, a->a, &one ); 1144 PLogFlops(a->nz); 1145 return 0; 1146 } 1147 1148 static int MatGetSubMatrices_SeqAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall, 1149 Mat **B) 1150 { 1151 int ierr,i; 1152 1153 if (scall == MAT_INITIAL_MATRIX) { 1154 *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B); 1155 } 1156 1157 for ( i=0; i<n; i++ ) { 1158 ierr = MatGetSubMatrix(A,irow[i],icol[i],scall,&(*B)[i]); CHKERRQ(ierr); 1159 } 1160 return 0; 1161 } 1162 1163 static int MatIncreaseOverlap_SeqAIJ(Mat A, int is_max, IS *is, int ov) 1164 { 1165 Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data; 1166 int shift, row, i,j,k,l,m,n, *idx,ierr, *nidx, isz, val; 1167 int start, end, *ai, *aj; 1168 char *table; 1169 shift = a->indexshift; 1170 m = a->m; 1171 ai = a->i; 1172 aj = a->j+shift; 1173 1174 if (ov < 0) SETERRQ(1,"MatIncreaseOverlap_SeqAIJ: illegal overlap value used"); 1175 1176 table = (char *) PetscMalloc((m/BITSPERBYTE +1)*sizeof(char)); CHKPTRQ(table); 1177 nidx = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(nidx); 1178 1179 for ( i=0; i<is_max; i++ ) { 1180 /* Initialise the two local arrays */ 1181 isz = 0; 1182 PetscMemzero(table,(m/BITSPERBYTE +1)*sizeof(char)); 1183 1184 /* Extract the indices, assume there can be duplicate entries */ 1185 ierr = ISGetIndices(is[i],&idx); CHKERRQ(ierr); 1186 ierr = ISGetLocalSize(is[i],&n); CHKERRQ(ierr); 1187 1188 /* Enter these into the temp arrays i.e mark table[row], enter row into new index */ 1189 for ( j=0; j<n ; ++j){ 1190 if(!BT_LOOKUP(table, idx[j])) { nidx[isz++] = idx[j];} 1191 } 1192 ierr = ISRestoreIndices(is[i],&idx); CHKERRQ(ierr); 1193 ierr = ISDestroy(is[i]); CHKERRQ(ierr); 1194 1195 k = 0; 1196 for ( j=0; j<ov; j++){ /* for each overlap*/ 1197 n = isz; 1198 for ( ; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */ 1199 row = nidx[k]; 1200 start = ai[row]; 1201 end = ai[row+1]; 1202 for ( l = start; l<end ; l++){ 1203 val = aj[l] + shift; 1204 if (!BT_LOOKUP(table,val)) {nidx[isz++] = val;} 1205 } 1206 } 1207 } 1208 ierr = ISCreateSeq(MPI_COMM_SELF, isz, nidx, (is+i)); CHKERRQ(ierr); 1209 } 1210 PetscFree(table); 1211 PetscFree(nidx); 1212 return 0; 1213 } 1214 1215 int MatPrintHelp_SeqAIJ(Mat A) 1216 { 1217 static int called = 0; 1218 MPI_Comm comm = A->comm; 1219 1220 if (called) return 0; else called = 1; 1221 MPIU_printf(comm," Options for MATSEQAIJ and MATMPIAIJ matrix formats (the defaults):\n"); 1222 MPIU_printf(comm," -mat_lu_pivotthreshold <threshold>\n"); 1223 MPIU_printf(comm," -mat_aij_oneindex - internal indices begin at 1 instead of the default 0.\n"); 1224 MPIU_printf(comm," -mat_aij_no_inode - Do not use inodes\n"); 1225 MPIU_printf(comm," -mat_aij_inode_limit <limit> - Set inode limit (max limit=5)\n"); 1226 #if defined(HAVE_ESSL) 1227 MPIU_printf(comm," -mat_aij_essl - Use IBM sparse LU factorization and solve.\n"); 1228 #endif 1229 return 0; 1230 } 1231 int MatEqual_SeqAIJ(Mat A,Mat B, int* flg); 1232 /* -------------------------------------------------------------------*/ 1233 static struct _MatOps MatOps = {MatSetValues_SeqAIJ, 1234 MatGetRow_SeqAIJ,MatRestoreRow_SeqAIJ, 1235 MatMult_SeqAIJ,MatMultAdd_SeqAIJ, 1236 MatMultTrans_SeqAIJ,MatMultTransAdd_SeqAIJ, 1237 MatSolve_SeqAIJ,MatSolveAdd_SeqAIJ, 1238 MatSolveTrans_SeqAIJ,MatSolveTransAdd_SeqAIJ, 1239 MatLUFactor_SeqAIJ,0, 1240 MatRelax_SeqAIJ, 1241 MatTranspose_SeqAIJ, 1242 MatGetInfo_SeqAIJ,MatEqual_SeqAIJ, 1243 MatGetDiagonal_SeqAIJ,MatDiagonalScale_SeqAIJ,MatNorm_SeqAIJ, 1244 0,MatAssemblyEnd_SeqAIJ, 1245 MatCompress_SeqAIJ, 1246 MatSetOption_SeqAIJ,MatZeroEntries_SeqAIJ,MatZeroRows_SeqAIJ, 1247 MatGetReordering_SeqAIJ, 1248 MatLUFactorSymbolic_SeqAIJ,MatLUFactorNumeric_SeqAIJ,0,0, 1249 MatGetSize_SeqAIJ,MatGetSize_SeqAIJ,MatGetOwnershipRange_SeqAIJ, 1250 MatILUFactorSymbolic_SeqAIJ,0, 1251 0,0,MatConvert_SeqAIJ, 1252 MatGetSubMatrix_SeqAIJ,0, 1253 MatConvertSameType_SeqAIJ,0,0, 1254 MatILUFactor_SeqAIJ,0,0, 1255 MatGetSubMatrices_SeqAIJ,MatIncreaseOverlap_SeqAIJ, 1256 MatGetValues_SeqAIJ,0, 1257 MatPrintHelp_SeqAIJ, 1258 MatScale_SeqAIJ}; 1259 1260 extern int MatUseSuperLU_SeqAIJ(Mat); 1261 extern int MatUseEssl_SeqAIJ(Mat); 1262 extern int MatUseDXML_SeqAIJ(Mat); 1263 1264 /*@C 1265 MatCreateSeqAIJ - Creates a sparse matrix in AIJ (compressed row) format 1266 (the default parallel PETSc format). For good matrix assembly performance 1267 the user should preallocate the matrix storage by setting the parameter nz 1268 (or nzz). By setting these parameters accurately, performance can be 1269 increased by more than a factor of 50. 1270 1271 Input Parameters: 1272 . comm - MPI communicator, set to MPI_COMM_SELF 1273 . m - number of rows 1274 . n - number of columns 1275 . nz - number of nonzeros per row (same for all rows) 1276 . nzz - number of nonzeros per row or null (possibly different for each row) 1277 1278 Output Parameter: 1279 . A - the matrix 1280 1281 Notes: 1282 The AIJ format (also called the Yale sparse matrix format or 1283 compressed row storage), is fully compatible with standard Fortran 77 1284 storage. That is, the stored row and column indices can begin at 1285 either one (as in Fortran) or zero. See the users manual for details. 1286 1287 Specify the preallocated storage with either nz or nnz (not both). 1288 Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory 1289 allocation. For additional details, see the users manual chapter on 1290 matrices and the file $(PETSC_DIR)/Performance. 1291 1292 By default, this format uses inodes (identical nodes) when possible, to 1293 improve numerical efficiency of Matrix vector products and solves. We 1294 search for consecutive rows with the same nonzero structure, thereby 1295 reusing matrix information to achieve increased efficiency. 1296 1297 Options Database Keys: 1298 $ -mat_aij_no_inode - Do not use inodes 1299 $ -mat_aij_inode_limit <limit> - Set inode limit. 1300 $ (max limit=5) 1301 $ -mat_aij_oneindex - Internally use indexing starting at 1 1302 $ rather than 0. Note: When calling MatSetValues(), 1303 $ the user still MUST index entries starting at 0! 1304 1305 .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues() 1306 @*/ 1307 int MatCreateSeqAIJ(MPI_Comm comm,int m,int n,int nz,int *nnz, Mat *A) 1308 { 1309 Mat B; 1310 Mat_SeqAIJ *b; 1311 int i,len,ierr, flg; 1312 1313 *A = 0; 1314 PetscHeaderCreate(B,_Mat,MAT_COOKIE,MATSEQAIJ,comm); 1315 PLogObjectCreate(B); 1316 B->data = (void *) (b = PetscNew(Mat_SeqAIJ)); CHKPTRQ(b); 1317 PetscMemcpy(&B->ops,&MatOps,sizeof(struct _MatOps)); 1318 B->destroy = MatDestroy_SeqAIJ; 1319 B->view = MatView_SeqAIJ; 1320 B->factor = 0; 1321 B->lupivotthreshold = 1.0; 1322 ierr = OptionsGetDouble(PETSC_NULL,"-mat_lu_pivotthreshold",&B->lupivotthreshold, 1323 &flg); CHKERRQ(ierr); 1324 b->ilu_preserve_row_sums = PETSC_FALSE; 1325 ierr = OptionsHasName(PETSC_NULL,"-pc_ilu_preserve_row_sums", 1326 (int*) &b->ilu_preserve_row_sums); CHKERRQ(ierr); 1327 b->row = 0; 1328 b->col = 0; 1329 b->indexshift = 0; 1330 b->reallocs = 0; 1331 ierr = OptionsHasName(PETSC_NULL,"-mat_aij_oneindex", &flg); CHKERRQ(ierr); 1332 if (flg) b->indexshift = -1; 1333 1334 b->m = m; 1335 b->n = n; 1336 b->imax = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(b->imax); 1337 if (nnz == PETSC_NULL) { 1338 if (nz == PETSC_DEFAULT) nz = 10; 1339 else if (nz <= 0) nz = 1; 1340 for ( i=0; i<m; i++ ) b->imax[i] = nz; 1341 nz = nz*m; 1342 } 1343 else { 1344 nz = 0; 1345 for ( i=0; i<m; i++ ) {b->imax[i] = nnz[i]; nz += nnz[i];} 1346 } 1347 1348 /* allocate the matrix space */ 1349 len = nz*(sizeof(int) + sizeof(Scalar)) + (b->m+1)*sizeof(int); 1350 b->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(b->a); 1351 b->j = (int *) (b->a + nz); 1352 PetscMemzero(b->j,nz*sizeof(int)); 1353 b->i = b->j + nz; 1354 b->singlemalloc = 1; 1355 1356 b->i[0] = -b->indexshift; 1357 for (i=1; i<m+1; i++) { 1358 b->i[i] = b->i[i-1] + b->imax[i-1]; 1359 } 1360 1361 /* b->ilen will count nonzeros in each row so far. */ 1362 b->ilen = (int *) PetscMalloc((m+1)*sizeof(int)); 1363 PLogObjectMemory(B,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ)); 1364 for ( i=0; i<b->m; i++ ) { b->ilen[i] = 0;} 1365 1366 b->nz = 0; 1367 b->maxnz = nz; 1368 b->sorted = 0; 1369 b->roworiented = 1; 1370 b->nonew = 0; 1371 b->diag = 0; 1372 b->solve_work = 0; 1373 b->spptr = 0; 1374 b->inode.node_count = 0; 1375 b->inode.size = 0; 1376 b->inode.limit = 5; 1377 b->inode.max_limit = 5; 1378 1379 *A = B; 1380 /* SuperLU is not currently supported through PETSc */ 1381 #if defined(HAVE_SUPERLU) 1382 ierr = OptionsHasName(PETSC_NULL,"-mat_aij_superlu", &flg); CHKERRQ(ierr); 1383 if (flg) { ierr = MatUseSuperLU_SeqAIJ(B); CHKERRQ(ierr); } 1384 #endif 1385 ierr = OptionsHasName(PETSC_NULL,"-mat_aij_essl", &flg); CHKERRQ(ierr); 1386 if (flg) { ierr = MatUseEssl_SeqAIJ(B); CHKERRQ(ierr); } 1387 ierr = OptionsHasName(PETSC_NULL,"-mat_aij_dxml", &flg); CHKERRQ(ierr); 1388 if (flg) { 1389 if (!b->indexshift) SETERRQ(1,"MatCreateSeqAIJ:need -mat_aij_oneindex with -mat_aij_dxml"); 1390 ierr = MatUseDXML_SeqAIJ(B); CHKERRQ(ierr); 1391 } 1392 ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr); 1393 if (flg) {ierr = MatPrintHelp(B); CHKERRQ(ierr); } 1394 return 0; 1395 } 1396 1397 int MatConvertSameType_SeqAIJ(Mat A,Mat *B,int cpvalues) 1398 { 1399 Mat C; 1400 Mat_SeqAIJ *c,*a = (Mat_SeqAIJ *) A->data; 1401 int i,len, m = a->m,shift = a->indexshift; 1402 1403 *B = 0; 1404 PetscHeaderCreate(C,_Mat,MAT_COOKIE,MATSEQAIJ,A->comm); 1405 PLogObjectCreate(C); 1406 C->data = (void *) (c = PetscNew(Mat_SeqAIJ)); CHKPTRQ(c); 1407 PetscMemcpy(&C->ops,&A->ops,sizeof(struct _MatOps)); 1408 C->destroy = MatDestroy_SeqAIJ; 1409 C->view = MatView_SeqAIJ; 1410 C->factor = A->factor; 1411 c->row = 0; 1412 c->col = 0; 1413 c->indexshift = shift; 1414 C->assembled = PETSC_TRUE; 1415 1416 c->m = a->m; 1417 c->n = a->n; 1418 1419 c->imax = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->imax); 1420 c->ilen = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(c->ilen); 1421 for ( i=0; i<m; i++ ) { 1422 c->imax[i] = a->imax[i]; 1423 c->ilen[i] = a->ilen[i]; 1424 } 1425 1426 /* allocate the matrix space */ 1427 c->singlemalloc = 1; 1428 len = (m+1)*sizeof(int)+(a->i[m])*(sizeof(Scalar)+sizeof(int)); 1429 c->a = (Scalar *) PetscMalloc( len ); CHKPTRQ(c->a); 1430 c->j = (int *) (c->a + a->i[m] + shift); 1431 c->i = c->j + a->i[m] + shift; 1432 PetscMemcpy(c->i,a->i,(m+1)*sizeof(int)); 1433 if (m > 0) { 1434 PetscMemcpy(c->j,a->j,(a->i[m]+shift)*sizeof(int)); 1435 if (cpvalues == COPY_VALUES) { 1436 PetscMemcpy(c->a,a->a,(a->i[m]+shift)*sizeof(Scalar)); 1437 } 1438 } 1439 1440 PLogObjectMemory(C,len+2*(m+1)*sizeof(int)+sizeof(struct _Mat)+sizeof(Mat_SeqAIJ)); 1441 c->sorted = a->sorted; 1442 c->roworiented = a->roworiented; 1443 c->nonew = a->nonew; 1444 c->ilu_preserve_row_sums = a->ilu_preserve_row_sums; 1445 1446 if (a->diag) { 1447 c->diag = (int *) PetscMalloc( (m+1)*sizeof(int) ); CHKPTRQ(c->diag); 1448 PLogObjectMemory(C,(m+1)*sizeof(int)); 1449 for ( i=0; i<m; i++ ) { 1450 c->diag[i] = a->diag[i]; 1451 } 1452 } 1453 else c->diag = 0; 1454 c->inode.limit = a->inode.limit; 1455 c->inode.max_limit = a->inode.max_limit; 1456 if (a->inode.size){ 1457 c->inode.size = (int *) PetscMalloc( m *sizeof(int) ); CHKPTRQ(c->inode.size); 1458 c->inode.node_count = a->inode.node_count; 1459 PetscMemcpy( c->inode.size, a->inode.size, m*sizeof(int)); 1460 } else { 1461 c->inode.size = 0; 1462 c->inode.node_count = 0; 1463 } 1464 c->nz = a->nz; 1465 c->maxnz = a->maxnz; 1466 c->solve_work = 0; 1467 c->spptr = 0; /* Dangerous -I'm throwing away a->spptr */ 1468 1469 *B = C; 1470 return 0; 1471 } 1472 1473 int MatLoad_SeqAIJ(Viewer bview,MatType type,Mat *A) 1474 { 1475 Mat_SeqAIJ *a; 1476 Mat B; 1477 int i, nz, ierr, fd, header[4],size,*rowlengths = 0,M,N,shift; 1478 PetscObject vobj = (PetscObject) bview; 1479 MPI_Comm comm = vobj->comm; 1480 1481 MPI_Comm_size(comm,&size); 1482 if (size > 1) SETERRQ(1,"MatLoad_SeqAIJ:view must have one processor"); 1483 ierr = ViewerFileGetDescriptor_Private(bview,&fd); CHKERRQ(ierr); 1484 ierr = SYRead(fd,header,4,SYINT); CHKERRQ(ierr); 1485 if (header[0] != MAT_COOKIE) SETERRQ(1,"MatLoad_SeqAIJ:not matrix object in file"); 1486 M = header[1]; N = header[2]; nz = header[3]; 1487 1488 /* read in row lengths */ 1489 rowlengths = (int*) PetscMalloc( M*sizeof(int) ); CHKPTRQ(rowlengths); 1490 ierr = SYRead(fd,rowlengths,M,SYINT); CHKERRQ(ierr); 1491 1492 /* create our matrix */ 1493 ierr = MatCreateSeqAIJ(comm,M,N,0,rowlengths,A); CHKERRQ(ierr); 1494 B = *A; 1495 a = (Mat_SeqAIJ *) B->data; 1496 shift = a->indexshift; 1497 1498 /* read in column indices and adjust for Fortran indexing*/ 1499 ierr = SYRead(fd,a->j,nz,SYINT); CHKERRQ(ierr); 1500 if (shift) { 1501 for ( i=0; i<nz; i++ ) { 1502 a->j[i] += 1; 1503 } 1504 } 1505 1506 /* read in nonzero values */ 1507 ierr = SYRead(fd,a->a,nz,SYSCALAR); CHKERRQ(ierr); 1508 1509 /* set matrix "i" values */ 1510 a->i[0] = -shift; 1511 for ( i=1; i<= M; i++ ) { 1512 a->i[i] = a->i[i-1] + rowlengths[i-1]; 1513 a->ilen[i-1] = rowlengths[i-1]; 1514 } 1515 PetscFree(rowlengths); 1516 1517 ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 1518 ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 1519 return 0; 1520 } 1521 1522 1523 1524 /* 1525 1526 flg =0 if error; 1527 flg =1 if both are equal; 1528 */ 1529 int MatEqual_SeqAIJ(Mat A,Mat B, int* flg) 1530 { 1531 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data, *b = (Mat_SeqAIJ *)B->data; 1532 1533 if (B->type !=MATSEQAIJ) SETERRQ(1,"MatEqual_SeqAIJ:Both matrices should be of type MATSEQAIJ"); 1534 1535 /* If the matrix dimensions are not equal, or no of nonzeros or shift */ 1536 if ((a->m != b->m ) || (a->n !=b->n) ||( a->nz != b->nz)|| 1537 (a->indexshift != b->indexshift)) { *flg =0 ; return 0; } 1538 1539 /* if the a->i are the same */ 1540 if(PetscMemcmp((char *)a->i, (char*)b->i, (a->n+1)*sizeof(int))) { *flg =0 ; return 0;} 1541 1542 /* if a->j are the same */ 1543 if(PetscMemcmp((char *)a->j, (char*)b->j, (a->nz)*sizeof(int))) { 1544 *flg =0 ; return 0; 1545 } 1546 1547 /* if a->j are the same */ 1548 if(PetscMemcmp((char *)a->a, (char*)b->a, (a->nz)*sizeof(Scalar))) { 1549 *flg =0 ; return 0; 1550 } 1551 *flg =1 ; 1552 return 0; 1553 1554 } 1555