Loading [MathJax]/extensions/tex2jax.js
OGS
VariableDependentNeumannBoundaryConditionLocalAssembler.h
Go to the documentation of this file.
1 
11 #pragma once
12 
14 #include "MeshLib/PropertyVector.h"
18 
19 namespace ProcessLib
20 {
22 {
27  // Used for mapping boundary nodes to bulk nodes.
29 };
30 
31 template <typename ShapeFunction, typename IntegrationMethod, int GlobalDim>
34  ShapeFunction, IntegrationMethod, GlobalDim>
35 {
37  ShapeFunction, IntegrationMethod, GlobalDim>;
40 
41 public:
45  MeshLib::Element const& e,
46  std::size_t const local_matrix_size,
47  bool const is_axially_symmetric,
48  unsigned const integration_order,
50  : Base(e, is_axially_symmetric, integration_order),
51  _data(data),
52  _local_matrix_size(local_matrix_size)
53  {
54  }
55 
56  void assemble(std::size_t const mesh_item_id,
57  NumLib::LocalToGlobalIndexMap const& dof_table_boundary,
58  double const t, std::vector<GlobalVector*> const& x,
59  int const process_id, GlobalMatrix& /*K*/, GlobalVector& b,
60  GlobalMatrix* /*Jac*/) override
61  {
63  _local_rhs.setZero();
64  // Get element nodes for the interpolation from nodes to
65  // integration point.
66  NodalVectorType const constant_node_values =
68  .template topRows<ShapeFunction::MeshElement::n_all_nodes>();
69  NodalVectorType const coefficient_current_variable_node_values =
72  .template topRows<ShapeFunction::MeshElement::n_all_nodes>();
73  NodalVectorType const coefficient_other_variable_node_values =
76  .template topRows<ShapeFunction::MeshElement::n_all_nodes>();
77  NodalVectorType const coefficient_mixed_variables_node_values =
80  .template topRows<ShapeFunction::MeshElement::n_all_nodes>();
81  unsigned const n_integration_points =
82  Base::_integration_method.getNumberOfPoints();
83 
84  auto const indices_current_variable =
85  NumLib::getIndices(mesh_item_id, dof_table_boundary);
86  auto const indices_other_variable = NumLib::getIndices(
88  std::vector<double> const local_current_variable =
89  x[process_id]->get(indices_current_variable);
90  std::vector<double> const local_other_variable =
91  x[process_id]->get(indices_other_variable);
92 
93  for (unsigned ip = 0; ip < n_integration_points; ip++)
94  {
95  auto const& n_and_weight = Base::_ns_and_weights[ip];
96  auto const& N = n_and_weight.N;
97  auto const& w = n_and_weight.weight;
98 
99  double current_variable_int_pt = 0.0;
100  double other_variable_int_pt = 0.0;
101 
102  NumLib::shapeFunctionInterpolate(local_current_variable, N,
103  current_variable_int_pt);
104  NumLib::shapeFunctionInterpolate(local_other_variable, N,
105  other_variable_int_pt);
106  NodalVectorType const neumann_node_values =
107  constant_node_values +
108  coefficient_current_variable_node_values *
109  current_variable_int_pt +
110  coefficient_other_variable_node_values * other_variable_int_pt +
111  coefficient_mixed_variables_node_values *
112  current_variable_int_pt * other_variable_int_pt;
113  _local_rhs.noalias() += N * neumann_node_values.dot(N) * w;
114  }
115 
116  b.add(indices_current_variable, _local_rhs);
117  }
118 
119 private:
121  std::size_t const _local_matrix_size;
122 };
123 
124 } // namespace ProcessLib
Global vector based on Eigen vector.
Definition: EigenVector.h:26
void add(IndexType rowId, double v)
add entry
Definition: EigenVector.h:80
std::vector< NAndWeight, Eigen::aligned_allocator< NAndWeight > > const _ns_and_weights
VariableDependentNeumannBoundaryConditionLocalAssembler(MeshLib::Element const &e, std::size_t const local_matrix_size, bool const is_axially_symmetric, unsigned const integration_order, VariableDependentNeumannBoundaryConditionData const &data)
void assemble(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table_boundary, double const t, std::vector< GlobalVector * > const &x, int const process_id, GlobalMatrix &, GlobalVector &b, GlobalMatrix *) override
void shapeFunctionInterpolate(const NodalValues &nodal_values, const ShapeMatrix &shape_matrix_N, double &interpolated_value, ScalarTypes &... interpolated_values)
Definition: Interpolation.h:79
std::vector< GlobalIndexType > getIndices(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
virtual Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > getNodalValuesOnElement(MeshLib::Element const &element, double const t) const
Returns a matrix of values for all nodes of the given element.
Definition: Parameter.h:163