OB.DAAC Logo
NASA Logo
Ocean Color Science Software

ocssw V2022
add_elements.f
Go to the documentation of this file.
1  subroutine add_elements(orbinit,cdrg,iyinit,idinit,secinit,
2  * igyr,igday,secst,tdifmax,irec)
3 
4 c $Header: /app/shared/RCS/irix-5.2/seawifsd/src/mops/mopsV4.1/swfnav/add_elements.f,v 1.1 1995 /01/17 23:01:58 seawifsd Exp seawifsd $
5 c $Log: add_elements.f,v $
6 c Revision 1.1 1995/01/17 23:01:58 seawifsd
7 c Initial revision
8 c
9 
10 c Calling Arguments:
11 c
12 c Name Type I/O Description
13 c -------- ---- --- -----------
14 c orbinit(6) R*8 I/O Orbital elements
15 c cdrg R*8 O Drag coefficient for ASAP model
16 c iyinit I*4 I/O Year of element epoch
17 c idinit I*4 I/O Day of element epoch
18 c secinit R*8 I/O Seconds of day of element epoch
19 c igyr I*4 I Year of data start time
20 c igday I*4 I Day of data start time
21 c secst R*8 I Seconds of data start time
22 c tdifmax R*8 I Maximum allowed time between element epoch
23 c and data start time
24 c irec I*4 I/O Record number for elements in file
25 c
26 c By: Frederick S. Patt, GSC, October 20, 1994
27 c
28 c Notes:
29 c
30 c Modification History:
31 c
32 c Added conditional compile for Sun OS record length differences. B. A.
33 c Franz, GSC, November 14, 1997.
34 c
35 
36  real*8 orbinit(6),secinit,secst,cdrg,ge,aj2,tdifmax
37  real*8 asap(6,1500),tsap(1500),tdif,orbupd(6)
38  integer*4 igyr,igday,jd,i,j,ilast,lun,jdl
39  integer*4 iyinit,idinit,irec,iddif,nstp,li,mi,iplt,ier
40  character*80 filnm
41  logical*4 init
42  data ge/3.9860050d5/,aj2/0.10826270e-02/
43  data li/30/,mi/15/,iplt/0/,lun/17/
44 
45  iddif = jd(igyr,1,igday) - jd(iyinit,1,idinit)
46  tdif = iddif*864.d2 + secst - secinit
47  print *,'ADD_ELEMENTS: time difference',tdif
48 
49 c Open elements file
50  filnm = '$ELEMENTS/elements.dat'
51  call filenv(filnm,filnm)
52  open(lun,file=filnm,status='old',access='direct',err=990,
53  * recl=128)
54 
55 c Continue while time difference between data start time and element epoch
56 c exceeds maximum allowed
57  dowhile(tdif.gt.tdifmax)
58 
59 c Integrate orbit using current elements
60 
61  nstp = tdif/60.d0 + 1
62  if (nstp.gt.1500) nstp = 1500
63  call asaps(li,mi,iplt,orbinit,iyinit,idinit,secinit,nstp,cdrg,
64  * tsap,asap)
65 
66 c Write any placeholder records at ascending node crossings
67  init = .true.
68  i = 2
69  dowhile(i.le.nstp)
70  if ((asap(3,i).gt.0.).and.(asap(3,i-1)).lt.0.) then
71  irec = irec + 1
72  ier = 1
73  i = i - 1
74  dowhile((ier.ne.0.).and.(i.le.nstp))
75  i = i + 1
76  call vec2mean(asap(1,i),ge,aj2,orbinit,orbupd,ier)
77  end do
78  call put_elements(lun,orbupd,cdrg,iyinit,idinit,tsap(i),
79  * irec,init)
80  ilast = i
81  end if
82  i = i + 1
83  end do
84 
85 c Update current elements and epoch (keep eccentricity from initial elements
86  orbinit(1) = orbupd(1)
87  do j=3,6
88  orbinit(j) = orbupd(j)
89  end do
90  idinit = idinit + tsap(ilast)/864.d2
91 
92 c Check for year end rollover
93  if (idinit.ge.365) then
94  jdl = jd(iyinit,1,idinit)
95  call jdate(jdl,iyinit,idinit)
96  end if
97  secinit = mod(tsap(ilast),864.d2)
98  iddif = jd(igyr,1,igday) - jd(iyinit,1,idinit)
99  tdif = iddif*864.d2 + secst - secinit
100 
101  end do
102  close(lun)
103  return
104 
105  990 print *, 'PUT_ELEMENTS: Error opening elements file'
106  close(lun)
107  return
108 
109  end
int jdate(int32_t julian, int32_t *year, int32_t *doy)
Definition: jdate.c:5
#define real
Definition: DbAlgOcean.cpp:26
subroutine add_elements(orbinit, cdrg, iyinit, idinit, secinit, igyr, igday, secst, tdifmax, irec)
Definition: add_elements.f:3
subroutine vec2mean(vec, ge, aj2, xinit, xmean, ier)
Definition: vec2mean.f:2
subroutine filenv(infil, outfil)
Definition: filenv.f:2
subroutine asaps(LI, MI, IPLT, ORBIN, IYR, IDAY, SEC, NSTP, CDRG, TVEC, XVEC)
Definition: asaps.f:3
Definition: jd.py:1
subroutine put_elements(lun, orbupd, cdrg, iyr, iday, sec, irec, init)
Definition: put_elements.f:2