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*327415f7SBarry Smith PetscFunctionBeginUser; 199566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 209566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); 219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 23c4762a1bSJed Brown n = 2*size; 24c4762a1bSJed Brown 259566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&C)); 269566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n)); 279566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C)); 289566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 29c4762a1bSJed Brown 30c4762a1bSJed Brown /* Create the matrix for the five point stencil, YET AGAIN */ 31c4762a1bSJed Brown for (i=0; i<m; i++) { 32c4762a1bSJed Brown for (j=2*rank; j<2*rank+2; j++) { 33c4762a1bSJed Brown v = -1.0; Ii = j + n*i; 349566063dSJacob Faibussowitsch if (i>0) {J = Ii - n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 359566063dSJacob Faibussowitsch if (i<m-1) {J = Ii + n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 369566063dSJacob Faibussowitsch if (j>0) {J = Ii - 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 379566063dSJacob Faibussowitsch if (j<n-1) {J = Ii + 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 389566063dSJacob Faibussowitsch v = 4.0; PetscCall(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES)); 39c4762a1bSJed Brown } 40c4762a1bSJed Brown } 41c4762a1bSJed Brown 42c4762a1bSJed Brown /* Add extra elements (to illustrate variants of MatGetInfo) */ 43c4762a1bSJed Brown Ii = n; J = n-2; v = 100.0; 449566063dSJacob Faibussowitsch PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES)); 45c4762a1bSJed Brown Ii = n-2; J = n; v = 100.0; 469566063dSJacob Faibussowitsch PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES)); 47c4762a1bSJed Brown 489566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 499566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 50c4762a1bSJed Brown 51c4762a1bSJed Brown /* Form vectors */ 529566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(C,&x,&y)); 539566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(x,&ldim)); 549566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(x,&low,&high)); 55c4762a1bSJed Brown for (i=0; i<ldim; i++) { 56c4762a1bSJed Brown iglobal = i + low; 57c4762a1bSJed Brown v = one*((PetscReal)i) + 100.0*rank; 589566063dSJacob Faibussowitsch PetscCall(VecSetValues(x,1,&iglobal,&v,INSERT_VALUES)); 59c4762a1bSJed Brown } 609566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(x)); 619566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(x)); 62c4762a1bSJed Brown 639566063dSJacob Faibussowitsch PetscCall(MatMult(C,x,y)); 64c4762a1bSJed Brown 659566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-view_info",&flg_info)); 66c4762a1bSJed Brown if (flg_info) { 679566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO)); 689566063dSJacob Faibussowitsch PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 69c4762a1bSJed Brown 709566063dSJacob Faibussowitsch PetscCall(MatGetInfo(C,MAT_GLOBAL_SUM,&info)); 719566063dSJacob Faibussowitsch PetscCall(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)); 729566063dSJacob Faibussowitsch PetscCall(MatGetInfo (C,MAT_GLOBAL_MAX,&info)); 739566063dSJacob Faibussowitsch PetscCall(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)); 74c4762a1bSJed Brown } 75c4762a1bSJed Brown 769566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-view_mat",&flg_mat)); 771baa6e33SBarry Smith if (flg_mat) PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 78c4762a1bSJed Brown 79c4762a1bSJed Brown /* Test MatCreateRedundantMatrix() */ 80c4762a1bSJed Brown nsubcomms = size; 819566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-nsubcomms",&nsubcomms,NULL)); 829566063dSJacob Faibussowitsch PetscCall(MatCreateRedundantMatrix(C,nsubcomms,MPI_COMM_NULL,MAT_INITIAL_MATRIX,&Credundant)); 839566063dSJacob Faibussowitsch PetscCall(MatCreateRedundantMatrix(C,nsubcomms,MPI_COMM_NULL,MAT_REUSE_MATRIX,&Credundant)); 84c4762a1bSJed Brown 859566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)Credundant,&subcomm)); 869566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(subcomm,&subsize)); 87c4762a1bSJed Brown 88c4762a1bSJed Brown if (subsize==2 && flg_mat) { 899566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_(subcomm),"\n[%d] Credundant:\n",rank)); 909566063dSJacob Faibussowitsch PetscCall(MatView(Credundant,PETSC_VIEWER_STDOUT_(subcomm))); 91c4762a1bSJed Brown } 929566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Credundant)); 93c4762a1bSJed Brown 94c4762a1bSJed Brown /* Test MatCreateRedundantMatrix() with user-provided subcomm */ 95c4762a1bSJed Brown { 96c4762a1bSJed Brown PetscSubcomm psubcomm; 97c4762a1bSJed Brown 989566063dSJacob Faibussowitsch PetscCall(PetscSubcommCreate(PETSC_COMM_WORLD,&psubcomm)); 999566063dSJacob Faibussowitsch PetscCall(PetscSubcommSetNumber(psubcomm,nsubcomms)); 1009566063dSJacob Faibussowitsch PetscCall(PetscSubcommSetType(psubcomm,PETSC_SUBCOMM_CONTIGUOUS)); 101c4762a1bSJed Brown /* enable runtime switch of psubcomm type, e.g., '-psubcomm_type interlaced */ 1029566063dSJacob Faibussowitsch PetscCall(PetscSubcommSetFromOptions(psubcomm)); 103c4762a1bSJed Brown 1049566063dSJacob Faibussowitsch PetscCall(MatCreateRedundantMatrix(C,nsubcomms,PetscSubcommChild(psubcomm),MAT_INITIAL_MATRIX,&Credundant)); 1059566063dSJacob Faibussowitsch PetscCall(MatCreateRedundantMatrix(C,nsubcomms,PetscSubcommChild(psubcomm),MAT_REUSE_MATRIX,&Credundant)); 106c4762a1bSJed Brown 1079566063dSJacob Faibussowitsch PetscCall(PetscSubcommDestroy(&psubcomm)); 1089566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Credundant)); 109c4762a1bSJed Brown } 110c4762a1bSJed Brown 1119566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1129566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 1139566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 1149566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 115b122ec5aSJacob Faibussowitsch return 0; 116c4762a1bSJed Brown } 117c4762a1bSJed Brown 118c4762a1bSJed Brown /*TEST 119c4762a1bSJed Brown 120c4762a1bSJed Brown test: 121c4762a1bSJed Brown nsize: 3 122c4762a1bSJed Brown args: -view_info 123c4762a1bSJed Brown 124c4762a1bSJed Brown test: 125c4762a1bSJed Brown suffix: 2 126c4762a1bSJed Brown nsize: 3 127c4762a1bSJed Brown args: -nsubcomms 2 -view_mat -psubcomm_type interlaced 128c4762a1bSJed Brown 129c4762a1bSJed Brown test: 130c4762a1bSJed Brown suffix: 3 131c4762a1bSJed Brown nsize: 3 132c4762a1bSJed Brown args: -nsubcomms 2 -view_mat -psubcomm_type contiguous 133c4762a1bSJed Brown 134c4762a1bSJed Brown test: 135c4762a1bSJed Brown suffix: 3_baij 136c4762a1bSJed Brown nsize: 3 137c4762a1bSJed Brown args: -mat_type baij -nsubcomms 2 -view_mat 138c4762a1bSJed Brown 139c4762a1bSJed Brown test: 140c4762a1bSJed Brown suffix: 3_sbaij 141c4762a1bSJed Brown nsize: 3 142c4762a1bSJed Brown args: -mat_type sbaij -nsubcomms 2 -view_mat 143c4762a1bSJed Brown 144c4762a1bSJed Brown test: 145c4762a1bSJed Brown suffix: 3_dense 146c4762a1bSJed Brown nsize: 3 147c4762a1bSJed Brown args: -mat_type dense -nsubcomms 2 -view_mat 148c4762a1bSJed Brown 149c4762a1bSJed Brown test: 150c4762a1bSJed Brown suffix: 4_baij 151c4762a1bSJed Brown nsize: 3 152c4762a1bSJed Brown args: -mat_type baij -nsubcomms 2 -view_mat -psubcomm_type interlaced 153c4762a1bSJed Brown 154c4762a1bSJed Brown test: 155c4762a1bSJed Brown suffix: 4_sbaij 156c4762a1bSJed Brown nsize: 3 157c4762a1bSJed Brown args: -mat_type sbaij -nsubcomms 2 -view_mat -psubcomm_type interlaced 158c4762a1bSJed Brown 159c4762a1bSJed Brown test: 160c4762a1bSJed Brown suffix: 4_dense 161c4762a1bSJed Brown nsize: 3 162c4762a1bSJed Brown args: -mat_type dense -nsubcomms 2 -view_mat -psubcomm_type interlaced 163c4762a1bSJed Brown 164c4762a1bSJed Brown TEST*/ 165