xref: /petsc/src/dm/tests/ex8.c (revision 58d68138c660dfb4e9f5b03334792cd4f2ffd7cc)
1 static char help[] = "Test parallel ruotines for GLVis\n\n";
2 
3 #include <petscdmshell.h>
4 #include <petsc/private/glvisvecimpl.h>
5 
6 PetscErrorCode VecView_Shell(Vec v, PetscViewer viewer) {
7   PetscViewerFormat format;
8   PetscBool         isglvis, isascii;
9 
10   PetscFunctionBegin;
11   PetscCall(PetscViewerGetFormat(viewer, &format));
12   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERGLVIS, &isglvis));
13   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
14   if (isglvis) {
15     DM dm;
16 
17     PetscCall(VecGetDM(v, &dm));
18     /* DMView() cannot be tested, as DMView_Shell defaults to VecView */
19     if (!dm) PetscFunctionReturn(0);
20     PetscCall(VecView_GLVis(v, viewer));
21   } else if (isascii) {
22     const char *name;
23     PetscInt    n;
24 
25     PetscCall(VecGetLocalSize(v, &n));
26     PetscCall(PetscObjectGetName((PetscObject)v, &name));
27     if (!PetscGlobalRank) { PetscCall(PetscViewerASCIIPrintf(viewer, "Hello from rank 0 -> vector name %s, size %" PetscInt_FMT "\n", name, n)); }
28   }
29   PetscFunctionReturn(0);
30 }
31 
32 PetscErrorCode DMSetUpGLVisViewer_Shell(PetscObject odm, PetscViewer viewer) {
33   DM          dm = (DM)odm;
34   Vec         V;
35   PetscInt    dim      = 2;
36   const char *fec_type = {"testme"};
37 
38   PetscFunctionBegin;
39   PetscCall(DMCreateGlobalVector(dm, &V));
40   PetscCall(PetscObjectSetName((PetscObject)V, "sample"));
41   PetscCall(PetscViewerGLVisSetFields(viewer, 1, &fec_type, &dim, NULL, (PetscObject *)&V, NULL, NULL));
42   PetscCall(VecDestroy(&V));
43   PetscFunctionReturn(0);
44 }
45 
46 int main(int argc, char **argv) {
47   DM  dm;
48   Vec v;
49 
50   PetscFunctionBeginUser;
51   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
52   PetscCall(DMShellCreate(PETSC_COMM_WORLD, &dm));
53   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMSetUpGLVisViewer_C", DMSetUpGLVisViewer_Shell));
54   PetscCall(VecCreateMPI(PETSC_COMM_WORLD, 1, PETSC_DECIDE, &v));
55   PetscCall(PetscObjectSetName((PetscObject)v, "seed"));
56   PetscCall(VecSetOperation(v, VECOP_VIEW, (void (*)(void))VecView_Shell));
57   PetscCall(DMShellSetGlobalVector(dm, v));
58   PetscCall(VecDestroy(&v));
59   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
60   PetscCall(DMGetGlobalVector(dm, &v));
61   PetscCall(VecViewFromOptions(v, NULL, "-vec_view"));
62   PetscCall(DMRestoreGlobalVector(dm, &v));
63   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMSetUpGLVisViewer_C", NULL));
64   PetscCall(DMDestroy(&dm));
65   PetscCall(PetscFinalize());
66   return 0;
67 }
68 
69 /*TEST
70 
71   test:
72     suffix: glvis_par
73     nsize: {{1 2}}
74     args: -dm_view glvis: -vec_view glvis:
75     output_file: output/ex8_glvis.out
76 
77 TEST*/
78