1c4762a1bSJed Brown static char help[]= "Test PetscSFFCompose when the ilocal arrays are not identity nor dense\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petsc.h> 4c4762a1bSJed Brown #include <petscsf.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown int main(int argc, char **argv) 7c4762a1bSJed Brown { 8c4762a1bSJed Brown PetscSF sfA, sfB, sfBA, sfAAm, sfBBm, sfAm, sfBm; 9c4762a1bSJed Brown PetscInt nrootsA, nleavesA, nrootsB, nleavesB; 10c4762a1bSJed Brown PetscInt *ilocalA, *ilocalB; 11c4762a1bSJed Brown PetscSFNode *iremoteA, *iremoteB; 12c4762a1bSJed Brown PetscMPIInt rank,size; 13c4762a1bSJed Brown PetscInt i,m,n,k,nl = 2,mA,mB,nldataA,nldataB; 14c4762a1bSJed Brown PetscInt *rdA,*rdB,*ldA,*ldB; 15c4762a1bSJed Brown PetscBool inverse = PETSC_FALSE; 16c4762a1bSJed Brown 17*327415f7SBarry Smith PetscFunctionBeginUser; 189566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,NULL,help)); 199566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-nl",&nl,NULL)); 209566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-explicit_inverse",&inverse,NULL)); 219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 23c4762a1bSJed Brown 249566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(PETSC_COMM_WORLD, &sfA)); 259566063dSJacob Faibussowitsch PetscCall(PetscSFCreate(PETSC_COMM_WORLD, &sfB)); 269566063dSJacob Faibussowitsch PetscCall(PetscSFSetFromOptions(sfA)); 279566063dSJacob Faibussowitsch PetscCall(PetscSFSetFromOptions(sfB)); 28c4762a1bSJed Brown 29c4762a1bSJed Brown n = 4*nl*size; 30c4762a1bSJed Brown m = 2*nl; 31c4762a1bSJed Brown k = nl; 32c4762a1bSJed Brown 33dd400576SPatrick Sanan nldataA = rank == 0 ? n : 0; 34c4762a1bSJed Brown nldataB = 3*nl; 35c4762a1bSJed Brown 36c4762a1bSJed Brown nrootsA = m; 37dd400576SPatrick Sanan nleavesA = rank == 0 ? size*m : 0; 38dd400576SPatrick Sanan nrootsB = rank == 0 ? n : 0; 39c4762a1bSJed Brown nleavesB = k; 40c4762a1bSJed Brown 419566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleavesA, &ilocalA)); 429566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleavesA, &iremoteA)); 439566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleavesB, &ilocalB)); 449566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nleavesB, &iremoteB)); 45c4762a1bSJed Brown 46c4762a1bSJed Brown /* sf A bcast is equivalent to a sparse gather on process 0 47c4762a1bSJed Brown process 0 receives data in the middle [nl,3*nl] of the leaf data array for A */ 48c4762a1bSJed Brown for (i = 0; i < nleavesA; i++) { 49c4762a1bSJed Brown iremoteA[i].rank = i/m; 50c4762a1bSJed Brown iremoteA[i].index = i%m; 51c4762a1bSJed Brown ilocalA[i] = nl + i/m * 4*nl + i%m; 52c4762a1bSJed Brown } 53c4762a1bSJed Brown 54c4762a1bSJed Brown /* sf B bcast is equivalent to a sparse scatter from process 0 55c4762a1bSJed Brown process 0 sends data from [nl,2*nl] of the leaf data array for A 56c4762a1bSJed Brown each process receives, in reverse order, in the middle [nl,2*nl] of the leaf data array for B */ 57c4762a1bSJed Brown for (i = 0; i < nleavesB; i++) { 58c4762a1bSJed Brown iremoteB[i].rank = 0; 59c4762a1bSJed Brown iremoteB[i].index = rank * 4*nl + nl + i%m; 60c4762a1bSJed Brown ilocalB[i] = 2*nl - i - 1; 61c4762a1bSJed Brown } 629566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(sfA, nrootsA, nleavesA, ilocalA, PETSC_OWN_POINTER, iremoteA, PETSC_OWN_POINTER)); 639566063dSJacob Faibussowitsch PetscCall(PetscSFSetGraph(sfB, nrootsB, nleavesB, ilocalB, PETSC_OWN_POINTER, iremoteB, PETSC_OWN_POINTER)); 649566063dSJacob Faibussowitsch PetscCall(PetscSFSetUp(sfA)); 659566063dSJacob Faibussowitsch PetscCall(PetscSFSetUp(sfB)); 669566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)sfA, "sfA")); 679566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)sfB, "sfB")); 689566063dSJacob Faibussowitsch PetscCall(PetscSFViewFromOptions(sfA, NULL, "-view")); 699566063dSJacob Faibussowitsch PetscCall(PetscSFViewFromOptions(sfB, NULL, "-view")); 70c4762a1bSJed Brown 719566063dSJacob Faibussowitsch PetscCall(PetscSFGetLeafRange(sfA, NULL, &mA)); 729566063dSJacob Faibussowitsch PetscCall(PetscSFGetLeafRange(sfB, NULL, &mB)); 739566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nrootsA, &rdA, nldataA, &ldA)); 749566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nrootsB, &rdB, nldataB, &ldB)); 75c4762a1bSJed Brown for (i = 0; i < nrootsA; i++) rdA[i] = m*rank + i; 76c4762a1bSJed Brown for (i = 0; i < nldataA; i++) ldA[i] = -1; 77c4762a1bSJed Brown for (i = 0; i < nldataB; i++) ldB[i] = -1; 78c4762a1bSJed Brown 799566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "BcastB(BcastA)\n")); 809566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "A: root data\n")); 819566063dSJacob Faibussowitsch PetscCall(PetscIntView(nrootsA, rdA, PETSC_VIEWER_STDOUT_WORLD)); 829566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sfA, MPIU_INT, rdA, ldA,MPI_REPLACE)); 839566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sfA, MPIU_INT, rdA, ldA,MPI_REPLACE)); 849566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "A: leaf data (all)\n")); 859566063dSJacob Faibussowitsch PetscCall(PetscIntView(nldataA, ldA, PETSC_VIEWER_STDOUT_WORLD)); 869566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sfB, MPIU_INT, ldA, ldB,MPI_REPLACE)); 879566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sfB, MPIU_INT, ldA, ldB,MPI_REPLACE)); 889566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "B: leaf data (all)\n")); 899566063dSJacob Faibussowitsch PetscCall(PetscIntView(nldataB, ldB, PETSC_VIEWER_STDOUT_WORLD)); 90c4762a1bSJed Brown 919566063dSJacob Faibussowitsch PetscCall(PetscSFCompose(sfA, sfB, &sfBA)); 929566063dSJacob Faibussowitsch PetscCall(PetscSFSetFromOptions(sfBA)); 939566063dSJacob Faibussowitsch PetscCall(PetscSFSetUp(sfBA)); 949566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)sfBA, "sfBA")); 959566063dSJacob Faibussowitsch PetscCall(PetscSFViewFromOptions(sfBA, NULL, "-view")); 96c4762a1bSJed Brown 97c4762a1bSJed Brown for (i = 0; i < nldataB; i++) ldB[i] = -1; 989566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "BcastBA\n")); 999566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "BA: root data\n")); 1009566063dSJacob Faibussowitsch PetscCall(PetscIntView(nrootsA, rdA, PETSC_VIEWER_STDOUT_WORLD)); 1019566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sfBA, MPIU_INT, rdA, ldB,MPI_REPLACE)); 1029566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd(sfBA, MPIU_INT, rdA, ldB,MPI_REPLACE)); 1039566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "BA: leaf data (all)\n")); 1049566063dSJacob Faibussowitsch PetscCall(PetscIntView(nldataB, ldB, PETSC_VIEWER_STDOUT_WORLD)); 105c4762a1bSJed Brown 1069566063dSJacob Faibussowitsch PetscCall(PetscSFCreateInverseSF(sfA, &sfAm)); 1079566063dSJacob Faibussowitsch PetscCall(PetscSFSetFromOptions(sfAm)); 1089566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)sfAm, "sfAm")); 1099566063dSJacob Faibussowitsch PetscCall(PetscSFViewFromOptions(sfAm, NULL, "-view")); 110c4762a1bSJed Brown 111c4762a1bSJed Brown if (!inverse) { 1129566063dSJacob Faibussowitsch PetscCall(PetscSFComposeInverse(sfA, sfA, &sfAAm)); 113c4762a1bSJed Brown } else { 1149566063dSJacob Faibussowitsch PetscCall(PetscSFCompose(sfA, sfAm, &sfAAm)); 115c4762a1bSJed Brown } 1169566063dSJacob Faibussowitsch PetscCall(PetscSFSetFromOptions(sfAAm)); 1179566063dSJacob Faibussowitsch PetscCall(PetscSFSetUp(sfAAm)); 1189566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)sfAAm, "sfAAm")); 1199566063dSJacob Faibussowitsch PetscCall(PetscSFViewFromOptions(sfAAm, NULL, "-view")); 120c4762a1bSJed Brown 1219566063dSJacob Faibussowitsch PetscCall(PetscSFCreateInverseSF(sfB, &sfBm)); 1229566063dSJacob Faibussowitsch PetscCall(PetscSFSetFromOptions(sfBm)); 1239566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)sfBm, "sfBm")); 1249566063dSJacob Faibussowitsch PetscCall(PetscSFViewFromOptions(sfBm, NULL, "-view")); 125c4762a1bSJed Brown 126c4762a1bSJed Brown if (!inverse) { 1279566063dSJacob Faibussowitsch PetscCall(PetscSFComposeInverse(sfB, sfB, &sfBBm)); 128c4762a1bSJed Brown } else { 1299566063dSJacob Faibussowitsch PetscCall(PetscSFCompose(sfB, sfBm, &sfBBm)); 130c4762a1bSJed Brown } 1319566063dSJacob Faibussowitsch PetscCall(PetscSFSetFromOptions(sfBBm)); 1329566063dSJacob Faibussowitsch PetscCall(PetscSFSetUp(sfBBm)); 1339566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)sfBBm, "sfBBm")); 1349566063dSJacob Faibussowitsch PetscCall(PetscSFViewFromOptions(sfBBm, NULL, "-view")); 135c4762a1bSJed Brown 1369566063dSJacob Faibussowitsch PetscCall(PetscFree2(rdA, ldA)); 1379566063dSJacob Faibussowitsch PetscCall(PetscFree2(rdB, ldB)); 138c4762a1bSJed Brown 1399566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfA)); 1409566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfB)); 1419566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfBA)); 1429566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfAm)); 1439566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfBm)); 1449566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfAAm)); 1459566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfBBm)); 146c4762a1bSJed Brown 1479566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 148b122ec5aSJacob Faibussowitsch return 0; 149c4762a1bSJed Brown } 150c4762a1bSJed Brown 151c4762a1bSJed Brown /*TEST 152c4762a1bSJed Brown 153c4762a1bSJed Brown test: 154c4762a1bSJed Brown suffix: 1 155c4762a1bSJed Brown args: -view -explicit_inverse {{0 1}} 156c4762a1bSJed Brown 157c4762a1bSJed Brown test: 158c4762a1bSJed Brown nsize: 7 159c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 160c4762a1bSJed Brown suffix: 2 161c4762a1bSJed Brown args: -view -nl 5 -explicit_inverse {{0 1}} 162c4762a1bSJed Brown 163c4762a1bSJed Brown # we cannot test for -sf_window_flavor dynamic because SFCompose with sparse leaves may change the root data pointer only locally, and this is not supported by the dynamic case 164c4762a1bSJed Brown test: 165c4762a1bSJed Brown nsize: 7 166c4762a1bSJed Brown suffix: 2_window 167c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 168c4762a1bSJed Brown output_file: output/ex5_2.out 169c4762a1bSJed Brown args: -view -nl 5 -explicit_inverse {{0 1}} -sf_type window -sf_window_sync {{fence lock active}} -sf_window_flavor {{create allocate}} 170dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW) 171c4762a1bSJed Brown 172c4762a1bSJed Brown # The nightly test suite with MPICH uses ch3:sock, which is broken when winsize == 0 in some of the processes 173c4762a1bSJed Brown test: 174c4762a1bSJed Brown nsize: 7 175c4762a1bSJed Brown suffix: 2_window_shared 176c4762a1bSJed Brown filter: grep -v "type" | grep -v "sort" 177c4762a1bSJed Brown output_file: output/ex5_2.out 178c4762a1bSJed Brown args: -view -nl 5 -explicit_inverse {{0 1}} -sf_type window -sf_window_sync {{fence lock active}} -sf_window_flavor shared 179dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY) !defined(PETSC_HAVE_MPICH_NUMVERSION) defined(PETSC_HAVE_MPI_ONE_SIDED) 180c4762a1bSJed Brown 181c4762a1bSJed Brown TEST*/ 182