|
| v0.14.0
|
Go to the documentation of this file.
13 using namespace MoFEM;
15 static char help[] =
"...\n\n";
40 if (
type == MBEDGE && side == getEdgeSideNumber()) {
42 MatrixDouble diff = getCoordsAtGaussPts() - getEdgeCoordsAtGaussPts();
44 const double eps = 1e-12;
45 if (norm_inf(diff) >
eps)
47 "coordinates at integration pts are different");
49 const size_t nb_dofs = data.
getN().size2() / 3;
50 const size_t nb_integration_pts = data.
getN().size1();
52 auto t_tangent = getFTensor1Direction();
56 if (getFEMethod()->nInTheLoop == 0)
61 dot_elem_data.resize(nb_integration_pts, nb_dofs,
false);
63 for (
size_t gg = 0; gg != nb_integration_pts; ++gg) {
64 for (
size_t bb = 0; bb != nb_dofs; ++bb) {
65 dot_elem_data(gg, bb) = t_tangent(
i) * t_hcurl_base(
i);
79 faceSideFe(m_field), elemData(elem_data) {
90 const size_t nb_dofs = data.
getN().size2() / 3;
91 const size_t nb_integration_pts = data.
getN().size1();
93 auto t_tangent = getFTensor1Direction();
96 elemData.
dotEdge.resize(nb_integration_pts, nb_dofs,
false);
97 elemData.
dotEleLeft.resize(nb_integration_pts, 0,
false);
98 elemData.
dotEleRight.resize(nb_integration_pts, 0,
false);
100 for (
size_t gg = 0; gg != nb_integration_pts; ++gg) {
101 for (
size_t bb = 0; bb != nb_dofs; ++bb) {
102 elemData.
dotEdge(gg, bb) = t_tangent(
i) * t_hcurl_base(
i);
107 CHKERR loopSideFaces(
"dFE", faceSideFe);
109 auto check_continuity_of_base = [&](
auto &vol_dot_data) {
112 if (vol_dot_data.size1() != elemData.
dotEdge.size1())
114 "Inconsistent number of integration points");
116 if (vol_dot_data.size2() != elemData.
dotEdge.size2())
118 "Inconsistent number of base functions");
119 const double eps = 1e-12;
120 for (
size_t gg = 0; gg != vol_dot_data.size1(); ++gg)
121 for (
size_t bb = 0; bb != vol_dot_data.size2(); ++bb) {
123 std::abs(vol_dot_data(gg, bb) - elemData.
dotEdge(gg, bb));
126 "Inconsistency (%d, %d) %3.4e != %3.4e", gg, bb,
127 vol_dot_data(gg, bb), elemData.
dotEdge(gg, bb));
145 int main(
int argc,
char *argv[]) {
160 auto core_log = logging::core::get();
166 DMType dm_name =
"DMMOFEM";
179 enum bases { AINSWORTH, DEMKOWICZ, LASTBASEOP };
180 const char *list_bases[] = {
"ainsworth",
"demkowicz"};
182 PetscInt choice_base_value = AINSWORTH;
184 LASTBASEOP, &choice_base_value, &flg);
185 if (flg == PETSC_TRUE) {
187 if (choice_base_value == AINSWORTH)
189 else if (choice_base_value == DEMKOWICZ)
197 auto base = get_base();
205 auto dm = simple_interface->
getDM();
209 boost::shared_ptr<EdgeEle> skeleton_fe =
210 boost::shared_ptr<EdgeEle>(
new EdgeEle(m_field));
213 skeleton_fe->getOpPtrVector().push_back(
#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)
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
UBlasMatrix< double > MatrixDouble
const std::string getSkeletonFEName() const
Get the Skeleton FE Name.
MoFEMErrorCode loadFile(const std::string options, const std::string mesh_file_name, LoadFileFunc loadFunc=defaultLoadFileFunc)
Load mesh file.
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
static boost::shared_ptr< SinkType > createSink(boost::shared_ptr< std::ostream > stream_ptr, std::string comm_filter)
Create a sink object.
Simple interface for fast problem set-up.
Deprecated interface functions.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
MoFEMErrorCode addSkeletonField(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_ZERO, int verb=-1)
Add field on skeleton.
DeprecatedCoreInterface Interface
MoFEMErrorCode getOptions()
get options
MoFEMErrorCode getDM(DM *dm)
Get DM.
#define CHKERR
Inline error check.
friend class UserDataOperator
implementation of Data Operators for Forces and Sources
MoFEMErrorCode addDomainField(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_ZERO, int verb=-1)
Add field on domain.
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
OpFaceSide(CommonData &elem_data)
static boost::shared_ptr< std::ostream > getStrmSelf()
Get the strm self object.
MoFEM::EdgeElementForcesAndSourcesCore EdgeEle
friend class UserDataOperator
FTensor::Index< 'i', SPACE_DIM > i
MoFEMErrorCode setFieldOrder(const std::string field_name, const int order, const Range *ents=NULL)
Set field order.
Add operators pushing bases from local to physical configuration.
default operator for EDGE element
int main(int argc, char *argv[])
SkeletonFE(MoFEM::Interface &m_field, CommonData &elem_data)
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
#define MOFEM_LOG(channel, severity)
Log.
#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....
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1N(FieldApproximationBase base)
Get base functions for Hdiv/Hcurl spaces.
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
boost::ptr_deque< UserDataOperator > & getOpPtrVector()
Use to push back operator for row operator.
@ HCURL
field with continuous tangents
PetscErrorCode PetscOptionsGetEList(PetscOptions *, const char pre[], const char name[], const char *const *list, PetscInt next, PetscInt *value, PetscBool *set)
default operator for Face element
FieldApproximationBase
approximation base
@ MOFEM_ATOM_TEST_INVALID
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
static LoggerType & setLog(const std::string channel)
Set ans resset chanel logger.
MoFEMErrorCode setUp(const PetscBool is_partitioned=PETSC_TRUE)
Setup problem.
PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Executes FEMethod for finite elements in DM.
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Base face element used to integrate on skeleton.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...