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 6d71ae5a4SJacob Faibussowitsch PetscErrorCode VecView_Shell(Vec v, PetscViewer viewer) 7d71ae5a4SJacob Faibussowitsch { 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 */ 203ba16761SJacob Faibussowitsch if (!dm) PetscFunctionReturn(PETSC_SUCCESS); 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)); 2848a46eb9SPierre Jolivet if (!PetscGlobalRank) PetscCall(PetscViewerASCIIPrintf(viewer, "Hello from rank 0 -> vector name %s, size %" PetscInt_FMT "\n", name, n)); 29c4762a1bSJed Brown } 303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 31c4762a1bSJed Brown } 32c4762a1bSJed Brown 33d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSetUpGLVisViewer_Shell(PetscObject odm, PetscViewer viewer) 34d71ae5a4SJacob Faibussowitsch { 35c4762a1bSJed Brown DM dm = (DM)odm; 36c4762a1bSJed Brown Vec V; 37c4762a1bSJed Brown PetscInt dim = 2; 38c4762a1bSJed Brown const char *fec_type = {"testme"}; 39c4762a1bSJed Brown 40c4762a1bSJed Brown PetscFunctionBegin; 419566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &V)); 429566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)V, "sample")); 439566063dSJacob Faibussowitsch PetscCall(PetscViewerGLVisSetFields(viewer, 1, &fec_type, &dim, NULL, (PetscObject *)&V, NULL, NULL)); 449566063dSJacob Faibussowitsch PetscCall(VecDestroy(&V)); 453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 46c4762a1bSJed Brown } 47c4762a1bSJed Brown 48d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 49d71ae5a4SJacob Faibussowitsch { 50c4762a1bSJed Brown DM dm; 51c4762a1bSJed Brown Vec v; 52c4762a1bSJed Brown 53327415f7SBarry Smith PetscFunctionBeginUser; 549566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 559566063dSJacob Faibussowitsch PetscCall(DMShellCreate(PETSC_COMM_WORLD, &dm)); 569566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMSetUpGLVisViewer_C", DMSetUpGLVisViewer_Shell)); 5777433607SBarry Smith PetscCall(VecCreateFromOptions(PETSC_COMM_WORLD, NULL, 1, 1, PETSC_DECIDE, &v)); 589566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)v, "seed")); 59*57d50842SBarry Smith PetscCall(VecSetOperation(v, VECOP_VIEW, (PetscErrorCodeFn *)VecView_Shell)); 609566063dSJacob Faibussowitsch PetscCall(DMShellSetGlobalVector(dm, v)); 619566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v)); 629566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 639566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &v)); 649566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(v, NULL, "-vec_view")); 659566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &v)); 662e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMSetUpGLVisViewer_C", NULL)); 679566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 689566063dSJacob 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