prismtrs_interp.F90

Go to the documentation of this file.
00001 !------------------------------------------------------------------------
00002 ! Copyright 2006-2010, CERFACS, Toulouse, France.
00003 ! Copyright 2006-2010, NEC High Performance Computing, Duesseldorf, Germany.
00004 ! All rights reserved. Use is subject to OASIS4 license terms.
00005 !-----------------------------------------------------------------------
00006 !BOP
00007 !
00008 ! !ROUTINE: PRISMTrs_Interp
00009 !
00010 ! !INTERFACE
00011 subroutine prismtrs_interp(id_exchange_id,       &
00012                            id_epio_id,           &
00013                            id_trans_out_size,    &
00014                            dda_trans_out,        &
00015                            id_trans_in_size,     &
00016                            dda_trans_in,         &
00017                            nbr_fields,           &
00018                            id_err)
00019 
00020 !
00021 ! !USES:
00022 !
00023   USE PRISMDrv, dummy_interface => PRISMTrs_Interp
00024 
00025   IMPLICIT NONE
00026 
00027 !
00028 ! !PARAMETERS:
00029 !
00030 ! exchange id
00031   INTEGER, INTENT (In)        :: id_exchange_id
00032 
00033 ! epio id
00034   INTEGER, INTENT (In)        :: id_epio_id
00035 
00036 ! size of the trans_out field
00037   INTEGER, INTENT (In)        :: id_trans_out_size
00038 
00039 ! trans_out field
00040   DOUBLE PRECISION, DIMENSION(id_trans_out_size), INTENT (IN) :: dda_trans_out
00041 
00042 ! size of the trans_in field
00043   INTEGER, INTENT (In)        :: id_trans_in_size
00044 
00045 ! nbr of bundle components
00046   INTEGER, INTENT (In)        :: nbr_fields
00047 
00048 !
00049 ! ! RETURN VALUE
00050 !
00051 ! trans_in field
00052   DOUBLE PRECISION, DIMENSION(id_trans_in_size), INTENT (Out) :: dda_trans_in
00053 
00054   INTEGER, INTENT (Out)       :: id_err   ! error value
00055 
00056 ! !DESCRIPTION
00057 ! Subroutine "PRISMTrs_Interp" manages the interpolation operations.
00058 !
00059 ! !REVISED HISTORY
00060 !   Date      Programmer   Description
00061 ! ----------  ----------   -----------
00062 ! 06/12/2002  D. Declat    Creation
00063 ! 06/03/2008  J. Charles   Modifications added for the use of bundle fields
00064 !
00065 ! EOP
00066 !----------------------------------------------------------------------
00067 ! $Id: prismtrs_interp.F90 2685 2010-10-28 14:05:10Z coquart $
00068 ! $Author: coquart $
00069 !----------------------------------------------------------------------
00070 !
00071 ! Local declarations
00072 !
00073   CHARACTER(LEN=len_cvs_string), SAVE  :: mycvs = 
00074      '$Id: prismtrs_interp.F90 2685 2010-10-28 14:05:10Z coquart $'
00075 
00076 !
00077 ! EPIOs
00078 !
00079 ! src and tgt sizes
00080   INTEGER                          :: il_src_size,  il_tgt_size
00081 !
00082 ! Interpolation arguments
00083 !
00084 ! method & type
00085   INTEGER, DIMENSION(3)                 :: ila_interp_method  
00086   INTEGER                               :: il_interp_type
00087 
00088   INTEGER                               :: il_max_links_map1 
00089   INTEGER                               :: il_index_size1
00090 
00091 ! number of neighbors on the source grid for each target point
00092   INTEGER                               :: il_nb_neighbors
00093 
00094 !
00095   INTEGER                               :: il_interp_id
00096   INTEGER                               :: il_method
00097 !
00098 ! Number of weights for 2D conservative remapping
00099   INTEGER                               :: il_num_wts = 3
00100   DOUBLE PRECISION                      :: dl_gaus_var
00101 
00102 ! Storage variables
00103   DOUBLE PRECISION, DIMENSION(:), ALLOCATABLE :: dla_temp_field
00104   REAl, DIMENSION(:), ALLOCATABLE             :: rla_temp_field
00105   INTEGER, DIMENSION(:), ALLOCATABLE          :: ila_temp_field
00106 
00107 !
00108 !     ... for error handling
00109 !
00110   INTEGER, PARAMETER              :: nerrp = 2
00111   INTEGER                         :: ierrp (nerrp)
00112 !
00113 ! ---------------------------------------------------------------------
00114 !
00115 #ifdef VERBOSE
00116   PRINT *, '| | | | | Enter PRISMTrs_Interp for epio ID = ', id_epio_id
00117   call psmile_flushstd
00118 #endif
00119 
00120 !
00121 #ifdef DEBUG
00122   IF ( Drv_Epios(id_epio_id)%tgt_size /= &
00123      Drv_Exchanges(id_exchange_id)%trans_in_field_size) THEN
00124       PRINT *,'|| tgt_size not equal to trans_in_field_size in prismtrs_interp'
00125       call psmile_abort
00126   ENDIF
00127   IF ( id_trans_in_size /= Drv_Epios(id_epio_id)%tgt_size*nbr_fields) THEN
00128       PRINT *,'|| id_trans_in_size not equal to tgt_size in prismtrs_interp'
00129       call psmile_abort
00130   ENDIF
00131 #endif
00132 ! 1. Get the information used in the interpolations operations
00133 !
00134 ! 1.1. Get the source and target sizes
00135   il_src_size = Drv_Epios(id_epio_id)%src_size
00136   il_tgt_size = Drv_Epios(id_epio_id)%tgt_size
00137 
00138 ! 1.5. Set the interpolation arguments
00139 !
00140 ! 1.5.1. Get the interpolation method
00141   il_interp_id = Drv_Exchanges(id_exchange_id)%interp_id
00142   ila_interp_method = Drv_Interps(il_interp_id)%interp_method
00143   il_interp_type = Drv_Interps(il_interp_id)%interp_type
00144 
00145   il_nb_neighbors = Drv_Interps(il_interp_id)%arg2(1)
00146   il_method       = Drv_Interps(il_interp_id)%arg5(1)
00147 
00148   IF (ila_interp_method(1) .eq. PSMILe_bilinear ) THEN
00149       IF (ila_interp_method(2) .eq. PSMILe_none) THEN
00150           il_nb_neighbors = 4
00151       ELSE
00152           il_nb_neighbors = 8
00153       END IF
00154   ELSE IF (ila_interp_method(1) .eq. PSMILe_trilinear) THEN
00155       il_nb_neighbors = 8
00156   ELSE IF ( ila_interp_method(1) .eq. PSMILe_bicubic ) THEN
00157       IF (ila_interp_method(2) .eq. PSMILe_none) THEN
00158           il_nb_neighbors = 16
00159       ELSE IF (ila_interp_method(2) .eq. PSMILe_linear) THEN
00160           il_nb_neighbors = 32
00161       ELSE
00162           il_nb_neighbors = 64
00163       END IF
00164   ELSE IF (ila_interp_method(1) .EQ. PSMILe_conserv2D) THEN
00165       il_nb_neighbors = 0
00166   ELSE
00167       if ( il_nb_neighbors < 1 ) then
00168          ierrp (1) = il_nb_neighbors
00169 
00170          id_err = PRISM_Error_Arg
00171          call psmile_error_common ( id_err, 'nb_neighbors', &
00172                              ierrp, 1, __FILE__, __LINE__ )
00173          return
00174       endif
00175   END IF
00176 !
00177 ! 2. Compute the weights
00178 !
00179   IF (Drv_Epios(id_epio_id)%weights_status /= PSMILe_trans_complete) THEN
00180 !
00181 ! 2.3. Allocate the pointer in the fields structure to store the weights
00182       IF (ila_interp_method(1) == PSMILe_conserv2D) THEN
00183           Drv_Epios(id_epio_id)%resize_increment = &
00184              0.1*MAX(il_tgt_size, il_src_size)
00185           Drv_Epios(id_epio_id)%max_links_map1 = 4* il_tgt_size
00186           il_max_links_map1 = Drv_Epios(id_epio_id)%max_links_map1
00187           ALLOCATE (Drv_Epios(id_epio_id)%grid1_add_map1(il_max_links_map1), &
00188              Drv_Epios(id_epio_id)%grid2_add_map1(il_max_links_map1),        &
00189              Drv_Epios(id_epio_id)%wts_map1(il_num_wts, il_max_links_map1),     &
00190              stat = id_err)
00191           IF ( id_err > 0 ) THEN
00192               ierrp (1) = id_err
00193               ierrp (2) = il_max_links_map1
00194               call psmile_error_common ( PRISM_Error_Alloc, 'wts_map1', &
00195                  ierrp, 2, __FILE__, __LINE__ )
00196           ENDIF
00197           Drv_Epios(id_epio_id)%grid1_add_map1(:) = zero
00198           Drv_Epios(id_epio_id)%grid2_add_map1(:) = zero
00199           Drv_Epios(id_epio_id)%wts_map1(:,:) = zero
00200           Drv_Epios(id_epio_id)%num_links_map1 = zero
00201       ELSE        
00202           ALLOCATE(Drv_Epios(id_epio_id)%weights(il_tgt_size,il_nb_neighbors),&
00203              stat = id_err)
00204           IF ( id_err > 0 ) THEN
00205               ierrp (1) = id_err
00206               ierrp (2) = il_tgt_size*il_nb_neighbors
00207               call psmile_error_common ( PRISM_Error_Alloc, 'weights', &
00208                  ierrp, 2, __FILE__, __LINE__ )
00209           ENDIF
00210       ENDIF
00211 
00212 ! 2.4. If first time, compute the weights, depending on the type and the method
00213 
00214       IF (il_interp_type .EQ. PSMILe_3D) THEN
00215 
00216           SELECT CASE (ila_interp_method(1))
00217 
00218           CASE (PSMILe_nnghbr3D)
00219 
00220               dl_gaus_var = Drv_Interps(il_interp_id)%arg8
00221               IF ( dl_gaus_var .EQ. PSMILe_dundef) THEN
00222 #ifdef VERBOSE
00223                   PRINT *, &
00224                      '| | | | | | Trs is asked to interp distwght 3d'
00225                   call psmile_flushstd
00226 #endif
00227                   IF (Drv_Epios(id_epio_id)%src_coord_type .EQ. PRISM_Real) THEN
00228                       call prismtrs_distwght_weight_3d(                     &
00229                          il_src_size,                                       &
00230                          DBLE(Drv_Epios(id_epio_id)%src_lat_pointer_real),  &
00231                          DBLE(Drv_Epios(id_epio_id)%src_lon_pointer_real),  &
00232                          DBLE(Drv_Epios(id_epio_id)%src_z_pointer_real),    &
00233                          Drv_Epios(id_epio_id)%src_mask_pointer,  & 
00234                          il_tgt_size,                                       &
00235                          DBLE(Drv_Epios(id_epio_id)%tgt_lat_pointer_real),  &
00236                          DBLE(Drv_Epios(id_epio_id)%tgt_lon_pointer_real),  &
00237                          DBLE(Drv_Epios(id_epio_id)%tgt_z_pointer_real),    &
00238                          Drv_Epios(id_epio_id)%tgt_mask_pointer,            &
00239                          il_nb_neighbors,                                   &
00240                          Drv_Epios(id_epio_id)%index_array,                 &
00241                          Drv_Epios(id_epio_id)%weights,  &
00242                          id_err)
00243                   ELSE IF (Drv_Epios(id_epio_id)%src_coord_type .EQ.   &
00244                      PRISM_Double_Precision) THEN
00245                       call prismtrs_distwght_weight_3d(                &
00246                          il_src_size,                                  &
00247                          Drv_Epios(id_epio_id)%src_lat_pointer_dble,   &
00248                          Drv_Epios(id_epio_id)%src_lon_pointer_dble,   &
00249                          Drv_Epios(id_epio_id)%src_z_pointer_dble,     &
00250                          Drv_Epios(id_epio_id)%src_mask_pointer,       & 
00251                          il_tgt_size,                                  &
00252                          Drv_Epios(id_epio_id)%tgt_lat_pointer_dble,   &
00253                          Drv_Epios(id_epio_id)%tgt_lon_pointer_dble,   &
00254                          Drv_Epios(id_epio_id)%tgt_z_pointer_dble,     &
00255                          Drv_Epios(id_epio_id)%tgt_mask_pointer,       &
00256                          il_nb_neighbors,                              &
00257                          Drv_Epios(id_epio_id)%index_array,            &
00258                          Drv_Epios(id_epio_id)%weights,                &
00259                          id_err)
00260                   ENDIF
00261               ELSE
00262 
00263 #ifdef VERBOSE
00264                   PRINT *, &
00265                      '| | | | | | Trs is asked to interp gauswght 3d'
00266                   call psmile_flushstd
00267 #endif
00268                   IF (Drv_Epios(id_epio_id)%src_coord_type .EQ. PRISM_Real) THEN
00269                       call prismtrs_gauswght_weight_3d(                     &
00270                          il_src_size,                                       &
00271                          DBLE(Drv_Epios(id_epio_id)%src_lat_pointer_real),  &
00272                          DBLE(Drv_Epios(id_epio_id)%src_lon_pointer_real),  &
00273                          DBLE(Drv_Epios(id_epio_id)%src_z_pointer_real),    &
00274                          Drv_Epios(id_epio_id)%src_mask_pointer,            &
00275                          il_tgt_size,                                       &
00276                          DBLE(Drv_Epios(id_epio_id)%tgt_lat_pointer_real),  &
00277                          DBLE(Drv_Epios(id_epio_id)%tgt_lon_pointer_real),  &
00278                          DBLE(Drv_Epios(id_epio_id)%tgt_z_pointer_real),    &
00279                          Drv_Epios(id_epio_id)%tgt_mask_pointer,            &
00280                          dl_gaus_var,                                       &
00281                          il_nb_neighbors,                                   &
00282                          Drv_Epios(id_epio_id)%index_array,                 &
00283                          Drv_Epios(id_epio_id)%weights,                     &
00284                          id_err)
00285                   ELSE IF (Drv_Epios(id_epio_id)%src_coord_type .EQ.        &
00286                      PRISM_Double_Precision) THEN
00287                       call prismtrs_gauswght_weight_3d(                     &
00288                          il_src_size,                                       &
00289                          Drv_Epios(id_epio_id)%src_lat_pointer_dble,        &
00290                          Drv_Epios(id_epio_id)%src_lon_pointer_dble,        &
00291                          Drv_Epios(id_epio_id)%src_z_pointer_dble ,         &
00292                          Drv_Epios(id_epio_id)%src_mask_pointer,            &
00293                          il_tgt_size,                                       &
00294                          Drv_Epios(id_epio_id)%tgt_lat_pointer_dble,        &
00295                          Drv_Epios(id_epio_id)%tgt_lon_pointer_dble,        &
00296                          Drv_Epios(id_epio_id)%tgt_z_pointer_dble,          &
00297                          Drv_Epios(id_epio_id)%tgt_mask_pointer,            &
00298                          dl_gaus_var,                                       &
00299                          il_nb_neighbors,                                   &
00300                          Drv_Epios(id_epio_id)%index_array,                 &
00301                          Drv_Epios(id_epio_id)%weights,                     &
00302                          id_err)
00303                   END IF
00304               ENDIF
00305           CASE (PSMILe_trilinear)
00306 #ifdef VERBOSE
00307               PRINT *, &
00308                  '| | | | | | Trs is asked to interp trilinear'
00309               call psmile_flushstd
00310 #endif
00311               IF (Drv_Epios(id_epio_id)%src_coord_type .EQ. PRISM_Real) THEN
00312                   call prismtrs_trilinear_weight(                       &
00313                      il_src_size,                                       &
00314                      DBLE(Drv_Epios(id_epio_id)%src_lat_pointer_real),  &
00315                      DBLE(Drv_Epios(id_epio_id)%src_lon_pointer_real),  &
00316                      DBLE(Drv_Epios(id_epio_id)%src_z_pointer_real),    &
00317                      il_tgt_size,                                       & 
00318                      DBLE(Drv_Epios(id_epio_id)%tgt_lat_pointer_real),  &
00319                      DBLE(Drv_Epios(id_epio_id)%tgt_lon_pointer_real),  &
00320                      DBLE(Drv_Epios(id_epio_id)%tgt_z_pointer_real),    &
00321                      Drv_Epios(id_epio_id)%tgt_mask_pointer,            &
00322                      Drv_Epios(id_epio_id)%index_array,                 &
00323                      Drv_Epios(id_epio_id)%weights,                     &
00324                      id_err)
00325               ELSE IF (Drv_Epios(id_epio_id)%src_coord_type .EQ. &
00326                  PRISM_Double_Precision) THEN
00327                   call prismtrs_trilinear_weight(                       &
00328                      il_src_size,                                       &
00329                      Drv_Epios(id_epio_id)%src_lat_pointer_dble,        &
00330                      Drv_Epios(id_epio_id)%src_lon_pointer_dble,        &
00331                      Drv_Epios(id_epio_id)%src_z_pointer_dble,          &
00332                      il_tgt_size,                                       & 
00333                      Drv_Epios(id_epio_id)%tgt_lat_pointer_dble,        &
00334                      Drv_Epios(id_epio_id)%tgt_lon_pointer_dble,        &
00335                      Drv_Epios(id_epio_id)%tgt_z_pointer_dble,          &
00336                      Drv_Epios(id_epio_id)%tgt_mask_pointer,            &
00337                      Drv_Epios(id_epio_id)%index_array,                 &
00338                      Drv_Epios(id_epio_id)%weights,                     &
00339                      id_err)
00340               ENDIF
00341 
00342           CASE DEFAULT
00343               PRINT *, &
00344                  '| | | | | | Trs is not inteligent enough to understand'
00345               PRINT *, '| | | | | | Unknown interpolation method ?'
00346               call psmile_flushstd
00347 
00348           END SELECT
00349 
00350       ELSE IF (il_interp_type .EQ. PSMILe_2D1D) THEN
00351           IF (ila_interp_method(1) .EQ. PSMILe_conserv2D) THEN
00352               !
00353               ! Calculate the 2d conservative remapping weights
00354               ! SVXXX: The logic used for the 2D conservative remapping
00355               ! should be used for the other interpolations too instead
00356               ! of the complicated loops after line 404
00357               !
00358               il_index_size1 = SUM(Drv_Epios(id_epio_id)%nbsrccells_pertgtpt)
00359               IF (Drv_Epios(id_epio_id)%src_coord_type .EQ. PRISM_Real) THEN
00360                   call prismtrs_conserv2d_weight (                      &
00361                      il_src_size,                                       &
00362                      il_tgt_size,                                       &
00363                      Drv_Epios(id_epio_id)%src_nbr_corner,              &
00364                      Drv_Epios(id_epio_id)%tgt_nbr_corner,              &
00365                      Drv_Epios(id_epio_id)%src_lonlatz_size,            &
00366                      il_index_size1,                                    &
00367                      DBLE(Drv_Epios(id_epio_id)%src_lat_pointer_real),  &
00368                      DBLE(Drv_Epios(id_epio_id)%src_lon_pointer_real),  &
00369                      DBLE(Drv_Epios(id_epio_id)%src_z_pointer_real),    &
00370                      DBLE(Drv_Epios(id_epio_id)%tgt_lat_pointer_real),  &
00371                      DBLE(Drv_Epios(id_epio_id)%tgt_lon_pointer_real),  &
00372                      DBLE(Drv_Epios(id_epio_id)%tgt_z_pointer_real),    &
00373                      Drv_Epios(id_epio_id)%nbsrccells_pertgtpt,         &
00374                      Drv_Epios(id_epio_id)%index_array,                 &
00375                      Drv_Epios(id_epio_id)%srcepio_add,                 &
00376                      id_epio_id,                                        &
00377                      il_interp_id,                                      &
00378                      1,                                                 &
00379                      il_num_wts,                                        &
00380                      id_err)
00381                ELSE IF (Drv_Epios(id_epio_id)%src_coord_type .EQ. &
00382                  PRISM_Double_Precision) THEN
00383                   call prismtrs_conserv2d_weight (                      &
00384                      il_src_size,                                       &
00385                      il_tgt_size,                                       &
00386                      Drv_Epios(id_epio_id)%src_nbr_corner,              &
00387                      Drv_Epios(id_epio_id)%tgt_nbr_corner,              &
00388                      Drv_Epios(id_epio_id)%src_lonlatz_size,            &
00389                      il_index_size1,                                    &
00390                      Drv_Epios(id_epio_id)%src_lat_pointer_dble,        &
00391                      Drv_Epios(id_epio_id)%src_lon_pointer_dble,        &
00392                      Drv_Epios(id_epio_id)%src_z_pointer_dble,          &
00393                      Drv_Epios(id_epio_id)%tgt_lat_pointer_dble,        &
00394                      Drv_Epios(id_epio_id)%tgt_lon_pointer_dble,        &
00395                      Drv_Epios(id_epio_id)%tgt_z_pointer_dble,          &
00396                      Drv_Epios(id_epio_id)%nbsrccells_pertgtpt,         &
00397                      Drv_Epios(id_epio_id)%index_array,                 &
00398                      Drv_Epios(id_epio_id)%srcepio_add,                 &
00399                      id_epio_id,                                        &
00400                      il_interp_id,                                      &
00401                      1,                                                 &
00402                      il_num_wts,                                        &
00403                      id_err)
00404               ENDIF
00405 
00406 
00407               IF (ila_interp_method(2) .NE. PSMILe_none) THEN
00408                   !SV modify weigths to take vertical interpolation into account
00409                   PRINT *, '2D1D with something else that PSMILe_none'
00410                   PRINT *, 'is not implemented for 2d conservatice remapping'
00411               ENDIF
00412 
00413           ELSE IF (ila_interp_method(2) .EQ. PSMILe_none) THEN
00414 
00415               SELECT CASE (ila_interp_method(1))
00416                   
00417               CASE (PSMILe_nnghbr2D)
00418 
00419                   dl_gaus_var = Drv_Interps(il_interp_id)%arg8
00420                   IF ( dl_gaus_var .EQ. PSMILe_dundef) THEN
00421 #ifdef VERBOSE
00422                       PRINT *, &
00423                          '| | | | | | Trs is asked to interp distwght 2d'
00424                       call psmile_flushstd
00425 #endif
00426                       IF (Drv_Epios(id_epio_id)%src_coord_type .EQ.           &
00427                          PRISM_Real) THEN
00428                           call prismtrs_distwght_weight_2d(                   &
00429                              il_src_size,                                     &
00430                              DBLE(Drv_Epios(id_epio_id)%src_lat_pointer_real),&
00431                              DBLE(Drv_Epios(id_epio_id)%src_lon_pointer_real),&
00432                              Drv_Epios(id_epio_id)%src_mask_pointer,          &
00433                              il_tgt_size,                                     &
00434                              DBLE(Drv_Epios(id_epio_id)%tgt_lat_pointer_real),&
00435                              DBLE(Drv_Epios(id_epio_id)%tgt_lon_pointer_real),&
00436                              Drv_Epios(id_epio_id)%tgt_mask_pointer,          &
00437                              il_nb_neighbors,                                 &
00438                              Drv_Epios(id_epio_id)%index_array,               &
00439                              Drv_Epios(id_epio_id)%weights,                   &
00440                              id_err)
00441                       ELSE IF (Drv_Epios(id_epio_id)%src_coord_type .EQ.      &
00442                          PRISM_Double_Precision) THEN
00443                           call prismtrs_distwght_weight_2d(                   &
00444                              il_src_size,                                     &
00445                              Drv_Epios(id_epio_id)%src_lat_pointer_dble,      &
00446                              Drv_Epios(id_epio_id)%src_lon_pointer_dble,      &
00447                              Drv_Epios(id_epio_id)%src_mask_pointer,          &
00448                              il_tgt_size,                                     &
00449                              Drv_Epios(id_epio_id)%tgt_lat_pointer_dble,      &
00450                              Drv_Epios(id_epio_id)%tgt_lon_pointer_dble,      &
00451                              Drv_Epios(id_epio_id)%tgt_mask_pointer,          &
00452                              il_nb_neighbors,                                 &
00453                              Drv_Epios(id_epio_id)%index_array,               &
00454                              Drv_Epios(id_epio_id)%weights,                   &
00455                              id_err)
00456                       ENDIF
00457                   ELSE
00458                       
00459 #ifdef VERBOSE
00460                       PRINT *, &
00461                          '| | | | | | Trs is asked to interp gauswght 2d'
00462                       call psmile_flushstd
00463 #endif
00464                       IF (Drv_Epios(id_epio_id)%src_coord_type .EQ. &
00465                          PRISM_Real) THEN
00466                           call prismtrs_gauswght_weight_2d(                   &
00467                              il_src_size,                                     &
00468                              DBLE(Drv_Epios(id_epio_id)%src_lat_pointer_real),&
00469                              DBLE(Drv_Epios(id_epio_id)%src_lon_pointer_real),&
00470                              Drv_Epios(id_epio_id)%src_mask_pointer,          &
00471                              il_tgt_size,                                     &
00472                              DBLE(Drv_Epios(id_epio_id)%tgt_lat_pointer_real),&
00473                              DBLE(Drv_Epios(id_epio_id)%tgt_lon_pointer_real),&
00474                              Drv_Epios(id_epio_id)%tgt_mask_pointer,          &
00475                              dl_gaus_var,                                     &
00476                              il_nb_neighbors,                                 &
00477                              Drv_Epios(id_epio_id)%index_array,               &
00478                              Drv_Epios(id_epio_id)%weights,                   &
00479                              id_err)
00480                       ELSE IF (Drv_Epios(id_epio_id)%src_coord_type .EQ. &
00481                          PRISM_Double_Precision) THEN
00482                           call prismtrs_gauswght_weight_2d(                   &
00483                              il_src_size,                                     &
00484                              Drv_Epios(id_epio_id)%src_lat_pointer_dble,      &
00485                              Drv_Epios(id_epio_id)%src_lon_pointer_dble,      &
00486                              Drv_Epios(id_epio_id)%src_mask_pointer,          &
00487                              il_tgt_size,                                     &
00488                              Drv_Epios(id_epio_id)%tgt_lat_pointer_dble,      &
00489                              Drv_Epios(id_epio_id)%tgt_lon_pointer_dble,      &
00490                              Drv_Epios(id_epio_id)%tgt_mask_pointer,          &
00491                              dl_gaus_var,                                     &
00492                              il_nb_neighbors,                                 &
00493                              Drv_Epios(id_epio_id)%index_array,               &
00494                              Drv_Epios(id_epio_id)%weights,                   &
00495                              id_err)
00496                       ENDIF
00497                   END IF
00498 
00499               CASE (PSMILe_bilinear)
00500                   
00501 #ifdef VERBOSE
00502                   PRINT *, &
00503                      '| | | | | | Trs is asked to interp bilinear'
00504                   call psmile_flushstd
00505 #endif
00506                   IF (Drv_Epios(id_epio_id)%src_coord_type .EQ. &
00507                      PRISM_Real) THEN
00508                       call prismtrs_bilinear_weight_2d(                       &
00509                          il_src_size,                                         &
00510                          DBLE(Drv_Epios(id_epio_id)%src_lat_pointer_real),    &
00511                          DBLE(Drv_Epios(id_epio_id)%src_lon_pointer_real),    &
00512                          il_tgt_size,                                         &
00513                          DBLE(Drv_Epios(id_epio_id)%tgt_lat_pointer_real),    &
00514                          DBLE(Drv_Epios(id_epio_id)%tgt_lon_pointer_real),    &
00515                          Drv_Epios(id_epio_id)%tgt_mask_pointer,              &
00516                          Drv_Epios(id_epio_id)%index_array,                   &
00517                          Drv_Epios(id_epio_id)%weights,                       &
00518                          id_err)
00519                   ELSE IF (Drv_Epios(id_epio_id)%src_coord_type .EQ.          &
00520                      PRISM_Double_Precision) THEN
00521                       call prismtrs_bilinear_weight_2d(                       &
00522                          il_src_size,                                         &
00523                          Drv_Epios(id_epio_id)%src_lat_pointer_dble,          &
00524                          Drv_Epios(id_epio_id)%src_lon_pointer_dble,          &
00525                          il_tgt_size,                                         &
00526                          Drv_Epios(id_epio_id)%tgt_lat_pointer_dble,          &
00527                          Drv_Epios(id_epio_id)%tgt_lon_pointer_dble,          &
00528                          Drv_Epios(id_epio_id)%tgt_mask_pointer,              &
00529                          Drv_Epios(id_epio_id)%index_array,                   &
00530                          Drv_Epios(id_epio_id)%weights,                       &
00531                          id_err)
00532                   ENDIF
00533               CASE (PSMILe_bicubic)
00534                   
00535 #ifdef VERBOSE
00536                   PRINT *, &
00537                      '| | | | | | Trs is asked to interp bicubic'
00538                   call psmile_flushstd
00539 #endif
00540 
00541                   IF ( il_method == PSMILe_gradient ) THEN
00542                       IF (Drv_Epios(id_epio_id)%src_coord_type .EQ. &
00543                          PRISM_Real) THEN
00544                           call prismtrs_bicubic_grad_2d(                      &
00545                              il_src_size,                                     &
00546                              DBLE(Drv_Epios(id_epio_id)%src_lat_pointer_real),&
00547                              DBLE(Drv_Epios(id_epio_id)%src_lon_pointer_real),&
00548                              Drv_Epios(id_epio_id)%src_mask_pointer,          &
00549                              il_tgt_size,                                     &
00550                              DBLE(Drv_Epios(id_epio_id)%tgt_lat_pointer_real),&
00551                              DBLE(Drv_Epios(id_epio_id)%tgt_lon_pointer_real),&
00552                              Drv_Epios(id_epio_id)%tgt_mask_pointer,          &
00553                              Drv_Epios(id_epio_id)%index_array,               &
00554                              Drv_Epios(id_epio_id)%weights,                   &
00555                              id_err)
00556                       ELSE IF (Drv_Epios(id_epio_id)%src_coord_type .EQ.      &
00557                          PRISM_Double_Precision) THEN
00558                           call prismtrs_bicubic_grad_2d(                      &
00559                              il_src_size,                                     &
00560                              Drv_Epios(id_epio_id)%src_lat_pointer_dble,      &
00561                              Drv_Epios(id_epio_id)%src_lon_pointer_dble,      &
00562                              Drv_Epios(id_epio_id)%src_mask_pointer,          &
00563                              il_tgt_size,                                     &
00564                              Drv_Epios(id_epio_id)%tgt_lat_pointer_dble,      &
00565                              Drv_Epios(id_epio_id)%tgt_lon_pointer_dble,      &
00566                              Drv_Epios(id_epio_id)%tgt_mask_pointer,          &
00567                              Drv_Epios(id_epio_id)%index_array,               &
00568                              Drv_Epios(id_epio_id)%weights,                   &
00569                              id_err)
00570                       ENDIF
00571                   ELSE
00572                       IF (Drv_Epios(id_epio_id)%src_coord_type .EQ.           &
00573                          PRISM_Real) THEN
00574                           call prismtrs_bicubic_weight_2d(                    &
00575                              il_src_size,                                     &
00576                              DBLE(Drv_Epios(id_epio_id)%src_lat_pointer_real),&
00577                              DBLE(Drv_Epios(id_epio_id)%src_lon_pointer_real),&
00578                              il_tgt_size,                                     &
00579                              DBLE(Drv_Epios(id_epio_id)%tgt_lat_pointer_real),&
00580                              DBLE(Drv_Epios(id_epio_id)%tgt_lon_pointer_real),&
00581                              Drv_Epios(id_epio_id)%tgt_mask_pointer,          &
00582                              Drv_Epios(id_epio_id)%index_array,               &
00583                              Drv_Epios(id_epio_id)%same_lat,                  &
00584                              Drv_Epios(id_epio_id)%weights,                   &
00585                              id_err)
00586                       ELSE IF (Drv_Epios(id_epio_id)%src_coord_type .EQ.      &
00587                          PRISM_Double_Precision) THEN
00588                           call prismtrs_bicubic_weight_2d(                    &
00589                              il_src_size,                                     &
00590                              Drv_Epios(id_epio_id)%src_lat_pointer_dble,      &
00591                              Drv_Epios(id_epio_id)%src_lon_pointer_dble,      &
00592                              il_tgt_size,                                     &
00593                              Drv_Epios(id_epio_id)%tgt_lat_pointer_dble,      &
00594                              Drv_Epios(id_epio_id)%tgt_lon_pointer_dble,      &
00595                              Drv_Epios(id_epio_id)%tgt_mask_pointer,          &
00596                              Drv_Epios(id_epio_id)%index_array,               &
00597                              Drv_Epios(id_epio_id)%same_lat,                  &
00598                              Drv_Epios(id_epio_id)%weights,                   &
00599                              id_err)
00600                       ENDIF
00601                   ENDIF
00602 
00603               CASE DEFAULT
00604 
00605                   PRINT *, &
00606                      '| | | | | | Trs is not inteligent enough to understand'
00607                   PRINT *, '| | | | | | Unknown interpolation method ?'
00608                   call psmile_flushstd
00609 
00610               END SELECT
00611 
00612           ELSE !matches IF (ila_interp_method(2) .EQ. PSMILe_none) THEN
00613                ! under ELSE IF (il_interp_type .eq. PSMILe_2D1D) THEN
00614                ! Here it is supposed 2D1D with linear for z
00615 
00616               SELECT CASE (ila_interp_method(1))
00617                   
00618               CASE (PSMILe_nnghbr2D)
00619 
00620                   dl_gaus_var = Drv_Interps(il_interp_id)%arg8
00621                   IF ( dl_gaus_var .eq. PSMILe_dundef) THEN
00622 #ifdef VERBOSE
00623                       PRINT *, &
00624                          '| | | | | | Trs is asked to interp distwght 2d first'
00625                       call psmile_flushstd
00626 #endif
00627                       IF (Drv_Epios(id_epio_id)%src_coord_type .EQ.           &
00628                          PRISM_Real) THEN
00629                           call prismtrs_distwght_weight_2d1d(                 &
00630                              il_src_size,                                     &
00631                              DBLE(Drv_Epios(id_epio_id)%src_lat_pointer_real),&
00632                              DBLE(Drv_Epios(id_epio_id)%src_lon_pointer_real),&
00633                              DBLE(Drv_Epios(id_epio_id)%src_z_pointer_real),  &
00634                              Drv_Epios(id_epio_id)%src_mask_pointer,          &
00635                              il_tgt_size,                                     &
00636                              DBLE(Drv_Epios(id_epio_id)%tgt_lat_pointer_real),&
00637                              DBLE(Drv_Epios(id_epio_id)%tgt_lon_pointer_real),&
00638                              DBLE(Drv_Epios(id_epio_id)%tgt_z_pointer_real),  &
00639                              Drv_Epios(id_epio_id)%tgt_mask_pointer,          &
00640                              il_nb_neighbors,                                 &
00641                              Drv_Epios(id_epio_id)%index_array,               &
00642                              Drv_Epios(id_epio_id)%weights,                   &
00643                              id_err)
00644                       ELSE IF (Drv_Epios(id_epio_id)%src_coord_type .EQ.      &
00645                          PRISM_Double_Precision) THEN
00646                           call prismtrs_distwght_weight_2d1d(                 &
00647                              il_src_size,                                     &
00648                              Drv_Epios(id_epio_id)%src_lat_pointer_dble,      &
00649                              Drv_Epios(id_epio_id)%src_lon_pointer_dble,      &
00650                              Drv_Epios(id_epio_id)%src_z_pointer_dble,        &
00651                              Drv_Epios(id_epio_id)%src_mask_pointer,          &
00652                              il_tgt_size,                                     &
00653                              Drv_Epios(id_epio_id)%tgt_lat_pointer_dble,      &
00654                              Drv_Epios(id_epio_id)%tgt_lon_pointer_dble,      &
00655                              Drv_Epios(id_epio_id)%tgt_z_pointer_dble,        &
00656                              Drv_Epios(id_epio_id)%tgt_mask_pointer,          &
00657                              il_nb_neighbors,                                 &
00658                              Drv_Epios(id_epio_id)%index_array,               &
00659                              Drv_Epios(id_epio_id)%weights,                   &
00660                              id_err) 
00661 
00662                       ENDIF
00663                   ELSE
00664                       
00665 #ifdef VERBOSE
00666                       PRINT *, &
00667                          '| | | | | | Trs is asked to interp gauswght 2d first'
00668                       call psmile_flushstd
00669 #endif
00670                       PRINT *, &
00671                          '| | | | | | gausswgt_2d1d not yet supported'
00672 
00673                       call psmile_abort
00674                   END IF
00675 
00676               CASE (PSMILe_bilinear)
00677                   
00678 #ifdef VERBOSE
00679                   PRINT *, &
00680                      '| | | | | | Trs is asked to interp bilinear first'
00681                   call psmile_flushstd
00682 #endif
00683                   IF (Drv_Epios(id_epio_id)%src_coord_type .EQ.               &
00684                      PRISM_Real) THEN
00685                       call prismtrs_bilinear_weight_2d1d(                     &
00686                          il_src_size,                                         &
00687                          DBLE(Drv_Epios(id_epio_id)%src_lat_pointer_real),    &
00688                          DBLE(Drv_Epios(id_epio_id)%src_lon_pointer_real),    &
00689                          DBLE(Drv_Epios(id_epio_id)%src_z_pointer_real),      &
00690                          il_tgt_size,                                         &
00691                          DBLE(Drv_Epios(id_epio_id)%tgt_lat_pointer_real),    &
00692                          DBLE(Drv_Epios(id_epio_id)%tgt_lon_pointer_real),    &
00693                          DBLE(Drv_Epios(id_epio_id)%tgt_z_pointer_real),      &
00694                          Drv_Epios(id_epio_id)%index_array,                   &
00695                          Drv_Epios(id_epio_id)%weights,                       &
00696                          id_err)
00697                   ELSE IF (Drv_Epios(id_epio_id)%src_coord_type .EQ.          &
00698                      PRISM_Double_Precision) THEN
00699                       call prismtrs_bilinear_weight_2d1d(                     &
00700                          il_src_size,                                         &
00701                          Drv_Epios(id_epio_id)%src_lat_pointer_dble,          &
00702                          Drv_Epios(id_epio_id)%src_lon_pointer_dble,          &
00703                          Drv_Epios(id_epio_id)%src_z_pointer_dble,            &
00704                          il_tgt_size,                                         &
00705                          Drv_Epios(id_epio_id)%tgt_lat_pointer_dble,          &
00706                          Drv_Epios(id_epio_id)%tgt_lon_pointer_dble,          &
00707                          Drv_Epios(id_epio_id)%tgt_z_pointer_dble,            &
00708                          Drv_Epios(id_epio_id)%index_array,                   &
00709                          Drv_Epios(id_epio_id)%weights,                       &
00710                          id_err)
00711                   ENDIF
00712               CASE (PSMILe_bicubic)
00713                   
00714 #ifdef VERBOSE
00715                   PRINT *, &
00716                      '| | | | | | Trs is asked to interp bicubic first'
00717                   call psmile_flushstd
00718 #endif
00719 
00720                   PRINT *, &
00721                      '| | | | | | bicubic_2d1d not yet supported'
00722 
00723                   call psmile_abort
00724 
00725               CASE DEFAULT
00726 
00727                   PRINT *, &
00728                      '| | | | | | Trs is not inteligent enough to understand'
00729                   PRINT *, '| | | | | | Unknown interpolation method ?'
00730                   call psmile_flushstd
00731 
00732               END SELECT
00733 
00734               SELECT CASE (ila_interp_method(2))
00735                   
00736               CASE DEFAULT
00737                   IF (Drv_Epios(id_epio_id)%src_coord_type .EQ. &
00738                      PRISM_Real) THEN
00739                       call prismtrs_linear_weight_for_2d1d(                   &
00740                          il_src_size,                                         &
00741                          DBLE(Drv_Epios(id_epio_id)%src_z_pointer_real),      &
00742                          il_tgt_size,                                         &
00743                          DBLE(Drv_Epios(id_epio_id)%tgt_z_pointer_real),      &
00744                          Drv_Epios(id_epio_id)%tgt_mask_pointer,              &
00745                          il_nb_neighbors,                                     &
00746                          Drv_Epios(id_epio_id)%index_array,                   &
00747                          Drv_Epios(id_epio_id)%weights,  &
00748                          id_err)
00749                   ELSE IF (Drv_Epios(id_epio_id)%src_coord_type .EQ. &
00750                      PRISM_Double_Precision) THEN
00751                       call prismtrs_linear_weight_for_2d1d(                   &
00752                          il_src_size,                                         &
00753                          Drv_Epios(id_epio_id)%src_z_pointer_dble,            &
00754                          il_tgt_size,                                         &
00755                          Drv_Epios(id_epio_id)%tgt_z_pointer_dble ,           &
00756                          Drv_Epios(id_epio_id)%tgt_mask_pointer,              &
00757                          il_nb_neighbors,                                     &
00758                          Drv_Epios(id_epio_id)%index_array,                   &
00759                          Drv_Epios(id_epio_id)%weights,  &
00760                          id_err)
00761                   ENDIF
00762               END SELECT
00763 
00764           END IF
00765           
00766       ELSE
00767 
00768           PRINT *, '********************************************************'
00769           PRINT *, 'Interpolation type not supported for weights computation'
00770           PRINT *, '********************************************************'
00771           call psmile_flushstd
00772 
00773           CALL MPI_Abort(MPI_COMM_WORLD, 1,id_err)
00774 
00775       END IF
00776 
00777       Drv_Epios(id_epio_id)%weights_status = PSMILe_trans_complete
00778 
00779   ENDIF
00780 !
00781 ! 3. Apply the weights
00782 !
00783   IF (ila_interp_method(1) .EQ. PSMILe_conserv2D) THEN
00784       IF (Drv_Epios(id_epio_id)%max_links_map1 >= 1) &
00785          call prismtrs_apply_weights_2dcons(      &
00786          Drv_Epios(id_epio_id)%src_size,          &
00787          dda_trans_out,                           &
00788          Drv_Epios(id_epio_id)%tgt_size,          &
00789          dda_trans_in,                            &
00790          Drv_Epios(id_epio_id)%tgt_mask_pointer,  &
00791          Drv_Epios(id_epio_id)%max_links_map1,    &
00792          Drv_Epios(id_epio_id)%grid1_add_map1,    &
00793          Drv_Epios(id_epio_id)%grid2_add_map1,    &
00794          il_num_wts,                              &
00795          Drv_Epios(id_epio_id)%wts_map1,          &
00796          nbr_fields,                              &
00797          id_err)
00798   ELSE
00799       IF ( il_method == PSMILe_gradient ) THEN
00800 
00801           call prismtrs_apply_grads(il_src_size,      &
00802              dda_trans_out,                           &
00803              Drv_Epios(id_epio_id)%src_mask_pointer,  &
00804              il_tgt_size,                             &
00805              dda_trans_in,                            &
00806              Drv_Epios(id_epio_id)%tgt_mask_pointer,  &
00807              il_nb_neighbors,                         &
00808              Drv_Epios(id_epio_id)%index_array(:,1:il_nb_neighbors), &
00809              Drv_Epios(id_epio_id)%weights,           &
00810              nbr_fields,                              &
00811              id_err)
00812 
00813       ELSE 
00814 
00815           call prismtrs_apply_weights(il_src_size,    &
00816              dda_trans_out,                           &
00817              il_tgt_size,                             &
00818              dda_trans_in,                            &
00819              Drv_Epios(id_epio_id)%tgt_mask_pointer,  &  
00820              il_nb_neighbors,                         &
00821              Drv_Epios(id_epio_id)%index_array(:,1:il_nb_neighbors), &
00822              Drv_Epios(id_epio_id)%weights,           &
00823              nbr_fields,                              &
00824              id_err)
00825 
00826       ENDIF
00827   ENDIF
00828  
00829 #ifdef VERBOSE
00830   PRINT *, '| | | | | Quit PRISMTrs_Interp'
00831   call psmile_flushstd
00832 #endif
00833    
00834 END SUBROUTINE PRISMTrs_Interp

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1