| 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, er, k ! 作業用配列
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 : | |||
| 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(16) :: tmpc
character(10) :: 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(16) :: tmpc
character(10) :: 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(16) :: tmpc
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 : | |||
| 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, er, k ! 作業用配列
character(10) :: 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