/[MITgcm]/MITgcm/pkg/ctrl/ctrl_set_pack_xz.F
ViewVC logotype

Diff of /MITgcm/pkg/ctrl/ctrl_set_pack_xz.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph | View Patch Patch

revision 1.2 by heimbach, Sat Jul 13 02:47:32 2002 UTC revision 1.3 by heimbach, Tue Jun 24 16:07:07 2003 UTC
# Line 3  Line 3 
3    
4    
5        subroutine ctrl_set_pack_xz(        subroutine ctrl_set_pack_xz(
6       &     cunit, ivartype, fname, masktype,       &     cunit, ivartype, fname, masktype,weighttype,
7       &     weightfld, lxxadxx, mythid)       &     weightfld, lxxadxx, mythid)
8    
9  c     ==================================================================  c     ==================================================================
# Line 13  c Line 13  c
13  c     o Compress the control vector such that only ocean points are  c     o Compress the control vector such that only ocean points are
14  c       written to file.  c       written to file.
15  c  c
16    c     o Open boundary packing finalized :
17    c          gebbie@mit.edu, 18-Mar-2003
18    c
19    c     changed: heimbach@mit.edu 17-Jun-2003
20    c              merged Armin's changes to replace write of
21    c              nr * globfld2d by 1 * globfld3d
22    c              (ad hoc fix to speed up global I/O)
23    c
24  c     ==================================================================  c     ==================================================================
25    
26        implicit none        implicit none
# Line 37  c     == routine arguments == Line 45  c     == routine arguments ==
45        integer ivartype        integer ivartype
46        character*( 80) fname        character*( 80) fname
47        character*(  9) masktype        character*(  9) masktype
48          character*( 80) weighttype
49        _RL     weightfld( nr,nobcs )        _RL     weightfld( nr,nobcs )
50        logical lxxadxx        logical lxxadxx
51        integer mythid        integer mythid
# Line 52  c     == local variables == Line 61  c     == local variables ==
61        integer i,j,k        integer i,j,k
62        integer ii        integer ii
63        integer il        integer il
64        integer irec,iobcs        integer irec,iobcs,nrec_nl
65        integer itlo,ithi        integer itlo,ithi
66        integer jtlo,jthi        integer jtlo,jthi
67        integer jmin,jmax        integer jmin,jmax
# Line 62  c     == local variables == Line 71  c     == local variables ==
71  cgg(  cgg(
72        integer igg        integer igg
73        _RL     gg        _RL     gg
74          character*(80) weightname
75  cgg)  cgg)
76    
77        _RL     cbuff      ( snx*nsx*npx*nsy*npy )        _RL     cbuff      ( snx*nsx*npx*nsy*npy )
       _RL     globmskxz  ( snx,nsx,npx,nsy,npy,nr )  
