| Class | file_operate |
| In: |
file_operate.f90
|
CReSS の計算結果出力ファイルデータの操作用モジュール
| Subroutine : | |||
| file_name : | character(*), intent(in)
| ||
| d2n : | integer, intent(in)
| ||
| d3n : | integer, intent(in)
| ||
| d2val : | character(*), intent(in)
| ||
| d3val : | character(*), intent(in)
| ||
| ad2val : | character(*), intent(in)
| ||
| ad3val : | character(*), intent(in)
| ||
| nx : | integer, intent(in)
| ||
| ny : | integer, intent(in)
| ||
| nz : | integer, intent(in)
| ||
| d2var : | real, dimension(nx,ny,d2n), intent(inout)
| ||
| d3var : | real, dimension(nx,ny,nz,d3n), intent(inout)
| ||
| specnz : | integer, intent(in), optional
| ||
| intvnz : | integer, intent(in), dimension(2), optional
|
読み込む変数の個数を指定することで, 不要な変数データはメモリに読み込まない ルーチン. (For CReSS) このルーチンを呼び出す前に必ず, val_counter で用意する変数の総数を求めておく.
subroutine auto_read_file( file_name, d2n, d3n, d2val, d3val, ad2val, ad3val, nx, ny, nz, d2var, d3var, specnz, intvnz )
! 読み込む変数の個数を指定することで, 不要な変数データはメモリに読み込まない
! ルーチン. (For CReSS)
! このルーチンを呼び出す前に必ず, val_counter で用意する変数の総数を求めておく.
implicit none
character(*), intent(in) :: file_name ! 読み出すファイル名
integer, intent(in) :: d2n ! d2 変数用に用意すべき配列総数 (変数の種類数)
integer, intent(in) :: d3n ! d3 変数用に用意すべき配列総数 (変数の種類数)
character(*), intent(in) :: d2val ! CReSS のダンプファイルの最初の d2 データ
character(*), intent(in) :: d3val ! CReSS のダンプファイルの最初の d3 データ
character(*), intent(in) :: ad2val ! CReSS のダンプファイルの後の d2 データ
character(*), intent(in) :: ad3val ! CReSS のダンプファイルの後の d3 データ
integer, intent(in) :: nx ! x 方向の要素数
integer, intent(in) :: ny ! y 方向の要素数
integer, intent(in) :: nz ! z 方向の要素数
integer, intent(in), optional :: specnz ! nz について指定高度のみ取り出す
! この場合, 必ず nz = 1 でなければならない.
integer, intent(in), dimension(2), optional :: intvnz
! nz について intvnz(1) から intvnz(2) までの高度を取り出す.
! この場合, nz はその変数の全要素数でなければならない.
real, dimension(nx,ny,d2n), intent(inout) :: d2var ! d2 変数用配列
real, dimension(nx,ny,nz,d3n), intent(inout) :: d3var ! d3 変数用配列
integer :: d2m, d3m, d2am, d3am, i, j, k, all_count, d2_count, d3_count
d2m=len_trim(d2val)
d3m=len_trim(d3val)
d2am=len_trim(ad2val)
d3am=len_trim(ad3val)
write(*,*) "len=", d2m, d3m, d2am, d3am
all_count=0
!-- ファイルの 1 レコード目から順に読むかどうかの判断を行って読む ---
!-- まず, 最初の d2 データについて
j=0
if(d2m>=1)then
do i=1,d2m
if(d2val(i:i)=='1')then
j=j+1
if(d2n<j)then
write(*,*) "*** ERROR ***"
write(*,*) "read array exceeds d2n. STOP."
stop
else
call read_file( file_name, nx, ny, i, d2var(:,:,j))
end if
end if
end do
end if
d2_count=j ! 最初の d2 データでいくら入ったか (読み飛ばしを入れない)
all_count=all_count+d2m ! 読み飛ばした分も含めて何レコード分読んだかのカウント
!-- 次, 最初の d3 データについて
j=0
if(d3m>=1)then
do i=1,d3m
if(d3val(i:i)=='1')then
j=j+1
if(d3n<j)then
write(*,*) "*** ERROR ***"
write(*,*) "read array exceeds d3n. STOP."
stop
else
if(present(specnz))then
call read_file( file_name, nx, ny, all_count+specnz, d3var(:,:,specnz,j) )
else
if(present(intvnz))then
do k=intvnz(1),intvnz(2)
call read_file( file_name, nx, ny, all_count+k, d3var(:,:,k,j) )
end do
else
do k=1,nz
call read_file( file_name, nx, ny, all_count+k, d3var(:,:,k,j) )
end do
end if
end if
end if
end if
all_count=all_count+nz ! if 文で読んでいなくても nz 行分カウント
end do
end if
d3_count=j ! 最初の d3 データでいくら入ったか (読み飛ばしを入れない)
!-- 続いて, 後の d2 データについて
j=d2_count
if(d2am>=1)then
do i=1,d2am
if(ad2val(i:i)=='1')then
j=j+1
if(d2n<j)then
write(*,*) "*** ERROR ***"
write(*,*) "read array exceeds d2n. STOP."
stop
else
call read_file( file_name, nx, ny, all_count+i, d2var(:,:,j))
end if
end if
end do
end if
all_count=all_count+d2am ! 読み飛ばした分も含めて何レコード分読んだかのカウント
!-- 最後に, 後の d3 データについて
j=d3_count
if(d3am>=1)then
do i=1,d3am
if(ad3val(i:i)=='1')then
j=j+1
if(d3n<j)then
write(*,*) "*** ERROR ***"
write(*,*) "read array exceeds d3n. STOP."
stop
else
if(present(specnz))then
call read_file( file_name, nx, ny, all_count+specnz, d3var(:,:,specnz,j) )
else
if(present(intvnz))then
do k=intvnz(1),intvnz(2)
call read_file( file_name, nx, ny, all_count+k, d3var(:,:,k,j) )
end do
else
do k=1,nz
call read_file( file_name, nx, ny, all_count+k, d3var(:,:,k,j) )
end do
end if
end if
end if
end if
all_count=all_count+nz ! if 文で読んでいなくても nz 行分カウント
end do
end if
end subroutine
| Function : | |||
| line_number_counter : | integer | ||
| fname : | character(*), intent(in)
|
テキストデータの行数を計算する関数
integer function line_number_counter( fname )
! テキストデータの行数を計算する関数
implicit none
character(*), intent(in) :: fname ! 読み取るテキストファイル
! integer, intent(inout) :: res
integer :: i, err
character(1) :: dummy
i=0
open(unit=10,file=trim(fname),iostat=err,status='old')
do while(err==0)
read(10,*,iostat=err) dummy
i=i+1
end do
close(unit=10,status='keep')
line_number_counter=i-1
return
end function
| Subroutine : | |||
| file_name : | character(*), intent(in)
| ||
| nx : | integer, intent(in)
| ||
| ny : | integer, intent(in)
| ||
| rec_num : | integer, intent(in)
| ||
| var(nx,ny) : | real, intent(inout)
|
出力結果読み取りルーチン 本ルーチンでは, ダイレクトアクセスを読み出す際, 1 変数のバイト数を 4 バイト と仮定して読み出すので, 4 バイト以外のファイルを読み出すときは, やはりソースファイルの書き換えが必要となる(要修正)
subroutine read_file( file_name, nx, ny, rec_num, var ) ! 出力結果読み取りルーチン
! 本ルーチンでは, ダイレクトアクセスを読み出す際, 1 変数のバイト数を 4 バイト
! と仮定して読み出すので, 4 バイト以外のファイルを読み出すときは,
! やはりソースファイルの書き換えが必要となる(要修正)
implicit none
integer, intent(in) :: nx ! データの x 方向の個数
integer, intent(in) :: ny ! データの y 方向の個数
integer, intent(in) :: rec_num ! 読み出すデータのレコード番号
character(*), intent(in) :: file_name ! 読み出すデータファイル名
real, intent(inout) :: var(nx,ny) ! 読み出すデータ
integer :: i, j ! 作業用配列
open(unit=11, file=file_name, access='direct', recl=4*nx*ny, status='old')
read(11,rec=rec_num) ((var(i,j),i=1,nx),j=1,ny)
close(unit=11, status='keep')
end subroutine read_file
| Subroutine : | |||
| file_name : | character(*), intent(in)
| ||
| nx : | integer, intent(in)
| ||
| ny : | integer, intent(in)
| ||
| nz : | integer, intent(in)
| ||
| rec_num : | integer, intent(in)
| ||
| var(nx,ny,nz) : | real, intent(inout)
|
read_file 3 次元版. 本ルーチンでは, ダイレクトアクセスを読み出す際, 1 変数のバイト数を 4 バイト と仮定して読み出すので, 4 バイト以外のファイルを読み出すときは, やはりソースファイルの書き換えが必要となる(要修正)
subroutine read_file_3d( file_name, nx, ny, nz, rec_num, var )
! read_file 3 次元版.
! 本ルーチンでは, ダイレクトアクセスを読み出す際, 1 変数のバイト数を 4 バイト
! と仮定して読み出すので, 4 バイト以外のファイルを読み出すときは,
! やはりソースファイルの書き換えが必要となる(要修正)
implicit none
integer, intent(in) :: nx ! データの x 方向の個数
integer, intent(in) :: ny ! データの y 方向の個数
integer, intent(in) :: nz ! データの z 方向の個数
integer, intent(in) :: rec_num ! 読み出しを開始するデータのレコード番号
character(*), intent(in) :: file_name ! 読み出すデータファイル名
real, intent(inout) :: var(nx,ny,nz) ! 読み出すデータ
integer :: i, j, k ! 作業用配列
open(unit=11, file=file_name, access='direct', recl=4*nx*ny, status='old')
do k=1,nz
read(11,rec=rec_num+k-1) ((var(i,j,k),i=1,nx),j=1,ny)
end do
close(unit=11, status='keep')
end subroutine read_file_3d
| Subroutine : | |||
| fname : | character(*), intent(in)
| ||
| val_name : | character(*), intent(in)
| ||
| dim_len : | integer, intent(in), dimension(4)
| ||
| intvz : | integer, intent(in), dimension(2)
| ||
| intvt : | integer, intent(in), dimension(2)
| ||
| val : | real, intent(inout),
dimension(dim_len(1),dim_len(2),intvz(1):intvz(2),intvt(1):intvt(2))
|
GrADS のコントロールファイルを読み込んで, val_name で指定した変数のデータを 取り出す.
subroutine read_file_grads( fname, val_name, dim_len, intvz, intvt, val )
! GrADS のコントロールファイルを読み込んで, val_name で指定した変数のデータを
! 取り出す.
implicit none
character(*), intent(in) :: fname ! GrADS のコントロールファイル
character(*), intent(in) :: val_name ! 読み出す変数名
integer, intent(in), dimension(4) :: dim_len ! 各次元要素の数 (x,y,z,t) の順
integer, intent(in), dimension(2) :: intvz ! 高度方向の読み出す要素区間
integer, intent(in), dimension(2) :: intvt ! 時間方向の読み出す要素区間
real, intent(inout), dimension(dim_len(1),dim_len(2),intvz(1):intvz(2),intvt(1):intvt(2)) :: val ! 読み出した値が格納される変数.
character(1000), allocatable, dimension(:) :: ctl_cval
character(100), allocatable, dimension(:,:) :: ctl_hval
integer :: nl, nx, ny, nz, nt, mon_num
integer :: i, j
character(1000) :: tmpc, ctmp, header, footer, iname
integer, dimension(10) :: cont_flag
integer :: t_fact, year, month, day, hour, counter_start, read_flag, counter
integer :: year_r, month_r, day_r, hour_r
logical, dimension(4) :: t_fact_flag, date_flag
nl=line_number_counter( trim(fname) )
nx=dim_len(1)
ny=dim_len(2)
nz=dim_len(3)
nt=dim_len(4)
allocate(ctl_cval(nl))
allocate(ctl_hval(10,nl))
ctl_cval=''
ctl_hval=''
tmpc=''
ctmp=''
header=''
footer=''
iname=''
!-- ctl file reading
call read_file_text( trim(fname), 1, nl, ctl_cval, forma='(a1000)' )
call read_file_text( trim(fname), 2, nl-1, ctl_hval(1:2,1:nl-1) )
!-- essential factor estimating
do i=1,nl
tmpc=trim(adjustl(ctl_cval(i)))
if(tmpc(1:4)=='vars')then
cont_flag(1)=i
end if
if(tmpc(1:7)=='endvars')then
cont_flag(2)=i
end if
if(tmpc(1:4)=='dset')then
cont_flag(4)=i
end if
if(tmpc(1:4)=='tdef')then
cont_flag(5)=i
end if
end do
call read_file_text( trim(fname), 4, 1, ctl_hval(1:4,cont_flag(5):cont_flag(5)), skip=cont_flag(5)-1 )
do i=cont_flag(1),cont_flag(2)
tmpc=trim(adjustl(ctl_hval(1,i)))
if(trim(tmpc)==trim(val_name))then
cont_flag(3)=i
end if
end do
!-- reading array estimate
read_flag=0
do i=cont_flag(1)+1,cont_flag(3)-1
read(ctl_hval(2,i),*) counter
if(counter==0)then
read_flag=read_flag+1
else
read_flag=read_flag+counter
end if
end do
!-- t_fact caluculate
t_fact_flag(1:4)=.false.
call read_file_text( trim(fname), 5, 1, ctl_hval(1:5,cont_flag(5)), skip=cont_flag(5)-1 )
tmpc=trim(ctl_hval(5,cont_flag(5)))
select case (tmpc(len_trim(adjustl(tmpc))-1:len_trim(adjustl(tmpc))))
case ('yr')
t_fact_flag(1)=.true.
read(tmpc(1:len_trim(adjustl(tmpc))-2),*) t_fact
case ('mn')
t_fact_flag(2)=.true.
read(tmpc(1:len_trim(adjustl(tmpc))-2),*) t_fact
case ('dy')
t_fact_flag(3)=.true.
read(tmpc(1:len_trim(adjustl(tmpc))-2),*) t_fact
case ('hr')
t_fact_flag(4)=.true.
read(tmpc(1:len_trim(adjustl(tmpc))-2),*) t_fact
end select
!-- file header estimate
counter=0
date_flag(1:4)=.false.
tmpc=ctl_hval(2,cont_flag(4))
do i=2,len_trim(adjustl(ctl_hval(2,cont_flag(4))))
if(tmpc(i:i)=='%')then
counter=counter+1
if(counter==1)then
counter_start=i
header(1:i-2)=tmpc(2:i-1)
end if
select case (tmpc(i+1:i+2))
case ('y4')
date_flag(1)=.true.
if(t_fact_flag(1).eqv..true.)then
footer=tmpc(i+3:len_trim(adjustl(ctl_hval(2,cont_flag(4)))))
end if
case ('m2')
date_flag(2)=.true.
if(t_fact_flag(1).eqv..true.)then
footer=tmpc(i+3:len_trim(adjustl(ctl_hval(2,cont_flag(4)))))
end if
case ('d2')
date_flag(3)=.true.
if(t_fact_flag(1).eqv..true.)then
footer=tmpc(i+3:len_trim(adjustl(ctl_hval(2,cont_flag(4)))))
end if
case ('h2')
date_flag(4)=.true.
if(t_fact_flag(1).eqv..true.)then
footer=tmpc(i+3:len_trim(adjustl(ctl_hval(2,cont_flag(4)))))
end if
end select
end if
end do
if(counter==0)then
tmpc=trim(adjustl(ctl_hval(2,cont_flag(4))))
header=trim(tmpc(2:len_trim(adjustl(ctl_hval(2,cont_flag(4))))))
footer=''
end if
!-- initial time reading
tmpc=trim(ctl_hval(4,cont_flag(5)))
read(tmpc(1:2),*) hour
select case (tmpc(3:3))
case (':')
read(tmpc(7:8),*) day
read(tmpc(12:15),*) year
mon_num=9 ! 月を表す文字が始まる位置
case ('Z')
read(tmpc(4:5),*) day
read(tmpc(9:12),*) year
mon_num=6 ! 月を表す文字が始まる位置
end select
select case (tmpc(mon_num:mon_num+2))
case('Jan')
month=1
case('Feb')
month=2
case('Mar')
month=3
case('Apr')
month=4
case('May')
month=5
case('Jun')
month=6
case('Jul')
month=7
case('Aug')
month=8
case('Sep')
month=9
case('Oct')
month=10
case('Nov')
month=11
case('Dec')
month=12
end select
year_r=year
month_r=month
day_r=day
hour_r=hour
!-- reading data
do i=intvt(1),intvt(2)
!-- default value setting
iname=trim(header)
do j=4,1,-1 ! 時間から足し合わせた方が都合がよいためループ逆回し
select case (j)
case (1)
if(t_fact_flag(1).eqv..true.)then
year=year_r+t_fact*(i-1)
end if
case (2)
if(t_fact_flag(2).eqv..true.)then
month=month_r+t_fact*(i-1)
end if
case (3)
if(t_fact_flag(3).eqv..true.)then
day=day_r+t_fact*(i-1)
end if
case (4)
if(t_fact_flag(4).eqv..true.)then
hour=hour_r+t_fact*(i-1)
end if
end select
do while (hour>24)
hour=hour-24
day=day+1
end do
do while (day>28)
do while (month>12)
year=year+1
month=month-12
end do
select case (month)
case (1)
if(day>31)then
month=month+1
day=day-31
else
exit
end if
case (3)
if(day>31)then
month=month+1
day=day-31
else
exit
end if
case (5)
if(day>31)then
month=month+1
day=day-31
else
exit
end if
case (7)
if(day>31)then
month=month+1
day=day-31
else
exit
end if
case (8)
if(day>31)then
month=month+1
day=day-31
else
exit
end if
case (10)
if(day>31)then
month=month+1
day=day-31
else
exit
end if
case (12)
if(day>31)then
month=month+1
day=day-31
else
exit
end if
case (4)
if(day>30)then
month=month+1
day=day-30
else
exit
end if
case (6)
if(day>30)then
month=month+1
day=day-30
else
exit
end if
case (9)
if(day>30)then
month=month+1
day=day-30
else
exit
end if
case (11)
if(day>30)then
month=month+1
day=day-30
else
exit
end if
case (2)
if(day>29.and.mod(year,4)==0)then
month=month+1
day=day-29
else
if(day>28.and.mod(year,4)/=0)then
month=month+1
day=day-28
else
exit
end if
end if
end select
end do
end do
do j=1,4
if(date_flag(j).eqv..true.)then
select case (j)
case (1)
write(ctmp,'(i4.4)') year
iname=trim(iname)//trim(adjustl(ctmp))
case (2)
write(ctmp,'(i2.2)') month
iname=trim(iname)//trim(adjustl(ctmp))
case (3)
write(ctmp,'(i2.2)') day
iname=trim(iname)//trim(adjustl(ctmp))
case (4)
write(ctmp,'(i2.2)') hour
iname=trim(iname)//trim(adjustl(ctmp))
end select
end if
end do
iname=trim(iname)//trim(footer)
do j=intvz(1),intvz(2)
call read_file( trim(iname), nx, ny, read_flag+j, val(:,:,j,i) )
end do
write(*,*) "finishing to read file", trim(iname), nx, ny, read_flag+intvz(1), read_flag+intvz(2), i
end do
end subroutine
| Subroutine : | |||
| fname : | character(*), intent(in)
| ||
| nx : | integer, intent(in)
| ||
| ny : | integer, intent(in)
| ||
| nz : | integer, intent(in)
| ||
| tstr : | integer, intent(in)
| ||
| tstep : | integer, intent(in)
| ||
| tnum : | integer, intent(in)
| ||
| val(nx,ny,nz,tnum) : | real, intent(inout) |
gtool3 形式のデータを読み出す.
subroutine read_file_gtool3( fname, nx, ny, nz, tstr, tstep, tnum, val )
! gtool3 形式のデータを読み出す.
implicit none
character(*), intent(in) :: fname ! 読み込むデータファイル
integer, intent(in) :: nx ! x 方向の配列数
integer, intent(in) :: ny ! y 方向の配列数
integer, intent(in) :: nz ! z 方向の配列数
integer, intent(in) :: tstr ! 読み込む時刻 (初期値から何ステップ目かで指定).
integer, intent(in) :: tstep ! 読み飛ばす時間間隔 (1 で毎時読み込み)
integer, intent(in) :: tnum ! 読み込む時間の数
real, intent(inout) :: val(nx,ny,nz,tnum)
character(16) :: header(64) ! ヘッダー情報読み飛ばし用
integer :: i, j, k, l, tmax, counter
real :: skip(nx,ny,nz)
tmax=tstr+tstep*(tnum-1)
open(unit=11,file=trim(fname),form='unformatted',access='sequential',status='old')
if(tstr/=1)then
do i=1,tstr-1
read(11) header
read(11) skip
end do
end if
counter=0
do i=tstr,tmax
read(11) header
read(11) skip
if(mod((i-tstr),tstep)==0)then
counter=counter+1
do l=1,nz
do k=1,ny
do j=1,nx
val(j,k,l,counter)=skip(j,k,l)
end do
end do
end do
end if
end do
close(unit=11)
end subroutine
| Subroutine : | |
| fname : | character(*), intent(in) |
| header(:) : | character(*), intent(in) |
| val(size(header)) : | character(16), intent(inout) |
gtool3 形式のデータのヘッダーから, 必要な変数名を入力し, その値を返す. 参照変数名は www.riam.kyushu-u.ac.jp/taikai/lab/others/Gtool/node111.html 参照. ただし, Fortran の予約語である unit, size については, uni, siz を使用.
Alias for read_file_gtool3_header_c
| Subroutine : | |
| fname : | character(*), intent(in) |
| header(:) : | character(*), intent(in) |
| val(size(header)) : | integer, intent(inout) |
gtool3 形式のデータのヘッダーから, 必要な変数名を入力し, その値を返す. 参照変数名は www.riam.kyushu-u.ac.jp/taikai/lab/others/Gtool/node111.html 参照. ただし, Fortran の予約語である unit, size については, uni, siz を使用.
Alias for read_file_gtool3_header_i
| Subroutine : | |
| fname : | character(*), intent(in) |
| header(:) : | character(*), intent(in) |
| val(size(header)) : | real, intent(inout) |
gtool3 形式のデータのヘッダーから, 必要な変数名を入力し, その値を返す. 参照変数名は www.riam.kyushu-u.ac.jp/taikai/lab/others/Gtool/node111.html 参照. ただし, Fortran の予約語である unit, size については, uni, siz を使用.
Alias for read_file_gtool3_header_r
| Subroutine : | |
| fname : | character(*), intent(in) |
| header(:) : | character(*), intent(in) |
| val(size(header)) : | character(16), intent(inout) |
gtool3 形式のデータのヘッダーから, 必要な変数名を入力し, その値を返す. 参照変数名は www.riam.kyushu-u.ac.jp/taikai/lab/others/Gtool/node111.html 参照. ただし, Fortran の予約語である unit, size については, uni, siz を使用.
subroutine read_file_gtool3_header_c( fname, header, val )
! gtool3 形式のデータのヘッダーから, 必要な変数名を入力し, その値を返す.
! 参照変数名は http://www.riam.kyushu-u.ac.jp/taikai/lab/others/Gtool/node111.html 参照.
! ただし, Fortran の予約語である unit, size については, uni, siz を使用.
implicit none
character(*), intent(in) :: fname
character(*), intent(in) :: header(:)
character(16), intent(inout) :: val(size(header))
integer, parameter :: vnum=43
character(16) :: raw_dat(64)
character(6) :: signc(vnum)
integer :: signi(vnum)
integer :: i, j
signc=(/'dset ', 'item ', 'edit1 ', 'edit2 ', 'edit3 ', 'edit4 ', 'edit5 ', 'edit6 ', 'edit7 ', 'edit8 ', 'titl1 ', 'titl2 ', 'uni ', 'ettl1 ', 'ettl2 ', 'ettl3 ', 'ettl4 ', 'ettl5 ', 'ettl6 ', 'ettl7 ', 'ettl8 ', 'date ', 'utim ', 'aitm1 ', 'aitm2 ', 'aitm3 ', 'dfmt ', 'coptn ', 'utim2 ', 'memo01', 'memo02', 'memo03', 'memo04', 'memo05', 'memo06', 'memo07', 'memo08', 'memo09', 'memo10', 'cdate ', 'csign ', 'mdate ', 'msign ' /)
signi=(/2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 14, 15, 16, 17, 18, 19, 20, 21, 22, 23, 24, 26, 27, 29, 32, 35, 38, 45, 49, 50, 51, 52, 53, 54, 55, 56, 57, 58, 59, 60, 61, 62, 63 /)
open(unit=11,file=trim(fname),form='unformatted',access='sequential',status='old')
read(11) raw_dat
close(unit=11)
do i=1,size(header)
do j=1,vnum
if(trim(header(i))==trim(signc(j)))then
val(i)=raw_dat(signi(j))
exit
end if
end do
end do
end subroutine
| Subroutine : | |
| fname : | character(*), intent(in) |
| header(:) : | character(*), intent(in) |
| val(size(header)) : | integer, intent(inout) |
gtool3 形式のデータのヘッダーから, 必要な変数名を入力し, その値を返す. 参照変数名は www.riam.kyushu-u.ac.jp/taikai/lab/others/Gtool/node111.html 参照. ただし, Fortran の予約語である unit, size については, uni, siz を使用.
subroutine read_file_gtool3_header_i( fname, header, val )
! gtool3 形式のデータのヘッダーから, 必要な変数名を入力し, その値を返す.
! 参照変数名は http://www.riam.kyushu-u.ac.jp/taikai/lab/others/Gtool/node111.html 参照.
! ただし, Fortran の予約語である unit, size については, uni, siz を使用.
implicit none
character(*), intent(in) :: fname
character(*), intent(in) :: header(:)
integer, intent(inout) :: val(size(header))
integer, parameter :: vnum=15
character(16) :: raw_dat(64)
character(6) :: signc(vnum)
integer :: signi(vnum)
integer :: i, j
signc=(/'idfm ', 'fnum ', 'dnum ', 'time ', 'tdur ', 'astr1 ', 'aend1 ', 'astr2 ', 'aend2 ', 'astr3 ', 'aend3 ', 'styp ', 'ioptn ', 'time2 ', 'size ' /)
signi=(/1, 12, 13, 25, 28, 30, 31, 33, 34, 36, 37, 44, 46, 48, 64 /)
open(unit=11,file=trim(fname),form='unformatted',access='sequential',status='old')
read(11) raw_dat
close(unit=11)
do i=1,size(header)
do j=1,vnum
if(trim(header(i))==trim(signc(j)))then
read(raw_dat(signi(j)),'(i16)') val(i)
exit
end if
end do
end do
end subroutine
| Subroutine : | |
| fname : | character(*), intent(in) |
| header(:) : | character(*), intent(in) |
| val(size(header)) : | real, intent(inout) |
gtool3 形式のデータのヘッダーから, 必要な変数名を入力し, その値を返す. 参照変数名は www.riam.kyushu-u.ac.jp/taikai/lab/others/Gtool/node111.html 参照. ただし, Fortran の予約語である unit, size については, uni, siz を使用.
subroutine read_file_gtool3_header_r( fname, header, val )
! gtool3 形式のデータのヘッダーから, 必要な変数名を入力し, その値を返す.
! 参照変数名は http://www.riam.kyushu-u.ac.jp/taikai/lab/others/Gtool/node111.html 参照.
! ただし, Fortran の予約語である unit, size については, uni, siz を使用.
implicit none
character(*), intent(in) :: fname
character(*), intent(in) :: header(:)
real, intent(inout) :: val(size(header))
integer, parameter :: vnum=6
character(16) :: raw_dat(64)
character(10) :: signc(vnum)
integer :: signi(vnum)
integer :: i, j
signc=(/'miss ', 'dmin ', 'dmax ', 'divs ', 'divl ', 'roptn ' /)
signi=(/39, 40, 41, 42, 43, 47 /)
open(unit=11,file=trim(fname),form='unformatted',access='sequential',status='old')
read(11) raw_dat
close(unit=11)
do i=1,size(header)
do j=1,vnum
if(trim(header(i))==trim(signc(j)))then
if(trim(header(i))=='roptn')then
read(raw_dat(signi(j)),'(e16.7)') val(i)
else
read(raw_dat(signi(j)),'(e15.7)') val(i)
end if
exit
end if
end do
end do
end subroutine
| Subroutine : | |||
| fname : | character(*), intent(in)
| ||
| nx : | integer, intent(in)
| ||
| ny : | integer, intent(in)
| ||
| val(nx,ny) : | character(*), intent(inout)
| ||
| skip : | integer, intent(in), optional
| ||
| forma : | character(*), intent(in), optional
|
テキスト形式のデータを読み込むルーチン forma と nx, ny の関係は 4x3 の i3 (空白込み) の行列データの場合, forma = ’(3a3)’ という form で読み込むなら, nx = 3, ny = 4 forma = ’(a)’ という form で読み込むなら, nx = 1, ny = 4 となる.
subroutine read_file_text( fname, nx, ny, val, skip, forma )
! テキスト形式のデータを読み込むルーチン
! forma と nx, ny の関係は 4x3 の i3 (空白込み) の行列データの場合,
! forma = '(3a3)' という form で読み込むなら, nx = 3, ny = 4
! forma = '(a)' という form で読み込むなら, nx = 1, ny = 4 となる.
implicit none
character(*), intent(in) :: fname ! 読み込むテキストファイル
integer, intent(in) :: nx ! 読み込むファイルの列数
integer, intent(in) :: ny ! 読み込むファイルの行数
character(*), intent(inout) :: val(nx,ny) ! 文字列データにも対応するため, 文字列で結果を返す.
integer, intent(in), optional :: skip ! 行頭何行を読み飛ばすか. default は 0.
character(*), intent(in), optional :: forma ! read のフォーマット '(f4.1)' など. デフォルトでは read のデフォルトフォーマット使用.
integer :: i, j
character(1) :: dummy
open(unit=10, file=trim(fname), status='old')
if(present(skip))then
do i=1,skip
read(10,*) dummy
end do
end if
if(present(forma))then
do j=1,ny
read(10,forma) (val(i,j),i=1,nx)
end do
else
do j=1,ny
read(10,*) (val(i,j),i=1,nx)
end do
end if
close(unit=10,status='keep')
end subroutine
| Subroutine : | |||
| fname : | character(*), intent(in)
| ||
| sst : | real, intent(inout), dimension(nx,ny)
|
mgdsst データを読み込むためのルーチン 読み込み用データは 1440 x 720 で固定. これは, mgdsst のフォーマットからくるものである.
subroutine read_mgdsst( fname, sst )
! mgdsst データを読み込むためのルーチン
! 読み込み用データは 1440 x 720 で固定.
! これは, mgdsst のフォーマットからくるものである.
implicit none
integer, parameter :: nx=1440 ! 経度方向の格子点数
integer, parameter :: ny=720 ! 緯度方向の格子点数
character(*), intent(in) :: fname ! 読み込む mgdsst ルーチン
real, intent(inout), dimension(nx,ny) :: sst ! 読み込んだ SST データ [K].
! ただし, 999.0 が未定義, 888.0 が海氷.
character(3), dimension(nx,ny) :: sstc
character(10) :: formatt
integer :: j, k
formatt='(1440a3)'
call read_file_text( trim(fname), nx, ny, sstc, skip=1, forma=trim(formatt))
! open(unit=12,file=trim(dat_list(i)),status='old')
! read(12,*) sstc(1,1)
do k=1,ny
! read(12,formatt) (sstc(j,k),j=1,nx)
do j=1,nx
if(sstc(j,k)=='888'.or.sstc(j,k)=='999')then
sst(j,ny-k+1)=999.0
else
read(sstc(j,k),*) sst(j,ny-k+1)
sst(j,ny-k+1)=sst(j,ny-k+1)/10.0+273.15
! write(*,*) j,k,sstc(j,k),sst(j,k)
end if
! tmp(j,ny-k+1)=sst(j,k)
! sst(j,k)=tmp(j,k)
end do
end do
end subroutine
| Subroutine : | |||
| fname : | character(*), intent(in)
| ||
| axis_name : | character(*), intent(in)
| ||
| x(:) : | real, intent(inout)
| ||
| y(:) : | real, optional, intent(inout)
|
GrADS のコントロールファイルから 1 次元方向の軸データを取得する. 現在, linear, levels タイプのみ対応.
subroutine request_axis_grads( fname, axis_name, x, y )
! GrADS のコントロールファイルから 1 次元方向の軸データを取得する.
! 現在, linear, levels タイプのみ対応.
implicit none
character(*), intent(in) :: fname ! GrADS のコントロールファイル
character(*), intent(in) :: axis_name ! 読み出す軸の名前.
! xdef, ydef, zdef の 3 パターン
real, intent(inout) :: x(:) ! 軸の座標値
real, optional, intent(inout) :: y(:) ! pdef の場合の y 軸の座標値
character(1000), allocatable, dimension(:) :: ctl_cval
character(100), allocatable, dimension(:) :: ctl_vval
character(100), allocatable, dimension(:,:) :: ctl_hval
integer :: nl, i, j, nx, ny
real :: vmin, dv
character(1000) :: tmpc
character(100), allocatable, dimension(:) :: dummy
nl=line_number_counter( trim(fname) )
nx=size(x)
if(present(y))then
if(trim(axis_name)=='pdef')then
ny=size(y)
else
write(*,*) "ERROR : request_axis_grads."
write(*,*) "If axis_name is pdef, argument y must need."
write(*,*) "STOP."
stop
end if
end if
allocate(ctl_cval(nl))
allocate(ctl_hval(2,nl))
allocate(ctl_vval(13))
allocate(dummy(nx+3))
ctl_cval=''
ctl_vval=''
ctl_hval=''
tmpc=''
dummy=''
call read_file_text( trim(fname), 1, nl, ctl_cval, forma='(a1000)' )
call read_file_text( trim(fname), 2, nl-1, ctl_hval(1:2,1:nl-1) )
do i=1,nl
tmpc=trim(adjustl(ctl_cval(i)))
if(tmpc(1:4)==trim(axis_name))then
if(tmpc(1:4)/='pdef')then
call read_file_text( trim(fname), 5, 1, ctl_vval(1:5), skip=i-1 )
if(trim(ctl_vval(3))=='linear')then
read(ctl_vval(4),*) vmin
read(ctl_vval(5),*) dv
x=(/((vmin+dv*(i-1)),i=1,nx)/)
else
if(trim(ctl_vval(3))=='levels')then
call read_file_text( trim(fname), nx+3, 1, dummy, skip=i-1 )
do j=4,nx+3
read(dummy(j),*) x(j-3)
end do
else
write(*,*) "ERROR (request_axis_grads) : ", trim(ctl_vval(3))
write(*,*) "Mode is neither 'levels' nor 'linear'."
write(*,*) "STOP."
stop
end if
end if
else
call read_file_text( trim(fname), 13, 1, ctl_vval(1:13), skip=i-1 )
read(ctl_vval(7),*) vmin
read(ctl_vval(12),*) dv
x=(/((vmin+dv*(i-1)),i=1,nx)/)
read(ctl_vval(8),*) vmin
read(ctl_vval(13),*) dv
y=(/((vmin+dv*(i-1)),i=1,ny)/)
end if
end if
end do
end subroutine
| Subroutine : | |||
| fname : | character(*), intent(in)
| ||
| dim_len : | integer, intent(inout), dimension(4)
|
GrADS のコントロールファイルから次元要素数を取得する.
subroutine request_dim_grads( fname, dim_len )
! GrADS のコントロールファイルから次元要素数を取得する.
implicit none
character(*), intent(in) :: fname ! GrADS のコントロールファイル
integer, intent(inout), dimension(4) :: dim_len ! x,y,z,t 方向の要素数
character(1000), allocatable, dimension(:) :: ctl_cval
character(100), allocatable, dimension(:,:) :: ctl_hval
integer, dimension(5) :: cont_flag
integer :: nl, i
character(1000) :: tmpc
logical :: pdeflag
nl=line_number_counter( trim(fname) )
allocate(ctl_cval(nl))
allocate(ctl_hval(10,nl))
pdeflag=.false.
ctl_cval=''
ctl_hval=''
tmpc=''
call read_file_text( trim(fname), 1, nl, ctl_cval, forma='(a1000)' )
call read_file_text( trim(fname), 2, nl-1, ctl_hval(1:2,1:nl-1) )
do i=1,nl
tmpc=trim(adjustl(ctl_cval(i)))
if(tmpc(1:4)=='xdef')then
cont_flag(1)=i
end if
if(tmpc(1:4)=='ydef')then
cont_flag(2)=i
end if
if(tmpc(1:4)=='zdef')then
cont_flag(3)=i
end if
if(tmpc(1:4)=='tdef')then
cont_flag(4)=i
end if
if(tmpc(1:4)=='pdef')then
cont_flag(5)=i
pdeflag=.true.
end if
end do
if(pdeflag.eqv..true.)then
read(ctl_hval(2,cont_flag(1)),*) dim_len(1)
read(ctl_hval(2,cont_flag(2)),*) dim_len(2)
read(ctl_hval(2,cont_flag(3)),*) dim_len(3)
read(ctl_hval(2,cont_flag(4)),*) dim_len(4)
else
read(ctl_hval(2,cont_flag(5)),*) dim_len(1)
read(ctl_hval(3,cont_flag(5)),*) dim_len(2)
read(ctl_hval(2,cont_flag(3)),*) dim_len(3)
read(ctl_hval(2,cont_flag(4)),*) dim_len(4)
end if
end subroutine
| Subroutine : | |||
| fname : | character(*), intent(in)
| ||
| var_name : | character(*), intent(in)
| ||
| nz : | integer, intent(inout)
|
GrADS のコントロールファイルから var_name で指定された変数の高度方向の要素数 を取得する.
subroutine request_vardim_grads( fname, var_name, nz )
! GrADS のコントロールファイルから var_name で指定された変数の高度方向の要素数
! を取得する.
implicit none
character(*), intent(in) :: fname ! GrADS のコントロールファイル
character(*), intent(in) :: var_name ! 要素数を取得する変数の名前.
integer, intent(inout) :: nz ! 要素数
character(1000), allocatable, dimension(:) :: ctl_cval, ctl_vval
character(100), allocatable, dimension(:,:) :: ctl_hval
integer :: nl, i, counter
character(1000) :: tmpc
nl=line_number_counter( trim(fname) )
allocate(ctl_cval(nl))
allocate(ctl_hval(2,nl))
allocate(ctl_vval(5))
ctl_vval=''
ctl_hval=''
ctl_cval=''
tmpc=''
call read_file_text( trim(fname), 1, nl, ctl_cval, forma='(a1000)' )
call read_file_text( trim(fname), 2, nl-1, ctl_hval(1:2,1:nl-1) )
counter=0
do i=1,nl
tmpc=trim(adjustl(ctl_hval(1,i)))
if(trim(tmpc)=='vars')then
counter=1
end if
if(trim(tmpc)==trim(var_name).and.counter==1)then
read(ctl_hval(2,i),*) nz
exit
end if
end do
end subroutine
| Subroutine : | |||
| file_name : | character(*), intent(in)
| ||
| nx : | integer, intent(in)
| ||
| ny : | integer, intent(in)
| ||
| rec_num : | integer, intent(in)
| ||
| var(nx,ny) : | real, intent(inout)
| ||
| mode : | character(*), optional, intent(in)
|
解析出力ルーチン 本ルーチンでは, ダイレクトアクセスを読み出す際, 1 変数のバイト数を 4 バイト と仮定して読み出すので, 4 バイト以外のファイルを読み出すときは, やはりソースファイルの書き換えが必要となる(要修正)
subroutine write_file( file_name, nx, ny, rec_num, var, mode ) ! 解析出力ルーチン
! 本ルーチンでは, ダイレクトアクセスを読み出す際, 1 変数のバイト数を 4 バイト
! と仮定して読み出すので, 4 バイト以外のファイルを読み出すときは,
! やはりソースファイルの書き換えが必要となる(要修正)
implicit none
integer, intent(in) :: nx ! データの x 方向の個数
integer, intent(in) :: ny ! データの y 方向の個数
integer, intent(in) :: rec_num ! 読み出すデータのレコード番号
character(*), intent(in) :: file_name ! 読み出すデータファイル名
real, intent(inout) :: var(nx,ny) ! 読み出すデータ
character(*), optional, intent(in) :: mode ! ファイルの書き出しオプション
integer :: i, j ! 作業用配列
character(10) :: cmode
cmode=''
if(present(mode))then
cmode=mode
else
cmode='unknown'
end if
open(unit=11, file=file_name, access='direct', recl=4*nx*ny, status=cmode)
write(11,rec=rec_num) ((var(i,j),i=1,nx),j=1,ny)
close(unit=11, status='keep')
end subroutine write_file
| Subroutine : | |||
| file_name : | character(*), intent(in)
| ||
| nx : | integer, intent(in)
| ||
| ny : | integer, intent(in)
| ||
| nz : | integer, intent(in)
| ||
| rec_num : | integer, intent(in)
| ||
| var(nx,ny,nz) : | real, intent(inout)
| ||
| mode : | character(*), optional, intent(in)
|
write_file 3 次元版. 本ルーチンでは, ダイレクトアクセスを読み出す際, 1 変数のバイト数を 4 バイト と仮定して読み出すので, 4 バイト以外のファイルを読み出すときは, やはりソースファイルの書き換えが必要となる(要修正)
subroutine write_file_3d( file_name, nx, ny, nz, rec_num, var, mode )
! write_file 3 次元版.
! 本ルーチンでは, ダイレクトアクセスを読み出す際, 1 変数のバイト数を 4 バイト
! と仮定して読み出すので, 4 バイト以外のファイルを読み出すときは,
! やはりソースファイルの書き換えが必要となる(要修正)
implicit none
integer, intent(in) :: nx ! データの x 方向の個数
integer, intent(in) :: ny ! データの y 方向の個数
integer, intent(in) :: nz ! データの z 方向の個数
integer, intent(in) :: rec_num ! 読み出しを開始するデータのレコード番号
character(*), intent(in) :: file_name ! 読み出すデータファイル名
real, intent(inout) :: var(nx,ny,nz) ! 読み出すデータ
character(*), optional, intent(in) :: mode ! ファイルの書き出しオプション
integer :: i, j, k ! 作業用配列
character(10) :: cmode
cmode=''
if(present(mode))then
cmode=mode
else
cmode='unknown'
end if
open(unit=11, file=file_name, access='direct', recl=4*nx*ny, status=cmode)
do k=1,nz
write(11,rec=rec_num+k-1) ((var(i,j,k),i=1,nx),j=1,ny)
end do
close(unit=11, status='keep')
end subroutine write_file_3d