IGModel-SW 1.0
|
00001 00009 module class_TestCase1 00010 00011 ! モジュール引用 ; Use statements 00012 ! 00013 00014 ! 種類型パラメータ 00015 ! Kind type parameter 00016 ! 00017 use dc_types, only: DP, & ! 倍精度実数型. Double precision. 00018 & TOKEN 00019 00020 ! 00021 ! 00022 ! 00023 use gtool_history 00024 00025 ! * IGMBaseLib * 00026 ! * 00027 00028 ! 数学 00029 ! Mathematic 00030 ! 00031 use igmcore_math, only: & 00032 & PI, toDegrees 00033 00034 ! 線形代数 00035 ! Linear algebra 00036 ! 00037 use igmcore_linear_algebra, only: & 00038 & rotateZ, rotateY 00039 00040 ! 座標変換 00041 ! Coordinate conversion 00042 ! 00043 use igmcore_coordinate_conversion, only: & 00044 & orth_to_geo_pos, geo_to_orth_pos 00045 00046 ! 正二十面体格子の管理クラス 00047 ! 00048 ! 00049 use IcGrid2D_FVM_Manager, only: & 00050 & IcGrid2D_FVM, & 00051 & get_glevel, get_IcRadius, & 00052 & RC_REGIONS_NUM 00053 00054 ! 正二十面体格子上に分布する物理場の管理クラス 00055 ! 00056 ! 00057 use Field_IcGrid2D_manager, only: & 00058 & Field_IcGrid2D, & 00059 & Field_IcGrid2D_Init, get_icgrid 00060 00061 ! 00062 ! 00063 ! 00064 use Field_Pattern_Builder, only: & 00065 & create_solid_rotation_field, create_cosine_bell_field 00066 00067 ! 正二十面体格子上の物理場における誤差ノルムチェック 00068 ! 00069 ! 00070 use Field_IcGrid2D_error_norm_check, only: & 00071 & numerical_error_norm1, numerical_error_norm2, & 00072 & numerical_error_norminfinity 00073 00074 ! 00075 ! 00076 ! 00077 use IcGrid_ncWriter_mod, only: & 00078 & IcGrid_ncWriter, IcGrid_ncWriter_Init, & 00079 & open_ncFile, ncdef_simulation_parameter, ncdef_GridData, ncdef_FieldData, end_ncdef, & 00080 & write_GridData, write_FieldData, close_ncFile, & 00081 & increase_recorde_counter 00082 00083 00084 ! * IGModel-SW **** 00085 ! * 00086 00087 ! パラメータ管理 00088 ! 00089 ! 00090 use param_manager, only: & 00091 & alpha, output_tick, integration_time, DelTime 00092 00093 ! 宣言文 ; Declaration statement 00094 ! 00095 implicit none 00096 private 00097 00098 ! 公開変数 00099 ! Public variables 00100 ! 00101 00104 real(DP), parameter :: h_0 = 1000.0d0 00105 00108 real(DP), parameter :: init_theta_c = 0.0d0 00109 00112 real(DP), parameter :: init_lambda_c = 3.0d0 * PI / 2.0d0 00113 00116 real(DP) :: cb_R 00117 00120 real(DP) :: u_0 00121 00122 ! クラス定義 00123 ! type definition 00124 ! 00125 00126 ! 00130 type, public :: TestCase1 00131 00132 ! 宣言文 ; Declaration statement 00133 ! 00134 00135 ! 非公開メンバ変数 00136 ! Private member variable 00137 ! 00138 type(Field_IcGrid2D) :: true_h 00139 ! f2003: type(Field_IcGrid2D), private :: true_h 00140 00141 type(Field_IcGrid2D) :: error_h 00142 ! f2003: type(Field_IcGrid2D), private :: error_h 00143 00144 integer :: errorh_ncVarID 00145 ! f2003: integer, private :: errorh_ncVarID 00146 00147 type(Field_IcGrid2D) :: solid_rot_field 00148 ! f2003: type(Field_IcGrid2D), private :: solid_rot_field 00149 00150 end type TestCase1 00151 00152 ! 00153 ! 00154 ! 00155 interface initialize_TestCase 00156 module procedure init_TestCase1 00157 end interface initialize_TestCase 00158 00159 interface finalize_TestCase 00160 module procedure finalize_TestCase1 00161 end interface finalize_TestCase 00162 00163 public :: initialize_TestCase, finalize_TestCase 00164 public :: set_initial_h, set_initial_hs, set_initial_v, timelevel_Updated 00165 00166 ! 非公開変数 00167 ! Private variables 00168 ! 00169 real(DP), private :: t_lambdaC 00170 real(DP), private :: t_thetaC 00171 00172 character(TOKEN) :: filename = 'error_norm.dat' 00173 type(IcGrid_ncWriter), save :: writer 00174 00175 contains 00176 00177 ! 00183 subroutine init_TestCase1(self, icgrid_ref ) 00184 00185 ! 宣言文 ; Declaration statement 00186 ! 00187 type(TestCase1), intent(inout) :: self 00188 type(IcGrid2D_FVM), intent(in) :: icgrid_ref 00189 00190 ! 作業変数 00191 ! Work variables 00192 ! 00193 real(DP) :: ic_radius 00194 real(DP) :: angular_speed 00195 00196 ! 実行文 ; Executable statement 00197 ! 00198 00199 ic_radius = get_IcRadius(icgrid_ref) 00200 00201 cb_R = ic_radius / 3.0d0 00202 u_0 = 2 * PI * ic_radius / ( 12.0d0 * 24.0d0 * 3600d0 ) 00203 angular_speed = 2 * PI / ( 12.0d0 * 24.0d0 * 3600d0 ) 00204 00205 call Field_IcGrid2D_Init(self%true_h, icgrid_ref, 'true_h', 1) 00206 call Field_IcGrid2D_Init(self%error_h, icgrid_ref, 'error_h', 1) 00207 call Field_IcGrid2D_Init(self%solid_rot_field, icgrid_ref, 'solid_rotation_field', 3) 00208 call create_solid_rotation_field(self%solid_rot_field, angular_speed, alpha) 00209 00210 t_lambdaC = init_lambda_c 00211 t_thetaC = init_theta_c 00212 00213 call HistoryCreate( & 00214 & file='error_norm.nc', title='error norm check', & 00215 & source=' ', institution=' ', & 00216 & dims=(/ 't' /), dimsizes=(/ 0 /), & 00217 & longnames=(/ 'time' /), units=(/ 's' /), & 00218 & origin=real(0d0), interval=real(output_tick) ) 00219 00220 call HistoryAddVariable( & 00221 & varname='norm1', dims=(/ 't' /), & 00222 & longname='norm1', units=' ', xtype='double' ) 00223 00224 call HistoryAddVariable( & 00225 & varname='norm2', dims=(/ 't' /), & 00226 & longname='norm2', units=' ', xtype='double' ) 00227 00228 call HistoryAddVariable( & 00229 & varname='norm_infinity', dims=(/ 't' /), & 00230 & longname='norm_infinity', units=' ', xtype='double' ) 00231 00232 call HistoryAddAttr('norm1', '+glevel', get_glevel(icgrid_ref)) 00233 call HistoryAddAttr('norm1', '+alpha', alpha) 00234 00235 ! 00236 ! 00237 call IcGrid_ncWriter_Init(writer, icgrid_ref) 00238 call open_ncFile(writer, 'h_error_check.nc') 00239 call ncdef_Simulation_Parameter(writer, integration_time, DelTime, output_tick) 00240 00241 ! 格子点情報の定義 00242 call ncdef_GridData(writer) 00243 self%errorh_ncVarID = ncdef_FieldData(writer, self%error_h) 00244 00245 ! 定義の終了. 00246 call end_ncdef(writer) 00247 00248 ! 格子情報を netCDF に書きこむ. 00249 ! 00250 call write_GridData(writer) 00251 00252 00253 end subroutine init_TestCase1 00254 00255 ! 00260 subroutine finalize_TestCase1(self) 00261 ! 宣言文 ; Declaration statement 00262 ! 00263 type(TestCase1), intent(inout) :: self 00264 00265 ! 実行文 ; Executable statement 00266 ! 00267 00268 ! 00269 call HistoryClose() 00270 ! 00271 call close_ncFile(writer) 00272 00273 end subroutine finalize_TestCase1 00274 00275 ! 00281 subroutine set_initial_v(self, init_v) 00282 ! 宣言文 ; Declaration statement 00283 ! 00284 type(TestCase1), intent(inout) :: self 00285 type(Field_IcGrid2D), intent(inout) :: init_v 00286 00287 ! 実行文 ; Executable statement 00288 ! 00289 00290 write(*,*) 'init_Testcase1_v' 00291 init_v%val(:,:,:,:) = self%solid_rot_field%val(:,:,:,:) 00292 00293 end subroutine set_initial_v 00294 00295 !! 00296 !! @param[in,out] self TestCase1 クラスのオブジェクトの参照. 00297 !! @param[in,out] init_h 00298 !! 00299 subroutine set_initial_h(self, init_h) 00300 ! 宣言文 ; Declaration statement 00301 ! 00302 type(TestCase1), intent(inout) :: self 00303 type(Field_IcGrid2D), intent(inout) :: init_h 00304 00305 ! 実行文 ; Executable statement 00306 ! 00307 00308 write(*,*) 'init_Testcase1_h' 00309 00310 ! 00311 call create_cosine_bell_field(init_h, init_lambda_c, init_theta_c, cb_R, h_0) 00312 !call eval_numcal_h_solution(0d0, self%true_h, init_h) 00313 00314 00315 end subroutine set_initial_h 00316 00317 ! 00323 subroutine set_initial_hs(self, init_hs) 00324 00325 ! 宣言文 ; Declaration statement 00326 ! 00327 type(TestCase1), intent(inout) :: self 00328 type(Field_IcGrid2D), intent(inout) :: init_hs 00329 00330 ! 実行文 ; Executable statement 00331 ! 00332 00333 write(*,*) 'Testcase1 module - initialize height surface field' 00334 init_hs%val(:,:,:,:) = 0.0d0 00335 00336 end subroutine set_initial_hs 00337 00338 ! 00347 subroutine timelevel_Updated(self, tstep, dt, v_n, h_n) 00348 00349 ! 宣言文 ; Declaration statement 00350 ! 00351 type(TestCase1), intent(inout) :: self 00352 integer, intent(in) :: tstep 00353 real(DP), intent(in) :: dt 00354 type(Field_IcGrid2D), intent(inout) :: v_n 00355 type(Field_IcGrid2D), intent(inout) :: h_n 00356 00357 ! 作業変数 00358 ! Work variables 00359 ! 00360 real(DP) :: new_bellC(3), orth_bellC_onXY(3) 00361 type(IcGrid2D_FVM), pointer :: icgrid 00362 real(DP) :: ic_radius 00363 00364 ! 実行文 ; Executable statement 00365 ! 00366 00367 ! 速度場を剛体回転場で上書きする. 00368 v_n%val(:,:,:,:) = self%solid_rot_field%val(:,:,:,:) 00369 00370 ! 00371 ! 00372 icgrid => get_icgrid(v_n) 00373 ic_radius = get_IcRadius(icgrid) 00374 00375 if ( tstep /= 0 ) then 00376 orth_bellC_onXY(:) = & 00377 & rotateZ( & 00378 & rotateY( geo_to_orth_pos( (/ t_lambdaC, t_thetaC, ic_radius /) ), -alpha ), & 00379 & 2.0d0*PI * dt / ( 12.0d0 * 24.0d0 * 3600d0 ) & 00380 & ) 00381 00382 new_bellC(:) = orth_to_geo_pos(rotateY( orth_bellC_onXY(:), alpha)) 00383 t_lambdaC = new_bellC(1) 00384 t_thetaC = new_bellC(2) 00385 ! write(*,*) '-----', toDegrees(t_lambdaC), toDegrees(t_thetaC) 00386 end if 00387 00388 if ( mod(tstep, int(output_tick / dt)) == 0 ) then 00389 ! 00390 call eval_numcal_h_solution( tstep*dt, self%true_h, h_n) 00391 00392 ! 00393 self%error_h%val(:,:,:,:) = self%true_h%val(:,:,:,:) - h_n%val(:,:,:,:) 00394 call write_FieldData(writer, self%errorh_ncVarID, self%error_h) 00395 call increase_recorde_counter(writer) 00396 end if 00397 00398 end subroutine timelevel_Updated 00399 00400 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00401 ! 00402 ! 非公開手続き 00403 ! 00404 !!!!!!!!!!!!!!!!!!!!!!!!! 00405 00406 ! 00412 subroutine eval_numcal_h_solution(t, true_h, h_n) 00413 00414 ! 宣言文 ; Declaration statement 00415 ! 00416 real(DP), intent(in) :: t 00417 type(Field_IcGrid2D), intent(inout) :: true_h 00418 type(Field_IcGrid2D), intent(in) :: h_n 00419 00420 ! 作業変数 00421 ! Work variables 00422 ! 00423 real(DP) norm1Val, norm2Val, normInfityVal 00424 00425 ! 実行文 ; Executable statement 00426 ! 00427 00428 ! 00429 call create_cosine_bell_field(true_h, t_lambdaC, t_thetaC, cb_R, h_0) 00430 00431 norm1Val = numerical_error_norm1(true_h, h_n) 00432 norm2Val = numerical_error_norm2(true_h, h_n) 00433 normInfityVal = numerical_error_norminfinity(true_h, h_n, maxmin_info_flag=.false.) 00434 ! write(*,*) 'h', true_h%val(10,33,33,:), h_n%val(10,33,33,:) 00435 00436 call HistoryPut('norm1', norm1Val) 00437 call HistoryPut('norm2', norm2Val) 00438 call HistoryPut('norm_infinity', normInfityVal) 00439 00440 end subroutine eval_numcal_h_solution 00441 00442 end module class_TestCase1