xref: /petsc/src/mat/tests/ex114.c (revision 4e879ede91dffd7f752e2ff7506bc58e6af186e3)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests MatGetRowMax(), MatGetRowMin(), MatGetRowMaxAbs()\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown #define M 5
7c4762a1bSJed Brown #define N 6
8c4762a1bSJed Brown 
9c4762a1bSJed Brown int main(int argc,char **args)
10c4762a1bSJed Brown {
11c4762a1bSJed Brown   Mat            A;
12fa213d2fSHong Zhang   Vec            min,max,maxabs,e;
13*4e879edeSHong Zhang   PetscInt       m,n,j,imin[M],imax[M],imaxabs[M],indices[N],row,testcase=0;
14c4762a1bSJed Brown   PetscScalar    values[N];
15c4762a1bSJed Brown   PetscErrorCode ierr;
161a254869SHong Zhang   PetscMPIInt    size,rank;
17fa213d2fSHong Zhang   PetscReal      enorm;
18c4762a1bSJed Brown 
19c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
20c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
211a254869SHong Zhang   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
221a254869SHong Zhang   ierr = PetscOptionsGetInt(NULL,NULL,"-testcase",&testcase,NULL);CHKERRQ(ierr);
23c4762a1bSJed Brown 
24c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
251a254869SHong Zhang   if (testcase == 1) { /* proc[0] holds entire A and other processes have no entry */
261a254869SHong Zhang     if (!rank) {
271a254869SHong Zhang       ierr = MatSetSizes(A,M,N,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
281a254869SHong Zhang     } else {
291a254869SHong Zhang       ierr = MatSetSizes(A,0,0,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
301a254869SHong Zhang     }
311a254869SHong Zhang     testcase = 0;
321a254869SHong Zhang   } else {
33c4762a1bSJed Brown     ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr);
341a254869SHong Zhang   }
35c4762a1bSJed Brown   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
36c4762a1bSJed Brown   ierr = MatSetUp(A);CHKERRQ(ierr);
371a254869SHong Zhang 
381a254869SHong Zhang   if (!rank) { /* proc[0] sets matrix A */
391a254869SHong Zhang     for (j=0; j<N; j++) indices[j] = j;
401a254869SHong Zhang     switch (testcase) {
411a254869SHong Zhang     case 1: /* see testcast 0 */
421a254869SHong Zhang       break;
431a254869SHong Zhang     case 2:
44c4762a1bSJed Brown       row = 0;
45f07e67edSHong Zhang       values[0]  = -2.0; values[1] = -2.0; values[2] = -2.0; values[3] = -4.0; values[4] = 1.0; values[5] = 1.0;
461a254869SHong Zhang       ierr = MatSetValues(A,1,&row,N,indices,values,INSERT_VALUES);CHKERRQ(ierr);
471a254869SHong Zhang       row = 2;
481a254869SHong Zhang       indices[0] = 0;    indices[1] = 3;    indices[2] = 5;
491a254869SHong Zhang       values[0]  = -2.0; values[1]  = -2.0; values[2]  = -2.0;
501a254869SHong Zhang       ierr = MatSetValues(A,1,&row,3,indices,values,INSERT_VALUES);CHKERRQ(ierr);
511a254869SHong Zhang       row = 3;
521a254869SHong Zhang       indices[0] = 0;    indices[1] = 1;    indices[2] = 4;
531a254869SHong Zhang       values[0]  = -2.0; values[1]  = -2.0; values[2]  = -2.0;
54c4762a1bSJed Brown       ierr = MatSetValues(A,1,&row,3,indices,values,INSERT_VALUES);CHKERRQ(ierr);
55c4762a1bSJed Brown       row = 4;
561a254869SHong Zhang       indices[0] = 0;    indices[1] = 1;    indices[2] = 2;
571a254869SHong Zhang       values[0]  = -2.0; values[1]  = -2.0; values[2]  = -2.0;
581a254869SHong Zhang       ierr = MatSetValues(A,1,&row,3,indices,values,INSERT_VALUES);CHKERRQ(ierr);
591a254869SHong Zhang       break;
601a254869SHong Zhang     case 3:
611a254869SHong Zhang       row = 0;
621a254869SHong Zhang       values[0]  = -2.0; values[1] = -2.0; values[2] = -2.0;
631a254869SHong Zhang       ierr = MatSetValues(A,1,&row,3,indices+1,values,INSERT_VALUES);CHKERRQ(ierr);
641a254869SHong Zhang       row = 1;
651a254869SHong Zhang       values[0]  = -2.0; values[1] = -2.0; values[2] = -2.0;
661a254869SHong Zhang       ierr = MatSetValues(A,1,&row,3,indices,values,INSERT_VALUES);CHKERRQ(ierr);
671a254869SHong Zhang       row = 2;
681a254869SHong Zhang       values[0]  = -2.0; values[1] = -2.0; values[2]  = -2.0;
691a254869SHong Zhang       ierr = MatSetValues(A,1,&row,3,indices,values,INSERT_VALUES);CHKERRQ(ierr);
701a254869SHong Zhang       row = 3;
711a254869SHong Zhang       values[0]  = -2.0; values[1] = -2.0; values[2]  = -2.0; values[3] = -1.0;
721a254869SHong Zhang       ierr = MatSetValues(A,1,&row,4,indices,values,INSERT_VALUES);CHKERRQ(ierr);
731a254869SHong Zhang       row = 4;
741a254869SHong Zhang       values[0]  = -2.0; values[1] = -2.0; values[2]  = -2.0; values[3] = -1.0;
751a254869SHong Zhang       ierr = MatSetValues(A,1,&row,4,indices,values,INSERT_VALUES);CHKERRQ(ierr);
761a254869SHong Zhang       break;
771a254869SHong Zhang 
781a254869SHong Zhang     default:
791a254869SHong Zhang       row  = 0;
801a254869SHong Zhang       values[0]  = -1.0; values[1] = 0.0; values[2] = 1.0; values[3] = 3.0; values[4] = 4.0; values[5] = -5.0;
811a254869SHong Zhang       ierr = MatSetValues(A,1,&row,N,indices,values,INSERT_VALUES);CHKERRQ(ierr);
821a254869SHong Zhang       row  = 1;
831a254869SHong Zhang       ierr = MatSetValues(A,1,&row,3,indices,values,INSERT_VALUES);CHKERRQ(ierr);
841a254869SHong Zhang       row  = 3;
85c4762a1bSJed Brown       ierr = MatSetValues(A,1,&row,1,indices+4,values+4,INSERT_VALUES);CHKERRQ(ierr);
86c4762a1bSJed Brown       row  = 4;
87c4762a1bSJed Brown       ierr = MatSetValues(A,1,&row,2,indices+4,values+4,INSERT_VALUES);CHKERRQ(ierr);
881a254869SHong Zhang     }
891a254869SHong Zhang   }
90c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
91c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
92c4762a1bSJed Brown   ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
93c4762a1bSJed Brown 
94c4762a1bSJed Brown   ierr = MatGetLocalSize(A, &m,&n);CHKERRQ(ierr);
95c4762a1bSJed Brown   ierr = VecCreate(PETSC_COMM_WORLD,&min);CHKERRQ(ierr);
96c4762a1bSJed Brown   ierr = VecSetSizes(min,m,PETSC_DECIDE);CHKERRQ(ierr);
97c4762a1bSJed Brown   ierr = VecSetFromOptions(min);CHKERRQ(ierr);
98c4762a1bSJed Brown   ierr = VecDuplicate(min,&max);CHKERRQ(ierr);
99c4762a1bSJed Brown   ierr = VecDuplicate(min,&maxabs);CHKERRQ(ierr);
100fa213d2fSHong Zhang   ierr = VecDuplicate(min,&e);CHKERRQ(ierr);
101c4762a1bSJed Brown 
102f07e67edSHong Zhang   /* MatGetRowMax() */
103f07e67edSHong Zhang   ierr = PetscPrintf(PETSC_COMM_WORLD,"\n MatGetRowMax\n");CHKERRQ(ierr);
1041a254869SHong Zhang   ierr = MatGetRowMax(A,max,NULL);CHKERRQ(ierr);
1051a254869SHong Zhang   ierr = MatGetRowMax(A,max,imax);CHKERRQ(ierr);
1061a254869SHong Zhang   ierr = VecView(max,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1071a254869SHong Zhang   ierr = VecGetLocalSize(max,&n);CHKERRQ(ierr);
1081a254869SHong Zhang   ierr = PetscIntView(n,imax,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
109fa213d2fSHong Zhang 
110f07e67edSHong Zhang   /* MatGetRowMin() */
111fa213d2fSHong Zhang   ierr = MatScale(A,-1.0);CHKERRQ(ierr);
112f07e67edSHong Zhang   ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
113f07e67edSHong Zhang   ierr = PetscPrintf(PETSC_COMM_WORLD,"\n MatGetRowMin\n");CHKERRQ(ierr);
114fa213d2fSHong Zhang   ierr = MatGetRowMin(A,min,NULL);CHKERRQ(ierr);
115fa213d2fSHong Zhang   ierr = MatGetRowMin(A,min,imin);CHKERRQ(ierr);
116fa213d2fSHong Zhang 
117fa213d2fSHong Zhang   ierr = VecWAXPY(e,1.0,max,min);CHKERRQ(ierr); /* e = max + min */
118fa213d2fSHong Zhang   ierr = VecNorm(e,NORM_INFINITY,&enorm);CHKERRQ(ierr);
119fa213d2fSHong Zhang   if (enorm > PETSC_MACHINE_EPSILON) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"max+min > PETSC_MACHINE_EPSILON ");
120f07e67edSHong Zhang   for (j = 0; j < n; j++) {
121*4e879edeSHong Zhang     if (imin[j] != imax[j]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"imin[%D] %D != imax %D",j,imin[j],imax[j]);
122f07e67edSHong Zhang   }
123fa213d2fSHong Zhang 
124f07e67edSHong Zhang   /* MatGetRowMaxAbs() */
125f07e67edSHong Zhang   ierr = PetscPrintf(PETSC_COMM_WORLD,"\n MatGetRowMaxAbs\n");CHKERRQ(ierr);
126475b8b61SHong Zhang   ierr = MatGetRowMaxAbs(A,maxabs,NULL);CHKERRQ(ierr);
127475b8b61SHong Zhang   ierr = MatGetRowMaxAbs(A,maxabs,imaxabs);CHKERRQ(ierr);
128475b8b61SHong Zhang   ierr = VecView(maxabs,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
129475b8b61SHong Zhang   ierr = PetscIntView(n,imaxabs,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1301a254869SHong Zhang 
131f07e67edSHong Zhang   /* MatGetRowMinAbs() */
132f07e67edSHong Zhang   ierr = PetscPrintf(PETSC_COMM_WORLD,"\n MatGetRowMinAbs\n");CHKERRQ(ierr);
133f07e67edSHong Zhang   ierr = MatGetRowMinAbs(A,min,NULL);CHKERRQ(ierr);
134f07e67edSHong Zhang   ierr = MatGetRowMinAbs(A,min,imin);CHKERRQ(ierr);
135c4762a1bSJed Brown   ierr = VecView(min,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
136f07e67edSHong Zhang   ierr = PetscIntView(n,imin,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
137f07e67edSHong Zhang 
138f07e67edSHong Zhang   if (size == 1) {
139f07e67edSHong Zhang     /* Test MatGetRowMax, MatGetRowMin and MatGetRowMaxAbs for SeqDense and MPIBAIJ matrix */
140*4e879edeSHong Zhang     Mat Adense;
141*4e879edeSHong Zhang     Vec max_d,maxabs_d;
142*4e879edeSHong Zhang     ierr = VecDuplicate(min,&max_d);CHKERRQ(ierr);
143*4e879edeSHong Zhang     ierr = VecDuplicate(min,&maxabs_d);CHKERRQ(ierr);
144f07e67edSHong Zhang 
145f07e67edSHong Zhang     ierr = MatScale(A,-1.0);CHKERRQ(ierr);
146*4e879edeSHong Zhang     ierr = MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&Adense);CHKERRQ(ierr);
147*4e879edeSHong Zhang     ierr = PetscPrintf(PETSC_COMM_WORLD,"MatGetRowMax for seqdense matrix\n");CHKERRQ(ierr);
148*4e879edeSHong Zhang     ierr = MatGetRowMax(Adense,max_d,imax);CHKERRQ(ierr);
149*4e879edeSHong Zhang 
150*4e879edeSHong Zhang     ierr = VecWAXPY(e,-1.0,max,max_d);CHKERRQ(ierr); /* e = -max + max_d */
151*4e879edeSHong Zhang     ierr = VecNorm(e,NORM_INFINITY,&enorm);CHKERRQ(ierr);
152*4e879edeSHong Zhang     if (enorm > PETSC_MACHINE_EPSILON) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"norm(-max + max_d) %g > PETSC_MACHINE_EPSILON",enorm);
153*4e879edeSHong Zhang 
154*4e879edeSHong Zhang     ierr = MatScale(Adense,-1.0);CHKERRQ(ierr);
155*4e879edeSHong Zhang     ierr = PetscPrintf(PETSC_COMM_WORLD,"MatGetRowMin for seqdense matrix\n");CHKERRQ(ierr);
156*4e879edeSHong Zhang     ierr = MatGetRowMin(Adense,min,imin);CHKERRQ(ierr);
157f07e67edSHong Zhang 
158f07e67edSHong Zhang     ierr = VecWAXPY(e,1.0,max,min);CHKERRQ(ierr); /* e = max + min */
159f07e67edSHong Zhang     ierr = VecNorm(e,NORM_INFINITY,&enorm);CHKERRQ(ierr);
160f07e67edSHong Zhang     if (enorm > PETSC_MACHINE_EPSILON) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"max+min > PETSC_MACHINE_EPSILON ");
161f07e67edSHong Zhang     for (j = 0; j < n; j++) {
162f07e67edSHong Zhang       if (imin[j] != imax[j]) {
163*4e879edeSHong Zhang         SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"imin[%D] %D != imax %D for seqdense matrix",j,imin[j],imax[j]);
164f07e67edSHong Zhang       }
165f07e67edSHong Zhang     }
166f07e67edSHong Zhang 
167*4e879edeSHong Zhang     ierr = PetscPrintf(PETSC_COMM_WORLD,"MatGetRowMaxAbs for seqdense matrix\n");CHKERRQ(ierr);
168*4e879edeSHong Zhang     ierr = MatGetRowMaxAbs(Adense,maxabs_d,imaxabs);CHKERRQ(ierr);
169*4e879edeSHong Zhang     ierr = VecWAXPY(e,-1.0,maxabs,maxabs_d);CHKERRQ(ierr); /* e = -maxabs + maxabs_d */
170*4e879edeSHong Zhang     ierr = VecNorm(e,NORM_INFINITY,&enorm);CHKERRQ(ierr);
171*4e879edeSHong Zhang     if (enorm > PETSC_MACHINE_EPSILON) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"norm(-maxabs + maxabs_d) %g > PETSC_MACHINE_EPSILON",enorm);
172c4762a1bSJed Brown 
173*4e879edeSHong Zhang     ierr = MatDestroy(&Adense);CHKERRQ(ierr);
174*4e879edeSHong Zhang     ierr = VecDestroy(&max_d);CHKERRQ(ierr);
175*4e879edeSHong Zhang     ierr = VecDestroy(&maxabs_d);CHKERRQ(ierr);
176c4762a1bSJed Brown   } else {
177c4762a1bSJed Brown     /* BAIJ */
178*4e879edeSHong Zhang     Mat      B;
179*4e879edeSHong Zhang     Vec      maxabsB;
180*4e879edeSHong Zhang     PetscInt imaxabsB[M];
181*4e879edeSHong Zhang     ierr = MatConvert(A,MATMPIBAIJ,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr);
182*4e879edeSHong Zhang     ierr = VecDuplicate(min,&maxabsB);CHKERRQ(ierr);
183*4e879edeSHong Zhang     ierr = MatGetRowMaxAbs(A,maxabsB,NULL);CHKERRQ(ierr);
184*4e879edeSHong Zhang     ierr = MatGetRowMaxAbs(A,maxabsB,imaxabsB);CHKERRQ(ierr);
185*4e879edeSHong Zhang     ierr = PetscPrintf(PETSC_COMM_WORLD,"\n MatGetRowMaxAbs for MPIBAIJ matrix\n");CHKERRQ(ierr);
186*4e879edeSHong Zhang     ierr = VecWAXPY(e,-1.0,maxabs,maxabsB);CHKERRQ(ierr); /* e = -maxabs + maxabsB */
187*4e879edeSHong Zhang     ierr = VecNorm(e,NORM_INFINITY,&enorm);CHKERRQ(ierr);
188*4e879edeSHong Zhang     if (enorm > PETSC_MACHINE_EPSILON) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"norm(-maxabs + maxabs_d) %g > PETSC_MACHINE_EPSILON",enorm);
189*4e879edeSHong Zhang 
190*4e879edeSHong Zhang     for (j = 0; j < n; j++) {
191*4e879edeSHong Zhang       if (imaxabs[j] != imaxabsB[j]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"imaxabs[%D] %D != imaxabsB %D",j,imin[j],imax[j]);
192c4762a1bSJed Brown     }
193*4e879edeSHong Zhang     ierr = MatDestroy(&B);CHKERRQ(ierr);
194*4e879edeSHong Zhang     ierr = VecDestroy(&maxabsB);CHKERRQ(ierr);
195c4762a1bSJed Brown   }
196c4762a1bSJed Brown 
197c4762a1bSJed Brown   ierr = VecDestroy(&min);CHKERRQ(ierr);
198c4762a1bSJed Brown   ierr = VecDestroy(&max);CHKERRQ(ierr);
199c4762a1bSJed Brown   ierr = VecDestroy(&maxabs);CHKERRQ(ierr);
200fa213d2fSHong Zhang   ierr = VecDestroy(&e);CHKERRQ(ierr);
201c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
202c4762a1bSJed Brown   ierr = PetscFinalize();
203c4762a1bSJed Brown   return ierr;
204c4762a1bSJed Brown }
205c4762a1bSJed Brown 
206c4762a1bSJed Brown /*TEST
207c4762a1bSJed Brown 
208c4762a1bSJed Brown    test:
209c4762a1bSJed Brown       output_file: output/ex114.out
210c4762a1bSJed Brown 
211c4762a1bSJed Brown    test:
212c4762a1bSJed Brown       suffix: 2
213*4e879edeSHong Zhang       args: -testcase 1
214*4e879edeSHong Zhang       output_file: output/ex114.out
215c4762a1bSJed Brown 
216c4762a1bSJed Brown    test:
217c4762a1bSJed Brown       suffix: 3
218*4e879edeSHong Zhang       args: -testcase 2
219*4e879edeSHong Zhang       output_file: output/ex114_3.out
220*4e879edeSHong Zhang 
221*4e879edeSHong Zhang    test:
222*4e879edeSHong Zhang       suffix: 4
223*4e879edeSHong Zhang       args: -testcase 3
224*4e879edeSHong Zhang       output_file: output/ex114_4.out
225*4e879edeSHong Zhang 
226*4e879edeSHong Zhang    test:
227*4e879edeSHong Zhang       suffix: 5
228*4e879edeSHong Zhang       nsize: 3
229*4e879edeSHong Zhang       args: -testcase 0
230*4e879edeSHong Zhang       output_file: output/ex114_5.out
231*4e879edeSHong Zhang 
232*4e879edeSHong Zhang    test:
233*4e879edeSHong Zhang       suffix: 6
234*4e879edeSHong Zhang       nsize: 3
235*4e879edeSHong Zhang       args: -testcase 1
236*4e879edeSHong Zhang       output_file: output/ex114_6.out
237*4e879edeSHong Zhang 
238*4e879edeSHong Zhang    test:
239*4e879edeSHong Zhang       suffix: 7
240*4e879edeSHong Zhang       nsize: 3
241*4e879edeSHong Zhang       args: -testcase 2
242*4e879edeSHong Zhang       output_file: output/ex114_7.out
243*4e879edeSHong Zhang 
244*4e879edeSHong Zhang    test:
245*4e879edeSHong Zhang       suffix: 8
246*4e879edeSHong Zhang       nsize: 3
247*4e879edeSHong Zhang       args: -testcase 3
248*4e879edeSHong Zhang       output_file: output/ex114_8.out
249c4762a1bSJed Brown 
250c4762a1bSJed Brown TEST*/
251