template<typename ShapeFunction, int GlobalDim>
class ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >
Definition at line 50 of file SteadyStateDiffusionFEM.h.
|
| LocalAssemblerData (MeshLib::Element const &element, std::size_t const, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, 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.
|
|
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 |
|
virtual void | setInitialConditions (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, 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_x_prev, 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_x_prev, 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_x_prev, int const process_id, 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_prev, 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, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id) |
|
void | postNonLinearSolver (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id) |
|
virtual Eigen::Vector3d | getFlux (MathLib::Point3d const &, double const, std::vector< std::vector< double > > const &) const |
| Fits to staggered scheme.
|
|
template<typename ShapeFunction , 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 131 of file SteadyStateDiffusionFEM.h.
134 {
135
136
137 double const dt = std::numeric_limits<double>::quiet_NaN();
138
139
140
141
142 auto const shape_matrices =
144 GlobalDim>(
146 std::array{p_local_coords})[0];
147
148
151 auto const& medium =
153
156 medium
157 .property(
159 .template value<double>(vars, pos, t, dt);
160 double pressure = 0.0;
163
166 .value(vars, pos, t, dt));
167
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());
172
173 return flux;
174 }
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, GlobalDim >::_element, ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >::_process_data, NumLib::computeShapeMatrices(), MaterialPropertyLib::diffusion, MaterialPropertyLib::formEigenTensor(), MeshLib::Element::getID(), MaterialPropertyLib::MaterialSpatialDistributionMap::getMedium(), MaterialPropertyLib::VariableArray::liquid_phase_pressure, ProcessLib::SteadyStateDiffusion::SteadyStateDiffusionData::media_map, MaterialPropertyLib::reference_temperature, ParameterLib::SpatialPosition::setElementID(), NumLib::detail::shapeFunctionInterpolate(), and MaterialPropertyLib::VariableArray::temperature.
template<typename ShapeFunction , int GlobalDim>
Implements ProcessLib::SteadyStateDiffusion::SteadyStateDiffusionLocalAssemblerInterface.
Definition at line 185 of file SteadyStateDiffusionFEM.h.
190 {
191
192
193 double const dt = std::numeric_limits<double>::quiet_NaN();
194
195 auto const n_integration_points =
197
198 int const process_id = 0;
199 auto const indices =
201 assert(!indices.empty());
202 auto const local_x = x[process_id]->get(indices);
203 auto const local_x_vec =
205 local_x, ShapeFunction::NPOINTS);
206
207 cache.clear();
209 Eigen::Matrix<double, GlobalDim, Eigen::Dynamic, Eigen::RowMajor>>(
210 cache, GlobalDim, n_integration_points);
211
214
215 auto const& medium =
217
220 medium
221 .property(
223 .template value<double>(vars, pos, t, dt);
224 double pressure = 0.0;
225 for (unsigned i = 0; i < n_integration_points; ++i)
226 {
227 pos.setIntegrationPoint(i);
229 pressure);
231
234 .value(vars, pos, t, dt));
235
236 cache_mat.col(i).noalias() =
238 }
239
240 return cache;
241 }
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)
References ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >::_element, ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >::_integration_method, ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >::_process_data, ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >::_shape_matrices, MathLib::createZeroedMatrix(), MaterialPropertyLib::diffusion, MaterialPropertyLib::formEigenTensor(), MeshLib::Element::getID(), NumLib::getIndices(), MaterialPropertyLib::MaterialSpatialDistributionMap::getMedium(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), MaterialPropertyLib::VariableArray::liquid_phase_pressure, ProcessLib::SteadyStateDiffusion::SteadyStateDiffusionData::media_map, MaterialPropertyLib::reference_temperature, ParameterLib::SpatialPosition::setElementID(), NumLib::detail::shapeFunctionInterpolate(), MaterialPropertyLib::VariableArray::temperature, and MathLib::toVector().