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

Go to the source code of this file.

Functions/Subroutines

subroutine interpolsubmodel (integerglob, doubleglob, value, coords, iselect, nselect, nodeface, tieset, istartset, iendset, ialset, ntie, entity)
 

Function/Subroutine Documentation

◆ interpolsubmodel()

subroutine interpolsubmodel ( integer, dimension(*), intent(in)  integerglob,
real*8, dimension(*), intent(in)  doubleglob,
real*8, dimension(*), intent(inout)  value,
real*8, dimension(3), intent(in)  coords,
integer, dimension(nselect), intent(in)  iselect,
integer, intent(in)  nselect,
integer, intent(in)  nodeface,
character*81, dimension(3,*), intent(in)  tieset,
integer, dimension(*), intent(in)  istartset,
integer, dimension(*), intent(in)  iendset,
integer, dimension(*), intent(in)  ialset,
integer, intent(in)  ntie,
character*1, intent(in)  entity 
)
22 !
23 ! interpolates for a node with coordinates in "coords" the
24 ! "nselect" values with relative positions in "iselect" within the
25 ! global mesh stored in integerglob and doubleglob. The fields
26 ! integerglob and doubleglob are created and filled in
27 ! getglobalresults.c.
28 
29 ! The domain of the global model within which
30 ! the interpolation takes place can be limited to an element
31 ! subset. To this end the submodel to which node "node" belongs is
32 ! determined. The submodels are stored in tieset(1..3,i), i=1..ntie.
33 ! tieset(1,i)(81:81)='S' if the tie is a submodel. In that case the
34 ! set number corresponding to the submodel boundary is stored in
35 ! tieset(2,i) and the set number corresponding to the global element
36 ! model in tieset(3,i). Notice that the submodel boundary can be
37 ! a nodal set or an element face set (the actual node and the actual
38 ! element face are stored in nodeface, respecively).
39 !
40  implicit none
41 !
42  character*1 entity
43  character*81 tieset(3,*)
44 !
45  integer integerglob(*),nselect,iselect(nselect),nodeface,
46  & istartset(*),iendset(*),ialset(*),ntie,i,islavset,iset,
47  & nlength,id,jfaces,nelems,nktet,netet,ne,nkon,nfaces,nfield,
48  & imastset,nterms,konl(20)
49 !
50  real*8 doubleglob(*),value(*),coords(3),ratio(20)
51 !
52  intent(in) integerglob,doubleglob,
53  & coords,iselect,nselect,nodeface,tieset,istartset,iendset,
54  & ialset,ntie,entity
55 !
56  intent(inout) value
57 !
58 ! if no global file was read, set results to zero
59 !
60  if(integerglob(1).eq.0) then
61  do i=1,nselect
62  value(i)=0.d0
63  enddo
64  return
65  endif
66 !
67 ! determining the submodel to which the entity "nodeface" belongs
68 !
69  islavset=0
70  do i=1,ntie
71  if(tieset(1,i)(81:81).ne.'S') cycle
72 !
73 ! check whether submodel is of the right kind (nodal or
74 ! element face)
75 !
76  if(tieset(2,i)(11:11).ne.entity) cycle
77  read(tieset(2,i)(1:10),'(i10)') iset
78  nlength=iendset(iset)-istartset(iset)+1
79  call nident(ialset(istartset(iset)),nodeface,nlength,id)
80  if(id.le.0) cycle
81  if(ialset(istartset(iset)+id-1).ne.nodeface) cycle
82 !
83 ! check whether slave set is of the right
84 !
85  islavset=iset
86  exit
87  enddo
88 !
89 ! check whether a submodel was found
90 !
91  if(islavset.eq.0) then
92  if(entity.eq.'N') then
93  write(*,*) '*ERROR in interpolsubmodel: node',nodeface
94  write(*,*) ' does not belong to any submodel'
95  call exit(201)
96  else
97  nelems=int(nodeface/10)
98  jfaces=nodeface-nelems*10
99  write(*,*) '*ERROR in interpolsubmodel: face',jfaces
100  write(*,*) ' of element',nelems
101  write(*,*) ' does not belong to any submodel'
102  call exit(201)
103  endif
104  endif
105 !
106 ! determining the global element set (if zero: all global elements
107 ! are taken)
108 !
109  read(tieset(3,i)(1:10),'(i10)') imastset
110 !
111 ! reading the number of nodes, tetrahedral interpolation elements,
112 ! global elements, connectivity size and number of faces
113 !
114  nktet=integerglob(1)
115  netet=integerglob(2)
116  ne=integerglob(3)
117  nkon=integerglob(4)
118  nfaces=integerglob(5)
119  nfield=13
120 !
121 ! perform the interpolation
122 !
123  call basis(doubleglob(1),doubleglob(netet+1),
124  & doubleglob(2*netet+1),
125  & doubleglob(3*netet+1),doubleglob(4*netet+1),
126  & doubleglob(5*netet+1),integerglob(6),integerglob(netet+6),
127  & integerglob(2*netet+6),doubleglob(6*netet+1),
128  & integerglob(3*netet+6),nktet,netet,
129  & doubleglob(4*nfaces+6*netet+1),nfield,
130  & doubleglob(13*nktet+4*nfaces+6*netet+1),
131  & integerglob(7*netet+6),integerglob(ne+7*netet+6),
132  & integerglob(2*ne+7*netet+6),integerglob(nkon+2*ne+7*netet+6),
133  & coords(1),coords(2),coords(3),value,ratio,iselect,nselect,
134  & istartset,iendset,ialset,imastset,
135  & integerglob(nkon+2*ne+8*netet+6),nterms,konl)
136 !
137  return
subroutine basis(x, y, z, xo, yo, zo, nx, ny, nz, planfa, ifatet, nktet, netet, field, nfield, cotet, kontyp, ipkon, kon, iparent, xp, yp, zp, value, ratio, iselect, nselect, istartset, iendset, ialset, imastset, ielemnr, nterms, konl)
Definition: basis.f:25
subroutine nident(x, px, n, id)
Definition: nident.f:26
Hosted by OpenAircraft.com, (Michigan UAV, LLC)