IGModel-SW 1.0
|
00001 00011 module conservation_analysis 00012 00013 ! モジュール引用 ; Use statement 00014 ! 00015 00016 !* gtool ** 00017 !* 00018 00019 ! 種類型パラメータ 00020 ! Kind type parameter 00021 ! 00022 use dc_types, only: DP ! 倍精度実数型. Double precision. 00023 00024 !* IGMBaseLib ** 00025 !* 00026 00027 ! 線形代数 00028 ! Linear algebra 00029 ! 00030 use igmcore_linear_algebra, only: & 00031 & dot 00032 00033 ! 正二十面体格子の管理 00034 ! Manager of icosahedral grid 00035 ! 00036 use IcGrid2D_FVM_Manager, only: & 00037 & IcGrid2D_FVM, & 00038 & get_EffSize_Min, get_EffSize_Max, & 00039 & RC_REGIONS_NUM 00040 00041 ! 正二十面格子上の物理場データの管理 00042 ! Manager of the physical field data on icosahedral grid 00043 ! 00044 use Field_IcGrid2D_Manager, only: & 00045 & Field_IcGrid2D, & 00046 & Field_IcGrid2D_Init, Field_IcGrid2D_Final, & 00047 & get_IcGrid, global_mean 00048 00049 ! 正二十面体格子上の物理場に作用する微分演算 00050 ! Some differntial operators acting on a physical field on icosahedral grid 00051 ! 00052 use Derivate_Field_IcGrid2D_Manager, only: & 00053 & Derivate_Field_IcGrid2D, & 00054 & vertical_curl_op 00055 00056 !* IGModel-SW ** 00057 !* 00058 00059 ! シミュレーションパラメータの管理 00060 ! Manager of the simulation parameters 00061 ! 00062 use param_manager, only: Grav 00063 00064 ! 宣言文 ; Declaration statement 00065 ! 00066 implicit none 00067 private 00068 00069 ! 公開手続き 00070 ! Public procedures 00071 ! 00072 public :: init_conservation_analysis, calc_total_energy, calc_potential_enstrophy 00073 00074 ! 非公開変数 00075 ! Private variables 00076 00079 real(DP) :: Ep0 00080 00081 contains 00082 00083 ! 00100 subroutine init_conservation_analysis(init_h, hs) 00101 ! 宣言文 ; Declaration statement 00102 ! 00103 type(Field_IcGrid2D), intent(in) :: init_h 00104 type(Field_IcGrid2D), intent(in) :: hs 00105 00106 ! 作業変数 00107 ! Work variables 00108 ! 00109 type(Field_IcGrid2D) :: tmp_Ep_field 00110 type(IcGrid2D_FVM), pointer :: icgrid 00111 integer :: rcID, i, j 00112 integer :: EMax, EMin 00113 00114 ! 実行文 ; Executable statements 00115 ! 00116 icgrid => get_IcGrid(init_h) 00117 EMin = get_EffSize_Min(icgrid) 00118 EMax = get_EffSize_Max(icgrid) 00119 00120 call Field_IcGrid2D_Init(tmp_Ep_field, icgrid, 'tmp_PE_field', 1) 00121 00122 do rcID=1, RC_REGIONS_NUM 00123 do j= EMin, EMax 00124 do i= EMin, EMax 00125 tmp_Ep_field%val(rcID,i,j,1) = & 00126 & 0.5d0 * Grav * & 00127 & ( init_h%val(rcID,i,j,1)**2 - hs%val(rcID,i,j,1)**2 ) 00128 00129 end do 00130 end do 00131 end do 00132 00133 Ep0 = global_mean(tmp_Ep_field) 00134 !write(*,*) 'Ep0=',Ep0 00135 call Field_IcGrid2D_Final(tmp_Ep_field) 00136 00137 end subroutine init_conservation_analysis 00138 00139 ! 00154 function calc_total_energy(h, hs, v) result(total_energy) 00155 ! 宣言文 ; Declaration statement 00156 ! 00157 type(Field_IcGrid2D), intent(in) :: h 00158 type(Field_IcGrid2D), intent(in) :: hs 00159 type(Field_IcGrid2D), intent(in) :: v 00160 real(DP) total_energy 00161 00162 ! 作業変数 00163 ! Work variables 00164 ! 00165 type(Field_IcGrid2D) :: tmp_TE_field 00166 type(IcGrid2D_FVM), pointer :: icgrid 00167 integer :: rcID, i, j 00168 integer :: EMax, EMin 00169 00170 ! 実行文 00171 ! 00172 icgrid => get_IcGrid(h) 00173 EMin = get_EffSize_Min(icgrid) 00174 EMax = get_EffSize_Max(icgrid) 00175 00176 call Field_IcGrid2D_Init(tmp_TE_field, icgrid, 'tmp_TE_field', 1) 00177 00178 do rcID=1, RC_REGIONS_NUM 00179 do j= EMin, EMax 00180 do i= EMin, EMax 00181 tmp_TE_field%val(rcID,i,j,1) = & 00182 & 0.5d0 * ( h%val(rcID,i,j,1) - hs%val(rcID,i,j,1) ) * dot(v%val(rcID,i,j,:), v%val(rcID,i,j,:)) & 00183 & + 0.5d0 * Grav * ( h%val(rcID,i,j,1)**2 - hs%val(rcID,i,j,1)**2 ) & 00184 & - Ep0 00185 end do 00186 end do 00187 end do 00188 00189 total_energy = global_mean(tmp_TE_field) 00190 call Field_IcGrid2D_Final(tmp_TE_field) 00191 00192 end function calc_total_energy 00193 00194 ! 00213 function calc_potential_enstrophy(h, hs, v, f, diff_eval) result(potential_enstrophy) 00214 ! 宣言文 ; Declaration statement 00215 ! 00216 type(Field_IcGrid2D), intent(in) :: h 00217 type(Field_IcGrid2D), intent(in) :: hs 00218 type(Field_IcGrid2D), intent(in) :: v 00219 type(Field_IcGrid2D), intent(in) :: f 00220 type(Derivate_Field_IcGrid2D), intent(inout) :: diff_eval 00221 real(DP) potential_enstrophy 00222 00223 ! 作業変数 00224 ! Work variables 00225 ! 00226 type(Field_IcGrid2D) :: tmp_PE_field 00227 type(Field_IcGrid2D) :: vorticity 00228 type(IcGrid2D_FVM), pointer :: icgrid 00229 integer :: rcID, i, j 00230 integer :: EMax, EMin 00231 00232 ! 実行文 00233 ! 00234 icgrid => get_IcGrid(h) 00235 EMin = get_EffSize_Min(icgrid) 00236 EMax = get_EffSize_Max(icgrid) 00237 00238 call Field_IcGrid2D_Init(tmp_PE_field, icgrid, 'tmp_PE_field', 1) 00239 call Field_IcGrid2D_Init(vorticity, icgrid, 'vorticity', 1) 00240 00241 ! 速度場に回転演算子を作用させて, 渦度場を計算する. 00242 ! 00243 call vertical_curl_op(diff_eval, v, vorticity ) 00244 00245 do rcID=1, RC_REGIONS_NUM 00246 do j= EMin, EMax 00247 do i= EMin, EMax 00248 tmp_PE_field%val(rcID,i,j,1) = & 00249 & 0.5d0 * ( vorticity%val(rcID,i,j,1) + f%val(rcID,i,j,1) )**2 & 00250 & / ( h%val(rcID,i,j,1) - hs%val(rcID,i,j,1) ) 00251 end do 00252 end do 00253 end do 00254 00255 potential_enstrophy = global_mean(tmp_PE_field) 00256 call Field_IcGrid2D_Final(tmp_PE_field) 00257 00258 end function calc_potential_enstrophy 00259 00260 end module conservation_analysis