SPEED
READ_FILEMATE.f90
Go to the documentation of this file.
1! Copyright (C) 2012 The SPEED FOUNDATION
2! Author: Ilario Mazzieri
3!
4! This file is part of SPEED.
5!
6! SPEED is free software; you can redistribute it and/or modify it
7! under the terms of the GNU Affero General Public License as
8! published by the Free Software Foundation, either version 3 of the
9! License, or (at your option) any later version.
10!
11! SPEED is distributed in the hope that it will be useful, but
12! WITHOUT ANY WARRANTY; without even the implied warranty of
13! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
14! Affero General Public License for more details.
15!
16! You should have received a copy of the GNU Affero General Public License
17! along with SPEED. If not, see <http://www.gnu.org/licenses/>.
18
135
136
137 subroutine read_filemate(filemate,nb_mate,char_mat,type_mat,trefm,lab_mat,Qua_S,Qua_P, &
138 nb_mate_nle,char_mat_nle,val_mat_nle,type_mat_nle,lab_mat_nle, &
139 nb_mate_rnd, lab_mat_rnd, &
140 nb_diriX,val_diriX,fnc_diriX,lab_diriX, &
141 nb_diriY,val_diriY,fnc_diriY,lab_diriY, &
142 nb_diriZ,val_diriZ,fnc_diriZ,lab_diriZ, &
143 nb_neuX,val_neuX,fnc_neuX,lab_neuX, &
144 nb_neuY,val_neuY,fnc_neuY,lab_neuY, &
145 nb_neuZ,val_neuZ,fnc_neuZ,lab_neuZ, &
146 nb_neuN,val_neuN,fnc_neuN,lab_neuN, &
147 nb_poiX,val_poiX,fnc_poiX, &
148 nb_poiY,val_poiY,fnc_poiY, &
149 nb_poiZ,val_poiZ,fnc_poiZ, &
150 ntX,valtX,ftX, &
151 ntY,valtY,ftY, &
152 ntZ,valtZ,ftZ, &
153 nb_plaX,val_plaX,fnc_plaX,lab_plaX, &
154 nb_plaY,val_plaY,fnc_plaY,lab_plaY, &
155 nb_plaZ,val_plaZ,fnc_plaZ,lab_plaZ, &
156 nb_forX,val_forX,fnc_forX, &
157 nb_forY,val_forY,fnc_forY, &
158 nb_forZ,val_forZ,fnc_forZ, &
159 nb_forF,val_forF,fnc_forF, &
160 nb_pre,val_pre,fnc_pre, &
161 nb_she,val_she,fnc_she, &
162 ntest,ftest, & !valftest,&
163 nb_abc,lab_abc, &
164 nb_dg,lab_dg,lab_dg_yn,lab_dg_link,&
165 lab_dg_frc, val_dg_frc, nb_frac,&
166 srcmodflag, szsism, &
167 nb_sism,val_sism,fnc_sism,lab_sism, &
168 nb_expl,val_expl,fnc_expl,lab_expl, &
169 nb_case,val_case,lab_case,tol_case, &
170 nb_nhee,val_nhe,tol_nhe, &
171 nb_fnc,type_fnc,ind_fnc,dat_fnc,lab_fnc,nb_fnc_data, &
172 fmax,fpeak,damping_val)
173
174
176! use speed_par, only: slip_type
177
178 implicit none
179
180 character*70 :: filemate,fileinput
181 character*100000 :: inline
182 character*4 :: keyword
183
184 integer*4 :: nb_mate,nb_fnc, nb_fnc_data, nb_frac, damping_val
185 integer*4 :: nb_mate_nle, nb_mate_rnd, nb_nhee
186 integer*4 :: nb_diriX,nb_diriY,nb_diriZ,nb_neuX,nb_neuY,nb_neuZ
187 integer*4 :: nb_neuN
188 integer*4 :: nb_poiX,nb_poiY,nb_poiZ,nb_forX,nb_forY,nb_forZ,nb_forF
189 integer*4 :: nb_plaX,nb_plaY,nb_plaZ,ntX,ntY,ntZ
190 integer*4 :: nb_abc, nb_dg
191 integer*4 :: nb_pre,nb_she
192 integer*4 :: isism,nb_sism
193 integer*4 :: iexpl,nb_expl
194 integer*4 :: icase,nb_case
195
196 integer*4 :: im,ifunc,idf,ndat_fnc,file_nd, inhee
197 integer*4 :: im_nle, im_rnd
198 integer*4 :: idX,idY,idZ,inX,inY,inZ, itX, itY, itZ
199 integer*4 :: inN, itest
200 integer*4 :: ipX,ipY,ipZ,ifX,ifY,ifZ,iff
201 integer*4 :: iplX,iplY,iplZ
202 integer*4 :: ipr,ish,iabc,idg
203 integer*4 :: ileft,iright
204 integer*4 :: i,j,dummy,status
205 integer*4 :: ntest, srcmodflag, szsism
206
207 integer*4, dimension(nb_diriX) :: fnc_diriX,lab_diriX
208 integer*4, dimension(nb_diriY) :: fnc_diriY,lab_diriY
209 integer*4, dimension(nb_diriZ) :: fnc_diriZ,lab_diriZ
210 integer*4, dimension(nb_neuX) :: fnc_neuX,lab_neuX
211 integer*4, dimension(nb_neuY) :: fnc_neuY,lab_neuY
212 integer*4, dimension(nb_neuZ) :: fnc_neuZ,lab_neuZ
213 integer*4, dimension(nb_neuN) :: fnc_neuN,lab_neuN
214 integer*4, dimension(ntest) :: ftest
215 integer*4, dimension(nb_forX) :: fnc_forX
216 integer*4, dimension(nb_forY) :: fnc_forY
217 integer*4, dimension(nb_forZ) :: fnc_forZ
218 integer*4, dimension(nb_forF) :: fnc_forF
219 integer*4, dimension(nb_pre) :: fnc_pre
220 integer*4, dimension(nb_she) :: fnc_she
221 integer*4, dimension(nb_abc) :: lab_abc
222 integer*4, dimension(nb_dg) :: lab_dg, lab_dg_yn, lab_dg_frc, lab_dg_link
223 integer*4, dimension(nb_sism) :: fnc_sism,lab_sism
224 integer*4, dimension(nb_expl) :: fnc_expl,lab_expl
225 integer*4, dimension(nb_case) :: lab_case
226 integer*4, dimension(nb_poiX) :: fnc_poiX !building_id_x
227 integer*4, dimension(nb_poiY) :: fnc_poiY !building_id_y
228 integer*4, dimension(nb_poiZ) :: fnc_poiZ !building_id_z
229 integer*4, dimension(ntX) :: ftX
230 integer*4, dimension(ntY) :: ftY
231 integer*4, dimension(ntZ) :: ftZ
232 integer*4, dimension(nb_plaX) :: fnc_plaX,lab_plaX
233 integer*4, dimension(nb_plaY) :: fnc_plaY,lab_plaY
234 integer*4, dimension(nb_plaZ) :: fnc_plaZ,lab_plaZ
235 integer*4, dimension(nb_mate) :: type_mat
236 integer*4, dimension(nb_mate) :: lab_mat
237 integer*4, dimension(nb_fnc) :: type_fnc
238 integer*4, dimension(nb_fnc) :: lab_fnc
239 integer*4, dimension(nb_fnc +1) :: ind_fnc
240 integer*4, dimension(nb_mate_nle) :: lab_mat_nle
241 integer*4, dimension(nb_mate_rnd) :: lab_mat_rnd
242
243 integer*4, dimension(nb_mate_nle) :: type_mat_nle
244 integer*4, dimension(nb_mate_nle,1) :: char_mat_nle
245 integer*4, dimension(nb_case) :: val_case
246 integer*4, dimension(nb_nhee) :: val_nhe
247
248 real*8 :: fmax
249 real*8 :: fpeak
250 real*8 :: rho, vs, vp
251
252 real*8, dimension(nb_case) :: tol_case
253 real*8, dimension(nb_nhee) :: tol_nhe
254 real*8, dimension(nb_fnc_data) :: dat_fnc
255 real*8, dimension(nb_mate) :: trefm, qua_s, qua_p
256
257 real*8, dimension(nb_diriX,4) :: val_dirix
258 real*8, dimension(nb_diriY,4) :: val_diriy
259 real*8, dimension(nb_diriZ,4) :: val_diriz
260 real*8, dimension(nb_neuX,4) :: val_neux
261 real*8, dimension(nb_neuY,4) :: val_neuy
262 real*8, dimension(nb_neuZ,4) :: val_neuz
263 real*8, dimension(nb_neuN,4) :: val_neun
264 real*8, dimension(nb_poiX,4) :: val_poix
265 real*8, dimension(nb_poiY,4) :: val_poiy
266 real*8, dimension(nb_poiZ,4) :: val_poiz
267
268 real*8, dimension(ntX,4) :: valtx
269 real*8, dimension(ntY,4) :: valty
270 real*8, dimension(ntZ,4) :: valtz
271 real*8, dimension(nb_plaX,1) :: val_plax
272 real*8, dimension(nb_plaY,1) :: val_play
273 real*8, dimension(nb_plaZ,1) :: val_plaz
274 real*8, dimension(nb_forX,4) :: val_forx
275 real*8, dimension(nb_forY,4) :: val_fory
276 real*8, dimension(nb_forZ,4) :: val_forz
277 real*8, dimension(nb_forF,10) :: val_forf
278 real*8, dimension(nb_pre,10) :: val_pre
279 real*8, dimension(nb_she,10) :: val_she
280 real*8, dimension(nb_sism,szsism) :: val_sism
281 real*8, dimension(nb_expl,20) :: val_expl
282 real*8, dimension(nb_mate_nle,1) :: val_mat_nle
283 real*8, dimension(nb_mate,4) :: char_mat
284 real*8, dimension(nb_dg,2) :: val_dg_frc
285
286 open(40,file=filemate)
287
288
289
290 im = 0; im_nle = 0; icase = 0 ; im_rnd = 0;
291 idx = 0; idy = 0; idz = 0
292 inx = 0; iny = 0; inz = 0; inn = 0;
293 ipx = 0; ipy = 0; ipz = 0
294 iplx = 0; iply = 0; iplz = 0
295 ifx = 0; ify = 0; ifz = 0
296 iff = 0; ipr = 0; ish = 0
297 iabc = 0; idg = 0
298 isism = 0; iexpl = 0; itz = 0;
299 itest = 0
300
301 ifunc = 0
302 fmax = 0.d0
303
304 inhee = 0;
305
306 !building_id_x = -1;
307 !building_id_y = -1;
308 !building_id_z = -1;
309
310 if (nb_fnc.gt.0) ind_fnc(1) = 1
311
312
313 do
314 read(40,'(A)',iostat = status) inline
315
316 if (status.ne.0) exit
317
318 keyword = inline(1:4)
319
320 ileft = 0
321 iright = len(inline)
322 do i = 1,iright
323 if (inline(i:i).eq.' ') exit
324 enddo
325 ileft = i
326
327 select case (keyword)
328
329 case('MATE')
330 im = im + 1
331! read(inline(ileft:iright),*) lab_mat(im),type_mat(im),&
332! char_mat(im,1),char_mat(im,2),char_mat(im,3), & !char_mat(im,4)
333! Qua_S(im), Qua_P(im)
334 read(inline(ileft:iright),*) lab_mat(im),type_mat(im),&
335 rho, vs, vp, & !char_mat(im,4)
336 qua_s(im), qua_p(im)
337
338 ! write(*,*) lab_mat(im),type_mat(im),&
339 ! rho, VS, VP, & !char_mat(im,4)
340 ! Qua_S(im), Qua_P(im)
341 ! read(*,*)
342
343 if(damping_val .eq. 2) then
344 qua_s(im) = 0.5d0*qua_s(im)
345 qua_p(im) = 0.5d0*qua_p(im)
346 endif
347
348
349 char_mat(im,1) = rho
350 !lambda
351 char_mat(im,2) = rho * (vp**2 - 2*vs**2)
352 !mu
353 char_mat(im,3) = rho * vs**2
354
355 case('MATN')
356 im_nle = im_nle + 1
357 read(inline(ileft:iright),*) lab_mat_nle(im_nle),type_mat_nle(im_nle),&
358 char_mat_nle(im_nle,1),val_mat_nle(im_nle,1)
359
360 case('MATR')
361 im_rnd = im_rnd + 1
362 read(inline(ileft:iright),*) lab_mat_rnd(im_rnd)
363
364 case('DIRX')
365 idx = idx + 1
366 read(inline(ileft:iright),*) lab_dirix(idx),fnc_dirix(idx),&
367 val_dirix(idx,1),val_dirix(idx,2),val_dirix(idx,3),val_dirix(idx,4)
368
369 case('DIRY')
370 idy = idy + 1
371 read(inline(ileft:iright),*) lab_diriy(idy),fnc_diriy(idy),&
372 val_diriy(idy,1),val_diriy(idy,2),val_diriy(idy,3),val_diriy(idy,4)
373
374 case('DIRZ')
375 idz = idz + 1
376 read(inline(ileft:iright),*) lab_diriz(idz),fnc_diriz(idz),&
377 val_diriz(idz,1),val_diriz(idz,2),val_diriz(idz,3),val_diriz(idz,4)
378
379 case('NEUX')
380 inx = inx + 1
381 read(inline(ileft:iright),*) lab_neux(inx),fnc_neux(inx),&
382 val_neux(inx,1),val_neux(inx,2),val_neux(inx,3),val_neux(inx,4)
383
384 case('NEUY')
385 iny = iny + 1
386 read(inline(ileft:iright),*) lab_neuy(iny),fnc_neuy(iny),&
387 val_neuy(iny,1),val_neuy(iny,2),val_neuy(iny,3),val_neuy(iny,4)
388
389 case('NEUZ')
390 inz = inz + 1
391 read(inline(ileft:iright),*) lab_neuz(inz),fnc_neuz(inz),&
392 val_neuz(inz,1),val_neuz(inz,2),val_neuz(inz,3),val_neuz(inz,4)
393
394 case('NEUN')
395 inn = inn + 1
396 read(inline(ileft:iright),*) lab_neun(inn),fnc_neun(inn),&
397 val_neun(inn,1),val_neun(inn,2),val_neun(inn,3),val_neun(inn,4)
398
399 case('PLOX')
400 ipx = ipx + 1
401 read(inline(ileft:iright),*) fnc_poix(ipx),&
402 val_poix(ipx,1),val_poix(ipx,2),val_poix(ipx,3),val_poix(ipx,4)
403
404 case('PLOY')
405 ipy = ipy + 1
406 read(inline(ileft:iright),*) fnc_poiy(ipy),&
407 val_poiy(ipy,1),val_poiy(ipy,2),val_poiy(ipy,3),val_poiy(ipy,4)
408
409 case('PLOZ')
410 ipz = ipz + 1
411 read(inline(ileft:iright),*) fnc_poiz(ipz),&
412 val_poiz(ipz,1),val_poiz(ipz,2),val_poiz(ipz,3),val_poiz(ipz,4)
413
414 !! For SSI (AH, SS)
415 !case('PLOD')
416 !ipX = ipX + 1; ipY = ipY + 1; ipZ = ipZ + 1
417 !read(inline(ileft:iright),*) fnc_poiX(ipX),val_poiX(ipX,4),building_id_x(ipX) !!! function id, value of the applied load, building id
418 !val_poiX(ipX,6)=1
419
420 !fnc_poiY(ipY)=fnc_poiX(ipX)
421 !val_poiY(ipY,4)=val_poiX(ipX,4); building_id_y(ipY)=building_id_x(ipX)
422 !val_poiY(ipY,6)=2
423
424 !fnc_poiZ(ipZ)=fnc_poiX(ipX)
425 !val_poiZ(ipZ,4)=val_poiX(ipX,4); building_id_z(ipZ)=building_id_x(ipX)
426 !val_poiZ(ipZ,6)=3
427 !!
428
429! case('TLOX')
430! itX = itX + 1
431! read(inline(ileft:iright),*) ftX(itX),valtX(itX,1),valtX(itX,4)
432
433! case('TLOY')
434! itY = itY + 1
435! read(inline(ileft:iright),*) ftY(itY),valtY(itY,1),valtY(itY,4)
436
437 case('TLOZ')
438 itz = itz + 1
439 ! function id, velocity of travelling load, amplitude wave
440 read(inline(ileft:iright),*) ftz(itz),valtz(itz,1),valtz(itz,2)
441
442 case('PLAX')
443 iplx = iplx + 1
444 read(inline(ileft:iright),*)fnc_plax(iplx),&
445 lab_plax(iplx),val_plax(iplx,1)
446
447 case('PLAY')
448 iply = iply + 1
449 read(inline(ileft:iright),*)fnc_play(iply),&
450 lab_play(iply),val_play(iply,1)
451
452 case('PLAZ')
453 iplz = iplz + 1
454 read(inline(ileft:iright),*)fnc_plaz(iplz),&
455 lab_plaz(iplz),val_plaz(iplz,1)
456
457 case('FORX')
458 ifx = ifx + 1
459 read(inline(ileft:iright),*) fnc_forx(ifx),&
460 val_forx(ifx,1),val_forx(ifx,2),val_forx(ifx,3),val_forx(ifx,4)
461
462 case('FORY')
463 ify = ify + 1
464 read(inline(ileft:iright),*) fnc_fory(ify),&
465 val_fory(ify,1),val_fory(ify,2),val_fory(ify,3),val_fory(ify,4)
466
467 case('FORZ')
468 ifz = ifz + 1
469 read(inline(ileft:iright),*) fnc_forz(ifz),&
470 val_forz(ifz,1),val_forz(ifz,2),val_forz(ifz,3),val_forz(ifz,4)
471
472 case('TEST')
473 itest = itest + 1
474 read(inline(ileft:iright),*) ftest(itest)
475
476 case('FORC')
477 iff = iff + 1
478 read(inline(ileft:iright),*) fnc_forf(iff),&
479 val_forf(iff,1),val_forf(iff,2),val_forf(iff,3),val_forf(iff,4),&
480 val_forf(iff,5),val_forf(iff,6),val_forf(iff,7),&
481 val_forf(iff,8),val_forf(iff,9),val_forf(iff,10)
482
483 case('PRES')
484 ipr = ipr + 1
485 read(inline(ileft:iright),*) fnc_pre(ipr),&
486 val_pre(ipr,1),val_pre(ipr,2),val_pre(ipr,3),val_pre(ipr,4),&
487 val_pre(ipr,5),val_pre(ipr,6),val_pre(ipr,7),&
488 val_pre(ipr,8),val_pre(ipr,9),val_pre(ipr,10)
489
490 case('SHEA')
491 ish = ish + 1
492 read(inline(ileft:iright),*) fnc_she(ish),&
493 val_she(ish,1),val_she(ish,2),val_she(ish,3),val_she(ish,4),&
494 val_she(ish,5),val_she(ish,6),val_she(ish,7),&
495 val_she(ish,8),val_she(ish,9),val_she(ish,10)
496
497 case('ABSO')
498 iabc = iabc + 1
499 read(inline(ileft:iright),*) lab_abc(iabc)
500
501 case('DGIC')
502 idg = idg + 1
503 if(nb_frac .eq. 0) then
504 read(inline(ileft:iright),*) lab_dg(idg), lab_dg_yn(idg)
505 else
506 read(inline(ileft:iright),*) lab_dg(idg), lab_dg_yn(idg), lab_dg_frc(idg), val_dg_frc(idg,1), val_dg_frc(idg,2)
507 endif
508
509 case('DGIX')
510 idg = idg+1
511 read(inline(ileft:iright),*) lab_dg(idg), lab_dg_yn(idg), lab_dg_link(idg)
512
513
514 case('SISM')
515 isism = isism + 1
516 if (srcmodflag.eq.0) then
517 read(inline(ileft:iright),*) fnc_sism(isism),&
518 lab_sism(isism),val_sism(isism,1),val_sism(isism,2),&
519 val_sism(isism,3),val_sism(isism,4),val_sism(isism,5),&
520 val_sism(isism,6),val_sism(isism,7),val_sism(isism,8),&
521 val_sism(isism,9),val_sism(isism,10),val_sism(isism,11),&
522 val_sism(isism,12),val_sism(isism,13),val_sism(isism,14),&
523 val_sism(isism,15),val_sism(isism,16),val_sism(isism,17),&
524 val_sism(isism,18),val_sism(isism,19),val_sism(isism,20),&
525 val_sism(isism,21)
526 elseif (srcmodflag.eq.1) then
527 read(inline(ileft:iright),*) fnc_sism(isism),&
528 lab_sism(isism),val_sism(isism,1),val_sism(isism,2),&
529 val_sism(isism,3),val_sism(isism,4),val_sism(isism,5),&
530 val_sism(isism,6),val_sism(isism,7),val_sism(isism,8),&
531 val_sism(isism,9),val_sism(isism,10),val_sism(isism,11),&
532 val_sism(isism,12),val_sism(isism,13),val_sism(isism,14),&
533 val_sism(isism,15)
534 endif
535
536 case('EXPL')
537 iexpl = iexpl + 1
538 read(inline(ileft:iright),*) fnc_expl(iexpl),&
539 lab_expl(iexpl),val_expl(iexpl,1),val_expl(iexpl,2),&
540 val_expl(iexpl,3),val_expl(iexpl,4),val_expl(iexpl,5),&
541 val_expl(iexpl,6),val_expl(iexpl,7),val_expl(iexpl,8),&
542 val_expl(iexpl,9),val_expl(iexpl,10),val_expl(iexpl,11),&
543 val_expl(iexpl,12),val_expl(iexpl,13),val_expl(iexpl,14),&
544 val_expl(iexpl,15),val_expl(iexpl,16),val_expl(iexpl,17),&
545 val_expl(iexpl,18),val_expl(iexpl,19),val_expl(iexpl,20)
546
547 case('CASE')
548 icase = icase + 1
549 read(inline(ileft:iright),*) &
550 !case tag, block id, tolerance
551 lab_case(icase),val_case(icase),tol_case(icase)
552 !lab_case,val_case,tol_case
553
554 case('NHEE')
555 inhee = inhee + 1
556 tol_nhe(inhee) = 0.0
557 read(inline(ileft:iright),*) val_nhe(inhee) !, tol_case(inhee)
558
559 case('FMAX') !!! frequency at which Qs and Qp are computed (for plane wave load)
560 read(inline(ileft:iright),*) fmax
561
562 case('FPEK') !!! peak frequency for non linear damping
563 read(inline(ileft:iright),*) fpeak
564
565! case('SLIP')
566! read(inline(ileft:iright),*) slip_type
567
568 case('FUNC')
569 ifunc = ifunc + 1
570 read(inline(ileft:iright),*) lab_fnc(ifunc),type_fnc(ifunc)
571
572 select case (type_fnc(ifunc))
573
574 case(0,32)
575 ind_fnc(ifunc +1) = ind_fnc(ifunc) + 0
576
577 case(1,2)
578 ind_fnc(ifunc +1) = ind_fnc(ifunc) + 2
579 read(inline(ileft:iright),*)dummy,dummy,&
580 (dat_fnc(j), j = ind_fnc(ifunc),ind_fnc(ifunc +1) -1)
581
582
583 case(3,30,31,33)
584 read(inline(ileft:iright),*)dummy,dummy,ndat_fnc,fileinput
585 ind_fnc(ifunc +1) = ind_fnc(ifunc) + 2*ndat_fnc
586
587 open(24,file=fileinput)
588 read(24,*) file_nd
589 if (ndat_fnc .ne. file_nd) then
590 write(*,*) .ne.'Error reading function ! ndat_fnc file_nd !'
591 write(*,*) 'Error reading function from file ', trim(fileinput), '!'
592 write(*,*) 'Line numbers not consistent with material file.'
593 call exit(exit_function_error)
594 endif
595 do j = 1, file_nd
596 i = ind_fnc(ifunc) + 2*(j -1)
597 read(24,*) dat_fnc(i), dat_fnc(i +1)
598 enddo
599 close(24)
600
601 case(4,5,6,7)
602 ind_fnc(ifunc +1) = ind_fnc(ifunc) + 2
603 read(inline(ileft:iright),*)dummy,dummy,&
604 (dat_fnc(j), j = ind_fnc(ifunc),ind_fnc(ifunc +1) -1)
605
606
607 case(8,9)
608 ind_fnc(ifunc +1) = ind_fnc(ifunc) + 1
609 read(inline(ileft:iright),*)dummy,dummy,&
610 (dat_fnc(j), j = ind_fnc(ifunc),ind_fnc(ifunc +1) -1)
611
612
613 !SIGMOIDAL FUNCTION
614 case(12)
615 ind_fnc(ifunc +1) = ind_fnc(ifunc) + 3
616 read(inline(ileft:iright),*)dummy,dummy,&
617 (dat_fnc(j), j = ind_fnc(ifunc),ind_fnc(ifunc +1) -1)
618
619 !GRENOBLE BENCHMARK
620 case(13)
621 ind_fnc(ifunc +1) = ind_fnc(ifunc) + 2
622 read(inline(ileft:iright),*)dummy,dummy,&
623 (dat_fnc(j), j = ind_fnc(ifunc),ind_fnc(ifunc +1) -1)
624
625 !SCEC BENCHMARK
626 case(14)
627 ind_fnc(ifunc +1) = ind_fnc(ifunc) + 2
628 read(inline(ileft:iright),*)dummy,dummy,&
629 (dat_fnc(j), j = ind_fnc(ifunc),ind_fnc(ifunc +1) -1)
630
631 !ERF FUNCTION
632 case(15)
633 ind_fnc(ifunc +1) = ind_fnc(ifunc) + 4
634 read(inline(ileft:iright),*)dummy,dummy,&
635 (dat_fnc(j), j = ind_fnc(ifunc),ind_fnc(ifunc +1) -1)
636
637 !ADD ARCHULETA FUNC
638 case(21)
639 ind_fnc(ifunc +1) = ind_fnc(ifunc) + 0
640
641
642 case(50,55)
643 ind_fnc(ifunc +1) = ind_fnc(ifunc) + 2
644 read(inline(ileft:iright),*)dummy,dummy,&
645 (dat_fnc(j), j = ind_fnc(ifunc),ind_fnc(ifunc +1) -1)
646
647 case(60,62)
648 read(inline(ileft:iright),*)dummy,dummy,ndat_fnc
649 ind_fnc(ifunc +1) = ind_fnc(ifunc) + 2*ndat_fnc
650 read(inline(ileft:iright),*)dummy,dummy,dummy,&
651 (dat_fnc(j), j = ind_fnc(ifunc),ind_fnc(ifunc +1) -1)
652
653 case(61,63)
654 read(inline(ileft:iright),*)dummy,dummy,ndat_fnc
655 ind_fnc(ifunc +1) = ind_fnc(ifunc) + 2*ndat_fnc
656 read(inline(ileft:iright),*)dummy,dummy,dummy,&
657 (dat_fnc(j), j = ind_fnc(ifunc),ind_fnc(ifunc +1) -1)
658
659 case(99)
660 ind_fnc(ifunc +1) = ind_fnc(ifunc) + 2
661 read(inline(ileft:iright),*)dummy,dummy,&
662 (dat_fnc(j), j = ind_fnc(ifunc),ind_fnc(ifunc +1) -1)
663
664 case(100)
665 ind_fnc(ifunc +1) = ind_fnc(ifunc) + 1
666 read(inline(ileft:iright),*) dummy,dummy,&
667 (dat_fnc(j), j = ind_fnc(ifunc),ind_fnc(ifunc +1) -1)
668
669 case(773)
670 read(inline(ileft:iright),*)dummy,dummy,ndat_fnc,fileinput
671 ind_fnc(ifunc +1) = ind_fnc(ifunc) + ndat_fnc
672
673 open(24,file=fileinput)
674 read(24,*) file_nd
675 if (ndat_fnc .ne. file_nd) then
676 write(*,*) .ne.'Error reading function ! ndat_fnc file_nd !'
677 write(*,*) 'Error reading function from file ', trim(fileinput), '!'
678 write(*,*) 'Line numbers not consistent with material file.'
679 call exit(exit_function_error)
680 endif
681 do j = 1, file_nd
682 i = ind_fnc(ifunc) + (j -1)
683 read(24,*) dat_fnc(i)
684 enddo
685 close(24)
686
687 !! AH
688 case(777) !!! reaction force from structure is used for the cases where this function is called AH
689 read(inline(ileft:iright),*) dummy, dummy, ndat_fnc, dat_fnc(ind_fnc(ifunc)) !, dat_fnc(ind_fnc(ifunc)+1)
690 ind_fnc(ifunc +1) = ind_fnc(ifunc) + 1 !mod SS
691
692 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!
693
694 end select
695
696 end select
697
698 enddo
699
700 close(40)
701 return
702
703 end subroutine read_filemate
704
subroutine read_filemate(filemate, nb_mate, char_mat, type_mat, trefm,
Reads and stores data from filemate (*.mate)
SPEED exit codes.
Definition MODULES.f90:25
integer, parameter exit_function_error
Definition MODULES.f90:45