OB.DAAC Logo
NASA Logo
Ocean Color Science Software

ocssw V2022
put_elements.f
Go to the documentation of this file.
1  subroutine put_elements(lun,orbupd,cdrg,iyr,iday,sec,irec,init)
2 
3 c $Header$
4 c $Log$
5 c
6 c Purpose: This routine writes orbital elements to a file. The updated
7 c elements for the current record and the initial elements
8 c for the following record are written.
9 c
10 c Calling Arguments:
11 c
12 c Name Type I/O Description
13 c -------- ---- --- -----------
14 c lun I*4 I Unit number for elements file
15 c orbupd(6) R*8 I Updated orbital elements
16 c cdrg R*8 I Drag coefficient for ASAP model
17 c iyr I*4 I Year of element epoch
18 c iday I*4 I Day of element epoch
19 c sec R*8 I Seconds of day of element epoch
20 c irec I*4 I Record number of elements in file
21 c init L*4 I Initial or updated elements flag
22 c
23 c By: Frederick S. Patt, GSC, December 23, 1993
24 c
25 c Notes:
26 c
27 c Modification History:
28 c
29 c Modified to write one set of elements per call, with argument specifying
30 c whether this is first or second set in record - F.S. Patt, 10/21/94.
31 c
32 c Modified to not overwrite updated (fitted) elements with initial elements.
33 c F. S. Patt, 1/28/96
34 c
35  implicit none
36 #include "nav_cnst.fin"
37 
38  real*8 orbupd(6),cdrg,sec
39  real*8 orb1(6),orb2(6),secorb,secin
40  integer*4 iyr,iday,iyorb,idorb,idin,irec,lun,spare
41  integer*4 jd,jdr,jdo,i,nrecs
42  logical forb,init
43 
44 c Read header record
45  read(lun,err=990,rec=1) nrecs
46 
47 c If initial elements for record, write integrated elements to file.
48  if (init) then
49 
50 c First check for existing record with fitted elements
51  forb = .false.
52  if (irec.le.nrecs) read (lun,rec=irec,err=990) iyorb,idorb,
53  * secorb,orb1,orb2,cdrg,forb,spare
54 
55 c If record doesn't exist or not previously fitted, write elements
56  if (.not.forb) then
57  do i=1,6
58  orb1(i) = orbupd(i)
59  orb2(i) = 0.d0
60  end do
61  idorb = iday + sec/864.d2
62  iyorb = iyr
63 
64 c Check for year end rollover
65  if (idorb.ge.365) then
66  jdr = jd(iyorb,1,idorb)
67  call jdate(jdr,iyorb,idorb)
68  end if
69  secorb = mod(sec,864.d2)
70  print *,'PUT_ELEMENTS:',irec,iyorb,idorb,secorb
71  print *,orb1,cdrg,forb
72  write (lun,rec=irec,err=991) iyorb,idorb,secorb,orb1,orb2,
73  * cdrg,forb,spare
74 
75 c Check to see if number of records in file needs to be updated
76  if (irec.gt.nrecs) then
77  nrecs = irec
78  write(lun,rec=1,err=991) nrecs
79  end if
80  end if
81 
82 c Else these are updated elements
83  else
84 
85 c Read record irec from file
86  read (lun,rec=irec,err=990) iyorb,idorb,secorb,orb1,orb2,
87  * cdrg,forb,spare
88  idin = iday + sec/864.d2
89  secin = mod(sec,864.d2)
90  jdo = jd(iyorb,1,idorb)
91  jdr = jd(iyr,1,idin)
92  if ((jdo.ne.jdr).or.(secorb.ne.secin)) then
93  print *,'PUT_ELEMENTS: Input time does not match record',
94  * irec,' time'
95  else
96 
97 c Load updated elements and write to file
98  do i=1,6
99  orb2(i) = orbupd(i)
100  end do
101  forb = .true.
102  print *,'PUT_ELEMENTS:',irec,iyorb,idorb,secorb
103  print *,orb2,cdrg,forb
104  write (lun,rec=irec,err=991) iyorb,idorb,secorb,orb1,orb2,
105  * cdrg,forb,spare
106  end if
107  end if
108 
109  return
110 
111  990 print *, 'PUT_ELEMENTS: Error reading elements file'
112  close(lun)
113 
114  return
115  991 print *, 'PUT_ELEMENTS: Error writing elements file'
116  close(lun)
117 
118  return
119  end
int jdate(int32_t julian, int32_t *year, int32_t *doy)
Definition: jdate.c:5
#define real
Definition: DbAlgOcean.cpp:26
Definition: jd.py:1
subroutine put_elements(lun, orbupd, cdrg, iyr, iday, sec, irec, init)
Definition: put_elements.f:2