xref: /petsc/src/dm/impls/da/grglvis.c (revision c6a7a37075f8bf8d34d92c4910d42445b7a3482d)
18135c375SStefano Zampini /* Routines to visualize DMDAs and fields through GLVis */
28135c375SStefano Zampini 
38135c375SStefano Zampini #include <petsc/private/dmdaimpl.h>
48135c375SStefano Zampini #include <petsc/private/glvisviewerimpl.h>
58135c375SStefano Zampini 
68135c375SStefano Zampini typedef struct {
746900a7fSStefano Zampini   PetscBool ll;
846900a7fSStefano Zampini } DMDAGhostedGLVisViewerCtx;
98135c375SStefano Zampini 
10d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDAGhostedDestroyGLVisViewerCtx_Private(void **vctx)
11d71ae5a4SJacob Faibussowitsch {
1246900a7fSStefano Zampini   PetscFunctionBegin;
139566063dSJacob Faibussowitsch   PetscCall(PetscFree(*vctx));
143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1546900a7fSStefano Zampini }
1646900a7fSStefano Zampini 
1746900a7fSStefano Zampini typedef struct {
1846900a7fSStefano Zampini   Vec xlocal;
1946900a7fSStefano Zampini } DMDAFieldGLVisViewerCtx;
2046900a7fSStefano Zampini 
21d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDAFieldDestroyGLVisViewerCtx_Private(void *vctx)
22d71ae5a4SJacob Faibussowitsch {
2346900a7fSStefano Zampini   DMDAFieldGLVisViewerCtx *ctx = (DMDAFieldGLVisViewerCtx *)vctx;
248135c375SStefano Zampini 
258135c375SStefano Zampini   PetscFunctionBegin;
269566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx->xlocal));
279566063dSJacob Faibussowitsch   PetscCall(PetscFree(vctx));
283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
298135c375SStefano Zampini }
308135c375SStefano Zampini 
3146900a7fSStefano Zampini /*
3246900a7fSStefano Zampini    dactx->ll is false -> all but the last proc per dimension claim the ghosted node on the right
3346900a7fSStefano Zampini    dactx->ll is true -> all but the first proc per dimension claim the ghosted node on the left
3446900a7fSStefano Zampini */
35d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDAGetNumElementsGhosted(DM da, PetscInt *nex, PetscInt *ney, PetscInt *nez)
36d71ae5a4SJacob Faibussowitsch {
3746900a7fSStefano Zampini   DMDAGhostedGLVisViewerCtx *dactx;
3846900a7fSStefano Zampini   PetscInt                   sx, sy, sz, ien, jen, ken;
398135c375SStefano Zampini 
408135c375SStefano Zampini   PetscFunctionBegin;
41aab5bcd8SJed Brown   /* Appease -Wmaybe-uninitialized */
42aab5bcd8SJed Brown   if (nex) *nex = -1;
43aab5bcd8SJed Brown   if (ney) *ney = -1;
44aab5bcd8SJed Brown   if (nez) *nez = -1;
459566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &sx, &sy, &sz, &ien, &jen, &ken));
469566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(da, &dactx));
4746900a7fSStefano Zampini   if (dactx->ll) {
4846900a7fSStefano Zampini     PetscInt dim;
4946900a7fSStefano Zampini 
509566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(da, &dim));
5146900a7fSStefano Zampini     if (!sx) ien--;
5246900a7fSStefano Zampini     if (!sy && dim > 1) jen--;
5346900a7fSStefano Zampini     if (!sz && dim > 2) ken--;
5446900a7fSStefano Zampini   } else {
5546900a7fSStefano Zampini     PetscInt M, N, P;
5646900a7fSStefano Zampini 
579566063dSJacob Faibussowitsch     PetscCall(DMDAGetInfo(da, NULL, &M, &N, &P, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
588135c375SStefano Zampini     if (sx + ien == M) ien--;
598135c375SStefano Zampini     if (sy + jen == N) jen--;
608135c375SStefano Zampini     if (sz + ken == P) ken--;
6146900a7fSStefano Zampini   }
628135c375SStefano Zampini   if (nex) *nex = ien;
638135c375SStefano Zampini   if (ney) *ney = jen;
648135c375SStefano Zampini   if (nez) *nez = ken;
653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
668135c375SStefano Zampini }
678135c375SStefano Zampini 
688135c375SStefano Zampini /* inherits number of vertices from DMDAGetNumElementsGhosted */
69d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDAGetNumVerticesGhosted(DM da, PetscInt *nvx, PetscInt *nvy, PetscInt *nvz)
70d71ae5a4SJacob Faibussowitsch {
7121414b21SStefano Zampini   PetscInt ien = 0, jen = 0, ken = 0, dim;
7221414b21SStefano Zampini   PetscInt tote;
738135c375SStefano Zampini 
748135c375SStefano Zampini   PetscFunctionBegin;
759566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(da, &dim));
769566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumElementsGhosted(da, &ien, &jen, &ken));
7721414b21SStefano Zampini   tote = ien * (dim > 1 ? jen : 1) * (dim > 2 ? ken : 1);
7821414b21SStefano Zampini   if (tote) {
798135c375SStefano Zampini     ien = ien + 1;
808135c375SStefano Zampini     jen = dim > 1 ? jen + 1 : 1;
818135c375SStefano Zampini     ken = dim > 2 ? ken + 1 : 1;
8221414b21SStefano Zampini   }
838135c375SStefano Zampini   if (nvx) *nvx = ien;
848135c375SStefano Zampini   if (nvy) *nvy = jen;
858135c375SStefano Zampini   if (nvz) *nvz = ken;
863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
878135c375SStefano Zampini }
888135c375SStefano Zampini 
89d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDASampleGLVisFields_Private(PetscObject oX, PetscInt nf, PetscObject oXf[], void *vctx)
90d71ae5a4SJacob Faibussowitsch {
918135c375SStefano Zampini   DM                         da;
9246900a7fSStefano Zampini   DMDAFieldGLVisViewerCtx   *ctx = (DMDAFieldGLVisViewerCtx *)vctx;
9346900a7fSStefano Zampini   DMDAGhostedGLVisViewerCtx *dactx;
948135c375SStefano Zampini   const PetscScalar         *array;
958135c375SStefano Zampini   PetscScalar              **arrayf;
968135c375SStefano Zampini   PetscInt                   i, f, ii, ien, jen, ken, ie, je, ke, bs, *bss;
978135c375SStefano Zampini   PetscInt                   sx, sy, sz, gsx, gsy, gsz, ist, jst, kst, gm, gn, gp;
988135c375SStefano Zampini 
998135c375SStefano Zampini   PetscFunctionBegin;
1009566063dSJacob Faibussowitsch   PetscCall(VecGetDM(ctx->xlocal, &da));
1017a8be351SBarry Smith   PetscCheck(da, PetscObjectComm(oX), PETSC_ERR_ARG_WRONG, "Vector not generated from a DMDA");
1029566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(da, &dactx));
1039566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(ctx->xlocal, &bs));
1049566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, (Vec)oX, INSERT_VALUES, ctx->xlocal));
1059566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, (Vec)oX, INSERT_VALUES, ctx->xlocal));
1069566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumVerticesGhosted(da, &ien, &jen, &ken));
1079566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gsx, &gsy, &gsz, &gm, &gn, &gp));
1089566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &sx, &sy, &sz, NULL, NULL, NULL));
10946900a7fSStefano Zampini   if (dactx->ll) {
11046900a7fSStefano Zampini     kst = jst = ist = 0;
11146900a7fSStefano Zampini   } else {
1128135c375SStefano Zampini     kst = gsz != sz ? 1 : 0;
1138135c375SStefano Zampini     jst = gsy != sy ? 1 : 0;
1148135c375SStefano Zampini     ist = gsx != sx ? 1 : 0;
11546900a7fSStefano Zampini   }
1169566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(nf, &arrayf, nf, &bss));
1179566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(ctx->xlocal, &array));
1188135c375SStefano Zampini   for (f = 0; f < nf; f++) {
1199566063dSJacob Faibussowitsch     PetscCall(VecGetBlockSize((Vec)oXf[f], &bss[f]));
1209566063dSJacob Faibussowitsch     PetscCall(VecGetArray((Vec)oXf[f], &arrayf[f]));
1218135c375SStefano Zampini   }
1228135c375SStefano Zampini   for (ke = kst, ii = 0; ke < kst + ken; ke++) {
1238135c375SStefano Zampini     for (je = jst; je < jst + jen; je++) {
1248135c375SStefano Zampini       for (ie = ist; ie < ist + ien; ie++) {
1258135c375SStefano Zampini         PetscInt cf, b;
1268135c375SStefano Zampini         i = ke * gm * gn + je * gm + ie;
1278135c375SStefano Zampini         for (f = 0, cf = 0; f < nf; f++)
1289371c9d4SSatish Balay           for (b = 0; b < bss[f]; b++) arrayf[f][bss[f] * ii + b] = array[i * bs + cf++];
1298135c375SStefano Zampini         ii++;
1308135c375SStefano Zampini       }
1318135c375SStefano Zampini     }
1328135c375SStefano Zampini   }
1339566063dSJacob Faibussowitsch   for (f = 0; f < nf; f++) PetscCall(VecRestoreArray((Vec)oXf[f], &arrayf[f]));
1349566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(ctx->xlocal, &array));
1359566063dSJacob Faibussowitsch   PetscCall(PetscFree2(arrayf, bss));
1363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1378135c375SStefano Zampini }
1388135c375SStefano Zampini 
139d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode DMSetUpGLVisViewer_DMDA(PetscObject oda, PetscViewer viewer)
140d71ae5a4SJacob Faibussowitsch {
14146900a7fSStefano Zampini   DM da = (DM)oda, daview;
1428135c375SStefano Zampini 
1438135c375SStefano Zampini   PetscFunctionBegin;
1449566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery(oda, "GLVisGraphicsDMDAGhosted", (PetscObject *)&daview));
14546900a7fSStefano Zampini   if (!daview) {
14646900a7fSStefano Zampini     DMDAGhostedGLVisViewerCtx *dactx;
14746900a7fSStefano Zampini     DM                         dacoord = NULL;
14846900a7fSStefano Zampini     Vec                        xcoor, xcoorl;
14946900a7fSStefano Zampini     PetscBool                  hashocoord = PETSC_FALSE;
15046900a7fSStefano Zampini     const PetscInt            *lx, *ly, *lz;
15146900a7fSStefano Zampini     PetscInt                   dim, M, N, P, m, n, p, dof, s, i;
15246900a7fSStefano Zampini 
1539566063dSJacob Faibussowitsch     PetscCall(PetscNew(&dactx));
154a7bda6aaSStefano Zampini     dactx->ll = PETSC_TRUE; /* default to match elements layout obtained by DMDAGetElements */
155d0609cedSBarry Smith     PetscOptionsBegin(PetscObjectComm(oda), oda->prefix, "GLVis PetscViewer DMDA Options", "PetscViewer");
1569566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-viewer_glvis_dm_da_ll", "Left-looking subdomain view", NULL, dactx->ll, &dactx->ll, NULL));
157d0609cedSBarry Smith     PetscOptionsEnd();
1588135c375SStefano Zampini     /* Create a properly ghosted DMDA to visualize the mesh and the fields associated with */
1599566063dSJacob Faibussowitsch     PetscCall(DMDAGetInfo(da, &dim, &M, &N, &P, &m, &n, &p, &dof, &s, NULL, NULL, NULL, NULL));
1609566063dSJacob Faibussowitsch     PetscCall(DMDAGetOwnershipRanges(da, &lx, &ly, &lz));
1619566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)da, "_glvis_mesh_coords", (PetscObject *)&xcoor));
16277eacf09SStefano Zampini     if (!xcoor) {
1639566063dSJacob Faibussowitsch       PetscCall(DMGetCoordinates(da, &xcoor));
16477eacf09SStefano Zampini     } else {
1658f07809aSStefano Zampini       hashocoord = PETSC_TRUE;
16677eacf09SStefano Zampini     }
1676aad120cSJose E. Roman     PetscCall(PetscInfo(da, "Creating auxiliary DMDA for managing GLVis graphics\n"));
1688135c375SStefano Zampini     switch (dim) {
1698135c375SStefano Zampini     case 1:
1709566063dSJacob Faibussowitsch       PetscCall(DMDACreate1d(PetscObjectComm((PetscObject)da), DM_BOUNDARY_NONE, M, dof, 1, lx, &daview));
17148a46eb9SPierre Jolivet       if (!hashocoord) PetscCall(DMDACreate1d(PetscObjectComm((PetscObject)da), DM_BOUNDARY_NONE, M, 1, 1, lx, &dacoord));
1728135c375SStefano Zampini       break;
1738135c375SStefano Zampini     case 2:
1749566063dSJacob Faibussowitsch       PetscCall(DMDACreate2d(PetscObjectComm((PetscObject)da), DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, M, N, m, n, dof, 1, lx, ly, &daview));
17548a46eb9SPierre Jolivet       if (!hashocoord) PetscCall(DMDACreate2d(PetscObjectComm((PetscObject)da), DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, M, N, m, n, 2, 1, lx, ly, &dacoord));
1768135c375SStefano Zampini       break;
1778135c375SStefano Zampini     case 3:
1789566063dSJacob Faibussowitsch       PetscCall(DMDACreate3d(PetscObjectComm((PetscObject)da), DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, M, N, P, m, n, p, dof, 1, lx, ly, lz, &daview));
17948a46eb9SPierre Jolivet       if (!hashocoord) PetscCall(DMDACreate3d(PetscObjectComm((PetscObject)da), DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, M, N, P, m, n, p, 3, 1, lx, ly, lz, &dacoord));
1808135c375SStefano Zampini       break;
181d71ae5a4SJacob Faibussowitsch     default:
182d71ae5a4SJacob Faibussowitsch       SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Unsupported dimension %" PetscInt_FMT, dim);
1838135c375SStefano Zampini     }
1849566063dSJacob Faibussowitsch     PetscCall(DMSetApplicationContext(daview, dactx));
1859566063dSJacob Faibussowitsch     PetscCall(DMSetApplicationContextDestroy(daview, DMDAGhostedDestroyGLVisViewerCtx_Private));
1869566063dSJacob Faibussowitsch     PetscCall(DMSetUp(daview));
1878135c375SStefano Zampini     if (!xcoor) {
1889566063dSJacob Faibussowitsch       PetscCall(DMDASetUniformCoordinates(daview, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
1899566063dSJacob Faibussowitsch       PetscCall(DMGetCoordinates(daview, &xcoor));
1908135c375SStefano Zampini     }
1918f07809aSStefano Zampini     if (dacoord) {
1929566063dSJacob Faibussowitsch       PetscCall(DMSetUp(dacoord));
1939566063dSJacob Faibussowitsch       PetscCall(DMCreateLocalVector(dacoord, &xcoorl));
1949566063dSJacob Faibussowitsch       PetscCall(DMGlobalToLocalBegin(dacoord, xcoor, INSERT_VALUES, xcoorl));
1959566063dSJacob Faibussowitsch       PetscCall(DMGlobalToLocalEnd(dacoord, xcoor, INSERT_VALUES, xcoorl));
1969566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&dacoord));
1978f07809aSStefano Zampini     } else {
1988f07809aSStefano Zampini       PetscInt    ien, jen, ken, nc, nl, cdof, deg;
1993e7cab22SLisandro Dalcin       char        fecmesh[64];
20021414b21SStefano Zampini       const char *name;
20121414b21SStefano Zampini       PetscBool   flg;
2028135c375SStefano Zampini 
2039566063dSJacob Faibussowitsch       PetscCall(DMDAGetNumElementsGhosted(daview, &ien, &jen, &ken));
2048f07809aSStefano Zampini       nc = ien * (jen > 0 ? jen : 1) * (ken > 0 ? ken : 1);
2058f07809aSStefano Zampini 
2069566063dSJacob Faibussowitsch       PetscCall(VecGetLocalSize(xcoor, &nl));
2071dca8a05SBarry Smith       PetscCheck(!nc || (nl % nc) == 0, PETSC_COMM_SELF, PETSC_ERR_SUP, "Incompatible local coordinate size %" PetscInt_FMT " and number of cells %" PetscInt_FMT, nl, nc);
2089566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(xcoor, &xcoorl));
2099566063dSJacob Faibussowitsch       PetscCall(VecCopy(xcoor, xcoorl));
2109566063dSJacob Faibussowitsch       PetscCall(VecSetDM(xcoorl, NULL));
2119566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetName((PetscObject)xcoor, &name));
2129566063dSJacob Faibussowitsch       PetscCall(PetscStrbeginswith(name, "FiniteElementCollection:", &flg));
21321414b21SStefano Zampini       if (!flg) {
21421414b21SStefano Zampini         deg = 0;
21521414b21SStefano Zampini         if (nc && nl) {
2168f07809aSStefano Zampini           cdof = nl / (nc * dim);
2178f07809aSStefano Zampini           deg  = 1;
2188f07809aSStefano Zampini           while (1) {
2198f07809aSStefano Zampini             PetscInt degd = 1;
2208f07809aSStefano Zampini             for (i = 0; i < dim; i++) degd *= (deg + 1);
22163a3b9bcSJacob Faibussowitsch             PetscCheck(degd <= cdof, PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell dofs %" PetscInt_FMT, cdof);
2228f07809aSStefano Zampini             if (degd == cdof) break;
2238f07809aSStefano Zampini             deg++;
2248f07809aSStefano Zampini           }
22521414b21SStefano Zampini         }
22663a3b9bcSJacob Faibussowitsch         PetscCall(PetscSNPrintf(fecmesh, sizeof(fecmesh), "FiniteElementCollection: L2_T1_%" PetscInt_FMT "D_P%" PetscInt_FMT, dim, deg));
2279566063dSJacob Faibussowitsch         PetscCall(PetscObjectSetName((PetscObject)xcoorl, fecmesh));
2281baa6e33SBarry Smith       } else PetscCall(PetscObjectSetName((PetscObject)xcoorl, name));
2298f07809aSStefano Zampini     }
2308135c375SStefano Zampini 
23146900a7fSStefano Zampini     /* xcoorl is composed with the ghosted DMDA, the ghosted coordinate DMDA (if present) is only available through this vector */
2329566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)daview, "GLVisGraphicsCoordsGhosted", (PetscObject)xcoorl));
2339566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)xcoorl));
2348135c375SStefano Zampini 
23577eacf09SStefano Zampini     /* daview is composed with the original DMDA */
2369566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose(oda, "GLVisGraphicsDMDAGhosted", (PetscObject)daview));
2379566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)daview));
23846900a7fSStefano Zampini   }
23977eacf09SStefano Zampini 
240c9493c35SStefano Zampini   /* customize the viewer if present */
241c9493c35SStefano Zampini   if (viewer) {
24246900a7fSStefano Zampini     DMDAFieldGLVisViewerCtx   *ctx;
24346900a7fSStefano Zampini     DMDAGhostedGLVisViewerCtx *dactx;
24446900a7fSStefano Zampini     char                       fec[64];
245ab166b77SStefano Zampini     Vec                        xlocal, *Ufield;
24646900a7fSStefano Zampini     const char               **dafieldname;
24746900a7fSStefano Zampini     char                     **fec_type, **fieldname;
24846900a7fSStefano Zampini     PetscInt                  *nlocal, *bss, *dims;
24946900a7fSStefano Zampini     PetscInt                   dim, M, N, P, dof, s, i, nf;
25046900a7fSStefano Zampini     PetscBool                  bsset;
2518f07809aSStefano Zampini 
2529566063dSJacob Faibussowitsch     PetscCall(DMDAGetInfo(daview, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
2539566063dSJacob Faibussowitsch     PetscCall(DMGetApplicationContext(daview, &dactx));
2549566063dSJacob Faibussowitsch     PetscCall(DMCreateLocalVector(daview, &xlocal));
2559566063dSJacob Faibussowitsch     PetscCall(DMDAGetFieldNames(da, (const char *const **)&dafieldname));
2569566063dSJacob Faibussowitsch     PetscCall(DMDAGetNumVerticesGhosted(daview, &M, &N, &P));
25763a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(fec, sizeof(fec), "FiniteElementCollection: H1_%" PetscInt_FMT "D_P1", dim));
2589566063dSJacob Faibussowitsch     PetscCall(PetscMalloc6(dof, &fec_type, dof, &nlocal, dof, &bss, dof, &dims, dof, &fieldname, dof, &Ufield));
2598135c375SStefano Zampini     for (i = 0; i < dof; i++) bss[i] = 1;
2608135c375SStefano Zampini     nf = dof;
2618135c375SStefano Zampini 
262d0609cedSBarry Smith     PetscOptionsBegin(PetscObjectComm(oda), oda->prefix, "GLVis PetscViewer DMDA Field options", "PetscViewer");
2639566063dSJacob Faibussowitsch     PetscCall(PetscOptionsIntArray("-viewer_glvis_dm_da_bs", "Block sizes for subfields; enable vector representation", NULL, bss, &nf, &bsset));
264d0609cedSBarry Smith     PetscOptionsEnd();
2658135c375SStefano Zampini     if (bsset) {
2668135c375SStefano Zampini       PetscInt t;
2678135c375SStefano Zampini       for (i = 0, t = 0; i < nf; i++) t += bss[i];
26863a3b9bcSJacob Faibussowitsch       PetscCheck(t == dof, PetscObjectComm(oda), PETSC_ERR_USER, "Sum of block sizes %" PetscInt_FMT " should equal %" PetscInt_FMT, t, dof);
2698135c375SStefano Zampini     } else nf = dof;
2708135c375SStefano Zampini 
2718135c375SStefano Zampini     for (i = 0, s = 0; i < nf; i++) {
2729566063dSJacob Faibussowitsch       PetscCall(PetscStrallocpy(fec, &fec_type[i]));
2738135c375SStefano Zampini       if (bss[i] == 1) {
2749566063dSJacob Faibussowitsch         PetscCall(PetscStrallocpy(dafieldname[s], &fieldname[i]));
2758135c375SStefano Zampini       } else {
276*c6a7a370SJeremy L Thompson         const char     prefix[]   = "Vector-";
277*c6a7a370SJeremy L Thompson         const size_t   prefix_len = PETSC_STATIC_ARRAY_LENGTH(prefix);
278*c6a7a370SJeremy L Thompson         const PetscInt bss_i      = bss[i];
279*c6a7a370SJeremy L Thompson         char          *fname_i    = NULL;
280*c6a7a370SJeremy L Thompson         size_t         tlen       = prefix_len, cur_len;
281*c6a7a370SJeremy L Thompson 
282*c6a7a370SJeremy L Thompson         for (PetscInt b = 0; b < bss_i; ++b) {
2838135c375SStefano Zampini           size_t len;
284*c6a7a370SJeremy L Thompson 
2859566063dSJacob Faibussowitsch           PetscCall(PetscStrlen(dafieldname[s + b], &len));
2868135c375SStefano Zampini           tlen += len + 1; /* field + "-" */
2878135c375SStefano Zampini         }
288*c6a7a370SJeremy L Thompson         PetscCall(PetscMalloc1(tlen, &fname_i));
289*c6a7a370SJeremy L Thompson         PetscCall(PetscArraycpy(fname_i, prefix, prefix_len - 1));
290*c6a7a370SJeremy L Thompson         cur_len = prefix_len;
291*c6a7a370SJeremy L Thompson         for (PetscInt b = 0; b < bss_i; ++b) {
292*c6a7a370SJeremy L Thompson           const char *dafname = dafieldname[s + b];
293*c6a7a370SJeremy L Thompson           size_t      len;
294*c6a7a370SJeremy L Thompson 
295*c6a7a370SJeremy L Thompson           PetscCall(PetscStrlen(dafname, &len));
296*c6a7a370SJeremy L Thompson           PetscCall(PetscStrncpy(fname_i + cur_len, dafname, len + 1));
297*c6a7a370SJeremy L Thompson           cur_len += len + 1;
298*c6a7a370SJeremy L Thompson           // don't add for final iteration of the loop
299*c6a7a370SJeremy L Thompson           if ((b + 1) < bss_i) fname_i[cur_len++] = '-';
3008135c375SStefano Zampini         }
301*c6a7a370SJeremy L Thompson         fname_i[cur_len] = '\0';
302*c6a7a370SJeremy L Thompson         fieldname[i]     = fname_i;
3038135c375SStefano Zampini       }
304bb77a09fSStefano Zampini       dims[i]   = dim;
3058135c375SStefano Zampini       nlocal[i] = M * N * P * bss[i];
3068135c375SStefano Zampini       s += bss[i];
3078135c375SStefano Zampini     }
3088135c375SStefano Zampini 
30946900a7fSStefano Zampini     /* the viewer context takes ownership of xlocal and destroys it in DMDAFieldDestroyGLVisViewerCtx_Private */
3109566063dSJacob Faibussowitsch     PetscCall(PetscNew(&ctx));
3118135c375SStefano Zampini     ctx->xlocal = xlocal;
3128135c375SStefano Zampini 
3134cac2994SStefano Zampini     /* create work vectors */
3144cac2994SStefano Zampini     for (i = 0; i < nf; i++) {
3159566063dSJacob Faibussowitsch       PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)da), nlocal[i], PETSC_DECIDE, &Ufield[i]));
3169566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)Ufield[i], fieldname[i]));
3179566063dSJacob Faibussowitsch       PetscCall(VecSetBlockSize(Ufield[i], bss[i]));
3189566063dSJacob Faibussowitsch       PetscCall(VecSetDM(Ufield[i], da));
3194cac2994SStefano Zampini     }
3204cac2994SStefano Zampini 
3219566063dSJacob Faibussowitsch     PetscCall(PetscViewerGLVisSetFields(viewer, nf, (const char **)fec_type, dims, DMDASampleGLVisFields_Private, (PetscObject *)Ufield, ctx, DMDAFieldDestroyGLVisViewerCtx_Private));
3228135c375SStefano Zampini     for (i = 0; i < nf; i++) {
3239566063dSJacob Faibussowitsch       PetscCall(PetscFree(fec_type[i]));
3249566063dSJacob Faibussowitsch       PetscCall(PetscFree(fieldname[i]));
3259566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&Ufield[i]));
3268135c375SStefano Zampini     }
3279566063dSJacob Faibussowitsch     PetscCall(PetscFree6(fec_type, nlocal, bss, dims, fieldname, Ufield));
328c9493c35SStefano Zampini   }
3293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3308135c375SStefano Zampini }
3318135c375SStefano Zampini 
332d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDAView_GLVis_ASCII(DM dm, PetscViewer viewer)
333d71ae5a4SJacob Faibussowitsch {
33477eacf09SStefano Zampini   DM                 da, cda;
3358135c375SStefano Zampini   Vec                xcoorl;
3363924b612SStefano Zampini   PetscMPIInt        size;
3378135c375SStefano Zampini   const PetscScalar *array;
3388135c375SStefano Zampini   PetscContainer     glvis_container;
33977eacf09SStefano Zampini   PetscInt           dim, sdim, i, vid[8], mid, cid, cdof;
34021414b21SStefano Zampini   PetscInt           sx, sy, sz, ie, je, ke, ien, jen, ken, nel;
3418135c375SStefano Zampini   PetscInt           gsx, gsy, gsz, gm, gn, gp, kst, jst, ist;
3428135c375SStefano Zampini   PetscBool          enabled = PETSC_TRUE, isascii;
3438f07809aSStefano Zampini   const char        *fmt;
3448135c375SStefano Zampini 
3458135c375SStefano Zampini   PetscFunctionBegin;
3468135c375SStefano Zampini   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3478135c375SStefano Zampini   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
3489566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
3497a8be351SBarry Smith   PetscCheck(isascii, PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Viewer must be of type VIEWERASCII");
3509566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)viewer), &size));
35108401ef6SPierre Jolivet   PetscCheck(size <= 1, PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Use single sequential viewers for parallel visualization");
3529566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
3538135c375SStefano Zampini 
3548135c375SStefano Zampini   /* get container: determines if a process visualizes is portion of the data or not */
3559566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)viewer, "_glvis_info_container", (PetscObject *)&glvis_container));
3567a8be351SBarry Smith   PetscCheck(glvis_container, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Missing GLVis container");
3578135c375SStefano Zampini   {
3588135c375SStefano Zampini     PetscViewerGLVisInfo glvis_info;
3599566063dSJacob Faibussowitsch     PetscCall(PetscContainerGetPointer(glvis_container, (void **)&glvis_info));
3608135c375SStefano Zampini     enabled = glvis_info->enabled;
36177eacf09SStefano Zampini     fmt     = glvis_info->fmt;
3628135c375SStefano Zampini   }
363c9493c35SStefano Zampini   /* this can happen if we are calling DMView outside of VecView_GLVis */
3649566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)dm, "GLVisGraphicsDMDAGhosted", (PetscObject *)&da));
3659566063dSJacob Faibussowitsch   if (!da) PetscCall(DMSetUpGLVisViewer_DMDA((PetscObject)dm, NULL));
3669566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)dm, "GLVisGraphicsDMDAGhosted", (PetscObject *)&da));
3677a8be351SBarry Smith   PetscCheck(da, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Missing GLVis ghosted DMDA");
3689566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(da, &sdim));
3698135c375SStefano Zampini 
3709566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "MFEM mesh v1.0\n"));
3719566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "\ndimension\n"));
37263a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT "\n", dim));
3738135c375SStefano Zampini 
3748135c375SStefano Zampini   if (!enabled) {
3759566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "\nelements\n"));
37663a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "0\n"));
3779566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "\nboundary\n"));
37863a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "0\n"));
3799566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "\nvertices\n"));
38063a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "0\n"));
38163a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT "\n", sdim));
3823ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
3838135c375SStefano Zampini   }
3848135c375SStefano Zampini 
3859566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumElementsGhosted(da, &ien, &jen, &ken));
38621414b21SStefano Zampini   nel = ien;
38721414b21SStefano Zampini   if (dim > 1) nel *= jen;
38821414b21SStefano Zampini   if (dim > 2) nel *= ken;
3899566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "\nelements\n"));
39063a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT "\n", nel));
3918135c375SStefano Zampini   switch (dim) {
3928135c375SStefano Zampini   case 1:
3938135c375SStefano Zampini     for (ie = 0; ie < ien; ie++) {
3948135c375SStefano Zampini       vid[0] = ie;
3958135c375SStefano Zampini       vid[1] = ie + 1;
3968135c375SStefano Zampini       mid    = 1; /* material id */
3978135c375SStefano Zampini       cid    = 1; /* segment */
39863a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n", mid, cid, vid[0], vid[1]));
3998135c375SStefano Zampini     }
4008135c375SStefano Zampini     break;
4018135c375SStefano Zampini   case 2:
4028135c375SStefano Zampini     for (je = 0; je < jen; je++) {
4038135c375SStefano Zampini       for (ie = 0; ie < ien; ie++) {
4048135c375SStefano Zampini         vid[0] = je * (ien + 1) + ie;
4058135c375SStefano Zampini         vid[1] = je * (ien + 1) + ie + 1;
4068135c375SStefano Zampini         vid[2] = (je + 1) * (ien + 1) + ie + 1;
4078135c375SStefano Zampini         vid[3] = (je + 1) * (ien + 1) + ie;
4088135c375SStefano Zampini         mid    = 1; /* material id */
4098135c375SStefano Zampini         cid    = 3; /* quad */
41063a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n", mid, cid, vid[0], vid[1], vid[2], vid[3]));
4118135c375SStefano Zampini       }
4128135c375SStefano Zampini     }
4138135c375SStefano Zampini     break;
4148135c375SStefano Zampini   case 3:
4158135c375SStefano Zampini     for (ke = 0; ke < ken; ke++) {
4168135c375SStefano Zampini       for (je = 0; je < jen; je++) {
4178135c375SStefano Zampini         for (ie = 0; ie < ien; ie++) {
4188135c375SStefano Zampini           vid[0] = ke * (jen + 1) * (ien + 1) + je * (ien + 1) + ie;
4198135c375SStefano Zampini           vid[1] = ke * (jen + 1) * (ien + 1) + je * (ien + 1) + ie + 1;
4208135c375SStefano Zampini           vid[2] = ke * (jen + 1) * (ien + 1) + (je + 1) * (ien + 1) + ie + 1;
4218135c375SStefano Zampini           vid[3] = ke * (jen + 1) * (ien + 1) + (je + 1) * (ien + 1) + ie;
4228135c375SStefano Zampini           vid[4] = (ke + 1) * (jen + 1) * (ien + 1) + je * (ien + 1) + ie;
4238135c375SStefano Zampini           vid[5] = (ke + 1) * (jen + 1) * (ien + 1) + je * (ien + 1) + ie + 1;
4248135c375SStefano Zampini           vid[6] = (ke + 1) * (jen + 1) * (ien + 1) + (je + 1) * (ien + 1) + ie + 1;
4258135c375SStefano Zampini           vid[7] = (ke + 1) * (jen + 1) * (ien + 1) + (je + 1) * (ien + 1) + ie;
4268135c375SStefano Zampini           mid    = 1; /* material id */
4278135c375SStefano Zampini           cid    = 5; /* hex */
42863a3b9bcSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n", mid, cid, vid[0], vid[1], vid[2], vid[3], vid[4], vid[5], vid[6], vid[7]));
4298135c375SStefano Zampini         }
4308135c375SStefano Zampini       }
4318135c375SStefano Zampini     }
4328135c375SStefano Zampini     break;
433d71ae5a4SJacob Faibussowitsch   default:
434d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Unsupported dimension %" PetscInt_FMT, dim);
4358135c375SStefano Zampini   }
4369566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "\nboundary\n"));
43763a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "0\n"));
4388135c375SStefano Zampini 
43977eacf09SStefano Zampini   /* vertex coordinates */
4409566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)da, "GLVisGraphicsCoordsGhosted", (PetscObject *)&xcoorl));
4417a8be351SBarry Smith   PetscCheck(xcoorl, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Missing GLVis ghosted coords");
4429566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumVerticesGhosted(da, &ien, &jen, &ken));
4439566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "\nvertices\n"));
44463a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT "\n", ien * jen * ken));
44521414b21SStefano Zampini   if (nel) {
4469566063dSJacob Faibussowitsch     PetscCall(VecGetDM(xcoorl, &cda));
4479566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(xcoorl, &array));
4488f07809aSStefano Zampini     if (!cda) { /* HO viz */
4498f07809aSStefano Zampini       const char *fecname;
4508f07809aSStefano Zampini       PetscInt    nc, nl;
4518f07809aSStefano Zampini 
4529566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetName((PetscObject)xcoorl, &fecname));
4539566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "nodes\n"));
4549566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "FiniteElementSpace\n"));
4559566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "%s\n", fecname));
45663a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "VDim: %" PetscInt_FMT "\n", sdim));
4579566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Ordering: 1\n\n")); /*Ordering::byVDIM*/
4588f07809aSStefano Zampini       /* L2 coordinates */
4599566063dSJacob Faibussowitsch       PetscCall(DMDAGetNumElementsGhosted(da, &ien, &jen, &ken));
4609566063dSJacob Faibussowitsch       PetscCall(VecGetLocalSize(xcoorl, &nl));
4618f07809aSStefano Zampini       nc   = ien * (jen > 0 ? jen : 1) * (ken > 0 ? ken : 1);
46221414b21SStefano Zampini       cdof = nc ? nl / nc : 0;
4638f07809aSStefano Zampini       if (!ien) ien++;
4648f07809aSStefano Zampini       if (!jen) jen++;
4658f07809aSStefano Zampini       if (!ken) ken++;
4668f07809aSStefano Zampini       ist = jst = kst = 0;
4678f07809aSStefano Zampini       gm              = ien;
4688f07809aSStefano Zampini       gn              = jen;
4698f07809aSStefano Zampini       gp              = ken;
47077eacf09SStefano Zampini     } else {
47146900a7fSStefano Zampini       DMDAGhostedGLVisViewerCtx *dactx;
47246900a7fSStefano Zampini 
4739566063dSJacob Faibussowitsch       PetscCall(DMGetApplicationContext(da, &dactx));
47463a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT "\n", sdim));
4758f07809aSStefano Zampini       cdof = sdim;
4769566063dSJacob Faibussowitsch       PetscCall(DMDAGetCorners(da, &sx, &sy, &sz, NULL, NULL, NULL));
4779566063dSJacob Faibussowitsch       PetscCall(DMDAGetGhostCorners(da, &gsx, &gsy, &gsz, &gm, &gn, &gp));
47846900a7fSStefano Zampini       if (dactx->ll) {
47946900a7fSStefano Zampini         kst = jst = ist = 0;
48046900a7fSStefano Zampini       } else {
4818135c375SStefano Zampini         kst = gsz != sz ? 1 : 0;
4828135c375SStefano Zampini         jst = gsy != sy ? 1 : 0;
4838135c375SStefano Zampini         ist = gsx != sx ? 1 : 0;
4848f07809aSStefano Zampini       }
48546900a7fSStefano Zampini     }
4868135c375SStefano Zampini     for (ke = kst; ke < kst + ken; ke++) {
4878135c375SStefano Zampini       for (je = jst; je < jst + jen; je++) {
4888135c375SStefano Zampini         for (ie = ist; ie < ist + ien; ie++) {
48977eacf09SStefano Zampini           PetscInt c;
4908135c375SStefano Zampini 
4918135c375SStefano Zampini           i = ke * gm * gn + je * gm + ie;
49277eacf09SStefano Zampini           for (c = 0; c < cdof / sdim; c++) {
49377eacf09SStefano Zampini             PetscInt d;
49448a46eb9SPierre Jolivet             for (d = 0; d < sdim; d++) PetscCall(PetscViewerASCIIPrintf(viewer, fmt, PetscRealPart(array[cdof * i + c * sdim + d])));
4959566063dSJacob Faibussowitsch             PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
49677eacf09SStefano Zampini           }
4978135c375SStefano Zampini         }
4988135c375SStefano Zampini       }
4998135c375SStefano Zampini     }
5009566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(xcoorl, &array));
50121414b21SStefano Zampini   }
5023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5038135c375SStefano Zampini }
5048135c375SStefano Zampini 
505d71ae5a4SJacob Faibussowitsch PetscErrorCode DMView_DA_GLVis(DM dm, PetscViewer viewer)
506d71ae5a4SJacob Faibussowitsch {
5078135c375SStefano Zampini   PetscFunctionBegin;
5089566063dSJacob Faibussowitsch   PetscCall(DMView_GLVis(dm, viewer, DMDAView_GLVis_ASCII));
5093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5108135c375SStefano Zampini }
511