#include "src/mat/impls/baij/seq/baij.h" #include "src/inline/spops.h" #include "src/inline/ilu.h" #include "petscbt.h" #undef __FUNCT__ #define __FUNCT__ "MatIncreaseOverlap_SeqBAIJ" PetscErrorCode MatIncreaseOverlap_SeqBAIJ(Mat A,int is_max,IS is[],int ov) { Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; PetscErrorCode ierr; int row,i,j,k,l,m,n,*idx,*nidx,isz,val,ival; int start,end,*ai,*aj,bs,*nidx2; PetscBT table; PetscFunctionBegin; m = a->mbs; ai = a->i; aj = a->j; bs = a->bs; if (ov < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified"); ierr = PetscBTCreate(m,table);CHKERRQ(ierr); ierr = PetscMalloc((m+1)*sizeof(int),&nidx);CHKERRQ(ierr); ierr = PetscMalloc((A->m+1)*sizeof(int),&nidx2);CHKERRQ(ierr); for (i=0; i=m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"index greater than mat-dim"); if(!PetscBTLookupSet(table,ival)) { nidx[isz++] = ival;} } ierr = ISRestoreIndices(is[i],&idx);CHKERRQ(ierr); ierr = ISDestroy(is[i]);CHKERRQ(ierr); k = 0; for (j=0; jdata,*c; PetscErrorCode ierr; int *smap,i,k,kstart,kend,oldcols = a->nbs,*lens; int row,mat_i,*mat_j,tcol,*mat_ilen; int *irow,*icol,nrows,ncols,*ssmap,bs=a->bs,bs2=a->bs2; int *aj = a->j,*ai = a->i; MatScalar *mat_a; Mat C; PetscTruth flag; PetscFunctionBegin; ierr = ISSorted(iscol,(PetscTruth*)&i);CHKERRQ(ierr); if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"IS is not sorted"); ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); ierr = PetscMalloc((1+oldcols)*sizeof(int),&smap);CHKERRQ(ierr); ssmap = smap; ierr = PetscMalloc((1+nrows)*sizeof(int),&lens);CHKERRQ(ierr); ierr = PetscMemzero(smap,oldcols*sizeof(int));CHKERRQ(ierr); for (i=0; iilen[irow[i]]; lens[i] = 0; for (k=kstart; kdata); if (c->mbs!=nrows || c->nbs!=ncols || c->bs!=bs) SETERRQ(PETSC_ERR_ARG_SIZ,"Submatrix wrong size"); ierr = PetscMemcmp(c->ilen,lens,c->mbs *sizeof(int),&flag);CHKERRQ(ierr); if (flag == PETSC_FALSE) { SETERRQ(PETSC_ERR_ARG_SIZ,"Cannot reuse matrix. wrong no of nonzeros"); } ierr = PetscMemzero(c->ilen,c->mbs*sizeof(int));CHKERRQ(ierr); C = *B; } else { ierr = MatCreate(A->comm,nrows*bs,ncols*bs,PETSC_DETERMINE,PETSC_DETERMINE,&C);CHKERRQ(ierr); ierr = MatSetType(C,A->type_name);CHKERRQ(ierr); ierr = MatSeqBAIJSetPreallocation(C,bs,0,lens);CHKERRQ(ierr); } c = (Mat_SeqBAIJ *)(C->data); for (i=0; iilen[row]; mat_i = c->i[i]; mat_j = c->j + mat_i; mat_a = c->a + mat_i*bs2; mat_ilen = c->ilen + i; for (k=kstart; kj[k]])) { *mat_j++ = tcol - 1; ierr = PetscMemcpy(mat_a,a->a+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); mat_a += bs2; (*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__ "MatGetSubMatrix_SeqBAIJ" PetscErrorCode MatGetSubMatrix_SeqBAIJ(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B) { Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; IS is1,is2; PetscErrorCode ierr; int *vary,*iary,*irow,*icol,nrows,ncols,i,bs=a->bs,count; PetscFunctionBegin; ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); /* Verify if the indices corespond to each element in a block and form the IS with compressed IS */ ierr = PetscMalloc(2*(a->mbs+1)*sizeof(int),&vary);CHKERRQ(ierr); iary = vary + a->mbs; ierr = PetscMemzero(vary,(a->mbs)*sizeof(int));CHKERRQ(ierr); for (i=0; imbs; i++) { if (vary[i]!=0 && vary[i]!=bs) SETERRQ(PETSC_ERR_ARG_SIZ,"Index set does not match blocks"); if (vary[i]==bs) iary[count++] = i; } ierr = ISCreateGeneral(PETSC_COMM_SELF,count,iary,&is1);CHKERRQ(ierr); ierr = PetscMemzero(vary,(a->mbs)*sizeof(int));CHKERRQ(ierr); for (i=0; imbs; i++) { if (vary[i]!=0 && vary[i]!=bs) SETERRQ(PETSC_ERR_PLIB,"Internal error in PETSc"); if (vary[i]==bs) iary[count++] = i; } ierr = ISCreateGeneral(PETSC_COMM_SELF,count,iary,&is2);CHKERRQ(ierr); ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); ierr = PetscFree(vary);CHKERRQ(ierr); ierr = MatGetSubMatrix_SeqBAIJ_Private(A,is1,is2,cs,scall,B);CHKERRQ(ierr); ISDestroy(is1); ISDestroy(is2); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatGetSubMatrices_SeqBAIJ" PetscErrorCode MatGetSubMatrices_SeqBAIJ(Mat A,int n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[]) { PetscErrorCode ierr; int i; PetscFunctionBegin; if (scall == MAT_INITIAL_MATRIX) { ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr); } for (i=0; idata; PetscScalar *x,*z,sum; MatScalar *v; PetscErrorCode ierr; int mbs=a->mbs,i,*idx,*ii,n; PetscFunctionBegin; ierr = VecGetArray(xx,&x);CHKERRQ(ierr); ierr = VecGetArray(zz,&z);CHKERRQ(ierr); idx = a->j; v = a->a; ii = a->i; for (i=0; inz - A->m); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatMult_SeqBAIJ_2" PetscErrorCode MatMult_SeqBAIJ_2(Mat A,Vec xx,Vec zz) { Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; PetscScalar *x,*z,*xb,sum1,sum2; PetscScalar x1,x2; MatScalar *v; PetscErrorCode ierr; int mbs=a->mbs,i,*idx,*ii,j,n; PetscFunctionBegin; ierr = VecGetArray(xx,&x);CHKERRQ(ierr); ierr = VecGetArray(zz,&z);CHKERRQ(ierr); idx = a->j; v = a->a; ii = a->i; for (i=0; inz - A->m); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatMult_SeqBAIJ_3" PetscErrorCode MatMult_SeqBAIJ_3(Mat A,Vec xx,Vec zz) { Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; PetscScalar *x,*z,*xb,sum1,sum2,sum3,x1,x2,x3; MatScalar *v; PetscErrorCode ierr; int mbs=a->mbs,i,*idx,*ii,j,n; #if defined(PETSC_HAVE_PRAGMA_DISJOINT) #pragma disjoint(*v,*z,*xb) #endif PetscFunctionBegin; ierr = VecGetArray(xx,&x);CHKERRQ(ierr); ierr = VecGetArray(zz,&z);CHKERRQ(ierr); idx = a->j; v = a->a; ii = a->i; for (i=0; inz - A->m); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatMult_SeqBAIJ_4" PetscErrorCode MatMult_SeqBAIJ_4(Mat A,Vec xx,Vec zz) { Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; PetscScalar *x,*z,*xb,sum1,sum2,sum3,sum4,x1,x2,x3,x4; MatScalar *v; PetscErrorCode ierr; int mbs=a->mbs,i,*idx,*ii,j,n; PetscFunctionBegin; ierr = VecGetArray(xx,&x);CHKERRQ(ierr); ierr = VecGetArray(zz,&z);CHKERRQ(ierr); idx = a->j; v = a->a; ii = a->i; for (i=0; inz - A->m); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatMult_SeqBAIJ_5" PetscErrorCode MatMult_SeqBAIJ_5(Mat A,Vec xx,Vec zz) { Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; PetscScalar sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5,*xb,*z,*x; MatScalar *v; PetscErrorCode ierr; int mbs=a->mbs,i,*idx,*ii,j,n; PetscFunctionBegin; ierr = VecGetArray(xx,&x);CHKERRQ(ierr); ierr = VecGetArray(zz,&z);CHKERRQ(ierr); idx = a->j; v = a->a; ii = a->i; for (i=0; inz - A->m); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatMult_SeqBAIJ_6" PetscErrorCode MatMult_SeqBAIJ_6(Mat A,Vec xx,Vec zz) { Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; PetscScalar *x,*z,*xb,sum1,sum2,sum3,sum4,sum5,sum6; PetscScalar x1,x2,x3,x4,x5,x6; MatScalar *v; PetscErrorCode ierr; int mbs=a->mbs,i,*idx,*ii,j,n; PetscFunctionBegin; ierr = VecGetArray(xx,&x);CHKERRQ(ierr); ierr = VecGetArray(zz,&z);CHKERRQ(ierr); idx = a->j; v = a->a; ii = a->i; for (i=0; inz - A->m); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatMult_SeqBAIJ_7" PetscErrorCode MatMult_SeqBAIJ_7(Mat A,Vec xx,Vec zz) { Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; PetscScalar *x,*z,*xb,sum1,sum2,sum3,sum4,sum5,sum6,sum7; PetscScalar x1,x2,x3,x4,x5,x6,x7; MatScalar *v; PetscErrorCode ierr; int mbs=a->mbs,i,*idx,*ii,j,n; PetscFunctionBegin; ierr = VecGetArray(xx,&x);CHKERRQ(ierr); ierr = VecGetArray(zz,&z);CHKERRQ(ierr); idx = a->j; v = a->a; ii = a->i; for (i=0; inz - A->m); PetscFunctionReturn(0); } /* This will not work with MatScalar == float because it calls the BLAS */ #undef __FUNCT__ #define __FUNCT__ "MatMult_SeqBAIJ_N" PetscErrorCode MatMult_SeqBAIJ_N(Mat A,Vec xx,Vec zz) { Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; PetscScalar *x,*z,*xb,*work,*workt; MatScalar *v; PetscErrorCode ierr; int mbs=a->mbs,i,*idx,*ii,bs=a->bs,j,n,bs2=a->bs2; int ncols,k; PetscFunctionBegin; ierr = VecGetArray(xx,&x);CHKERRQ(ierr); ierr = VecGetArray(zz,&z);CHKERRQ(ierr); idx = a->j; v = a->a; ii = a->i; if (!a->mult_work) { k = PetscMax(A->m,A->n); ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr); } work = a->mult_work; for (i=0; inz*bs2 - A->m); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatMultAdd_SeqBAIJ_1" PetscErrorCode MatMultAdd_SeqBAIJ_1(Mat A,Vec xx,Vec yy,Vec zz) { Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; PetscScalar *x,*y,*z,sum; MatScalar *v; PetscErrorCode ierr; int mbs=a->mbs,i,*idx,*ii,n; 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; } idx = a->j; v = a->a; ii = a->i; for (i=0; inz); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatMultAdd_SeqBAIJ_2" PetscErrorCode MatMultAdd_SeqBAIJ_2(Mat A,Vec xx,Vec yy,Vec zz) { Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; PetscScalar *x,*y,*z,*xb,sum1,sum2; PetscScalar x1,x2; MatScalar *v; PetscErrorCode ierr; int mbs=a->mbs,i,*idx,*ii,j,n; 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; } idx = a->j; v = a->a; ii = a->i; for (i=0; inz); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatMultAdd_SeqBAIJ_3" PetscErrorCode MatMultAdd_SeqBAIJ_3(Mat A,Vec xx,Vec yy,Vec zz) { Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; PetscScalar *x,*y,*z,*xb,sum1,sum2,sum3,x1,x2,x3; MatScalar *v; PetscErrorCode ierr; int mbs=a->mbs,i,*idx,*ii,j,n; 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; } idx = a->j; v = a->a; ii = a->i; for (i=0; inz); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatMultAdd_SeqBAIJ_4" PetscErrorCode MatMultAdd_SeqBAIJ_4(Mat A,Vec xx,Vec yy,Vec zz) { Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; PetscScalar *x,*y,*z,*xb,sum1,sum2,sum3,sum4,x1,x2,x3,x4; MatScalar *v; PetscErrorCode ierr; int mbs=a->mbs,i,*idx,*ii; int j,n; 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; } idx = a->j; v = a->a; ii = a->i; for (i=0; inz); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatMultAdd_SeqBAIJ_5" PetscErrorCode MatMultAdd_SeqBAIJ_5(Mat A,Vec xx,Vec yy,Vec zz) { Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; PetscScalar *x,*y,*z,*xb,sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5; MatScalar *v; PetscErrorCode ierr; int mbs=a->mbs,i,*idx,*ii,j,n; 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; } idx = a->j; v = a->a; ii = a->i; for (i=0; inz); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatMultAdd_SeqBAIJ_6" PetscErrorCode MatMultAdd_SeqBAIJ_6(Mat A,Vec xx,Vec yy,Vec zz) { Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; PetscScalar *x,*y,*z,*xb,sum1,sum2,sum3,sum4,sum5,sum6; PetscScalar x1,x2,x3,x4,x5,x6; MatScalar *v; PetscErrorCode ierr; int mbs=a->mbs,i,*idx,*ii,j,n; 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; } idx = a->j; v = a->a; ii = a->i; for (i=0; inz); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatMultAdd_SeqBAIJ_7" PetscErrorCode MatMultAdd_SeqBAIJ_7(Mat A,Vec xx,Vec yy,Vec zz) { Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; PetscScalar *x,*y,*z,*xb,sum1,sum2,sum3,sum4,sum5,sum6,sum7; PetscScalar x1,x2,x3,x4,x5,x6,x7; MatScalar *v; PetscErrorCode ierr; int mbs=a->mbs,i,*idx,*ii,j,n; 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; } idx = a->j; v = a->a; ii = a->i; for (i=0; inz); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatMultAdd_SeqBAIJ_N" PetscErrorCode MatMultAdd_SeqBAIJ_N(Mat A,Vec xx,Vec yy,Vec zz) { Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; PetscScalar *x,*z,*xb,*work,*workt,*y; MatScalar *v; PetscErrorCode ierr; int mbs=a->mbs,i,*idx,*ii,bs=a->bs,j,n,bs2=a->bs2; int ncols,k; PetscFunctionBegin; ierr = VecGetArray(xx,&x);CHKERRQ(ierr); ierr = VecGetArray(zz,&z);CHKERRQ(ierr); if (zz != yy) { ierr = VecGetArray(yy,&y);CHKERRQ(ierr); ierr = PetscMemcpy(z,y,yy->n*sizeof(PetscScalar));CHKERRQ(ierr); ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); } idx = a->j; v = a->a; ii = a->i; if (!a->mult_work) { k = PetscMax(A->m,A->n); ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr); } work = a->mult_work; for (i=0; inz*bs2); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatMultTranspose_SeqBAIJ" PetscErrorCode MatMultTranspose_SeqBAIJ(Mat A,Vec xx,Vec zz) { Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; PetscScalar *xg,*zg,*zb,zero = 0.0; PetscScalar *x,*z,*xb,x1,x2,x3,x4,x5,x6,x7; MatScalar *v; int mbs=a->mbs,i,*idx,*ii,*ai=a->i,rval; PetscErrorCode ierr; int bs=a->bs,j,n,bs2=a->bs2,*ib; PetscFunctionBegin; ierr = VecSet(&zero,zz);CHKERRQ(ierr); ierr = VecGetArray(xx,&xg);CHKERRQ(ierr); x = xg; ierr = VecGetArray(zz,&zg);CHKERRQ(ierr); z = zg; idx = a->j; v = a->a; ii = a->i; xb = x; switch (bs) { case 1: for (i=0; imult_work) { k = PetscMax(A->m,A->n); ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr); } work = a->mult_work; for (i=0; inz*a->bs2 - A->n); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatMultTransposeAdd_SeqBAIJ" PetscErrorCode MatMultTransposeAdd_SeqBAIJ(Mat A,Vec xx,Vec yy,Vec zz) { Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; PetscScalar *xg,*zg,*zb,*x,*z,*xb,x1,x2,x3,x4,x5; MatScalar *v; PetscErrorCode ierr; int mbs=a->mbs,i,*idx,*ii,*ai=a->i,rval,bs=a->bs,j,n,bs2=a->bs2,*ib; PetscFunctionBegin; if (yy != zz) { ierr = VecCopy(yy,zz);CHKERRQ(ierr); } ierr = VecGetArray(xx,&xg);CHKERRQ(ierr); x = xg; ierr = VecGetArray(zz,&zg);CHKERRQ(ierr); z = zg; idx = a->j; v = a->a; ii = a->i; xb = x; switch (bs) { case 1: for (i=0; imult_work) { k = PetscMax(A->m,A->n); ierr = PetscMalloc((k+1)*sizeof(PetscScalar),&a->mult_work);CHKERRQ(ierr); } work = a->mult_work; for (i=0; inz*a->bs2); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatScale_SeqBAIJ" PetscErrorCode MatScale_SeqBAIJ(const PetscScalar *alpha,Mat inA) { Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)inA->data; int totalnz = a->bs2*a->nz; #if defined(PETSC_USE_MAT_SINGLE) int i; #else int one = 1; #endif PetscFunctionBegin; #if defined(PETSC_USE_MAT_SINGLE) for (i=0; ia[i] *= *alpha; #else BLscal_(&totalnz,(PetscScalar*)alpha,a->a,&one); #endif PetscLogFlops(totalnz); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatNorm_SeqBAIJ" PetscErrorCode MatNorm_SeqBAIJ(Mat A,NormType type,PetscReal *norm) { Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; MatScalar *v = a->a; PetscReal sum = 0.0; int i,j,k,bs = a->bs,nz=a->nz,bs2=a->bs2,k1; PetscFunctionBegin; if (type == NORM_FROBENIUS) { for (i=0; i< bs2*nz; i++) { #if defined(PETSC_USE_COMPLEX) sum += PetscRealPart(PetscConj(*v)*(*v)); v++; #else sum += (*v)*(*v); v++; #endif } *norm = sqrt(sum); } else if (type == NORM_INFINITY) { /* maximum row sum */ *norm = 0.0; for (k=0; kmbs; j++) { v = a->a + bs2*a->i[j] + k; sum = 0.0; for (i=0; ii[j+1]-a->i[j]; i++) { for (k1=0; k1 *norm) *norm = sum; } } } else { SETERRQ(PETSC_ERR_SUP,"No support for this norm yet"); } PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatEqual_SeqBAIJ" PetscErrorCode MatEqual_SeqBAIJ(Mat A,Mat B,PetscTruth* flg) { Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data,*b = (Mat_SeqBAIJ *)B->data; PetscErrorCode ierr; PetscFunctionBegin; /* If the matrix/block dimensions are not equal, or no of nonzeros or shift */ if ((A->m != B->m) || (A->n != B->n) || (a->bs != b->bs)|| (a->nz != b->nz)) { *flg = PETSC_FALSE; PetscFunctionReturn(0); } /* if the a->i are the same */ ierr = PetscMemcmp(a->i,b->i,(a->mbs+1)*sizeof(int),flg);CHKERRQ(ierr); if (*flg == PETSC_FALSE) { PetscFunctionReturn(0); } /* if a->j are the same */ ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(int),flg);CHKERRQ(ierr); if (*flg == PETSC_FALSE) { PetscFunctionReturn(0); } /* if a->a are the same */ ierr = PetscMemcmp(a->a,b->a,(a->nz)*(a->bs)*(a->bs)*sizeof(PetscScalar),flg);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatGetDiagonal_SeqBAIJ" PetscErrorCode MatGetDiagonal_SeqBAIJ(Mat A,Vec v) { Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; PetscErrorCode ierr; int i,j,k,n,row,bs,*ai,*aj,ambs,bs2; PetscScalar *x,zero = 0.0; MatScalar *aa,*aa_j; PetscFunctionBegin; if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); bs = a->bs; aa = a->a; ai = a->i; aj = a->j; ambs = a->mbs; bs2 = a->bs2; ierr = VecSet(&zero,v);CHKERRQ(ierr); ierr = VecGetArray(v,&x);CHKERRQ(ierr); ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr); if (n != A->m) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); for (i=0; idata; PetscScalar *l,*r,x,*li,*ri; MatScalar *aa,*v; PetscErrorCode ierr; int i,j,k,lm,rn,M,m,n,*ai,*aj,mbs,tmp,bs,bs2; PetscFunctionBegin; ai = a->i; aj = a->j; aa = a->a; m = A->m; n = A->n; bs = a->bs; mbs = a->mbs; bs2 = a->bs2; if (ll) { ierr = VecGetArray(ll,&l);CHKERRQ(ierr); ierr = VecGetLocalSize(ll,&lm);CHKERRQ(ierr); if (lm != m) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vector wrong length"); for (i=0; inz); } if (rr) { ierr = VecGetArray(rr,&r);CHKERRQ(ierr); ierr = VecGetLocalSize(rr,&rn);CHKERRQ(ierr); if (rn != n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vector wrong length"); for (i=0; inz); } PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatGetInfo_SeqBAIJ" PetscErrorCode MatGetInfo_SeqBAIJ(Mat A,MatInfoType flag,MatInfo *info) { Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; PetscFunctionBegin; info->rows_global = (double)A->m; info->columns_global = (double)A->n; info->rows_local = (double)A->m; info->columns_local = (double)A->n; info->block_size = a->bs2; info->nz_allocated = a->maxnz; info->nz_used = a->bs2*a->nz; info->nz_unneeded = (double)(info->nz_allocated - info->nz_used); info->assemblies = A->num_ass; info->mallocs = 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__ "MatZeroEntries_SeqBAIJ" PetscErrorCode MatZeroEntries_SeqBAIJ(Mat A) { Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; PetscErrorCode ierr; PetscFunctionBegin; ierr = PetscMemzero(a->a,a->bs2*a->i[a->mbs]*sizeof(MatScalar));CHKERRQ(ierr); PetscFunctionReturn(0); }