32template <
typename BcOrStData,
typename ShapeFunction,
33 typename LowerOrderShapeFunction,
int GlobalDim>
47 bool is_axially_symmetric,
54 GlobalDim>(e, is_axially_symmetric,
61 void assemble(std::size_t
const boundary_element_id,
67 Eigen::MatrixXd
const primary_variables_mat =
ProcessLib::
70 dof_table_boundary, x);
73 auto const& indices_this_component =
74 dof_table_boundary(boundary_element_id,
77 std::vector<GlobalIndexType> indices_all_components_for_Jac;
79 auto const num_dof_this_component = indices_this_component.
size();
81 Eigen::VectorXd local_rhs =
82 Eigen::VectorXd::Zero(num_dof_this_component);
83 Eigen::MatrixXd local_Jac;
87 indices_all_components_for_Jac =
90 auto const num_dof_all_components =
91 indices_all_components_for_Jac.size();
93 local_Jac = Eigen::MatrixXd::Zero(num_dof_this_component,
94 num_dof_all_components);
98 auto const num_comp_total =
100 std::vector<double> prim_vars_data(num_comp_total);
104 unsigned const num_integration_points =
109 for (
unsigned ip = 0; ip < num_integration_points; ip++)
115 primary_variables_mat,
119 auto const [flag, flux, dFlux] =
120 bc_or_st_data.getFlagAndFluxAndDFlux(t, coords, prim_vars_data);
131 if (
static_cast<int>(dFlux.size()) != num_comp_total)
137 "The Python BC or ST must return the derivative of the "
138 "flux w.r.t. each component of each primary variable. "
139 "{:d} components expected. {:d} components returned from "
141 "The expected number of components depends on multiple "
142 "factors. In the case of vectorial primary variables (such "
143 "as displacement), it will depend on the space dimension. "
144 "In the case of staggered coupling, the number of "
145 "components can be different for each staggered process "
146 "and different from monolithic coupling.",
147 num_comp_total, dFlux.size());
156 b.
add(indices_this_component, local_rhs);
161 indices_this_component, indices_all_components_for_Jac};
162 Jac->
add(rci, local_Jac);
185 Eigen::VectorXd& local_rhs)
const
187 auto const w = ns_and_weight.
weight();
190 local_rhs.noalias() += N.transpose() * (flux * w);
206 std::vector<double>
const& dFlux,
208 Eigen::MatrixXd& local_Jac)
const
210 auto const& pv_refs =
212 auto const weight = ns_and_weight.
weight();
215 auto const N_this = ns_and_weight.
N(this_shp_fct_order);
217 Eigen::Index left = 0;
219 for (
auto pv_ref : pv_refs)
221 auto const& pv = pv_ref.get();
222 auto const other_num_comp = pv.getNumberOfGlobalComponents();
223 auto const other_shp_fct_order = pv.getShapeFunctionOrder();
224 auto const N_other = ns_and_weight.
N(other_shp_fct_order);
225 auto const width = N_other.size();
227 for (
decltype(+other_num_comp) other_comp = 0;
228 other_comp < other_num_comp;
233 local_Jac.middleCols(left, width).noalias() -=
234 N_this.transpose() * (dFlux[other_comp] * weight) * N_other;
240 assert(left == local_Jac.cols());
int add(IndexType row, IndexType col, double val)
Global vector based on Eigen vector.
void add(IndexType rowId, double v)
add entry
std::size_t getID() const
Get id of the mesh.
unsigned getNumberOfPoints() const
int getNumberOfGlobalComponents() const
Template class for isoparametric elements.
std::array< double, 3 > interpolateCoordinates(typename ShapeMatrices::ShapeType const &N) const
Interpolates the coordinates of the element with the given shape matrix.
Eigen::Map< const Vector > toVector(std::vector< double > const &data, Eigen::VectorXd::Index size)
Creates an Eigen mapped vector from the given data vector.
std::vector< GlobalIndexType > getIndices(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
NumLib::TemplateIsoparametric< ShapeFunction, ShapeMatricesType > createIsoparametricFiniteElement(MeshLib::Element const &e)
auto computeNsAndWeights(MeshLib::Element const &element, bool const is_axially_symmetric, IntegrationMethod const &integration_method)
Eigen::MatrixXd collectDofsToMatrix(MeshLib::Element const &element, std::size_t const mesh_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x)
void interpolate(Eigen::MatrixXd const &primary_variables_mat, std::vector< std::reference_wrapper< ProcessVariable > > const &pv_refs, NsAndWeight const &ns_and_weight, Eigen::Ref< Eigen::VectorXd > interpolated_primary_variables)
void assembleLocalRhs(double const flux, typename Traits::NsAndWeight const &ns_and_weight, Eigen::VectorXd &local_rhs) const
void assemble(std::size_t const boundary_element_id, NumLib::LocalToGlobalIndexMap const &dof_table_boundary, double const t, GlobalVector const &x, GlobalVector &b, GlobalMatrix *const Jac) const
BcAndStLocalAssemblerImpl(MeshLib::Element const &e, NumLib::GenericIntegrationMethod const &integration_method_, bool is_axially_symmetric, BcOrStData const &data)
std::array< double, 3 > interpolateCoords(typename Traits::NsAndWeight const &ns_and_weight, NumLib::TemplateIsoparametric< ShapeFunction, ShapeMatrixPolicy > const &fe) const
void assembleLocalJacobian(std::vector< double > const &dFlux, typename Traits::NsAndWeight const &ns_and_weight, Eigen::MatrixXd &local_Jac) const
typename Traits::ShapeMatrixPolicy ShapeMatrixPolicy
MeshLib::Element const & element
BcOrStData const & bc_or_st_data
NumLib::GenericIntegrationMethod const & integration_method
std::vector< typename Traits::NsAndWeight > const nss_and_weights
MeshLib::Mesh const & bc_or_st_mesh
The domain, where this BC or ST will be applied.
std::vector< std::reference_wrapper< ProcessVariable > > const & all_process_variables_for_this_process
int const global_component_id
unsigned const shape_function_order
ShapeMatrix const & NHigherOrder() const
Eigen::Ref< const Eigen::RowVectorXd > N(unsigned order) const
Collects common type aliases needed when working with NsAndWeight.
ShapeMatrixPolicyType< ShapeFunction, GlobalDim > ShapeMatrixPolicy