|
| v0.14.0
|
Go to the source code of this file.
|
struct | MyOp< OP > |
| Operator used to check consistency between local coordinates and global cooridnates for integrated points and evaluated points. More...
|
|
|
int | main (int argc, char *argv[]) |
|
◆ main()
int main |
( |
int |
argc, |
|
|
char * |
argv[] |
|
) |
| |
- Examples
- field_evaluator.cpp.
Definition at line 77 of file field_evaluator.cpp.
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);
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);
◆ get_rule
◆ help
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
RuleHookFun getRuleHook
Hook to get rule.
virtual int get_comm_rank() const =0
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.
Simple interface for fast problem set-up.
Deprecated interface functions.
DeprecatedCoreInterface Interface
MoFEMErrorCode getOptions()
get options
MoFEMErrorCode buildTree3D(boost::shared_ptr< SetPtsData > spd_ptr, const std::string finite_element)
Build spatial tree.
MoFEMErrorCode getDM(DM *dm)
Get DM.
#define CHKERR
Inline error check.
friend class UserDataOperator
MoFEMErrorCode buildTree2D(boost::shared_ptr< SetPtsData > spd_ptr, const std::string finite_element)
Build spatial tree.
default operator for TRI element
Field evaluator interface.
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.
friend class UserDataOperator
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.
const std::string getDomainFEName() const
Get the Domain FE Name.
Volume finite element base.
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 setFieldOrder(const std::string field_name, const int order, const Range *ents=NULL)
Set field order.
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
#define CATCH_ERRORS
Catch errors.
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Operator used to check consistency between local coordinates and global cooridnates for integrated po...
@ MOFEM_DATA_INCONSISTENCY
PetscErrorCode DMMoFEMGetProblemPtr(DM dm, const MoFEM::Problem **problem_ptr)
Get pointer to problem data structure.
int getDim() const
Get the problem dimension.
keeps basic data about problem
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.
MoFEMErrorCode setUp(const PetscBool is_partitioned=PETSC_TRUE)
Setup problem.