FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
mechgauss.f90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
2 ! Copyright (c) 2019 FrontISTR Commons
3 ! This software is released under the MIT License, see LICENSE.txt
4 !-------------------------------------------------------------------------------
6 module mmechgauss
7  use hecmw_util
8  use mmaterial
9  implicit none
10 
11  ! ----------------------------------------------------------------------------
14  type(tmaterial), pointer :: pmaterial => null()
15  real(kind=kreal) :: strain(6)
16  real(kind=kreal) :: stress(6)
17  integer, pointer :: istatus(:) =>null()
18  real(kind=kreal), pointer :: fstatus(:) => null()
19  real(kind=kreal) :: plstrain
20  real(kind=kreal) :: strain_bak(6)
21  real(kind=kreal) :: stress_bak(6)
22  real(kind=kreal) :: nqm(12)
23  real(kind=kreal) :: strain_out(6)
24  real(kind=kreal) :: stress_out(6)
25  real(kind=kreal) :: strain_energy
26  real(kind=kreal) :: strain_energy_bak
27  real(kind=kreal) :: plpotential
28  end type
29 
30  ! ----------------------------------------------------------------------------
32  type telement
33  integer :: etype
34  integer :: iset
35  real(kind=kreal), pointer :: equiforces(:) => null()
36  type(tgaussstatus), pointer :: gausses(:) => null()
37  real(kind=kreal), pointer :: aux(:,:) => null()
38  end type
39 
40 contains
41 
43  subroutine fstr_init_gauss( gauss )
45  type( tgaussstatus ), intent(inout) :: gauss
46  integer :: n
47  gauss%strain=0.d0; gauss%stress=0.d0
48  gauss%strain_bak=0.d0; gauss%stress_bak=0.d0
49  gauss%strain_out=0.d0; gauss%stress_out=0.d0
50  gauss%plstrain =0.d0
51  gauss%nqm =0.d0
52  gauss%strain_energy =0.d0
53  gauss%strain_energy_bak =0.d0
54  if( gauss%pMaterial%mtype==usermaterial ) then
55  if( gauss%pMaterial%nfstatus> 0 ) then
56  allocate( gauss%fstatus(gauss%pMaterial%nfstatus) )
57  gauss%fstatus(:) = 0.d0
58  endif
59  else if( iselastoplastic(gauss%pMaterial%mtype) ) then
60  allocate( gauss%istatus(1) ) ! 0:elastic 1:plastic
61  if( getyieldfunction( gauss%pMaterial%mtype )==3 ) then ! user defined
62  n = uelastoplasticnumstatus( gauss%pMaterial%variables )
63  if( n>0 ) allocate( gauss%fstatus(n) )
64  elseif( iskinematicharden( gauss%pMaterial%mtype ) ) then
65  allocate( gauss%fstatus(7+6) ) ! plastic strain, back stress
66  else
67  allocate( gauss%fstatus(2) ) ! plastic strain
68  endif
69  gauss%istatus = 0
70  gauss%fstatus = 0.d0
71  else if( isviscoelastic(gauss%pMaterial%mtype) ) then
72  n = fetch_tablerow( mc_viscoelastic, gauss%pMaterial%dict )
73  if( n>0 ) then
74  allocate( gauss%fstatus(12*n+6) ) ! visco stress components
75  gauss%fstatus = 0.d0
76  else
77  stop "Viscoelastic properties not defined"
78  endif
79  else if( gauss%pMaterial%mtype==norton ) then
80  allocate( gauss%fstatus(2) ) ! effective stress, effective viscoplastic strain
81  gauss%fstatus = 0.d0
82  gauss%plstrain = 0.d0
83  endif
84  end subroutine fstr_init_gauss
85 
87  subroutine fstr_finalize_gauss( gauss )
88  type( tgaussstatus ), intent(inout) :: gauss
89  if( associated( gauss%istatus ) ) deallocate( gauss%istatus )
90  if( associated( gauss%fstatus ) ) deallocate( gauss%fstatus )
91  end subroutine
92 
94  subroutine fstr_copy_gauss( gauss1, gauss2 )
95  type( tgaussstatus ), intent(in) :: gauss1
96  type( tgaussstatus ), intent(inout) :: gauss2
97 
98  gauss2%strain = gauss1%strain
99  gauss2%stress = gauss1%stress
100  gauss2%strain_bak = gauss1%strain_bak
101  gauss2%stress_bak = gauss1%stress_bak
102  gauss2%plstrain = gauss1%plstrain
103  gauss2%strain_energy = gauss1%strain_energy
104 
105  if( associated(gauss1%istatus) .and. associated(gauss2%istatus) ) then
106  gauss2%istatus = gauss1%istatus
107  end if
108  if( associated(gauss1%fstatus) .and. associated(gauss2%fstatus) ) then
109  gauss2%fstatus = gauss1%fstatus
110  end if
111  end subroutine fstr_copy_gauss
112 
113 
114 end module
115 
116 
117 
mmaterial::mc_viscoelastic
character(len=dict_key_length) mc_viscoelastic
Definition: material.f90:145
mmaterial::iskinematicharden
logical function iskinematicharden(mtype)
If it is a kinematic hardening material?
Definition: material.f90:343
muyield
This subroutine read in used-defined material properties tangent.
Definition: uyield.f90:7
mmechgauss::fstr_finalize_gauss
subroutine fstr_finalize_gauss(gauss)
Finializer.
Definition: mechgauss.f90:88
mmechgauss::fstr_copy_gauss
subroutine fstr_copy_gauss(gauss1, gauss2)
Copy.
Definition: mechgauss.f90:95
mmaterial::usermaterial
integer(kind=kint), parameter usermaterial
Definition: material.f90:63
mmaterial::iselastoplastic
logical function iselastoplastic(mtype)
If it is an elastoplastic material?
Definition: material.f90:361
mmechgauss
This modules defines a structure to record history dependent parameter in static analysis.
Definition: mechgauss.f90:6
hecmw_util
I/O and Utility.
Definition: hecmw_util_f.F90:7
muyield::uelastoplasticnumstatus
integer(kind=kint) function, public uelastoplasticnumstatus(matl)
This function returns the number of real state variables.
Definition: uyield.f90:25
hecmw_util::kreal
integer(kind=4), parameter kreal
Definition: hecmw_util_f.F90:16
mmechgauss::fstr_init_gauss
subroutine fstr_init_gauss(gauss)
Initializer.
Definition: mechgauss.f90:44
mmaterial::getyieldfunction
integer function getyieldfunction(mtype)
Get type of yield function.
Definition: material.f90:319
mmechgauss::tgaussstatus
All data should be recorded in every quadrature points.
Definition: mechgauss.f90:13
mmaterial::tmaterial
Structure to manage all material related data.
Definition: material.f90:167
mmaterial::isviscoelastic
logical function isviscoelastic(mtype)
If it is an viscoelastic material?
Definition: material.f90:379
mmaterial
This module summarizes all information of material properties.
Definition: material.f90:6
mmaterial::norton
integer(kind=kint), parameter norton
Definition: material.f90:78
mmechgauss::telement
All data should be recorded in every elements.
Definition: mechgauss.f90:32