Back to home page

MITgcm

 
 

    


File indexing completed on 2018-03-02 18:45:26 UTC

view on githubraw file Latest commit add29e06 on 2018-01-31 20:35:05 UTC
fad027e390 Ed H*0001 #include "FIZHI_OPTIONS.h"
0815b638b3 Jean*0002       subroutine fizhi_tendency_apply_u(
                0003      U                        gU_arr,
                0004      I                        iMin,iMax,jMin,jMax, kLev, bi, bj,
                0005      I                        myTime, myIter, myThid )
fad027e390 Ed H*0006 C=======================================================================
                0007 C Routine: fizhi_tendency_apply_u
                0008 C     Interpolate tendencies from physics grid to dynamics grid and
                0009 C     add fizhi tendency terms to U tendency.
0815b638b3 Jean*0010 C
                0011 C INPUT:
fad027e390 Ed H*0012 C     iMin - Working range of tile for applying forcing.
                0013 C     iMax
                0014 C     jMin
                0015 C     jMax
                0016 C     kLev
                0017 C
                0018 C Notes: Routine works for one level at a time
                0019 C        Assumes that U and V tendencies are already on C-Grid
                0020 C=======================================================================
                0021       implicit none
                0022 
                0023 #include "SIZE.h"
                0024 #include "GRID.h"
                0025 #include "EEPARAMS.h"
                0026 #include "DYNVARS.h"
                0027 #include "fizhi_SIZE.h"
                0028 #include "fizhi_land_SIZE.h"
                0029 #include "fizhi_coms.h"
                0030 
