module_initialize_heldsuarez.F

References to this file elsewhere.
1 !IDEAL:MODEL_LAYER:INITIALIZATION
2 
3 !  This MODULE holds the routines which are used to perform various initializations
4 !  for the individual domains.  
5 
6 !-----------------------------------------------------------------------
7 
8 MODULE module_initialize_ideal
9 
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
21 
22 
23 CONTAINS
24 
25 
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.
34 
35    SUBROUTINE init_domain ( grid )
36 
37    IMPLICIT NONE
38 
39    !  Input data.
40    TYPE (domain), POINTER :: grid 
41    !  Local data.
42    INTEGER                :: dyn_opt 
43    INTEGER :: idum1, idum2
44 
45    CALL nl_get_dyn_opt( 1, dyn_opt )
46    
47    CALL set_scalar_indices_from_config ( head_grid%id , idum1, idum2 )
48 
49    IF (      dyn_opt .eq. 1 &
50         .or. dyn_opt .eq. 2 &
51         .or. dyn_opt .eq. 3 &
52                                        ) THEN
53      CALL init_domain_rk( grid &
54 !
55 #include <em_actual_new_args.inc>
56 !
57                         )
58 
59    ELSE
60      WRITE(0,*)' init_domain: unknown or unimplemented dyn_opt = ',dyn_opt
61      call wrf_error_fatal ( ' init_domain: unknown or unimplemented dyn_opt ' )
62    ENDIF
63 
64    END SUBROUTINE init_domain
65 
66 !-------------------------------------------------------------------
67 
68    SUBROUTINE init_domain_rk ( grid &
69 !
70 # include <em_dummy_new_args.inc>
71 !
72 )
73    IMPLICIT NONE
74 
75    !  Input data.
76    TYPE (domain), POINTER :: grid
77 
78 # include <em_dummy_decl.inc>
79 
80    TYPE (grid_config_rec_type)              :: config_flags
81 
82    !  Local data
83    INTEGER                             ::                       &
84                                   ids, ide, jds, jde, kds, kde, &
85                                   ims, ime, jms, jme, kms, kme, &
86                                   its, ite, jts, jte, kts, kte, &
87                                   i, j, k
88 
89    INTEGER :: nxx, nyy, ig, jg, im, error
90 
91    REAL :: dlam, dphi, vlat, tperturb
92    REAL :: p_surf, p_level, pd_surf, qvf1, qvf2, qvf
93    REAL :: thtmp, ptmp, temp(3), cof1, cof2
94 
95    INTEGER :: icm,jcm
96 
97 #ifdef DM_PARALLEL
98 #    include <em_data_calls.inc>
99 #endif
100 
101    SELECT CASE ( model_data_order )
102          CASE ( DATA_ORDER_ZXY )
103    kds = grid%sd31 ; kde = grid%ed31 ;
104    ids = grid%sd32 ; ide = grid%ed32 ;
105    jds = grid%sd33 ; jde = grid%ed33 ;
106 
107    kms = grid%sm31 ; kme = grid%em31 ;
108    ims = grid%sm32 ; ime = grid%em32 ;
109    jms = grid%sm33 ; jme = grid%em33 ;
110 
111    kts = grid%sp31 ; kte = grid%ep31 ;   ! note that tile is entire patch
112    its = grid%sp32 ; ite = grid%ep32 ;   ! note that tile is entire patch
113    jts = grid%sp33 ; jte = grid%ep33 ;   ! note that tile is entire patch
114          CASE ( DATA_ORDER_XYZ )
115    ids = grid%sd31 ; ide = grid%ed31 ;
116    jds = grid%sd32 ; jde = grid%ed32 ;
117    kds = grid%sd33 ; kde = grid%ed33 ;
118 
119    ims = grid%sm31 ; ime = grid%em31 ;
120    jms = grid%sm32 ; jme = grid%em32 ;
121    kms = grid%sm33 ; kme = grid%em33 ;
122 
123    its = grid%sp31 ; ite = grid%ep31 ;   ! note that tile is entire patch
124    jts = grid%sp32 ; jte = grid%ep32 ;   ! note that tile is entire patch
125    kts = grid%sp33 ; kte = grid%ep33 ;   ! note that tile is entire patch
126          CASE ( DATA_ORDER_XZY )
127    ids = grid%sd31 ; ide = grid%ed31 ;
128    kds = grid%sd32 ; kde = grid%ed32 ;
129    jds = grid%sd33 ; jde = grid%ed33 ;
130 
131    ims = grid%sm31 ; ime = grid%em31 ;
132    kms = grid%sm32 ; kme = grid%em32 ;
133    jms = grid%sm33 ; jme = grid%em33 ;
134 
135    its = grid%sp31 ; ite = grid%ep31 ;   ! note that tile is entire patch
136    kts = grid%sp32 ; kte = grid%ep32 ;   ! note that tile is entire patch
137    jts = grid%sp33 ; jte = grid%ep33 ;   ! note that tile is entire patch
138 
139    END SELECT
140 
141    CALL model_to_grid_config_rec ( grid%id , model_config_rec , config_flags )
142 
143 ! here we check to see if the boundary conditions are set properly
144 
145    CALL boundary_condition_check( config_flags, bdyzone, error, grid%id )
146 
147     grid%itimestep=0
148    grid%step_number = 0
149 
150 
151 #ifdef DM_PARALLEL
152    CALL wrf_dm_bcast_bytes( icm , IWORDSIZE )
153    CALL wrf_dm_bcast_bytes( jcm , IWORDSIZE )
154 #endif
155 
156    ! Initialize 2D surface arrays
157 
158    nxx = ide-ids ! Don't include u-stagger
159    nyy = jde-jds ! Don't include v-stagger
160    dphi = 180./REAL(nyy)
161    dlam = 360./REAL(nxx)
162 
163    DO j = jts, jte
164    DO i = its, ite
165       ! ig is the I index in the global (domain) span of the array.
166       ! jg is the J index in the global (domain) span of the array.
167       ig = i - ids + 1  ! ids is not necessarily 1
168       jg = j - jds + 1  ! jds is not necessarily 1
169 
170       grid%xlat(i,j)  = (REAL(jg)-0.5)*dphi-90.
171       grid%xlong(i,j) = (REAL(ig)-0.5)*dlam-180.
172       vlat       = grid%xlat(i,j) - 0.5*dphi
173 
174       grid%clat(i,j) = grid%xlat(i,j)
175       grid%clong(i,j) = grid%xlong(i,j)
176 
177       grid%msftx(i,j) = 1./COS(grid%xlat(i,j)*degrad)
178       grid%msfty(i,j) = 1.
179       grid%msfux(i,j) = 1./COS(grid%xlat(i,j)*degrad)
180       grid%msfuy(i,j) = 1.
181       grid%e(i,j)     = 2*EOMEG*COS(grid%xlat(i,j)*degrad)
182       grid%f(i,j)     = 2*EOMEG*SIN(grid%xlat(i,j)*degrad)
183 
184       ! The following two are the cosine and sine of the rotation
185       ! of projection.  Simple cylindrical is *simple* ... no rotation!
186       grid%sina(i,j) = 0.
187       grid%cosa(i,j) = 1.
188 
189    END DO
190    END DO
191 
192 !   DO j = max(jds+1,jts), min(jde-1,jte)
193    DO j = jts, jte
194    DO i = its, ite
195       vlat       = grid%xlat(i,j) - 0.5*dphi
196       grid%msfvx(i,j) = 1./COS(vlat*degrad)
197       grid%msfvy(i,j) = 1.
198       grid%msfvx_inv(i,j) = 1./grid%msfvx(i,j)
199    END DO
200    END DO
201 
202    IF(jts == jds) THEN
203    DO i = its, ite
204       grid%msfvx(i,jts) = 00.
205       grid%msfvx_inv(i,jts) = 0.
206    END DO
207    END IF
208 
209    IF(jte == jde) THEN
210    DO i = its, ite
211       grid%msfvx(i,jte) = 00.
212       grid%msfvx_inv(i,jte) = 0.
213    END DO
214    END IF
215 
216    DO j=jts,jte
217      vlat       = grid%xlat(its,j) - 0.5*dphi
218      write(6,*) j,vlat,grid%msfvx(its,j),grid%msfvx_inv(its,j)
219    ENDDO
220 
221    DO j=jts,jte
222    DO i=its,ite
223       grid%ht(i,j)       = 0.
224 
225       grid%albedo(i,j)   = 0.
226       grid%thc(i,j)      = 1000.
227       grid%znt(i,j)      = 0.01
228       grid%emiss(i,j)    = 1.
229       grid%ivgtyp(i,j)   = 1
230       grid%lu_index(i,j) = REAL(ivgtyp(i,j))
231       grid%xland(i,j)    = 1.
232       grid%mavail(i,j)   = 0.
233    END DO
234    END DO
235 
236    grid%dx = dlam*degrad/reradius
237    grid%dy = dphi*degrad/reradius
238    grid%rdx = 1./grid%dx
239    grid%rdy = 1./grid%dy
240 
241    !WRITE(*,*) ''
242    !WRITE(*,'(A,1PG14.6,A,1PG14.6)') ' For the namelist: dx =',grid%dx,', dy =',grid%dy
243 
244    CALL nl_set_mminlu(1,'    ')
245    grid%iswater = 0
246    grid%cen_lat = 0.
247    grid%cen_lon = 0.
248    grid%truelat1 = 0.
249    grid%truelat2 = 0.
250    grid%moad_cen_lat = 0.
251    grid%stand_lon    = 0.
252    ! Apparently, map projection 0 is "none" which actually turns out to be
253    ! a regular grid of latitudes and longitudes, the simple cylindrical projection
254    grid%map_proj = 0
255 
256    DO k = kds, kde
257       grid%em_znw(k) = 1. - REAL(k-kds)/REAL(kde-kds)
258    END DO
259 
260    DO k=1, kde-1
261     grid%em_dnw(k) = grid%em_znw(k+1) - grid%em_znw(k)
262     grid%em_rdnw(k) = 1./grid%em_dnw(k)
263     grid%em_znu(k) = 0.5*(grid%em_znw(k+1)+grid%em_znw(k))
264    ENDDO
265    DO k=2, kde-1
266     grid%em_dn(k) = 0.5*(grid%em_dnw(k)+grid%em_dnw(k-1))
267     grid%em_rdn(k) = 1./grid%em_dn(k)
268     grid%em_fnp(k) = .5* grid%em_dnw(k  )/grid%em_dn(k)
269     grid%em_fnm(k) = .5* grid%em_dnw(k-1)/grid%em_dn(k)
270    ENDDO
271 
272    cof1 = (2.*grid%em_dn(2)+grid%em_dn(3))/(grid%em_dn(2)+grid%em_dn(3))*grid%em_dnw(1)/grid%em_dn(2) 
273    cof2 =     grid%em_dn(2)        /(grid%em_dn(2)+grid%em_dn(3))*grid%em_dnw(1)/grid%em_dn(3) 
274    grid%cf1  = grid%em_fnp(2) + cof1
275    grid%cf2  = grid%em_fnm(2) - cof1 - cof2
276    grid%cf3  = cof2       
277 
278    grid%cfn  = (.5*grid%em_dnw(kde-1)+grid%em_dn(kde-1))/grid%em_dn(kde-1)
279    grid%cfn1 = -.5*grid%em_dnw(kde-1)/grid%em_dn(kde-1)
280 
281 
282    ! Need to add perturbations to initial profile.  Set up random number
283    ! seed here.
284    CALL random_seed
285 
286    ! General assumption from here after is that the initial temperature
287    ! profile is isothermal at a value of T0, and the initial winds are
288    ! all 0.
289 
290    ! find ptop for the desired ztop (ztop is input from the namelist)
291    grid%p_top =  p0 * EXP(-(g*config_flags%ztop)/(r_d*T0))
292 
293 
294    ! Values of geopotential (base, perturbation, and at p0) at the surface
295    DO j = jts, jte
296    DO i = its, ite
297       grid%em_phb(i,1,j) = grid%ht(i,j)*g
298       grid%em_php(i,1,j) = 0.         ! This is perturbation geopotential
299                               ! Since this is an initial condition, there
300                               ! should be no perturbation!
301       grid%em_ph0(i,1,j) = grid%ht(i,j)*g
302    ENDDO
303    ENDDO
304 
305 
306    DO J = jts, jte
307    DO I = its, ite
308 
309       p_surf = p0 * EXP(-(g*grid%em_phb(i,1,j)/g)/(r_d*T0))
310       grid%em_mub(i,j) = p_surf-grid%p_top
311 
312       ! given p (coordinate), calculate theta and compute 1/rho from equation
313       ! of state
314 
315       DO K = kts, kte-1
316          p_level = grid%em_znu(k)*(p_surf - grid%p_top) + grid%p_top
317          grid%em_pb(i,k,j) = p_level
318 
319          grid%em_t_init(i,k,j) = T0*(p0/p_level)**rcp
320          grid%em_t_init(i,k,j) = grid%em_t_init(i,k,j) - t0
321 
322          grid%em_alb(i,k,j)=(r_d/p1000mb)*(grid%em_t_init(i,k,j)+t0)*(grid%em_pb(i,k,j)/p1000mb)**cvpm
323       END DO
324 
325       ! calculate hydrostatic balance (alternatively we could interpolate
326       ! the geopotential from the sounding, but this assures that the base
327       ! state is in exact hydrostatic balance with respect to the model eqns.
328 
329       DO k = kts+1, kte
330          grid%em_phb(i,k,j) = grid%em_phb(i,k-1,j) - grid%em_dnw(k-1)*grid%em_mub(i,j)*grid%em_alb(i,k-1,j)
331       ENDDO
332 
333    ENDDO
334    ENDDO
335 
336    DO im = PARAM_FIRST_SCALAR, num_moist
337    DO J = jts, jte
338    DO K = kts, kte-1
339    DO I = its, ite
340       grid%moist(i,k,j,im) = 0.
341    END DO
342    END DO
343    END DO
344    END DO
345 
346    ! Now calculate the full (hydrostatically-balanced) state for each column
347    ! We will include moisture
348    DO J = jts, jte
349    DO I = its, ite
350 
351       ! At this point p_top is already set. find the DRY mass in the column
352       pd_surf = p0 * EXP(-(g*grid%em_phb(i,1,j)/g)/(r_d*T0))
353 
354       ! compute the perturbation mass (mu/mu_1/mu_2) and the full mass
355       grid%em_mu_1(i,j) = pd_surf-grid%p_top - grid%em_mub(i,j)
356       grid%em_mu_2(i,j) = grid%em_mu_1(i,j)
357       grid%em_mu0(i,j)  = grid%em_mu_1(i,j) + grid%em_mub(i,j)
358 
359       ! given the dry pressure and coordinate system, calculate the
360       ! perturbation potential temperature (t/t_1/t_2)
361 
362       DO k = kds, kde-1
363          p_level = grid%em_znu(k)*(pd_surf - grid%p_top) + grid%p_top
364          grid%em_t_1(i,k,j) = T0*(p0/p_level)**rcp
365          ! Add a small perturbation to initial isothermal profile
366          CALL random_number(tperturb)
367          grid%em_t_1(i,k,j)=grid%em_t_1(i,k,j)*(1.0+0.004*(tperturb-0.5))
368          grid%em_t_1(i,k,j) = grid%em_t_1(i,k,j)-t0
369          grid%em_t_2(i,k,j) = grid%em_t_1(i,k,j)
370       END DO
371 
372 
373       ! integrate the hydrostatic equation (from the RHS of the bigstep
374       ! vertical momentum equation) down from the top to get p.
375       ! first from the top of the model to the top pressure
376 
377       k = kte-1  ! top level
378 
379       qvf1 = 0.5*(grid%moist(i,k,j,P_QV)+grid%moist(i,k,j,P_QV))
380       qvf2 = 1./(1.+qvf1)
381       qvf1 = qvf1*qvf2
382 
383       ! grid%em_p(i,k,j) = - 0.5*grid%em_mu_1(i,j)/grid%em_rdnw(k)
384       grid%em_p(i,k,j) = - 0.5*(grid%em_mu_1(i,j)+qvf1*grid%em_mub(i,j))/grid%em_rdnw(k)/qvf2
385       qvf = 1. + rvovrd*grid%moist(i,k,j,P_QV)
386       grid%em_alt(i,k,j) = (r_d/p1000mb)*(grid%em_t_1(i,k,j)+t0)*qvf* &
387                   (((grid%em_p(i,k,j)+grid%em_pb(i,k,j))/p1000mb)**cvpm)
388       grid%em_al(i,k,j) = grid%em_alt(i,k,j) - grid%em_alb(i,k,j)
389 
390       !  down the column
391 
392       do k=kte-2,kts,-1
393          qvf1 = 0.5*(grid%moist(i,k,j,P_QV)+grid%moist(i,k+1,j,P_QV))
394          qvf2 = 1./(1.+qvf1)
395          qvf1 = qvf1*qvf2
396          grid%em_p(i,k,j) = grid%em_p(i,k+1,j) - (grid%em_mu_1(i,j) + qvf1*grid%em_mub(i,j))/qvf2/grid%em_rdn(k+1)
397          qvf = 1. + rvovrd*grid%moist(i,k,j,P_QV)
398          grid%em_alt(i,k,j) = (r_d/p1000mb)*(grid%em_t_1(i,k,j)+t0)*qvf* &
399                      (((grid%em_p(i,k,j)+grid%em_pb(i,k,j))/p1000mb)**cvpm)
400          grid%em_al(i,k,j) = grid%em_alt(i,k,j) - grid%em_alb(i,k,j)
401       enddo
402 
403       ! this is the hydrostatic equation used in the model after the
404       ! small timesteps.  In the model, al (inverse density)
405       ! is computed from the geopotential.
406 
407       grid%em_ph_1(i,1,j) = 0.
408       DO k  = kts+1,kte
409          grid%em_ph_1(i,k,j) = grid%em_ph_1(i,k-1,j) - (1./grid%em_rdnw(k-1))*(  &
410                       (grid%em_mub(i,j)+grid%em_mu_1(i,j))*grid%em_al(i,k-1,j)+  &
411                        grid%em_mu_1(i,j)*grid%em_alb(i,k-1,j)  )
412 
413          grid%em_ph_2(i,k,j) = grid%em_ph_1(i,k,j)
414          grid%em_ph0(i,k,j) = grid%em_ph_1(i,k,j) + grid%em_phb(i,k,j)
415       ENDDO
416 
417    END DO
418    END DO
419 
420 
421 
422    ! Now set U & V
423    DO J = jts, jte
424    DO K = kts, kte-1
425    DO I = its, ite
426       grid%em_u_1(i,k,j) = 0.
427       grid%em_u_2(i,k,j) = 0.
428       grid%em_v_1(i,k,j) = 0.
429       grid%em_v_2(i,k,j) = 0.
430    END DO
431    END DO
432    END DO
433 
434 
435    DO j=jts, jte
436    DO k=kds, kde
437    DO i=its, ite
438       grid%em_ww(i,k,j)  = 0.
439       grid%em_w_1(i,k,j) = 0.
440       grid%em_w_2(i,k,j) = 0.
441       grid%h_diabatic(i,k,j) = 0.
442    END DO
443    END DO
444    END DO
445 
446   
447    DO k=kts,kte
448       grid%em_t_base(k)  = grid%em_t_init(its,k,jts)
449       grid%qv_base(k) = 0.
450       grid%u_base(k)  = 0.
451       grid%v_base(k)  = 0.
452    END DO
453 
454    ! One subsurface layer: infinite slab at constant temperature below
455    ! the surface.  Surface temperature is an infinitely thin "skin" on
456    ! top of a half-infinite slab.  The temperature of both the skin and
457    ! the slab are determined from the initial nearest-surface-air-layer
458    ! temperature.
459    DO J = jts, MIN(jte, jde-1)
460    DO I = its, MIN(ite, ide-1)
461       thtmp   = grid%em_t_2(i,1,j)+t0
462       ptmp    = grid%em_p(i,1,j)+grid%em_pb(i,1,j)
463       temp(1) = thtmp * (ptmp/p1000mb)**rcp
464       thtmp   = grid%em_t_2(i,2,j)+t0
465       ptmp    = grid%em_p(i,2,j)+grid%em_pb(i,2,j)
466       temp(2) = thtmp * (ptmp/p1000mb)**rcp
467       thtmp   = grid%em_t_2(i,3,j)+t0
468       ptmp    = grid%em_p(i,3,j)+grid%em_pb(i,3,j)
469       temp(3) = thtmp * (ptmp/p1000mb)**rcp
470       grid%tsk(I,J)=cf1*temp(1)+cf2*temp(2)+cf3*temp(3)
471       grid%tmn(I,J)=grid%tsk(I,J)-0.5
472    END DO
473    END DO
474 
475    RETURN
476 
477  END SUBROUTINE init_domain_rk
478 
479 !---------------------------------------------------------------------
480 
481  SUBROUTINE init_module_initialize
482  END SUBROUTINE init_module_initialize
483 
484 !---------------------------------------------------------------------
485 
486 END MODULE module_initialize_ideal