v0.14.0
Loading...
Searching...
No Matches
Public Member Functions | Public Attributes | List of all members
MyTransport Struct Reference

Application of mix transport data structure. More...

Inheritance diagram for MyTransport:
[legend]
Collaboration diagram for MyTransport:
[legend]

Public Member Functions

 MyTransport (MoFEM::Interface &m_field)
 
MoFEMErrorCode getSource (EntityHandle ent, const double x, const double y, const double z, double &flux)
 set source term More...
 
MoFEMErrorCode getBcOnValues (const EntityHandle ent, const double x, const double y, const double z, double &value)
 
MoFEMErrorCode getBcOnFluxes (const EntityHandle ent, const double x, const double y, const double z, double &flux)
 essential (Neumann) boundary condition (set fluxes) More...
 
 MyTransport (MoFEM::Interface &m_field, BcFluxMap &bc_flux_map)
 
MoFEMErrorCode getSource (EntityHandle ent, const double x, const double y, const double z, double &flux)
 set source term More...
 
MoFEMErrorCode getBcOnValues (const EntityHandle ent, const double x, const double y, const double z, double &value)
 natural (Dirihlet) boundary conditions (set values) More...
 
MoFEMErrorCode getBcOnFluxes (const EntityHandle ent, const double x, const double y, const double z, double &flux)
 essential (Neumann) boundary condition (set fluxes) More...
 
MoFEMErrorCode addBoundaryElements (BitRefLevel &ref_level)
 set-up boundary conditions More...
 
MoFEMErrorCode refineMesh (MixTransportElement &ufe, const int nb_levels, const int order)
 Refine mesh. More...
 
MoFEMErrorCode squashBits ()
 Squash bits of entities. More...
 
MoFEMErrorCode updateMeshsetsFieldsAndElements (const int nb_levels)
 update meshsets with new entities after mesh refinement More...
 
- Public Member Functions inherited from MixTransport::MixTransportElement
 MixTransportElement (MoFEM::Interface &m_field)
 construction of this data structure More...
 
virtual ~MixTransportElement ()
 destructor More...
 
MoFEMErrorCode getDirichletBCIndices (IS *is)
 get dof indices where essential boundary conditions are applied More...
 
virtual MoFEMErrorCode getSource (const EntityHandle ent, const double x, const double y, const double z, double &flux)
 set source term More...
 
virtual MoFEMErrorCode getResistivity (const EntityHandle ent, const double x, const double y, const double z, MatrixDouble3by3 &inv_k)
 natural (Dirichlet) boundary conditions (set values) More...
 
virtual MoFEMErrorCode getBcOnValues (const EntityHandle ent, const int gg, const double x, const double y, const double z, double &value)
 evaluate natural (Dirichlet) boundary conditions More...
 
virtual MoFEMErrorCode getBcOnFluxes (const EntityHandle ent, const double x, const double y, const double z, double &flux)
 essential (Neumann) boundary condition (set fluxes) More...
 
MoFEMErrorCode addFields (const std::string &values, const std::string &fluxes, const int order)
 Add fields to database. More...
 
MoFEMErrorCode addFiniteElements (const std::string &fluxes_name, const std::string &values_name)
 add finite elements More...
 
MoFEMErrorCode buildProblem (BitRefLevel &ref_level)
 Build problem. More...
 
MoFEMErrorCode postProc (const string out_file)
 Post process results. More...
 
MoFEMErrorCode createMatrices ()
 create matrices More...
 
MoFEMErrorCode solveLinearProblem ()
 solve problem More...
 
MoFEMErrorCode calculateResidual ()
 calculate residual More...
 
MoFEMErrorCode evaluateError ()
 Calculate error on elements. More...
 
MoFEMErrorCode destroyMatrices ()
 destroy matrices More...
 

Public Attributes

BcFluxMapbcFluxMap
 
EntityHandle lastEnt
 
double lastFlux
 
- Public Attributes inherited from MixTransport::MixTransportElement
MoFEM::InterfacemField
 
MyVolumeFE feVol
 Instance of volume element. More...
 
MyTriFE feTri
 Instance of surface element. More...
 
VectorDouble valuesAtGaussPts
 values at integration points on element More...
 
MatrixDouble valuesGradientAtGaussPts
 gradients at integration points on element More...
 
VectorDouble divergenceAtGaussPts
 divergence at integration points on element More...
 
MatrixDouble fluxesAtGaussPts
 fluxes at integration points on element More...
 
