OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
scinit3.F File Reference
#include "implicit_f.inc"
#include "mvsiz_p.inc"
#include "com04_c.inc"
#include "param_c.inc"
#include "scr12_c.inc"
#include "scr17_c.inc"
#include "scry_c.inc"
#include "vect01_c.inc"

Go to the source code of this file.

Functions/Subroutines

subroutine scinit3 (elbuf_str, mas, ixs, pm, x, mss, detonators, geo, veul, ale_connectivity, iparg, dtelem, sigi, nel, skew, igeo, stifn, partsav, v, iparts, ipart, sigsp, nsigi, msnf, mssf, ipm, iuser, nsigs, volnod, bvolnod, vns, bns, wma, ptsol, bufmat, mcp, mcps, temp, npf, tf, mssa, strsglob, straglob, orthoglob, fail_ini, iloadp, facload, rnoise, perturb, glob_therm)
subroutine svalue0 (rho, vol, off, sig, eint, dtx, rhog, volg, offg, sigg, eintg, dtxg, nel)
subroutine sczero3 (rhog, sigg, eintg, nel)
subroutine sdlensh (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)
subroutine clsys3 (rx, ry, rz, sx, sy, sz, vq, det, nel, mvsiz)

Function/Subroutine Documentation

◆ clsys3()

subroutine clsys3 ( intent(in) rx,
intent(in) ry,
intent(in) rz,
intent(in) sx,
intent(in) sy,
intent(in) sz,
intent(out) vq,
intent(out) det,
integer, intent(in) nel,
integer, intent(in) mvsiz )

Definition at line 714 of file scinit3.F.

715C-----------------------------------------------
716C I m p l i c i t T y p e s
717C-----------------------------------------------
718#include "implicit_f.inc"
719C-----------------------------------------------
720C D u m m y A r g u m e n t s
721C-----------------------------------------------
722 INTEGER ,INTENT(IN) :: NEL,MVSIZ
723C
724 my_real,DIMENSION(MVSIZ),INTENT(IN) ::
725 . rx , ry , rz,
726 . sx , sy , sz
727 my_real,DIMENSION(MVSIZ),INTENT(OUT) :: det
728 my_real,DIMENSION(MVSIZ,3,3),INTENT(OUT) :: vq
729C-----------------------------------------------
730C C o m m o n B l o c k s
731C-----------------------------------------------
732C---------+---------+---+---+--------------------------------------------
733C VAR | SIZE |TYP| RW| DEFINITION
734C---------+---------+---+---+--------------------------------------------
735C NEL | 1 | I | R | SIZE OF THE ELEMENT GROUP NEL
736C RX | NEL | R | R | X-of covariant vector g1
737C RY | NEL | R | R | Y-of covariant vector g1
738C RZ | NEL | R | R | Z-of covariant vector g1
739C SX | NEL | R | R | X-of covariant vector g2
740C SY | NEL | R | R | Y-of covariant vector g2
741C SZ | NEL | R | R | Z-of covariant vector g2
742C VQ |3*3*NEL | R | W | Local elem sys bases
743C DET | NEL | R | W | det of g1 ^ g2
744C---------+---------+---+---+--------------------------------------------
745C-----------------------------------------------
746C L o c a l V a r i a b l e s
747C-----------------------------------------------
748 INTEGER I
749C
750 my_real
751 . e1x(nel), e1y(nel), e1z(nel),
752 . e2x(nel), e2y(nel), e2z(nel),
753 . e3x(nel), e3y(nel), e3z(nel),
754 . c1,c2,cc,c1c1,c2c2,c1_1,c2_1
755C~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
756 DO i=1,nel
757C---------E3------------
758 e3x(i) = ry(i) * sz(i) - rz(i) * sy(i)
759 e3y(i) = rz(i) * sx(i) - rx(i) * sz(i)
760 e3z(i) = rx(i) * sy(i) - ry(i) * sx(i)
761 det(i) = sqrt(e3x(i)*e3x(i) + e3y(i)*e3y(i) + e3z(i)*e3z(i))
762C ----- EM20=1.0E-20
763 det(i) = max(em20,det(i))
764 e3x(i) = e3x(i) / det(i)
765 e3y(i) = e3y(i) / det(i)
766 e3z(i) = e3z(i) / det(i)
767 ENDDO
768C~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
769 DO i=1,nel
770 c1c1 = rx(i)*rx(i) + ry(i)*ry(i) + rz(i)*rz(i)
771 c2c2 = sx(i)*sx(i) + sy(i)*sy(i) + sz(i)*sz(i)
772C ----- ZERO=0., ONE=1.0
773 IF(c1c1 /= zero) THEN
774 c2_1 = sqrt(c2c2/max(em20,c1c1))
775 c1_1 = one
776 ELSEIF(c2c2 /= zero)THEN
777 c2_1 = one
778 c1_1 = sqrt(c1c1/max(em20,c2c2))
779 END IF
780 e1x(i) = rx(i)*c2_1+(sy(i)*e3z(i)-sz(i)*e3y(i))*c1_1
781 e1y(i) = ry(i)*c2_1+(sz(i)*e3x(i)-sx(i)*e3z(i))*c1_1
782 e1z(i) = rz(i)*c2_1+(sx(i)*e3y(i)-sy(i)*e3x(i))*c1_1
783 ENDDO
784C
785 DO i=1,nel
786 c1 = sqrt(e1x(i)*e1x(i) + e1y(i)*e1y(i) + e1z(i)*e1z(i))
787 IF ( c1 /= zero) c1 = one / max(em20,c1)
788 e1x(i) = e1x(i)*c1
789 e1y(i) = e1y(i)*c1
790 e1z(i) = e1z(i)*c1
791 e2x(i) = e3y(i) * e1z(i) - e3z(i) * e1y(i)
792 e2y(i) = e3z(i) * e1x(i) - e3x(i) * e1z(i)
793 e2z(i) = e3x(i) * e1y(i) - e3y(i) * e1x(i)
794 ENDDO
795 DO i=1,nel
796 vq(i,1,1)=e1x(i)
797 vq(i,2,1)=e1y(i)
798 vq(i,3,1)=e1z(i)
799 vq(i,1,2)=e2x(i)
800 vq(i,2,2)=e2y(i)
801 vq(i,3,2)=e2z(i)
802 vq(i,1,3)=e3x(i)
803 vq(i,2,3)=e3y(i)
804 vq(i,3,3)=e3z(i)
805 ENDDO
806C-----------
807 RETURN
#define my_real
Definition cppsort.cpp:32
#define max(a, b)
Definition macros.h:21

