xref: /petsc/src/dm/tests/ex8.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerGetFormat(viewer,&format));
13*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERGLVIS,&isglvis));
14*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
15c4762a1bSJed Brown   if (isglvis) {
16c4762a1bSJed Brown     DM dm;
17c4762a1bSJed Brown 
18*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetDM(v,&dm));
19c4762a1bSJed Brown     /* DMView() cannot be tested, as DMView_Shell defaults to VecView */
20c4762a1bSJed Brown     if (!dm) PetscFunctionReturn(0);
21*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView_GLVis(v,viewer));
22c4762a1bSJed Brown   } else if (isascii) {
23c4762a1bSJed Brown     const char* name;
24c4762a1bSJed Brown     PetscInt    n;
25c4762a1bSJed Brown 
26*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetLocalSize(v,&n));
27*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectGetName((PetscObject)v,&name));
28c4762a1bSJed Brown     if (!PetscGlobalRank) {
29*5f80ce2aSJacob Faibussowitsch       CHKERRQ(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*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(dm,&V));
44*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)V,"sample"));
45*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerGLVisSetFields(viewer,1,&fec_type,&dim,NULL,(PetscObject*)&V,NULL,NULL));
46*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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   PetscErrorCode ierr;
55c4762a1bSJed Brown 
56c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
57*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMShellCreate(PETSC_COMM_WORLD,&dm));
58*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)dm,"DMSetUpGLVisViewer_C",DMSetUpGLVisViewer_Shell));
59*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateMPI(PETSC_COMM_WORLD,1,PETSC_DECIDE,&v));
60*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject)v,"seed"));
61*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetOperation(v,VECOP_VIEW,(void (*)(void))VecView_Shell));
62*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMShellSetGlobalVector(dm,v));
63*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&v));
64*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMViewFromOptions(dm,NULL,"-dm_view"));
65*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetGlobalVector(dm,&v));
66*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecViewFromOptions(v,NULL,"-vec_view"));
67*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreGlobalVector(dm,&v));
68*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dm));
69c4762a1bSJed Brown   ierr = PetscFinalize();
70c4762a1bSJed Brown   return ierr;
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