OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
s8sforc3.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!|| s8sforc3 ../engine/source/elements/solid/solide8s/s8sforc3.F
25!||--- called by ------------------------------------------------------
26!|| forint ../engine/source/elements/forint.F
27!||--- calls -----------------------------------------------------
28!|| ancmsg ../engine/source/output/message/message.F
29!|| basis8 ../engine/source/elements/solid/solide8/basis8.F
30!|| csmall3 ../engine/source/elements/solid/solide/csmall3.F
31!|| degenes8 ../engine/source/elements/solid/solide/degenes8.F
32!|| mmain ../engine/source/materials/mat_share/mmain.F90
33!|| s8edefoc3 ../engine/source/elements/solid/solide8e/s8edefoc3.F
34!|| s8edefw3 ../engine/source/elements/solid/solide8e/s8edefw3.F
35!|| s8ederi_2 ../engine/source/elements/solid/solide8e/s8ederi_2.F
36!|| s8ederic3 ../engine/source/elements/solid/solide8e/s8ederic3.F
37!|| s8edericm3 ../engine/source/elements/solid/solide8e/s8edericm3.F
38!|| s8ederig3 ../engine/source/elements/solid/solide8e/s8ederig3.F
39!|| s8ederipr3 ../engine/source/elements/solid/solide8e/s8ederipr3.F
40!|| s8efmoy3 ../engine/source/elements/solid/solide8e/s8efmoy3.F
41!|| s8ejacip3 ../engine/source/elements/solid/solide8e/s8ejacip3.F
42!|| s8eprst_ini ../engine/source/elements/solid/solide8e/s8eprst_ini.F
43!|| s8sbdefo3 ../engine/source/elements/solid/solide8s/s8sbdefo3.F
44!|| s8sderi3 ../engine/source/elements/solid/solide8s/s8sderi3.F
45!|| s8sfint3 ../engine/source/elements/solid/solide8s/s8sfint3.F
46!|| s8sfint3_crimp ../engine/source/elements/solid/solide8s/s8sfint3_crimp.F
47!|| s8xref_imp ../engine/source/elements/solid/solide8s/s8xref_imp.F
48!|| s8zderims3 ../engine/source/elements/solid/solide8z/s8zderims3.F
49!|| s8zfint_reg ../engine/source/elements/solid/solide8z/s8zfint_reg.F
50!|| s8zfintp3 ../engine/source/elements/solid/solide8z/s8zfintp3.f
51!|| s8zzero3 ../engine/source/elements/solid/solide8z/s8zzero3.F
52!|| sbilan ../engine/source/elements/solid/solide/sbilan.F
53!|| scumu3 ../engine/source/elements/solid/solide/scumu3.F
54!|| scumu3p ../engine/source/elements/solid/solide/scumu3p.F
55!|| sdlen_dege ../engine/source/elements/solid/solide/sdlen_dege.F
56!|| sfillopt ../engine/source/elements/solid/solide/sfillopt.F
57!|| smallb3 ../engine/source/elements/solid/solide/smallb3.F
58!|| smallg3 ../engine/source/elements/solid/solide/smallg3.f
59!|| srbilan ../engine/source/elements/solid/solide/srbilan.F
60!|| srcoor3_imp ../engine/source/elements/solid/solide8s/srcoor3_imp.F
61!|| srho3 ../engine/source/elements/solid/solide/srho3.F
62!|| sstra3 ../engine/source/elements/solid/solide/sstra3.F
63!|| startime ../engine/source/system/timer_mod.F90
64!|| stoptime ../engine/source/system/timer_mod.F90
65!||--- uses -----------------------------------------------------
66!|| ale_connectivity_mod ../common_source/modules/ale/ale_connectivity_mod.F
67!|| dt_mod ../engine/source/modules/dt_mod.F
68!|| elbufdef_mod ../common_source/modules/mat_elem/elbufdef_mod.F90
69!|| glob_therm_mod ../common_source/modules/mat_elem/glob_therm_mod.F90
70!|| mat_elem_mod ../common_source/modules/mat_elem/mat_elem_mod.F90
71!|| message_mod ../engine/share/message_module/message_mod.F
72!|| mmain_mod ../engine/source/materials/mat_share/mmain.F90
73!|| nlocal_reg_mod ../common_source/modules/nlocal_reg_mod.F
74!|| output_mod ../common_source/modules/output/output_mod.F90
75!|| sensor_mod ../common_source/modules/sensor_mod.F90
76!|| table_mod ../engine/share/modules/table_mod.F
77!|| timer_mod ../engine/source/system/timer_mod.F90
78!||====================================================================
79 SUBROUTINE s8sforc3(TIMERS, OUTPUT, ELBUF_TAB,NG ,
80 1 PM ,GEO ,IXS ,X ,
81 2 A ,V ,MS ,W ,FLUX ,
82 3 FLU1 ,VEUL ,FV ,ALE_CONNECT ,IPARG ,
83 4 TF ,NPF ,BUFMAT ,PARTSAV ,NLOC_DMG,
84 5 DT2T ,NELTST ,ITYPTST ,STIFN ,FSKY ,
85 6 IADS ,OFFSET ,EANI ,IPARTS ,ICP ,
86 7 F11 ,F21 ,F31 ,F12 ,F22 ,
87 8 F32 ,F13 ,F23 ,F33 ,F14 ,
88 9 F24 ,F34 ,F15 ,F25 ,F35 ,
89 A F16 ,F26 ,F36 ,F17 ,F27 ,
90 B F37 ,F18 ,F28 ,F38 ,NEL ,
91 C NVC ,IPM ,ITASK ,ISTRAIN ,
92 D TEMP ,FTHE ,FTHESKY ,IEXPAN ,GRESAV ,
93 E GRTH ,IGRTH ,MSSA ,DMELS ,TABLE ,
94 F IGEO ,XDP ,VOLN ,CONDN ,CONDNSKY,
95 G D ,IOUTPRT,MAT_ELEM,H3D_STRAIN,SNPC,STF,
96 H SBUFMAT ,SVIS,NSVOIS,IDTMINS,IRESP,IDTMIN,MAXFUNC,
97 I USERL_AVAIL,DT ,GLOB_THERM,SENSORS)
98C-----------------------------------------------
99C M o d u l e s
100C-----------------------------------------------
101 USE timer_mod
102 USE output_mod, only : output_
103 USE mmain_mod
104 USE table_mod
105 USE mat_elem_mod
106 USE message_mod
109 USE elbufdef_mod
110 USE dt_mod
111 use glob_therm_mod
112 USE sensor_mod
113C-----------------------------------------------
114C I m p l i c i t T y p e s
115C-----------------------------------------------
116#include "implicit_f.inc"
117C-----------------------------------------------
118C G l o b a l P a r a m e t e r s
119C-----------------------------------------------
120#include "mvsiz_p.inc"
121C-----------------------------------------------
122C C o m m o n B l o c k s
123C-----------------------------------------------
124#include "com01_c.inc"
125#include "com04_c.inc"
126#include "com08_c.inc"
127#include "scr07_c.inc"
128#include "vect01_c.inc"
129#include "parit_c.inc"
130#include "param_c.inc"
131#include "timeri_c.inc"
132#include "impl1_c.inc"
133#include "scr17_c.inc"
134C-----------------------------------------------
135C D u m m y A r g u m e n t s
136C-----------------------------------------------
137 TYPE(timer_) , INTENT(INOUT) :: TIMERS
138 TYPE(OUTPUT_), INTENT(INOUT) :: OUTPUT
139 INTEGER, INTENT(IN) :: SNPC
140 INTEGER, INTENT(IN) :: STF
141 INTEGER, INTENT(IN) :: SBUFMAT
142 INTEGER, INTENT(IN) :: NSVOIS
143 INTEGER, INTENT(IN) :: IDTMINS
144 INTEGER ,INTENT(IN) :: IRESP
145 integer,dimension(102) :: IDTMIN
146 INTEGER ,INTENT(IN) :: MAXFUNC
147 INTEGER, INTENT(IN) :: USERL_AVAIL
148
149 INTEGER IXS(NIXS,*), IPARG(NPARG,NGROUP),NPF(*),
150 . IADS(8,*),IPARTS(*), IPM(NPROPMI,*),GRTH(*),IGRTH(*),IGEO(*),
151 . IOUTPRT
152C
153 INTEGER NELTST,ITYPTST,OFFSET,NEL,ICP,
154 . ICSIG, NVC ,ITASK, ISTRAIN, IEXPAN,NG,H3D_STRAIN
155
156 double precision
157 . xdp(3,*)
158 my_real
159 . dt2t
160C
161 my_real
162 . pm(npropm,*), geo(*), x(*), a(*), v(*), ms(*), w(*),
163 . flux(6,*),flu1(*), veul(*), fv(*), tf(*),
164 . partsav(*),stifn(*), fsky(*),eani(*),bufmat(*),
165 . f11(mvsiz),f21(mvsiz),f31(mvsiz),
166 . f12(mvsiz),f22(mvsiz),f32(mvsiz),
167 . f13(mvsiz),f23(mvsiz),f33(mvsiz),
168 . f14(mvsiz),f24(mvsiz),f34(mvsiz),
169 . f15(mvsiz),f25(mvsiz),f35(mvsiz),
170 . f16(mvsiz),f26(mvsiz),f36(mvsiz),
171 . f17(mvsiz),f27(mvsiz),f37(mvsiz),
172 . f18(mvsiz),f28(mvsiz),f38(mvsiz),
173 . temp(*),fthe(*),fthesky(*),gresav(*), mssa(*), dmels(*), voln(mvsiz),
174 . condn(*),condnsky(*),d(*)
175 my_real, DIMENSION(MVSIZ,6), INTENT(INOUT) :: svis
176 TYPE(ttable) TABLE(*)
177 TYPE (ELBUF_STRUCT_), TARGET, DIMENSION(NGROUP) :: ELBUF_TAB
178 TYPE (NLOCAL_STR_) , TARGET :: NLOC_DMG
179 TYPE(t_ale_connectivity), INTENT(IN) :: ALE_CONNECT
180 TYPE (MAT_ELEM_) ,INTENT(INOUT) :: MAT_ELEM
181 TYPE (DT_), INTENT(IN) :: DT
182 type (glob_therm_) ,intent(inout) :: glob_therm
183 type (sensors_),INTENT(INOUT) :: SENSORS
184C-----------------------------------------------
185C L o c a l V a r i a b l e s
186C-----------------------------------------------
187 INTEGER LCO, NF1, IFLAG, I,
188 . ILAY,IP,IR,IS,IT,MPT,NPTR,NPTS,NPTT,ICR,ICS,ICT,j
189 INTEGER IBID,IBIDV(1), MX, IFVM22
190 INTEGER MXT(MVSIZ),NGL(MVSIZ),NGEO(MVSIZ)
191 my_real
192 . VD2(MVSIZ) , DVOL(MVSIZ),DELTAX(MVSIZ),
193 . vis(mvsiz) , qvis(mvsiz), cxx(mvsiz) ,
194 . s1(mvsiz) , s2(mvsiz) , s3(mvsiz) ,
195 . s4(mvsiz) , s5(mvsiz) , s6(mvsiz) ,
196 . dxx(mvsiz) , dyy(mvsiz) , dzz(mvsiz) ,
197 . d4(mvsiz) , d5(mvsiz) , d6(mvsiz) ,
198 . ajc1(mvsiz) , ajc2(mvsiz) , ajc3(mvsiz) ,
199 . ajc4(mvsiz) , ajc5(mvsiz) , ajc6(mvsiz) ,
200 . ajc7(mvsiz) , ajc8(mvsiz) , ajc9(mvsiz) ,
201 . aj1(mvsiz) , aj2(mvsiz) , aj3(mvsiz) ,
202 . aj4(mvsiz) , aj5(mvsiz) , aj6(mvsiz) ,
203 . ajp1(mvsiz,8) , ajp2(mvsiz,8) , ajp3(mvsiz,8) ,
204 . ajp4(mvsiz,8) , ajp5(mvsiz,8) , ajp6(mvsiz,8) ,
205 . ajp7(mvsiz,8) , ajp8(mvsiz,8) , ajp9(mvsiz,8) ,
206 . vdx(mvsiz) , vdy(mvsiz) , vdz(mvsiz),ssp_eq(mvsiz),aire(mvsiz),
207 . e0(mvsiz),c1,fac(mvsiz) ,pr(8,8),ps(8,8),pt(8,8)
208C-----
209C Variables utilisees en argument par les materiaux.
210 my_real
211 . sti(mvsiz),stin(mvsiz),gama(mvsiz,6),
212 . wxx(mvsiz) , wyy(mvsiz) , wzz(mvsiz),
213 . wxx0(mvsiz) , wyy0(mvsiz) , wzz0(mvsiz),
214 . conde(mvsiz), conden(mvsiz),divde(mvsiz)
215C Variables utilisees en argument par les materiaux si SPH uniquement.
216 my_real
217 . muvoid(mvsiz)
218 INTEGER IOFFS, N,ITET
219 my_real
220 . OFFS(MVSIZ),DSV(MVSIZ),SDV(MVSIZ)
221 INTEGER SZ_R1_FREE
222C-----
223C Variables utilisees dans les routines solides uniquement (en arguments).
224 INTEGER NC1(MVSIZ), NC2(MVSIZ), NC3(MVSIZ), NC4(MVSIZ),
225 . NC5(MVSIZ), NC6(MVSIZ), NC7(MVSIZ), NC8(MVSIZ),
226 . MAT(MVSIZ)
227
228 DOUBLE PRECISION
229 . XD1(MVSIZ), XD2(MVSIZ), XD3(MVSIZ), XD4(MVSIZ),
230 . XD5(MVSIZ), XD6(MVSIZ), XD7(MVSIZ), XD8(MVSIZ),
231 . YD1(MVSIZ), YD2(MVSIZ), YD3(MVSIZ), YD4(MVSIZ),
232 . YD5(MVSIZ), YD6(MVSIZ), YD7(MVSIZ), YD8(MVSIZ),
233 . ZD1(MVSIZ), ZD2(MVSIZ), ZD3(MVSIZ), ZD4(MVSIZ),
234 . zd5(mvsiz), zd6(mvsiz), zd7(mvsiz), zd8(mvsiz),
235 . x0(mvsiz,8),y0(mvsiz,8),z0(mvsiz,8)
236
237 my_real
238 . off(mvsiz) ,offg(mvsiz) , rhoo(mvsiz),
239 . x1(mvsiz), x2(mvsiz), x3(mvsiz), x4(mvsiz),
240 . x5(mvsiz), x6(mvsiz), x7(mvsiz), x8(mvsiz),
241 . y1(mvsiz), y2(mvsiz), y3(mvsiz), y4(mvsiz),
242 . y5(mvsiz), y6(mvsiz), y7(mvsiz), y8(mvsiz),
243 . z1(mvsiz), z2(mvsiz), z3(mvsiz), z4(mvsiz),
244 . z5(mvsiz), z6(mvsiz), z7(mvsiz), z8(mvsiz),
245 . vx1(mvsiz),vx2(mvsiz),vx3(mvsiz),vx4(mvsiz),
246 . vx5(mvsiz),vx6(mvsiz),vx7(mvsiz),vx8(mvsiz),
247 . vy1(mvsiz),vy2(mvsiz),vy3(mvsiz),vy4(mvsiz),
248 . vy5(mvsiz),vy6(mvsiz),vy7(mvsiz),vy8(mvsiz),
249 . vz1(mvsiz),vz2(mvsiz),vz3(mvsiz),vz4(mvsiz),
250 . vz5(mvsiz),vz6(mvsiz),vz7(mvsiz),vz8(mvsiz),
251 . pxc1(mvsiz),pxc2(mvsiz),pxc3(mvsiz),pxc4(mvsiz),
252 . pyc1(mvsiz),pyc2(mvsiz),pyc3(mvsiz),pyc4(mvsiz),
253 . pzc1(mvsiz),pzc2(mvsiz),pzc3(mvsiz),pzc4(mvsiz),
254 . vdx1(mvsiz),vdx2(mvsiz),vdx3(mvsiz),vdx4(mvsiz),
255 . vdx5(mvsiz),vdx6(mvsiz),vdx7(mvsiz),vdx8(mvsiz),
256 . vdy1(mvsiz),vdy2(mvsiz),vdy3(mvsiz),vdy4(mvsiz),
257 . vdy5(mvsiz),vdy6(mvsiz),vdy7(mvsiz),vdy8(mvsiz),
258 . vdz1(mvsiz),vdz2(mvsiz),vdz3(mvsiz),vdz4(mvsiz),
259 . vdz5(mvsiz),vdz6(mvsiz),vdz7(mvsiz),vdz8(mvsiz),
260 . vgxa(mvsiz),vgya(mvsiz),vgza(mvsiz), vga2(mvsiz),
261 . hx(mvsiz,4), hy(mvsiz,4), hz(mvsiz,4),offg0(mvsiz),
262 . xgxa(mvsiz),xgya(mvsiz),xgza(mvsiz),
263 . xgxya(mvsiz),xgyza(mvsiz),xgzxa(mvsiz),
264 . xgxa2(mvsiz),xgya2(mvsiz),xgza2(mvsiz)
265 my_real
266 . dxy(mvsiz),dyx(mvsiz),
267 . dyz(mvsiz),dzy(mvsiz),
268 . dzx(mvsiz),dxz(mvsiz)
269 my_real
270 . r11(mvsiz),r12(mvsiz),r13(mvsiz),
271 . r21(mvsiz),r22(mvsiz),r23(mvsiz),
272 . r31(mvsiz),r32(mvsiz),r33(mvsiz)
273 my_real
274 . wi,smax(mvsiz),volg(mvsiz),pp(mvsiz),bid(mvsiz)
275 my_real
276 . sigy(mvsiz), et(mvsiz), nu(mvsiz),nu1(mvsiz),
277 . r1_free(mvsiz),r3_free(mvsiz),r4_free(mvsiz)
278 my_real
279 . vx0(mvsiz,8),vy0(mvsiz,8),vz0(mvsiz,8),
280 . tempel(mvsiz),them(mvsiz,8),die(mvsiz)
281 INTEGER NNPT,IDEG(MVSIZ)
282 PARAMETER (NNPT = 8)
283 my_real
284 . mfxx(mvsiz,nnpt),mfxy(mvsiz,nnpt),mfyx(mvsiz,nnpt),
285 . mfyy(mvsiz,nnpt),mfyz(mvsiz,nnpt),mfzy(mvsiz,nnpt),
286 . mfzz(mvsiz,nnpt),mfzx(mvsiz,nnpt),mfxz(mvsiz,nnpt),
287 . bxx(mvsiz,nnpt),byy(mvsiz,nnpt),bzz(mvsiz,nnpt),
288 . bxy(mvsiz,nnpt),byz(mvsiz,nnpt),bxz(mvsiz,nnpt),
289 . ni(8),detf0(mvsiz),jfac(mvsiz,nnpt)
290C-----
291 TYPE(g_bufel_) ,POINTER :: GBUF
292 TYPE(L_BUFEL_) ,POINTER :: LBUF
293c-----------------------------------------------------
294 my_real
295 . bxy1(mvsiz),bxy2(mvsiz),bxy3(mvsiz),bxy4(mvsiz),
296 . bxy5(mvsiz),bxy6(mvsiz),bxy7(mvsiz),bxy8(mvsiz),
297 . byx1(mvsiz),byx2(mvsiz),byx3(mvsiz),byx4(mvsiz),
298 . byx5(mvsiz),byx6(mvsiz),byx7(mvsiz),byx8(mvsiz),
299 . bxz1(mvsiz),bxz2(mvsiz),bxz3(mvsiz),bxz4(mvsiz),
300 . bxz5(mvsiz),bxz6(mvsiz),bxz7(mvsiz),bxz8(mvsiz),
301 . bzx1(mvsiz),bzx2(mvsiz),bzx3(mvsiz),bzx4(mvsiz),
302 . bzx5(mvsiz),bzx6(mvsiz),bzx7(mvsiz),bzx8(mvsiz),
303 . byz1(mvsiz),byz2(mvsiz),byz3(mvsiz),byz4(mvsiz),
304 . byz5(mvsiz),byz6(mvsiz),byz7(mvsiz),byz8(mvsiz),
305 . bzy1(mvsiz),bzy2(mvsiz),bzy3(mvsiz),bzy4(mvsiz),
306 . bzy5(mvsiz),bzy6(mvsiz),bzy7(mvsiz),bzy8(mvsiz),
307 . bxx1(mvsiz),bxx2(mvsiz),bxx3(mvsiz),bxx4(mvsiz),
308 . bxx5(mvsiz),bxx6(mvsiz),bxx7(mvsiz),bxx8(mvsiz),
309 . byy1(mvsiz),byy2(mvsiz),byy3(mvsiz),byy4(mvsiz),
310 . byy5(mvsiz),byy6(mvsiz),byy7(mvsiz),byy8(mvsiz),
311 . bzz1(mvsiz),bzz2(mvsiz),bzz3(mvsiz),bzz4(mvsiz),
312 . bzz5(mvsiz),bzz6(mvsiz),bzz7(mvsiz),bzz8(mvsiz)
313C
314 INTEGER NNEGA,INDEX(MVSIZ),ISELECT,iel,ipr
315 my_real
316 . PXY1(MVSIZ),PXY2(MVSIZ),PXY3(MVSIZ),PXY4(MVSIZ),
317 . PXY5(MVSIZ),PXY6(MVSIZ),PXY7(MVSIZ),PXY8(MVSIZ),
318 . PYX1(MVSIZ),PYX2(MVSIZ),PYX3(MVSIZ),PYX4(MVSIZ),
319 . PYX5(MVSIZ),PYX6(MVSIZ),PYX7(MVSIZ),PYX8(MVSIZ),
320 . PXZ1(MVSIZ),PXZ2(MVSIZ),PXZ3(MVSIZ),PXZ4(MVSIZ),
321 . PXZ5(MVSIZ),PXZ6(MVSIZ),PXZ7(MVSIZ),PXZ8(MVSIZ),
322 . PZX1(MVSIZ),PZX2(MVSIZ),PZX3(MVSIZ),PZX4(MVSIZ),
323 . PZX5(MVSIZ),PZX6(MVSIZ),PZX7(MVSIZ),PZX8(MVSIZ),
324 . PYZ1(MVSIZ),PYZ2(MVSIZ),PYZ3(MVSIZ),PYZ4(MVSIZ),
325 . PYZ5(MVSIZ),PYZ6(MVSIZ),PYZ7(MVSIZ),PYZ8(MVSIZ),
326 . PZY1(MVSIZ),PZY2(MVSIZ),PZY3(MVSIZ),PZY4(MVSIZ),
327 . PZY5(MVSIZ),PZY6(MVSIZ),PZY7(MVSIZ),PZY8(MVSIZ),
328 . AJI1(MVSIZ,NNPT) , AJI2(MVSIZ,NNPT) , AJI3(MVSIZ,NNPT) ,
329 . AJI4(MVSIZ,NNPT) , AJI5(MVSIZ,NNPT) , AJI6(MVSIZ,NNPT) ,
330 . AJI7(MVSIZ,NNPT) , AJI8(MVSIZ,NNPT) , AJI9(MVSIZ,NNPT),
331 . PX1(MVSIZ,NNPT),PX2(MVSIZ,NNPT),PX3(MVSIZ,NNPT),PX4(MVSIZ,NNPT),
332 . px5(mvsiz,nnpt),px6(mvsiz,nnpt),px7(mvsiz,nnpt),px8(mvsiz,nnpt),
333 . py1(mvsiz,nnpt),py2(mvsiz,nnpt),py3(mvsiz,nnpt),py4(mvsiz,nnpt),
334 . py5(mvsiz,nnpt),py6(mvsiz,nnpt),py7(mvsiz,nnpt),py8(mvsiz,nnpt),
335 . pz1(mvsiz,nnpt),pz2(mvsiz,nnpt),pz3(mvsiz,nnpt),pz4(mvsiz,nnpt),
336 . pz5(mvsiz,nnpt),pz6(mvsiz,nnpt),pz7(mvsiz,nnpt),pz8(mvsiz,nnpt),
337 . p0xy1(mvsiz,2),p0xy2(mvsiz,2),p0xy3(mvsiz,2),p0xy4(mvsiz,2),
338 . p0xy5(mvsiz,2),p0xy6(mvsiz,2),p0xy7(mvsiz,2),p0xy8(mvsiz,2),
339 . p0yx1(mvsiz,2),p0yx2(mvsiz,2),p0yx3(mvsiz,2),p0yx4(mvsiz,2),
340 . p0yx5(mvsiz,2),p0yx6(mvsiz,2),p0yx7(mvsiz,2),p0yx8(mvsiz,2),
341 . p0xz1(mvsiz,2),p0xz2(mvsiz,2),p0xz3(mvsiz,2),p0xz4(mvsiz,2),
342 . p0xz5(mvsiz,2),p0xz6(mvsiz,2),p0xz7(mvsiz,2),p0xz8(mvsiz,2),
343 . p0zx1(mvsiz,2),p0zx2(mvsiz,2),p0zx3(mvsiz,2),p0zx4(mvsiz,2),
344 . p0zx5(mvsiz,2),p0zx6(mvsiz,2),p0zx7(mvsiz,2),p0zx8(mvsiz,2),
345 . p0yz1(mvsiz,2),p0yz2(mvsiz,2),p0yz3(mvsiz,2),p0yz4(mvsiz,2),
346 . p0yz5(mvsiz,2),p0yz6(mvsiz,2),p0yz7(mvsiz,2),p0yz8(mvsiz,2),
347 . p0zy1(mvsiz,2),p0zy2(mvsiz,2),p0zy3(mvsiz,2),p0zy4(mvsiz,2),
348 . p0zy5(mvsiz,2),p0zy6(mvsiz,2),p0zy7(mvsiz,2),p0zy8(mvsiz,2),
349 . px1h1(mvsiz), px1h2(mvsiz), px1h3(mvsiz), px1h4(mvsiz),
350 . px2h1(mvsiz), px2h2(mvsiz), px2h3(mvsiz), px2h4(mvsiz),
351 . px3h1(mvsiz), px3h2(mvsiz), px3h3(mvsiz), px3h4(mvsiz),
352 . px4h1(mvsiz), px4h2(mvsiz), px4h3(mvsiz), px4h4(mvsiz),
353 . jr_1(mvsiz),js_1(mvsiz),jt_1(mvsiz),
354 . bb(6,24,mvsiz),amu(mvsiz)
355 double precision
356 . volp(mvsiz,nnpt)
357C----- Variables utilisees pour le calcul de puissance iteree
358C-----
359C Variables utilis es pour le non-local
360 my_real,
361 . DIMENSION(:,:), ALLOCATABLE :: var_reg
362 INTEGER :: INLOC, L_NLOC, INOD(8), IPOS(8), IMAT
363 my_real,
364 . DIMENSION(:), POINTER :: DNL
365 my_real
366 . H(8),PS2(8),PR2(8),PT2(8), ZR,ZS,ZT
367C
368 DOUBLE PRECISION
369 . ULX1(MVSIZ), ULX2(MVSIZ), ULX3(MVSIZ), ULX4(MVSIZ),
370 . ULX5(MVSIZ), ULX6(MVSIZ), ULX7(MVSIZ), ULX8(MVSIZ),
371 . ULY1(MVSIZ), ULY2(MVSIZ), ULY3(MVSIZ), ULY4(MVSIZ),
372 . ULY5(MVSIZ), ULY6(MVSIZ), ULY7(MVSIZ), ULY8(MVSIZ),
373 . ulz1(mvsiz), ulz2(mvsiz), ulz3(mvsiz), ulz4(mvsiz),
374 . ulz5(mvsiz), ulz6(mvsiz), ulz7(mvsiz), ulz8(mvsiz),
375 . dn_x(mvsiz,8),dn_y(mvsiz,8),dn_z(mvsiz,8) ,
376 . trm(nel,24,24),dn_r(8),dn_s(8),dn_t(8),invj(9,mvsiz)
377 my_real
378 . a11(mvsiz), a12(mvsiz), a13(mvsiz),
379 . a21(mvsiz), a22(mvsiz), a23(mvsiz),
380 . a31(mvsiz), a32(mvsiz), a33(mvsiz)
381 INTEGER SZ_IX
382C
383 my_real
384 . w_gauss(9,9),a_gauss(9,9)
385 DATA w_gauss /
386 1 2. ,0. ,0. ,
387 1 0. ,0. ,0. ,
388 1 0. ,0. ,0. ,
389 2 1. ,1. ,0. ,
390 2 0. ,0. ,0. ,
391 2 0. ,0. ,0. ,
392 3 0.555555555555556,0.888888888888889,0.555555555555556,
393 3 0. ,0. ,0. ,
394 3 0. ,0. ,0. ,
395 4 0.347854845137454,0.652145154862546,0.652145154862546,
396 4 0.347854845137454,0. ,0. ,
397 4 0. ,0. ,0. ,
398 5 0.236926885056189,0.478628670499366,0.568888888888889,
399 5 0.478628670499366,0.236926885056189,0. ,
400 5 0. ,0. ,0. ,
401 6 0.171324492379170,0.360761573048139,0.467913934572691,
402 6 0.467913934572691,0.360761573048139,0.171324492379170,
403 6 0. ,0. ,0. ,
404 7 0.129484966168870,0.279705391489277,0.381830050505119,
405 7 0.417959183673469,0.381830050505119,0.279705391489277,
406 7 0.129484966168870,0. ,0. ,
407 8 0.101228536290376,0.222381034453374,0.313706645877887,
408 8 0.362683783378362,0.362683783378362,0.313706645877887,
409 8 0.222381034453374,0.101228536290376,0. ,
410 9 0.081274388361574,0.180648160694857,0.260610696402935,
411 9 0.312347077040003,0.330239355001260,0.312347077040003,
412 9 0.260610696402935,0.180648160694857,0.081274388361574/
413 DATA a_gauss /
414 1 0. ,0. ,0. ,
415 1 0. ,0. ,0. ,
416 1 0. ,0. ,0. ,
417 2 -.577350269189626,0.577350269189626,0. ,
418 2 0. ,0. ,0. ,
419 2 0. ,0. ,0. ,
420 3 -.774596669241483,0. ,0.774596669241483,
421 3 0. ,0. ,0. ,
422 3 0. ,0. ,0. ,
423 4 -.861136311594053,-.339981043584856,0.339981043584856,
424 4 0.861136311594053,0. ,0. ,
425 4 0. ,0. ,0. ,
426 5 -.906179845938664,-.538469310105683,0. ,
427 5 0.538469310105683,0.906179845938664,0. ,
428 5 0. ,0. ,0. ,
429 6 -.932469514203152,-.661209386466265,-.238619186083197,
430 6 0.238619186083197,0.661209386466265,0.932469514203152,
431 6 0. ,0. ,0. ,
432 7 -.949107912342759,-.741531185599394,-.405845151377397,
433 7 0. ,0.405845151377397,0.741531185599394,
434 7 0.949107912342759,0. ,0. ,
435 8 -.960289856497536,-.796666477413627,-.525532409916329,
436 8 -.183434642495650,0.183434642495650,0.525532409916329,
437 8 0.796666477413627,0.960289856497536,0. ,
438 9 -.968160239507626,-.836031107326636,-.613371432700590,
439 9 -.324253423403809,0. ,0.324253423403809,
440 9 0.613371432700590,0.836031107326636,0.968160239507626/
441C
442c Flag Bolt Preloading
443 INTEGER IBOLTP,NBPRELD
444 my_real,
445 . DIMENSION(:), POINTER :: BPRELD
446C-----------------------------------------------
447C S o u r c e L i n e s
448C-----------------------------------------------
449 NPTR = elbuf_tab(ng)%NPTR
450 npts = elbuf_tab(ng)%NPTS
451 nptt = elbuf_tab(ng)%NPTT
452 gbuf => elbuf_tab(ng)%GBUF
453C
454 iboltp = iparg(72,ng)
455 inloc = iparg(78,ng)
456 ALLOCATE(var_reg(nel,nptr*npts*nptt))
457 nbpreld = gbuf%G_BPRELD
458 bpreld =>gbuf%BPRELD(1:nbpreld*nel)
459 tempel(:) = zero
460 bid(:) = zero
461 sz_r1_free=mvsiz
462 sz_ix=numelq+numels+nsvois
463c
464 ilay = 1
465 ibid = 0
466 ibidv = 0
467 nf1=nft+1
468C------------------------shear locking is alleviated with ANS for shear
469 iselect=0
470C-----------
471C GATHER NODAL VARIABLES AND COMPUTE INTRINSIC ROTATION.
472C-----------
473C GATHER NODAL VARIABLES AND COMPUTE INTRINSIC ROTATION.
474 CALL srcoor3_imp(x,ixs(1,nf1),v,w,gbuf%GAMA,gama,
475 . x1, x2, x3, x4, x5, x6, x7, x8,
476 . y1, y2, y3, y4, y5, y6, y7, y8,
477 . z1, z2, z3, z4, z5, z6, z7, z8,
478 . vx1, vx2, vx3, vx4, vx5, vx6, vx7, vx8,
479 . vy1, vy2, vy3, vy4, vy5, vy6, vy7, vy8,
480 . vz1, vz2, vz3, vz4, vz5, vz6, vz7, vz8,
481 . vd2,vis,gbuf%OFF,offg,gbuf%SMSTR,gbuf%RHO,rhoo, gbuf%COR_FR,
482 . nc1,nc2,nc3,nc4,nc5,nc6,nc7,nc8,
483 . ngl,mxt,ngeo,ioutprt, vgxa, vgya, vgza, vga2,
484 . xd1, xd2, xd3, xd4, xd5, xd6, xd7, xd8,
485 . yd1, yd2, yd3, yd4, yd5, yd6, yd7, yd8,
486 . zd1, zd2, zd3, zd4, zd5, zd6, zd7, zd8,
487 . xdp, x0 , y0 , z0 , nel, trm, gbuf%COR_XR,
488 . ulx1 ,ulx2 ,ulx3 ,ulx4 ,ulx5 ,ulx6 ,ulx7 ,ulx8 ,
489 . uly1 ,uly2 ,uly3 ,uly4 ,uly5 ,uly6 ,uly7 ,uly8 ,
490 . ulz1 ,ulz2 ,ulz3 ,ulz4 ,ulz5 ,ulz6 ,ulz7 ,ulz8 ,
491 . xgxa ,xgya ,xgza, xgxa2,xgya2,xgza2,xgxya,xgyza,xgzxa,iparg(1,ng))
492c---
493C-----------
494C GATHER NODAL VARIABLES FOR TOTAL STRAIN CASE.
495C-----------
496C-----------
497 ioffs=0
498 DO i=1,nel
499 offs(i)=ep20
500 deltax(i)=ep30
501 ideg(i)=0
502 ENDDO
503 IF(jthe < 0) them(1:nel,1:8) =zero
504
505 CALL s8eprst_ini(pr ,ps ,pt )
506C-------now [B](Pij) is always in global system w/ ISMSTR10
507C-----------------------------
508C----------JACOBIEN-inverse- in case of Vol<=0-------
509 IF (ismstr/=11) THEN
510 CALL s8ederic3(
511 1 offg, volg, ngl, xd1,
512 2 xd2, xd3, xd4, xd5,
513 3 xd6, xd7, xd8, yd1,
514 4 yd2, yd3, yd4, yd5,
515 5 yd6, yd7, yd8, zd1,
516 6 zd2, zd3, zd4, zd5,
517 7 zd6, zd7, zd8, pxc1,
518 8 pxc2, pxc3, pxc4, pyc1,
519 9 pyc2, pyc3, pyc4, pzc1,
520 a pzc2, pzc3, pzc4, px1h1,
521 b px1h2, px1h3, px1h4, px2h1,
522 c px2h2, px2h3, px2h4, px3h1,
523 d px3h2, px3h3, px3h4, px4h1,
524 e px4h2, px4h3, px4h4, hx,
525 f hy, hz, jr_1, js_1,
526 g jt_1, ajc1, ajc2, ajc3,
527 h ajc4, ajc5, ajc6, ajc7,
528 i ajc8, ajc9, smax, gbuf%SMSTR,
529 j gbuf%OFF, nel, ismstr, jlag)
530 CALL s8ejacip3(
531 1 hx, hy, hz, ajc1,
532 2 ajc2, ajc3, ajc4, ajc5,
533 3 ajc6, ajc7, ajc8, ajc9,
534 4 ajp1, ajp2, ajp3, ajp4,
535 5 ajp5, ajp6, ajp7, ajp8,
536 6 ajp9, nel)
537 IF (idts6>0) THEN
538 CALL sdlen_dege(
539 1 volg, deltax, x1, x2,
540 2 x3, x4, x5, x6,
541 3 x7, x8, y1, y2,
542 4 y3, y4, y5, y6,
543 5 y7, y8, z1, z2,
544 6 z3, z4, z5, z6,
545 7 z7, z8, ixs(1,nf1),ideg,
546 8 nel)
547 END IF
548 END IF !(ISMSTR/=11) THEN
549C----------JACOBIEN-inverse- in case of Vol<=0-------
550 nnega = 0
551C-----------Begin integrating points-----
552c
553 DO ir=1,nptr
554 DO is=1,npts
555 DO it=1,nptt
556 ip = ir + ( (is-1) + (it-1)*npts )*nptr
557 wi = w_gauss(ir,nptr)*w_gauss(is,npts)*w_gauss(it,nptt)
558C
559 CALL s8ederipr3(
560 1 gbuf%OFF, volp(1,ip),ngl, wi,
561 2 ajp1(1,ip),ajp2(1,ip),ajp3(1,ip),ajp4(1,ip),
562 3 ajp5(1,ip),ajp6(1,ip),ajp7(1,ip),ajp8(1,ip),
563 4 ajp9(1,ip),aji1(1,ip),aji2(1,ip),aji3(1,ip),
564 5 aji4(1,ip),aji5(1,ip),aji6(1,ip),aji7(1,ip),
565 6 aji8(1,ip),aji9(1,ip),nnega, index,
566 7 ip, nel)
567 ENDDO ! IT=1,NPTT
568 ENDDO ! IS=1,NPTS
569 ENDDO ! IR=1,NPTR
570C
571 IF (nnega>0) THEN
572 CALL s8edericm3(
573 1 volg, ngl, xd1, xd2,
574 2 xd3, xd4, xd5, xd6,
575 3 xd7, xd8, yd1, yd2,
576 4 yd3, yd4, yd5, yd6,
577 5 yd7, yd8, zd1, zd2,
578 6 zd3, zd4, zd5, zd6,
579 7 zd7, zd8, pxc1, pxc2,
580 8 pxc3, pxc4, pyc1, pyc2,
581 9 pyc3, pyc4, pzc1, pzc2,
582 a pzc3, pzc4, hx, hy,
583 b hz, ajc1, ajc2, ajc3,
584 c ajc4, ajc5, ajc6, ajc7,
585 d ajc8, ajc9, smax, gbuf%SMSTR,
586 e gbuf%OFF, nnega, index, nel,
587 f ismstr)
588 DO ir=1,nptr
589 DO is=1,npts
590 DO it=1,nptt
591 ip = ir + ( (is-1) + (it-1)*npts )*nptr
592 wi = w_gauss(ir,nptr)*w_gauss(is,npts)*w_gauss(it,nptt)
593C
594 CALL s8zderims3(volp(1,ip),
595 . a_gauss(ir,nptr),a_gauss(is,npts),a_gauss(it,nptt),wi,
596 . hx, hy, hz,
597 . ajc1,ajc2,ajc3,
598 . ajc4,ajc5,ajc6,
599 . ajc7,ajc8,ajc9,
600 . ajp1(1,ip),ajp2(1,ip),ajp3(1,ip),
601 . ajp4(1,ip),ajp5(1,ip),ajp6(1,ip),
602 . ajp7(1,ip),ajp8(1,ip),ajp9(1,ip),
603 . aji1(1,ip),aji2(1,ip),aji3(1,ip),
604 . aji4(1,ip),aji5(1,ip),aji6(1,ip),
605 . aji7(1,ip),aji8(1,ip),aji9(1,ip),nnega,index)
606 ENDDO ! IT=1,NPTT
607 ENDDO ! IS=1,NPTS
608 ENDDO ! IR=1,NPTR
609 IF (ineg_v ==0) THEN
610C--- /NEGAVOL/STOP
611 CALL ancmsg(msgid=280,anmode=aninfo)
612 mstop = 1
613 END IF !(INEG_V /=0) THEN
614 END IF
615C --------------------------
616C --- UPDATE REF CONFIGURATION (possible future change to small strain option)
617C --- ! Total strain option doesn't change the Ref CONF.
618C --------------------------
619C
620 IF (icp==1.OR.icp==2) THEN
621 CALL s8edefoc3(
622 1 pxc1, pxc2, pxc3, pxc4,
623 2 pyc1, pyc2, pyc3, pyc4,
624 3 pzc1, pzc2, pzc3, pzc4,
625 4 vx1, vx2, vx3, vx4,
626 5 vx5, vx6, vx7, vx8,
627 6 vy1, vy2, vy3, vy4,
628 7 vy5, vy6, vy7, vy8,
629 8 vz1, vz2, vz3, vz4,
630 9 vz5, vz6, vz7, vz8,
631 a dsv, nel)
632C-----don't do it w/ degenerated elm, but limit to law28 after QA
633 IF (idts6==0) CALL degenes8(
634 1 ixs(1,nf1),ideg, nel)
635 IF (icp==1.AND.mtn==28) THEN
636 DO i=1,nel
637 IF (ideg(i)>0) ideg(i) = ideg(i) + 10
638 ENDDO
639 END IF !(ICP==1.AND.MTN==28) THEN
640 ENDIF
641C-----------
642C INITIALIZATION---Before increment-----
643C-----------
644 CALL s8zzero3(
645 1 f11, f21, f31, f12,
646 2 f22, f32, f13, f23,
647 3 f33, f14, f24, f34,
648 4 f15, f25, f35, f16,
649 5 f26, f36, f17, f27,
650 6 f37, f18, f28, f38,
651 7 gbuf%SIG, gbuf%EINT, gbuf%RHO, gbuf%QVIS,
652 8 gbuf%PLA, gbuf%EPSD, stin, pp,
653 9 gbuf%G_PLA, gbuf%G_EPSD,iexpan, gbuf%EINTTH,
654 a nel, conden)
655C-------------------------------------------
656C COMPUTE AVERAGE TEMPERATURE IN ELEMENT (for output)
657C-------------------------------------------
658 IF(jthe < 0 ) THEN
659 DO i=1,nel
660 IF(gbuf%OFF(i) == zero) cycle
661 tempel(i) = one_over_8 *( temp(nc1(i)) + temp(nc2(i))
662 . + temp(nc3(i)) + temp(nc4(i))
663 . + temp(nc5(i)) + temp(nc6(i))
664 . + temp(nc7(i)) + temp(nc8(i)))
665 gbuf%TEMP(i) = tempel(i)
666 ENDDO
667 ENDIF
668c-------------------------------------------
669c COMPUTE Regularized non local variable in Gauss point
670c-------------------------------------------
671 IF (inloc > 0) THEN
672 l_nloc = nloc_dmg%L_NLOC
673 dnl => nloc_dmg%DNL(1:l_nloc) ! DNL = non local variable increment
674 DO ir=1,nptr
675 DO is=1,npts
676 DO it=1,nptt
677 zr = a_gauss(ir,nptr)
678 zs = a_gauss(is,npts)
679 zt = a_gauss(it,nptt)
680 ip = ir + ( (is-1) + (it-1)*npts )*nptr
681 CALL basis8 (zr,zs,zt,h,pr2,ps2,pt2)
682 DO i=1,nel
683 inod(1) = nloc_dmg%IDXI(nc1(i))
684 inod(2) = nloc_dmg%IDXI(nc2(i))
685 inod(3) = nloc_dmg%IDXI(nc3(i))
686 inod(4) = nloc_dmg%IDXI(nc4(i))
687 inod(5) = nloc_dmg%IDXI(nc5(i))
688 inod(6) = nloc_dmg%IDXI(nc6(i))
689 inod(7) = nloc_dmg%IDXI(nc7(i))
690 inod(8) = nloc_dmg%IDXI(nc8(i))
691 DO j = 1, 8
692 ipos(j) = nloc_dmg%POSI(inod(j))
693 ENDDO
694 var_reg(i,ip) = h(1)*dnl(ipos(1)) + h(2)*dnl(ipos(2)) + h(3)*dnl(ipos(3))
695 . + h(4)*dnl(ipos(4)) + h(5)*dnl(ipos(5)) + h(6)*dnl(ipos(6))
696 . + h(7)*dnl(ipos(7)) + h(8)*dnl(ipos(8))
697 ENDDO
698 ENDDO
699 ENDDO
700 ENDDO
701 ENDIF
702C---------[Bj] first---
703 DO ir=1,nptr
704 DO is=1,npts
705 DO it=1,nptt
706C-----------
707 ip = ir + ( (is-1) + (it-1)*npts )*nptr
708C
709 CALL s8ederig3(
710 1 px1(1,ip), px2(1,ip), px3(1,ip), px4(1,ip),
711 2 px5(1,ip), px6(1,ip), px7(1,ip), px8(1,ip),
712 3 py1(1,ip), py2(1,ip), py3(1,ip), py4(1,ip),
713 4 py5(1,ip), py6(1,ip), py7(1,ip), py8(1,ip),
714 5 pz1(1,ip), pz2(1,ip), pz3(1,ip), pz4(1,ip),
715 6 pz5(1,ip), pz6(1,ip), pz7(1,ip), pz8(1,ip),
716 7 aji1(1,ip),aji2(1,ip),aji3(1,ip),aji4(1,ip),
717 8 aji5(1,ip),aji6(1,ip),aji7(1,ip),aji8(1,ip),
718 9 aji9(1,ip),pr(1,ip), ps(1,ip), pt(1,ip),
719 a nel)
720 ENDDO ! IT=1,NPTT
721 ENDDO ! IS=1,NPTS
722 ENDDO ! IR=1,NPTR
723C-----------Begin integrating points-----
724 DO is=1,npts
725 DO ir=1,nptr
726 DO it=1,nptt
727 lbuf => elbuf_tab(ng)%BUFLY(ilay)%LBUF(ir,is,it)
728C-----avoid multi-print
729 IF (ioffs == 1)THEN
730 DO i=1,nel
731 IF (offs(i)<=two) lbuf%OFF(i)=offs(i)
732 ENDDO
733 END IF
734C-----------
735 ip = ir + ( (is-1) + (it-1)*npts )*nptr
736 wi = w_gauss(ir,nptr)*w_gauss(is,npts)*w_gauss(it,nptt)
737C
738 CALL s8ederi_2(
739 1 offg, off, voln, a_gauss(ir,nptr),
740 2 a_gauss(is,npts),a_gauss(it,nptt),wi, px1(1,ip),
741 3 px2(1,ip), px3(1,ip), px4(1,ip), px5(1,ip),
742 4 px6(1,ip), px7(1,ip), px8(1,ip), py1(1,ip),
743 5 py2(1,ip), py3(1,ip), py4(1,ip), py5(1,ip),
744 6 py6(1,ip), py7(1,ip), py8(1,ip), pz1(1,ip),
745 7 pz2(1,ip), pz3(1,ip), pz4(1,ip), pz5(1,ip),
746 8 pz6(1,ip), pz7(1,ip), pz8(1,ip), pxc1,
747 9 pxc2, pxc3, pxc4, pyc1,
748 a pyc2, pyc3, pyc4, pzc1,
749 b pzc2, pzc3, pzc4, bxy1,
750 c bxy2, bxy3, bxy4, bxy5,
751 d bxy6, bxy7, bxy8, byx1,
752 e byx2, byx3, byx4, byx5,
753 f byx6, byx7, byx8, bxz1,
754 g bxz2, bxz3, bxz4, bxz5,
755 h bxz6, bxz7, bxz8, bzx1,
756 i bzx2, bzx3, bzx4, bzx5,
757 j bzx6, bzx7, bzx8, byz1,
758 k byz2, byz3, byz4, byz5,
759 l byz6, byz7, byz8, bzy1,
760 m bzy2, bzy3, bzy4, bzy5,
761 n bzy6, bzy7, bzy8, bxx1,
762 o bxx2, bxx3, bxx4, bxx5,
763 p bxx6, bxx7, bxx8, byy1,
764 q byy2, byy3, byy4, byy5,
765 r byy6, byy7, byy8, bzz1,
766 s bzz2, bzz3, bzz4, bzz5,
767 t bzz6, bzz7, bzz8, ajp4(1,ip),
768 u ajp5(1,ip), ajp6(1,ip), ajp7(1,ip), ajp8(1,ip),
769 v ajp9(1,ip), aj1, aj2, aj3,
770 w aj4, aj5, aj6, smax,
771 x deltax, icp, ideg, nu,
772 y volp(1,ip), nel)
773C
774 CALL s8sderi3(
775 1 offg, off, volp(1,ip), ngl,
776 2 a_gauss(it,nptt),a_gauss(ir,nptr),a_gauss(is,npts),wi,
777 3 xd1, xd2, xd3, xd4,
778 4 xd5, xd6, xd7, xd8,
779 5 yd1, yd2, yd3, yd4,
780 6 yd5, yd6, yd7, yd8,
781 7 zd1, zd2, zd3, zd4,
782 8 zd5, zd6, zd7, zd8,
783 9 a11, a12, a13, a21,
784 a a22, a23, a31, a32,
785 b a33, dn_r, dn_s, dn_t,
786 c invj, dn_x, dn_y, dn_z,
787 d voln, nel)
788 CALL s8sbdefo3(
789 1 ulx1, ulx2, ulx3, ulx4,
790 2 ulx5, ulx6, ulx7, ulx8,
791 3 uly1, uly2, uly3, uly4,
792 4 uly5, uly6, uly7, uly8,
793 5 ulz1, ulz2, ulz3, ulz4,
794 6 ulz5, ulz6, ulz7, ulz8,
795 7 xd1, xd2, xd3, xd4,
796 8 xd5, xd6, xd7, xd8,
797 9 yd1, yd2, yd3, yd4,
798 a yd5, yd6, yd7, yd8,
799 b zd1, zd2, zd3, zd4,
800 c zd5, zd6, zd7, zd8,
801 d invj, a_gauss(it,nptt),a_gauss(ir,nptr),a_gauss(is,npts),
802 e a11, a12, a13, a21,
803 f a22, a23, a31, a32,
804 g a33, dn_r, dn_s, dn_t,
805 h bb, dxx, dxy, dxz,
806 i dyx, dyy, dyz, dzx,
807 j dzy, dzz, d4, d5,
808 k d6, wxx, wyy, wzz,
809 l lbuf%VOL, off, lbuf%EINT, gbuf%OFF,
810 m dsv, icp, fac, sdv,
811 n iselect, ideg, nel, ismstr)
812c
813 DO i=1,nel
814 rhoo(i) = lbuf%RHO(i)
815 ENDDO
816 divde(1:nel) = dt1*(dxx(1:nel)+ dyy(1:nel)+ dzz(1:nel))+sdv(1:nel)
817 CALL srho3(
818 1 pm, lbuf%VOL, lbuf%RHO, lbuf%EINT,
819 2 divde, flux(1,nf1),flu1(nf1), voln,
820 3 dvol, ngl, mxt, off,
821 4 0, gbuf%TAG22, volp(1,ip), lbuf%VOL0DP,
822 5 amu, gbuf%OFF, nel, mtn,
823 6 jale, ismstr, jeul, jlag)
824C
825C-----------------------------
826C EXTRACT STRESSES + SMALL STRAIN
827C-----------------------------
828 CALL csmall3(lbuf%SIG,s1,s2,s3,s4,s5,s6,
829 . gbuf%OFF,off,nel)
830C
831C for heat transfert
832C
833C------------------------------------------------------
834C CALCUL DES CONTRAINTES SUIVANT LOIS CONSTITUTIVES
835C TIME STEP EST CALCULE ICI DT2T
836C------------------------------------------------------
837 IF ((itask==0).AND.(imon_mat==1)) CALL startime(timers,35)
838c
839c
840 CALL mmain(timers,output,
841 1 elbuf_tab, ng, pm, geo,
842 2 ale_connect, ixs, iparg,
843 3 v, tf, npf, bufmat,
844 4 sti, x, dt2t, neltst,
845 5 ityptst, offset, nel, w,
846 6 off, ngeo, mxt, ngl,
847 7 voln, vd2, dvol, deltax,
848 8 vis, qvis, cxx, s1,
849 9 s2, s3, s4, s5,
850 a s6, dxx, dyy, dzz,
851 b d4, d5, d6, wxx,
852 c wyy, wzz, aj1, aj2,
853 d aj3, aj4, aj5, aj6,
854 e vdx, vdy, vdz, muvoid,
855 f ssp_eq, aire, sigy, et,
856 g r1_free, lbuf%PLA, r3_free, amu,
857 h mfxx(1,ip), mfxy(1,ip), mfxz(1,ip), mfyx(1,ip),
858 i mfyy(1,ip), mfyz(1,ip), mfzx(1,ip), mfzy(1,ip),
859 j mfzz(1,ip), ipm, gama, bid,
860 k bid, bid, bid, bid,
861 l bid, bid, istrain, tempel,
862 m die, iexpan, ilay, mssa,
863 n dmels, ir, is, it,
864 o table, bid, bid, bid,
865 p bid, iparg(1,ng), igeo, conde,
866 q itask, nloc_dmg, var_reg(1,ip),mat_elem,
867 r h3d_strain, jplasol, jsph, sz_r1_free,
868 s snpc, stf, sbufmat, glob_therm,
869 * svis, sz_ix, iresp,
870 t n2d, th_strain, ngroup, tt,
871 . dt1, ntable, numelq, nummat,
872 . numgeo, numnod, numels,
873 . idel7nok, idtmin, maxfunc,
874 . imon_mat, userl_avail, impl_s,
875 . idyna, dt, bid ,sensors)
876C
877 IF (istrain == 1)THEN
878 CALL sstra3(
879 1 dxx, dyy, dzz, d4,
880 2 d5, d6, lbuf%STRA,wxx,
881 3 wyy, wzz, off, nel,
882 4 jcvt)
883 ENDIF
884 IF ((itask==0).AND.(imon_mat==1)) CALL stoptime(timers,35)
885C----------------------------
886C INTERNAL FORCES
887C----------------------------
888 CALL s8sfint3(
889 1 lbuf%SIG, f11, f21, f31,
890 2 f12, f22, f32, f13,
891 3 f23, f33, f14, f24,
892 4 f34, f15, f25, f35,
893 5 f16, f26, f36, f17,
894 6 f27, f37, f18, f28,
895 7 f38, dn_x, dn_y, dn_z,
896 8 bb, voln, qvis, icp,
897 9 jfac(1,ip),nel, iselect, ideg,
898 a ismstr, svis)
899 CALL s8efmoy3(
900 1 lbuf%SIG, voln, qvis, pp,
901 2 lbuf%EINT, lbuf%RHO, lbuf%QVIS, lbuf%PLA,
902 3 lbuf%EPSD, gbuf%EPSD, gbuf%SIG, gbuf%EINT,
903 4 gbuf%RHO, gbuf%QVIS, gbuf%PLA, volg,
904 5 sti, stin, icp, off,
905 6 lbuf%VOL, gbuf%VOL, gbuf%G_PLA, gbuf%G_EPSD,
906 7 lbuf%EINTTH,gbuf%EINTTH,iexpan, nel,
907 8 conde, conden,svis,glob_therm%NODADT_THERM,
908 9 gbuf%WPLA, lbuf%WPLA, gbuf%G_WPLA)
909c-------------------------
910c Virtual internal forces of regularized non local ddl
911c--------------------------
912 IF (inloc > 0) THEN
913 imat = mxt(1)
914 zr = a_gauss(ir,nptr)
915 zs = a_gauss(is,npts)
916 zt = a_gauss(it,nptt)
917 CALL basis8 (zr,zs,zt,h,pr2,ps2,pt2)
918 CALL s8zfint_reg(
919 1 nloc_dmg, var_reg(1,ip),nel, gbuf%OFF,
920 2 voln, nc1, nc2, nc3,
921 3 nc4, nc5, nc6, nc7,
922 4 nc8, px1, px2, px3,
923 5 px4, px5, px6, px7,
924 6 px8, py1, py2, py3,
925 7 py4, py5, py6, py7,
926 8 py8, pz1, pz2, pz3,
927 9 pz4, pz5, pz6, pz7,
928 a pz8, imat, h, wi,
929 b ip, itask, dt2t, gbuf%VOL,
930 c nft)
931 ENDIF
932 DO i=1,nel
933 offg(i)=min(offg(i),off(i))
934 IF (lbuf%OFF(i)>one .AND. gbuf%OFF(i) == one) THEN
935 offs(i)=min(lbuf%OFF(i),offs(i))
936 ioffs =1
937 END IF
938 ENDDO
939C-------------------------
940c finite element heat transfert
941C--------------------------
942C
943 ENDDO ! IT=1,NPTT
944 ENDDO ! IS=1,NPTS
945 ENDDO ! IR=1,NPTR
946C-----------End integrating points-----
947c
948 IF (ioffs == 1) THEN
949 DO i=1,nel
950 IF (offs(i)<=two) gbuf%OFF(i)=offs(i)
951 END DO
952c-------------------
953 DO ir=1,nptr
954 DO is=1,npts
955 DO it=1,nptt
956 lbuf => elbuf_tab(ng)%BUFLY(ilay)%LBUF(ir,is,it)
957 ip = ir + ( (is-1) + (it-1)*npts )*nptr
958 DO i=1,nel
959 IF (gbuf%OFF(i) > one) lbuf%OFF(i)=gbuf%OFF(i)
960 END DO
961 END DO
962 END DO
963 END DO
964 END IF ! IOFFS == 1
965 IF(icp==1.AND.ismstr/=10.AND.ismstr/=12) THEN
966 CALL s8zfintp3(
967 1 pxc1, pxc2, pxc3, pxc4,
968 2 pyc1, pyc2, pyc3, pyc4,
969 3 pzc1, pzc2, pzc3, pzc4,
970 4 f11, f21, f31, f12,
971 5 f22, f32, f13, f23,
972 6 f33, f14, f24, f34,
973 7 f15, f25, f35, f16,
974 8 f26, f36, f17, f27,
975 9 f37, f18, f28, f38,
976 a volg, pp, ideg, nel)
977 ENDIF
978C-----------------------------
979C SMALL STRAIN
980C-----------------------------
981 CALL smallb3(
982 1 gbuf%OFF,offg, nel, ismstr)
983 IF (ismstr==11.OR.(jcvt == 0 .AND. ismstr>0)) THEN
984 CALL s8edefw3(
985 1 pxc1, pxc2, pxc3, pxc4,
986 2 pyc1, pyc2, pyc3, pyc4,
987 3 pzc1, pzc2, pzc3, pzc4,
988 4 vx1, vx2, vx3, vx4,
989 5 vx5, vx6, vx7, vx8,
990 6 vy1, vy2, vy3, vy4,
991 7 vy5, vy6, vy7, vy8,
992 8 vz1, vz2, vz3, vz4,
993 9 vz5, vz6, vz7, vz8,
994 a wxx0, wyy0, wzz0, nel)
995 CALL smallg3(
996 1 gbuf%SMSTR,offg, wxx0, wyy0,
997 2 wzz0, r11, r12, r13,
998 3 r21, r22, r23, r31,
999 4 r32, r33, nel, ismstr,
1000 5 jcvt)
1001 END IF
1002 IF (inconv==1.AND.tt>zero ) THEN
1003 CALL s8xref_imp(
1004 1 gbuf%OFF, gbuf%COR_XR,x1, x2,
1005 2 x3, x4, x5, x6,
1006 3 x7, x8, y1, y2,
1007 4 y3, y4, y5, y6,
1008 5 y7, y8, z1, z2,
1009 6 z3, z4, z5, z6,
1010 7 z7, z8, gbuf%COR_FR,nel)
1011 END IF
1012
1013C--------------------------
1014C BILANS PAR MATERIAU
1015C--------------------------
1016 iflag=mod(ncycle,ncpri)
1017 IF(ioutprt>0)THEN
1018 IF(jcvt == 0) THEN
1019 ifvm22 = 0
1020 CALL sbilan(partsav,gbuf%EINT,gbuf%RHO,gbuf%RK,gbuf%VOL,
1021 . vx1, vx2, vx3, vx4, vx5, vx6, vx7, vx8,
1022 . vy1, vy2, vy3, vy4, vy5, vy6, vy7, vy8,
1023 . vz1, vz2, vz3, vz4, vz5, vz6, vz7, vz8,
1024 . volg,iparts,gresav,grth,igrth,gbuf%OFF,
1025 . iexpan,gbuf%EINTTH,gbuf%FILL,gbuf%MOM,nel,ifvm22,
1026 . x1, x2, x3, x4, x5, x6, x7, x8,
1027 . y1, y2, y3, y4, y5, y6, y7, y8,
1028 . z1, z2, z3, z4, z5, z6, z7, z8,
1029 . itask,iparg(1,ng),sensors,gbuf%G_WPLA,gbuf%WPLA)
1030 ELSE
1031 CALL srbilan(partsav,gbuf%EINT,gbuf%RHO,gbuf%RK,gbuf%VOL,
1032 . vgxa, vgya, vgza, vga2, volg,iparts,
1033 . gresav,grth,igrth,gbuf%OFF,iexpan,gbuf%EINTTH,
1034 . gbuf%FILL, xgxa, xgya, xgza,
1035 . xgxa2,xgya2,xgza2,xgxya,xgyza,xgzxa,itask,iparg(1,ng),sensors,
1036 . nel,gbuf%G_WPLA,gbuf%WPLA)
1037 ENDIF
1038 ENDIF
1039C----------------------------
1040C CONVECTE --> GLOBAL.
1041C----------------------------
1042 CALL s8sfint3_crimp(trm,gbuf%COR_NF,gbuf%COR_FR,
1043 . f11,f21,f31,f12,f22,f32,f13,f23,f33,f14,f24,f34,
1044 . f15,f25,f35,f16,f26,f36,f17,f27,f37,f18,f28,f38,
1045 . nel )
1046C----------------------------
1047 IF(nfilsol/=0) CALL sfillopt(
1048 1 gbuf%FILL,stin, f11, f21,
1049 2 f31, f12, f22, f32,
1050 3 f13, f23, f33, f14,
1051 4 f24, f34, f15, f25,
1052 5 f35, f16, f26, f36,
1053 6 f17, f27, f37, f18,
1054 7 f28, f38, nel)
1055C----------------------------
1056 IF(iparit == 0)THEN
1057 CALL scumu3(
1058 1 gbuf%OFF,a, nc1, nc2,
1059 2 nc3, nc4, nc5, nc6,
1060 3 nc7, nc8, stifn, stin,
1061 4 f11, f21, f31, f12,
1062 5 f22, f32, f13, f23,
1063 6 f33, f14, f24, f34,
1064 7 f15, f25, f35, f16,
1065 8 f26, f36, f17, f27,
1066 9 f37, f18, f28, f38,
1067 a nvc, bid, bid, bid,
1068 b bid, bid, bid, bid,
1069 c bid, bid, bid, bid,
1070 d bid, bid, bid, bid,
1071 e bid, bid, bid, bid,
1072 f bid, bid, bid, bid,
1073 g bid, bid, bid, bid,
1074 h them, fthe, condn, conden,
1075 i nel, jthe, isrot, ipartsph,glob_therm%NODADT_THERM)
1076 ELSE
1077 CALL scumu3p(
1078 1 gbuf%OFF,stin, fsky, fsky,
1079 2 iads, f11, f21, f31,
1080 3 f12, f22, f32, f13,
1081 4 f23, f33, f14, f24,
1082 5 f34, f15, f25, f35,
1083 6 f16, f26, f36, f17,
1084 7 f27, f37, f18, f28,
1085 8 f38, nc1, nc2, nc3,
1086 9 nc4, nc5, nc6, nc7,
1087 a nc8, bid, bid, bid,
1088 b bid, bid, bid, bid,
1089 c bid, bid, bid, bid,
1090 d bid, bid, bid, bid,
1091 e bid, bid, bid, bid,
1092 f bid, bid, bid, bid,
1093 g bid, bid, bid, bid,
1094 h them, fthesky, condnsky,conden,
1095 i nel, nft, jthe, isrot,
1096 j ipartsph,glob_therm%NODADT_THERM)
1097 ENDIF
1098C-----------
1099 IF (ALLOCATED(var_reg)) DEALLOCATE(var_reg)
1100 RETURN
1101 END
1102
#define my_real
Definition cppsort.cpp:32
subroutine csmall3(sig, s1, s2, s3, s4, s5, s6, offg, off, nel)
Definition csmall3.F:35
subroutine degenes8(ixs, idege, nel)
Definition degenes8.F:37
#define min(a, b)
Definition macros.h:20
subroutine mmain(pm, elbuf_str, ix, nix, x, geo, iparg, nel, skew, bufmat, ipart, ipartel, nummat, matparam, imat, ipm, ngl, pid, npf, tf, mfxx, mfxy, mfxz, mfyx, mfyy, mfyz, mfzx, mfzy, mfzz, rx, ry, rz, sx, sy, sz, gama, voln, dvol, s1, s2, s3, s4, s5, s6, dxx, dyy, dzz, d4, d5, d6, wxx, wyy, wzz)
Definition mmain.F:43
subroutine s8edefoc3(pxc1, pxc2, pxc3, pxc4, pyc1, pyc2, pyc3, pyc4, pzc1, pzc2, pzc3, pzc4, vx1, vx2, vx3, vx4, vx5, vx6, vx7, vx8, vy1, vy2, vy3, vy4, vy5, vy6, vy7, vy8, vz1, vz2, vz3, vz4, vz5, vz6, vz7, vz8, dvca, nel)
Definition s8edefoc3.F:40
subroutine s8edefw3(px1, px2, px3, px4, py1, py2, py3, py4, pz1, pz2, pz3, pz4, vx1, vx2, vx3, vx4, vx5, vx6, vx7, vx8, vy1, vy2, vy3, vy4, vy5, vy6, vy7, vy8, vz1, vz2, vz3, vz4, vz5, vz6, vz7, vz8, wxx, wyy, wzz, nel)
Definition s8edefw3.F:41
subroutine s8ederi_2(offg, off, vol, ksi, eta, zeta, wi, px1, px2, px3, px4, px5, px6, px7, px8, py1, py2, py3, py4, py5, py6, py7, py8, pz1, pz2, pz3, pz4, pz5, pz6, pz7, pz8, pxc1, pxc2, pxc3, pxc4, pyc1, pyc2, pyc3, pyc4, pzc1, pzc2, pzc3, pzc4, bxy1, bxy2, bxy3, bxy4, bxy5, bxy6, bxy7, bxy8, byx1, byx2, byx3, byx4, byx5, byx6, byx7, byx8, bxz1, bxz2, bxz3, bxz4, bxz5, bxz6, bxz7, bxz8, bzx1, bzx2, bzx3, bzx4, bzx5, bzx6, bzx7, bzx8, byz1, byz2, byz3, byz4, byz5, byz6, byz7, byz8, bzy1, bzy2, bzy3, bzy4, bzy5, bzy6, bzy7, bzy8, bxx1, bxx2, bxx3, bxx4, bxx5, bxx6, bxx7, bxx8, byy1, byy2, byy3, byy4, byy5, byy6, byy7, byy8, bzz1, bzz2, bzz3, bzz4, bzz5, bzz6, bzz7, bzz8, aj4, aj5, aj6, aj7, aj8, aj9, rx, ry, rz, sx, sy, sz, smax, deltax, icp, ideg, nu, volp, nel)
Definition s8ederi_2.F:66
subroutine s8ederic3(off, det, ngl, x1, x2, x3, x4, x5, x6, x7, x8, y1, y2, y3, y4, y5, y6, y7, y8, z1, z2, z3, z4, z5, z6, z7, z8, px1, px2, px3, px4, py1, py2, py3, py4, pz1, pz2, pz3, pz4, px1h1, px1h2, px1h3, px1h4, px2h1, px2h2, px2h3, px2h4, px3h1, px3h2, px3h3, px3h4, px4h1, px4h2, px4h3, px4h4, hx, hy, hz, jr_1, js_1, jt_1, aj1, aj2, aj3, aj4, aj5, aj6, aj7, aj8, aj9, smax, sav, offg, nel, ismstr, jlag)
Definition s8ederic3.F:55
subroutine s8edericm3(det, ngl, x1, x2, x3, x4, x5, x6, x7, x8, y1, y2, y3, y4, y5, y6, y7, y8, z1, z2, z3, z4, z5, z6, z7, z8, px1, px2, px3, px4, py1, py2, py3, py4, pz1, pz2, pz3, pz4, hx, hy, hz, aj1, aj2, aj3, aj4, aj5, aj6, aj7, aj8, aj9, smax, sav, offg, nnega, index, nel, ismstr)
Definition s8edericm3.F:49
subroutine s8ederig3(px1, px2, px3, px4, px5, px6, px7, px8, py1, py2, py3, py4, py5, py6, py7, py8, pz1, pz2, pz3, pz4, pz5, pz6, pz7, pz8, aji1, aji2, aji3, aji4, aji5, aji6, aji7, aji8, aji9, pr, ps, pt, nel)
Definition s8ederig3.F:40
subroutine s8ederipr3(offg, voldp, ngl, wi, aj1, aj2, aj3, aj4, aj5, aj6, aj7, aj8, aj9, aji1, aji2, aji3, aji4, aji5, aji6, aji7, aji8, aji9, nnega, index, ipt, nel)
Definition s8ederipr3.F:39
subroutine s8efmoy3(sigor, vol, qvis, pp, eint, rho, q, defp, epsd, epsdm, sigm, eintm, rhom, qm, defpm, volg, sti, stin, icp, off, vol0, vol0g, g_pla, g_epsd, eintth, eintthm, iexpan, nel, conde, conden, svis, nodadt_therm, g_wpla, l_wpla, g_wpla_flag)
Definition s8efmoy3.F:40
subroutine s8eprst_ini(pr, ps, pt)
Definition s8eprst_ini.F:30
subroutine s8sbdefo3(ulx1, ulx2, ulx3, ulx4, ulx5, ulx6, ulx7, ulx8, uly1, uly2, uly3, uly4, uly5, uly6, uly7, uly8, ulz1, ulz2, ulz3, ulz4, ulz5, ulz6, ulz7, ulz8, x1, x2, x3, x4, x5, x6, x7, x8, y1, y2, y3, y4, y5, y6, y7, y8, z1, z2, z3, z4, z5, z6, z7, z8, invj, ksi, eta, zeta, a11, a12, a13, a21, a22, a23, a31, a32, a33, dn_r, dn_s, dn_t, bb, dxx, dxy, dxz, dyx, dyy, dyz, dzx, dzy, dzz, d4, d5, d6, wxx, wyy, wzz, volo, off, eint, offs, dsv, icp, fac, sdv, i_sh, idege, nel, ismstr)
Definition s8sbdefo3.F:52
subroutine s8sderi3(offg, off, voldp, ngl, ksi, eta, zeta, wi, x1, x2, x3, x4, x5, x6, x7, x8, y1, y2, y3, y4, y5, y6, y7, y8, z1, z2, z3, z4, z5, z6, z7, z8, a11, a12, a13, a21, a22, a23, a31, a32, a33, dn_r, dn_s, dn_t, invj, dn_x, dn_y, dn_z, voln, nel)
Definition s8sderi3.F:43
subroutine s8sfint3(sig, f11, f21, f31, f12, f22, f32, f13, f23, f33, f14, f24, f34, f15, f25, f35, f16, f26, f36, f17, f27, f37, f18, f28, f38, dn_x, dn_y, dn_z, bb, vol, qvis, icp, jfac, nel, i_sh, idege, ismstr, svis)
Definition s8sfint3.F:39
subroutine s8sfint3_crimp(trm, qf, r, f11, f21, f31, f12, f22, f32, f13, f23, f33, f14, f24, f34, f15, f25, f35, f16, f26, f36, f17, f27, f37, f18, f28, f38, nel)
subroutine s8sforc3(timers, output, elbuf_tab, ng, pm, geo, ixs, x, a, v, ms, w, flux, flu1, veul, fv, ale_connect, iparg, tf, npf, bufmat, partsav, nloc_dmg, dt2t, neltst, ityptst, stifn, fsky, iads, offset, eani, iparts, icp, f11, f21, f31, f12, f22, f32, f13, f23, f33, f14, f24, f34, f15, f25, f35, f16, f26, f36, f17, f27, f37, f18, f28, f38, nel, nvc, ipm, itask, istrain, temp, fthe, fthesky, iexpan, gresav, grth, igrth, mssa, dmels, table, igeo, xdp, voln, condn, condnsky, d, ioutprt, mat_elem, h3d_strain, snpc, stf, sbufmat, svis, nsvois, idtmins, iresp, idtmin, maxfunc, userl_avail, dt, glob_therm, sensors)
Definition s8sforc3.F:98
subroutine s8xref_imp(offg, xref, xd1, xd2, xd3, xd4, xd5, xd6, xd7, xd8, yd1, yd2, yd3, yd4, yd5, yd6, yd7, yd8, zd1, zd2, zd3, zd4, zd5, zd6, zd7, zd8, r, nel)
Definition s8xref_imp.F:38
subroutine s8zderims3(voldp, ksi, eta, zeta, wi, hx, hy, hz, cj1, cj2, cj3, cj4, cj5, cj6, cj7, cj8, cj9, jac1, jac2, jac3, jac4, jac5, jac6, jac7, jac8, jac9, jaci1, jaci2, jaci3, jaci4, jaci5, jaci6, jaci7, jaci8, jaci9, nnega, index)
Definition s8zderims3.F:42
subroutine s8zfint_reg(nloc_dmg, var_reg, nel, offg, vol, nc1, nc2, nc3, nc4, nc5, nc6, nc7, nc8, px1, px2, px3, px4, px5, px6, px7, px8, py1, py2, py3, py4, py5, py6, py7, py8, pz1, pz2, pz3, pz4, pz5, pz6, pz7, pz8, imat, h, wi, ip, itask, dt2t, vol0, nft)
Definition s8zfint_reg.F:45
subroutine s8zfintp3(px1, px2, px3, px4, py1, py2, py3, py4, pz1, pz2, pz3, pz4, f11, f21, f31, f12, f22, f32, f13, f23, f33, f14, f24, f34, f15, f25, f35, f16, f26, f36, f17, f27, f37, f18, f28, f38, vol, pp, idege, nel)
Definition s8zfintp3.F:42
subroutine s8zzero3(fx1, fy1, fz1, fx2, fy2, fz2, fx3, fy3, fz3, fx4, fy4, fz4, fx5, fy5, fz5, fx6, fy6, fz6, fx7, fy7, fz7, fx8, fy8, fz8, sigm, eintm, rhom, qm, defpm, epsdm, stin, pp, g_pla, g_epsd, iexpan, eintthm, nel, conden)
Definition s8zzero3.F:42
subroutine sbilan(partsav, eint, rho, rk, vol, vx1, vx2, vx3, vx4, vx5, vx6, vx7, vx8, vy1, vy2, vy3, vy4, vy5, vy6, vy7, vy8, vz1, vz2, vz3, vz4, vz5, vz6, vz7, vz8, vnew, iparts, gresav, grth, igrth, off, iexpan, eintth, fill, mom, nel, ifvm22, x1, x2, x3, x4, x5, x6, x7, x8, y1, y2, y3, y4, y5, y6, y7, y8, z1, z2, z3, z4, z5, z6, z7, z8, itask, iparg, sensors, g_wpla, wpla)
Definition sbilan.F:47
subroutine scumu3(offg, e, nc1, nc2, nc3, nc4, nc5, nc6, nc7, nc8, stifn, sti, f11, f21, f31, f12, f22, f32, f13, f23, f33, f14, f24, f34, f15, f25, f35, f16, f26, f36, f17, f27, f37, f18, f28, f38, nvc, ar, fr_wave, fr_wav, mx1, my1, mz1, mx2, my2, mz2, mx3, my3, mz3, mx4, my4, mz4, mx5, my5, mz5, mx6, my6, mz6, mx7, my7, mz7, mx8, my8, mz8, them, fthe, condn, conde, nel, jthe, isrot, ipartsph, nodadt_therm)
Definition scumu3.F:55
subroutine scumu3p(offg, sti, fsky, fskyv, iads, f11, f21, f31, f12, f22, f32, f13, f23, f33, f14, f24, f34, f15, f25, f35, f16, f26, f36, f17, f27, f37, f18, f28, f38, nc1, nc2, nc3, nc4, nc5, nc6, nc7, nc8, ar, fr_wave, fr_wav, mx1, my1, mz1, mx2, my2, mz2, mx3, my3, mz3, mx4, my4, mz4, mx5, my5, mz5, mx6, my6, mz6, mx7, my7, mz7, mx8, my8, mz8, them, fthesky, condnsky, conde, nel, nft, jthe, isrot, ipartsph, nodadt_therm)
Definition scumu3p.F:56
subroutine sdlen_dege(volg, lat, x1, x2, x3, x4, x5, x6, x7, x8, y1, y2, y3, y4, y5, y6, y7, y8, z1, z2, z3, z4, z5, z6, z7, z8, ixs, idege, nel)
Definition sdlen_dege.F:45
subroutine sfillopt(fill, sti, f11, f21, f31, f12, f22, f32, f13, f23, f33, f14, f24, f34, f15, f25, f35, f16, f26, f36, f17, f27, f37, f18, f28, f38, nel)
Definition sfillopt.F:43
subroutine smallb3(offg, off, nel, ismstr)
Definition smallb3.F:44
subroutine smallg3(sav, offg, wxx, wyy, wzz, r11, r12, r13, r21, r22, r23, r31, r32, r33, nel, ismstr, jcvt)
Definition smallg3.F:36
subroutine srbilan(partsav, eint, rho, rk, vol, vxa, vya, vza, va2, vnew, iparts, gresav, grth, igrth, off, iexpan, eintth, fill, xx, yy, zz, xx2, yy2, zz2, xy, yz, zx, itask, iparg, sensors, nel, g_wpla, wpla)
Definition srbilan.F:45
subroutine srcoor3_imp(x, ixs, v, w, gama0, gama, x1, x2, x3, x4, x5, x6, x7, x8, y1, y2, y3, y4, y5, y6, y7, y8, z1, z2, z3, z4, z5, z6, z7, z8, vx1, vx2, vx3, vx4, vx5, vx6, vx7, vx8, vy1, vy2, vy3, vy4, vy5, vy6, vy7, vy8, vz1, vz2, vz3, vz4, vz5, vz6, vz7, vz8, vd2, vis, offg, off, sav, rho, rhoo, r, nc1, nc2, nc3, nc4, nc5, nc6, nc7, nc8, ngl, mxt, ngeo, ioutprt, vgax, vgay, vgaz, vga2, xd1, xd2, xd3, xd4, xd5, xd6, xd7, xd8, yd1, yd2, yd3, yd4, yd5, yd6, yd7, yd8, zd1, zd2, zd3, zd4, zd5, zd6, zd7, zd8, xdp, x0, y0, z0, nel, trm, xref, ulx1, ulx2, ulx3, ulx4, ulx5, ulx6, ulx7, ulx8, uly1, uly2, uly3, uly4, uly5, uly6, uly7, uly8, ulz1, ulz2, ulz3, ulz4, ulz5, ulz6, ulz7, ulz8, xgax, xgay, xgaz, xgxa2, xgya2, xgza2, xgxya, xgyza, xgzxa, iparg)
Definition srcoor3_imp.F:53
subroutine sstra3(dxx, dyy, dzz, d4, d5, d6, strain, wxx, wyy, wzz, off, nel, jcvt)
Definition sstra3.F:46
subroutine basis8(r, s, t, h, pr, ps, pt)
Definition basis8.F:30
subroutine s8ejacip3(hx, hy, hz, cj1, cj2, cj3, cj4, cj5, cj6, cj7, cj8, cj9, jac1, jac2, jac3, jac4, jac5, jac6, jac7, jac8, jac9)
Definition s8zderi3.F:2180
subroutine srho3(pm, volo, rhon, eint, dxx, dyy, dzz, voln, dvol, mat)
Definition srho3.F:31
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 startime(event, itask)
Definition timer.F:93
subroutine stoptime(event, itask)
Definition timer.F:135