◆ scinit3()

subroutine scinit3 ( type(elbuf_struct_), target elbuf_str,
mas,
integer, dimension(nixs,*) ixs,
pm,
x,
mss,
type(detonators_struct_) detonators,
geo,
veul,
type(t_ale_connectivity), intent(inout) ale_connectivity,
integer, dimension(*) iparg,
dtelem,
sigi,
integer nel,
skew,
integer, dimension(npropgi,*) igeo,
stifn,
partsav,
v,
integer, dimension(*) iparts,
integer, dimension(lipart1,*) ipart,
sigsp,
integer nsigi,
msnf,
mssf,
integer, dimension(npropmi,*) ipm,
integer iuser,
integer nsigs,
volnod,
bvolnod,
vns,
bns,
wma,
integer, dimension(*) ptsol,
bufmat,
mcp,
mcps,
temp,
integer, dimension(*) npf,
tf,
mssa,
integer, dimension(*) strsglob,
integer, dimension(*) straglob,
integer, dimension(*) orthoglob,
integer, dimension(*) fail_ini,
integer, dimension(sizloadp,*), intent(in) iloadp,
dimension(lfacload,*), intent(in) facload,
rnoise,
integer, dimension(nperturb) perturb,
type (glob_therm_), intent(in) glob_therm )

Definition at line 48 of file scinit3.F.

