1 | function add_ini_no3(ininame,grdname,oaname,cycle,NO3min); |
---|
2 | %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
---|
3 | % |
---|
4 | % function [longrd,latgrd,no3]=add_ini_no3(ininame,grdname,... |
---|
5 | % seas_datafile,ann_datafile,... |
---|
6 | % cycle); |
---|
7 | % |
---|
8 | % Add nitrate (mMol N m-3) in a ROMS initial file. |
---|
9 | % take seasonal data for the upper levels and annual data for the |
---|
10 | % lower levels. |
---|
11 | % do a temporal interpolation to have values at initial |
---|
12 | % time. |
---|
13 | % |
---|
14 | % input: |
---|
15 | % |
---|
16 | % ininame : roms initial file to process (netcdf) |
---|
17 | % grdname : roms grid file (netcdf) |
---|
18 | % seas_datafile : regular longitude - latitude - z seasonal data |
---|
19 | % file used for the upper levels (netcdf) |
---|
20 | % ann_datafile : regular longitude - latitude - z annual data |
---|
21 | % file used for the lower levels (netcdf) |
---|
22 | % cycle : time length (days) of climatology cycle (ex:360 for |
---|
23 | % annual cycle) - 0 if no cycle. |
---|
24 | % |
---|
25 | % output: |
---|
26 | % |
---|
27 | % [longrd,latgrd,no3] : surface field to plot (as an illustration) |
---|
28 | % |
---|
29 | % Further Information: |
---|
30 | % http://www.brest.ird.fr/Roms_tools/ |
---|
31 | % |
---|
32 | % This file is part of ROMSTOOLS |
---|
33 | % |
---|
34 | % ROMSTOOLS is free software; you can redistribute it and/or modify |
---|
35 | % it under the terms of the GNU General Public License as published |
---|
36 | % by the Free Software Foundation; either version 2 of the License, |
---|
37 | % or (at your option) any later version. |
---|
38 | % |
---|
39 | % ROMSTOOLS is distributed in the hope that it will be useful, but |
---|
40 | % WITHOUT ANY WARRANTY; without even the implied warranty of |
---|
41 | % MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the |
---|
42 | % GNU General Public License for more details. |
---|
43 | % |
---|
44 | % You should have received a copy of the GNU General Public License |
---|
45 | % along with this program; if not, write to the Free Software |
---|
46 | % Foundation, Inc., 59 Temple Place, Suite 330, Boston, |
---|
47 | % MA 02111-1307 USA |
---|
48 | % |
---|
49 | % Copyright (c) 2001-2006 by Pierrick Penven |
---|
50 | % e-mail:Pierrick.Penven@ird.fr |
---|
51 | % |
---|
52 | %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% |
---|
53 | % |
---|
54 | disp('Add_ini_no3: creating variable and attribute') |
---|
55 | % |
---|
56 | % open the grid file |
---|
57 | % |
---|
58 | nc=netcdf(grdname); |
---|
59 | h=nc{'h'}(:); |
---|
60 | close(nc); |
---|
61 | % |
---|
62 | % open the initial file |
---|
63 | % |
---|
64 | nc=netcdf(ininame,'write'); |
---|
65 | theta_s = nc{'theta_s'}(:); |
---|
66 | if isempty(theta_s) |
---|
67 | disp('Restart file') |
---|
68 | theta_s=nc.theta_s(:); |
---|
69 | theta_b=nc.theta_b(:); |
---|
70 | hc=nc.hc(:); |
---|
71 | else |
---|
72 | theta_b = nc{'theta_b'}(:); |
---|
73 | hc = nc{'hc'}(:); |
---|
74 | end |
---|
75 | N = length(nc('s_rho')); |
---|
76 | % |
---|
77 | % open the oa file |
---|
78 | % |
---|
79 | noa=netcdf(oaname); |
---|
80 | z=-noa{'Zno3'}(:); |
---|
81 | oatime=noa{'no3_time'}(:); |
---|
82 | tlen=length(oatime); |
---|
83 | % |
---|
84 | % Get the sigma depths |
---|
85 | % |
---|
86 | zroms=zlevs(h,0.*h,theta_s,theta_b,hc,N,'r'); |
---|
87 | zmin=min(min(min(zroms))); |
---|
88 | zmax=max(max(max(zroms))); |
---|
89 | % |
---|
90 | % Check if the min z level is below the min sigma level |
---|
91 | % (if not add a deep layer) |
---|
92 | % |
---|
93 | %addsurf=max(z)<zmax; |
---|
94 | addsurf=1; |
---|
95 | addbot=min(z)>zmin; |
---|
96 | if addsurf |
---|
97 | z=[100;z]; |
---|
98 | end |
---|
99 | if addbot |
---|
100 | z=[z;-100000]; |
---|
101 | end |
---|
102 | Nz=min(find(z<zmin)); |
---|
103 | z=z(1:Nz); |
---|
104 | % |
---|
105 | % read in the initial file |
---|
106 | % |
---|
107 | scrum_time = nc{'scrum_time'}(:); |
---|
108 | scrum_time = scrum_time / (24*3600); |
---|
109 | tinilen = length(scrum_time); |
---|
110 | redef(nc); |
---|
111 | nc{'NO3'} = ncdouble('time','s_rho','eta_rho','xi_rho'); |
---|
112 | nc{'NO3'}.long_name = ncchar('Nitrate'); |
---|
113 | nc{'NO3'}.long_name = 'Nitrate'; |
---|
114 | nc{'NO3'}.units = ncchar('mMol N m-3'); |
---|
115 | nc{'NO3'}.units = 'mMol N m-3'; |
---|
116 | nc{'NO3'}.fields = ncchar('NO3, scalar, series'); |
---|
117 | nc{'NO3'}.fields = 'NO3, scalar, series'; |
---|
118 | endef(nc); |
---|
119 | % |
---|
120 | % loop on initial time |
---|
121 | % |
---|
122 | for l=1:tinilen |
---|
123 | disp(['time index: ',num2str(l),' of total: ',num2str(tinilen)]) |
---|
124 | % |
---|
125 | % get data time indices and weights for temporal interpolation |
---|
126 | % |
---|
127 | if cycle~=0 |
---|
128 | modeltime=mod(scrum_time(l),cycle); |
---|
129 | else |
---|
130 | modeltime=scrum_time; |
---|
131 | end |
---|
132 | l1=find(modeltime==oatime); |
---|
133 | if isempty(l1) |
---|
134 | disp('temporal interpolation') |
---|
135 | l1=max(find(oatime<modeltime)); |
---|
136 | time1=oatime(l1); |
---|
137 | if isempty(l1) |
---|
138 | if cycle~=0 |
---|
139 | l1=tlen; |
---|
140 | time1=oatime(l1)-cycle; |
---|
141 | else |
---|
142 | error('No previous time in the dataset') |
---|
143 | end |
---|
144 | end |
---|
145 | l2=min(find(oatime>modeltime)); |
---|
146 | time2=oatime(l2); |
---|
147 | if isempty(l2) |
---|
148 | if cycle~=0 |
---|
149 | l2=1; |
---|
150 | time2=oatime(l2)+cycle; |
---|
151 | else |
---|
152 | error('No posterious time in the dataset') |
---|
153 | end |
---|
154 | end |
---|
155 | disp(['Initialisation time: ',num2str(modeltime),... |
---|
156 | ' - Time 1: ',num2str(time1),... |
---|
157 | ' - Time 2: ',num2str(time2)]) |
---|
158 | cff1=(modeltime-time2)/(time1-time2); |
---|
159 | cff2=(time1-modeltime)/(time1-time2); |
---|
160 | else |
---|
161 | cff1=1; |
---|
162 | l2=l1; |
---|
163 | cff2=0; |
---|
164 | end |
---|
165 | % |
---|
166 | % interpole the seasonal dataset on the horizontal roms grid |
---|
167 | % |
---|
168 | disp(['Add_ini_no3: vertical interpolation']) |
---|
169 | var=squeeze(noa{'NO3'}(l1,:,:,:)); |
---|
170 | if addsurf |
---|
171 | var=cat(1,var(1,:,:),var); |
---|
172 | end |
---|
173 | if addbot |
---|
174 | var=cat(1,var,var(end,:,:)); |
---|
175 | end |
---|
176 | var2=squeeze(noa{'NO3'}(l2,:,:,:)); |
---|
177 | if addsurf |
---|
178 | var2=cat(1,var2(1,:,:),var2); |
---|
179 | end |
---|
180 | if addbot |
---|
181 | var2=cat(1,var2,var2(end,:,:)); |
---|
182 | end |
---|
183 | var=cff1*var + cff2*var2; |
---|
184 | zeta = squeeze(nc{'zeta'}(l,:,:)); |
---|
185 | zroms=zlevs(h,zeta,theta_s,theta_b,hc,N,'r'); |
---|
186 | nc{'NO3'}(l,:,:,:)=ztosigma(flipdim(var,1),zroms,flipud(z)); |
---|
187 | end |
---|
188 | close(noa); |
---|
189 | % |
---|
190 | % Remove low values for oligotrophic areas |
---|
191 | % |
---|
192 | for l=1:tinilen |
---|
193 | NO3=nc{'NO3'}(l,:,:,:); |
---|
194 | NO3(NO3<NO3min)=0; |
---|
195 | nc{'NO3'}(l,:,:,:)=NO3; |
---|
196 | end |
---|
197 | close(nc); |
---|
198 | return |
---|