/[MITgcm]/MITgcm_contrib/gael/matlab_class/gcmfaces_diags/diags_set_A.m
ViewVC logotype

Annotation of /MITgcm_contrib/gael/matlab_class/gcmfaces_diags/diags_set_A.m

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.3 - (hide annotations) (download)
Tue Dec 11 02:35:07 2012 UTC (12 years, 7 months ago) by gforget
Branch: MAIN
Changes since 1.2: +11 -3 lines
diags_grid_parms.m      add cs510 params and params switch
diags_select.m          in case mygrid.memoryLimit=1, load the
                        stuff that was not saved to diags_grid_parms.mat
diags_set_A.m           restrict list of basins to global ocean when not llc90

1 gforget 1.1
2     if userStep==1;%diags to be computed
3     listDiags='fldBAR fldTRANSPORTS gloOV gloOVbolus gloOVres gloMT_H gloMT_FW gloMT_SLT';
4 gforget 1.3 if sum([90 1170]~=mygrid.ioSize)==0;
5     listDiags=[listDiags ' atlOV atlOVbolus atlOVres atlMT_H atlMT_FW atlMT_SLT'];
6     listDiags=[listDiags ' pacindOV pacindOVbolus pacindOVres pacindMT_H pacindMT_FW pacindMT_SLT'];
7     end;
8 gforget 1.1 elseif userStep==2;%input files and variables
9     listFlds={'UVELMASS','VVELMASS','GM_PsiX','GM_PsiY'};
10     listFlds={listFlds{:},'ADVx_TH ','ADVy_TH ','DFxE_TH ','DFyE_TH '};
11     listFlds={listFlds{:},'ADVx_SLT','ADVy_SLT','DFxE_SLT','DFyE_SLT'};
12     listFldsNames=deblank(listFlds);
13     listFiles={'trsp_3d_set1','trsp_3d_set2'};
14     listSubdirs={[dirModel 'diags/TRSP/']};
15     elseif userStep==3;%computational part;
16     %mask fields:
17     fldU=UVELMASS.*mygrid.mskW; fldV=VVELMASS.*mygrid.mskS;
18     if ~isempty(whos('ADVx_TH'));%for backward compatibility
19     if size(ADVx_TH{1},3)>1; %assume full 3D fields
20     mskW=mygrid.mskW; mskS=mygrid.mskS;
21     else; %assume vertically integrated (2D fields)
22     mskW=mygrid.mskW(:,:,1); mskS=mygrid.mskS(:,:,1);
23     end;
24     ADVx_TH=ADVx_TH.*mskW; ADVy_TH=ADVy_TH.*mskS;
25     ADVx_SLT=ADVx_SLT.*mskW; ADVy_SLT=ADVy_SLT.*mskS;
26     DFxE_TH=DFxE_TH.*mskW; DFyE_TH=DFyE_TH.*mskS;
27     DFxE_SLT=DFxE_SLT.*mskW; DFyE_SLT=DFyE_SLT.*mskS;
28     end;
29     if ~isempty(whos('GM_PsiX'));%for backward compatibility
30     [fldUbolus,fldVbolus,fldWbolus]=calc_bolus(GM_PsiX,GM_PsiY);
31     fldUbolus=fldUbolus.*mygrid.mskW; fldVbolus=fldVbolus.*mygrid.mskS;
32     fldUres=fldU+fldUbolus; fldVres=fldV+fldVbolus;
33     end;
34    
35     %compute barotropic stream function:
36     [fldBAR]=calc_barostream(fldU,fldV);
37     %compute transports along transects:
38 gforget 1.2 [fldTRANSPORTS]=1e-6*calc_transports(fldU,fldV,mygrid.LINES_MASKS,{'dh','dz'});
39 gforget 1.3
40 gforget 1.1 %compute overturning stream functions:
41 gforget 1.3 if sum([90 1170]~=mygrid.ioSize)==0;
42     listBasins=[1:3];
43     else;
44     listBasins=1;
45     end;
46     for bb=listBasins;
47 gforget 1.1 %mask : global, atlantic or Pac+Ind
48     if bb==1; mskC=mygrid.mskC; mskC(:)=1;
49     elseif bb==2; mskC=v4_basin({'atlExt'}); mskC=mk3D(mskC,fldU);
50     elseif bb==3; mskC=v4_basin({'pacExt','indExt'}); mskC=mk3D(mskC,fldU);
51     end;
52     %note: while mskC is a basin mask for tracer points, it can be applied to U/V below
53     %compute overturning: eulerian contribution
54     [fldOV]=calc_overturn(fldU.*mskC,fldV.*mskC);
55     if ~isempty(whos('GM_PsiX'));%for backward compatibility
56     %compute overturning: eddy contribution
57     [fldOVbolus]=calc_overturn(fldUbolus.*mskC,fldVbolus.*mskC);
58     %compute overturning: residual overturn
59     [fldOVres]=calc_overturn(fldUres.*mskC,fldVres.*mskC);
60     else;
61     fldOVbolus=NaN*fldOV; fldOVres=NaN*fldOV;
62     end;
63    
64     if ~isempty(whos('ADVx_TH'));%for backward compatibility
65     %compute meridional heat transports:
66     tmpU=(ADVx_TH+DFxE_TH).*mskC(:,:,1:size(ADVx_TH{1},3));
67     tmpV=(ADVy_TH+DFyE_TH).*mskC(:,:,1:size(ADVx_TH{1},3));
68     [fldMT_H]=1e-15*4e6*calc_MeridionalTransport(tmpU,tmpV,0);
69     %compute meridional fresh water transports:
70     %... using the virtual salt flux formula:
71     %[fldMT_FW]=1e-6/35*calc_MeridionalTransport(ADVx_SLT+DFxE_SLT,ADVy_SLT+DFyE_SLT,0);
72     %[fldMT_FW]=1e-6/35*calc_MeridionalTransport(ADVx_SLT,ADVy_SLT,0);
73     %... using the real freshwater flux formula:
74     [fldMT_FW]=1e-6*calc_MeridionalTransport(fldU.*mskC,fldV.*mskC,1);
75     %compute meridional salt transports:
76     tmpU=(ADVx_SLT+DFxE_SLT).*mskC(:,:,1:size(ADVx_TH{1},3));
77     tmpV=(ADVy_SLT+DFyE_SLT).*mskC(:,:,1:size(ADVx_TH{1},3));
78     [fldMT_SLT]=1e-6*calc_MeridionalTransport(tmpU,tmpV,0);
79     else;
80     fldMT_H=NaN*mygrid.LATS; fldMT_FW=NaN*mygrid.LATS; fldMT_SLT=NaN*mygrid.LATS;
81     end;
82    
83     %store to global, atlantic or Pac+Ind arrays:
84     if bb==1;
85     gloMT_H=fldMT_H; gloMT_FW=fldMT_FW; gloMT_SLT=fldMT_SLT;
86     gloOV=fldOV; gloOVbolus=fldOVbolus; gloOVres=fldOVres;
87     elseif bb==2;
88     kk=find(mygrid.LATS<-30|mygrid.LATS>58);
89     fldMT_H(kk,:)=NaN; fldMT_FW(kk,:)=NaN; fldMT_SLT(kk,:)=NaN;
90     fldOV(kk,:)=NaN; fldOVbolus(kk,:)=NaN; fldOVres(kk,:)=NaN;
91     atlMT_H=fldMT_H; atlMT_FW=fldMT_FW; atlMT_SLT=fldMT_SLT;
92     atlOV=fldOV; atlOVbolus=fldOVbolus; atlOVres=fldOVres;
93     elseif bb==3;
94     kk=find(mygrid.LATS<-30|mygrid.LATS>58);
95     fldMT_H(kk,:)=NaN; fldMT_FW(kk,:)=NaN; fldMT_SLT(kk,:)=NaN;
96     fldOV(kk,:)=NaN; fldOVbolus(kk,:)=NaN; fldOVres(kk,:)=NaN;
97     pacindMT_H=fldMT_H; pacindMT_FW=fldMT_FW; pacindMT_SLT=fldMT_SLT;
98     pacindOV=fldOV; pacindOVbolus=fldOVbolus; pacindOVres=fldOVres;
99     end;
100     end;
101     end;

  ViewVC Help
Powered by ViewVC 1.1.22