v0.14.0
Public Member Functions | Private Attributes | List of all members
FreeSurfaceOps::OpRhsU Struct Reference

Rhs for U. More...

#include <tutorials/vec-5/src/FreeSurfaceOps.hpp>

Inheritance diagram for FreeSurfaceOps::OpRhsU:
[legend]
Collaboration diagram for FreeSurfaceOps::OpRhsU:
[legend]

Public Member Functions

 OpRhsU (const std::string field_name, boost::shared_ptr< MatrixDouble > dot_u_ptr, boost::shared_ptr< MatrixDouble > u_ptr, boost::shared_ptr< MatrixDouble > grad_u_ptr, boost::shared_ptr< VectorDouble > h_ptr, boost::shared_ptr< MatrixDouble > grad_h_ptr, boost::shared_ptr< VectorDouble > g_ptr, boost::shared_ptr< VectorDouble > p_ptr)
 
MoFEMErrorCode iNtegrate (EntitiesFieldData::EntData &data)
 
- Public Member Functions inherited from MoFEM::OpBaseImpl< A, EleOp >
 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...
 

Private Attributes

boost::shared_ptr< MatrixDouble > dotUPtr
 
boost::shared_ptr< MatrixDouble > uPtr
 
boost::shared_ptr< MatrixDouble > gradUPtr
 
boost::shared_ptr< VectorDouble > hPtr
 
boost::shared_ptr< MatrixDouble > gradHPtr
 
boost::shared_ptr< VectorDouble > gPtr
 
boost::shared_ptr< VectorDouble > pPtr
 

Additional Inherited Members

- Public Types inherited from MoFEM::OpBaseImpl< A, EleOp >
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)>
 
- Public Attributes inherited from MoFEM::OpBaseImpl< A, EleOp >
TimeFun timeScalingFun
 assumes that time variable is set More...
 
FEFun feScalingFun
 assumes that time variable is set More...
 
boost::shared_ptr< RangeentsPtr
 Entities on which element is run. More...
 
- Static Public Attributes inherited from MoFEM::OpBaseImpl< A, EleOp >
static MatSetValuesHook matSetValuesHook
 
- Protected Member Functions inherited from MoFEM::OpBaseImpl< A, EleOp >
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...
 
- Protected Attributes inherited from MoFEM::OpBaseImpl< A, EleOp >
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...
 

Detailed Description

Rhs for U.

Examples
free_surface.cpp.

Definition at line 356 of file FreeSurfaceOps.hpp.

Constructor & Destructor Documentation

◆ OpRhsU()

FreeSurfaceOps::OpRhsU::OpRhsU ( const std::string  field_name,
boost::shared_ptr< MatrixDouble >  dot_u_ptr,
boost::shared_ptr< MatrixDouble >  u_ptr,
boost::shared_ptr< MatrixDouble >  grad_u_ptr,
boost::shared_ptr< VectorDouble >  h_ptr,
boost::shared_ptr< MatrixDouble >  grad_h_ptr,
boost::shared_ptr< VectorDouble >  g_ptr,
boost::shared_ptr< VectorDouble >  p_ptr 
)
inline

Definition at line 358 of file FreeSurfaceOps.hpp.

366  : AssemblyDomainEleOp(field_name, field_name, AssemblyDomainEleOp::OPROW),
367  dotUPtr(dot_u_ptr), uPtr(u_ptr), gradUPtr(grad_u_ptr), hPtr(h_ptr),
368  gradHPtr(grad_h_ptr), gPtr(g_ptr), pPtr(p_ptr) {}

Member Function Documentation

◆ iNtegrate()

MoFEMErrorCode FreeSurfaceOps::OpRhsU::iNtegrate ( EntitiesFieldData::EntData data)
inline

Definition at line 370 of file FreeSurfaceOps.hpp.

