1!
2!     CalculiX - A 3-dimensional finite element program
3!     Copyright (C) 1998-2021 Guido Dhondt
4!
5!     This program is free software; you can redistribute it and/or
6!     modify it under the terms of the GNU General Public License as
7!     published by the Free Software Foundation(version 2);
8!
9!
10!     This program is distributed in the hope that it will be useful,
11!     but WITHOUT ANY WARRANTY; without even the implied warranty of
12!     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13!     GNU General Public License for more details.
14!
15!     You should have received a copy of the GNU General Public License
16!     along with this program; if not, write to the Free Software
17!     Foundation, Inc., 675 Mass Ave, Cambridge, MA 02139, USA.
18!
19!     check wether the property array contains values or IDs
20!     if first gas iteration (iin =0) then
21!     for each property of each fluid element
22!     if IDs: the IDs contained in prop array are stored in prop_store array;
23!     the IDs are interpreted and values are set/stored in prop array
24!     for each property of each fluid element
25!     else (after convergence)
26!     if IDs: the IDs stored in prop_store array are copied back in
27!     prop array
28!
29!     2016.05.19 Y. Muller
30!
31      subroutine checkinputvaluesnet(ieg,nflow,prop,ielprop,lakon)
32!
33      implicit none
34!
35      character*8 lakon(*)
36!
37      integer ieg(*),nflow,i,ielprop(*),index,nelem
38!
39      real*8 prop(*)
40!
41      do i=1,nflow
42        nelem=ieg(i)
43        index=ielprop(nelem)
44        if(index.lt.0) cycle
45!
46!     modifying the prop array (formerly done in fluidsections.f)
47!
48        if(lakon(nelem)(2:8).eq.'REBRJI2') then
49          if(1.d0-(prop(index+5)+prop(index+6))/
50     &         prop(index+4).gt.0.01d0)then
51            write(*,*) '*ERROR in checkinputvaluesnet:'
52            write(*,*) '       in element type RESTRICTOR
53     &BRANCH JOINT IDELCHIK2'
54            write(*,*) '       A0 ist not equal to A1+A2'
55            write(*,*) '       element number: ',nelem
56            call exit(201)
57          endif
58!
59        elseif((lakon(nelem)(2:3).eq.'OR').and.
60     &         (lakon(nelem)(2:5).ne.'ORC1').and.
61     &         (lakon(nelem)(2:5).ne.'ORB1').and.
62     &         (lakon(nelem)(2:5).ne.'ORB2').and.
63     &         (lakon(nelem)(2:5).ne.'ORBT').and.
64     &         (lakon(nelem)(2:5).ne.'ORPN').and.
65     &         (lakon(nelem)(2:5).ne.'ORFL'))then
66          if(prop(index+2).lt.0.d0) then
67            write(*,*) '*ERROR in checkinputvaluesnet: diameter '
68            write(*,*) '       of the orifice is not positive'
69            write(*,*) '       element number: ',nelem
70            call exit(201)
71          endif
72          if(prop(index+3).lt.0.d0) then
73            write(*,*) '*ERROR in checkinputvaluesnet:'
74            write(*,*) '        length of the orifice is ',
75     &           'not positive'
76            write(*,*) '       element number: ',nelem
77            call exit(201)
78          endif
79          if((lakon(nelem)(2:5).ne.'ORC1').and.
80     &         (lakon(nelem)(2:5).ne.'ORBF').and.
81     &         (lakon(nelem)(2:5).ne.'ORRG').and.
82     &         (lakon(nelem)(2:5).ne.'ORSG').and.
83     &         (lakon(nelem)(2:5).ne.'ORGA').and.
84     &         (lakon(nelem)(2:5).ne.'ORBO'))then
85            if((prop(index+4).gt.1.d-20).and.
86     &           (prop(index+5).gt.1.d-20)) then
87              write(*,*)
88     &             '*ERROR in checkinputvaluesnet:'
89              write(*,*) 'either the radius of ',
90     &             'the orifice must be zero or the chamfer angle'
91              write(*,*) '       element number: ',nelem
92              call exit(201)
93            endif
94          endif
95          if(prop(index+4)/prop(index+2).lt.0.d0) then
96            write(*,*) '*ERROR in checkinputvaluesnet: '
97            write(*,*)
98     &           '        r/d of an orifice must not be negative'
99            write(*,*) '       element number: ',nelem
100            call exit(201)
101          endif
102          if(prop(index+4)/prop(index+2).gt.0.82d0) then
103            write(*,*)
104     &           '*WARNING in checkinputvaluesnet: '
105            write(*,*) '        r/d of an orifice ',
106     &           'must not exceed 0.82'
107            write(*,*) '         element number: ',nelem
108          endif
109          if(prop(index+5).lt.0.d0) then
110            write(*,*) '*ERROR in checkinputvaluesnet:'
111            write(*,*) '       the chamfer angle of an ',
112     &           'orifice must not be negative'
113            write(*,*) '       element number: ',nelem
114            call exit(201)
115          endif
116          if(prop(index+5).gt.90.d0) then
117            write(*,*)
118     &           '*ERROR in checkinputvaluesnet: '
119            write(*,*) '       the chamfer angle of an ',
120     &           'orifice must not exceed 90�'
121            write(*,*) '       element number: ',nelem
122            call exit(201)
123          endif
124          if(prop(index+6).lt.0.d0) then
125            write(*,*) '*ERROR in checkinputvaluesnet:'
126            write(*,*) '      d/D (orifice)must not be negative'
127            write(*,*) '       element number: ',nelem
128            call exit(201)
129          endif
130          if(prop(index+6).gt.0.5d0) then
131            write(*,*) '*ERROR in checkinputvaluesnet:'
132            write(*,*) '       d/D (orifice)must not exceed 0.5'
133            write(*,*) '       element number: ',nelem
134            call exit(201)
135          endif
136          if(prop(index+3)/prop(index+2).lt.0.d0) then
137            write(*,*)
138     &           '*ERROR in checkinputvaluesnet:'
139            write(*,*) '        L/d of an orifice ',
140     &           'must not be negative'
141            write(*,*) '       element number: ',nelem
142            call exit(201)
143          endif
144          if(lakon(nelem)(4:4).eq.'P') then
145            if(prop(index+3)/prop(index+2).gt.2.d0) then
146              write(*,*)
147     &             '*WARNING in checkinputvaluesnet: '
148              write(*,*) '       L/d of an orifice ',
149     &             'with Parker must not exceed 2'
150              write(*,*) '       element number: ',nelem
151!     call exit(201)
152            endif
153          endif
154        elseif((lakon(nelem)(2:5).eq.'ORBT'))then
155          if(prop(index+2).lt.0.d0) then
156            write(*,*) '*ERROR in checkinputvaluesnet:'
157            write(*,*) '        ps1/pt1 (bleedtapping) ',
158     &           'must not be negative'
159            write(*,*) '       element number: ',nelem
160            call exit(201)
161          endif
162          if(prop(index+2).gt.1.d0) then
163            write(*,*) '*ERROR in checkinputvaluesnet: '
164            write(*,*) '       ps1/pt1 (bleed tapping) ',
165     &           'must not exceed 1'
166            write(*,*) '       element number: ',nelem
167            call exit(201)
168          endif
169        elseif((lakon(nelem)(2:5).eq.'ORPN'))then
170          if(prop(index+2).lt.0.d0) then
171            write(*,*) '*ERROR in checkinputvaluesnet: '
172            write(*,*) '       theta (preswirlnozzle) ',
173     &           'must not be negative'
174            write(*,*) '       element number: ',nelem
175            call exit(201)
176          endif
177          if(prop(index+2).gt.90.d0) then
178            write(*,*) '*ERROR in checkinputvaluesnet:'
179            write(*,*) '       theta (preswirl nozzle) ',
180     &           'must not exceed 90�'
181            write(*,*) '       element number: ',nelem
182            call exit(201)
183          endif
184          if(prop(index+3).lt.0.d0) then
185            write(*,*) '*ERROR in checkinputvaluesnet: '
186            write(*,*) '       k_phi (preswirl nozzle) ',
187     &           'must not be negative'
188            write(*,*) '       element number: ',nelem
189            call exit(201)
190          endif
191          if(prop(index+3).gt.1.05d0) then
192            write(*,*) '*ERROR in checkinputvaluesnet: '
193            write(*,*) '       k_phi (preswirlnozzle) ',
194     &           'must not exceed 1.05'
195            write(*,*) '       element number: ',nelem
196            call exit(201)
197          endif
198        elseif((lakon(nelem)(2:5).eq.'ORBG'))then
199          if(prop(index+1).lt.0.d0) then
200            write(*,*) '*ERROR in checkinputvaluesnet: '
201            write(*,*) '       section area is not positive'
202            write(*,*) '       element number: ',nelem
203            call exit(201)
204          endif
205          if(prop(index+2).lt.0.d0 .or.
206     &         prop(index+2).ge.1.d0) then
207            write(*,*) '*ERROR in checkinputvaluesnet: '
208            write(*,*) '       using Bragg Method'
209            write(*,*) '       Cd by crtitical pressure ratio '
210            write(*,*) '       *FLUID SECTIONS position 2'
211            write(*,*) '       0 < Cd _crit < 1'
212            write(*,*) '       element number: ',nelem
213            call exit(201)
214          endif
215        elseif(lakon(nelem)(2:5).eq.'ORB1')then
216          if(dabs(1.d0-prop(index+4)/50.d0).lt.0.00001d0.and.
217     &         dabs(1.d0-prop(index+5)/0.00015d0).gt.0.00001d0)then
218            write(*,*) '*ERROR in checkinputvaluesnet: '
219            write(*,*) '       For a brush seal with a bristle'
220            write(*,*) '       density of 50 the bristle diameter'
221            write(*,*) '       has to be 0.00015m!'
222            write(*,*) '       For element number:', nelem
223            write(*,*) '       the value is:',prop(index+5)
224            call exit(201)
225          elseif(dabs(1.d0-prop(index+4)/100.d0).lt.0.00001d0.and.
226     &           dabs(1.d0-prop(index+5)/0.00007d0).gt.0.00001d0)then
227            write(*,*) '*ERROR in checkinputvaluesnet: '
228            write(*,*) '       For a brush seal with a bristle'
229            write(*,*) '       density of 100 the bristle diameter'
230            write(*,*) '       has to be 0.00007m!'
231            write(*,*) '       For element number:', nelem
232            write(*,*) '       the value is:',prop(index+5)
233            call exit(201)
234          elseif(dabs(1.d0-prop(index+4)/140.d0).lt.0.00001d0.and.
235     &           dabs(1.d0-prop(index+5)/0.0001d0).gt.0.00001d0)then
236            write(*,*) '*ERROR in checkinputvaluesnet: '
237            write(*,*) '       For a brush seal with a bristle'
238            write(*,*) '       density of 140 the bristle diameter'
239            write(*,*) '       has to be 0.00010m!'
240            write(*,*) '       For element number:', nelem
241            write(*,*) '       the value is:',prop(index+5)
242            call exit(201)
243          elseif(dabs(1.d0-prop(index+4)/200.d0).lt.0.00001d0.and.
244     &           dabs(1.d0-prop(index+5)/0.00007d0).gt.0.00001d0)then
245            write(*,*) '*ERROR in checkinputvaluesnet: '
246            write(*,*) '       For a brush seal with a bristle'
247            write(*,*) '       density of 200 the bristle diameter'
248            write(*,*) '       has to be 0.00007m!'
249            write(*,*) '       For element number:', nelem
250            write(*,*) '       the value is:',prop(index+5)
251            call exit(201)
252          endif
253          if(dabs(1.d0-prop(index+4)/50.d0).gt.0.00001d0.and.
254     &         dabs(1.d0-prop(index+4)/100.d0).gt.0.00001d0.and.
255     &         dabs(1.d0-prop(index+4)/140.d0).gt.0.00001d0.and.
256     &         dabs(1.d0-prop(index+4)/200.d0).gt.0.00001d0)then
257            write(*,*) '*ERROR in checkinputvaluesnet: '
258            write(*,*) '       Only a brsitle density of'
259            write(*,*) '       50, 100, 140 and 200 is supported'
260            write(*,*) '       with corresponding bristle diameter'
261            write(*,*) '       For element number:', nelem
262            write(*,*) '       the value is:',prop(index+4)
263            call exit(201)
264          endif
265!
266        elseif(lakon(nelem)(2:5).eq.'ORB2')then
267          if(dabs(1.d0-prop(index+3)/50.d0).lt.0.00001d0.and.
268     &         dabs(1.d0-prop(index+4)/0.00015d0).gt.0.00001d0)then
269            write(*,*) '*ERROR in checkinputvaluesnet: '
270            write(*,*) '       For a brush seal with a bristle'
271            write(*,*) '       density of 50 the bristle diameter'
272            write(*,*) '       has to be 0.00015m!'
273            write(*,*) '       For element number:', nelem
274            write(*,*) '       the value is:',prop(index+5)
275            call exit(201)
276          elseif(dabs(1.d0-prop(index+3)/100.d0).lt.0.00001d0.and.
277     &           dabs(1.d0-prop(index+4)/0.00007d0).gt.0.00001d0)then
278            write(*,*) '*ERROR in checkinputvaluesnet: '
279            write(*,*) '       For a brush seal with a bristle'
280            write(*,*) '       density of 100 the bristle diameter'
281            write(*,*) '       has to be 0.00007m!'
282            write(*,*) '       For element number:', nelem
283            write(*,*) '       the value is:',prop(index+5)
284            call exit(201)
285          elseif(dabs(1.d0-prop(index+3)/140.d0).lt.0.00001d0.and.
286     &           dabs(1.d0-prop(index+4)/0.0001d0).gt.0.00001d0)then
287            write(*,*) '*ERROR in checkinputvaluesnet: '
288            write(*,*) '       For a brush seal with a bristle'
289            write(*,*) '       density of 140 the bristle diameter'
290            write(*,*) '       has to be 0.00010m!'
291            write(*,*) '       For element number:', nelem
292            write(*,*) '       the value is:',prop(index+5)
293            call exit(201)
294          elseif(dabs(1.d0-prop(index+3)/200.d0).lt.0.00001d0.and.
295     &           dabs(1.d0-prop(index+4)/0.00007d0).gt.0.00001d0)then
296            write(*,*) '*ERROR in checkinputvaluesnet: '
297            write(*,*) '       For a brush seal with a bristle'
298            write(*,*) '       density of 200 the bristle diameter'
299            write(*,*) '       has to be 0.00007m!'
300            write(*,*) '       For element number:', nelem
301            write(*,*) '       the value is:',prop(index+5)
302            call exit(201)
303          endif
304          if(dabs(1.d0-prop(index+3)/50.d0).gt.0.00001d0.and.
305     &         dabs(1.d0-prop(index+3)/100.d0).gt.0.00001d0.and.
306     &         dabs(1.d0-prop(index+3)/140.d0).gt.0.00001d0.and.
307     &         dabs(1.d0-prop(index+3)/200.d0).gt.0.00001d0)then
308            write(*,*) '*ERROR in checkinputvaluesnet: '
309            write(*,*) '       Only a brsitle density of'
310            write(*,*) '       50, 100, 140 and 200 is supported'
311            write(*,*) '       with corresponding bristle diameter'
312            write(*,*) '       For element number:', nelem
313            write(*,*) '       the value is:',prop(index+4)
314            call exit(201)
315          endif
316!
317        elseif((lakon(nelem)(2:4).eq.'LAB').and.
318     &         (lakon(nelem)(2:5).ne.'LABF').and.
319     &         (lakon(nelem)(2:5).ne.'LABD')) then
320          if((prop(index+1).gt.1000.d0)
321     &         .or.(prop(index+1).lt.0.d0)) then
322            write(*,*)
323     &           '*ERROR in checkinputvaluesnet: ',
324     &           'the selected pitch t'
325            write(*,*) '       for the labyrinth is not correct'
326            write(*,*) '       0<=t<=1000 mm'
327            write(*,*) '       element number: ',nelem
328            call exit(201)
329          endif
330          if((prop(index+2).gt.100.d0)
331     &         .or.(prop(index+2).lt.0.d0)) then
332            write(*,*)
333     &           '*ERROR in checkinputvaluesnet: ',
334     &           'the selected gap s'
335            write(*,*) '       for the labyrinth is not correct'
336            write(*,*) '       0<=s<=100 mm'
337            write(*,*) '       element number: ',nelem
338            call exit(201)
339          endif
340          if((prop(index+3).gt.5000.d0)
341     &         .or.(prop(index+3).lt.0.d0)) then
342            write(*,*) '*ERROR in checkinputvaluesnet:'
343            write(*,*) '       the selected diameter d'
344            write(*,*) '       for the labyrinth is not correct'
345            write(*,*) '       0<=d<=5000'
346            write(*,*) '       element number: ',nelem
347            call exit(201)
348          endif
349          if((prop(index+5).gt.9.d0)
350     &         .or.(prop(index+5).lt.0.d0)) then
351            write(*,*) '*ERROR in checkinputvaluesnet:'
352            write(*,*)'       the selected spike number n'
353            write(*,*) '       for the labyrinth is not correct'
354            write(*,*) '       0<=n<=9'
355            write(*,*) prop(index+4)
356            write(*,*) '       element number: ',nelem
357            call exit(201)
358          endif
359          if((prop(index+5).gt.100.d0)
360     &         .or.(prop(index+5).lt.0.d0)) then
361            write(*,*) '*ERROR in checkinputvaluesnet:'
362            write(*,*) '       the selected spike breadth'
363            write(*,*) '       for the labyrinth is not correct'
364            write(*,*) '       0<=b<=100 mm'
365            write(*,*) '       element number: ',nelem
366            call exit(201)
367          endif
368          if((prop(index+6).gt.9.d0)
369     &         .or.(prop(index+6).lt.0.d0)) then
370            write(*,*) '*ERROR in checkinputvaluesnet:'
371            write(*,*) '       the selected spike height'
372            write(*,*) '       for the labyrinth is not correct'
373            write(*,*) '       0<=b<=20 mm'
374            write(*,*) '       element number: ',nelem
375            call exit(201)
376          endif
377          if((prop(index+7).gt.4.d0)
378     &         .or.(prop(index+7).lt.0.d0)) then
379            write(*,*) '*ERROR in checkinputvaluesnet:'
380            write(*,*) '       the selected Honeycomb cell width'
381            write(*,*) '       for the labyrinth is not correct'
382            write(*,*) '       0<=L<=4 mm'
383            write(*,*) '       element number: ',nelem
384            call exit(201)
385          endif
386          if((prop(index+8).gt.5.d0)
387     &         .or.(prop(index+8).lt.0.d0)) then
388            write(*,*) '*ERROR in checkinputvaluesnet:'
389            write(*,*) '       the selected edge radius'
390            write(*,*) '       for the labyrinth is not correct'
391            write(*,*) '       0<=r<=5 mm'
392            write(*,*) '       element number: ',nelem
393            call exit(201)
394          endif
395          if((prop(index+9).gt.100.d0)
396     &         .or.(prop(index+9).lt.0.d0)) then
397            write(*,*) '*ERROR in checkinputvaluesnet:'
398            write(*,*) '       the selected position of the spike'
399            write(*,*) '       for the labyrinth is not correct'
400            write(*,*) '       0<=X<=0.1 mm'
401            write(*,*) '       element number: ',nelem
402            call exit(201)
403          endif
404          if((prop(index+10).gt.100.d0)
405     &         .or.(prop(index+10).lt.0.d0)) then
406            write(*,*) '*ERROR in checkinputvaluesnet:'
407            write(*,*) '       the selected height of the step'
408            write(*,*) '       for the labyrinth is not correct'
409            write(*,*) '       0<=Hst<=0.1 mm'
410            write(*,*) '       element number: ',nelem
411            call exit(201)
412          endif
413!
414        elseif((lakon(nelem)(2:6).eq.'CARBS'))then
415          if(lakon(nelem)(2:8).eq.'CARBSGE') then
416            if(prop(index+1).le.0.d0) then
417              write(*,*) '*ERROR in checkinputvaluesnet:'
418              write(*,*) '       the selected diameter of ',
419     &             'the carbon seal'
420              write(*,*)
421     &             '       has been defined as less or equal to 0'
422              write(*,*) '       element number: ',nelem
423              call exit(201)
424            endif
425          else
426            if((prop(index+1).le.0.d0)
427     &           .or.(prop(index+2).le.0.d0)
428     &           .or.(prop(index+3).le.0.d0)) then
429              write(*,*) '*ERROR in checkinputvaluesnet:'
430              write(*,*) '       the selected diameter'
431              write(*,*) '       or the selected length'
432              write(*,*)
433     &             '       or the selected gap of the carbon seal'
434              write(*,*)
435     &             '       has been defined as less or equal to 0'
436              write(*,*) '       element number: ',nelem
437              call exit(201)
438            endif
439          endif
440!
441        elseif((lakon(nelem)(2:4).eq.'RCVL'))then
442          if(prop(index+2).lt.(prop(index+3))) then
443            write(*,*) '*ERROR in checkinputvaluesnet: '
444            write(*,*) '       element TYPE=ROTATING CAVITY ',
445     &           '(Radial inflow)'
446            write(*,*) '       the specified upstream radius is ',
447     &           'smaller than'
448            write(*,*) '       the specified downstream radius!'
449            write(*,*) '       Please check the element definition.'
450            write(*,*) '       element number: ',nelem
451            call exit(201)
452          endif
453        elseif((lakon(nelem)(2:4).eq.'RCVN'))then
454          if(prop(index+1).lt.(prop(index+2))) then
455            write(*,*) '*ERROR in checkinputvaluesnet: '
456            write(*,*) '       element TYPE=ROTATING CAVITY ',
457     &           '(Radial inflow)'
458            write(*,*) '       the specified upstream radius is ',
459     &           'smaller than'
460            write(*,*) '       the specified downstream radius!'
461            write(*,*) '       Please check the element definition.'
462            write(*,*) '       element number: ',nelem
463            call exit(201)
464          endif
465!
466        elseif(lakon(nelem)(2:3).eq.'RE') then
467          if(((prop(index+1).le.0.d0)
468     &         .or.(prop(index+2).le.0.d0)
469     &         .or.(prop(index+3).le.0.d0))
470     &         .and.(lakon(nelem)(2:5).ne.'REBR')
471     &         .and.(lakon(nelem)(2:5).ne.'REEX')
472     &         .and.(lakon(nelem)(2:7).ne.'REWAOR')
473     &         .and.(lakon(nelem)(2:5).ne.'REEN'))then
474            write(*,*) '*ERROR in checkinputvaluesnet:'
475            write(*,*) '       A1,A2 or Dh less or equal 0'
476            write(*,*) '       element number: ',nelem
477            call exit(201)
478!
479          elseif((lakon(nelem)(2:5).eq.'REEL'))then
480            if(prop(index+1).ge.(prop(index+2))) then
481              write(*,*) '*ERROR in checkinputvaluesnet:'
482              write(*,*)
483     &             '       Section A1 is greater than section A2'
484              write(*,*) '       element number: ',nelem
485              call exit(201)
486            endif
487!
488          elseif((lakon(nelem)(2:5).eq.'RECO'))then
489            if(prop(index+1).lt.(prop(index+2))) then
490              write(*,*) '*ERROR in checkinputvaluesnet:'
491              write(*,*)
492     &             '       Section A2 is greater than section A1'
493              write(*,*) '       element number: ',nelem
494              call exit(201)
495            endif
496          endif
497!
498          if((lakon(nelem)(2:5).eq.'REBR'))then
499            if((prop(index+1).le.0.d0)
500     &           .or.(prop(index+2).le.0.d0)
501     &           .or.(prop(index+3).le.0.d0))then
502              write(*,*) '*ERROR in checkinputvaluesnet:'
503              write(*,*) '       trying to define a branch '
504              write(*,*) '       all three elements must be ',
505     &             'different from 0'
506              write(*,*) '       element number: ',nelem
507              call exit(201)
508!
509            elseif((prop(index+4).le.0.d0)
510     &             .or.(prop(index+5).le.0.d0)
511     &             .or.(prop(index+6).le.0.d0))then
512              write(*,*) '*ERROR in checkinputvaluesnet:'
513              write(*,*) '       trying to define a branch '
514              write(*,*) '       all sections must be positive'
515              write(*,*) '       element number: ',nelem
516              call exit(201)
517!
518            elseif((prop(index+7).lt.0)
519     &             .or.(prop(index+8).lt.0))then
520              write(*,*) '*ERROR in checkinputvaluesnet:'
521              write(*,*) '       trying to define a branch '
522              write(*,*) '       alpha1 & 2 cannot be negative'
523              write(*,*) '       element number: ',nelem
524              call exit(201)
525!
526            elseif((prop(index+7).gt.90)
527     &             .or.(prop(index+8).gt.90)) then
528              write(*,*) '*ERROR in checkinputvaluesnet:'
529              write(*,*) '       trying to define a branch '
530              write(*,*) '       alpha1 & 2 cannot greater than ',
531     &             '90 gegrees'
532              write(*,*) '       element number: ',nelem
533              call exit(201)
534!
535            elseif((lakon(nelem)(6:8).eq.'SI1')
536     &             .or.(lakon(nelem)(6:8).eq.'JI1')
537     &             .or.(lakon(nelem)(6:8).eq.'JI2')) then
538              if(prop(index+7).gt.0) then
539                write(*,*) '*ERROR in checkinputvaluesnet:'
540                write(*,*) '       trying to define a branch '
541                write(*,*)
542     &               '       Type IDELCHIK JOINT1 or SPLIT1&2'
543                write(*,*) '       alpha1 must be 0 degrees'
544                write(*,*) '       element number: ',nelem
545                call exit(201)
546              endif
547!
548            elseif((lakon(nelem)(6:8).eq.'SI1').or.
549     &             (lakon(nelem)(6:8).eq.'JI1'))then
550              if(prop(index+4).ne.(prop(index+5)))then
551                write(*,*) '*ERROR in checkinputvaluesnet:'
552                write(*,*) '       trying to define a branch '
553                write(*,*) '       Type IDELCHIK SPLIT1 or JOINT1 '
554                write(*,*) '       A1=A0'
555                write(*,*) '       element number: ',nelem
556                call exit(201)
557              endif
558!
559            elseif(lakon(nelem)(6:8).eq.'JI2') then
560              if((prop(index+5)+(prop(index+6)))
561     &             .ne.prop(index+4))then
562                write(*,*) '*ERROR in checkinputvaluesnet:'
563                write(*,*) '       trying to define a branch '
564                write(*,*) '       Type IDELCHIK JOINT2 '
565                write(*,*) '       A1+A2 must be equal to A0'
566                write(*,*) '       element number: ',nelem
567                call exit(201)
568              endif
569            endif
570          endif
571!
572!     General Vortex
573!
574        elseif((lakon(nelem)(2:3).eq.'VO'))then
575!
576!     inner and outer radius less or equal to 0
577!
578          if((prop(index+1).le.0.d0) .or.
579     &         (prop(index+2).le.0.d0)) then
580            write(*,*)'*ERROR in checkinputvaluesnet:'
581            write(*,*)'       trying to define a VORTEX'
582            write(*,*)'       R1 and R2 must be positive'
583            write(*,*) '       element number: ',nelem
584            call exit(201)
585!
586          elseif(prop(index+3).le.0.d0) then
587            write(*,*)'*ERROR in checkinputvaluesnet:'
588            write(*,*)'       trying to define a VORTEX'
589            write(*,*)'       eta must be different positive'
590            write(*,*) '       element number: ',nelem
591            call exit(201)
592          endif
593!
594!     FREE VORTEX
595!
596          if((lakon(nelem)(2:5).eq.'VOFR'))then
597!
598!     the swirl comes from another upstream  element
599!
600            if(prop(index+6).ne.0d0) then
601!
602!     the rotation speed must be 0
603!
604              if(prop(index+7).ne.0.d0) then
605                write(*,*)'*ERROR in checkinputvaluesnet:'
606                write(*,*)'       trying to define a FREE VORTEX'
607                write(*,*)
608     &               '       rotation speed and upstream element'
609                write(*,*)'       cannot be simultaneously used '
610                write(*,*) '       element number: ',nelem
611                call exit(201)
612              endif
613            endif
614!
615!     FORCED VORTEX
616!
617          elseif((lakon(nelem)(2:5).eq.'VOFO'))then
618!
619!     Core swirl ratio must be defined and positive
620!
621            if(prop(index+4).le.0.d0) then
622              write(*,*)'*ERROR in checkinputvaluesnet:'
623              write(*,*)'       trying to define a FORCED VORTEX'
624              write(*,*)'       Core swirl ratio Kr is strictly ',
625     &             'positive'
626              write(*,*) '       element number: ',nelem
627              call exit(201)
628            endif
629            if(prop(index+4).gt.1.d0) then
630              write(*,*)'*WARNING in checkinputvaluesnet:'
631              write(*,*)'       trying to define a FORCED VORTEX'
632              write(*,*)'       Core swirl ratio Kr is ',
633     &             'greater than 1'
634              write(*,*) '       element number: ',nelem
635!     call exit(201)
636            endif
637!
638!     Rotation speed must be defined and positive
639!
640            if(prop(index+5).le.0.d0) then
641              write(*,*)'*ERROR in checkinputvaluesnet:'
642              write(*,*)'       trying to define a FORCED VORTEX'
643              write(*,*)
644     &             '       Rotation speed n is strictly positive'
645              write(*,*) '       element number: ',nelem
646              call exit(201)
647            endif
648          endif
649!
650!     Absolute/relative system
651!
652        elseif((lakon(nelem)(2:4).eq.'ATR').or.
653     &         (lakon(nelem)(2:4).eq.'RTA')) then
654          if(prop(index+1).le.0.d0) then
655            write(*,*)'*ERROR in checkinputvaluesnet:'
656            write(*,*)'       trying to define an element'
657            write(*,*)'       TYPE= ABSOLUTE TO RELATIVE or'
658            write(*,*)'       TYPE= RELATIVE TO ABSOLUTE'
659            write(*,*)'       Rotation speed is strictly positive'
660            write(*,*)'       element number: ',nelem
661            call exit(201)
662          elseif(prop(index+3).ne.0.d0) then
663            if(prop(index+2).ne.0.d0) then
664              write(*,*)'*ERROR in checkinputvaluesnet:'
665              write(*,*)'       trying to define an element'
666              write(*,*)'       TYPE= ABSOLUTE TO RELATIVE or'
667              write(*,*)'       TYPE= RELATIVE TO ABSOLUTE'
668              write(*,*)'       reference element has been provided'
669              write(*,*)'       but tangential velocity ',
670     &             'is already defined'
671              write(*,*)'       element number: ',nelem
672              call exit(201)
673            endif
674          endif
675!
676!     Air Valve
677!
678        elseif(lakon(nelem)(2:5).eq.'AVLV')then
679!
680          if((dabs(prop(index+1)-1.d0).gt.1.d-8).and.
681     &         (dabs(prop(index+1)-1.25d0).gt.1.d-8).and.
682     &         (dabs(prop(index+1)-1.5d0).gt.1.d-8).and.
683     &         (dabs(prop(index+1)-2.d0).gt.1.d-8).and.
684     &         (dabs(prop(index+1)-2.5d0).gt.1.d-8).and.
685     &         (dabs(prop(index+1)-3.d0).gt.1.d-8).and.
686     &         (dabs(prop(index+1)-4.d0).gt.1.d-8).and.
687     &         (dabs(prop(index+1)-6.d0).gt.1.d-8).and.
688     &         (dabs(prop(index+1)-6.5d0).gt.1.d-8).and.
689     &         (dabs(prop(index+1)-8.d0).gt.1.d-8).and.
690     &         nint(prop(index+3)).eq.1)then
691!
692            write(*,*) '*ERROR in air_valve: '
693            write(*,*) '  Specified diameter is: ',
694     &           prop(index+1)
695            write(*,*) '  inch(s) diameter should be a value'
696            write(*,*) '  in inch in1, 1.25, 1.5, 2, 2.5, 3,'
697            write(*,*) '  4, 6, 6.5, 8,'
698            call exit(201)
699!
700          endif
701!
702          if((prop(index+1).gt.90d0).or.
703     &         (prop(index+1).le.0)) then
704!
705            write(*,*) '*ERROR in air_valve: '
706            write(*,*) '  valve opening should be a value'
707            write(*,*) '  comprised between 0 grad and 90 grad'
708            call exit(201)
709!
710          endif
711!
712        elseif(lakon(nelem)(2:5).eq.'GAPF') then
713          if(prop(index+5).eq.0.d0) then
714            write(*,*) '*ERROR in checkinputvaluesnet: form factor'
715            write(*,*) '       is equal to zero'
716            write(*,*) '       element number:',nelem
717            call exit(201)
718          endif
719!
720        elseif((lakon(nelem)(2:5).ne.'LIPU').and.
721     &         (lakon(nelem).ne.'       ').and.
722     &         (lakon(nelem)(2:4).ne.'LAB').and.
723     &         (lakon(nelem)(2:5).ne.'CHAR').and.
724     &         (lakon(nelem)(2:6).ne.'CARBS'))then
725          if((prop(index+1).lt.0.d0)) then
726            write(*,*) '*ERROR in checkinputvaluesnet: section area'
727            write(*,*) '       is not positive'
728            write(*,*) '       element number: ',nelem
729            call exit(201)
730          endif
731!
732        endif
733!
734      enddo
735!
736      return
737      end
738