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
DMDAGhostedDestroyGLVisViewerCtx_Private(PetscCtxRt vctx)10*2a8381b2SBarry Smith static PetscErrorCode DMDAGhostedDestroyGLVisViewerCtx_Private(PetscCtxRt vctx)
11d71ae5a4SJacob Faibussowitsch {
1246900a7fSStefano Zampini PetscFunctionBegin;
13*2a8381b2SBarry Smith PetscCall(PetscFree(*(void **)vctx));
143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1546900a7fSStefano Zampini }
1646900a7fSStefano Zampini
1746900a7fSStefano Zampini typedef struct {
1846900a7fSStefano Zampini Vec xlocal;
1946900a7fSStefano Zampini } DMDAFieldGLVisViewerCtx;
2046900a7fSStefano Zampini
DMDAFieldDestroyGLVisViewerCtx_Private(PetscCtxRt vctx)21*2a8381b2SBarry Smith static PetscErrorCode DMDAFieldDestroyGLVisViewerCtx_Private(PetscCtxRt vctx)
22d71ae5a4SJacob Faibussowitsch {
23e6aa7a3bSBarry Smith DMDAFieldGLVisViewerCtx *ctx = *(DMDAFieldGLVisViewerCtx **)vctx;
248135c375SStefano Zampini
258135c375SStefano Zampini PetscFunctionBegin;
269566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->xlocal));
27e6aa7a3bSBarry Smith PetscCall(PetscFree(ctx));
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 */
DMDAGetNumElementsGhosted(DM da,PetscInt * nex,PetscInt * ney,PetscInt * nez)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 */
DMDAGetNumVerticesGhosted(DM da,PetscInt * nvx,PetscInt * nvy,PetscInt * nvz)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
DMDASampleGLVisFields_Private(PetscObject oX,PetscInt nf,PetscObject oXf[],void * vctx)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;
9658b7e2c1SStefano Zampini PetscInt i, f, ii, ien, jen, ken, ie, je, ke, bs, *bss, *nfs;
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 }
11658b7e2c1SStefano Zampini PetscCall(PetscMalloc3(nf, &arrayf, nf, &bss, nf, &nfs));
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]));
12158b7e2c1SStefano Zampini PetscCall(VecGetLocalSize((Vec)oXf[f], &nfs[f]));
1228135c375SStefano Zampini }
1238135c375SStefano Zampini for (ke = kst, ii = 0; ke < kst + ken; ke++) {
1248135c375SStefano Zampini for (je = jst; je < jst + jen; je++) {
1258135c375SStefano Zampini for (ie = ist; ie < ist + ien; ie++) {
1268135c375SStefano Zampini PetscInt cf, b;
1278135c375SStefano Zampini i = ke * gm * gn + je * gm + ie;
12858b7e2c1SStefano Zampini for (f = 0, cf = 0; f < nf; f++) {
12958b7e2c1SStefano Zampini if (!nfs[f]) continue;
1309371c9d4SSatish Balay for (b = 0; b < bss[f]; b++) arrayf[f][bss[f] * ii + b] = array[i * bs + cf++];
13158b7e2c1SStefano Zampini }
1328135c375SStefano Zampini ii++;
1338135c375SStefano Zampini }
1348135c375SStefano Zampini }
1358135c375SStefano Zampini }
1369566063dSJacob Faibussowitsch for (f = 0; f < nf; f++) PetscCall(VecRestoreArray((Vec)oXf[f], &arrayf[f]));
1379566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(ctx->xlocal, &array));
13858b7e2c1SStefano Zampini PetscCall(PetscFree3(arrayf, bss, nfs));
1393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1408135c375SStefano Zampini }
1418135c375SStefano Zampini
DMSetUpGLVisViewer_DMDA(PetscObject oda,PetscViewer viewer)142ba38deedSJacob Faibussowitsch PetscErrorCode DMSetUpGLVisViewer_DMDA(PetscObject oda, PetscViewer viewer)
143d71ae5a4SJacob Faibussowitsch {
14446900a7fSStefano Zampini DM da = (DM)oda, daview;
1458135c375SStefano Zampini
1468135c375SStefano Zampini PetscFunctionBegin;
1479566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery(oda, "GLVisGraphicsDMDAGhosted", (PetscObject *)&daview));
14846900a7fSStefano Zampini if (!daview) {
14946900a7fSStefano Zampini DMDAGhostedGLVisViewerCtx *dactx;
15046900a7fSStefano Zampini DM dacoord = NULL;
15146900a7fSStefano Zampini Vec xcoor, xcoorl;
15246900a7fSStefano Zampini PetscBool hashocoord = PETSC_FALSE;
15346900a7fSStefano Zampini const PetscInt *lx, *ly, *lz;
15446900a7fSStefano Zampini PetscInt dim, M, N, P, m, n, p, dof, s, i;
15546900a7fSStefano Zampini
1569566063dSJacob Faibussowitsch PetscCall(PetscNew(&dactx));
157a7bda6aaSStefano Zampini dactx->ll = PETSC_TRUE; /* default to match elements layout obtained by DMDAGetElements */
158d0609cedSBarry Smith PetscOptionsBegin(PetscObjectComm(oda), oda->prefix, "GLVis PetscViewer DMDA Options", "PetscViewer");
1599566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-viewer_glvis_dm_da_ll", "Left-looking subdomain view", NULL, dactx->ll, &dactx->ll, NULL));
160d0609cedSBarry Smith PetscOptionsEnd();
1618135c375SStefano Zampini /* Create a properly ghosted DMDA to visualize the mesh and the fields associated with */
1629566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, &M, &N, &P, &m, &n, &p, &dof, &s, NULL, NULL, NULL, NULL));
1639566063dSJacob Faibussowitsch PetscCall(DMDAGetOwnershipRanges(da, &lx, &ly, &lz));
1649566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)da, "_glvis_mesh_coords", (PetscObject *)&xcoor));
16577eacf09SStefano Zampini if (!xcoor) {
1669566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(da, &xcoor));
16777eacf09SStefano Zampini } else {
1688f07809aSStefano Zampini hashocoord = PETSC_TRUE;
16977eacf09SStefano Zampini }
1706aad120cSJose E. Roman PetscCall(PetscInfo(da, "Creating auxiliary DMDA for managing GLVis graphics\n"));
1718135c375SStefano Zampini switch (dim) {
1728135c375SStefano Zampini case 1:
1739566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(PetscObjectComm((PetscObject)da), DM_BOUNDARY_NONE, M, dof, 1, lx, &daview));
17448a46eb9SPierre Jolivet if (!hashocoord) PetscCall(DMDACreate1d(PetscObjectComm((PetscObject)da), DM_BOUNDARY_NONE, M, 1, 1, lx, &dacoord));
1758135c375SStefano Zampini break;
1768135c375SStefano Zampini case 2:
1779566063dSJacob Faibussowitsch PetscCall(DMDACreate2d(PetscObjectComm((PetscObject)da), DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, M, N, m, n, dof, 1, lx, ly, &daview));
17848a46eb9SPierre 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));
1798135c375SStefano Zampini break;
1808135c375SStefano Zampini case 3:
1819566063dSJacob 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));
18248a46eb9SPierre 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));
1838135c375SStefano Zampini break;
184d71ae5a4SJacob Faibussowitsch default:
185d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Unsupported dimension %" PetscInt_FMT, dim);
1868135c375SStefano Zampini }
1879566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(daview, dactx));
1889566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContextDestroy(daview, DMDAGhostedDestroyGLVisViewerCtx_Private));
1899566063dSJacob Faibussowitsch PetscCall(DMSetUp(daview));
1908135c375SStefano Zampini if (!xcoor) {
1919566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(daview, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
1929566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(daview, &xcoor));
1938135c375SStefano Zampini }
1948f07809aSStefano Zampini if (dacoord) {
1959566063dSJacob Faibussowitsch PetscCall(DMSetUp(dacoord));
1969566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(dacoord, &xcoorl));
1979566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dacoord, xcoor, INSERT_VALUES, xcoorl));
1989566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dacoord, xcoor, INSERT_VALUES, xcoorl));
1999566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dacoord));
2008f07809aSStefano Zampini } else {
2018f07809aSStefano Zampini PetscInt ien, jen, ken, nc, nl, cdof, deg;
2023e7cab22SLisandro Dalcin char fecmesh[64];
20321414b21SStefano Zampini const char *name;
20421414b21SStefano Zampini PetscBool flg;
2058135c375SStefano Zampini
2069566063dSJacob Faibussowitsch PetscCall(DMDAGetNumElementsGhosted(daview, &ien, &jen, &ken));
2078f07809aSStefano Zampini nc = ien * (jen > 0 ? jen : 1) * (ken > 0 ? ken : 1);
2088f07809aSStefano Zampini
2099566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(xcoor, &nl));
2101dca8a05SBarry 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);
2119566063dSJacob Faibussowitsch PetscCall(VecDuplicate(xcoor, &xcoorl));
2129566063dSJacob Faibussowitsch PetscCall(VecCopy(xcoor, xcoorl));
2139566063dSJacob Faibussowitsch PetscCall(VecSetDM(xcoorl, NULL));
2149566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)xcoor, &name));
2159566063dSJacob Faibussowitsch PetscCall(PetscStrbeginswith(name, "FiniteElementCollection:", &flg));
21621414b21SStefano Zampini if (!flg) {
21721414b21SStefano Zampini deg = 0;
21821414b21SStefano Zampini if (nc && nl) {
2198f07809aSStefano Zampini cdof = nl / (nc * dim);
2208f07809aSStefano Zampini deg = 1;
2218f07809aSStefano Zampini while (1) {
2228f07809aSStefano Zampini PetscInt degd = 1;
2238f07809aSStefano Zampini for (i = 0; i < dim; i++) degd *= (deg + 1);
22463a3b9bcSJacob Faibussowitsch PetscCheck(degd <= cdof, PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell dofs %" PetscInt_FMT, cdof);
2258f07809aSStefano Zampini if (degd == cdof) break;
2268f07809aSStefano Zampini deg++;
2278f07809aSStefano Zampini }
22821414b21SStefano Zampini }
22963a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(fecmesh, sizeof(fecmesh), "FiniteElementCollection: L2_T1_%" PetscInt_FMT "D_P%" PetscInt_FMT, dim, deg));
2309566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)xcoorl, fecmesh));
2311baa6e33SBarry Smith } else PetscCall(PetscObjectSetName((PetscObject)xcoorl, name));
2328f07809aSStefano Zampini }
2338135c375SStefano Zampini
23446900a7fSStefano Zampini /* xcoorl is composed with the ghosted DMDA, the ghosted coordinate DMDA (if present) is only available through this vector */
2359566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)daview, "GLVisGraphicsCoordsGhosted", (PetscObject)xcoorl));
2369566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)xcoorl));
2378135c375SStefano Zampini
23877eacf09SStefano Zampini /* daview is composed with the original DMDA */
2399566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose(oda, "GLVisGraphicsDMDAGhosted", (PetscObject)daview));
2409566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)daview));
24146900a7fSStefano Zampini }
24277eacf09SStefano Zampini
243c9493c35SStefano Zampini /* customize the viewer if present */
244c9493c35SStefano Zampini if (viewer) {
24546900a7fSStefano Zampini DMDAFieldGLVisViewerCtx *ctx;
24646900a7fSStefano Zampini DMDAGhostedGLVisViewerCtx *dactx;
24746900a7fSStefano Zampini char fec[64];
248ab166b77SStefano Zampini Vec xlocal, *Ufield;
24946900a7fSStefano Zampini const char **dafieldname;
25046900a7fSStefano Zampini char **fec_type, **fieldname;
25146900a7fSStefano Zampini PetscInt *nlocal, *bss, *dims;
25246900a7fSStefano Zampini PetscInt dim, M, N, P, dof, s, i, nf;
25346900a7fSStefano Zampini PetscBool bsset;
2548f07809aSStefano Zampini
2559566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daview, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
2569566063dSJacob Faibussowitsch PetscCall(DMGetApplicationContext(daview, &dactx));
2579566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(daview, &xlocal));
2589566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldNames(da, (const char *const **)&dafieldname));
2599566063dSJacob Faibussowitsch PetscCall(DMDAGetNumVerticesGhosted(daview, &M, &N, &P));
26063a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(fec, sizeof(fec), "FiniteElementCollection: H1_%" PetscInt_FMT "D_P1", dim));
2619566063dSJacob Faibussowitsch PetscCall(PetscMalloc6(dof, &fec_type, dof, &nlocal, dof, &bss, dof, &dims, dof, &fieldname, dof, &Ufield));
2628135c375SStefano Zampini for (i = 0; i < dof; i++) bss[i] = 1;
2638135c375SStefano Zampini nf = dof;
2648135c375SStefano Zampini
265d0609cedSBarry Smith PetscOptionsBegin(PetscObjectComm(oda), oda->prefix, "GLVis PetscViewer DMDA Field options", "PetscViewer");
2669566063dSJacob Faibussowitsch PetscCall(PetscOptionsIntArray("-viewer_glvis_dm_da_bs", "Block sizes for subfields; enable vector representation", NULL, bss, &nf, &bsset));
267d0609cedSBarry Smith PetscOptionsEnd();
2688135c375SStefano Zampini if (bsset) {
2698135c375SStefano Zampini PetscInt t;
2708135c375SStefano Zampini for (i = 0, t = 0; i < nf; i++) t += bss[i];
27163a3b9bcSJacob Faibussowitsch PetscCheck(t == dof, PetscObjectComm(oda), PETSC_ERR_USER, "Sum of block sizes %" PetscInt_FMT " should equal %" PetscInt_FMT, t, dof);
2728135c375SStefano Zampini } else nf = dof;
2738135c375SStefano Zampini
2748135c375SStefano Zampini for (i = 0, s = 0; i < nf; i++) {
2759566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(fec, &fec_type[i]));
2768135c375SStefano Zampini if (bss[i] == 1) {
2779566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(dafieldname[s], &fieldname[i]));
2788135c375SStefano Zampini } else {
279c6a7a370SJeremy L Thompson const char prefix[] = "Vector-";
280c6a7a370SJeremy L Thompson const size_t prefix_len = PETSC_STATIC_ARRAY_LENGTH(prefix);
281c6a7a370SJeremy L Thompson const PetscInt bss_i = bss[i];
282c6a7a370SJeremy L Thompson char *fname_i = NULL;
283c6a7a370SJeremy L Thompson size_t tlen = prefix_len, cur_len;
284c6a7a370SJeremy L Thompson
285c6a7a370SJeremy L Thompson for (PetscInt b = 0; b < bss_i; ++b) {
2868135c375SStefano Zampini size_t len;
287c6a7a370SJeremy L Thompson
2889566063dSJacob Faibussowitsch PetscCall(PetscStrlen(dafieldname[s + b], &len));
2898135c375SStefano Zampini tlen += len + 1; /* field + "-" */
2908135c375SStefano Zampini }
291c6a7a370SJeremy L Thompson PetscCall(PetscMalloc1(tlen, &fname_i));
292c6a7a370SJeremy L Thompson PetscCall(PetscArraycpy(fname_i, prefix, prefix_len - 1));
293c6a7a370SJeremy L Thompson cur_len = prefix_len;
294c6a7a370SJeremy L Thompson for (PetscInt b = 0; b < bss_i; ++b) {
295c6a7a370SJeremy L Thompson const char *dafname = dafieldname[s + b];
296c6a7a370SJeremy L Thompson size_t len;
297c6a7a370SJeremy L Thompson
298c6a7a370SJeremy L Thompson PetscCall(PetscStrlen(dafname, &len));
299363da2dcSJacob Faibussowitsch PetscCall(PetscArraycpy(fname_i + cur_len, dafname, len + 1));
300c6a7a370SJeremy L Thompson cur_len += len + 1;
301c6a7a370SJeremy L Thompson // don't add for final iteration of the loop
302c6a7a370SJeremy L Thompson if ((b + 1) < bss_i) fname_i[cur_len++] = '-';
3038135c375SStefano Zampini }
304c6a7a370SJeremy L Thompson fname_i[cur_len] = '\0';
305c6a7a370SJeremy L Thompson fieldname[i] = fname_i;
3068135c375SStefano Zampini }
307bb77a09fSStefano Zampini dims[i] = dim;
3088135c375SStefano Zampini nlocal[i] = M * N * P * bss[i];
3098135c375SStefano Zampini s += bss[i];
3108135c375SStefano Zampini }
3118135c375SStefano Zampini
31246900a7fSStefano Zampini /* the viewer context takes ownership of xlocal and destroys it in DMDAFieldDestroyGLVisViewerCtx_Private */
3139566063dSJacob Faibussowitsch PetscCall(PetscNew(&ctx));
3148135c375SStefano Zampini ctx->xlocal = xlocal;
3158135c375SStefano Zampini
3164cac2994SStefano Zampini /* create work vectors */
3174cac2994SStefano Zampini for (i = 0; i < nf; i++) {
3189566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)da), nlocal[i], PETSC_DECIDE, &Ufield[i]));
3199566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)Ufield[i], fieldname[i]));
32058b7e2c1SStefano Zampini PetscCall(VecSetBlockSize(Ufield[i], PetscMax(bss[i], 1)));
3219566063dSJacob Faibussowitsch PetscCall(VecSetDM(Ufield[i], da));
3224cac2994SStefano Zampini }
3234cac2994SStefano Zampini
3249566063dSJacob Faibussowitsch PetscCall(PetscViewerGLVisSetFields(viewer, nf, (const char **)fec_type, dims, DMDASampleGLVisFields_Private, (PetscObject *)Ufield, ctx, DMDAFieldDestroyGLVisViewerCtx_Private));
3258135c375SStefano Zampini for (i = 0; i < nf; i++) {
3269566063dSJacob Faibussowitsch PetscCall(PetscFree(fec_type[i]));
3279566063dSJacob Faibussowitsch PetscCall(PetscFree(fieldname[i]));
3289566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Ufield[i]));
3298135c375SStefano Zampini }
3309566063dSJacob Faibussowitsch PetscCall(PetscFree6(fec_type, nlocal, bss, dims, fieldname, Ufield));
331c9493c35SStefano Zampini }
3323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3338135c375SStefano Zampini }
3348135c375SStefano Zampini
DMDAView_GLVis_ASCII(DM dm,PetscViewer viewer)335d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDAView_GLVis_ASCII(DM dm, PetscViewer viewer)
336d71ae5a4SJacob Faibussowitsch {
33777eacf09SStefano Zampini DM da, cda;
3388135c375SStefano Zampini Vec xcoorl;
3393924b612SStefano Zampini PetscMPIInt size;
3408135c375SStefano Zampini const PetscScalar *array;
3418135c375SStefano Zampini PetscContainer glvis_container;
34277eacf09SStefano Zampini PetscInt dim, sdim, i, vid[8], mid, cid, cdof;
34321414b21SStefano Zampini PetscInt sx, sy, sz, ie, je, ke, ien, jen, ken, nel;
3448135c375SStefano Zampini PetscInt gsx, gsy, gsz, gm, gn, gp, kst, jst, ist;
3458135c375SStefano Zampini PetscBool enabled = PETSC_TRUE, isascii;
3468f07809aSStefano Zampini const char *fmt;
3478135c375SStefano Zampini
3488135c375SStefano Zampini PetscFunctionBegin;
3498135c375SStefano Zampini PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3508135c375SStefano Zampini PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
3519566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
3527a8be351SBarry Smith PetscCheck(isascii, PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Viewer must be of type VIEWERASCII");
3539566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)viewer), &size));
35408401ef6SPierre Jolivet PetscCheck(size <= 1, PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Use single sequential viewers for parallel visualization");
3559566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim));
3568135c375SStefano Zampini
3578135c375SStefano Zampini /* get container: determines if a process visualizes is portion of the data or not */
3589566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)viewer, "_glvis_info_container", (PetscObject *)&glvis_container));
3597a8be351SBarry Smith PetscCheck(glvis_container, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Missing GLVis container");
3608135c375SStefano Zampini {
3618135c375SStefano Zampini PetscViewerGLVisInfo glvis_info;
362*2a8381b2SBarry Smith PetscCall(PetscContainerGetPointer(glvis_container, &glvis_info));
3638135c375SStefano Zampini enabled = glvis_info->enabled;
36477eacf09SStefano Zampini fmt = glvis_info->fmt;
3658135c375SStefano Zampini }
366c9493c35SStefano Zampini /* this can happen if we are calling DMView outside of VecView_GLVis */
3679566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)dm, "GLVisGraphicsDMDAGhosted", (PetscObject *)&da));
3689566063dSJacob Faibussowitsch if (!da) PetscCall(DMSetUpGLVisViewer_DMDA((PetscObject)dm, NULL));
3699566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)dm, "GLVisGraphicsDMDAGhosted", (PetscObject *)&da));
3707a8be351SBarry Smith PetscCheck(da, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Missing GLVis ghosted DMDA");
3719566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(da, &sdim));
3728135c375SStefano Zampini
3739566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "MFEM mesh v1.0\n"));
3749566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\ndimension\n"));
37563a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT "\n", dim));
3768135c375SStefano Zampini
3778135c375SStefano Zampini if (!enabled) {
3789566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\nelements\n"));
37963a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "0\n"));
3809566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\nboundary\n"));
38163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "0\n"));
3829566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\nvertices\n"));
38363a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "0\n"));
38463a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT "\n", sdim));
3853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3868135c375SStefano Zampini }
3878135c375SStefano Zampini
3889566063dSJacob Faibussowitsch PetscCall(DMDAGetNumElementsGhosted(da, &ien, &jen, &ken));
38921414b21SStefano Zampini nel = ien;
39021414b21SStefano Zampini if (dim > 1) nel *= jen;
39121414b21SStefano Zampini if (dim > 2) nel *= ken;
3929566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\nelements\n"));
39363a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT "\n", nel));
3948135c375SStefano Zampini switch (dim) {
3958135c375SStefano Zampini case 1:
3968135c375SStefano Zampini for (ie = 0; ie < ien; ie++) {
3978135c375SStefano Zampini vid[0] = ie;
3988135c375SStefano Zampini vid[1] = ie + 1;
3998135c375SStefano Zampini mid = 1; /* material id */
4008135c375SStefano Zampini cid = 1; /* segment */
40163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n", mid, cid, vid[0], vid[1]));
4028135c375SStefano Zampini }
4038135c375SStefano Zampini break;
4048135c375SStefano Zampini case 2:
4058135c375SStefano Zampini for (je = 0; je < jen; je++) {
4068135c375SStefano Zampini for (ie = 0; ie < ien; ie++) {
4078135c375SStefano Zampini vid[0] = je * (ien + 1) + ie;
4088135c375SStefano Zampini vid[1] = je * (ien + 1) + ie + 1;
4098135c375SStefano Zampini vid[2] = (je + 1) * (ien + 1) + ie + 1;
4108135c375SStefano Zampini vid[3] = (je + 1) * (ien + 1) + ie;
4118135c375SStefano Zampini mid = 1; /* material id */
4128135c375SStefano Zampini cid = 3; /* quad */
41363a3b9bcSJacob 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]));
4148135c375SStefano Zampini }
4158135c375SStefano Zampini }
4168135c375SStefano Zampini break;
4178135c375SStefano Zampini case 3:
4188135c375SStefano Zampini for (ke = 0; ke < ken; ke++) {
4198135c375SStefano Zampini for (je = 0; je < jen; je++) {
4208135c375SStefano Zampini for (ie = 0; ie < ien; ie++) {
4218135c375SStefano Zampini vid[0] = ke * (jen + 1) * (ien + 1) + je * (ien + 1) + ie;
4228135c375SStefano Zampini vid[1] = ke * (jen + 1) * (ien + 1) + je * (ien + 1) + ie + 1;
4238135c375SStefano Zampini vid[2] = ke * (jen + 1) * (ien + 1) + (je + 1) * (ien + 1) + ie + 1;
4248135c375SStefano Zampini vid[3] = ke * (jen + 1) * (ien + 1) + (je + 1) * (ien + 1) + ie;
4258135c375SStefano Zampini vid[4] = (ke + 1) * (jen + 1) * (ien + 1) + je * (ien + 1) + ie;
4268135c375SStefano Zampini vid[5] = (ke + 1) * (jen + 1) * (ien + 1) + je * (ien + 1) + ie + 1;
4278135c375SStefano Zampini vid[6] = (ke + 1) * (jen + 1) * (ien + 1) + (je + 1) * (ien + 1) + ie + 1;
4288135c375SStefano Zampini vid[7] = (ke + 1) * (jen + 1) * (ien + 1) + (je + 1) * (ien + 1) + ie;
4298135c375SStefano Zampini mid = 1; /* material id */
4308135c375SStefano Zampini cid = 5; /* hex */
43163a3b9bcSJacob 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]));
4328135c375SStefano Zampini }
4338135c375SStefano Zampini }
4348135c375SStefano Zampini }
4358135c375SStefano Zampini break;
436d71ae5a4SJacob Faibussowitsch default:
437d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_SUP, "Unsupported dimension %" PetscInt_FMT, dim);
4388135c375SStefano Zampini }
4399566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\nboundary\n"));
44063a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "0\n"));
4418135c375SStefano Zampini
44277eacf09SStefano Zampini /* vertex coordinates */
4439566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)da, "GLVisGraphicsCoordsGhosted", (PetscObject *)&xcoorl));
4447a8be351SBarry Smith PetscCheck(xcoorl, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Missing GLVis ghosted coords");
4459566063dSJacob Faibussowitsch PetscCall(DMDAGetNumVerticesGhosted(da, &ien, &jen, &ken));
4469566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\nvertices\n"));
44763a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT "\n", ien * jen * ken));
44821414b21SStefano Zampini if (nel) {
4499566063dSJacob Faibussowitsch PetscCall(VecGetDM(xcoorl, &cda));
4509566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(xcoorl, &array));
4518f07809aSStefano Zampini if (!cda) { /* HO viz */
4528f07809aSStefano Zampini const char *fecname;
4538f07809aSStefano Zampini PetscInt nc, nl;
4548f07809aSStefano Zampini
4559566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)xcoorl, &fecname));
4569566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "nodes\n"));
4579566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "FiniteElementSpace\n"));
4589566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%s\n", fecname));
45963a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "VDim: %" PetscInt_FMT "\n", sdim));
4609566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Ordering: 1\n\n")); /*Ordering::byVDIM*/
4618f07809aSStefano Zampini /* L2 coordinates */
4629566063dSJacob Faibussowitsch PetscCall(DMDAGetNumElementsGhosted(da, &ien, &jen, &ken));
4639566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(xcoorl, &nl));
4648f07809aSStefano Zampini nc = ien * (jen > 0 ? jen : 1) * (ken > 0 ? ken : 1);
46521414b21SStefano Zampini cdof = nc ? nl / nc : 0;
4668f07809aSStefano Zampini if (!ien) ien++;
4678f07809aSStefano Zampini if (!jen) jen++;
4688f07809aSStefano Zampini if (!ken) ken++;
4698f07809aSStefano Zampini ist = jst = kst = 0;
4708f07809aSStefano Zampini gm = ien;
4718f07809aSStefano Zampini gn = jen;
4728f07809aSStefano Zampini gp = ken;
47377eacf09SStefano Zampini } else {
47446900a7fSStefano Zampini DMDAGhostedGLVisViewerCtx *dactx;
47546900a7fSStefano Zampini
4769566063dSJacob Faibussowitsch PetscCall(DMGetApplicationContext(da, &dactx));
47763a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT "\n", sdim));
4788f07809aSStefano Zampini cdof = sdim;
4799566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &sx, &sy, &sz, NULL, NULL, NULL));
4809566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, &gsx, &gsy, &gsz, &gm, &gn, &gp));
48146900a7fSStefano Zampini if (dactx->ll) {
48246900a7fSStefano Zampini kst = jst = ist = 0;
48346900a7fSStefano Zampini } else {
4848135c375SStefano Zampini kst = gsz != sz ? 1 : 0;
4858135c375SStefano Zampini jst = gsy != sy ? 1 : 0;
4868135c375SStefano Zampini ist = gsx != sx ? 1 : 0;
4878f07809aSStefano Zampini }
48846900a7fSStefano Zampini }
4898135c375SStefano Zampini for (ke = kst; ke < kst + ken; ke++) {
4908135c375SStefano Zampini for (je = jst; je < jst + jen; je++) {
4918135c375SStefano Zampini for (ie = ist; ie < ist + ien; ie++) {
49277eacf09SStefano Zampini PetscInt c;
4938135c375SStefano Zampini
4948135c375SStefano Zampini i = ke * gm * gn + je * gm + ie;
49577eacf09SStefano Zampini for (c = 0; c < cdof / sdim; c++) {
49677eacf09SStefano Zampini PetscInt d;
49748a46eb9SPierre Jolivet for (d = 0; d < sdim; d++) PetscCall(PetscViewerASCIIPrintf(viewer, fmt, PetscRealPart(array[cdof * i + c * sdim + d])));
4989566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
49977eacf09SStefano Zampini }
5008135c375SStefano Zampini }
5018135c375SStefano Zampini }
5028135c375SStefano Zampini }
5039566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(xcoorl, &array));
50421414b21SStefano Zampini }
5053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
5068135c375SStefano Zampini }
5078135c375SStefano Zampini
DMView_DA_GLVis(DM dm,PetscViewer viewer)508d71ae5a4SJacob Faibussowitsch PetscErrorCode DMView_DA_GLVis(DM dm, PetscViewer viewer)
509d71ae5a4SJacob Faibussowitsch {
5108135c375SStefano Zampini PetscFunctionBegin;
5119566063dSJacob Faibussowitsch PetscCall(DMView_GLVis(dm, viewer, DMDAView_GLVis_ASCII));
5123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
5138135c375SStefano Zampini }
514