Merge remote-tracking branch 'origin/release-v4.6.1'
[WRF.git] / hydro / Routing / module_UDMAP.F90
blob0b72a8eb0ce6bfbe84fc69e7e96fcd5322bacdb5
1 !  Program Name:
2 !  Author(s)/Contact(s):
3 !  Abstract:
4 !  History Log:
6 !  Usage:
7 !  Parameters: <Specify typical arguments passed>
8 !  Input Files:
9 !        <list file names and briefly describe the data they include>
10 !  Output Files:
11 !        <list file names and briefly describe the information they include>
13 !  Condition codes:
14 !        <list exit condition or error codes returned >
15 !        If appropriate, descriptive troubleshooting instructions or
16 !        likely causes for failures could be mentioned here with the
17 !        appropriate error code
19 !  User controllable options: <if applicable>
21 ! This subrouting includs the data structure and tools used for NHDPlus network mapping.
22 module module_UDMAP
24   use config_base, only: nlst
25 #ifdef MPP_LAND
26   use module_mpp_land, only: my_id, local_startx_rt, local_starty_rt,  &
27        local_endx_rt,local_endy_rt, left_id, right_id, down_id, up_id, mpp_collect_1d_int_mem, &
28        IO_id , numprocs
29   use module_mpp_land, only: mpp_land_bcast_int, mpp_land_bcast_real8_1d, mpp_land_bcast_int1, mpp_land_bcast_int8
31   use module_mpp_land, only: sum_int1d, global_rt_nx, global_rt_ny, write_IO_rt_int, MPP_LAND_COM_INTEGER
33   use MODULE_mpp_ReachLS, only : updatelinkv, ReachLS_write_io, com_write1dInt, &
34        com_decomp1dInt, pack_decomp_int, pack_decomp_real8, &
35        com_decomp1dint8, pack_decomp_int8, com_write1dInt8
36   use module_hydro_stop, only:HYDRO_stop
37   use hashtable
39   use iso_fortran_env, only: int64
40 #endif
42   implicit none
44 #ifndef MPP_LAND
45     integer, parameter :: numprocs=1
46 #endif
48 #include <netcdf.inc>
50 type userDefineMapping
51     integer, allocatable, dimension(:) :: grid_i, grid_j
52     real, allocatable, dimension(:) :: weight, nodeArea, cellArea
53     integer :: ngrids
54     integer(kind=int64) :: myid
55 !   for bucket model definition
56     real, allocatable, dimension(:) :: cellWeight
57     integer, allocatable, dimension(:) :: cell_i, cell_j
58     integer :: ncell
59 end type userDefineMapping
61 TYPE ( userDefineMapping ), allocatable, DIMENSION (:) :: LUDRSL
63 integer(kind=int64), allocatable, dimension(:) :: bufid
64 real*8 , allocatable, dimension(:) :: bufw
65 integer :: LNUMRSL  ! number of local links
66 integer :: ter_rt_flag
67 real*8, allocatable, dimension(:) :: basns_area
68 integer :: gnpid, lnsize
69 integer, allocatable, dimension(:) :: bufi,bufj
71 contains
72     subroutine UDMP_ini(nlinksl,ixrt,jxrt,rtmask, OVRTSWCRT, SUBRTSWCRT,cell_area)
73 !This is the driver for user defined mapping file funciton application.
74         integer :: ixrt, jxrt, OVRTSWCRT, SUBRTSWCRT, nlinksl
75         integer, intent(in), dimension(ixrt,jxrt):: rtmask
76         integer :: npid    !local variable.
77         real,dimension(:,:) :: cell_area
78         ter_rt_flag = 0
79         if(OVRTSWCRT .eq. 1 .or. SUBRTSWCRT .eq. 1) then
80             ter_rt_flag = 1
81         endif
82         call readUDMP(ixrt,jxrt,npid,nlinksl)
83         call UDMP2LOCAL(npid,ixrt,jxrt,rtmask, ter_rt_flag)
84         call getUDMP_area(cell_area)
85     end subroutine UDMP_ini
87     subroutine readUDMP(ixrt,jxrt,npid, nlinksl)
88         implicit none
89         integer :: i,j,Ndata, did, Npid, nlinksl, k, m, kk
90         integer,allocatable,dimension(:) ::  nprocs_map, lnsizes, istart
91         integer(kind=int64), allocatable, dimension(:) :: g1bufid, gbufid, linkid, bufid_tmp, bufidflag
92         integer :: ix_bufid, ii, ixrt,jxrt
93         integer, allocatable, dimension(:) :: gbufi,gbufj,bufsize
94         real*8 , allocatable, dimension(:) :: gbufw
96         did = 1
97         call get_dimension(trim(nlst(did)%UDMAP_FILE), ndata, npid)
99 #ifdef MPP_LAND
100         gnpid = npid
101         allocate (lnsizes(numprocs))
102         if(my_id .eq. io_id) then
103            allocate (istart(numprocs))
104            allocate (nprocs_map(ndata))
105            allocate(gbufi(ndata))
106            allocate(gbufj(ndata))
107            call get1d_int(trim(nlst(did)%UDMAP_FILE),"i_index",gbufi)
108            call get1d_int(trim(nlst(did)%UDMAP_FILE),"j_index",gbufj)
109         endif
110            call get_nprocs_map(ixrt,jxrt,gbufi,gbufj,nprocs_map,ndata)
112         if(my_id .eq. io_id) then
113            lnsizes = 0
114            do i =1 , ndata
115                if(nprocs_map(i) .gt. 0) then
116                   lnsizes(nprocs_map(i)) = lnsizes(nprocs_map(i)) + 1
117                endif
118            enddo
119         endif
120         call mpp_land_bcast_int(numprocs,lnsizes)
122      if(my_id .eq. io_id ) then
123         kk = 0
124         do i = 1, numprocs
125            kk = kk + lnsizes(i)
126         end do
127      end if
129       if(my_id .eq. IO_id) then
130           ii = 1
131           do i = 1, numprocs
132              istart(i) = ii
133              if(lnsizes(i) .gt. 0) then
134                 ii = lnsizes(i) + ii
135              else
136                 istart(i) = -999
137              endif
138           end do
139       endif
141       if(lnsizes(my_id+1) .gt. 0)  allocate(bufi(lnsizes(my_id+1) ))
142       call pack_decomp_int(gbufi, ndata, nprocs_map, lnsizes, istart,bufi)
143       if(my_id .eq. io_id) then
144            if(allocated(gbufi))  deallocate(gbufi)
145       endif
148       if(lnsizes(my_id+1) .gt. 0) allocate(bufj(lnsizes(my_id+1) ))
149       call pack_decomp_int(gbufj, ndata, nprocs_map, lnsizes, istart,bufj)
150       if(my_id .eq. io_id)  then
151          if(allocated(gbufj)) deallocate(gbufj)
152       endif
155 ! check bufid
156 !      check  polyid and linkid
157         allocate(linkid(nlinksl))
158         if(my_id .eq. io_id) then
159             call get1d_int64(trim(nlst(did)%route_link_f),"link",linkid)
160             allocate(gbufid(npid))
161             call get1d_int64(trim(nlst(did)%UDMAP_FILE),"polyid",gbufid)
162         endif
163 #ifdef MPP_LAND
164        if(nlinksl .gt. 0) then
165           call mpp_land_bcast_int8(nlinksl,linkid)
166        endif
167        call com_decomp1dInt8(gbufid,npid,bufid_tmp,ix_bufid)
168 #endif
169        if(ix_bufid .gt. 0) then
170           allocate(bufidflag(ix_bufid))
171           bufidflag = -999
172        endif
174        ! The following loops are replaced by a hashtable-based algorithm
175        !  do i = 1, ix_bufid
176        !           do j = 1, nlinksl
177        !                if(bufid_tmp(i) .eq. linkid(j)) then
178        !                   bufidflag(i) = bufid_tmp(i)
179        !                   goto 102
180        !                endif
181        !           end do
182        ! 102       continue
183        !        end do
185        block
186          type(hash_t) :: hash_table
187          integer(kind=int64) :: val,it
188          logical :: found
190          call hash_table%set_all_idx(bufid_tmp, ix_bufid)
191          do it = 1, nlinksl
192             call hash_table%get(linkid(it), val, found)
193             if(found .eqv. .true.) then
194                bufidflag(val) = bufid_tmp(val)
195             end if
196          end do
197          call hash_table%clear()
198        end block
200 #ifdef MPP_LAND
201       call com_write1dInt8(bufidflag,ix_bufid,gbufid,npid)
202 #endif
203       if(ix_bufid .gt. 0) then
204           if(allocated(bufidflag)) deallocate(bufidflag)
205           if(allocated(bufid_tmp)) deallocate(bufid_tmp)
206       endif
207       if(allocated(linkid)) deallocate(linkid)
208       if(my_id .eq. io_id) then
209           allocate(bufsize(npid))
210           allocate(g1bufid(ndata))
211           call get1d_int(trim(nlst(did)%UDMAP_FILE),"overlaps",bufsize)
212           g1bufid = -999
213           i = 1
214           do k = 1, npid
215                do j = 1, bufsize(k)
216                  g1bufid(i) = gbufid(k)
217                  i = i + 1
218                end do
219           enddo
220           if(allocated(bufsize))  deallocate(bufsize)
221       endif
224       if(my_id .eq. io_id) then
225            if(allocated(gbufid)) deallocate(gbufid)
226       endif
229       if(lnsizes(my_id+1) .gt. 0) allocate(bufid(lnsizes(my_id+1) ))
230       call pack_decomp_int8(g1bufid, ndata, nprocs_map, lnsizes, istart,bufid)
231       if(my_id .eq. io_id) then
232             if(allocated(g1bufid)) deallocate(g1bufid)
233       endif
236       if(my_id .eq. io_id) then
237           allocate(gbufw(ndata))
238           call get1d_real8(trim(nlst(did)%UDMAP_FILE),"regridweight",gbufw)
239       endif
240       if(lnsizes(my_id+1) .gt. 0) allocate(bufw(lnsizes(my_id+1) ))
241       call pack_decomp_real8(gbufw, ndata, nprocs_map, lnsizes, istart,bufw)
242       if(my_id .eq. io_id) then
243           if(allocated(gbufw))     deallocate(gbufw)
244       endif
247         if(my_id .eq. io_id) then
248            if(allocated(nprocs_map)) deallocate (nprocs_map)
249            if(allocated(istart)) deallocate (istart)
250         endif
251         lnsize = lnsizes(my_id + 1)
252         if(allocated(lnsizes)) deallocate(lnsizes)
253 #else
254        call hydro_stop("FATAL ERROR in UDMP : sequential not defined.")
255 #endif
257     end subroutine readUDMP
259     subroutine UDMP2LOCAL(npid,ix,jx,rtmask, ter_rt_flag)
260         implicit none
261         integer :: i,j,k, ngrids, ix,jx, starti,startj, endi,endj, ii,jj, npid, kk
262         integer, intent(in), dimension(ix,jx) :: rtmask
263         integer, dimension(lnsize) :: lndflag,gridflag , tmpgridflag
264         integer :: ter_rt_flag, m, c
267 !   find ngrids is 0 so that we need to mapping from subsurface runoff.
268 #ifdef MPP_LAND
269         if(left_id .ge. 0) then
270            starti = local_startx_rt  + 1
271         else
272            starti = local_startx_rt
273         endif
274         if(down_id .ge. 0) then
275            startj = local_starty_rt  + 1
276         else
277            startj = local_starty_rt
278         endif
279         if(right_id .ge. 0) then
280            endi = local_startx_rt + ix -2
281         else
282            endi = local_startx_rt + ix -1
283         endif
284         if(up_id .ge. 0) then
285            endj = local_starty_rt + jx -2
286         else
287            endj = local_starty_rt + jx -1
288         endif
289 #else
290         starti = 1
291         startj = 1
292         endi = ix
293         endj = jx
294 #endif
295         gridflag = 0
296         lndflag = 0
298 #ifdef MPP_LAND
299         k = 0
300         do i = 1, lnsize
301            if(bufid(i) .gt. 0) then
302                 if(bufi(i) .ge. starti .and. bufj(i) .ge. startj .and. &
303                     bufi(i) .le. endi   .and. bufj(i) .le. endj) then
304                     if(k .eq. 0) then
305                        k = 1
306                     else
307                        if(bufid(i) .ne. bufid(i-1)) k = k + 1
308                     endif
309                     lndflag(k) = lndflag(k) + 1
310                     if(ter_rt_flag .eq. 1) then
311                         if(rtmask(bufi(i)-local_startx_rt+1,bufj(i)-local_starty_rt+1) .ge. 0) then
312                              gridflag(k) = gridflag(k) + 1
313                         endif
314                     endif
315                  endif
316            endif
317         end do
319 ! decide how many mapping land grids on current domain
320 !       tmpgridflag = gridflag
321 #ifdef MPP_LAND
322 !       call mpp_collect_1d_int_mem(npid,tmpgridflag)
323 #endif
325 ! decide how many user defined links on current domain
326         kk = k
327         LNUMRSL = 0
328         do k = 1, lnsize
329            if(lndflag(k) .gt. 0) LNUMRSL = LNUMRSL + 1
330         enddo
333         if(LNUMRSL .gt. 0) then
335                allocate(LUDRSL(LNUMRSL))
336                allocate( basns_area(LNUMRSL) )
337         else
338 ! When MPI is performed,for every subdomain in each process, all the links 
339 ! are listed and if there is no link in the subdomain then it is calling 
340 ! cleanBuf (memory cleaning purposes), this used to print a warning
341 ! that is not necessary for the user to see it, therefore it is been commented out here
342 !               write(6,*) "Warning: no routing links found."
343                call cleanBuf()
344                return
345         endif
347         kk = 0
348         do k = 1, lnsize
349            if( bufid(k) .ge. 0 ) then
350              if (bufi(k) .ge. starti .and. bufj(k) .ge. startj .and. &
351                     bufi(k) .le. endi   .and. bufj(k) .le. endj ) then
352                  if(kk .eq. 0) then
353                        kk = 1
354                  else
355                        if(bufid(k) .ne. bufid(k-1)) kk = kk + 1
356                  endif
357                  LUDRSL(kk)%myid = bufid(k)
358                  LUDRSL(kk)%ngrids = -999
359                  if(gridflag(kk) .gt. 0) then
360                    LUDRSL(kk)%ngrids = gridflag(kk)
361                    if(.not. allocated(LUDRSL(kk)%weight) ) then
362                          allocate( LUDRSL(kk)%weight(LUDRSL(kk)%ngrids ))
363                          allocate( LUDRSL(kk)%grid_i(LUDRSL(kk)%ngrids ))
364                          allocate( LUDRSL(kk)%grid_j(LUDRSL(kk)%ngrids ))
365                          allocate( LUDRSL(kk)%nodeArea(LUDRSL(kk)%ngrids ))
366                    endif
367                  endif
368 !  define bucket variables
369                  LUDRSL(kk)%ncell = lndflag(kk)
370                  if(.not. allocated(LUDRSL(kk)%cellweight) ) then
371                      allocate( LUDRSL(kk)%cellweight(LUDRSL(kk)%ncell))
372                      allocate( LUDRSL(kk)%cell_i(LUDRSL(kk)%ncell))
373                      allocate( LUDRSL(kk)%cell_j(LUDRSL(kk)%ncell))
374                      allocate( LUDRSL(kk)%cellArea(LUDRSL(kk)%ncell))
375                  endif
376              endif
377            endif
378         enddo
381 ! maping grid_i, grid_j and weight
382         kk = 0
383         m  = 1
384         c  = 1
385         do i = 1, lnsize
386                if( (bufid(i) .ge. 0)  ) then
387                    if(bufi(i) .ge. starti .and. bufj(i) .ge. startj .and. &
388                       bufi(i) .le. endi   .and. bufj(i) .le. endj) then
389                       if(kk .eq. 0) then
390                          kk = 1
391                       else
392                          if(bufid(i) .ne. bufid(i-1)) then
393                              kk = kk + 1
394                              m  = 1
395                              c  = 1
396                          endif
397                       endif
399                       if(LUDRSL(kk)%ngrids .gt. 0) then
400                           if(rtmask(bufi(i)-local_startx_rt+1,bufj(i)-local_starty_rt+1) .ge. 0) then
401                              LUDRSL(kk)%grid_i(m) = bufi(i) - local_startx_rt+1
402                              LUDRSL(kk)%grid_j(m) = bufj(i) - local_starty_rt+1
403                              LUDRSL(kk)%weight(m) = bufw(i)
404                              m  = m  + 1
405                           endif
406                       endif
407 !! begin define bucket variables
408                           LUDRSL(kk)%cell_i(c) = bufi(i) - local_startx_rt+1
409                           LUDRSL(kk)%cell_j(c) = bufj(i) - local_starty_rt+1
410                           LUDRSL(kk)%cellWeight(c) = bufw(i)
411                           c  = c  + 1
412 !! end define bucket variables
413                    endif
414                 endif
415         end do
417         call cleanBuf()
419 #else
420         call hydro_stop("FATAL ERROR in UDMP: Sequential not work.")
421 #endif
423     end subroutine UDMP2LOCAL
425     subroutine cleanBuf()
426         if(allocated(bufi))  deallocate(bufi)
427         if(allocated(bufj))  deallocate(bufj)
428         if(allocated(bufw))  deallocate(bufw)
429         if(allocated(bufid))  deallocate(bufid)
430     end subroutine cleanBuf
432      subroutine get_dimension(fileName, ndata,npid)
433             implicit none
434             character(len=*) fileName
435             integer ncid , iret, ndata,npid, dimid
436 #ifdef MPP_LAND
437             if(my_id .eq. IO_id) then
438 #endif
439                iret = nf_open(fileName, NF_NOWRITE, ncid)
440                if (iret /= 0) then
441                   write(*,'("FATAL ERROR: Problem opening mapping file: ''", A, "''")') &
442                        trim(fileName)
443                   call hydro_stop("In get_dimension() - Problem opening mapping file.")
444                endif
446                iret = nf_inq_dimid(ncid, "polyid", dimid)
448                if (iret /= 0) then
449                   print*, "nf_inq_dimid:  polyid"
450                   call hydro_stop("In get_dimension() - nf_inq_dimid:  polyid")
451                endif
453                iret = nf_inq_dimlen(ncid, dimid, npid)
455                iret = nf_inq_dimid(ncid, "data", dimid)
456                if (iret /= 0) then
457                           print*, "nf_inq_dimid:  data"
458                           call hydro_stop("In get_file_dimension() - nf_inq_dimid:  data")
459                endif
461                iret = nf_inq_dimlen(ncid, dimid, ndata)
462                iret = nf_close(ncid)
463 #ifdef MPP_LAND
464             endif
465             call mpp_land_bcast_int1(ndata)
466             call mpp_land_bcast_int1(npid)
467 #endif
468             return
469      end subroutine get_dimension
471        subroutine get1d_real8(fileName,var_name,out_buff)
472           implicit none
473           integer :: ivar, iret,varid,ncid
474           real*8 out_buff(:)
475           character(len=*), intent(in) :: var_name
476           character(len=*), intent(in) :: fileName
478           iret = nf_open(trim(fileName), NF_NOWRITE, ncid)
479           if (iret .ne. 0) then
480             print*,"failed to open the netcdf file: ",trim(fileName)
481             call hydro_stop("In get1d_real8() - failed to open the netcdf file.")
482             return
483           endif
484           ivar = nf_inq_varid(ncid,trim(var_name),  varid)
485           if(ivar .ne. 0) then
486                write(6,*) "Read Variable Error file: ",trim(fileName)
487                write(6,*) "Read Error: could not find ",trim(var_name)
488                call hydro_stop("In get1d_real8() - failed to read netcdf varialbe name. ")
489           end if
490           iret = nf_get_var_double(ncid, varid, out_buff)
491           iret = nf_close(ncid)
492       end subroutine get1d_real8
494        subroutine get1d_int(fileName,var_name,out_buff)
495           implicit none
496           integer :: ivar, iret,varid,ncid
497           integer out_buff(:)
498           character(len=*), intent(in) :: var_name
499           character(len=*), intent(in) :: fileName
501           iret = nf_open(trim(fileName), NF_NOWRITE, ncid)
502           if (iret .ne. 0) then
503             print*,"FATAL ERROR: Failed to open the netcdf file: ",trim(fileName)
504             call hydro_stop("In get1d_int() -  Failed to open the netcdf file")
505             return
506           endif
507           ivar = nf_inq_varid(ncid,trim(var_name),  varid)
508           if(ivar .ne. 0) then
509                write(6,*) "Read Variable Error file: ",trim(fileName)
510                write(6,*) "Read Error: could not find ",trim(var_name)
511                call hydro_stop("In get1d_int() - failed to read netcdf variable name.")
512           end if
513           iret = nf_get_var_int(ncid, varid, out_buff)
514           iret = nf_close(ncid)
515       end subroutine get1d_int
517     subroutine get1d_int64(fileName,var_name,out_buff)
518         implicit none
519         integer :: ivar, iret,varid,ncid
520         integer(kind=int64) out_buff(:)
521         character(len=*), intent(in) :: var_name
522         character(len=*), intent(in) :: fileName
524         iret = nf_open(trim(fileName), NF_NOWRITE, ncid)
525         if (iret .ne. 0) then
526             print*,"FATAL ERROR: Failed to open the netcdf file: ",trim(fileName)
527             call hydro_stop("In get1d_int() -  Failed to open the netcdf file")
528             return
529         endif
530         ivar = nf_inq_varid(ncid,trim(var_name),  varid)
531         if(ivar .ne. 0) then
532             write(6,*) "Read Variable Error file: ",trim(fileName)
533             write(6,*) "Read Error: could not find ",trim(var_name)
534             call hydro_stop("In get1d_int() - failed to read netcdf variable name.")
535         end if
536         iret = nf_get_var_int64(ncid, varid, out_buff)
537         iret = nf_close(ncid)
538     end subroutine get1d_int64
540       subroutine getUDMP_area(cell_area)
541          implicit none
542          integer i,j,k, m
543          real, dimension(:,:) :: cell_area
544          do k  = 1, LNUMRSL
545             if(LUDRSL(k)%ngrids .gt. 0) then
546                 do m = 1, LUDRSL(k)%ngrids
547                     LUDRSL(k)%nodeArea(m) = cell_area(LUDRSL(k)%grid_i(m),LUDRSL(k)%grid_j(m)) 
548                 enddo
549             endif
550                 do m = 1, LUDRSL(k)%ncell
551                     LUDRSL(k)%cellArea(m) = cell_area(LUDRSL(k)%cell_i(m),LUDRSL(k)%cell_j(m)) 
552                 enddo
554             basns_area(k) = 0
555             do m = 1, LUDRSL(k)%ncell
556                     basns_area(k) = basns_area(k) + &
557                           cell_area(LUDRSL(k)%cell_i(m),LUDRSL(k)%cell_j(m)) * LUDRSL(k)%cellWeight(m) 
558             enddo
560          end do
561       end subroutine getUDMP_area
563       subroutine get_basn_area_nhd(inOut)
564          implicit none
565          real, dimension(:) :: inOut
566          real, dimension(gnpid) :: buf
567 #ifdef MPP_LAND
568          call updateLinkV(basns_area, inOut)
569 #else
570          inOut = basns_area
571 #endif
574       end subroutine get_basn_area_nhd
576       subroutine get_nprocs_map(ix,jx,bufi,bufj,nprocs_map,ndata)
577           implicit none
578           integer,dimension(:)  :: bufi, bufj,nprocs_map
579 !          integer, allocatable, dimension(:) ::  lbufi,lbufj, lmap
580           integer  :: ndata, lsize, ix,jx
581           integer, dimension(ix,jx) :: mask
582           integer, allocatable,dimension(:,:) :: gmask
584         integer :: i,j,k, starti,startj, endi,endj, ii,jj, npid, kk
585 #ifdef MPP_LAND
587           mask = my_id + 1
588           if(my_id .eq. IO_id) allocate(gmask(global_rt_nx, global_rt_ny))
590           call MPP_LAND_COM_INTEGER(mask,IX,JX,99)
591           call write_IO_rt_int(mask, gmask)
593           if(my_id .eq. io_id ) then
594              nprocs_map = -999
595              do i = 1, ndata
596                   if( (bufi(i) .gt. 0 .and. bufi(i) .le. global_rt_nx) .and.  &
597                      (bufj(i) .gt. 0 .and. bufj(i) .le. global_rt_ny) ) then
598                      nprocs_map(i) = gmask(bufi(i), bufj(i))
599                      if( gmask(bufi(i), bufj(i)) .lt. 0) then
600                          write(6,*) "mapping error in gmask : ", bufi(i) ,bufj(i)
601                      endif
602                   else
603                       write(6,*) "no mapping for i,j : ", bufi(i) ,bufj(i)
604                   endif
605              end do
607              if(allocated(gmask)) deallocate(gmask)
608           endif
609 #else
610         call hydro_stop("FATAL ERROR in UDMP: Sequential not work.")
611 #endif
614       end subroutine get_nprocs_map
617 end module module_UDMAP