diff --git a/ProcessLib/CalculateSurfaceFlux/CalculateSurfaceFlux.cpp b/ProcessLib/CalculateSurfaceFlux/CalculateSurfaceFlux.cpp index 2116437a3edbf9d2b34feaff60ee38314a0a3371..8d06ccea54f6ef87742ad7692fc68ee6040a8948 100644 --- a/ProcessLib/CalculateSurfaceFlux/CalculateSurfaceFlux.cpp +++ b/ProcessLib/CalculateSurfaceFlux/CalculateSurfaceFlux.cpp @@ -53,16 +53,20 @@ CalculateSurfaceFlux::CalculateSurfaceFlux( *bulk_element_ids, *bulk_face_ids); } -void CalculateSurfaceFlux::integrate(GlobalVector const& x, - MeshLib::PropertyVector<double>& balance, - double const t, - Process const& bulk_process) +void CalculateSurfaceFlux::integrate( + GlobalVector const& x, + MeshLib::PropertyVector<double>& balance, + double const t, + MeshLib::Mesh const& bulk_mesh, + std::function<Eigen::Vector3d(std::size_t const, MathLib::Point3d const&, + double const, GlobalVector const&)> const& + getFlux) { DBUG("Integrate CalculateSurfaceFlux."); GlobalExecutor::executeMemberOnDereferenced( &CalculateSurfaceFluxLocalAssemblerInterface::integrate, - _local_assemblers, x, balance, t, bulk_process); + _local_assemblers, x, balance, t, bulk_mesh, getFlux); } } // namespace ProcessLib diff --git a/ProcessLib/CalculateSurfaceFlux/CalculateSurfaceFlux.h b/ProcessLib/CalculateSurfaceFlux/CalculateSurfaceFlux.h index 983b33816bde2a8e7d0a4ecf7183ff76fb26c1c9..732418bdd1daecf85235a2d79a5a8a4ab923b872 100644 --- a/ProcessLib/CalculateSurfaceFlux/CalculateSurfaceFlux.h +++ b/ProcessLib/CalculateSurfaceFlux/CalculateSurfaceFlux.h @@ -30,12 +30,18 @@ public: /// @param x The global solution the intergration values will be fetched of. /// @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 - /// integration. + /// @param bulk_mesh Stores a reference to the bulk mesh that is needed to + /// fetch the information for the integration over the surface mesh. + /// @param getFlux function that calculates the flux in the integration + /// points of the face elements of the bulk element that belongs to the + /// surface. void integrate(GlobalVector const& x, MeshLib::PropertyVector<double>& balance, double const t, - Process const& bulk_process); + MeshLib::Mesh const& bulk_mesh, + std::function<Eigen::Vector3d( + std::size_t const, MathLib::Point3d const&, double const, + GlobalVector const&)> const& getFlux); private: // the local assemblers for each element of the mesh diff --git a/ProcessLib/CalculateSurfaceFlux/CalculateSurfaceFluxLocalAssembler.h b/ProcessLib/CalculateSurfaceFlux/CalculateSurfaceFluxLocalAssembler.h index 6634fff48800e5692243ba94e46d602501fb1543..398a346369867076ea6da3e384bfac278c12c42a 100644 --- a/ProcessLib/CalculateSurfaceFlux/CalculateSurfaceFluxLocalAssembler.h +++ b/ProcessLib/CalculateSurfaceFlux/CalculateSurfaceFluxLocalAssembler.h @@ -28,11 +28,14 @@ class CalculateSurfaceFluxLocalAssemblerInterface public: virtual ~CalculateSurfaceFluxLocalAssemblerInterface() = default; - virtual void integrate(std::size_t element_id, + virtual void integrate(std::size_t const element_id, GlobalVector const& x, MeshLib::PropertyVector<double>& balance, double const t, - Process const& bulk_process) = 0; + MeshLib::Mesh const& bulk_mesh, + std::function<Eigen::Vector3d( + std::size_t const, MathLib::Point3d const&, + double const, GlobalVector const&)>) = 0; }; template <typename ShapeFunction, typename IntegrationMethod, @@ -103,13 +106,20 @@ public: /// into, where balance has the same number of entries like mesh elements /// exists. /// @param t The integration is performed at the time t. - /// @param bulk_process Reference to the bulk process is used to compute the - /// flux. - void integrate(std::size_t element_id, + /// @param bulk_mesh Stores a reference to the bulk mesh that is needed to + /// fetch the information for the integration over the surface mesh. + /// @param getFlux function that calculates the flux in the integration + /// points of the face elements of the bulk element that belongs to the + /// surface. + void integrate(std::size_t const element_id, GlobalVector const& x, MeshLib::PropertyVector<double>& balance, double const t, - Process const& bulk_process) override + MeshLib::Mesh const& bulk_mesh, + std::function<Eigen::Vector3d( + std::size_t const, MathLib::Point3d const&, double const, + GlobalVector const&)> + getFlux) override { auto surface_element_normal = MeshLib::FaceRule::getSurfaceNormal(&_surface_element); @@ -128,9 +138,9 @@ public: auto const& wp = _integration_method.getWeightedPoint(ip); auto const bulk_element_point = MeshLib::getBulkElementPoint( - bulk_process.getMesh(), _bulk_element_id, _bulk_face_id, wp); - auto const bulk_flux = bulk_process.getFlux( - _bulk_element_id, bulk_element_point, t, x); + bulk_mesh, _bulk_element_id, _bulk_face_id, wp); + auto const bulk_flux = + getFlux(_bulk_element_id, bulk_element_point, t, x); for (int component_id(0); component_id < balance.getNumberOfComponents(); diff --git a/ProcessLib/GroundwaterFlow/GroundwaterFlowProcess.h b/ProcessLib/GroundwaterFlow/GroundwaterFlowProcess.h index c11bb4cdba7eead17d7ca362b63511ddbff49017..199d3248ac3941085d4dced8630a3d440242597f 100644 --- a/ProcessLib/GroundwaterFlow/GroundwaterFlowProcess.h +++ b/ProcessLib/GroundwaterFlow/GroundwaterFlowProcess.h @@ -88,8 +88,13 @@ public: _balance_mesh->getProperties() .template getPropertyVector<double>(_balance_pv_name); - balance.integrate(x, *balance_pv, t, *this); - // post: surface_mesh has vectorial element property + balance.integrate(x, *balance_pv, t, _mesh, + [this](std::size_t const element_id, + MathLib::Point3d const& pnt, + double const t, GlobalVector const& x) { + return getFlux(element_id, pnt, t, x); + }); + // post: surface_mesh has scalar element property // TODO output, if output classes are ready this has to be // changed diff --git a/ProcessLib/HT/HTProcess.cpp b/ProcessLib/HT/HTProcess.cpp index 765328721d65ca7257c6d59a7300c81be5c4826a..97787d9410dc54237059726edbe4abd38e241ff6 100644 --- a/ProcessLib/HT/HTProcess.cpp +++ b/ProcessLib/HT/HTProcess.cpp @@ -288,7 +288,12 @@ void HTProcess::postTimestepConcreteProcess(GlobalVector const& x, getProcessVariables(process_id)[0].get().getNumberOfComponents(), _integration_order); - balance.integrate(x, *balance_pv, t, *this); + balance.integrate( + x, *balance_pv, t, _mesh, + [this](std::size_t const element_id, MathLib::Point3d const& pnt, + double const t, GlobalVector const& x) { + return getFlux(element_id, pnt, t, x); + }); // post: surface_mesh has scalar element property // TODO output, if output classes are ready this has to be