OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
write_matparam.F File Reference
#include "implicit_f.inc"

Go to the source code of this file.

Functions/Subroutines

subroutine write_matparam (mat_elem, len)

Function/Subroutine Documentation

◆ write_matparam()

subroutine write_matparam ( type(mat_elem_), intent(in) mat_elem,
integer, intent(inout) len )

Definition at line 38 of file write_matparam.F.

39C-----------------------------------------------
40C M o d u l e s
41C-----------------------------------------------
42 USE mat_elem_mod
44 use write_therpmaram_mod
45 use write_ale_rezoning_param_mod , only : write_ale_rezoning_param
46! use write_eosparam_mod , only : write_eosparam
47C-----------------------------------------------
48C I m p l i c i t T y p e s
49C-----------------------------------------------
50#include "implicit_f.inc"
51C-----------------------------------------------
52C D u m m y A r g u m e n t s
53C-----------------------------------------------
54 INTEGER ,INTENT(INOUT) :: LEN
55 TYPE(MAT_ELEM_) ,INTENT(IN) :: MAT_ELEM
56C-----------------------------------------------
57C L o c a l V a r i a b l e s
58C-----------------------------------------------
59 INTEGER :: I,NUMMAT,IMAT,NUPARAM,NIPARAM,NUMTABL,NFAIL
60 INTEGER :: IAD,NFIX,NFIXR,LENI,LENR,NMOD,MOD
61 INTEGER ,DIMENSION(NCHARTITLE) :: NAME
62 INTEGER ,DIMENSION(:) ,ALLOCATABLE :: IBUF
63 INTEGER :: NBSUBMAT, OLDFORMAT
64 my_real ,DIMENSION(:) ,ALLOCATABLE :: rbuf
65C=======================================================================
66 nfix = 24
67 nummat = mat_elem%NUMMAT
68 len = 0
69 leni = nfix*nummat
70 ALLOCATE (ibuf(leni) )
71
72 ! write integer parameters and flags
73 iad = 0
74 DO imat = 1,nummat
75 ibuf(iad+ 1) = mat_elem%MAT_PARAM(imat)%ILAW
76 ibuf(iad+ 2) = mat_elem%MAT_PARAM(imat)%MAT_ID
77 ibuf(iad+ 3) = mat_elem%MAT_PARAM(imat)%NUPARAM
78 ibuf(iad+ 4) = mat_elem%MAT_PARAM(imat)%NIPARAM
79 ibuf(iad+ 5) = mat_elem%MAT_PARAM(imat)%NFUNC
80 ibuf(iad+ 6) = mat_elem%MAT_PARAM(imat)%NTABLE
81 ibuf(iad+ 7) = mat_elem%MAT_PARAM(imat)%NSUBMAT
82 ibuf(iad+ 8) = mat_elem%MAT_PARAM(imat)%NFAIL
83 ibuf(iad+ 9) = mat_elem%MAT_PARAM(imat)%IVISC
84 ibuf(iad+10) = mat_elem%MAT_PARAM(imat)%IEOS
85 ibuf(iad+11) = mat_elem%MAT_PARAM(imat)%ITHERM
86 ibuf(iad+12) = mat_elem%MAT_PARAM(imat)%IEXPAN
87 ibuf(iad+13) = mat_elem%MAT_PARAM(imat)%IALE
88 ibuf(iad+14) = mat_elem%MAT_PARAM(imat)%ITURB
89 ibuf(iad+15) = mat_elem%MAT_PARAM(imat)%HEAT_FLAG
90 ibuf(iad+16) = mat_elem%MAT_PARAM(imat)%COMPRESSIBILITY
91 ibuf(iad+17) = mat_elem%MAT_PARAM(imat)%SMSTR
92 ibuf(iad+18) = mat_elem%MAT_PARAM(imat)%STRAIN_FORMULATION
93 ibuf(iad+19) = mat_elem%MAT_PARAM(imat)%IPRES
94 ibuf(iad+20) = mat_elem%MAT_PARAM(imat)%ORTHOTROPY
95 ibuf(iad+21) = mat_elem%MAT_PARAM(imat)%NLOC
96 ibuf(iad+22) = mat_elem%MAT_PARAM(imat)%IFAILWAVE
97 ibuf(iad+23) = mat_elem%MAT_PARAM(imat)%IXFEM
98 ibuf(iad+24) = mat_elem%MAT_PARAM(imat)%NMOD
99
100 iad = iad + nfix
101 END DO
102
103 CALL write_i_c(ibuf,leni)
104 len = len + leni
105 DEALLOCATE(ibuf)
106
107 ! write real parameters
108 nfixr = 9
109 lenr = nfixr*nummat
110 iad = 0
111 ALLOCATE (rbuf(lenr) )
112 DO imat = 1,nummat
113 rbuf(iad + 1) = mat_elem%MAT_PARAM(imat)%RHO
114 rbuf(iad + 2) = mat_elem%MAT_PARAM(imat)%RHO0
115 rbuf(iad + 3) = mat_elem%MAT_PARAM(imat)%YOUNG
116 rbuf(iad + 4) = mat_elem%MAT_PARAM(imat)%BULK
117 rbuf(iad + 5) = mat_elem%MAT_PARAM(imat)%SHEAR
118 rbuf(iad + 6) = mat_elem%MAT_PARAM(imat)%NU
119 rbuf(iad + 7) = mat_elem%MAT_PARAM(imat)%STIFF_CONTACT
120 rbuf(iad + 8) = mat_elem%MAT_PARAM(imat)%STIFF_HGLASS
121 rbuf(iad + 9) = mat_elem%MAT_PARAM(imat)%STIFF_TSTEP
122 iad = iad + nfixr
123 END DO
124 CALL write_db(rbuf ,lenr)
125 len = len + lenr
126 DEALLOCATE(rbuf)
127
128 ! write material title
129 DO imat = 1,nummat
130 IF(imat < nummat) THEN
131 DO i=1,nchartitle
132 name(i) = ichar(mat_elem%MAT_PARAM(imat)%TITLE(i:i))
133 END DO
134 ELSE
135 name = 0
136 ENDIF
137 CALL write_c_c(name,nchartitle)
138 END DO
139
140 ! write material parameter array
141 DO imat = 1,nummat
142 nuparam = mat_elem%MAT_PARAM(imat)%NUPARAM
143 niparam = mat_elem%MAT_PARAM(imat)%NIPARAM
144 IF (nuparam > 0) THEN
145 CALL write_db(mat_elem%MAT_PARAM(imat)%UPARAM ,nuparam)
146 END IF
147 IF (niparam > 0) THEN
148 CALL write_i_c(mat_elem%MAT_PARAM(imat)%IPARAM ,niparam)
149 END IF
150 len = len + nuparam + niparam
151 END DO
152
153 ! write function tables 4D
154 DO imat = 1,nummat
155 numtabl = mat_elem%MAT_PARAM(imat)%NTABLE
156 IF (numtabl > 0) THEN
157 CALL write_mat_table(mat_elem%MAT_PARAM(imat)%TABLE, numtabl)
158 len = len + leni + lenr
159 END IF
160 END DO
161
162 ! write viscosity model parameters
163 DO imat = 1,nummat
164 IF (mat_elem%MAT_PARAM(imat)%IVISC > 0) THEN
165 CALL write_viscparam(mat_elem%MAT_PARAM(imat)%VISC,len)
166 END IF
167 END DO
168c
169c write thermal parameters
170c
171 DO imat = 1,nummat
172 CALL write_thermparam(mat_elem%MAT_PARAM(imat)%THERM,len)
173 END DO
174c
175c write parameters of failure models per material
176c
177 DO imat = 1,nummat
178 nfail = mat_elem%MAT_PARAM(imat)%NFAIL
179 IF (nfail > 0) THEN
180 DO i = 1,nfail
181 CALL write_failparam(mat_elem%MAT_PARAM(imat)%FAIL(i),len)
182 END DO
183 END IF
184 END DO
185
186 ! write damage modes
187 DO imat = 1,nummat
188 nmod = mat_elem%MAT_PARAM(imat)%NMOD
189 IF (nmod > 0) THEN
190 DO mod = 1,nmod
191 DO i=1,nchartitle
192 name(i) = ichar(mat_elem%MAT_PARAM(imat)%MODE(mod)(i:i))
193 END DO
194 CALL write_c_c(name,nchartitle)
195 ENDDO
196 ENDIF
197 ENDDO
198
199 ! write eos model parameters
200 DO imat = 1,nummat
201 IF (mat_elem%MAT_PARAM(imat)%IEOS > 0) THEN
202 CALL write_eosparam(mat_elem%MAT_PARAM(imat)%EOS)
203 END IF
204 END DO
205
206
207 ! Write multimaterial buffer
208 ALLOCATE (ibuf(21) ) !< mex is 21 submaterials
209 ALLOCATE (rbuf(21) )
210 DO imat = 1,nummat
211 nbsubmat = mat_elem%MAT_PARAM(imat)%MULTIMAT%NB
212 oldformat = mat_elem%MAT_PARAM(imat)%MULTIMAT%old_data_format
213 ibuf(1) = nbsubmat
214 ibuf(2) = oldformat
215 CALL write_i_c(ibuf, 2)
216 IF(nbsubmat > 0)THEN
217 ibuf(1:nbsubmat) = mat_elem%MAT_PARAM(imat)%MULTIMAT%mid
218 rbuf(1:nbsubmat) = mat_elem%MAT_PARAM(imat)%MULTIMAT%vfrac
219 CALL write_i_c(ibuf, nbsubmat)
220 CALL write_db(rbuf, nbsubmat)
221 ! write eos model parameters
222 IF(oldformat == 1)THEN
223 !embedded EOS parameters
224 DO i = 1,nbsubmat
225 CALL write_eosparam(mat_elem%MAT_PARAM(imat)%MULTIMAT%EOS(i))
226 END DO
227 ENDIF
228 ENDIF
229 ENDDO
230 DEALLOCATE(ibuf)
231 DEALLOCATE(rbuf)
232
233
234
235
236 ! write ALE rezoning parameters
237 DO imat = 1,nummat
238 CALL write_ale_rezoning_param(mat_elem%MAT_PARAM(imat)%REZON)
239 END DO
240
241!-----------
242 RETURN
#define my_real
Definition cppsort.cpp:32
integer, parameter nchartitle
subroutine write_failparam(fail, len)
subroutine write_viscparam(visc, len)
subroutine write_mat_table(table, numtabl)
subroutine write_db(a, n)
Definition write_db.F:142
void write_i_c(int *w, int *len)
void write_c_c(int *w, int *len)