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