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

Go to the source code of this file.

Functions/Subroutines

subroutine user_network_element_p0 (node1, node2, nodem, nelem, lakon, kon, ipkon, nactdog, identity, ielprop, prop, iflag, v, xflow, f, nodef, idirf, df, cp, R, physcon, dvi, numf, set, co, vold, mi, ttime, time, iaxial)
 

Function/Subroutine Documentation

◆ user_network_element_p0()

subroutine user_network_element_p0 ( integer, intent(in)  node1,
integer, intent(in)  node2,
integer, intent(in)  nodem,
integer, intent(in)  nelem,
character*8, dimension(*), intent(in)  lakon,
integer, dimension(*), intent(in)  kon,
integer, dimension(*), intent(in)  ipkon,
integer, dimension(0:3,*), intent(in)  nactdog,
logical, intent(inout)  identity,
integer, dimension(*), intent(in)  ielprop,
real*8, dimension(*), intent(in)  prop,
integer, intent(in)  iflag,
real*8, dimension(0:mi(2),*), intent(in)  v,
real*8, intent(inout)  xflow,
real*8, intent(inout)  f,
integer, dimension(*), intent(inout)  nodef,
integer, dimension(*), intent(inout)  idirf,
real*8, dimension(*), intent(inout)  df,
real*8, intent(in)  cp,
real*8, intent(in)  R,
real*8, dimension(*), intent(in)  physcon,
real*8, intent(in)  dvi,
integer, intent(inout)  numf,
character*81, dimension(*), intent(in)  set,
real*8, dimension(3,*), intent(in)  co,
real*8, dimension(0:mi(2),*), intent(in)  vold,
integer, dimension(*), intent(in)  mi,
real*8, intent(in)  ttime,
real*8, intent(in)  time,
integer, intent(in)  iaxial 
)
23 !
24 ! user subroutine user_network_element
25 !
26 ! skeleton file
27 !
28 ! INPUT:
29 !
30 ! node1 first node in element topology
31 ! node2 third node in element topology
32 ! nodem second node in element topology (middle node)
33 ! nelem element number
34 ! lakon(i) label of element i
35 ! kon connectivity list of all elements; the topology
36 ! of element i starts at kon(ipkon(i))
37 ! ipkon(i) pointer of element i into list kon
38 ! nactdog(j,i) global degree of freedom in the network
39 ! equation system of local dof j (0,1 or 2 for
40 ! networks) of node i. If nactdog(j,i)=0 the
41 ! variable is known
42 ! ielprop(i) pointer for element i into field prop. The
43 ! properties for element i start at
44 ! prop(ielprop(i)+1,...)
45 ! prop field of all properties
46 ! iflag indicates what information should be returned
47 ! by the routine:
48 ! 0: identity
49 ! 1: xflow
50 ! 2: numf, nodef, idirf, f, df
51 ! 3: none
52 ! v(0..mi(2),i) values at node i in the current network
53 ! iteration (0=total temperature,
54 ! 1=mass flow, 2=total pressure for network nodes,
55 ! 0=temperature, 1..3=displacements for structural
56 ! nodes)
57 ! cp specific heat at constant pressure corresponding
58 ! to a mean static temperature across the element
59 ! R specific gas constant
60 ! physcon(1..) physical constants (e.g. physcon(1) is absolute
61 ! zero in the unit systemof the user; cf. the
62 ! user's manual for the other entries)
63 ! dvi dynamical viscosity corresponding
64 ! to a mean static temperature across the element
65 ! set(i) set name corresponding to set i
66 ! co(1..3,i) coordinates of node i in the global system
67 ! vold(0..mi(2),i) values at node i at the start of the current network
68 ! iterations (0=total temperature,
69 ! 1=mass flow, 2=total pressure for network nodes,
70 ! 0=temperature, 1..3=displacements for structural
71 ! nodes)
72 ! mi(2) max degree of freedom per node (max over all
73 ! nodes) in fields like v(0:mi(2))...; the other
74 ! values of mi are not relevant here
75 ! ttime total time at the end of the current
76 ! thermo-mechanical increment. To reach the end
77 ! of this increment several thermo-mechanical
78 ! iterations are performed. For each of these
79 ! iterations a loop of network iterations is
80 ! performed
81 ! time step time a the end of the current thermo-
82 ! mechanical increment
83 ! iaxial number of times the current structure fits into
84 ! 360 degrees
85 !
86 !
87 ! OUTPUT:
88 !
89 ! identity if .true. the user_network_element routine is
90 ! not needed (all variables known)
91 ! xflow mass flow
92 ! f value of the element equation
93 ! nodef nodes corresponding to the variables in the
94 ! element equation
95 ! idirf degrees of freedom corresponding to the variables
96 ! in the element equation
97 ! df derivatives of the element equation w.r.t. its
98 ! variables
99 ! numf number of variables in the element equation
100 !
101 ! NOTE: to convert total temperature into static temperatures
102 ! subroutine
103 ! call ts_calc(xflow,Tt1,pt1,kappa,r,A,T1,icase)
104 ! may be used (cf. user_netowrk_element_p1.f for an example
105 ! of its use).
106 !
107  implicit none
108 !
109  logical identity
110  character*8 lakon(*)
111  character*81 set(*)
112 !
113  integer nelem,nactdog(0:3,*),node1,node2,nodem,numf,
114  & ielprop(*),nodef(*),idirf(*),iflag,ipkon(*),kon(*),
115  & iaxial,mi(*),inv,index,icase
116 !
117  real*8 prop(*),v(0:mi(2),*),xflow,f,df(*),r,cp,physcon(*),dvi,
118  & co(3,*),vold(0:mi(2),*),ttime,time,xmach,xflow_oil,tt1,tt2,
119  & reynolds,pt1,pt2,kappa,c,a,d,pi,t1,t2
120 !
121  intent(in) node1,node2,nodem,nelem,lakon,
122  & kon,ipkon,nactdog,ielprop,prop,iflag,v,
123  & cp,r,physcon,dvi,set,co,vold,mi,ttime,
124  & time,iaxial
125 !
126  intent(inout) xflow,f,df,identity,nodef,idirf,numf
127 !
128  if(iflag.eq.0) then
129 !
130 ! called by envtemp.f:
131 !
132 ! check whether element equation is needed (this is the
133 ! case if identity=.false.
134 !
135 ! identity=?
136 !
137  elseif(iflag.eq.1)then
138 !
139 ! called by initialnet.f:
140 !
141 ! calculation of the mass flow if everything else is known
142 !
143 ! xflow=?
144 !
145  elseif(iflag.eq.2)then
146 !
147 ! called by resultnet.f and mafillnet.f
148 !
149 ! numf=?
150 ! nodef(1...numf)=?
151 ! idirf(1...numf)=?
152 ! f=?
153 ! df(1...numf)=?
154 !
155  elseif(iflag.eq.3) then
156 !
157 ! called by flowoutput.f
158 !
159 ! storage in the .net-file
160 ! this is an example
161 !
162  write(1,*) ''
163  write(1,55) ' from node ',node1,
164  & ' to node ', node2,' : air massflow rate = '
165  &,inv*xflow,' ',
166  & ', oil massflow rate = ',xflow_oil,' '
167 !
168  if(inv.eq.1) then
169  write(1,56)' Inlet node ',node1,' : Tt1 = ',tt1,
170  & ' , Ts1 = ',t1,' , Pt1 = ',pt1, ' '
171 !
172  write(1,*)' Element ',nelem,lakon(nelem)
173  write(1,57)' dyn.visc = ',dvi,' , Re = '
174  & ,reynolds,', M = ',xmach
175 !
176  write(1,56)' Outlet node ',node2,' : Tt2 = ',tt2,
177  & ' , Ts2 = ',t2,' , Pt2 = ',pt2,' '
178 !
179  else if(inv.eq.-1) then
180  write(1,56)' Inlet node ',node2,': Tt1 = ',tt1,
181  & ' , Ts1 = ',t1,' , Pt1 = ',pt1, ' '
182  &
183  write(1,*)' element R ',set(numf)(1:30)
184  write(1,57)' dyn.visc. = ',dvi,' , Re ='
185  & ,reynolds,', M = ',xmach
186 !
187  write(1,56)' Outlet node ',node1,': Tt2 = ',tt2,
188  & ' , Ts2 = ',t2,' , Pt2 = ',pt2, ' '
189  endif
190 !
191  endif
192 !
193  55 format(1x,a,i6,a,i6,a,e11.4,a,a,e11.4,a)
194  56 format(1x,a,i6,a,e11.4,a,e11.4,a,e11.4,a)
195  57 format(1x,a,e11.4,a,e11.4,a,e11.4)
196 !
197 ! following lines are needed because the element equations
198 ! usually specify the mass flow for the complete cross section
199 !
200 ! in CalculiX, axisymmetric elements are expanded into 3D
201 ! using a sector of 360°/iaxial. Therefore, the mass flow and
202 ! the derivative of f w.r.t. the mass flow have to be adjusted
203 ! appropriately.
204 !
205  xflow=xflow/iaxial
206 !
207 ! only if the mass flow is an active degree of freedom:
208 !
209 ! df(mass_flow_dof)=df(mass flow dof)*iaxial
210 !
211  return
subroutine df(x, u, uprime, rpar, nev)
Definition: subspace.f:133
static double * e11
Definition: radflowload.c:42
Hosted by OpenAircraft.com, (Michigan UAV, LLC)