From 314208810c5ea04bf43c5e55ccee45eccd009e43 Mon Sep 17 00:00:00 2001 From: Wenqing Wang <wenqing.wang@ufz.de> Date: Wed, 10 Jan 2018 14:12:31 +0100 Subject: [PATCH] [HM] Documented on new members of HM process. --- ProcessLib/HydroMechanics/HydroMechanicsFEM.h | 53 +++++++++++++++++++ .../HydroMechanicsProcess-impl.h | 19 +++---- .../HydroMechanics/HydroMechanicsProcess.h | 17 +++++- 3 files changed, 79 insertions(+), 10 deletions(-) diff --git a/ProcessLib/HydroMechanics/HydroMechanicsFEM.h b/ProcessLib/HydroMechanics/HydroMechanicsFEM.h index e8f318c01ff..df32f26f290 100644 --- a/ProcessLib/HydroMechanics/HydroMechanicsFEM.h +++ b/ProcessLib/HydroMechanics/HydroMechanicsFEM.h @@ -326,6 +326,31 @@ private: return cache; } + /** + * Assemble local matrices and vectors arise from the linearized discretized + * weak form of the residual of the momentum balance equation, + * \f[ + * \nabla (\sigma - \alpha_b p \mathrm{I}) = f + * \f] + * where \f$ \sigma\f$ is the effective stress tensor, \f$p\f$ is the pore + * pressure, \f$\alpha_b\f$ is the Biots constant, \f$\mathrm{I}\f$ is the + * identity tensor, and \f$f\f$ is the body force. + * + * @param t Time + * @param local_xdot Nodal values of \f$\dot{x}\f$ of an element. + * @param dxdot_dx Value of \f$\dot{x} \cdot dx\f$. + * @param dx_dx Value of \f$ x \cdot dx\f$. + * @param local_M_data Mass matrix of an element, which takes the form of + * \f$ \inta N^T N\mathrm{d}\Omega\f$. Not used. + * @param local_K_data Lappacian matrix of an element, which takes the + * form of \f$ \int (\nabla N)^T K \nabla N\mathrm{d}\Omega\f$. + * Not used. + * @param local_b_data Right hand side vector of an element. + * @param local_Jac_data Element Jacobian matrix from the Newton-Raphson + * method. + * @param local_coupled_solutions Nodal values of solutions of the coupled + * processes of an element. + */ void assembleWithJacobianForDeformationEquations( double const t, std::vector<double> const& local_xdot, const double dxdot_dx, const double dx_dx, @@ -333,6 +358,34 @@ private: std::vector<double>& local_b_data, std::vector<double>& local_Jac_data, LocalCoupledSolutions const& local_coupled_solutions); + /** + * Assemble local matrices and vectors arise from the linearized discretized + * weak form of the residual of the mass balance equation of single phase + * flow, + * \f[ + * \alpha \cdot{p} - \nabla (K (\nabla p + \rho g \nabla z) + + * \alpha_b \nabla \cdot \dot{u} = Q + * \f] + * where \f$ alpha\f$ is a coefficient may depend on storage or the fluid + * density change, \f$ \rho\f$ is the fluid density, \f$\g\f$ is the + * gravitational acceleration, \f$z\f$ is the vertical coordinate, \f$u\f$ + * is the displacement, and \f$Q\f$ is the source/sink term. + * + * @param t Time + * @param local_xdot Nodal values of \f$\dot{x}\f$ of an element. + * @param dxdot_dx Value of \f$\dot{x} \cdot dx\f$. + * @param dx_dx Value of \f$ x \cdot dx\f$. + * @param local_M_data Mass matrix of an element, which takes the form of + * \f$ \inta N^T N\mathrm{d}\Omega\f$. Not used. + * @param local_K_data Lappacian matrix of an element, which takes the + * form of \f$ \int (\nabla N)^T K \nabla N\mathrm{d}\Omega\f$. + * Not used. + * @param local_b_data Right hand side vector of an element. + * @param local_Jac_data Element Jacobian matrix from the Newton-Raphson + * method. + * @param local_coupled_solutions Nodal values of solutions of the coupled + * processes of an element. + */ void assembleWithJacobianForPressureEquations( double const t, std::vector<double> const& local_xdot, const double dxdot_dx, const double dx_dx, diff --git a/ProcessLib/HydroMechanics/HydroMechanicsProcess-impl.h b/ProcessLib/HydroMechanics/HydroMechanicsProcess-impl.h index 87a07436760..3c00ac7124b 100644 --- a/ProcessLib/HydroMechanics/HydroMechanicsProcess-impl.h +++ b/ProcessLib/HydroMechanics/HydroMechanicsProcess-impl.h @@ -53,19 +53,18 @@ bool HydroMechanicsProcess<DisplacementDim>::isLinear() const template <int DisplacementDim> MathLib::MatrixSpecifications HydroMechanicsProcess<DisplacementDim>::getMatrixSpecifications( - const int equation_id) const + const int process_id) const { - // For the staggered scheme, equation_id == 1 indicates that the matrix - // specifications are for the momentum balance equation (deformation). - if (_use_monolithic_scheme || equation_id == 1) + // For the monolithic scheme or the M process (deformation) in the staggered + // scheme. + if (_use_monolithic_scheme || process_id == 1) { auto const& l = *_local_to_global_index_map; return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(), &l.getGhostIndices(), &this->_sparsity_pattern}; } - // For staggered scheme and the mass conservation balance equation - // (pressure). + // For staggered scheme and H process (pressure). auto const& l = *_local_to_global_index_map_with_base_nodes; return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(), &l.getGhostIndices(), &_sparsity_pattern_with_linear_element}; @@ -158,7 +157,7 @@ void HydroMechanicsProcess<DisplacementDim>::initializeConcreteProcess( ProcessLib::HydroMechanics::createLocalAssemblers< DisplacementDim, HydroMechanicsLocalAssembler>( mesh.getDimension(), mesh.getElements(), dof_table, - // use displacment process variable for shapefunction order + // use displacement process variable to set shape function order getProcessVariables(mechinical_process_id)[deformation_variable_id] .get() .getShapeFunctionOrder(), @@ -253,10 +252,12 @@ void HydroMechanicsProcess<DisplacementDim>::assembleConcreteProcess( const double t, GlobalVector const& x, GlobalMatrix& M, GlobalMatrix& K, GlobalVector& b) { - // Not available because that only the Newton-Raphson method is available - // for HydroMechanics DBUG("Assemble the equations for HydroMechanics"); + // Note: This assembly function is for the Picard nonlinear solver. Since + // only the Newton-Raphson method is employed to simulate coupled HM + // processes in this class, this function is actually not used so far. + // Call global assembler for each local assembly item. GlobalExecutor::executeMemberDereferenced( _global_assembler, &VectorMatrixAssembler::assemble, _local_assemblers, diff --git a/ProcessLib/HydroMechanics/HydroMechanicsProcess.h b/ProcessLib/HydroMechanics/HydroMechanicsProcess.h index 716a1290c00..97e2593c38e 100644 --- a/ProcessLib/HydroMechanics/HydroMechanicsProcess.h +++ b/ProcessLib/HydroMechanics/HydroMechanicsProcess.h @@ -48,8 +48,20 @@ public: bool isLinear() const override; //! @} + /** + * Get the size and the sparse pattern of the global matrix in order to + * create the global matrices and vectors for the system equations of this + * process. + * + * @param process_id Process ID. If the monolithic scheme is applied, + * process_id = 0. For the staggered scheme, + * process_id = 0 represents the + * hydraulic (H) process, while process_id = 1 + * represents the mechanical (M) process. + * @return Matrix specifications including size and sparse pattern. + */ MathLib::MatrixSpecifications getMatrixSpecifications( - const int equation_id) const override; + const int process_id) const override; private: void constructDofTable() override; @@ -105,6 +117,9 @@ private: /// Solutions of the previous time step std::array<std::unique_ptr<GlobalVector>, 2> _xs_previous_timestep; + /** + * @copydoc ProcessLib::Process::getDOFTableForExtrapolatorData() + */ std::tuple<NumLib::LocalToGlobalIndexMap*, bool> getDOFTableForExtrapolatorData() const override; -- GitLab