xref: /petsc/src/mat/tests/ex9.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] = "Tests MPI parallel matrix creation. Test MatCreateRedundantMatrix() \n\n";
3*c4762a1bSJed Brown 
4*c4762a1bSJed Brown #include <petscmat.h>
5*c4762a1bSJed Brown 
6*c4762a1bSJed Brown int main(int argc,char **args)
7*c4762a1bSJed Brown {
8*c4762a1bSJed Brown   Mat            C,Credundant;
9*c4762a1bSJed Brown   MatInfo        info;
10*c4762a1bSJed Brown   PetscMPIInt    rank,size,subsize;
11*c4762a1bSJed Brown   PetscInt       i,j,m = 3,n = 2,low,high,iglobal;
12*c4762a1bSJed Brown   PetscInt       Ii,J,ldim,nsubcomms;
13*c4762a1bSJed Brown   PetscErrorCode ierr;
14*c4762a1bSJed Brown   PetscBool      flg_info,flg_mat;
15*c4762a1bSJed Brown   PetscScalar    v,one = 1.0;
16*c4762a1bSJed Brown   Vec            x,y;
17*c4762a1bSJed Brown   MPI_Comm       subcomm;
18*c4762a1bSJed Brown 
19*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
20*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);CHKERRQ(ierr);
21*c4762a1bSJed Brown   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
22*c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
23*c4762a1bSJed Brown   n    = 2*size;
24*c4762a1bSJed Brown 
25*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&C);CHKERRQ(ierr);
26*c4762a1bSJed Brown   ierr = MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);CHKERRQ(ierr);
27*c4762a1bSJed Brown   ierr = MatSetFromOptions(C);CHKERRQ(ierr);
28*c4762a1bSJed Brown   ierr = MatSetUp(C);CHKERRQ(ierr);
29*c4762a1bSJed Brown 
30*c4762a1bSJed Brown   /* Create the matrix for the five point stencil, YET AGAIN */
31*c4762a1bSJed Brown   for (i=0; i<m; i++) {
32*c4762a1bSJed Brown     for (j=2*rank; j<2*rank+2; j++) {
33*c4762a1bSJed Brown       v = -1.0;  Ii = j + n*i;
34*c4762a1bSJed Brown       if (i>0)   {J = Ii - n; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
35*c4762a1bSJed Brown       if (i<m-1) {J = Ii + n; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
36*c4762a1bSJed Brown       if (j>0)   {J = Ii - 1; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
37*c4762a1bSJed Brown       if (j<n-1) {J = Ii + 1; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
38*c4762a1bSJed Brown       v = 4.0; ierr = MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES);CHKERRQ(ierr);
39*c4762a1bSJed Brown     }
40*c4762a1bSJed Brown   }
41*c4762a1bSJed Brown 
42*c4762a1bSJed Brown   /* Add extra elements (to illustrate variants of MatGetInfo) */
43*c4762a1bSJed Brown   Ii   = n; J = n-2; v = 100.0;
44*c4762a1bSJed Brown   ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);
45*c4762a1bSJed Brown   Ii   = n-2; J = n; v = 100.0;
46*c4762a1bSJed Brown   ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);
47*c4762a1bSJed Brown 
48*c4762a1bSJed Brown   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
49*c4762a1bSJed Brown   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
50*c4762a1bSJed Brown 
51*c4762a1bSJed Brown   /* Form vectors */
52*c4762a1bSJed Brown   ierr = MatCreateVecs(C,&x,&y);CHKERRQ(ierr);
53*c4762a1bSJed Brown   ierr = VecGetLocalSize(x,&ldim);CHKERRQ(ierr);
54*c4762a1bSJed Brown   ierr = VecGetOwnershipRange(x,&low,&high);CHKERRQ(ierr);
55*c4762a1bSJed Brown   for (i=0; i<ldim; i++) {
56*c4762a1bSJed Brown     iglobal = i + low;
57*c4762a1bSJed Brown     v       = one*((PetscReal)i) + 100.0*rank;
58*c4762a1bSJed Brown     ierr    = VecSetValues(x,1,&iglobal,&v,INSERT_VALUES);CHKERRQ(ierr);
59*c4762a1bSJed Brown   }
60*c4762a1bSJed Brown   ierr = VecAssemblyBegin(x);CHKERRQ(ierr);
61*c4762a1bSJed Brown   ierr = VecAssemblyEnd(x);CHKERRQ(ierr);
62*c4762a1bSJed Brown 
63*c4762a1bSJed Brown   ierr = MatMult(C,x,y);CHKERRQ(ierr);
64*c4762a1bSJed Brown 
65*c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-view_info",&flg_info);CHKERRQ(ierr);
66*c4762a1bSJed Brown   if (flg_info)  {
67*c4762a1bSJed Brown     ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
68*c4762a1bSJed Brown     ierr = MatView(C,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
69*c4762a1bSJed Brown 
70*c4762a1bSJed Brown     ierr = MatGetInfo(C,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr);
71*c4762a1bSJed Brown     ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"matrix information (global sums):\nnonzeros = %D, allocated nonzeros = %D\n",(PetscInt)info.nz_used,(PetscInt)info.nz_allocated);CHKERRQ(ierr);
72*c4762a1bSJed Brown     ierr = MatGetInfo (C,MAT_GLOBAL_MAX,&info);CHKERRQ(ierr);
73*c4762a1bSJed Brown     ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"matrix information (global max):\nnonzeros = %D, allocated nonzeros = %D\n",(PetscInt)info.nz_used,(PetscInt)info.nz_allocated);CHKERRQ(ierr);
74*c4762a1bSJed Brown   }
75*c4762a1bSJed Brown 
76*c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-view_mat",&flg_mat);CHKERRQ(ierr);
77*c4762a1bSJed Brown   if (flg_mat) {
78*c4762a1bSJed Brown     ierr = MatView(C,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
79*c4762a1bSJed Brown   }
80*c4762a1bSJed Brown 
81*c4762a1bSJed Brown   /* Test MatCreateRedundantMatrix() */
82*c4762a1bSJed Brown   nsubcomms = size;
83*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-nsubcomms",&nsubcomms,NULL);CHKERRQ(ierr);
84*c4762a1bSJed Brown   ierr = MatCreateRedundantMatrix(C,nsubcomms,MPI_COMM_NULL,MAT_INITIAL_MATRIX,&Credundant);CHKERRQ(ierr);
85*c4762a1bSJed Brown   ierr = MatCreateRedundantMatrix(C,nsubcomms,MPI_COMM_NULL,MAT_REUSE_MATRIX,&Credundant);CHKERRQ(ierr);
86*c4762a1bSJed Brown 
87*c4762a1bSJed Brown   ierr = PetscObjectGetComm((PetscObject)Credundant,&subcomm);CHKERRQ(ierr);
88*c4762a1bSJed Brown   ierr = MPI_Comm_size(subcomm,&subsize);CHKERRQ(ierr);
89*c4762a1bSJed Brown 
90*c4762a1bSJed Brown   if (subsize==2 && flg_mat) {
91*c4762a1bSJed Brown     ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_(subcomm),"\n[%d] Credundant:\n",rank);CHKERRQ(ierr);
92*c4762a1bSJed Brown     ierr = MatView(Credundant,PETSC_VIEWER_STDOUT_(subcomm));CHKERRQ(ierr);
93*c4762a1bSJed Brown   }
94*c4762a1bSJed Brown   ierr = MatDestroy(&Credundant);CHKERRQ(ierr);
95*c4762a1bSJed Brown 
96*c4762a1bSJed Brown   /* Test MatCreateRedundantMatrix() with user-provided subcomm */
97*c4762a1bSJed Brown   {
98*c4762a1bSJed Brown     PetscSubcomm psubcomm;
99*c4762a1bSJed Brown 
100*c4762a1bSJed Brown     ierr = PetscSubcommCreate(PETSC_COMM_WORLD,&psubcomm);CHKERRQ(ierr);
101*c4762a1bSJed Brown     ierr = PetscSubcommSetNumber(psubcomm,nsubcomms);CHKERRQ(ierr);
102*c4762a1bSJed Brown     ierr = PetscSubcommSetType(psubcomm,PETSC_SUBCOMM_CONTIGUOUS);CHKERRQ(ierr);
103*c4762a1bSJed Brown     /* enable runtime switch of psubcomm type, e.g., '-psubcomm_type interlaced */
104*c4762a1bSJed Brown     ierr = PetscSubcommSetFromOptions(psubcomm);CHKERRQ(ierr);
105*c4762a1bSJed Brown 
106*c4762a1bSJed Brown     ierr = MatCreateRedundantMatrix(C,nsubcomms,PetscSubcommChild(psubcomm),MAT_INITIAL_MATRIX,&Credundant);CHKERRQ(ierr);
107*c4762a1bSJed Brown     ierr = MatCreateRedundantMatrix(C,nsubcomms,PetscSubcommChild(psubcomm),MAT_REUSE_MATRIX,&Credundant);CHKERRQ(ierr);
108*c4762a1bSJed Brown 
109*c4762a1bSJed Brown     ierr = PetscSubcommDestroy(&psubcomm);CHKERRQ(ierr);
110*c4762a1bSJed Brown     ierr = MatDestroy(&Credundant);CHKERRQ(ierr);
111*c4762a1bSJed Brown   }
112*c4762a1bSJed Brown 
113*c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);
114*c4762a1bSJed Brown   ierr = VecDestroy(&y);CHKERRQ(ierr);
115*c4762a1bSJed Brown   ierr = MatDestroy(&C);CHKERRQ(ierr);
116*c4762a1bSJed Brown   ierr = PetscFinalize();
117*c4762a1bSJed Brown   return ierr;
118*c4762a1bSJed Brown }
119*c4762a1bSJed Brown 
120*c4762a1bSJed Brown 
121*c4762a1bSJed Brown /*TEST
122*c4762a1bSJed Brown 
123*c4762a1bSJed Brown    test:
124*c4762a1bSJed Brown       nsize: 3
125*c4762a1bSJed Brown       args: -view_info
126*c4762a1bSJed Brown 
127*c4762a1bSJed Brown    test:
128*c4762a1bSJed Brown       suffix: 2
129*c4762a1bSJed Brown       nsize: 3
130*c4762a1bSJed Brown       args: -nsubcomms 2 -view_mat -psubcomm_type interlaced
131*c4762a1bSJed Brown 
132*c4762a1bSJed Brown    test:
133*c4762a1bSJed Brown       suffix: 3
134*c4762a1bSJed Brown       nsize: 3
135*c4762a1bSJed Brown       args: -nsubcomms 2 -view_mat -psubcomm_type contiguous
136*c4762a1bSJed Brown 
137*c4762a1bSJed Brown    test:
138*c4762a1bSJed Brown       suffix: 3_baij
139*c4762a1bSJed Brown       nsize: 3
140*c4762a1bSJed Brown       args: -mat_type baij -nsubcomms 2 -view_mat
141*c4762a1bSJed Brown 
142*c4762a1bSJed Brown    test:
143*c4762a1bSJed Brown       suffix: 3_sbaij
144*c4762a1bSJed Brown       nsize: 3
145*c4762a1bSJed Brown       args: -mat_type sbaij -nsubcomms 2 -view_mat
146*c4762a1bSJed Brown 
147*c4762a1bSJed Brown    test:
148*c4762a1bSJed Brown       suffix: 3_dense
149*c4762a1bSJed Brown       nsize: 3
150*c4762a1bSJed Brown       args: -mat_type dense -nsubcomms 2 -view_mat
151*c4762a1bSJed Brown 
152*c4762a1bSJed Brown    test:
153*c4762a1bSJed Brown       suffix: 4_baij
154*c4762a1bSJed Brown       nsize: 3
155*c4762a1bSJed Brown       args: -mat_type baij -nsubcomms 2 -view_mat -psubcomm_type interlaced
156*c4762a1bSJed Brown 
157*c4762a1bSJed Brown    test:
158*c4762a1bSJed Brown       suffix: 4_sbaij
159*c4762a1bSJed Brown       nsize: 3
160*c4762a1bSJed Brown       args: -mat_type sbaij -nsubcomms 2 -view_mat -psubcomm_type interlaced
161*c4762a1bSJed Brown 
162*c4762a1bSJed Brown    test:
163*c4762a1bSJed Brown       suffix: 4_dense
164*c4762a1bSJed Brown       nsize: 3
165*c4762a1bSJed Brown       args: -mat_type dense -nsubcomms 2 -view_mat -psubcomm_type interlaced
166*c4762a1bSJed Brown 
167*c4762a1bSJed Brown TEST*/
168