Commit 37b993a2825a0b226450d548e8ef9527e2b01614

Authored by Dany Dumont
1 parent 31a3ff4d
Exists in master

Nettoyage et resolution du bogue netcdf

nml/parameter.nml
... ... @@ -20,7 +20,7 @@
20 20  
21 21 Tm =4.5
22 22 Hs =1.5
23   -disp =0
  23 +disp =1
24 24 cont =1
25 25 /
26 26  
... ... @@ -37,8 +37,8 @@ cont =1
37 37 !--------------------------------------------------------------------
38 38 &model_parameter
39 39  
40   -nx =1000
41   -dx =5
  40 +nx =100
  41 +dx =10
42 42 cfl =1.0
43 43 name_sim ='test'
44 44 root = 'output/'
... ... @@ -91,7 +91,7 @@ swell_Hs =0.09
91 91 X_ice =50
92 92 cice =1
93 93 ice_thick=1
94   -hice =0
  94 +hice =0.8
95 95 hmax =0.8
96 96 Xh =250
97 97 strain_crit=3e-5
... ... @@ -131,7 +131,7 @@ nr =25
131 131 !--------------------------------------------------------------------
132 132 &itd_parameters
133 133  
134   -itd_scheme =1
  134 +itd_scheme =0
135 135 mu_itd =0.5
136 136 mincat_h =0.1
137 137 maxcat_h =10
... ...
src/attenuation.f90
1   -
2   -!____________________________________________________________________
3   -!
  1 +!--------------------------------------------------------------------
4 2 !DESCRIPTION: In this routine, the attenuation coefficient is
5 3 ! calculated according Kohout and Meylan (2008) and the
6 4 ! spectrum attenuation is computed.
7   -!
8   -!____________________________________________________________________
  5 +!--------------------------------------------------------------------
9 6  
10 7 !INTERFACE:
11 8 subroutine attenuation
12 9  
13   -!USES:
14   -use parameters
15   -
16   -implicit none
17   -
18   -double precision :: q11,q12,q13,q14,q15,q21,q22,q23,q24,q25
19   -double precision :: q31,q32,q33,q34,q35
20   -double precision, allocatable :: p1(:),p2(:),p3(:)
21   -double precision, allocatable :: att(:)
22   -
23   -allocate(p1 (nf))
24   -allocate(p2 (nf))
25   -allocate(p3 (nf))
26   -allocate(att(nf))
27   -
28   -
29   -q11 = -0.00000777
30   -q12 = 0.00032080
31   -q13 = -0.00437542
32   -q14 = 0.02047559
33   -q15 = -0.01356537
34   -
35   -q21 = 0.00003635
36   -q22 = -0.00153484
37   -q23 = 0.02121709
38   -q24 = -0.09289399
39   -q25 = -0.03693082
40   -
41   -q31 = -0.00004509
42   -q32 = 0.00214484
43   -q33 = -0.03663425
44   -q34 = 0.26065369
45   -q35 = -0.62474085
46   -
47   -p1 = q11*T**4 + q12*T**3 + q13*T**2 + q14*T + q15
48   -p2 = q21*T**4 + q22*T**3 + q23*T**2 + q24*T + q25
49   -p3 = q31*T**4 + q32*T**3 + q33*T**2 + q34*T + q35
50   -
51   -alpha(i,:)=-1*(p1*h(i)**2 + p2*h(i) + p3)
52   -
53   -do ii=1,nf
54   - if ( alpha(i,ii).lt.0 )then
55   - alpha(i,ii)=0d0
56   - end if
57   -end do
58   -
59   -att=C_ice(i)*alpha(i,:)/(Dave(i)+3e-14)
60   -!S_ice(n,i,1:nf)=-att*E(n,i,1:nf)
61   -S_ice(i,:,n)=-att*E(i,:,n)
62   -E(i,:,n)=E(i,:,n)*exp(-att*Cg*dt)
  10 + !USES:
  11 + use parameters
  12 +
  13 + implicit none
  14 +
  15 + double precision :: q11,q12,q13,q14,q15,q21,q22,q23,q24,q25
  16 + double precision :: q31,q32,q33,q34,q35
  17 + double precision, allocatable :: p1(:),p2(:),p3(:)
  18 + double precision, allocatable :: att(:)
  19 +
  20 + allocate(p1 (nf))
  21 + allocate(p2 (nf))
  22 + allocate(p3 (nf))
  23 + allocate(att(nf))
  24 +
  25 + q11 = -0.00000777
  26 + q12 = 0.00032080
  27 + q13 = -0.00437542
  28 + q14 = 0.02047559
  29 + q15 = -0.01356537
  30 +
  31 + q21 = 0.00003635
  32 + q22 = -0.00153484
  33 + q23 = 0.02121709
  34 + q24 = -0.09289399
  35 + q25 = -0.03693082
  36 +
  37 + q31 = -0.00004509
  38 + q32 = 0.00214484
  39 + q33 = -0.03663425
  40 + q34 = 0.26065369
  41 + q35 = -0.62474085
  42 +
  43 + p1 = q11*T**4 + q12*T**3 + q13*T**2 + q14*T + q15
  44 + p2 = q21*T**4 + q22*T**3 + q23*T**2 + q24*T + q25
  45 + p3 = q31*T**4 + q32*T**3 + q33*T**2 + q34*T + q35
  46 +
  47 + alpha(i,:)=-1*(p1*h(i)**2 + p2*h(i) + p3)
  48 +
  49 + do ii=1,nf
  50 + if ( alpha(i,ii).lt.0 )then
  51 + alpha(i,ii)=0d0
  52 + end if
  53 + end do
  54 +
  55 + att=C_ice(i)*alpha(i,:)/(Dave(i)+3e-14)
  56 + !S_ice(n,i,1:nf)=-att*E(n,i,1:nf)
  57 + S_ice(i,:,n)=-att*E(i,:,n)
  58 + E(i,:,n)=E(i,:,n)*exp(-att*Cg*dt)
