OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
sboltlaw.F File Reference
#include "implicit_f.inc"
#include "mvsiz_p.inc"
#include "com08_c.inc"
#include "param_c.inc"
#include "scr05_c.inc"

Go to the source code of this file.

Functions/Subroutines

subroutine sboltlaw (pm, sig, mat, d1, d2, d3, d4, d5, d6, nel, rho, bpreld, eint, qold, vol0, stifn, dt2t, neltst, ityptst, offg, geo, pid, mumax, ngl, ssp, dvol, aire, vnew, vd2, deltax, vis, pnew, psh, qnew, ssp_eq, sold1, sold2, sold3, sold4, sold5, sold6, mssa, dmels, conde, amu, vol_avg, dtel, g_dt, off, ipm, rhoref, rhosp, vol0dp, ismstr, jsph, jtur, ity, jthe, jsms, npg, glob_therm)

Function/Subroutine Documentation

◆ sboltlaw()

subroutine sboltlaw ( pm,
sig,
integer, dimension(*) mat,
d1,
d2,
d3,
d4,
d5,
d6,
integer nel,
rho,
bpreld,
eint,
qold,
vol0,
stifn,
dt2t,
integer neltst,
integer ityptst,
offg,
geo,
integer, dimension(*) pid,
mumax,
integer, dimension(*) ngl,
ssp,
dvol,
aire,
vnew,
vd2,
deltax,
vis,
pnew,
psh,
qnew,
ssp_eq,
sold1,
sold2,
sold3,
sold4,
sold5,
sold6,
mssa,
dmels,
conde,
amu,
vol_avg,
dtel,
integer g_dt,
off,
integer, dimension(npropmi,*) ipm,
rhoref,
rhosp,
double precision, dimension(*) vol0dp,
integer, intent(in) ismstr,
integer, intent(in) jsph,
integer, intent(in) jtur,
integer, intent(in) ity,
integer, intent(in) jthe,
integer, intent(in) jsms,
integer, intent(in) npg,
type (glob_therm_), intent(inout) glob_therm )

Definition at line 33 of file sboltlaw.F.

