!== Sample program for gt4_history/gt4f90io ! ! Authors:: Yasuhiro MORIKAWA, Eizi TOYODA ! Version:: $Id: histtest.f90,v 1.9 2007/08/22 18:34:28 morikawa Exp $ ! Tag Name:: $Name: gt4f90io-20070827 $ ! Copyright:: Copyright (C) GFD Dennou Club, 2000-2005. All rights reserved. ! License:: See COPYRIGHT[link:../COPYRIGHT] ! program histtest ! ! Test Program for "HistoryCreate", "HistoryAddVariable", ! "HistoryPut", "HistorySetTime", "HistoryClose". ! use dc_types, only: STRING use dc_trace, only: SetDebug use gt4_history, only: GT_HISTORY, HistoryCreate, HistoryAddVariable, & & HistoryClose, HistoryCopy, HistoryPut, & & HistorySetTime, HistoryAddAttr use dc_string, only: StoA, CPrintf implicit none integer:: i character(STRING) :: char type(GT_HISTORY) :: gt4hist_test_5, gt4hist_test_7, gt4hist_test_8 real(8) :: Zeta(4,8) continue !----- デバッグモードへ ----- call SetDebug !----------------------------------------------------------------- ! HistoryCreate1 を用いた基本出力テスト ! (Conventions と gt_version は自動出力) !----------------------------------------------------------------- call HistoryCreate(file='xhisttest/xhisttest1.nc', & & title='gt4_history test 1', & & source='gt4f90io/Fortran library test kit', & & institution='GFD Dennou Club', & & dims=(/'time'/), & & dimsizes=(/0/), & & longnames=(/'time'/), & & units=(/'s @ 2001-12-16T10:00'/), & & origin=0.0, interval=0.0, & & xtypes=(/'real'/)) call HistoryAddVariable('u', dims=(/'time'/), longname='any quantity', & units='non-dimensional') call HistoryAddVariable('scalar', dims=(/''/), longname='scalar quantity', & units='1') call HistoryPut('scalar', 999.9) do, i = 1, 24 call HistoryPut('time', real(i)) call HistoryPut('u', real(i) * 10.0) enddo call HistoryClose !----------------------------------------------------------------- ! HistoryCreate1 を用いた出力テスト ! 時間は gt4_history 内部の TimeGoAhead にて進める。 ! ここでは時間を気にせずただ HistoryPut を呼んでみる。 ! この場合、HistoryCreate に渡す origin と interval が使われる。 ! (gt_version のみを手動で出力) !----------------------------------------------------------------- call HistoryCreate(file='xhisttest/xhisttest2.nc', & & title='gt4_history test 2', & & source='gt4f90io/Fortran library test kit', & & institution='GFD Dennou Club', & & gt_version="4.0", & & dims=(/'x ', 'time'/), & & dimsizes=(/3, 0/), & & longnames=(/'eastward length', 'time '/), & & units=(/'m','s'/), & & origin=100.0, interval=2.5, & & xtypes=(/'real', 'real'/)) call HistoryPut('x', (/100.0, 200.0, 300.0/)) call HistoryAddVariable('u', dims=(/'x ', 'time'/), & & longname='any quantity', units='m/s') do, i = 1, 3 call HistoryPut('u', (/1.0*i, 2.0*i, 3.0*i/)) enddo call HistoryClose !----------------------------------------------------------------- ! HistoryCreate1 を用いた出力テスト ! 時間は HistoryPut にて与える。 ! (Conventions と gt_version は手動で出力) ! (origin と interval は自動設定) !----------------------------------------------------------------- call HistoryCreate(file='xhisttest/xhisttest3.nc', & & title='gt4_history test 3', & & source='gt4f90io/Fortran library test kit', & & institution='GFD DennouClub', & & conventions="http://www.gfd-dennou.org/library/gtool4/conventions/", & & gt_version="4.1", & & dims=(/'x ', 'time'/), & & dimsizes=(/3, 0/), & & longnames=(/'eastward length', 'time '/), & & units=(/'m','s'/), & & xtypes=(/'float', 'float'/) ) call HistoryAddVariable('u', dims=(/'x ', 'time'/), & & longname='foo quantity', units='m/s') call HistoryAddVariable('v', dims=(/'x ', 'time'/), & & longname='foo quantity', units='m/s', xtype='int') call HistoryPut('x', (/0.0, 1.0, 2.0/)) call HistoryPut('time', 0.0) call HistoryPut('u', (/0.0, 0.1, 0.2/)) call HistoryPut('v', (/10, 11, 12/)) call HistoryPut('time', 1.0) call HistoryPut('u', (/2.0, 2.1, 2.2/)) call HistoryPut('time', 2.0) call HistoryPut('u', (/3.0, 3.1, 3.2/)) call HistoryPut('v', (/40, 41, 42/)) call HistoryClose !----------------------------------------------------------------- ! HistoryCreate1 を用いた出力テスト ! 時間は HistorySetTime にて与える。 ! (Conventions のみを手動で出力) !----------------------------------------------------------------- call HistoryCreate(file='xhisttest/xhisttest4.nc', & & title='gt4_history test 4', & & source='gt4f90io/Fortran library test kit', & & institution='GFD DennouClub', & & conventions="Another certain Conventions", & & dims=(/'x ', 'time'/), & & dimsizes=(/3, 0/), & & longnames=(/'eastward length', 'time '/), & & units=(/'m','s'/), & & origin=0.0, interval=0.0, & & xtypes=(/'float', 'float'/) ) call HistoryAddVariable('u', dims=(/'x ', 'time'/), & & longname='foo quantity', units='m/s') call HistoryAddVariable('v', dims=(/'x ', 'time'/), & & longname='foo quantity', units='m/s') call HistoryPut('x', (/0.0, 1.0, 2.0/)) call HistorySetTime(0.0) call HistoryPut('u', (/0.0, 0.1, 0.2/)) call HistoryPut('v', (/1.0, 1.1, 1.2/)) call HistorySetTime(1.0) call HistoryPut('u', (/2.0, 2.1, 2.2/)) call HistorySetTime(2.0) call HistoryPut('u', (/3.0, 3.1, 3.2/)) call HistorySetTime(1.0) call HistoryPut('v', (/4.0, 4.1, 4.2/)) call HistorySetTime(2.0) call HistoryPut('v', (/5.0, 5.1, 5.2/)) call HistoryClose !----------------------------------------------------------------- ! * HistoryCreate1 を用いた出力テスト ! * 明示的に引数キーワードを与えずに出力 ! GT_HISTORY 変数を利用 ! ! * HistoryCopy によるコピーのテスト (GT_HISTORY 変数を利用した場合) !----------------------------------------------------------------- call HistoryCreate('xhisttest/xhisttest5.nc', & & 'gt4_history test 5', & & 'gt4f90io/Fortran library test kit', & & 'GFD Dennou Club', & & (/'x ', 'time'/), & & (/3, 0/), & & (/'eastward length', 'time '/), & & (/'m','s'/), & & 100.0, 2.5, & & (/'real', 'real'/), & & gt4hist_test_5, & & 'http://www.gfd-dennou.org/library/gtool4/conventions/', '4.1') call HistoryPut('x', (/100.0, 200.0, 300.0/), gt4hist_test_5) call HistoryAddVariable('u', dims=(/'x ', 'time'/), & & longname='any quantity', units='m/s', history=gt4hist_test_5) do, i = 1, 3 call HistoryPut('u', (/1.0*i, 2.0*i, 3.0*i/), gt4hist_test_5) enddo call HistoryCopy(gt4hist_test_7, 'xhisttest/xhisttest7.nc', & & gt4hist_test_5, & & title='gt4_history test 7', origin=0.0, interval=0.0) call HistoryAddVariable('u', dims=(/'x ', 'time'/), & & longname='any quantity', units='m/s', history=gt4hist_test_7) do, i = 1, 3 call HistoryPut('u', (/1.0*i, 2.0*i, 3.0*i/), gt4hist_test_7) enddo call HistoryClose(gt4hist_test_7) call HistoryClose(gt4hist_test_5) !----------------------------------------------------------------- ! * HistoryAddAttr を用いた出力テスト ! * HistoryCopy によるコピーのテスト (GT_HISTORY 変数を利用しない場合) !----------------------------------------------------------------- call HistoryCreate(file='xhisttest/xhisttest6.nc', & & title='gt4_history test 6', & & source='gt4f90io/Fortran library test kit', & & institution='GFD Dennou Club', & & dims=(/'lon','lat','t '/), & & dimsizes=(/4,8,0/), & & longnames=(/'longitude','latitude ','time '/), & & units=(/'deg.','deg.','sec.'/), & & origin=111.0, interval=12.3 ) call HistoryAddAttr('lon', 'topology', 'circular') call HistoryAddAttr('lon', 'modulo', 360.0) call HistoryPut( 'lon', (/0., 90., 180., 270./) ) call HistoryPut( 'lat', & & (/-70.0, -50.0, -30.0, -10.0, 10.0, 30.0, 50.0, 70.0/) ) call HistoryAddVariable( & ! 変数定義 & varname='zeta', dims=(/'lon','lat','t '/) , & & longname='vorticity', units='1/s', xtype='double' ) Zeta(:,:) = 1.0 call HistoryPut('zeta', Zeta) Zeta(:,:) = 2.0 call HistoryPut('zeta', Zeta) call HistoryCopy(gt4hist_test_8, 'xhisttest/xhisttest8.nc', & & title='gt4_history test 8', origin=13.13, interval=12.12) call HistoryAddVariable('zeta', dims=(/'lon ','lat ', 't '/), & & longname='any quantity', units='m/s', history=gt4hist_test_8) do, i = 1, 3 call HistoryPut('zeta', Zeta * i, gt4hist_test_8) enddo call HistoryClose(gt4hist_test_8) call HistoryClose !----------------------------------------------------------------- ! HistoryPut の range オプションの動作確認 !----------------------------------------------------------------- call HistoryCreate(file='xhisttest/xhisttest9.nc', & & title='gt4_history test 9', & & source='gt4f90io/Fortran library test kit', & & institution='GFD DennouClub', & & dims=StoA('x', 'l', 'time'), & & dimsizes=(/3, 5, 0/), & & longnames=StoA('eastward length', 'any parameter', 'time'), & & units=StoA('m','1','s'), & & origin=0.0, interval=0.0 ,& & xtypes=StoA('float', 'int', 'float')) call HistoryAddVariable('u', dims=StoA('x', 'l', 'time'), & & longname='foo quantity', units='m/s') call HistoryPut('x', (/10.0, 20.0, 30.0/)) call HistoryPut('l', (/1.0, 2.0, 3.0, 4.0, 5.0/)) call HistorySetTime(0.0) do i = 1, 5 char = CPrintf("l=%d", i=(/i/)) call HistoryPut('u', & & (/1.0 * real(i), 1.1 * real(i), 1.2 * real(i)/), range=char) end do call HistorySetTime(1.0) do i = 1, 5 char = CPrintf("l=%d", i=(/i/)) call HistoryPut('u', & & (/2.0 * real(i), 2.1 * real(i), 2.2 * real(i)/), range=char) end do call HistorySetTime(2.0) do i = 1, 5 char = CPrintf("l=%d", i=(/i/)) call HistoryPut('u', & & (/3.0 * real(i), 3.1 * real(i), 3.2 * real(i)/), range=char) end do call HistoryClose end program histtest