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

Last change on this file since 3532 was 3214, checked in by jbclement, 10 months ago

PEM:

  • It is now possible to set the number of initial PCM calls independently of the number of "inter-PEM" PCM calls. It is useful to get a stable situation for the start of the simulation.
  • Correction of a bug: 'reshape_XIOS_output' treated the first two PCM runs instead of the last two. The numbering has been adapted in "launch_pem.sh" to get it right.
  • PEM outputs ("diagpem.nc" and "diagsoil_pem.nc") has been improved: 'Time' now evolves according to 'dt_pem' (usually year by year) and 'ecritpem' is a full variable of "time_evol_mod.F90".

JBC

File size: 9.9 KB
Line 
1PROGRAM reshape_XIOS_output
2
3!=======================================================================
4!
5! Purpose: Read XIOS files, and convert them into the correct PCM grid
6!          XIOS  longitudes start at -180 but stop before -180 (not duplicated)
7!          We basically add the last point, and complete the XIOS file. Looped
8!          over the two PCM runs
9!
10! Authors: RV & LL
11!=======================================================================
12use netcdf
13
14implicit none
15
16integer                               :: state, ncid1, ncid2, nDims, nVars, nGlobalAtts, unlimDimID
17integer                               :: i, include_parents, cstat
18integer, dimension(:),    allocatable :: dimids, varids, dimids_2, varids_2, dimid_var
19real, dimension(:),       allocatable :: tempvalues_1d, values_1d
20real, dimension(:,:),     allocatable :: tempvalues_2d, values_2d
21real, dimension(:,:,:),   allocatable :: tempvalues_3d, values_3d
22real, dimension(:,:,:,:), allocatable :: tempvalues_4d, values_4d
23character(1)                          :: str
24character(30)                         :: name_, namevar
25integer                               :: xtype_var, len_, len_1, len_2, len_lat, len_lon, len_time, len_soil
26integer                               :: dimid_lon, dimid_lat, dimid_time, dimid_soil, dimid_2, numdims, numatts, numyear
27logical                               :: yes
28
29do numyear = 1,2
30    write(*,*) 'numyear',numyear
31    write(str(1:1),'(i1.1)') numyear
32
33    state = nf90_open(path = "data2reshape_Y"//str//".nc", mode = nf90_nowrite, ncid = ncid1)
34    if (state /= nf90_noerr) call handle_err(state)
35
36    inquire(file = 'data_PCM_Y'//str//'.nc', exist = yes)
37    if (yes) then
38        call execute_command_line('rm data_PCM_Y'//str//'.nc',cmdstat = cstat)
39        if (cstat > 0) then
40            error stop 'Command exection failed!'
41        else if (cstat < 0) then
42            error stop 'Command execution not supported!'
43        endif
44    endif
45    state = nf90_create(path = "data_PCM_Y"//str//".nc", cmode=or(nf90_noclobber,nf90_64bit_offset), ncid = ncid2)
46    if (state /= nf90_noerr) call handle_err(state)
47
48    state = nf90_inquire(ncid1, ndims, nvars, nglobalatts, unlimdimid)
49    if (state /= nf90_noerr) call handle_err(state)
50
51    allocate(dimids(ndims))
52    allocate(varids(nvars))
53
54    allocate(dimids_2(ndims))
55    allocate(varids_2(nvars))
56
57    state = nf90_inq_dimids(ncid1, ndims, dimids, include_parents)
58    if (state /= nf90_noerr) call handle_err(state)
59    state = nf90_inq_varids(ncid1, nvars, varids)
60    if (state /= nf90_noerr) call handle_err(state)
61
62    do i = 1,ndims
63        state = nf90_inquire_dimension(ncid1, dimids(i), name_, len_)
64        if (state /= nf90_noerr) call handle_err(state)
65        if (name_ == "lon" .or. name_ == "longitude") then
66            dimid_lon = dimids(i)
67            len_lon = len_
68            len_ = len_ + 1
69        else if (name_ == "lat".or. name_ == "latitude") then
70            dimid_lat=dimids(i)
71            len_lat=len_
72        else if (name_ == "time_counter".or. name_ ==  "Time") then
73            dimid_time=dimids(i)
74            len_time=len_
75        else if (name_ == "soil_layers".or. name_ ==  "subsurface_layers") then
76            dimid_soil=dimids(i)
77            len_soil = len_
78        endif
79        state = nf90_def_dim(ncid2, name_,len_,dimid_2)
80        if (state /= nf90_noerr) call handle_err(state)
81        dimids_2(i) = dimid_2
82    enddo
83
84    do i = 1,nvars
85        state = nf90_inquire_variable(ncid1, varids(i),name = namevar,xtype = xtype_var,ndims = numdims,natts = numatts)
86        write(*,*) "namevar00= ", namevar
87        if (state /= nf90_noerr) call handle_err(state)
88        allocate(dimid_var(numdims))
89        state = nf90_inquire_variable(ncid1,varids(i),name = namevar,xtype = xtype_var,ndims = numdims,dimids = dimid_var,natts = numatts)
90        if (state /= nf90_noerr) call handle_err(state)
91        if (numdims == 1) then
92            if (namevar == "lon") then
93                allocate(tempvalues_1d(len_lon))
94                allocate(values_1d(len_lon + 1))
95                state = nf90_get_var(ncid1,varids(i),tempvalues_1d)
96                if (state /= nf90_noerr) call handle_err(state)
97                state = nf90_def_var(ncid2,namevar,xtype_var, dimid_var, varids_2(i))
98                if (state /= nf90_noerr) call handle_err(state)
99                values_1d(1:len_lon) = tempvalues_1d(:)
100                values_1d(len_lon + 1) = values_1d(1)
101                state = nf90_enddef(ncid2)
102                if (state /= nf90_noerr) call handle_err(state)
103                state = nf90_put_var(ncid2, varids_2(i), values_1d)
104                if (state /= nf90_noerr) call handle_err(state)
105                state = nf90_redef(ncid2)
106                if (state /= nf90_noerr) call handle_err(state)
107                deallocate(tempvalues_1d)
108                deallocate(values_1d)
109            else
110                state = nf90_inquire_dimension(ncid1,dimid_var(1),name_,len_)
111                if (state /= nf90_noerr) call handle_err(state)
112                allocate(tempvalues_1d(len_))
113                state = nf90_get_var(ncid1,varids(i),tempvalues_1d)
114                if (state /= nf90_noerr) call handle_err(state)
115                state = nf90_def_var(ncid2,namevar,xtype_var,dimid_var,varids_2(i))
116                if (state /= nf90_noerr) call handle_err(state)
117                state = nf90_enddef(ncid2)
118                if (state /= nf90_noerr) call handle_err(state)
119                state = nf90_put_var(ncid2, varids_2(i), tempvalues_1d)
120                if (state /= nf90_noerr) call handle_err(state)
121                state = nf90_redef(ncid2)
122                if (state /= nf90_noerr) call handle_err(state)
123                deallocate(tempvalues_1d)
124            endif
125        else if (numdims == 2) then
126            if (namevar == "area") then
127                allocate(tempvalues_2d(len_lon,len_lat))
128                allocate(values_2d(len_lon + 1,len_lat))
129                state = nf90_get_var(ncid1,varids(i),tempvalues_2d)
130                if (state /= nf90_noerr) call handle_err(state)
131                state = nf90_def_var(ncid2,namevar,xtype_var,dimid_var,varids_2(i))
132                if (state /= nf90_noerr) call handle_err(state)
133                values_2d(1:len_lon,:) = tempvalues_2d(:,:)
134                values_2d(len_lon+1,:) = values_2d(1,:)
135                state = nf90_enddef(ncid2)
136                if (state /= nf90_noerr) call handle_err(state)
137                state = nf90_put_var(ncid2,varids_2(i),values_2d)
138                if (state /= nf90_noerr) call handle_err(state)
139                state = nf90_redef(ncid2)
140                if (state /= nf90_noerr) call handle_err(state)
141                deallocate(tempvalues_2d)
142                deallocate(values_2d)
143            else
144                state = nf90_inquire_dimension(ncid1,dimid_var(1),name_,len_1)
145                if (state /= nf90_noerr) call handle_err(state)
146                state = nf90_inquire_dimension(ncid1,dimid_var(2),name_,len_2)
147                if (state /= nf90_noerr) call handle_err(state)
148                allocate(tempvalues_2d(len_1,len_2))
149                state = nf90_get_var(ncid1, varids(i), tempvalues_2d)
150                if (state /= nf90_noerr) call handle_err(state)
151                state = nf90_def_var(ncid2,namevar,xtype_var,dimid_var,varids_2(i))
152                if (state /= nf90_noerr) call handle_err(state)
153                state = nf90_enddef(ncid2)
154                if (state /= nf90_noerr) call handle_err(state)
155                state = nf90_put_var(ncid2, varids_2(i), tempvalues_2d)
156                if (state /= nf90_noerr) call handle_err(state)
157                state = nf90_redef(ncid2)
158                if (state /= nf90_noerr) call handle_err(state)
159                deallocate(tempvalues_2d)
160            endif
161        else if (numdims == 3) then
162            allocate(tempvalues_3d(len_lon,len_lat,len_time))
163            allocate(values_3d(len_lon + 1,len_lat,len_time))
164            state = nf90_get_var(ncid1,varids(i),tempvalues_3d)
165            if (state /= nf90_noerr) call handle_err(state)
166            state = nf90_def_var(ncid2,namevar,xtype_var,dimid_var,varids_2(i))
167            if (state /= nf90_noerr) call handle_err(state)
168            values_3d(1:len_lon,:,:) = tempvalues_3d(:,:,:)
169            values_3d(len_lon+1,:,:) = values_3d(1,:,:)
170            state = nf90_enddef(ncid2)
171            if (state /= nf90_noerr) call handle_err(state)
172            state = nf90_put_var(ncid2, varids_2(i), values_3d)
173            if (state /= nf90_noerr) call handle_err(state)
174            state = nf90_redef(ncid2)
175            if (state /= nf90_noerr) call handle_err(state)
176            deallocate(tempvalues_3d)
177            deallocate(values_3d)
178        else if (numdims == 4) then
179            allocate(tempvalues_4d(len_lon,len_lat,len_soil,len_time))
180            allocate(values_4d(len_lon+1,len_lat,len_soil,len_time))
181            state = nf90_get_var(ncid1, varids(i), tempvalues_4d)
182            if (state /= nf90_noerr) call handle_err(state)
183            state = nf90_def_var(ncid2, namevar, xtype_var, dimid_var, varids_2(i))
184            if (state /= nf90_noerr) call handle_err(state)
185            state = nf90_enddef(ncid2)
186            values_4d(1:len_lon,:,:,:) = tempvalues_4d(:,:,:,:)
187            values_4d(len_lon+1,:,:,:) = values_4d(1,:,:,:)
188            if (state /= nf90_noerr) call handle_err(state)
189            state = nf90_put_var(ncid2, varids_2(i), values_4d)
190            if (state /= nf90_noerr) call handle_err(state)
191            state = nf90_redef(ncid2)
192            if (state /= nf90_noerr) call handle_err(state)
193            deallocate(tempvalues_4d)
194            deallocate(values_4d)
195        endif
196        deallocate(dimid_var)
197    enddo
198
199    state = nf90_enddef(ncid2)
200    if (state /= nf90_noerr) call handle_err(state)
201    state = nf90_close(ncid1)
202    if (state /= nf90_noerr) call handle_err(state)
203    state = nf90_close(ncid2)
204    if (state /= nf90_noerr) call handle_err(state)
205
206    deallocate(dimids,varids,dimids_2,varids_2)
207enddo
208
209END PROGRAM reshape_XIOS_output
Note: See TracBrowser for help on using the repository browser.