prismtrs_store_link_cnsrv.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 ! !ROUTINE: PRISMTrs_store_link_cnsrv
00007 !
00008       subroutine prismtrs_store_link_cnsrv(add1, add2, lcl_weights, id_action, id_epio_id, num_wts)
00009 
00010 !-----------------------------------------------------------------------
00011 !
00012 !     input variables
00013 !
00014 !-----------------------------------------------------------------------
00015       USE PRISMDrv, dummy_INTERFACE => PRISMTrs_store_link_cnsrv
00016       IMPLICIT NONE
00017       integer, intent(in) :: 
00018               add1,   ! address on grid1
00019               add2,   ! address on grid2
00020               id_action,  ! 0 if first call, 2 if last call, 1 otherwise
00021               id_epio_id, 
00022               num_wts ! Number of weights for 2D conservative remapping
00023 
00024 
00025       DOUBLE PRECISION, dimension(2*num_wts), intent(in) :: 
00026               lcl_weights ! array of remapping lcl_weights for this link
00027 !
00028 ! !DESCRIPTION
00029 !
00030 ! SUBROUTINE PRISMTrs_store_link_cnsrv stores the address and weight for 
00031 ! this link in the appropriate address and weight arrays and resizes those
00032 ! arrays if necessary. This routine comes from the SCRIP 1.4 library 
00033 ! available from Los Alamos National Laboratory.
00034 !
00035 ! !REVISED HISTORY
00036 !   Date      Programmer   Description
00037 ! ----------  ----------   -----------
00038 ! 21/03/2007  S. Valcke    Duplicated from the SCRIP library
00039 !
00040 ! EOP
00041 !----------------------------------------------------------------------
00042 ! $Id: prismtrs_store_link_cnsrv.F90 2082 2009-10-21 13:31:19Z hanke $
00043 ! $Author: hanke $
00044 !----------------------------------------------------------------------
00045 !
00046 !-----------------------------------------------------------------------
00047 !
00048 !     local variables
00049 !
00050       CHARACTER(LEN=len_cvs_string), SAVE  :: mycvs = 
00051      '$Id: prismtrs_store_link_cnsrv.F90 2082 2009-10-21 13:31:19Z hanke $'
00052 
00053       integer :: nlink, min_link, max_link ! link index
00054       integer :: il_num_links_map1
00055       integer, dimension(:,:), allocatable, save :: 
00056               link_add1,   ! min,max link add to restrict search
00057               link_add2   ! min,max link add to restrict search
00058 
00059 #ifdef DEBUGX
00060   PRINT *, '| | | | | | | Enter PRISMTrs_store_link_cnsrv'
00061   call psmile_flushstd
00062 #endif
00063 !-----------------------------------------------------------------------
00064 !
00065 !     if all lcl_weights are zero, do not bother storing the link
00066 !
00067 !-----------------------------------------------------------------------
00068 
00069       !For oasis4, we keep next line as fracnnei is not an option in oasis4
00070 #ifdef DEBUGX
00071       WRITE(*,*) 'id_action, lcl_weights =', id_action, lcl_weights(:)
00072       call psmile_flushstd
00073 #endif
00074       if (all(lcl_weights == zero)) return
00075       ! If store_link_cnsrv is called only for deallocation, do it
00076       IF (id_action == 2) THEN
00077           DEALLOCATE(link_add1, link_add2)
00078           RETURN
00079       ENDIF
00080 !-----------------------------------------------------------------------
00081 !
00082 !     restrict the range of links to search for existing links
00083 !
00084 !-----------------------------------------------------------------------
00085       ! if (first_call) then
00086       if (id_action == 0) then
00087         ! if (.not. first_conserv) then
00088         !  deallocate(link_add1, link_add2)
00089         ! endif
00090         allocate(link_add1(2,Drv_Epios(id_epio_id)%src_size), &
00091            link_add2(2,Drv_Epios(id_epio_id)%tgt_size))
00092         link_add1 = 0
00093         link_add2 = 0
00094         ! first_call = .false.
00095         min_link = 1
00096         max_link = 0
00097       else
00098         min_link = min(link_add1(1,add1),link_add2(1,add2))
00099         max_link = max(link_add1(2,add1),link_add2(2,add2))
00100         if (min_link == 0) then
00101           min_link = 1
00102           max_link = 0
00103         endif
00104       endif
00105 
00106 !-----------------------------------------------------------------------
00107 !
00108 !     if the link already exists, add the weight to the current weight
00109 !     arrays
00110 !
00111 !-----------------------------------------------------------------------
00112       do nlink=min_link,max_link
00113         if (add1 ==  Drv_Epios(id_epio_id)%grid1_add_map1(nlink)) then
00114         if (add2 ==  Drv_Epios(id_epio_id)%grid2_add_map1(nlink)) then
00115 
00116           Drv_Epios(id_epio_id)%wts_map1(:,nlink) = &
00117              Drv_Epios(id_epio_id)%wts_map1(:,nlink) + lcl_weights(1:num_wts)
00118           return
00119 
00120         endif
00121         endif
00122       end do
00123 
00124 !-----------------------------------------------------------------------
00125 !
00126 !     if the link does not yet exist, increment number of links and 
00127 !     check to see if remap arrays need to be increased to accomodate 
00128 !     the new link.  then store the link.
00129 !
00130 !-----------------------------------------------------------------------
00131       il_num_links_map1 = Drv_Epios(id_epio_id)%num_links_map1 + 1
00132       if (il_num_links_map1 > Drv_Epios(id_epio_id)%max_links_map1)  &
00133          call prismtrs_resize_remap_vars(Drv_Epios(id_epio_id)%resize_increment, id_epio_id, num_wts)
00134       Drv_Epios(id_epio_id)%grid1_add_map1(il_num_links_map1) = add1
00135       Drv_Epios(id_epio_id)%grid2_add_map1(il_num_links_map1) = add2
00136       Drv_Epios(id_epio_id)%wts_map1    (:,il_num_links_map1) = lcl_weights(1:num_wts)
00137 
00138       if (link_add1(1,add1) == 0) link_add1(1,add1) = il_num_links_map1
00139       if (link_add2(1,add2) == 0) link_add2(1,add2) = il_num_links_map1
00140       link_add1(2,add1) = il_num_links_map1
00141       link_add2(2,add2) = il_num_links_map1
00142 
00143       Drv_Epios(id_epio_id)%num_links_map1  = il_num_links_map1
00144 
00145 #ifdef DEBUGX
00146   PRINT *, '| | | | | | | Quit PRISMTrs_store_link_cnsrv'
00147   call psmile_flushstd
00148 #endif
00149 !-----------------------------------------------------------------------
00150 
00151       end subroutine PRISMTrs_store_link_cnsrv
00152 

Generated on 18 Mar 2011 for Oasis4 by  doxygen 1.6.1