1!
2!     CalculiX - A 3-dimensional finite element program
3!              Copyright (C) 1998-2021 Guido Dhondt
4!
5!     This program is free software; you can redistribute it and/or
6!     modify it under the terms of the GNU General Public License as
7!     published by the Free Software Foundation(version 2);
8!
9!
10!     This program is distributed in the hope that it will be useful,
11!     but WITHOUT ANY WARRANTY; without even the implied warranty of
12!     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13!     GNU General Public License for more details.
14!
15!     You should have received a copy of the GNU General Public License
16!     along with this program; if not, write to the Free Software
17!     Foundation, Inc., 675 Mass Ave, Cambridge, MA 02139, USA.
18!
19      subroutine statics(inpc,textpart,nmethod,iperturb,isolver,istep,
20     &  istat,n,tinc,tper,tmin,tmax,idrct,iline,ipol,inl,ipoinp,inp,
21     &  ithermal,cs,ics,tieset,istartset,
22     &  iendset,ialset,ipompc,nodempc,coefmpc,nmpc,nmpc_,ikmpc,
23     &  ilmpc,mpcfree,mcs,set,nset,labmpc,ipoinpc,iexpl,nef,ttime,
24     &  iaxial,nelcon,nmat,tincf,ier)
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!             8: pastix
35!
36!      iexpl==0:  structure:implicit, fluid:incompressible
37!
38      implicit none
39!
40      logical timereset
41!
42      character*1 inpc(*)
43      character*20 labmpc(*),solver
44      character*81 set(*),tieset(3,*)
45      character*132 textpart(16)
46!
47      integer nmethod,iperturb(*),isolver,istep,istat,n,key,i,idrct,
48     &  iline,ipol,inl,ipoinp(2,*),inp(3,*),ithermal(*),ics(*),iexpl,
49     &  istartset(*),iendset(*),ialset(*),ipompc(*),nodempc(3,*),
50     &  nmpc,nmpc_,ikmpc(*),ilmpc(*),mpcfree,nset,mcs,ipoinpc(0:*),
51     &  nef,iaxial,nelcon(2,*),nmat,ier
52!
53      real*8 tinc,tper,tmin,tmax,cs(17,*),coefmpc(*),ttime,tincf
54!
55      idrct=0
56      tmin=0.d0
57      tmax=0.d0
58      timereset=.false.
59!
60c      if((iperturb.eq.1).and.(istep.ge.1)) then
61c         write(*,*) '*ERROR reading *STATIC: perturbation analysis is'
62c         write(*,*) '       not provided in a *STATIC step. Perform'
63c         write(*,*) '       a genuine nonlinear geometric calculation'
64c         write(*,*) '       instead (parameter NLGEOM)'
65c         ier=1
66c         return
67c      endif
68!
69      if(istep.lt.1) then
70         write(*,*) '*ERROR reading *STATIC: *STATIC can only be used'
71         write(*,*) '       within a STEP'
72         ier=1
73         return
74      endif
75c!
76c!     no creep allowed in a *STATIC step
77c!
78c      do i=1,nmat
79c         if(nelcon(1,i).eq.-52) nelcon(1,i)=-51
80c      enddo
81!
82!     no heat transfer analysis
83!
84      if(ithermal(1).gt.1) then
85         ithermal(1)=1
86      endif
87!
88!     default solver
89!
90      solver='                    '
91      if(isolver.eq.0) then
92         solver(1:7)='SPOOLES'
93      elseif(isolver.eq.2) then
94         solver(1:16)='ITERATIVESCALING'
95      elseif(isolver.eq.3) then
96         solver(1:17)='ITERATIVECHOLESKY'
97      elseif(isolver.eq.4) then
98         solver(1:3)='SGI'
99      elseif(isolver.eq.5) then
100         solver(1:5)='TAUCS'
101      elseif(isolver.eq.7) then
102         solver(1:7)='PARDISO'
103      elseif(isolver.eq.8) then
104         solver(1:6)='PASTIX'
105      endif
106!
107      do i=2,n
108         if(textpart(i)(1:7).eq.'SOLVER=') then
109            read(textpart(i)(8:27),'(a20)') solver
110         elseif((textpart(i)(1:6).eq.'DIRECT').and.
111     &          (textpart(i)(1:9).ne.'DIRECT=NO')) then
112            idrct=1
113         elseif(textpart(i)(1:9).eq.'TIMERESET') then
114            timereset=.true.
115         elseif(textpart(i)(1:17).eq.'TOTALTIMEATSTART=') then
116            read(textpart(i)(18:37),'(f20.0)',iostat=istat) ttime
117         else
118            write(*,*)
119     &        '*WARNING reading *STATIC: parameter not recognized:'
120            write(*,*) '         ',
121     &                 textpart(i)(1:index(textpart(i),' ')-1)
122            call inputwarning(inpc,ipoinpc,iline,
123     &"*STATIC%")
124         endif
125      enddo
126!
127      if(solver(1:7).eq.'SPOOLES') then
128         isolver=0
129      elseif(solver(1:16).eq.'ITERATIVESCALING') then
130         isolver=2
131      elseif(solver(1:17).eq.'ITERATIVECHOLESKY') then
132         isolver=3
133      elseif(solver(1:3).eq.'SGI') then
134         isolver=4
135      elseif(solver(1:5).eq.'TAUCS') then
136         isolver=5
137      elseif(solver(1:7).eq.'PARDISO') then
138         isolver=7
139      elseif(solver(1:6).eq.'PASTIX') then
140         isolver=8
141      else
142         write(*,*) '*WARNING reading *STATIC: unknown solver;'
143         write(*,*) '         the default solver is used'
144      endif
145!
146      nmethod=1
147!
148!     check for nodes on a cyclic symmetry axis
149!
150      if((mcs.eq.0).or.(iaxial.eq.180)) then
151         call getnewline(inpc,textpart,istat,n,key,iline,ipol,inl,
152     &        ipoinp,inp,ipoinpc)
153      else
154         n=3
155         textpart(2)='NMIN=0
156     &
157     &                '
158         textpart(3)='NMAX=0
159     &
160     &                '
161         nmethod=2
162         call selectcyclicsymmetrymodess(inpc,textpart,cs,ics,tieset,
163     &        istartset,
164     &        iendset,ialset,ipompc,nodempc,coefmpc,nmpc,nmpc_,ikmpc,
165     &        ilmpc,mpcfree,mcs,set,nset,labmpc,istep,istat,n,iline,
166     &        ipol,inl,ipoinp,inp,nmethod,key,ipoinpc)
167         nmethod=1
168         do i=1,mcs
169            cs(2,i)=-0.5d0
170            cs(3,i)=-0.5d0
171         enddo
172      endif
173!
174      if((istat.lt.0).or.(key.eq.1)) then
175         if((iperturb(1).ge.2).or.(nef.gt.0)) then
176            write(*,*) '*WARNING reading *STATIC: a nonlinear analysis i
177     &s requested'
178            write(*,*) '         but no time increment nor step is speci
179     &fied'
180            write(*,*) '         the defaults (1,1) are used'
181            write(*,*)
182            tinc=1.d0
183            tper=1.d0
184            tmin=1.d-5
185            tmax=1.d+30
186            tincf=-1.d0
187c            tincf=1.d-2
188         else
189            tper=1.d0
190         endif
191         if(timereset)ttime=ttime-tper
192         return
193      endif
194!
195      read(textpart(1)(1:20),'(f20.0)',iostat=istat) tinc
196      if(istat.gt.0) then
197         call inputerror(inpc,ipoinpc,iline,
198     &        "*STATIC%",ier)
199         return
200      endif
201      read(textpart(2)(1:20),'(f20.0)',iostat=istat) tper
202      if(istat.gt.0) then
203         call inputerror(inpc,ipoinpc,iline,
204     &        "*STATIC%",ier)
205         return
206      endif
207      read(textpart(3)(1:20),'(f20.0)',iostat=istat) tmin
208      if(istat.gt.0) then
209         call inputerror(inpc,ipoinpc,iline,
210     &        "*STATIC%",ier)
211         return
212      endif
213      read(textpart(4)(1:20),'(f20.0)',iostat=istat) tmax
214      if(istat.gt.0) then
215         call inputerror(inpc,ipoinpc,iline,
216     &        "*STATIC%",ier)
217         return
218      endif
219      read(textpart(5)(1:20),'(f20.0)',iostat=istat) tincf
220      if(istat.gt.0) then
221         call inputerror(inpc,ipoinpc,iline,
222     &        "*STATIC%",ier)
223         return
224      endif
225!
226      if(tper.lt.0.d0) then
227         write(*,*) '*ERROR reading *STATIC: step size is negative'
228         ier=1
229         return
230      elseif(tper.le.0.d0) then
231         tper=1.d0
232      endif
233      if(tinc.lt.0.d0) then
234         write(*,*) '*ERROR reading *STATIC: initial increment size is n
235     &egative'
236         ier=1
237         return
238      elseif(tinc.le.0.d0) then
239         tinc=tper
240      endif
241      if(tinc.gt.tper) then
242         write(*,*) '*ERROR reading *STATIC: initial increment size exce
243     &eds step size'
244         ier=1
245         return
246      endif
247!
248      if(idrct.ne.1) then
249         if(dabs(tmin).lt.1.d-6*tper) then
250            write(*,*) '*WARNING reading *STATIC:'
251            write(*,*) '         the minimum increment ',tmin
252            write(*,*) '         is smaller then 1.e-6 times the '
253            write(*,*) '         step time;'
254            write(*,*) '         the minimum increment is changed'
255            write(*,*) '         to ',min(tinc,1.d-6*tper)
256            write(*,*) '         which is the minimum of the initial'
257            write(*,*)
258     &         '         increment time and 1.e-6 times the step time'
259            tmin=min(tinc,1.d-6*tper)
260         endif
261         if(dabs(tmax).lt.1.d-10) then
262            tmax=1.d+30
263         endif
264         if(tinc.gt.dabs(tmax)) then
265            write(*,*) '*WARNING reading *STATIC:'
266            write(*,*) '         the initial increment ',tinc
267            write(*,*) '         exceeds the maximum increment ',
268     &          tmax
269            write(*,*) '         the initial increment is reduced'
270            write(*,*) '         to the maximum value'
271            tinc=dabs(tmax)
272         endif
273      endif
274!
275      if(timereset)ttime=ttime-tper
276!
277      call getnewline(inpc,textpart,istat,n,key,iline,ipol,inl,
278     &     ipoinp,inp,ipoinpc)
279!
280      return
281      end
282
283