Skip to content
Snippets Groups Projects
Commit 9b47442c authored by Christoph Lehmann's avatar Christoph Lehmann
Browse files

[PL] made Process Jacobian-assembly-ready

parent 11752654
No related branches found
No related tags found
No related merge requests found
...@@ -8,43 +8,21 @@ ...@@ -8,43 +8,21 @@
*/ */
#include "LocalAssemblerInterface.h" #include "LocalAssemblerInterface.h"
#include <cassert>
#include "NumLib/DOF/DOFTableUtil.h" #include "NumLib/DOF/DOFTableUtil.h"
namespace ProcessLib namespace ProcessLib
{ {
void LocalAssemblerInterface::assemble( void LocalAssemblerInterface::assembleWithJacobian(
const std::size_t mesh_item_id, double const /*t*/, std::vector<double> const& /*local_x*/,
const NumLib::LocalToGlobalIndexMap& dof_table, const double t, std::vector<double> const& /*local_xdot*/, const double /*dxdot_dx*/,
const GlobalVector& x, GlobalMatrix& M, GlobalMatrix& K, GlobalVector& b) const double /*dx_dx*/, std::vector<double>& /*local_M_data*/,
{ std::vector<double>& /*local_K_data*/,
auto const indices = NumLib::getIndices(mesh_item_id, dof_table); std::vector<double>& /*local_b_data*/,
auto const local_x = x.get(indices); std::vector<double>& /*local_Jac_data*/)
auto const r_c_indices =
NumLib::LocalToGlobalIndexMap::RowColumnIndices(indices, indices);
assembleConcrete(t, local_x, r_c_indices, M, K, b);
}
void LocalAssemblerInterface::assembleJacobian(
const std::size_t mesh_item_id,
const NumLib::LocalToGlobalIndexMap& dof_table, const double t,
const GlobalVector& x, GlobalMatrix& Jac)
{
auto const indices = NumLib::getIndices(mesh_item_id, dof_table);
auto const local_x = x.get(indices);
auto const r_c_indices =
NumLib::LocalToGlobalIndexMap::RowColumnIndices(indices, indices);
assembleJacobianConcrete(t, local_x, r_c_indices, Jac);
}
void LocalAssemblerInterface::assembleJacobianConcrete(
double const /*t*/, std::vector<double> const& /*local_x*/,
NumLib::LocalToGlobalIndexMap::RowColumnIndices const& /*indices*/,
GlobalMatrix& /*Jac*/)
{ {
OGS_FATAL( OGS_FATAL(
"The assembleJacobian() function is not implemented in the local " "The assembleWithJacobian() function is not implemented in the local "
"assembler."); "assembler.");
} }
......
...@@ -10,9 +10,13 @@ ...@@ -10,9 +10,13 @@
#ifndef PROCESSLIB_LOCALASSEMBLERINTERFACE_H #ifndef PROCESSLIB_LOCALASSEMBLERINTERFACE_H
#define PROCESSLIB_LOCALASSEMBLERINTERFACE_H #define PROCESSLIB_LOCALASSEMBLERINTERFACE_H
#include "NumLib/DOF/LocalToGlobalIndexMap.h"
#include "NumLib/NumericsConfig.h" #include "NumLib/NumericsConfig.h"
namespace NumLib
{
class LocalToGlobalIndexMap;
} // NumLib
namespace ProcessLib namespace ProcessLib
{ {
...@@ -26,15 +30,19 @@ class LocalAssemblerInterface ...@@ -26,15 +30,19 @@ class LocalAssemblerInterface
public: public:
virtual ~LocalAssemblerInterface() = default; virtual ~LocalAssemblerInterface() = default;
void assemble(std::size_t const mesh_item_id, virtual void assemble(
NumLib::LocalToGlobalIndexMap const& dof_table, double const t, std::vector<double> const& local_x,
double const t, GlobalVector const& x, GlobalMatrix& M, std::vector<double>& local_M_data, std::vector<double>& local_K_data,
GlobalMatrix& K, GlobalVector& b); std::vector<double>& local_b_data) = 0;
void assembleJacobian(std::size_t const mesh_item_id, virtual void assembleWithJacobian(double const t,
NumLib::LocalToGlobalIndexMap const& dof_table, std::vector<double> const& local_x,
double const t, GlobalVector const& x, std::vector<double> const& local_xdot,
GlobalMatrix& Jac); 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 preTimestep(std::size_t const mesh_item_id, virtual void preTimestep(std::size_t const mesh_item_id,
NumLib::LocalToGlobalIndexMap const& dof_table, NumLib::LocalToGlobalIndexMap const& dof_table,
...@@ -45,17 +53,7 @@ public: ...@@ -45,17 +53,7 @@ public:
NumLib::LocalToGlobalIndexMap const& dof_table, NumLib::LocalToGlobalIndexMap const& dof_table,
GlobalVector const& x); GlobalVector const& x);
protected: private:
virtual void assembleConcrete(
double const t, std::vector<double> const& local_x,
NumLib::LocalToGlobalIndexMap::RowColumnIndices const& indices,
GlobalMatrix& M, GlobalMatrix& K, GlobalVector& b) = 0;
virtual void assembleJacobianConcrete(
double const t, std::vector<double> const& local_x,
NumLib::LocalToGlobalIndexMap::RowColumnIndices const& indices,
GlobalMatrix& Jac);
virtual void preTimestepConcrete(std::vector<double> const& /*local_x*/, virtual void preTimestepConcrete(std::vector<double> const& /*local_x*/,
double const /*t*/, double const /*dt*/) double const /*t*/, double const /*dt*/)
{ {
......
...@@ -20,6 +20,7 @@ namespace ProcessLib ...@@ -20,6 +20,7 @@ namespace ProcessLib
{ {
Process::Process( Process::Process(
MeshLib::Mesh& mesh, MeshLib::Mesh& mesh,
std::unique_ptr<ProcessLib::AbstractJacobianAssembler>&& jacobian_assembler,
std::vector<std::unique_ptr<ParameterBase>> const& parameters, std::vector<std::unique_ptr<ParameterBase>> const& parameters,
std::vector<std::reference_wrapper<ProcessVariable>>&& process_variables, std::vector<std::reference_wrapper<ProcessVariable>>&& process_variables,
SecondaryVariableCollection&& secondary_variables, SecondaryVariableCollection&& secondary_variables,
...@@ -27,6 +28,7 @@ Process::Process( ...@@ -27,6 +28,7 @@ Process::Process(
: _mesh(mesh), : _mesh(mesh),
_secondary_variables(std::move(secondary_variables)), _secondary_variables(std::move(secondary_variables)),
_named_function_caller(std::move(named_function_caller)), _named_function_caller(std::move(named_function_caller)),
_global_assembler(std::move(jacobian_assembler)),
_process_variables(std::move(process_variables)), _process_variables(std::move(process_variables)),
_boundary_conditions(parameters) _boundary_conditions(parameters)
{ {
...@@ -115,43 +117,16 @@ void Process::assemble(const double t, GlobalVector const& x, GlobalMatrix& M, ...@@ -115,43 +117,16 @@ void Process::assemble(const double t, GlobalVector const& x, GlobalMatrix& M,
_boundary_conditions.applyNaturalBC(t, x, K, b); _boundary_conditions.applyNaturalBC(t, x, K, b);
} }
void Process::assembleJacobian(const double t, GlobalVector const& x, void Process::assembleWithJacobian(const double t, GlobalVector const& x,
GlobalVector const& xdot, const double dxdot_dx, GlobalVector const& xdot,
GlobalMatrix const& M, const double dx_dx, const double dxdot_dx, const double dx_dx,
GlobalMatrix const& K, GlobalMatrix& Jac) GlobalMatrix& M, GlobalMatrix& K,
GlobalVector& b, GlobalMatrix& Jac)
{ {
assembleJacobianConcreteProcess(t, x, xdot, dxdot_dx, M, dx_dx, K, Jac); assembleWithJacobianConcreteProcess(t, x, xdot, dxdot_dx, dx_dx, M, K, b, Jac);
// TODO In this method one could check if the user wants to use an
// analytical or a numerical Jacobian. Then the right
// assembleJacobianConcreteProcess() method will be chosen.
// Additionally in the default implementation of said method one
// could provide a fallback to a numerical Jacobian. However, that
// would be in a sense implicit behaviour and it might be better to
// just abort, as is currently the case.
// In order to implement the Jacobian assembly entirely, in addition
// to the operator() in VectorMatrixAssembler there has to be a method
// that dispatches the Jacobian assembly.
// Similarly, then the NeumannBC specialization of VectorMatrixAssembler
// probably can be merged into the main class s.t. one has only one
// type of VectorMatrixAssembler (for each equation type) with the
// three methods assemble(), assembleJacobian() and assembleNeumannBC().
// That list can be extended, e.g. by methods for the assembly of
// source terms.
// UPDATE: Probably it is better to keep a separate NeumannBC version of the
// VectoMatrixAssembler since that will work for all kinds of processes.
}
void Process::assembleJacobianConcreteProcess( // TODO apply BCs to Jacobian.
const double /*t*/, GlobalVector const& /*x*/, GlobalVector const& /*xdot*/, _boundary_conditions.applyNaturalBC(t, x, K, b);
const double /*dxdot_dx*/, GlobalMatrix const& /*M*/,
const double /*dx_dx*/, GlobalMatrix const& /*K*/, GlobalMatrix& /*Jac*/)
{
OGS_FATAL(
"The concrete implementation of this Process did not override the"
" assembleJacobianConcreteProcess() method."
" Hence, no analytical Jacobian is provided for this process"
" and the Newton-Raphson method cannot be used to solve it.");
} }
void Process::constructDofTable() void Process::constructDofTable()
......
...@@ -40,6 +40,8 @@ public: ...@@ -40,6 +40,8 @@ public:
Process( Process(
MeshLib::Mesh& mesh, MeshLib::Mesh& mesh,
std::unique_ptr<ProcessLib::AbstractJacobianAssembler>&&
jacobian_assembler,
std::vector<std::unique_ptr<ParameterBase>> const& parameters, std::vector<std::unique_ptr<ParameterBase>> const& parameters,
std::vector<std::reference_wrapper<ProcessVariable>>&& std::vector<std::reference_wrapper<ProcessVariable>>&&
process_variables, process_variables,
...@@ -68,14 +70,14 @@ public: ...@@ -68,14 +70,14 @@ public:
void assemble(const double t, GlobalVector const& x, GlobalMatrix& M, void assemble(const double t, GlobalVector const& x, GlobalMatrix& M,
GlobalMatrix& K, GlobalVector& b) override final; GlobalMatrix& K, GlobalVector& b) override final;
void assembleJacobian(const double t, GlobalVector const& x, void assembleWithJacobian(const double t, GlobalVector const& x,
GlobalVector const& xdot, const double dxdot_dx, GlobalVector const& xdot, const double dxdot_dx,
GlobalMatrix const& M, const double dx_dx, const double dx_dx, GlobalMatrix& M,
GlobalMatrix const& K, GlobalMatrix& K, GlobalVector& b,
GlobalMatrix& Jac) override final; GlobalMatrix& Jac) override final;
std::vector<NumLib::IndexValueVector<GlobalIndexType>> const* std::vector<NumLib::IndexValueVector<GlobalIndexType>> const*
getKnownSolutions( getKnownSolutions(
double const t) const override final double const t) const override final
{ {
return _boundary_conditions.getKnownSolutions(t); return _boundary_conditions.getKnownSolutions(t);
...@@ -121,6 +123,12 @@ private: ...@@ -121,6 +123,12 @@ private:
GlobalMatrix& M, GlobalMatrix& K, GlobalMatrix& M, GlobalMatrix& K,
GlobalVector& b) = 0; 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;
virtual void preTimestepConcreteProcess(GlobalVector const& /*x*/, virtual void preTimestepConcreteProcess(GlobalVector const& /*x*/,
const double /*t*/, const double /*t*/,
const double /*delta_t*/) const double /*delta_t*/)
...@@ -169,6 +177,8 @@ protected: ...@@ -169,6 +177,8 @@ protected:
_cached_secondary_variables; _cached_secondary_variables;
SecondaryVariableContext _secondary_variable_context; SecondaryVariableContext _secondary_variable_context;
VectorMatrixAssembler _global_assembler;
private: private:
unsigned const _integration_order = 2; unsigned const _integration_order = 2;
GlobalSparsityPattern _sparsity_pattern; GlobalSparsityPattern _sparsity_pattern;
......
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