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

Go to the source code of this file.

Functions/Subroutines

subroutine boundarys (inpc, textpart, set, istartset, iendset, ialset, nset, nodeboun, ndirboun, xboun, nboun, nboun_, nk, iamboun, amname, nam, ipompc, nodempc, coefmpc, nmpc, nmpc_, mpcfree, inotr, trab, ntrans, ikboun, ilboun, ikmpc, ilmpc, nk_, co, labmpc, boun_flag, typeboun, istat, n, iline, ipol, inl, ipoinp, inp, nam_, namtot_, namta, amta, nmethod, iperturb, iaxial, ipoinpc, vold, mi)
 

Function/Subroutine Documentation

◆ boundarys()

subroutine boundarys ( character*1, dimension(*)  inpc,
character*132, dimension(16)  textpart,
character*81, dimension(*)  set,
integer, dimension(*)  istartset,
integer, dimension(*)  iendset,
integer, dimension(*)  ialset,
integer  nset,
integer, dimension(*)  nodeboun,
integer, dimension(*)  ndirboun,
real*8, dimension(*)  xboun,
integer  nboun,
integer  nboun_,
integer  nk,
integer, dimension(*)  iamboun,
character*80, dimension(*)  amname,
integer  nam,
integer, dimension(*)  ipompc,
integer, dimension(3,*)  nodempc,
real*8, dimension(*)  coefmpc,
integer  nmpc,
integer  nmpc_,
integer  mpcfree,
integer, dimension(2,*)  inotr,
real*8, dimension(7,*)  trab,
integer  ntrans,
integer, dimension(*)  ikboun,
integer, dimension(*)  ilboun,
integer, dimension(*)  ikmpc,
integer, dimension(*)  ilmpc,
integer  nk_,
real*8, dimension(3,*)  co,
character*20, dimension(*)  labmpc,
logical  boun_flag,
character*1, dimension(*)  typeboun,
integer  istat,
integer  n,
integer  iline,
integer  ipol,
integer  inl,
integer, dimension(2,*)  ipoinp,
integer, dimension(3,*)  inp,
integer  nam_,
integer  namtot_,
integer, dimension(3,*)  namta,
real*8, dimension(2,*)  amta,
integer  nmethod,
integer  iperturb,
integer  iaxial,
integer, dimension(0:*)  ipoinpc,
real*8, dimension(0:mi(2),*)  vold,
integer, dimension(*)  mi 
)
26 !
27 ! reading the input deck: *BOUNDARY
28 !
29  implicit none
30 !
31  logical boun_flag,user,massflowrate,fixed,submodel
32 !
33  character*1 typeboun(*),type,inpc(*)
34  character*20 labmpc(*),label
35  character*80 amname(*),amplitude
36  character*81 set(*),noset
37  character*132 textpart(16)
38 !
39  integer istartset(*),iendset(*),ialset(*),nodeboun(*),
40  & ndirboun(*),
41  & nset,nboun,nboun_,istat,n,i,j,k,l,ibounstart,ibounend,
42  & key,nk,iamboun(*),nam,iamplitude,ipompc(*),nodempc(3,*),
43  & nmpc,nmpc_,mpcfree,inotr(2,*),ikboun(*),ilboun(*),ikmpc(*),
44  & ilmpc(*),nmpcold,id,idof,index1,ntrans,nk_,ipos,m,node,is,ie,
45  & iline,ipol,inl,ipoinp(2,*),inp(3,*),nam_,namtot,namtot_,
46  & namta(3,*),idelay,nmethod,iperturb,lc,iaxial,ipoinpc(0:*),
47  & ktrue,mi(*),iglobstep
48 !
49  real*8 xboun(*),bounval,coefmpc(*),trab(7,*),co(3,*),amta(2,*),
50  & vold(0:mi(2),*)
51 !
52  type='B'
53  label=' '
54  iamplitude=0
55  idelay=0
56  user=.false.
57  massflowrate=.false.
58  fixed=.false.
59  lc=1
60  iglobstep=0
61  submodel=.false.
62 !
63  do i=2,n
64  if((textpart(i)(1:6).eq.'OP=NEW').and.(.not.boun_flag)) then
65 !
66 ! spc's in nonglobal coordinates result in mpc's
67 ! removing these mpc's
68 ! necessary and sufficient condition for a MPC to be removed:
69 ! - on the dependent side a node "a" corresponding to SPC "b"
70 ! (no matter which DOF); SPC "b" is applied in direction "c" of
71 ! node "a" and corresponds to node "d" to account for the
72 ! inhomogeneous term
73 ! - on the independent side a term for node "d" in direction "c".
74 !
75  if(ntrans.gt.0) then
76  nmpcold=nmpc
77  do j=1,nk
78  if(inotr(2,j).gt.0) then
79  do k=1,3
80  idof=8*(inotr(2,j)-1)+k
81  call nident(ikboun,idof,nboun,id)
82  if(id.gt.0) then
83  if(ikboun(id).eq.idof) then
84 !
85 ! if a SPC is defined in direction k for a node j for which a
86 ! local coordinate system applies, then the coordinate system
87 ! number is stored in inotr(1,j) and the additional node
88 ! for the inhomogeneous term is stored in inotr(2,j). The
89 ! SPC DOF is (inotr(2,j)-1)*3+k, however, the independent
90 ! MPC DOF is (j-1)*3+l, where l can be different from k,
91 ! since (j-1)*3+k might already be taken by another MPC, or
92 ! the coefficient for this direction might be zero.
93 !
94  loop: do l=1,3
95  idof=8*(j-1)+l
96  call nident(ikmpc,idof,nmpc,id)
97  if(id.gt.0) then
98  if(ikmpc(id).eq.idof) then
99  index1=ipompc(ilmpc(id))
100  if(index1.eq.0) cycle
101  do
102  if((nodempc(1,index1).eq.
103  & inotr(2,j)).and.
104  & (nodempc(2,index1).eq.k))
105  & then
106  nodempc(3,index1)=mpcfree
107  mpcfree=ipompc(ilmpc(id))
108  ipompc(ilmpc(id))=0
109  do m=id,nmpc-1
110  ikmpc(m)=ikmpc(m+1)
111  ilmpc(m)=ilmpc(m+1)
112  enddo
113  ikmpc(nmpc)=0
114  ilmpc(nmpc)=0
115  nmpc=nmpc-1
116  exit
117  endif
118  index1=nodempc(3,index1)
119  if(index1.eq.0) exit
120  enddo
121  endif
122  endif
123  enddo loop
124 !
125  endif
126  endif
127  enddo
128  endif
129  enddo
130 !
131 ! getting rid of the superfluous lines in ipompc and labmpc
132 !
133  k=0
134  do j=1,nmpcold
135  if(ipompc(j).ne.0) then
136  k=k+1
137  ipompc(k)=ipompc(j)
138  labmpc(k)=labmpc(j)
139  index1=ipompc(j)
140  idof=8*(nodempc(1,index1)-1)+nodempc(2,index1)
141  call nident(ikmpc,idof,nmpc,id)
142  if(id.eq.0) then
143  write(*,*) '*ERROR reading *BOUNDARY'
144  call exit(201)
145  elseif(ikmpc(id).ne.idof) then
146  write(*,*) '*ERROR reading *BOUNDARY'
147  call exit(201)
148  endif
149  ilmpc(id)=k
150  endif
151  enddo
152  endif
153 !
154 ! removing the boundary conditions defined by a *BOUNDARY
155 ! statement
156 !
157  loop1: do
158  if(nboun.gt.0) then
159  do j=1,nboun
160  if(typeboun(j).eq.'B') then
161  node=nodeboun(j)
162  is=ndirboun(j)
163  ie=ndirboun(j)
164  call bounrem(node,is,j,nodeboun,ndirboun,xboun,
165  & nboun,iamboun,nam,ikboun,ilboun,typeboun)
166  cycle loop1
167  endif
168  enddo
169  exit
170  endif
171  exit
172  enddo loop1
173  elseif(textpart(i)(1:10).eq.'AMPLITUDE=') then
174  read(textpart(i)(11:90),'(a80)') amplitude
175  do j=nam,1,-1
176  if(amname(j).eq.amplitude) then
177  iamplitude=j
178  exit
179  endif
180  enddo
181  if(j.eq.0) then
182  write(*,*)
183  & '*ERROR reading *BOUNDARY: nonexistent amplitude'
184  write(*,*) ' '
185  call inputerror(inpc,ipoinpc,iline,
186  &"*BOUNDARY%")
187  call exit(201)
188  endif
189  iamplitude=j
190  elseif(textpart(i)(1:10).eq.'TIMEDELAY=') THEN
191  if(idelay.ne.0) then
192  write(*,*) '*ERROR reading *BOUNDARY: the parameter TIME'
193  write(*,*) ' DELAY is used twice in the same'
194  write(*,*) ' keyword; '
195  call inputerror(inpc,ipoinpc,iline,
196  &"*BOUNDARY%")
197  call exit(201)
198  else
199  idelay=1
200  endif
201  nam=nam+1
202  if(nam.gt.nam_) then
203  write(*,*) '*ERROR reading *BOUNDARY: increase nam_'
204  call exit(201)
205  endif
206  amname(nam)='
207  & '
208  if(iamplitude.eq.0) then
209  write(*,*) '*ERROR reading *BOUNDARY: time delay must be'
210  write(*,*) ' preceded by the amplitude parameter'
211  call exit(201)
212  endif
213  namta(3,nam)=sign(iamplitude,namta(3,iamplitude))
214  iamplitude=nam
215  if(nam.eq.1) then
216  namtot=0
217  else
218  namtot=namta(2,nam-1)
219  endif
220  namtot=namtot+1
221  if(namtot.gt.namtot_) then
222  write(*,*) '*ERROR boundaries: increase namtot_'
223  call exit(201)
224  endif
225  namta(1,nam)=namtot
226  namta(2,nam)=namtot
227  read(textpart(i)(11:30),'(f20.0)',iostat=istat)
228  & amta(1,namtot)
229  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
230  &"*BOUNDARY%")
231  elseif(textpart(i)(1:9).eq.'LOADCASE=') then
232  read(textpart(i)(10:19),'(i10)',iostat=istat) lc
233  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
234  &"*BOUNDARY%")
235  if(nmethod.ne.5) then
236  write(*,*) '*ERROR reading *BOUNDARY: the parameter LOAD'
237  write(*,*) ' CASE is only allowed in STEADY STATE'
238  write(*,*) ' DYNAMICS calculations'
239  call exit(201)
240  endif
241  elseif(textpart(i)(1:4).eq.'USER') then
242  user=.true.
243  elseif(textpart(i)(1:8).eq.'MASSFLOW') then
244  massflowrate=.true.
245  elseif(textpart(i)(1:5).eq.'FIXED') then
246  fixed=.true.
247  elseif(textpart(i)(1:8).eq.'SUBMODEL') then
248  submodel=.true.
249  elseif(textpart(i)(1:5).eq.'STEP=') then
250  read(textpart(i)(6:15),'(i10)',iostat=istat) iglobstep
251  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
252  &"*BOUNDARY%")
253  else
254  write(*,*)
255  & '*WARNING reading *BOUNDARY: parameter not recognized:'
256  write(*,*) ' ',
257  & textpart(i)(1:index(textpart(i),' ')-1)
258  call inputwarning(inpc,ipoinpc,iline,
259  &"*BOUNDARY%")
260  endif
261  enddo
262 !
263 ! check whether global step was specified for submodel
264 !
265  if((submodel).and.(iglobstep.eq.0)) then
266  write(*,*) '*ERROR reading *BOUNDARY: no global step'
267  write(*,*) ' step specified for the submodel'
268  call inputerror(inpc,ipoinpc,iline,
269  &"*BOUNDARY%")
270  endif
271 !
272 ! storing the step for submodels in iamboun
273 !
274  if(submodel) then
275  if(iamplitude.ne.0) then
276  write(*,*) '*WARNING reading *BOUNDARY:'
277  write(*,*) ' no amplitude definition is allowed'
278  write(*,*) ' in combination with a submodel'
279  endif
280  iamplitude=iglobstep
281  endif
282 !
283  if(user.and.(iamplitude.ne.0)) then
284  write(*,*) '*WARNING reading *BOUNDARY: '
285  write(*,*) ' no amplitude definition is allowed'
286  write(*,*) ' for boundary conditions defined by a'
287  write(*,*) ' user routine'
288  iamplitude=0
289  endif
290 !
291  do
292  call getnewline(inpc,textpart,istat,n,key,iline,ipol,inl,
293  & ipoinp,inp,ipoinpc)
294  if((istat.lt.0).or.(key.eq.1)) return
295 !
296  read(textpart(2)(1:10),'(i10)',iostat=istat) ibounstart
297  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
298  &"*BOUNDARY%")
299 !
300  if(textpart(3)(1:1).eq.' ') then
301  ibounend=ibounstart
302  else
303  read(textpart(3)(1:10),'(i10)',iostat=istat) ibounend
304  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
305  &"*BOUNDARY%")
306  endif
307 !
308  if(textpart(4)(1:1).eq.' ') then
309  bounval=0.d0
310  else
311  read(textpart(4)(1:20),'(f20.0)',iostat=istat) bounval
312  if(istat.gt.0) call inputerror(inpc,ipoinpc,iline,
313  &"*BOUNDARY%")
314  endif
315  if((massflowrate).and.(iaxial.eq.180)) bounval=bounval/iaxial
316 !
317 ! dummy boundary condition consisting of the first primes
318 !
319  if(user) bounval=1.2357111317d0
320  if(submodel) bounval=1.9232931374d0
321 !
322  read(textpart(1)(1:10),'(i10)',iostat=istat) l
323  if(istat.eq.0) then
324  if((l.gt.nk).or.(l.le.0)) then
325  write(*,*) '*ERROR reading *BOUNDARY:'
326  write(*,*) ' node ',l,' is not defined'
327  call exit(201)
328  endif
329  if(submodel) then
330  if(ntrans.gt.0) then
331  if(inotr(1,l).gt.0) then
332  write(*,*) '*ERROR reading *BOUNDARY: in submodel'
333  write(*,*) ' node',l,' a local coordinate'
334  write(*,*) ' system was defined. This is not'
335  write(*,*) ' allowed'
336  call exit(201)
337  endif
338  endif
339  endif
340  ktrue=l
341  if(lc.ne.1) l=l+nk
342  call bounadd(l,ibounstart,ibounend,bounval,
343  & nodeboun,ndirboun,xboun,nboun,nboun_,
344  & iamboun,iamplitude,nam,ipompc,nodempc,
345  & coefmpc,nmpc,nmpc_,mpcfree,inotr,trab,
346  & ntrans,ikboun,ilboun,ikmpc,ilmpc,co,nk,nk_,labmpc,
347  & type,typeboun,nmethod,iperturb,fixed,vold,ktrue,
348  & mi,label)
349  else
350  read(textpart(1)(1:80),'(a80)',iostat=istat) noset
351  noset(81:81)=' '
352  ipos=index(noset,' ')
353  noset(ipos:ipos)='N'
354  do i=1,nset
355  if(set(i).eq.noset) exit
356  enddo
357  if(i.gt.nset) then
358  noset(ipos:ipos)=' '
359  write(*,*) '*ERROR reading *BOUNDARY: node set ',noset
360  write(*,*) ' has not yet been defined. '
361  call inputerror(inpc,ipoinpc,iline,
362  &"*BOUNDARY%")
363  call exit(201)
364  endif
365  do j=istartset(i),iendset(i)
366  if(ialset(j).gt.0) then
367  k=ialset(j)
368  if(submodel) then
369  if(ntrans.gt.0) then
370  if(inotr(1,k).gt.0) then
371  write(*,*)
372  & '*ERROR reading *BOUNDARY: in submodel'
373  write(*,*) ' node',k,
374  & ' a local coordinate'
375  write(*,*)
376  & ' system was defined. This is not'
377  write(*,*) ' allowed'
378  call exit(201)
379  endif
380  endif
381  endif
382  ktrue=k
383  if(lc.ne.1) k=k+nk
384  call bounadd(k,ibounstart,ibounend,bounval,
385  & nodeboun,ndirboun,xboun,nboun,nboun_,
386  & iamboun,iamplitude,nam,ipompc,nodempc,
387  & coefmpc,nmpc,nmpc_,mpcfree,inotr,trab,
388  & ntrans,ikboun,ilboun,ikmpc,ilmpc,co,nk,nk_,labmpc,
389  & type,typeboun,nmethod,iperturb,fixed,vold,ktrue,
390  & mi,label)
391  else
392  k=ialset(j-2)
393  do
394  k=k-ialset(j)
395  if(k.ge.ialset(j-1)) exit
396  if(submodel) then
397  if(ntrans.gt.0) then
398  if(inotr(1,k).gt.0) then
399  write(*,*)
400  & '*ERROR reading *BOUNDARY: in submodel'
401  write(*,*) ' node',k,
402  & ' a local coordinate'
403  write(*,*)
404  & ' system was defined. This is not'
405  write(*,*) ' allowed'
406  call exit(201)
407  endif
408  endif
409  endif
410  ktrue=k
411  if(lc.ne.1) k=k+nk
412  call bounadd(k,ibounstart,ibounend,bounval,
413  & nodeboun,ndirboun,xboun,nboun,nboun_,
414  & iamboun,iamplitude,nam,ipompc,nodempc,
415  & coefmpc,nmpc,nmpc_,mpcfree,inotr,trab,
416  & ntrans,ikboun,ilboun,ikmpc,ilmpc,co,nk,nk_,
417  & labmpc,type,typeboun,nmethod,iperturb,fixed,
418  & vold,ktrue,mi,label)
419  enddo
420  endif
421  enddo
422  endif
423  enddo
424 !
425  return
subroutine inputwarning(inpc, ipoinpc, iline, text)
Definition: inputwarning.f:20
subroutine bounrem(node, is, iboun, nodeboun, ndirboun, xboun, nboun, iamboun, nam, ikboun, ilboun, typeboun)
Definition: bounrem.f:21
subroutine getnewline(inpc, textpart, istat, n, key, iline, ipol, inl, ipoinp, inp, ipoinpc)
Definition: getnewline.f:21
subroutine nident(x, px, n, id)
Definition: nident.f:26
subroutine bounadd(node, is, ie, val, nodeboun, ndirboun, xboun, nboun, nboun_, iamboun, iamplitude, nam, ipompc, nodempc, coefmpc, nmpc, nmpc_, mpcfree, inotr, trab, ntrans, ikboun, ilboun, ikmpc, ilmpc, co, nk, nk_, labmpc, type, typeboun, nmethod, iperturb, fixed, vold, nodetrue, mi, label)
Definition: bounadd.f:24
subroutine inputerror(inpc, ipoinpc, iline, text)
Definition: inputerror.f:20
Hosted by OpenAircraft.com, (Michigan UAV, LLC)