00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 subroutine psmile_read_4d_int(unit &
00012 ,var &
00013 ,domain &
00014 ,a &
00015 ,a_shape &
00016 ,v_shape &
00017 ,itime &
00018 ,time_used &
00019 ,id_blockid &
00020 ,block_used &
00021 ,domain_used &
00022 ,ierror)
00023
00024
00025
00026 use psmile
00027 #ifdef __PSMILE_WITH_IO
00028 implicit none
00029 include 'prism.inc'
00030
00031
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 integer , 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 integer, intent(in) :: itime
00043 logical,intent(in) :: time_used
00044 integer,intent(in) :: id_blockid
00045 logical,intent(in) :: block_used
00046 logical,intent(in) :: domain_used
00047
00048
00049
00050 integer, intent(out) :: ierror
00051
00052
00053
00054 double precision,allocatable :: atmp(:,:,:,:)
00055 integer :: ierrp(2)
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066
00067
00068
00069
00070
00071
00072 character(len=len_cvs_string),save :: mcvs =
00073 '$Id: psmile_read_4d_int.F90 2325 2010-04-21 15:00:07Z valcke $'
00074
00075 ierror=0
00076
00077 #ifdef VERBOSE
00078 print*,trim(ch_id),' : psmile_read_4d_dble: start'
00079 print*,trim(ch_id),' : psmile_read_4d_dble: 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)))
00088 if ( ierror /= 0 ) then
00089 ierrp (1) = 1
00090 ierror = PRISM_Error_Alloc
00091 call psmile_error ( ierror, 'atmp', ierrp, 1, __FILE__, __LINE__ )
00092 return
00093 endif
00094
00095 if(domain_used) then
00096 if(block_used) then
00097 if(id_blockid.le.0) then
00098 ierror = PRISM_Error_Internal
00099 call psmile_error ( ierror, 'id_blockid <= 0! ', &
00100 ierrp, 0, __FILE__, __LINE__ )
00101 endif
00102 if(time_used) then
00103 call mpp_read(unit,var,domain,atmp(v_shape(1,1):v_shape(2,1) &
00104 ,v_shape(1,2):v_shape(2,2) &
00105 ,v_shape(1,3):v_shape(2,3) &
00106 ,v_shape(1,4):v_shape(2,4)) &
00107 ,tindex=itime,blockid=id_blockid)
00108 else
00109 call mpp_read( unit,var,domain,atmp(v_shape(1,1):v_shape(2,1) &
00110 ,v_shape(1,2):v_shape(2,2) &
00111 ,v_shape(1,3):v_shape(2,3) &
00112 ,v_shape(1,4):v_shape(2,4)) &
00113 , blockid=id_blockid)
00114 endif
00115 else
00116 if(time_used) then
00117 call mpp_read(unit,var,domain,atmp(v_shape(1,1):v_shape(2,1) &
00118 ,v_shape(1,2):v_shape(2,2) &
00119 ,v_shape(1,3):v_shape(2,3) &
00120 ,v_shape(1,4):v_shape(2,4)) &
00121 ,tindex=itime)
00122 else
00123 call mpp_read( unit,var,domain,atmp(v_shape(1,1):v_shape(2,1) &
00124 ,v_shape(1,2):v_shape(2,2) &
00125 ,v_shape(1,3):v_shape(2,3) &
00126 ,v_shape(1,4):v_shape(2,4)))
00127 endif
00128 endif
00129
00130 else
00131
00132 if(block_used) then
00133 if(id_blockid.le.0) then
00134 ierror = PRISM_Error_Internal
00135 call psmile_error ( ierror, 'id_blockid <= 0! ', &
00136 ierrp, 0, __FILE__, __LINE__ )
00137 endif
00138 if(time_used) then
00139 call mpp_read(unit,var,atmp(v_shape(1,1):v_shape(2,1) &
00140 ,v_shape(1,2):v_shape(2,2) &
00141 ,v_shape(1,3):v_shape(2,3) &
00142 ,v_shape(1,4):v_shape(2,4)) &
00143 ,tindex=itime,blockid=id_blockid)
00144 else
00145 call mpp_read( unit,var,atmp(v_shape(1,1):v_shape(2,1) &
00146 ,v_shape(1,2):v_shape(2,2) &
00147 ,v_shape(1,3):v_shape(2,3) &
00148 ,v_shape(1,4):v_shape(2,4)) &
00149 , blockid=id_blockid)
00150 endif
00151 else
00152 if(time_used) then
00153 call mpp_read(unit,var,atmp(v_shape(1,1):v_shape(2,1) &
00154 ,v_shape(1,2):v_shape(2,2) &
00155 ,v_shape(1,3):v_shape(2,3) &
00156 ,v_shape(1,4):v_shape(2,4)) &
00157 ,tindex=itime)
00158 else
00159 call mpp_read( unit,var,atmp(v_shape(1,1):v_shape(2,1) &
00160 ,v_shape(1,2):v_shape(2,2) &
00161 ,v_shape(1,3):v_shape(2,3) &
00162 ,v_shape(1,4):v_shape(2,4)))
00163 endif
00164 endif
00165 endif
00166
00167
00168 a(v_shape(1,1):v_shape(2,1) &
00169 ,v_shape(1,2):v_shape(2,2) &
00170 ,v_shape(1,3):v_shape(2,3) &
00171 ,v_shape(1,4):v_shape(2,4)) = &
00172 atmp(v_shape(1,1):v_shape(2,1) &
00173 ,v_shape(1,2):v_shape(2,2) &
00174 ,v_shape(1,3):v_shape(2,3) &
00175 ,v_shape(1,4):v_shape(2,4))
00176
00177
00178 deallocate(atmp,stat=ierror)
00179
00180 if ( ierror /= 0 ) then
00181 ierrp (1) = 1
00182 ierror = PRISM_Error_Alloc
00183 call psmile_error ( ierror, 'deallocate(atmp)', ierrp, 1, __FILE__ &
00184 , __LINE__ )
00185 return
00186 endif
00187
00188 #ifdef VERBOSE
00189 print*,trim(ch_id),' : psmile_read_4d_int: end'
00190 call psmile_flushstd
00191
00192 #endif
00193 #endif
00194 end subroutine psmile_read_4d_int