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

Go to the source code of this file.

Functions/Subroutines

subroutine bounaddf (iface, is, ie, val, nodeboun, ndirboun, xboun, nboun, nboun_, iamboun, iamplitude, nam, ipompc, nodempc, coefmpc, nmpc, nmpc_, mpcfree, trab, ntrans, ikboun, ilboun, ikmpc, ilmpc, co, nk, nk_, labmpc, type, typeboun, nmethod, iperturb, vold, mi, nelemload, sideload, xload, nload, lakon, ipkon, kon)
 

Function/Subroutine Documentation

◆ bounaddf()

subroutine bounaddf ( integer  iface,
integer  is,
integer  ie,
real*8  val,
integer, dimension(*)  nodeboun,
integer, dimension(*)  ndirboun,
real*8, dimension(*)  xboun,
integer  nboun,
integer  nboun_,
integer, dimension(*)  iamboun,
integer  iamplitude,
integer  nam,
integer, dimension(*)  ipompc,
integer, dimension(3,*)  nodempc,
real*8, dimension(*)  coefmpc,
integer  nmpc,
integer  nmpc_,
integer  mpcfree,
real*8, dimension(7,*)  trab,
integer  ntrans,
integer, dimension(*)  ikboun,
integer, dimension(*)  ilboun,
integer, dimension(*)  ikmpc,
integer, dimension(*)  ilmpc,
real*8, dimension(3,*)  co,
integer  nk,
integer  nk_,
character*20, dimension(*)  labmpc,
character*1  type,
character*1, dimension(*)  typeboun,
integer  nmethod,
integer  iperturb,
real*8, dimension(0:mi(2),*)  vold,
integer, dimension(*)  mi,
integer, dimension(2,*)  nelemload,
character*20, dimension(*)  sideload,
real*8, dimension(2,*)  xload,
integer  nload,
character*8, dimension(*)  lakon,
integer, dimension(*)  ipkon,
integer, dimension(*)  kon 
)
25 !
26 ! adds a facial boundary condition to the data base
27 !
28 ! for facial boundary conditions: dof=-[8*(face-1)+idof]
29 ! (where face = 10*element+local_face_number)
30 ! for nodal boundary conditions: dof=8*(node-1)+idof
31 !
32  implicit none
33 !
34  character*1 type,typeboun(*)
35  character*8 lakon(*)
36  character*20 labmpc(*),label,sideload(*)
37 !
38  integer nodeboun(*),ndirboun(*),is,ie,nboun,nboun_,i,j,
39  & iamboun(*),iamplitude,nam,ipompc(*),nodempc(3,*),nmpc,nmpc_,
40  & mpcfree,ntrans,ikboun(*),ilboun(*),ikmpc(*),ipkon(*),indexe,
41  & ilmpc(*),itr,idof,newnode,number,id,idofnew,idnew,nk,nk_,
42  & mpcfreenew,nmethod,iperturb,ii,mi(*),three,kflag,
43  & iy(3),inumber,iface,nload,nelemload(2,*),nopes,kon(*),nope,
44  & nelem,loadid,ifacel,ifaceq(8,6),ifacet(6,4),ifacew(8,5)
45 !
46  real*8 xboun(*),val,coefmpc(*),trab(7,*),a(3,3),co(3,*),cg(3),
47  & vold(0:mi(2),*),dx(3),xload(2,*)
48 !
49  data ifaceq /4,3,2,1,11,10,9,12,
50  & 5,6,7,8,13,14,15,16,
51  & 1,2,6,5,9,18,13,17,
52  & 2,3,7,6,10,19,14,18,
53  & 3,4,8,7,11,20,15,19,
54  & 4,1,5,8,12,17,16,20/
55  data ifacet /1,3,2,7,6,5,
56  & 1,2,4,5,9,8,
57  & 2,3,4,6,10,9,
58  & 1,4,3,8,10,7/
59  data ifacew /1,3,2,9,8,7,0,0,
60  & 4,5,6,10,11,12,0,0,
61  & 1,2,5,4,7,14,10,13,
62  & 2,3,6,5,8,15,11,14,
63  & 4,6,3,1,12,15,9,13/
64 !
65  if(ntrans.le.0) then
66  itr=0
67  else
68  nelem=int(iface/10.d0)
69  ifacel=iface-10*nelem
70  label(1:20)='T '
71  write(label(2:2),'(i1)') ifacel
72  call identifytransform(nelem,label,nelemload,sideload,
73  & nload,loadid)
74  if(loadid.eq.0) then
75  itr=0
76  else
77  itr=nelemload(2,loadid)
78  endif
79  endif
80 !
81  loop: do ii=is,ie
82  if((itr.eq.0).or.(ii.eq.0).or.(ii.eq.11).or.(ii.eq.8)) then
83 !
84 ! no transformation applies: simple SPC
85 !
86  if(ii.le.3) then
87  i=ii
88  elseif(ii.eq.4) then
89  write(*,*) '*ERROR in bounaddf: a boundary condition'
90  write(*,*) ' on DOF 4 is not allowed'
91  call exit(201)
92  elseif(ii.eq.5) then
93  write(*,*) '*ERROR in bounaddf: a boundary condition'
94  write(*,*) ' on DOF 5 is not allowed'
95  call exit(201)
96  elseif(ii.eq.6) then
97  write(*,*) '*ERROR in bounaddf: a boundary condition'
98  write(*,*) ' on DOF 6 is not allowed'
99  call exit(201)
100  elseif(ii.eq.8) then
101  i=4
102  elseif(ii.eq.11) then
103  i=0
104  else
105  write(*,*) '*ERROR in bounadd: unknown DOF: ',
106  & ii
107  call exit(201)
108  endif
109  idof=-(8*(iface-1)+i)
110  call nident(ikboun,idof,nboun,id)
111  if(id.gt.0) then
112  if(ikboun(id).eq.idof) then
113  j=ilboun(id)
114  if(typeboun(j).ne.type) cycle loop
115  xboun(j)=val
116  if(nam.gt.0) iamboun(j)=iamplitude
117  cycle loop
118  endif
119  endif
120  nboun=nboun+1
121 c write(*,*) 'bounaddf boun',iface,i,val,idof
122  if(nboun.gt.nboun_) then
123  write(*,*) '*ERROR in bounadd: increase nboun_'
124  call exit(201)
125  endif
126  if((nmethod.eq.4).and.(iperturb.le.1)) then
127  write(*,*) '*ERROR in bounadd: in a modal dynamic step'
128  write(*,*) ' new SPCs are not allowed'
129  call exit(201)
130  endif
131  nodeboun(nboun)=iface
132  ndirboun(nboun)=i
133  xboun(nboun)=val
134  typeboun(nboun)=type
135  if(nam.gt.0) iamboun(nboun)=iamplitude
136 !
137 ! updating ikboun and ilboun
138 !
139  do j=nboun,id+2,-1
140  ikboun(j)=ikboun(j-1)
141  ilboun(j)=ilboun(j-1)
142  enddo
143  ikboun(id+1)=idof
144  ilboun(id+1)=nboun
145  else
146 !
147 ! transformation applies: SPC is MPC in global carthesian
148 ! coordinates
149 !
150 ! number of nodes belonging to the face
151 !
152  if(lakon(nelem)(4:4).eq.'8') then
153  nope=8
154  nopes=4
155  elseif(lakon(nelem)(4:4).eq.'4') then
156  nope=4
157  nopes=3
158  elseif(lakon(nelem)(4:4).eq.'6') then
159  nope=6
160  if(ifacel.le.2) then
161  nopes=3
162  else
163  nopes=4
164  endif
165  endif
166 !
167 ! determining the center of gravity
168 !
169  do j=1,3
170  cg(j)=0.d0
171  enddo
172 !
173  indexe=ipkon(nelem)
174  if(nope.eq.8) then
175  do i=1,nopes
176  do j=1,3
177  cg(j)=cg(j)+co(j,kon(indexe+ifaceq(i,ifacel)))
178  enddo
179  enddo
180  elseif(nope.eq.4) then
181  do i=1,nopes
182  do j=1,3
183  cg(j)=cg(j)+co(j,kon(indexe+ifacet(i,ifacel)))
184  enddo
185  enddo
186  else
187  do i=1,nopes
188  do j=1,3
189  cg(j)=cg(j)+co(j,kon(indexe+ifacew(i,ifacel)))
190  enddo
191  enddo
192  endif
193  do j=1,3
194  cg(j)=cg(j)/nopes
195  enddo
196 !
197 ! determining the transformation coefficients at the center
198 ! of gravity
199 !
200  call transformatrix(trab(1,itr),cg,a)
201  if(ii.le.3) then
202  i=ii
203  elseif(ii.eq.4) then
204  write(*,*) '*ERROR in bounaddf: a boundary condition'
205  write(*,*) ' on DOF 4 is not allowed'
206  call exit(201)
207  elseif(ii.eq.5) then
208  write(*,*) '*ERROR in bounaddf: a boundary condition'
209  write(*,*) ' on DOF 5 is not allowed'
210  call exit(201)
211  elseif(ii.eq.6) then
212  write(*,*) '*ERROR in bounaddf: a boundary condition'
213  write(*,*) ' on DOF 6 is not allowed'
214  call exit(201)
215  elseif(ii.eq.8) then
216  i=4
217  elseif(ii.eq.11) then
218  i=0
219  else
220  write(*,*) '*ERROR in bounadd: unknown DOF: ',
221  & ii
222  call exit(201)
223  endif
224  if(int(xload(1,loadid)).ne.0) then
225  newnode=int(xload(1,loadid))
226  idofnew=8*(newnode-1)+i
227  call nident(ikboun,idofnew,nboun,idnew)
228  if(idnew.gt.0) then
229  if(ikboun(idnew).eq.idofnew) then
230  j=ilboun(idnew)
231  if(typeboun(j).ne.type) cycle
232  xboun(j)=val
233  if(nam.gt.0) iamboun(j)=iamplitude
234  cycle
235  endif
236  endif
237  else
238 !
239 ! new node is generated for the inhomogeneous MPC term
240 !
241  if((nmethod.eq.4).and.(iperturb.le.1)) then
242  write(*,*)'*ERROR in bounadd: in a modal dynamic step'
243  write(*,*) ' new SPCs are not allowed'
244  call exit(201)
245  endif
246  nk=nk+1
247  if(nk.gt.nk_) then
248  write(*,*) '*ERROR in bounadd: increase nk_'
249  call exit(201)
250  endif
251  newnode=nk
252  xload(1,loadid)=newnode+0.5d0
253  idofnew=8*(newnode-1)+i
254  idnew=nboun
255  endif
256 !
257 ! new mpc
258 !
259  iy(1)=1
260  iy(2)=2
261  iy(3)=3
262  dx(1)=dabs(a(1,i))
263  dx(2)=dabs(a(2,i))
264  dx(3)=dabs(a(3,i))
265  three=3
266  kflag=-2
267  call dsort(dx,iy,three,kflag)
268  do inumber=1,3
269  number=iy(inumber)
270  idof=-(8*(iface-1)+number)
271  call nident(ikmpc,idof,nmpc,id)
272  if(id.ne.0) then
273  if(ikmpc(id).eq.idof) cycle
274  endif
275  if(dabs(a(number,i)).lt.1.d-5) cycle
276  nmpc=nmpc+1
277  if(nmpc.gt.nmpc_) then
278  write(*,*) '*ERROR in bounadd: increase nmpc_'
279  call exit(201)
280  endif
281  labmpc(nmpc)='FLUIDSPC '
282 c write(*,*) nmpc,labmpc(nmpc),'bounaddf'
283  ipompc(nmpc)=mpcfree
284  do j=nmpc,id+2,-1
285  ikmpc(j)=ikmpc(j-1)
286  ilmpc(j)=ilmpc(j-1)
287  enddo
288  ikmpc(id+1)=idof
289  ilmpc(id+1)=nmpc
290  exit
291  enddo
292 !
293  inumber=inumber-1
294  do j=1,3
295  inumber=inumber+1
296  if(inumber.gt.3) inumber=1
297  number=iy(inumber)
298  if(dabs(a(number,i)).lt.1.d-30) cycle
299  nodempc(1,mpcfree)=iface
300  nodempc(2,mpcfree)=number
301  coefmpc(mpcfree)=a(number,i)
302  mpcfree=nodempc(3,mpcfree)
303  if(mpcfree.eq.0) then
304  write(*,*) '*ERROR in bounadd: increase memmpc_'
305  call exit(201)
306  endif
307  enddo
308 !
309 ! storage of the boundary condition (faster than
310 ! storage of the new node); the negative of the
311 ! condition number is stored as tag for a SPC
312 !
313  nodempc(1,mpcfree)=-(nboun+1)
314  nodempc(2,mpcfree)=i
315  coefmpc(mpcfree)=-1.d0
316  mpcfreenew=nodempc(3,mpcfree)
317  if(mpcfreenew.eq.0) then
318  write(*,*) '*ERROR in bounadd: increase nmpc_'
319  call exit(201)
320  endif
321  nodempc(3,mpcfree)=0
322  mpcfree=mpcfreenew
323 !
324 ! nonhomogeneous term
325 !
326  nboun=nboun+1
327 c write(*,*) 'bounaddf inhom',newnode,i,val,idofnew
328  if(nboun.gt.nboun_) then
329  write(*,*) '*ERROR in bounadd: increase nboun_'
330  call exit(201)
331  endif
332  nodeboun(nboun)=newnode
333  ndirboun(nboun)=i
334  xboun(nboun)=val
335  typeboun(nboun)=type
336  if(nam.gt.0) iamboun(nboun)=iamplitude
337 !
338 ! updating ikboun and ilboun
339 !
340  do j=nboun,idnew+2,-1
341  ikboun(j)=ikboun(j-1)
342  ilboun(j)=ilboun(j-1)
343  enddo
344  ikboun(idnew+1)=idofnew
345  ilboun(idnew+1)=nboun
346 !
347  endif
348  enddo loop
349 !
350  return
subroutine identifytransform(nelement, label, nelemload, sideload, nload, loadid)
Definition: identifytransform.f:21
subroutine transformatrix(xab, p, a)
Definition: transformatrix.f:20
subroutine nident(x, px, n, id)
Definition: nident.f:26
subroutine dsort(dx, iy, n, kflag)
Definition: dsort.f:6
Hosted by OpenAircraft.com, (Michigan UAV, LLC)