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   IF (add2 == 4848) THEN
00072       WRITE(88,*), 'Enter PRISMTrs_store_link_cnsrv'
00073       WRITE(88,*) 'add1, add2, id_action, lcl_weights =', add1,add2,id_action, lcl_weights(:)
00074   ENDIF
00075       call psmile_flushstd
00076 #endif
00077       ! If store_link_cnsrv is called only for deallocation, do it
00078       IF (id_action == 2) THEN
00079           DEALLOCATE(link_add1, link_add2)
00080           RETURN
00081       ENDIF
00082 !-----------------------------------------------------------------------
00083 !
00084 !     restrict the range of links to search for existing links
00085 !
00086 !-----------------------------------------------------------------------
00087       ! if (first_call) then
00088       if (id_action == 0) then
00089         ! if (.not. first_conserv) then
00090         !  deallocate(link_add1, link_add2)
00091         ! endif
00092         allocate(link_add1(2,Drv_Epios(id_epio_id)%src_size), &
00093            link_add2(2,Drv_Epios(id_epio_id)%tgt_size))
00094         link_add1 = 0
00095         link_add2 = 0
00096         ! first_call = .false.
00097         min_link = 1
00098         max_link = 0
00099 
00100         if (all(lcl_weights == zero)) return
00101 
00102       else
00103 
00104         if (all(lcl_weights == zero)) return
00105 
00106         min_link = min(link_add1(1,add1),link_add2(1,add2))
00107         max_link = max(link_add1(2,add1),link_add2(2,add2))
00108         if (min_link == 0) then
00109           min_link = 1
00110           max_link = 0
00111         endif
00112       endif
00113 
00114 !-----------------------------------------------------------------------
00115 !
00116 !     if the link already exists, add the weight to the current weight
00117 !     arrays
00118 !
00119 !-----------------------------------------------------------------------
00120       do nlink=min_link,max_link
00121         if (add1 ==  Drv_Epios(id_epio_id)%grid1_add_map1(nlink)) then
00122         if (add2 ==  Drv_Epios(id_epio_id)%grid2_add_map1(nlink)) then
00123 
00124           Drv_Epios(id_epio_id)%wts_map1(:,nlink) = &
00125              Drv_Epios(id_epio_id)%wts_map1(:,nlink) + lcl_weights(1:num_wts)
00126           return
00127 
00128         endif
00129         endif
00130       end do
00131 
00132 !-----------------------------------------------------------------------
00133 !
00134 !     if the link does not yet exist, increment number of links and 
00135 !     check to see if remap arrays need to be increased to accomodate 
00136 !     the new link.  then store the link.
00137 !
00138 !-----------------------------------------------------------------------
00139       il_num_links_map1 = Drv_Epios(id_epio_id)%num_links_map1 + 1
00140       if (il_num_links_map1 > Drv_Epios(id_epio_id)%max_links_map1)  &
00141          call prismtrs_resize_remap_vars(Drv_Epios(id_epio_id)%resize_increment, id_epio_id, num_wts)
00142       Drv_Epios(id_epio_id)%grid1_add_map1(il_num_links_map1) = add1
00143       Drv_Epios(id_epio_id)%grid2_add_map1(il_num_links_map1) = add2
00144       Drv_Epios(id_epio_id)%wts_map1    (:,il_num_links_map1) = lcl_weights(1:num_wts)
00145 
00146       if (link_add1(1,add1) == 0) link_add1(1,add1) = il_num_links_map1
00147       if (link_add2(1,add2) == 0) link_add2(1,add2) = il_num_links_map1
00148       link_add1(2,add1) = il_num_links_map1
00149       link_add2(2,add2) = il_num_links_map1
00150 
00151       Drv_Epios(id_epio_id)%num_links_map1  = il_num_links_map1
00152 
00153 #ifdef DEBUGX
00154   PRINT *, '| | | | | | | Quit PRISMTrs_store_link_cnsrv'
00155   call psmile_flushstd
00156 #endif
00157 !-----------------------------------------------------------------------
00158 
00159       end subroutine PRISMTrs_store_link_cnsrv
00160 

Generated on 1 Dec 2011 for Oasis4 by  doxygen 1.6.1