59C-----------------------------------------------
60C D e s c r i p t i o n
61C Initialize 8-nodes thick shell HQEPH, co-rotational formulation
62C-----------------------------------------------
63C M o d u l e s
64C-----------------------------------------------
65 USE message_mod
66 USE elbufdef_mod
70 use glob_therm_mod
71C-----------------------------------------------
72C I m p l i c i t T y p e s
73C-----------------------------------------------
74#include "implicit_f.inc"
75C-----------------------------------------------
76C G l o b a l P a r a m e t e r s
77C-----------------------------------------------
78#include "mvsiz_p.inc"
79C-----------------------------------------------
80C C o m m o n B l o c k s
81C-----------------------------------------------
82#include "com04_c.inc"
83#include "param_c.inc"
84#include "scr12_c.inc"
85#include "scr17_c.inc"
86#include "scry_c.inc"
87#include "vect01_c.inc"
88C-----------------------------------------------
89C D u m m y A r g u m e n t s
90C-----------------------------------------------
91 INTEGER NEL,NSIGI, IUSER, NSIGS
92 INTEGER IXS(NIXS,*),IPARG(*),IPARTS(*),IPART(LIPART1,*),
93 . IPM(NPROPMI,*),PTSOL(*), NPF(*),IGEO(NPROPGI,*),
94 . STRSGLOB(*),STRAGLOB(*),ORTHOGLOB(*),FAIL_INI(*),PERTURB(NPERTURB)
96 . mas(*), pm(npropm,*), x(*), geo(npropg,*),
97 . veul(lveul,*), dtelem(*),sigi(nsigs,*),skew(lskew,*),stifn(*),
98 . partsav(20,*), v(*), mss(8,*),
99 . sigsp(nsigi, *),msnf(*), mssf(8,*), wma(*),
100 . volnod(*), bvolnod(*), vns(8,*), bns(8,*), bufmat(*),
101 . mcp(*),mcps(8,*),temp(*), tf(*), mssa(*),rnoise(nperturb,*)
102 TYPE(ELBUF_STRUCT_), TARGET :: ELBUF_STR
103 INTEGER,INTENT(IN) :: ILOADP(SIZLOADP,*)
104 my_real,INTENT(IN) :: facload(lfacload,*)
105 TYPE(DETONATORS_STRUCT_)::DETONATORS
106 TYPE(t_ale_connectivity), INTENT(INOUT) :: ALE_CONNECTIVITY
107 type (glob_therm_) ,intent(in) :: glob_therm
108C-----------------------------------------------
109C L o c a l V a r i a b l e s
110C-----------------------------------------------
111 INTEGER NF1, IBID, I, NLAY,IGTYP,NLYMAX,IS,NUVAR,IREP,NCC,JHBE,
112 . IDEF, IP, IPANG, IPTHK, IPPOS, IPMAT,IG,IM,MTN0,IPID1,
113 . NPTR,NPTS,NPTT,L_PLA,L_SIGB
114 INTEGER IX1(MVSIZ), IX2(MVSIZ), IX3(MVSIZ), IX4(MVSIZ),
115 . IX5(MVSIZ), IX6(MVSIZ), IX7(MVSIZ), IX8(MVSIZ)
116 INTEGER MAT(MVSIZ), PID(MVSIZ), NGL(MVSIZ) , MAT0(MVSIZ)
117 CHARACTER(LEN=NCHARTITLE)::TITR1
118 my_real
119 . bid, fv, sti, zi,wi
120 my_real
121 . v8loc(51,mvsiz),volu(mvsiz),dtx(mvsiz),vzl(mvsiz),vzq(mvsiz),
122 . x1(mvsiz),x2(mvsiz),x3(mvsiz),x4(mvsiz),x5(mvsiz),x6(mvsiz),
123 . x7(mvsiz),x8(mvsiz),y1(mvsiz),y2(mvsiz),y3(mvsiz),y4(mvsiz),
124 . y5(mvsiz),y6(mvsiz),y7(mvsiz),y8(mvsiz),z1(mvsiz),z2(mvsiz),
125 . z3(mvsiz),z4(mvsiz),z5(mvsiz),z6(mvsiz),z7(mvsiz),z8(mvsiz),
126 . rx(mvsiz) ,ry(mvsiz) ,rz(mvsiz) ,sx(mvsiz) ,
127 . sy(mvsiz) ,sz(mvsiz) ,tx(mvsiz) ,ty(mvsiz) ,
128 . tz(mvsiz) ,e1x(mvsiz),e1y(mvsiz),e1z(mvsiz),e2x(mvsiz),
129 . e2y(mvsiz),e2z(mvsiz),e3x(mvsiz),e3y(mvsiz),e3z(mvsiz),
130 . f1x(mvsiz) ,f1y(mvsiz) ,f1z(mvsiz) ,
131 . f2x(mvsiz) ,f2y(mvsiz) ,f2z(mvsiz) ,gama(6,mvsiz),
132 . rhocp(mvsiz) ,temp0(mvsiz),angle(mvsiz),dtx0(mvsiz),
133 . deltax(mvsiz), aire(mvsiz),llsh(mvsiz)
134 my_real, DIMENSION(8,MVSIZ) :: bid8mvsiz
135 my_real, DIMENSION(MVSIZ) :: bidmvsiz
136 my_real :: tempel(nel)
137C--------------------------------
138C-----------------------------------------------
139 TYPE(G_BUFEL_) ,POINTER :: GBUF
140 TYPE(BUF_LAY_) ,POINTER :: BUFLY
141 TYPE(L_BUFEL_) ,POINTER :: LBUF
142 TYPE(BUF_MAT_) ,POINTER :: MBUF
143C-----------------------------------------------
144 my_real
145 . w_gauss(9,9),a_gauss(9,9)
146C-----------------------------------------------
147 DATA w_gauss /
148 1 2. ,0. ,0. ,
149 1 0. ,0. ,0. ,
150 1 0. ,0. ,0. ,
151 2 1. ,1. ,0. ,
152 2 0. ,0. ,0. ,
153 2 0. ,0. ,0. ,
154 3 0.555555555555556,0.888888888888889,0.555555555555556,
155 3 0. ,0. ,0. ,
156 3 0. ,0. ,0. ,
157 4 0.347854845137454,0.652145154862546,0.652145154862546,
158 4 0.347854845137454,0. ,0. ,
159 4 0. ,0. ,0. ,
160 5 0.236926885056189,0.478628670499366,0.568888888888889,
161 5 0.478628670499366,0.236926885056189,0. ,
162 5 0. ,0. ,0. ,
163 6 0.171324492379170,0.360761573048139,0.467913934572691,
164 6 0.467913934572691,0.360761573048139,0.171324492379170,
165 6 0. ,0. ,0. ,
166 7 0.129484966168870,0.279705391489277,0.381830050505119,
167 7 0.417959183673469,0.381830050505119,0.279705391489277,
168 7 0.129484966168870,0. ,0. ,
169 8 0.101228536290376,0.222381034453374,0.313706645877887,
170 8 0.362683783378362,0.362683783378362,0.313706645877887,
171 8 0.222381034453374,0.101228536290376,0. ,
172 9 0.081274388361574,0.180648160694857,0.260610696402935,
173 9 0.312347077040003,0.330239355001260,0.312347077040003,
174 9 0.260610696402935,0.180648160694857,0.081274388361574/
175 DATA a_gauss /
176 1 0. ,0. ,0. ,
177 1 0. ,0. ,0. ,
178 1 0. ,0. ,0. ,
179 2 -.577350269189626,0.577350269189626,0. ,
180 2 0. ,0. ,0. ,
181 2 0. ,0. ,0. ,
182 3 -.774596669241483,0. ,0.774596669241483,
183 3 0. ,0. ,0. ,
184 3 0. ,0. ,0. ,
185 4 -.861136311594053,-.339981043584856,0.339981043584856,
186 4 0.861136311594053,0. ,0. ,
187 4 0. ,0. ,0. ,
188 5 -.906179845938664,-.538469310105683,0. ,
189 5 0.538469310105683,0.906179845938664,0. ,
190 5 0. ,0. ,0. ,
191 6 -.932469514203152,-.661209386466265,-.238619186083197,
192 6 0.238619186083197,0.661209386466265,0.932469514203152,
193 6 0. ,0. ,0. ,
194 7 -.949107912342759,-.741531185599394,-.405845151377397,
195 7 0. ,0.405845151377397,0.741531185599394,
196 7 0.949107912342759,0. ,0. ,
197 8 -.960289856497536,-.796666477413627,-.525532409916329,
198 8 -.183434642495650,0.183434642495650,0.525532409916329,
199 8 0.796666477413627,0.960289856497536,0. ,
200 9 -.968160239507626,-.836031107326636,-.613371432700590,
201 9 -.324253423403809,0. ,0.324253423403809,
202 9 0.613371432700590,0.836031107326636,0.968160239507626/
203C-----------------------------------------------
204C S o u r c e L i n e s
205C=======================================================================
206 gbuf => elbuf_str%GBUF
207 lbuf => elbuf_str%BUFLY(1)%LBUF(1,1,1)
208 mbuf => elbuf_str%BUFLY(1)%MAT(1,1,1)
209 bufly => elbuf_str%BUFLY(1)
210 nptr = elbuf_str%NPTR
211 npts = elbuf_str%NPTS
212 nptt = elbuf_str%NPTT
213 nlay = elbuf_str%NLAY
214C
215 jeul = iparg(11)
216 irep = iparg(35)
217 igtyp = iparg(38)
218 jhbe = iparg(23)
219 nf1=nft+1
220 IF (jcvt==1.AND.isorth/=0) jcvt=2
221 IF (igtyp /= 22) isorth = 0
222 ibid = 0
223 idef = 0
224C
225 DO i=1,nel
226 rhocp(i) = pm(69,ixs(1,nft+i))
227 temp0(i) = pm(79,ixs(1,nft+i))
228 ENDDO
229 CALL sccoor3(x ,ixs(1,nf1),geo ,mat ,pid ,ngl ,
230 . ix1 ,ix2 ,ix3 ,ix4 ,ix5 ,ix6 ,ix7 ,ix8 ,
231 . x1 ,x2 ,x3 ,x4 ,x5 ,x6 ,x7 ,x8 ,
232 . y1 ,y2 ,y3 ,y4 ,y5 ,y6 ,y7 ,y8 ,
233 . z1 ,z2 ,z3 ,z4 ,z5 ,z6 ,z7 ,z8 ,
234 . rx ,ry ,rz ,sx ,sy ,sz ,tx ,ty ,tz ,
235 . e1x ,e1y ,e1z ,e2x ,e2y ,e2z ,e3x ,e3y ,e3z ,
236 . f1x ,f1y ,f1z ,f2x ,f2y ,f2z ,temp0, temp,glob_therm%NINTEMP)
237 IF (igtyp == 21 .OR. igtyp == 22) THEN
238 DO i=1,nel
239 angle(i) = geo(1,pid(i))
240 END DO
241 CALL scmorth3(pid ,geo ,igeo ,skew ,irep ,gbuf%GAMA ,
242 . rx ,ry ,rz ,sx ,sy ,sz ,tx ,ty ,tz ,
243 . e1x ,e1y ,e1z ,e2x ,e2y ,e2z ,e3x ,e3y ,e3z ,
244 . ngl ,angle,nsigi,sigsp,nsigs,sigi ,ixs(1,nf1) ,1 ,
245 . orthoglob,ptsol,nel)
246 IF (igtyp == 22) THEN
247 nlymax= 200
248 ipang = 200
249 ipthk = ipang+nlymax
250 ippos = ipthk+nlymax
251 ipmat = 100
252 ig=pid(1)
253 mtn0=mtn
254 DO i=1,nel
255 mat0(i) = mat(i)
256 dtx0(i) = ep20
257 ENDDO
258 END IF
259 END IF
260 CALL scderi3(nel,gbuf%VOL,jeul,veul(1,nf1),geo ,
261 . vzl ,vzq ,ngl ,pid ,
262 . x1 ,x2 ,x3 ,x4 ,x5 ,x6 ,x7 ,x8 ,
263 . y1 ,y2 ,y3 ,y4 ,y5 ,y6 ,y7 ,y8 ,
264 . z1 ,z2 ,z3 ,z4 ,z5 ,z6 ,z7 ,z8 , volu)
265C
266 CALL sdlen3(x1 ,x2 ,x3 ,x4 ,x5 ,x6 ,x7 ,x8 ,
267 . y1 ,y2 ,y3 ,y4 ,y5 ,y6 ,y7 ,y8 ,
268 . z1 ,z2 ,z3 ,z4 ,z5 ,z6 ,z7 ,z8 ,
269 . deltax, volu)
270 IF (idttsh > 0) THEN
271 CALL sdlensh(nel,llsh,
272 . x1, x2, x3, x4, x5, x6, x7, x8,
273 . y1, y2, y3, y4, y5, y6, y7, y8,
274 . z1, z2, z3, z4, z5, z6, z7, z8)
275 DO i=1,nel
276 IF (gbuf%IDT_TSH(i)>0)
277 . deltax(i)=max(llsh(i),deltax(i))
278 ENDDO
279 END IF
280C
281 IF (jthe == 0 .and. glob_therm%NINTEMP > 0) THEN
282 DO i=1,nel
283 tempel(i) = one_over_8 *(temp(ixs(2,i)) + temp(ixs(3,i))
284 . + temp(ixs(4,i)) + temp(ixs(5,i))
285 . + temp(ixs(6,i)) + temp(ixs(7,i))
286 . + temp(ixs(8,i)) + temp(ixs(9,i)))
287 ENDDO
288 ELSE
289 tempel(1:nel) = temp0(1:nel)
290 END IF
291!
292 ip=0
293 CALL matini(pm ,ixs ,nixs ,x ,
294 . geo ,ale_connectivity ,detonators ,iparg ,
295 . sigi ,nel ,skew ,igeo ,
296 . ipart ,iparts ,
297 . mat ,ipm ,nsigs ,numsol ,ptsol ,
298 . ip ,ngl ,npf ,tf ,bufmat ,
299 . gbuf ,lbuf ,mbuf ,elbuf_str ,iloadp ,
300 . facload, deltax ,tempel )
301C
302 IF (igtyp == 22) CALL sczero3(gbuf%RHO,gbuf%SIG,gbuf%EINT,nel)
303C----------------------------------------
304C Thermal initialization
305 IF(jthe /=0) CALL atheri(mat,pm,gbuf%TEMP)
306C------------------------
307C Loop on integration points
308 DO is=1,nlay
309C
310 lbuf => elbuf_str%BUFLY(is)%LBUF(1,1,1)
311 mbuf => elbuf_str%BUFLY(is)%MAT(1,1,1)
312 l_pla = elbuf_str%BUFLY(is)%L_PLA
313 l_sigb= elbuf_str%BUFLY(is)%L_SIGB
314C
315 IF (igtyp == 22) THEN
316 zi = geo(ippos+is,ig)
317 wi = geo(ipthk+is,ig)
318 im=igeo(ipmat+is,ig)
319 mtn=nint(pm(19,im))
320 DO i=1,nel
321 mat(i)=im
322 angle(i) = geo(ipang+is,pid(i))
323 ENDDO
324 ELSE
325 zi = a_gauss(is,nlay)
326 wi = w_gauss(is,nlay)
327 ENDIF
328 DO i=1,nel
329 lbuf%VOL0DP(i) = half*wi*(gbuf%VOL(i)+vzl(i)*zi)
330 lbuf%VOL(i) = lbuf%VOL0DP(i)
331 ENDDO
332 IF (igtyp == 22)
333 . CALL scmorth3(pid ,geo ,igeo ,skew ,irep ,lbuf%GAMA ,
334 . rx ,ry ,rz ,sx ,sy ,sz ,tx ,ty ,tz ,
335 . e1x ,e1y ,e1z ,e2x ,e2y ,e2z ,e3x ,e3y ,e3z ,
336 . ngl ,angle,nsigi,sigsp,nsigs,sigi ,ixs(1,nf1) ,is ,
337 . orthoglob,ptsol,nel)
338C
339 IF (jthe == 0 .and. glob_therm%NINTEMP > 0) THEN
340 DO i=1,nel
341 tempel(i) = one_over_8 *(temp(ixs(2,i)) + temp(ixs(3,i))
342 . + temp(ixs(4,i)) + temp(ixs(5,i))
343 . + temp(ixs(6,i)) + temp(ixs(7,i))
344 . + temp(ixs(8,i)) + temp(ixs(9,i)))
345 ENDDO
346 ELSE
347 tempel(1:nel) = temp0(1:nel)
348 END IF
349!
350 CALL matini(pm ,ixs ,nixs ,x ,
351 . geo ,ale_connectivity ,detonators ,iparg ,
352 . sigi ,nel ,skew ,igeo ,
353 . ipart ,iparts,
354 . mat ,ipm ,nsigs ,numsol ,ptsol,
355 . is ,ngl ,npf ,tf ,bufmat,
356 . gbuf ,lbuf ,mbuf ,elbuf_str,iloadp,
357 . facload, deltax ,tempel )
358c
359 nuvar = elbuf_str%BUFLY(is)%NVAR_MAT
360 IF(mtn>=28)THEN
361 idef =1
362 ELSE
363 IF(mtn == 14 .OR. mtn == 12)THEN
364 idef =1
365 ELSEIF(mtn == 24)THEN
366 idef =1
367 ELSEIF(istrain == 1)THEN
368 IF(mtn == 1)THEN
369 idef =1
370 ELSEIF(mtn == 2)THEN
371 idef =1
372 ELSEIF(mtn == 4)THEN
373 idef =1
374 ELSEIF(mtn == 3.OR.mtn == 6.OR.mtn == 10.OR.
375 . mtn == 21.OR.mtn == 22.OR.mtn == 23.
376 . or.mtn == 49)THEN
377 idef =1
378 ENDIF
379 ENDIF
380 ENDIF
381c
382 CALL sigin20b(lbuf%SIG,
383 . pm, lbuf%VOL,sigsp,sigi,lbuf%EINT,lbuf%RHO,mbuf%VAR ,
384 . lbuf%STRA,ixs ,nixs,nsigi, is, nuvar,nel,iuser,idef,
385 . nsigs,strsglob,straglob,jhbe,igtyp,x,gbuf%GAMA,
386 . mat ,lbuf%PLA,l_pla,ptsol,lbuf%SIGB,l_sigb,ipm ,
387 . bufmat,lbuf%VOL0DP)
388 IF(igtyp == 22) THEN
389C
390 aire(:) = zero
391 CALL dtmain(geo ,pm ,ipm ,pid ,mat ,fv ,
392 . lbuf%EINT ,lbuf%TEMP ,lbuf%DELTAX ,lbuf%RK ,lbuf%RE ,bufmat, deltax, aire,
393 . volu, dtx , igeo ,igtyp)
394C Average density, stresses, ...
395 CALL svalue0(
396 . lbuf%RHO,lbuf%VOL,lbuf%OFF,lbuf%SIG,lbuf%EINT,dtx,
397 . gbuf%RHO,gbuf%VOL,gbuf%OFF,gbuf%SIG,gbuf%EINT,dtx0,
398 . nel )
399 ENDIF
400 ENDDO ! IS=1,NLAY
401C----------------
402 IF(igtyp == 22) THEN
403 mtn=mtn0
404 DO i=1,nel
405 mat(i)=mat0(i)
406 ENDDO
407 ENDIF
408C----------------------------------------
409C Masses initialization
410 bid8mvsiz(1:8,1:mvsiz) = zero
411 bidmvsiz(1:mvsiz) = zero
412 CALL smass3(
413 . gbuf%RHO ,mas ,partsav ,x ,v ,
414 . iparts(nf1),mss(1,nf1),volu ,
415 . msnf ,mssf(1,nf1),bid ,
416 . bid ,bid8mvsiz ,wma ,rhocp ,mcp ,
417 . mcps(1,nf1),mssa ,bidmvsiz ,bidmvsiz ,gbuf%FILL,
418 . ix1, ix2, ix3, ix4, ix5, ix6, ix7, ix8)
419C----------------------------------------
420C Failure model initialization
421 CALL failini(elbuf_str,nptr,npts,nptt,nlay,
422 . ipm,sigsp,nsigi,fail_ini ,
423 . sigi,nsigs,ixs,nixs,ptsol,rnoise,perturb,bufmat)
424C------------------------------------------
425C Assemble nodal volumes and moduli for interface stiffness
426C Warning : IX1, IX2 ... IX8 <=> NC(MVSIZ,8)
427 IF(i7stifs/=0)THEN
428 ncc=8
429 CALL sbulk3(volu ,ix1 ,ncc,mat,pm ,
430 2 volnod,bvolnod,vns(1,nf1),bns(1,nf1),bid,
431 3 bid ,gbuf%FILL)
432 ENDIF
433C------------------------------------------
434C Element time step
435 aire(:) = zero
436 CALL dtmain(geo ,pm ,ipm ,pid ,mat ,fv ,
437 . lbuf%EINT ,lbuf%TEMP ,lbuf%DELTAX ,lbuf%RK ,lbuf%RE ,bufmat, deltax, aire,
438 . volu, dtx , igeo ,igtyp)
439C
440 IF(igtyp == 22) THEN
441 DO i=1,nel
442 dtx(i)=dtx0(i)
443 ENDDO
444 ENDIF
445C
446 DO i=1,nel
447 IF (ixs(10,i+nft) /= 0) THEN
448 IF (igtyp < 20 .OR. igtyp > 22) THEN
449 ipid1=ixs(nixs-1,i+nft)
450 CALL fretitl2(titr1,igeo(npropgi-ltitr+1,ipid1),ltitr)
451 CALL ancmsg(msgid=226,
452 . msgtype=msgerror,
453 . anmode=aninfo_blind_1,
454 . i1=igeo(1,ipid1),
455 . c1=titr1,
456 . i2=igtyp)
457 ENDIF
458 ENDIF
459 dtelem(nft+i)=dtx(i)
460C STI = 0.25 * RHO * VOL / (DT*DT)
461 sti = fourth * gbuf%FILL(i) * gbuf%RHO(i) * volu(i) /
462 . max(em20,dtx(i)*dtx(i))
463 stifn(ixs(2,i+nft))=stifn(ixs(2,i+nft))+sti
464 stifn(ixs(3,i+nft))=stifn(ixs(3,i+nft))+sti
465 stifn(ixs(4,i+nft))=stifn(ixs(4,i+nft))+sti
466 stifn(ixs(5,i+nft))=stifn(ixs(5,i+nft))+sti
467 stifn(ixs(6,i+nft))=stifn(ixs(6,i+nft))+sti
468 stifn(ixs(7,i+nft))=stifn(ixs(7,i+nft))+sti
469 stifn(ixs(8,i+nft))=stifn(ixs(8,i+nft))+sti
470 stifn(ixs(9,i+nft))=stifn(ixs(9,i+nft))+sti
471 END DO
472C-----------
473 RETURN
subroutine atheri(mat, pm, temp)
Definition atheri.F:41
subroutine dtmain(geo, pm, ipm, pid, mat, fv, eint, temp, deltax, rk, re, bufmat, ddeltax, aire, vol, dtx, igeo, igtyp)
Definition dtmain.F:67
subroutine failini(elbuf_str, nptr, npts, nptt, nlay, ipm, sigsp, nsigi, fail_ini, sigi, nsigs, ix, nix, pt, rnoise, perturb, mat_param)
Definition failini.F:43
subroutine matini(pm, ix, nix, x, geo, ale_connectivity, detonators, iparg, sigi, nel, skew, igeo, ipart, ipartel, mat, ipm, nsig, nums, pt, ipt, ngl, npf, tf, bufmat, gbuf, lbuf, mbuf, elbuf_str, iloadp, facload, ddeltax, tempel)
Definition matini.F:81
integer, parameter nchartitle
subroutine sigin20b(sig, pm, vol, sigsp, sigi, eint, rho, uvar, eps, ix, nix, nsigi, ipt, nuvar, nel, iuser, idef, nsigs, strsglob, straglob, jhbe, igtyp, x, bufgama, mat, epsp, l_pla, pt, sigb, l_sigb, ipm, bufmat, voldp)
Definition s20mass3.F:350
subroutine sbulk3(volu, nc, nnc, mat, pm, volnod, bvolnod, vns, bns, vnsx, bnsx, fill)
Definition sbulk3.F:42
subroutine sccoor3(x, ixs, geo, mxt, ngeo, ngl, ix1, ix2, ix3, ix4, ix5, ix6, ix7, ix8, x1, x2, x3, x4, x5, x6, x7, x8, y1, y2, y3, y4, y5, y6, y7, y8, z1, z2, z3, z4, z5, z6, z7, z8, rx, ry, rz, sx, sy, sz, tx, ty, tz, r11, r21, r31, r12, r22, r32, r13, r23, r33, f1x, f1y, f1z, f2x, f2y, f2z, temp0, temp, nintemp)
Definition sccoor3.F:43
subroutine sczero3(rhog, sigg, eintg, nel)
Definition scinit3.F:532
subroutine sdlensh(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)
Definition scinit3.F:571
subroutine svalue0(rho, vol, off, sig, eint, dtx, rhog, volg, offg, sigg, eintg, dtxg, nel)
Definition scinit3.F:487
subroutine scmorth3(pid, geo, igeo, skew, irep, gama, rx, ry, rz, sx, sy, sz, tx, ty, tz, e1x, e1y, e1z, e2x, e2y, e2z, e3x, e3y, e3z, ngl, angle, nsigi, sigsp, nsigs, sigi, ixs, ilay, orthoglob, pt, nel)
Definition scmorth3.F:40
subroutine smass3(rho, ms, partsav, x, v, ipart, mss, volu, msnf, mssf, in, vr, ins, wma, rhocp, mcp, mcps, mssa, rhof, frac, fill, nc1, nc2, nc3, nc4, nc5, nc6, nc7, nc8)
Definition smass3.F:44
subroutine sdlen3(x1, x2, x3, x4, x5, x6, x7, x8, y1, y2, y3, y4, y5, y6, y7, y8, z1, z2, z3, z4, z5, z6, z7, z8, deltax, voln)
Definition sdlen3.F:41
subroutine scderi3(nel, vol, jeul, veul, geo, vzl, vzq, ngl, ngeo, x1, x2, x3, x4, x5, x6, x7, x8, y1, y2, y3, y4, y5, y6, y7, y8, z1, z2, z3, z4, z5, z6, z7, z8, det)
Definition scderi3.F:37
subroutine ancmsg(msgid, msgtype, anmode, i1, i2, i3, i4, i5, i6, i7, i8, i9, i10, i11, i12, i13, i14, i15, i16, i17, i18, i19, i20, r1, r2, r3, r4, r5, r6, r7, r8, r9, c1, c2, c3, c4, c5, c6, c7, c8, c9, prmode)
Definition message.F:889
subroutine fretitl2(titr, iasc, l)
Definition freform.F:804

