prismtrs_conserv2d_weight.F90

Go to the documentation of this file.
00001 !------------------------------------------------------------------------
00002 ! Copyright 2006-2010, CERFACS, Toulouse, France.
00003 ! All rights reserved. Use is subject to OASIS4 license terms.
00004 !-----------------------------------------------------------------------
00005 !BOP
00006 
00007 ! !ROUTINE: PRISMTrs_Conserv2D_weight
00008 !
00009 ! !INTERFACE
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   ! !USES:
00027   !
00028   USE PRISMDrv, dummy_interface => PRISMTrs_Conserv2D_weight
00029   !
00030   IMPLICIT NONE
00031   !
00032   ! !PARAMETERS:
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     ! EPIO Id
00048   INTEGER :: id_interp_id   ! Interpolation Id
00049   INTEGER :: id_idim        ! Dimension concerned with current interpolation
00050   INTEGER :: id_num_wts     ! Number of weights for 2D conservative remapping
00051   INTEGER :: id_err
00052   !
00053   ! !DESCRIPTION
00054   ! Subroutine "PRISMTrs_Conserv2D_weight" sets up the environment
00055   ! needed to call PRISMTrs_remap_conserv that computes the weights with a 
00056   ! 2d conservative remapping method. 
00057   ! This routine is called by prismtrs_interp.F90.
00058   !
00059   ! !REVISED HISTORY
00060   !   Date      Programmer   Description
00061   ! ----------  ----------   -----------
00062   ! 11/12/2006  S. Valcke     Creation
00063   !
00064   ! EOP
00065   !----------------------------------------------------------------------
00066   ! $Id: prismtrs_conserv2d_weight.F90 2082 2009-10-21 13:31:19Z hanke $
00067   ! $Author: hanke $
00068   !----------------------------------------------------------------------
00069   !
00070   ! Local declarations
00071   !
00072   INTEGER           :: stat, nc_scpid, il_loc, il_l
00073   CHARACTER (len=1) :: cl_charepio
00074   CHARACTER (len=8) :: crmpfile 
00075   CHARACTER (len=12):: cweight          ! string for weights
00076   CHARACTER (len=11):: csrcadd, cdstadd ! string for source/dest grid addresses
00077   CHARACTER (len=13):: cdstare, cdstfra ! string for destination grid area/frac
00078   LOGICAL           :: lcalc, first_call
00079 !
00080 #ifdef VERBOSE
00081   PRINT *, '| | | | | Enter PRISMTrs_Conserv2D_weight'
00082   PRINT *, '            '
00083   PRINT *, ' Calculation of SCRIP remapping matrix for CONSERV2D '
00084   call psmile_flushstd 
00085 #endif 
00086   !
00087   !-----------------------------------------------------------------------
00088   !
00089   !  In T, lats/lons are in radians (converted from degrees when received)
00090   !  Here, convert longitudes to 0,2pi interval
00091   !
00092   !-----------------------------------------------------------------------
00093   !
00094   WHERE (dda_src_lon .GT. dble_pi2)  dda_src_lon = dda_src_lon - dble_pi2
00095   WHERE (dda_src_lon .LT. zero) dda_src_lon = dda_src_lon + dble_pi2
00096   WHERE (dda_tgt_lon .GT. dble_pi2)  dda_tgt_lon = dda_tgt_lon - dble_pi2
00097   WHERE (dda_tgt_lon .LT. zero) dda_tgt_lon = dda_tgt_lon + dble_pi2
00098   !
00099   !-----------------------------------------------------------------------
00100   !
00101   !  Here, make sure input latitude range is within the machine values
00102   !     for +/- pi/2 
00103   !
00104   !-----------------------------------------------------------------------
00105   !
00106   WHERE (dda_src_lat >  dble_pih) dda_src_lat =  dble_pih
00107   WHERE (dda_src_lat < -dble_pih) dda_src_lat = -dble_pih
00108   WHERE (dda_tgt_lat >  dble_pih) dda_tgt_lat =  dble_pih
00109   WHERE (dda_tgt_lat < -dble_pih) dda_tgt_lat = -dble_pih
00110   !
00111   call prismtrs_remap_conserv (                            &
00112      id_epiosrc_size, id_epiotgt_size,                     &
00113      id_src_nbr_corner, id_tgt_nbr_corner,                 &
00114      id_src_lonlatz_size,                                  &
00115      id_index_size1,                                       &
00116      dda_src_lat, dda_src_lon,                             &
00117      dda_tgt_lat, dda_tgt_lon,                             &
00118      ida_nbsrccells_pertgtpt,                              &
00119      ida_index_array,                                      &
00120      ida_srcepio_add,                                      &
00121      id_epio_id,                                           &
00122      id_interp_id,                                         &
00123      id_idim,                                              &
00124      id_num_wts,                                           &
00125      id_err)
00126   !
00127   IF (Drv_Epios(id_epio_id)%num_links_map1 >= 1) &
00128      call prismtrs_sort_add(         &
00129      Drv_Epios(id_epio_id)%grid2_add_map1,      &
00130      Drv_Epios(id_epio_id)%grid1_add_map1,                 &
00131      Drv_Epios(id_epio_id)%wts_map1,                       &
00132      Drv_Epios(id_epio_id)%num_links_map1, id_num_wts)
00133   !
00134   ! No call to fracnnei (option not supported yet)
00135   !
00136   IF (Drv_Epios(id_epio_id)%num_links_map1 /=  &
00137      Drv_Epios(id_epio_id)%max_links_map1) THEN
00138       il_loc = Drv_Epios(id_epio_id)%num_links_map1 - &
00139          Drv_Epios(id_epio_id)%max_links_map1
00140       call prismtrs_resize_remap_vars(il_loc, id_epio_id, id_num_wts)
00141   ENDIF
00142 #ifdef DEBUGX
00143       DO il_l=1,Drv_Epios(id_epio_id)%num_links_map1
00144         il_loc = 7202
00145         IF (Drv_Epios(id_epio_id)%grid2_add_map1(il_l) == il_loc) THEN
00146             PRINT *, 'Source and weight end PRISMTrs_Conserv2D_weight for target point ', il_loc
00147             PRINT *, Drv_Epios(id_epio_id)%grid1_add_map1(il_l), Drv_Epios(id_epio_id)%wts_map1(1,il_l)
00148             call psmile_flushstd
00149         ENDIF
00150       ENDDO
00151 #endif
00152 #ifdef VERBOSE
00153   PRINT *, '| | | | | Quit PRISMTrs_Conserv2D_weight'
00154   PRINT *, '            '
00155   call psmile_flushstd 
00156 #endif 
00157   !
00158   RETURN
00159 END SUBROUTINE PRISMTrs_Conserv2D_weight

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1