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

Go to the source code of this file.

Functions/Subroutines

subroutine calc_ider_cross_split (df, pt1, Tt1, xflow1, xflow2, pt2, Tt2, ichan_num, A1, A2, A_s, dh1, dh2, alpha, zeta_fac, kappa, R, ider, iflag)
 

Function/Subroutine Documentation

◆ calc_ider_cross_split()

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