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 Mat seqmat,mpimat; 7c4762a1bSJed Brown PetscMPIInt rank; 8c4762a1bSJed Brown PetscScalar value[3],*vals; 9c4762a1bSJed Brown PetscInt i,col[3],n=5,bs=1; 10c4762a1bSJed Brown 11*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help)); 125f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 135f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL)); 14c4762a1bSJed Brown 15c4762a1bSJed Brown /* Create seqaij matrices of size (n+rank) by n */ 165f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_SELF,&seqmat)); 175f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(seqmat,(n+rank)*bs,PETSC_DECIDE,PETSC_DECIDE,n*bs)); 185f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(seqmat)); 195f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJSetPreallocation(seqmat,3*bs,NULL)); 205f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqBAIJSetPreallocation(seqmat,bs,3,NULL)); 21c4762a1bSJed Brown 22c4762a1bSJed Brown if (bs == 1) { 23c4762a1bSJed Brown value[0] = -1.0; value[1] = 2.0; value[2] = -1.0; 24c4762a1bSJed Brown for (i=1; i<n-1; i++) { 25c4762a1bSJed Brown col[0] = i-1; col[1] = i; col[2] = i+1; 265f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(seqmat,1,&i,3,col,value,INSERT_VALUES)); 27c4762a1bSJed Brown } 28c4762a1bSJed Brown i = n - 1; col[0] = n - 2; col[1] = n - 1; 295f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(seqmat,1,&i,2,col,value,INSERT_VALUES)); 30c4762a1bSJed Brown 31c4762a1bSJed Brown i = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0; 325f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(seqmat,1,&i,2,col,value,INSERT_VALUES)); 33c4762a1bSJed Brown } else { 34c4762a1bSJed Brown PetscInt *rows,*cols,j; 355f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc3(bs*bs,&vals,bs,&rows,bs,&cols)); 36c4762a1bSJed Brown /* diagonal blocks */ 37c4762a1bSJed Brown for (i=0; i<bs*bs; i++) vals[i] = 2.0; 38c4762a1bSJed Brown for (i=0; i<n*bs; i+=bs) { 39c4762a1bSJed Brown for (j=0; j<bs; j++) {rows[j] = i+j; cols[j] = i+j;} 405f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(seqmat,bs,rows,bs,cols,vals,INSERT_VALUES)); 41c4762a1bSJed Brown } 42c4762a1bSJed Brown /* off-diagonal blocks */ 43c4762a1bSJed Brown for (i=0; i<bs*bs; i++) vals[i] = -1.0; 44c4762a1bSJed Brown for (i=0; i<(n-1)*bs; i+=bs) { 45c4762a1bSJed Brown for (j=0; j<bs; j++) {rows[j] = i+j; cols[j] = i+bs+j;} 465f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(seqmat,bs,rows,bs,cols,vals,INSERT_VALUES)); 47c4762a1bSJed Brown } 48c4762a1bSJed Brown 495f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree3(vals,rows,cols)); 50c4762a1bSJed Brown } 515f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(seqmat,MAT_FINAL_ASSEMBLY)); 525f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(seqmat,MAT_FINAL_ASSEMBLY)); 53dd400576SPatrick Sanan if (rank == 0) { 545f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"[%d] seqmat:\n",rank)); 555f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(seqmat,PETSC_VIEWER_STDOUT_SELF)); 56c4762a1bSJed Brown } 57c4762a1bSJed Brown 585f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateMPIMatConcatenateSeqMat(PETSC_COMM_WORLD,seqmat,PETSC_DECIDE,MAT_INITIAL_MATRIX,&mpimat)); 595f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateMPIMatConcatenateSeqMat(PETSC_COMM_WORLD,seqmat,PETSC_DECIDE,MAT_REUSE_MATRIX,&mpimat)); 605f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(mpimat,PETSC_VIEWER_STDOUT_WORLD)); 61c4762a1bSJed Brown 625f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&seqmat)); 635f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&mpimat)); 64*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 65*b122ec5aSJacob Faibussowitsch return 0; 66c4762a1bSJed Brown } 67c4762a1bSJed Brown 68c4762a1bSJed Brown /*TEST 69c4762a1bSJed Brown 70c4762a1bSJed Brown test: 71c4762a1bSJed Brown nsize: 3 72c4762a1bSJed Brown 73c4762a1bSJed Brown test: 74c4762a1bSJed Brown suffix: 2 75c4762a1bSJed Brown nsize: 3 76c4762a1bSJed Brown args: -mat_type baij 77c4762a1bSJed Brown 78c4762a1bSJed Brown test: 79c4762a1bSJed Brown suffix: 3 80c4762a1bSJed Brown nsize: 3 81c4762a1bSJed Brown args: -mat_type baij -bs 2 82c4762a1bSJed Brown 83c4762a1bSJed Brown TEST*/ 84