xref: /petsc/src/vec/is/sf/tutorials/ex3.c (revision a2dece7af36d571dbc9b13292ed2989a59daf513)
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 
main(int argc,char ** argv)7d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
8d71ae5a4SJacob Faibussowitsch {
9c4762a1bSJed Brown   PetscSF      sf;
10c4762a1bSJed Brown   Vec          A, Aout;
11c4762a1bSJed Brown   PetscScalar *bufA;
12c4762a1bSJed Brown   PetscScalar *bufAout;
13c4762a1bSJed Brown   PetscMPIInt  rank, size;
14c4762a1bSJed Brown   PetscInt     nroots, nleaves;
15c4762a1bSJed Brown   PetscInt     i;
16c4762a1bSJed Brown   PetscInt    *ilocal;
17c4762a1bSJed Brown   PetscSFNode *iremote;
18c4762a1bSJed Brown   PetscBool    test_dupped_type;
19c4762a1bSJed Brown   MPI_Datatype contig;
20c4762a1bSJed Brown 
21327415f7SBarry Smith   PetscFunctionBeginUser;
229566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
239566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
249566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
25c4762a1bSJed Brown 
26*bd158744SPierre Jolivet   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only coded for one MPI process");
27c4762a1bSJed Brown 
28d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, "", "PetscSF type freeing options", "none");
29c4762a1bSJed Brown   test_dupped_type = PETSC_FALSE;
309566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_dupped_type", "Test dupped input type", "", test_dupped_type, &test_dupped_type, NULL));
31d0609cedSBarry Smith   PetscOptionsEnd();
32c4762a1bSJed Brown 
339566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(PETSC_COMM_WORLD, &sf));
349566063dSJacob Faibussowitsch   PetscCall(PetscSFSetFromOptions(sf));
35c4762a1bSJed Brown 
36c4762a1bSJed Brown   nleaves = 1;
37c4762a1bSJed Brown   nroots  = 1;
389566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nleaves, &ilocal));
39c4762a1bSJed Brown 
40ad540459SPierre Jolivet   for (i = 0; i < nleaves; i++) ilocal[i] = i;
41c4762a1bSJed Brown 
429566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nleaves, &iremote));
43c4762a1bSJed Brown   iremote[0].rank  = 0;
44c4762a1bSJed Brown   iremote[0].index = 0;
459566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(sf, nroots, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER));
469566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sf));
479566063dSJacob Faibussowitsch   PetscCall(PetscSFView(sf, PETSC_VIEWER_STDOUT_WORLD));
489566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &A));
499566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(A, 4, PETSC_DETERMINE));
509566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(A));
519566063dSJacob Faibussowitsch   PetscCall(VecSetUp(A));
52c4762a1bSJed Brown 
539566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(A, &Aout));
549566063dSJacob Faibussowitsch   PetscCall(VecGetArray(A, &bufA));
55ad540459SPierre Jolivet   for (i = 0; i < 4; i++) bufA[i] = (PetscScalar)i;
569566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(A, &bufA));
57c4762a1bSJed Brown 
589566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(A, (const PetscScalar **)&bufA));
599566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Aout, &bufAout));
60c4762a1bSJed Brown 
619566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_contiguous(4, MPIU_SCALAR, &contig));
629566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_commit(&contig));
63c4762a1bSJed Brown 
64c4762a1bSJed Brown   if (test_dupped_type) {
65c4762a1bSJed Brown     MPI_Datatype tmp;
669566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_dup(contig, &tmp));
679566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_free(&contig));
68c4762a1bSJed Brown     contig = tmp;
69c4762a1bSJed Brown   }
708e54d7e8SToby Isaac   PetscCall(PetscSFRegisterPersistent(sf, contig, bufA, bufAout));
71c4762a1bSJed Brown   for (i = 0; i < 10000; i++) {
729566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(sf, contig, bufA, bufAout, MPI_REPLACE));
739566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(sf, contig, bufA, bufAout, MPI_REPLACE));
74c4762a1bSJed Brown   }
758e54d7e8SToby Isaac   PetscCall(PetscSFDeregisterPersistent(sf, contig, bufA, bufAout));
769566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(A, (const PetscScalar **)&bufA));
779566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Aout, &bufAout));
78c4762a1bSJed Brown 
799566063dSJacob Faibussowitsch   PetscCall(VecView(Aout, PETSC_VIEWER_STDOUT_WORLD));
809566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&A));
819566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Aout));
829566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sf));
839566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_free(&contig));
849566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
85b122ec5aSJacob Faibussowitsch   return 0;
86c4762a1bSJed Brown }
87c4762a1bSJed Brown 
88c4762a1bSJed Brown /*TEST
89c4762a1bSJed Brown 
90c4762a1bSJed Brown    test:
91c4762a1bSJed Brown       suffix: basic
92c4762a1bSJed Brown       args: -sf_type basic
93c4762a1bSJed Brown 
94c4762a1bSJed Brown    test:
95c4762a1bSJed Brown       suffix: basic_dupped
96c4762a1bSJed Brown       args: -test_dupped_type -sf_type basic
97c4762a1bSJed Brown 
98c4762a1bSJed Brown    test:
99c4762a1bSJed Brown       suffix: window
100c4762a1bSJed Brown       filter: grep -v "type" | grep -v "sort"
101c4762a1bSJed Brown       args: -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create allocate dynamic}}
102dfd57a17SPierre Jolivet       requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
103c4762a1bSJed Brown 
104c4762a1bSJed Brown    test:
105c4762a1bSJed Brown       suffix: window_dupped
106c4762a1bSJed Brown       filter: grep -v "type" | grep -v "sort"
107c4762a1bSJed Brown       args: -test_dupped_type -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor {{create allocate dynamic}}
108dfd57a17SPierre Jolivet       requires: defined(PETSC_HAVE_MPI_ONE_SIDED) defined(PETSC_HAVE_MPI_FEATURE_DYNAMIC_WINDOW)
109c4762a1bSJed Brown 
110c4762a1bSJed Brown    test:
111c4762a1bSJed Brown       suffix: window_shared
112c4762a1bSJed Brown       output_file: output/ex3_window.out
113c4762a1bSJed Brown       filter: grep -v "type" | grep -v "sort"
114c4762a1bSJed Brown       args: -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor shared
1150f85934cSJunchao Zhang       requires: defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY) defined(PETSC_HAVE_MPI_ONE_SIDED) !defined(PETSC_HAVE_I_MPI)
116c4762a1bSJed Brown 
117c4762a1bSJed Brown    test:
118c4762a1bSJed Brown       suffix: window_dupped_shared
119c4762a1bSJed Brown       output_file: output/ex3_window_dupped.out
120c4762a1bSJed Brown       filter: grep -v "type" | grep -v "sort"
121c4762a1bSJed Brown       args: -test_dupped_type -sf_type window -sf_window_sync {{fence active lock}} -sf_window_flavor shared
1220f85934cSJunchao Zhang       requires: defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY) defined(PETSC_HAVE_MPI_ONE_SIDED) !defined(PETSC_HAVE_I_MPI)
123c4762a1bSJed Brown 
124c4762a1bSJed Brown TEST*/
125