Calculate volume by integrating volume elements and using divergence theorem by integrating surface elements.
static char help[] =
"...\n\n";
struct OpVolume :
public VolumeElementForcesAndSourcesCore::UserDataOperator {
: VolumeElementForcesAndSourcesCore::UserDataOperator(
field_name, OPROW),
if (type != MBVERTEX)
const int nb_int_pts = getGaussPts().size2();
auto t_w = getFTensor0IntegrationWeight();
double vol = 0;
for (int gg = 0; gg != nb_int_pts; gg++) {
++t_w;
}
}
EntityType col_type,
}
};
struct OpFace :
public FaceElementForcesAndSourcesCore::UserDataOperator {
if (type != MBVERTEX)
double vol = 0;
for (int gg = 0; gg != nb_int_pts; gg++) {
vol += (t_coords(
i) * t_normal(
i)) * t_w;
++t_normal;
++t_w;
++t_coords;
}
vol /= 6;
}
EntityType col_type,
}
};
};
};
int main(
int argc,
char *argv[]) {
auto core_log = logging::core::get();
core_log->add_sink(
LogManager::createSink(LogManager::getStrmWorld(), "TEST"));
LogManager::setLog("TEST");
core_log->add_sink(
LogManager::createSink(LogManager::getStrmSync(), "TESTSYNC"));
LogManager::setLog("TESTSYNC");
try {
moab::Core moab_core;
moab::Interface &moab = moab_core;
DMType dm_name = "DMMOFEM";
{
auto dm = simple_interface->
getDM();
auto domain_fe =
boost::make_shared<VolumeElementForcesAndSourcesCore>(m_field);
auto boundary_fe =
boost::make_shared<FaceElementForcesAndSourcesCore>(m_field);
domain_fe->getRuleHook =
VolRule();
int ghosts[] = {0};
auto material_grad_mat = boost::make_shared<MatrixDouble>();
auto material_det_vec = boost::make_shared<VectorDouble>();
domain_fe->getOpPtrVector().push_back(
material_grad_mat));
domain_fe->getOpPtrVector().push_back(
domain_fe->getOpPtrVector().push_back(
domain_fe->getOpPtrVector().push_back(
domain_fe->getOpPtrVector().push_back(
new OpVolume(
"MESH_NODE_POSITIONS", vol));
boundary_fe->getOpPtrVector().push_back(
boundary_fe->getOpPtrVector().push_back(
new OpFace(
"MESH_NODE_POSITIONS", surf_vol));
domain_fe);
boundary_fe);
auto skeleton_fe = boost::make_shared<FEMethod>();
skeleton_fe->x = x;
skeleton_fe->x_t = x_t;
skeleton_fe->x_tt = x_tt;
skeleton_fe->preProcessHook = [&]() {
if (f != skeleton_fe->ts_F)
if (A != skeleton_fe->ts_A)
if (
B != skeleton_fe->ts_B)
if (x != skeleton_fe->ts_u)
if (x_t != skeleton_fe->ts_u_t)
if (x_tt != skeleton_fe->ts_u_tt)
};
skeleton_fe->postProcessHook = []() { return 0; };
skeleton_fe->operatorHook = []() { return 0; };
skeleton_fe);
CHKERR VecGhostUpdateBegin(vol, ADD_VALUES, SCATTER_REVERSE);
CHKERR VecGhostUpdateEnd(vol, ADD_VALUES, SCATTER_REVERSE);
CHKERR VecAssemblyBegin(surf_vol);
CHKERR VecAssemblyEnd(surf_vol);
CHKERR VecGhostUpdateBegin(surf_vol, ADD_VALUES, SCATTER_REVERSE);
CHKERR VecGhostUpdateEnd(surf_vol, ADD_VALUES, SCATTER_REVERSE);
double *a_vol;
CHKERR VecGetArray(vol, &a_vol);
double *a_surf_vol;
CHKERR VecGetArray(surf_vol, &a_surf_vol);
MOFEM_LOG(
"TEST", Sev::inform) <<
"Volume = " << a_vol[0];
MOFEM_LOG(
"TEST", Sev::inform) <<
"Surf Volume = " << a_surf_vol[0];
if (fabs(a_vol[0] - a_surf_vol[0]) > 1e-12) {
}
CHKERR VecRestoreArray(vol, &a_vol);
CHKERR VecRestoreArray(vol, &a_surf_vol);
}
using ForcesAndSourcesCore::ForcesAndSourcesCore;
const auto type = numeredEntFiniteElementPtr->getEntType();
if (type != MBENTITYSET) {
"Expected entity set as a finite element");
}
}
};
auto fe_ptr = boost::make_shared<MeshsetFE>(m_field);
struct OpMeshset : ForcesAndSourcesCore::UserDataOperator {
using OP = ForcesAndSourcesCore::UserDataOperator;
using OP::OP;
}
MoFEMErrorCode doWork(
int row_side,
int col_side, EntityType row_type,
EntityType col_type,
}
};
fe_ptr->getOpPtrVector().push_back(
new OpMeshset("GLOBAL", OpMeshset::OPROW));
fe_ptr->getOpPtrVector().push_back(
new OpMeshset("GLOBAL", "GLOBAL", OpMeshset::OPROWCOL));
}
}
return 0;
}
ForcesAndSourcesCore::UserDataOperator UserDataOperator
#define MOFEM_LOG_SEVERITY_SYNC(comm, severity)
Synchronise "SYNC" on curtain severity level.
#define CATCH_ERRORS
Catch errors.
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base nme:nme847.
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ NOFIELD
scalar or vector of scalars describe (no true field)
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_ATOM_TEST_INVALID
@ MOFEM_DATA_INCONSISTENCY
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
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.
auto createDMMatrix(DM dm)
Get smart matrix from DM.
PetscErrorCode DMoFEMLoopFiniteElementsUpAndLowRank(DM dm, const char fe_name[], MoFEM::FEMethod *method, int low_rank, int up_rank, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Executes FEMethod for finite elements in DM.
#define MOFEM_LOG(channel, severity)
Log.
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
virtual MoFEMErrorCode loop_dofs(const Problem *problem_ptr, const std::string &field_name, RowColData rc, DofMethod &method, int lower_rank, int upper_rank, int verb=DEFAULT_VERBOSITY)=0
Make a loop over dofs.
FTensor::Index< 'i', SPACE_DIM > i
const double v
phase velocity of light in medium (cm/ns)
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
implementation of Data Operators for Forces and Sources
SmartPetscObj< Vec > vectorDuplicate(Vec vec)
Create duplicate vector of smart vector.
auto createGhostVector(MPI_Comm comm, PetscInt n, PetscInt N, PetscInt nghost, const PetscInt ghosts[])
Create smart ghost vector.
constexpr auto field_name
Set integration rule to boundary elements.
int operator()(int, int, int) const
virtual int get_comm_size() const =0
virtual moab::Interface & get_moab()=0
virtual MPI_Comm & get_comm() const =0
virtual int get_comm_rank() const =0
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)
const VectorInt & getIndices() const
Get global indices of dofs on entity.
auto getFTensor1NormalsAtGaussPts()
get normal at integration points
auto getFTensor1CoordsAtGaussPts()
Get coordinates at integration points assuming linear geometry.
auto getFTensor0IntegrationWeight()
Get integration weights.
@ OPROW
operator doWork function is executed on FE rows
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element
structure to get information form mofem into EntitiesFieldData
Calculate HO coordinates at gauss points.
Set inverse jacobian to base functions.
Projection of edge entities with one mid-node on hierarchical basis.
Simple interface for fast problem set-up.
MoFEMErrorCode addDomainField(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_ZERO, int verb=-1)
Add field on domain.
const std::string getBoundaryFEName() const
Get the Boundary FE Name.
MoFEMErrorCode loadFile(const std::string options, const std::string mesh_file_name, LoadFileFunc loadFunc=defaultLoadFileFunc)
Load mesh file.
MoFEMErrorCode addMeshsetField(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_ZERO, int verb=-1)
Add meshset field.
const std::string getSkeletonFEName() const
Get the Skeleton FE Name.
MoFEMErrorCode addBoundaryField(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_ZERO, int verb=-1)
Add field on boundary.
MoFEMErrorCode getOptions()
get options
MoFEMErrorCode getDM(DM *dm)
Get DM.
std::string & getMeshsetFEName()
Get the Skeleton FE Name.
MoFEMErrorCode setFieldOrder(const std::string field_name, const int order, const Range *ents=NULL)
Set field order.
auto & getMeshsetFiniteElementEntities()
Get the Domain Fields.
MoFEMErrorCode setUp(const PetscBool is_partitioned=PETSC_TRUE)
Setup problem.
const std::string getDomainFEName() const
Get the Domain FE Name.
MoFEMErrorCode addSkeletonField(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_ZERO, int verb=-1)
Add field on skeleton.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
OpFace(const Range &skin_ents, const std::vector< unsigned char > &marker)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
OpVolume(const std::string &field_name, Vec vol)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Set integration rule to volume elements.
int operator()(int, int, int) const