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

Go to the source code of this file.

Functions/Subroutines

subroutine generateeminterfaces (istartset, iendset, ialset, iactive, ipkon, lakon, kon, ikmpc, nmpc, nafaces)
 

Function/Subroutine Documentation

◆ generateeminterfaces()

subroutine generateeminterfaces ( integer, dimension(*)  istartset,
integer, dimension(*)  iendset,
integer, dimension(*)  ialset,
integer, dimension(3)  iactive,
integer, dimension(*)  ipkon,
character*8, dimension(*)  lakon,
integer, dimension(*)  kon,
integer, dimension(*)  ikmpc,
integer  nmpc,
integer  nafaces 
)
21 !
22 ! determines the interfaces in between the phi-domain and any
23 ! other domain, i.e.
24 !
25 ! faces belonging to the phi domain and adjacent to the A- or
26 ! A-V-domain
27 ! faces belonging to the A-domain adjacent to the phi-domain
28 ! faces belonging to the A-V-domain adjacent to the phi-domain
29 !
30  implicit none
31 !
32  character*8 lakon(*)
33 !
34  integer istartset(*),iendset(*),ialset(*),iactive(3),
35  & ipkon(*),kon(*),iset,ntotal,nope,nopes,indexe,nodef(9),
36  & ifaceq(8,6),ifacet(6,4),ifacew1(4,5),ifacew2(8,5),m,
37  & node,idof,jdof(3),idummy,nlength,kflag,ikmpc(*),nmpc,id,
38  & iface,nelem,jface,i,j,k,nafaces
39 !
40  data kflag /1/
41 !
42 ! nodes per face for hex elements
43 !
44  data ifaceq /4,3,2,1,11,10,9,12,
45  & 5,6,7,8,13,14,15,16,
46  & 1,2,6,5,9,18,13,17,
47  & 2,3,7,6,10,19,14,18,
48  & 3,4,8,7,11,20,15,19,
49  & 4,1,5,8,12,17,16,20/
50 !
51 ! nodes per face for tet elements
52 !
53  data ifacet /1,3,2,7,6,5,
54  & 1,2,4,5,9,8,
55  & 2,3,4,6,10,9,
56  & 1,4,3,8,10,7/
57 !
58 ! nodes per face for linear wedge elements
59 !
60  data ifacew1 /1,3,2,0,
61  & 4,5,6,0,
62  & 1,2,5,4,
63  & 2,3,6,5,
64  & 3,1,4,6/
65 !
66 ! nodes per face for quadratic wedge elements
67 !
68  data ifacew2 /1,3,2,9,8,7,0,0,
69  & 4,5,6,10,11,12,0,0,
70  & 1,2,5,4,7,14,10,13,
71  & 2,3,6,5,8,15,11,14,
72  & 3,1,4,6,9,13,12,15/
73 !
74 ! The sets iactive(1),iactive(2) and iactive(3) contain all
75 ! external surfaces of the phi-domain, the A-V-domain and the
76 ! A-domain, respectively. In what follows these surfaces are
77 ! shrunk to those faces in common between the phi-domain and the
78 ! others, i.e. iactive(1) contains the faces of the phi-domain
79 ! adjacent to the A-V-domain and A-domain, iactive(2) contains the
80 ! faces of the A-V-domain adjacent to the phi-domain and iactive(3)
81 ! contains the faces of the A-domain adjacent to the phi-domain
82 !
83 ! 1) determining the faces of the phi-domain in all nodes of which a
84 ! MPC for the x-component of A exists.
85 ! 2) determining the faces of the A-V-domain in all nodes of which
86 ! a MPC for phi exists (dof 5)
87 ! 3) determining the faces of the A-domain in all nodes of which
88 ! a MPC for phi exists (dof 5)
89 !
90  jdof(1)=1
91  jdof(2)=5
92  jdof(3)=5
93 !
94  do m=1,3
95  if(iactive(m).gt.0) then
96  iset=iactive(m)
97  ntotal=0
98  loop1: do j=istartset(iset),iendset(iset)
99  iface=ialset(j)
100  nelem=int(iface/10.d0)
101  jface=iface-10*nelem
102 !
103 ! nodes belonging to the element (nope) and to the
104 ! face (nopes)
105 !
106  if(lakon(nelem)(4:5).eq.'8R') then
107  nopes=4
108  nope=8
109  elseif(lakon(nelem)(4:4).eq.'8') then
110  nopes=4
111  nope=8
112  elseif(lakon(nelem)(4:6).eq.'20R') then
113  nopes=8
114  nope=20
115  elseif(lakon(nelem)(4:5).eq.'20') then
116  nopes=8
117  nope=20
118  elseif(lakon(nelem)(4:5).eq.'10') then
119  nopes=6
120  nope=10
121  elseif(lakon(nelem)(4:4).eq.'4') then
122  nopes=3
123  nope=4
124 !
125 ! treatment of wedge faces
126 !
127  elseif(lakon(nelem)(4:4).eq.'6') then
128  nope=6
129  if(jface.le.2) then
130  nopes=3
131  else
132  nopes=4
133  endif
134  elseif(lakon(nelem)(4:5).eq.'15') then
135  nope=15
136  if(jface.le.2) then
137  nopes=6
138  else
139  nopes=8
140  endif
141  endif
142 !
143 ! nodes belonging to the face
144 !
145  indexe=ipkon(nelem)
146  if((nope.eq.20).or.(nope.eq.8)) then
147  do k=1,nopes
148  nodef(k)=kon(indexe+ifaceq(k,jface))
149  enddo
150  elseif((nope.eq.10).or.(nope.eq.4)) then
151  do k=1,nopes
152  nodef(k)=kon(indexe+ifacet(k,jface))
153  enddo
154  elseif(nope.eq.15) then
155  do k=1,nopes
156  nodef(k)=kon(indexe+ifacew2(k,jface))
157  enddo
158  else
159  do k=1,nopes
160  nodef(k)=kon(indexe+ifacew1(k,jface))
161  enddo
162  endif
163 !
164  do i=1,nopes
165  node=nodef(i)
166  idof=8*(node-1)+jdof(m)
167  call nident(ikmpc,idof,nmpc,id)
168  if(id.gt.0) then
169  if(ikmpc(id).eq.idof) cycle
170  endif
171  cycle loop1
172  enddo
173 !
174  ialset(istartset(iset)+ntotal)=iface
175  ntotal=ntotal+1
176  cycle
177  enddo loop1
178  iendset(iset)=istartset(iset)+ntotal-1
179  endif
180  enddo
181 !
182 ! sorting the sets
183 !
184  do m=1,3
185  if(iactive(m).eq.0) cycle
186  nlength=iendset(iactive(m))-istartset(iactive(m))+1
187  call isortii(ialset(istartset(iactive(m))),idummy,
188  & nlength,kflag)
189  enddo
190 !
191 ! nafaces is the total number of faces belonging to the A-V-domain
192 ! or A-domain and adjacent to the phi-domain. The nodes on these
193 ! faces are subject to A.n=0
194 !
195  if((iactive(2).eq.0).and.(iactive(3).eq.0)) then
196  nafaces=0
197  elseif(iactive(2).eq.0) then
198  nafaces=iendset(iactive(3))-istartset(iactive(3))+1
199  elseif(iactive(3).eq.0) then
200  nafaces=iendset(iactive(2))-istartset(iactive(2))+1
201  else
202  nafaces=iendset(iactive(2))-istartset(iactive(2))+1
203  & +iendset(iactive(3))-istartset(iactive(3))+1
204  endif
205 !
206  return
subroutine isortii(ix, iy, n, kflag)
Definition: isortii.f:6
subroutine nident(x, px, n, id)
Definition: nident.f:26
Hosted by OpenAircraft.com, (Michigan UAV, LLC)