158b5cd2aSSatish Balay static char help[] = "Tests PetscGlobalMinMax\n\n";
258b5cd2aSSatish Balay
358b5cd2aSSatish Balay #include <petscsys.h>
458b5cd2aSSatish Balay
main(int argc,char ** argv)5d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
6d71ae5a4SJacob Faibussowitsch {
758b5cd2aSSatish Balay PetscMPIInt size, rank;
858b5cd2aSSatish Balay PetscInt li[2], gi[2] = {-1, -1};
958b5cd2aSSatish Balay PetscReal lr[2], gr[2] = {-1., -1.};
1058b5cd2aSSatish Balay
11327415f7SBarry Smith PetscFunctionBeginUser;
12c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
139566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
149566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
1558b5cd2aSSatish Balay
1658b5cd2aSSatish Balay li[0] = 4 + rank;
1758b5cd2aSSatish Balay li[1] = -3 + size - rank;
189566063dSJacob Faibussowitsch PetscCall(PetscGlobalMinMaxInt(PETSC_COMM_WORLD, li, gi));
199566063dSJacob Faibussowitsch if (gi[0] != 4 || gi[1] != -3 + size) PetscCall(PetscPrintf(PETSC_COMM_SELF, "1) Error MIN/MAX %" PetscInt_FMT " %" PetscInt_FMT "\n", gi[0], gi[1]));
209566063dSJacob Faibussowitsch PetscCall(PetscGlobalMinMaxInt(PETSC_COMM_WORLD, li, li));
219566063dSJacob Faibussowitsch if (li[0] != gi[0] || li[1] != gi[1]) PetscCall(PetscPrintf(PETSC_COMM_SELF, "2) Error MIN/MAX %" PetscInt_FMT " %" PetscInt_FMT "\n", li[0], li[1]));
2258b5cd2aSSatish Balay
2358b5cd2aSSatish Balay if (rank == 0) {
241690c2aeSBarry Smith li[0] = PETSC_INT_MAX;
251690c2aeSBarry Smith li[1] = PETSC_INT_MIN;
2658b5cd2aSSatish Balay } else if (rank == 1) {
271690c2aeSBarry Smith li[0] = PETSC_INT_MIN;
281690c2aeSBarry Smith li[1] = PETSC_INT_MAX;
2958b5cd2aSSatish Balay }
3058b5cd2aSSatish Balay
319566063dSJacob Faibussowitsch PetscCall(PetscGlobalMinMaxInt(PETSC_COMM_WORLD, li, gi));
329566063dSJacob Faibussowitsch if (gi[0] > li[0] || gi[1] < li[1]) PetscCall(PetscPrintf(PETSC_COMM_SELF, "3) Error MIN/MAX %" PetscInt_FMT " %" PetscInt_FMT "\n", gi[0], gi[1]));
3358b5cd2aSSatish Balay
3458b5cd2aSSatish Balay lr[0] = 4.0 + rank;
3558b5cd2aSSatish Balay lr[1] = -3.0 + size - rank;
369566063dSJacob Faibussowitsch PetscCall(PetscGlobalMinMaxReal(PETSC_COMM_WORLD, lr, gr));
379566063dSJacob Faibussowitsch if (gr[0] != 4.0 || gr[1] != -3.0 + size) PetscCall(PetscPrintf(PETSC_COMM_SELF, "4) Error MIN/MAX %g %g\n", (double)gr[0], (double)gr[1]));
389566063dSJacob Faibussowitsch PetscCall(PetscGlobalMinMaxReal(PETSC_COMM_WORLD, lr, lr));
399566063dSJacob Faibussowitsch if (lr[0] != gr[0] || lr[1] != gr[1]) PetscCall(PetscPrintf(PETSC_COMM_SELF, "5) Error MIN/MAX %g %g\n", (double)lr[0], (double)li[1]));
4058b5cd2aSSatish Balay
419566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
42b122ec5aSJacob Faibussowitsch return 0;
4358b5cd2aSSatish Balay }
4458b5cd2aSSatish Balay
4558b5cd2aSSatish Balay /*TEST
4658b5cd2aSSatish Balay
4758b5cd2aSSatish Balay test:
48*3886731fSPierre Jolivet output_file: output/empty.out
4958b5cd2aSSatish Balay
5058b5cd2aSSatish Balay test:
5158b5cd2aSSatish Balay suffix: 2
52*3886731fSPierre Jolivet output_file: output/empty.out
5358b5cd2aSSatish Balay nsize: 2
5458b5cd2aSSatish Balay
5558b5cd2aSSatish Balay TEST*/
56