![]() |
FrontISTR
5.7.0
Large-scale structural analysis program with finit element method
|
This module provides functions on nonlinear analysis. More...
Functions/Subroutines | |
subroutine | fstr_quasi_newton (cstep, hecMESH, hecMAT, fstrSOLID, fstrPARAM, restrt_step_num, sub_step, ctime, dtime) |
This subroutine solve nonlinear solid mechanics problems by Newton-Raphson method. More... | |
subroutine | fstr_calc_direction_lbfgs (hecMESH, g_prev, s_k, y_k, rho_k, z_k, n_mem) |
subroutine | fstr_addbc_to_direction_vector (z_k, hecMESH, fstrSOLID, cstep) |
subroutine | fstr_apply_alpha0 (hecMESH, hecMAT, fstrSOLID, ctime, tincr, iter, cstep, dtime, fstrPARAM, z_k, h_prime, pot) |
subroutine | fstr_apply_alpha (hecMESH, hecMAT, fstrSOLID, ctime, tincr, iter, cstep, dtime, fstrPARAM, z_k, alpha, h_prime, pot) |
subroutine | fstr_init_line_search_range (hecMESH, hecMAT, fstrSOLID, ctime, tincr, iter, cstep, dtime, fstrPARAM, z_k, h_prime_0, pot_0, alpha_S, h_prime_S, pot_S, alpha_E, h_prime_E, pot_E) |
subroutine | fstr_update_range (hecMESH, hecMAT, fstrSOLID, ctime, tincr, iter, cstep, dtime, fstrPARAM, z_k, alpha_S, h_prime_S, pot_S, alpha_E, h_prime_E, pot_E, alpha_c, h_prime_c, pot_c, h_prime_0, pot_0) |
subroutine | fstr_get_new_range_with_potential (hecMESH, hecMAT, fstrSOLID, ctime, tincr, iter, cstep, dtime, fstrPARAM, z_k, pot_0, alpha_a, h_prime_a, pot_a, alpha_b, h_prime_b, pot_b, alpha_c, h_prime_c, pot_c, alpha_a_bar, h_prime_a_bar, pot_a_bar, alpha_b_bar, h_prime_b_bar, pot_b_bar) |
subroutine | fstr_line_search_along_direction (hecMESH, hecMAT, fstrSOLID, ctime, tincr, iter, cstep, dtime, fstrPARAM, z_k) |
subroutine | fstr_get_secant (a, Fa, b, Fb, c) |
logical function | fstr_wolfe_condition (alpha_c, h_prime_c, pot_c, h_prime_0, pot_0) |
logical function | fstr_approx_wolfe_condition (alpha_c, h_prime_c, pot_c, h_prime_0, pot_0) |
Variables | |
real(kind=kreal), parameter | c_line_search =2.0 |
real(kind=kreal), parameter | psi0_line_search =1.0d-2 |
real(kind=kreal), parameter | delta_wolfe =0.2 |
real(kind=kreal), parameter | sigma_wolfe =0.9 |
real(kind=kreal), parameter | eps_wolfe =1.0d-3 |
real(kind=kreal), parameter | omega_wolfe =0.001 |
real(kind=kreal), parameter | delta_approx_wolfe =0.7 |
real(kind=kreal) | c_wolfe |
real(kind=kreal) | q_wolfe |
integer, parameter | n_mem_max =10 |
integer(kind=kint), parameter | maxiter_ls = 10 |
integer(kind=kint), parameter | pot_type =1 |
This module provides functions on nonlinear analysis.
subroutine m_fstr_quasinewton::fstr_addbc_to_direction_vector | ( | real(kind=kreal), dimension(:) | z_k, |
type (hecmwst_local_mesh) | hecMESH, | ||
type (fstr_solid) | fstrSOLID, | ||
integer(kind=kint) | cstep | ||
) |
hecmesh | hecmw mesh |
fstrsolid | fstr_solid |
Definition at line 231 of file fstr_solve_QuasiNewton.f90.
subroutine m_fstr_quasinewton::fstr_apply_alpha | ( | type (hecmwst_local_mesh) | hecMESH, |
type (hecmwst_matrix) | hecMAT, | ||
type (fstr_solid) | fstrSOLID, | ||
real(kind=kreal), intent(in) | ctime, | ||
real(kind=kreal), intent(in) | tincr, | ||
integer(kind=kint) | iter, | ||
integer, intent(in) | cstep, | ||
real(kind=kreal), intent(in) | dtime, | ||
type (fstr_param) | fstrPARAM, | ||
real(kind=kreal), dimension(:), intent(in) | z_k, | ||
real(kind=kreal), intent(in) | alpha, | ||
real(kind=kreal) | h_prime, | ||
real(kind=kreal) | pot | ||
) |
hecmesh | hecmw mesh | |
hecmat | hecmw matrix | |
fstrsolid | fstr_solid | |
[in] | ctime | current time |
[in] | cstep | current loading step |
[in] | dtime | time increment |
fstrparam | type fstr_param |
Definition at line 282 of file fstr_solve_QuasiNewton.f90.
subroutine m_fstr_quasinewton::fstr_apply_alpha0 | ( | type (hecmwst_local_mesh) | hecMESH, |
type (hecmwst_matrix) | hecMAT, | ||
type (fstr_solid) | fstrSOLID, | ||
real(kind=kreal), intent(in) | ctime, | ||
real(kind=kreal), intent(in) | tincr, | ||
integer(kind=kint) | iter, | ||
integer, intent(in) | cstep, | ||
real(kind=kreal), intent(in) | dtime, | ||
type (fstr_param) | fstrPARAM, | ||
real(kind=kreal), dimension(:), intent(in) | z_k, | ||
real(kind=kreal) | h_prime, | ||
real(kind=kreal) | pot | ||
) |
hecmesh | hecmw mesh | |
hecmat | hecmw matrix | |
fstrsolid | fstr_solid | |
[in] | ctime | current time |
[in] | cstep | current loading step |
[in] | dtime | time increment |
fstrparam | type fstr_param |
Definition at line 263 of file fstr_solve_QuasiNewton.f90.
logical function m_fstr_quasinewton::fstr_approx_wolfe_condition | ( | real(kind=kreal), intent(in) | alpha_c, |
real(kind=kreal), intent(in) | h_prime_c, | ||
real(kind=kreal), intent(in) | pot_c, | ||
real(kind=kreal), intent(in) | h_prime_0, | ||
real(kind=kreal), intent(in) | pot_0 | ||
) |
Definition at line 635 of file fstr_solve_QuasiNewton.f90.
subroutine m_fstr_quasinewton::fstr_calc_direction_lbfgs | ( | type (hecmwst_local_mesh) | hecMESH, |
real(kind=kreal), dimension(:) | g_prev, | ||
real(kind=kreal), dimension(:,:) | s_k, | ||
real(kind=kreal), dimension(:,:) | y_k, | ||
real(kind=kreal), dimension(:) | rho_k, | ||
real(kind=kreal), dimension(:) | z_k, | ||
integer | n_mem | ||
) |
hecmesh | hecmw mesh |
Definition at line 175 of file fstr_solve_QuasiNewton.f90.
subroutine m_fstr_quasinewton::fstr_get_new_range_with_potential | ( | type (hecmwst_local_mesh) | hecMESH, |
type (hecmwst_matrix) | hecMAT, | ||
type (fstr_solid) | fstrSOLID, | ||
real(kind=kreal), intent(in) | ctime, | ||
real(kind=kreal), intent(in) | tincr, | ||
integer(kind=kint) | iter, | ||
integer, intent(in) | cstep, | ||
real(kind=kreal), intent(in) | dtime, | ||
type (fstr_param) | fstrPARAM, | ||
real(kind=kreal), dimension(:), intent(in) | z_k, | ||
real(kind=kreal), intent(in) | pot_0, | ||
real(kind=kreal), intent(in) | alpha_a, | ||
real(kind=kreal), intent(in) | h_prime_a, | ||
real(kind=kreal), intent(in) | pot_a, | ||
real(kind=kreal), intent(in) | alpha_b, | ||
real(kind=kreal), intent(in) | h_prime_b, | ||
real(kind=kreal), intent(in) | pot_b, | ||
real(kind=kreal), intent(in) | alpha_c, | ||
real(kind=kreal), intent(in) | h_prime_c, | ||
real(kind=kreal), intent(in) | pot_c, | ||
real(kind=kreal) | alpha_a_bar, | ||
real(kind=kreal) | h_prime_a_bar, | ||
real(kind=kreal) | pot_a_bar, | ||
real(kind=kreal) | alpha_b_bar, | ||
real(kind=kreal) | h_prime_b_bar, | ||
real(kind=kreal) | pot_b_bar | ||
) |
hecmesh | hecmw mesh | |
hecmat | hecmw matrix | |
fstrsolid | fstr_solid | |
[in] | ctime | current time |
[in] | cstep | current loading step |
[in] | dtime | time increment |
fstrparam | type fstr_param |
Definition at line 442 of file fstr_solve_QuasiNewton.f90.
subroutine m_fstr_quasinewton::fstr_get_secant | ( | real(kind=kreal), intent(in) | a, |
real(kind=kreal), intent(in) | Fa, | ||
real(kind=kreal), intent(in) | b, | ||
real(kind=kreal), intent(in) | Fb, | ||
real(kind=kreal), intent(out) | c | ||
) |
Definition at line 607 of file fstr_solve_QuasiNewton.f90.
subroutine m_fstr_quasinewton::fstr_init_line_search_range | ( | type (hecmwst_local_mesh) | hecMESH, |
type (hecmwst_matrix) | hecMAT, | ||
type (fstr_solid) | fstrSOLID, | ||
real(kind=kreal), intent(in) | ctime, | ||
real(kind=kreal), intent(in) | tincr, | ||
integer(kind=kint) | iter, | ||
integer, intent(in) | cstep, | ||
real(kind=kreal), intent(in) | dtime, | ||
type (fstr_param) | fstrPARAM, | ||
real(kind=kreal), dimension(:), intent(in) | z_k, | ||
real(kind=kreal), intent(in) | h_prime_0, | ||
real(kind=kreal), intent(in) | pot_0, | ||
real(kind=kreal) | alpha_S, | ||
real(kind=kreal) | h_prime_S, | ||
real(kind=kreal) | pot_S, | ||
real(kind=kreal) | alpha_E, | ||
real(kind=kreal) | h_prime_E, | ||
real(kind=kreal) | pot_E | ||
) |
hecmesh | hecmw mesh | |
hecmat | hecmw matrix | |
fstrsolid | fstr_solid | |
[in] | ctime | current time |
[in] | cstep | current loading step |
[in] | dtime | time increment |
fstrparam | type fstr_param |
Definition at line 304 of file fstr_solve_QuasiNewton.f90.
subroutine m_fstr_quasinewton::fstr_line_search_along_direction | ( | type (hecmwst_local_mesh) | hecMESH, |
type (hecmwst_matrix) | hecMAT, | ||
type (fstr_solid) | fstrSOLID, | ||
real(kind=kreal), intent(in) | ctime, | ||
real(kind=kreal), intent(in) | tincr, | ||
integer(kind=kint) | iter, | ||
integer, intent(in) | cstep, | ||
real(kind=kreal), intent(in) | dtime, | ||
type (fstr_param) | fstrPARAM, | ||
real(kind=kreal), dimension(:) | z_k | ||
) |
hecmesh | hecmw mesh | |
hecmat | hecmw matrix | |
fstrsolid | fstr_solid | |
[in] | ctime | current time |
[in] | cstep | current loading step |
[in] | dtime | time increment |
fstrparam | type fstr_param |
Definition at line 539 of file fstr_solve_QuasiNewton.f90.
subroutine m_fstr_quasinewton::fstr_quasi_newton | ( | integer, intent(in) | cstep, |
type (hecmwst_local_mesh) | hecMESH, | ||
type (hecmwst_matrix) | hecMAT, | ||
type (fstr_solid) | fstrSOLID, | ||
type (fstr_param) | fstrPARAM, | ||
integer(kind=kint) | restrt_step_num, | ||
integer, intent(in) | sub_step, | ||
real(kind=kreal), intent(in) | ctime, | ||
real(kind=kreal), intent(in) | dtime | ||
) |
This subroutine solve nonlinear solid mechanics problems by Newton-Raphson method.
[in] | cstep | current loading step |
hecmesh | hecmw mesh | |
hecmat | hecmw matrix | |
fstrsolid | fstr_solid | |
[in] | sub_step | substep number of current loading step |
[in] | ctime | current time |
[in] | dtime | time increment |
fstrparam | type fstr_param |
Definition at line 27 of file fstr_solve_QuasiNewton.f90.
subroutine m_fstr_quasinewton::fstr_update_range | ( | type (hecmwst_local_mesh) | hecMESH, |
type (hecmwst_matrix) | hecMAT, | ||
type (fstr_solid) | fstrSOLID, | ||
real(kind=kreal), intent(in) | ctime, | ||
real(kind=kreal), intent(in) | tincr, | ||
integer(kind=kint) | iter, | ||
integer, intent(in) | cstep, | ||
real(kind=kreal), intent(in) | dtime, | ||
type (fstr_param) | fstrPARAM, | ||
real(kind=kreal), dimension(:) | z_k, | ||
real(kind=kreal) | alpha_S, | ||
real(kind=kreal) | h_prime_S, | ||
real(kind=kreal) | pot_S, | ||
real(kind=kreal) | alpha_E, | ||
real(kind=kreal) | h_prime_E, | ||
real(kind=kreal) | pot_E, | ||
real(kind=kreal) | alpha_c, | ||
real(kind=kreal) | h_prime_c, | ||
real(kind=kreal) | pot_c, | ||
real(kind=kreal) | h_prime_0, | ||
real(kind=kreal) | pot_0 | ||
) |
hecmesh | hecmw mesh | |
hecmat | hecmw matrix | |
fstrsolid | fstr_solid | |
[in] | ctime | current time |
[in] | cstep | current loading step |
[in] | dtime | time increment |
fstrparam | type fstr_param |
Definition at line 369 of file fstr_solve_QuasiNewton.f90.
logical function m_fstr_quasinewton::fstr_wolfe_condition | ( | real(kind=kreal), intent(in) | alpha_c, |
real(kind=kreal), intent(in) | h_prime_c, | ||
real(kind=kreal), intent(in) | pot_c, | ||
real(kind=kreal), intent(in) | h_prime_0, | ||
real(kind=kreal), intent(in) | pot_0 | ||
) |
Definition at line 618 of file fstr_solve_QuasiNewton.f90.
real(kind=kreal), parameter m_fstr_quasinewton::c_line_search =2.0 |
Definition at line 13 of file fstr_solve_QuasiNewton.f90.
real(kind=kreal) m_fstr_quasinewton::c_wolfe |
Definition at line 16 of file fstr_solve_QuasiNewton.f90.
real(kind=kreal), parameter m_fstr_quasinewton::delta_approx_wolfe =0.7 |
Definition at line 15 of file fstr_solve_QuasiNewton.f90.
real(kind=kreal), parameter m_fstr_quasinewton::delta_wolfe =0.2 |
Definition at line 14 of file fstr_solve_QuasiNewton.f90.
real(kind=kreal), parameter m_fstr_quasinewton::eps_wolfe =1.0d-3 |
Definition at line 14 of file fstr_solve_QuasiNewton.f90.
integer(kind=kint), parameter m_fstr_quasinewton::maxiter_ls = 10 |
Definition at line 18 of file fstr_solve_QuasiNewton.f90.
integer, parameter m_fstr_quasinewton::n_mem_max =10 |
Definition at line 17 of file fstr_solve_QuasiNewton.f90.
real(kind=kreal), parameter m_fstr_quasinewton::omega_wolfe =0.001 |
Definition at line 15 of file fstr_solve_QuasiNewton.f90.
integer(kind=kint), parameter m_fstr_quasinewton::pot_type =1 |
Definition at line 19 of file fstr_solve_QuasiNewton.f90.
real(kind=kreal), parameter m_fstr_quasinewton::psi0_line_search =1.0d-2 |
Definition at line 13 of file fstr_solve_QuasiNewton.f90.
real(kind=kreal) m_fstr_quasinewton::q_wolfe |
Definition at line 16 of file fstr_solve_QuasiNewton.f90.
real(kind=kreal), parameter m_fstr_quasinewton::sigma_wolfe =0.9 |
Definition at line 14 of file fstr_solve_QuasiNewton.f90.