FrontISTR  5.8.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  integer :: elemact_flag
39  real(kind=kreal) :: elemact_coeff
40  end type
41 
42 contains
43 
45  subroutine fstr_init_gauss( gauss )
47  type( tgaussstatus ), intent(inout) :: gauss
48  integer :: n
49  gauss%strain=0.d0; gauss%stress=0.d0
50  gauss%strain_bak=0.d0; gauss%stress_bak=0.d0
51  gauss%strain_out=0.d0; gauss%stress_out=0.d0
52  gauss%plstrain =0.d0
53  gauss%nqm =0.d0
54  gauss%strain_energy =0.d0
55  gauss%strain_energy_bak =0.d0
56  if( gauss%pMaterial%mtype==usermaterial ) then
57  if( gauss%pMaterial%nfstatus> 0 ) then
58  allocate( gauss%fstatus(gauss%pMaterial%nfstatus) )
59  gauss%fstatus(:) = 0.d0
60  endif
61  else if( iselastoplastic(gauss%pMaterial%mtype) ) then
62  allocate( gauss%istatus(1) ) ! 0:elastic 1:plastic
63  if( getyieldfunction( gauss%pMaterial%mtype )==3 ) then ! user defined
64  n = uelastoplasticnumstatus( gauss%pMaterial%variables )
65  if( n>0 ) allocate( gauss%fstatus(n) )
66  elseif( iskinematicharden( gauss%pMaterial%mtype ) ) then
67  allocate( gauss%fstatus(7+6) ) ! plastic strain, back stress
68  else
69  allocate( gauss%fstatus(2) ) ! plastic strain
70  endif
71  gauss%istatus = 0
72  gauss%fstatus = 0.d0
73  else if( isviscoelastic(gauss%pMaterial%mtype) ) then
74  n = fetch_tablerow( mc_viscoelastic, gauss%pMaterial%dict )
75  if( n>0 ) then
76  allocate( gauss%fstatus(12*n+6) ) ! visco stress components
77  gauss%fstatus = 0.d0
78  else
79  stop "Viscoelastic properties not defined"
80  endif
81  else if( gauss%pMaterial%mtype==norton ) then
82  allocate( gauss%fstatus(2) ) ! effective stress, effective viscoplastic strain
83  gauss%fstatus = 0.d0
84  gauss%plstrain = 0.d0
85  endif
86  end subroutine fstr_init_gauss
87 
89  subroutine fstr_finalize_gauss( gauss )
90  type( tgaussstatus ), intent(inout) :: gauss
91  if( associated( gauss%istatus ) ) deallocate( gauss%istatus )
92  if( associated( gauss%fstatus ) ) deallocate( gauss%fstatus )
93  end subroutine
94 
96  subroutine fstr_copy_gauss( gauss1, gauss2 )
97  type( tgaussstatus ), intent(in) :: gauss1
98  type( tgaussstatus ), intent(inout) :: gauss2
99 
100  gauss2%strain = gauss1%strain
101  gauss2%stress = gauss1%stress
102  gauss2%strain_bak = gauss1%strain_bak
103  gauss2%stress_bak = gauss1%stress_bak
104  gauss2%plstrain = gauss1%plstrain
105  gauss2%strain_energy = gauss1%strain_energy
106 
107  if( associated(gauss1%istatus) .and. associated(gauss2%istatus) ) then
108  gauss2%istatus = gauss1%istatus
109  end if
110  if( associated(gauss1%fstatus) .and. associated(gauss2%fstatus) ) then
111  gauss2%fstatus = gauss1%fstatus
112  end if
113  end subroutine fstr_copy_gauss
114 
115 
116 end module
117 
118 
119 
I/O and Utility.
Definition: hecmw_util_f.F90:7
integer(kind=4), parameter kreal
This module summarizes all information of material properties.
Definition: material.f90:6
integer function getyieldfunction(mtype)
Get type of yield function.
Definition: material.f90:319
character(len=dict_key_length) mc_viscoelastic
Definition: material.f90:145
integer(kind=kint), parameter norton
Definition: material.f90:78
logical function iskinematicharden(mtype)
If it is a kinematic hardening material?
Definition: material.f90:343
integer(kind=kint), parameter usermaterial
Definition: material.f90:63
logical function isviscoelastic(mtype)
If it is an viscoelastic material?
Definition: material.f90:379
logical function iselastoplastic(mtype)
If it is an elastoplastic material?
Definition: material.f90:361
This modules defines a structure to record history dependent parameter in static analysis.
Definition: mechgauss.f90:6
subroutine fstr_init_gauss(gauss)
Initializer.
Definition: mechgauss.f90:46
subroutine fstr_finalize_gauss(gauss)
Finializer.
Definition: mechgauss.f90:90
subroutine fstr_copy_gauss(gauss1, gauss2)
Copy.
Definition: mechgauss.f90:97
This subroutine read in used-defined material properties tangent.
Definition: uyield.f90:7
integer(kind=kint) function, public uelastoplasticnumstatus(matl)
This function returns the number of real state variables.
Definition: uyield.f90:25
Structure to manage all material related data.
Definition: material.f90:167
All data should be recorded in every elements.
Definition: mechgauss.f90:32
All data should be recorded in every quadrature points.
Definition: mechgauss.f90:13