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

Go to the source code of this file.

Functions/Subroutines

subroutine dflux (flux, sol, kstep, kinc, time, noel, npt, coords, jltyp, temp, press, loadtype, area, vold, co, lakonl, konl, ipompc, nodempc, coefmpc, nmpc, ikmpc, ilmpc, iscale, mi, sti, xstateini, xstate, nstate_, dtime)
 

Function/Subroutine Documentation

◆ dflux()

subroutine dflux ( real*8, dimension(2), intent(out)  flux,
real*8, intent(in)  sol,
integer, intent(in)  kstep,
integer, intent(in)  kinc,
real*8, dimension(2), intent(in)  time,
integer, intent(in)  noel,
integer, intent(in)  npt,
real*8, dimension(3), intent(in)  coords,
integer, intent(in)  jltyp,
real*8, intent(in)  temp,
real*8, intent(in)  press,
character*20, intent(in)  loadtype,
real*8, intent(in)  area,
real*8, dimension(0:mi(2),*), intent(in)  vold,
real*8, dimension(3,*), intent(in)  co,
character*8, intent(in)  lakonl,
integer, dimension(20), intent(in)  konl,
integer, dimension(*), intent(in)  ipompc,
integer, dimension(3,*), intent(in)  nodempc,
real*8, dimension(*), intent(in)  coefmpc,
integer, intent(in)  nmpc,
integer, dimension(*), intent(in)  ikmpc,
integer, dimension(*), intent(in)  ilmpc,
integer, intent(out)  iscale,
integer, dimension(*), intent(in)  mi,
real*8, dimension(6,mi(1),*), intent(in)  sti,
real*8, dimension(nstate_,mi(1),*), intent(in)  xstateini,
real*8, dimension(nstate_,mi(1),*), intent(in)  xstate,
integer, intent(in)  nstate_,
real*8, intent(in)  dtime 
)
23 !
24 ! user subroutine dflux
25 !
26 !
27 ! INPUT:
28 !
29 ! sol current temperature value
30 ! kstep step number
31 ! kinc increment number
32 ! time(1) current step time
33 ! time(2) current total time
34 ! noel element number
35 ! npt integration point number
36 ! coords(1..3) global coordinates of the integration point
37 ! jltyp loading face kode:
38 ! 1 = body flux
39 ! 11 = face 1
40 ! 12 = face 2
41 ! 13 = face 3
42 ! 14 = face 4
43 ! 15 = face 5
44 ! 16 = face 6
45 ! temp currently not used
46 ! press currently not used
47 ! loadtype load type label
48 ! area for surface flux: area covered by the
49 ! integration point
50 ! for body flux: volume covered by the
51 ! integration point
52 ! vold(0..4,1..nk) solution field in all nodes
53 ! 0: temperature
54 ! 1: displacement in global x-direction
55 ! 2: displacement in global y-direction
56 ! 3: displacement in global z-direction
57 ! 4: static pressure
58 ! co(3,1..nk) coordinates of all nodes
59 ! 1: coordinate in global x-direction
60 ! 2: coordinate in global y-direction
61 ! 3: coordinate in global z-direction
62 ! lakonl element label
63 ! konl(1..20) nodes belonging to the element
64 ! ipompc(1..nmpc)) ipompc(i) points to the first term of
65 ! MPC i in field nodempc
66 ! nodempc(1,*) node number of a MPC term
67 ! nodempc(2,*) coordinate direction of a MPC term
68 ! nodempc(3,*) if not 0: points towards the next term
69 ! of the MPC in field nodempc
70 ! if 0: MPC definition is finished
71 ! coefmpc(*) coefficient of a MPC term
72 ! nmpc number of MPC's
73 ! ikmpc(1..nmpc) ordered global degrees of freedom of the MPC's
74 ! the global degree of freedom is
75 ! 8*(node-1)+direction of the dependent term of
76 ! the MPC (direction = 0: temperature;
77 ! 1-3: displacements; 4: static pressure;
78 ! 5-7: rotations)
79 ! ilmpc(1..nmpc) ilmpc(i) is the MPC number corresponding
80 ! to the reference number in ikmpc(i)
81 ! mi(1) max # of integration points per element (max
82 ! over all elements)
83 ! mi(2) max degree of freedomm per node (max over all
84 ! nodes) in fields like v(0:mi(2))...
85 ! sti(i,j,k) actual Cauchy stress component i at integration
86 ! point j in element k. The components are
87 ! in the order xx,yy,zz,xy,xz,yz
88 ! xstateini(i,j,k) value of the state variable i at integration
89 ! point j in element k at the beginning of the
90 ! present increment
91 ! xstateini(i,j,k) value of the state variable i at integration
92 ! point j in element k at the end of the
93 ! present increment
94 ! nstate_ number of state variables
95 ! dtime time length of the increment
96 !
97 !
98 ! OUTPUT:
99 !
100 ! flux(1) magnitude of the flux
101 ! flux(2) not used; please do NOT assign any value
102 ! iscale determines whether the flux has to be
103 ! scaled for increments smaller than the
104 ! step time in static calculations
105 ! 0: no scaling
106 ! 1: scaling (default)
107 !
108  implicit none
109 !
110  character*8 lakonl
111  character*20 loadtype
112 !
113  integer kstep,kinc,noel,npt,jltyp,konl(20),ipompc(*),nstate_,i,
114  & nodempc(3,*),nmpc,ikmpc(*),ilmpc(*),node,idof,id,iscale,mi(*)
115 !
116  real*8 flux(2),time(2),coords(3),sol,temp,press,vold(0:mi(2),*),
117  & area,co(3,*),coefmpc(*),sti(6,mi(1),*),xstate(nstate_,mi(1),*),
118  & xstateini(nstate_,mi(1),*),dtime
119 !
120  intent(in) sol,kstep,kinc,time,noel,npt,coords,
121  & jltyp,temp,press,loadtype,area,vold,co,lakonl,konl,
122  & ipompc,nodempc,coefmpc,nmpc,ikmpc,ilmpc,mi,sti,
123  & xstateini,xstate,nstate_,dtime
124 !
125  intent(out) flux,iscale
126 !
127 ! Start of your own code.
128 !
129 ! example: heat generation due to plasticity for an isotropic
130 ! material
131 !
132 ! the plastic strain tensor is stored in state variables 2...7
133 !
134 !
135  flux(1)=(sti(1,npt,noel)*
136  & (xstate(2,npt,noel)-xstateini(2,npt,noel))+
137  & sti(2,npt,noel)*
138  & (xstate(3,npt,noel)-xstateini(3,npt,noel))+
139  & sti(3,npt,noel)*
140  & (xstate(4,npt,noel)-xstateini(4,npt,noel))+
141  & 2.d0*(sti(4,npt,noel)*
142  & (xstate(5,npt,noel)-xstateini(5,npt,noel))+
143  & sti(5,npt,noel)*
144  & (xstate(6,npt,noel)-xstateini(6,npt,noel))+
145  & sti(6,npt,noel)*
146  & (xstate(7,npt,noel)-xstateini(7,npt,noel))))/
147  & dtime
148 !
149  return
subroutine flux(node1, node2, nodem, nelem, lakon, kon, ipkon, nactdog, identity, ielprop, prop, kflag, v, xflow, f, nodef, idirf, df, cp, R, rho, physcon, g, co, dvi, numf, vold, set, shcon, nshcon, rhcon, nrhcon, ntmat_, mi, ider, ttime, time, iaxial)
Definition: flux.f:24
Hosted by OpenAircraft.com, (Michigan UAV, LLC)