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

Go to the source code of this file.

Functions/Subroutines

subroutine extrapolate_u1 (yi, yn, ipkon, inum, kon, lakon, nfield, nk, ne, mi, ndim, orab, ielorien, co, iorienloc, cflag, vold, force, ielmat, thicke, ielprop, prop, i)
 

Function/Subroutine Documentation

◆ extrapolate_u1()

subroutine extrapolate_u1 ( real*8, dimension(ndim,mi(1),*), intent(in)  yi,
real*8, dimension(nfield,*), intent(inout)  yn,
integer, dimension(*), intent(in)  ipkon,
integer, dimension(*), intent(in)  inum,
integer, dimension(*), intent(in)  kon,
character*8, dimension(*), intent(in)  lakon,
integer, intent(in)  nfield,
integer, intent(in)  nk,
integer, intent(in)  ne,
integer, dimension(*), intent(in)  mi,
integer, intent(in)  ndim,
real*8, dimension(7,*), intent(in)  orab,
integer, dimension(mi(3),*), intent(in)  ielorien,
real*8, dimension(3,*), intent(in)  co,
integer, intent(in)  iorienloc,
character*1, intent(in)  cflag,
real*8, dimension(0:mi(2),*), intent(in)  vold,
logical, intent(in)  force,
integer, dimension(mi(3),*), intent(in)  ielmat,
real*8, dimension(mi(3),*), intent(in)  thicke,
integer, dimension(*), intent(in)  ielprop,
real*8, dimension(*), intent(in)  prop,
integer, intent(in)  i 
)
22 !
23 ! extrapolates field values at the integration points to the
24 ! nodes for user element i of type u1
25 !
26 !
27 ! INPUT:
28 !
29 ! yi(ndim,mi(1),*) value of the variables at the integration
30 ! points
31 ! ipkon(i) points to the location in field kon preceding
32 ! the topology of element i
33 ! inum(i) < 0: node i is a network node
34 ! > 0: node i is a structural node
35 ! =0: node i is not used
36 ! kon(*) contains the topology of all elements. The
37 ! topology of element i starts at kon(ipkon(i)+1)
38 ! and continues until all nodes are covered. The
39 ! number of nodes depends on the element label
40 ! lakon(i) contains the label of element i
41 ! nfield number of variables to be extrapolated
42 ! nk maximum node number in the mesh
43 ! ne maximum element number in the mesh
44 ! ndim number of variables in the integration point
45 ! field to be extrapolated
46 ! orab(7,*) description of all local coordinate systems.
47 ! (cf. List of variables and their meaning in the
48 ! User's manual)
49 ! ielorien(i) orientation in element i
50 ! co(1..3,i) global coordinates of node i
51 ! iorienloc 0: extrapolated variables requested in global
52 ! coordinates
53 ! 1: extrapolated variables requested in local
54 ! coordinates
55 ! cflag (char*1) I: interpolate 3D results onto 1D/2D
56 ! E: store extrapolated 1D/2D results
57 ! M: store 1D section forces
58 ! blank: any other case
59 ! vold(j,i) value of variable j in node i at the end
60 ! of the previous iteration
61 ! force logical variable; if true the values to
62 ! be extrapolated are force values; important for
63 ! interpolation from 3D expanded structures on the
64 ! original 1D/2D structure: forces across the
65 ! expansion have to be summed, not interpolated
66 ! ielmat(i) material of element i
67 ! thicke(j,i) thickness of layer j in node i
68 ! ielprop(i) properties for element i are stored in
69 ! prop(ielprop(i)+1),prop(ielprop(i)+2),....
70 ! (number of properties depends on the type of
71 ! element)
72 ! prop property field
73 ! i number of the element for which the extrapolation
74 ! is to be performed
75 !
76 ! OUTPUT:
77 !
78 ! yn(nfield,*) value of the variables at the nodes
79 !
80  implicit none
81 !
82  logical force
83 !
84  character*1 cflag
85  character*8 lakon(*)
86 !
87  integer ipkon(*),inum(*),kon(*),mi(*),ne,nfield,nk,i,ndim,
88  & iorienloc,ielorien(mi(3),*),ielmat(mi(3),*),ielprop(*)
89 !
90  real*8 yi(ndim,mi(1),*),yn(nfield,*),orab(7,*),co(3,*),prop(*),
91  & vold(0:mi(2),*),thicke(mi(3),*)
92 !
93  intent(in) yi,ipkon,inum,kon,lakon,nfield,nk,
94  & ne,mi,ndim,orab,ielorien,co,iorienloc,cflag,
95  & vold,force,ielmat,thicke,ielprop,prop,i
96 !
97  intent(inout) yn
98 !
99 ! START OF THIS SUBROUTINE
100 !
101  integer indexe,j,k,node
102 !
103  if(iorienloc.ne.0) then
104  write(*,*) '*ERROR in extrapolate_u1'
105  write(*,*) ' no local orientation for variables'
106  write(*,*) ' belonging to this type of element'
107  write(*,*) ' allowed'
108  call exit(201)
109  endif
110 !
111  if(nfield.eq.6) then
112  indexe=ipkon(i)
113  do j=1,2
114  node=kon(indexe+j)
115  do k=1,nfield
116  yn(k,node)=yn(k,node)+yi(k,1,i)
117  enddo
118  enddo
119  else
120  write(*,*) '*ERROR in extrapolate_u1'
121  write(*,*) ' extropolation for element of type u1'
122  write(*,*) ' is only coded for fields with 6'
123  write(*,*) ' entries'
124  call exit(201)
125  endif
126 !
127 ! END OF THIS SUBROUTINE
128 !
129  return
Hosted by OpenAircraft.com, (Michigan UAV, LLC)