Merge branch 'master' into devel
[wrffire.git] / wrfv2_fire / dyn_em / module_initialize_heldsuarez.F
blob6c2c36479d0101b3036c8db4c5a01c7d8130867c
1 !IDEAL:MODEL_LAYER:INITIALIZATION
3 !  This MODULE holds the routines which are used to perform various initializations
4 !  for the individual domains.  
6 !-----------------------------------------------------------------------
8 MODULE module_initialize_ideal
10    USE module_domain             ! frame/module_domain.F
11    USE module_io_domain          ! share
12    USE module_state_description  ! frame
13    USE module_model_constants    ! share
14    USE module_bc                 ! share
15    USE module_timing             ! frame
16    USE module_configure          ! frame
17    USE module_init_utilities     ! dyn_em
18 #ifdef DM_PARALLEL
19    USE module_dm
20 #endif
23 CONTAINS
26 !-------------------------------------------------------------------
27 ! this is a wrapper for the solver-specific init_domain routines.
28 ! Also dereferences the grid variables and passes them down as arguments.
29 ! This is crucial, since the lower level routines may do message passing
30 ! and this will get fouled up on machines that insist on passing down
31 ! copies of assumed-shape arrays (by passing down as arguments, the 
32 ! data are treated as assumed-size -- ie. f77 -- arrays and the copying
33 ! business is avoided).  Fie on the F90 designers.  Fie and a pox.
35    SUBROUTINE init_domain ( grid )
37    IMPLICIT NONE
39    !  Input data.
40    TYPE (domain), POINTER :: grid 
41    !  Local data.
42    INTEGER :: idum1, idum2
44    CALL set_scalar_indices_from_config ( head_grid%id , idum1, idum2 )
46      CALL init_domain_rk( grid &
48 #include <actual_new_args.inc>
50                         )
52    END SUBROUTINE init_domain
54 !-------------------------------------------------------------------
56    SUBROUTINE init_domain_rk ( grid &
58 # include <dummy_new_args.inc>
61    IMPLICIT NONE
63    !  Input data.
64    TYPE (domain), POINTER :: grid
66 # include <dummy_decl.inc>
68    TYPE (grid_config_rec_type)              :: config_flags
70    !  Local data
71    INTEGER                             ::                       &
72                                   ids, ide, jds, jde, kds, kde, &
73                                   ims, ime, jms, jme, kms, kme, &
74                                   its, ite, jts, jte, kts, kte, &
75                                   i, j, k
77    INTEGER :: nxx, nyy, ig, jg, im, error
79    REAL :: dlam, dphi, vlat, tperturb
80    REAL :: p_surf, p_level, pd_surf, qvf1, qvf2, qvf
81    REAL :: thtmp, ptmp, temp(3), cof1, cof2
83    INTEGER :: icm,jcm
85    SELECT CASE ( model_data_order )
86          CASE ( DATA_ORDER_ZXY )
87    kds = grid%sd31 ; kde = grid%ed31 ;
88    ids = grid%sd32 ; ide = grid%ed32 ;
89    jds = grid%sd33 ; jde = grid%ed33 ;
91    kms = grid%sm31 ; kme = grid%em31 ;
92    ims = grid%sm32 ; ime = grid%em32 ;
93    jms = grid%sm33 ; jme = grid%em33 ;
95    kts = grid%sp31 ; kte = grid%ep31 ;   ! note that tile is entire patch
96    its = grid%sp32 ; ite = grid%ep32 ;   ! note that tile is entire patch
97    jts = grid%sp33 ; jte = grid%ep33 ;   ! note that tile is entire patch
98          CASE ( DATA_ORDER_XYZ )
99    ids = grid%sd31 ; ide = grid%ed31 ;
100    jds = grid%sd32 ; jde = grid%ed32 ;
101    kds = grid%sd33 ; kde = grid%ed33 ;
103    ims = grid%sm31 ; ime = grid%em31 ;
104    jms = grid%sm32 ; jme = grid%em32 ;
105    kms = grid%sm33 ; kme = grid%em33 ;
107    its = grid%sp31 ; ite = grid%ep31 ;   ! note that tile is entire patch
108    jts = grid%sp32 ; jte = grid%ep32 ;   ! note that tile is entire patch
109    kts = grid%sp33 ; kte = grid%ep33 ;   ! note that tile is entire patch
110          CASE ( DATA_ORDER_XZY )
111    ids = grid%sd31 ; ide = grid%ed31 ;
112    kds = grid%sd32 ; kde = grid%ed32 ;
113    jds = grid%sd33 ; jde = grid%ed33 ;
115    ims = grid%sm31 ; ime = grid%em31 ;
116    kms = grid%sm32 ; kme = grid%em32 ;
117    jms = grid%sm33 ; jme = grid%em33 ;
119    its = grid%sp31 ; ite = grid%ep31 ;   ! note that tile is entire patch
120    kts = grid%sp32 ; kte = grid%ep32 ;   ! note that tile is entire patch
121    jts = grid%sp33 ; jte = grid%ep33 ;   ! note that tile is entire patch
123    END SELECT
125    CALL model_to_grid_config_rec ( grid%id , model_config_rec , config_flags )
127 ! here we check to see if the boundary conditions are set properly
129    CALL boundary_condition_check( config_flags, bdyzone, error, grid%id )
131     grid%itimestep=0
132    grid%step_number = 0
135 #ifdef DM_PARALLEL
136    CALL wrf_dm_bcast_bytes( icm , IWORDSIZE )
137    CALL wrf_dm_bcast_bytes( jcm , IWORDSIZE )
138 #endif
140    ! Initialize 2D surface arrays
142    nxx = ide-ids ! Don't include u-stagger
143    nyy = jde-jds ! Don't include v-stagger
144    dphi = 180./REAL(nyy)
145    dlam = 360./REAL(nxx)
147    DO j = jts, jte
148    DO i = its, ite
149       ! ig is the I index in the global (domain) span of the array.
150       ! jg is the J index in the global (domain) span of the array.
151       ig = i - ids + 1  ! ids is not necessarily 1
152       jg = j - jds + 1  ! jds is not necessarily 1
154       grid%xlat(i,j)  = (REAL(jg)-0.5)*dphi-90.
155       grid%xlong(i,j) = (REAL(ig)-0.5)*dlam-180.
156       vlat       = grid%xlat(i,j) - 0.5*dphi
158       grid%clat(i,j) = grid%xlat(i,j)
160       grid%msftx(i,j) = 1./COS(grid%xlat(i,j)*degrad)
161       grid%msfty(i,j) = 1.
162       grid%msfux(i,j) = 1./COS(grid%xlat(i,j)*degrad)
163       grid%msfuy(i,j) = 1.
164       grid%e(i,j)     = 2*EOMEG*COS(grid%xlat(i,j)*degrad)
165       grid%f(i,j)     = 2*EOMEG*SIN(grid%xlat(i,j)*degrad)
167       ! The following two are the cosine and sine of the rotation
168       ! of projection.  Simple cylindrical is *simple* ... no rotation!
169       grid%sina(i,j) = 0.
170       grid%cosa(i,j) = 1.
172    END DO
173    END DO
175 !   DO j = max(jds+1,jts), min(jde-1,jte)
176    DO j = jts, jte
177    DO i = its, ite
178       vlat       = grid%xlat(i,j) - 0.5*dphi
179       grid%msfvx(i,j) = 1./COS(vlat*degrad)
180       grid%msfvy(i,j) = 1.
181       grid%msfvx_inv(i,j) = 1./grid%msfvx(i,j)
182    END DO
183    END DO
185    IF(jts == jds) THEN
186    DO i = its, ite
187       grid%msfvx(i,jts) = 00.
188       grid%msfvx_inv(i,jts) = 0.
189    END DO
190    END IF
192    IF(jte == jde) THEN
193    DO i = its, ite
194       grid%msfvx(i,jte) = 00.
195       grid%msfvx_inv(i,jte) = 0.
196    END DO
197    END IF
199    DO j=jts,jte
200      vlat       = grid%xlat(its,j) - 0.5*dphi
201      write(6,*) j,vlat,grid%msfvx(its,j),grid%msfvx_inv(its,j)
202    ENDDO
204    DO j=jts,jte
205    DO i=its,ite
206       grid%ht(i,j)       = 0.
208       grid%albedo(i,j)   = 0.
209       grid%thc(i,j)      = 1000.
210       grid%znt(i,j)      = 0.01
211       grid%emiss(i,j)    = 1.
212       grid%ivgtyp(i,j)   = 1
213       grid%lu_index(i,j) = REAL(ivgtyp(i,j))
214       grid%xland(i,j)    = 1.
215       grid%mavail(i,j)   = 0.
216    END DO
217    END DO
219    grid%dx = dlam*degrad/reradius
220    grid%dy = dphi*degrad/reradius
221    grid%rdx = 1./grid%dx
222    grid%rdy = 1./grid%dy
224    !WRITE(*,*) ''
225    !WRITE(*,'(A,1PG14.6,A,1PG14.6)') ' For the namelist: dx =',grid%dx,', dy =',grid%dy
227    CALL nl_set_mminlu(1,'    ')
228    grid%iswater = 0
229    grid%cen_lat = 0.
230    grid%cen_lon = 0.
231    grid%truelat1 = 0.
232    grid%truelat2 = 0.
233    grid%moad_cen_lat = 0.
234    grid%stand_lon    = 0.
235    grid%pole_lon = 0.
236    grid%pole_lat = 90.
237    ! Apparently, map projection 0 is "none" which actually turns out to be
238    ! a regular grid of latitudes and longitudes, the simple cylindrical projection
239    grid%map_proj = 0
241    DO k = kds, kde
242       grid%znw(k) = 1. - REAL(k-kds)/REAL(kde-kds)
243    END DO
245    DO k=1, kde-1
246     grid%dnw(k) = grid%znw(k+1) - grid%znw(k)
247     grid%rdnw(k) = 1./grid%dnw(k)
248     grid%znu(k) = 0.5*(grid%znw(k+1)+grid%znw(k))
249    ENDDO
250    DO k=2, kde-1
251     grid%dn(k) = 0.5*(grid%dnw(k)+grid%dnw(k-1))
252     grid%rdn(k) = 1./grid%dn(k)
253     grid%fnp(k) = .5* grid%dnw(k  )/grid%dn(k)
254     grid%fnm(k) = .5* grid%dnw(k-1)/grid%dn(k)
255    ENDDO
257    cof1 = (2.*grid%dn(2)+grid%dn(3))/(grid%dn(2)+grid%dn(3))*grid%dnw(1)/grid%dn(2) 
258    cof2 =     grid%dn(2)        /(grid%dn(2)+grid%dn(3))*grid%dnw(1)/grid%dn(3) 
259    grid%cf1  = grid%fnp(2) + cof1
260    grid%cf2  = grid%fnm(2) - cof1 - cof2
261    grid%cf3  = cof2       
263    grid%cfn  = (.5*grid%dnw(kde-1)+grid%dn(kde-1))/grid%dn(kde-1)
264    grid%cfn1 = -.5*grid%dnw(kde-1)/grid%dn(kde-1)
267    ! Need to add perturbations to initial profile.  Set up random number
268    ! seed here.
269    CALL random_seed
271    ! General assumption from here after is that the initial temperature
272    ! profile is isothermal at a value of T0, and the initial winds are
273    ! all 0.
275    ! find ptop for the desired ztop (ztop is input from the namelist)
276    grid%p_top =  p0 * EXP(-(g*config_flags%ztop)/(r_d*T0))
279    ! Values of geopotential (base, perturbation, and at p0) at the surface
280    DO j = jts, jte
281    DO i = its, ite
282       grid%phb(i,1,j) = grid%ht(i,j)*g
283       grid%php(i,1,j) = 0.         ! This is perturbation geopotential
284                               ! Since this is an initial condition, there
285                               ! should be no perturbation!
286       grid%ph0(i,1,j) = grid%ht(i,j)*g
287    ENDDO
288    ENDDO
291    DO J = jts, jte
292    DO I = its, ite
294       p_surf = p0 * EXP(-(g*grid%phb(i,1,j)/g)/(r_d*T0))
295       grid%mub(i,j) = p_surf-grid%p_top
297       ! given p (coordinate), calculate theta and compute 1/rho from equation
298       ! of state
300       DO K = kts, kte-1
301          p_level = grid%znu(k)*(p_surf - grid%p_top) + grid%p_top
302          grid%pb(i,k,j) = p_level
304          grid%t_init(i,k,j) = T0*(p0/p_level)**rcp
305          grid%t_init(i,k,j) = grid%t_init(i,k,j) - t0
307          grid%alb(i,k,j)=(r_d/p1000mb)*(grid%t_init(i,k,j)+t0)*(grid%pb(i,k,j)/p1000mb)**cvpm
308       END DO
310       ! calculate hydrostatic balance (alternatively we could interpolate
311       ! the geopotential from the sounding, but this assures that the base
312       ! state is in exact hydrostatic balance with respect to the model eqns.
314       DO k = kts+1, kte
315          grid%phb(i,k,j) = grid%phb(i,k-1,j) - grid%dnw(k-1)*grid%mub(i,j)*grid%alb(i,k-1,j)
316       ENDDO
318    ENDDO
319    ENDDO
321    DO im = PARAM_FIRST_SCALAR, num_moist
322    DO J = jts, jte
323    DO K = kts, kte-1
324    DO I = its, ite
325       grid%moist(i,k,j,im) = 0.
326    END DO
327    END DO
328    END DO
329    END DO
331    ! Now calculate the full (hydrostatically-balanced) state for each column
332    ! We will include moisture
333    DO J = jts, jte
334    DO I = its, ite
336       ! At this point p_top is already set. find the DRY mass in the column
337       pd_surf = p0 * EXP(-(g*grid%phb(i,1,j)/g)/(r_d*T0))
339       ! compute the perturbation mass (mu/mu_1/mu_2) and the full mass
340       grid%mu_1(i,j) = pd_surf-grid%p_top - grid%mub(i,j)
341       grid%mu_2(i,j) = grid%mu_1(i,j)
342       grid%mu0(i,j)  = grid%mu_1(i,j) + grid%mub(i,j)
344       ! given the dry pressure and coordinate system, calculate the
345       ! perturbation potential temperature (t/t_1/t_2)
347       DO k = kds, kde-1
348          p_level = grid%znu(k)*(pd_surf - grid%p_top) + grid%p_top
349          grid%t_1(i,k,j) = T0*(p0/p_level)**rcp
350          ! Add a small perturbation to initial isothermal profile
351          CALL random_number(tperturb)
352          grid%t_1(i,k,j)=grid%t_1(i,k,j)*(1.0+0.004*(tperturb-0.5))
353          grid%t_1(i,k,j) = grid%t_1(i,k,j)-t0
354          grid%t_2(i,k,j) = grid%t_1(i,k,j)
355       END DO
358       ! integrate the hydrostatic equation (from the RHS of the bigstep
359       ! vertical momentum equation) down from the top to get p.
360       ! first from the top of the model to the top pressure
362       k = kte-1  ! top level
364       qvf1 = 0.5*(grid%moist(i,k,j,P_QV)+grid%moist(i,k,j,P_QV))
365       qvf2 = 1./(1.+qvf1)
366       qvf1 = qvf1*qvf2
368       ! grid%p(i,k,j) = - 0.5*grid%mu_1(i,j)/grid%rdnw(k)
369       grid%p(i,k,j) = - 0.5*(grid%mu_1(i,j)+qvf1*grid%mub(i,j))/grid%rdnw(k)/qvf2
370       qvf = 1. + rvovrd*grid%moist(i,k,j,P_QV)
371       grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* &
372                   (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm)
373       grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j)
375       !  down the column
377       do k=kte-2,kts,-1
378          qvf1 = 0.5*(grid%moist(i,k,j,P_QV)+grid%moist(i,k+1,j,P_QV))
379          qvf2 = 1./(1.+qvf1)
380          qvf1 = qvf1*qvf2
381          grid%p(i,k,j) = grid%p(i,k+1,j) - (grid%mu_1(i,j) + qvf1*grid%mub(i,j))/qvf2/grid%rdn(k+1)
382          qvf = 1. + rvovrd*grid%moist(i,k,j,P_QV)
383          grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* &
384                      (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm)
385          grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j)
386       enddo
388       ! this is the hydrostatic equation used in the model after the
389       ! small timesteps.  In the model, al (inverse density)
390       ! is computed from the geopotential.
392       grid%ph_1(i,1,j) = 0.
393       DO k  = kts+1,kte
394          grid%ph_1(i,k,j) = grid%ph_1(i,k-1,j) - (1./grid%rdnw(k-1))*(  &
395                       (grid%mub(i,j)+grid%mu_1(i,j))*grid%al(i,k-1,j)+  &
396                        grid%mu_1(i,j)*grid%alb(i,k-1,j)  )
398          grid%ph_2(i,k,j) = grid%ph_1(i,k,j)
399          grid%ph0(i,k,j) = grid%ph_1(i,k,j) + grid%phb(i,k,j)
400       ENDDO
402    END DO
403    END DO
407    ! Now set U & V
408    DO J = jts, jte
409    DO K = kts, kte-1
410    DO I = its, ite
411       grid%u_1(i,k,j) = 0.
412       grid%u_2(i,k,j) = 0.
413       grid%v_1(i,k,j) = 0.
414       grid%v_2(i,k,j) = 0.
415    END DO
416    END DO
417    END DO
420    DO j=jts, jte
421    DO k=kds, kde
422    DO i=its, ite
423       grid%ww(i,k,j)  = 0.
424       grid%w_1(i,k,j) = 0.
425       grid%w_2(i,k,j) = 0.
426       grid%h_diabatic(i,k,j) = 0.
427    END DO
428    END DO
429    END DO
431   
432    DO k=kts,kte
433       grid%t_base(k)  = grid%t_init(its,k,jts)
434       grid%qv_base(k) = 0.
435       grid%u_base(k)  = 0.
436       grid%v_base(k)  = 0.
437    END DO
439    ! One subsurface layer: infinite slab at constant temperature below
440    ! the surface.  Surface temperature is an infinitely thin "skin" on
441    ! top of a half-infinite slab.  The temperature of both the skin and
442    ! the slab are determined from the initial nearest-surface-air-layer
443    ! temperature.
444    DO J = jts, MIN(jte, jde-1)
445    DO I = its, MIN(ite, ide-1)
446       thtmp   = grid%t_2(i,1,j)+t0
447       ptmp    = grid%p(i,1,j)+grid%pb(i,1,j)
448       temp(1) = thtmp * (ptmp/p1000mb)**rcp
449       thtmp   = grid%t_2(i,2,j)+t0
450       ptmp    = grid%p(i,2,j)+grid%pb(i,2,j)
451       temp(2) = thtmp * (ptmp/p1000mb)**rcp
452       thtmp   = grid%t_2(i,3,j)+t0
453       ptmp    = grid%p(i,3,j)+grid%pb(i,3,j)
454       temp(3) = thtmp * (ptmp/p1000mb)**rcp
455       grid%tsk(I,J)=cf1*temp(1)+cf2*temp(2)+cf3*temp(3)
456       grid%tmn(I,J)=grid%tsk(I,J)-0.5
457    END DO
458    END DO
460    RETURN
462  END SUBROUTINE init_domain_rk
464 !---------------------------------------------------------------------
466  SUBROUTINE init_module_initialize
467  END SUBROUTINE init_module_initialize
469 !---------------------------------------------------------------------
471 END MODULE module_initialize_ideal