1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests MPI parallel matrix creation. Test MatCreateRedundantMatrix() \n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown int main(int argc,char **args) 7c4762a1bSJed Brown { 8c4762a1bSJed Brown Mat C,Credundant; 9c4762a1bSJed Brown MatInfo info; 10c4762a1bSJed Brown PetscMPIInt rank,size,subsize; 11c4762a1bSJed Brown PetscInt i,j,m = 3,n = 2,low,high,iglobal; 12c4762a1bSJed Brown PetscInt Ii,J,ldim,nsubcomms; 13c4762a1bSJed Brown PetscBool flg_info,flg_mat; 14c4762a1bSJed Brown PetscScalar v,one = 1.0; 15c4762a1bSJed Brown Vec x,y; 16c4762a1bSJed Brown MPI_Comm subcomm; 17c4762a1bSJed Brown 18*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 195f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); 205f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 215f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 22c4762a1bSJed Brown n = 2*size; 23c4762a1bSJed Brown 245f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C)); 255f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n)); 265f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(C)); 275f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(C)); 28c4762a1bSJed Brown 29c4762a1bSJed Brown /* Create the matrix for the five point stencil, YET AGAIN */ 30c4762a1bSJed Brown for (i=0; i<m; i++) { 31c4762a1bSJed Brown for (j=2*rank; j<2*rank+2; j++) { 32c4762a1bSJed Brown v = -1.0; Ii = j + n*i; 335f80ce2aSJacob Faibussowitsch if (i>0) {J = Ii - n; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 345f80ce2aSJacob Faibussowitsch if (i<m-1) {J = Ii + n; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 355f80ce2aSJacob Faibussowitsch if (j>0) {J = Ii - 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 365f80ce2aSJacob Faibussowitsch if (j<n-1) {J = Ii + 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 375f80ce2aSJacob Faibussowitsch v = 4.0; CHKERRQ(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES)); 38c4762a1bSJed Brown } 39c4762a1bSJed Brown } 40c4762a1bSJed Brown 41c4762a1bSJed Brown /* Add extra elements (to illustrate variants of MatGetInfo) */ 42c4762a1bSJed Brown Ii = n; J = n-2; v = 100.0; 435f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES)); 44c4762a1bSJed Brown Ii = n-2; J = n; v = 100.0; 455f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES)); 46c4762a1bSJed Brown 475f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 485f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 49c4762a1bSJed Brown 50c4762a1bSJed Brown /* Form vectors */ 515f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(C,&x,&y)); 525f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(x,&ldim)); 535f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetOwnershipRange(x,&low,&high)); 54c4762a1bSJed Brown for (i=0; i<ldim; i++) { 55c4762a1bSJed Brown iglobal = i + low; 56c4762a1bSJed Brown v = one*((PetscReal)i) + 100.0*rank; 575f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetValues(x,1,&iglobal,&v,INSERT_VALUES)); 58c4762a1bSJed Brown } 595f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyBegin(x)); 605f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyEnd(x)); 61c4762a1bSJed Brown 625f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(C,x,y)); 63c4762a1bSJed Brown 645f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-view_info",&flg_info)); 65c4762a1bSJed Brown if (flg_info) { 665f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO)); 675f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 68c4762a1bSJed Brown 695f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetInfo(C,MAT_GLOBAL_SUM,&info)); 705f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"matrix information (global sums):\nnonzeros = %" PetscInt_FMT ", allocated nonzeros = %" PetscInt_FMT "\n",(PetscInt)info.nz_used,(PetscInt)info.nz_allocated)); 715f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetInfo (C,MAT_GLOBAL_MAX,&info)); 725f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"matrix information (global max):\nnonzeros = %" PetscInt_FMT ", allocated nonzeros = %" PetscInt_FMT "\n",(PetscInt)info.nz_used,(PetscInt)info.nz_allocated)); 73c4762a1bSJed Brown } 74c4762a1bSJed Brown 755f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-view_mat",&flg_mat)); 76c4762a1bSJed Brown if (flg_mat) { 775f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 78c4762a1bSJed Brown } 79c4762a1bSJed Brown 80c4762a1bSJed Brown /* Test MatCreateRedundantMatrix() */ 81c4762a1bSJed Brown nsubcomms = size; 825f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-nsubcomms",&nsubcomms,NULL)); 835f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateRedundantMatrix(C,nsubcomms,MPI_COMM_NULL,MAT_INITIAL_MATRIX,&Credundant)); 845f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateRedundantMatrix(C,nsubcomms,MPI_COMM_NULL,MAT_REUSE_MATRIX,&Credundant)); 85c4762a1bSJed Brown 865f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetComm((PetscObject)Credundant,&subcomm)); 875f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(subcomm,&subsize)); 88c4762a1bSJed Brown 89c4762a1bSJed Brown if (subsize==2 && flg_mat) { 905f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_(subcomm),"\n[%d] Credundant:\n",rank)); 915f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(Credundant,PETSC_VIEWER_STDOUT_(subcomm))); 92c4762a1bSJed Brown } 935f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Credundant)); 94c4762a1bSJed Brown 95c4762a1bSJed Brown /* Test MatCreateRedundantMatrix() with user-provided subcomm */ 96c4762a1bSJed Brown { 97c4762a1bSJed Brown PetscSubcomm psubcomm; 98c4762a1bSJed Brown 995f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSubcommCreate(PETSC_COMM_WORLD,&psubcomm)); 1005f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSubcommSetNumber(psubcomm,nsubcomms)); 1015f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSubcommSetType(psubcomm,PETSC_SUBCOMM_CONTIGUOUS)); 102c4762a1bSJed Brown /* enable runtime switch of psubcomm type, e.g., '-psubcomm_type interlaced */ 1035f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSubcommSetFromOptions(psubcomm)); 104c4762a1bSJed Brown 1055f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateRedundantMatrix(C,nsubcomms,PetscSubcommChild(psubcomm),MAT_INITIAL_MATRIX,&Credundant)); 1065f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateRedundantMatrix(C,nsubcomms,PetscSubcommChild(psubcomm),MAT_REUSE_MATRIX,&Credundant)); 107c4762a1bSJed Brown 1085f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSubcommDestroy(&psubcomm)); 1095f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Credundant)); 110c4762a1bSJed Brown } 111c4762a1bSJed Brown 1125f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 1135f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&y)); 1145f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 115*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 116*b122ec5aSJacob Faibussowitsch return 0; 117c4762a1bSJed Brown } 118c4762a1bSJed Brown 119c4762a1bSJed Brown /*TEST 120c4762a1bSJed Brown 121c4762a1bSJed Brown test: 122c4762a1bSJed Brown nsize: 3 123c4762a1bSJed Brown args: -view_info 124c4762a1bSJed Brown 125c4762a1bSJed Brown test: 126c4762a1bSJed Brown suffix: 2 127c4762a1bSJed Brown nsize: 3 128c4762a1bSJed Brown args: -nsubcomms 2 -view_mat -psubcomm_type interlaced 129c4762a1bSJed Brown 130c4762a1bSJed Brown test: 131c4762a1bSJed Brown suffix: 3 132c4762a1bSJed Brown nsize: 3 133c4762a1bSJed Brown args: -nsubcomms 2 -view_mat -psubcomm_type contiguous 134c4762a1bSJed Brown 135c4762a1bSJed Brown test: 136c4762a1bSJed Brown suffix: 3_baij 137c4762a1bSJed Brown nsize: 3 138c4762a1bSJed Brown args: -mat_type baij -nsubcomms 2 -view_mat 139c4762a1bSJed Brown 140c4762a1bSJed Brown test: 141c4762a1bSJed Brown suffix: 3_sbaij 142c4762a1bSJed Brown nsize: 3 143c4762a1bSJed Brown args: -mat_type sbaij -nsubcomms 2 -view_mat 144c4762a1bSJed Brown 145c4762a1bSJed Brown test: 146c4762a1bSJed Brown suffix: 3_dense 147c4762a1bSJed Brown nsize: 3 148c4762a1bSJed Brown args: -mat_type dense -nsubcomms 2 -view_mat 149c4762a1bSJed Brown 150c4762a1bSJed Brown test: 151c4762a1bSJed Brown suffix: 4_baij 152c4762a1bSJed Brown nsize: 3 153c4762a1bSJed Brown args: -mat_type baij -nsubcomms 2 -view_mat -psubcomm_type interlaced 154c4762a1bSJed Brown 155c4762a1bSJed Brown test: 156c4762a1bSJed Brown suffix: 4_sbaij 157c4762a1bSJed Brown nsize: 3 158c4762a1bSJed Brown args: -mat_type sbaij -nsubcomms 2 -view_mat -psubcomm_type interlaced 159c4762a1bSJed Brown 160c4762a1bSJed Brown test: 161c4762a1bSJed Brown suffix: 4_dense 162c4762a1bSJed Brown nsize: 3 163c4762a1bSJed Brown args: -mat_type dense -nsubcomms 2 -view_mat -psubcomm_type interlaced 164c4762a1bSJed Brown 165c4762a1bSJed Brown TEST*/ 166