370  {
372 
373  const double vol = getMeasure();
374  auto t_dot_u = getFTensor1FromMat<U_FIELD_DIM>(*dotUPtr);
375  auto t_u = getFTensor1FromMat<U_FIELD_DIM>(*uPtr);
376  auto t_p = getFTensor0FromVec(*pPtr);
377  auto t_grad_u = getFTensor2FromMat<U_FIELD_DIM, SPACE_DIM>(*gradUPtr);
378  auto t_h = getFTensor0FromVec(*hPtr);
379  auto t_grad_h = getFTensor1FromMat<SPACE_DIM>(*gradHPtr);
380  auto t_g = getFTensor0FromVec(*gPtr);
381  auto t_coords = getFTensor1CoordsAtGaussPts();
382 
383  auto t_base = data.getFTensor0N();
384  auto t_diff_base = data.getFTensor1DiffN<SPACE_DIM>();
385 
386  auto t_w = getFTensor0IntegrationWeight();
387 
396 
397  t_buoyancy(i) = 0;
398  t_gravity(i) = 0;
399 
400  for (int gg = 0; gg != nbIntegrationPts; gg++) {
401 
402  const double r = t_coords(0);
403  const double alpha = t_w * vol * cylindrical(r);
404 
405  const double rho = phase_function(t_h, rho_diff, rho_ave);
406  const double mu = phase_function(t_h, mu_diff, mu_ave);
407 
408  auto t_D = get_D(2 * mu);
409 
410  t_inertia_force(i) = (rho * alpha) * (t_dot_u(i));
411  // t_buoyancy(SPACE_DIM - 1) = -(alpha * rho * a0) * t_h;
412  t_gravity(SPACE_DIM - 1) = (alpha * rho * a0);
413  t_phase_force(i) = -alpha * kappa * t_g * t_grad_h(i);
414  t_convection(i) = (rho * alpha) * (t_u(j) * t_grad_u(i, j));
415 
416  t_stress(i, j) =
417  alpha * (t_D(i, j, k, l) * t_grad_u(k, l) + t_kd(i, j) * t_p);
418 
419  auto t_nf = getFTensor1FromArray<U_FIELD_DIM, U_FIELD_DIM>(locF);
420 
421  t_forces(i) = t_inertia_force(i) + t_buoyancy(i) + t_gravity(i) +
422  t_convection(i) + t_phase_force(i);
423 
424  int bb = 0;
425  for (; bb != nbRows / U_FIELD_DIM; ++bb) {
426 
427  t_nf(i) += t_base * t_forces(i);
428  t_nf(i) += t_diff_base(j) * t_stress(i, j);
429 
430  if (coord_type == CYLINDRICAL) {
431  t_nf(0) += (t_base * (alpha / t_coords(0))) * (2 * mu * t_u(0) + t_p);
432  }
433 
434  ++t_base;
435  ++t_diff_base;
436  ++t_nf;
437  }
438 
439  for (; bb < nbRowBaseFunctions; ++bb) {
440  ++t_diff_base;
441  ++t_base;
442  }
443 
444  ++t_dot_u;
445  ++t_u;
446  ++t_grad_u;
447  ++t_h;
448  ++t_grad_h;
449  ++t_g;
450  ++t_p;
451 
452  ++t_w;
453  ++t_coords;
454  }
455 
457  }

Member Data Documentation

◆ dotUPtr

boost::shared_ptr<MatrixDouble> FreeSurfaceOps::OpRhsU::dotUPtr
private

Definition at line 460 of file FreeSurfaceOps.hpp.

◆ gPtr

boost::shared_ptr<VectorDouble> FreeSurfaceOps::OpRhsU::gPtr
private

Definition at line 465 of file FreeSurfaceOps.hpp.

◆ gradHPtr

boost::shared_ptr<MatrixDouble> FreeSurfaceOps::OpRhsU::gradHPtr
private

Definition at line 464 of file FreeSurfaceOps.hpp.

◆ gradUPtr

boost::shared_ptr<MatrixDouble> FreeSurfaceOps::OpRhsU::gradUPtr
private

Definition at line 462 of file FreeSurfaceOps.hpp.

◆ hPtr

boost::shared_ptr<VectorDouble> FreeSurfaceOps::OpRhsU::hPtr
private

Definition at line 463 of file FreeSurfaceOps.hpp.

◆ pPtr

boost::shared_ptr<VectorDouble> FreeSurfaceOps::OpRhsU::pPtr
private

Definition at line 466 of file FreeSurfaceOps.hpp.

◆ uPtr

boost::shared_ptr<MatrixDouble> FreeSurfaceOps::OpRhsU::uPtr
private

