Skip to content
Snippets Groups Projects
Commit dd8c5ca1 authored by Tom Fischer's avatar Tom Fischer
Browse files

[PL] Time dependency for the flux calculation.

parent 3f11cb84
No related branches found
No related tags found
No related merge requests found
...@@ -55,13 +55,14 @@ CalculateSurfaceFlux::CalculateSurfaceFlux( ...@@ -55,13 +55,14 @@ CalculateSurfaceFlux::CalculateSurfaceFlux(
void CalculateSurfaceFlux::integrate(GlobalVector const& x, void CalculateSurfaceFlux::integrate(GlobalVector const& x,
MeshLib::PropertyVector<double>& balance, MeshLib::PropertyVector<double>& balance,
double const t,
Process const& bulk_process) Process const& bulk_process)
{ {
DBUG("Integrate CalculateSurfaceFlux."); DBUG("Integrate CalculateSurfaceFlux.");
GlobalExecutor::executeMemberOnDereferenced( GlobalExecutor::executeMemberOnDereferenced(
&CalculateSurfaceFluxLocalAssemblerInterface::integrate, &CalculateSurfaceFluxLocalAssemblerInterface::integrate,
_local_assemblers, x, balance, bulk_process); _local_assemblers, x, balance, t, bulk_process);
} }
} // namespace ProcessLib } // namespace ProcessLib
...@@ -29,10 +29,12 @@ public: ...@@ -29,10 +29,12 @@ public:
/// Executes for each element of the mesh the local intergration procedure. /// Executes for each element of the mesh the local intergration procedure.
/// @param x The global solution the intergration values will be fetched of. /// @param x The global solution the intergration values will be fetched of.
/// @param balance The vector the integration results will be stored in. /// @param balance The vector the integration results will be stored in.
/// @param t The balance will be computed at the time t.
/// @param bulk_process Stores the variable that is used for the /// @param bulk_process Stores the variable that is used for the
/// integration. /// integration.
void integrate(GlobalVector const& x, void integrate(GlobalVector const& x,
MeshLib::PropertyVector<double>& balance, MeshLib::PropertyVector<double>& balance,
double const t,
Process const& bulk_process); Process const& bulk_process);
private: private:
......
...@@ -31,6 +31,7 @@ public: ...@@ -31,6 +31,7 @@ public:
virtual void integrate(std::size_t element_id, virtual void integrate(std::size_t element_id,
GlobalVector const& x, GlobalVector const& x,
MeshLib::PropertyVector<double>& balance, MeshLib::PropertyVector<double>& balance,
double const t,
Process const& bulk_process) = 0; Process const& bulk_process) = 0;
}; };
...@@ -101,11 +102,13 @@ public: ...@@ -101,11 +102,13 @@ public:
/// @param balance PropertyVector the integration result will be stored /// @param balance PropertyVector the integration result will be stored
/// into, where balance has the same number of entries like mesh elements /// into, where balance has the same number of entries like mesh elements
/// exists. /// exists.
/// @param t The integration is performed at the time t.
/// @param bulk_process Reference to the bulk process is used to compute the /// @param bulk_process Reference to the bulk process is used to compute the
/// flux. /// flux.
void integrate(std::size_t element_id, void integrate(std::size_t element_id,
GlobalVector const& x, GlobalVector const& x,
MeshLib::PropertyVector<double>& balance, MeshLib::PropertyVector<double>& balance,
double const t,
Process const& bulk_process) override Process const& bulk_process) override
{ {
auto surface_element_normal = auto surface_element_normal =
...@@ -126,8 +129,8 @@ public: ...@@ -126,8 +129,8 @@ public:
auto const bulk_element_point = MeshLib::getBulkElementPoint( auto const bulk_element_point = MeshLib::getBulkElementPoint(
bulk_process.getMesh(), _bulk_element_id, _bulk_face_id, wp); bulk_process.getMesh(), _bulk_element_id, _bulk_face_id, wp);
auto const bulk_flux = auto const bulk_flux = bulk_process.getFlux(
bulk_process.getFlux(_bulk_element_id, bulk_element_point, x); _bulk_element_id, bulk_element_point, t, x);
for (int component_id(0); for (int component_id(0);
component_id < balance.getNumberOfComponents(); component_id < balance.getNumberOfComponents();
......
...@@ -104,10 +104,9 @@ public: ...@@ -104,10 +104,9 @@ public:
/// Computes the flux in the point \c p_local_coords that is given in local /// Computes the flux in the point \c p_local_coords that is given in local
/// coordinates using the values from \c local_x. /// coordinates using the values from \c local_x.
// TODO add time dependency Eigen::Vector3d getFlux(MathLib::Point3d const& p_local_coords,
Eigen::Vector3d getFlux( double const t,
MathLib::Point3d const& p_local_coords, std::vector<double> const& local_x) const override
std::vector<double> const& local_x) const override
{ {
// eval dNdx and invJ at p // eval dNdx and invJ at p
using FemType = using FemType =
...@@ -127,8 +126,6 @@ public: ...@@ -127,8 +126,6 @@ public:
// fetch hydraulic conductivity // fetch hydraulic conductivity
SpatialPosition pos; SpatialPosition pos;
pos.setElementID(_element.getID()); pos.setElementID(_element.getID());
// TODO remove follwing line if time dependency is implemented
double const t = 0.0;
auto const k = _process_data.hydraulic_conductivity(t, pos)[0]; auto const k = _process_data.hydraulic_conductivity(t, pos)[0];
Eigen::Vector3d flux; Eigen::Vector3d flux;
......
...@@ -48,8 +48,9 @@ public: ...@@ -48,8 +48,9 @@ public:
//! @} //! @}
Eigen::Vector3d getFlux(std::size_t element_id, Eigen::Vector3d getFlux(std::size_t element_id,
MathLib::Point3d const& p, MathLib::Point3d const& p,
GlobalVector const& x) const override double const t,
GlobalVector const& x) const override
{ {
// fetch local_x from primary variable // fetch local_x from primary variable
std::vector<GlobalIndexType> indices_cache; std::vector<GlobalIndexType> indices_cache;
...@@ -57,11 +58,11 @@ public: ...@@ -57,11 +58,11 @@ public:
element_id, *_local_to_global_index_map, indices_cache); element_id, *_local_to_global_index_map, indices_cache);
std::vector<double> local_x(x.get(r_c_indices.rows)); std::vector<double> local_x(x.get(r_c_indices.rows));
return _local_assemblers[element_id]->getFlux(p, local_x); return _local_assemblers[element_id]->getFlux(p, t, local_x);
} }
void postTimestepConcreteProcess(GlobalVector const& x, void postTimestepConcreteProcess(GlobalVector const& x,
const double /*t*/, const double t,
const double /*delta_t*/, const double /*delta_t*/,
int const process_id) override int const process_id) override
{ {
...@@ -89,7 +90,7 @@ public: ...@@ -89,7 +90,7 @@ public:
_balance_mesh->getProperties() _balance_mesh->getProperties()
.template getPropertyVector<double>(_balance_pv_name); .template getPropertyVector<double>(_balance_pv_name);
balance.integrate(x, *balance_pv, *this); balance.integrate(x, *balance_pv, t, *this);
// post: surface_mesh has vectorial element property // post: surface_mesh has vectorial element property
// TODO output, if output classes are ready this has to be // TODO output, if output classes are ready this has to be
......
...@@ -90,9 +90,9 @@ public: ...@@ -90,9 +90,9 @@ public:
/// Computes the flux in the point \c pnt_local_coords that is given in /// Computes the flux in the point \c pnt_local_coords that is given in
/// local coordinates using the values from \c local_x. /// local coordinates using the values from \c local_x.
// TODO add time dependency
Eigen::Vector3d getFlux( Eigen::Vector3d getFlux(
MathLib::Point3d const& pnt_local_coords, MathLib::Point3d const& pnt_local_coords,
double const t,
std::vector<double> const& local_x) const override std::vector<double> const& local_x) const override
{ {
// eval dNdx and invJ at given point // eval dNdx and invJ at given point
...@@ -123,8 +123,6 @@ public: ...@@ -123,8 +123,6 @@ public:
vars[static_cast<int>( vars[static_cast<int>(
MaterialLib::Fluid::PropertyVariableType::p)]); MaterialLib::Fluid::PropertyVariableType::p)]);
// TODO remove following line if time dependency is implemented
double const t = 0.0;
auto const K = auto const K =
this->_material_properties.porous_media_properties this->_material_properties.porous_media_properties
.getIntrinsicPermeability(t, pos) .getIntrinsicPermeability(t, pos)
......
...@@ -61,6 +61,7 @@ public: ...@@ -61,6 +61,7 @@ public:
virtual Eigen::Vector3d getFlux( virtual Eigen::Vector3d getFlux(
MathLib::Point3d const& pnt_local_coords, MathLib::Point3d const& pnt_local_coords,
double const t,
std::vector<double> const& local_x) const = 0; std::vector<double> const& local_x) const = 0;
protected: protected:
......
...@@ -242,6 +242,7 @@ void HTProcess::setCoupledSolutionsOfPreviousTimeStep() ...@@ -242,6 +242,7 @@ void HTProcess::setCoupledSolutionsOfPreviousTimeStep()
Eigen::Vector3d HTProcess::getFlux(std::size_t element_id, Eigen::Vector3d HTProcess::getFlux(std::size_t element_id,
MathLib::Point3d const& p, MathLib::Point3d const& p,
double const t,
GlobalVector const& x) const GlobalVector const& x) const
{ {
// fetch local_x from primary variable // fetch local_x from primary variable
...@@ -250,7 +251,7 @@ Eigen::Vector3d HTProcess::getFlux(std::size_t element_id, ...@@ -250,7 +251,7 @@ Eigen::Vector3d HTProcess::getFlux(std::size_t element_id,
element_id, *_local_to_global_index_map, indices_cache); element_id, *_local_to_global_index_map, indices_cache);
std::vector<double> local_x(x.get(r_c_indices.rows)); std::vector<double> local_x(x.get(r_c_indices.rows));
return _local_assemblers[element_id]->getFlux(p, local_x); return _local_assemblers[element_id]->getFlux(p, t, local_x);
} }
// this is almost a copy of the implemention in the GroundwaterFlow // this is almost a copy of the implemention in the GroundwaterFlow
...@@ -284,7 +285,7 @@ void HTProcess::postTimestepConcreteProcess(GlobalVector const& x, ...@@ -284,7 +285,7 @@ void HTProcess::postTimestepConcreteProcess(GlobalVector const& x,
getProcessVariables(process_id)[0].get().getNumberOfComponents(), getProcessVariables(process_id)[0].get().getNumberOfComponents(),
_integration_order); _integration_order);
balance.integrate(x, *balance_pv, *this); balance.integrate(x, *balance_pv, t, *this);
// post: surface_mesh has scalar element property // post: surface_mesh has scalar element property
// TODO output, if output classes are ready this has to be // TODO output, if output classes are ready this has to be
......
...@@ -70,8 +70,9 @@ public: ...@@ -70,8 +70,9 @@ public:
//! @} //! @}
Eigen::Vector3d getFlux(std::size_t element_id, Eigen::Vector3d getFlux(std::size_t element_id,
MathLib::Point3d const& p, MathLib::Point3d const& p,
GlobalVector const& x) const override; double const t,
GlobalVector const& x) const override;
void setCoupledTermForTheStaggeredSchemeToLocalAssemblers() override; void setCoupledTermForTheStaggeredSchemeToLocalAssemblers() override;
......
...@@ -87,6 +87,7 @@ public: ...@@ -87,6 +87,7 @@ public:
/// coordinates using the values from \c local_x. /// coordinates using the values from \c local_x.
virtual Eigen::Vector3d getFlux( virtual Eigen::Vector3d getFlux(
MathLib::Point3d const& /*p_local_coords*/, MathLib::Point3d const& /*p_local_coords*/,
double const /*t*/,
std::vector<double> const& /*local_x*/) const std::vector<double> const& /*local_x*/) const
{ {
return Eigen::Vector3d{}; return Eigen::Vector3d{};
......
...@@ -138,6 +138,7 @@ public: ...@@ -138,6 +138,7 @@ public:
virtual Eigen::Vector3d getFlux(std::size_t /*element_id*/, virtual Eigen::Vector3d getFlux(std::size_t /*element_id*/,
MathLib::Point3d const& /*p*/, MathLib::Point3d const& /*p*/,
double const /*t*/,
GlobalVector const& /*x*/) const GlobalVector const& /*x*/) const
{ {
return Eigen::Vector3d{}; return Eigen::Vector3d{};
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment