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