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 PetscErrorCode ierr; 14c4762a1bSJed Brown PetscBool flg_info,flg_mat; 15c4762a1bSJed Brown PetscScalar v,one = 1.0; 16c4762a1bSJed Brown Vec x,y; 17c4762a1bSJed Brown MPI_Comm subcomm; 18c4762a1bSJed Brown 19c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 20*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); 21*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 22*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 23c4762a1bSJed Brown n = 2*size; 24c4762a1bSJed Brown 25*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C)); 26*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n)); 27*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(C)); 28*5f80ce2aSJacob Faibussowitsch CHKERRQ(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; 34*5f80ce2aSJacob Faibussowitsch if (i>0) {J = Ii - n; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 35*5f80ce2aSJacob Faibussowitsch if (i<m-1) {J = Ii + n; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 36*5f80ce2aSJacob Faibussowitsch if (j>0) {J = Ii - 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 37*5f80ce2aSJacob Faibussowitsch if (j<n-1) {J = Ii + 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 38*5f80ce2aSJacob Faibussowitsch v = 4.0; CHKERRQ(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; 44*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES)); 45c4762a1bSJed Brown Ii = n-2; J = n; v = 100.0; 46*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES)); 47c4762a1bSJed Brown 48*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 49*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 50c4762a1bSJed Brown 51c4762a1bSJed Brown /* Form vectors */ 52*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(C,&x,&y)); 53*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(x,&ldim)); 54*5f80ce2aSJacob Faibussowitsch CHKERRQ(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; 58*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetValues(x,1,&iglobal,&v,INSERT_VALUES)); 59c4762a1bSJed Brown } 60*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyBegin(x)); 61*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyEnd(x)); 62c4762a1bSJed Brown 63*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(C,x,y)); 64c4762a1bSJed Brown 65*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-view_info",&flg_info)); 66c4762a1bSJed Brown if (flg_info) { 67*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO)); 68*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 69c4762a1bSJed Brown 70*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetInfo(C,MAT_GLOBAL_SUM,&info)); 71*5f80ce2aSJacob 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)); 72*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetInfo (C,MAT_GLOBAL_MAX,&info)); 73*5f80ce2aSJacob 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)); 74c4762a1bSJed Brown } 75c4762a1bSJed Brown 76*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-view_mat",&flg_mat)); 77c4762a1bSJed Brown if (flg_mat) { 78*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 79c4762a1bSJed Brown } 80c4762a1bSJed Brown 81c4762a1bSJed Brown /* Test MatCreateRedundantMatrix() */ 82c4762a1bSJed Brown nsubcomms = size; 83*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-nsubcomms",&nsubcomms,NULL)); 84*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateRedundantMatrix(C,nsubcomms,MPI_COMM_NULL,MAT_INITIAL_MATRIX,&Credundant)); 85*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateRedundantMatrix(C,nsubcomms,MPI_COMM_NULL,MAT_REUSE_MATRIX,&Credundant)); 86c4762a1bSJed Brown 87*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetComm((PetscObject)Credundant,&subcomm)); 88*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(subcomm,&subsize)); 89c4762a1bSJed Brown 90c4762a1bSJed Brown if (subsize==2 && flg_mat) { 91*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_(subcomm),"\n[%d] Credundant:\n",rank)); 92*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(Credundant,PETSC_VIEWER_STDOUT_(subcomm))); 93c4762a1bSJed Brown } 94*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Credundant)); 95c4762a1bSJed Brown 96c4762a1bSJed Brown /* Test MatCreateRedundantMatrix() with user-provided subcomm */ 97c4762a1bSJed Brown { 98c4762a1bSJed Brown PetscSubcomm psubcomm; 99c4762a1bSJed Brown 100*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSubcommCreate(PETSC_COMM_WORLD,&psubcomm)); 101*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSubcommSetNumber(psubcomm,nsubcomms)); 102*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSubcommSetType(psubcomm,PETSC_SUBCOMM_CONTIGUOUS)); 103c4762a1bSJed Brown /* enable runtime switch of psubcomm type, e.g., '-psubcomm_type interlaced */ 104*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSubcommSetFromOptions(psubcomm)); 105c4762a1bSJed Brown 106*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateRedundantMatrix(C,nsubcomms,PetscSubcommChild(psubcomm),MAT_INITIAL_MATRIX,&Credundant)); 107*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateRedundantMatrix(C,nsubcomms,PetscSubcommChild(psubcomm),MAT_REUSE_MATRIX,&Credundant)); 108c4762a1bSJed Brown 109*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSubcommDestroy(&psubcomm)); 110*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Credundant)); 111c4762a1bSJed Brown } 112c4762a1bSJed Brown 113*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 114*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&y)); 115*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 116c4762a1bSJed Brown ierr = PetscFinalize(); 117c4762a1bSJed Brown return ierr; 118c4762a1bSJed Brown } 119c4762a1bSJed Brown 120c4762a1bSJed Brown /*TEST 121c4762a1bSJed Brown 122c4762a1bSJed Brown test: 123c4762a1bSJed Brown nsize: 3 124c4762a1bSJed Brown args: -view_info 125c4762a1bSJed Brown 126c4762a1bSJed Brown test: 127c4762a1bSJed Brown suffix: 2 128c4762a1bSJed Brown nsize: 3 129c4762a1bSJed Brown args: -nsubcomms 2 -view_mat -psubcomm_type interlaced 130c4762a1bSJed Brown 131c4762a1bSJed Brown test: 132c4762a1bSJed Brown suffix: 3 133c4762a1bSJed Brown nsize: 3 134c4762a1bSJed Brown args: -nsubcomms 2 -view_mat -psubcomm_type contiguous 135c4762a1bSJed Brown 136c4762a1bSJed Brown test: 137c4762a1bSJed Brown suffix: 3_baij 138c4762a1bSJed Brown nsize: 3 139c4762a1bSJed Brown args: -mat_type baij -nsubcomms 2 -view_mat 140c4762a1bSJed Brown 141c4762a1bSJed Brown test: 142c4762a1bSJed Brown suffix: 3_sbaij 143c4762a1bSJed Brown nsize: 3 144c4762a1bSJed Brown args: -mat_type sbaij -nsubcomms 2 -view_mat 145c4762a1bSJed Brown 146c4762a1bSJed Brown test: 147c4762a1bSJed Brown suffix: 3_dense 148c4762a1bSJed Brown nsize: 3 149c4762a1bSJed Brown args: -mat_type dense -nsubcomms 2 -view_mat 150c4762a1bSJed Brown 151c4762a1bSJed Brown test: 152c4762a1bSJed Brown suffix: 4_baij 153c4762a1bSJed Brown nsize: 3 154c4762a1bSJed Brown args: -mat_type baij -nsubcomms 2 -view_mat -psubcomm_type interlaced 155c4762a1bSJed Brown 156c4762a1bSJed Brown test: 157c4762a1bSJed Brown suffix: 4_sbaij 158c4762a1bSJed Brown nsize: 3 159c4762a1bSJed Brown args: -mat_type sbaij -nsubcomms 2 -view_mat -psubcomm_type interlaced 160c4762a1bSJed Brown 161c4762a1bSJed Brown test: 162c4762a1bSJed Brown suffix: 4_dense 163c4762a1bSJed Brown nsize: 3 164c4762a1bSJed Brown args: -mat_type dense -nsubcomms 2 -view_mat -psubcomm_type interlaced 165c4762a1bSJed Brown 166c4762a1bSJed Brown TEST*/ 167