#include <tutorials/vec-5/src/FreeSurfaceOps.hpp>
|
| OpRhsH (const std::string field_name, boost::shared_ptr< MatrixDouble > u_ptr, boost::shared_ptr< VectorDouble > dot_h_ptr, boost::shared_ptr< VectorDouble > h_ptr, boost::shared_ptr< MatrixDouble > grad_h_ptr, boost::shared_ptr< MatrixDouble > grad_g_ptr) |
|
MoFEMErrorCode | iNtegrate (EntitiesFieldData::EntData &data) |
|
| OpBaseImpl (const std::string row_field_name, const std::string col_field_name, const OpType type, boost::shared_ptr< Range > ents_ptr=nullptr) |
|
MoFEMErrorCode | doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, EntData &row_data, EntData &col_data) |
| Do calculations for the left hand side. More...
|
|
MoFEMErrorCode | doWork (int row_side, EntityType row_type, EntData &row_data) |
| Do calculations for the right hand side. More...
|
|
|
boost::shared_ptr< MatrixDouble > | uPtr |
|
boost::shared_ptr< VectorDouble > | dotHPtr |
|
boost::shared_ptr< VectorDouble > | hPtr |
|
boost::shared_ptr< MatrixDouble > | gradHPtr |
|
boost::shared_ptr< MatrixDouble > | gradGPtr |
|
|
using | OpType = typename EleOp::OpType |
|
using | EntData = EntitiesFieldData::EntData |
|
using | MatSetValuesHook = boost::function< MoFEMErrorCode(ForcesAndSourcesCore::UserDataOperator *op_ptr, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, MatrixDouble &m)> |
|
TimeFun | timeScalingFun |
| assumes that time variable is set More...
|
|
FEFun | feScalingFun |
| assumes that time variable is set More...
|
|
boost::shared_ptr< Range > | entsPtr |
| Entities on which element is run. More...
|
|
static MatSetValuesHook | matSetValuesHook |
|
template<int DIM> |
FTensor::Tensor1< FTensor::PackPtr< double *, DIM >, DIM > | getNf () |
|
template<int DIM> |
FTensor::Tensor2< FTensor::PackPtr< double *, DIM >, DIM, DIM > | getLocMat (const int rr) |
|
virtual MoFEMErrorCode | iNtegrate (EntData &row_data, EntData &col_data) |
| Integrate grad-grad operator. More...
|
|
virtual MoFEMErrorCode | aSsemble (EntData &row_data, EntData &col_data, const bool trans) |
|
virtual MoFEMErrorCode | iNtegrate (EntData &data) |
| Class dedicated to integrate operator. More...
|
|
virtual MoFEMErrorCode | aSsemble (EntData &data) |
|
virtual size_t | getNbOfBaseFunctions (EntitiesFieldData::EntData &data) |
| Get number of base functions. More...
|
|
int | nbRows |
| number of dofs on rows More...
|
|
int | nbCols |
| number if dof on column More...
|
|
int | nbIntegrationPts |
| number of integration points More...
|
|
int | nbRowBaseFunctions |
| number or row base functions More...
|
|
int | rowSide |
| row side number More...
|
|
int | colSide |
| column side number More...
|
|
EntityType | rowType |
| row type More...
|
|
EntityType | colType |
| column type More...
|
|
bool | assembleTranspose |
|
bool | onlyTranspose |
|
MatrixDouble | locMat |
| local entity block matrix More...
|
|
MatrixDouble | locMatTranspose |
| local entity block matrix More...
|
|
VectorDouble | locF |
| local entity vector More...
|
|
template<bool I>
struct FreeSurfaceOps::OpRhsH< I >
- Examples
- free_surface.cpp.
Definition at line 796 of file FreeSurfaceOps.hpp.
◆ OpRhsH()
template<bool I>
FreeSurfaceOps::OpRhsH< I >::OpRhsH |
( |
const std::string |
field_name, |
|
|
boost::shared_ptr< MatrixDouble > |
u_ptr, |
|
|
boost::shared_ptr< VectorDouble > |
dot_h_ptr, |
|
|
boost::shared_ptr< VectorDouble > |
h_ptr, |
|
|
boost::shared_ptr< MatrixDouble > |
grad_h_ptr, |
|
|
boost::shared_ptr< MatrixDouble > |
grad_g_ptr |
|
) |
| |
|
inline |
◆ iNtegrate()
Definition at line 807 of file FreeSurfaceOps.hpp.
810 const double vol = getMeasure();
811 auto t_w = getFTensor0IntegrationWeight();
812 auto t_coords = getFTensor1CoordsAtGaussPts();
813 auto t_base = data.getFTensor0N();
814 auto t_diff_base = data.getFTensor1DiffN<
SPACE_DIM>();
817 if (data.getDiffN().size1() != data.getN().size1())
819 if (data.getDiffN().size2() != data.getN().size2() *
SPACE_DIM) {
822 MOFEM_LOG(
"SELF", Sev::error) << data.getN();
823 MOFEM_LOG(
"SELF", Sev::error) << data.getDiffN();
831 auto t_grad_g = getFTensor1FromMat<SPACE_DIM>(*
gradGPtr);
835 const double r = t_coords(0);
838 const double set_h =
init_h(t_coords(0), t_coords(1), t_coords(2));
839 const double m =
get_M(set_h) * alpha;
842 for (; bb !=
nbRows; ++bb) {
843 locF[bb] += (t_base * alpha) * (t_h - set_h);
844 locF[bb] += (t_diff_base(
i) *
m) * t_grad_g(
i);
865 auto t_u = getFTensor1FromMat<U_FIELD_DIM>(*
uPtr);
866 auto t_grad_h = getFTensor1FromMat<SPACE_DIM>(*
gradHPtr);
867 auto t_grad_g = getFTensor1FromMat<SPACE_DIM>(*
gradGPtr);
871 const double r = t_coords(0);
874 const double m =
get_M(t_h) * alpha;
877 for (; bb !=
nbRows; ++bb) {
878 locF[bb] += (t_base * alpha) * (t_dot_h);
879 locF[bb] += (t_base * alpha) * (t_grad_h(
i) * t_u(
i));
880 locF[bb] += (t_diff_base(
i) * t_grad_g(
i)) *
m;
◆ dotHPtr
◆ gradGPtr
◆ gradHPtr
◆ hPtr
◆ uPtr
The documentation for this struct was generated from the following file: