xref: /petsc/src/vec/is/sf/tutorials/ex3.c (revision 11486bccf1aeea1ca5536228f99d437b39bdaca6)
1 static const char help[] = "Test freeing of MPI types in PetscSF\n\n";
2 
3 #include <petscvec.h>
4 #include <petscsf.h>
5 #include <petscviewer.h>
6 
7 int main(int argc, char **argv)
8 {
9   PetscSF     sf;
10   Vec         A,Aout;
11   PetscScalar *bufA;
12   PetscScalar *bufAout;
13   PetscMPIInt rank, size;
14   PetscInt    nroots, nleaves;
15   PetscInt    i;
16   PetscInt    *ilocal;
17   PetscSFNode *iremote;
18   PetscBool   test_dupped_type;
19   MPI_Datatype contig;
20 
21   PetscFunctionBeginUser;
22   PetscCall(PetscInitialize(&argc,&argv,NULL,help));
23   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
24   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
25 
26   PetscCheck(size == 1,PETSC_COMM_WORLD, PETSC_ERR_USER, "Only coded for one MPI process");
27 
28   PetscOptionsBegin(PETSC_COMM_WORLD,"","PetscSF type freeing options","none");
29   test_dupped_type = PETSC_FALSE;
30   PetscCall(PetscOptionsBool("-test_dupped_type", "Test dupped input type","",test_dupped_type,&test_dupped_type,NULL));
31   PetscOptionsEnd();
32 
33   PetscCall(PetscSFCreate(PETSC_COMM_WORLD,&sf));
34   PetscCall(PetscSFSetFromOptions(sf));
35 
36   nleaves = 1;
37   nroots = 1;
38   PetscCall(PetscMalloc1(nleaves,&ilocal));
39 
40   for (i = 0; i<nleaves; i++) {
41     ilocal[i] = i;
42   }
43 
44   PetscCall(PetscMalloc1(nleaves,&iremote));
45   iremote[0].rank = 0;
46   iremote[0].index = 0;
47   PetscCall(PetscSFSetGraph(sf,nroots,nleaves,ilocal,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER));
48   PetscCall(PetscSFSetUp(sf));
49   PetscCall(PetscSFView(sf,PETSC_VIEWER_STDOUT_WORLD));
50   PetscCall(VecCreate(PETSC_COMM_WORLD,&A));
51   PetscCall(VecSetSizes(A,4,PETSC_DETERMINE));
52   PetscCall(VecSetFromOptions(A));
53   PetscCall(VecSetUp(A));
54 
55   PetscCall(VecDuplicate(A,&Aout));
56   PetscCall(VecGetArray(A,&bufA));
57   for (i=0; i<4; i++) {
58     bufA[i] = (PetscScalar)i;
59   }
60   PetscCall(VecRestoreArray(A,&bufA));
61 
62   PetscCall(VecGetArrayRead(A,(const PetscScalar**)&bufA));
63   PetscCall(VecGetArray(Aout,&bufAout));
64 
65   PetscCallMPI(MPI_Type_contiguous(4, MPIU_SCALAR, &contig));
66   PetscCallMPI(MPI_Type_commit(&contig));
67 
68   if (test_dupped_type) {
69     MPI_Datatype tmp;
70     PetscCallMPI(MPI_Type_dup(contig, &tmp));
71     PetscCallMPI(MPI_Type_free(&contig));
72     contig = tmp;
73   }
74   for (i=0;i<10000;i++) {
75     PetscCall(PetscSFBcastBegin(sf,contig,bufA,bufAout,MPI_REPLACE));
76     PetscCall(PetscSFBcastEnd(sf,contig,bufA,bufAout,MPI_REPLACE));
77   }
78   PetscCall(VecRestoreArrayRead(A,(const PetscScalar**)&bufA));
79   PetscCall(VecRestoreArray(Aout,&bufAout));
80 
81   PetscCall(VecView(Aout,PETSC_VIEWER_STDOUT_WORLD));
82   PetscCall(VecDestroy(&A));
83   PetscCall(VecDestroy(&Aout));
84   PetscCall(PetscSFDestroy(&sf));
85   PetscCallMPI(MPI_Type_free(&contig));
86   PetscCall(PetscFinalize());
87   return 0;
88 }
89 
90 /*TEST
91 
92    test:
93       suffix: basic
94       args: -sf_type basic
95 
96    test:
97       suffix: basic_dupped
98       args: -test_dupped_type -sf_type basic
99 
100    test:
101       suffix: window
102       filter: grep -v "type" | grep -v "sort"
103       args: -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create allocate dynamic}}
104       requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
105 
106    test:
107       suffix: window_dupped
108       filter: grep -v "type" | grep -v "sort"
109       args: -test_dupped_type -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create allocate dynamic}}
110       requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
111 
112    test:
113       suffix: window_shared
114       output_file: output/ex3_window.out
115       filter: grep -v "type" | grep -v "sort"
116       args: -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor shared
117       requires: defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY) defined(PETSC_HAVE_MPI_ONE_SIDED) !defined(PETSC_HAVE_I_MPI_NUMVERSION)
118 
119    test:
120       suffix: window_dupped_shared
121       output_file: output/ex3_window_dupped.out
122       filter: grep -v "type" | grep -v "sort"
123       args: -test_dupped_type -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor shared
124       requires: defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY) defined(PETSC_HAVE_MPI_ONE_SIDED) !defined(PETSC_HAVE_I_MPI_NUMVERSION)
125 
126 TEST*/
127