Merge remote-tracking branch 'origin/release-v4.5.2'
[WRF.git] / dyn_em / module_initialize_heldsuarez.F
blobc63dc23ce46efe7e601a32081fe8fde50b95f72a
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    :: B1, B2, B3, B4, B5
81    REAL :: p_surf, p_level, pd_surf, qvf1, qvf2, qvf
82    REAL :: thtmp, ptmp, temp(3), cof1, cof2
84    INTEGER :: icm,jcm
86    SELECT CASE ( model_data_order )
87          CASE ( DATA_ORDER_ZXY )
88    kds = grid%sd31 ; kde = grid%ed31 ;
89    ids = grid%sd32 ; ide = grid%ed32 ;
90    jds = grid%sd33 ; jde = grid%ed33 ;
92    kms = grid%sm31 ; kme = grid%em31 ;
93    ims = grid%sm32 ; ime = grid%em32 ;
94    jms = grid%sm33 ; jme = grid%em33 ;
96    kts = grid%sp31 ; kte = grid%ep31 ;   ! note that tile is entire patch
97    its = grid%sp32 ; ite = grid%ep32 ;   ! note that tile is entire patch
98    jts = grid%sp33 ; jte = grid%ep33 ;   ! note that tile is entire patch
99          CASE ( DATA_ORDER_XYZ )
100    ids = grid%sd31 ; ide = grid%ed31 ;
101    jds = grid%sd32 ; jde = grid%ed32 ;
102    kds = grid%sd33 ; kde = grid%ed33 ;
104    ims = grid%sm31 ; ime = grid%em31 ;
105    jms = grid%sm32 ; jme = grid%em32 ;
106    kms = grid%sm33 ; kme = grid%em33 ;
108    its = grid%sp31 ; ite = grid%ep31 ;   ! note that tile is entire patch
109    jts = grid%sp32 ; jte = grid%ep32 ;   ! note that tile is entire patch
110    kts = grid%sp33 ; kte = grid%ep33 ;   ! note that tile is entire patch
111          CASE ( DATA_ORDER_XZY )
112    ids = grid%sd31 ; ide = grid%ed31 ;
113    kds = grid%sd32 ; kde = grid%ed32 ;
114    jds = grid%sd33 ; jde = grid%ed33 ;
116    ims = grid%sm31 ; ime = grid%em31 ;
117    kms = grid%sm32 ; kme = grid%em32 ;
118    jms = grid%sm33 ; jme = grid%em33 ;
120    its = grid%sp31 ; ite = grid%ep31 ;   ! note that tile is entire patch
121    kts = grid%sp32 ; kte = grid%ep32 ;   ! note that tile is entire patch
122    jts = grid%sp33 ; jte = grid%ep33 ;   ! note that tile is entire patch
124    END SELECT
126    CALL model_to_grid_config_rec ( grid%id , model_config_rec , config_flags )
128 ! here we check to see if the boundary conditions are set properly
130    CALL boundary_condition_check( config_flags, bdyzone, error, grid%id )
132     grid%itimestep=0
133    grid%step_number = 0
136 #ifdef DM_PARALLEL
137    CALL wrf_dm_bcast_bytes( icm , IWORDSIZE )
138    CALL wrf_dm_bcast_bytes( jcm , IWORDSIZE )
139 #endif
141    ! Initialize 2D surface arrays
143    nxx = ide-ids ! Don't include u-stagger
144    nyy = jde-jds ! Don't include v-stagger
145    dphi = 180./REAL(nyy)
146    dlam = 360./REAL(nxx)
148    DO j = jts, jte
149    DO i = its, ite
150       ! ig is the I index in the global (domain) span of the array.
151       ! jg is the J index in the global (domain) span of the array.
152       ig = i - ids + 1  ! ids is not necessarily 1
153       jg = j - jds + 1  ! jds is not necessarily 1
155       grid%xlat(i,j)  = (REAL(jg)-0.5)*dphi-90.
156       grid%xlong(i,j) = (REAL(ig)-0.5)*dlam-180.
157       vlat       = grid%xlat(i,j) - 0.5*dphi
159       grid%clat(i,j) = grid%xlat(i,j)
161       grid%msftx(i,j) = 1./COS(grid%xlat(i,j)*degrad)
162       grid%msfty(i,j) = 1.
163       grid%msfux(i,j) = 1./COS(grid%xlat(i,j)*degrad)
164       grid%msfuy(i,j) = 1.
165       grid%e(i,j)     = 2*EOMEG*COS(grid%xlat(i,j)*degrad)
166       grid%f(i,j)     = 2*EOMEG*SIN(grid%xlat(i,j)*degrad)
168       ! The following two are the cosine and sine of the rotation
169       ! of projection.  Simple cylindrical is *simple* ... no rotation!
170       grid%sina(i,j) = 0.
171       grid%cosa(i,j) = 1.
173    END DO
174    END DO
176 !   DO j = max(jds+1,jts), min(jde-1,jte)
177    DO j = jts, jte
178    DO i = its, ite
179       vlat       = grid%xlat(i,j) - 0.5*dphi
180       grid%msfvx(i,j) = 1./COS(vlat*degrad)
181       grid%msfvy(i,j) = 1.
182       grid%msfvx_inv(i,j) = 1./grid%msfvx(i,j)
183    END DO
184    END DO
186    IF(jts == jds) THEN
187    DO i = its, ite
188       grid%msfvx(i,jts) = 00.
189       grid%msfvx_inv(i,jts) = 0.
190    END DO
191    END IF
193    IF(jte == jde) THEN
194    DO i = its, ite
195       grid%msfvx(i,jte) = 00.
196       grid%msfvx_inv(i,jte) = 0.
197    END DO
198    END IF
200    DO j=jts,jte
201      vlat       = grid%xlat(its,j) - 0.5*dphi
202      write(6,*) j,vlat,grid%msfvx(its,j),grid%msfvx_inv(its,j)
203    ENDDO
205    DO j=jts,jte
206    DO i=its,ite
207       grid%ht(i,j)       = 0.
209       grid%albedo(i,j)   = 0.
210       grid%thc(i,j)      = 1000.
211       grid%znt(i,j)      = 0.01
212       grid%emiss(i,j)    = 1.
213       grid%ivgtyp(i,j)   = 1
214       grid%lu_index(i,j) = REAL(ivgtyp(i,j))
215       grid%xland(i,j)    = 1.
216       grid%mavail(i,j)   = 0.
217    END DO
218    END DO
220    grid%dx = dlam*degrad/reradius
221    grid%dy = dphi*degrad/reradius
222    grid%rdx = 1./grid%dx
223    grid%rdy = 1./grid%dy
225    !WRITE(*,*) ''
226    !WRITE(*,'(A,1PG14.6,A,1PG14.6)') ' For the namelist: dx =',grid%dx,', dy =',grid%dy
228    CALL nl_set_mminlu(1,'    ')
229    grid%iswater = 0
230    grid%cen_lat = 0.
231    grid%cen_lon = 0.
232    grid%truelat1 = 0.
233    grid%truelat2 = 0.
234    grid%moad_cen_lat = 0.
235    grid%stand_lon    = 0.
236    grid%pole_lon = 0.
237    grid%pole_lat = 90.
238    ! Apparently, map projection 0 is "none" which actually turns out to be
239    ! a regular grid of latitudes and longitudes, the simple cylindrical projection
240    grid%map_proj = 0
242    DO k = kds, kde
243       grid%znw(k) = 1. - REAL(k-kds)/REAL(kde-kds)
244    END DO
246    DO k=1, kde-1
247     grid%dnw(k) = grid%znw(k+1) - grid%znw(k)
248     grid%rdnw(k) = 1./grid%dnw(k)
249     grid%znu(k) = 0.5*(grid%znw(k+1)+grid%znw(k))
250    ENDDO
252    DO k=2, kde-1
253     grid%dn(k) = 0.5*(grid%dnw(k)+grid%dnw(k-1))
254     grid%rdn(k) = 1./grid%dn(k)
255     grid%fnp(k) = .5* grid%dnw(k  )/grid%dn(k)
256     grid%fnm(k) = .5* grid%dnw(k-1)/grid%dn(k)
257    ENDDO
259    cof1 = (2.*grid%dn(2)+grid%dn(3))/(grid%dn(2)+grid%dn(3))*grid%dnw(1)/grid%dn(2) 
260    cof2 =     grid%dn(2)        /(grid%dn(2)+grid%dn(3))*grid%dnw(1)/grid%dn(3) 
261    grid%cf1  = grid%fnp(2) + cof1
262    grid%cf2  = grid%fnm(2) - cof1 - cof2
263    grid%cf3  = cof2       
265    grid%cfn  = (.5*grid%dnw(kde-1)+grid%dn(kde-1))/grid%dn(kde-1)
266    grid%cfn1 = -.5*grid%dnw(kde-1)/grid%dn(kde-1)
269    ! Need to add perturbations to initial profile.  Set up random number
270    ! seed here.
271    CALL random_seed
273    ! General assumption from here after is that the initial temperature
274    ! profile is isothermal at a value of T0, and the initial winds are
275    ! all 0.
277    ! find ptop for the desired ztop (ztop is input from the namelist)
278    grid%p_top =  p0 * EXP(-(g*config_flags%ztop)/(r_d*T0))
280    !  For hybrid coord
282    DO k=kts, kte
283       IF      ( config_flags%hybrid_opt .EQ. 0 ) THEN
284          grid%c3f(k) = grid%znw(k)
285       ELSE IF ( config_flags%hybrid_opt .EQ. 1 ) THEN
286          grid%c3f(k) = grid%znw(k)
287       ELSE IF ( config_flags%hybrid_opt .EQ. 2 ) THEN
288          B1 = 2. * grid%etac**2 * ( 1. - grid%etac )
289          B2 = -grid%etac * ( 4. - 3. * grid%etac - grid%etac**3 )
290          B3 = 2. * ( 1. - grid%etac**3 )
291          B4 = - ( 1. - grid%etac**2 )
292          B5 = (1.-grid%etac)**4
293          grid%c3f(k) = ( B1 + B2*grid%znw(k) + B3*grid%znw(k)**2 + B4*grid%znw(k)**3 ) / B5
294          IF ( grid%znw(k) .LT. grid%etac ) THEN
295             grid%c3f(k) = 0.
296          END IF
297          IF ( k .EQ. kds ) THEN
298             grid%c3f(k) = 1.
299          ELSE IF ( k .EQ. kde ) THEN
300             grid%c3f(k) = 0.
301          END IF
302       ELSE IF ( config_flags%hybrid_opt .EQ. 3 ) THEN
303          grid%c3f(k) = grid%znw(k)*sin(0.5*3.14159*grid%znw(k))**2
304          IF      ( k .EQ. kds ) THEN
305             grid%c3f(k) = 1.
306          ELSE IF ( k .EQ. kds ) THEN
307             grid%c3f(kde) = 0.
308          END IF
309       ELSE
310          CALL wrf_message     ( 'ERROR: --- hybrid_opt' )
311          CALL wrf_message     ( 'ERROR: --- hybrid_opt=0    ==> Standard WRF terrain-following coordinate' )
312          CALL wrf_message     ( 'ERROR: --- hybrid_opt=1    ==> Standard WRF terrain-following coordinate, hybrid c1, c2, c3, c4' )
313          CALL wrf_message     ( 'ERROR: --- hybrid_opt=2    ==> Hybrid, Klemp polynomial' )
314          CALL wrf_message     ( 'ERROR: --- hybrid_opt=3    ==> Hybrid, sin^2' )
315          CALL wrf_error_fatal ( 'ERROR: --- Invalid option' )
316       END IF
317    END DO
319    !  c4 is a function of c3 and eta.
321    DO k=1, kde
322       grid%c4f(k) = ( grid%znw(k) - grid%c3f(k) ) * ( p1000mb - grid%p_top )
323    ENDDO
324        
325    !  Now on half levels, just add up and divide by 2 (for c3h).  Use (eta-c3)*(p00-pt) for c4 on half levels.
326        
327    DO k=1, kde-1
328       grid%znu(k) = ( grid%znw(k+1) + grid%znw(k) ) * 0.5
329       grid%c3h(k) = ( grid%c3f(k+1) + grid%c3f(k) ) * 0.5
330       grid%c4h(k) = ( grid%znu(k) - grid%c3h(k) ) * ( p1000mb - grid%p_top )
331    ENDDO
332        
333    !  c1 = d(B)/d(eta).  We define c1f as c1 on FULL levels.  For a vertical difference,
334    !  we need to use B and eta on half levels.  The k-loop ends up referring to the
335    !  full levels, neglecting the top and bottom.
336        
337    DO k=kds+1, kde-1
338       grid%c1f(k) = ( grid%c3h(k) - grid%c3h(k-1) ) / ( grid%znu(k) - grid%znu(k-1) )
339    ENDDO
340        
341    !  The boundary conditions to get the coefficients:
342    !  1) At k=kts: define d(B)/d(eta) = 1.  This gives us the same value of B and d(B)/d(eta)
343    !     when doing the sigma-only B=eta.
344    !  2) At k=kte: define d(B)/d(eta) = 0.  The curve B SMOOTHLY goes to zero, and at the very
345    !     top, B continues to SMOOTHLY go to zero.  Note that for almost all cases of non B=eta,
346    !     B is ALREADY=ZERO at the top, so this is a reasonable BC to assume.
347        
348    grid%c1f(kds) = 1.
349    IF ( ( config_flags%hybrid_opt .EQ. 0 ) .OR. ( config_flags%hybrid_opt .EQ. 1 ) ) THEN
350       grid%c1f(kde) = 1.
351    ELSE
352       grid%c1f(kde) = 0.
353    END IF
354        
355    !  c2 = ( 1. - c1(k) ) * (p00 - pt).  There is no vertical differencing, so we can do the
356    !  full kds to kde looping.
357        
358    DO k=kds, kde
359       grid%c2f(k) = ( 1. - grid%c1f(k) ) * ( p1000mb - grid%p_top )
360    ENDDO
361        
362    !  Now on half levels for c1 and c2.  The c1h will result from the full level c3 and full
363    !  level eta differences.  The c2 value use the half level c1(k).
364        
365    DO k=1, kde-1
366       grid%c1h(k) = ( grid%c3f(k+1) - grid%c3f(k) ) / ( grid%znw(k+1) - grid%znw(k) )
367       grid%c2h(k) = ( 1. - grid%c1h(k) ) * ( p1000mb - grid%p_top )
368    ENDDO
370    ! Values of geopotential (base, perturbation, and at p0) at the surface
371    DO j = jts, jte
372    DO i = its, ite
373       grid%phb(i,1,j) = grid%ht(i,j)*g
374       grid%php(i,1,j) = 0.         ! This is perturbation geopotential
375                               ! Since this is an initial condition, there
376                               ! should be no perturbation!
377       grid%ph0(i,1,j) = grid%ht(i,j)*g
378    ENDDO
379    ENDDO
382    DO J = jts, jte
383    DO I = its, ite
385       p_surf = p0 * EXP(-(g*grid%phb(i,1,j)/g)/(r_d*T0))
386       grid%mub(i,j) = p_surf-grid%p_top
388       ! given p (coordinate), calculate theta and compute 1/rho from equation
389       ! of state
391       DO K = kts, kte-1
392          p_level = grid%c3h(k)*(p_surf - grid%p_top)+grid%c4h(k) + grid%p_top
393          grid%pb(i,k,j) = p_level
395          grid%t_init(i,k,j) = T0*(p0/p_level)**rcp
396          grid%t_init(i,k,j) = grid%t_init(i,k,j) - t0
398          grid%alb(i,k,j)=(r_d/p1000mb)*(grid%t_init(i,k,j)+t0)*(grid%pb(i,k,j)/p1000mb)**cvpm
399       END DO
401       ! calculate hydrostatic balance (alternatively we could interpolate
402       ! the geopotential from the sounding, but this assures that the base
403       ! state is in exact hydrostatic balance with respect to the model eqns.
405       DO k = kts+1, kte
406          grid%phb(i,k,j) = grid%phb(i,k-1,j) - grid%dnw(k-1)*(grid%c1h(k-1)*grid%mub(i,j)+grid%c2h(k-1))*grid%alb(i,k-1,j)
407       ENDDO
409    ENDDO
410    ENDDO
412    DO im = PARAM_FIRST_SCALAR, num_moist
413    DO J = jts, jte
414    DO K = kts, kte-1
415    DO I = its, ite
416       grid%moist(i,k,j,im) = 0.
417    END DO
418    END DO
419    END DO
420    END DO
422    ! Now calculate the full (hydrostatically-balanced) state for each column
423    ! We will include moisture
424    DO J = jts, jte
425    DO I = its, ite
427       ! At this point p_top is already set. find the DRY mass in the column
428       pd_surf = p0 * EXP(-(g*grid%phb(i,1,j)/g)/(r_d*T0))
430       ! compute the perturbation mass (mu/mu_1/mu_2) and the full mass
431       grid%mu_1(i,j) = pd_surf-grid%p_top - grid%mub(i,j)
432       grid%mu_2(i,j) = grid%mu_1(i,j)
433       grid%mu0(i,j)  = grid%mu_1(i,j) + grid%mub(i,j)
435       ! given the dry pressure and coordinate system, calculate the
436       ! perturbation potential temperature (t/t_1/t_2)
438       DO k = kds, kde-1
439          p_level = grid%c3h(k)*(pd_surf - grid%p_top)+grid%c4h(k) + grid%p_top
440          grid%t_1(i,k,j) = T0*(p0/p_level)**rcp
441          ! Add a small perturbation to initial isothermal profile
442          CALL random_number(tperturb)
443          grid%t_1(i,k,j)=grid%t_1(i,k,j)*(1.0+0.004*(tperturb-0.5))
444          grid%t_1(i,k,j) = grid%t_1(i,k,j)-t0
445          grid%t_2(i,k,j) = grid%t_1(i,k,j)
446       END DO
449       ! integrate the hydrostatic equation (from the RHS of the bigstep
450       ! vertical momentum equation) down from the top to get p.
451       ! first from the top of the model to the top pressure
453       k = kte-1  ! top level
455       qvf1 = 0.5*(grid%moist(i,k,j,P_QV)+grid%moist(i,k,j,P_QV))
456       qvf2 = 1./(1.+qvf1)
457       qvf1 = qvf1*qvf2
459       grid%p(i,k,j) = - 0.5*((grid%c1f(k+1)*grid%mu_1(i,j))+qvf1*(grid%c1f(k+1)*grid%mub(i,j)+grid%c2f(k+1)))/grid%rdnw(k)/qvf2
460       qvf = 1. + rvovrd*grid%moist(i,k,j,P_QV)
461       grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* &
462                   (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm)
463       grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j)
465       !  down the column
467       do k=kte-2,kts,-1
468          qvf1 = 0.5*(grid%moist(i,k,j,P_QV)+grid%moist(i,k+1,j,P_QV))
469          qvf2 = 1./(1.+qvf1)
470          qvf1 = qvf1*qvf2
471          grid%p(i,k,j) = grid%p(i,k+1,j) - ((grid%c1f(k+1)*grid%mu_1(i,j)) + qvf1*(grid%c1f(k+1)*grid%mub(i,j)+grid%c2f(k+1)))/qvf2/grid%rdn(k+1)
472          qvf = 1. + rvovrd*grid%moist(i,k,j,P_QV)
473          grid%alt(i,k,j) = (r_d/p1000mb)*(grid%t_1(i,k,j)+t0)*qvf* &
474                      (((grid%p(i,k,j)+grid%pb(i,k,j))/p1000mb)**cvpm)
475          grid%al(i,k,j) = grid%alt(i,k,j) - grid%alb(i,k,j)
476       enddo
478       ! this is the hydrostatic equation used in the model after the
479       ! small timesteps.  In the model, al (inverse density)
480       ! is computed from the geopotential.
482       grid%ph_1(i,1,j) = 0.
483       DO k  = kts+1,kte
484          grid%ph_1(i,k,j) = grid%ph_1(i,k-1,j) - (grid%dnw(k-1))*(       &    
485                       ((grid%c1h(k-1)*grid%mub(i,j)+grid%c2h(k-1))+(grid%c1h(k-1)*grid%mu_1(i,j)))*grid%al(i,k-1,j)+ &
486                        (grid%c1h(k-1)*grid%mu_1(i,j))*grid%alb(i,k-1,j)  )
488          grid%ph_2(i,k,j) = grid%ph_1(i,k,j)
489          grid%ph0(i,k,j) = grid%ph_1(i,k,j) + grid%phb(i,k,j)
490       ENDDO
492    END DO
493    END DO
497    ! Now set U & V
498    DO J = jts, jte
499    DO K = kts, kte-1
500    DO I = its, ite
501       grid%u_1(i,k,j) = 0.
502       grid%u_2(i,k,j) = 0.
503       grid%v_1(i,k,j) = 0.
504       grid%v_2(i,k,j) = 0.
505    END DO
506    END DO
507    END DO
510    DO j=jts, jte
511    DO k=kds, kde
512    DO i=its, ite
513       grid%ww(i,k,j)  = 0.
514       grid%w_1(i,k,j) = 0.
515       grid%w_2(i,k,j) = 0.
516       grid%h_diabatic(i,k,j) = 0.
517    END DO
518    END DO
519    END DO
521   
522    DO k=kts,kte
523       grid%t_base(k)  = grid%t_init(its,k,jts)
524       grid%qv_base(k) = 0.
525       grid%u_base(k)  = 0.
526       grid%v_base(k)  = 0.
527    END DO
529    ! One subsurface layer: infinite slab at constant temperature below
530    ! the surface.  Surface temperature is an infinitely thin "skin" on
531    ! top of a half-infinite slab.  The temperature of both the skin and
532    ! the slab are determined from the initial nearest-surface-air-layer
533    ! temperature.
534    DO J = jts, MIN(jte, jde-1)
535    DO I = its, MIN(ite, ide-1)
536       thtmp   = grid%t_2(i,1,j)+t0
537       ptmp    = grid%p(i,1,j)+grid%pb(i,1,j)
538       temp(1) = thtmp * (ptmp/p1000mb)**rcp
539       thtmp   = grid%t_2(i,2,j)+t0
540       ptmp    = grid%p(i,2,j)+grid%pb(i,2,j)
541       temp(2) = thtmp * (ptmp/p1000mb)**rcp
542       thtmp   = grid%t_2(i,3,j)+t0
543       ptmp    = grid%p(i,3,j)+grid%pb(i,3,j)
544       temp(3) = thtmp * (ptmp/p1000mb)**rcp
545       grid%tsk(I,J)=cf1*temp(1)+cf2*temp(2)+cf3*temp(3)
546       grid%tmn(I,J)=grid%tsk(I,J)-0.5
547    END DO
548    END DO
550       !  Save the dry perturbation potential temperature.
552       DO j = jts, min(jde-1,jte)
553          DO k = kts, kte
554             DO i = its, min(ide-1,ite)
555                grid%th_phy_m_t0(i,k,j) = grid%t_2(i,k,j)
556             END DO
557          END DO
558       END DO
560    !  Turn dry potential temperature into moist potential temperature
561       !  at the very end of this routine
562       !  This field will be in the model IC and and used to construct the 
563       !  BC file.
565       IF ( ( config_flags%use_theta_m .EQ. 1 ) .AND. (P_Qv .GE. PARAM_FIRST_SCALAR) ) THEN
566       DO J  = jts, min(jde-1,jte)
567          DO K = kts, kte-1
568             DO I = its, min(ide-1,ite)
569                grid%t_2(i,k,j) = ( grid%t_2(i,k,j) + T0 ) * (1. + (R_v/R_d) * moist(i,k,j,p_qv)) - T0
570             END DO
571          END DO
572       END DO
573       ENDIF
575    RETURN
577  END SUBROUTINE init_domain_rk
579 !---------------------------------------------------------------------
581  SUBROUTINE init_module_initialize
582  END SUBROUTINE init_module_initialize
584 !---------------------------------------------------------------------
586 END MODULE module_initialize_ideal