00001 module mrg_x2i_mct
00002
00003 use shr_kind_mod, only: r8 => shr_kind_r8, cl => shr_kind_cl
00004 use mct_mod
00005 use seq_flds_mod
00006 use seq_flds_indices
00007 use seq_comm_mct
00008 use seq_cdata_mod
00009 use seq_infodata_mod
00010
00011 implicit none
00012 save
00013 private
00014
00015
00016
00017
00018
00019 public :: mrg_x2i_init_mct
00020 public :: mrg_x2i_run_mct
00021 public :: mrg_x2i_final_mct
00022
00023
00024
00025
00026
00027
00028 contains
00029
00030
00031 subroutine mrg_x2i_init_mct( cdata_i, a2x_i, o2x_i )
00032
00033 type(seq_cdata), intent(in) :: cdata_i
00034 type(mct_aVect), intent(inout) :: a2x_i
00035 type(mct_aVect), intent(inout) :: o2x_i
00036
00037 type(mct_gsMap), pointer :: gsMap_ice
00038 integer :: mpicom
00039
00040
00041 call seq_cdata_setptrs(cdata_i,gsMap=gsMap_ice,mpicom=mpicom)
00042
00043
00044
00045 call MCT_aVect_init(a2x_i, rList=seq_flds_a2x_fields, &
00046 lsize=mct_gsMap_lsize(gsMap_ice, mpicom))
00047 call MCT_aVect_zero(a2x_i)
00048
00049
00050
00051 call MCT_aVect_init(o2x_i, rList=seq_flds_o2x_fields, &
00052 lsize=mct_gsMap_lsize(gsMap_ice, mpicom))
00053 call MCT_aVect_zero(o2x_i)
00054
00055 end subroutine mrg_x2i_init_mct
00056
00057
00058
00059 subroutine mrg_x2i_run_mct( cdata_i, a2x_i, o2x_i, x2i_i )
00060
00061
00062
00063
00064
00065 type(seq_cdata),intent(in) :: cdata_i
00066 type(mct_aVect),intent(in) :: a2x_i
00067 type(mct_aVect),intent(in) :: o2x_i
00068 type(mct_aVect),intent(out):: x2i_i
00069
00070
00071
00072 logical :: usevector
00073 integer :: i
00074 real(r8):: flux_epbalfact
00075 character(len=cl) :: flux_epbal
00076 type(seq_infodata_type),pointer :: infodata
00077
00078
00079 character(*),parameter :: F01 = "('(mrg_x2i_run_mct) ',a,3e11.3,a,f9.6)"
00080 character(*),parameter :: subName = '(mrg_x2i_run_mct) '
00081
00082
00083
00084
00085 #ifdef CPP_VECTOR
00086 usevector = .true.
00087 #else
00088 usevector = .false.
00089 #endif
00090
00091 call seq_cdata_setptrs(cdata_i,infodata=infodata)
00092
00093 call mct_aVect_copy(aVin=o2x_i, aVout=x2i_i, vector=usevector)
00094 call mct_aVect_copy(aVin=a2x_i, aVout=x2i_i, vector=usevector)
00095
00096
00097
00098 call seq_infodata_GetData(infodata, flux_epbal=flux_epbal)
00099 flux_epbalfact = 1.0_r8
00100 if (trim(flux_epbal) == 'ocn') then
00101 call seq_infodata_GetData(infodata, precip_fact = flux_epbalfact)
00102 flux_epbalfact = flux_epbalfact * 1.0e-6_r8
00103 end if
00104 if (flux_epbalfact <= 0.0_R8) then
00105 if (loglevel > 0) write(logunit,F01) 'WARNING: factor from ocn = ',flux_epbalfact
00106 if (loglevel > 0) write(logunit,F01) 'WARNING: resetting flux_epbalfact to 1.0'
00107 flux_epbalfact = 1.0_r8
00108 end if
00109
00110
00111 do i = 1,mct_aVect_lsize(x2i_i)
00112 x2i_i%rAttr(index_x2i_Faxa_rain,i) = a2x_i%rAttr(index_a2x_Faxa_rainc,i) + &
00113 a2x_i%rAttr(index_a2x_Faxa_rainl,i)
00114 x2i_i%rAttr(index_x2i_Faxa_snow,i) = a2x_i%rAttr(index_a2x_Faxa_snowc,i) + &
00115 a2x_i%rAttr(index_a2x_Faxa_snowl,i)
00116
00117
00118
00119
00120 x2i_i%rAttr(index_x2i_Faxa_rain,i) = x2i_i%rAttr(index_x2i_Faxa_rain,i) * flux_epbalfact
00121 x2i_i%rAttr(index_x2i_Faxa_snow,i) = x2i_i%rAttr(index_x2i_Faxa_snow,i) * flux_epbalfact
00122
00123 end do
00124
00125 end subroutine mrg_x2i_run_mct
00126
00127
00128
00129 subroutine mrg_x2i_final_mct
00130
00131
00132
00133 end subroutine mrg_x2i_final_mct
00134
00135 end module mrg_x2i_mct
00136
00137