set< int > bcIndices
 
std::map< int, BlockDatasetOfBlocks
 maps block set id with appropriate BlockData More...
 
Vec D
 
Vec D0
 
Vec F
 
Mat Aij
 
map< double, EntityHandleerrorMap
 
double sumErrorFlux
 
double sumErrorDiv
 
double sumErrorJump
 

Detailed Description

Application of mix transport data structure.

define sources and other stuff

MixTransportElement is a class collecting functions, operators and data for mix implementation of transport element. See there to learn how elements are created or how operators look like.

Some methods in MixTransportElement are abstract, f.e. user need to implement own source therm.

MixTransportElement is a class collecting functions, operators and data for mix implementation of transport element. See there to learn how elements are created or how operators look like.

Some methods in MixTransportElement are abstract, f.e. user need to implement own source therm.

Definition at line 24 of file mix_transport.cpp.

Constructor & Destructor Documentation

◆ MyTransport() [1/2]

MyTransport::MyTransport ( MoFEM::Interface m_field)
inline

Definition at line 26 of file mix_transport.cpp.

26: MixTransportElement(m_field){};

◆ MyTransport() [2/2]

MyTransport::MyTransport ( MoFEM::Interface m_field,
BcFluxMap bc_flux_map 
)
inline

Definition at line 53 of file h_adaptive_transport.cpp.

54 : MixTransportElement(m_field), bcFluxMap(bc_flux_map), lastEnt(0),
55 lastFlux(0) {}
EntityHandle lastEnt
BcFluxMap & bcFluxMap

Member Function Documentation

◆ addBoundaryElements()

MoFEMErrorCode MyTransport::addBoundaryElements ( BitRefLevel ref_level)
inline

set-up boundary conditions

Parameters
ref_levelmesh refinement level
Note
It is assumed that user would like to something non-standard with boundary conditions, have a own type of data structures to pass to functions calculating values and fluxes on boundary. For example BcFluxMap. That way this function is implemented here not in generic class MixTransportElement.
Returns
error code

Definition at line 129 of file h_adaptive_transport.cpp.

129 {
131 Range tets;
132 CHKERR mField.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
133 ref_level, BitRefLevel().set(), MBTET, tets);
134 Skinner skin(&mField.get_moab());
135 Range skin_faces; // skin faces from 3d ents
136 CHKERR skin.find_skin(0, tets, false, skin_faces);
137 // note: what is essential (dirichlet) is natural (neumann) for mix-FE
138 // compared to classical FE
139 Range natural_bc;
141 mField, NODESET | TEMPERATURESET, it)) {
142 Range tris;
143 CHKERR it->getMeshsetIdEntitiesByDimension(mField.get_moab(), 2, tris,
144 true);
145 natural_bc.insert(tris.begin(), tris.end());
146 }
148 mField, SIDESET | HEATFLUXSET, it)) {
149 HeatFluxCubitBcData mydata;
150 CHKERR it->getBcDataStructure(mydata);
151 if (mydata.data.flag1 == 1) {
152 Range tris;
153 CHKERR it->getMeshsetIdEntitiesByDimension(mField.get_moab(), 2, tris,
154 true);
155 bcFluxMap[it->getMeshsetId()].eNts = tris;
156 bcFluxMap[it->getMeshsetId()].fLux = mydata.data.value1;
157 // cerr << bcFluxMap[it->getMeshsetId()].eNts << endl;
158 // cerr << bcFluxMap[it->getMeshsetId()].fLux << endl;
159 }
160 }
161 Range essential_bc = subtract(skin_faces, natural_bc);
162 Range bit_tris;
163 CHKERR mField.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
164 ref_level, BitRefLevel().set(), MBTRI, bit_tris);
165 essential_bc = intersect(bit_tris, essential_bc);
166 natural_bc = intersect(bit_tris, natural_bc);
168 "MIX_BCFLUX");
170 "MIX_BCVALUE");
171 // CHKERR
172 // mField.add_ents_to_finite_element_by_type(skin_faces,MBTRI,"MIX_BCVALUE");
174 }
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
@ TEMPERATURESET
Definition: definitions.h:155
@ HEATFLUXSET
Definition: definitions.h:156
@ NODESET
Definition: definitions.h:146
@ SIDESET
Definition: definitions.h:147
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
#define CHKERR
Inline error check.
Definition: definitions.h:535
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
#define _IT_CUBITMESHSETS_BY_BCDATA_TYPE_FOR_LOOP_(MESHSET_MANAGER, CUBITBCTYPE, IT)
Iterator that loops over a specific Cubit MeshSet in a moFEM field.
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:40
Managing BitRefLevels.
virtual moab::Interface & get_moab()=0
Definition of the heat flux bc data structure.
Definition: BCData.hpp:422
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.

