static char help[] = "Tests cholesky, icc factorization and solve on sequential aij, baij and sbaij matrices. \n"; #include int main(int argc,char **args) { Vec x,y,b; Mat A; /* linear system matrix */ Mat sA,sC; /* symmetric part of the matrices */ PetscInt n,mbs=16,bs=1,nz=3,prob=1,i,j,col[3],block, row,Ii,J,n1,lvl; PetscErrorCode ierr; PetscMPIInt size; PetscReal norm2; PetscScalar neg_one = -1.0,four=4.0,value[3]; IS perm,cperm; PetscRandom rdm; PetscBool reorder = PETSC_FALSE,displ = PETSC_FALSE; MatFactorInfo factinfo; PetscBool equal; PetscBool TestAIJ = PETSC_FALSE,TestBAIJ = PETSC_TRUE; PetscInt TestShift=0; ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); ierr = PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL);CHKERRQ(ierr); ierr = PetscOptionsGetInt(NULL,NULL,"-mbs",&mbs,NULL);CHKERRQ(ierr); ierr = PetscOptionsGetBool(NULL,NULL,"-reorder",&reorder,NULL);CHKERRQ(ierr); ierr = PetscOptionsGetBool(NULL,NULL,"-testaij",&TestAIJ,NULL);CHKERRQ(ierr); ierr = PetscOptionsGetInt(NULL,NULL,"-testShift",&TestShift,NULL);CHKERRQ(ierr); ierr = PetscOptionsGetBool(NULL,NULL,"-displ",&displ,NULL);CHKERRQ(ierr); n = mbs*bs; if (TestAIJ) { /* A is in aij format */ ierr = MatCreateSeqAIJ(PETSC_COMM_WORLD,n,n,nz,NULL,&A);CHKERRQ(ierr); TestBAIJ = PETSC_FALSE; } else { /* A is in baij format */ ierr =MatCreateSeqBAIJ(PETSC_COMM_WORLD,bs,n,n,nz,NULL,&A);CHKERRQ(ierr); TestAIJ = PETSC_FALSE; } /* Assemble matrix */ if (bs == 1) { ierr = PetscOptionsGetInt(NULL,NULL,"-test_problem",&prob,NULL);CHKERRQ(ierr); if (prob == 1) { /* tridiagonal matrix */ value[0] = -1.0; value[1] = 2.0; value[2] = -1.0; for (i=1; i0) { J = Ii - n1; ierr = MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES);CHKERRQ(ierr); } if (i0) { J = Ii - 1; ierr = MatSetValues(A,1,&Ii,1,&J,&neg_one,INSERT_VALUES);CHKERRQ(ierr); } if (j 1 */ for (block=0; block