Ignore:
Timestamp:
Jul 23, 2024, 8:22:55 AM (2 months ago)
Author:
abarral
Message:

Handle DEBUG_IO in lmdz_cppkeys_wrapper.F90
Transform some files .F -> .[fF]90
[ne compile pas à cause de writefield_u non défini - en attente de réponse Laurent]

File:
1 moved

Legend:

Unmodified
Added
Removed
  • LMDZ6/branches/Amaury_dev/libf/dyn3dmem/caldyn_loc.f90

    r5100 r5101  
    1 
    21! $Id: $
    32
    4 #undef DEBUG_IO
    5 !#define DEBUG_IO
     3SUBROUTINE caldyn_loc &
     4        (itau, ucov, vcov, teta, ps, masse, pk, pkf, phis, &
     5        phi, conser, du, dv, dteta, dp, w, pbaru, pbarv, time)
     6  USE parallel_lmdz
     7  USE Write_Field_loc
     8  USE caldyn_mod, ONLY: vcont, ucont, ang, p, massebx, masseby, &
     9          vorpot, ecin, bern, massebxy, convm
     10  USE comvert_mod, ONLY: ap, bp
     11  USE lmdz_cppkeys_wrapper, ONLY: CPPKEY_DEBUGIO
    612
    7       SUBROUTINE caldyn_loc
    8      $ (itau,ucov,vcov,teta,ps,masse,pk,pkf,phis ,
    9      $  phi,conser,du,dv,dteta,dp,w,pbaru,pbarv,time )
    10       USE parallel_lmdz
    11       USE Write_Field_loc
    12       USE caldyn_mod, ONLY: vcont, ucont, ang, p, massebx, masseby,
    13      &                      vorpot, ecin, bern, massebxy, convm
    14       USE comvert_mod, ONLY: ap, bp
    15      
    16       IMPLICIT NONE
     13  IMPLICIT NONE
    1714
    18 !=======================================================================
     15  !=======================================================================
    1916
    20 !  Auteur :  P. Le Van
     17  !  Auteur :  P. Le Van
    2118
    22 !   Objet:
    23 !   ------
     19  !   Objet:
     20  !   ------
    2421
    25 !   Calcul des tendances dynamiques.
     22  !   Calcul des tendances dynamiques.
    2623
    27 ! Modif 04/93 F.Forget
    28 !=======================================================================
     24  ! Modif 04/93 F.Forget
     25  !=======================================================================
    2926
    30 !-----------------------------------------------------------------------
    31 !   0. Declarations:
    32 !   ----------------
     27  !-----------------------------------------------------------------------
     28  !   0. Declarations:
     29  !   ----------------
    3330
    34       include "dimensions.h"
    35       include "paramet.h"
    36       include "comgeom.h"
     31  include "dimensions.h"
     32  include "paramet.h"
     33  include "comgeom.h"
    3734
    38 !   Arguments:
    39 !   ----------
     35  !   Arguments:
     36  !   ----------
    4037
    41       LOGICAL,INTENT(IN) :: conser ! triggers printing some diagnostics ! not used
    42       INTEGER,INTENT(IN) :: itau ! time step index ! not used
    43       REAL,INTENT(IN) :: vcov(ijb_v:ije_v,llm) ! covariant meridional wind
    44       REAL,INTENT(IN) :: ucov(ijb_u:ije_u,llm) ! covariant zonal wind
    45       REAL,INTENT(IN) :: teta(ijb_u:ije_u,llm) ! potential temperature
    46       REAL,INTENT(IN) :: ps(ijb_u:ije_u) ! surface pressure
    47       REAL,INTENT(IN) :: phis(ijb_u:ije_u) ! geopotential at the surface
    48       REAL,INTENT(IN) :: pk(iip1,jjb_u:jje_u,llm) ! Exner at mid-layer
    49       REAL,INTENT(IN) :: pkf(ijb_u:ije_u,llm) ! filtered Exner
    50       REAL,INTENT(IN) :: phi(ijb_u:ije_u,llm) ! geopotential
    51       REAL,INTENT(OUT) :: masse(ijb_u:ije_u,llm) ! air mass
    52       REAL,INTENT(OUT) :: dv(ijb_v:ije_v,llm) ! tendency on vcov
    53       REAL,INTENT(OUT) :: du(ijb_u:ije_u,llm) ! tendency on ucov
    54       REAL,INTENT(OUT) :: dteta(ijb_u:ije_u,llm) ! tenddency on teta
    55       REAL,INTENT(OUT) :: dp(ijb_u:ije_u) ! tendency on ps
    56       REAL,INTENT(OUT) :: w(ijb_u:ije_u,llm) ! vertical velocity
    57       REAL,INTENT(OUT) :: pbaru(ijb_u:ije_u,llm) ! mass flux in the zonal direction
    58       REAL,INTENT(OUT) :: pbarv(ijb_v:ije_v,llm) ! mass flux in the meridional direction
    59       REAL,INTENT(IN) :: time ! current time
     38  LOGICAL, INTENT(IN) :: conser ! triggers printing some diagnostics ! not used
     39  INTEGER, INTENT(IN) :: itau ! time step index ! not used
     40  REAL, INTENT(IN) :: vcov(ijb_v:ije_v, llm) ! covariant meridional wind
     41  REAL, INTENT(IN) :: ucov(ijb_u:ije_u, llm) ! covariant zonal wind
     42  REAL, INTENT(IN) :: teta(ijb_u:ije_u, llm) ! potential temperature
     43  REAL, INTENT(IN) :: ps(ijb_u:ije_u) ! surface pressure
     44  REAL, INTENT(IN) :: phis(ijb_u:ije_u) ! geopotential at the surface
     45  REAL, INTENT(IN) :: pk(iip1, jjb_u:jje_u, llm) ! Exner at mid-layer
     46  REAL, INTENT(IN) :: pkf(ijb_u:ije_u, llm) ! filtered Exner
     47  REAL, INTENT(IN) :: phi(ijb_u:ije_u, llm) ! geopotential
     48  REAL, INTENT(OUT) :: masse(ijb_u:ije_u, llm) ! air mass
     49  REAL, INTENT(OUT) :: dv(ijb_v:ije_v, llm) ! tendency on vcov
     50  REAL, INTENT(OUT) :: du(ijb_u:ije_u, llm) ! tendency on ucov
     51  REAL, INTENT(OUT) :: dteta(ijb_u:ije_u, llm) ! tenddency on teta
     52  REAL, INTENT(OUT) :: dp(ijb_u:ije_u) ! tendency on ps
     53  REAL, INTENT(OUT) :: w(ijb_u:ije_u, llm) ! vertical velocity
     54  REAL, INTENT(OUT) :: pbaru(ijb_u:ije_u, llm) ! mass flux in the zonal direction
     55  REAL, INTENT(OUT) :: pbarv(ijb_v:ije_v, llm) ! mass flux in the meridional direction
     56  REAL, INTENT(IN) :: time ! current time
    6057
    61 !   Local:
    62 !   ------
     58  !   Local:
     59  !   ------
    6360
    64       INTEGER   ij,l,ijb,ije,ierr
     61  INTEGER :: ij, l, ijb, ije, ierr
    6562
    6663
    67 !-----------------------------------------------------------------------
    68 !   Compute dynamical tendencies:
    69 !--------------------------------
     64  !-----------------------------------------------------------------------
     65  !   Compute dynamical tendencies:
     66  !--------------------------------
    7067
    71       ! compute contravariant winds ucont() and vcont
    72       CALL covcont_loc  ( llm    , ucov    , vcov , ucont, vcont     )
    73       ! compute pressure p()
    74       CALL pression_loc ( ip1jmp1, ap      , bp   ,  ps  , p         )
    75 cym      CALL psextbar (   ps   , psexbarxy                          )
    76 c$OMP BARRIER
    77       ! compute mass in each atmospheric mesh: masse()
    78       CALL massdair_loc (    p   , masse                             )
    79       ! compute X and Y-averages of mass, massebx() and masseby()
    80       CALL massbar_loc  (   masse, massebx , masseby                 )
    81       ! compute XY-average of mass, massebxy()
    82       call massbarxy_loc(   masse, massebxy                          )
    83       ! compute mass fluxes pbaru() and pbarv()
    84       CALL flumass_loc  ( massebx, masseby,vcont,ucont,pbaru,pbarv   )
    85       ! compute dteta() , horizontal converging flux of theta
    86       CALL dteta1_loc   (   teta , pbaru   , pbarv, dteta            )
    87       ! compute convm(), horizontal converging flux of mass
    88       CALL convmas1_loc  (   pbaru, pbarv   , convm                  )
    89 c$OMP BARRIER     
    90       CALL convmas2_loc  (   convm                      )
    91 c$OMP BARRIER
    92 #ifdef DEBUG_IO
    93       call WriteField_u('ucont',ucont)
    94       call WriteField_v('vcont',vcont)
    95       call WriteField_u('p',p)
    96       call WriteField_u('masse',masse)
    97       call WriteField_u('massebx',massebx)
    98       call WriteField_v('masseby',masseby)
    99       call WriteField_v('massebxy',massebxy)
    100       call WriteField_u('pbaru',pbaru)
    101       call WriteField_v('pbarv',pbarv)
    102       call WriteField_u('dteta',dteta)
    103       call WriteField_u('convm',convm)
    104 #endif     
     68  ! ! compute contravariant winds ucont() and vcont
     69  CALL covcont_loc  (llm, ucov, vcov, ucont, vcont)
     70  ! ! compute pressure p()
     71  CALL pression_loc (ip1jmp1, ap, bp, ps, p)
     72  !ym      CALL psextbar (   ps   , psexbarxy                          )
     73  !$OMP BARRIER
     74  ! ! compute mass in each atmospheric mesh: masse()
     75  CALL massdair_loc (p, masse)
     76  ! ! compute X and Y-averages of mass, massebx() and masseby()
     77  CALL massbar_loc  (masse, massebx, masseby)
     78  ! ! compute XY-average of mass, massebxy()
     79  CALL massbarxy_loc(masse, massebxy)
     80  ! ! compute mass fluxes pbaru() and pbarv()
     81  CALL flumass_loc  (massebx, masseby, vcont, ucont, pbaru, pbarv)
     82  ! ! compute dteta() , horizontal converging flux of theta
     83  CALL dteta1_loc   (teta, pbaru, pbarv, dteta)
     84  ! ! compute convm(), horizontal converging flux of mass
     85  CALL convmas1_loc  (pbaru, pbarv, convm)
     86  !$OMP BARRIER
     87  CALL convmas2_loc  (convm)
     88  !$OMP BARRIER
     89  IF (CPPKEY_DEBUGIO) THEN
     90    CALL WriteField_u('ucont', ucont)
     91    CALL WriteField_v('vcont', vcont)
     92    CALL WriteField_u('p', p)
     93    CALL WriteField_u('masse', masse)
     94    CALL WriteField_u('massebx', massebx)
     95    CALL WriteField_v('masseby', masseby)
     96    CALL WriteField_v('massebxy', massebxy)
     97    CALL WriteField_u('pbaru', pbaru)
     98    CALL WriteField_v('pbarv', pbarv)
     99    CALL WriteField_u('dteta', dteta)
     100    CALL WriteField_u('convm', convm)
     101  END IF
    105102
    106 c$OMP BARRIER
    107 c$OMP MASTER
    108       ijb=ij_begin
    109       ije=ij_end
    110       ! compute pressure variation due to mass convergence
    111       DO ij =ijb, ije
    112          dp( ij ) = convm( ij,1 ) / airesurg( ij )
    113       ENDDO
    114 c$OMP END MASTER
    115 c$OMP BARRIER
    116      
    117       ! compute vertical velocity w()
    118       CALL vitvert_loc ( convm  , w                                )
    119       ! compute potential vorticity vorpot()
    120       CALL tourpot_loc ( vcov   , ucov  , massebxy  , vorpot       )
    121       ! compute rotation induced du() and dv()
    122       CALL dudv1_loc   ( vorpot , pbaru , pbarv     , du     , dv  )
     103  !$OMP BARRIER
     104  !$OMP MASTER
     105  ijb = ij_begin
     106  ije = ij_end
     107  ! ! compute pressure variation due to mass convergence
     108  DO ij = ijb, ije
     109    dp(ij) = convm(ij, 1) / airesurg(ij)
     110  ENDDO
     111  !$OMP END MASTER
     112  !$OMP BARRIER
    123113
    124 #ifdef DEBUG_IO     
    125       call WriteField_u('w',w)
    126       call WriteField_v('vorpot',vorpot)
    127       call WriteField_u('du',du)
    128       call WriteField_v('dv',dv)
    129 #endif     
    130      
    131       ! compute kinetic energy ecin()
    132       CALL enercin_loc ( vcov   , ucov  , vcont   , ucont  , ecin  )
    133       ! compute Bernouilli function bern()
    134       CALL bernoui_loc ( ip1jmp1, llm   , phi       , ecin   , bern)
    135       ! compute and add du() and dv() contributions from Bernouilli and pressure
    136       CALL dudv2_loc   ( teta   , pkf   , bern      , du     , dv  )
     114  ! ! compute vertical velocity w()
     115  CALL vitvert_loc (convm, w)
     116  ! ! compute potential vorticity vorpot()
     117  CALL tourpot_loc (vcov, ucov, massebxy, vorpot)
     118  ! ! compute rotation induced du() and dv()
     119  CALL dudv1_loc   (vorpot, pbaru, pbarv, du, dv)
    137120
    138 #ifdef DEBUG_IO
    139       call WriteField_u('ecin',ecin)
    140       call WriteField_u('bern',bern)
    141       call WriteField_u('du',du)
    142       call WriteField_v('dv',dv)
    143       call WriteField_u('pkf',pkf)
    144 #endif
    145      
    146       ijb=ij_begin-iip1
    147       ije=ij_end+iip1
    148      
    149       if (pole_nord) ijb=ij_begin
    150       if (pole_sud) ije=ij_end
     121  IF (CPPKEY_DEBUGIO) THEN
     122    CALL WriteField_u('w', w)
     123    CALL WriteField_v('vorpot', vorpot)
     124    CALL WriteField_u('du', du)
     125    CALL WriteField_v('dv', dv)
     126  END IF
    151127
    152 c$OMP DO SCHEDULE(STATIC,OMP_CHUNK)       
    153       DO l=1,llm
    154          DO ij=ijb,ije
    155             ang(ij,l) = ucov(ij,l) + constang(ij)
    156         ENDDO
    157       ENDDO
    158 c$OMP END DO
     128  ! ! compute kinetic energy ecin()
     129  CALL enercin_loc (vcov, ucov, vcont, ucont, ecin)
     130  ! ! compute Bernouilli function bern()
     131  CALL bernoui_loc (ip1jmp1, llm, phi, ecin, bern)
     132  ! ! compute and add du() and dv() contributions from Bernouilli and pressure
     133  CALL dudv2_loc   (teta, pkf, bern, du, dv)
    159134
    160       ! compute vertical advection contributions to du(), dv() and dteta()
    161       CALL advect_new_loc(ang,vcov,teta,w,massebx,masseby,du,dv,dteta)
     135  IF (CPPKEY_DEBUGIO) THEN
     136    CALL WriteField_u('ecin', ecin)
     137    CALL WriteField_u('bern', bern)
     138    CALL WriteField_u('du', du)
     139    CALL WriteField_v('dv', dv)
     140    CALL WriteField_u('pkf', pkf)
     141  END IF
    162142
    163 C  WARNING probleme de peridocite de dv sur les PC/linux. Pb d'arrondi
    164 C          probablement. Observe sur le code compile avec pgf90 3.0-1
    165       ijb=ij_begin
    166       ije=ij_end
    167       if (pole_sud) ije=ij_end-iip1
     143  ijb = ij_begin - iip1
     144  ije = ij_end + iip1
    168145
    169 c$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
    170       DO l = 1, llm
    171          DO ij = ijb, ije, iip1
    172            IF( dv(ij,l)/=dv(ij+iim,l) )  THEN
    173 c         PRINT *,'!!!ATTENTION!!! probleme de periodicite sur vcov', 
    174 c    ,   ' dans caldyn'
    175 c         PRINT *,' l,  ij = ', l, ij, ij+iim,dv(ij+iim,l),dv(ij,l)
    176           dv(ij+iim,l) = dv(ij,l)
    177           endif
    178          enddo
    179       enddo
    180 c$OMP END DO NOWAIT     
     146  if (pole_nord) ijb = ij_begin
     147  if (pole_sud) ije = ij_end
    181148
    182 ! Ehouarn: NB: output of control variables not implemented...
     149  !$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
     150  DO l = 1, llm
     151    DO ij = ijb, ije
     152      ang(ij, l) = ucov(ij, l) + constang(ij)
     153    ENDDO
     154  ENDDO
     155  !$OMP END DO
    183156
    184       RETURN
    185       END
     157  ! ! compute vertical advection contributions to du(), dv() and dteta()
     158  CALL advect_new_loc(ang, vcov, teta, w, massebx, masseby, du, dv, dteta)
     159
     160  !  WARNING probleme de peridocite de dv sur les PC/linux. Pb d'arrondi
     161  ! probablement. Observe sur le code compile avec pgf90 3.0-1
     162  ijb = ij_begin
     163  ije = ij_end
     164  if (pole_sud) ije = ij_end - iip1
     165
     166  !$OMP DO SCHEDULE(STATIC,OMP_CHUNK)
     167  DO l = 1, llm
     168    DO ij = ijb, ije, iip1
     169      IF(dv(ij, l)/=dv(ij + iim, l))  THEN
     170        ! PRINT *,'!!!ATTENTION!!! probleme de periodicite sur vcov',
     171        !    ,   ' dans caldyn'
     172        ! PRINT *,' l,  ij = ', l, ij, ij+iim,dv(ij+iim,l),dv(ij,l)
     173        dv(ij + iim, l) = dv(ij, l)
     174      endif
     175    enddo
     176  enddo
     177  !$OMP END DO NOWAIT
     178
     179  ! Ehouarn: NB: output of control variables not implemented...
     180
     181  RETURN
     182END SUBROUTINE caldyn_loc
Note: See TracChangeset for help on using the changeset viewer.