OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
sdlensh14.F File Reference
#include "implicit_f.inc"
#include "mvsiz_p.inc"

Go to the source code of this file.

Functions/Subroutines

subroutine sdlensh14 (nel, llsh, x1, x2, x3, x4, x5, x6, x7, x8, y1, y2, y3, y4, y5, y6, y7, y8, z1, z2, z3, z4, z5, z6, z7, z8, ics, idt1sol)

Function/Subroutine Documentation

◆ sdlensh14()

subroutine sdlensh14 ( integer, intent(in) nel,
intent(out) llsh,
intent(in) x1,
intent(in) x2,
intent(in) x3,
intent(in) x4,
intent(in) x5,
intent(in) x6,
intent(in) x7,
intent(in) x8,
intent(in) y1,
intent(in) y2,
intent(in) y3,
intent(in) y4,
intent(in) y5,
intent(in) y6,
intent(in) y7,
intent(in) y8,
intent(in) z1,
intent(in) z2,
intent(in) z3,
intent(in) z4,
intent(in) z5,
intent(in) z6,
intent(in) z7,
intent(in) z8,
integer, intent(in) ics,
integer, intent(in) idt1sol )

Definition at line 30 of file sdlensh14.F.

38C-----------------------------------------------
39C I m p l i c i t T y p e s
40C-----------------------------------------------
41#include "implicit_f.inc"
42C-----------------------------------------------
43C G l o b a l P a r a m e t e r s
44C-----------------------------------------------
45#include "mvsiz_p.inc"
46C-----------------------------------------------
47C D u m m y A r g u m e n t s
48C-----------------------------------------------
49 INTEGER, INTENT(IN) :: NEL,ICS,IDT1SOL
50 my_real, DIMENSION(MVSIZ) , INTENT(OUT) :: llsh
51 my_real, DIMENSION(MVSIZ) , INTENT(IN) ::
52 . x1, x2, x3, x4, x5, x6, x7, x8,
53 . y1, y2, y3, y4, y5, y6, y7, y8,
54 . z1, z2, z3, z4, z5, z6, z7, z8
55C-----------------------------------------------
56C L o c a l V a r i a b l e s
57C-----------------------------------------------
58 INTEGER I, J, N
60 . rx(mvsiz),ry(mvsiz),rz(mvsiz),sx(mvsiz),sy(mvsiz),sz(mvsiz),
61 . vq(3,3,mvsiz), lxyz0(3),deta1(mvsiz),xx,yy,zz,
62 . xl2(mvsiz),xl3(mvsiz),xl4(mvsiz),yl2(mvsiz),
63 . yl3(mvsiz),yl4(mvsiz),zl1(mvsiz),
64 . xn(mvsiz,4) , yn(mvsiz,4) , zn(mvsiz,4) ,area(mvsiz)
66 . al1,al2,ll(mvsiz),corel(2,4)
68 . x13,x24,y13,y24,l13,l24,c1,c2,thkly,posly,
69 . fac,visce,rx1,ry1,sx1,sy1,s1,fac1,fac2,faci,fac11,facdt
70C=======================================================================
71 SELECT CASE(ics)
72 CASE (1)
73 DO i=1,nel
74 xn(i,1) = half*(x1(i)+x4(i))
75 yn(i,1) = half*(y1(i)+y4(i))
76 zn(i,1) = half*(z1(i)+z4(i))
77 xn(i,2) = half*(x2(i)+x3(i))
78 yn(i,2) = half*(y2(i)+y3(i))
79 zn(i,2) = half*(z2(i)+z3(i))
80 xn(i,3) = half*(x6(i)+x7(i))
81 yn(i,3) = half*(y6(i)+y7(i))
82 zn(i,3) = half*(z6(i)+z7(i))
83 xn(i,4) = half*(x5(i)+x8(i))
84 yn(i,4) = half*(y5(i)+y8(i))
85 zn(i,4) = half*(z5(i)+z8(i))
86 ENDDO
87 CASE (10)
88 DO i=1,nel
89 xn(i,1) = half*(x1(i)+x5(i))
90 yn(i,1) = half*(y1(i)+y5(i))
91 zn(i,1) = half*(z1(i)+z5(i))
92 xn(i,2) = half*(x2(i)+x6(i))
93 yn(i,2) = half*(y2(i)+y6(i))
94 zn(i,2) = half*(z2(i)+z6(i))
95 xn(i,3) = half*(x3(i)+x7(i))
96 yn(i,3) = half*(y3(i)+y7(i))
97 zn(i,3) = half*(z3(i)+z7(i))
98 xn(i,4) = half*(x4(i)+x8(i))
99 yn(i,4) = half*(y4(i)+y8(i))
100 zn(i,4) = half*(z4(i)+z8(i))
101 ENDDO
102 CASE (100)
103 DO i=1,nel
104 xn(i,1) = half*(x1(i)+x2(i))
105 yn(i,1) = half*(y1(i)+y2(i))
106 zn(i,1) = half*(z1(i)+z2(i))
107 xn(i,2) = half*(x5(i)+x6(i))
108 yn(i,2) = half*(y5(i)+y6(i))
109 zn(i,2) = half*(z5(i)+z6(i))
110 xn(i,3) = half*(x8(i)+x7(i))
111 yn(i,3) = half*(y8(i)+y7(i))
112 zn(i,3) = half*(z8(i)+z7(i))
113 xn(i,4) = half*(x4(i)+x3(i))
114 yn(i,4) = half*(y4(i)+y3(i))
115 zn(i,4) = half*(z4(i)+z3(i))
116 ENDDO
117 END SELECT
118C------g1,g2 :
119 DO i=1,nel
120 rx(i)=xn(i,2)+xn(i,3)-xn(i,1)-xn(i,4)
121 ry(i)=yn(i,2)+yn(i,3)-yn(i,1)-yn(i,4)
122 rz(i)=zn(i,2)+zn(i,3)-zn(i,1)-zn(i,4)
123 sx(i)=xn(i,3)+xn(i,4)-xn(i,1)-xn(i,2)
124 sy(i)=yn(i,3)+yn(i,4)-yn(i,1)-yn(i,2)
125 sz(i)=zn(i,3)+zn(i,4)-zn(i,1)-zn(i,2)
126 ENDDO
127C------Local elem. base:
128 CALL clsys3(rx, ry, rz, sx, sy, sz,
129 . vq, deta1,nel,mvsiz)
130C------ Global -> Local Coordinate FOURTH=0.25 ;
131 DO i=1,nel
132 lxyz0(1)=fourth*(xn(i,1)+xn(i,2)+xn(i,3)+xn(i,4))
133 lxyz0(2)=fourth*(yn(i,1)+yn(i,2)+yn(i,3)+yn(i,4))
134 lxyz0(3)=fourth*(zn(i,1)+zn(i,2)+zn(i,3)+zn(i,4))
135 xx=xn(i,2)-xn(i,1)
136 yy=yn(i,2)-yn(i,1)
137 zz=zn(i,2)-zn(i,1)
138 xl2(i)=vq(1,1,i)*xx+vq(2,1,i)*yy+vq(3,1,i)*zz
139 yl2(i)=vq(1,2,i)*xx+vq(2,2,i)*yy+vq(3,2,i)*zz
140 xx=xn(i,2)-lxyz0(1)
141 yy=yn(i,2)-lxyz0(2)
142 zz=zn(i,2)-lxyz0(3)
143 zl1(i)=vq(1,3,i)*xx+vq(2,3,i)*yy+vq(3,3,i)*zz
144C
145 xx=xn(i,3)-xn(i,1)
146 yy=yn(i,3)-yn(i,1)
147 zz=zn(i,3)-zn(i,1)
148 xl3(i)=vq(1,1,i)*xx+vq(2,1,i)*yy+vq(3,1,i)*zz
149 yl3(i)=vq(1,2,i)*xx+vq(2,2,i)*yy+vq(3,2,i)*zz
150C
151 xx=xn(i,4)-xn(i,1)
152 yy=yn(i,4)-yn(i,1)
153 zz=zn(i,4)-zn(i,1)
154 xl4(i)=vq(1,1,i)*xx+vq(2,1,i)*yy+vq(3,1,i)*zz
155 yl4(i)=vq(1,2,i)*xx+vq(2,2,i)*yy+vq(3,2,i)*zz
156 area(i)=fourth*deta1(i)
157 ENDDO
158 fac = two
159 facdt = five_over_4
160C-------same than QBAT
161 IF (idt1sol>0) facdt =four_over_3
162C---- compute COREL(2,4) mean surface and area
163 DO i=1,nel
164 lxyz0(1)=fourth*(xl2(i)+xl3(i)+xl4(i))
165 lxyz0(2)=fourth*(yl2(i)+yl3(i)+yl4(i))
166 corel(1,1)=-lxyz0(1)
167 corel(1,2)=xl2(i)-lxyz0(1)
168 corel(1,3)=xl3(i)-lxyz0(1)
169 corel(1,4)=xl4(i)-lxyz0(1)
170 corel(2,1)=-lxyz0(2)
171 corel(2,2)=yl2(i)-lxyz0(2)
172 corel(2,3)=yl3(i)-lxyz0(2)
173 corel(2,4)=yl4(i)-lxyz0(2)
174 x13=(corel(1,1)-corel(1,3))*half
175 x24=(corel(1,2)-corel(1,4))*half
176 y13=(corel(2,1)-corel(2,3))*half
177 y24=(corel(2,2)-corel(2,4))*half
178C
179 l13=x13*x13+y13*y13
180 l24=x24*x24+y24*y24
181 al1=max(l13,l24)
182 c1 =corel(1,2)*corel(2,4)-corel(2,2)*corel(1,4)
183 c2 =corel(1,1)*corel(2,3)-corel(2,1)*corel(1,3)
184 al2 =max(abs(c1),abs(c2))/area(i)
185 rx1=x24-x13
186 ry1=y24-y13
187 sx1=-x24-x13
188 sy1=-y24-y13
189 c1=sqrt(rx1*rx1+ry1*ry1)
190 c2=sqrt(sx1*sx1+sy1*sy1)
191 s1=fourth*(max(c1,c2)/min(c1,c2)-one)
192 fac1=min(half,s1)+one
193 fac2=area(i)/(c1*c2)
194 fac2=3.413*max(zero,fac2-0.7071)
195 fac2=0.78+0.22*fac2*fac2*fac2
196 faci=two*fac1*fac2
197 s1 = sqrt(faci*(facdt+al2)*al1)
198 s1 = max(s1,em20)
199 llsh(i) = area(i)/s1
200 ENDDO
201C
202 RETURN
#define my_real
Definition cppsort.cpp:32
subroutine area(d1, x, x2, y, y2, eint, stif0)
#define min(a, b)
Definition macros.h:20
#define max(a, b)
Definition macros.h:21
subroutine clsys3(rx, ry, rz, sx, sy, sz, vq, det, nel, mvsiz)
Definition scinit3.F:715