File: characteristic.f

package info (click to toggle)
calculix-ccx 2.11-1
  • links: PTS, VCS
  • area: main
  • in suites: buster, stretch
  • size: 10,188 kB
  • sloc: fortran: 115,312; ansic: 34,480; sh: 374; makefile: 35; perl: 15
file content (239 lines) | stat: -rw-r--r-- 7,066 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
!     
!     CalculiX - A 3-dimensional finite element program
!     Copyright (C) 1998-2015 Guido Dhondt
!     
!     This program is free software; you can redistribute it and/or
!     modify it under the terms of the GNU General Public License as
!     published by the Free Software Foundation(version 2);
!     
!     
!     This program is distributed in the hope that it will be useful,
!     but WITHOUT ANY WARRANTY; without even the implied warranty of 
!     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 
!     GNU General Public License for more details.
!     
!     You should have received a copy of the GNU General Public License
!     along with this program; if not, write to the Free Software
!     Foundation, Inc., 675 Mass Ave, Cambridge, MA 02139, USA.
!     
      subroutine  characteristic(node1,node2,nodem,nelem,lakon,
     &     kon,ipkon,
     &     nactdog,identity,ielprop,prop,iflag,v,xflow,f,
     &     nodef,idirf,df,cp,r,physcon,dvi,numf,set,
     &     mi,ttime,time,iaxial)
!     
!     This subroutine is used to enables the processing of empiric 
!     given under the form
!     massflow*dsqrt(T1)/Pt1=f((Pt1-Pt2)/Pt1) and T1=T2
!     characteristics the subroutine proceeds using 
!     linear interpolation to estimate the values for the whole characteristic
!     note that the characteristic is implicitly containing the point (0,0)
!
!     author: Yannick Muller
!     
      implicit none
!     
      logical identity
!
      character*8 lakon(*)
      character*81 set(*)
!     
      integer nelem,nactdog(0:3,*),node1,node2,nodem,kon(*),ipkon(*),
     &     ielprop(*),nodef(4),idirf(4),index,iflag,
     &     inv,id,numf,npu,i,mi(*),iaxial
!     
      real*8 prop(*),v(0:mi(2),*),xflow,f,df(4),cp,r,dvi,
     &     p1,p2,physcon(*),ttime,time,xmach,kappa,
     &     xpu(100),ypu(100),Qred,p1mp2zp1,T1,scal,T2
!     
      if (iflag.eq.0) then
         identity=.true.
!     
         if(nactdog(2,node1).ne.0)then
            identity=.false.
         elseif(nactdog(2,node2).ne.0)then
            identity=.false.
         elseif(nactdog(1,nodem).ne.0)then
            identity=.false.
         endif
!     
      elseif ((iflag.eq.1).or.(iflag.eq.2)) then
!     
         index=ielprop(nelem)
!     
         npu=nint(prop(index+2))
         scal=prop(index+1)

         do i=1, 100
            xpu(i)=0
            ypu(i)=0
         enddo
!
         do i=1,npu
            xpu(i)=prop(index+2*i+1)
            ypu(i)=prop(index+2*i+2)
         enddo
!     
         p1=v(2,node1)
         p2=v(2,node2)       
!     
         if(p1.ge.p2) then
            inv=1
            T1=v(0,node1)+physcon(1)
         else
            inv=-1
            p1=v(2,node2)
            p2=v(2,node1)
            T1=v(0,node2)+physcon(1)
         endif
!     
         p1mp2zp1=(P1-P2)/P1
!     
         if(iflag.eq.1) then
            
            call ident(xpu,p1mp2zp1,npu,id)
            if(id.eq.0) then
               Qred=scal*ypu(1)
               xflow=inv*Qred*P1/dsqrt(T1)
            elseif(id.ge.npu) then
               Qred=scal*ypu(npu)
               xflow=inv*Qred*P1/dsqrt(T1)
            else
               Qred=scal*(ypu(id)+(ypu(id+1)-ypu(id))
     &             *(p1mp2zp1-xpu(id))/(xpu(id+1)-xpu(id)))
               xflow=inv*Qred*P1/dsqrt(T1)
            endif
!
         elseif (iflag.eq.2) then
            numf=4
!     
            p1=v(2,node1)
            p2=v(2,node2) 
            xflow=v(1,nodem)*iaxial
!     
            if (p1.ge.p2) then
