IGModel-SW 1.0

src/IGModelSW_main.f90

Go to the documentation of this file.
00001 !################################################################################
00002 !
00072 
00081 program IGModelSW_main  
00082 
00083   ! モジュール引用 ; Use statements
00084   !
00085   
00086   ! * gtool *
00087   ! *
00088 
00089   ! 種類型パラメータ
00090   ! Kind type parameter
00091   !
00092   use dc_types, only: DP, &    ! 倍精度実数型. Double precision. 
00093     &                 TOKEN, & !
00094     &                 STRING
00095 
00096   ! メッセージ出力
00097   ! Output messgae
00098   !
00099   use dc_message, only: MessageNotify
00100 
00101   ! * IGMBaseLib *
00102   ! *
00103 
00104   ! 正二十面体格子の管理クラス
00105   ! Manager class of icoshedral grid
00106   !
00107   use IcGrid2D_FVM_Manager, only: &
00108     & IcGrid2D_FVM, &
00109     & IcGrid2D_FVM_Init
00110 
00111   ! 正二十面体格子座標データのローダー
00112   !
00113   !
00114   use IcGrid_ncGridLoader, only: &
00115     & load_gridData
00116 
00117   ! 正二十面体格子上の物理場に作用する微分演算を提供
00118   ! Providing some differential operators acting on the physical field on icosahedral grid
00119   !
00120   use Derivate_Field_IcGrid2D_Manager, only: &
00121     & Derivate_Field_IcGrid2D
00122 
00123   ! * IGModel-SW *
00124   ! *
00125 
00126   ! 物理場データの管理
00127   ! Manager of the fields 
00128   !
00129   use field_manager, only: &
00130     & xy_Vel_TL_list, xy_Height_TL_list, &
00131     & DVelDt_TL_list, DHeightDt_TL_list, &
00132     & xy_Htopo, &
00133     & TL_N, TL_Nplus1, &
00134     & TL_DDT_N, TL_DDT_Nminus1, TL_DDT_Nminus2, &
00135     & init_field_manager, update_timelevel
00136   
00137   ! 初期場の生成
00138   ! Generation of the initial field. 
00139   !
00140   !use initial_field_builder, only: &
00141   !  & init_initial_field_builder, &
00142   !  & set_initial_fields
00143 
00144   ! シミュレーションパラメータの管理
00145   ! Manager of simulation parameter
00146   !
00147   use param_manager, only: &
00148     & read_namelistFile, &
00149     & gridData_ncFile, &
00150     & div_level, ic_radius, &
00151     & DelTime, integration_time, &
00152     & output_data_ncFile, output_tick, &
00153     & data_output_flag
00154 
00155   ! 浅水方程式ソルバ
00156   ! Solve shallow water equations
00157   !
00158   use sw_equation_solver, only: &
00159     & init_sw_equation_solver,  &
00160     & temporal_integration
00161 
00162   ! 物理場データの NetCDF ファイルへの出力 
00163   ! Output some physical field data into a NetCDF file
00164   !
00165   use output_field_data, only: &
00166     & init_output_field_data, &
00167     & finalize_output_field_data, &
00168     & output_field
00169 
00170   ! テストケース
00171   ! Test case
00172   !
00173   use class_TestCase6, only: &
00174     & TestCase => TestCase6, &
00175     & initialize_TestCase, finalize_TestCase, &
00176     & set_initial_v, set_initial_h, set_initial_hs, &
00177     & timelevel_updated
00178 
00179 
00180   ! 宣言文 ; Declaration statements
00181   !
00182   implicit none
00183 
00186   type(IcGrid2D_FVM) :: icgrid
00187 
00190   character(STRING) :: nml_filename = 'config.nml' 
00191   
00194   integer :: tstep = 0
00195  
00198   integer :: end_tstep
00199 
00202   integer :: output_tstep
00203 
00206   integer :: ten_percent_tstep
00207 
00210   type(TestCase) :: test_case
00211 
00212 
00213   !********************************************************
00214   !
00215   ! 実行文 ; Executable statements
00216   !
00217   !********************************************************
00218 
00219   !#
00220   !# 下準備
00221   !# Preparation
00222 
00223   ! コマンドライン引数の解析を行う. 
00224   ! Analize the comand line arguments. 
00225   call analize_cmdline_arg()
00226 
00227   ! ネームリストファイルを読み込んで, 各シミュレーションパラメータを設定する. 
00228   ! Read the namelist file and configure each parameter. 
00229   call read_namelistFile(nml_filename)  ! (in)
00230 
00231   ! シミュレーション終了時間, データ書き出しの時間間隔を, 
00232   ! 時間のステップ数に換算する. 
00233   end_tstep = integration_time / DelTime 
00234   ten_percent_tstep = 0.1d0 * end_tstep
00235   output_tstep = output_tick / DelTime
00236 
00237   ! 格子生成
00238   ! Grid generation. 
00239 
00240   ! 正二十面体格子管理クラスのインスタンスを初期化する. 
00241   ! また, 既存のデータファイルから格子座標データを読み込み, IcGrid2D_FVM クラスのオブジェクトに格子点情報を設定する.
00242   call IcGrid2D_FVM_Init(  &
00243     & icgrid, &               ! (inout)
00244     & div_level, ic_radius &  ! (in)
00245     & )
00246 
00247   call load_gridData( &
00248     & gridData_ncFile, &  ! (in)
00249     & icgrid &            ! (inout)
00250     & )
00251 
00252   ! 各種モジュールの初期化
00253   ! Initialize each of the modules. 
00254   !
00255 
00256   ! field_manager モジュールの初期化
00257   call init_field_manager(icgrid)       ! (in)
00258 
00259   ! sw_equation_solver モジュールの初期化
00260   call init_sw_equation_solver(icgrid)  ! (in)
00261 
00262   ! output_field_data モジュールの初期化
00263   call init_output_field_data( &
00264     & output_data_ncFile, icgrid, &        ! (in) 
00265     & integration_time, DelTime, output_tick )  ! (in)
00266 
00267   ! テストケースクラスのオブジェクトの初期化  
00268   call initialize_TestCase(test_case, icgrid)
00269 
00270 ! initial_field_builder モジュールの初期化
00271 !  call init_initial_field_builder(icgrid, test_case)
00272 
00273   ! 速度および高度の初期場を設定する. 
00274   call set_initial_v(test_case, xy_Vel_TL_list(TL_N))
00275   call set_initial_h(test_case, xy_Height_TL_list(TL_N))
00276   call set_initial_hs(test_case, xy_Htopo)
00277 
00278 !  call set_initial_fields( &
00279 !    & v_tl_list(TL_N), h_tl_list(TL_N), h_s )  ! (inout)
00280 
00281   ! 初期場を netCDF ファイルに書き出す. 
00282   call output_field( &
00283     & 0, xy_Vel_TL_list(TL_N), xy_Height_TL_list(TL_N), xy_Htopo )  ! (in) 
00284   
00285 
00286   !# 時間積分ループ
00287   !# The loop for temporal integration. 
00288   !#
00289 
00290   ! Test_Case モジュールに対して, 時間積分タイムレベルの更新を行う.
00291   call timelevel_Updated( &
00292     & test_case,                       &  ! (inout)
00293     & tstep, DelTime,                       &  ! (in)
00294     & xy_Vel_TL_list(TL_N), xy_Height_TL_list(TL_N) &  ! (inout)
00295      & )
00296 
00297   ! 時間積分ループの開始 
00298   ! Starts the loop for temporal integration.
00299 
00300   do tstep=1, end_tstep
00301     
00302     ! 時間積分が全体の 10 % 分完了したごとに, 進捗状況を標準出力する. 
00303     if( .not. ten_percent_tstep == 0 ) then
00304       if ( mod(tstep, ten_percent_tstep ) == 0 ) then
00305        
00306         call MessageNotify( "M", "IGModel-SW", "-------- %d/%d", i=(/ tstep, end_tstep /) )    
00307       end if
00308     end if
00309 
00310     ! 浅水方程式系の時間積分 
00311     ! Performs time integration of shallow water equations. 
00312 
00313     call temporal_integration( &
00314       &    tstep, DelTime, icgrid )     ! (in)
00315 
00316     ! タイムレベルの更新
00317     ! Updates the time level. 
00318     
00319     ! field_manager モジュールに対して, タイムレベルの更新を行う. 
00320     ! Updates the time level for field_manager manager. 
00321     call update_timeLevel()
00322 
00323     ! Test_Case モジュールに対して, タイムレベルの更新を行う. 
00324     ! Updates the time level for Test_Case manager. 
00325     call timelevel_Updated( &
00326       &   test_case,                       &  ! (inout)
00327       &   tstep, DelTime,                       &  ! (in)
00328       &   xy_Vel_TL_list(TL_N), xy_Height_TL_list(TL_N) &  ! (inout)
00329       & )
00330 
00331     ! データ書き出し
00332     ! Outputs data
00333 
00334     ! タイムレベル更新後のタイムレベル n の物理場データを netCDF ファイルに書き出す. 
00335     if( data_output_flag .and. mod(tstep, output_tstep) == 0 ) then
00336 
00337       call output_field( &
00338         &    tstep, xy_Vel_TL_list(TL_N), xy_Height_TL_list(TL_N) &  ! (in)
00339         & )
00340 
00341     end if
00342 
00343   end do
00344 
00345   !#
00346   !# 後片付け
00347   !# Finalization
00348         
00349   ! テストケースを提供する構造型の変数の最終化.
00350   ! Finalizes the variable of derived type providing the test case.
00351   call finalize_TestCase(test_case)
00352 
00353   ! output_field_data モジュールの最終化. 
00354   ! Finalizes output_field_data module. 
00355   call finalize_output_field_data()
00356 
00357 contains
00358 
00359 !
00370 subroutine analize_cmdline_arg()
00371   ! モジュール引用 ; Use statement 
00372   !
00373 
00374   ! 文字列型
00375   ! String
00376   !
00377   use dc_string, only: StoA
00378 
00379   ! コマンドライン引数の解析
00380   ! Analization of command-line arguments
00381   !
00382   use dc_args, only: ARGS, &
00383     & DCArgsOpen, DCArgsClose, DCArgsOption, &
00384     & DCArgsPutLine, DCArgsDebug, DCArgsHelp, DCArgsStrict 
00385 
00386   ! 作業変数
00387   ! Work variables
00388   !
00389 
00390   type(ARGS) :: arg
00391   logical :: OPT_ncfname
00392   character(STRING) :: VAL_ncfname
00393 
00394   ! 実行文; Executable statement
00395   !
00396 
00397   call DCArgsOpen(arg)
00398   call DCArgsOption(arg, StoA('-C', '--config'), OPT_ncfname, VAL_ncfname, &
00399     &    help="Specify namelist filename for the simulation parameters.[defalut value is " // trim(nml_filename) // "]." )
00400   
00401   call DCArgsDebug(arg)
00402   call DCArgsHelp(arg)
00403   call DCArgsStrict(arg)
00404 
00405   if( OPT_ncfname ) then
00406     nml_filename = VAL_ncfname
00407   end if
00408 
00409 end subroutine analize_cmdline_arg
00410 
00411 end program IGModelSW_main
 All Classes Namespaces Files Functions Variables