static char help[] =
"...\n\n";
static double fUn(const double x, const double y, double z) {
for (
int i = 0;
i <= o; ++
i) {
r += pow(x,
i) * pow(y,
j);
}
}
}
double z) {
for (
int i = 0;
i <= o; ++
i) {
r(0) +=
i > 0 ?
i * pow(x,
i - 1) * pow(y,
j) : 0;
r(1) +=
j > 0 ?
j * pow(x,
i) * pow(y,
j - 1) : 0;
}
}
}
}
};
static double fUn(const double x, const double y, double z) {
for (
int i = 0;
i <= o; ++
i) {
for (
int j = 0;
j <= o -
i;
j++) {
r += pow(x,
i) * pow(y,
j) * pow(z,
k);
}
}
}
}
}
double z) {
for (
int i = 0;
i <= o; ++
i) {
for (
int j = 0;
j <= o -
i;
j++) {
r(0) +=
i > 0 ?
i * pow(x,
i - 1) * pow(y,
j) * pow(z,
k) : 0;
r(1) +=
j > 0 ?
j * pow(x,
i) * pow(y,
j - 1) * pow(z,
k) : 0;
r(2) +=
k > 0 ?
k * pow(x,
i) * pow(y,
j) * pow(z,
k - 1) : 0;
}
}
}
}
}
};
OpValsDiffVals(VectorDouble &vals, MatrixDouble &diff_vals,
bool check_grads)
const int nb_gauss_pts = getGaussPts().size2();
if (type == MBVERTEX) {
vAls.resize(nb_gauss_pts,
false);
}
if (nb_dofs) {
<< "Type " << moab::CN::EntityTypeName(type) << " side " << side;
for (int gg = 0; gg != nb_gauss_pts; gg++) {
for (int bb = 0; bb != nb_dofs; bb++) {
t_vals += t_base_fun * t_data;
++t_base_fun;
++t_data;
}
++t_vals;
}
for (int gg = 0; gg != nb_gauss_pts; gg++) {
for (int bb = 0; bb != nb_dofs; bb++) {
t_diff_vals(
i) += t_diff_base_fun(
i) * t_data;
++t_diff_base_fun;
++t_data;
}
if (!std::isnormal(t_diff_vals(d)))
++t_diff_vals;
}
}
}
}
};
boost::shared_ptr<VectorDouble>
ptrVals;
boost::shared_ptr<VectorDouble> &ptr_vals,
boost::shared_ptr<MatrixDouble> &ptr_diff_vals,
bool check_grads)
}
const int nb_gauss_pts = data.
getN().size1();
for (int gg = 0; gg != nb_gauss_pts; gg++) {
double err_val;
err_val = std::abs(t_vals - t_ptr_vals);
MOFEM_LOG(
"AT", Sev::noisy) <<
"Val op error " << err_val;
"Wrong value from operator %4.3e", err_val);
err_val = std::fabs(delta_val * delta_val);
MOFEM_LOG(
"AT", Sev::verbose) << err_val <<
" : " << t_vals;
err_val);
++t_vals;
++t_ptr_vals;
}
for (int gg = 0; gg != nb_gauss_pts; gg++) {
double err_diff_val;
t_delta_diff_val(
i) = t_diff_vals(
i) - t_ptr_diff_vals(
i);
err_diff_val = sqrt(t_delta_diff_val(
i) * t_delta_diff_val(
i));
MOFEM_LOG(
"AT", Sev::noisy) <<
"Diff val op error " << err_diff_val;
"Wrong derivatives from operator %4.3e", err_diff_val);
t_delta_diff_val(
i) = t_diff_vals(
i) - t_diff_anal(
i);
err_diff_val = sqrt(t_delta_diff_val(
i) * t_delta_diff_val(
i));
<< "Diff val " << err_diff_val << " : "
<< sqrt(t_diff_vals(
i) * t_diff_vals(
i)) <<
" : "
<< t_diff_vals(0) << " (" << t_diff_anal(0) << ") "
<< t_diff_vals(1) << " (" << t_diff_anal(1) << ") "
<< t_diff_vals(2) << " (" << t_diff_anal(2) << ")";
else
<< "Diff val " << err_diff_val << " : "
<< sqrt(t_diff_vals(
i) * t_diff_vals(
i)) <<
" : "
<< t_diff_vals(0) << " (" << t_diff_anal(0) << ") "
<< t_diff_vals(1) << " (" << t_diff_anal(1) << ")";
MOFEM_LOG(
"AT", Sev::verbose) <<
"Diff val error " << err_diff_val;
"Wrong derivative of value %4.3e %4.3e", err_diff_val,
t_diff_anal.l2());
++t_diff_vals;
++t_ptr_diff_vals;
}
}
}
};
int main(
int argc,
char *argv[]) {
try {
DMType dm_name = "DMMOFEM";
moab::Core mb_instance;
moab::Interface &moab = mb_instance;
auto core_log = logging::core::get();
core_log->add_sink(
LogManager::createSink(LogManager::getStrmWorld(), "AT"));
LogManager::setLog("AT");
simple->getAddBoundaryFE() =
true;
enum bases {
AINSWORTH,
AINSWORTH_LOBATTO,
DEMKOWICZ,
BERNSTEIN,
LASBASETOP
};
const char *list_bases[] = {"ainsworth", "ainsworth_lobatto", "demkowicz",
"bernstein"};
PetscBool flg;
PetscInt choice_base_value = AINSWORTH;
LASBASETOP, &choice_base_value, &flg);
if (flg != PETSC_TRUE)
if (choice_base_value == AINSWORTH)
if (choice_base_value == AINSWORTH_LOBATTO)
else if (choice_base_value == DEMKOWICZ)
else if (choice_base_value == BERNSTEIN)
enum spaces { H1SPACE, L2SPACE, LASBASETSPACE };
const char *list_spaces[] = {"h1", "l2"};
PetscInt choice_space_value = H1SPACE;
LASBASETSPACE, &choice_space_value, &flg);
if (flg != PETSC_TRUE)
if (choice_space_value == H1SPACE)
else if (choice_space_value == L2SPACE)
PETSC_NULL);
CHKERR moab.get_entities_by_dimension(0, 1, edges);
CHKERR moab.get_entities_by_dimension(0, 2, faces);
if (choice_base_value != BERNSTEIN) {
for (auto e : edges) {
rise_order.insert(e);
}
}
for (auto f : faces) {
rise_order.insert(f);
}
}
for (auto e : rise_order) {
rise_twice.insert(e);
}
}
}
auto volume_adj = [](moab::Interface &moab,
const Field &field,
std::vector<EntityHandle> &adjacency) {
switch (field.getSpace()) {
CHKERR moab.get_connectivity(&fe_ent, 1, adjacency,
true);
CHKERR moab.get_adjacencies(&fe_ent, 1, 1,
false, adjacency,
moab::Interface::UNION);
CHKERR moab.get_adjacencies(&fe_ent, 1, 2,
false, adjacency,
moab::Interface::UNION);
adjacency.push_back(fe_ent);
for (auto ent : adjacency)
fe.getSideNumberPtr(ent);
break;
CHKERR moab.get_entities_by_handle(field.getMeshset(), adjacency,
false);
for (auto ent : adjacency) {
.insert(
boost::shared_ptr<SideNumber>(
new SideNumber(ent, -1, 0, 0)));
}
} break;
default:
"this field is not implemented for TRI finite element");
}
};
simple->getDomainFEName(), MBTET, volume_adj);
simple->getDomainFEName(), MBHEX, volume_adj);
VectorDouble vals;
MatrixDouble diff_vals;
auto assemble_matrices_and_vectors = [&]() {
pipeline_mng->getOpDomainLhsPipeline(), {space});
pipeline_mng->getOpDomainRhsPipeline(), {space});
pipeline_mng->getOpDomainLhsPipeline().push_back(
pipeline_mng->getOpDomainLhsPipeline().push_back(
new OpMass(
"FIELD1", "FIELD1", [](double, double, double) { return 1.; }));
pipeline_mng->getOpDomainLhsPipeline().push_back(
pipeline_mng->getOpDomainRhsPipeline().push_back(
new OpSource("FIELD1", ApproxFunctions::fUn));
return 2 * p_data + 1;
};
};
auto solve_problem = [&] {
auto solver = pipeline_mng->createKSP();
CHKERR KSPSetFromOptions(solver);
CHKERR VecGhostUpdateBegin(
D, INSERT_VALUES, SCATTER_FORWARD);
CHKERR VecGhostUpdateEnd(
D, INSERT_VALUES, SCATTER_FORWARD);
};
auto check_solution = [&] {
auto ptr_values = boost::make_shared<VectorDouble>();
auto ptr_diff_vals = boost::make_shared<MatrixDouble>();
pipeline_mng->getDomainLhsFE().reset();
pipeline_mng->getOpDomainRhsPipeline().clear();
pipeline_mng->getOpDomainRhsPipeline(), {space});
pipeline_mng->getOpDomainRhsPipeline().push_back(
pipeline_mng->getOpDomainRhsPipeline().push_back(
pipeline_mng->getOpDomainRhsPipeline().push_back(
ptr_diff_vals));
vals, diff_vals, ptr_values, ptr_diff_vals, true));
CHKERR pipeline_mng->loopFiniteElements();
};
auto post_proc = [&] {
auto get_domain_post_proc_fe = [&](auto post_proc_mesh_ptr) {
auto post_proc_fe =
boost::make_shared<PostProcBrokenMeshInMoabBaseCont<DomainEle>>(
m_field, post_proc_mesh_ptr);
post_proc_fe->getOpPtrVector(), {space});
auto ptr_values = boost::make_shared<VectorDouble>();
auto ptr_diff_vals = boost::make_shared<MatrixDouble>();
post_proc_fe->getOpPtrVector().push_back(
post_proc_fe->getOpPtrVector().push_back(
ptr_diff_vals));
post_proc_fe->getOpPtrVector().push_back(
post_proc_fe->getPostProcMesh(), post_proc_fe->getMapGaussPts(),
{{"FIELD1", ptr_values}},
{{"FIELD1_GRAD", ptr_diff_vals}},
{}, {})
);
return post_proc_fe;
};
auto get_bdy_post_proc_fe = [&](auto post_proc_mesh_ptr) {
auto bdy_post_proc_fe =
boost::make_shared<PostProcBrokenMeshInMoabBaseCont<BoundaryEle>>(
m_field, post_proc_mesh_ptr);
boost::make_shared<
ForcesAndSourcesCore::UserDataOperator::AdjCache>());
auto ptr_values = boost::make_shared<VectorDouble>();
auto ptr_diff_vals = boost::make_shared<MatrixDouble>();
op_loop_side->getOpPtrVector(), {space});
op_loop_side->getOpPtrVector().push_back(
op_loop_side->getOpPtrVector().push_back(
ptr_diff_vals));
bdy_post_proc_fe->getOpPtrVector().push_back(op_loop_side);
bdy_post_proc_fe->getOpPtrVector().push_back(
bdy_post_proc_fe->getPostProcMesh(),
bdy_post_proc_fe->getMapGaussPts(),
{{"FIELD1", ptr_values}},
{{"FIELD1_GRAD", ptr_diff_vals}},
{}, {})
);
return bdy_post_proc_fe;
};
auto post_proc_mesh_ptr = boost::make_shared<moab::Core>();
auto post_proc_begin =
boost::make_shared<PostProcBrokenMeshInMoabBaseBegin>(
m_field, post_proc_mesh_ptr);
auto post_proc_end =
boost::make_shared<PostProcBrokenMeshInMoabBaseEnd>(
m_field, post_proc_mesh_ptr);
auto domain_post_proc_fe = get_domain_post_proc_fe(post_proc_mesh_ptr);
auto bdy_post_proc_fe = get_bdy_post_proc_fe(post_proc_mesh_ptr);
domain_post_proc_fe);
bdy_post_proc_fe);
CHKERR post_proc_end->writeFile(
"out_post_proc.h5m");
};
CHKERR assemble_matrices_and_vectors();
}
}
void simple(double P1[], double P2[], double P3[], double c[], const int N)
ElementsAndOps< SPACE_DIM >::BoundaryEle BoundaryEle
ElementsAndOps< SPACE_DIM >::DomainEle DomainEle
[Define dimension]
#define CATCH_ERRORS
Catch errors.
FieldApproximationBase
approximation base
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base nme:nme847.
@ AINSWORTH_BERNSTEIN_BEZIER_BASE
FieldSpace
approximation spaces
@ L2
field with C-1 continuity
@ NOFIELD
scalar or vector of scalars describe (no true field)
@ HCURL
field with continuous tangents
@ HDIV
field with continuous normal traction
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_ATOM_TEST_INVALID
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
PetscErrorCode DMoFEMPostProcessFiniteElements(DM dm, MoFEM::FEMethod *method)
execute finite element method for each element in dm (problem)
PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode)
set local (or ghosted) vector values on mesh for partition only
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Executes FEMethod for finite elements in DM.
auto createDMVector(DM dm)
Get smart vector from DM.
PetscErrorCode DMoFEMPreProcessFiniteElements(DM dm, MoFEM::FEMethod *method)
execute finite element method for each element in dm (problem)
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.
virtual MoFEMErrorCode modify_finite_element_adjacency_table(const std::string &fe_name, const EntityType type, ElementAdjacencyFunct function)=0
modify finite element table, only for advanced user
#define MOFEM_LOG(channel, severity)
Log.
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
FTensor::Index< 'i', SPACE_DIM > i
const double v
phase velocity of light in medium (cm/ns)
FTensor::Index< 'j', 3 > j
FTensor::Index< 'k', 3 > k
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
implementation of Data Operators for Forces and Sources
OpSchurAssembleBase * createOpSchurAssembleEnd(std::vector< std::string > fields_name, std::vector< boost::shared_ptr< Range > > field_ents, SmartPetscObj< AO > ao, SmartPetscObj< Mat > schur, bool sym_schur, bool symm_op)
Construct a new Op Schur Assemble End object.
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
SmartPetscObj< Vec > vectorDuplicate(Vec vec)
Create duplicate vector of smart vector.
FTensor::Tensor1< FTensor::PackPtr< T *, S >, Tensor_Dim > getFTensor1FromMat(ublas::matrix< T, L, A > &data)
Get tensor rank 1 (vector) form data matrix.
static auto getFTensor0FromVec(ublas::vector< T, A > &data)
Get tensor rank 0 (scalar) form data vector.
PetscErrorCode PetscOptionsGetEList(PetscOptions *, const char pre[], const char name[], const char *const *list, PetscInt next, PetscInt *value, PetscBool *set)
OpSchurAssembleBase * createOpSchurAssembleBegin()
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
#define EXECUTABLE_DIMENSION
static constexpr int approx_order
PipelineManager::ElementsAndOpsByDim< SPACE_DIM >::FaceSideEle EleOnSide
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, SPACE_DIM > OpMass
[Only used with Hooke equation (linear material model)]
static FTensor::Tensor2< double, BASE_DIM, SPACE_DIM > diffFun(const double x, const double y)
static FTensor::Tensor1< double, BASE_DIM > fUn(const double x, const double y, double z)
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.
std::array< bool, MBMAXTYPE > doEntities
If true operator is executed for entity.
Deprecated interface functions.
Finite element data for entity.
Data on single entity (This is passed as argument to DataOperator::doWork)
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1DiffN(const FieldApproximationBase base)
Get derivatives of base functions.
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
MatrixDouble & getDiffN(const FieldApproximationBase base)
get derivatives of base functions
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0FieldData()
Resturn scalar files as a FTensor of rank 0.
const VectorInt & getIndices() const
Get global indices of dofs on entity.
VectorDouble & getCoords()
get triangle coordinates
Provide data structure for (tensor) field approximation.
MatrixDouble & getCoordsAtGaussPts()
Gauss points and weight, matrix (nb. of points x 3)
@ OPROW
operator doWork function is executed on FE rows
Get field gradients at integration pts for scalar filed rank 0, i.e. vector field.
Get value at integration points for scalar field.
Element used to execute operators on side of the element.
Post post-proc data at points from hash maps.
PipelineManager interface.
keeps information about side number for the finite element
Simple interface for fast problem set-up.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
FTensor::Index< 'i', 3 > i
boost::shared_ptr< MatrixDouble > ptrVals
const bool checkGradients
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
boost::shared_ptr< MatrixDouble > ptrDiffVals
OpCheckValsDiffVals(boost::shared_ptr< MatrixDouble > ptr_vals, boost::shared_ptr< VectorDouble > ptr_div, boost::shared_ptr< MatrixDouble > ptr_grad, boost::shared_ptr< MatrixDouble > ptr_hess)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
const bool checkGradients
OpValsDiffVals(VectorDouble &vals, MatrixDouble &diff_vals, bool check_grads)
FTensor::Index< 'i', SPACE_DIM > i