prismtrs_resize_remap_vars.F90

Go to the documentation of this file.
00001 !------------------------------------------------------------------------
00002 ! Copyright 2007-2010, CERFACS, Toulouse, France.
00003 ! All rights reserved. Use is subject to OASIS4 license terms.
00004 !-----------------------------------------------------------------------
00005 !BOP
00006 !
00007 ! !ROUTINE: PRISMTrs_resize_remap_vars
00008 !
00009 ! !INTERFACE!
00010       subroutine prismtrs_resize_remap_vars(increment, id_epio_id, num_wts)
00011 
00012 !-----------------------------------------------------------------------
00013 !
00014 !     input variables
00015 !
00016 !-----------------------------------------------------------------------
00017       USE PRISMDrv, dummy_INTERFACE => PRISMTrs_resize_remap_vars
00018       IMPLICIT NONE
00019       integer, intent(in) :: 
00020            increment,   ! the number of links to add(subtract) to arrays
00021            id_epio_id, 
00022            num_wts ! Number of weights for 2D conservative remapping
00023 !
00024 ! !DESCRIPTION
00025 !
00026 ! SUBROUTINE PRISMTrs_resize_remap_vars resizes remapping arrays by 
00027 ! increasing (decreasing) the max_links by increment. This routine
00028 ! comes from the SCRIP 1.4 library available from Los Alamos National 
00029 ! Laboratory.
00030 !
00031 ! !REVISED HISTORY
00032 !   Date      Programmer   Description
00033 ! ----------  ----------   -----------
00034 ! 20/03/2007  S. Valcke    Duplicated from the SCRIP library
00035 !
00036 ! EOP
00037 !----------------------------------------------------------------------
00038 ! $Id: prismtrs_resize_remap_vars.F90 2082 2009-10-21 13:31:19Z hanke $
00039 ! $Author: hanke $
00040 !----------------------------------------------------------------------
00041 !
00042 !-----------------------------------------------------------------------
00043 !
00044 !     local variables
00045 !
00046   CHARACTER(LEN=len_cvs_string), SAVE  :: mycvs = 
00047      '$Id: prismtrs_resize_remap_vars.F90 2082 2009-10-21 13:31:19Z hanke $'
00048 !-----------------------------------------------------------------------
00049 
00050       integer :: 
00051          ierr,       ! error flag
00052          mxlinks,    ! size of link arrays
00053          il_max_links_map1
00054       integer, dimension(:), allocatable :: 
00055          add1_tmp,  ! temp array for resizing address arrays
00056          add2_tmp  ! temp array for resizing address arrays
00057 
00058       double precision, dimension(:,:), allocatable :: 
00059          wts_tmp   ! temp array for resizing weight arrays
00060 
00061 #ifdef VERBOSE
00062   PRINT *, '| | | | | | | Enter PRISMTrs_resize_remap_vars'
00063   call psmile_flushstd
00064 #endif
00065 !-----------------------------------------------------------------------
00066 !
00067 !     resize map 1 arrays if required.
00068 !
00069 !-----------------------------------------------------------------------
00070 
00071         !***
00072         !*** allocate temporaries to hold original values
00073         !***
00074 
00075         mxlinks = size(Drv_Epios(id_epio_id)%grid1_add_map1)
00076         allocate (add1_tmp(mxlinks), add2_tmp(mxlinks), &
00077                   wts_tmp(num_wts,mxlinks))
00078 
00079         add1_tmp = Drv_Epios(id_epio_id)%grid1_add_map1
00080         add2_tmp = Drv_Epios(id_epio_id)%grid2_add_map1
00081         wts_tmp  = Drv_Epios(id_epio_id)%wts_map1
00082         
00083         !***
00084         !*** deallocate originals and increment max_links then
00085         !*** reallocate arrays at new size
00086         !***
00087 
00088         deallocate (Drv_Epios(id_epio_id)%grid1_add_map1, &
00089            Drv_Epios(id_epio_id)%grid2_add_map1, &
00090            Drv_Epios(id_epio_id)%wts_map1)
00091         Drv_Epios(id_epio_id)%max_links_map1 = mxlinks + increment
00092         il_max_links_map1 = Drv_Epios(id_epio_id)%max_links_map1
00093         allocate (Drv_Epios(id_epio_id)%grid1_add_map1(il_max_links_map1), &
00094                   Drv_Epios(id_epio_id)%grid2_add_map1(il_max_links_map1), &
00095                   Drv_Epios(id_epio_id)%wts_map1(num_wts,il_max_links_map1))
00096 
00097         !***
00098         !*** restore original values from temp arrays and
00099         !*** deallocate temps
00100         !***
00101 
00102         mxlinks = min(mxlinks, il_max_links_map1)
00103         Drv_Epios(id_epio_id)%grid1_add_map1(1:mxlinks) = add1_tmp (1:mxlinks)
00104         Drv_Epios(id_epio_id)%grid2_add_map1(1:mxlinks) = add2_tmp (1:mxlinks)
00105         Drv_Epios(id_epio_id)%wts_map1    (:,1:mxlinks) = wts_tmp(:,1:mxlinks)
00106         deallocate(add1_tmp, add2_tmp, wts_tmp)
00107 
00108 #ifdef VERBOSE
00109   PRINT *, '| | | | | | | Quit PRISMTrs_resize_remap_vars'
00110   call psmile_flushstd
00111 #endif
00112 !-----------------------------------------------------------------------
00113 
00114       end subroutine prismtrs_resize_remap_vars
00115 !----------------------------------------------------------------------- 

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1