v0.12.1
Classes | Typedefs | Functions | Variables
prism_elements_from_surface.cpp File Reference
#include <MoFEM.hpp>

Go to the source code of this file.

Classes

struct  CoordsAndHandle
 
struct  PrismOp
 
struct  PrismFE
 
struct  Op< OP >
 
struct  TriFE
 
struct  QuadFE
 
struct  EdgeFE
 

Typedefs

typedef multi_index_container< CoordsAndHandle, indexed_by< hashed_unique< composite_key< CoordsAndHandle, member< CoordsAndHandle, int, &CoordsAndHandle::x >, member< CoordsAndHandle, int, &CoordsAndHandle::y >, member< CoordsAndHandle, int, &CoordsAndHandle::z > > > > > MapCoords
 

Functions

int main (int argc, char *argv[])
 

Variables

static char help [] = "...\n\n"
 
static constexpr int precision_exponent = 5
 
static constexpr int number_of_prisms_layers = 18
 
static constexpr double delta
 
static constexpr std::array< double, 3 > d3 = {0, 0, 0}
 
static constexpr std::array< double, 3 > d4 = {0, 0, delta}
 

Typedef Documentation

◆ MapCoords

typedef multi_index_container< CoordsAndHandle, indexed_by< hashed_unique<composite_key< CoordsAndHandle, member<CoordsAndHandle, int, &CoordsAndHandle::x>, member<CoordsAndHandle, int, &CoordsAndHandle::y>, member<CoordsAndHandle, int, &CoordsAndHandle::z> > > > > MapCoords

Function Documentation

◆ main()

int main ( int  argc,
char *  argv[] 
)
Examples
prism_elements_from_surface.cpp.

Definition at line 143 of file prism_elements_from_surface.cpp.

