Changeset 3473 for trunk/LMDZ.PLUTO/libf/phypluto
- Timestamp:
- Oct 24, 2024, 1:10:47 PM (4 weeks ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.PLUTO/libf/phypluto/spreadglacier_simple.F90
r3412 r3473 1 1 subroutine spreadglacier_simple(ngrid,nq,pqsurf,timstep) 2 2 3 3 use geometry_mod, only: latitude, longitude, cell_area 4 4 use radinc_h, only : naerkind … … 12 12 ! ------- 13 13 ! Spreading the glacier in a circular crater using a speed of flow 14 ! 14 ! 15 15 ! Inputs 16 ! ------ 16 ! ------ 17 17 ! 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 ! 20 20 ! Outputs 21 21 ! ------- 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 ! 24 24 ! Authors 25 ! ------- 25 ! ------- 26 26 ! Tanguy Bertrand (2015) 27 ! 27 ! 28 28 !================================================================== 29 29 … … 35 35 INTEGER ngrid, nq, ig, ig2 36 36 REAL :: pqsurf(ngrid,nq) 37 REAL timstep 37 REAL timstep 38 38 !----------------------------------------------------------------------- 39 39 ! Local variables … … 48 48 LOGICAL firstcall 49 49 DATA firstcall/.true./ 50 INTEGER indglac(ngrid)51 INTEGER nglac50 INTEGER,SAVE,ALLOCATABLE :: indglac(:) 51 INTEGER, SAVE :: nglac 52 52 53 53 !---------------- INPUT ------------------------------------------------ 54 ! Defined for Tombaugh crater : 54 ! Defined for Tombaugh crater : 55 55 56 56 dist=300. ! reference radius km … … 62 62 !--------------- FIRSTCALL : define crater index ----------------------- 63 63 IF (firstcall) then 64 ALLOCATE(indglac(ngrid)) 65 nglac=0 64 66 indglac(:)=0 65 67 DO ig=1,ngrid … … 68 70 indglac(nglac)=ig 69 71 ENDIF 70 ENDDO 71 !write(*,*) 'inglac=',indglac 72 write(*,*) 'nglac=',nglac 73 74 firstcall=.false. 72 ENDDO 73 firstcall=.false. 75 74 76 75 ENDIF 77 78 76 !--------------- redistribution ------------------------------------- 79 77 DO ig=1,nglac … … 83 81 ! and close to the reservoir 84 82 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 86 84 !write(*,*) 'loop=',latitude(indglac(ig2))*180/3.14,longitude(indglac(ig2))*180/3.14 87 85 !.and.(phisfi(indglac(ig2)).le.phisfi(indglac(ig)))) THEN … … 98 96 ENDIF 99 97 ENDDO 100 ENDIF 98 ENDIF 101 99 ENDDO 102 100 … … 117 115 c = 2 * atan2( sqrt(a), sqrt(1-a) ) 118 116 dist_pluto2 = rad * c 119 return 117 return 120 118 end function dist_pluto2 121 119
Note: See TracChangeset
for help on using the changeset viewer.