63 59  
64 60 end subroutine attenuation
... ...
src/ice_fracture.f90
  1 +!--------------------------------------------------------------------
  2 +!DESCRIPTION: This routine does the wave-induced floe breaking by
  3 +! modifying the floe size distribution.
  4 +!--------------------------------------------------------------------
  5 +
  6 +
  7 +! There is a major issue with as the integral of the distribution
  8 +! is not consistent with the ice area and the redistribution is not
  9 +! performed conservatively.
1 10  
2 11 subroutine ice_fracture
3 12  
... ... @@ -6,20 +15,18 @@ subroutine ice_fracture
6 15 !local parameters
7 16 implicit none
8 17  
9   - integer :: Nband=500
  18 + !integer :: Nband=500
10 19 integer :: freq1
11 20 integer :: freq2
12 21 integer :: ninterp=1
13 22  
14   - double precision :: Y=5.5e9
15   - double precision :: TT=13
16   - double precision :: poisson=0.3
  23 + !double precision :: TT=13
17 24 double precision :: m_0
18 25 double precision :: m_0_stress
19 26 double precision :: m_0_amp
20 27 double precision :: m_2_amp
21   - double precision :: Tw
22   - double precision :: kw
  28 + !double precision :: Tw
  29 + !double precision :: kw
23 30 double precision :: Lmin
24 31  
25 32 double precision, allocatable :: P(:)
... ... @@ -70,9 +77,13 @@ subroutine ice_fracture
70 77  
71 78 do iii=1,nh
72 79 ! Lmin according Mellor 1986
73   - Lmin=pi*0.5*((5e6*middle_h_cat(iii)**3)/(3*10*(1-0.3**2)))**0.25
74   -
75   - W=(9.81*kice(iii,:))/omega**2
  80 + !Lmin=pi*0.5*((Y*middle_h_cat(iii)**3)/ &
  81 + ! (3*rhoi*g*(1-nu**2)))**0.25
  82 + ! Lmin corrected from Mellor (1986)
  83 + Lmin=pi*0.25*((Y*middle_h_cat(iii)**3)/ &
  84 + (3*rhoi*g*(1-nu**2)))**0.25
  85 +
  86 + W=(g*kice(iii,:))/omega**2
76 87 wavelength=2*pi/kice(iii,:)
77 88  
78 89 strain=middle_h_cat(iii)*0.5*kice(iii,:)**2*W
... ... @@ -82,8 +93,8 @@ subroutine ice_fracture
82 93  
83 94 do ii=1,nr
84 95  
85   - if ( 2d0*floe_cat(ii+1).lt.wavelength(nf) .or. &
86   - 2d0*floe_cat(ii+1).gt.wavelength(1)) then
  96 + if ( 2d0*floe_cat(ii+1).lt.wavelength(nf) .or. &
  97 + 2d0*floe_cat(ii+1).gt.wavelength( 1) ) then
