OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
hm_read_fail_biquad.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!|| hm_read_fail_biquad ../starter/source/materials/fail/biquad/hm_read_fail_biquad.F
25!||--- called by ------------------------------------------------------
26!|| hm_read_fail ../starter/source/materials/fail/hm_read_fail.F
27!||--- calls -----------------------------------------------------
28!|| ancmsg ../starter/source/output/message/message.F
29!|| biquad_coefficients ../starter/source/materials/fail/biquad/biquad_coefficients.F
30!|| hm_get_floatv ../starter/source/devtools/hm_reader/hm_get_floatv.F
31!|| hm_get_floatv_dim ../starter/source/devtools/hm_reader/hm_get_floatv_dim.F
32!|| hm_get_intv ../starter/source/devtools/hm_reader/hm_get_intv.F
33!|| hm_option_is_encrypted ../starter/source/devtools/hm_reader/hm_option_is_encrypted.F
34!||--- uses -----------------------------------------------------
35!|| hm_option_read_mod ../starter/share/modules1/hm_option_read_mod.F
36!|| message_mod ../starter/share/message_module/message_mod.F
37!|| submodel_mod ../starter/share/modules1/submodel_mod.F
38!||====================================================================
40 . FAIL ,MAT_ID ,FAIL_ID ,IRUPT ,
41 . TITR ,LSUBMODEL,UNITAB )
42C-----------------------------------------------
43c ROUTINE DESCRIPTION :
44c Read read double parabolic analytical failure model parameters
45c-----------------------------------------------
46c UVAR1 = damage due to instability (triax between 1/3 and 2/3)
47c UVAR2 = integration point
48c UVAR3-8 = perturbated parameter
49c UVAR3 (if perturbation is not used) or UVAR9 (if used) = initial element length
50C-----------------------------------------------
51C M o d u l e s
52C-----------------------------------------------
53 USE fail_param_mod
54 USE unitab_mod
55 USE message_mod
56 USE submodel_mod
59C-----------------------------------------------
60C I m p l i c i t T y p e s
61C-----------------------------------------------
62#include "implicit_f.inc"
63C-----------------------------------------------
64C C o m m o n B l o c k s
65C-----------------------------------------------
66#include "com04_c.inc"
67#include "units_c.inc"
68C-----------------------------------------------
69C D u m m y A r g u m e n t s
70C-----------------------------------------------
71 INTEGER ,INTENT(IN) :: FAIL_ID ! failure model ID
72 INTEGER ,INTENT(IN) :: MAT_ID ! material law ID
73 INTEGER ,INTENT(IN) :: IRUPT ! failure model type number
74 CHARACTER(LEN=NCHARTITLE) ,INTENT(IN) :: TITR ! material model title
75 TYPE(unit_type_) ,INTENT(IN) :: UNITAB ! table of input units
76 TYPE(submodel_data),INTENT(IN) :: LSUBMODEL(*) ! submodel table
77 TYPE(fail_param_) ,INTENT(INOUT) :: FAIL ! failure model data structure
78C-----------------------------------------------
79C L o c a l V a r i a b l e s
80C-----------------------------------------------
81 INTEGER :: MFLAG,SFLAG,REG_FUNC,ICOUP,NFUNC,NUVAR,NUPARAM,FAILIP
82 my_real :: C1,C2,C3,C4,C5,E1,E2,E3,E4,PTHK,INST,REF_LEN,REF_SIZ_UNIT
83 my_real :: x_1(2)
84 my_real :: x_2(3)
85 my_real :: xmin,ymin,dcrit,exp
86 LOGICAL :: IS_AVAILABLE,IS_ENCRYPTED
87C=======================================================================
88 is_encrypted = .false.
89 is_available = .false.
90C--------------------------------------------------
91C EXTRACT DATA (IS OPTION CRYPTED)
92C--------------------------------------------------
93 CALL hm_option_is_encrypted(is_encrypted)
94C--------------------------------------------------
95C EXTRACT INPUT DATA
96C--------------------------------------------------
97card1
98 CALL hm_get_floatv ('C1' ,c1 ,is_available,lsubmodel,unitab)
99 CALL hm_get_floatv ('C2' ,c2 ,is_available,lsubmodel,unitab)
100 CALL hm_get_floatv ('C3' ,c3 ,is_available,lsubmodel,unitab)
101 CALL hm_get_floatv ('C4' ,c4 ,is_available,lsubmodel,unitab)
102 CALL hm_get_floatv ('C5' ,c5 ,is_available,lsubmodel,unitab)
103c
104card2 damage accumulation parametars
105c
106 CALL hm_get_intv ('FAILIP' ,failip ,is_available,lsubmodel)
107 IF (failip == 0) failip = 1
108 CALL hm_get_floatv ('P_thickfail' ,pthk ,is_available,lsubmodel,unitab)
109 CALL hm_get_intv ('M_Flag' ,mflag ,is_available,lsubmodel)
110 CALL hm_get_intv ('S_Flag' ,sflag ,is_available,lsubmodel)
111 CALL hm_get_floatv ('Inst_start' ,inst ,is_available,lsubmodel,unitab)
112 CALL hm_get_intv ('fct_IDel' ,reg_func ,is_available,lsubmodel)
113 CALL hm_get_floatv ('EI_ref' ,ref_len ,is_available,lsubmodel,unitab)
114 IF (reg_func > 0 .AND. ref_len == zero) THEN
115 CALL hm_get_floatv_dim('EI_ref' ,ref_siz_unit,is_available, lsubmodel, unitab)
116 ref_len = one*ref_siz_unit
117 ENDIF
118c---------------------------------------------------
119c Optional input
120c---------------------------------------------------
121 CALL hm_get_floatv ('R1' ,e1 ,is_available,lsubmodel,unitab)
122 CALL hm_get_floatv ('R2' ,e2 ,is_available,lsubmodel,unitab)
123 CALL hm_get_floatv ('R4' ,e3 ,is_available,lsubmodel,unitab)
124 CALL hm_get_floatv ('R5' ,e4 ,is_available,lsubmodel,unitab)
125c---------------------------------------------------
126c
127c---------------------------------------------------
128c Stress softening input
129c---------------------------------------------------
130 CALL hm_get_intv ('ICOUP' ,icoup ,is_available,lsubmodel)
131 CALL hm_get_floatv ('DCRIT' ,dcrit ,is_available,lsubmodel,unitab)
132 CALL hm_get_floatv ('EXP' ,exp ,is_available,lsubmodel,unitab)
133c---------------------------------------------------
134c
135 ! Check Pthickfail parameter
136 pthk = min(pthk, one)
137 pthk = max(pthk,-one)
138 IF (pthk == zero) pthk = em06
139c
140 ! Check SFLAG formulation parameter
141 IF (sflag == 0) sflag = 2
142c
143c---------------------------------------------------
144c pre definition for user-input data when only
145c tension test data are provided
146c---------------------------------------------------
147 IF (c3 == zero) THEN
148 SELECT CASE (mflag)
149c
150 CASE (1) ! Mild Seel
151 c3 = 0.6
152 CASE (2) ! HSS Seel light e-Body DP600
153 c3 = 0.5
154 CASE (3) ! UHSS Seel light_eBody Boron
155 c3 = 0.12
156 CASE (4) ! Aluminum light_eBody AA5182
157 c3 = 0.3
158 CASE (5) ! Aluminum light_eBody AA6082-T6
159 c3 = 0.17
160 CASE (6) ! Plastic light_eBody PA6GF30
161 c3 = 0.1
162 CASE (7) ! Plastic light_eBody PP T40
163 c3 = 0.11
164 CASE DEFAULT
165 c3 = .6
166c
167 END SELECT
168 ENDIF
169c---------------------------------------------------
170c
171 CALL biquad_coefficients(c1,c2,c3,c4,c5,mflag,x_1,x_2,e1,e2,e3,e4)
172c
173 ! Check necking instability parameters
174 IF (sflag == 3 .AND. inst <= zero) THEN
175 CALL ancmsg(msgid=3042, msgtype=msgwarning, anmode=aninfo_blind_1,
176 . i1=mat_id,
177 . c1=titr)
178 sflag = 2
179 ELSEIF (sflag == 3 .AND. inst >= c4) THEN
180 CALL ancmsg(msgid=3043, msgtype=msgwarning, anmode=aninfo_blind_1,
181 . i1=mat_id,
182 . c1=titr)
183 sflag = 2
184 ENDIF
185c
186 ! Check stress softening parameters
187 dcrit = min(dcrit,one)
188 dcrit = max(dcrit,zero)
189 exp = abs(exp)
190 IF (exp == zero) exp = one
191 IF (dcrit /= zero .AND. icoup == 0) icoup = 1
192 IF (sflag /= 3 .AND. icoup == 2) THEN
193 CALL ancmsg(msgid=3044, msgtype=msgwarning, anmode=aninfo_blind_1,
194 . i1=mat_id,
195 . c1=titr)
196 icoup = 0
197 ENDIF
198c
199 ! Check if minimum of first parabola is negative
200 IF(x_1(2) .ne. zero) then
201 xmin = -x_1(1)/(two*x_1(2))
202 ymin = x_1(2)*(xmin**2) + x_1(1)*xmin + c2
203 IF (ymin < zero) THEN
204 CALL ancmsg(msgid=3004, msgtype=msgwarning, anmode=aninfo_blind_1,
205 . i1=mat_id,
206 . c1=titr)
207 ENDIF
208 ENDIF
209 ! Check if minimum of second parabola is negative
210 IF (sflag == 1 .and. x_2(3) .ne. 0) THEN
211 xmin = -x_2(2)/(two*x_2(3))
212 ymin = x_2(3)*(xmin**2) + x_2(2)*xmin + x_2(1)
213 IF (ymin < zero) THEN
214 CALL ancmsg(msgid=3005, msgtype=msgwarning, anmode=aninfo_blind_1,
215 . i1=mat_id,
216 . c1=titr)
217 ENDIF
218 ENDIF
219c
220c---------------------------------------------------
221 nuparam = 17
222 IF (reg_func == 0) THEN
223 nfunc = 0
224 nuvar = 2
225 IF (nperturb /= 0) nuvar = 8
226 ELSE
227 nfunc = 1
228 nuvar = 3
229 IF (nperturb /= 0) nuvar = 9
230 ENDIF
231c-------------------------
232 fail%KEYWORD = 'BIQUAD'
233 fail%IRUPT = irupt
234 fail%FAIL_ID = fail_id
235 fail%NUPARAM = nuparam
236 fail%NIPARAM = 0
237 fail%NUVAR = nuvar
238 fail%NFUNC = nfunc
239 fail%NTABLE = 0
240 fail%NMOD = 0
241c
242 fail%PTHK = pthk
243c
244 ALLOCATE (fail%UPARAM(fail%NUPARAM))
245 ALLOCATE (fail%IPARAM(fail%NIPARAM))
246 ALLOCATE (fail%IFUNC (fail%NFUNC))
247 ALLOCATE (fail%TABLE (fail%NTABLE))
248
249 IF (nfunc == 1) fail%IFUNC(1) = reg_func
250
251 fail%UPARAM(1) = c2
252 fail%UPARAM(2) = x_1(1)
253 fail%UPARAM(3) = x_1(2)
254 fail%UPARAM(4) = x_2(1)
255 fail%UPARAM(5) = x_2(2)
256 fail%UPARAM(6) = x_2(3)
257 fail%UPARAM(7) = pthk
258 fail%UPARAM(8) = 0
259 fail%UPARAM(9) = c3
260 fail%UPARAM(10) = mflag
261 fail%UPARAM(11) = sflag
262 fail%UPARAM(12) = inst
263 fail%UPARAM(13) = ref_len
264 fail%UPARAM(14) = icoup
265 fail%UPARAM(15) = dcrit
266 fail%UPARAM(16) = exp
267 fail%UPARAM(17)= failip
268c---------------------------------------------------
269 IF (is_encrypted)THEN
270 WRITE(iout,'(5X,A,//)')'CONFIDENTIAL DATA'
271 ELSE
272 WRITE(iout,1000)
273 IF (mflag /= 0) WRITE(iout, 1100) mflag
274 WRITE(iout,1200) c1,c2,c3,c4,c5
275 WRITE(iout,1300) x_1(2),x_1(1),c2
276 WRITE(iout,1400) x_2(3),x_2(2),x_2(1)
277 WRITE(iout,1500) sflag
278 IF (sflag == 3) WRITE(iout,1600) inst
279 IF (reg_func > 0) WRITE(iout, 1700) reg_func,ref_len
280 IF (icoup > 0) THEN
281 WRITE(iout,1800) icoup,dcrit,exp
282 ENDIF
283 WRITE(iout, 1900) pthk,failip
284 WRITE(iout, 2000)
285 ENDIF
286c-----------------------------------------------------------------------
287 RETURN
288c-----------------------------------------------------------------------
289 1000 FORMAT(
290 & 5x,'-----------------------------------------------',/,
291 & 5x,' BIQUADRATIC FAILURE MODEL ',/,
292 & 5x,'-----------------------------------------------',/)
293 1100 FORMAT(
294 & 5x,'MATERIAL PARAMETER SELECTOR M-FLAG. . . . . . .=',i10/,
295 & 5x,' = 1 : MILD STEEL ',/,
296 & 5x,' = 2 : HSS STEEL ',/,
297 & 5x,' = 3 : UHSS STEEL ',/,
298 & 5x,' = 4 : ALUMINUM AA5182 ',/,
299 & 5x,' = 5 : ALUMINUM AA6082-T6 ',/,
300 & 5x,' = 6 : PLASTIC PA6GF30 ',/,
301 & 5x,' = 7 : PLASTIC PP T40 ',/,
302 & 5x,' = 99: USER DEFINED STRAIN RATIO ',/)
303 1200 FORMAT(
304 & 5x,'PLASTIC STRAINS AT FAILURE: ',/,
305 & 5x,'--------------------------- ',/,
306 & 5x,'C1 (SIMPLE COMPRESSION). . . . . . . . . . . .=',1pg20.13,/
307 & 5x,'C2 (SHEAR) . . . . . . . . . . . . . . . . . .=',1pg20.13,/
308 & 5x,'C3 (SIMPLE TENSION). . . . . . . . . . . . . .=',1pg20.13,/
309 & 5x,'C4 (PLANE STRAIN). . . . . . . . . . . . . . .=',1pg20.13,/
310 & 5x,'C5 (BIAXIAL TENSION) . . . . . . . . . . . . .=',1pg20.13,/)
311 1300 FORMAT(
312 & 5x,'COEFFICIENTS OF FIRST PARABOLA: ',/,
313 & 5x,'------------------------------ ',/,
314 & 5x,'A. . . . . . . . . . . . . . . . . . . . . . .=',1pg20.13,/
315 & 5x,'B. . . . . . . . . . . . . . . . . . . . . . .=',1pg20.13,/
316 & 5x,'C. . . . . . . . . . . . . . . . . . . . . . .=',1pg20.13,/)
317 1400 FORMAT(
318 & 5x,'COEFFICIENTS OF SECOND PARABOLA: ',/,
319 & 5x,'-------------------------------- ',/,
320 & 5x,'D. . . . . . . . . . . . . . . . . . . . . . .=',1pg20.13,/
321 & 5x,'E. . . . . . . . . . . . . . . . . . . . . . .=',1pg20.13,/
322 & 5x,'F. . . . . . . . . . . . . . . . . . . . . . .=',1pg20.13,/)
323 1500 FORMAT(
324 & 5x,'specific behavior flag s-flag. . . . . . . . .=',I10,/
325 & 5X,' = 1: two quadratic functions ',/,
326 & 5X,' = 2: plane strain VALUE is the global minimum',/,
327 & 5X,' = 3: plane strain is global minimum + instability necking(shells only)',/)
328 1600 FORMAT(
329 & 5X,'instability strain(shells only) . . . . . . .=',1PG20.13,/)
330 1700 FORMAT(
331 & 5X,'element length regularization: ',/,
332 & 5X,'------------------------------ ',/,
333 & 5X,'regularization FUNCTION id . . . . . . . . . .=',I10,/
334 & 5X,'referenze element length . . . . . . . . . . .=',1PG20.13,/)
335 1800 FORMAT(
336 & 5X,'stress softening: ',/,
337 & 5X,'----------------- ',/,
338 & 5X,'coupling method flag icoup . . . . . . . . . .=',I10,/,
339 & 5X,' = 1: classical coupling using critical damage',/,
340 & 5X,' = 2: necking instability coupling (SHELLS ONLY)',/,
341 & 5X,'damage critical value dcrit (IF ICOUP = 1) . .=',1PG20.13,/,
342 & 5X,'stress softening exponent exp. . . . . . . . .=',1PG20.13,/)
343 1900 FORMAT(
344 & 5X,'element deletion: ',/,
345 & 5X,'----------------- ',/,
346 & 5X,'shell element deletion parameter pthickfail. .=',1PG20.13,/
347 & 5X,' > 0.0: fraction of failed thickness ',/,
348 & 5X,' < 0.0: fraction of failed intg. points ',/,
349 & 5X,'number of failed intg. points prior to elem deletion .=',I10/)
350 2000 FORMAT(
351 & 5X,'-----------------------------------------------',/)
352c-----------------------------------------------------------------------
353 END
subroutine biquad_coefficients(c1, c2, c3, c4, c5, l, x_1, x_2, e1, e2, e3, e4)
end diagonal values have been computed in the(sparse) matrix id.SOL
subroutine hm_get_floatv(name, rval, is_available, lsubmodel, unitab)
subroutine hm_get_floatv_dim(name, dim_fac, is_available, lsubmodel, unitab)
subroutine hm_get_intv(name, ival, is_available, lsubmodel)
subroutine hm_option_is_encrypted(is_encrypted)
subroutine hm_read_fail_biquad(fail, mat_id, fail_id, irupt, titr, lsubmodel, unitab)
#define min(a, b)
Definition macros.h:20
#define max(a, b)
Definition macros.h:21
initmumps id
integer, parameter nchartitle
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:895