|
| v0.14.0
|
Go to the documentation of this file.
5 namespace bio = boost::iostreams;
11 static char help[] =
"...\n\n";
13 static const double face_coords[4][9] = {{0, 0, 0, 1, 0, 0, 0, 0, 1},
15 {1, 0, 0, 0, 1, 0, 0, 0, 1},
17 {0, 0, 0, 0, 1, 0, 0, 0, 1},
19 {0, 0, 0, 1, 0, 0, 0, 1, 0}};
21 int main(
int argc,
char *argv[]) {
27 enum bases { HDIV_AINSWORTH, HDIV_DEMKOWICZ,
LASTOP };
29 const char *list[] = {
"hdiv_ainsworth",
"hdiv_demkowicz"};
32 PetscInt choise_value = HDIV_AINSWORTH;
35 if (flg != PETSC_TRUE) {
42 MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
45 #if PETSC_VERSION_GE(3, 6, 4)
52 if (flg != PETSC_TRUE) {
53 SETERRQ(PETSC_COMM_SELF, 1,
"*** ERROR -my_file (MESH FILE NEEDED)");
66 CHKERR moab.create_meshset(MESHSET_SET, meshset_level0);
73 switch (choise_value) {
92 "MESH_NODE_POSITIONS");
98 "MESH_NODE_POSITIONS");
104 "MESH_NODE_POSITIONS");
131 CHKERR skin.find_skin(0, tets,
false, skin_faces);
138 faces = subtract(faces, skin_faces);
177 CHKERR VecSetRandom(
v, PETSC_NULL);
179 "TEST_PROBLEM",
ROW,
v, INSERT_VALUES, SCATTER_REVERSE);
182 typedef tee_device<std::ostream, std::ofstream>
TeeDevice;
184 std::ofstream ofs(
"forces_and_sources_hdiv_continuity_check.txt");
197 m_field(m_field), tH(
th) {}
208 boost::shared_ptr<const NumeredEntFiniteElement> mofem_fe =
209 getNumeredEntFiniteElementPtr();
212 .find(boost::make_tuple(
type, side))
215 int sense = side_table.get<1>()
216 .find(boost::make_tuple(
type, side))
221 int nb_dofs = data.
getN().size2() / 3;
222 for (
int dd = 0;
dd != nb_dofs;
dd++) {
223 for (
int ddd = 0; ddd != 3; ddd++)
230 (
const void **)&t_ptr);
231 for (
int dd = 0;
dd < 3;
dd++)
232 t_ptr[
dd] += sense *
t[
dd];
243 int getRule(
int order) {
return -1; };
252 gaussPts.resize(4, 4);
253 for (
int ff = 0; ff < 4; ff++) {
254 for (
int dd = 0;
dd < 3;
dd++) {
265 struct OpFacesSkinFluxes
276 m_field(m_field), tH1(th1), tH2(th2), mySplit(my_split) {}
285 EntityHandle face = getNumeredEntFiniteElementPtr()->getEnt();
289 (
const void **)&t_ptr);
292 (
const void **)&tn_ptr);
294 *tn_ptr = getNormalsAtGaussPts()(0, 0) * t_ptr[0] +
295 getNormalsAtGaussPts()(0, 1) * t_ptr[1] +
296 getNormalsAtGaussPts()(0, 2) * t_ptr[2];
298 int nb_dofs = data.
getN().size2() / 3;
300 for (;
dd < nb_dofs;
dd++) {
301 *tn_ptr += -getNormalsAtGaussPts()(0, 0) *
303 getNormalsAtGaussPts()(0, 1) * data.
getN()(0, 3 *
dd + 1) *
305 getNormalsAtGaussPts()(0, 2) * data.
getN()(0, 3 *
dd + 2) *
309 const double eps = 1e-8;
310 if (fabs(*tn_ptr) >
eps) {
312 "HDiv continuity failed %6.4e", *tn_ptr);
315 mySplit.precision(5);
317 mySplit << face <<
" " << fabs(*tn_ptr) << std::endl;
334 m_field(m_field), tH1(th1), tH2(th2), mySplit(my_split) {}
343 EntityHandle face = getNumeredEntFiniteElementPtr()->getEnt();
347 (
const void **)&t_ptr);
350 (
const void **)&tn_ptr);
352 *tn_ptr = getNormalsAtGaussPts()(0, 0) * t_ptr[0] +
353 getNormalsAtGaussPts()(0, 1) * t_ptr[1] +
354 getNormalsAtGaussPts()(0, 2) * t_ptr[2];
356 const double eps = 1e-8;
357 if (fabs(*tn_ptr) >
eps) {
359 "HDiv continuity failed %6.4e", *tn_ptr);
362 mySplit.precision(5);
364 mySplit << face <<
" " << fabs(*tn_ptr) << std::endl;
374 int getRule(
int order) {
return -1; };
379 gaussPts.resize(3, 1);
388 MyTetFE tet_fe(m_field);
389 MyTriFE tri_fe(m_field);
390 MyTriFE skin_fe(m_field);
393 double def_val[] = {0, 0, 0};
394 CHKERR moab.tag_get_handle(
"T", 3, MB_TYPE_DOUBLE, th1,
395 MB_TAG_CREAT | MB_TAG_SPARSE, &def_val);
397 auto material_grad_mat = boost::make_shared<MatrixDouble>();
398 auto material_det_vec = boost::make_shared<VectorDouble>();
399 auto material_inv_grad_mat = boost::make_shared<MatrixDouble>();
402 "MESH_NODE_POSITIONS", material_grad_mat));
404 material_grad_mat, material_det_vec, material_inv_grad_mat));
405 tet_fe.getOpPtrVector().push_back(
new OpSetHOWeights(material_det_vec));
407 HDIV, material_det_vec, material_grad_mat));
408 tet_fe.getOpPtrVector().push_back(
410 tet_fe.getOpPtrVector().push_back(
new OpTetFluxes(m_field, th1));
413 CHKERR moab.tag_get_handle(
"TN", 1, MB_TYPE_DOUBLE, th2,
414 MB_TAG_CREAT | MB_TAG_SPARSE, &def_val);
415 tri_fe.getOpPtrVector().push_back(
417 tri_fe.getOpPtrVector().push_back(
419 tri_fe.getOpPtrVector().push_back(
422 skin_fe.getOpPtrVector().push_back(
424 skin_fe.getOpPtrVector().push_back(
426 skin_fe.getOpPtrVector().push_back(
427 new OpFacesSkinFluxes(m_field, th1, th2, my_split));
429 for (Range::iterator fit = faces.begin(); fit != faces.end(); fit++) {
430 CHKERR moab.tag_set_data(th1, &*fit, 1, &def_val);
431 CHKERR moab.tag_set_data(th2, &*fit, 1, &def_val);
435 my_split <<
"internal\n";
437 my_split <<
"skin\n";
441 CHKERR moab.create_meshset(MESHSET_SET, meshset);
444 CHKERR moab.write_file(
"out.vtk",
"VTK",
"", &meshset, 1);
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Data on single entity (This is passed as argument to DataOperator::doWork)
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.
virtual MoFEMErrorCode loop_dofs(const Problem *problem_ptr, const std::string &field_name, RowColData rc, DofMethod &method, int lower_rank, int upper_rank, int verb=DEFAULT_VERBOSITY)=0
Make a loop over dofs.
MoFEMErrorCode buildProblem(const std::string name, const bool square_matrix, int verb=VERBOSE)
build problem data structures
Problem manager is used to build and partition problems.
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
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
UBlasMatrix< double > MatrixDouble
const VectorDouble & getFieldData() const
get dofs values
Projection of edge entities with one mid-node on hierarchical basis.
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
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.
Deprecated interface functions.
DeprecatedCoreInterface Interface
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 CHKERR
Inline error check.
friend class UserDataOperator
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
virtual moab::Interface & get_moab()=0
implementation of Data Operators for Forces and Sources
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_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
default operator for TRI element
multi_index_container< boost::shared_ptr< SideNumber >, indexed_by< ordered_unique< member< SideNumber, EntityHandle, &SideNumber::ent > >, ordered_non_unique< composite_key< SideNumber, const_mem_fun< SideNumber, EntityType, &SideNumber::getEntType >, member< SideNumber, signed char, &SideNumber::side_number > > > > > SideNumber_multiIndex
SideNumber_multiIndex for SideNumber.
friend class UserDataOperator
Calculate normals at Gauss points of triangle element.
constexpr double t
plate stiffness
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_field)=0
set finite element field data
Volume finite element base.
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_ref_level_add_bit(const std::string &name_problem, const BitRefLevel &bit)=0
add ref level to problem
Vector manager is used to create vectors \mofem_vectors.
Get field gradients at integration pts for scalar filed rank 0, i.e. vector field.
int main(int argc, char *argv[])
const MatrixAdaptor getVectorN(const FieldApproximationBase base, const int gg)
get Hdiv of base functions at Gauss pts
static const double face_coords[4][9]
const double v
phase velocity of light in medium (cm/ns)
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
const Tensor2_symmetric_Expr< const ddTensor0< T, Dim, i, j >, typename promote< T, double >::V, Dim, i, j > dd(const Tensor0< T * > &a, const Index< i, Dim > index1, const Index< j, Dim > index2, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
#define CATCH_ERRORS
Catch errors.
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
MoFEMErrorCode partitionSimpleProblem(const std::string name, int verb=VERBOSE)
partition problem dofs
Set inverse jacobian to base functions.
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
const FTensor::Tensor2< T, Dim, Dim > Vec
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
PetscErrorCode PetscOptionsGetEList(PetscOptions *, const char pre[], const char name[], const char *const *list, PetscInt next, PetscInt *value, PetscBool *set)
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
UBlasVector< double > VectorDouble
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
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
transform local reference derivatives of shape function to global derivatives if higher order geometr...
@ MOFEM_ATOM_TEST_INVALID
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
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.
MoFEMErrorCode partitionGhostDofs(const std::string name, int verb=VERBOSE)
determine ghost nodes
virtual MoFEMErrorCode add_problem(const std::string &name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add problem.
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ HDIV
field with continuous normal traction
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.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ OPROW
operator doWork function is executed on FE rows