Actual source code: ex1.cxx
petsc-3.7.3 2016-08-01
1: static char help[] = "Simple MOAB example\n\n";
3: #include <petscdmmoab.h>
4: #include "moab/ScdInterface.hpp"
6: typedef struct {
7: DM dm; /* DM implementation using the MOAB interface */
8: PetscLogEvent createMeshEvent;
9: /* Domain and mesh definition */
10: PetscInt dim;
11: char filename[PETSC_MAX_PATH_LEN];
12: char tagname[PETSC_MAX_PATH_LEN];
13: } AppCtx;
17: PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
18: {
20: PetscBool flg;
23: PetscStrcpy(options->filename, "");
24: PetscStrcpy(options->tagname, "petsc_tag");
25: options->dim = -1;
27: PetscOptionsBegin(comm, "", "MOAB example options", "DMMOAB");
28: PetscOptionsInt("-dim", "The topological mesh dimension", "ex1.c", options->dim, &options->dim, NULL);
29: PetscOptionsString("-filename", "The file containing the mesh", "ex1.c", options->filename, options->filename, sizeof(options->filename), NULL);
30: PetscOptionsString("-tagname", "The tag name from which to create a vector", "ex1.c", options->tagname, options->tagname, sizeof(options->tagname), &flg);
31: PetscOptionsEnd();
33: PetscLogEventRegister("CreateMesh", DM_CLASSID, &options->createMeshEvent);
34: return(0);
35: }
39: PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
40: {
41: moab::Interface *iface=NULL;
42: moab::ParallelComm *pcomm=NULL;
43: moab::Tag tag=NULL;
44: moab::Tag ltog_tag=NULL;
45: moab::Range range;
46: PetscInt tagsize;
47: moab::ErrorCode merr;
51: PetscLogEventBegin(user->createMeshEvent,0,0,0,0);
52: DMMoabCreateMoab(comm, iface, pcomm, <og_tag, &range, dm);
53: std::cout << "Created DMMoab using DMMoabCreateMoab." << std::endl;
54: DMMoabGetInterface(*dm, &iface);
56: // load file and get entities of requested or max dimension
57: if (strlen(user->filename) > 0) {
58: merr = iface->load_file(user->filename);MBERRNM(merr);
59: std::cout << "Read mesh from file " << user->filename << std::endl;
60: }
61: else {
62: // make a simple structured mesh
63: moab::ScdInterface *scdi;
64: merr = iface->query_interface(scdi);
66: moab::ScdBox *box;
67: merr = scdi->construct_box (moab::HomCoord(0,0,0), moab::HomCoord(5,5,5), NULL, 0, box);MBERRNM(merr);
68: user->dim = 3;
69: merr = iface->set_dimension(user->dim);MBERRNM(merr);
70: std::cout << "Created structured 5x5x5 mesh." << std::endl;
71: }
72: if (-1 == user->dim) {
73: moab::Range tmp_range;
74: merr = iface->get_entities_by_handle(0, tmp_range);MBERRNM(merr);
75: if (tmp_range.empty()) {
76: MBERRNM(moab::MB_FAILURE);
77: }
78: user->dim = iface->dimension_from_handle(*tmp_range.rbegin());
79: }
80: merr = iface->get_entities_by_dimension(0, user->dim, range);MBERRNM(merr);
81: DMMoabSetLocalVertices(*dm, &range);
83: // get the requested tag and create if necessary
84: std::cout << "Creating tag with name: " << user->tagname << ";\n";
85: merr = iface->tag_get_handle(user->tagname, 1, moab::MB_TYPE_DOUBLE, tag, moab::MB_TAG_CREAT | moab::MB_TAG_DENSE);MBERRNM(merr);
86: {
87: // initialize new tag with gids
88: std::vector<double> tag_vals(range.size());
89: merr = iface->tag_get_data(ltog_tag, range, tag_vals.data());MBERRNM(merr); // read them into ints
90: double *dval = tag_vals.data(); int *ival = reinterpret_cast<int*>(dval); // "stretch" them into doubles, from the end
91: for (int i = tag_vals.size()-1; i >= 0; i--) dval[i] = ival[i];
92: merr = iface->tag_set_data(tag, range, tag_vals.data());MBERRNM(merr); // write them into doubles
93: }
94: merr = iface->tag_get_length(tag, tagsize);MBERRNM(merr);
96: DMSetUp(*dm);
98: // create the dmmoab and initialize its data
99: PetscObjectSetName((PetscObject) *dm, "MOAB mesh");
100: PetscLogEventEnd(user->createMeshEvent,0,0,0,0);
101: user->dm = *dm;
102: return(0);
103: }
107: int main(int argc, char **argv)
108: {
109: AppCtx user; /* user-defined work context */
110: moab::ErrorCode merr;
112: Vec vec;
113: moab::Interface* mbImpl=NULL;
114: moab::Tag datatag=NULL;
116: PetscInitialize(&argc, &argv, NULL, help);
117: ProcessOptions(PETSC_COMM_WORLD, &user);
119: CreateMesh(PETSC_COMM_WORLD, &user, &user.dm); /* create the MOAB dm and the mesh */
121: DMMoabGetInterface(user.dm, &mbImpl);
122: merr = mbImpl->tag_get_handle(user.tagname, datatag);MBERRNM(merr);
123: DMMoabCreateVector(user.dm, datatag, NULL, PETSC_TRUE, PETSC_FALSE,
124: &vec); /* create a vec from user-input tag */
126: std::cout << "Created VecMoab from existing tag." << std::endl;
127: VecDestroy(&vec);
128: std::cout << "Destroyed VecMoab." << std::endl;
129: DMDestroy(&user.dm);
130: std::cout << "Destroyed DMMoab." << std::endl;
131: PetscFinalize();
132: return 0;
133: }