IGMBaseLib 1.0

src/core/grid/STDGCGrid_Builder.f90

Go to the documentation of this file.
00001 
00024 module STDGCGrid_Builder
00025 
00026   ! モジュール引用 ; Use statements
00027   !
00028 
00029   ! 種類型パラメタ
00030   ! Kind type parameter
00031   !
00032   use dc_types, only: DP         ! 倍精度実数型. Double precision.
00033 
00034   ! 数学定数
00035   ! Mathematic
00036   !
00037   use igmcore_math, only: PI   ! 円周率
00038 
00039   ! 線形代数
00040   ! Linear algebra
00041   !
00042   use igmcore_linear_algebra, only: &
00043     & vec_length, vec_normarize, &
00044     & cross
00045 
00046   ! 球面三角法
00047   ! Spherical trigonometry
00048   !
00049   use igmcore_spherical_trigonometry, only: geodesic_arc_length
00050 
00051   ! 正二十面体格子データの管理
00052   ! Manager of Icosahedral grid data
00053   !
00054   use IcGrid2D_FVM_Manager, only: &
00055     & IcGrid2D_FVM, &
00056     & paste_margin_width, &
00057     & get_EffSize_Min, get_EffSize_Max, get_IcRadius, &
00058     & RC_REGIONS_NUM
00059 
00060 
00061   ! STD-grid の構築
00062   ! STD-grid construction
00063   !
00064   use STDGrid_Builder, only: &
00065     & calc_STDGrid_GP_vertex, calc_CV_vertex, calc_CV_Area, &
00066     & project_on_sphere
00067 
00068   ! 宣言文 ; Declaration statements
00069   !
00070   implicit none
00071   private
00072 
00073   ! 公開手続き
00074   ! Public procedures
00075   !
00076   public :: construct_STDGCGrid, move_GP_into_CVGravPt
00077   public :: calc_CV_vertex, calc_CV_Area, calc_STDGrid_GP_vertex, project_on_sphere ! Cascade
00078 
00079 contains
00080 
00091 subroutine construct_STDGCGrid( &
00092   & icgrid &  ! (inout)
00093   & )
00094 
00095   ! 宣言文 ; Declaration statements
00096   !
00097   type(IcGrid2D_FVM), intent(inout) :: icgrid
00098 
00099   ! 実行文 ; Executable statements
00100   !
00101 
00102   ! 最初に STD-GC-grid の格子点を生成する.
00103   call calc_STDGrid_GP_vertex(icgrid)
00104 
00105   ! コントロールボリュームの頂点の座標を計算する.
00106   call calc_CV_vertex(icgrid)
00107 
00108   ! 各格子点に付随するコントロールボリュームの重心に, 格子点を移動させる.
00109   call move_GP_into_CVGravPt(icgrid)
00110 
00111   ! 各格子点に付随するコントロールボリュームの面積を計算する.
00112   call calc_CV_Area(icgrid)
00113 
00114 end subroutine construct_STDGCGrid
00115 
00116 !!!!!!!!!!!!!!!!!!!!!!!!!!
00117 !
00118 ! モジュール外非公開手続き
00119 !
00120 !!!!!!!!!!!!!!!!!!!!!!!!!!
00121 
00122 !
00133 subroutine move_GP_into_CVGravPt( &
00134   & icgrid &  ! (inout)
00135   & )
00136 
00137   ! 宣言文 ; Declaration statemnts
00138   !
00139   type(IcGrid2D_FVM), intent(inout) :: icgrid
00140 
00141   ! 作業変数
00142   ! Work variables
00143   !
00144   integer :: rcID, i, j
00145   integer :: EMin, EMax
00146   real(DP) :: ic_radius
00147   real(DP) :: mod_GP(3)
00148   integer :: CV_num, cvID
00149 
00150   ! 実行文 ; Executable statements
00151   !
00152 
00153   ! 格子点座標を保持する配列の境界インデックスに関する情報を取得する.
00154   EMin = get_EffSize_Min(icgrid)
00155   EMax = get_EffSize_Max(icgrid)
00156   ic_radius = get_IcRadius(icgrid)
00157 
00158   do rcID=1, RC_REGIONS_NUM
00159     do j=EMin, EMax
00160       do i=EMin, EMax
00161 
00162         ! インデックスに応じて, 格子点に付随するコントールボリュームの頂点数を設定.
00163         if ( ( i==EMin .and. j==EMin ) .or. ( i==EMax .and. j== EMin ) .or. &
00164           &  ( i==EMin .and. j== EMax ) .or. ( i==EMax .and. j== EMax ) ) then
00165           CV_num = 5
00166         else
00167           CV_num = 6
00168         end if
00169 
00170         ! 格子点の位置を, 近接の 3 個のコントールボリューム頂点の重心に移動させる.
00171         mod_GP(:) = 0.0d0
00172         do cvID=1, CV_num
00173           mod_GP(:) = mod_GP(:) + icgrid%rcs_CV(rcID,i,j,cvID, :)
00174         end do
00175 
00176         icgrid%rcs_AGrid(rcID,i,j,:) = mod_GP(:) / CV_num
00177 
00178         ! 現在, コントールボリュームの 3 頂点で作られる平面三角形上に存在する格子点を,
00179         ! 球面上に投影する.
00180         call project_on_sphere(ic_radius, icgrid%rcs_AGrid(rcID,i,j,:))
00181 
00182       end do
00183     end do
00184   end do
00185 
00186   ! のりしろを貼り合わせる
00187   call paste_margin_width(icgrid)
00188 
00189 end subroutine move_GP_into_CVGravPt
00190 
00191 end module STDGCGrid_Builder
 All Classes Namespaces Files Functions Variables