source: trunk/LMDZ.COMMON/libf/evolution/reshape_XIOS_output.F90 @ 2852

Last change on this file since 2852 was 2850, checked in by llange, 3 years ago

PEM
Specific commit for the reshape routines (I forgot to add it before commiting...)

Reshape routine adds one element in longitude as XIOS output have one element less compared to the longitude grid of the physics (point -180 and 180)

LL

File size: 8.3 KB
Line 
1program reshape_XIOS_output
2    use netcdf
3    implicit none
4    integer :: status, ncid, ncid1, ncid2
5    integer :: nDims, nVars, nGlobalAtts, unlimDimID
6    integer i,j
7
8    integer :: include_parents
9
10    integer, dimension(:),allocatable :: dimids
11    integer, dimension(:),allocatable :: varids
12
13    integer, dimension(:),allocatable :: dimids_2
14    integer, dimension(:),allocatable :: varids_2
15
16    integer, dimension(:),allocatable :: dimid_var
17
18    real, dimension(:), allocatable :: tempvalues_1d
19    real, dimension(:), allocatable :: values_1d
20
21    real, dimension(:,:), allocatable :: tempvalues_2d
22    real, dimension(:,:), allocatable :: values_2d
23
24    real, dimension(:,:,:), allocatable :: tempvalues_3d
25    real, dimension(:,:,:), allocatable :: values_3d
26
27    real, dimension(:,:,:,:), allocatable :: tempvalues_4d
28    real, dimension(:,:,:,:), allocatable :: values_4d
29
30  character*1 str2
31  character*30 :: name_
32  character*30 :: namevar
33  integer  :: xtype_var
34  integer :: len_
35  integer :: len_1,len_2
36  integer :: len_lat, len_lon, len_time, len_soil
37  integer :: dimid_lon, dimid_lat, dimid_time, dimid_soil
38  integer :: dimid_2
39  integer :: numdims
40  integer :: numatts
41  integer :: numyear
42
43
44
45
46
47DO numyear=1, 2
48write(*,*) 'numyear',numyear
49write(str2(1:1),'(i1.1)') numyear
50!nf90_open                 ! open existing netCDF dataset
51!integer :: ncid, status
52!...
53status = nf90_open(path = "Xdiurnalave"//str2//".nc", mode = nf90_nowrite, ncid = ncid1)
54if(status /= nf90_noerr) call handle_err(status)
55
56status = nf90_create(path = "Xdiurnalave"//str2//".nc_new", cmode=or(nf90_noclobber,nf90_64bit_offset), ncid = ncid2)
57if(status /= nf90_noerr) call handle_err(status)
58
59status = nf90_inquire(ncid1, ndims, nvars, nglobalatts, unlimdimid)
60if(status /= nf90_noerr) call handle_err(status)
61
62allocate(dimids(ndims))
63allocate(varids(nvars))
64
65allocate(dimids_2(ndims))
66allocate(varids_2(nvars))
67
68status = nf90_inq_dimids(ncid1, ndims, dimids, include_parents)
69if(status /= nf90_noerr) call handle_err(status)
70status = nf90_inq_varids(ncid1, nvars, varids)
71if(status /= nf90_noerr) call handle_err(status)
72
73do i=1,ndims
74  status = nf90_inquire_dimension(ncid1, dimids(i), name_, len_)
75  if(status /= nf90_noerr) call handle_err(status)
76  if(name_.eq."lon")  then
77     dimid_lon=dimids(i)
78     len_lon=len_
79     len_=len_+1
80  elseif(name_.eq."lat") then
81     dimid_lat=dimids(i)
82     len_lat=len_
83  elseif(name_.eq."time_counter") then
84     dimid_time=dimids(i)
85     len_time=len_
86  elseif(name_.eq."soil_layers") then
87     dimid_soil=dimids(i)
88     len_soil=len_
89  endif
90  status = nf90_def_dim(ncid2, name_, len_, dimid_2)
91  if(status /= nf90_noerr) call handle_err(status)
92  dimids_2(i)=dimid_2
93enddo
94
95
96
97allocate(tempvalues_3d(len_lon,len_lat,len_time))
98allocate(values_3d(len_lon+1,len_lat,len_time))
99
100allocate(tempvalues_4d(len_lon,len_lat,len_soil,len_time))
101allocate(values_4d(len_lon+1,len_lat,len_soil,len_time))
102
103
104do i=1,nvars
105  status = nf90_inquire_variable(ncid1, varids(i), name=namevar, xtype=xtype_var, ndims = numdims,natts = numatts)
106      print *, "namevar00= ", namevar
107  if(status /= nf90_noerr) call handle_err(status)
108  allocate(dimid_var(numdims))
109  status = nf90_inquire_variable(ncid1, varids(i), name=namevar, xtype=xtype_var, ndims = numdims, dimids=dimid_var, natts = numatts)
110  if(status /= nf90_noerr) call handle_err(status)
111  if(numdims.eq.1) then
112    if(namevar.eq."lon") then
113      allocate(tempvalues_1d(len_lon))
114      allocate(values_1d(len_lon+1))
115      status = nf90_get_var(ncid1, varids(i), tempvalues_1d)
116      if(status /= nf90_noerr) call handle_err(status)
117      status = nf90_def_var(ncid2, namevar, xtype_var, dimid_var, varids_2(i))
118      if(status /= nf90_noerr) call handle_err(status)
119      values_1d(1:len_lon)=tempvalues_1d(:)
120      values_1d(len_lon+1)=values_1d(1)
121      status = nf90_enddef(ncid2)
122      if(status /= nf90_noerr) call handle_err(status)
123      status = nf90_put_var(ncid2, varids_2(i), values_1d)
124      if(status /= nf90_noerr) call handle_err(status)
125      status = nf90_redef(ncid2)
126      if(status /= nf90_noerr) call handle_err(status)
127      deallocate(tempvalues_1d)
128      deallocate(values_1d) 
129    else
130      status = nf90_inquire_dimension(ncid1, dimid_var(1), name_, len_)
131      if(status /= nf90_noerr) call handle_err(status)
132      allocate(tempvalues_1d(len_))
133      status = nf90_get_var(ncid1, varids(i), tempvalues_1d)
134      if(status /= nf90_noerr) call handle_err(status)
135      status = nf90_def_var(ncid2, namevar, xtype_var, dimid_var, varids_2(i))
136      if(status /= nf90_noerr) call handle_err(status)
137      status = nf90_enddef(ncid2)
138      if(status /= nf90_noerr) call handle_err(status)
139      status = nf90_put_var(ncid2, varids_2(i), tempvalues_1d) 
140      if(status /= nf90_noerr) call handle_err(status)
141      status = nf90_redef(ncid2)
142      if(status /= nf90_noerr) call handle_err(status)
143      deallocate(tempvalues_1d)   
144    endif
145  elseif(numdims.eq.2) then
146    if(namevar.eq."area") then
147      allocate(tempvalues_2d(len_lon,len_lat))
148      allocate(values_2d(len_lon+1,len_lat))     
149      status = nf90_get_var(ncid1, varids(i), tempvalues_2d)
150      if(status /= nf90_noerr) call handle_err(status)
151      status = nf90_def_var(ncid2, namevar, xtype_var, dimid_var, varids_2(i))
152      if(status /= nf90_noerr) call handle_err(status)
153      values_2d(1:len_lon,:)=tempvalues_2d(:,:)
154      values_2d(len_lon+1,:)=values_2d(1,:)
155      status = nf90_enddef(ncid2)
156      if(status /= nf90_noerr) call handle_err(status)
157      status = nf90_put_var(ncid2, varids_2(i), values_2d)   
158      if(status /= nf90_noerr) call handle_err(status)
159      status = nf90_redef(ncid2)
160      if(status /= nf90_noerr) call handle_err(status)
161      deallocate(tempvalues_2d)
162      deallocate(values_2d)
163    else
164      status = nf90_inquire_dimension(ncid1, dimid_var(1), name_, len_1)
165      if(status /= nf90_noerr) call handle_err(status)
166      status = nf90_inquire_dimension(ncid1, dimid_var(2), name_, len_2)
167      if(status /= nf90_noerr) call handle_err(status)
168      allocate(tempvalues_2d(len_1,len_2))
169      status = nf90_get_var(ncid1, varids(i), tempvalues_2d)
170      if(status /= nf90_noerr) call handle_err(status)
171      status = nf90_def_var(ncid2, namevar, xtype_var, dimid_var, varids_2(i))
172      if(status /= nf90_noerr) call handle_err(status)
173      status = nf90_enddef(ncid2)
174      if(status /= nf90_noerr) call handle_err(status)
175      status = nf90_put_var(ncid2, varids_2(i), tempvalues_2d)
176      if(status /= nf90_noerr) call handle_err(status)
177      status = nf90_redef(ncid2)
178      if(status /= nf90_noerr) call handle_err(status)
179      deallocate(tempvalues_2d)
180    endif
181  elseif(numdims.eq.3) then
182      status = nf90_get_var(ncid1, varids(i), tempvalues_3d)
183      if(status /= nf90_noerr) call handle_err(status)
184      status = nf90_def_var(ncid2, namevar, xtype_var, dimid_var, varids_2(i))
185      if(status /= nf90_noerr) call handle_err(status)
186      values_3d(1:len_lon,:,:)=tempvalues_3d(:,:,:)
187      values_3d(len_lon+1,:,:)=values_3d(1,:,:)
188      status = nf90_enddef(ncid2)
189      if(status /= nf90_noerr) call handle_err(status)
190      status = nf90_put_var(ncid2, varids_2(i), values_3d)
191      if(status /= nf90_noerr) call handle_err(status)
192      status = nf90_redef(ncid2)
193      if(status /= nf90_noerr) call handle_err(status)
194  elseif(numdims.eq.4) then
195      status = nf90_get_var(ncid1, varids(i), tempvalues_4d)
196      if(status /= nf90_noerr) call handle_err(status)
197      status = nf90_def_var(ncid2, namevar, xtype_var, dimid_var, varids_2(i))
198      if(status /= nf90_noerr) call handle_err(status)
199      status = nf90_enddef(ncid2)
200      values_4d(1:len_lon,:,:,:)=tempvalues_4d(:,:,:,:)
201      values_4d(len_lon+1,:,:,:)=values_4d(1,:,:,:)
202      if(status /= nf90_noerr) call handle_err(status)
203      status = nf90_put_var(ncid2, varids_2(i), values_4d)
204      if(status /= nf90_noerr) call handle_err(status)
205      status = nf90_redef(ncid2)
206      if(status /= nf90_noerr) call handle_err(status)
207  endif
208
209  deallocate(dimid_var)
210enddo
211
212status = nf90_enddef(ncid2)
213if(status /= nf90_noerr) call handle_err(status)
214status = nf90_close(ncid1)
215if(status /= nf90_noerr) call handle_err(status)
216status = nf90_close(ncid2)
217if(status /= nf90_noerr) call handle_err(status)
218
219
220deallocate(dimids)
221deallocate(varids)
222deallocate(dimids_2)
223deallocate(varids_2)
224deallocate(tempvalues_3d)
225deallocate(values_3d)
226deallocate(tempvalues_4d)
227deallocate(values_4d)
228
229
230
231enddo
232
233end program reshape_XIOS_output
234
Note: See TracBrowser for help on using the repository browser.