16static char help[] =
"...\n\n";
18struct OpVolume :
public VolumeElementForcesAndSourcesCore::UserDataOperator {
29 const int nb_int_pts = getGaussPts().size2();
31 auto t_w = getFTensor0IntegrationWeight();
32 double v = getMeasure();
34 for (
int gg = 0; gg != nb_int_pts; gg++) {
38 CHKERR VecSetValue(
vOl, 0, vol, ADD_VALUES);
51struct OpFace :
public FaceElementForcesAndSourcesCore::UserDataOperator {
68 for (
int gg = 0; gg != nb_int_pts; gg++) {
69 vol += (t_coords(
i) * t_normal(
i)) * t_w;
75 CHKERR VecSetValue(
vOl, 0, vol, ADD_VALUES);
95int main(
int argc,
char *argv[]) {
100 auto core_log = logging::core::get();
113 moab::Core moab_core;
114 moab::Interface &moab = moab_core;
121 DMType dm_name =
"DMMOFEM";
150 CHKERR m_field.
get_moab().get_entities_by_type(0, MBVERTEX, verts);
164 auto dm = simple_interface->
getDM();
167 boost::make_shared<VolumeElementForcesAndSourcesCore>(m_field);
169 boost::make_shared<FaceElementForcesAndSourcesCore>(m_field);
172 domain_fe->getRuleHook =
VolRule();
173 boundary_fe->getRuleHook =
FaceRule();
183 auto material_grad_mat = boost::make_shared<MatrixDouble>();
184 auto material_det_vec = boost::make_shared<VectorDouble>();
186 domain_fe->getOpPtrVector().push_back(
189 domain_fe->getOpPtrVector().push_back(
191 domain_fe->getOpPtrVector().push_back(
193 domain_fe->getOpPtrVector().push_back(
195 domain_fe->getOpPtrVector().push_back(
196 new OpVolume(
"MESH_NODE_POSITIONS", vol));
198 boundary_fe->getOpPtrVector().push_back(
200 boundary_fe->getOpPtrVector().push_back(
201 new OpFace(
"MESH_NODE_POSITIONS", surf_vol));
209 auto skeleton_fe = boost::make_shared<FEMethod>();
220 skeleton_fe->x_t = x_t;
221 skeleton_fe->x_tt = x_tt;
223 skeleton_fe->preProcessHook = [&]() {
225 if (f != skeleton_fe->ts_F)
227 if (A != skeleton_fe->ts_A)
229 if (
B != skeleton_fe->ts_B)
231 if (x != skeleton_fe->ts_u)
233 if (x_t != skeleton_fe->ts_u_t)
235 if (x_tt != skeleton_fe->ts_u_tt)
240 skeleton_fe->postProcessHook = []() {
return 0; };
241 skeleton_fe->operatorHook = []() {
return 0; };
247 CHKERR VecAssemblyBegin(vol);
248 CHKERR VecAssemblyEnd(vol);
249 CHKERR VecGhostUpdateBegin(vol, ADD_VALUES, SCATTER_REVERSE);
250 CHKERR VecGhostUpdateEnd(vol, ADD_VALUES, SCATTER_REVERSE);
251 CHKERR VecAssemblyBegin(surf_vol);
252 CHKERR VecAssemblyEnd(surf_vol);
253 CHKERR VecGhostUpdateBegin(surf_vol, ADD_VALUES, SCATTER_REVERSE);
254 CHKERR VecGhostUpdateEnd(surf_vol, ADD_VALUES, SCATTER_REVERSE);
257 CHKERR VecGetArray(vol, &a_vol);
259 CHKERR VecGetArray(surf_vol, &a_surf_vol);
260 MOFEM_LOG(
"TEST", Sev::inform) <<
"Volume = " << a_vol[0];
261 MOFEM_LOG(
"TEST", Sev::inform) <<
"Surf Volume = " << a_surf_vol[0];
262 if (fabs(a_vol[0] - a_surf_vol[0]) > 1e-12) {
265 CHKERR VecRestoreArray(vol, &a_vol);
266 CHKERR VecRestoreArray(vol, &a_surf_vol);
270 using ForcesAndSourcesCore::ForcesAndSourcesCore;
273 const auto type = numeredEntFiniteElementPtr->getEntType();
274 if (type != MBENTITYSET) {
276 "Expected entity set as a finite element");
285 auto fe_ptr = boost::make_shared<MeshsetFE>(m_field);
287 struct OpMeshset : ForcesAndSourcesCore::UserDataOperator {
288 using OP = ForcesAndSourcesCore::UserDataOperator;
298 MoFEMErrorCode doWork(
int row_side,
int col_side, EntityType row_type,
309 fe_ptr->getOpPtrVector().push_back(
310 new OpMeshset(
"GLOBAL", OpMeshset::OPROW));
311 fe_ptr->getOpPtrVector().push_back(
312 new OpMeshset(
"GLOBAL",
"GLOBAL", OpMeshset::OPROWCOL));
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.
static LoggerType & setLog(const std::string channel)
Set ans resset chanel logger.
#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
MoFEMErrorCode loopOverOperators()
Iterate user data operators.
virtual MoFEMErrorCode operator()()
function is run for every finite element
virtual MoFEMErrorCode preProcess()
function is run at the beginning of loop
virtual MoFEMErrorCode postProcess()
function is run at the end of loop
static boost::shared_ptr< SinkType > createSink(boost::shared_ptr< std::ostream > stream_ptr, std::string comm_filter)
Create a sink object.
static boost::shared_ptr< std::ostream > getStrmWorld()
Get the strm world object.
static boost::shared_ptr< std::ostream > getStrmSync()
Get the strm sync object.
Calculate HO coordinates at gauss points.
Get field gradients at integration pts for scalar filed rank 0, i.e. vector field.
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.
Volume finite element base.
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
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.
OpFace(const std::string &field_name, Vec vol)
MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
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