17static char help[] =
"...\n\n";
21 static inline double fun(
double x,
double y) {
24 for (
int i = 0;
i <=
o; ++
i) {
27 r += pow(x,
i) * pow(y,
j);
37 for (
int i = 0;
i <=
o; ++
i) {
40 r[0] +=
i > 0 ?
i * pow(x,
i - 1) * pow(y,
j) : 0;
41 r[1] +=
j > 0 ?
j * pow(x,
i) * pow(y,
j - 1) : 0;
51 QuadOpCheck(boost::shared_ptr<VectorDouble> &field_vals,
52 boost::shared_ptr<MatrixDouble> &diff_field_vals);
83int main(
int argc,
char *argv[]) {
97 const char *list_bases[] = {
"ainsworth",
"ainsworth_labatto",
"demkowicz",
100 PetscInt choice_base_value = AINSWORTH;
102 LASBASETOP, &choice_base_value, &flg);
104 if (flg != PETSC_TRUE)
107 if (choice_base_value == AINSWORTH)
109 if (choice_base_value == AINSWORTH_LOBATTO)
111 else if (choice_base_value == DEMKOWICZ)
113 else if (choice_base_value == BERNSTEIN)
116 enum spaces { H1SPACE, L2SPACE, LASBASETSPACE };
117 const char *list_spaces[] = {
"h1",
"l2"};
118 PetscInt choice_space_value = H1SPACE;
120 LASBASETSPACE, &choice_space_value, &flg);
121 if (flg != PETSC_TRUE)
124 if (choice_space_value == H1SPACE)
126 else if (choice_space_value == L2SPACE)
129 moab::Core mb_instance;
130 moab::Interface &moab = mb_instance;
132 std::array<double, 12> one_quad_coords = {0, 0, 0,
140 std::array<EntityHandle, 4> one_quad_nodes;
141 for (
int n = 0;
n != 4; ++
n)
142 CHKERR moab.create_vertex(&one_quad_coords[3 *
n], one_quad_nodes[
n]);
144 CHKERR moab.create_element(MBQUAD, one_quad_nodes.data(), 4, one_quad);
145 Range one_quad_range;
146 one_quad_range.insert(one_quad);
147 Range one_quad_adj_ents;
148 CHKERR moab.get_adjacencies(one_quad_range, 1,
true, one_quad_adj_ents,
149 moab::Interface::UNION);
202 ->createMPIAIJWithArrays<PetscGlobalIdx_mi_tag>(
"TEST_PROBLEM",
A);
210 auto rule = [&](int, int,
int p) {
return 2 * (
p + 1); };
212 auto assemble_matrices_and_vectors = [&]() {
216 auto jac_ptr = boost::make_shared<MatrixDouble>();
217 auto inv_jac_ptr = boost::make_shared<MatrixDouble>();
218 auto det_ptr = boost::make_shared<VectorDouble>();
234 CHKERR MatAssemblyBegin(
A, MAT_FINAL_ASSEMBLY);
235 CHKERR MatAssemblyEnd(
A, MAT_FINAL_ASSEMBLY);
239 auto solve_problem = [&] {
241 auto solver =
createKSP(PETSC_COMM_WORLD);
242 CHKERR KSPSetOperators(solver,
A,
A);
243 CHKERR KSPSetFromOptions(solver);
246 CHKERR VecGhostUpdateBegin(
D, INSERT_VALUES, SCATTER_FORWARD);
247 CHKERR VecGhostUpdateEnd(
D, INSERT_VALUES, SCATTER_FORWARD);
249 "TEST_PROBLEM",
COL,
D, INSERT_VALUES, SCATTER_REVERSE);
253 auto check_solution = [&] {
257 auto field_vals_ptr = boost::make_shared<VectorDouble>();
258 auto diff_field_vals_ptr = boost::make_shared<MatrixDouble>();
259 auto jac_ptr = boost::make_shared<MatrixDouble>();
260 auto inv_jac_ptr = boost::make_shared<MatrixDouble>();
261 auto det_ptr = boost::make_shared<VectorDouble>();
274 "FIELD1", diff_field_vals_ptr, space ==
L2 ? MBQUAD : MBVERTEX));
276 new QuadOpCheck(field_vals_ptr, diff_field_vals_ptr));
281 CHKERR assemble_matrices_and_vectors();
292 boost::shared_ptr<MatrixDouble> &diff_field_vals)
294 fieldVals(field_vals), diffFieldVals(diff_field_vals) {}
300 if (type == MBQUAD) {
301 const int nb_gauss_pts = data.
getN().size1();
303 for (
int gg = 0; gg != nb_gauss_pts; ++gg) {
305 constexpr double eps = 1e-6;
307 std::cout << f - (*fieldVals)[gg] << std::endl;
311 "Wrong value %d : %6.4e != %6.4e", gg, f, (*
fieldVals)[gg]);
314 for (
auto d : {0, 1})
316 std::cout << std::endl;
317 for (
auto d : {0, 1})
320 "Wrong derivative value (%d) %6.4e != %6.4e", diff_f[d],
339 const int nb_gauss_pts = data.
getN().size1();
346 for (
int gg = 0; gg != nb_gauss_pts; ++gg) {
348 double v =
a * t_w * f;
349 double *val = &*nf.begin();
350 for (
int bb = 0; bb != nb_dofs; ++bb) {
365 :
OpEle(
"FIELD1",
"FIELD1",
378 const int row_nb_dofs = row_data.
getIndices().size();
379 const int col_nb_dofs = col_data.
getIndices().size();
380 if (row_nb_dofs && col_nb_dofs) {
381 const int nb_gauss_pts = row_data.
getN().size1();
386 double *row_base_ptr = &*row_data.
getN().data().begin();
387 double *col_base_ptr = &*col_data.
getN().data().begin();
388 for (
int gg = 0; gg != nb_gauss_pts; ++gg) {
390 cblas_dger(CblasRowMajor, row_nb_dofs, col_nb_dofs,
v, row_base_ptr, 1,
391 col_base_ptr, 1, &*
m.data().begin(), col_nb_dofs);
392 row_base_ptr += row_nb_dofs;
393 col_base_ptr += col_nb_dofs;
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
const double v
phase velocity of light in medium (cm/ns)
FTensor::Index< 'j', 3 > j
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
VectorBoundedArray< double, 3 > VectorDouble3
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
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 VectorDouble3 diff_fun(double x, double y, double z)
static VectorDouble3 diff_fun(double x, double y)
static double fun(double x, double y)
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.
bool sYmm
If true assume that matrix is symmetric structure.
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
double getMeasure()
get measure of element
auto getFTensor1CoordsAtGaussPts()
Get coordinates at integration points assuming linear geometry.
auto getFTensor0IntegrationWeight()
Get integration weights.
structure to get information form mofem into EntitiesFieldData
boost::ptr_deque< UserDataOperator > & getOpPtrVector()
Use to push back operator for row operator.
RuleHookFun getRuleHook
Hook to get rule.
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.