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

Go to the source code of this file.

Functions/Subroutines

subroutine heattransfers (inpc, textpart, nmethod, iperturb, isolver, istep, istat, n, tinc, tper, tmin, tmax, idrct, ithermal, iline, ipol, inl, ipoinp, inp, alpha, mei, fei, ipoinpc, ctrl, ttime)
 

Function/Subroutine Documentation

◆ heattransfers()

subroutine heattransfers ( character*1, dimension(*)  inpc,
character*132, dimension(16)  textpart,
integer  nmethod,
integer  iperturb,
integer  isolver,
integer  istep,
integer  istat,
integer  n,
real*8  tinc,
real*8  tper,
real*8  tmin,
real*8  tmax,
integer  idrct,
integer  ithermal,
integer  iline,
integer  ipol,
integer  inl,
integer, dimension(2,*)  ipoinp,
integer, dimension(3,*)  inp,
real*8  alpha,
integer, dimension(4)  mei,
real*8, dimension(3)  fei,
integer, dimension(0:*)  ipoinpc,
real*8, dimension(*)  ctrl,
real*8  ttime 
)
22 !
23 ! reading the input deck: *HEAT TRANSFER
24 !
25 ! isolver=0: SPOOLES
26 ! 2: iterative solver with diagonal scaling
27 ! 3: iterative solver with Cholesky preconditioning
28 ! 4: sgi solver
29 ! 5: TAUCS
30 ! 7: pardiso
31 !
32  implicit none
33 !
34  logical timereset
35 !
36  character*1 inpc(*)
37  character*20 solver
38  character*132 textpart(16)
39 !
40  integer nmethod,iperturb,isolver,istep,istat,n,key,i,idrct,nev,
41  & ithermal,iline,ipol,inl,ipoinp(2,*),inp(3,*),mei(4),ncv,mxiter,
42  & ipoinpc(0:*),idirect
43 !
44  real*8 tinc,tper,tmin,tmax,alpha,fei(3),tol,fmin,fmax,ctrl(*),
45  & ttime
46 !
47  tmin=0.d0
48  tmax=0.d0
49  nmethod=4
50  alpha=0.d0
51  mei(4)=0
52  timereset=.false.
53 !
54 ! defaults for fmin and fmax
55 !
56  fmin=-1.d0
57  fmax=-1.d0
58 !
59  if(iperturb.eq.0) then
60  iperturb=2
61  elseif((iperturb.eq.1).and.(istep.gt.1)) then
62  write(*,*) '*ERROR in heattransfers: perturbation analysis is'
63  write(*,*) ' not provided in a *HEAT TRANSFER step.'
64  call exit(201)
65  endif
66 !
67  if(istep.lt.1) then
68  write(*,*) '*ERROR in heattransfers: *HEAT TRANSFER can only'
69  write(*,*) ' be used within a STEP'
70  call exit(201)
71  endif
72 !
73 ! default solver
74 !
75  solver=' '
76  if(isolver.eq.0) then
77  solver(1:7)='SPOOLES'
78  elseif(isolver.eq.2) then
79  solver(1:16)='ITERATIVESCALING'
80  elseif(isolver.eq.3) then
81  solver(1:17)='ITERATIVECHOLESKY'
82  elseif(isolver.eq.4) then
83  solver(1:3)='SGI'
84  elseif(isolver.eq.5) then
85  solver(1:5)='TAUCS'
86  elseif(isolver.eq.7) then
87  solver(1:7)='PARDISO'
88  endif
89 !
90  idirect=2
91  do i=2,n
92  if(textpart(i)(1:7).eq.'SOLVER=') then
93  read(textpart(i)(8:27),'(a20)') solver
94  elseif((textpart(i)(1:6).eq.'DIRECT').and.
95  & (textpart(i)(1:9).ne.'DIRECT=NO')) then
96  idirect=1
97  elseif(textpart(i)(1:9).eq.'DIRECT=NO') then
98  idirect=0
99  elseif(textpart(i)(1:11).eq.'STEADYSTATE') then
100  nmethod=1
101  elseif(textpart(i)(1:9).eq.'FREQUENCY') then
102  nmethod=2
103  elseif(textpart(i)(1:12).eq.'MODALDYNAMIC') then
104  iperturb=0
105  elseif(textpart(i)(1:11).eq.'STORAGE=YES') then
106  mei(4)=1
107  elseif(textpart(i)(1:7).eq.'DELTMX=') then
108  read(textpart(i)(8:27),'(f20.0)',iostat=istat) ctrl(27)
109  elseif(textpart(i)(1:9).eq.'TIMERESET') then
110  timereset=.true.
111  elseif(textpart(i)(1:17).eq.'TOTALTIMEATSTART=') then
112  read(textpart(i)(18:37),'(f20.0)',iostat=istat) ttime
113  else
114  write(*,*)
115  & '*WARNING in heattransfers: parameter not recognized:'
116  write(*,*) ' ',
117  & textpart(i)(1:index(textpart(i),' ')-1)
118  call inputwarning(inpc,ipoinpc,iline,
119  &"*HEAT TRANSFER%")
120  endif
121  enddo
122  if(nmethod.eq.1) ctrl(27)=1.d30
123 !
124 ! default for modal dynamic calculations is DIRECT,
125 ! for static or dynamic calculations DIRECT=NO
126 !
127  if(iperturb.eq.0) then
128  idrct=1
129  if(idirect.eq.0)idrct=0
130  else
131  idrct=0
132  if(idirect.eq.1)idrct=1
133  endif
134 !
135  if((ithermal.eq.0).and.(nmethod.ne.1).and.
136  & (nmethod.ne.2).and.(iperturb.ne.0)) then
137  write(*,*) '*ERROR in heattransfers: please define initial '
138  write(*,*) ' conditions for the temperature'
139  call exit(201)
140  else
141  ithermal=2
142  endif
143 !
144  if((nmethod.ne.2).and.(iperturb.ne.0)) then
145 !
146 ! static or dynamic thermal analysis
147 !
148  if(solver(1:7).eq.'SPOOLES') then
149  isolver=0
150  elseif(solver(1:16).eq.'ITERATIVESCALING') then
151  isolver=2
152  elseif(solver(1:17).eq.'ITERATIVECHOLESKY') then
153  isolver=3
154  elseif(solver(1:3).eq.'SGI') then
155  isolver=4
156  elseif(solver(1:5).eq.'TAUCS') then
157  isolver=5
158  elseif(solver(1:7).eq.'PARDISO') then
159  isolver=7
160  else
161  write(*,*) '*WARNING in heattransfers: unknown solver;'
162  write(*,*) ' the default solver is used'
163  endif
164 !
165  call getnewline(inpc,textpart,istat,n,key,iline,ipol,inl,
166  & ipoinp,inp,ipoinpc)
167  if((istat.lt.0).or.(key.eq.1)) then
168  if(iperturb.ge.2) then
169  write(*,*) '*WARNING in heattransfers: a nonlinear geomet
170  &ric analysis is requested'
171  write(*,*) ' but no time increment nor step is sp
172  &ecified'
173  write(*,*) ' the defaults (1,1) are used'
174  tinc=1.d0
175  tper=1.d0
176  tmin=1.d-5
177  tmax=1.d+30
178  endif
179  if(timereset)ttime=ttime-tper
180  return
181  endif
182 !
183  read(textpart(1)(1:20),'(f20.0)',iostat=istat) tinc
184  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
185  &"*HEAT TRANSFER%")
186  read(textpart(2)(1:20),'(f20.0)',iostat=istat) tper
187  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
188  &"*HEAT TRANSFER%")
189  read(textpart(3)(1:20),'(f20.0)',iostat=istat) tmin
190  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
191  &"*HEAT TRANSFER%")
192  read(textpart(4)(1:20),'(f20.0)',iostat=istat) tmax
193  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
194  &"*HEAT TRANSFER%")
195 !
196  if(tinc.le.0.d0) then
197  write(*,*) '*ERROR in heattransfers: initial increment size
198  &is negative'
199  endif
200  if(tper.le.0.d0) then
201  write(*,*) '*ERROR in heattransfers: step size is negative'
202  endif
203  if(tinc.gt.tper) then
204  write(*,*) '*ERROR in heattransfers: initial increment size
205  &exceeds step size'
206  endif
207 !
208  if(idrct.ne.1) then
209 c if(dabs(tmin).lt.1.d-10) then
210 c tmin=min(tinc,1.d-5*tper)
211  if(dabs(tmin).lt.1.d-6*tper) then
212  tmin=min(tinc,1.d-6*tper)
213  endif
214  if(dabs(tmax).lt.1.d-10) then
215  tmax=1.d+30
216  endif
217  endif
218  elseif(nmethod.eq.2) then
219 !
220 ! thermal eigenmode analysis
221 !
222  iperturb=0
223 !
224  call getnewline(inpc,textpart,istat,n,key,iline,ipol,inl,
225  & ipoinp,inp,ipoinpc)
226  if((istat.lt.0).or.(key.eq.1)) then
227  write(*,*)'*ERROR in heattransfers: definition not complete'
228  write(*,*) ' '
229  call inputerror(inpc,ipoinpc,iline,
230  &"*HEAT TRANSFER%")
231  call exit(201)
232  endif
233  read(textpart(1)(1:10),'(i10)',iostat=istat) nev
234  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
235  &"*HEAT TRANSFER%")
236  if(nev.le.0) then
237  write(*,*) '*ERROR in frequencies: less than 1 eigenvalue re
238  &quested'
239  call exit(201)
240  endif
241  tol=1.d-2
242  ncv=4*nev
243  ncv=ncv+nev
244  mxiter=1000
245  if(textpart(2)(1:1).ne.' ') then
246  read(textpart(2)(1:20),'(f20.0)',iostat=istat) fmin
247  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
248  &"*HEAT TRANSFER%")
249  endif
250  if(textpart(3)(1:1).ne.' ') then
251  read(textpart(3)(1:20),'(f20.0)',iostat=istat) fmax
252  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
253  &"*HEAT TRANSFER%")
254  endif
255 !
256  mei(1)=nev
257  mei(2)=ncv
258  mei(3)=mxiter
259  fei(1)=tol
260  fei(2)=fmin
261  fei(3)=fmax
262  else
263 !
264 ! modal dynamic analysis for variables which satisfy the
265 ! Helmholtz equation
266 !
267  call getnewline(inpc,textpart,istat,n,key,iline,ipol,inl,
268  & ipoinp,inp,ipoinpc)
269  if((istat.lt.0).or.(key.eq.1)) then
270  write(*,*)'*ERROR in heattransfers: definition not complete'
271  write(*,*) ' '
272  call inputerror(inpc,ipoinpc,iline,
273  &"*HEAT TRANSFER%")
274  call exit(201)
275  endif
276  read(textpart(1)(1:20),'(f20.0)',iostat=istat) tinc
277  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
278  &"*HEAT TRANSFER%")
279  read(textpart(2)(1:20),'(f20.0)',iostat=istat) tper
280  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
281  &"*HEAT TRANSFER%")
282  endif
283 !
284  if(timereset)ttime=ttime-tper
285 !
286  call getnewline(inpc,textpart,istat,n,key,iline,ipol,inl,
287  & ipoinp,inp,ipoinpc)
288 !
289  return
subroutine inputwarning(inpc, ipoinpc, iline, text)
Definition: inputwarning.f:20
#define min(a, b)
Definition: cascade.c:31
subroutine getnewline(inpc, textpart, istat, n, key, iline, ipol, inl, ipoinp, inp, ipoinpc)
Definition: getnewline.f:21
subroutine inputerror(inpc, ipoinpc, iline, text)
Definition: inputerror.f:20
Hosted by OpenAircraft.com, (Michigan UAV, LLC)