/*! * \file TFEL/Material/StationaryHeatTransfer.hxx * \brief this file implements the StationaryHeatTransfer Behaviour. * File generated by tfel version 3.3.0 * \author Thomas Helfer * \date 15 / 02 / 2019 */ #ifndef LIB_TFELMATERIAL_STATIONARYHEATTRANSFER_HXX #define LIB_TFELMATERIAL_STATIONARYHEATTRANSFER_HXX #include #include #include #include #include #include"TFEL/Raise.hxx" #include"TFEL/PhysicalConstants.hxx" #include"TFEL/Config/TFELConfig.hxx" #include"TFEL/Config/TFELTypes.hxx" #include"TFEL/Metaprogramming/StaticAssert.hxx" #include"TFEL/TypeTraits/IsFundamentalNumericType.hxx" #include"TFEL/TypeTraits/IsReal.hxx" #include"TFEL/Math/General/IEEE754.hxx" #include"TFEL/Math/Vector/TVectorView.hxx" #include"TFEL/Math/Matrix/TMatrixView.hxx" #include"TFEL/Material/MaterialException.hxx" #include"TFEL/Material/MechanicalBehaviour.hxx" #include"TFEL/Material/MechanicalBehaviourTraits.hxx" #include"TFEL/Material/OutOfBoundsPolicy.hxx" #include"TFEL/Material/BoundsCheck.hxx" #include"TFEL/Material/IsotropicPlasticity.hxx" #include"TFEL/Material/Lame.hxx" #include"TFEL/Material/Hosford1972YieldCriterion.hxx" #include"TFEL/Material/StationaryHeatTransferBehaviourData.hxx" #include"TFEL/Material/StationaryHeatTransferIntegrationData.hxx" #include "MFront/GenericBehaviour/State.hxx" #include "MFront/GenericBehaviour/BehaviourData.hxx" namespace tfel{ namespace material{ struct StationaryHeatTransferParametersInitializer { static StationaryHeatTransferParametersInitializer& get(); double A; double B; double minimal_time_step_scaling_factor; double maximal_time_step_scaling_factor; void set(const char* const,const double); /*! * \brief convert a string to double * \param[in] p : parameter * \param[in] v : value */ static double getDouble(const std::string&,const std::string&); private : StationaryHeatTransferParametersInitializer(); StationaryHeatTransferParametersInitializer(const StationaryHeatTransferParametersInitializer&); StationaryHeatTransferParametersInitializer& operator=(const StationaryHeatTransferParametersInitializer&); /*! * \brief read the parameters from the given file * \param[out] pi : parameters initializer * \param[in] fn : file name */ static void readParameters(StationaryHeatTransferParametersInitializer&,const char* const); }; //! \brief forward declaration template class StationaryHeatTransfer; //! \brief forward declaration template std::ostream& operator <<(std::ostream&,const StationaryHeatTransfer&); /*! * \class StationaryHeatTransfer * \brief This class implements the StationaryHeatTransfer behaviour. * \param hypothesis, modelling hypothesis. * \param Type, numerical type. * \author Thomas Helfer * \date 15 / 02 / 2019 */ template class StationaryHeatTransfer final : public MechanicalBehaviour, public StationaryHeatTransferBehaviourData, public StationaryHeatTransferIntegrationData { static constexpr unsigned short N = ModellingHypothesisToSpaceDimension::value; TFEL_STATIC_ASSERT(N==1||N==2||N==3); TFEL_STATIC_ASSERT(tfel::typetraits::IsFundamentalNumericType::cond); TFEL_STATIC_ASSERT(tfel::typetraits::IsReal::cond); friend std::ostream& operator<< <>(std::ostream&,const StationaryHeatTransfer&); static constexpr unsigned short TVectorSize = N; typedef tfel::math::StensorDimeToSize StensorDimeToSize; static constexpr unsigned short StensorSize = StensorDimeToSize::value; typedef tfel::math::TensorDimeToSize TensorDimeToSize; static constexpr unsigned short TensorSize = TensorDimeToSize::value; using ushort = unsigned short; using Types = tfel::config::Types; using real = typename Types::real; using time = typename Types::time; using length = typename Types::length; using frequency = typename Types::frequency; using stress = typename Types::stress; using strain = typename Types::strain; using strainrate = typename Types::strainrate; using stressrate = typename Types::stressrate; using temperature = typename Types::temperature; using thermalexpansion = typename Types::thermalexpansion; using thermalconductivity = typename Types::thermalconductivity; using massdensity = typename Types::massdensity; using energydensity = typename Types::energydensity; using TVector = typename Types::TVector; using Stensor = typename Types::Stensor; using Stensor4 = typename Types::Stensor4; using FrequencyStensor = typename Types::FrequencyStensor; using ForceTVector = typename Types::ForceTVector; using StressStensor = typename Types::StressStensor; using StressRateStensor = typename Types::StressRateStensor; using DisplacementTVector = typename Types::DisplacementTVector; using StrainStensor = typename Types::StrainStensor; using StrainRateStensor = typename Types::StrainRateStensor; using StiffnessTensor = typename Types::StiffnessTensor; using Tensor = typename Types::Tensor; using FrequencyTensor = typename Types::FrequencyTensor; using StressTensor = typename Types::StressTensor; using ThermalExpansionCoefficientTensor = typename Types::ThermalExpansionCoefficientTensor; using DeformationGradientTensor = typename Types::DeformationGradientTensor; using DeformationGradientRateTensor = typename Types::DeformationGradientRateTensor; using TemperatureGradient = typename Types::TemperatureGradient; using HeatFlux = typename Types::HeatFlux; using TangentOperator = tfel::math::tvector<(TVectorSize)*(TVectorSize)+(TVectorSize)*(1),real>; using PhysicalConstants = tfel::PhysicalConstants; public : typedef StationaryHeatTransferBehaviourData BehaviourData; typedef StationaryHeatTransferIntegrationData IntegrationData; typedef typename MechanicalBehaviour::SMFlag SMFlag; typedef typename MechanicalBehaviour::SMType SMType; using MechanicalBehaviour::ELASTIC; using MechanicalBehaviour::SECANTOPERATOR; using MechanicalBehaviour::TANGENTOPERATOR; using MechanicalBehaviour::CONSISTENTTANGENTOPERATOR; using MechanicalBehaviour::NOSTIFFNESSREQUESTED; using IntegrationResult = typename MechanicalBehaviour::IntegrationResult; using MechanicalBehaviour::SUCCESS; using MechanicalBehaviour::FAILURE; using MechanicalBehaviour::UNRELIABLE_RESULTS; private : #line 19 "StationaryHeatTransfer.mfront" thermalconductivity k; #line 16 "StationaryHeatTransfer.mfront" real A; #line 17 "StationaryHeatTransfer.mfront" real B; real minimal_time_step_scaling_factor; real maximal_time_step_scaling_factor; //! Tangent operator; TangentOperator Dt; tfel::math::TMatrixView dj_ddgT; tfel::math::TVectorView dj_ddT; /*! * \brief Update internal variables at end of integration */ void updateIntegrationVariables() {} /*! * \brief Update internal variables at end of integration */ void updateStateVariables() {} /*! * \brief Update auxiliary state variables at end of integration */ void updateAuxiliaryStateVariables() {} //! \brief Default constructor (disabled) StationaryHeatTransfer() =delete ; //! \brief Copy constructor (disabled) StationaryHeatTransfer(const StationaryHeatTransfer&) = delete; //! \brief Assignement operator (disabled) StationaryHeatTransfer& operator = (const StationaryHeatTransfer&) = delete; public: /*! * \brief Constructor */ StationaryHeatTransfer(const StationaryHeatTransferBehaviourData& src1, const StationaryHeatTransferIntegrationData& src2) : StationaryHeatTransferBehaviourData(src1), StationaryHeatTransferIntegrationData(src2), dj_ddgT(Dt.begin()), dj_ddT(Dt.begin()+TVectorSize*TVectorSize) { using namespace std; using namespace tfel::math; using std::vector; this->A = StationaryHeatTransferParametersInitializer::get().A; this->B = StationaryHeatTransferParametersInitializer::get().B; this->minimal_time_step_scaling_factor = StationaryHeatTransferParametersInitializer::get().minimal_time_step_scaling_factor; this->maximal_time_step_scaling_factor = StationaryHeatTransferParametersInitializer::get().maximal_time_step_scaling_factor; } /* * \brief constructor for the Generic interface * \param[in] mgb_d: behaviour data */ StationaryHeatTransfer(const mfront::gb::BehaviourData& mgb_d) : StationaryHeatTransferBehaviourData(mgb_d), StationaryHeatTransferIntegrationData(mgb_d), dj_ddgT(Dt.begin()), dj_ddT(Dt.begin()+TVectorSize*TVectorSize) { using namespace std; using namespace tfel::math; using std::vector; this->A = StationaryHeatTransferParametersInitializer::get().A; this->B = StationaryHeatTransferParametersInitializer::get().B; this->minimal_time_step_scaling_factor = StationaryHeatTransferParametersInitializer::get().minimal_time_step_scaling_factor; this->maximal_time_step_scaling_factor = StationaryHeatTransferParametersInitializer::get().maximal_time_step_scaling_factor; tfel::fsalgo::copy::exe(mgb_d.s0.gradients,this->gT.begin()); tfel::fsalgo::transform::exe(mgb_d.s1.gradients,mgb_d.s0.gradients,this->dgT.begin(),std::minus()); tfel::fsalgo::copy::exe(mgb_d.s0.thermodynamic_forces,this->j.begin()); } /*! * \ brief initialize the behaviour with user code */ void initialize(){ using namespace std; using namespace tfel::math; using std::vector; } /*! * \brief set the policy for "out of bounds" conditions */ void setOutOfBoundsPolicy(const OutOfBoundsPolicy policy_value){ this->policy = policy_value; } // end of setOutOfBoundsPolicy /*! * \return the modelling hypothesis */ constexpr ModellingHypothesis::Hypothesis getModellingHypothesis() const{ return hypothesis; } // end of getModellingHypothesis /*! * \brief check bounds */ void checkBounds() const{ } // end of checkBounds IntegrationResult computePredictionOperator(const SMFlag,const SMType) override{ tfel::raise("StationaryHeatTransfer::computePredictionOperator: " "unsupported prediction operator flag"); } real getMinimalTimeStepScalingFactor() const override{ return this->minimal_time_step_scaling_factor; } std::pair computeAPrioriTimeStepScalingFactor(const real current_time_step_scaling_factor) const override{ const auto time_scaling_factor = this->computeAPrioriTimeStepScalingFactorII(); return {time_scaling_factor.first, std::min(std::min(std::max(time_scaling_factor.second, this->minimal_time_step_scaling_factor), this->maximal_time_step_scaling_factor), current_time_step_scaling_factor)}; } /*! * \brief Integrate behaviour over the time step */ IntegrationResult integrate(const SMFlag smflag, const SMType smt) override{ using namespace std; using namespace tfel::math; raise_if(smflag!=MechanicalBehaviour::STANDARDTANGENTOPERATOR, "invalid tangent operator flag"); bool computeTangentOperator_ = smt!=NOSTIFFNESSREQUESTED; #line 22 "StationaryHeatTransfer.mfront" const auto T_ = this->T + this->dT; #line 23 "StationaryHeatTransfer.mfront" this->k = 1 / (this->A + this->B * T_); #line 24 "StationaryHeatTransfer.mfront" this->j = this->k * (this->gT + this->dgT); this->updateIntegrationVariables(); this->updateStateVariables(); this->updateAuxiliaryStateVariables(); if(computeTangentOperator_){ if(!this->computeConsistentTangentOperator(smt)){ return MechanicalBehaviour::FAILURE; } } return MechanicalBehaviour::SUCCESS; } std::pair computeAPosterioriTimeStepScalingFactor(const real current_time_step_scaling_factor) const override{ const auto time_scaling_factor = this->computeAPosterioriTimeStepScalingFactorII(); return {time_scaling_factor.first, std::min(std::min(std::max(time_scaling_factor.second, this->minimal_time_step_scaling_factor), this->maximal_time_step_scaling_factor), current_time_step_scaling_factor)}; } /*! * \brief Update the internal energy at end of the time step * \param[in] Psi_s: internal energy at end of the time step */ void computeInternalEnergy(real& Psi_s) const { Psi_s=0; } /*! * \brief Update the dissipated energy at end of the time step * \param[in] Psi_d: dissipated energy at end of the time step */ void computeDissipatedEnergy(real& Psi_d) const { Psi_d=0; } bool computeConsistentTangentOperator(const SMType smt){ using namespace std; using namespace tfel::math; using std::vector; #line 28 "StationaryHeatTransfer.mfront" this->dj_ddgT = this->k * tmatrix::Id(); #line 29 "StationaryHeatTransfer.mfront" dj_ddT = -this->B * this->k * this->k * (this->gT + this->dgT); return true; } const TangentOperator& getTangentOperator() const{ return this->Dt; } void updateExternalStateVariables(){ this->gT += this->dgT; this->T += this->dT; } //! ~StationaryHeatTransfer() override = default; private: std::pair computeAPrioriTimeStepScalingFactorII() const{ return {true,this->maximal_time_step_scaling_factor}; } std::pair computeAPosterioriTimeStepScalingFactorII() const{ return {true,this->maximal_time_step_scaling_factor}; } //! policy for treating out of bounds conditions OutOfBoundsPolicy policy = None; }; // end of StationaryHeatTransfer class template std::ostream& operator <<(std::ostream& os,const StationaryHeatTransfer& b) { os << "gT : " << b.gT << '\n'; os << "ΔgT : " << b.dgT << '\n'; os << "j : " << b.j << '\n'; os << "Δt : " << b.dt << '\n'; os << "T : " << b.T << '\n'; os << "ΔT : " << b.dT << '\n'; os << "k : " << b.k << '\n'; os << "A : " << b.A << '\n'; os << "B : " << b.B << '\n'; os << "minimal_time_step_scaling_factor : " << b.minimal_time_step_scaling_factor << '\n'; os << "maximal_time_step_scaling_factor : " << b.maximal_time_step_scaling_factor << '\n'; return os; } /*! * Partial specialisation for StationaryHeatTransfer. */ template class MechanicalBehaviourTraits > { static constexpr unsigned short N = ModellingHypothesisToSpaceDimension::value; static constexpr unsigned short TVectorSize = N; typedef tfel::math::StensorDimeToSize StensorDimeToSize; static constexpr unsigned short StensorSize = StensorDimeToSize::value; typedef tfel::math::TensorDimeToSize TensorDimeToSize; static constexpr unsigned short TensorSize = TensorDimeToSize::value; public: static constexpr bool is_defined = true; static constexpr bool use_quantities = false; static constexpr bool hasStressFreeExpansion = false; static constexpr bool handlesThermalExpansion = false; static constexpr unsigned short dimension = N; typedef Type NumType; static constexpr unsigned short material_properties_nb = 0; static constexpr unsigned short internal_variables_nb = 0; static constexpr unsigned short external_variables_nb = 1; static constexpr unsigned short external_variables_nb2 = 0; static constexpr bool hasConsistentTangentOperator = true; static constexpr bool isConsistentTangentOperatorSymmetric = false; static constexpr bool hasPredictionOperator = false; static constexpr bool hasAPrioriTimeStepScalingFactor = false; static constexpr bool hasComputeInternalEnergy = false; static constexpr bool hasComputeDissipatedEnergy = false; /*! * \return the name of the class. */ static const char* getName(){ return "StationaryHeatTransfer"; } }; /*! * Partial specialisation for StationaryHeatTransfer. */ template class MechanicalBehaviourTraits > { public: static constexpr bool is_defined = false; static constexpr bool use_quantities = false; static constexpr bool hasStressFreeExpansion = false; static constexpr bool handlesThermalExpansion = false; static constexpr unsigned short dimension = 0u; typedef Type NumType; static constexpr unsigned short material_properties_nb = 0; static constexpr unsigned short internal_variables_nb = 0; static constexpr unsigned short external_variables_nb = 0; static constexpr unsigned short external_variables_nb2 = 0; static constexpr bool hasConsistentTangentOperator = false; static constexpr bool isConsistentTangentOperatorSymmetric = false; static constexpr bool hasPredictionOperator = false; static constexpr bool hasAPrioriTimeStepScalingFactor = false; static constexpr bool hasComputeInternalEnergy = false; static constexpr bool hasComputeDissipatedEnergy = false; /*! * \return the name of the class. */ static const char* getName(){ return "StationaryHeatTransfer"; } }; /*! * Partial specialisation for StationaryHeatTransfer. */ template class MechanicalBehaviourTraits > { public: static constexpr bool is_defined = false; static constexpr bool use_quantities = false; static constexpr bool hasStressFreeExpansion = false; static constexpr bool handlesThermalExpansion = false; static constexpr unsigned short dimension = 0u; typedef Type NumType; static constexpr unsigned short material_properties_nb = 0; static constexpr unsigned short internal_variables_nb = 0; static constexpr unsigned short external_variables_nb = 0; static constexpr unsigned short external_variables_nb2 = 0; static constexpr bool hasConsistentTangentOperator = false; static constexpr bool isConsistentTangentOperatorSymmetric = false; static constexpr bool hasPredictionOperator = false; static constexpr bool hasAPrioriTimeStepScalingFactor = false; static constexpr bool hasComputeInternalEnergy = false; static constexpr bool hasComputeDissipatedEnergy = false; /*! * \return the name of the class. */ static const char* getName(){ return "StationaryHeatTransfer"; } }; } // end of namespace material } // end of namespace tfel #endif /* LIB_TFELMATERIAL_STATIONARYHEATTRANSFER_HXX */