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