OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
sdlensh.F
Go to the documentation of this file.
1Copyright> OpenRadioss
2Copyright> Copyright (C) 1986-2025 Altair Engineering Inc.
3Copyright>
4Copyright> This program is free software: you can redistribute it and/or modify
5Copyright> it under the terms of the GNU Affero General Public License as published by
6Copyright> the Free Software Foundation, either version 3 of the License, or
7Copyright> (at your option) any later version.
8Copyright>
9Copyright> This program is distributed in the hope that it will be useful,
10Copyright> but WITHOUT ANY WARRANTY; without even the implied warranty of
11Copyright> MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12Copyright> GNU Affero General Public License for more details.
13Copyright>
14Copyright> You should have received a copy of the GNU Affero General Public License
15Copyright> along with this program. If not, see <https://www.gnu.org/licenses/>.
16Copyright>
17Copyright>
18Copyright> Commercial Alternative: Altair Radioss Software
19Copyright>
20Copyright> As an alternative to this open-source version, Altair also offers Altair Radioss
21Copyright> software under a commercial license. Contact Altair to discuss further if the
22Copyright> commercial version may interest you: https://www.altair.com/radioss/.
23!||====================================================================
24!|| sdlensh ../engine/source/elements/thickshell/solidec/sdlensh.F
25!||--- called by ------------------------------------------------------
26!|| scforc3 ../engine/source/elements/thickshell/solidec/scforc3.F
27!||--- calls -----------------------------------------------------
28!|| clsys3 ../engine/source/output/h3d/h3d_results/h3d_shell_tensor.f
29!||====================================================================
30 SUBROUTINE sdlensh(
31 1 VOLN, LLSH, AREA, X1, X2,
32 2 X3, X4, X5, X6,
33 3 X7, X8, Y1, Y2,
34 4 Y3, Y4, Y5, Y6,
35 5 Y7, Y8, Z1, Z2,
36 6 Z3, Z4, Z5, Z6,
37 7 Z7, Z8, NEL)
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 C o m m o n B l o c k s
48C-----------------------------------------------
49#include "scr17_c.inc"
50C-----------------------------------------------
51C D u m m y A r g u m e n t s
52C-----------------------------------------------
53 INTEGER, INTENT(IN) :: NEL
54 my_real, DIMENSION(MVSIZ) , INTENT(OUT) :: AREA
55 my_real
56 . VOLN(*),LLSH(*),
57 . X1(*), X2(*), X3(*), X4(*), X5(*), X6(*), X7(*), X8(*),
58 . Y1(*), Y2(*), Y3(*), Y4(*), Y5(*), Y6(*), Y7(*), Y8(*),
59 . Z1(*), Z2(*), Z3(*), Z4(*), Z5(*), Z6(*), Z7(*), Z8(*)
60C-----------------------------------------------
61C L o c a l V a r i a b l e s
62C-----------------------------------------------
63 INTEGER I, J, N
64 my_real
65 . RX(MVSIZ),RY(MVSIZ),RZ(MVSIZ),SX(MVSIZ),SY(MVSIZ),SZ(MVSIZ),
66 . VQ(3,3,MVSIZ), LXYZ0(3),DETA1(MVSIZ),XX,YY,ZZ,
67 . XL2(MVSIZ),XL3(MVSIZ),XL4(MVSIZ),YL2(MVSIZ),
68 . YL3(MVSIZ),YL4(MVSIZ),ZL1(MVSIZ),
69 . XN(4,MVSIZ) , YN(4,MVSIZ) , ZN(4,MVSIZ)
70 my_real
71 . al1,al2,ll(mvsiz),corel(2,4)
72 my_real
73 . x13,x24,y13,y24,l13,l24,c1,c2,thkly,posly,
74 . fac,visce,rx1,ry1,sx1,sy1,s1,fac1,fac2,faci,fac11,facdt
75C=======================================================================
76 DO i=1,nel
77 xn(1,i) = half*(x1(i)+x5(i))
78 yn(1,i) = half*(y1(i)+y5(i))
79 zn(1,i) = half*(z1(i)+z5(i))
80 xn(2,i) = half*(x2(i)+x6(i))
81 yn(2,i) = half*(y2(i)+y6(i))
82 zn(2,i) = half*(z2(i)+z6(i))
83 xn(3,i) = half*(x3(i)+x7(i))
84 yn(3,i) = half*(y3(i)+y7(i))
85 zn(3,i) = half*(z3(i)+z7(i))
86 xn(4,i) = half*(x4(i)+x8(i))
87 yn(4,i) = half*(y4(i)+y8(i))
88 zn(4,i) = half*(z4(i)+z8(i))
89 ENDDO
90C------g1,g2 :
91 DO i=1,nel
92 rx(i)=xn(2,i)+xn(3,i)-xn(1,i)-xn(4,i)
93 ry(i)=yn(2,i)+yn(3,i)-yn(1,i)-yn(4,i)
94 rz(i)=zn(2,i)+zn(3,i)-zn(1,i)-zn(4,i)
95 sx(i)=xn(3,i)+xn(4,i)-xn(1,i)-xn(2,i)
96 sy(i)=yn(3,i)+yn(4,i)-yn(1,i)-yn(2,i)
97 sz(i)=zn(3,i)+zn(4,i)-zn(1,i)-zn(2,i)
98 ENDDO
99C------Local elem. base:
100 CALL clsys3(rx, ry, rz, sx, sy, sz,
101 . vq, deta1,nel)
102C------ Global -> Local Coordinate FOURTH=0.25 ;
103 DO i=1,nel
104 lxyz0(1)=fourth*(xn(1,i)+xn(2,i)+xn(3,i)+xn(4,i))
105 lxyz0(2)=fourth*(yn(1,i)+yn(2,i)+yn(3,i)+yn(4,i))
106 lxyz0(3)=fourth*(zn(1,i)+zn(2,i)+zn(3,i)+zn(4,i))
107 xx=xn(2,i)-xn(1,i)
108 yy=yn(2,i)-yn(1,i)
109 zz=zn(2,i)-zn(1,i)
110 xl2(i)=vq(1,1,i)*xx+vq(2,1,i)*yy+vq(3,1,i)*zz
111 yl2(i)=vq(1,2,i)*xx+vq(2,2,i)*yy+vq(3,2,i)*zz
112 xx=xn(2,i)-lxyz0(1)
113 yy=yn(2,i)-lxyz0(2)
114 zz=zn(2,i)-lxyz0(3)
115 zl1(i)=vq(1,3,i)*xx+vq(2,3,i)*yy+vq(3,3,i)*zz
116C
117 xx=xn(3,i)-xn(1,i)
118 yy=yn(3,i)-yn(1,i)
119 zz=zn(3,i)-zn(1,i)
120 xl3(i)=vq(1,1,i)*xx+vq(2,1,i)*yy+vq(3,1,i)*zz
121 yl3(i)=vq(1,2,i)*xx+vq(2,2,i)*yy+vq(3,2,i)*zz
122C
123 xx=xn(4,i)-xn(1,i)
124 yy=yn(4,i)-yn(1,i)
125 zz=zn(4,i)-zn(1,i)
126 xl4(i)=vq(1,1,i)*xx+vq(2,1,i)*yy+vq(3,1,i)*zz
127 yl4(i)=vq(1,2,i)*xx+vq(2,2,i)*yy+vq(3,2,i)*zz
128 area(i)=fourth*deta1(i)
129 ENDDO
130 fac = two
131 facdt = five_over_4
132C-------same than QBAT
133 IF (idt1sol>0) facdt =four_over_3
134C---- compute COREL(2,4) mean surface and area
135 DO i=1,nel
136 lxyz0(1)=fourth*(xl2(i)+xl3(i)+xl4(i))
137 lxyz0(2)=fourth*(yl2(i)+yl3(i)+yl4(i))
138 corel(1,1)=-lxyz0(1)
139 corel(1,2)=xl2(i)-lxyz0(1)
140 corel(1,3)=xl3(i)-lxyz0(1)
141 corel(1,4)=xl4(i)-lxyz0(1)
142 corel(2,1)=-lxyz0(2)
143 corel(2,2)=yl2(i)-lxyz0(2)
144 corel(2,3)=yl3(i)-lxyz0(2)
145 corel(2,4)=yl4(i)-lxyz0(2)
146 x13=(corel(1,1)-corel(1,3))*half
147 x24=(corel(1,2)-corel(1,4))*half
148 y13=(corel(2,1)-corel(2,3))*half
149 y24=(corel(2,2)-corel(2,4))*half
150C
151 l13=x13*x13+y13*y13
152 l24=x24*x24+y24*y24
153 al1=max(l13,l24)
154 c1 =corel(1,2)*corel(2,4)-corel(2,2)*corel(1,4)
155 c2 =corel(1,1)*corel(2,3)-corel(2,1)*corel(1,3)
156 al2 =max(abs(c1),abs(c2))/area(i)
157 rx1=x24-x13
158 ry1=y24-y13
159 sx1=-x24-x13
160 sy1=-y24-y13
161 c1=sqrt(rx1*rx1+ry1*ry1)
162 c2=sqrt(sx1*sx1+sy1*sy1)
163 s1=fourth*(max(c1,c2)/min(c1,c2)-one)
164 fac1=min(half,s1)+one
165 fac2=area(i)/(c1*c2)
166 fac2=3.413*max(zero,fac2-0.7071)
167 fac2=0.78+0.22*fac2*fac2*fac2
168 faci=two*fac1*fac2
169 s1 = sqrt(faci*(facdt+al2)*al1)
170 s1 = max(s1,em20)
171 llsh(i) = area(i)/s1
172 ENDDO
173C
174 RETURN
175 END
subroutine h3d_shell_tensor(elbuf_tab, shell_tensor, iparg, itens, invert, nelcut, el2fa, nbf, tens, epsdot, iadp, nbf_l, nbpart, iadg, x, ixc, igeo, ixtg, ipm, stack, id_elem, ity_elem, info1, info2, is_written_shell, ipartc, iparttg, layer_input, ipt_input, ply_input, gauss_input, iuvar_input, h3d_part, keyword, d, id, bufmat, mat_param, geo, drape_sh4n, drape_sh3n, drapeg)
#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
subroutine sdlensh(voln, llsh, area, x1, x2, x3, x4, x5, x6, x7, x8, y1, y2, y3, y4, y5, y6, y7, y8, z1, z2, z3, z4, z5, z6, z7, z8, nel)
Definition sdlensh.F:38