xref: /petsc/src/sys/tests/ex58.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
158b5cd2aSSatish Balay 
258b5cd2aSSatish Balay static char help[] = "Tests PetscGlobalMinMax\n\n";
358b5cd2aSSatish Balay 
458b5cd2aSSatish Balay #include <petscsys.h>
558b5cd2aSSatish Balay 
658b5cd2aSSatish Balay int main(int argc,char **argv)
758b5cd2aSSatish Balay {
858b5cd2aSSatish Balay   PetscErrorCode ierr;
958b5cd2aSSatish Balay   PetscMPIInt    size,rank;
1058b5cd2aSSatish Balay   PetscInt       li[2],gi[2] = {-1, -1};
1158b5cd2aSSatish Balay   PetscReal      lr[2],gr[2] = {-1., -1.};
1258b5cd2aSSatish Balay 
1358b5cd2aSSatish Balay   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
14*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
15*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
1658b5cd2aSSatish Balay 
1758b5cd2aSSatish Balay   li[0] = 4 + rank;
1858b5cd2aSSatish Balay   li[1] = -3 + size - rank;
19*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscGlobalMinMaxInt(PETSC_COMM_WORLD,li,gi));
20*5f80ce2aSJacob Faibussowitsch   if (gi[0] != 4 || gi[1] != -3+size) CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"1) Error MIN/MAX %" PetscInt_FMT " %" PetscInt_FMT "\n",gi[0],gi[1]));
21*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscGlobalMinMaxInt(PETSC_COMM_WORLD,li,li));
22*5f80ce2aSJacob Faibussowitsch   if (li[0] != gi[0] || li[1] != gi[1]) CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"2) Error MIN/MAX %" PetscInt_FMT " %" PetscInt_FMT "\n",li[0],li[1]));
2358b5cd2aSSatish Balay 
2458b5cd2aSSatish Balay   if (rank == 0) {
2558b5cd2aSSatish Balay     li[0] = PETSC_MAX_INT;
2658b5cd2aSSatish Balay     li[1] = PETSC_MIN_INT;
2758b5cd2aSSatish Balay   } else if (rank == 1) {
2858b5cd2aSSatish Balay     li[0] = PETSC_MIN_INT;
2958b5cd2aSSatish Balay     li[1] = PETSC_MAX_INT;
3058b5cd2aSSatish Balay   }
3158b5cd2aSSatish Balay 
32*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscGlobalMinMaxInt(PETSC_COMM_WORLD,li,gi));
33*5f80ce2aSJacob Faibussowitsch   if (gi[0] > li[0] || gi[1] < li[1]) CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"3) Error MIN/MAX %" PetscInt_FMT " %" PetscInt_FMT "\n",gi[0],gi[1]));
3458b5cd2aSSatish Balay 
3558b5cd2aSSatish Balay   lr[0] = 4.0 + rank;
3658b5cd2aSSatish Balay   lr[1] = -3.0 + size - rank;
37*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscGlobalMinMaxReal(PETSC_COMM_WORLD,lr,gr));
38*5f80ce2aSJacob Faibussowitsch   if (gr[0] != 4.0 || gr[1] != -3.0+size) CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"4) Error MIN/MAX %g %g\n",(double)gr[0],(double)gr[1]));
39*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscGlobalMinMaxReal(PETSC_COMM_WORLD,lr,lr));
40*5f80ce2aSJacob Faibussowitsch   if (lr[0] != gr[0] || lr[1] != gr[1]) CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"5) Error MIN/MAX %g %g\n",(double)lr[0],(double)li[1]));
4158b5cd2aSSatish Balay 
4258b5cd2aSSatish Balay   ierr = PetscFinalize();
4358b5cd2aSSatish Balay   return ierr;
4458b5cd2aSSatish Balay }
4558b5cd2aSSatish Balay 
4658b5cd2aSSatish Balay /*TEST
4758b5cd2aSSatish Balay 
4858b5cd2aSSatish Balay    test:
4958b5cd2aSSatish Balay      output_file: output/ex58_1.out
5058b5cd2aSSatish Balay 
5158b5cd2aSSatish Balay    test:
5258b5cd2aSSatish Balay      suffix: 2
5358b5cd2aSSatish Balay      output_file: output/ex58_1.out
5458b5cd2aSSatish Balay      nsize: 2
5558b5cd2aSSatish Balay 
5658b5cd2aSSatish Balay TEST*/
57