Definition at line 461 of file FreeSurfaceOps.hpp.


The documentation for this struct was generated from the following file:
get_D
auto get_D
Definition: free_surface.cpp:252
rho_ave
double rho_ave
Definition: free_surface.cpp:178
MoFEM::OpBaseImpl::nbIntegrationPts
int nbIntegrationPts
number of integration points
Definition: FormsIntegrators.hpp:238
FTensor::Tensor1
Definition: Tensor1_value.hpp:8
FreeSurfaceOps::OpRhsU::hPtr
boost::shared_ptr< VectorDouble > hPtr
Definition: FreeSurfaceOps.hpp:463
rho
double rho
Definition: plastic.cpp:140
FreeSurfaceOps::OpRhsU::gPtr
boost::shared_ptr< VectorDouble > gPtr
Definition: FreeSurfaceOps.hpp:465
FTensor::Tensor2_symmetric
Definition: Tensor2_symmetric_value.hpp:13
sdf.r
int r
Definition: sdf.py:8
FreeSurfaceOps::OpRhsU::gradHPtr
boost::shared_ptr< MatrixDouble > gradHPtr
Definition: FreeSurfaceOps.hpp:464
cylindrical
auto cylindrical
Definition: free_surface.cpp:187
MoFEM::OpBaseImpl::locF
VectorDouble locF
local entity vector
Definition: FormsIntegrators.hpp:251
U_FIELD_DIM
constexpr int U_FIELD_DIM
Definition: free_surface.cpp:79
SPACE_DIM
constexpr int SPACE_DIM
Definition: child_and_parent.cpp:16
rho_diff
double rho_diff
Definition: free_surface.cpp:179
MoFEM::getFTensor0FromVec
static auto getFTensor0FromVec(ublas::vector< T, A > &data)
Get tensor rank 0 (scalar) form data vector.
Definition: Templates.hpp:135
FreeSurfaceOps::OpRhsU::gradUPtr
boost::shared_ptr< MatrixDouble > gradUPtr
Definition: FreeSurfaceOps.hpp:462
phase_function
auto phase_function
Definition: free_surface.cpp:211
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
t_kd
constexpr auto t_kd
Definition: free_surface.cpp:137
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
mu_diff
double mu_diff
Definition: free_surface.cpp:181
MoFEM::OpBaseImpl::nbRowBaseFunctions
int nbRowBaseFunctions
number or row base functions
Definition: FormsIntegrators.hpp:239
kappa
double kappa
Definition: free_surface.cpp:183
MoFEM::OpBaseImpl::nbRows
int nbRows
number of dofs on rows
Definition: FormsIntegrators.hpp:236
CYLINDRICAL
@ CYLINDRICAL
Definition: definitions.h:130
mu_ave
double mu_ave
Definition: free_surface.cpp:180
FreeSurfaceOps::OpRhsU::pPtr
boost::shared_ptr< VectorDouble > pPtr
Definition: FreeSurfaceOps.hpp:466
j
FTensor::Index< 'j', 3 > j
Definition: matrix_function.cpp:19
coord_type
constexpr CoordinateTypes coord_type
Definition: incompressible_elasticity.cpp:19
mu
double mu
Definition: dynamic_first_order_con_law.cpp:98
a0
constexpr double a0
Definition: hcurl_check_approx_in_2d.cpp:37
FreeSurfaceOps::OpRhsU::uPtr
boost::shared_ptr< MatrixDouble > uPtr
Definition: FreeSurfaceOps.hpp:461
FTensor::Kronecker_Delta_symmetric
Kronecker Delta class symmetric.
Definition: Kronecker_Delta.hpp:49
AssemblyDomainEleOp
FormsIntegrators< DomainEleOp >::Assembly< A >::OpBase AssemblyDomainEleOp
Definition: tensor_divergence_operator.cpp:59
k
FTensor::Index< 'k', 3 > k
Definition: matrix_function.cpp:20
FreeSurfaceOps::OpRhsU::dotUPtr
boost::shared_ptr< MatrixDouble > dotUPtr
Definition: FreeSurfaceOps.hpp:460
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
l
FTensor::Index< 'l', 3 > l
Definition: matrix_function.cpp:21