template<typename ShapeFunction, typename IntegrationMethod, int GlobalDim>
class ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >
Definition at line 49 of file SteadyStateDiffusionFEM.h.
|
| LocalAssemblerData (MeshLib::Element const &element, std::size_t const, bool is_axially_symmetric, unsigned const integration_order, SteadyStateDiffusionData const &process_data) |
|
void | assemble (double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &, std::vector< double > &, std::vector< double > &local_K_data, std::vector< double > &) override |
|
Eigen::Vector3d | getFlux (MathLib::Point3d const &p_local_coords, double const t, std::vector< double > const &local_x) const override |
|
Eigen::Map< const Eigen::RowVectorXd > | getShapeMatrix (const unsigned integration_point) const override |
| Provides the shape matrix at the given integration point. More...
|
|
std::vector< double > const & | getIntPtDarcyVelocity (const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override |
|
virtual | ~LocalAssemblerInterface ()=default |
|
void | setInitialConditions (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, double const t, bool const use_monolithic_scheme, int const process_id) |
|
virtual void | initialize (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table) |
|
virtual void | preAssemble (double const, double const, std::vector< double > const &) |
|
virtual void | assembleForStaggeredScheme (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_xdot, int const process_id, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data) |
|
virtual void | assembleWithJacobian (double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_xdot, const double dxdot_dx, const double dx_dx, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data) |
|
virtual void | assembleWithJacobianForStaggeredScheme (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_xdot, const double dxdot_dx, const double dx_dx, int const process_id, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data) |
|
virtual void | computeSecondaryVariable (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_dot, int const process_id) |
|
virtual void | preTimestep (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, double const t, double const delta_t) |
|
virtual void | postTimestep (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, double const dt) |
|
void | postNonLinearSolver (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, GlobalVector const &xdot, double const t, double const dt, bool const use_monolithic_scheme, int const process_id) |
|
virtual std::vector< double > | interpolateNodalValuesToIntegrationPoints (std::vector< double > const &) |
|
virtual Eigen::Vector3d | getFlux (MathLib::Point3d const &, double const, std::vector< std::vector< double >> const &) const |
| Fits to staggered scheme. More...
|
|
template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
Computes the flux in the point p_local_coords
that is given in local coordinates using the values from local_x
.
Reimplemented from ProcessLib::LocalAssemblerInterface.
Definition at line 130 of file SteadyStateDiffusionFEM.h.
136 double const dt = std::numeric_limits<double>::quiet_NaN();
141 auto const shape_matrices =
145 std::array{p_local_coords})[0];
158 .template value<double>(vars, pos, t, dt);
159 double pressure = 0.0;
164 auto const k = MaterialPropertyLib::formEigenTensor<GlobalDim>(
166 .value(vars, pos, t, dt));
168 Eigen::Vector3d flux(0.0, 0.0, 0.0);
169 flux.head<GlobalDim>() =
170 -k * shape_matrices.dNdx *
171 Eigen::Map<const NodalVectorType>(local_x.data(), local_x.size());
std::vector< typename ShapeMatricesType::ShapeMatrices, Eigen::aligned_allocator< typename ShapeMatricesType::ShapeMatrices > > computeShapeMatrices(MeshLib::Element const &e, bool const is_axially_symmetric, PointContainer const &points)
References ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_element, ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_process_data, NumLib::computeShapeMatrices(), MaterialPropertyLib::diffusion, MeshLib::Element::getID(), ProcessLib::SteadyStateDiffusion::SteadyStateDiffusionData::media_map, MaterialPropertyLib::phase_pressure, MaterialPropertyLib::reference_temperature, ParameterLib::SpatialPosition::setElementID(), NumLib::shapeFunctionInterpolate(), and MaterialPropertyLib::temperature.
template<typename ShapeFunction , typename IntegrationMethod , int GlobalDim>
Implements ProcessLib::SteadyStateDiffusion::SteadyStateDiffusionLocalAssemblerInterface.
Definition at line 185 of file SteadyStateDiffusionFEM.h.
193 double const dt = std::numeric_limits<double>::quiet_NaN();
195 auto const n_integration_points =
198 int const process_id = 0;
201 assert(!indices.empty());
202 auto const local_x = x[process_id]->get(indices);
203 auto const local_x_vec =
204 MathLib::toVector<Eigen::Matrix<double, ShapeFunction::NPOINTS, 1>>(
205 local_x, ShapeFunction::NPOINTS);
209 Eigen::Matrix<double, GlobalDim, Eigen::Dynamic, Eigen::RowMajor>>(
210 cache, GlobalDim, n_integration_points);
223 .template value<double>(vars, pos, t, dt);
224 double pressure = 0.0;
225 for (
unsigned i = 0; i < n_integration_points; ++i)
227 pos.setIntegrationPoint(i);
230 vars[
static_cast<int>(
233 auto const k = MaterialPropertyLib::formEigenTensor<GlobalDim>(
235 .value(vars, pos, t, dt));
237 cache_mat.col(i).noalias() =
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
std::vector< GlobalIndexType > getIndices(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
References ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_element, ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_integration_method, ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_process_data, ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_shape_matrices, MathLib::createZeroedMatrix(), MaterialPropertyLib::diffusion, MeshLib::Element::getID(), NumLib::getIndices(), ProcessLib::SteadyStateDiffusion::SteadyStateDiffusionData::media_map, MaterialPropertyLib::phase_pressure, MaterialPropertyLib::reference_temperature, ParameterLib::SpatialPosition::setElementID(), NumLib::shapeFunctionInterpolate(), and MaterialPropertyLib::temperature.