149c89c76SBlaise Bourdin #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/
249c89c76SBlaise Bourdin
349c89c76SBlaise Bourdin #include <netcdf.h>
449c89c76SBlaise Bourdin #include <exodusII.h>
549c89c76SBlaise Bourdin
649c89c76SBlaise Bourdin #include <petsc/private/viewerimpl.h>
749c89c76SBlaise Bourdin #include <petsc/private/viewerexodusiiimpl.h>
849c89c76SBlaise Bourdin /*@C
949c89c76SBlaise Bourdin PETSC_VIEWER_EXODUSII_ - Creates an `PETSCVIEWEREXODUSII` `PetscViewer` shared by all processors in a communicator.
1049c89c76SBlaise Bourdin
1149c89c76SBlaise Bourdin Collective; No Fortran Support
1249c89c76SBlaise Bourdin
1349c89c76SBlaise Bourdin Input Parameter:
1449c89c76SBlaise Bourdin . comm - the MPI communicator to share the `PETSCVIEWEREXODUSII` `PetscViewer`
1549c89c76SBlaise Bourdin
1649c89c76SBlaise Bourdin Level: intermediate
1749c89c76SBlaise Bourdin
1849c89c76SBlaise Bourdin Note:
1949c89c76SBlaise Bourdin Unlike almost all other PETSc routines, `PETSC_VIEWER_EXODUSII_()` does not return
2049c89c76SBlaise Bourdin an error code. The GLVIS PetscViewer is usually used in the form
21b44f4de4SBarry Smith .vb
22b44f4de4SBarry Smith XXXView(XXX object, PETSC_VIEWER_EXODUSII_(comm));
23b44f4de4SBarry Smith .ve
2449c89c76SBlaise Bourdin
2549c89c76SBlaise Bourdin .seealso: `PETSCVIEWEREXODUSII`, `PetscViewer`, `PetscViewer`, `PetscViewerExodusIIOpen()`, `PetscViewerType`, `PetscViewerCreate()`, `PetscViewerDestroy()`
2649c89c76SBlaise Bourdin @*/
PETSC_VIEWER_EXODUSII_(MPI_Comm comm)2749c89c76SBlaise Bourdin PetscViewer PETSC_VIEWER_EXODUSII_(MPI_Comm comm)
2849c89c76SBlaise Bourdin {
2949c89c76SBlaise Bourdin PetscViewer viewer;
3049c89c76SBlaise Bourdin
3149c89c76SBlaise Bourdin PetscFunctionBegin;
3249c89c76SBlaise Bourdin PetscCallNull(PetscViewerExodusIIOpen(comm, "mesh.exo", FILE_MODE_WRITE, &viewer));
3349c89c76SBlaise Bourdin PetscCallNull(PetscObjectRegisterDestroy((PetscObject)viewer));
3449c89c76SBlaise Bourdin PetscFunctionReturn(viewer);
3549c89c76SBlaise Bourdin }
3649c89c76SBlaise Bourdin
PetscViewerView_ExodusII(PetscViewer v,PetscViewer viewer)3749c89c76SBlaise Bourdin static PetscErrorCode PetscViewerView_ExodusII(PetscViewer v, PetscViewer viewer)
3849c89c76SBlaise Bourdin {
3949c89c76SBlaise Bourdin PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)v->data;
4049c89c76SBlaise Bourdin
4149c89c76SBlaise Bourdin PetscFunctionBegin;
4249c89c76SBlaise Bourdin if (exo->filename) PetscCall(PetscViewerASCIIPrintf(viewer, "Filename: %s\n", exo->filename));
430a5cf5d8SBlaise Bourdin if (exo->exoid) PetscCall(PetscViewerASCIIPrintf(viewer, "exoid: %" PetscExodusIIInt_FMT "\n", exo->exoid));
4449c89c76SBlaise Bourdin if (exo->btype) PetscCall(PetscViewerASCIIPrintf(viewer, "IO Mode: %d\n", exo->btype));
450a5cf5d8SBlaise Bourdin if (exo->order) PetscCall(PetscViewerASCIIPrintf(viewer, "Mesh order: %" PetscInt_FMT "\n", exo->order));
460a5cf5d8SBlaise Bourdin PetscCall(PetscViewerASCIIPrintf(viewer, "Number of nodal variables: %" PetscExodusIIInt_FMT "\n", exo->numNodalVariables));
4749c89c76SBlaise Bourdin for (int i = 0; i < exo->numNodalVariables; i++) PetscCall(PetscViewerASCIIPrintf(viewer, " %d: %s\n", i, exo->nodalVariableNames[i]));
480a5cf5d8SBlaise Bourdin PetscCall(PetscViewerASCIIPrintf(viewer, "Number of zonal variables: %" PetscExodusIIInt_FMT "\n", exo->numZonalVariables));
4949c89c76SBlaise Bourdin for (int i = 0; i < exo->numZonalVariables; i++) PetscCall(PetscViewerASCIIPrintf(viewer, " %d: %s\n", i, exo->zonalVariableNames[i]));
5049c89c76SBlaise Bourdin PetscFunctionReturn(PETSC_SUCCESS);
5149c89c76SBlaise Bourdin }
5249c89c76SBlaise Bourdin
PetscViewerFlush_ExodusII(PetscViewer v)5349c89c76SBlaise Bourdin static PetscErrorCode PetscViewerFlush_ExodusII(PetscViewer v)
5449c89c76SBlaise Bourdin {
5549c89c76SBlaise Bourdin PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)v->data;
5649c89c76SBlaise Bourdin
5749c89c76SBlaise Bourdin PetscFunctionBegin;
5849c89c76SBlaise Bourdin if (exo->exoid >= 0) PetscCallExternal(ex_update, exo->exoid);
5949c89c76SBlaise Bourdin PetscFunctionReturn(PETSC_SUCCESS);
6049c89c76SBlaise Bourdin }
6149c89c76SBlaise Bourdin
PetscViewerSetFromOptions_ExodusII(PetscViewer v,PetscOptionItems PetscOptionsObject)62ce78bad3SBarry Smith static PetscErrorCode PetscViewerSetFromOptions_ExodusII(PetscViewer v, PetscOptionItems PetscOptionsObject)
6349c89c76SBlaise Bourdin {
6449c89c76SBlaise Bourdin PetscFunctionBegin;
6549c89c76SBlaise Bourdin PetscOptionsHeadBegin(PetscOptionsObject, "ExodusII PetscViewer Options");
6649c89c76SBlaise Bourdin PetscOptionsHeadEnd();
6749c89c76SBlaise Bourdin PetscFunctionReturn(PETSC_SUCCESS);
6849c89c76SBlaise Bourdin }
6949c89c76SBlaise Bourdin
PetscViewerDestroy_ExodusII(PetscViewer viewer)7049c89c76SBlaise Bourdin static PetscErrorCode PetscViewerDestroy_ExodusII(PetscViewer viewer)
7149c89c76SBlaise Bourdin {
7249c89c76SBlaise Bourdin PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data;
7349c89c76SBlaise Bourdin
7449c89c76SBlaise Bourdin PetscFunctionBegin;
7549c89c76SBlaise Bourdin if (exo->exoid >= 0) PetscCallExternal(ex_close, exo->exoid);
765341eb2eSStefano Zampini for (PetscInt i = 0; i < exo->numZonalVariables; i++) PetscCall(PetscFree(exo->zonalVariableNames[i]));
7749c89c76SBlaise Bourdin PetscCall(PetscFree(exo->zonalVariableNames));
785341eb2eSStefano Zampini for (PetscInt i = 0; i < exo->numNodalVariables; i++) PetscCall(PetscFree(exo->nodalVariableNames[i]));
7949c89c76SBlaise Bourdin PetscCall(PetscFree(exo->nodalVariableNames));
8049c89c76SBlaise Bourdin PetscCall(PetscFree(exo->filename));
8149c89c76SBlaise Bourdin PetscCall(PetscFree(exo));
8249c89c76SBlaise Bourdin PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileSetName_C", NULL));
8349c89c76SBlaise Bourdin PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileGetName_C", NULL));
8449c89c76SBlaise Bourdin PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileSetMode_C", NULL));
8549c89c76SBlaise Bourdin PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileGetMode_C", NULL));
8649c89c76SBlaise Bourdin PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerGetId_C", NULL));
8749c89c76SBlaise Bourdin PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerGetOrder_C", NULL));
8849c89c76SBlaise Bourdin PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerSetOrder_C", NULL));
8949c89c76SBlaise Bourdin PetscFunctionReturn(PETSC_SUCCESS);
9049c89c76SBlaise Bourdin }
9149c89c76SBlaise Bourdin
PetscViewerFileSetName_ExodusII(PetscViewer viewer,const char name[])9249c89c76SBlaise Bourdin static PetscErrorCode PetscViewerFileSetName_ExodusII(PetscViewer viewer, const char name[])
9349c89c76SBlaise Bourdin {
9449c89c76SBlaise Bourdin PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data;
9549c89c76SBlaise Bourdin PetscMPIInt rank;
960a5cf5d8SBlaise Bourdin PetscExodusIIInt CPU_word_size, IO_word_size, EXO_mode;
9749c89c76SBlaise Bourdin MPI_Info mpi_info = MPI_INFO_NULL;
980a5cf5d8SBlaise Bourdin PetscExodusIIFloat EXO_version;
9949c89c76SBlaise Bourdin
10049c89c76SBlaise Bourdin PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
10149c89c76SBlaise Bourdin CPU_word_size = sizeof(PetscReal);
10249c89c76SBlaise Bourdin IO_word_size = sizeof(PetscReal);
10349c89c76SBlaise Bourdin
10449c89c76SBlaise Bourdin PetscFunctionBegin;
10549c89c76SBlaise Bourdin if (exo->exoid >= 0) {
10649c89c76SBlaise Bourdin PetscCallExternal(ex_close, exo->exoid);
10749c89c76SBlaise Bourdin exo->exoid = -1;
10849c89c76SBlaise Bourdin }
10949c89c76SBlaise Bourdin if (exo->filename) PetscCall(PetscFree(exo->filename));
11049c89c76SBlaise Bourdin PetscCall(PetscStrallocpy(name, &exo->filename));
11149c89c76SBlaise Bourdin switch (exo->btype) {
11249c89c76SBlaise Bourdin case FILE_MODE_READ:
11349c89c76SBlaise Bourdin EXO_mode = EX_READ;
11449c89c76SBlaise Bourdin break;
11549c89c76SBlaise Bourdin case FILE_MODE_APPEND:
11649c89c76SBlaise Bourdin case FILE_MODE_UPDATE:
11749c89c76SBlaise Bourdin case FILE_MODE_APPEND_UPDATE:
11849c89c76SBlaise Bourdin /* Will fail if the file does not already exist */
11949c89c76SBlaise Bourdin EXO_mode = EX_WRITE;
12049c89c76SBlaise Bourdin break;
12149c89c76SBlaise Bourdin case FILE_MODE_WRITE:
12249c89c76SBlaise Bourdin /*
12349c89c76SBlaise Bourdin exodus only allows writing geometry upon file creation, so we will let DMView create the file.
12449c89c76SBlaise Bourdin */
12549c89c76SBlaise Bourdin PetscFunctionReturn(PETSC_SUCCESS);
12649c89c76SBlaise Bourdin break;
12749c89c76SBlaise Bourdin default:
12849c89c76SBlaise Bourdin SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
12949c89c76SBlaise Bourdin }
13049c89c76SBlaise Bourdin #if defined(PETSC_USE_64BIT_INDICES)
13149c89c76SBlaise Bourdin EXO_mode += EX_ALL_INT64_API;
13249c89c76SBlaise Bourdin #endif
13349c89c76SBlaise Bourdin exo->exoid = ex_open_par(name, EXO_mode, &CPU_word_size, &IO_word_size, &EXO_version, PETSC_COMM_WORLD, mpi_info);
13449c89c76SBlaise Bourdin PetscCheck(exo->exoid >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_open_par failed for %s", name);
13549c89c76SBlaise Bourdin PetscFunctionReturn(PETSC_SUCCESS);
13649c89c76SBlaise Bourdin }
13749c89c76SBlaise Bourdin
PetscViewerFileGetName_ExodusII(PetscViewer viewer,const char * name[])13849c89c76SBlaise Bourdin static PetscErrorCode PetscViewerFileGetName_ExodusII(PetscViewer viewer, const char *name[])
13949c89c76SBlaise Bourdin {
14049c89c76SBlaise Bourdin PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data;
14149c89c76SBlaise Bourdin
14249c89c76SBlaise Bourdin PetscFunctionBegin;
14349c89c76SBlaise Bourdin *name = exo->filename;
14449c89c76SBlaise Bourdin PetscFunctionReturn(PETSC_SUCCESS);
14549c89c76SBlaise Bourdin }
14649c89c76SBlaise Bourdin
PetscViewerFileSetMode_ExodusII(PetscViewer viewer,PetscFileMode type)14749c89c76SBlaise Bourdin static PetscErrorCode PetscViewerFileSetMode_ExodusII(PetscViewer viewer, PetscFileMode type)
14849c89c76SBlaise Bourdin {
14949c89c76SBlaise Bourdin PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data;
15049c89c76SBlaise Bourdin
15149c89c76SBlaise Bourdin PetscFunctionBegin;
15249c89c76SBlaise Bourdin exo->btype = type;
15349c89c76SBlaise Bourdin PetscFunctionReturn(PETSC_SUCCESS);
15449c89c76SBlaise Bourdin }
15549c89c76SBlaise Bourdin
PetscViewerFileGetMode_ExodusII(PetscViewer viewer,PetscFileMode * type)15649c89c76SBlaise Bourdin static PetscErrorCode PetscViewerFileGetMode_ExodusII(PetscViewer viewer, PetscFileMode *type)
15749c89c76SBlaise Bourdin {
15849c89c76SBlaise Bourdin PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data;
15949c89c76SBlaise Bourdin
16049c89c76SBlaise Bourdin PetscFunctionBegin;
16149c89c76SBlaise Bourdin *type = exo->btype;
16249c89c76SBlaise Bourdin PetscFunctionReturn(PETSC_SUCCESS);
16349c89c76SBlaise Bourdin }
16449c89c76SBlaise Bourdin
PetscViewerExodusIIGetId_ExodusII(PetscViewer viewer,PetscExodusIIInt * exoid)1650a5cf5d8SBlaise Bourdin static PetscErrorCode PetscViewerExodusIIGetId_ExodusII(PetscViewer viewer, PetscExodusIIInt *exoid)
16649c89c76SBlaise Bourdin {
16749c89c76SBlaise Bourdin PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data;
16849c89c76SBlaise Bourdin
16949c89c76SBlaise Bourdin PetscFunctionBegin;
17049c89c76SBlaise Bourdin *exoid = exo->exoid;
17149c89c76SBlaise Bourdin PetscFunctionReturn(PETSC_SUCCESS);
17249c89c76SBlaise Bourdin }
17349c89c76SBlaise Bourdin
PetscViewerExodusIIGetOrder_ExodusII(PetscViewer viewer,PetscInt * order)1740a5cf5d8SBlaise Bourdin static PetscErrorCode PetscViewerExodusIIGetOrder_ExodusII(PetscViewer viewer, PetscInt *order)
17549c89c76SBlaise Bourdin {
17649c89c76SBlaise Bourdin PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data;
17749c89c76SBlaise Bourdin
17849c89c76SBlaise Bourdin PetscFunctionBegin;
17949c89c76SBlaise Bourdin *order = exo->order;
18049c89c76SBlaise Bourdin PetscFunctionReturn(PETSC_SUCCESS);
18149c89c76SBlaise Bourdin }
18249c89c76SBlaise Bourdin
PetscViewerExodusIISetOrder_ExodusII(PetscViewer viewer,PetscInt order)1830a5cf5d8SBlaise Bourdin static PetscErrorCode PetscViewerExodusIISetOrder_ExodusII(PetscViewer viewer, PetscInt order)
18449c89c76SBlaise Bourdin {
18549c89c76SBlaise Bourdin PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data;
18649c89c76SBlaise Bourdin
18749c89c76SBlaise Bourdin PetscFunctionBegin;
18849c89c76SBlaise Bourdin exo->order = order;
18949c89c76SBlaise Bourdin PetscFunctionReturn(PETSC_SUCCESS);
19049c89c76SBlaise Bourdin }
19149c89c76SBlaise Bourdin
19249c89c76SBlaise Bourdin /*@
193caff39ffSPierre Jolivet PetscViewerExodusIISetZonalVariable - Sets the number of zonal variables in an ExodusII file
19449c89c76SBlaise Bourdin
19549c89c76SBlaise Bourdin Collective;
19649c89c76SBlaise Bourdin
19749c89c76SBlaise Bourdin Input Parameters:
19849c89c76SBlaise Bourdin + viewer - a `PetscViewer` of type `PETSCVIEWEREXODUSII`
199caff39ffSPierre Jolivet - num - the number of zonal variables in the ExodusII file
20049c89c76SBlaise Bourdin
20149c89c76SBlaise Bourdin Level: intermediate
20249c89c76SBlaise Bourdin
20349c89c76SBlaise Bourdin Notes:
204caff39ffSPierre Jolivet The ExodusII API does not allow changing the number of variables in a file so this function will return an error
20549c89c76SBlaise Bourdin if called twice, called on a read-only file, or called on file for which the number of variables has already been specified
20649c89c76SBlaise Bourdin
20749c89c76SBlaise Bourdin .seealso: `PETSCVIEWEREXODUSII`, `PetscViewer`, `PetscViewerCreate()`, `PetscViewerDestroy()`, `PetscViewerExodusIIOpen()`, `PetscViewerSetType()`, `PetscViewerType`, `PetscViewerExodusIIGetZonalVariable()`
20849c89c76SBlaise Bourdin @*/
PetscViewerExodusIISetZonalVariable(PetscViewer viewer,PetscExodusIIInt num)2090a5cf5d8SBlaise Bourdin PetscErrorCode PetscViewerExodusIISetZonalVariable(PetscViewer viewer, PetscExodusIIInt num)
21049c89c76SBlaise Bourdin {
21149c89c76SBlaise Bourdin PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data;
21249c89c76SBlaise Bourdin MPI_Comm comm;
2130a5cf5d8SBlaise Bourdin PetscExodusIIInt exoid = -1;
21449c89c76SBlaise Bourdin
21549c89c76SBlaise Bourdin PetscFunctionBegin;
21649c89c76SBlaise Bourdin PetscCall(PetscObjectGetComm((PetscObject)viewer, &comm));
21749c89c76SBlaise Bourdin PetscCheck(exo->numZonalVariables == -1, comm, PETSC_ERR_SUP, "The number of zonal variables has already been set to %d and cannot be overwritten", exo->numZonalVariables);
21849c89c76SBlaise Bourdin PetscCheck((exo->btype != FILE_MODE_READ) && (exo->btype != FILE_MODE_UNDEFINED), comm, PETSC_ERR_FILE_WRITE, "Cannot set the number of variables because the file is not writable");
21949c89c76SBlaise Bourdin
22049c89c76SBlaise Bourdin exo->numZonalVariables = num;
22149c89c76SBlaise Bourdin PetscCall(PetscMalloc1(num, &exo->zonalVariableNames));
222ac530a7eSPierre Jolivet for (int i = 0; i < num; i++) exo->zonalVariableNames[i] = NULL;
22349c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIIGetId(viewer, &exoid));
22449c89c76SBlaise Bourdin PetscCallExternal(ex_put_variable_param, exoid, EX_ELEM_BLOCK, num);
22549c89c76SBlaise Bourdin PetscFunctionReturn(PETSC_SUCCESS);
22649c89c76SBlaise Bourdin }
22749c89c76SBlaise Bourdin
22849c89c76SBlaise Bourdin /*@
229caff39ffSPierre Jolivet PetscViewerExodusIISetNodalVariable - Sets the number of nodal variables in an ExodusII file
23049c89c76SBlaise Bourdin
23149c89c76SBlaise Bourdin Collective;
23249c89c76SBlaise Bourdin
23349c89c76SBlaise Bourdin Input Parameters:
23449c89c76SBlaise Bourdin + viewer - a `PetscViewer` of type `PETSCVIEWEREXODUSII`
235caff39ffSPierre Jolivet - num - the number of nodal variables in the ExodusII file
23649c89c76SBlaise Bourdin
23749c89c76SBlaise Bourdin Level: intermediate
23849c89c76SBlaise Bourdin
23949c89c76SBlaise Bourdin Notes:
240caff39ffSPierre Jolivet The ExodusII API does not allow changing the number of variables in a file so this function will return an error
24149c89c76SBlaise Bourdin if called twice, called on a read-only file, or called on file for which the number of variables has already been specified
24249c89c76SBlaise Bourdin
24349c89c76SBlaise Bourdin .seealso: `PETSCVIEWEREXODUSII`, `PetscViewer`, `PetscViewerCreate()`, `PetscViewerDestroy()`, `PetscViewerExodusIIOpen()`, `PetscViewerSetType()`, `PetscViewerType`, `PetscViewerExodusIIGetNodalVariable()`
24449c89c76SBlaise Bourdin @*/
PetscViewerExodusIISetNodalVariable(PetscViewer viewer,PetscExodusIIInt num)2450a5cf5d8SBlaise Bourdin PetscErrorCode PetscViewerExodusIISetNodalVariable(PetscViewer viewer, PetscExodusIIInt num)
24649c89c76SBlaise Bourdin {
24749c89c76SBlaise Bourdin PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data;
24849c89c76SBlaise Bourdin MPI_Comm comm;
2490a5cf5d8SBlaise Bourdin PetscExodusIIInt exoid = -1;
25049c89c76SBlaise Bourdin
25149c89c76SBlaise Bourdin PetscFunctionBegin;
25249c89c76SBlaise Bourdin PetscCall(PetscObjectGetComm((PetscObject)viewer, &comm));
25349c89c76SBlaise Bourdin PetscCheck(exo->numNodalVariables == -1, comm, PETSC_ERR_SUP, "The number of nodal variables has already been set to %d and cannot be overwritten", exo->numNodalVariables);
25449c89c76SBlaise Bourdin PetscCheck((exo->btype != FILE_MODE_READ) && (exo->btype != FILE_MODE_UNDEFINED), comm, PETSC_ERR_FILE_WRITE, "Cannot set the number of variables because the file is not writable");
25549c89c76SBlaise Bourdin
25649c89c76SBlaise Bourdin exo->numNodalVariables = num;
25749c89c76SBlaise Bourdin PetscCall(PetscMalloc1(num, &exo->nodalVariableNames));
258ac530a7eSPierre Jolivet for (int i = 0; i < num; i++) exo->nodalVariableNames[i] = NULL;
25949c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIIGetId(viewer, &exoid));
26049c89c76SBlaise Bourdin PetscCallExternal(ex_put_variable_param, exoid, EX_NODAL, num);
26149c89c76SBlaise Bourdin PetscFunctionReturn(PETSC_SUCCESS);
26249c89c76SBlaise Bourdin }
26349c89c76SBlaise Bourdin
26449c89c76SBlaise Bourdin /*@
265caff39ffSPierre Jolivet PetscViewerExodusIIGetZonalVariable - Gets the number of zonal variables in an ExodusII file
26649c89c76SBlaise Bourdin
26749c89c76SBlaise Bourdin Collective
26849c89c76SBlaise Bourdin
26949c89c76SBlaise Bourdin Input Parameters:
27049c89c76SBlaise Bourdin . viewer - a `PetscViewer` of type `PETSCVIEWEREXODUSII`
27149c89c76SBlaise Bourdin
272ce78bad3SBarry Smith Output Parameter:
273caff39ffSPierre Jolivet . num - the number variables in the ExodusII file
27449c89c76SBlaise Bourdin
27549c89c76SBlaise Bourdin Level: intermediate
27649c89c76SBlaise Bourdin
27749c89c76SBlaise Bourdin Notes:
278caff39ffSPierre Jolivet The number of variables in the ExodusII file is cached in the viewer
27949c89c76SBlaise Bourdin
28049c89c76SBlaise Bourdin .seealso: `PETSCVIEWEREXODUSII`, `PetscViewer`, `PetscViewerCreate()`, `PetscViewerDestroy()`, `PetscViewerExodusIIOpen()`, `PetscViewerSetType()`, `PetscViewerType`, `PetscViewerExodusIIsetZonalVariable()`
28149c89c76SBlaise Bourdin @*/
PetscViewerExodusIIGetZonalVariable(PetscViewer viewer,PetscExodusIIInt * num)2820a5cf5d8SBlaise Bourdin PetscErrorCode PetscViewerExodusIIGetZonalVariable(PetscViewer viewer, PetscExodusIIInt *num)
28349c89c76SBlaise Bourdin {
28449c89c76SBlaise Bourdin PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data;
28549c89c76SBlaise Bourdin MPI_Comm comm;
2860a5cf5d8SBlaise Bourdin PetscExodusIIInt exoid = -1;
28749c89c76SBlaise Bourdin
28849c89c76SBlaise Bourdin PetscFunctionBegin;
28949c89c76SBlaise Bourdin if (exo->numZonalVariables > -1) {
29049c89c76SBlaise Bourdin *num = exo->numZonalVariables;
29149c89c76SBlaise Bourdin } else {
29249c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIIGetId(viewer, &exoid));
29349c89c76SBlaise Bourdin PetscCall(PetscObjectGetComm((PetscObject)viewer, &comm));
29449c89c76SBlaise Bourdin PetscCheck(exoid > 0, comm, PETSC_ERR_FILE_OPEN, "Exodus file is not open");
29549c89c76SBlaise Bourdin PetscCallExternal(ex_get_variable_param, exoid, EX_ELEM_BLOCK, num);
29649c89c76SBlaise Bourdin exo->numZonalVariables = *num;
29749c89c76SBlaise Bourdin PetscCall(PetscMalloc1(*num, &exo->zonalVariableNames));
298ac530a7eSPierre Jolivet for (int i = 0; i < *num; i++) exo->zonalVariableNames[i] = NULL;
29949c89c76SBlaise Bourdin }
30049c89c76SBlaise Bourdin PetscFunctionReturn(PETSC_SUCCESS);
30149c89c76SBlaise Bourdin }
30249c89c76SBlaise Bourdin
30349c89c76SBlaise Bourdin /*@
304caff39ffSPierre Jolivet PetscViewerExodusIIGetNodalVariable - Gets the number of nodal variables in an ExodusII file
30549c89c76SBlaise Bourdin
30649c89c76SBlaise Bourdin Collective
30749c89c76SBlaise Bourdin
30849c89c76SBlaise Bourdin Input Parameters:
30949c89c76SBlaise Bourdin . viewer - a `PetscViewer` of type `PETSCVIEWEREXODUSII`
31049c89c76SBlaise Bourdin
311ce78bad3SBarry Smith Output Parameter:
312caff39ffSPierre Jolivet . num - the number variables in the ExodusII file
31349c89c76SBlaise Bourdin
31449c89c76SBlaise Bourdin Level: intermediate
31549c89c76SBlaise Bourdin
31649c89c76SBlaise Bourdin Notes:
31749c89c76SBlaise Bourdin This function gets the number of nodal variables and saves it in the address of num.
31849c89c76SBlaise Bourdin
31949c89c76SBlaise Bourdin .seealso: `PETSCVIEWEREXODUSII`, `PetscViewer`, `PetscViewerCreate()`, `PetscViewerDestroy()`, `PetscViewerExodusIIOpen()`, `PetscViewerSetType()`, `PetscViewerType`, `PetscViewerExodusIISetNodalVariable()`
32049c89c76SBlaise Bourdin @*/
PetscViewerExodusIIGetNodalVariable(PetscViewer viewer,PetscExodusIIInt * num)3210a5cf5d8SBlaise Bourdin PetscErrorCode PetscViewerExodusIIGetNodalVariable(PetscViewer viewer, PetscExodusIIInt *num)
32249c89c76SBlaise Bourdin {
32349c89c76SBlaise Bourdin PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data;
32449c89c76SBlaise Bourdin MPI_Comm comm;
3250a5cf5d8SBlaise Bourdin PetscExodusIIInt exoid = -1;
32649c89c76SBlaise Bourdin
32749c89c76SBlaise Bourdin PetscFunctionBegin;
32849c89c76SBlaise Bourdin if (exo->numNodalVariables > -1) {
32949c89c76SBlaise Bourdin *num = exo->numNodalVariables;
33049c89c76SBlaise Bourdin } else {
33149c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIIGetId(viewer, &exoid));
33249c89c76SBlaise Bourdin PetscCall(PetscObjectGetComm((PetscObject)viewer, &comm));
33349c89c76SBlaise Bourdin PetscCheck(exoid > 0, comm, PETSC_ERR_FILE_OPEN, "Exodus file is not open");
33449c89c76SBlaise Bourdin PetscCallExternal(ex_get_variable_param, exoid, EX_NODAL, num);
33549c89c76SBlaise Bourdin exo->numNodalVariables = *num;
33649c89c76SBlaise Bourdin PetscCall(PetscMalloc1(*num, &exo->nodalVariableNames));
337ac530a7eSPierre Jolivet for (int i = 0; i < *num; i++) exo->nodalVariableNames[i] = NULL;
33849c89c76SBlaise Bourdin }
33949c89c76SBlaise Bourdin PetscFunctionReturn(PETSC_SUCCESS);
34049c89c76SBlaise Bourdin }
34149c89c76SBlaise Bourdin
34249c89c76SBlaise Bourdin /*@
34349c89c76SBlaise Bourdin PetscViewerExodusIISetZonalVariableName - Sets the name of a zonal variable.
34449c89c76SBlaise Bourdin
34549c89c76SBlaise Bourdin Collective;
34649c89c76SBlaise Bourdin
34749c89c76SBlaise Bourdin Input Parameters:
34849c89c76SBlaise Bourdin + viewer - a `PetscViewer` of type `PETSCVIEWEREXODUSII`
34949c89c76SBlaise Bourdin . idx - the index for which you want to save the name
35049c89c76SBlaise Bourdin - name - string containing the name characters
35149c89c76SBlaise Bourdin
35249c89c76SBlaise Bourdin Level: intermediate
35349c89c76SBlaise Bourdin
35449c89c76SBlaise Bourdin .seealso: `PETSCVIEWEREXODUSII`, `PetscViewer`, `PetscViewerCreate()`, `PetscViewerDestroy()`, `PetscViewerExodusIIOpen()`, `PetscViewerSetType()`, `PetscViewerType`, `PetscViewerExodusIIGetZonalVariableName()`
35549c89c76SBlaise Bourdin @*/
PetscViewerExodusIISetZonalVariableName(PetscViewer viewer,PetscExodusIIInt idx,const char name[])3560a5cf5d8SBlaise Bourdin PetscErrorCode PetscViewerExodusIISetZonalVariableName(PetscViewer viewer, PetscExodusIIInt idx, const char name[])
35749c89c76SBlaise Bourdin {
35849c89c76SBlaise Bourdin PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data;
35949c89c76SBlaise Bourdin
36049c89c76SBlaise Bourdin PetscFunctionBegin;
36149c89c76SBlaise Bourdin PetscCheck((idx >= 0) && (idx < exo->numZonalVariables), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Variable index out of range. Was PetscViewerExodusIISetZonalVariable called?");
36249c89c76SBlaise Bourdin PetscCall(PetscStrallocpy(name, (char **)&exo->zonalVariableNames[idx]));
36349c89c76SBlaise Bourdin PetscCallExternal(ex_put_variable_name, exo->exoid, EX_ELEM_BLOCK, idx + 1, name);
36449c89c76SBlaise Bourdin PetscFunctionReturn(PETSC_SUCCESS);
36549c89c76SBlaise Bourdin }
36649c89c76SBlaise Bourdin
36749c89c76SBlaise Bourdin /*@
36849c89c76SBlaise Bourdin PetscViewerExodusIISetNodalVariableName - Sets the name of a nodal variable.
36949c89c76SBlaise Bourdin
37049c89c76SBlaise Bourdin Collective;
37149c89c76SBlaise Bourdin
37249c89c76SBlaise Bourdin Input Parameters:
37349c89c76SBlaise Bourdin + viewer - a `PetscViewer` of type `PETSCVIEWEREXODUSII`
37449c89c76SBlaise Bourdin . idx - the index for which you want to save the name
37549c89c76SBlaise Bourdin - name - string containing the name characters
37649c89c76SBlaise Bourdin
37749c89c76SBlaise Bourdin Level: intermediate
37849c89c76SBlaise Bourdin
37949c89c76SBlaise Bourdin .seealso: `PETSCVIEWEREXODUSII`, `PetscViewer`, `PetscViewerCreate()`, `PetscViewerDestroy()`, `PetscViewerExodusIIOpen()`, `PetscViewerSetType()`, `PetscViewerType`, `PetscViewerExodusIIGetNodalVariableName()`
38049c89c76SBlaise Bourdin @*/
PetscViewerExodusIISetNodalVariableName(PetscViewer viewer,PetscExodusIIInt idx,const char name[])3810a5cf5d8SBlaise Bourdin PetscErrorCode PetscViewerExodusIISetNodalVariableName(PetscViewer viewer, PetscExodusIIInt idx, const char name[])
38249c89c76SBlaise Bourdin {
38349c89c76SBlaise Bourdin PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data;
38449c89c76SBlaise Bourdin
38549c89c76SBlaise Bourdin PetscFunctionBegin;
38649c89c76SBlaise Bourdin PetscCheck((idx >= 0) && (idx < exo->numNodalVariables), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Variable index out of range. Was PetscViewerExodusIISetNodalVariable called?");
38749c89c76SBlaise Bourdin PetscCall(PetscStrallocpy(name, (char **)&exo->nodalVariableNames[idx]));
38849c89c76SBlaise Bourdin PetscCallExternal(ex_put_variable_name, exo->exoid, EX_NODAL, idx + 1, name);
38949c89c76SBlaise Bourdin PetscFunctionReturn(PETSC_SUCCESS);
39049c89c76SBlaise Bourdin }
39149c89c76SBlaise Bourdin
39249c89c76SBlaise Bourdin /*@
39349c89c76SBlaise Bourdin PetscViewerExodusIIGetZonalVariableName - Gets the name of a zonal variable.
39449c89c76SBlaise Bourdin
39549c89c76SBlaise Bourdin Collective;
39649c89c76SBlaise Bourdin
39749c89c76SBlaise Bourdin Input Parameters:
39849c89c76SBlaise Bourdin + viewer - a `PetscViewer` of type `PETSCVIEWEREXODUSII`
39949c89c76SBlaise Bourdin - idx - the index for which you want to get the name
40049c89c76SBlaise Bourdin
40149c89c76SBlaise Bourdin Output Parameter:
40249c89c76SBlaise Bourdin . name - pointer to the string containing the name characters
40349c89c76SBlaise Bourdin
40449c89c76SBlaise Bourdin Level: intermediate
40549c89c76SBlaise Bourdin
40649c89c76SBlaise Bourdin .seealso: `PETSCVIEWEREXODUSII`, `PetscViewer`, `PetscViewerCreate()`, `PetscViewerDestroy()`, `PetscViewerExodusIIOpen()`, `PetscViewerSetType()`, `PetscViewerType`, `PetscViewerExodusIISetZonalVariableName()`
40749c89c76SBlaise Bourdin @*/
PetscViewerExodusIIGetZonalVariableName(PetscViewer viewer,PetscExodusIIInt idx,const char * name[])4080a5cf5d8SBlaise Bourdin PetscErrorCode PetscViewerExodusIIGetZonalVariableName(PetscViewer viewer, PetscExodusIIInt idx, const char *name[])
40949c89c76SBlaise Bourdin {
41049c89c76SBlaise Bourdin PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data;
4110a5cf5d8SBlaise Bourdin PetscExodusIIInt exoid = -1;
41249c89c76SBlaise Bourdin char tmpName[MAX_NAME_LENGTH + 1];
41349c89c76SBlaise Bourdin
41449c89c76SBlaise Bourdin PetscFunctionBegin;
41549c89c76SBlaise Bourdin PetscCheck(idx >= 0 && idx < exo->numZonalVariables, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Variable index out of range. Was PetscViewerExodusIISetZonalVariable called?");
41649c89c76SBlaise Bourdin if (!exo->zonalVariableNames[idx]) {
41749c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIIGetId(viewer, &exoid));
41849c89c76SBlaise Bourdin PetscCallExternal(ex_get_variable_name, exoid, EX_ELEM_BLOCK, idx + 1, tmpName);
41949c89c76SBlaise Bourdin PetscCall(PetscStrallocpy(tmpName, (char **)&exo->zonalVariableNames[idx]));
42049c89c76SBlaise Bourdin }
42149c89c76SBlaise Bourdin *name = exo->zonalVariableNames[idx];
42249c89c76SBlaise Bourdin PetscFunctionReturn(PETSC_SUCCESS);
42349c89c76SBlaise Bourdin }
42449c89c76SBlaise Bourdin
42549c89c76SBlaise Bourdin /*@
42649c89c76SBlaise Bourdin PetscViewerExodusIIGetNodalVariableName - Gets the name of a nodal variable.
42749c89c76SBlaise Bourdin
42849c89c76SBlaise Bourdin Collective;
42949c89c76SBlaise Bourdin
43049c89c76SBlaise Bourdin Input Parameters:
43149c89c76SBlaise Bourdin + viewer - a `PetscViewer` of type `PETSCVIEWEREXODUSII`
43249c89c76SBlaise Bourdin - idx - the index for which you want to save the name
43349c89c76SBlaise Bourdin
43449c89c76SBlaise Bourdin Output Parameter:
43549c89c76SBlaise Bourdin . name - string array containing name characters
43649c89c76SBlaise Bourdin
43749c89c76SBlaise Bourdin Level: intermediate
43849c89c76SBlaise Bourdin
43949c89c76SBlaise Bourdin .seealso: `PETSCVIEWEREXODUSII`, `PetscViewer`, `PetscViewerCreate()`, `PetscViewerDestroy()`, `PetscViewerExodusIIOpen()`, `PetscViewerSetType()`, `PetscViewerType`, `PetscViewerExodusIISetNodalVariableName()`
44049c89c76SBlaise Bourdin @*/
PetscViewerExodusIIGetNodalVariableName(PetscViewer viewer,PetscExodusIIInt idx,const char * name[])4410a5cf5d8SBlaise Bourdin PetscErrorCode PetscViewerExodusIIGetNodalVariableName(PetscViewer viewer, PetscExodusIIInt idx, const char *name[])
44249c89c76SBlaise Bourdin {
44349c89c76SBlaise Bourdin PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data;
4440a5cf5d8SBlaise Bourdin PetscExodusIIInt exoid = -1;
44549c89c76SBlaise Bourdin char tmpName[MAX_NAME_LENGTH + 1];
44649c89c76SBlaise Bourdin
44749c89c76SBlaise Bourdin PetscFunctionBegin;
44849c89c76SBlaise Bourdin PetscCheck((idx >= 0) && (idx < exo->numNodalVariables), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Variable index out of range. Was PetscViewerExodusIISetNodalVariable called?");
44949c89c76SBlaise Bourdin if (!exo->nodalVariableNames[idx]) {
45049c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIIGetId(viewer, &exoid));
45149c89c76SBlaise Bourdin PetscCallExternal(ex_get_variable_name, exoid, EX_NODAL, idx + 1, tmpName);
45249c89c76SBlaise Bourdin PetscCall(PetscStrallocpy(tmpName, (char **)&exo->nodalVariableNames[idx]));
45349c89c76SBlaise Bourdin }
45449c89c76SBlaise Bourdin *name = exo->nodalVariableNames[idx];
45549c89c76SBlaise Bourdin PetscFunctionReturn(PETSC_SUCCESS);
45649c89c76SBlaise Bourdin }
45749c89c76SBlaise Bourdin
45849c89c76SBlaise Bourdin /*@C
45949c89c76SBlaise Bourdin PetscViewerExodusIISetZonalVariableNames - Sets the names of all nodal variables
46049c89c76SBlaise Bourdin
46149c89c76SBlaise Bourdin Collective; No Fortran Support
46249c89c76SBlaise Bourdin
46349c89c76SBlaise Bourdin Input Parameters:
46449c89c76SBlaise Bourdin + viewer - a `PetscViewer` of type `PETSCVIEWEREXODUSII`
465ce78bad3SBarry Smith - names - an array of string names to be set, the strings are copied into the `PetscViewer`
46649c89c76SBlaise Bourdin
46749c89c76SBlaise Bourdin Level: intermediate
46849c89c76SBlaise Bourdin
46949c89c76SBlaise Bourdin Notes:
47049c89c76SBlaise Bourdin This function allows users to set multiple zonal variable names at a time.
47149c89c76SBlaise Bourdin
47249c89c76SBlaise Bourdin .seealso: `PETSCVIEWEREXODUSII`, `PetscViewer`, `PetscViewerCreate()`, `PetscViewerDestroy()`, `PetscViewerExodusIIOpen()`, `PetscViewerSetType()`, `PetscViewerType`, `PetscViewerExodusIIGetZonalVariableNames()`
47349c89c76SBlaise Bourdin @*/
PetscViewerExodusIISetZonalVariableNames(PetscViewer viewer,const char * const names[])474ce78bad3SBarry Smith PetscErrorCode PetscViewerExodusIISetZonalVariableNames(PetscViewer viewer, const char *const names[])
47549c89c76SBlaise Bourdin {
4760a5cf5d8SBlaise Bourdin PetscExodusIIInt numNames;
4770a5cf5d8SBlaise Bourdin PetscExodusIIInt exoid = -1;
47849c89c76SBlaise Bourdin PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data;
47949c89c76SBlaise Bourdin
48049c89c76SBlaise Bourdin PetscFunctionBegin;
48149c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIIGetZonalVariable(viewer, &numNames));
48249c89c76SBlaise Bourdin PetscCheck(numNames >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of zonal variables not set. Was PetscViewerExodusIISetZonalVariable called?");
48349c89c76SBlaise Bourdin
48449c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIIGetId(viewer, &exoid));
48549c89c76SBlaise Bourdin for (int i = 0; i < numNames; i++) {
48649c89c76SBlaise Bourdin PetscCall(PetscStrallocpy(names[i], &exo->zonalVariableNames[i]));
48749c89c76SBlaise Bourdin PetscCallExternal(ex_put_variable_name, exoid, EX_ELEM_BLOCK, i + 1, names[i]);
48849c89c76SBlaise Bourdin }
48949c89c76SBlaise Bourdin PetscFunctionReturn(PETSC_SUCCESS);
49049c89c76SBlaise Bourdin }
49149c89c76SBlaise Bourdin
49249c89c76SBlaise Bourdin /*@C
49349c89c76SBlaise Bourdin PetscViewerExodusIISetNodalVariableNames - Sets the names of all nodal variables.
49449c89c76SBlaise Bourdin
49549c89c76SBlaise Bourdin Collective; No Fortran Support
49649c89c76SBlaise Bourdin
49749c89c76SBlaise Bourdin Input Parameters:
49849c89c76SBlaise Bourdin + viewer - a `PetscViewer` of type `PETSCVIEWEREXODUSII`
499ce78bad3SBarry Smith - names - an array of string names to be set, the strings are copied into the `PetscViewer`
50049c89c76SBlaise Bourdin
50149c89c76SBlaise Bourdin Level: intermediate
50249c89c76SBlaise Bourdin
50349c89c76SBlaise Bourdin Notes:
50449c89c76SBlaise Bourdin This function allows users to set multiple nodal variable names at a time.
50549c89c76SBlaise Bourdin
50649c89c76SBlaise Bourdin .seealso: `PETSCVIEWEREXODUSII`, `PetscViewer`, `PetscViewerCreate()`, `PetscViewerDestroy()`, `PetscViewerExodusIIOpen()`, `PetscViewerSetType()`, `PetscViewerType`, `PetscViewerExodusIIGetNodalVariableNames()`
50749c89c76SBlaise Bourdin @*/
PetscViewerExodusIISetNodalVariableNames(PetscViewer viewer,const char * const names[])508ce78bad3SBarry Smith PetscErrorCode PetscViewerExodusIISetNodalVariableNames(PetscViewer viewer, const char *const names[])
50949c89c76SBlaise Bourdin {
5100a5cf5d8SBlaise Bourdin PetscExodusIIInt numNames;
5110a5cf5d8SBlaise Bourdin PetscExodusIIInt exoid = -1;
51249c89c76SBlaise Bourdin PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data;
51349c89c76SBlaise Bourdin
51449c89c76SBlaise Bourdin PetscFunctionBegin;
51549c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIIGetNodalVariable(viewer, &numNames));
51649c89c76SBlaise Bourdin PetscCheck(numNames >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of nodal variables not set. Was PetscViewerExodusIISetNodalVariable called?");
51749c89c76SBlaise Bourdin
51849c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIIGetId(viewer, &exoid));
51949c89c76SBlaise Bourdin for (int i = 0; i < numNames; i++) {
52049c89c76SBlaise Bourdin PetscCall(PetscStrallocpy(names[i], &exo->nodalVariableNames[i]));
52149c89c76SBlaise Bourdin PetscCallExternal(ex_put_variable_name, exoid, EX_NODAL, i + 1, names[i]);
52249c89c76SBlaise Bourdin }
52349c89c76SBlaise Bourdin PetscFunctionReturn(PETSC_SUCCESS);
52449c89c76SBlaise Bourdin }
52549c89c76SBlaise Bourdin
52649c89c76SBlaise Bourdin /*@C
52749c89c76SBlaise Bourdin PetscViewerExodusIIGetZonalVariableNames - Gets the names of all zonal variables.
52849c89c76SBlaise Bourdin
52949c89c76SBlaise Bourdin Collective; No Fortran Support
53049c89c76SBlaise Bourdin
53149c89c76SBlaise Bourdin Input Parameters:
53249c89c76SBlaise Bourdin + viewer - a `PetscViewer` of type `PETSCVIEWEREXODUSII`
53349c89c76SBlaise Bourdin - numVars - the number of zonal variable names to retrieve
53449c89c76SBlaise Bourdin
535ce78bad3SBarry Smith Output Parameter:
536ce78bad3SBarry Smith . varNames - returns an array of char pointers where the zonal variable names are
53749c89c76SBlaise Bourdin
53849c89c76SBlaise Bourdin Level: intermediate
53949c89c76SBlaise Bourdin
54049c89c76SBlaise Bourdin Notes:
54149c89c76SBlaise Bourdin This function returns a borrowed pointer which should not be freed.
54249c89c76SBlaise Bourdin
54349c89c76SBlaise Bourdin .seealso: `PETSCVIEWEREXODUSII`, `PetscViewer`, `PetscViewerCreate()`, `PetscViewerDestroy()`, `PetscViewerExodusIIOpen()`, `PetscViewerSetType()`, `PetscViewerType`, `PetscViewerExodusIISetZonalVariableNames()`
54449c89c76SBlaise Bourdin @*/
PetscViewerExodusIIGetZonalVariableNames(PetscViewer viewer,PetscExodusIIInt * numVars,const char * const * varNames[])545ce78bad3SBarry Smith PetscErrorCode PetscViewerExodusIIGetZonalVariableNames(PetscViewer viewer, PetscExodusIIInt *numVars, const char *const *varNames[])
54649c89c76SBlaise Bourdin {
54749c89c76SBlaise Bourdin PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data;
5480a5cf5d8SBlaise Bourdin PetscExodusIIInt idx;
54949c89c76SBlaise Bourdin char tmpName[MAX_NAME_LENGTH + 1];
5500a5cf5d8SBlaise Bourdin PetscExodusIIInt exoid = -1;
55149c89c76SBlaise Bourdin
55249c89c76SBlaise Bourdin PetscFunctionBegin;
55349c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIIGetZonalVariable(viewer, numVars));
55449c89c76SBlaise Bourdin /*
55549c89c76SBlaise Bourdin Cache variable names if necessary
55649c89c76SBlaise Bourdin */
55749c89c76SBlaise Bourdin for (idx = 0; idx < *numVars; idx++) {
55849c89c76SBlaise Bourdin if (!exo->zonalVariableNames[idx]) {
55949c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIIGetId(viewer, &exoid));
56049c89c76SBlaise Bourdin PetscCallExternal(ex_get_variable_name, exoid, EX_ELEM_BLOCK, idx + 1, tmpName);
56149c89c76SBlaise Bourdin PetscCall(PetscStrallocpy(tmpName, (char **)&exo->zonalVariableNames[idx]));
56249c89c76SBlaise Bourdin }
56349c89c76SBlaise Bourdin }
564ce78bad3SBarry Smith *varNames = (const char *const *)exo->zonalVariableNames;
56549c89c76SBlaise Bourdin PetscFunctionReturn(PETSC_SUCCESS);
56649c89c76SBlaise Bourdin }
56749c89c76SBlaise Bourdin
56849c89c76SBlaise Bourdin /*@C
56949c89c76SBlaise Bourdin PetscViewerExodusIIGetNodalVariableNames - Gets the names of all nodal variables.
57049c89c76SBlaise Bourdin
57149c89c76SBlaise Bourdin Collective; No Fortran Support
57249c89c76SBlaise Bourdin
57349c89c76SBlaise Bourdin Input Parameters:
57449c89c76SBlaise Bourdin + viewer - a `PetscViewer` of type `PETSCVIEWEREXODUSII`
57549c89c76SBlaise Bourdin - numVars - the number of nodal variable names to retrieve
57649c89c76SBlaise Bourdin
577ce78bad3SBarry Smith Output Parameter:
578ce78bad3SBarry Smith . varNames - returns an array of char pointers where the nodal variable names are
57949c89c76SBlaise Bourdin
58049c89c76SBlaise Bourdin Level: intermediate
58149c89c76SBlaise Bourdin
58249c89c76SBlaise Bourdin Notes:
58349c89c76SBlaise Bourdin This function returns a borrowed pointer which should not be freed.
58449c89c76SBlaise Bourdin
58549c89c76SBlaise Bourdin .seealso: `PETSCVIEWEREXODUSII`, `PetscViewer`, `PetscViewerCreate()`, `PetscViewerDestroy()`, `PetscViewerExodusIIOpen()`, `PetscViewerSetType()`, `PetscViewerType`, `PetscViewerExodusIISetNodalVariableNames()`
58649c89c76SBlaise Bourdin @*/
PetscViewerExodusIIGetNodalVariableNames(PetscViewer viewer,PetscExodusIIInt * numVars,const char * const * varNames[])587ce78bad3SBarry Smith PetscErrorCode PetscViewerExodusIIGetNodalVariableNames(PetscViewer viewer, PetscExodusIIInt *numVars, const char *const *varNames[])
58849c89c76SBlaise Bourdin {
58949c89c76SBlaise Bourdin PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data;
5900a5cf5d8SBlaise Bourdin PetscExodusIIInt idx;
59149c89c76SBlaise Bourdin char tmpName[MAX_NAME_LENGTH + 1];
5920a5cf5d8SBlaise Bourdin PetscExodusIIInt exoid = -1;
59349c89c76SBlaise Bourdin
59449c89c76SBlaise Bourdin PetscFunctionBegin;
59549c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIIGetNodalVariable(viewer, numVars));
59649c89c76SBlaise Bourdin /*
59749c89c76SBlaise Bourdin Cache variable names if necessary
59849c89c76SBlaise Bourdin */
59949c89c76SBlaise Bourdin for (idx = 0; idx < *numVars; idx++) {
60049c89c76SBlaise Bourdin if (!exo->nodalVariableNames[idx]) {
60149c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIIGetId(viewer, &exoid));
60249c89c76SBlaise Bourdin PetscCallExternal(ex_get_variable_name, exoid, EX_NODAL, idx + 1, tmpName);
60349c89c76SBlaise Bourdin PetscCall(PetscStrallocpy(tmpName, (char **)&exo->nodalVariableNames[idx]));
60449c89c76SBlaise Bourdin }
60549c89c76SBlaise Bourdin }
606ce78bad3SBarry Smith *varNames = (const char *const *)exo->nodalVariableNames;
60749c89c76SBlaise Bourdin PetscFunctionReturn(PETSC_SUCCESS);
60849c89c76SBlaise Bourdin }
60949c89c76SBlaise Bourdin
61049c89c76SBlaise Bourdin /*MC
61149c89c76SBlaise Bourdin PETSCVIEWEREXODUSII - A viewer that writes to an Exodus II file
61249c89c76SBlaise Bourdin
61349c89c76SBlaise Bourdin Level: beginner
61449c89c76SBlaise Bourdin
61549c89c76SBlaise Bourdin .seealso: `PetscViewerExodusIIOpen()`, `PetscViewerCreate()`, `PETSCVIEWERBINARY`, `PETSCVIEWERHDF5`, `DMView()`,
61649c89c76SBlaise Bourdin `PetscViewerFileSetName()`, `PetscViewerFileSetMode()`, `PetscViewerFormat`, `PetscViewerType`, `PetscViewerSetType()`
61749c89c76SBlaise Bourdin M*/
PetscViewerCreate_ExodusII(PetscViewer v)61849c89c76SBlaise Bourdin PETSC_EXTERN PetscErrorCode PetscViewerCreate_ExodusII(PetscViewer v)
61949c89c76SBlaise Bourdin {
62049c89c76SBlaise Bourdin PetscViewer_ExodusII *exo;
62149c89c76SBlaise Bourdin
62249c89c76SBlaise Bourdin PetscFunctionBegin;
62349c89c76SBlaise Bourdin PetscCall(PetscNew(&exo));
62449c89c76SBlaise Bourdin
62549c89c76SBlaise Bourdin v->data = (void *)exo;
62649c89c76SBlaise Bourdin v->ops->destroy = PetscViewerDestroy_ExodusII;
62749c89c76SBlaise Bourdin v->ops->setfromoptions = PetscViewerSetFromOptions_ExodusII;
62849c89c76SBlaise Bourdin v->ops->view = PetscViewerView_ExodusII;
62949c89c76SBlaise Bourdin v->ops->flush = PetscViewerFlush_ExodusII;
63049c89c76SBlaise Bourdin exo->btype = FILE_MODE_UNDEFINED;
63149c89c76SBlaise Bourdin exo->filename = 0;
63249c89c76SBlaise Bourdin exo->exoid = -1;
63349c89c76SBlaise Bourdin exo->numNodalVariables = -1;
63449c89c76SBlaise Bourdin exo->numZonalVariables = -1;
63549c89c76SBlaise Bourdin exo->nodalVariableNames = NULL;
63649c89c76SBlaise Bourdin exo->zonalVariableNames = NULL;
63749c89c76SBlaise Bourdin
63849c89c76SBlaise Bourdin PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetName_C", PetscViewerFileSetName_ExodusII));
63949c89c76SBlaise Bourdin PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetName_C", PetscViewerFileGetName_ExodusII));
64049c89c76SBlaise Bourdin PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetMode_C", PetscViewerFileSetMode_ExodusII));
64149c89c76SBlaise Bourdin PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetMode_C", PetscViewerFileGetMode_ExodusII));
64249c89c76SBlaise Bourdin PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerGetId_C", PetscViewerExodusIIGetId_ExodusII));
64349c89c76SBlaise Bourdin PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerSetOrder_C", PetscViewerExodusIISetOrder_ExodusII));
64449c89c76SBlaise Bourdin PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerGetOrder_C", PetscViewerExodusIIGetOrder_ExodusII));
64549c89c76SBlaise Bourdin PetscFunctionReturn(PETSC_SUCCESS);
64649c89c76SBlaise Bourdin }
64749c89c76SBlaise Bourdin
64849c89c76SBlaise Bourdin /*@
649caff39ffSPierre Jolivet PetscViewerExodusIIGetNodalVariableIndex - return the location of a nodal variable in an ExodusII file given its name
65049c89c76SBlaise Bourdin
65149c89c76SBlaise Bourdin Collective
65249c89c76SBlaise Bourdin
65349c89c76SBlaise Bourdin Input Parameters:
65449c89c76SBlaise Bourdin + viewer - a `PetscViewer` of type `PETSCVIEWEREXODUSII`
65549c89c76SBlaise Bourdin - name - the name of the result
65649c89c76SBlaise Bourdin
65749c89c76SBlaise Bourdin Output Parameter:
65849c89c76SBlaise Bourdin . varIndex - the location of the variable in the exodus file or -1 if the variable is not found
65949c89c76SBlaise Bourdin
66049c89c76SBlaise Bourdin Level: beginner
66149c89c76SBlaise Bourdin
66249c89c76SBlaise Bourdin Notes:
66349c89c76SBlaise Bourdin The exodus variable index is obtained by comparing the name argument to the
66449c89c76SBlaise Bourdin names of zonal variables declared in the exodus file. For instance if name is "V"
66549c89c76SBlaise Bourdin the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
66649c89c76SBlaise Bourdin amongst all variables of type obj_type.
66749c89c76SBlaise Bourdin
66849c89c76SBlaise Bourdin .seealso: `PetscViewerExodusIISetNodalVariable()`, `PetscViewerExodusIIGetNodalVariable()`, `PetscViewerExodusIISetNodalVariableName()`, `PetscViewerExodusIISetNodalVariableNames()`, `PetscViewerExodusIIGetNodalVariableName()`, `PetscViewerExodusIIGetNodalVariableNames()`
66949c89c76SBlaise Bourdin @*/
PetscViewerExodusIIGetNodalVariableIndex(PetscViewer viewer,const char name[],PetscExodusIIInt * varIndex)6700a5cf5d8SBlaise Bourdin PetscErrorCode PetscViewerExodusIIGetNodalVariableIndex(PetscViewer viewer, const char name[], PetscExodusIIInt *varIndex)
67149c89c76SBlaise Bourdin {
6720a5cf5d8SBlaise Bourdin PetscExodusIIInt num_vars = 0, i, j;
67349c89c76SBlaise Bourdin char ext_name[MAX_STR_LENGTH + 1];
674ce78bad3SBarry Smith const char *const *var_names;
67549c89c76SBlaise Bourdin const int num_suffix = 5;
67649c89c76SBlaise Bourdin char *suffix[5];
67749c89c76SBlaise Bourdin PetscBool flg;
67849c89c76SBlaise Bourdin
67949c89c76SBlaise Bourdin PetscFunctionBegin;
68049c89c76SBlaise Bourdin suffix[0] = (char *)"";
68149c89c76SBlaise Bourdin suffix[1] = (char *)"_X";
68249c89c76SBlaise Bourdin suffix[2] = (char *)"_XX";
68349c89c76SBlaise Bourdin suffix[3] = (char *)"_1";
68449c89c76SBlaise Bourdin suffix[4] = (char *)"_11";
68549c89c76SBlaise Bourdin *varIndex = -1;
68649c89c76SBlaise Bourdin
68749c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIIGetNodalVariableNames(viewer, &num_vars, &var_names));
68849c89c76SBlaise Bourdin for (i = 0; i < num_vars; ++i) {
68949c89c76SBlaise Bourdin for (j = 0; j < num_suffix; ++j) {
69049c89c76SBlaise Bourdin PetscCall(PetscStrncpy(ext_name, name, MAX_STR_LENGTH));
69149c89c76SBlaise Bourdin PetscCall(PetscStrlcat(ext_name, suffix[j], MAX_STR_LENGTH));
69249c89c76SBlaise Bourdin PetscCall(PetscStrcasecmp(ext_name, var_names[i], &flg));
69349c89c76SBlaise Bourdin if (flg) *varIndex = i;
69449c89c76SBlaise Bourdin }
69549c89c76SBlaise Bourdin if (flg) break;
69649c89c76SBlaise Bourdin }
69749c89c76SBlaise Bourdin PetscFunctionReturn(PETSC_SUCCESS);
69849c89c76SBlaise Bourdin }
69949c89c76SBlaise Bourdin
70049c89c76SBlaise Bourdin /*@
701caff39ffSPierre Jolivet PetscViewerExodusIIGetZonalVariableIndex - return the location of a zonal variable in an ExodusII file given its name
70249c89c76SBlaise Bourdin
70349c89c76SBlaise Bourdin Collective
70449c89c76SBlaise Bourdin
70549c89c76SBlaise Bourdin Input Parameters:
70649c89c76SBlaise Bourdin + viewer - a `PetscViewer` of type `PETSCVIEWEREXODUSII`
70749c89c76SBlaise Bourdin - name - the name of the result
70849c89c76SBlaise Bourdin
70949c89c76SBlaise Bourdin Output Parameter:
71049c89c76SBlaise Bourdin . varIndex - the location of the variable in the exodus file or -1 if the variable is not found
71149c89c76SBlaise Bourdin
71249c89c76SBlaise Bourdin Level: beginner
71349c89c76SBlaise Bourdin
71449c89c76SBlaise Bourdin Notes:
71549c89c76SBlaise Bourdin The exodus variable index is obtained by comparing the name argument to the
71649c89c76SBlaise Bourdin names of zonal variables declared in the exodus file. For instance if name is "V"
71749c89c76SBlaise Bourdin the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
71849c89c76SBlaise Bourdin amongst all variables of type obj_type.
71949c89c76SBlaise Bourdin
72049c89c76SBlaise Bourdin .seealso: `PetscViewerExodusIISetNodalVariable()`, `PetscViewerExodusIIGetNodalVariable()`, `PetscViewerExodusIISetNodalVariableName()`, `PetscViewerExodusIISetNodalVariableNames()`, `PetscViewerExodusIIGetNodalVariableName()`, `PetscViewerExodusIIGetNodalVariableNames()`
72149c89c76SBlaise Bourdin @*/
PetscViewerExodusIIGetZonalVariableIndex(PetscViewer viewer,const char name[],int * varIndex)72249c89c76SBlaise Bourdin PetscErrorCode PetscViewerExodusIIGetZonalVariableIndex(PetscViewer viewer, const char name[], int *varIndex)
72349c89c76SBlaise Bourdin {
7240a5cf5d8SBlaise Bourdin PetscExodusIIInt num_vars = 0, i, j;
72549c89c76SBlaise Bourdin char ext_name[MAX_STR_LENGTH + 1];
726ce78bad3SBarry Smith const char *const *var_names;
72749c89c76SBlaise Bourdin const int num_suffix = 5;
72849c89c76SBlaise Bourdin char *suffix[5];
72949c89c76SBlaise Bourdin PetscBool flg;
73049c89c76SBlaise Bourdin
73149c89c76SBlaise Bourdin PetscFunctionBegin;
73249c89c76SBlaise Bourdin suffix[0] = (char *)"";
73349c89c76SBlaise Bourdin suffix[1] = (char *)"_X";
73449c89c76SBlaise Bourdin suffix[2] = (char *)"_XX";
73549c89c76SBlaise Bourdin suffix[3] = (char *)"_1";
73649c89c76SBlaise Bourdin suffix[4] = (char *)"_11";
73749c89c76SBlaise Bourdin *varIndex = -1;
73849c89c76SBlaise Bourdin
73949c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIIGetZonalVariableNames(viewer, &num_vars, &var_names));
74049c89c76SBlaise Bourdin for (i = 0; i < num_vars; ++i) {
74149c89c76SBlaise Bourdin for (j = 0; j < num_suffix; ++j) {
74249c89c76SBlaise Bourdin PetscCall(PetscStrncpy(ext_name, name, MAX_STR_LENGTH));
74349c89c76SBlaise Bourdin PetscCall(PetscStrlcat(ext_name, suffix[j], MAX_STR_LENGTH));
74449c89c76SBlaise Bourdin PetscCall(PetscStrcasecmp(ext_name, var_names[i], &flg));
74549c89c76SBlaise Bourdin if (flg) *varIndex = i;
74649c89c76SBlaise Bourdin }
74749c89c76SBlaise Bourdin if (flg) break;
74849c89c76SBlaise Bourdin }
74949c89c76SBlaise Bourdin PetscFunctionReturn(PETSC_SUCCESS);
75049c89c76SBlaise Bourdin }
75149c89c76SBlaise Bourdin
DMView_PlexExodusII(DM dm,PetscViewer viewer)75249c89c76SBlaise Bourdin PetscErrorCode DMView_PlexExodusII(DM dm, PetscViewer viewer)
75349c89c76SBlaise Bourdin {
75449c89c76SBlaise Bourdin enum ElemType {
75549c89c76SBlaise Bourdin SEGMENT,
75649c89c76SBlaise Bourdin TRI,
75749c89c76SBlaise Bourdin QUAD,
75849c89c76SBlaise Bourdin TET,
75949c89c76SBlaise Bourdin HEX
76049c89c76SBlaise Bourdin };
76149c89c76SBlaise Bourdin MPI_Comm comm;
7620a5cf5d8SBlaise Bourdin PetscInt degree; /* the order of the mesh */
76349c89c76SBlaise Bourdin /* Connectivity Variables */
76449c89c76SBlaise Bourdin PetscInt cellsNotInConnectivity;
76549c89c76SBlaise Bourdin /* Cell Sets */
76649c89c76SBlaise Bourdin DMLabel csLabel;
76749c89c76SBlaise Bourdin IS csIS;
76849c89c76SBlaise Bourdin const PetscInt *csIdx;
76949c89c76SBlaise Bourdin PetscInt num_cs, cs;
77049c89c76SBlaise Bourdin enum ElemType *type;
77149c89c76SBlaise Bourdin PetscBool hasLabel;
77249c89c76SBlaise Bourdin /* Coordinate Variables */
77349c89c76SBlaise Bourdin DM cdm;
77449c89c76SBlaise Bourdin PetscSection coordSection;
77549c89c76SBlaise Bourdin Vec coord;
77649c89c76SBlaise Bourdin PetscInt **nodes;
77749c89c76SBlaise Bourdin PetscInt depth, d, dim, skipCells = 0;
77849c89c76SBlaise Bourdin PetscInt pStart, pEnd, p, cStart, cEnd, numCells, vStart, vEnd, numVertices, eStart, eEnd, numEdges, fStart, fEnd, numFaces, numNodes;
77949c89c76SBlaise Bourdin PetscInt num_vs, num_fs;
78049c89c76SBlaise Bourdin PetscMPIInt rank, size;
78149c89c76SBlaise Bourdin const char *dmName;
78249c89c76SBlaise Bourdin PetscInt nodesLineP1[4] = {2, 0, 0, 0};
78349c89c76SBlaise Bourdin PetscInt nodesLineP2[4] = {2, 0, 0, 1};
78449c89c76SBlaise Bourdin PetscInt nodesTriP1[4] = {3, 0, 0, 0};
78549c89c76SBlaise Bourdin PetscInt nodesTriP2[4] = {3, 3, 0, 0};
78649c89c76SBlaise Bourdin PetscInt nodesQuadP1[4] = {4, 0, 0, 0};
78749c89c76SBlaise Bourdin PetscInt nodesQuadP2[4] = {4, 4, 0, 1};
78849c89c76SBlaise Bourdin PetscInt nodesTetP1[4] = {4, 0, 0, 0};
78949c89c76SBlaise Bourdin PetscInt nodesTetP2[4] = {4, 6, 0, 0};
79049c89c76SBlaise Bourdin PetscInt nodesHexP1[4] = {8, 0, 0, 0};
79149c89c76SBlaise Bourdin PetscInt nodesHexP2[4] = {8, 12, 6, 1};
7920a5cf5d8SBlaise Bourdin PetscExodusIIInt CPU_word_size, IO_word_size, EXO_mode;
7930a5cf5d8SBlaise Bourdin PetscExodusIIFloat EXO_version;
79449c89c76SBlaise Bourdin
79549c89c76SBlaise Bourdin PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data;
79649c89c76SBlaise Bourdin
79749c89c76SBlaise Bourdin PetscFunctionBegin;
79849c89c76SBlaise Bourdin PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
79949c89c76SBlaise Bourdin PetscCallMPI(MPI_Comm_rank(comm, &rank));
80049c89c76SBlaise Bourdin PetscCallMPI(MPI_Comm_size(comm, &size));
80149c89c76SBlaise Bourdin
80249c89c76SBlaise Bourdin /*
80349c89c76SBlaise Bourdin Creating coordSection is a collective operation so we do it somewhat out of sequence
80449c89c76SBlaise Bourdin */
80549c89c76SBlaise Bourdin PetscCall(PetscSectionCreate(comm, &coordSection));
80649c89c76SBlaise Bourdin PetscCall(DMGetCoordinatesLocalSetUp(dm));
80749c89c76SBlaise Bourdin /*
80849c89c76SBlaise Bourdin Check that all points are on rank 0 since we don't know how to save distributed DM in exodus format
80949c89c76SBlaise Bourdin */
81049c89c76SBlaise Bourdin PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
81149c89c76SBlaise Bourdin PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd));
81249c89c76SBlaise Bourdin PetscCall(DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd));
81349c89c76SBlaise Bourdin PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
81449c89c76SBlaise Bourdin numCells = cEnd - cStart;
81549c89c76SBlaise Bourdin numEdges = eEnd - eStart;
81649c89c76SBlaise Bourdin numVertices = vEnd - vStart;
817caff39ffSPierre Jolivet PetscCheck(!(rank && (numCells || numEdges || numVertices)), PETSC_COMM_SELF, PETSC_ERR_SUP, "Writing distributed DM in ExodusII format not supported");
81849c89c76SBlaise Bourdin if (rank == 0) {
81949c89c76SBlaise Bourdin switch (exo->btype) {
82049c89c76SBlaise Bourdin case FILE_MODE_READ:
82149c89c76SBlaise Bourdin case FILE_MODE_APPEND:
82249c89c76SBlaise Bourdin case FILE_MODE_UPDATE:
82349c89c76SBlaise Bourdin case FILE_MODE_APPEND_UPDATE:
824caff39ffSPierre Jolivet /* ExodusII does not allow writing geometry to an existing file */
82549c89c76SBlaise Bourdin SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "cannot add geometry to existing file %s", exo->filename);
82649c89c76SBlaise Bourdin case FILE_MODE_WRITE:
82749c89c76SBlaise Bourdin /* Create an empty file if one already exists*/
82849c89c76SBlaise Bourdin EXO_mode = EX_CLOBBER;
82949c89c76SBlaise Bourdin #if defined(PETSC_USE_64BIT_INDICES)
83049c89c76SBlaise Bourdin EXO_mode += EX_ALL_INT64_API;
83149c89c76SBlaise Bourdin #endif
83249c89c76SBlaise Bourdin CPU_word_size = sizeof(PetscReal);
83349c89c76SBlaise Bourdin IO_word_size = sizeof(PetscReal);
83449c89c76SBlaise Bourdin exo->exoid = ex_create(exo->filename, EXO_mode, &CPU_word_size, &IO_word_size);
83549c89c76SBlaise Bourdin PetscCheck(exo->exoid >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_create failed for %s", exo->filename);
83649c89c76SBlaise Bourdin
83749c89c76SBlaise Bourdin break;
83849c89c76SBlaise Bourdin default:
83949c89c76SBlaise Bourdin SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
84049c89c76SBlaise Bourdin }
84149c89c76SBlaise Bourdin
84249c89c76SBlaise Bourdin /* --- Get DM info --- */
84349c89c76SBlaise Bourdin PetscCall(PetscObjectGetName((PetscObject)dm, &dmName));
84449c89c76SBlaise Bourdin PetscCall(DMPlexGetDepth(dm, &depth));
84549c89c76SBlaise Bourdin PetscCall(DMGetDimension(dm, &dim));
84649c89c76SBlaise Bourdin PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
84749c89c76SBlaise Bourdin if (depth == 3) {
84849c89c76SBlaise Bourdin numFaces = fEnd - fStart;
84949c89c76SBlaise Bourdin } else {
85049c89c76SBlaise Bourdin numFaces = 0;
85149c89c76SBlaise Bourdin }
85249c89c76SBlaise Bourdin PetscCall(DMGetLabelSize(dm, "Cell Sets", &num_cs));
85349c89c76SBlaise Bourdin PetscCall(DMGetLabelSize(dm, "Vertex Sets", &num_vs));
85449c89c76SBlaise Bourdin PetscCall(DMGetLabelSize(dm, "Face Sets", &num_fs));
85549c89c76SBlaise Bourdin PetscCall(DMGetCoordinatesLocal(dm, &coord));
85649c89c76SBlaise Bourdin PetscCall(DMGetCoordinateDM(dm, &cdm));
85749c89c76SBlaise Bourdin if (num_cs > 0) {
85849c89c76SBlaise Bourdin PetscCall(DMGetLabel(dm, "Cell Sets", &csLabel));
85949c89c76SBlaise Bourdin PetscCall(DMLabelGetValueIS(csLabel, &csIS));
86049c89c76SBlaise Bourdin PetscCall(ISGetIndices(csIS, &csIdx));
86149c89c76SBlaise Bourdin }
86249c89c76SBlaise Bourdin PetscCall(PetscMalloc1(num_cs, &nodes));
86349c89c76SBlaise Bourdin /* Set element type for each block and compute total number of nodes */
86449c89c76SBlaise Bourdin PetscCall(PetscMalloc1(num_cs, &type));
86549c89c76SBlaise Bourdin numNodes = numVertices;
86649c89c76SBlaise Bourdin
86749c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIIGetOrder(viewer, °ree));
86849c89c76SBlaise Bourdin if (degree == 2) numNodes += numEdges;
86949c89c76SBlaise Bourdin cellsNotInConnectivity = numCells;
87049c89c76SBlaise Bourdin for (cs = 0; cs < num_cs; ++cs) {
87149c89c76SBlaise Bourdin IS stratumIS;
87249c89c76SBlaise Bourdin const PetscInt *cells;
87349c89c76SBlaise Bourdin PetscScalar *xyz = NULL;
87449c89c76SBlaise Bourdin PetscInt csSize, closureSize;
87549c89c76SBlaise Bourdin
87649c89c76SBlaise Bourdin PetscCall(DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS));
87749c89c76SBlaise Bourdin PetscCall(ISGetIndices(stratumIS, &cells));
87849c89c76SBlaise Bourdin PetscCall(ISGetSize(stratumIS, &csSize));
87949c89c76SBlaise Bourdin PetscCall(DMPlexVecGetClosure(cdm, NULL, coord, cells[0], &closureSize, &xyz));
88049c89c76SBlaise Bourdin switch (dim) {
88149c89c76SBlaise Bourdin case 1:
882966bd95aSPierre Jolivet PetscCheck(closureSize == 2 * dim, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of vertices %" PetscInt_FMT " in dimension %" PetscInt_FMT " has no ExodusII type", closureSize / dim, dim);
88349c89c76SBlaise Bourdin type[cs] = SEGMENT;
884fac68d15SPierre Jolivet break;
88549c89c76SBlaise Bourdin case 2:
88649c89c76SBlaise Bourdin if (closureSize == 3 * dim) {
88749c89c76SBlaise Bourdin type[cs] = TRI;
88849c89c76SBlaise Bourdin } else if (closureSize == 4 * dim) {
88949c89c76SBlaise Bourdin type[cs] = QUAD;
89049c89c76SBlaise Bourdin } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of vertices %" PetscInt_FMT " in dimension %" PetscInt_FMT " has no ExodusII type", closureSize / dim, dim);
89149c89c76SBlaise Bourdin break;
89249c89c76SBlaise Bourdin case 3:
89349c89c76SBlaise Bourdin if (closureSize == 4 * dim) {
89449c89c76SBlaise Bourdin type[cs] = TET;
89549c89c76SBlaise Bourdin } else if (closureSize == 8 * dim) {
89649c89c76SBlaise Bourdin type[cs] = HEX;
89749c89c76SBlaise Bourdin } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of vertices %" PetscInt_FMT " in dimension %" PetscInt_FMT " has no ExodusII type", closureSize / dim, dim);
89849c89c76SBlaise Bourdin break;
89949c89c76SBlaise Bourdin default:
90049c89c76SBlaise Bourdin SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %" PetscInt_FMT " not handled by ExodusII viewer", dim);
90149c89c76SBlaise Bourdin }
90249c89c76SBlaise Bourdin if ((degree == 2) && (type[cs] == SEGMENT)) numNodes += csSize;
90349c89c76SBlaise Bourdin if ((degree == 2) && (type[cs] == QUAD)) numNodes += csSize;
90449c89c76SBlaise Bourdin if ((degree == 2) && (type[cs] == HEX)) {
90549c89c76SBlaise Bourdin numNodes += csSize;
90649c89c76SBlaise Bourdin numNodes += numFaces;
90749c89c76SBlaise Bourdin }
90849c89c76SBlaise Bourdin PetscCall(DMPlexVecRestoreClosure(cdm, NULL, coord, cells[0], &closureSize, &xyz));
90949c89c76SBlaise Bourdin /* Set nodes and Element type */
91049c89c76SBlaise Bourdin if (type[cs] == SEGMENT) {
91149c89c76SBlaise Bourdin if (degree == 1) nodes[cs] = nodesLineP1;
91249c89c76SBlaise Bourdin else if (degree == 2) nodes[cs] = nodesLineP2;
91349c89c76SBlaise Bourdin } else if (type[cs] == TRI) {
91449c89c76SBlaise Bourdin if (degree == 1) nodes[cs] = nodesTriP1;
91549c89c76SBlaise Bourdin else if (degree == 2) nodes[cs] = nodesTriP2;
91649c89c76SBlaise Bourdin } else if (type[cs] == QUAD) {
91749c89c76SBlaise Bourdin if (degree == 1) nodes[cs] = nodesQuadP1;
91849c89c76SBlaise Bourdin else if (degree == 2) nodes[cs] = nodesQuadP2;
91949c89c76SBlaise Bourdin } else if (type[cs] == TET) {
92049c89c76SBlaise Bourdin if (degree == 1) nodes[cs] = nodesTetP1;
92149c89c76SBlaise Bourdin else if (degree == 2) nodes[cs] = nodesTetP2;
92249c89c76SBlaise Bourdin } else if (type[cs] == HEX) {
92349c89c76SBlaise Bourdin if (degree == 1) nodes[cs] = nodesHexP1;
92449c89c76SBlaise Bourdin else if (degree == 2) nodes[cs] = nodesHexP2;
92549c89c76SBlaise Bourdin }
92649c89c76SBlaise Bourdin /* Compute the number of cells not in the connectivity table */
92749c89c76SBlaise Bourdin cellsNotInConnectivity -= nodes[cs][3] * csSize;
92849c89c76SBlaise Bourdin
92949c89c76SBlaise Bourdin PetscCall(ISRestoreIndices(stratumIS, &cells));
93049c89c76SBlaise Bourdin PetscCall(ISDestroy(&stratumIS));
93149c89c76SBlaise Bourdin }
93249c89c76SBlaise Bourdin if (num_cs) PetscCallExternal(ex_put_init, exo->exoid, dmName, dim, numNodes, numCells, num_cs, num_vs, num_fs);
93349c89c76SBlaise Bourdin /* --- Connectivity --- */
93449c89c76SBlaise Bourdin for (cs = 0; cs < num_cs; ++cs) {
93549c89c76SBlaise Bourdin IS stratumIS;
93649c89c76SBlaise Bourdin const PetscInt *cells;
93749c89c76SBlaise Bourdin PetscInt *connect, off = 0;
93849c89c76SBlaise Bourdin PetscInt edgesInClosure = 0, facesInClosure = 0, verticesInClosure = 0;
93949c89c76SBlaise Bourdin PetscInt csSize, c, connectSize, closureSize;
94049c89c76SBlaise Bourdin char *elem_type = NULL;
94149c89c76SBlaise Bourdin char elem_type_bar2[] = "BAR2", elem_type_bar3[] = "BAR3";
94249c89c76SBlaise Bourdin char elem_type_tri3[] = "TRI3", elem_type_quad4[] = "QUAD4";
94349c89c76SBlaise Bourdin char elem_type_tri6[] = "TRI6", elem_type_quad9[] = "QUAD9";
94449c89c76SBlaise Bourdin char elem_type_tet4[] = "TET4", elem_type_hex8[] = "HEX8";
94549c89c76SBlaise Bourdin char elem_type_tet10[] = "TET10", elem_type_hex27[] = "HEX27";
94649c89c76SBlaise Bourdin
94749c89c76SBlaise Bourdin PetscCall(DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS));
94849c89c76SBlaise Bourdin PetscCall(ISGetIndices(stratumIS, &cells));
94949c89c76SBlaise Bourdin PetscCall(ISGetSize(stratumIS, &csSize));
95049c89c76SBlaise Bourdin /* Set Element type */
95149c89c76SBlaise Bourdin if (type[cs] == SEGMENT) {
95249c89c76SBlaise Bourdin if (degree == 1) elem_type = elem_type_bar2;
95349c89c76SBlaise Bourdin else if (degree == 2) elem_type = elem_type_bar3;
95449c89c76SBlaise Bourdin } else if (type[cs] == TRI) {
95549c89c76SBlaise Bourdin if (degree == 1) elem_type = elem_type_tri3;
95649c89c76SBlaise Bourdin else if (degree == 2) elem_type = elem_type_tri6;
95749c89c76SBlaise Bourdin } else if (type[cs] == QUAD) {
95849c89c76SBlaise Bourdin if (degree == 1) elem_type = elem_type_quad4;
95949c89c76SBlaise Bourdin else if (degree == 2) elem_type = elem_type_quad9;
96049c89c76SBlaise Bourdin } else if (type[cs] == TET) {
96149c89c76SBlaise Bourdin if (degree == 1) elem_type = elem_type_tet4;
96249c89c76SBlaise Bourdin else if (degree == 2) elem_type = elem_type_tet10;
96349c89c76SBlaise Bourdin } else if (type[cs] == HEX) {
96449c89c76SBlaise Bourdin if (degree == 1) elem_type = elem_type_hex8;
96549c89c76SBlaise Bourdin else if (degree == 2) elem_type = elem_type_hex27;
96649c89c76SBlaise Bourdin }
96749c89c76SBlaise Bourdin connectSize = nodes[cs][0] + nodes[cs][1] + nodes[cs][2] + nodes[cs][3];
96849c89c76SBlaise Bourdin PetscCall(PetscMalloc1(PetscMax(27, connectSize) * csSize, &connect));
96949c89c76SBlaise Bourdin PetscCallExternal(ex_put_block, exo->exoid, EX_ELEM_BLOCK, csIdx[cs], elem_type, csSize, connectSize, 0, 0, 1);
97049c89c76SBlaise Bourdin /* Find number of vertices, edges, and faces in the closure */
97149c89c76SBlaise Bourdin verticesInClosure = nodes[cs][0];
97249c89c76SBlaise Bourdin if (depth > 1) {
97349c89c76SBlaise Bourdin if (dim == 2) {
97449c89c76SBlaise Bourdin PetscCall(DMPlexGetConeSize(dm, cells[0], &edgesInClosure));
97549c89c76SBlaise Bourdin } else if (dim == 3) {
97649c89c76SBlaise Bourdin PetscInt *closure = NULL;
97749c89c76SBlaise Bourdin
97849c89c76SBlaise Bourdin PetscCall(DMPlexGetConeSize(dm, cells[0], &facesInClosure));
97949c89c76SBlaise Bourdin PetscCall(DMPlexGetTransitiveClosure(dm, cells[0], PETSC_TRUE, &closureSize, &closure));
98049c89c76SBlaise Bourdin edgesInClosure = closureSize - facesInClosure - 1 - verticesInClosure;
98149c89c76SBlaise Bourdin PetscCall(DMPlexRestoreTransitiveClosure(dm, cells[0], PETSC_TRUE, &closureSize, &closure));
98249c89c76SBlaise Bourdin }
98349c89c76SBlaise Bourdin }
98449c89c76SBlaise Bourdin /* Get connectivity for each cell */
98549c89c76SBlaise Bourdin for (c = 0; c < csSize; ++c) {
98649c89c76SBlaise Bourdin PetscInt *closure = NULL;
98749c89c76SBlaise Bourdin PetscInt temp, i;
98849c89c76SBlaise Bourdin
98949c89c76SBlaise Bourdin PetscCall(DMPlexGetTransitiveClosure(dm, cells[c], PETSC_TRUE, &closureSize, &closure));
99049c89c76SBlaise Bourdin for (i = 0; i < connectSize; ++i) {
99149c89c76SBlaise Bourdin if (i < nodes[cs][0]) { /* Vertices */
99249c89c76SBlaise Bourdin connect[i + off] = closure[(i + edgesInClosure + facesInClosure + 1) * 2] + 1;
99349c89c76SBlaise Bourdin connect[i + off] -= cellsNotInConnectivity;
99449c89c76SBlaise Bourdin } else if (i < nodes[cs][0] + nodes[cs][1]) { /* Edges */
99549c89c76SBlaise Bourdin connect[i + off] = closure[(i - verticesInClosure + facesInClosure + 1) * 2] + 1;
99649c89c76SBlaise Bourdin if (nodes[cs][2] == 0) connect[i + off] -= numFaces;
99749c89c76SBlaise Bourdin connect[i + off] -= cellsNotInConnectivity;
99849c89c76SBlaise Bourdin } else if (i < nodes[cs][0] + nodes[cs][1] + nodes[cs][3]) { /* Cells */
99949c89c76SBlaise Bourdin connect[i + off] = closure[0] + 1;
100049c89c76SBlaise Bourdin connect[i + off] -= skipCells;
100149c89c76SBlaise Bourdin } else if (i < nodes[cs][0] + nodes[cs][1] + nodes[cs][3] + nodes[cs][2]) { /* Faces */
100249c89c76SBlaise Bourdin connect[i + off] = closure[(i - edgesInClosure - verticesInClosure) * 2] + 1;
100349c89c76SBlaise Bourdin connect[i + off] -= cellsNotInConnectivity;
100449c89c76SBlaise Bourdin } else {
100549c89c76SBlaise Bourdin connect[i + off] = -1;
100649c89c76SBlaise Bourdin }
100749c89c76SBlaise Bourdin }
100849c89c76SBlaise Bourdin /* Tetrahedra are inverted */
100949c89c76SBlaise Bourdin if (type[cs] == TET) {
101049c89c76SBlaise Bourdin temp = connect[0 + off];
101149c89c76SBlaise Bourdin connect[0 + off] = connect[1 + off];
101249c89c76SBlaise Bourdin connect[1 + off] = temp;
101349c89c76SBlaise Bourdin if (degree == 2) {
101449c89c76SBlaise Bourdin temp = connect[5 + off];
101549c89c76SBlaise Bourdin connect[5 + off] = connect[6 + off];
101649c89c76SBlaise Bourdin connect[6 + off] = temp;
101749c89c76SBlaise Bourdin temp = connect[7 + off];
101849c89c76SBlaise Bourdin connect[7 + off] = connect[8 + off];
101949c89c76SBlaise Bourdin connect[8 + off] = temp;
102049c89c76SBlaise Bourdin }
102149c89c76SBlaise Bourdin }
102249c89c76SBlaise Bourdin /* Hexahedra are inverted */
102349c89c76SBlaise Bourdin if (type[cs] == HEX) {
102449c89c76SBlaise Bourdin temp = connect[1 + off];
102549c89c76SBlaise Bourdin connect[1 + off] = connect[3 + off];
102649c89c76SBlaise Bourdin connect[3 + off] = temp;
102749c89c76SBlaise Bourdin if (degree == 2) {
102849c89c76SBlaise Bourdin temp = connect[8 + off];
102949c89c76SBlaise Bourdin connect[8 + off] = connect[11 + off];
103049c89c76SBlaise Bourdin connect[11 + off] = temp;
103149c89c76SBlaise Bourdin temp = connect[9 + off];
103249c89c76SBlaise Bourdin connect[9 + off] = connect[10 + off];
103349c89c76SBlaise Bourdin connect[10 + off] = temp;
103449c89c76SBlaise Bourdin temp = connect[16 + off];
103549c89c76SBlaise Bourdin connect[16 + off] = connect[17 + off];
103649c89c76SBlaise Bourdin connect[17 + off] = temp;
103749c89c76SBlaise Bourdin temp = connect[18 + off];
103849c89c76SBlaise Bourdin connect[18 + off] = connect[19 + off];
103949c89c76SBlaise Bourdin connect[19 + off] = temp;
104049c89c76SBlaise Bourdin
104149c89c76SBlaise Bourdin temp = connect[12 + off];
104249c89c76SBlaise Bourdin connect[12 + off] = connect[16 + off];
104349c89c76SBlaise Bourdin connect[16 + off] = temp;
104449c89c76SBlaise Bourdin temp = connect[13 + off];
104549c89c76SBlaise Bourdin connect[13 + off] = connect[17 + off];
104649c89c76SBlaise Bourdin connect[17 + off] = temp;
104749c89c76SBlaise Bourdin temp = connect[14 + off];
104849c89c76SBlaise Bourdin connect[14 + off] = connect[18 + off];
104949c89c76SBlaise Bourdin connect[18 + off] = temp;
105049c89c76SBlaise Bourdin temp = connect[15 + off];
105149c89c76SBlaise Bourdin connect[15 + off] = connect[19 + off];
105249c89c76SBlaise Bourdin connect[19 + off] = temp;
105349c89c76SBlaise Bourdin
105449c89c76SBlaise Bourdin temp = connect[23 + off];
105549c89c76SBlaise Bourdin connect[23 + off] = connect[26 + off];
105649c89c76SBlaise Bourdin connect[26 + off] = temp;
105749c89c76SBlaise Bourdin temp = connect[24 + off];
105849c89c76SBlaise Bourdin connect[24 + off] = connect[25 + off];
105949c89c76SBlaise Bourdin connect[25 + off] = temp;
106049c89c76SBlaise Bourdin temp = connect[25 + off];
106149c89c76SBlaise Bourdin connect[25 + off] = connect[26 + off];
106249c89c76SBlaise Bourdin connect[26 + off] = temp;
106349c89c76SBlaise Bourdin }
106449c89c76SBlaise Bourdin }
106549c89c76SBlaise Bourdin off += connectSize;
106649c89c76SBlaise Bourdin PetscCall(DMPlexRestoreTransitiveClosure(dm, cells[c], PETSC_TRUE, &closureSize, &closure));
106749c89c76SBlaise Bourdin }
106849c89c76SBlaise Bourdin PetscCallExternal(ex_put_conn, exo->exoid, EX_ELEM_BLOCK, csIdx[cs], connect, 0, 0);
106949c89c76SBlaise Bourdin skipCells += (nodes[cs][3] == 0) * csSize;
107049c89c76SBlaise Bourdin PetscCall(PetscFree(connect));
107149c89c76SBlaise Bourdin PetscCall(ISRestoreIndices(stratumIS, &cells));
107249c89c76SBlaise Bourdin PetscCall(ISDestroy(&stratumIS));
107349c89c76SBlaise Bourdin }
107449c89c76SBlaise Bourdin PetscCall(PetscFree(type));
107549c89c76SBlaise Bourdin /* --- Coordinates --- */
107649c89c76SBlaise Bourdin PetscCall(PetscSectionSetChart(coordSection, pStart, pEnd));
107749c89c76SBlaise Bourdin if (num_cs) {
107849c89c76SBlaise Bourdin for (d = 0; d < depth; ++d) {
107949c89c76SBlaise Bourdin PetscCall(DMPlexGetDepthStratum(dm, d, &pStart, &pEnd));
108049c89c76SBlaise Bourdin for (p = pStart; p < pEnd; ++p) PetscCall(PetscSectionSetDof(coordSection, p, nodes[0][d] > 0));
108149c89c76SBlaise Bourdin }
108249c89c76SBlaise Bourdin }
108349c89c76SBlaise Bourdin for (cs = 0; cs < num_cs; ++cs) {
108449c89c76SBlaise Bourdin IS stratumIS;
108549c89c76SBlaise Bourdin const PetscInt *cells;
108649c89c76SBlaise Bourdin PetscInt csSize, c;
108749c89c76SBlaise Bourdin
108849c89c76SBlaise Bourdin PetscCall(DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS));
108949c89c76SBlaise Bourdin PetscCall(ISGetIndices(stratumIS, &cells));
109049c89c76SBlaise Bourdin PetscCall(ISGetSize(stratumIS, &csSize));
109149c89c76SBlaise Bourdin for (c = 0; c < csSize; ++c) PetscCall(PetscSectionSetDof(coordSection, cells[c], nodes[cs][3] > 0));
109249c89c76SBlaise Bourdin PetscCall(ISRestoreIndices(stratumIS, &cells));
109349c89c76SBlaise Bourdin PetscCall(ISDestroy(&stratumIS));
109449c89c76SBlaise Bourdin }
109549c89c76SBlaise Bourdin if (num_cs) {
109649c89c76SBlaise Bourdin PetscCall(ISRestoreIndices(csIS, &csIdx));
109749c89c76SBlaise Bourdin PetscCall(ISDestroy(&csIS));
109849c89c76SBlaise Bourdin }
109949c89c76SBlaise Bourdin PetscCall(PetscFree(nodes));
110049c89c76SBlaise Bourdin PetscCall(PetscSectionSetUp(coordSection));
110149c89c76SBlaise Bourdin if (numNodes) {
110249c89c76SBlaise Bourdin const char *coordNames[3] = {"x", "y", "z"};
110349c89c76SBlaise Bourdin PetscScalar *closure, *cval;
110449c89c76SBlaise Bourdin PetscReal *coords;
110549c89c76SBlaise Bourdin PetscInt hasDof, n = 0;
110649c89c76SBlaise Bourdin
110749c89c76SBlaise Bourdin /* There can't be more than 24 values in the closure of a point for the coord coordSection */
110849c89c76SBlaise Bourdin PetscCall(PetscCalloc3(numNodes * 3, &coords, dim, &cval, 24, &closure));
110949c89c76SBlaise Bourdin PetscCall(DMGetCoordinatesLocalNoncollective(dm, &coord));
111049c89c76SBlaise Bourdin PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
111149c89c76SBlaise Bourdin for (p = pStart; p < pEnd; ++p) {
111249c89c76SBlaise Bourdin PetscCall(PetscSectionGetDof(coordSection, p, &hasDof));
111349c89c76SBlaise Bourdin if (hasDof) {
111449c89c76SBlaise Bourdin PetscInt closureSize = 24, j;
111549c89c76SBlaise Bourdin
111649c89c76SBlaise Bourdin PetscCall(DMPlexVecGetClosure(cdm, NULL, coord, p, &closureSize, &closure));
111749c89c76SBlaise Bourdin for (d = 0; d < dim; ++d) {
111849c89c76SBlaise Bourdin cval[d] = 0.0;
111949c89c76SBlaise Bourdin for (j = 0; j < closureSize / dim; j++) cval[d] += closure[j * dim + d];
112049c89c76SBlaise Bourdin coords[d * numNodes + n] = PetscRealPart(cval[d]) * dim / closureSize;
112149c89c76SBlaise Bourdin }
112249c89c76SBlaise Bourdin ++n;
112349c89c76SBlaise Bourdin }
112449c89c76SBlaise Bourdin }
112549c89c76SBlaise Bourdin PetscCallExternal(ex_put_coord, exo->exoid, &coords[0 * numNodes], &coords[1 * numNodes], &coords[2 * numNodes]);
112649c89c76SBlaise Bourdin PetscCall(PetscFree3(coords, cval, closure));
112749c89c76SBlaise Bourdin PetscCallExternal(ex_put_coord_names, exo->exoid, (char **)coordNames);
112849c89c76SBlaise Bourdin }
112949c89c76SBlaise Bourdin
113049c89c76SBlaise Bourdin /* --- Node Sets/Vertex Sets --- */
113149c89c76SBlaise Bourdin PetscCall(DMHasLabel(dm, "Vertex Sets", &hasLabel));
113249c89c76SBlaise Bourdin if (hasLabel) {
113349c89c76SBlaise Bourdin PetscInt i, vs, vsSize;
113449c89c76SBlaise Bourdin const PetscInt *vsIdx, *vertices;
113549c89c76SBlaise Bourdin PetscInt *nodeList;
113649c89c76SBlaise Bourdin IS vsIS, stratumIS;
113749c89c76SBlaise Bourdin DMLabel vsLabel;
113849c89c76SBlaise Bourdin PetscCall(DMGetLabel(dm, "Vertex Sets", &vsLabel));
113949c89c76SBlaise Bourdin PetscCall(DMLabelGetValueIS(vsLabel, &vsIS));
114049c89c76SBlaise Bourdin PetscCall(ISGetIndices(vsIS, &vsIdx));
114149c89c76SBlaise Bourdin for (vs = 0; vs < num_vs; ++vs) {
114249c89c76SBlaise Bourdin PetscCall(DMLabelGetStratumIS(vsLabel, vsIdx[vs], &stratumIS));
114349c89c76SBlaise Bourdin PetscCall(ISGetIndices(stratumIS, &vertices));
114449c89c76SBlaise Bourdin PetscCall(ISGetSize(stratumIS, &vsSize));
114549c89c76SBlaise Bourdin PetscCall(PetscMalloc1(vsSize, &nodeList));
114649c89c76SBlaise Bourdin for (i = 0; i < vsSize; ++i) nodeList[i] = vertices[i] - skipCells + 1;
114749c89c76SBlaise Bourdin PetscCallExternal(ex_put_set_param, exo->exoid, EX_NODE_SET, vsIdx[vs], vsSize, 0);
114849c89c76SBlaise Bourdin PetscCallExternal(ex_put_set, exo->exoid, EX_NODE_SET, vsIdx[vs], nodeList, NULL);
114949c89c76SBlaise Bourdin PetscCall(ISRestoreIndices(stratumIS, &vertices));
115049c89c76SBlaise Bourdin PetscCall(ISDestroy(&stratumIS));
115149c89c76SBlaise Bourdin PetscCall(PetscFree(nodeList));
115249c89c76SBlaise Bourdin }
115349c89c76SBlaise Bourdin PetscCall(ISRestoreIndices(vsIS, &vsIdx));
115449c89c76SBlaise Bourdin PetscCall(ISDestroy(&vsIS));
115549c89c76SBlaise Bourdin }
115649c89c76SBlaise Bourdin /* --- Side Sets/Face Sets --- */
115749c89c76SBlaise Bourdin PetscCall(DMHasLabel(dm, "Face Sets", &hasLabel));
115849c89c76SBlaise Bourdin if (hasLabel) {
115949c89c76SBlaise Bourdin PetscInt i, j, fs, fsSize;
116049c89c76SBlaise Bourdin const PetscInt *fsIdx, *faces;
116149c89c76SBlaise Bourdin IS fsIS, stratumIS;
116249c89c76SBlaise Bourdin DMLabel fsLabel;
116349c89c76SBlaise Bourdin PetscInt numPoints, *points;
116449c89c76SBlaise Bourdin PetscInt elem_list_size = 0;
116549c89c76SBlaise Bourdin PetscInt *elem_list, *elem_ind, *side_list;
116649c89c76SBlaise Bourdin
116749c89c76SBlaise Bourdin PetscCall(DMGetLabel(dm, "Face Sets", &fsLabel));
116849c89c76SBlaise Bourdin /* Compute size of Node List and Element List */
116949c89c76SBlaise Bourdin PetscCall(DMLabelGetValueIS(fsLabel, &fsIS));
117049c89c76SBlaise Bourdin PetscCall(ISGetIndices(fsIS, &fsIdx));
117149c89c76SBlaise Bourdin for (fs = 0; fs < num_fs; ++fs) {
117249c89c76SBlaise Bourdin PetscCall(DMLabelGetStratumIS(fsLabel, fsIdx[fs], &stratumIS));
117349c89c76SBlaise Bourdin PetscCall(ISGetSize(stratumIS, &fsSize));
117449c89c76SBlaise Bourdin elem_list_size += fsSize;
117549c89c76SBlaise Bourdin PetscCall(ISDestroy(&stratumIS));
117649c89c76SBlaise Bourdin }
117749c89c76SBlaise Bourdin if (num_fs) {
117849c89c76SBlaise Bourdin PetscCall(PetscMalloc3(num_fs, &elem_ind, elem_list_size, &elem_list, elem_list_size, &side_list));
117949c89c76SBlaise Bourdin elem_ind[0] = 0;
118049c89c76SBlaise Bourdin for (fs = 0; fs < num_fs; ++fs) {
118149c89c76SBlaise Bourdin PetscCall(DMLabelGetStratumIS(fsLabel, fsIdx[fs], &stratumIS));
118249c89c76SBlaise Bourdin PetscCall(ISGetIndices(stratumIS, &faces));
118349c89c76SBlaise Bourdin PetscCall(ISGetSize(stratumIS, &fsSize));
118449c89c76SBlaise Bourdin /* Set Parameters */
118549c89c76SBlaise Bourdin PetscCallExternal(ex_put_set_param, exo->exoid, EX_SIDE_SET, fsIdx[fs], fsSize, 0);
118649c89c76SBlaise Bourdin /* Indices */
118749c89c76SBlaise Bourdin if (fs < num_fs - 1) elem_ind[fs + 1] = elem_ind[fs] + fsSize;
118849c89c76SBlaise Bourdin
118949c89c76SBlaise Bourdin for (i = 0; i < fsSize; ++i) {
119049c89c76SBlaise Bourdin /* Element List */
119149c89c76SBlaise Bourdin points = NULL;
119249c89c76SBlaise Bourdin PetscCall(DMPlexGetTransitiveClosure(dm, faces[i], PETSC_FALSE, &numPoints, &points));
119349c89c76SBlaise Bourdin elem_list[elem_ind[fs] + i] = points[2] + 1;
119449c89c76SBlaise Bourdin PetscCall(DMPlexRestoreTransitiveClosure(dm, faces[i], PETSC_FALSE, &numPoints, &points));
119549c89c76SBlaise Bourdin
119649c89c76SBlaise Bourdin /* Side List */
119749c89c76SBlaise Bourdin points = NULL;
119849c89c76SBlaise Bourdin PetscCall(DMPlexGetTransitiveClosure(dm, elem_list[elem_ind[fs] + i] - 1, PETSC_TRUE, &numPoints, &points));
119949c89c76SBlaise Bourdin for (j = 1; j < numPoints; ++j) {
120049c89c76SBlaise Bourdin if (points[j * 2] == faces[i]) break;
120149c89c76SBlaise Bourdin }
120249c89c76SBlaise Bourdin /* Convert HEX sides */
120349c89c76SBlaise Bourdin if (numPoints == 27) {
120449c89c76SBlaise Bourdin if (j == 1) {
120549c89c76SBlaise Bourdin j = 5;
120649c89c76SBlaise Bourdin } else if (j == 2) {
120749c89c76SBlaise Bourdin j = 6;
120849c89c76SBlaise Bourdin } else if (j == 3) {
120949c89c76SBlaise Bourdin j = 1;
121049c89c76SBlaise Bourdin } else if (j == 4) {
121149c89c76SBlaise Bourdin j = 3;
121249c89c76SBlaise Bourdin } else if (j == 5) {
121349c89c76SBlaise Bourdin j = 2;
121449c89c76SBlaise Bourdin } else if (j == 6) {
121549c89c76SBlaise Bourdin j = 4;
121649c89c76SBlaise Bourdin }
121749c89c76SBlaise Bourdin }
121849c89c76SBlaise Bourdin /* Convert TET sides */
121949c89c76SBlaise Bourdin if (numPoints == 15) {
122049c89c76SBlaise Bourdin --j;
122149c89c76SBlaise Bourdin if (j == 0) j = 4;
122249c89c76SBlaise Bourdin }
122349c89c76SBlaise Bourdin side_list[elem_ind[fs] + i] = j;
122449c89c76SBlaise Bourdin PetscCall(DMPlexRestoreTransitiveClosure(dm, elem_list[elem_ind[fs] + i] - 1, PETSC_TRUE, &numPoints, &points));
122549c89c76SBlaise Bourdin }
122649c89c76SBlaise Bourdin PetscCall(ISRestoreIndices(stratumIS, &faces));
122749c89c76SBlaise Bourdin PetscCall(ISDestroy(&stratumIS));
122849c89c76SBlaise Bourdin }
122949c89c76SBlaise Bourdin PetscCall(ISRestoreIndices(fsIS, &fsIdx));
123049c89c76SBlaise Bourdin PetscCall(ISDestroy(&fsIS));
123149c89c76SBlaise Bourdin
123249c89c76SBlaise Bourdin /* Put side sets */
123349c89c76SBlaise Bourdin for (fs = 0; fs < num_fs; ++fs) PetscCallExternal(ex_put_set, exo->exoid, EX_SIDE_SET, fsIdx[fs], &elem_list[elem_ind[fs]], &side_list[elem_ind[fs]]);
123449c89c76SBlaise Bourdin PetscCall(PetscFree3(elem_ind, elem_list, side_list));
123549c89c76SBlaise Bourdin }
123649c89c76SBlaise Bourdin }
123749c89c76SBlaise Bourdin /*
123849c89c76SBlaise Bourdin close the exodus file
123949c89c76SBlaise Bourdin */
124049c89c76SBlaise Bourdin ex_close(exo->exoid);
124149c89c76SBlaise Bourdin exo->exoid = -1;
124249c89c76SBlaise Bourdin }
124349c89c76SBlaise Bourdin PetscCall(PetscSectionDestroy(&coordSection));
124449c89c76SBlaise Bourdin
124549c89c76SBlaise Bourdin /*
124649c89c76SBlaise Bourdin reopen the file in parallel
124749c89c76SBlaise Bourdin */
124849c89c76SBlaise Bourdin EXO_mode = EX_WRITE;
124949c89c76SBlaise Bourdin #if defined(PETSC_USE_64BIT_INDICES)
125049c89c76SBlaise Bourdin EXO_mode += EX_ALL_INT64_API;
125149c89c76SBlaise Bourdin #endif
125249c89c76SBlaise Bourdin CPU_word_size = sizeof(PetscReal);
125349c89c76SBlaise Bourdin IO_word_size = sizeof(PetscReal);
125449c89c76SBlaise Bourdin exo->exoid = ex_open_par(exo->filename, EXO_mode, &CPU_word_size, &IO_word_size, &EXO_version, comm, MPI_INFO_NULL);
125549c89c76SBlaise Bourdin PetscCheck(exo->exoid >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_open_par failed for %s", exo->filename);
125649c89c76SBlaise Bourdin PetscFunctionReturn(PETSC_SUCCESS);
125749c89c76SBlaise Bourdin }
125849c89c76SBlaise Bourdin
12590a5cf5d8SBlaise Bourdin static PetscErrorCode VecViewPlex_ExodusII_Nodal_Internal(Vec v, PetscExodusIIInt exoid, PetscExodusIIInt step, PetscExodusIIInt offset);
12600a5cf5d8SBlaise Bourdin static PetscErrorCode VecViewPlex_ExodusII_Zonal_Internal(Vec v, PetscExodusIIInt exoid, PetscExodusIIInt step, PetscExodusIIInt offset);
12610a5cf5d8SBlaise Bourdin static PetscErrorCode VecLoadPlex_ExodusII_Nodal_Internal(Vec v, PetscExodusIIInt exoid, PetscExodusIIInt step, PetscExodusIIInt offset);
12620a5cf5d8SBlaise Bourdin static PetscErrorCode VecLoadPlex_ExodusII_Zonal_Internal(Vec v, PetscExodusIIInt exoid, PetscExodusIIInt step, PetscExodusIIInt offset);
126349c89c76SBlaise Bourdin
VecView_PlexExodusII_Internal(Vec v,PetscViewer viewer)126449c89c76SBlaise Bourdin PetscErrorCode VecView_PlexExodusII_Internal(Vec v, PetscViewer viewer)
126549c89c76SBlaise Bourdin {
126649c89c76SBlaise Bourdin DM dm;
126749c89c76SBlaise Bourdin MPI_Comm comm;
126849c89c76SBlaise Bourdin PetscMPIInt rank;
126949c89c76SBlaise Bourdin
12700a5cf5d8SBlaise Bourdin PetscExodusIIInt exoid, offsetN = -1, offsetZ = -1;
127149c89c76SBlaise Bourdin const char *vecname;
127249c89c76SBlaise Bourdin PetscInt step;
127349c89c76SBlaise Bourdin
127449c89c76SBlaise Bourdin PetscFunctionBegin;
127549c89c76SBlaise Bourdin PetscCall(PetscObjectGetComm((PetscObject)v, &comm));
127649c89c76SBlaise Bourdin PetscCallMPI(MPI_Comm_rank(comm, &rank));
127749c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIIGetId(viewer, &exoid));
127849c89c76SBlaise Bourdin PetscCall(VecGetDM(v, &dm));
127949c89c76SBlaise Bourdin PetscCall(PetscObjectGetName((PetscObject)v, &vecname));
128049c89c76SBlaise Bourdin
128149c89c76SBlaise Bourdin PetscCall(DMGetOutputSequenceNumber(dm, &step, NULL));
128249c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIIGetNodalVariableIndex(viewer, vecname, &offsetN));
128349c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIIGetZonalVariableIndex(viewer, vecname, &offsetZ));
128449c89c76SBlaise Bourdin PetscCheck(!(offsetN >= 0 && offsetZ >= 0), comm, PETSC_ERR_FILE_UNEXPECTED, "Found both nodal and zonal variable %s in exodus file. ", vecname);
128549c89c76SBlaise Bourdin if (offsetN >= 0) {
12860a5cf5d8SBlaise Bourdin PetscCall(VecViewPlex_ExodusII_Nodal_Internal(v, exoid, (PetscExodusIIInt)step + 1, offsetN + 1));
128749c89c76SBlaise Bourdin } else if (offsetZ >= 0) {
12880a5cf5d8SBlaise Bourdin PetscCall(VecViewPlex_ExodusII_Zonal_Internal(v, exoid, (PetscExodusIIInt)step + 1, offsetZ + 1));
128949c89c76SBlaise Bourdin } else SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "Could not find nodal or zonal variable %s in exodus file. ", vecname);
129049c89c76SBlaise Bourdin PetscFunctionReturn(PETSC_SUCCESS);
129149c89c76SBlaise Bourdin }
129249c89c76SBlaise Bourdin
VecLoad_PlexExodusII_Internal(Vec v,PetscViewer viewer)129349c89c76SBlaise Bourdin PetscErrorCode VecLoad_PlexExodusII_Internal(Vec v, PetscViewer viewer)
129449c89c76SBlaise Bourdin {
129549c89c76SBlaise Bourdin DM dm;
129649c89c76SBlaise Bourdin MPI_Comm comm;
129749c89c76SBlaise Bourdin PetscMPIInt rank;
129849c89c76SBlaise Bourdin
12990a5cf5d8SBlaise Bourdin PetscExodusIIInt exoid, offsetN = 0, offsetZ = 0;
130049c89c76SBlaise Bourdin const char *vecname;
130149c89c76SBlaise Bourdin PetscInt step;
130249c89c76SBlaise Bourdin
130349c89c76SBlaise Bourdin PetscFunctionBegin;
130449c89c76SBlaise Bourdin PetscCall(PetscObjectGetComm((PetscObject)v, &comm));
130549c89c76SBlaise Bourdin PetscCallMPI(MPI_Comm_rank(comm, &rank));
130649c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIIGetId(viewer, &exoid));
130749c89c76SBlaise Bourdin PetscCall(VecGetDM(v, &dm));
130849c89c76SBlaise Bourdin PetscCall(PetscObjectGetName((PetscObject)v, &vecname));
130949c89c76SBlaise Bourdin
131049c89c76SBlaise Bourdin PetscCall(DMGetOutputSequenceNumber(dm, &step, NULL));
131149c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIIGetNodalVariableIndex(viewer, vecname, &offsetN));
131249c89c76SBlaise Bourdin PetscCall(PetscViewerExodusIIGetZonalVariableIndex(viewer, vecname, &offsetZ));
131349c89c76SBlaise Bourdin PetscCheck(!(offsetN >= 0 && offsetZ >= 0), comm, PETSC_ERR_FILE_UNEXPECTED, "Found both nodal and zonal variable %s in exodus file. ", vecname);
13140a5cf5d8SBlaise Bourdin if (offsetN >= 0) PetscCall(VecLoadPlex_ExodusII_Nodal_Internal(v, exoid, (PetscExodusIIInt)step + 1, offsetN + 1));
13150a5cf5d8SBlaise Bourdin else if (offsetZ >= 0) PetscCall(VecLoadPlex_ExodusII_Zonal_Internal(v, exoid, (PetscExodusIIInt)step + 1, offsetZ + 1));
131649c89c76SBlaise Bourdin else SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "Could not find nodal or zonal variable %s in exodus file. ", vecname);
131749c89c76SBlaise Bourdin PetscFunctionReturn(PETSC_SUCCESS);
131849c89c76SBlaise Bourdin }
131949c89c76SBlaise Bourdin
VecViewPlex_ExodusII_Nodal_Internal(Vec v,PetscExodusIIInt exoid,PetscExodusIIInt step,PetscExodusIIInt offset)13200a5cf5d8SBlaise Bourdin static PetscErrorCode VecViewPlex_ExodusII_Nodal_Internal(Vec v, PetscExodusIIInt exoid, PetscExodusIIInt step, PetscExodusIIInt offset)
132149c89c76SBlaise Bourdin {
132249c89c76SBlaise Bourdin MPI_Comm comm;
132349c89c76SBlaise Bourdin PetscMPIInt size;
132449c89c76SBlaise Bourdin DM dm;
132549c89c76SBlaise Bourdin Vec vNatural, vComp;
132649c89c76SBlaise Bourdin const PetscScalar *varray;
132749c89c76SBlaise Bourdin PetscInt xs, xe, bs;
132849c89c76SBlaise Bourdin PetscBool useNatural;
132949c89c76SBlaise Bourdin
133049c89c76SBlaise Bourdin PetscFunctionBegin;
133149c89c76SBlaise Bourdin PetscCall(PetscObjectGetComm((PetscObject)v, &comm));
133249c89c76SBlaise Bourdin PetscCallMPI(MPI_Comm_size(comm, &size));
133349c89c76SBlaise Bourdin PetscCall(VecGetDM(v, &dm));
133449c89c76SBlaise Bourdin PetscCall(DMGetUseNatural(dm, &useNatural));
133549c89c76SBlaise Bourdin useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE;
133649c89c76SBlaise Bourdin if (useNatural) {
133749c89c76SBlaise Bourdin PetscCall(DMPlexCreateNaturalVector(dm, &vNatural));
133849c89c76SBlaise Bourdin PetscCall(DMPlexGlobalToNaturalBegin(dm, v, vNatural));
133949c89c76SBlaise Bourdin PetscCall(DMPlexGlobalToNaturalEnd(dm, v, vNatural));
134049c89c76SBlaise Bourdin } else {
134149c89c76SBlaise Bourdin vNatural = v;
134249c89c76SBlaise Bourdin }
134349c89c76SBlaise Bourdin
134449c89c76SBlaise Bourdin /* Write local chunk of the result in the exodus file
134549c89c76SBlaise Bourdin exodus stores each component of a vector-valued field as a separate variable.
134649c89c76SBlaise Bourdin We assume that they are stored sequentially */
134749c89c76SBlaise Bourdin PetscCall(VecGetOwnershipRange(vNatural, &xs, &xe));
134849c89c76SBlaise Bourdin PetscCall(VecGetBlockSize(vNatural, &bs));
134949c89c76SBlaise Bourdin if (bs == 1) {
135049c89c76SBlaise Bourdin PetscCall(VecGetArrayRead(vNatural, &varray));
135149c89c76SBlaise Bourdin PetscCallExternal(ex_put_partial_var, exoid, step, EX_NODAL, offset, 1, xs + 1, xe - xs, varray);
135249c89c76SBlaise Bourdin PetscCall(VecRestoreArrayRead(vNatural, &varray));
135349c89c76SBlaise Bourdin } else {
135449c89c76SBlaise Bourdin IS compIS;
135549c89c76SBlaise Bourdin PetscInt c;
135649c89c76SBlaise Bourdin
135749c89c76SBlaise Bourdin PetscCall(ISCreateStride(comm, (xe - xs) / bs, xs, bs, &compIS));
135849c89c76SBlaise Bourdin for (c = 0; c < bs; ++c) {
135949c89c76SBlaise Bourdin PetscCall(ISStrideSetStride(compIS, (xe - xs) / bs, xs + c, bs));
136049c89c76SBlaise Bourdin PetscCall(VecGetSubVector(vNatural, compIS, &vComp));
136149c89c76SBlaise Bourdin PetscCall(VecGetArrayRead(vComp, &varray));
136249c89c76SBlaise Bourdin PetscCallExternal(ex_put_partial_var, exoid, step, EX_NODAL, offset + c, 1, xs / bs + 1, (xe - xs) / bs, varray);
136349c89c76SBlaise Bourdin PetscCall(VecRestoreArrayRead(vComp, &varray));
136449c89c76SBlaise Bourdin PetscCall(VecRestoreSubVector(vNatural, compIS, &vComp));
136549c89c76SBlaise Bourdin }
136649c89c76SBlaise Bourdin PetscCall(ISDestroy(&compIS));
136749c89c76SBlaise Bourdin }
136849c89c76SBlaise Bourdin if (useNatural) PetscCall(VecDestroy(&vNatural));
136949c89c76SBlaise Bourdin PetscFunctionReturn(PETSC_SUCCESS);
137049c89c76SBlaise Bourdin }
137149c89c76SBlaise Bourdin
VecLoadPlex_ExodusII_Nodal_Internal(Vec v,PetscExodusIIInt exoid,PetscExodusIIInt step,PetscExodusIIInt offset)13720a5cf5d8SBlaise Bourdin static PetscErrorCode VecLoadPlex_ExodusII_Nodal_Internal(Vec v, PetscExodusIIInt exoid, PetscExodusIIInt step, PetscExodusIIInt offset)
137349c89c76SBlaise Bourdin {
137449c89c76SBlaise Bourdin MPI_Comm comm;
137549c89c76SBlaise Bourdin PetscMPIInt size;
137649c89c76SBlaise Bourdin DM dm;
137749c89c76SBlaise Bourdin Vec vNatural, vComp;
137849c89c76SBlaise Bourdin PetscScalar *varray;
137949c89c76SBlaise Bourdin PetscInt xs, xe, bs;
138049c89c76SBlaise Bourdin PetscBool useNatural;
138149c89c76SBlaise Bourdin
138249c89c76SBlaise Bourdin PetscFunctionBegin;
138349c89c76SBlaise Bourdin PetscCall(PetscObjectGetComm((PetscObject)v, &comm));
138449c89c76SBlaise Bourdin PetscCallMPI(MPI_Comm_size(comm, &size));
138549c89c76SBlaise Bourdin PetscCall(VecGetDM(v, &dm));
138649c89c76SBlaise Bourdin PetscCall(DMGetUseNatural(dm, &useNatural));
138749c89c76SBlaise Bourdin useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE;
138849c89c76SBlaise Bourdin if (useNatural) PetscCall(DMPlexCreateNaturalVector(dm, &vNatural));
138949c89c76SBlaise Bourdin else vNatural = v;
139049c89c76SBlaise Bourdin
139149c89c76SBlaise Bourdin /* Read local chunk from the file */
139249c89c76SBlaise Bourdin PetscCall(VecGetOwnershipRange(vNatural, &xs, &xe));
139349c89c76SBlaise Bourdin PetscCall(VecGetBlockSize(vNatural, &bs));
139449c89c76SBlaise Bourdin if (bs == 1) {
139549c89c76SBlaise Bourdin PetscCall(VecGetArray(vNatural, &varray));
139649c89c76SBlaise Bourdin PetscCallExternal(ex_get_partial_var, exoid, step, EX_NODAL, offset, 1, xs + 1, xe - xs, varray);
139749c89c76SBlaise Bourdin PetscCall(VecRestoreArray(vNatural, &varray));
139849c89c76SBlaise Bourdin } else {
139949c89c76SBlaise Bourdin IS compIS;
140049c89c76SBlaise Bourdin PetscInt c;
140149c89c76SBlaise Bourdin
140249c89c76SBlaise Bourdin PetscCall(ISCreateStride(comm, (xe - xs) / bs, xs, bs, &compIS));
140349c89c76SBlaise Bourdin for (c = 0; c < bs; ++c) {
140449c89c76SBlaise Bourdin PetscCall(ISStrideSetStride(compIS, (xe - xs) / bs, xs + c, bs));
140549c89c76SBlaise Bourdin PetscCall(VecGetSubVector(vNatural, compIS, &vComp));
140649c89c76SBlaise Bourdin PetscCall(VecGetArray(vComp, &varray));
140749c89c76SBlaise Bourdin PetscCallExternal(ex_get_partial_var, exoid, step, EX_NODAL, offset + c, 1, xs / bs + 1, (xe - xs) / bs, varray);
140849c89c76SBlaise Bourdin PetscCall(VecRestoreArray(vComp, &varray));
140949c89c76SBlaise Bourdin PetscCall(VecRestoreSubVector(vNatural, compIS, &vComp));
141049c89c76SBlaise Bourdin }
141149c89c76SBlaise Bourdin PetscCall(ISDestroy(&compIS));
141249c89c76SBlaise Bourdin }
141349c89c76SBlaise Bourdin if (useNatural) {
141449c89c76SBlaise Bourdin PetscCall(DMPlexNaturalToGlobalBegin(dm, vNatural, v));
141549c89c76SBlaise Bourdin PetscCall(DMPlexNaturalToGlobalEnd(dm, vNatural, v));
141649c89c76SBlaise Bourdin PetscCall(VecDestroy(&vNatural));
141749c89c76SBlaise Bourdin }
141849c89c76SBlaise Bourdin PetscFunctionReturn(PETSC_SUCCESS);
141949c89c76SBlaise Bourdin }
142049c89c76SBlaise Bourdin
VecViewPlex_ExodusII_Zonal_Internal(Vec v,PetscExodusIIInt exoid,PetscExodusIIInt step,PetscExodusIIInt offset)14210a5cf5d8SBlaise Bourdin static PetscErrorCode VecViewPlex_ExodusII_Zonal_Internal(Vec v, PetscExodusIIInt exoid, PetscExodusIIInt step, PetscExodusIIInt offset)
142249c89c76SBlaise Bourdin {
142349c89c76SBlaise Bourdin MPI_Comm comm;
142449c89c76SBlaise Bourdin PetscMPIInt size;
142549c89c76SBlaise Bourdin DM dm;
142649c89c76SBlaise Bourdin Vec vNatural, vComp;
142749c89c76SBlaise Bourdin const PetscScalar *varray;
142849c89c76SBlaise Bourdin PetscInt xs, xe, bs;
142949c89c76SBlaise Bourdin PetscBool useNatural;
143049c89c76SBlaise Bourdin IS compIS;
143149c89c76SBlaise Bourdin PetscInt *csSize, *csID;
14320a5cf5d8SBlaise Bourdin PetscExodusIIInt numCS, set, csxs = 0;
143349c89c76SBlaise Bourdin
143449c89c76SBlaise Bourdin PetscFunctionBegin;
143549c89c76SBlaise Bourdin PetscCall(PetscObjectGetComm((PetscObject)v, &comm));
143649c89c76SBlaise Bourdin PetscCallMPI(MPI_Comm_size(comm, &size));
143749c89c76SBlaise Bourdin PetscCall(VecGetDM(v, &dm));
143849c89c76SBlaise Bourdin PetscCall(DMGetUseNatural(dm, &useNatural));
143949c89c76SBlaise Bourdin useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE;
144049c89c76SBlaise Bourdin if (useNatural) {
144149c89c76SBlaise Bourdin PetscCall(DMPlexCreateNaturalVector(dm, &vNatural));
144249c89c76SBlaise Bourdin PetscCall(DMPlexGlobalToNaturalBegin(dm, v, vNatural));
144349c89c76SBlaise Bourdin PetscCall(DMPlexGlobalToNaturalEnd(dm, v, vNatural));
144449c89c76SBlaise Bourdin } else {
144549c89c76SBlaise Bourdin vNatural = v;
144649c89c76SBlaise Bourdin }
144749c89c76SBlaise Bourdin
144849c89c76SBlaise Bourdin /* Write local chunk of the result in the exodus file
144949c89c76SBlaise Bourdin exodus stores each component of a vector-valued field as a separate variable.
145049c89c76SBlaise Bourdin We assume that they are stored sequentially
145149c89c76SBlaise Bourdin Zonal variables are accessed one element block at a time, so we loop through the cell sets,
145249c89c76SBlaise Bourdin but once the vector has been reordered to natural size, we cannot use the label information
145349c89c76SBlaise Bourdin to figure out what to save where. */
14548a1c49cdSMatthew G. Knepley numCS = (PetscExodusIIInt)ex_inquire_int(exoid, EX_INQ_ELEM_BLK); // This is an int64_t
145549c89c76SBlaise Bourdin PetscCall(PetscMalloc2(numCS, &csID, numCS, &csSize));
145649c89c76SBlaise Bourdin PetscCallExternal(ex_get_ids, exoid, EX_ELEM_BLOCK, csID);
145749c89c76SBlaise Bourdin for (set = 0; set < numCS; ++set) {
145849c89c76SBlaise Bourdin ex_block block;
145949c89c76SBlaise Bourdin
146049c89c76SBlaise Bourdin block.id = csID[set];
146149c89c76SBlaise Bourdin block.type = EX_ELEM_BLOCK;
146249c89c76SBlaise Bourdin PetscCallExternal(ex_get_block_param, exoid, &block);
146362e4c2b1SMatthew G. Knepley PetscCall(PetscIntCast(block.num_entry, &csSize[set])); // This is an int64_t
146449c89c76SBlaise Bourdin }
146549c89c76SBlaise Bourdin PetscCall(VecGetOwnershipRange(vNatural, &xs, &xe));
146649c89c76SBlaise Bourdin PetscCall(VecGetBlockSize(vNatural, &bs));
146749c89c76SBlaise Bourdin if (bs > 1) PetscCall(ISCreateStride(comm, (xe - xs) / bs, xs, bs, &compIS));
146849c89c76SBlaise Bourdin for (set = 0; set < numCS; set++) {
146949c89c76SBlaise Bourdin PetscInt csLocalSize, c;
147049c89c76SBlaise Bourdin
147149c89c76SBlaise Bourdin /* range of indices for set setID[set]: csxs:csxs + csSize[set]-1
147249c89c76SBlaise Bourdin local slice of zonal values: xs/bs,xm/bs-1
147349c89c76SBlaise Bourdin intersection: max(xs/bs,csxs),min(xm/bs-1,csxs + csSize[set]-1) */
147449c89c76SBlaise Bourdin csLocalSize = PetscMax(0, PetscMin(xe / bs, csxs + csSize[set]) - PetscMax(xs / bs, csxs));
147549c89c76SBlaise Bourdin if (bs == 1) {
147649c89c76SBlaise Bourdin PetscCall(VecGetArrayRead(vNatural, &varray));
147749c89c76SBlaise Bourdin PetscCallExternal(ex_put_partial_var, exoid, step, EX_ELEM_BLOCK, offset, csID[set], PetscMax(xs - csxs, 0) + 1, csLocalSize, &varray[PetscMax(0, csxs - xs)]);
147849c89c76SBlaise Bourdin PetscCall(VecRestoreArrayRead(vNatural, &varray));
147949c89c76SBlaise Bourdin } else {
148049c89c76SBlaise Bourdin for (c = 0; c < bs; ++c) {
148149c89c76SBlaise Bourdin PetscCall(ISStrideSetStride(compIS, (xe - xs) / bs, xs + c, bs));
148249c89c76SBlaise Bourdin PetscCall(VecGetSubVector(vNatural, compIS, &vComp));
148349c89c76SBlaise Bourdin PetscCall(VecGetArrayRead(vComp, &varray));
148449c89c76SBlaise Bourdin PetscCallExternal(ex_put_partial_var, exoid, step, EX_ELEM_BLOCK, offset + c, csID[set], PetscMax(xs / bs - csxs, 0) + 1, csLocalSize, &varray[PetscMax(0, csxs - xs / bs)]);
148549c89c76SBlaise Bourdin PetscCall(VecRestoreArrayRead(vComp, &varray));
148649c89c76SBlaise Bourdin PetscCall(VecRestoreSubVector(vNatural, compIS, &vComp));
148749c89c76SBlaise Bourdin }
148849c89c76SBlaise Bourdin }
148949c89c76SBlaise Bourdin csxs += csSize[set];
149049c89c76SBlaise Bourdin }
149149c89c76SBlaise Bourdin PetscCall(PetscFree2(csID, csSize));
149249c89c76SBlaise Bourdin if (bs > 1) PetscCall(ISDestroy(&compIS));
149349c89c76SBlaise Bourdin if (useNatural) PetscCall(VecDestroy(&vNatural));
149449c89c76SBlaise Bourdin PetscFunctionReturn(PETSC_SUCCESS);
149549c89c76SBlaise Bourdin }
149649c89c76SBlaise Bourdin
VecLoadPlex_ExodusII_Zonal_Internal(Vec v,PetscExodusIIInt exoid,PetscExodusIIInt step,PetscExodusIIInt offset)14970a5cf5d8SBlaise Bourdin static PetscErrorCode VecLoadPlex_ExodusII_Zonal_Internal(Vec v, PetscExodusIIInt exoid, PetscExodusIIInt step, PetscExodusIIInt offset)
149849c89c76SBlaise Bourdin {
149949c89c76SBlaise Bourdin MPI_Comm comm;
150049c89c76SBlaise Bourdin PetscMPIInt size;
150149c89c76SBlaise Bourdin DM dm;
150249c89c76SBlaise Bourdin Vec vNatural, vComp;
150349c89c76SBlaise Bourdin PetscScalar *varray;
150449c89c76SBlaise Bourdin PetscInt xs, xe, bs;
150549c89c76SBlaise Bourdin PetscBool useNatural;
150649c89c76SBlaise Bourdin IS compIS;
150749c89c76SBlaise Bourdin PetscInt *csSize, *csID;
15080a5cf5d8SBlaise Bourdin PetscExodusIIInt numCS, set, csxs = 0;
150949c89c76SBlaise Bourdin
151049c89c76SBlaise Bourdin PetscFunctionBegin;
151149c89c76SBlaise Bourdin PetscCall(PetscObjectGetComm((PetscObject)v, &comm));
151249c89c76SBlaise Bourdin PetscCallMPI(MPI_Comm_size(comm, &size));
151349c89c76SBlaise Bourdin PetscCall(VecGetDM(v, &dm));
151449c89c76SBlaise Bourdin PetscCall(DMGetUseNatural(dm, &useNatural));
151549c89c76SBlaise Bourdin useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE;
151649c89c76SBlaise Bourdin if (useNatural) PetscCall(DMPlexCreateNaturalVector(dm, &vNatural));
151749c89c76SBlaise Bourdin else vNatural = v;
151849c89c76SBlaise Bourdin
151949c89c76SBlaise Bourdin /* Read local chunk of the result in the exodus file
152049c89c76SBlaise Bourdin exodus stores each component of a vector-valued field as a separate variable.
152149c89c76SBlaise Bourdin We assume that they are stored sequentially
152249c89c76SBlaise Bourdin Zonal variables are accessed one element block at a time, so we loop through the cell sets,
152349c89c76SBlaise Bourdin but once the vector has been reordered to natural size, we cannot use the label information
152449c89c76SBlaise Bourdin to figure out what to save where. */
15258a1c49cdSMatthew G. Knepley numCS = (PetscExodusIIInt)ex_inquire_int(exoid, EX_INQ_ELEM_BLK); // This is an int64_t
152649c89c76SBlaise Bourdin PetscCall(PetscMalloc2(numCS, &csID, numCS, &csSize));
152749c89c76SBlaise Bourdin PetscCallExternal(ex_get_ids, exoid, EX_ELEM_BLOCK, csID);
152849c89c76SBlaise Bourdin for (set = 0; set < numCS; ++set) {
152949c89c76SBlaise Bourdin ex_block block;
153049c89c76SBlaise Bourdin
153149c89c76SBlaise Bourdin block.id = csID[set];
153249c89c76SBlaise Bourdin block.type = EX_ELEM_BLOCK;
153349c89c76SBlaise Bourdin PetscCallExternal(ex_get_block_param, exoid, &block);
153462e4c2b1SMatthew G. Knepley PetscCall(PetscIntCast(block.num_entry, &csSize[set])); // This is an int64_t
153549c89c76SBlaise Bourdin }
153649c89c76SBlaise Bourdin PetscCall(VecGetOwnershipRange(vNatural, &xs, &xe));
153749c89c76SBlaise Bourdin PetscCall(VecGetBlockSize(vNatural, &bs));
153849c89c76SBlaise Bourdin if (bs > 1) PetscCall(ISCreateStride(comm, (xe - xs) / bs, xs, bs, &compIS));
153949c89c76SBlaise Bourdin for (set = 0; set < numCS; ++set) {
154049c89c76SBlaise Bourdin PetscInt csLocalSize, c;
154149c89c76SBlaise Bourdin
154249c89c76SBlaise Bourdin /* range of indices for set setID[set]: csxs:csxs + csSize[set]-1
154349c89c76SBlaise Bourdin local slice of zonal values: xs/bs,xm/bs-1
154449c89c76SBlaise Bourdin intersection: max(xs/bs,csxs),min(xm/bs-1,csxs + csSize[set]-1) */
154549c89c76SBlaise Bourdin csLocalSize = PetscMax(0, PetscMin(xe / bs, csxs + csSize[set]) - PetscMax(xs / bs, csxs));
154649c89c76SBlaise Bourdin if (bs == 1) {
154749c89c76SBlaise Bourdin PetscCall(VecGetArray(vNatural, &varray));
154849c89c76SBlaise Bourdin PetscCallExternal(ex_get_partial_var, exoid, step, EX_ELEM_BLOCK, offset, csID[set], PetscMax(xs - csxs, 0) + 1, csLocalSize, &varray[PetscMax(0, csxs - xs)]);
154949c89c76SBlaise Bourdin PetscCall(VecRestoreArray(vNatural, &varray));
155049c89c76SBlaise Bourdin } else {
155149c89c76SBlaise Bourdin for (c = 0; c < bs; ++c) {
155249c89c76SBlaise Bourdin PetscCall(ISStrideSetStride(compIS, (xe - xs) / bs, xs + c, bs));
155349c89c76SBlaise Bourdin PetscCall(VecGetSubVector(vNatural, compIS, &vComp));
155449c89c76SBlaise Bourdin PetscCall(VecGetArray(vComp, &varray));
155549c89c76SBlaise Bourdin PetscCallExternal(ex_get_partial_var, exoid, step, EX_ELEM_BLOCK, offset + c, csID[set], PetscMax(xs / bs - csxs, 0) + 1, csLocalSize, &varray[PetscMax(0, csxs - xs / bs)]);
155649c89c76SBlaise Bourdin PetscCall(VecRestoreArray(vComp, &varray));
155749c89c76SBlaise Bourdin PetscCall(VecRestoreSubVector(vNatural, compIS, &vComp));
155849c89c76SBlaise Bourdin }
155949c89c76SBlaise Bourdin }
156049c89c76SBlaise Bourdin csxs += csSize[set];
156149c89c76SBlaise Bourdin }
156249c89c76SBlaise Bourdin PetscCall(PetscFree2(csID, csSize));
156349c89c76SBlaise Bourdin if (bs > 1) PetscCall(ISDestroy(&compIS));
156449c89c76SBlaise Bourdin if (useNatural) {
156549c89c76SBlaise Bourdin PetscCall(DMPlexNaturalToGlobalBegin(dm, vNatural, v));
156649c89c76SBlaise Bourdin PetscCall(DMPlexNaturalToGlobalEnd(dm, vNatural, v));
156749c89c76SBlaise Bourdin PetscCall(VecDestroy(&vNatural));
156849c89c76SBlaise Bourdin }
156949c89c76SBlaise Bourdin PetscFunctionReturn(PETSC_SUCCESS);
157049c89c76SBlaise Bourdin }
157149c89c76SBlaise Bourdin
ExodusGetCellType_Internal(const char * elem_type,DMPolytopeType * ct)157249c89c76SBlaise Bourdin static PetscErrorCode ExodusGetCellType_Internal(const char *elem_type, DMPolytopeType *ct)
157349c89c76SBlaise Bourdin {
157449c89c76SBlaise Bourdin PetscBool flg;
157549c89c76SBlaise Bourdin
157649c89c76SBlaise Bourdin PetscFunctionBegin;
157749c89c76SBlaise Bourdin *ct = DM_POLYTOPE_UNKNOWN;
157849c89c76SBlaise Bourdin PetscCall(PetscStrcmp(elem_type, "BAR2", &flg));
157949c89c76SBlaise Bourdin if (flg) {
158049c89c76SBlaise Bourdin *ct = DM_POLYTOPE_SEGMENT;
158149c89c76SBlaise Bourdin goto done;
158249c89c76SBlaise Bourdin }
158349c89c76SBlaise Bourdin PetscCall(PetscStrcmp(elem_type, "BAR3", &flg));
158449c89c76SBlaise Bourdin if (flg) {
158549c89c76SBlaise Bourdin *ct = DM_POLYTOPE_SEGMENT;
158649c89c76SBlaise Bourdin goto done;
158749c89c76SBlaise Bourdin }
158849c89c76SBlaise Bourdin PetscCall(PetscStrcmp(elem_type, "TRI", &flg));
158949c89c76SBlaise Bourdin if (flg) {
159049c89c76SBlaise Bourdin *ct = DM_POLYTOPE_TRIANGLE;
159149c89c76SBlaise Bourdin goto done;
159249c89c76SBlaise Bourdin }
159349c89c76SBlaise Bourdin PetscCall(PetscStrcmp(elem_type, "TRI3", &flg));
159449c89c76SBlaise Bourdin if (flg) {
159549c89c76SBlaise Bourdin *ct = DM_POLYTOPE_TRIANGLE;
159649c89c76SBlaise Bourdin goto done;
159749c89c76SBlaise Bourdin }
159849c89c76SBlaise Bourdin PetscCall(PetscStrcmp(elem_type, "QUAD", &flg));
159949c89c76SBlaise Bourdin if (flg) {
160049c89c76SBlaise Bourdin *ct = DM_POLYTOPE_QUADRILATERAL;
160149c89c76SBlaise Bourdin goto done;
160249c89c76SBlaise Bourdin }
160349c89c76SBlaise Bourdin PetscCall(PetscStrcmp(elem_type, "QUAD4", &flg));
160449c89c76SBlaise Bourdin if (flg) {
160549c89c76SBlaise Bourdin *ct = DM_POLYTOPE_QUADRILATERAL;
160649c89c76SBlaise Bourdin goto done;
160749c89c76SBlaise Bourdin }
160849c89c76SBlaise Bourdin PetscCall(PetscStrcmp(elem_type, "SHELL4", &flg));
160949c89c76SBlaise Bourdin if (flg) {
161049c89c76SBlaise Bourdin *ct = DM_POLYTOPE_QUADRILATERAL;
161149c89c76SBlaise Bourdin goto done;
161249c89c76SBlaise Bourdin }
161349c89c76SBlaise Bourdin PetscCall(PetscStrcmp(elem_type, "TETRA", &flg));
161449c89c76SBlaise Bourdin if (flg) {
161549c89c76SBlaise Bourdin *ct = DM_POLYTOPE_TETRAHEDRON;
161649c89c76SBlaise Bourdin goto done;
161749c89c76SBlaise Bourdin }
161849c89c76SBlaise Bourdin PetscCall(PetscStrcmp(elem_type, "TET4", &flg));
161949c89c76SBlaise Bourdin if (flg) {
162049c89c76SBlaise Bourdin *ct = DM_POLYTOPE_TETRAHEDRON;
162149c89c76SBlaise Bourdin goto done;
162249c89c76SBlaise Bourdin }
162349c89c76SBlaise Bourdin PetscCall(PetscStrcmp(elem_type, "WEDGE", &flg));
162449c89c76SBlaise Bourdin if (flg) {
162549c89c76SBlaise Bourdin *ct = DM_POLYTOPE_TRI_PRISM;
162649c89c76SBlaise Bourdin goto done;
162749c89c76SBlaise Bourdin }
162849c89c76SBlaise Bourdin PetscCall(PetscStrcmp(elem_type, "HEX", &flg));
162949c89c76SBlaise Bourdin if (flg) {
163049c89c76SBlaise Bourdin *ct = DM_POLYTOPE_HEXAHEDRON;
163149c89c76SBlaise Bourdin goto done;
163249c89c76SBlaise Bourdin }
163349c89c76SBlaise Bourdin PetscCall(PetscStrcmp(elem_type, "HEX8", &flg));
163449c89c76SBlaise Bourdin if (flg) {
163549c89c76SBlaise Bourdin *ct = DM_POLYTOPE_HEXAHEDRON;
163649c89c76SBlaise Bourdin goto done;
163749c89c76SBlaise Bourdin }
163849c89c76SBlaise Bourdin PetscCall(PetscStrcmp(elem_type, "HEXAHEDRON", &flg));
163949c89c76SBlaise Bourdin if (flg) {
164049c89c76SBlaise Bourdin *ct = DM_POLYTOPE_HEXAHEDRON;
164149c89c76SBlaise Bourdin goto done;
164249c89c76SBlaise Bourdin }
164349c89c76SBlaise Bourdin SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unrecognized element type %s", elem_type);
164449c89c76SBlaise Bourdin done:
164549c89c76SBlaise Bourdin PetscFunctionReturn(PETSC_SUCCESS);
164649c89c76SBlaise Bourdin }
164749c89c76SBlaise Bourdin
164849c89c76SBlaise Bourdin /*@
164949c89c76SBlaise Bourdin DMPlexCreateExodus - Create a `DMPLEX` mesh from an ExodusII file ID.
165049c89c76SBlaise Bourdin
165149c89c76SBlaise Bourdin Collective
165249c89c76SBlaise Bourdin
165349c89c76SBlaise Bourdin Input Parameters:
165449c89c76SBlaise Bourdin + comm - The MPI communicator
165549c89c76SBlaise Bourdin . exoid - The ExodusII id associated with a exodus file and obtained using ex_open
165649c89c76SBlaise Bourdin - interpolate - Create faces and edges in the mesh
165749c89c76SBlaise Bourdin
165849c89c76SBlaise Bourdin Output Parameter:
165949c89c76SBlaise Bourdin . dm - The `DM` object representing the mesh
166049c89c76SBlaise Bourdin
166149c89c76SBlaise Bourdin Level: beginner
166249c89c76SBlaise Bourdin
166349c89c76SBlaise Bourdin .seealso: [](ch_unstructured), `DM`, `PETSCVIEWEREXODUSII`, `DMPLEX`, `DMCreate()`
166449c89c76SBlaise Bourdin @*/
DMPlexCreateExodus(MPI_Comm comm,PetscExodusIIInt exoid,PetscBool interpolate,DM * dm)16650a5cf5d8SBlaise Bourdin PetscErrorCode DMPlexCreateExodus(MPI_Comm comm, PetscExodusIIInt exoid, PetscBool interpolate, DM *dm)
166649c89c76SBlaise Bourdin {
166749c89c76SBlaise Bourdin PetscMPIInt num_proc, rank;
166849c89c76SBlaise Bourdin DMLabel cellSets = NULL, faceSets = NULL, vertSets = NULL;
166949c89c76SBlaise Bourdin PetscSection coordSection;
167049c89c76SBlaise Bourdin Vec coordinates;
167149c89c76SBlaise Bourdin PetscScalar *coords;
167249c89c76SBlaise Bourdin PetscInt coordSize, v;
167349c89c76SBlaise Bourdin /* Read from ex_get_init() */
167449c89c76SBlaise Bourdin char title[PETSC_MAX_PATH_LEN + 1];
167549c89c76SBlaise Bourdin int dim = 0, dimEmbed = 0, numVertices = 0, numCells = 0;
167649c89c76SBlaise Bourdin int num_cs = 0, num_vs = 0, num_fs = 0;
167749c89c76SBlaise Bourdin
167849c89c76SBlaise Bourdin PetscFunctionBegin;
167949c89c76SBlaise Bourdin PetscCallMPI(MPI_Comm_rank(comm, &rank));
168049c89c76SBlaise Bourdin PetscCallMPI(MPI_Comm_size(comm, &num_proc));
168149c89c76SBlaise Bourdin PetscCall(DMCreate(comm, dm));
168249c89c76SBlaise Bourdin PetscCall(DMSetType(*dm, DMPLEX));
168349c89c76SBlaise Bourdin /* Open EXODUS II file and read basic information on rank 0, then broadcast to all processors */
168449c89c76SBlaise Bourdin if (rank == 0) {
168549c89c76SBlaise Bourdin PetscCall(PetscMemzero(title, PETSC_MAX_PATH_LEN + 1));
168649c89c76SBlaise Bourdin PetscCallExternal(ex_get_init, exoid, title, &dimEmbed, &numVertices, &numCells, &num_cs, &num_vs, &num_fs);
168749c89c76SBlaise Bourdin PetscCheck(num_cs, PETSC_COMM_SELF, PETSC_ERR_SUP, "Exodus file does not contain any cell set");
168849c89c76SBlaise Bourdin }
168949c89c76SBlaise Bourdin PetscCallMPI(MPI_Bcast(title, PETSC_MAX_PATH_LEN + 1, MPI_CHAR, 0, comm));
169049c89c76SBlaise Bourdin PetscCallMPI(MPI_Bcast(&dim, 1, MPI_INT, 0, comm));
169149c89c76SBlaise Bourdin PetscCall(PetscObjectSetName((PetscObject)*dm, title));
169249c89c76SBlaise Bourdin PetscCall(DMPlexSetChart(*dm, 0, numCells + numVertices));
169349c89c76SBlaise Bourdin /* We do not want this label automatically computed, instead we compute it here */
169449c89c76SBlaise Bourdin PetscCall(DMCreateLabel(*dm, "celltype"));
169549c89c76SBlaise Bourdin
169649c89c76SBlaise Bourdin /* Read cell sets information */
169749c89c76SBlaise Bourdin if (rank == 0) {
169849c89c76SBlaise Bourdin PetscInt *cone;
169949c89c76SBlaise Bourdin int c, cs, ncs, c_loc, v, v_loc;
170049c89c76SBlaise Bourdin /* Read from ex_get_elem_blk_ids() */
170149c89c76SBlaise Bourdin int *cs_id, *cs_order;
170249c89c76SBlaise Bourdin /* Read from ex_get_elem_block() */
170349c89c76SBlaise Bourdin char buffer[PETSC_MAX_PATH_LEN + 1];
170449c89c76SBlaise Bourdin int num_cell_in_set, num_vertex_per_cell, num_hybrid, num_attr;
170549c89c76SBlaise Bourdin /* Read from ex_get_elem_conn() */
170649c89c76SBlaise Bourdin int *cs_connect;
170749c89c76SBlaise Bourdin
170849c89c76SBlaise Bourdin /* Get cell sets IDs */
170949c89c76SBlaise Bourdin PetscCall(PetscMalloc2(num_cs, &cs_id, num_cs, &cs_order));
171049c89c76SBlaise Bourdin PetscCallExternal(ex_get_ids, exoid, EX_ELEM_BLOCK, cs_id);
171149c89c76SBlaise Bourdin /* Read the cell set connectivity table and build mesh topology
171249c89c76SBlaise Bourdin EXO standard requires that cells in cell sets be numbered sequentially and be pairwise disjoint. */
171349c89c76SBlaise Bourdin /* Check for a hybrid mesh */
171449c89c76SBlaise Bourdin for (cs = 0, num_hybrid = 0; cs < num_cs; ++cs) {
171549c89c76SBlaise Bourdin DMPolytopeType ct;
171649c89c76SBlaise Bourdin char elem_type[PETSC_MAX_PATH_LEN];
171749c89c76SBlaise Bourdin
171849c89c76SBlaise Bourdin PetscCall(PetscArrayzero(elem_type, sizeof(elem_type)));
171949c89c76SBlaise Bourdin PetscCallExternal(ex_get_elem_type, exoid, cs_id[cs], elem_type);
172049c89c76SBlaise Bourdin PetscCall(ExodusGetCellType_Internal(elem_type, &ct));
172149c89c76SBlaise Bourdin dim = PetscMax(dim, DMPolytopeTypeGetDim(ct));
172249c89c76SBlaise Bourdin PetscCallExternal(ex_get_block, exoid, EX_ELEM_BLOCK, cs_id[cs], buffer, &num_cell_in_set, &num_vertex_per_cell, 0, 0, &num_attr);
172349c89c76SBlaise Bourdin switch (ct) {
172449c89c76SBlaise Bourdin case DM_POLYTOPE_TRI_PRISM:
172549c89c76SBlaise Bourdin cs_order[cs] = cs;
172649c89c76SBlaise Bourdin ++num_hybrid;
172749c89c76SBlaise Bourdin break;
172849c89c76SBlaise Bourdin default:
172949c89c76SBlaise Bourdin for (c = cs; c > cs - num_hybrid; --c) cs_order[c] = cs_order[c - 1];
173049c89c76SBlaise Bourdin cs_order[cs - num_hybrid] = cs;
173149c89c76SBlaise Bourdin }
173249c89c76SBlaise Bourdin }
173349c89c76SBlaise Bourdin /* First set sizes */
173449c89c76SBlaise Bourdin for (ncs = 0, c = 0; ncs < num_cs; ++ncs) {
173549c89c76SBlaise Bourdin DMPolytopeType ct;
173649c89c76SBlaise Bourdin char elem_type[PETSC_MAX_PATH_LEN];
173749c89c76SBlaise Bourdin const PetscInt cs = cs_order[ncs];
173849c89c76SBlaise Bourdin
173949c89c76SBlaise Bourdin PetscCall(PetscArrayzero(elem_type, sizeof(elem_type)));
174049c89c76SBlaise Bourdin PetscCallExternal(ex_get_elem_type, exoid, cs_id[cs], elem_type);
174149c89c76SBlaise Bourdin PetscCall(ExodusGetCellType_Internal(elem_type, &ct));
174249c89c76SBlaise Bourdin PetscCallExternal(ex_get_block, exoid, EX_ELEM_BLOCK, cs_id[cs], buffer, &num_cell_in_set, &num_vertex_per_cell, 0, 0, &num_attr);
174349c89c76SBlaise Bourdin for (c_loc = 0; c_loc < num_cell_in_set; ++c_loc, ++c) {
174449c89c76SBlaise Bourdin PetscCall(DMPlexSetConeSize(*dm, c, num_vertex_per_cell));
174549c89c76SBlaise Bourdin PetscCall(DMPlexSetCellType(*dm, c, ct));
174649c89c76SBlaise Bourdin }
174749c89c76SBlaise Bourdin }
174849c89c76SBlaise Bourdin for (v = numCells; v < numCells + numVertices; ++v) PetscCall(DMPlexSetCellType(*dm, v, DM_POLYTOPE_POINT));
174949c89c76SBlaise Bourdin PetscCall(DMSetUp(*dm));
175049c89c76SBlaise Bourdin for (ncs = 0, c = 0; ncs < num_cs; ++ncs) {
175149c89c76SBlaise Bourdin const PetscInt cs = cs_order[ncs];
175249c89c76SBlaise Bourdin PetscCallExternal(ex_get_block, exoid, EX_ELEM_BLOCK, cs_id[cs], buffer, &num_cell_in_set, &num_vertex_per_cell, 0, 0, &num_attr);
175349c89c76SBlaise Bourdin PetscCall(PetscMalloc2(num_vertex_per_cell * num_cell_in_set, &cs_connect, num_vertex_per_cell, &cone));
175449c89c76SBlaise Bourdin PetscCallExternal(ex_get_conn, exoid, EX_ELEM_BLOCK, cs_id[cs], cs_connect, NULL, NULL);
175549c89c76SBlaise Bourdin /* EXO uses Fortran-based indexing, DMPlex uses C-style and numbers cell first then vertices. */
175649c89c76SBlaise Bourdin for (c_loc = 0, v = 0; c_loc < num_cell_in_set; ++c_loc, ++c) {
175749c89c76SBlaise Bourdin DMPolytopeType ct;
175849c89c76SBlaise Bourdin
175949c89c76SBlaise Bourdin for (v_loc = 0; v_loc < num_vertex_per_cell; ++v_loc, ++v) cone[v_loc] = cs_connect[v] + numCells - 1;
176049c89c76SBlaise Bourdin PetscCall(DMPlexGetCellType(*dm, c, &ct));
176149c89c76SBlaise Bourdin PetscCall(DMPlexInvertCell(ct, cone));
176249c89c76SBlaise Bourdin PetscCall(DMPlexSetCone(*dm, c, cone));
176349c89c76SBlaise Bourdin PetscCall(DMSetLabelValue_Fast(*dm, &cellSets, "Cell Sets", c, cs_id[cs]));
176449c89c76SBlaise Bourdin }
176549c89c76SBlaise Bourdin PetscCall(PetscFree2(cs_connect, cone));
176649c89c76SBlaise Bourdin }
176749c89c76SBlaise Bourdin PetscCall(PetscFree2(cs_id, cs_order));
176849c89c76SBlaise Bourdin }
176949c89c76SBlaise Bourdin {
177049c89c76SBlaise Bourdin PetscInt ints[] = {dim, dimEmbed};
177149c89c76SBlaise Bourdin
177249c89c76SBlaise Bourdin PetscCallMPI(MPI_Bcast(ints, 2, MPIU_INT, 0, comm));
177349c89c76SBlaise Bourdin PetscCall(DMSetDimension(*dm, ints[0]));
177449c89c76SBlaise Bourdin PetscCall(DMSetCoordinateDim(*dm, ints[1]));
177549c89c76SBlaise Bourdin dim = ints[0];
177649c89c76SBlaise Bourdin dimEmbed = ints[1];
177749c89c76SBlaise Bourdin }
177849c89c76SBlaise Bourdin PetscCall(DMPlexSymmetrize(*dm));
177949c89c76SBlaise Bourdin PetscCall(DMPlexStratify(*dm));
178049c89c76SBlaise Bourdin if (interpolate) {
178149c89c76SBlaise Bourdin DM idm;
178249c89c76SBlaise Bourdin
178349c89c76SBlaise Bourdin PetscCall(DMPlexInterpolate(*dm, &idm));
178449c89c76SBlaise Bourdin PetscCall(DMDestroy(dm));
178549c89c76SBlaise Bourdin *dm = idm;
178649c89c76SBlaise Bourdin }
178749c89c76SBlaise Bourdin
178849c89c76SBlaise Bourdin /* Create vertex set label */
178949c89c76SBlaise Bourdin if (rank == 0 && (num_vs > 0)) {
179049c89c76SBlaise Bourdin int vs, v;
179149c89c76SBlaise Bourdin /* Read from ex_get_node_set_ids() */
179249c89c76SBlaise Bourdin int *vs_id;
179349c89c76SBlaise Bourdin /* Read from ex_get_node_set_param() */
179449c89c76SBlaise Bourdin int num_vertex_in_set;
179549c89c76SBlaise Bourdin /* Read from ex_get_node_set() */
179649c89c76SBlaise Bourdin int *vs_vertex_list;
179749c89c76SBlaise Bourdin
179849c89c76SBlaise Bourdin /* Get vertex set ids */
179949c89c76SBlaise Bourdin PetscCall(PetscMalloc1(num_vs, &vs_id));
180049c89c76SBlaise Bourdin PetscCallExternal(ex_get_ids, exoid, EX_NODE_SET, vs_id);
180149c89c76SBlaise Bourdin for (vs = 0; vs < num_vs; ++vs) {
180249c89c76SBlaise Bourdin PetscCallExternal(ex_get_set_param, exoid, EX_NODE_SET, vs_id[vs], &num_vertex_in_set, NULL);
180349c89c76SBlaise Bourdin PetscCall(PetscMalloc1(num_vertex_in_set, &vs_vertex_list));
180449c89c76SBlaise Bourdin PetscCallExternal(ex_get_set, exoid, EX_NODE_SET, vs_id[vs], vs_vertex_list, NULL);
180549c89c76SBlaise Bourdin for (v = 0; v < num_vertex_in_set; ++v) PetscCall(DMSetLabelValue_Fast(*dm, &vertSets, "Vertex Sets", vs_vertex_list[v] + numCells - 1, vs_id[vs]));
180649c89c76SBlaise Bourdin PetscCall(PetscFree(vs_vertex_list));
180749c89c76SBlaise Bourdin }
180849c89c76SBlaise Bourdin PetscCall(PetscFree(vs_id));
180949c89c76SBlaise Bourdin }
181049c89c76SBlaise Bourdin /* Read coordinates */
181149c89c76SBlaise Bourdin PetscCall(DMGetCoordinateSection(*dm, &coordSection));
181249c89c76SBlaise Bourdin PetscCall(PetscSectionSetNumFields(coordSection, 1));
181349c89c76SBlaise Bourdin PetscCall(PetscSectionSetFieldComponents(coordSection, 0, dimEmbed));
181449c89c76SBlaise Bourdin PetscCall(PetscSectionSetChart(coordSection, numCells, numCells + numVertices));
181549c89c76SBlaise Bourdin for (v = numCells; v < numCells + numVertices; ++v) {
181649c89c76SBlaise Bourdin PetscCall(PetscSectionSetDof(coordSection, v, dimEmbed));
181749c89c76SBlaise Bourdin PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, dimEmbed));
181849c89c76SBlaise Bourdin }
181949c89c76SBlaise Bourdin PetscCall(PetscSectionSetUp(coordSection));
182049c89c76SBlaise Bourdin PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize));
182149c89c76SBlaise Bourdin PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates));
182249c89c76SBlaise Bourdin PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates"));
182349c89c76SBlaise Bourdin PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE));
182449c89c76SBlaise Bourdin PetscCall(VecSetBlockSize(coordinates, dimEmbed));
182549c89c76SBlaise Bourdin PetscCall(VecSetType(coordinates, VECSTANDARD));
182649c89c76SBlaise Bourdin PetscCall(VecGetArray(coordinates, &coords));
182749c89c76SBlaise Bourdin if (rank == 0) {
182849c89c76SBlaise Bourdin PetscReal *x, *y, *z;
182949c89c76SBlaise Bourdin
183049c89c76SBlaise Bourdin PetscCall(PetscMalloc3(numVertices, &x, numVertices, &y, numVertices, &z));
183149c89c76SBlaise Bourdin PetscCallExternal(ex_get_coord, exoid, x, y, z);
183249c89c76SBlaise Bourdin if (dimEmbed > 0) {
183349c89c76SBlaise Bourdin for (v = 0; v < numVertices; ++v) coords[v * dimEmbed + 0] = x[v];
183449c89c76SBlaise Bourdin }
183549c89c76SBlaise Bourdin if (dimEmbed > 1) {
183649c89c76SBlaise Bourdin for (v = 0; v < numVertices; ++v) coords[v * dimEmbed + 1] = y[v];
183749c89c76SBlaise Bourdin }
183849c89c76SBlaise Bourdin if (dimEmbed > 2) {
183949c89c76SBlaise Bourdin for (v = 0; v < numVertices; ++v) coords[v * dimEmbed + 2] = z[v];
184049c89c76SBlaise Bourdin }
184149c89c76SBlaise Bourdin PetscCall(PetscFree3(x, y, z));
184249c89c76SBlaise Bourdin }
184349c89c76SBlaise Bourdin PetscCall(VecRestoreArray(coordinates, &coords));
184449c89c76SBlaise Bourdin PetscCall(DMSetCoordinatesLocal(*dm, coordinates));
184549c89c76SBlaise Bourdin PetscCall(VecDestroy(&coordinates));
184649c89c76SBlaise Bourdin
184749c89c76SBlaise Bourdin /* Create side set label */
184849c89c76SBlaise Bourdin if (rank == 0 && interpolate && (num_fs > 0)) {
184949c89c76SBlaise Bourdin int fs, f, voff;
185049c89c76SBlaise Bourdin /* Read from ex_get_side_set_ids() */
185149c89c76SBlaise Bourdin int *fs_id;
185249c89c76SBlaise Bourdin /* Read from ex_get_side_set_param() */
185349c89c76SBlaise Bourdin int num_side_in_set;
185449c89c76SBlaise Bourdin /* Read from ex_get_side_set_node_list() */
185549c89c76SBlaise Bourdin int *fs_vertex_count_list, *fs_vertex_list, *fs_side_list;
185649c89c76SBlaise Bourdin /* Read side set labels */
185749c89c76SBlaise Bourdin char fs_name[MAX_STR_LENGTH + 1];
185849c89c76SBlaise Bourdin size_t fs_name_len;
185949c89c76SBlaise Bourdin
186049c89c76SBlaise Bourdin /* Get side set ids */
186149c89c76SBlaise Bourdin PetscCall(PetscMalloc1(num_fs, &fs_id));
186249c89c76SBlaise Bourdin PetscCallExternal(ex_get_ids, exoid, EX_SIDE_SET, fs_id);
186349c89c76SBlaise Bourdin // Ids 1 and 2 are reserved by ExodusII for indicating things in 3D
186449c89c76SBlaise Bourdin for (fs = 0; fs < num_fs; ++fs) {
186549c89c76SBlaise Bourdin PetscCallExternal(ex_get_set_param, exoid, EX_SIDE_SET, fs_id[fs], &num_side_in_set, NULL);
186649c89c76SBlaise Bourdin PetscCall(PetscMalloc3(num_side_in_set, &fs_vertex_count_list, num_side_in_set * 4, &fs_vertex_list, num_side_in_set, &fs_side_list));
186749c89c76SBlaise Bourdin PetscCallExternal(ex_get_side_set_node_list, exoid, fs_id[fs], fs_vertex_count_list, fs_vertex_list);
186849c89c76SBlaise Bourdin PetscCallExternal(ex_get_set, exoid, EX_SIDE_SET, fs_id[fs], NULL, fs_side_list);
186949c89c76SBlaise Bourdin
187049c89c76SBlaise Bourdin /* Get the specific name associated with this side set ID. */
187149c89c76SBlaise Bourdin int fs_name_err = ex_get_name(exoid, EX_SIDE_SET, fs_id[fs], fs_name);
187249c89c76SBlaise Bourdin if (!fs_name_err) {
187349c89c76SBlaise Bourdin PetscCall(PetscStrlen(fs_name, &fs_name_len));
187449c89c76SBlaise Bourdin if (fs_name_len == 0) PetscCall(PetscStrncpy(fs_name, "Face Sets", MAX_STR_LENGTH));
187549c89c76SBlaise Bourdin }
187649c89c76SBlaise Bourdin for (f = 0, voff = 0; f < num_side_in_set; ++f) {
187749c89c76SBlaise Bourdin const PetscInt *faces = NULL;
187849c89c76SBlaise Bourdin PetscInt faceSize = fs_vertex_count_list[f], numFaces;
187949c89c76SBlaise Bourdin PetscInt faceVertices[4], v;
188049c89c76SBlaise Bourdin
188149c89c76SBlaise Bourdin PetscCheck(faceSize <= 4, comm, PETSC_ERR_ARG_WRONG, "ExodusII side cannot have %" PetscInt_FMT " > 4 vertices", faceSize);
188249c89c76SBlaise Bourdin for (v = 0; v < faceSize; ++v, ++voff) faceVertices[v] = fs_vertex_list[voff] + numCells - 1;
188349c89c76SBlaise Bourdin PetscCall(DMPlexGetFullJoin(*dm, faceSize, faceVertices, &numFaces, &faces));
188449c89c76SBlaise Bourdin PetscCheck(numFaces == 1, comm, PETSC_ERR_ARG_WRONG, "Invalid ExodusII side %d in set %d maps to %" PetscInt_FMT " faces", f, fs, numFaces);
188549c89c76SBlaise Bourdin PetscCheck(dim == 1 || faces[0] >= numCells + numVertices, comm, PETSC_ERR_ARG_WRONG, "Invalid ExodusII side %d in set %d maps to point %" PetscInt_FMT " which is not a face", f, fs, faces[0]);
188649c89c76SBlaise Bourdin PetscCall(DMSetLabelValue_Fast(*dm, &faceSets, "Face Sets", faces[0], fs_id[fs]));
188749c89c76SBlaise Bourdin /* Only add the label if one has been detected for this side set. */
188849c89c76SBlaise Bourdin if (!fs_name_err) PetscCall(DMSetLabelValue(*dm, fs_name, faces[0], fs_id[fs]));
188949c89c76SBlaise Bourdin PetscCall(DMPlexRestoreJoin(*dm, faceSize, faceVertices, &numFaces, &faces));
189049c89c76SBlaise Bourdin }
189149c89c76SBlaise Bourdin PetscCall(PetscFree3(fs_vertex_count_list, fs_vertex_list, fs_side_list));
189249c89c76SBlaise Bourdin }
189349c89c76SBlaise Bourdin PetscCall(PetscFree(fs_id));
189449c89c76SBlaise Bourdin }
189549c89c76SBlaise Bourdin
189649c89c76SBlaise Bourdin { /* Create Cell/Face/Vertex Sets labels at all processes */
189749c89c76SBlaise Bourdin enum {
189849c89c76SBlaise Bourdin n = 3
189949c89c76SBlaise Bourdin };
190049c89c76SBlaise Bourdin PetscBool flag[n];
190149c89c76SBlaise Bourdin
190249c89c76SBlaise Bourdin flag[0] = cellSets ? PETSC_TRUE : PETSC_FALSE;
190349c89c76SBlaise Bourdin flag[1] = faceSets ? PETSC_TRUE : PETSC_FALSE;
190449c89c76SBlaise Bourdin flag[2] = vertSets ? PETSC_TRUE : PETSC_FALSE;
1905*5440e5dcSBarry Smith PetscCallMPI(MPI_Bcast(flag, n, MPI_C_BOOL, 0, comm));
190649c89c76SBlaise Bourdin if (flag[0]) PetscCall(DMCreateLabel(*dm, "Cell Sets"));
190749c89c76SBlaise Bourdin if (flag[1]) PetscCall(DMCreateLabel(*dm, "Face Sets"));
190849c89c76SBlaise Bourdin if (flag[2]) PetscCall(DMCreateLabel(*dm, "Vertex Sets"));
190949c89c76SBlaise Bourdin }
191049c89c76SBlaise Bourdin PetscFunctionReturn(PETSC_SUCCESS);
191149c89c76SBlaise Bourdin }
1912