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