78        _RL     globfldxz  ( snx,nsx,npx,nsy,npy,nr )        _RL     globfldxz  ( snx,nsx,npx,nsy,npy,nr )
79          _RL     globfld3d  ( snx,nsx,npx,sny,nsy,npy,nr )
80          _RL     globmskxz  ( snx,nsx,npx,nsy,npy,nr,nobcs )
81    #ifdef CTRL_PACK_PRECISE
82          _RL     weightfldxz( snx,nsx,npx,nsy,npy,nr,nobcs )
83    #endif
84    
85  c     == external ==  c     == external ==
86    
# Line 96  c     Initialise temporary file Line 110  c     Initialise temporary file
110                    do bi = itlo,ithi                    do bi = itlo,ithi
111                       do i = imin,imax                       do i = imin,imax
112                          globfldxz(i,bi,ip,bj,jp,k) = 0. _d 0                          globfldxz(i,bi,ip,bj,jp,k) = 0. _d 0
113                          globmskxz(i,bi,ip,bj,jp,k) = 0. _d 0                          do iobcs=1,nobcs
114                               globmskxz(i,bi,ip,bj,jp,k,iobcs) = 0. _d 0
115                            enddo
116                         enddo
117                      enddo
118                   enddo
119                enddo
120             enddo
121          enddo
122    c     Initialise temporary file
123          do k = 1,nr
124             do jp = 1,nPy
125                do bj = jtlo,jthi
126                   do j = jmin,jmax
127                      do ip = 1,nPx
128                         do bi = itlo,ithi
129                            do i = imin,imax
130                               globfld3d(i,bi,ip,j,bj,jp,k) = 0. _d 0
131                            enddo
132                       enddo                       enddo
133                    enddo                    enddo
134                 enddo                 enddo
# Line 107  c     Initialise temporary file Line 139  c     Initialise temporary file
139  c--   Only the master thread will do I/O.  c--   Only the master thread will do I/O.
140        _BEGIN_MASTER( mythid )        _BEGIN_MASTER( mythid )
141    
142        do irec = 1, ncvarrecs(ivartype)        do iobcs = 1, nobcs
143             call MDSREADFIELD_XZ_GL(
144         &        masktype, ctrlprec, 'RL',
145         &        Nr, globmskxz(1,1,1,1,1,1,iobcs), iobcs, mythid)
146    #ifdef CTRL_PACK_PRECISE
147             il=ilnblnk( weighttype)
148             write(weightname(1:80),'(80a)') ' '
149             write(weightname(1:80),'(a)') weighttype(1:il)
150             call MDSREADFIELD_XZ_GL(
151         &        weightname, ctrlprec, 'RL',
152         &        Nr, weightfldxz(1,1,1,1,1,1,iobcs), iobcs, mythid)
153    CGG   One special exception: barotropic velocity should be nondimensionalized
154    cgg   differently. Probably introduce new variable.
155             if (iobcs .eq. 3 .or. iobcs .eq. 4) then
156                k = 1
157                do jp = 1,nPy
158                   do bj = jtlo,jthi
159                      do ip = 1,nPx
160                         do bi = itlo,ithi
161                            do i = imin,imax
162                               weightfldxz(i,bi,ip,bj,jp,k,iobcs) = wbaro
163                            enddo
164                         enddo
165                      enddo
166                   enddo
167                enddo
168             endif
169    #endif
170          enddo
171    
172          nrec_nl=int(ncvarrecs(ivartype)/sny)
173          do irec = 1, nrec_nl
174             call MDSREADFIELD_3D_GL( fname, ctrlprec, 'RL',
175         &        nr, globfld3d, irec, mythid)
176             do j=1,sny
177                iobcs= mod((irec-1)*sny+j-1,nobcs)+1
178    
179    CGG   One special exception: barotropic velocity should be nondimensionalized
180    cgg   differently. Probably introduce new variable.
181                if (iobcs .eq. 3 .or. iobcs .eq. 4) then
182                   k = 1
183                   do jp = 1,nPy
184                      do bj = jtlo,jthi
185                         do ip = 1,nPx
186                            do bi = itlo,ithi
187                               do i = imin,imax
188    #ifdef NO_CONTROL_BAROTROPIC_VELOCITY
189                                  if (.not. lxxadxx) then
190    cgg    Get rid of any sensitivity to barotropic velocity.
191                                     globfld3d(i,bi,ip,j,bj,jp,k) = 0.
192                                  endif
193    #endif
194                               enddo
195                            enddo
196                         enddo
197                      enddo
198                   enddo
199                endif
200    
201                write(cunit) ncvarindex(ivartype)
202                write(cunit) 1
203                write(cunit) 1
204                do k = 1,nr
205                 cbuffindex = 0
206                 do jp = 1,nPy
207                  do bj = jtlo,jthi
208                   do ip = 1,nPx
209                    do bi = itlo,ithi
210                     do i = imin,imax
211                      if (globmskxz(i,bi,ip,bj,jp,k,iobcs)  .ne. 0. ) then
212                         cbuffindex = cbuffindex + 1
213    #ifdef ALLOW_NONDIMENSIONAL_CONTROL_IO
214                         if (lxxadxx) then
215                            cbuff(cbuffindex) =
216         &                       globfld3d(i,bi,ip,j,bj,jp,k) *
217    # ifdef CTRL_PACK_PRECISE
218         &                       sqrt(weightfldxz(i,bi,ip,bj,jp,k,iobcs))
219    # else
220         &                       sqrt(weightfld(k,iobcs))
221    # endif
222                         else
223                            cbuff(cbuffindex) =
224         &                       globfld3d(i,bi,ip,j,bj,jp,k) /
225    # ifdef CTRL_PACK_PRECISE
226         &                       sqrt(weightfldxz(i,bi,ip,bj,jp,k,iobcs))
227    # else
228         &                       sqrt(weightfld(k,iobcs))
229    # endif
230                         endif
231    #else /* ALLOW_NONDIMENSIONAL_CONTROL_IO undef */
232                         cbuff(cbuffindex) = globfld3d(i,bi,ip,j,bj,jp,k)
233    #endif /* ALLOW_NONDIMENSIONAL_CONTROL_IO */
234                      endif
235                     enddo
236                    enddo
237                   enddo
238                  enddo
239                 enddo
240    c           --> check cbuffindex.
241                 if ( cbuffindex .gt. 0) then
242                    write(cunit) cbuffindex
243                    write(cunit) k
244                    write(cunit) (cbuff(ii), ii=1,cbuffindex)
245                 endif
246    c     -- end of k loop --
247                enddo
248    c     -- end of j loop --
249             enddo
250    c     -- end of irec loop --
251          enddo
252    
253          do irec = nrec_nl*sny+1, ncvarrecs(ivartype)
254  cgg       do iobcs = 1, nobcs  cgg       do iobcs = 1, nobcs
255  cgg    Need to solve for what iobcs would have been.  cgg    Need to solve for what iobcs would have been.
256            gg   = (irec-1)/nobcs           iobcs= mod(irec-1,nobcs)+1
           igg  = int(gg)  
           iobcs = irec - igg*nobcs  
   
          call MDSREADFIELD_XZ_GL(  
      &        masktype, ctrlprec, 'RL',  
      &        Nr, globmskxz, iobcs, mythid)  
