OB.DAAC Logo
NASA Logo
Ocean Color Science Software

ocssw V2022
get_elements.f
Go to the documentation of this file.
1  subroutine get_elements(iyr,iday,sec,secend,fit,orbout,cdrg,
2  * iyorb,idorb,secorb,irec,ier)
3 
4 c $Header$
5 c $Log$
6 c
7 c Purpose: This routine reads orbital elements for use by the ASAP orbit
8 c integrator from a file, according to the specified time.
9 c The first element set is read which precedes the specified
10 c time. The file contains both fitted and unfitted (to GPS
11 c data) element sets.
12 c
13 c Calling Arguments:
14 c
15 c Name Type I/O Description
16 c -------- ---- --- -----------
17 c iyr I*4 I Year of date for requested elements
18 c iday I*4 I Day of date for requested elements
19 c sec R*8 I Seconds of day for requested elements
20 c secend R*8 I Seconds of day at end of interval
21 c fit L*4 I/O Fitting flag for elements
22 c =.true., use fitted elements if available
23 c =.false., use unfitted elements only
24 c output value reflects elements used.
25 c orbout(6) R*8 O Orbital elements
26 c cdrg R*8 O Drag coefficient for ASAP model
27 c iyorb I*4 O Year of element epoch
28 c idorb I*4 O Day of element epoch
29 c secorb R*8 O Seconds of day of element epoch
30 c irec I*4 O Record number of elements in file
31 c ier I*4 O Error code: =0, success
32 c =1, error
33 c
34 c By: Frederick S. Patt, GSC, December 23, 1993
35 c
36 c Notes:
37 c
38 c Modification History:
39 c
40 c Added check for non-fitted elements within time span of data, which will
41 c cause fit flag to be set to false. F.S. Patt, GSC, October 26, 1994.
42 c
43 c Fixed a typo in the comparison of orbital element epochs to input times.
44 c F.S. Patt, GSC, August 15, 1996.
45 c
46 c Added conditional compile for Sun OS record length differences. B. A.
47 c Franz, GSC, November 14, 1997.
48 c
49 c Added record number in direct access read for Sun OS compatibility.
50 c B. A. Franz, November 14, 1997.
51 c
52 c Added check for last record in file more than one orbit before end
53 c of data. F.S. Patt, SAIC GSC, February 12, 1998
54 c
55 c Fixed a bug introduced in the previous change, which caused the fit
56 c flag in the elements.dat record not to be checked if only the last record
57 c was read. F.S. Patt, SAIC GSC, March 10, 1998.
58 
59  implicit none
60 #include "nav_cnst.fin"
61 
62  real*8 orbout(6),orb1(6),orb2(6),cdrg,sec,secend
63  real*8 secorb,tdif,s,cd,tdifen,oneorb
64  integer*4 iyr,iday,iyorb,idorb,irec,ier,lun,spare
65  integer*4 jd,jdr,jdo,i,nrecs,iy,id
66  character*256 filnm
67  logical fit,forb
68  data oneorb/5940.d0/
69  data lun/17/
70 
71 c Open elements file
72  filnm = '$ELEMENTS/elements.dat'
73  call filenv(filnm,filnm)
74  open(lun,file=filnm,status='old',access='direct',err=990,
75  * recl=128, convert='big_endian')
76 
77 c Read header record
78  read(lun,rec=1,err=990) nrecs
79 
80  jdr = jd(iyr,1,iday)
81 
82 c Read last record and check if more than one orbit before end of data
83  irec = nrecs
84  write(*,*) irec
85  read (lun,rec=irec,err=990) iy,id,s,orb1,orb2,cd,forb,spare
86  write(*,*) iy,id,s,orb1,orb2,cd,forb
87  jdo = jd(iy,1,id)
88  tdif = (jdo-jdr)*864.d2 + s - sec
89  tdifen = tdif + sec - secend
90  if (tdifen.lt.(-oneorb)) fit = .false.
91 
92 c Read backwards through file to find most recent elements set which
93 c precedes start of data
94  dowhile(tdif.gt.0.)
95 
96  irec = irec - 1
97  if (irec.le.1) then
98  print *,'GET_ELEMENTS:',
99  *' No elements in file preceding requested time'
100  go to 990
101  end if
102  read (lun,rec=irec,err=990) iy,id,s,orb1,orb2,cd,forb,spare
103  jdo = jd(iy,1,id)
104  tdif = (jdo-jdr)*864.d2 + s - sec
105 
106 c Check if unfitted elements in interval
107  tdifen = tdif + sec - secend
108  if (tdifen.lt.0.) fit = (fit.and.forb)
109  end do
110 
111 c Set fit flag using input and record flags
112  fit = (fit.and.forb)
113 
114 c Select element set from record according to logical flags
115  if (fit) then
116  do i=1,6
117  orbout(i) = orb2(i)
118  end do
119  else
120  do i=1,6
121  orbout(i) = orb1(i)
122  end do
123  end if
124  iyorb = iy
125  idorb = id
126  secorb = s
127  cdrg = cd
128  print *,'GET_ELEMENTS:',iyorb,idorb,secorb,fit
129  print *,orbout
130  ier = 0
131  close(lun)
132 
133  return
134 
135  990 print *, 'Error reading elements file'
136  ier = 1
137  close(lun)
138 
139  return
140  end
#define real
Definition: DbAlgOcean.cpp:26
subroutine filenv(infil, outfil)
Definition: filenv.f:2
Definition: jd.py:1
subroutine get_elements(iyr, iday, sec, secend, fit, orbout, cdrg, iyorb, idorb, secorb, irec, ier)
Definition: get_elements.f:3