From 89b5a79a126c36dfa1bf5d46ba5011470dc79502 Mon Sep 17 00:00:00 2001 From: Dmitri Naumov <dmitri.naumov@ufz.de> Date: Fri, 27 Sep 2019 13:00:06 +0200 Subject: [PATCH] [PL] Update staggered processes. The flat vector of staggered solution simplifies the data access. Fixing the location of phasefield/temperature/displacement subvectors in phase field processes. --- .../ComponentTransportFEM.h | 19 ++-- ProcessLib/HT/HTFEM.h | 11 ++- ProcessLib/HT/MonolithicHTFEM.h | 10 +- ProcessLib/HT/StaggeredHTFEM-impl.h | 57 ++++++------ .../HydroMechanics/HydroMechanicsFEM-impl.h | 39 ++++---- ProcessLib/PhaseField/PhaseFieldFEM-impl.h | 92 ++++++------------- ProcessLib/PhaseField/PhaseFieldFEM.h | 22 +++-- .../ThermoMechanicalPhaseFieldFEM-impl.h | 62 +++++-------- .../ThermoMechanicalPhaseFieldFEM.h | 26 ++++-- .../ThermoMechanics/ThermoMechanicsFEM-impl.h | 25 ++--- 10 files changed, 165 insertions(+), 198 deletions(-) diff --git a/ProcessLib/ComponentTransport/ComponentTransportFEM.h b/ProcessLib/ComponentTransport/ComponentTransportFEM.h index ea06d4b2469..581e50b67db 100644 --- a/ProcessLib/ComponentTransport/ComponentTransportFEM.h +++ b/ProcessLib/ComponentTransport/ComponentTransportFEM.h @@ -424,10 +424,9 @@ public: LocalCoupledSolutions const& coupled_xs) { auto local_p = Eigen::Map<const NodalVectorType>( - coupled_xs.local_coupled_xs[hydraulic_process_id].data(), - pressure_size); + &coupled_xs.local_coupled_xs[pressure_index], pressure_size); auto local_C = Eigen::Map<const NodalVectorType>( - coupled_xs.local_coupled_xs[first_transport_process_id].data(), + &coupled_xs.local_coupled_xs[first_concentration_index], concentration_size); auto local_C0 = Eigen::Map<const NodalVectorType>( &coupled_xs.local_coupled_xs0[first_concentration_index], @@ -535,11 +534,12 @@ public: LocalCoupledSolutions const& coupled_xs, int const transport_process_id) { auto local_C = Eigen::Map<const NodalVectorType>( - coupled_xs.local_coupled_xs[transport_process_id].data(), + &coupled_xs.local_coupled_xs[first_concentration_index + + (transport_process_id - 1) * + concentration_size], concentration_size); auto local_p = Eigen::Map<const NodalVectorType>( - coupled_xs.local_coupled_xs[hydraulic_process_id].data(), - pressure_size); + &coupled_xs.local_coupled_xs[pressure_index], pressure_size); auto local_p0 = Eigen::Map<const NodalVectorType>( &coupled_xs.local_coupled_xs0[pressure_index], pressure_size); @@ -733,9 +733,10 @@ public: auto const local_xs = getCurrentLocalSolutions( *(this->_coupled_solutions), indices_of_all_coupled_processes); - auto const local_p = MathLib::toVector(local_xs[hydraulic_process_id]); - auto const local_C = - MathLib::toVector(local_xs[first_transport_process_id]); + auto const local_p = Eigen::Map<const NodalVectorType>( + &local_xs[pressure_index], pressure_size); + auto const local_C = Eigen::Map<const NodalVectorType>( + &local_xs[first_concentration_index], concentration_size); return calculateIntPtDarcyVelocity(t, local_p, local_C, cache); } diff --git a/ProcessLib/HT/HTFEM.h b/ProcessLib/HT/HTFEM.h index 52d85f23ffb..a748da5cdb6 100644 --- a/ProcessLib/HT/HTFEM.h +++ b/ProcessLib/HT/HTFEM.h @@ -251,9 +251,16 @@ protected: } std::vector<double> const& getIntPtDarcyVelocityLocal( - const double t, std::vector<double> const& local_p, - std::vector<double> const& local_T, std::vector<double>& cache) const + const double t, std::vector<double> const& local_x, + std::vector<double>& cache) const { + std::vector<double> local_p{ + local_x.data() + pressure_index, + local_x.data() + pressure_index + pressure_size}; + std::vector<double> local_T{ + local_x.data() + temperature_index, + local_x.data() + temperature_index + temperature_size}; + auto const n_integration_points = _integration_method.getNumberOfPoints(); diff --git a/ProcessLib/HT/MonolithicHTFEM.h b/ProcessLib/HT/MonolithicHTFEM.h index 67423f05ff7..d1b5db88c6f 100644 --- a/ProcessLib/HT/MonolithicHTFEM.h +++ b/ProcessLib/HT/MonolithicHTFEM.h @@ -211,15 +211,9 @@ public: auto const indices = NumLib::getIndices(this->_element.getID(), dof_table); assert(!indices.empty()); - auto local_x = current_solution.get(indices); + auto const& local_x = current_solution.get(indices); - std::vector<double> local_p( - std::make_move_iterator(local_x.begin() + local_x.size() / 2), - std::make_move_iterator(local_x.end())); - // only T is kept in local_x - local_x.erase(local_x.begin() + local_x.size() / 2, local_x.end()); - - return this->getIntPtDarcyVelocityLocal(t, local_p, local_x, cache); + return this->getIntPtDarcyVelocityLocal(t, local_x, cache); } private: diff --git a/ProcessLib/HT/StaggeredHTFEM-impl.h b/ProcessLib/HT/StaggeredHTFEM-impl.h index c0ed2f84435..ba1c705d620 100644 --- a/ProcessLib/HT/StaggeredHTFEM-impl.h +++ b/ProcessLib/HT/StaggeredHTFEM-impl.h @@ -51,25 +51,25 @@ void StaggeredHTFEM<ShapeFunction, IntegrationMethod, GlobalDim>:: std::vector<double>& local_b_data, LocalCoupledSolutions const& coupled_xs) { - auto const& local_p = coupled_xs.local_coupled_xs[_hydraulic_process_id]; - auto const local_matrix_size = local_p.size(); - // This assertion is valid only if all nodal d.o.f. use the same shape - // matrices. - assert(local_matrix_size == ShapeFunction::NPOINTS); - - auto const& local_T1 = - coupled_xs.local_coupled_xs[_heat_transport_process_id]; + auto const local_p = Eigen::Map< + typename ShapeMatricesType::template VectorType<pressure_size> const>( + &coupled_xs.local_coupled_xs[pressure_index], pressure_size); + + auto const local_T1 = + Eigen::Map<typename ShapeMatricesType::template VectorType< + temperature_size> const>( + &coupled_xs.local_coupled_xs[temperature_index], temperature_size); auto const local_T0 = Eigen::Map<typename ShapeMatricesType::template VectorType< temperature_size> const>( &coupled_xs.local_coupled_xs0[temperature_index], temperature_size); auto local_M = MathLib::createZeroedMatrix<LocalMatrixType>( - local_M_data, local_matrix_size, local_matrix_size); + local_M_data, pressure_size, pressure_size); auto local_K = MathLib::createZeroedMatrix<LocalMatrixType>( - local_K_data, local_matrix_size, local_matrix_size); - auto local_b = MathLib::createZeroedVector<LocalVectorType>( - local_b_data, local_matrix_size); + local_K_data, pressure_size, pressure_size); + auto local_b = MathLib::createZeroedVector<LocalVectorType>(local_b_data, + pressure_size); ParameterLib::SpatialPosition pos; pos.setElementID(this->_element.getID()); @@ -188,22 +188,19 @@ void StaggeredHTFEM<ShapeFunction, IntegrationMethod, GlobalDim>:: std::vector<double>& /*local_b_data*/, LocalCoupledSolutions const& coupled_xs) { - auto const& local_p = coupled_xs.local_coupled_xs[_hydraulic_process_id]; - auto const local_matrix_size = local_p.size(); - // This assertion is valid only if all nodal d.o.f. use the same shape - // matrices. - assert(local_matrix_size == ShapeFunction::NPOINTS); + auto const local_p = Eigen::Map< + typename ShapeMatricesType::template VectorType<pressure_size> const>( + &coupled_xs.local_coupled_xs[pressure_index], pressure_size); - auto local_p_Eigen_type = - Eigen::Map<const NodalVectorType>(&local_p[0], local_matrix_size); - - auto const& local_T1 = - coupled_xs.local_coupled_xs[_heat_transport_process_id]; + auto const local_T1 = + Eigen::Map<typename ShapeMatricesType::template VectorType< + temperature_size> const>( + &coupled_xs.local_coupled_xs[temperature_index], temperature_size); auto local_M = MathLib::createZeroedMatrix<LocalMatrixType>( - local_M_data, local_matrix_size, local_matrix_size); + local_M_data, temperature_size, temperature_size); auto local_K = MathLib::createZeroedMatrix<LocalMatrixType>( - local_K_data, local_matrix_size, local_matrix_size); + local_K_data, temperature_size, temperature_size); ParameterLib::SpatialPosition pos; pos.setElementID(this->_element.getID()); @@ -276,9 +273,9 @@ void StaggeredHTFEM<ShapeFunction, IntegrationMethod, GlobalDim>:: intrinsic_permeability / viscosity; GlobalDimVectorType const velocity = process_data.has_gravity - ? GlobalDimVectorType(-K_over_mu * (dNdx * local_p_Eigen_type - - fluid_density * b)) - : GlobalDimVectorType(-K_over_mu * dNdx * local_p_Eigen_type); + ? GlobalDimVectorType(-K_over_mu * + (dNdx * local_p - fluid_density * b)) + : GlobalDimVectorType(-K_over_mu * dNdx * local_p); GlobalDimMatrixType const thermal_conductivity_dispersivity = this->getThermalConductivityDispersivity( @@ -305,12 +302,10 @@ StaggeredHTFEM<ShapeFunction, IntegrationMethod, GlobalDim>:: assert(!indices.empty()); std::vector<std::vector<GlobalIndexType>> indices_of_all_coupled_processes = {indices, indices}; - auto const local_xs = getCurrentLocalSolutions( + auto const& local_xs = getCurrentLocalSolutions( *(this->_coupled_solutions), indices_of_all_coupled_processes); - return this->getIntPtDarcyVelocityLocal( - t, local_xs[_hydraulic_process_id], - local_xs[_heat_transport_process_id], cache); + return this->getIntPtDarcyVelocityLocal(t, local_xs, cache); } } // namespace HT } // namespace ProcessLib diff --git a/ProcessLib/HydroMechanics/HydroMechanicsFEM-impl.h b/ProcessLib/HydroMechanics/HydroMechanicsFEM-impl.h index 49372e36f3a..5102fd192d5 100644 --- a/ProcessLib/HydroMechanics/HydroMechanicsFEM-impl.h +++ b/ProcessLib/HydroMechanics/HydroMechanicsFEM-impl.h @@ -357,25 +357,24 @@ void HydroMechanicsLocalAssembler<ShapeFunctionDisplacement, std::vector<double>& local_b_data, std::vector<double>& local_Jac_data, const LocalCoupledSolutions& local_coupled_solutions) { - auto const& local_p = local_coupled_solutions.local_coupled_xs[0]; - auto const& local_u = local_coupled_solutions.local_coupled_xs[1]; - assert(local_p.size() == pressure_size); - assert(local_u.size() == displacement_size); - - auto const local_matrix_size = local_p.size(); auto local_rhs = MathLib::createZeroedVector<typename ShapeMatricesTypeDisplacement:: template VectorType<pressure_size>>( - local_b_data, local_matrix_size); + local_b_data, pressure_size); ParameterLib::SpatialPosition pos; pos.setElementID(this->_element.getID()); - auto p = Eigen::Map<typename ShapeMatricesTypePressure::template VectorType< - pressure_size> const>(local_p.data(), pressure_size); - auto u = + auto const p = + Eigen::Map<typename ShapeMatricesTypePressure::template VectorType< + pressure_size> const>( + &local_coupled_solutions.local_coupled_xs[pressure_index], + pressure_size); + auto const u = Eigen::Map<typename ShapeMatricesTypeDisplacement::template VectorType< - displacement_size> const>(local_u.data(), displacement_size); + displacement_size> const>( + &local_coupled_solutions.local_coupled_xs[displacement_index], + displacement_size); auto p_dot = Eigen::Map<typename ShapeMatricesTypePressure::template VectorType< @@ -469,14 +468,16 @@ void HydroMechanicsLocalAssembler<ShapeFunctionDisplacement, std::vector<double>& local_b_data, std::vector<double>& local_Jac_data, const LocalCoupledSolutions& local_coupled_solutions) { - auto const& local_p = local_coupled_solutions.local_coupled_xs[0]; - auto const& local_u = local_coupled_solutions.local_coupled_xs[1]; - assert(local_p.size() == pressure_size); - assert(local_u.size() == displacement_size); - - auto u = + auto const p = + Eigen::Map<typename ShapeMatricesTypePressure::template VectorType< + pressure_size> const>( + &local_coupled_solutions.local_coupled_xs[pressure_index], + pressure_size); + auto const u = Eigen::Map<typename ShapeMatricesTypeDisplacement::template VectorType< - displacement_size> const>(local_u.data(), displacement_size); + displacement_size> const>( + &local_coupled_solutions.local_coupled_xs[displacement_index], + displacement_size); auto local_Jac = MathLib::createZeroedMatrix< typename ShapeMatricesTypeDisplacement::template MatrixType< @@ -534,7 +535,7 @@ void HydroMechanicsLocalAssembler<ShapeFunctionDisplacement, local_Jac.noalias() += B.transpose() * C * B * w; double p_at_xi = 0.; - NumLib::shapeFunctionInterpolate(local_p, N_p, p_at_xi); + NumLib::shapeFunctionInterpolate(p, N_p, p_at_xi); double const rho = rho_sr * (1 - porosity) + porosity * rho_fr; local_rhs.noalias() -= diff --git a/ProcessLib/PhaseField/PhaseFieldFEM-impl.h b/ProcessLib/PhaseField/PhaseFieldFEM-impl.h index 5098caa42ea..e6ace4fba0a 100644 --- a/ProcessLib/PhaseField/PhaseFieldFEM-impl.h +++ b/ProcessLib/PhaseField/PhaseFieldFEM-impl.h @@ -55,31 +55,25 @@ void PhaseFieldLocalAssembler<ShapeFunction, IntegrationMethod, std::vector<double>& local_b_data, std::vector<double>& local_Jac_data, LocalCoupledSolutions const& local_coupled_solutions) { - using DeformationVector = - typename ShapeMatricesType::template VectorType<displacement_size>; using DeformationMatrix = typename ShapeMatricesType::template MatrixType<displacement_size, displacement_size>; - using PhaseFieldVector = - typename ShapeMatricesType::template VectorType<phasefield_size>; - auto const& local_u = local_coupled_solutions.local_coupled_xs[0]; - auto const& local_d = local_coupled_solutions.local_coupled_xs[1]; - assert(local_u.size() == displacement_size); - assert(local_d.size() == phasefield_size); + assert(local_coupled_solutions.local_coupled_xs.size() == + phasefield_size + displacement_size); - auto const local_matrix_size = local_u.size(); - auto d = - Eigen::Map<PhaseFieldVector const>(local_d.data(), phasefield_size); - - auto u = - Eigen::Map<DeformationVector const>(local_u.data(), displacement_size); + auto const d = Eigen::Map<PhaseFieldVector const>( + &local_coupled_solutions.local_coupled_xs[phasefield_index], + phasefield_size); + auto const u = Eigen::Map<DeformationVector const>( + &local_coupled_solutions.local_coupled_xs[displacement_index], + displacement_size); auto local_Jac = MathLib::createZeroedMatrix<DeformationMatrix>( - local_Jac_data, local_matrix_size, local_matrix_size); + local_Jac_data, displacement_size, displacement_size); auto local_rhs = MathLib::createZeroedVector<DeformationVector>( - local_b_data, local_matrix_size); + local_b_data, displacement_size); ParameterLib::SpatialPosition x_position; x_position.setElementID(_element.getID()); @@ -158,29 +152,17 @@ void PhaseFieldLocalAssembler<ShapeFunction, IntegrationMethod, std::vector<double>& local_b_data, std::vector<double>& local_Jac_data, LocalCoupledSolutions const& local_coupled_solutions) { - using DeformationVector = - typename ShapeMatricesType::template VectorType<displacement_size>; - using PhaseFieldVector = - typename ShapeMatricesType::template VectorType<phasefield_size>; - using PhaseFieldMatrix = - typename ShapeMatricesType::template MatrixType<phasefield_size, - phasefield_size>; - - auto const& local_u = local_coupled_solutions.local_coupled_xs[0]; - auto const& local_d = local_coupled_solutions.local_coupled_xs[1]; - assert(local_u.size() == displacement_size); - assert(local_d.size() == phasefield_size); - - auto const local_matrix_size = local_d.size(); - auto d = - Eigen::Map<PhaseFieldVector const>(local_d.data(), phasefield_size); - auto u = - Eigen::Map<DeformationVector const>(local_u.data(), displacement_size); + auto const d = Eigen::Map<PhaseFieldVector const>( + &local_coupled_solutions.local_coupled_xs[phasefield_index], + phasefield_size); + auto const u = Eigen::Map<DeformationVector const>( + &local_coupled_solutions.local_coupled_xs[displacement_index], + displacement_size); auto local_Jac = MathLib::createZeroedMatrix<PhaseFieldMatrix>( - local_Jac_data, local_matrix_size, local_matrix_size); + local_Jac_data, phasefield_size, phasefield_size); auto local_rhs = MathLib::createZeroedVector<PhaseFieldVector>( - local_b_data, local_matrix_size); + local_b_data, phasefield_size); ParameterLib::SpatialPosition x_position; x_position.setElementID(_element.getID()); @@ -281,20 +263,12 @@ void PhaseFieldLocalAssembler<ShapeFunction, IntegrationMethod, auto local_coupled_xs = getCurrentLocalSolutions(*cpl_xs, indices_of_processes); - assert(local_coupled_xs.size() == 2); - - auto const& local_u = local_coupled_xs[0]; - auto const& local_d = local_coupled_xs[1]; - - assert(local_u.size() == displacement_size); - assert(local_d.size() == phasefield_size); - - auto d = Eigen::Map< - typename ShapeMatricesType::template VectorType<phasefield_size> const>( - local_d.data(), phasefield_size); + assert(local_coupled_xs.size() == displacement_size + phasefield_size); - auto u = Eigen::Map<typename ShapeMatricesType::template VectorType< - displacement_size> const>(local_u.data(), displacement_size); + auto const d = Eigen::Map<PhaseFieldVector const>( + &local_coupled_xs[phasefield_index], phasefield_size); + auto const u = Eigen::Map<DeformationVector const>( + &local_coupled_xs[displacement_index], displacement_size); ParameterLib::SpatialPosition x_position; x_position.setElementID(_element.getID()); @@ -346,22 +320,14 @@ void PhaseFieldLocalAssembler<ShapeFunction, IntegrationMethod, return NumLib::getIndices(mesh_item_id, dof_table); }); - auto local_coupled_xs = + auto const local_coupled_xs = getCurrentLocalSolutions(*cpl_xs, indices_of_processes); - assert(local_coupled_xs.size() == 2); - - auto const& local_u = local_coupled_xs[0]; - auto const& local_d = local_coupled_xs[1]; - - assert(local_u.size() == displacement_size); - assert(local_d.size() == phasefield_size); - - auto d = Eigen::Map< - typename ShapeMatricesType::template VectorType<phasefield_size> const>( - local_d.data(), phasefield_size); + assert(local_coupled_xs.size() == displacement_size + phasefield_size); - auto u = Eigen::Map<typename ShapeMatricesType::template VectorType< - displacement_size> const>(local_u.data(), displacement_size); + auto const d = Eigen::Map<PhaseFieldVector const>( + &local_coupled_xs[phasefield_index], phasefield_size); + auto const u = Eigen::Map<DeformationVector const>( + &local_coupled_xs[displacement_index], displacement_size); ParameterLib::SpatialPosition x_position; x_position.setElementID(_element.getID()); diff --git a/ProcessLib/PhaseField/PhaseFieldFEM.h b/ProcessLib/PhaseField/PhaseFieldFEM.h index 9727fd3346c..14bd3e546c1 100644 --- a/ProcessLib/PhaseField/PhaseFieldFEM.h +++ b/ProcessLib/PhaseField/PhaseFieldFEM.h @@ -107,6 +107,14 @@ template <typename ShapeFunction, typename IntegrationMethod, int DisplacementDim> class PhaseFieldLocalAssembler : public PhaseFieldLocalAssemblerInterface { +private: + static constexpr int displacement_index = 0; + static constexpr int displacement_size = + ShapeFunction::NPOINTS * DisplacementDim; + static constexpr int phasefield_index = + displacement_index + displacement_size; + static constexpr int phasefield_size = ShapeFunction::NPOINTS; + public: using ShapeMatricesType = ShapeMatrixPolicyType<ShapeFunction, DisplacementDim>; @@ -117,6 +125,14 @@ public: using NodalForceVectorType = typename BMatricesType::NodalForceVectorType; + using DeformationVector = + typename ShapeMatricesType::template VectorType<displacement_size>; + using PhaseFieldVector = + typename ShapeMatricesType::template VectorType<phasefield_size>; + using PhaseFieldMatrix = + typename ShapeMatricesType::template MatrixType<phasefield_size, + phasefield_size>; + PhaseFieldLocalAssembler(PhaseFieldLocalAssembler const&) = delete; PhaseFieldLocalAssembler(PhaseFieldLocalAssembler&&) = delete; @@ -337,12 +353,6 @@ private: MeshLib::Element const& _element; SecondaryData<typename ShapeMatrices::ShapeType> _secondary_data; bool const _is_axially_symmetric; - - static const int phasefield_index = 0; - static const int phasefield_size = ShapeFunction::NPOINTS; - static const int displacement_index = ShapeFunction::NPOINTS; - static const int displacement_size = - ShapeFunction::NPOINTS * DisplacementDim; }; } // namespace PhaseField diff --git a/ProcessLib/ThermoMechanicalPhaseField/ThermoMechanicalPhaseFieldFEM-impl.h b/ProcessLib/ThermoMechanicalPhaseField/ThermoMechanicalPhaseFieldFEM-impl.h index 5e2cdf6b12a..32584dc29b4 100644 --- a/ProcessLib/ThermoMechanicalPhaseField/ThermoMechanicalPhaseFieldFEM-impl.h +++ b/ProcessLib/ThermoMechanicalPhaseField/ThermoMechanicalPhaseFieldFEM-impl.h @@ -64,25 +64,18 @@ void ThermoMechanicalPhaseFieldLocalAssembler<ShapeFunction, IntegrationMethod, std::vector<double>& local_b_data, std::vector<double>& local_Jac_data, LocalCoupledSolutions const& local_coupled_solutions) { - auto const& local_d = - local_coupled_solutions.local_coupled_xs[_phase_field_process_id]; - auto const& local_u = - local_coupled_solutions.local_coupled_xs[_mechanics_related_process_id]; - auto const& local_T = - local_coupled_solutions.local_coupled_xs[_heat_conduction_process_id]; - assert(local_T.size() == temperature_size); - assert(local_d.size() == phasefield_size); - assert(local_u.size() == displacement_size); - - auto d = Eigen::Map< - typename ShapeMatricesType::template VectorType<phasefield_size> const>( - local_d.data(), phasefield_size); - - auto u = Eigen::Map<typename ShapeMatricesType::template VectorType< - displacement_size> const>(local_u.data(), displacement_size); - - auto T = Eigen::Map<typename ShapeMatricesType::template VectorType< - temperature_size> const>(local_T.data(), temperature_size); + assert(local_coupled_solutions.local_coupled_xs.size() == + phasefield_size + displacement_size + temperature_size); + + auto const d = Eigen::Map<PhaseFieldVector const>( + &local_coupled_solutions.local_coupled_xs[phasefield_index], + phasefield_size); + auto const u = Eigen::Map<DeformationVector const>( + &local_coupled_solutions.local_coupled_xs[displacement_index], + displacement_size); + auto const T = Eigen::Map<TemperatureVector const>( + &local_coupled_solutions.local_coupled_xs[temperature_index], + temperature_size); auto local_Jac = MathLib::createZeroedMatrix< typename ShapeMatricesType::template MatrixType<displacement_size, @@ -172,19 +165,15 @@ void ThermoMechanicalPhaseFieldLocalAssembler<ShapeFunction, IntegrationMethod, std::vector<double>& local_b_data, std::vector<double>& local_Jac_data, LocalCoupledSolutions const& local_coupled_solutions) { - auto const& local_d = - local_coupled_solutions.local_coupled_xs[_phase_field_process_id]; - auto const& local_T = - local_coupled_solutions.local_coupled_xs[_heat_conduction_process_id]; - assert(local_T.size() == temperature_size); - assert(local_d.size() == phasefield_size); + assert(local_coupled_solutions.local_coupled_xs.size() == + phasefield_size + displacement_size + temperature_size); - auto d = Eigen::Map< - typename ShapeMatricesType::template VectorType<phasefield_size> const>( - local_d.data(), phasefield_size); - - auto T = Eigen::Map<typename ShapeMatricesType::template VectorType< - temperature_size> const>(local_T.data(), temperature_size); + auto const d = Eigen::Map<PhaseFieldVector const>( + &local_coupled_solutions.local_coupled_xs[phasefield_index], + phasefield_size); + auto const T = Eigen::Map<TemperatureVector const>( + &local_coupled_solutions.local_coupled_xs[temperature_index], + temperature_size); auto T_dot = Eigen::Map<typename ShapeMatricesType::template VectorType< temperature_size> const>(local_xdot.data(), temperature_size); @@ -276,13 +265,12 @@ void ThermoMechanicalPhaseFieldLocalAssembler<ShapeFunction, IntegrationMethod, std::vector<double>& local_b_data, std::vector<double>& local_Jac_data, LocalCoupledSolutions const& local_coupled_solutions) { - auto const& local_d = - local_coupled_solutions.local_coupled_xs[_phase_field_process_id]; - assert(local_d.size() == phasefield_size); + assert(local_coupled_solutions.local_coupled_xs.size() == + phasefield_size + displacement_size + temperature_size); - auto d = Eigen::Map< - typename ShapeMatricesType::template VectorType<phasefield_size> const>( - local_d.data(), phasefield_size); + auto const d = Eigen::Map<PhaseFieldVector const>( + &local_coupled_solutions.local_coupled_xs[phasefield_index], + phasefield_size); auto local_Jac = MathLib::createZeroedMatrix< typename ShapeMatricesType::template MatrixType<phasefield_size, diff --git a/ProcessLib/ThermoMechanicalPhaseField/ThermoMechanicalPhaseFieldFEM.h b/ProcessLib/ThermoMechanicalPhaseField/ThermoMechanicalPhaseFieldFEM.h index fff7bf8ecaf..e3dcf0e5bcf 100644 --- a/ProcessLib/ThermoMechanicalPhaseField/ThermoMechanicalPhaseFieldFEM.h +++ b/ProcessLib/ThermoMechanicalPhaseField/ThermoMechanicalPhaseFieldFEM.h @@ -111,6 +111,17 @@ template <typename ShapeFunction, typename IntegrationMethod, class ThermoMechanicalPhaseFieldLocalAssembler : public ThermoMechanicalPhaseFieldLocalAssemblerInterface { +private: + static constexpr int temperature_index = 0; + static constexpr int temperature_size = ShapeFunction::NPOINTS; + static constexpr int displacement_index = + temperature_index + temperature_size; + static constexpr int displacement_size = + ShapeFunction::NPOINTS * DisplacementDim; + static constexpr int phasefield_index = + displacement_index + displacement_size; + static constexpr int phasefield_size = ShapeFunction::NPOINTS; + public: using ShapeMatricesType = ShapeMatrixPolicyType<ShapeFunction, DisplacementDim>; @@ -123,6 +134,13 @@ public: using GlobalDimVectorType = typename ShapeMatricesType::GlobalDimVectorType; + using TemperatureVector = + typename ShapeMatricesType::template VectorType<temperature_size>; + using DeformationVector = + typename ShapeMatricesType::template VectorType<displacement_size>; + using PhaseFieldVector = + typename ShapeMatricesType::template VectorType<phasefield_size>; + ThermoMechanicalPhaseFieldLocalAssembler( ThermoMechanicalPhaseFieldLocalAssembler const&) = delete; ThermoMechanicalPhaseFieldLocalAssembler( @@ -371,14 +389,6 @@ private: SecondaryData<typename ShapeMatrices::ShapeType> _secondary_data; bool const _is_axially_symmetric; - static const int temperature_index = 0; - static const int temperature_size = ShapeFunction::NPOINTS; - static const int phasefield_index = ShapeFunction::NPOINTS; - static const int phasefield_size = ShapeFunction::NPOINTS; - static const int displacement_index = 2 * ShapeFunction::NPOINTS; - static const int displacement_size = - ShapeFunction::NPOINTS * DisplacementDim; - /// ID of the processes that contains mechanical process. int const _mechanics_related_process_id; diff --git a/ProcessLib/ThermoMechanics/ThermoMechanicsFEM-impl.h b/ProcessLib/ThermoMechanics/ThermoMechanicsFEM-impl.h index b18a550ecbd..7e51944ba75 100644 --- a/ProcessLib/ThermoMechanics/ThermoMechanicsFEM-impl.h +++ b/ProcessLib/ThermoMechanics/ThermoMechanicsFEM-impl.h @@ -338,13 +338,11 @@ void ThermoMechanicsLocalAssembler<ShapeFunction, IntegrationMethod, std::vector<double>& local_b_data, std::vector<double>& local_Jac_data, const LocalCoupledSolutions& local_coupled_solutions) { - auto const& local_T_vector = - local_coupled_solutions - .local_coupled_xs[_process_data.heat_conduction_process_id]; - assert(local_T_vector.size() == temperature_size); auto const local_T = Eigen::Map<typename ShapeMatricesType::template VectorType< - temperature_size> const>(local_T_vector.data(), temperature_size); + temperature_size> const>( + &local_coupled_solutions.local_coupled_xs[temperature_index], + temperature_size); auto const local_T0 = Eigen::Map<typename ShapeMatricesType::template VectorType< @@ -352,12 +350,10 @@ void ThermoMechanicsLocalAssembler<ShapeFunction, IntegrationMethod, &local_coupled_solutions.local_coupled_xs0[temperature_index], temperature_size); - auto const& local_u_vector = - local_coupled_solutions - .local_coupled_xs[_process_data.mechanics_process_id]; - assert(local_u_vector.size() == displacement_size); auto const u = Eigen::Map<typename ShapeMatricesType::template VectorType< - displacement_size> const>(local_u_vector.data(), displacement_size); + displacement_size> const>( + &local_coupled_solutions.local_coupled_xs[displacement_index], + displacement_size); auto local_Jac = MathLib::createZeroedMatrix< typename ShapeMatricesType::template MatrixType<displacement_size, @@ -476,19 +472,18 @@ void ThermoMechanicsLocalAssembler<ShapeFunction, IntegrationMethod, std::vector<double>& local_b_data, std::vector<double>& local_Jac_data, const LocalCoupledSolutions& local_coupled_solutions) { - auto const& local_T_vector = - local_coupled_solutions - .local_coupled_xs[_process_data.heat_conduction_process_id]; - assert(local_T_vector.size() == temperature_size); auto const local_T = Eigen::Map<typename ShapeMatricesType::template VectorType< - temperature_size> const>(local_T_vector.data(), temperature_size); + temperature_size> const>( + &local_coupled_solutions.local_coupled_xs[temperature_index], + temperature_size); auto const local_T0 = Eigen::Map<typename ShapeMatricesType::template VectorType< temperature_size> const>( &local_coupled_solutions.local_coupled_xs0[temperature_index], temperature_size); + auto const local_dT = local_T - local_T0; auto local_Jac = MathLib::createZeroedMatrix< -- GitLab