257    
258           call MDSREADFIELD_XZ_GL( fname, ctrlprec, 'RL',           call MDSREADFIELD_XZ_GL( fname, ctrlprec, 'RL',
259       &        nr, globfldxz, irec, mythid)       &        nr, globfldxz, irec, mythid)
260    
261    CGG   One special exception: barotropic velocity should be nondimensionalized
262    cgg   differently. Probably introduce new variable.
263             if (iobcs .eq. 3 .or. iobcs .eq. 4) then
264                k = 1
265                do jp = 1,nPy
266                   do bj = jtlo,jthi
267                      do ip = 1,nPx
268                         do bi = itlo,ithi
269                            do i = imin,imax
270    #ifdef NO_CONTROL_BAROTROPIC_VELOCITY
271                               if (.not. lxxadxx) then
272    cgg    Get rid of any sensitivity to barotropic velocity.
273                                  globfldxz(i,bi,ip,bj,jp,k) = 0.
274                               endif
275    #endif
276                            enddo
277                         enddo
278                      enddo
279                   enddo
280                enddo
281             endif
282    
283           write(cunit) ncvarindex(ivartype)           write(cunit) ncvarindex(ivartype)
284           write(cunit) 1           write(cunit) 1
285           write(cunit) 1           write(cunit) 1
# Line 131  cgg    Need to solve for what iobcs woul Line 290  cgg    Need to solve for what iobcs woul
290                do ip = 1,nPx                do ip = 1,nPx
291                 do bi = itlo,ithi                 do bi = itlo,ithi
292                  do i = imin,imax                  do i = imin,imax
293                   if (globmskxz(i,bi,ip,bj,jp,k)  .ne. 0. ) then                   if (globmskxz(i,bi,ip,bj,jp,k,iobcs)  .ne. 0. ) then
294                       cbuffindex = cbuffindex + 1                       cbuffindex = cbuffindex + 1
295  #ifdef ALLOW_NONDIMENSIONAL_CONTROL_IO  #ifdef ALLOW_NONDIMENSIONAL_CONTROL_IO
296                       if (lxxadxx) then                       if (lxxadxx) then
297                          cbuff(cbuffindex) =                          cbuff(cbuffindex) =
298       &                       globfldxz(i,bi,ip,bj,jp,k) *       &                       globfldxz(i,bi,ip,bj,jp,k) *
299    # ifdef CTRL_PACK_PRECISE
300         &                       sqrt(weightfldxz(i,bi,ip,bj,jp,k,iobcs))
301    # else
302       &                       sqrt(weightfld(k,iobcs))       &                       sqrt(weightfld(k,iobcs))
303    # endif
304                       else                       else
305                          cbuff(cbuffindex) =                          cbuff(cbuffindex) =
306       &                       globfldxz(i,bi,ip,bj,jp,k) /       &                       globfldxz(i,bi,ip,bj,jp,k) /
307    # ifdef CTRL_PACK_PRECISE
308         &                       sqrt(weightfldxz(i,bi,ip,bj,jp,k,iobcs))
309    # else
310       &                       sqrt(weightfld(k,iobcs))       &                       sqrt(weightfld(k,iobcs))
311    # endif
312                       endif                       endif
313  #else  #else
314                       cbuff(cbuffindex) = globfldxz(i,bi,ip,bj,jp,k)                       cbuff(cbuffindex) = globfldxz(i,bi,ip,bj,jp,k)
# Line 158  c           --> check cbuffindex. Line 325  c           --> check cbuffindex.
325                 write(cunit) k                 write(cunit) k
326                 write(cunit) (cbuff(ii), ii=1,cbuffindex)                 write(cunit) (cbuff(ii), ii=1,cbuffindex)
327              endif              endif
          enddo  
328  c  c
329    c     -- end of k loop --
330             enddo
331  c     -- end of iobcs loop --  c     -- end of iobcs loop --
332  cgg       enddo  cgg       enddo
333  c     -- end of irec loop --  c     -- end of irec loop --

Legend:
Removed from v.1.2  
changed lines
  Added in v.1.3

  ViewVC Help
Powered by ViewVC 1.1.22