v0.13.0
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

using EdgeEle = EdgeElementForcesAndSourcesCore
 
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

◆ EdgeEle

Definition at line 35 of file prism_elements_from_surface.cpp.

◆ 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 144 of file prism_elements_from_surface.cpp.

144  {
145 
146  MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
147 
148  try {
149 
150  moab::Core mb_instance;
151  moab::Interface &moab = mb_instance;
152  int rank;
153  MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
154 
155  // Read parameters from line command
156  PetscBool flg = PETSC_TRUE;
157  char mesh_file_name[255];
158  CHKERR PetscOptionsGetString(PETSC_NULL, "", "-my_file", mesh_file_name,
159  255, &flg);
160  if (flg != PETSC_TRUE)
161  SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA,
162  "error -my_file (MESH FILE NEEDED)");
163 
164  // Read mesh to MOAB
165  const char *option;
166  option = ""; //"PARALLEL=BCAST;";//;DEBUG_IO";
167  CHKERR moab.load_file(mesh_file_name, 0, option);
168 
169  Range tris;
170  CHKERR moab.get_entities_by_type(0, MBTRI, tris, false);
171  Range tris_verts;
172  CHKERR moab.get_connectivity(tris, tris_verts);
173  MatrixDouble tri_coords(tris_verts.size(), 3);
174  CHKERR moab.get_coords(tris_verts, &tri_coords(0, 0));
175 
176  MoFEM::Core core(moab);
177  MoFEM::Interface &m_field = core;
178 
179  std::array<Range, 3> edge_block;
180  for (auto b : {1, 2, 3})
182  b, BLOCKSET, 1, edge_block[b - 1]);
183 
184  PrismsFromSurfaceInterface *prisms_from_surface_interface;
185  CHKERR m_field.getInterface(prisms_from_surface_interface);
186 
187  Range prisms;
188  CHKERR prisms_from_surface_interface->createPrisms(tris, prisms);
189  prisms_from_surface_interface->setThickness(prisms, d3.data(), d4.data());
190  Range add_prims_layer;
191  Range extrude_prisms = prisms;
192 
193  for (int ll = 1; ll != number_of_prisms_layers; ++ll) {
194  prisms_from_surface_interface->createdVertices.clear();
195  CHKERR prisms_from_surface_interface->createPrismsFromPrisms(
196  extrude_prisms, false, add_prims_layer);
197  prisms_from_surface_interface->setThickness(add_prims_layer, d3.data(),
198  d4.data());
199  extrude_prisms = add_prims_layer;
200  prisms.merge(add_prims_layer);
201  add_prims_layer.clear();
202  }
203 
204  EntityHandle meshset;
205  CHKERR moab.create_meshset(MESHSET_SET | MESHSET_TRACK_OWNER, meshset);
206  CHKERR moab.add_entities(meshset, prisms);
207 
208  MapCoords map_coords;
209  Range verts;
210  CHKERR moab.get_connectivity(prisms, verts);
211  MatrixDouble coords(verts.size(), 3);
212  CHKERR moab.get_coords(verts, &coords(0, 0));
213 
214  for (size_t v = 0; v != verts.size(); ++v)
215  map_coords.insert(CoordsAndHandle(&coords(v, 0), verts[v]));
216 
217  EntityHandle one_prism_meshset;
218  CHKERR moab.create_meshset(MESHSET_SET | MESHSET_TRACK_OWNER,
219  one_prism_meshset);
220 
221  std::array<double, 18> one_prism_coords = {0, 0, 0, 1, 0, 0, 0, 1, 0,
222  0, 0, 1, 1, 0, 1, 0, 1, 1};
223  std::array<EntityHandle, 6> one_prism_nodes;
224  for (int n = 0; n != 6; ++n)
225  CHKERR moab.create_vertex(&one_prism_coords[3 * n], one_prism_nodes[n]);
226  EntityHandle one_prism;
227  CHKERR m_field.get_moab().create_element(MBPRISM, one_prism_nodes.data(), 6,
228  one_prism);
229  Range one_prism_range;
230  one_prism_range.insert(one_prism);
231  CHKERR moab.add_entities(one_prism_meshset, one_prism_range);
232  Range one_prism_verts;
233  CHKERR moab.get_connectivity(one_prism_range, one_prism_verts);
234  CHKERR moab.add_entities(one_prism_meshset, one_prism_verts);
235  Range one_prism_adj_ents;
236  for (int d = 1; d != 3; ++d)
237  CHKERR moab.get_adjacencies(one_prism_range, d, true, one_prism_adj_ents,
238  moab::Interface::UNION);
239  CHKERR moab.add_entities(one_prism_meshset, one_prism_adj_ents);
240 
241  BitRefLevel bit_level0;
242  bit_level0.set(0);
243  CHKERR m_field.getInterface<BitRefManager>()->setEntitiesBitRefLevel(
244  one_prism_range, bit_level0);
245  CHKERR prisms_from_surface_interface->seedPrismsEntities(one_prism_range,
246  bit_level0);
247 
248  // Fields
249  CHKERR m_field.add_field("FIELD1", H1, AINSWORTH_LEGENDRE_BASE, 1);
250  CHKERR m_field.add_ents_to_field_by_type(one_prism_meshset, MBPRISM,
251  "FIELD1", VERY_NOISY);
252 
253  CHKERR m_field.set_field_order(one_prism_meshset, MBVERTEX, "FIELD1", 1);
254  CHKERR m_field.set_field_order(one_prism_meshset, MBEDGE, "FIELD1", 5,
255  VERY_NOISY);
256  CHKERR m_field.set_field_order(one_prism_meshset, MBTRI, "FIELD1", 5);
257  CHKERR m_field.set_field_order(one_prism_meshset, MBQUAD, "FIELD1", 5,
258  VERY_NOISY);
259  CHKERR m_field.set_field_order(one_prism_meshset, MBPRISM, "FIELD1", 7,
260  VERY_NOISY);
261  CHKERR m_field.build_fields(VERY_NOISY);
262 
263  // FE
264  CHKERR m_field.add_finite_element("PRISM");
265  CHKERR m_field.add_finite_element("EDGE");
266  CHKERR m_field.add_finite_element("TRI");
267  CHKERR m_field.add_finite_element("QUAD");
268 
269  // Define rows/cols and element data
270  CHKERR m_field.modify_finite_element_add_field_row("PRISM", "FIELD1");
271  CHKERR m_field.modify_finite_element_add_field_col("PRISM", "FIELD1");
272  CHKERR m_field.modify_finite_element_add_field_data("PRISM", "FIELD1");
273  // Define rows/cols and element data
274  CHKERR m_field.modify_finite_element_add_field_row("EDGE", "FIELD1");
275  CHKERR m_field.modify_finite_element_add_field_col("EDGE", "FIELD1");
276  CHKERR m_field.modify_finite_element_add_field_data("EDGE", "FIELD1");
277  // Define rows/cols and element data
278  CHKERR m_field.modify_finite_element_add_field_row("TRI", "FIELD1");
279  CHKERR m_field.modify_finite_element_add_field_col("TRI", "FIELD1");
280  CHKERR m_field.modify_finite_element_add_field_data("TRI", "FIELD1");
281  // Define rows/cols and element data
282  CHKERR m_field.modify_finite_element_add_field_row("QUAD", "FIELD1");
283  CHKERR m_field.modify_finite_element_add_field_col("QUAD", "FIELD1");
284  CHKERR m_field.modify_finite_element_add_field_data("QUAD", "FIELD1");
285 
286  CHKERR m_field.add_ents_to_finite_element_by_type(one_prism_meshset, MBEDGE,
287  "EDGE");
288  CHKERR m_field.add_ents_to_finite_element_by_type(one_prism_meshset, MBTRI,
289  "TRI");
290  CHKERR m_field.add_ents_to_finite_element_by_type(one_prism_meshset, MBQUAD,
291  "QUAD");
292  CHKERR m_field.add_ents_to_finite_element_by_type(one_prism_meshset,
293  MBPRISM, "PRISM");
294 
295  // build finite elemnts
296  CHKERR m_field.build_finite_elements();
297  // //build adjacencies
298  CHKERR m_field.build_adjacencies(bit_level0);
299 
300  // Problem
301  CHKERR m_field.add_problem("TEST_PROBLEM");
302 
303  // set finite elements for problem
304  CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM", "PRISM");
305  CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM", "EDGE");
306  CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM", "TRI");
307  CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM", "QUAD");
308  // set refinement level for problem
309  CHKERR m_field.modify_problem_ref_level_add_bit("TEST_PROBLEM", bit_level0);
310 
311  // build problem
312  ProblemsManager *prb_mng_ptr;
313  CHKERR m_field.getInterface(prb_mng_ptr);
314  CHKERR prb_mng_ptr->buildProblem("TEST_PROBLEM", true);
315  // partition
316  CHKERR prb_mng_ptr->partitionSimpleProblem("TEST_PROBLEM");
317  CHKERR prb_mng_ptr->partitionFiniteElements("TEST_PROBLEM");
318  // what are ghost nodes, see Petsc Manual
319  CHKERR prb_mng_ptr->partitionGhostDofs("TEST_PROBLEM");
320 
321  PrismFE fe_prism(m_field, tri_coords);
322  fe_prism.getOpPtrVector().push_back(new PrismOp(moab, map_coords));
323  CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "PRISM", fe_prism);
324 
325  EdgeFE fe_edge(m_field, edge_block, one_prism);
326  fe_edge.getOpPtrVector().push_back(
328  moab, map_coords, one_prism));
329  CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "EDGE", fe_edge);
330 
331  TriFE fe_tri(m_field, tri_coords, one_prism);
332  fe_tri.getOpPtrVector().push_back(
334  moab, map_coords, one_prism));
335  CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "TRI", fe_tri);
336 
337  QuadFE fe_quad(m_field, edge_block, one_prism);
338  fe_quad.getOpPtrVector().push_back(
340  moab, map_coords, one_prism));
341  CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "QUAD", fe_quad);
342 
343  CHKERR moab.write_file("prism_mesh.vtk", "VTK", "", &meshset, 1);
344  CHKERR moab.write_file("one_prism_mesh.vtk", "VTK", "", &one_prism_meshset,
345  1);
346  }
347  CATCH_ERRORS;
348 
350  return 0;
351 }
static Index< 'n', 3 > n
@ VERY_NOISY
Definition: definitions.h:225
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:73
@ H1
continuous field
Definition: definitions.h:98
@ BLOCKSET
Definition: definitions.h:161
@ MOFEM_INVALID_DATA
Definition: definitions.h:49
#define CHKERR
Inline error check.
Definition: definitions.h:548
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.
virtual MoFEMErrorCode loop_finite_elements(const std::string &problem_name, const std::string &fe_name, FEMethod &method, boost::shared_ptr< NumeredEntFiniteElement_multiIndex > fe_ptr=nullptr, MoFEMTypes bh=MF_EXIST, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr(), int verb=DEFAULT_VERBOSITY)=0
Make a loop over finite elements.
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
double v
phase velocity of light in medium (cm/ns)
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
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:88
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:51
CoreTmp< 0 > Core
Definition: Core.hpp:1096
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
DeprecatedCoreInterface Interface
Definition: Interface.hpp:1965
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:92
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
Definition: Core.cpp:85
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:125
Deprecated interface functions.
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(IFACE *&iface) const
Get interface refernce 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