source: XIOS/dev/dev_ym/XIOS_COUPLING/src/test/generic_testcase.f90 @ 2239

Last change on this file since 2239 was 2239, checked in by ymipsl, 3 years ago

implement cyclic distribution for servers.
in param.def :
proc_distribution = contiguous/cyclic

YM

  • Property svn:executable set to *
File size: 89.7 KB
Line 
1PROGRAM generic_testcase
2  USE xios
3  USE mod_wait
4  IMPLICIT NONE
5  INCLUDE "mpif.h"
6  INTEGER,PARAMETER :: unit=10
7  INTEGER,PARAMETER :: len_str=255
8  INTEGER :: ierr
9
10  INTEGER :: nb_proc_atm
11  INTEGER :: nb_proc_oce
12  INTEGER :: nb_proc_surf
13  INTEGER :: nb_proc_server
14  CHARACTER(LEN=len_str)  :: proc_distribution
15  CHARACTER(LEN=len_str)  :: start_date
16  CHARACTER(LEN=len_str)  :: duration
17  DOUBLE PRECISION        :: sypd
18  NAMELIST /params_run/    nb_proc_atm, nb_proc_surf, nb_proc_oce, & 
19                           duration, sypd, start_date, proc_distribution
20
21
22  TYPE tmodel_params
23    CHARACTER(len_str)  :: timestep=""
24    CHARACTER(len_str)  :: domain=""
25    DOUBLE PRECISION    :: domain_proc_frac
26    INTEGER             :: domain_proc_n
27    CHARACTER(len_str)  :: axis=""
28    DOUBLE PRECISION    :: axis_proc_frac
29    INTEGER             :: axis_proc_n
30    INTEGER             :: ensemble_proc_n
31    INTEGER             :: ni 
32    INTEGER             :: nj 
33    INTEGER             :: nlev
34    CHARACTER(len_str)  :: init_field2D=""
35    DOUBLE PRECISION    :: pressure_factor
36    LOGICAL             :: domain_mask 
37    CHARACTER(len_str)  :: domain_mask_type=""
38    LOGICAL             :: scalar_mask
39    CHARACTER(len_str) :: scalar_mask_type=""
40    LOGICAL :: axis_mask
41    LOGICAL :: mask3d
42    INTEGER             :: field_sub_freq 
43    INTEGER             :: field_sub_offset
44  END TYPE  tmodel_params
45
46
47
48  CHARACTER(len_str)  :: timestep_atm
49  CHARACTER(len_str)  :: domain_atm 
50  CHARACTER(len_str)  :: axis_atm 
51  INTEGER             :: ni_atm 
52  INTEGER             :: nj_atm 
53  INTEGER             :: nlev_atm
54  CHARACTER(len_str)  :: init_field2D_atm 
55 
56 
57  LOGICAL :: i_am_atm, i_am_oce, i_am_surf, i_am_server
58  INTEGER :: rank, size_loc
59  INTEGER :: nb_procs(4)
60  LOGICAL :: who_i_am(4)
61   
62  OPEN(unit=unit, file='param.def',status='old',iostat=ierr)
63  nb_proc_atm=1
64  nb_proc_oce=0 
65  nb_proc_surf=0
66  duration="1d"
67  sypd=10000
68  start_date="2018-01-01"
69  proc_distribution="contiguous"
70  READ (unit, nml=params_run)
71  CLOSE(unit)
72!!! MPI Initialization
73
74  CALL MPI_INIT(ierr)
75  CALL init_wait
76  CALL MPI_COMM_RANK(MPI_COMM_WORLD,rank,ierr)
77  CALL MPI_COMM_SIZE(MPI_COMM_WORLD,size_loc,ierr)
78  if (size_loc < nb_proc_atm + nb_proc_oce + nb_proc_surf) then
79     PRINT *,"inconsistency : size=",size_loc," nb_proc_atm=",nb_proc_atm,&
80          " nb_proc_oce=",nb_proc_oce," nb_proc_surf=",nb_proc_surf
81     CALL MPI_ABORT()
82  endif
83
84  nb_proc_server=size_loc-(nb_proc_atm + nb_proc_oce + nb_proc_surf)
85  i_am_atm=.FALSE. ; i_am_oce=.FALSE. ; i_am_surf=.FALSE. ; i_am_server=.FALSE. 
86
87  IF (proc_distribution=="contiguous") THEN
88
89   ! First ranks are dedicated to Atm, next to Ocean, next to Surf, last to Server
90
91    IF (rank < nb_proc_atm) THEN
92      i_am_atm=.TRUE. 
93    ELSE IF (rank < nb_proc_atm+nb_proc_oce) THEN
94      i_am_oce=.TRUE. 
95    ELSE IF (rank < nb_proc_atm+nb_proc_oce+nb_proc_surf) THEN
96      i_am_surf=.TRUE. 
97    ELSE
98      i_am_server=.TRUE. 
99    ENDIF ;
100  ELSE IF (proc_distribution=="cyclic") THEN
101   ! server ranks are dispatched all along the partition
102    nb_procs(1)=nb_proc_atm ; nb_procs(2)=nb_proc_oce ; nb_procs(3)=nb_proc_surf ; nb_procs(4)=nb_proc_server
103    CALL compute_cyclic_distribution(rank, nb_procs, who_i_am )
104    i_am_atm=who_i_am(1)
105    i_am_oce=who_i_am(2)
106    i_am_surf=who_i_am(3)
107    i_am_server=who_i_am(4)
108  ENDIF
109 
110 
111  IF (i_am_server) THEN
112    CALL xios_init_server()
113 
114  ELSE
115
116    IF (i_am_atm) CALL model("atm")
117    IF (i_am_oce) CALL model("oce")
118    IF (i_am_surf) CALL model("surf")
119 
120    CALL xios_finalize()
121  ENDIF
122
123  CALL MPI_FINALIZE(ierr)
124 
125CONTAINS
126
127  SUBROUTINE model(model_id)
128  IMPLICIT NONE
129    CHARACTER(LEN=*) :: model_id
130    TYPE(tmodel_params) :: params, other_params
131    TYPE(xios_context) :: ctx_hdl
132    INTEGER :: comm
133    TYPE(xios_date)  :: cdate, edate
134    TYPE(xios_duration)  :: dtime
135    INTEGER :: rank
136    INTEGER :: size_loc
137    DOUBLE PRECISION :: timestep_in_seconds, simulated_seconds_per_seconds, elapsed_per_timestep
138    INTEGER :: ts
139
140    DOUBLE PRECISION, POINTER :: lon(:)
141    DOUBLE PRECISION, POINTER :: lat(:)
142    DOUBLE PRECISION, POINTER :: axis_value(:)
143    LOGICAL, POINTER          :: domain_mask(:)
144    INTEGER, POINTER          :: domain_index(:)
145    LOGICAL, POINTER          :: axis_mask(:)
146    LOGICAL                   :: scalar_mask
147    INTEGER, POINTER          :: axis_index(:)
148    DOUBLE PRECISION, POINTER :: ensemble_value(:)
149   
150    DOUBLE PRECISION, POINTER :: X_lon(:)
151    DOUBLE PRECISION, POINTER :: X_lat(:)
152    LOGICAL, POINTER          :: X_mask(:)
153    INTEGER, POINTER          :: X_index(:)
154
155    DOUBLE PRECISION, POINTER :: Y_lon(:)
156    DOUBLE PRECISION, POINTER :: Y_lat(:)
157    LOGICAL, POINTER          :: Y_mask(:)
158    INTEGER, POINTER          :: Y_index(:)
159               
160    DOUBLE PRECISION, POINTER :: field2D_init(:)
161    DOUBLE PRECISION, POINTER :: fieldX_init(:)
162    DOUBLE PRECISION, POINTER :: fieldY_init(:)
163    DOUBLE PRECISION, POINTER :: fieldXY_init(:,:)
164
165    DOUBLE PRECISION, POINTER :: field2D(:), other_field2D(:)
166    DOUBLE PRECISION          :: field0D, other_field0D
167    DOUBLE PRECISION, POINTER :: field_X(:), other_field_X(:)
168    DOUBLE PRECISION, POINTER :: field_Y(:), other_field_Y(:)
169    DOUBLE PRECISION, POINTER :: field_Z(:), other_field_Z(:)
170    DOUBLE PRECISION, POINTER :: field_XY(:,:), other_field_XY(:,:)
171    DOUBLE PRECISION, POINTER :: field_XYZ(:,:,:), other_field_XYZ(:,:,:)
172    DOUBLE PRECISION, POINTER :: field_XZ(:,:), other_field_XZ(:,:)
173    DOUBLE PRECISION, POINTER :: field_YZ(:,:), other_field_YZ(:,:)
174   
175    DOUBLE PRECISION, POINTER :: field2D_sub(:), other_field2D_sub(:)
176    DOUBLE PRECISION, POINTER :: field3D(:,:), other_field3D(:,:)
177    DOUBLE PRECISION, POINTER :: field3D_sub(:,:), other_field3D_sub(:,:)
178    DOUBLE PRECISION, POINTER :: field3D_recv(:,:), other_field3D_recv(:,:)
179    DOUBLE PRECISION, POINTER :: pressure(:,:), other_pressure(:,:)
180
181    DOUBLE PRECISION, POINTER :: field2D_W(:,:), other_field2D_W(:,:)
182    DOUBLE PRECISION, POINTER :: field0D_W(:), other_field0D_W(:)
183    DOUBLE PRECISION, POINTER :: field_XW(:,:), other_field_XW(:,:)
184    DOUBLE PRECISION, POINTER :: field_YW(:,:), other_field_YW(:,:)
185    DOUBLE PRECISION, POINTER :: field_ZW(:,:), other_field_ZW(:,:)
186    DOUBLE PRECISION, POINTER :: field_XYW(:,:,:), other_field_XYW(:,:,:)
187    DOUBLE PRECISION, POINTER :: field_XYZW(:,:,:,:), other_field_XYZW(:,:,:,:)
188    DOUBLE PRECISION, POINTER :: field_XZW(:,:,:), other_field_XZW(:,:,:)
189    DOUBLE PRECISION, POINTER :: field_YZW(:,:,:), other_field_YZW(:,:,:)
190   
191    DOUBLE PRECISION, POINTER :: field2D_sub_W(:,:), other_field2D_sub_W(:,:)
192    DOUBLE PRECISION, POINTER :: field3D_W(:,:,:), other_field3D_W(:,:,:)
193    DOUBLE PRECISION, POINTER :: field3D_sub_W(:,:,:), other_field3D_sub_W(:,:,:)
194    DOUBLE PRECISION, POINTER :: field3D_recv_W(:,:,:), other_field3D_recv_W(:,:,:)
195    DOUBLE PRECISION, POINTER :: pressure_W(:,:,:), other_pressure_W(:,:,:)
196
197
198
199   
200    TYPE(xios_duration) :: freq_op_recv, freq_offset_recv
201    INTEGER :: ni,nj,nk
202    INTEGER :: i,j,k,xy,x,y,z,w
203    DOUBLE PRECISION :: scale,dist
204    LOGICAL :: ok
205    INTEGER :: ierr     
206
207    LOGICAL :: ok_field0D, ok_field2D, ok_field3D, ok_pressure, ok_field2D_sub, ok_field3D_sub,ok_field3D_recv, ok_field3D_send
208    LOGICAL :: ok_field_X, ok_field_Y, ok_field_XY, ok_field_Z, ok_field_XYZ, ok_field_XZ, ok_field_YZ
209    LOGICAL :: ok_field0D_W, ok_field2D_W, ok_field3D_W, ok_pressure_W, ok_field2D_sub_W, ok_field3D_sub_W,ok_field3D_recv_W, ok_field3D_send_W
210    LOGICAL :: ok_field_XW, ok_field_YW, ok_field_XYW, ok_field_ZW, ok_field_XYZW, ok_field_XZW, ok_field_YZW
211
212    LOGICAL :: ok_other_field0D, ok_other_field2D, ok_other_field3D, ok_other_pressure, ok_other_field2D_sub, ok_other_field3D_sub,ok_other_field3D_recv, ok_other_field3D_send
213    LOGICAL :: ok_other_field_X, ok_other_field_Y, ok_other_field_XY, ok_other_field_Z, ok_other_field_XYZ, ok_other_field_XZ, ok_other_field_YZ
214    LOGICAL :: ok_other_field0D_W, ok_other_field2D_W, ok_other_field3D_W, ok_other_pressure_W, ok_other_field2D_sub_W, ok_other_field3D_sub_W,ok_other_field3D_recv_W, ok_other_field3D_send_W
215    LOGICAL :: ok_other_field_XW, ok_other_field_YW, ok_other_field_XYW, ok_other_field_ZW, ok_other_field_XYZW, ok_other_field_XZW, ok_other_field_YZW
216
217    TYPE(xios_domain)  :: domain_handle
218   
219      !! XIOS Initialization (get the local communicator)
220    CALL xios_initialize(trim(model_id),return_comm=comm)
221    CALL MPI_COMM_RANK(comm,rank,ierr)
222    CALL MPI_COMM_SIZE(comm,size_loc,ierr)
223    CALL xios_context_initialize(trim(model_id),comm)
224    CALL xios_get_handle(trim(model_id),ctx_hdl)
225    CALL xios_set_current_context(ctx_hdl)
226   
227    CALL init_model_params('',params)
228    CALL init_model_params('other_',other_params)
229   !!! Definition du découpage domain/axis
230
231   
232   !!! Definition de la start date et du timestep
233    CALL xios_set_start_date(xios_date_convert_from_string(start_date//" 00:00:00"))
234    dtime=xios_duration_convert_from_string(TRIM(params%timestep))
235    CALL xios_set_timestep(timestep=dtime)
236
237     !!! Calcul de temps elaps par seconde pour respecter le SYPD (hyp : pas de delai d'I/O)
238    CALL xios_get_start_date(cdate)
239    edate=cdate+xios_duration_convert_from_string(TRIM(duration))
240    timestep_in_seconds=xios_date_convert_to_seconds(cdate+dtime) - xios_date_convert_to_seconds(cdate)
241    simulated_seconds_per_seconds=sypd * 365 
242    elapsed_per_timestep=timestep_in_seconds/simulated_seconds_per_seconds ! in seconds
243
244
245!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
246!!!             standard initialisation of domain, axis, grid, field                               !
247!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
248
249    CALL init_domain("domain", comm, params, ni, nj, lon, lat, domain_mask, domain_index, &
250                                                     X_lon, X_lat, X_mask, X_index,       &
251                                                     Y_lon, Y_lat, Y_mask, Y_index)
252                                                     
253    CALL init_axis("axis", comm, params, axis_value, axis_mask, axis_index)
254    CALL init_scalar("scalar", comm, params,  scalar_mask)
255    CALL init_ensemble("ensemble", comm, params, ensemble_value)
256
257    CALL set_mask3d("grid3D",params, ni, nj, lon, lat , axis_value)
258     !!! Fin de la definition du contexte
259
260    ok_field3D_recv=xios_is_valid_field("field3D_recv").AND.xios_is_valid_field("field3D_resend") ;
261    IF (ok_field3D_recv) THEN
262      CALL xios_is_defined_field_attr("field3D_recv",freq_op=ok)
263      IF (ok) THEN
264        CALL xios_get_field_attr("field3D_recv",freq_op=freq_op_recv)
265       ELSE
266         freq_op_recv%timestep=1
267         CALL xios_set_field_attr("field3D_recv",freq_op=freq_op_recv)
268      ENDIF
269      CALL xios_is_defined_field_attr("field3D_recv",freq_offset=ok)
270      IF (ok) THEN
271        CALL xios_get_field_attr("field3D_recv",freq_offset=freq_offset_recv)
272      ELSE
273         freq_offset_recv%timestep=0
274         CALL xios_set_field_attr("field3D_recv",freq_op=freq_op_recv)
275      ENDIF
276      CALL xios_set_field_attr("field3D_resend",freq_op=freq_op_recv,freq_offset=freq_offset_recv)
277    ENDIF
278   
279    ok_field3D_recv_W=xios_is_valid_field("field3D_recv_W").AND.xios_is_valid_field("field3D_resend_W") ;
280    IF (ok_field3D_recv_W) THEN
281      CALL xios_is_defined_field_attr("field3D_recv_W",freq_op=ok)
282      IF (ok) THEN
283        CALL xios_get_field_attr("field3D_recv_W",freq_op=freq_op_recv)
284       ELSE
285         freq_op_recv%timestep=1
286         CALL xios_set_field_attr("field3D_recv_W",freq_op=freq_op_recv)
287      ENDIF
288      CALL xios_is_defined_field_attr("field3D_recv_W",freq_offset=ok)
289      IF (ok) THEN
290        CALL xios_get_field_attr("field3D_recv_W",freq_offset=freq_offset_recv)
291      ELSE
292         freq_offset_recv%timestep=0
293         CALL xios_set_field_attr("field3D_recv_W",freq_op=freq_op_recv)
294      ENDIF
295      CALL xios_set_field_attr("field3D_resend_W",freq_op=freq_op_recv,freq_offset=freq_offset_recv)
296    ENDIF
297
298
299    CALL init_field2D(comm, params, lon, lat, domain_mask, field2D_init, &
300                      x_lon, x_lat, x_mask, fieldX_init, y_lon, y_lat, y_mask, fieldY_init, fieldXY_init)
301
302    xy=size(domain_index)
303    x=size(X_index)
304    y=size(Y_index)
305    z=size(axis_index)
306    w=size(ensemble_value)
307
308    ALLOCATE(field2D(0:xy-1))
309    ALLOCATE(field3D(0:xy-1,0:z-1))
310    ALLOCATE(pressure(0:xy-1,0:z-1))
311    ALLOCATE(field3D_recv(0:xy-1,0:z-1))
312    ALLOCATE(field_Z(0:z-1))
313    ALLOCATE(field_X(0:x-1))
314    ALLOCATE(field_Y(0:y-1))
315    ALLOCATE(field_XY(0:x-1,0:y-1))
316    ALLOCATE(field_XYZ(0:x-1,0:y-1,0:z-1))
317    ALLOCATE(field_XZ(0:x-1,0:z-1))
318    ALLOCATE(field_YZ(0:y-1,0:z-1))
319
320    ALLOCATE(field2D_W(0:xy-1,0:w-1))
321    ALLOCATE(field3D_W(0:xy-1,0:z-1,0:w-1))
322    ALLOCATE(pressure_W(0:xy-1,0:z-1,0:w-1))
323    ALLOCATE(field3D_recv_W(0:xy-1,0:z-1,0:w-1))
324    ALLOCATE(field_ZW(0:z-1,0:w-1))
325    ALLOCATE(field_XW(0:x-1,0:w-1))
326    ALLOCATE(field_YW(0:y-1,0:w-1))
327    ALLOCATE(field_XYW(0:x-1,0:y-1,0:w-1))
328    ALLOCATE(field_XYZW(0:x-1,0:y-1,0:z-1,0:w-1))
329    ALLOCATE(field_XZW(0:x-1,0:z-1,0:w-1))
330    ALLOCATE(field_YZW(0:y-1,0:z-1,0:w-1))
331    ALLOCATE(field0D_W(0:w-1))
332   
333
334
335    DO i=0,xy-1
336      IF (domain_index(i)/=-1) THEN
337        field2D(i)=field2D_init(domain_index(i))
338      ENDIF
339    ENDDO
340
341    DO i=0,x-1
342      IF (X_index(i)/=-1) THEN
343        field_X(i)=fieldX_init(X_index(i))
344      ENDIF
345    ENDDO
346
347    DO i=0,y-1
348      IF (Y_index(i)/=-1) THEN
349        field_Y(i)=fieldY_init(Y_index(i))
350      ENDIF
351    ENDDO
352
353    DO j=0,y-1
354      DO i=0,x-1
355        IF (Y_index(j)/=-1 .AND. X_index(i)/=-1) THEN
356          field_XY(i,j)=fieldXY_init(X_index(i),Y_index(j))
357        ENDIF
358      ENDDO
359    ENDDO
360
361    DO i=0,z-1
362      IF (axis_index(i)/=-1) THEN
363        k=axis_index(i)
364        IF (axis_mask(k)) THEN
365          field_Z(i)=axis_value(axis_index(k))
366          field3D(:,i)=field2D(:)*axis_value(axis_index(k))
367          field_XYZ(:,:,i)=field_XY(:,:)*axis_value(axis_index(k))
368          field_XZ(:,i)=field_X(:)*axis_value(axis_index(k))
369          field_YZ(:,i)=field_Y(:)*axis_value(axis_index(k))
370        ENDIF
371      ENDIF
372    ENDDO
373         
374    pressure=1e20
375    DO j=0,z-1
376      pressure(:,j)=axis_value(j)*100000 ; ! Pa
377      DO i=0,xy-1
378        IF (domain_index(i)/=-1) THEN
379          k=domain_index(i)
380          IF (domain_mask(k)) THEN
381            dist=sqrt((lat(k)/90.)**2+(lon(k)/180.)**2) ;
382            pressure(i,j)=pressure(i,j)*(1+params%pressure_factor*exp(-(4*dist)**2))
383          ENDIF
384        ENDIF
385      ENDDO
386    ENDDO
387   
388    field0D=1
389
390
391    field2D_W(:,0) = field2D(:)*(1+0.1*ensemble_value(0))
392    field3D_W(:,:,0)= field3D(:,:)*(1+0.1*ensemble_value(0))
393    pressure_W(:,:,0) = pressure(:,:)*(1+0.1*ensemble_value(0))
394    field_XW(:,0) = field_X(:)*(1+0.1*ensemble_value(0))
395    field_YW(:,0) = field_Y(:)*(1+0.1*ensemble_value(0))
396    field_XYW(:,:,0) = field_XY(:,:)*(1+0.1*ensemble_value(0))
397    field_ZW(:,0) = field_Z(:)*(1+0.1*ensemble_value(0))
398    field_XYZW(:,:,:,0) = field_XYZ(:,:,:)*(1+0.1*ensemble_value(0))
399    field_XZW(:,:,0) = field_XZ(:,:)*(1+0.1*ensemble_value(0))
400    field_YZW(:,:,0) = field_YZ(:,:)*(1+0.1*ensemble_value(0))
401    field0D_W(0) = field0D*(1+0.1*ensemble_value(0))
402   
403    ok_field2D=xios_is_valid_field("field2D") ;
404    ok_field3D=xios_is_valid_field("field3D") ;
405    ok_pressure=xios_is_valid_field("pressure") ;
406    ok_field2D_sub=xios_is_valid_field("field2D_sub") ;
407    ok_field3D_sub=xios_is_valid_field("field3D_sub") ;
408    ok_field_X=xios_is_valid_field("field_X") ;
409    ok_field_Y=xios_is_valid_field("field_Y") ;
410    ok_field_XY=xios_is_valid_field("field_XY") ;
411    ok_field_Z=xios_is_valid_field("field_Z") ; 
412    ok_field_XYZ=xios_is_valid_field("field_XYZ") ;
413    ok_field_XZ=xios_is_valid_field("field_XZ") ;
414    ok_field_YZ=xios_is_valid_field("field_YZ") ;
415    ok_field0D=xios_is_valid_field("field0D") ;
416
417    ok_field2D_W=xios_is_valid_field("field2D_W") ;
418    ok_field3D_W=xios_is_valid_field("field3D_W") ;
419    ok_pressure_W=xios_is_valid_field("pressure_W") ;
420    ok_field2D_sub_W=xios_is_valid_field("field2D_sub_W") ;
421    ok_field3D_sub_W=xios_is_valid_field("field3D_sub_W") ;
422    ok_field_XW=xios_is_valid_field("field_XW") ;
423    ok_field_YW=xios_is_valid_field("field_YW") ;
424    ok_field_XYW=xios_is_valid_field("field_XYW") ;
425    ok_field_ZW=xios_is_valid_field("field_ZW") ; 
426    ok_field_XYZW=xios_is_valid_field("field_XYZW") ;
427    ok_field_XZW=xios_is_valid_field("field_XZW") ;
428    ok_field_YZW=xios_is_valid_field("field_YZW") ;
429    ok_field0D_W=xios_is_valid_field("field0D_W") ;
430
431
432
433    IF (ASSOCIATED(lon)) DEALLOCATE(lon)
434    IF (ASSOCIATED(lat)) DEALLOCATE(lat)
435    IF (ASSOCIATED(axis_value)) DEALLOCATE(axis_value)
436    IF (ASSOCIATED(domain_mask)) DEALLOCATE(domain_mask)
437    IF (ASSOCIATED(domain_index)) DEALLOCATE(domain_index)
438    IF (ASSOCIATED(axis_mask)) DEALLOCATE(axis_mask)
439    IF (ASSOCIATED(axis_index)) DEALLOCATE(axis_index)
440    IF (ASSOCIATED(ensemble_value)) DEALLOCATE(ensemble_value)
441    IF (ASSOCIATED(X_lon)) DEALLOCATE(X_lon)
442    IF (ASSOCIATED(X_lat)) DEALLOCATE(X_lat)
443    IF (ASSOCIATED(X_mask)) DEALLOCATE(X_mask)
444    IF (ASSOCIATED(X_index)) DEALLOCATE(X_index)
445    IF (ASSOCIATED(Y_lon)) DEALLOCATE(Y_lon)
446    IF (ASSOCIATED(Y_lat)) DEALLOCATE(Y_lat)
447    IF (ASSOCIATED(Y_mask)) DEALLOCATE(Y_mask)
448    IF (ASSOCIATED(Y_index)) DEALLOCATE(Y_index)
449    IF (ASSOCIATED(field2D_init)) DEALLOCATE(field2D_init)
450    IF (ASSOCIATED(fieldX_init)) DEALLOCATE(fieldX_init)
451    IF (ASSOCIATED(fieldY_init)) DEALLOCATE(fieldY_init)
452    IF (ASSOCIATED(fieldXY_init)) DEALLOCATE(fieldXY_init)
453
454
455!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
456!              duplicated section for other_                    !
457!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
458
459    CALL init_domain("other_domain", comm, other_params, ni, nj, lon, lat, domain_mask, domain_index, &
460                                                     X_lon, X_lat, X_mask, X_index,       &
461                                                     Y_lon, Y_lat, Y_mask, Y_index)
462                                                     
463    CALL init_axis("other_axis", comm, other_params, axis_value, axis_mask, axis_index)
464    CALL init_scalar("other_scalar", comm, params, scalar_mask)
465    CALL init_ensemble("other_ensemble", comm, other_params, ensemble_value)
466
467    CALL set_mask3d("other_grid3D",other_params, ni, nj, lon, lat , axis_value)
468     !!! Fin de la definition du contexte
469
470    ok_other_field3D_recv=xios_is_valid_field("other_field3D_recv").AND.xios_is_valid_field("other_field3D_resend") ;
471    IF (ok_other_field3D_recv) THEN
472      CALL xios_is_defined_field_attr("other_field3D_recv",freq_op=ok)
473      IF (ok) THEN
474        CALL xios_get_field_attr("other_field3D_recv",freq_op=freq_op_recv)
475       ELSE
476         freq_op_recv%timestep=1
477         CALL xios_set_field_attr("other_field3D_recv",freq_op=freq_op_recv)
478      ENDIF
479      CALL xios_is_defined_field_attr("other_field3D_recv",freq_offset=ok)
480      IF (ok) THEN
481        CALL xios_get_field_attr("other_field3D_recv",freq_offset=freq_offset_recv)
482      ELSE
483         freq_offset_recv%timestep=0
484         CALL xios_set_field_attr("other_field3D_recv",freq_op=freq_op_recv)
485      ENDIF
486      CALL xios_set_field_attr("other_field3D_resend",freq_op=freq_op_recv,freq_offset=freq_offset_recv)
487    ENDIF
488   
489    ok_other_field3D_recv_W=xios_is_valid_field("other_field3D_recv_W").AND.xios_is_valid_field("other_field3D_resend_W") ;
490    IF (ok_other_field3D_recv_W) THEN
491      CALL xios_is_defined_field_attr("other_field3D_recv_W",freq_op=ok)
492      IF (ok) THEN
493        CALL xios_get_field_attr("other_field3D_recv_W",freq_op=freq_op_recv)
494       ELSE
495         freq_op_recv%timestep=1
496         CALL xios_set_field_attr("other_field3D_recv_W",freq_op=freq_op_recv)
497      ENDIF
498      CALL xios_is_defined_field_attr("other_field3D_recv_W",freq_offset=ok)
499      IF (ok) THEN
500        CALL xios_get_field_attr("other_field3D_recv_W",freq_offset=freq_offset_recv)
501      ELSE
502         freq_offset_recv%timestep=0
503         CALL xios_set_field_attr("other_field3D_recv_W",freq_op=freq_op_recv)
504      ENDIF
505      CALL xios_set_field_attr("other_field3D_resend_W",freq_op=freq_op_recv,freq_offset=freq_offset_recv)
506    ENDIF
507
508
509    CALL init_field2D(comm, other_params, lon, lat, domain_mask, field2D_init, &
510                      x_lon, x_lat, x_mask, fieldX_init, y_lon, y_lat, y_mask, fieldY_init, fieldXY_init)
511
512    xy=size(domain_index)
513    x=size(X_index)
514    y=size(Y_index)
515    z=size(axis_index)
516    w=size(ensemble_value)
517
518    ALLOCATE(other_field2D(0:xy-1))
519    ALLOCATE(other_field3D(0:xy-1,0:z-1))
520    ALLOCATE(other_pressure(0:xy-1,0:z-1))
521    ALLOCATE(other_field3D_recv(0:xy-1,0:z-1))
522    ALLOCATE(other_field_Z(0:z-1))
523    ALLOCATE(other_field_X(0:x-1))
524    ALLOCATE(other_field_Y(0:y-1))
525    ALLOCATE(other_field_XY(0:x-1,0:y-1))
526    ALLOCATE(other_field_XYZ(0:x-1,0:y-1,0:z-1))
527    ALLOCATE(other_field_XZ(0:x-1,0:z-1))
528    ALLOCATE(other_field_YZ(0:y-1,0:z-1))
529
530    ALLOCATE(other_field2D_W(0:xy-1,0:w-1))
531    ALLOCATE(other_field3D_W(0:xy-1,0:z-1,0:w-1))
532    ALLOCATE(other_pressure_W(0:xy-1,0:z-1,0:w-1))
533    ALLOCATE(other_field3D_recv_W(0:xy-1,0:z-1,0:w-1))
534    ALLOCATE(other_field_ZW(0:z-1,0:w-1))
535    ALLOCATE(other_field_XW(0:x-1,0:w-1))
536    ALLOCATE(other_field_YW(0:y-1,0:w-1))
537    ALLOCATE(other_field_XYW(0:x-1,0:y-1,0:w-1))
538    ALLOCATE(other_field_XYZW(0:x-1,0:y-1,0:z-1,0:w-1))
539    ALLOCATE(other_field_XZW(0:x-1,0:z-1,0:w-1))
540    ALLOCATE(other_field_YZW(0:y-1,0:z-1,0:w-1))
541    ALLOCATE(other_field0D_W(0:w-1))
542   
543
544
545    DO i=0,xy-1
546      IF (domain_index(i)/=-1) THEN
547        other_field2D(i)=field2D_init(domain_index(i))
548      ENDIF
549    ENDDO
550
551    DO i=0,x-1
552      IF (X_index(i)/=-1) THEN
553        other_field_X(i)=fieldX_init(X_index(i))
554      ENDIF
555    ENDDO
556
557    DO i=0,y-1
558      IF (Y_index(i)/=-1) THEN
559        other_field_Y(i)=fieldY_init(Y_index(i))
560      ENDIF
561    ENDDO
562
563    DO j=0,y-1
564      DO i=0,x-1
565        IF (Y_index(j)/=-1 .AND. X_index(i)/=-1) THEN
566          other_field_XY(i,j)=fieldXY_init(X_index(i),Y_index(j))
567        ENDIF
568      ENDDO
569    ENDDO
570
571    DO i=0,z-1
572      IF (axis_index(i)/=-1) THEN
573        k=axis_index(i)
574        IF (axis_mask(k)) THEN
575          other_field_Z(i)=axis_value(axis_index(k))
576          other_field3D(:,i)=other_field2D(:)*axis_value(axis_index(k))
577          other_field_XYZ(:,:,i)=other_field_XY(:,:)*axis_value(axis_index(k))
578          other_field_XZ(:,i)=other_field_X(:)*axis_value(axis_index(k))
579          other_field_YZ(:,i)=other_field_Y(:)*axis_value(axis_index(k))
580        ENDIF
581      ENDIF
582    ENDDO
583         
584    other_pressure=1e20
585    DO j=0,z-1
586      other_pressure(:,j)=axis_value(j)*100000 ; ! Pa
587      DO i=0,xy-1
588        IF (domain_index(i)/=-1) THEN
589          k=domain_index(i)
590          IF (domain_mask(k)) THEN
591            dist=sqrt((lat(k)/90.)**2+(lon(k)/180.)**2) ;
592            other_pressure(i,j)=other_pressure(i,j)*(1+other_params%pressure_factor*exp(-(4*dist)**2))
593          ENDIF
594        ENDIF
595      ENDDO
596    ENDDO
597   
598    other_field0D = 1
599
600
601    other_field2D_W(:,0) = other_field2D(:)*(1+0.1*ensemble_value(0))
602    other_field3D_W(:,:,0)= other_field3D(:,:)*(1+0.1*ensemble_value(0))
603    other_pressure_W(:,:,0) = other_pressure(:,:)*(1+0.1*ensemble_value(0))
604    other_field_XW(:,0) = other_field_X(:)*(1+0.1*ensemble_value(0))
605    other_field_YW(:,0) = other_field_Y(:)*(1+0.1*ensemble_value(0))
606    other_field_XYW(:,:,0) = other_field_XY(:,:)*(1+0.1*ensemble_value(0))
607    other_field_ZW(:,0) = other_field_Z(:)*(1+0.1*ensemble_value(0))
608    other_field_XYZW(:,:,:,0) = other_field_XYZ(:,:,:)*(1+0.1*ensemble_value(0))
609    other_field_XZW(:,:,0) = other_field_XZ(:,:)*(1+0.1*ensemble_value(0))
610    other_field_YZW(:,:,0) = other_field_YZ(:,:)*(1+0.1*ensemble_value(0))
611    other_field0D_W(0) = other_field0D*(1+0.1*ensemble_value(0))
612   
613   
614    ok_other_field2D=xios_is_valid_field("other_field2D") ;
615    ok_other_field3D=xios_is_valid_field("other_field3D") ;
616    ok_other_pressure=xios_is_valid_field("other_pressure") ;
617    ok_other_field2D_sub=xios_is_valid_field("other_field2D_sub") ;
618    ok_other_field3D_sub=xios_is_valid_field("other_field3D_sub") ;
619    ok_other_field_X=xios_is_valid_field("other_field_X") ;
620    ok_other_field_Y=xios_is_valid_field("other_field_Y") ;
621    ok_other_field_XY=xios_is_valid_field("other_field_XY") ;
622    ok_other_field_Z=xios_is_valid_field("other_field_Z") ; 
623    ok_other_field_XYZ=xios_is_valid_field("other_field_XYZ") ;
624    ok_other_field_XZ=xios_is_valid_field("other_field_XZ") ;
625    ok_other_field_YZ=xios_is_valid_field("other_field_YZ") ;
626    ok_other_field0D=xios_is_valid_field("other_field0D") ;
627
628    ok_other_field2D_W=xios_is_valid_field("other_field2D_W") ;
629    ok_other_field3D_W=xios_is_valid_field("other_field3D_W") ;
630    ok_other_pressure_W=xios_is_valid_field("other_pressure_W") ;
631    ok_other_field2D_sub_W=xios_is_valid_field("other_field2D_sub_W") ;
632    ok_other_field3D_sub_W=xios_is_valid_field("other_field3D_sub_W") ;
633    ok_other_field_XW=xios_is_valid_field("other_field_XW") ;
634    ok_other_field_YW=xios_is_valid_field("other_field_YW") ;
635    ok_other_field_XYZW=xios_is_valid_field("other_field_XYW") ;
636    ok_other_field_ZW=xios_is_valid_field("other_field_ZW") ; 
637    ok_other_field_XYZW=xios_is_valid_field("other_field_XYZW") ;
638    ok_other_field_XZW=xios_is_valid_field("other_field_XZW") ;
639    ok_other_field_YZW=xios_is_valid_field("other_field_YZW") ;
640    ok_other_field0D_W=xios_is_valid_field("other_field0D_W") ;
641
642
643!!!!!!!!!!!!!!!!!!!!! end of duplicated section   !!!!!!!!!!!!!!!!!!!!!!!!!!
644
645    ts=1
646    cdate=cdate+dtime
647    CALL xios_close_context_definition()
648    CALL xios_set_current_context(trim(model_id))
649!    CALL xios_get_handle("field3D::domain",domain_handle)
650   
651    DO while ( cdate <= edate )
652 
653      !!! Mise a jour du pas de temps
654      CALL xios_update_calendar(ts)
655
656      IF (ok_field2D) CALL xios_send_field("field2D",field2D)
657      IF (ok_field3D) CALL xios_send_field("field3D",field3D)
658      IF (ok_pressure) CALL xios_send_field("pressure",pressure)
659      IF (ok_field_X) CALL xios_send_field("field_X",field_X)
660      IF (ok_field_Y) CALL xios_send_field("field_Y",field_Y)
661      IF (ok_field_XY) CALL xios_send_field("field_XY",field_XY)
662      IF (ok_field_Z) CALL xios_send_field("field_Z",field_Z)
663      IF (ok_field_XYZ) CALL xios_send_field("field_XYZ",field_XYZ)
664      IF (ok_field_XZ) CALL xios_send_field("field_XZ",field_XZ)
665      IF (ok_field_YZ) CALL xios_send_field("field_YZ",field_YZ)
666      IF (ok_field0D)  CALL xios_send_field("field0D",field0D)
667     
668      IF ( MOD(params%field_sub_offset+ts-1,params%field_sub_freq)==0) THEN
669        IF (ok_field2D_sub) CALL xios_send_field("field2D_sub",field2D*10)
670        IF (ok_field3D_sub) CALL xios_send_field("field3D_sub",field3D*10)
671      ENDIF
672
673      IF (ok_field3D_recv) THEN
674        IF (xios_field_is_active("field3D_resend",.TRUE.)) THEN
675          CALL xios_recv_field("field3D_recv",field3D_recv)
676          CALL xios_send_field("field3D_resend",field3D_recv)         
677         ENDIF
678      ENDIF
679
680! ensemble
681      IF (ok_field2D_W) CALL xios_send_field("field2D_W",field2D_W)
682      IF (ok_field3D_W) CALL xios_send_field("field3D_W",field3D_W)
683      IF (ok_pressure_W) CALL xios_send_field("pressure_W",pressure_W)
684      IF (ok_field_XW) CALL xios_send_field("field_XW",field_XW)
685      IF (ok_field_YW) CALL xios_send_field("field_YW",field_YW)
686      IF (ok_field_XYW) CALL xios_send_field("field_XYW",field_XYW)
687      IF (ok_field_ZW) CALL xios_send_field("field_ZW",field_ZW)
688      IF (ok_field_XYZW) CALL xios_send_field("field_XYZW",field_XYZW)
689      IF (ok_field_XZW) CALL xios_send_field("field_XZW",field_XZW)
690      IF (ok_field_YZW) CALL xios_send_field("field_YZW",field_YZW)
691      IF (ok_field0D_W) CALL xios_send_field("field0D_W",field0D_W)
692     
693      IF ( MOD(params%field_sub_offset+ts-1,params%field_sub_freq)==0) THEN
694        IF (ok_field2D_sub_W) CALL xios_send_field("field2D_sub_W",field2D_W*10)
695        IF (ok_field3D_sub_W) CALL xios_send_field("field3D_sub_W",field3D_W*10)
696      ENDIF
697
698      IF (ok_field3D_recv_W) THEN
699        IF (xios_field_is_active("field3D_resend_W",.TRUE.)) THEN
700          CALL xios_recv_field("field3D_recv_W",field3D_recv_W)
701          CALL xios_send_field("field3D_resend_W",field3D_recv_W)         
702         ENDIF
703      ENDIF
704
705      field2D=field2D+1
706      field3D=field3D+1
707      field0D=field0D+1
708
709!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
710!              duplicated section for other_                    !
711!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
712
713
714      IF (ok_other_field2D) CALL xios_send_field("other_field2D", other_field2D)
715      IF (ok_other_field3D) CALL xios_send_field("other_field3D", other_field3D)
716      IF (ok_other_pressure) CALL xios_send_field("other_pressure", other_pressure)
717      IF (ok_other_field_X) CALL xios_send_field("other_field_X", other_field_X)
718      IF (ok_other_field_Y) CALL xios_send_field("other_field_Y", other_field_Y)
719      IF (ok_other_field_XY) CALL xios_send_field("other_field_XY", other_field_XY)
720      IF (ok_other_field_Z) CALL xios_send_field("other_field_Z", other_field_Z)
721      IF (ok_other_field_XY) CALL xios_send_field("other_field_XYZ", other_field_XYZ)
722      IF (ok_other_field_XY) CALL xios_send_field("other_field_XZ", other_field_XZ)
723      IF (ok_other_field_XY) CALL xios_send_field("other_field_YZ", other_field_YZ)
724      IF (ok_other_field0D)  CALL xios_send_field("other_field0D", other_field0D)
725     
726      IF ( MOD(other_params%field_sub_offset+ts-1,other_params%field_sub_freq)==0) THEN
727        IF (ok_other_field2D_sub) CALL xios_send_field("other_field2D_sub",other_field2D*10)
728        IF (ok_other_field3D_sub) CALL xios_send_field("other_field3D_sub",other_field3D*10)
729      ENDIF
730
731      IF (ok_other_field3D_recv) THEN
732        IF (xios_field_is_active("other_field3D_resend",.TRUE.)) THEN
733          CALL xios_recv_field("other_field3D_recv",other_field3D_recv)
734          CALL xios_send_field("other_field3D_resend",other_field3D_recv)         
735         ENDIF
736      ENDIF
737
738! ensemble
739      IF (ok_other_field2D_W) CALL xios_send_field("other_field2D_W",other_field2D_W)
740      IF (ok_other_field3D_W) CALL xios_send_field("other_field3D_W",other_field3D_W)
741      IF (ok_other_pressure_W) CALL xios_send_field("other_pressure_W",other_pressure_W)
742      IF (ok_other_field_XW)  CALL xios_send_field("other_field_XW",other_field_XW)
743      IF (ok_other_field_YW)  CALL xios_send_field("other_field_YW",other_field_YW)
744      IF (ok_other_field_XYW) CALL xios_send_field("other_field_XYW",other_field_XYW)
745      IF (ok_other_field_ZW)  CALL xios_send_field("other_field_ZW",other_field_ZW)
746      IF (ok_other_field_XYW) CALL xios_send_field("other_field_XYZW",other_field_XYZW)
747      IF (ok_other_field_XYW) CALL xios_send_field("other_field_XZW",other_field_XZW)
748      IF (ok_other_field_XYW) CALL xios_send_field("other_field_YZW",other_field_YZW)
749      IF (ok_other_field0D_W)   CALL xios_send_field("other_field0D_W",other_field0D_W)
750     
751      IF ( MOD(other_params%field_sub_offset+ts-1,other_params%field_sub_freq)==0) THEN
752        IF (ok_other_field2D_sub_W) CALL xios_send_field("other_field2D_sub_W",other_field2D_W*10)
753        IF (ok_other_field3D_sub_W) CALL xios_send_field("other_field3D_sub_W",other_field3D_W*10)
754      ENDIF
755
756      IF (ok_other_field3D_recv_W) THEN
757        IF (xios_field_is_active("other_field3D_resend_W",.TRUE.)) THEN
758          CALL xios_recv_field("other_field3D_recv_W",other_field3D_recv_W)
759          CALL xios_send_field("other_field3D_resend_W",other_field3D_recv_W)         
760         ENDIF
761      ENDIF
762
763      other_field2D=other_field2D+1
764      other_field3D=other_field3D+1
765      other_field0D=other_field0D+1
766
767
768!!!!!!!!!!!!!!!!!!!!!! end of duplicated section !!!!!!!!!!!!!!!!!
769
770
771
772      !! Boucle d'attente
773      CALL wait_us(int(elapsed_per_timestep*1.e6))   ! micro-secondes
774      cdate=cdate+dtime
775      ts=ts+1
776    ENDDO
777   
778    CALL xios_context_finalize()
779    CALL MPI_COMM_FREE(comm, ierr)
780   
781  END SUBROUTINE model
782
783  SUBROUTINE init_model_params(prefix,params)
784  IMPLICIT NONE
785    CHARACTER(LEN=*) :: prefix
786    TYPE(tmodel_params) :: params
787
788
789    IF (.NOT. xios_getvar(prefix//"timestep", params%timestep) )         params%timestep="1h"
790    IF (.NOT. xios_getvar(prefix//"domain", params%domain) )             params%domain="lmdz"
791    IF (.NOT. xios_getvar(prefix//"domain_mask", params%domain_mask) )   params%domain_mask=.FALSE.
792    IF (.NOT. xios_getvar(prefix//"domain_mask_type", params%domain_mask_type) )   params%domain_mask_type="cross"
793    IF (.NOT. xios_getvar(prefix//"scalar_mask", params%scalar_mask) )   params%scalar_mask=.FALSE.
794    IF (.NOT. xios_getvar(prefix//"scalar_mask_type", params%scalar_mask_type) )   params%scalar_mask_type="none"
795    IF (.NOT. xios_getvar(prefix//"axis", params%axis) )                 params%axis="pressure"
796    IF (.NOT. xios_getvar(prefix//"axis_mask", params%axis_mask) )       params%axis_mask=.FALSE.
797    IF (.NOT. xios_getvar(prefix//"ni", params%ni) )                     params%ni=36
798    IF (.NOT. xios_getvar(prefix//"nj", params%nj) )                     params%nj=18
799    IF (.NOT. xios_getvar(prefix//"nlev", params%nlev) )                 params%nlev=10
800    IF (.NOT. xios_getvar(prefix//"init_field2D", params%init_field2D) ) params%init_field2D="academic"
801    IF (.NOT. xios_getvar(prefix//"pressure_factor", params%pressure_factor) ) params%pressure_factor=0.
802    IF (.NOT. xios_getvar(prefix//"mask3d", params%mask3d) ) params%mask3d=.FALSE.
803    IF (.NOT. xios_getvar(prefix//"field_sub_freq", params%field_sub_freq) ) params%field_sub_freq=1
804    IF (.NOT. xios_getvar(prefix//"field_sub_offset", params%field_sub_offset) ) params%field_sub_offset=0
805
806    IF (.NOT. xios_getvar(prefix//"domain_proc_n", params%domain_proc_n)) params%domain_proc_n=0
807    IF (.NOT. xios_getvar(prefix//"axis_proc_n", params%axis_proc_n)) params%axis_proc_n=0
808    IF (.NOT. xios_getvar(prefix//"ensemble_proc_n", params%ensemble_proc_n)) params%ensemble_proc_n=1
809    IF ((.NOT. xios_getvar(prefix//"domain_proc_frac", params%domain_proc_frac)) .AND.  &
810       (.NOT. xios_getvar("prefix//axis_proc_frac", params%axis_proc_frac))) THEN
811       params%domain_proc_frac=1.0
812       params%axis_proc_frac=0.0
813    ELSE IF (.NOT. xios_getvar(prefix//"domain_proc_frac", params%domain_proc_frac)) THEN
814      params%domain_proc_frac=0.0
815    ELSE IF (.NOT. xios_getvar(prefix//"axis_proc_frac", params%axis_proc_frac)) THEN
816      params%axis_proc_frac=0.0
817    ENDIF
818   
819  END SUBROUTINE init_model_params
820
821 
822  SUBROUTINE init_domain(domain_id, comm, params, return_ni, return_nj,             &
823                         return_lon, return_lat, return_mask, return_index,         &
824                         return_X_lon,return_X_lat, return_X_mask, return_X_index,  &
825                         return_Y_lon,return_Y_lat, return_Y_mask, return_Y_index)
826  IMPLICIT NONE
827    CHARACTER(LEN=*),INTENT(IN) :: domain_id
828    TYPE(tmodel_params),INTENT(IN) :: params
829    INTEGER,INTENT(IN) :: comm
830    INTEGER,INTENT(OUT) :: return_ni
831    INTEGER,INTENT(OUT) :: return_nj
832    DOUBLE PRECISION, POINTER :: return_lon(:)
833    DOUBLE PRECISION, POINTER :: return_lat(:)
834    LOGICAL, POINTER :: return_mask(:)
835    INTEGER, POINTER :: return_index(:)
836    DOUBLE PRECISION, POINTER :: return_X_lon(:)
837    DOUBLE PRECISION, POINTER :: return_X_lat(:)
838    LOGICAL, POINTER :: return_X_mask(:)
839    INTEGER, POINTER :: return_X_index(:)
840    DOUBLE PRECISION, POINTER :: return_Y_lon(:)
841    DOUBLE PRECISION, POINTER :: return_Y_lat(:)
842    LOGICAL, POINTER :: return_Y_mask(:)
843    INTEGER, POINTER :: return_Y_index(:)
844       
845    IF (params%domain=="lmdz") THEN
846      CALL init_domain_lmdz(domain_id,comm,params, return_ni, return_nj,               &
847                            return_lon, return_lat, return_mask, return_index,         &
848                            return_X_lon,return_X_lat, return_X_mask, return_X_index,  &
849                            return_Y_lon,return_Y_lat, return_Y_mask, return_Y_index)
850
851    ELSE IF (params%domain=="orchidee") THEN
852     CALL init_domain_orchidee(domain_id,comm,params, return_ni, return_nj,               &
853                               return_lon, return_lat, return_mask, return_index,         &
854                               return_X_lon,return_X_lat, return_X_mask, return_X_index,  &
855                               return_Y_lon,return_Y_lat, return_Y_mask, return_Y_index)
856
857    ELSE IF (params%domain=="nemo") THEN
858       CALL init_domain_nemo(domain_id,comm,params, return_ni, return_nj,               &
859                             return_lon, return_lat, return_mask, return_index,         &
860                             return_X_lon,return_X_lat, return_X_mask, return_X_index,  &
861                             return_Y_lon,return_Y_lat, return_Y_mask, return_Y_index)
862
863     ELSE IF (params%domain=="dynamico") THEN
864       CALL init_domain_dynamico(domain_id,comm,params, return_ni, return_nj,               &
865                                 return_lon, return_lat, return_mask, return_index,         &
866                                 return_X_lon,return_X_lat, return_X_mask, return_X_index,  &
867                                 return_Y_lon,return_Y_lat, return_Y_mask, return_Y_index)
868
869     ELSE IF (params%domain=="gaussian") THEN
870       CALL init_domain_gaussian(domain_id,comm,params, return_ni, return_nj,               &
871                                 return_lon, return_lat, return_mask, return_index,         &
872                                 return_X_lon,return_X_lat, return_X_mask, return_X_index,  &
873                                 return_Y_lon,return_Y_lat, return_Y_mask, return_Y_index)
874    ENDIF
875   
876  END SUBROUTINE init_domain
877
878  SUBROUTINE init_field2D(comm,params, lon, lat, mask, return_field,        &
879                                       x_lon, x_lat, x_mask, return_fieldx, &
880                                       y_lon, y_lat, y_mask, return_fieldy, return_fieldXY)
881  IMPLICIT NONE
882    TYPE(tmodel_params) :: params
883    TYPE(xios_context) :: ctx_hdl
884    INTEGER :: comm
885    DOUBLE PRECISION, POINTER :: lon(:)
886    DOUBLE PRECISION, POINTER :: lat(:)
887    LOGICAL, POINTER :: mask(:)
888    DOUBLE PRECISION, POINTER :: return_field(:)
889
890    DOUBLE PRECISION, POINTER :: X_lon(:)
891    DOUBLE PRECISION, POINTER :: X_lat(:)
892    LOGICAL, POINTER :: X_mask(:)
893    DOUBLE PRECISION, POINTER :: return_fieldX(:)
894
895    DOUBLE PRECISION, POINTER :: Y_lon(:)
896    DOUBLE PRECISION, POINTER :: Y_lat(:)
897    LOGICAL, POINTER :: Y_mask(:)
898    DOUBLE PRECISION, POINTER :: return_fieldY(:)
899    DOUBLE PRECISION, POINTER :: return_fieldXY(:,:)
900       
901    IF (params%init_field2D=="academic") THEN
902      CALL init_field2D_academic(comm,params, lon, lat, mask, return_field, X_lon, X_lat, X_mask, return_fieldX, &
903                                Y_lon, Y_lat, Y_mask, return_fieldY, return_fieldXY)
904    ELSE IF (params%init_field2D=="constant") THEN
905      CALL init_field2D_constant(comm,params, lon, lat, mask, return_field, X_lon, X_lat, X_mask, return_fieldX, &
906                                 Y_lon, Y_lat, Y_mask, return_fieldY, return_fieldXY)
907    ELSE IF (params%init_field2D=="rank") THEN
908      CALL init_field2D_rank(comm,params, lon, lat, mask, return_field, X_lon, X_lat, X_mask, &
909                             return_fieldX, Y_lon, Y_lat, Y_mask, return_fieldY, return_fieldXY)
910    ENDIF
911   
912  END SUBROUTINE init_field2D
913
914 
915  SUBROUTINE init_axis(axis_id,comm,params, return_value, return_mask, return_index)
916  IMPLICIT NONE
917    CHARACTER(LEN=*) :: axis_id
918    TYPE(tmodel_params) :: params
919    TYPE(xios_context) :: ctx_hdl
920    INTEGER :: comm
921    DOUBLE PRECISION, POINTER :: return_value(:)
922    LOGICAL, POINTER          :: return_mask(:)
923    INTEGER, POINTER          :: return_index(:)
924   
925    IF (params%axis=="pressure") THEN
926      CALL init_axis_pressure(axis_id,comm,params, return_value, return_mask, return_index)
927    ENDIF
928   
929   END SUBROUTINE init_axis
930   
931  SUBROUTINE init_ensemble(ensemble_id,comm,params, return_value)
932  IMPLICIT NONE
933    CHARACTER(LEN=*) :: ensemble_id
934    TYPE(tmodel_params) :: params
935    TYPE(xios_context) :: ctx_hdl
936    INTEGER :: comm
937    DOUBLE PRECISION, POINTER :: return_value(:)
938    INTEGER :: domain_proc_rank, domain_proc_size   
939    INTEGER :: axis_proc_rank, axis_proc_size
940    INTEGER :: ensemble_proc_size, ensemble_proc_rank
941
942    CALL get_decomposition(comm, params, domain_proc_size, domain_proc_rank, axis_proc_size, axis_proc_rank, ensemble_proc_size, ensemble_proc_rank)
943    ALLOCATE(return_value(0:0))
944    return_value(0)=ensemble_proc_rank
945
946    IF (xios_is_valid_axis(TRIM(ensemble_id))) THEN
947      CALL xios_set_axis_attr(ensemble_id, n_glo=ensemble_proc_size, begin=ensemble_proc_rank, n=1, value=return_value(:), unit='-')   
948    ENDIF
949   
950   END SUBROUTINE init_ensemble
951   
952
953  SUBROUTINE init_domain_gaussian(domain_id, comm, params, return_ni, return_nj,               &
954                                  return_lon,return_lat,return_mask, return_index,             &
955                                  return_X_lon,return_X_lat, return_X_mask, return_X_index,    &
956                                  return_Y_lon,return_Y_lat, return_Y_mask, return_Y_index)
957  IMPLICIT NONE
958    CHARACTER(LEN=*) :: domain_id
959    TYPE(tmodel_params) :: params
960    TYPE(xios_context) :: ctx_hdl
961    INTEGER :: comm
962    INTEGER :: return_ni
963    INTEGER :: return_nj
964    DOUBLE PRECISION, POINTER :: return_lon(:)
965    DOUBLE PRECISION, POINTER :: return_lat(:)
966    LOGICAL, POINTER :: return_mask(:)
967    INTEGER, POINTER :: return_index(:)
968    DOUBLE PRECISION, POINTER :: return_X_lon(:)
969    DOUBLE PRECISION, POINTER :: return_X_lat(:)
970    LOGICAL, POINTER :: return_X_mask(:)
971    INTEGER, POINTER :: return_X_index(:)
972    DOUBLE PRECISION, POINTER :: return_Y_lon(:)
973    DOUBLE PRECISION, POINTER :: return_Y_lat(:)
974    LOGICAL, POINTER :: return_Y_mask(:)
975    INTEGER, POINTER :: return_Y_index(:)
976    INTEGER :: domain_proc_rank, domain_proc_size   
977    INTEGER :: axis_proc_rank, axis_proc_size
978    INTEGER :: ensemble_proc_size, ensemble_proc_rank
979
980    INTEGER :: mpi_rank, mpi_size
981    INTEGER ::  ierr
982   
983!  INTEGER, PARAMETER :: nlon=60
984!  INTEGER, PARAMETER :: nlat=30
985!  INTEGER,PARAMETER :: ni_glo=100
986!  INTEGER,PARAMETER :: nj_glo=100
987!  INTEGER,PARAMETER :: llm=5
988!  DOUBLE PRECISION  :: lval(llm)=1
989!  TYPE(xios_field) :: field_hdl
990!  TYPE(xios_fieldgroup) :: fieldgroup_hdl
991!  TYPE(xios_file) :: file_hdl
992!  LOGICAL :: ok
993
994!  DOUBLE PRECISION,ALLOCATABLE :: lon_glo(:),lat_glo(:)
995!  DOUBLE PRECISION,ALLOCATABLE :: bounds_lon_glo(:,:),bounds_lat_glo(:,:)
996!  DOUBLE PRECISION,ALLOCATABLE :: field_A_glo(:,:)
997!  INTEGER,ALLOCATABLE :: i_index_glo(:)
998!  INTEGER,ALLOCATABLE :: i_index(:)
999!  LOGICAL,ALLOCATABLE :: mask_glo(:),mask(:)
1000!  INTEGER,ALLOCATABLE :: n_local(:),local_neighbor(:,:)
1001!  DOUBLE PRECISION,ALLOCATABLE :: lon(:),lat(:),field_A_srf(:,:), lonvalue(:) ;
1002!  DOUBLE PRECISION,ALLOCATABLE :: bounds_lon(:,:),bounds_lat(:,:) ;
1003!  INTEGER :: ni,ibegin,iend,nj,jbegin,jend
1004!  INTEGER :: i,j,l,ts,n, nbMax
1005!  INTEGER :: ncell_glo,ncell,ind
1006!  REAL :: ilon,ilat
1007!  DOUBLE PRECISION, PARAMETER :: Pi=3.14159265359
1008!  INTEGER :: list_ind(nlon,nlat)
1009!  INTEGER :: rank,j1,j2,np,ncell_x
1010!  INTEGER :: data_n_index
1011!  INTEGER,ALLOCATABLE :: data_i_index(:)
1012!  DOUBLE PRECISION,ALLOCATABLE :: field_A_compressed(:,:)
1013
1014
1015    INTEGER :: nlon, nlat
1016    DOUBLE PRECISION,ALLOCATABLE :: lon_glo(:),lat_glo(:)
1017    DOUBLE PRECISION,ALLOCATABLE :: bounds_lon_glo(:,:),bounds_lat_glo(:,:)
1018    INTEGER,ALLOCATABLE :: i_index_glo(:)
1019    INTEGER,ALLOCATABLE :: i_index(:)
1020    LOGICAL,ALLOCATABLE :: mask_glo(:),mask(:)
1021    DOUBLE PRECISION,ALLOCATABLE :: lon(:),lat(:),field_A_srf(:,:), lonvalue(:) ;
1022    DOUBLE PRECISION,ALLOCATABLE :: bounds_lon(:,:),bounds_lat(:,:) ;
1023
1024    INTEGER :: ni,ibegin,iend,nj,jbegin,jend
1025    INTEGER :: i,j,l,ts,n, nbMax
1026    INTEGER :: ncell_glo,ncell,ind
1027    REAL :: ilon,ilat
1028    DOUBLE PRECISION, PARAMETER :: Pi=3.14159265359
1029    INTEGER,ALLOCATABLE :: list_ind(:,:)
1030    INTEGER :: rank,j1,j2,np,ncell_x
1031    INTEGER :: data_n_index
1032    INTEGER,ALLOCATABLE :: data_i_index(:)
1033
1034
1035
1036    CALL MPI_COMM_RANK(comm,mpi_rank,ierr)
1037    CALL MPI_COMM_SIZE(comm,mpi_size,ierr)
1038
1039    CALL get_decomposition(comm, params, domain_proc_size, domain_proc_rank, axis_proc_size, axis_proc_rank, ensemble_proc_size, ensemble_proc_rank)
1040
1041    mpi_rank=domain_proc_rank
1042    mpi_size=domain_proc_size
1043    nlon=params%ni
1044    nlat=params%nj
1045
1046   
1047    ncell_glo=0
1048    DO j=1,nlat
1049      n =  NINT(COS(Pi/2-(j-0.5)*PI/nlat)*nlon)
1050      IF (n<8) n=8
1051      ncell_glo=ncell_glo+n
1052    ENDDO
1053
1054    ALLOCATE(lon_glo(ncell_glo))
1055    ALLOCATE(lat_glo(ncell_glo))
1056    ALLOCATE(bounds_lon_glo(4,ncell_glo))
1057    ALLOCATE(bounds_lat_glo(4,ncell_glo))
1058    ALLOCATE(i_index_glo(ncell_glo))
1059    ALLOCATE(mask_glo(ncell_glo))
1060    ALLOCATE(list_ind(nlon,nlat))
1061   
1062    ind=0
1063    DO j=1,nlat
1064      n = NINT(COS(Pi/2-(j-0.5)*PI/nlat)*nlon)
1065      if (j==1) PRINT*,"--- ",n
1066      if (j==nlat) PRINT*,"--- ",n
1067      IF (n<8) n=8
1068
1069      DO i=1,n
1070        ind=ind+1
1071        list_ind(i,j)=ind
1072        ilon=i-0.5
1073        ilat=j-0.5
1074
1075        lat_glo(ind)= 90-(ilat*180./nlat)
1076        lon_glo(ind)= (ilon*360./n)-180.
1077
1078
1079        bounds_lat_glo(1,ind)= 90-((ilat-0.5)*180./nlat)
1080        bounds_lon_glo(1,ind)=((ilon-0.5)*360./n)-180.
1081
1082        bounds_lat_glo(2,ind)= 90-((ilat-0.5)*180./nlat)
1083        bounds_lon_glo(2,ind)=((ilon+0.5)*360./n)-180.
1084
1085        bounds_lat_glo(3,ind)= 90-((ilat+0.5)*180./nlat)
1086        bounds_lon_glo(3,ind)=((ilon+0.5)*360./n)-180.
1087
1088        bounds_lat_glo(4,ind)= 90-((ilat+0.5)*180./nlat)
1089        bounds_lon_glo(4,ind)=((ilon-0.5)*360./n)-180.
1090
1091      ENDDO
1092    ENDDO
1093
1094!  mpi_size=32
1095    rank=(mpi_size-1)/2
1096    ncell_x=sqrt(ncell_glo*1./mpi_size)
1097
1098    j1=nlat/2
1099    DO WHILE(rank>=0)
1100      j2=MAX(j1-ncell_x+1,1)
1101      j=(j1+j2)/2
1102      n=NINT(COS(Pi/2-(j-0.5)*PI/nlat)*nlon)
1103      np = MIN(n/ncell_x,rank+1) ;
1104      if (j2==1) np=rank+1
1105 
1106      PRINT *,"domain ",j2,j1,rank,np ;
1107      DO j=j2,j1
1108        n=NINT(COS(Pi/2-(j-0.5)*PI/nlat)*nlon)
1109        IF (n<8) n=8
1110        DO i=1,n
1111          ind=list_ind(i,j)
1112          IF ( (i-1) < MOD(n,np)*(n/np+1)) THEN 
1113            i_index_glo(ind) = rank - (i-1)/(n/np+1)
1114          ELSE
1115            i_index_glo(ind) = rank-(MOD(n,np)+ (i-1-MOD(n,np)*(n/np+1))/(n/np))
1116          ENDIF
1117        ENDDO
1118      ENDDO
1119      rank=rank-np
1120      j1=j2-1
1121    ENDDO
1122
1123    rank=(mpi_size-1)/2+1
1124    ncell_x=sqrt(ncell_glo*1./mpi_size)
1125
1126    j1=nlat/2+1
1127    DO WHILE(rank<=mpi_size-1)
1128      j2=MIN(j1+ncell_x-1,nlat)
1129      j=(j1+j2)/2
1130      n=NINT(COS(Pi/2-(j-0.5)*PI/nlat)*nlon)
1131      np = MIN(n/ncell_x,mpi_size-rank) ;
1132      if (j2==nlat) np=mpi_size-rank
1133
1134      PRINT *,"domain ",j2,j1,rank,np ;
1135      DO j=j1,j2
1136        n=NINT(COS(Pi/2-(j-0.5)*PI/nlat)*nlon)
1137        IF (n<8) n=8
1138        DO i=1,n
1139          ind=list_ind(i,j)
1140          IF ( (i-1) < MOD(n,np)*(n/np+1)) THEN
1141            i_index_glo(ind) = rank + (i-1)/(n/np+1)
1142          ELSE
1143            i_index_glo(ind) = rank+(MOD(n,np)+ (i-1-MOD(n,np)*(n/np+1))/(n/np))
1144          ENDIF
1145        ENDDO
1146      ENDDO
1147      rank=rank+np
1148      j1=j2+1
1149   ENDDO
1150
1151    ncell=0
1152    DO ind=1,ncell_glo
1153      IF (i_index_glo(ind)==mpi_rank) ncell=ncell+1
1154    ENDDO
1155    ALLOCATE(i_index(ncell))
1156    ALLOCATE(lon(ncell))
1157    ALLOCATE(lat(ncell))
1158    ALLOCATE(bounds_lon(4,ncell))
1159    ALLOCATE(bounds_lat(4,ncell))
1160!    ALLOCATE(field_A_srf(ncell,llm))
1161    ALLOCATE(mask(ncell))
1162    ncell=0
1163    data_n_index=0
1164    DO ind=1,ncell_glo
1165      IF (i_index_glo(ind)==mpi_rank) THEN
1166        ncell=ncell+1
1167        i_index(ncell)=ind-1
1168        lon(ncell)=lon_glo(ind)
1169        lat(ncell)=lat_glo(ind)
1170        bounds_lon(:,ncell)=bounds_lon_glo(:,ind)
1171        bounds_lat(:,ncell)=bounds_lat_glo(:,ind)
1172        field_A_srf(ncell,:)=i_index_glo(ind)
1173        IF (MOD(ind,8)>=0 .AND. MOD(ind,8)<2) THEN
1174!        mask(ncell)=.FALSE.
1175
1176         mask(ncell)=.TRUE.
1177          data_n_index=data_n_index+1
1178        ELSE
1179          mask(ncell)=.TRUE.
1180          data_n_index=data_n_index+1
1181        ENDIF
1182      ENDIF
1183    ENDDO
1184
1185!    ALLOCATE(field_A_compressed(data_n_index,llm))
1186    ALLOCATE(data_i_index(data_n_index))
1187    data_n_index=0
1188    DO ind=1,ncell
1189      IF (mask(ind)) THEN
1190        data_n_index=data_n_index+1
1191        data_i_index(data_n_index)=ind-1
1192!        field_A_compressed(data_n_index,:)=field_A_srf(ind,:)
1193      ENDIF
1194    ENDDO
1195
1196    ALLOCATE(return_lon(0:ncell-1))
1197    ALLOCATE(return_lat(0:ncell-1))
1198    ALLOCATE(return_mask(0:ncell-1))
1199    ALLOCATE(return_index(0:data_n_index-1))
1200    return_lon(0:ncell-1)=lon(1:ncell)   
1201    return_lat(0:ncell-1)=lat(1:ncell)
1202    return_index(0:data_n_index-1)= data_i_index(1:data_n_index) 
1203    CALL set_domain_mask(params, lon, lat, return_mask)
1204
1205
1206    ALLOCATE(return_X_lon(0:ncell-1))
1207    ALLOCATE(return_X_lat(0:ncell-1))
1208    ALLOCATE(return_X_mask(0:ncell-1))
1209    ALLOCATE(return_X_index(0:ncell-1))
1210   
1211    return_X_lon=return_lon
1212    return_X_lat=return_lat
1213    return_X_index=return_index
1214    CALL set_domain_mask(params, return_X_lon,return_X_lat, return_X_mask)
1215
1216
1217    ALLOCATE(return_Y_lon(0:0))
1218    ALLOCATE(return_Y_lat(0:0))
1219    ALLOCATE(return_Y_mask(0:0))
1220    ALLOCATE(return_Y_index(0:0))
1221
1222    return_Y_lon(0)=lon_glo(ncell_glo/2)
1223    return_Y_lat(0)=lat_glo(ncell_glo/2)
1224    CALL set_domain_mask(params, return_Y_lon,return_Y_lat, return_Y_mask)
1225    return_Y_index(0)=0
1226
1227    return_ni=ncell
1228    return_nj=1
1229
1230
1231    IF (xios_is_valid_domain(TRIM(domain_id))) THEN
1232      CALL xios_set_domain_attr(TRIM(domain_id), type="unstructured", ni_glo=ncell_glo, ni=ncell, ibegin=0, i_index=i_index)
1233      CALL xios_set_domain_attr(TRIM(domain_id), data_dim=1, data_ni=data_n_index, data_i_index=data_i_index, mask_1d=return_mask)
1234      CALL xios_set_domain_attr(TRIM(domain_id), lonvalue_1D=lon, latvalue_1D=lat, nvertex=4, bounds_lon_1D=bounds_lon, bounds_lat_1D=bounds_lat)
1235    ENDIF
1236
1237    IF (xios_is_valid_axis(TRIM(domain_id)//"_X")) THEN
1238      CALL xios_set_axis_attr(TRIM(domain_id)//"_X", n_glo=ncell_glo, n=ncell, begin=0, index=i_index,value=return_X_lon)
1239      CALL xios_set_axis_attr(TRIM(domain_id)//"_X", data_n=data_n_index, data_index=data_i_index, mask=return_X_mask)
1240    ENDIF
1241   
1242    IF (xios_is_valid_axis(TRIM(domain_id)//"_Y")) THEN   
1243      CALL xios_set_axis_attr(TRIM(domain_id)//"_Y", n_glo=1, begin=0, n=1, value=return_Y_lat, mask=return_Y_mask)
1244    ENDIF
1245
1246
1247
1248  END SUBROUTINE init_domain_gaussian
1249
1250
1251
1252  SUBROUTINE init_domain_dynamico(domain_id, comm, params, return_ni, return_nj,               &
1253                              return_lon,return_lat,return_mask, return_index,             &
1254                              return_X_lon,return_X_lat, return_X_mask, return_X_index,    &
1255                              return_Y_lon,return_Y_lat, return_Y_mask, return_Y_index)
1256  IMPLICIT NONE
1257    CHARACTER(LEN=*) :: domain_id
1258    TYPE(tmodel_params) :: params
1259    TYPE(xios_context) :: ctx_hdl
1260    INTEGER :: comm
1261    INTEGER :: return_ni
1262    INTEGER :: return_nj
1263    DOUBLE PRECISION, POINTER :: return_lon(:)
1264    DOUBLE PRECISION, POINTER :: return_lat(:)
1265    LOGICAL, POINTER :: return_mask(:)
1266    INTEGER, POINTER :: return_index(:)
1267    DOUBLE PRECISION, POINTER :: return_X_lon(:)
1268    DOUBLE PRECISION, POINTER :: return_X_lat(:)
1269    LOGICAL, POINTER :: return_X_mask(:)
1270    INTEGER, POINTER :: return_X_index(:)
1271    DOUBLE PRECISION, POINTER :: return_Y_lon(:)
1272    DOUBLE PRECISION, POINTER :: return_Y_lat(:)
1273    LOGICAL, POINTER :: return_Y_mask(:)
1274    INTEGER, POINTER :: return_Y_index(:)
1275    INTEGER :: domain_proc_rank, domain_proc_size   
1276    INTEGER :: axis_proc_rank, axis_proc_size
1277    INTEGER :: ensemble_proc_size, ensemble_proc_rank
1278
1279    INTEGER :: mpi_rank, mpi_size
1280    INTEGER ::  ierr
1281    INTEGER :: ni,nj,ni_glo,nj_glo,nvertex
1282    INTEGER :: ibegin,jbegin
1283    INTEGER :: nbp,nbp_glo, offset
1284    DOUBLE PRECISION, ALLOCATABLE :: lon_glo(:), lat_glo(:), lon(:), lat(:)
1285    DOUBLE PRECISION, ALLOCATABLE :: bounds_lon_glo(:,:), bounds_lat_glo(:,:), bounds_lon(:,:), bounds_lat(:,:)
1286    LOGICAL,ALLOCATABLE :: mask(:)
1287    LOGICAL,ALLOCATABLE :: dom_mask(:)
1288    INTEGER :: i,j
1289    INTEGER,ALLOCATABLE :: ibegin_all(:), ni_all(:)
1290       
1291    CALL MPI_COMM_RANK(comm,mpi_rank,ierr)
1292    CALL MPI_COMM_SIZE(comm,mpi_size,ierr)
1293
1294    CALL get_decomposition(comm, params, domain_proc_size, domain_proc_rank, axis_proc_size, axis_proc_rank, ensemble_proc_size, ensemble_proc_rank)
1295
1296    CALL xios_get_current_context(ctx_hdl)
1297   
1298    CALL xios_context_initialize("grid_dynamico",comm)
1299    CALL xios_close_context_definition()
1300    CALL xios_get_domain_attr("Ai::",ni_glo=ni_glo,nj_glo=nj_glo,ni=ni,nj=nj, ibegin=ibegin, jbegin=jbegin ,nvertex=nvertex)
1301    ALLOCATE(lon(ni),lat(ni),bounds_lon(nvertex,ni),bounds_lat(nvertex,ni))
1302    CALL xios_get_domain_attr("Ai::", lonvalue_1d=lon, latvalue_1d=lat, bounds_lon_1d=bounds_lon, bounds_lat_1d=bounds_lat)
1303    CALL xios_context_finalize
1304
1305    CALL xios_set_current_context(ctx_hdl)
1306
1307   
1308    ALLOCATE(ibegin_all(mpi_size))
1309    ALLOCATE(ni_all(mpi_size))
1310    ALLOCATE(lon_glo(0:ni_glo-1))
1311    ALLOCATE(lat_glo(0:ni_glo-1))
1312    ALLOCATE(bounds_lon_glo(nvertex,0:ni_glo-1))
1313    ALLOCATE(bounds_lat_glo(nvertex,0:ni_glo-1))
1314   
1315    CALL MPI_Allgather(ibegin, 1, MPI_INTEGER, ibegin_all, 1, MPI_INTEGER, comm, ierr)
1316    CALL MPI_Allgather(ni, 1, MPI_INTEGER, ni_all, 1, MPI_INTEGER, comm, ierr)
1317
1318    CALL MPI_AllgatherV(lon, ni, MPI_REAL8, lon_glo, ni_all, ibegin_all, MPI_REAL8, comm, ierr) 
1319    CALL MPI_AllgatherV(lat, ni, MPI_REAL8, lat_glo, ni_all, ibegin_all, MPI_REAL8, comm, ierr) 
1320    CALL MPI_AllgatherV(bounds_lon, ni*nvertex, MPI_REAL8, bounds_lon_glo, ni_all*nvertex, ibegin_all*nvertex, MPI_REAL8, comm, ierr) 
1321    CALL MPI_AllgatherV(bounds_lat, ni*nvertex, MPI_REAL8, bounds_lat_glo, ni_all*nvertex, ibegin_all*nvertex, MPI_REAL8, comm, ierr) 
1322
1323   
1324    nbp_glo=ni_glo
1325    nbp=nbp_glo/domain_proc_size
1326    IF (domain_proc_rank < MOD(nbp_glo,domain_proc_size) ) THEN
1327      nbp=nbp+1
1328      offset=nbp*domain_proc_rank
1329    ELSE
1330      offset=nbp*domain_proc_rank + MOD(nbp_glo, domain_proc_size)
1331    ENDIF
1332
1333    DEALLOCATE(lat,lon,bounds_lon,bounds_lat)
1334    ALLOCATE(lat(0:nbp-1))
1335    ALLOCATE(lon(0:nbp-1))
1336    ALLOCATE(bounds_lon(nvertex,0:nbp-1))
1337    ALLOCATE(bounds_lat(nvertex,0:nbp-1))
1338    ALLOCATE(return_lon(0:nbp-1))
1339    ALLOCATE(return_lat(0:nbp-1))
1340    ALLOCATE(return_index(0:nbp-1))
1341
1342    DO i=0,nbp-1
1343      lat(i)=lat_glo(i+offset)
1344      lon(i)=lon_glo(i+offset)
1345      bounds_lon(:,i) = bounds_lon_glo(:,i+offset) 
1346      bounds_lat(:,i) = bounds_lat_glo(:,i+offset)
1347      return_index(i)=i 
1348    ENDDO
1349    return_lon=lon
1350    return_lat=lat
1351   
1352    ALLOCATE(return_mask(0:nbp-1))
1353    CALL set_domain_mask(params, lon, lat, return_mask)
1354
1355
1356    ALLOCATE(return_X_lon(0:nbp-1))
1357    ALLOCATE(return_X_lat(0:nbp-1))
1358    ALLOCATE(return_X_mask(0:nbp-1))
1359    ALLOCATE(return_X_index(0:nbp-1))
1360   
1361    DO i=0,nbp-1
1362      return_X_lon(i)=lon_glo(i+offset)
1363      return_X_lat(i)=lat_glo(i+offset)
1364    ENDDO
1365    CALL set_domain_mask(params, return_X_lon,return_X_lat, return_X_mask)
1366    DO i=0,nbp-1
1367      return_X_index(i)=i
1368    ENDDO
1369
1370    ALLOCATE(return_Y_lon(0:0))
1371    ALLOCATE(return_Y_lat(0:0))
1372    ALLOCATE(return_Y_mask(0:0))
1373    ALLOCATE(return_Y_index(0:0))
1374
1375    return_Y_lon(0)=lon_glo(nbp_glo/2)
1376    return_Y_lat(0)=lat_glo(nbp_glo/2)
1377    CALL set_domain_mask(params, return_Y_lon,return_Y_lat, return_Y_mask)
1378    return_Y_index(0)=0
1379 
1380
1381   
1382    return_ni=nbp
1383    return_nj=1
1384
1385    IF (xios_is_valid_domain(TRIM(domain_id))) THEN
1386      CALL xios_set_domain_attr(TRIM(domain_id), type="unstructured", ni_glo=ni_glo, ibegin=offset, ni=nbp, nvertex=nvertex)
1387      CALL xios_set_domain_attr(TRIM(domain_id), data_dim=1, lonvalue_1d=lon, latvalue_1d=lat, bounds_lon_1d=bounds_lon, bounds_lat_1d=bounds_lat, mask_1d=return_mask)
1388    ENDIF   
1389
1390    IF (xios_is_valid_axis(TRIM(domain_id)//"_X")) THEN
1391      CALL xios_set_axis_attr(TRIM(domain_id)//"_X", n_glo=ni_glo, begin=offset, n=nbp, value=return_X_lon)
1392    ENDIF
1393
1394    IF (xios_is_valid_axis(TRIM(domain_id)//"_Y")) THEN   
1395      CALL xios_set_axis_attr(TRIM(domain_id)//"_Y", n_glo=1, begin=0, n=1, value=return_Y_lat, mask=return_Y_mask)
1396    ENDIF
1397
1398   END SUBROUTINE init_domain_dynamico
1399   
1400
1401
1402
1403  SUBROUTINE init_domain_lmdz(domain_id, comm, params, return_ni, return_nj,               &
1404                              return_lon,return_lat,return_mask, return_index,             &
1405                              return_X_lon,return_X_lat, return_X_mask, return_X_index,    &
1406                              return_Y_lon,return_Y_lat, return_Y_mask, return_Y_index)
1407  IMPLICIT NONE
1408    CHARACTER(LEN=*) :: domain_id
1409    TYPE(tmodel_params) :: params
1410    TYPE(xios_context) :: ctx_hdl
1411    INTEGER :: comm
1412    INTEGER :: return_ni
1413    INTEGER :: return_nj
1414    DOUBLE PRECISION, POINTER :: return_lon(:)
1415    DOUBLE PRECISION, POINTER :: return_lat(:)
1416    LOGICAL, POINTER :: return_mask(:)
1417    INTEGER, POINTER :: return_index(:)
1418    DOUBLE PRECISION, POINTER :: return_X_lon(:)
1419    DOUBLE PRECISION, POINTER :: return_X_lat(:)
1420    LOGICAL, POINTER :: return_X_mask(:)
1421    INTEGER, POINTER :: return_X_index(:)
1422    DOUBLE PRECISION, POINTER :: return_Y_lon(:)
1423    DOUBLE PRECISION, POINTER :: return_Y_lat(:)
1424    LOGICAL, POINTER :: return_Y_mask(:)
1425    INTEGER, POINTER :: return_Y_index(:)
1426    INTEGER :: domain_proc_rank, domain_proc_size   
1427    INTEGER :: axis_proc_rank, axis_proc_size
1428    INTEGER :: ensemble_proc_size, ensemble_proc_rank
1429
1430    INTEGER :: mpi_rank, mpi_size
1431    INTEGER ::  ierr
1432    INTEGER :: ni,nj,ni_glo,nj_glo
1433    INTEGER :: ibegin,jbegin
1434    INTEGER :: nbp,nbp_glo, offset
1435    DOUBLE PRECISION, ALLOCATABLE :: lon_glo(:), lat_glo(:), lon(:), lat(:)
1436    LOGICAL,ALLOCATABLE :: mask(:)
1437    LOGICAL,ALLOCATABLE :: dom_mask(:)
1438    INTEGER :: i,j
1439       
1440    CALL MPI_COMM_RANK(comm,mpi_rank,ierr)
1441    CALL MPI_COMM_SIZE(comm,mpi_size,ierr)
1442
1443    CALL get_decomposition(comm, params, domain_proc_size, domain_proc_rank, axis_proc_size, axis_proc_rank, ensemble_proc_size, ensemble_proc_rank)
1444    ni_glo=params%ni
1445    nj_glo=params%nj
1446    nbp_glo=ni_glo*nj_glo
1447    nbp=nbp_glo/domain_proc_size
1448    IF (domain_proc_rank < MOD(nbp_glo,domain_proc_size) ) THEN
1449     nbp=nbp+1
1450     offset=nbp*domain_proc_rank
1451    ELSE
1452      offset=nbp*domain_proc_rank + MOD(nbp_glo, domain_proc_size)
1453    ENDIF
1454   
1455   
1456    ibegin=0 ; ni=ni_glo
1457    jbegin=offset / ni_glo
1458   
1459    nj = (offset+nbp-1) / ni_glo - jbegin + 1 
1460   
1461    offset=offset-jbegin*ni
1462   
1463    ALLOCATE(lon(0:ni-1), lat(0:nj-1), mask(0:ni*nj-1), dom_mask(0:ni*nj-1))
1464    mask(:)=.FALSE.
1465    mask(offset:offset+nbp-1)=.TRUE.
1466   
1467    ALLOCATE(lon_glo(0:ni_glo-1), lat_glo(0:nj_glo-1))
1468   
1469    DO i=0,ni_glo-1
1470      lon_glo(i)=-180+(i+0.5)*(360./ni_glo)
1471    ENDDO
1472
1473    DO j=0,nj_glo-1
1474      lat_glo(j)=-90+(j+0.5)*(180./nj_glo)
1475    ENDDO
1476     
1477    lon(:)=lon_glo(:)
1478    lat(:)=lat_glo(jbegin:jbegin+nj-1)
1479
1480    ALLOCATE(return_lon(0:ni*nj-1))
1481    ALLOCATE(return_lat(0:ni*nj-1))
1482    ALLOCATE(return_mask(0:ni*nj-1))
1483    ALLOCATE(return_index(0:ni*nj-1))
1484
1485    ALLOCATE(return_X_lon(0:ni-1))
1486    ALLOCATE(return_X_lat(0:ni-1))
1487    ALLOCATE(return_X_mask(0:ni-1))
1488    ALLOCATE(return_X_index(0:ni-1))
1489
1490    ALLOCATE(return_Y_lon(0:nj-1))
1491    ALLOCATE(return_Y_lat(0:nj-1))
1492    ALLOCATE(return_Y_mask(0:nj-1))
1493    ALLOCATE(return_Y_index(0:nj-1))
1494
1495    DO j=0,nj-1
1496      DO i=0,ni-1
1497        return_lon(i+j*ni)=lon(i)
1498        return_lat(i+j*ni)=lat(j)
1499        return_index(i+j*ni)=i+j*ni
1500      ENDDO
1501    ENDDO
1502
1503    CALL set_domain_mask(params, return_lon,return_lat, dom_mask)
1504   
1505    return_mask = mask .AND. dom_mask
1506
1507    return_X_lat(:)=lat_glo(nj_glo/2)
1508    return_X_lon(:)=lon_glo(:)
1509    CALL set_domain_mask(params, return_X_lon,return_X_lat, return_X_mask)
1510    DO i=0,ni-1
1511      return_X_index(i)=i
1512    ENDDO
1513
1514    return_Y_lon(:)=lon_glo(ni_glo/2)
1515    return_Y_lat(:)=lat_glo(jbegin:jbegin+nj-1)
1516    CALL set_domain_mask(params, return_Y_lon,return_Y_lat, return_Y_mask)
1517
1518    DO j=0,nj-1
1519      return_Y_index(j)=j
1520    ENDDO
1521
1522    return_ni=ni
1523    return_nj=nj
1524
1525    IF (xios_is_valid_domain(TRIM(domain_id))) THEN
1526      CALL xios_set_domain_attr(TRIM(domain_id), type="rectilinear", ni_glo=ni_glo, ibegin=ibegin, ni=ni, nj_glo=nj_glo, jbegin=jbegin, nj=nj)
1527      CALL xios_set_domain_attr(TRIM(domain_id), data_dim=2, lonvalue_1d=lon, latvalue_1d=lat, mask_1d=return_mask)
1528    ENDIF
1529   
1530    IF (xios_is_valid_axis(TRIM(domain_id)//"_X")) THEN
1531      CALL xios_set_axis_attr(TRIM(domain_id)//"_X", n_glo=ni_glo, begin=ibegin, n=ni, value=return_X_lon, mask=return_X_mask)
1532    ENDIF
1533
1534    IF (xios_is_valid_axis(TRIM(domain_id)//"_Y")) THEN   
1535      CALL xios_set_axis_attr(TRIM(domain_id)//"_Y", n_glo=nj_glo, begin=jbegin, n=nj, value=return_Y_lat, mask=return_Y_mask)
1536    ENDIF
1537
1538  END SUBROUTINE init_domain_lmdz
1539
1540
1541  SUBROUTINE init_domain_orchidee(domain_id, comm, params, return_ni, return_nj,               &
1542                                  return_lon,return_lat,return_mask, return_index,             &
1543                                  return_X_lon,return_X_lat, return_X_mask, return_X_index,    &
1544                                  return_Y_lon,return_Y_lat, return_Y_mask, return_Y_index)
1545  IMPLICIT NONE
1546    CHARACTER(LEN=*) :: domain_id
1547    TYPE(tmodel_params) :: params
1548    TYPE(xios_context) :: ctx_hdl
1549    INTEGER :: comm
1550    INTEGER :: return_ni
1551    INTEGER :: return_nj
1552    DOUBLE PRECISION, POINTER :: return_lon(:)
1553    DOUBLE PRECISION, POINTER :: return_lat(:)
1554    LOGICAL, POINTER :: return_mask(:)
1555    INTEGER, POINTER :: return_index(:)
1556    DOUBLE PRECISION, POINTER :: return_X_lon(:)
1557    DOUBLE PRECISION, POINTER :: return_X_lat(:)
1558    LOGICAL, POINTER :: return_X_mask(:)
1559    INTEGER, POINTER :: return_X_index(:)
1560    DOUBLE PRECISION, POINTER :: return_Y_lon(:)
1561    DOUBLE PRECISION, POINTER :: return_Y_lat(:)
1562    LOGICAL, POINTER :: return_Y_mask(:)
1563    INTEGER, POINTER :: return_Y_index(:)
1564    INTEGER :: domain_proc_rank, domain_proc_size   
1565    INTEGER :: axis_proc_rank, axis_proc_size
1566    INTEGER :: ensemble_proc_size, ensemble_proc_rank
1567
1568    INTEGER :: mpi_rank, mpi_size
1569    INTEGER ::  ierr
1570    INTEGER :: ni,nj,ni_glo,nj_glo
1571    INTEGER :: ibegin,jbegin
1572    INTEGER :: nbp,nbp_glo, offset
1573    DOUBLE PRECISION, ALLOCATABLE :: lon_glo(:), lat_glo(:), lon(:), lat(:)
1574    LOGICAL,ALLOCATABLE :: mask(:)
1575    INTEGER :: i,j, ij, pos, i_glo, j_glo, ij_begin, ij_end
1576       
1577    CALL MPI_COMM_RANK(comm,mpi_rank,ierr)
1578    CALL MPI_COMM_SIZE(comm,mpi_size,ierr)
1579
1580    CALL get_decomposition(comm, params, domain_proc_size, domain_proc_rank, axis_proc_size, axis_proc_rank, ensemble_proc_size, ensemble_proc_rank)
1581    ni_glo=params%ni
1582    nj_glo=params%nj
1583    nbp_glo=ni_glo*nj_glo
1584    nbp=nbp_glo/domain_proc_size
1585    IF (domain_proc_rank < MOD(nbp_glo,domain_proc_size) ) THEN
1586     nbp=nbp+1
1587     offset=nbp*domain_proc_rank
1588    ELSE
1589      offset=nbp*domain_proc_rank + MOD(nbp_glo, domain_proc_size)
1590    ENDIF
1591   
1592   
1593    ibegin=0 ; ni=ni_glo
1594    jbegin=offset / ni_glo
1595   
1596    nj = (offset+nbp-1) / ni_glo - jbegin + 1 
1597   
1598    offset=offset-jbegin*ni
1599
1600    ij_begin=offset+jbegin*ni
1601    ij_end=ij_begin+nbp-1
1602   
1603    ALLOCATE(lon(0:ni-1), lat(0:nj-1), mask(0:ni*nj-1))
1604    mask(:)=.FALSE.
1605    mask(offset:offset+nbp-1)=.TRUE.
1606   
1607    ALLOCATE(lon_glo(0:ni_glo-1), lat_glo(0:nj_glo-1))
1608   
1609    DO i=0,ni_glo-1
1610      lon_glo(i)=-180+(i+0.5)*(360./ni_glo)
1611    ENDDO
1612
1613    DO j=0,nj_glo-1
1614      lat_glo(j)=-90+(j+0.5)*(180./nj_glo)
1615    ENDDO
1616     
1617    lon(:)=lon_glo(:)
1618    lat(:)=lat_glo(jbegin:jbegin+nj-1)
1619
1620    ALLOCATE(return_lon(0:ni*nj-1))
1621    ALLOCATE(return_lat(0:ni*nj-1))
1622    ALLOCATE(return_mask(0:ni*nj-1))
1623
1624    ALLOCATE(return_X_lon(0:ni-1))
1625    ALLOCATE(return_X_lat(0:ni-1))
1626    ALLOCATE(return_X_mask(0:ni-1))
1627
1628    ALLOCATE(return_Y_lon(0:nj-1))
1629    ALLOCATE(return_Y_lat(0:nj-1))
1630    ALLOCATE(return_Y_mask(0:nj-1))
1631
1632
1633     DO j=0,nj-1
1634      DO i=0,ni-1
1635        return_lon(i+j*ni)=lon(i)
1636        return_lat(i+j*ni)=lat(j)
1637      ENDDO
1638    ENDDO
1639
1640    pos=0
1641    DO j=0,nj-1
1642      DO i=0,ni-1
1643        i_glo=i
1644        j_glo=jbegin+j
1645        ij=j_glo*ni_glo+i_glo
1646        IF ( ij>=ij_begin .AND. ij<=ij_end) THEN
1647           IF ((MOD(j_glo,3)==1 .OR. MOD(j_glo,3)==2) .AND. (MOD(i_glo,5)==3 .OR. MOD(i_glo,5)==4)) CYCLE
1648           pos=pos+1
1649        ENDIF
1650      ENDDO
1651   ENDDO
1652
1653   ALLOCATE(return_index(0:pos-1))
1654   return_index(:)=-1
1655   pos=0
1656    DO j=0,nj-1
1657      DO i=0,ni-1
1658        i_glo=i
1659        j_glo=jbegin+j
1660        ij=j_glo*ni_glo+i_glo
1661        IF ( ij>=ij_begin .AND. ij<=ij_end) THEN
1662           IF ((MOD(j_glo,3)==1 .OR. MOD(j_glo,3)==2) .AND. (MOD(i_glo,5)==3 .OR. MOD(i_glo,5)==4)) CYCLE
1663            return_index(pos)=i+j*ni
1664           pos=pos+1
1665        ENDIF
1666      ENDDO
1667   ENDDO     
1668
1669    CALL set_domain_mask(params, return_lon,return_lat, mask)
1670
1671    return_mask=mask
1672
1673    return_X_lat(:)=lat_glo(nj_glo/2)
1674    return_X_lon(:)=lon_glo(:)
1675    CALL set_domain_mask(params, return_X_lon,return_X_lat, return_X_mask)
1676
1677   pos=0
1678    DO i=0,ni-1
1679      i_glo=i
1680      j_glo=nj_glo/2
1681      IF ((MOD(j_glo,3)==1 .OR. MOD(j_glo,3)==2) .AND. (MOD(i_glo,5)==3 .OR. MOD(i_glo,5)==4)) CYCLE
1682      pos=pos+1
1683    ENDDO
1684
1685    ALLOCATE(return_X_index(0:pos-1))
1686    return_X_index(:)=-1
1687    pos=0
1688    DO i=0,ni-1
1689      i_glo=i
1690      j_glo=nj_glo/2
1691      IF ((MOD(j_glo,3)==1 .OR. MOD(j_glo,3)==2) .AND. (MOD(i_glo,5)==3 .OR. MOD(i_glo,5)==4)) CYCLE
1692      return_X_index(pos)=i
1693      pos=pos+1
1694    ENDDO 
1695
1696
1697   
1698    return_Y_lon(:)=lon_glo(ni_glo/2)
1699    return_Y_lat(:)=lat_glo(jbegin:jbegin+nj-1)
1700    CALL set_domain_mask(params, return_Y_lon,return_Y_lat, return_Y_mask)
1701
1702    pos=0
1703    DO j=0,nj-1
1704      i_glo=ni_glo/2
1705      j_glo=j+jbegin
1706      IF ((MOD(j_glo,3)==1 .OR. MOD(j_glo,3)==2) .AND. (MOD(i_glo,5)==3 .OR. MOD(i_glo,5)==4)) CYCLE
1707      pos=pos+1
1708    ENDDO
1709
1710    ALLOCATE(return_Y_index(0:pos-1))
1711    return_Y_index=-1
1712    pos=0
1713    DO j=0,nj-1
1714      i_glo=ni_glo/2
1715      j_glo=j+jbegin
1716      IF ((MOD(j_glo,3)==1 .OR. MOD(j_glo,3)==2) .AND. (MOD(i_glo,5)==3 .OR. MOD(i_glo,5)==4)) CYCLE
1717      return_Y_index(pos)=j
1718      pos=pos+1
1719    ENDDO
1720
1721    return_ni=ni
1722    return_nj=nj
1723
1724    IF (xios_is_valid_domain(TRIM(domain_id))) THEN
1725      CALL xios_set_domain_attr(TRIM(domain_id), type="rectilinear", ni_glo=ni_glo, ibegin=ibegin, ni=ni, nj_glo=nj_glo, jbegin=jbegin, nj=nj)
1726      CALL xios_set_domain_attr(TRIM(domain_id), data_dim=1, data_ni=size(return_index), data_i_index=return_index)
1727      CALL xios_set_domain_attr(TRIM(domain_id), lonvalue_1d=lon, latvalue_1d=lat, mask_1d=mask)
1728    ENDIF
1729   
1730    IF (xios_is_valid_axis(TRIM(domain_id)//"_X")) THEN
1731      CALL xios_set_axis_attr(TRIM(domain_id)//"_X", n_glo=ni_glo, begin=ibegin, n=ni, value=return_X_lon, data_n=size(return_X_index), data_index=return_X_index)
1732    ENDIF
1733
1734    IF (xios_is_valid_axis(TRIM(domain_id)//"_Y")) THEN   
1735      CALL xios_set_axis_attr(TRIM(domain_id)//"_Y", n_glo=nj_glo, begin=jbegin, n=nj, value=return_Y_lat, data_n=size(return_Y_index), data_index=return_Y_index)
1736    ENDIF
1737
1738  END SUBROUTINE init_domain_orchidee
1739
1740
1741
1742  SUBROUTINE init_domain_nemo(domain_id, comm, params, return_ni, return_nj,               &
1743                              return_lon,return_lat,return_mask, return_index,             &
1744                              return_X_lon,return_X_lat, return_X_mask, return_X_index,    &
1745                              return_Y_lon,return_Y_lat, return_Y_mask, return_Y_index)
1746  IMPLICIT NONE
1747    CHARACTER(LEN=*) :: domain_id
1748    TYPE(tmodel_params) :: params
1749    TYPE(xios_context) :: ctx_hdl
1750    INTEGER :: comm
1751    INTEGER :: return_ni
1752    INTEGER :: return_nj
1753    DOUBLE PRECISION, POINTER :: return_lon(:)
1754    DOUBLE PRECISION, POINTER :: return_lat(:)
1755    LOGICAL, POINTER :: return_mask(:)
1756    INTEGER, POINTER :: return_index(:)
1757    DOUBLE PRECISION, POINTER :: return_X_lon(:)
1758    DOUBLE PRECISION, POINTER :: return_X_lat(:)
1759    LOGICAL, POINTER :: return_X_mask(:)
1760    INTEGER, POINTER :: return_X_index(:)
1761    DOUBLE PRECISION, POINTER :: return_Y_lon(:)
1762    DOUBLE PRECISION, POINTER :: return_Y_lat(:)
1763    LOGICAL, POINTER :: return_Y_mask(:)
1764    INTEGER, POINTER :: return_Y_index(:)
1765    INTEGER :: domain_proc_rank, domain_proc_size   
1766    INTEGER :: axis_proc_rank, axis_proc_size
1767    INTEGER :: ensemble_proc_size, ensemble_proc_rank
1768
1769    INTEGER :: mpi_rank, mpi_size
1770    INTEGER ::  ierr
1771    INTEGER :: ni,nj,ni_glo,nj_glo
1772    INTEGER :: ibegin,jbegin
1773    INTEGER :: offset_i, offset_j
1774    DOUBLE PRECISION, ALLOCATABLE :: lon_glo(:), lat_glo(:), bounds_lon_glo(:,:), bounds_lat_glo(:,:)
1775    DOUBLE PRECISION, ALLOCATABLE :: lon(:,:), lat(:,:), bounds_lon(:,:,:), bounds_lat(:,:,:) 
1776    LOGICAL,ALLOCATABLE :: mask(:)
1777    INTEGER :: i,j, ij, n, rank
1778    INTEGER :: nproc_i, nproc_j, nholes, pos_hole
1779    INTEGER,ALLOCATABLE :: size_i(:), begin_i(:), size_j(:), begin_j(:)
1780       
1781    CALL MPI_COMM_RANK(comm,mpi_rank,ierr)
1782    CALL MPI_COMM_SIZE(comm,mpi_size,ierr)
1783
1784    CALL get_decomposition(comm, params, domain_proc_size, domain_proc_rank, axis_proc_size, axis_proc_rank, ensemble_proc_size, ensemble_proc_rank)
1785    ni_glo=params%ni
1786    nj_glo=params%nj
1787
1788    n=INT(SQRT(mpi_size*1.))
1789    nproc_i=n
1790    nproc_j=n
1791    IF ( n*n == mpi_size) THEN
1792    ! do nothing
1793    ELSE IF ( (n+1)*n < mpi_size) THEN
1794      nproc_i=nproc_i+1
1795      nproc_j=nproc_j+1
1796    ELSE
1797      nproc_i=nproc_i+1
1798    ENDIF
1799   
1800    nholes=nproc_i*nproc_j-mpi_size
1801
1802    ALLOCATE(size_i(0:nproc_i-1))
1803    ALLOCATE(begin_i(0:nproc_i-1))
1804    DO i=0,nproc_i-1
1805      size_i(i)=ni_glo/nproc_i
1806      IF (i<MOD(ni_glo,nproc_i)) size_i(i)=size_i(i)+1
1807      IF (i==0) THEN
1808        begin_i(i)=0
1809      ELSE
1810        begin_i(i)=begin_i(i-1)+size_i(i-1)
1811      ENDIF
1812    ENDDO
1813   
1814    ALLOCATE(size_j(0:nproc_j-1))
1815    ALLOCATE(begin_j(0:nproc_j-1))
1816    DO j=0,nproc_i-1
1817      size_j(j)=nj_glo/nproc_j
1818      IF (j<MOD(nj_glo,nproc_j)) size_j(j)=size_j(j)+1
1819      IF (j==0) THEN
1820        begin_j(j)=0
1821      ELSE
1822        begin_j(j)=begin_j(j-1)+size_j(j-1)
1823      ENDIF
1824    ENDDO
1825
1826
1827    pos_hole=0
1828    rank=0
1829    DO j=0,nproc_j-1
1830      DO i=0,nproc_i-1
1831
1832        ij = j*nproc_i + i
1833        IF (pos_hole<nholes) THEN
1834          IF ( MOD(ij,(nproc_i*nproc_j/nholes)) == 0) THEN
1835            pos_hole=pos_hole+1
1836            CYCLE
1837          ENDIF
1838        ENDIF
1839       
1840        IF (mpi_rank==rank) THEN
1841          ibegin = begin_i(i)
1842          ni = size_i(i)
1843          jbegin = begin_j(j)
1844          nj = size_j(j)
1845        ENDIF
1846        rank=rank+1
1847      ENDDO
1848    ENDDO
1849
1850    ALLOCATE(lon_glo(0:ni_glo-1), lat_glo(0:nj_glo-1))
1851    ALLOCATE(bounds_lon_glo(4,0:ni_glo-1), bounds_lat_glo(4,0:nj_glo-1))
1852   
1853    DO i=0,ni_glo-1
1854      lon_glo(i)=-180+(i+0.5)*(360./ni_glo)
1855      bounds_lon_glo(1,i)= -180+(i)*(360./ni_glo)
1856      bounds_lon_glo(2,i)= -180+(i+1)*(360./ni_glo)
1857    ENDDO
1858
1859    DO j=0,nj_glo-1
1860      lat_glo(j)=-90+(j+0.5)*(180./nj_glo)
1861      bounds_lat_glo(1,j)= -90+(j)*(180./nj_glo)
1862      bounds_lat_glo(2,j)= -90+(j+1)*(180./nj_glo)
1863    ENDDO
1864
1865    offset_i=2    ! halo of 2 on i
1866    offset_j=1    ! halo of 1 on j
1867
1868    ALLOCATE(lon(0:ni-1,0:nj-1))
1869    ALLOCATE(lat(0:ni-1,0:nj-1))
1870    ALLOCATE(bounds_lon(4,0:ni-1,0:nj-1))
1871    ALLOCATE(bounds_lat(4,0:ni-1,0:nj-1))
1872!    ALLOCATE(mask(0:ni-1,0:nj-1))
1873   
1874    ALLOCATE(return_lon(0:ni*nj-1))
1875    ALLOCATE(return_lat(0:ni*nj-1))
1876    ALLOCATE(return_mask(0:ni*nj-1))
1877    ALLOCATE(return_index(0:(ni+2*offset_i)*(nj+2*offset_j)-1))
1878
1879    ALLOCATE(return_X_lon(0:ni-1))
1880    ALLOCATE(return_X_lat(0:ni-1))
1881    ALLOCATE(return_X_mask(0:ni-1))
1882    ALLOCATE(return_X_index(0:ni-1))
1883
1884    ALLOCATE(return_Y_lon(0:nj-1))
1885    ALLOCATE(return_Y_lat(0:nj-1))
1886    ALLOCATE(return_Y_mask(0:nj-1))
1887    ALLOCATE(return_Y_index(0:nj-1))
1888   
1889    return_index=-1 
1890    DO j=0,nj-1
1891      DO i=0,ni-1
1892        ij=j*ni+i
1893        return_lon(ij)=lon_glo(ibegin+i)
1894        return_lat(ij)=lat_glo(jbegin+j)
1895        lon(i,j)=return_lon(ij)
1896        lat(i,j)=return_lat(ij)
1897        bounds_lon(1,i,j)=bounds_lon_glo(2,ibegin+i)
1898        bounds_lon(2,i,j)=bounds_lon_glo(1,ibegin+i)
1899        bounds_lon(3,i,j)=bounds_lon_glo(1,ibegin+i)
1900        bounds_lon(4,i,j)=bounds_lon_glo(2,ibegin+i)
1901        bounds_lat(1,i,j)=bounds_lat_glo(1,jbegin+j)
1902        bounds_lat(2,i,j)=bounds_lat_glo(1,jbegin+j)
1903        bounds_lat(3,i,j)=bounds_lat_glo(2,jbegin+j)
1904        bounds_lat(4,i,j)=bounds_lat_glo(2,jbegin+j)
1905
1906        ij=(j+offset_j)*(ni+2*offset_i)+i+offset_i
1907        return_index(ij)=i+j*ni
1908      ENDDO
1909    ENDDO
1910
1911    CALL set_domain_mask(params, return_lon,return_lat, return_mask)
1912
1913    ALLOCATE(return_X_lon(0:ni-1))
1914    ALLOCATE(return_X_lat(0:ni-1))
1915    ALLOCATE(return_X_mask(0:ni-1))
1916    ALLOCATE(return_X_index(0:ni+2*offset_i-1))
1917
1918    return_X_lat(:)=lat_glo(nj_glo/2)
1919    return_X_lon(:)=lon_glo(ibegin:ibegin+ni-1)
1920    CALL set_domain_mask(params, return_X_lon,return_X_lat, return_X_mask)
1921
1922    return_X_index(:)=-1
1923    DO i=0,ni-1
1924      return_X_index(offset_i+i)=i
1925    ENDDO
1926
1927    ALLOCATE(return_Y_lon(0:nj-1))
1928    ALLOCATE(return_Y_lat(0:nj-1))
1929    ALLOCATE(return_Y_mask(0:nj-1))
1930    ALLOCATE(return_Y_index(0:nj+2*offset_j-1))
1931
1932    return_Y_lat(:)=lat_glo(jbegin:jbegin+nj-1)
1933    return_Y_lon(:)=lon_glo(ni_glo/2)
1934    CALL set_domain_mask(params, return_Y_lon,return_Y_lat, return_Y_mask)
1935
1936    return_Y_index(:)=-1
1937    DO j=0,nj-1
1938      return_Y_index(offset_j+j)=j
1939    ENDDO
1940
1941    return_ni=ni
1942    return_nj=nj
1943
1944    IF (xios_is_valid_domain(TRIM(domain_id))) THEN
1945      CALL xios_set_domain_attr(TRIM(domain_id), type="curvilinear", data_dim=2)
1946      CALL xios_set_domain_attr(TRIM(domain_id), ni_glo=ni_glo, ibegin=ibegin, ni=ni, data_ibegin=-offset_i, data_ni=ni+2*offset_i)
1947      CALL xios_set_domain_attr(TRIM(domain_id), nj_glo=nj_glo, jbegin=jbegin, nj=nj, data_jbegin=-offset_j, data_nj=nj+2*offset_j)
1948      CALL xios_set_domain_attr(TRIM(domain_id), data_dim=2, lonvalue_2d=lon, latvalue_2d=lat, mask_1d=return_mask)
1949      CALL xios_set_domain_attr(TRIM(domain_id), bounds_lon_2d=bounds_lon, bounds_lat_2d=bounds_lat, nvertex=4)
1950    ENDIF
1951
1952   
1953    IF (xios_is_valid_axis(TRIM(domain_id)//"_X")) THEN
1954      CALL xios_set_axis_attr(TRIM(domain_id)//"_X", n_glo=ni_glo, begin=ibegin, n=ni, data_begin=-offset_i, data_n=ni+2*offset_i, value=return_X_lon, mask=return_X_mask)
1955!      CALL xios_set_axis_attr(TRIM(domain_id)//"_X", n_glo=ni_glo, begin=ibegin, n=ni, data_index=return_X_index, data_n=ni+2*offset_i, value=return_X_lon)
1956    ENDIF
1957
1958    IF (xios_is_valid_axis(TRIM(domain_id)//"_Y")) THEN   
1959      CALL xios_set_axis_attr(TRIM(domain_id)//"_Y", n_glo=nj_glo, begin=jbegin, n=nj, data_begin=-offset_j, data_n=nj+2*offset_j, value=return_Y_lat, mask=return_Y_mask)
1960!      CALL xios_set_axis_attr(TRIM(domain_id)//"_Y", n_glo=nj_glo, begin=jbegin, n=nj, data_index=return_Y_index, data_n=nj+2*offset_j, value=return_Y_lat, mask=return_Y_mask)
1961    ENDIF
1962
1963  END SUBROUTINE init_domain_nemo
1964
1965
1966
1967
1968 
1969   SUBROUTINE set_domain_mask(params, lon,lat, mask)
1970   IMPLICIT NONE
1971     TYPE(tmodel_params) :: params
1972     DOUBLE PRECISION  :: lon(:)
1973     DOUBLE PRECISION  :: lat(:)
1974     LOGICAL           :: mask(:)
1975     INTEGER :: i,x
1976
1977     mask(:)=.TRUE.
1978     IF (params%domain_mask) THEN
1979       IF (params%domain_mask_type=="cross") THEN
1980         WHERE (lon(:)-2*lat(:)>-10 .AND. lon(:)-2*lat(:) <10) mask(:)=.FALSE. 
1981         WHERE (2*lat(:)+lon(:)>-10 .AND. 2*lat(:)+lon(:)<10) mask(:)=.FALSE.
1982       ELSE IF (params%domain_mask_type=="latitude_band") THEN
1983         WHERE (lat(:)>-30 .AND. lat(:)<30) mask(:)=.FALSE.
1984      ENDIF
1985     ENDIF
1986
1987  END SUBROUTINE set_domain_mask
1988     
1989
1990  SUBROUTINE set_mask3d(grid_id, params, ni, nj, lon,lat, axis_value)
1991  IMPLICIT NONE
1992    CHARACTER(LEN=*) :: grid_id
1993    TYPE(tmodel_params) :: params
1994    INTEGER :: comm
1995    INTEGER :: ni
1996    INTEGER :: nj
1997    DOUBLE PRECISION, POINTER :: lon(:)
1998    DOUBLE PRECISION, POINTER :: lat(:)
1999    INTEGER, POINTER          :: domain_index(:)   
2000    DOUBLE PRECISION, POINTER :: axis_value(:)
2001    INTEGER, POINTER          :: axis_index(:)
2002    INTEGER ::i,j,ij,k,nk
2003    LOGICAL, ALLOCATABLE :: mask3d(:,:,:)
2004    DOUBLE PRECISION :: r
2005
2006    nk=size(axis_value)
2007
2008    ALLOCATE(mask3d(0:ni-1,0:nj-1,0:nk-1))
2009
2010    mask3d=.TRUE.
2011    DO k=0,nk-1
2012      DO j=0,nj-1
2013        DO i=0,ni-1
2014          ij=j*ni+i
2015          r=sqrt((lon(ij)/2)**2 + lat(ij)**2) / ((nk-k)*1./nk) 
2016          if (r < 60) mask3d(i,j,k)=.FALSE.
2017        ENDDO
2018      ENDDO
2019    ENDDO
2020
2021    IF (params%mask3d) CALL  xios_set_grid_attr(grid_id, mask_3d=mask3d)
2022
2023  END SUBROUTINE set_mask3d
2024   
2025  SUBROUTINE init_axis_pressure(axis_id,comm,params, return_value, return_mask, return_index)
2026  IMPLICIT NONE
2027    CHARACTER(LEN=*) :: axis_id
2028    TYPE(tmodel_params) :: params
2029    INTEGER :: comm
2030    DOUBLE PRECISION, POINTER :: return_value(:)
2031    LOGICAL, POINTER          :: return_mask(:)
2032    INTEGER, POINTER          :: return_index(:)
2033   
2034    INTEGER :: nlev_glo
2035    INTEGER :: nlev, begin, end
2036    DOUBLE PRECISION, ALLOCATABLE :: value_glo(:) 
2037    DOUBLE PRECISION, ALLOCATABLE :: bounds_value_glo(:,:) 
2038    DOUBLE PRECISION, ALLOCATABLE :: value(:) 
2039    DOUBLE PRECISION, ALLOCATABLE :: bounds_value(:,:) 
2040    DOUBLE PRECISION :: dp
2041    INTEGER :: i
2042
2043    INTEGER :: domain_proc_rank, domain_proc_size   
2044    INTEGER :: axis_proc_rank, axis_proc_size
2045    INTEGER :: ensemble_proc_size, ensemble_proc_rank
2046
2047    CALL get_decomposition(comm, params, domain_proc_size, domain_proc_rank, axis_proc_size, axis_proc_rank, ensemble_proc_size, ensemble_proc_rank)     
2048
2049    nlev_glo=params%nlev
2050   
2051    ALLOCATE(value_glo(0:nlev_glo-1), bounds_value_glo(2,0:nlev_glo-1) )
2052   
2053    dp=(1-0.1)/(nlev_glo-1)
2054    DO i=0,nlev_glo-1
2055     value_glo(i)=1-i*dp
2056    ENDDO
2057   
2058    bounds_value_glo(2,0)=value_glo(0)-(value_glo(1)-value_glo(0))/2
2059    DO i=1,nlev_glo-1
2060     bounds_value_glo(2,i)=(value_glo(i-1)+value_glo(i)) /2
2061    ENDDO
2062
2063    DO i=0,nlev_glo-2
2064     bounds_value_glo(1,i)=(value_glo(i)+value_glo(i+1)) /2
2065    ENDDO
2066    bounds_value_glo(1,nlev_glo-1)=value_glo(nlev_glo-1)-(value_glo(nlev_glo-2)-value_glo(nlev_glo-1))/2
2067
2068    nlev=nlev_glo/axis_proc_size
2069    IF (axis_proc_rank < MOD(nlev_glo,axis_proc_size)) THEN
2070      nlev=nlev+1
2071      begin= axis_proc_rank*nlev
2072      end=begin+nlev-1
2073    ELSE
2074      begin=MOD(nlev_glo,axis_proc_size)*(nlev+1) + (axis_proc_rank-MOD(nlev_glo,axis_proc_size))*nlev
2075      end=begin+nlev-1
2076    ENDIF
2077
2078    ALLOCATE(value(0:nlev-1), bounds_value(2,0:nlev-1) )
2079    value(:)=value_glo(begin:end)
2080    bounds_value(:,:)=bounds_value_glo(:,begin:end)
2081   
2082    ALLOCATE(return_value(0:nlev-1))
2083    ALLOCATE(return_mask(0:nlev-1))
2084    return_value=value
2085    return_mask=.TRUE.
2086    CALL set_axis_mask(params,value,return_mask)   
2087    CALL xios_set_axis_attr(axis_id, n_glo=nlev_glo, begin=begin, n=nlev, value=value*100000, mask=return_mask, bounds=bounds_value*100000, unit='Pa', positive='up')   
2088   
2089
2090    ALLOCATE(return_index(0:nlev-1))
2091
2092    DO i=0,nlev-1
2093      return_index(i)=i
2094    ENDDO 
2095
2096
2097  END SUBROUTINE init_axis_pressure
2098
2099  SUBROUTINE set_axis_mask(params, value, mask)
2100  IMPLICIT NONE
2101    TYPE(tmodel_params) :: params
2102     DOUBLE PRECISION  :: value(:)
2103     LOGICAL           :: mask(:)
2104     INTEGER :: i,x
2105     
2106     mask(:)=.TRUE.
2107     x=size(mask)
2108     IF (params%axis_mask) THEN
2109       DO i=0,x-1
2110         IF (MOD(i,3)==0) mask(i+1)=.FALSE.
2111         IF (MOD(i,4)==0) mask(i+1)=.FALSE.
2112       ENDDO
2113     ENDIF
2114
2115  END SUBROUTINE set_axis_mask 
2116
2117  SUBROUTINE init_scalar(scalar_id, comm, params,  return_mask)
2118  IMPLICIT NONE
2119     CHARACTER(LEN=*) :: scalar_id
2120     TYPE(tmodel_params) :: params
2121     INTEGER             :: comm
2122     LOGICAL             :: return_mask
2123     DOUBLE PRECISION    :: value =10.
2124
2125    CALL set_scalar_mask(comm, params, return_mask)   
2126    CALL xios_set_scalar_attr(scalar_id, value=value, mask=return_mask) 
2127
2128  END SUBROUTINE init_scalar
2129
2130  SUBROUTINE set_scalar_mask(comm, params, mask)
2131   IMPLICIT NONE
2132     TYPE(tmodel_params) :: params
2133     INTEGER             :: comm
2134     LOGICAL             :: mask
2135     INTEGER             :: ierr,rank
2136
2137     mask=.TRUE.
2138     IF (params%scalar_mask) THEN
2139       IF (params%scalar_mask_type=="none") THEN
2140         mask=.TRUE.
2141       ELSE IF (params%scalar_mask_type=="full") THEN
2142         mask=.FALSE.
2143       ELSE IF (params%scalar_mask_type=="root") THEN
2144         CALL MPI_COMM_RANK(comm,rank,ierr)
2145         mask = (rank==0)
2146       ELSE IF (params%scalar_mask_type=="sparse") THEN
2147         CALL MPI_COMM_RANK(comm,rank,ierr)
2148         mask = (MOD(rank,2)==0)
2149      ENDIF
2150     ENDIF
2151
2152  END SUBROUTINE set_scalar_mask
2153
2154  SUBROUTINE init_field2D_academic(comm,params, lon, lat, mask, return_field,            &
2155                                   X_lon, X_lat, X_mask, return_fieldX,                  &
2156                                   Y_lon, Y_lat, Y_mask, return_fieldY, return_fieldXY)
2157  IMPLICIT NONE
2158    TYPE(tmodel_params) :: params
2159    INTEGER :: comm
2160    DOUBLE PRECISION, POINTER :: lon(:)
2161    DOUBLE PRECISION, POINTER :: lat(:)
2162    LOGICAL, POINTER :: mask(:)
2163    DOUBLE PRECISION, POINTER :: return_field(:)
2164
2165    DOUBLE PRECISION, POINTER :: X_lon(:)
2166    DOUBLE PRECISION, POINTER :: X_lat(:)
2167    LOGICAL, POINTER :: X_mask(:)
2168    DOUBLE PRECISION, POINTER :: return_fieldX(:)
2169
2170    DOUBLE PRECISION, POINTER :: Y_lon(:)
2171    DOUBLE PRECISION, POINTER :: Y_lat(:)
2172    LOGICAL, POINTER :: Y_mask(:)
2173    DOUBLE PRECISION, POINTER :: return_fieldY(:)
2174    DOUBLE PRECISION, POINTER :: return_fieldXY(:,:)
2175   
2176    DOUBLE PRECISION, PARAMETER :: coef=2., dp_pi=3.14159265359
2177    DOUBLE PRECISION :: dp_length, dp_conv
2178    INTEGER :: i,j,x,y,xy
2179   
2180     ! Parameters for analytical function
2181    dp_length= 1.2*dp_pi
2182    dp_conv=dp_pi/180.
2183   
2184    xy=size(mask)
2185    x=size(X_mask)
2186    y=size(Y_mask)
2187
2188    ALLOCATE(return_field(0:xy-1))
2189   
2190    DO i=0,xy-1
2191      IF (mask(i)) THEN
2192         return_field(i)=(coef-SIN(dp_pi*(ACOS(COS(lat(i)*dp_conv)*&
2193                   COS(lon(i)*dp_conv))/dp_length)))
2194      ENDIF
2195    ENDDO       
2196
2197
2198    ALLOCATE(return_fieldX(0:x-1))
2199   
2200    DO i=0,x-1
2201      IF (X_mask(i)) THEN
2202         return_fieldX(i)=(coef-SIN(dp_pi*(ACOS(COS(X_lat(i)*dp_conv)*&
2203                            COS(X_lon(i)*dp_conv))/dp_length)))
2204      ENDIF
2205    ENDDO           
2206
2207
2208    ALLOCATE(return_fieldY(0:y-1))
2209   
2210    DO i=0,y-1
2211      IF (Y_mask(i)) THEN
2212         return_fieldY(i)=(coef-SIN(dp_pi*(ACOS(COS(Y_lat(i)*dp_conv)*&
2213                            COS(Y_lon(i)*dp_conv))/dp_length)))
2214      ENDIF
2215    ENDDO
2216
2217    ALLOCATE(return_fieldXY(0:x-1,0:y-1))
2218   
2219    DO j=0,y-1
2220      DO i=0,x-1
2221        IF (Y_mask(j) .AND. X_mask(i)) THEN
2222           return_fieldXY(i,j)=(coef-SIN(dp_pi*(ACOS(COS(Y_lat(j)*dp_conv)*&
2223                                COS(X_lon(i)*dp_conv))/dp_length)))
2224        ENDIF
2225      ENDDO
2226    ENDDO   
2227       
2228  END SUBROUTINE init_field2D_academic
2229
2230
2231  SUBROUTINE init_field2D_constant(comm,params, lon, lat, mask, return_field,             &
2232                                   X_lon, X_lat, X_mask, return_fieldX,                   &
2233                                   Y_lon, Y_lat, Y_mask, return_fieldY, return_fieldXY)
2234  IMPLICIT NONE
2235    TYPE(tmodel_params) :: params
2236    INTEGER :: comm
2237    DOUBLE PRECISION, POINTER :: lon(:)
2238    DOUBLE PRECISION, POINTER :: lat(:)
2239    LOGICAL, POINTER :: mask(:)
2240    DOUBLE PRECISION, POINTER :: return_field(:)
2241   
2242    DOUBLE PRECISION, POINTER :: X_lon(:)
2243    DOUBLE PRECISION, POINTER :: X_lat(:)
2244    LOGICAL, POINTER :: X_mask(:)
2245    DOUBLE PRECISION, POINTER :: return_fieldX(:)
2246
2247    DOUBLE PRECISION, POINTER :: Y_lon(:)
2248    DOUBLE PRECISION, POINTER :: Y_lat(:)
2249    LOGICAL, POINTER :: Y_mask(:)
2250    DOUBLE PRECISION, POINTER :: return_fieldY(:)
2251    DOUBLE PRECISION, POINTER :: return_fieldXY(:,:)
2252    INTEGER :: x,y,xy
2253   
2254    xy=size(mask)
2255    x=size(X_mask)
2256    y=size(Y_mask)
2257       
2258    ALLOCATE(return_field(0:xy-1))
2259    return_field(:)=1
2260
2261    ALLOCATE(return_fieldX(0:x-1))
2262    return_fieldX=1
2263
2264    ALLOCATE(return_fieldY(0:y-1))
2265    return_fieldY=1
2266
2267    ALLOCATE(return_fieldXY(0:x-1,0:y-1))
2268    return_fieldXY=1
2269   
2270  END SUBROUTINE init_field2D_constant
2271
2272  SUBROUTINE init_field2D_rank(comm,params, lon, lat, mask, return_field,         &
2273                               X_lon, X_lat, X_mask, return_fieldX,               &
2274                               Y_lon, Y_lat, Y_mask, return_fieldY, return_fieldXY)
2275
2276  IMPLICIT NONE
2277    TYPE(tmodel_params) :: params
2278    INTEGER :: comm
2279    DOUBLE PRECISION, POINTER :: lon(:)
2280    DOUBLE PRECISION, POINTER :: lat(:)
2281    LOGICAL, POINTER :: mask(:)
2282    DOUBLE PRECISION, POINTER :: return_field(:)
2283
2284    DOUBLE PRECISION, POINTER :: X_lon(:)
2285    DOUBLE PRECISION, POINTER :: X_lat(:)
2286    LOGICAL, POINTER :: X_mask(:)
2287    DOUBLE PRECISION, POINTER :: return_fieldX(:)
2288
2289    DOUBLE PRECISION, POINTER :: Y_lon(:)
2290    DOUBLE PRECISION, POINTER :: Y_lat(:)
2291    LOGICAL, POINTER :: Y_mask(:)
2292    DOUBLE PRECISION, POINTER :: return_fieldY(:)
2293    DOUBLE PRECISION, POINTER :: return_fieldXY(:,:)   
2294    INTEGER ::  rank,ierr
2295    INTEGER :: x,y,xy
2296
2297    CALL MPI_COMM_RANK(comm,rank,ierr)
2298
2299   
2300   
2301    xy=size(mask)
2302    x=size(X_mask)
2303    y=size(Y_mask)
2304       
2305    ALLOCATE(return_field(0:xy-1))
2306    return_field(:)=rank
2307
2308    ALLOCATE(return_fieldX(0:x-1))
2309    return_fieldX=rank
2310
2311    ALLOCATE(return_fieldY(0:y-1))
2312    return_fieldY=rank
2313
2314    ALLOCATE(return_fieldXY(0:x-1,0:y-1))
2315    return_fieldXY=rank
2316
2317  END SUBROUTINE init_field2D_rank
2318
2319
2320
2321   SUBROUTINE get_decomposition(comm, params, domain_proc_size, domain_proc_rank, axis_proc_size, axis_proc_rank, ensemble_proc_size, ensemble_proc_rank)
2322   IMPLICIT NONE
2323     INTEGER,INTENT(IN) :: comm
2324     TYPE(tmodel_params) :: params
2325     INTEGER,INTENT(OUT) :: domain_proc_size
2326     INTEGER,INTENT(OUT) :: domain_proc_rank
2327     INTEGER,INTENT(OUT) :: axis_proc_size
2328     INTEGER,INTENT(OUT) :: axis_proc_rank
2329     INTEGER,INTENT(OUT) :: ensemble_proc_size
2330     INTEGER,INTENT(OUT) :: ensemble_proc_rank
2331
2332     INTEGER :: mpi_rank,mpi_size,rank, ensemble_number
2333     INTEGER :: ierr
2334     INTEGER :: n_domain,n_axis, n_ensemble, min_dist, new_dist, best
2335     INTEGER :: axis_ind, domain_ind
2336     LOGICAL :: axis_layer
2337
2338     CALL MPI_COMM_RANK(comm,mpi_rank,ierr)
2339     CALL MPI_COMM_SIZE(comm,mpi_size,ierr)
2340
2341     n_ensemble = params%ensemble_proc_n
2342     ensemble_proc_size = mpi_size / n_ensemble
2343     IF (  mpi_rank < MOD(mpi_size,n_ensemble) * (ensemble_proc_size+1) ) THEN
2344       ensemble_proc_size=ensemble_proc_size+1
2345       ensemble_proc_rank = MOD(mpi_rank,ensemble_proc_size)
2346       ensemble_number = mpi_rank / ensemble_proc_size
2347     ELSE
2348       ensemble_number = MOD(mpi_size,n_ensemble)
2349       rank =  mpi_rank - MOD(mpi_size,n_ensemble) * (ensemble_proc_size+1)
2350       ensemble_proc_rank= MOD(rank,ensemble_proc_size)
2351       ensemble_number = ensemble_number + rank / ensemble_proc_size
2352     ENDIF
2353
2354     mpi_size=ensemble_proc_size
2355     mpi_rank=ensemble_proc_rank
2356     ensemble_proc_size = n_ensemble
2357     ensemble_proc_rank = ensemble_number
2358
2359    IF (params%axis_proc_n > 0 ) THEN
2360      n_axis=params%axis_proc_n
2361      n_domain = mpi_size / n_axis
2362      axis_layer=.TRUE.
2363    ELSE IF (params%domain_proc_n > 0 ) THEN
2364      n_domain=params%domain_proc_n
2365      n_axis = mpi_size / n_domain
2366      axis_layer=.FALSE.
2367    ELSE
2368      IF (params%axis_proc_frac==0) THEN
2369         params%axis_proc_frac=1
2370         params%domain_proc_frac=mpi_size
2371      ELSE IF (params%domain_proc_frac==0) THEN
2372         params%domain_proc_frac=1
2373         params%axis_proc_frac=mpi_size
2374      ENDIF       
2375   
2376      n_domain = INT(sqrt(params%domain_proc_frac * mpi_size/ params%axis_proc_frac))
2377      n_axis =   INT(sqrt(params%axis_proc_frac * mpi_size/ params%domain_proc_frac))
2378
2379
2380      min_dist= mpi_size - n_domain*n_axis
2381      best=0
2382   
2383      new_dist = mpi_size -(n_domain+1)*n_axis
2384      IF (new_dist < min_dist .AND. new_dist >= 0 ) THEN
2385         min_dist=new_dist
2386         best=1
2387      ENDIF
2388
2389      new_dist=mpi_size-n_domain*(n_axis+1)
2390      IF (new_dist < min_dist .AND. new_dist >= 0 ) THEN
2391         min_dist=new_dist
2392         best=2
2393      ENDIF
2394
2395      IF (best==0) THEN
2396      ELSE IF (best==1) THEN
2397        n_domain=n_domain+1
2398      ELSE IF (best==2) THEN
2399        n_axis=n_axis+1
2400      ENDIF
2401
2402      IF ( MOD(mpi_size,n_axis) <= MOD(mpi_size,n_domain)) axis_layer=.TRUE.
2403
2404    ENDIF
2405   
2406    IF ( axis_layer) THEN
2407      !!! n_axis layer
2408      IF (mpi_rank < MOD(mpi_size,n_axis)*(n_domain+1)) THEN
2409        axis_proc_rank=mpi_rank/(n_domain+1)
2410        domain_proc_rank=MOD(mpi_rank,n_domain+1)
2411        axis_proc_size=n_axis
2412        domain_proc_size=n_domain+1
2413      ELSE
2414        rank=mpi_rank-MOD(mpi_size,n_axis)*(n_domain+1)
2415        axis_proc_size=n_axis
2416        axis_proc_rank=MOD(mpi_size,n_axis)+rank/n_domain
2417        domain_proc_rank=MOD(rank,n_domain)
2418        domain_proc_size=n_domain
2419      ENDIF
2420    ELSE
2421      !!! n_domain column
2422      IF (mpi_rank < MOD(mpi_size,n_domain)*(n_axis+1)) THEN
2423        domain_proc_rank=mpi_rank/(n_axis+1)
2424        axis_proc_rank=MOD(mpi_rank,n_axis+1)
2425        domain_proc_size=n_domain
2426        axis_proc_size=n_axis+1
2427      ELSE
2428        rank=mpi_rank-MOD(mpi_size,n_domain)*(n_axis+1)
2429        domain_proc_size=n_domain
2430        domain_proc_rank=MOD(mpi_size,n_domain)+rank/n_axis
2431        axis_proc_rank=MOD(rank,n_axis)
2432        axis_proc_size=n_axis
2433      ENDIF
2434    ENDIF 
2435
2436
2437  END SUBROUTINE get_decomposition
2438
2439  SUBROUTINE compute_cyclic_distribution(rank, nb_procs, who_i_am)
2440  IMPLICIT NONE
2441    INTEGER,INTENT(IN)  :: rank
2442    INTEGER,INTENT(IN)  :: nb_procs(:)
2443    LOGICAL,INTENT(OUT) :: who_i_am(:)
2444    INTEGER             :: start(SIZE(nb_procs))
2445    INTEGER :: nb_comp 
2446    INTEGER :: nb_client
2447    INTEGER :: nb_server
2448    INTEGER :: current
2449    INTEGER :: current_client
2450    INTEGER :: i,j,k,nq,q,r
2451
2452    who_i_am=.FALSE.
2453    nb_comp = SIZE(nb_procs)
2454    nb_client = SUM(nb_procs(1:nb_comp-1))
2455    nb_server = nb_procs(nb_comp)
2456   
2457    start(1)=0
2458    DO k=2,nb_comp
2459     start(k)=start(k-1)+nb_procs(k-1)
2460    ENDDO
2461
2462    current=0
2463    current_client=0
2464
2465    IF (nb_client >= nb_server) THEN
2466      q=nb_client/nb_server
2467      r=MOD(nb_client,nb_server)
2468       
2469      DO i=1,nb_server
2470        IF (i<=r) THEN ; nq=q+1 ; ELSE ; nq=q ; ENDIF
2471        DO j=1,nq
2472          IF (current==rank) THEN
2473            DO k=1,nb_comp-1
2474              IF (current_client >= start(k) .AND. current_client<start(k+1)) THEN
2475                who_i_am(k)=.TRUE.
2476                RETURN
2477              ENDIF 
2478            ENDDO
2479          ENDIF
2480          current=current+1
2481          current_client=current_client+1
2482        ENDDO
2483
2484        IF (current==rank) THEN
2485          who_i_am(nb_comp)=.TRUE.
2486          RETURN
2487        ENDIF
2488        current=current+1
2489      ENDDO
2490 
2491    ELSE
2492      q=nb_server/nb_client
2493      r=MOD(nb_server,nb_client)
2494
2495      DO i=1,nb_client
2496        IF (i<=r) THEN ; nq=q+1 ; ELSE ; nq=q ; ENDIF
2497        DO j=1,nq
2498          IF (current==rank) THEN
2499            who_i_am(nb_comp)=.TRUE.
2500            RETURN
2501          ENDIF
2502          current=current+1
2503        ENDDO
2504         
2505        IF (current==rank) THEN
2506          DO k=1,nb_comp-1
2507            IF (current_client >= start(k) .AND. current_client<start(k+1)) THEN
2508              who_i_am(k)=.TRUE.
2509              RETURN
2510            ENDIF 
2511          ENDDO
2512        ENDIF
2513        current=current+1
2514        current_client=current_client+1
2515      ENDDO
2516    ENDIF
2517
2518  END SUBROUTINE compute_cyclic_distribution   
2519
2520END PROGRAM generic_testcase
2521
2522
Note: See TracBrowser for help on using the repository browser.