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

Go to the source code of this file.

Functions/Subroutines

subroutine applympc (nface, ielfa, is, ie, ifabou, ipompc, vfa, coefmpc, nodempc, ipnei, neifa, labmpc, xbounact, nactdoh, ifaext, nfaext)
 

Function/Subroutine Documentation

◆ applympc()

subroutine applympc ( integer, intent(in)  nface,
integer, dimension(4,*), intent(in)  ielfa,
integer, intent(in)  is,
integer, intent(in)  ie,
integer, dimension(*), intent(in)  ifabou,
integer, dimension(*), intent(in)  ipompc,
real*8, dimension(0:7,*), intent(inout)  vfa,
real*8, dimension(*), intent(in)  coefmpc,
integer, dimension(3,*), intent(in)  nodempc,
integer, dimension(*), intent(in)  ipnei,
integer, dimension(*), intent(in)  neifa,
character*20, dimension(*), intent(in)  labmpc,
real*8, dimension(*), intent(in)  xbounact,
integer, dimension(*), intent(in)  nactdoh,
integer, dimension(*), intent(in)  ifaext,
integer, intent(in)  nfaext 
)
21 !
22 ! applies MPC's to the faces
23 !
24  implicit none
25 !
26  character*20 labmpc(*)
27 !
28  integer i,j,nface,ielfa(4,*),ipointer,is,ie,ifabou(*),mpc,
29  & ipompc(*),index,iel,iface,nodempc(3,*),ipnei(*),neifa(*),
30  & nactdoh(*),ielorig,ifaext(*),nfaext,k
31 !
32  real*8 coefmpc(*),denominator,vfa(0:7,*),sum,xbounact(*),
33  & coefnorm
34 !
35  intent(in) nface,ielfa,is,ie,ifabou,ipompc,coefmpc,
36  & nodempc,ipnei,neifa,labmpc,xbounact,nactdoh,ifaext,nfaext
37 !
38  intent(inout) vfa
39 !
40 c$omp parallel default(none)
41 c$omp& shared(nfaext,ifaext,ielfa,ifabou,ipompc,nodempc,coefmpc,
42 c$omp& xbounact,nactdoh,neifa,ipnei,vfa,is,ie)
43 c$omp& private(k,i,ipointer,mpc,index,sum,coefnorm,ielorig,iel,iface)
44 c$omp do
45  do k=1,nfaext
46  i=ifaext(k)
47  if(ielfa(2,i).ge.0) cycle
48  ipointer=-ielfa(2,i)
49  do j=is,ie
50  if(ifabou(ipointer+j).ge.0) cycle
51  mpc=-ifabou(ipointer+j)
52 !
53  index=ipompc(mpc)
54  sum=0.d0
55  coefnorm=0.d0
56  do
57  if(index.eq.0) exit
58  if(nodempc(1,index).lt.0) then
59 !
60 ! a negative number refers to a boundary
61 ! condition (fields nodeboun, ndirboun..)
62 ! resulting from a SPC in local coordinates
63 !
64  sum=sum+coefmpc(index)*xbounact(-nodempc(1,index))
65  else
66 !
67 ! face term
68 !
69  ielorig=int(nodempc(1,index)/10.d0)
70  iel=nactdoh(ielorig)
71  iface=nodempc(1,index)-10*ielorig
72  sum=sum+coefmpc(index)
73  & *vfa(nodempc(2,index),neifa(ipnei(iel)+iface))
74  coefnorm=coefnorm+coefmpc(index)**2
75  endif
76  index=nodempc(3,index)
77  enddo
78 !
79 ! distribute the sum across all terms which are not
80 ! fixed by boundary conditions
81 !
82  index=ipompc(mpc)
83  do
84  if(index.eq.0) exit
85  if(nodempc(1,index).gt.0) then
86  ielorig=int(nodempc(1,index)/10.d0)
87  iel=nactdoh(ielorig)
88  iface=nodempc(1,index)-10*ielorig
89  vfa(nodempc(2,index),neifa(ipnei(iel)+iface))=
90  & vfa(nodempc(2,index),neifa(ipnei(iel)+iface))-
91  & sum*coefmpc(index)/coefnorm
92  endif
93  index=nodempc(3,index)
94  enddo
95  enddo
96  enddo
97 c$omp end do
98 c$omp end parallel
99 !
100  return
Hosted by OpenAircraft.com, (Michigan UAV, LLC)