◆ getBcOnFluxes() [1/2]

MoFEMErrorCode MyTransport::getBcOnFluxes ( const EntityHandle  ent,
const double  x,
const double  y,
const double  z,
double flux 
)
inlinevirtual

essential (Neumann) boundary condition (set fluxes)

Parameters
enthandle to finite element entity
xcoord
ycoord
zcoord
fluxreference to flux which is set by function
Returns
[description]

Reimplemented from MixTransport::MixTransportElement.

Definition at line 43 of file mix_transport.cpp.

44 {
46 flux = 0;
48 }
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440

◆ getBcOnFluxes() [2/2]

MoFEMErrorCode MyTransport::getBcOnFluxes ( const EntityHandle  ent,
const double  x,
const double  y,
const double  z,
double flux 
)
inlinevirtual

essential (Neumann) boundary condition (set fluxes)

Parameters
enthandle to finite element entity
xcoord
ycoord
zcoord
fluxreference to flux which is set by function
Returns
[description]

Reimplemented from MixTransport::MixTransportElement.

Definition at line 98 of file h_adaptive_transport.cpp.

99 {
101 if (lastEnt == ent) {
102 flux = lastFlux;
103 } else {
104 flux = 0;
105 for (BcFluxMap::iterator mit = bcFluxMap.begin(); mit != bcFluxMap.end();
106 mit++) {
107 Range &tris = mit->second.eNts;
108 if (tris.find(ent) != tris.end()) {
109 flux = mit->second.fLux;
110 }
111 }
112 lastEnt = ent;
113 lastFlux = flux;
114 }
116 }

◆ getBcOnValues() [1/2]

MoFEMErrorCode MyTransport::getBcOnValues ( const EntityHandle  ent,
const double  x,
const double  y,
const double  z,
double value 
)
inline

Definition at line 36 of file mix_transport.cpp.

37 {
39 value = 1;
41 }

◆ getBcOnValues() [2/2]

MoFEMErrorCode MyTransport::getBcOnValues ( const EntityHandle  ent,
const double  x,
const double  y,
const double  z,
double value 
)
inline

natural (Dirihlet) boundary conditions (set values)

Parameters
enthandle to finite element entity
xcoord
ycoord
zcoord
valuereference to value set by function
Returns
error code

Definition at line 82 of file h_adaptive_transport.cpp.

83 {
85 value = 0;
87 }

◆ getSource() [1/2]

MoFEMErrorCode MyTransport::getSource ( EntityHandle  ent,
const double  x,
const double  y,
const double  z,
double flux 
)
inlinevirtual

set source term

Parameters
enthandle to entity on which function is evaluated
xcoord
ycoord
zcoord
fluxreference to source term set by function
Returns
error code

Reimplemented from MixTransport::MixTransportElement.

Definition at line 28 of file mix_transport.cpp.

29 {
31 // double d = std::sqrt(x*x+y*y+z*z);
32 flux = 1; //-pow(d,5./4.);
34 }

◆ getSource() [2/2]

MoFEMErrorCode MyTransport::getSource ( EntityHandle  ent,
const double  x,
const double  y,
const double  z,
double flux 
)
inlinevirtual

set source term

Parameters
enthandle to entity on which function is evaluated
xcoord
ycoord
zcoord
fluxreference to source term set by function
Returns
error code

Reimplemented from MixTransport::MixTransportElement.

Definition at line 66 of file h_adaptive_transport.cpp.

67 {
69 flux = 0;
71 }

◆ refineMesh()

MoFEMErrorCode MyTransport::refineMesh ( MixTransportElement ufe,
const int  nb_levels,
const int  order 
)
inline

Refine mesh.

Parameters
ufegeneral data structure
nb_levelsnumber of refinement levels
orderset order of approximation
Returns
errpr code

Refinement of could result in distorted mesh, for example, imagine when you have two levels of non-uniform refinement. Some tetrahedra on the mesh at first refinement instance are only refined by splitting subset of edges on it. Then refined child tetrahedra usually will have worse quality than quality of parent element. Refining such element in subsequent mesh refinement, potentially will deteriorate elements quality even worse. To prevent that adding new refinement level, recreate whole hierarchy of meshes.

