/** * \file * \copyright * Copyright (c) 2012-2021, OpenGeoSys Community (http://www.opengeosys.org) * Distributed under a Modified BSD License. * See accompanying file LICENSE.txt or * http://www.opengeosys.org/project/license * */ #pragma once #include #include #include "IntegrationPointData.h" #include "LocalAssemblerInterface.h" #include "MaterialLib/PhysicalConstant.h" #include "MaterialLib/SolidModels/LinearElasticIsotropic.h" #include "MathLib/KelvinVector.h" #include "MathLib/LinAlg/Eigen/EigenMapTools.h" #include "NumLib/DOF/DOFTableUtil.h" #include "NumLib/Fem/InitShapeMatrices.h" #include "NumLib/Fem/ShapeMatrixPolicy.h" #include "ParameterLib/Parameter.h" #include "ProcessLib/Deformation/BMatrixPolicy.h" #include "ProcessLib/Deformation/LinearBMatrix.h" #include "ProcessLib/LocalAssemblerTraits.h" #include "ProcessLib/Utils/SetOrGetIntegrationPointData.h" #include "ThermoHydroMechanicsProcessData.h" namespace ProcessLib { namespace ThermoHydroMechanics { namespace MPL = MaterialPropertyLib; /// Used by for extrapolation of the integration point values. It is ordered /// (and stored) by integration points. template struct SecondaryData { std::vector> N_u; }; template class ThermoHydroMechanicsLocalAssembler : public LocalAssemblerInterface { public: using ShapeMatricesTypeDisplacement = ShapeMatrixPolicyType; // Types for pressure. using ShapeMatricesTypePressure = ShapeMatrixPolicyType; using GlobalDimMatrixType = typename ShapeMatricesTypePressure::GlobalDimMatrixType; using GlobalDimVectorType = typename ShapeMatricesTypePressure::GlobalDimVectorType; static int const KelvinVectorSize = MathLib::KelvinVector::KelvinVectorDimensions::value; using Invariants = MathLib::KelvinVector::Invariants; using SymmetricTensor = Eigen::Matrix; ThermoHydroMechanicsLocalAssembler( ThermoHydroMechanicsLocalAssembler const&) = delete; ThermoHydroMechanicsLocalAssembler(ThermoHydroMechanicsLocalAssembler&&) = delete; ThermoHydroMechanicsLocalAssembler( MeshLib::Element const& e, std::size_t const /*local_matrix_size*/, bool const is_axially_symmetric, unsigned const integration_order, ThermoHydroMechanicsProcessData& process_data); void assemble(double const /*t*/, double const /*dt*/, std::vector const& /*local_x*/, std::vector const& /*local_xdot*/, std::vector& /*local_M_data*/, std::vector& /*local_K_data*/, std::vector& /*local_rhs_data*/) override { OGS_FATAL( "ThermoHydroMechanicsLocalAssembler: assembly without Jacobian is " "not implemented."); } void assembleWithJacobian(double const t, double const dt, std::vector const& local_x, std::vector const& local_xdot, const double /*dxdot_dx*/, const double /*dx_dx*/, std::vector& /*local_M_data*/, std::vector& /*local_K_data*/, std::vector& local_rhs_data, std::vector& local_Jac_data) override; void initializeConcrete() override { unsigned const n_integration_points = _integration_method.getNumberOfPoints(); for (unsigned ip = 0; ip < n_integration_points; ip++) { _ip_data[ip].pushBackState(); } } void postTimestepConcrete(Eigen::VectorXd const& /*local_x*/, double const /*t*/, double const /*dt*/) override { unsigned const n_integration_points = _integration_method.getNumberOfPoints(); for (unsigned ip = 0; ip < n_integration_points; ip++) { _ip_data[ip].pushBackState(); } } void computeSecondaryVariableConcrete( double const t, double const dt, Eigen::VectorXd const& local_x, Eigen::VectorXd const& local_x_dot) override; void postNonLinearSolverConcrete(std::vector const& local_x, std::vector const& local_xdot, double const t, double const dt, bool const use_monolithic_scheme, int const process_id) override; Eigen::Map getShapeMatrix( const unsigned integration_point) const override { auto const& N_u = _secondary_data.N_u[integration_point]; // assumes N is stored contiguously in memory return Eigen::Map(N_u.data(), N_u.size()); } std::vector const& getIntPtDarcyVelocity( const double t, std::vector const& x, std::vector const& dof_table, std::vector& cache) const override; private: std::size_t setSigma(double const* values) { return ProcessLib::setIntegrationPointKelvinVectorData( values, _ip_data, &IpData::sigma_eff); } // TODO (naumov) This method is same as getIntPtSigma but for arguments and // the ordering of the cache_mat. // There should be only one. std::vector getSigma() const override { return ProcessLib::getIntegrationPointKelvinVectorData( _ip_data, &IpData::sigma_eff); } std::vector const& getIntPtSigma( const double /*t*/, std::vector const& /*x*/, std::vector const& /*dof_table*/, std::vector& cache) const override { return ProcessLib::getIntegrationPointKelvinVectorData( _ip_data, &IpData::sigma_eff, cache); } virtual std::vector const& getIntPtEpsilon( const double /*t*/, std::vector const& /*x*/, std::vector const& /*dof_table*/, std::vector& cache) const override { return ProcessLib::getIntegrationPointKelvinVectorData( _ip_data, &IpData::eps, cache); } private: ThermoHydroMechanicsProcessData& _process_data; using BMatricesType = BMatrixPolicyType; using IpData = IntegrationPointData; std::vector> _ip_data; IntegrationMethod _integration_method; MeshLib::Element const& _element; bool const _is_axially_symmetric; SecondaryData< typename ShapeMatricesTypeDisplacement::ShapeMatrices::ShapeType> _secondary_data; // The shape function of pressure has the same form with the shape function // of temperature static const int temperature_index = 0; static const int temperature_size = ShapeFunctionPressure::NPOINTS; static const int pressure_index = ShapeFunctionPressure::NPOINTS; static const int pressure_size = ShapeFunctionPressure::NPOINTS; static const int displacement_index = ShapeFunctionPressure::NPOINTS * 2; static const int displacement_size = ShapeFunctionDisplacement::NPOINTS * DisplacementDim; }; } // namespace ThermoHydroMechanics } // namespace ProcessLib #include "ThermoHydroMechanicsFEM-impl.h"