143  {
144 
145  MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
146 
147  try {
148 
149  moab::Core mb_instance;
150  moab::Interface &moab = mb_instance;
151  int rank;
152  MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
153 
154  // Read parameters from line command
155  PetscBool flg = PETSC_TRUE;
156  char mesh_file_name[255];
157  CHKERR PetscOptionsGetString(PETSC_NULL, "", "-my_file", mesh_file_name,
158  255, &flg);
159  if (flg != PETSC_TRUE)
160  SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA,
161  "error -my_file (MESH FILE NEEDED)");
162 
163  // Read mesh to MOAB
164  const char *option;
165  option = ""; //"PARALLEL=BCAST;";//;DEBUG_IO";
166  CHKERR moab.load_file(mesh_file_name, 0, option);
167 
168  Range tris;
169  CHKERR moab.get_entities_by_type(0, MBTRI, tris, false);
170  Range tris_verts;
171  CHKERR moab.get_connectivity(tris, tris_verts);
172  MatrixDouble tri_coords(tris_verts.size(), 3);
173  CHKERR moab.get_coords(tris_verts, &tri_coords(0, 0));
174 
175  MoFEM::Core core(moab);
176  MoFEM::Interface &m_field = core;
177 
178  std::array<Range, 3> edge_block;
179  for (auto b : {1, 2, 3})
181  b, BLOCKSET, 1, edge_block[b - 1]);
182 
183  PrismsFromSurfaceInterface *prisms_from_surface_interface;
184  CHKERR m_field.getInterface(prisms_from_surface_interface);
185 
186  Range prisms;
187  CHKERR prisms_from_surface_interface->createPrisms(tris, prisms);
188  prisms_from_surface_interface->setThickness(prisms, d3.data(), d4.data());
189  Range add_prims_layer;
190  Range extrude_prisms = prisms;
191 
192  for (int ll = 1; ll != number_of_prisms_layers; ++ll) {
193  prisms_from_surface_interface->createdVertices.clear();
194  CHKERR prisms_from_surface_interface->createPrismsFromPrisms(
195  extrude_prisms, false, add_prims_layer);
196  prisms_from_surface_interface->setThickness(add_prims_layer, d3.data(),
197  d4.data());
198  extrude_prisms = add_prims_layer;
199  prisms.merge(add_prims_layer);
200  add_prims_layer.clear();
201  }
202 
203  EntityHandle meshset;
204  CHKERR moab.create_meshset(MESHSET_SET | MESHSET_TRACK_OWNER, meshset);
205  CHKERR moab.add_entities(meshset, prisms);
206 
207  MapCoords map_coords;
208  Range verts;
209  CHKERR moab.get_connectivity(prisms, verts);
210  MatrixDouble coords(verts.size(), 3);
211  CHKERR moab.get_coords(verts, &coords(0, 0));
212 
213  for (size_t v = 0; v != verts.size(); ++v)
214  map_coords.insert(CoordsAndHandle(&coords(v, 0), verts[v]));
215 
216  EntityHandle one_prism_meshset;
217  CHKERR moab.create_meshset(MESHSET_SET | MESHSET_TRACK_OWNER,
218  one_prism_meshset);
219 
220  std::array<double, 18> one_prism_coords = {0, 0, 0, 1, 0, 0, 0, 1, 0,
221  0, 0, 1, 1, 0, 1, 0, 1, 1};
222  std::array<EntityHandle, 6> one_prism_nodes;
223  for (int n = 0; n != 6; ++n)
224  CHKERR moab.create_vertex(&one_prism_coords[3 * n], one_prism_nodes[n]);
225  EntityHandle one_prism;
226  CHKERR m_field.get_moab().create_element(MBPRISM, one_prism_nodes.data(), 6,
227  one_prism);
228  Range one_prism_range;
229  one_prism_range.insert(one_prism);
230  CHKERR moab.add_entities(one_prism_meshset, one_prism_range);
231  Range one_prism_verts;
232  CHKERR moab.get_connectivity(one_prism_range, one_prism_verts);
233  CHKERR moab.add_entities(one_prism_meshset, one_prism_verts);
234  Range one_prism_adj_ents;
235  for (int d = 1; d != 3; ++d)
236  CHKERR moab.get_adjacencies(one_prism_range, d, true, one_prism_adj_ents,
237  moab::Interface::UNION);
238  CHKERR moab.add_entities(one_prism_meshset, one_prism_adj_ents);
239 
240  BitRefLevel bit_level0;
241  bit_level0.set(0);
242  CHKERR m_field.getInterface<BitRefManager>()->setEntitiesBitRefLevel(
243  one_prism_range, bit_level0);
244  CHKERR prisms_from_surface_interface->seedPrismsEntities(one_prism_range,
245  bit_level0);
246 
247  // Fields
248  CHKERR m_field.add_field("FIELD1", H1, AINSWORTH_LEGENDRE_BASE, 1);
249  CHKERR m_field.add_ents_to_field_by_type(one_prism_meshset, MBPRISM,
250  "FIELD1", VERY_NOISY);
251 
252  CHKERR m_field.set_field_order(one_prism_meshset, MBVERTEX, "FIELD1", 1);
253  CHKERR m_field.set_field_order(one_prism_meshset, MBEDGE, "FIELD1", 5,
254  VERY_NOISY);
255  CHKERR m_field.set_field_order(one_prism_meshset, MBTRI, "FIELD1", 5);
256  CHKERR m_field.set_field_order(one_prism_meshset, MBQUAD, "FIELD1", 5,
257  VERY_NOISY);
258  CHKERR m_field.set_field_order(one_prism_meshset, MBPRISM, "FIELD1", 7,
259  VERY_NOISY);
260  CHKERR m_field.build_fields(VERY_NOISY);
261 
262  // FE
263  CHKERR m_field.add_finite_element("PRISM");
264  CHKERR m_field.add_finite_element("EDGE");
265  CHKERR m_field.add_finite_element("TRI");
266  CHKERR m_field.add_finite_element("QUAD");
267 
268  // Define rows/cols and element data
269  CHKERR m_field.modify_finite_element_add_field_row("PRISM", "FIELD1");
270  CHKERR m_field.modify_finite_element_add_field_col("PRISM", "FIELD1");
271  CHKERR m_field.modify_finite_element_add_field_data("PRISM", "FIELD1");
272  // Define rows/cols and element data
273  CHKERR m_field.modify_finite_element_add_field_row("EDGE", "FIELD1");
274  CHKERR m_field.modify_finite_element_add_field_col("EDGE", "FIELD1");
275  CHKERR m_field.modify_finite_element_add_field_data("EDGE", "FIELD1");
276  // Define rows/cols and element data
277  CHKERR m_field.modify_finite_element_add_field_row("TRI", "FIELD1");
278  CHKERR m_field.modify_finite_element_add_field_col("TRI", "FIELD1");
279  CHKERR m_field.modify_finite_element_add_field_data("TRI", "FIELD1");
280  // Define rows/cols and element data
281  CHKERR m_field.modify_finite_element_add_field_row("QUAD", "FIELD1");
282  CHKERR m_field.modify_finite_element_add_field_col("QUAD", "FIELD1");
283  CHKERR m_field.modify_finite_element_add_field_data("QUAD", "FIELD1");
284 
285  CHKERR m_field.add_ents_to_finite_element_by_type(one_prism_meshset, MBEDGE,
286  "EDGE");
287  CHKERR m_field.add_ents_to_finite_element_by_type(one_prism_meshset, MBTRI,
288  "TRI");
289  CHKERR m_field.add_ents_to_finite_element_by_type(one_prism_meshset, MBQUAD,
290  "QUAD");
291  CHKERR m_field.add_ents_to_finite_element_by_type(one_prism_meshset,
292  MBPRISM, "PRISM");
293 
294  // build finite elemnts
295  CHKERR m_field.build_finite_elements();
296  // //build adjacencies
297  CHKERR m_field.build_adjacencies(bit_level0);
298 
299  // Problem
300  CHKERR m_field.add_problem("TEST_PROBLEM");
301 
302  // set finite elements for problem
303  CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM", "PRISM");
304  CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM", "EDGE");
305  CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM", "TRI");
306  CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM", "QUAD");
307  // set refinement level for problem
308  CHKERR m_field.modify_problem_ref_level_add_bit("TEST_PROBLEM", bit_level0);
309 
310  // build problem
311  ProblemsManager *prb_mng_ptr;
312  CHKERR m_field.getInterface(prb_mng_ptr);
313  CHKERR prb_mng_ptr->buildProblem("TEST_PROBLEM", true);
314  // partition
315  CHKERR prb_mng_ptr->partitionSimpleProblem("TEST_PROBLEM");
316  CHKERR prb_mng_ptr->partitionFiniteElements("TEST_PROBLEM");
317  // what are ghost nodes, see Petsc Manual
318  CHKERR prb_mng_ptr->partitionGhostDofs("TEST_PROBLEM");
319 
320  PrismFE fe_prism(m_field, tri_coords);
321  fe_prism.getOpPtrVector().push_back(new PrismOp(moab, map_coords));
322  CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "PRISM", fe_prism);
323 
324  EdgeFE fe_edge(m_field, edge_block, one_prism);
325  fe_edge.getOpPtrVector().push_back(
327  moab, map_coords, one_prism));
328  CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "EDGE", fe_edge);
329 
330  TriFE fe_tri(m_field, tri_coords, one_prism);
331  fe_tri.getOpPtrVector().push_back(
333  moab, map_coords, one_prism));
334  CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "TRI", fe_tri);
335 
336  QuadFE fe_quad(m_field, edge_block, one_prism);
337  fe_quad.getOpPtrVector().push_back(
339  moab, map_coords, one_prism));
340  CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "QUAD", fe_quad);
341 
342  CHKERR moab.write_file("prism_mesh.vtk", "VTK", "", &meshset, 1);
343  CHKERR moab.write_file("one_prism_mesh.vtk", "VTK", "", &one_prism_meshset,
344  1);
345  }
346  CATCH_ERRORS;
347 
349  return 0;
350 }
static Index< 'n', 3 > n
@ VERY_NOISY
Definition: definitions.h:287
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:447
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:153
@ H1
continuous field
Definition: definitions.h:178
@ BLOCKSET
Definition: definitions.h:223
@ MOFEM_INVALID_DATA
Definition: definitions.h:129
#define CHKERR
Inline error check.
Definition: definitions.h:610
virtual MoFEMErrorCode modify_finite_element_add_field_row(const std::string &fe_name, const std::string &name_row)=0
set field row which finite element use
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
virtual MoFEMErrorCode add_ents_to_finite_element_by_type(const EntityHandle entities, const EntityType type, const std::string &name, const bool recursive=true)=0
add entities to finite element
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string &name_filed)=0
set finite element field data
virtual MoFEMErrorCode modify_finite_element_add_field_col(const std::string &fe_name, const std::string &name_row)=0
set field col which finite element use
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
virtual MoFEMErrorCode set_field_order(const EntityHandle meshset, const EntityType type, const std::string &name, const ApproximationOrder order, int verb=DEFAULT_VERBOSITY)=0
Set order approximation of the entities in the field.
virtual MoFEMErrorCode add_ents_to_field_by_type(const Range &ents, const EntityType type, const std::string &name, int verb=DEFAULT_VERBOSITY)=0
Add entities to field meshset.
MoFEMErrorCode getEntitiesByDimension(const int ms_id, const unsigned int cubit_bc_type, const int dimension, Range &entities, const bool recursive=true) const
get entities from CUBIT/meshset of a particular entity dimension
MoFEMErrorCode partitionGhostDofs(const std::string name, int verb=VERBOSE)
determine ghost nodes
MoFEMErrorCode partitionSimpleProblem(const std::string name, int verb=VERBOSE)
partition problem dofs
MoFEMErrorCode buildProblem(const std::string name, const bool square_matrix, int verb=VERBOSE)
build problem data structures
MoFEMErrorCode partitionFiniteElements(const std::string name, bool part_from_moab=false, int low_proc=-1, int hi_proc=-1, int verb=VERBOSE)
partition finite elements
virtual MoFEMErrorCode modify_problem_add_finite_element(const std::string &name_problem, const std::string &fe_name)=0
add finite element to problem, this add entities assigned to finite element to a particular problem
virtual MoFEMErrorCode add_problem(const std::string &name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add problem.
virtual MoFEMErrorCode modify_problem_ref_level_add_bit(const std::string &name_problem, const BitRefLevel &bit)=0
add ref level to problem
char mesh_file_name[255]
const Tensor1_Expr< const dTensor0< T, Dim, i >, typename promote< T, double >::V, Dim, i > d(const Tensor0< T * > &a, const Index< i, Dim > index, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
Definition: dTensor0.hpp:27
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:51
ublas::matrix< double, ublas::row_major, DoubleAllocator > MatrixDouble
Definition: Types.hpp:76
CoreTmp< 0 > Core
Definition: Core.hpp:1124
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
DeprecatedCoreInterface Interface
Definition: Interface.hpp:1954
static constexpr int number_of_prisms_layers
static char help[]
static constexpr std::array< double, 3 > d4
multi_index_container< CoordsAndHandle, indexed_by< hashed_unique< composite_key< CoordsAndHandle, member< CoordsAndHandle, int, &CoordsAndHandle::x >, member< CoordsAndHandle, int, &CoordsAndHandle::y >, member< CoordsAndHandle, int, &CoordsAndHandle::z > > > > > MapCoords
static constexpr std::array< double, 3 > d3
Managing BitRefLevels.
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
virtual MoFEMErrorCode add_field(const std::string &name, const FieldSpace space, const FieldApproximationBase base, const FieldCoefficientsNumber nb_of_coefficients, const TagType tag_type=MB_TAG_SPARSE, const enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add field.
virtual moab::Interface & get_moab()=0
Core (interface) class.
Definition: Core.hpp:93
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
Definition: Core.cpp:84
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:124
Deprecated interface functions.
DEPRECATED MoFEMErrorCode loop_finite_elements(const Problem *problem_ptr, const std::string &fe_name, FEMethod &method, int lower_rank, int upper_rank, MoFEMTypes bh, int verb=DEFAULT_VERBOSITY)
Interface for managing meshsets containing materials and boundary conditions.
std::map< EntityHandle, EntityHandle > createdVertices
MoFEMErrorCode createPrisms(const Range &ents, const SwapType swap_type, Range &prisms, int verb=-1)
Make prisms from triangles.
MoFEMErrorCode setThickness(const Range &prisms, const double director3[], const double director4[])
MoFEMErrorCode seedPrismsEntities(Range &prisms, const BitRefLevel &bit, int verb=-1)
Seed prism entities by bit level.
MoFEMErrorCode createPrismsFromPrisms(const Range &prisms, bool from_down, Range &out_prisms, int verb=-1)
Make prisms by extruding top or bottom prisms.
Problem manager is used to build and partition problems.
MoFEMErrorCode getInterface(const MOFEMuuid &uuid, IFACE *&iface) const
Get interface by uuid and return reference to pointer of interface.

Variable Documentation

◆ d3

constexpr std::array<double, 3> d3 = {0, 0, 0}
staticconstexpr

◆ d4

constexpr std::array<double, 3> d4 = {0, 0, delta}
staticconstexpr

◆ delta

constexpr double delta
staticconstexpr

◆ help

char help[] = "...\n\n"
static

◆ number_of_prisms_layers

constexpr int number_of_prisms_layers = 18
staticconstexpr

◆ precision_exponent

constexpr int precision_exponent = 5
staticconstexpr