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

Go to the source code of this file.

Functions/Subroutines

subroutine electromagneticss (inpc, textpart, nmethod, iperturb, isolver, istep, istat, n, tinc, tper, tmin, tmax, idrct, ithermal, iline, ipol, inl, ipoinp, inp, alpha, mei, fei, ipoinpc, ctrl, ttime)
 

Function/Subroutine Documentation

◆ electromagneticss()

subroutine electromagneticss ( 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  ithermal,
integer  iline,
integer  ipol,
integer  inl,
integer, dimension(2,*)  ipoinp,
integer, dimension(3,*)  inp,
real*8  alpha,
integer, dimension(4)  mei,
real*8, dimension(3)  fei,
integer, dimension(0:*)  ipoinpc,
real*8, dimension(*)  ctrl,
real*8  ttime 
)
22 !
23 ! reading the input deck: *ELECTROMAGNETICS
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  implicit none
33 !
34  logical timereset
35 !
36  character*1 inpc(*)
37  character*20 solver
38  character*132 textpart(16)
39 !
40  integer nmethod,iperturb,isolver,istep,istat,n,key,i,idrct,nev,
41  & ithermal,iline,ipol,inl,ipoinp(2,*),inp(3,*),mei(4),ncv,mxiter,
42  & ipoinpc(0:*),idirect,iheat
43 !
44  real*8 tinc,tper,tmin,tmax,alpha,fei(3),tol,fmin,fmax,ctrl(*),
45  & ttime
46 !
47  tmin=0.d0
48  tmax=0.d0
49  nmethod=9
50  alpha=0.d0
51  mei(4)=0
52  timereset=.false.
53 !
54 ! defaults for fmin and fmax
55 !
56  fmin=-1.d0
57  fmax=-1.d0
58 !
59  if(iperturb.eq.0) then
60  iperturb=2
61  elseif((iperturb.eq.1).and.(istep.gt.1)) then
62  write(*,*)
63  & '*ERROR reading *ELECTROMAGNETICS: perturbation analysis is'
64  write(*,*) ' not provided in a *HEAT TRANSFER step.'
65  call exit(201)
66  endif
67 !
68  if(istep.lt.1) then
69  write(*,*)
70  & '*ERROR reading *ELECTROMAGNETICS: *HEAT TRANSFER can only'
71  write(*,*) ' be used within a STEP'
72  call exit(201)
73  endif
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  idirect=2
93  iheat=1
94 !
95  do i=2,n
96  if(textpart(i)(1:7).eq.'SOLVER=') then
97  read(textpart(i)(8:27),'(a20)') solver
98  elseif((textpart(i)(1:6).eq.'DIRECT').and.
99  & (textpart(i)(1:9).ne.'DIRECT=NO')) then
100  idirect=1
101  elseif(textpart(i)(1:9).eq.'DIRECT=NO') then
102  idirect=0
103  elseif(textpart(i)(1:14).eq.'MAGNETOSTATICS') then
104  nmethod=8
105  elseif(textpart(i)(1:9).eq.'FREQUENCY') then
106  nmethod=10
107  elseif(textpart(i)(1:7).eq.'DELTMX=') then
108  read(textpart(i)(8:27),'(f20.0)',iostat=istat) ctrl(27)
109  elseif(textpart(i)(1:9).eq.'TIMERESET') then
110  timereset=.true.
111  elseif(textpart(i)(1:17).eq.'TOTALTIMEATSTART=') then
112  read(textpart(i)(18:37),'(f20.0)',iostat=istat) ttime
113  elseif(textpart(i)(1:14).eq.'NOHEATTRANSFER') then
114  iheat=0
115  else
116  write(*,*)
117  & '*WARNING reading *ELECTROMAGNETICS: parameter not recognized:'
118  write(*,*) ' ',
119  & textpart(i)(1:index(textpart(i),' ')-1)
120  call inputwarning(inpc,ipoinpc,iline,
121  &"*ELECTROMAGNETICS%")
122  endif
123  enddo
124  if(nmethod.eq.8) ctrl(27)=1.d30
125 !
126 ! default for modal dynamic calculations is DIRECT,
127 ! for static or dynamic calculations DIRECT=NO
128 !
129  if(iperturb.eq.0) then
130  idrct=1
131  if(idirect.eq.0)idrct=0
132  else
133  idrct=0
134  if(idirect.eq.1)idrct=1
135  endif
136 !
137  if(nmethod.eq.9) then
138  if(iheat.eq.1) then
139  if(ithermal.eq.0) then
140  write(*,*)
141  & '*ERROR reading *ELECTROMAGNETICS: please define initial '
142  write(*,*) ' conditions for the temperature'
143  call exit(201)
144  endif
145  ithermal=3
146  endif
147  endif
148 !
149  if((nmethod.ne.10).and.(iperturb.ne.0)) then
150 !
151 ! static or dynamic thermal analysis
152 !
153  if(solver(1:7).eq.'SPOOLES') then
154  isolver=0
155  elseif(solver(1:16).eq.'ITERATIVESCALING') then
156  isolver=2
157  elseif(solver(1:17).eq.'ITERATIVECHOLESKY') then
158  isolver=3
159  elseif(solver(1:3).eq.'SGI') then
160  isolver=4
161  elseif(solver(1:5).eq.'TAUCS') then
162  isolver=5
163  elseif(solver(1:7).eq.'PARDISO') then
164  isolver=7
165  else
166  write(*,*)
167  & '*WARNING reading *ELECTROMAGNETICS: unknown solver;'
168  write(*,*) ' the default solver is used'
169  endif
170 !
171  call getnewline(inpc,textpart,istat,n,key,iline,ipol,inl,
172  & ipoinp,inp,ipoinpc)
173  if((istat.lt.0).or.(key.eq.1)) then
174  if(iperturb.ge.2) then
175  write(*,*) '*WARNING reading *ELECTROMAGNETICS: a nonline
176  &ar geometric analysis is requested'
177  write(*,*) ' but no time increment nor step is sp
178  &ecified'
179  write(*,*) ' the defaults (1,1) are used'
180  tinc=1.d0
181  tper=1.d0
182  tmin=1.d-5
183  tmax=1.d+30
184  endif
185  if(timereset)ttime=ttime-tper
186  return
187  endif
188 !
189  read(textpart(1)(1:20),'(f20.0)',iostat=istat) tinc
190  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
191  &"*ELECTROMAGNETICS%")
192  read(textpart(2)(1:20),'(f20.0)',iostat=istat) tper
193  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
194  &"*ELECTROMAGNETICS%")
195  read(textpart(3)(1:20),'(f20.0)',iostat=istat) tmin
196  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
197  &"*ELECTROMAGNETICS%")
198  read(textpart(4)(1:20),'(f20.0)',iostat=istat) tmax
199  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
200  &"*ELECTROMAGNETICS%")
201 !
202  if(tinc.le.0.d0) then
203  write(*,*) '*ERROR reading *ELECTROMAGNETICS: initial increm
204  &ent size is negative'
205  endif
206  if(tper.le.0.d0) then
207  write(*,*)
208  & '*ERROR reading *ELECTROMAGNETICS: step size is negative'
209  endif
210  if(tinc.gt.tper) then
211  write(*,*) '*ERROR reading *ELECTROMAGNETICS: initial increm
212  &ent size exceeds step size'
213  endif
214 !
215  if(idrct.ne.1) then
216 c if(dabs(tmin).lt.1.d-10) then
217 c tmin=min(tinc,1.d-5*tper)
218  if(dabs(tmin).lt.1.d-6*tper) then
219  tmin=min(tinc,1.d-6*tper)
220  endif
221  if(dabs(tmax).lt.1.d-10) then
222  tmax=1.d+30
223  endif
224  endif
225  elseif(nmethod.eq.10) then
226 !
227 ! thermal eigenmode analysis
228 !
229  iperturb=0
230 !
231  call getnewline(inpc,textpart,istat,n,key,iline,ipol,inl,
232  & ipoinp,inp,ipoinpc)
233  if((istat.lt.0).or.(key.eq.1)) then
234  write(*,*)
235  & '*ERROR reading *ELECTROMAGNETICS: definition not complete'
236  write(*,*) ' '
237  call inputerror(inpc,ipoinpc,iline,
238  &"*ELECTROMAGNETICS%")
239  call exit(201)
240  endif
241  read(textpart(1)(1:10),'(i10)',iostat=istat) nev
242  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
243  &"*ELECTROMAGNETICS%")
244  if(nev.le.0) then
245  write(*,*) '*ERROR reading *ELECTROMAGNETICS: less than 1 ei
246  &genvalue requested'
247  call exit(201)
248  endif
249  tol=1.d-2
250  ncv=4*nev
251  ncv=ncv+nev
252  mxiter=1000
253  if(textpart(2)(1:1).ne.' ') then
254  read(textpart(2)(1:20),'(f20.0)',iostat=istat) fmin
255  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
256  &"*ELECTROMAGNETICS%")
257  endif
258  if(textpart(3)(1:1).ne.' ') then
259  read(textpart(3)(1:20),'(f20.0)',iostat=istat) fmax
260  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
261  &"*ELECTROMAGNETICS%")
262  endif
263 !
264  mei(1)=nev
265  mei(2)=ncv
266  mei(3)=mxiter
267  fei(1)=tol
268  fei(2)=fmin
269  fei(3)=fmax
270  endif
271 !
272  if(timereset)ttime=ttime-tper
273 !
274  call getnewline(inpc,textpart,istat,n,key,iline,ipol,inl,
275  & ipoinp,inp,ipoinpc)
276 !
277  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)