source: trunk/LMDZ.MARS/util/startarchive2icosa/start_archive2icosa.f90 @ 2813

Last change on this file since 2813 was 2532, checked in by adelavois, 4 years ago

update of start_archive2icosa tool

File size: 22.3 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=7 ! number of tracers
26  INTEGER :: src_nt=1 ! number of time steps
27  INTEGER :: src_nsoil=18 ! number of soil layers
28  DOUBLE PRECISION,ALLOCATABLE :: lev_values(:) ! vertical axis
29  DOUBLE PRECISION,ALLOCATABLE :: lev_p1_values(:) ! vertical axis
30  DOUBLE PRECISION,ALLOCATABLE :: nq_values(:) ! tracer # axis
31  DOUBLE PRECISION,ALLOCATABLE :: soil_layers_values(:) ! soil axis
32  DOUBLE PRECISION,ALLOCATABLE :: src_lon(:) ! mesh center coordinate
33  DOUBLE PRECISION,ALLOCATABLE :: src_lat(:)
34  DOUBLE PRECISION,ALLOCATABLE :: src_ap(:)
35  DOUBLE PRECISION,ALLOCATABLE :: src_bp(:)
36  DOUBLE PRECISION,ALLOCATABLE :: src_controle(:)
37  DOUBLE PRECISION,ALLOCATABLE :: src_field_2D(:,:)
38  DOUBLE PRECISION,ALLOCATABLE :: src_pk(:,:)
39  DOUBLE PRECISION,ALLOCATABLE :: src_field_3D(:,:,:)
40  DOUBLE PRECISION,ALLOCATABLE :: src_pressure(:,:,:)
41  DOUBLE PRECISION,ALLOCATABLE :: src_theta_rhodz(:,:,:)
42  DOUBLE PRECISION,ALLOCATABLE :: src_topo_lon(:) ! mesh center coordinate
43  DOUBLE PRECISION,ALLOCATABLE :: src_topo_lat(:)
44  DOUBLE PRECISION,ALLOCATABLE :: src_topo(:,:)
45  DOUBLE PRECISION,ALLOCATABLE :: src_field_3D_soil(:,:,:) ! 3D grid in soil, for tsoil et inertia
46  DOUBLE PRECISION,ALLOCATABLE :: src_soil_layers(:) ! soil_depth
47  DOUBLE PRECISION,ALLOCATABLE :: src_field_3D_p1(:,:,:) !q2 field
48 
49  CHARACTER(LEN=*),PARAMETER :: src_file="start_archive_nc4.nc"
50  CHARACTER(LEN=*),PARAMETER :: src_topo_file="surface_nc4.nc"
51  CHARACTER(LEN=*),PARAMETER :: dst_coord_file="start_icosa_ref.nc"
52  CHARACTER(LEN=*),PARAMETER :: src_controle_file="startphy_icosa_ref.nc"
53  CHARACTER(LEN=*),PARAMETER :: output_start_file="start_icosa_prefinalize.nc"
54  CHARACTER(LEN=*),PARAMETER :: output_startfi_file="startfi_prefinalize.nc"
55  DOUBLE PRECISION,ALLOCATABLE :: dst_lon(:),dst_lat(:)
56  DOUBLE PRECISION,ALLOCATABLE :: dst_boundslon(:,:) ! mesh corner coordinates
57  DOUBLE PRECISION,ALLOCATABLE :: dst_boundslat(:,:)
58  INTEGER :: dst_ibegin !, dst_iend
59  INTEGER :: dst_ni, dst_ni_glo
60  INTEGER :: dst_nvertex
61  INTEGER :: ncid
62  INTEGER :: dimids(4)
63  INTEGER :: varid
64 
65  INTEGER :: div, remain
66  INTEGER :: ts ! time step #
67  DOUBLE PRECISION,PARAMETER :: pi=acos(-1.d0)
68  DOUBLE PRECISION :: gravity,kappa,preff
69! Tracers
70  CHARACTER(LEN=11)          :: i_trac,format_string
71  INTEGER                    :: i
72
73!!! MPI Initialization
74  CALL MPI_INIT(ierr)
75  CALL init_wait
76
77!!! XIOS Initialization (get the local communicator)
78  CALL xios_initialize(id,return_comm=comm)
79! get local rank of MPI process
80  CALL MPI_COMM_RANK(comm,rank,ierr)
81! get total number of MPI processes
82  CALL MPI_COMM_SIZE(comm,size,ierr)
83
84!!! Open files and load sizes and coordinates
85  ierr=NF90_OPEN(src_topo_file, NF90_NOWRITE, ncid)
86  ierr=NF90_INQ_VARID(ncid,"zMOL",varid)
87  ierr=NF90_INQUIRE_VARIABLE(ncid, varid,dimids=dimids)
88  write(*,*) "rank=",rank,"dimids=",dimids
89  ierr=NF90_INQUIRE_DIMENSION(ncid, dimids(1), len=src_topo_ni_glo)
90  ierr=NF90_INQUIRE_DIMENSION(ncid, dimids(2), len=src_topo_nj_glo)
91  write(*,*) "rank=",rank," src_topo_ni_glo=",src_topo_ni_glo ! longitude
92  write(*,*) "rank=",rank," src_topo_nj_glo=",src_topo_nj_glo ! latitude
93
94! assume domain splitup with MPI only along latitudes
95  src_topo_ni=src_topo_ni_glo
96  src_topo_ibegin=0
97  src_topo_iend=src_topo_ibegin+src_ni-1
98  write(*,*) "rank=",rank," src_topo_ni=",src_topo_ni
99 
100  src_topo_jbegin=0
101  DO n=0,size-1
102    src_topo_nj=src_topo_nj_glo/size
103    IF (n<MOD(src_topo_nj_glo,size)) src_topo_nj=src_topo_nj+1
104    IF (n==rank) exit
105    src_topo_jbegin=src_topo_jbegin+src_topo_nj
106  ENDDO
107  src_topo_jend=src_topo_jbegin+src_topo_nj-1
108  write(*,*) "rank=",rank," src_topo_nj=",src_topo_nj, &
109             " src_topo_jbegin=",src_topo_jbegin
110
111  ALLOCATE(src_topo_lon(src_topo_ni))
112  ALLOCATE(src_topo_lat(src_topo_nj))
113  ALLOCATE(src_topo(src_topo_ni,src_topo_nj))
114
115! load src_topo_lon and src_topo_lat
116  ierr=NF90_INQ_VARID(ncid,"longitude",varid)
117  ierr=NF90_GET_VAR(ncid,varid, src_topo_lon, &
118                    start=(/src_topo_ibegin+1/),count=(/src_topo_ni/))
119  WRITE(*,*) rank,":src_topo_lon(1:2)=",src_topo_lon(1:2)
120  ierr=NF90_INQ_VARID(ncid,"latitude",varid)
121  ierr=NF90_GET_VAR(ncid,varid, src_topo_lat, &
122                    start=(/src_topo_jbegin+1/),count=(/src_topo_nj/))
123  WRITE(*,*) rank,":src_topo_lat(1:2)=",src_topo_lat(1:2)
124
125! from start_archive.nc file
126  ierr=NF90_OPEN(src_file, NF90_NOWRITE, ncid)
127  ierr=NF90_INQ_VARID(ncid,"temp",varid)
128  ierr=NF90_INQUIRE_VARIABLE(ncid, varid,dimids=dimids)
129  write(*,*) "rank=",rank,"dimids=",dimids
130  ierr=NF90_INQUIRE_DIMENSION(ncid, dimids(1), len=src_ni_glo)
131  ierr=NF90_INQUIRE_DIMENSION(ncid, dimids(2), len=src_nj_glo)
132  ierr=NF90_INQUIRE_DIMENSION(ncid, dimids(3), len=src_nlev)
133  write(*,*) "rank=",rank," src_ni_glo=",src_ni_glo ! longitude
134  write(*,*) "rank=",rank," src_nj_glo=",src_nj_glo ! latitude
135  write(*,*) "rank=",rank," src_nlev=",src_nlev ! number of vertical layers
136  write(*,*) "rank=",rank," src_nq=",src_nq ! number of tracers
137! soil_depth with tsoil variable
138  ierr=NF90_INQ_VARID(ncid,"tsoil",varid)
139  ierr=NF90_INQUIRE_VARIABLE(ncid, varid,dimids=dimids)
140  write(*,*) "rank=",rank,"dimids tsoil =",dimids
141  ierr=NF90_INQUIRE_DIMENSION(ncid, dimids(3), len=src_nsoil)
142  write(*,*) "rank=",rank," src_nsoil=",src_nsoil ! number of vertical layers
143
144! assume domain splitup with MPI only along latitudes
145  src_ni=src_ni_glo
146  src_ibegin=0
147  src_iend=src_ibegin+src_ni-1
148  write(*,*) "rank=",rank," src_ni=",src_ni
149 
150  src_jbegin=0
151  DO n=0,size-1
152    src_nj=src_nj_glo/size
153    IF (n<MOD(src_nj_glo,size)) src_nj=src_nj+1
154    IF (n==rank) exit
155    src_jbegin=src_jbegin+src_nj
156  ENDDO
157  src_jend=src_jbegin+src_nj-1
158  write(*,*) "rank=",rank," src_nj=",src_nj," src_jbegin=",src_jbegin
159
160  ALLOCATE(src_lon(src_ni))
161  ALLOCATE(src_lat(src_nj))
162  ALLOCATE(src_field_2D(src_ni,src_nj))
163  ALLOCATE(src_pk(src_ni,src_nj))
164  ALLOCATE(src_field_3D(src_ni,src_nj,src_nlev))
165  ALLOCATE(src_field_3D_p1(src_ni,src_nj,src_nlev+1))
166  ALLOCATE(src_field_3D_soil(src_ni,src_nj,src_nsoil))
167!  ALLOCATE(src_field_soil_layers(src_nsoil))
168  ALLOCATE(src_pressure(src_ni,src_nj,src_nlev+1))
169  ALLOCATE(src_theta_rhodz(src_ni,src_nj,src_nlev))
170
171! load src_lon and src_lat
172  ierr=NF90_INQ_VARID(ncid,"rlonv",varid)
173  ierr=NF90_GET_VAR(ncid,varid, src_lon, &
174                    start=(/src_ibegin+1/),count=(/src_ni/))
175! convert rad to deg
176  src_lon(1:src_ni)=src_lon(1:src_ni)*(180.d0/pi)
177  WRITE(*,*) rank,":src_lon=",src_lon
178  ierr=NF90_INQ_VARID(ncid,"rlatu",varid)
179  ierr=NF90_GET_VAR(ncid,varid, src_lat, &
180                    start=(/src_jbegin+1/),count=(/src_nj/))
181! convert rad to deg
182  src_lat(1:src_nj)=src_lat(1:src_nj)*(180.d0/pi)
183  WRITE(*,*) rank,":src_lat=",src_lat
184
185! load ap, bp and controle
186  ALLOCATE(src_ap(src_nlev+1),src_bp(src_nlev+1))
187  ierr=NF90_INQ_VARID(ncid,"ap",varid)
188  ierr=NF90_GET_VAR(ncid,varid,src_ap)
189  WRITE(*,*) rank,":src_ap(1:5)=",src_ap(1:5)
190  ierr=NF90_INQ_VARID(ncid,"bp",varid)
191  ierr=NF90_GET_VAR(ncid,varid,src_bp)
192  WRITE(*,*) rank,":src_bp(1:5)=",src_bp(1:5)
193 
194! controle is taken in startphy_icosa_ref as start_archive_nc4 is too old
195  ierr=NF90_OPEN(src_controle_file, NF90_NOWRITE, ncid)
196  ALLOCATE(src_controle(100))
197  ierr=NF90_INQ_VARID(ncid,"controle",varid)
198  ierr=NF90_GET_VAR(ncid,varid,src_controle)
199  ! day_ini set to 0 as startphy_icosa_ref is a restart
200  src_controle(3)=0
201  WRITE(*,*) rank,":src_controle(1:5)=",src_controle(1:5)
202  gravity=src_controle(7)
203  WRITE(*,*) rank,":gravity=",gravity
204  kappa=src_controle(9)
205  WRITE(*,*) rank,":kappa=",kappa
206!  preff=src_controle(18)
207! AD: Warning preff is hardcoded because not in controle of restartfi
208  preff=610.
209  WRITE(*,*) rank,":preff=",preff
210
211! destination coordinates
212  ierr=NF90_OPEN(dst_coord_file, NF90_NOWRITE, ncid)
213  ierr=NF90_INQ_VARID(ncid,"bounds_lon_mesh",varid)
214  ierr=NF90_INQUIRE_VARIABLE(ncid, varid,dimids=dimids)
215  ierr=NF90_INQUIRE_DIMENSION(ncid, dimids(1), len=dst_nvertex)
216  ierr=NF90_INQUIRE_DIMENSION(ncid, dimids(2), len=dst_ni_glo)
217  write(*,*) "rank=",rank," dst_nvertex=",dst_nvertex ! vertex
218  write(*,*) "rank=",rank," dst_ni_glo=",dst_ni_glo ! vertex boundaries
219
220! evenly split into MPI domains
221  div    = dst_ni_glo/size
222  remain = MOD( dst_ni_glo, size )
223  IF (rank < remain) THEN
224    dst_ni=div+1 ;
225    dst_ibegin=rank*(div+1) ;
226  ELSE
227    dst_ni=div ;
228    dst_ibegin= remain * (div+1) + (rank-remain) * div ;
229  ENDIF
230  write(*,*) "rank=",rank," dst_ni=",dst_ni
231
232  ALLOCATE(dst_lon(dst_ni))
233  ALLOCATE(dst_lat(dst_ni))
234  ALLOCATE(dst_boundslon(dst_nvertex,dst_ni))
235  ALLOCATE(dst_boundslat(dst_nvertex,dst_ni))
236
237! load dst_lon, dst_lat, dst_boundslon and dst_boundslat
238  ierr=NF90_INQ_VARID(ncid,"lon_mesh",varid)
239  ierr=NF90_GET_VAR(ncid,varid, dst_lon, &
240                    start=(/dst_ibegin+1/), count=(/dst_ni/))
241  WRITE(*,*) rank,":dst_lon(1:5)=",dst_lon(1:5)
242  ierr=NF90_INQ_VARID(ncid,"lat_mesh",varid)
243  ierr=NF90_GET_VAR(ncid,varid, dst_lat, &
244                    start=(/dst_ibegin+1/), count=(/dst_ni/))
245  WRITE(*,*) rank,":dst_lat(1:5)=",dst_lat(1:5)
246  ierr=NF90_INQ_VARID(ncid,"bounds_lon_mesh",varid)
247  ierr=NF90_GET_VAR(ncid,varid,dst_boundslon, &
248                    start=(/1,dst_ibegin+1/), count=(/dst_nvertex,dst_ni/))
249  WRITE(*,*) rank,":dst_boundslon(:,1:2)=",dst_boundslon(:,1:2)
250  ierr=NF90_INQ_VARID(ncid,"bounds_lat_mesh",varid)
251  ierr=NF90_GET_VAR(ncid,varid, dst_boundslat, &
252                    start=(/1,dst_ibegin+1/), count=(/dst_nvertex,dst_ni/))
253  WRITE(*,*) rank,":dst_boundslat(:,1:2)=",dst_boundslat(:,1:2)
254
255
256! Initialize XIOS context
257  WRITE(*,*) rank,":CALL xios_context_initialize()"
258  CALL xios_context_initialize("test",comm)
259  CALL xios_get_handle("test",ctx_hdl)
260  CALL xios_set_current_context(ctx_hdl)
261
262! Set XIOS calendar and timestep
263  CALL xios_get_calendar_type(calendar_type)
264  WRITE(*,*) rank,":calendar_type = ", calendar_type
265  dtime%second = 3600
266  CALL xios_set_timestep(dtime)
267
268! Set axes
269  ! vertical atm axis
270  ALLOCATE(lev_values(src_nlev))
271  lev_values=(/ (l,l=1,src_nlev) /)
272  CALL xios_set_axis_attr("lev",n_glo=src_nlev,value=lev_values)
273  ! lev+1 for q2
274  ALLOCATE(lev_p1_values(src_nlev+1))
275  lev_p1_values=(/ (l,l=1,src_nlev+1) /)
276  CALL xios_set_axis_attr("lev_p1",n_glo=src_nlev+1,value=lev_p1_values)
277  ! tracers axis
278  ALLOCATE(nq_values(src_nq))
279  nq_values=(/(l,l=1,src_nq)/)
280  CALL xios_set_axis_attr("nq",n_glo=src_nq,value=nq_values)
281  ! soil layers axis
282  ALLOCATE(soil_layers_values(src_nsoil))
283  soil_layers_values=(/ (l,l=1,src_nsoil) /)
284  CALL xios_set_axis_attr("soil_layers",n_glo=src_nsoil,value=soil_layers_values)
285
286! Set domains
287  CALL xios_set_domain_attr("src_domain_regular", &
288                            ni_glo=src_ni_glo, nj_glo=src_nj_glo, &
289                            ibegin=src_ibegin, ni=src_ni, &
290                            jbegin=src_jbegin, nj=src_nj, &
291                            type='rectilinear')
292  CALL xios_set_domain_attr("src_domain_regular", &
293                             data_dim=2, &
294                             data_ibegin=0, data_ni=src_ni, &
295                             data_jbegin=0, data_nj=src_nj)
296  CALL xios_set_domain_attr("src_domain_regular", &
297                            lonvalue_1D=src_lon, &
298                            latvalue_1D=src_lat)
299
300  CALL xios_set_domain_attr("src_topo_domain_regular", &
301                            ni_glo=src_topo_ni_glo, nj_glo=src_topo_nj_glo, &
302                            ibegin=src_topo_ibegin, ni=src_topo_ni, &
303                            jbegin=src_topo_jbegin, nj=src_topo_nj, &
304                            type='rectilinear')
305  CALL xios_set_domain_attr("src_topo_domain_regular", &
306                             data_dim=2, &
307                             data_ibegin=0, data_ni=src_topo_ni, &
308                             data_jbegin=0, data_nj=src_topo_nj)
309  CALL xios_set_domain_attr("src_topo_domain_regular", &
310                            lonvalue_1D=src_topo_lon, &
311                            latvalue_1D=src_topo_lat)
312
313  CALL xios_set_domain_attr("src_domain_regular_clean", &
314                            ni_glo=src_ni_glo-1, nj_glo=src_nj_glo, &
315                            ibegin=src_ibegin, ni=src_ni-1, &
316                            jbegin=src_jbegin, nj=src_nj, &
317                            type='rectilinear')
318  CALL xios_set_domain_attr("src_domain_regular_clean", &
319                             data_dim=2, &
320                             data_ibegin=0, data_ni=src_ni-1, &
321                             data_jbegin=0, data_nj=src_nj)
322  CALL xios_set_domain_attr("src_domain_regular_clean", &
323                            lonvalue_1D=src_lon(1:src_ni-1), &
324                            latvalue_1D=src_lat)
325
326  CALL xios_set_domain_attr("dst_domain_unstructured", &
327                            ni_glo=dst_ni_glo, &
328                            ibegin=dst_ibegin, &
329                            ni=dst_ni, &
330                            type="unstructured")
331  CALL xios_set_domain_attr("dst_domain_unstructured", &
332                            lonvalue_1D=dst_lon, &
333                            latvalue_1D=dst_lat, &
334                            bounds_lon_1D=dst_boundslon, &
335                            bounds_lat_1D=dst_boundslat, &
336                            nvertex=dst_nvertex)
337
338! Finalize XIOS context definition
339  WRITE(*,*) rank,":CALL xios_close_context_definition()"
340  CALL xios_close_context_definition()
341  CALL xios_get_handle("test",ctx_hdl)
342  CALL xios_set_current_context(ctx_hdl)
343
344! Temporal loop
345  DO ts=1,src_nt
346    WRITE(*,*) rank,":ts=",ts
347    ! Update calendar
348    CALL xios_update_calendar(ts)
349
350    ! Topography
351    CALL xios_recv_field("zMOL",src_topo)
352    WRITE(*,*) rank,":topo(1:2,1:3)=",src_topo(1:2,1:3)
353    ! Send surface geopotential
354    CALL xios_send_field("topo",src_topo(:,:)*gravity*1000.)
355
356    ! Surface pressure:
357    !! get data using XIOS:
358    CALL xios_recv_field("src_ps",src_field_2D)
359    WRITE(*,*) rank,":src_ps(1:2,1:3)=",src_field_2D(1:2,1:3)
360    !! write data using XIOS
361    CALL xios_send_field("ps_clean",src_field_2D(1:src_ni-1,1:src_nj))
362
363    ! compute inter-layer pressures
364    DO l=1,src_nlev+1
365      src_pressure(:,:,l) = src_ap(l)+src_bp(l)*src_field_2D(:,:)
366    ENDDO
367   
368    ! surface temperature:
369    CALL xios_recv_field("src_tsurf",src_field_2D)
370    WRITE(*,*) rank,":src_tsurf(1:2,1:3)=",src_field_2D(1:2,1:3)
371    CALL xios_send_field("tsurf_clean",src_field_2D(1:src_ni-1,1:src_nj))
372
373    ! co2 ice coverage:
374    CALL xios_recv_field("src_co2ice",src_field_2D)
375    WRITE(*,*) rank,":src_co2ice(1:2,1:3)=",src_field_2D(1:2,1:3)
376    CALL xios_send_field("co2ice_clean",src_field_2D(1:src_ni-1,1:src_nj))
377
378    ! emissivity:
379    CALL xios_recv_field("src_emis",src_field_2D)
380    WRITE(*,*) rank,":src_emis(1:2,1:3)=",src_field_2D(1:2,1:3)
381    CALL xios_send_field("emis_clean",src_field_2D(1:src_ni-1,1:src_nj))
382
383    ! q2 : q2surf for first layer of q2, the rest is q2atm
384    CALL xios_recv_field("src_q2surf",src_field_2D)
385    WRITE(*,*) rank,":src_q2(1:2,1:3)=",src_field_2D(1:2,1:3)
386    CALL xios_recv_field("src_q2atm",src_field_3D)
387    src_field_3D_p1(:,:,1)=src_field_2D(:,:)
388    src_field_3D_p1(:,:,2:src_nlev+1)=src_field_3D(:,:,:)
389
390    CALL xios_send_field("q2_clean",src_field_3D_p1(1:src_ni-1,1:src_nj,1:src_nlev+1))
391   
392    ! Temperature:
393    CALL xios_recv_field("src_temp",src_field_3D)
394    ! compute theta_rhodz
395    DO l=1,src_nlev
396      src_pk(:,:)=((.5/preff)*(src_pressure(:,:,l)+src_pressure(:,:,l+1)))**kappa
397      src_theta_rhodz(:,:,l) = src_field_3D(:,:,l) * &
398      ((src_pressure(:,:,l)-src_pressure(:,:,l+1))/gravity)/src_pk(:,:)
399    ENDDO
400    CALL xios_send_field("theta_rhodz_clean", &
401                         src_theta_rhodz(1:src_ni-1,1:src_nj,1:src_nlev))
402
403    ! zonal wind
404    CALL xios_recv_field("src_u",src_field_3D)
405    CALL xios_send_field("u_clean", &
406                          src_field_3D(1:src_ni-1,1:src_nj,1:src_nlev))
407   
408    ! meridional wind
409    CALL xios_recv_field("src_v",src_field_3D)
410    CALL xios_send_field("v_clean", &
411                          src_field_3D(1:src_ni-1,1:src_nj,1:src_nlev))
412    ! tracers
413    DO i=1,src_nq
414       if (i < 10) then
415           format_string = "(A8,I1)"
416       else
417           format_string = "(A8,I2)"
418       endif
419       write (i_trac, format_string) "src_trac", i
420       CALL xios_recv_field(i_trac,src_field_3D(:,:,:))
421       if (i==1) THEN
422          CALL xios_send_field("co2_clean", &
423                        src_field_3D(1:src_ni-1,1:src_nj,1:src_nlev))
424          CALL xios_recv_field("src_co2_surf",src_field_2D)
425          CALL xios_send_field("co2_surf_clean", &
426                        src_field_2D(1:src_ni-1,1:src_nj))
427       elseif (i==2) THEN
428          CALL xios_send_field("dust_number_clean", &
429                        src_field_3D(1:src_ni-1,1:src_nj,1:src_nlev))
430          CALL xios_recv_field("src_dust_number_surf",src_field_2D)
431          CALL xios_send_field("dust_number_surf_clean", &
432                        src_field_2D(1:src_ni-1,1:src_nj))
433       elseif (i==3) THEN
434          CALL xios_send_field("dust_mass_clean", &
435                        src_field_3D(1:src_ni-1,1:src_nj,1:src_nlev))
436          CALL xios_recv_field("src_dust_mass_surf",src_field_2D)
437          CALL xios_send_field("dust_mass_surf_clean", &
438                        src_field_2D(1:src_ni-1,1:src_nj))
439       elseif (i==4) THEN
440          CALL xios_send_field("ccn_number_clean", &
441                        src_field_3D(1:src_ni-1,1:src_nj,1:src_nlev))
442          CALL xios_recv_field("src_ccn_number_surf",src_field_2D)
443          CALL xios_send_field("ccn_number_surf_clean", &
444                        src_field_2D(1:src_ni-1,1:src_nj))
445       elseif (i==5) THEN
446          CALL xios_send_field("ccn_mass_clean", &
447                        src_field_3D(1:src_ni-1,1:src_nj,1:src_nlev))
448          CALL xios_recv_field("src_ccn_mass_surf",src_field_2D)
449          CALL xios_send_field("ccn_mass_surf_clean", &
450                        src_field_2D(1:src_ni-1,1:src_nj))
451       elseif (i==6) THEN
452          CALL xios_send_field("h2o_ice_clean", &
453                        src_field_3D(1:src_ni-1,1:src_nj,1:src_nlev))
454          CALL xios_recv_field("src_h2o_ice_surf",src_field_2D)
455          CALL xios_send_field("h2o_ice_surf_clean", &
456                        src_field_2D(1:src_ni-1,1:src_nj))
457       elseif (i==7) THEN
458          CALL xios_send_field("h2o_vap_clean", &
459                        src_field_3D(1:src_ni-1,1:src_nj,1:src_nlev))
460          CALL xios_recv_field("src_h2o_vap_surf",src_field_2D)
461          CALL xios_send_field("h2o_vap_surf_clean", &
462                        src_field_2D(1:src_ni-1,1:src_nj))
463       ENDIF
464    ENDDO
465
466    ! soil temperature
467    CALL xios_recv_field("src_tsoil",src_field_3D_soil)
468    CALL xios_send_field("tsoil_clean", &
469                          src_field_3D_soil(1:src_ni-1,1:src_nj,1:src_nsoil))
470    ! soil thermal intertia
471    CALL xios_recv_field("src_inertiedat",src_field_3D_soil)
472    CALL xios_send_field("inertiedat_clean", &
473                          src_field_3D_soil(1:src_ni-1,1:src_nj,1:src_nsoil))
474
475    ! ZMEA "zmea Orographie sous-maille"
476    CALL xios_recv_field("src_ZMEA",src_field_2D)
477    WRITE(*,*) rank,":src_ZMEA(1:2,1:3)=",src_field_2D(1:2,1:3)
478    CALL xios_send_field("ZMEA_clean",src_field_2D(1:src_ni-1,1:src_nj))
479   
480    ! ZSTD "zstd Orographie sous-maille"
481    CALL xios_recv_field("src_ZSTD",src_field_2D)
482    WRITE(*,*) rank,":src_ZSTD(1:2,1:3)=",src_field_2D(1:2,1:3)
483    CALL xios_send_field("ZSTD_clean",src_field_2D(1:src_ni-1,1:src_nj))
484   
485    ! ZSIG "zsig Orographie sous-maille"
486    CALL xios_recv_field("src_ZSIG",src_field_2D)
487    WRITE(*,*) rank,":src_ZSIG(1:2,1:3)=",src_field_2D(1:2,1:3)
488    CALL xios_send_field("ZSIG_clean",src_field_2D(1:src_ni-1,1:src_nj))
489   
490    ! ZGAM "zgam Orographie sous-maille"
491    CALL xios_recv_field("src_ZGAM",src_field_2D)
492    WRITE(*,*) rank,":src_ZGAM(1:2,1:3)=",src_field_2D(1:2,1:3)
493    CALL xios_send_field("ZGAM_clean",src_field_2D(1:src_ni-1,1:src_nj))
494   
495    ! ZTHE "zthe Orographie sous-maille"
496    CALL xios_recv_field("src_ZTHE",src_field_2D)
497    WRITE(*,*) rank,":src_ZTHE(1:2,1:3)=",src_field_2D(1:2,1:3)
498    CALL xios_send_field("ZTHE_clean",src_field_2D(1:src_ni-1,1:src_nj))
499
500    ! albedodat "albedo"
501    CALL xios_recv_field("src_albedodat",src_field_2D)
502    WRITE(*,*) rank,":src_albedodat(1:2,1:3)=",src_field_2D(1:2,1:3)
503    CALL xios_send_field("albedodat_clean",src_field_2D(1:src_ni-1,1:src_nj))
504
505    ! z0 "surface roughness"
506    CALL xios_recv_field("src_z0",src_field_2D)
507    WRITE(*,*) rank,":src_z0(1:2,1:3)=",src_field_2D(1:2,1:3)
508    CALL xios_send_field("z0_clean",src_field_2D(1:src_ni-1,1:src_nj))
509
510
511  ENDDO ! of DO ts=1,src_nt
512 
513!! Finalize
514  write(*,*) rank,":Finalize: call xios_context_finalize"
515  CALL xios_context_finalize()
516
517  write(*,*) rank,":Finalize: call MPI_COMM_FREE"
518  CALL MPI_COMM_FREE(comm, ierr)
519  write(*,*) rank,":Finalize: call xios_finalize"
520  CALL xios_finalize()
521
522  if (rank==0) then
523    ! add a couple of things in the "startphy_icosa.nc" file
524    write(*,*) rank,"Write controle() to startphy_icosa.nc"
525    ierr=NF90_OPEN(output_startfi_file,NF90_WRITE,ncid)
526    ierr=NF90_REDEF(ncid) ! switch to define mode
527    ierr=NF90_DEF_DIM(ncid,"index",100,dimids(1))
528    ierr=NF90_DEF_VAR(ncid,"controle",NF90_DOUBLE,dimids(1),varid)
529    ierr=NF90_ENDDEF(ncid) ! switch out of define mode
530    ierr=NF90_PUT_VAR(ncid,varid,src_controle(1:100))
531    if (ierr.ne.NF90_NOERR) then
532      write(*,*) "NetCDF Error:",NF90_STRERROR(ierr)
533    endif
534    ierr=NF90_CLOSE(ncid)
535    ! add a couple of things in the "start_icosa.nc" file
536    ierr=NF90_OPEN(output_start_file,NF90_WRITE,ncid)
537    ierr=NF90_REDEF(ncid)
538    ierr=NF90_DEF_DIM(ncid,"nvertex_u",2,dimids(1))
539    ierr=NF90_DEF_VAR(ncid,"iteration",NF90_FLOAT,varid)
540    ierr=NF90_ENDDEF(ncid)
541    if (ierr.ne.NF90_NOERR) then
542      write(*,*) "NetCDF Error:",NF90_STRERROR(ierr)
543    endif
544    ierr=NF90_PUT_VAR(ncid,varid,0) ! set "iteration" value to 0
545    ierr=NF90_CLOSE(ncid)
546   
547  endif ! of if (rank==0)
548
549  write(*,*) rank,":Finalize: call MPI_FINALIZE"
550  CALL MPI_FINALIZE(ierr)
551
552  write(*,*) rank,":my_remap: all is well that ends well!"
553
554END PROGRAM start_archive2icosa
Note: See TracBrowser for help on using the repository browser.