xref: /petsc/src/mat/tests/ex9.c (revision 9566063d113dddea24716c546802770db7481bc0)
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*9566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
19*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
20*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
21*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
22c4762a1bSJed Brown   n    = 2*size;
23c4762a1bSJed Brown 
24*9566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
25*9566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n));
26*9566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(C));
27*9566063dSJacob 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;
33*9566063dSJacob Faibussowitsch       if (i>0)   {J = Ii - n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
34*9566063dSJacob Faibussowitsch       if (i<m-1) {J = Ii + n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
35*9566063dSJacob Faibussowitsch       if (j>0)   {J = Ii - 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
36*9566063dSJacob Faibussowitsch       if (j<n-1) {J = Ii + 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
37*9566063dSJacob 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;
43*9566063dSJacob Faibussowitsch   PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));
44c4762a1bSJed Brown   Ii   = n-2; J = n; v = 100.0;
45*9566063dSJacob Faibussowitsch   PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));
46c4762a1bSJed Brown 
47*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
48*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
49c4762a1bSJed Brown 
50c4762a1bSJed Brown   /* Form vectors */
51*9566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(C,&x,&y));
52*9566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(x,&ldim));
53*9566063dSJacob 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;
57*9566063dSJacob Faibussowitsch     PetscCall(VecSetValues(x,1,&iglobal,&v,INSERT_VALUES));
58c4762a1bSJed Brown   }
59*9566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(x));
60*9566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(x));
61c4762a1bSJed Brown 
62*9566063dSJacob Faibussowitsch   PetscCall(MatMult(C,x,y));
63c4762a1bSJed Brown 
64*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL,"-view_info",&flg_info));
65c4762a1bSJed Brown   if (flg_info)  {
66*9566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO));
67*9566063dSJacob Faibussowitsch     PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
68c4762a1bSJed Brown 
69*9566063dSJacob Faibussowitsch     PetscCall(MatGetInfo(C,MAT_GLOBAL_SUM,&info));
70*9566063dSJacob 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));
71*9566063dSJacob Faibussowitsch     PetscCall(MatGetInfo (C,MAT_GLOBAL_MAX,&info));
72*9566063dSJacob 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 
75*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL,"-view_mat",&flg_mat));
76c4762a1bSJed Brown   if (flg_mat) {
77*9566063dSJacob Faibussowitsch     PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
78c4762a1bSJed Brown   }
79c4762a1bSJed Brown 
80c4762a1bSJed Brown   /* Test MatCreateRedundantMatrix() */
81c4762a1bSJed Brown   nsubcomms = size;
82*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-nsubcomms",&nsubcomms,NULL));
83*9566063dSJacob Faibussowitsch   PetscCall(MatCreateRedundantMatrix(C,nsubcomms,MPI_COMM_NULL,MAT_INITIAL_MATRIX,&Credundant));
84*9566063dSJacob Faibussowitsch   PetscCall(MatCreateRedundantMatrix(C,nsubcomms,MPI_COMM_NULL,MAT_REUSE_MATRIX,&Credundant));
85c4762a1bSJed Brown 
86*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)Credundant,&subcomm));
87*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(subcomm,&subsize));
88c4762a1bSJed Brown 
89c4762a1bSJed Brown   if (subsize==2 && flg_mat) {
90*9566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_(subcomm),"\n[%d] Credundant:\n",rank));
91*9566063dSJacob Faibussowitsch     PetscCall(MatView(Credundant,PETSC_VIEWER_STDOUT_(subcomm)));
92c4762a1bSJed Brown   }
93*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Credundant));
94c4762a1bSJed Brown 
95c4762a1bSJed Brown   /* Test MatCreateRedundantMatrix() with user-provided subcomm */
96c4762a1bSJed Brown   {
97c4762a1bSJed Brown     PetscSubcomm psubcomm;
98c4762a1bSJed Brown 
99*9566063dSJacob Faibussowitsch     PetscCall(PetscSubcommCreate(PETSC_COMM_WORLD,&psubcomm));
100*9566063dSJacob Faibussowitsch     PetscCall(PetscSubcommSetNumber(psubcomm,nsubcomms));
101*9566063dSJacob Faibussowitsch     PetscCall(PetscSubcommSetType(psubcomm,PETSC_SUBCOMM_CONTIGUOUS));
102c4762a1bSJed Brown     /* enable runtime switch of psubcomm type, e.g., '-psubcomm_type interlaced */
103*9566063dSJacob Faibussowitsch     PetscCall(PetscSubcommSetFromOptions(psubcomm));
104c4762a1bSJed Brown 
105*9566063dSJacob Faibussowitsch     PetscCall(MatCreateRedundantMatrix(C,nsubcomms,PetscSubcommChild(psubcomm),MAT_INITIAL_MATRIX,&Credundant));
106*9566063dSJacob Faibussowitsch     PetscCall(MatCreateRedundantMatrix(C,nsubcomms,PetscSubcommChild(psubcomm),MAT_REUSE_MATRIX,&Credundant));
107c4762a1bSJed Brown 
108*9566063dSJacob Faibussowitsch     PetscCall(PetscSubcommDestroy(&psubcomm));
109*9566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Credundant));
110c4762a1bSJed Brown   }
111c4762a1bSJed Brown 
112*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
113*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&y));
114*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
115*9566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
116b122ec5aSJacob 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