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