- Timestamp:
- 2015-07-09T12:14:37+02:00 (9 years ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
branches/UKMO/dev_r5107_hadgem3_cplseq/NEMOGCM/NEMO/OPA_SRC/SBC/sbccpl.F90
r5491 r5572 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 25 26 #if defined key_lim3 26 USE par_ice ! ice parameters27 27 USE ice ! ice variables 28 28 #endif … … 33 33 USE cpl_oasis3 ! OASIS3 coupling 34 34 USE geo2ocean ! 35 USE oce , ONLY : tsn, un, vn 35 USE oce , ONLY : tsn, un, vn, sshn, ub, vb, sshb, fraqsr_1lev 36 36 USE albedo ! 37 37 USE in_out_manager ! I/O manager … … 41 41 USE timing ! Timing 42 42 USE lbclnk ! ocean lateral boundary conditions (or mpp link) 43 USE eosbn2 44 USE sbcrnf , ONLY : l_rnfcpl 43 45 #if defined key_cpl_carbon_cycle 44 46 USE p4zflx, ONLY : oce_co2 … … 47 49 USE ice_domain_size, only: ncat 48 50 #endif 51 #if defined key_lim3 52 USE limthd_dh ! for CALL lim_thd_snwblow 53 #endif 54 49 55 IMPLICIT NONE 50 56 PRIVATE 51 !EM XIOS-OASIS-MCT compliance 57 52 58 PUBLIC sbc_cpl_init ! routine called by sbcmod.F90 53 59 PUBLIC sbc_cpl_rcv ! routine called by sbc_ice_lim(_2).F90 … … 90 96 INTEGER, PARAMETER :: jpr_topm = 32 ! topmeltn 91 97 INTEGER, PARAMETER :: jpr_botm = 33 ! botmeltn 92 INTEGER, PARAMETER :: jprcv = 33 ! total number of fields received 93 94 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 95 110 INTEGER, PARAMETER :: jps_toce = 2 ! ocean temperature 96 111 INTEGER, PARAMETER :: jps_tice = 3 ! ice temperature … … 107 122 INTEGER, PARAMETER :: jps_ivz1 = 14 ! 108 123 INTEGER, PARAMETER :: jps_co2 = 15 109 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 110 138 111 139 ! !!** namelist namsbc_cpl ** … … 126 154 LOGICAL :: ln_usecplmask ! use a coupling mask file to merge data received from several models 127 155 ! -> file cplmask.nc with the float variable called cplmask (jpi,jpj,nn_cplmodel) 128 129 REAL(wp), ALLOCATABLE, DIMENSION(:,:,:) :: xcplmask130 131 156 TYPE :: DYNARR 132 157 REAL(wp), POINTER, DIMENSION(:,:,:) :: z3 … … 140 165 141 166 !! Substitution 167 # include "domzgr_substitute.h90" 142 168 # include "vectopt_loop_substitute.h90" 143 169 !!---------------------------------------------------------------------- … … 162 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) 163 189 #endif 164 !ALLOCATE( xcplmask(jpi,jpj, nn_cplmodel) , STAT=ierr(3) )190 !ALLOCATE( xcplmask(jpi,jpj,0:nn_cplmodel) , STAT=ierr(3) ) 165 191 ! Hardwire only two models as nn_cplmodel has not been read in 166 192 ! from the namelist yet. 167 ALLOCATE( xcplmask(jpi,jpj, 2) , STAT=ierr(3) )193 ALLOCATE( xcplmask(jpi,jpj,0:2) , STAT=ierr(3) ) 168 194 ! 169 195 sbc_cpl_alloc = MAXVAL( ierr ) … … 186 212 !! * initialise the OASIS coupler 187 213 !!---------------------------------------------------------------------- 188 INTEGER, INTENT(in) :: k_ice ! ice management in the sbc (=0/1/2/3)214 INTEGER, INTENT(in) :: k_ice ! ice management in the sbc (=0/1/2/3) 189 215 !! 190 216 INTEGER :: jn ! dummy loop index … … 220 246 WRITE(numout,*)'sbc_cpl_init : namsbc_cpl namelist ' 221 247 WRITE(numout,*)'~~~~~~~~~~~~' 248 ENDIF 249 IF( lwp .AND. ln_cpl ) THEN ! control print 222 250 WRITE(numout,*)' received fields (mutiple ice categogies)' 223 251 WRITE(numout,*)' 10m wind module = ', TRIM(sn_rcv_w10m%cldes ), ' (', TRIM(sn_rcv_w10m%clcat ), ')' … … 363 391 srcv(jpr_oemp)%clname = 'OOEvaMPr' ! ocean water budget = ocean Evap - ocean precip 364 392 SELECT CASE( TRIM( sn_rcv_emp%cldes ) ) 393 CASE( 'none' ) ! nothing to do 365 394 CASE( 'oce only' ) ; srcv( jpr_oemp )%laction = .TRUE. 366 395 CASE( 'conservative' ) … … 374 403 ! ! Runoffs & Calving ! 375 404 ! ! ------------------------- ! 376 srcv(jpr_rnf )%clname = 'O_Runoff' ; IF( TRIM( sn_rcv_rnf%cldes ) == 'coupled' ) srcv(jpr_rnf)%laction = .TRUE. 377 ! This isn't right - really just want ln_rnf_emp changed 378 ! IF( TRIM( sn_rcv_rnf%cldes ) == 'climato' ) THEN ; ln_rnf = .TRUE. 379 ! ELSE ; ln_rnf = .FALSE. 380 ! ENDIF 405 srcv(jpr_rnf )%clname = 'O_Runoff' 406 IF( TRIM( sn_rcv_rnf%cldes ) == 'coupled' ) THEN 407 srcv(jpr_rnf)%laction = .TRUE. 408 l_rnfcpl = .TRUE. ! -> no need to read runoffs in sbcrnf 409 ln_rnf = nn_components /= jp_iam_sas ! -> force to go through sbcrnf if not sas 410 IF(lwp) WRITE(numout,*) 411 IF(lwp) WRITE(numout,*) ' runoffs received from oasis -> force ln_rnf = ', ln_rnf 412 ENDIF 413 ! 381 414 srcv(jpr_cal )%clname = 'OCalving' ; IF( TRIM( sn_rcv_cal%cldes ) == 'coupled' ) srcv(jpr_cal)%laction = .TRUE. 382 415 … … 388 421 srcv(jpr_qnsmix)%clname = 'O_QnsMix' 389 422 SELECT CASE( TRIM( sn_rcv_qns%cldes ) ) 423 CASE( 'none' ) ! nothing to do 390 424 CASE( 'oce only' ) ; srcv( jpr_qnsoce )%laction = .TRUE. 391 425 CASE( 'conservative' ) ; srcv( (/jpr_qnsice, jpr_qnsmix/) )%laction = .TRUE. … … 403 437 srcv(jpr_qsrmix)%clname = 'O_QsrMix' 404 438 SELECT CASE( TRIM( sn_rcv_qsr%cldes ) ) 439 CASE( 'none' ) ! nothing to do 405 440 CASE( 'oce only' ) ; srcv( jpr_qsroce )%laction = .TRUE. 406 441 CASE( 'conservative' ) ; srcv( (/jpr_qsrice, jpr_qsrmix/) )%laction = .TRUE. … … 418 453 ! 419 454 ! non solar sensitivity mandatory for LIM ice model 420 IF( TRIM( sn_rcv_dqnsdt%cldes ) == 'none' .AND. k_ice /= 0 .AND. k_ice /= 4 ) &455 IF( TRIM( sn_rcv_dqnsdt%cldes ) == 'none' .AND. k_ice /= 0 .AND. k_ice /= 4 .AND. nn_components /= jp_iam_sas ) & 421 456 CALL ctl_stop( 'sbc_cpl_init: sn_rcv_dqnsdt%cldes must be coupled in namsbc_cpl namelist' ) 422 457 ! non solar sensitivity mandatory for mixed oce-ice solar radiation coupling technique … … 451 486 srcv(jpr_topm:jpr_botm)%laction = .TRUE. 452 487 ENDIF 453 454 ! Allocate all parts of frcv used for received fields 488 ! ! ------------------------------- ! 489 ! ! OPA-SAS coupling - rcv by opa ! 490 ! ! ------------------------------- ! 491 srcv(jpr_sflx)%clname = 'O_SFLX' 492 srcv(jpr_fice)%clname = 'RIceFrc' 493 ! 494 IF( nn_components == jp_iam_opa ) THEN ! OPA coupled to SAS via OASIS: force received field by OPA (sent by SAS) 495 srcv(:)%laction = .FALSE. ! force default definition in case of opa <-> sas coupling 496 srcv(:)%clgrid = 'T' ! force default definition in case of opa <-> sas coupling 497 srcv(:)%nsgn = 1. ! force default definition in case of opa <-> sas coupling 498 srcv( (/jpr_qsroce, jpr_qnsoce, jpr_oemp, jpr_sflx, jpr_fice, jpr_otx1, jpr_oty1, jpr_taum/) )%laction = .TRUE. 499 srcv(jpr_otx1)%clgrid = 'U' ! oce components given at U-point 500 srcv(jpr_oty1)%clgrid = 'V' ! and V-point 501 ! Vectors: change of sign at north fold ONLY if on the local grid 502 srcv( (/jpr_otx1,jpr_oty1/) )%nsgn = -1. 503 sn_rcv_tau%clvgrd = 'U,V' 504 sn_rcv_tau%clvor = 'local grid' 505 sn_rcv_tau%clvref = 'spherical' 506 sn_rcv_emp%cldes = 'oce only' 507 ! 508 IF(lwp) THEN ! control print 509 WRITE(numout,*) 510 WRITE(numout,*)' Special conditions for SAS-OPA coupling ' 511 WRITE(numout,*)' OPA component ' 512 WRITE(numout,*) 513 WRITE(numout,*)' received fields from SAS component ' 514 WRITE(numout,*)' ice cover ' 515 WRITE(numout,*)' oce only EMP ' 516 WRITE(numout,*)' salt flux ' 517 WRITE(numout,*)' mixed oce-ice solar flux ' 518 WRITE(numout,*)' mixed oce-ice non solar flux ' 519 WRITE(numout,*)' wind stress U,V on local grid and sperical coordinates ' 520 WRITE(numout,*)' wind stress module' 521 WRITE(numout,*) 522 ENDIF 523 ENDIF 524 ! ! -------------------------------- ! 525 ! ! OPA-SAS coupling - rcv by sas ! 526 ! ! -------------------------------- ! 527 srcv(jpr_toce )%clname = 'I_SSTSST' 528 srcv(jpr_soce )%clname = 'I_SSSal' 529 srcv(jpr_ocx1 )%clname = 'I_OCurx1' 530 srcv(jpr_ocy1 )%clname = 'I_OCury1' 531 srcv(jpr_ssh )%clname = 'I_SSHght' 532 srcv(jpr_e3t1st)%clname = 'I_E3T1st' 533 srcv(jpr_fraqsr)%clname = 'I_FraQsr' 534 ! 535 IF( nn_components == jp_iam_sas ) THEN 536 IF( .NOT. ln_cpl ) srcv(:)%laction = .FALSE. ! force default definition in case of opa <-> sas coupling 537 IF( .NOT. ln_cpl ) srcv(:)%clgrid = 'T' ! force default definition in case of opa <-> sas coupling 538 IF( .NOT. ln_cpl ) srcv(:)%nsgn = 1. ! force default definition in case of opa <-> sas coupling 539 srcv( (/jpr_toce, jpr_soce, jpr_ssh, jpr_fraqsr, jpr_ocx1, jpr_ocy1/) )%laction = .TRUE. 540 srcv( jpr_e3t1st )%laction = lk_vvl 541 srcv(jpr_ocx1)%clgrid = 'U' ! oce components given at U-point 542 srcv(jpr_ocy1)%clgrid = 'V' ! and V-point 543 ! Vectors: change of sign at north fold ONLY if on the local grid 544 srcv(jpr_ocx1:jpr_ocy1)%nsgn = -1. 545 ! Change first letter to couple with atmosphere if already coupled OPA 546 ! this is nedeed as each variable name used in the namcouple must be unique: 547 ! for example O_Runoff received by OPA from SAS and therefore O_Runoff received by SAS from the Atmosphere 548 DO jn = 1, jprcv 549 IF ( srcv(jn)%clname(1:1) == "O" ) srcv(jn)%clname = "S"//srcv(jn)%clname(2:LEN(srcv(jn)%clname)) 550 END DO 551 ! 552 IF(lwp) THEN ! control print 553 WRITE(numout,*) 554 WRITE(numout,*)' Special conditions for SAS-OPA coupling ' 555 WRITE(numout,*)' SAS component ' 556 WRITE(numout,*) 557 IF( .NOT. ln_cpl ) THEN 558 WRITE(numout,*)' received fields from OPA component ' 559 ELSE 560 WRITE(numout,*)' Additional received fields from OPA component : ' 561 ENDIF 562 WRITE(numout,*)' sea surface temperature (Celcius) ' 563 WRITE(numout,*)' sea surface salinity ' 564 WRITE(numout,*)' surface currents ' 565 WRITE(numout,*)' sea surface height ' 566 WRITE(numout,*)' thickness of first ocean T level ' 567 WRITE(numout,*)' fraction of solar net radiation absorbed in the first ocean level' 568 WRITE(numout,*) 569 ENDIF 570 ENDIF 571 572 ! =================================================== ! 573 ! Allocate all parts of frcv used for received fields ! 574 ! =================================================== ! 455 575 DO jn = 1, jprcv 456 576 IF ( srcv(jn)%laction ) ALLOCATE( frcv(jn)%z3(jpi,jpj,srcv(jn)%nct) ) … … 458 578 ! Allocate taum part of frcv which is used even when not received as coupling field 459 579 IF ( .NOT. srcv(jpr_taum)%laction ) ALLOCATE( frcv(jpr_taum)%z3(jpi,jpj,srcv(jpr_taum)%nct) ) 580 ! Allocate w10m part of frcv which is used even when not received as coupling field 581 IF ( .NOT. srcv(jpr_w10m)%laction ) ALLOCATE( frcv(jpr_w10m)%z3(jpi,jpj,srcv(jpr_w10m)%nct) ) 582 ! Allocate jpr_otx1 part of frcv which is used even when not received as coupling field 583 IF ( .NOT. srcv(jpr_otx1)%laction ) ALLOCATE( frcv(jpr_otx1)%z3(jpi,jpj,srcv(jpr_otx1)%nct) ) 584 IF ( .NOT. srcv(jpr_oty1)%laction ) ALLOCATE( frcv(jpr_oty1)%z3(jpi,jpj,srcv(jpr_oty1)%nct) ) 460 585 ! Allocate itx1 and ity1 as they are used in sbc_cpl_ice_tau even if srcv(jpr_itx1)%laction = .FALSE. 461 586 IF( k_ice /= 0 ) THEN … … 481 606 ssnd(jps_tmix)%clname = 'O_TepMix' 482 607 SELECT CASE( TRIM( sn_snd_temp%cldes ) ) 483 CASE( 'none' ) ! nothing to do484 CASE( 'oce only' ) ; ssnd( jps_toce)%laction = .TRUE.485 CASE( ' weighted oce and ice' )608 CASE( 'none' ) ! nothing to do 609 CASE( 'oce only' ) ; ssnd( jps_toce )%laction = .TRUE. 610 CASE( 'oce and ice' , 'weighted oce and ice' ) 486 611 ssnd( (/jps_toce, jps_tice/) )%laction = .TRUE. 487 612 IF ( TRIM( sn_snd_temp%clcat ) == 'yes' ) ssnd(jps_tice)%nct = jpl 488 CASE( 'mixed oce-ice' ) ; ssnd( jps_tmix)%laction = .TRUE.613 CASE( 'mixed oce-ice' ) ; ssnd( jps_tmix )%laction = .TRUE. 489 614 CASE default ; CALL ctl_stop( 'sbc_cpl_init: wrong definition of sn_snd_temp%cldes' ) 490 615 END SELECT 491 616 492 617 ! ! ------------------------- ! 493 618 ! ! Albedo ! … … 496 621 ssnd(jps_albmix)%clname = 'O_AlbMix' 497 622 SELECT CASE( TRIM( sn_snd_alb%cldes ) ) 498 CASE( 'none' )! nothing to do499 CASE( ' weighted ice' ) ;ssnd(jps_albice)%laction = .TRUE.500 CASE( 'mixed oce-ice' ) ;ssnd(jps_albmix)%laction = .TRUE.623 CASE( 'none' ) ! nothing to do 624 CASE( 'ice' , 'weighted ice' ) ; ssnd(jps_albice)%laction = .TRUE. 625 CASE( 'mixed oce-ice' ) ; ssnd(jps_albmix)%laction = .TRUE. 501 626 CASE default ; CALL ctl_stop( 'sbc_cpl_init: wrong definition of sn_snd_alb%cldes' ) 502 627 END SELECT … … 522 647 IF ( TRIM( sn_snd_thick%clcat ) == 'yes' ) ssnd(jps_fice)%nct = jpl 523 648 ENDIF 524 649 525 650 SELECT CASE ( TRIM( sn_snd_thick%cldes ) ) 526 651 CASE( 'none' ) ! nothing to do … … 529 654 IF ( TRIM( sn_snd_thick%clcat ) == 'yes' ) THEN 530 655 ssnd(jps_hice:jps_hsnw)%nct = jpl 531 ELSE532 IF ( jpl > 1 ) THEN533 CALL ctl_stop( 'sbc_cpl_init: use weighted ice and snow option for sn_snd_thick%cldes if not exchanging category fields' )534 ENDIF535 656 ENDIF 536 657 CASE ( 'weighted ice and snow' ) … … 571 692 ! ! ------------------------- ! 572 693 ssnd(jps_co2)%clname = 'O_CO2FLX' ; IF( TRIM(sn_snd_co2%cldes) == 'coupled' ) ssnd(jps_co2 )%laction = .TRUE. 694 695 ! ! ------------------------------- ! 696 ! ! OPA-SAS coupling - snd by opa ! 697 ! ! ------------------------------- ! 698 ssnd(jps_ssh )%clname = 'O_SSHght' 699 ssnd(jps_soce )%clname = 'O_SSSal' 700 ssnd(jps_e3t1st)%clname = 'O_E3T1st' 701 ssnd(jps_fraqsr)%clname = 'O_FraQsr' 702 ! 703 IF( nn_components == jp_iam_opa ) THEN 704 ssnd(:)%laction = .FALSE. ! force default definition in case of opa <-> sas coupling 705 ssnd( (/jps_toce, jps_soce, jps_ssh, jps_fraqsr, jps_ocx1, jps_ocy1/) )%laction = .TRUE. 706 ssnd( jps_e3t1st )%laction = lk_vvl 707 ! vector definition: not used but cleaner... 708 ssnd(jps_ocx1)%clgrid = 'U' ! oce components given at U-point 709 ssnd(jps_ocy1)%clgrid = 'V' ! and V-point 710 sn_snd_crt%clvgrd = 'U,V' 711 sn_snd_crt%clvor = 'local grid' 712 sn_snd_crt%clvref = 'spherical' 713 ! 714 IF(lwp) THEN ! control print 715 WRITE(numout,*) 716 WRITE(numout,*)' sent fields to SAS component ' 717 WRITE(numout,*)' sea surface temperature (T before, Celcius) ' 718 WRITE(numout,*)' sea surface salinity ' 719 WRITE(numout,*)' surface currents U,V on local grid and spherical coordinates' 720 WRITE(numout,*)' sea surface height ' 721 WRITE(numout,*)' thickness of first ocean T level ' 722 WRITE(numout,*)' fraction of solar net radiation absorbed in the first ocean level' 723 WRITE(numout,*) 724 ENDIF 725 ENDIF 726 ! ! ------------------------------- ! 727 ! ! OPA-SAS coupling - snd by sas ! 728 ! ! ------------------------------- ! 729 ssnd(jps_sflx )%clname = 'I_SFLX' 730 ssnd(jps_fice2 )%clname = 'IIceFrc' 731 ssnd(jps_qsroce)%clname = 'I_QsrOce' 732 ssnd(jps_qnsoce)%clname = 'I_QnsOce' 733 ssnd(jps_oemp )%clname = 'IOEvaMPr' 734 ssnd(jps_otx1 )%clname = 'I_OTaux1' 735 ssnd(jps_oty1 )%clname = 'I_OTauy1' 736 ssnd(jps_rnf )%clname = 'I_Runoff' 737 ssnd(jps_taum )%clname = 'I_TauMod' 738 ! 739 IF( nn_components == jp_iam_sas ) THEN 740 IF( .NOT. ln_cpl ) ssnd(:)%laction = .FALSE. ! force default definition in case of opa <-> sas coupling 741 ssnd( (/jps_qsroce, jps_qnsoce, jps_oemp, jps_fice2, jps_sflx, jps_otx1, jps_oty1, jps_taum/) )%laction = .TRUE. 742 ! 743 ! Change first letter to couple with atmosphere if already coupled with sea_ice 744 ! this is nedeed as each variable name used in the namcouple must be unique: 745 ! for example O_SSTSST sent by OPA to SAS and therefore S_SSTSST sent by SAS to the Atmosphere 746 DO jn = 1, jpsnd 747 IF ( ssnd(jn)%clname(1:1) == "O" ) ssnd(jn)%clname = "S"//ssnd(jn)%clname(2:LEN(ssnd(jn)%clname)) 748 END DO 749 ! 750 IF(lwp) THEN ! control print 751 WRITE(numout,*) 752 IF( .NOT. ln_cpl ) THEN 753 WRITE(numout,*)' sent fields to OPA component ' 754 ELSE 755 WRITE(numout,*)' Additional sent fields to OPA component : ' 756 ENDIF 757 WRITE(numout,*)' ice cover ' 758 WRITE(numout,*)' oce only EMP ' 759 WRITE(numout,*)' salt flux ' 760 WRITE(numout,*)' mixed oce-ice solar flux ' 761 WRITE(numout,*)' mixed oce-ice non solar flux ' 762 WRITE(numout,*)' wind stress U,V components' 763 WRITE(numout,*)' wind stress module' 764 ENDIF 765 ENDIF 766 573 767 ! 574 768 ! ================================ ! … … 576 770 ! ================================ ! 577 771 578 CALL cpl_define(jprcv, jpsnd,nn_cplmodel) 772 CALL cpl_define(jprcv, jpsnd, nn_cplmodel) 773 579 774 IF (ln_usecplmask) THEN 580 775 xcplmask(:,:,:) = 0. … … 586 781 xcplmask(:,:,:) = 1. 587 782 ENDIF 588 ! 589 IF( ln_dm2dc .AND. ( cpl_freq( jpr_qsroce ) + cpl_freq( jpr_qsrmix ) /= 86400 ) ) & 783 xcplmask(:,:,0) = 1. - SUM( xcplmask(:,:,1:nn_cplmodel), dim = 3 ) 784 ! 785 ncpl_qsr_freq = cpl_freq( 'O_QsrOce' ) + cpl_freq( 'O_QsrMix' ) + cpl_freq( 'I_QsrOce' ) + cpl_freq( 'I_QsrMix' ) 786 IF( ln_dm2dc .AND. ln_cpl .AND. ncpl_qsr_freq /= 86400 ) & 590 787 & CALL ctl_stop( 'sbc_cpl_init: diurnal cycle reconstruction (ln_dm2dc) needs daily couping for solar radiation' ) 788 ncpl_qsr_freq = 86400 / ncpl_qsr_freq 591 789 592 790 CALL wrk_dealloc( jpi,jpj, zacs, zaos ) … … 642 840 !! emp upward mass flux [evap. - precip. (- runoffs) (- calving)] (ocean only case) 643 841 !!---------------------------------------------------------------------- 644 INTEGER, INTENT(in) :: kt ! ocean model time step index 645 INTEGER, INTENT(in) :: k_fsbc ! frequency of sbc (-> ice model) computation 646 INTEGER, INTENT(in) :: k_ice ! ice management in the sbc (=0/1/2/3) 647 !! 648 LOGICAL :: llnewtx, llnewtau ! update wind stress components and module?? 842 INTEGER, INTENT(in) :: kt ! ocean model time step index 843 INTEGER, INTENT(in) :: k_fsbc ! frequency of sbc (-> ice model) computation 844 INTEGER, INTENT(in) :: k_ice ! ice management in the sbc (=0/1/2/3) 845 846 !! 847 LOGICAL :: llnewtx, llnewtau ! update wind stress components and module?? 649 848 INTEGER :: ji, jj, jn ! dummy loop indices 650 849 INTEGER :: isec ! number of seconds since nit000 (assuming rdttra did not change since nit000) … … 654 853 REAL(wp) :: zcdrag = 1.5e-3 ! drag coefficient 655 854 REAL(wp) :: zzx, zzy ! temporary variables 656 REAL(wp), POINTER, DIMENSION(:,:) :: ztx, zty 855 REAL(wp), POINTER, DIMENSION(:,:) :: ztx, zty, zmsk, zemp, zqns, zqsr 657 856 !!---------------------------------------------------------------------- 658 857 ! 659 858 IF( nn_timing == 1 ) CALL timing_start('sbc_cpl_rcv') 660 859 ! 661 CALL wrk_alloc( jpi,jpj, ztx, zty ) 662 ! ! Receive all the atmos. fields (including ice information) 663 isec = ( kt - nit000 ) * NINT( rdttra(1) ) ! date of exchanges 664 DO jn = 1, jprcv ! received fields sent by the atmosphere 665 IF( srcv(jn)%laction ) CALL cpl_rcv( jn, isec, frcv(jn)%z3, xcplmask, nrcvinfo(jn) ) 860 CALL wrk_alloc( jpi,jpj, ztx, zty, zmsk, zemp, zqns, zqsr ) 861 ! 862 IF( ln_mixcpl ) zmsk(:,:) = 1. - xcplmask(:,:,0) 863 ! 864 ! ! ======================================================= ! 865 ! ! Receive all the atmos. fields (including ice information) 866 ! ! ======================================================= ! 867 isec = ( kt - nit000 ) * NINT( rdttra(1) ) ! date of exchanges 868 DO jn = 1, jprcv ! received fields sent by the atmosphere 869 IF( srcv(jn)%laction ) CALL cpl_rcv( jn, isec, frcv(jn)%z3, xcplmask(:,:,1:nn_cplmodel), nrcvinfo(jn) ) 666 870 END DO 667 871 … … 723 927 ! 724 928 ENDIF 725 726 929 ! ! ========================= ! 727 930 ! ! wind stress module ! (taum) … … 752 955 ENDIF 753 956 ENDIF 754 957 ! 755 958 ! ! ========================= ! 756 959 ! ! 10 m wind speed ! (wndm) … … 765 968 !CDIR NOVERRCHK 766 969 DO ji = 1, jpi 767 wndm(ji,jj) = SQRT( frcv(jpr_taum)%z3(ji,jj,1) * zcoef )970 frcv(jpr_w10m)%z3(ji,jj,1) = SQRT( frcv(jpr_taum)%z3(ji,jj,1) * zcoef ) 768 971 END DO 769 972 END DO 770 973 ENDIF 771 ELSE772 IF ( nrcvinfo(jpr_w10m) == OASIS_Rcv ) wndm(:,:) = frcv(jpr_w10m)%z3(:,:,1)773 974 ENDIF 774 975 … … 777 978 IF( MOD( kt-1, k_fsbc ) == 0 ) THEN 778 979 ! 779 utau(:,:) = frcv(jpr_otx1)%z3(:,:,1) 780 vtau(:,:) = frcv(jpr_oty1)%z3(:,:,1) 781 taum(:,:) = frcv(jpr_taum)%z3(:,:,1) 980 IF( ln_mixcpl ) THEN 981 utau(:,:) = utau(:,:) * xcplmask(:,:,0) + frcv(jpr_otx1)%z3(:,:,1) * zmsk(:,:) 982 vtau(:,:) = vtau(:,:) * xcplmask(:,:,0) + frcv(jpr_oty1)%z3(:,:,1) * zmsk(:,:) 983 taum(:,:) = taum(:,:) * xcplmask(:,:,0) + frcv(jpr_taum)%z3(:,:,1) * zmsk(:,:) 984 wndm(:,:) = wndm(:,:) * xcplmask(:,:,0) + frcv(jpr_w10m)%z3(:,:,1) * zmsk(:,:) 985 ELSE 986 utau(:,:) = frcv(jpr_otx1)%z3(:,:,1) 987 vtau(:,:) = frcv(jpr_oty1)%z3(:,:,1) 988 taum(:,:) = frcv(jpr_taum)%z3(:,:,1) 989 wndm(:,:) = frcv(jpr_w10m)%z3(:,:,1) 990 ENDIF 782 991 CALL iom_put( "taum_oce", taum ) ! output wind stress module 783 992 ! … … 785 994 786 995 #if defined key_cpl_carbon_cycle 787 ! ! atmosph. CO2 (ppm) 996 ! ! ================== ! 997 ! ! atmosph. CO2 (ppm) ! 998 ! ! ================== ! 788 999 IF( srcv(jpr_co2)%laction ) atm_co2(:,:) = frcv(jpr_co2)%z3(:,:,1) 789 1000 #endif 790 1001 1002 ! Fields received by SAS when OASIS coupling 1003 ! (arrays no more filled at sbcssm stage) 1004 ! ! ================== ! 1005 ! ! SSS ! 1006 ! ! ================== ! 1007 IF( srcv(jpr_soce)%laction ) THEN ! received by sas in case of opa <-> sas coupling 1008 sss_m(:,:) = frcv(jpr_soce)%z3(:,:,1) 1009 CALL iom_put( 'sss_m', sss_m ) 1010 ENDIF 1011 ! 1012 ! ! ================== ! 1013 ! ! SST ! 1014 ! ! ================== ! 1015 IF( srcv(jpr_toce)%laction ) THEN ! received by sas in case of opa <-> sas coupling 1016 sst_m(:,:) = frcv(jpr_toce)%z3(:,:,1) 1017 IF( srcv(jpr_soce)%laction .AND. ln_useCT ) THEN ! make sure that sst_m is the potential temperature 1018 sst_m(:,:) = eos_pt_from_ct( sst_m(:,:), sss_m(:,:) ) 1019 ENDIF 1020 ENDIF 1021 ! ! ================== ! 1022 ! ! SSH ! 1023 ! ! ================== ! 1024 IF( srcv(jpr_ssh )%laction ) THEN ! received by sas in case of opa <-> sas coupling 1025 ssh_m(:,:) = frcv(jpr_ssh )%z3(:,:,1) 1026 CALL iom_put( 'ssh_m', ssh_m ) 1027 ENDIF 1028 ! ! ================== ! 1029 ! ! surface currents ! 1030 ! ! ================== ! 1031 IF( srcv(jpr_ocx1)%laction ) THEN ! received by sas in case of opa <-> sas coupling 1032 ssu_m(:,:) = frcv(jpr_ocx1)%z3(:,:,1) 1033 ub (:,:,1) = ssu_m(:,:) ! will be used in sbcice_lim in the call of lim_sbc_tau 1034 CALL iom_put( 'ssu_m', ssu_m ) 1035 ENDIF 1036 IF( srcv(jpr_ocy1)%laction ) THEN 1037 ssv_m(:,:) = frcv(jpr_ocy1)%z3(:,:,1) 1038 vb (:,:,1) = ssv_m(:,:) ! will be used in sbcice_lim in the call of lim_sbc_tau 1039 CALL iom_put( 'ssv_m', ssv_m ) 1040 ENDIF 1041 ! ! ======================== ! 1042 ! ! first T level thickness ! 1043 ! ! ======================== ! 1044 IF( srcv(jpr_e3t1st )%laction ) THEN ! received by sas in case of opa <-> sas coupling 1045 e3t_m(:,:) = frcv(jpr_e3t1st )%z3(:,:,1) 1046 CALL iom_put( 'e3t_m', e3t_m(:,:) ) 1047 ENDIF 1048 ! ! ================================ ! 1049 ! ! fraction of solar net radiation ! 1050 ! ! ================================ ! 1051 IF( srcv(jpr_fraqsr)%laction ) THEN ! received by sas in case of opa <-> sas coupling 1052 frq_m(:,:) = frcv(jpr_fraqsr)%z3(:,:,1) 1053 CALL iom_put( 'frq_m', frq_m ) 1054 ENDIF 1055 791 1056 ! ! ========================= ! 792 IF( k_ice <= 1 ) THEN! heat & freshwater fluxes ! (Ocean only case)1057 IF( k_ice <= 1 .AND. MOD( kt-1, k_fsbc ) == 0 ) THEN ! heat & freshwater fluxes ! (Ocean only case) 793 1058 ! ! ========================= ! 794 1059 ! 795 1060 ! ! total freshwater fluxes over the ocean (emp) 796 SELECT CASE( TRIM( sn_rcv_emp%cldes ) ) ! evaporation - precipitation 797 CASE( 'conservative' ) 798 emp(:,:) = frcv(jpr_tevp)%z3(:,:,1) - ( frcv(jpr_rain)%z3(:,:,1) + frcv(jpr_snow)%z3(:,:,1) ) 799 CASE( 'oce only', 'oce and ice' ) 800 emp(:,:) = frcv(jpr_oemp)%z3(:,:,1) 801 CASE default 802 CALL ctl_stop( 'sbc_cpl_rcv: wrong definition of sn_rcv_emp%cldes' ) 803 END SELECT 1061 IF( srcv(jpr_oemp)%laction .OR. srcv(jpr_rain)%laction ) THEN 1062 SELECT CASE( TRIM( sn_rcv_emp%cldes ) ) ! evaporation - precipitation 1063 CASE( 'conservative' ) 1064 zemp(:,:) = frcv(jpr_tevp)%z3(:,:,1) - ( frcv(jpr_rain)%z3(:,:,1) + frcv(jpr_snow)%z3(:,:,1) ) 1065 CASE( 'oce only', 'oce and ice' ) 1066 zemp(:,:) = frcv(jpr_oemp)%z3(:,:,1) 1067 CASE default 1068 CALL ctl_stop( 'sbc_cpl_rcv: wrong definition of sn_rcv_emp%cldes' ) 1069 END SELECT 1070 ELSE 1071 zemp(:,:) = 0._wp 1072 ENDIF 804 1073 ! 805 1074 ! ! runoffs and calving (added in emp) 806 IF( srcv(jpr_rnf)%laction ) emp(:,:) = emp(:,:) - frcv(jpr_rnf)%z3(:,:,1) 807 IF( srcv(jpr_cal)%laction ) emp(:,:) = emp(:,:) - frcv(jpr_cal)%z3(:,:,1) 808 ! 809 !!gm : this seems to be internal cooking, not sure to need that in a generic interface 810 !!gm at least should be optional... 811 !! IF( TRIM( sn_rcv_rnf%cldes ) == 'coupled' ) THEN ! add to the total freshwater budget 812 !! ! remove negative runoff 813 !! zcumulpos = SUM( MAX( frcv(jpr_rnf)%z3(:,:,1), 0.e0 ) * e1t(:,:) * e2t(:,:) * tmask_i(:,:) ) 814 !! zcumulneg = SUM( MIN( frcv(jpr_rnf)%z3(:,:,1), 0.e0 ) * e1t(:,:) * e2t(:,:) * tmask_i(:,:) ) 815 !! IF( lk_mpp ) CALL mpp_sum( zcumulpos ) ! sum over the global domain 816 !! IF( lk_mpp ) CALL mpp_sum( zcumulneg ) 817 !! IF( zcumulpos /= 0. ) THEN ! distribute negative runoff on positive runoff grid points 818 !! zcumulneg = 1.e0 + zcumulneg / zcumulpos 819 !! frcv(jpr_rnf)%z3(:,:,1) = MAX( frcv(jpr_rnf)%z3(:,:,1), 0.e0 ) * zcumulneg 820 !! ENDIF 821 !! ! add runoff to e-p 822 !! emp(:,:) = emp(:,:) - frcv(jpr_rnf)%z3(:,:,1) 823 !! ENDIF 824 !!gm end of internal cooking 1075 IF( srcv(jpr_rnf)%laction ) rnf(:,:) = frcv(jpr_rnf)%z3(:,:,1) 1076 IF( srcv(jpr_cal)%laction ) zemp(:,:) = zemp(:,:) - frcv(jpr_cal)%z3(:,:,1) 1077 1078 IF( ln_mixcpl ) THEN ; emp(:,:) = emp(:,:) * xcplmask(:,:,0) + zemp(:,:) * zmsk(:,:) 1079 ELSE ; emp(:,:) = zemp(:,:) 1080 ENDIF 825 1081 ! 826 1082 ! ! non solar heat flux over the ocean (qns) 827 IF( srcv(jpr_qnsoce)%laction ) qns(:,:) = frcv(jpr_qnsoce)%z3(:,:,1) 828 IF( srcv(jpr_qnsmix)%laction ) qns(:,:) = frcv(jpr_qnsmix)%z3(:,:,1) 1083 IF( srcv(jpr_qnsoce)%laction ) THEN ; zqns(:,:) = frcv(jpr_qnsoce)%z3(:,:,1) 1084 ELSE IF( srcv(jpr_qnsmix)%laction ) THEN ; zqns(:,:) = frcv(jpr_qnsmix)%z3(:,:,1) 1085 ELSE ; zqns(:,:) = 0._wp 1086 END IF 829 1087 ! update qns over the free ocean with: 830 qns(:,:) = qns(:,:) - emp(:,:) * sst_m(:,:) * rcp ! remove heat content due to mass flux (assumed to be at SST) 831 IF( srcv(jpr_snow )%laction ) THEN 832 qns(:,:) = qns(:,:) - frcv(jpr_snow)%z3(:,:,1) * lfus ! energy for melting solid precipitation over the free ocean 1088 IF( nn_components /= jp_iam_opa ) THEN 1089 zqns(:,:) = zqns(:,:) - zemp(:,:) * sst_m(:,:) * rcp ! remove heat content due to mass flux (assumed to be at SST) 1090 IF( srcv(jpr_snow )%laction ) THEN 1091 zqns(:,:) = zqns(:,:) - frcv(jpr_snow)%z3(:,:,1) * lfus ! energy for melting solid precipitation over the free ocean 1092 ENDIF 1093 ENDIF 1094 IF( ln_mixcpl ) THEN ; qns(:,:) = qns(:,:) * xcplmask(:,:,0) + zqns(:,:) * zmsk(:,:) 1095 ELSE ; qns(:,:) = zqns(:,:) 833 1096 ENDIF 834 1097 835 1098 ! ! solar flux over the ocean (qsr) 836 IF( srcv(jpr_qsroce)%laction ) qsr(:,:) = frcv(jpr_qsroce)%z3(:,:,1) 837 IF( srcv(jpr_qsrmix)%laction ) qsr(:,:) = frcv(jpr_qsrmix)%z3(:,:,1) 838 IF( ln_dm2dc ) qsr(:,:) = sbc_dcy( qsr ) ! modify qsr to include the diurnal cycle 1099 IF ( srcv(jpr_qsroce)%laction ) THEN ; zqsr(:,:) = frcv(jpr_qsroce)%z3(:,:,1) 1100 ELSE IF( srcv(jpr_qsrmix)%laction ) then ; zqsr(:,:) = frcv(jpr_qsrmix)%z3(:,:,1) 1101 ELSE ; zqsr(:,:) = 0._wp 1102 ENDIF 1103 IF( ln_dm2dc .AND. ln_cpl ) zqsr(:,:) = sbc_dcy( zqsr ) ! modify qsr to include the diurnal cycle 1104 IF( ln_mixcpl ) THEN ; qsr(:,:) = qsr(:,:) * xcplmask(:,:,0) + zqsr(:,:) * zmsk(:,:) 1105 ELSE ; qsr(:,:) = zqsr(:,:) 1106 ENDIF 839 1107 ! 840 841 ENDIF 842 ! 843 CALL wrk_dealloc( jpi,jpj, ztx, zty ) 1108 ! salt flux over the ocean (received by opa in case of opa <-> sas coupling) 1109 IF( srcv(jpr_sflx )%laction ) sfx(:,:) = frcv(jpr_sflx )%z3(:,:,1) 1110 ! Ice cover (received by opa in case of opa <-> sas coupling) 1111 IF( srcv(jpr_fice )%laction ) fr_i(:,:) = frcv(jpr_fice )%z3(:,:,1) 1112 ! 1113 1114 ENDIF 1115 ! 1116 CALL wrk_dealloc( jpi,jpj, ztx, zty, zmsk, zemp, zqns, zqsr ) 844 1117 ! 845 1118 IF( nn_timing == 1 ) CALL timing_stop('sbc_cpl_rcv') … … 938 1211 ! 939 1212 ENDIF 940 941 1213 ! ! ======================= ! 942 1214 ! ! put on ice grid ! … … 1060 1332 1061 1333 1062 SUBROUTINE sbc_cpl_ice_flx( p_frld , palbi , psst , pist)1334 SUBROUTINE sbc_cpl_ice_flx( p_frld, palbi, psst, pist ) 1063 1335 !!---------------------------------------------------------------------- 1064 1336 !! *** ROUTINE sbc_cpl_ice_flx *** … … 1102 1374 REAL(wp), INTENT(in ), DIMENSION(:,:) :: p_frld ! lead fraction [0 to 1] 1103 1375 ! optional arguments, used only in 'mixed oce-ice' case 1104 REAL(wp), INTENT(in ), DIMENSION(:,:,:), OPTIONAL :: palbi ! all skies ice albedo 1105 REAL(wp), INTENT(in ), DIMENSION(:,: ), OPTIONAL :: psst ! sea surface temperature [Celsius] 1106 REAL(wp), INTENT(in ), DIMENSION(:,:,:), OPTIONAL :: pist ! ice surface temperature [Kelvin] 1107 ! 1108 INTEGER :: jl ! dummy loop index 1109 REAL(wp), POINTER, DIMENSION(:,:) :: zcptn, ztmp, zicefr 1376 REAL(wp), INTENT(in ), DIMENSION(:,:,:), OPTIONAL :: palbi ! all skies ice albedo 1377 REAL(wp), INTENT(in ), DIMENSION(:,: ), OPTIONAL :: psst ! sea surface temperature [Celsius] 1378 REAL(wp), INTENT(in ), DIMENSION(:,:,:), OPTIONAL :: pist ! ice surface temperature [Kelvin] 1379 ! 1380 INTEGER :: jl ! dummy loop index 1381 REAL(wp), POINTER, DIMENSION(:,: ) :: zcptn, ztmp, zicefr, zmsk 1382 REAL(wp), POINTER, DIMENSION(:,: ) :: zemp_tot, zemp_ice, zsprecip, ztprecip, zqns_tot, zqsr_tot 1383 REAL(wp), POINTER, DIMENSION(:,:,:) :: zqns_ice, zqsr_ice, zdqns_ice 1384 REAL(wp), POINTER, DIMENSION(:,: ) :: zevap, zsnw, zqns_oce, zqsr_oce, zqprec_ice, zqemp_oce ! for LIM3 1110 1385 !!---------------------------------------------------------------------- 1111 1386 ! 1112 1387 IF( nn_timing == 1 ) CALL timing_start('sbc_cpl_ice_flx') 1113 1388 ! 1114 CALL wrk_alloc( jpi,jpj, zcptn, ztmp, zicefr ) 1115 1389 CALL wrk_alloc( jpi,jpj, zcptn, ztmp, zicefr, zmsk, zemp_tot, zemp_ice, zsprecip, ztprecip, zqns_tot, zqsr_tot ) 1390 CALL wrk_alloc( jpi,jpj,jpl, zqns_ice, zqsr_ice, zdqns_ice ) 1391 1392 IF( ln_mixcpl ) zmsk(:,:) = 1. - xcplmask(:,:,0) 1116 1393 zicefr(:,:) = 1.- p_frld(:,:) 1117 1394 zcptn(:,:) = rcp * sst_m(:,:) … … 1121 1398 ! ! ========================= ! 1122 1399 ! 1123 ! ! total Precipitations - total Evaporation (emp_tot) 1124 ! ! solid precipitation - sublimation (emp_ice) 1125 ! ! solid Precipitation (sprecip) 1400 ! ! total Precipitation - total Evaporation (emp_tot) 1401 ! ! solid precipitation - sublimation (emp_ice) 1402 ! ! solid Precipitation (sprecip) 1403 ! ! liquid + solid Precipitation (tprecip) 1126 1404 SELECT CASE( TRIM( sn_rcv_emp%cldes ) ) 1127 1405 CASE( 'conservative' ) ! received fields: jpr_rain, jpr_snow, jpr_ievp, jpr_tevp 1128 sprecip(:,:) = frcv(jpr_snow)%z3(:,:,1)! May need to ensure positive here1129 tprecip(:,:) = frcv(jpr_rain)%z3(:,:,1) + sprecip (:,:)! May need to ensure positive here1130 emp_tot(:,:) = frcv(jpr_tevp)%z3(:,:,1) -tprecip(:,:)1131 emp_ice(:,:) = frcv(jpr_ievp)%z3(:,:,1) - frcv(jpr_snow)%z3(:,:,1)1406 zsprecip(:,:) = frcv(jpr_snow)%z3(:,:,1) ! May need to ensure positive here 1407 ztprecip(:,:) = frcv(jpr_rain)%z3(:,:,1) + zsprecip(:,:) ! May need to ensure positive here 1408 zemp_tot(:,:) = frcv(jpr_tevp)%z3(:,:,1) - ztprecip(:,:) 1409 zemp_ice(:,:) = frcv(jpr_ievp)%z3(:,:,1) - frcv(jpr_snow)%z3(:,:,1) 1132 1410 CALL iom_put( 'rain' , frcv(jpr_rain)%z3(:,:,1) ) ! liquid precipitation 1133 1411 IF( iom_use('hflx_rain_cea') ) & … … 1140 1418 CALL iom_put( 'hflx_evap_cea', ztmp(:,:) * zcptn(:,:) ) ! heat flux from from evap (cell average) 1141 1419 CASE( 'oce and ice' ) ! received fields: jpr_sbpr, jpr_semp, jpr_oemp, jpr_ievp 1142 emp_tot(:,:) = p_frld(:,:) * frcv(jpr_oemp)%z3(:,:,1) + zicefr(:,:) * frcv(jpr_sbpr)%z3(:,:,1) 1143 emp_ice(:,:) = frcv(jpr_semp)%z3(:,:,1) 1144 sprecip(:,:) = - frcv(jpr_semp)%z3(:,:,1) + frcv(jpr_ievp)%z3(:,:,1) 1420 zemp_tot(:,:) = p_frld(:,:) * frcv(jpr_oemp)%z3(:,:,1) + zicefr(:,:) * frcv(jpr_sbpr)%z3(:,:,1) 1421 zemp_ice(:,:) = frcv(jpr_semp)%z3(:,:,1) 1422 zsprecip(:,:) = frcv(jpr_ievp)%z3(:,:,1) - frcv(jpr_semp)%z3(:,:,1) 1423 ztprecip(:,:) = frcv(jpr_semp)%z3(:,:,1) - frcv(jpr_sbpr)%z3(:,:,1) + zsprecip(:,:) 1145 1424 END SELECT 1425 1426 IF( iom_use('subl_ai_cea') ) & 1427 CALL iom_put( 'subl_ai_cea', frcv(jpr_ievp)%z3(:,:,1) * zicefr(:,:) ) ! Sublimation over sea-ice (cell average) 1428 ! 1429 ! ! runoffs and calving (put in emp_tot) 1430 IF( srcv(jpr_rnf)%laction ) rnf(:,:) = frcv(jpr_rnf)%z3(:,:,1) 1431 IF( srcv(jpr_cal)%laction ) THEN 1432 zemp_tot(:,:) = zemp_tot(:,:) - frcv(jpr_cal)%z3(:,:,1) 1433 CALL iom_put( 'calving_cea', frcv(jpr_cal)%z3(:,:,1) ) 1434 ENDIF 1435 1436 IF( ln_mixcpl ) THEN 1437 emp_tot(:,:) = emp_tot(:,:) * xcplmask(:,:,0) + zemp_tot(:,:) * zmsk(:,:) 1438 emp_ice(:,:) = emp_ice(:,:) * xcplmask(:,:,0) + zemp_ice(:,:) * zmsk(:,:) 1439 sprecip(:,:) = sprecip(:,:) * xcplmask(:,:,0) + zsprecip(:,:) * zmsk(:,:) 1440 tprecip(:,:) = tprecip(:,:) * xcplmask(:,:,0) + ztprecip(:,:) * zmsk(:,:) 1441 ELSE 1442 emp_tot(:,:) = zemp_tot(:,:) 1443 emp_ice(:,:) = zemp_ice(:,:) 1444 sprecip(:,:) = zsprecip(:,:) 1445 tprecip(:,:) = ztprecip(:,:) 1446 ENDIF 1146 1447 1147 1448 CALL iom_put( 'snowpre' , sprecip ) ! Snow … … 1150 1451 IF( iom_use('snow_ai_cea') ) & 1151 1452 CALL iom_put( 'snow_ai_cea', sprecip(:,:) * zicefr(:,:) ) ! Snow over sea-ice (cell average) 1152 IF( iom_use('subl_ai_cea') ) &1153 CALL iom_put( 'subl_ai_cea', frcv(jpr_ievp)%z3(:,:,1) * zicefr(:,:) ) ! Sublimation over sea-ice (cell average)1154 !1155 ! ! runoffs and calving (put in emp_tot)1156 IF( srcv(jpr_rnf)%laction ) THEN1157 emp_tot(:,:) = emp_tot(:,:) - frcv(jpr_rnf)%z3(:,:,1)1158 CALL iom_put( 'runoffs' , frcv(jpr_rnf)%z3(:,:,1) ) ! rivers1159 IF( iom_use('hflx_rnf_cea') ) &1160 CALL iom_put( 'hflx_rnf_cea' , frcv(jpr_rnf)%z3(:,:,1) * zcptn(:,:) ) ! heat flux from rivers1161 ENDIF1162 IF( srcv(jpr_cal)%laction ) THEN1163 emp_tot(:,:) = emp_tot(:,:) - frcv(jpr_cal)%z3(:,:,1)1164 CALL iom_put( 'calving', frcv(jpr_cal)%z3(:,:,1) )1165 ENDIF1166 !1167 !!gm : this seems to be internal cooking, not sure to need that in a generic interface1168 !!gm at least should be optional...1169 !! ! remove negative runoff ! sum over the global domain1170 !! zcumulpos = SUM( MAX( frcv(jpr_rnf)%z3(:,:,1), 0.e0 ) * e1t(:,:) * e2t(:,:) * tmask_i(:,:) )1171 !! zcumulneg = SUM( MIN( frcv(jpr_rnf)%z3(:,:,1), 0.e0 ) * e1t(:,:) * e2t(:,:) * tmask_i(:,:) )1172 !! IF( lk_mpp ) CALL mpp_sum( zcumulpos )1173 !! IF( lk_mpp ) CALL mpp_sum( zcumulneg )1174 !! IF( zcumulpos /= 0. ) THEN ! distribute negative runoff on positive runoff grid points1175 !! zcumulneg = 1.e0 + zcumulneg / zcumulpos1176 !! frcv(jpr_rnf)%z3(:,:,1) = MAX( frcv(jpr_rnf)%z3(:,:,1), 0.e0 ) * zcumulneg1177 !! ENDIF1178 !! emp_tot(:,:) = emp_tot(:,:) - frcv(jpr_rnf)%z3(:,:,1) ! add runoff to e-p1179 !!1180 !!gm end of internal cooking1181 1453 1182 1454 ! ! ========================= ! … … 1184 1456 ! ! ========================= ! 1185 1457 CASE( 'oce only' ) ! the required field is directly provided 1186 qns_tot(:,: ) = frcv(jpr_qnsoce)%z3(:,:,1)1458 zqns_tot(:,: ) = frcv(jpr_qnsoce)%z3(:,:,1) 1187 1459 CASE( 'conservative' ) ! the required fields are directly provided 1188 qns_tot(:,: ) = frcv(jpr_qnsmix)%z3(:,:,1)1460 zqns_tot(:,: ) = frcv(jpr_qnsmix)%z3(:,:,1) 1189 1461 IF ( TRIM(sn_rcv_qns%clcat) == 'yes' ) THEN 1190 qns_ice(:,:,1:jpl) = frcv(jpr_qnsice)%z3(:,:,1:jpl)1462 zqns_ice(:,:,1:jpl) = frcv(jpr_qnsice)%z3(:,:,1:jpl) 1191 1463 ELSE 1192 1464 ! Set all category values equal for the moment 1193 1465 DO jl=1,jpl 1194 qns_ice(:,:,jl) = frcv(jpr_qnsice)%z3(:,:,1)1466 zqns_ice(:,:,jl) = frcv(jpr_qnsice)%z3(:,:,1) 1195 1467 ENDDO 1196 1468 ENDIF 1197 1469 CASE( 'oce and ice' ) ! the total flux is computed from ocean and ice fluxes 1198 qns_tot(:,: ) = p_frld(:,:) * frcv(jpr_qnsoce)%z3(:,:,1)1470 zqns_tot(:,: ) = p_frld(:,:) * frcv(jpr_qnsoce)%z3(:,:,1) 1199 1471 IF ( TRIM(sn_rcv_qns%clcat) == 'yes' ) THEN 1200 1472 DO jl=1,jpl 1201 qns_tot(:,: ) =qns_tot(:,:) + a_i(:,:,jl) * frcv(jpr_qnsice)%z3(:,:,jl)1202 qns_ice(:,:,jl) = frcv(jpr_qnsice)%z3(:,:,jl)1473 zqns_tot(:,: ) = zqns_tot(:,:) + a_i(:,:,jl) * frcv(jpr_qnsice)%z3(:,:,jl) 1474 zqns_ice(:,:,jl) = frcv(jpr_qnsice)%z3(:,:,jl) 1203 1475 ENDDO 1204 1476 ELSE 1477 qns_tot(:,: ) = qns_tot(:,:) + zicefr(:,:) * frcv(jpr_qnsice)%z3(:,:,1) 1205 1478 DO jl=1,jpl 1206 qns_tot(:,: ) =qns_tot(:,:) + zicefr(:,:) * frcv(jpr_qnsice)%z3(:,:,1)1207 qns_ice(:,:,jl) = frcv(jpr_qnsice)%z3(:,:,1)1479 zqns_tot(:,: ) = zqns_tot(:,:) + zicefr(:,:) * frcv(jpr_qnsice)%z3(:,:,1) 1480 zqns_ice(:,:,jl) = frcv(jpr_qnsice)%z3(:,:,1) 1208 1481 ENDDO 1209 1482 ENDIF 1210 1483 CASE( 'mixed oce-ice' ) ! the ice flux is cumputed from the total flux, the SST and ice informations 1211 1484 ! ** NEED TO SORT OUT HOW THIS SHOULD WORK IN THE MULTI-CATEGORY CASE - CURRENTLY NOT ALLOWED WHEN INTERFACE INITIALISED ** 1212 qns_tot(:,: ) = frcv(jpr_qnsmix)%z3(:,:,1)1213 qns_ice(:,:,1) = frcv(jpr_qnsmix)%z3(:,:,1) &1485 zqns_tot(:,: ) = frcv(jpr_qnsmix)%z3(:,:,1) 1486 zqns_ice(:,:,1) = frcv(jpr_qnsmix)%z3(:,:,1) & 1214 1487 & + frcv(jpr_dqnsdt)%z3(:,:,1) * ( pist(:,:,1) - ( (rt0 + psst(:,: ) ) * p_frld(:,:) & 1215 1488 & + pist(:,:,1) * zicefr(:,:) ) ) 1216 1489 END SELECT 1217 ztmp(:,:) = p_frld(:,:) * sprecip(:,:) * lfus1218 qns_tot(:,:) = qns_tot(:,:) & ! qns_tot update over free ocean with:1219 & - ztmp(:,:) & ! remove the latent heat flux of solid precip. melting1220 & - ( emp_tot(:,:) & ! remove the heat content of mass flux (assumed to be at SST)1221 & - emp_ice(:,:) * zicefr(:,:) ) * zcptn(:,:)1222 IF( iom_use('hflx_snow_cea') ) &1223 CALL iom_put( 'hflx_snow_cea', ztmp + sprecip(:,:) * zcptn(:,:) ) ! heat flux from snow (cell average)1224 1490 !!gm 1225 !! currently it is taken into account in leads budget but not in the qns_tot, and thus not in1491 !! currently it is taken into account in leads budget but not in the zqns_tot, and thus not in 1226 1492 !! the flux that enter the ocean.... 1227 1493 !! moreover 1 - it is not diagnose anywhere.... … … 1232 1498 IF( srcv(jpr_cal)%laction ) THEN ! Iceberg melting 1233 1499 ztmp(:,:) = frcv(jpr_cal)%z3(:,:,1) * lfus ! add the latent heat of iceberg melting 1234 qns_tot(:,:) =qns_tot(:,:) - ztmp(:,:)1500 zqns_tot(:,:) = zqns_tot(:,:) - ztmp(:,:) 1235 1501 IF( iom_use('hflx_cal_cea') ) & 1236 1502 CALL iom_put( 'hflx_cal_cea', ztmp + frcv(jpr_cal)%z3(:,:,1) * zcptn(:,:) ) ! heat flux from calving 1237 1503 ENDIF 1504 1505 ztmp(:,:) = p_frld(:,:) * zsprecip(:,:) * lfus 1506 IF( iom_use('hflx_snow_cea') ) CALL iom_put( 'hflx_snow_cea', ztmp + sprecip(:,:) * zcptn(:,:) ) ! heat flux from snow (cell average) 1507 1508 #if defined key_lim3 1509 CALL wrk_alloc( jpi,jpj, zevap, zsnw, zqns_oce, zqprec_ice, zqemp_oce ) 1510 1511 ! --- evaporation --- ! 1512 ! clem: evap_ice is set to 0 for LIM3 since we still do not know what to do with sublimation 1513 ! the problem is: the atm. imposes both mass evaporation and heat removed from the snow/ice 1514 ! but it is incoherent WITH the ice model 1515 DO jl=1,jpl 1516 evap_ice(:,:,jl) = 0._wp ! should be: frcv(jpr_ievp)%z3(:,:,1) 1517 ENDDO 1518 zevap(:,:) = zemp_tot(:,:) + ztprecip(:,:) ! evaporation over ocean 1519 1520 ! --- evaporation minus precipitation --- ! 1521 emp_oce(:,:) = emp_tot(:,:) - emp_ice(:,:) 1522 1523 ! --- non solar flux over ocean --- ! 1524 ! note: p_frld cannot be = 0 since we limit the ice concentration to amax 1525 zqns_oce = 0._wp 1526 WHERE( p_frld /= 0._wp ) zqns_oce(:,:) = ( zqns_tot(:,:) - SUM( a_i * zqns_ice, dim=3 ) ) / p_frld(:,:) 1527 1528 ! --- heat flux associated with emp --- ! 1529 zsnw(:,:) = 0._wp 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 #else 1563 1564 ! clem: this formulation is certainly wrong... but better than it was... 1565 zqns_tot(:,:) = zqns_tot(:,:) & ! zqns_tot update over free ocean with: 1566 & - ztmp(:,:) & ! remove the latent heat flux of solid precip. melting 1567 & - ( zemp_tot(:,:) & ! remove the heat content of mass flux (assumed to be at SST) 1568 & - zemp_ice(:,:) * zicefr(:,:) ) * zcptn(:,:) 1569 1570 IF( ln_mixcpl ) THEN 1571 qns_tot(:,:) = qns(:,:) * p_frld(:,:) + SUM( qns_ice(:,:,:) * a_i(:,:,:), dim=3 ) ! total flux from blk 1572 qns_tot(:,:) = qns_tot(:,:) * xcplmask(:,:,0) + zqns_tot(:,:)* zmsk(:,:) 1573 DO jl=1,jpl 1574 qns_ice(:,:,jl) = qns_ice(:,:,jl) * xcplmask(:,:,0) + zqns_ice(:,:,jl)* zmsk(:,:) 1575 ENDDO 1576 ELSE 1577 qns_tot(:,: ) = zqns_tot(:,: ) 1578 qns_ice(:,:,:) = zqns_ice(:,:,:) 1579 ENDIF 1580 1581 #endif 1238 1582 1239 1583 ! ! ========================= ! … … 1241 1585 ! ! ========================= ! 1242 1586 CASE( 'oce only' ) 1243 qsr_tot(:,: ) = MAX( 0._wp , frcv(jpr_qsroce)%z3(:,:,1) )1587 zqsr_tot(:,: ) = MAX( 0._wp , frcv(jpr_qsroce)%z3(:,:,1) ) 1244 1588 CASE( 'conservative' ) 1245 qsr_tot(:,: ) = frcv(jpr_qsrmix)%z3(:,:,1)1589 zqsr_tot(:,: ) = frcv(jpr_qsrmix)%z3(:,:,1) 1246 1590 IF ( TRIM(sn_rcv_qsr%clcat) == 'yes' ) THEN 1247 qsr_ice(:,:,1:jpl) = frcv(jpr_qsrice)%z3(:,:,1:jpl)1591 zqsr_ice(:,:,1:jpl) = frcv(jpr_qsrice)%z3(:,:,1:jpl) 1248 1592 ELSE 1249 1593 ! Set all category values equal for the moment 1250 1594 DO jl=1,jpl 1251 qsr_ice(:,:,jl) = frcv(jpr_qsrice)%z3(:,:,1)1595 zqsr_ice(:,:,jl) = frcv(jpr_qsrice)%z3(:,:,1) 1252 1596 ENDDO 1253 1597 ENDIF 1254 qsr_tot(:,: ) = frcv(jpr_qsrmix)%z3(:,:,1)1255 qsr_ice(:,:,1) = frcv(jpr_qsrice)%z3(:,:,1)1598 zqsr_tot(:,: ) = frcv(jpr_qsrmix)%z3(:,:,1) 1599 zqsr_ice(:,:,1) = frcv(jpr_qsrice)%z3(:,:,1) 1256 1600 CASE( 'oce and ice' ) 1257 qsr_tot(:,: ) = p_frld(:,:) * frcv(jpr_qsroce)%z3(:,:,1)1601 zqsr_tot(:,: ) = p_frld(:,:) * frcv(jpr_qsroce)%z3(:,:,1) 1258 1602 IF ( TRIM(sn_rcv_qsr%clcat) == 'yes' ) THEN 1259 1603 DO jl=1,jpl 1260 qsr_tot(:,: ) =qsr_tot(:,:) + a_i(:,:,jl) * frcv(jpr_qsrice)%z3(:,:,jl)1261 qsr_ice(:,:,jl) = frcv(jpr_qsrice)%z3(:,:,jl)1604 zqsr_tot(:,: ) = zqsr_tot(:,:) + a_i(:,:,jl) * frcv(jpr_qsrice)%z3(:,:,jl) 1605 zqsr_ice(:,:,jl) = frcv(jpr_qsrice)%z3(:,:,jl) 1262 1606 ENDDO 1263 1607 ELSE 1608 qsr_tot(:,: ) = qsr_tot(:,:) + zicefr(:,:) * frcv(jpr_qsrice)%z3(:,:,1) 1264 1609 DO jl=1,jpl 1265 qsr_tot(:,: ) =qsr_tot(:,:) + zicefr(:,:) * frcv(jpr_qsrice)%z3(:,:,1)1266 qsr_ice(:,:,jl) = frcv(jpr_qsrice)%z3(:,:,1)1610 zqsr_tot(:,: ) = zqsr_tot(:,:) + zicefr(:,:) * frcv(jpr_qsrice)%z3(:,:,1) 1611 zqsr_ice(:,:,jl) = frcv(jpr_qsrice)%z3(:,:,1) 1267 1612 ENDDO 1268 1613 ENDIF 1269 1614 CASE( 'mixed oce-ice' ) 1270 qsr_tot(:,: ) = frcv(jpr_qsrmix)%z3(:,:,1)1615 zqsr_tot(:,: ) = frcv(jpr_qsrmix)%z3(:,:,1) 1271 1616 ! ** NEED TO SORT OUT HOW THIS SHOULD WORK IN THE MULTI-CATEGORY CASE - CURRENTLY NOT ALLOWED WHEN INTERFACE INITIALISED ** 1272 1617 ! Create solar heat flux over ice using incoming solar heat flux and albedos 1273 1618 ! ( see OASIS3 user guide, 5th edition, p39 ) 1274 qsr_ice(:,:,1) = frcv(jpr_qsrmix)%z3(:,:,1) * ( 1.- palbi(:,:,1) ) &1619 zqsr_ice(:,:,1) = frcv(jpr_qsrmix)%z3(:,:,1) * ( 1.- palbi(:,:,1) ) & 1275 1620 & / ( 1.- ( albedo_oce_mix(:,: ) * p_frld(:,:) & 1276 1621 & + palbi (:,:,1) * zicefr(:,:) ) ) 1277 1622 END SELECT 1278 IF( ln_dm2dc ) THEN ! modify qsr to include the diurnal cycle1279 qsr_tot(:,: ) = sbc_dcy(qsr_tot(:,: ) )1623 IF( ln_dm2dc .AND. ln_cpl ) THEN ! modify qsr to include the diurnal cycle 1624 zqsr_tot(:,: ) = sbc_dcy( zqsr_tot(:,: ) ) 1280 1625 DO jl=1,jpl 1281 qsr_ice(:,:,jl) = sbc_dcy(qsr_ice(:,:,jl) )1626 zqsr_ice(:,:,jl) = sbc_dcy( zqsr_ice(:,:,jl) ) 1282 1627 ENDDO 1628 ENDIF 1629 1630 #if defined key_lim3 1631 CALL wrk_alloc( jpi,jpj, zqsr_oce ) 1632 ! --- solar flux over ocean --- ! 1633 ! note: p_frld cannot be = 0 since we limit the ice concentration to amax 1634 zqsr_oce = 0._wp 1635 WHERE( p_frld /= 0._wp ) zqsr_oce(:,:) = ( zqsr_tot(:,:) - SUM( a_i * zqsr_ice, dim=3 ) ) / p_frld(:,:) 1636 1637 IF( ln_mixcpl ) THEN ; qsr_oce(:,:) = qsr_oce(:,:) * xcplmask(:,:,0) + zqsr_oce(:,:)* zmsk(:,:) 1638 ELSE ; qsr_oce(:,:) = zqsr_oce(:,:) ; ENDIF 1639 1640 CALL wrk_dealloc( jpi,jpj, zqsr_oce ) 1641 #endif 1642 1643 IF( ln_mixcpl ) THEN 1644 qsr_tot(:,:) = qsr(:,:) * p_frld(:,:) + SUM( qsr_ice(:,:,:) * a_i(:,:,:), dim=3 ) ! total flux from blk 1645 qsr_tot(:,:) = qsr_tot(:,:) * xcplmask(:,:,0) + zqsr_tot(:,:)* zmsk(:,:) 1646 DO jl=1,jpl 1647 qsr_ice(:,:,jl) = qsr_ice(:,:,jl) * xcplmask(:,:,0) + zqsr_ice(:,:,jl)* zmsk(:,:) 1648 ENDDO 1649 ELSE 1650 qsr_tot(:,: ) = zqsr_tot(:,: ) 1651 qsr_ice(:,:,:) = zqsr_ice(:,:,:) 1283 1652 ENDIF 1284 1653 … … 1288 1657 CASE ('coupled') 1289 1658 IF ( TRIM(sn_rcv_dqnsdt%clcat) == 'yes' ) THEN 1290 dqns_ice(:,:,1:jpl) = frcv(jpr_dqnsdt)%z3(:,:,1:jpl)1659 zdqns_ice(:,:,1:jpl) = frcv(jpr_dqnsdt)%z3(:,:,1:jpl) 1291 1660 ELSE 1292 1661 ! Set all category values equal for the moment 1293 1662 DO jl=1,jpl 1294 dqns_ice(:,:,jl) = frcv(jpr_dqnsdt)%z3(:,:,1)1663 zdqns_ice(:,:,jl) = frcv(jpr_dqnsdt)%z3(:,:,1) 1295 1664 ENDDO 1296 1665 ENDIF 1297 1666 END SELECT 1298 1667 1668 IF( ln_mixcpl ) THEN 1669 DO jl=1,jpl 1670 dqns_ice(:,:,jl) = dqns_ice(:,:,jl) * xcplmask(:,:,0) + zdqns_ice(:,:,jl) * zmsk(:,:) 1671 ENDDO 1672 ELSE 1673 dqns_ice(:,:,:) = zdqns_ice(:,:,:) 1674 ENDIF 1675 1299 1676 ! ! ========================= ! 1300 1677 SELECT CASE( TRIM( sn_rcv_iceflx%cldes ) ) ! topmelt and botmelt ! … … 1312 1689 fr2_i0(:,:) = ( 0.82 * ( 1.0 - cldf_ice ) + 0.65 * cldf_ice ) 1313 1690 1314 CALL wrk_dealloc( jpi,jpj, zcptn, ztmp, zicefr ) 1691 CALL wrk_dealloc( jpi,jpj, zcptn, ztmp, zicefr, zmsk, zemp_tot, zemp_ice, zsprecip, ztprecip, zqns_tot, zqsr_tot ) 1692 CALL wrk_dealloc( jpi,jpj,jpl, zqns_ice, zqsr_ice, zdqns_ice ) 1315 1693 ! 1316 1694 IF( nn_timing == 1 ) CALL timing_stop('sbc_cpl_ice_flx') … … 1332 1710 INTEGER :: ji, jj, jl ! dummy loop indices 1333 1711 INTEGER :: isec, info ! local integer 1712 REAL(wp) :: zumax, zvmax 1334 1713 REAL(wp), POINTER, DIMENSION(:,:) :: zfr_l, ztmp1, ztmp2, zotx1, zoty1, zotz1, zitx1, zity1, zitz1 1335 1714 REAL(wp), POINTER, DIMENSION(:,:,:) :: ztmp3, ztmp4 … … 1348 1727 ! ! ------------------------- ! 1349 1728 IF( ssnd(jps_toce)%laction .OR. ssnd(jps_tice)%laction .OR. ssnd(jps_tmix)%laction ) THEN 1350 SELECT CASE( sn_snd_temp%cldes) 1351 CASE( 'oce only' ) ; ztmp1(:,:) = tsn(:,:,1,jp_tem) + rt0 1352 CASE( 'weighted oce and ice' ) ; ztmp1(:,:) = ( tsn(:,:,1,jp_tem) + rt0 ) * zfr_l(:,:) 1353 SELECT CASE( sn_snd_temp%clcat ) 1354 CASE( 'yes' ) 1355 ztmp3(:,:,1:jpl) = tn_ice(:,:,1:jpl) * a_i(:,:,1:jpl) 1356 CASE( 'no' ) 1357 ztmp3(:,:,:) = 0.0 1729 1730 IF ( nn_components == jp_iam_opa ) THEN 1731 ztmp1(:,:) = tsn(:,:,1,jp_tem) ! send temperature as it is (potential or conservative) -> use of ln_useCT on the received part 1732 ELSE 1733 ! we must send the surface potential temperature 1734 IF( ln_useCT ) THEN ; ztmp1(:,:) = eos_pt_from_ct( tsn(:,:,1,jp_tem), tsn(:,:,1,jp_sal) ) 1735 ELSE ; ztmp1(:,:) = tsn(:,:,1,jp_tem) 1736 ENDIF 1737 ! 1738 SELECT CASE( sn_snd_temp%cldes) 1739 CASE( 'oce only' ) ; ztmp1(:,:) = ztmp1(:,:) + rt0 1740 CASE( 'oce and ice' ) ; ztmp1(:,:) = ztmp1(:,:) + rt0 1741 SELECT CASE( sn_snd_temp%clcat ) 1742 CASE( 'yes' ) 1743 ztmp3(:,:,1:jpl) = tn_ice(:,:,1:jpl) 1744 CASE( 'no' ) 1745 WHERE( SUM( a_i, dim=3 ) /= 0. ) 1746 ztmp3(:,:,1) = SUM( tn_ice * a_i, dim=3 ) / SUM( a_i, dim=3 ) 1747 ELSEWHERE 1748 ztmp3(:,:,1) = rt0 ! TODO: Is freezing point a good default? (Maybe SST is better?) 1749 END WHERE 1750 CASE default ; CALL ctl_stop( 'sbc_cpl_snd: wrong definition of sn_snd_temp%clcat' ) 1751 END SELECT 1752 CASE( 'weighted oce and ice' ) ; ztmp1(:,:) = ( ztmp1(:,:) + rt0 ) * zfr_l(:,:) 1753 SELECT CASE( sn_snd_temp%clcat ) 1754 CASE( 'yes' ) 1755 ztmp3(:,:,1:jpl) = tn_ice(:,:,1:jpl) * a_i(:,:,1:jpl) 1756 CASE( 'no' ) 1757 ztmp3(:,:,:) = 0.0 1758 DO jl=1,jpl 1759 ztmp3(:,:,1) = ztmp3(:,:,1) + tn_ice(:,:,jl) * a_i(:,:,jl) 1760 ENDDO 1761 CASE default ; CALL ctl_stop( 'sbc_cpl_snd: wrong definition of sn_snd_temp%clcat' ) 1762 END SELECT 1763 CASE( 'mixed oce-ice' ) 1764 ztmp1(:,:) = ( ztmp1(:,:) + rt0 ) * zfr_l(:,:) 1358 1765 DO jl=1,jpl 1359 ztmp 3(:,:,1) = ztmp3(:,:,1) + tn_ice(:,:,jl) * a_i(:,:,jl)1766 ztmp1(:,:) = ztmp1(:,:) + tn_ice(:,:,jl) * a_i(:,:,jl) 1360 1767 ENDDO 1361 CASE default ; CALL ctl_stop( 'sbc_cpl_snd: wrong definition of sn_snd_temp%clcat' )1768 CASE default ; CALL ctl_stop( 'sbc_cpl_snd: wrong definition of sn_snd_temp%cldes' ) 1362 1769 END SELECT 1363 CASE( 'mixed oce-ice' ) 1364 ztmp1(:,:) = ( tsn(:,:,1,jp_tem) + rt0 ) * zfr_l(:,:) 1365 DO jl=1,jpl 1366 ztmp1(:,:) = ztmp1(:,:) + tn_ice(:,:,jl) * a_i(:,:,jl) 1367 ENDDO 1368 CASE default ; CALL ctl_stop( 'sbc_cpl_snd: wrong definition of sn_snd_temp%cldes' ) 1369 END SELECT 1770 ENDIF 1370 1771 IF( ssnd(jps_toce)%laction ) CALL cpl_snd( jps_toce, isec, RESHAPE ( ztmp1, (/jpi,jpj,1/) ), info ) 1371 1772 IF( ssnd(jps_tice)%laction ) CALL cpl_snd( jps_tice, isec, ztmp3, info ) … … 1376 1777 ! ! ------------------------- ! 1377 1778 IF( ssnd(jps_albice)%laction ) THEN ! ice 1378 ztmp3(:,:,1:jpl) = alb_ice(:,:,1:jpl) * a_i(:,:,1:jpl) 1779 SELECT CASE( sn_snd_alb%cldes ) 1780 CASE( 'ice' ) ; ztmp3(:,:,1:jpl) = alb_ice(:,:,1:jpl) 1781 CASE( 'weighted ice' ) ; ztmp3(:,:,1:jpl) = alb_ice(:,:,1:jpl) * a_i(:,:,1:jpl) 1782 CASE default ; CALL ctl_stop( 'sbc_cpl_snd: wrong definition of sn_snd_alb%cldes' ) 1783 END SELECT 1379 1784 CALL cpl_snd( jps_albice, isec, ztmp3, info ) 1380 1785 ENDIF … … 1389 1794 ! ! Ice fraction & Thickness ! 1390 1795 ! ! ------------------------- ! 1391 ! Send ice fraction field 1796 ! Send ice fraction field to atmosphere 1392 1797 IF( ssnd(jps_fice)%laction ) THEN 1393 1798 SELECT CASE( sn_snd_thick%clcat ) … … 1396 1801 CASE default ; CALL ctl_stop( 'sbc_cpl_snd: wrong definition of sn_snd_thick%clcat' ) 1397 1802 END SELECT 1398 CALL cpl_snd( jps_fice, isec, ztmp3, info ) 1803 IF( ssnd(jps_fice)%laction ) CALL cpl_snd( jps_fice, isec, ztmp3, info ) 1804 ENDIF 1805 1806 ! Send ice fraction field to OPA (sent by SAS in SAS-OPA coupling) 1807 IF( ssnd(jps_fice2)%laction ) THEN 1808 ztmp3(:,:,1) = fr_i(:,:) 1809 IF( ssnd(jps_fice2)%laction ) CALL cpl_snd( jps_fice2, isec, ztmp3, info ) 1399 1810 ENDIF 1400 1811 … … 1417 1828 END SELECT 1418 1829 CASE( 'ice and snow' ) 1419 ztmp3(:,:,1:jpl) = ht_i(:,:,1:jpl) 1420 ztmp4(:,:,1:jpl) = ht_s(:,:,1:jpl) 1830 SELECT CASE( sn_snd_thick%clcat ) 1831 CASE( 'yes' ) 1832 ztmp3(:,:,1:jpl) = ht_i(:,:,1:jpl) 1833 ztmp4(:,:,1:jpl) = ht_s(:,:,1:jpl) 1834 CASE( 'no' ) 1835 WHERE( SUM( a_i, dim=3 ) /= 0. ) 1836 ztmp3(:,:,1) = SUM( ht_i * a_i, dim=3 ) / SUM( a_i, dim=3 ) 1837 ztmp4(:,:,1) = SUM( ht_s * a_i, dim=3 ) / SUM( a_i, dim=3 ) 1838 ELSEWHERE 1839 ztmp3(:,:,1) = 0. 1840 ztmp4(:,:,1) = 0. 1841 END WHERE 1842 CASE default ; CALL ctl_stop( 'sbc_cpl_snd: wrong definition of sn_snd_thick%clcat' ) 1843 END SELECT 1421 1844 CASE default ; CALL ctl_stop( 'sbc_cpl_snd: wrong definition of sn_snd_thick%cldes' ) 1422 1845 END SELECT … … 1444 1867 ! i-1 i i 1445 1868 ! i i+1 (for I) 1446 SELECT CASE( TRIM( sn_snd_crt%cldes ) ) 1447 CASE( 'oce only' ) ! C-grid ==> T 1448 DO jj = 2, jpjm1 1449 DO ji = fs_2, fs_jpim1 ! vector opt. 1450 zotx1(ji,jj) = 0.5 * ( un(ji,jj,1) + un(ji-1,jj ,1) ) 1451 zoty1(ji,jj) = 0.5 * ( vn(ji,jj,1) + vn(ji ,jj-1,1) ) 1452 END DO 1453 END DO 1454 CASE( 'weighted oce and ice' ) 1455 SELECT CASE ( cp_ice_msh ) 1456 CASE( 'C' ) ! Ocean and Ice on C-grid ==> T 1869 IF( nn_components == jp_iam_opa ) THEN 1870 zotx1(:,:) = un(:,:,1) 1871 zoty1(:,:) = vn(:,:,1) 1872 ELSE 1873 SELECT CASE( TRIM( sn_snd_crt%cldes ) ) 1874 CASE( 'oce only' ) ! C-grid ==> T 1457 1875 DO jj = 2, jpjm1 1458 1876 DO ji = fs_2, fs_jpim1 ! vector opt. 1459 zotx1(ji,jj) = 0.5 * ( un (ji,jj,1) + un (ji-1,jj ,1) ) * zfr_l(ji,jj) 1460 zoty1(ji,jj) = 0.5 * ( vn (ji,jj,1) + vn (ji ,jj-1,1) ) * zfr_l(ji,jj) 1461 zitx1(ji,jj) = 0.5 * ( u_ice(ji,jj ) + u_ice(ji-1,jj ) ) * fr_i(ji,jj) 1462 zity1(ji,jj) = 0.5 * ( v_ice(ji,jj ) + v_ice(ji ,jj-1 ) ) * fr_i(ji,jj) 1877 zotx1(ji,jj) = 0.5 * ( un(ji,jj,1) + un(ji-1,jj ,1) ) 1878 zoty1(ji,jj) = 0.5 * ( vn(ji,jj,1) + vn(ji ,jj-1,1) ) 1463 1879 END DO 1464 1880 END DO 1465 CASE( 'I' ) ! Ocean on C grid, Ice on I-point (B-grid) ==> T 1466 DO jj = 2, jpjm1 1467 DO ji = 2, jpim1 ! NO vector opt. 1468 zotx1(ji,jj) = 0.5 * ( un(ji,jj,1) + un(ji-1,jj ,1) ) * zfr_l(ji,jj) 1469 zoty1(ji,jj) = 0.5 * ( vn(ji,jj,1) + vn(ji ,jj-1,1) ) * zfr_l(ji,jj) 1470 zitx1(ji,jj) = 0.25 * ( u_ice(ji+1,jj+1) + u_ice(ji,jj+1) & 1471 & + u_ice(ji+1,jj ) + u_ice(ji,jj ) ) * fr_i(ji,jj) 1472 zity1(ji,jj) = 0.25 * ( v_ice(ji+1,jj+1) + v_ice(ji,jj+1) & 1473 & + v_ice(ji+1,jj ) + v_ice(ji,jj ) ) * fr_i(ji,jj) 1881 CASE( 'weighted oce and ice' ) 1882 SELECT CASE ( cp_ice_msh ) 1883 CASE( 'C' ) ! Ocean and Ice on C-grid ==> T 1884 DO jj = 2, jpjm1 1885 DO ji = fs_2, fs_jpim1 ! vector opt. 1886 zotx1(ji,jj) = 0.5 * ( un (ji,jj,1) + un (ji-1,jj ,1) ) * zfr_l(ji,jj) 1887 zoty1(ji,jj) = 0.5 * ( vn (ji,jj,1) + vn (ji ,jj-1,1) ) * zfr_l(ji,jj) 1888 zitx1(ji,jj) = 0.5 * ( u_ice(ji,jj ) + u_ice(ji-1,jj ) ) * fr_i(ji,jj) 1889 zity1(ji,jj) = 0.5 * ( v_ice(ji,jj ) + v_ice(ji ,jj-1 ) ) * fr_i(ji,jj) 1890 END DO 1474 1891 END DO 1475 END DO1476 CASE( 'F' ) ! Ocean on C grid, Ice on F-point (B-grid) ==> T1477 DO jj = 2, jpjm11478 DO ji = 2, jpim1 ! NO vector opt.1479 zotx1(ji,jj) = 0.5 * ( un(ji,jj,1) + un(ji-1,jj,1) ) * zfr_l(ji,jj)1480 zoty1(ji,jj) = 0.5 * ( vn(ji,jj,1) + vn(ji ,jj-1,1) ) * zfr_l(ji,jj)1481 zitx1(ji,jj) = 0.25 * ( u_ice(ji-1,jj-1) + u_ice(ji,jj-1) &1482 & + u_ice(ji-1,jj ) + u_ice(ji,jj ) ) * fr_i(ji,jj)1483 zity1(ji,jj) = 0.25 * ( v_ice(ji-1,jj-1) + v_ice(ji,jj-1) &1484 & + v_ice(ji-1,jj ) + v_ice(ji,jj ) ) * fr_i(ji,jj)1892 CASE( 'I' ) ! Ocean on C grid, Ice on I-point (B-grid) ==> T 1893 DO jj = 2, jpjm1 1894 DO ji = 2, jpim1 ! NO vector opt. 1895 zotx1(ji,jj) = 0.5 * ( un(ji,jj,1) + un(ji-1,jj ,1) ) * zfr_l(ji,jj) 1896 zoty1(ji,jj) = 0.5 * ( vn(ji,jj,1) + vn(ji ,jj-1,1) ) * zfr_l(ji,jj) 1897 zitx1(ji,jj) = 0.25 * ( u_ice(ji+1,jj+1) + u_ice(ji,jj+1) & 1898 & + u_ice(ji+1,jj ) + u_ice(ji,jj ) ) * fr_i(ji,jj) 1899 zity1(ji,jj) = 0.25 * ( v_ice(ji+1,jj+1) + v_ice(ji,jj+1) & 1900 & + v_ice(ji+1,jj ) + v_ice(ji,jj ) ) * fr_i(ji,jj) 1901 END DO 1485 1902 END DO 1486 END DO 1903 CASE( 'F' ) ! Ocean on C grid, Ice on F-point (B-grid) ==> T 1904 DO jj = 2, jpjm1 1905 DO ji = 2, jpim1 ! NO vector opt. 1906 zotx1(ji,jj) = 0.5 * ( un(ji,jj,1) + un(ji-1,jj ,1) ) * zfr_l(ji,jj) 1907 zoty1(ji,jj) = 0.5 * ( vn(ji,jj,1) + vn(ji ,jj-1,1) ) * zfr_l(ji,jj) 1908 zitx1(ji,jj) = 0.25 * ( u_ice(ji-1,jj-1) + u_ice(ji,jj-1) & 1909 & + u_ice(ji-1,jj ) + u_ice(ji,jj ) ) * fr_i(ji,jj) 1910 zity1(ji,jj) = 0.25 * ( v_ice(ji-1,jj-1) + v_ice(ji,jj-1) & 1911 & + v_ice(ji-1,jj ) + v_ice(ji,jj ) ) * fr_i(ji,jj) 1912 END DO 1913 END DO 1914 END SELECT 1915 CALL lbc_lnk( zitx1, 'T', -1. ) ; CALL lbc_lnk( zity1, 'T', -1. ) 1916 CASE( 'mixed oce-ice' ) 1917 SELECT CASE ( cp_ice_msh ) 1918 CASE( 'C' ) ! Ocean and Ice on C-grid ==> T 1919 DO jj = 2, jpjm1 1920 DO ji = fs_2, fs_jpim1 ! vector opt. 1921 zotx1(ji,jj) = 0.5 * ( un (ji,jj,1) + un (ji-1,jj ,1) ) * zfr_l(ji,jj) & 1922 & + 0.5 * ( u_ice(ji,jj ) + u_ice(ji-1,jj ) ) * fr_i(ji,jj) 1923 zoty1(ji,jj) = 0.5 * ( vn (ji,jj,1) + vn (ji ,jj-1,1) ) * zfr_l(ji,jj) & 1924 & + 0.5 * ( v_ice(ji,jj ) + v_ice(ji ,jj-1 ) ) * fr_i(ji,jj) 1925 END DO 1926 END DO 1927 CASE( 'I' ) ! Ocean on C grid, Ice on I-point (B-grid) ==> T 1928 DO jj = 2, jpjm1 1929 DO ji = 2, jpim1 ! NO vector opt. 1930 zotx1(ji,jj) = 0.5 * ( un(ji,jj,1) + un(ji-1,jj ,1) ) * zfr_l(ji,jj) & 1931 & + 0.25 * ( u_ice(ji+1,jj+1) + u_ice(ji,jj+1) & 1932 & + u_ice(ji+1,jj ) + u_ice(ji,jj ) ) * fr_i(ji,jj) 1933 zoty1(ji,jj) = 0.5 * ( vn(ji,jj,1) + vn(ji ,jj-1,1) ) * zfr_l(ji,jj) & 1934 & + 0.25 * ( v_ice(ji+1,jj+1) + v_ice(ji,jj+1) & 1935 & + v_ice(ji+1,jj ) + v_ice(ji,jj ) ) * fr_i(ji,jj) 1936 END DO 1937 END DO 1938 CASE( 'F' ) ! Ocean on C grid, Ice on F-point (B-grid) ==> T 1939 DO jj = 2, jpjm1 1940 DO ji = 2, jpim1 ! NO vector opt. 1941 zotx1(ji,jj) = 0.5 * ( un(ji,jj,1) + un(ji-1,jj ,1) ) * zfr_l(ji,jj) & 1942 & + 0.25 * ( u_ice(ji-1,jj-1) + u_ice(ji,jj-1) & 1943 & + u_ice(ji-1,jj ) + u_ice(ji,jj ) ) * fr_i(ji,jj) 1944 zoty1(ji,jj) = 0.5 * ( vn(ji,jj,1) + vn(ji ,jj-1,1) ) * zfr_l(ji,jj) & 1945 & + 0.25 * ( v_ice(ji-1,jj-1) + v_ice(ji,jj-1) & 1946 & + v_ice(ji-1,jj ) + v_ice(ji,jj ) ) * fr_i(ji,jj) 1947 END DO 1948 END DO 1949 END SELECT 1487 1950 END SELECT 1488 CALL lbc_lnk( zitx1, 'T', -1. ) ; CALL lbc_lnk( zity1, 'T', -1. ) 1489 CASE( 'mixed oce-ice' ) 1490 SELECT CASE ( cp_ice_msh ) 1491 CASE( 'C' ) ! Ocean and Ice on C-grid ==> T 1492 DO jj = 2, jpjm1 1493 DO ji = fs_2, fs_jpim1 ! vector opt. 1494 zotx1(ji,jj) = 0.5 * ( un (ji,jj,1) + un (ji-1,jj ,1) ) * zfr_l(ji,jj) & 1495 & + 0.5 * ( u_ice(ji,jj ) + u_ice(ji-1,jj ) ) * fr_i(ji,jj) 1496 zoty1(ji,jj) = 0.5 * ( vn (ji,jj,1) + vn (ji ,jj-1,1) ) * zfr_l(ji,jj) & 1497 & + 0.5 * ( v_ice(ji,jj ) + v_ice(ji ,jj-1 ) ) * fr_i(ji,jj) 1498 END DO 1499 END DO 1500 CASE( 'I' ) ! Ocean on C grid, Ice on I-point (B-grid) ==> T 1501 DO jj = 2, jpjm1 1502 DO ji = 2, jpim1 ! NO vector opt. 1503 zotx1(ji,jj) = 0.5 * ( un(ji,jj,1) + un(ji-1,jj ,1) ) * zfr_l(ji,jj) & 1504 & + 0.25 * ( u_ice(ji+1,jj+1) + u_ice(ji,jj+1) & 1505 & + u_ice(ji+1,jj ) + u_ice(ji,jj ) ) * fr_i(ji,jj) 1506 zoty1(ji,jj) = 0.5 * ( vn(ji,jj,1) + vn(ji ,jj-1,1) ) * zfr_l(ji,jj) & 1507 & + 0.25 * ( v_ice(ji+1,jj+1) + v_ice(ji,jj+1) & 1508 & + v_ice(ji+1,jj ) + v_ice(ji,jj ) ) * fr_i(ji,jj) 1509 END DO 1510 END DO 1511 CASE( 'F' ) ! Ocean on C grid, Ice on F-point (B-grid) ==> T 1512 DO jj = 2, jpjm1 1513 DO ji = 2, jpim1 ! NO vector opt. 1514 zotx1(ji,jj) = 0.5 * ( un(ji,jj,1) + un(ji-1,jj ,1) ) * zfr_l(ji,jj) & 1515 & + 0.25 * ( u_ice(ji-1,jj-1) + u_ice(ji,jj-1) & 1516 & + u_ice(ji-1,jj ) + u_ice(ji,jj ) ) * fr_i(ji,jj) 1517 zoty1(ji,jj) = 0.5 * ( vn(ji,jj,1) + vn(ji ,jj-1,1) ) * zfr_l(ji,jj) & 1518 & + 0.25 * ( v_ice(ji-1,jj-1) + v_ice(ji,jj-1) & 1519 & + v_ice(ji-1,jj ) + v_ice(ji,jj ) ) * fr_i(ji,jj) 1520 END DO 1521 END DO 1522 END SELECT 1523 END SELECT 1524 CALL lbc_lnk( zotx1, ssnd(jps_ocx1)%clgrid, -1. ) ; CALL lbc_lnk( zoty1, ssnd(jps_ocy1)%clgrid, -1. ) 1951 CALL lbc_lnk( zotx1, ssnd(jps_ocx1)%clgrid, -1. ) ; CALL lbc_lnk( zoty1, ssnd(jps_ocy1)%clgrid, -1. ) 1952 ! 1953 ENDIF 1525 1954 ! 1526 1955 ! … … 1562 1991 ENDIF 1563 1992 ! 1993 ! 1994 ! Fields sent by OPA to SAS when doing OPA<->SAS coupling 1995 ! ! SSH 1996 IF( ssnd(jps_ssh )%laction ) THEN 1997 ! ! removed inverse barometer ssh when Patm 1998 ! forcing is used (for sea-ice dynamics) 1999 IF( ln_apr_dyn ) THEN ; ztmp1(:,:) = sshb(:,:) - 0.5 * ( ssh_ib(:,:) + ssh_ibb(:,:) ) 2000 ELSE ; ztmp1(:,:) = sshn(:,:) 2001 ENDIF 2002 CALL cpl_snd( jps_ssh , isec, RESHAPE ( ztmp1 , (/jpi,jpj,1/) ), info ) 2003 2004 ENDIF 2005 ! ! SSS 2006 IF( ssnd(jps_soce )%laction ) THEN 2007 CALL cpl_snd( jps_soce , isec, RESHAPE ( tsn(:,:,1,jp_sal), (/jpi,jpj,1/) ), info ) 2008 ENDIF 2009 ! ! first T level thickness 2010 IF( ssnd(jps_e3t1st )%laction ) THEN 2011 CALL cpl_snd( jps_e3t1st, isec, RESHAPE ( fse3t_n(:,:,1) , (/jpi,jpj,1/) ), info ) 2012 ENDIF 2013 ! ! Qsr fraction 2014 IF( ssnd(jps_fraqsr)%laction ) THEN 2015 CALL cpl_snd( jps_fraqsr, isec, RESHAPE ( fraqsr_1lev(:,:) , (/jpi,jpj,1/) ), info ) 2016 ENDIF 2017 ! 2018 ! Fields sent by SAS to OPA when OASIS coupling 2019 ! ! Solar heat flux 2020 IF( ssnd(jps_qsroce)%laction ) CALL cpl_snd( jps_qsroce, isec, RESHAPE ( qsr , (/jpi,jpj,1/) ), info ) 2021 IF( ssnd(jps_qnsoce)%laction ) CALL cpl_snd( jps_qnsoce, isec, RESHAPE ( qns , (/jpi,jpj,1/) ), info ) 2022 IF( ssnd(jps_oemp )%laction ) CALL cpl_snd( jps_oemp , isec, RESHAPE ( emp , (/jpi,jpj,1/) ), info ) 2023 IF( ssnd(jps_sflx )%laction ) CALL cpl_snd( jps_sflx , isec, RESHAPE ( sfx , (/jpi,jpj,1/) ), info ) 2024 IF( ssnd(jps_otx1 )%laction ) CALL cpl_snd( jps_otx1 , isec, RESHAPE ( utau, (/jpi,jpj,1/) ), info ) 2025 IF( ssnd(jps_oty1 )%laction ) CALL cpl_snd( jps_oty1 , isec, RESHAPE ( vtau, (/jpi,jpj,1/) ), info ) 2026 IF( ssnd(jps_rnf )%laction ) CALL cpl_snd( jps_rnf , isec, RESHAPE ( rnf , (/jpi,jpj,1/) ), info ) 2027 IF( ssnd(jps_taum )%laction ) CALL cpl_snd( jps_taum , isec, RESHAPE ( taum, (/jpi,jpj,1/) ), info ) 2028 1564 2029 CALL wrk_dealloc( jpi,jpj, zfr_l, ztmp1, ztmp2, zotx1, zoty1, zotz1, zitx1, zity1, zitz1 ) 1565 2030 CALL wrk_dealloc( jpi,jpj,jpl, ztmp3, ztmp4 )
Note: See TracChangeset
for help on using the changeset viewer.