/[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.1 - (hide annotations) (download)
Thu Nov 1 19:07:17 2012 UTC (12 years, 8 months ago) by gforget
Branch: MAIN
- more modular, and suposedly more user-friendly, version
  of ecco_v4/comp_driver.m and basic_diags_ecco.m

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

  ViewVC Help
Powered by ViewVC 1.1.22