12static char help[] =
"...\n\n";
17 static inline double fun(
double x,
double y,
double z) {
20 for (
int i = 0;
i <=
o; ++
i) {
21 for (
int j = 0;
j <= (
o -
i); ++
j) {
24 r += pow(x,
i) * pow(y,
j) * pow(z,
k);
36 for (
int i = 0;
i <=
o; ++
i) {
37 for (
int j = 0;
j <= (
o -
i); ++
j) {
40 r[0] +=
i > 0 ?
i * pow(x,
i - 1) * pow(y,
j) * pow(z,
k) : 0;
41 r[1] +=
j > 0 ?
j * pow(x,
i) * pow(y,
j - 1) * pow(z,
k) : 0;
42 r[2] +=
k > 0 ?
k * pow(x,
i) * pow(y,
j) * pow(z,
k - 1) : 0;
54 PrismOpCheck(boost::shared_ptr<VectorDouble> &field_vals,
55 boost::shared_ptr<MatrixDouble> &diff_field_vals);
91 FatPrismElementForcesAndSourcesCore;
96int main(
int argc,
char *argv[]) {
102 moab::Core mb_instance;
103 moab::Interface &moab = mb_instance;
105 ParallelComm *pcomm = ParallelComm::get_pcomm(&moab,
MYPCOMM_INDEX);
106 auto moab_comm_wrap =
107 boost::make_shared<WrapMPIComm>(PETSC_COMM_WORLD,
false);
109 pcomm =
new ParallelComm(&moab, moab_comm_wrap->get_comm());
111 std::array<double, 18> one_prism_coords = {0, 0, 0, 1, 0, 0, 0, 1, 0,
112 0, 0, 1, 1, 0, 1, 0, 1, 1};
113 std::array<EntityHandle, 6> one_prism_nodes;
114 for (
int n = 0;
n != 6; ++
n)
115 CHKERR moab.create_vertex(&one_prism_coords[3 *
n], one_prism_nodes[
n]);
117 CHKERR moab.create_element(MBPRISM, one_prism_nodes.data(), 6, one_prism);
118 Range one_prism_range;
119 one_prism_range.insert(one_prism);
120 Range one_prism_adj_ents;
121 for (
int d = 1; d != 3; ++d)
122 CHKERR moab.get_adjacencies(one_prism_range, d,
true, one_prism_adj_ents,
123 moab::Interface::UNION);
133 one_prism_range, bit_level0);
183 ->createMPIAIJWithArrays<PetscGlobalIdx_mi_tag>(
"TEST_PROBLEM",
A);
191 auto assemble_matrices_and_vectors = [&]() {
208 CHKERR MatAssemblyBegin(
A, MAT_FINAL_ASSEMBLY);
209 CHKERR MatAssemblyEnd(
A, MAT_FINAL_ASSEMBLY);
213 auto solve_problem = [&] {
215 auto solver =
createKSP(PETSC_COMM_WORLD);
216 CHKERR KSPSetOperators(solver,
A,
A);
217 CHKERR KSPSetFromOptions(solver);
220 CHKERR VecGhostUpdateBegin(
D, INSERT_VALUES, SCATTER_FORWARD);
221 CHKERR VecGhostUpdateEnd(
D, INSERT_VALUES, SCATTER_FORWARD);
223 "TEST_PROBLEM",
COL,
D, INSERT_VALUES, SCATTER_REVERSE);
227 auto check_solution = [&] {
230 boost::shared_ptr<VectorDouble> field_vals_ptr(
new VectorDouble());
231 boost::shared_ptr<MatrixDouble> diff_field_vals_ptr(
new MatrixDouble());
247 CHKERR assemble_matrices_and_vectors();
258 boost::shared_ptr<MatrixDouble> &diff_field_vals)
261 fieldVals(field_vals), diffFieldVals(diff_field_vals) {}
266 if (type == MBVERTEX) {
267 const int nb_gauss_pts = data.
getN().size2();
269 for (
int gg = 0; gg != nb_gauss_pts; ++gg) {
274 std::cout << f - (*fieldVals)[gg] <<
" : ";
275 for (
auto d : {0, 1, 2})
277 std::cout << std::endl;
279 constexpr double eps = 1e-6;
283 "Wrong value %6.4e != %6.4e (%6.4e)", f, (*
fieldVals)[gg],
284 f - (*fieldVals)[gg]);
286 for (
auto d : {0, 1, 2})
290 "Wrong diff value %6.4e != %6.4e (%6.4e)", diff_f[d],
308 const int nb_dofs = data.
getN().size2();
310 const int nb_gauss_pts = data.
getN().size1();
317 for (
int gg = 0; gg != nb_gauss_pts; ++gg) {
319 double v = t_w * vol * f;
320 double *val = &*nf.begin();
321 for (
int bb = 0; bb != nb_dofs; ++bb) {
347 const int row_nb_dofs = row_data.
getN().size2();
348 const int col_nb_dofs = col_data.
getN().size2();
349 if (row_nb_dofs && col_nb_dofs) {
350 const int nb_gauss_pts = row_data.
getN().size1();
355 double *row_base_ptr = &*row_data.
getN().data().begin();
356 double *col_base_ptr = &*col_data.
getN().data().begin();
357 for (
int gg = 0; gg != nb_gauss_pts; ++gg) {
359 double v = t_w * vol;
360 cblas_dger(CblasRowMajor, row_nb_dofs, col_nb_dofs,
v, row_base_ptr, 1,
361 col_base_ptr, 1, &*
m.data().begin(), col_nb_dofs);
363 row_base_ptr += row_nb_dofs;
364 col_base_ptr += col_nb_dofs;
ForcesAndSourcesCore::UserDataOperator UserDataOperator
#define CATCH_ERRORS
Catch errors.
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
#define MYPCOMM_INDEX
default communicator number PCOMM
#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
FTensor::Index< 'k', 3 > k
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.
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.
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....
default operator for Flat Prism 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.
Matrix manager is used to build and partition problems.
Calculate inverse of jacobian for face element.
Get field gradients at integration pts for scalar filed rank 0, i.e. vector field.
Get value at integration points for scalar field.
Operator for fat prism element updating integration weights in the volume.
Transform local reference derivatives of shape functions to global derivatives.
merge node from two bit levels
MoFEMErrorCode seedPrismsEntities(Range &prisms, const BitRefLevel &bit, int verb=-1)
Seed prism entities by bit level.
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.
double getVolume() const
element volume (linear geometry)
int getRuleThroughThickness(int order)
int getRuleTrianglesOnly(int order)
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 > diffFieldVals
boost::shared_ptr< VectorDouble > fieldVals
PrismOpCheck(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.
PrismOpLhs(SmartPetscObj< Mat > &a)
PrismOpRhs(SmartPetscObj< Vec > &f)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.