xref: /petsc/src/mat/tests/ex9.c (revision 1baa6e3354bfe224b33a0c158f482508792a8a6e)
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 
189566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
199566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
209566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
219566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
22c4762a1bSJed Brown   n    = 2*size;
23c4762a1bSJed Brown 
249566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
259566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n));
269566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(C));
279566063dSJacob Faibussowitsch   PetscCall(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;
339566063dSJacob Faibussowitsch       if (i>0)   {J = Ii - n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
349566063dSJacob Faibussowitsch       if (i<m-1) {J = Ii + n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
359566063dSJacob Faibussowitsch       if (j>0)   {J = Ii - 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
369566063dSJacob Faibussowitsch       if (j<n-1) {J = Ii + 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
379566063dSJacob Faibussowitsch       v = 4.0; PetscCall(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;
439566063dSJacob Faibussowitsch   PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));
44c4762a1bSJed Brown   Ii   = n-2; J = n; v = 100.0;
459566063dSJacob Faibussowitsch   PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));
46c4762a1bSJed Brown 
479566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
489566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
49c4762a1bSJed Brown 
50c4762a1bSJed Brown   /* Form vectors */
519566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(C,&x,&y));
529566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(x,&ldim));
539566063dSJacob Faibussowitsch   PetscCall(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;
579566063dSJacob Faibussowitsch     PetscCall(VecSetValues(x,1,&iglobal,&v,INSERT_VALUES));
58c4762a1bSJed Brown   }
599566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(x));
609566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(x));
61c4762a1bSJed Brown 
629566063dSJacob Faibussowitsch   PetscCall(MatMult(C,x,y));
63c4762a1bSJed Brown 
649566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL,"-view_info",&flg_info));
65c4762a1bSJed Brown   if (flg_info)  {
669566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO));
679566063dSJacob Faibussowitsch     PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
68c4762a1bSJed Brown 
699566063dSJacob Faibussowitsch     PetscCall(MatGetInfo(C,MAT_GLOBAL_SUM,&info));
709566063dSJacob 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));
719566063dSJacob Faibussowitsch     PetscCall(MatGetInfo (C,MAT_GLOBAL_MAX,&info));
729566063dSJacob 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));
73c4762a1bSJed Brown   }
74c4762a1bSJed Brown 
759566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL,"-view_mat",&flg_mat));
76*1baa6e33SBarry Smith   if (flg_mat) PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
77c4762a1bSJed Brown 
78c4762a1bSJed Brown   /* Test MatCreateRedundantMatrix() */
79c4762a1bSJed Brown   nsubcomms = size;
809566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-nsubcomms",&nsubcomms,NULL));
819566063dSJacob Faibussowitsch   PetscCall(MatCreateRedundantMatrix(C,nsubcomms,MPI_COMM_NULL,MAT_INITIAL_MATRIX,&Credundant));
829566063dSJacob Faibussowitsch   PetscCall(MatCreateRedundantMatrix(C,nsubcomms,MPI_COMM_NULL,MAT_REUSE_MATRIX,&Credundant));
83c4762a1bSJed Brown 
849566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Credundant,&subcomm));
859566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(subcomm,&subsize));
86c4762a1bSJed Brown 
87c4762a1bSJed Brown   if (subsize==2 && flg_mat) {
889566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_(subcomm),"\n[%d] Credundant:\n",rank));
899566063dSJacob Faibussowitsch     PetscCall(MatView(Credundant,PETSC_VIEWER_STDOUT_(subcomm)));
90c4762a1bSJed Brown   }
919566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Credundant));
92c4762a1bSJed Brown 
93c4762a1bSJed Brown   /* Test MatCreateRedundantMatrix() with user-provided subcomm */
94c4762a1bSJed Brown   {
95c4762a1bSJed Brown     PetscSubcomm psubcomm;
96c4762a1bSJed Brown 
979566063dSJacob Faibussowitsch     PetscCall(PetscSubcommCreate(PETSC_COMM_WORLD,&psubcomm));
989566063dSJacob Faibussowitsch     PetscCall(PetscSubcommSetNumber(psubcomm,nsubcomms));
999566063dSJacob Faibussowitsch     PetscCall(PetscSubcommSetType(psubcomm,PETSC_SUBCOMM_CONTIGUOUS));
100c4762a1bSJed Brown     /* enable runtime switch of psubcomm type, e.g., '-psubcomm_type interlaced */
1019566063dSJacob Faibussowitsch     PetscCall(PetscSubcommSetFromOptions(psubcomm));
102c4762a1bSJed Brown 
1039566063dSJacob Faibussowitsch     PetscCall(MatCreateRedundantMatrix(C,nsubcomms,PetscSubcommChild(psubcomm),MAT_INITIAL_MATRIX,&Credundant));
1049566063dSJacob Faibussowitsch     PetscCall(MatCreateRedundantMatrix(C,nsubcomms,PetscSubcommChild(psubcomm),MAT_REUSE_MATRIX,&Credundant));
105c4762a1bSJed Brown 
1069566063dSJacob Faibussowitsch     PetscCall(PetscSubcommDestroy(&psubcomm));
1079566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Credundant));
108c4762a1bSJed Brown   }
109c4762a1bSJed Brown 
1109566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1119566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
1129566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
1139566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
114b122ec5aSJacob Faibussowitsch   return 0;
115c4762a1bSJed Brown }
116c4762a1bSJed Brown 
117c4762a1bSJed Brown /*TEST
118c4762a1bSJed Brown 
119c4762a1bSJed Brown    test:
120c4762a1bSJed Brown       nsize: 3
121c4762a1bSJed Brown       args: -view_info
122c4762a1bSJed Brown 
123c4762a1bSJed Brown    test:
124c4762a1bSJed Brown       suffix: 2
125c4762a1bSJed Brown       nsize: 3
126c4762a1bSJed Brown       args: -nsubcomms 2 -view_mat -psubcomm_type interlaced
127c4762a1bSJed Brown 
128c4762a1bSJed Brown    test:
129c4762a1bSJed Brown       suffix: 3
130c4762a1bSJed Brown       nsize: 3
131c4762a1bSJed Brown       args: -nsubcomms 2 -view_mat -psubcomm_type contiguous
132c4762a1bSJed Brown 
133c4762a1bSJed Brown    test:
134c4762a1bSJed Brown       suffix: 3_baij
135c4762a1bSJed Brown       nsize: 3
136c4762a1bSJed Brown       args: -mat_type baij -nsubcomms 2 -view_mat
137c4762a1bSJed Brown 
138c4762a1bSJed Brown    test:
139c4762a1bSJed Brown       suffix: 3_sbaij
140c4762a1bSJed Brown       nsize: 3
141c4762a1bSJed Brown       args: -mat_type sbaij -nsubcomms 2 -view_mat
142c4762a1bSJed Brown 
143c4762a1bSJed Brown    test:
144c4762a1bSJed Brown       suffix: 3_dense
145c4762a1bSJed Brown       nsize: 3
146c4762a1bSJed Brown       args: -mat_type dense -nsubcomms 2 -view_mat
147c4762a1bSJed Brown 
148c4762a1bSJed Brown    test:
149c4762a1bSJed Brown       suffix: 4_baij
150c4762a1bSJed Brown       nsize: 3
151c4762a1bSJed Brown       args: -mat_type baij -nsubcomms 2 -view_mat -psubcomm_type interlaced
152c4762a1bSJed Brown 
153c4762a1bSJed Brown    test:
154c4762a1bSJed Brown       suffix: 4_sbaij
155c4762a1bSJed Brown       nsize: 3
156c4762a1bSJed Brown       args: -mat_type sbaij -nsubcomms 2 -view_mat -psubcomm_type interlaced
157c4762a1bSJed Brown 
158c4762a1bSJed Brown    test:
159c4762a1bSJed Brown       suffix: 4_dense
160c4762a1bSJed Brown       nsize: 3
161c4762a1bSJed Brown       args: -mat_type dense -nsubcomms 2 -view_mat -psubcomm_type interlaced
162c4762a1bSJed Brown 
163c4762a1bSJed Brown TEST*/
164