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