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

Go to the source code of this file.

Functions/Subroutines

subroutine create_iau6 (nef, ipnei, neiel, jq, irow, nzs, iau6, lakonf)
 

Function/Subroutine Documentation

◆ create_iau6()

subroutine create_iau6 ( integer, intent(in)  nef,
integer, dimension(*), intent(in)  ipnei,
integer, dimension(*), intent(in)  neiel,
integer, dimension(*), intent(in)  jq,
integer, dimension(*), intent(in)  irow,
integer, intent(in)  nzs,
integer, dimension(6,*), intent(inout)  iau6,
character*8, dimension(*)  lakonf 
)
20 !
21 ! sets up a field of pointers iau6(j,i) for neighbor j of
22 ! element (cell) i (for CFD-applications) into field auv,aup,aut
23 !
24  implicit none
25 !
26  character*8 lakonf(*)
27 !
28  integer nef,ipnei(*),neiel(*),jq(*),irow(*),nzs,iau6(6,*),
29  & numfaces,id,i,j,iel,indexf
30 !
31  intent(in) nef,ipnei,neiel,jq,irow,nzs
32 !
33  intent(inout) iau6
34 !
35  do i=1,nef
36  indexf=ipnei(i)
37 !
38  do j=1,ipnei(i+1)-ipnei(i)
39  indexf=indexf+1
40  iel=neiel(indexf)
41  if(iel.eq.0) cycle
42  if(i.gt.iel) then
43  call nident(irow(jq(iel)),i,jq(iel+1)-jq(iel),id)
44  iau6(j,i)=jq(iel)+id-1
45  elseif(i.lt.iel) then
46  call nident(irow(jq(i)),iel,jq(i+1)-jq(i),id)
47  iau6(j,i)=nzs+jq(i)+id-1
48  endif
49  enddo
50  enddo
51 !
52  return
subroutine nident(x, px, n, id)
Definition: nident.f:26
Hosted by OpenAircraft.com, (Michigan UAV, LLC)