forked from athurnherr/LDEO_IX
-
Notifications
You must be signed in to change notification settings - Fork 0
/
loadrdi.m
1779 lines (1526 loc) · 45.9 KB
/
loadrdi.m
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
function [d,p,de]=loadrdi(f,p)
% function [d,p,de]=loadrdi(f,p)
% LADCP-2 software
%
% Martin Visbeck and Gerd Krahmann , LDEO April-2000
% read RDI raw data
% thanks to Christian Mertens who contributed most of the functions
%
% Added pg_save variable which saves percent good data, 6/2003
% Added added up-looking beam coordiante data, 6/2004
%
%======================================================================
% L O A D R D I . M
% doc: Fri Jun 18 18:21:56 2004
% dlm: Tue Jun 25 10:46:50 2024
% (c) 2004 ladcp@
% uE-Info: 55 74 NIL 0 0 72 2 2 8 NIL ofnI
%======================================================================
% CHANGES BY ANT
% Jun 18, 2004: - added p.mask_dn_bins, p.mask_up_bins (later moved)
% Jun 21, 2004: - clarified large-velocity warning message
% Jun 22, 2004: - changed large-velocity warning to check only
% central hour of cast, as large velocities were found
% to occur commonly in uplooker near beginning and end
% of cast, but not in middle. The max number of allowed
% large velocities before warning is issued was reduced
% from 25 to zero, on the other hand.
% - removed logging messages that are not useful in general
% (e.g. number of bytes per ensemble)
% - removed `removed' messages if no ensembles were removed
% Jun 24, 2004: - BUG: new large-velocity warning code failed for casts
% shorter than one hour
% Jul 4, 2004: - BUG: new large-velocity warning code failed because
% in some cases l.tim contains NaNs
% Jul 21, 2004: - moved bin masking to [edit_data.m]
% Nov 17, 2007: - bad error message in b2earth()
% - removed a lot of commented-out code from b2earth()
% Nov 18, 2007: - added code for p.allow_3beam_solutions, p.ignore_beam
% Jan 7, 2009: - tightened use of exist()
% Oct 28, 2009: - modified p.ignore_beam for dual-head systems
% Jun 27, 2011: - l.blen removed because bin lenght can be different for UL & DL
% - apparently unused z-variable commented out
% Jun 30, 2011: - buggy bin-remapping disabled
% Aug 18, 2011: - added comment to coord-transformation code (gimbal pitch)
% Jun 24, 2013: - blen re-added but separately for DL/UL
% - added separate nbin, blnk, dist for DL/UL to p struct
% Jan 23, 2015: - made updown() bomb when UL file is not found
% Apr 15, 2015: - modified ambiguity-velocity warning as suggested by Diana Cardoso
% May 27, 2015: - clarified time-related warnings
% Feb 23, 2016: - clarified header id error message
% Jan 27, 2020: - moved magdev call further back, where start time is known
% (with GK's new magdev code, this howto is miraculously correct again - I think)
% Jun 25, 2024: - BUG: RDI BT data were used when processing only UL data
% p=setdefv(p,'pg_save',[1 2 3 4]);
% Default =3 for loadctd_whoi.
p=setdefv(p,'drot',nan);
% how many db should the last bin be below bin 1
p=setdefv(p,'ts_signal_min',-5);
p=setdefv(p,'ignore_beam',[nan nan]);
if nargin<2, p.name='unknown'; end
if existf(f,'ladcpdo')==0
error([' need file name f.ladcpdo !!!!! '])
end
if ~exist(f.ladcpdo,'file') | length(f.ladcpdo)==sum(f.ladcpdo==' ')
error([' can not find ADCP data file : ',f.ladcpdo])
end
f=setdefv(f,'ladcpup',' ');
disp('LOADRDI:')
%
% set some defaults, in case they have not been set yet
%
p=setdefv(p,'pglim',0);
p=setdefv(p,'elim',0.5);
p=setdefv(p,'vlim',2.5);
p=setdefv(p,'wlim',0.2);
p=setdefv(p,'timoff',0);
p=setdefv(p,'drot',NaN);
p=setdefv(p,'orig',0);
p=setdefv(p,'weighbin1',1);
p=setdefv(p,'tiltmax',[22 4]);
p=setdefv(p,'name','unknown');
p=setdefv(p,'maxbinrange',0);
p=setdefv(p,'ts_save',0);
p=setdefv(p,'cm_save',0);
p=setdefv(p,'pg_save',3);
p=setdefv(p,'xmv_min',0);
%
% load RDI data
%
[l,message,le]=updown(f.ladcpdo,f.ladcpup,p.pglim,p.elim,...
p.maxbinrange,p.ts_save,p.cm_save,p.pg_save,p);
p.warn=l.warn;
% store original data
%
if p.orig
d.l = l;
end
% check for data set
%
de=0;
if le==1
disp(' !!! NO DATA !!! ')
disp(message)
de=1;
return
end
%
% remember which bin come from which instrument (up-down) configuration
% get instrument configuration
%
d.izd = 1:length(l.zd);
d.zd = l.zd;
d.bbadcp = l.bbadcp;
p.serial_cpu_d = l.serial_cpu_d;
p.nping_total = l.npng_d*l.nens_d;
p.instid(1) = prod(p.serial_cpu_d+1)+sum(p.serial_cpu_d);
p.blen_d = l.blen_d;
p.nbin_d = l.nbin_d;
p.blnk_d = l.blnk_d;
p.dist_d = l.dist_d;
[dummy,d.down]=rditype(f.ladcpdo);
if d.down.Up
warn=(' up looking instrument detected in do-file');
p.warnp(size(p.warnp,1)+1,1:length(warn))=warn;
disp(warn)
d.zd=-d.zd;
end
if existf(l,'zu')
d.izu = fliplr(1:length(l.zu));
d.izd = d.izd+length(d.izu);
d.zu = l.zu;
p.serial_cpu_u = l.serial_cpu_u;
p.instid(2) = prod(p.serial_cpu_u+1)+sum(p.serial_cpu_u);
p.nping_total(2) = l.npng_u*l.nens_u;
p.blen_u = l.blen_u;
p.nbin_u = l.nbin_u;
p.blnk_u = l.blnk_u;
p.dist_u = l.dist_u;
[dummy,d.up]=rditype(f.ladcpup);
if d.up.Up==0
warn=(' down looking instrument detected in up-file');
p.warnp(size(p.warnp,1)+1,1:length(warn))=warn;
disp(warn)
d.zu=-d.zu;
end
else
d.izu=[];
d.zu=[];
end
%
% apply w velocity threshold
%
d.wrange=5;
izr=d.izd(1:d.wrange);
w(d.izd,:)=meshgrid(medianan(l.w(izr,:),1),d.izd,2);
if existf(l,'zu')
izr=[izr,d.izu(1:d.wrange)];
w(d.izu,:)=meshgrid(medianan(l.w(d.izu(1:d.wrange),:),1),d.izu,2);
end
p=setdefv(p,'wizr',izr);
% to normal velocity data
j = find(abs(l.w-w) > p.wlim);
l.u(j) = NaN;
l.v(j) = NaN;
l.w(j) = NaN;
if p.orig
d.l.problem(j) = d.l.problem(j)+1;
end
% Estimate single ping velocity error from std(W)
nmax=min([length(d.izd),6]);
sw=stdnan(l.w(d.izd(2:nmax),:));
ii=find(sw>0); sw=medianan(sw(ii));
d.down.Single_Ping_Err=sw/tan(d.down.Beam_angle*pi/180)/...
sqrt(d.down.Pings_per_Ensemble);
p.beamangle=d.down.Beam_angle;
if existf(l,'zu')
nmax=min([length(d.izu),6]);
sw=stdnan(l.w(d.izu(2:nmax),:));
ii=find(sw>0); sw=medianan(sw(ii));
d.up.Single_Ping_Err=sw/tan(d.up.Beam_angle*pi/180)/...
sqrt(d.up.Pings_per_Ensemble);
end
% to bottom track velocity data
j = find(abs(l.wb-w(d.izd(1),:)) > p.wlim);
l.ub(j) = NaN;
l.vb(j) = NaN;
l.wb(j) = NaN;
if p.orig
d.l.problemb(j) = d.l.problemb(j)+1;
end
% Horizontal Velocity limit
%
% to normal velocity data
vel=sqrt(l.u.^2+l.v.^2);
j = find(vel > p.vlim);
l.u(j) = NaN;
l.v(j) = NaN;
l.w(j) = NaN;
if length(j) > 0
disp(sprintf(' removed %d values because of horizontal speed > %g m/s',length(j),p.vlim));
end
% only warn if large velocities occur during middle hour of cast; this
% excludes near-surface effects when large velocities can be common.
% However, reduce number of allowed large velocities before warning is
% issued from 25 to 10.
nens = size(l.u,2);
enstime = 86400 * (max(l.tim(1,:))-min(l.tim(1,:))) / nens;
jstart = floor(nens/2-1800/enstime);
jend = ceil(nens/2+1800/enstime);
if (jstart < 1), jstart = 1; end
if (jend > length(vel)), jend = length(vel); end
jj = find(vel(jstart:jend) > p.vlim);
skipnens = 1200 / enstime;
[j1,j2] = find(vel > p.vlim);
jj = find(j2>skipnens & j2<size(l.u,2)-skipnens);
if length(jj)>10
warn = sprintf('** found %d (%.1f%% of total) velocity measurements > %g m/s',...
length(jj),length(jj)/length(vel)*100,p.vlim);
disp(warn);
if length(jj)>100
p.warn(size(p.warn,1)+1,1:length(warn))=warn;
end
disp('** WARNING check ambiguity velocity setting in CMD-file ** ')
end
if p.orig
d.l.problem(j) = d.l.problem(j)+1;
end
% to bottom track velocity data
vel=sqrt(l.ub.^2+l.vb.^2);
j = find(vel > p.vlim);
l.ub(j) = NaN;
l.vb(j) = NaN;
l.wb(j) = NaN;
if p.orig
d.l.problemb(j) = d.l.problemb(j)+1;
end
%
% apply a time offset, if given
%
if p.timoff~=0
disp([' WARNING adjusted ADCP time by ',num2str(p.timoff),' days']),
l.tim=l.tim+p.timoff;
end
% cut time range and apply a time offset, if given
if existf(p,'time_start')==0
disp(' using whole profile since no start time was given')
it=find(isfinite(l.tim(1,:)));
p.time_start=gregoria(l.tim(1,it(1)));
p.time_end=gregoria(l.tim(1,it(end)));
d.time_jul=l.tim(1,it);
else
% fix time for NB-ADCP data
if l.bbadcp==0
dum=gregoria(l.tim(1));
l.tim=l.tim-julian(dum)+julian([p.time_start(1) dum(2:end)]);
disp(' adjust year for NB-ADCP using given start time ')
end
d.time_jul=l.tim(1,:);
it=find(julian(p.time_start)<=d.time_jul & julian(p.time_end)>=d.time_jul);
it2=find(julian(p.time_start)>d.time_jul | julian(p.time_end)<d.time_jul);
if p.orig
d.l.problem(:,it2) = d.l.problem(:,it2)+10;
d.l.problemb(:,it2) = d.l.problemb(:,it2)+10;
end
d.time_jul = d.time_jul(it);
disp([' extracting ',int2str(length(it)),' ensembles as profile'])
end
%
% check whether the given time ranges have lead to a profile
%
if length(it)==0
disp(' ')
disp(' given time range resulted in empty extracted array!')
disp(' ')
disp(' check times')
disp([' given start time :'])
disp(p.time_start)
disp([' internal LADCP start time :'])
disp(gregoria(l.tim(1,1)))
disp([' given end time :'])
disp(p.time_end)
disp([' internal LADCP end time :'])
disp(gregoria(maxnan(l.tim(1,:))))
disp('!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!')
disp('Will try to use all data instead')
disp('!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!')
pause(2)
it=find(isfinite(l.tim(1,:)));
d.time_jul=l.tim(1,it);
% p.time_start=gregoria(l.tim(it(1)));
% p.time_end=gregoria(l.tim(it(end)));
end
if existf(p,'poss')
if isfinite(sum(p.poss))
drot=magdev(p.poss(1)+p.poss(2)/60, p.poss(3)+p.poss(4)/60,0,p.time_start(1));
if existf(p,'drot')
if isfinite(p.drot)
disp([' found drot:',num2str(p.drot),' should be ',num2str(drot)])
else
p.drot=drot;
end
else
p.drot=drot;
end
end
end
%
% rotate for magnetic deviation
%
d.soundc=0;
if isfinite(p.drot)
[d.ru,d.rv]=uvrot(l.u(:,it),l.v(:,it),p.drot);
[ub,vb]=uvrot(l.ub(it),l.vb(it),p.drot);
disp(' apply magnetic deviation, rotate bottom track and water velocities')
else
d.ru=l.u(:,it);
d.rv=l.v(:,it);
ub=l.ub(it);
vb=l.vb(it);
end
%
% save all bottom track data together and remove bad data
%
d.bvel=[ub',vb',l.wb(it)',l.eb(it)'];
ii=find(d.bvel<-30);
d.bvel(ii)=NaN;
d.hbot=l.hb(it);
d.hbot4=l.hb4(:,it);
p.btrk_used=l.btrk_used;
if existf(l,'hs')==1
d.hsurf=l.hs(it);
end
%
% extract the profile from all recorded data
%
d.firstlastindx = [it(1),it(end)];
d.rw=l.w(:,it);
d.re=l.e(:,it);
d.ts=l.ts(:,it);
if p.ts_save(1)~=0
d.ts_all_d=l.ts_all_d(it,:,:);
if size(l.pit,1)==2
d.ts_all_u=l.ts_all_u(it,:,:);
end
end
if p.cm_save(1)~=0
d.cm_all_d=l.cm_all_d(it,:,:);
if size(l.pit,1)==2
d.cm_all_u=l.cm_all_u(it,:,:);
end
end
if p.pg_save(1)~=0
d.pg_all_d=l.pg_all_d(it,:,:);
if size(l.pit,1)==2
d.pg_all_u=l.pg_all_u(it,:,:);
end
end
d.hdg=l.hdg(:,it);
d.xmc=l.xmc(:,it);
d.xmv=l.xmv(:,it);
d.tint=l.tint(:,it);
d.sv=l.sv(:,it);
d.temp=l.t(:,it);
d.weight=l.cm(:,it);
d.weight=d.weight./medianan(maxnan(d.weight));
d.pit = l.pit(:,it);
d.rol = l.rol(:,it);
%d.tilt=sqrt(l.pit(1,it).^2 + l.rol(1,it).^2);
% more accurate calculation
d.tilt=real(asin(sqrt(sin(l.pit(1,it)/180*pi).^2 +...
sin(l.rol(1,it)/180*pi).^2)))/pi*180;
% compute tilt difference
rold=mean(abs(diff([0,d.rol(1,:);d.rol(1,:),0]'))');
pitd=mean(abs(diff([0,d.pit(1,:);d.pit(1,:),0]'))');
d.tiltd=sqrt(rold.^2+pitd.^2);
% reduce weight for strong tilt difference
ii=find(d.tilt>p.tiltmax(1));
if length(ii) > 0
disp([' removed ',num2str(length(ii)),...
' profiles due to tilt > ',num2str(p.tiltmax(1)) ' degrees'])
end
d.weight(:,ii)=NaN;
if length(ii)>length(d.tilt)*0.1;
warn=([' ',int2str(length(ii)*100/length(d.tilt)),...
'% tilt > ',int2str(p.tiltmax(1)),' ']);
disp(warn)
p.warn(size(p.warn,1)+1,1:length(warn))=warn;
p.warnp(size(p.warnp,1)+1,1:length(warn))=warn;
end
ii=find(d.tiltd>p.tiltmax(2));
if length(ii) > 0
disp([' removed ',num2str(length(ii)),...
' profiles due to tilt derivative > ',num2str(p.tiltmax(2)) ' degrees'])
end
d.weight(:,ii)=NaN;
% reduce weight for strong echos possibly from crosstalk or bottom
d.tsw=d.weight*0+1;
for i=1:size(d.tsw,1)
tsmed=median(d.ts(i,:));
ts=d.ts(i,:)-tsmed;
ii=find(ts>0);
if length(ii)>0
d.weight(i,ii)=d.weight(i,ii).*(1-(ts(ii)/max(ts)).^1.5);
d.tsw(i,ii)=(1-(ts(ii)/max(ts)).^1.5);
end
end
% save transmit current volt and internal temperature
for j=1:size(d.xmv,1)
p.xmc(j)=medianan(d.xmc(j,:),size(d.xmc,2)/4);
p.xmv(j)=medianan(d.xmv(j,:),size(d.xmv,2)/4);
p.tint(j)=medianan(d.tint(j,:),size(d.tint,2)/4);
end
% warn if battery low
if p.xmv(1)<p.xmv_min
warn=([' median Xmit-volt ',num2str(p.xmv(1)),' < ',...
num2str(p.xmv_min),' BATTERY weak ? ']);
disp(warn)
p.warn(size(p.warn,1)+1,1:length(warn))=warn;
end
% check for dead instrument (no pings just listen)
%
% remove suspect ensembles
drw=medianan(abs(diff(d.rw(d.izd,:))));
dru=medianan(abs(diff(d.rv(d.izd,:))));
drv=medianan(abs(diff(d.ru(d.izd,:))));
nbad=find(abs(drw)<0.005 & abs(dru)<0.005 & abs(dru)<0.005);
if length(nbad) > 0.2*length(it)
warn=([' down looker ',int2str(length(nbad)),' ensembles ',...
' have no flow gradient. ']);
disp(warn)
p.warn(size(p.warn,1)+1,1:length(warn))=warn;
warn=(['DOWN LOOKER NOT PINGING ?']);
disp(warn)
p.warn(size(p.warn,1)+1,1:length(warn))=warn;
disp(' WARNING WARNING ')
end
if length(nbad)>0
d.weight(d.izd,nbad)=nan;
disp([' removed ',int2str(length(nbad)),...
' suspect non pinging? (low velocity gradient) ensembles from down-looker'])
end
if length(d.izu)>1
% remove suspect ensembles
drw=medianan(abs(diff(d.rw(d.izu,:))));
dru=medianan(abs(diff(d.rv(d.izu,:))));
drv=medianan(abs(diff(d.ru(d.izu,:))));
nbad=find(abs(drw)<0.005 & abs(dru)<0.005 & abs(dru)<0.005);
if length(nbad) > 0.2*length(it)
warn=([' up looker ',int2str(length(nbad)),' ensembles ',...
' have no flow gradient.']);
disp(warn)
p.warn(size(p.warn,1)+1,1:length(warn))=warn;
warn=(['UP LOOKER NOT PINGING ?']);
disp(warn)
p.warn(size(p.warn,1)+1,1:length(warn))=warn;
disp(' WARNING WARNING ')
end
if length(nbad)>0
d.weight(d.izu,nbad)=nan;
disp([' removed ',int2str(length(nbad)),...
' suspect non pinging? (low velocity gradient) ensembles from up-looker'])
end
end
%
% reduce certainty in bin 1
%
idb1=d.izd(1);
if length(idb1)==1
d.weight(idb1,:)=d.weight(idb1,:)*p.weighbin1;
end
if length(d.izu)>1
iub1=find(d.izu==1);
if length(iub1)==1
d.weight(iub1,:)=d.weight(iub1,:)*p.weighbin1;
end
end
if p.weighbin1~=1
disp([' multiply weight of bin 1 by ',num2str(p.weighbin1)])
end
%
% prepare array
%
d.izm=d.weight+NaN;
%
% save mean correlation and echo amp profile
%
if length(size(l.tsd_m))==2
d.tsd_m=reshape(l.tsd_m,length(l.tsd_m)/4,4);
else
d.tsd_m=squeeze(l.tsd_m);
end
if length(size(l.cmd_m))==2
d.cmd_m=reshape(l.cmd_m,length(l.cmd_m)/4,4);
else
d.cmd_m=squeeze(l.cmd_m);
end
t=d.tsd_m(1,:);
if abs(min(t)-median(t))>15
disp('!!!! WARNING one beam might be broken !!!!!!!!!')
[m,it]=min(t);
warn=([' broken down looking beam ',int2str(it)]);
disp(warn)
p.warn(size(p.warn,1)+1,1:length(warn))=warn;
elseif abs(min(t)-median(t))>5
disp('WARNING one beam might be weak !')
[m,it]=min(t);
warn=([' weak down looking beam ',int2str(it)]);
disp(warn)
p.warn(size(p.warn,1)+1,1:length(warn))=warn;
end
if existf(l,'tsu_m')
if length(size(l.tsu_m))==2
d.tsu_m=reshape(l.tsu_m,length(l.tsu_m)/4,4);
else
d.tsu_m=squeeze(l.tsu_m);
end
if length(size(l.cmu_m))==2
d.cmu_m=reshape(l.cmu_m,length(l.cmu_m)/4,4);
else
d.cmu_m=squeeze(l.cmu_m);
end
t=d.tsu_m(1,:);
if abs(min(t)-median(t))>10
disp('!!!! WARNING one beam might be broken !!!!!!!!!')
[m,it]=min(t);
warn=([' broken up looking beam ',int2str(it)]);
disp(warn)
p.warn(size(p.warn,1)+1,1:length(warn))=warn;
elseif abs(min(t)-median(t))>5
disp('WARNING one beam might be weak ! ')
[m,it]=min(t);
warn=([' weak up looking beam ',int2str(it)]);
disp(warn)
p.warn(size(p.warn,1)+1,1:length(warn))=warn;
end
end
[p.nbins,p.nt]=size(d.ru);
% check for outlier within the whole data set
[d,p]=outlier(d,p);
%-----------------------------------------------------------
function [l,message,le] = updown(fdown,fup,pglim,elim,...
bmax,tssave,cmsave,pgsave,p)
%UPDOWN Load and merge upward and downward looking ADCP raw data.
% L = UPDOWN('filedown','fileup') reads ADCP raw data from the specified
% files. L is a structure array with the following fields:
%
% blen: bin length %%% ANT: REMOVED 2011/06/28 BECAUSE IT CAN BE DIFFERENT FOR UL/DL
% nbin: number of bins
% blnk: blank after transmit
% dist: distance of bin 1 from transducer
% tim: time axis
% pit: pitch
% rol: roll
% hdg: heading
% s: salinity
% t: temperature
% sv: sound velocity
% u: east velocity
% v: north velocity
% w: vertical velocity
% e: error velocity
% ts: target strength
% cm: correlation
% hb: bottom track distance
% ub: bottom track east velocity
% vb: bottom track north velocity
% wb: bottom track vertical velocity
% eb: bottom track error velocity
%
% [L,MESSAGE] = UPDOWN('filedown','fileup') returns a system dependent error
% message if the opening of 'filedown' is not successful. In this case -1 is
% returned for L.
% Christian Mertens, IfM Kiel
% default editing parameters
if nargin < 5
bmax = 0;
end
if nargin < 4
elim = 0.5;
end
if nargin < 3
pglim = 0.3;
end
% fixed leader
f.nbin = 1; % number of depth cells
f.npng = 2; % pings per ensemble
f.blen = 3; % depth cell length
f.blnk = 4; % blank after transmit
f.dist = 5; % distance to the middle of the first depth cell
f.plen = 6; % transmit pulse length
f.serial = '7:14'; % serial number of CPU board
% variable leader
v.tim = 1; % true time (Julian days)
v.pit = 2; % pitch
v.rol = 3; % roll
v.hdg = 4; % heading
v.t = 5; % temperature
v.s = 6; % salinity
v.sv = 7; % sound velocity
v.xmc = 8; % transmit current
v.xmv = 9; % transmit volt
v.tint = 10; % internal temperature
% load downward looking ADCP
[fid,message] = fopen(fdown,'r','l');
le=0;
if fid == -1
le = 1;
message = sprintf('%s: %s',fdown,message);
disp(' LOADRDI problem with down looking RDI file ')
disp(message)
error('terminate LADCP processing')
end
disp([' loading down-data ',fdown])
% check if BB data
if isbb(fid)
[fd,vd,veld,cmd,ead,pgd,btd] = rdread(fid);
l.bbadcp=1;
else
fclose(fid);
fid = fopen(fdown,'r','b');
[fd,vd,veld,swd,ead,pgd] = nbread(fid);
cmd=ead*0+100;
btd = NaN*ones(size(vd,1),1,16);
ok = double(prod(veld,3)~=sum(veld,3));
ii=find(ok==0);
ok(ii)=NaN;
disp([' removed ',int2str(length(ii)),...
' values because of 0 in nbdata'])
for k = 1:4
veld(:,:,k) = veld(:,:,k).*ok;
end
l.bbadcp=0;
end
% check for beam coordinates
[dummy,dd]=rditype(fdown);
if dd.Coordinates==0
disp(' DETECTED BEAM coordinates: rotating to EARTH coordinates')
veld=b2earth(veld,vd,dd,p,p.ignore_beam(1));
end
%
fclose(fid);
l1=size(veld,1);
l2=size(veld,2);
disp([' read ',int2str(l1),' ensembles with ',int2str(l2),' bins each'])
% remove extra (?) bottom track dimension
btd = squeeze(btd);
if ndims(btd)>2
warning(' removal of extra bottom track dimension failed !!!')
end
% median echoamplitude and correlation
%ead = targs(mean(ead,3)',z(:))';
ead_m=medianan(ead);
cmd_m=medianan(cmd);
if tssave(1)~=0
ead_all=ead(:,:,tssave);
end
if cmsave(1)~=0
cmd_all=cmd(:,:,cmsave);
end
if pgsave(1)~=0
pgd_all=pgd(:,:,pgsave);
end
ead=median(ead,3);
cmd=median(cmd,3);
% prepare vector containing info why a value has been discarded
% b is the same for the bottom track data
%
% 1st digit : w deviates more than p.wlim from median w (checked later)
% 2nd digit : out of time range (checked later)
% 3rd digit : below percent good threshold
% 4th digit : above error velocity threshold
% 5th digit : 3 beam solution
% 6th digit : no vel
%
%GK
dproblem = repmat(0,[size(veld,1),size(veld,2)]);
problemb = repmat(0,[size(btd,1),size(btd,2)]);
% grep 3-beam solution and no velocities at all
i = find( isnan(veld(:,:,4)) & ~isnan(veld(:,:,3)) );
dproblem(i) = dproblem(i) + 10000;
i = find( isnan(veld(:,:,1)) );
dproblem(i) = dproblem(i) + 100000;
l.warn=('LADCP WARNINGS');
l.btrk_used = 0;
if sum(isfinite(btd(:)))>0
[dummy,db]=rditype(fdown);
if db.Up
warn=(' Warning: ignoring RDI bottom track data from upward-looking instrument');
disp(warn)
l.warn(size(l.warn,1)+1,1:length(warn))=warn;
else
l.btrk_used = 1;
end
end
% transform to earth coordinates
if l.btrk_used == 1
[dummy,db]=rditype(fdown);
if db.Coordinates==0
for i=1:4
velb(:,1,i)=btd(:,4+i);
velb(:,2,i)=NaN;
end
disp(' DETECTED BEAM bottom track coordinates!')
db.use_binremap=0;
velb=b2earth(velb,vd,db,p,p.ignore_beam(1));
for i=1:4
btd(:,4+i)=velb(:,1,i);
end
end
end
% apply percent-good threshold
pgd = pgd(:,:,4);
i = pgd < pglim;
pgd(i) = NaN;
pgd(~i) = 1;
if length(find(i)) > 0
disp(sprintf(' removed %d downlooker values because of percent good < %g',...
length(find(i)),pglim));
end
for k = 1:4
veld(:,:,k) = veld(:,:,k).*pgd;
end
dproblem(i) = dproblem(i) + 100;
% load upward looking ADCP
up = nargin>1;
if strcmp(fup,'') || strcmp(fup,' ')
up = 0;
end
if up
fid = fopen(fup,'r','l');
if fid == -1
error(sprintf('%s: no such file or directory',fup));
end
if length(bmax)<2, bmax(2)=bmax(1); end
disp([' loading up-data ',fup])
[fu,vu,velu,cmu,eau,pgu,btu] = rdread(fid);
% check for beam coordinates
[dummy,du]=rditype(fup);
if du.Coordinates==0
disp(' DETECTED BEAM coordinates: rotating to EARTH coordinates')
velu=b2earth(velu,vu,du,p,p.ignore_beam(2));
end
%
l1=size(velu,1);
l2=size(velu,2);
btu = squeeze(btu);
disp([' read ',int2str(l1),' ensemble and ',int2str(l2),' bins '])
% median echoamplitude and correlation
%eau = targs(mean(eau,3)',z(:))';
eau_m=medianan(eau);
cmu_m=medianan(cmu);
if tssave(1)~=0
eau_all=eau(:,:,tssave);
end
if cmsave(1)~=0
cmu_all=cmu(:,:,cmsave);
end
if pgsave(1)~=0
pgu_all=pgu(:,:,pgsave);
end
eau=median(eau,3);
cmu=median(cmu,3);
fclose(fid);
% prepare vector containing info why a value has been discarded
% b is the same for the bottom track data
uproblem = repmat(0,[size(velu,1),size(velu,2)]);
% grep 3-beam solution and no velocities at all
i = find( velu(:,:,4) == 0 & velu(:,:,3)~=0 );
uproblem(i) = uproblem(i) + 10000;
i = find( velu(:,:,1) == 0 );
uproblem(i) = uproblem(i) + 100000;
% apply percent-good threshold
pgu = pgu(:,:,4);
i = pgu < pglim;
if length(find(i)) > 0
disp(sprintf(' removed %d uplooker values because of percent good < %g',...
length(find(i)),pglim));
end
pgu(i) = NaN;
pgu(~i) = 1;
for k = 1:4
velu(:,:,k) = velu(:,:,k).*pgu;
end
uproblem(i) = uproblem(i) + 100;
end
% distance vectors
%%% z = fd(f.dist) + fd(f.blen)*([1:fd(f.nbin)] - 1); % unused?! ANT 2011/06/28
if bmax(1)>0, fd(f.nbin)=bmax(1); end
idb=1:fd(f.nbin);
if up
if bmax(2)>0, fu(f.nbin)=bmax(2); end
iub=1:fu(f.nbin);
% check if ping rate is the same for both instruments
timd = vd(:,1,v.tim)';
timu = vu(:,1,v.tim)';
pingconst=abs(maxnan(diff(timd))+maxnan(-diff(timd))) > (0.05/24/3600);
pingdiff=abs(medianan(diff(timd))-medianan(diff(timu))) > (0.05/24/3600);
ii=min(length(timd),length(timu));
timediff=abs((timd(1)-timd(ii)-(timu(1)-timu(ii)))) > (5/24/3600);
if pingdiff | timediff | pingconst
disp(' WARNING WARNING WARNING WARNING WARNING WARNING WARNING WARNING ')
if pingconst
warn=(' Warning: non-constant ping rate in downlooker data (staggered pinging?)');
disp(warn)
disp([' min down ping rate :',num2str(-24*3600*maxnan(-diff(timd))),...
' max down ping rate :',num2str(24*3600*maxnan(diff(timd)))])
end
if pingdiff
warn=(' Warning: mean ping rates differ in downlooker/uplooker data ');
disp(warn)
l.warn(size(l.warn,1)+1,1:length(warn))=warn;
disp([' mean down ping rate :',num2str(24*3600*meannan(diff(timd))),...
' mean up ping rate :',num2str(24*3600*meannan(diff(timu)))])
end
if timediff
warn=(' Warning: cast duration differs in downlooker/uplooker data ');
disp(warn)
l.warn(size(l.warn,1)+1,1:length(warn))=warn;
disp([' down dt for common ping number:',num2str((timd(ii)-timd(1))*24),...
' up dt :',num2str((timu(ii)-timu(1))*24),' hours '])
end
iu=1:length(timd);
ii=find(iu>length(timu));
iu(ii)=length(timu);
disp(' find best time match of up-looking ADCP to down looking ADCP')
for i=find(isfinite(timd))
[m,iu(i)]=min(abs(timu-timd(i)));
end
ilast=min(length(iu),length(timu));
disp([' up instrument is different by ',num2str(iu(ilast)-ilast),...
' ensembles'])
disp(' WARNING WARNING WARNING WARNING WARNING WARNING WARNING WARNING ')
id=1:length(timd);
else
id=1:length(timd);
iu=1:length(timu);
end
% find best lag to match up vertical velocity
wu=squeeze(velu(iu,:,3));
wd=squeeze(veld(id,:,3));
wb2u=medianan(wu');
wb2d=medianan(wd');
maxlag=20;
[lag,iiu,id,co]=bestlag(wb2u,wb2d,maxlag);
disp([' try to shift timeseries by lag: ',num2str(lag),...
' correlation: ',num2str(co)])
if abs(lag)==maxlag | co<0.9,
disp(' best lag not obvious! use time to match up-down looking ADCP')
id=1:length(timd);
iu=1:length(timd);
ii=find(iu>length(timu));
iu(ii)=length(timu);
for i=find(isfinite(timd))
[m,iu(i)]=min(abs(timu-timd(i)));
end
lag=mean(iu-id);
disp([' mean lag is ',num2str(lag),' ensembles']);
else
disp([' shift ADCP timeseries by ',num2str(lag),' ensembles']);
iu=iu(iiu);
end
disp([' number of joint ensembles is : ',num2str(length(iu))]);
% merge upward and downward
l.zu=[0:(fu(f.nbin)-1)]*fu(f.blen)+fu(f.dist);
l.zd=[0:(fd(f.nbin)-1)]*fd(f.blen)+fd(f.dist);
eval(['l.serial_cpu_u=fu(',f.serial,');']);
eval(['l.serial_cpu_d=fd(',f.serial,');']);
l.npng_u = fu(f.npng);
l.npng_d = fd(f.npng);
l.nens_u = size(vu,1);
l.nens_d = size(vd,1);
l.blen_u = fu(f.blen);
l.blen_d = fd(f.blen);
l.nbin_u = fu(f.nbin);
l.nbin_d = fd(f.nbin);
l.blnk_u = fu(f.blnk);
l.blnk_d = fd(f.blnk);
l.dist_u = fu(f.dist);
l.dist_d = fd(f.dist);
l.tim = [vd(id,1,v.tim),vu(iu,1,v.tim)]';
l.pit = [vd(id,1,v.pit),vu(iu,1,v.pit)]';
l.rol = [vd(id,1,v.rol),vu(iu,1,v.rol)]';
l.hdg = [vd(id,1,v.hdg),vu(iu,1,v.hdg)]';
l.s = [vd(id,1,v.s),vu(iu,1,v.s)]';
l.t = [vd(id,1,v.t),vu(iu,1,v.t)]';
l.sv = [vd(id,1,v.sv),vu(iu,1,v.sv)]';
l.xmc = [vd(id,1,v.xmc),vu(iu,1,v.xmc)]';
l.xmv = [vd(id,1,v.xmv),vu(iu,1,v.xmv)]';
l.tint = [vd(id,1,v.tint),vu(iu,1,v.tint)]';
l.u = [fliplr(velu(iu,iub,1)) veld(id,idb,1)]';
l.v = [fliplr(velu(iu,iub,2)) veld(id,idb,2)]';
l.w = [fliplr(velu(iu,iub,3)) veld(id,idb,3)]';
l.e = [fliplr(velu(iu,iub,4)) veld(id,idb,4)]';
l.ts = [fliplr(eau(iu,iub)) ead(id,idb)]';
l.cm = [fliplr(cmu(iu,iub)) cmd(id,idb)]';
% No reason to keep this since pgu and pgd don't mean much anymore