#define PETSCMAT_DLL /* Defines the basic matrix operations for the AIJ (compressed row) matrix storage format. */ #include "src/mat/impls/aij/seq/aij.h" /*I "petscmat.h" I*/ #include "src/inline/spops.h" #include "src/inline/dot.h" #include "petscbt.h" #undef __FUNCT__ #define __FUNCT__ "MatDiagonalSet_SeqAIJ" PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalSet_SeqAIJ(Mat Y,Vec D,InsertMode is) { PetscErrorCode ierr; Mat_SeqAIJ *aij = (Mat_SeqAIJ*) Y->data; PetscInt i,*diag, m = Y->rmap.n; PetscScalar *v,*aa = aij->a; PetscTruth missing; PetscFunctionBegin; if (Y->assembled) { ierr = MatMissingDiagonal_SeqAIJ(Y,&missing,PETSC_NULL);CHKERRQ(ierr); if (!missing) { diag = aij->diag; ierr = VecGetArray(D,&v);CHKERRQ(ierr); if (is == INSERT_VALUES) { for (i=0; idata; PetscErrorCode ierr; PetscInt i,ishift; PetscFunctionBegin; *m = A->rmap.n; if (!ia) PetscFunctionReturn(0); ishift = 0; if (symmetric && !A->structurally_symmetric) { ierr = MatToSymmetricIJ_SeqAIJ(A->rmap.n,a->i,a->j,ishift,oshift,ia,ja);CHKERRQ(ierr); } else if (oshift == 1) { PetscInt nz = a->i[A->rmap.n]; /* malloc space and add 1 to i and j indices */ ierr = PetscMalloc((A->rmap.n+1)*sizeof(PetscInt),ia);CHKERRQ(ierr); ierr = PetscMalloc((nz+1)*sizeof(PetscInt),ja);CHKERRQ(ierr); for (i=0; ij[i] + 1; for (i=0; irmap.n+1; i++) (*ia)[i] = a->i[i] + 1; } else { *ia = a->i; *ja = a->j; } PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatRestoreRowIJ_SeqAIJ" PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscInt *n,PetscInt *ia[],PetscInt *ja[],PetscTruth *done) { PetscErrorCode ierr; PetscFunctionBegin; if (!ia) PetscFunctionReturn(0); if ((symmetric && !A->structurally_symmetric) || oshift == 1) { ierr = PetscFree(*ia);CHKERRQ(ierr); ierr = PetscFree(*ja);CHKERRQ(ierr); } PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatGetColumnIJ_SeqAIJ" PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat A,PetscInt oshift,PetscTruth symmetric,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscTruth *done) { Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; PetscErrorCode ierr; PetscInt i,*collengths,*cia,*cja,n = A->cmap.n,m = A->rmap.n; PetscInt nz = a->i[m],row,*jj,mr,col; PetscFunctionBegin; *nn = n; if (!ia) PetscFunctionReturn(0); if (symmetric) { ierr = MatToSymmetricIJ_SeqAIJ(A->rmap.n,a->i,a->j,0,oshift,ia,ja);CHKERRQ(ierr); } else { ierr = PetscMalloc((n+1)*sizeof(PetscInt),&collengths);CHKERRQ(ierr); ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr); ierr = PetscMalloc((n+1)*sizeof(PetscInt),&cia);CHKERRQ(ierr); ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cja);CHKERRQ(ierr); jj = a->j; for (i=0; ij; for (row=0; rowi[row+1] - a->i[row]; for (i=0; idata; PetscInt *ai = a->i; PetscErrorCode ierr; PetscFunctionBegin; ierr = PetscMemcpy(a->a+ai[row],v,(ai[row+1]-ai[row])*sizeof(PetscScalar));CHKERRQ(ierr); PetscFunctionReturn(0); } #define CHUNKSIZE 15 #undef __FUNCT__ #define __FUNCT__ "MatSetValues_SeqAIJ" PetscErrorCode MatSetValues_SeqAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],const PetscScalar v[],InsertMode is) { Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; PetscInt *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N; PetscInt *imax = a->imax,*ai = a->i,*ailen = a->ilen; PetscErrorCode ierr; PetscInt *aj = a->j,nonew = a->nonew,lastcol = -1; PetscScalar *ap,value,*aa = a->a; PetscTruth ignorezeroentries = ((a->ignorezeroentries && is == ADD_VALUES) ? PETSC_TRUE:PETSC_FALSE); PetscTruth roworiented = a->roworiented; PetscFunctionBegin; for (k=0; k= A->rmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap.n-1); #endif rp = aj + ai[row]; ap = aa + ai[row]; rmax = imax[row]; nrow = ailen[row]; low = 0; high = nrow; for (l=0; l= A->cmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->cmap.n-1); #endif col = in[l]; if (roworiented) { value = v[l + k*n]; } else { value = v[k + l*m]; } if (value == 0.0 && ignorezeroentries && (is == ADD_VALUES)) continue; if (col <= lastcol) low = 0; else high = nrow; lastcol = col; while (high-low > 5) { t = (low+high)/2; if (rp[t] > col) high = t; else low = t; } for (i=low; i col) break; if (rp[i] == col) { if (is == ADD_VALUES) ap[i] += value; else ap[i] = value; goto noinsert; } } if (value == 0.0 && ignorezeroentries) goto noinsert; if (nonew == 1) goto noinsert; if (nonew == -1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero at (%D,%D) in the matrix",row,col); MatSeqXAIJReallocateAIJ(A,A->rmap.n,1,nrow,row,col,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar); N = nrow++ - 1; a->nz++; high++; /* shift up all the later entries in this row */ for (ii=N; ii>=i; ii--) { rp[ii+1] = rp[ii]; ap[ii+1] = ap[ii]; } rp[i] = col; ap[i] = value; noinsert:; low = i + 1; } ailen[row] = nrow; } A->same_nonzero = PETSC_FALSE; PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatGetValues_SeqAIJ" PetscErrorCode MatGetValues_SeqAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt n,const PetscInt in[],PetscScalar v[]) { Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; PetscInt *rp,k,low,high,t,row,nrow,i,col,l,*aj = a->j; PetscInt *ai = a->i,*ailen = a->ilen; PetscScalar *ap,*aa = a->a,zero = 0.0; PetscFunctionBegin; for (k=0; k= A->rmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",row,A->rmap.n-1); rp = aj + ai[row]; ap = aa + ai[row]; nrow = ailen[row]; for (l=0; l= A->cmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",in[l],A->cmap.n-1); col = in[l] ; high = nrow; low = 0; /* assume unsorted */ while (high-low > 5) { t = (low+high)/2; if (rp[t] > col) high = t; else low = t; } for (i=low; i col) break; if (rp[i] == col) { *v++ = ap[i]; goto finished; } } *v++ = zero; finished:; } } PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatView_SeqAIJ_Binary" PetscErrorCode MatView_SeqAIJ_Binary(Mat A,PetscViewer viewer) { Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; PetscErrorCode ierr; PetscInt i,*col_lens; int fd; PetscFunctionBegin; ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); ierr = PetscMalloc((4+A->rmap.n)*sizeof(PetscInt),&col_lens);CHKERRQ(ierr); col_lens[0] = MAT_FILE_COOKIE; col_lens[1] = A->rmap.n; col_lens[2] = A->cmap.n; col_lens[3] = a->nz; /* store lengths of each row and write (including header) to file */ for (i=0; irmap.n; i++) { col_lens[4+i] = a->i[i+1] - a->i[i]; } ierr = PetscBinaryWrite(fd,col_lens,4+A->rmap.n,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr); ierr = PetscFree(col_lens);CHKERRQ(ierr); /* store column indices (zero start index) */ ierr = PetscBinaryWrite(fd,a->j,a->nz,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr); /* store nonzero values */ ierr = PetscBinaryWrite(fd,a->a,a->nz,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr); PetscFunctionReturn(0); } EXTERN PetscErrorCode MatSeqAIJFactorInfo_Matlab(Mat,PetscViewer); #undef __FUNCT__ #define __FUNCT__ "MatView_SeqAIJ_ASCII" PetscErrorCode MatView_SeqAIJ_ASCII(Mat A,PetscViewer viewer) { Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; PetscErrorCode ierr; PetscInt i,j,m = A->rmap.n,shift=0; const char *name; PetscViewerFormat format; PetscFunctionBegin; ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); if (format == PETSC_VIEWER_ASCII_MATLAB) { PetscInt nofinalvalue = 0; if ((a->i[m] == a->i[m-1]) || (a->j[a->nz-1] != A->cmap.n-!shift)) { nofinalvalue = 1; } ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",m,A->cmap.n);CHKERRQ(ierr); ierr = PetscViewerASCIIPrintf(viewer,"%% Nonzeros = %D \n",a->nz);CHKERRQ(ierr); ierr = PetscViewerASCIIPrintf(viewer,"zzz = zeros(%D,3);\n",a->nz+nofinalvalue);CHKERRQ(ierr); ierr = PetscViewerASCIIPrintf(viewer,"zzz = [\n");CHKERRQ(ierr); for (i=0; ii[i]+shift; ji[i+1]+shift; j++) { #if defined(PETSC_USE_COMPLEX) ierr = PetscViewerASCIIPrintf(viewer,"%D %D %18.16e + %18.16ei \n",i+1,a->j[j]+!shift,PetscRealPart(a->a[j]),PetscImaginaryPart(a->a[j]));CHKERRQ(ierr); #else ierr = PetscViewerASCIIPrintf(viewer,"%D %D %18.16e\n",i+1,a->j[j]+!shift,a->a[j]);CHKERRQ(ierr); #endif } } if (nofinalvalue) { ierr = PetscViewerASCIIPrintf(viewer,"%D %D %18.16e\n",m,A->cmap.n,0.0);CHKERRQ(ierr); } ierr = PetscViewerASCIIPrintf(viewer,"];\n %s = spconvert(zzz);\n",name);CHKERRQ(ierr); ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); } else if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) { PetscFunctionReturn(0); } else if (format == PETSC_VIEWER_ASCII_COMMON) { ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); for (i=0; ii[i]+shift; ji[i+1]+shift; j++) { #if defined(PETSC_USE_COMPLEX) if (PetscImaginaryPart(a->a[j]) > 0.0 && PetscRealPart(a->a[j]) != 0.0) { ierr = PetscViewerASCIIPrintf(viewer," (%D, %G + %G i)",a->j[j]+shift,PetscRealPart(a->a[j]),PetscImaginaryPart(a->a[j]));CHKERRQ(ierr); } else if (PetscImaginaryPart(a->a[j]) < 0.0 && PetscRealPart(a->a[j]) != 0.0) { ierr = PetscViewerASCIIPrintf(viewer," (%D, %G - %G i)",a->j[j]+shift,PetscRealPart(a->a[j]),-PetscImaginaryPart(a->a[j]));CHKERRQ(ierr); } else if (PetscRealPart(a->a[j]) != 0.0) { ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",a->j[j]+shift,PetscRealPart(a->a[j]));CHKERRQ(ierr); } #else if (a->a[j] != 0.0) {ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",a->j[j]+shift,a->a[j]);CHKERRQ(ierr);} #endif } ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); } ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); } else if (format == PETSC_VIEWER_ASCII_SYMMODU) { PetscInt nzd=0,fshift=1,*sptr; ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); ierr = PetscMalloc((m+1)*sizeof(PetscInt),&sptr);CHKERRQ(ierr); for (i=0; ii[i]+shift; ji[i+1]+shift; j++) { if (a->j[j] >= i) { #if defined(PETSC_USE_COMPLEX) if (PetscImaginaryPart(a->a[j]) != 0.0 || PetscRealPart(a->a[j]) != 0.0) nzd++; #else if (a->a[j] != 0.0) nzd++; #endif } } } sptr[m] = nzd+1; ierr = PetscViewerASCIIPrintf(viewer," %D %D\n\n",m,nzd);CHKERRQ(ierr); for (i=0; ii[i]+shift; ji[i+1]+shift; j++) { if (a->j[j] >= i) {ierr = PetscViewerASCIIPrintf(viewer," %D ",a->j[j]+fshift);CHKERRQ(ierr);} } ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); } ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); for (i=0; ii[i]+shift; ji[i+1]+shift; j++) { if (a->j[j] >= i) { #if defined(PETSC_USE_COMPLEX) if (PetscImaginaryPart(a->a[j]) != 0.0 || PetscRealPart(a->a[j]) != 0.0) { ierr = PetscViewerASCIIPrintf(viewer," %18.16e %18.16e ",PetscRealPart(a->a[j]),PetscImaginaryPart(a->a[j]));CHKERRQ(ierr); } #else if (a->a[j] != 0.0) {ierr = PetscViewerASCIIPrintf(viewer," %18.16e ",a->a[j]);CHKERRQ(ierr);} #endif } } ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); } ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); } else if (format == PETSC_VIEWER_ASCII_DENSE) { PetscInt cnt = 0,jcnt; PetscScalar value; ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); for (i=0; icmap.n; j++) { if (jcnt < a->i[i+1]-a->i[i] && j == a->j[cnt]) { value = a->a[cnt++]; jcnt++; } else { value = 0.0; } #if defined(PETSC_USE_COMPLEX) ierr = PetscViewerASCIIPrintf(viewer," %7.5e+%7.5e i ",PetscRealPart(value),PetscImaginaryPart(value));CHKERRQ(ierr); #else ierr = PetscViewerASCIIPrintf(viewer," %7.5e ",value);CHKERRQ(ierr); #endif } ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); } ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); } else { ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); for (i=0; ii[i]+shift; ji[i+1]+shift; j++) { #if defined(PETSC_USE_COMPLEX) if (PetscImaginaryPart(a->a[j]) > 0.0) { ierr = PetscViewerASCIIPrintf(viewer," (%D, %G + %G i)",a->j[j]+shift,PetscRealPart(a->a[j]),PetscImaginaryPart(a->a[j]));CHKERRQ(ierr); } else if (PetscImaginaryPart(a->a[j]) < 0.0) { ierr = PetscViewerASCIIPrintf(viewer," (%D, %G - %G i)",a->j[j]+shift,PetscRealPart(a->a[j]),-PetscImaginaryPart(a->a[j]));CHKERRQ(ierr); } else { ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",a->j[j]+shift,PetscRealPart(a->a[j]));CHKERRQ(ierr); } #else ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",a->j[j]+shift,a->a[j]);CHKERRQ(ierr); #endif } ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); } ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); } ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatView_SeqAIJ_Draw_Zoom" PetscErrorCode MatView_SeqAIJ_Draw_Zoom(PetscDraw draw,void *Aa) { Mat A = (Mat) Aa; Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; PetscErrorCode ierr; PetscInt i,j,m = A->rmap.n,color; PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r,maxv = 0.0; PetscViewer viewer; PetscViewerFormat format; PetscFunctionBegin; ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); /* loop over matrix elements drawing boxes */ if (format != PETSC_VIEWER_DRAW_CONTOUR) { /* Blue for negative, Cyan for zero and Red for positive */ color = PETSC_DRAW_BLUE; for (i=0; ii[i]; ji[i+1]; j++) { x_l = a->j[j] ; x_r = x_l + 1.0; #if defined(PETSC_USE_COMPLEX) if (PetscRealPart(a->a[j]) >= 0.) continue; #else if (a->a[j] >= 0.) continue; #endif ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); } } color = PETSC_DRAW_CYAN; for (i=0; ii[i]; ji[i+1]; j++) { x_l = a->j[j]; x_r = x_l + 1.0; if (a->a[j] != 0.) continue; ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); } } color = PETSC_DRAW_RED; for (i=0; ii[i]; ji[i+1]; j++) { x_l = a->j[j]; x_r = x_l + 1.0; #if defined(PETSC_USE_COMPLEX) if (PetscRealPart(a->a[j]) <= 0.) continue; #else if (a->a[j] <= 0.) continue; #endif ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); } } } else { /* use contour shading to indicate magnitude of values */ /* first determine max of all nonzero values */ PetscInt nz = a->nz,count; PetscDraw popup; PetscReal scale; for (i=0; ia[i]) > maxv) maxv = PetscAbsScalar(a->a[i]); } scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv; ierr = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr); if (popup) {ierr = PetscDrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);} count = 0; for (i=0; ii[i]; ji[i+1]; j++) { x_l = a->j[j]; x_r = x_l + 1.0; color = PETSC_DRAW_BASIC_COLORS + (PetscInt)(scale*PetscAbsScalar(a->a[count])); ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); count++; } } } PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatView_SeqAIJ_Draw" PetscErrorCode MatView_SeqAIJ_Draw(Mat A,PetscViewer viewer) { PetscErrorCode ierr; PetscDraw draw; PetscReal xr,yr,xl,yl,h,w; PetscTruth isnull; PetscFunctionBegin; ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0); ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); xr = A->cmap.n; yr = A->rmap.n; h = yr/10.0; w = xr/10.0; xr += w; yr += h; xl = -w; yl = -h; ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); ierr = PetscDrawZoom(draw,MatView_SeqAIJ_Draw_Zoom,A);CHKERRQ(ierr); ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatView_SeqAIJ" PetscErrorCode MatView_SeqAIJ(Mat A,PetscViewer viewer) { PetscErrorCode ierr; PetscTruth iascii,isbinary,isdraw; PetscFunctionBegin; ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr); ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); if (iascii) { ierr = MatView_SeqAIJ_ASCII(A,viewer);CHKERRQ(ierr); } else if (isbinary) { ierr = MatView_SeqAIJ_Binary(A,viewer);CHKERRQ(ierr); } else if (isdraw) { ierr = MatView_SeqAIJ_Draw(A,viewer);CHKERRQ(ierr); } else { SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by SeqAIJ matrices",((PetscObject)viewer)->type_name); } ierr = MatView_Inode(A,viewer);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatAssemblyEnd_SeqAIJ" PetscErrorCode MatAssemblyEnd_SeqAIJ(Mat A,MatAssemblyType mode) { Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; PetscErrorCode ierr; PetscInt fshift = 0,i,j,*ai = a->i,*aj = a->j,*imax = a->imax; PetscInt m = A->rmap.n,*ip,N,*ailen = a->ilen,rmax = 0; PetscScalar *aa = a->a,*ap; PetscReal ratio=0.6; PetscFunctionBegin; if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(0); if (m) rmax = ailen[0]; /* determine row with most nonzeros */ for (i=1; inz = ai[m]; ierr = MatMarkDiagonal_SeqAIJ(A);CHKERRQ(ierr); ierr = PetscInfo4(A,"Matrix size: %D X %D; storage space: %D unneeded,%D used\n",m,A->cmap.n,fshift,a->nz);CHKERRQ(ierr); ierr = PetscInfo1(A,"Number of mallocs during MatSetValues() is %D\n",a->reallocs);CHKERRQ(ierr); ierr = PetscInfo1(A,"Maximum nonzeros in any row is %D\n",rmax);CHKERRQ(ierr); a->reallocs = 0; A->info.nz_unneeded = (double)fshift; a->rmax = rmax; /* check for zero rows. If found a large number of zero rows, use CompressedRow functions */ ierr = Mat_CheckCompressedRow(A,&a->compressedrow,a->i,m,ratio);CHKERRQ(ierr); A->same_nonzero = PETSC_TRUE; ierr = MatAssemblyEnd_Inode(A,mode);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatRealPart_SeqAIJ" PetscErrorCode MatRealPart_SeqAIJ(Mat A) { Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; PetscInt i,nz = a->nz; PetscScalar *aa = a->a; PetscFunctionBegin; for (i=0; idata; PetscInt i,nz = a->nz; PetscScalar *aa = a->a; PetscFunctionBegin; for (i=0; idata; PetscErrorCode ierr; PetscFunctionBegin; ierr = PetscMemzero(a->a,(a->i[A->rmap.n])*sizeof(PetscScalar));CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatDestroy_SeqAIJ" PetscErrorCode MatDestroy_SeqAIJ(Mat A) { Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; PetscErrorCode ierr; PetscFunctionBegin; #if defined(PETSC_USE_LOG) PetscLogObjectState((PetscObject)A,"Rows=%D, Cols=%D, NZ=%D",A->rmap.n,A->cmap.n,a->nz); #endif ierr = MatSeqXAIJFreeAIJ(A,&a->a,&a->j,&a->i);CHKERRQ(ierr); if (a->row) { ierr = ISDestroy(a->row);CHKERRQ(ierr); } if (a->col) { ierr = ISDestroy(a->col);CHKERRQ(ierr); } ierr = PetscFree(a->diag);CHKERRQ(ierr); ierr = PetscFree2(a->imax,a->ilen);CHKERRQ(ierr); ierr = PetscFree(a->idiag);CHKERRQ(ierr); ierr = PetscFree(a->solve_work);CHKERRQ(ierr); if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);} ierr = PetscFree(a->saved_values);CHKERRQ(ierr); if (a->coloring) {ierr = ISColoringDestroy(a->coloring);CHKERRQ(ierr);} ierr = PetscFree(a->xtoy);CHKERRQ(ierr); if (a->compressedrow.use){ierr = PetscFree(a->compressedrow.i);} ierr = MatDestroy_Inode(A);CHKERRQ(ierr); ierr = PetscFree(a);CHKERRQ(ierr); ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr); ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSeqAIJSetColumnIndices_C","",PETSC_NULL);CHKERRQ(ierr); ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatStoreValues_C","",PETSC_NULL);CHKERRQ(ierr); ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatRetrieveValues_C","",PETSC_NULL);CHKERRQ(ierr); ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatConvert_seqaij_seqsbaij_C","",PETSC_NULL);CHKERRQ(ierr); ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatConvert_seqaij_seqbaij_C","",PETSC_NULL);CHKERRQ(ierr); ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatConvert_seqaij_seqcsrperm_C","",PETSC_NULL);CHKERRQ(ierr); ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatIsTranspose_C","",PETSC_NULL);CHKERRQ(ierr); ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSeqAIJSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr); ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSeqAIJSetPreallocationCSR_C","",PETSC_NULL);CHKERRQ(ierr); ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatReorderForNonzeroDiagonal_C","",PETSC_NULL);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatCompress_SeqAIJ" PetscErrorCode MatCompress_SeqAIJ(Mat A) { PetscFunctionBegin; PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatSetOption_SeqAIJ" PetscErrorCode MatSetOption_SeqAIJ(Mat A,MatOption op) { Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; PetscErrorCode ierr; PetscFunctionBegin; switch (op) { case MAT_ROW_ORIENTED: a->roworiented = PETSC_TRUE; break; case MAT_KEEP_ZEROED_ROWS: a->keepzeroedrows = PETSC_TRUE; break; case MAT_COLUMN_ORIENTED: a->roworiented = PETSC_FALSE; break; case MAT_COLUMNS_SORTED: a->sorted = PETSC_TRUE; break; case MAT_COLUMNS_UNSORTED: a->sorted = PETSC_FALSE; break; case MAT_NO_NEW_NONZERO_LOCATIONS: a->nonew = 1; break; case MAT_NEW_NONZERO_LOCATION_ERR: a->nonew = -1; break; case MAT_NEW_NONZERO_ALLOCATION_ERR: a->nonew = -2; break; case MAT_YES_NEW_NONZERO_LOCATIONS: a->nonew = 0; break; case MAT_IGNORE_ZERO_ENTRIES: a->ignorezeroentries = PETSC_TRUE; break; case MAT_USE_COMPRESSEDROW: a->compressedrow.use = PETSC_TRUE; break; case MAT_DO_NOT_USE_COMPRESSEDROW: a->compressedrow.use = PETSC_FALSE; break; case MAT_ROWS_SORTED: case MAT_ROWS_UNSORTED: case MAT_YES_NEW_DIAGONALS: case MAT_IGNORE_OFF_PROC_ENTRIES: case MAT_USE_HASH_TABLE: ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr); break; case MAT_NO_NEW_DIAGONALS: SETERRQ(PETSC_ERR_SUP,"MAT_NO_NEW_DIAGONALS"); default: break; } ierr = MatSetOption_Inode(A,op);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatGetDiagonal_SeqAIJ" PetscErrorCode MatGetDiagonal_SeqAIJ(Mat A,Vec v) { Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; PetscErrorCode ierr; PetscInt i,j,n; PetscScalar *x,zero = 0.0; PetscFunctionBegin; ierr = VecSet(v,zero);CHKERRQ(ierr); ierr = VecGetArray(v,&x);CHKERRQ(ierr); ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr); if (n != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); for (i=0; irmap.n; i++) { for (j=a->i[i]; ji[i+1]; j++) { if (a->j[j] == i) { x[i] = a->a[j]; break; } } } ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatMultTransposeAdd_SeqAIJ" PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A,Vec xx,Vec zz,Vec yy) { Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; PetscScalar *x,*y; PetscErrorCode ierr; PetscInt m = A->rmap.n; #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ) PetscScalar *v,alpha; PetscInt n,i,*idx,*ii,*ridx=PETSC_NULL; Mat_CompressedRow cprow = a->compressedrow; PetscTruth usecprow = cprow.use; #endif PetscFunctionBegin; if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} ierr = VecGetArray(xx,&x);CHKERRQ(ierr); ierr = VecGetArray(yy,&y);CHKERRQ(ierr); #if defined(PETSC_USE_FORTRAN_KERNEL_MULTTRANSPOSEAIJ) fortranmulttransposeaddaij_(&m,x,a->i,a->j,a->a,y); #else if (usecprow){ m = cprow.nrows; ii = cprow.i; ridx = cprow.rindex; } else { ii = a->i; } for (i=0; ij + ii[i] ; v = a->a + ii[i] ; n = ii[i+1] - ii[i]; if (usecprow){ alpha = x[ridx[i]]; } else { alpha = x[i]; } while (n-->0) {y[*idx++] += alpha * *v++;} } #endif ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr); ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatMultTranspose_SeqAIJ" PetscErrorCode MatMultTranspose_SeqAIJ(Mat A,Vec xx,Vec yy) { PetscScalar zero = 0.0; PetscErrorCode ierr; PetscFunctionBegin; ierr = VecSet(yy,zero);CHKERRQ(ierr); ierr = MatMultTransposeAdd_SeqAIJ(A,xx,yy,yy);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatMult_SeqAIJ" PetscErrorCode MatMult_SeqAIJ(Mat A,Vec xx,Vec yy) { Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; PetscScalar *x,*y,*aa; PetscErrorCode ierr; PetscInt m=A->rmap.n,*aj,*ii; PetscInt n,i,j,*ridx=PETSC_NULL; PetscScalar sum; PetscTruth usecprow=a->compressedrow.use; #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTAIJ) PetscInt jrow; #endif #if defined(PETSC_HAVE_PRAGMA_DISJOINT) #pragma disjoint(*x,*y,*aa) #endif PetscFunctionBegin; ierr = VecGetArray(xx,&x);CHKERRQ(ierr); ierr = VecGetArray(yy,&y);CHKERRQ(ierr); aj = a->j; aa = a->a; ii = a->i; if (usecprow){ /* use compressed row format */ m = a->compressedrow.nrows; ii = a->compressedrow.i; ridx = a->compressedrow.rindex; for (i=0; ij + ii[i]; aa = a->a + ii[i]; sum = 0.0; for (j=0; jnz - m);CHKERRQ(ierr); ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatMultAdd_SeqAIJ" PetscErrorCode MatMultAdd_SeqAIJ(Mat A,Vec xx,Vec yy,Vec zz) { Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; PetscScalar *x,*y,*z,*aa; PetscErrorCode ierr; PetscInt m = A->rmap.n,*aj,*ii; #if !defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ) PetscInt n,i,jrow,j,*ridx=PETSC_NULL; PetscScalar sum; PetscTruth usecprow=a->compressedrow.use; #endif PetscFunctionBegin; ierr = VecGetArray(xx,&x);CHKERRQ(ierr); ierr = VecGetArray(yy,&y);CHKERRQ(ierr); if (zz != yy) { ierr = VecGetArray(zz,&z);CHKERRQ(ierr); } else { z = y; } aj = a->j; aa = a->a; ii = a->i; #if defined(PETSC_USE_FORTRAN_KERNEL_MULTADDAIJ) fortranmultaddaij_(&m,x,ii,aj,aa,y,z); #else if (usecprow){ /* use compressed row format */ if (zz != yy){ ierr = PetscMemcpy(z,y,m*sizeof(PetscScalar));CHKERRQ(ierr); } m = a->compressedrow.nrows; ii = a->compressedrow.i; ridx = a->compressedrow.rindex; for (i=0; ij + ii[i]; aa = a->a + ii[i]; sum = y[*ridx]; for (j=0; jnz);CHKERRQ(ierr); ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); if (zz != yy) { ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); } PetscFunctionReturn(0); } /* Adds diagonal pointers to sparse matrix structure. */ #undef __FUNCT__ #define __FUNCT__ "MatMarkDiagonal_SeqAIJ" PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat A) { Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; PetscErrorCode ierr; PetscInt i,j,m = A->rmap.n; PetscFunctionBegin; if (!a->diag) { ierr = PetscMalloc(m*sizeof(PetscInt),&a->diag);CHKERRQ(ierr); } for (i=0; irmap.n; i++) { a->diag[i] = a->i[i+1]; for (j=a->i[i]; ji[i+1]; j++) { if (a->j[j] == i) { a->diag[i] = j; break; } } } PetscFunctionReturn(0); } /* Checks for missing diagonals */ #undef __FUNCT__ #define __FUNCT__ "MatMissingDiagonal_SeqAIJ" PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat A,PetscTruth *missing,PetscInt *d) { Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; PetscInt *diag,*jj = a->j,i; PetscFunctionBegin; *missing = PETSC_FALSE; if (A->rmap.n > 0 && !jj) { *missing = PETSC_TRUE; if (d) *d = 0; PetscInfo(A,"Matrix has no entries therefor is missing diagonal"); } else { diag = a->diag; for (i=0; irmap.n; i++) { if (jj[diag[i]] != i) { *missing = PETSC_TRUE; if (d) *d = i; PetscInfo1(A,"Matrix is missing diagonal number %D",i); } } } PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatRelax_SeqAIJ" PetscErrorCode MatRelax_SeqAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx) { Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; PetscScalar *x,d,*xs,sum,*t,scale,*idiag=0,*mdiag; const PetscScalar *v = a->a, *b, *bs,*xb, *ts; PetscErrorCode ierr; PetscInt n = A->cmap.n,m = A->rmap.n,i; const PetscInt *idx,*diag; PetscFunctionBegin; its = its*lits; if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits); diag = a->diag; if (!a->idiag) { ierr = PetscMalloc(3*m*sizeof(PetscScalar),&a->idiag);CHKERRQ(ierr); a->ssor = a->idiag + m; mdiag = a->ssor + m; v = a->a; /* this is wrong when fshift omega changes each iteration */ if (omega == 1.0 && !fshift) { for (i=0; iidiag[i] = 1.0/v[diag[i]]; } ierr = PetscLogFlops(m);CHKERRQ(ierr); } else { for (i=0; iidiag[i] = omega/(fshift + v[diag[i]]); } ierr = PetscLogFlops(2*m);CHKERRQ(ierr); } } t = a->ssor; idiag = a->idiag; mdiag = a->idiag + 2*m; ierr = VecGetArray(xx,&x);CHKERRQ(ierr); if (xx != bb) { ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); } else { b = x; } /* We count flops by assuming the upper triangular and lower triangular parts have the same number of nonzeros */ xs = x; if (flag == SOR_APPLY_UPPER) { /* apply (U + D/omega) to the vector */ bs = b; for (i=0; ia[diag[i]]; n = a->i[i+1] - diag[i] - 1; idx = a->j + diag[i] + 1; v = a->a + diag[i] + 1; sum = b[i]*d/omega; SPARSEDENSEDOT(sum,bs,v,idx,n); x[i] = sum; } ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); if (bb != xx) {ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);} ierr = PetscLogFlops(a->nz);CHKERRQ(ierr); PetscFunctionReturn(0); } /* Let A = L + U + D; where L is lower trianglar, U is upper triangular, E is diagonal; This routine applies (L + E)^{-1} A (U + E)^{-1} to a vector efficiently using Eisenstat's trick. This is for the case of SSOR preconditioner, so E is D/omega where omega is the relaxation factor. */ if (flag == SOR_APPLY_LOWER) { SETERRQ(PETSC_ERR_SUP,"SOR_APPLY_LOWER is not implemented"); } else if (flag & SOR_EISENSTAT) { /* Let A = L + U + D; where L is lower trianglar, U is upper triangular, E is diagonal; This routine applies (L + E)^{-1} A (U + E)^{-1} to a vector efficiently using Eisenstat's trick. This is for the case of SSOR preconditioner, so E is D/omega where omega is the relaxation factor. */ scale = (2.0/omega) - 1.0; /* x = (E + U)^{-1} b */ for (i=m-1; i>=0; i--) { n = a->i[i+1] - diag[i] - 1; idx = a->j + diag[i] + 1; v = a->a + diag[i] + 1; sum = b[i]; SPARSEDENSEMDOT(sum,xs,v,idx,n); x[i] = sum*idiag[i]; } /* t = b - (2*E - D)x */ v = a->a; for (i=0; idiag; for (i=0; ii[i]; idx = a->j + a->i[i]; v = a->a + a->i[i]; sum = t[i]; SPARSEDENSEMDOT(sum,ts,v,idx,n); t[i] = sum*idiag[i]; /* x = x + t */ x[i] += t[i]; } ierr = PetscLogFlops(6*m-1 + 2*a->nz);CHKERRQ(ierr); ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); if (bb != xx) {ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);} PetscFunctionReturn(0); } if (flag & SOR_ZERO_INITIAL_GUESS) { if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ #if defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ) fortranrelaxaijforwardzero_(&m,&omega,x,a->i,a->j,(PetscInt*)diag,idiag,a->a,(void*)b); #else for (i=0; ii[i]; idx = a->j + a->i[i]; v = a->a + a->i[i]; sum = b[i]; SPARSEDENSEMDOT(sum,xs,v,idx,n); x[i] = sum*idiag[i]; } #endif xb = x; ierr = PetscLogFlops(a->nz);CHKERRQ(ierr); } else xb = b; if ((flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) && (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP)) { for (i=0; ii,a->j,(PetscInt*)diag,idiag,a->a,(void*)xb); #else for (i=m-1; i>=0; i--) { n = a->i[i+1] - diag[i] - 1; idx = a->j + diag[i] + 1; v = a->a + diag[i] + 1; sum = xb[i]; SPARSEDENSEMDOT(sum,xs,v,idx,n); x[i] = sum*idiag[i]; } #endif ierr = PetscLogFlops(a->nz);CHKERRQ(ierr); } its--; } while (its--) { if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ #if defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ) fortranrelaxaijforward_(&m,&omega,x,a->i,a->j,(PetscInt*)diag,a->a,(void*)b); #else for (i=0; ii[i+1] - a->i[i]; idx = a->j + a->i[i]; v = a->a + a->i[i]; sum = b[i]; SPARSEDENSEMDOT(sum,xs,v,idx,n); x[i] = (1. - omega)*x[i] + (sum + mdiag[i]*x[i])*idiag[i]; } #endif ierr = PetscLogFlops(a->nz);CHKERRQ(ierr); } if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ #if defined(PETSC_USE_FORTRAN_KERNEL_RELAXAIJ) fortranrelaxaijbackward_(&m,&omega,x,a->i,a->j,(PetscInt*)diag,a->a,(void*)b); #else for (i=m-1; i>=0; i--) { n = a->i[i+1] - a->i[i]; idx = a->j + a->i[i]; v = a->a + a->i[i]; sum = b[i]; SPARSEDENSEMDOT(sum,xs,v,idx,n); x[i] = (1. - omega)*x[i] + (sum + mdiag[i]*x[i])*idiag[i]; } #endif ierr = PetscLogFlops(a->nz);CHKERRQ(ierr); } } ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); if (bb != xx) {ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr);} PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatGetInfo_SeqAIJ" PetscErrorCode MatGetInfo_SeqAIJ(Mat A,MatInfoType flag,MatInfo *info) { Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; PetscFunctionBegin; info->rows_global = (double)A->rmap.n; info->columns_global = (double)A->cmap.n; info->rows_local = (double)A->rmap.n; info->columns_local = (double)A->cmap.n; info->block_size = 1.0; info->nz_allocated = (double)a->maxnz; info->nz_used = (double)a->nz; info->nz_unneeded = (double)(a->maxnz - a->nz); info->assemblies = (double)A->num_ass; info->mallocs = (double)a->reallocs; info->memory = A->mem; if (A->factor) { info->fill_ratio_given = A->info.fill_ratio_given; info->fill_ratio_needed = A->info.fill_ratio_needed; info->factor_mallocs = A->info.factor_mallocs; } else { info->fill_ratio_given = 0; info->fill_ratio_needed = 0; info->factor_mallocs = 0; } PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatZeroRows_SeqAIJ" PetscErrorCode MatZeroRows_SeqAIJ(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag) { Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; PetscInt i,m = A->rmap.n - 1,d; PetscErrorCode ierr; PetscTruth missing; PetscFunctionBegin; if (a->keepzeroedrows) { for (i=0; i m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %D out of range", rows[i]); ierr = PetscMemzero(&a->a[a->i[rows[i]]],a->ilen[rows[i]]*sizeof(PetscScalar));CHKERRQ(ierr); } if (diag != 0.0) { ierr = MatMissingDiagonal_SeqAIJ(A,&missing,&d);CHKERRQ(ierr); if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry in row %D",d); for (i=0; ia[a->diag[rows[i]]] = diag; } } A->same_nonzero = PETSC_TRUE; } else { if (diag != 0.0) { for (i=0; i m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %D out of range", rows[i]); if (a->ilen[rows[i]] > 0) { a->ilen[rows[i]] = 1; a->a[a->i[rows[i]]] = diag; a->j[a->i[rows[i]]] = rows[i]; } else { /* in case row was completely empty */ ierr = MatSetValues_SeqAIJ(A,1,&rows[i],1,&rows[i],&diag,INSERT_VALUES);CHKERRQ(ierr); } } } else { for (i=0; i m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"row %D out of range", rows[i]); a->ilen[rows[i]] = 0; } } A->same_nonzero = PETSC_FALSE; } ierr = MatAssemblyEnd_SeqAIJ(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatGetRow_SeqAIJ" PetscErrorCode MatGetRow_SeqAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) { Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; PetscInt *itmp; PetscFunctionBegin; if (row < 0 || row >= A->rmap.n) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Row %D out of range",row); *nz = a->i[row+1] - a->i[row]; if (v) *v = a->a + a->i[row]; if (idx) { itmp = a->j + a->i[row]; if (*nz) { *idx = itmp; } else *idx = 0; } PetscFunctionReturn(0); } /* remove this function? */ #undef __FUNCT__ #define __FUNCT__ "MatRestoreRow_SeqAIJ" PetscErrorCode MatRestoreRow_SeqAIJ(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) { PetscFunctionBegin; PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatNorm_SeqAIJ" PetscErrorCode MatNorm_SeqAIJ(Mat A,NormType type,PetscReal *nrm) { Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; PetscScalar *v = a->a; PetscReal sum = 0.0; PetscErrorCode ierr; PetscInt i,j; PetscFunctionBegin; if (type == NORM_FROBENIUS) { for (i=0; inz; i++) { #if defined(PETSC_USE_COMPLEX) sum += PetscRealPart(PetscConj(*v)*(*v)); v++; #else sum += (*v)*(*v); v++; #endif } *nrm = sqrt(sum); } else if (type == NORM_1) { PetscReal *tmp; PetscInt *jj = a->j; ierr = PetscMalloc((A->cmap.n+1)*sizeof(PetscReal),&tmp);CHKERRQ(ierr); ierr = PetscMemzero(tmp,A->cmap.n*sizeof(PetscReal));CHKERRQ(ierr); *nrm = 0.0; for (j=0; jnz; j++) { tmp[*jj++] += PetscAbsScalar(*v); v++; } for (j=0; jcmap.n; j++) { if (tmp[j] > *nrm) *nrm = tmp[j]; } ierr = PetscFree(tmp);CHKERRQ(ierr); } else if (type == NORM_INFINITY) { *nrm = 0.0; for (j=0; jrmap.n; j++) { v = a->a + a->i[j]; sum = 0.0; for (i=0; ii[j+1]-a->i[j]; i++) { sum += PetscAbsScalar(*v); v++; } if (sum > *nrm) *nrm = sum; } } else { SETERRQ(PETSC_ERR_SUP,"No support for two norm"); } PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatTranspose_SeqAIJ" PetscErrorCode MatTranspose_SeqAIJ(Mat A,Mat *B) { Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; Mat C; PetscErrorCode ierr; PetscInt i,*aj = a->j,*ai = a->i,m = A->rmap.n,len,*col; PetscScalar *array = a->a; PetscFunctionBegin; if (!B && m != A->cmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Square matrix only for in-place"); ierr = PetscMalloc((1+A->cmap.n)*sizeof(PetscInt),&col);CHKERRQ(ierr); ierr = PetscMemzero(col,(1+A->cmap.n)*sizeof(PetscInt));CHKERRQ(ierr); for (i=0; icomm,&C);CHKERRQ(ierr); ierr = MatSetSizes(C,A->cmap.n,m,A->cmap.n,m);CHKERRQ(ierr); ierr = MatSetType(C,A->type_name);CHKERRQ(ierr); ierr = MatSeqAIJSetPreallocation_SeqAIJ(C,0,col);CHKERRQ(ierr); ierr = PetscFree(col);CHKERRQ(ierr); for (i=0; idata,*bij = (Mat_SeqAIJ*) A->data; PetscInt *adx,*bdx,*aii,*bii,*aptr,*bptr; PetscScalar *va,*vb; PetscErrorCode ierr; PetscInt ma,na,mb,nb, i; PetscFunctionBegin; bij = (Mat_SeqAIJ *) B->data; ierr = MatGetSize(A,&ma,&na);CHKERRQ(ierr); ierr = MatGetSize(B,&mb,&nb);CHKERRQ(ierr); if (ma!=nb || na!=mb){ *f = PETSC_FALSE; PetscFunctionReturn(0); } aii = aij->i; bii = bij->i; adx = aij->j; bdx = bij->j; va = aij->a; vb = bij->a; ierr = PetscMalloc(ma*sizeof(PetscInt),&aptr);CHKERRQ(ierr); ierr = PetscMalloc(mb*sizeof(PetscInt),&bptr);CHKERRQ(ierr); for (i=0; i tol) { *f = PETSC_FALSE; goto done; } else { aptr[i]++; if (B || i!=idc) bptr[idc]++; } } } done: ierr = PetscFree(aptr);CHKERRQ(ierr); if (B) { ierr = PetscFree(bptr);CHKERRQ(ierr); } PetscFunctionReturn(0); } EXTERN_C_END #undef __FUNCT__ #define __FUNCT__ "MatIsSymmetric_SeqAIJ" PetscErrorCode MatIsSymmetric_SeqAIJ(Mat A,PetscReal tol,PetscTruth *f) { PetscErrorCode ierr; PetscFunctionBegin; ierr = MatIsTranspose_SeqAIJ(A,A,tol,f);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatDiagonalScale_SeqAIJ" PetscErrorCode MatDiagonalScale_SeqAIJ(Mat A,Vec ll,Vec rr) { Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; PetscScalar *l,*r,x,*v; PetscErrorCode ierr; PetscInt i,j,m = A->rmap.n,n = A->cmap.n,M,nz = a->nz,*jj; PetscFunctionBegin; if (ll) { /* The local size is used so that VecMPI can be passed to this routine by MatDiagonalScale_MPIAIJ */ ierr = VecGetLocalSize(ll,&m);CHKERRQ(ierr); if (m != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector wrong length"); ierr = VecGetArray(ll,&l);CHKERRQ(ierr); v = a->a; for (i=0; ii[i+1] - a->i[i]; for (j=0; jcmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vector wrong length"); ierr = VecGetArray(rr,&r);CHKERRQ(ierr); v = a->a; jj = a->j; for (i=0; idata,*c; PetscErrorCode ierr; PetscInt *smap,i,k,kstart,kend,oldcols = A->cmap.n,*lens; PetscInt row,mat_i,*mat_j,tcol,first,step,*mat_ilen,sum,lensi; PetscInt *irow,*icol,nrows,ncols; PetscInt *starts,*j_new,*i_new,*aj = a->j,*ai = a->i,ii,*ailen = a->ilen; PetscScalar *a_new,*mat_a; Mat C; PetscTruth stride; PetscFunctionBegin; ierr = ISSorted(isrow,(PetscTruth*)&i);CHKERRQ(ierr); if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"ISrow is not sorted"); ierr = ISSorted(iscol,(PetscTruth*)&i);CHKERRQ(ierr); if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"IScol is not sorted"); ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); ierr = ISStrideGetInfo(iscol,&first,&step);CHKERRQ(ierr); ierr = ISStride(iscol,&stride);CHKERRQ(ierr); if (stride && step == 1) { /* special case of contiguous rows */ ierr = PetscMalloc((2*nrows+1)*sizeof(PetscInt),&lens);CHKERRQ(ierr); starts = lens + nrows; /* loop over new rows determining lens and starting points */ for (i=0; i= first) { starts[i] = k; break; } } sum = 0; while (k < kend) { if (aj[k++] >= first+ncols) break; sum++; } lens[i] = sum; } /* create submatrix */ if (scall == MAT_REUSE_MATRIX) { PetscInt n_cols,n_rows; ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr); if (n_rows != nrows || n_cols != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Reused submatrix wrong size"); ierr = MatZeroEntries(*B);CHKERRQ(ierr); C = *B; } else { ierr = MatCreate(A->comm,&C);CHKERRQ(ierr); ierr = MatSetSizes(C,nrows,ncols,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); ierr = MatSetType(C,A->type_name);CHKERRQ(ierr); ierr = MatSeqAIJSetPreallocation_SeqAIJ(C,0,lens);CHKERRQ(ierr); } c = (Mat_SeqAIJ*)C->data; /* loop over rows inserting into submatrix */ a_new = c->a; j_new = c->j; i_new = c->i; for (i=0; ia + starts[i],lensi*sizeof(PetscScalar));CHKERRQ(ierr); a_new += lensi; i_new[i+1] = i_new[i] + lensi; c->ilen[i] = lensi; } ierr = PetscFree(lens);CHKERRQ(ierr); } else { ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); ierr = PetscMalloc((1+oldcols)*sizeof(PetscInt),&smap);CHKERRQ(ierr); ierr = PetscMalloc((1+nrows)*sizeof(PetscInt),&lens);CHKERRQ(ierr); ierr = PetscMemzero(smap,oldcols*sizeof(PetscInt));CHKERRQ(ierr); for (i=0; iilen[irow[i]]; lens[i] = 0; for (k=kstart; kdata); if ((*B)->rmap.n != nrows || (*B)->cmap.n != ncols) SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong size"); ierr = PetscMemcmp(c->ilen,lens,(*B)->rmap.n*sizeof(PetscInt),&equal);CHKERRQ(ierr); if (!equal) { SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros"); } ierr = PetscMemzero(c->ilen,(*B)->rmap.n*sizeof(PetscInt));CHKERRQ(ierr); C = *B; } else { ierr = MatCreate(A->comm,&C);CHKERRQ(ierr); ierr = MatSetSizes(C,nrows,ncols,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); ierr = MatSetType(C,A->type_name);CHKERRQ(ierr); ierr = MatSeqAIJSetPreallocation_SeqAIJ(C,0,lens);CHKERRQ(ierr); } c = (Mat_SeqAIJ *)(C->data); for (i=0; iilen[row]; mat_i = c->i[i]; mat_j = c->j + mat_i; mat_a = c->a + mat_i; mat_ilen = c->ilen + i; for (k=kstart; kj[k]])) { *mat_j++ = tcol - 1; *mat_a++ = a->a[k]; (*mat_ilen)++; } } } /* Free work space */ ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); ierr = PetscFree(smap);CHKERRQ(ierr); ierr = PetscFree(lens);CHKERRQ(ierr); } ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); *B = C; PetscFunctionReturn(0); } /* */ #undef __FUNCT__ #define __FUNCT__ "MatILUFactor_SeqAIJ" PetscErrorCode MatILUFactor_SeqAIJ(Mat inA,IS row,IS col,MatFactorInfo *info) { Mat_SeqAIJ *a = (Mat_SeqAIJ*)inA->data; PetscErrorCode ierr; Mat outA; PetscTruth row_identity,col_identity; PetscFunctionBegin; if (info->levels != 0) SETERRQ(PETSC_ERR_SUP,"Only levels=0 supported for in-place ilu"); ierr = ISIdentity(row,&row_identity);CHKERRQ(ierr); ierr = ISIdentity(col,&col_identity);CHKERRQ(ierr); outA = inA; inA->factor = FACTOR_LU; a->row = row; a->col = col; ierr = PetscObjectReference((PetscObject)row);CHKERRQ(ierr); ierr = PetscObjectReference((PetscObject)col);CHKERRQ(ierr); /* Create the inverse permutation so that it can be used in MatLUFactorNumeric() */ if (a->icol) {ierr = ISDestroy(a->icol);CHKERRQ(ierr);} /* need to remove old one */ ierr = ISInvertPermutation(col,PETSC_DECIDE,&a->icol);CHKERRQ(ierr); ierr = PetscLogObjectParent(inA,a->icol);CHKERRQ(ierr); if (!a->solve_work) { /* this matrix may have been factored before */ ierr = PetscMalloc((inA->rmap.n+1)*sizeof(PetscScalar),&a->solve_work);CHKERRQ(ierr); } ierr = MatMarkDiagonal_SeqAIJ(inA);CHKERRQ(ierr); if (row_identity && col_identity) { ierr = MatLUFactorNumeric_SeqAIJ(inA,info,&outA);CHKERRQ(ierr); } else { ierr = MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(inA,info,&outA);CHKERRQ(ierr); } PetscFunctionReturn(0); } #include "petscblaslapack.h" #undef __FUNCT__ #define __FUNCT__ "MatScale_SeqAIJ" PetscErrorCode MatScale_SeqAIJ(Mat inA,PetscScalar alpha) { Mat_SeqAIJ *a = (Mat_SeqAIJ*)inA->data; PetscBLASInt bnz = (PetscBLASInt)a->nz,one = 1; PetscScalar oalpha = alpha; PetscErrorCode ierr; PetscFunctionBegin; BLASscal_(&bnz,&oalpha,a->a,&one); ierr = PetscLogFlops(a->nz);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatGetSubMatrices_SeqAIJ" PetscErrorCode MatGetSubMatrices_SeqAIJ(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[]) { PetscErrorCode ierr; PetscInt i; PetscFunctionBegin; if (scall == MAT_INITIAL_MATRIX) { ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr); } for (i=0; idata; PetscErrorCode ierr; PetscInt row,i,j,k,l,m,n,*idx,*nidx,isz,val; PetscInt start,end,*ai,*aj; PetscBT table; PetscFunctionBegin; m = A->rmap.n; ai = a->i; aj = a->j; if (ov < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"illegal negative overlap value used"); ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nidx);CHKERRQ(ierr); ierr = PetscBTCreate(m,table);CHKERRQ(ierr); for (i=0; idata; PetscErrorCode ierr; PetscInt i,nz,m = A->rmap.n,n = A->cmap.n,*col; PetscInt *row,*cnew,j,*lens; IS icolp,irowp; PetscInt *cwork; PetscScalar *vwork; PetscFunctionBegin; ierr = ISInvertPermutation(rowp,PETSC_DECIDE,&irowp);CHKERRQ(ierr); ierr = ISGetIndices(irowp,&row);CHKERRQ(ierr); ierr = ISInvertPermutation(colp,PETSC_DECIDE,&icolp);CHKERRQ(ierr); ierr = ISGetIndices(icolp,&col);CHKERRQ(ierr); /* determine lengths of permuted rows */ ierr = PetscMalloc((m+1)*sizeof(PetscInt),&lens);CHKERRQ(ierr); for (i=0; ii[i+1] - a->i[i]; } ierr = MatCreate(A->comm,B);CHKERRQ(ierr); ierr = MatSetSizes(*B,m,n,m,n);CHKERRQ(ierr); ierr = MatSetType(*B,A->type_name);CHKERRQ(ierr); ierr = MatSeqAIJSetPreallocation_SeqAIJ(*B,0,lens);CHKERRQ(ierr); ierr = PetscFree(lens);CHKERRQ(ierr); ierr = PetscMalloc(n*sizeof(PetscInt),&cnew);CHKERRQ(ierr); for (i=0; iassembled = PETSC_FALSE; ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); ierr = ISRestoreIndices(irowp,&row);CHKERRQ(ierr); ierr = ISRestoreIndices(icolp,&col);CHKERRQ(ierr); ierr = ISDestroy(irowp);CHKERRQ(ierr); ierr = ISDestroy(icolp);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatCopy_SeqAIJ" PetscErrorCode MatCopy_SeqAIJ(Mat A,Mat B,MatStructure str) { PetscErrorCode ierr; PetscFunctionBegin; /* If the two matrices have the same copy implementation, use fast copy. */ if (str == SAME_NONZERO_PATTERN && (A->ops->copy == B->ops->copy)) { Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data; if (a->i[A->rmap.n] != b->i[B->rmap.n]) { SETERRQ(PETSC_ERR_ARG_INCOMP,"Number of nonzeros in two matrices are different"); } ierr = PetscMemcpy(b->a,a->a,(a->i[A->rmap.n])*sizeof(PetscScalar));CHKERRQ(ierr); } else { ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr); } PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatSetUpPreallocation_SeqAIJ" PetscErrorCode MatSetUpPreallocation_SeqAIJ(Mat A) { PetscErrorCode ierr; PetscFunctionBegin; ierr = MatSeqAIJSetPreallocation_SeqAIJ(A,PETSC_DEFAULT,0);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatGetArray_SeqAIJ" PetscErrorCode MatGetArray_SeqAIJ(Mat A,PetscScalar *array[]) { Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; PetscFunctionBegin; *array = a->a; PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatRestoreArray_SeqAIJ" PetscErrorCode MatRestoreArray_SeqAIJ(Mat A,PetscScalar *array[]) { PetscFunctionBegin; PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatFDColoringApply_SeqAIJ" PetscErrorCode MatFDColoringApply_SeqAIJ(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx) { PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void *))coloring->f; PetscErrorCode ierr; PetscInt k,N,start,end,l,row,col,srow,**vscaleforrow,m1,m2; PetscScalar dx,*y,*xx,*w3_array; PetscScalar *vscale_array; PetscReal epsilon = coloring->error_rel,umin = coloring->umin; Vec w1,w2,w3; void *fctx = coloring->fctx; PetscTruth flg; PetscFunctionBegin; if (!coloring->w1) { ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr); ierr = PetscLogObjectParent(coloring,coloring->w1);CHKERRQ(ierr); ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr); ierr = PetscLogObjectParent(coloring,coloring->w2);CHKERRQ(ierr); ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr); ierr = PetscLogObjectParent(coloring,coloring->w3);CHKERRQ(ierr); } w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3; ierr = MatSetUnfactored(J);CHKERRQ(ierr); ierr = PetscOptionsHasName(coloring->prefix,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr); if (flg) { ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr); } else { PetscTruth assembled; ierr = MatAssembled(J,&assembled);CHKERRQ(ierr); if (assembled) { ierr = MatZeroEntries(J);CHKERRQ(ierr); } } ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr); ierr = VecGetSize(x1,&N);CHKERRQ(ierr); /* This is a horrible, horrible, hack. See DMMGComputeJacobian_Multigrid() it inproperly sets coloring->F for the coarser grids from the finest */ if (coloring->F) { ierr = VecGetLocalSize(coloring->F,&m1);CHKERRQ(ierr); ierr = VecGetLocalSize(w1,&m2);CHKERRQ(ierr); if (m1 != m2) { coloring->F = 0; } } if (coloring->F) { w1 = coloring->F; coloring->F = 0; } else { ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr); ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); } /* Compute all the scale factors and share with other processors */ ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start; ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start; for (k=0; kncolors; k++) { /* Loop over each column associated with color adding the perturbation to the vector w3. */ for (l=0; lncolumns[k]; l++) { col = coloring->columns[k][l]; /* column of the matrix we are probing for */ dx = xx[col]; if (dx == 0.0) dx = 1.0; #if !defined(PETSC_USE_COMPLEX) if (dx < umin && dx >= 0.0) dx = umin; else if (dx < 0.0 && dx > -umin) dx = -umin; #else if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; #endif dx *= epsilon; vscale_array[col] = 1.0/dx; } } vscale_array = vscale_array + start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); /* ierr = VecView(coloring->vscale,PETSC_VIEWER_STDOUT_WORLD); ierr = VecView(x1,PETSC_VIEWER_STDOUT_WORLD);*/ if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow; else vscaleforrow = coloring->columnsforrow; ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); /* Loop over each color */ for (k=0; kncolors; k++) { coloring->currentcolor = k; ierr = VecCopy(x1,w3);CHKERRQ(ierr); ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start; /* Loop over each column associated with color adding the perturbation to the vector w3. */ for (l=0; lncolumns[k]; l++) { col = coloring->columns[k][l]; /* column of the matrix we are probing for */ dx = xx[col]; if (dx == 0.0) dx = 1.0; #if !defined(PETSC_USE_COMPLEX) if (dx < umin && dx >= 0.0) dx = umin; else if (dx < 0.0 && dx > -umin) dx = -umin; #else if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0) dx = umin; else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin; #endif dx *= epsilon; if (!PetscAbsScalar(dx)) SETERRQ(PETSC_ERR_PLIB,"Computed 0 differencing parameter"); w3_array[col] += dx; } w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr); /* Evaluate function at x1 + dx (here dx is a vector of perturbations) */ ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr); ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr); /* Loop over rows of vector, putting results into Jacobian matrix */ ierr = VecGetArray(w2,&y);CHKERRQ(ierr); for (l=0; lnrows[k]; l++) { row = coloring->rows[k][l]; col = coloring->columnsforrow[k][l]; y[row] *= vscale_array[vscaleforrow[k][l]]; srow = row + start; ierr = MatSetValues_SeqAIJ(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr); } ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr); } coloring->currentcolor = k; ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr); xx = xx + start; ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr); ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); PetscFunctionReturn(0); } #include "petscblaslapack.h" #undef __FUNCT__ #define __FUNCT__ "MatAXPY_SeqAIJ" PetscErrorCode MatAXPY_SeqAIJ(Mat Y,PetscScalar a,Mat X,MatStructure str) { PetscErrorCode ierr; PetscInt i; Mat_SeqAIJ *x = (Mat_SeqAIJ *)X->data,*y = (Mat_SeqAIJ *)Y->data; PetscBLASInt one=1,bnz = (PetscBLASInt)x->nz; PetscFunctionBegin; if (str == SAME_NONZERO_PATTERN) { PetscScalar alpha = a; BLASaxpy_(&bnz,&alpha,x->a,&one,y->a,&one); } else if (str == SUBSET_NONZERO_PATTERN) { /* nonzeros of X is a subset of Y's */ if (y->xtoy && y->XtoY != X) { ierr = PetscFree(y->xtoy);CHKERRQ(ierr); ierr = MatDestroy(y->XtoY);CHKERRQ(ierr); } if (!y->xtoy) { /* get xtoy */ ierr = MatAXPYGetxtoy_Private(X->rmap.n,x->i,x->j,PETSC_NULL, y->i,y->j,PETSC_NULL, &y->xtoy);CHKERRQ(ierr); y->XtoY = X; } for (i=0; inz; i++) y->a[y->xtoy[i]] += a*(x->a[i]); ierr = PetscInfo3(0,"ratio of nnz(X)/nnz(Y): %d/%d = %G\n",x->nz,y->nz,(PetscReal)(x->nz)/y->nz);CHKERRQ(ierr); } else { ierr = MatAXPY_Basic(Y,a,X,str);CHKERRQ(ierr); } PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatSetBlockSize_SeqAIJ" PetscErrorCode MatSetBlockSize_SeqAIJ(Mat A,PetscInt bs) { PetscFunctionBegin; PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatConjugate_SeqAIJ" PetscErrorCode PETSCMAT_DLLEXPORT MatConjugate_SeqAIJ(Mat mat) { #if defined(PETSC_USE_COMPLEX) Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data; PetscInt i,nz; PetscScalar *a; PetscFunctionBegin; nz = aij->nz; a = aij->a; for (i=0; idata; PetscErrorCode ierr; PetscInt i,j,m = A->rmap.n,*ai,*aj,ncols,n; PetscReal atmp; PetscScalar *x,zero = 0.0; MatScalar *aa; PetscFunctionBegin; if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); aa = a->a; ai = a->i; aj = a->j; ierr = VecSet(v,zero);CHKERRQ(ierr); ierr = VecGetArray(v,&x);CHKERRQ(ierr); ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr); if (n != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); for (i=0; idata; PetscInt i,nz,n; PetscFunctionBegin; nz = aij->maxnz; n = mat->cmap.n; for (i=0; ij[i] = indices[i]; } aij->nz = nz; for (i=0; iilen[i] = aij->imax[i]; } PetscFunctionReturn(0); } EXTERN_C_END #undef __FUNCT__ #define __FUNCT__ "MatSeqAIJSetColumnIndices" /*@ MatSeqAIJSetColumnIndices - Set the column indices for all the rows in the matrix. Input Parameters: + mat - the SeqAIJ matrix - indices - the column indices Level: advanced Notes: This can be called if you have precomputed the nonzero structure of the matrix and want to provide it to the matrix object to improve the performance of the MatSetValues() operation. You MUST have set the correct numbers of nonzeros per row in the call to MatCreateSeqAIJ(), and the columns indices MUST be sorted. MUST be called before any calls to MatSetValues(); The indices should start with zero, not one. @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetColumnIndices(Mat mat,PetscInt *indices) { PetscErrorCode ierr,(*f)(Mat,PetscInt *); PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); PetscValidPointer(indices,2); ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSeqAIJSetColumnIndices_C",(void (**)(void))&f);CHKERRQ(ierr); if (f) { ierr = (*f)(mat,indices);CHKERRQ(ierr); } else { SETERRQ(PETSC_ERR_SUP,"Wrong type of matrix to set column indices"); } PetscFunctionReturn(0); } /* ----------------------------------------------------------------------------------------*/ EXTERN_C_BEGIN #undef __FUNCT__ #define __FUNCT__ "MatStoreValues_SeqAIJ" PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues_SeqAIJ(Mat mat) { Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data; PetscErrorCode ierr; size_t nz = aij->i[mat->rmap.n]; PetscFunctionBegin; if (aij->nonew != 1) { SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); } /* allocate space for values if not already there */ if (!aij->saved_values) { ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&aij->saved_values);CHKERRQ(ierr); } /* copy values over */ ierr = PetscMemcpy(aij->saved_values,aij->a,nz*sizeof(PetscScalar));CHKERRQ(ierr); PetscFunctionReturn(0); } EXTERN_C_END #undef __FUNCT__ #define __FUNCT__ "MatStoreValues" /*@ MatStoreValues - Stashes a copy of the matrix values; this allows, for example, reuse of the linear part of a Jacobian, while recomputing the nonlinear portion. Collect on Mat Input Parameters: . mat - the matrix (currently only AIJ matrices support this option) Level: advanced Common Usage, with SNESSolve(): $ Create Jacobian matrix $ Set linear terms into matrix $ Apply boundary conditions to matrix, at this time matrix must have $ final nonzero structure (i.e. setting the nonlinear terms and applying $ boundary conditions again will not change the nonzero structure $ ierr = MatSetOption(mat,MAT_NO_NEW_NONZERO_LOCATIONS); $ ierr = MatStoreValues(mat); $ Call SNESSetJacobian() with matrix $ In your Jacobian routine $ ierr = MatRetrieveValues(mat); $ Set nonlinear terms in matrix Common Usage without SNESSolve(), i.e. when you handle nonlinear solve yourself: $ // build linear portion of Jacobian $ ierr = MatSetOption(mat,MAT_NO_NEW_NONZERO_LOCATIONS); $ ierr = MatStoreValues(mat); $ loop over nonlinear iterations $ ierr = MatRetrieveValues(mat); $ // call MatSetValues(mat,...) to set nonliner portion of Jacobian $ // call MatAssemblyBegin/End() on matrix $ Solve linear system with Jacobian $ endloop Notes: Matrix must already be assemblied before calling this routine Must set the matrix option MatSetOption(mat,MAT_NO_NEW_NONZERO_LOCATIONS); before calling this routine. When this is called multiple times it overwrites the previous set of stored values and does not allocated additional space. .seealso: MatRetrieveValues() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues(Mat mat) { PetscErrorCode ierr,(*f)(Mat); PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); ierr = PetscObjectQueryFunction((PetscObject)mat,"MatStoreValues_C",(void (**)(void))&f);CHKERRQ(ierr); if (f) { ierr = (*f)(mat);CHKERRQ(ierr); } else { SETERRQ(PETSC_ERR_SUP,"Wrong type of matrix to store values"); } PetscFunctionReturn(0); } EXTERN_C_BEGIN #undef __FUNCT__ #define __FUNCT__ "MatRetrieveValues_SeqAIJ" PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues_SeqAIJ(Mat mat) { Mat_SeqAIJ *aij = (Mat_SeqAIJ *)mat->data; PetscErrorCode ierr; PetscInt nz = aij->i[mat->rmap.n]; PetscFunctionBegin; if (aij->nonew != 1) { SETERRQ(PETSC_ERR_ORDER,"Must call MatSetOption(A,MAT_NO_NEW_NONZERO_LOCATIONS);first"); } if (!aij->saved_values) { SETERRQ(PETSC_ERR_ORDER,"Must call MatStoreValues(A);first"); } /* copy values over */ ierr = PetscMemcpy(aij->a,aij->saved_values,nz*sizeof(PetscScalar));CHKERRQ(ierr); PetscFunctionReturn(0); } EXTERN_C_END #undef __FUNCT__ #define __FUNCT__ "MatRetrieveValues" /*@ MatRetrieveValues - Retrieves the copy of the matrix values; this allows, for example, reuse of the linear part of a Jacobian, while recomputing the nonlinear portion. Collect on Mat Input Parameters: . mat - the matrix (currently on AIJ matrices support this option) Level: advanced .seealso: MatStoreValues() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues(Mat mat) { PetscErrorCode ierr,(*f)(Mat); PetscFunctionBegin; PetscValidHeaderSpecific(mat,MAT_COOKIE,1); if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); ierr = PetscObjectQueryFunction((PetscObject)mat,"MatRetrieveValues_C",(void (**)(void))&f);CHKERRQ(ierr); if (f) { ierr = (*f)(mat);CHKERRQ(ierr); } else { SETERRQ(PETSC_ERR_SUP,"Wrong type of matrix to retrieve values"); } PetscFunctionReturn(0); } /* --------------------------------------------------------------------------------*/ #undef __FUNCT__ #define __FUNCT__ "MatCreateSeqAIJ" /*@C MatCreateSeqAIJ - Creates a sparse matrix in AIJ (compressed row) format (the default parallel PETSc format). For good matrix assembly performance the user should preallocate the matrix storage by setting the parameter nz (or the array nnz). By setting these parameters accurately, performance during matrix assembly can be increased by more than a factor of 50. Collective on MPI_Comm Input Parameters: + comm - MPI communicator, set to PETSC_COMM_SELF . m - number of rows . n - number of columns . nz - number of nonzeros per row (same for all rows) - nnz - array containing the number of nonzeros in the various rows (possibly different for each row) or PETSC_NULL Output Parameter: . A - the matrix Notes: If nnz is given then nz is ignored The AIJ format (also called the Yale sparse matrix format or compressed row storage), is fully compatible with standard Fortran 77 storage. That is, the stored row and column indices can begin at either one (as in Fortran) or zero. See the users' manual for details. Specify the preallocated storage with either nz or nnz (not both). Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory allocation. For large problems you MUST preallocate memory or you will get TERRIBLE performance, see the users' manual chapter on matrices. By default, this format uses inodes (identical nodes) when possible, to improve numerical efficiency of matrix-vector products and solves. We search for consecutive rows with the same nonzero structure, thereby reusing matrix information to achieve increased efficiency. Options Database Keys: + -mat_no_inode - Do not use inodes . -mat_inode_limit - Sets inode limit (max limit=5) - -mat_aij_oneindex - Internally use indexing starting at 1 rather than 0. Note that when calling MatSetValues(), the user still MUST index entries starting at 0! Level: intermediate .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues(), MatSeqAIJSetColumnIndices(), MatCreateSeqAIJWithArrays() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJ(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A) { PetscErrorCode ierr; PetscFunctionBegin; ierr = MatCreate(comm,A);CHKERRQ(ierr); ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr); ierr = MatSetType(*A,MATSEQAIJ);CHKERRQ(ierr); ierr = MatSeqAIJSetPreallocation_SeqAIJ(*A,nz,(PetscInt*)nnz);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatSeqAIJSetPreallocation" /*@C MatSeqAIJSetPreallocation - For good matrix assembly performance the user should preallocate the matrix storage by setting the parameter nz (or the array nnz). By setting these parameters accurately, performance during matrix assembly can be increased by more than a factor of 50. Collective on MPI_Comm Input Parameters: + B - The matrix . nz - number of nonzeros per row (same for all rows) - nnz - array containing the number of nonzeros in the various rows (possibly different for each row) or PETSC_NULL Notes: If nnz is given then nz is ignored The AIJ format (also called the Yale sparse matrix format or compressed row storage), is fully compatible with standard Fortran 77 storage. That is, the stored row and column indices can begin at either one (as in Fortran) or zero. See the users' manual for details. Specify the preallocated storage with either nz or nnz (not both). Set nz=PETSC_DEFAULT and nnz=PETSC_NULL for PETSc to control dynamic memory allocation. For large problems you MUST preallocate memory or you will get TERRIBLE performance, see the users' manual chapter on matrices. Developers: Use nz of MAT_SKIP_ALLOCATION to not allocate any space for the matrix entries or columns indices By default, this format uses inodes (identical nodes) when possible, to improve numerical efficiency of matrix-vector products and solves. We search for consecutive rows with the same nonzero structure, thereby reusing matrix information to achieve increased efficiency. Options Database Keys: + -mat_no_inode - Do not use inodes . -mat_inode_limit - Sets inode limit (max limit=5) - -mat_aij_oneindex - Internally use indexing starting at 1 rather than 0. Note that when calling MatSetValues(), the user still MUST index entries starting at 0! Level: intermediate .seealso: MatCreate(), MatCreateMPIAIJ(), MatSetValues(), MatSeqAIJSetColumnIndices(), MatCreateSeqAIJWithArrays() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocation(Mat B,PetscInt nz,const PetscInt nnz[]) { PetscErrorCode ierr,(*f)(Mat,PetscInt,const PetscInt[]); PetscFunctionBegin; ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqAIJSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); if (f) { ierr = (*f)(B,nz,nnz);CHKERRQ(ierr); } PetscFunctionReturn(0); } EXTERN_C_BEGIN #undef __FUNCT__ #define __FUNCT__ "MatSeqAIJSetPreallocation_SeqAIJ" PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocation_SeqAIJ(Mat B,PetscInt nz,PetscInt *nnz) { Mat_SeqAIJ *b; PetscTruth skipallocation = PETSC_FALSE; PetscErrorCode ierr; PetscInt i; PetscFunctionBegin; if (nz == MAT_SKIP_ALLOCATION) { skipallocation = PETSC_TRUE; nz = 0; } B->rmap.bs = B->cmap.bs = 1; ierr = PetscMapInitialize(B->comm,&B->rmap);CHKERRQ(ierr); ierr = PetscMapInitialize(B->comm,&B->cmap);CHKERRQ(ierr); if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 5; if (nz < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"nz cannot be less than 0: value %d",nz); if (nnz) { for (i=0; irmap.n; i++) { if (nnz[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be less than 0: local row %d value %d",i,nnz[i]); if (nnz[i] > B->cmap.n) SETERRQ3(PETSC_ERR_ARG_OUTOFRANGE,"nnz cannot be greater than row length: local row %d value %d rowlength %d",i,nnz[i],B->cmap.n); } } B->preallocated = PETSC_TRUE; b = (Mat_SeqAIJ*)B->data; if (!skipallocation) { ierr = PetscMalloc2(B->rmap.n,PetscInt,&b->imax,B->rmap.n,PetscInt,&b->ilen);CHKERRQ(ierr); if (!nnz) { if (nz == PETSC_DEFAULT || nz == PETSC_DECIDE) nz = 10; else if (nz <= 0) nz = 1; for (i=0; irmap.n; i++) b->imax[i] = nz; nz = nz*B->rmap.n; } else { nz = 0; for (i=0; irmap.n; i++) {b->imax[i] = nnz[i]; nz += nnz[i];} } /* b->ilen will count nonzeros in each row so far. */ for (i=0; irmap.n; i++) { b->ilen[i] = 0;} /* allocate the matrix space */ ierr = PetscMalloc3(nz,PetscScalar,&b->a,nz,PetscInt,&b->j,B->rmap.n+1,PetscInt,&b->i);CHKERRQ(ierr); b->i[0] = 0; for (i=1; irmap.n+1; i++) { b->i[i] = b->i[i-1] + b->imax[i-1]; } b->singlemalloc = PETSC_TRUE; b->free_a = PETSC_TRUE; b->free_ij = PETSC_TRUE; } else { b->free_a = PETSC_FALSE; b->free_ij = PETSC_FALSE; } b->nz = 0; b->maxnz = nz; B->info.nz_unneeded = (double)b->maxnz; PetscFunctionReturn(0); } EXTERN_C_END #undef __FUNCT__ #define __FUNCT__ "MatSeqAIJSetPreallocationCSR" /*@C MatSeqAIJSetPreallocationCSR - Allocates memory for a sparse sequential matrix in AIJ format. Input Parameters: + B - the matrix . i - the indices into j for the start of each row (starts with zero) . j - the column indices for each row (starts with zero) these must be sorted for each row - v - optional values in the matrix Contributed by: Lisandro Dalchin Level: developer .keywords: matrix, aij, compressed row, sparse, sequential .seealso: MatCreate(), MatCreateSeqAIJ(), MatSetValues(), MatSeqAIJSetPreallocation(), MatCreateSeqAIJ(), SeqAIJ @*/ PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat B,const PetscInt i[],const PetscInt j[],const PetscScalar v[]) { PetscErrorCode (*f)(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]); PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(B,MAT_COOKIE,1); ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqAIJSetPreallocationCSR_C",(void (**)(void))&f);CHKERRQ(ierr); if (f) { ierr = (*f)(B,i,j,v);CHKERRQ(ierr); } PetscFunctionReturn(0); } EXTERN_C_BEGIN #undef __FUNCT__ #define __FUNCT__ "MatSeqAIJSetPreallocationCSR_SeqAIJ" PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocationCSR_SeqAIJ(Mat B,const PetscInt Ii[],const PetscInt J[],const PetscScalar v[]) { PetscInt i; PetscInt m,n; PetscInt nz; PetscInt *nnz, nz_max = 0; PetscScalar *values; PetscErrorCode ierr; PetscFunctionBegin; ierr = MatGetSize(B, &m, &n);CHKERRQ(ierr); if (Ii[0]) { SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Ii[0] must be 0 it is %D", Ii[0]); } ierr = PetscMalloc((m+1) * sizeof(PetscInt), &nnz);CHKERRQ(ierr); for(i = 0; i < m; i++) { nz = Ii[i+1]- Ii[i]; nz_max = PetscMax(nz_max, nz); if (nz < 0) { SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE, "Local row %D has a negative number of columns %D", i, nnz); } nnz[i] = nz; } ierr = MatSeqAIJSetPreallocation(B, 0, nnz);CHKERRQ(ierr); ierr = PetscFree(nnz);CHKERRQ(ierr); if (v) { values = (PetscScalar*) v; } else { ierr = PetscMalloc((nz_max+1)*sizeof(PetscScalar), &values);CHKERRQ(ierr); ierr = PetscMemzero(values, nz_max*sizeof(PetscScalar));CHKERRQ(ierr); } ierr = MatSetOption(B,MAT_COLUMNS_SORTED);CHKERRQ(ierr); for(i = 0; i < m; i++) { nz = Ii[i+1] - Ii[i]; ierr = MatSetValues_SeqAIJ(B, 1, &i, nz, J+Ii[i], values + (v ? Ii[i] : 0), INSERT_VALUES);CHKERRQ(ierr); } ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); ierr = MatSetOption(B,MAT_COLUMNS_UNSORTED);CHKERRQ(ierr); if (!v) { ierr = PetscFree(values);CHKERRQ(ierr); } PetscFunctionReturn(0); } EXTERN_C_END /*MC MATSEQAIJ - MATSEQAIJ = "seqaij" - A matrix type to be used for sequential sparse matrices, based on compressed sparse row format. Options Database Keys: . -mat_type seqaij - sets the matrix type to "seqaij" during a call to MatSetFromOptions() Level: beginner .seealso: MatCreateSeqAIJ(), MatSetFromOptions(), MatSetType(), MatCreate(), MatType M*/ EXTERN_C_BEGIN extern PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SeqCRL(Mat,MatType,MatReuse,Mat*); EXTERN_C_END EXTERN_C_BEGIN #undef __FUNCT__ #define __FUNCT__ "MatCreate_SeqAIJ" PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SeqAIJ(Mat B) { Mat_SeqAIJ *b; PetscErrorCode ierr; PetscMPIInt size; PetscFunctionBegin; ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); if (size > 1) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Comm must be of size 1"); ierr = PetscNew(Mat_SeqAIJ,&b);CHKERRQ(ierr); B->data = (void*)b; ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); B->factor = 0; B->mapping = 0; b->row = 0; b->col = 0; b->icol = 0; b->reallocs = 0; b->sorted = PETSC_FALSE; b->ignorezeroentries = PETSC_FALSE; b->roworiented = PETSC_TRUE; b->nonew = 0; b->diag = 0; b->solve_work = 0; B->spptr = 0; b->saved_values = 0; b->idiag = 0; b->ssor = 0; b->keepzeroedrows = PETSC_FALSE; b->xtoy = 0; b->XtoY = 0; b->compressedrow.use = PETSC_FALSE; b->compressedrow.nrows = B->rmap.n; b->compressedrow.i = PETSC_NULL; b->compressedrow.rindex = PETSC_NULL; b->compressedrow.checked = PETSC_FALSE; B->same_nonzero = PETSC_FALSE; ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);CHKERRQ(ierr); ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqAIJSetColumnIndices_C", "MatSeqAIJSetColumnIndices_SeqAIJ", MatSeqAIJSetColumnIndices_SeqAIJ);CHKERRQ(ierr); ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatStoreValues_C", "MatStoreValues_SeqAIJ", MatStoreValues_SeqAIJ);CHKERRQ(ierr); ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatRetrieveValues_C", "MatRetrieveValues_SeqAIJ", MatRetrieveValues_SeqAIJ);CHKERRQ(ierr); ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_seqsbaij_C", "MatConvert_SeqAIJ_SeqSBAIJ", MatConvert_SeqAIJ_SeqSBAIJ);CHKERRQ(ierr); ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_seqbaij_C", "MatConvert_SeqAIJ_SeqBAIJ", MatConvert_SeqAIJ_SeqBAIJ);CHKERRQ(ierr); ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_seqcsrperm_C", "MatConvert_SeqAIJ_SeqCSRPERM", MatConvert_SeqAIJ_SeqCSRPERM);CHKERRQ(ierr); ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqaij_seqcrl_C", "MatConvert_SeqAIJ_SeqCRL", MatConvert_SeqAIJ_SeqCRL);CHKERRQ(ierr); ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatIsTranspose_C", "MatIsTranspose_SeqAIJ", MatIsTranspose_SeqAIJ);CHKERRQ(ierr); ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqAIJSetPreallocation_C", "MatSeqAIJSetPreallocation_SeqAIJ", MatSeqAIJSetPreallocation_SeqAIJ);CHKERRQ(ierr); ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqAIJSetPreallocationCSR_C", "MatSeqAIJSetPreallocationCSR_SeqAIJ", MatSeqAIJSetPreallocationCSR_SeqAIJ);CHKERRQ(ierr); ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatReorderForNonzeroDiagonal_C", "MatReorderForNonzeroDiagonal_SeqAIJ", MatReorderForNonzeroDiagonal_SeqAIJ);CHKERRQ(ierr); ierr = MatCreate_Inode(B);CHKERRQ(ierr); ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQAIJ);CHKERRQ(ierr); PetscFunctionReturn(0); } EXTERN_C_END #undef __FUNCT__ #define __FUNCT__ "MatDuplicate_SeqAIJ" PetscErrorCode MatDuplicate_SeqAIJ(Mat A,MatDuplicateOption cpvalues,Mat *B) { Mat C; Mat_SeqAIJ *c,*a = (Mat_SeqAIJ*)A->data; PetscErrorCode ierr; PetscInt i,m = A->rmap.n; PetscFunctionBegin; *B = 0; ierr = MatCreate(A->comm,&C);CHKERRQ(ierr); ierr = MatSetSizes(C,A->rmap.n,A->cmap.n,A->rmap.n,A->cmap.n);CHKERRQ(ierr); ierr = MatSetType(C,A->type_name);CHKERRQ(ierr); ierr = PetscMemcpy(C->ops,A->ops,sizeof(struct _MatOps));CHKERRQ(ierr); ierr = PetscMapCopy(A->comm,&A->rmap,&C->rmap);CHKERRQ(ierr); ierr = PetscMapCopy(A->comm,&A->cmap,&C->cmap);CHKERRQ(ierr); c = (Mat_SeqAIJ*)C->data; C->factor = A->factor; c->row = 0; c->col = 0; c->icol = 0; c->reallocs = 0; C->assembled = PETSC_TRUE; ierr = PetscMalloc2(m,PetscInt,&c->imax,m,PetscInt,&c->ilen);CHKERRQ(ierr); for (i=0; iimax[i] = a->imax[i]; c->ilen[i] = a->ilen[i]; } /* allocate the matrix space */ ierr = PetscMalloc3(a->i[m],PetscScalar,&c->a,a->i[m],PetscInt,&c->j,m+1,PetscInt,&c->i);CHKERRQ(ierr); c->singlemalloc = PETSC_TRUE; ierr = PetscMemcpy(c->i,a->i,(m+1)*sizeof(PetscInt));CHKERRQ(ierr); if (m > 0) { ierr = PetscMemcpy(c->j,a->j,(a->i[m])*sizeof(PetscInt));CHKERRQ(ierr); if (cpvalues == MAT_COPY_VALUES) { ierr = PetscMemcpy(c->a,a->a,(a->i[m])*sizeof(PetscScalar));CHKERRQ(ierr); } else { ierr = PetscMemzero(c->a,(a->i[m])*sizeof(PetscScalar));CHKERRQ(ierr); } } c->sorted = a->sorted; c->ignorezeroentries = a->ignorezeroentries; c->roworiented = a->roworiented; c->nonew = a->nonew; if (a->diag) { ierr = PetscMalloc((m+1)*sizeof(PetscInt),&c->diag);CHKERRQ(ierr); ierr = PetscLogObjectMemory(C,(m+1)*sizeof(PetscInt));CHKERRQ(ierr); for (i=0; idiag[i] = a->diag[i]; } } else c->diag = 0; c->solve_work = 0; c->saved_values = 0; c->idiag = 0; c->ssor = 0; c->keepzeroedrows = a->keepzeroedrows; c->free_a = PETSC_TRUE; c->free_ij = PETSC_TRUE; c->xtoy = 0; c->XtoY = 0; c->nz = a->nz; c->maxnz = a->maxnz; C->preallocated = PETSC_TRUE; c->compressedrow.use = a->compressedrow.use; c->compressedrow.nrows = a->compressedrow.nrows; c->compressedrow.checked = a->compressedrow.checked; if ( a->compressedrow.checked && a->compressedrow.use){ i = a->compressedrow.nrows; ierr = PetscMalloc((2*i+1)*sizeof(PetscInt),&c->compressedrow.i);CHKERRQ(ierr); c->compressedrow.rindex = c->compressedrow.i + i + 1; ierr = PetscMemcpy(c->compressedrow.i,a->compressedrow.i,(i+1)*sizeof(PetscInt));CHKERRQ(ierr); ierr = PetscMemcpy(c->compressedrow.rindex,a->compressedrow.rindex,i*sizeof(PetscInt));CHKERRQ(ierr); } else { c->compressedrow.use = PETSC_FALSE; c->compressedrow.i = PETSC_NULL; c->compressedrow.rindex = PETSC_NULL; } C->same_nonzero = A->same_nonzero; ierr = MatDuplicate_Inode(A,cpvalues,&C);CHKERRQ(ierr); *B = C; ierr = PetscFListDuplicate(A->qlist,&C->qlist);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatLoad_SeqAIJ" PetscErrorCode MatLoad_SeqAIJ(PetscViewer viewer, MatType type,Mat *A) { Mat_SeqAIJ *a; Mat B; PetscErrorCode ierr; PetscInt i,sum,nz,header[4],*rowlengths = 0,M,N; int fd; PetscMPIInt size; MPI_Comm comm; PetscFunctionBegin; ierr = PetscObjectGetComm((PetscObject)viewer,&comm);CHKERRQ(ierr); ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); if (size > 1) SETERRQ(PETSC_ERR_ARG_SIZ,"view must have one processor"); ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"not matrix object in file"); M = header[1]; N = header[2]; nz = header[3]; if (nz < 0) { SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Matrix stored in special format on disk,cannot load as SeqAIJ"); } /* read in row lengths */ ierr = PetscMalloc(M*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr); ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); /* check if sum of rowlengths is same as nz */ for (i=0,sum=0; i< M; i++) sum +=rowlengths[i]; if (sum != nz) SETERRQ2(PETSC_ERR_FILE_READ,"Inconsistant matrix data in file. no-nonzeros = %d, sum-row-lengths = %d\n",nz,sum); /* create our matrix */ ierr = MatCreate(comm,&B);CHKERRQ(ierr); ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr); ierr = MatSetType(B,type);CHKERRQ(ierr); ierr = MatSeqAIJSetPreallocation_SeqAIJ(B,0,rowlengths);CHKERRQ(ierr); a = (Mat_SeqAIJ*)B->data; ierr = PetscBinaryRead(fd,a->j,nz,PETSC_INT);CHKERRQ(ierr); /* read in nonzero values */ ierr = PetscBinaryRead(fd,a->a,nz,PETSC_SCALAR);CHKERRQ(ierr); /* set matrix "i" values */ a->i[0] = 0; for (i=1; i<= M; i++) { a->i[i] = a->i[i-1] + rowlengths[i-1]; a->ilen[i-1] = rowlengths[i-1]; } ierr = PetscFree(rowlengths);CHKERRQ(ierr); ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); *A = B; PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatEqual_SeqAIJ" PetscErrorCode MatEqual_SeqAIJ(Mat A,Mat B,PetscTruth* flg) { Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data,*b = (Mat_SeqAIJ *)B->data; PetscErrorCode ierr; PetscFunctionBegin; /* If the matrix dimensions are not equal,or no of nonzeros */ if ((A->rmap.n != B->rmap.n) || (A->cmap.n != B->cmap.n) ||(a->nz != b->nz)) { *flg = PETSC_FALSE; PetscFunctionReturn(0); } /* if the a->i are the same */ ierr = PetscMemcmp(a->i,b->i,(A->rmap.n+1)*sizeof(PetscInt),flg);CHKERRQ(ierr); if (!*flg) PetscFunctionReturn(0); /* if a->j are the same */ ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),flg);CHKERRQ(ierr); if (!*flg) PetscFunctionReturn(0); /* if a->a are the same */ ierr = PetscMemcmp(a->a,b->a,(a->nz)*sizeof(PetscScalar),flg);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatCreateSeqAIJWithArrays" /*@ MatCreateSeqAIJWithArrays - Creates an sequential AIJ matrix using matrix elements (in CSR format) provided by the user. Collective on MPI_Comm Input Parameters: + comm - must be an MPI communicator of size 1 . m - number of rows . n - number of columns . i - row indices . j - column indices - a - matrix values Output Parameter: . mat - the matrix Level: intermediate Notes: The i, j, and a arrays are not copied by this routine, the user must free these arrays once the matrix is destroyed You cannot set new nonzero locations into this matrix, that will generate an error. The i and j indices are 0 based .seealso: MatCreate(), MatCreateMPIAIJ(), MatCreateSeqAIJ(), MatCreateMPIAIJWithArrays(), MatMPIAIJSetPreallocationCSR() @*/ PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJWithArrays(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt* i,PetscInt*j,PetscScalar *a,Mat *mat) { PetscErrorCode ierr; PetscInt ii; Mat_SeqAIJ *aij; PetscFunctionBegin; if (i[0]) { SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"i (row indices) must start with 0"); } ierr = MatCreate(comm,mat);CHKERRQ(ierr); ierr = MatSetSizes(*mat,m,n,m,n);CHKERRQ(ierr); ierr = MatSetType(*mat,MATSEQAIJ);CHKERRQ(ierr); ierr = MatSeqAIJSetPreallocation_SeqAIJ(*mat,MAT_SKIP_ALLOCATION,0);CHKERRQ(ierr); aij = (Mat_SeqAIJ*)(*mat)->data; ierr = PetscMalloc2(m,PetscInt,&aij->imax,m,PetscInt,&aij->ilen);CHKERRQ(ierr); aij->i = i; aij->j = j; aij->a = a; aij->singlemalloc = PETSC_FALSE; aij->nonew = -1; /*this indicates that inserting a new value in the matrix that generates a new nonzero is an error*/ aij->free_a = PETSC_FALSE; aij->free_ij = PETSC_FALSE; for (ii=0; iiilen[ii] = aij->imax[ii] = i[ii+1] - i[ii]; #if defined(PETSC_USE_DEBUG) if (i[ii+1] - i[ii] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Negative row length in i (row indices) row = %d length = %d",ii,i[ii+1] - i[ii]); #endif } #if defined(PETSC_USE_DEBUG) for (ii=0; iii[m]; ii++) { if (j[ii] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Negative column index at location = %d index = %d",ii,j[ii]); if (j[ii] > n - 1) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column index to large at location = %d index = %d",ii,j[ii]); } #endif ierr = MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); ierr = MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatSetColoring_SeqAIJ" PetscErrorCode MatSetColoring_SeqAIJ(Mat A,ISColoring coloring) { PetscErrorCode ierr; Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; PetscFunctionBegin; if (coloring->ctype == IS_COLORING_GLOBAL) { ierr = ISColoringReference(coloring);CHKERRQ(ierr); a->coloring = coloring; } else if (coloring->ctype == IS_COLORING_GHOSTED) { PetscInt i,*larray; ISColoring ocoloring; ISColoringValue *colors; /* set coloring for diagonal portion */ ierr = PetscMalloc((A->cmap.n+1)*sizeof(PetscInt),&larray);CHKERRQ(ierr); for (i=0; icmap.n; i++) { larray[i] = i; } ierr = ISGlobalToLocalMappingApply(A->mapping,IS_GTOLM_MASK,A->cmap.n,larray,PETSC_NULL,larray);CHKERRQ(ierr); ierr = PetscMalloc((A->cmap.n+1)*sizeof(ISColoringValue),&colors);CHKERRQ(ierr); for (i=0; icmap.n; i++) { colors[i] = coloring->colors[larray[i]]; } ierr = PetscFree(larray);CHKERRQ(ierr); ierr = ISColoringCreate(PETSC_COMM_SELF,coloring->n,A->cmap.n,colors,&ocoloring);CHKERRQ(ierr); a->coloring = ocoloring; } PetscFunctionReturn(0); } #if defined(PETSC_HAVE_ADIC) EXTERN_C_BEGIN #include "adic/ad_utils.h" EXTERN_C_END #undef __FUNCT__ #define __FUNCT__ "MatSetValuesAdic_SeqAIJ" PetscErrorCode MatSetValuesAdic_SeqAIJ(Mat A,void *advalues) { Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; PetscInt m = A->rmap.n,*ii = a->i,*jj = a->j,nz,i,j,nlen; PetscScalar *v = a->a,*values = ((PetscScalar*)advalues)+1; ISColoringValue *color; PetscFunctionBegin; if (!a->coloring) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Coloring not set for matrix"); nlen = PetscADGetDerivTypeSize()/sizeof(PetscScalar); color = a->coloring->colors; /* loop over rows */ for (i=0; idata; PetscInt m = A->rmap.n,*ii = a->i,*jj = a->j,nz,i,j; PetscScalar *v = a->a,*values = (PetscScalar *)advalues; ISColoringValue *color; PetscFunctionBegin; if (!a->coloring) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Coloring not set for matrix"); color = a->coloring->colors; /* loop over rows */ for (i=0; icomm,ierr) #undef SETERRQ2 #define SETERRQ2(ierr,b,c,d) CHKERRABORT(A->comm,ierr) EXTERN_C_BEGIN #undef __FUNCT__ #define __FUNCT__ "matsetvaluesseqaij_" void PETSC_STDCALL matsetvaluesseqaij_(Mat *AA,PetscInt *mm,const PetscInt im[],PetscInt *nn,const PetscInt in[],const PetscScalar v[],InsertMode *isis, PetscErrorCode *_ierr) { Mat A = *AA; PetscInt m = *mm, n = *nn; InsertMode is = *isis; Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; PetscInt *rp,k,low,high,t,ii,row,nrow,i,col,l,rmax,N; PetscInt *imax,*ai,*ailen; PetscErrorCode ierr; PetscInt *aj,nonew = a->nonew,lastcol = -1; PetscScalar *ap,value,*aa; PetscTruth ignorezeroentries = ((a->ignorezeroentries && is == ADD_VALUES) ? PETSC_TRUE:PETSC_FALSE); PetscTruth roworiented = a->roworiented; PetscFunctionBegin; MatPreallocated(A); imax = a->imax; ai = a->i; ailen = a->ilen; aj = a->j; aa = a->a; for (k=0; k= A->rmap.n) SETERRABORT(A->comm,PETSC_ERR_ARG_OUTOFRANGE,"Row too large"); #endif rp = aj + ai[row]; ap = aa + ai[row]; rmax = imax[row]; nrow = ailen[row]; low = 0; high = nrow; for (l=0; l= A->cmap.n) SETERRABORT(A->comm,PETSC_ERR_ARG_OUTOFRANGE,"Column too large"); #endif col = in[l]; if (roworiented) { value = v[l + k*n]; } else { value = v[k + l*m]; } if (value == 0.0 && ignorezeroentries && (is == ADD_VALUES)) continue; if (col <= lastcol) low = 0; else high = nrow; lastcol = col; while (high-low > 5) { t = (low+high)/2; if (rp[t] > col) high = t; else low = t; } for (i=low; i col) break; if (rp[i] == col) { if (is == ADD_VALUES) ap[i] += value; else ap[i] = value; goto noinsert; } } if (value == 0.0 && ignorezeroentries) goto noinsert; if (nonew == 1) goto noinsert; if (nonew == -1) SETERRABORT(A->comm,PETSC_ERR_ARG_OUTOFRANGE,"Inserting a new nonzero in the matrix"); MatSeqXAIJReallocateAIJ(A,A->rmap.n,1,nrow,row,col,rmax,aa,ai,aj,rp,ap,imax,nonew,MatScalar); N = nrow++ - 1; a->nz++; high++; /* shift up all the later entries in this row */ for (ii=N; ii>=i; ii--) { rp[ii+1] = rp[ii]; ap[ii+1] = ap[ii]; } rp[i] = col; ap[i] = value; noinsert:; low = i + 1; } ailen[row] = nrow; } A->same_nonzero = PETSC_FALSE; PetscFunctionReturnVoid(); } EXTERN_C_END