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

Go to the source code of this file.

Functions/Subroutines

subroutine storeresidual (nactdof, b, fn, filab, ithermal, nk, sti, stn, ipkon, inum, kon, lakon, ne, mi, orab, ielorien, co, itg, ntg, vold, ielmat, thicke, ielprop, prop)
 

Function/Subroutine Documentation

◆ storeresidual()

subroutine storeresidual ( integer, dimension(0:mi(2),*)  nactdof,
real*8, dimension(*)  b,
real*8, dimension(0:mi(2),*)  fn,
character*87, dimension(*)  filab,
integer, dimension(2)  ithermal,
integer  nk,
real*8, dimension(6,mi(1),*)  sti,
real*8, dimension(6,*)  stn,
integer, dimension(*)  ipkon,
integer, dimension(*)  inum,
integer, dimension(*)  kon,
character*8, dimension(*)  lakon,
integer  ne,
integer, dimension(*)  mi,
real*8, dimension(7,*)  orab,
integer  ielorien,
real*8, dimension(3,*)  co,
integer, dimension(*)  itg,
integer  ntg,
real*8, dimension(0:mi(2),*)  vold,
integer, dimension(mi(3),*)  ielmat,
real*8, dimension(mi(3),*)  thicke,
integer, dimension(*)  ielprop,
real*8, dimension(*)  prop 
)
22 !
23 ! This routine is called in case of divergence:
24 ! stores the residual forces in fn and changes the
25 ! file storage labels so that the independent
26 ! variables (displacements and/or temperatures) and
27 ! the corresponding residual forces are stored in the
28 ! frd file
29 !
30  implicit none
31 !
32  logical force
33 !
34  character*1 cflag
35  character*8 lakon(*)
36  character*87 filab(*)
37 !
38  integer mi(*),nactdof(0:mi(2),*),ithermal(2),i,j,nk,
39  & nfield,ndim,iorienglob,icfdout,ielmat(mi(3),*),
40  & ipkon(*),inum(*),kon(*),ielprop(*),
41  & ne,ielorien,itg(*),ntg,mt,nlabel
42 !
43  real*8 b(*),fn(0:mi(2),*),sti(6,mi(1),*),stn(6,*),orab(7,*),
44  & co(3,*),vold(0:mi(2),*),thicke(mi(3),*),prop(*)
45 !
46  mt=mi(2)+1
47 !
48  nlabel=47
49 !
50 ! storing the residual forces in field fn
51 !
52  do i=1,nk
53  do j=0,mi(2)
54  if(nactdof(j,i).gt.0) then
55  fn(j,i)=b(nactdof(j,i))
56  else
57  fn(j,i)=0.d0
58  endif
59  enddo
60  enddo
61 !
62 ! adapting the storage labels
63 !
64  do i=1,nlabel
65  filab(i)(1:4)=' '
66  enddo
67 !
68  if(ithermal(1).ne.2) then
69  filab(1)(1:3)='U '
70  filab(5)(1:4)='RF '
71  else
72  filab(1)(1:4)=' '
73  filab(5)(1:4)=' '
74  endif
75 !
76  if(ithermal(1).gt.1) then
77  filab(2)(1:4)='NT '
78  filab(10)(1:4)='RFL '
79  filab(14)(1:4)='TT '
80  filab(15)(1:4)='MF '
81  filab(16)(1:4)='TP '
82  filab(17)(1:4)='ST '
83  else
84  filab(2)(1:4)=' '
85  filab(10)(1:4)=' '
86  filab(14)(1:4)=' '
87  filab(15)(1:4)=' '
88  filab(16)(1:4)=' '
89  filab(17)(1:4)=' '
90  endif
91 !
92 ! calculating inum
93 !
94  nfield=0
95  ndim=0
96  iorienglob=0
97  cflag=filab(1)(5:5)
98  icfdout=0
99  force=.false.
100  call extrapolate(sti,stn,ipkon,inum,kon,lakon,nfield,nk,
101  & ne,mi(1),ndim,orab,ielorien,co,iorienglob,cflag,
102  & vold,force,ielmat,thicke,ielprop,prop)
103 !
104  if(ithermal(1).gt.1) then
105  call networkextrapolate(vold,ipkon,inum,kon,lakon,ne,mi)
106  endif
107 !
108 ! interpolation for 1d/2d elements
109 !
110  if(filab(1)(5:5).eq.'I') then
111  nfield=mt
112  cflag=filab(1)(5:5)
113  force=.false.
114  call map3dto1d2d(vold,ipkon,inum,kon,lakon,nfield,nk,
115  & ne,cflag,co,vold,force,mi)
116  endif
117 !
118 ! marking gas nodes by multiplying inum by -1
119 !
120  do i=1,ntg
121  inum(itg(i))=-inum(itg(i))
122  enddo
123 !
124  return
subroutine networkextrapolate(v, ipkon, inum, kon, lakon, ne, mi)
Definition: networkextrapolate.f:21
subroutine extrapolate(yi, yn, ipkon, inum, kon, lakon, nfield, nk, ne, mi, ndim, orab, ielorien, co, iorienloc, cflag, vold, force, ielmat, thicke, ielprop, prop)
Definition: extrapolate.f:22
subroutine map3dto1d2d(yn, ipkon, inum, kon, lakon, nfield, nk, ne, cflag, co, vold, force, mi)
Definition: map3dto1d2d.f:21
Hosted by OpenAircraft.com, (Michigan UAV, LLC)