psmile_write_4d_log.F90

Go to the documentation of this file.
00001 !-----------------------------------------------------------------------
00002 ! Copyright 2006-2010, SGI Germany, Munich, Germany.
00003 ! All rights reserved. Use is subject to OASIS4 license terms.
00004 !-----------------------------------------------------------------------
00005 !BOP
00006 !
00007 ! !ROUTINE:   psmile_write_2d_dble
00008 !
00009 ! !INTERFACE:
00010 
00011       subroutine psmile_write_4d_log(unit                 &
00012                                      ,var                 &
00013                                      ,domain              &
00014                                      ,a                   &
00015                                      ,a_shape             &
00016                                      ,v_shape             &
00017                                      ,time                &
00018                                      ,time_used           &
00019                                      ,id_blockid          &
00020                                      ,block_used          &
00021                                      ,ierror)
00022 !
00023 ! !USES:
00024 !
00025 
00026       use psmile
00027 #ifdef __PSMILE_WITH_IO
00028       implicit none
00029       include 'prism.inc'
00030 !
00031 ! !INPUT PARAMETERS:
00032 !
00033       integer,intent(in) :: unit
00034       Type(fieldtype),intent(in) :: var
00035       Type(domain2D),intent(inout)  :: domain
00036       integer,intent(in) :: a_shape(2,*)
00037       integer,intent(in) :: v_shape(2,*)
00038       logical            , intent(inout) :: a(a_shape(1,1):a_shape(2,1)  
00039                                        ,a_shape(1,2):a_shape(2,2)  
00040                                        ,a_shape(1,3):a_shape(2,3)  
00041                                        ,a_shape(1,4):a_shape(2,4)) 
00042       double precision, intent(in) :: time
00043       logical,intent(in) :: time_used
00044       integer,intent(in) :: id_blockid
00045       logical,intent(in) :: block_used
00046 !
00047 ! !OUTPUT PARAMETERS:
00048 !
00049       integer, intent(out) :: ierror
00050 !
00051 ! !LOCAL VARIABLES:
00052 !
00053       integer :: ierrp(2)
00054       integer :: i,j,k,l
00055       double precision,allocatable    :: atmp(:,:,:,:)
00056 !
00057 ! !DESCRIPTION:
00058 !
00059 ! Writes a logical partitioned 4d array to a file.
00060 ! Subblocks are supported.
00061 !
00062 ! !REVISION HISTORY:
00063 !
00064 !   Date      Programmer   Description
00065 ! ----------  ----------   -----------
00066 !  09.12.03   R. Vogelsang created
00067 !  6.05.04    R. Vogelsang bugfixes a=atmp
00068 !
00069 !EOP
00070 !----------------------------------------------------------------------
00071 
00072       character(len=len_cvs_string),save :: mcvs = 
00073 '$Id: psmile_write_4d_log.F90 2325 2010-04-21 15:00:07Z valcke $'
00074       
00075       ierror=0
00076 
00077 #ifdef VERBOSE
00078       print*,trim(ch_id),' :   psmile_write_4d_log: start'
00079       print*,trim(ch_id),' :   psmile_write_4d_log: size ',size(a)
00080       call psmile_flushstd
00081 
00082 #endif
00083 
00084       allocate( atmp(a_shape(1,1):a_shape(2,1)  &
00085                     ,a_shape(1,2):a_shape(2,2)  &
00086                     ,a_shape(1,3):a_shape(2,3)  &
00087                     ,a_shape(1,4):a_shape(2,4)),stat=ierror)
00088 
00089       if ( ierror /= 0 ) then
00090          ierrp (1) = 1
00091          ierror = PRISM_Error_Alloc
00092          call psmile_error ( ierror, 'atmp', ierrp, 1, __FILE__, __LINE__ )
00093          return
00094       endif
00095 
00096 !
00097 !     Workaround. One can not write directly logical as integers.
00098 !     Internal representation of logicals are compiler dependent.
00099 !
00100       do l=a_shape(1,4),a_shape(2,4)
00101         do k=a_shape(1,3),a_shape(2,3)
00102           do j=a_shape(1,2),a_shape(2,2)
00103             do i=a_shape(1,1),a_shape(2,1)
00104               if(a(i,j,k,l)) then
00105                 atmp(i,j,k,l)=1
00106               else
00107                 atmp(i,j,k,l)=0
00108               endif
00109             enddo
00110           enddo
00111         enddo
00112       enddo
00113 
00114       if(block_used) then
00115       if(id_blockid.le.0) then
00116          ierror = PRISM_Error_Internal
00117          call psmile_error ( ierror, 'id_blockid <= 0! ', &
00118                              ierrp, 0, __FILE__, __LINE__ )
00119       endif
00120       if(time_used) then
00121         call mpp_write(unit,var,domain,atmp(v_shape(1,1):v_shape(2,1) &
00122                                              ,v_shape(1,2):v_shape(2,2) &
00123                                              ,v_shape(1,3):v_shape(2,3) &
00124                                              ,v_shape(1,4):v_shape(2,4)) &
00125                                       ,tstamp=time,blockid=id_blockid)
00126       else
00127         call mpp_write( unit,var,domain,atmp(v_shape(1,1):v_shape(2,1) &
00128                                               ,v_shape(1,2):v_shape(2,2) &
00129                                              ,v_shape(1,3):v_shape(2,3) &
00130                                               ,v_shape(1,4):v_shape(2,4)) &
00131                                               , blockid=id_blockid)
00132       endif
00133       else
00134       if(time_used) then
00135         call mpp_write(unit,var,domain,atmp(v_shape(1,1):v_shape(2,1) &
00136                                              ,v_shape(1,2):v_shape(2,2) &
00137                                              ,v_shape(1,3):v_shape(2,3) &
00138                                              ,v_shape(1,4):v_shape(2,4)) &
00139                                       ,tstamp=time)
00140       else
00141         call mpp_write( unit,var,domain,atmp(v_shape(1,1):v_shape(2,1) &
00142                                               ,v_shape(1,2):v_shape(2,2) &
00143                                              ,v_shape(1,3):v_shape(2,3) &
00144                                               ,v_shape(1,4):v_shape(2,4)))
00145       endif
00146       endif
00147 
00148       deallocate(atmp,stat=ierror)
00149 
00150       if ( ierror /= 0 ) then
00151          ierrp (1) = 1
00152          ierror = PRISM_Error_Alloc
00153          call psmile_error ( ierror, 'deallocate(atmp)', ierrp, 1, __FILE__ &
00154                            , __LINE__ )
00155          return
00156       endif
00157 #ifdef __PSMILE_IO_SYNC
00158       call mpp_flush(unit)
00159 #endif
00160 
00161 
00162 
00163 #ifdef VERBOSE
00164       print*,trim(ch_id),' :   psmile_write_4d_log: end'
00165       call psmile_flushstd
00166 
00167 #endif
00168 #endif
00169       end subroutine psmile_write_4d_log

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1