!     
               inv=1
               xflow=v(1,nodem)*iaxial
               T1=v(0,node1)+physcon(1)
               nodef(1)=node1
               nodef(2)=node1
               nodef(3)=nodem
               nodef(4)=node2
!     
            else 
!     
               inv=-1
               p1=v(2,node2)
               p2=v(2,node1) 
               T1=v(0,node2)+physcon(1)
               xflow=-v(1,nodem)*iaxial
               nodef(1)=node2
               nodef(2)=node2
               nodef(3)=nodem
               nodef(4)=node1
            endif
!     
            idirf(1)=2
            idirf(2)=0
            idirf(3)=1
            idirf(4)=2
!     
            df(2)=xflow/(2.d0*P1*dsqrt(T1))
            df(3)=inv*dsqrt(T1)/P1
!     
            call ident(xpu,p1mp2zp1,npu,id)
!     
            if(id.eq.0) then
               f=dabs(xflow)/P1*dsqrt(T1)-scal*ypu(1)
               df(4)=0.01d0
               df(1)=-xflow*dsqrt(T1)/P1**2
!     
            elseif(id.ge.npu) then
               f=dabs(xflow)/P1*dsqrt(T1)-scal*ypu(npu)
               df(4)=0.01d0
               df(1)=-xflow*dsqrt(T1)/P1**2
!    
            else
               f=dabs(xflow)/P1*dsqrt(T1)-(scal*ypu(id)
     &             +scal*(ypu(id+1)-ypu(id))
     &              *(p1mp2zp1-xpu(id))/(xpu(id+1)-xpu(id)))
!
               df(4)=scal*(ypu(id+1)-ypu(id))/(xpu(id+1)-xpu(id))*1/p1
!
               df(1)=-xflow*dsqrt(T1)/P1**2-(P2/P1**2)
     &              *(scal*(ypu(id+1)-ypu(id))/(xpu(id+1)-xpu(id)))
            endif
         endif

      elseif(iflag.eq.3)  then
         p1=v(2,node1)
         p2=v(2,node2) 
         xflow=v(1,nodem)*iaxial
         kappa=(cp/(cp-r))
         xmach=dsqrt(((p1/p2)**((kappa-1.d0)/kappa)-1.d0)*2.d0/
     &          (kappa-1.d0))
!     
         if (p1.ge.p2) then
!     
            inv=1
            xflow=v(1,nodem)*iaxial
            T1=v(0,node1)+physcon(1)
            T2=v(0,node2)+physcon(1)
            nodef(1)=node1
            nodef(2)=node1
            nodef(3)=nodem
            nodef(4)=node2
!     
         else 
!     
            inv=-1
            p1=v(2,node2)
            p2=v(2,node1) 
            T1=v(0,node2)+physcon(1)
            T2=v(0,node1)+physcon(1)
            xflow=-v(1,nodem)*iaxial
            nodef(1)=node2
            nodef(2)=node2
            nodef(3)=nodem
            nodef(4)=node1
         endif
!
         write(1,*) ''
         write(1,55) ' from node',node1,
     &        ' to node', node2,':   air massflow rate=',xflow
!
 55      FORMAT(1X,A,I6,A,I6,A,e11.4,A,A,e11.4,A)
!
         if(inv.eq.1) then
            write(1,56)'       Inlet node  ',node1,':   Tt1=',T1,
     &           ', Ts1=',T1,', Pt1=',P1
         
            write(1,*)'             Element ',nelem,lakon(nelem)
            write(1,57) 'M = ',xmach
!
            write(1,56)'       Outlet node ',node2,':   Tt2=',T2,
     &           ', Ts2=',T2,', Pt2=',P2
!     
         else if(inv.eq.-1) then
            write(1,56)'       Inlet node  ',node2,':    Tt1=',T1,
     &           ', Ts1=',T1,', Pt1=',P1
     &          
            write(1,*)'             Element ',nelem,lakon(nelem)
            write(1,57) 'M = ',xmach
!
            write(1,56)'       Outlet node ',node1,':    Tt2=',T2,
     &           ', Ts2=',T2,', Pt2=',P2
!               
         endif
!      
 56      FORMAT(1X,A,I6,A,e11.4,A,e11.4,A,e11.4,A)
 57      format(40x,a,e11.4)
!     
      endif
!     
      xflow=xflow/iaxial
      df(3)=df(3)*iaxial
!     
      return
      end