0815b638b3 Jean*0031       _RL     gU_arr(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
                0032       INTEGER iMin, iMax, jMin, jMax
                0033       INTEGER kLev, bi, bj
                0034       _RL     myTime
                0035       INTEGER myIter
                0036       INTEGER myThid
fad027e390 Ed H*0037 
                0038       integer i, j
                0039 
                0040       do j=jMin,jMax
                0041        do i=iMin,iMax
0815b638b3 Jean*0042         gU_arr(i,j) = gU_arr(i,j) +
fad027e390 Ed H*0043      .      maskW(i,j,kLev,bi,bj) * guphy(i,j,kLev,bi,bj)
                0044        enddo
                0045       enddo
                0046 
                0047       return
                0048       end
0815b638b3 Jean*0049       subroutine fizhi_tendency_apply_v(
                0050      U                        gV_arr,
                0051      I                        iMin,iMax,jMin,jMax, kLev, bi, bj,
                0052      I                        myTime, myIter, myThid )
fad027e390 Ed H*0053 C=======================================================================
                0054 C Routine: fizhi_tendency_apply_v
                0055 C     Interpolate tendencies from physics grid to dynamics grid and
                0056 C     add fizhi tendency terms to V tendency.
0815b638b3 Jean*0057 C
                0058 C INPUT:
fad027e390 Ed H*0059 C     iMin - Working range of tile for applying forcing.
                0060 C     iMax
                0061 C     jMin
                0062 C     jMax
                0063 C     kLev
                0064 C
                0065 C Notes: Routine works for one level at a time
                0066 C        Assumes that U and V tendencies are already on C-Grid
                0067 C=======================================================================
                0068       implicit none
                0069 
                0070 #include "SIZE.h"
                0071 #include "GRID.h"
                0072 #include "EEPARAMS.h"
                0073 #include "DYNVARS.h"
                0074 #include "fizhi_SIZE.h"
                0075 #include "fizhi_land_SIZE.h"
                0076 #include "fizhi_coms.h"
                0077 
0815b638b3 Jean*0078       _RL     gV_arr(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
                0079       INTEGER iMin, iMax, jMin, jMax
                0080       INTEGER kLev, bi, bj
                0081       _RL     myTime
                0082       INTEGER myIter
                0083       INTEGER myThid
fad027e390 Ed H*0084 
                0085       integer i, j
                0086 
                0087       do j=jMin,jMax
                0088        do i=iMin,iMax
0815b638b3 Jean*0089         gV_arr(i,j) = gV_arr(i,j) +
fad027e390 Ed H*0090      .      maskS(i,j,kLev,bi,bj) * gvphy(i,j,kLev,bi,bj)
                0091        enddo
                0092       enddo
                0093 
                0094       return
                0095       end
0815b638b3 Jean*0096       subroutine fizhi_tendency_apply_t(
                0097      U                        gT_arr,
                0098      I                        iMin,iMax,jMin,jMax, kLev, bi, bj,
                0099      I                        myTime, myIter, myThid )
fad027e390 Ed H*0100 C=======================================================================
                0101 C Routine: fizhi_tendency_apply_t
                0102 C     Interpolate tendencies from physics grid to dynamics grid and
                0103 C     add fizhi tendency terms to T (theta) tendency.
0815b638b3 Jean*0104 C
                0105 C INPUT:
fad027e390 Ed H*0106 C     iMin - Working range of tile for applying forcing.
                0107 C     iMax
                0108 C     jMin
                0109 C     jMax
                0110 C     kLev
                0111 C
                0112 C Notes: Routine works for one level at a time
                0113 C=======================================================================
                0114       implicit none
                0115 
                0116 #include "SIZE.h"
                0117 #include "GRID.h"
                0118 #include "EEPARAMS.h"
                0119 #include "DYNVARS.h"
                0120 #include "fizhi_SIZE.h"
                0121 #include "fizhi_land_SIZE.h"
                0122 #include "fizhi_coms.h"
                0123 
0815b638b3 Jean*0124       _RL     gT_arr(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
                0125       INTEGER iMin, iMax, jMin, jMax
                0126       INTEGER kLev, bi, bj
                0127       _RL     myTime
                0128       INTEGER myIter
                0129       INTEGER myThid
fad027e390 Ed H*0130 
                0131       integer i, j
                0132 
                0133       do j=jMin,jMax
                0134        do i=iMin,iMax
0815b638b3 Jean*0135         gT_arr(i,j) = maskC(i,j,kLev,bi,bj)
                0136      .       *( gT_arr(i,j) + gthphy(i,j,kLev,bi,bj) )
fad027e390 Ed H*0137        enddo
                0138       enddo
                0139 
                0140       return
                0141       end
0815b638b3 Jean*0142       subroutine fizhi_tendency_apply_s(
                0143      U                        gS_arr,
                0144      I                        iMin,iMax,jMin,jMax, kLev, bi, bj,
                0145      I                        myTime, myIter, myThid )
fad027e390 Ed H*0146 C=======================================================================
                0147 C Routine: fizhi_tendency_apply_s
                0148 C     Interpolate tendencies from physics grid to dynamics grid and
                0149 C     add fizhi tendency terms to S tendency.
0815b638b3 Jean*0150 C
                0151 C INPUT:
fad027e390 Ed H*0152 C     iMin - Working range of tile for applying forcing.
                0153 C     iMax
                0154 C     jMin
                0155 C     jMax
                0156 C     kLev
                0157 C
                0158 C Notes: Routine works for one level at a time
                0159 C=======================================================================
                0160       implicit none
                0161 
                0162 #include "SIZE.h"
                0163 #include "GRID.h"
                0164 #include "EEPARAMS.h"
                0165 #include "DYNVARS.h"
                0166 #include "fizhi_SIZE.h"
                0167 #include "fizhi_land_SIZE.h"
                0168 #include "fizhi_coms.h"
                0169 
0815b638b3 Jean*0170       _RL     gS_arr(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
                0171       INTEGER iMin, iMax, jMin, jMax
                0172       INTEGER kLev, bi, bj
                0173       _RL     myTime
                0174       INTEGER myIter
                0175       INTEGER myThid
fad027e390 Ed H*0176 
                0177       integer i, j
                0178 
                0179       do j=jMin,jMax
                0180        do i=iMin,iMax
0815b638b3 Jean*0181         gS_arr(i,j) = maskC(i,j,kLev,bi,bj)
                0182      .       *( gS_arr(i,j) + gsphy(i,j,kLev,bi,bj) )
fad027e390 Ed H*0183        enddo
                0184       enddo
                0185 
                0186       return
                0187       end