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

Last change on this file since 3093 was 3076, checked in by jbclement, 14 months ago

PEM:
Big cleaning/improvements of the PEM:

  • Conversion of "abort_pem.F" and "soil_settings_PEM.F" into Fortran 90;
  • Transformation of every PEM subroutines into module;
  • Rewriting of many subroutines with modern Fortran syntax;
  • Correction of a bug in "pem.F90" when calling 'recomp_tend_co2_slope'. The arguments were given in disorder and emissivity was missing;
  • Update of "launch_pem.sh" in deftank.

JBC

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