[ English | Japanese ] [ Gtool5 Documents | Gtool5 Tutorial ]
1 !----------------------------------------------------------------------
2 ! Copyright (C) 2001--2008 SPMODEL Development Group. All rights reserved.
3 !----------------------------------------------------------------------
4 ! Sample program for gtool_history/gtool5 and ISPACK 2002/08/21 S.Takehiro
5 ! 2004/01/26 M.Odaka
6 ! 2009/02/27 Y.Morikawa
7 !
8 ! Solving a linear 2-D shallow water system on a sphere
9 ! with an isolated mountain
10 ! du/dt + u0/a\cos\phi du/d\lambda + v/a du0/d\phi - u0v\tan\phi/a
11 ! - 2\Omega\sin\phi v = -g/a\cos\phi dh/d\lambda - \nu\lapla^4 u,
12 ! dv/dt + u0/a\cos\phi dv/d\lambda + 2 u0 u \tan\phi/a
13 ! + 2\Omega\sin\phi u = -g/a dh/d\phi - \nu\lapla^4 v,
14 ! dh/dt + ( 1/\cos\phi d( (H+h0)u+(h-ht)u0 )/d\lambda
15 ! + 1/\cos\phi d( (H+h0)v\cos\phi)/d\phi ) = - \nu\lapla^4 h.
16 !
17 ! A setup is similar to the experiment of Grose and Hoskins (1979)
18 ! with a superrotating rigid-rotation zonal wind profile.
19 !
20 program sp_topo_gtauto
21
22 use w_module
23 use gtool5
24 use gtool_historyauto ! モジュール指定
25 implicit none
26
27 !---- 空間解像度設定 ----
28 integer, parameter :: im=64, jm=32 ! 格子点の設定(X,Y)
29 integer, parameter :: nm=21
30
31 !---- 変数 ----
32 real(8) :: xy_U(0:im-1,jm) ! 格子点データ(速度経度成分)
33 real(8) :: xy_V(0:im-1,jm) ! 格子点データ(速度緯度成分)
34 real(8) :: xy_H(0:im-1,jm) ! 格子点データ(変位)
:
:
44 real(8) :: xy_Zeta(0:im-1,jm) ! 格子点データ(渦度)
45 real(8) :: xy_Htopo(0:im-1,jm) ! 格子点データ(地形)
46
47 !---- 時間積分パラメター ----
48 real(8):: dt ! 時間間隔 [s]
49 type(DC_DIFFTIME):: deltime, dispint ! 時間間隔, 出力時間間隔
50 type(DC_DIFFTIME):: endtime ! 計算終了時間
51 type(DC_DIFFTIME):: ct ! 現在時刻
:
:
68 !------ NAMELIST ファイル名 ------
69 character(*), parameter:: nmlfile = 'sp_topo_gtauto.nml'
:
:
75 !---- 時間積分パラメターの設定 ----
76 dt = 100. ! 時間間隔 [s]
77 call DCDiffTimeCreate(deltime, dt, 'sec') ! 時間間隔
78 dispint = deltime * 500 ! 出力時間間隔
79 endtime = deltime * 10000 ! 計算終了時間
80 call DCDiffTimeCreate(ct, 0., 'sec') ! 現在時刻
81
82 !---------------- 座標値の設定 ---------------------
83 call w_Initial(nm,im,jm) ! ISPACK初期化
84
85 !------------------- 初期値設定 ----------------------
86 xy_U0 = U0*cos(xy_Lat)
87 xy_H0 = ( Omega*R0*U0/(2*Grav) + U0**2/(4*Grav) )*cos(2*xy_Lat)
88
89 xy_U = 0 ; xy_V = 0 ; xy_H = 0
90
91 w_U0 = w_xy(xy_U0) !; w_H0 = w_xy(xy_H0)
92 w_U = w_xy(xy_U) ; w_V = w_xy(xy_V) ; w_H = w_xy(xy_H)
:
:
107 !------------------- ヒストリー初期設定 ----------------------
108 call output_gtool5_init ! ヒストリー初期化
109 call output_gtool5 ! 初期値出力
110
111 !------------------- 時間積分 ----------------------
112 do while ( ct <= endtime )
113 if ( mod(ct, dispint) == 0 ) then
114 write(6,*) 'it = ', int( ct / deltime )
115 end if
116
117 ct = ct + deltime ! 時刻の進行
118
119 w_U = ( w_U &
120 + dt * w_xy( - xy_U0 * xy_GradLon_w(w_U) / R0 &
121 - xy_V * xy_GradLat_w(w_U0) / R0 &
122 + xy_U0 * xy_V * tan(xy_Lat) / R0 &
123 + 2 * Omega * sin(xy_Lat) * xy_V &
124 - Grav * xy_GradLon_w(w_H)/ R0 ) &
125 )/(1+Nu*(-rn(:,1)/R0**2)**(ndiff/2)*dt)
:
:
144 xy_H = xy_w(w_H)
145
146 call output_gtool5 ! 出力
147 enddo
148
149 call output_gtool5_close ! ファイルのクローズ
150 stop
151
152 contains
153 subroutine output_gtool5_init
154 write(6,'(a)',advance='NO') ' Input NAMELIST file: '
155 read (5,'(a)') nmlfile
156
157 call HistoryAutoCreate( & ! ヒストリー作成
158 title='Shallow water equation on a sphere', &
159 source='Sample program of gtool_historyauto/gtool5', &
160 institution='GFD_Dennou Club davis/spmodel project', &
161 dims=(/'lon','lat','t '/), dimsizes=(/im,jm,0/), &
162 longnames=(/'longitude','latitude ','time '/), &
163 units=(/'degree_east ','degree_north','sec. '/), &
164 origin=ct, interval=dispint, terminus=endtime, &
165 namelist_filename=nmlfile )
166
167 call HistoryAutoPutAxis('lon',x_Lon*180/pi) ! 変数出力
168 call HistoryAutoAddAttr('lon','topology','circular') ! 周期属性
169 call HistoryAutoAddAttr('lon','modulo',360.0) ! 周期属性
170 call HistoryAutoPutAxis('lat',y_Lat*180/pi) ! 変数出力
171
172 call HistoryAutoAddVariable( & ! 変数定義
173 varname='h', dims=(/'lon','lat','t '/), &
174 longname='surface displacement ', units='m')
175
176 call HistoryAutoAddVariable( & ! 変数定義
177 varname='u', dims=(/'lon','lat','t '/), &
178 longname='velocity(longitude) ', units='m/s')
179
180 call HistoryAutoAddVariable( & ! 変数定義
181 varname='v', dims=(/'lon','lat','t '/), &
182 longname='velocity(latitude) ', units='m/s')
183
184 call HistoryAutoAddVariable( & ! 変数定義
185 varname='zeta', dims=(/'lon','lat','t '/), &
186 longname='vorticity', units='1/s')
187
188 end subroutine output_gtool5_init
189
190 subroutine output_gtool5
191 call HistoryAutoPut(ct, 'u', xy_U)
192 call HistoryAutoPut(ct, 'v', xy_V)
193 call HistoryAutoPut(ct, 'h', xy_H)
194 xy_Zeta = xy_w(w_Divlon_xy(xy_V) - w_Divlat_xy(xy_U))/r0
195 call HistoryAutoPut(ct, 'zeta', xy_Zeta)
196 end subroutine output_gtool5
197
198 subroutine output_gtool5_close
199 call HistoryAutoClose
200 end subroutine output_gtool5_close
201
202 end program sp_topo_gtauto
netcdf u {
dimensions:
lon = 32 ;
lat = 16 ;
t = UNLIMITED ; // (21 currently)
variables:
float lon(lon) ;
lon:long_name = "longitude" ;
lon:units = "degree_east" ;
lon:topology = "circular" ;
lon:modulo = 360.f ;
float lat(lat) ;
lat:long_name = "latitude" ;
lat:units = "degree_north" ;
float t(t) ;
t:long_name = "time" ;
t:units = "sec." ;
float u(t, lat, lon) ;
u:long_name = "velocity(longitude)" ;
u:units = "m/s" ;
// global attributes:
:Conventions = "http://www.gfd-dennou.org/library/gtool4/conventions/" ;
:
data:
t = 0, 50000, 100000, 150000, 200000, 250000, 300000, 350000, 400000,
450000, 500000, 550000, 600000, 650000, 700000, 750000, 800000, 850000,
900000, 950000, 1000000 ;
}
47 !---- 時間積分パラメター ----
48 real(8):: dt ! 時間間隔 [s]
49 type(DC_DIFFTIME):: deltime, dispint ! 時間間隔, 出力時間間隔
:
:
75 !---- 時間積分パラメターの設定 ----
76 dt = 100. ! 時間間隔 [s]
77 call DCDiffTimeCreate(deltime, dt, 'sec') ! 時間間隔
78 dispint = deltime * 500 ! 出力時間間隔
:
:
157 call HistoryAutoCreate( & ! ヒストリー作成
158 title='Shallow water equation on a sphere', &
159 source='Sample program of gtool_historyauto/gtool5', &
160 institution='GFD_Dennou Club davis/spmodel project', &
161 dims=(/'lon','lat','t '/), dimsizes=(/im,jm,0/), &
162 longnames=(/'longitude','latitude ','time '/), &
163 units=(/'degree_east ','degree_north','sec. '/), &
164 origin=ct, interval=dispint, terminus=endtime, &
165 namelist_filename=nmlfile )
166
netcdf u {
dimensions:
lon = 32 ;
lat = 16 ;
t = UNLIMITED ; // (21 currently)
variables:
:
float t(t) ;
t:long_name = "time" ;
t:units = "day" ;
float u(t, lat, lon) ;
u:long_name = "velocity(longitude)" ;
u:units = "m/s" ;
// global attributes:
:Conventions = "http://www.gfd-dennou.org/library/gtool4/conventions/" ;
:
data:
t = 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11 ;
}
>ool_historyauto_nml
IntValue = 1., ! 出力間隔の数値
IntUnit = 'day' ! 出力間隔の単位
/
>ool_historyauto_nml
Name = 'u, v, h, zeta' ! 出力変数
/
:
float t(t) ;
t:long_name = "time" ;
t:units = "hour" ;
float u(t, lat, lon) ;
u:long_name = "velocity(longitude)" ;
u:units = "m/s" ;
:
data:
t = 0, 12, 24, 36, 48, 60, 72, 84, 96, 108, 120, 132, 144, 156, 168, 180,
192, 204, 216, 228, 240, 252, 264, 276 ;
}
:
float t(t) ;
t:long_name = "time" ;
t:units = "day" ;
float zeta(t, lat, lon) ;
zeta:long_name = "vorticity" ;
zeta:units = "1/s" ;
:
data:
t = 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11 ;
}
!
! データ出力の全体設定
!
>ool_historyauto_nml
IntValue = 1., ! 出力間隔の数値
IntUnit = 'day', ! 出力間隔の単位
/
!
! データ出力の個別設定
!
>ool_historyauto_nml
Name = 'u, v' ! 出力変数
IntValue = 12., ! 出力間隔の数値
IntUnit = 'hour', ! 出力間隔の単位
/
>ool_historyauto_nml
Name = 'zeta' ! 出力変数
/
$Id: gtauto_first.rd,v 1.3 2009-03-01 02:39:15 morikawa Exp $
Gtool4 Development Group / GFD Dennou Staff