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

Go to the source code of this file.

Functions/Subroutines

subroutine createtiedsurfs (nodface, ipoface, set, istartset, iendset, ialset, tieset, inomat, ne, ipkon, lakon, kon, ntie, tietol, nalset, nk, nset, iactive)
 

Function/Subroutine Documentation

◆ createtiedsurfs()

subroutine createtiedsurfs ( integer, dimension(5,*)  nodface,
integer, dimension(*)  ipoface,
character*81, dimension(*)  set,
integer, dimension(*)  istartset,
integer, dimension(*)  iendset,
integer, dimension(*)  ialset,
character*81, dimension(3,*)  tieset,
integer, dimension(*)  inomat,
integer  ne,
integer, dimension(*)  ipkon,
character*8, dimension(*)  lakon,
integer, dimension(*)  kon,
integer  ntie,
real*8, dimension(3,*)  tietol,
integer  nalset,
integer  nk,
integer  nset,
integer, dimension(3)  iactive 
)
22 !
23  implicit none
24 !
25 ! creates ties for the surfaces in between domains of an
26 ! electromagnetic calculation
27 !
28  character*8 lakon(*)
29  character*81 set(*),tieset(3,*)
30 !
31  integer ipkon(*),kon(*),ne,nodface(5,*),ipoface(*),istartset(*),
32  & iendset(*),ialset(*),inomat(*),ithree,ifour,ifaceq(8,6),
33  & ifacet(6,4),ifacew(8,5),ifree,ifreenew,index,indexold,kflag,
34  & i,j,k,iactive(3),ntie,nodes(4),iaux,nalset,nk,nset,indexe
35 !
36  real*8 tietol(3,*)
37 !
38 ! nodes belonging to the element faces
39 !
40  data ifaceq /4,3,2,1,11,10,9,12,
41  & 5,6,7,8,13,14,15,16,
42  & 1,2,6,5,9,18,13,17,
43  & 2,3,7,6,10,19,14,18,
44  & 3,4,8,7,11,20,15,19,
45  & 4,1,5,8,12,17,16,20/
46  data ifacet /1,3,2,7,6,5,
47  & 1,2,4,5,9,8,
48  & 2,3,4,6,10,9,
49  & 1,4,3,8,10,7/
50  data ifacew /1,3,2,9,8,7,0,0,
51  & 4,5,6,10,11,12,0,0,
52  & 1,2,5,4,7,14,10,13,
53  & 2,3,6,5,8,15,11,14,
54  & 4,6,3,1,12,15,9,13/
55 !
56  kflag=1
57  ithree=3
58  ifour=4
59 !
60 ! determining the external element faces of the electromagnetic mesh
61 ! the faces are catalogued by the three lowes nodes numbers
62 ! in ascending order. ipoface(i) points to a face for which
63 ! node i is the lowest node and nodface(1,ipoface(i)) and
64 ! nodface(2,ipoface(i)) are the next lower ones.
65 ! nodface(3,ipoface(i)) contains the element number,
66 ! nodface(4,ipoface(i)) the face number and nodface(5,ipoface(i))
67 ! is a pointer to the next surface for which node i is the
68 ! lowest node; if there are no more such surfaces the pointer
69 ! has the value zero
70 ! An external element face is one which belongs to one element
71 ! only
72 !
73  ifree=1
74  do i=1,6*ne-1
75  nodface(5,i)=i+1
76  enddo
77  do i=1,ne
78  if(ipkon(i).lt.0) cycle
79  if(lakon(i)(1:3).ne.'C3D') cycle
80  indexe=ipkon(i)
81  if((lakon(i)(4:4).eq.'2').or.(lakon(i)(4:4).eq.'8')) then
82  do j=1,6
83  do k=1,4
84  nodes(k)=kon(indexe+ifaceq(k,j))
85  enddo
86  call isortii(nodes,iaux,ifour,kflag)
87  indexold=0
88  index=ipoface(nodes(1))
89  do
90 !
91 ! adding a surface which has not been
92 ! catalogued so far
93 !
94  if(index.eq.0) then
95  ifreenew=nodface(5,ifree)
96  nodface(1,ifree)=nodes(2)
97  nodface(2,ifree)=nodes(3)
98  nodface(3,ifree)=i
99  nodface(4,ifree)=j
100  nodface(5,ifree)=ipoface(nodes(1))
101  ipoface(nodes(1))=ifree
102  ifree=ifreenew
103  exit
104  endif
105 !
106 ! removing a surface which has already
107 ! been catalogued
108 !
109  if((nodface(1,index).eq.nodes(2)).and.
110  & (nodface(2,index).eq.nodes(3))) then
111  if(indexold.eq.0) then
112  ipoface(nodes(1))=nodface(5,index)
113  else
114  nodface(5,indexold)=nodface(5,index)
115  endif
116  nodface(5,index)=ifree
117  ifree=index
118  exit
119  endif
120  indexold=index
121  index=nodface(5,index)
122  enddo
123  enddo
124  elseif((lakon(i)(4:4).eq.'4').or.(lakon(i)(4:5).eq.'10')) then
125  do j=1,4
126  do k=1,3
127  nodes(k)=kon(indexe+ifacet(k,j))
128  enddo
129  call isortii(nodes,iaux,ithree,kflag)
130  indexold=0
131  index=ipoface(nodes(1))
132  do
133 !
134 ! adding a surface which has not been
135 ! catalogued so far
136 !
137  if(index.eq.0) then
138  ifreenew=nodface(5,ifree)
139  nodface(1,ifree)=nodes(2)
140  nodface(2,ifree)=nodes(3)
141  nodface(3,ifree)=i
142  nodface(4,ifree)=j
143  nodface(5,ifree)=ipoface(nodes(1))
144  ipoface(nodes(1))=ifree
145  ifree=ifreenew
146  exit
147  endif
148 !
149 ! removing a surface which has already
150 ! been catalogued
151 !
152  if((nodface(1,index).eq.nodes(2)).and.
153  & (nodface(2,index).eq.nodes(3))) then
154  if(indexold.eq.0) then
155  ipoface(nodes(1))=nodface(5,index)
156  else
157  nodface(5,indexold)=nodface(5,index)
158  endif
159  nodface(5,index)=ifree
160  ifree=index
161  exit
162  endif
163  indexold=index
164  index=nodface(5,index)
165  enddo
166  enddo
167  else
168  do j=1,5
169  if(j.le.2) then
170  do k=1,3
171  nodes(k)=kon(indexe+ifacew(k,j))
172  enddo
173  call isortii(nodes,iaux,ithree,kflag)
174  else
175  do k=1,4
176  nodes(k)=kon(indexe+ifacew(k,j))
177  enddo
178  call isortii(nodes,iaux,ifour,kflag)
179  endif
180  indexold=0
181  index=ipoface(nodes(1))
182  do
183 !
184 ! adding a surface which has not been
185 ! catalogued so far
186 !
187  if(index.eq.0) then
188  ifreenew=nodface(5,ifree)
189  nodface(1,ifree)=nodes(2)
190  nodface(2,ifree)=nodes(3)
191  nodface(3,ifree)=i
192  nodface(4,ifree)=j
193  nodface(5,ifree)=ipoface(nodes(1))
194  ipoface(nodes(1))=ifree
195  ifree=ifreenew
196  exit
197  endif
198 !
199 ! removing a surface which has already
200 ! been catalogued
201 !
202  if((nodface(1,index).eq.nodes(2)).and.
203  & (nodface(2,index).eq.nodes(3))) then
204  if(indexold.eq.0) then
205  ipoface(nodes(1))=nodface(5,index)
206  else
207  nodface(5,indexold)=nodface(5,index)
208  endif
209  nodface(5,index)=ifree
210  ifree=index
211  exit
212  endif
213  indexold=index
214  index=nodface(5,index)
215  enddo
216  enddo
217  endif
218  enddo
219 !
220 ! boundary of domain 1 (phi-domain)
221 !
222  nset=nset+1
223  istartset(nset)=nalset+1
224  do i=1,nk
225  if(ipoface(i).eq.0) cycle
226  if(inomat(i).ne.1) cycle
227  index=ipoface(i)
228  do
229  if(index.eq.0) exit
230  nalset=nalset+1
231  ialset(nalset)=10*nodface(3,index)+nodface(4,index)
232  index=nodface(5,index)
233  enddo
234  enddo
235  if(istartset(nset).gt.nalset) then
236  iactive(1)=0
237  nset=nset-1
238  else
239  iactive(1)=nset
240  set(nset)(1:22)='ELECTROMAGNETICSZONE1T'
241  do i=23,81
242  set(nset)(i:i)=' '
243  enddo
244  iendset(nset)=nalset
245  endif
246 !
247 ! boundary of domain 2 (A-V-domain)
248 !
249  nset=nset+1
250  istartset(nset)=nalset+1
251  do i=1,nk
252  if(ipoface(i).eq.0) cycle
253  if(inomat(i).ne.2) cycle
254  index=ipoface(i)
255  do
256  if(index.eq.0) exit
257  nalset=nalset+1
258  ialset(nalset)=10*nodface(3,index)+nodface(4,index)
259  index=nodface(5,index)
260  enddo
261  enddo
262  if(istartset(nset).gt.nalset) then
263  iactive(2)=0
264  nset=nset-1
265  else
266  iactive(2)=nset
267  set(nset)(1:22)='ELECTROMAGNETICSZONE2T'
268  do i=23,81
269  set(nset)(i:i)=' '
270  enddo
271  iendset(nset)=nalset
272  endif
273 !
274 ! boundary of domain 3 (A-domain)
275 !
276  nset=nset+1
277  istartset(nset)=nalset+1
278  do i=1,nk
279  if(ipoface(i).eq.0) cycle
280  if(inomat(i).ne.3) cycle
281  index=ipoface(i)
282  do
283  if(index.eq.0) exit
284  nalset=nalset+1
285  ialset(nalset)=10*nodface(3,index)+nodface(4,index)
286  index=nodface(5,index)
287  enddo
288  enddo
289  if(istartset(nset).gt.nalset) then
290  iactive(3)=0
291  nset=nset-1
292  else
293  iactive(3)=nset
294  set(nset)(1:22)='ELECTROMAGNETICSZONE3T'
295  do i=23,81
296  set(nset)(i:i)=' '
297  enddo
298  iendset(nset)=nalset
299  endif
300 !
301 ! create contact ties between the domains
302 !
303  do i=1,3
304  do j=1,3
305  if((i.eq.j).or.((i.eq.3).and.(j.eq.2))) cycle
306  if(iactive(i)*iactive(j).gt.0) then
307  ntie=ntie+1
308  tietol(1,ntie)=0.d0
309  tietol(2,ntie)=0.d0
310  write(tieset(1,ntie)(1:1),'(i1)')i
311  write(tieset(1,ntie)(2:2),'(i1)')j
312  do k=3,80
313  tieset(1,ntie)(k:k)=' '
314  enddo
315  tieset(1,ntie)(81:81)='E'
316 !
317 ! slave set
318 ! this set does not have to be identified as nodal
319 ! or facial at this stage
320 !
321  tieset(2,ntie)(1:20)='ELECTROMAGNETICSZONE'
322  write(tieset(2,ntie)(21:21),'(i1)')i
323  do k=22,81
324  tieset(2,ntie)(k:k)=' '
325  enddo
326 !
327 ! master set
328 !
329  tieset(3,ntie)(1:22)='ELECTROMAGNETICSZONE T'
330  write(tieset(3,ntie)(21:21),'(i1)')j
331  do k=23,81
332  tieset(3,ntie)(k:k)=' '
333  enddo
334  endif
335  enddo
336  enddo
337 !
338  return
subroutine isortii(ix, iy, n, kflag)
Definition: isortii.f:6
subroutine nodes(inpc, textpart, co, nk, nk_, set, istartset, iendset, ialset, nset, nset_, nalset, nalset_, istep, istat, n, iline, ipol, inl, ipoinp, inp, ipoinpc)
Definition: nodes.f:22
Hosted by OpenAircraft.com, (Michigan UAV, LLC)