xref: /petsc/src/vec/is/sf/tests/ex10.c (revision 2eab0914758aeb8383261f911e918804bd166472)
1*14357523SPierre Jolivet static char help[] = "Test PetscSFCompose() against some corner cases \n\n";
20ea77edaSksagiyam 
30ea77edaSksagiyam #include <petscsf.h>
40ea77edaSksagiyam 
main(int argc,char ** argv)5d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
6d71ae5a4SJacob Faibussowitsch {
70ea77edaSksagiyam   PetscMPIInt  size;
80ea77edaSksagiyam   PetscSF      sfA0, sfA1, sfA2, sfB;
90ea77edaSksagiyam   PetscInt     nroots, nleaves;
100ea77edaSksagiyam   PetscInt    *ilocalA0, *ilocalA1, *ilocalA2, *ilocalB;
110ea77edaSksagiyam   PetscSFNode *iremoteA0, *iremoteA1, *iremoteA2, *iremoteB;
120ea77edaSksagiyam 
130ea77edaSksagiyam   PetscFunctionBeginUser;
140ea77edaSksagiyam   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
150ea77edaSksagiyam   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
16bd158744SPierre Jolivet   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only coded for one MPI process");
170ea77edaSksagiyam   PetscCall(PetscSFCreate(PETSC_COMM_WORLD, &sfA0));
180ea77edaSksagiyam   PetscCall(PetscSFCreate(PETSC_COMM_WORLD, &sfA1));
190ea77edaSksagiyam   PetscCall(PetscSFCreate(PETSC_COMM_WORLD, &sfA2));
200ea77edaSksagiyam   PetscCall(PetscSFCreate(PETSC_COMM_WORLD, &sfB));
210ea77edaSksagiyam   /* sfA0 */
220ea77edaSksagiyam   nroots  = 1;
230ea77edaSksagiyam   nleaves = 0;
240ea77edaSksagiyam   PetscCall(PetscMalloc1(nleaves, &ilocalA0));
250ea77edaSksagiyam   PetscCall(PetscMalloc1(nleaves, &iremoteA0));
260ea77edaSksagiyam   PetscCall(PetscSFSetGraph(sfA0, nroots, nleaves, ilocalA0, PETSC_OWN_POINTER, iremoteA0, PETSC_OWN_POINTER));
270ea77edaSksagiyam   PetscCall(PetscSFSetUp(sfA0));
280ea77edaSksagiyam   PetscCall(PetscObjectSetName((PetscObject)sfA0, "sfA0"));
290ea77edaSksagiyam   PetscCall(PetscSFView(sfA0, NULL));
300ea77edaSksagiyam   /* sfA1 */
310ea77edaSksagiyam   nroots  = 1;
320ea77edaSksagiyam   nleaves = 1;
330ea77edaSksagiyam   PetscCall(PetscMalloc1(nleaves, &ilocalA1));
340ea77edaSksagiyam   PetscCall(PetscMalloc1(nleaves, &iremoteA1));
350ea77edaSksagiyam   ilocalA1[0]        = 1;
360ea77edaSksagiyam   iremoteA1[0].rank  = 0;
370ea77edaSksagiyam   iremoteA1[0].index = 0;
380ea77edaSksagiyam   PetscCall(PetscSFSetGraph(sfA1, nroots, nleaves, ilocalA1, PETSC_OWN_POINTER, iremoteA1, PETSC_OWN_POINTER));
390ea77edaSksagiyam   PetscCall(PetscSFSetUp(sfA1));
400ea77edaSksagiyam   PetscCall(PetscObjectSetName((PetscObject)sfA1, "sfA1"));
410ea77edaSksagiyam   PetscCall(PetscSFView(sfA1, NULL));
420ea77edaSksagiyam   /* sfA2 */
430ea77edaSksagiyam   nroots  = 1;
440ea77edaSksagiyam   nleaves = 1;
450ea77edaSksagiyam   PetscCall(PetscMalloc1(nleaves, &ilocalA2));
460ea77edaSksagiyam   PetscCall(PetscMalloc1(nleaves, &iremoteA2));
470ea77edaSksagiyam   ilocalA2[0]        = 0;
480ea77edaSksagiyam   iremoteA2[0].rank  = 0;
490ea77edaSksagiyam   iremoteA2[0].index = 0;
500ea77edaSksagiyam   PetscCall(PetscSFSetGraph(sfA2, nroots, nleaves, ilocalA2, PETSC_OWN_POINTER, iremoteA2, PETSC_OWN_POINTER));
510ea77edaSksagiyam   PetscCall(PetscSFSetUp(sfA2));
520ea77edaSksagiyam   PetscCall(PetscObjectSetName((PetscObject)sfA2, "sfA2"));
530ea77edaSksagiyam   PetscCall(PetscSFView(sfA2, NULL));
540ea77edaSksagiyam   /* sfB */
550ea77edaSksagiyam   nroots  = 2;
560ea77edaSksagiyam   nleaves = 2;
570ea77edaSksagiyam   PetscCall(PetscMalloc1(nleaves, &ilocalB));
580ea77edaSksagiyam   PetscCall(PetscMalloc1(nleaves, &iremoteB));
590ea77edaSksagiyam   ilocalB[0]        = 100;
600ea77edaSksagiyam   iremoteB[0].rank  = 0;
610ea77edaSksagiyam   iremoteB[0].index = 0;
620ea77edaSksagiyam   ilocalB[1]        = 101;
630ea77edaSksagiyam   iremoteB[1].rank  = 0;
640ea77edaSksagiyam   iremoteB[1].index = 1;
650ea77edaSksagiyam   PetscCall(PetscSFSetGraph(sfB, nroots, nleaves, ilocalB, PETSC_OWN_POINTER, iremoteB, PETSC_OWN_POINTER));
660ea77edaSksagiyam   PetscCall(PetscSFSetUp(sfB));
670ea77edaSksagiyam   PetscCall(PetscObjectSetName((PetscObject)sfB, "sfB"));
680ea77edaSksagiyam   PetscCall(PetscSFView(sfB, NULL));
690ea77edaSksagiyam   /* Test 0 */
700ea77edaSksagiyam   {
710ea77edaSksagiyam     PetscSF sfC;
720ea77edaSksagiyam 
730ea77edaSksagiyam     PetscCall(PetscSFCompose(sfA0, sfB, &sfC));
740ea77edaSksagiyam     PetscCall(PetscObjectSetName((PetscObject)sfC, "PetscSFCompose(sfA0, sfB)"));
750ea77edaSksagiyam     PetscCall(PetscSFView(sfC, NULL));
760ea77edaSksagiyam     PetscCall(PetscSFDestroy(&sfC));
770ea77edaSksagiyam   }
780ea77edaSksagiyam   /* Test 1 */
790ea77edaSksagiyam   {
800ea77edaSksagiyam     PetscSF sfC;
810ea77edaSksagiyam 
820ea77edaSksagiyam     PetscCall(PetscSFCompose(sfA1, sfB, &sfC));
830ea77edaSksagiyam     PetscCall(PetscObjectSetName((PetscObject)sfC, "PetscSFCompose(sfA1, sfB)"));
840ea77edaSksagiyam     PetscCall(PetscSFView(sfC, NULL));
850ea77edaSksagiyam     PetscCall(PetscSFDestroy(&sfC));
860ea77edaSksagiyam   }
870ea77edaSksagiyam   /* Test 2 */
880ea77edaSksagiyam   {
890ea77edaSksagiyam     PetscSF sfC;
900ea77edaSksagiyam 
910ea77edaSksagiyam     PetscCall(PetscSFCompose(sfA2, sfB, &sfC));
920ea77edaSksagiyam     PetscCall(PetscObjectSetName((PetscObject)sfC, "PetscSFCompose(sfA2, sfB)"));
930ea77edaSksagiyam     PetscCall(PetscSFView(sfC, NULL));
940ea77edaSksagiyam     PetscCall(PetscSFDestroy(&sfC));
950ea77edaSksagiyam   }
960ea77edaSksagiyam   PetscCall(PetscSFDestroy(&sfA0));
970ea77edaSksagiyam   PetscCall(PetscSFDestroy(&sfA1));
980ea77edaSksagiyam   PetscCall(PetscSFDestroy(&sfA2));
990ea77edaSksagiyam   PetscCall(PetscSFDestroy(&sfB));
1000ea77edaSksagiyam   PetscCall(PetscFinalize());
1010ea77edaSksagiyam   return 0;
1020ea77edaSksagiyam }
1030ea77edaSksagiyam 
1040ea77edaSksagiyam /*TEST
1050ea77edaSksagiyam 
1060ea77edaSksagiyam    test:
1070ea77edaSksagiyam      suffix: 0
1080ea77edaSksagiyam 
1090ea77edaSksagiyam TEST*/
110