|
||||
File indexing completed on 2018-03-02 18:41:31 UTC
view on githubraw file Latest commit add29e06 on 2018-01-31 20:35:05 UTCcc36b35673 Jean*0001 #include "GRIDALT_OPTIONS.h" 0002 8aa94f3b31 Andr*0003 subroutine phys2dyn(qphy,pephy,im1,im2,jm1,jm2,lmphy,Nsx,Nsy, 0004 . idim1,idim2,jdim1,jdim2,bi,bj,pedyn,Lbot,lmdyn,nlperdyn,qdyn) 0005 C*********************************************************************** 0006 C Purpose: 0007 C To interpolate an arbitrary quantity from the 'dynamics' eta (pstar) cc36b35673 Jean*0008 C grid to the higher resolution physics grid 8aa94f3b31 Andr*0009 C Algorithm: 0010 C Routine works one layer (edge to edge pressure) at a time. 0011 C Physics -> Dynamics computes the physics layer mean value, 616b995d4d Andr*0012 C weighted by dp**kappa (interp1) or by dp. 8aa94f3b31 Andr*0013 C 0014 C Input: 0015 C qphy..... [im,jm,lmphy] Arbitrary Quantity on Input Grid 0016 C pephy.... [im,jm,lmphy+1] Pressures at bottom edges of input levels 0017 C im1,2 ... Limits for Longitude Dimension of Input 0018 C jm1,2 ... Limits for Latitude Dimension of Input 0019 C lmphy.... Vertical Dimension of Input 0020 C Nsx...... Number of processes in x-direction 0021 C Nsy...... Number of processes in y-direction 0022 C idim1,2.. Beginning and ending i-values to calculate 0023 C jdim1,2.. Beginning and ending j-values to calculate 0024 C bi....... Index of process number in x-direction 0025 C bj....... Index of process number in x-direction 0026 C pedyn.... [im,jm,lmdyn+1] Pressures at bottom edges of output levels 0027 C lmdyn.... Vertical Dimension of Output 0028 C nlperdyn. Mapping Array-Highest Physics level in each dynmics level 0029 C 0030 C Output: 0031 C qdyn..... [im,jm,lmdyn] Quantity at output grid (physics grid) 0032 C 0033 C Notes: 0034 C 1) This algorithm assumes that the output (physics) grid levels 0035 C fit exactly into the input (dynamics) grid levels 0036 C*********************************************************************** 0037 implicit none 0038 0039 integer im1, im2, jm1, jm2, lmdyn, lmphy, Nsx, Nsy 0040 integer idim1, idim2, jdim1, jdim2, bi, bj 0041 _RL qphy(im1:im2,jm1:jm2,lmphy,Nsx,Nsy) 0042 _RL pedyn(im1:im2,jm1:jm2,lmdyn+1,Nsx,Nsy) 0043 _RL pephy(im1:im2,jm1:jm2,lmphy+1,Nsx,Nsy) 0044 integer nlperdyn(im1:im2,jm1:jm2,lmdyn,Nsx,Nsy) 0045 _RL qdyn(im1:im2,jm1:jm2,lmdyn,Nsx,Nsy) 0046 integer Lbot(im1:im2,jm1:jm2,Nsx,Nsy) 0047 be5ec7d59d Andr*0048 integer i,j,L,Lout1,Lout1p1,Lout2,Lphy 4d198cba86 Jean*0049 _RL dpkephy, dpkedyn, sum 8aa94f3b31 Andr*0050 4d198cba86 Jean*0051 cinterp1 _RL kappa 0052 #ifdef ALLOW_FIZHI 0053 cinterp1 _RL getcon 0054 #else 0055 cinterp1 #include 'SIZE.h' 0056 cinterp1 #include 'EEPARAMS.h' 0057 cinterp1 #include 'PARAMS.h' 0058 #endif 0059 0060 #ifdef ALLOW_FIZHI 0061 cinterp1 kappa = getcon('KAPPA') 0062 #else 0063 cinterp1 kappa = atm_kappa 0064 #endif 8aa94f3b31 Andr*0065 0066 c do loop for all dynamics (output) levels 8bb9f0d9db Andr*0067 do L = 1,lmdyn 0068 c do loop for all grid points 0069 do j = jdim1,jdim2 0070 do i = idim1,idim2 8aa94f3b31 Andr*0071 qdyn(i,j,L,bi,bj) = 0. 0072 c Check to make sure we are above ground - otherwise do nothing 0073 if(L.ge.Lbot(i,j,bi,bj))then 0074 if(L.eq.Lbot(i,j,bi,bj)) then 0075 Lout1 = 0 0076 else 0077 Lout1 = nlperdyn(i,j,L-1,bi,bj) 0078 endif 0079 Lout2 = nlperdyn(i,j,L,bi,bj) 0080 c do loop for all physics levels contained in this dynamics level 616b995d4d Andr*0081 cinterp1 dpkedyn = (pedyn(i,j,L,bi,bj)**kappa)- 0082 cinterp1 (pedyn(i,j,L+1,bi,bj)**kappa) 0083 dpkedyn = pedyn(i,j,L,bi,bj)-pedyn(i,j,L+1,bi,bj) 8bb9f0d9db Andr*0084 sum = 0. be5ec7d59d Andr*0085 Lout1p1 = Lout1+1 0086 do Lphy = Lout1p1,Lout2 616b995d4d Andr*0087 cinterp1 dpkephy = (pephy(i,j,Lphy,bi,bj)**kappa)- 0088 cinterp1 (pephy(i,j,Lphy+1,bi,bj)**kappa) 0089 dpkephy = pephy(i,j,Lphy,bi,bj)-pephy(i,j,Lphy+1,bi,bj) 8bb9f0d9db Andr*0090 sum=sum+qphy(i,j,Lphy,bi,bj)*(dpkephy/dpkedyn) 8aa94f3b31 Andr*0091 enddo 8bb9f0d9db Andr*0092 qdyn(i,j,L,bi,bj) = sum 8aa94f3b31 Andr*0093 endif 0094 enddo 0095 enddo 0096 enddo 0097 0098 return 0099 end
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated from https://github.com/MITgcm/MITgcm by the 2.2.1-MITgcm-0.1 LXR engine. The LXR team |