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

Go to the source code of this file.

Functions/Subroutines

subroutine spcmatch (xboun, nodeboun, ndirboun, nboun, xbounold, nodebounold, ndirbounold, nbounold, ikboun, ilboun, vold, reorder, nreorder, mi)
 

Function/Subroutine Documentation

◆ spcmatch()

subroutine spcmatch ( real*8, dimension(*)  xboun,
integer, dimension(*)  nodeboun,
integer, dimension(*)  ndirboun,
integer  nboun,
real*8, dimension(*)  xbounold,
integer, dimension(*)  nodebounold,
integer, dimension(*)  ndirbounold,
integer  nbounold,
integer, dimension(*)  ikboun,
integer, dimension(*)  ilboun,
real*8, dimension(0:mi(2),*)  vold,
real*8, dimension(*)  reorder,
integer, dimension(*)  nreorder,
integer, dimension(*)  mi 
)
22 !
23 ! matches SPC boundary conditions of one step with those of
24 ! the previous step
25 !
26  implicit none
27 !
28  integer nodeboun(*),ndirboun(*),nboun,nodebounold(*),ilboun(*),
29  & ndirbounold(*),nbounold,i,kflag,idof,id,nreorder(*),ikboun(*),
30  & mi(*)
31 !
32  real*8 xboun(*),xbounold(*),vold(0:mi(2),*),reorder(*)
33 !
34  kflag=2
35 !
36  do i=1,nboun
37  nreorder(i)=0
38  enddo
39 !
40  do i=1,nbounold
41  idof=8*(nodebounold(i)-1)+ndirbounold(i)
42  if(nboun.gt.0) then
43  call nident(ikboun,idof,nboun,id)
44  else
45  id=0
46  endif
47  if((id.gt.0).and.(ikboun(id).eq.idof)) then
48  reorder(ilboun(id))=xbounold(i)
49  nreorder(ilboun(id))=1
50  endif
51  enddo
52 !
53  do i=1,nboun
54  if(nreorder(i).eq.0) then
55  if(ndirboun(i).gt.4) then
56  reorder(i)=0.d0
57  else
58  reorder(i)=vold(ndirboun(i),nodeboun(i))
59  endif
60  endif
61  enddo
62 !
63  do i=1,nboun
64  xbounold(i)=reorder(i)
65  enddo
66 !
67  return
subroutine nident(x, px, n, id)
Definition: nident.f:26
Hosted by OpenAircraft.com, (Michigan UAV, LLC)