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

Last change on this file since 2965 was 2963, checked in by romain.vande, 3 years ago

Mars PEM :

Adapt PEM to the subslope PCM configuration, it is now fully compatible.

Create a PEM folder in deftank that contains:

run_pem1: a bash file that runs chained simulation of PEM as well as running a parameterizable number of PCM simulation in between.

It also takes care of reshaping XIOS output as well as renaming outputs etc… in the spirit of run_month1.

It is written for Irene machine and the header needs to be adapted for other machines.

run_PEM.def: A text file that shows the possible parameters to choose before a PEM simulation.

It should be included at the end of run.def just like callphys.def

ob_ex_lsp.asc: An ascii file containing the obliquity, eccentricity, ls_peri data from Laskar in Martian year.
README: A txt file explaining the content of the folder

Adapt field_def_physics_mars.xml to consider the case with 7 subslopes in the PCM.
Change context_lmdz_physics.xml to be able to output the file needed by the PEM.

Correct a few other minor bugs.

RV

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