Note on subsequent improvement could include refinement of tetrahedra from different levels, including initial mesh. So refinement two could split elements created during refinement one and also split elements from an initial mesh.

That adding the new refinement level creates refinement hierarchy of meshes from a scratch, not adding to existing one.

Entities from previous hierarchy are used in that process, but bit levels on those entities are squashed.

Definition at line 203 of file h_adaptive_transport.cpp.

204 {
205 MeshRefinement *refine_ptr;
207 // get refined edges having child vertex
208 auto ref_ents_ptr = mField.get_ref_ents();
209 typedef RefEntity_multiIndex::index<
210 Composite_EntType_and_ParentEntType_mi_tag>::type RefEntsByComposite;
211 const RefEntsByComposite &ref_ents =
213 RefEntsByComposite::iterator rit, hi_rit;
214 rit = ref_ents.lower_bound(boost::make_tuple(MBVERTEX, MBEDGE));
215 hi_rit = ref_ents.upper_bound(boost::make_tuple(MBVERTEX, MBEDGE));
216 Range refined_edges;
217 // thist loop is over vertices which parent is edge
218 for (; rit != hi_rit; rit++) {
219 refined_edges.insert((*rit)->getParentEnt()); // get parent edge
220 }
221 // get tets which has large error
222 Range tets_to_refine;
223 const double max_error = ufe.errorMap.rbegin()->first;
224 // int size = ((double)5/6)*ufe.errorMap.size();
225 for (map<double, EntityHandle>::iterator mit = ufe.errorMap.begin();
226 mit != ufe.errorMap.end(); mit++) {
227 // cerr << mit->first << " " << mit->second << endl;
228 // if((size--)>0) continue;
229 if (mit->first < 0.25 * max_error)
230 continue;
231 tets_to_refine.insert(mit->second);
232 }
233 Range tets_to_refine_edges;
234 CHKERR mField.get_moab().get_adjacencies(
235 tets_to_refine, 1, false, tets_to_refine_edges, moab::Interface::UNION);
236 refined_edges.merge(tets_to_refine_edges);
237 CHKERR mField.getInterface(refine_ptr);
238 for (int ll = 0; ll != nb_levels; ll++) {
239 Range edges;
240 CHKERR mField.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
241 BitRefLevel().set(ll), BitRefLevel().set(), MBEDGE, edges);
242 edges = intersect(edges, refined_edges);
243 // add edges to refine at current level edges (some of the where refined
244 // before)
246 edges, BitRefLevel().set(ll + 1));
247 // get tets at current level
248 Range tets;
249 CHKERR mField.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
250 BitRefLevel().set(ll), BitRefLevel().set(), MBTET, tets);
251 CHKERR refine_ptr->refineTets(tets, BitRefLevel().set(ll + 1));
253 }
254
255 // update fields and elements
256 EntityHandle ref_meshset;
257 CHKERR mField.get_moab().create_meshset(MESHSET_SET, ref_meshset);
258 {
259 // cerr << BitRefLevel().set(nb_levels) << endl;
260 CHKERR mField.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
261 BitRefLevel().set(nb_levels), BitRefLevel().set(), MBTET,
262 ref_meshset);
263
264 Range ref_tets;
265 CHKERR mField.get_moab().get_entities_by_type(ref_meshset, MBTET,
266 ref_tets);
267
268 // add entities to field
269 CHKERR mField.add_ents_to_field_by_type(ref_meshset, MBTET, "FLUXES");
270 CHKERR mField.add_ents_to_field_by_type(ref_meshset, MBTET, "VALUES");
271 CHKERR mField.set_field_order(0, MBTET, "FLUXES", order + 1);
272 CHKERR mField.set_field_order(0, MBTRI, "FLUXES", order + 1);
273 CHKERR mField.set_field_order(0, MBTET, "VALUES", order);
274
275 // add entities to skeleton
276 Range ref_tris;
277 CHKERR mField.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
278 BitRefLevel().set(nb_levels), BitRefLevel().set(), MBTRI, ref_tris);
280 "MIX_SKELETON");
281
282 // add entities to finite elements
285 Mat_Thermal temp_data;
286 CHKERR it->getAttributeDataStructure(temp_data);
287 setOfBlocks[it->getMeshsetId()].cOnductivity =
288 temp_data.data.Conductivity;
289 setOfBlocks[it->getMeshsetId()].cApacity = temp_data.data.HeatCapacity;
290 CHKERR mField.get_moab().get_entities_by_type(
291 it->meshset, MBTET, setOfBlocks[it->getMeshsetId()].tEts, true);
292 setOfBlocks[it->getMeshsetId()].tEts =
293 intersect(ref_tets, setOfBlocks[it->getMeshsetId()].tEts);
295 setOfBlocks[it->getMeshsetId()].tEts, MBTET, "MIX");
296 }
297 }
298 CHKERR mField.get_moab().delete_entities(&ref_meshset, 1);
300 }
@ MAT_THERMALSET
block name is "MAT_THERMAL"
Definition: definitions.h:161
@ BLOCKSET
Definition: definitions.h:148
virtual const RefEntity_multiIndex * get_ref_ents() const =0
Get the ref ents object.
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.
constexpr int nb_levels
Definition: level_set.cpp:41
map< double, EntityHandle > errorMap
std::map< int, BlockData > setOfBlocks
maps block set id with appropriate BlockData
Thermal material data structure.
Mesh refinement interface.
MoFEMErrorCode refineTets(const EntityHandle meshset, const BitRefLevel &bit, int verb=QUIET, const bool debug=false)
refine TET in the meshset
MoFEMErrorCode addVerticesInTheMiddleOfEdges(const EntityHandle meshset, const BitRefLevel &bit, const bool recursive=false, int verb=QUIET, EntityHandle start_v=0)
make vertices in the middle of edges in meshset and add them to refinement levels defined by bit
MoFEMErrorCode updateMeshsetsFieldsAndElements(const int nb_levels)
update meshsets with new entities after mesh refinement

