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

[PL] TES w/o VetMatAsm

parent 2fbdd7c5
No related branches found
No related tags found
No related merge requests found
...@@ -132,13 +132,16 @@ TESLocalAssembler< ...@@ -132,13 +132,16 @@ TESLocalAssembler<
template <typename ShapeFunction_, typename IntegrationMethod_, template <typename ShapeFunction_, typename IntegrationMethod_,
unsigned GlobalDim> unsigned GlobalDim>
void TESLocalAssembler<ShapeFunction_, IntegrationMethod_, void TESLocalAssembler<ShapeFunction_, IntegrationMethod_, GlobalDim>::assemble(
GlobalDim>::assemble(const double /*t*/, std::size_t const id, NumLib::LocalToGlobalIndexMap const& dof_table,
std::vector<double> const& local_x) double const /*t*/, GlobalVector const& x, GlobalMatrix& M, GlobalMatrix& K,
GlobalVector& b)
{ {
_local_M.setZero(); _local_M.setZero();
_local_K.setZero(); _local_K.setZero();
_local_b.setZero(); _local_b.setZero();
auto const indices = NumLib::detail::getIndices(id, dof_table);
auto const local_x = NumLib::detail::getLocalNodalDOFs(x, indices);
IntegrationMethod_ integration_method(_integration_order); IntegrationMethod_ integration_method(_integration_order);
unsigned const n_integration_points = integration_method.getNumberOfPoints(); unsigned const n_integration_points = integration_method.getNumberOfPoints();
...@@ -182,18 +185,12 @@ void TESLocalAssembler<ShapeFunction_, IntegrationMethod_, ...@@ -182,18 +185,12 @@ void TESLocalAssembler<ShapeFunction_, IntegrationMethod_,
ogs5OutVec(_local_b); ogs5OutVec(_local_b);
std::printf("\n"); std::printf("\n");
} }
}
template <typename ShapeFunction_, typename IntegrationMethod_, auto const r_c_indices =
unsigned GlobalDim> NumLib::LocalToGlobalIndexMap::RowColumnIndices(indices, indices);
void TESLocalAssembler<ShapeFunction_, IntegrationMethod_, GlobalDim>:: M.add(r_c_indices, _local_M);
addToGlobal( K.add(r_c_indices, _local_K);
NumLib::LocalToGlobalIndexMap::RowColumnIndices const& indices, b.add(indices, _local_b);
GlobalMatrix& M, GlobalMatrix& K, GlobalVector& b) const
{
M.add(indices, _local_M);
K.add(indices, _local_K);
b.add(indices.rows, _local_b);
} }
template <typename ShapeFunction_, typename IntegrationMethod_, template <typename ShapeFunction_, typename IntegrationMethod_,
......
...@@ -23,18 +23,12 @@ namespace ProcessLib ...@@ -23,18 +23,12 @@ namespace ProcessLib
namespace TES namespace TES
{ {
class TESLocalAssemblerInterface class TESLocalAssemblerInterface
: public NumLib::Extrapolatable<TESIntPtVariables> : public ProcessLib::LocalAssemblerInterface,
public NumLib::Extrapolatable<TESIntPtVariables>
{ {
public: public:
virtual ~TESLocalAssemblerInterface() = default; virtual ~TESLocalAssemblerInterface() = default;
virtual void assemble(double const t,
std::vector<double> const& local_x) = 0;
virtual void addToGlobal(
NumLib::LocalToGlobalIndexMap::RowColumnIndices const&,
GlobalMatrix& M, GlobalMatrix& K, GlobalVector& b) const = 0;
virtual bool checkBounds(std::vector<double> const& local_x, virtual bool checkBounds(std::vector<double> const& local_x,
std::vector<double> const& local_x_prev_ts) = 0; std::vector<double> const& local_x_prev_ts) = 0;
}; };
...@@ -54,11 +48,11 @@ public: ...@@ -54,11 +48,11 @@ public:
unsigned const integration_order, unsigned const integration_order,
AssemblyParams const& asm_params); AssemblyParams const& asm_params);
void assemble(double const t, std::vector<double> const& local_x) override; void assemble(std::size_t const id,
NumLib::LocalToGlobalIndexMap const& dof_table,
void addToGlobal( double const t, GlobalVector const& x,
NumLib::LocalToGlobalIndexMap::RowColumnIndices const& indices, GlobalMatrix& M, GlobalMatrix& K,
GlobalMatrix& M, GlobalMatrix& K, GlobalVector& b) const override; GlobalVector& b) override;
Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix( Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
const unsigned integration_point) const override const unsigned integration_point) const override
......
...@@ -171,9 +171,6 @@ void TESProcess::initializeConcreteProcess( ...@@ -171,9 +171,6 @@ void TESProcess::initializeConcreteProcess(
NumLib::LocalToGlobalIndexMap const& dof_table, NumLib::LocalToGlobalIndexMap const& dof_table,
MeshLib::Mesh const& mesh, unsigned const integration_order) MeshLib::Mesh const& mesh, unsigned const integration_order)
{ {
DBUG("Create global assembler.");
_global_assembler.reset(new GlobalAssembler(dof_table));
ProcessLib::createLocalAssemblers<TESLocalAssembler>( ProcessLib::createLocalAssemblers<TESLocalAssembler>(
mesh.getDimension(), mesh.getElements(), dof_table, integration_order, mesh.getDimension(), mesh.getElements(), dof_table, integration_order,
_local_assemblers, _assembly_params); _local_assemblers, _assembly_params);
...@@ -249,9 +246,9 @@ void TESProcess::assembleConcreteProcess(const double t, ...@@ -249,9 +246,9 @@ void TESProcess::assembleConcreteProcess(const double t,
DBUG("Assemble TESProcess."); DBUG("Assemble TESProcess.");
// Call global assembler for each local assembly item. // Call global assembler for each local assembly item.
GlobalExecutor::executeMemberDereferenced(*_global_assembler, GlobalExecutor::executeMemberOnDereferenced(
&GlobalAssembler::assemble, &TESLocalAssemblerInterface::assemble, _local_assemblers,
_local_assemblers, t, x, M, K, b); *_local_to_global_index_map, t, x, M, K, b);
} }
void TESProcess::preTimestep(GlobalVector const& x, const double t, void TESProcess::preTimestep(GlobalVector const& x, const double t,
......
...@@ -52,10 +52,6 @@ public: ...@@ -52,10 +52,6 @@ public:
bool isLinear() const override { return false; } bool isLinear() const override { return false; }
private: private:
using GlobalAssembler = NumLib::VectorMatrixAssembler<
TESLocalAssemblerInterface,
NumLib::ODESystemTag::FirstOrderImplicitQuasilinear>;
using ExtrapolatorInterface = using ExtrapolatorInterface =
NumLib::Extrapolator<TESIntPtVariables, TESLocalAssemblerInterface>; NumLib::Extrapolator<TESIntPtVariables, TESLocalAssemblerInterface>;
using ExtrapolatorImplementation = using ExtrapolatorImplementation =
...@@ -85,7 +81,6 @@ private: ...@@ -85,7 +81,6 @@ private:
NumLib::LocalToGlobalIndexMap const& dof_table, NumLib::LocalToGlobalIndexMap const& dof_table,
std::unique_ptr<GlobalVector>& result_cache); std::unique_ptr<GlobalVector>& result_cache);
std::unique_ptr<GlobalAssembler> _global_assembler;
std::vector<std::unique_ptr<TESLocalAssemblerInterface>> _local_assemblers; std::vector<std::unique_ptr<TESLocalAssemblerInterface>> _local_assemblers;
AssemblyParams _assembly_params; AssemblyParams _assembly_params;
......
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