xref: /petsc/src/dm/tests/ex8.c (revision 2e956fe4fc852fabc23b437482e1fb7b77fddb0d)
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;
129566063dSJacob Faibussowitsch   PetscCall(PetscViewerGetFormat(viewer,&format));
139566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERGLVIS,&isglvis));
149566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
15c4762a1bSJed Brown   if (isglvis) {
16c4762a1bSJed Brown     DM dm;
17c4762a1bSJed Brown 
189566063dSJacob Faibussowitsch     PetscCall(VecGetDM(v,&dm));
19c4762a1bSJed Brown     /* DMView() cannot be tested, as DMView_Shell defaults to VecView */
20c4762a1bSJed Brown     if (!dm) PetscFunctionReturn(0);
219566063dSJacob Faibussowitsch     PetscCall(VecView_GLVis(v,viewer));
22c4762a1bSJed Brown   } else if (isascii) {
23c4762a1bSJed Brown     const char* name;
24c4762a1bSJed Brown     PetscInt    n;
25c4762a1bSJed Brown 
269566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(v,&n));
279566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetName((PetscObject)v,&name));
28c4762a1bSJed Brown     if (!PetscGlobalRank) {
2963a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"Hello from rank 0 -> vector name %s, size %" PetscInt_FMT "\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;
439566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dm,&V));
449566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)V,"sample"));
459566063dSJacob Faibussowitsch   PetscCall(PetscViewerGLVisSetFields(viewer,1,&fec_type,&dim,NULL,(PetscObject*)&V,NULL,NULL));
469566063dSJacob 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 
559566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,NULL,help));
569566063dSJacob Faibussowitsch   PetscCall(DMShellCreate(PETSC_COMM_WORLD,&dm));
579566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)dm,"DMSetUpGLVisViewer_C",DMSetUpGLVisViewer_Shell));
589566063dSJacob Faibussowitsch   PetscCall(VecCreateMPI(PETSC_COMM_WORLD,1,PETSC_DECIDE,&v));
599566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)v,"seed"));
609566063dSJacob Faibussowitsch   PetscCall(VecSetOperation(v,VECOP_VIEW,(void (*)(void))VecView_Shell));
619566063dSJacob Faibussowitsch   PetscCall(DMShellSetGlobalVector(dm,v));
629566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&v));
639566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm,NULL,"-dm_view"));
649566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dm,&v));
659566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(v,NULL,"-vec_view"));
669566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dm,&v));
67*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)dm,"DMSetUpGLVisViewer_C",NULL));
689566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
699566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
70b122ec5aSJacob Faibussowitsch   return 0;
71c4762a1bSJed Brown }
72c4762a1bSJed Brown 
73c4762a1bSJed Brown /*TEST
74c4762a1bSJed Brown 
75c4762a1bSJed Brown   test:
76c4762a1bSJed Brown     suffix: glvis_par
77c4762a1bSJed Brown     nsize: {{1 2}}
78c4762a1bSJed Brown     args: -dm_view glvis: -vec_view glvis:
79c4762a1bSJed Brown     output_file: output/ex8_glvis.out
80c4762a1bSJed Brown 
81c4762a1bSJed Brown TEST*/
82