source: trunk/LMDZ.GENERIC/libf/phystd/calcenergy_kcm.F90 @ 1243

Last change on this file since 1243 was 787, checked in by aslmd, 12 years ago

LMDZ.GENERIC. (Sorry for long text but this is a quite major commit)

Paving the path for parallel computations. And moving towards a more flexible code.

Automatic allocation is used within all routines in phystd. No further mention to ngridmx and nqmx.

  1. ngridmx and nqmx are still used in LMDZ.GENERIC in the dyn3d part
  2. if the LMDZ4/LMDZ5 dynamical core is used, there is no more fixed dimensions ngridmx and nqmx --> a fully flexible parallel implementation is now possible (e.g. no need to recompile when changing numbers of processors)

The important stuff :

  • Compilation checked with ifort. OK with and without debug mode. No errors. Checked for: gcm, newstart, rcm1d, kcm1d
  • RUN GCM: Running an Earth test case. Comparison with previous revision --> debug mode : perfect match. bit by bit (diff command). checked with plots --> O1 mode : close match (checked with plots) --> O2 mode : sometimes up to 0.5 K departure.... BUT in this new version O2 and O1 are quite close while in previous version O1 and O2 differed by about, well, typically 0.5 K (pictures available on request)
  • RUN NEWSTART : perfect match (bit-by-bit) in either debug or normal mode.
  • RUN RCM1D : perfect match in normal mode.
  • RUN KCM1D : not tested (I don't know what is the use of kcm1d)

List of main changes :

  • Additional arguments to some subroutines (ngrid and nq)
  • F77 include strategy is obsolete and replaced by F90 module strategy In this new strategy arrays are allocatable and allocated once at first use This has to be done for all common featuring arrays defined with ngridmx or nqmx

surfdat.h >> surfdat_h.F90
tracer.h >> tracer_h.F90
comsaison.h >> comsaison_h.F90
comgeomfi.h >> comgeomfi_h.F90
comsoil.h >> comsoil_h.F90
comdiurn.h >> comdiurn_h.F90
fisice.h >> DELETED. was not used. probably a fossil.
watercap.h >> DELETED. variable put in surfdat_h.F90

  • F77 'save' strategy is obsolete and replaced by F90 'allocatable save' strategy (see previous point and e.g. new version of physiq.F90)
  • Suppressing any mention to advtrac.h which is a common in the dynamics and needs nqmx This was easily solved by adding an argument with tracer names, coming from the dynamics This is probably not a definitive solution, ... but this allows for generic physics to work easily with either LMDZ.GENERIC or LMDZ dynamical cores
  • Removing consistency tests between nq and nqmx ; and ngrid and ngridmx. No use now!
  • Adaptation of rcm1d, kcm1d, newstart given above-mentioned changes

A note on phyetat0 and soil_setting:

  • Now written so that a slice of horizontal size 'ngrid' starting at grid point 'cursor' is read in startfi.nc 'cursor' is defined in dimphys.h and initialized by inifis (or in newstart) this is useful for parallel computations. default behavior is the usual one : sequential runs, cursor is 1, size ngrid is the whole global domain

A note on an additional change :

  • nueffrad is now an argument to callcorrk as is the case for reffrad both are saved in physiq this is for consistency and lisibility (previously nueffrad was saved in callcorrk) ... but there is a call to a function which modifies nueffrad in physiq ... previously this was not modifying nueffrad (although it was quite cumbersome to detect this) ... to be conservative I kept this behaviour and highlighted it with an array nueffrad_dummy ... I added a comment because someone might want to change this
  • Property svn:executable set to *
File size: 1.7 KB
RevLine 
[305]1subroutine calcenergy_kcm(Tsurf,T,Play,Plev,Qsurf,Q,muvar,Eatmtot)
2
3
4use params_h, only : Rc
5use watercommon_h, only : mH2O
6implicit none
7
8!     ----------------------------------------------------------------
9!     Purpose: Calculate total energy of the (steam) atmosphere
10!     Authour: R. Wordsworth (2011)
11!     ----------------------------------------------------------------
12
13#include "dimensions.h"
14#include "dimphys.h"
15#include "comcstfi.h"
16!#include "callkeys.h"
17
18
19
20  ! inputs
21  real Tsurf,T(1:nlayermx)
22  real Play(1:nlayermx),Plev(1:nlayermx+1)
23  real Qsurf,Q(1:nlayermx)
[787]24  real muvar(1,nlayermx+1)
[305]25
26  ! internal
27  integer il
28  real m_n,m_v,cp_n,cp_v,cp_a,Eatmlay
29  real s_c,rho_v,L
30  double precision p_v, s_v, nul
31  real VMR(1:nlayermx)
32
33  ! output
34  real Eatmtot
35
36  ! functions
37  double precision cp_neutral
38
39  m_n = mugaz/1000.
40  m_v = mH2O/1000.
41
42
43  do il=1,nlayermx
44     VMR(il)=Q(il)*(muvar(1,il+1)/mH2O)
45  end do
46
47
48  Eatmtot = 0.0
49  do il=1,nlayermx
50
51
52     cp_n = cp_neutral(dble(T(il)))
53     cp_v = (32.24+1.923d-3*T(il) + 1.055d-5*T(il)**2-3.511d-9*T(il)**3)/m_v
54     ! / by m_v makes it per kg not per mol
55
56     call psat_H2O(dble(T(il)),p_v)
57     p_v   = p_v*1d6
58     rho_v = real(p_v)*m_v/(real(Rc)*T(il))
59
60     call therm(dble(T(il)),dble(rho_v)*1d-3,nul,nul,nul,nul,nul,nul,nul,&
61                nul,nul,nul,s_v,nul)
62
63     s_c = 2.06 * log(T(il)/273.15)
64     s_v = s_v * 1d3
65     s_c = s_c * 1d3
66     L   = (real(s_v)-s_c)*T(il)
67
68     cp_a    = (1-VMR(il))*cp_n + VMR(il)*cp_v
69     !cp_a    = (1-Q(il))*cp_n + Q(il)*cp_v
70
71     Eatmlay = ( cp_a*T(il) + L*VMR(il) ) * (Plev(il)-Plev(il+1))/g
72
73     Eatmtot = Eatmtot + Eatmlay
74  enddo
75
76  print*,'WARNING: E-calc currently assumes H2O saturation!'
77
78end subroutine calcenergy_kcm
79
Note: See TracBrowser for help on using the repository browser.