◆ squashBits()

MoFEMErrorCode MyTransport::squashBits ( )
inline

Squash bits of entities.

Information about hierarchy of meshsets is lost, but entities are not deleted from the mesh. After squash entities bits, new hierarchy can be created.

Returns
error code

Definition at line 310 of file h_adaptive_transport.cpp.

310 {
312 BitRefLevel all_but_0;
313 all_but_0.set(0);
314 all_but_0.flip();
315 BitRefLevel garbage_bit;
316 garbage_bit.set(BITREFLEVEL_SIZE - 1); // Garbage level
317 auto ref_ents_ptr = mField.get_ref_ents();
318 RefEntity_multiIndex::iterator mit = ref_ents_ptr->begin();
319 for (; mit != ref_ents_ptr->end(); mit++) {
320 if (mit->get()->getEntType() == MBENTITYSET)
321 continue;
322 BitRefLevel bit = mit->get()->getBitRefLevel();
323 if ((all_but_0 & bit) == bit) {
324 *(const_cast<RefEntity *>(mit->get())->getBitRefLevelPtr()) =
325 garbage_bit;
326 } else {
327 *(const_cast<RefEntity *>(mit->get())->getBitRefLevelPtr()) =
328 BitRefLevel().set(0);
329 }
330 }
332 }
#define BITREFLEVEL_SIZE
max number of refinements
Definition: definitions.h:219
auto bit
set bit
Struct keeps handle to refined handle.

◆ updateMeshsetsFieldsAndElements()

MoFEMErrorCode MyTransport::updateMeshsetsFieldsAndElements ( const int  nb_levels)
inline

update meshsets with new entities after mesh refinement

Parameters
nb_levelsnb_levels
orderappropriate order
Returns
error code

Definition at line 340 of file h_adaptive_transport.cpp.

340 {
341 BitRefLevel ref_level;
343 ref_level.set(nb_levels);
345 EntityHandle meshset = it->meshset;
347 ->updateMeshsetByEntitiesChildren(meshset, ref_level, meshset,
348 MBMAXTYPE, true);
349 }
351 }
#define _IT_CUBITMESHSETS_FOR_LOOP_(MESHSET_MANAGER, IT)
Iterator that loops over all the Cubit MeshSets in a moFEM field.

Member Data Documentation

◆ bcFluxMap

BcFluxMap& MyTransport::bcFluxMap

Definition at line 49 of file h_adaptive_transport.cpp.

◆ lastEnt

EntityHandle MyTransport::lastEnt

Definition at line 50 of file h_adaptive_transport.cpp.

◆ lastFlux

double MyTransport::lastFlux

Definition at line 51 of file h_adaptive_transport.cpp.


The documentation for this struct was generated from the following files: