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

Go to the source code of this file.

Functions/Subroutines

subroutine viscos (inpc, textpart, nmethod, iperturb, isolver, istep, istat, n, tinc, tper, tmin, tmax, idrct, iline, ipol, inl, ipoinp, inp, ipoinpc, nelcon, nmat, ncmat_, ctrl)
 

Function/Subroutine Documentation

◆ viscos()

subroutine viscos ( 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, dimension(0:*)  ipoinpc,
integer, dimension(2,*)  nelcon,
integer  nmat,
integer  ncmat_,
real*8, dimension(*)  ctrl 
)
22 !
23 ! reading the input deck: *VISCO (provided for compatibility
24 ! reasons with ABAQUS)
25 !
26 ! isolver=0: SPOOLES
27 ! 2: iterative solver with diagonal scaling
28 ! 3: iterative solver with Cholesky preconditioning
29 ! 4: sgi solver
30 ! 5: TAUCS
31 ! 7: pardiso
32 !
33  implicit none
34 !
35  character*1 inpc(*)
36  character*20 solver
37  character*132 textpart(16)
38 !
39  integer nmethod,iperturb,isolver,istep,istat,n,key,i,idrct,
40  & iline,ipol,inl,ipoinp(2,*),inp(3,*),ipoinpc(0:*),nelcon(2,*),
41  & nmat,ncmat_
42 !
43  real*8 tinc,tper,tmin,tmax,ctrl(*)
44 !
45  idrct=0
46  tmin=0.d0
47  tmax=0.d0
48 !
49  if((iperturb.eq.1).and.(istep.gt.1)) then
50  write(*,*) '*ERROR reading *VISCO: perturbation analysis is'
51  write(*,*) ' not provided in a *VISCO step. Perform'
52  write(*,*) ' a genuine nonlinear geometric calculation'
53  write(*,*) ' instead (parameter NLGEOM)'
54  call exit(201)
55  endif
56 !
57  if(istep.lt.1) then
58  write(*,*) '*ERROR reading *VISCO: *VISCO can only be used'
59  write(*,*) ' within a STEP'
60  call exit(201)
61  endif
62 !
63 ! default solver
64 !
65  solver=' '
66  if(isolver.eq.0) then
67  solver(1:7)='SPOOLES'
68  elseif(isolver.eq.2) then
69  solver(1:16)='ITERATIVESCALING'
70  elseif(isolver.eq.3) then
71  solver(1:17)='ITERATIVECHOLESKY'
72  elseif(isolver.eq.4) then
73  solver(1:3)='SGI'
74  elseif(isolver.eq.5) then
75  solver(1:5)='TAUCS'
76  elseif(isolver.eq.7) then
77  solver(1:7)='PARDISO'
78  endif
79 !
80  do i=2,n
81  if(textpart(i)(1:7).eq.'SOLVER=') then
82  read(textpart(i)(8:27),'(a20)') solver
83  elseif((textpart(i)(1:6).eq.'DIRECT').and.
84  & (textpart(i)(1:9).ne.'DIRECT=NO')) then
85  idrct=1
86  elseif(textpart(i)(1:6).eq.'CETOL=') then
87  read(textpart(i)(7:26),'(f20.0)',iostat=istat) ctrl(40)
88  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
89  &"*VISCO%")
90  else
91  write(*,*)
92  & '*WARNING reading *VISCO: parameter not recognized:'
93  write(*,*) ' ',
94  & textpart(i)(1:index(textpart(i),' ')-1)
95  call inputwarning(inpc,ipoinpc,iline,
96  &"*VISCO%")
97  endif
98  enddo
99 !
100  if(ctrl(40).le.0.d0) then
101  write(*,*) '*ERROR reading *VISCO:'
102  write(*,*) ' no strictly positive value'
103  write(*,*) ' for the parameter CETOL given'
104  call exit(201)
105  endif
106 !
107  if(solver(1:7).eq.'SPOOLES') then
108  isolver=0
109  elseif(solver(1:16).eq.'ITERATIVESCALING') then
110  isolver=2
111  elseif(solver(1:17).eq.'ITERATIVECHOLESKY') then
112  isolver=3
113  elseif(solver(1:3).eq.'SGI') then
114  isolver=4
115  elseif(solver(1:5).eq.'TAUCS') then
116  isolver=5
117  elseif(solver(1:7).eq.'PARDISO') then
118  isolver=7
119  else
120  write(*,*) '*WARNING reading *VISCO: unknown solver;'
121  write(*,*) ' the default solver is used'
122  endif
123 !
124  nmethod=-1
125 !
126  call getnewline(inpc,textpart,istat,n,key,iline,ipol,inl,
127  & ipoinp,inp,ipoinpc)
128  if((istat.lt.0).or.(key.eq.1)) then
129  if(iperturb.ge.2) then
130  write(*,*) '*WARNING reading *VISCO: a nonlinear analysis is
131  & requested'
132  write(*,*) ' but no time increment nor step is speci
133  &fied'
134  write(*,*) ' the defaults (1,1) are used'
135  tinc=1.d0
136  tper=1.d0
137  tmin=1.d-5
138  tmax=1.d+30
139  endif
140  return
141  endif
142 !
143  read(textpart(1)(1:20),'(f20.0)',iostat=istat) tinc
144  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
145  &"*VISCO%")
146  read(textpart(2)(1:20),'(f20.0)',iostat=istat) tper
147  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
148  &"*VISCO%")
149  read(textpart(3)(1:20),'(f20.0)',iostat=istat) tmin
150  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
151  &"*VISCO%")
152  read(textpart(4)(1:20),'(f20.0)',iostat=istat) tmax
153  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
154  &"*VISCO%")
155 !
156  if(tinc.le.0.d0) then
157  write(*,*) '*ERROR reading *VISCO: initial increment size is ne
158  &gative'
159  endif
160  if(tper.le.0.d0) then
161  write(*,*) '*ERROR reading *VISCO: step size is negative'
162  endif
163  if(tinc.gt.tper) then
164  write(*,*) '*ERROR reading *VISCO: initial increment size excee
165  &ds step size'
166  endif
167 !
168  if(idrct.ne.1) then
169  if(dabs(tmin).lt.1.d-6*tper) then
170  tmin=min(tinc,1.d-6*tper)
171  endif
172  if(dabs(tmax).lt.1.d-10) then
173  tmax=1.d+30
174  if(tinc.gt.dabs(tmax)) then
175  write(*,*) '*WARNING reading *VISCO:'
176  write(*,*) ' the initial increment ',tinc
177  write(*,*) ' exceeds the maximum increment ',
178  & tmax
179  write(*,*) ' the initial increment is reduced'
180  write(*,*) ' to the maximum value'
181  tinc=dabs(tmax)
182  endif
183  endif
184  endif
185 !
186  call getnewline(inpc,textpart,istat,n,key,iline,ipol,inl,
187  & ipoinp,inp,ipoinpc)
188 !
189  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)