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