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

Go to the source code of this file.

Functions/Subroutines

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

Function/Subroutine Documentation

◆ dynamics()

subroutine dynamics ( character*1, dimension(*)  inpc,
character*132, dimension(16)  textpart,
integer  nmethod,
integer  iperturb,
real*8  tinc,
real*8  tper,
real*8  tmin,
real*8  tmax,
integer  idrct,
real*8  alpha,
integer  iexpl,
integer  isolver,
integer  istep,
integer  istat,
integer  n,
integer  iline,
integer  ipol,
integer  inl,
integer, dimension(2,*)  ipoinp,
integer, dimension(3,*)  inp,
integer  ithermal,
integer, dimension(0:*)  ipoinpc,
integer  cfd,
real*8, dimension(*)  ctrl,
real*8  tincf,
integer  nener 
)
22 !
23 ! reading the input deck: *DYNAMIC
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==2: structure:explicit, fluid:incompressible
34 !
35  implicit none
36 !
37  character*1 inpc(*)
38  character*20 solver
39  character*132 textpart(16)
40 !
41  integer nmethod,istep,istat,n,key,i,iperturb,idrct,iexpl,
42  & isolver,iline,ipol,inl,ipoinp(2,*),inp(3,*),ithermal,
43  & ipoinpc(0:*),cfd,nener
44 !
45  real*8 tinc,tper,tmin,tmax,alpha,ctrl(*),tincf
46 !
47  if(istep.lt.1) then
48  write(*,*) '*ERROR reading *DYNAMIC: *DYNAMIC can only'
49  write(*,*) ' be used within a STEP'
50  call exit(201)
51  endif
52 !
53 ! default is implicit
54 !
55  iexpl=0
56 !
57 ! no heat transfer analysis
58 !
59  if(ithermal.gt.1) then
60  ithermal=1
61  endif
62 !
63 ! only nonlinear analysis allowed for this procedure
64 !
65  if(iperturb.lt.2) iperturb=2
66 !
67 ! default values
68 !
69  idrct=0
70  alpha=-0.05d0
71  tmin=0.d0
72  tmax=0.d0
73  tincf=1.d-2
74 !
75 ! default solver
76 !
77  solver=' '
78  if(isolver.eq.0) then
79  solver(1:7)='SPOOLES'
80  elseif(isolver.eq.2) then
81  solver(1:16)='ITERATIVESCALING'
82  elseif(isolver.eq.3) then
83  solver(1:17)='ITERATIVECHOLESKY'
84  elseif(isolver.eq.4) then
85  solver(1:3)='SGI'
86  elseif(isolver.eq.5) then
87  solver(1:5)='TAUCS'
88  elseif(isolver.eq.7) then
89  solver(1:7)='PARDISO'
90  endif
91 !
92  do i=2,n
93  if(textpart(i)(1:6).eq.'ALPHA=') then
94  read(textpart(i)(7:26),'(f20.0)',iostat=istat) alpha
95  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
96  &"*DYNAMIC%")
97  if(alpha.lt.-1.d0/3.d0) then
98  write(*,*) '*WARNING reading *DYNAMIC: alpha is smaller'
99  write(*,*) ' than -1/3 and is reset to -1/3'
100  alpha=-1.d0/3.d0
101  elseif(alpha.gt.0.d0) then
102  write(*,*) '*WARNING reading *DYNAMIC: alpha is greater'
103  write(*,*) ' than 0 and is reset to 0'
104  alpha=0.d0
105  endif
106  elseif(textpart(i)(1:8).eq.'EXPLICIT') then
107  iexpl=2
108  elseif((textpart(i)(1:6).eq.'DIRECT').and.
109  & (textpart(i)(1:9).ne.'DIRECT=NO')) then
110  idrct=1
111  elseif(textpart(i)(1:7).eq.'SOLVER=') then
112  read(textpart(i)(8:27),'(a20)') solver
113  else
114  write(*,*)
115  & '*WARNING reading *DYNAMIC: parameter not recognized:'
116  write(*,*) ' ',
117  & textpart(i)(1:index(textpart(i),' ')-1)
118  call inputwarning(inpc,ipoinpc,iline,
119  &"*DYNAMIC%")
120  endif
121  enddo
122 !
123  if(solver(1:7).eq.'SPOOLES') then
124  isolver=0
125  elseif(solver(1:16).eq.'ITERATIVESCALING') then
126  isolver=2
127  elseif(solver(1:17).eq.'ITERATIVECHOLESKY') then
128  isolver=3
129  elseif(solver(1:3).eq.'SGI') then
130  isolver=4
131  elseif(solver(1:5).eq.'TAUCS') then
132  isolver=5
133  elseif(solver(1:7).eq.'PARDISO') then
134  isolver=7
135  else
136  write(*,*) '*WARNING reading *DYNAMIC: unknown solver;'
137  write(*,*) ' the default solver is used'
138  endif
139 !
140  call getnewline(inpc,textpart,istat,n,key,iline,ipol,inl,
141  & ipoinp,inp,ipoinpc)
142  if((istat.lt.0).or.(key.eq.1)) then
143  if((iperturb.ge.2).or.(cfd.eq.1)) then
144  write(*,*)'*WARNING reading *DYNAMIC: a nonlinear analysis i
145  &s requested'
146  write(*,*) ' but no time increment nor step is speci
147  &fied'
148  write(*,*) ' the defaults (1,1) are used'
149  tinc=1.d0
150  tper=1.d0
151  tmin=1.d-5
152  tmax=1.d+30
153  tincf=1.d-2
154  endif
155  nmethod=4
156  return
157  endif
158 !
159  read(textpart(1)(1:20),'(f20.0)',iostat=istat) tinc
160  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
161  &"*DYNAMIC%")
162  read(textpart(2)(1:20),'(f20.0)',iostat=istat) tper
163  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
164  &"*DYNAMIC%")
165  read(textpart(3)(1:20),'(f20.0)',iostat=istat) tmin
166  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
167  &"*DYNAMIC%")
168  read(textpart(4)(1:20),'(f20.0)',iostat=istat) tmax
169  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
170  &"*DYNAMIC%")
171  read(textpart(4)(1:20),'(f20.0)',iostat=istat) tincf
172  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
173  &"*DYNAMIC%")
174 !
175  if(tinc.le.0.d0) then
176  write(*,*)'*ERROR reading *DYNAMIC: initial increment size is n
177  &egative'
178  endif
179  if(tper.le.0.d0) then
180  write(*,*) '*ERROR reading *DYNAMIC: step size is negative'
181  endif
182  if(tinc.gt.tper) then
183  write(*,*)'*ERROR reading *DYNAMIC: initial increment size exce
184  &eds step size'
185  endif
186  if((cfd.eq.1).and.(tincf.le.0.d0)) then
187  write(*,*) '*WARNING reading *DYNAMIC: initial CFD increment si
188  &ze is zero or negative; the default of 0.01 is taken'
189  tincf=1.d-2
190  endif
191 !
192  if(idrct.ne.1) then
193  if(dabs(tmin).lt.1.d-10*tper) then
194  tmin=min(tinc,1.e-10*tper)
195  endif
196  if(dabs(tmax).lt.1.d-10) then
197  tmax=1.d+30
198  endif
199  if(tinc.gt.dabs(tmax)) then
200  write(*,*) '*WARNING reading *DYNAMIC:'
201  write(*,*) ' the initial increment ',tinc
202  write(*,*) ' exceeds the maximum increment ',
203  & tmax
204  write(*,*) ' the initial increment is reduced'
205  write(*,*) ' to the maximum value'
206  tinc=dabs(tmax)
207  endif
208  endif
209 !
210 ! 10 cutbacks allowed for dynamics (because of contact)
211 !
212  ctrl(8)=10.5d0
213 !
214  nmethod=4
215  nener=1
216 !
217  call getnewline(inpc,textpart,istat,n,key,iline,ipol,inl,
218  & ipoinp,inp,ipoinpc)
219 !
220  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)