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

Go to the source code of this file.

Functions/Subroutines

subroutine statics (inpc, textpart, nmethod, iperturb, isolver, istep, istat, n, tinc, tper, tmin, tmax, idrct, iline, ipol, inl, ipoinp, inp, ithermal, cs, ics, tieset, istartset, iendset, ialset, ipompc, nodempc, coefmpc, nmpc, nmpc_, ikmpc, ilmpc, mpcfree, mcs, set, nset, labmpc, ipoinpc, iexpl, cfd, ttime, iaxial, nelcon, nmat, tincf)
 

Function/Subroutine Documentation

◆ statics()

subroutine statics ( 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  iline,
integer  ipol,
integer  inl,
integer, dimension(2,*)  ipoinp,
integer, dimension(3,*)  inp,
integer  ithermal,
real*8, dimension(17,*)  cs,
integer, dimension(*)  ics,
character*81, dimension(3,*)  tieset,
integer, dimension(*)  istartset,
integer, dimension(*)  iendset,
integer, dimension(*)  ialset,
integer, dimension(*)  ipompc,
integer, dimension(3,*)  nodempc,
real*8, dimension(*)  coefmpc,
integer  nmpc,
integer  nmpc_,
integer, dimension(*)  ikmpc,
integer, dimension(*)  ilmpc,
integer  mpcfree,
integer  mcs,
character*81, dimension(*)  set,
integer  nset,
character*20, dimension(*)  labmpc,
integer, dimension(0:*)  ipoinpc,
integer  iexpl,
integer  cfd,
real*8  ttime,
integer  iaxial,
integer, dimension(2,*)  nelcon,
integer  nmat,
real*8  tincf 
)
25 !
26 ! reading the input deck: *STATIC
27 !
28 ! isolver=0: SPOOLES
29 ! 2: iterative solver with diagonal scaling
30 ! 3: iterative solver with Cholesky preconditioning
31 ! 4: sgi solver
32 ! 5: TAUCS
33 ! 7: pardiso
34 !
35 ! iexpl==0: structure:implicit, fluid:incompressible
36 !
37  implicit none
38 !
39  logical timereset
40 !
41  character*1 inpc(*)
42  character*20 labmpc(*),solver
43  character*81 set(*),tieset(3,*)
44  character*132 textpart(16)
45 !
46  integer nmethod,iperturb,isolver,istep,istat,n,key,i,idrct,
47  & iline,ipol,inl,ipoinp(2,*),inp(3,*),ithermal,ics(*),iexpl,
48  & istartset(*),iendset(*),ialset(*),ipompc(*),nodempc(3,*),
49  & nmpc,nmpc_,ikmpc(*),ilmpc(*),mpcfree,nset,mcs,ipoinpc(0:*),
50  & cfd,iaxial,nelcon(2,*),nmat
51 !
52  real*8 tinc,tper,tmin,tmax,cs(17,*),coefmpc(*),ttime,tincf
53 !
54  idrct=0
55  tmin=0.d0
56  tmax=0.d0
57  timereset=.false.
58 !
59  if((iperturb.eq.1).and.(istep.ge.1)) then
60  write(*,*) '*ERROR reading *STATIC: perturbation analysis is'
61  write(*,*) ' not provided in a *STATIC step. Perform'
62  write(*,*) ' a genuine nonlinear geometric calculation'
63  write(*,*) ' instead (parameter NLGEOM)'
64  call exit(201)
65  endif
66 !
67  if(istep.lt.1) then
68  write(*,*) '*ERROR reading *STATIC: *STATIC can only be used'
69  write(*,*) ' within a STEP'
70  call exit(201)
71  endif
72 c!
73 c! no creep allowed in a *STATIC step
74 c!
75 c do i=1,nmat
76 c if(nelcon(1,i).eq.-52) nelcon(1,i)=-51
77 c enddo
78 !
79 ! no heat transfer analysis
80 !
81  if(ithermal.gt.1) then
82  ithermal=1
83  endif
84 !
85 ! default solver
86 !
87  solver=' '
88  if(isolver.eq.0) then
89  solver(1:7)='SPOOLES'
90  elseif(isolver.eq.2) then
91  solver(1:16)='ITERATIVESCALING'
92  elseif(isolver.eq.3) then
93  solver(1:17)='ITERATIVECHOLESKY'
94  elseif(isolver.eq.4) then
95  solver(1:3)='SGI'
96  elseif(isolver.eq.5) then
97  solver(1:5)='TAUCS'
98  elseif(isolver.eq.7) then
99  solver(1:7)='PARDISO'
100  endif
101 !
102  do i=2,n
103  if(textpart(i)(1:7).eq.'SOLVER=') then
104  read(textpart(i)(8:27),'(a20)') solver
105  elseif((textpart(i)(1:6).eq.'DIRECT').and.
106  & (textpart(i)(1:9).ne.'DIRECT=NO')) then
107  idrct=1
108  elseif(textpart(i)(1:9).eq.'TIMERESET') then
109  timereset=.true.
110  elseif(textpart(i)(1:17).eq.'TOTALTIMEATSTART=') then
111  read(textpart(i)(18:37),'(f20.0)',iostat=istat) ttime
112  else
113  write(*,*)
114  & '*WARNING reading *STATIC: parameter not recognized:'
115  write(*,*) ' ',
116  & textpart(i)(1:index(textpart(i),' ')-1)
117  call inputwarning(inpc,ipoinpc,iline,
118  &"*STATIC%")
119  endif
120  enddo
121 !
122  if(solver(1:7).eq.'SPOOLES') then
123  isolver=0
124  elseif(solver(1:16).eq.'ITERATIVESCALING') then
125  isolver=2
126  elseif(solver(1:17).eq.'ITERATIVECHOLESKY') then
127  isolver=3
128  elseif(solver(1:3).eq.'SGI') then
129  isolver=4
130  elseif(solver(1:5).eq.'TAUCS') then
131  isolver=5
132  elseif(solver(1:7).eq.'PARDISO') then
133  isolver=7
134  else
135  write(*,*) '*WARNING reading *STATIC: unknown solver;'
136  write(*,*) ' the default solver is used'
137  endif
138 !
139  nmethod=1
140 !
141 ! check for nodes on a cyclic symmetry axis
142 !
143  if((mcs.eq.0).or.(iaxial.eq.180)) then
144  call getnewline(inpc,textpart,istat,n,key,iline,ipol,inl,
145  & ipoinp,inp,ipoinpc)
146  else
147  n=3
148  textpart(2)='NMIN=0
149  &
150  & '
151  textpart(3)='NMAX=0
152  &
153  & '
154  nmethod=2
155  call selectcyclicsymmetrymodess(inpc,textpart,cs,ics,tieset,
156  & istartset,
157  & iendset,ialset,ipompc,nodempc,coefmpc,nmpc,nmpc_,ikmpc,
158  & ilmpc,mpcfree,mcs,set,nset,labmpc,istep,istat,n,iline,
159  & ipol,inl,ipoinp,inp,nmethod,key,ipoinpc)
160  nmethod=1
161  do i=1,mcs
162  cs(2,i)=-0.5
163  cs(3,i)=-0.5
164  enddo
165  endif
166 !
167  if((istat.lt.0).or.(key.eq.1)) then
168  if((iperturb.ge.2).or.(cfd.eq.1)) then
169  write(*,*) '*WARNING reading *STATIC: a nonlinear analysis i
170  &s requested'
171  write(*,*) ' but no time increment nor step is speci
172  &fied'
173  write(*,*) ' the defaults (1,1) are used'
174  write(*,*)
175  tinc=1.d0
176  tper=1.d0
177  tmin=1.d-5
178  tmax=1.d+30
179  tincf=-1.d0
180 c tincf=1.d-2
181  else
182  tper=1.d0
183  endif
184  if(timereset)ttime=ttime-tper
185  return
186  endif
187 !
188  read(textpart(1)(1:20),'(f20.0)',iostat=istat) tinc
189  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
190  &"*STATIC%")
191  read(textpart(2)(1:20),'(f20.0)',iostat=istat) tper
192  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
193  &"*STATIC%")
194  read(textpart(3)(1:20),'(f20.0)',iostat=istat) tmin
195  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
196  &"*STATIC%")
197  read(textpart(4)(1:20),'(f20.0)',iostat=istat) tmax
198  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
199  &"*STATIC%")
200  read(textpart(5)(1:20),'(f20.0)',iostat=istat) tincf
201  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
202  &"*STATIC%")
203 !
204  if(tper.lt.0.d0) then
205  write(*,*) '*ERROR reading *STATIC: step size is negative'
206  call exit(201)
207  elseif(tper.le.0.d0) then
208  tper=1.d0
209  endif
210  if(tinc.lt.0.d0) then
211  write(*,*) '*ERROR reading *STATIC: initial increment size is n
212  &egative'
213  call exit(201)
214  elseif(tinc.le.0.d0) then
215  tinc=tper
216  endif
217  if(tinc.gt.tper) then
218  write(*,*) '*ERROR reading *STATIC: initial increment size exce
219  &eds step size'
220  call exit(201)
221  endif
222 !
223  if(idrct.ne.1) then
224  if(dabs(tmin).lt.1.d-6*tper) then
225  tmin=min(tinc,1.d-6*tper)
226  endif
227  if(dabs(tmax).lt.1.d-10) then
228  tmax=1.d+30
229  endif
230  if(tinc.gt.dabs(tmax)) then
231  write(*,*) '*WARNING reading *STATIC:'
232  write(*,*) ' the initial increment ',tinc
233  write(*,*) ' exceeds the maximum increment ',
234  & tmax
235  write(*,*) ' the initial increment is reduced'
236  write(*,*) ' to the maximum value'
237  tinc=dabs(tmax)
238  endif
239  endif
240 !
241  if(timereset)ttime=ttime-tper
242 !
243  call getnewline(inpc,textpart,istat,n,key,iline,ipol,inl,
244  & ipoinp,inp,ipoinpc)
245 !
246  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 selectcyclicsymmetrymodess(inpc, textpart, cs, ics, tieset, istartset, iendset, ialset, ipompc, nodempc, coefmpc, nmpc, nmpc_, ikmpc, ilmpc, mpcfree, mcs, set, nset, labmpc, istep, istat, n, iline, ipol, inl, ipoinp, inp, nmethod, key, ipoinpc)
Definition: selectcyclicsymmetrymodess.f:24
subroutine inputerror(inpc, ipoinpc, iline, text)
Definition: inputerror.f:20
Hosted by OpenAircraft.com, (Michigan UAV, LLC)