[Date Prev][Date Next][Thread Prev][Thread Next][Date Index][Thread Index]
[Getfem-commits] (no subject)
From: |
Tetsuo Koyama |
Subject: |
[Getfem-commits] (no subject) |
Date: |
Sun, 3 Nov 2019 03:01:21 -0500 (EST) |
branch: devel-tetsuo-hht_method
commit bf0c8e15947826c88d86641a698511f35c42a08c
Author: Tetsuo Koyama <address@hidden>
Date: Fri Nov 1 04:46:20 2019 +0900
ADD: Section "The implicit Hilber Hughes Taylor (HHT) scheme"
---
.../source/userdoc/model_time_integration.rst | 36 ++++++++++++++++++++++
interface/tests/python/demo_wave_equation.py | 2 +-
2 files changed, 37 insertions(+), 1 deletion(-)
diff --git a/doc/sphinx/source/userdoc/model_time_integration.rst
b/doc/sphinx/source/userdoc/model_time_integration.rst
index 20963a8..5180564 100644
--- a/doc/sphinx/source/userdoc/model_time_integration.rst
+++ b/doc/sphinx/source/userdoc/model_time_integration.rst
@@ -258,6 +258,42 @@ The addition of this scheme to a variable is to be done
thanks to::
add_Houbolt_scheme(model &md, const std::string &varname);
+The implicit Hilber Hughes Taylor (HHT) scheme
+**********************************************
+
+For a problem which reads
+
+.. math::
+
+ \left(\left(1+\alpha\right)K+\left(1+\alpha\right)\dfrac{\gamma}{\beta\Delta
t}C+\dfrac{1}{\beta\Delta
t^{2}}M\right)U^{n}=F^{n}+M\left(\left(\dfrac{1}{2\beta}-1\right)A^{n-1}+\dfrac{1}{\beta\Delta
t}V^{n-1}+\dfrac{1}{\beta\Delta
t^{2}}U^{n-1}\right)+\left(1+\alpha\right)C\left[\left(\dfrac{\gamma}{2\beta}-1\right)\Delta
tA^{n-1}+\left(\dfrac{\gamma}{\beta}-1\right)V^{n-1}+\dfrac{\gamma}{\beta\Delta
t}U^{n-1}\right]+\alpha KU^{n-1}+\alpha CV^{n-1}
+
+where :math:`dt` means a time step, :math:`M` the matrix in term of "Dot2_u",
:math:`C` the matrix in term of "Dot_u" and :math:`K` the matrix in term of "u".
+
+The affine dependences are thus given by::
+
+ Dot_u =
(1-gamma/beta)*Previous_Dot_u+(1-gamma/(2*beta))*dt*Previous_Dot2_u+gamma/(beta*dt)*(u-Previous_u)
+ Dot2_u =
(1-1/(2*beta))*Previous_Dot2_u-1/(beta*dt)*Previous_Dot_u+1/(beta*dt**2)*(u-Previous_u)
+
+When aplying this scheme to a variable "u" of the model, the following affine
dependent variables are added to the model::
+
+ "Dot_u", "Dot2_u"
+
+which represent the first and second order time derivative of the variable and
can be used in some brick definition.
+
+The following data are also added::
+
+ "Previous_u", "Previous_Dot_u", "Previous_Dot2_u"
+
+which correspond to the values of "u", "Dot_u" and "Dot2_u" at the previous
time step.
+
+Before the first solve, the data "Previous_u" and "Previous_Dot_u"
(corresponding to :math:`U^0` in the example) have to be initialized. The data
"Previous_Dot2_u" is to be given or precomputed (see
:ref:`precomp_time_der_section` and except for :math:`\alpha = 0, \beta = 1/2,
\gamma = 1`).
+
+
+The addition of this scheme to a variable is to be done thanks to::
+
+ add_Hilber_Hughes_Taylor_scheme(model &md, const std::string &varname,
+ scalar_type alpha, scalar_type beta,
scalar_type gamma);
+
Transient terms
***************
diff --git a/interface/tests/python/demo_wave_equation.py
b/interface/tests/python/demo_wave_equation.py
index 3069a9a..a56505e 100644
--- a/interface/tests/python/demo_wave_equation.py
+++ b/interface/tests/python/demo_wave_equation.py
@@ -75,7 +75,7 @@ gamma = 0.5;
md.add_Newmark_scheme('u', beta, gamma)
md.add_Houbolt_scheme('u1')
-md.add_Hilber_Hughes_Taylor_scheme('u2', 'F', alpha, beta, gamma)
+md.add_Hilber_Hughes_Taylor_scheme('u2', alpha, beta, gamma)
md.add_mass_brick(mim, 'Dot2_u')
md.add_mass_brick(mim, 'Dot2_u1')
md.set_time_step(dt)