source: LMDZ6/branches/Amaury_dev/libf/dyn3d/groupe.F90 @ 5218

Last change on this file since 5218 was 5159, checked in by abarral, 5 months ago

Put dimensions.h and paramet.h into modules

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 2.9 KB
RevLine 
[524]1! $Header$
[5099]2
[5158]3MODULE lmdz_groupe
[5159]4  USE lmdz_paramet
[5158]5  IMPLICIT NONE; PRIVATE
6  PUBLIC groupe
[524]7
[5158]8CONTAINS
[524]9
[5158]10  SUBROUTINE groupe(pbaru, pbarv, pbarum, pbarvm, wm)
[524]11
[5158]12    USE comconst_mod, ONLY: ngroup
13    USE lmdz_ssum_scopy, ONLY: scopy
14    USE lmdz_comgeom2
[524]15
[5159]16  USE lmdz_dimensions, ONLY: iim, jjm, llm, ndm
[5158]17    IMPLICIT NONE
[524]18
[5158]19    !   sous-programme servant a fitlrer les champs de flux de masse aux
20    !   poles en "regroupant" les mailles 2 par 2 puis 4 par 4 etc. au fur
21    !   et a mesure qu'on se rapproche du pole.
[5159]22
[5158]23    !   en entree: pbaru et pbarv
[5159]24
[5158]25    !   en sortie:  pbarum,pbarvm et wm.
[5159]26
[5158]27    !   remarque, le wm est recalcule a partir des pbaru pbarv et on n'a donc
28    !   pas besoin de w en entree.
[524]29
30
[5159]31
32
[5158]33    ! integer ngroup
34    ! parameter (ngroup=3)
[524]35
[5158]36    REAL :: pbaru(iip1, jjp1, llm), pbarv(iip1, jjm, llm)
[524]37
[5158]38    REAL :: pbarum(iip1, jjp1, llm), pbarvm(iip1, jjm, llm)
39    REAL :: wm(iip1, jjp1, llm)
[524]40
[5158]41    REAL :: zconvm(iip1, jjp1, llm), zconvmm(iip1, jjp1, llm)
42    ! (gfortran) Warning: Array ‘zconvm’ at (1) is larger than limit set by ‘-fmax-stack-var-size=’, moved from stack to static storage.
43    ! This makes the procedure unsafe when called recursively, or concurrently from multiple threads.
[524]44
[5158]45    REAL :: uu
[524]46
[5158]47    INTEGER :: i, j, l
[1674]48
[5158]49    LOGICAL :: firstcall, groupe_ok
50    save firstcall, groupe_ok
[524]51
[5158]52    data firstcall/.TRUE./
53    data groupe_ok/.TRUE./
[1674]54
[5158]55    IF (iim==1) THEN
56      groupe_ok = .FALSE.
57    ENDIF
[524]58
[5158]59    IF (firstcall) THEN
60      IF (groupe_ok) THEN
61        IF(mod(iim, 2**ngroup)/=0) &
62                CALL abort_gcm('groupe', 'probleme du nombre de point', 1)
63      endif
64      firstcall = .FALSE.
65    ENDIF
[524]66
67
[5158]68    !   Champs 1D
[524]69
[5158]70    CALL convflu(pbaru, pbarv, llm, zconvm)
[5103]71
[5158]72    CALL scopy(ijp1llm, zconvm, 1, zconvmm, 1)
73    CALL scopy(ijmllm, pbarv, 1, pbarvm, 1)
74
75    IF (groupe_ok) THEN
76      CALL groupeun(jjp1, llm, zconvmm)
77      CALL groupeun(jjm, llm, pbarvm)
78
79      !   Champs 3D
80      DO l = 1, llm
81        DO j = 2, jjm
82          uu = pbaru(iim, j, l)
83          DO i = 1, iim
84            uu = uu + pbarvm(i, j, l) - pbarvm(i, j - 1, l) - zconvmm(i, j, l)
85            pbarum(i, j, l) = uu
86            ! zconvm(i,j,l ) =  xflu(i-1,j,l)-xflu(i,j,l)+
87            !    *                      yflu(i,j,l)-yflu(i,j-1,l)
88          enddo
89          pbarum(iip1, j, l) = pbarum(1, j, l)
[5103]90        enddo
[524]91      enddo
92
[5158]93    else
94      pbarum(:, :, :) = pbaru(:, :, :)
95      pbarvm(:, :, :) = pbarv(:, :, :)
96    ENDIF
[1674]97
[5158]98    !    integration de la convergence de masse de haut  en bas ......
99    DO l = 1, llm
100      DO j = 1, jjp1
101        DO i = 1, iip1
102          zconvmm(i, j, l) = zconvmm(i, j, l)
103        enddo
[524]104      enddo
[5103]105    enddo
[5158]106    DO  l = llm - 1, 1, -1
107      DO j = 1, jjp1
108        DO i = 1, iip1
109          zconvmm(i, j, l) = zconvmm(i, j, l) + zconvmm(i, j, l + 1)
110        enddo
[524]111      enddo
[5103]112    enddo
[524]113
[5158]114    CALL vitvert(zconvmm, wm)
[524]115
[5158]116  END SUBROUTINE  groupe
[524]117
[5158]118END MODULE lmdz_groupe
Note: See TracBrowser for help on using the repository browser.