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