14static char help[] =
"testing field evaluator\n\n";
21template <
typename OP>
struct MyOp :
public OP {
25 MyOp(std::array<double, 12> &eval_points)
33 if (type == MBVERTEX) {
35 MOFEM_LOG(
"SELF", Sev::inform) <<
"FE " << OP::getFEEntityHandle();
37 MOFEM_LOG(
"SELF", Sev::inform) <<
"Integration pts" << std::endl;
38 MOFEM_LOG(
"SELF", Sev::inform) << OP::getGaussPts() << endl;
40 MOFEM_LOG(
"SELF", Sev::inform) <<
"Global coordinates " << endl;
41 MOFEM_LOG(
"SELF", Sev::inform) << OP::getCoordsAtGaussPts() << std::endl;
43 for (
int gg = 0; gg != OP::getCoordsAtGaussPts().size1(); ++gg) {
44 int pt_number = OP::getGaussPts()(OP::getGaussPts().size1() - 1, gg);
46 MOFEM_LOG(
"SELF", Sev::inform) <<
"gg " << gg << std::endl;
47 MOFEM_LOG(
"SELF", Sev::inform) <<
"pt " << pt_number << std::endl;
49 ublas::matrix_row<MatrixDouble> coord_at_gauss_pt(
50 OP::getCoordsAtGaussPts(), gg);
51 ublas::matrix_row<MatrixShallowArrayAdaptor<double>> eval_coord(
54 MOFEM_LOG(
"SELF", Sev::inform) <<
"coord_at_gauss_pt ";
55 MOFEM_LOG(
"SELF", Sev::inform) << coord_at_gauss_pt << std::endl;
57 MOFEM_LOG(
"SELF", Sev::inform) <<
"eval_coord ";
58 MOFEM_LOG(
"SELF", Sev::inform) << eval_coord << std::endl;
60 double error = norm_2(coord_at_gauss_pt - eval_coord);
64 "Difference at %d error = %3.4e", pt_number, error);
73auto get_rule = [](
int order_row,
int order_col,
int order_data) {
77int main(
int argc,
char *argv[]) {
84 moab::Core mb_instance;
85 moab::Interface &moab = mb_instance;
91 DMType dm_name =
"DMMOFEM";
115 auto dm = simple_interface->
getDM();
119 if (simple_interface->
getDim() == 3) {
121 std::array<double, 3> point = {0, 0, 0};
122 const double dist = 0.1;
124 std::array<double, 12> eval_points = {0.0, 0.0, 0.0, 0.0, 0.01, 0.0,
125 -0.1, 0.0, 0.0, 0.1, 0.0, 0.0};
135 if (
auto fe_method = data->feMethodPtr.lock()) {
137 fe_method->getOpPtrVector().push_back(
new MyOp<VolOp>(eval_points));
140 "Pointer to element does not exists");
146 data->setEvalPoints(eval_points.data(), eval_points.size() / 3);
150 &point[0], dist, prb_ptr->
getName(),
155 if (simple_interface->
getDim() == 2) {
157 std::array<double, 3> point = {0, 0, 0};
158 const double dist = 0.1;
160 std::array<double, 12> eval_points = {
174 if (
auto fe_method = data->feMethodPtr.lock()) {
176 fe_method->getOpPtrVector().push_back(
new MyOp<FaceOp>(eval_points));
179 "Pointer to element does not exists");
185 data->setEvalPoints(eval_points.data(), eval_points.size() / 3);
189 &point[0], dist, prb_ptr->
getName(),
#define CATCH_ERRORS
Catch errors.
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base nme:nme847.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_ATOM_TEST_INVALID
@ MOFEM_DATA_INCONSISTENCY
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
PetscErrorCode DMMoFEMGetProblemPtr(DM dm, const MoFEM::Problem **problem_ptr)
Get pointer to problem data structure.
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
#define MOFEM_LOG(channel, severity)
Log.
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
ublas::matrix< double, ublas::row_major, ublas::shallow_array_adaptor< double > > MatrixShallowArrayAdaptor
implementation of Data Operators for Forces and Sources
auto getMatrixAdaptor(T1 ptr, const size_t n, const size_t m)
Get Matrix adaptor.
virtual int get_comm_rank() const =0
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Deprecated interface functions.
Data on single entity (This is passed as argument to DataOperator::doWork)
default operator for TRI element
friend class UserDataOperator
Field evaluator interface.
MoFEMErrorCode buildTree3D(boost::shared_ptr< SetPtsData > spd_ptr, const std::string finite_element)
Build spatial tree.
MoFEMErrorCode buildTree2D(boost::shared_ptr< SetPtsData > spd_ptr, const std::string finite_element)
Build spatial tree.
MoFEMErrorCode evalFEAtThePoint2D(const double *const point, const double distance, const std::string problem, const std::string finite_element, boost::shared_ptr< SetPtsData > data_ptr, int lower_rank, int upper_rank, boost::shared_ptr< CacheTuple > cache_ptr, MoFEMTypes bh=MF_EXIST, VERBOSITY_LEVELS verb=QUIET)
Evaluate field at artbitray position.
boost::shared_ptr< SPD > getData(const double *ptr=nullptr, const int nb_eval_points=0, const double eps=1e-12, VERBOSITY_LEVELS verb=QUIET)
Get the Data object.
MoFEMErrorCode evalFEAtThePoint3D(const double *const point, const double distance, const std::string problem, const std::string finite_element, boost::shared_ptr< SetPtsData > data_ptr, int lower_rank, int upper_rank, boost::shared_ptr< CacheTuple > cache_ptr, MoFEMTypes bh=MF_EXIST, VERBOSITY_LEVELS verb=QUIET)
Evaluate field at artbitray position.
@ OPROW
operator doWork function is executed on FE rows
keeps basic data about problem
Simple interface for fast problem set-up.
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.
int getDim() const
Get the problem dimension.
MoFEMErrorCode loadFile(const std::string options, const std::string mesh_file_name, LoadFileFunc loadFunc=defaultLoadFileFunc)
Load mesh file.
MoFEMErrorCode getOptions()
get options
MoFEMErrorCode getDM(DM *dm)
Get DM.
MoFEMErrorCode setFieldOrder(const std::string field_name, const int order, const Range *ents=NULL)
Set field order.
MoFEMErrorCode setUp(const PetscBool is_partitioned=PETSC_TRUE)
Setup problem.
const std::string getDomainFEName() const
Get the Domain FE Name.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Volume finite element base.
friend class UserDataOperator
Operator used to check consistency between local coordinates and global cooridnates for integrated po...
MatrixShallowArrayAdaptor< double > evalPoints
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
MyOp(std::array< double, 12 > &eval_points)