Path: | anvargetnum.f90 |
Last Update: | Mon Mar 20 20:58:38 JST 2006 |
Authors: | Yasuhiro MORIKAWA, Eizi TOYODA |
Version: | $Id: anvargetnum.f90,v 1.4 2006/03/20 11:58:38 morikawa Exp $ |
Tag Name: | $Name: gt4f90io-20060818 $ |
Copyright: | Copyright (C) GFD Dennou Club, 2000-2005. All rights reserved. |
License: | See COPYRIGHT |
以下のサブルーチン、関数は an_generic から gtdata_generic#Get として提供されます。
Subroutine : | |
var : | type(AN_VARIABLE), intent(in) |
start(:) : | integer, intent(in) |
cnt(:) : | integer, intent(in) |
stride(:) : | integer, intent(in) |
imap(:) : | integer, intent(in) |
siz : | integer, intent(in) |
value(siz) : | real(DP), intent(out) |
iostat : | integer, intent(out) |
subroutine ANVarGetDouble(var, start, cnt, stride, imap, siz, value, iostat) use an_types, only: AN_VARIABLE use an_vartable, only: AN_VARIABLE_ENTRY, vtable_lookup use netcdf_f77, only: nf_noerr, nf_einval, nf_get_varm_Double, nf_get_var1_Double use dc_types, only: DP use dc_trace, only: BeginSub, EndSub, DbgMessage implicit none type(AN_VARIABLE), intent(in):: var integer, intent(in):: start(:) integer, intent(in):: cnt(:) integer, intent(in):: stride(:) integer, intent(in):: imap(:) integer, intent(in):: siz real(DP), intent(out):: value(siz) integer, intent(out):: iostat integer:: nd, ipos, i type(AN_VARIABLE_ENTRY):: ent integer, allocatable:: istart(:), istride(:), iimap(:) continue call BeginSub('ANVarGetDouble', fmt='varmap=%d, start=%*d, cnt=%*d, stride=%*d, imap=%*d siz=%d', i=(/var%id, start(:), cnt(:), stride(:), imap(:), siz/), n=(/size(start), size(cnt), size(stride), size(imap)/)) iostat = vtable_lookup(var, ent) if (iostat /= nf_noerr) goto 999 ! --- nd check --- nd = 0 if (associated(ent%dimids)) nd = size(ent%dimids) if (min(size(start), size(cnt), size(stride), size(imap)) < nd) then iostat = nf_einval goto 999 endif if (nd == 0) then iostat = nf_get_var1_Double(ent%fileid, ent%varid, start, value(1)) goto 999 endif ! --- stride kakikae buffer --- allocate(istart(nd), istride(nd), iimap(nd)) istart(1:nd) = start(1:nd) istride(1:nd) = stride(1:nd) iimap(1:nd) = imap(1:nd) ipos = 1 ! --- do read --- if (ent%varid <= 0 .or. count(cnt(1:nd) == 1) >= 0) then call BeginSub('fake_map_get') call fake_map_get call EndSub('fake_map_get', 'iostat=%d', i=(/iostat/)) else ! negative stride is not allowed for netcdf do, i = 1, nd if (stride(i) > 0) cycle ipos = ipos + (cnt(i) - 1) * imap(i) istart(i) = start(i) + (cnt(i) - 1) * stride(i) istride(i) = -stride(i) iimap(i) = -imap(i) call DbgMessage('dim %d negate: stride->%d start->%d map->%d', i=(/i, istride(i), istart(i), iimap(i)/)) enddo iostat = nf_get_varm_Double(ent%fileid, ent%varid, istart, cnt, istride, iimap, value(ipos)) endif deallocate(istart, istride, iimap) 999 continue call EndSub('ANVarGetDouble', 'iostat=%d', i=(/iostat/)) return contains subroutine fake_map_get use netcdf_f77, only: NF_EINDEFINE, nf_enddef, nf_redef, nf_get_var1_Double implicit none integer:: ofs(nd), here(nd) integer:: j, result continue iostat = nf_noerr ofs(1:nd) = 0 do j = ipos + dot_product(ofs(1:nd), imap(1:nd)) here(1:nd) = istart(1:nd) + ofs(1:nd) * istride(1:nd) if (j < 1 .or. j > siz) then iostat = nf_einval call DbgMessage('nf_get_var1_Double(ncid=%d, varid=%d,& i=(/ent%fileid, ent%varid, here(1:nd), j/), n=(/nd/)) return endif if (ent%varid == 0) then value(j) = j iostat = nf_noerr else iostat = nf_get_var1_Double(ent%fileid, ent%varid, here(1), value(j)) if (iostat == NF_EINDEFINE) then result = nf_enddef(ent%fileid) iostat = nf_get_var1_Double(ent%fileid, ent%varid, here(1), value(j)) result = nf_redef(ent%fileid) end if endif if (iostat /= nf_noerr) return ofs(1) = ofs(1) + 1 do, j = 1, nd - 1 if (ofs(j) < cnt(j)) exit ofs(j) = 0 ofs(j + 1) = ofs(j + 1) + 1 enddo if (ofs(nd) >= cnt(nd)) exit enddo end subroutine fake_map_get end subroutine ANVarGetDouble
Subroutine : | |
var : | type(AN_VARIABLE), intent(in) |
start(:) : | integer, intent(in) |
cnt(:) : | integer, intent(in) |
stride(:) : | integer, intent(in) |
imap(:) : | integer, intent(in) |
siz : | integer, intent(in) |
value(siz) : | real, intent(out) |
iostat : | integer, intent(out) |
subroutine ANVarGetReal(var, start, cnt, stride, imap, siz, value, iostat) use an_types, only: AN_VARIABLE use an_vartable, only: AN_VARIABLE_ENTRY, vtable_lookup use netcdf_f77, only: nf_noerr, nf_einval, nf_get_varm_Real, nf_get_var1_Real use dc_types, only: DP use dc_trace, only: BeginSub, EndSub, DbgMessage implicit none type(AN_VARIABLE), intent(in):: var integer, intent(in):: start(:) integer, intent(in):: cnt(:) integer, intent(in):: stride(:) integer, intent(in):: imap(:) integer, intent(in):: siz real, intent(out):: value(siz) integer, intent(out):: iostat integer:: nd, ipos, i type(AN_VARIABLE_ENTRY):: ent integer, allocatable:: istart(:), istride(:), iimap(:) continue call BeginSub('ANVarGetReal', fmt='varmap=%d, start=%*d, cnt=%*d, stride=%*d, imap=%*d siz=%d', i=(/var%id, start(:), cnt(:), stride(:), imap(:), siz/), n=(/size(start), size(cnt), size(stride), size(imap)/)) iostat = vtable_lookup(var, ent) if (iostat /= nf_noerr) goto 999 ! --- nd check --- nd = 0 if (associated(ent%dimids)) nd = size(ent%dimids) if (min(size(start), size(cnt), size(stride), size(imap)) < nd) then iostat = nf_einval goto 999 endif if (nd == 0) then iostat = nf_get_var1_Real(ent%fileid, ent%varid, start, value(1)) goto 999 endif ! --- stride kakikae buffer --- allocate(istart(nd), istride(nd), iimap(nd)) istart(1:nd) = start(1:nd) istride(1:nd) = stride(1:nd) iimap(1:nd) = imap(1:nd) ipos = 1 ! --- do read --- if (ent%varid <= 0 .or. count(cnt(1:nd) == 1) >= 0) then call BeginSub('fake_map_get') call fake_map_get call EndSub('fake_map_get', 'iostat=%d', i=(/iostat/)) else ! negative stride is not allowed for netcdf do, i = 1, nd if (stride(i) > 0) cycle ipos = ipos + (cnt(i) - 1) * imap(i) istart(i) = start(i) + (cnt(i) - 1) * stride(i) istride(i) = -stride(i) iimap(i) = -imap(i) call DbgMessage('dim %d negate: stride->%d start->%d map->%d', i=(/i, istride(i), istart(i), iimap(i)/)) enddo iostat = nf_get_varm_Real(ent%fileid, ent%varid, istart, cnt, istride, iimap, value(ipos)) endif deallocate(istart, istride, iimap) 999 continue call EndSub('ANVarGetReal', 'iostat=%d', i=(/iostat/)) return contains subroutine fake_map_get use netcdf_f77, only: NF_EINDEFINE, nf_enddef, nf_redef, nf_get_var1_Real implicit none integer:: ofs(nd), here(nd) integer:: j, result continue iostat = nf_noerr ofs(1:nd) = 0 do j = ipos + dot_product(ofs(1:nd), imap(1:nd)) here(1:nd) = istart(1:nd) + ofs(1:nd) * istride(1:nd) if (j < 1 .or. j > siz) then iostat = nf_einval call DbgMessage('nf_get_var1_Real(ncid=%d, varid=%d,& i=(/ent%fileid, ent%varid, here(1:nd), j/), n=(/nd/)) return endif if (ent%varid == 0) then value(j) = j iostat = nf_noerr else iostat = nf_get_var1_Real(ent%fileid, ent%varid, here(1), value(j)) if (iostat == NF_EINDEFINE) then result = nf_enddef(ent%fileid) iostat = nf_get_var1_Real(ent%fileid, ent%varid, here(1), value(j)) result = nf_redef(ent%fileid) end if endif if (iostat /= nf_noerr) return ofs(1) = ofs(1) + 1 do, j = 1, nd - 1 if (ofs(j) < cnt(j)) exit ofs(j) = 0 ofs(j + 1) = ofs(j + 1) + 1 enddo if (ofs(nd) >= cnt(nd)) exit enddo end subroutine fake_map_get end subroutine ANVarGetReal
Subroutine : |
subroutine fake_map_get use netcdf_f77, only: NF_EINDEFINE, nf_enddef, nf_redef, nf_get_var1_Real implicit none integer:: ofs(nd), here(nd) integer:: j, result continue iostat = nf_noerr ofs(1:nd) = 0 do j = ipos + dot_product(ofs(1:nd), imap(1:nd)) here(1:nd) = istart(1:nd) + ofs(1:nd) * istride(1:nd) if (j < 1 .or. j > siz) then iostat = nf_einval call DbgMessage('nf_get_var1_Real(ncid=%d, varid=%d,& i=(/ent%fileid, ent%varid, here(1:nd), j/), n=(/nd/)) return endif if (ent%varid == 0) then value(j) = j iostat = nf_noerr else iostat = nf_get_var1_Real(ent%fileid, ent%varid, here(1), value(j)) if (iostat == NF_EINDEFINE) then result = nf_enddef(ent%fileid) iostat = nf_get_var1_Real(ent%fileid, ent%varid, here(1), value(j)) result = nf_redef(ent%fileid) end if endif if (iostat /= nf_noerr) return ofs(1) = ofs(1) + 1 do, j = 1, nd - 1 if (ofs(j) < cnt(j)) exit ofs(j) = 0 ofs(j + 1) = ofs(j + 1) + 1 enddo if (ofs(nd) >= cnt(nd)) exit enddo end subroutine fake_map_get
Subroutine : |
subroutine fake_map_get use netcdf_f77, only: NF_EINDEFINE, nf_enddef, nf_redef, nf_get_var1_Double implicit none integer:: ofs(nd), here(nd) integer:: j, result continue iostat = nf_noerr ofs(1:nd) = 0 do j = ipos + dot_product(ofs(1:nd), imap(1:nd)) here(1:nd) = istart(1:nd) + ofs(1:nd) * istride(1:nd) if (j < 1 .or. j > siz) then iostat = nf_einval call DbgMessage('nf_get_var1_Double(ncid=%d, varid=%d,& i=(/ent%fileid, ent%varid, here(1:nd), j/), n=(/nd/)) return endif if (ent%varid == 0) then value(j) = j iostat = nf_noerr else iostat = nf_get_var1_Double(ent%fileid, ent%varid, here(1), value(j)) if (iostat == NF_EINDEFINE) then result = nf_enddef(ent%fileid) iostat = nf_get_var1_Double(ent%fileid, ent%varid, here(1), value(j)) result = nf_redef(ent%fileid) end if endif if (iostat /= nf_noerr) return ofs(1) = ofs(1) + 1 do, j = 1, nd - 1 if (ofs(j) < cnt(j)) exit ofs(j) = 0 ofs(j + 1) = ofs(j + 1) + 1 enddo if (ofs(nd) >= cnt(nd)) exit enddo end subroutine fake_map_get