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

Go to the source code of this file.

Functions/Subroutines

subroutine complete_hel_cyclic (nef, bv, hel, adv, auv, jq, irow, ipnei, neiel, ifatie, c, lakonf, neifa, nzs)
 

Function/Subroutine Documentation

◆ complete_hel_cyclic()

subroutine complete_hel_cyclic ( integer  nef,
real*8, dimension(nef,3)  bv,
real*8, dimension(3,*)  hel,
real*8, dimension(*)  adv,
real*8, dimension(*)  auv,
integer, dimension(*)  jq,
integer, dimension(*)  irow,
integer, dimension(*)  ipnei,
integer, dimension(1)  neiel,
integer, dimension(*)  ifatie,
real*8, dimension(3,3)  c,
character*8, dimension(*)  lakonf,
integer, dimension(*)  neifa,
integer  nzs 
)
27 !
28  implicit none
29 !
30  character*8 lakonf(*)
31 !
32  integer irow(*),nef,nzs,j,k,l,jdof1,jq(*),ifa,neifa(*),
33  & indexf,ipnei(*),neiel(1),ifatie(*),iel,i
34 !
35  real*8 hel(3,*),bv(nef,3),auv(*),adv(*),c(3,3)
36 !
37 ! off-diagonal terms
38 !
39 c$omp parallel default(none)
40 c$omp& shared(nef,ipnei,neiel,neifa,ifatie,hel,auv,bv,c)
41 c$omp& private(i,indexf,iel,ifa,k)
42 c$omp do
43  do i=1,nef
44  do indexf=ipnei(i)+1,ipnei(i+1)
45 !
46  iel=neiel(indexf)
47  if(iel.eq.0) cycle
48  ifa=neifa(indexf)
49 !
50  if(ifatie(ifa).eq.0) then
51  do k=1,3
52  hel(k,i)=hel(k,i)-auv(indexf)*bv(iel,k)
53  enddo
54  elseif(ifatie(ifa).gt.0) then
55 c if(i.gt.iel) then
56  do k=1,3
57  hel(k,i)=hel(k,i)-auv(indexf)*
58  & (c(k,1)*bv(iel,1)+c(k,2)*bv(iel,2)
59  & +c(k,3)*bv(iel,3))
60  enddo
61 c else
62 c do k=1,3
63 c hel(k,i)=hel(k,i)-auv(indexf)*
64 c & (c(1,k)*bv(iel,1)+c(2,k)*bv(iel,2)
65 c & +c(3,k)*bv(iel,3))
66 c enddo
67 c endif
68  else
69 c if(i.gt.iel) then
70  do k=1,3
71  hel(k,i)=hel(k,i)-auv(indexf)*
72  & (c(1,k)*bv(iel,1)+c(2,k)*bv(iel,2)
73  & +c(3,k)*bv(iel,3))
74  enddo
75 c else
76 c do k=1,3
77 c hel(k,i)=hel(k,i)-auv(indexf)*
78 c & (c(k,1)*bv(iel,1)+c(k,2)*bv(iel,2)
79 c & +c(k,3)*bv(iel,3))
80 c enddo
81 c endif
82  endif
83  enddo
84  enddo
85 c$omp end do
86 c$omp end parallel
87 !
88  return
Hosted by OpenAircraft.com, (Michigan UAV, LLC)