Skip to content

Commit

Permalink
[PL/LF] Use shape function from matrix cache
Browse files Browse the repository at this point in the history
  • Loading branch information
TomFischer committed Mar 27, 2024
1 parent 5d344ea commit 2b5594e
Show file tree
Hide file tree
Showing 2 changed files with 17 additions and 10 deletions.
12 changes: 9 additions & 3 deletions ProcessLib/LiquidFlow/LiquidFlowLocalAssembler-impl.h
Original file line number Diff line number Diff line change
Expand Up @@ -143,12 +143,15 @@ void LiquidFlowLocalAssembler<ShapeFunction, GlobalDim>::
_process_data.element_rotation_matrices[_element.getID()].transpose() *
_process_data.specific_body_force;

auto const& N = _shape_matrix_cache
.NsHigherOrder<typename ShapeFunction::MeshElement>();

for (unsigned ip = 0; ip < n_integration_points; ip++)
{
auto const& ip_data = _ip_data[ip];

double p = 0.;
NumLib::shapeFunctionInterpolate(local_x, ip_data.N, p);
NumLib::shapeFunctionInterpolate(local_x, N[ip], p);
vars.liquid_phase_pressure = p;

// Compute density:
Expand All @@ -173,7 +176,7 @@ void LiquidFlowLocalAssembler<ShapeFunction, GlobalDim>::
// Assemble mass matrix, M
local_M.noalias() +=
(porosity * ddensity_dpressure / fluid_density + storage) *
ip_data.N.transpose() * ip_data.N * ip_data.integration_weight;
N[ip].transpose() * N[ip] * ip_data.integration_weight;

// Compute viscosity:
auto const viscosity =
Expand Down Expand Up @@ -292,11 +295,14 @@ void LiquidFlowLocalAssembler<ShapeFunction, GlobalDim>::
_process_data.element_rotation_matrices[_element.getID()].transpose() *
_process_data.specific_body_force;

auto const& N = _shape_matrix_cache
.NsHigherOrder<typename ShapeFunction::MeshElement>();

for (unsigned ip = 0; ip < n_integration_points; ip++)
{
auto const& ip_data = _ip_data[ip];
double p = 0.;
NumLib::shapeFunctionInterpolate(local_x, ip_data.N, p);
NumLib::shapeFunctionInterpolate(local_x, N[ip], p);
vars.liquid_phase_pressure = p;

// Compute density:
Expand Down
15 changes: 8 additions & 7 deletions ProcessLib/LiquidFlow/LiquidFlowLocalAssembler.h
Original file line number Diff line number Diff line change
Expand Up @@ -37,14 +37,12 @@ namespace LiquidFlow
template <typename NodalRowVectorType, typename GlobalDimNodalMatrixType>
struct IntegrationPointData final
{
explicit IntegrationPointData(NodalRowVectorType const& N_,
GlobalDimNodalMatrixType const& dNdx_,
explicit IntegrationPointData(GlobalDimNodalMatrixType const& dNdx_,
double const& integration_weight_)
: N(N_), dNdx(dNdx_), integration_weight(integration_weight_)
: dNdx(dNdx_), integration_weight(integration_weight_)
{
}

NodalRowVectorType const N;
GlobalDimNodalMatrixType const dNdx;
double const integration_weight;

Expand Down Expand Up @@ -115,7 +113,7 @@ class LiquidFlowLocalAssembler : public LiquidFlowLocalAssemblerInterface
for (unsigned ip = 0; ip < n_integration_points; ip++)
{
_ip_data.emplace_back(
shape_matrices[ip].N, shape_matrices[ip].dNdx,
shape_matrices[ip].dNdx,
_integration_method.getWeightedPoint(ip).getWeight() *
shape_matrices[ip].integralMeasure *
shape_matrices[ip].detJ * aperture_size);
Expand All @@ -138,10 +136,13 @@ class LiquidFlowLocalAssembler : public LiquidFlowLocalAssemblerInterface
Eigen::Map<const Eigen::RowVectorXd> getShapeMatrix(
const unsigned integration_point) const override
{
auto const& N = _ip_data[integration_point].N;
auto const& N =
_shape_matrix_cache
.NsHigherOrder<typename ShapeFunction::MeshElement>();

// assumes N is stored contiguously in memory
return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
return Eigen::Map<const Eigen::RowVectorXd>(
N[integration_point].data(), N[integration_point].size());
}

std::vector<double> const& getIntPtDarcyVelocity(
Expand Down

0 comments on commit 2b5594e

Please sign in to comment.