Ignore:
Timestamp:
Oct 24, 2024, 1:10:47 PM (4 weeks ago)
Author:
afalco
Message:

Pluto PCM: fixed glacier simple initializations bug with gfortran.
AF

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.PLUTO/libf/phypluto/spreadglacier_simple.F90

    r3412 r3473  
    11      subroutine spreadglacier_simple(ngrid,nq,pqsurf,timstep)
    2          
     2
    33      use geometry_mod, only: latitude, longitude, cell_area
    44      use radinc_h, only : naerkind
     
    1212!     -------
    1313!     Spreading the glacier in a circular crater using a speed of flow
    14 ! 
     14!
    1515!     Inputs
    16 !     ------ 
     16!     ------
    1717!     ngrid                 Number of vertical columns
    18 !     pqsurf(ngrid,nq)      N2 ice tracer on surface kg/m2 
    19 !     
     18!     pqsurf(ngrid,nq)      N2 ice tracer on surface kg/m2
     19!
    2020!     Outputs
    2121!     -------
    22 !     pqsurf(ngrid,nq)      new value for N2 ice tracer on surface kg/m2 
    23 !     
     22!     pqsurf(ngrid,nq)      new value for N2 ice tracer on surface kg/m2
     23!
    2424!     Authors
    25 !     ------- 
     25!     -------
    2626!     Tanguy Bertrand (2015)
    27 !     
     27!
    2828!==================================================================
    2929
     
    3535      INTEGER ngrid, nq, ig, ig2
    3636      REAL :: pqsurf(ngrid,nq)
    37       REAL timstep 
     37      REAL timstep
    3838!-----------------------------------------------------------------------
    3939!     Local variables
     
    4848      LOGICAL firstcall
    4949      DATA firstcall/.true./
    50       INTEGER indglac(ngrid)
    51       INTEGER nglac
     50      INTEGER,SAVE,ALLOCATABLE :: indglac(:)
     51      INTEGER, SAVE :: nglac
    5252
    5353!---------------- INPUT ------------------------------------------------
    54 ! Defined for  Tombaugh crater : 
     54! Defined for  Tombaugh crater :
    5555
    5656      dist=300.  ! reference radius km
     
    6262!--------------- FIRSTCALL : define crater index -----------------------
    6363      IF (firstcall) then
     64         ALLOCATE(indglac(ngrid))
     65         nglac=0
    6466         indglac(:)=0
    6567         DO ig=1,ngrid
     
    6870               indglac(nglac)=ig
    6971            ENDIF
    70          ENDDO
    71          !write(*,*) 'inglac=',indglac
    72          write(*,*) 'nglac=',nglac
    73  
    74          firstcall=.false.
     72         ENDDO
     73         firstcall=.false.
    7574
    7675      ENDIF
    77 
    7876!--------------- redistribution -------------------------------------
    7977      DO ig=1,nglac
     
    8381         ! and close to the reservoir
    8482         DO ig2=1,nglac
    85            IF ((pqsurf(indglac(ig2),igcm_n2).lt.pqsurf(indglac(ig),igcm_n2)*0.8)) THEN   
     83           IF ((pqsurf(indglac(ig2),igcm_n2).lt.pqsurf(indglac(ig),igcm_n2)*0.8)) THEN
    8684             !write(*,*) 'loop=',latitude(indglac(ig2))*180/3.14,longitude(indglac(ig2))*180/3.14
    8785             !.and.(phisfi(indglac(ig2)).le.phisfi(indglac(ig)))) THEN
     
    9896           ENDIF
    9997         ENDDO
    100         ENDIF 
     98        ENDIF
    10199      ENDDO
    102100
     
    117115      c = 2 * atan2( sqrt(a), sqrt(1-a) )
    118116      dist_pluto2 = rad * c
    119       return 
     117      return
    120118      end function dist_pluto2
    121119
Note: See TracChangeset for help on using the changeset viewer.