xref: /petsc/src/dm/impls/da/grglvis.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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 
109371c9d4SSatish Balay static PetscErrorCode DMDAGhostedDestroyGLVisViewerCtx_Private(void **vctx) {
1146900a7fSStefano Zampini   PetscFunctionBegin;
129566063dSJacob Faibussowitsch   PetscCall(PetscFree(*vctx));
1346900a7fSStefano Zampini   PetscFunctionReturn(0);
1446900a7fSStefano Zampini }
1546900a7fSStefano Zampini 
1646900a7fSStefano Zampini typedef struct {
1746900a7fSStefano Zampini   Vec xlocal;
1846900a7fSStefano Zampini } DMDAFieldGLVisViewerCtx;
1946900a7fSStefano Zampini 
209371c9d4SSatish Balay static PetscErrorCode DMDAFieldDestroyGLVisViewerCtx_Private(void *vctx) {
2146900a7fSStefano Zampini   DMDAFieldGLVisViewerCtx *ctx = (DMDAFieldGLVisViewerCtx *)vctx;
228135c375SStefano Zampini 
238135c375SStefano Zampini   PetscFunctionBegin;
249566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx->xlocal));
259566063dSJacob Faibussowitsch   PetscCall(PetscFree(vctx));
268135c375SStefano Zampini   PetscFunctionReturn(0);
278135c375SStefano Zampini }
288135c375SStefano Zampini 
2946900a7fSStefano Zampini /*
3046900a7fSStefano Zampini    dactx->ll is false -> all but the last proc per dimension claim the ghosted node on the right
3146900a7fSStefano Zampini    dactx->ll is true -> all but the first proc per dimension claim the ghosted node on the left
3246900a7fSStefano Zampini */
339371c9d4SSatish Balay static PetscErrorCode DMDAGetNumElementsGhosted(DM da, PetscInt *nex, PetscInt *ney, PetscInt *nez) {
3446900a7fSStefano Zampini   DMDAGhostedGLVisViewerCtx *dactx;
3546900a7fSStefano Zampini   PetscInt                   sx, sy, sz, ien, jen, ken;
368135c375SStefano Zampini 
378135c375SStefano Zampini   PetscFunctionBegin;
38aab5bcd8SJed Brown   /* Appease -Wmaybe-uninitialized */
39aab5bcd8SJed Brown   if (nex) *nex = -1;
40aab5bcd8SJed Brown   if (ney) *ney = -1;
41aab5bcd8SJed Brown   if (nez) *nez = -1;
429566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &sx, &sy, &sz, &ien, &jen, &ken));
439566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(da, &dactx));
4446900a7fSStefano Zampini   if (dactx->ll) {
4546900a7fSStefano Zampini     PetscInt dim;
4646900a7fSStefano Zampini 
479566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(da, &dim));
4846900a7fSStefano Zampini     if (!sx) ien--;
4946900a7fSStefano Zampini     if (!sy && dim > 1) jen--;
5046900a7fSStefano Zampini     if (!sz && dim > 2) ken--;
5146900a7fSStefano Zampini   } else {
5246900a7fSStefano Zampini     PetscInt M, N, P;
5346900a7fSStefano Zampini 
549566063dSJacob Faibussowitsch     PetscCall(DMDAGetInfo(da, NULL, &M, &N, &P, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
558135c375SStefano Zampini     if (sx + ien == M) ien--;
568135c375SStefano Zampini     if (sy + jen == N) jen--;
578135c375SStefano Zampini     if (sz + ken == P) ken--;
5846900a7fSStefano Zampini   }
598135c375SStefano Zampini   if (nex) *nex = ien;
608135c375SStefano Zampini   if (ney) *ney = jen;
618135c375SStefano Zampini   if (nez) *nez = ken;
628135c375SStefano Zampini   PetscFunctionReturn(0);
638135c375SStefano Zampini }
648135c375SStefano Zampini 
658135c375SStefano Zampini /* inherits number of vertices from DMDAGetNumElementsGhosted */
669371c9d4SSatish Balay static PetscErrorCode DMDAGetNumVerticesGhosted(DM da, PetscInt *nvx, PetscInt *nvy, PetscInt *nvz) {
6721414b21SStefano Zampini   PetscInt ien = 0, jen = 0, ken = 0, dim;
6821414b21SStefano Zampini   PetscInt tote;
698135c375SStefano Zampini 
708135c375SStefano Zampini   PetscFunctionBegin;
719566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(da, &dim));
729566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumElementsGhosted(da, &ien, &jen, &ken));
7321414b21SStefano Zampini   tote = ien * (dim > 1 ? jen : 1) * (dim > 2 ? ken : 1);
7421414b21SStefano Zampini   if (tote) {
758135c375SStefano Zampini     ien = ien + 1;
768135c375SStefano Zampini     jen = dim > 1 ? jen + 1 : 1;
778135c375SStefano Zampini     ken = dim > 2 ? ken + 1 : 1;
7821414b21SStefano Zampini   }
798135c375SStefano Zampini   if (nvx) *nvx = ien;
808135c375SStefano Zampini   if (nvy) *nvy = jen;
818135c375SStefano Zampini   if (nvz) *nvz = ken;
828135c375SStefano Zampini   PetscFunctionReturn(0);
838135c375SStefano Zampini }
848135c375SStefano Zampini 
859371c9d4SSatish Balay static PetscErrorCode DMDASampleGLVisFields_Private(PetscObject oX, PetscInt nf, PetscObject oXf[], void *vctx) {
868135c375SStefano Zampini   DM                         da;
8746900a7fSStefano Zampini   DMDAFieldGLVisViewerCtx   *ctx = (DMDAFieldGLVisViewerCtx *)vctx;
8846900a7fSStefano Zampini   DMDAGhostedGLVisViewerCtx *dactx;
898135c375SStefano Zampini   const PetscScalar         *array;
908135c375SStefano Zampini   PetscScalar              **arrayf;
918135c375SStefano Zampini   PetscInt                   i, f, ii, ien, jen, ken, ie, je, ke, bs, *bss;
928135c375SStefano Zampini   PetscInt                   sx, sy, sz, gsx, gsy, gsz, ist, jst, kst, gm, gn, gp;
938135c375SStefano Zampini 
948135c375SStefano Zampini   PetscFunctionBegin;
959566063dSJacob Faibussowitsch   PetscCall(VecGetDM(ctx->xlocal, &da));
967a8be351SBarry Smith   PetscCheck(da, PetscObjectComm(oX), PETSC_ERR_ARG_WRONG, "Vector not generated from a DMDA");
979566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(da, &dactx));
989566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(ctx->xlocal, &bs));
999566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, (Vec)oX, INSERT_VALUES, ctx->xlocal));
1009566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, (Vec)oX, INSERT_VALUES, ctx->xlocal));
1019566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumVerticesGhosted(da, &ien, &jen, &ken));
1029566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gsx, &gsy, &gsz, &gm, &gn, &gp));
1039566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &sx, &sy, &sz, NULL, NULL, NULL));
10446900a7fSStefano Zampini   if (dactx->ll) {
10546900a7fSStefano Zampini     kst = jst = ist = 0;
10646900a7fSStefano Zampini   } else {
1078135c375SStefano Zampini     kst = gsz != sz ? 1 : 0;
1088135c375SStefano Zampini     jst = gsy != sy ? 1 : 0;
1098135c375SStefano Zampini     ist = gsx != sx ? 1 : 0;
11046900a7fSStefano Zampini   }
1119566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(nf, &arrayf, nf, &bss));
1129566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(ctx->xlocal, &array));
1138135c375SStefano Zampini   for (f = 0; f < nf; f++) {
1149566063dSJacob Faibussowitsch     PetscCall(VecGetBlockSize((Vec)oXf[f], &bss[f]));
1159566063dSJacob Faibussowitsch     PetscCall(VecGetArray((Vec)oXf[f], &arrayf[f]));
1168135c375SStefano Zampini   }
1178135c375SStefano Zampini   for (ke = kst, ii = 0; ke < kst + ken; ke++) {
1188135c375SStefano Zampini     for (je = jst; je < jst + jen; je++) {
1198135c375SStefano Zampini       for (ie = ist; ie < ist + ien; ie++) {
1208135c375SStefano Zampini         PetscInt cf, b;
1218135c375SStefano Zampini         i = ke * gm * gn + je * gm + ie;
1228135c375SStefano Zampini         for (f = 0, cf = 0; f < nf; f++)
1239371c9d4SSatish Balay           for (b = 0; b < bss[f]; b++) arrayf[f][bss[f] * ii + b] = array[i * bs + cf++];
1248135c375SStefano Zampini         ii++;
1258135c375SStefano Zampini       }
1268135c375SStefano Zampini     }
1278135c375SStefano Zampini   }
1289566063dSJacob Faibussowitsch   for (f = 0; f < nf; f++) PetscCall(VecRestoreArray((Vec)oXf[f], &arrayf[f]));
1299566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(ctx->xlocal, &array));
1309566063dSJacob Faibussowitsch   PetscCall(PetscFree2(arrayf, bss));
1318135c375SStefano Zampini   PetscFunctionReturn(0);
1328135c375SStefano Zampini }
1338135c375SStefano Zampini 
1349371c9d4SSatish Balay PETSC_INTERN PetscErrorCode DMSetUpGLVisViewer_DMDA(PetscObject oda, PetscViewer viewer) {
13546900a7fSStefano Zampini   DM da = (DM)oda, daview;
1368135c375SStefano Zampini 
1378135c375SStefano Zampini   PetscFunctionBegin;
1389566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery(oda, "GLVisGraphicsDMDAGhosted", (PetscObject *)&daview));
13946900a7fSStefano Zampini   if (!daview) {
14046900a7fSStefano Zampini     DMDAGhostedGLVisViewerCtx *dactx;
14146900a7fSStefano Zampini     DM                         dacoord = NULL;
14246900a7fSStefano Zampini     Vec                        xcoor, xcoorl;
14346900a7fSStefano Zampini     PetscBool                  hashocoord = PETSC_FALSE;
14446900a7fSStefano Zampini     const PetscInt            *lx, *ly, *lz;
14546900a7fSStefano Zampini     PetscInt                   dim, M, N, P, m, n, p, dof, s, i;
14646900a7fSStefano Zampini 
1479566063dSJacob Faibussowitsch     PetscCall(PetscNew(&dactx));
148a7bda6aaSStefano Zampini     dactx->ll = PETSC_TRUE; /* default to match elements layout obtained by DMDAGetElements */
149d0609cedSBarry Smith     PetscOptionsBegin(PetscObjectComm(oda), oda->prefix, "GLVis PetscViewer DMDA Options", "PetscViewer");
1509566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-viewer_glvis_dm_da_ll", "Left-looking subdomain view", NULL, dactx->ll, &dactx->ll, NULL));
151d0609cedSBarry Smith     PetscOptionsEnd();
1528135c375SStefano Zampini     /* Create a properly ghosted DMDA to visualize the mesh and the fields associated with */
1539566063dSJacob Faibussowitsch     PetscCall(DMDAGetInfo(da, &dim, &M, &N, &P, &m, &n, &p, &dof, &s, NULL, NULL, NULL, NULL));
1549566063dSJacob Faibussowitsch     PetscCall(DMDAGetOwnershipRanges(da, &lx, &ly, &lz));
1559566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)da, "_glvis_mesh_coords", (PetscObject *)&xcoor));
15677eacf09SStefano Zampini     if (!xcoor) {
1579566063dSJacob Faibussowitsch       PetscCall(DMGetCoordinates(da, &xcoor));
15877eacf09SStefano Zampini     } else {
1598f07809aSStefano Zampini       hashocoord = PETSC_TRUE;
16077eacf09SStefano Zampini     }
1616aad120cSJose E. Roman     PetscCall(PetscInfo(da, "Creating auxiliary DMDA for managing GLVis graphics\n"));
1628135c375SStefano Zampini     switch (dim) {
1638135c375SStefano Zampini     case 1:
1649566063dSJacob Faibussowitsch       PetscCall(DMDACreate1d(PetscObjectComm((PetscObject)da), DM_BOUNDARY_NONE, M, dof, 1, lx, &daview));
165*48a46eb9SPierre Jolivet       if (!hashocoord) PetscCall(DMDACreate1d(PetscObjectComm((PetscObject)da), DM_BOUNDARY_NONE, M, 1, 1, lx, &dacoord));
1668135c375SStefano Zampini       break;
1678135c375SStefano Zampini     case 2:
1689566063dSJacob Faibussowitsch       PetscCall(DMDACreate2d(PetscObjectComm((PetscObject)da), DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, M, N, m, n, dof, 1, lx, ly, &daview));
169*48a46eb9SPierre 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));
1708135c375SStefano Zampini       break;
1718135c375SStefano Zampini     case 3:
1729566063dSJacob 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));
173*48a46eb9SPierre 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));
1748135c375SStefano Zampini       break;
1759371c9d4SSatish Balay     default: SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Unsupported dimension %" PetscInt_FMT, dim);
1768135c375SStefano Zampini     }
1779566063dSJacob Faibussowitsch     PetscCall(DMSetApplicationContext(daview, dactx));
1789566063dSJacob Faibussowitsch     PetscCall(DMSetApplicationContextDestroy(daview, DMDAGhostedDestroyGLVisViewerCtx_Private));
1799566063dSJacob Faibussowitsch     PetscCall(DMSetUp(daview));
1808135c375SStefano Zampini     if (!xcoor) {
1819566063dSJacob Faibussowitsch       PetscCall(DMDASetUniformCoordinates(daview, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
1829566063dSJacob Faibussowitsch       PetscCall(DMGetCoordinates(daview, &xcoor));
1838135c375SStefano Zampini     }
1848f07809aSStefano Zampini     if (dacoord) {
1859566063dSJacob Faibussowitsch       PetscCall(DMSetUp(dacoord));
1869566063dSJacob Faibussowitsch       PetscCall(DMCreateLocalVector(dacoord, &xcoorl));
1879566063dSJacob Faibussowitsch       PetscCall(DMGlobalToLocalBegin(dacoord, xcoor, INSERT_VALUES, xcoorl));
1889566063dSJacob Faibussowitsch       PetscCall(DMGlobalToLocalEnd(dacoord, xcoor, INSERT_VALUES, xcoorl));
1899566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&dacoord));
1908f07809aSStefano Zampini     } else {
1918f07809aSStefano Zampini       PetscInt    ien, jen, ken, nc, nl, cdof, deg;
1923e7cab22SLisandro Dalcin       char        fecmesh[64];
19321414b21SStefano Zampini       const char *name;
19421414b21SStefano Zampini       PetscBool   flg;
1958135c375SStefano Zampini 
1969566063dSJacob Faibussowitsch       PetscCall(DMDAGetNumElementsGhosted(daview, &ien, &jen, &ken));
1978f07809aSStefano Zampini       nc = ien * (jen > 0 ? jen : 1) * (ken > 0 ? ken : 1);
1988f07809aSStefano Zampini 
1999566063dSJacob Faibussowitsch       PetscCall(VecGetLocalSize(xcoor, &nl));
2001dca8a05SBarry 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);
2019566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(xcoor, &xcoorl));
2029566063dSJacob Faibussowitsch       PetscCall(VecCopy(xcoor, xcoorl));
2039566063dSJacob Faibussowitsch       PetscCall(VecSetDM(xcoorl, NULL));
2049566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetName((PetscObject)xcoor, &name));
2059566063dSJacob Faibussowitsch       PetscCall(PetscStrbeginswith(name, "FiniteElementCollection:", &flg));
20621414b21SStefano Zampini       if (!flg) {
20721414b21SStefano Zampini         deg = 0;
20821414b21SStefano Zampini         if (nc && nl) {
2098f07809aSStefano Zampini           cdof = nl / (nc * dim);
2108f07809aSStefano Zampini           deg  = 1;
2118f07809aSStefano Zampini           while (1) {
2128f07809aSStefano Zampini             PetscInt degd = 1;
2138f07809aSStefano Zampini             for (i = 0; i < dim; i++) degd *= (deg + 1);
21463a3b9bcSJacob Faibussowitsch             PetscCheck(degd <= cdof, PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell dofs %" PetscInt_FMT, cdof);
2158f07809aSStefano Zampini             if (degd == cdof) break;
2168f07809aSStefano Zampini             deg++;
2178f07809aSStefano Zampini           }
21821414b21SStefano Zampini         }
21963a3b9bcSJacob Faibussowitsch         PetscCall(PetscSNPrintf(fecmesh, sizeof(fecmesh), "FiniteElementCollection: L2_T1_%" PetscInt_FMT "D_P%" PetscInt_FMT, dim, deg));
2209566063dSJacob Faibussowitsch         PetscCall(PetscObjectSetName((PetscObject)xcoorl, fecmesh));
2211baa6e33SBarry Smith       } else PetscCall(PetscObjectSetName((PetscObject)xcoorl, name));
2228f07809aSStefano Zampini     }
2238135c375SStefano Zampini 
22446900a7fSStefano Zampini     /* xcoorl is composed with the ghosted DMDA, the ghosted coordinate DMDA (if present) is only available through this vector */
2259566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)daview, "GLVisGraphicsCoordsGhosted", (PetscObject)xcoorl));
2269566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)xcoorl));
2278135c375SStefano Zampini 
22877eacf09SStefano Zampini     /* daview is composed with the original DMDA */
2299566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose(oda, "GLVisGraphicsDMDAGhosted", (PetscObject)daview));
2309566063dSJacob Faibussowitsch     PetscCall(PetscObjectDereference((PetscObject)daview));
23146900a7fSStefano Zampini   }
23277eacf09SStefano Zampini 
233c9493c35SStefano Zampini   /* customize the viewer if present */
234c9493c35SStefano Zampini   if (viewer) {
23546900a7fSStefano Zampini     DMDAFieldGLVisViewerCtx   *ctx;
23646900a7fSStefano Zampini     DMDAGhostedGLVisViewerCtx *dactx;
23746900a7fSStefano Zampini     char                       fec[64];
238ab166b77SStefano Zampini     Vec                        xlocal, *Ufield;
23946900a7fSStefano Zampini     const char               **dafieldname;
24046900a7fSStefano Zampini     char                     **fec_type, **fieldname;
24146900a7fSStefano Zampini     PetscInt                  *nlocal, *bss, *dims;
24246900a7fSStefano Zampini     PetscInt                   dim, M, N, P, dof, s, i, nf;
24346900a7fSStefano Zampini     PetscBool                  bsset;
2448f07809aSStefano Zampini 
2459566063dSJacob Faibussowitsch     PetscCall(DMDAGetInfo(daview, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
2469566063dSJacob Faibussowitsch     PetscCall(DMGetApplicationContext(daview, &dactx));
2479566063dSJacob Faibussowitsch     PetscCall(DMCreateLocalVector(daview, &xlocal));
2489566063dSJacob Faibussowitsch     PetscCall(DMDAGetFieldNames(da, (const char *const **)&dafieldname));
2499566063dSJacob Faibussowitsch     PetscCall(DMDAGetNumVerticesGhosted(daview, &M, &N, &P));
25063a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(fec, sizeof(fec), "FiniteElementCollection: H1_%" PetscInt_FMT "D_P1", dim));
2519566063dSJacob Faibussowitsch     PetscCall(PetscMalloc6(dof, &fec_type, dof, &nlocal, dof, &bss, dof, &dims, dof, &fieldname, dof, &Ufield));
2528135c375SStefano Zampini     for (i = 0; i < dof; i++) bss[i] = 1;
2538135c375SStefano Zampini     nf = dof;
2548135c375SStefano Zampini 
255d0609cedSBarry Smith     PetscOptionsBegin(PetscObjectComm(oda), oda->prefix, "GLVis PetscViewer DMDA Field options", "PetscViewer");
2569566063dSJacob Faibussowitsch     PetscCall(PetscOptionsIntArray("-viewer_glvis_dm_da_bs", "Block sizes for subfields; enable vector representation", NULL, bss, &nf, &bsset));
257d0609cedSBarry Smith     PetscOptionsEnd();
2588135c375SStefano Zampini     if (bsset) {
2598135c375SStefano Zampini       PetscInt t;
2608135c375SStefano Zampini       for (i = 0, t = 0; i < nf; i++) t += bss[i];
26163a3b9bcSJacob Faibussowitsch       PetscCheck(t == dof, PetscObjectComm(oda), PETSC_ERR_USER, "Sum of block sizes %" PetscInt_FMT " should equal %" PetscInt_FMT, t, dof);
2628135c375SStefano Zampini     } else nf = dof;
2638135c375SStefano Zampini 
2648135c375SStefano Zampini     for (i = 0, s = 0; i < nf; i++) {
2659566063dSJacob Faibussowitsch       PetscCall(PetscStrallocpy(fec, &fec_type[i]));
2668135c375SStefano Zampini       if (bss[i] == 1) {
2679566063dSJacob Faibussowitsch         PetscCall(PetscStrallocpy(dafieldname[s], &fieldname[i]));
2688135c375SStefano Zampini       } else {
2698135c375SStefano Zampini         PetscInt b;
2708135c375SStefano Zampini         size_t   tlen = 9; /* "Vector-" + end */
2718135c375SStefano Zampini         for (b = 0; b < bss[i]; b++) {
2728135c375SStefano Zampini           size_t len;
2739566063dSJacob Faibussowitsch           PetscCall(PetscStrlen(dafieldname[s + b], &len));
2748135c375SStefano Zampini           tlen += len + 1; /* field + "-" */
2758135c375SStefano Zampini         }
2769566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(tlen, &fieldname[i]));
2779566063dSJacob Faibussowitsch         PetscCall(PetscStrcpy(fieldname[i], "Vector-"));
2788135c375SStefano Zampini         for (b = 0; b < bss[i] - 1; b++) {
2799566063dSJacob Faibussowitsch           PetscCall(PetscStrcat(fieldname[i], dafieldname[s + b]));
2809566063dSJacob Faibussowitsch           PetscCall(PetscStrcat(fieldname[i], "-"));
2818135c375SStefano Zampini         }
2829566063dSJacob Faibussowitsch         PetscCall(PetscStrcat(fieldname[i], dafieldname[s + b]));
2838135c375SStefano Zampini       }
284bb77a09fSStefano Zampini       dims[i]   = dim;
2858135c375SStefano Zampini       nlocal[i] = M * N * P * bss[i];
2868135c375SStefano Zampini       s += bss[i];
2878135c375SStefano Zampini     }
2888135c375SStefano Zampini 
28946900a7fSStefano Zampini     /* the viewer context takes ownership of xlocal and destroys it in DMDAFieldDestroyGLVisViewerCtx_Private */
2909566063dSJacob Faibussowitsch     PetscCall(PetscNew(&ctx));
2918135c375SStefano Zampini     ctx->xlocal = xlocal;
2928135c375SStefano Zampini 
2934cac2994SStefano Zampini     /* create work vectors */
2944cac2994SStefano Zampini     for (i = 0; i < nf; i++) {
2959566063dSJacob Faibussowitsch       PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)da), nlocal[i], PETSC_DECIDE, &Ufield[i]));
2969566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)Ufield[i], fieldname[i]));
2979566063dSJacob Faibussowitsch       PetscCall(VecSetBlockSize(Ufield[i], bss[i]));
2989566063dSJacob Faibussowitsch       PetscCall(VecSetDM(Ufield[i], da));
2994cac2994SStefano Zampini     }
3004cac2994SStefano Zampini 
3019566063dSJacob Faibussowitsch     PetscCall(PetscViewerGLVisSetFields(viewer, nf, (const char **)fec_type, dims, DMDASampleGLVisFields_Private, (PetscObject *)Ufield, ctx, DMDAFieldDestroyGLVisViewerCtx_Private));
3028135c375SStefano Zampini     for (i = 0; i < nf; i++) {
3039566063dSJacob Faibussowitsch       PetscCall(PetscFree(fec_type[i]));
3049566063dSJacob Faibussowitsch       PetscCall(PetscFree(fieldname[i]));
3059566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&Ufield[i]));
3068135c375SStefano Zampini     }
3079566063dSJacob Faibussowitsch     PetscCall(PetscFree6(fec_type, nlocal, bss, dims, fieldname, Ufield));
308c9493c35SStefano Zampini   }
3098135c375SStefano Zampini   PetscFunctionReturn(0);
3108135c375SStefano Zampini }
3118135c375SStefano Zampini 
3129371c9d4SSatish Balay static PetscErrorCode DMDAView_GLVis_ASCII(DM dm, PetscViewer viewer) {
31377eacf09SStefano Zampini   DM                 da, cda;
3148135c375SStefano Zampini   Vec                xcoorl;
3153924b612SStefano Zampini   PetscMPIInt        size;
3168135c375SStefano Zampini   const PetscScalar *array;
3178135c375SStefano Zampini   PetscContainer     glvis_container;
31877eacf09SStefano Zampini   PetscInt           dim, sdim, i, vid[8], mid, cid, cdof;
31921414b21SStefano Zampini   PetscInt           sx, sy, sz, ie, je, ke, ien, jen, ken, nel;
3208135c375SStefano Zampini   PetscInt           gsx, gsy, gsz, gm, gn, gp, kst, jst, ist;
3218135c375SStefano Zampini   PetscBool          enabled = PETSC_TRUE, isascii;
3228f07809aSStefano Zampini   const char        *fmt;
3238135c375SStefano Zampini 
3248135c375SStefano Zampini   PetscFunctionBegin;
3258135c375SStefano Zampini   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3268135c375SStefano Zampini   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
3279566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
3287a8be351SBarry Smith   PetscCheck(isascii, PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Viewer must be of type VIEWERASCII");
3299566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)viewer), &size));
33008401ef6SPierre Jolivet   PetscCheck(size <= 1, PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Use single sequential viewers for parallel visualization");
3319566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
3328135c375SStefano Zampini 
3338135c375SStefano Zampini   /* get container: determines if a process visualizes is portion of the data or not */
3349566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)viewer, "_glvis_info_container", (PetscObject *)&glvis_container));
3357a8be351SBarry Smith   PetscCheck(glvis_container, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Missing GLVis container");
3368135c375SStefano Zampini   {
3378135c375SStefano Zampini     PetscViewerGLVisInfo glvis_info;
3389566063dSJacob Faibussowitsch     PetscCall(PetscContainerGetPointer(glvis_container, (void **)&glvis_info));
3398135c375SStefano Zampini     enabled = glvis_info->enabled;
34077eacf09SStefano Zampini     fmt     = glvis_info->fmt;
3418135c375SStefano Zampini   }
342c9493c35SStefano Zampini   /* this can happen if we are calling DMView outside of VecView_GLVis */
3439566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)dm, "GLVisGraphicsDMDAGhosted", (PetscObject *)&da));
3449566063dSJacob Faibussowitsch   if (!da) PetscCall(DMSetUpGLVisViewer_DMDA((PetscObject)dm, NULL));
3459566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)dm, "GLVisGraphicsDMDAGhosted", (PetscObject *)&da));
3467a8be351SBarry Smith   PetscCheck(da, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Missing GLVis ghosted DMDA");
3479566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(da, &sdim));
3488135c375SStefano Zampini 
3499566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "MFEM mesh v1.0\n"));
3509566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "\ndimension\n"));
35163a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT "\n", dim));
3528135c375SStefano Zampini 
3538135c375SStefano Zampini   if (!enabled) {
3549566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "\nelements\n"));
35563a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "0\n"));
3569566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "\nboundary\n"));
35763a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "0\n"));
3589566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "\nvertices\n"));
35963a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "0\n"));
36063a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT "\n", sdim));
3618135c375SStefano Zampini     PetscFunctionReturn(0);
3628135c375SStefano Zampini   }
3638135c375SStefano Zampini 
3649566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumElementsGhosted(da, &ien, &jen, &ken));
36521414b21SStefano Zampini   nel = ien;
36621414b21SStefano Zampini   if (dim > 1) nel *= jen;
36721414b21SStefano Zampini   if (dim > 2) nel *= ken;
3689566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "\nelements\n"));
36963a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT "\n", nel));
3708135c375SStefano Zampini   switch (dim) {
3718135c375SStefano Zampini   case 1:
3728135c375SStefano Zampini     for (ie = 0; ie < ien; ie++) {
3738135c375SStefano Zampini       vid[0] = ie;
3748135c375SStefano Zampini       vid[1] = ie + 1;
3758135c375SStefano Zampini       mid    = 1; /* material id */
3768135c375SStefano Zampini       cid    = 1; /* segment */
37763a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n", mid, cid, vid[0], vid[1]));
3788135c375SStefano Zampini     }
3798135c375SStefano Zampini     break;
3808135c375SStefano Zampini   case 2:
3818135c375SStefano Zampini     for (je = 0; je < jen; je++) {
3828135c375SStefano Zampini       for (ie = 0; ie < ien; ie++) {
3838135c375SStefano Zampini         vid[0] = je * (ien + 1) + ie;
3848135c375SStefano Zampini         vid[1] = je * (ien + 1) + ie + 1;
3858135c375SStefano Zampini         vid[2] = (je + 1) * (ien + 1) + ie + 1;
3868135c375SStefano Zampini         vid[3] = (je + 1) * (ien + 1) + ie;
3878135c375SStefano Zampini         mid    = 1; /* material id */
3888135c375SStefano Zampini         cid    = 3; /* quad */
38963a3b9bcSJacob 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]));
3908135c375SStefano Zampini       }
3918135c375SStefano Zampini     }
3928135c375SStefano Zampini     break;
3938135c375SStefano Zampini   case 3:
3948135c375SStefano Zampini     for (ke = 0; ke < ken; ke++) {
3958135c375SStefano Zampini       for (je = 0; je < jen; je++) {
3968135c375SStefano Zampini         for (ie = 0; ie < ien; ie++) {
3978135c375SStefano Zampini           vid[0] = ke * (jen + 1) * (ien + 1) + je * (ien + 1) + ie;
3988135c375SStefano Zampini           vid[1] = ke * (jen + 1) * (ien + 1) + je * (ien + 1) + ie + 1;
3998135c375SStefano Zampini           vid[2] = ke * (jen + 1) * (ien + 1) + (je + 1) * (ien + 1) + ie + 1;
4008135c375SStefano Zampini           vid[3] = ke * (jen + 1) * (ien + 1) + (je + 1) * (ien + 1) + ie;
4018135c375SStefano Zampini           vid[4] = (ke + 1) * (jen + 1) * (ien + 1) + je * (ien + 1) + ie;
4028135c375SStefano Zampini           vid[5] = (ke + 1) * (jen + 1) * (ien + 1) + je * (ien + 1) + ie + 1;
4038135c375SStefano Zampini           vid[6] = (ke + 1) * (jen + 1) * (ien + 1) + (je + 1) * (ien + 1) + ie + 1;
4048135c375SStefano Zampini           vid[7] = (ke + 1) * (jen + 1) * (ien + 1) + (je + 1) * (ien + 1) + ie;
4058135c375SStefano Zampini           mid    = 1; /* material id */
4068135c375SStefano Zampini           cid    = 5; /* hex */
40763a3b9bcSJacob 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]));
4088135c375SStefano Zampini         }
4098135c375SStefano Zampini       }
4108135c375SStefano Zampini     }
4118135c375SStefano Zampini     break;
4129371c9d4SSatish Balay   default: SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Unsupported dimension %" PetscInt_FMT, dim);
4138135c375SStefano Zampini   }
4149566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "\nboundary\n"));
41563a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "0\n"));
4168135c375SStefano Zampini 
41777eacf09SStefano Zampini   /* vertex coordinates */
4189566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)da, "GLVisGraphicsCoordsGhosted", (PetscObject *)&xcoorl));
4197a8be351SBarry Smith   PetscCheck(xcoorl, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Missing GLVis ghosted coords");
4209566063dSJacob Faibussowitsch   PetscCall(DMDAGetNumVerticesGhosted(da, &ien, &jen, &ken));
4219566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "\nvertices\n"));
42263a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT "\n", ien * jen * ken));
42321414b21SStefano Zampini   if (nel) {
4249566063dSJacob Faibussowitsch     PetscCall(VecGetDM(xcoorl, &cda));
4259566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(xcoorl, &array));
4268f07809aSStefano Zampini     if (!cda) { /* HO viz */
4278f07809aSStefano Zampini       const char *fecname;
4288f07809aSStefano Zampini       PetscInt    nc, nl;
4298f07809aSStefano Zampini 
4309566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetName((PetscObject)xcoorl, &fecname));
4319566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "nodes\n"));
4329566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "FiniteElementSpace\n"));
4339566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "%s\n", fecname));
43463a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "VDim: %" PetscInt_FMT "\n", sdim));
4359566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Ordering: 1\n\n")); /*Ordering::byVDIM*/
4368f07809aSStefano Zampini       /* L2 coordinates */
4379566063dSJacob Faibussowitsch       PetscCall(DMDAGetNumElementsGhosted(da, &ien, &jen, &ken));
4389566063dSJacob Faibussowitsch       PetscCall(VecGetLocalSize(xcoorl, &nl));
4398f07809aSStefano Zampini       nc   = ien * (jen > 0 ? jen : 1) * (ken > 0 ? ken : 1);
44021414b21SStefano Zampini       cdof = nc ? nl / nc : 0;
4418f07809aSStefano Zampini       if (!ien) ien++;
4428f07809aSStefano Zampini       if (!jen) jen++;
4438f07809aSStefano Zampini       if (!ken) ken++;
4448f07809aSStefano Zampini       ist = jst = kst = 0;
4458f07809aSStefano Zampini       gm              = ien;
4468f07809aSStefano Zampini       gn              = jen;
4478f07809aSStefano Zampini       gp              = ken;
44877eacf09SStefano Zampini     } else {
44946900a7fSStefano Zampini       DMDAGhostedGLVisViewerCtx *dactx;
45046900a7fSStefano Zampini 
4519566063dSJacob Faibussowitsch       PetscCall(DMGetApplicationContext(da, &dactx));
45263a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT "\n", sdim));
4538f07809aSStefano Zampini       cdof = sdim;
4549566063dSJacob Faibussowitsch       PetscCall(DMDAGetCorners(da, &sx, &sy, &sz, NULL, NULL, NULL));
4559566063dSJacob Faibussowitsch       PetscCall(DMDAGetGhostCorners(da, &gsx, &gsy, &gsz, &gm, &gn, &gp));
45646900a7fSStefano Zampini       if (dactx->ll) {
45746900a7fSStefano Zampini         kst = jst = ist = 0;
45846900a7fSStefano Zampini       } else {
4598135c375SStefano Zampini         kst = gsz != sz ? 1 : 0;
4608135c375SStefano Zampini         jst = gsy != sy ? 1 : 0;
4618135c375SStefano Zampini         ist = gsx != sx ? 1 : 0;
4628f07809aSStefano Zampini       }
46346900a7fSStefano Zampini     }
4648135c375SStefano Zampini     for (ke = kst; ke < kst + ken; ke++) {
4658135c375SStefano Zampini       for (je = jst; je < jst + jen; je++) {
4668135c375SStefano Zampini         for (ie = ist; ie < ist + ien; ie++) {
46777eacf09SStefano Zampini           PetscInt c;
4688135c375SStefano Zampini 
4698135c375SStefano Zampini           i = ke * gm * gn + je * gm + ie;
47077eacf09SStefano Zampini           for (c = 0; c < cdof / sdim; c++) {
47177eacf09SStefano Zampini             PetscInt d;
472*48a46eb9SPierre Jolivet             for (d = 0; d < sdim; d++) PetscCall(PetscViewerASCIIPrintf(viewer, fmt, PetscRealPart(array[cdof * i + c * sdim + d])));
4739566063dSJacob Faibussowitsch             PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
47477eacf09SStefano Zampini           }
4758135c375SStefano Zampini         }
4768135c375SStefano Zampini       }
4778135c375SStefano Zampini     }
4789566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(xcoorl, &array));
47921414b21SStefano Zampini   }
4808135c375SStefano Zampini   PetscFunctionReturn(0);
4818135c375SStefano Zampini }
4828135c375SStefano Zampini 
4839371c9d4SSatish Balay PetscErrorCode DMView_DA_GLVis(DM dm, PetscViewer viewer) {
4848135c375SStefano Zampini   PetscFunctionBegin;
4859566063dSJacob Faibussowitsch   PetscCall(DMView_GLVis(dm, viewer, DMDAView_GLVis_ASCII));
4868135c375SStefano Zampini   PetscFunctionReturn(0);
4878135c375SStefano Zampini }
488