v0.14.0
Loading...
Searching...
No Matches
Classes | Functions | Variables
mix_transport.cpp File Reference
#include <BasicFiniteElements.hpp>
#include <MixTransportElement.hpp>

Go to the source code of this file.

Classes

struct  MyTransport
 Application of mix transport data structure. More...
 

Functions

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

Variables

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

Function Documentation

◆ main()

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

Definition at line 51 of file mix_transport.cpp.

51 {
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();
119 CHKERR ufe.solveLinearProblem();
120 CHKERR ufe.calculateResidual();
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 }
136
138
139 return 0;
140}
static const double eps
#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:273
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:372
#define MYPCOMM_INDEX
default communicator number PCOMM
Definition: definitions.h:215
#define CHKERRG(n)
Check error code of MoFEM/MOAB/PETSc function.
Definition: definitions.h:483
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h:40
#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:251
#define CHKERR
Inline error check.
Definition: definitions.h:535
constexpr int order
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
char mesh_file_name[255]
static char help[]
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:40
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
Managing BitRefLevels.
Core (interface) class.
Definition: Core.hpp:82
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
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
Deprecated interface functions.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
Application of mix transport data structure.

Variable Documentation

◆ help

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

Definition at line 12 of file mix_transport.cpp.