Checking approximation functions for quad.
static char help[] =
"...\n\n";
static inline double fun(
double x,
double y) {
double r = 1;
for (
int i = 0;
i <=
o; ++
i) {
r += pow(x,
i) * pow(y,
j);
}
}
return r;
}
r.clear();
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;
}
}
}
return r;
}
};
QuadOpCheck(boost::shared_ptr<VectorDouble> &field_vals,
boost::shared_ptr<MatrixDouble> &diff_field_vals);
private:
};
private:
};
private:
};
int main(
int argc,
char *argv[]) {
try {
enum bases {
AINSWORTH,
AINSWORTH_LOBATTO,
DEMKOWICZ,
BERNSTEIN,
LASBASETOP
};
const char *list_bases[] = {"ainsworth", "ainsworth_labatto", "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)
moab::Core mb_instance;
moab::Interface &moab = mb_instance;
std::array<double, 12> one_quad_coords = {0, 0, 0,
2, 0, 0,
1, 1, 0,
0, 1, 0};
std::array<EntityHandle, 4> one_quad_nodes;
for (
int n = 0;
n != 4; ++
n)
CHKERR moab.create_vertex(&one_quad_coords[3 *
n], one_quad_nodes[
n]);
CHKERR moab.create_element(MBQUAD, one_quad_nodes.data(), 4, one_quad);
one_quad_range.insert(one_quad);
CHKERR moab.get_adjacencies(one_quad_range, 1,
true, one_quad_adj_ents,
moab::Interface::UNION);
0, 2, bit_level0);
->createMPIAIJWithArrays<PetscGlobalIdx_mi_tag>(
"TEST_PROBLEM",
A);
auto rule = [&](
int,
int,
int p) {
return 2 * (
p + 1); };
auto assemble_matrices_and_vectors = [&]() {
fe.getRuleHook = rule;
auto jac_ptr = boost::make_shared<MatrixDouble>();
auto inv_jac_ptr = boost::make_shared<MatrixDouble>();
auto det_ptr = boost::make_shared<VectorDouble>();
fe.getOpPtrVector().push_back(
fe.getOpPtrVector().push_back(
fe.getOpPtrVector().push_back(
CHKERR MatAssemblyBegin(
A, MAT_FINAL_ASSEMBLY);
CHKERR MatAssemblyEnd(
A, MAT_FINAL_ASSEMBLY);
};
auto solve_problem = [&] {
CHKERR KSPSetFromOptions(solver);
CHKERR VecGhostUpdateBegin(
D, INSERT_VALUES, SCATTER_FORWARD);
CHKERR VecGhostUpdateEnd(
D, INSERT_VALUES, SCATTER_FORWARD);
"TEST_PROBLEM",
COL,
D, INSERT_VALUES, SCATTER_REVERSE);
};
auto check_solution = [&] {
fe.getRuleHook = rule;
auto field_vals_ptr = boost::make_shared<VectorDouble>();
auto diff_field_vals_ptr = boost::make_shared<MatrixDouble>();
auto jac_ptr = boost::make_shared<MatrixDouble>();
auto inv_jac_ptr = boost::make_shared<MatrixDouble>();
auto det_ptr = boost::make_shared<VectorDouble>();
fe.getOpPtrVector().push_back(
fe.getOpPtrVector().push_back(
fe.getOpPtrVector().push_back(
fe.getOpPtrVector().push_back(
"FIELD1", diff_field_vals_ptr, space ==
L2 ? MBQUAD : MBVERTEX));
fe.getOpPtrVector().push_back(
};
CHKERR assemble_matrices_and_vectors();
}
return 0;
}
boost::shared_ptr<MatrixDouble> &diff_field_vals)
fieldVals(field_vals), diffFieldVals(diff_field_vals) {}
if (type == MBQUAD) {
const int nb_gauss_pts = data.
getN().size1();
for (int gg = 0; gg != nb_gauss_pts; ++gg) {
constexpr double eps = 1e-6;
std::cout <<
f - (*fieldVals)[gg] << std::endl;
"Wrong value %d : %6.4e != %6.4e", gg, f, (*
fieldVals)[gg]);
for (auto d : {0, 1})
std::cout << std::endl;
for (auto d : {0, 1})
"Wrong derivative value (%d) %6.4e != %6.4e", diff_f[d],
++t_coords;
}
}
}
if (nb_dofs) {
const int nb_gauss_pts = data.
getN().size1();
nf.clear();
for (int gg = 0; gg != nb_gauss_pts; ++gg) {
double *val = &*nf.begin();
for (int bb = 0; bb != nb_dofs; ++bb) {
++t_base;
++val;
}
++t_coords;
++t_w;
}
}
}
:
OpEle(
"FIELD1",
"FIELD1",
sYmm = false;
}
const int row_nb_dofs = row_data.
getIndices().size();
const int col_nb_dofs = col_data.
getIndices().size();
if (row_nb_dofs && col_nb_dofs) {
const int nb_gauss_pts = row_data.
getN().size1();
double *row_base_ptr = &*row_data.
getN().data().begin();
double *col_base_ptr = &*col_data.
getN().data().begin();
for (int gg = 0; gg != nb_gauss_pts; ++gg) {
cblas_dger(CblasRowMajor, row_nb_dofs, col_nb_dofs,
v, row_base_ptr, 1,
col_base_ptr, 1, &*
m.data().begin(), col_nb_dofs);
row_base_ptr += row_nb_dofs;
col_base_ptr += col_nb_dofs;
++t_w;
}
}
}
ForcesAndSourcesCore::UserDataOperator UserDataOperator
#define CATCH_ERRORS
Catch errors.
FieldApproximationBase
approximation base
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
@ AINSWORTH_BERNSTEIN_BEZIER_BASE
FieldSpace
approximation spaces
@ L2
field with C-1 continuity
#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.
FTensor::Index< 'n', SPACE_DIM > n
FTensor::Index< 'm', SPACE_DIM > m
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
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
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
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_filed)=0
set finite element field data
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 MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
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.
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.
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 partitionGhostDofs(const std::string name, int verb=VERBOSE)
determine ghost nodes
MoFEMErrorCode partitionSimpleProblem(const std::string name, int verb=VERBOSE)
partition problem dofs
MoFEMErrorCode buildProblem(const std::string name, const bool square_matrix, int verb=VERBOSE)
build problem data structures
MoFEMErrorCode partitionFiniteElements(const std::string name, bool part_from_moab=false, int low_proc=-1, int hi_proc=-1, int verb=VERBOSE)
partition finite elements
virtual MoFEMErrorCode add_problem(const std::string &name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add problem.
virtual MoFEMErrorCode modify_problem_ref_level_add_bit(const std::string &name_problem, const BitRefLevel &bit)=0
add ref level to problem
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
FTensor::Index< 'i', SPACE_DIM > i
auto fun
Function to approximate.
auto diff_fun
Function derivative.
const double v
phase velocity of light in medium (cm/ns)
FTensor::Index< 'j', 3 > j
const Tensor1_Expr< const dTensor0< T, Dim, i >, typename promote< T, double >::V, Dim, i > d(const Tensor0< T * > &a, const Index< i, Dim > index, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
VectorBoundedArray< double, 3 > VectorDouble3
UBlasMatrix< double > MatrixDouble
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
UBlasVector< double > VectorDouble
implementation of Data Operators for Forces and Sources
auto createKSP(MPI_Comm comm)
MoFEMErrorCode MatSetValues(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const double *ptr, InsertMode iora)
Assemble PETSc matrix.
PetscErrorCode PetscOptionsGetEList(PetscOptions *, const char pre[], const char name[], const char *const *list, PetscInt next, PetscInt *value, PetscBool *set)
MoFEMErrorCode VecSetValues(Vec V, const EntitiesFieldData::EntData &data, const double *ptr, InsertMode iora)
Assemble PETSc vector.
static constexpr int approx_order
static constexpr int approx_order
static VectorDouble3 diff_fun(double x, double y, double z)
static double fun(double x, double y, double z)
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
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.
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)
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
const VectorInt & getIndices() const
Get global indices of dofs on entity.
default operator for TRI element
auto getFTensor1CoordsAtGaussPts()
Get coordinates at integration points assuming linear geometry.
auto getFTensor0IntegrationWeight()
Get integration weights.
double getMeasure() const
get measure of element
structure to get information form mofem into EntitiesFieldData
Matrix manager is used to build and partition problems.
Get field gradients at integration pts for scalar filed rank 0, i.e. vector field.
Get value at integration points for scalar field.
Set inverse jacobian to base functions.
Modify integration weights on face to take in account higher-order geometry.
Problem manager is used to build and partition problems.
intrusive_ptr for managing petsc objects
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
Vector manager is used to create vectors \mofem_vectors.
boost::shared_ptr< VectorDouble > fieldVals
boost::shared_ptr< MatrixDouble > diffFieldVals
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
QuadOpCheck(boost::shared_ptr< VectorDouble > &field_vals, boost::shared_ptr< MatrixDouble > &diff_field_vals)
MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
Operator for bi-linear form, usually to calculate values on left hand side.
QuadOpLhs(SmartPetscObj< Mat > &a)
QuadOpRhs(SmartPetscObj< Vec > &f)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.