source: Roms_tools/Preprocessing_tools/add_ini_no3.m @ 2

Last change on this file since 2 was 1, checked in by cholod, 13 years ago

import Roms_Agrif

File size: 5.1 KB
Line 
1function 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%
54disp('Add_ini_no3: creating variable and attribute')
55%
56% open the grid file 
57%
58nc=netcdf(grdname);
59h=nc{'h'}(:);
60close(nc);
61%
62% open the initial file 
63%
64nc=netcdf(ininame,'write');
65theta_s = nc{'theta_s'}(:);
66if isempty(theta_s)
67  disp('Restart file')
68  theta_s=nc.theta_s(:);
69  theta_b=nc.theta_b(:);
70  hc=nc.hc(:);
71else
72  theta_b =  nc{'theta_b'}(:);
73  hc  =  nc{'hc'}(:);
74end
75N =  length(nc('s_rho'));
76%
77% open the oa file 
78%
79noa=netcdf(oaname);
80z=-noa{'Zno3'}(:);
81oatime=noa{'no3_time'}(:);
82tlen=length(oatime);
83%
84% Get the sigma depths
85%
86zroms=zlevs(h,0.*h,theta_s,theta_b,hc,N,'r');
87zmin=min(min(min(zroms)));
88zmax=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;
94addsurf=1;
95addbot=min(z)>zmin;
96if addsurf
97 z=[100;z];
98end
99if addbot
100 z=[z;-100000];
101end
102Nz=min(find(z<zmin));
103z=z(1:Nz);
104%
105% read in the initial file 
106%
107scrum_time = nc{'scrum_time'}(:);
108scrum_time = scrum_time / (24*3600);
109tinilen = length(scrum_time);
110redef(nc);
111nc{'NO3'} = ncdouble('time','s_rho','eta_rho','xi_rho');
112nc{'NO3'}.long_name = ncchar('Nitrate');
113nc{'NO3'}.long_name = 'Nitrate';
114nc{'NO3'}.units = ncchar('mMol N m-3');
115nc{'NO3'}.units = 'mMol N m-3';
116nc{'NO3'}.fields = ncchar('NO3, scalar, series');
117nc{'NO3'}.fields = 'NO3, scalar, series';
118endef(nc);
119%
120%  loop on initial time
121%
122for 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));
187end
188close(noa);
189%
190% Remove low values for oligotrophic areas
191%
192for l=1:tinilen
193  NO3=nc{'NO3'}(l,:,:,:);
194  NO3(NO3<NO3min)=0;
195  nc{'NO3'}(l,:,:,:)=NO3;
196end
197close(nc);
198return
Note: See TracBrowser for help on using the repository browser.