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

Go to the source code of this file.

Functions/Subroutines

subroutine checkinputvaluesnet (ieg, nflow, prop, ielprop, lakon)
 

Function/Subroutine Documentation

◆ checkinputvaluesnet()

subroutine checkinputvaluesnet ( integer, dimension(*)  ieg,
integer  nflow,
real*8, dimension(*)  prop,
integer, dimension(*)  ielprop,
character*8, dimension(*)  lakon 
)
31 !
32  implicit none
33 !
34  character*8 lakon(*)
35 !
36  integer ieg(*),nflow,i,ielprop(*),index,nelem
37 !
38  real*8 prop(*)
39 !
40  do i=1,nflow
41  nelem=ieg(i)
42  index=ielprop(nelem)
43  if(index.lt.0) cycle
44 !
45 ! modifying the prop array (formerly done in fluidsections.f)
46 !
47  if(lakon(nelem)(2:8).eq.'REBRJI2') then
48  if(1.d0-(prop(index+5)+prop(index+6))/
49  & prop(index+4).gt.0.01d0)then
50  write(*,*) '*ERROR in checkinputvaluesnet:'
51  write(*,*) ' in element type RESTRICTOR
52  & BRANCH JOINT IDELCHIK2'
53  write(*,*) ' A0 ist not equal to A1+A2'
54  write(*,*) ' element number: ',nelem
55  call exit(201)
56  endif
57 !
58  elseif((lakon(nelem)(2:3).eq.'OR').and.
59  & (lakon(nelem)(2:5).ne.'ORC1').and.
60  & (lakon(nelem)(2:5).ne.'ORBT').and.
61  & (lakon(nelem)(2:5).ne.'ORPN').and.
62  & (lakon(nelem)(2:5).ne.'ORFL'))then
63  if(prop(index+2).lt.0.d0) then
64  write(*,*) '*ERROR in checkinputvaluesnet: diameter '
65  write(*,*) ' of the orifice is not positive'
66  write(*,*) ' element number: ',nelem
67  call exit(201)
68  endif
69  if(prop(index+3).lt.0.d0) then
70  write(*,*) '*ERROR in checkinputvaluesnet:'
71  write(*,*) ' length of the orifice is
72  & not positive'
73  write(*,*) ' element number: ',nelem
74  call exit(201)
75  endif
76  if((lakon(nelem)(2:5).ne.'ORC1').and.
77  & (lakon(nelem)(2:5).ne.'ORBF').and.
78  & (lakon(nelem)(2:5).ne.'ORRG').and.
79  & (lakon(nelem)(2:5).ne.'ORSG').and.
80  & (lakon(nelem)(2:5).ne.'ORGA').and.
81  & (lakon(nelem)(2:5).ne.'ORBO'))then
82  if((prop(index+4).gt.1.d-20).and.
83  & (prop(index+5).gt.1.d-20)) then
84  write(*,*)
85  & '*ERROR in checkinputvaluesnet:'
86  write(*,*) 'either the radius of
87  & the orifice must be zero or the chamfer angle'
88  write(*,*) ' element number: ',nelem
89  call exit(201)
90  endif
91  endif
92  if(prop(index+4)/prop(index+2).lt.0.d0) then
93  write(*,*) '*ERROR in checkinputvaluesnet: '
94  write(*,*)
95  & ' r/d of an orifice must not be negative'
96  write(*,*) ' element number: ',nelem
97  call exit(201)
98  endif
99  if(prop(index+4)/prop(index+2).gt.0.82d0) then
100  write(*,*)
101  & '*ERROR in checkinputvaluesnet: '
102  write(*,*) ' r/d of an orifice
103  & must not not exceed 0.82'
104  write(*,*) ' element number: ',nelem
105  call exit(201)
106  endif
107  if(prop(index+5).lt.0.d0) then
108  write(*,*) '*ERROR in checkinputvaluesnet:'
109  write(*,*) ' the chamfer angle of an
110  & orifice must not be negative'
111  write(*,*) ' element number: ',nelem
112  call exit(201)
113  endif
114  if(prop(index+5).gt.90.d0) then
115  write(*,*)
116  & '*ERROR in checkinputvaluesnet: '
117  write(*,*) ' the chamfer angle of an
118  & orifice must not exceed 90°'
119  write(*,*) ' element number: ',nelem
120  call exit(201)
121  endif
122  if(prop(index+6).lt.0.d0) then
123  write(*,*) '*ERROR in checkinputvaluesnet:'
124  write(*,*) ' d/D (orifice)must not be negative'
125  write(*,*) ' element number: ',nelem
126  call exit(201)
127  endif
128  if(prop(index+6).gt.0.5d0) then
129  write(*,*) '*ERROR in checkinputvaluesnet:'
130  write(*,*) ' d/D (orifice)must not exceed 0.5'
131  write(*,*) ' element number: ',nelem
132  call exit(201)
133  endif
134  if(prop(index+3)/prop(index+2).lt.0.d0) then
135  write(*,*)
136  & '*ERROR in checkinputvaluesnet:'
137  write(*,*) ' L/d of an orifice
138  & must not be negative'
139  write(*,*) ' element number: ',nelem
140  call exit(201)
141  endif
142  if(lakon(nelem)(4:4).eq.'P') then
143  if(prop(index+3)/prop(index+2).gt.2.d0) then
144  write(*,*)
145  & '*ERROR in checkinputvaluesnet: '
146  write(*,*) ' L/d of an orifice
147  & with Parker must not exceed 2'
148  write(*,*) ' element number: ',nelem
149  call exit(201)
150  endif
151  endif
152  elseif((lakon(nelem)(2:5).eq.'ORBT'))then
153  if(prop(index+2).lt.0.d0) then
154  write(*,*) '*ERROR in checkinputvaluesnet:'
155  write(*,*) ' ps1/pt1 (bleedtapping)
156  & must not be negative'
157  write(*,*) ' element number: ',nelem
158  call exit(201)
159  endif
160  if(prop(index+2).gt.1.d0) then
161  write(*,*) '*ERROR in checkinputvaluesnet: '
162  write(*,*) ' ps1/pt1 (bleed tapping)
163  & must not exceed 1'
164  write(*,*) ' element number: ',nelem
165  call exit(201)
166  endif
167  elseif((lakon(nelem)(2:5).eq.'ORPN'))then
168  if(prop(index+2).lt.0.d0) then
169  write(*,*) '*ERROR in checkinputvaluesnet: '
170  write(*,*) ' theta (preswirlnozzle)
171  & Nmust not be negative'
172  write(*,*) ' element number: ',nelem
173  call exit(201)
174  endif
175  if(prop(index+2).gt.90.d0) then
176  write(*,*) '*ERROR in checkinputvaluesnet:'
177  write(*,*) ' theta (preswirl nozzle)
178  & must not exceed 90°'
179  write(*,*) ' element number: ',nelem
180  call exit(201)
181  endif
182  if(prop(index+3).lt.0.d0) then
183  write(*,*) '*ERROR in checkinputvaluesnet: '
184  write(*,*) ' k_phi (preswirl nozzle)
185  & must not be negative'
186  write(*,*) ' element number: ',nelem
187  call exit(201)
188  endif
189  if(prop(index+3).gt.1.05d0) then
190  write(*,*) '*ERROR in checkinputvaluesnet: '
191  write(*,*) ' k_phi (preswirlnozzle)
192  & must not exceed 1.05'
193  write(*,*) ' element number: ',nelem
194  call exit(201)
195  endif
196  elseif((lakon(nelem)(2:5).eq.'ORBG'))then
197  if(prop(index+1).lt.0.d0) then
198  write(*,*) '*ERROR in checkinputvaluesnet: '
199  write(*,*) ' section area is not positive'
200  write(*,*) ' element number: ',nelem
201  call exit(201)
202  endif
203  if(prop(index+2).lt.0.d0 .or.
204  & prop(index+2).ge.1.d0) then
205  write(*,*) '*ERROR in checkinputvaluesnet: '
206  write(*,*) ' using Bragg Method'
207  write(*,*) ' Cd by crtitical pressure ratio '
208  write(*,*) ' *FLUID SECTIONS position 2'
209  write(*,*) ' 0 < Cd _crit < 1'
210  write(*,*) ' element number: ',nelem
211  call exit(201)
212  endif
213 !
214  elseif((lakon(nelem)(2:4).eq.'LAB').and.
215  & (lakon(nelem)(2:5).ne.'LABF').and.
216  & (lakon(nelem)(2:5).ne.'LABD')) then
217  if((prop(index+1).gt.1000.d0)
218  & .or.(prop(index+1).lt.0.d0)) then
219  write(*,*)
220  & '*ERROR in checkinputvaluesnet:
221  & the selected pitch t'
222  write(*,*) ' for the labyrinth is not correct'
223  write(*,*) ' 0<=t<=1000 mm'
224  write(*,*) ' element number: ',nelem
225  call exit(201)
226  endif
227  if((prop(index+2).gt.100d0)
228  & .or.(prop(index+2).lt.0.d0)) then
229  write(*,*)
230  & '*ERROR in checkinputvaluesnet:
231  & the selected gap s'
232  write(*,*) ' for the labyrinth is not correct'
233  write(*,*) ' 0<=s<=100 mm'
234  write(*,*) ' element number: ',nelem
235  call exit(201)
236  endif
237  if((prop(index+3).gt.5000.d0)
238  & .or.(prop(index+3).lt.0.d0)) then
239  write(*,*) '*ERROR in checkinputvaluesnet:'
240  write(*,*) ' the selected diameter d'
241  write(*,*) ' for the labyrinth is not correct'
242  write(*,*) ' 0<=d<=5000'
243  write(*,*) ' element number: ',nelem
244  call exit(201)
245  endif
246  if((prop(index+5).gt.9.d0)
247  & .or.(prop(index+5).lt.0.d0)) then
248  write(*,*) '*ERROR in checkinputvaluesnet:'
249  write(*,*)' the selected spike number n'
250  write(*,*) ' for the labyrinth is not correct'
251  write(*,*) ' 0<=n<=9'
252  write(*,*) prop(index+4)
253  write(*,*) ' element number: ',nelem
254  call exit(201)
255  endif
256  if((prop(index+5).gt.100.d0)
257  & .or.(prop(index+5).lt.0.d0)) then
258  write(*,*) '*ERROR in checkinputvaluesnet:'
259  write(*,*) ' the selected spike breadth'
260  write(*,*) ' for the labyrinth is not correct'
261  write(*,*) ' 0<=b<=100 mm'
262  write(*,*) ' element number: ',nelem
263  call exit(201)
264  endif
265  if((prop(index+6).gt.9.d0)
266  & .or.(prop(index+6).lt.0.d0)) then
267  write(*,*) '*ERROR in checkinputvaluesnet:'
268  write(*,*) ' the selected spike height'
269  write(*,*) ' for the labyrinth is not correct'
270  write(*,*) ' 0<=b<=20 mm'
271  write(*,*) ' element number: ',nelem
272  call exit(201)
273  endif
274  if((prop(index+7).gt.4.d0)
275  & .or.(prop(index+7).lt.0.d0)) then
276  write(*,*) '*ERROR in checkinputvaluesnet:'
277  write(*,*) ' the selected Honeycomb cell width'
278  write(*,*) ' for the labyrinth is not correct'
279  write(*,*) ' 0<=L<=4 mm'
280  write(*,*) ' element number: ',nelem
281  call exit(201)
282  endif
283  if((prop(index+8).gt.5.d0)
284  & .or.(prop(index+8).lt.0.d0)) then
285  write(*,*) '*ERROR in checkinputvaluesnet:'
286  write(*,*) ' the selected edge radius'
287  write(*,*) ' for the labyrinth is not correct'
288  write(*,*) ' 0<=r<=5 mm'
289  write(*,*) ' element number: ',nelem
290  call exit(201)
291  endif
292  if((prop(index+9).gt.100.d0)
293  & .or.(prop(index+9).lt.0.d0)) then
294  write(*,*) '*ERROR in checkinputvaluesnet:'
295  write(*,*) ' the selected position of the spike'
296  write(*,*) ' for the labyrinth is not correct'
297  write(*,*) ' 0<=X<=0.1 mm'
298  write(*,*) ' element number: ',nelem
299  call exit(201)
300  endif
301  if((prop(index+10).gt.100.d0)
302  & .or.(prop(index+10).lt.0.d0)) then
303  write(*,*) '*ERROR in checkinputvaluesnet:'
304  write(*,*) ' the selected height of the step'
305  write(*,*) ' for the labyrinth is not correct'
306  write(*,*) ' 0<=Hst<=0.1 mm'
307  write(*,*) ' element number: ',nelem
308  call exit(201)
309  endif
310 !
311  elseif((lakon(nelem)(2:6).eq.'CARBS'))then
312  if(lakon(nelem)(2:8).eq.'CARBSGE') then
313  if(prop(index+1).le.0.d0) then
314  write(*,*) '*ERROR in checkinputvaluesnet:'
315  write(*,*) ' the selected diameter of
316  & the carbon seal'
317  write(*,*)
318  & ' has been defined as less or equal to 0'
319  write(*,*) ' element number: ',nelem
320  call exit(201)
321  endif
322  else
323  if((prop(index+1).le.0.d0)
324  & .or.(prop(index+2).le.0.d0)
325  & .or.(prop(index+3).le.0.d0)) then
326  write(*,*) '*ERROR in checkinputvaluesnet:'
327  write(*,*) ' the selected diameter'
328  write(*,*) ' or the selected length'
329  write(*,*)
330  & ' or the selected gap of the carbon seal'
331  write(*,*)
332  & ' has been defined as less or equal to 0'
333  write(*,*) ' element number: ',nelem
334  call exit(201)
335  endif
336  endif
337 !
338  elseif((lakon(nelem)(2:4).eq.'RCVL'))then
339  if(prop(index+2).lt.(prop(index+3))) then
340  write(*,*) '*ERROR in checkinputvaluesnet: '
341  write(*,*) ' element TYPE=ROTATING CAVITY
342  & (Radial inflow)'
343  write(*,*) ' the specified upstream radius is smaller
344  & than'
345  write(*,*) ' the specified downstream radius!'
346  write(*,*) ' Please check the element definition.'
347  write(*,*) ' element number: ',nelem
348  call exit(201)
349  endif
350  elseif((lakon(nelem)(2:4).eq.'RCVN'))then
351  if(prop(index+1).lt.(prop(index+2))) then
352  write(*,*) '*ERROR in checkinputvaluesnet: '
353  write(*,*) ' element TYPE=ROTATING CAVITY
354  & (Radial inflow)'
355  write(*,*) ' the specified upstream radius is smaller
356  & than'
357  write(*,*) ' the specified downstream radius!'
358  write(*,*) ' Please check the element definition.'
359  write(*,*) ' element number: ',nelem
360  call exit(201)
361  endif
362 !
363  elseif(lakon(nelem)(2:3).eq.'RE') then
364  if(((prop(index+1).le.0.d0)
365  & .or.(prop(index+2).le.0.d0)
366  & .or.(prop(index+3).le.0.d0))
367  & .and.(lakon(nelem)(2:5).ne.'REBR')
368  & .and.(lakon(nelem)(2:5).ne.'REEX')
369  & .and.(lakon(nelem)(2:7).ne.'REWAOR')
370  & .and.(lakon(nelem)(2:5).ne.'REEN'))then
371  write(*,*) '*ERROR in checkinputvaluesnet:'
372  write(*,*) ' A1,A2 or Dh less or equal 0'
373  write(*,*) ' element number: ',nelem
374  call exit(201)
375 !
376  elseif((lakon(nelem)(2:5).eq.'REEL'))then
377  if(prop(index+1).ge.(prop(index+2))) then
378  write(*,*) '*ERROR in checkinputvaluesnet:'
379  write(*,*)
380  & ' Section A1 is greater than section A2'
381  write(*,*) ' element number: ',nelem
382  call exit(201)
383  endif
384 !
385  elseif((lakon(nelem)(2:5).eq.'RECO'))then
386  if(prop(index+1).lt.(prop(index+2))) then
387  write(*,*) '*ERROR in checkinputvaluesnet:'
388  write(*,*)
389  & ' Section A2 is greater than section A1'
390  write(*,*) ' element number: ',nelem
391  call exit(201)
392  endif
393  endif
394 !
395  if((lakon(nelem)(2:5).eq.'REBR'))then
396  if((prop(index+1).le.0.d0)
397  & .or.(prop(index+2).le.0.d0)
398  & .or.(prop(index+3).le.0.d0))then
399  write(*,*) '*ERROR in checkinputvaluesnet:'
400  write(*,*) ' trying to define a branch '
401  write(*,*) ' all three elements must be different
402  & from 0'
403  write(*,*) ' element number: ',nelem
404  call exit(201)
405 !
406  elseif((prop(index+4).le.0.d0)
407  & .or.(prop(index+5).le.0.d0)
408  & .or.(prop(index+6).le.0.d0))then
409  write(*,*) '*ERROR in checkinputvaluesnet:'
410  write(*,*) ' trying to define a branch '
411  write(*,*) ' all sections must be positive'
412  write(*,*) ' element number: ',nelem
413  call exit(201)
414 !
415  elseif((prop(index+7).lt.0)
416  & .or.(prop(index+8).lt.0))then
417  write(*,*) '*ERROR in checkinputvaluesnet:'
418  write(*,*) ' trying to define a branch '
419  write(*,*) ' alpha1 & 2 cannot be negative'
420  write(*,*) ' element number: ',nelem
421  call exit(201)
422 !
423  elseif((prop(index+7).gt.90)
424  & .or.(prop(index+8).gt.90)) then
425  write(*,*) '*ERROR in checkinputvaluesnet:'
426  write(*,*) ' trying to define a branch '
427  write(*,*) ' alpha1 & 2 cannot greater than
428  & 90 gegrees'
429  write(*,*) ' element number: ',nelem
430  call exit(201)
431 !
432  elseif((lakon(nelem)(6:8).eq.'SI1')
433  & .or.(lakon(nelem)(6:8).eq.'JI1')
434  & .or.(lakon(nelem)(6:8).eq.'JI2')) then
435  if(prop(index+7).gt.0) then
436  write(*,*) '*ERROR in checkinputvaluesnet:'
437  write(*,*) ' trying to define a branch '
438  write(*,*)
439  & ' Type IDELCHIK JOINT1 or SPLIT1&2'
440  write(*,*) ' alpha1 must be 0 degrees'
441  write(*,*) ' element number: ',nelem
442  call exit(201)
443  endif
444 !
445  elseif((lakon(nelem)(6:8).eq.'SI1').or.
446  & (lakon(nelem)(6:8).eq.'JI1'))then
447  if(prop(index+4).ne.(prop(index+5)))then
448  write(*,*) '*ERROR in checkinputvaluesnet:'
449  write(*,*) ' trying to define a branch '
450  write(*,*) ' Type IDELCHIK SPLIT1 or JOINT1 '
451  write(*,*) ' A1=A0'
452  write(*,*) ' element number: ',nelem
453  call exit(201)
454  endif
455 !
456  elseif(lakon(nelem)(6:8).eq.'JI2') then
457  if((prop(index+5)+(prop(index+6)))
458  & .ne.prop(index+4))then
459  write(*,*) '*ERROR in checkinputvaluesnet:'
460  write(*,*) ' trying to define a branch '
461  write(*,*) ' Type IDELCHIK JOINT2 '
462  write(*,*) ' A1+A2 must be equal to A0'
463  write(*,*) ' element number: ',nelem
464  call exit(201)
465  endif
466  endif
467  endif
468 !
469 ! General Vortex
470 !
471  elseif((lakon(nelem)(2:3).eq.'VO'))then
472 !
473 ! inner and outer radius less or equal to 0
474 !
475  if((prop(index+1).le.0.d0) .or.
476  & (prop(index+2).le.0.d0)) then
477  write(*,*)'*ERROR in checkinputvaluesnet:'
478  write(*,*)' trying to define a VORTEX'
479  write(*,*)' R1 and R2 must be positive'
480  write(*,*) ' element number: ',nelem
481  call exit(201)
482 !
483  elseif(prop(index+3).le.0.d0) then
484  write(*,*)'*ERROR in checkinputvaluesnet:'
485  write(*,*)' trying to define a VORTEX'
486  write(*,*)' eta must be different positive'
487  write(*,*) ' element number: ',nelem
488  call exit(201)
489  endif
490 !
491 ! FREE VORTEX
492 !
493  if((lakon(nelem)(2:5).eq.'VOFR'))then
494 !
495 ! the swirl comes from another upstream element
496 !
497  if(prop(index+6).ne.0d0) then
498 !
499 ! the rotation speed must be 0
500 !
501  if(prop(index+7).ne.0.d0) then
502  write(*,*)'*ERROR in checkinputvaluesnet:'
503  write(*,*)' trying to define a FREE VORTEX'
504  write(*,*)
505  & ' rotation speed and upstream element'
506  write(*,*)' cannot be simultaneously used '
507  write(*,*) ' element number: ',nelem
508  call exit(201)
509  endif
510  endif
511 !
512 ! FORCED VORTEX
513 !
514  elseif((lakon(nelem)(2:5).eq.'VOFO'))then
515 !
516 ! Core swirl ratio must be defined and positive
517 !
518  if(prop(index+4).le.0.d0) then
519  write(*,*)'*ERROR in checkinputvaluesnet:'
520  write(*,*)' trying to define a FORCED VORTEX'
521  write(*,*)' Core swirl ratio Kr is strictly
522  & positive'
523  write(*,*) ' element number: ',nelem
524  call exit(201)
525  endif
526  if(prop(index+4).gt.1.d0) then
527  write(*,*)'*ERROR in checkinputvaluesnet:'
528  write(*,*)' trying to define a FORCED VORTEX'
529  write(*,*)' Core swirl ratio Kr is
530  & greater than 1'
531  write(*,*) ' element number: ',nelem
532  call exit(201)
533  endif
534 !
535 ! Rotation speed must be defined and positive
536 !
537  if(prop(index+5).le.0.d0) then
538  write(*,*)'*ERROR in checkinputvaluesnet:'
539  write(*,*)' trying to define a FORCED VORTEX'
540  write(*,*)
541  & ' Rotation speed n is strictly positive'
542  write(*,*) ' element number: ',nelem
543  call exit(201)
544  endif
545  endif
546 !
547 ! Absolute/relative system
548 !
549  elseif((lakon(nelem)(2:4).eq.'ATR').or.
550  & (lakon(nelem)(2:4).eq.'RTA')) then
551  if(prop(index+1).le.0.d0) then
552  write(*,*)'*ERROR in checkinputvaluesnet:'
553  write(*,*)' trying to define an element'
554  write(*,*)' TYPE= ABSOLUTE TO RELATIVE or'
555  write(*,*)' TYPE= RELATIVE TO ABSOLUTE'
556  write(*,*)' Rotation speed is strictly positive'
557  write(*,*)' element number: ',nelem
558  call exit(201)
559  elseif(prop(index+3).ne.0.d0) then
560  if(prop(index+2).ne.0.d0) then
561  write(*,*)'*ERROR in checkinputvaluesnet:'
562  write(*,*)' trying to define an element'
563  write(*,*)' TYPE= ABSOLUTE TO RELATIVE or'
564  write(*,*)' TYPE= RELATIVE TO ABSOLUTE'
565  write(*,*)' reference element has been provided'
566  write(*,*)' but tengential velocity
567  & is already defined'
568  write(*,*)' element number: ',nelem
569  call exit(201)
570  endif
571  endif
572  elseif((lakon(nelem)(2:5).ne.'LIPU').and.
573  & (lakon(nelem).ne.' ').and.
574  & (lakon(nelem)(2:4).ne.'LAB').and.
575  & (lakon(nelem)(2:5).ne.'CHAR').and.
576  & (lakon(nelem)(2:6).ne.'CARBS'))then
577  if((prop(index+1).lt.0.d0)) then
578  write(*,*) '*ERROR in checkinputvaluesnet: section area'
579  write(*,*) ' is not positive'
580  write(*,*) ' element number: ',nelem
581  call exit(201)
582  endif
583  endif
584 !
585  enddo
586 !
587  return
Hosted by OpenAircraft.com, (Michigan UAV, LLC)