CalculiX  2.13
A Free Software Three-Dimensional Structural Finite Element Program
anisomaxwavspd.f File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine anisomaxwavspd (elas, rho, iorth, wavspd)
 

Function/Subroutine Documentation

◆ anisomaxwavspd()

subroutine anisomaxwavspd ( real*8, dimension(21), intent(inout)  elas,
real*8, intent(in)  rho,
integer, intent(in)  iorth,
real*8, intent(inout)  wavspd 
)
20 !
21 ! Calculates the propagation wave speed in a material, up to its 21
22 ! constants. Subroutine for calcmatwavsps.f
23 
24 ! Based on the procedure in:
25 ! C. Lane. The Development of a 2D Ultrasonic Array Inspection
26 ! for Single Crystal Turbine Blades.
27 ! Switzerland: Springer International Publishing, 2014.
28 !
29 ! CARLO MONJARAZ TEC (CMT)
30 !
31 ! INPUT:
32 !
33 ! elas: double(21) - The elasticity vector, containing 21 entries.
34 ! Non used are zero. If material is orthotropic, values are
35 ! rearranged to match indexes from anisotropic
36 ! material card.
37 !
38 ! rho: double - Density of the material
39 !
40 ! iorth: INTEGER - if the value is 1 : material is iorthtropic
41 ! for other vaules: material is anisotropic
42 !
43 ! OUTPUT:
44 !
45 ! wavspd
46 !
47  implicit none
48 !
49  integer i,j,k,im,imin,jm,jmin,iorth
50 !
51  real*8 elas(21),c(3,3,3,3),rho,xi(-1:1,-1:1),et(-1:1,-1:1),
52  & wavspd,d1,distmin,a
53 !
54  intent(in) rho,iorth
55 !
56  intent(inout) elas,wavspd
57 !
58  write(*,*)'++cMT: calculating max. speed in ANISOTROPIC...'
59 !
60 !--------IF IORTHTROPIC-----------------------------
61  if(iorth.eq.1)then
62 !
63  elas(10)=elas(9)
64  elas(15)=elas(8)
65  elas(21)=elas(7)
66  elas(9)=0.d0
67  elas(8)=0.d0
68  elas(7)=0.d0
69 !
70  endif
71 !
72 !--------FIlling c voigt Matrix-----------------------------
73 !
74  call anisotropic(elas,c)
75 !
76  d1=1.d0
77 !
78  xi(0,0)=0.d0
79  et(0,0)=0.d0
80  call inversewavspd(xi(0,0),et(0,0),c,rho,a)
81  distmin=a
82  imin=0
83  jmin=0
84 !
85  do k=1,8
86 !
87 ! initialisation
88 !
89  d1=d1/10.d0
90 !
91  do i=-1,1
92  do j=-1,1
93  if((i.eq.0).and.(j.eq.0)) cycle
94 !
95  xi(i,j)=xi(0,0)+i*d1
96  et(i,j)=et(0,0)+j*d1
97 !
98 ! check whether inside the (-1,1)x(-1,1) domain
99 !
100  if((xi(i,j).le.1.d0).and.
101  & (xi(i,j).ge.-1.d0).and.
102  & (et(i,j).le.1.d0).and.
103  & (et(i,j).ge.-1.d0)) then
104  call inversewavspd(xi(i,j),et(i,j),c,rho,a)
105 !
106 ! checking for smallest initial distance
107 !
108  if(a.lt.distmin) then
109  distmin=a
110  imin=i
111  jmin=j
112  endif
113  endif
114 !
115  enddo
116  enddo
117 !
118 ! minimizing the distance from the face to the node
119 !
120  do
121 !
122 ! exit if minimum found
123 !
124  if((imin.eq.0).and.(jmin.eq.0)) exit
125 !
126 ! new center of 3x3 matrix
127 !
128  xi(0,0)=xi(imin,jmin)
129  et(0,0)=et(imin,jmin)
130 !
131  im=imin
132  jm=jmin
133 !
134  imin=0
135  jmin=0
136 !
137  do i=-1,1
138  do j=-1,1
139  if((i+im.lt.-1).or.(i+im.gt.1).or.
140  & (j+jm.lt.-1).or.(j+jm.gt.1)) then
141 !
142  xi(i,j)=xi(0,0)+i*d1
143  et(i,j)=et(0,0)+j*d1
144 !
145 ! check whether inside the (-1,1)x(-1,1) domain
146 !
147  if((xi(i,j).le.1.d0).and.
148  & (xi(i,j).ge.-1.d0).and.
149  & (et(i,j).le.1.d0).and.
150  & (et(i,j).ge.-1.d0)) then
151  call inversewavspd(xi(i,j),et(i,j),c,rho,a)
152 !
153 ! check for new minimum
154 !
155  if(a.lt.distmin) then
156  distmin=a
157  imin=i
158  jmin=j
159  endif
160  endif
161 !
162  endif
163  enddo
164  enddo
165  enddo
166  enddo
167 !
168  call inversewavspd(xi(0,0),et(0,0),c,rho,a)
169 !
170  wavspd=1.d0/a
171 !
172  return
subroutine inversewavspd(xi, et, c, rho, a)
Definition: inversewavspd.f:20
subroutine anisotropic(anisol, anisox)
Definition: anisotropic.f:20
Hosted by OpenAircraft.com, (Michigan UAV, LLC)