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

Annotation of /MITgcm_contrib/gael/matlab_class/gcmfaces_diags/diags_set_D.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     %select kBudget:
3     if ~isempty(setDiagsParams);
4     kBudget=setDiagsParams{1};
5     else;
6     kBudget=1;
7     end;
8    
9     %preliminary tests and pre-processing for budget computations
10     test1=isempty(dir([dirModel 'diags/BUDG/budg2d_snap_set1*']));
11     test2=isempty(dir([dirMat 'BUDG/rate_budg2d_snap_set1*']));
12    
13     if (strcmp(setDiags,'D')&test1&test2);
14     fprintf('\n abort : global and regional budgets, due to missing \n');
15     fprintf(['\n ' dirModel 'diags/BUDG/budg2d_snap_set1* \n']);
16     return;
17     end;
18    
19     if (strcmp(setDiags,'D')&test2);
20     fprintf('\n abort : global and regional budgets, due to missing \n');
21     fprintf(['\n ' dirModel 'diags/BUDG/rate_budg2d_snap_set1* \n']);
22     return;
23     end;
24    
25     if userStep==1;%diags to be computed
26     listDiags=['glo_vol_ocn glo_vol_tot glo_vol_ice glo_bp'];
27     listDiags=[listDiags ' north_vol_ocn north_vol_tot north_vol_ice north_bp'];
28     listDiags=[listDiags ' south_vol_ocn south_vol_tot south_vol_ice south_bp'];
29     listDiags=[listDiags ' glo_heat_ocn glo_heat_tot glo_heat_ice'];
30     listDiags=[listDiags ' north_heat_ocn north_heat_tot north_heat_ice'];
31     listDiags=[listDiags ' south_heat_ocn south_heat_tot south_heat_ice'];
32     listDiags=[listDiags ' glo_salt_ocn glo_salt_tot glo_salt_ice'];
33     listDiags=[listDiags ' north_salt_ocn north_salt_tot north_salt_ice'];
34     listDiags=[listDiags ' south_salt_ocn south_salt_tot south_salt_ice'];
35     elseif userStep==2;%input files and variables
36     listFlds={ 'ETAN','SIheff','SIhsnow','THETA ','SALT ','PHIBOT'};
37     listFlds={listFlds{:},'SIatmFW ','oceFWflx','SItflux','TFLUX','SFLUX','oceSPflx','SRELAX'};
38     listFlds={listFlds{:},'oceQnet ','SIatmQnt','SIaaflux','SIsnPrcp','SIacSubl'};
39     listFlds={listFlds{:},'TRELAX','WTHMASS','WSLTMASS','oceSflux','oceQsw','oceSPtnd'};
40     if kBudget>1;
41     listFlds={listFlds{:},'ADVr_TH','DFrE_TH','DFrI_TH','ADVr_SLT','DFrE_SLT','DFrI_SLT','WVELMASS'};
42     end;
43     listFlds={listFlds{:},'UVELMASS','VVELMASS','AB_gT','AB_gS'};
44     listFlds={listFlds{:},'ADVx_TH ','ADVy_TH ','DFxE_TH ','DFyE_TH '};
45     listFlds={listFlds{:},'ADVx_SLT','ADVy_SLT','DFxE_SLT','DFyE_SLT'};
46     listFlds={listFlds{:},'ADVxHEFF','ADVyHEFF','DFxEHEFF','DFyEHEFF'};
47     listFlds={listFlds{:},'ADVxSNOW','ADVySNOW','DFxESNOW','DFyESNOW'};
48     listFldsNames=deblank(listFlds);
49     %
50     listFiles={'rate_budg2d_snap_set1','budg2d_hflux_set1','budg2d_zflux_set1','budg2d_zflux_set2'};
51     if kBudget==1;
52     listFiles={listFiles{:},'rate_budg2d_snap_set2','budg2d_hflux_set2'};
53     else;
54     tmp1=sprintf('rate_budg2d_snap_set3_%02i',kBudget);
55     tmp2=sprintf('budg2d_zflux_set3_%02i',kBudget);
56     tmp3=sprintf('budg2d_hflux_set3_%02i',kBudget);
57     listFiles={listFiles{:},tmp1,tmp2,tmp3};
58     end;
59     listSubdirs={[dirMat 'BUDG/' ],[dirModel 'diags/BUDG/']};
60     elseif userStep==3;%computational part;
61    
62     %override default file name:
63     %---------------------------
64     tmp1=setDiags;
65     if kBudget>1;
66     tmp1=sprintf('D%02i',kBudget);
67     end;
68     fileMat=['diags_set_' tmp1 '_' num2str(tt) '.mat'];
69    
70     %fill in optional fields:
71     %------------------------
72     if isempty(who('TRELAX')); TRELAX=0*mygrid.XC; end;
73     if isempty(who('SRELAX')); SRELAX=0*mygrid.XC; end;
74     if isempty(who('AB_gT')); AB_gT=0*mygrid.XC; end;
75     if isempty(who('AB_gS')); AB_gS=0*mygrid.XC; end;
76     if isempty(who('oceSPtnd')); oceSPtnd=0*mygrid.XC; end;
77     if isempty(who('oceSPflx')); oceSPflx=0*mygrid.XC; end;
78     if isempty(who('PHIBOT')); PHIBOT=0*mygrid.XC; end;
79    
80     %=======MASS=========
81    
82     %compute mapped budget:
83     %----------------------
84    
85     %mass = myparms.rhoconst * sea level
86     contOCN=ETAN*myparms.rhoconst;
87     contICE=(SIheff*myparms.rhoi+SIhsnow*myparms.rhosn);
88     %for deep ocean layer :
89     if kBudget>1&myparms.useNLFS<2;
90     contOCN=0;
91     elseif kBudget>1;%rstar case
92     tmp1=mk3D(mygrid.DRF,mygrid.hFacC).*mygrid.hFacC;
93     tmp2=sum(tmp1(:,:,kBudget:length(mygrid.RC)),3)./mygrid.Depth;
94     contOCN=tmp2.*ETAN*myparms.rhoconst;
95     end;
96     %
97     contTOT=contOCN+contICE;
98     %vertical divergence (air-sea fluxes or vertical advection)
99     zdivOCN=oceFWflx;
100     zdivICE=SIatmFW-oceFWflx;
101     %in virtual salt flux we omit :
102     if ~myparms.useRFWF; zdivOCN=0*zdivOCN; end;
103     %for deep ocean layer :
104     if kBudget>1; zdivOCN=-WVELMASS*myparms.rhoconst; end;
105     %
106     zdivTOT=zdivOCN+zdivICE;
107     %horizontal divergence (advection and ice diffusion)
108     hdivOCN=myparms.rhoconst*calc_UV_div(UVELMASS,VVELMASS,{'dh'}); %for 2D, already vertically integrated, fields
109     tmpU=(myparms.rhoi*DFxEHEFF+myparms.rhosn*DFxESNOW+myparms.rhoi*ADVxHEFF+myparms.rhosn*ADVxSNOW);
110     tmpV=(myparms.rhoi*DFyEHEFF+myparms.rhosn*DFyESNOW+myparms.rhoi*ADVyHEFF+myparms.rhosn*ADVySNOW);
111     hdivICE=calc_UV_div(tmpU,tmpV); %no dh needed here
112     hdivTOT=hdivOCN+hdivICE;
113    
114     %bottom pressure for comparison:
115     bp=myparms.rhoconst/9.81*PHIBOT;
116    
117     %compute global integrals:
118     %-------------------------
119     msk=mygrid.mskC(:,:,kBudget);
120     glo_vol_tot=calc_budget_mean_mask(contTOT,zdivTOT,hdivTOT,msk);
121     glo_vol_ocn=calc_budget_mean_mask(contOCN,zdivOCN,hdivOCN,msk);
122     glo_vol_ice=calc_budget_mean_mask(contICE,zdivICE,hdivICE,msk);
123     glo_bp=nansum(bp.*msk.*mygrid.RAC)/nansum(msk.*mygrid.RAC);
124    
125     %compute northern hemisphere integrals:
126     msk=mygrid.mskC(:,:,kBudget).*(mygrid.YC>0);
127     north_vol_tot=calc_budget_mean_mask(contTOT,zdivTOT,hdivTOT,msk);
128     north_vol_ocn=calc_budget_mean_mask(contOCN,zdivOCN,hdivOCN,msk);
129     north_vol_ice=calc_budget_mean_mask(contICE,zdivICE,hdivICE,msk);
130     north_bp=nansum(bp.*msk.*mygrid.RAC)/nansum(msk.*mygrid.RAC);
131    
132     %and southern hemisphere integrals:
133     msk=mygrid.mskC(:,:,kBudget).*(mygrid.YC<=0);
134     south_vol_tot=calc_budget_mean_mask(contTOT,zdivTOT,hdivTOT,msk);
135     south_vol_ocn=calc_budget_mean_mask(contOCN,zdivOCN,hdivOCN,msk);
136     south_vol_ice=calc_budget_mean_mask(contICE,zdivICE,hdivICE,msk);
137     south_bp=nansum(bp.*msk.*mygrid.RAC)/nansum(msk.*mygrid.RAC);
138    
139     %=======HEAT=======
140    
141     contOCN=myparms.rcp*THETA-myparms.rcp*AB_gT;
142     contICE=-myparms.flami*(SIheff*myparms.rhoi+SIhsnow*myparms.rhosn);
143     contTOT=contOCN+contICE;
144     %vertical divergence (air-sea fluxes or vertical adv/dif)
145     zdivOCN=TFLUX;
146     zdivICE=-(SItflux+TFLUX-TRELAX);
147     %in linear surface we omit :
148     if ~myparms.useNLFS; zdivOCN=zdivOCN-myparms.rcp*WTHMASS; end;
149     %in virtual salt flux we omit :
150     if ~myparms.useRFWF|~myparms.useNLFS; zdivICE=zdivICE+SIaaflux; end;
151     %working approach for real fresh water (?) and virtual salt flux
152     if 0; zdivICE=-oceQnet-SIatmQnt-myparms.flami*(SIsnPrcp-SIacSubl); end;
153     %for deep ocean layer :
154     if kBudget>1;
155     zdivOCN=-(ADVr_TH+DFrE_TH+DFrI_TH)./mygrid.RAC*myparms.rcp;
156     dd=mygrid.RF(kBudget); msk=mygrid.mskC(:,:,kBudget);
157     swfrac=0.62*exp(dd/0.6)+(1-0.62)*exp(dd/20);
158     if dd<-200; swfrac=0; end;
159     zdivOCN=zdivOCN+swfrac*oceQsw;%.*msk;
160     end;
161     %
162     zdivTOT=zdivOCN+zdivICE;
163     %horizontal divergence (advection and diffusion)
164     tmpU=myparms.rcp*(ADVx_TH+DFxE_TH); tmpV=myparms.rcp*(ADVy_TH+DFyE_TH);
165     hdivOCN=calc_UV_div(tmpU,tmpV);
166     tmpU=-myparms.flami*(myparms.rhoi*DFxEHEFF+myparms.rhosn*DFxESNOW+myparms.rhoi*ADVxHEFF+myparms.rhosn*ADVxSNOW);
167     tmpV=-myparms.flami*(myparms.rhoi*DFyEHEFF+myparms.rhosn*DFyESNOW+myparms.rhoi*ADVyHEFF+myparms.rhosn*ADVySNOW);
168     hdivICE=calc_UV_div(tmpU,tmpV); %no dh needed here
169     hdivTOT=hdivOCN+hdivICE;
170    
171     %compute global integrals:
172     %-------------------------
173     msk=mygrid.mskC(:,:,kBudget);
174     glo_heat_tot=calc_budget_mean_mask(contTOT,zdivTOT,hdivTOT,msk);
175     glo_heat_ocn=calc_budget_mean_mask(contOCN,zdivOCN,hdivOCN,msk);
176     glo_heat_ice=calc_budget_mean_mask(contICE,zdivICE,hdivICE,msk);
177    
178     %compute northern hemisphere integrals:
179     msk=mygrid.mskC(:,:,kBudget).*(mygrid.YC>0);
180     north_heat_tot=calc_budget_mean_mask(contTOT,zdivTOT,hdivTOT,msk);
181     north_heat_ocn=calc_budget_mean_mask(contOCN,zdivOCN,hdivOCN,msk);
182     north_heat_ice=calc_budget_mean_mask(contICE,zdivICE,hdivICE,msk);
183    
184     %and southern hemisphere integrals:
185     msk=mygrid.mskC(:,:,kBudget).*(mygrid.YC<=0);
186     south_heat_tot=calc_budget_mean_mask(contTOT,zdivTOT,hdivTOT,msk);
187     south_heat_ocn=calc_budget_mean_mask(contOCN,zdivOCN,hdivOCN,msk);
188     south_heat_ice=calc_budget_mean_mask(contICE,zdivICE,hdivICE,msk);
189    
190     %=======SALT=======
191    
192     contOCN=myparms.rhoconst*SALT-myparms.rhoconst*AB_gS;
193     contICE=myparms.SIsal0*myparms.rhoi*SIheff;
194     contTOT=contOCN+contICE;
195     %vertical divergence (air-sea fluxes or vertical adv/dif)
196     zdivOCN=SFLUX+oceSPflx;
197     zdivICE=-zdivOCN+SRELAX;
198     %in linear surface we omit :
199     if ~myparms.useNLFS; zdivOCN=zdivOCN-myparms.rhoconst*WSLTMASS; end;
200     %working approach for real fresh water (?) and virtual salt flux
201     if ~myparms.useRFWF|~myparms.useNLFS; zdivICE=-oceSflux; end;
202     %for deep ocean layer :
203     if kBudget>1;
204     zdivOCN=-(ADVr_SLT+DFrE_SLT+DFrI_SLT)./mygrid.RAC*myparms.rhoconst;
205     zdivOCN=zdivOCN+oceSPtnd;%.*msk;
206     end;
207     zdivTOT=zdivOCN+zdivICE;
208     %horizontal divergence (advection and diffusion)
209     tmpU=myparms.rhoconst*(ADVx_SLT+DFxE_SLT); tmpV=myparms.rhoconst*(ADVy_SLT+DFyE_SLT);
210     hdivOCN=calc_UV_div(tmpU,tmpV);
211     tmpU=myparms.SIsal0*(myparms.rhoi*DFxEHEFF+myparms.rhoi*ADVxHEFF);
212     tmpV=myparms.SIsal0*(myparms.rhoi*DFyEHEFF+myparms.rhoi*ADVyHEFF);
213     hdivICE=calc_UV_div(tmpU,tmpV); %no dh needed here
214     hdivTOT=hdivOCN+hdivICE;
215    
216     %compute global integrals:
217     %-------------------------
218     msk=mygrid.mskC(:,:,kBudget);
219     glo_salt_tot=calc_budget_mean_mask(contTOT,zdivTOT,hdivTOT,msk);
220     glo_salt_ocn=calc_budget_mean_mask(contOCN,zdivOCN,hdivOCN,msk);
221     glo_salt_ice=calc_budget_mean_mask(contICE,zdivICE,hdivICE,msk);
222    
223     %compute northern hemisphere integrals:
224     msk=mygrid.mskC(:,:,kBudget).*(mygrid.YC>0);
225     north_salt_tot=calc_budget_mean_mask(contTOT,zdivTOT,hdivTOT,msk);
226     north_salt_ocn=calc_budget_mean_mask(contOCN,zdivOCN,hdivOCN,msk);
227     north_salt_ice=calc_budget_mean_mask(contICE,zdivICE,hdivICE,msk);
228    
229     %and southern hemisphere integrals:
230     msk=mygrid.mskC(:,:,kBudget).*(mygrid.YC<=0);
231     south_salt_tot=calc_budget_mean_mask(contTOT,zdivTOT,hdivTOT,msk);
232     south_salt_ocn=calc_budget_mean_mask(contOCN,zdivOCN,hdivOCN,msk);
233     south_salt_ice=calc_budget_mean_mask(contICE,zdivICE,hdivICE,msk);
234     end;

  ViewVC Help
Powered by ViewVC 1.1.22