New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
Changeset 7646 for trunk/NEMOGCM/NEMO/TOP_SRC/PISCES/P4Z/p4zflx.F90 – NEMO

Ignore:
Timestamp:
2017-02-06T10:25:03+01:00 (7 years ago)
Author:
timgraham
Message:

Merge of dev_merge_2016 into trunk. UPDATE TO ARCHFILES NEEDED for XIOS2.
LIM_SRC_s/limrhg.F90 to follow in next commit due to change of kind (I'm unable to do it in this commit).
Merged using the following steps:

1) svn merge --reintegrate svn+ssh://forge.ipsl.jussieu.fr/ipsl/forge/projets/nemo/svn/trunk .
2) Resolve minor conflicts in sette.sh and namelist_cfg for ORCA2LIM3 (due to a change in trunk after branch was created)
3) svn commit
4) svn switch svn+ssh://forge.ipsl.jussieu.fr/ipsl/forge/projets/nemo/svn/trunk
5) svn merge svn+ssh://forge.ipsl.jussieu.fr/ipsl/forge/projets/nemo/svn/branches/2016/dev_merge_2016 .
6) At this stage I checked out a clean copy of the branch to compare against what is about to be committed to the trunk.
6) svn commit #Commit code to the trunk

In this commit I have also reverted a change to Fcheck_archfile.sh which was causing problems on the Paris machine.

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/NEMOGCM/NEMO/TOP_SRC/PISCES/P4Z/p4zflx.F90

    r6962 r7646  
    1111   !!                  !  2011-02  (J. Simeon, J. Orr) Include total atm P correction  
    1212   !!---------------------------------------------------------------------- 
    13 #if defined key_pisces 
    14    !!---------------------------------------------------------------------- 
    15    !!   'key_pisces'                                       PISCES bio-model 
    16    !!---------------------------------------------------------------------- 
    1713   !!   p4z_flx       :   CALCULATES GAS EXCHANGE AND CHEMISTRY AT SEA SURFACE 
    1814   !!   p4z_flx_init  :   Read the namelist 
     
    2622   USE iom                          !  I/O manager 
    2723   USE fldread                      !  read input fields 
    28 #if defined key_cpl_carbon_cycle 
    29    USE sbc_oce, ONLY :  atm_co2     !  atmospheric pCO2                
    30 #endif 
    3124 
    3225   IMPLICIT NONE 
     
    4841 
    4942   !                               !!* nampisatm namelist (Atmospheric PRessure) * 
    50    LOGICAL, PUBLIC ::   ln_presatm  !: ref. pressure: global mean Patm (F) or a constant (F) 
    51  
    52    REAL(wp) , ALLOCATABLE, SAVE, DIMENSION(:,:)  ::  patm      ! atmospheric pressure at kt                 [N/m2] 
    53    TYPE(FLD), ALLOCATABLE,       DIMENSION(:)    ::  sf_patm   ! structure of input fields (file informations, fields read) 
    54  
    55  
    56    REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: oce_co2   !: ocean carbon flux  
     43   LOGICAL, PUBLIC ::   ln_presatm     !: ref. pressure: global mean Patm (F) or a constant (F) 
     44   LOGICAL, PUBLIC ::   ln_presatmco2  !: accounting for spatial atm CO2 in the compuation of carbon flux (T) or not (F) 
     45 
     46   REAL(wp) , ALLOCATABLE, SAVE, DIMENSION(:,:) ::  patm      ! atmospheric pressure at kt                 [N/m2] 
     47   TYPE(FLD), ALLOCATABLE,       DIMENSION(:)   ::  sf_patm   ! structure of input fields (file informations, fields read) 
     48   TYPE(FLD), ALLOCATABLE,       DIMENSION(:)   ::  sf_atmco2 ! structure of input fields (file informations, fields read) 
     49 
    5750   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: satmco2   !: atmospheric pco2  
    5851 
     
    7467      !! ** Method  :  
    7568      !!              - Include total atm P correction via Esbensen & Kushnir (1981)  
    76       !!              - Pressure correction NOT done for key_cpl_carbon_cycle 
    7769      !!              - Remove Wanninkhof chemical enhancement; 
    7870      !!              - Add option for time-interpolation of atcco2.txt   
     
    8577      REAL(wp) ::   zfld, zflu, zfld16, zflu16, zfact 
    8678      REAL(wp) ::   zvapsw, zsal, zfco2, zxc2, xCO2approx, ztkel, zfugcoeff 
    87       REAL(wp) ::   zph, zah2, zbot, zdic, zalk, zsch_o2, zalka, zsch_co2 
     79      REAL(wp) ::   zph, zdic, zsch_o2, zsch_co2 
    8880      REAL(wp) ::   zyr_dec, zdco2dt 
    8981      CHARACTER (len=25) :: charout 
     
    10092      !     IS USED TO COMPUTE AIR-SEA FLUX OF CO2 
    10193 
    102       IF( kt /= nit000 .AND. knt == 1 ) CALL p4z_patm( kt )    ! Get sea-level pressure (E&K [1981] climatology) for use in flux calcs 
    103  
    104       IF( ln_co2int ) THEN  
     94      IF( kt /= nit000 .AND. .NOT.l_co2cpl .AND. knt == 1 ) CALL p4z_patm( kt )    ! Get sea-level pressure (E&K [1981] climatology) for use in flux calcs 
     95 
     96      IF( ln_co2int .AND. .NOT.ln_presatmco2 .AND. .NOT.l_co2cpl ) THEN  
    10597         ! Linear temporal interpolation  of atmospheric pco2.  atcco2.txt has annual values. 
    10698         ! Caveats: First column of .txt must be in years, decimal  years preferably.  
     
    116108      ENDIF 
    117109 
    118 #if defined key_cpl_carbon_cycle 
    119       satmco2(:,:) = atm_co2(:,:) 
    120 #endif 
    121  
    122       DO jm = 1, 10 
    123          DO jj = 1, jpj 
    124             DO ji = 1, jpi 
    125  
    126                ! DUMMY VARIABLES FOR DIC, H+, AND BORATE 
    127                zbot  = borat(ji,jj,1) 
    128                zfact = rhop(ji,jj,1) / 1000. + rtrn 
    129                zdic  = trb(ji,jj,1,jpdic) / zfact 
    130                zph   = MAX( hi(ji,jj,1), 1.e-10 ) / zfact 
    131                zalka = trb(ji,jj,1,jptal) / zfact 
    132  
    133                ! CALCULATE [ALK]([CO3--], [HCO3-]) 
    134                zalk  = zalka - (  akw3(ji,jj,1) / zph - zph / aphscale(ji,jj,1)    & 
    135                &       + zbot / ( 1.+ zph / akb3(ji,jj,1) )  ) 
    136  
    137                ! CALCULATE [H+] AND [H2CO3] 
    138                zah2   = SQRT(  (zdic-zalk)**2 + 4.* ( zalk * ak23(ji,jj,1)   & 
    139                   &                                        / ak13(ji,jj,1) ) * ( 2.* zdic - zalk )  ) 
    140                zah2   = 0.5 * ak13(ji,jj,1) / zalk * ( ( zdic - zalk ) + zah2 ) 
    141                zh2co3(ji,jj) = ( 2.* zdic - zalk ) / ( 2.+ ak13(ji,jj,1) / zah2 ) * zfact 
    142                hi(ji,jj,1)   = zah2 * zfact 
    143             END DO 
     110      IF( l_co2cpl )   satmco2(:,:) = atm_co2(:,:) 
     111 
     112      DO jj = 1, jpj 
     113         DO ji = 1, jpi 
     114            ! DUMMY VARIABLES FOR DIC, H+, AND BORATE 
     115            zfact = rhop(ji,jj,1) / 1000. + rtrn 
     116            zdic  = trb(ji,jj,1,jpdic) 
     117            zph   = MAX( hi(ji,jj,1), 1.e-10 ) / zfact 
     118            ! CALCULATE [H2CO3] 
     119            zh2co3(ji,jj) = zdic/(1. + ak13(ji,jj,1)/zph + ak13(ji,jj,1)*ak23(ji,jj,1)/zph**2) 
    144120         END DO 
    145121      END DO 
    146  
    147122 
    148123      ! -------------- 
     
    167142            zkgwan = 0.251 * zws 
    168143            zkgwan = zkgwan * xconv * ( 1.- fr_i(ji,jj) ) * tmask(ji,jj,1) 
    169 # if defined key_degrad 
    170             zkgwan = zkgwan * facvol(ji,jj,1) 
    171 #endif  
    172144            ! compute gas exchange for CO2 and O2 
    173145            zkgco2(ji,jj) = zkgwan * SQRT( 660./ zsch_co2 ) 
     
    176148      END DO 
    177149 
     150 
    178151      DO jj = 1, jpj 
    179152         DO ji = 1, jpi 
    180             ztkel     = tsn(ji,jj,1,jp_tem) + 273.15 
    181             zsal      = tsn(ji,jj,1,jp_sal) + ( 1.- tmask(ji,jj,1) ) * 35. 
     153            ztkel = tempis(ji,jj,1) + 273.15 
     154            zsal  = salinprac(ji,jj,1) + ( 1.- tmask(ji,jj,1) ) * 35. 
    182155            zvapsw    = EXP(24.4543 - 67.4509*(100.0/ztkel) - 4.8489*LOG(ztkel/100) - 0.000544*zsal) 
    183156            zpco2atm(ji,jj) = satmco2(ji,jj) * ( patm(ji,jj) - zvapsw ) 
     
    232205         ENDIF 
    233206         IF( iom_use( "Dpo2" ) )  THEN 
    234            zw2d(:,:) = ( atcox * patm(:,:) - atcox * trn(:,:,1,jpoxy) / ( chemo2(:,:,1) + rtrn ) ) * tmask(:,:,1) 
     207           zw2d(:,:) = ( atcox * patm(:,:) - atcox * trb(:,:,1,jpoxy) / ( chemo2(:,:,1) + rtrn ) ) * tmask(:,:,1) 
    235208           CALL iom_put( "Dpo2"  , zw2d ) 
    236209         ENDIF 
     
    239212         ! 
    240213         CALL wrk_dealloc( jpi, jpj, zw2d ) 
    241       ELSE 
    242          IF( ln_diatrc ) THEN 
    243             trc2d(:,:,jp_pcs0_2d    ) = oce_co2(:,:) / e1e2t(:,:) * rfact2r  
    244             trc2d(:,:,jp_pcs0_2d + 1) = zoflx(:,:) * 1000 * tmask(:,:,1)  
    245             trc2d(:,:,jp_pcs0_2d + 2) = zkgco2(:,:) * tmask(:,:,1)  
    246             trc2d(:,:,jp_pcs0_2d + 3) = ( zpco2atm(:,:) - zh2co3(:,:) / ( chemc(:,:,1) + rtrn ) ) * tmask(:,:,1)  
    247          ENDIF 
    248214      ENDIF 
    249215      ! 
     
    287253         WRITE(numout,*) ' ' 
    288254      ENDIF 
    289       IF( .NOT.ln_co2int ) THEN 
     255     IF( .NOT.ln_co2int .AND. .NOT.ln_presatmco2 ) THEN 
    290256         IF(lwp) THEN                         ! control print 
    291257            WRITE(numout,*) '    Constant Atmospheric pCO2 value  atcco2    =', atcco2 
     
    293259         ENDIF 
    294260         satmco2(:,:)  = atcco2      ! Initialisation of atmospheric pco2 
    295       ELSE 
     261      ELSEIF( ln_co2int .AND. .NOT.ln_presatmco2 ) THEN 
    296262         IF(lwp)  THEN 
    297263            WRITE(numout,*) '    Atmospheric pCO2 value  from file clname      =', TRIM( clname ) 
     
    315281         END DO 
    316282         CLOSE(numco2) 
    317       ENDIF 
     283      ELSEIF( .NOT.ln_co2int .AND. ln_presatmco2 ) THEN 
     284         IF(lwp)  THEN 
     285            WRITE(numout,*) '    Spatialized Atmospheric pCO2 from an external file' 
     286            WRITE(numout,*) ' ' 
     287         ENDIF 
     288      ELSE 
     289         IF(lwp)  THEN 
     290            WRITE(numout,*) '    Spatialized Atmospheric pCO2 from an external file' 
     291            WRITE(numout,*) ' ' 
     292         ENDIF 
     293      ENDIF 
     294 
    318295      ! 
    319296      oce_co2(:,:)  = 0._wp                ! Initialization of Flux of Carbon 
     
    341318      CHARACTER(len=100) ::  cn_dir   ! Root directory for location of ssr files 
    342319      TYPE(FLD_N)        ::  sn_patm  ! informations about the fields to be read 
    343       !! 
    344       NAMELIST/nampisatm/ ln_presatm, sn_patm, cn_dir 
     320      TYPE(FLD_N)        ::  sn_atmco2 ! informations about the fields to be read 
     321      !! 
     322      NAMELIST/nampisatm/ ln_presatm, ln_presatmco2, sn_patm, sn_atmco2, cn_dir 
    345323 
    346324      !                                         ! ----------------------- ! 
     
    361339            WRITE(numout,*) '   Namelist nampisatm : Atmospheric Pressure as external forcing' 
    362340            WRITE(numout,*) '      constant atmopsheric pressure (F) or from a file (T)  ln_presatm = ', ln_presatm 
     341            WRITE(numout,*) '      spatial atmopsheric CO2 for flux calcs  ln_presatmco2 = ', ln_presatmco2 
    363342            WRITE(numout,*) 
    364343         ENDIF 
     
    373352         ENDIF 
    374353         !                                          
     354         IF( ln_presatmco2 ) THEN 
     355            ALLOCATE( sf_atmco2(1), STAT=ierr )           !* allocate and fill sf_atmco2 (forcing structure) with sn_atmco2 
     356            IF( ierr > 0 )   CALL ctl_stop( 'STOP', 'p4z_flx: unable to allocate sf_atmco2 structure' ) 
     357            ! 
     358            CALL fld_fill( sf_atmco2, (/ sn_atmco2 /), cn_dir, 'p4z_flx', 'Atmospheric co2 partial pressure ', 'nampisatm' ) 
     359                                   ALLOCATE( sf_atmco2(1)%fnow(jpi,jpj,1)   ) 
     360            IF( sn_atmco2%ln_tint )  ALLOCATE( sf_atmco2(1)%fdta(jpi,jpj,1,2) ) 
     361         ENDIF 
     362         ! 
    375363         IF( .NOT.ln_presatm )   patm(:,:) = 1.e0    ! Initialize patm if no reading from a file 
    376364         ! 
     
    382370      ENDIF 
    383371      ! 
     372      IF( ln_presatmco2 ) THEN 
     373         CALL fld_read( kt, 1, sf_atmco2 )               !* input atmco2 provided at kt + 1/2 
     374         satmco2(:,:) = sf_atmco2(1)%fnow(:,:,1)                        ! atmospheric pressure 
     375      ELSE 
     376         satmco2(:,:) = atcco2    ! Initialize atmco2 if no reading from a file 
     377      ENDIF 
     378      ! 
    384379   END SUBROUTINE p4z_patm 
    385380 
     381 
    386382   INTEGER FUNCTION p4z_flx_alloc() 
    387383      !!---------------------------------------------------------------------- 
    388384      !!                     ***  ROUTINE p4z_flx_alloc  *** 
    389385      !!---------------------------------------------------------------------- 
    390       ALLOCATE( oce_co2(jpi,jpj), satmco2(jpi,jpj), patm(jpi,jpj), STAT=p4z_flx_alloc ) 
     386      ALLOCATE( satmco2(jpi,jpj), patm(jpi,jpj), STAT=p4z_flx_alloc ) 
    391387      ! 
    392388      IF( p4z_flx_alloc /= 0 )   CALL ctl_warn('p4z_flx_alloc : failed to allocate arrays') 
    393389      ! 
    394390   END FUNCTION p4z_flx_alloc 
    395  
    396 #else 
    397    !!====================================================================== 
    398    !!  Dummy module :                                   No PISCES bio-model 
    399    !!====================================================================== 
    400 CONTAINS 
    401    SUBROUTINE p4z_flx( kt )                   ! Empty routine 
    402       INTEGER, INTENT( in ) ::   kt 
    403       WRITE(*,*) 'p4z_flx: You should not have seen this print! error?', kt 
    404    END SUBROUTINE p4z_flx 
    405 #endif  
    406391 
    407392   !!====================================================================== 
Note: See TracChangeset for help on using the changeset viewer.