xref: /petsc/src/dm/tests/ex8.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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 
69371c9d4SSatish Balay PetscErrorCode VecView_Shell(Vec v, PetscViewer viewer) {
7c4762a1bSJed Brown   PetscViewerFormat format;
8c4762a1bSJed Brown   PetscBool         isglvis, isascii;
9c4762a1bSJed Brown 
10c4762a1bSJed Brown   PetscFunctionBegin;
119566063dSJacob Faibussowitsch   PetscCall(PetscViewerGetFormat(viewer, &format));
129566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERGLVIS, &isglvis));
139566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
14c4762a1bSJed Brown   if (isglvis) {
15c4762a1bSJed Brown     DM dm;
16c4762a1bSJed Brown 
179566063dSJacob Faibussowitsch     PetscCall(VecGetDM(v, &dm));
18c4762a1bSJed Brown     /* DMView() cannot be tested, as DMView_Shell defaults to VecView */
19c4762a1bSJed Brown     if (!dm) PetscFunctionReturn(0);
209566063dSJacob Faibussowitsch     PetscCall(VecView_GLVis(v, viewer));
21c4762a1bSJed Brown   } else if (isascii) {
22c4762a1bSJed Brown     const char *name;
23c4762a1bSJed Brown     PetscInt    n;
24c4762a1bSJed Brown 
259566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(v, &n));
269566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetName((PetscObject)v, &name));
27*48a46eb9SPierre Jolivet     if (!PetscGlobalRank) PetscCall(PetscViewerASCIIPrintf(viewer, "Hello from rank 0 -> vector name %s, size %" PetscInt_FMT "\n", name, n));
28c4762a1bSJed Brown   }
29c4762a1bSJed Brown   PetscFunctionReturn(0);
30c4762a1bSJed Brown }
31c4762a1bSJed Brown 
329371c9d4SSatish Balay PetscErrorCode DMSetUpGLVisViewer_Shell(PetscObject odm, PetscViewer viewer) {
33c4762a1bSJed Brown   DM          dm = (DM)odm;
34c4762a1bSJed Brown   Vec         V;
35c4762a1bSJed Brown   PetscInt    dim      = 2;
36c4762a1bSJed Brown   const char *fec_type = {"testme"};
37c4762a1bSJed Brown 
38c4762a1bSJed Brown   PetscFunctionBegin;
399566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dm, &V));
409566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)V, "sample"));
419566063dSJacob Faibussowitsch   PetscCall(PetscViewerGLVisSetFields(viewer, 1, &fec_type, &dim, NULL, (PetscObject *)&V, NULL, NULL));
429566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&V));
43c4762a1bSJed Brown   PetscFunctionReturn(0);
44c4762a1bSJed Brown }
45c4762a1bSJed Brown 
469371c9d4SSatish Balay int main(int argc, char **argv) {
47c4762a1bSJed Brown   DM  dm;
48c4762a1bSJed Brown   Vec v;
49c4762a1bSJed Brown 
50327415f7SBarry Smith   PetscFunctionBeginUser;
519566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
529566063dSJacob Faibussowitsch   PetscCall(DMShellCreate(PETSC_COMM_WORLD, &dm));
539566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMSetUpGLVisViewer_C", DMSetUpGLVisViewer_Shell));
549566063dSJacob Faibussowitsch   PetscCall(VecCreateMPI(PETSC_COMM_WORLD, 1, PETSC_DECIDE, &v));
559566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)v, "seed"));
569566063dSJacob Faibussowitsch   PetscCall(VecSetOperation(v, VECOP_VIEW, (void (*)(void))VecView_Shell));
579566063dSJacob Faibussowitsch   PetscCall(DMShellSetGlobalVector(dm, v));
589566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&v));
599566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
609566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dm, &v));
619566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(v, NULL, "-vec_view"));
629566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dm, &v));
632e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMSetUpGLVisViewer_C", NULL));
649566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
659566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
66b122ec5aSJacob Faibussowitsch   return 0;
67c4762a1bSJed Brown }
68c4762a1bSJed Brown 
69c4762a1bSJed Brown /*TEST
70c4762a1bSJed Brown 
71c4762a1bSJed Brown   test:
72c4762a1bSJed Brown     suffix: glvis_par
73c4762a1bSJed Brown     nsize: {{1 2}}
74c4762a1bSJed Brown     args: -dm_view glvis: -vec_view glvis:
75c4762a1bSJed Brown     output_file: output/ex8_glvis.out
76c4762a1bSJed Brown 
77c4762a1bSJed Brown TEST*/
78