Changeset 5630 for branches/UKMO/dev_r5107_restart_func_and_date/NEMOGCM/NEMO/OPA_SRC/SBC/sbccpl.F90
- Timestamp:
- 2015-07-23T18:05:51+02:00 (9 years ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
branches/UKMO/dev_r5107_restart_func_and_date/NEMOGCM/NEMO/OPA_SRC/SBC/sbccpl.F90
r5500 r5630 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 … … 477 603 ssnd(jps_tmix)%clname = 'O_TepMix' 478 604 SELECT CASE( TRIM( sn_snd_temp%cldes ) ) 479 CASE( 'none' ) ! nothing to do480 CASE( 'oce only' ) ; ssnd( jps_toce)%laction = .TRUE.481 CASE( ' weighted oce and ice' )605 CASE( 'none' ) ! nothing to do 606 CASE( 'oce only' ) ; ssnd( jps_toce )%laction = .TRUE. 607 CASE( 'oce and ice' , 'weighted oce and ice' ) 482 608 ssnd( (/jps_toce, jps_tice/) )%laction = .TRUE. 483 609 IF ( TRIM( sn_snd_temp%clcat ) == 'yes' ) ssnd(jps_tice)%nct = jpl 484 CASE( 'mixed oce-ice' ) ; ssnd( jps_tmix)%laction = .TRUE.610 CASE( 'mixed oce-ice' ) ; ssnd( jps_tmix )%laction = .TRUE. 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 ! … … 492 618 ssnd(jps_albmix)%clname = 'O_AlbMix' 493 619 SELECT CASE( TRIM( sn_snd_alb%cldes ) ) 494 CASE( 'none' )! nothing to do495 CASE( ' weighted ice' ) ;ssnd(jps_albice)%laction = .TRUE.496 CASE( 'mixed oce-ice' ) ;ssnd(jps_albmix)%laction = .TRUE.620 CASE( 'none' ) ! nothing to do 621 CASE( 'ice' , 'weighted ice' ) ; ssnd(jps_albice)%laction = .TRUE. 622 CASE( 'mixed oce-ice' ) ; ssnd(jps_albmix)%laction = .TRUE. 497 623 CASE default ; CALL ctl_stop( 'sbc_cpl_init: wrong definition of sn_snd_alb%cldes' ) 498 624 END SELECT … … 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 … … 525 651 IF ( TRIM( sn_snd_thick%clcat ) == 'yes' ) THEN 526 652 ssnd(jps_hice:jps_hsnw)%nct = jpl 527 ELSE528 IF ( jpl > 1 ) THEN529 CALL ctl_stop( 'sbc_cpl_init: use weighted ice and snow option for sn_snd_thick%cldes if not exchanging category fields' )530 ENDIF531 653 ENDIF 532 654 CASE ( 'weighted ice and snow' ) … … 567 689 ! ! ------------------------- ! 568 690 ssnd(jps_co2)%clname = 'O_CO2FLX' ; IF( TRIM(sn_snd_co2%cldes) == 'coupled' ) ssnd(jps_co2 )%laction = .TRUE. 691 692 ! ! ------------------------------- ! 693 ! ! OPA-SAS coupling - snd by opa ! 694 ! ! ------------------------------- ! 695 ssnd(jps_ssh )%clname = 'O_SSHght' 696 ssnd(jps_soce )%clname = 'O_SSSal' 697 ssnd(jps_e3t1st)%clname = 'O_E3T1st' 698 ssnd(jps_fraqsr)%clname = 'O_FraQsr' 699 ! 700 IF( nn_components == jp_iam_opa ) THEN 701 ssnd(:)%laction = .FALSE. ! force default definition in case of opa <-> sas coupling 702 ssnd( (/jps_toce, jps_soce, jps_ssh, jps_fraqsr, jps_ocx1, jps_ocy1/) )%laction = .TRUE. 703 ssnd( jps_e3t1st )%laction = lk_vvl 704 ! vector definition: not used but cleaner... 705 ssnd(jps_ocx1)%clgrid = 'U' ! oce components given at U-point 706 ssnd(jps_ocy1)%clgrid = 'V' ! and V-point 707 sn_snd_crt%clvgrd = 'U,V' 708 sn_snd_crt%clvor = 'local grid' 709 sn_snd_crt%clvref = 'spherical' 710 ! 711 IF(lwp) THEN ! control print 712 WRITE(numout,*) 713 WRITE(numout,*)' sent fields to SAS component ' 714 WRITE(numout,*)' sea surface temperature (T before, Celcius) ' 715 WRITE(numout,*)' sea surface salinity ' 716 WRITE(numout,*)' surface currents U,V on local grid and spherical coordinates' 717 WRITE(numout,*)' sea surface height ' 718 WRITE(numout,*)' thickness of first ocean T level ' 719 WRITE(numout,*)' fraction of solar net radiation absorbed in the first ocean level' 720 WRITE(numout,*) 721 ENDIF 722 ENDIF 723 ! ! ------------------------------- ! 724 ! ! OPA-SAS coupling - snd by sas ! 725 ! ! ------------------------------- ! 726 ssnd(jps_sflx )%clname = 'I_SFLX' 727 ssnd(jps_fice2 )%clname = 'IIceFrc' 728 ssnd(jps_qsroce)%clname = 'I_QsrOce' 729 ssnd(jps_qnsoce)%clname = 'I_QnsOce' 730 ssnd(jps_oemp )%clname = 'IOEvaMPr' 731 ssnd(jps_otx1 )%clname = 'I_OTaux1' 732 ssnd(jps_oty1 )%clname = 'I_OTauy1' 733 ssnd(jps_rnf )%clname = 'I_Runoff' 734 ssnd(jps_taum )%clname = 'I_TauMod' 735 ! 736 IF( nn_components == jp_iam_sas ) THEN 737 IF( .NOT. ln_cpl ) ssnd(:)%laction = .FALSE. ! force default definition in case of opa <-> sas coupling 738 ssnd( (/jps_qsroce, jps_qnsoce, jps_oemp, jps_fice2, jps_sflx, jps_otx1, jps_oty1, jps_taum/) )%laction = .TRUE. 739 ! 740 ! Change first letter to couple with atmosphere if already coupled with sea_ice 741 ! this is nedeed as each variable name used in the namcouple must be unique: 742 ! for example O_SSTSST sent by OPA to SAS and therefore S_SSTSST sent by SAS to the Atmosphere 743 DO jn = 1, jpsnd 744 IF ( ssnd(jn)%clname(1:1) == "O" ) ssnd(jn)%clname = "S"//ssnd(jn)%clname(2:LEN(ssnd(jn)%clname)) 745 END DO 746 ! 747 IF(lwp) THEN ! control print 748 WRITE(numout,*) 749 IF( .NOT. ln_cpl ) THEN 750 WRITE(numout,*)' sent fields to OPA component ' 751 ELSE 752 WRITE(numout,*)' Additional sent fields to OPA component : ' 753 ENDIF 754 WRITE(numout,*)' ice cover ' 755 WRITE(numout,*)' oce only EMP ' 756 WRITE(numout,*)' salt flux ' 757 WRITE(numout,*)' mixed oce-ice solar flux ' 758 WRITE(numout,*)' mixed oce-ice non solar flux ' 759 WRITE(numout,*)' wind stress U,V components' 760 WRITE(numout,*)' wind stress module' 761 ENDIF 762 ENDIF 763 569 764 ! 570 765 ! ================================ ! … … 572 767 ! ================================ ! 573 768 574 CALL cpl_define(jprcv, jpsnd,nn_cplmodel) 769 CALL cpl_define(jprcv, jpsnd, nn_cplmodel) 770 575 771 IF (ln_usecplmask) THEN 576 772 xcplmask(:,:,:) = 0. … … 582 778 xcplmask(:,:,:) = 1. 583 779 ENDIF 584 ! 585 IF( ln_dm2dc .AND. ( cpl_freq( jpr_qsroce ) + cpl_freq( jpr_qsrmix ) /= 86400 ) ) & 780 xcplmask(:,:,0) = 1. - SUM( xcplmask(:,:,1:nn_cplmodel), dim = 3 ) 781 ! 782 ncpl_qsr_freq = cpl_freq( 'O_QsrOce' ) + cpl_freq( 'O_QsrMix' ) + cpl_freq( 'I_QsrOce' ) + cpl_freq( 'I_QsrMix' ) 783 IF( ln_dm2dc .AND. ln_cpl .AND. ncpl_qsr_freq /= 86400 ) & 586 784 & CALL ctl_stop( 'sbc_cpl_init: diurnal cycle reconstruction (ln_dm2dc) needs daily couping for solar radiation' ) 785 ncpl_qsr_freq = 86400 / ncpl_qsr_freq 587 786 588 787 CALL wrk_dealloc( jpi,jpj, zacs, zaos ) … … 638 837 !! emp upward mass flux [evap. - precip. (- runoffs) (- calving)] (ocean only case) 639 838 !!---------------------------------------------------------------------- 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?? 839 INTEGER, INTENT(in) :: kt ! ocean model time step index 840 INTEGER, INTENT(in) :: k_fsbc ! frequency of sbc (-> ice model) computation 841 INTEGER, INTENT(in) :: k_ice ! ice management in the sbc (=0/1/2/3) 842 843 !! 844 LOGICAL :: llnewtx, llnewtau ! update wind stress components and module?? 645 845 INTEGER :: ji, jj, jn ! dummy loop indices 646 846 INTEGER :: isec ! number of seconds since nit000 (assuming rdttra did not change since nit000) … … 650 850 REAL(wp) :: zcdrag = 1.5e-3 ! drag coefficient 651 851 REAL(wp) :: zzx, zzy ! temporary variables 652 REAL(wp), POINTER, DIMENSION(:,:) :: ztx, zty 852 REAL(wp), POINTER, DIMENSION(:,:) :: ztx, zty, zmsk, zemp, zqns, zqsr 653 853 !!---------------------------------------------------------------------- 654 854 ! 655 855 IF( nn_timing == 1 ) CALL timing_start('sbc_cpl_rcv') 656 856 ! 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) ) 857 CALL wrk_alloc( jpi,jpj, ztx, zty, zmsk, zemp, zqns, zqsr ) 858 ! 859 IF( ln_mixcpl ) zmsk(:,:) = 1. - xcplmask(:,:,0) 860 ! 861 ! ! ======================================================= ! 862 ! ! Receive all the atmos. fields (including ice information) 863 ! ! ======================================================= ! 864 isec = ( kt - nit000 ) * NINT( rdttra(1) ) ! date of exchanges 865 DO jn = 1, jprcv ! received fields sent by the atmosphere 866 IF( srcv(jn)%laction ) CALL cpl_rcv( jn, isec, frcv(jn)%z3, xcplmask(:,:,1:nn_cplmodel), nrcvinfo(jn) ) 662 867 END DO 663 868 … … 719 924 ! 720 925 ENDIF 721 722 926 ! ! ========================= ! 723 927 ! ! wind stress module ! (taum) … … 748 952 ENDIF 749 953 ENDIF 750 954 ! 751 955 ! ! ========================= ! 752 956 ! ! 10 m wind speed ! (wndm) … … 761 965 !CDIR NOVERRCHK 762 966 DO ji = 1, jpi 763 wndm(ji,jj) = SQRT( frcv(jpr_taum)%z3(ji,jj,1) * zcoef )967 frcv(jpr_w10m)%z3(ji,jj,1) = SQRT( frcv(jpr_taum)%z3(ji,jj,1) * zcoef ) 764 968 END DO 765 969 END DO 766 970 ENDIF 767 ELSE768 IF ( nrcvinfo(jpr_w10m) == OASIS_Rcv ) wndm(:,:) = frcv(jpr_w10m)%z3(:,:,1)769 971 ENDIF 770 972 … … 773 975 IF( MOD( kt-1, k_fsbc ) == 0 ) THEN 774 976 ! 775 utau(:,:) = frcv(jpr_otx1)%z3(:,:,1) 776 vtau(:,:) = frcv(jpr_oty1)%z3(:,:,1) 777 taum(:,:) = frcv(jpr_taum)%z3(:,:,1) 977 IF( ln_mixcpl ) THEN 978 utau(:,:) = utau(:,:) * xcplmask(:,:,0) + frcv(jpr_otx1)%z3(:,:,1) * zmsk(:,:) 979 vtau(:,:) = vtau(:,:) * xcplmask(:,:,0) + frcv(jpr_oty1)%z3(:,:,1) * zmsk(:,:) 980 taum(:,:) = taum(:,:) * xcplmask(:,:,0) + frcv(jpr_taum)%z3(:,:,1) * zmsk(:,:) 981 wndm(:,:) = wndm(:,:) * xcplmask(:,:,0) + frcv(jpr_w10m)%z3(:,:,1) * zmsk(:,:) 982 ELSE 983 utau(:,:) = frcv(jpr_otx1)%z3(:,:,1) 984 vtau(:,:) = frcv(jpr_oty1)%z3(:,:,1) 985 taum(:,:) = frcv(jpr_taum)%z3(:,:,1) 986 wndm(:,:) = frcv(jpr_w10m)%z3(:,:,1) 987 ENDIF 778 988 CALL iom_put( "taum_oce", taum ) ! output wind stress module 779 989 ! … … 781 991 782 992 #if defined key_cpl_carbon_cycle 783 ! ! atmosph. CO2 (ppm) 993 ! ! ================== ! 994 ! ! atmosph. CO2 (ppm) ! 995 ! ! ================== ! 784 996 IF( srcv(jpr_co2)%laction ) atm_co2(:,:) = frcv(jpr_co2)%z3(:,:,1) 785 997 #endif 786 998 999 ! Fields received by SAS when OASIS coupling 1000 ! (arrays no more filled at sbcssm stage) 1001 ! ! ================== ! 1002 ! ! SSS ! 1003 ! ! ================== ! 1004 IF( srcv(jpr_soce)%laction ) THEN ! received by sas in case of opa <-> sas coupling 1005 sss_m(:,:) = frcv(jpr_soce)%z3(:,:,1) 1006 CALL iom_put( 'sss_m', sss_m ) 1007 ENDIF 1008 ! 1009 ! ! ================== ! 1010 ! ! SST ! 1011 ! ! ================== ! 1012 IF( srcv(jpr_toce)%laction ) THEN ! received by sas in case of opa <-> sas coupling 1013 sst_m(:,:) = frcv(jpr_toce)%z3(:,:,1) 1014 IF( srcv(jpr_soce)%laction .AND. ln_useCT ) THEN ! make sure that sst_m is the potential temperature 1015 sst_m(:,:) = eos_pt_from_ct( sst_m(:,:), sss_m(:,:) ) 1016 ENDIF 1017 ENDIF 1018 ! ! ================== ! 1019 ! ! SSH ! 1020 ! ! ================== ! 1021 IF( srcv(jpr_ssh )%laction ) THEN ! received by sas in case of opa <-> sas coupling 1022 ssh_m(:,:) = frcv(jpr_ssh )%z3(:,:,1) 1023 CALL iom_put( 'ssh_m', ssh_m ) 1024 ENDIF 1025 ! ! ================== ! 1026 ! ! surface currents ! 1027 ! ! ================== ! 1028 IF( srcv(jpr_ocx1)%laction ) THEN ! received by sas in case of opa <-> sas coupling 1029 ssu_m(:,:) = frcv(jpr_ocx1)%z3(:,:,1) 1030 ub (:,:,1) = ssu_m(:,:) ! will be used in sbcice_lim in the call of lim_sbc_tau 1031 CALL iom_put( 'ssu_m', ssu_m ) 1032 ENDIF 1033 IF( srcv(jpr_ocy1)%laction ) THEN 1034 ssv_m(:,:) = frcv(jpr_ocy1)%z3(:,:,1) 1035 vb (:,:,1) = ssv_m(:,:) ! will be used in sbcice_lim in the call of lim_sbc_tau 1036 CALL iom_put( 'ssv_m', ssv_m ) 1037 ENDIF 1038 ! ! ======================== ! 1039 ! ! first T level thickness ! 1040 ! ! ======================== ! 1041 IF( srcv(jpr_e3t1st )%laction ) THEN ! received by sas in case of opa <-> sas coupling 1042 e3t_m(:,:) = frcv(jpr_e3t1st )%z3(:,:,1) 1043 CALL iom_put( 'e3t_m', e3t_m(:,:) ) 1044 ENDIF 1045 ! ! ================================ ! 1046 ! ! fraction of solar net radiation ! 1047 ! ! ================================ ! 1048 IF( srcv(jpr_fraqsr)%laction ) THEN ! received by sas in case of opa <-> sas coupling 1049 frq_m(:,:) = frcv(jpr_fraqsr)%z3(:,:,1) 1050 CALL iom_put( 'frq_m', frq_m ) 1051 ENDIF 1052 787 1053 ! ! ========================= ! 788 IF( k_ice <= 1 ) THEN! heat & freshwater fluxes ! (Ocean only case)1054 IF( k_ice <= 1 .AND. MOD( kt-1, k_fsbc ) == 0 ) THEN ! heat & freshwater fluxes ! (Ocean only case) 789 1055 ! ! ========================= ! 790 1056 ! 791 1057 ! ! 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 1058 IF( srcv(jpr_oemp)%laction .OR. srcv(jpr_rain)%laction ) THEN 1059 SELECT CASE( TRIM( sn_rcv_emp%cldes ) ) ! evaporation - precipitation 1060 CASE( 'conservative' ) 1061 zemp(:,:) = frcv(jpr_tevp)%z3(:,:,1) - ( frcv(jpr_rain)%z3(:,:,1) + frcv(jpr_snow)%z3(:,:,1) ) 1062 CASE( 'oce only', 'oce and ice' ) 1063 zemp(:,:) = frcv(jpr_oemp)%z3(:,:,1) 1064 CASE default 1065 CALL ctl_stop( 'sbc_cpl_rcv: wrong definition of sn_rcv_emp%cldes' ) 1066 END SELECT 1067 ELSE 1068 zemp(:,:) = 0._wp 1069 ENDIF 800 1070 ! 801 1071 ! ! 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 1072 IF( srcv(jpr_rnf)%laction ) rnf(:,:) = frcv(jpr_rnf)%z3(:,:,1) 1073 IF( srcv(jpr_cal)%laction ) zemp(:,:) = zemp(:,:) - frcv(jpr_cal)%z3(:,:,1) 1074 1075 IF( ln_mixcpl ) THEN ; emp(:,:) = emp(:,:) * xcplmask(:,:,0) + zemp(:,:) * zmsk(:,:) 1076 ELSE ; emp(:,:) = zemp(:,:) 1077 ENDIF 821 1078 ! 822 1079 ! ! 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) 1080 IF( srcv(jpr_qnsoce)%laction ) THEN ; zqns(:,:) = frcv(jpr_qnsoce)%z3(:,:,1) 1081 ELSE IF( srcv(jpr_qnsmix)%laction ) THEN ; zqns(:,:) = frcv(jpr_qnsmix)%z3(:,:,1) 1082 ELSE ; zqns(:,:) = 0._wp 1083 END IF 825 1084 ! 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 1085 IF( nn_components /= jp_iam_opa ) THEN 1086 zqns(:,:) = zqns(:,:) - zemp(:,:) * sst_m(:,:) * rcp ! remove heat content due to mass flux (assumed to be at SST) 1087 IF( srcv(jpr_snow )%laction ) THEN 1088 zqns(:,:) = zqns(:,:) - frcv(jpr_snow)%z3(:,:,1) * lfus ! energy for melting solid precipitation over the free ocean 1089 ENDIF 1090 ENDIF 1091 IF( ln_mixcpl ) THEN ; qns(:,:) = qns(:,:) * xcplmask(:,:,0) + zqns(:,:) * zmsk(:,:) 1092 ELSE ; qns(:,:) = zqns(:,:) 829 1093 ENDIF 830 1094 831 1095 ! ! 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 1096 IF ( srcv(jpr_qsroce)%laction ) THEN ; zqsr(:,:) = frcv(jpr_qsroce)%z3(:,:,1) 1097 ELSE IF( srcv(jpr_qsrmix)%laction ) then ; zqsr(:,:) = frcv(jpr_qsrmix)%z3(:,:,1) 1098 ELSE ; zqsr(:,:) = 0._wp 1099 ENDIF 1100 IF( ln_dm2dc .AND. ln_cpl ) zqsr(:,:) = sbc_dcy( zqsr ) ! modify qsr to include the diurnal cycle 1101 IF( ln_mixcpl ) THEN ; qsr(:,:) = qsr(:,:) * xcplmask(:,:,0) + zqsr(:,:) * zmsk(:,:) 1102 ELSE ; qsr(:,:) = zqsr(:,:) 1103 ENDIF 835 1104 ! 836 837 ENDIF 838 ! 839 CALL wrk_dealloc( jpi,jpj, ztx, zty ) 1105 ! salt flux over the ocean (received by opa in case of opa <-> sas coupling) 1106 IF( srcv(jpr_sflx )%laction ) sfx(:,:) = frcv(jpr_sflx )%z3(:,:,1) 1107 ! Ice cover (received by opa in case of opa <-> sas coupling) 1108 IF( srcv(jpr_fice )%laction ) fr_i(:,:) = frcv(jpr_fice )%z3(:,:,1) 1109 ! 1110 1111 ENDIF 1112 ! 1113 CALL wrk_dealloc( jpi,jpj, ztx, zty, zmsk, zemp, zqns, zqsr ) 840 1114 ! 841 1115 IF( nn_timing == 1 ) CALL timing_stop('sbc_cpl_rcv') … … 934 1208 ! 935 1209 ENDIF 936 937 1210 ! ! ======================= ! 938 1211 ! ! put on ice grid ! … … 1056 1329 1057 1330 1058 SUBROUTINE sbc_cpl_ice_flx( p_frld , palbi , psst , pist)1331 SUBROUTINE sbc_cpl_ice_flx( p_frld, palbi, psst, pist ) 1059 1332 !!---------------------------------------------------------------------- 1060 1333 !! *** ROUTINE sbc_cpl_ice_flx *** … … 1098 1371 REAL(wp), INTENT(in ), DIMENSION(:,:) :: p_frld ! lead fraction [0 to 1] 1099 1372 ! 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 1373 REAL(wp), INTENT(in ), DIMENSION(:,:,:), OPTIONAL :: palbi ! all skies ice albedo 1374 REAL(wp), INTENT(in ), DIMENSION(:,: ), OPTIONAL :: psst ! sea surface temperature [Celsius] 1375 REAL(wp), INTENT(in ), DIMENSION(:,:,:), OPTIONAL :: pist ! ice surface temperature [Kelvin] 1376 ! 1377 INTEGER :: jl ! dummy loop index 1378 REAL(wp), POINTER, DIMENSION(:,: ) :: zcptn, ztmp, zicefr, zmsk 1379 REAL(wp), POINTER, DIMENSION(:,: ) :: zemp_tot, zemp_ice, zsprecip, ztprecip, zqns_tot, zqsr_tot 1380 REAL(wp), POINTER, DIMENSION(:,:,:) :: zqns_ice, zqsr_ice, zdqns_ice 1381 REAL(wp), POINTER, DIMENSION(:,: ) :: zevap, zsnw, zqns_oce, zqsr_oce, zqprec_ice, zqemp_oce ! for LIM3 1106 1382 !!---------------------------------------------------------------------- 1107 1383 ! 1108 1384 IF( nn_timing == 1 ) CALL timing_start('sbc_cpl_ice_flx') 1109 1385 ! 1110 CALL wrk_alloc( jpi,jpj, zcptn, ztmp, zicefr ) 1111 1386 CALL wrk_alloc( jpi,jpj, zcptn, ztmp, zicefr, zmsk, zemp_tot, zemp_ice, zsprecip, ztprecip, zqns_tot, zqsr_tot ) 1387 CALL wrk_alloc( jpi,jpj,jpl, zqns_ice, zqsr_ice, zdqns_ice ) 1388 1389 IF( ln_mixcpl ) zmsk(:,:) = 1. - xcplmask(:,:,0) 1112 1390 zicefr(:,:) = 1.- p_frld(:,:) 1113 1391 zcptn(:,:) = rcp * sst_m(:,:) … … 1117 1395 ! ! ========================= ! 1118 1396 ! 1119 ! ! total Precipitations - total Evaporation (emp_tot) 1120 ! ! solid precipitation - sublimation (emp_ice) 1121 ! ! solid Precipitation (sprecip) 1397 ! ! total Precipitation - total Evaporation (emp_tot) 1398 ! ! solid precipitation - sublimation (emp_ice) 1399 ! ! solid Precipitation (sprecip) 1400 ! ! liquid + solid Precipitation (tprecip) 1122 1401 SELECT CASE( TRIM( sn_rcv_emp%cldes ) ) 1123 1402 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)1403 zsprecip(:,:) = frcv(jpr_snow)%z3(:,:,1) ! May need to ensure positive here 1404 ztprecip(:,:) = frcv(jpr_rain)%z3(:,:,1) + zsprecip(:,:) ! May need to ensure positive here 1405 zemp_tot(:,:) = frcv(jpr_tevp)%z3(:,:,1) - ztprecip(:,:) 1406 zemp_ice(:,:) = frcv(jpr_ievp)%z3(:,:,1) - frcv(jpr_snow)%z3(:,:,1) 1128 1407 CALL iom_put( 'rain' , frcv(jpr_rain)%z3(:,:,1) ) ! liquid precipitation 1129 1408 IF( iom_use('hflx_rain_cea') ) & … … 1136 1415 CALL iom_put( 'hflx_evap_cea', ztmp(:,:) * zcptn(:,:) ) ! heat flux from from evap (cell average) 1137 1416 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) 1417 zemp_tot(:,:) = p_frld(:,:) * frcv(jpr_oemp)%z3(:,:,1) + zicefr(:,:) * frcv(jpr_sbpr)%z3(:,:,1) 1418 zemp_ice(:,:) = frcv(jpr_semp)%z3(:,:,1) 1419 zsprecip(:,:) = frcv(jpr_ievp)%z3(:,:,1) - frcv(jpr_semp)%z3(:,:,1) 1420 ztprecip(:,:) = frcv(jpr_semp)%z3(:,:,1) - frcv(jpr_sbpr)%z3(:,:,1) + zsprecip(:,:) 1141 1421 END SELECT 1422 1423 IF( iom_use('subl_ai_cea') ) & 1424 CALL iom_put( 'subl_ai_cea', frcv(jpr_ievp)%z3(:,:,1) * zicefr(:,:) ) ! Sublimation over sea-ice (cell average) 1425 ! 1426 ! ! runoffs and calving (put in emp_tot) 1427 IF( srcv(jpr_rnf)%laction ) rnf(:,:) = frcv(jpr_rnf)%z3(:,:,1) 1428 IF( srcv(jpr_cal)%laction ) THEN 1429 zemp_tot(:,:) = zemp_tot(:,:) - frcv(jpr_cal)%z3(:,:,1) 1430 CALL iom_put( 'calving_cea', frcv(jpr_cal)%z3(:,:,1) ) 1431 ENDIF 1432 1433 IF( ln_mixcpl ) THEN 1434 emp_tot(:,:) = emp_tot(:,:) * xcplmask(:,:,0) + zemp_tot(:,:) * zmsk(:,:) 1435 emp_ice(:,:) = emp_ice(:,:) * xcplmask(:,:,0) + zemp_ice(:,:) * zmsk(:,:) 1436 sprecip(:,:) = sprecip(:,:) * xcplmask(:,:,0) + zsprecip(:,:) * zmsk(:,:) 1437 tprecip(:,:) = tprecip(:,:) * xcplmask(:,:,0) + ztprecip(:,:) * zmsk(:,:) 1438 ELSE 1439 emp_tot(:,:) = zemp_tot(:,:) 1440 emp_ice(:,:) = zemp_ice(:,:) 1441 sprecip(:,:) = zsprecip(:,:) 1442 tprecip(:,:) = ztprecip(:,:) 1443 ENDIF 1142 1444 1143 1445 CALL iom_put( 'snowpre' , sprecip ) ! Snow … … 1146 1448 IF( iom_use('snow_ai_cea') ) & 1147 1449 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', 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 1450 1178 1451 ! ! ========================= ! … … 1180 1453 ! ! ========================= ! 1181 1454 CASE( 'oce only' ) ! the required field is directly provided 1182 qns_tot(:,: ) = frcv(jpr_qnsoce)%z3(:,:,1)1455 zqns_tot(:,: ) = frcv(jpr_qnsoce)%z3(:,:,1) 1183 1456 CASE( 'conservative' ) ! the required fields are directly provided 1184 qns_tot(:,: ) = frcv(jpr_qnsmix)%z3(:,:,1)1457 zqns_tot(:,: ) = frcv(jpr_qnsmix)%z3(:,:,1) 1185 1458 IF ( TRIM(sn_rcv_qns%clcat) == 'yes' ) THEN 1186 qns_ice(:,:,1:jpl) = frcv(jpr_qnsice)%z3(:,:,1:jpl)1459 zqns_ice(:,:,1:jpl) = frcv(jpr_qnsice)%z3(:,:,1:jpl) 1187 1460 ELSE 1188 1461 ! Set all category values equal for the moment 1189 1462 DO jl=1,jpl 1190 qns_ice(:,:,jl) = frcv(jpr_qnsice)%z3(:,:,1)1463 zqns_ice(:,:,jl) = frcv(jpr_qnsice)%z3(:,:,1) 1191 1464 ENDDO 1192 1465 ENDIF 1193 1466 CASE( 'oce and ice' ) ! the total flux is computed from ocean and ice fluxes 1194 qns_tot(:,: ) = p_frld(:,:) * frcv(jpr_qnsoce)%z3(:,:,1)1467 zqns_tot(:,: ) = p_frld(:,:) * frcv(jpr_qnsoce)%z3(:,:,1) 1195 1468 IF ( TRIM(sn_rcv_qns%clcat) == 'yes' ) THEN 1196 1469 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)1470 zqns_tot(:,: ) = zqns_tot(:,:) + a_i(:,:,jl) * frcv(jpr_qnsice)%z3(:,:,jl) 1471 zqns_ice(:,:,jl) = frcv(jpr_qnsice)%z3(:,:,jl) 1199 1472 ENDDO 1200 1473 ELSE 1201 1474 qns_tot(:,: ) = qns_tot(:,:) + zicefr(:,:) * frcv(jpr_qnsice)%z3(:,:,1) 1202 1475 DO jl=1,jpl 1203 qns_ice(:,:,jl) = frcv(jpr_qnsice)%z3(:,:,1) 1476 zqns_tot(:,: ) = zqns_tot(:,:) + zicefr(:,:) * frcv(jpr_qnsice)%z3(:,:,1) 1477 zqns_ice(:,:,jl) = frcv(jpr_qnsice)%z3(:,:,1) 1204 1478 ENDDO 1205 1479 ENDIF 1206 1480 CASE( 'mixed oce-ice' ) ! the ice flux is cumputed from the total flux, the SST and ice informations 1207 1481 ! ** 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) &1482 zqns_tot(:,: ) = frcv(jpr_qnsmix)%z3(:,:,1) 1483 zqns_ice(:,:,1) = frcv(jpr_qnsmix)%z3(:,:,1) & 1210 1484 & + frcv(jpr_dqnsdt)%z3(:,:,1) * ( pist(:,:,1) - ( (rt0 + psst(:,: ) ) * p_frld(:,:) & 1211 1485 & + pist(:,:,1) * zicefr(:,:) ) ) 1212 1486 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 1487 !!gm 1221 !! currently it is taken into account in leads budget but not in the qns_tot, and thus not in1488 !! currently it is taken into account in leads budget but not in the zqns_tot, and thus not in 1222 1489 !! the flux that enter the ocean.... 1223 1490 !! moreover 1 - it is not diagnose anywhere.... … … 1228 1495 IF( srcv(jpr_cal)%laction ) THEN ! Iceberg melting 1229 1496 ztmp(:,:) = frcv(jpr_cal)%z3(:,:,1) * lfus ! add the latent heat of iceberg melting 1230 qns_tot(:,:) =qns_tot(:,:) - ztmp(:,:)1497 zqns_tot(:,:) = zqns_tot(:,:) - ztmp(:,:) 1231 1498 IF( iom_use('hflx_cal_cea') ) & 1232 1499 CALL iom_put( 'hflx_cal_cea', ztmp + frcv(jpr_cal)%z3(:,:,1) * zcptn(:,:) ) ! heat flux from calving 1233 1500 ENDIF 1501 1502 ztmp(:,:) = p_frld(:,:) * zsprecip(:,:) * lfus 1503 IF( iom_use('hflx_snow_cea') ) CALL iom_put( 'hflx_snow_cea', ztmp + sprecip(:,:) * zcptn(:,:) ) ! heat flux from snow (cell average) 1504 1505 #if defined key_lim3 1506 CALL wrk_alloc( jpi,jpj, zevap, zsnw, zqns_oce, zqprec_ice, zqemp_oce ) 1507 1508 ! --- evaporation --- ! 1509 ! clem: evap_ice is set to 0 for LIM3 since we still do not know what to do with sublimation 1510 ! the problem is: the atm. imposes both mass evaporation and heat removed from the snow/ice 1511 ! but it is incoherent WITH the ice model 1512 DO jl=1,jpl 1513 evap_ice(:,:,jl) = 0._wp ! should be: frcv(jpr_ievp)%z3(:,:,1) 1514 ENDDO 1515 zevap(:,:) = zemp_tot(:,:) + ztprecip(:,:) ! evaporation over ocean 1516 1517 ! --- evaporation minus precipitation --- ! 1518 emp_oce(:,:) = emp_tot(:,:) - emp_ice(:,:) 1519 1520 ! --- non solar flux over ocean --- ! 1521 ! note: p_frld cannot be = 0 since we limit the ice concentration to amax 1522 zqns_oce = 0._wp 1523 WHERE( p_frld /= 0._wp ) zqns_oce(:,:) = ( zqns_tot(:,:) - SUM( a_i * zqns_ice, dim=3 ) ) / p_frld(:,:) 1524 1525 ! --- heat flux associated with emp --- ! 1526 zsnw(:,:) = 0._wp 1527 CALL lim_thd_snwblow( p_frld, zsnw ) ! snow distribution over ice after wind blowing 1528 zqemp_oce(:,:) = - zevap(:,:) * p_frld(:,:) * zcptn(:,:) & ! evap 1529 & + ( ztprecip(:,:) - zsprecip(:,:) ) * zcptn(:,:) & ! liquid precip 1530 & + zsprecip(:,:) * ( 1._wp - zsnw ) * ( zcptn(:,:) - lfus ) ! solid precip over ocean 1531 qemp_ice(:,:) = - frcv(jpr_ievp)%z3(:,:,1) * zicefr(:,:) * zcptn(:,:) & ! ice evap 1532 & + zsprecip(:,:) * zsnw * ( zcptn(:,:) - lfus ) ! solid precip over ice 1533 1534 ! --- heat content of precip over ice in J/m3 (to be used in 1D-thermo) --- ! 1535 zqprec_ice(:,:) = rhosn * ( zcptn(:,:) - lfus ) 1536 1537 ! --- total non solar flux --- ! 1538 zqns_tot(:,:) = zqns_tot(:,:) + qemp_ice(:,:) + zqemp_oce(:,:) 1539 1540 ! --- in case both coupled/forced are active, we must mix values --- ! 1541 IF( ln_mixcpl ) THEN 1542 qns_tot(:,:) = qns_tot(:,:) * xcplmask(:,:,0) + zqns_tot(:,:)* zmsk(:,:) 1543 qns_oce(:,:) = qns_oce(:,:) * xcplmask(:,:,0) + zqns_oce(:,:)* zmsk(:,:) 1544 DO jl=1,jpl 1545 qns_ice(:,:,jl) = qns_ice(:,:,jl) * xcplmask(:,:,0) + zqns_ice(:,:,jl)* zmsk(:,:) 1546 ENDDO 1547 qprec_ice(:,:) = qprec_ice(:,:) * xcplmask(:,:,0) + zqprec_ice(:,:)* zmsk(:,:) 1548 qemp_oce (:,:) = qemp_oce(:,:) * xcplmask(:,:,0) + zqemp_oce(:,:)* zmsk(:,:) 1549 !!clem evap_ice(:,:) = evap_ice(:,:) * xcplmask(:,:,0) 1550 ELSE 1551 qns_tot (:,: ) = zqns_tot (:,: ) 1552 qns_oce (:,: ) = zqns_oce (:,: ) 1553 qns_ice (:,:,:) = zqns_ice (:,:,:) 1554 qprec_ice(:,:) = zqprec_ice(:,:) 1555 qemp_oce (:,:) = zqemp_oce (:,:) 1556 ENDIF 1557 1558 CALL wrk_dealloc( jpi,jpj, zevap, zsnw, zqns_oce, zqprec_ice, zqemp_oce ) 1559 #else 1560 1561 ! clem: this formulation is certainly wrong... but better than it was... 1562 zqns_tot(:,:) = zqns_tot(:,:) & ! zqns_tot update over free ocean with: 1563 & - ztmp(:,:) & ! remove the latent heat flux of solid precip. melting 1564 & - ( zemp_tot(:,:) & ! remove the heat content of mass flux (assumed to be at SST) 1565 & - zemp_ice(:,:) * zicefr(:,:) ) * zcptn(:,:) 1566 1567 IF( ln_mixcpl ) THEN 1568 qns_tot(:,:) = qns(:,:) * p_frld(:,:) + SUM( qns_ice(:,:,:) * a_i(:,:,:), dim=3 ) ! total flux from blk 1569 qns_tot(:,:) = qns_tot(:,:) * xcplmask(:,:,0) + zqns_tot(:,:)* zmsk(:,:) 1570 DO jl=1,jpl 1571 qns_ice(:,:,jl) = qns_ice(:,:,jl) * xcplmask(:,:,0) + zqns_ice(:,:,jl)* zmsk(:,:) 1572 ENDDO 1573 ELSE 1574 qns_tot(:,: ) = zqns_tot(:,: ) 1575 qns_ice(:,:,:) = zqns_ice(:,:,:) 1576 ENDIF 1577 1578 #endif 1234 1579 1235 1580 ! ! ========================= ! … … 1237 1582 ! ! ========================= ! 1238 1583 CASE( 'oce only' ) 1239 qsr_tot(:,: ) = MAX( 0._wp , frcv(jpr_qsroce)%z3(:,:,1) )1584 zqsr_tot(:,: ) = MAX( 0._wp , frcv(jpr_qsroce)%z3(:,:,1) ) 1240 1585 CASE( 'conservative' ) 1241 qsr_tot(:,: ) = frcv(jpr_qsrmix)%z3(:,:,1)1586 zqsr_tot(:,: ) = frcv(jpr_qsrmix)%z3(:,:,1) 1242 1587 IF ( TRIM(sn_rcv_qsr%clcat) == 'yes' ) THEN 1243 qsr_ice(:,:,1:jpl) = frcv(jpr_qsrice)%z3(:,:,1:jpl)1588 zqsr_ice(:,:,1:jpl) = frcv(jpr_qsrice)%z3(:,:,1:jpl) 1244 1589 ELSE 1245 1590 ! Set all category values equal for the moment 1246 1591 DO jl=1,jpl 1247 qsr_ice(:,:,jl) = frcv(jpr_qsrice)%z3(:,:,1)1592 zqsr_ice(:,:,jl) = frcv(jpr_qsrice)%z3(:,:,1) 1248 1593 ENDDO 1249 1594 ENDIF 1250 qsr_tot(:,: ) = frcv(jpr_qsrmix)%z3(:,:,1)1251 qsr_ice(:,:,1) = frcv(jpr_qsrice)%z3(:,:,1)1595 zqsr_tot(:,: ) = frcv(jpr_qsrmix)%z3(:,:,1) 1596 zqsr_ice(:,:,1) = frcv(jpr_qsrice)%z3(:,:,1) 1252 1597 CASE( 'oce and ice' ) 1253 qsr_tot(:,: ) = p_frld(:,:) * frcv(jpr_qsroce)%z3(:,:,1)1598 zqsr_tot(:,: ) = p_frld(:,:) * frcv(jpr_qsroce)%z3(:,:,1) 1254 1599 IF ( TRIM(sn_rcv_qsr%clcat) == 'yes' ) THEN 1255 1600 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)1601 zqsr_tot(:,: ) = zqsr_tot(:,:) + a_i(:,:,jl) * frcv(jpr_qsrice)%z3(:,:,jl) 1602 zqsr_ice(:,:,jl) = frcv(jpr_qsrice)%z3(:,:,jl) 1258 1603 ENDDO 1259 1604 ELSE 1260 1605 qsr_tot(:,: ) = qsr_tot(:,:) + zicefr(:,:) * frcv(jpr_qsrice)%z3(:,:,1) 1261 1606 DO jl=1,jpl 1262 qsr_ice(:,:,jl) = frcv(jpr_qsrice)%z3(:,:,1) 1607 zqsr_tot(:,: ) = zqsr_tot(:,:) + zicefr(:,:) * frcv(jpr_qsrice)%z3(:,:,1) 1608 zqsr_ice(:,:,jl) = frcv(jpr_qsrice)%z3(:,:,1) 1263 1609 ENDDO 1264 1610 ENDIF 1265 1611 CASE( 'mixed oce-ice' ) 1266 qsr_tot(:,: ) = frcv(jpr_qsrmix)%z3(:,:,1)1612 zqsr_tot(:,: ) = frcv(jpr_qsrmix)%z3(:,:,1) 1267 1613 ! ** NEED TO SORT OUT HOW THIS SHOULD WORK IN THE MULTI-CATEGORY CASE - CURRENTLY NOT ALLOWED WHEN INTERFACE INITIALISED ** 1268 1614 ! Create solar heat flux over ice using incoming solar heat flux and albedos 1269 1615 ! ( see OASIS3 user guide, 5th edition, p39 ) 1270 qsr_ice(:,:,1) = frcv(jpr_qsrmix)%z3(:,:,1) * ( 1.- palbi(:,:,1) ) &1616 zqsr_ice(:,:,1) = frcv(jpr_qsrmix)%z3(:,:,1) * ( 1.- palbi(:,:,1) ) & 1271 1617 & / ( 1.- ( albedo_oce_mix(:,: ) * p_frld(:,:) & 1272 1618 & + palbi (:,:,1) * zicefr(:,:) ) ) 1273 1619 END SELECT 1274 IF( ln_dm2dc ) THEN ! modify qsr to include the diurnal cycle1275 qsr_tot(:,: ) = sbc_dcy(qsr_tot(:,: ) )1620 IF( ln_dm2dc .AND. ln_cpl ) THEN ! modify qsr to include the diurnal cycle 1621 zqsr_tot(:,: ) = sbc_dcy( zqsr_tot(:,: ) ) 1276 1622 DO jl=1,jpl 1277 qsr_ice(:,:,jl) = sbc_dcy(qsr_ice(:,:,jl) )1623 zqsr_ice(:,:,jl) = sbc_dcy( zqsr_ice(:,:,jl) ) 1278 1624 ENDDO 1625 ENDIF 1626 1627 #if defined key_lim3 1628 CALL wrk_alloc( jpi,jpj, zqsr_oce ) 1629 ! --- solar flux over ocean --- ! 1630 ! note: p_frld cannot be = 0 since we limit the ice concentration to amax 1631 zqsr_oce = 0._wp 1632 WHERE( p_frld /= 0._wp ) zqsr_oce(:,:) = ( zqsr_tot(:,:) - SUM( a_i * zqsr_ice, dim=3 ) ) / p_frld(:,:) 1633 1634 IF( ln_mixcpl ) THEN ; qsr_oce(:,:) = qsr_oce(:,:) * xcplmask(:,:,0) + zqsr_oce(:,:)* zmsk(:,:) 1635 ELSE ; qsr_oce(:,:) = zqsr_oce(:,:) ; ENDIF 1636 1637 CALL wrk_dealloc( jpi,jpj, zqsr_oce ) 1638 #endif 1639 1640 IF( ln_mixcpl ) THEN 1641 qsr_tot(:,:) = qsr(:,:) * p_frld(:,:) + SUM( qsr_ice(:,:,:) * a_i(:,:,:), dim=3 ) ! total flux from blk 1642 qsr_tot(:,:) = qsr_tot(:,:) * xcplmask(:,:,0) + zqsr_tot(:,:)* zmsk(:,:) 1643 DO jl=1,jpl 1644 qsr_ice(:,:,jl) = qsr_ice(:,:,jl) * xcplmask(:,:,0) + zqsr_ice(:,:,jl)* zmsk(:,:) 1645 ENDDO 1646 ELSE 1647 qsr_tot(:,: ) = zqsr_tot(:,: ) 1648 qsr_ice(:,:,:) = zqsr_ice(:,:,:) 1279 1649 ENDIF 1280 1650 … … 1284 1654 CASE ('coupled') 1285 1655 IF ( TRIM(sn_rcv_dqnsdt%clcat) == 'yes' ) THEN 1286 dqns_ice(:,:,1:jpl) = frcv(jpr_dqnsdt)%z3(:,:,1:jpl)1656 zdqns_ice(:,:,1:jpl) = frcv(jpr_dqnsdt)%z3(:,:,1:jpl) 1287 1657 ELSE 1288 1658 ! Set all category values equal for the moment 1289 1659 DO jl=1,jpl 1290 dqns_ice(:,:,jl) = frcv(jpr_dqnsdt)%z3(:,:,1)1660 zdqns_ice(:,:,jl) = frcv(jpr_dqnsdt)%z3(:,:,1) 1291 1661 ENDDO 1292 1662 ENDIF 1293 1663 END SELECT 1294 1664 1665 IF( ln_mixcpl ) THEN 1666 DO jl=1,jpl 1667 dqns_ice(:,:,jl) = dqns_ice(:,:,jl) * xcplmask(:,:,0) + zdqns_ice(:,:,jl) * zmsk(:,:) 1668 ENDDO 1669 ELSE 1670 dqns_ice(:,:,:) = zdqns_ice(:,:,:) 1671 ENDIF 1672 1295 1673 ! ! ========================= ! 1296 1674 SELECT CASE( TRIM( sn_rcv_iceflx%cldes ) ) ! topmelt and botmelt ! … … 1308 1686 fr2_i0(:,:) = ( 0.82 * ( 1.0 - cldf_ice ) + 0.65 * cldf_ice ) 1309 1687 1310 CALL wrk_dealloc( jpi,jpj, zcptn, ztmp, zicefr ) 1688 CALL wrk_dealloc( jpi,jpj, zcptn, ztmp, zicefr, zmsk, zemp_tot, zemp_ice, zsprecip, ztprecip, zqns_tot, zqsr_tot ) 1689 CALL wrk_dealloc( jpi,jpj,jpl, zqns_ice, zqsr_ice, zdqns_ice ) 1311 1690 ! 1312 1691 IF( nn_timing == 1 ) CALL timing_stop('sbc_cpl_ice_flx') … … 1328 1707 INTEGER :: ji, jj, jl ! dummy loop indices 1329 1708 INTEGER :: isec, info ! local integer 1709 REAL(wp) :: zumax, zvmax 1330 1710 REAL(wp), POINTER, DIMENSION(:,:) :: zfr_l, ztmp1, ztmp2, zotx1, zoty1, zotz1, zitx1, zity1, zitz1 1331 1711 REAL(wp), POINTER, DIMENSION(:,:,:) :: ztmp3, ztmp4 … … 1344 1724 ! ! ------------------------- ! 1345 1725 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 1726 1727 IF ( nn_components == jp_iam_opa ) THEN 1728 ztmp1(:,:) = tsn(:,:,1,jp_tem) ! send temperature as it is (potential or conservative) -> use of ln_useCT on the received part 1729 ELSE 1730 ! we must send the surface potential temperature 1731 IF( ln_useCT ) THEN ; ztmp1(:,:) = eos_pt_from_ct( tsn(:,:,1,jp_tem), tsn(:,:,1,jp_sal) ) 1732 ELSE ; ztmp1(:,:) = tsn(:,:,1,jp_tem) 1733 ENDIF 1734 ! 1735 SELECT CASE( sn_snd_temp%cldes) 1736 CASE( 'oce only' ) ; ztmp1(:,:) = ztmp1(:,:) + rt0 1737 CASE( 'oce and ice' ) ; ztmp1(:,:) = ztmp1(:,:) + rt0 1738 SELECT CASE( sn_snd_temp%clcat ) 1739 CASE( 'yes' ) 1740 ztmp3(:,:,1:jpl) = tn_ice(:,:,1:jpl) 1741 CASE( 'no' ) 1742 WHERE( SUM( a_i, dim=3 ) /= 0. ) 1743 ztmp3(:,:,1) = SUM( tn_ice * a_i, dim=3 ) / SUM( a_i, dim=3 ) 1744 ELSEWHERE 1745 ztmp3(:,:,1) = rt0 ! TODO: Is freezing point a good default? (Maybe SST is better?) 1746 END WHERE 1747 CASE default ; CALL ctl_stop( 'sbc_cpl_snd: wrong definition of sn_snd_temp%clcat' ) 1748 END SELECT 1749 CASE( 'weighted oce and ice' ) ; ztmp1(:,:) = ( ztmp1(:,:) + rt0 ) * zfr_l(:,:) 1750 SELECT CASE( sn_snd_temp%clcat ) 1751 CASE( 'yes' ) 1752 ztmp3(:,:,1:jpl) = tn_ice(:,:,1:jpl) * a_i(:,:,1:jpl) 1753 CASE( 'no' ) 1754 ztmp3(:,:,:) = 0.0 1755 DO jl=1,jpl 1756 ztmp3(:,:,1) = ztmp3(:,:,1) + tn_ice(:,:,jl) * a_i(:,:,jl) 1757 ENDDO 1758 CASE default ; CALL ctl_stop( 'sbc_cpl_snd: wrong definition of sn_snd_temp%clcat' ) 1759 END SELECT 1760 CASE( 'mixed oce-ice' ) 1761 ztmp1(:,:) = ( ztmp1(:,:) + rt0 ) * zfr_l(:,:) 1354 1762 DO jl=1,jpl 1355 ztmp 3(:,:,1) = ztmp3(:,:,1) + tn_ice(:,:,jl) * a_i(:,:,jl)1763 ztmp1(:,:) = ztmp1(:,:) + tn_ice(:,:,jl) * a_i(:,:,jl) 1356 1764 ENDDO 1357 CASE default ; CALL ctl_stop( 'sbc_cpl_snd: wrong definition of sn_snd_temp%clcat' )1765 CASE default ; CALL ctl_stop( 'sbc_cpl_snd: wrong definition of sn_snd_temp%cldes' ) 1358 1766 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 1767 ENDIF 1366 1768 IF( ssnd(jps_toce)%laction ) CALL cpl_snd( jps_toce, isec, RESHAPE ( ztmp1, (/jpi,jpj,1/) ), info ) 1367 1769 IF( ssnd(jps_tice)%laction ) CALL cpl_snd( jps_tice, isec, ztmp3, info ) … … 1372 1774 ! ! ------------------------- ! 1373 1775 IF( ssnd(jps_albice)%laction ) THEN ! ice 1374 ztmp3(:,:,1:jpl) = alb_ice(:,:,1:jpl) * a_i(:,:,1:jpl) 1776 SELECT CASE( sn_snd_alb%cldes ) 1777 CASE( 'ice' ) ; ztmp3(:,:,1:jpl) = alb_ice(:,:,1:jpl) 1778 CASE( 'weighted ice' ) ; ztmp3(:,:,1:jpl) = alb_ice(:,:,1:jpl) * a_i(:,:,1:jpl) 1779 CASE default ; CALL ctl_stop( 'sbc_cpl_snd: wrong definition of sn_snd_alb%cldes' ) 1780 END SELECT 1375 1781 CALL cpl_snd( jps_albice, isec, ztmp3, info ) 1376 1782 ENDIF … … 1385 1791 ! ! Ice fraction & Thickness ! 1386 1792 ! ! ------------------------- ! 1387 ! Send ice fraction field 1793 ! Send ice fraction field to atmosphere 1388 1794 IF( ssnd(jps_fice)%laction ) THEN 1389 1795 SELECT CASE( sn_snd_thick%clcat ) … … 1392 1798 CASE default ; CALL ctl_stop( 'sbc_cpl_snd: wrong definition of sn_snd_thick%clcat' ) 1393 1799 END SELECT 1394 CALL cpl_snd( jps_fice, isec, ztmp3, info ) 1800 IF( ssnd(jps_fice)%laction ) CALL cpl_snd( jps_fice, isec, ztmp3, info ) 1801 ENDIF 1802 1803 ! Send ice fraction field to OPA (sent by SAS in SAS-OPA coupling) 1804 IF( ssnd(jps_fice2)%laction ) THEN 1805 ztmp3(:,:,1) = fr_i(:,:) 1806 IF( ssnd(jps_fice2)%laction ) CALL cpl_snd( jps_fice2, isec, ztmp3, info ) 1395 1807 ENDIF 1396 1808 … … 1413 1825 END SELECT 1414 1826 CASE( 'ice and snow' ) 1415 ztmp3(:,:,1:jpl) = ht_i(:,:,1:jpl) 1416 ztmp4(:,:,1:jpl) = ht_s(:,:,1:jpl) 1827 SELECT CASE( sn_snd_thick%clcat ) 1828 CASE( 'yes' ) 1829 ztmp3(:,:,1:jpl) = ht_i(:,:,1:jpl) 1830 ztmp4(:,:,1:jpl) = ht_s(:,:,1:jpl) 1831 CASE( 'no' ) 1832 WHERE( SUM( a_i, dim=3 ) /= 0. ) 1833 ztmp3(:,:,1) = SUM( ht_i * a_i, dim=3 ) / SUM( a_i, dim=3 ) 1834 ztmp4(:,:,1) = SUM( ht_s * a_i, dim=3 ) / SUM( a_i, dim=3 ) 1835 ELSEWHERE 1836 ztmp3(:,:,1) = 0. 1837 ztmp4(:,:,1) = 0. 1838 END WHERE 1839 CASE default ; CALL ctl_stop( 'sbc_cpl_snd: wrong definition of sn_snd_thick%clcat' ) 1840 END SELECT 1417 1841 CASE default ; CALL ctl_stop( 'sbc_cpl_snd: wrong definition of sn_snd_thick%cldes' ) 1418 1842 END SELECT … … 1440 1864 ! i-1 i i 1441 1865 ! 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 1866 IF( nn_components == jp_iam_opa ) THEN 1867 zotx1(:,:) = un(:,:,1) 1868 zoty1(:,:) = vn(:,:,1) 1869 ELSE 1870 SELECT CASE( TRIM( sn_snd_crt%cldes ) ) 1871 CASE( 'oce only' ) ! C-grid ==> T 1453 1872 DO jj = 2, jpjm1 1454 1873 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) 1874 zotx1(ji,jj) = 0.5 * ( un(ji,jj,1) + un(ji-1,jj ,1) ) 1875 zoty1(ji,jj) = 0.5 * ( vn(ji,jj,1) + vn(ji ,jj-1,1) ) 1459 1876 END DO 1460 1877 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) 1878 CASE( 'weighted oce and ice' ) 1879 SELECT CASE ( cp_ice_msh ) 1880 CASE( 'C' ) ! Ocean and Ice on C-grid ==> T 1881 DO jj = 2, jpjm1 1882 DO ji = fs_2, fs_jpim1 ! vector opt. 1883 zotx1(ji,jj) = 0.5 * ( un (ji,jj,1) + un (ji-1,jj ,1) ) * zfr_l(ji,jj) 1884 zoty1(ji,jj) = 0.5 * ( vn (ji,jj,1) + vn (ji ,jj-1,1) ) * zfr_l(ji,jj) 1885 zitx1(ji,jj) = 0.5 * ( u_ice(ji,jj ) + u_ice(ji-1,jj ) ) * fr_i(ji,jj) 1886 zity1(ji,jj) = 0.5 * ( v_ice(ji,jj ) + v_ice(ji ,jj-1 ) ) * fr_i(ji,jj) 1887 END DO 1470 1888 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)1889 CASE( 'I' ) ! Ocean on C grid, Ice on I-point (B-grid) ==> T 1890 DO jj = 2, jpjm1 1891 DO ji = 2, jpim1 ! NO vector opt. 1892 zotx1(ji,jj) = 0.5 * ( un(ji,jj,1) + un(ji-1,jj ,1) ) * zfr_l(ji,jj) 1893 zoty1(ji,jj) = 0.5 * ( vn(ji,jj,1) + vn(ji ,jj-1,1) ) * zfr_l(ji,jj) 1894 zitx1(ji,jj) = 0.25 * ( u_ice(ji+1,jj+1) + u_ice(ji,jj+1) & 1895 & + u_ice(ji+1,jj ) + u_ice(ji,jj ) ) * fr_i(ji,jj) 1896 zity1(ji,jj) = 0.25 * ( v_ice(ji+1,jj+1) + v_ice(ji,jj+1) & 1897 & + v_ice(ji+1,jj ) + v_ice(ji,jj ) ) * fr_i(ji,jj) 1898 END DO 1481 1899 END DO 1482 END DO 1900 CASE( 'F' ) ! Ocean on C grid, Ice on F-point (B-grid) ==> T 1901 DO jj = 2, jpjm1 1902 DO ji = 2, jpim1 ! NO vector opt. 1903 zotx1(ji,jj) = 0.5 * ( un(ji,jj,1) + un(ji-1,jj ,1) ) * zfr_l(ji,jj) 1904 zoty1(ji,jj) = 0.5 * ( vn(ji,jj,1) + vn(ji ,jj-1,1) ) * zfr_l(ji,jj) 1905 zitx1(ji,jj) = 0.25 * ( u_ice(ji-1,jj-1) + u_ice(ji,jj-1) & 1906 & + u_ice(ji-1,jj ) + u_ice(ji,jj ) ) * fr_i(ji,jj) 1907 zity1(ji,jj) = 0.25 * ( v_ice(ji-1,jj-1) + v_ice(ji,jj-1) & 1908 & + v_ice(ji-1,jj ) + v_ice(ji,jj ) ) * fr_i(ji,jj) 1909 END DO 1910 END DO 1911 END SELECT 1912 CALL lbc_lnk( zitx1, 'T', -1. ) ; CALL lbc_lnk( zity1, 'T', -1. ) 1913 CASE( 'mixed oce-ice' ) 1914 SELECT CASE ( cp_ice_msh ) 1915 CASE( 'C' ) ! Ocean and Ice on C-grid ==> T 1916 DO jj = 2, jpjm1 1917 DO ji = fs_2, fs_jpim1 ! vector opt. 1918 zotx1(ji,jj) = 0.5 * ( un (ji,jj,1) + un (ji-1,jj ,1) ) * zfr_l(ji,jj) & 1919 & + 0.5 * ( u_ice(ji,jj ) + u_ice(ji-1,jj ) ) * fr_i(ji,jj) 1920 zoty1(ji,jj) = 0.5 * ( vn (ji,jj,1) + vn (ji ,jj-1,1) ) * zfr_l(ji,jj) & 1921 & + 0.5 * ( v_ice(ji,jj ) + v_ice(ji ,jj-1 ) ) * fr_i(ji,jj) 1922 END DO 1923 END DO 1924 CASE( 'I' ) ! Ocean on C grid, Ice on I-point (B-grid) ==> T 1925 DO jj = 2, jpjm1 1926 DO ji = 2, jpim1 ! NO vector opt. 1927 zotx1(ji,jj) = 0.5 * ( un(ji,jj,1) + un(ji-1,jj ,1) ) * zfr_l(ji,jj) & 1928 & + 0.25 * ( u_ice(ji+1,jj+1) + u_ice(ji,jj+1) & 1929 & + u_ice(ji+1,jj ) + u_ice(ji,jj ) ) * fr_i(ji,jj) 1930 zoty1(ji,jj) = 0.5 * ( vn(ji,jj,1) + vn(ji ,jj-1,1) ) * zfr_l(ji,jj) & 1931 & + 0.25 * ( v_ice(ji+1,jj+1) + v_ice(ji,jj+1) & 1932 & + v_ice(ji+1,jj ) + v_ice(ji,jj ) ) * fr_i(ji,jj) 1933 END DO 1934 END DO 1935 CASE( 'F' ) ! Ocean on C grid, Ice on F-point (B-grid) ==> T 1936 DO jj = 2, jpjm1 1937 DO ji = 2, jpim1 ! NO vector opt. 1938 zotx1(ji,jj) = 0.5 * ( un(ji,jj,1) + un(ji-1,jj ,1) ) * zfr_l(ji,jj) & 1939 & + 0.25 * ( u_ice(ji-1,jj-1) + u_ice(ji,jj-1) & 1940 & + u_ice(ji-1,jj ) + u_ice(ji,jj ) ) * fr_i(ji,jj) 1941 zoty1(ji,jj) = 0.5 * ( vn(ji,jj,1) + vn(ji ,jj-1,1) ) * zfr_l(ji,jj) & 1942 & + 0.25 * ( v_ice(ji-1,jj-1) + v_ice(ji,jj-1) & 1943 & + v_ice(ji-1,jj ) + v_ice(ji,jj ) ) * fr_i(ji,jj) 1944 END DO 1945 END DO 1946 END SELECT 1483 1947 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. ) 1948 CALL lbc_lnk( zotx1, ssnd(jps_ocx1)%clgrid, -1. ) ; CALL lbc_lnk( zoty1, ssnd(jps_ocy1)%clgrid, -1. ) 1949 ! 1950 ENDIF 1521 1951 ! 1522 1952 ! … … 1558 1988 ENDIF 1559 1989 ! 1990 ! 1991 ! Fields sent by OPA to SAS when doing OPA<->SAS coupling 1992 ! ! SSH 1993 IF( ssnd(jps_ssh )%laction ) THEN 1994 ! ! removed inverse barometer ssh when Patm 1995 ! forcing is used (for sea-ice dynamics) 1996 IF( ln_apr_dyn ) THEN ; ztmp1(:,:) = sshb(:,:) - 0.5 * ( ssh_ib(:,:) + ssh_ibb(:,:) ) 1997 ELSE ; ztmp1(:,:) = sshn(:,:) 1998 ENDIF 1999 CALL cpl_snd( jps_ssh , isec, RESHAPE ( ztmp1 , (/jpi,jpj,1/) ), info ) 2000 2001 ENDIF 2002 ! ! SSS 2003 IF( ssnd(jps_soce )%laction ) THEN 2004 CALL cpl_snd( jps_soce , isec, RESHAPE ( tsn(:,:,1,jp_sal), (/jpi,jpj,1/) ), info ) 2005 ENDIF 2006 ! ! first T level thickness 2007 IF( ssnd(jps_e3t1st )%laction ) THEN 2008 CALL cpl_snd( jps_e3t1st, isec, RESHAPE ( fse3t_n(:,:,1) , (/jpi,jpj,1/) ), info ) 2009 ENDIF 2010 ! ! Qsr fraction 2011 IF( ssnd(jps_fraqsr)%laction ) THEN 2012 CALL cpl_snd( jps_fraqsr, isec, RESHAPE ( fraqsr_1lev(:,:) , (/jpi,jpj,1/) ), info ) 2013 ENDIF 2014 ! 2015 ! Fields sent by SAS to OPA when OASIS coupling 2016 ! ! Solar heat flux 2017 IF( ssnd(jps_qsroce)%laction ) CALL cpl_snd( jps_qsroce, isec, RESHAPE ( qsr , (/jpi,jpj,1/) ), info ) 2018 IF( ssnd(jps_qnsoce)%laction ) CALL cpl_snd( jps_qnsoce, isec, RESHAPE ( qns , (/jpi,jpj,1/) ), info ) 2019 IF( ssnd(jps_oemp )%laction ) CALL cpl_snd( jps_oemp , isec, RESHAPE ( emp , (/jpi,jpj,1/) ), info ) 2020 IF( ssnd(jps_sflx )%laction ) CALL cpl_snd( jps_sflx , isec, RESHAPE ( sfx , (/jpi,jpj,1/) ), info ) 2021 IF( ssnd(jps_otx1 )%laction ) CALL cpl_snd( jps_otx1 , isec, RESHAPE ( utau, (/jpi,jpj,1/) ), info ) 2022 IF( ssnd(jps_oty1 )%laction ) CALL cpl_snd( jps_oty1 , isec, RESHAPE ( vtau, (/jpi,jpj,1/) ), info ) 2023 IF( ssnd(jps_rnf )%laction ) CALL cpl_snd( jps_rnf , isec, RESHAPE ( rnf , (/jpi,jpj,1/) ), info ) 2024 IF( ssnd(jps_taum )%laction ) CALL cpl_snd( jps_taum , isec, RESHAPE ( taum, (/jpi,jpj,1/) ), info ) 2025 1560 2026 CALL wrk_dealloc( jpi,jpj, zfr_l, ztmp1, ztmp2, zotx1, zoty1, zotz1, zitx1, zity1, zitz1 ) 1561 2027 CALL wrk_dealloc( jpi,jpj,jpl, ztmp3, ztmp4 )
Note: See TracChangeset
for help on using the changeset viewer.