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