|
| v0.14.0
|
Go to the documentation of this file.
21 namespace bio = boost::iostreams;
23 using bio::tee_device;
25 using namespace MoFEM;
27 static char help[] =
"...\n\n";
60 boost::shared_ptr<BaseFunctionCtx> ctx_ptr) {
65 int nb_gauss_pts = pts.size2();
70 if (pts.size1() < 3) {
72 "Wrong dimension of pts, should be at least 3 rows with "
78 CHKERR getValueHdivForCGGBubble(pts);
99 "Wrong base, should be USER_BASE");
106 int nb_gauss_pts = pts.size2();
109 shapeFun.resize(nb_gauss_pts, 4,
false);
111 &pts(2, 0), nb_gauss_pts);
113 double diff_shape_fun[12];
120 double *diff_phi_f[4];
123 for (
int ff = 0; ff != 4; ff++) {
125 int order = volume_order > face_order ? volume_order : face_order;
131 phi_f[ff] = &*data.
dataOnEntities[MBTRI][ff].getN(base).data().begin();
138 diff_shape_fun, phi_f[ff], diff_phi_f[ff], nb_gauss_pts, 4);
147 double *phi_v = &*data.
dataOnEntities[MBTET][0].getN(base).data().begin();
151 volume_order, &*shapeFun.data().begin(), diff_shape_fun, p_f, phi_f,
152 diff_phi_f, phi_v, diff_phi_v, nb_gauss_pts);
156 for (
int ff = 0; ff != 4; ff++) {
168 int main(
int argc,
char *argv[]) {
181 PetscBool flg = PETSC_TRUE;
183 #if PETSC_VERSION_GE(3, 6, 4)
190 if (flg != PETSC_TRUE) {
192 "*** ERROR -my_file (MESH FILE NEEDED)");
219 auto field_order_table =
220 const_cast<Field *
>(field_ptr)->getFieldOrderTable();
223 auto get_cgg_bubble_order_zero = [](
int p) {
return 0; };
225 auto get_cgg_bubble_order_face = [](
int p) {
228 auto get_cgg_bubble_order_tet = [](
int p) {
231 field_order_table[MBVERTEX] = get_cgg_bubble_order_zero;
232 field_order_table[MBEDGE] = get_cgg_bubble_order_zero;
233 field_order_table[MBTRI] = get_cgg_bubble_order_face;
234 field_order_table[MBTET] = get_cgg_bubble_order_tet;
235 const_cast<Field *
>(field_ptr)->rebuildDofsOrderMap();
238 for(
auto d = 0;
d!=10; ++
d) {
239 MOFEM_LOG(
"WORLD", Sev::noisy) <<
"dof " << dof_order_map[
d];
297 typedef tee_device<std::ostream, std::ofstream>
TeeDevice;
300 std::ofstream ofs(
"forces_and_sources_testing_users_base.txt");
312 MyOp1(
const std::string &row_filed,
const std::string &col_field,
315 row_filed, col_field,
type),
316 my_split(_my_split) {
326 my_split << rowFieldName << endl;
327 my_split <<
"side: " << side <<
" type: " <<
type << std::endl;
328 my_split << data << endl;
329 my_split << data.
getN() << endl;
334 MoFEMErrorCode doWork(
int row_side,
int col_side, EntityType row_type,
343 my_split << rowFieldName <<
" : " << colFieldName << endl;
344 my_split <<
"row side: " << row_side <<
" row_type: " << row_type
346 my_split <<
"col side: " << col_side <<
" col_type: " << col_type
348 my_split << row_data.
getIndices().size() <<
" : "
363 new MyOp1(
"FILED_CGG",
"FILED_CGG", my_split,
366 new MyOp1(
"FILED_CGG",
"FILED_RT", my_split,
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce 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.
MoFEMErrorCode Hdiv_Demkowicz_Face_MBTET_ON_FACE(int *faces_nodes, int p, double *N, double *diffN, double *phi_f, double *diff_phi_f, int gdim, int nb)
const std::array< ApproximationOrder, MAX_DOFS_ON_ENTITY > & getDofOrderMap(const EntityType type) const
get hash-map relating dof index on entity with its order
Problem manager is used to build and partition problems.
Class used to pass element data to calculate base functions on tet,triangle,edge.
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
virtual const Field * get_field_structure(const std::string &name, enum MoFEMTypes bh=MF_EXIST) const =0
get field structure
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
UBlasMatrix< double > MatrixDouble
MoFEMErrorCode Hdiv_Demkowicz_Interior_MBTET(int p, double *N, double *diffN, int p_face[], double *phi_f[4], double *diff_phi_f[4], double *phi_v, double *diff_phi_v, int gdim)
Base class if inherited used to calculate base functions.
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
@ OPROWCOL
operator doWork is executed on FE rows &columns
MatrixInt facesNodes
nodes on finite element faces
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.
Provide data structure for (tensor) field approximation.
@ USER_BASE
user implemented approximation base
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.
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
auto & getUserPolynomialBase()
Get the User Polynomial Base object.
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
MoFEMErrorCode getValue(MatrixDouble &pts, boost::shared_ptr< BaseFunctionCtx > ctx_ptr)
Calculate base functions at intergeneration points.
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
#define NBFACETRI_DEMKOWICZ_HDIV(P)
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.
friend class UserDataOperator
const VectorInt & getIndices() const
Get global indices of dofs on entity.
#define NBVOLUMETET_DEMKOWICZ_HDIV(P)
Volume finite element base.
virtual MoFEMErrorCode modify_problem_ref_level_add_bit(const std::string &name_problem, const BitRefLevel &bit)=0
add ref level to problem
const FieldApproximationBase bAse
base class for all interface classes
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....
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
Return interface to this class when one ask for for tetrahedron, otherisw return interface class for ...
Implementation of H-curl base function.
boost::ptr_deque< UserDataOperator > & getOpPtrVector()
Use to push back operator for row operator.
@ MOFEM_DATA_INCONSISTENCY
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_filed)=0
set finite element field data
FieldApproximationBase
approximation base
std::array< boost::ptr_vector< EntData >, MBMAXTYPE > dataOnEntities
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.
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 ...
MoFEMErrorCode getValueHdivForCGGBubble(MatrixDouble &pts)
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.
data structure for finite element entity
Class used to calculate base functions at integration points.
int main(int argc, char *argv[])
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()
EntPolynomialBaseCtx * cTx
@ HDIV
field with continuous normal traction
#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