source: lmdz_wrf/trunk/WRFV3/phys/module_fr_sfire_util.F

Last change on this file was 1, checked in by lfita, 10 years ago
  • -- --- Opening of the WRF+LMDZ coupling repository --- -- -

WRF: version v3.3
LMDZ: version v1818

More details in:

File size: 41.4 KB
Line 
1!
2! This module contains general purpose utilities and WRF wrappers because some
3! applications require this module to run standalone. No physics here.
4! Some are dependent on WRF indexing scheme. Some violate WRF conventions but these
5! are not called from the WRF dependent code. Some are not called at all.
6!
7
8module module_fr_sfire_util
9
10! various method selection parameters
11! 1. add the parameter and its static default here
12! optional:
13! 2. add copy from config_flags in module_fr_sfire_driver%%set_flags
14! 3. add a line in Registry.EM to define the variable and set default value
15! 4. add the variable and value in namelist.input
16
17
18integer,save::          &
19 fire_print_msg=1,      &  ! print SFIRE progress
20 fire_print_file=1,     &  ! write many files by write_array_m; compile with DEBUG_OUT, do not run in parallel
21 fuel_left_method=1,    &  ! 1=simple, 2=exact in linear case
22 fuel_left_irl=2,       &  ! refinement for fuel calculation, must be even
23 fuel_left_jrl=2,       &
24 boundary_guard=-1,     &  ! crash if fire gets this many cells to domain boundary, -1=off
25 fire_grows_only=1,     &  ! fire can spread out only (level set functions may not increase)
26 fire_upwinding=3,      &  ! upwind normal spread: 1=standard, 2=godunov, 3=eno, 4=sethian
27 fire_upwind_split=0,   &  ! 1=upwind advection separately from normal direction spread
28 fire_test_steps=0,     &  ! 0=no fire, 1=normal, >1 = do specified number of steps and terminate (testing only)
29 fire_topo_from_atm=1,  &  ! 0 = expect ZSF set correctly on entry, 1 = populate by interploating from atmosphere
30 fire_advection=0          ! 0 = fire spread from normal wind/slope (CAWFE), 1 = full speed projected
31 
32
33real, save::            &
34 fire_const_time=-1,    &  ! time from ignition to start constant heat output  <0=never
35 fire_const_grnhfx=-1., &  ! if both >=0, the constant heat flux to output then
36 fire_const_grnqfx=-1., &  ! if both >=0, the constant heat flux to output then
37 fire_atm_feedback=1. , &  ! 1 = normal, 0. = one way coupling atmosphere -> fire only
38 fire_back_weight=0.5,  &  ! RK parameter, 0 = Euler method, 0.5 = Heun, 1 = fake backward Euler
39 fire_viscosity=0.4,    &  ! artificial viscosity
40 fire_lfn_ext_up=1.        ! 0.=extend level set function at boundary by reflection, 1.=always up
41
42integer, parameter:: REAL_SUM=10, REAL_MAX=20, RNRM_SUM=30, RNRM_MAX=40
43
44contains
45
46!
47!****************
48!
49subroutine crash(s)
50use module_wrf_error
51implicit none
52character(len=*), intent(in)::s
53character(len=128)msg
54msg='crash: '//s
55call message(msg,level=0)
56!$OMP CRITICAL(SFIRE_MESSAGE_CRIT)
57call wrf_error_fatal(msg)
58!$OMP END CRITICAL(SFIRE_MESSAGE_CRIT)
59end subroutine crash
60
61!
62!****************
63!
64
65subroutine message(s,level)
66use module_wrf_error
67#ifdef _OPENMP
68use OMP_LIB
69#endif
70implicit none
71!*** arguments
72character(len=*), intent(in)::s
73integer,intent(in),optional::level
74!*** local
75character(len=128)::msg
76character(len=118)::t
77integer m,mlevel
78logical op
79!*** executable
80if(present(level))then
81    mlevel=level
82else
83    mlevel=2  ! default message level
84endif
85if(fire_print_msg.ge.mlevel)then
86      m=0
87!$OMP CRITICAL(SFIRE_MESSAGE_CRIT)
88#ifdef _OPENMP
89      m=omp_get_thread_num()
90      t=s
91      write(msg,'(a6,i3,a1,a118)')'SFIRE:',m,':',t
92#else
93      msg='SFIRE:'//s
94#endif
95      call wrf_message(msg)
96!      flush(6)
97!      flush(0)
98!$OMP END CRITICAL(SFIRE_MESSAGE_CRIT)
99endif
100end subroutine message
101
102!
103!****************
104!
105
106
107integer function open_text_file(filename,rw)
108implicit none
109character(len=*),intent(in):: filename,rw
110!$ integer, external:: OMP_GET_THREAD_NUM
111character(len=128):: msg
112character(len=1)::act
113integer::iounit,ierr
114logical::op
115
116!$  if (OMP_GET_THREAD_NUM() .ne. 0)then
117!$     call crash('open_input_text_file: called from parallel loop')
118!$  endif
119
120    do iounit=19,99
121       inquire(iounit,opened=op)
122       if(.not.op)goto 1
123    enddo
124    call crash('open_text_file: Cannot find any available I/O unit')
1251   continue
126    act=rw(1:1)
127    select case (act)
128        case ('r','R')
129            OPEN(iounit, FILE=filename,FORM='FORMATTED',STATUS='OLD',ACTION='READ',IOSTAT=ierr)
130        case ('w','W')
131            OPEN(iounit, FILE=filename,FORM='FORMATTED',STATUS='UNKNOWN',ACTION='WRITE',IOSTAT=ierr)
132        case default
133            write(msg,*)'open_text_file: bad mode ',trim(rw),' for file ',trim(filename)
134    end select
135   
136    if(ierr.ne.0)then
137        write(msg,*)'open_text_file: Cannot open file ',filename
138        call crash(msg)
139    endif
140    open_text_file=iounit
141
142end function open_text_file
143
144!
145!****************
146!
147
148
149subroutine set_ideal_coord( dxf,dyf, &
150                ifds,ifde,jfds,jfde,  &
151                ifms,ifme,jfms,jfme,  &
152                ifts,ifte,jfts,jfte,  &
153                fxlong,fxlat          &
154            )
155implicit none
156! arguments
157real, intent(in)::dxf,dyf
158integer, intent(in):: &
159                ifds,ifde,jfds,jfde,  &
160                ifms,ifme,jfms,jfme,  &
161                ifts,ifte,jfts,jfte
162real, intent(out),dimension(ifms:ifme,jfms:jfme)::fxlong,fxlat
163! local
164integer::i,j
165                ! set fake  coordinates, in m
166                do j=jfts,jfte
167                    do i=ifts,ifte
168                        ! uniform mesh, lower left domain corner is (0,0)
169                        fxlong(i,j)=(i-ifds+0.5)*dxf
170                        fxlat (i,j)=(j-jfds+0.5)*dyf
171                    enddo
172                enddo
173end subroutine set_ideal_coord
174
175!
176!****************
177!
178
179
180subroutine continue_at_boundary(ix,iy,bias, & ! do x direction or y direction
181    ims,ime,jms,jme, &                ! memory dims
182    ids,ide,jds,jde, &                ! domain dims
183    ips,ipe,jps,jpe, &                ! patch dims
184    its,ite,jts,jte, &                ! tile dims
185    itso,iteo,jtso,jteo, &            ! tile dims where set
186    lfn)                              ! array
187implicit none
188!*** description
189! extend array by one beyond the domain by linear continuation
190!*** arguments
191integer, intent(in)::ix,iy              ! not 0 = do x or y (1 or 2) direction
192real,intent(in)::bias                   ! 0=none, 1.=max
193integer, intent(in)::ims,ime,jms,jme, &                ! memory dims
194    ids,ide,jds,jde, &                ! domain dims
195    ips,ipe,jps,jpe, &                ! patch dims
196    its,ite,jts,jte                   ! tile dims
197integer, intent(out)::itso,jtso,iteo,jteo ! where set
198real,intent(inout),dimension(ims:ime,jms:jme)::lfn
199!*** local
200integer i,j
201character(len=128)::msg
202integer::its1,ite1,jts1,jte1
203integer,parameter::halo=1           ! width of halo region to update
204!*** executable
205! check if there is space for the extension
206call check_mesh_2dim(its-1,ite+1,jts-1,jte+1,ims,ime,jms,jme)
207! for dislay only
208itso=its
209jtso=jts
210iteo=ite
211jteo=jte
212! go halo width beyond if at patch boundary but not at domain boudnary
213! assume we have halo need to compute the value we do not have
214! the next thread that would conveniently computer the extended values at patch corners
215! besides halo may not transfer values outside of the domain
216!
217its1=its
218jts1=jts
219ite1=ite
220jte1=jte
221if(its.eq.ips.and..not.its.eq.ids)its1=its-halo
222if(jts.eq.jps.and..not.jts.eq.jds)jts1=jts-halo
223if(ite.eq.ipe.and..not.ite.eq.ide)ite1=ite+halo
224if(jte.eq.jpe.and..not.jte.eq.jde)jte1=jte+halo
225!$OMP CRITICAL(SFIRE_UTIL_CRIT)
226write(msg,'(a,2i5,a,f5.2)')'continue_at_boundary: directions',ix,iy,' bias ',bias
227call message(msg)
228!$OMP END CRITICAL(SFIRE_UTIL_CRIT)
229if(ix.ne.0)then
230    if(its.eq.ids)then
231        do j=jts1,jte1
232            lfn(ids-1,j)=EX(lfn(ids,j),lfn(ids+1,j))
233        enddo
234        itso=ids-1
235    endif
236    if(ite.eq.ide)then
237        do j=jts1,jte1
238            lfn(ide+1,j)=EX(lfn(ide,j),lfn(ide-1,j))
239        enddo
240        iteo=ide+1
241    endif
242!$OMP CRITICAL(SFIRE_UTIL_CRIT)
243    write(msg,'(8(a,i5))')'continue_at_boundary: x:',its,':',ite,',',jts,':',jte,' ->',itso,':',iteo,',',jts1,':',jte1
244    call message(msg)
245!$OMP END CRITICAL(SFIRE_UTIL_CRIT)
246endif
247if(iy.ne.0)then
248    if(jts.eq.jds)then
249        do i=its1,ite1
250            lfn(i,jds-1)=EX(lfn(i,jds),lfn(i,jds+1))
251        enddo
252        jtso=jds-1
253    endif
254    if(jte.eq.jde)then
255        do i=its1,ite1
256            lfn(i,jde+1)=EX(lfn(i,jde),lfn(i,jde-1))
257        enddo
258        jteo=jde+1
259    endif
260!$OMP CRITICAL(SFIRE_UTIL_CRIT)
261    write(msg,'(8(a,i5))')'continue_at_boundary: y:',its,':',ite,',',jts,':',jte,' ->',its1,':',ite1,',',jtso,':',jteo
262!$OMP END CRITICAL(SFIRE_UTIL_CRIT)
263    call message(msg)
264endif
265! corners of the domain
266if(ix.ne.0.and.iy.ne.0)then
267    if(its.eq.ids.and.jts.eq.jds)lfn(ids-1,jds-1)=EX(lfn(ids,jds),lfn(ids+1,jds+1))
268    if(its.eq.ids.and.jte.eq.jde)lfn(ids-1,jde+1)=EX(lfn(ids,jde),lfn(ids+1,jde-1))
269    if(ite.eq.ide.and.jts.eq.jds)lfn(ide+1,jds-1)=EX(lfn(ide,jds),lfn(ide-1,jds+1))
270    if(ite.eq.ide.and.jte.eq.jde)lfn(ide+1,jde+1)=EX(lfn(ide,jde),lfn(ide-1,jde-1))
271endif
272return
273contains
274real function EX(a,b)
275!*** statement function
276real a,b
277EX=(1.-bias)*(2.*a-b)+bias*max(2.*a-b,a,b)   ! extrapolation, max quarded
278end function EX
279end subroutine continue_at_boundary
280
281!
282!*****************************
283!
284
285subroutine check_mesh_2dim(ids,ide,jds,jde,ims,ime,jms,jme)
286implicit none
287integer, intent(in)::ids,ide,jds,jde,ims,ime,jms,jme
288character(len=128)msg
289if(ids<ims.or.ide>ime.or.jds<jms.or.jde>jme)then
290!$OMP CRITICAL(SFIRE_UTIL_CRIT)
291    write(msg,*)'mesh dimensions:  ',ids,ide,jds,jde
292    call message(msg)
293    write(msg,*)'memory dimensions:',ims,ime,jms,jme
294!$OMP END CRITICAL(SFIRE_UTIL_CRIT)
295    call message(msg)
296    call crash('check_mesh_2dim: memory dimensions too small')
297endif
298end subroutine check_mesh_2dim
299
300!
301!****************
302!
303
304subroutine check_mesh_3dim(ids,ide,kds,kde,jds,jde,ims,ime,kms,kme,jms,jme)
305integer, intent(in)::ids,ide,jds,jde,ims,ime,jms,jme,kds,kde,kms,kme
306if(ids<ims.or.ide>ime.or.jds<jms.or.jde>jme.or.kds<kms.or.kde>kme) then
307    call crash('memory dimensions too small')
308endif
309end subroutine check_mesh_3dim
310
311!
312!****************
313!
314
315subroutine sum_2d_cells(       &
316       ims2,ime2,jms2,jme2,    &
317       its2,ite2,jts2,jte2,    &
318       v2,                     &  ! input       
319       ims1,ime1,jms1,jme1,    &
320       its1,ite1,jts1,jte1,    &
321       v1)                        ! output
322implicit none
323
324!*** purpose
325! sum cell values in mesh2 to cell values of coarser mesh1
326
327!*** arguments
328! the dimensions are in cells, not nodes!
329
330integer, intent(in)::its1,ite1,jts1,jte1,ims1,ime1,jms1,jme1
331real, intent(out)::v1(ims1:ime1,jms1:jme1)
332integer, intent(in)::its2,ite2,jts2,jte2,ims2,ime2,jms2,jme2
333real, intent(in)::v2(ims2:ime2,jms2:jme2)
334
335!*** local
336integer:: i1,i2,j1,j2,ir,jr,isz1,isz2,jsz1,jsz2,ioff,joff,ibase,jbase
337real t
338character(len=128)msg
339
340!*** executable
341
342!check mesh dimensions and domain dimensions
343call check_mesh_2dim(its1,ite1,jts1,jte1,ims1,ime1,jms1,jme1)
344call check_mesh_2dim(its2,ite2,jts2,jte2,ims2,ime2,jms2,jme2)
345
346! compute mesh sizes
347isz1 = ite1-its1+1
348jsz1 = jte1-jts1+1
349isz2 = ite2-its2+1
350jsz2 = jte2-jts2+1
351
352! check mesh sizes
353if(isz1.le.0.or.jsz1.le.0.or.isz2.le.0.or.jsz2.le.0)then
354    call message('all mesh sizes must be positive')
355    goto 9
356endif
357
358! compute mesh ratios
359ir=isz2/isz1
360jr=jsz2/jsz1
361
362if(isz2.ne.isz1*ir .or. jsz2.ne.jsz1*jr)then
363    call message('input mesh size must be multiple of output mesh size')
364    goto 9
365endif
366
367
368! v1 = sum(v2)
369do j1=jts1,jte1
370    jbase=jts2+jr*(j1-jts1)
371    do i1=its1,ite1
372       ibase=its2+ir*(i1-its1)
373       t=0.
374       do joff=0,jr-1
375           j2=joff+jbase
376           do ioff=0,ir-1
377               i2=ioff+ibase
378               t=t+v2(i2,j2)
379           enddo
380       enddo
381       v1(i1,j1)=t
382    enddo
383enddo
384
385return
386
3879 continue
388!$OMP CRITICAL(SFIRE_UTIL_CRIT)
389write(msg,91)its2,ite2,jts2,jte2,ims2,ime2,jms2,jme2
390call message(msg)
391write(msg,91)its1,ite1,jts1,jte1,ims1,ime1,jms1,jme1
392call message(msg)
393write(msg,92)'input  mesh size:',isz2,jsz2
394call message(msg)
39591 format('dimensions: ',8i8)
396write(msg,92)'output mesh size:',isz1,jsz1
397call message(msg)
39892 format(a,2i8)
399!$OMP END CRITICAL(SFIRE_UTIL_CRIT)
400call crash('module_fr_spread_util:sum_mesh_cells: bad mesh sizes')
401
402end subroutine sum_2d_cells
403
404
405
406! module_fr_sfire_util%%interpolate_2d
407subroutine interpolate_2d(  &
408    ims2,ime2,jms2,jme2, & ! array coarse grid
409    its2,ite2,jts2,jte2, & ! dimensions coarse grid
410    ims1,ime1,jms1,jme1, & ! array coarse grid
411    its1,ite1,jts1,jte1, & ! dimensions fine grid
412    ir,jr,               & ! refinement ratio
413    rip2,rjp2,rip1,rjp1, & ! (rip2,rjp2) on grid 2 lines up with (rip1,rjp1) on grid 1
414    v2, &                  ! in coarse grid 
415    v1  )                  ! out fine grid
416implicit none
417
418!*** purpose
419! interpolate nodal values in mesh2 to nodal values in mesh1
420! interpolation runs over the mesh2 region its2:ite2,jts2:jte2
421! only the part of mesh 1 in the region its1:ite1,jts1:jte1 is modified
422
423!*** arguments
424
425integer, intent(in)::its1,ite1,jts1,jte1,ims1,ime1,jms1,jme1
426integer, intent(in)::its2,ite2,jts2,jte2,ims2,ime2,jms2,jme2
427integer, intent(in)::ir,jr
428real,intent(in):: rjp1,rip1,rjp2,rip2
429real, intent(out)::v1(ims1:ime1,jms1:jme1)
430real, intent(in)::v2(ims2:ime2,jms2:jme2)
431
432!*** local
433integer:: i1,i2,j1,j2,is,ie,js,je
434real:: tx,ty,rx,ry
435real:: rio,rjo
436intrinsic::ceiling,floor
437
438!*** executable
439
440!check mesh dimensions and domain dimensions
441call check_mesh_2dim(its1,ite1,jts1,jte1,ims1,ime1,jms1,jme1)
442call check_mesh_2dim(its2,ite2,jts2,jte2,ims2,ime2,jms2,jme2)
443
444! compute mesh ratios
445rx=1./ir
446ry=1./jr
447
448do j2=jts2,jte2-1               ! loop over mesh 2 cells
449    rjo=rjp1+jr*(j2-rjp2)       ! mesh 1 coordinate of the mesh 2 patch start
450    js=max(jts1,ceiling(rjo))   ! lower bound of mesh 1 patch for this mesh 2 cell
451    je=min(jte1,floor(rjo)+jr)  ! upper bound of mesh 1 patch for this mesh 2 cell
452    do i2=its2,ite2-1
453        rio=rip1+ir*(i2-rip2)
454        is=max(its1,ceiling(rio))
455        ie=min(ite1,floor(rio)+ir)
456        do j1=js,je
457            ty = (j1-rjo)*ry
458            do i1=is,ie
459                ! in case mesh 1 node lies on the boundary of several mesh 2 cells
460                ! the result will be written multiple times with the same value
461                ! up to a rounding error
462                tx = (i1-rio)*rx
463                !print *,'coarse ',i2,j2,'to',i2+1,j2+1,' fine ',is,js,' to ',ie,je
464                v1(i1,j1)=                     &
465                      (1-tx)*(1-ty)*v2(i2,j2)  &
466                 +    (1-tx)*ty  *v2(i2,j2+1)  &
467                 +      tx*(1-ty)*v2(i2+1,j2)  &
468                 +        tx*ty  *v2(i2+1,j2+1) 
469                !print *,'coarse ',i2,j2,' fine ',i1,j1, ' offset ',io,jo,' weights ',tx,ty, &
470                ! 'in ',v2(i2,j2),v2(i2,j2+1),v2(i2+1,j2),v2(i2+1,j2+1),' out ',v1(i1,j1)
471                !write(*,'(a,2i5,a,2f8.2,a,4f8.2,a,2i5,a,f8.2)') &
472                !'coarse ',i2,j2,' coord',rio,rjo,' val',v2(i2,j2),v2(i2,j2+1),v2(i2+1,j2),v2(i2+1,j2+1),&
473                !' fine ',i1,j1,' out ',v1(i1,j1)
474           enddo
475       enddo
476    enddo
477enddo
478
479end subroutine interpolate_2d
480
481!
482!****************
483!
484
485subroutine interpolate_2d_cells2cells(              &
486      ids2,ide2,jds2,jde2,ims2,ime2,jms2,jme2,v2, &  ! in 
487      ids1,ide1,jds1,jde1,ims1,ime1,jms1,jme1,v1  )  ! out
488implicit none
489
490!*** purpose
491! interpolate nodal values in mesh2 to nodal values in mesh1
492! input mesh 2 is coarse output mesh 1 is fine
493
494!*** arguments
495
496integer, intent(in)::ids1,ide1,jds1,jde1,ims1,ime1,jms1,jme1
497real, intent(out)::v1(ims1:ime1,jms1:jme1)
498integer, intent(in)::ids2,ide2,jds2,jde2,ims2,ime2,jms2,jme2
499real, intent(in)::v2(ims2:ime2,jms2:jme2)
500
501! Example with mesh ratio=4. | = cell boundary,  x = cell center
502!
503!  mesh2   |-------x-------|-------x-------|
504!  mesh1   |-x-|-x-|-x-|-x-|-x-|-x-|-x-|-x-|
505!
506
507!*** local
508integer:: ir,jr,isz1,isz2,jsz1,jsz2,ip,jp,ih,jh
509character(len=128)msg
510
511!*** executable
512
513!check mesh dimensions and domain dimensions
514call check_mesh_2dim(ids1,ide1,jds1,jde1,ims1,ime1,jms1,jme1)
515call check_mesh_2dim(ids2,ide2,jds2,jde2,ims2,ime2,jms2,jme2)
516
517! compute mesh sizes
518isz1 = ide1-ids1+1
519jsz1 = jde1-jds1+1
520isz2 = ide2-ids2+1
521jsz2 = jde2-jds2+1
522
523! check mesh sizes
524if(isz1.le.0.or.jsz1.le.0.or.isz2.le.0.or.jsz2.le.0)goto 9
525if(mod(isz1,isz2).ne.0.or.mod(jsz1,jsz2).ne.0)goto 9
526
527! compute mesh ratios
528ir=isz1/isz2
529jr=jsz1/jsz2
530!
531!  mesh2   |-------x-------|-------x-------|
532!  mesh1   |-x-|-x-|-x-|-x-|-x-|-x-|-x-|-x-|
533
534!  mesh2   |-----x-----|-----x-----|  rx=3
535!  mesh1   |-x-|-x-|-x-|-x-|-x-|-x-|
536!  i2            1   1   1   2
537!  i1        1   2   3   4   5
538!  ioff          0   1   2   0
539!  tx            0  1/3 2/3
540
541!  mesh2   |---x---|---x---| rx=2
542!  mesh1   |-x-|-x-|-x-|-x-|
543!  i2            1   1   2 
544!  i1            2   3   4
545!  ioff          0   1   2   
546!  tx           1/4 3/4
547
548
549! offset of the last node in the 1st half of the cell
550ih=ir/2
551jh=jr/2
552! 0 if coarse cell center coincides with fine, 1 if not
553ip=mod(ir+1,2)
554jp=mod(jr+1,2)
555
556call interpolate_2d_w(ip,jp,ih,jh,ir,jr,              &
557      ids2,ide2,jds2,jde2,ims2,ime2,jms2,jme2,v2, &  ! in 
558      ids1,ide1,jds1,jde1,ims1,ime1,jms1,jme1,v1  )  ! out
559
560return
561
5629 continue
563!$OMP CRITICAL(SFIRE_UTIL_CRIT)
564write(msg,91)ids2,ide2,jds2,jde2,ims2,ime2,jms2,jme2
565call message(msg)
566write(msg,91)ids1,ide1,jds1,jde1,ims1,ime1,jms1,jme1
567call message(msg)
568write(msg,92)'input  mesh size:',isz2,jsz2
569call message(msg)
57091 format('dimensions: ',8i8)
571write(msg,92)'output mesh size:',isz1,jsz1
572call message(msg)
57392 format(a,2i8)
574call crash("module_fr_sfire_util:interpolate_2dmesh_cells: bad mesh sizes")
575!$OMP END CRITICAL(SFIRE_UTIL_CRIT)
576end subroutine interpolate_2d_cells2cells
577
578!
579!****************
580!
581
582subroutine interpolate_2d_cells2nodes(              &
583      ids2,ide2,jds2,jde2,ims2,ime2,jms2,jme2,v2, &  ! in 
584      ids1,ide1,jds1,jde1,ims1,ime1,jms1,jme1,v1  )  ! out
585implicit none
586
587!*** purpose
588! interpolate nodal values in mesh2 to nodal values in mesh1
589! input mesh 2 is coarse output mesh 1 is fine
590
591!*** arguments
592
593integer, intent(in)::ids1,ide1,jds1,jde1,ims1,ime1,jms1,jme1
594real, intent(out)::v1(ims1:ime1,jms1:jme1)
595integer, intent(in)::ids2,ide2,jds2,jde2,ims2,ime2,jms2,jme2
596real, intent(in)::v2(ims2:ime2,jms2:jme2)
597
598! Example with mesh ratio=4. | = cell boundary,  x = cell center
599!
600!  mesh2   |-------x-------|-------x-------|
601!  mesh1   x-|-x-|-x-|-x-|-x-|-x-|-x-|-x-|-x
602!
603
604!*** local
605integer:: ir,jr,isz1,isz2,jsz1,jsz2,ip,jp,ih,jh
606character(len=128)msg
607
608!*** executable
609
610!check mesh dimensions and domain dimensions
611call check_mesh_2dim(ids1,ide1+1,jds1,jde1+1,ims1,ime1,jms1,jme1)
612call check_mesh_2dim(ids2,ide2,jds2,jde2,ims2,ime2,jms2,jme2)
613
614! compute mesh sizes
615isz1 = ide1-ids1+1
616jsz1 = jde1-jds1+1
617isz2 = ide2-ids2+1
618jsz2 = jde2-jds2+1
619
620! check mesh sizes
621if(isz1.le.0.or.jsz1.le.0.or.isz2.le.0.or.jsz2.le.0)goto 9
622if(mod(isz1,isz2).ne.0.or.mod(jsz1,jsz2).ne.0)goto 9
623
624! compute mesh ratios
625ir=isz1/isz2
626jr=jsz1/jsz2
627!
628!  mesh2   |-------x-------|-------x-------|
629!  mesh1   x-|-x-|-x-|-x-|-x-|-x-|-x-|-x-|-x
630
631!  mesh2   |-----x-----|-----x-----|  rx=3
632!  mesh1   x-|-x-|-x-|-x-|-x-|-x-|-x
633
634!  mesh2   |---x---|---x---| rx=2
635!  mesh1   x-|-x-|-x-|-x-|-x
636
637! offset of the last node in the 1st half of the cell
638ih=(ir+1)/2
639jh=(jr+1)/2
640! 0 if coarse cell center coincides with fine, 1 if not
641ip=mod(ir,2)
642jp=mod(jr,2)
643
644
645call interpolate_2d_w(ip,jp,ih,jh,ir,jr,              &
646      ids2,ide2,jds2,jde2,ims2,ime2,jms2,jme2,v2, &  ! in 
647      ids1,ide1+1,jds1,jde1+1,ims1,ime1,jms1,jme1,v1  )  ! out
648
649
650return
6519 continue
652!$OMP CRITICAL(SFIRE_UTIL_CRIT)
653write(msg,91)ids2,ide2,jds2,jde2,ims2,ime2,jms2,jme2
654call message(msg)
655write(msg,91)ids1,ide1,jds1,jde1,ims1,ime1,jms1,jme1
656call message(msg)
657write(msg,92)'input  mesh size:',isz2,jsz2
658call message(msg)
65991 format('dimensions: ',8i8)
660write(msg,92)'output mesh size:',isz1,jsz1
661call message(msg)
66292 format(a,2i8)
663call crash("module_fr_sfire_util:interpolate_2d_cells2nodes: bad mesh sizes")
664!$OMP END CRITICAL(SFIRE_UTIL_CRIT)
665end subroutine interpolate_2d_cells2nodes
666!
667!****************
668!
669
670subroutine interpolate_2d_w(ip,jp,ih,jh,ir,jr,             &
671      ids2,ide2,jds2,jde2,ims2,ime2,jms2,jme2,v2, &  ! in 
672      ids1,ide1,jds1,jde1,ims1,ime1,jms1,jme1,v1  )  ! out
673implicit none
674!*** EXCEPTION: THIS SUBROUTINE IS NEITHER CELL NOR NODE BASED.
675
676integer, intent(in)::ip,jp,ih,jh,ir,jr
677integer, intent(in)::ids1,ide1,jds1,jde1,ims1,ime1,jms1,jme1
678real, intent(out)::v1(ims1:ime1,jms1:jme1)
679integer, intent(in)::ids2,ide2,jds2,jde2,ims2,ime2,jms2,jme2
680real, intent(in)::v2(ims2:ime2,jms2:jme2)
681
682real:: tx,ty,rx,ry,half,xoff,yoff
683integer:: i1,i2,j1,j2,ioff,joff
684parameter(half=0.5)
685
686rx=ir
687ry=jr
688
689xoff = ip*half
690yoff = jp*half
691
692! the inside, ids1+ih:ide1-ih,jds1+jh:jde1-jh
693do j2=jds2,jde2-1     ! interpolate from nodes j2 and j2+1
694    do i2=ids2,ide2-1
695        do ioff=0,ir-ip
696            do joff=0,jr-jp
697                ! compute fine mesh index
698                i1=ioff+(ih+ids1)+ir*(i2-ids2)
699                j1=joff+(jh+jds1)+jr*(j2-jds2)
700                ! weights
701                tx = (ioff+xoff)/rx
702                ty = (joff+yoff)/ry
703                ! interpolation
704                v1(i1,j1)=                     &
705                      (1-tx)*(1-ty)*v2(i2,j2)  &
706                 +    (1-tx)*ty  *v2(i2,j2+1)  &
707                 +      tx*(1-ty)*v2(i2+1,j2)  &
708                 +        tx*ty  *v2(i2+1,j2+1) 
709                !write(*,'(3(a,2i5),a,2f7.4)')'coarse ',i2,j2,' fine ',i1,j1, &
710                ! ' offset ',ioff,joff,' weights ',tx,ty
711                !write(*,'(a,4f7.4,a,f7.4)')'in ',v2(i2,j2),v2(i2,j2+1),v2(i2+1,j2), &
712                !  v2(i2+1,j2+1),' out ',v1(i1,j1)
713           enddo
714       enddo
715    enddo
716enddo
717
718! extend to the boundary strips from the nearest known
719do ioff=0,ih-1  ! top and bottom strips
720    do j2=jds2,jde2-1
721        do joff=0,jr-jp
722           j1=joff+(jh+jds1)+jr*(j2-jds2)
723           ! weights
724           ty = (joff+yoff)/ry
725           ! interpolation
726           v1(ids1+ioff,j1)=(1-ty)*v2(ids2,j2)+ty*v2(ids2,j2+1)
727           v1(ide1-ioff,j1)=(1-ty)*v2(ide2,j2)+ty*v2(ide2,j2+1)
728       enddo
729    enddo
730enddo
731do joff=0,jh-1  ! left and right strips
732    do i2=ids2,ide2-1
733        do ioff=0,ir-ip
734           i1=ioff+(ih+ids1)+ir*(i2-ids2)
735           ! weights
736           tx = (ioff+xoff)/rx
737           ! interpolation
738           v1(i1,jds1+joff)=(1-tx)*v2(i2,jds2)+tx*v2(i2+1,jds2)
739           v1(i1,jde1-joff)=(1-tx)*v2(i2,jde2)+tx*v2(i2+1,jde2)
740       enddo
741    enddo
742enddo
743! extend to the 4 corner squares from the nearest known
744do ioff=0,ih-1 
745    do joff=0,jh-1
746        v1(ids1+ioff,jds1+joff)=v2(ids2,jds2)
747        v1(ide1-ioff,jds1+joff)=v2(ide2,jds2)
748        v1(ids1+ioff,jde1-joff)=v2(ids2,jde2)
749        v1(ide1-ioff,jde1-joff)=v2(ide2,jde2)
750    enddo
751enddo         
752end subroutine interpolate_2d_w 
753
754!
755!****************
756!
757               
758real function interp(ids,ide,jds,jde,ims,ime,jms,jme,x,y,v)
759implicit none
760!*** purpose
761! general interpolation in a rectangular
762
763!*** arguments
764
765integer, intent(in)::ids,ide,jds,jde,ims,ime,jms,jme
766real, intent(in)::x,y,v(ims:ime,jms:jme)
767! the mesh is cell based so the used dimension of v is ids:ide+1,jds:jde+1
768
769!*** calls
770intrinsic floor,min,max
771
772!*** local
773integer i,j
774real tx,ty
775
776! executable
777
778! indices of the lower left corner of the cell in the mesh that contains (x,y)
779i = floor(x)
780i=max(min(i,ide),ids)
781j = floor(y)
782j=max(min(j,jde),jds)
783
784! the leftover
785tx = x - real(i)
786ty = y - real(j)
787
788! interpolate the values
789interp = &
790                    (1-tx)*(1-ty)*v(i,j)    &
791                 +    tx*(1-ty)  *v(i+1,j)  &
792                 +    (1-tx)*ty  *v(i,j+1)  &
793                 +        tx*ty  *v(i+1,j+1) 
794
795!print *,'x,y=',x,y,'i1,i2=',i1,j1,'tx,ty=',tx,ty,' interp=',interp
796end function interp
797
798subroutine meshdiffc_2d(ids, ide, jds,jde , &    ! mesh area used (in cells, end +1)
799                   ims1,ime1,jms1,jme1, &       ! memory dimensiuons
800                   dx,dy,               &       ! mesh spacing
801                   lfn,                 &       ! input
802                   diffCx,diffCy) ! output
803implicit none
804
805!*** purpose
806! central differences on a 2d mesh
807
808!*** arguments
809
810integer, intent(in)::ids,ide,jds,jde,ims1,ime1,jms1,jme1
811real, intent(in):: dx,dy
812real, intent(in), dimension(ims1:ime1,jms1:jme1):: lfn
813real, intent(out), dimension(ims1:ime1,jms1:jme1):: diffCx,diffCy
814
815!*** local
816integer:: i,j
817real, dimension(ims1:ime1,jms1:jme1):: diffLx,diffRx,diffLy,diffRy
818
819! get one-sided differences; dumb but had that already...
820call meshdiff_2d(ids, ide, jds,jde , &    ! mesh area used (in cells, end +1)
821                   ims1,ime1,jms1,jme1, &       ! dimensions of lfn
822                   dx,dy,               &       ! mesh spacing
823                   lfn,                 &       ! input
824                   diffLx,diffRx,diffLy,diffRy) ! output
825
826! make into central
827do j=jds,jde+1
828    do i=ids,ide+1
829        diffCx(i,j)=0.5*(diffLx(i,j) + diffRx(i,j))
830        diffCy(i,j)=0.5*(diffLy(i,j) + diffRy(i,j))
831    enddo
832enddo
833end subroutine meshdiffc_2d
834
835subroutine meshdiff_2d(ids, ide, jds,jde , &    ! mesh area used (in cells, end +1)
836                   ims1,ime1,jms1,jme1, &       ! dimensions of lfn
837                   dx,dy,               &       ! mesh spacing
838                   lfn,                 &       ! input
839                   diffLx,diffRx,diffLy,diffRy) ! output
840implicit none
841
842!*** purpose
843! one-sided differences on a 2d mesh
844
845!*** arguments
846
847integer, intent(in)::ids,ide,jds,jde,ims1,ime1,jms1,jme1
848real, intent(in):: dx,dy
849real, intent(in), dimension(ims1:ime1,jms1:jme1):: lfn
850real, intent(out), dimension(ims1:ime1,jms1:jme1):: diffLx,diffRx,diffLy,diffRy
851
852!*** local
853integer:: i,j
854real:: tmpx,tmpy
855
856!*** executable
857
858    call check_mesh_2dim(ids,ide+1,jds,jde+1,ims1,ime1,jms1,jme1)
859 
860    ! the bulk of the work
861    do j=jds,jde
862        do i=ids,ide
863            tmpx = (lfn(i+1,j)-lfn(i,j))/dx
864            diffLx(i+1,j) = tmpx
865            diffRx(i,j)   = tmpx
866            tmpy = (lfn(i,j+1)-lfn(i,j))/dy
867            diffLy(i,j+1) = tmpy
868            diffRy(i,j)   = tmpy
869        enddo
870        ! missing values - put there the other one
871        diffLx(ids,j)  = diffLx(ids+1,j)
872        diffRx(ide+1,j)= diffRx(ide,j)
873    enddo
874    ! cleanup
875    ! j=jde+1 from above loop
876    do i=ids,ide
877        tmpx = (lfn(i+1,j)-lfn(i,j))/dx
878        diffLx(i+1,j) = tmpx
879        diffRx(i,j)   = tmpx
880    enddo
881    ! i=ide+1 from above loop
882    do j=jds,jde
883        tmpy = (lfn(i,j+1)-lfn(i,j))/dy
884        diffLy(i,j+1) = tmpy
885        diffRy(i,j)   = tmpy
886    enddo
887    ! missing values - put there the other one
888    ! j=jde+1 from above loop, j=jds:jde done before in main bulk loop
889    diffLx(ids,j)   = diffLx(ids+1,j)
890    diffRx(ide+1,j) = diffRx(ide,j)
891    do i=ids,ide+1
892        diffLy(i,jds)   = diffLy(i,jds+1)
893        diffRy(i,jde+1) = diffRy(i,jde)
894    enddo   
895
896end subroutine meshdiff_2d
897
898
899
900
901real pure function sum_2darray( its,ite,jts,jte,               &
902                                ims,ime,jms,jme,               &
903                                a)
904integer, intent(in)::its,ite,jts,jte,ims,ime,jms,jme
905real, intent(in)::a(ims:ime,jms:jme)
906!*** local
907integer:: i,j
908real:: t
909t=0.
910do j=jts,jte
911    do i=its,ite
912        t=t+a(i,j)
913    enddo
914enddo
915sum_2darray = t
916end function sum_2darray
917
918real pure function max_2darray( its,ite,jts,jte,               &
919                                ims,ime,jms,jme,               &
920                                a)
921integer, intent(in)::its,ite,jts,jte,ims,ime,jms,jme
922real, intent(in)::a(ims:ime,jms:jme)
923!*** local
924integer:: i,j
925real:: t
926t=0.
927do j=jts,jte
928    do i=its,ite
929        t=max(t,a(i,j))
930    enddo
931enddo
932max_2darray = t
933end function max_2darray
934
935subroutine print_2d_stats_vec(ips,ipe,jps,jpe, &
936                         ims,ime,jms,jme, &
937                         ax,ay,name)
938implicit none
939integer, intent(in)::ips,ipe,jps,jpe,ims,ime,jms,jme
940real, intent(in), dimension(ims:ime,jms:jme)::ax,ay
941character(len=*),intent(in)::name
942integer:: i,j
943real:: t
944real:: avg_a,max_a,min_a
945character(len=25)::id
946id=name
947call print_2d_stats(ips,ipe,jps,jpe, &
948                         ims,ime,jms,jme, &
949                         ax,id//'/x ')
950call print_2d_stats(ips,ipe,jps,jpe, &
951                         ims,ime,jms,jme, &
952                         ay,id//'/y ')
953avg_a=0
954max_a=-huge(max_a)
955min_a= huge(min_a)
956do j=jps,jpe
957    do i=ips,ipe
958        t=sqrt(ax(i,j)**2+ay(i,j)**2)
959        max_a=max(max_a,t)
960        min_a=min(min_a,t)
961        avg_a=avg_a+t
962    enddo
963enddo
964avg_a = avg_a/((ipe-ips+1)*(jpe-jps+1))
965call print_stat_line(id//'/sz',ips,ipe,jps,jpe,min_a,max_a,avg_a)
966end subroutine print_2d_stats_vec
967
968
969subroutine print_stat_line(name,ips,ipe,jps,jpe,min_a,max_a,avg_a)
970!*** encapsulate line with statistics
971implicit none
972!*** arguments
973integer, intent(in)::ips,ipe,jps,jpe
974character(len=*),intent(in)::name
975real,intent(in)::min_a,max_a,avg_a
976!*** local
977character(len=128)::msg
978character(len=24)::id
979if(fire_print_msg.eq.0)return
980id=name
981!$OMP CRITICAL(SFIRE_UTIL_CRIT)
982write(msg,'(a,4i4,3g11.3)')id,ips,ipe,jps,jpe,min_a,max_a,avg_a
983!$OMP END CRITICAL(SFIRE_UTIL_CRIT)
984call message(msg)
985if(.not.avg_a.eq.avg_a)call crash('NaN detected')
986end subroutine print_stat_line
987
988
989subroutine print_3d_stats(ips,ipe,kps,kpe,jps,jpe, &
990                         ims,ime,kms,kme,jms,jme, &
991                         a,name)
992implicit none
993integer, intent(in)::ips,ipe,jps,jpe,ims,ime,jms,jme,kms,kme,kps,kpe
994real, intent(in)::a(ims:ime,kms:kme,jms:jme)
995character(len=*),intent(in)::name
996integer:: i,j,k
997real:: avg_a,max_a,min_a,t,aa,bb
998character(len=128)::msg
999! if(fire_print_msg.eq.0)return
1000! check for nans in any case
1001bb=0.
1002do j=jps,jpe
1003  do k=kps,kpe
1004    do i=ips,ipe
1005       bb=bb+a(i,k,j)
1006    enddo
1007  enddo
1008enddo
1009if(bb.eq.bb.and.fire_print_msg.eq.0)return
1010avg_a=0
1011max_a=-huge(max_a)
1012min_a= huge(min_a)
1013t=huge(t)
1014do j=jps,jpe
1015  do k=kps,kpe
1016    do i=ips,ipe
1017        aa=a(i,k,j)
1018        if(aa.ne.aa.or..not.aa.le.t.or..not.aa.ge.-t)goto 9
1019        max_a=max(max_a,aa)
1020        min_a=min(min_a,aa)
1021        avg_a=avg_a+aa
1022    enddo
1023  enddo
1024enddo
1025if(bb.ne.bb)goto 10 ! should never happen
1026if(fire_print_msg.eq.0)return
1027avg_a = avg_a/((ipe-ips+1)*(jpe-jps+1)*(kpe-kps+1))
1028call print_stat_line(name,ips,ipe,jps,jpe,min_a,max_a,avg_a)
1029return
10309 continue
1031!$OMP CRITICAL(SFIRE_UTIL_CRIT)
1032write(msg,1)name,i,k,j,aa
1033call message(msg)
10341 format(a30,'(',i6,',',i6,',',i6,') = ',g13.5)
1035!$OMP END CRITICAL(SFIRE_UTIL_CRIT)
1036call print_stat_line(name,ips,ipe,jps,jpe,aa,aa,aa)
1037if(aa.ne.aa)goto 10
1038msg='Invalid floating point number detected in '//name
1039call crash(msg)
104010 msg='NaN detected in '//name
1041call crash(msg)
1042end subroutine print_3d_stats
1043
1044subroutine print_2d_stats(ips,ipe,jps,jpe, &
1045                         ims,ime,jms,jme, &
1046                         a,name)
1047implicit none
1048integer, intent(in)::ips,ipe,jps,jpe,ims,ime,jms,jme
1049real, intent(in)::a(ims:ime,jms:jme)
1050character(len=*),intent(in)::name
1051!!character(len=128)::msg
1052!if(fire_print_msg.eq.0)return
1053call print_3d_stats(ips,ipe,1,1,jps,jpe, &
1054                         ims,ime,1,1,jms,jme, &
1055                         a,name)
1056!!write(msg,'(2a,z16)')name,' address =',loc(a)
1057!!call message(msg)
1058end subroutine print_2d_stats
1059
1060real pure function avg_2darray( its,ite,jts,jte,               &
1061                                ims,ime,jms,jme,               &
1062                                a)
1063integer, intent(in)::its,ite,jts,jte,ims,ime,jms,jme
1064real, intent(in)::a(ims:ime,jms:jme)
1065!*** local
1066!*** executable
1067avg_2darray = sum_2darray( its,ite,jts,jte,               &
1068                           ims,ime,jms,jme,               &
1069                           a)/((ite-its+1)*(jte-jts+1))
1070end function avg_2darray
1071
1072real pure function avg_2darray_vec( its,ite,jts,jte,           &
1073                                ims,ime,jms,jme,               &
1074                                ax,ay)
1075integer, intent(in)::its,ite,jts,jte,ims,ime,jms,jme
1076real, intent(in), dimension(ims:ime,jms:jme):: ax,ay
1077!*** local
1078integer:: i,j
1079real:: t
1080t=0.
1081do j=jts,jte
1082    do i=its,ite
1083        t=t+sqrt(ax(i,j)**2+ay(i,j)**2)
1084    enddo
1085enddo
1086t = t/((ite-its+1)*(jte-jts+1))
1087avg_2darray_vec = t
1088end function avg_2darray_vec
1089
1090
1091subroutine print_array(its,ite,jts,jte,           &
1092                         ims,ime,jms,jme,               &
1093                         a,name,id)
1094! debug
1095!*** arguments
1096integer, intent(in)::its,ite,jts,jte,ims,ime,jms,jme,id
1097real, intent(in), dimension(ims:ime,jms:jme):: a
1098character(len=*),intent(in)::name
1099!****
1100integer i,j
1101character(len=128)::msg
1102!****
1103!$OMP CRITICAL(SFIRE_UTIL_CRIT)
1104write(msg,*)name,' start ',id,' dim ',its,ite,jts,jte
1105call message(msg)
1106do j=jts,jte
1107    do i=its,ite
1108         write(msg,*)i,j,a(i,j)
1109         call message(msg)
1110    enddo
1111enddo
1112write(msg,*)name,' end ',id
1113call message(msg)
1114!$OMP END CRITICAL(SFIRE_UTIL_CRIT)
1115end subroutine print_array
1116
1117subroutine write_array_m(its,ite,jts,jte,           &
1118                         ims,ime,jms,jme,               &
1119                         a,name,id)
1120! debug
1121!*** arguments
1122integer, intent(in)::its,ite,jts,jte,ims,ime,jms,jme,id
1123real, intent(in), dimension(ims:ime,jms:jme):: a
1124character(len=*),intent(in)::name
1125!****
1126call write_array_m3(its,ite,1,1,jts,jte,           &
1127                         ims,ime,1,1,jms,jme,               &
1128                         a,name,id)
1129end subroutine write_array_m
1130
1131
1132subroutine write_array_m3(its,ite,kts,kte,jts,jte,           &
1133                         ims,ime,kms,kme,jms,jme,               &
1134                         a,name,id)
1135
1136implicit none
1137! debug
1138!*** arguments
1139integer, intent(in)::its,ite,jts,jte,ims,ime,jms,jme,kts,kte,kms,kme,id
1140real, intent(in), dimension(ims:ime,kms:kme,jms:jme):: a
1141character(len=*),intent(in)::name
1142!****
1143integer i,j,k,iu,ilen,myproc,nprocs
1144logical op
1145character(len=128)::fname,msg
1146!****
1147if(fire_print_file.eq.0.or.id.le.0)return
1148call check_mesh_2dim(its,ite,jts,jte,ims,ime,jms,jme)
1149call wrf_get_nproc (nprocs)
1150call wrf_get_myproc(myproc)
1151
1152!$OMP CRITICAL(SFIRE_UTIL_CRIT)
1153if(nprocs.eq.1)then
1154    write(fname,3)name,'_',id,'.txt'
1155else
1156    write(fname,4)name,'_',id,'.',myproc,'.txt'
1157endif
1158
1159iu=0
1160do i=6,99
1161    inquire(unit=i,opened=op)
1162    if(.not.op.and.iu.le.0)iu=i
1163enddo
1164if(iu.gt.0)open(iu,file=trim(fname),form='formatted',status='unknown')
1165
1166if(iu.le.0)call crash('write_array_m: cannot find available fortran unit')
1167
1168write(iu,1)real(its)
1169write(iu,1)real(ite)
1170write(iu,1)real(jts)
1171write(iu,1)real(jte)
1172write(iu,1)real(kts)
1173write(iu,1)real(kte)
1174write(iu,1)(((a(i,k,j),i=its,ite),j=jts,jte),k=kts,kte)
1175close(iu)
1176write(msg,2)name,'(',its,':',ite,',',jts,':',jte,',', &
1177kts,':',kte,') -> ',trim(fname)
1178!$OMP END CRITICAL(SFIRE_UTIL_CRIT)
1179call message(msg)
1180return
1181
11821 format(e20.12)
11832 format(2a,3(i5,a,i5,a),2a)
11843 format(a,a,i8.8,a)
11854 format(a,a,i8.8,a,i4.4,a)
1186
1187
1188end subroutine write_array_m3
1189!
1190!***
1191!
1192subroutine read_array_2d_integer(filename,ia,its,ite,jts,jte,ims,ime,jms,jme)
1193!*** arguments
1194integer, intent(in)::its,ite,jts,jte,ims,ime,jms,jme
1195integer, intent(out), dimension(ims:ime,jms:jme):: ia
1196character(len=*),intent(in)::filename
1197!*** local
1198real, allocatable, dimension(:,:)::a
1199integer::i,j
1200real::r
1201character(len=256)msg
1202!*** executable
1203allocate(a(ims:ime,jms:jme))
1204call read_array_2d_real(filename,a,its,ite,jts,jte,ims,ime,jms,jme)
1205do j=jts,jte
1206    do i=its,ite
1207        ia(i,j)=a(i,j)
1208        r=ia(i,j)
1209        if(r.ne.a(i,j))then
1210           write(6,*)'Reading file ',trim(filename)
1211           write(6,*)'value at position ',i,j,' cannot be converted to integer'
1212           write(6,*)'read ',a(i,j),' converted as',k,' converted as ',r
1213           msg=trim(filename)//' is not integer file'
1214           call crash(msg)
1215        endif
1216    enddo
1217enddo
1218end subroutine read_array_2d_integer
1219
1220       
1221
1222subroutine read_array_2d_real(filename,a,its,ite,jts,jte,ims,ime,jms,jme)
1223#ifdef _OPENMP
1224use OMP_LIB
1225#endif
1226implicit none
1227!*** purpose: read a 2D array from a text file
1228! file format: line 1: array dimensions ni,nj
1229!              following lines: one row of a at a time
1230!              each row may be split between several lines
1231!*** arguments
1232integer, intent(in)::its,ite,jts,jte,ims,ime,jms,jme
1233real, intent(out), dimension(ims:ime,jms:jme):: a
1234character(len=*),intent(in)::filename
1235!*** local
1236integer i,j,ni,nj,mi,mj,nprocs,myproc,mythread,iu
1237logical op
1238character(len=128)::fname,msg
1239!*** executable
1240
1241call wrf_get_nproc (nprocs)
1242call wrf_get_myproc( myproc )
1243mythread=0
1244#ifdef _OPENMP
1245    mythread=omp_get_thread_num()
1246#endif
1247if(nprocs.ne.1.or.myproc.ne.0.or.mythread.ne.0) &
1248   call crash('read_array_2d: parallel execution not supported')
1249
1250! print line
1251mi=ite-its+1
1252mj=jte-jts+1
1253write(msg,2)'reading array size ',mi,mj,' from file ',trim(filename)
12542 format(a,2i6,2a)
1255call message(msg)
1256
1257! check array index overflow
1258call check_mesh_2dim(its,ite,jts,jte,ims,ime,jms,jme)
1259
1260! find available unit
1261iu=0
1262do i=11,99
1263    inquire(unit=i,opened=op)
1264    if(.not.op.and.iu.le.0)iu=i
1265enddo
1266if(iu.le.0)call crash('read_array_2d: cannot find available fortran unit')
1267
1268if(iu.gt.0)open(iu,file=filename,form='formatted',status='old',err=9)
1269rewind(iu,err=9)
1270
1271read(iu,*,err=10)ni,nj
1272if(ni.ne.mi.or.nj.ne.mj)then
1273    write(msg,'(a,2i6,a,2i6)')'Array dimensions',ni,nj,' in the input file should be ',mi,mj
1274    call message(msg)
1275    goto 10
1276endif
1277do i=its,ite
1278   read(iu,*,err=10)(a(i,j),j=jts,jte)
1279enddo
1280close(iu,err=11)
1281return
1282
12839  msg='Error opening file '//trim(filename)
1284call crash(msg)
128510 msg='Error reading file '//trim(filename)
1286call crash(msg)
128711 msg='Error closing file '//trim(filename)
1288call crash(msg)
1289end subroutine read_array_2d_real
1290!
1291!***
1292!
1293
1294! general conditional expression
1295pure integer function ifval(l,i,j)
1296implicit none
1297logical, intent(in)::l
1298integer, intent(in)::i,j
1299if(l)then
1300        ifval=i
1301else
1302        ifval=j
1303endif
1304end function ifval
1305
1306! function to go beyond domain boundary if tile is next to it
1307pure integer function snode(t,d,i)
1308implicit none
1309integer, intent(in)::t,d,i
1310if(t.ne.d)then
1311    snode=t
1312else
1313    snode=t+i
1314endif
1315end function snode
1316
1317subroutine print_chsum( id, &
1318    ims,ime,kms,kme,jms,jme, &                ! memory dims
1319    ids,ide,kds,kde,jds,jde, &                ! domain dims
1320    ips,ipe,kps,kpe,jps,jpe, &                ! patch or tile dims
1321    istag,kstag,jstag,       &
1322    a,name)                     
1323
1324#ifdef DM_PARALLEL
1325    USE module_dm , only : wrf_dm_bxor_integer
1326#endif
1327
1328integer, intent(in):: id, &
1329    ims,ime,kms,kme,jms,jme, &                ! memory dims
1330    ids,ide,kds,kde,jds,jde, &                ! domain dims
1331    ips,ipe,kps,kpe,jps,jpe, &                ! patch dims
1332    istag,kstag,jstag
1333real, intent(in),dimension(ims:ime,kms:kme,jms:jme)::a
1334character(len=*)::name
1335
1336!$ external, logical:: omp_in_parallel
1337!$ external, integer:: omp_get_thread_num
1338
1339!*** local
1340integer::lsum
1341integer::i,j,k,n,ipe1,jpe1,kpe1,iel,thread,is,js,ks
1342integer, save::psum,gsum
1343real::rel
1344equivalence(rel,iel)
1345character(len=256)msg
1346
1347if(fire_print_msg.le.0)return
1348
1349ipe1=ifval(ipe.eq.ide.and.istag.ne.0,ipe+1,ipe)
1350kpe1=ifval(kpe.eq.kde.and.kstag.ne.0,kpe+1,kpe)
1351jpe1=ifval(jpe.eq.jde.and.jstag.ne.0,jpe+1,jpe)
1352is=ifval(istag.ne.0,1,0)
1353ks=ifval(kstag.ne.0,1,0)
1354js=ifval(jstag.ne.0,1,0)
1355
1356lsum=0
1357do j=jps,jpe1
1358  do k=kps,kpe1
1359    do i=ips,ipe1
1360      rel=a(i,k,j)
1361      lsum=ieor(lsum,iel)
1362    enddo
1363  enddo
1364enddo
1365
1366! get process sum over all threads
1367thread=0
1368!$ thread=omp_get_thread_num()
1369if(thread.eq.0)psum=0
1370!$OMP BARRIER
1371!$OMP CRITICAL(CHSUM)
1372psum=ieor(psum,lsum)
1373!$OMP END CRITICAL(CHSUM)
1374!$OMP BARRIER
1375
1376! get global sum over all processes
1377if(thread.eq.0)then
1378#ifdef DM_PARALLEL
1379    gsum = wrf_dm_bxor_integer ( psum )
1380#else
1381    gsum = psum
1382#endif
1383    write(msg,1)id,name,ids,ide+is,kds,kde+ks,jds,jde+js,gsum
13841   format(i6,1x,a10,' dims',6i5,' chsum ',z8.8)
1385    call message(msg)
1386endif
1387
1388end subroutine print_chsum
1389
1390
1391real function fun_real(fun,  &
1392    ims,ime,kms,kme,jms,jme, &                ! memory dims
1393    ids,ide,kds,kde,jds,jde, &                ! domain dims
1394    ips,ipe,kps,kpe,jps,jpe, &                ! patch or tile dims
1395    istag,kstag,jstag,       &                ! staggering
1396    a,b)                     
1397
1398#ifdef DM_PARALLEL
1399    USE module_dm , only : wrf_dm_sum_real , wrf_dm_max_real
1400#endif
1401
1402integer, intent(in)::  fun, &
1403    ims,ime,kms,kme,jms,jme, &                ! memory dims
1404    ids,ide,kds,kde,jds,jde, &                ! domain dims
1405    ips,ipe,kps,kpe,jps,jpe, &                ! patch dims
1406    istag,kstag,jstag                         ! staggering
1407real, intent(in),dimension(ims:ime,kms:kme,jms:jme)::a,b
1408
1409!*** local
1410real::lsum,void
1411integer::i,j,k,n,ipe1,jpe1,kpe1,iel,thread,is,js,ks
1412real, save::psum,gsum
1413real::rel
1414logical:: dosum,domax
1415character(len=256)msg
1416
1417ipe1=ifval(ipe.eq.ide.and.istag.ne.0,ipe+1,ipe)
1418kpe1=ifval(kpe.eq.kde.and.kstag.ne.0,kpe+1,kpe)
1419jpe1=ifval(jpe.eq.jde.and.jstag.ne.0,jpe+1,jpe)
1420is=ifval(istag.ne.0,1,0)
1421ks=ifval(kstag.ne.0,1,0)
1422js=ifval(jstag.ne.0,1,0)
1423
1424if(fun.eq.REAL_SUM)then
1425  void=0.
1426  lsum=void
1427  do j=jps,jpe1
1428    do k=kps,kpe1
1429      do i=ips,ipe1
1430        lsum=lsum+a(i,k,j)
1431      enddo
1432    enddo
1433  enddo
1434elseif(fun.eq.RNRM_SUM)then
1435  void=0.
1436  lsum=void
1437  do j=jps,jpe1
1438    do k=kps,kpe1
1439      do i=ips,ipe1
1440        lsum=lsum+sqrt(a(i,k,j)*a(i,k,j)+b(i,k,j)*b(i,k,j))
1441      enddo
1442    enddo
1443  enddo
1444elseif(fun.eq.REAL_MAX)then
1445  void=-huge(lsum)
1446  lsum=void
1447  do j=jps,jpe1
1448    do k=kps,kpe1
1449      do i=ips,ipe1
1450        lsum=max(lsum,a(i,k,j))
1451      enddo
1452    enddo
1453  enddo
1454elseif(fun.eq.RNRM_MAX)then
1455  void=0.
1456  lsum=void
1457  do j=jps,jpe1
1458    do k=kps,kpe1
1459      do i=ips,ipe1
1460        lsum=max(lsum,sqrt(a(i,k,j)*a(i,k,j)+b(i,k,j)*b(i,k,j)))
1461      enddo
1462    enddo
1463  enddo
1464else
1465  call crash('fun_real: bad fun')
1466endif
1467
1468if(lsum.ne.lsum)call crash('fun_real: NaN detected')
1469
1470dosum=fun.eq.REAL_SUM.or.fun.eq.RNRM_SUM
1471domax=fun.eq.REAL_MAX.or.fun.eq.RNRM_MAX
1472
1473! get process sum over all threads
1474!$OMP SINGLE
1475! only one thread should write to shared variable
1476psum=void
1477!$OMP END SINGLE
1478!$OMP BARRIER
1479! now all threads know psum
1480
1481!$OMP CRITICAL(RDSUM)
1482! each thread adds its own lsum
1483if(dosum)psum=psum+lsum
1484if(domax)psum=max(psum,lsum)
1485!$OMP END CRITICAL(RDSUM)
1486
1487! wait till all theads are done
1488!$OMP BARRIER
1489
1490! get global sum over all processes
1491!$OMP SINGLE
1492! only one threads will do the mpi communication
1493#ifdef DM_PARALLEL
1494    if(dosum) gsum = wrf_dm_sum_real ( psum )
1495    if(domax) gsum = wrf_dm_max_real ( psum )
1496#else
1497    gsum = psum
1498#endif
1499if(gsum.ne.gsum)call crash('fun_real: NaN detected')
1500!$OMP END SINGLE
1501
1502!$OMP BARRIER
1503! now gsum is known to all threads
1504
1505fun_real=gsum
1506
1507end function fun_real
1508
1509end module module_fr_sfire_util
Note: See TracBrowser for help on using the repository browser.