1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220
|
!
! CalculiX - A 3-dimensional finite element program
! Copyright (C) 1998-2015 Guido Dhondt
!
! This program is free software; you can redistribute it and/or
! modify it under the terms of the GNU General Public License as
! published by the Free Software Foundation(version 2);
!
!
! This program is distributed in the hope that it will be useful,
! but WITHOUT ANY WARRANTY; without even the implied warranty of
! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
! GNU General Public License for more details.
!
! You should have received a copy of the GNU General Public License
! along with this program; if not, write to the Free Software
! Foundation, Inc., 675 Mass Ave, Cambridge, MA 02139, USA.
!
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)
!
! reading the input deck: *CFD
!
! isolver=0: SPOOLES
! 2: iterative solver with diagonal scaling
! 3: iterative solver with Cholesky preconditioning
! 4: sgi solver
! 5: TAUCS
! 7: pardiso
!
! iexpl==0: structure:implicit, fluid:incompressible
! iexpl==1: structure:implicit, fluid:compressible
!
implicit none
!
logical timereset
!
character*1 inpc(*)
character*20 solver
character*132 textpart(16)
!
integer nmethod,iperturb,isolver,istep,istat,n,key,i,idrct,
& ithermal,iline,ipol,inl,ipoinp(2,*),inp(3,*),ipoinpc(0:*),
& iexpl
!
real*8 tinc,tper,tmin,tmax,alpha,ctrl(*),tincf,ttime
!
idrct=0
alpha=-0.05d0
tmin=0.d0
tmax=0.d0
tincf=-1.d0
nmethod=4
timereset=.false.
!
if(iperturb.eq.0) then
iperturb=2
elseif((iperturb.eq.1).and.(istep.gt.1)) then
write(*,*) '*ERROR reading *CFD: perturbation analysis is'
write(*,*) ' not provided in a *HEAT TRANSFER step.'
call exit(201)
endif
!
if(istep.lt.1) then
write(*,*) '*ERROR reading *CFD: *HEAT TRANSFER can only '
write(*,*) ' be used within a STEP'
call exit(201)
endif
!
! default solver
!
solver=' '
if(isolver.eq.0) then
solver(1:7)='SPOOLES'
elseif(isolver.eq.2) then
solver(1:16)='ITERATIVESCALING'
elseif(isolver.eq.3) then
solver(1:17)='ITERATIVECHOLESKY'
elseif(isolver.eq.4) then
solver(1:3)='SGI'
elseif(isolver.eq.5) then
solver(1:5)='TAUCS'
elseif(isolver.eq.7) then
solver(1:7)='PARDISO'
endif
!
do i=2,n
if(textpart(i)(1:7).eq.'SOLVER=') then
read(textpart(i)(8:27),'(a20)') solver
elseif(textpart(i)(1:12).eq.'COMPRESSIBLE') then
iexpl=1
elseif(textpart(i)(1:11).eq.'STEADYSTATE') then
nmethod=1
elseif(textpart(i)(1:9).eq.'TIMERESET') then
timereset=.true.
elseif(textpart(i)(1:17).eq.'TOTALTIMEATSTART=') then
read(textpart(i)(18:37),'(f20.0)',iostat=istat) ttime
else
write(*,*)
& '*WARNING reading *CFD: parameter not recognized:'
write(*,*) ' ',
& textpart(i)(1:index(textpart(i),' ')-1)
call inputwarning(inpc,ipoinpc,iline,
&"*CFD%")
endif
enddo
if(nmethod.eq.1) ctrl(27)=1.d30
!
if((ithermal.eq.0).and.(iexpl.eq.1)) then
write(*,*) '*ERROR reading *CFD: please define initial '
write(*,*) ' conditions for the temperature'
call exit(201)
elseif(ithermal.gt.0) then
ithermal=3
endif
!
if(solver(1:7).eq.'SPOOLES') then
isolver=0
elseif(solver(1:16).eq.'ITERATIVESCALING') then
isolver=2
elseif(solver(1:17).eq.'ITERATIVECHOLESKY') then
isolver=3
elseif(solver(1:3).eq.'SGI') then
isolver=4
elseif(solver(1:5).eq.'TAUCS') then
isolver=5
elseif(solver(1:7).eq.'PARDISO') then
isolver=7
else
write(*,*) '*WARNING reading *CFD: unknown solver;'
write(*,*) ' the default solver is used'
endif
!
call getnewline(inpc,textpart,istat,n,key,iline,ipol,inl,
& ipoinp,inp,ipoinpc)
if((istat.lt.0).or.(key.eq.1)) then
if(nmethod.eq.4) then
!
! transient cfd
!
write(*,*) '*WARNING reading *CFD: a nonlinear geometric
& analysis is requested'
write(*,*) ' but no time increment nor step is speci
&fied'
write(*,*) ' the defaults (1,1) are used'
tinc=1.d0
tper=1.d0
tmin=1.d-5
tmax=1.d+30
tincf=-1.d0
elseif(nmethod.eq.1) then
!
! steady state cfd: initial increment time and step
! time are set to large values since only the first
! increment is calculated
!
tinc=1.d+30
tper=1.d+30
tmin=1.d-5
tmax=1.d+30
tincf=-1.d0
endif
if(timereset)ttime=ttime-tper
return
endif
!
read(textpart(1)(1:20),'(f20.0)',iostat=istat) tinc
if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
&"*CFD%")
read(textpart(2)(1:20),'(f20.0)',iostat=istat) tper
if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
&"*CFD%")
read(textpart(3)(1:20),'(f20.0)',iostat=istat) tmin
if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
&"*CFD%")
read(textpart(4)(1:20),'(f20.0)',iostat=istat) tmax
if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
&"*CFD%")
read(textpart(5)(1:20),'(f20.0)',iostat=istat) tincf
if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
&"*CFD%")
!
if(tinc.le.0.d0) then
write(*,*) '*ERROR reading *CFD: initial increment size is
&negative'
endif
if(tper.le.0.d0) then
write(*,*) '*ERROR reading *CFD: step size is negative'
endif
if(tinc.gt.tper) then
write(*,*) '*ERROR reading *CFD: initial increment size exc
&eeds step size'
endif
!
if(idrct.ne.1) then
if(dabs(tmin).lt.1.d-6*tper) then
tmin=min(tinc,1.d-6*tper)
endif
if(dabs(tmax).lt.1.d-10) then
tmax=1.d+30
endif
endif
!
if(timereset)ttime=ttime-tper
!
call getnewline(inpc,textpart,istat,n,key,iline,ipol,inl,
& ipoinp,inp,ipoinpc)
!
return
end
|