From c42a65ea769b2d2e882daed3ec40b4dfe3891d87 Mon Sep 17 00:00:00 2001 From: Dmitri Naumov <dmitri.naumov@ufz.de> Date: Wed, 18 Sep 2019 20:03:17 +0200 Subject: [PATCH] [PL] Pass dt to assembleWithJacobian. Exemplary on SmallDeformationProcess, storage of dt is no longer required and the assembly can be repeated without dependency on the t/dt values. --- NumLib/ODESolver/ODESystem.h | 3 ++- NumLib/ODESolver/TimeDiscretization.h | 9 +++++++ NumLib/ODESolver/TimeDiscretizedODESystem.cpp | 5 ++-- ProcessLib/AbstractJacobianAssembler.h | 16 +++++++------ ProcessLib/AnalyticalJacobianAssembler.cpp | 4 ++-- ProcessLib/AnalyticalJacobianAssembler.h | 16 +++++++------ .../CentralDifferencesJacobianAssembler.cpp | 2 +- .../CentralDifferencesJacobianAssembler.h | 16 +++++++------ .../CompareJacobiansJacobianAssembler.cpp | 6 ++--- .../CompareJacobiansJacobianAssembler.h | 2 +- ProcessLib/LocalAssemblerInterface.cpp | 3 ++- ProcessLib/LocalAssemblerInterface.h | 2 +- ProcessLib/Process.cpp | 7 +++--- ProcessLib/Process.h | 15 ++++++------ .../SmallDeformation/SmallDeformationFEM.h | 6 ++--- .../SmallDeformationProcess.cpp | 24 +++++-------------- .../SmallDeformationProcess.h | 11 ++++----- .../SmallDeformationProcessData.h | 2 -- ProcessLib/VectorMatrixAssembler.cpp | 8 +++---- ProcessLib/VectorMatrixAssembler.h | 6 ++--- 20 files changed, 83 insertions(+), 80 deletions(-) diff --git a/NumLib/ODESolver/ODESystem.h b/NumLib/ODESolver/ODESystem.h index 233ea92da8f..1f7437ebb34 100644 --- a/NumLib/ODESolver/ODESystem.h +++ b/NumLib/ODESolver/ODESystem.h @@ -125,7 +125,8 @@ public: * \mathtt{Jac} \f$. * \endparblock */ - virtual void assembleWithJacobian(const double t, GlobalVector const& x, + virtual void assembleWithJacobian(const double t, double const dt, + GlobalVector const& x, GlobalVector const& xdot, const double dxdot_dx, const double dx_dx, GlobalMatrix& M, GlobalMatrix& K, diff --git a/NumLib/ODESolver/TimeDiscretization.h b/NumLib/ODESolver/TimeDiscretization.h index e65f20459dc..4033a5f232f 100644 --- a/NumLib/ODESolver/TimeDiscretization.h +++ b/NumLib/ODESolver/TimeDiscretization.h @@ -166,6 +166,10 @@ public: //! assembled. virtual double getCurrentTime() const = 0; + //! Returns \f$ \Delta t_C \f$, i.e., the time at which the equation will be + //! assembled. + virtual double getCurrentTimeIncrement() const = 0; + //! Returns \f$ \hat x \f$, i.e. the discretized approximation of \f$ \dot x //! \f$. void getXdot(GlobalVector const& x_at_new_timestep, GlobalVector& xdot) const @@ -286,6 +290,7 @@ public: } double getCurrentTime() const override { return _t; } + double getCurrentTimeIncrement() const override { return _delta_t; } double getNewXWeight() const override { return 1.0 / _delta_t; } void getWeightedOldX(GlobalVector& y) const override { @@ -349,6 +354,8 @@ public: return _t_old; // forward Euler does assembly at the preceding timestep } + double getCurrentTimeIncrement() const override { return _delta_t; } + GlobalVector const& getCurrentX( const GlobalVector& /*x_at_new_timestep*/) const override { @@ -430,6 +437,7 @@ public: } double getCurrentTime() const override { return _t; } + double getCurrentTimeIncrement() const override { return _delta_t; } double getNewXWeight() const override { return 1.0 / _delta_t; } void getWeightedOldX(GlobalVector& y) const override { @@ -509,6 +517,7 @@ public: } double getCurrentTime() const override { return _t; } + double getCurrentTimeIncrement() const override { return _delta_t; } double getNewXWeight() const override; diff --git a/NumLib/ODESolver/TimeDiscretizedODESystem.cpp b/NumLib/ODESolver/TimeDiscretizedODESystem.cpp index 0b6dc1741e9..6fd0a0b847f 100644 --- a/NumLib/ODESolver/TimeDiscretizedODESystem.cpp +++ b/NumLib/ODESolver/TimeDiscretizedODESystem.cpp @@ -76,6 +76,7 @@ void TimeDiscretizedODESystem< namespace LinAlg = MathLib::LinAlg; auto const t = _time_disc.getCurrentTime(); + auto const dt = _time_disc.getCurrentTimeIncrement(); auto const& x_curr = _time_disc.getCurrentX(x_new_timestep); auto const dxdot_dx = _time_disc.getNewXWeight(); auto const dx_dx = _time_disc.getDxDx(); @@ -89,8 +90,8 @@ void TimeDiscretizedODESystem< _Jac->setZero(); _ode.preAssemble(t, x_curr); - _ode.assembleWithJacobian(t, x_curr, xdot, dxdot_dx, dx_dx, *_M, *_K, *_b, - *_Jac); + _ode.assembleWithJacobian(t, dt, x_curr, xdot, dxdot_dx, dx_dx, *_M, *_K, + *_b, *_Jac); LinAlg::finalizeAssembly(*_M); LinAlg::finalizeAssembly(*_K); diff --git a/ProcessLib/AbstractJacobianAssembler.h b/ProcessLib/AbstractJacobianAssembler.h index 54a31db8452..7eb946c1157 100644 --- a/ProcessLib/AbstractJacobianAssembler.h +++ b/ProcessLib/AbstractJacobianAssembler.h @@ -24,13 +24,15 @@ class AbstractJacobianAssembler public: //! Assembles the Jacobian, the matrices \f$M\f$ and \f$K\f$, and the vector //! \f$b\f$. - virtual void assembleWithJacobian( - LocalAssemblerInterface& local_assembler, double const t, - 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) = 0; + virtual void assembleWithJacobian(LocalAssemblerInterface& local_assembler, + 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) = 0; //! Assembles the Jacobian, the matrices \f$M\f$ and \f$K\f$, and the vector //! \f$b\f$ with coupling. diff --git a/ProcessLib/AnalyticalJacobianAssembler.cpp b/ProcessLib/AnalyticalJacobianAssembler.cpp index 8d9fd3d0495..5bf99eca048 100644 --- a/ProcessLib/AnalyticalJacobianAssembler.cpp +++ b/ProcessLib/AnalyticalJacobianAssembler.cpp @@ -15,13 +15,13 @@ namespace ProcessLib { void AnalyticalJacobianAssembler::assembleWithJacobian( - LocalAssemblerInterface& local_assembler, double const t, + LocalAssemblerInterface& local_assembler, 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) { - local_assembler.assembleWithJacobian(t, local_x, local_xdot, dxdot_dx, + local_assembler.assembleWithJacobian(t, dt, local_x, local_xdot, dxdot_dx, dx_dx, local_M_data, local_K_data, local_b_data, local_Jac_data); } diff --git a/ProcessLib/AnalyticalJacobianAssembler.h b/ProcessLib/AnalyticalJacobianAssembler.h index 3df92c3839d..e71d15dc886 100644 --- a/ProcessLib/AnalyticalJacobianAssembler.h +++ b/ProcessLib/AnalyticalJacobianAssembler.h @@ -31,13 +31,15 @@ public: //! \f$b\f$. //! In this implementation the call is only forwarded to the respective //! method of the given \c local_assembler. - void assembleWithJacobian( - LocalAssemblerInterface& local_assembler, double const t, - 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) override; + void assembleWithJacobian(LocalAssemblerInterface& local_assembler, + 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) override; void assembleWithJacobianForStaggeredScheme( LocalAssemblerInterface& local_assembler, diff --git a/ProcessLib/CentralDifferencesJacobianAssembler.cpp b/ProcessLib/CentralDifferencesJacobianAssembler.cpp index c335fa99933..3a93d3496c2 100644 --- a/ProcessLib/CentralDifferencesJacobianAssembler.cpp +++ b/ProcessLib/CentralDifferencesJacobianAssembler.cpp @@ -26,7 +26,7 @@ CentralDifferencesJacobianAssembler::CentralDifferencesJacobianAssembler( } void CentralDifferencesJacobianAssembler::assembleWithJacobian( - LocalAssemblerInterface& local_assembler, const double t, + LocalAssemblerInterface& local_assembler, const double t, double const dt, const std::vector<double>& local_x_data, const std::vector<double>& local_xdot_data, const double dxdot_dx, const double dx_dx, std::vector<double>& local_M_data, diff --git a/ProcessLib/CentralDifferencesJacobianAssembler.h b/ProcessLib/CentralDifferencesJacobianAssembler.h index 27b68a30684..2d5736fd18b 100644 --- a/ProcessLib/CentralDifferencesJacobianAssembler.h +++ b/ProcessLib/CentralDifferencesJacobianAssembler.h @@ -50,13 +50,15 @@ public: //! //! \attention It is assumed that the local vectors and matrices are ordered //! by component. - void assembleWithJacobian( - LocalAssemblerInterface& local_assembler, double const t, - 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) override; + void assembleWithJacobian(LocalAssemblerInterface& local_assembler, + 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) override; private: std::vector<double> const _absolute_epsilons; diff --git a/ProcessLib/CompareJacobiansJacobianAssembler.cpp b/ProcessLib/CompareJacobiansJacobianAssembler.cpp index 0955881bc42..816f60130c0 100644 --- a/ProcessLib/CompareJacobiansJacobianAssembler.cpp +++ b/ProcessLib/CompareJacobiansJacobianAssembler.cpp @@ -125,7 +125,7 @@ const std::string msg_fatal = namespace ProcessLib { void CompareJacobiansJacobianAssembler::assembleWithJacobian( - LocalAssemblerInterface& local_assembler, double const t, + LocalAssemblerInterface& local_assembler, 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, @@ -147,7 +147,7 @@ void CompareJacobiansJacobianAssembler::assembleWithJacobian( // First assembly -- the one whose results will be added to the global // equation system finally. - _asm1->assembleWithJacobian(local_assembler, t, local_x, local_xdot, + _asm1->assembleWithJacobian(local_assembler, t, dt, local_x, local_xdot, dxdot_dx, dx_dx, local_M_data, local_K_data, local_b_data, local_Jac_data); @@ -161,7 +161,7 @@ void CompareJacobiansJacobianAssembler::assembleWithJacobian( std::vector<double> local_Jac_data2; // Second assembly -- used for checking only. - _asm2->assembleWithJacobian(local_assembler, t, local_x, local_xdot, + _asm2->assembleWithJacobian(local_assembler, t, dt, local_x, local_xdot, dxdot_dx, dx_dx, local_M_data2, local_K_data2, local_b_data2, local_Jac_data2); diff --git a/ProcessLib/CompareJacobiansJacobianAssembler.h b/ProcessLib/CompareJacobiansJacobianAssembler.h index e4170e1287f..3c429372b3f 100644 --- a/ProcessLib/CompareJacobiansJacobianAssembler.h +++ b/ProcessLib/CompareJacobiansJacobianAssembler.h @@ -52,7 +52,7 @@ public: } void assembleWithJacobian(LocalAssemblerInterface& local_assembler, - double const t, + 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, diff --git a/ProcessLib/LocalAssemblerInterface.cpp b/ProcessLib/LocalAssemblerInterface.cpp index a7f7d059506..4c35edd8728 100644 --- a/ProcessLib/LocalAssemblerInterface.cpp +++ b/ProcessLib/LocalAssemblerInterface.cpp @@ -39,7 +39,8 @@ void LocalAssemblerInterface::assembleForStaggeredScheme( } void LocalAssemblerInterface::assembleWithJacobian( - double const /*t*/, std::vector<double> const& /*local_x*/, + 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*/, diff --git a/ProcessLib/LocalAssemblerInterface.h b/ProcessLib/LocalAssemblerInterface.h index b19832a86dd..dcac7e9e4aa 100644 --- a/ProcessLib/LocalAssemblerInterface.h +++ b/ProcessLib/LocalAssemblerInterface.h @@ -58,7 +58,7 @@ public: std::vector<double>& local_b_data, LocalCoupledSolutions const& coupled_solutions); - virtual void assembleWithJacobian(double const t, + 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, diff --git a/ProcessLib/Process.cpp b/ProcessLib/Process.cpp index f131c7d716b..f4c04257a40 100644 --- a/ProcessLib/Process.cpp +++ b/ProcessLib/Process.cpp @@ -210,7 +210,8 @@ void Process::assemble(const double t, GlobalVector const& x, GlobalMatrix& M, _source_term_collections[pcs_id].integrate(t, x, b, nullptr); } -void Process::assembleWithJacobian(const double t, GlobalVector const& x, +void Process::assembleWithJacobian(const double t, double const dt, + GlobalVector const& x, GlobalVector const& xdot, const double dxdot_dx, const double dx_dx, GlobalMatrix& M, GlobalMatrix& K, @@ -219,8 +220,8 @@ void Process::assembleWithJacobian(const double t, GlobalVector const& x, MathLib::LinAlg::setLocalAccessibleVector(x); MathLib::LinAlg::setLocalAccessibleVector(xdot); - assembleWithJacobianConcreteProcess(t, x, xdot, dxdot_dx, dx_dx, M, K, b, - Jac); + assembleWithJacobianConcreteProcess(t, dt, x, xdot, dxdot_dx, dx_dx, M, K, + b, Jac); // TODO: apply BCs to Jacobian. const auto pcs_id = diff --git a/ProcessLib/Process.h b/ProcessLib/Process.h index 3de497f10bf..22ce9b9a881 100644 --- a/ProcessLib/Process.h +++ b/ProcessLib/Process.h @@ -101,10 +101,10 @@ public: void assemble(const double t, GlobalVector const& x, GlobalMatrix& M, GlobalMatrix& K, GlobalVector& b) final; - void assembleWithJacobian(const double t, GlobalVector const& x, - GlobalVector const& xdot, const double dxdot_dx, - const double dx_dx, GlobalMatrix& M, - GlobalMatrix& K, GlobalVector& b, + void assembleWithJacobian(const double t, double const dt, + GlobalVector const& x, GlobalVector const& xdot, + const double dxdot_dx, const double dx_dx, + GlobalMatrix& M, GlobalMatrix& K, GlobalVector& b, GlobalMatrix& Jac) final; std::vector<NumLib::IndexValueVector<GlobalIndexType>> const* @@ -194,9 +194,10 @@ private: GlobalVector& b) = 0; virtual void assembleWithJacobianConcreteProcess( - const double t, GlobalVector const& x, GlobalVector const& xdot, - const double dxdot_dx, const double dx_dx, GlobalMatrix& M, - GlobalMatrix& K, GlobalVector& b, GlobalMatrix& Jac) = 0; + const double t, double const dt, GlobalVector const& x, + GlobalVector const& xdot, const double dxdot_dx, const double dx_dx, + GlobalMatrix& M, GlobalMatrix& K, GlobalVector& b, + GlobalMatrix& Jac) = 0; virtual void preTimestepConcreteProcess(GlobalVector const& /*x*/, const double /*t*/, diff --git a/ProcessLib/SmallDeformation/SmallDeformationFEM.h b/ProcessLib/SmallDeformation/SmallDeformationFEM.h index 82bad1bd52f..0b4988ef9df 100644 --- a/ProcessLib/SmallDeformation/SmallDeformationFEM.h +++ b/ProcessLib/SmallDeformation/SmallDeformationFEM.h @@ -230,7 +230,7 @@ public: "implemented."); } - void assembleWithJacobian(double const t, + 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*/, @@ -295,8 +295,8 @@ public: local_x.data(), ShapeFunction::NPOINTS * DisplacementDim); auto&& solution = _ip_data[ip].solid_material.integrateStress( - t, x_position, _process_data.dt, eps_prev, eps, sigma_prev, - *state, _process_data.reference_temperature); + t, x_position, dt, eps_prev, eps, sigma_prev, *state, + _process_data.reference_temperature); if (!solution) { diff --git a/ProcessLib/SmallDeformation/SmallDeformationProcess.cpp b/ProcessLib/SmallDeformation/SmallDeformationProcess.cpp index f884f4c9698..beca0fe087c 100644 --- a/ProcessLib/SmallDeformation/SmallDeformationProcess.cpp +++ b/ProcessLib/SmallDeformation/SmallDeformationProcess.cpp @@ -260,12 +260,10 @@ void SmallDeformationProcess<DisplacementDim>::assembleConcreteProcess( template <int DisplacementDim> void SmallDeformationProcess<DisplacementDim>:: - assembleWithJacobianConcreteProcess(const double t, GlobalVector const& x, - GlobalVector const& xdot, - const double dxdot_dx, - const double dx_dx, GlobalMatrix& M, - GlobalMatrix& K, GlobalVector& b, - GlobalMatrix& Jac) + assembleWithJacobianConcreteProcess( + const double t, double const dt, GlobalVector const& x, + GlobalVector const& xdot, const double dxdot_dx, const double dx_dx, + GlobalMatrix& M, GlobalMatrix& K, GlobalVector& b, GlobalMatrix& Jac) { DBUG("AssembleWithJacobian SmallDeformationProcess."); @@ -277,23 +275,13 @@ void SmallDeformationProcess<DisplacementDim>:: // Call global assembler for each local assembly item. GlobalExecutor::executeSelectedMemberDereferenced( _global_assembler, &VectorMatrixAssembler::assembleWithJacobian, - _local_assemblers, pv.getActiveElementIDs(), dof_table, t, x, - xdot, dxdot_dx, dx_dx, M, K, b, Jac, _coupled_solutions); + _local_assemblers, pv.getActiveElementIDs(), dof_table, t, dt, x, xdot, + dxdot_dx, dx_dx, M, K, b, Jac, _coupled_solutions); transformVariableFromGlobalVector(b, 0, *_local_to_global_index_map, *_nodal_forces, std::negate<double>()); } -template <int DisplacementDim> -void SmallDeformationProcess<DisplacementDim>::preTimestepConcreteProcess( - GlobalVector const& /*x*/, double const /*t*/, double const dt, - const int /*process_id*/) -{ - DBUG("PreTimestep SmallDeformationProcess."); - - _process_data.dt = dt; -} - template <int DisplacementDim> void SmallDeformationProcess<DisplacementDim>::postTimestepConcreteProcess( GlobalVector const& x, const double t, const double dt, diff --git a/ProcessLib/SmallDeformation/SmallDeformationProcess.h b/ProcessLib/SmallDeformation/SmallDeformationProcess.h index 9529ec70a94..8cb3a8ef46d 100644 --- a/ProcessLib/SmallDeformation/SmallDeformationProcess.h +++ b/ProcessLib/SmallDeformation/SmallDeformationProcess.h @@ -55,13 +55,10 @@ private: GlobalVector& b) override; void assembleWithJacobianConcreteProcess( - const double t, GlobalVector const& x, GlobalVector const& xdot, - const double dxdot_dx, const double dx_dx, GlobalMatrix& M, - GlobalMatrix& K, GlobalVector& b, GlobalMatrix& Jac) override; - - void preTimestepConcreteProcess( - GlobalVector const& x, double const t, double const dt, - const int process_id) override; + const double t, double const dt, GlobalVector const& x, + GlobalVector const& xdot, const double dxdot_dx, const double dx_dx, + GlobalMatrix& M, GlobalMatrix& K, GlobalVector& b, + GlobalMatrix& Jac) override; void postTimestepConcreteProcess(GlobalVector const& x, const double t, const double delta_t, diff --git a/ProcessLib/SmallDeformation/SmallDeformationProcessData.h b/ProcessLib/SmallDeformation/SmallDeformationProcessData.h index f0dcd2cc555..8b5a21e98db 100644 --- a/ProcessLib/SmallDeformation/SmallDeformationProcessData.h +++ b/ProcessLib/SmallDeformation/SmallDeformationProcessData.h @@ -53,8 +53,6 @@ struct SmallDeformationProcessData double const reference_temperature = std::numeric_limits<double>::quiet_NaN(); - - double dt = std::numeric_limits<double>::quiet_NaN(); }; } // namespace SmallDeformation diff --git a/ProcessLib/VectorMatrixAssembler.cpp b/ProcessLib/VectorMatrixAssembler.cpp index efdda6826e9..45abf2db3bb 100644 --- a/ProcessLib/VectorMatrixAssembler.cpp +++ b/ProcessLib/VectorMatrixAssembler.cpp @@ -108,9 +108,9 @@ void VectorMatrixAssembler::assembleWithJacobian( std::size_t const mesh_item_id, LocalAssemblerInterface& local_assembler, std::vector<std::reference_wrapper<NumLib::LocalToGlobalIndexMap>> const& dof_tables, - const double t, GlobalVector const& x, GlobalVector const& xdot, - const double dxdot_dx, const double dx_dx, GlobalMatrix& M, GlobalMatrix& K, - GlobalVector& b, GlobalMatrix& Jac, + const double t, double const dt, GlobalVector const& x, + GlobalVector const& xdot, const double dxdot_dx, const double dx_dx, + GlobalMatrix& M, GlobalMatrix& K, GlobalVector& b, GlobalMatrix& Jac, CoupledSolutionsForStaggeredScheme const* const cpl_xs) { std::vector<std::vector<GlobalIndexType>> indices_of_processes; @@ -135,7 +135,7 @@ void VectorMatrixAssembler::assembleWithJacobian( { auto const local_x = x.get(indices); _jacobian_assembler->assembleWithJacobian( - local_assembler, t, local_x, local_xdot, dxdot_dx, dx_dx, + local_assembler, t, dt, local_x, local_xdot, dxdot_dx, dx_dx, _local_M_data, _local_K_data, _local_b_data, _local_Jac_data); } else diff --git a/ProcessLib/VectorMatrixAssembler.h b/ProcessLib/VectorMatrixAssembler.h index 62da2ad0ec1..677e91d846c 100644 --- a/ProcessLib/VectorMatrixAssembler.h +++ b/ProcessLib/VectorMatrixAssembler.h @@ -59,9 +59,9 @@ public: std::vector< std::reference_wrapper<NumLib::LocalToGlobalIndexMap>> const& dof_tables, - const double t, GlobalVector const& x, GlobalVector const& xdot, - const double dxdot_dx, const double dx_dx, GlobalMatrix& M, - GlobalMatrix& K, GlobalVector& b, GlobalMatrix& Jac, + const double t, double const dt, GlobalVector const& x, + GlobalVector const& xdot, const double dxdot_dx, const double dx_dx, + GlobalMatrix& M, GlobalMatrix& K, GlobalVector& b, GlobalMatrix& Jac, CoupledSolutionsForStaggeredScheme const* const cpl_xs); private: -- GitLab