00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 subroutine psmile_unfolded_check_real( tgt_cell_index, & !< (in)
00012 unfolded_check, &
00013 grid_type, &
00014 visited, &
00015 nbr_tgt_corners, &
00016 tgt_corners, &
00017 grid_shape, &
00018 corner_shape_3d, &
00019 nbr_src_corners, &
00020 src_corner_x, &
00021 src_corner_y, &
00022 stacksize, &
00023 stack, &
00024 stack_index, &
00025 ierror )
00026
00027
00028
00029 use PRISM_constants
00030
00031 use PSMILe
00032
00033 implicit none
00034
00035
00036
00037 Integer, Intent (In) :: tgt_cell_index
00038
00039
00040
00041 Integer, Intent (In) :: unfolded_check
00042
00043
00044
00045
00046
00047 Integer, Intent (In) :: grid_type
00048
00049
00050
00051 Integer, Intent (In) :: grid_shape(2,ndim_3d)
00052
00053
00054
00055 Integer, Intent (InOut) :: visited(grid_shape(1,1):grid_shape(2,1),
00056 grid_shape(1,2):grid_shape(2,2))
00057
00058
00059
00060
00061 Integer, Intent (In) :: nbr_tgt_corners
00062
00063
00064
00065 Type(line_real), Intent (In) :: tgt_corners(nbr_tgt_corners)
00066
00067
00068
00069 Integer, Intent (In) :: nbr_src_corners
00070
00071
00072
00073 Integer, Intent (In) :: corner_shape_3d(2,ndim_3d,ndim_3d)
00074
00075
00076
00077 Real, Intent (In) :: src_corner_x (corner_shape_3d (1,1,1):corner_shape_3d(2,1,1),
00078 corner_shape_3d (1,2,1):corner_shape_3d(2,2,1),
00079 corner_shape_3d (1,3,1):corner_shape_3d(2,3,1),
00080 nbr_src_corners )
00081
00082
00083
00084 Real, Intent (In) :: src_corner_y (corner_shape_3d (1,1,2):corner_shape_3d(2,1,2),
00085 corner_shape_3d (1,2,2):corner_shape_3d(2,2,2),
00086 corner_shape_3d (1,3,2):corner_shape_3d(2,3,2),
00087 nbr_src_corners )
00088
00089
00090
00091 Integer, Intent (In) :: stacksize
00092
00093
00094
00095 Type(memo), Intent (InOut) :: stack(stacksize)
00096
00097
00098
00099 Integer, Intent (InOut) :: stack_index
00100
00101
00102
00103
00104
00105
00106 Integer, Intent (Out) :: ierror
00107
00108
00109
00110
00111
00112
00113
00114
00115
00116
00117
00118
00119 Type (line_real), Allocatable :: src(:)
00120
00121 Logical :: overlap
00122 Integer :: j1, j2
00123 Integer :: ii, jj
00124
00125
00126
00127 Integer, Parameter :: nerrp = 2
00128 Integer :: ierrp (nerrp)
00129
00130
00131
00132
00133
00134
00135
00136
00137
00138
00139
00140
00141
00142
00143
00144
00145
00146
00147
00148 Character(len=len_cvs_string), save :: mycvs =
00149 '$Id$'
00150
00151
00152
00153
00154
00155 #ifdef VERBOSE
00156 print 9990, trim(ch_id)
00157
00158 call psmile_flushstd
00159 #endif /* VERBOSE */
00160
00161 ierror = 0
00162
00163 jj = grid_shape(unfolded_check,2)
00164
00165 allocate(src(nbr_src_corners), STAT=ierror)
00166
00167 if ( ierror > 0 ) then
00168 ierrp (1) = ierror
00169 ierrp (2) = 8
00170 ierror = PRISM_Error_Alloc
00171 call psmile_error ( ierror, 'src struct', &
00172 ierrp, 2, __FILE__, __LINE__ )
00173 return
00174 endif
00175
00176 do ii = grid_shape(1,1), grid_shape(2,1)
00177
00178 if ( visited(ii,jj) == tgt_cell_index ) cycle
00179
00180 Select Case ( grid_type )
00181
00182 Case ( PRISM_Irrlonlat_regvrt )
00183
00184 do j1 = 1, nbr_src_corners
00185 j2 = mod(j1,nbr_tgt_corners)+1
00186 src(j1)%p1%x = src_corner_x(ii,jj,1,j1)
00187 src(j1)%p1%y = src_corner_y(ii,jj,1,j1)
00188 src(j1)%p2%x = src_corner_x(ii,jj,1,j2)
00189 src(j1)%p2%y = src_corner_y(ii,jj,1,j2)
00190 enddo
00191
00192 Case ( PRISM_Reglonlatvrt )
00193
00194 do j1 = 1, nbr_src_corners
00195 j2 = mod(j1,nbr_tgt_corners)+1
00196 src(j1)%p1%x = src_corner_x(ii, 1,1,j1)
00197 src(j1)%p1%y = src_corner_y( 1,jj,1,j1)
00198 src(j1)%p2%x = src_corner_x(ii, 1,1,j2)
00199 src(j1)%p2%y = src_corner_y( 1,jj,1,j2)
00200 enddo
00201
00202 Case Default
00203
00204 ierrp (1) = grid_type
00205 ierror = PRISM_Error_Internal
00206
00207 call psmile_error ( ierror, 'unsupported grid generation type', &
00208 ierrp, 1, __FILE__, __LINE__ )
00209
00210 End Select
00211
00212 call psmile_overlap_real ( nbr_src_corners, nbr_tgt_corners, src, tgt_corners, overlap )
00213
00214 visited(ii,jj) = tgt_cell_index
00215
00216 stack_index = 1
00217 stack(stack_index)%dir = 0
00218
00219 stack(stack_index)%overlap = overlap
00220
00221 if ( .not. overlap ) then
00222
00223 cycle
00224
00225 else
00226
00227 stack(stack_index)%i = ii
00228 stack(stack_index)%j = jj
00229
00230 exit
00231
00232 endif
00233
00234 enddo
00235
00236 if ( ii > grid_shape(2,1) ) then
00237
00238 stack_index = 0
00239 endif
00240
00241 deallocate(src, STAT=ierror)
00242
00243 if ( ierror > 0 ) then
00244 ierrp (1) = ierror
00245 ierrp (2) = 8
00246 ierror = PRISM_Error_Dealloc
00247 call psmile_error ( ierror, 'src struct', &
00248 ierrp, 2, __FILE__, __LINE__ )
00249 return
00250 endif
00251
00252 #ifdef VERBOSE
00253 print 9980, trim(ch_id)
00254
00255 call psmile_flushstd
00256 #endif /* VERBOSE */
00257
00258
00259
00260 9990 format (1x, a, ': psmile_unfolded_check_real')
00261 9980 format (1x, a, ': psmile_unfolded_check_real: eof')
00262
00263 end subroutine psmile_unfolded_check_real