IGModel-SW 1.0

src/continious_equation.f90

Go to the documentation of this file.
00001 
00012 module continious_equation
00013 
00014   ! モジュール引用 ; Use statements
00015   !
00016 
00017   !* gtool ****
00018   !*
00019 
00020   ! 種類型パラメータ
00021   ! Kind type parameter
00022   !
00023   use dc_types, only: DP  ! 倍精度実数型. Double precision. 
00024 
00025   !* IGMBaseLib ****
00026   !
00027 
00028   ! 正二十面体格子上の物理場の管理
00029   ! Manager of a physical field on icosahedral grid
00030   !
00031   use Field_IcGrid2D_Manager, only: &
00032     & Field_IcGrid2D, &
00033     & Field_IcGrid2D_Init, &
00034     & get_icgrid, paste_field2D_margin
00035 
00036   ! 正二十面体格子の管理
00037   ! Manager of icoshedral grid
00038   !
00039   use IcGrid2D_FVM_Manager, only: &
00040     & IcGrid2D_FVM, &
00041     & get_EffSize_Min, get_EffSize_Max, &
00042     & RC_REGIONS_NUM
00043 
00044   ! 正二十面体格子上の物理場に作用する微分演算子の提供
00045   ! Providing some differential operators acting on the physical field on icosahedral grid
00046   !
00047   use Derivate_Field_IcGrid2D_Manager, only: &
00048     & Derivate_Field_IcGrid2D, &
00049     & divergence_op
00050 
00051   !* IGModel-SW ****
00052   !*
00053 
00054   ! シミュレーションパラメータの管理
00055   ! Manager of simulation parameter
00056   !
00057   use param_manager, only: &
00058     & Grav, Omega, alpha
00059   
00060   ! OpenMP
00061   !
00062   !
00063   use omp_lib
00064 
00065   ! 宣言文 ; Declaration statements
00066   !
00067   implicit none
00068   private
00069 
00070   ! 公開手続き
00071   ! Public statements
00072   !
00073   public :: continious_equation_Init, calc_continious_eq_dhdt
00074 
00075   ! 非公開変数
00076   ! Pricate statements
00077   !
00078 
00081   type(Field_IcGrid2D), save :: mass_flux
00082 
00083 contains
00084 
00085 !
00096 subroutine continious_equation_Init( &
00097   & icgrid                           &  ! (inout)
00098   & )
00099   
00100   ! 宣言文 ; Declaration statements
00101   !
00102   type(IcGrid2D_FVM), intent(inout) :: icgrid
00103  
00104   ! 実行文 ; Executable statement
00105   !
00106 
00107   !  質量フラックスを管理する構造型 Field_IcGrid2D の変数を初期化する.
00108   call Field_IcGrid2D_Init(self=mass_flux, icgrid=icgrid, name='mass_flux', rank=3)
00109  
00110   
00111 end subroutine continious_equation_Init
00112 
00113 !
00148 subroutine calc_continious_eq_dhdt( &
00149   & DHeightDtN,                    &  ! (inout)
00150   & xy_VelN, xy_HeightN, xy_Htopo, &  ! (in)
00151   & diff_eval                      &  ! (inout)
00152   & )
00153  
00154   ! 宣言文 ; Declaration statements
00155   !
00156   type(Field_IcGrid2D), intent(inout) :: DHeightDtN
00157   type(Field_IcGrid2D), intent(in) :: xy_VelN
00158   type(Field_IcGrid2D), intent(in) :: xy_HeightN
00159   type(Field_IcGrid2D), intent(in) :: xy_Htopo
00160   type(Derivate_Field_IcGrid2D), intent(inout) :: diff_eval
00161 
00162 
00163   ! 作業変数
00164   ! Work variables
00165   !
00166   integer :: EMax, EMin
00167   integer :: rcID,i,j
00168   real(DP) :: s_pos(3)
00169   real(DP) :: f
00170   type(IcGrid2D_FVM), pointer :: icgrid
00171 
00172   ! 実行文 ; Executable statements
00173   !
00174 
00175   icgrid => get_icgrid(xy_velN)
00176   EMin = get_EffSize_Min(icgrid)
00177   EMax = get_EffSize_Max(icgrid)
00178 
00179   ! 質量フラックスを計算する.
00180   ! ただし, ここではマイナス符号を付加した \f$ -h\Dvect{v}\f$ を計算することにする. 
00181   ! Calculates the mass flux.
00182   ! Note that \f$ -h\Dvect{v}\f$ which is mass flux muitiplied by -1 is actually calculated here.
00183 
00184   do rcID=1, RC_REGIONS_NUM
00185     !$omp parallel do
00186     do j=EMin, EMax
00187       do i=EMin, EMax
00188 
00189         mass_flux%val(rcID,i,j,1:3) = &
00190           & - ( xy_HeightN%val(rcID,i,j,1) - xy_Htopo%val(rcID,i,j, 1) ) * xy_VelN%val(rcID,i,j,1:3)
00191 
00192       end do
00193     end do
00194   end do
00195 
00196   ! ノリシロの貼りあわせ
00197   call paste_field2D_margin(mass_flux)
00198   
00199   ! 質量フラックスの収束を計算し, その結果を高度の時間変化率とする.
00200   ! Calculates the convergence of mass flux,
00201   ! and the result becomes the rate of change of height wit respect to time.
00202 
00203   call divergence_op(diff_eval, mass_flux, DHeightDtN)
00204   
00205 end subroutine calc_continious_eq_dhdt
00206 
00207 end module continious_equation
 All Classes Namespaces Files Functions Variables