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

Go to the source code of this file.

Functions/Subroutines

subroutine calc_ider_tee (df, pt1, Tt1, xflow1, xflow2, pt2, Tt2, A1, A2, zeta_fac, kappa, R, ider, iflag)
 

Function/Subroutine Documentation

◆ calc_ider_tee()

subroutine calc_ider_tee ( real*8, dimension(6)  df,
real*8  pt1,
real*8  Tt1,
real*8  xflow1,
real*8  xflow2,
real*8  pt2,
real*8  Tt2,
real*8  A1,
real*8  A2,
real*8  zeta_fac,
real*8  kappa,
real*8  R,
integer  ider,
integer  iflag 
)
22 !
23  implicit none
24 !
25  integer ider,iflag
26 !
27  real*8
28  &df(6),
29  &pt1,
30  &pt2,
31  &tt1,
32  &tt2,
33  &xflow1,
34  &xflow2,
35  &a1,
36  &a2,
37  &kappa,
38  &r,
40 ! Accuracy of the numerical derivatives
41  &eps,
42 ! Step for the derivative
43  &h,
44  &f0,
45  &zeta_fac
46 !
47  eps = 1.0e-8
48 !
49  f0 = calc_residual_tee(pt1,tt1,xflow1,xflow2,pt2,
50  &tt2,a1,a2,zeta_fac,kappa,r,ider,iflag)
51 !
52  h = eps*dabs(pt1)
53  if(h.eq.0)then
54  h = eps
55  endif
56  df(1) = (calc_residual_tee(pt1+h,tt1,xflow1,xflow2,pt2,
57  &tt2,a1,a2,zeta_fac,kappa,r,ider,iflag)-f0)/h
58 !
59  h = eps*dabs(tt1)
60  if(h.eq.0)then
61  h = eps
62  endif
63  df(2) = (calc_residual_tee(pt1,tt1+h,xflow1,xflow2,pt2,
64  &tt2,a1,a2,zeta_fac,kappa,r,ider,iflag)-f0)/h
65 !
66  h = eps*dabs(xflow1)
67  if(h.eq.0)then
68  h = eps
69  endif
70  df(3) = (calc_residual_tee(pt1,tt1,xflow1+h,xflow2,pt2,
71  &tt2,a1,a2,zeta_fac,kappa,r,ider,iflag)-f0)/h
72 !
73  h = eps*dabs(xflow2)
74  if(h.eq.0)then
75  h = eps
76  endif
77  df(4) = (calc_residual_tee(pt1,tt1,xflow1,xflow2+h,pt2,
78  &tt2,a1,a2,zeta_fac,kappa,r,ider,iflag)-f0)/h
79 !
80  h = eps*dabs(pt2)
81  if(h.eq.0)then
82  h = eps
83  endif
84  df(5) = (calc_residual_tee(pt1,tt1,xflow1,xflow2,pt2+h,
85  &tt2,a1,a2,zeta_fac,kappa,r,ider,iflag)-f0)/h
86 !
87  h = eps*dabs(tt2)
88  if(h.eq.0)then
89  h = eps
90  endif
91  df(6) = (calc_residual_tee(pt1,tt1,xflow1,xflow2,pt2,
92  &tt2+h,a1,a2,zeta_fac,kappa,r,ider,iflag)-f0)/h
93 !
94  return
subroutine df(x, u, uprime, rpar, nev)
Definition: subspace.f:133
real *8 function calc_residual_tee(pt1, Tt1, xflow1, xflow2, pt2, Tt2, A1, A2, zeta_fac, kappa, R, ider, iflag)
Definition: calc_residual_tee.f:20
Hosted by OpenAircraft.com, (Michigan UAV, LLC)