FrontISTR  5.8.0
Large-scale structural analysis program with finit element method
m_element_activation.f90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
2 ! Copyright (c) 2016 The University of Tokyo
3 ! This software is released under the MIT License, see LICENSE.txt
4 !-------------------------------------------------------------------------------
6 module m_elemact
7  use hecmw
8  use mmechgauss
9 
10  implicit none
11 
12  type telemact
14  integer(kind=kint) :: elemact_egrp_tot
15  integer(kind=kint), pointer :: elemact_egrp_grpid (:) =>null()
16  integer(kind=kint), pointer :: elemact_egrp_id (:) =>null()
17  integer(kind=kint), pointer :: elemact_egrp_amp (:) =>null()
18  real(kind=kreal), pointer :: elemact_egrp_eps(:) =>null()
19  integer(kind=kint), pointer :: elemact_egrp_depends(:) =>null()
20  real(kind=kreal), pointer :: elemact_egrp_ts_lower(:) =>null()
21  real(kind=kreal), pointer :: elemact_egrp_ts_upper(:) =>null()
22  integer(kind=kint), pointer :: elemact_egrp_state (:) =>null() ! Element activation state
23  end type
24 
25  ! Element status flags
26  integer, parameter :: kelact_undefined = -1 ! Undefined state
27  integer, parameter :: kelact_active = 0 ! Element is active in calculation
28  integer, parameter :: kelact_inactive = 1 ! Element is inactive (deactivated) in calculation
29 
30  ! Dependency type flags
31  integer, parameter :: kelactd_none = 1 ! No dependencies
32  integer, parameter :: kelactd_stress = 2 ! Stress dependent
33  integer, parameter :: kelactd_strain = 3 ! Strain dependent
34 
35 contains
36 
37 
38  subroutine print_elemact_info( elact )
39  type(telemact), intent(in) :: elact
40 
41  integer(kind=kint) :: i
42 
43  write(*,'(A,I0)') 'ELEMACT_egrp_tot: ', elact%ELEMACT_egrp_tot
44  do i=1,elact%ELEMACT_egrp_tot
45  write(*,*) 'DUMMY : ',i
46  write(*,*) 'GRPID : ',elact%ELEMACT_egrp_GRPID
47  write(*,*) 'EGRPID : ',elact%ELEMACT_egrp_ID
48  write(*,*) 'AMPLITUDE: ',elact%ELEMACT_egrp_amp
49  write(*,*) 'EPSILON : ',elact%ELEMACT_egrp_eps
50  write(*,*) 'DEPENDS : ',elact%ELEMACT_egrp_depends
51  write(*,*) 'THLOWER : ',elact%ELEMACT_egrp_ts_lower
52  write(*,*) 'THUPPER : ',elact%ELEMACT_egrp_ts_upper
53  end do
54  end subroutine
55 
56 
57  subroutine stf_dummy( ndof, nn, ecoord, u, stiff, element )
58  use mmechgauss
59  integer(kind=kint), intent(in) :: ndof
60  integer(kind=kint), intent(in) :: nn
61  real(kind=kreal), intent(in) :: ecoord(3,nn)
62  real(kind=kreal), intent(in) :: u(3,nn)
63  type(telement), intent(inout) :: element
64  real(kind=kreal), intent(out) :: stiff(:,:)
65 
66  integer(kind=kint) :: i,j,k,m,n
67  real(kind=kreal) :: dnn, coeff, xmax(3), xmin(3), dl
68 
69  !get average element size
70  xmax(1:3) = ecoord(1:3,1)
71  xmin(1:3) = ecoord(1:3,1)
72  do i=2,nn
73  do k=1,ndof
74  if( ecoord(k,i) > xmax(k) ) xmax(k) = ecoord(k,i)
75  if( ecoord(k,i) < xmin(k) ) xmin(k) = ecoord(k,i)
76  end do
77  end do
78  dl = 0.33333333333d0*((xmax(1)-xmin(1))+(xmax(2)-xmin(2))+(xmax(3)-xmin(3)))
79  if( dl < 1.d-8 ) dl = 1.d0
80 
81  coeff = element%elemact_coeff/dl
82  dnn = 1.d0/dble(nn)
83  stiff(:,:) = 0.d0
84  do i=1,nn
85  m = ndof*(i-1)
86 
87  do j=1,nn
88  n = ndof*(j-1)
89  do k=1,ndof
90  stiff(m+k,n+k) = stiff(m+k,n+k)-coeff*dnn
91  end do
92  end do
93 
94  do k=1,ndof
95  stiff(m+k,m+k) = stiff(m+k,m+k)+coeff
96  end do
97  end do
98 
99  end subroutine
100 
101 
102  subroutine update_dummy( ndof, nn, ecoord, u, du, qf, element )
103  use mmechgauss
104  integer(kind=kint), intent(in) :: ndof
105  integer(kind=kint), intent(in) :: nn
106  real(kind=kreal), intent(in) :: ecoord(3,nn)
107  real(kind=kreal), intent(in) :: u(3,nn)
108  real(kind=kreal), intent(in) :: du(3,nn)
109  type(telement), intent(inout) :: element
110  real(kind=kreal), intent(out) :: qf(:)
111 
112  integer(kind=kint) :: i,k,m
113  real(kind=kreal) :: dnn, coeff, aveu(3), xmax(3), xmin(3), dl
114 
115  !get average element size
116  xmax(1:3) = ecoord(1:3,1)
117  xmin(1:3) = ecoord(1:3,1)
118  do i=2,nn
119  do k=1,ndof
120  if( ecoord(k,i) > xmax(k) ) xmax(k) = ecoord(k,i)
121  if( ecoord(k,i) < xmin(k) ) xmin(k) = ecoord(k,i)
122  end do
123  end do
124  dl = 0.33333333333d0*((xmax(1)-xmin(1))+(xmax(2)-xmin(2))+(xmax(3)-xmin(3)))
125  if( dl < 1.d-8 ) dl = 1.d0
126 
127  coeff = element%elemact_coeff/dl
128  dnn = 1.d0/dble(nn)
129  aveu(:) = 0.d0
130  do i=1,nn
131  do k=1,ndof
132  aveu(k) = aveu(k) + u(k,i) + du(k,i)
133  end do
134  end do
135  aveu(:) = aveu(:)*dnn
136 
137  qf(:) = 0.d0
138  do i=1,nn
139  m = ndof*(i-1)
140 
141  do k=1,ndof
142  qf(m+k) = coeff*(u(k,i)+du(k,i)-aveu(k))
143  end do
144  end do
145 
146  do i=1,size(element%gausses)
147  element%gausses(i)%strain = 0.d0
148  element%gausses(i)%stress = 0.d0
149  element%gausses(i)%strain_out = 0.d0
150  element%gausses(i)%stress_out = 0.d0
151  if( associated(element%gausses(i)%istatus) ) element%gausses(i)%istatus = 0
152  if( associated(element%gausses(i)%fstatus) ) element%gausses(i)%fstatus = 0.d0
153  end do
154 
155  end subroutine
156 
157 
158 end module m_elemact
Definition: hecmw.f90:6
This module defined elemact data and function.
integer, parameter kelactd_none
integer, parameter kelact_active
subroutine update_dummy(ndof, nn, ecoord, u, du, qf, element)
integer, parameter kelactd_stress
integer, parameter kelact_undefined
integer, parameter kelact_inactive
integer, parameter kelactd_strain
subroutine print_elemact_info(elact)
subroutine stf_dummy(ndof, nn, ecoord, u, stiff, element)
This modules defines a structure to record history dependent parameter in static analysis.
Definition: mechgauss.f90:6
All data should be recorded in every elements.
Definition: mechgauss.f90:32