87 98 P_cat(ii)=0d0
88 99 else
89 100 do while ( wavelength(freq1).lt.2d0*floe_cat(ii+1) .or. &
... ...
src/initialization.f90
... ... @@ -42,6 +42,7 @@ subroutine initialization
42 42 end do
43 43  
44 44 write(*,*) ' * dt = ',dt,'s'
  45 + write(*,*) ' * dx = ',dx,'m'
45 46 if ( disp.eq.1 ) then
46 47 write(*,*) ' * dispersion is on'
47 48 else
... ... @@ -49,7 +50,7 @@ subroutine initialization
49 50 end if
50 51  
51 52 if ( cont.eq.1 ) then
52   - write(*,*) ' * wave energy is continuously injected'
  53 + write(*,*) ' * wave energy continuously injected'
53 54 else
54 55 write(*,*) ' * only a single wave energy pulse'
55 56 end if
... ... @@ -62,6 +63,8 @@ subroutine initialization
62 63  
63 64 if ( init_spec.eq.1 ) then
64 65  
  66 + write(*,*) ' * Jonswap spectrum'
  67 +
65 68 ! JONSWAP spectrum
66 69 do i=1,nf
67 70 if ( freq(i).le.freq_s ) then
... ... @@ -72,18 +75,22 @@ subroutine initialization
72 75 end do
73 76  
74 77 Gf=gamma_s**(exp((-(freq-freq_s)**2)/(2*sigma_s**2*freq_s**2)))
75   - PM=alpha_s*Hs**2*(freq_s**4/freq**5)*exp(-beta_s*(freq_s/freq)**4)
  78 + PM=alpha_s*Hs**2*(freq_s**4/freq**5)* &
  79 + exp(-beta_s*(freq_s/freq)**4)
76 80 Ei=Gf*PM
77   -
  81 +
78 82 else if ( init_spec.eq.0 ) then
79 83  
  84 + write(*,*) ' * Bretschneider spectrum'
  85 +
80 86 ! Bretschneider spectrum
81   - Ei=(1.25*Hs**2*(1/freq)**5)/(8*pi*Tm**4)* &
  87 + Ei=(1.25*Hs**2*(1/freq)**5)/(8*pi*Tm**4)* &
82 88 exp(-1.25*((1/freq)/Tm)**4)
83 89  
84 90 else
85 91  
86   - Ei=1/(0.01*sqrt(2*pi))*exp(-(omega-2*pi/swell_T)**2/(2*0.01**2))
  92 + Ei=1/(0.01*sqrt(2*pi))* &
  93 + exp(-(omega-2*pi/swell_T)**2/(2*0.01**2))
87 94 Ei=(swell_Hs/4d0)**2*Ei/(sum(Ei*domega))
88 95  
89 96 end if
... ... @@ -94,12 +101,12 @@ subroutine initialization
94 101  
95 102 S_ice(1,:,1)=0
96 103 x1=floor(X_ice/dx) !find in which grid bin is the ice edge
97   - C_ice( 1:x1)=0 !ice concentration is 0 before ice edge!
  104 + C_ice( 1:x1)=0 !ice concentration is 0 before ice edge
98 105 C_ice(x1:nx)=cice !ice concentration in the transect
99 106  
100   - Dave( 1:x1)=0 !initalize mean floe size before ice edge
101   - Dmax( 1:x1)=0 !initialize max floe size before ice edge
102   - Dave(x1:nx)=maxfloe !initalize mean floe size in ice transect
  107 + Dave( 1:x1)=0 !initialize mean floe size before ice edge
  108 + Dmax( 1:x1)=0 !initialize max floe size before ice edge
  109 + Dave(x1:nx)=maxfloe !initialize mean floe size in ice transect
103 110 Dmax(x1:nx)=maxfloe !initialize max floe size in ice transect
104 111  
105 112 if ( ice_thick.eq.0 ) then
... ... @@ -118,10 +125,10 @@ subroutine initialization
118 125  
119 126 !--------------------- FLOE SIZE DISTRIBUTION -----------------------
120 127  
121   - res=abs(minfloe-maxfloe)/nedge
  128 + dr=abs(minfloe - maxfloe)/nedge
122 129 floe_cat(1)=minfloe
123 130 do i=2,nedge
124   - floe_cat(i)=floe_cat(1)+i*res
  131 + floe_cat(i)=floe_cat(1)+i*dr
125 132 end do
126 133  
127 134 do i=1,nr
... ... @@ -134,18 +141,17 @@ subroutine initialization
134 141 end do
135 142  
136 143 Dave(x1:nx)=middle_floe_cat(nr)
137   - !initialize max floe size in ice transect
138 144  
139 145 !--------------------- ICE THICKNESS DISTRIBUTION -------------------
140 146  
141   - resh=abs(mincat_h-maxcat_h)/(nedge_h)
  147 + dh=abs(mincat_h - maxcat_h)/nedge_h
142 148 h_cat(1)=mincat_h
143 149 do i=2,nedge_h
144   - h_cat(i)=h_cat(1)+i*resh
  150 + h_cat(i)=h_cat(1) + i*dh
145 151 end do
146 152  
147 153 do i=1,nh
148   - middle_h_cat(i)=(h_cat(i)+h_cat(i+1))*0.5
  154 + middle_h_cat(i)=0.5*(h_cat(i) + h_cat(i+1))
149 155 end do
150 156  
151 157 itd(:,:)=0d0
... ... @@ -170,7 +176,6 @@ subroutine initialization
170 176  
171 177 end if
172 178  
173   -
174 179 !--------------------- STATISTICS -----------------------------------
175 180  
176 181 h_sign(1,1)=Hs
... ...
src/main.f90
1   -!
2   -!____________________________________________________________________
3   -!
  1 +!--------------------------------------------------------------------
4 2 !DESCRIPTION: This is the main program of WIM. This routine merely
5 3 ! calls other subroutines and do the main time loop.
6 4 ! It also contains the subroutine progress which display
7 5 ! a progress bar in the terminal while the model is
8 6 ! running.
9   -!____________________________________________________________________
  7 +!--------------------------------------------------------------------
10 8  
11 9 !INTERFACE:
12 10 program wim2
... ... @@ -113,7 +111,7 @@ subroutine message(info)
113 111 endif
114 112  
115 113 if ( info.eq.3 ) then
116   - write(*,*)' Allocate memory... '
  114 + write(*,*)' Allocate memory ... '
117 115 endif
118 116  
119 117 if ( info.eq.4 ) then
... ... @@ -129,12 +127,9 @@ subroutine message(info)
129 127 endif
130 128  
131 129 if ( info.eq.7 ) then
132   - write(*,*)' done!'
  130 + write(*,*)' done!'
133 131 endif
134 132  
135 133 end subroutine message
136   -!____________________________________________________________________
137   -
138 134  
139 135 end program wim2
140   -
... ...
src/parameters.f90
1   -
2   -!____________________________________________________________________
3   -!
  1 +!--------------------------------------------------------------------
4 2 !DESCRIPTION: In this routine, the shared parameters for other
5 3 ! subroutines are declared. It also contains the
6 4 ! subroutine read_namelist which read parameters value
7 5 ! from the namelist and the subroutine array_allocation
8 6 ! for memory allocation.
9   -!____________________________________________________________________
  7 +!--------------------------------------------------------------------
10 8  
11   -!INTERFACE:
12 9 module parameters
13 10  
14 11 implicit none
15 12  
16 13 !NAME OF OUTPUTS FILES
17   -character(len=100) :: root='output/' ! Output directory
18   -character(len=100) :: name_sim='test'! Name of the simulation
19   -character(len=100) :: namefile
20   -character(len=100) :: dirout
  14 +character(len=100) :: root='output/' ! Output directory
  15 +character(len=100) :: name_sim='test'! Name of the simulation
  16 +character(len=100) :: namefile
  17 +character(len=100) :: dirout
21 18  
22 19 !DUMMIES
23   -integer :: i,ii,iii
24   -integer :: j,jj,jjj
25   -integer :: n
  20 +integer :: i,ii,iii
  21 +integer :: j,jj,jjj
  22 +integer :: n
26 23  
27   -!GLOBAL PARAMETERS
28   -double precision :: g=9.81 ! Gravitationnal acceleration
29   -double precision :: pi=3.1416 ! Pi
  24 +! physical and useful quantities
  25 +real :: g=9.81 ! Gravitationnal acc. [m s-2]
  26 +real :: Y=5.5e9 ! Young's modulus [Pa]
  27 +real :: nu=0.3 ! Poisson's ratio [-]
  28 +real :: rhoi=925.0 ! Sea ice density [kg m-3]
  29 +double precision :: pi=3.141592653 ! Pi
30 30  
31 31 !SPECTRUM
32   -integer :: nf=60 ! Nb of frequency bins
33   -integer :: init_spec=1
34   -
35   -double precision :: alpha_s=8.1e-3
36   -double precision :: beta_s=1.25
37   -double precision :: Tmin=2.5 ! Minimum period
38   -double precision :: Tmax=20 ! Maximum period
39   -double precision :: Tm=6 ! Peak period
40   -double precision :: Hs=1 ! Significant Height
41   -double precision :: gamma_s=3.3
42   -double precision :: freq_s ! Peak frequency
43   -double precision :: domega ! Frequency increment
44   -double precision :: swell_T=8
45   -double precision :: swell_Hs=1
  32 +integer :: nf=60 ! Nb of frequency bins
  33 +integer :: init_spec=1
  34 +
  35 +double precision :: alpha_s=8.1e-3 ! Jonswap parameter
  36 +double precision :: beta_s=1.25 ! Jonswap parameter
  37 +double precision :: gamma_s=3.3 ! Jonswap parameter
  38 +double precision :: Tmin=2.5 ! Minimum period
  39 +double precision :: Tmax=20 ! Maximum period
  40 +double precision :: Tm=6 ! Peak period [s]
  41 +double precision :: Hs=1 ! Significant wave height [m]
  42 +double precision :: freq_s ! Peak frequency
  43 +double precision :: domega ! Frequency increment
  44 +double precision :: swell_T=8
  45 +double precision :: swell_Hs=1
46 46  
47 47 double precision, allocatable :: E(:,:,:) ! Energy spectrum
48 48 double precision, allocatable :: Ei(:) ! Initial spectrum
... ... @@ -67,7 +67,7 @@ integer :: nx=10 ! Nb of spatial bins
67 67 integer :: nt ! Duration of the simulation
68 68  
69 69 real :: dt ! Time step
70   -double precision :: dx=500 ! Grid resolution
  70 +real :: dx=500 ! Grid resolution
71 71 double precision :: cfl=1 ! Courant
72 72  
73 73 double precision, allocatable :: x_axis(:) ! Vector of grid width
... ... @@ -106,7 +106,7 @@ integer :: nedge
106 106  
107 107 double precision :: minfloe=5
108 108 double precision :: maxfloe=500
109   -double precision :: res
  109 +double precision :: dr
110 110  
111 111 double precision, allocatable :: fsd(:,:,:,:)
112 112 double precision, allocatable :: floe_cat(:)
... ... @@ -119,8 +119,8 @@ integer :: nedge_h
119 119  
120 120 double precision :: mu_itd=1.2
121 121 double precision :: mincat_h=0.1
122   -double precision :: maxcat_h=10
123   -double precision :: resh
  122 +double precision :: maxcat_h=5.0
  123 +double precision :: dh
124 124  
125 125 double precision, allocatable :: itd(:,:)
126 126 double precision, allocatable :: middle_h_cat(:)
... ... @@ -135,7 +135,7 @@ contains
135 135  
136 136 subroutine read_namelist
137 137  
138   - namelist /spectrum_parameters/ nf,alpha_s,beta_s, &
  138 + namelist /spectrum_parameters/ nf,alpha_s,beta_s, &
139 139 Tmin,Tmax,gamma_s,init_spec, &
140 140 swell_T,swell_Hs
141 141  
... ... @@ -195,7 +195,7 @@ subroutine array_allocation
195 195 allocate(C_ice (nx))
196 196 allocate(x_axis (nx))
197 197  
198   - domega=(2*pi/Tmin-2*pi/Tmax)/(nf)
  198 + domega=(2*pi/Tmin - 2*pi/Tmax)/nf
199 199 dwl=(g*Tmax**2/(2*pi)-g*Tmin**2/(2*pi))/(nf)
200 200 ii=1
201 201 do i=0,nf
... ... @@ -222,11 +222,11 @@ subroutine array_allocation
222 222 nt=ceiling(nx/minval(CN))
223 223  
224 224 allocate(E (nx,nf,nt))
225   - allocate(Dmax (nx))
226   - allocate(Dave (nx))
227   - allocate(time (nt))
  225 + allocate(Dmax (nx ))
  226 + allocate(Dave (nx ))
  227 + allocate(time (nt ))
228 228 allocate(S_ice (nx,nf,nt))
229   - allocate(h_sign(nx,nt))
  229 + allocate(h_sign(nx,nt ))
230 230  
231 231 nedge=nr+1
232 232  
... ...
src/write_output.f90
... ... @@ -22,11 +22,13 @@ subroutine write_output
22 22 include 'netcdf.inc'
23 23  
24 24 integer :: ncid
25   - integer :: omID
26 25 integer :: xID
  26 + integer :: fID
  27 + integer :: rID
  28 + integer :: hID
27 29 integer :: tID
28 30 integer :: EID
29   - integer :: om_varID
  31 + integer :: f_varID
30 32 integer :: x_varID
31 33 integer :: t_varID
32 34 integer :: stat
... ... @@ -34,12 +36,10 @@ subroutine write_output
34 36 integer :: Dave_varID
35 37 integer :: c_varID
36 38 integer :: h_varID
37   - integer :: floeID
38 39 integer :: fsdID
39   - integer :: floe_varID
  40 + integer :: r_varID
40 41 integer :: SrcID
41 42 integer :: HsID
42   - integer :: hID
43 43 integer :: hcat_varID
44 44 integer :: itdID
45 45  
... ... @@ -47,7 +47,7 @@ subroutine write_output
47 47 call handle_err(stat)
48 48  
49 49 ! dimensions
50   - stat=nf90_def_dim(ncid,"omega",nf,omID)
  50 + stat=nf90_def_dim(ncid,"omega",nf,fID)
51 51 call handle_err(stat)
52 52  
53 53 stat=nf90_def_dim(ncid,"x" ,nx,xID)
... ... @@ -56,14 +56,14 @@ subroutine write_output
56 56 stat=nf90_def_dim(ncid,"time",nt,tID)
57 57 call handle_err(stat)
58 58  
59   - stat=nf90_def_dim(ncid,"dcat",nr,floeID)
  59 + stat=nf90_def_dim(ncid,"dcat",nr,rID)
60 60 call handle_err(stat)
61 61  
62 62 stat=nf90_def_dim(ncid,"hcat",nh,hID)
63 63 call handle_err(stat)
64 64  
65 65 ! coordinates
66   - stat=nf90_def_var(ncid,"omega",nf90_double,omID,om_varID)
  66 + stat=nf90_def_var(ncid,"omega",nf90_double,fID,f_varID)
67 67 call handle_err(stat)
68 68  
69 69 stat=nf90_def_var(ncid,"x",nf90_double,xID,x_varID)
... ... @@ -75,15 +75,15 @@ subroutine write_output
75 75 stat=nf90_def_var(ncid,"hcat",nf90_double,hID,hcat_varID)
76 76 call handle_err(stat)
77 77  
78   - stat=nf90_def_var(ncid,"dcat",nf90_double,floeID,floe_varID)
  78 + stat=nf90_def_var(ncid,"dcat",nf90_double,rID,r_varID)
79 79 call handle_err(stat)
80 80  
81 81 ! variables
82 82 ! (x, f, t)
83   - stat=nf90_def_var(ncid,"E",nf90_double,(/xID,omID,tID/),EID)
  83 + stat=nf90_def_var(ncid,"E",nf90_double,(/xID,fID,tID/),EID)
84 84 call handle_err(stat)
85 85  
86   - stat=nf90_def_var(ncid,"Sice",nf90_double,(/xID,omID,tID/),SrcID)
  86 + stat=nf90_def_var(ncid,"Sice",nf90_double,(/xID,fID,tID/),SrcID)
87 87 call handle_err(stat)
88 88  
89 89 ! (x)
... ... @@ -100,7 +100,7 @@ subroutine write_output
100 100 call handle_err(stat)
101 101  
102 102 ! (x, r, h, t)
103   - stat=nf90_def_var(ncid,"fsd",nf90_double,(/xID,floeID,hID,tID/),fsdID)
  103 + stat=nf90_def_var(ncid,"fsd",nf90_double,(/xID,rID,hID,tID/),fsdID)
104 104 call handle_err(stat)
105 105  
106 106 ! (x, t)
... ... @@ -113,37 +113,58 @@ subroutine write_output
113 113  
114 114 ! set attributes
115 115 stat=nf_put_att_text(ncid,xID,'units',10,'km')
  116 + call handle_err(stat)
116 117 stat=nf_put_att_text(ncid,xID,'axis',1,'X')
117   - stat=nf_put_att_text(ncid,xID,'long_name',7,'Distance')
  118 + call handle_err(stat)
  119 + stat=nf_put_att_text(ncid,xID,'long_name',10,'Distance')
  120 + call handle_err(stat)
118 121  
119   - stat=nf_put_att_text(ncid,omID,'units',5,'Hz')
120   - stat=nf_put_att_text(ncid,omID,'long_name',7,'Frequency')
  122 + stat=nf_put_att_text(ncid,rID,'units',5,'m')
  123 + call handle_err(stat)
  124 + stat=nf_put_att_text(ncid,rID,'long_name',10,'Floe size')
  125 + call handle_err(stat)
  126 +
  127 + stat=nf_put_att_text(ncid,fID,'units',5,'Hz')
  128 + call handle_err(stat)
  129 + stat=nf_put_att_text(ncid,fID,'long_name',10,'Frequency')
  130 + call handle_err(stat)
121 131  
122 132 stat=nf_put_att_text(ncid,hID,'units',5,'m')
  133 + call handle_err(stat)
123 134 stat=nf_put_att_text(ncid,hID,'axis',1,'Z')
124   - stat=nf_put_att_text(ncid,hID,'long_name',7,'Ice thickness')
  135 + call handle_err(stat)
  136 + stat=nf_put_att_text(ncid,hID,'long_name',20,'Ice thickness')
  137 + call handle_err(stat)
125 138  
126 139 stat=nf_put_att_text(ncid,tID,'units',7,'hours')
  140 + call handle_err(stat)
127 141 stat=nf_put_att_text(ncid,tID,'axis',1,'T')
  142 + call handle_err(stat)
128 143 stat=nf_put_att_text(ncid,tID,'long_name',4,'Time')
  144 + call handle_err(stat)
129 145  
130 146 stat=nf_put_att_text(ncid,HsID,'units',1,'m')
  147 + call handle_err(stat)
131 148 stat=nf_put_att_text(ncid,HsID,'long_name',23, &
132   - 'Significant Wave Heigth')
  149 + 'Significant Wave Height')
  150 + call handle_err(stat)
