IGModel-SW 1.0

src/field_manager.f90

Go to the documentation of this file.
00001 
00011 module field_manager
00012 
00013   ! モジュール引用 ; Use statements
00014   !
00015 
00016   !* gtool **
00017   !* 
00018 
00019   ! 種類型パラメータ
00020   ! Kind type parameter 
00021   !
00022   use dc_types, only: TOKEN  ! 倍精度実数型. Double precision.
00023 
00024   !* IGMBaseLib **
00025   !*
00026 
00027   ! 正二十面体格子を管理する管理するクラス
00028   ! Manager class of icosahedral grid
00029   !
00030   use IcGrid2D_FVM_Manager, only: &
00031     & IcGrid2D_FVM
00032 
00033   ! 正二十面体格子上の物理場データを管理するクラス
00034   ! Manager class of a physical field data on an icosahedral grid
00035   !
00036   use Field_IcGrid2D_Manager, only: &
00037     & Field_IcGrid2D, &
00038     & Field_IcGrid2D_Init, paste_field2D_margin
00039   
00040   ! 正二十面体格子上の物理場に対して微分演算を提供するクラス
00041   ! Class providing some differential operators acting on a physcial field on icosahedral grid
00042   !
00043   use Derivate_Field_IcGrid2D_Manager, only: &
00044     & Derivate_Field_IcGrid2D, Derivate_Field_IcGrid2D_Init
00045 
00046   ! 典型的な物理場パターンの生成
00047   ! Generation of the typical physical fields pattern
00048   !
00049   use Field_Pattern_Builder, only: &
00050     & create_planetaryVorticity_field
00051 
00052   !
00053   !
00054   !
00055   use IcGrid_ncWriter_mod, only: &
00056     & IcGrid_ncWriter, &
00057     & IcGrid_ncWriter_Init
00058 
00059   !* IGModel-SW
00060   !*
00061 
00062   ! シミュレーションパラメータの管理
00063   ! Manager of simulation parameters
00064   !
00065   use param_manager, only: &
00066     & alpha, Omega
00067 
00068   ! 宣言文 ; Declaration statements
00069   !
00070   implicit none
00071   private
00072   
00073   ! 公開変数
00074   ! Public variables
00075   !
00076 
00077   ! 各時間レベルを表す配列インデックス管理する変数の宣言.
00078   ! 時間微分を含む物理場変数は, そうでない変数に対して分別的に時間レベルを管理する.
00079   !
00080 
00083   integer, public :: TL_Nplus1 = 1
00084   
00087   integer, public :: TL_N      = 2
00088 
00091   integer, parameter, public :: TL_NUMS = 2
00092 
00095   integer, public :: TL_DDT_N       = 1
00096   
00099   integer, public :: TL_DDT_Nminus1 = 2
00100   
00103   integer, public :: TL_DDT_Nminus2 = 3
00104 
00107   integer, parameter, public :: TL_DDT_NUMS = 3
00108 
00109   ! 正二十面格子上の物理場データを管理する構造型 Field_IcGrid2D の変数の宣言.
00110   !
00111 
00114   type(Field_IcGrid2D), public, save :: xy_Coli
00115   
00118   type(Field_IcGrid2D), public, save :: xy_Htopo
00119 
00122   type(Field_IcGrid2D), public, save :: xy_Vel_TL_list(TL_NUMS)
00123   
00126   type(Field_IcGrid2D), public, save :: xy_Height_TL_list(TL_NUMS)
00127 
00130   type(Field_IcGrid2D), public, save :: DVelDt_TL_list(TL_DDT_NUMS)
00131  
00134   type(Field_IcGrid2D), public, save :: DHeightDt_TL_list(TL_DDT_NUMS)
00135   
00139   type(Derivate_Field_IcGrid2D), public :: diff_eval
00140 
00141 
00142   ! 公開手続き
00143   ! Public procedure
00144   !
00145   public :: init_field_manager, update_timeLevel
00146 
00147 contains
00148 
00149 !
00160 subroutine init_field_manager( &
00161   & icgrid &
00162   & )
00163 
00164   ! 宣言文 ; Declaration statements
00165   !
00166   type(IcGrid2D_FVM), intent(in) :: icgrid
00167  
00168   ! 作業変数
00169   ! Work variables
00170   !
00171   integer :: tlID
00172   
00173   ! 実行文 ; Executable statements
00174   !
00175 
00176 !  write(*,*) 'init field_manager'
00177 
00178   ! 海底地形の高度場を保持する Field_IcGrid2D クラスのインスタンスの初期化.
00179   call Field_IcGrid2D_Init(self=xy_Htopo, icgrid=icgrid, name='h_s', rank=1)
00180 
00181   ! 惑星渦度場を保持する Field_IcGrid2D クラスのインスタンスの初期化.
00182   call Field_IcGrid2D_Init(self=xy_Coli, icgrid=icgrid, name='planetary_vorticity', rank=1)
00183   call create_planetaryVorticity_field(xy_Coli, Omega, alpha)
00184 
00185   ! タイムレベル n+1,n に対する速度・高度場のデータを保持する Field_IcGrid2D クラスのインスタンスの初期化.
00186   do tlID=1, TL_NUMS
00187 
00188     call Field_IcGrid2D_Init( &
00189       & self=xy_Vel_TL_list(tlID), icgrid=icgrid, &
00190       & name='v', rank=3, long_name='velocity', units='m/s' )
00191 
00192     call Field_IcGrid2D_Init( &
00193       & self=xy_Height_TL_list(tlID), icgrid=icgrid, &
00194       & name='h', rank=1, long_name='surface_height', units='m' )
00195 
00196   end do
00197 
00198   ! タイムレベル n, n-1, n-2 の速度・高度場を局所時間微分した場のデータを保持する, 
00199   ! Field_IcGrid2D クラスのインスタンスの初期化する.
00200   do tlID=1, TL_DDT_NUMS
00201     call Field_IcGrid2D_Init( &
00202       & self=DVelDt_TL_list(tlID), icgrid=icgrid, &
00203       & name='dvdt', rank=3 )
00204 
00205     call Field_IcGrid2D_Init( &
00206       & self=DHeightDt_TL_list(tlID), icgrid=icgrid, &
00207       & name='dhdt', rank=1 )
00208 
00209   end do
00210 
00211   ! 微分演算を提供するクラスのインスタンスを初期化する.
00212   call Derivate_Field_IcGrid2D_Init(diff_eval, icgrid)
00213 
00214 end subroutine init_field_manager
00215 
00216 ! 
00241 subroutine update_timelevel()
00242 
00243   ! 作業変数 
00244   ! Work variables
00245   !
00246   integer :: tlID
00247 
00248   ! 実行文 ; Executable statement
00249   !
00250   
00251   ! タイムレベル n+ 1 の高度・速度場
00252   call paste_field2D_margin(xy_Vel_TL_list(TL_Nplus1))
00253   call paste_field2D_margin(xy_Height_TL_list(TL_Nplus1))
00254 
00255   !
00256   TL_Nplus1 = get_next_TLid(TL_Nplus1, TL_NUMS)
00257   TL_N = get_next_TLid(TL_N, TL_NUMS)
00258 
00259   TL_DDT_N = get_next_TLid( TL_DDT_N, TL_DDT_NUMS)
00260   TL_DDT_Nminus1 = get_next_TLid( TL_DDT_Nminus1, TL_DDT_NUMS)
00261   TL_DDT_Nminus2 = get_next_TLid( TL_DDT_Nminus2, TL_DDT_NUMS)
00262   
00263   !write(*,*) 'dt_n', TL_DDT_N
00264   !write(*,*) 'dt_n-1', TL_DDT_Nminus1
00265   !write(*,*) 'dt_n-2', TL_DDT_Nminus2
00266 end subroutine update_timelevel
00267 
00268 !!!!!!!!!!!!!!!!!!!!!!!!!!!!
00269 ! 
00270 ! 非公開手続き
00271 !
00272 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00273 
00274 !
00282 function get_next_TLid( &
00283   & target_TL, TimeLevelNums  &
00284   & ) result( NextTLid )
00285   
00286   ! 宣言文 ; Declaration statement
00287   !
00288   integer, intent(in) :: target_TL
00289   integer, intent(in) :: TimeLevelNums
00290   integer nextTLid
00291   
00292   ! 実行文 ; Executable statement
00293   !
00294   
00295   nextTLid = target_TL - 1
00296 
00297   if ( nextTLid == 0 ) then
00298     nextTLid = TimeLevelNums
00299   end if
00300 end function get_next_TLid
00301 
00302 end module field_manager
 All Classes Namespaces Files Functions Variables