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

Go to the source code of this file.

Functions/Subroutines

subroutine tempload_em (xforcold, xforc, xforcact, iamforc, nforc, xloadold, xload, xloadact, iamload, nload, ibody, xbody, nbody, xbodyold, xbodyact, t1old, t1, t1act, iamt1, nk, amta, namta, nam, ampli, time, reltime, ttime, dtime, ithermal, nmethod, xbounold, xboun, xbounact, iamboun, nboun, nodeboun, ndirboun, nodeforc, ndirforc, istep, iinc, co, vold, itg, ntg, amname, ikboun, ilboun, nelemload, sideload, mi, ntrans, trab, inotr, veold, integerglob, doubleglob, tieset, istartset, iendset, ialset, ntie, nmpc, ipompc, ikmpc, ilmpc, nodempc, coefmpc, h0scale, inomat, ipobody, iponoel, inoel)
 

Function/Subroutine Documentation

◆ tempload_em()

subroutine tempload_em ( real*8, dimension(*)  xforcold,
real*8, dimension(*)  xforc,
real*8, dimension(*)  xforcact,
integer, dimension(*)  iamforc,
integer  nforc,
real*8, dimension(2,*)  xloadold,
real*8, dimension(2,*)  xload,
real*8, dimension(2,*)  xloadact,
integer, dimension(2,*)  iamload,
integer  nload,
integer, dimension(3,*)  ibody,
real*8, dimension(7,*)  xbody,
integer  nbody,
real*8, dimension(7,*)  xbodyold,
real*8, dimension(7,*)  xbodyact,
real*8, dimension(*)  t1old,
real*8, dimension(*)  t1,
real*8, dimension(*)  t1act,
integer, dimension(*)  iamt1,
integer  nk,
real*8, dimension(2,*)  amta,
integer, dimension(3,*)  namta,
integer  nam,
real*8, dimension(*)  ampli,
real*8  time,
real*8  reltime,
real*8  ttime,
real*8  dtime,
integer  ithermal,
integer  nmethod,
real*8, dimension(*)  xbounold,
real*8, dimension(*)  xboun,
real*8, dimension(*)  xbounact,
integer, dimension(*)  iamboun,
integer  nboun,
integer, dimension(*)  nodeboun,
integer, dimension(*)  ndirboun,
integer, dimension(2,*)  nodeforc,
integer, dimension(*)  ndirforc,
integer  istep,
integer  iinc,
real*8, dimension(3,*)  co,
real*8, dimension(0:mi(2),*)  vold,
integer, dimension(*)  itg,
integer  ntg,
character*80, dimension(*)  amname,
integer, dimension(*)  ikboun,
integer, dimension(*)  ilboun,
integer, dimension(2,*)  nelemload,
character*20, dimension(*)  sideload,
integer, dimension(*)  mi,
integer  ntrans,
real*8, dimension(7,*)  trab,
integer, dimension(2,*)  inotr,
real*8, dimension(0:mi(2),*)  veold,
integer, dimension(*)  integerglob,
real*8, dimension(*)  doubleglob,
character*81, dimension(3,*)  tieset,
integer, dimension(*)  istartset,
integer, dimension(*)  iendset,
integer, dimension(*)  ialset,
integer  ntie,
integer  nmpc,
integer, dimension(*)  ipompc,
integer, dimension(*)  ikmpc,
integer, dimension(*)  ilmpc,
integer, dimension(3,*)  nodempc,
real*8, dimension(*)  coefmpc,
real*8  h0scale,
integer, dimension(*)  inomat,
integer, dimension(2,*)  ipobody,
integer, dimension(*)  iponoel,
integer, dimension(2,*)  inoel 
)
29 !
30 ! calculates the loading at a given time
31 !
32  implicit none
33 !
34  logical gasnode
35 !
36  character*1 entity
37  character*20 sideload(*)
38  character*80 amname(*)
39  character*81 tieset(3,*)
40 !
41  integer iamforc(*),iamload(2,*),iamt1(*),nelemload(2,*),
42  & nam,i,istart,iend,id,nforc,nload,nk,namta(3,*),ithermal,
43  & nmethod,iamt1i,iamboun(*),nboun,iamforci,iambouni,
44  & iamloadi1,iamloadi2,ibody(3,*),itg(*),ntg,idof,one,
45  & nbody,iambodyi,nodeboun(*),ndirboun(*),nodeforc(2,*),
46  & ndirforc(*),istep,iinc,msecpt,node,j,ikboun(*),ilboun(*),
47  & ipresboun,mi(*),ntrans,inotr(2,*),idummy,integerglob(*),
48  & istartset(*),iendset(*),ialset(*),ntie,iselect(1),
49  & nmpc,ikmpc(*),ilmpc(*),nodempc(3,*),k,ist,index,ipompc(*),
50  & inomat(*),ipobody(2,*),iponoel(*),inoel(2,*)
51 !
52  real*8 xforc(*),xforcact(*),xload(2,*),xloadact(2,*),
53  & t1(*),t1act(*),amta(2,*),ampli(*),time,fixed_temp,
54  & xforcold(*),xloadold(2,*),t1old(*),reltime,coefmpc(*),
55  & xbounold(*),xboun(*),xbounact(*),ttime,dtime,reftime,
56  & xbody(7,*),xbodyold(7,*),xbodyact(7,*),co(3,*),
57  & vold(0:mi(2),*),abqtime(2),coords(3),trab(7,*),
58  & veold(0:mi(2),*),ddummy,doubleglob(*),h0scale
59 !
60  data msecpt /1/
61  data one /1/
62 !
63  h0scale=1.d0
64 !
65 ! if an amplitude is active, the loading is scaled according to
66 ! the actual time. If no amplitude is active, then the load is
67 ! - scaled according to the relative time for a static step
68 ! - applied as a step loading for a dynamic step
69 !
70 ! calculating all amplitude values for the current time
71 !
72  do i=1,nam
73  if(namta(3,i).lt.0) then
74  reftime=ttime+time
75  else
76  reftime=time
77  endif
78  if(abs(namta(3,i)).ne.i) then
79  reftime=reftime-amta(1,namta(1,i))
80  istart=namta(1,abs(namta(3,i)))
81  iend=namta(2,abs(namta(3,i)))
82  if(istart.eq.0) then
83  call uamplitude(reftime,amname(namta(3,i)),ampli(i))
84  cycle
85  endif
86  else
87  istart=namta(1,i)
88  iend=namta(2,i)
89  if(istart.eq.0) then
90  call uamplitude(reftime,amname(i),ampli(i))
91  cycle
92  endif
93  endif
94  call identamta(amta,reftime,istart,iend,id)
95  if(id.lt.istart) then
96  ampli(i)=amta(2,istart)
97  elseif(id.eq.iend) then
98  ampli(i)=amta(2,iend)
99  else
100  ampli(i)=amta(2,id)+(amta(2,id+1)-amta(2,id))
101  & *(reftime-amta(1,id))/(amta(1,id+1)-amta(1,id))
102  endif
103  enddo
104 !
105 ! scaling the boundary conditions
106 !
107  do i=1,nboun
108  if((xboun(i).lt.1.2357111318d0).and.
109  & (xboun(i).gt.1.2357111316d0)) then
110 !
111 ! user subroutine for boundary conditions
112 !
113  node=nodeboun(i)
114 !
115 ! check whether node is a gasnode
116 !
117  gasnode=.false.
118  call nident(itg,node,ntg,id)
119  if(id.gt.0) then
120  if(itg(id).eq.node) then
121  gasnode=.true.
122  endif
123  endif
124 !
125  abqtime(1)=time
126  abqtime(2)=ttime+time
127 !
128 ! a gasnode cannot move (displacement DOFs are used
129 ! for other purposes, e.g. mass flow and pressure)
130 !
131  if(gasnode) then
132  do j=1,3
133  coords(j)=co(j,node)
134  enddo
135  else
136  do j=1,3
137  coords(j)=co(j,node)+vold(j,node)
138  enddo
139  endif
140 !
141  if(ndirboun(i).eq.0) then
142  call utemp(xbounact(i),msecpt,istep,iinc,abqtime,node,
143  & coords,vold,mi,iponoel,inoel,
144  & ipobody,xbodyact,ibody)
145  else
146  call uboun(xbounact(i),istep,iinc,abqtime,node,
147  & ndirboun(i),coords,vold,mi,iponoel,inoel,
148  & ipobody,xbodyact,ibody)
149  endif
150  cycle
151  endif
152  if((xboun(i).lt.1.9232931375d0).and.
153  & (xboun(i).gt.1.9232931373d0)) then
154 !
155 ! boundary conditions for submodel
156 !
157  node=nodeboun(i)
158 !
159 ! check whether node is a gasnode
160 !
161  gasnode=.false.
162  call nident(itg,node,ntg,id)
163  if(id.gt.0) then
164  if(itg(id).eq.node) then
165  gasnode=.true.
166  endif
167  endif
168 !
169 ! for the interpolation of submodels the undeformed
170 ! geometry is taken
171 !
172  do j=1,3
173  coords(j)=co(j,node)
174  enddo
175 !
176  entity='N'
177  one=1
178  iselect(1)=ndirboun(i)+1
179  call interpolsubmodel(integerglob,doubleglob,xbounact(i),
180  & coords,iselect,one,node,tieset,istartset,iendset,
181  & ialset,ntie,entity)
182 c write(*,*) 'tempload ',node,ndirboun(i),xbounact(i)
183 !
184  if(nmethod.eq.1) then
185  xbounact(i)=xbounold(i)+
186  & (xbounact(i)-xbounold(i))*reltime
187  endif
188  cycle
189  endif
190 !
191  if(nam.gt.0) then
192  iambouni=iamboun(i)
193  else
194  iambouni=0
195  endif
196  if(iambouni.gt.0) then
197  xbounact(i)=xboun(i)*ampli(iambouni)
198 !
199 ! scaling of the voltage across a coil in an
200 ! electromagnetic calculation
201 !
202  node=nodeboun(i)
203  if(inomat(node).eq.0) h0scale=ampli(iambouni)
204 !
205  elseif(nmethod.eq.1) then
206  xbounact(i)=xbounold(i)+
207  & (xboun(i)-xbounold(i))*reltime
208  else
209  xbounact(i)=xboun(i)
210  endif
211  enddo
212 !
213 ! scaling the loading
214 !
215  do i=1,nforc
216  if(ndirforc(i).eq.0) then
217  if((xforc(i).lt.1.2357111318d0).and.
218  & (xforc(i).gt.1.2357111316d0)) then
219 !
220 ! user subroutine for the concentrated heat flux
221 !
222  node=nodeforc(1,i)
223 !
224 ! check whether node is a gasnode
225 !
226  gasnode=.false.
227  call nident(itg,node,ntg,id)
228  if(id.gt.0) then
229  if(itg(id).eq.node) then
230  gasnode=.true.
231  endif
232  endif
233 !
234  abqtime(1)=time
235  abqtime(2)=ttime+time
236 !
237 ! a gasnode cannot move (displacement DOFs are used
238 ! for other purposes, e.g. mass flow and pressure)
239 !
240  if(gasnode) then
241  do j=1,3
242  coords(j)=co(j,node)
243  enddo
244  else
245  do j=1,3
246  coords(j)=co(j,node)+vold(j,node)
247  enddo
248  endif
249 !
250  call cflux(xforcact(i),msecpt,istep,iinc,abqtime,node,
251  & coords,vold,mi)
252  cycle
253  endif
254  else
255  if((xforc(i).lt.1.2357111318d0).and.
256  & (xforc(i).gt.1.2357111316d0)) then
257 !
258 ! user subroutine for the concentrated load
259 !
260  node=nodeforc(1,i)
261 !
262  abqtime(1)=time
263  abqtime(2)=ttime+time
264 !
265  do j=1,3
266  coords(j)=co(j,node)+vold(j,node)
267  enddo
268 !
269  call cload(xforcact(i),istep,iinc,abqtime,node,
270  & ndirforc(i),coords,vold,mi,ntrans,trab,inotr,veold)
271  cycle
272  elseif((xforc(i).lt.1.9232931375d0).and.
273  & (xforc(i).gt.1.9232931373d0)) then
274 !
275 ! boundary conditions for submodel
276 !
277  node=nodeforc(1,i)
278 !
279 ! for the interpolation of submodels the undeformed
280 ! geometry is taken
281 !
282  do j=1,3
283  coords(j)=co(j,node)
284  enddo
285 !
286  entity='N'
287  one=1
288  iselect(1)=ndirforc(i)+10
289  call interpolsubmodel(integerglob,doubleglob,xforcact(i),
290  & coords,iselect,one,node,tieset,istartset,iendset,
291  & ialset,ntie,entity)
292 !
293  if(nmethod.eq.1) then
294  xforcact(i)=xforcold(i)+
295  & (xforcact(i)-xforcold(i))*reltime
296  endif
297  cycle
298  endif
299  endif
300 !
301  if(nam.gt.0) then
302  iamforci=iamforc(i)
303  else
304  iamforci=0
305  endif
306  if(iamforci.gt.0) then
307  xforcact(i)=xforc(i)*ampli(iamforci)
308  elseif(nmethod.eq.1) then
309  xforcact(i)=xforcold(i)+
310  & (xforc(i)-xforcold(i))*reltime
311  else
312  xforcact(i)=xforc(i)
313  endif
314  enddo
315 !
316  do i=1,nload
317  ipresboun=0
318 !
319 ! check for pressure boundary conditions
320 !
321  if(sideload(i)(3:4).eq.'NP') then
322  node=nelemload(2,i)
323  idof=8*(node-1)+2
324  call nident(ikboun,idof,nboun,id)
325  if(id.gt.0) then
326  if(ikboun(id).eq.idof) then
327  ipresboun=1
328  xloadact(1,i)=xbounact(ilboun(id))
329  endif
330  endif
331  endif
332 !
333  if(ipresboun.eq.0) then
334  if(nam.gt.0) then
335  iamloadi1=iamload(1,i)
336  iamloadi2=iamload(2,i)
337  else
338  iamloadi1=0
339  iamloadi2=0
340  endif
341  if(iamloadi1.gt.0) then
342  xloadact(1,i)=xload(1,i)*ampli(iamloadi1)
343  elseif(nmethod.eq.1) then
344  xloadact(1,i)=xloadold(1,i)+
345  & (xload(1,i)-xloadold(1,i))*reltime
346  else
347  xloadact(1,i)=xload(1,i)
348  endif
349  if(iamloadi2.gt.0) then
350  xloadact(2,i)=xload(2,i)*ampli(iamloadi2)
351 c elseif(nmethod.eq.1) then
352 c xloadact(2,i)=xload(2,i)
353  else
354  xloadact(2,i)=xload(2,i)
355  endif
356  endif
357  enddo
358 !
359  do i=1,nbody
360  if(nam.gt.0) then
361  iambodyi=ibody(2,i)
362  else
363  iambodyi=0
364  endif
365  if(iambodyi.gt.0) then
366  xbodyact(1,i)=xbody(1,i)*ampli(iambodyi)
367  elseif(nmethod.eq.1) then
368  xbodyact(1,i)=xbodyold(1,i)+
369  & (xbody(1,i)-xbodyold(1,i))*reltime
370  else
371  xbodyact(1,i)=xbody(1,i)
372  endif
373  enddo
374 !
375 ! scaling the temperatures
376 !
377  if(ithermal.eq.1) then
378  do i=1,nk
379  if((t1(i).lt.1.2357111318d0).and.
380  & (t1(i).gt.1.2357111316d0)) then
381 !
382  abqtime(1)=time
383  abqtime(2)=ttime+time
384 !
385  do j=1,3
386  coords(j)=co(j,i)+vold(j,i)
387  enddo
388  call utemp(t1act(i),msecpt,istep,iinc,abqtime,i,
389  & coords,vold,mi,iponoel,inoel,
390  & ipobody,xbodyact,ibody)
391  cycle
392  endif
393 !
394  if((t1(i).lt.1.9232931375d0).and.
395  & (t1(i).gt.1.9232931373d0)) then
396 !
397 ! for the interpolation of submodels the undeformed
398 ! geometry is taken
399 !
400  do j=1,3
401  coords(j)=co(j,i)
402  enddo
403 !
404  entity='N'
405  one=1
406  iselect(1)=1
407  call interpolsubmodel(integerglob,doubleglob,t1act(i),
408  & coords,iselect,one,i,tieset,istartset,iendset,
409  & ialset,ntie,entity)
410 !
411  if(nmethod.eq.1) then
412  t1act(i)=t1old(i)+(t1act(i)-t1old(i))*reltime
413  endif
414  cycle
415  endif
416 !
417  if(nam.gt.0) then
418  iamt1i=iamt1(i)
419  else
420  iamt1i=0
421  endif
422  if(iamt1i.gt.0) then
423  t1act(i)=t1(i)*ampli(iamt1i)
424  elseif(nmethod.eq.1) then
425  t1act(i)=t1old(i)+(t1(i)-t1old(i))*reltime
426  else
427  t1act(i)=t1(i)
428  endif
429  enddo
430 !
431 ! taking temperature MPC's into account
432 !
433  do j=1,nmpc
434  k=mod(ikmpc(j),8)
435  if(k.ne.0) cycle
436  i=ilmpc(j)
437  ist=ipompc(i)
438  node=nodempc(1,ist)
439  index=nodempc(3,ist)
440  fixed_temp=0.d0
441  if(index.ne.0) then
442  do
443  fixed_temp=fixed_temp-
444  & coefmpc(index)*t1act(nodempc(1,index))
445  index=nodempc(3,index)
446  if(index.eq.0) exit
447  enddo
448  endif
449  t1act(node)=fixed_temp/coefmpc(ist)
450  enddo
451  endif
452 c write(*,*) 'nboun'
453 c do i=1,nboun
454 c write(*,'(i7,1x,e11.4,1x,e11.4)') i,xbounact(i),xboun(i)
455 c enddo
456 c write(*,*) 'nforc'
457 c do i=1,nforc
458 c write(*,'(i7,1x,e11.4,1x,e11.4)') i,xforcact(i),xforc(i)
459 c enddo
460 c write(*,*) 'nload'
461 c do i=1,nload
462 c write(*,'(i7,1x,e11.4,1x,e11.4)') i,xloadact(1,i),xload(1,i)
463 c enddo
464 !
465  return
subroutine uamplitude(time, name, amplitude)
Definition: uamplitude.f:20
subroutine cflux(flux, msecpt, kstep, kinc, time, node, coords, vold, mi)
Definition: cflux.f:21
subroutine utemp(temp, msecpt, kstep, kinc, time, node, coords, vold, mi, iponoel, inoel, ipobody, xbody, ibody)
Definition: utemp.f:21
subroutine identamta(amta, reftime, istart, iend, id)
Definition: identamta.f:26
subroutine nident(x, px, n, id)
Definition: nident.f:26
subroutine uboun(boun, kstep, kinc, time, node, idof, coords, vold, mi, iponoel, inoel, ipobody, xbody, ibody)
Definition: uboun.f:21
subroutine cload(xload, kstep, kinc, time, node, idof, coords, vold, mi, ntrans, trab, inotr, veold)
Definition: cload.f:21
subroutine interpolsubmodel(integerglob, doubleglob, value, coords, iselect, nselect, nodeface, tieset, istartset, iendset, ialset, ntie, entity)
Definition: interpolsubmodel.f:22
Hosted by OpenAircraft.com, (Michigan UAV, LLC)