source: trunk/LMDZ.VENUS/Tools/startarchive2icosa/start_archive2icosa.f90 @ 3594

Last change on this file since 3594 was 3392, checked in by emillour, 7 months ago

Venus PCM:
Add some utilities to generate Venus PCM DYNAMICO start files
from a lon-lat start_archive.nc file
EM

File size: 20.2 KB
Line 
1PROGRAM start_archive2icosa
2
3  USE xios
4  USE mod_wait
5  USE netcdf
6 
7  IMPLICIT NONE
8  INCLUDE "mpif.h"
9  INTEGER :: rank
10  INTEGER :: size
11  INTEGER :: ierr
12
13  CHARACTER(len=*),PARAMETER :: id="client"
14  INTEGER :: comm
15  TYPE(xios_duration) :: dtime
16  CHARACTER(len=15) :: calendar_type
17  TYPE(xios_context) :: ctx_hdl
18
19  INTEGER :: n,l
20  INTEGER :: src_ibegin, src_iend, src_topo_ibegin, src_topo_iend
21  INTEGER :: src_jbegin, src_jend, src_topo_jbegin, src_topo_jend
22  INTEGER :: src_ni, src_ni_glo, src_topo_ni, src_topo_ni_glo
23  INTEGER :: src_nj, src_nj_glo, src_topo_nj, src_topo_nj_glo
24  INTEGER :: src_nlev ! number of vertical layers
25  INTEGER :: src_nq=1 ! number of tracers
26  INTEGER :: src_nt=1 ! number of time steps
27  DOUBLE PRECISION,ALLOCATABLE :: lev_values(:) ! vertical axis
28  DOUBLE PRECISION,ALLOCATABLE :: nq_values(:) ! tracer # axis
29  DOUBLE PRECISION,ALLOCATABLE :: src_lon(:) ! mesh center coordinate
30  DOUBLE PRECISION,ALLOCATABLE :: src_lat(:)
31  DOUBLE PRECISION,ALLOCATABLE :: src_ap(:)
32  DOUBLE PRECISION,ALLOCATABLE :: src_bp(:)
33  DOUBLE PRECISION,ALLOCATABLE :: src_controle(:)
34  DOUBLE PRECISION,ALLOCATABLE :: src_field_2D(:,:)
35  DOUBLE PRECISION,ALLOCATABLE :: src_pk(:,:)
36  DOUBLE PRECISION,ALLOCATABLE :: src_field_3D(:,:,:)
37  DOUBLE PRECISION,ALLOCATABLE :: src_pressure(:,:,:)
38  DOUBLE PRECISION,ALLOCATABLE :: src_theta_rhodz(:,:,:)
39  DOUBLE PRECISION,ALLOCATABLE :: src_tracers(:,:,:,:)
40  DOUBLE PRECISION,ALLOCATABLE :: src_topo_lon(:) ! mesh center coordinate
41  DOUBLE PRECISION,ALLOCATABLE :: src_topo_lat(:)
42  DOUBLE PRECISION,ALLOCATABLE :: src_topo(:,:)
43 
44  CHARACTER(LEN=*),PARAMETER :: src_file="start_archive_nc4.nc"
45  CHARACTER(LEN=*),PARAMETER :: src_topo_file="topo_Thomas_inv_nc4.nc"
46!  CHARACTER(LEN=*),PARAMETER :: src_topo_file="vt1x1inv_nc4.nc"
47  CHARACTER(LEN=*),PARAMETER :: dst_coord_file="start_icosa_ref.nc"
48  DOUBLE PRECISION,ALLOCATABLE :: dst_lon(:),dst_lat(:)
49  DOUBLE PRECISION,ALLOCATABLE :: dst_boundslon(:,:) ! mesh corner coordinates
50  DOUBLE PRECISION,ALLOCATABLE :: dst_boundslat(:,:)
51  INTEGER :: dst_ibegin !, dst_iend
52  INTEGER :: dst_ni, dst_ni_glo
53  INTEGER :: dst_nvertex
54  INTEGER :: ncid
55  INTEGER :: dimids(4)
56  INTEGER :: varid
57 
58  INTEGER :: div, remain
59  INTEGER :: ts ! time step #
60  DOUBLE PRECISION,PARAMETER :: pi=acos(-1.d0)
61  DOUBLE PRECISION :: gravity,kappa,preff
62
63!!! MPI Initialization
64  CALL MPI_INIT(ierr)
65  CALL init_wait
66
67!!! XIOS Initialization (get the local communicator)
68  CALL xios_initialize(id,return_comm=comm)
69! get local rank of MPI process
70  CALL MPI_COMM_RANK(comm,rank,ierr)
71! get total number of MPI processes
72  CALL MPI_COMM_SIZE(comm,size,ierr)
73
74!!! Open files and load sizes and coordinates
75  ierr=NF90_OPEN(src_topo_file, NF90_NOWRITE, ncid)
76  ierr=NF90_INQ_VARID(ncid,"RELIEF",varid)
77  ierr=NF90_INQUIRE_VARIABLE(ncid, varid,dimids=dimids)
78  write(*,*) "rank=",rank,"dimids=",dimids
79  ierr=NF90_INQUIRE_DIMENSION(ncid, dimids(1), len=src_topo_ni_glo)
80  ierr=NF90_INQUIRE_DIMENSION(ncid, dimids(2), len=src_topo_nj_glo)
81  write(*,*) "rank=",rank," src_topo_ni_glo=",src_topo_ni_glo ! longitude
82  write(*,*) "rank=",rank," src_topo_nj_glo=",src_topo_nj_glo ! latitude
83
84! assume domain splitup with MPI only along latitudes
85  src_topo_ni=src_topo_ni_glo
86  src_topo_ibegin=0
87  src_topo_iend=src_topo_ibegin+src_ni-1
88  write(*,*) "rank=",rank," src_topo_ni=",src_topo_ni
89 
90  src_topo_jbegin=0
91  DO n=0,size-1
92    src_topo_nj=src_topo_nj_glo/size
93    IF (n<MOD(src_topo_nj_glo,size)) src_topo_nj=src_topo_nj+1
94    IF (n==rank) exit
95    src_topo_jbegin=src_topo_jbegin+src_topo_nj
96  ENDDO
97  src_topo_jend=src_topo_jbegin+src_topo_nj-1
98  write(*,*) "rank=",rank," src_topo_nj=",src_topo_nj, &
99             " src_topo_jbegin=",src_topo_jbegin
100
101  ALLOCATE(src_topo_lon(src_topo_ni))
102  ALLOCATE(src_topo_lat(src_topo_nj))
103  ALLOCATE(src_topo(src_topo_ni,src_topo_nj))
104
105! load src_topo_lon and src_topo_lat
106  ierr=NF90_INQ_VARID(ncid,"longitude",varid)
107  ierr=NF90_GET_VAR(ncid,varid, src_topo_lon, &
108                    start=(/src_topo_ibegin+1/),count=(/src_topo_ni/))
109  WRITE(*,*) rank,":src_topo_lon(1:2)=",src_topo_lon(1:2)
110  ierr=NF90_INQ_VARID(ncid,"latitude",varid)
111  ierr=NF90_GET_VAR(ncid,varid, src_topo_lat, &
112                    start=(/src_topo_jbegin+1/),count=(/src_topo_nj/))
113  WRITE(*,*) rank,":src_topo_lat(1:2)=",src_topo_lat(1:2)
114
115! from start_archive.nc file
116  ierr=NF90_OPEN(src_file, NF90_NOWRITE, ncid)
117  ierr=NF90_INQ_VARID(ncid,"temp",varid)
118  ierr=NF90_INQUIRE_VARIABLE(ncid, varid,dimids=dimids)
119  write(*,*) "rank=",rank,"dimids=",dimids
120  ierr=NF90_INQUIRE_DIMENSION(ncid, dimids(1), len=src_ni_glo)
121  ierr=NF90_INQUIRE_DIMENSION(ncid, dimids(2), len=src_nj_glo)
122  ierr=NF90_INQUIRE_DIMENSION(ncid, dimids(3), len=src_nlev)
123  write(*,*) "rank=",rank," src_ni_glo=",src_ni_glo ! longitude
124  write(*,*) "rank=",rank," src_nj_glo=",src_nj_glo ! latitude
125  write(*,*) "rank=",rank," src_nlev=",src_nlev ! number of vertical layers
126  write(*,*) "rank=",rank," src_nq=",src_nq ! number of tracers
127
128! assume domain splitup with MPI only along latitudes
129  src_ni=src_ni_glo
130  src_ibegin=0
131  src_iend=src_ibegin+src_ni-1
132  write(*,*) "rank=",rank," src_ni=",src_ni
133 
134  src_jbegin=0
135  DO n=0,size-1
136    src_nj=src_nj_glo/size
137    IF (n<MOD(src_nj_glo,size)) src_nj=src_nj+1
138    IF (n==rank) exit
139    src_jbegin=src_jbegin+src_nj
140  ENDDO
141  src_jend=src_jbegin+src_nj-1
142  write(*,*) "rank=",rank," src_nj=",src_nj," src_jbegin=",src_jbegin
143
144  ALLOCATE(src_lon(src_ni))
145  ALLOCATE(src_lat(src_nj))
146  ALLOCATE(src_field_2D(src_ni,src_nj))
147  ALLOCATE(src_pk(src_ni,src_nj))
148  ALLOCATE(src_field_3D(src_ni,src_nj,src_nlev))
149  ALLOCATE(src_pressure(src_ni,src_nj,src_nlev+1))
150  ALLOCATE(src_theta_rhodz(src_ni,src_nj,src_nlev))
151  ALLOCATE(src_tracers(src_ni,src_nj,src_nlev,src_nq))
152
153! load src_lon and src_lat
154  ierr=NF90_INQ_VARID(ncid,"rlonv",varid)
155  ierr=NF90_GET_VAR(ncid,varid, src_lon, &
156                    start=(/src_ibegin+1/),count=(/src_ni/))
157! convert rad to deg
158  src_lon(1:src_ni)=src_lon(1:src_ni)*(180.d0/pi)
159  WRITE(*,*) rank,":src_lon=",src_lon
160  ierr=NF90_INQ_VARID(ncid,"rlatu",varid)
161  ierr=NF90_GET_VAR(ncid,varid, src_lat, &
162                    start=(/src_jbegin+1/),count=(/src_nj/))
163! convert rad to deg
164  src_lat(1:src_nj)=src_lat(1:src_nj)*(180.d0/pi)
165  WRITE(*,*) rank,":src_lat=",src_lat
166
167! load ap, bp and controle
168  ALLOCATE(src_ap(src_nlev+1),src_bp(src_nlev+1),src_controle(200))
169  ierr=NF90_INQ_VARID(ncid,"ap",varid)
170  ierr=NF90_GET_VAR(ncid,varid,src_ap)
171  WRITE(*,*) rank,":src_ap(1:5)=",src_ap(1:5)
172  ierr=NF90_INQ_VARID(ncid,"bp",varid)
173  ierr=NF90_GET_VAR(ncid,varid,src_bp)
174  WRITE(*,*) rank,":src_bp(1:5)=",src_bp(1:5)
175  ierr=NF90_INQ_VARID(ncid,"controle",varid)
176  ierr=NF90_GET_VAR(ncid,varid,src_controle)
177  gravity=src_controle(8)
178  WRITE(*,*) rank,":gravity=",gravity
179  kappa=src_controle(10)
180  WRITE(*,*) rank,":kappa=",kappa
181  preff=src_controle(19)
182  WRITE(*,*) rank,":preff=",preff
183
184! destination coordinates
185  ierr=NF90_OPEN(dst_coord_file, NF90_NOWRITE, ncid)
186  ierr=NF90_INQ_VARID(ncid,"bounds_lon_mesh",varid)
187  ierr=NF90_INQUIRE_VARIABLE(ncid, varid,dimids=dimids)
188  ierr=NF90_INQUIRE_DIMENSION(ncid, dimids(1), len=dst_nvertex)
189  ierr=NF90_INQUIRE_DIMENSION(ncid, dimids(2), len=dst_ni_glo)
190  write(*,*) "rank=",rank," dst_nvertex=",dst_nvertex ! vertex
191  write(*,*) "rank=",rank," dst_ni_glo=",dst_ni_glo ! vertex boundaries
192
193! evenly split into MPI domains
194  div    = dst_ni_glo/size
195  remain = MOD( dst_ni_glo, size )
196  IF (rank < remain) THEN
197    dst_ni=div+1 ;
198    dst_ibegin=rank*(div+1) ;
199  ELSE
200    dst_ni=div ;
201    dst_ibegin= remain * (div+1) + (rank-remain) * div ;
202  ENDIF
203  write(*,*) "rank=",rank," dst_ni=",dst_ni
204
205  ALLOCATE(dst_lon(dst_ni))
206  ALLOCATE(dst_lat(dst_ni))
207  ALLOCATE(dst_boundslon(dst_nvertex,dst_ni))
208  ALLOCATE(dst_boundslat(dst_nvertex,dst_ni))
209
210! load dst_lon, dst_lat, dst_boundslon and dst_boundslat
211  ierr=NF90_INQ_VARID(ncid,"lon_mesh",varid)
212  ierr=NF90_GET_VAR(ncid,varid, dst_lon, &
213                    start=(/dst_ibegin+1/), count=(/dst_ni/))
214  WRITE(*,*) rank,":dst_lon(1:5)=",dst_lon(1:5)
215  ierr=NF90_INQ_VARID(ncid,"lat_mesh",varid)
216  ierr=NF90_GET_VAR(ncid,varid, dst_lat, &
217                    start=(/dst_ibegin+1/), count=(/dst_ni/))
218  WRITE(*,*) rank,":dst_lat(1:5)=",dst_lat(1:5)
219  ierr=NF90_INQ_VARID(ncid,"bounds_lon_mesh",varid)
220  ierr=NF90_GET_VAR(ncid,varid,dst_boundslon, &
221                    start=(/1,dst_ibegin+1/), count=(/dst_nvertex,dst_ni/))
222  WRITE(*,*) rank,":dst_boundslon(:,1:2)=",dst_boundslon(:,1:2)
223  ierr=NF90_INQ_VARID(ncid,"bounds_lat_mesh",varid)
224  ierr=NF90_GET_VAR(ncid,varid, dst_boundslat, &
225                    start=(/1,dst_ibegin+1/), count=(/dst_nvertex,dst_ni/))
226  WRITE(*,*) rank,":dst_boundslat(:,1:2)=",dst_boundslat(:,1:2)
227
228
229! Initialize XIOS context
230  WRITE(*,*) rank,":CALL xios_context_initialize()"
231  CALL xios_context_initialize("test",comm)
232  CALL xios_get_handle("test",ctx_hdl)
233  CALL xios_set_current_context(ctx_hdl)
234
235! Set XIOS calendar and timestep
236  CALL xios_get_calendar_type(calendar_type)
237  WRITE(*,*) rank,":calendar_type = ", calendar_type
238  dtime%second = 3600
239  CALL xios_set_timestep(dtime)
240
241! Set axes
242  ALLOCATE(lev_values(src_nlev))
243  lev_values=(/ (l,l=1,src_nlev) /)
244  write(*,*) rank,":lev_values()=",lev_values
245  CALL xios_set_axis_attr("lev",n_glo=src_nlev,value=lev_values)
246  ALLOCATE(nq_values(src_nq))
247  nq_values=(/(l,l=1,src_nq)/)
248  CALL xios_set_axis_attr("nq",n_glo=src_nq,value=nq_values)
249
250! Set domains
251  CALL xios_set_domain_attr("src_domain_regular", &
252                            ni_glo=src_ni_glo, nj_glo=src_nj_glo, &
253                            ibegin=src_ibegin, ni=src_ni, &
254                            jbegin=src_jbegin, nj=src_nj, &
255                            type='rectilinear')
256  CALL xios_set_domain_attr("src_domain_regular", &
257                             data_dim=2, &
258                             data_ibegin=0, data_ni=src_ni, &
259                             data_jbegin=0, data_nj=src_nj)
260  CALL xios_set_domain_attr("src_domain_regular", &
261                            lonvalue_1D=src_lon, &
262                            latvalue_1D=src_lat)
263
264  CALL xios_set_domain_attr("src_topo_domain_regular", &
265                            ni_glo=src_topo_ni_glo, nj_glo=src_topo_nj_glo, &
266                            ibegin=src_topo_ibegin, ni=src_topo_ni, &
267                            jbegin=src_topo_jbegin, nj=src_topo_nj, &
268                            type='rectilinear')
269  CALL xios_set_domain_attr("src_topo_domain_regular", &
270                             data_dim=2, &
271                             data_ibegin=0, data_ni=src_topo_ni, &
272                             data_jbegin=0, data_nj=src_topo_nj)
273  CALL xios_set_domain_attr("src_topo_domain_regular", &
274                            lonvalue_1D=src_topo_lon, &
275                            latvalue_1D=src_topo_lat)
276
277  CALL xios_set_domain_attr("src_domain_regular_clean", &
278                            ni_glo=src_ni_glo-1, nj_glo=src_nj_glo, &
279                            ibegin=src_ibegin, ni=src_ni-1, &
280                            jbegin=src_jbegin, nj=src_nj, &
281                            type='rectilinear')
282  CALL xios_set_domain_attr("src_domain_regular_clean", &
283                             data_dim=2, &
284                             data_ibegin=0, data_ni=src_ni-1, &
285                             data_jbegin=0, data_nj=src_nj)
286  CALL xios_set_domain_attr("src_domain_regular_clean", &
287                            lonvalue_1D=src_lon(1:src_ni-1), &
288                            latvalue_1D=src_lat)
289
290  CALL xios_set_domain_attr("dst_domain_unstructured", &
291                            ni_glo=dst_ni_glo, &
292                            ibegin=dst_ibegin, &
293                            ni=dst_ni, &
294                            type="unstructured")
295  CALL xios_set_domain_attr("dst_domain_unstructured", &
296                            lonvalue_1D=dst_lon, &
297                            latvalue_1D=dst_lat, &
298                            bounds_lon_1D=dst_boundslon, &
299                            bounds_lat_1D=dst_boundslat, &
300                            nvertex=dst_nvertex)
301
302! Finalize XIOS context definition
303  WRITE(*,*) rank,":CALL xios_close_context_definition()"
304  CALL xios_close_context_definition()
305  CALL xios_get_handle("test",ctx_hdl)
306  CALL xios_set_current_context(ctx_hdl)
307
308! Temporal loop
309  DO ts=1,src_nt
310    WRITE(*,*) rank,":ts=",ts
311    ! Update calendar
312    CALL xios_update_calendar(ts)
313
314    ! Topography
315    CALL xios_recv_field("RELIEF",src_topo)
316    WRITE(*,*) rank,":topo(1:2,1:3)=",src_topo(1:2,1:3)
317    ! Send surface geopotential
318    CALL xios_send_field("topo",src_topo(:,:)*gravity)
319
320    ! Surface pressure:
321    !! get data using XIOS:
322    CALL xios_recv_field("src_ps",src_field_2D)
323    WRITE(*,*) rank,":src_ps(1:2,1:3)=",src_field_2D(1:2,1:3)
324    !! write data using XIOS
325    CALL xios_send_field("ps_clean",src_field_2D(1:src_ni-1,1:src_nj))
326
327    ! compute inter-layer pressures
328    DO l=1,src_nlev+1
329      src_pressure(:,:,l) = src_ap(l)+src_bp(l)*src_field_2D(:,:)
330    ENDDO
331   
332    ! surface temperature:
333    CALL xios_recv_field("src_tsurf",src_field_2D)
334    WRITE(*,*) rank,":src_tsurf(1:2,1:3)=",src_field_2D(1:2,1:3)
335    CALL xios_send_field("tsurf_clean",src_field_2D(1:src_ni-1,1:src_nj))
336   
337    ! Temperature:
338    CALL xios_recv_field("src_temp",src_field_3D)
339    ! compute theta_rhodz
340    DO l=1,src_nlev
341      src_pk(:,:)=((.5/preff)*(src_pressure(:,:,l)+src_pressure(:,:,l+1)))**kappa
342      src_theta_rhodz(:,:,l) = src_field_3D(:,:,l) * &
343      ((src_pressure(:,:,l)-src_pressure(:,:,l+1))/gravity)/src_pk(:,:)
344    ENDDO
345    CALL xios_send_field("theta_rhodz_clean", &
346                         src_theta_rhodz(1:src_ni-1,1:src_nj,1:src_nlev))
347
348    ! zonal wind
349    CALL xios_recv_field("src_u",src_field_3D)
350    CALL xios_send_field("u_clean", &
351                          src_field_3D(1:src_ni-1,1:src_nj,1:src_nlev))
352   
353    ! meridional wind
354    CALL xios_recv_field("src_v",src_field_3D)
355    CALL xios_send_field("v_clean", &
356                          src_field_3D(1:src_ni-1,1:src_nj,1:src_nlev))
357    ! tracers
358    !CALL xios_recv_field("src_tracer001",src_tracers(:,:,:,1))
359    CALL xios_recv_field("src_co2",src_tracers(:,:,:,1))
360    CALL xios_send_field("tracers_clean", &
361                          src_tracers(1:src_ni-1,1:src_nj,1:src_nlev,1:src_nq))
362
363    ! subsurface temperatures
364    CALL xios_recv_field("src_Tsoil01",src_field_2D)
365    WRITE(*,*) rank,":src_Tsoil01(1:2,1:3)=",src_field_2D(1:2,1:3)
366    CALL xios_send_field("Tsoil01_clean",src_field_2D(1:src_ni-1,1:src_nj))
367    CALL xios_recv_field("src_Tsoil02",src_field_2D)
368    CALL xios_send_field("Tsoil02_clean",src_field_2D(1:src_ni-1,1:src_nj))
369    CALL xios_recv_field("src_Tsoil03",src_field_2D)
370    CALL xios_send_field("Tsoil03_clean",src_field_2D(1:src_ni-1,1:src_nj))
371    CALL xios_recv_field("src_Tsoil04",src_field_2D)
372    CALL xios_send_field("Tsoil04_clean",src_field_2D(1:src_ni-1,1:src_nj))
373    CALL xios_recv_field("src_Tsoil05",src_field_2D)
374    CALL xios_send_field("Tsoil05_clean",src_field_2D(1:src_ni-1,1:src_nj))
375    CALL xios_recv_field("src_Tsoil06",src_field_2D)
376    CALL xios_send_field("Tsoil06_clean",src_field_2D(1:src_ni-1,1:src_nj))
377    CALL xios_recv_field("src_Tsoil07",src_field_2D)
378    CALL xios_send_field("Tsoil07_clean",src_field_2D(1:src_ni-1,1:src_nj))
379    CALL xios_recv_field("src_Tsoil08",src_field_2D)
380    CALL xios_send_field("Tsoil08_clean",src_field_2D(1:src_ni-1,1:src_nj))
381    CALL xios_recv_field("src_Tsoil09",src_field_2D)
382    CALL xios_send_field("Tsoil09_clean",src_field_2D(1:src_ni-1,1:src_nj))
383    CALL xios_recv_field("src_Tsoil10",src_field_2D)
384    CALL xios_send_field("Tsoil10_clean",src_field_2D(1:src_ni-1,1:src_nj))
385    CALL xios_recv_field("src_Tsoil11",src_field_2D)
386    CALL xios_send_field("Tsoil11_clean",src_field_2D(1:src_ni-1,1:src_nj))
387
388    ! Albedo
389    CALL xios_recv_field("src_albe",src_field_2D)
390    WRITE(*,*) rank,":src_albe(1:2,1:3)=",src_field_2D(1:2,1:3)
391    CALL xios_send_field("ALBE_clean",src_field_2D(1:src_ni-1,1:src_nj))
392   
393    ! SW flux at the surface
394    CALL xios_recv_field("src_solsw",src_field_2D)
395    WRITE(*,*) rank,":src_solsw(1:2,1:3)=",src_field_2D(1:2,1:3)
396    CALL xios_send_field("solsw_clean",src_field_2D(1:src_ni-1,1:src_nj))
397   
398    ! LW flux at the surface
399    CALL xios_recv_field("src_sollw",src_field_2D)
400    WRITE(*,*) rank,":src_sollw(1:2,1:3)=",src_field_2D(1:2,1:3)
401    CALL xios_send_field("sollw_clean",src_field_2D(1:src_ni-1,1:src_nj))
402   
403    ! fder "Derive de flux"
404    CALL xios_recv_field("src_fder",src_field_2D)
405    WRITE(*,*) rank,":src_fder(1:2,1:3)=",src_field_2D(1:2,1:3)
406    CALL xios_send_field("fder_clean",src_field_2D(1:src_ni-1,1:src_nj))
407   
408    ! dlw "Derivee flux IR"
409    CALL xios_recv_field("src_dlw",src_field_2D)
410    WRITE(*,*) rank,":src_dlw(1:2,1:3)=",src_field_2D(1:2,1:3)
411    CALL xios_send_field("dlw_clean",src_field_2D(1:src_ni-1,1:src_nj))
412   
413    ! sollwdown "Flux IR vers le bas a la surface"
414    CALL xios_recv_field("src_sollwdown",src_field_2D)
415    WRITE(*,*) rank,":src_sollwdown(1:2,1:3)=",src_field_2D(1:2,1:3)
416    CALL xios_send_field("sollwdown_clean",src_field_2D(1:src_ni-1,1:src_nj))
417
418    ! RADS "Net flux at surface"
419    CALL xios_recv_field("src_RADS",src_field_2D)
420    WRITE(*,*) rank,":src_RADS(1:2,1:3)=",src_field_2D(1:2,1:3)
421    CALL xios_send_field("RADS_clean",src_field_2D(1:src_ni-1,1:src_nj))
422
423    ! ZMEA "zmea Orographie sous-maille"
424    CALL xios_recv_field("src_ZMEA",src_field_2D)
425    WRITE(*,*) rank,":src_ZMEA(1:2,1:3)=",src_field_2D(1:2,1:3)
426    CALL xios_send_field("ZMEA_clean",src_field_2D(1:src_ni-1,1:src_nj))
427   
428    ! ZSTD "zstd Orographie sous-maille"
429    CALL xios_recv_field("src_ZSTD",src_field_2D)
430    WRITE(*,*) rank,":src_ZSTD(1:2,1:3)=",src_field_2D(1:2,1:3)
431    CALL xios_send_field("ZSTD_clean",src_field_2D(1:src_ni-1,1:src_nj))
432   
433    ! ZSIG "zsig Orographie sous-maille"
434    CALL xios_recv_field("src_ZSIG",src_field_2D)
435    WRITE(*,*) rank,":src_ZSIG(1:2,1:3)=",src_field_2D(1:2,1:3)
436    CALL xios_send_field("ZSIG_clean",src_field_2D(1:src_ni-1,1:src_nj))
437   
438    ! ZGAM "zgam Orographie sous-maille"
439    CALL xios_recv_field("src_ZGAM",src_field_2D)
440    WRITE(*,*) rank,":src_ZGAM(1:2,1:3)=",src_field_2D(1:2,1:3)
441    CALL xios_send_field("ZGAM_clean",src_field_2D(1:src_ni-1,1:src_nj))
442   
443    ! ZTHE "zthe Orographie sous-maille"
444    CALL xios_recv_field("src_ZTHE",src_field_2D)
445    WRITE(*,*) rank,":src_ZTHE(1:2,1:3)=",src_field_2D(1:2,1:3)
446    CALL xios_send_field("ZTHE_clean",src_field_2D(1:src_ni-1,1:src_nj))
447   
448    ! ZPIC "zpic Orographie sous-maille"
449    CALL xios_recv_field("src_ZPIC",src_field_2D)
450    WRITE(*,*) rank,":src_ZPIC(1:2,1:3)=",src_field_2D(1:2,1:3)
451    CALL xios_send_field("ZPIC_clean",src_field_2D(1:src_ni-1,1:src_nj))
452   
453    ! ZVAL "zval Orographie sous-maille"
454    CALL xios_recv_field("src_ZVAL",src_field_2D)
455    WRITE(*,*) rank,":src_ZVAL(1:2,1:3)=",src_field_2D(1:2,1:3)
456    CALL xios_send_field("ZVAL_clean",src_field_2D(1:src_ni-1,1:src_nj))
457   
458    ! TANCIEN => not necessary
459  ENDDO ! of DO ts=1,src_nt
460 
461!! Finalize
462  write(*,*) rank,":Finalize: call xios_context_finalize"
463  CALL xios_context_finalize()
464
465  write(*,*) rank,":Finalize: call MPI_COMM_FREE"
466  CALL MPI_COMM_FREE(comm, ierr)
467
468  write(*,*) rank,":Finalize: call xios_finalize"
469  CALL xios_finalize()
470
471  if (rank==0) then
472    ! add a couple of things in the "startphy_icosa.nc" file
473    write(*,*) rank,"Write controle() to startphy_icosa.nc"
474    ierr=NF90_OPEN("startphy_icosa.nc",NF90_WRITE,ncid)
475    ierr=NF90_REDEF(ncid) ! switch to define mode
476    ierr=NF90_DEF_DIM(ncid,"index",100,dimids(1))
477    ierr=NF90_DEF_VAR(ncid,"controle",NF90_DOUBLE,dimids(1),varid)
478    ierr=NF90_ENDDEF(ncid) ! switch out of define mode
479    ierr=NF90_PUT_VAR(ncid,varid,src_controle(101:200))
480    if (ierr.ne.NF90_NOERR) then
481      write(*,*) "NetCDF Error:",NF90_STRERROR(ierr)
482    endif
483    ierr=NF90_CLOSE(ncid)
484   
485    ! add a couple of things in the "start_icosa.nc" file
486    ierr=NF90_OPEN("start_icosa.nc",NF90_WRITE,ncid)
487    ierr=NF90_REDEF(ncid)
488    ierr=NF90_DEF_DIM(ncid,"nvertex_u",2,dimids(1))
489    ierr=NF90_DEF_VAR(ncid,"iteration",NF90_FLOAT,varid)
490    ierr=NF90_ENDDEF(ncid)
491    if (ierr.ne.NF90_NOERR) then
492      write(*,*) "NetCDF Error:",NF90_STRERROR(ierr)
493    endif
494    ierr=NF90_PUT_VAR(ncid,varid,0) ! set "iteration" value to 0
495    ierr=NF90_CLOSE(ncid)
496   
497  endif ! of if (rank==0)
498
499  write(*,*) rank,":Finalize: call MPI_FINALIZE"
500  CALL MPI_FINALIZE(ierr)
501
502  write(*,*) rank,":my_remap: all is well that ends well!"
503
504END PROGRAM start_archive2icosa
Note: See TracBrowser for help on using the repository browser.