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*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 14*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL)); 15c4762a1bSJed Brown 16c4762a1bSJed Brown /* Create seqaij matrices of size (n+rank) by n */ 17*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_SELF,&seqmat)); 18*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(seqmat,(n+rank)*bs,PETSC_DECIDE,PETSC_DECIDE,n*bs)); 19*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(seqmat)); 20*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJSetPreallocation(seqmat,3*bs,NULL)); 21*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqBAIJSetPreallocation(seqmat,bs,3,NULL)); 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; 27*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(seqmat,1,&i,3,col,value,INSERT_VALUES)); 28c4762a1bSJed Brown } 29c4762a1bSJed Brown i = n - 1; col[0] = n - 2; col[1] = n - 1; 30*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(seqmat,1,&i,2,col,value,INSERT_VALUES)); 31c4762a1bSJed Brown 32c4762a1bSJed Brown i = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0; 33*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(seqmat,1,&i,2,col,value,INSERT_VALUES)); 34c4762a1bSJed Brown } else { 35c4762a1bSJed Brown PetscInt *rows,*cols,j; 36*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc3(bs*bs,&vals,bs,&rows,bs,&cols)); 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;} 41*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(seqmat,bs,rows,bs,cols,vals,INSERT_VALUES)); 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;} 47*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(seqmat,bs,rows,bs,cols,vals,INSERT_VALUES)); 48c4762a1bSJed Brown } 49c4762a1bSJed Brown 50*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree3(vals,rows,cols)); 51c4762a1bSJed Brown } 52*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(seqmat,MAT_FINAL_ASSEMBLY)); 53*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(seqmat,MAT_FINAL_ASSEMBLY)); 54dd400576SPatrick Sanan if (rank == 0) { 55*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"[%d] seqmat:\n",rank)); 56*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(seqmat,PETSC_VIEWER_STDOUT_SELF)); 57c4762a1bSJed Brown } 58c4762a1bSJed Brown 59*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateMPIMatConcatenateSeqMat(PETSC_COMM_WORLD,seqmat,PETSC_DECIDE,MAT_INITIAL_MATRIX,&mpimat)); 60*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateMPIMatConcatenateSeqMat(PETSC_COMM_WORLD,seqmat,PETSC_DECIDE,MAT_REUSE_MATRIX,&mpimat)); 61*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(mpimat,PETSC_VIEWER_STDOUT_WORLD)); 62c4762a1bSJed Brown 63*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&seqmat)); 64*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&mpimat)); 65c4762a1bSJed Brown ierr = PetscFinalize(); 66c4762a1bSJed Brown return ierr; 67c4762a1bSJed Brown } 68c4762a1bSJed Brown 69c4762a1bSJed Brown /*TEST 70c4762a1bSJed Brown 71c4762a1bSJed Brown test: 72c4762a1bSJed Brown nsize: 3 73c4762a1bSJed Brown 74c4762a1bSJed Brown test: 75c4762a1bSJed Brown suffix: 2 76c4762a1bSJed Brown nsize: 3 77c4762a1bSJed Brown args: -mat_type baij 78c4762a1bSJed Brown 79c4762a1bSJed Brown test: 80c4762a1bSJed Brown suffix: 3 81c4762a1bSJed Brown nsize: 3 82c4762a1bSJed Brown args: -mat_type baij -bs 2 83c4762a1bSJed Brown 84c4762a1bSJed Brown TEST*/ 85