Changeset 5407 for trunk/NEMOGCM/NEMO/OPA_SRC/SBC/sbccpl.F90
- Timestamp:
- 2015-06-11T21:13:22+02:00 (9 years ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/NEMOGCM/NEMO/OPA_SRC/SBC/sbccpl.F90
r5363 r5407 21 21 USE sbc_oce ! Surface boundary condition: ocean fields 22 22 USE sbc_ice ! Surface boundary condition: ice fields 23 USE sbcapr 23 24 USE sbcdcy ! surface boundary condition: diurnal cycle 24 25 USE phycst ! physical constants … … 32 33 USE cpl_oasis3 ! OASIS3 coupling 33 34 USE geo2ocean ! 34 USE oce , ONLY : tsn, un, vn 35 USE oce , ONLY : tsn, un, vn, sshn, ub, vb, sshb, fraqsr_1lev 35 36 USE albedo ! 36 37 USE in_out_manager ! I/O manager … … 40 41 USE timing ! Timing 41 42 USE lbclnk ! ocean lateral boundary conditions (or mpp link) 43 USE eosbn2 44 USE sbcrnf , ONLY : l_rnfcpl 42 45 #if defined key_cpl_carbon_cycle 43 46 USE p4zflx, ONLY : oce_co2 … … 46 49 USE ice_domain_size, only: ncat 47 50 #endif 51 #if defined key_lim3 52 USE limthd_dh ! for CALL lim_thd_snwblow 53 #endif 54 48 55 IMPLICIT NONE 49 56 PRIVATE 50 !EM XIOS-OASIS-MCT compliance 57 51 58 PUBLIC sbc_cpl_init ! routine called by sbcmod.F90 52 59 PUBLIC sbc_cpl_rcv ! routine called by sbc_ice_lim(_2).F90 … … 89 96 INTEGER, PARAMETER :: jpr_topm = 32 ! topmeltn 90 97 INTEGER, PARAMETER :: jpr_botm = 33 ! botmeltn 91 INTEGER, PARAMETER :: jprcv = 33 ! total number of fields received 92 93 INTEGER, PARAMETER :: jps_fice = 1 ! ice fraction 98 INTEGER, PARAMETER :: jpr_sflx = 34 ! salt flux 99 INTEGER, PARAMETER :: jpr_toce = 35 ! ocean temperature 100 INTEGER, PARAMETER :: jpr_soce = 36 ! ocean salinity 101 INTEGER, PARAMETER :: jpr_ocx1 = 37 ! ocean current on grid 1 102 INTEGER, PARAMETER :: jpr_ocy1 = 38 ! 103 INTEGER, PARAMETER :: jpr_ssh = 39 ! sea surface height 104 INTEGER, PARAMETER :: jpr_fice = 40 ! ice fraction 105 INTEGER, PARAMETER :: jpr_e3t1st = 41 ! first T level thickness 106 INTEGER, PARAMETER :: jpr_fraqsr = 42 ! fraction of solar net radiation absorbed in the first ocean level 107 INTEGER, PARAMETER :: jprcv = 42 ! total number of fields received 108 109 INTEGER, PARAMETER :: jps_fice = 1 ! ice fraction sent to the atmosphere 94 110 INTEGER, PARAMETER :: jps_toce = 2 ! ocean temperature 95 111 INTEGER, PARAMETER :: jps_tice = 3 ! ice temperature … … 106 122 INTEGER, PARAMETER :: jps_ivz1 = 14 ! 107 123 INTEGER, PARAMETER :: jps_co2 = 15 108 INTEGER, PARAMETER :: jpsnd = 15 ! total number of fields sended 124 INTEGER, PARAMETER :: jps_soce = 16 ! ocean salinity 125 INTEGER, PARAMETER :: jps_ssh = 17 ! sea surface height 126 INTEGER, PARAMETER :: jps_qsroce = 18 ! Qsr above the ocean 127 INTEGER, PARAMETER :: jps_qnsoce = 19 ! Qns above the ocean 128 INTEGER, PARAMETER :: jps_oemp = 20 ! ocean freshwater budget (evap - precip) 129 INTEGER, PARAMETER :: jps_sflx = 21 ! salt flux 130 INTEGER, PARAMETER :: jps_otx1 = 22 ! 2 atmosphere-ocean stress components on grid 1 131 INTEGER, PARAMETER :: jps_oty1 = 23 ! 132 INTEGER, PARAMETER :: jps_rnf = 24 ! runoffs 133 INTEGER, PARAMETER :: jps_taum = 25 ! wind stress module 134 INTEGER, PARAMETER :: jps_fice2 = 26 ! ice fraction sent to OPA (by SAS when doing SAS-OPA coupling) 135 INTEGER, PARAMETER :: jps_e3t1st = 27 ! first level depth (vvl) 136 INTEGER, PARAMETER :: jps_fraqsr = 28 ! fraction of solar net radiation absorbed in the first ocean level 137 INTEGER, PARAMETER :: jpsnd = 28 ! total number of fields sended 109 138 110 139 ! !!** namelist namsbc_cpl ** … … 125 154 LOGICAL :: ln_usecplmask ! use a coupling mask file to merge data received from several models 126 155 ! -> file cplmask.nc with the float variable called cplmask (jpi,jpj,nn_cplmodel) 127 128 REAL(wp), ALLOCATABLE, DIMENSION(:,:,:) :: xcplmask129 130 156 TYPE :: DYNARR 131 157 REAL(wp), POINTER, DIMENSION(:,:,:) :: z3 … … 139 165 140 166 !! Substitution 167 # include "domzgr_substitute.h90" 141 168 # include "vectopt_loop_substitute.h90" 142 169 !!---------------------------------------------------------------------- … … 161 188 ALLOCATE( a_i(jpi,jpj,1) , STAT=ierr(2) ) ! used in sbcice_if.F90 (done here as there is no sbc_ice_if_init) 162 189 #endif 163 ALLOCATE( xcplmask(jpi,jpj, nn_cplmodel) , STAT=ierr(3) )190 ALLOCATE( xcplmask(jpi,jpj,0:nn_cplmodel) , STAT=ierr(3) ) 164 191 ! 165 192 sbc_cpl_alloc = MAXVAL( ierr ) … … 182 209 !! * initialise the OASIS coupler 183 210 !!---------------------------------------------------------------------- 184 INTEGER, INTENT(in) :: k_ice ! ice management in the sbc (=0/1/2/3)211 INTEGER, INTENT(in) :: k_ice ! ice management in the sbc (=0/1/2/3) 185 212 !! 186 213 INTEGER :: jn ! dummy loop index … … 216 243 WRITE(numout,*)'sbc_cpl_init : namsbc_cpl namelist ' 217 244 WRITE(numout,*)'~~~~~~~~~~~~' 245 ENDIF 246 IF( lwp .AND. ln_cpl ) THEN ! control print 218 247 WRITE(numout,*)' received fields (mutiple ice categogies)' 219 248 WRITE(numout,*)' 10m wind module = ', TRIM(sn_rcv_w10m%cldes ), ' (', TRIM(sn_rcv_w10m%clcat ), ')' … … 359 388 srcv(jpr_oemp)%clname = 'OOEvaMPr' ! ocean water budget = ocean Evap - ocean precip 360 389 SELECT CASE( TRIM( sn_rcv_emp%cldes ) ) 390 CASE( 'none' ) ! nothing to do 361 391 CASE( 'oce only' ) ; srcv( jpr_oemp )%laction = .TRUE. 362 392 CASE( 'conservative' ) … … 370 400 ! ! Runoffs & Calving ! 371 401 ! ! ------------------------- ! 372 srcv(jpr_rnf )%clname = 'O_Runoff' ; IF( TRIM( sn_rcv_rnf%cldes ) == 'coupled' ) srcv(jpr_rnf)%laction = .TRUE. 373 ! This isn't right - really just want ln_rnf_emp changed 374 ! IF( TRIM( sn_rcv_rnf%cldes ) == 'climato' ) THEN ; ln_rnf = .TRUE. 375 ! ELSE ; ln_rnf = .FALSE. 376 ! ENDIF 402 srcv(jpr_rnf )%clname = 'O_Runoff' 403 IF( TRIM( sn_rcv_rnf%cldes ) == 'coupled' ) THEN 404 srcv(jpr_rnf)%laction = .TRUE. 405 l_rnfcpl = .TRUE. ! -> no need to read runoffs in sbcrnf 406 ln_rnf = nn_components /= jp_iam_sas ! -> force to go through sbcrnf if not sas 407 IF(lwp) WRITE(numout,*) 408 IF(lwp) WRITE(numout,*) ' runoffs received from oasis -> force ln_rnf = ', ln_rnf 409 ENDIF 410 ! 377 411 srcv(jpr_cal )%clname = 'OCalving' ; IF( TRIM( sn_rcv_cal%cldes ) == 'coupled' ) srcv(jpr_cal)%laction = .TRUE. 378 412 … … 384 418 srcv(jpr_qnsmix)%clname = 'O_QnsMix' 385 419 SELECT CASE( TRIM( sn_rcv_qns%cldes ) ) 420 CASE( 'none' ) ! nothing to do 386 421 CASE( 'oce only' ) ; srcv( jpr_qnsoce )%laction = .TRUE. 387 422 CASE( 'conservative' ) ; srcv( (/jpr_qnsice, jpr_qnsmix/) )%laction = .TRUE. … … 399 434 srcv(jpr_qsrmix)%clname = 'O_QsrMix' 400 435 SELECT CASE( TRIM( sn_rcv_qsr%cldes ) ) 436 CASE( 'none' ) ! nothing to do 401 437 CASE( 'oce only' ) ; srcv( jpr_qsroce )%laction = .TRUE. 402 438 CASE( 'conservative' ) ; srcv( (/jpr_qsrice, jpr_qsrmix/) )%laction = .TRUE. … … 414 450 ! 415 451 ! non solar sensitivity mandatory for LIM ice model 416 IF( TRIM( sn_rcv_dqnsdt%cldes ) == 'none' .AND. k_ice /= 0 .AND. k_ice /= 4 ) &452 IF( TRIM( sn_rcv_dqnsdt%cldes ) == 'none' .AND. k_ice /= 0 .AND. k_ice /= 4 .AND. nn_components /= jp_iam_sas ) & 417 453 CALL ctl_stop( 'sbc_cpl_init: sn_rcv_dqnsdt%cldes must be coupled in namsbc_cpl namelist' ) 418 454 ! non solar sensitivity mandatory for mixed oce-ice solar radiation coupling technique … … 447 483 srcv(jpr_topm:jpr_botm)%laction = .TRUE. 448 484 ENDIF 449 450 ! Allocate all parts of frcv used for received fields 485 ! ! ------------------------------- ! 486 ! ! OPA-SAS coupling - rcv by opa ! 487 ! ! ------------------------------- ! 488 srcv(jpr_sflx)%clname = 'O_SFLX' 489 srcv(jpr_fice)%clname = 'RIceFrc' 490 ! 491 IF( nn_components == jp_iam_opa ) THEN ! OPA coupled to SAS via OASIS: force received field by OPA (sent by SAS) 492 srcv(:)%laction = .FALSE. ! force default definition in case of opa <-> sas coupling 493 srcv(:)%clgrid = 'T' ! force default definition in case of opa <-> sas coupling 494 srcv(:)%nsgn = 1. ! force default definition in case of opa <-> sas coupling 495 srcv( (/jpr_qsroce, jpr_qnsoce, jpr_oemp, jpr_sflx, jpr_fice, jpr_otx1, jpr_oty1, jpr_taum/) )%laction = .TRUE. 496 srcv(jpr_otx1)%clgrid = 'U' ! oce components given at U-point 497 srcv(jpr_oty1)%clgrid = 'V' ! and V-point 498 ! Vectors: change of sign at north fold ONLY if on the local grid 499 srcv( (/jpr_otx1,jpr_oty1/) )%nsgn = -1. 500 sn_rcv_tau%clvgrd = 'U,V' 501 sn_rcv_tau%clvor = 'local grid' 502 sn_rcv_tau%clvref = 'spherical' 503 sn_rcv_emp%cldes = 'oce only' 504 ! 505 IF(lwp) THEN ! control print 506 WRITE(numout,*) 507 WRITE(numout,*)' Special conditions for SAS-OPA coupling ' 508 WRITE(numout,*)' OPA component ' 509 WRITE(numout,*) 510 WRITE(numout,*)' received fields from SAS component ' 511 WRITE(numout,*)' ice cover ' 512 WRITE(numout,*)' oce only EMP ' 513 WRITE(numout,*)' salt flux ' 514 WRITE(numout,*)' mixed oce-ice solar flux ' 515 WRITE(numout,*)' mixed oce-ice non solar flux ' 516 WRITE(numout,*)' wind stress U,V on local grid and sperical coordinates ' 517 WRITE(numout,*)' wind stress module' 518 WRITE(numout,*) 519 ENDIF 520 ENDIF 521 ! ! -------------------------------- ! 522 ! ! OPA-SAS coupling - rcv by sas ! 523 ! ! -------------------------------- ! 524 srcv(jpr_toce )%clname = 'I_SSTSST' 525 srcv(jpr_soce )%clname = 'I_SSSal' 526 srcv(jpr_ocx1 )%clname = 'I_OCurx1' 527 srcv(jpr_ocy1 )%clname = 'I_OCury1' 528 srcv(jpr_ssh )%clname = 'I_SSHght' 529 srcv(jpr_e3t1st)%clname = 'I_E3T1st' 530 srcv(jpr_fraqsr)%clname = 'I_FraQsr' 531 ! 532 IF( nn_components == jp_iam_sas ) THEN 533 IF( .NOT. ln_cpl ) srcv(:)%laction = .FALSE. ! force default definition in case of opa <-> sas coupling 534 IF( .NOT. ln_cpl ) srcv(:)%clgrid = 'T' ! force default definition in case of opa <-> sas coupling 535 IF( .NOT. ln_cpl ) srcv(:)%nsgn = 1. ! force default definition in case of opa <-> sas coupling 536 srcv( (/jpr_toce, jpr_soce, jpr_ssh, jpr_fraqsr, jpr_ocx1, jpr_ocy1/) )%laction = .TRUE. 537 srcv( jpr_e3t1st )%laction = lk_vvl 538 srcv(jpr_ocx1)%clgrid = 'U' ! oce components given at U-point 539 srcv(jpr_ocy1)%clgrid = 'V' ! and V-point 540 ! Vectors: change of sign at north fold ONLY if on the local grid 541 srcv(jpr_ocx1:jpr_ocy1)%nsgn = -1. 542 ! Change first letter to couple with atmosphere if already coupled OPA 543 ! this is nedeed as each variable name used in the namcouple must be unique: 544 ! for example O_Runoff received by OPA from SAS and therefore O_Runoff received by SAS from the Atmosphere 545 DO jn = 1, jprcv 546 IF ( srcv(jn)%clname(1:1) == "O" ) srcv(jn)%clname = "S"//srcv(jn)%clname(2:LEN(srcv(jn)%clname)) 547 END DO 548 ! 549 IF(lwp) THEN ! control print 550 WRITE(numout,*) 551 WRITE(numout,*)' Special conditions for SAS-OPA coupling ' 552 WRITE(numout,*)' SAS component ' 553 WRITE(numout,*) 554 IF( .NOT. ln_cpl ) THEN 555 WRITE(numout,*)' received fields from OPA component ' 556 ELSE 557 WRITE(numout,*)' Additional received fields from OPA component : ' 558 ENDIF 559 WRITE(numout,*)' sea surface temperature (Celcius) ' 560 WRITE(numout,*)' sea surface salinity ' 561 WRITE(numout,*)' surface currents ' 562 WRITE(numout,*)' sea surface height ' 563 WRITE(numout,*)' thickness of first ocean T level ' 564 WRITE(numout,*)' fraction of solar net radiation absorbed in the first ocean level' 565 WRITE(numout,*) 566 ENDIF 567 ENDIF 568 569 ! =================================================== ! 570 ! Allocate all parts of frcv used for received fields ! 571 ! =================================================== ! 451 572 DO jn = 1, jprcv 452 573 IF ( srcv(jn)%laction ) ALLOCATE( frcv(jn)%z3(jpi,jpj,srcv(jn)%nct) ) … … 454 575 ! Allocate taum part of frcv which is used even when not received as coupling field 455 576 IF ( .NOT. srcv(jpr_taum)%laction ) ALLOCATE( frcv(jpr_taum)%z3(jpi,jpj,srcv(jpr_taum)%nct) ) 577 ! Allocate w10m part of frcv which is used even when not received as coupling field 578 IF ( .NOT. srcv(jpr_w10m)%laction ) ALLOCATE( frcv(jpr_w10m)%z3(jpi,jpj,srcv(jpr_w10m)%nct) ) 579 ! Allocate jpr_otx1 part of frcv which is used even when not received as coupling field 580 IF ( .NOT. srcv(jpr_otx1)%laction ) ALLOCATE( frcv(jpr_otx1)%z3(jpi,jpj,srcv(jpr_otx1)%nct) ) 581 IF ( .NOT. srcv(jpr_oty1)%laction ) ALLOCATE( frcv(jpr_oty1)%z3(jpi,jpj,srcv(jpr_oty1)%nct) ) 456 582 ! Allocate itx1 and ity1 as they are used in sbc_cpl_ice_tau even if srcv(jpr_itx1)%laction = .FALSE. 457 583 IF( k_ice /= 0 ) THEN … … 485 611 CASE default ; CALL ctl_stop( 'sbc_cpl_init: wrong definition of sn_snd_temp%cldes' ) 486 612 END SELECT 487 613 488 614 ! ! ------------------------- ! 489 615 ! ! Albedo ! … … 518 644 IF ( TRIM( sn_snd_thick%clcat ) == 'yes' ) ssnd(jps_fice)%nct = jpl 519 645 ENDIF 520 646 521 647 SELECT CASE ( TRIM( sn_snd_thick%cldes ) ) 522 648 CASE( 'none' ) ! nothing to do … … 567 693 ! ! ------------------------- ! 568 694 ssnd(jps_co2)%clname = 'O_CO2FLX' ; IF( TRIM(sn_snd_co2%cldes) == 'coupled' ) ssnd(jps_co2 )%laction = .TRUE. 695 696 ! ! ------------------------------- ! 697 ! ! OPA-SAS coupling - snd by opa ! 698 ! ! ------------------------------- ! 699 ssnd(jps_ssh )%clname = 'O_SSHght' 700 ssnd(jps_soce )%clname = 'O_SSSal' 701 ssnd(jps_e3t1st)%clname = 'O_E3T1st' 702 ssnd(jps_fraqsr)%clname = 'O_FraQsr' 703 ! 704 IF( nn_components == jp_iam_opa ) THEN 705 ssnd(:)%laction = .FALSE. ! force default definition in case of opa <-> sas coupling 706 ssnd( (/jps_toce, jps_soce, jps_ssh, jps_fraqsr, jps_ocx1, jps_ocy1/) )%laction = .TRUE. 707 ssnd( jps_e3t1st )%laction = lk_vvl 708 ! vector definition: not used but cleaner... 709 ssnd(jps_ocx1)%clgrid = 'U' ! oce components given at U-point 710 ssnd(jps_ocy1)%clgrid = 'V' ! and V-point 711 sn_snd_crt%clvgrd = 'U,V' 712 sn_snd_crt%clvor = 'local grid' 713 sn_snd_crt%clvref = 'spherical' 714 ! 715 IF(lwp) THEN ! control print 716 WRITE(numout,*) 717 WRITE(numout,*)' sent fields to SAS component ' 718 WRITE(numout,*)' sea surface temperature (T before, Celcius) ' 719 WRITE(numout,*)' sea surface salinity ' 720 WRITE(numout,*)' surface currents U,V on local grid and spherical coordinates' 721 WRITE(numout,*)' sea surface height ' 722 WRITE(numout,*)' thickness of first ocean T level ' 723 WRITE(numout,*)' fraction of solar net radiation absorbed in the first ocean level' 724 WRITE(numout,*) 725 ENDIF 726 ENDIF 727 ! ! ------------------------------- ! 728 ! ! OPA-SAS coupling - snd by sas ! 729 ! ! ------------------------------- ! 730 ssnd(jps_sflx )%clname = 'I_SFLX' 731 ssnd(jps_fice2 )%clname = 'IIceFrc' 732 ssnd(jps_qsroce)%clname = 'I_QsrOce' 733 ssnd(jps_qnsoce)%clname = 'I_QnsOce' 734 ssnd(jps_oemp )%clname = 'IOEvaMPr' 735 ssnd(jps_otx1 )%clname = 'I_OTaux1' 736 ssnd(jps_oty1 )%clname = 'I_OTauy1' 737 ssnd(jps_rnf )%clname = 'I_Runoff' 738 ssnd(jps_taum )%clname = 'I_TauMod' 739 ! 740 IF( nn_components == jp_iam_sas ) THEN 741 IF( .NOT. ln_cpl ) ssnd(:)%laction = .FALSE. ! force default definition in case of opa <-> sas coupling 742 ssnd( (/jps_qsroce, jps_qnsoce, jps_oemp, jps_fice2, jps_sflx, jps_otx1, jps_oty1, jps_taum/) )%laction = .TRUE. 743 ! 744 ! Change first letter to couple with atmosphere if already coupled with sea_ice 745 ! this is nedeed as each variable name used in the namcouple must be unique: 746 ! for example O_SSTSST sent by OPA to SAS and therefore S_SSTSST sent by SAS to the Atmosphere 747 DO jn = 1, jpsnd 748 IF ( ssnd(jn)%clname(1:1) == "O" ) ssnd(jn)%clname = "S"//ssnd(jn)%clname(2:LEN(ssnd(jn)%clname)) 749 END DO 750 ! 751 IF(lwp) THEN ! control print 752 WRITE(numout,*) 753 IF( .NOT. ln_cpl ) THEN 754 WRITE(numout,*)' sent fields to OPA component ' 755 ELSE 756 WRITE(numout,*)' Additional sent fields to OPA component : ' 757 ENDIF 758 WRITE(numout,*)' ice cover ' 759 WRITE(numout,*)' oce only EMP ' 760 WRITE(numout,*)' salt flux ' 761 WRITE(numout,*)' mixed oce-ice solar flux ' 762 WRITE(numout,*)' mixed oce-ice non solar flux ' 763 WRITE(numout,*)' wind stress U,V components' 764 WRITE(numout,*)' wind stress module' 765 ENDIF 766 ENDIF 767 569 768 ! 570 769 ! ================================ ! … … 572 771 ! ================================ ! 573 772 574 CALL cpl_define(jprcv, jpsnd,nn_cplmodel) 773 CALL cpl_define(jprcv, jpsnd, nn_cplmodel) 774 575 775 IF (ln_usecplmask) THEN 576 776 xcplmask(:,:,:) = 0. … … 582 782 xcplmask(:,:,:) = 1. 583 783 ENDIF 584 ! 585 IF( ln_dm2dc .AND. ( cpl_freq( jpr_qsroce ) + cpl_freq( jpr_qsrmix ) /= 86400 ) ) & 784 xcplmask(:,:,0) = 1. - SUM( xcplmask(:,:,1:nn_cplmodel), dim = 3 ) 785 ! 786 ncpl_qsr_freq = cpl_freq( 'O_QsrOce' ) + cpl_freq( 'O_QsrMix' ) + cpl_freq( 'S_QsrOce' ) + cpl_freq( 'S_QsrMix' ) 787 IF( ln_dm2dc .AND. ln_cpl .AND. ncpl_qsr_freq /= 86400 ) & 586 788 & CALL ctl_stop( 'sbc_cpl_init: diurnal cycle reconstruction (ln_dm2dc) needs daily couping for solar radiation' ) 789 ncpl_qsr_freq = 86400 / ncpl_qsr_freq 587 790 588 791 CALL wrk_dealloc( jpi,jpj, zacs, zaos ) … … 638 841 !! emp upward mass flux [evap. - precip. (- runoffs) (- calving)] (ocean only case) 639 842 !!---------------------------------------------------------------------- 640 INTEGER, INTENT(in) :: kt ! ocean model time step index 641 INTEGER, INTENT(in) :: k_fsbc ! frequency of sbc (-> ice model) computation 642 INTEGER, INTENT(in) :: k_ice ! ice management in the sbc (=0/1/2/3) 643 !! 644 LOGICAL :: llnewtx, llnewtau ! update wind stress components and module?? 843 INTEGER, INTENT(in) :: kt ! ocean model time step index 844 INTEGER, INTENT(in) :: k_fsbc ! frequency of sbc (-> ice model) computation 845 INTEGER, INTENT(in) :: k_ice ! ice management in the sbc (=0/1/2/3) 846 847 !! 848 LOGICAL :: llnewtx, llnewtau ! update wind stress components and module?? 645 849 INTEGER :: ji, jj, jn ! dummy loop indices 646 850 INTEGER :: isec ! number of seconds since nit000 (assuming rdttra did not change since nit000) … … 650 854 REAL(wp) :: zcdrag = 1.5e-3 ! drag coefficient 651 855 REAL(wp) :: zzx, zzy ! temporary variables 652 REAL(wp), POINTER, DIMENSION(:,:) :: ztx, zty 856 REAL(wp), POINTER, DIMENSION(:,:) :: ztx, zty, zmsk, zemp, zqns, zqsr 653 857 !!---------------------------------------------------------------------- 654 858 ! 655 859 IF( nn_timing == 1 ) CALL timing_start('sbc_cpl_rcv') 656 860 ! 657 CALL wrk_alloc( jpi,jpj, ztx, zty ) 658 ! ! Receive all the atmos. fields (including ice information) 659 isec = ( kt - nit000 ) * NINT( rdttra(1) ) ! date of exchanges 660 DO jn = 1, jprcv ! received fields sent by the atmosphere 661 IF( srcv(jn)%laction ) CALL cpl_rcv( jn, isec, frcv(jn)%z3, xcplmask, nrcvinfo(jn) ) 861 CALL wrk_alloc( jpi,jpj, ztx, zty, zmsk, zemp, zqns, zqsr ) 862 ! 863 IF( ln_mixcpl ) zmsk(:,:) = 1. - xcplmask(:,:,0) 864 ! 865 ! ! ======================================================= ! 866 ! ! Receive all the atmos. fields (including ice information) 867 ! ! ======================================================= ! 868 isec = ( kt - nit000 ) * NINT( rdttra(1) ) ! date of exchanges 869 DO jn = 1, jprcv ! received fields sent by the atmosphere 870 IF( srcv(jn)%laction ) CALL cpl_rcv( jn, isec, frcv(jn)%z3, xcplmask(:,:,1:nn_cplmodel), nrcvinfo(jn) ) 662 871 END DO 663 872 … … 719 928 ! 720 929 ENDIF 721 722 930 ! ! ========================= ! 723 931 ! ! wind stress module ! (taum) … … 748 956 ENDIF 749 957 ENDIF 750 958 ! 751 959 ! ! ========================= ! 752 960 ! ! 10 m wind speed ! (wndm) … … 761 969 !CDIR NOVERRCHK 762 970 DO ji = 1, jpi 763 wndm(ji,jj) = SQRT( frcv(jpr_taum)%z3(ji,jj,1) * zcoef )971 frcv(jpr_w10m)%z3(ji,jj,1) = SQRT( frcv(jpr_taum)%z3(ji,jj,1) * zcoef ) 764 972 END DO 765 973 END DO 766 974 ENDIF 767 ELSE768 IF ( nrcvinfo(jpr_w10m) == OASIS_Rcv ) wndm(:,:) = frcv(jpr_w10m)%z3(:,:,1)769 975 ENDIF 770 976 … … 773 979 IF( MOD( kt-1, k_fsbc ) == 0 ) THEN 774 980 ! 775 utau(:,:) = frcv(jpr_otx1)%z3(:,:,1) 776 vtau(:,:) = frcv(jpr_oty1)%z3(:,:,1) 777 taum(:,:) = frcv(jpr_taum)%z3(:,:,1) 981 IF( ln_mixcpl ) THEN 982 utau(:,:) = utau(:,:) * xcplmask(:,:,0) + frcv(jpr_otx1)%z3(:,:,1) * zmsk(:,:) 983 vtau(:,:) = vtau(:,:) * xcplmask(:,:,0) + frcv(jpr_oty1)%z3(:,:,1) * zmsk(:,:) 984 taum(:,:) = taum(:,:) * xcplmask(:,:,0) + frcv(jpr_taum)%z3(:,:,1) * zmsk(:,:) 985 wndm(:,:) = wndm(:,:) * xcplmask(:,:,0) + frcv(jpr_w10m)%z3(:,:,1) * zmsk(:,:) 986 ELSE 987 utau(:,:) = frcv(jpr_otx1)%z3(:,:,1) 988 vtau(:,:) = frcv(jpr_oty1)%z3(:,:,1) 989 taum(:,:) = frcv(jpr_taum)%z3(:,:,1) 990 wndm(:,:) = frcv(jpr_w10m)%z3(:,:,1) 991 ENDIF 778 992 CALL iom_put( "taum_oce", taum ) ! output wind stress module 779 993 ! … … 781 995 782 996 #if defined key_cpl_carbon_cycle 783 ! ! atmosph. CO2 (ppm) 997 ! ! ================== ! 998 ! ! atmosph. CO2 (ppm) ! 999 ! ! ================== ! 784 1000 IF( srcv(jpr_co2)%laction ) atm_co2(:,:) = frcv(jpr_co2)%z3(:,:,1) 785 1001 #endif 786 1002 1003 ! Fields received by SAS when OASIS coupling 1004 ! (arrays no more filled at sbcssm stage) 1005 ! ! ================== ! 1006 ! ! SSS ! 1007 ! ! ================== ! 1008 IF( srcv(jpr_soce)%laction ) THEN ! received by sas in case of opa <-> sas coupling 1009 sss_m(:,:) = frcv(jpr_soce)%z3(:,:,1) 1010 CALL iom_put( 'sss_m', sss_m ) 1011 ENDIF 1012 ! 1013 ! ! ================== ! 1014 ! ! SST ! 1015 ! ! ================== ! 1016 IF( srcv(jpr_toce)%laction ) THEN ! received by sas in case of opa <-> sas coupling 1017 sst_m(:,:) = frcv(jpr_toce)%z3(:,:,1) 1018 IF( srcv(jpr_soce)%laction .AND. ln_useCT ) THEN ! make sure that sst_m is the potential temperature 1019 sst_m(:,:) = eos_pt_from_ct( sst_m(:,:), sss_m(:,:) ) 1020 ENDIF 1021 ENDIF 1022 ! ! ================== ! 1023 ! ! SSH ! 1024 ! ! ================== ! 1025 IF( srcv(jpr_ssh )%laction ) THEN ! received by sas in case of opa <-> sas coupling 1026 ssh_m(:,:) = frcv(jpr_ssh )%z3(:,:,1) 1027 CALL iom_put( 'ssh_m', ssh_m ) 1028 ENDIF 1029 ! ! ================== ! 1030 ! ! surface currents ! 1031 ! ! ================== ! 1032 IF( srcv(jpr_ocx1)%laction ) THEN ! received by sas in case of opa <-> sas coupling 1033 ssu_m(:,:) = frcv(jpr_ocx1)%z3(:,:,1) 1034 ub (:,:,1) = ssu_m(:,:) ! will be used in sbcice_lim in the call of lim_sbc_tau 1035 CALL iom_put( 'ssu_m', ssu_m ) 1036 ENDIF 1037 IF( srcv(jpr_ocy1)%laction ) THEN 1038 ssv_m(:,:) = frcv(jpr_ocy1)%z3(:,:,1) 1039 vb (:,:,1) = ssv_m(:,:) ! will be used in sbcice_lim in the call of lim_sbc_tau 1040 CALL iom_put( 'ssv_m', ssv_m ) 1041 ENDIF 1042 ! ! ======================== ! 1043 ! ! first T level thickness ! 1044 ! ! ======================== ! 1045 IF( srcv(jpr_e3t1st )%laction ) THEN ! received by sas in case of opa <-> sas coupling 1046 e3t_m(:,:) = frcv(jpr_e3t1st )%z3(:,:,1) 1047 CALL iom_put( 'e3t_m', e3t_m(:,:) ) 1048 ENDIF 1049 ! ! ================================ ! 1050 ! ! fraction of solar net radiation ! 1051 ! ! ================================ ! 1052 IF( srcv(jpr_fraqsr)%laction ) THEN ! received by sas in case of opa <-> sas coupling 1053 frq_m(:,:) = frcv(jpr_fraqsr)%z3(:,:,1) 1054 CALL iom_put( 'frq_m', frq_m ) 1055 ENDIF 1056 787 1057 ! ! ========================= ! 788 IF( k_ice <= 1 ) THEN! heat & freshwater fluxes ! (Ocean only case)1058 IF( k_ice <= 1 .AND. MOD( kt-1, k_fsbc ) == 0 ) THEN ! heat & freshwater fluxes ! (Ocean only case) 789 1059 ! ! ========================= ! 790 1060 ! 791 1061 ! ! total freshwater fluxes over the ocean (emp) 792 SELECT CASE( TRIM( sn_rcv_emp%cldes ) ) ! evaporation - precipitation 793 CASE( 'conservative' ) 794 emp(:,:) = frcv(jpr_tevp)%z3(:,:,1) - ( frcv(jpr_rain)%z3(:,:,1) + frcv(jpr_snow)%z3(:,:,1) ) 795 CASE( 'oce only', 'oce and ice' ) 796 emp(:,:) = frcv(jpr_oemp)%z3(:,:,1) 797 CASE default 798 CALL ctl_stop( 'sbc_cpl_rcv: wrong definition of sn_rcv_emp%cldes' ) 799 END SELECT 1062 IF( srcv(jpr_oemp)%laction .OR. srcv(jpr_rain)%laction ) THEN 1063 SELECT CASE( TRIM( sn_rcv_emp%cldes ) ) ! evaporation - precipitation 1064 CASE( 'conservative' ) 1065 zemp(:,:) = frcv(jpr_tevp)%z3(:,:,1) - ( frcv(jpr_rain)%z3(:,:,1) + frcv(jpr_snow)%z3(:,:,1) ) 1066 CASE( 'oce only', 'oce and ice' ) 1067 zemp(:,:) = frcv(jpr_oemp)%z3(:,:,1) 1068 CASE default 1069 CALL ctl_stop( 'sbc_cpl_rcv: wrong definition of sn_rcv_emp%cldes' ) 1070 END SELECT 1071 ELSE 1072 zemp(:,:) = 0._wp 1073 ENDIF 800 1074 ! 801 1075 ! ! runoffs and calving (added in emp) 802 IF( srcv(jpr_rnf)%laction ) emp(:,:) = emp(:,:) - frcv(jpr_rnf)%z3(:,:,1) 803 IF( srcv(jpr_cal)%laction ) emp(:,:) = emp(:,:) - frcv(jpr_cal)%z3(:,:,1) 804 ! 805 !!gm : this seems to be internal cooking, not sure to need that in a generic interface 806 !!gm at least should be optional... 807 !! IF( TRIM( sn_rcv_rnf%cldes ) == 'coupled' ) THEN ! add to the total freshwater budget 808 !! ! remove negative runoff 809 !! zcumulpos = SUM( MAX( frcv(jpr_rnf)%z3(:,:,1), 0.e0 ) * e1t(:,:) * e2t(:,:) * tmask_i(:,:) ) 810 !! zcumulneg = SUM( MIN( frcv(jpr_rnf)%z3(:,:,1), 0.e0 ) * e1t(:,:) * e2t(:,:) * tmask_i(:,:) ) 811 !! IF( lk_mpp ) CALL mpp_sum( zcumulpos ) ! sum over the global domain 812 !! IF( lk_mpp ) CALL mpp_sum( zcumulneg ) 813 !! IF( zcumulpos /= 0. ) THEN ! distribute negative runoff on positive runoff grid points 814 !! zcumulneg = 1.e0 + zcumulneg / zcumulpos 815 !! frcv(jpr_rnf)%z3(:,:,1) = MAX( frcv(jpr_rnf)%z3(:,:,1), 0.e0 ) * zcumulneg 816 !! ENDIF 817 !! ! add runoff to e-p 818 !! emp(:,:) = emp(:,:) - frcv(jpr_rnf)%z3(:,:,1) 819 !! ENDIF 820 !!gm end of internal cooking 1076 IF( srcv(jpr_rnf)%laction ) rnf(:,:) = frcv(jpr_rnf)%z3(:,:,1) 1077 IF( srcv(jpr_cal)%laction ) zemp(:,:) = zemp(:,:) - frcv(jpr_cal)%z3(:,:,1) 1078 1079 IF( ln_mixcpl ) THEN ; emp(:,:) = emp(:,:) * xcplmask(:,:,0) + zemp(:,:) * zmsk(:,:) 1080 ELSE ; emp(:,:) = zemp(:,:) 1081 ENDIF 821 1082 ! 822 1083 ! ! non solar heat flux over the ocean (qns) 823 IF( srcv(jpr_qnsoce)%laction ) qns(:,:) = frcv(jpr_qnsoce)%z3(:,:,1) 824 IF( srcv(jpr_qnsmix)%laction ) qns(:,:) = frcv(jpr_qnsmix)%z3(:,:,1) 1084 IF( srcv(jpr_qnsoce)%laction ) THEN ; zqns(:,:) = frcv(jpr_qnsoce)%z3(:,:,1) 1085 ELSE IF( srcv(jpr_qnsmix)%laction ) THEN ; zqns(:,:) = frcv(jpr_qnsmix)%z3(:,:,1) 1086 ELSE ; zqns(:,:) = 0._wp 1087 END IF 825 1088 ! update qns over the free ocean with: 826 qns(:,:) = qns(:,:) - emp(:,:) * sst_m(:,:) * rcp ! remove heat content due to mass flux (assumed to be at SST) 827 IF( srcv(jpr_snow )%laction ) THEN 828 qns(:,:) = qns(:,:) - frcv(jpr_snow)%z3(:,:,1) * lfus ! energy for melting solid precipitation over the free ocean 1089 IF( nn_components /= jp_iam_opa ) THEN 1090 zqns(:,:) = zqns(:,:) - zemp(:,:) * sst_m(:,:) * rcp ! remove heat content due to mass flux (assumed to be at SST) 1091 IF( srcv(jpr_snow )%laction ) THEN 1092 zqns(:,:) = zqns(:,:) - frcv(jpr_snow)%z3(:,:,1) * lfus ! energy for melting solid precipitation over the free ocean 1093 ENDIF 1094 ENDIF 1095 IF( ln_mixcpl ) THEN ; qns(:,:) = qns(:,:) * xcplmask(:,:,0) + zqns(:,:) * zmsk(:,:) 1096 ELSE ; qns(:,:) = zqns(:,:) 829 1097 ENDIF 830 1098 831 1099 ! ! solar flux over the ocean (qsr) 832 IF( srcv(jpr_qsroce)%laction ) qsr(:,:) = frcv(jpr_qsroce)%z3(:,:,1) 833 IF( srcv(jpr_qsrmix)%laction ) qsr(:,:) = frcv(jpr_qsrmix)%z3(:,:,1) 834 IF( ln_dm2dc ) qsr(:,:) = sbc_dcy( qsr ) ! modify qsr to include the diurnal cycle 1100 IF ( srcv(jpr_qsroce)%laction ) THEN ; zqsr(:,:) = frcv(jpr_qsroce)%z3(:,:,1) 1101 ELSE IF( srcv(jpr_qsrmix)%laction ) then ; zqsr(:,:) = frcv(jpr_qsrmix)%z3(:,:,1) 1102 ELSE ; zqsr(:,:) = 0._wp 1103 ENDIF 1104 IF( ln_dm2dc .AND. ln_cpl ) zqsr(:,:) = sbc_dcy( zqsr ) ! modify qsr to include the diurnal cycle 1105 IF( ln_mixcpl ) THEN ; qsr(:,:) = qsr(:,:) * xcplmask(:,:,0) + zqsr(:,:) * zmsk(:,:) 1106 ELSE ; qsr(:,:) = zqsr(:,:) 1107 ENDIF 835 1108 ! 836 837 ENDIF 838 ! 839 CALL wrk_dealloc( jpi,jpj, ztx, zty ) 1109 ! salt flux over the ocean (received by opa in case of opa <-> sas coupling) 1110 IF( srcv(jpr_sflx )%laction ) sfx(:,:) = frcv(jpr_sflx )%z3(:,:,1) 1111 ! Ice cover (received by opa in case of opa <-> sas coupling) 1112 IF( srcv(jpr_fice )%laction ) fr_i(:,:) = frcv(jpr_fice )%z3(:,:,1) 1113 ! 1114 1115 ENDIF 1116 ! 1117 CALL wrk_dealloc( jpi,jpj, ztx, zty, zmsk, zemp, zqns, zqsr ) 840 1118 ! 841 1119 IF( nn_timing == 1 ) CALL timing_stop('sbc_cpl_rcv') … … 934 1212 ! 935 1213 ENDIF 936 937 1214 ! ! ======================= ! 938 1215 ! ! put on ice grid ! … … 1056 1333 1057 1334 1058 SUBROUTINE sbc_cpl_ice_flx( p_frld , palbi , psst , pist)1335 SUBROUTINE sbc_cpl_ice_flx( p_frld, palbi, psst, pist ) 1059 1336 !!---------------------------------------------------------------------- 1060 1337 !! *** ROUTINE sbc_cpl_ice_flx *** … … 1098 1375 REAL(wp), INTENT(in ), DIMENSION(:,:) :: p_frld ! lead fraction [0 to 1] 1099 1376 ! optional arguments, used only in 'mixed oce-ice' case 1100 REAL(wp), INTENT(in ), DIMENSION(:,:,:), OPTIONAL :: palbi ! all skies ice albedo 1101 REAL(wp), INTENT(in ), DIMENSION(:,: ), OPTIONAL :: psst ! sea surface temperature [Celsius] 1102 REAL(wp), INTENT(in ), DIMENSION(:,:,:), OPTIONAL :: pist ! ice surface temperature [Kelvin] 1103 ! 1104 INTEGER :: jl ! dummy loop index 1105 REAL(wp), POINTER, DIMENSION(:,:) :: zcptn, ztmp, zicefr 1377 REAL(wp), INTENT(in ), DIMENSION(:,:,:), OPTIONAL :: palbi ! all skies ice albedo 1378 REAL(wp), INTENT(in ), DIMENSION(:,: ), OPTIONAL :: psst ! sea surface temperature [Celsius] 1379 REAL(wp), INTENT(in ), DIMENSION(:,:,:), OPTIONAL :: pist ! ice surface temperature [Kelvin] 1380 ! 1381 INTEGER :: jl ! dummy loop index 1382 REAL(wp), POINTER, DIMENSION(:,: ) :: zcptn, ztmp, zicefr, zmsk 1383 REAL(wp), POINTER, DIMENSION(:,: ) :: zemp_tot, zemp_ice, zsprecip, ztprecip, zqns_tot, zqsr_tot 1384 REAL(wp), POINTER, DIMENSION(:,:,:) :: zqns_ice, zqsr_ice, zdqns_ice 1385 REAL(wp), POINTER, DIMENSION(:,: ) :: zevap, zsnw, zqns_oce, zqprec_ice, zqemp_oce ! for LIM3 1106 1386 !!---------------------------------------------------------------------- 1107 1387 ! 1108 1388 IF( nn_timing == 1 ) CALL timing_start('sbc_cpl_ice_flx') 1109 1389 ! 1110 CALL wrk_alloc( jpi,jpj, zcptn, ztmp, zicefr ) 1111 1390 CALL wrk_alloc( jpi,jpj, zcptn, ztmp, zicefr, zmsk, zemp_tot, zemp_ice, zsprecip, ztprecip, zqns_tot, zqsr_tot ) 1391 CALL wrk_alloc( jpi,jpj,jpl, zqns_ice, zqsr_ice, zdqns_ice ) 1392 1393 IF( ln_mixcpl ) zmsk(:,:) = 1. - xcplmask(:,:,0) 1112 1394 zicefr(:,:) = 1.- p_frld(:,:) 1113 1395 zcptn(:,:) = rcp * sst_m(:,:) … … 1117 1399 ! ! ========================= ! 1118 1400 ! 1119 ! ! total Precipitations - total Evaporation (emp_tot) 1120 ! ! solid precipitation - sublimation (emp_ice) 1121 ! ! solid Precipitation (sprecip) 1401 ! ! total Precipitation - total Evaporation (emp_tot) 1402 ! ! solid precipitation - sublimation (emp_ice) 1403 ! ! solid Precipitation (sprecip) 1404 ! ! liquid + solid Precipitation (tprecip) 1122 1405 SELECT CASE( TRIM( sn_rcv_emp%cldes ) ) 1123 1406 CASE( 'conservative' ) ! received fields: jpr_rain, jpr_snow, jpr_ievp, jpr_tevp 1124 sprecip(:,:) = frcv(jpr_snow)%z3(:,:,1)! May need to ensure positive here1125 tprecip(:,:) = frcv(jpr_rain)%z3(:,:,1) + sprecip (:,:)! May need to ensure positive here1126 emp_tot(:,:) = frcv(jpr_tevp)%z3(:,:,1) -tprecip(:,:)1127 emp_ice(:,:) = frcv(jpr_ievp)%z3(:,:,1) - frcv(jpr_snow)%z3(:,:,1)1407 zsprecip(:,:) = frcv(jpr_snow)%z3(:,:,1) ! May need to ensure positive here 1408 ztprecip(:,:) = frcv(jpr_rain)%z3(:,:,1) + zsprecip(:,:) ! May need to ensure positive here 1409 zemp_tot(:,:) = frcv(jpr_tevp)%z3(:,:,1) - ztprecip(:,:) 1410 zemp_ice(:,:) = frcv(jpr_ievp)%z3(:,:,1) - frcv(jpr_snow)%z3(:,:,1) 1128 1411 CALL iom_put( 'rain' , frcv(jpr_rain)%z3(:,:,1) ) ! liquid precipitation 1129 1412 IF( iom_use('hflx_rain_cea') ) & … … 1136 1419 CALL iom_put( 'hflx_evap_cea', ztmp(:,:) * zcptn(:,:) ) ! heat flux from from evap (cell average) 1137 1420 CASE( 'oce and ice' ) ! received fields: jpr_sbpr, jpr_semp, jpr_oemp, jpr_ievp 1138 emp_tot(:,:) = p_frld(:,:) * frcv(jpr_oemp)%z3(:,:,1) + zicefr(:,:) * frcv(jpr_sbpr)%z3(:,:,1) 1139 emp_ice(:,:) = frcv(jpr_semp)%z3(:,:,1) 1140 sprecip(:,:) = - frcv(jpr_semp)%z3(:,:,1) + frcv(jpr_ievp)%z3(:,:,1) 1421 zemp_tot(:,:) = p_frld(:,:) * frcv(jpr_oemp)%z3(:,:,1) + zicefr(:,:) * frcv(jpr_sbpr)%z3(:,:,1) 1422 zemp_ice(:,:) = frcv(jpr_semp)%z3(:,:,1) 1423 zsprecip(:,:) = frcv(jpr_ievp)%z3(:,:,1) - frcv(jpr_semp)%z3(:,:,1) 1424 ztprecip(:,:) = frcv(jpr_semp)%z3(:,:,1) - frcv(jpr_sbpr)%z3(:,:,1) + zsprecip(:,:) 1141 1425 END SELECT 1426 1427 IF( iom_use('subl_ai_cea') ) & 1428 CALL iom_put( 'subl_ai_cea', frcv(jpr_ievp)%z3(:,:,1) * zicefr(:,:) ) ! Sublimation over sea-ice (cell average) 1429 ! 1430 ! ! runoffs and calving (put in emp_tot) 1431 IF( srcv(jpr_rnf)%laction ) rnf(:,:) = frcv(jpr_rnf)%z3(:,:,1) 1432 IF( srcv(jpr_cal)%laction ) THEN 1433 zemp_tot(:,:) = zemp_tot(:,:) - frcv(jpr_cal)%z3(:,:,1) 1434 CALL iom_put( 'calving_cea', frcv(jpr_cal)%z3(:,:,1) ) 1435 ENDIF 1436 1437 IF( ln_mixcpl ) THEN 1438 emp_tot(:,:) = emp_tot(:,:) * xcplmask(:,:,0) + zemp_tot(:,:) * zmsk(:,:) 1439 emp_ice(:,:) = emp_ice(:,:) * xcplmask(:,:,0) + zemp_ice(:,:) * zmsk(:,:) 1440 sprecip(:,:) = sprecip(:,:) * xcplmask(:,:,0) + zsprecip(:,:) * zmsk(:,:) 1441 tprecip(:,:) = tprecip(:,:) * xcplmask(:,:,0) + ztprecip(:,:) * zmsk(:,:) 1442 ELSE 1443 emp_tot(:,:) = zemp_tot(:,:) 1444 emp_ice(:,:) = zemp_ice(:,:) 1445 sprecip(:,:) = zsprecip(:,:) 1446 tprecip(:,:) = ztprecip(:,:) 1447 ENDIF 1142 1448 1143 1449 CALL iom_put( 'snowpre' , sprecip ) ! Snow … … 1146 1452 IF( iom_use('snow_ai_cea') ) & 1147 1453 CALL iom_put( 'snow_ai_cea', sprecip(:,:) * zicefr(:,:) ) ! Snow over sea-ice (cell average) 1148 IF( iom_use('subl_ai_cea') ) &1149 CALL iom_put( 'subl_ai_cea', frcv(jpr_ievp)%z3(:,:,1) * zicefr(:,:) ) ! Sublimation over sea-ice (cell average)1150 !1151 ! ! runoffs and calving (put in emp_tot)1152 IF( srcv(jpr_rnf)%laction ) THEN1153 emp_tot(:,:) = emp_tot(:,:) - frcv(jpr_rnf)%z3(:,:,1)1154 CALL iom_put( 'runoffs' , frcv(jpr_rnf)%z3(:,:,1) ) ! rivers1155 IF( iom_use('hflx_rnf_cea') ) &1156 CALL iom_put( 'hflx_rnf_cea' , frcv(jpr_rnf)%z3(:,:,1) * zcptn(:,:) ) ! heat flux from rivers1157 ENDIF1158 IF( srcv(jpr_cal)%laction ) THEN1159 emp_tot(:,:) = emp_tot(:,:) - frcv(jpr_cal)%z3(:,:,1)1160 CALL iom_put( 'calving_cea', frcv(jpr_cal)%z3(:,:,1) )1161 ENDIF1162 !1163 !!gm : this seems to be internal cooking, not sure to need that in a generic interface1164 !!gm at least should be optional...1165 !! ! remove negative runoff ! sum over the global domain1166 !! zcumulpos = SUM( MAX( frcv(jpr_rnf)%z3(:,:,1), 0.e0 ) * e1t(:,:) * e2t(:,:) * tmask_i(:,:) )1167 !! zcumulneg = SUM( MIN( frcv(jpr_rnf)%z3(:,:,1), 0.e0 ) * e1t(:,:) * e2t(:,:) * tmask_i(:,:) )1168 !! IF( lk_mpp ) CALL mpp_sum( zcumulpos )1169 !! IF( lk_mpp ) CALL mpp_sum( zcumulneg )1170 !! IF( zcumulpos /= 0. ) THEN ! distribute negative runoff on positive runoff grid points1171 !! zcumulneg = 1.e0 + zcumulneg / zcumulpos1172 !! frcv(jpr_rnf)%z3(:,:,1) = MAX( frcv(jpr_rnf)%z3(:,:,1), 0.e0 ) * zcumulneg1173 !! ENDIF1174 !! emp_tot(:,:) = emp_tot(:,:) - frcv(jpr_rnf)%z3(:,:,1) ! add runoff to e-p1175 !!1176 !!gm end of internal cooking1177 1454 1178 1455 ! ! ========================= ! … … 1180 1457 ! ! ========================= ! 1181 1458 CASE( 'oce only' ) ! the required field is directly provided 1182 qns_tot(:,: ) = frcv(jpr_qnsoce)%z3(:,:,1)1459 zqns_tot(:,: ) = frcv(jpr_qnsoce)%z3(:,:,1) 1183 1460 CASE( 'conservative' ) ! the required fields are directly provided 1184 qns_tot(:,: ) = frcv(jpr_qnsmix)%z3(:,:,1)1461 zqns_tot(:,: ) = frcv(jpr_qnsmix)%z3(:,:,1) 1185 1462 IF ( TRIM(sn_rcv_qns%clcat) == 'yes' ) THEN 1186 qns_ice(:,:,1:jpl) = frcv(jpr_qnsice)%z3(:,:,1:jpl)1463 zqns_ice(:,:,1:jpl) = frcv(jpr_qnsice)%z3(:,:,1:jpl) 1187 1464 ELSE 1188 1465 ! Set all category values equal for the moment 1189 1466 DO jl=1,jpl 1190 qns_ice(:,:,jl) = frcv(jpr_qnsice)%z3(:,:,1)1467 zqns_ice(:,:,jl) = frcv(jpr_qnsice)%z3(:,:,1) 1191 1468 ENDDO 1192 1469 ENDIF 1193 1470 CASE( 'oce and ice' ) ! the total flux is computed from ocean and ice fluxes 1194 qns_tot(:,: ) = p_frld(:,:) * frcv(jpr_qnsoce)%z3(:,:,1)1471 zqns_tot(:,: ) = p_frld(:,:) * frcv(jpr_qnsoce)%z3(:,:,1) 1195 1472 IF ( TRIM(sn_rcv_qns%clcat) == 'yes' ) THEN 1196 1473 DO jl=1,jpl 1197 qns_tot(:,: ) =qns_tot(:,:) + a_i(:,:,jl) * frcv(jpr_qnsice)%z3(:,:,jl)1198 qns_ice(:,:,jl) = frcv(jpr_qnsice)%z3(:,:,jl)1474 zqns_tot(:,: ) = zqns_tot(:,:) + a_i(:,:,jl) * frcv(jpr_qnsice)%z3(:,:,jl) 1475 zqns_ice(:,:,jl) = frcv(jpr_qnsice)%z3(:,:,jl) 1199 1476 ENDDO 1200 1477 ELSE 1201 1478 qns_tot(:,: ) = qns_tot(:,:) + zicefr(:,:) * frcv(jpr_qnsice)%z3(:,:,1) 1202 1479 DO jl=1,jpl 1203 qns_ice(:,:,jl) = frcv(jpr_qnsice)%z3(:,:,1) 1480 zqns_tot(:,: ) = zqns_tot(:,:) + zicefr(:,:) * frcv(jpr_qnsice)%z3(:,:,1) 1481 zqns_ice(:,:,jl) = frcv(jpr_qnsice)%z3(:,:,1) 1204 1482 ENDDO 1205 1483 ENDIF 1206 1484 CASE( 'mixed oce-ice' ) ! the ice flux is cumputed from the total flux, the SST and ice informations 1207 1485 ! ** NEED TO SORT OUT HOW THIS SHOULD WORK IN THE MULTI-CATEGORY CASE - CURRENTLY NOT ALLOWED WHEN INTERFACE INITIALISED ** 1208 qns_tot(:,: ) = frcv(jpr_qnsmix)%z3(:,:,1)1209 qns_ice(:,:,1) = frcv(jpr_qnsmix)%z3(:,:,1) &1486 zqns_tot(:,: ) = frcv(jpr_qnsmix)%z3(:,:,1) 1487 zqns_ice(:,:,1) = frcv(jpr_qnsmix)%z3(:,:,1) & 1210 1488 & + frcv(jpr_dqnsdt)%z3(:,:,1) * ( pist(:,:,1) - ( (rt0 + psst(:,: ) ) * p_frld(:,:) & 1211 1489 & + pist(:,:,1) * zicefr(:,:) ) ) 1212 1490 END SELECT 1213 ztmp(:,:) = p_frld(:,:) * sprecip(:,:) * lfus1214 qns_tot(:,:) = qns_tot(:,:) & ! qns_tot update over free ocean with:1215 & - ztmp(:,:) & ! remove the latent heat flux of solid precip. melting1216 & - ( emp_tot(:,:) & ! remove the heat content of mass flux (assumed to be at SST)1217 & - emp_ice(:,:) * zicefr(:,:) ) * zcptn(:,:)1218 IF( iom_use('hflx_snow_cea') ) &1219 CALL iom_put( 'hflx_snow_cea', ztmp + sprecip(:,:) * zcptn(:,:) ) ! heat flux from snow (cell average)1220 1491 !!gm 1221 !! currently it is taken into account in leads budget but not in the qns_tot, and thus not in1492 !! currently it is taken into account in leads budget but not in the zqns_tot, and thus not in 1222 1493 !! the flux that enter the ocean.... 1223 1494 !! moreover 1 - it is not diagnose anywhere.... … … 1228 1499 IF( srcv(jpr_cal)%laction ) THEN ! Iceberg melting 1229 1500 ztmp(:,:) = frcv(jpr_cal)%z3(:,:,1) * lfus ! add the latent heat of iceberg melting 1230 qns_tot(:,:) =qns_tot(:,:) - ztmp(:,:)1501 zqns_tot(:,:) = zqns_tot(:,:) - ztmp(:,:) 1231 1502 IF( iom_use('hflx_cal_cea') ) & 1232 1503 CALL iom_put( 'hflx_cal_cea', ztmp + frcv(jpr_cal)%z3(:,:,1) * zcptn(:,:) ) ! heat flux from calving 1233 1504 ENDIF 1505 1506 ztmp(:,:) = p_frld(:,:) * zsprecip(:,:) * lfus 1507 IF( iom_use('hflx_snow_cea') ) CALL iom_put( 'hflx_snow_cea', ztmp + sprecip(:,:) * zcptn(:,:) ) ! heat flux from snow (cell average) 1508 1509 #if defined key_lim3 1510 CALL wrk_alloc( jpi,jpj, zevap, zsnw, zqns_oce, zqprec_ice, zqemp_oce ) 1511 1512 ! --- evaporation --- ! 1513 ! clem: evap_ice is set to 0 for LIM3 since we still do not know what to do with sublimation 1514 ! the problem is: the atm. imposes both mass evaporation and heat removed from the snow/ice 1515 ! but it is incoherent WITH the ice model 1516 DO jl=1,jpl 1517 evap_ice(:,:,jl) = 0._wp ! should be: frcv(jpr_ievp)%z3(:,:,1) 1518 ENDDO 1519 zevap(:,:) = zemp_tot(:,:) + ztprecip(:,:) ! evaporation over ocean 1520 1521 ! --- evaporation minus precipitation --- ! 1522 emp_oce(:,:) = emp_tot(:,:) - emp_ice(:,:) 1523 1524 ! --- non solar flux over ocean --- ! 1525 ! note: p_frld cannot be = 0 since we limit the ice concentration to amax 1526 zqns_oce = 0._wp 1527 WHERE( p_frld /= 0._wp ) zqns_oce(:,:) = ( zqns_tot(:,:) - SUM( a_i * zqns_ice, dim=3 ) ) / p_frld(:,:) 1528 1529 ! --- heat flux associated with emp --- ! 1530 CALL lim_thd_snwblow( p_frld, zsnw ) ! snow distribution over ice after wind blowing 1531 zqemp_oce(:,:) = - zevap(:,:) * p_frld(:,:) * zcptn(:,:) & ! evap 1532 & + ( ztprecip(:,:) - zsprecip(:,:) ) * zcptn(:,:) & ! liquid precip 1533 & + zsprecip(:,:) * ( 1._wp - zsnw ) * ( zcptn(:,:) - lfus ) ! solid precip over ocean 1534 qemp_ice(:,:) = - frcv(jpr_ievp)%z3(:,:,1) * zicefr(:,:) * zcptn(:,:) & ! ice evap 1535 & + zsprecip(:,:) * zsnw * ( zcptn(:,:) - lfus ) ! solid precip over ice 1536 1537 ! --- heat content of precip over ice in J/m3 (to be used in 1D-thermo) --- ! 1538 zqprec_ice(:,:) = rhosn * ( zcptn(:,:) - lfus ) 1539 1540 ! --- total non solar flux --- ! 1541 zqns_tot(:,:) = zqns_tot(:,:) + qemp_ice(:,:) + zqemp_oce(:,:) 1542 1543 ! --- in case both coupled/forced are active, we must mix values --- ! 1544 IF( ln_mixcpl ) THEN 1545 qns_tot(:,:) = qns_tot(:,:) * xcplmask(:,:,0) + zqns_tot(:,:)* zmsk(:,:) 1546 qns_oce(:,:) = qns_oce(:,:) * xcplmask(:,:,0) + zqns_oce(:,:)* zmsk(:,:) 1547 DO jl=1,jpl 1548 qns_ice(:,:,jl) = qns_ice(:,:,jl) * xcplmask(:,:,0) + zqns_ice(:,:,jl)* zmsk(:,:) 1549 ENDDO 1550 qprec_ice(:,:) = qprec_ice(:,:) * xcplmask(:,:,0) + zqprec_ice(:,:)* zmsk(:,:) 1551 qemp_oce (:,:) = qemp_oce(:,:) * xcplmask(:,:,0) + zqemp_oce(:,:)* zmsk(:,:) 1552 !!clem evap_ice(:,:) = evap_ice(:,:) * xcplmask(:,:,0) 1553 ELSE 1554 qns_tot (:,: ) = zqns_tot (:,: ) 1555 qns_oce (:,: ) = zqns_oce (:,: ) 1556 qns_ice (:,:,:) = zqns_ice (:,:,:) 1557 qprec_ice(:,:) = zqprec_ice(:,:) 1558 qemp_oce (:,:) = zqemp_oce (:,:) 1559 ENDIF 1560 1561 CALL wrk_dealloc( jpi,jpj, zevap, zsnw, zqns_oce, zqprec_ice, zqemp_oce ) 1562 1563 #else 1564 1565 ! clem: this formulation is certainly wrong... but better than it was... 1566 zqns_tot(:,:) = zqns_tot(:,:) & ! zqns_tot update over free ocean with: 1567 & - ztmp(:,:) & ! remove the latent heat flux of solid precip. melting 1568 & - ( zemp_tot(:,:) & ! remove the heat content of mass flux (assumed to be at SST) 1569 & - zemp_ice(:,:) * zicefr(:,:) ) * zcptn(:,:) 1570 1571 IF( ln_mixcpl ) THEN 1572 qns_tot(:,:) = qns(:,:) * p_frld(:,:) + SUM( qns_ice(:,:,:) * a_i(:,:,:), dim=3 ) ! total flux from blk 1573 qns_tot(:,:) = qns_tot(:,:) * xcplmask(:,:,0) + zqns_tot(:,:)* zmsk(:,:) 1574 DO jl=1,jpl 1575 qns_ice(:,:,jl) = qns_ice(:,:,jl) * xcplmask(:,:,0) + zqns_ice(:,:,jl)* zmsk(:,:) 1576 ENDDO 1577 ELSE 1578 qns_tot(:,: ) = zqns_tot(:,: ) 1579 qns_ice(:,:,:) = zqns_ice(:,:,:) 1580 ENDIF 1581 1582 #endif 1234 1583 1235 1584 ! ! ========================= ! … … 1237 1586 ! ! ========================= ! 1238 1587 CASE( 'oce only' ) 1239 qsr_tot(:,: ) = MAX( 0._wp , frcv(jpr_qsroce)%z3(:,:,1) )1588 zqsr_tot(:,: ) = MAX( 0._wp , frcv(jpr_qsroce)%z3(:,:,1) ) 1240 1589 CASE( 'conservative' ) 1241 qsr_tot(:,: ) = frcv(jpr_qsrmix)%z3(:,:,1)1590 zqsr_tot(:,: ) = frcv(jpr_qsrmix)%z3(:,:,1) 1242 1591 IF ( TRIM(sn_rcv_qsr%clcat) == 'yes' ) THEN 1243 qsr_ice(:,:,1:jpl) = frcv(jpr_qsrice)%z3(:,:,1:jpl)1592 zqsr_ice(:,:,1:jpl) = frcv(jpr_qsrice)%z3(:,:,1:jpl) 1244 1593 ELSE 1245 1594 ! Set all category values equal for the moment 1246 1595 DO jl=1,jpl 1247 qsr_ice(:,:,jl) = frcv(jpr_qsrice)%z3(:,:,1)1596 zqsr_ice(:,:,jl) = frcv(jpr_qsrice)%z3(:,:,1) 1248 1597 ENDDO 1249 1598 ENDIF 1250 qsr_tot(:,: ) = frcv(jpr_qsrmix)%z3(:,:,1)1251 qsr_ice(:,:,1) = frcv(jpr_qsrice)%z3(:,:,1)1599 zqsr_tot(:,: ) = frcv(jpr_qsrmix)%z3(:,:,1) 1600 zqsr_ice(:,:,1) = frcv(jpr_qsrice)%z3(:,:,1) 1252 1601 CASE( 'oce and ice' ) 1253 qsr_tot(:,: ) = p_frld(:,:) * frcv(jpr_qsroce)%z3(:,:,1)1602 zqsr_tot(:,: ) = p_frld(:,:) * frcv(jpr_qsroce)%z3(:,:,1) 1254 1603 IF ( TRIM(sn_rcv_qsr%clcat) == 'yes' ) THEN 1255 1604 DO jl=1,jpl 1256 qsr_tot(:,: ) =qsr_tot(:,:) + a_i(:,:,jl) * frcv(jpr_qsrice)%z3(:,:,jl)1257 qsr_ice(:,:,jl) = frcv(jpr_qsrice)%z3(:,:,jl)1605 zqsr_tot(:,: ) = zqsr_tot(:,:) + a_i(:,:,jl) * frcv(jpr_qsrice)%z3(:,:,jl) 1606 zqsr_ice(:,:,jl) = frcv(jpr_qsrice)%z3(:,:,jl) 1258 1607 ENDDO 1259 1608 ELSE 1260 1609 qsr_tot(:,: ) = qsr_tot(:,:) + zicefr(:,:) * frcv(jpr_qsrice)%z3(:,:,1) 1261 1610 DO jl=1,jpl 1262 qsr_ice(:,:,jl) = frcv(jpr_qsrice)%z3(:,:,1) 1611 zqsr_tot(:,: ) = zqsr_tot(:,:) + zicefr(:,:) * frcv(jpr_qsrice)%z3(:,:,1) 1612 zqsr_ice(:,:,jl) = frcv(jpr_qsrice)%z3(:,:,1) 1263 1613 ENDDO 1264 1614 ENDIF 1265 1615 CASE( 'mixed oce-ice' ) 1266 qsr_tot(:,: ) = frcv(jpr_qsrmix)%z3(:,:,1)1616 zqsr_tot(:,: ) = frcv(jpr_qsrmix)%z3(:,:,1) 1267 1617 ! ** NEED TO SORT OUT HOW THIS SHOULD WORK IN THE MULTI-CATEGORY CASE - CURRENTLY NOT ALLOWED WHEN INTERFACE INITIALISED ** 1268 1618 ! Create solar heat flux over ice using incoming solar heat flux and albedos 1269 1619 ! ( see OASIS3 user guide, 5th edition, p39 ) 1270 qsr_ice(:,:,1) = frcv(jpr_qsrmix)%z3(:,:,1) * ( 1.- palbi(:,:,1) ) &1620 zqsr_ice(:,:,1) = frcv(jpr_qsrmix)%z3(:,:,1) * ( 1.- palbi(:,:,1) ) & 1271 1621 & / ( 1.- ( albedo_oce_mix(:,: ) * p_frld(:,:) & 1272 1622 & + palbi (:,:,1) * zicefr(:,:) ) ) 1273 1623 END SELECT 1274 IF( ln_dm2dc ) THEN ! modify qsr to include the diurnal cycle1275 qsr_tot(:,: ) = sbc_dcy(qsr_tot(:,: ) )1624 IF( ln_dm2dc .AND. ln_cpl ) THEN ! modify qsr to include the diurnal cycle 1625 zqsr_tot(:,: ) = sbc_dcy( zqsr_tot(:,: ) ) 1276 1626 DO jl=1,jpl 1277 qsr_ice(:,:,jl) = sbc_dcy(qsr_ice(:,:,jl) )1627 zqsr_ice(:,:,jl) = sbc_dcy( zqsr_ice(:,:,jl) ) 1278 1628 ENDDO 1629 ENDIF 1630 1631 IF( ln_mixcpl ) THEN 1632 qsr_tot(:,:) = qsr(:,:) * p_frld(:,:) + SUM( qsr_ice(:,:,:) * a_i(:,:,:), dim=3 ) ! total flux from blk 1633 qsr_tot(:,:) = qsr_tot(:,:) * xcplmask(:,:,0) + zqsr_tot(:,:)* zmsk(:,:) 1634 DO jl=1,jpl 1635 qsr_ice(:,:,jl) = qsr_ice(:,:,jl) * xcplmask(:,:,0) + zqsr_ice(:,:,jl)* zmsk(:,:) 1636 ENDDO 1637 ELSE 1638 qsr_tot(:,: ) = zqsr_tot(:,: ) 1639 qsr_ice(:,:,:) = zqsr_ice(:,:,:) 1279 1640 ENDIF 1280 1641 … … 1284 1645 CASE ('coupled') 1285 1646 IF ( TRIM(sn_rcv_dqnsdt%clcat) == 'yes' ) THEN 1286 dqns_ice(:,:,1:jpl) = frcv(jpr_dqnsdt)%z3(:,:,1:jpl)1647 zdqns_ice(:,:,1:jpl) = frcv(jpr_dqnsdt)%z3(:,:,1:jpl) 1287 1648 ELSE 1288 1649 ! Set all category values equal for the moment 1289 1650 DO jl=1,jpl 1290 dqns_ice(:,:,jl) = frcv(jpr_dqnsdt)%z3(:,:,1)1651 zdqns_ice(:,:,jl) = frcv(jpr_dqnsdt)%z3(:,:,1) 1291 1652 ENDDO 1292 1653 ENDIF 1293 1654 END SELECT 1294 1655 1656 IF( ln_mixcpl ) THEN 1657 DO jl=1,jpl 1658 dqns_ice(:,:,jl) = dqns_ice(:,:,jl) * xcplmask(:,:,0) + zdqns_ice(:,:,jl) * zmsk(:,:) 1659 ENDDO 1660 ELSE 1661 dqns_ice(:,:,:) = zdqns_ice(:,:,:) 1662 ENDIF 1663 1295 1664 ! ! ========================= ! 1296 1665 SELECT CASE( TRIM( sn_rcv_iceflx%cldes ) ) ! topmelt and botmelt ! … … 1308 1677 fr2_i0(:,:) = ( 0.82 * ( 1.0 - cldf_ice ) + 0.65 * cldf_ice ) 1309 1678 1310 CALL wrk_dealloc( jpi,jpj, zcptn, ztmp, zicefr ) 1679 CALL wrk_dealloc( jpi,jpj, zcptn, ztmp, zicefr, zmsk, zemp_tot, zemp_ice, zsprecip, ztprecip, zqns_tot, zqsr_tot ) 1680 CALL wrk_dealloc( jpi,jpj,jpl, zqns_ice, zqsr_ice, zdqns_ice ) 1311 1681 ! 1312 1682 IF( nn_timing == 1 ) CALL timing_stop('sbc_cpl_ice_flx') … … 1328 1698 INTEGER :: ji, jj, jl ! dummy loop indices 1329 1699 INTEGER :: isec, info ! local integer 1700 REAL(wp) :: zumax, zvmax 1330 1701 REAL(wp), POINTER, DIMENSION(:,:) :: zfr_l, ztmp1, ztmp2, zotx1, zoty1, zotz1, zitx1, zity1, zitz1 1331 1702 REAL(wp), POINTER, DIMENSION(:,:,:) :: ztmp3, ztmp4 … … 1344 1715 ! ! ------------------------- ! 1345 1716 IF( ssnd(jps_toce)%laction .OR. ssnd(jps_tice)%laction .OR. ssnd(jps_tmix)%laction ) THEN 1346 SELECT CASE( sn_snd_temp%cldes) 1347 CASE( 'oce only' ) ; ztmp1(:,:) = tsn(:,:,1,jp_tem) + rt0 1348 CASE( 'weighted oce and ice' ) ; ztmp1(:,:) = ( tsn(:,:,1,jp_tem) + rt0 ) * zfr_l(:,:) 1349 SELECT CASE( sn_snd_temp%clcat ) 1350 CASE( 'yes' ) 1351 ztmp3(:,:,1:jpl) = tn_ice(:,:,1:jpl) * a_i(:,:,1:jpl) 1352 CASE( 'no' ) 1353 ztmp3(:,:,:) = 0.0 1717 1718 IF ( nn_components == jp_iam_opa ) THEN 1719 ztmp1(:,:) = tsn(:,:,1,jp_tem) ! send temperature as it is (potential or conservative) -> use of ln_useCT on the received part 1720 ELSE 1721 ! we must send the surface potential temperature 1722 IF( ln_useCT ) THEN ; ztmp1(:,:) = eos_pt_from_ct( tsn(:,:,1,jp_tem), tsn(:,:,1,jp_sal) ) 1723 ELSE ; ztmp1(:,:) = tsn(:,:,1,jp_tem) 1724 ENDIF 1725 ! 1726 SELECT CASE( sn_snd_temp%cldes) 1727 CASE( 'oce only' ) ; ztmp1(:,:) = ztmp1(:,:) + rt0 1728 CASE( 'weighted oce and ice' ) ; ztmp1(:,:) = ( ztmp1(:,:) + rt0 ) * zfr_l(:,:) 1729 SELECT CASE( sn_snd_temp%clcat ) 1730 CASE( 'yes' ) 1731 ztmp3(:,:,1:jpl) = tn_ice(:,:,1:jpl) * a_i(:,:,1:jpl) 1732 CASE( 'no' ) 1733 ztmp3(:,:,:) = 0.0 1734 DO jl=1,jpl 1735 ztmp3(:,:,1) = ztmp3(:,:,1) + tn_ice(:,:,jl) * a_i(:,:,jl) 1736 ENDDO 1737 CASE default ; CALL ctl_stop( 'sbc_cpl_snd: wrong definition of sn_snd_temp%clcat' ) 1738 END SELECT 1739 CASE( 'mixed oce-ice' ) 1740 ztmp1(:,:) = ( ztmp1(:,:) + rt0 ) * zfr_l(:,:) 1354 1741 DO jl=1,jpl 1355 ztmp 3(:,:,1) = ztmp3(:,:,1) + tn_ice(:,:,jl) * a_i(:,:,jl)1742 ztmp1(:,:) = ztmp1(:,:) + tn_ice(:,:,jl) * a_i(:,:,jl) 1356 1743 ENDDO 1357 CASE default ; CALL ctl_stop( 'sbc_cpl_snd: wrong definition of sn_snd_temp%clcat' )1744 CASE default ; CALL ctl_stop( 'sbc_cpl_snd: wrong definition of sn_snd_temp%cldes' ) 1358 1745 END SELECT 1359 CASE( 'mixed oce-ice' ) 1360 ztmp1(:,:) = ( tsn(:,:,1,jp_tem) + rt0 ) * zfr_l(:,:) 1361 DO jl=1,jpl 1362 ztmp1(:,:) = ztmp1(:,:) + tn_ice(:,:,jl) * a_i(:,:,jl) 1363 ENDDO 1364 CASE default ; CALL ctl_stop( 'sbc_cpl_snd: wrong definition of sn_snd_temp%cldes' ) 1365 END SELECT 1746 ENDIF 1366 1747 IF( ssnd(jps_toce)%laction ) CALL cpl_snd( jps_toce, isec, RESHAPE ( ztmp1, (/jpi,jpj,1/) ), info ) 1367 1748 IF( ssnd(jps_tice)%laction ) CALL cpl_snd( jps_tice, isec, ztmp3, info ) … … 1385 1766 ! ! Ice fraction & Thickness ! 1386 1767 ! ! ------------------------- ! 1387 ! Send ice fraction field 1768 ! Send ice fraction field to atmosphere 1388 1769 IF( ssnd(jps_fice)%laction ) THEN 1389 1770 SELECT CASE( sn_snd_thick%clcat ) … … 1392 1773 CASE default ; CALL ctl_stop( 'sbc_cpl_snd: wrong definition of sn_snd_thick%clcat' ) 1393 1774 END SELECT 1394 CALL cpl_snd( jps_fice, isec, ztmp3, info ) 1775 IF( ssnd(jps_fice)%laction ) CALL cpl_snd( jps_fice, isec, ztmp3, info ) 1776 ENDIF 1777 1778 ! Send ice fraction field to OPA (sent by SAS in SAS-OPA coupling) 1779 IF( ssnd(jps_fice2)%laction ) THEN 1780 ztmp3(:,:,1) = fr_i(:,:) 1781 IF( ssnd(jps_fice2)%laction ) CALL cpl_snd( jps_fice2, isec, ztmp3, info ) 1395 1782 ENDIF 1396 1783 … … 1440 1827 ! i-1 i i 1441 1828 ! i i+1 (for I) 1442 SELECT CASE( TRIM( sn_snd_crt%cldes ) ) 1443 CASE( 'oce only' ) ! C-grid ==> T 1444 DO jj = 2, jpjm1 1445 DO ji = fs_2, fs_jpim1 ! vector opt. 1446 zotx1(ji,jj) = 0.5 * ( un(ji,jj,1) + un(ji-1,jj ,1) ) 1447 zoty1(ji,jj) = 0.5 * ( vn(ji,jj,1) + vn(ji ,jj-1,1) ) 1448 END DO 1449 END DO 1450 CASE( 'weighted oce and ice' ) 1451 SELECT CASE ( cp_ice_msh ) 1452 CASE( 'C' ) ! Ocean and Ice on C-grid ==> T 1829 IF( nn_components == jp_iam_opa ) THEN 1830 zotx1(:,:) = un(:,:,1) 1831 zoty1(:,:) = vn(:,:,1) 1832 ELSE 1833 SELECT CASE( TRIM( sn_snd_crt%cldes ) ) 1834 CASE( 'oce only' ) ! C-grid ==> T 1453 1835 DO jj = 2, jpjm1 1454 1836 DO ji = fs_2, fs_jpim1 ! vector opt. 1455 zotx1(ji,jj) = 0.5 * ( un (ji,jj,1) + un (ji-1,jj ,1) ) * zfr_l(ji,jj) 1456 zoty1(ji,jj) = 0.5 * ( vn (ji,jj,1) + vn (ji ,jj-1,1) ) * zfr_l(ji,jj) 1457 zitx1(ji,jj) = 0.5 * ( u_ice(ji,jj ) + u_ice(ji-1,jj ) ) * fr_i(ji,jj) 1458 zity1(ji,jj) = 0.5 * ( v_ice(ji,jj ) + v_ice(ji ,jj-1 ) ) * fr_i(ji,jj) 1837 zotx1(ji,jj) = 0.5 * ( un(ji,jj,1) + un(ji-1,jj ,1) ) 1838 zoty1(ji,jj) = 0.5 * ( vn(ji,jj,1) + vn(ji ,jj-1,1) ) 1459 1839 END DO 1460 1840 END DO 1461 CASE( 'I' ) ! Ocean on C grid, Ice on I-point (B-grid) ==> T 1462 DO jj = 2, jpjm1 1463 DO ji = 2, jpim1 ! NO vector opt. 1464 zotx1(ji,jj) = 0.5 * ( un(ji,jj,1) + un(ji-1,jj ,1) ) * zfr_l(ji,jj) 1465 zoty1(ji,jj) = 0.5 * ( vn(ji,jj,1) + vn(ji ,jj-1,1) ) * zfr_l(ji,jj) 1466 zitx1(ji,jj) = 0.25 * ( u_ice(ji+1,jj+1) + u_ice(ji,jj+1) & 1467 & + u_ice(ji+1,jj ) + u_ice(ji,jj ) ) * fr_i(ji,jj) 1468 zity1(ji,jj) = 0.25 * ( v_ice(ji+1,jj+1) + v_ice(ji,jj+1) & 1469 & + v_ice(ji+1,jj ) + v_ice(ji,jj ) ) * fr_i(ji,jj) 1841 CASE( 'weighted oce and ice' ) 1842 SELECT CASE ( cp_ice_msh ) 1843 CASE( 'C' ) ! Ocean and Ice on C-grid ==> T 1844 DO jj = 2, jpjm1 1845 DO ji = fs_2, fs_jpim1 ! vector opt. 1846 zotx1(ji,jj) = 0.5 * ( un (ji,jj,1) + un (ji-1,jj ,1) ) * zfr_l(ji,jj) 1847 zoty1(ji,jj) = 0.5 * ( vn (ji,jj,1) + vn (ji ,jj-1,1) ) * zfr_l(ji,jj) 1848 zitx1(ji,jj) = 0.5 * ( u_ice(ji,jj ) + u_ice(ji-1,jj ) ) * fr_i(ji,jj) 1849 zity1(ji,jj) = 0.5 * ( v_ice(ji,jj ) + v_ice(ji ,jj-1 ) ) * fr_i(ji,jj) 1850 END DO 1470 1851 END DO 1471 END DO1472 CASE( 'F' ) ! Ocean on C grid, Ice on F-point (B-grid) ==> T1473 DO jj = 2, jpjm11474 DO ji = 2, jpim1 ! NO vector opt.1475 zotx1(ji,jj) = 0.5 * ( un(ji,jj,1) + un(ji-1,jj,1) ) * zfr_l(ji,jj)1476 zoty1(ji,jj) = 0.5 * ( vn(ji,jj,1) + vn(ji ,jj-1,1) ) * zfr_l(ji,jj)1477 zitx1(ji,jj) = 0.25 * ( u_ice(ji-1,jj-1) + u_ice(ji,jj-1) &1478 & + u_ice(ji-1,jj ) + u_ice(ji,jj ) ) * fr_i(ji,jj)1479 zity1(ji,jj) = 0.25 * ( v_ice(ji-1,jj-1) + v_ice(ji,jj-1) &1480 & + v_ice(ji-1,jj ) + v_ice(ji,jj ) ) * fr_i(ji,jj)1852 CASE( 'I' ) ! Ocean on C grid, Ice on I-point (B-grid) ==> T 1853 DO jj = 2, jpjm1 1854 DO ji = 2, jpim1 ! NO vector opt. 1855 zotx1(ji,jj) = 0.5 * ( un(ji,jj,1) + un(ji-1,jj ,1) ) * zfr_l(ji,jj) 1856 zoty1(ji,jj) = 0.5 * ( vn(ji,jj,1) + vn(ji ,jj-1,1) ) * zfr_l(ji,jj) 1857 zitx1(ji,jj) = 0.25 * ( u_ice(ji+1,jj+1) + u_ice(ji,jj+1) & 1858 & + u_ice(ji+1,jj ) + u_ice(ji,jj ) ) * fr_i(ji,jj) 1859 zity1(ji,jj) = 0.25 * ( v_ice(ji+1,jj+1) + v_ice(ji,jj+1) & 1860 & + v_ice(ji+1,jj ) + v_ice(ji,jj ) ) * fr_i(ji,jj) 1861 END DO 1481 1862 END DO 1482 END DO 1863 CASE( 'F' ) ! Ocean on C grid, Ice on F-point (B-grid) ==> T 1864 DO jj = 2, jpjm1 1865 DO ji = 2, jpim1 ! NO vector opt. 1866 zotx1(ji,jj) = 0.5 * ( un(ji,jj,1) + un(ji-1,jj ,1) ) * zfr_l(ji,jj) 1867 zoty1(ji,jj) = 0.5 * ( vn(ji,jj,1) + vn(ji ,jj-1,1) ) * zfr_l(ji,jj) 1868 zitx1(ji,jj) = 0.25 * ( u_ice(ji-1,jj-1) + u_ice(ji,jj-1) & 1869 & + u_ice(ji-1,jj ) + u_ice(ji,jj ) ) * fr_i(ji,jj) 1870 zity1(ji,jj) = 0.25 * ( v_ice(ji-1,jj-1) + v_ice(ji,jj-1) & 1871 & + v_ice(ji-1,jj ) + v_ice(ji,jj ) ) * fr_i(ji,jj) 1872 END DO 1873 END DO 1874 END SELECT 1875 CALL lbc_lnk( zitx1, 'T', -1. ) ; CALL lbc_lnk( zity1, 'T', -1. ) 1876 CASE( 'mixed oce-ice' ) 1877 SELECT CASE ( cp_ice_msh ) 1878 CASE( 'C' ) ! Ocean and Ice on C-grid ==> T 1879 DO jj = 2, jpjm1 1880 DO ji = fs_2, fs_jpim1 ! vector opt. 1881 zotx1(ji,jj) = 0.5 * ( un (ji,jj,1) + un (ji-1,jj ,1) ) * zfr_l(ji,jj) & 1882 & + 0.5 * ( u_ice(ji,jj ) + u_ice(ji-1,jj ) ) * fr_i(ji,jj) 1883 zoty1(ji,jj) = 0.5 * ( vn (ji,jj,1) + vn (ji ,jj-1,1) ) * zfr_l(ji,jj) & 1884 & + 0.5 * ( v_ice(ji,jj ) + v_ice(ji ,jj-1 ) ) * fr_i(ji,jj) 1885 END DO 1886 END DO 1887 CASE( 'I' ) ! Ocean on C grid, Ice on I-point (B-grid) ==> T 1888 DO jj = 2, jpjm1 1889 DO ji = 2, jpim1 ! NO vector opt. 1890 zotx1(ji,jj) = 0.5 * ( un(ji,jj,1) + un(ji-1,jj ,1) ) * zfr_l(ji,jj) & 1891 & + 0.25 * ( u_ice(ji+1,jj+1) + u_ice(ji,jj+1) & 1892 & + u_ice(ji+1,jj ) + u_ice(ji,jj ) ) * fr_i(ji,jj) 1893 zoty1(ji,jj) = 0.5 * ( vn(ji,jj,1) + vn(ji ,jj-1,1) ) * zfr_l(ji,jj) & 1894 & + 0.25 * ( v_ice(ji+1,jj+1) + v_ice(ji,jj+1) & 1895 & + v_ice(ji+1,jj ) + v_ice(ji,jj ) ) * fr_i(ji,jj) 1896 END DO 1897 END DO 1898 CASE( 'F' ) ! Ocean on C grid, Ice on F-point (B-grid) ==> T 1899 DO jj = 2, jpjm1 1900 DO ji = 2, jpim1 ! NO vector opt. 1901 zotx1(ji,jj) = 0.5 * ( un(ji,jj,1) + un(ji-1,jj ,1) ) * zfr_l(ji,jj) & 1902 & + 0.25 * ( u_ice(ji-1,jj-1) + u_ice(ji,jj-1) & 1903 & + u_ice(ji-1,jj ) + u_ice(ji,jj ) ) * fr_i(ji,jj) 1904 zoty1(ji,jj) = 0.5 * ( vn(ji,jj,1) + vn(ji ,jj-1,1) ) * zfr_l(ji,jj) & 1905 & + 0.25 * ( v_ice(ji-1,jj-1) + v_ice(ji,jj-1) & 1906 & + v_ice(ji-1,jj ) + v_ice(ji,jj ) ) * fr_i(ji,jj) 1907 END DO 1908 END DO 1909 END SELECT 1483 1910 END SELECT 1484 CALL lbc_lnk( zitx1, 'T', -1. ) ; CALL lbc_lnk( zity1, 'T', -1. ) 1485 CASE( 'mixed oce-ice' ) 1486 SELECT CASE ( cp_ice_msh ) 1487 CASE( 'C' ) ! Ocean and Ice on C-grid ==> T 1488 DO jj = 2, jpjm1 1489 DO ji = fs_2, fs_jpim1 ! vector opt. 1490 zotx1(ji,jj) = 0.5 * ( un (ji,jj,1) + un (ji-1,jj ,1) ) * zfr_l(ji,jj) & 1491 & + 0.5 * ( u_ice(ji,jj ) + u_ice(ji-1,jj ) ) * fr_i(ji,jj) 1492 zoty1(ji,jj) = 0.5 * ( vn (ji,jj,1) + vn (ji ,jj-1,1) ) * zfr_l(ji,jj) & 1493 & + 0.5 * ( v_ice(ji,jj ) + v_ice(ji ,jj-1 ) ) * fr_i(ji,jj) 1494 END DO 1495 END DO 1496 CASE( 'I' ) ! Ocean on C grid, Ice on I-point (B-grid) ==> T 1497 DO jj = 2, jpjm1 1498 DO ji = 2, jpim1 ! NO vector opt. 1499 zotx1(ji,jj) = 0.5 * ( un(ji,jj,1) + un(ji-1,jj ,1) ) * zfr_l(ji,jj) & 1500 & + 0.25 * ( u_ice(ji+1,jj+1) + u_ice(ji,jj+1) & 1501 & + u_ice(ji+1,jj ) + u_ice(ji,jj ) ) * fr_i(ji,jj) 1502 zoty1(ji,jj) = 0.5 * ( vn(ji,jj,1) + vn(ji ,jj-1,1) ) * zfr_l(ji,jj) & 1503 & + 0.25 * ( v_ice(ji+1,jj+1) + v_ice(ji,jj+1) & 1504 & + v_ice(ji+1,jj ) + v_ice(ji,jj ) ) * fr_i(ji,jj) 1505 END DO 1506 END DO 1507 CASE( 'F' ) ! Ocean on C grid, Ice on F-point (B-grid) ==> T 1508 DO jj = 2, jpjm1 1509 DO ji = 2, jpim1 ! NO vector opt. 1510 zotx1(ji,jj) = 0.5 * ( un(ji,jj,1) + un(ji-1,jj ,1) ) * zfr_l(ji,jj) & 1511 & + 0.25 * ( u_ice(ji-1,jj-1) + u_ice(ji,jj-1) & 1512 & + u_ice(ji-1,jj ) + u_ice(ji,jj ) ) * fr_i(ji,jj) 1513 zoty1(ji,jj) = 0.5 * ( vn(ji,jj,1) + vn(ji ,jj-1,1) ) * zfr_l(ji,jj) & 1514 & + 0.25 * ( v_ice(ji-1,jj-1) + v_ice(ji,jj-1) & 1515 & + v_ice(ji-1,jj ) + v_ice(ji,jj ) ) * fr_i(ji,jj) 1516 END DO 1517 END DO 1518 END SELECT 1519 END SELECT 1520 CALL lbc_lnk( zotx1, ssnd(jps_ocx1)%clgrid, -1. ) ; CALL lbc_lnk( zoty1, ssnd(jps_ocy1)%clgrid, -1. ) 1911 CALL lbc_lnk( zotx1, ssnd(jps_ocx1)%clgrid, -1. ) ; CALL lbc_lnk( zoty1, ssnd(jps_ocy1)%clgrid, -1. ) 1912 ! 1913 ENDIF 1521 1914 ! 1522 1915 ! … … 1558 1951 ENDIF 1559 1952 ! 1953 ! 1954 ! Fields sent by OPA to SAS when doing OPA<->SAS coupling 1955 ! ! SSH 1956 IF( ssnd(jps_ssh )%laction ) THEN 1957 ! ! removed inverse barometer ssh when Patm 1958 ! forcing is used (for sea-ice dynamics) 1959 IF( ln_apr_dyn ) THEN ; ztmp1(:,:) = sshb(:,:) - 0.5 * ( ssh_ib(:,:) + ssh_ibb(:,:) ) 1960 ELSE ; ztmp1(:,:) = sshn(:,:) 1961 ENDIF 1962 CALL cpl_snd( jps_ssh , isec, RESHAPE ( ztmp1 , (/jpi,jpj,1/) ), info ) 1963 1964 ENDIF 1965 ! ! SSS 1966 IF( ssnd(jps_soce )%laction ) THEN 1967 CALL cpl_snd( jps_soce , isec, RESHAPE ( tsn(:,:,1,jp_sal), (/jpi,jpj,1/) ), info ) 1968 ENDIF 1969 ! ! first T level thickness 1970 IF( ssnd(jps_e3t1st )%laction ) THEN 1971 CALL cpl_snd( jps_e3t1st, isec, RESHAPE ( fse3t_n(:,:,1) , (/jpi,jpj,1/) ), info ) 1972 ENDIF 1973 ! ! Qsr fraction 1974 IF( ssnd(jps_fraqsr)%laction ) THEN 1975 CALL cpl_snd( jps_fraqsr, isec, RESHAPE ( fraqsr_1lev(:,:) , (/jpi,jpj,1/) ), info ) 1976 ENDIF 1977 ! 1978 ! Fields sent by SAS to OPA when OASIS coupling 1979 ! ! Solar heat flux 1980 IF( ssnd(jps_qsroce)%laction ) CALL cpl_snd( jps_qsroce, isec, RESHAPE ( qsr , (/jpi,jpj,1/) ), info ) 1981 IF( ssnd(jps_qnsoce)%laction ) CALL cpl_snd( jps_qnsoce, isec, RESHAPE ( qns , (/jpi,jpj,1/) ), info ) 1982 IF( ssnd(jps_oemp )%laction ) CALL cpl_snd( jps_oemp , isec, RESHAPE ( emp , (/jpi,jpj,1/) ), info ) 1983 IF( ssnd(jps_sflx )%laction ) CALL cpl_snd( jps_sflx , isec, RESHAPE ( sfx , (/jpi,jpj,1/) ), info ) 1984 IF( ssnd(jps_otx1 )%laction ) CALL cpl_snd( jps_otx1 , isec, RESHAPE ( utau, (/jpi,jpj,1/) ), info ) 1985 IF( ssnd(jps_oty1 )%laction ) CALL cpl_snd( jps_oty1 , isec, RESHAPE ( vtau, (/jpi,jpj,1/) ), info ) 1986 IF( ssnd(jps_rnf )%laction ) CALL cpl_snd( jps_rnf , isec, RESHAPE ( rnf , (/jpi,jpj,1/) ), info ) 1987 IF( ssnd(jps_taum )%laction ) CALL cpl_snd( jps_taum , isec, RESHAPE ( taum, (/jpi,jpj,1/) ), info ) 1988 1560 1989 CALL wrk_dealloc( jpi,jpj, zfr_l, ztmp1, ztmp2, zotx1, zoty1, zotz1, zitx1, zity1, zitz1 ) 1561 1990 CALL wrk_dealloc( jpi,jpj,jpl, ztmp3, ztmp4 )
Note: See TracChangeset
for help on using the changeset viewer.