source: trunk/LMDZ.PLUTO/libf/phypluto/testconservfast.F90 @ 3453

Last change on this file since 3453 was 3412, checked in by afalco, 10 months ago

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

File size: 2.6 KB
Line 
1      subroutine testconservfast(ngrid,nlayer,nq,ptimestep, &
2          pplev,zdq,zdqs,car1,car2)
3
4      use comgeomfi_h         
5      use comcstfi_mod, only: pi, g
6      use geometry_mod, only: latitude, longitude, cell_area
7      implicit none
8
9!==================================================================
10!     Purpose
11!     -------
12!     Test conservation of tracers
13
14!     Inputs
15!     ------
16!     ngrid                 Number of vertical columns
17!     nlayer                Number of layers
18!     pplev(ngrid,nlayer+1) Pressure levels
19!     zdq
20!     
21!     Outputs
22!     -------
23!     
24!     Both
25!     ----
26!
27!     Authors
28!     -------
29!     Tanguy Bertrand 2015
30!     
31!==================================================================
32
33#include "dimensions.h"
34
35!-----------------------------------------------------------------------
36!     Arguments
37
38      INTEGER ngrid, nlayer, nq
39      REAL zdq(ngrid)
40      REAL zdqs(ngrid)
41      REAL ptimestep
42      REAL pplev(ngrid,nlayer+1)
43      character(len=3) :: car1     
44      character(len=7) :: car2     
45
46      ! LOCAL VARIABLES
47      INTEGER l,ig,iq
48      REAL masse
49
50      ! OUTPUT
51      REAL dWtot   
52      REAL dWtots
53      REAL nconsMAX
54      INTEGER myig
55      REAL vdifcncons(ngrid)
56!-----------------------------------------------------------------------
57
58      dWtot=0.0
59      dWtots=0.0
60      nconsMAX=0.0
61      do ig = 1, ngrid
62         vdifcncons(ig)=0.0
63
64        ! sum of atmospheric mass : kg
65         dWtot = dWtot + zdq(ig)*ptimestep*cell_area(ig)*pplev(ig,1)/g
66        ! for each column, total mass lost per sec : kg(tracer) / m2 /s
67         vdifcncons(ig)=vdifcncons(ig) + pplev(ig,1)/g*zdq(ig)
68 
69         dWtots = dWtots + zdqs(ig)*ptimestep*cell_area(ig)
70         vdifcncons(ig)=vdifcncons(ig)+zdqs(ig)
71
72             
73             ! vdifcncons is the total amount of material that appear or
74             ! disapear per second in the routine
75             ! it is the non conservative factor
76
77         if(vdifcncons(ig).gt.nconsMAX)then
78              nconsMAX=vdifcncons(ig)
79              myig=ig
80         endif
81
82      enddo
83
84      dWtot  = dWtot/totarea
85      dWtots = dWtots/totarea
86      print*,'-------------------------------------------'
87      print*,'In ',car2,' atmospheric ',car1,' change=',dWtot,' kg m-2'
88      print*,'In ',car2,' surface ',car1,' change  =',dWtots,' kg m-2'
89      print*,'--> non-cons factor     =',dWtot+dWtots,' kg m-2'
90      print*,'--> MAX non-cons factor =',nconsMAX,' kg m-2 s-1'
91      IF (nconsMAX.gt.0.) then
92        print*,'--> obtained at lat/lon=',latitude(myig)*180./pi,longitude(myig)*180./pi
93      ENDIF
94      end subroutine testconservfast
95
Note: See TracBrowser for help on using the repository browser.