00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 subroutine psmile_mg_get_cyclic ( grid_id, range, tol, ierror)
00012
00013
00014
00015 use PRISM_constants
00016
00017 use PSMILe, dummy_interface => PSMILe_MG_get_cyclic
00018
00019 implicit none
00020
00021
00022
00023 integer, Intent (In) :: grid_id
00024
00025
00026
00027 Integer, Intent (In) :: range (2, ndim_3d)
00028
00029
00030
00031 Double precision, Intent (In) :: tol
00032
00033
00034
00035
00036
00037 integer, Intent (Out) :: ierror
00038
00039
00040
00041
00042
00043
00044
00045 Integer :: i
00046
00047 Type(Grid), Pointer :: grid_info
00048
00049
00050
00051 Integer :: nbr_lats
00052 #if 0
00053 Type (Corner_Block), Pointer :: corner_pointer
00054
00055
00056
00057 integer, parameter :: nerrp = 1
00058 integer :: ierrp (nerrp)
00059 #endif
00060
00061
00062
00063
00064
00065
00066
00067
00068
00069
00070
00071
00072
00073
00074
00075
00076
00077
00078 Character(len=len_cvs_string), save :: mycvs =
00079 '$Id: psmile_mg_get_cyclic.F90 2325 2010-04-21 15:00:07Z valcke $'
00080
00081
00082
00083 #ifdef VERBOSE
00084 print 9990, trim(ch_id), grid_id
00085
00086 call psmile_flushstd
00087 #endif /* VERBOSE */
00088
00089
00090
00091
00092
00093
00094
00095
00096 #if 1
00097 ierror = 0
00098 grid_info => Grids(grid_id)
00099
00100 grid_info%cyclic = grid_info%periodic == PSMILE_true
00101
00102 if (grid_info%grid_type /= PRISM_Gaussreduced_regvrt) then
00103 if (Associated (grid_info%partition)) then
00104
00105
00106
00107 do i = 1, ndim_3d
00108 grid_info%cyclic (i) = grid_info%cyclic (i) .and. &
00109 grid_info%len_periodic(i) == &
00110 grid_info%grid_shape(2,i)-grid_info%grid_shape(1,i)+1
00111 end do
00112 else
00113
00114 endif
00115 else
00116
00117
00118
00119 nbr_lats = size (grid_info%extent(:,1))
00120
00121 if (grid_info%cyclic (1)) then
00122
00123
00124
00125
00126
00127
00128
00129
00130
00131
00132 do i = 1, nbr_lats
00133
00134
00135
00136
00137 if (grid_info%extent(i, 1) == &
00138 grid_info%nbr_points_per_lat(grid_info%l2g(i,1)) ) exit
00139 end do
00140
00141 grid_info%cyclic (1) = i <= nbr_lats
00142 endif
00143
00144 if (grid_info%cyclic (2)) then
00145
00146 grid_info%cyclic (2) = nbr_lats == grid_info%nbr_latitudes
00147 endif
00148
00149 if (grid_info%cyclic (3)) then
00150
00151 grid_info%cyclic (3) = grid_info%len_periodic(3) == &
00152 grid_info%grid_shape(2,3)-grid_info%grid_shape(1,3)+1
00153 endif
00154
00155 endif
00156 #else
00157
00158 Grids(grid_id)%cyclic = .false.
00159
00160
00161
00162 corner_pointer => Grids(grid_id)%corner_pointer
00163
00164 if (corner_pointer%corner_datatype == MPI_REAL) then
00165
00166
00167
00168 call psmile_mg_get_cyclic_real ( grid_id, range, &
00169 real(tol), ierror)
00170
00171 else if (corner_pointer%corner_datatype == MPI_DOUBLE_PRECISION) then
00172
00173
00174
00175 call psmile_mg_get_cyclic_dble ( grid_id, range, &
00176 tol, ierror)
00177
00178 #if defined ( PRISM_QUAD_TYPE )
00179 else if (corner_pointer%corner_datatype == MPI_REAL16) then
00180
00181
00182
00183 call psmile_mg_get_cyclic_quad ( grid_id, range, &
00184 real(tol, 16), ierror)
00185 #endif
00186
00187 else
00188
00189
00190
00191 ierrp (1) = corner_pointer%corner_datatype
00192 ierror = PRISM_Error_Internal
00193 call psmile_error ( ierror, 'unsupported data type', &
00194 ierrp, 1, __FILE__, __LINE__ )
00195 endif
00196 #endif
00197
00198
00199
00200 #ifdef VERBOSE
00201 print 9980, trim(ch_id), ierror, Grids(grid_id)%cyclic
00202
00203 call psmile_flushstd
00204 #endif /* VERBOSE */
00205
00206
00207
00208
00209
00210 #ifdef VERBOSE
00211
00212 9990 format (1x, a, ': psmile_mg_get_cyclic: grid_id', i3)
00213 9980 format (1x, a, ': psmile_mg_get_cyclic: eof ierror =', i3, '; cyclic:', 3l3)
00214
00215 #endif /* VERBOSE */
00216
00217
00218 end subroutine PSMILe_MG_get_cyclic