source: trunk/LMDZ.PLUTO/libf/phypluto/spreadglacier_simple.F90 @ 3421

Last change on this file since 3421 was 3412, checked in by afalco, 16 months ago

Pluto PCM: Imported glaciers & conserv mass routines from pluto.old.
AF

File size: 4.1 KB
Line 
1      subroutine spreadglacier_simple(ngrid,nq,pqsurf,timstep)
2         
3      use geometry_mod, only: latitude, longitude, cell_area
4      use radinc_h, only : naerkind
5      use surfdat_h, only: phisfi
6      use tracer_h, only: igcm_n2
7
8      implicit none
9
10!==================================================================
11!     Purpose
12!     -------
13!     Spreading the glacier in a circular crater using a speed of flow
14
15!     Inputs
16!     ------
17!     ngrid                 Number of vertical columns
18!     pqsurf(ngrid,nq)      N2 ice tracer on surface kg/m2
19!     
20!     Outputs
21!     -------
22!     pqsurf(ngrid,nq)      new value for N2 ice tracer on surface kg/m2
23!     
24!     Authors
25!     -------
26!     Tanguy Bertrand (2015)
27!     
28!==================================================================
29
30#include "dimensions.h"
31
32!-----------------------------------------------------------------------
33!     Arguments
34
35      INTEGER ngrid, nq, ig, ig2
36      REAL :: pqsurf(ngrid,nq)
37      REAL timstep
38!-----------------------------------------------------------------------
39!     Local variables
40
41      REAL :: zdqsurf(ngrid)
42      REAL phiref
43      REAL dist
44      REAL locdist
45      REAL tau
46      REAL speed
47      REAL dist_pluto2  ! function
48      LOGICAL firstcall
49      DATA firstcall/.true./
50      INTEGER indglac(ngrid)
51      INTEGER nglac
52
53!---------------- INPUT ------------------------------------------------
54! Defined for  Tombaugh crater :
55
56      dist=300.  ! reference radius km
57      phiref=-1000. ! geopotential reference
58!      tau=1.e10
59      speed=6.4e-7   ! glacier speed km/s
60!      tau=1. ! direct redistribution glacier
61
62!--------------- FIRSTCALL : define crater index -----------------------
63      IF (firstcall) then
64         indglac(:)=0
65         DO ig=1,ngrid
66            IF (phisfi(ig).lt.phiref) THEN
67               nglac=nglac+1
68               indglac(nglac)=ig
69            ENDIF
70         ENDDO
71         !write(*,*) 'inglac=',indglac
72         write(*,*) 'nglac=',nglac
73 
74         firstcall=.false.
75
76      ENDIF
77
78!--------------- redistribution -------------------------------------
79      DO ig=1,nglac
80        ! redistribution if large amount of ice only : detection of reservoirs
81        IF (pqsurf(indglac(ig),igcm_n2).gt.50.) THEN
82         ! looking where to redistribute : if less ice, lower altitude
83         ! and close to the reservoir
84         DO ig2=1,nglac
85           IF ((pqsurf(indglac(ig2),igcm_n2).lt.pqsurf(indglac(ig),igcm_n2)*0.8)) THEN   
86             !write(*,*) 'loop=',latitude(indglac(ig2))*180/3.14,longitude(indglac(ig2))*180/3.14
87             !.and.(phisfi(indglac(ig2)).le.phisfi(indglac(ig)))) THEN
88             locdist=dist_pluto2(latitude(indglac(ig2)),longitude(indglac(ig2)),latitude(indglac(ig)),longitude(indglac(ig)))
89             !write(*,*) 'dist=',locdist
90             IF (locdist.lt.dist) THEN
91              !write(*,*) 'calcul=',latitude(indglac(ig2))*180/3.14,longitude(indglac(ig2))*180/3.14,latitude(indglac(ig))*180/3.14,longitude(indglac(ig))*180/3.14
92              tau=locdist/speed  ! characteristic time for redistribution
93              zdqsurf(indglac(ig2))=pqsurf(indglac(ig),igcm_n2)*(1.-exp(-timstep/tau))/timstep*cell_area(indglac(ig))/cell_area(indglac(ig2))
94              zdqsurf(indglac(ig))=pqsurf(indglac(ig),igcm_n2)*(exp(-timstep/tau)-1.)/timstep
95              pqsurf(indglac(ig2),igcm_n2)=pqsurf(indglac(ig2),igcm_n2)+zdqsurf(indglac(ig2))*timstep
96              pqsurf(indglac(ig),igcm_n2)=pqsurf(indglac(ig),igcm_n2)+zdqsurf(indglac(ig))*timstep
97             ENDIF
98           ENDIF
99         ENDDO
100        ENDIF
101      ENDDO
102
103      end subroutine spreadglacier_simple
104
105      real function dist_pluto2(lat1,lon1,lat2,lon2)
106      implicit none
107      !#include "comcstfi.h"
108      real, intent(in) ::  lat1,lon1,lat2,lon2
109      real dlon,dlat
110      real a,c,d
111      real rad
112      rad=1187. ! planetary radius (km) 1187 km sur Pluton
113
114      dlon = lon2 - lon1
115      dlat = lat2 - lat1
116      a = (sin(dlat/2))**2 + cos(lat1) * cos(lat2) * (sin(dlon/2))**2
117      c = 2 * atan2( sqrt(a), sqrt(1-a) )
118      dist_pluto2 = rad * c
119      return
120      end function dist_pluto2
121
122
Note: See TracBrowser for help on using the repository browser.