Due to the lapse in federal government funding, NASA is not updating this website. We sincerely regret this inconvenience.
NASA Logo
Ocean Color Science Software

ocssw V2022
attcmp.f
Go to the documentation of this file.
1  subroutine attcmp(nlines, gaclac, tlm, pos, vel, timref,
2  1 time, navctl, navqc, attxfm, navblk, tiltpr, tiltfl)
3 c
4 c attcmp(nlines, gaclac, tlm, pos, vel, timref, time, navctl,
5 c navqc, attxfm, navblk)
6 c
7 c Purpose: Compute attitude from the spacecraft sensor data
8 c
9 c Calling Arguments:
10 c
11 c Name Type I/O Description
12 c -------- ---- --- -----------
13 c nlines I*4 I number of scan lines to process
14 c gaclac I*4 I flag for GAC or LAC data type,
15 c tlm struct I telemetry block data
16 c pos R*4 I size 3 by nlines orbit position in km
17 c vel R*4 I size 3 by nlines orbit velocity in km/sec
18 c timref R*8 I size 3 reference time at start line
19 c of data: year, day, sec
20 c time R*8 I array of time in seconds relative to
21 c timref for every scan line
22 c navctl struct I controls for processing and sensor,
23 c offset data
24 c navqc struct I navigation quality control parameters
25 c attxfm R*4 I size 3 by 3 by nlines attitude transform
26 c matricies
27 c navblk struct O navigation block containing navigation
28 c information
29 c tiltpr R*4 O size nlines processed tilt data
30 c for each line and for two tilt measurements
31 c tiltfl I*4 O flags for goodness of tilts 0 - good,
32 c 1 - bad
33 c
34 c By: W. Robinson, GSC, 1 Apr 93
35 c
36 c Notes:
37 c
38 c Modification History:
39 c
40 c Added rs to call to l_sun. F.S. Patt, GSC, August 15, 1996.
41 c
42 c Modified to compute sen_mat and scan_ell using zero tilt if tilt flags
43 c are set. F. S. Patt, GSC, October 9, 1996.
44 c
45 c Minor modifications for compatibility with Sun OS. B. A. Franz, GSC,
46 c November 14, 1997.
47 c
48 c Added initial index to some array references for Sun OS compatibility.
49 c B. A. Franz, GSC, November 14, 1997.
50 c
51 c Added check for large yaw angle (indicating large uncertainty) and set
52 c nav flag. F. S. Patt, SAIC GSC, December 11, 1997.
53 c
54 c Modified to work with new (Kalman filter) version of attdet.f.
55 c F. S. Patt, SAIC GSC, June 2, 1998.
56 c
57 c Modified to add two-pass (forward and backward) filtering.
58 c F. S. Patt, SAIC GSC, June 4, 1998.
59 c
60 c Fixed bug which caused incorrect initialization for flagged frames at
61 c start of scene. F. S. Patt, SAIC GSC, June 16, 1998.
62 c
63 c Added test for final yaw covariance greater than a limit, indicating
64 c poor geometry. F. S. Patt, SAIC GSC, September 24, 1998.
65 c
66 c Fixed bug which caused the last scan line not to be checked for covariance
67 c greater that the limit. F. S. Patt, SAIC GSC, December 10, 1998.
68 c
69 c Fixed to compute nadbodfl from widphfl before writing to file.
70 c F. S. Patt, SAIC GSC, February 5, 1999.
71 c
72 c Modified to set nflag(8) (navwarn) for nflag(5) = 1 (high yaw uncertainty)
73 c or nflag(7) = 2 (tilt change). F. S. Patt, SAIC, October 10, 2001.
74 c
75 c Modified to recompute yaw covariance from forward and backward passes
76 c before checking against limit. F. S. Patt, SAIC, March 21, 2003.
77 c
78 c Modified check on input flags to check for valid orbit vector only, in order
79 c to allow attitude determination even with flagged time (e.g., time shift).
80 c F. S. Patt, SAIC, October 30, 2007.
81 
82 
83  implicit none
84 #include "tlm_str.fin"
85 #include "navctl_s.fin"
86 #include "navqc_s.fin"
87 #include "navblk_s.fin"
88 c
89  type(tlm_struct) tlm
90  type(navctl_struct) navctl
91  type(navqc_struct) navqc
92  type(navblk_struct) navblk(maxlin)
93  type(navblk_struct) navblk1
94 c
95  integer*4 nlines, gaclac, isens
96  real*4 pos(3,nlines), vel(3,nlines), attxfm(3,3,nlines)
97  real*8 timref(3), time(nlines)
98 c
99  integer*4 ilin, i, l, sunbodfl(maxlin), nadbodfl(maxlin), year,
100  1 day, delday, tiltfl(maxlin), widphfl(maxlin,2)
101  real*4 sun_bod(3,maxlin), nad_bod(3,maxlin), tiltpr(maxlin),
102  1 widphse(2,maxlin,2), yaw_lim, covar(3), delta_t,
103  1 att_prev(3), covpass1(3,maxlin), maxcov
104  real*8 sec, rs
105  logical first
106  data yaw_lim/6.0/, att_prev/3*0.0/, maxcov/1.0e-6/
107 c
108 c start, compute the array of spacecraft sun vectors for each line
109 c
110  call suncomp( tlm, time, nlines, gaclac, navqc, navctl,
111  1 sun_bod, sunbodfl )
112 c
113 c debug printout of sun body vectors
114  if (navctl%lvdbug.gt.1) then
115  open(unit=21,file='sun_bod.out',status = 'unknown')
116  do ilin = 1, nlines
117  write(21,1200)(sun_bod(isens,ilin),isens=1,3),
118  1 sunbodfl(ilin)
119  1200 format(2x,3f12.7,i7)
120  end do
121  close(unit = 21)
122  end if
123 c
124 c compute the array of spacecraft nadir vectors for each line
125 c
126  call earcomp(tlm, time, pos, vel, nlines, gaclac, navqc,
127  1 navctl, nad_bod, nadbodfl, widphse, widphfl )
128 c
129 c prepare the tilt values reported
130 c
131  call tiltcomp( nlines, tlm, timref, time, gaclac, navqc,
132  1 tiltpr, tiltfl )
133 c
134 c (future) compute the array of spacecraft magnetic vectors
135 c
136 c call magcomp(
137 c
138 c Loop through lines and compute the attitude
139 c First pass of Kalman filter
140 c
141  first = .true.
142  do ilin = 1, nlines
143 c
144 c compute the model sun vector
145 c (get actual year, day, sec from ref time and time array)
146 c
147  year = timref(1)
148  day = timref(2)
149  sec = timref(3)
150  delday = 0
151  call ydsadd( year, day, sec, delday, time(ilin) )
152 c
153  call l_sun( year, day, sec, navblk(ilin)%sun_ref(1), rs )
154 c
155 c access the local vertical from the attitude transform matrix
156 c
157  navblk(ilin)%l_vert(1) = attxfm(1,1,ilin)
158  navblk(ilin)%l_vert(2) = attxfm(1,2,ilin)
159  navblk(ilin)%l_vert(3) = attxfm(1,3,ilin)
160 c
161 c determine attitude from the model and measured vectors
162 c
163  if( (widphfl(ilin,1).ne.0) .and. (widphfl(ilin,2).ne.0)) then
164  navblk(ilin)%nflag(4) = 1
165  end if
166 c
167  if( sunbodfl(ilin) .ne. 0 ) then
168  navblk(ilin)%nflag(3) = 1
169  end if
170 c
171  if( navblk(ilin)%nflag(2) .eq. 0 ) then
172 c
173  if (first) then
174  first = .false.
175  delta_t = 0.0
176  else
177  delta_t = time(ilin) - time(ilin-1)
178  end if
179 
180  call attdet(navblk(ilin), navctl, nad_bod(1,ilin),
181  1 attxfm(1,1,ilin), pos(1,ilin), vel(1,ilin),
182  1 widphse(1,ilin,1), widphfl(ilin,1), widphse(1,ilin,2),
183  1 widphfl(ilin,2), sun_bod(1,ilin), sunbodfl(ilin),
184  1 att_prev, delta_t, covar)
185 
186 c Store attitude for next line and save covariance for second pass
187  do i=1,3
188  att_prev(i) = navblk(ilin)%att_ang(i)
189  covpass1(i,ilin) = covar(i)
190  end do
191 
192  end if
193  end do
194 
195 c Check last scan for yaw angle or covariance outside of limit
196  if ((covar(1).ge.maxcov).or.
197  1 (abs(navblk(nlines)%att_ang(1)).gt.yaw_lim)) then
198  navblk(nlines)%nflag(8) = 1
199  navblk(nlines)%nflag(5) = 1
200  end if
201 
202 c Now make a second pass backwards
203  do l=1,nlines-1
204  ilin = nlines - l
205 
206 c Already have sun_ref and flags from first pass
207 
208  if( navblk(ilin)%nflag(2) .eq. 0 ) then
209 
210  do i=1,3
211  navblk1%sun_ref(i) = navblk(ilin)%sun_ref(i)
212  end do
213  delta_t = time(ilin) - time(ilin+1)
214 
215  call attdet(navblk1, navctl, nad_bod(1,ilin),
216  1 attxfm(1,1,ilin), pos(1,ilin), vel(1,ilin),
217  1 widphse(1,ilin,1), widphfl(ilin,1), widphse(1,ilin,2),
218  1 widphfl(ilin,2), sun_bod(1,ilin), sunbodfl(ilin),
219  1 att_prev, delta_t, covar)
220 
221  do i=1,3
222  att_prev(i) = navblk1%att_ang(i)
223  end do
224 
225 c Compute weighted average of attitude from two passes
226  do i=1,3
227  navblk(ilin)%att_ang(i) =
228  1 ( navblk(ilin)%att_ang(i)*covar(i) +
229  1 navblk1%att_ang(i)*covpass1(i,ilin)) /
230  1 (covar(i) + covpass1(i,ilin))
231  covpass1(i,ilin) = 1./(1./covpass1(i,ilin) +
232  1 1./covar(i))
233  end do
234 
235 c write (65,*) (covpass1(i,ilin),i=1,3)
236 
237 c Check for yaw angle or covariance outside of limit
238  if ((covpass1(1,ilin).ge.maxcov).or.
239  1 (abs(navblk(ilin)%att_ang(1)).gt.yaw_lim)) then
240  navblk(ilin)%nflag(8) = 1
241  navblk(ilin)%nflag(5) = 1
242  end if
243 
244  end if
245  end do
246 c
247 c debug printout of nadir body vectors
248  if (navctl%lvdbug.gt.1) then
249  open(unit=21,file='nad_bod.out',status = 'unknown')
250  do ilin = 1, nlines
251  if ((widphfl(ilin,1).ne.0) .and. (widphfl(ilin,2).ne.0)) then
252  nadbodfl(ilin) = 1
253  else
254  nadbodfl(ilin) = 0
255  end if
256  write(21,1200) (nad_bod(isens,ilin),isens=1,3),
257  1 nadbodfl(ilin)
258  end do
259  close(unit = 21)
260  end if
261 c
262 c add scan ellipse coefficients and tilt flags to the nav block
263 c
264  do ilin = 1, nlines
265 
266  navblk(ilin)%nflag(7) = tiltfl(ilin)
267  if ( tiltfl(ilin) .eq. 1) then
268  navblk(ilin)%nflag(1) = 1
269  else if ( tiltfl(ilin) .eq. 2) then
270  navblk(ilin)%nflag(8) = 1
271  end if
272 
273  call ellxfm( attxfm(1,1,ilin), navblk(ilin)%att_ang(1),
274  1 tiltpr(ilin), pos(1,ilin), navctl,
275  1 navblk(ilin)%sen_mat(1,1), navblk(ilin)%scan_ell(1) )
276 
277  end do
278 c
279 c and return
280 c
281  return
282  end
subroutine suncomp(tlm, time, nlines, gaclac, navqc,
Definition: suncomp.f:2
subroutine attdet(navblk, navctl, nad_bod, attxfm,
Definition: attdet.f:2
Definition: l_sun.py:1
#define real
Definition: DbAlgOcean.cpp:26
subroutine ydsadd(iy, id, sec, deld, delsec)
Definition: ydsadd.f:2
subroutine earcomp(tlm, time, pos, vel, nlines, gaclac, navqc,
Definition: earcomp.f:2
subroutine tiltcomp(nlines, tlm, timref, time, gaclac, navqc, tiltpr, tiltfl)
Definition: tiltcomp.f:3
subroutine ellxfm(attxfm, att_ang, tilt, p, navctl, smat, coef)
Definition: ellxfm.f:2
subroutine attcmp(nlines, gaclac, tlm, pos, vel, timref, time, navctl, navqc, attxfm, navblk, tiltpr, tiltfl)
Definition: attcmp.f:3
#define abs(a)
Definition: misc.h:90