◆ sczero3()

subroutine sczero3 ( rhog,
sigg,
eintg,
integer nel )

Definition at line 531 of file scinit3.F.

532C-----------------------------------------------
533C I m p l i c i t T y p e s
534C-----------------------------------------------
535#include "implicit_f.inc"
536C-----------------------------------------------
537C D u m m y A r g u m e n t s
538C-----------------------------------------------
539 INTEGER NEL
540 my_real
541 . sigg(nel,6),eintg(*),rhog(*)
542C-----------------------------------------------
543C L o c a l V a r i a b l e s
544C-----------------------------------------------
545 INTEGER I
546C
547 DO i=1,nel
548 sigg(i,1) = zero
549 sigg(i,2) = zero
550 sigg(i,3) = zero
551 sigg(i,4) = zero
552 sigg(i,5) = zero
553 sigg(i,6) = zero
554 rhog(i) = zero
555 eintg(i) = zero
556 ENDDO
557C
558 RETURN

◆ sdlensh()

subroutine sdlensh ( integer 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 )

Definition at line 567 of file scinit3.F.

571C-----------------------------------------------
572C I m p l i c i t T y p e s
573C-----------------------------------------------
574#include "implicit_f.inc"
575C-----------------------------------------------
576C G l o b a l P a r a m e t e r s
577C-----------------------------------------------
578#include "mvsiz_p.inc"
579C-----------------------------------------------
580C C o m m o n B l o c k s
581C-----------------------------------------------
582#include "scr17_c.inc"
583C-----------------------------------------------
584C D u m m y A r g u m e n t s
585C-----------------------------------------------
586 INTEGER :: NEL
587 my_real,DIMENSION(MVSIZ),INTENT(IN) ::
588 . x1, x2, x3, x4, x5, x6, x7, x8,
589 . y1, y2, y3, y4, y5, y6, y7, y8,
590 . z1, z2, z3, z4, z5, z6, z7, z8
591 my_real,DIMENSION(MVSIZ),INTENT(OUT) :: llsh
592C-----------------------------------------------
593C L o c a l V a r i a b l e s
594C-----------------------------------------------
595 INTEGER I, J, N
596 my_real
597 . rx(mvsiz),ry(mvsiz),rz(mvsiz),sx(mvsiz),sy(mvsiz),sz(mvsiz),
598 . vq(mvsiz,3,3), lxyz0(3),deta1(mvsiz),xx,yy,zz,
599 . xl2(mvsiz),xl3(mvsiz),xl4(mvsiz),yl2(mvsiz),
600 . yl3(mvsiz),yl4(mvsiz),zl1(mvsiz),area(mvsiz),
601 . xn(mvsiz,4) , yn(mvsiz,4) , zn(mvsiz,4)
602 my_real
603 . al1,al2,ll(mvsiz),corel(2,4)
604 my_real
605 . x13,x24,y13,y24,l13,l24,c1,c2,thkly,posly,
606 . fac,visce,rx1,ry1,sx1,sy1,s1,fac1,fac2,faci,fac11,facdt
607C=======================================================================
608 DO i=1,nel
609 xn(i,1) = half*(x1(i)+x5(i))
610 yn(i,1) = half*(y1(i)+y5(i))
611 zn(i,1) = half*(z1(i)+z5(i))
612 xn(i,2) = half*(x2(i)+x6(i))
613 yn(i,2) = half*(y2(i)+y6(i))
614 zn(i,2) = half*(z2(i)+z6(i))
615 xn(i,3) = half*(x3(i)+x7(i))
616 yn(i,3) = half*(y3(i)+y7(i))
617 zn(i,3) = half*(z3(i)+z7(i))
618 xn(i,4) = half*(x4(i)+x8(i))
619 yn(i,4) = half*(y4(i)+y8(i))
620 zn(i,4) = half*(z4(i)+z8(i))
621 ENDDO
622C------g1,g2 :
623 DO i=1,nel
624 rx(i)=xn(i,2)+xn(i,3)-xn(i,1)-xn(i,4)
625 ry(i)=yn(i,2)+yn(i,3)-yn(i,1)-yn(i,4)
626 rz(i)=zn(i,2)+zn(i,3)-zn(i,1)-zn(i,4)
627 sx(i)=xn(i,3)+xn(i,4)-xn(i,1)-xn(i,2)
628 sy(i)=yn(i,3)+yn(i,4)-yn(i,1)-yn(i,2)
629 sz(i)=zn(i,3)+zn(i,4)-zn(i,1)-zn(i,2)
630 ENDDO
631C------Local elem. base:
632 CALL clsys3(rx, ry, rz, sx, sy, sz,
633 . vq, deta1,nel ,mvsiz)
634C------ Global -> Local Coordinate FOURTH=0.25 ;
635 DO i=1,nel
636 lxyz0(1)=fourth*(xn(i,1)+xn(i,2)+xn(i,3)+xn(i,4))
637 lxyz0(2)=fourth*(yn(i,1)+yn(i,2)+yn(i,3)+yn(i,4))
638 lxyz0(3)=fourth*(zn(i,1)+zn(i,2)+zn(i,3)+zn(i,4))
639 xx=xn(i,2)-xn(i,1)
640 yy=yn(i,2)-yn(i,1)
641 zz=zn(i,2)-zn(i,1)
642 xl2(i)=vq(i,1,1)*xx+vq(i,2,1)*yy+vq(i,3,1)*zz
643 yl2(i)=vq(i,1,2)*xx+vq(i,2,2)*yy+vq(i,3,2)*zz
644 xx=xn(i,2)-lxyz0(1)
645 yy=yn(i,2)-lxyz0(2)
646 zz=zn(i,2)-lxyz0(3)
647 zl1(i)=vq(i,1,3)*xx+vq(i,2,3)*yy+vq(i,3,3)*zz
648C
649 xx=xn(i,3)-xn(i,1)
650 yy=yn(i,3)-yn(i,1)
651 zz=zn(i,3)-zn(i,1)
652 xl3(i)=vq(i,1,1)*xx+vq(i,2,1)*yy+vq(i,3,1)*zz
653 yl3(i)=vq(i,1,2)*xx+vq(i,2,2)*yy+vq(i,3,2)*zz
654C
655 xx=xn(i,4)-xn(i,1)
656 yy=yn(i,4)-yn(i,1)
657 zz=zn(i,4)-zn(i,1)
658 xl4(i)=vq(i,1,1)*xx+vq(i,2,1)*yy+vq(i,3,1)*zz
659 yl4(i)=vq(i,1,2)*xx+vq(i,2,2)*yy+vq(i,3,2)*zz
660 area(i)=fourth*deta1(i)
661 ENDDO
662 fac = two
663 facdt = five_over_4
664C-------same as QBAT
665 IF (idt1sol>0) facdt =four_over_3
666C---- compute COREL(2,4) mean surface and area
667 DO i=1,nel
668 lxyz0(1)=fourth*(xl2(i)+xl3(i)+xl4(i))
669 lxyz0(2)=fourth*(yl2(i)+yl3(i)+yl4(i))
670 corel(1,1)=-lxyz0(1)
671 corel(1,2)=xl2(i)-lxyz0(1)
672 corel(1,3)=xl3(i)-lxyz0(1)
673 corel(1,4)=xl4(i)-lxyz0(1)
674 corel(2,1)=-lxyz0(2)
675 corel(2,2)=yl2(i)-lxyz0(2)
676 corel(2,3)=yl3(i)-lxyz0(2)
677 corel(2,4)=yl4(i)-lxyz0(2)
678 x13=(corel(1,1)-corel(1,3))*half
679 x24=(corel(1,2)-corel(1,4))*half
680 y13=(corel(2,1)-corel(2,3))*half
681 y24=(corel(2,2)-corel(2,4))*half
682C
683 l13=x13*x13+y13*y13
684 l24=x24*x24+y24*y24
685 al1=max(l13,l24)
686 c1 =corel(1,2)*corel(2,4)-corel(2,2)*corel(1,4)
687 c2 =corel(1,1)*corel(2,3)-corel(2,1)*corel(1,3)
688 al2 =max(abs(c1),abs(c2))/area(i)
689 rx1=x24-x13
690 ry1=y24-y13
691 sx1=-x24-x13
692 sy1=-y24-y13
693 c1=sqrt(rx1*rx1+ry1*ry1)
694 c2=sqrt(sx1*sx1+sy1*sy1)
695 s1=fourth*(max(c1,c2)/min(c1,c2)-one)
696 fac1=min(half,s1)+one
697 fac2=area(i)/(c1*c2)
698 fac2=3.413*max(zero,fac2-0.7071)
699 fac2=0.78+0.22*fac2*fac2*fac2
700 faci=two*fac1*fac2
701 s1 = sqrt(faci*(facdt+al2)*al1)
702 s1 = max(s1,em20)
703 llsh(i) = area(i)/s1
704 ENDDO
705C
706 RETURN
subroutine area(d1, x, x2, y, y2, eint, stif0)
#define min(a, b)
Definition macros.h:20
subroutine clsys3(rx, ry, rz, sx, sy, sz, vq, det, nel, mvsiz)
Definition scinit3.F:715

