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

Go to the source code of this file.

Functions/Subroutines

subroutine dload (f, kstep, kinc, time, noel, npt, layer, kspt, coords, jltyp, loadtype, vold, co, lakonl, konl, ipompc, nodempc, coefmpc, nmpc, ikmpc, ilmpc, iscale, veold, rho, amat, mi)
 

Function/Subroutine Documentation

◆ dload()

subroutine dload ( real*8, intent(inout)  f,
integer, intent(in)  kstep,
integer, intent(in)  kinc,
real*8, dimension(2), intent(in)  time,
integer, intent(in)  noel,
integer, intent(in)  npt,
integer, intent(in)  layer,
integer, intent(in)  kspt,
real*8, dimension(3), intent(in)  coords,
integer, intent(in)  jltyp,
character*20, intent(in)  loadtype,
real*8, dimension(0:mi(2),*), intent(in)  vold,
real*8, dimension(3,*), intent(in)  co,
character*8, intent(in)  lakonl,
integer, dimension(20), intent(in)  konl,
integer, dimension(*), intent(in)  ipompc,
integer, dimension(3,*), intent(in)  nodempc,
real*8, dimension(*), intent(in)  coefmpc,
integer, intent(in)  nmpc,
integer, dimension(*), intent(in)  ikmpc,
integer, dimension(*), intent(in)  ilmpc,
integer, intent(inout)  iscale,
real*8, dimension(0:mi(2),*), intent(in)  veold,
real*8, intent(in)  rho,
character*80, intent(in)  amat,
integer, dimension(*), intent(in)  mi 
)
23 !
24 ! user subroutine dload
25 !
26 !
27 ! INPUT:
28 !
29 ! kstep step number
30 ! kinc increment number
31 ! time(1) current step time
32 ! time(2) current total time
33 ! noel element number
34 ! npt integration point number
35 ! layer currently not used
36 ! kspt currently not used
37 ! coords(1..3) global coordinates of the integration point
38 ! jltyp loading face kode:
39 ! 21 = face 1
40 ! 22 = face 2
41 ! 23 = face 3
42 ! 24 = face 4
43 ! 25 = face 5
44 ! 26 = face 6
45 ! loadtype load type label
46 ! vold(0..4,1..nk) solution field in all nodes
47 ! 0: temperature
48 ! 1: displacement in global x-direction
49 ! 2: displacement in global y-direction
50 ! 3: displacement in global z-direction
51 ! 4: static pressure
52 ! veold(0..3,1..nk) derivative of the solution field w.r.t.
53 ! time in all nodes
54 ! 0: temperature rate
55 ! 1: velocity in global x-direction
56 ! 2: velocity in global y-direction
57 ! 3: velocity in global z-direction
58 ! co(3,1..nk) coordinates of all nodes
59 ! 1: coordinate in global x-direction
60 ! 2: coordinate in global y-direction
61 ! 3: coordinate in global z-direction
62 ! lakonl element label
63 ! konl(1..20) nodes belonging to the element
64 ! ipompc(1..nmpc)) ipompc(i) points to the first term of
65 ! MPC i in field nodempc
66 ! nodempc(1,*) node number of a MPC term
67 ! nodempc(2,*) coordinate direction of a MPC term
68 ! nodempc(3,*) if not 0: points towards the next term
69 ! of the MPC in field nodempc
70 ! if 0: MPC definition is finished
71 ! coefmpc(*) coefficient of a MPC term
72 ! nmpc number of MPC's
73 ! ikmpc(1..nmpc) ordered global degrees of freedom of the MPC's
74 ! the global degree of freedom is
75 ! 8*(node-1)+direction of the dependent term of
76 ! the MPC (direction = 0: temperature;
77 ! 1-3: displacements; 4: static pressure;
78 ! 5-7: rotations)
79 ! ilmpc(1..nmpc) ilmpc(i) is the MPC number corresponding
80 ! to the reference number in ikmpc(i)
81 ! rho local density
82 ! amat material name
83 ! mi(1) max # of integration points per element (max
84 ! over all elements)
85 ! mi(2) max degree of freedomm per node (max over all
86 ! nodes) in fields like v(0:mi(2))...
87 !
88 ! OUTPUT:
89 !
90 ! f magnitude of the distributed load
91 ! iscale determines whether the flux has to be
92 ! scaled for increments smaller than the
93 ! step time in static calculations
94 ! 0: no scaling
95 ! 1: scaling (default)
96 !
97  implicit none
98 !
99  character*8 lakonl
100  character*20 loadtype
101  character*80 amat
102 !
103  integer kstep,kinc,noel,npt,jltyp,layer,kspt,konl(20),iscale,
104  & mi(*)
105 !
106  real*8 f,time(2),coords(3),vold(0:mi(2),*),co(3,*),rho
107 !
108  intent(in) kstep,kinc,time,noel,npt,layer,kspt,
109  & coords,jltyp,loadtype,vold,co,lakonl,konl,
110  & ipompc,nodempc,coefmpc,nmpc,ikmpc,ilmpc,veold,
111  & rho,amat,mi
112 !
113  intent(inout) f,iscale
114 !
115 ! the code starting here up to the end of the file serves as
116 ! an example for combined mechanical-lubrication problems.
117 ! Please replace it by your own code for your concrete application.
118 !
119  integer ifaceq(8,6),ifacet(6,4),ifacew(8,5),ig,nelem,nopes,
120  & iflag,i,j,nope,ipompc(*),nodempc(3,*),nmpc,ikmpc(*),ilmpc(*),
121  & node,idof,id
122 !
123  real*8 xl2(3,8),pres(8),xi,et,xsj2(3),xs2(3,7),shp2(7,8),
124  & coefmpc(*),veold(0:mi(2),*)
125 !
126  include "gauss.f"
127 !
128  ifaceq=reshape((/4,3,2,1,11,10,9,12,
129  & 5,6,7,8,13,14,15,16,
130  & 1,2,6,5,9,18,13,17,
131  & 2,3,7,6,10,19,14,18,
132  & 3,4,8,7,11,20,15,19,
133  & 4,1,5,8,12,17,16,20/),(/8,6/))
134  ifacet=reshape((/1,3,2,7,6,5,
135  & 1,2,4,5,9,8,
136  & 2,3,4,6,10,9,
137  & 1,4,3,8,10,7/),(/6,4/))
138  ifacew=reshape((/1,3,2,9,8,7,0,0,
139  & 4,5,6,10,11,12,0,0,
140  & 1,2,5,4,7,14,10,13,
141  & 2,3,6,5,8,15,11,14,
142  & 4,6,3,1,12,15,9,13/),(/8,5/))
143  iflag=2
144 !
145  nelem=noel
146  ig=jltyp-20
147 !
148  if(lakonl(4:4).eq.'2') then
149  nope=20
150  nopes=8
151  elseif(lakonl(4:4).eq.'8') then
152  nope=8
153  nopes=4
154  elseif(lakonl(4:5).eq.'10') then
155  nope=10
156  nopes=6
157  elseif(lakonl(4:4).eq.'4') then
158  nope=4
159  nopes=3
160  elseif(lakonl(4:5).eq.'15') then
161  nope=15
162  elseif(lakonl(4:4).eq.'6') then
163  nope=6
164  endif
165 !
166 ! treatment of wedge faces
167 !
168  if(lakonl(4:4).eq.'6') then
169  if(ig.le.2) then
170  nopes=3
171  else
172  nopes=4
173  endif
174  endif
175  if(lakonl(4:5).eq.'15') then
176  if(ig.le.2) then
177  nopes=6
178  else
179  nopes=8
180  endif
181  endif
182 !
183  do i=1,nopes
184  do j=1,3
185  xl2(j,i)=0.d0
186  enddo
187  enddo
188 !
189  if((nope.eq.20).or.(nope.eq.8)) then
190  do i=1,nopes
191  node=konl(ifaceq(i,ig))
192  idof=8*(node-1)
193  call nident(ikmpc,idof,nmpc,id)
194  if((id.eq.0).or.(ikmpc(id).ne.idof)) then
195  write(*,*) '*ERROR in dload: node ',node
196  write(*,*) ' is not connected to the oil film'
197  call exit(201)
198  endif
199  node=nodempc(1,nodempc(3,ipompc(ilmpc(id))))
200  pres(i)=vold(0,node)
201  enddo
202  elseif((nope.eq.10).or.(nope.eq.4)) then
203  do i=1,nopes
204  node=konl(ifacet(i,ig))
205  node=konl(ifaceq(i,ig))
206  idof=8*(node-1)
207  call nident(ikmpc,idof,nmpc,id)
208  if((id.eq.0).or.(ikmpc(id).ne.idof)) then
209  write(*,*) '*ERROR in dload: node ',node
210  write(*,*) ' is not connected to the oil film'
211  call exit(201)
212  endif
213  node=nodempc(1,nodempc(3,ipompc(ilmpc(id))))
214  pres(i)=vold(0,node)
215  enddo
216  else
217  do i=1,nopes
218  node=konl(ifacew(i,ig))
219  node=konl(ifaceq(i,ig))
220  idof=8*(node-1)
221  call nident(ikmpc,idof,nmpc,id)
222  if((id.eq.0).or.(ikmpc(id).ne.idof)) then
223  write(*,*) '*ERROR in dload: node ',node
224  write(*,*) ' is not connected to the oil film'
225  call exit(201)
226  endif
227  node=nodempc(1,nodempc(3,ipompc(ilmpc(id))))
228  pres(i)=vold(0,node)
229  enddo
230  endif
231 !
232  i=npt
233 !
234  if((lakonl(4:5).eq.'8R').or.
235  & ((lakonl(4:4).eq.'6').and.(nopes.eq.4))) then
236  xi=gauss2d1(1,i)
237  et=gauss2d1(2,i)
238  elseif((lakonl(4:4).eq.'8').or.
239  & (lakonl(4:6).eq.'20R').or.
240  & ((lakonl(4:5).eq.'15').and.(nopes.eq.8))) then
241  xi=gauss2d2(1,i)
242  et=gauss2d2(2,i)
243  elseif(lakonl(4:4).eq.'2') then
244  xi=gauss2d3(1,i)
245  et=gauss2d3(2,i)
246  elseif((lakonl(4:5).eq.'10').or.
247  & ((lakonl(4:5).eq.'15').and.(nopes.eq.6))) then
248  xi=gauss2d5(1,i)
249  et=gauss2d5(2,i)
250  elseif((lakonl(4:4).eq.'4').or.
251  & ((lakonl(4:4).eq.'6').and.(nopes.eq.3))) then
252  xi=gauss2d4(1,i)
253  et=gauss2d4(2,i)
254  endif
255 !
256  if(nopes.eq.8) then
257  call shape8q(xi,et,xl2,xsj2,xs2,shp2,iflag)
258  elseif(nopes.eq.4) then
259  call shape4q(xi,et,xl2,xsj2,xs2,shp2,iflag)
260  elseif(nopes.eq.6) then
261  call shape6tri(xi,et,xl2,xsj2,xs2,shp2,iflag)
262  else
263  call shape3tri(xi,et,xl2,xsj2,xs2,shp2,iflag)
264  endif
265 !
266 ! determining the pressure
267 !
268  f=0.d0
269  do j=1,nopes
270  f=f+pres(j)*shp2(4,j)
271  enddo
272 !
273  iscale=0
274 !
275  return
subroutine shape8q(xi, et, xl, xsj, xs, shp, iflag)
Definition: shape8q.f:20
subroutine shape3tri(xi, et, xl, xsj, xs, shp, iflag)
Definition: shape3tri.f:20
subroutine shape4q(xi, et, xl, xsj, xs, shp, iflag)
Definition: shape4q.f:20
subroutine nident(x, px, n, id)
Definition: nident.f:26
subroutine shape6tri(xi, et, xl, xsj, xs, shp, iflag)
Definition: shape6tri.f:20
Hosted by OpenAircraft.com, (Michigan UAV, LLC)