133 151  
134 152 stat=nf_put_att_text(ncid,EID,'units',10,'m2 s rad-1')
  153 + call handle_err(stat)
135 154 stat=nf_put_att_text(ncid,EID,'long_name',50, &
136 155 'Sea surface wave variance spectral density')
  156 + call handle_err(stat)
137 157  
138 158 stat=nf_put_att_text(ncid,fsdID,'units',10,'m m-1')
139 159 stat=nf_put_att_text(ncid,fsdID,'long_name',50, &
140 160 'Floe size and thickness distribution')
  161 +
141 162 ! end file definition
142 163 stat=nf90_enddef(ncid)
143 164 call handle_err(stat)
144 165  
145 166 ! populate variables
146   - stat=nf90_put_var(ncid,om_varID,omega)
  167 + stat=nf90_put_var(ncid,f_varID,omega)
147 168 call handle_err(stat)
148 169  
149 170 stat=nf90_put_var(ncid,EID,E)
... ... @@ -173,7 +194,7 @@ subroutine write_output
173 194 stat=nf90_put_var(ncid,fsdID,fsd)
174 195 call handle_err(stat)
175 196  
176   - stat=nf90_put_var(ncid,floe_varID,middle_floe_cat)
  197 + stat=nf90_put_var(ncid,r_varID,middle_floe_cat)
177 198 call handle_err(stat)
178 199  
179 200 stat=nf90_put_var(ncid,hcat_varID,middle_h_cat)
... ...