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

Go to the source code of this file.

Functions/Subroutines

subroutine cfds (inpc, textpart, nmethod, iperturb, isolver, istep, istat, n, tinc, tper, tmin, tmax, idrct, ithermal, iline, ipol, inl, ipoinp, inp, ipoinpc, alpha, ctrl, iexpl, tincf, ttime, physcon)
 

Function/Subroutine Documentation

◆ cfds()

subroutine cfds ( 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,
integer, dimension(0:*)  ipoinpc,
real*8  alpha,
real*8, dimension(*)  ctrl,
integer  iexpl,
real*8  tincf,
real*8  ttime,
real*8, dimension(*)  physcon 
)
22 !
23 ! reading the input deck: *CFD
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 ! iexpl==0: structure:implicit, fluid:incompressible
33 ! iexpl==1: structure:implicit, fluid:compressible
34 !
35  implicit none
36 !
37  logical timereset
38 !
39  character*1 inpc(*)
40  character*20 solver
41  character*132 textpart(16)
42 !
43  integer nmethod,iperturb,isolver,istep,istat,n,key,i,idrct,
44  & ithermal,iline,ipol,inl,ipoinp(2,*),inp(3,*),ipoinpc(0:*),
45  & iexpl
46 !
47  real*8 tinc,tper,tmin,tmax,alpha,ctrl(*),tincf,ttime,physcon(*)
48 !
49  idrct=0
50  alpha=-0.05d0
51  tmin=0.d0
52  tmax=0.d0
53  tincf=-1.d0
54  nmethod=4
55  timereset=.false.
56  physcon(9)=0.5d0
57 !
58  if(iperturb.eq.0) then
59  iperturb=2
60  elseif((iperturb.eq.1).and.(istep.gt.1)) then
61  write(*,*) '*ERROR reading *CFD: perturbation analysis is'
62  write(*,*) ' not provided in a *HEAT TRANSFER step.'
63  call exit(201)
64  endif
65 !
66  if(istep.lt.1) then
67  write(*,*) '*ERROR reading *CFD: *HEAT TRANSFER can only '
68  write(*,*) ' be used within a STEP'
69  call exit(201)
70  endif
71 !
72 ! default solver
73 !
74  solver=' '
75  if(isolver.eq.0) then
76  solver(1:7)='SPOOLES'
77  elseif(isolver.eq.2) then
78  solver(1:16)='ITERATIVESCALING'
79  elseif(isolver.eq.3) then
80  solver(1:17)='ITERATIVECHOLESKY'
81  elseif(isolver.eq.4) then
82  solver(1:3)='SGI'
83  elseif(isolver.eq.5) then
84  solver(1:5)='TAUCS'
85  elseif(isolver.eq.7) then
86  solver(1:7)='PARDISO'
87  endif
88 !
89  do i=2,n
90  if(textpart(i)(1:7).eq.'SOLVER=') then
91  read(textpart(i)(8:27),'(a20)') solver
92  elseif(textpart(i)(1:12).eq.'COMPRESSIBLE') then
93  iexpl=1
94  elseif(textpart(i)(1:11).eq.'STEADYSTATE') then
95  nmethod=1
96  elseif(textpart(i)(1:9).eq.'TIMERESET') then
97  timereset=.true.
98  elseif(textpart(i)(1:17).eq.'TOTALTIMEATSTART=') then
99  read(textpart(i)(18:37),'(f20.0)',iostat=istat) ttime
100  elseif(textpart(i)(1:16).eq.'TURBULENCEMODEL=') then
101 !
102 ! turbulence model
103 !
104  if(textpart(i)(17:25).eq.'NONE') then
105  physcon(9)=0.5d0
106  elseif(textpart(i)(17:25).eq.'K-EPSILON') then
107  physcon(9)=1.5d0
108  elseif(textpart(i)(17:23).eq.'K-OMEGA') then
109  physcon(9)=2.5d0
110  elseif(textpart(i)(17:19).eq.'SST') then
111  physcon(9)=3.5d0
112  endif
113  else
114  write(*,*)
115  & '*WARNING reading *CFD: parameter not recognized:'
116  write(*,*) ' ',
117  & textpart(i)(1:index(textpart(i),' ')-1)
118  call inputwarning(inpc,ipoinpc,iline,
119  &"*CFD%")
120  endif
121  enddo
122  if(nmethod.eq.1) ctrl(27)=1.d30
123 !
124  if((ithermal.eq.0).and.(iexpl.eq.1)) then
125  write(*,*) '*ERROR reading *CFD: please define initial '
126  write(*,*) ' conditions for the temperature'
127  call exit(201)
128  elseif(ithermal.gt.0) then
129  ithermal=3
130  endif
131 !
132  if(solver(1:7).eq.'SPOOLES') then
133  isolver=0
134  elseif(solver(1:16).eq.'ITERATIVESCALING') then
135  isolver=2
136  elseif(solver(1:17).eq.'ITERATIVECHOLESKY') then
137  isolver=3
138  elseif(solver(1:3).eq.'SGI') then
139  isolver=4
140  elseif(solver(1:5).eq.'TAUCS') then
141  isolver=5
142  elseif(solver(1:7).eq.'PARDISO') then
143  isolver=7
144  else
145  write(*,*) '*WARNING reading *CFD: unknown solver;'
146  write(*,*) ' the default solver is used'
147  endif
148 !
149  call getnewline(inpc,textpart,istat,n,key,iline,ipol,inl,
150  & ipoinp,inp,ipoinpc)
151  if((istat.lt.0).or.(key.eq.1)) then
152  if(nmethod.eq.4) then
153 !
154 ! transient cfd
155 !
156  write(*,*) '*WARNING reading *CFD: a nonlinear geometric
157  & analysis is requested'
158  write(*,*) ' but no time increment nor step is speci
159  &fied'
160  write(*,*) ' the defaults (1,1) are used'
161  tinc=1.d0
162  tper=1.d0
163  tmin=1.d-5
164  tmax=1.d+30
165  tincf=-1.d0
166  elseif(nmethod.eq.1) then
167 !
168 ! steady state cfd: initial increment time and step
169 ! time are set to large values since only the first
170 ! increment is calculated
171 !
172  tinc=1.d+30
173  tper=1.d+30
174  tmin=1.d-5
175  tmax=1.d+30
176  tincf=-1.d0
177  endif
178  if(timereset)ttime=ttime-tper
179  return
180  endif
181 !
182  read(textpart(1)(1:20),'(f20.0)',iostat=istat) tinc
183  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
184  &"*CFD%")
185  read(textpart(2)(1:20),'(f20.0)',iostat=istat) tper
186  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
187  &"*CFD%")
188  read(textpart(3)(1:20),'(f20.0)',iostat=istat) tmin
189  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
190  &"*CFD%")
191  read(textpart(4)(1:20),'(f20.0)',iostat=istat) tmax
192  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
193  &"*CFD%")
194  read(textpart(5)(1:20),'(f20.0)',iostat=istat) tincf
195  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
196  &"*CFD%")
197 !
198  if(tinc.le.0.d0) then
199  write(*,*) '*ERROR reading *CFD: initial increment size is
200  &negative'
201  endif
202  if(tper.le.0.d0) then
203  write(*,*) '*ERROR reading *CFD: step size is negative'
204  endif
205  if(tinc.gt.tper) then
206  write(*,*) '*ERROR reading *CFD: initial increment size exc
207  &eeds step size'
208  endif
209 !
210  if(idrct.ne.1) then
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  if(tinc.gt.dabs(tmax)) then
218  write(*,*) '*WARNING reading *CFD:'
219  write(*,*) ' the initial increment ',tinc
220  write(*,*) ' exceeds the maximum increment ',
221  & tmax
222  write(*,*) ' the initial increment is reduced'
223  write(*,*) ' to the maximum value'
224  tinc=dabs(tmax)
225  endif
226  endif
227 !
228  if(timereset)ttime=ttime-tper
229 !
230  call getnewline(inpc,textpart,istat,n,key,iline,ipol,inl,
231  & ipoinp,inp,ipoinpc)
232 !
233  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)