Class | Arare4fileio |
In: |
../src/io/arare4fileio.f90
|
リスタート用の場の情報を netCDF ファイルに出力するためのルーチン
Subroutine : | |||
xyz_PressBZ(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out)
| ||
xyz_ExnerBZ(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out)
| ||
xyz_TempBZ(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out)
| ||
xyz_PTempBZ(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out)
| ||
xyz_DensBZ(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out)
| ||
xyz_VelSoundBZ(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out)
| ||
xyzf_QMixBZ(imin:imax,jmin:jmax,kmin:kmax, ncmax) : | real(DP), intent(out)
| ||
xyz_EffMolWtBZ(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out)
|
subroutine Arare4fileio_BZ_Get( xyz_PressBZ, xyz_ExnerBZ, xyz_TempBZ, xyz_PTempBZ, xyz_DensBZ, xyz_VelSoundBZ, xyzf_QMixBZ, xyz_EffMolWtBZ ) ! !リスタートファイルから情報取得 ! !暗黙の型宣言禁止 implicit none !変数定義 real(DP) :: Var3D(1:nx,1,1:nz) real(DP), intent(out) :: xyz_ExnerBZ (imin:imax,jmin:jmax,kmin:kmax) !基本場のエクスナー関数 real(DP), intent(out) :: xyz_DensBZ (imin:imax,jmin:jmax,kmin:kmax) !基本場の密度 real(DP), intent(out) :: xyz_PTempBZ (imin:imax,jmin:jmax,kmin:kmax) !基本場の温位 real(DP), intent(out) :: xyz_VelSoundBZ (imin:imax,jmin:jmax,kmin:kmax) !基本場の音速 real(DP), intent(out) :: xyz_PressBZ (imin:imax,jmin:jmax,kmin:kmax) !基本場の圧力 real(DP), intent(out) :: xyz_TempBZ (imin:imax,jmin:jmax,kmin:kmax) !基本場の温度 real(DP), intent(out) :: xyzf_QMixBZ (imin:imax,jmin:jmax,kmin:kmax, ncmax) !基本場の混合比 real(DP), intent(out) :: xyz_EffMolWtBZ (imin:imax,jmin:jmax,kmin:kmax) !基本場の分子量効果 character(STRING) :: name !変数名 character(STRING) :: InputFile integer :: f !------------------------------------------------------------- ! 基本場の取得 ! InputFile = trim(Arare4Prefix) // "BasicZ.nc" name = "DensBasicZ" call HistoryGet( InputFile, name, Var3D(:,1,:), flag_mpi_split = FLAG_LIB_MPI ) xyz_DensBZ(1:nx,1,1:nz) = Var3D(1:nx,1,1:nz) call SetMargin_xyz( xyz_DensBZ ) name = "ExnerBasicZ" call HistoryGet( InputFile, name, Var3D(:,1,:), flag_mpi_split = FLAG_LIB_MPI ) xyz_ExnerBZ(1:nx,1,1:nz) = Var3D(1:nx,1,1:nz) call SetMargin_xyz( xyz_ExnerBZ ) name = "PotTempBasicZ" call HistoryGet( InputFile, name, Var3D(:,1,:), flag_mpi_split = FLAG_LIB_MPI ) xyz_PTempBZ(1:nx,1,1:nz) = Var3D(1:nx,1,1:nz) call SetMargin_xyz( xyz_PTempBZ ) name = "VelSoundBasicZ" call HistoryGet( InputFile, name, Var3D(:,1,:), flag_mpi_split = FLAG_LIB_MPI ) xyz_VelSoundBZ(1:nx,1,1:nz) = Var3D(1:nx,1,1:nz) call SetMargin_xyz( xyz_VelSoundBZ ) name = "TempBasicZ" call HistoryGet( InputFile, name, Var3D(:,1,:), flag_mpi_split = FLAG_LIB_MPI ) xyz_TempBZ(1:nx,1,1:nz) = Var3D(1:nx,1,1:nz) call SetMargin_xyz( xyz_TempBZ ) name = "PressBasicZ" call HistoryGet( InputFile, name, Var3D(:,1,:), flag_mpi_split = FLAG_LIB_MPI ) xyz_PressBZ(1:nx,1,1:nz) = Var3D(1:nx,1,1:nz) call SetMargin_xyz( xyz_PressBZ) ! name = "MixRtBasicZ" ! call HistoryGet( InputFile, name, Var4D(:,1,:,:), flag_mpi_split = FLAG_LIB_MPI ) ! xyzf_QMixBZ(1:nx,1,1:nz,1:ncmax) = Var4D(1:nx,1,1:nz,1:ncmax) do f = 1, ncmax name = trim(SpcWetSymbol(f))//"BasicZ" call HistoryGet( InputFile, name, Var3D, flag_mpi_split = FLAG_LIB_MPI ) xyzf_QMixBZ(1:nx,1,1:nz,f) = Var3D(1:nx,1,1:nz) end do call SetMargin_xyzf( xyzf_QMixBZ ) name = "EffMolWtBasicZ" call HistoryGet( InputFile, name, Var3D(:,1,:), flag_mpi_split = FLAG_LIB_MPI ) xyz_EffMolWtBZ(1:nx,1,1:nz) = Var3D(1:nx,1,1:nz) call SetMargin_xyz( xyz_EffMolWtBZ ) end subroutine Arare4fileio_BZ_Get
Subroutine : |
リスタートファイルの書き出し
暗黙の型宣言禁止
This procedure input/output NAMELIST#arare4fileio_nml .
subroutine Arare4fileio_Init ! !リスタートファイルの書き出し ! !暗黙の型宣言禁止 implicit none !変数 integer :: unit !装置番号 !NAMELIST から情報を取得 NAMELIST /arare4fileio_nml/ Arare4Prefix call FileOpen(unit, file=namelist_filename, mode='r') read(unit, NML=arare4fileio_nml) close(unit) !確認 if (myrank == 0) then call MessageNotify( "M", "arare4FileIO_init", "Arare4Prefix = %c", c1=trim(Arare4Prefix) ) end if end subroutine Arare4fileio_Init
Subroutine : | |||
xyz_PressBZ(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out)
| ||
xyz_ExnerBZ(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out)
| ||
xyz_TempBZ(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out)
| ||
xyz_PTempBZ(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out)
| ||
xyz_DensBZ(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out)
| ||
xyz_VelSoundBZ(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out)
| ||
xyzf_QMixBZ(imin:imax,jmin:jmax,kmin:kmax, ncmax) : | real(DP), intent(out)
| ||
xyz_EffMolWtBZ(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out)
|
subroutine Arare4fileio_MMC_BZ_Get( xyz_PressBZ, xyz_ExnerBZ, xyz_TempBZ, xyz_PTempBZ, xyz_DensBZ, xyz_VelSoundBZ, xyzf_QMixBZ, xyz_EffMolWtBZ ) ! !リスタートファイルから情報取得 ! !暗黙の型宣言禁止 implicit none !変数定義 real(DP) :: Var3D(1:nx,1,1:nz) ! real(DP) :: Var4D(1:nx,1,1:nz,1:ncmax) real(DP), intent(out) :: xyz_ExnerBZ (imin:imax,jmin:jmax,kmin:kmax) !基本場のエクスナー関数 real(DP), intent(out) :: xyz_DensBZ (imin:imax,jmin:jmax,kmin:kmax) !基本場の密度 real(DP), intent(out) :: xyz_PTempBZ (imin:imax,jmin:jmax,kmin:kmax) !基本場の温位 real(DP), intent(out) :: xyz_VelSoundBZ (imin:imax,jmin:jmax,kmin:kmax) !基本場の音速 real(DP), intent(out) :: xyz_PressBZ (imin:imax,jmin:jmax,kmin:kmax) !基本場の圧力 real(DP), intent(out) :: xyz_TempBZ (imin:imax,jmin:jmax,kmin:kmax) !基本場の温度 real(DP), intent(out) :: xyzf_QMixBZ (imin:imax,jmin:jmax,kmin:kmax, ncmax) !基本場の混合比 real(DP), intent(out) :: xyz_EffMolWtBZ (imin:imax,jmin:jmax,kmin:kmax) !基本場の分子量効果 character(STRING) :: name !変数名 character(STRING) :: InputFile !------------------------------------------------------------- ! 基本場の取得 ! InputFile = trim(Arare4Prefix) // "BasicZ.nc" name = "DensBasicZ" call HistoryGet( InputFile, name, Var3D(:,1,:), flag_mpi_split = FLAG_LIB_MPI ) xyz_DensBZ(1:nx,1,1:nz) = Var3D(1:nx,1,1:nz) call SetMargin_xyz( xyz_DensBZ ) name = "ExnerBasicZ" call HistoryGet( InputFile, name, Var3D(:,1,:), flag_mpi_split = FLAG_LIB_MPI ) xyz_ExnerBZ(1:nx,1,1:nz) = Var3D(1:nx,1,1:nz) call SetMargin_xyz( xyz_ExnerBZ ) name = "PotTempBasicZ" call HistoryGet( InputFile, name, Var3D(:,1,:), flag_mpi_split = FLAG_LIB_MPI ) xyz_PTempBZ(1:nx,1,1:nz) = Var3D(1:nx,1,1:nz) call SetMargin_xyz( xyz_PTempBZ ) name = "VelSoundBasicZ" call HistoryGet( InputFile, name, Var3D(:,1,:), flag_mpi_split = FLAG_LIB_MPI ) xyz_VelSoundBZ(1:nx,1,1:nz) = Var3D(1:nx,1,1:nz) call SetMargin_xyz( xyz_VelSoundBZ ) name = "TempBasicZ" call HistoryGet( InputFile, name, Var3D(:,1,:), flag_mpi_split = FLAG_LIB_MPI ) xyz_TempBZ(1:nx,1,1:nz) = Var3D(1:nx,1,1:nz) call SetMargin_xyz( xyz_TempBZ ) name = "PressBasicZ" call HistoryGet( InputFile, name, Var3D(:,1,:), flag_mpi_split = FLAG_LIB_MPI ) xyz_PressBZ(1:nx,1,1:nz) = Var3D(1:nx,1,1:nz) call SetMargin_xyz( xyz_PressBZ) ! name = "MixRtBasicZ" ! call HistoryGet( InputFile, name, Var4D(:,1,:,:), flag_mpi_split = FLAG_LIB_MPI ) ! xyzf_QMixBZ(1:nx,1,1:nz) = Var4D(1:nx,1,1:nz) xyzf_QMixBZ = 0.0d0 name = "EffMolWtBasicZ" call HistoryGet( InputFile, name, Var3D(:,1,:), flag_mpi_split = FLAG_LIB_MPI ) xyz_EffMolWtBZ(1:nx,1,1:nz) = Var3D(1:nx,1,1:nz) call SetMargin_xyz( xyz_EffMolWtBZ ) end subroutine Arare4fileio_MMC_BZ_Get
Subroutine : | |
pyz_VelX(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out) |
xqz_VelY(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out) |
xyr_VelZ(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out) |
xyz_PTemp(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out) |
xyz_Exner(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out) |
xyzf_QMix(imin:imax,jmin:jmax,kmin:kmax,ncmax) : | real(DP), intent(out) |
xyz_Km(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out) |
xyz_Kh(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out) |
xyz_CDens(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out) |
(out)
リスタートファイルから情報取得
subroutine Arare4fileio_MMC_Var_Get( pyz_VelX, xqz_VelY, xyr_VelZ, xyz_PTemp, xyz_Exner, xyzf_QMix, xyz_Km, xyz_Kh, xyz_CDens ) ! (out) ! !リスタートファイルから情報取得 ! !暗黙の型宣言禁止 implicit none !変数定義 real(DP) :: Var3D(1:nx, 1, 1:nz) ! real(DP) :: Var4D(1:nx, 1, 1:nz, 1:ncmax) real(DP), intent(out) :: pyz_VelX (imin:imax,jmin:jmax,kmin:kmax) real(DP), intent(out) :: xqz_VelY (imin:imax,jmin:jmax,kmin:kmax) real(DP), intent(out) :: xyr_VelZ (imin:imax,jmin:jmax,kmin:kmax) real(DP), intent(out) :: xyz_Exner (imin:imax,jmin:jmax,kmin:kmax) real(DP), intent(out) :: xyz_PTemp (imin:imax,jmin:jmax,kmin:kmax) real(DP), intent(out) :: xyz_Km (imin:imax,jmin:jmax,kmin:kmax) real(DP), intent(out) :: xyz_Kh (imin:imax,jmin:jmax,kmin:kmax) real(DP), intent(out) :: xyz_CDens (imin:imax,jmin:jmax,kmin:kmax) real(DP), intent(out) :: xyzf_QMix (imin:imax,jmin:jmax,kmin:kmax,ncmax) real(DP) :: RestartTime character(STRING) :: name !変数名 character(STRING) :: name1 character(STRING) :: Time character(STRING) :: InputFile !------------------------------------------------------------- ! Get a Value from netCDF File !------------------------------------------------------------- RestartTime = 0.0d0 time = 't=' // trim(toChar( RestartTime )) name = "PotTemp" name1 = "PotTempDist" InputFile = trim(Arare4Prefix) // trim(name) // ".nc" call HistoryGet( InputFile, name1, Var3D, range=Time, flag_mpi_split = FLAG_LIB_MPI ) xyz_PTemp(1:nx,1,1:nz) = Var3D(1:nx,1,1:nz) call SetMargin_xyz( xyz_PTemp ) name = "VelX" InputFile = trim(Arare4Prefix) // trim(name) // ".nc" call HistoryGet( InputFile, name, Var3D, range=Time, flag_mpi_split = FLAG_LIB_MPI ) pyz_VelX(1:nx,1,1:nz) = Var3D(1:nx,1,1:nz) call SetMargin_pyz( pyz_VelX ) name = "VelY" xqz_VelY = 0.0d0 name = "VelZ" InputFile = trim(Arare4Prefix) // trim(name) // ".nc" call HistoryGet( InputFile, name, Var3D, range=Time, flag_mpi_split = FLAG_LIB_MPI ) xyr_VelZ(1:nx,1,1:nz) = Var3D(1:nx,1,1:nz) call SetMargin_xyr( xyr_VelZ ) name = "Exner" InputFile = trim(Arare4Prefix) // trim(name) // ".nc" call HistoryGet( InputFile, name, Var3D, range=Time, flag_mpi_split = FLAG_LIB_MPI ) xyz_Exner(1:nx,1,1:nz) = Var3D(1:nx,1,1:nz) call SetMargin_xyz( xyz_Exner ) name = "Km" InputFile = trim(Arare4Prefix) // trim(name) // ".nc" call HistoryGet( InputFile, name, Var3D, range=Time, flag_mpi_split = FLAG_LIB_MPI ) xyz_Km(1:nx,1,1:nz) = Var3D(1:nx,1,1:nz) call SetMargin_xyz( xyz_Km ) name = "Kh" InputFile = trim(Arare4Prefix) // trim(name) // ".nc" call HistoryGet( InputFile, name, Var3D, range=Time, flag_mpi_split = FLAG_LIB_MPI ) xyz_Kh(1:nx,1,1:nz) = Var3D(1:nx,1,1:nz) call SetMargin_xyz( xyz_Kh ) name = "DensCloud" InputFile = trim(Arare4Prefix) // trim(name) // ".nc" call HistoryGet( InputFile, name, Var3D, range=Time, flag_mpi_split = FLAG_LIB_MPI ) xyz_CDens(1:nx,1,1:nz) = Var3D(1:nx,1,1:nz) call SetMargin_xyz( xyz_CDens ) ! name = "MixRt" ! InputFile = trim(Arare4Prefix) // trim(name) // ".nc" ! call HistoryGet( InputFile, name, Var4D, range=Time, flag_mpi_split = FLAG_LIB_MPI ) ! xyzf_QMix(1:nx,1,1:nz,1:ncmax) = Var4D(1:nx,1,1:nz,1:ncmax) xyzf_QMix(1:nx,1,1:nz,1:ncmax) = 0.0d0 end subroutine Arare4fileio_MMC_Var_Get
Subroutine : | |
pyz_VelX(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out) |
xqz_VelY(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out) |
xyr_VelZ(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out) |
xyz_PTemp(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out) |
xyz_Exner(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out) |
xyzf_QMix(imin:imax,jmin:jmax,kmin:kmax,ncmax) : | real(DP), intent(out) |
xyz_Km(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out) |
xyz_Kh(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out) |
xyz_CDens(imin:imax,jmin:jmax,kmin:kmax) : | real(DP), intent(out) |
(out)
リスタートファイルから情報取得
subroutine Arare4fileio_Var_Get( pyz_VelX, xqz_VelY, xyr_VelZ, xyz_PTemp, xyz_Exner, xyzf_QMix, xyz_Km, xyz_Kh, xyz_CDens ) ! (out) ! !リスタートファイルから情報取得 ! !暗黙の型宣言禁止 implicit none !変数定義 real(DP) :: Var3D(1:nx, 1, 1:nz) real(DP), intent(out) :: pyz_VelX (imin:imax,jmin:jmax,kmin:kmax) real(DP), intent(out) :: xqz_VelY (imin:imax,jmin:jmax,kmin:kmax) real(DP), intent(out) :: xyr_VelZ (imin:imax,jmin:jmax,kmin:kmax) real(DP), intent(out) :: xyz_Exner (imin:imax,jmin:jmax,kmin:kmax) real(DP), intent(out) :: xyz_PTemp (imin:imax,jmin:jmax,kmin:kmax) real(DP), intent(out) :: xyz_Km (imin:imax,jmin:jmax,kmin:kmax) real(DP), intent(out) :: xyz_Kh (imin:imax,jmin:jmax,kmin:kmax) real(DP), intent(out) :: xyz_CDens (imin:imax,jmin:jmax,kmin:kmax) real(DP), intent(out) :: xyzf_QMix (imin:imax,jmin:jmax,kmin:kmax,ncmax) real(DP) :: RestartTime integer :: f character(STRING) :: name !変数名 ! character(STRING) :: name1 character(STRING) :: Time character(STRING) :: InputFile !------------------------------------------------------------- ! Get a Value from netCDF File !------------------------------------------------------------- RestartTime = 0.0d0 time = 't=' // trim(toChar( RestartTime )) name = "PotTemp" ! name1 = "PotTempDist" InputFile = trim(Arare4Prefix) // trim(name) // ".nc" call HistoryGet( InputFile, name, Var3D, range=Time, flag_mpi_split = FLAG_LIB_MPI ) xyz_PTemp(1:nx,1,1:nz) = Var3D(1:nx,1,1:nz) call SetMargin_xyz( xyz_PTemp ) name = "VelX" InputFile = trim(Arare4Prefix) // trim(name) // ".nc" call HistoryGet( InputFile, name, Var3D, range=Time, flag_mpi_split = FLAG_LIB_MPI ) pyz_VelX(1:nx,1,1:nz) = Var3D(1:nx,1,1:nz) call SetMargin_pyz( pyz_VelX ) name = "VelY" xqz_VelY = 0.0d0 name = "VelZ" InputFile = trim(Arare4Prefix) // trim(name) // ".nc" call HistoryGet( InputFile, name, Var3D, range=Time, flag_mpi_split = FLAG_LIB_MPI ) xyr_VelZ(1:nx,1,1:nz) = Var3D(1:nx,1,1:nz) call SetMargin_xyr( xyr_VelZ ) name = "Exner" InputFile = trim(Arare4Prefix) // trim(name) // ".nc" call HistoryGet( InputFile, name, Var3D, range=Time, flag_mpi_split = FLAG_LIB_MPI ) xyz_Exner(1:nx,1,1:nz) = Var3D(1:nx,1,1:nz) call SetMargin_xyz( xyz_Exner ) name = "Km" InputFile = trim(Arare4Prefix) // trim(name) // ".nc" call HistoryGet( InputFile, name, Var3D, range=Time, flag_mpi_split = FLAG_LIB_MPI ) xyz_Km(1:nx,1,1:nz) = Var3D(1:nx,1,1:nz) call SetMargin_xyz( xyz_Km ) name = "Kh" InputFile = trim(Arare4Prefix) // trim(name) // ".nc" call HistoryGet( InputFile, name, Var3D, range=Time, flag_mpi_split = FLAG_LIB_MPI ) xyz_Kh(1:nx,1,1:nz) = Var3D(1:nx,1,1:nz) call SetMargin_xyz( xyz_Kh ) ! name = "MixRt" ! InputFile = trim(Arare4Prefix) // trim(name) // ".nc" ! call HistoryGet( InputFile, name, Var4D, range=Time, flag_mpi_split = FLAG_LIB_MPI ) ! xyzf_QMix(1:nx,1,1:nz,1:ncmax) = Var4D(1:nx,1,1:nz,1:ncmax) do f = 1, ncmax name = trim(SpcWetSymbol(f)) InputFile = trim(Arare4Prefix) // trim(SpcWetSymbol(f)) // ".nc" call HistoryGet( InputFile, name, Var3D, range=Time, flag_mpi_split = FLAG_LIB_MPI ) xyzf_QMix(1:nx,1,1:nz,f) = Var3D(1:nx,1,1:nz) end do call SetMargin_xyzf( xyzf_QMix ) name = "CloudDensity" xyz_CDens = 0.0d0 end subroutine Arare4fileio_Var_Get