1 static char help[] = "Testing MatCreateMPIAIJSumSeqAIJ().\n\n"; 2 3 #include <petscmat.h> 4 5 int main(int argc,char **argv) 6 { 7 Mat A,B; 8 MatScalar a[1],alpha; 9 PetscMPIInt size,rank; 10 PetscInt m,n,i,col, prid; 11 PetscErrorCode ierr; 12 13 ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 14 ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 15 ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 16 prid = size; 17 ierr = PetscOptionsGetInt(NULL,NULL,"-prid",&prid,NULL);CHKERRQ(ierr); 18 19 m = n = 10*size; 20 ierr = MatCreate(PETSC_COMM_SELF,&A);CHKERRQ(ierr); 21 ierr = MatSetSizes(A,PETSC_DETERMINE,PETSC_DETERMINE,m,n);CHKERRQ(ierr); 22 ierr = MatSetType(A,MATSEQAIJ);CHKERRQ(ierr); 23 ierr = MatSetUp(A);CHKERRQ(ierr); 24 25 a[0] = rank+1; 26 for (i=0; i<m-rank; i++) { 27 col = i+rank; 28 ierr = MatSetValues(A,1,&i,1,&col,a,INSERT_VALUES);CHKERRQ(ierr); 29 } 30 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 31 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 32 33 if (rank == prid) { 34 ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] A: \n",rank);CHKERRQ(ierr); 35 ierr = MatView(A,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 36 } 37 38 /* Test MatCreateMPIAIJSumSeqAIJ */ 39 ierr = MatCreateMPIAIJSumSeqAIJ(PETSC_COMM_WORLD,A,PETSC_DECIDE,PETSC_DECIDE,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr); 40 41 /* Test MAT_REUSE_MATRIX */ 42 alpha = 0.1; 43 for (i=0; i<3; i++) { 44 ierr = MatScale(A,alpha);CHKERRQ(ierr); 45 ierr = MatCreateMPIAIJSumSeqAIJ(PETSC_COMM_WORLD,A,PETSC_DECIDE,PETSC_DECIDE,MAT_REUSE_MATRIX,&B);CHKERRQ(ierr); 46 } 47 ierr = MatView(B, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 48 ierr = MatDestroy(&B);CHKERRQ(ierr); 49 ierr = MatDestroy(&A);CHKERRQ(ierr); 50 ierr = PetscFinalize(); 51 return ierr; 52 } 53 54 55 /*TEST 56 57 test: 58 nsize: 3 59 filter: grep -v "MPI processes" 60 61 test: 62 suffix: 2 63 filter: grep -v "MPI processes" 64 65 TEST*/ 66