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