[629] | 1 | module iophy |
---|
| 2 | |
---|
| 3 | REAL,private,allocatable,dimension(:,:),save :: tmp_tab2d |
---|
| 4 | REAL,private,allocatable,dimension(:,:,:),save :: tmp_tab3d |
---|
| 5 | INTEGER,private,allocatable,dimension(:),save :: ndex2d |
---|
| 6 | INTEGER,private,allocatable,dimension(:),save :: ndex3d |
---|
| 7 | REAL,private,allocatable,dimension(:),save :: io_lat |
---|
| 8 | REAL,private,allocatable,dimension(:),save :: io_lon |
---|
| 9 | |
---|
| 10 | INTERFACE histwrite_phy |
---|
| 11 | MODULE PROCEDURE histwrite2d_phy,histwrite3d_phy |
---|
| 12 | END INTERFACE |
---|
| 13 | |
---|
| 14 | contains |
---|
| 15 | |
---|
| 16 | subroutine init_iophy(lat,lon) |
---|
| 17 | use dimphy |
---|
| 18 | implicit none |
---|
| 19 | include 'dimensions90.h' |
---|
| 20 | real,dimension(iim),intent(in) :: lon |
---|
| 21 | real,dimension(jjm+1),intent(in) :: lat |
---|
| 22 | |
---|
| 23 | allocate(io_lat(jjm+1)) |
---|
| 24 | io_lat(:)=lat(:) |
---|
| 25 | allocate(io_lon(iim)) |
---|
| 26 | io_lon(:)=lon(:) |
---|
| 27 | allocate(tmp_tab2d(iim,jjphy_nb)) |
---|
| 28 | allocate(tmp_tab3d(iim,jjphy_nb,klev)) |
---|
| 29 | allocate(ndex2d(iim*jjphy_nb)) |
---|
| 30 | allocate(ndex3d(iim*jjphy_nb*klev)) |
---|
| 31 | ndex2d(:)=0 |
---|
| 32 | ndex3d(:)=0 |
---|
| 33 | |
---|
| 34 | end subroutine init_iophy |
---|
| 35 | |
---|
| 36 | subroutine histbeg_phy(name,itau0,zjulian,dtime,nhori,nid_day) |
---|
| 37 | use dimphy |
---|
| 38 | use ioipsl |
---|
| 39 | use write_field |
---|
| 40 | implicit none |
---|
| 41 | include 'dimensions90.h' |
---|
| 42 | |
---|
| 43 | character*(*), intent(IN) :: name |
---|
| 44 | integer, intent(in) :: itau0 |
---|
| 45 | real,intent(in) :: zjulian |
---|
| 46 | real,intent(in) :: dtime |
---|
| 47 | integer,intent(out) :: nhori |
---|
| 48 | integer,intent(out) :: nid_day |
---|
| 49 | |
---|
| 50 | if (monocpu) then |
---|
| 51 | call histbeg(name,iim,io_lon, jjphy_nb,io_lat(jjphy_begin:jjphy_end), & |
---|
| 52 | 1,iim,1,jjphy_nb,itau0, zjulian, dtime, nhori, nid_day) |
---|
| 53 | else |
---|
| 54 | call histbeg(name//'_'//trim(int2str(phy_rank)),iim,io_lon, jjphy_nb,io_lat(jjphy_begin:jjphy_end), & |
---|
| 55 | 1,iim,1,jjphy_nb,itau0, zjulian, dtime, nhori, nid_day) |
---|
| 56 | endif |
---|
| 57 | |
---|
| 58 | end subroutine histbeg_phy |
---|
| 59 | |
---|
| 60 | subroutine histwrite2d_phy(nid,name,itau,field) |
---|
| 61 | use dimphy |
---|
| 62 | use ioipsl |
---|
| 63 | implicit none |
---|
| 64 | include 'dimensions90.h' |
---|
| 65 | |
---|
| 66 | integer,intent(in) :: nid |
---|
| 67 | character*(*), intent(IN) :: name |
---|
| 68 | integer, intent(in) :: itau |
---|
| 69 | real,dimension(klon),intent(in) :: field |
---|
| 70 | |
---|
| 71 | CALL phy2dyn(field,tmp_tab2d,1) |
---|
| 72 | CALL histwrite(nid,name,itau,tmp_tab2d,iim*jjphy_nb,ndex2d) |
---|
| 73 | |
---|
| 74 | end subroutine histwrite2d_phy |
---|
| 75 | |
---|
| 76 | subroutine histwrite3d_phy(nid,name,itau,field) |
---|
| 77 | use dimphy |
---|
| 78 | use ioipsl |
---|
| 79 | implicit none |
---|
| 80 | include 'dimensions90.h' |
---|
| 81 | |
---|
| 82 | integer,intent(in) :: nid |
---|
| 83 | character*(*), intent(IN) :: name |
---|
| 84 | integer, intent(in) :: itau |
---|
| 85 | real,dimension(klon,klev),intent(in) :: field |
---|
| 86 | |
---|
| 87 | CALL phy2dyn(field,tmp_tab3d,klev) |
---|
| 88 | CALL histwrite(nid,name,itau,tmp_tab3d,iim*jjphy_nb*klev,ndex3d) |
---|
| 89 | |
---|
| 90 | end subroutine histwrite3d_phy |
---|
| 91 | |
---|
| 92 | |
---|
| 93 | subroutine phy2dyn(field_phy,field_dyn,nlev) |
---|
| 94 | use dimphy |
---|
| 95 | implicit none |
---|
| 96 | include 'dimensions90.h' |
---|
| 97 | |
---|
| 98 | real,dimension(klon,nlev),intent(in) :: field_phy |
---|
| 99 | real,dimension(iim,jjphy_nb,nlev),intent(out) :: field_dyn |
---|
| 100 | integer,intent(in) :: nlev |
---|
| 101 | |
---|
| 102 | integer :: next |
---|
| 103 | integer :: j,l |
---|
| 104 | |
---|
| 105 | do l=1,nlev |
---|
| 106 | |
---|
| 107 | if (jjphy_begin==jjphy_end) then |
---|
| 108 | field_dyn(:,1,l)=0. |
---|
| 109 | field_dyn(iiphy_begin:iiphy_end,1,l)=field_phy(1:klon,l) |
---|
| 110 | else |
---|
| 111 | |
---|
| 112 | if (jjphy_begin==1) then |
---|
| 113 | field_dyn(:,1,l)=field_phy(1,l) |
---|
| 114 | next=2 |
---|
| 115 | else |
---|
| 116 | field_dyn(:,1,l)=0. |
---|
| 117 | next=iim-iiphy_begin+1 |
---|
| 118 | field_dyn(iiphy_begin:iim,1,l)=field_phy(1:next-1,l) |
---|
| 119 | endif |
---|
| 120 | |
---|
| 121 | do j=2,jjphy_nb-1 |
---|
| 122 | field_dyn(:,j,l)=field_phy(next:next+iim-1,l) |
---|
| 123 | next=next+iim |
---|
| 124 | enddo |
---|
| 125 | |
---|
| 126 | if (jjphy_end==jjm+1) then |
---|
| 127 | field_dyn(:,jjphy_nb,l)=field_phy(klon,l) |
---|
| 128 | else |
---|
| 129 | field_dyn(:,jjphy_nb,l)=0. |
---|
| 130 | field_dyn(1:iiphy_end,jjphy_nb,l)=field_phy(next:next+iiphy_end-1,l) |
---|
| 131 | endif |
---|
| 132 | |
---|
| 133 | endif |
---|
| 134 | |
---|
| 135 | enddo |
---|
| 136 | |
---|
| 137 | end subroutine phy2dyn |
---|
| 138 | |
---|
| 139 | |
---|
| 140 | end module iophy |
---|