CalculiX  2.13
A Free Software Three-Dimensional Structural Finite Element Program
fluidsections.f File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine fluidsections (inpc, textpart, set, istartset, iendset, ialset, nset, ielmat, matname, nmat, irstrt, istep, istat, n, iline, ipol, inl, ipoinp, inp, lakon, ielprop, nprop, nprop_, prop, ipoinpc, mi)
 

Function/Subroutine Documentation

◆ fluidsections()

subroutine fluidsections ( character*1, dimension(*)  inpc,
character*132, dimension(16)  textpart,
character*81, dimension(*)  set,
integer, dimension(*)  istartset,
integer, dimension(*)  iendset,
integer, dimension(*)  ialset,
integer  nset,
integer, dimension(mi(3),*)  ielmat,
character*80, dimension(*)  matname,
integer  nmat,
integer  irstrt,
integer  istep,
integer  istat,
integer  n,
integer  iline,
integer  ipol,
integer  inl,
integer, dimension(2,*)  ipoinp,
integer, dimension(3,*)  inp,
character*8, dimension(*)  lakon,
integer, dimension(*)  ielprop,
integer  nprop,
integer  nprop_,
real*8, dimension(*)  prop,
integer, dimension(0:*)  ipoinpc,
integer, dimension(*)  mi 
)
25 !
26 ! reading the input deck: *FLUID SECTION
27 !
28  implicit none
29 !
30  logical liquid,manning
31 !
32  character*1 inpc(*)
33  character*7 elname
34  character*8 lakon(*)
35  character*80 matname(*),material,typename,typename_oil
36  character*81 set(*),elset
37  character*132 textpart(16)
38 !
39  integer mi(*),istartset(*),iendset(*),ialset(*),nconstants,
40  & ielmat(mi(3),*),irstrt,nset,nmat,ndprop,npropstart,ndpropread,
41  & istep,istat,n,key,i,j,k,imaterial,ipos,lprop,ipoinpc(0:*),
42  & iline,ipol,inl,ipoinp(2,*),inp(3,*),ielprop(*),nprop,nprop_,
43  & nodea,nodeb,noil _mat,npu,nfix,nstart,iset
44 !
45  real*8 prop(*)
46 
47  noil_mat=0
48  liquid=.false.
49  manning=.false.
50 !
51  if((istep.gt.0).and.(irstrt.ge.0)) then
52  write(*,*)
53  & '*ERROR reading *FLUID SECTION: *FLUID SECTION should'
54  write(*,*) ' be placed before all step definitions'
55  call exit(201)
56  endif
57 !
58  typename='
59  & '
60 !
61  do i=2,n
62  if(textpart(i)(1:9).eq.'MATERIAL=') then
63  material=textpart(i)(10:89)
64  elseif(textpart(i)(1:6).eq.'ELSET=') then
65  elset=textpart(i)(7:86)
66  elset(81:81)=' '
67  ipos=index(elset,' ')
68  elset(ipos:ipos)='E'
69 !
70  elseif(textpart(i)(1:5).eq.'TYPE=') then
71  read(textpart(i)(6:85),'(a80)',iostat=istat) typename
72  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
73  &"*FLUID SECTION%")
74 !
75  elseif(textpart(i)(1:4).eq.'OIL=') then
76  read(textpart(i)(5:85),'(a80)',iostat=istat) typename_oil
77  do j=1, nmat
78  if(matname(j)(1:80).eq. typename_oil(1:80)) then
79  noil_mat=j
80  exit
81  endif
82  enddo
83  if(j.gt.nmat) then
84  write(*,*)
85  & '*ERROR reading *FLUID SECTION: no oil with the'
86  write(*,*) ' name',typename_oil,'has been defined'
87  call exit(201)
88  endif
89 !
90  elseif(textpart(i)(1:6).eq.'LIQUID') then
91  liquid=.true.
92  elseif(textpart(i)(1:7).eq.'MANNING') then
93  manning=.true.
94  elseif(textpart(i)(1:10).eq.'CONSTANTS=') then
95  read(textpart(i)(11:20),'(i10)',iostat=istat) nconstants
96  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
97  &"*FLUID SECTION%")
98  else
99  write(*,*)
100  & '*WARNING reading *FLUID SECTION: parameter not recognized:'
101  write(*,*) ' ',
102  & textpart(i)(1:index(textpart(i),' ')-1)
103  call inputwarning(inpc,ipoinpc,iline,
104  &"*FLUID SECTION%")
105  endif
106  enddo
107 !
108 ! types of sections
109 !
110  if(typename(1:18).eq.'ABSOLUTETORELATIVE') then
111  elname='ATR '
112  ndprop=3
113 !
114 ! version of Yavor Dobrev
115 !
116  elseif(typename(1:10).eq.'ACCTUBEONE') then
117  elname ='ACCTUBO'
118  ndprop=19
119 !
120 ! version of Stefanie Asenbauer
121 !
122  elseif(typename(1:7).eq.'ACCTUBE') then
123  elname ='ACCTUBE'
124  ndprop=24
125 !
126  elseif(typename(1:12).eq.'BLEEDTAPPING') then
127  elname='ORBT '
128  ndprop=41
129 !
130 ! elements "BLINDLINK" with boundary massflows
131 ! (pressures and temperatures unknown)
132 !
133  elseif(typename(1:9).eq.'BLINDLINK') then
134  elname='BCMF '
135  ndprop=0
136 !
137  elseif(typename(1:6).eq.'BRANCH') then
138  if(typename(7:13).eq.'JOINTGE')then
139  elname='REBRJG '
140  ndprop=16
141  elseif(typename(7:20).eq.'JOINTIDELCHIK1')then
142  elname='REBRJI1'
143  ndprop=16
144  elseif(typename(7:20).eq.'JOINTIDELCHIK2')then
145  elname='REBRJI2'
146  ndprop=16
147  elseif(typename(7:13).eq.'SPLITGE')then
148  elname='REBRSG '
149  ndprop=16
150  elseif(typename(7:20).eq.'SPLITIDELCHIK1')then
151  elname='REBRSI1'
152  ndprop=16
153  elseif(typename(7:20).eq.'SPLITIDELCHIK2')then
154  elname='REBRSI2'
155  ndprop=16
156  endif
157 !
158  elseif(typename(1:10).eq.'CARBONSEAL') then
159  ndprop=4
160  if(typename(11:12).eq.'GE') then
161  elname='CARBSGE'
162  else
163  elname='CARBS '
164  endif
165 !
166  elseif(typename(1:12).eq.'CHANNELINOUT')then
167  elname='LICHIO '
168  ndprop=0
169 !
170  elseif(typename(1:7).eq.'CHANNEL') then
171  if(typename(8:17).eq.'SLUICEGATE') then
172  elname='LICHSG '
173  ndprop=7
174  elseif(typename(8:20).eq.'SLUICEOPENING') then
175  elname='LICHSO '
176  ndprop=6
177  elseif(typename(8:16).eq.'WEIRCREST') then
178  elname='LICHWE '
179  ndprop=7
180  elseif(typename(8:16).eq.'WEIRSLOPE') then
181  elname='LICHWO '
182  ndprop=6
183  elseif(typename(8:17).eq.'STRAIGHT') then
184  elname='LICH '
185  ndprop=7
186  elseif(typename(8:16).eq.'RESERVOIR') then
187  elname='LICHRE '
188  ndprop=6
189  elseif(typename(8:25).eq.'DISCONTINUOUSSLOPE') then
190  elname='LICHDS '
191  ndprop=9
192  elseif(typename(8:27).eq.'DISCONTINUOUSOPENING') then
193  elname='LICHDO '
194  ndprop=6
195  elseif(typename(8:18).eq.'CONTRACTION') then
196  elname='LICHCO '
197  ndprop=6
198  elseif(typename(8:18).eq.'ENLARGEMENT') then
199  elname='LICHEL '
200  ndprop=6
201  elseif(typename(8:11).eq.'STEP') then
202  elname='LICHST '
203  ndprop=6
204  elseif(typename(8:11).eq.'DROP') then
205  elname='LICHDR '
206  ndprop=6
207  else
208  write(*,*) '*ERROR reading *FLUID SECTION:'
209  write(*,*) ' unknown channel section'
210  call inputerror(inpc,ipoinpc,iline,
211  &"*FLUID SECTION%")
212  endif
213 !
214  elseif(typename(1:14).eq.'CHARACTERISTIC') then
215  ndprop=23
216  elname='CHAR '
217 !
218  elseif(typename(1:10).eq.'CROSSSPLIT')then
219  elname='CROSPL'
220  ndprop=12
221 !
222  elseif (typename(1:18).eq.'FREECONVECTIONFLOW') then
223  elname='FCVF '
224  ndprop=6
225 !
226  elseif (typename(1:16).eq.'FREEDISCPUMPFLOW') then
227  elname='FDPF '
228  ndprop=7
229 !
230  elseif(typename(1:12).eq.'GASPIPEFANNO') then
231 !
232 ! gaspipe version Fanno(friction and oil)
233 !
234  if(typename(13:21).eq.'ADIABATIC') then
235  if(typename(22:27).eq.'ALBERS') then
236  elname='GAPFAA '
237  ndprop=9
238  elseif(typename(22:28).eq.'FRIEDEL') then
239  elname='GAPFAF '
240  ndprop=9
241  elseif(typename(22:35).eq.'FLEXIBLERADIUS') then
242  elname='GAPFAFR'
243  ndprop=9
244  elseif(typename(22:44).eq.'FLEXIBLERADIUSANDLENGTH') then
245  elname='GAPFARL'
246  ndprop=9
247  else
248  elname='GAPFA '
249  ndprop=9
250  endif
251 
252  elseif(typename(13:22).eq.'ISOTHERMAL') then
253  if(typename(23:28).eq.'ALBERS') then
254  elname='GAPFIA '
255  ndprop=9
256  elseif(typename(23:29).eq.'FRIEDEL') then
257  elname='GAPFIF '
258  ndprop=9
259  elseif(typename(22:35).eq.'FLEXIBLERADIUS') then
260  elname='GAPFIFR'
261  ndprop=9
262  elseif(typename(22:44).eq.'FLEXIBLERADIUSANDLENGTH') then
263  elname='GAPFIRL'
264  ndprop=9
265  else
266  elname='GAPFI '
267  ndprop=9
268  endif
269  endif
270 !
271 ! inlet for the general case (flows, pressures and
272 ! temperatures unknown)
273 !
274  elseif(typename(1:5).eq.'INLET') then
275  elname='INLT '
276  ndprop=0
277 !
278  elseif(typename(1:5).eq.'INOUT')then
279  elname='IO '
280  ndprop=0
281 !
282  elseif(typename(1:9).eq.'LABYRINTH') then
283  if(typename(10:17).eq.'FLEXIBLE') then
284  ndprop=23
285  if(typename(18:23).eq.'SINGLE') then
286  elname='LABFSN '
287  elseif(typename(18:25).eq.'STRAIGHT') then
288  elname='LABFSR '
289  elseif(typename(18:24).eq.'STEPPED') then
290  elname='LABFSP '
291  endif
292  elseif(typename(10:14).eq.'DUMMY') then
293  ndprop=1
294  elname='LABD '
295  else
296  ndprop=15
297  if(typename(10:15).eq.'SINGLE') then
298  elname='LABSN '
299  elseif(typename(10:17).eq.'STRAIGHT') then
300  elname='LABSR '
301  elseif(typename(10:16).eq.'STEPPED') then
302  elname='LABSP '
303  endif
304  endif
305 !
306  elseif(typename(1:10).eq.'LIQUIDPUMP') then
307 C ndprop=20
308  ndprop=19
309  elname='LIPU '
310 !
311  elseif (typename(1:15).eq.'MASSFLOWPERCENT') then
312  elname='MFPC '
313  ndprop=11
314 !
315  elseif(typename(1:8).eq.'MOEHRING') then
316  if(typename(9:19).eq.'CENTRIFUGAL') then
317  elname='MRGF '
318  elseif(typename(9:19).eq.'CENTRIPETAL') then
319  elname='MRGP '
320  endif
321  ndprop=11
322 !
323  else if (typename(1:6).eq.'O-PUMP') then
324  elname='LDOP '
325  ndprop=3
326 !
327  elseif(typename(1:7).eq.'ORIFICE') then
328  ndprop=9
329  if(typename(8:13).eq.'OWN_AL') then
330  elname='ORMA '
331  elseif(typename(8:12).eq.'PK_AL') then
332  elname='ORPA '
333  elseif(typename(8:12).eq.'MS_MS') then
334  elname='ORMM '
335  elseif(typename(8:12).eq.'PK_MS') then
336  elname='ORPM '
337  elseif(typename(8:12).eq.'BRAGG') then
338  elname='ORBG '
339  elseif(typename(8:14).eq.'RINGGAP') then
340  elname='ORRG '
341  elseif(typename(8:17).eq.'SEGMENTGAP') then
342  elname='ORSG '
343  elseif(typename(8:11).eq.'BORE') then
344  elname='ORBO '
345  elseif(typename(8:22).eq.'GEOMETRICALAREA') then
346  elname='ORGA '
347  elseif(typename(8:20).eq.'BRUSHFLEXIBLE') then
348  elname='ORBF '
349  else
350  elname='ORC1 '
351  endif
352 !
353 ! inlet for the general case (flows, pressures and
354 ! temperatures unknown)
355 !
356  elseif(typename(1:6).eq.'OUTLET') then
357  elname='OUTLT '
358  ndprop=0
359 !
360  elseif(typename(1:9).eq.'PIPEINOUT')then
361  elname='LIPIO '
362  ndprop=0
363 !
364 ! liquid pipe
365 !
366  elseif(typename(1:4).eq.'PIPE') then
367  if(typename(5:8).eq.'BEND') then
368  elname='LIPIBE '
369  ndprop=4
370  elseif(typename(5:10).eq.'BRANCH') then
371  elname='LIPIBR '
372  ndprop=6
373  elseif(typename(5:15).eq.'CONTRACTION') then
374  elname='LIPICO '
375  ndprop=2
376  elseif(typename(5:13).eq.'DIAPHRAGM') then
377  elname='LIPIDI '
378  ndprop=2
379  elseif(typename(5:15).eq.'ENLARGEMENT') then
380  elname='LIPIEL '
381  ndprop=2
382  elseif(typename(5:12).eq.'ENTRANCE') then
383  elname='LIPIEN '
384  ndprop=2
385  elseif(typename(5:13).eq.'GATEVALVE') then
386  elname='LIPIGV '
387  ndprop=2
388  elseif(typename(5:11).eq.'MANNING') then
389  if(typename(12:19).eq.'FLEXIBLE') then
390  elname='LIPIMAF'
391  ndprop=3
392  else
393  elname='LIPIMA '
394  ndprop=3
395  endif
396  elseif(typename(5:19).eq.'WHITE-COLEBROOK') then
397  if(typename(20:27).eq.'FLEXIBLE') then
398  elname='LIPIWCF'
399  ndprop=5
400  else
401  elname='LIPIWC '
402  ndprop=5
403  endif
404  endif
405 !
406  elseif(typename(1:14).eq.'PRESWIRLNOZZLE') then
407  elname='ORPN '
408  ndprop=41
409 !
410  elseif(typename(1:13).eq.'RADIALOUTFLOW') then
411  if(typename(14:24).eq.'RADIALINLET')then
412  elname='RORI '
413 
414  else if(typename(14:23).eq.'AXIALINLET')then
415  elname='ROAI '
416  endif
417  ndprop=7
418 !
419  elseif(typename(1:18).eq.'RELATIVETOABSOLUTE') then
420  elname='RTA '
421  ndprop=3
422 !
423  elseif(typename(1:10).eq.'RESTRICTOR') then
424  if(typename(11:15).eq.'USER') then
425  elname='REUS '
426  ndprop=10
427 !
428  elseif(typename(11:15).eq.'BRUSH') then
429  elname='RBSF'
430  ndprop=10
431  elseif(typename(11:29).eq.'LONGORIFICEIDELCHIK') then
432  elname='RELOID '
433  ndprop=10
434 !
435  elseif(typename(11:21).eq.'WALLORIFICE') then
436  elname='REWAOR '
437  ndprop=10
438  elseif(typename(11:18).eq.'ENTRANCE') then
439  elname='REEN '
440  ndprop=10
441  elseif(typename(11:21).eq.'ENLARGEMENT') then
442  elname='REEL '
443  ndprop=10
444  elseif(typename(11:21).eq.'CONTRACTION') then
445  elname='RECO '
446  ndprop=10
447  elseif(typename(11:23).eq.'BENDIDELCIRC') then
448  elname='REBEIDC'
449  ndprop=10
450  elseif(typename(11:23).eq.'BENDIDELRECT') then
451  elname='REBEIDR'
452  ndprop=10
453  elseif(typename(11:20).eq.'BENDMILLER') then
454  elname='REBEMI '
455  ndprop=10
456  elseif(typename(11:17).eq.'BENDOWN') then
457  elname='REBEMA '
458  ndprop=10
459  elseif(typename(11:14).eq.'EXIT') then
460  elname='REEX '
461  ndprop=10
462  elseif(typename(11:33).eq.'LONGORIFICELICHTAROWICZ') then
463  elname='RELOLI '
464  ndprop=10
465  endif
466 !
467  elseif((typename(1:7).eq.'RIMSEAL').and.
468  & (typename(1:15).ne.'RIMSEALFLEXIBLE')) then
469  elname='RIMS '
470  ndprop=6
471 !
472  elseif(typename(1:15).eq.'RIMSEALFLEXIBLE') then
473  elname='RIMFLEX'
474  ndprop=14
475 !
476  elseif(typename(1:14).eq.'ROTATINGCAVITY') then
477  if(typename(15:20).eq.'LINEAR')then
478  elname='RCVL '
479  else if(typename(15:23).eq.'NONLINEAR')then
480  elname='RCVN '
481  endif
482  ndprop=6
483 !
484  else if (typename(1:6).eq.'S-PUMP') then
485  elname='SPUMP '
486  ndprop=6
487 !
488  elseif(typename(1:6).eq.'VORTEX') then
489  if(typename(7:10).eq.'FREE') then
490  elname='VOFR '
491  ndprop=10
492  elseif(typename(7:12).eq.'FORCED') then
493  elname='VOFO '
494  ndprop=10
495  endif
496 !
497  elseif(typename(1:1).eq.' ') then
498  elname=' '
499  ndprop=0
500  elseif(typename(1:1).eq.'U') then
501  elname(1:7)=typename(1:7)
502  ndprop=nconstants
503  else
504  write(*,*) '*ERROR reading *FLUID SECTION: ',typename
505  write(*,*) ' is an unknown fluid section type'
506  call exit(201)
507  endif
508 !
509 ! check for the existence of the set and the material
510 !
511  do i=1,nmat
512  if(matname(i).eq.material) exit
513  enddo
514  if(i.gt.nmat) then
515  write(*,*)
516  & '*ERROR reading *FLUID SECTION: nonexistent material'
517  write(*,*) ' '
518  call inputerror(inpc,ipoinpc,iline,
519  &"*FLUID SECTION%")
520  call exit(201)
521  endif
522  imaterial=i
523 !
524  do i=1,nset
525  if(set(i).eq.elset) exit
526  enddo
527  if(i.gt.nset) then
528  elset(ipos:ipos)=' '
529  write(*,*) '*ERROR reading *FLUID SECTION: element set ',elset
530  write(*,*) ' has not yet been defined. '
531  call inputerror(inpc,ipoinpc,iline,
532  &"*FLUID SECTION%")
533  call exit(201)
534  endif
535  iset=i
536 !
537  npropstart=nprop
538 !
539 ! reading the element properties
540 !
541 ! liquid pump / gas characteristic /
542 ! preswirl nozzle / bleed tapping
543 ! these network elements contain curves with a variable
544 ! amount of data points
545 !
546  if((elname(1:4).eq.'LIPU').or.(elname(1:4).eq.'CHAR').or.
547  & (elname(1:4).eq.'ORPN').or.(elname(1:4).eq.'ORBT')) then
548 !
549 ! determine the number of fixed inputs (before entering
550 ! a curve with arbitrary many data points)
551 !
552  if(elname(1:4).eq.'LIPU') then
553  nfix=1
554  elseif(elname(1:4).eq.'CHAR') then
555  nfix=4
556  elseif(elname(1:4).eq.'ORPN') then
557  nfix=6
558  elseif(elname(1:4).eq.'ORBT') then
559  nfix=4
560  endif
561  ndprop=0
562  npu=0
563 c do
564 c call getnewline(inpc,textpart,istat,n,key,iline,ipol,inl,
565 c & ipoinp,inp,ipoinpc)
566 c if((istat.lt.0).or.(key.eq.1)) exit
567 c if(ndprop.eq.0) then
568 c do j=1,nfix
569 c ndprop=ndprop+1
570 c read(textpart(j),'(f40.0)',iostat=istat)
571 c & prop(nprop+ndprop)
572 c if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
573 c &"*FLUID SECTION%")
574 c enddo
575 c nstart=nfix+1
576 c!
577 c! adding the point (0.,0.) for characteristics
578 c!
579 c if(elname(1:4).eq.'CHAR') then
580 c prop(nprop+3)=0.d0
581 c prop(nprop+4)=0.d0
582 c npu=2
583 c nfix=2
584 c endif
585 c else
586 c nstart=1
587 c endif
588 c do j=nstart,n
589 c npu=npu+1
590 c ndprop=ndprop+1
591 !
592 ! npu is the number of data points, i.e. the number of
593 ! abscissa-values + the number of ordinate values
594 !
595  do
596  call getnewline(inpc,textpart,istat,n,key,iline,ipol,inl,
597  & ipoinp,inp,ipoinpc)
598  if((istat.lt.0).or.(key.eq.1)) exit
599  do j=1,n
600  ndprop=ndprop+1
601  if(ndprop.eq.nfix) then
602  if(elname(1:4).eq.'CHAR') then
603  prop(nprop+3)=0.d0
604  prop(nprop+4)=0.d0
605  npu=2
606  nfix=2
607  cycle
608  endif
609  elseif(ndprop.gt.nfix) then
610  npu=npu+1
611  endif
612 !
613 ! check whether number of data points is smaller than allowed
614 !
615  if(elname(1:4).eq.'LIPU') then
616  if(ndprop.gt.19) then
617  write(*,*) 'ERROR reading *FLUID SECTION: more'
618  write(*,*) ' than 9 pairs were defined for'
619  write(*,*) ' fluid section type LIQUID PUMP'
620  call exit(201)
621  endif
622  elseif(elname(1:4).eq.'CHAR') then
623  if(ndprop.gt.23) then
624  write(*,*) 'ERROR reading *FLUID SECTION: more'
625  write(*,*) ' than 9 pairs were defined for'
626  write(*,*)
627  & ' fluid section type CHARACTERISTIC'
628  call exit(201)
629  endif
630  elseif(elname(1:4).eq.'ORPN') then
631  if(ndprop.gt.41) then
632  write(*,*) 'ERROR reading *FLUID SECTION: more'
633  write(*,*) ' than 17 pairs were defined for'
634  write(*,*)
635  & ' fluid section type PRESWIRL NOZZLE'
636  call exit(201)
637  endif
638  elseif(elname(1:4).eq.'ORBT') then
639  if(ndprop.gt.41) then
640  write(*,*) 'ERROR reading *FLUID SECTION: more'
641  write(*,*) ' than 18 pairs were defined for'
642  write(*,*) ' fluid section type BLEED TAPPING'
643  call exit(201)
644  endif
645  endif
646 !
647  read(textpart(j),'(f40.0)',iostat=istat)
648  & prop(nprop+ndprop)
649  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
650  &"*FLUID SECTION%")
651  enddo
652  enddo
653 !
654 ! filling up to nfix if less values were specified
655 !
656  if(ndprop.lt.nfix) then
657  do j=ndprop+1,nfix
658  prop(nprop+j)=0.d0
659  enddo
660  ndprop=nfix
661  endif
662 !
663 ! check whether data points are paired
664 !
665  if(2*(npu/2).ne.npu) then
666  write(*,*) '*WARNING reading *FLUID SECTION: less y data'
667  write(*,*) ' points x data points for fluid'
668  write(*,*) ' section type',elname(1:4)
669  npu=npu-1
670 c call exit(201)
671  endif
672 !
673  prop(nprop+nfix)=npu/2
674 !
675 ! check whether data points are in appropriate order
676 ! (only for liquid pump so far)
677 !
678  if(elname(1:4).eq.'LIPU') then
679  do j=2,npu/2,2
680  if(prop(nprop+nfix+2*j-1).le.prop(nprop+nfix+2*j-3)) then
681  write(*,*) '*ERROR reading *FLUID SECTION: volumetric'
682  write(*,*) ' flow data points must be in'
683  write(*,*) ' strict increasing order'
684  call exit(201)
685  elseif(prop(nprop+nfix+2*j).gt.prop(nprop+nfix+2*j-2))
686  & then
687  write(*,*) '*ERROR reading *FLUID SECTION: total'
688  write(*,*) ' head data points must be '
689  write(*,*) ' decreasing'
690  call exit(201)
691  endif
692  enddo
693  endif
694 !
695  nprop=nprop+ndprop+1
696 !
697  elseif(elname.eq.'ACCTUBE') then
698 !
699 ! reading the element properties
700 !
701 ! acc-tube (proprietary)
702 !
703 ! Read the first 20 elements
704 !
705  lprop=0
706  ndpropread=ndprop
707 !
708  do j=1,(ndpropread-1)/8+1
709  call getnewline(inpc,textpart,istat,n,key,iline,ipol,inl,
710  & ipoinp,inp,ipoinpc)
711  if((istat.lt.0).or.(key.eq.1)) exit
712  do k=1,8
713  lprop=lprop+1
714  if(lprop.gt.ndpropread) exit
715  read(textpart(k),'(f40.0)',iostat=istat)
716  & prop(nprop+lprop)
717 
718  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
719  &"*FLUID SECTION%")
720 
721 ! If 20 elements have been read, check how many more
722 ! are to be read
723  if (lprop.eq.20) then
724  ndpropread = lprop +
725  & (prop(nprop+19)+prop(nprop+20))*2
726  ndprop = ndpropread
727  endif
728  enddo
729  enddo
730 !
731 ! Until now 3 lines have been read
732 ! If necessary read the rest
733 !
734  if(ndpropread.gt.lprop) then
735  do j=3,(ndpropread-1)/8+1
736  call getnewline(inpc,textpart,istat,n,key,iline,ipol,inl,
737  & ipoinp,inp,ipoinpc)
738  if((istat.lt.0).or.(key.eq.1)) exit
739  do k=1,8
740  lprop=lprop+1
741  if(lprop.gt.ndpropread) exit
742  read(textpart(k),'(f40.0)',iostat=istat)
743  & prop(nprop+lprop)
744 
745  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
746  &"*FLUID SECTION%")
747  enddo
748  enddo
749  endif
750 !
751  nprop=nprop+ndprop
752 !
753  call getnewline(inpc,textpart,istat,n,key,iline,ipol,inl,
754  & ipoinp,inp,ipoinpc)
755 !
756  elseif(ndprop.gt.0) then
757 !
758 ! reading the element properties
759 !
760 ! general case
761 !
762  lprop=0
763 !
764  do
765  call getnewline(inpc,textpart,istat,n,key,iline,ipol,inl,
766  & ipoinp,inp,ipoinpc)
767  if((istat.lt.0).or.(key.eq.1)) exit
768  do k=1,n
769  lprop=lprop+1
770  if(lprop.gt.ndprop) exit
771  read(textpart(k),'(f20.0)',iostat=istat)
772  & prop(nprop+lprop)
773  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
774  &"*FLUID SECTION%")
775  enddo
776  enddo
777  nprop=nprop+ndprop
778  else
779  call getnewline(inpc,textpart,istat,n,key,iline,ipol,inl,
780  & ipoinp,inp,ipoinpc)
781  endif
782 !
783  if(nprop.gt.nprop_) then
784  write(*,*) '*ERROR reading *FLUID SECTION: increase nprop_'
785  call exit(201)
786  endif
787 !
788 ! complementing the properties of restrictor elements
789 !
790  if((elname(1:6).eq.'REBEMI').or.
791  & (elname(1:6).eq.'REBEMA')) then
792  prop(npropstart+7)=noil_mat
793 !
794  elseif(elname(1:7).eq.'REBEIDC') then
795  prop(npropstart+7)=noil_mat
796 !
797  elseif(elname(1:7).eq.'REBEIDR') then
798  prop(npropstart+9)=noil_mat
799 !
800  elseif(elname(1:4).eq.'REEX') then
801 !
802  prop(npropstart+6)=noil_mat
803 !
804  elseif(elname(1:6).eq.'REWAOR') then
805 !
806  prop(npropstart+6)=noil_mat
807 !
808  elseif(elname(1:4).eq.'REEN') then
809 !
810  prop(npropstart+6)=noil_mat
811 !
812 ! zeta (loss coefficient for an entry)
813 !
814  elseif(elname(1:7).eq.'REBRJI1') then
815  prop(npropstart+11)=noil_mat
816 !
817  elseif(elname(1:7).eq.'REBRJI2') then
818  prop(npropstart+11)=noil_mat
819  elseif(elname(1:6).eq.'REBRJG') then
820  prop(npropstart+11)=noil_mat
821  elseif(elname(1:6).eq.'REBRSG') then
822  prop(npropstart+11)=noil_mat
823  elseif(elname(1:7).eq.'REBRSI1') then
824  prop(npropstart+15)=noil_mat
825  elseif(elname(1:7).eq.'REBRSI2') then
826  prop(npropstart+13)=noil_mat
827  endif
828 !
829 ! k_oil for restrictors type USER, ENTRENCE, LONG ORIFICE IDELCHIK,
830 ! EXIT, ORIFICE IN A WALL, LONG ORIFICE LICHTAROWICZ
831 !
832  if((elname(1:4).eq.'REUS').or.
833  & (elname(1:6).eq.'RELOID').or.
834  & (elname(1:6).eq.'RELOLI')) then
835 
836  prop(npropstart+6)= noil_mat
837 !
838 ! k_oil for restrictors type ENLARGEMENT
839 !
840  elseif (elname(1:4).eq.'REEL') then
841  prop(npropstart+5)= noil_mat
842 !
843 ! k_oil for restrictors type CONTRACTION, BEN MILLER, BEND MA
844 ! (BEND IDELCHICK has already been treated)
845 !
846  elseif(elname(1:4).eq.'RECO') then
847  prop(npropstart+7)= noil_mat
848 !
849 ! k_oil for pipe elements
850 !
851  elseif(elname(1:3).eq.'GAP') then
852  prop(npropstart+7)=noil_mat
853  elseif(elname(1:4).eq.'LICH') then
854 !
855 ! if Manning: change sign of 5th entry as marker
856 !
857  if(manning) then
858  if((elname(5:6).eq.'SO').or.
859  & (elname(5:6).eq.'WO').or.
860  & (elname(5:6).eq.' ').or.
861  & (elname(5:6).eq.'RE').or.
862  & (elname(5:6).eq.'DS').or.
863  & (elname(5:6).eq.'DO')) then
864  prop(npropstart+5)=-prop(npropstart+5)
865  endif
866  endif
867  endif
868 !
869 ! assigning the elements of the set the appropriate material
870 ! and property pointer
871 !
872  do j=istartset(iset),iendset(iset)
873  if(ialset(j).gt.0) then
874  if(lakon(ialset(j))(1:1).ne.'D') then
875  write(*,*) '*ERROR reading *FLUID SECTION: element ',
876  & ialset(j),' is no fluid element.'
877  call exit(201)
878  endif
879  lakon(ialset(j))(2:8)=elname(1:7)
880 !
881 ! gas type elements used for liquids are labeled with LP
882 !
883  if((liquid.and.(lakon(ialset(j))(2:3).eq.'RE')).or.
884  & (liquid.and.(lakon(ialset(j))(2:3).eq.'OR')))
885  & lakon(ialset(j))(2:3)='LP'
886 !
887  if(liquid.and.(lakon(ialset(j))(2:3).eq.'VO')) then
888  lakon(ialset(j))(2:3)='LP'
889  if(lakon(ialset(j))(4:5).eq.'FR') then
890  lakon(ialset(j))(4:5)='VF'
891  else if(lakon(ialset(j))(4:5).eq.'FO') then
892  lakon(ialset(j))(4:5)='VS'
893  endif
894  endif
895 !
896  if(liquid.and.((lakon(ialset(j))(4:5).eq.'BG').or.
897  & (lakon(ialset(j))(4:5).eq.'BT').or.
898  & (lakon(ialset(j))(4:5).eq.'MA').or.
899  & (lakon(ialset(j))(4:5).eq.'MM').or.
900  & (lakon(ialset(j))(4:5).eq.'PA').or.
901  & (lakon(ialset(j))(4:5).eq.'PM').or.
902  & (lakon(ialset(j))(4:5).eq.'PN'))) then
903  write(*,*) ''
904  write(*,*) '*ERROR reading *FLUID SECTION: element ',k,
905  &' is no valid incompressible orifice element.'
906  write(*,*) ' Please change element type',
907  &' and/or definition.'
908  write(*,*) ' Calculation stopped.'
909  call exit(201)
910  endif
911 !
912  ielmat(1,ialset(j))=imaterial
913  if(ndprop.gt.0) ielprop(ialset(j))=npropstart
914  else
915  k=ialset(j-2)
916  do
917  k=k-ialset(j)
918  if(k.ge.ialset(j-1)) exit
919  if(lakon(k)(1:1).ne.'D') then
920  write(*,*) '*ERROR reading *FLUID SECTION: element ',
921  & k,' is no fluid element.'
922  call exit(201)
923  endif
924  lakon(k)(2:8)=elname(1:7)
925 !
926 ! gas type elements used for liquids are labeled with LP
927 !
928  if((liquid.and.(lakon(ialset(j))(2:3).eq.'RE')).or.
929  & (liquid.and.(lakon(ialset(j))(2:3).eq.'OR')))
930  & lakon(ialset(j))(2:3)='LP'
931 !
932  if(liquid.and.(lakon(ialset(j))(2:3).eq.'VO')) then
933  lakon(ialset(j))(2:3)='LP'
934  if(lakon(ialset(j))(4:5).eq.'FR') then
935  lakon(ialset(j))(4:5)='VF'
936  else if(lakon(ialset(j))(4:5).eq.'FO') then
937  lakon(ialset(j))(4:5)='VS'
938  endif
939  endif
940 !
941  if(liquid.and.((lakon(ialset(j))(4:5).eq.'BG').or.
942  & (lakon(ialset(j))(4:5).eq.'BT').or.
943  & (lakon(ialset(j))(4:5).eq.'MA').or.
944  & (lakon(ialset(j))(4:5).eq.'MM').or.
945  & (lakon(ialset(j))(4:5).eq.'PA').or.
946  & (lakon(ialset(j))(4:5).eq.'PM').or.
947  & (lakon(ialset(j))(4:5).eq.'PN'))) then
948  write(*,*) ''
949  write(*,*)'*ERROR reading *FLUID SECTION: element ',k,
950  & ' is no valid incompressible orifice element.'
951  write(*,*) ' Please change element type ',
952  &'and/or definition.'
953  write(*,*) ' Calculation stopped.'
954  call exit(201)
955  endif
956  ielmat(1,k)=imaterial
957  if(ndprop.gt.0) ielprop(k)=npropstart
958  enddo
959  endif
960  enddo
961 !
962  return
subroutine inputwarning(inpc, ipoinpc, iline, text)
Definition: inputwarning.f:20
subroutine getnewline(inpc, textpart, istat, n, key, iline, ipol, inl, ipoinp, inp, ipoinpc)
Definition: getnewline.f:21
subroutine inputerror(inpc, ipoinpc, iline, text)
Definition: inputerror.f:20
Hosted by OpenAircraft.com, (Michigan UAV, LLC)