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

Go to the source code of this file.

Functions/Subroutines

subroutine envtemp (itg, ieg, ntg, ntr, sideload, nelemload, ipkon, kon, lakon, ielmat, ne, nload, kontri, ntri, nloadtr, nflow, ndirboun, nactdog, nodeboun, nacteq, nboun, ielprop, prop, nteq, v, network, physcon, shcon, ntmat_, co, vold, set, nshcon, rhcon, nrhcon, mi, nmpc, nodempc, ipompc, labmpc, ikboun, nasym, ttime, time, iaxial)
 

Function/Subroutine Documentation

◆ envtemp()

subroutine envtemp ( integer, dimension(*)  itg,
integer, dimension(*)  ieg,
integer  ntg,
integer  ntr,
character*20, dimension(*)  sideload,
integer, dimension(2,*)  nelemload,
integer, dimension(*)  ipkon,
integer, dimension(*)  kon,
character*8, dimension(*)  lakon,
integer, dimension(mi(3),*)  ielmat,
integer  ne,
integer  nload,
integer, dimension(4,*)  kontri,
integer  ntri,
integer, dimension(*)  nloadtr,
integer  nflow,
integer, dimension(*)  ndirboun,
integer, dimension(0:3,*)  nactdog,
integer, dimension(*)  nodeboun,
integer, dimension(0:3,*)  nacteq,
integer  nboun,
integer, dimension(*)  ielprop,
real*8, dimension(*)  prop,
integer  nteq,
real*8, dimension(0:mi(2),*)  v,
integer  network,
real*8, dimension(*)  physcon,
real*8, dimension(0:3,ntmat_,*)  shcon,
integer  ntmat_,
real*8, dimension(3,*)  co,
real*8, dimension(0:mi(2),*)  vold,
character*81, dimension(*)  set,
integer, dimension(*)  nshcon,
real*8, dimension(*)  rhcon,
integer, dimension(*)  nrhcon,
integer, dimension(*)  mi,
integer  nmpc,
integer, dimension(3,*)  nodempc,
integer, dimension(*)  ipompc,
character*20, dimension(*)  labmpc,
integer, dimension(*)  ikboun,
integer  nasym,
real*8  ttime,
real*8  time,
integer  iaxial 
)
25 !
26 ! determines the number of gas temperatures and radiation
27 ! temperatures
28 !
29  implicit none
30 !
31  logical identity,walltemp,temperaturebc,pressurebc,massflowbcall,
32  & pressurebcall
33 !
34  character*8 lakon(*)
35  character*20 labmpc(*)
36  character*20 sideload(*)
37  character*81 set(*)
38 !
39  integer itg(*),ntg,ntr,nelemload(2,*),ipkon(*),network,mi(*),
40  & kon(*),ielmat(mi(3),*),ne,i,j,k,l,index,id,node,nload,
41  & ifaceq(8,6),ider,nasym,indexe,iaxial,networkmpcs,
42  & ifacet(6,4),ifacew(8,5),kontri3(3,1),kontri4(3,2),
43  & kontri6(3,4),kontri8(3,6),kontri(4,*),ntri,
44  & konf(8),nloadtr(*),nelem,nope,nopes,ig,nflow,ieg(*),
45  & ndirboun(*),nactdog(0:3,*),nboun,nodeboun(*),ntmat_,
46  & idir,ntq,nteq,nacteq(0:3,*),node1,node2,nodem,
47  & ielprop(*),idirf(8),iflag,imat,numf,nrhcon(*),nshcon(*),
48  & nmpc,nodempc(3,*),ipompc(*),ikboun(*)
49 !
50  real*8 prop(*),f,xflow,nodef(8),df(8),v(0:mi(2),*),g(3),
51  & cp,r,physcon(*),shcon(0:3,ntmat_,*),rho,ttime,time,
52  & co(3,*),dvi,vold(0:mi(2),*),rhcon(*)
53 !
54  data ifaceq /4,3,2,1,11,10,9,12,
55  & 5,6,7,8,13,14,15,16,
56  & 1,2,6,5,9,18,13,17,
57  & 2,3,7,6,10,19,14,18,
58  & 3,4,8,7,11,20,15,19,
59  & 4,1,5,8,12,17,16,20/
60  data ifacet /1,3,2,7,6,5,
61  & 1,2,4,5,9,8,
62  & 2,3,4,6,10,9,
63  & 1,4,3,8,10,7/
64  data ifacew /1,3,2,9,8,7,0,0,
65  & 4,5,6,10,11,12,0,0,
66  & 1,2,5,4,7,14,10,13,
67  & 2,3,6,5,8,15,11,14,
68  & 4,6,3,1,12,15,9,13/
69  data kontri3 /1,2,3/
70  data kontri4 /1,2,4,2,3,4/
71  data kontri6 /1,4,6,4,5,6,4,2,5,6,5,3/
72  data kontri8 /1,5,8,8,5,7,8,7,4,5,2,6,5,6,7,7,6,3/
73 !
74  ntg=0
75  ntr=0
76  ntri=0
77 !
78  walltemp=.false.
79  temperaturebc=.false.
80  pressurebc=.false.
81  massflowbcall=.true.
82  pressurebcall=.true.
83 !
84  networkmpcs=0
85 !
86 ! ordering the gas temperature nodes and counting them
87 ! counting the radiation temperatures
88 !
89  do i=1,nload
90  if(sideload(i)(3:4).eq.'FC') then
91  walltemp=.true.
92  call nident(itg,nelemload(2,i),ntg,id)
93  if(id.gt.0) then
94  if(itg(id).eq.nelemload(2,i)) then
95  nactdog(0,nelemload(2,i))=1
96  cycle
97  endif
98  endif
99  ntg=ntg+1
100  do j=ntg,id+2,-1
101  itg(j)=itg(j-1)
102  enddo
103  itg(id+1)=nelemload(2,i)
104  nactdog(0,nelemload(2,i))=1
105 !
106  elseif(sideload(i)(3:4).eq.'NP') then
107  call nident(itg,nelemload(2,i),ntg,id)
108  if(id.gt.0) then
109  if(itg(id).eq.nelemload(2,i)) then
110  nactdog(2,nelemload(2,i))=1
111  cycle
112  endif
113  endif
114  ntg=ntg+1
115  do j=ntg,id+2,-1
116  itg(j)=itg(j-1)
117  enddo
118  itg(id+1)=nelemload(2,i)
119  nactdog(2,nelemload(2,i))=1
120 !
121  elseif(sideload(i)(3:4).eq.'CR') then
122  ntr=ntr+1
123  nelem=nelemload(1,i)
124  read(sideload(i)(2:2),'(i1)') ig
125 !
126 ! number of nodes in the face
127 !
128  if(lakon(nelem)(4:4).eq.'2') then
129  nope=20
130  nopes=8
131  elseif(lakon(nelem)(4:4).eq.'8') then
132  nope=8
133  nopes=4
134  elseif(lakon(nelem)(4:5).eq.'10') then
135  nope=10
136  nopes=6
137  elseif(lakon(nelem)(4:4).eq.'4') then
138  nope=4
139  nopes=3
140  elseif(lakon(nelem)(4:4).eq.'6') then
141  nope=6
142  if(ig.le.2) then
143  nopes=3
144  else
145  nopes=4
146  endif
147  elseif(lakon(nelem)(4:5).eq.'15') then
148  nope=15
149  if(ig.le.2) then
150  nopes=6
151  else
152  nopes=8
153  endif
154  endif
155 !
156 ! nodes in the face
157 !
158  if((nope.eq.20).or.(nope.eq.8)) then
159  do k=1,nopes
160  konf(k)=kon(ipkon(nelem)+ifaceq(k,ig))
161  enddo
162  elseif((nope.eq.10).or.(nope.eq.4)) then
163  do k=1,nopes
164  konf(k)=kon(ipkon(nelem)+ifacet(k,ig))
165  enddo
166  else
167  do k=1,nopes
168  konf(k)=kon(ipkon(nelem)+ifacew(k,ig))
169  enddo
170  endif
171 !
172 ! triangulation of the face
173 !
174  nloadtr(ntr)=i
175  if((lakon(nelem)(4:4).eq.'2').or.
176  & ((lakon(nelem)(4:5).eq.'15').and.(ig.gt.2))) then
177  do k=1,6
178  ntri=ntri+1
179  do l=1,3
180  kontri(l,ntri)=konf(kontri8(l,k))
181  enddo
182  kontri(4,ntri)=ntr
183  enddo
184  elseif((lakon(nelem)(4:4).eq.'8').or.
185  & ((lakon(nelem)(4:4).eq.'6').and.(ig.gt.2))) then
186  do k=1,2
187  ntri=ntri+1
188  do l=1,3
189  kontri(l,ntri)=konf(kontri4(l,k))
190  enddo
191  kontri(4,ntri)=ntr
192  enddo
193  elseif((lakon(nelem)(4:5).eq.'10').or.
194  & ((lakon(nelem)(4:5).eq.'15').and.(ig.le.2))) then
195  do k=1,4
196  ntri=ntri+1
197  do l=1,3
198  kontri(l,ntri)=konf(kontri6(l,k))
199  enddo
200  kontri(4,ntri)=ntr
201  enddo
202  elseif((lakon(nelem)(4:4).eq.'4').or.
203  & ((lakon(nelem)(4:4).eq.'6').and.(ig.le.2))) then
204  do k=1,1
205  ntri=ntri+1
206  do l=1,3
207  kontri(l,ntri)=konf(kontri3(l,k))
208  enddo
209  kontri(4,ntri)=ntr
210  enddo
211  endif
212  endif
213  enddo
214 !
215 ! storing the gas elements in a dedicated array
216 !
217  nflow=0
218 !
219  do i=1,ne
220  if(lakon(i)(1:1).eq.'D') then
221  if((lakon(i)(2:2).ne.' ').and.(network.ne.1)) then
222  nflow=nflow+1
223  ieg(nflow)=i
224  else
225  nasym=1
226 !
227 ! removing gas nodes belonging to 'D '-elements
228 ! in which a 'FC'-film condition was applied from the
229 ! itg vector
230 !
231  indexe=ipkon(i)
232  do j=1,3,2
233  node=kon(indexe+j)
234  call nident(itg,node,ntg,id)
235  if(id.gt.0) then
236  if(itg(id).eq.node) then
237  ntg=ntg-1
238  do k=id,ntg
239  itg(k)=itg(k+1)
240  enddo
241  nactdog(0,node)=0
242  nactdog(2,node)=0
243  endif
244  endif
245  enddo
246 !
247  endif
248  endif
249 !
250 ! removing gas nodes belonging to advective elements
251 ! (last node in the element topology)
252 ! these are usually gas nodes not belonging to any
253 ! "D"-element
254 !
255  if(lakon(i)(1:7).eq.'ESPRNGF') then
256  read(lakon(i)(8:8),'(i1)') nopes
257  nope=nopes+1
258  node=kon(ipkon(i)+nope)
259 !
260  call nident(itg,node,ntg,id)
261  if(id.gt.0) then
262  if(itg(id).eq.node) then
263  ntg=ntg-1
264  do k=id,ntg
265  itg(k)=itg(k+1)
266  enddo
267  nactdog(0,node)=0
268  nactdog(2,node)=0
269  endif
270  endif
271  endif
272  enddo
273 !
274 ! mass flux nodes are also taken as unknowns in the
275 ! gas temperature system; determining the active
276 ! degrees of freedom
277 !
278 ! first node of the flow element
279 !
280  do i=1,nflow
281  index=ipkon(ieg(i))
282  node=kon(index+1)
283  if (node.eq.0) cycle
284  call nident(itg,node,ntg,id)
285  if(id.gt.0) then
286  if(itg(id).eq.node) then
287 !
288 ! upstream depth of SO,WO and DO is known
289 !
290  if((lakon(ieg(i))(4:7).eq.'CHSO').or.
291  & (lakon(ieg(i))(4:7).eq.'CHWO').or.
292  & (lakon(ieg(i))(4:7).eq.'CHDO')) cycle
293  nactdog(0,node)=1
294  nactdog(2,node)=1
295  cycle
296  endif
297  endif
298  ntg=ntg+1
299  do j=ntg,id+2,-1
300  itg(j)=itg(j-1)
301  enddo
302  itg(id+1)=node
303 !
304 ! upstream depth of SO,WO and DO is known
305 !
306  if((lakon(ieg(i))(4:7).eq.'CHSO').or.
307  & (lakon(ieg(i))(4:7).eq.'CHWO').or.
308  & (lakon(ieg(i))(4:7).eq.'CHDO')) cycle
309  nactdog(0,node)=1
310  nactdog(2,node)=1
311  enddo
312 !
313 ! middle node of the flow element :flux
314 !
315  do i=1,nflow
316  index=ipkon(ieg(i))
317  node=kon(index+2)
318  call nident(itg,node,ntg,id)
319  if(id.gt.0) then
320  if(itg(id).eq.node) cycle
321  endif
322  ntg=ntg+1
323  do j=ntg,id+2,-1
324  itg(j)=itg(j-1)
325  enddo
326  itg(id+1)=node
327  nactdog(1,node)=1
328 !
329 ! variable geometric property
330 !
331  if(lakon(ieg(i))(6:7).eq.'GV') then
332  index=ielprop(ieg(i))
333  if(prop(index+2).le.0.d0) nactdog(3,node)=1
334  elseif((lakon(ieg(i))(4:7).eq.'CHSG').or.
335  & (lakon(ieg(i))(4:7).eq.'CHWE').or.
336  & (lakon(ieg(i))(4:7).eq.'CHDS')) then
337  nactdog(3,node)=1
338  elseif(lakon(ieg(i))(2:7).eq.'ACCTUB') then
339 !
340  index=ielprop(ieg(i))
341  if(prop(index+1).eq.2) then
342 ! Interval factor unknown,setting a DOF for geometry
343  nactdog(3,node)=1
344  elseif(prop(index+1).eq.3) then
345 ! Hole diameter factor unknown,setting a DOF for geometry
346  nactdog(3,node)=1
347  endif
348  endif
349  enddo
350 !
351 ! third node of the flow element
352 !
353  do i=1,nflow
354  index=ipkon(ieg(i))
355  node=kon(index+3)
356  if (node.eq.0) cycle
357  call nident(itg,node,ntg,id)
358  if(id.gt.0) then
359  if(itg(id).eq.node) then
360 !
361 ! downstream depth of SG,WE and DS is known
362 !
363  if((lakon(ieg(i))(4:7).eq.'CHSG').or.
364  & (lakon(ieg(i))(4:7).eq.'CHWE').or.
365  & (lakon(ieg(i))(4:7).eq.'CHDS')) cycle
366  nactdog(0,node)=1
367  nactdog(2,node)=1
368  cycle
369  endif
370  endif
371  ntg=ntg+1
372  do j=ntg,id+2,-1
373  itg(j)=itg(j-1)
374  enddo
375  itg(id+1)=node
376 !
377 ! downstream depth of SG,WE and DS is known
378 !
379  if((lakon(ieg(i))(4:7).eq.'CHSG').or.
380  & (lakon(ieg(i))(4:7).eq.'CHWE').or.
381  & (lakon(ieg(i))(4:7).eq.'CHDS')) cycle
382  nactdog(0,node)=1
383  nactdog(2,node)=1
384  enddo
385 !
386 ! tagging the network MPC's
387 !
388  do i=1,nmpc
389 !
390 ! check whether network MPC
391 !
392  index=ipompc(i)
393  do
394  node=nodempc(1,index)
395  call nident(itg,node,ntg,id)
396  if(id.gt.0) then
397  if(itg(id).eq.node) then
398  labmpc(i)(1:7)='NETWORK'
399  networkmpcs=1
400  exit
401  endif
402  endif
403  index=nodempc(3,index)
404  if(index.eq.0) exit
405  enddo
406  enddo
407 !
408 ! taking the MPC network nodes into account
409 !
410  if(networkmpcs.eq.1) then
411  do i=1,nmpc
412  if(labmpc(i)(1:7).ne.'NETWORK') cycle
413 !
414  index=ipompc(i)
415  do
416  node=nodempc(1,index)
417  call nident(itg,node,ntg,id)
418  if(id.gt.0) then
419  if(itg(id).eq.node) then
420  nactdog(nodempc(2,index),node)=1
421  index=nodempc(3,index)
422  if(index.eq.0) then
423  exit
424  else
425  cycle
426  endif
427  endif
428  endif
429 !
430 ! adding a node to itg
431 !
432  ntg=ntg+1
433  do j=ntg,id+2,-1
434  itg(j)=itg(j-1)
435  enddo
436  itg(id+1)=node
437  nactdog(nodempc(2,index),node)=1
438  index=nodempc(3,index)
439  if(index.eq.0) exit
440  enddo
441  enddo
442  endif
443 !
444 ! subtracting the SPC conditions
445 !
446  do i=1,nboun
447  node=nodeboun(i)
448  call nident(itg,node,ntg,id)
449  if (id.gt.0) then
450  if (itg(id).eq.node) then
451  idir=ndirboun(i)
452  nactdog(idir,node)=0
453  if(idir.eq.0) then
454  temperaturebc=.true.
455  elseif(idir.eq.2) then
456  pressurebc=.true.
457  endif
458  endif
459  endif
460  enddo
461 !
462 ! determining the active equations
463 !
464 ! element contributions
465 !
466  iflag=0
467  do i=1,nflow
468  nelem=ieg(i)
469  index=ipkon(nelem)
470  node1=kon(index+1)
471  nodem=kon(index+2)
472  node2=kon(index+3)
473 !
474 ! "end of network" element ---X---O
475 ! 1 m 2
476 !
477  if(node1.eq.0)then
478  if ((nactdog(1,nodem).ne.0).or.(nactdog(3,nodem).ne.0))then
479  nacteq(1,node2)=1 ! mass equation
480  endif
481  if (nactdog(0,node2).ne.0)then
482  nacteq(0,node2)=1 ! energy equation
483  endif
484 !
485 ! "end of network" element node O---X---
486 ! 1 m 2
487  elseif (node2.eq.0) then
488  if ((nactdog(1,nodem).ne.0).or.(nactdog(3,nodem).ne.0))then
489  nacteq(1,node1)=1 ! mass equation
490  endif
491  if (nactdog(0,node1).ne.0)then
492  nacteq(0,node1)=1 ! energy equation
493  endif
494 !
495 ! "flow element" O---X---O
496 ! 1 m 2
497 !
498  else
499  if((nactdog(1,nodem).ne.0).or.(nactdog(3,nodem).ne.0)) then
500  nacteq(1,node2)=1 ! mass equation
501  nacteq(1,node1)=1 ! mass equation
502  endif
503 !
504  ider=0
505  call flux(node1,node2,nodem,nelem,lakon,kon,ipkon,
506  & nactdog,identity,ielprop,prop,iflag,v,xflow,f,
507  & nodef,idirf,df,cp,r,rho,physcon,g,co,dvi,numf,
508  & vold,set,shcon,nshcon,rhcon,nrhcon,ntmat_,mi,ider,
509  & ttime,time,iaxial)
510 !
511  if (.not.identity) then
512  nacteq(2,nodem)=1 ! momentum equation
513  endif
514 !
515  if (nactdog(0,node1).ne.0)then
516  nacteq(0,node1)=1 ! energy equation
517  endif
518 !
519  if (nactdog(0,node2).ne.0)then
520  nacteq(0,node2)=1 ! energy equation
521  endif
522  endif
523  enddo
524 !
525 ! wall convective contributions
526 !
527  do i=1, nload
528  if((sideload(i)(3:4)).eq.'FC')then
529  node=nelemload(2,i)
530  if (nactdog(0,node).ne.0) then
531  nacteq(0,node)=1
532  endif
533  endif
534  enddo
535 !
536 ! removing the energy equation from those end nodes for which
537 ! the temperature constitutes the first term in a network MPC
538 !
539  if(networkmpcs.eq.1) then
540  do i=1,nmpc
541  if(labmpc(i)(1:7).ne.'NETWORK') cycle
542  index=ipompc(i)
543  idir=nodempc(2,index)
544  if(idir.eq.0) then
545  node=nodempc(1,index)
546  call nident(itg,node,ntg,id)
547  if(id.gt.0) then
548  if(itg(id).eq.node) then
549  nacteq(0,node)=0
550  endif
551  endif
552  endif
553  enddo
554  endif
555 !
556 ! check whether all mass flow is known
557 !
558  do i=1,ntg
559  node=itg(i)
560  if((nactdog(1,node).ne.0).or.(nactdog(3,node).ne.0)) then
561  massflowbcall=.false.
562  exit
563  endif
564  enddo
565 !
566 ! check whether all pressures are known
567 !
568  do i=1,ntg
569  node=itg(i)
570  if(nactdog(2,node).ne.0) then
571  pressurebcall=.false.
572  exit
573  endif
574  enddo
575 !
576 ! check for special cases
577 !
578 c network=1
579 !
580  if(massflowbcall.and.((.not.pressurebc).or.(pressurebcall))) then
581 !
582 ! purely thermal (only set to 1 if D-type elements are present)
583 !
584  if(network.gt.1) network=2
585  do i=1,nflow
586  nelem=ieg(i)
587  index=ipkon(nelem)
588  node1=kon(index+1)
589  nodem=kon(index+2)
590  node2=kon(index+3)
591  nacteq(2,nodem)=0
592  if(node1.ne.0) nactdog(2,node1)=0
593  if(node2.ne.0) nactdog(2,node2)=0
594  enddo
595  elseif((.not.temperaturebc).and.(.not.walltemp)) then
596 !
597 ! pure liquid dynamics
598 !
599  write(*,*) '*INFO in envtemp: no thermal boundary conditions'
600  write(*,*) ' detected; the network is considered to be'
601  write(*,*) ' athermal and no gas temperatures will be'
602  write(*,*) ' calculated'
603  network=4
604  do i=1,ntg
605  node=itg(i)
606  nactdog(0,node)=0
607  nacteq(0,node)=0
608  enddo
609  elseif((.not.temperaturebc).and.walltemp) then
610  write(*,*) '*ERROR in envtemp: at least one temperature'
611  write(*,*) ' boundary condition must be given'
612  call exit(201)
613  elseif(.not.pressurebc) then
614  write(*,*) '*ERROR in envtemp: at least one pressure'
615  write(*,*) ' boundary condition must be given'
616  call exit(201)
617  endif
618 !
619 ! check whether a specific gas constant was defined for all fluid
620 ! elements (except for purely thermal calculations)
621 !
622  if(network.gt.2) then
623  do i=1,nflow
624  nelem=ieg(i)
625  if((lakon(nelem)(2:3).eq.'LI').or.
626  & (lakon(nelem)(2:3).eq.'LP').or.
627  & (lakon(nelem)(2:3).eq.' ')) cycle
628  imat=ielmat(1,nelem)
629  r=shcon(3,1,imat)
630  if(r.lt.1.d-10) then
631  write(*,*)'*ERROR in envtemp: specific gas',
632  & 'constant is close to zero'
633  call exit(201)
634  endif
635  enddo
636  endif
637 !
638 ! numbering the active equations
639 !
640  nteq=0
641  do i=1,ntg
642  node=itg(i)
643  do j=0,2
644  if (nacteq(j,node).ne.0) then
645  nteq=nteq+1
646  nacteq(j,node)=nteq
647  endif
648  enddo
649 ! write(30,*) 'unknowns ',node,(nactdog(j,node),j=0,3)
650  enddo
651 ! do i=1,ntg
652 ! node=itg(i)
653 ! write(30,*) 'equations',node,(nacteq(j,node),j=0,2)
654 ! enddo
655 !
656 ! taking network MPC's into account
657 !
658  if(networkmpcs.eq.1) then
659  do i=1,nmpc
660  if(labmpc(i)(1:7).eq.'NETWORK') nteq=nteq+1
661  enddo
662  endif
663 !
664 ! numbering the active degrees of freedom
665 !
666  ntq=0
667  do i=1,ntg
668  node=itg(i)
669  do j=0,3
670  if (nactdog(j,node).ne.0) then
671  ntq=ntq+1
672  nactdog(j,node)=ntq
673  endif
674  enddo
675  enddo
676 c
677 c open(30,file='dummy',status='unknown')
678 c write(30,*) 'nactdog'
679 c do i=1,ntg
680 c write(30,*) itg(i),(nactdog(j,itg(i)),j=0,3)
681 c enddo
682 c
683 c write(30,*) ''
684 c write(30,*) 'nacteq'
685 c do i=1,ntg
686 c write(30,*) itg(i),(nacteq(j,itg(i)),j=0,3)
687 c enddo
688 c close(30)
689 !
690  if(ntq.ne.nteq) then
691  write(*,*) '*ERROR in envtemp:'
692  write(*,*) '*****number of network equations is not equal to'
693  write(*,*) ' number of active degrees of freedom*****'
694  write(*,*) ' # of network equations = ',nteq
695  write(*,*) ' # of active degrees of freedom= ',ntq
696  call exit(201)
697  endif
698 !
699 ! for isothermal gas pipes the energy equation in the
700 ! topologically downstream node is replaced by an equation
701 ! expressing the equality of the static temperature at both
702 ! ends of the pipe. To this end these downstream nodes are
703 ! referring in nacteq(3,*) to the topologically upstream node
704 !
705 ! if the temperature in the downstream node is a boundary
706 ! condition (i.e. the energy equation is not built), the
707 ! energy equation in the upstream node is replaced. In that
708 ! case nacteq(3,upstreamnode) refers to the downstream node.
709 !
710 ! if both nodes are boundary conditions, nothing is done
711 !
712  do i=1,nflow
713  nelem=ieg(i)
714  if((lakon(nelem)(1:4).eq."DGAP")
715  & .and.(lakon(nelem)(6:6).eq."I")) then
716 !
717  index=ipkon(nelem)
718  node1=kon(index+1)
719  node2=kon(index+3)
720  if((node1.eq.0).or.(node2.eq.0)) cycle
721 !
722  if(nacteq(0,node2).ne.0) then
723  nacteq(3,node2)=node1
724  elseif(nacteq(0,node1).ne.0) then
725  nacteq(3,node1)=node2
726  endif
727  endif
728  enddo
729 !
730  return
subroutine df(x, u, uprime, rpar, nev)
Definition: subspace.f:133
subroutine networkmpcs(inpc, textpart, ipompc, nodempc, coefmpc, nmpc, nmpc_, mpcfree, nk, ikmpc, ilmpc, labmpc, istep, istat, n, iline, ipol, inl, ipoinp, inp, ipoinpc)
Definition: networkmpcs.f:22
subroutine flux(node1, node2, nodem, nelem, lakon, kon, ipkon, nactdog, identity, ielprop, prop, kflag, v, xflow, f, nodef, idirf, df, cp, R, rho, physcon, g, co, dvi, numf, vold, set, shcon, nshcon, rhcon, nrhcon, ntmat_, mi, ider, ttime, time, iaxial)
Definition: flux.f:24
subroutine nident(x, px, n, id)
Definition: nident.f:26
Hosted by OpenAircraft.com, (Michigan UAV, LLC)