source: branches/iLoveclim/SOURCES/Hemin40_files/output_hemin40_mod.f90 @ 146

Last change on this file since 146 was 146, checked in by aquiquet, 7 years ago

Grisli-iLoveclim branch: merged to trunk at revision 145

File size: 20.3 KB
Line 
1module  output_hemin40_mod
2
3       USE module3D_phy
4       use bilan_eau_mod
5
6implicit none
7
8real ::  bmean                        !
9real ::  accmean                      ! accumulation moyenne
10real ::  ablmean                      ! ablation moyenne
11real ::  calvmean                     ! moyenne calving
12real ::  ablbordmean                  !
13real ::  ablatotmean
14real ::  bmeltmean                    ! moyenne fusion basale
15real ::  tbmean                       ! temperature basale moyenne
16real ::  tbdotmean                    ! moyenne variation / temps temperature basale
17real ::  vsmean                       ! vitesse de surface moyenne
18!real ::  vsdotmean                    ! moyenne variation / temps vitesse de surface
19real ::  uzsmean   !!!! utilise ?     ! vitesse verticale de surface moyenne
20real ::  uzsdotmean                   ! moyenne variation / temps vitesse verticale de surface
21real ::  uzkmean                      ! moyenne vitesse verticale de surface
22real ::  hdotmean                     ! moyenne derivee / temps de H
23real ::  bdotmean                     ! moyenne bedrock derive / temps
24!real ::  pf0mean                      ! moyenne de PF0
25!real ::  pf1mean                      ! moyenne de PF1
26!real ::  evapmean                     ! moyenne de EVAP
27!real ::  pwmean                       ! moyenne PW
28!real ::  pdfmean                      ! moyenne PDF
29logical,dimension(nx,ny,13) :: mask_cal ! masque regions calotte
30! afq for CONSEAU REAL, dimension(nx,ny) :: old_H_dtt      ! Epaisseur de glace au pas de temps precedent
31integer, dimension(nx,ny) :: write_mask
32real :: ilc_units                     !afq, unit change from Grisli to iLoveclim
33
34CONTAINS
35!_________________________________________________________________________
36subroutine init_outshort
37
38!ndisp sorite courte tous les ndisp
39NDISP=100
40mask_cal(:,:,:)=.false.
41mask_cal(:,:,1)=.true.
42! afq for CONSEAU old_H_dtt(:,:)=H(:,:) ! afq missing in mab?
43! calcul d'un masque pour les regions des calottes
44do j=1,ny
45   do i=1,nx
46! -- Laurentide
47      IF( (((xlong(i,j).ge.190).AND.(xlong(i,j).lt.210)).AND. &
48           ((ylat(i,j).ge.50).AND.(ylat(i,j).le.70))).OR.      &
49           (((xlong(i,j).ge.210).AND.(xlong(i,j).lt.220)).AND. &
50           ((ylat(i,j).ge.55).AND.(ylat(i,j).le.75))).OR.      &
51           (((xlong(i,j).ge.220).AND.(xlong(i,j).lt.250)).AND. & 
52           ((ylat(i,j).ge.40).AND.(ylat(i,j).le.85))).OR.      &
53           (((xlong(i,j).ge.250).AND.(xlong(i,j).lt.290)).AND. &
54           ((ylat(i,j).ge.35).AND.(ylat(i,j).le.85))).OR.      &
55           (((xlong(i,j).ge.290).AND.(xlong(i,j).lt.300)).AND. &
56           ((ylat(i,j).ge.35).AND.(ylat(i,j).le.75))).OR.      &
57           (((xlong(i,j).ge.300).AND.(xlong(i,j).le.310)).AND. &
58           ((ylat(i,j).ge.35).AND.(ylat(i,j).lt.60))) ) mask_cal(i,j,2)=.true.
59! -- Groenland
60
61! -- Labrador Sector
62        IF((((xlong(i,j).ge.285).AND.(xlong(i,j).le.300)).AND.  &
63            ((ylat(i,j).ge.35).AND.(ylat(i,j).le.70))).OR.      &
64            (((xlong(i,j).ge.300).AND.(xlong(i,j).lt.310)).AND. &
65            ((ylat(i,j).ge.35).AND.(ylat(i,j).lt.60)))) mask_cal(i,j,5)=.true.
66! -- Keewatin Sector
67        IF(  ((xlong(i,j).gt.240).AND.(xlong(i,j).lt.285)).AND. &
68             ((ylat(i,j).ge.35).AND.(ylat(i,j).le.70)) ) mask_cal(i,j,6)=.true.
69! -- Innuitian Ice Sheet
70        IF( ((xlong(i,j).gt.230).AND.(xlong(i,j).lt.290)).AND. &
71            ((ylat(i,j).gt.70).AND.(ylat(i,j).le.85))) mask_cal(i,j,7)=.true.
72! -- Cordilleran Ice Sheet
73        IF( (((xlong(i,j).ge.190).AND.(xlong(i,j).lt.210)).AND. &
74            ((ylat(i,j).ge.50).AND.(ylat(i,j).le.70))).OR.      &
75            (((xlong(i,j).ge.210).AND.(xlong(i,j).lt.220)).AND. &
76            ((ylat(i,j).ge.55).AND.(ylat(i,j).le.75))).OR.      &
77            (((xlong(i,j).ge.210).AND.(xlong(i,j).lt.220)).AND. &
78             ((ylat(i,j).ge.55).AND.(ylat(i,j).le.75))).OR.     &
79            (((xlong(i,j).ge.220).AND.(xlong(i,j).le.240)).AND. &
80             ((ylat(i,j).ge.40).AND.(ylat(i,j).le.70))) ) mask_cal(i,j,8)=.true.
81! -- Groenland
82        IF( (((xlong(i,j).ge.290).AND.(xlong(i,j).le.310)).AND. &
83            ((ylat(i,j).ge.75).AND.(ylat(i,j).le.85))).OR.      &
84            (((xlong(i,j).ge.300).AND.(xlong(i,j).lt.310)).AND. &
85            ((ylat(i,j).ge.60).AND.(ylat(i,j).le.75))).OR.      &
86!           (((xlong(i,j).ge.310).AND.(xlong(i,j).le.350)).AND. &
87!           ((ylat(i,j).ge.40).AND.(ylat(i,j).le.85))) ) THEN
88            (((xlong(i,j).ge.310).AND.(xlong(i,j).le.345)).AND. &
89            ((ylat(i,j).ge.54).AND.(ylat(i,j).le.85))) ) mask_cal(i,j,3)=.true.
90! -- Fennoscandie
91        IF( (((xlong(i,j).ge.345).AND.(xlong(i,j).lt.360)).AND. & !British
92            ((ylat(i,j).ge.50).AND.(ylat(i,j).le.60))).OR.      &
93            (((xlong(i,j).ge.0).AND.(xlong(i,j).le.40)).AND.    & !Scand
94            ((ylat(i,j).ge.50).AND.(ylat(i,j).le.70))).OR.      &
95!            (((xlong(i,j).ge.50).AND.(xlong(i,j).le.120)).AND.  &
96            (((xlong(i,j).ge.10).AND.(xlong(i,j).le.122)).AND.  & !BK
97            ((ylat(i,j).ge.65).AND.(ylat(i,j).le.85))) ) mask_cal(i,j,4)=.true.
98!! -- Scandinavian Ice Sheet Pauline
99        IF( ((xlong(i,j).gt.0.).AND.(xlong(i,j).le.40)).AND. &
100             ((ylat(i,j).ge.50).AND.(ylat(i,j).le.70))) mask_cal(i,j,9)=.true.
101! -- Barents Ice Sheet Pauline
102        IF ( (((xlong(i,j).gt.10.).AND.(xlong(i,j).le.60.)).AND. &
103             ((ylat(i,j).gt.70).AND.(ylat(i,j).le.85.))).OR.    &
104             (((xlong(i,j).gt.40.).AND.(xlong(i,j).le.60.)).AND. &
105              ((ylat(i,j).ge.65.).AND.(ylat(i,j).lt.70.)))) mask_cal(i,j,10)=.true.
106! -- Kara Ice Sheet Pauline
107        IF( (((xlong(i,j).gt.60.).AND.(xlong(i,j).le.122)).AND. &
108             ((ylat(i,j).ge.70).AND.(ylat(i,j).le.85))).OR.    &
109             (((xlong(i,j).gt.60.).AND.(xlong(i,j).le.80.)).AND. &
110              ((ylat(i,j).ge.65.).AND.(ylat(i,j).lt.70.))))  mask_cal(i,j,11)=.true.
111! -- British Ice Sheet Pauline
112        IF( ((xlong(i,j).gt.350.).AND.(xlong(i,j).le.360)).AND. &
113             ((ylat(i,j).ge.50).AND.(ylat(i,j).le.60))) mask_cal(i,j,13)=.true.
114! -- Siberian Ice Sheet
115        IF( ((xlong(i,j).gt.122.).AND.(xlong(i,j).lt.190)).AND. &
116             ((ylat(i,j).ge.60).AND.(ylat(i,j).le.85))) mask_cal(i,j,12)=.true.
117       enddo
118    enddo
119!        open (UNIT=7777,file='mask_hemin40-02.dat')
120!   do j=1,ny
121!        do i=1,nx
122!        if (mask_cal(i,j,2)) then
123!                write_mask(i,j)=2
124!        else
125!                write_mask(i,j)=0
126!        endif   
127!        write(7777,*) i,j,write_mask(i,j)
128!        enddo     
129!   enddo
130!   close(UNIT=7777)
131!   stop
132
133! ecriture entete avec les colonnes du fichier short :
134WRITE(num_ritz, '(a91)')"! hemin40    : time,isvol,deltavol,inp,isvolf,inf,isacc,isabl,ISABLBORD,ABLATOT,ISCALV,ISBM"
135WRITE(num_ritz, '(a90)')"! colonne    :  1    2      3       4    5     6    7     8      9        10      11    12"
136WRITE(num_ritz,'(a157)')"! calotte    : Laurentide Groenland Fennoscandie LabradorSector KeewatinSector InnuitianIS CordilleranIS ScandinavianIS BarentsIS KaraIS SiberianIS BritishIS"
137WRITE(num_ritz,'(a153)')"! isvol      :    13        26          39            52            65             78           91            104         117      130      143       156"
138WRITE(num_ritz,'(a153)')"! inp        :    14        27          40            53            66             79           92            105         118      131      144       157"
139WRITE(num_ritz,'(a153)')"! isvolf     :    15        28          41            54            67             80           93            106         119      132      145       158"
140WRITE(num_ritz,'(a153)')"! inf        :    16        29          42            55            68             81           94            107         120      133      146       159"
141WRITE(num_ritz,'(a153)')"! isacc      :    17        30          43            56            69             82           95            108         121      134      147       160"
142WRITE(num_ritz,'(a153)')"! isabl      :    18        31          44            57            70             83           96            109         122      135      148       161"
143WRITE(num_ritz,'(a153)')"! isablbord  :    19        32          45            58            71             84           97            110         123      136      149       162"
144WRITE(num_ritz,'(a153)')"! ablatot    :    20        33          46            59            72             85           98            111         124      137      150       163"
145WRITE(num_ritz,'(a153)')"! iscalv     :    21        34          47            60            73             86           99            112         125      138      151       164"
146WRITE(num_ritz,'(a153)')"! isbm       :    22        35          48            61            74             87          100            113         126      139      152       165"
147WRITE(num_ritz,'(a153)')"! itjjamean_ :    23        36          49            62            75             88          101            114         127      140      153       166"
148WRITE(num_ritz,'(a153)')"! hmean_     :    24        37          50            63            76             89          102            115         128      141      154       167"
149WRITE(num_ritz,'(a153)')"! Hmax_      :    25        38          51            64            77             90          103            116         129      142      155       168"
150WRITE(num_ritz,*)
151end subroutine init_outshort
152
153!_________________________________________________________________________
154subroutine shortoutput
155
156! 1_initialization
157!------------------
158      integer KK
159!     integer inpl, INPG, INPF
160!     integer inplab,inpkeew,inpinn,inpcord
161!     integer inpscaN,inpbar,inpkara,inpsib,inpbrit
162     integer inp(13) ! Surface posee (nb de noeuds)
163     integer inf(13) ! surface flottante
164
165!Variables pour sommer     
166      real isvol(13),isvolf(13) ! volume posé et flottants
167      REAL ISCALV(13),ISACC(13),ISBM(13),ISABL(13)
168      REAL ISABLBORD(13),ABLATOT(13),TACC(13),TBM(13)
169      REAL ITJJA(13)
170!      REAL ABLAMEAN
171!moyennes utilisées en output
172      REAL HMAX_(13) , HMEAN_(13) 
173      REAL BMEAN_(13) , ACCMEAN_(13) , ABLMEAN_(13) , ABLBORDMEAN_(13)
174      REAL ABLATOTMEAN_(13) , CALVMEAN_(13) , ITJJAMEAN_(13)
175
176
177
178!     REAL ITJJAMEAN_L, ITJJAMEAN_G, ITJJAMEAN_F ...
179!     REAL ITJJAMEAN_LAB,ITJJAMEAN_KEEW,ITJJAMEAN_INN,ITJJAMEAN_CORD
180!     REAL ITJJAMEAN_SCAN,ITJJAMEAN_BAR,ITJJAMEAN_KARA,ITJJAMEAN_SIB,ITJJAMEAN_BRIT
181
182      REAL HDOTMEAN_G
183!      REAL ABLA(NX,NY)
184      REAL DELTAVOL
185! afq for CONSEAU      REAL, dimension(nx,ny) :: delta_H_dtt
186
187      ilc_units = dx*dy/DICE
188
189!     open(unit=4145,file='reg_output_nord.dat')
190
191      BMELTMEAN=0.
192      TBMEAN=0.
193      TBDOTMEAN=0.
194      VSMEAN=0.
195      UZSMEAN=0.
196      UZSDOTMEAN=0.
197      UZKMEAN=0.
198      HDOTMEAN=0.
199      HDOTMEAN_G = 0.
200      BDOTMEAN=0.
201
202! afq for CONSEAU      runof_oc(:,:)    = 0.
203! afq for CONSEAU      bmelt_oc(:,:)    = 0.
204! afq for CONSEAU      delta_H_dtt(:,:) = 0.
205     
206      DO kk = 1,13
207        INP(kk) = 0
208        INF(kk) = 0
209        ISVOL(kk) = 0.
210        ISVOLF(kk) = 0.
211        ISBM(kk) = 0.
212        ISACC(kk) = 0.
213        ISCALV(kk) = 0.
214        ISABL(kk) = 0.
215        ISABLBORD(kk) = 0.
216        ABLATOT(kk)= 0.
217        TACC(kk) = 0.
218        TBM(kk) = 0.
219        ITJJA(kk) = 0.
220       
221! nouveau tof mai 2009
222!        where (mask_cal(:,:,kk).and.H(:,:).gt.2..and.flot(:,:)) ISVOLF(kk) = ISVOLF(kk) + H(I,J)
223        isvolf(kk) = sum(H(:,:),mask=(mask_cal(:,:,kk).and.(H(:,:).gt.2.).and.flot(:,:)))
224        INF(kk) = count(mask_cal(:,:,kk).and.(H(:,:).gt.2.).and.flot(:,:))
225
226        INP(kk) = count(mask_cal(:,:,kk).and.(H(:,:).gt.2.).and..not.flot(:,:))
227        isvol(kk) = sum(H(:,:),mask=(mask_cal(:,:,kk).and.(H(:,:).gt.2.).and..not.flot(:,:)))
228        isacc(kk) = sum(Acc(:,:),mask=(mask_cal(:,:,kk).and.(H(:,:).gt.2.).and..not.flot(:,:)))
229        isabl(kk) = sum(Abl(:,:),mask=(mask_cal(:,:,kk).and.(H(:,:).gt.2.).and..not.flot(:,:)))
230        itjja(kk) = sum(Tjuly(:,:),mask=(mask_cal(:,:,kk).and.(H(:,:).gt.2.).and..not.flot(:,:)))
231        Hmax_(kk) = maxval(H(:,:),mask=(mask_cal(:,:,kk).and.(H(:,:).gt.2.).and..not.flot(:,:)))
232        Tacc(kk) = sum(Acc(:,:),mask=(mask_cal(:,:,kk).and.(H(:,:).gt.2.)))
233        isablbord(kk) = sum(ablbord(:,:),mask=(mask_cal(:,:,kk)))
234        iscalv(kk) =  sum(calv(:,:),mask=(mask_cal(:,:,kk).and.(H(:,:).gt.2.)))/dt
235
236        ablatot(kk) = isabl(kk) + isablbord(kk)
237        isbm(kk) = isacc(kk)+isabl(kk)+iscalv(kk)+isablbord(kk)
238        Tbm(kk) = Tacc(kk)+isabl(kk)+iscalv(kk)+isablbord(kk)
239
240
241     enddo
242
243! afq for CONSEAU! afq -- thickness variation
244! afq for CONSEAU     delta_H_dtt(:,:) = old_H_dtt(:,:) - H(:,:)
245! afq for CONSEAU! afq -- continental runoff     
246! afq for CONSEAU!        this goes to ECBilt, will be sum up, then routed to CLIO
247! afq for CONSEAU!     where (mask_cal(:,:,1).and.((H(:,:).gt.1.).or.(old_H_dtt(:,:).gt.1.)).and..not.flot(:,:))
248! afq for CONSEAU!     afq I think we do not need the H condition:
249! afq for CONSEAU     where  (mask_cal(:,:,1).and..not.flot(:,:))
250! afq for CONSEAU        runof_oc(:,:) = runof_oc(:,:) + DBLE(delta_H_dtt(:,:))*DX*DY
251! afq for CONSEAU     endwhere
252! afq for CONSEAU! afq -- shelf runoff
253! afq for CONSEAU!        this goes to bilinear interpolation, then put to CLIO
254! afq for CONSEAU!        we have to make sure that we are not considering a "calved" point
255! afq for CONSEAU     where (mask_cal(:,:,1).and.(H(:,:).gt.1.).and.flot(:,:).and.(calv(:,:).eq.0.))
256! afq for CONSEAU        bmelt_oc(:,:) = bmelt_oc(:,:) + DBLE(delta_H_dtt(:,:))*DX*DY
257! afq for CONSEAU     endwhere
258
259! afq for CONSEAU! afq -- old_h_dtt has to be reset
260! afq for CONSEAU     old_H_dtt(:,:)  = H(:,:)
261
262! afq for CONSEAU ! afq -- the accumulated calving need to be reset, BUT not here because of coupling freq.   
263! afq for CONSEAU     calvin_GRIS(:,:) = calvin_GRIS(:,:)+calvingGRIS(:,:)
264! afq for CONSEAU    calvingGRIS(:,:) = 0.
265     
266
267 
268do K=1,13
269
270! == Les moyennes     
271      IF(INP(K).ne.0)THEN
272          HMEAN_(K) = ISVOL(K)   /INP(K) ! /(DX*DY*INP(K))
273          BMEAN_(K) = ISBM(K)    /INP(K)
274          ACCMEAN_(K) = ISACC(K) /INP(K)
275          ABLMEAN_(K) = ISABL(K) /INP(K)
276          ABLBORDMEAN_(K) = ISABLBORD(K) /INP(K)
277          CALVMEAN_(K)    = ISCALV(K)    /INP(K)
278          ABLATOTMEAN_(K) = ABLATOT(K)   /INP(K)
279          ITJJAMEAN_(K)   = ITJJA(K)     /INP(K)
280       ENDIF
281! == Les volmes intergrées (3D)
282      ISVOL(K) = ISVOL(K)*DX*DY     
283              isacc(k)=isacc(k)*dx*dy
284              isabl(k)=isabl(k)*dx*dy
285              isablbord(k)=isablbord(k)*dx*dy
286              ablatot(k)=ablatot(k)*dx*dy
287              iscalv(k)=iscalv(k)*dx*dy
288              isbm(k)=isbm(k)*dx*dy
289              tacc(k)=tacc(k)*dx*dy
290              tbm(k)=tbm(k)*dx*dy
291enddo
292!nom_table='abl'
293!call printtable_r(abl,nom_table)
294
295!====================================================
296! -- Hemin40                : : 1   
297
298      if (NP.ne.0) then 
299!        HMEAN=VOL/NP
300!        VOL=VOL*DX*DY
301!        ISVOL(1) = VOL
302!        BMEAN=BMEAN/NP
303!        ACCMEAN=ACCMEAN/NP
304!        ITJJAMEAN = ITJJA(1)/NP
305!        ABLMEAN=BMEAN-ACCMEAN
306!!!     ABLMEAN = ABLMEAN/NP
307!        CALVMEAN=CALVMEAN/NP
308!!!     BMELTMEAN=BMELTMEAN/NP
309!!!     ABLBORDMEAN=ABLBORDMEAN/NP
310!!!     ABLAMEAN = ABLAMEAN/NP
311        TBMEAN=TBMEAN/NP
312        TBDOTMEAN=TBDOTMEAN/NP
313        VSMEAN=VSMEAN/NP
314!       VSDOTMEAN=VSDOTMEAN/NP
315        UZSMEAN=UZSMEAN/NP
316        UZSDOTMEAN=UZSDOTMEAN/NP
317        UZKMEAN=UZKMEAN/NP
318        DELTAVOL = HDOTMEAN*DX*DY
319        HDOTMEAN=HDOTMEAN/NP 
320    ELSE
321        TBMEAN=0.
322        TBDOTMEAN=0.
323        VSMEAN=0.
324        UZSMEAN=0.
325        UZSDOTMEAN=0.
326        UZKMEAN=0.
327        DELTAVOL = 0.
328        HDOTMEAN=0.
329      endif
330
331      BDOTMEAN=BDOTMEAN/NX/NY 
332
333! afq -- iLOVECLIM water conservation:
334!        units: GRISLI is in m i.e. / yr
335!               iLOVECLIM expects m3 w.e. /yr
336      trendWAC        = diff_H*ilc_units !or water_bilan
337      smbWAC(:,:)     = (-ablbord_dtt(:,:)+ Bm_dtt(:,:))*ilc_units/dtt
338      bmeltWAC(:,:)   = -Bmelt_dtt(:,:)*ilc_units/dtt
339      calvingWAC(:,:) = Calv_dtt(:,:)*ilc_units/dtt
340!check->
341      if (abs(1-diff_H/(water_bilan+1e-20)).gt.0.01) then
342         write(*,*) "Water not conserved in GRISLI", diff_H, water_bilan, abs(1-diff_H/(water_bilan+1e-20))
343      endif
344
345! 2_writing outputs
346!------------------
347
348!         WRITE(35,904)TIME,VOL,DELTAVOL,NP,ISACC(1),ISABL(1),      &
349!          ISABLBORD(1),ABLATOT(1),ISCALV(1),ISBM(1),               &
350!          ISVOL(2),INPL,ISACC(2),                                  &
351!          ISABL(2),ISABLBORD(2),ABLATOT(2),ISCALV(2),ISBM(2),      &
352!          ISVOL(3),INPG,ISACC(3),ISABL(3),ISABLBORD(3),            &
353!          ABLATOT(3),ISCALV(3),ISBM(3),ISVOL(4),INPF,ISACC(4),     &
354!          ISABL(4),ISABLBORD(4),ABLATOT(4),ISCALV(4),ISBM(4),      &
355!          nint(HMEAN),nint(HMAX),BMEAN,TBMEAN,nint(VSMEAN),        &
356!          TBDOTMEAN,HDOTMEAN,BDOTMEAN,BMELTMEAN,NPAB,    &
357!!          TBDOTMEAN,VSDOTMEAN,HDOTMEAN,BDOTMEAN,BMELTMEAN,NPAB,    &
358!          NPCAL,HDOTMEAN_G,TACC(1),TBM(1),TACC(2),TBM(2),          &
359!          TACC(3),TBM(3),TACC(4),TBM(4),                           &
360!           ITJJAMEAN,ITJJAMEAN_L,ITJJAMEAN_G,ITJJAMEAN_F,          &
361!           ISVOL(5),INPLAB,ISACC(5),ABLATOT(5),ITJJAMEAN_LAB,      &
362!           ISVOL(6),INPKEEW,ISACC(6),ABLATOT(6),ITJJAMEAN_KEEW,    &
363!           ISVOL(7),INPINN,ISACC(7),ABLATOT(7),ITJJAMEAN_INN,      &
364!           ISVOL(8),INPCORD,ISACC(8),ABLATOT(8),ITJJAMEAN_CORD,    &
365!           ISVOL(9),INPSCAN,ISACC(9),ABLATOT(9),ITJJAMEAN_SCAN,    &
366!           ISVOL(10),INPBAR,ISACC(10),ABLATOT(10),ITJJAMEAN_BAR,   &
367!           ISVOL(11),INPKARA,ISACC(11),ABLATOT(11),ITJJAMEAN_KARA, &     
368!           ISVOL(12),INPSIB,ISACC(12),ABLATOT(12),ITJJAMEAN_SIB,   &     
369!           ISVOL(13),INPBRIT,ISACC(13),ABLATOT(13),ITJJAMEAN_BRIT, &
370!           HMEAN,HMEAN_L,HMEAN_G,HMEAN_F,HMEAN_LAB,HMEAN_KEEW,     &
371!           HMEAN_INN,HMEAN_CORD,HMEAN_SCAN,HMEAN_BAR,HMEAN_KARA,   &
372!           HMEAN_SIB,HMEAN_BRIT,                                   &
373!           HMAX,HMAX_L,HMAX_G,HMAX_F,HMAX_LAB,HMAX_KEEW,HMAX_INN,  &
374!           HMAX_CORD, HMAX_SCAN,HMAX_BAR,HMAX_KARA,HMAX_SIB,HMAX_BRIT     
375
376
377
378
379!print*,'========================='
380!print*,'ecrit short, time = ',time,dt
381!print*,'hdot ISABLBORD',hdot(30:40,100),ISABLBORD(2)
382
383          write(num_ritz,905)   time,isvol(1),deltavol,inp(1),isvolf(1),inf(1),isacc(1),isabl(1),    &   ! 8
384              ISABLBORD(1),ABLATOT(1),ISCALV(1),ISBM(1),     &                                     ! 4
385  isvol(2),inp(2),isvolf(2),inf(2),isacc(2),isabl(2),isablbord(2),ablatot(2),iscalv(2),isbm(2),itjjamean_(2),hmean_(2),Hmax_(2), & ! 13
386  isvol(3),inp(3),isvolf(3),inf(3),isacc(3),isabl(3),isablbord(3),ablatot(3),iscalv(3),isbm(3),itjjamean_(3),hmean_(3),hmax_(3), & ! 13
387  isvol(4),inp(4),isvolf(4),inf(4),isacc(4),isabl(4),isablbord(4),ablatot(4),iscalv(4),isbm(4),itjjamean_(4),hmean_(4),hmax_(4), & ! 13
388  isvol(5),inp(5),isvolf(5),inf(5),isacc(5),isabl(5),isablbord(5),ablatot(5),iscalv(5),isbm(5),itjjamean_(5),hmean_(5),hmax_(5), & ! 13
389  isvol(6),inp(6),isvolf(6),inf(6),isacc(6),isabl(6),isablbord(6),ablatot(6),iscalv(6),isbm(6),itjjamean_(6),hmean_(6),hmax_(6), & ! 13
390  isvol(7),inp(7),isvolf(7),inf(7),isacc(7),isabl(7),isablbord(7),ablatot(7),iscalv(7),isbm(7),itjjamean_(7),hmean_(7),hmax_(7), & ! 13
391  isvol(8),inp(8),isvolf(8),inf(8),isacc(8),isabl(8),isablbord(8),ablatot(8),iscalv(8),isbm(8),itjjamean_(8),hmean_(8),hmax_(8), & ! 13
392  isvol(9),inp(9),isvolf(9),inf(9),isacc(9),isabl(9),isablbord(9),ablatot(9),iscalv(9),isbm(9),itjjamean_(9),hmean_(9),hmax_(9), & ! 13
393  isvol(10),inp(10),isvolf(10),inf(10),isacc(10),isabl(10),isablbord(10),ablatot(10),iscalv(10),isbm(10),itjjamean_(10),hmean_(10),hmax_(10), & ! 13
394  isvol(11),inp(11),isvolf(11),inf(11),isacc(11),isabl(11),isablbord(11),ablatot(11),iscalv(11),isbm(11),itjjamean_(11),hmean_(11),hmax_(11), & ! 13
395  isvol(12),inp(12),isvolf(12),inf(12),isacc(12),isabl(12),isablbord(12),ablatot(12),iscalv(12),isbm(12),itjjamean_(12),hmean_(12),hmax_(12), & ! 13
396  isvol(13),inp(13),isvolf(13),inf(13),isacc(13),isabl(13),isablbord(13),ablatot(13),iscalv(13),isbm(13),itjjamean_(13),hmean_(13),hmax_(13), & ! 13
397  diff_H, water_bilan, sum(calv_dtt(:,:))/dtt, sum(ablbord_dtt(:,:))/dtt, sum(Bm_dtt(:,:))/dtt, sum(bmelt_dtt(:,:))/dtt !6
398
399 
400! !! format 900 faux, a reprendre
401
402905 format(f9.1,1x, e11.4,1x,e11.5, 1x,i5, 1x,e10.4,1x,i5 , 6(1x,e12.5), &
403         12(2(1x,e10.4,1x,i5), 9(1x,e11.4) ) , 6(1x,e15.8)) 
404 
405   
406940   format('%%%% ',a,'   time=',f8.0,' %%%%')
407end subroutine shortoutput
408
409
410end module output_hemin40_mod
Note: See TracBrowser for help on using the repository browser.