v0.14.0
mix_transport.cpp
Go to the documentation of this file.
1 
2 
3 #include <MoFEM.hpp>
5 using namespace MoFEM;
6 using namespace MixTransport;
7 
8 namespace bio = boost::iostreams;
9 using bio::stream;
10 using bio::tee_device;
11 
12 static char help[] = "...\n\n";
13 
14 /** define sources and other stuff
15  *
16  * MixTransportElement is a class collecting functions, operators and
17  * data for mix implementation of transport element. See there to
18  * learn how elements are created or how operators look like.
19  *
20  * Some methods in MixTransportElement are abstract, f.e. user need to
21  * implement own source therm.
22  *
23  */
25 
27 
28  MoFEMErrorCode getSource(EntityHandle ent, const double x, const double y,
29  const double z, double &flux) {
31  // double d = std::sqrt(x*x+y*y+z*z);
32  flux = 1; //-pow(d,5./4.);
34  }
35 
36  MoFEMErrorCode getBcOnValues(const EntityHandle ent, const double x,
37  const double y, const double z, double &value) {
39  value = 1;
41  }
42 
43  MoFEMErrorCode getBcOnFluxes(const EntityHandle ent, const double x,
44  const double y, const double z, double &flux) {
46  flux = 0;
48  }
49 };
50 
51 int main(int argc, char *argv[]) {
52 
53  MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
54 
55  try {
56 
57  moab::Core mb_instance;
58  moab::Interface &moab = mb_instance;
59  int rank;
60  MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
61 
62  PetscBool flg = PETSC_TRUE;
63  char mesh_file_name[255];
64  CHKERR PetscOptionsGetString(PETSC_NULL, PETSC_NULL, "-my_file",
65  mesh_file_name, 255, &flg);
66  if (flg != PETSC_TRUE) {
67  SETERRQ(PETSC_COMM_SELF, 1, "*** ERROR -my_file (MESH FILE NEEDED)");
68  }
69 
70  ParallelComm *pcomm = ParallelComm::get_pcomm(&moab, MYPCOMM_INDEX);
71  if (pcomm == NULL)
72  pcomm = new ParallelComm(&moab, PETSC_COMM_WORLD);
73 
74  const char *option;
75  option = ""; //"PARALLEL=BCAST;";//;DEBUG_IO";
77  CHKERR moab.load_file(mesh_file_name, 0, option);
79 
80  // Create MoFEM (Joseph) database
81  MoFEM::Core core(moab);
82  MoFEM::Interface &m_field = core;
83 
84  // set entities bit level
85  BitRefLevel ref_level;
86  ref_level.set(0);
87  CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
88  0, 3, ref_level);
89 
90  // set app. order
91  // see Hierarchic Finite Element Bases on Unstructured Tetrahedral Meshes
92  // (Mark Ainsworth & Joe Coyle)
93  PetscInt order;
94  CHKERR PetscOptionsGetInt(PETSC_NULL, PETSC_NULL, "-my_order", &order,
95  &flg);
96  if (flg != PETSC_TRUE) {
97  order = 2;
98  }
99 
100  // finite elements
101  MyTransport ufe(m_field);
102 
103  CHKERR ufe.addFields("VALUES", "FLUXES", order);
104  CHKERR ufe.addFiniteElements("FLUXES", "VALUES");
105 
106  // Set boundary conditions
107  Range tets;
108  ierr = m_field.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
109  BitRefLevel().set(0), BitRefLevel().set(), MBTET, tets);
110  CHKERRG(ierr);
111  Skinner skin(&moab);
112  Range skin_faces; // skin faces from 3d ents
113  CHKERR skin.find_skin(0, tets, false, skin_faces);
114  CHKERR m_field.add_ents_to_finite_element_by_type(skin_faces, MBTRI,
115  "MIX_BCVALUE");
116 
117  CHKERR ufe.buildProblem(ref_level);
118  CHKERR ufe.createMatrices();
121  CHKERR ufe.evaluateError();
122 
123  double nrm2_F;
124  CHKERR VecNorm(ufe.F, NORM_2, &nrm2_F);
125  // PetscPrintf(PETSC_COMM_WORLD,"nrm2_F = %6.4e\n",nrm2_F);
126  const double eps = 1e-8;
127  if (nrm2_F > eps) {
128  SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
129  "problem with residual");
130  }
131 
132  CHKERR ufe.destroyMatrices();
133 
134  }
135  CATCH_ERRORS;
136 
138 
139  return 0;
140 }
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
BARRIER_PCOMM_RANK_START
#define BARRIER_PCOMM_RANK_START(PCMB)
set barrier start Run code in sequence, starting from process 0, and ends on last process.
Definition: definitions.h:264
MixTransport::MixTransportElement::addFiniteElements
MoFEMErrorCode addFiniteElements(const std::string &fluxes_name, const std::string &values_name)
add finite elements
Definition: MixTransportElement.hpp:213
MYPCOMM_INDEX
#define MYPCOMM_INDEX
default communicator number PCOMM
Definition: definitions.h:228
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
EntityHandle
MixTransport::MixTransportElement::calculateResidual
MoFEMErrorCode calculateResidual()
calculate residual
Definition: MixTransportElement.hpp:617
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
MyTransport::getBcOnFluxes
MoFEMErrorCode getBcOnFluxes(const EntityHandle ent, const double x, const double y, const double z, double &flux)
essential (Neumann) boundary condition (set fluxes)
Definition: mix_transport.cpp:43
MoFEM.hpp
help
static char help[]
Definition: mix_transport.cpp:12
MoFEM::CoreTmp< 0 >::Finalize
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
MyTransport::getBcOnValues
MoFEMErrorCode getBcOnValues(const EntityHandle ent, const double x, const double y, const double z, double &value)
Definition: mix_transport.cpp:36
order
constexpr int order
Definition: dg_projection.cpp:18
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:2010
MixTransport::MixTransportElement::createMatrices
MoFEMErrorCode createMatrices()
create matrices
Definition: MixTransportElement.hpp:479
MoFEM::CoreInterface::add_ents_to_finite_element_by_type
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
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
MixTransport::MixTransportElement::solveLinearProblem
MoFEMErrorCode solveLinearProblem()
solve problem
Definition: MixTransportElement.hpp:493
MixTransport::MixTransportElement::F
Vec F
Definition: MixTransportElement.hpp:475
MixTransport::MixTransportElement
Mix transport problem.
Definition: MixTransportElement.hpp:34
MyTransport
Application of mix transport data structure.
Definition: mix_transport.cpp:24
BARRIER_PCOMM_RANK_END
#define BARRIER_PCOMM_RANK_END(PCMB)
set barrier start Run code in sequence, starting from process 0, and ends on last process.
Definition: definitions.h:286
mesh_file_name
char mesh_file_name[255]
Definition: mesh_smoothing.cpp:22
main
int main(int argc, char *argv[])
Definition: mix_transport.cpp:51
Range
MoFEM::CoreTmp< 0 >::Initialize
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
Definition: Core.cpp:72
MixTransport::MixTransportElement::buildProblem
MoFEMErrorCode buildProblem(BitRefLevel &ref_level)
Build problem.
Definition: MixTransportElement.hpp:296
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
MixTransport
Definition: MaterialUnsaturatedFlow.hpp:11
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1148
eps
static const double eps
Definition: check_base_functions_derivatives_on_tet.cpp:11
MixTransportElement.hpp
Mix implementation of transport element.
MixTransport::MixTransportElement::addFields
MoFEMErrorCode addFields(const std::string &values, const std::string &fluxes, const int order)
Add fields to database.
Definition: MixTransportElement.hpp:194
MoFEM::Exceptions::ierr
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76
MyTransport::MyTransport
MyTransport(MoFEM::Interface &m_field)
Definition: mix_transport.cpp:26
MoFEM::PetscOptionsGetString
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
Definition: DeprecatedPetsc.hpp:172
MoFEM::BitRefManager
Managing BitRefLevels.
Definition: BitRefManager.hpp:21
MoFEM::Types::BitRefLevel
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:40
MOFEM_ATOM_TEST_INVALID
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h:40
MixTransport::MixTransportElement::destroyMatrices
MoFEMErrorCode destroyMatrices()
destroy matrices
Definition: MixTransportElement.hpp:705
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:453
MixTransport::MixTransportElement::evaluateError
MoFEMErrorCode evaluateError()
Calculate error on elements.
Definition: MixTransportElement.hpp:675
MyTransport::getSource
MoFEMErrorCode getSource(EntityHandle ent, const double x, const double y, const double z, double &flux)
set source term
Definition: mix_transport.cpp:28
MoFEM::PetscOptionsGetInt
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
Definition: DeprecatedPetsc.hpp:142
CHKERRG
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:496