1 2 static char help[] = "Tests using MatShift() to create a constant diagonal matrix\n\n"; 3 4 #include <petscmat.h> 5 6 int main(int argc,char **argv) 7 { 8 Mat A,F; 9 MatFactorInfo info; 10 PetscErrorCode ierr; 11 PetscInt m = 10; 12 IS perm; 13 PetscMPIInt size; 14 PetscBool issbaij; 15 16 ierr = PetscInitialize(&argc,&argv,(char*) 0,help);if (ierr) return ierr; 17 ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 18 19 ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 20 ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,m);CHKERRQ(ierr); 21 ierr = MatSetFromOptions(A);CHKERRQ(ierr); 22 ierr = MatSetUp(A);CHKERRQ(ierr); 23 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 24 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 25 26 ierr = MatShift(A,1.0);CHKERRQ(ierr); 27 28 ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr); 29 if (size == 1 && !issbaij) { 30 ierr = MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_LU,&F);CHKERRQ(ierr); 31 ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr); 32 ierr = ISCreateStride(PETSC_COMM_SELF,m,0,1,&perm);CHKERRQ(ierr); 33 ierr = MatLUFactorSymbolic(F,A,perm,perm,&info);CHKERRQ(ierr); 34 ierr = MatLUFactorNumeric(F,A,&info);CHKERRQ(ierr); 35 ierr = MatDestroy(&F);CHKERRQ(ierr); 36 ierr = ISDestroy(&perm);CHKERRQ(ierr); 37 } 38 ierr = MatDestroy(&A);CHKERRQ(ierr); 39 ierr = PetscFinalize(); 40 return ierr; 41 } 42 43 /*TEST 44 45 test: 46 requires: define(PETSC_USE_INFO) 47 args: -info 48 filter: grep malloc | sort -b 49 50 test: 51 suffix: 2 52 nsize: 2 53 requires: define(PETSC_USE_INFO) 54 args: -info ex182info 55 filter: grep -h malloc "ex182info.0" | sort -b 56 57 test: 58 suffix: 3 59 requires: define(PETSC_USE_INFO) 60 args: -info -mat_type baij 61 filter: grep malloc | sort -b 62 63 test: 64 suffix: 4 65 nsize: 2 66 requires: define(PETSC_USE_INFO) 67 args: -info ex182info -mat_type baij 68 filter: grep -h malloc "ex182info.1" | sort -b 69 70 test: 71 suffix: 5 72 requires: define(PETSC_USE_INFO) 73 args: -info -mat_type sbaij 74 filter: grep malloc | sort -b 75 76 test: 77 suffix: 6 78 nsize: 2 79 requires: define(PETSC_USE_INFO) 80 args: -info ex182info -mat_type sbaij 81 filter: grep -h malloc "ex182info.0" | sort -b 82 83 TEST*/ 84