00001
00002
00003
00004
00005
00006
00007
00008
00009
00010 subroutine prismtrs_conserv2d_weight( &
00011 id_epiosrc_size, id_epiotgt_size, &
00012 id_src_nbr_corner, id_tgt_nbr_corner, &
00013 id_src_lonlatz_size, &
00014 id_index_size1, &
00015 dda_src_lat, dda_src_lon, dda_src_z, &
00016 dda_tgt_lat, dda_tgt_lon, dda_tgt_z, &
00017 ida_nbsrccells_pertgtpt, &
00018 ida_index_array, &
00019 ida_srcepio_add, &
00020 id_epio_id, &
00021 id_interp_id, &
00022 id_idim, &
00023 id_num_wts, &
00024 id_err)
00025
00026
00027
00028 USE PRISMDrv, dummy_interface => PRISMTrs_Conserv2D_weight
00029
00030 IMPLICIT NONE
00031
00032
00033
00034 INTEGER :: id_epiosrc_size, id_epiotgt_size
00035 INTEGER :: id_src_nbr_corner, id_tgt_nbr_corner
00036 INTEGER :: id_src_lonlatz_size
00037 INTEGER :: id_index_size1
00038 DOUBLE PRECISION, DIMENSION(id_src_lonlatz_size) :: dda_src_lat
00039 DOUBLE PRECISION, DIMENSION(id_src_lonlatz_size) :: dda_src_lon
00040 DOUBLE PRECISION, DIMENSION(id_src_lonlatz_size) :: dda_src_z
00041 DOUBLE PRECISION, DIMENSION(id_epiotgt_size,id_tgt_nbr_corner) :: dda_tgt_lat
00042 DOUBLE PRECISION, DIMENSION(id_epiotgt_size,id_tgt_nbr_corner) :: dda_tgt_lon
00043 DOUBLE PRECISION, DIMENSION(id_epiotgt_size,id_tgt_nbr_corner) :: dda_tgt_z
00044 INTEGER, DIMENSION(id_epiotgt_size) :: ida_nbsrccells_pertgtpt
00045 INTEGER, DIMENSION(id_index_size1,id_src_nbr_corner) :: ida_index_array
00046 INTEGER, DIMENSION(id_index_size1) :: ida_srcepio_add
00047 INTEGER :: id_epio_id
00048 INTEGER :: id_interp_id
00049 INTEGER :: id_idim
00050 INTEGER :: id_num_wts
00051 INTEGER :: id_err
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066
00067
00068
00069
00070
00071
00072 INTEGER :: stat, nc_scpid, il_loc, il_l, ib
00073 CHARACTER (len=1) :: cl_charepio
00074 CHARACTER (len=8) :: crmpfile
00075 CHARACTER (len=12):: cweight
00076 CHARACTER (len=11):: csrcadd, cdstadd
00077 CHARACTER (len=13):: cdstare, cdstfra
00078 LOGICAL :: lcalc, first_call
00079 DOUBLE PRECISION :: common_grid_range(2)
00080 DOUBLE PRECISION :: shiftSize
00081
00082 #ifdef VERBOSE
00083 PRINT *, '| | | | | Enter PRISMTrs_Conserv2D_weight'
00084 PRINT *, ' '
00085 PRINT *, ' Calculation of SCRIP remapping matrix for CONSERV2D '
00086 call psmile_flushstd
00087 #endif
00088
00089
00090
00091
00092
00093
00094
00095 common_grid_range(1) = 0.
00096 common_grid_range(2) = dble_pi2
00097 shiftSize = common_grid_range(2) - common_grid_range(1)
00098
00099 DO ib = 1, id_src_lonlatz_size
00100 DO WHILE (dda_src_lon(ib) < common_grid_range(1))
00101 dda_src_lon(ib) = dda_src_lon(ib) + shiftSize
00102 ENDDO
00103 DO WHILE (dda_src_lon(ib) >= common_grid_range(2))
00104 dda_src_lon(ib) = dda_src_lon(ib) - shiftSize
00105 ENDDO
00106 ENDDO
00107 DO ib = 1, id_epiotgt_size
00108 DO il_l=1,id_tgt_nbr_corner
00109 DO WHILE (dda_tgt_lon(ib,il_l) < common_grid_range(1))
00110 dda_tgt_lon(ib,il_l) = dda_tgt_lon(ib,il_l) + shiftSize
00111 ENDDO
00112 DO WHILE (dda_tgt_lon(ib,il_l) >= common_grid_range(2))
00113 dda_tgt_lon(ib,il_l) = dda_tgt_lon(ib,il_l) - shiftSize
00114 ENDDO
00115 ENDDO
00116 ENDDO
00117
00118
00119
00120
00121
00122
00123
00124
00125 DO ib = 1, id_src_lonlatz_size
00126 IF ( (dda_src_lat(ib) > dble_pih) .OR. (dda_src_lat(ib) < -dble_pih)) THEN
00127 WRITE(*,*) 'Problem :Latitudes of source are greater than +/-90 degrees'
00128 CALL psmile_abort
00129 ENDIF
00130 ENDDO
00131 DO ib = 1, id_epiotgt_size
00132 DO il_l=1,id_tgt_nbr_corner
00133 IF ( (dda_tgt_lat(ib,il_l) > dble_pih) .OR. (dda_tgt_lat(ib,il_l) < -dble_pih) ) THEN
00134 WRITE(*,*) 'Problem :Latitudes of target are greater than +/-90 degrees'
00135 CALL psmile_abort
00136 ENDIF
00137 ENDDO
00138 ENDDO
00139
00140 call prismtrs_remap_conserv ( &
00141 id_epiosrc_size, id_epiotgt_size, &
00142 id_src_nbr_corner, id_tgt_nbr_corner, &
00143 id_src_lonlatz_size, &
00144 id_index_size1, &
00145 dda_src_lat, dda_src_lon, &
00146 dda_tgt_lat, dda_tgt_lon, &
00147 ida_nbsrccells_pertgtpt, &
00148 ida_index_array, &
00149 ida_srcepio_add, &
00150 id_epio_id, &
00151 id_interp_id, &
00152 id_idim, &
00153 id_num_wts, &
00154 id_err)
00155
00156 IF (Drv_Epios(id_epio_id)%num_links_map1 >= 1) &
00157 call prismtrs_sort_add( &
00158 Drv_Epios(id_epio_id)%grid2_add_map1, &
00159 Drv_Epios(id_epio_id)%grid1_add_map1, &
00160 Drv_Epios(id_epio_id)%wts_map1, &
00161 Drv_Epios(id_epio_id)%num_links_map1, id_num_wts)
00162
00163
00164
00165 IF (Drv_Epios(id_epio_id)%num_links_map1 /= &
00166 Drv_Epios(id_epio_id)%max_links_map1) THEN
00167 il_loc = Drv_Epios(id_epio_id)%num_links_map1 - &
00168 Drv_Epios(id_epio_id)%max_links_map1
00169 call prismtrs_resize_remap_vars(il_loc, id_epio_id, id_num_wts)
00170 ENDIF
00171 #ifdef DEBUGX
00172 DO il_l=1,Drv_Epios(id_epio_id)%num_links_map1
00173 il_loc = 4848
00174 IF (Drv_Epios(id_epio_id)%grid2_add_map1(il_l) == il_loc) THEN
00175 WRITE(88,*) ''
00176 WRITE(88,*), 'Source and weight end PRISMTrs_Conserv2D_weight for target point ', il_loc
00177 WRITE(88,*), Drv_Epios(id_epio_id)%grid1_add_map1(il_l), Drv_Epios(id_epio_id)%wts_map1(1,il_l)
00178 call psmile_flushstd
00179 ENDIF
00180 ENDDO
00181 #endif
00182 #ifdef VERBOSE
00183 PRINT *, '| | | | | Quit PRISMTrs_Conserv2D_weight'
00184 PRINT *, ' '
00185 call psmile_flushstd
00186 #endif
00187
00188 RETURN
00189 END SUBROUTINE PRISMTrs_Conserv2D_weight