◆ svalue0()

subroutine svalue0 ( rho,
vol,
off,
sig,
eint,
dtx,
rhog,
volg,
offg,
sigg,
eintg,
dtxg,
integer nel )

Definition at line 484 of file scinit3.F.

487C
488C-----------------------------------------------
489C I m p l i c i t T y p e s
490C-----------------------------------------------
491#include "implicit_f.inc"
492C-----------------------------------------------
493C D u m m y A r g u m e n t s
494C-----------------------------------------------
495 INTEGER NEL
496 my_real
497 . rho(*), vol(*),sig(nel,6),eint(*),off(*),dtx(*),
498 . sigg(nel,6),eintg(*),rhog(*),offg(*),volg(*),dtxg(*)
499C-----------------------------------------------
500C L o c a l V a r i a b l e s
501C-----------------------------------------------
502 INTEGER I
503 my_real
504 . fac
505C
506 DO i=1,nel
507 fac = off(i)*vol(i)/volg(i)
508 sigg(i,1) = sigg(i,1) + fac * sig(i,1)
509 sigg(i,2) = sigg(i,2) + fac * sig(i,2)
510 sigg(i,3) = sigg(i,3) + fac * sig(i,3)
511 sigg(i,4) = sigg(i,4) + fac * sig(i,4)
512 sigg(i,5) = sigg(i,5) + fac * sig(i,5)
513 sigg(i,6) = sigg(i,6) + fac * sig(i,6)
514 rhog(i) = rhog(i) + fac * rho(i)
515 eintg(i) = eintg(i) + fac * eint(i)
516 dtxg(i) = min(dtxg(i),dtx(i))
517 ENDDO
518C
519 RETURN