v0.14.0
field_blas.cpp

test field blas interface

/** \file field_blas.cpp
* \example field_blas.cpp
* \brief test field blas interface
*
* \ingroup mofem_field_algebra
*/
#include <MoFEM.hpp>
using namespace MoFEM;
static char help[] = "...\n\n";
int main(int argc, char *argv[]) {
MoFEM::Core::Initialize(&argc, &argv, PETSC_NULL, help);
try {
moab::Core mb_instance;
moab::Interface &moab = mb_instance;
int rank;
MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
// Reade parameters from line command
PetscBool flg = PETSC_TRUE;
char mesh_file_name[255];
#if PETSC_VERSION_GE(3, 6, 4)
CHKERR PetscOptionsGetString(PETSC_NULL, "", "-my_file", mesh_file_name,
255, &flg);
#else
CHKERR PetscOptionsGetString(PETSC_NULL, PETSC_NULL, "-my_file",
mesh_file_name, 255, &flg);
#endif
if (flg != PETSC_TRUE) {
SETERRQ(PETSC_COMM_SELF, 1, "*** ERROR -my_file (MESH FILE NEEDED)");
}
int order = 1;
#if PETSC_VERSION_GE(3, 6, 4)
CHKERR PetscOptionsGetInt(PETSC_NULL, "", "-my_order", &order, PETSC_NULL);
#else
CHKERR PetscOptionsGetInt(PETSC_NULL, PETSC_NULL, "-my_order", &order,
PETSC_NULL);
#endif
// Read mesh to MOAB
const char *option;
option = "";
CHKERR moab.load_file(mesh_file_name, 0, option);
// Create MoFEM database
// Note: Is MoFEM::CoreTmp<1> for testing purposes only
MoFEM::CoreTmp<-1> core(moab);
MoFEM::Interface &m_field = core;
CHKERR m_field.getInterface<BitRefManager>()->setBitRefLevelByDim(
0, 3, BitRefLevel().set(0));
// Create fields, add entities to field and set approximation order
CHKERR m_field.add_field("FIELD_A", H1, AINSWORTH_LEGENDRE_BASE, 3,
MB_TAG_DENSE);
CHKERR m_field.add_field("FIELD_B", H1, AINSWORTH_LEGENDRE_BASE, 3,
MB_TAG_DENSE);
CHKERR m_field.add_field("FIELD_C", H1, AINSWORTH_LEGENDRE_BASE, 1,
MB_TAG_DENSE);
CHKERR m_field.add_ents_to_field_by_type(0, MBTET, "FIELD_A");
CHKERR m_field.add_ents_to_field_by_type(0, MBTET, "FIELD_B");
CHKERR m_field.add_ents_to_field_by_type(0, MBTET, "FIELD_C");
CHKERR m_field.set_field_order(0, MBTET, "FIELD_A", order + 1);
CHKERR m_field.set_field_order(0, MBTRI, "FIELD_A", order + 1);
CHKERR m_field.set_field_order(0, MBEDGE, "FIELD_A", order + 1);
CHKERR m_field.set_field_order(0, MBVERTEX, "FIELD_A", 1);
CHKERR m_field.set_field_order(0, MBTET, "FIELD_B", order);
CHKERR m_field.set_field_order(0, MBTRI, "FIELD_B", order);
CHKERR m_field.set_field_order(0, MBEDGE, "FIELD_B", order);
CHKERR m_field.set_field_order(0, MBVERTEX, "FIELD_B", 1);
CHKERR m_field.set_field_order(0, MBTET, "FIELD_C", order);
CHKERR m_field.set_field_order(0, MBTRI, "FIELD_C", order);
CHKERR m_field.set_field_order(0, MBEDGE, "FIELD_C", order);
CHKERR m_field.set_field_order(0, MBVERTEX, "FIELD_C", 1);
// build field
CHKERR m_field.build_fields();
// get access to field agebra interface
auto fb = m_field.getInterface<FieldBlas>();
auto check_axpy = [&]() {
// set value to field
CHKERR fb->setField(+1, MBVERTEX, "FIELD_A");
CHKERR fb->setField(-2, MBVERTEX, "FIELD_B");
// FIELD_A = FIELD_A + 0.5 * FIELD_B
CHKERR fb->fieldAxpy(+0.5, "FIELD_B", "FIELD_A");
// FIELD_B *= -0.5
CHKERR fb->fieldScale(-0.5, "FIELD_B");
auto dofs_ptr = m_field.get_dofs();
for (auto dit : *dofs_ptr) {
auto check = [&](const std::string name, const double expected) {
if (dit->getName() == name) {
cout << name << " " << dit->getFieldData() << " " << expected
<< endl;
if (dit->getFieldData() != expected)
SETERRQ2(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
"Wrong DOF value 0 != %4.3e for %s", dit->getFieldData(),
boost::lexical_cast<std::string>(*dit).c_str());
}
};
CHKERR check("FIELD_A", 0);
CHKERR check("FIELD_B", 1);
}
};
auto check_set_vertex = [&]() {
auto set_distance = [&](VectorAdaptor &&field_data, double *xcoord,
double *ycoord, double *zcoord) {
xcoord, ycoord, zcoord);
field_data[0] = sqrt(t_coord(i) * t_coord(i));
};
CHKERR fb->setVertexDofs(set_distance, "FIELD_A");
// Test set veryex
struct TestMethod : EntityMethod {
TestMethod(moab::Interface &moab) : EntityMethod(), moab(moab) {}
MoFEMErrorCode preProcess() { return 0; }
MoFEMErrorCode operator()() {
if (entPtr->getEntType() == MBVERTEX) {
EntityHandle v = entPtr->getEnt();
VectorDouble3 coords(3);
CHKERR moab.get_coords(&v, 1, &coords[0]);
if (std::abs(entPtr->getEntFieldData()[0] - norm_2(coords)) >
std::numeric_limits<double>::epsilon())
SETERRQ2(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
"Wrong vals %3.4e != %3.4e",
entPtr->getEntFieldData()[0], norm_2(coords));
}
}
MoFEMErrorCode postProcess() { return 0; }
private:
};
TestMethod method(moab);
// checking if all is ok
CHKERR m_field.loop_entities("FIELD_A", method);
};
auto check_lambda = [&]() {
Range meshset_ents;
m_field, SIDESET | PRESSURESET, it)) {
if (it->getMeshsetId() != 1)
continue;
Range ents, nodes;
CHKERR it->getMeshsetIdEntitiesByType(m_field.get_moab(), MBTRI, ents,
true);
CHKERR m_field.get_moab().get_connectivity(ents, nodes, true);
meshset_ents.merge(nodes);
}
// set value to field
CHKERR fb->setField(+1, MBVERTEX, "FIELD_A");
CHKERR fb->setField(-1, MBVERTEX, "FIELD_B");
CHKERR fb->setField(1.5, MBVERTEX, "FIELD_C");
auto field_axpy = [&](const double val_y, const double val_x) {
return 0.5 * val_x + val_y;
};
auto field_scale = [&](const double val) { return -0.5 * val; };
// note that y_ent is first
// FIXME: this can be confusing?
auto vector_times_scalar_field =
[&](boost::shared_ptr<FieldEntity> ent_ptr_y,
boost::shared_ptr<FieldEntity> ent_ptr_x) {
auto x_data = ent_ptr_x->getEntFieldData();
auto y_data = ent_ptr_y->getEntFieldData();
// const auto size_x = x_data.size(); // scalar
const auto size_y = y_data.size(); // vector
for (size_t dd = 0; dd != size_y; ++dd)
y_data[dd] *= x_data[0];
// y_data *= x_data[0]; //would work as well
};
Range ents; // TODO: create test for subentities
// FIELD_A = FIELD_A + 0.5 * FIELD_B
CHKERR fb->fieldLambdaOnValues(field_axpy, "FIELD_B", "FIELD_A", true);
// CHKERR fb->fieldAxpy(+0.5, "FIELD_B", "FIELD_A");
// FIELD_B *= -0.5
// CHKERR fb->fieldScale(-0.5, "FIELD_B");
CHKERR fb->fieldLambdaOnValues(field_scale, "FIELD_B", &meshset_ents);
// FIELD_B(i) *= FIELD_C
CHKERR fb->fieldLambdaOnEntities(vector_times_scalar_field, "FIELD_C",
"FIELD_B", true, &meshset_ents);
auto dofs_ptr = m_field.get_dofs();
constexpr double eps = std::numeric_limits<double>::epsilon();
for (auto dit : *dofs_ptr) {
auto check = [&](const std::string name, const double expected) {
if (dit->getName() == name && dit->getEntType() == MBVERTEX) {
cout << name << " " << dit->getFieldData() << " " << expected
<< endl;
if (abs(dit->getFieldData() - expected) > eps)
SETERRQ3(PETSC_COMM_WORLD, MOFEM_ATOM_TEST_INVALID,
"Wrong DOF value %4.3e != %4.3e for %s",
dit->getFieldData(), expected,
boost::lexical_cast<std::string>(*dit).c_str());
}
};
CHKERR check("FIELD_A", 0.5);
if (meshset_ents.find(dit->getEnt()) != meshset_ents.end()) {
CHKERR check("FIELD_B", 0.75);
CHKERR check("FIELD_C", 1.5);
}
}
};
CHKERR check_axpy();
CHKERR check_set_vertex();
CHKERR check_lambda();
}
return 0;
}
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
SIDESET
@ SIDESET
Definition: definitions.h:160
H1
@ H1
continuous field
Definition: definitions.h:85
MoFEM::Types::VectorDouble3
VectorBoundedArray< double, 3 > VectorDouble3
Definition: Types.hpp:92
FTensor::Tensor1
Definition: Tensor1_value.hpp:8
EntityHandle
MoFEM::EntityMethod
Data structure to exchange data between mofem and User Loop Methods on entities.
Definition: LoopMethods.hpp:471
PRESSURESET
@ PRESSURESET
Definition: definitions.h:165
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
_IT_CUBITMESHSETS_BY_BCDATA_TYPE_FOR_LOOP_
#define _IT_CUBITMESHSETS_BY_BCDATA_TYPE_FOR_LOOP_(MESHSET_MANAGER, CUBITBCTYPE, IT)
Iterator that loops over a specific Cubit MeshSet in a moFEM field.
Definition: MeshsetsManager.hpp:49
MoFEM.hpp
MoFEM::CoreTmp< 0 >::Finalize
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
MoFEM::CoreInterface::add_ents_to_field_by_type
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.
order
constexpr int order
Definition: dg_projection.cpp:18
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
MoFEM::CoreInterface::get_dofs
virtual const DofEntity_multiIndex * get_dofs() const =0
Get the dofs object.
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:2010
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
MoFEM::CoreTmp
Definition: Core.hpp:36
MoFEM::Types::VectorAdaptor
VectorShallowArrayAdaptor< double > VectorAdaptor
Definition: Types.hpp:115
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
FTensor::Index< 'i', 3 >
mesh_file_name
char mesh_file_name[255]
Definition: mesh_smoothing.cpp:22
v
const double v
phase velocity of light in medium (cm/ns)
Definition: initial_diffusion.cpp:40
Range
MoFEM::CoreTmp< 0 >::Initialize
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
Definition: Core.cpp:72
FTensor::dd
const Tensor2_symmetric_Expr< const ddTensor0< T, Dim, i, j >, typename promote< T, double >::V, Dim, i, j > dd(const Tensor0< T * > &a, const Index< i, Dim > index1, const Index< j, Dim > index2, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
Definition: ddTensor0.hpp:33
help
static char help[]
Definition: field_blas.cpp:14
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:385
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1148
eps
static const double eps
Definition: check_base_functions_derivatives_on_tet.cpp:11
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
MoFEM::PetscOptionsGetString
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
Definition: DeprecatedPetsc.hpp:172
MoFEM::CoreInterface::build_fields
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
main
int main(int argc, char *argv[])
Definition: field_blas.cpp:16
MoFEM::BitRefManager
Managing BitRefLevels.
Definition: BitRefManager.hpp:21
MoFEM::Types::BitRefLevel
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:40
MOFEM_ATOM_TEST_INVALID
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h:40
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:453
MoFEM::CoreInterface::set_field_order
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.
MoFEM::PetscOptionsGetInt
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
Definition: DeprecatedPetsc.hpp:142
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
MoFEM::FieldBlas
Basic algebra on fields.
Definition: FieldBlas.hpp:21
MoFEM::CoreInterface::add_field
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.
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
MoFEM::CoreInterface::loop_entities
virtual MoFEMErrorCode loop_entities(const std::string field_name, EntityMethod &method, Range const *const ents=nullptr, int verb=DEFAULT_VERBOSITY)=0
Loop over field entities.