source: CPL/oasis3/trunk/src/lib/mpp_io/include/mpp_write_2Ddecomp.h

Last change on this file was 1677, checked in by aclsce, 12 years ago

Imported oasis3 (tag ipslcm5a) from cvs server to svn server (igcmg project).

File size: 4.3 KB
Line 
1subroutine MPP_WRITE_2DDECOMP_1D_( unit, field, domain, data, tstamp )
2      integer, intent(in) :: unit
3      type(fieldtype), intent(in) :: field
4      type(domain2D), intent(inout) :: domain
5      MPP_TYPE_, intent(inout) :: data(:)
6      real(DOUBLE_KIND), intent(in), optional :: tstamp
7      integer::il_xbegin,il_xend,il_ybegin,il_yend,ij
8      MPP_TYPE_,allocatable :: data3D(:,:,:)
9      call mpp_get_compute_domain(domain &
10                                ,xbegin=il_xbegin,xend=il_xend &
11                                ,ybegin=il_ybegin,yend=il_yend)
12      allocate(data3D(1:il_xend-il_xbegin+1,1:il_yend-il_ybegin+1,1:1))
13      data3D = RESHAPE( data, SHAPE(data3D) )
14      call mpp_write( unit, field, domain, data3D, tstamp )
15      data   = RESHAPE( data3D, SHAPE(data) )
16      deallocate(data3D)
17
18      return
19    end subroutine MPP_WRITE_2DDECOMP_1D_
20
21    subroutine MPP_WRITE_2DDECOMP_2D_( unit, field, domain, data, tstamp )
22      integer, intent(in) :: unit
23      type(fieldtype), intent(in) :: field
24      type(domain2D), intent(inout) :: domain
25      MPP_TYPE_, intent(inout) :: data(:,:)
26      real(DOUBLE_KIND), intent(in), optional :: tstamp
27      MPP_TYPE_ :: data3D(size(data,1),size(data,2),1)
28#ifdef use_CRI_pointers
29      pointer( ptr, data3D )
30      ptr = LOC(data)
31      call mpp_write( unit, field, domain, data3D, tstamp )
32#else
33      data3D = RESHAPE( data, SHAPE(data3D) )
34      call mpp_write( unit, field, domain, data3D, tstamp )
35      data   = RESHAPE( data3D, SHAPE(data) )
36#endif
37      return
38    end subroutine MPP_WRITE_2DDECOMP_2D_
39
40    subroutine MPP_WRITE_2DDECOMP_3D_( unit, field, domain, data, tstamp )
41!mpp_write writes <data> which has the domain decomposition <domain>
42      integer, intent(in) :: unit
43      type(fieldtype), intent(in) :: field
44      type(domain2D), intent(inout) :: domain !must have intent(out) as well because active domain might be reset
45      MPP_TYPE_, intent(inout) :: data(:,:,:)
46      real(DOUBLE_KIND), intent(in), optional :: tstamp
47!cdata is used to store compute domain as contiguous data
48!gdata is used to globalize data for multi-PE single-threaded I/O
49      MPP_TYPE_, allocatable, dimension(:,:,:) :: cdata, gdata
50!NEW: data may be on compute OR data domain
51      logical :: data_has_halos, halos_are_global, x_is_global, y_is_global
52      integer :: is, ie, js, je, isd, ied, jsd, jed, isg, ieg, jsg, jeg
53
54      if( .NOT.module_is_initialized )call mpp_error( FATAL, 'MPP_WRITE: must first call mpp_io_init.' )
55      if( .NOT.mpp_file(unit)%opened )call mpp_error( FATAL, 'MPP_WRITE: invalid unit number.' )
56      call mpp_get_compute_domain( domain, is,  ie,  js,  je  )
57      call mpp_get_data_domain   ( domain, isd, ied, jsd, jed, x_is_global=x_is_global, y_is_global=y_is_global )
58      call mpp_get_global_domain ( domain, isg, ieg, jsg, jeg )
59      if( size(data,1).EQ.ie-is+1 .AND. size(data,2).EQ.je-js+1 )then
60          data_has_halos = .FALSE.
61      else if( size(data,1).EQ.ied-isd+1 .AND. size(data,2).EQ.jed-jsd+1 )then
62          data_has_halos = .TRUE.
63      else
64          call mpp_error( FATAL, 'MPP_WRITE: data must be either on compute domain or data domain.' )
65      end if
66      halos_are_global = x_is_global .AND. y_is_global
67      if( npes.GT.1 .AND. mpp_file(unit)%threading.EQ.MPP_SINGLE )then
68          if( halos_are_global )then
69              call mpp_update_domains( data, domain )
70!all non-0 PEs have passed their data to PE 0 and may now exit
71              if( pe.NE.0 )return
72              call write_record( unit, field, size(data), data, tstamp )
73          else
74!put field onto global domain
75              allocate( gdata(isg:ieg,jsg:jeg,size(data,3)) )
76              call mpp_global_field( domain, data, gdata )
77!all non-0 PEs have passed their data to PE 0 and may now exit
78              if( pe.NE.0 )return
79              call write_record( unit, field, size(gdata), gdata, tstamp )
80          end if
81      else if( data_has_halos )then
82!store compute domain as contiguous data and pass to write_record
83          allocate( cdata(is:ie,js:je,size(data,3)) )
84          cdata(:,:,:) = data(is-isd+1:ie-isd+1,js-jsd+1:je-jsd+1,:)
85          call write_record( unit, field, size(cdata), cdata, tstamp, domain )
86      else
87!data is already contiguous
88          call write_record( unit, field, size(data), data, tstamp, domain )
89      end if
90      return
91    end subroutine MPP_WRITE_2DDECOMP_3D_
Note: See TracBrowser for help on using the repository browser.