50C-----------------------------------------------
51C M o d u l e s
52C-----------------------------------------------
53 use glob_therm_mod
54C-----------------------------------------------
55C I m p l i c i t T y p e s
56C-----------------------------------------------
57#include "implicit_f.inc"
58C-----------------------------------------------
59C G l o b a l P a r a m e t e r s
60C-----------------------------------------------
61#include "mvsiz_p.inc"
62C-----------------------------------------------
63C C o m m o n B l o c k s
64C-----------------------------------------------
65#include "com08_c.inc"
66#include "param_c.inc"
67#include "scr05_c.inc"
68C-----------------------------------------------
69C D u m m y A r g u m e n t s
70C-----------------------------------------------
71 INTEGER, INTENT(IN) :: JSMS
72 INTEGER, INTENT(IN) :: JTUR
73 INTEGER, INTENT(IN) :: ITY
74 INTEGER, INTENT(IN) :: JTHE
75 INTEGER, INTENT(IN) :: ISMSTR
76 INTEGER, INTENT(IN) :: JSPH,NPG
77C
78 INTEGER NEL
80 . pm(npropm,*), sig(nel,6), rho(mvsiz),
81 . d1(*), d2(*), d3(*), d4(*),d5(*), d6(*),
82 . bpreld(nel,*), off(*)
83 INTEGER NELTST,ITYPTST,PID(*),G_DT
84 INTEGER MAT(*),NGL(*),IPM(NPROPMI,*)
86 . dt2t
87
89 . eint(*), qold(*),
90 . vol0(*),stifn(*), offg(*),geo(npropg,*),mumax(*)
92 . vnew(*), vd2(*), deltax(*), ssp(*), aire(*), vis(*),
93 . psh(*), pnew(*),qnew(*) ,ssp_eq(*), dvol(*),
94 . sold1(*), sold2(*), sold3(*), sold4(*), sold5(*), sold6(*),
95 . mssa(*), dmels(*),conde(*),amu(*),vol_avg(*),dtel(*),
96 . rhoref(*) ,rhosp(*)
97 DOUBLE PRECISION VOL0DP(*)
98 type (glob_therm_) ,intent(inout) :: glob_therm
99C-----------------------------------------------
100C L o c a l V a r i a b l e s
101C-----------------------------------------------
102 INTEGER I, MX,IBID !, J
103 my_real
104 . reduc ,reduc1 , t1, t2 , ee, nu, ggt, treps, p, tp
105 my_real
106 . c1(mvsiz), lmbd(mvsiz), gg(mvsiz),
107 . rho0(mvsiz), rho0_1,c1_1
108C-----------------------------------------------
109 my_real
110 . df,dpdm,facq0,
111 . e1, e2, e3, e4, e5, e6, einc,
112 . bid1, bid2, bid3, dta, ym, dpdmp
113C-----------------------------------------------
114 facq0 = one
115 reduc1 = em04 !EM03 !ONE
116
117 DO i=1,nel
118 t1 = bpreld(i,1)+zep4*(bpreld(i,2)-bpreld(i,1))
119 t2 = bpreld(i,1)+zep7*(bpreld(i,2)-bpreld(i,1))
120 reduc = reduc1
121 tp = tt-t1
122 IF(tp > zero ) THEN
123 reduc = min(reduc1*(one-tp/(t2-t1))+tp/(t2-t1), one)
124 ENDIF
125C
126 mx = mat(i)
127 rho0_1 = pm(1, mx)
128 nu = pm(21, mx)
129 ee = pm(20, mx)*reduc
130 gg(i) = pm(22, mx)*reduc
131 c1_1 = pm(32, mx)
132 lmbd(i) = dt1*ee*nu/((one+nu)*(one-two*nu))
133 c1(i) = c1_1*reduc
134 rho0(i) =rho0_1
135C
136 !SSP(I)=SQRT((ONEP333*GG(I)+C1(I))/RHO0(I))
137 ssp(i)=sqrt((onep333*pm(22, mx)+c1_1)/rho0(i))
138 ENDDO
139C
140 IF (jsph==0)THEN
141 CALL mqviscb(
142 1 pm, off, rho, bid1,
143 2 bid2, ssp, bid3, stifn,
144 3 dt2t, neltst, ityptst, aire,
145 4 offg, geo, pid, vnew,
146 5 vd2, deltax, vis, d1,
147 6 d2, d3, pnew, psh,
148 7 mat, ngl, qnew, ssp_eq,
149 8 vol0, mssa, dmels, ibid,
150 9 facq0, conde, dtel, g_dt,
151 a ipm, rhoref, rhosp, nel,
152 b ity, ismstr, jtur, jthe,
153 c jsms, npg , glob_therm)
154 ELSE
155 CALL mdtsph(
156 1 pm, off, rho, bid1,
157 2 bid2, bid3, stifn, dt2t,
158 3 neltst, ityptst, offg, geo,
159 4 pid, mumax, ssp, vnew,
160 5 vd2, deltax, vis, d1,
161 6 d2, d3, pnew, psh,
162 7 mat, ngl, qnew, ssp_eq,
163 8 g_dt, dtel, nel, ity,
164 9 jtur, jthe)
165 ENDIF
166
167 dta =half*dt1
168C
169 DO i=1,nel
170 pnew(i) = c1(i)*amu(i) !
171 treps = d1(i)+d2(i)+d3(i)
172 ggt = dt1*gg(i)
173 sig(i,1) = sig(i,1)+two*ggt*d1(i)+lmbd(i)*treps
174 sig(i,2) = sig(i,2)+two*ggt*d2(i)+lmbd(i)*treps
175 sig(i,3) = sig(i,3)+two*ggt*d3(i)+lmbd(i)*treps
176 sig(i,4) = sig(i,4)+ggt*d4(i)
177 sig(i,5) = sig(i,5)+ggt*d5(i)
178 sig(i,6) = sig(i,6)+ggt*d6(i)
179 qold(i) = qnew(i)
180 ENDDO
181 !write(*,'("***** sboltlaw ")')
182 !DO I=1,NEL
183 ! write(*,'(i5,6e10.3)')i,sig(i,1:6)
184 !enddo
185C
186C - To avoid discontinuous element pressure when switching material law
187C cf Law1, Law2, etc <=> Pressure is computed using total formulation P = C1*(rho/rho0-1)
188 IF(ismstr == 1)THEN
189 IF(iresp==0)THEN
190 DO i=1,nel
191 IF(tt >= (bpreld(i,2)-two*dt1)) THEN
192 p = -third*(sig(i,1)+sig(i,2)+sig(i,3))
193 rho(i)=rho0_1*(one+p/c1_1) ! rho will be used in next cycles to compute
194 ! rho_new = rho_old - rho0 dV
195 ! AMU = rho/rho0 - 1 in double precision
196 ENDIF
197 ENDDO
198 ELSE ! Single precision
199 DO i=1,nel
200 IF(tt >= (bpreld(i,2)-two*dt1)) THEN
201 p = -third*(sig(i,1)+sig(i,2)+sig(i,3))
202 rho(i)=rho0_1*(one+p/c1_1)
203 vol0dp(i) = vol0(i)*(one+p/c1_1) ! VOL0DP will be used in next cycle to compute
204 ! AMU = VOL0DP/VOLDP - 1 in single precision, VOLDP=V0 in small strain
205 ENDIF
206 ENDDO
207 END IF
208 ELSEIF(ismstr==2)THEN
209 IF(iresp==0)THEN
210 DO i=1,nel
211 IF(offg(i) > one)THEN
212 IF(tt >= (bpreld(i,2)-two*dt1)) THEN
213 p = -third*(sig(i,1)+sig(i,2)+sig(i,3))
214 rho(i)=rho0_1*(one+p/c1_1) ! rho will be used in next cycles to compute
215 ! rho_new = rho_old - rho0 dV
216 ! AMU = rho/rho0 - 1 in double precision
217 ENDIF
218 ELSE
219 p = -third*(sig(i,1)+sig(i,2)+sig(i,3))
220 vol0(i) = vnew(i)*(one+p/c1_1)
221 rho(i) = rho0_1*vol0(i)/vnew(i) ! cf rho = rho0 * V0 / V at each time step
222 ! AMU = rho/rho0 - 1 in double precision
223
224 END IF
225 ENDDO
226 ELSE ! Single precision
227 DO i=1,nel
228 IF(offg(i) > one)THEN
229 IF(tt >= (bpreld(i,2)-two*dt1)) THEN
230 p = -third*(sig(i,1)+sig(i,2)+sig(i,3))
231 rho(i) = rho0_1*(one+p/c1_1)
232 vol0dp(i) = vol0(i)*(one+p/c1_1) ! VOL0DP will be used in next cycle to compute
233 ! AMU = VOL0DP/VOLDP - 1 in single precision, VOLDP=V0 in small strain
234 ENDIF
235 ELSE
236 p = -third*(sig(i,1)+sig(i,2)+sig(i,3))
237 vol0(i) = vnew(i)*(one+p/c1_1)
238 vol0dp(i) = vol0(i) ! cf rho = rho0 * VOL0DP / VOLDP
239 ! AMU = VOL0DP/VOLDP - 1 in single precision, VOLDP = V in large strain
240 rho(i) = rho0_1*vol0(i)/vnew(i)
241 END IF
242 ENDDO
243 END IF
244 ELSEIF(ismstr<=4)THEN
245 IF(iresp==0)THEN
246 DO i=1,nel
247 p = -third*(sig(i,1)+sig(i,2)+sig(i,3))
248 vol0(i)=vnew(i)*(one+p/c1_1)
249 rho(i)= rho0_1*vol0(i)/vnew(i) ! cf rho = rho0 * V0 / V at each time step
250 ! AMU = rho/rho0 - 1 in double precision
251 ENDDO
252 ELSE ! Single precision
253 DO i=1,nel
254 p = -third*(sig(i,1)+sig(i,2)+sig(i,3))
255 vol0(i) = vnew(i)*(one+p/c1_1)
256 rho(i) = rho0_1*vol0(i)/vnew(i)
257 vol0dp(i) = vol0(i) ! cf rho = rho0 * VOL0DP / VOLDP
258 ! AMU = VOL0DP/VOLDP - 1 in single precision, VOLDP = V in large strain
259 ENDDO
260 END IF
261 END IF
262C
263 DO i=1,nel
264 e1=d1(i)*(sold1(i)+sig(i,1))
265 e2=d2(i)*(sold2(i)+sig(i,2))
266 e3=d3(i)*(sold3(i)+sig(i,3))
267 e4=d4(i)*(sold4(i)+sig(i,4))
268 e5=d5(i)*(sold5(i)+sig(i,5))
269 e6=d6(i)*(sold6(i)+sig(i,6))
270 einc= zero !VOL_AVG(I)*(E1+E2+E3+E4+E5+E6)*DTA - HALF*DVOL(I)*(QOLD(I)+QNEW(I))
271 eint(i)=(eint(i)+einc*off(i)) / max(em15,vol0(i))
272 ENDDO
273C
274 RETURN
#define my_real
Definition cppsort.cpp:32
subroutine dtel(ssp, pm, geo, pid, mat, rho0, vis, deltax, aire, vol, dtx)
Definition dtel.F:46
#define min(a, b)
Definition macros.h:20
#define max(a, b)
Definition macros.h:21
subroutine mdtsph(pm, off, rho, rk, t, re, sti, dt2t, neltst, ityptst, offg, geo, pid, mumax, ssp, vol, vd2, deltax, vis, d1, d2, d3, pnew, psh, mat, ngl, qvis, ssp_eq, g_dt, dtsph, nel, ity, jtur, jthe)
Definition mdtsph.F:47
subroutine mqviscb(pm, off, rho, rk, temp, ssp, re, sti, dt2t, neltst, ityptst, aire, offg, geo, pid, vol, vd2, deltax, vis, d1, d2, d3, pnew, psh, mat, ngl, qvis, ssp_eq, vol0, mssa, dmels, igeo, facq0, conde, dtel, g_dt, ipm, rhoref, rhosp, nel, ity, ismstr, jtur, jthe, jsms, npg, glob_therm)
Definition mqviscb.F:56