1*c4762a1bSJed Brown static char help[] = "Testing MatCreateMPIMatConcatenateSeqMat().\n\n"; 2*c4762a1bSJed Brown 3*c4762a1bSJed Brown #include <petscmat.h> 4*c4762a1bSJed Brown int main(int argc,char **argv) 5*c4762a1bSJed Brown { 6*c4762a1bSJed Brown PetscErrorCode ierr; 7*c4762a1bSJed Brown Mat seqmat,mpimat; 8*c4762a1bSJed Brown PetscMPIInt rank; 9*c4762a1bSJed Brown PetscScalar value[3],*vals; 10*c4762a1bSJed Brown PetscInt i,col[3],n=5,bs=1; 11*c4762a1bSJed Brown 12*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 13*c4762a1bSJed Brown ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 14*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL);CHKERRQ(ierr); 15*c4762a1bSJed Brown 16*c4762a1bSJed Brown /* Create seqaij matrices of size (n+rank) by n */ 17*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_SELF,&seqmat);CHKERRQ(ierr); 18*c4762a1bSJed Brown ierr = MatSetSizes(seqmat,(n+rank)*bs,PETSC_DECIDE,PETSC_DECIDE,n*bs);CHKERRQ(ierr); 19*c4762a1bSJed Brown ierr = MatSetFromOptions(seqmat);CHKERRQ(ierr); 20*c4762a1bSJed Brown ierr = MatSeqAIJSetPreallocation(seqmat,3*bs,NULL);CHKERRQ(ierr); 21*c4762a1bSJed Brown ierr = MatSeqBAIJSetPreallocation(seqmat,bs,3,NULL);CHKERRQ(ierr); 22*c4762a1bSJed Brown 23*c4762a1bSJed Brown if (bs == 1) { 24*c4762a1bSJed Brown value[0] = -1.0; value[1] = 2.0; value[2] = -1.0; 25*c4762a1bSJed Brown for (i=1; i<n-1; i++) { 26*c4762a1bSJed Brown col[0] = i-1; col[1] = i; col[2] = i+1; 27*c4762a1bSJed Brown ierr = MatSetValues(seqmat,1,&i,3,col,value,INSERT_VALUES);CHKERRQ(ierr); 28*c4762a1bSJed Brown } 29*c4762a1bSJed Brown i = n - 1; col[0] = n - 2; col[1] = n - 1; 30*c4762a1bSJed Brown ierr = MatSetValues(seqmat,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr); 31*c4762a1bSJed Brown 32*c4762a1bSJed Brown i = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0; 33*c4762a1bSJed Brown ierr = MatSetValues(seqmat,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr); 34*c4762a1bSJed Brown } else { 35*c4762a1bSJed Brown PetscInt *rows,*cols,j; 36*c4762a1bSJed Brown ierr = PetscMalloc3(bs*bs,&vals,bs,&rows,bs,&cols);CHKERRQ(ierr); 37*c4762a1bSJed Brown /* diagonal blocks */ 38*c4762a1bSJed Brown for (i=0; i<bs*bs; i++) vals[i] = 2.0; 39*c4762a1bSJed Brown for (i=0; i<n*bs; i+=bs) { 40*c4762a1bSJed Brown for (j=0; j<bs; j++) {rows[j] = i+j; cols[j] = i+j;} 41*c4762a1bSJed Brown ierr = MatSetValues(seqmat,bs,rows,bs,cols,vals,INSERT_VALUES);CHKERRQ(ierr); 42*c4762a1bSJed Brown } 43*c4762a1bSJed Brown /* off-diagonal blocks */ 44*c4762a1bSJed Brown for (i=0; i<bs*bs; i++) vals[i] = -1.0; 45*c4762a1bSJed Brown for (i=0; i<(n-1)*bs; i+=bs) { 46*c4762a1bSJed Brown for (j=0; j<bs; j++) {rows[j] = i+j; cols[j] = i+bs+j;} 47*c4762a1bSJed Brown ierr = MatSetValues(seqmat,bs,rows,bs,cols,vals,INSERT_VALUES);CHKERRQ(ierr); 48*c4762a1bSJed Brown } 49*c4762a1bSJed Brown 50*c4762a1bSJed Brown ierr = PetscFree3(vals,rows,cols);CHKERRQ(ierr); 51*c4762a1bSJed Brown } 52*c4762a1bSJed Brown ierr = MatAssemblyBegin(seqmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 53*c4762a1bSJed Brown ierr = MatAssemblyEnd(seqmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 54*c4762a1bSJed Brown if (!rank) { 55*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] seqmat:\n",rank);CHKERRQ(ierr); 56*c4762a1bSJed Brown ierr = MatView(seqmat,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 57*c4762a1bSJed Brown } 58*c4762a1bSJed Brown 59*c4762a1bSJed Brown ierr = MatCreateMPIMatConcatenateSeqMat(PETSC_COMM_WORLD,seqmat,PETSC_DECIDE,MAT_INITIAL_MATRIX,&mpimat);CHKERRQ(ierr); 60*c4762a1bSJed Brown ierr = MatCreateMPIMatConcatenateSeqMat(PETSC_COMM_WORLD,seqmat,PETSC_DECIDE,MAT_REUSE_MATRIX,&mpimat);CHKERRQ(ierr); 61*c4762a1bSJed Brown ierr = MatView(mpimat,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 62*c4762a1bSJed Brown 63*c4762a1bSJed Brown ierr = MatDestroy(&seqmat);CHKERRQ(ierr); 64*c4762a1bSJed Brown ierr = MatDestroy(&mpimat);CHKERRQ(ierr); 65*c4762a1bSJed Brown ierr = PetscFinalize(); 66*c4762a1bSJed Brown return ierr; 67*c4762a1bSJed Brown } 68*c4762a1bSJed Brown 69*c4762a1bSJed Brown 70*c4762a1bSJed Brown 71*c4762a1bSJed Brown /*TEST 72*c4762a1bSJed Brown 73*c4762a1bSJed Brown test: 74*c4762a1bSJed Brown nsize: 3 75*c4762a1bSJed Brown 76*c4762a1bSJed Brown test: 77*c4762a1bSJed Brown suffix: 2 78*c4762a1bSJed Brown nsize: 3 79*c4762a1bSJed Brown args: -mat_type baij 80*c4762a1bSJed Brown 81*c4762a1bSJed Brown test: 82*c4762a1bSJed Brown suffix: 3 83*c4762a1bSJed Brown nsize: 3 84*c4762a1bSJed Brown args: -mat_type baij -bs 2 85*c4762a1bSJed Brown 86*c4762a1bSJed Brown TEST*/ 87