WAVEWATCH III  beta 0.0.1
wmesmfmd Module Reference

National Unified Prediction Capability (NUOPC) based Earth System Modeling Framework (ESMF) interface module for multi-grid wave model. More...

Functions/Subroutines

subroutine, public setservices (gcomp, rc)
 Wave model ESMF set services. More...
 
subroutine setexport (gcomp, rc)
 Set export fields from internal data structures. More...
 
subroutine createexpmesh (gcomp, rc)
 Create ESMF mesh (unstructured) for export fields. More...
 
subroutine setupimpbmsk (bmskField, impField, missingVal, rc)
 Setup background blending mask field for an import field. More...
 
integer function fieldindex (fnameList, fname, rc)
 Return index associated with field name. More...
 
subroutine calcstokes3d (a, usxField, usyField, rc)
 Calculate 3D Stokes drift current for export. More...
 
subroutine readfromfile (idfld, fldwx, fldwy, time0, timen, rc)
 Read input file to fill unmapped point for regional applications. More...
 

Detailed Description

National Unified Prediction Capability (NUOPC) based Earth System Modeling Framework (ESMF) interface module for multi-grid wave model.

All module variables and types are scoped private by default. The private module variables and types are not listed in this section.

Author
T. J. Campell
J. Meixner
A. J. van der Westhuysen
Date
09-Aug-2017

Function/Subroutine Documentation

◆ calcstokes3d()

subroutine wmesmfmd::calcstokes3d ( real, dimension(nth,nk,0:nseal)  a,
type(esmf_field)  usxField,
type(esmf_field)  usyField,
integer  rc 
)

Calculate 3D Stokes drift current for export.

Parameters
aInput spectra (in par list to change shape)
usxField3D SDC eastward-component export field
usyField3D SDC northward-component export field
rcReturn code
Author
T. J. Campbell
Date
09-Aug-2017

Definition at line 7232 of file wmesmfmd.F90.

7232  !/
7233  !/ +-----------------------------------+
7234  !/ | WAVEWATCH III NOAA/NCEP |
7235  !/ | T. J. Campbell, NRL |
7236  !/ | FORTRAN 90 |
7237  !/ | Last update : 09-Aug-2017 |
7238  !/ +-----------------------------------+
7239  !/
7240  !/ 09-Aug-2017 : Origination. ( version 6.03 )
7241  !/
7242  ! 1. Purpose :
7243  !
7244  ! Calculate 3D Stokes drift current for export
7245  !
7246  ! 2. Method :
7247  !
7248  ! Kenyon, K. E. (1969), J. Geophys. R., Vol 74, No 28, p 6991
7249  !
7250  ! U_vec(z)
7251  ! //
7252  ! = 2 g || ( F(f,theta) k_vec/C cosh(2k(D+z))/sinh(2kD) ) dsig dtheta
7253  ! //
7254  !
7255  ! //
7256  ! = || (Ac(k,theta) sig^2 k_vec/Cg cosh(2k(D+z))/sinh^2(kD) ) dsig dtheta
7257  ! //
7258  !
7259  ! Where:
7260  ! Ac(k,theta) = wave action density
7261  ! k_vec = k*[cos(theta),sin(theta)]
7262  ! D = depth
7263  ! z = height (0 = mean sea level)
7264  !
7265  ! In deep water (kD large): cosh(2k(D+z))/sinh^2(kD) --> 2*exp(2kz)
7266  !
7267  ! 3. Parameters :
7268  !
7269  ! Parameter list
7270  ! ----------------------------------------------------------------
7271  ! a Real I Input spectra (in par list to change shape)
7272  ! usxField Type I/O 3D SDC eastward-component export field
7273  ! usyField Type I/O 3D SDC northward-component export field
7274  ! rc Int O Return code
7275  ! ----------------------------------------------------------------
7276  !
7277  ! 4. Subroutines used :
7278  !
7279  ! Name Type Module Description
7280  ! ----------------------------------------------------------------
7281  ! NONE
7282  ! ----------------------------------------------------------------
7283  !
7284  ! 5. Called by :
7285  !
7286  ! 6. Error messages :
7287  !
7288  ! 7. Remarks :
7289  !
7290  ! 8. Structure :
7291  !
7292  ! 9. Switches :
7293  !
7294  ! 10. Source code :
7295  !
7296  !/ ------------------------------------------------------------------- /
7297  !/
7298  !/ ------------------------------------------------------------------- /
7299  !/ Parameter list
7300  !/
7301  implicit none
7302  real :: a(nth,nk,0:nseal)
7303  type(ESMF_Field) :: usxField
7304  type(ESMF_Field) :: usyField
7305  integer :: rc
7306  !/
7307  !/ ------------------------------------------------------------------- /
7308  !/ Local parameters
7309  !/
7310  real(8), parameter :: zero = 0.0
7311  real(8), parameter :: half = 0.5
7312  real(8), parameter :: one = 1.0
7313  real(8), parameter :: two = 2.0
7314  ! kdmin = 1e-7: sinh(kdmin)**2 ~ 1e-14
7315  real(8), parameter :: kdmin = 1e-7
7316  ! kdmax = 18.0: cosh(2*(kdmax+kz))/sinh(kdmax)**2 - 2*exp(2*kz) < 1e-14
7317  real(8), parameter :: kdmax = 18.0
7318  ! kdmin & kdmax settings used in w3iogomd
7319  real(8), parameter :: kdmin_us3d = 1e-3
7320  real(8), parameter :: kdmax_us3d = 6.0
7321  integer :: isea, jsea, ik, ith, iz
7322  real(8) :: depth
7323  real(8) :: akx, aky, kd, kz, fac1, fac2, fac3
7324  real(8) :: uzx(nz), uzy(nz)
7325  real(8), allocatable :: fack(:)
7326  type(ESMF_Field) :: usxnField, usynField
7327  real(ESMF_KIND_RX), pointer :: usxn(:,:), usyn(:,:)
7328  integer, save :: timeSlice = 1
7329  ! Need this workaround to deal with ESMF_FieldCreate not setting up the
7330  ! Fortran arrays with the ungridded dimension as the first array dimension
7331 #define ESMF_ARBSEQ_WORKAROUND
7332 #ifdef ESMF_ARBSEQ_WORKAROUND
7333  type(ESMF_DistGrid) :: natDistGrid
7334  type(ESMF_Array) :: usxnArray, usynArray
7335 #endif
7336  !
7337  ! -------------------------------------------------------------------- /
7338  !
7339  rc = esmf_success
7340 
7341 #ifdef ESMF_ARBSEQ_WORKAROUND
7342  call esmf_gridget( natgrid, distgrid=natdistgrid, rc=rc )
7343  if (esmf_logfounderror(rc, passthru)) return
7344  usxnarray = esmf_arraycreate( natdistgrid, esmf_typekind_rx, &
7345  distgridtoarraymap=(/2/), undistlbound=(/1/), undistubound=(/nz/), rc=rc )
7346  if (esmf_logfounderror(rc, passthru)) return
7347  usxnfield = esmf_fieldcreate( natgrid, usxnarray, &
7348  gridtofieldmap=(/2,3/), ungriddedlbound=(/1/), ungriddedubound=(/nz/), &
7349  staggerloc=natstaggerloc, rc=rc )
7350  if (esmf_logfounderror(rc, passthru)) return
7351  usynarray = esmf_arraycreate( natdistgrid, esmf_typekind_rx, &
7352  distgridtoarraymap=(/2/), undistlbound=(/1/), undistubound=(/nz/), rc=rc )
7353  if (esmf_logfounderror(rc, passthru)) return
7354  usynfield = esmf_fieldcreate( natgrid, usynarray, &
7355  gridtofieldmap=(/2,3/), ungriddedlbound=(/1/), ungriddedubound=(/nz/), &
7356  staggerloc=natstaggerloc, rc=rc )
7357  if (esmf_logfounderror(rc, passthru)) return
7358 #else
7359  usxnfield = esmf_fieldcreate( natgrid, natarrayspec3d, &
7360  gridtofieldmap=(/2,3/), ungriddedlbound=(/1/), ungriddedubound=(/nz/), &
7361  staggerloc=natstaggerloc, rc=rc )
7362  if (esmf_logfounderror(rc, passthru)) return
7363  usynfield = esmf_fieldcreate( natgrid, natarrayspec3d, &
7364  gridtofieldmap=(/2,3/), ungriddedlbound=(/1/), ungriddedubound=(/nz/), &
7365  staggerloc=natstaggerloc, rc=rc )
7366  if (esmf_logfounderror(rc, passthru)) return
7367 #endif
7368 
7369  call fieldfill( usxnfield, zerovalue, rc=rc )
7370  if (esmf_logfounderror(rc, passthru)) return
7371  call fieldfill( usynfield, zerovalue, rc=rc )
7372  if (esmf_logfounderror(rc, passthru)) return
7373 
7374  if ( natgridislocal ) then
7375 
7376  call esmf_fieldget( usxnfield, farrayptr=usxn, rc=rc )
7377  if (esmf_logfounderror(rc, passthru)) return
7378  call esmf_fieldget( usynfield, farrayptr=usyn, rc=rc )
7379  if (esmf_logfounderror(rc, passthru)) return
7380 
7381  allocate( fack(1:nk) )
7382  fack(1:nk) = dden(1:nk) * sig(1:nk)
7383 
7384  jsea_loop: do jsea = 1,nseal
7385 #ifdef W3_DIST
7386  isea = iaproc + (jsea-1)*naproc
7387 #endif
7388 #ifdef W3_SHRD
7389  isea = jsea
7390 #endif
7391  if ( dw(isea).le.zero ) cycle jsea_loop
7392  depth = max(dmin,dw(isea))
7393  uzx(:) = zero
7394  uzy(:) = zero
7395 #ifdef USE_W3OUTG_FOR_EXPORT
7396  ik_loop: do ik = us3df(2),us3df(3)
7397  fac1 = tpiinv*dsii(ik)
7398  kd = max(kdmin_us3d,wn(ik,isea)*dw(isea))
7399  iz_loop: do iz = 1,nz
7400  if ( dw(isea)+zl(iz).le.zero ) cycle iz_loop
7401  kz = wn(ik,isea)*zl(iz)
7402  if ( kd .lt. kdmax_us3d ) then
7403  fac2 = fac1*cosh(two*max(zero,kd+kz))/cosh(two*kd)
7404  else
7405  fac2 = fac1*exp(two*kz)
7406  endif
7407  uzx(iz) = uzx(iz) + us3d(jsea,ik )*fac2
7408  uzy(iz) = uzy(iz) + us3d(jsea,nk+ik)*fac2
7409  enddo iz_loop
7410  enddo ik_loop
7411 #else
7412  ik_loop: do ik = 1,nk
7413  akx = zero
7414  aky = zero
7415  ith_loop: do ith = 1,nth
7416  akx = akx + a(ith,ik,jsea)*ecos(ith)
7417  aky = aky + a(ith,ik,jsea)*esin(ith)
7418  enddo ith_loop
7419  fac1 = fack(ik)*wn(ik,isea)/cg(ik,isea)
7420  kd = max(kdmin,wn(ik,isea)*depth)
7421  if ( kd .lt. kdmax ) then
7422  fac2 = fac1/sinh(kd)**2
7423  else
7424  fac2 = fac1*two
7425  endif
7426  akx = akx*fac2
7427  aky = aky*fac2
7428  iz_loop: do iz = 1,nz
7429  if ( depth+zl(iz).le.zero ) cycle iz_loop
7430  kz = wn(ik,isea)*zl(iz)
7431  if ( kd .lt. kdmax ) then
7432  fac3 = cosh(two*max(zero,kd+kz))
7433  else
7434  fac3 = exp(two*kz)
7435  endif
7436  uzx(iz) = uzx(iz) + akx*fac3
7437  uzy(iz) = uzy(iz) + aky*fac3
7438  enddo iz_loop
7439  enddo ik_loop
7440 #endif
7441  usxn(:,jsea) = uzx(:)
7442  usyn(:,jsea) = uzy(:)
7443  enddo jsea_loop
7444 
7445  deallocate( fack )
7446 
7447  endif !natGridIsLocal
7448 
7449  call esmf_fieldredist( usxnfield, usxfield, n2erh, rc=rc )
7450  if (esmf_logfounderror(rc, passthru)) return
7451  call esmf_fieldredist( usynfield, usyfield, n2erh, rc=rc )
7452  if (esmf_logfounderror(rc, passthru)) return
7453 
7454 #ifdef ESMF_ARBSEQ_WORKAROUND
7455  call esmf_arraydestroy( usxnarray, rc=rc )
7456  if (esmf_logfounderror(rc, passthru)) return
7457  call esmf_arraydestroy( usynarray, rc=rc )
7458  if (esmf_logfounderror(rc, passthru)) return
7459 #endif
7460  call esmf_fielddestroy( usxnfield, rc=rc )
7461  if (esmf_logfounderror(rc, passthru)) return
7462  call esmf_fielddestroy( usynfield, rc=rc )
7463  if (esmf_logfounderror(rc, passthru)) return
7464 
7465 #ifdef TEST_WMESMFMD_STOKES3D
7466  call esmf_fieldwrite( usxfield, "wmesmfmd_stokes3d_usx.nc", &
7467  overwrite=.true., timeslice=timeslice, rc=rc )
7468  if (esmf_logfounderror(rc, passthru)) return
7469  call esmf_fieldwrite( usyfield, "wmesmfmd_stokes3d_usy.nc", &
7470  overwrite=.true., timeslice=timeslice, rc=rc )
7471  if (esmf_logfounderror(rc, passthru)) return
7472  timeslice = timeslice + 1
7473 #endif
7474  !/
7475  !/ End of CalcStokes3D ----------------------------------------------- /
7476  !/

References w3iogomd::calc_u3stokes(), w3adatmd::cg, w3gdatmd::dden, w3gdatmd::dmin, w3gdatmd::dsii, w3adatmd::dw, w3gdatmd::ecos, w3gdatmd::esin, w3odatmd::iaproc, w3odatmd::naproc, w3gdatmd::sig, constants::tpiinv, w3adatmd::us3d, w3gdatmd::us3df, w3adatmd::ussp, and w3adatmd::wn.

◆ createexpmesh()

subroutine wmesmfmd::createexpmesh ( type(esmf_gridcomp)  gcomp,
integer, intent(out)  rc 
)

Create ESMF mesh (unstructured) for export fields.

Create an ESMF Mesh for export using the unstructured mesh description in W3GDATMD. At present, this export mesh is not domain decomposed, but instead is defined on PET 0 only. (In future, when the unstructured mesh will run on domain decomposition, we will use that decomposition.)

Since the internal parallel data is currently stored accross grid points in a "card deck" fashion, we will define an intermediate native grid, as is done for regular/curvilinear grids, and perform an ESMF regrid to the export mesh. This code segment is taken from T. J. Campbell, and modified to 1D, because the internal data structure for unstructred meshes is an array with dimensions [NX,NY=1].

Parameters
gcompGridded component
[out]rcReturn code
Author
A. J. van der Westhuysen
Date
28-Feb-2018

Definition at line 4554 of file wmesmfmd.F90.

4554  !/
4555  !/ +-----------------------------------+
4556  !/ | WAVEWATCH III NOAA/NCEP |
4557  !/ | A. J. van der Westhuysen |
4558  !/ | FORTRAN 90 |
4559  !/ | Last update : 28-FEB-2018 |
4560  !/ +-----------------------------------+
4561  !/
4562  !/ 28-Feb-2018 : Origination. ( version 6.06 )
4563  !/
4564  ! 1. Purpose :
4565  !
4566  ! Create ESMF mesh (unstructured) for export fields
4567  !
4568  ! 2. Method :
4569  !
4570  ! Create an ESMF Mesh for export using the unstructured mesh description
4571  ! in W3GDATMD. At present, this export mesh is not domain decomposed,
4572  ! but instead is defined on PET 0 only. (In future, when the unstructured
4573  ! mesh will run on domain decomposition, we will use that decomposition.)
4574  !
4575  ! Since the internal parallel data is currently stored accross grid points
4576  ! in a "card deck" fashion, we will define an intermediate native grid, as
4577  ! is done for regular/curvilinear grids, and perform an ESMF regrid to the
4578  ! export mesh. This code segment is taken from T. J. Campbell, and
4579  ! modified to 1D, because the internal data structure for unstructred
4580  ! meshes is an array with dimensions [NX,NY=1].
4581  !
4582  ! 3. Parameters :
4583  !
4584  ! Parameter list
4585  ! ----------------------------------------------------------------
4586  ! gcomp Type I/O Gridded component
4587  ! rc Int. O Return code
4588  ! ----------------------------------------------------------------
4589  !
4590  ! 4. Subroutines used :
4591  !
4592  ! Name Type Module Description
4593  ! ----------------------------------------------------------------
4594  ! NONE
4595  ! ----------------------------------------------------------------
4596  !
4597  ! 5. Called by :
4598  !
4599  ! 6. Error messages :
4600  !
4601  ! 7. Remarks :
4602  !
4603  ! 8. Structure :
4604  !
4605  ! 9. Switches :
4606  !
4607  ! !/SHRD Switch for shared / distributed memory architecture.
4608  ! !/DIST Id.
4609  !
4610  ! 10. Source code :
4611  !
4612  !/ ------------------------------------------------------------------- /
4613  !/
4614 #ifdef W3_PDLIB
4615  use yownodepool, only: npa, iplg, nodes_global
4616  use yowelementpool, only: ne, ielg, ine
4617 #endif
4618  !/
4619  !/ ------------------------------------------------------------------- /
4620  !/ Parameter list
4621  !/
4622  implicit none
4623  type(ESMF_GridComp) :: gcomp
4624  integer,intent(out) :: rc
4625  !/
4626  !/ ------------------------------------------------------------------- /
4627  !/ Local parameters
4628  !/
4629  character(ESMF_MAXSTR) :: cname
4630  character(128) :: msg
4631  integer :: nproc, nxproc, nyproc, n, nfac, irp
4632  real :: gr, rp, pr, diff
4633  integer, parameter :: lde = 0
4634  integer :: ldecnt
4635  integer :: i, j, pos, ix, iy, isea, jsea, iproc
4636  integer :: elb(2), eub(2)
4637  integer(ESMF_KIND_I4), pointer :: iptr(:,:)
4638  real(ESMF_KIND_RX), pointer :: rptrx(:,:), rptry(:,:)
4639  real(ESMF_KIND_RX), pointer :: rptr(:,:)
4640  integer :: arbIndexCount
4641  integer, allocatable :: arbIndexList(:,:)
4642  type(ESMF_Field) :: nField, eField
4643  type(ESMF_Field) :: tmpField
4644  integer(ESMF_KIND_I4), allocatable :: nodeIds(:)
4645  real(ESMF_KIND_R8), allocatable :: nodeCoords(:)
4646  integer(ESMF_KIND_I4), allocatable :: nodeOwners(:)
4647  integer(ESMF_KIND_I4), allocatable :: elemIds(:)
4648  integer(ESMF_KIND_I4), allocatable :: elemTypes(:)
4649  integer(ESMF_KIND_I4), allocatable :: elemConn(:)
4650  !
4651  ! -------------------------------------------------------------------- /
4652  ! Prep
4653  !
4654  rc = esmf_success
4655  if ( noactiveexpfields ) return
4656  call esmf_gridcompget(gcomp, name=cname, rc=rc)
4657  if (esmf_logfounderror(rc, passthru)) return
4658  if (verbosity.gt.0) call esmf_logwrite(trim(cname)// &
4659  ': entered CreateExpMesh', esmf_logmsg_info)
4660  !
4661  ! -------------------------------------------------------------------- /
4662  ! Set flag to indicate that this processor has local native grid storage
4663  !
4664  natgridislocal = iaproc .gt. 0 .and. iaproc .le. naproc
4665 
4666  ! 1. Setup
4667  !
4668  ! 1.a Set grid pointers
4669  !
4670  expgridid = 1 !TODO: only export from grid 1
4671  call w3setg ( expgridid, mdse, mdst )
4672  call w3setw ( expgridid, mdse, mdst )
4673  call w3seta ( expgridid, mdse, mdst )
4674  call w3seti ( expgridid, mdse, mdst )
4675  call w3seto ( expgridid, mdse, mdst )
4676  call wmsetm ( expgridid, mdse, mdst )
4677  natgridid = expgridid
4678  nproc = naproc
4679  !
4680  ! 1.b Set arraySpec, staggerLoc, and indexFlag for native fields.
4681  ! NOTE: For unstructured meshes the native grid is a 1D array (NY=1)
4682  !
4683  call esmf_arrayspecset( natarrayspec1d, rank=1, &
4684  typekind=esmf_typekind_rx, rc=rc )
4685  if (esmf_logfounderror(rc, passthru)) return
4686  natstaggerloc = esmf_staggerloc_center
4687  natindexflag = esmf_index_delocal
4688  !
4689  ! -------------------------------------------------------------------- /
4690  ! 2. Create ESMF mesh for export
4691  !
4692  ! 2.a Create ESMF export mesh
4693  !
4694  ! Allocate and fill the node id array.
4695 #ifdef W3_PDLIB
4696  if ( lpdlib .EQV. .false. ) then
4697 #endif
4698  allocate(nodeids(nx))
4699  do i = 1,nx
4700  nodeids(i)=i
4701  enddo
4702 #ifdef W3_PDLIB
4703  else
4704  ! -------------------------------------------------------------------
4705  ! ESMF Definition: The global id's of the nodes resident on this processor
4706  ! -------------------------------------------------------------------
4707  ! Allocate global node ids, including ghost nodes (npa=np+ng)
4708  allocate(nodeids(npa))
4709  do i = 1,npa
4710  nodeids(i)=iplg(i)
4711  enddo
4712  endif
4713  !
4714  ! call ESMF_LogWrite(trim(cname)//': In CreateExpMesh, nodeIds=', &
4715  ! ESMF_LOGMSG_INFO)
4716  ! do i = 1,npa
4717  ! write(msg,*) trim(cname)//': nodeIds(i)',i, &
4718  ! ' ',nodeIds(i)
4719  ! call ESMF_LogWrite(trim(msg), ESMF_LOGMSG_INFO)
4720  ! enddo
4721 #endif
4722 
4723  ! call ESMF_LogWrite(trim(cname)//': In CreateExpMesh, nodeIds=', &
4724  ! ESMF_LOGMSG_INFO)
4725  ! do i = 1,NX
4726  ! write(msg,*) trim(cname)//': ',i, &
4727  ! ' ',nodeIds(i)
4728  ! call ESMF_LogWrite(trim(msg), ESMF_LOGMSG_INFO)
4729  ! enddo
4730 
4731  ! Allocate and fill node coordinate array.
4732  ! Since this is a 2D Mesh the size is 2x the
4733  ! number of nodes.
4734 #ifdef W3_PDLIB
4735  if ( lpdlib .EQV. .false. ) then
4736 #endif
4737  allocate(nodecoords(2*nx))
4738  do i = 1,nx
4739  do j = 1,2
4740  pos=2*(i-1)+j
4741  if (j == 1) then
4742  nodecoords(pos) = xgrd(1,i)
4743  else
4744  nodecoords(pos) = ygrd(1,i)
4745  endif
4746  enddo
4747  enddo
4748 #ifdef W3_PDLIB
4749  else
4750  ! -------------------------------------------------------------------
4751  ! ESMF Definition: Physical coordinates of the nodes
4752  ! -------------------------------------------------------------------
4753  allocate(nodecoords(2*npa))
4754  do i = 1,npa
4755  do j = 1,2
4756  pos=2*(i-1)+j
4757  if ( j == 1) then
4758  nodecoords(pos) = xgrd(1,iplg(i))
4759  else
4760  nodecoords(pos) = ygrd(1,iplg(i))
4761  endif
4762  enddo
4763  enddo
4764  endif
4765  !
4766  ! call ESMF_LogWrite(trim(cname)//': In CreateExpMesh, nodeCoords=', &
4767  ! ESMF_LOGMSG_INFO)
4768  ! do i = 1,(2*npa)
4769  ! write(msg,*) trim(cname)//': nodeCoords(i)',i, &
4770  ! ' ',nodeCoords(i)
4771  ! call ESMF_LogWrite(trim(msg), ESMF_LOGMSG_INFO)
4772  ! enddo
4773 #endif
4774 
4775  ! call ESMF_LogWrite(trim(cname)//': In CreateExpMesh, nodeCoords=', &
4776  ! ESMF_LOGMSG_INFO)
4777  ! do i = 1,(2*NX)
4778  ! write(msg,*) trim(cname)//': ',i, &
4779  ! ' ',nodeCoords(i)
4780  ! call ESMF_LogWrite(trim(msg), ESMF_LOGMSG_INFO)
4781  ! enddo
4782 
4783  ! Allocate and fill the node owner array.
4784 #ifdef W3_PDLIB
4785  if ( lpdlib .EQV. .false. ) then
4786 #endif
4787  allocate(nodeowners(nx))
4788  nodeowners=0 ! TODO: For now, export everything via PET 0
4789 #ifdef W3_PDLIB
4790  else
4791  ! -------------------------------------------------------------------
4792  ! ESMF Definition: Processor that owns the node
4793  ! -------------------------------------------------------------------
4794  allocate(nodeowners(npa))
4795  nodeowners=nodes_global(iplg(1:npa))%domainID-1
4796  endif
4797  !
4798  ! call ESMF_LogWrite(trim(cname)//': In CreateExpMesh, nodeOwners=', &
4799  ! ESMF_LOGMSG_INFO)
4800  ! do i = 1,npa
4801  ! write(msg,*) trim(cname)//': nodeOwners(i)',i, &
4802  ! ' ',nodeOwners(i)
4803  ! call ESMF_LogWrite(trim(msg), ESMF_LOGMSG_INFO)
4804  ! enddo
4805 #endif
4806 
4807  ! call ESMF_LogWrite(trim(cname)//': In CreateExpMesh, nodeOwners=', &
4808  ! ESMF_LOGMSG_INFO)
4809  ! do i = 1,NX
4810  ! write(msg,*) trim(cname)//': ',i, &
4811  ! ' ',nodeOwners(i)
4812  ! call ESMF_LogWrite(trim(msg), ESMF_LOGMSG_INFO)
4813  ! enddo
4814 
4815  ! Allocate and fill the element id array.
4816 #ifdef W3_PDLIB
4817  if ( lpdlib .EQV. .false. ) then
4818 #endif
4819  allocate(elemids(ntri))
4820  do i = 1,ntri
4821  elemids(i)=i
4822  enddo
4823 #ifdef W3_PDLIB
4824  else
4825  ! -------------------------------------------------------------------
4826  ! ESMF Definition: The global id's of the elements resident on this processor
4827  ! -------------------------------------------------------------------
4828  allocate(elemids(ne))
4829  do i = 1,ne
4830  elemids(i)=ielg(i)
4831  enddo
4832  endif
4833  !
4834  ! call ESMF_LogWrite(trim(cname)//': In CreateExpMesh, elemIds=', &
4835  ! ESMF_LOGMSG_INFO)
4836  ! do i = 1,ne
4837  ! write(msg,*) trim(cname)//': elemIds(i)',i, &
4838  ! ' ',elemIds(i)
4839  ! call ESMF_LogWrite(trim(msg), ESMF_LOGMSG_INFO)
4840  ! enddo
4841 #endif
4842 
4843  ! call ESMF_LogWrite(trim(cname)//': In CreateExpMesh, elemIds=', &
4844  ! ESMF_LOGMSG_INFO)
4845  ! do i = 1,NTRI
4846  ! write(msg,*) trim(cname)//': ',i, &
4847  ! ' ',elemIds(i)
4848  ! call ESMF_LogWrite(trim(msg), ESMF_LOGMSG_INFO)
4849  ! enddo
4850 
4851  ! Allocate and fill the element topology type array.
4852 #ifdef W3_PDLIB
4853  if ( lpdlib .EQV. .false. ) then
4854 #endif
4855  allocate(elemtypes(ntri))
4856  do i = 1,ntri
4857  elemtypes(i)=esmf_meshelemtype_tri
4858  enddo
4859 #ifdef W3_PDLIB
4860  else
4861  ! -------------------------------------------------------------------
4862  ! ESMF Definition: Topology of the given element (one of ESMF_MeshElement)
4863  ! -------------------------------------------------------------------
4864  allocate(elemtypes(ne))
4865  do i = 1,ne
4866  elemtypes(i)=esmf_meshelemtype_tri
4867  enddo
4868  endif
4869  !
4870  ! call ESMF_LogWrite(trim(cname)//': In CreateExpMesh, elemTypes=', &
4871  ! ESMF_LOGMSG_INFO)
4872  ! do i = 1,ne
4873  ! write(msg,*) trim(cname)//': elemTypes(i)',i, &
4874  ! ' ',elemTypes(i)
4875  ! call ESMF_LogWrite(trim(msg), ESMF_LOGMSG_INFO)
4876  ! enddo
4877 #endif
4878 
4879  ! call ESMF_LogWrite(trim(cname)//': In CreateExpMesh, elemTypes=', &
4880  ! ESMF_LOGMSG_INFO)
4881  ! do i = 1,NTRI
4882  ! write(msg,*) trim(cname)//': ',i, &
4883  ! ' ',elemTypes(i)
4884  ! call ESMF_LogWrite(trim(msg), ESMF_LOGMSG_INFO)
4885  ! enddo
4886 
4887  ! Allocate and fill the element connection type array.
4888 #ifdef W3_PDLIB
4889  if ( lpdlib .EQV. .false. ) then
4890 #endif
4891  allocate(elemconn(3*ntri))
4892  do i = 1,ntri
4893  do j = 1,3
4894  pos=3*(i-1)+j
4895  elemconn(pos)=trigp(j,i)
4896  enddo
4897  enddo
4898 #ifdef W3_PDLIB
4899  else
4900  ! -------------------------------------------------------------------
4901  ! ESMF Definition: Connectivity table. The number of entries should
4902  ! be equal to the number of nodes in the given topology. The indices
4903  ! should be the local index (1 based) into the array of nodes that
4904  ! was declared with MeshAddNodes.
4905  ! -------------------------------------------------------------------
4906  ! > INE is local element array. it stores the local node IDs
4907  ! > first index from 1 to 3.
4908  ! > second index from 1 to ne.
4909  allocate(elemconn(3*ne))
4910  do i = 1,ne
4911  do j = 1,3
4912  pos=3*(i-1)+j
4913  elemconn(pos)=ine(j,i)
4914  enddo
4915  enddo
4916  endif
4917  !
4918  ! call ESMF_LogWrite(trim(cname)//': In CreateExpMesh, elemConn=', &
4919  ! ESMF_LOGMSG_INFO)
4920  ! do i = 1,(3*ne)
4921  ! write(msg,*) trim(cname)//': elemConn(i)',i, &
4922  ! ' ',elemConn(i)
4923  ! call ESMF_LogWrite(trim(msg), ESMF_LOGMSG_INFO)
4924  ! enddo
4925 #endif
4926 
4927  ! call ESMF_LogWrite(trim(cname)//': In CreateExpMesh, elemConn=', &
4928  ! ESMF_LOGMSG_INFO)
4929  ! do i = 1,(3*NTRI)
4930  ! write(msg,*) trim(cname)//': ',i, &
4931  ! ' ',elemConn(i)
4932  ! call ESMF_LogWrite(trim(msg), ESMF_LOGMSG_INFO)
4933  ! enddo
4934 
4935  expmesh = esmf_meshcreate( parametricdim=2,spatialdim=2, &
4936  nodeids=nodeids, nodecoords=nodecoords, &
4937  nodeowners=nodeowners, elementids=elemids,&
4938  elementtypes=elemtypes, elementconn=elemconn, &
4939  rc=rc )
4940  if (esmf_logfounderror(rc, passthru)) return
4941 
4942  deallocate(nodeids)
4943  deallocate(nodecoords)
4944  deallocate(nodeowners)
4945  deallocate(elemids)
4946  deallocate(elemtypes)
4947  deallocate(elemconn)
4948  !
4949  ! Set flag to indicate that this processor has local export grid storage
4950  !
4951  !AW if (lpet == 0) then
4952  !AW call ESMF_GridGet( expMesh, localDECount=ldecnt, rc=rc )
4953  !AW if (ESMF_LogFoundError(rc, PASSTHRU)) return
4954  !AW expGridIsLocal = ldecnt.gt.0
4955  !AW endif
4956 
4957  !
4958  ! -------------------------------------------------------------------- /
4959  ! 3. Create ESMF grid with arbitrary domain decomposition to match
4960  ! the native domain decomposition of the non-excluded points
4961  ! Note that the native grid layout is dim1=Y, dim2=X
4962  ! Note that coordinates and mask are not needed since this
4963  ! grid is only used to define fields for a redist operation
4964  !
4965  ! 3.a Set flag to indicate that this processor has local native grid storage
4966  !
4967  natgridislocal = iaproc .gt. 0 .and. iaproc .le. naproc
4968 
4969 #ifdef W3_PDLIB
4970  if ( lpdlib .EQV. .false. ) then
4971 #endif
4972  !
4973  ! 3.b Setup arbitrary sequence index list
4974  !
4975  do ipass = 1,2
4976  if (ipass.eq.2) then
4977  allocate (arbindexlist(arbindexcount,2), stat=rc)
4978  if (esmf_logfoundallocerror(rc, passthru)) return
4979  endif
4980  arbindexcount = 0
4981  ! list local native grid non-excluded points
4982  if ( natgridislocal ) then
4983  do jsea = 1,nseal
4984 #ifdef W3_DIST
4985  isea = iaproc + (jsea-1)*naproc
4986 #endif
4987 #ifdef W3_SHRD
4988  isea = jsea
4989 #endif
4990  arbindexcount = arbindexcount+1
4991  if (ipass.eq.2) then
4992  ix = mapsf(isea,1)
4993  iy = mapsf(isea,2)
4994  ! native grid layout: dim1=Y, dim2=X
4995  arbindexlist(arbindexcount,1) = iy
4996  arbindexlist(arbindexcount,2) = ix
4997  endif
4998  enddo
4999  endif
5000  enddo !ipass
5001  !
5002  ! 3.c Create ESMF native grid
5003  !
5004  natgrid = esmf_gridcreatenoperidim( &
5005  minindex=(/ 1, 1/), &
5006  maxindex=(/ny,nx/), &
5007  coorddep1=(/esmf_dim_arb,esmf_dim_arb/), &
5008  coorddep2=(/esmf_dim_arb,esmf_dim_arb/), &
5009  arbindexcount=arbindexcount, &
5010  arbindexlist=arbindexlist, &
5011  coordtypekind=esmf_typekind_rx, &
5012  coordsys=esmf_coordsys_sph_deg, &
5013  name=trim(cname)//"_native_grid", rc=rc )
5014  if (esmf_logfounderror(rc, passthru)) return
5015  !
5016  ! 3.d Deallocate arbitrary sequence index list
5017  !
5018  deallocate (arbindexlist, stat=rc)
5019  if (esmf_logfounddeallocerror(rc, passthru)) return
5020  !
5021  ! -------------------------------------------------------------------- /
5022  ! 4. Create route handle for redist between native grid domain
5023  ! decomposition and the export grid domain decomposition
5024  !
5025  ! 4.a Create temporary fields
5026  !
5027  nfield = esmf_fieldcreate( natgrid, natarrayspec1d, &
5028  staggerloc=natstaggerloc, rc=rc )
5029  if (esmf_logfounderror(rc, passthru)) return
5030  efield = esmf_fieldcreate(expmesh, typekind=esmf_typekind_rx, rc=rc)
5031  if (esmf_logfounderror(rc, passthru)) return
5032  !
5033  ! 4.b Store route handle
5034  !
5035  call esmf_fieldrediststore( nfield, efield, n2erh, rc=rc )
5036  if (esmf_logfounderror(rc, passthru)) return
5037  !
5038  ! 4.c Clean up
5039  !
5040  call esmf_fielddestroy( nfield, rc=rc )
5041  if (esmf_logfounderror(rc, passthru)) return
5042  call esmf_fielddestroy( efield, rc=rc )
5043  if (esmf_logfounderror(rc, passthru)) return
5044 
5045 #ifdef W3_PDLIB
5046  endif
5047 #endif
5048 
5049  call esmf_logwrite(trim(cname)//': In CreateExpMesh, created expMesh', &
5050  esmf_logmsg_info)
5051 
5052  rc = esmf_success
5053  if (verbosity.gt.0) call esmf_logwrite(trim(cname)// &
5054  ': leaving CreateExpMesh', esmf_logmsg_info)
5055  !/
5056  !/ End of CreateExpMesh ---------------------------------------------- /
5057  !/

References w3odatmd::iaproc, yowelementpool::ielg, yowelementpool::ine, w3adatmd::ipass, yownodepool::iplg, constants::lpdlib, w3gdatmd::mapsf, wmmdatmd::mdse, wmmdatmd::mdst, w3odatmd::naproc, yowelementpool::ne, yownodepool::nodes_global, yownodepool::npa, w3gdatmd::nseal, w3gdatmd::ntri, w3gdatmd::nx, w3gdatmd::ny, w3gdatmd::trigp, w3adatmd::w3seta(), w3gdatmd::w3setg(), w3idatmd::w3seti(), w3odatmd::w3seto(), w3wdatmd::w3setw(), wmmdatmd::wmsetm(), w3gdatmd::xgrd, and w3gdatmd::ygrd.

◆ fieldindex()

integer function wmesmfmd::fieldindex ( character (len=*), dimension(:)  fnameList,
character (len=*)  fname,
integer  rc 
)

Return index associated with field name.

Parameters
fnameListArray of field names
fnameField name
rcReturn code
Returns
indx Returned index of fname
Author
T. J. Campbell
Date
20-Jan-2017

Definition at line 5857 of file wmesmfmd.F90.

5857  !/
5858  !/ +-----------------------------------+
5859  !/ | WAVEWATCH III NOAA/NCEP |
5860  !/ | T. J. Campbell, NRL |
5861  !/ | FORTRAN 90 |
5862  !/ | Last update : 20-Jan-2017 |
5863  !/ +-----------------------------------+
5864  !/
5865  !/ 20-Jan-2017 : Origination. ( version 6.02 )
5866  !/
5867  ! 1. Purpose :
5868  !
5869  ! Return index associated with field name
5870  !
5871  ! 2. Method :
5872  !
5873  ! 3. Parameters :
5874  !
5875  ! Parameter list
5876  ! ----------------------------------------------------------------
5877  ! fnameList StrA I Array of field names
5878  ! fname Str I Field name
5879  ! rc Int. O Return code
5880  ! indx Int I Returned index of fname
5881  ! ----------------------------------------------------------------
5882  !
5883  ! 4. Subroutines used :
5884  !
5885  ! Name Type Module Description
5886  ! ----------------------------------------------------------------
5887  ! NONE
5888  ! ----------------------------------------------------------------
5889  !
5890  ! 5. Called by :
5891  !
5892  ! 6. Error messages :
5893  !
5894  ! 7. Remarks :
5895  !
5896  ! 8. Structure :
5897  !
5898  ! 9. Switches :
5899  !
5900  ! 10. Source code :
5901  !
5902  !/ ------------------------------------------------------------------- /
5903  !/
5904  !/ ------------------------------------------------------------------- /
5905  !/ Parameter list
5906  !/
5907  implicit none
5908  character (len=*) :: fnameList(:)
5909  character (len=*) :: fname
5910  integer :: rc
5911  integer :: indx
5912  !/
5913  !/ ------------------------------------------------------------------- /
5914  !/ Local parameters
5915  !/
5916  integer :: i, check
5917  !
5918  ! -------------------------------------------------------------------- /
5919  ! Find name in fnameList that matches fname
5920  !
5921  check = lbound(fnamelist,1) - 1
5922  indx = check
5923  do i = lbound(fnamelist,1),ubound(fnamelist,1)
5924  if ( trim(fnamelist(i)).eq.trim(fname) ) then
5925  indx = i
5926  exit
5927  endif
5928  enddo
5929  if ( indx.eq.check ) then
5930  call esmf_logseterror(esmf_failure, rctoreturn=rc, &
5931  msg='FieldIndex: input name ('//fname//') not in list')
5932  endif
5933  !/
5934  !/ End of FieldIndex ------------------------------------------------- /
5935  !/

References w3adatmd::aba, w3adatmd::cg, w3adatmd::charn, check(), w3gdatmd::clgtype, w3gdatmd::dden, w3gdatmd::dmin, w3adatmd::dw, constants::dwat, w3gdatmd::ec2, w3gdatmd::ecos, w3gdatmd::es2, w3gdatmd::esc, w3gdatmd::esin, file(), w3gdatmd::fte, constants::grav, w3gdatmd::gtype, w3odatmd::iaproc, yownodepool::iplg, constants::lpdlib, w3gdatmd::mapsf, w3gdatmd::mapsta, w3odatmd::naproc, w3gdatmd::nk, yownodepool::np, w3gdatmd::nseal, w3gdatmd::nth, w3gdatmd::nx, w3gdatmd::ny, w3gdatmd::rlgtype, w3gdatmd::sig, w3gdatmd::sig2, w3adatmd::sxx, w3adatmd::sxy, w3adatmd::syy, constants::tpi, w3adatmd::u10, w3adatmd::u10d, w3adatmd::uba, w3adatmd::ubd, w3gdatmd::ungtype, w3wdatmd::ust, w3wdatmd::va, w3src3md::w3spr3(), w3src4md::w3spr4(), w3adatmd::wn, and w3adatmd::wnmean.

◆ readfromfile()

subroutine wmesmfmd::readfromfile ( character(len=3), intent(inout)  idfld,
type(esmf_field), intent(inout)  fldwx,
type(esmf_field), intent(inout)  fldwy,
integer, dimension(2), intent(in)  time0,
integer, dimension(2), intent(in)  timen,
integer, intent(inout), optional  rc 
)

Read input file to fill unmapped point for regional applications.

Parameters
[in,out]idfldField name
[in,out]fldwx2D eastward-component of field
[in,out]fldwy2D northward-component of field
[in]time0Time stamp for current time
[in]timenTime stamp for end time
[in,out]rcReturn code
Author
U. Turuncoglu
Date
18-May-2021

Definition at line 7709 of file wmesmfmd.F90.

7709  !/
7710  !/ +-----------------------------------+
7711  !/ | WAVEWATCH III NOAA/NCEP |
7712  !/ | U. Turuncoglu |
7713  !/ | FORTRAN 90 |
7714  !/ | Last update : 18-May-2021 |
7715  !/ +-----------------------------------+
7716  !/
7717  !/ 18-May-2021 : Origination. ( version 7.13 )
7718  !/
7719  ! 1. Purpose :
7720  !
7721  ! Read input file to fill unmapped point for regional applications
7722  !
7723  ! 2. Method :
7724  !
7725  ! 3. Parameters :
7726  !
7727  ! Parameter list
7728  ! ----------------------------------------------------------------
7729  ! idfld Str I/O Field name
7730  ! fldwx Type I/O 2D eastward-component of field
7731  ! fldwy Type I/O 2D northward-component of field
7732  ! time0 Int I Time stamp for current time
7733  ! timen Int I Time stamp for end time
7734  ! rc Int I/O Return code
7735  ! ----------------------------------------------------------------
7736  !
7737  ! 4. Subroutines used :
7738  !
7739  ! Name Type Module Description
7740  ! ----------------------------------------------------------------
7741  ! NONE
7742  ! ----------------------------------------------------------------
7743  !
7744  ! 5. Called by :
7745  !
7746  ! 6. Error messages :
7747  !
7748  ! 7. Remarks :
7749  !
7750  ! 8. Structure :
7751  !
7752  ! 9. Switches :
7753  !
7754  ! 10. Source code :
7755  !
7756  !/ ------------------------------------------------------------------- /
7757  !/
7758  USE w3fldsmd, ONLY: w3fldo, w3fldg
7759  USE wmunitmd, ONLY: wmuget, wmuset
7760  IMPLICIT NONE
7761  !/ ------------------------------------------------------------------- /
7762  !/ Parameter list
7763  !/
7764  character(len=3), intent(inout) :: idfld
7765  type(ESMF_Field), intent(inout) :: fldwx
7766  type(ESMF_Field), intent(inout) :: fldwy
7767  integer, intent(in) :: time0(2)
7768  integer, intent(in) :: timen(2)
7769  integer, intent(inout), optional :: rc
7770  !/
7771  !/ ------------------------------------------------------------------- /
7772  !/ Local parameters
7773  !/
7774  integer :: ierr, tw0l(2), twnl(2), lb(2), ub(2)
7775  real :: wx0l(nx,ny), wy0l(nx,ny)
7776  real :: wxnl(nx,ny), wynl(nx,ny)
7777  real :: dt0l(nx,ny), dtnl(nx,ny)
7778  real(ESMF_KIND_RX), pointer :: dptr(:,:)
7779  integer :: mdse = 6
7780  integer :: mdst = 10
7781  integer, save :: mdsf
7782  character(256) :: logmsg
7783  logical :: flagsc = .false.
7784  integer, parameter :: lde = 0
7785  logical, save :: firstCall = .true.
7786  character(len=13) :: tsstr
7787  character(len=3) :: tsfld
7788  integer :: nxt, nyt, gtypet, filler(3), tideflag
7789 #if defined(TEST_WMESMFMD) || defined(TEST_WMESMFMD_READFROMFILE)
7790  integer, save :: timeSlice = 1
7791 #endif
7792  !
7793  ! -------------------------------------------------------------------- /
7794  !
7795  rc = esmf_success
7796 
7797  if (firstcall) then
7798  ! assign unit number for input file
7799  call wmuget(mdse, mdst, mdsf, 'INP')
7800  call wmuset(mdse, mdst, mdsf, .true., desc='Input data file')
7801 
7802  ! open file
7803  call w3fldo('READ', idfld, mdsf, mdst, mdse, nx, ny, gtype, ierr)
7804  if (ierr.ne.0) then
7805  write(logmsg,*) "Error in opening "//idfld//", iostat = ", ierr
7806  call esmf_logwrite(trim(logmsg), esmf_logmsg_error)
7807  rc = esmf_failure
7808  return
7809  endif
7810 
7811  firstcall = .false.
7812  end if
7813 
7814  ! init variables
7815  wx0l = 0.0
7816  wy0l = 0.0
7817  dt0l = 0.0
7818  wxnl = 0.0
7819  wynl = 0.0
7820  dtnl = 0.0
7821 
7822  ! need to rewind to the begining of the file to access
7823  ! data of requested date correctly
7824  rewind(mdsf)
7825 
7826  ! read header information
7827  ! this was inside of w3fldo call but since we are opening file
7828  ! once and rewinding, the header need to be read
7829  read(mdsf, iostat=ierr) tsstr, tsfld, nxt, nyt, &
7830  gtypet, filler(1:2), tideflag
7831 
7832  ! read input
7833  call w3fldg('READ', idfld, mdsf, mdst, mdse, nx, ny, &
7834  nx, ny, time0, timen, tw0l, wx0l, wy0l, dt0l, twnl, &
7835  wxnl, wynl, dtnl, ierr, flagsc)
7836 
7837  ! fill fields with data belong to current time
7838  if ( impgridislocal ) then
7839  call esmf_fieldget(fldwx, localde=lde, farrayptr=dptr, &
7840  exclusivelbound=lb, exclusiveubound=ub, rc=rc)
7841  if (esmf_logfounderror(rc, passthru)) return
7842  dptr(lb(1):ub(1),lb(2):ub(2)) = wx0l(lb(1):ub(1),lb(2):ub(2))
7843  if (associated(dptr)) nullify(dptr)
7844  call esmf_fieldget(fldwy, localde=lde, farrayptr=dptr, &
7845  exclusivelbound=lb, exclusiveubound=ub, rc=rc)
7846  if (esmf_logfounderror(rc, passthru)) return
7847  dptr(lb(1):ub(1),lb(2):ub(2)) = wy0l(lb(1):ub(1),lb(2):ub(2))
7848  if (associated(dptr)) nullify(dptr)
7849  end if
7850 
7851 #if defined(TEST_WMESMFMD) || defined(TEST_WMESMFMD_READFROMFILE)
7852  write(logmsg,*) 'time0 = ', time0(1), time0(2)
7853  call esmf_logwrite(trim(logmsg), esmf_logmsg_info, rc=rc)
7854  if (esmf_logfounderror(rc, passthru)) return
7855  write(logmsg,*) 'timen = ', timen(1), timen(2)
7856  call esmf_logwrite(trim(logmsg), esmf_logmsg_info, rc=rc)
7857  if (esmf_logfounderror(rc, passthru)) return
7858  write(logmsg,*) 'tw0 = ', tw0l(1), tw0l(2)
7859  call esmf_logwrite(trim(logmsg), esmf_logmsg_info, rc=rc)
7860  if (esmf_logfounderror(rc, passthru)) return
7861  write(logmsg,*) 'twn = ', twnl(1), twnl(2)
7862  call esmf_logwrite(trim(logmsg), esmf_logmsg_info, rc=rc)
7863  if (esmf_logfounderror(rc, passthru)) return
7864  write(logmsg,*) 'wx0 min, max = ', minval(wx0l), maxval(wx0l)
7865  call esmf_logwrite(trim(logmsg), esmf_logmsg_info, rc=rc)
7866  if (esmf_logfounderror(rc, passthru)) return
7867  write(logmsg,*) 'wy0 min, max = ', minval(wy0l), maxval(wy0l)
7868  call esmf_logwrite(trim(logmsg), esmf_logmsg_info, rc=rc)
7869  if (esmf_logfounderror(rc, passthru)) return
7870  write(logmsg,*) 'wxn min, max = ', minval(wxnl), maxval(wxnl)
7871  call esmf_logwrite(trim(logmsg), esmf_logmsg_info, rc=rc)
7872  if (esmf_logfounderror(rc, passthru)) return
7873  write(logmsg,*) 'wyn min, max = ', minval(wynl), maxval(wynl)
7874  call esmf_logwrite(trim(logmsg), esmf_logmsg_info, rc=rc)
7875  if (esmf_logfounderror(rc, passthru)) return
7876  call esmf_fieldwrite( fldwx, &
7877  "wmesmfmd_read_wx0.nc", &
7878  overwrite=.true., timeslice=timeslice, rc=rc )
7879  if (esmf_logfounderror(rc, passthru)) return
7880  call esmf_fieldwrite( fldwy, &
7881  "wmesmfmd_read_wy0.nc", &
7882  overwrite=.true., timeslice=timeslice, rc=rc )
7883  if (esmf_logfounderror(rc, passthru)) return
7884  timeslice = timeslice + 1
7885 #endif
7886  !/
7887  !/ End of ReadFromFile ------------------------------------------- /
7888  !/

References w3gdatmd::gtype, w3fldsmd::w3fldg(), w3fldsmd::w3fldo(), wmunitmd::wmuget(), and wmunitmd::wmuset().

◆ setexport()

subroutine wmesmfmd::setexport ( type(esmf_gridcomp)  gcomp,
integer, intent(out)  rc 
)

Set export fields from internal data structures.

Parameters
gcompGridded component
[out]rcReturn code
Author
T. J. Campbell
Date
09-Aug-2017

Definition at line 2821 of file wmesmfmd.F90.

2821  !/
2822  !/ +-----------------------------------+
2823  !/ | WAVEWATCH III NOAA/NCEP |
2824  !/ | T. J. Campbell, NRL |
2825  !/ | FORTRAN 90 |
2826  !/ | Last update : 09-Aug-2017 |
2827  !/ +-----------------------------------+
2828  !/
2829  !/ 20-Jan-2017 : Origination. ( version 6.02 )
2830  !/ 09-Aug-2017 : Add ocean forcing export fields ( version 6.03 )
2831  !/
2832  ! 1. Purpose :
2833  !
2834  ! Set export fields from internal data structures
2835  !
2836  ! 2. Method :
2837  !
2838  ! 3. Parameters :
2839  !
2840  ! Parameter list
2841  ! ----------------------------------------------------------------
2842  ! gcomp Type I/O Gridded component
2843  ! rc Int. O Return code
2844  ! ----------------------------------------------------------------
2845  !
2846  ! 4. Subroutines used :
2847  !
2848  ! Name Type Module Description
2849  ! ----------------------------------------------------------------
2850  ! NONE
2851  ! ----------------------------------------------------------------
2852  !
2853  ! 5. Called by :
2854  !
2855  ! 6. Error messages :
2856  !
2857  ! 7. Remarks :
2858  !
2859  ! 8. Structure :
2860  !
2861  ! 9. Switches :
2862  !
2863  ! 10. Source code :
2864  !
2865  !/ ------------------------------------------------------------------- /
2866  !/
2867  !/ ------------------------------------------------------------------- /
2868  !/ Parameter list
2869  !/
2870  implicit none
2871  type(ESMF_GridComp) :: gcomp
2872  integer,intent(out) :: rc
2873  !/
2874  !/ ------------------------------------------------------------------- /
2875  !/ Local parameters
2876  !/
2877  character(ESMF_MAXSTR) :: cname
2878  integer, parameter :: iwt=8
2879  real(8) :: wstime, wftime
2880  integer :: i1, i2, i3, i4, i5, i6
2881  logical :: flpart = .false., floutg = .false., floutg2 = .true.
2882 #if defined(TEST_WMESMFMD) || defined(TEST_WMESMFMD_SETEXPORT)
2883  type(ESMF_State) :: dumpState
2884  integer, save :: timeSlice = 1
2885 #endif
2886  real(ESMF_KIND_R8), pointer :: farrayptr(:,:)
2887  !
2888  ! -------------------------------------------------------------------- /
2889  ! Prep
2890  !
2891  rc = esmf_success
2892  if ( noactiveexpfields ) return
2893  call esmf_vmwtime(wstime)
2894  call esmf_gridcompget(gcomp, name=cname, rc=rc)
2895  if (esmf_logfounderror(rc, passthru)) return
2896  if (verbosity.gt.0) call esmf_logwrite(trim(cname)// &
2897  ': entered SetExport', esmf_logmsg_info)
2898  !
2899  ! -------------------------------------------------------------------- /
2900  ! Setup
2901  !
2902  call w3setg ( expgridid, mdse, mdst )
2903  call w3setw ( expgridid, mdse, mdst )
2904  call w3seta ( expgridid, mdse, mdst )
2905  call w3seti ( expgridid, mdse, mdst )
2906  call w3seto ( expgridid, mdse, mdst )
2907  call wmsetm ( expgridid, mdse, mdst )
2908 #ifdef USE_W3OUTG_FOR_EXPORT
2909  if ( natgridislocal ) call w3outg( va, flpart, floutg, floutg2 )
2910 #endif
2911  !
2912  ! -------------------------------------------------------------------- /
2913  ! Charnock
2914  !
2915  i1 = fieldindex( expfieldname, 'charno', rc )
2916  if (esmf_logfounderror(rc, passthru)) return
2917  if ( expfieldactive(i1) ) then
2918  call calccharnk( expfield(i1), rc )
2919  if (esmf_logfounderror(rc, passthru)) return
2920  endif
2921  !
2922  ! -------------------------------------------------------------------- /
2923  ! Surface Roughness
2924  !
2925  i1 = fieldindex( expfieldname, 'z0rlen', rc )
2926  if (esmf_logfounderror(rc, passthru)) return
2927  if ( expfieldactive(i1) ) then
2928  call calcroughl( expfield(i1), rc )
2929  if (esmf_logfounderror(rc, passthru)) return
2930  endif
2931  !
2932  ! -------------------------------------------------------------------- /
2933  ! Stokes Drift 3D
2934  !
2935  i1 = fieldindex( expfieldname, 'uscurr', rc )
2936  if (esmf_logfounderror(rc, passthru)) return
2937  i2 = fieldindex( expfieldname, 'vscurr', rc )
2938  if (esmf_logfounderror(rc, passthru)) return
2939  if ( expfieldactive(i1) .and. &
2940  expfieldactive(i2) ) then
2941  call calcstokes3d( va, expfield(i1), expfield(i2), rc )
2942  if (esmf_logfounderror(rc, passthru)) return
2943  endif
2944  !
2945  ! -------------------------------------------------------------------- /
2946  ! Partitioned Stokes Drift 3 2D fields
2947  !
2948  i1 = fieldindex( expfieldname, 'x1pstk', rc )
2949  if (esmf_logfounderror(rc, passthru)) return
2950  i2 = fieldindex( expfieldname, 'y1pstk', rc )
2951  if (esmf_logfounderror(rc, passthru)) return
2952  i3 = fieldindex( expfieldname, 'x2pstk', rc )
2953  if (esmf_logfounderror(rc, passthru)) return
2954  i4 = fieldindex( expfieldname, 'y2pstk', rc )
2955  if (esmf_logfounderror(rc, passthru)) return
2956  i5 = fieldindex( expfieldname, 'x3pstk', rc )
2957  if (esmf_logfounderror(rc, passthru)) return
2958  i6 = fieldindex( expfieldname, 'y3pstk', rc )
2959  if (esmf_logfounderror(rc, passthru)) return
2960  if ( expfieldactive(i1) .and. &
2961  expfieldactive(i2) .and. &
2962  expfieldactive(i3) .and. &
2963  expfieldactive(i4) .and. &
2964  expfieldactive(i5) .and. &
2965  expfieldactive(i6) ) then
2966  call calcpstokes( va, expfield(i1), expfield(i2), expfield(i3), &
2967  expfield(i4), expfield(i5), expfield(i6), rc )
2968  if (esmf_logfounderror(rc, passthru)) return
2969  endif
2970  !
2971  ! -------------------------------------------------------------------- /
2972  ! Bottom Currents
2973  !
2974  i1 = fieldindex( expfieldname, 'wbcuru', rc )
2975  if (esmf_logfounderror(rc, passthru)) return
2976  i2 = fieldindex( expfieldname, 'wbcurv', rc )
2977  if (esmf_logfounderror(rc, passthru)) return
2978  i3 = fieldindex( expfieldname, 'wbcurp', rc )
2979  if (esmf_logfounderror(rc, passthru)) return
2980  if ( expfieldactive(i1) .and. &
2981  expfieldactive(i2) .and. &
2982  expfieldactive(i3) ) then
2983  call calcbotcur( va, expfield(i1), expfield(i2), expfield(i3), rc )
2984  if (esmf_logfounderror(rc, passthru)) return
2985  endif
2986  !
2987  ! -------------------------------------------------------------------- /
2988  ! Radiation stresses 2D
2989  !
2990  i1 = fieldindex( expfieldname, 'wavsuu', rc )
2991  if (esmf_logfounderror(rc, passthru)) return
2992  i2 = fieldindex( expfieldname, 'wavsuv', rc )
2993  if (esmf_logfounderror(rc, passthru)) return
2994  i3 = fieldindex( expfieldname, 'wavsvv', rc )
2995  if (esmf_logfounderror(rc, passthru)) return
2996  if ( expfieldactive(i1) .and. &
2997  expfieldactive(i2) .and. &
2998  expfieldactive(i3) ) then
2999  call calcradstr2d( va, expfield(i1), expfield(i2), expfield(i3), rc )
3000  if (esmf_logfounderror(rc, passthru)) return
3001  endif
3002  !
3003  ! -------------------------------------------------------------------- /
3004  ! cpl_scalars - grid sizes
3005  !
3006  if (med_present) then
3007  i1 = fieldindex( expfieldname, trim(flds_scalar_name), rc=rc )
3008  if (esmf_logfounderror(rc, passthru)) return
3009  if ( expfieldactive(i1) ) then
3010  call esmf_fieldget(expfield(i1), farrayptr=farrayptr, rc=rc)
3011  if (esmf_logfounderror(rc, passthru)) return
3012  if (flds_scalar_index_nx > 0 .and. flds_scalar_index_nx < flds_scalar_num) then
3013  farrayptr(flds_scalar_index_nx,1) = dble(nx)
3014  endif
3015  if (flds_scalar_index_ny > 0 .and. flds_scalar_index_ny < flds_scalar_num) then
3016  farrayptr(flds_scalar_index_ny,1) = dble(ny)
3017  endif
3018  endif
3019  endif
3020  !
3021  ! -------------------------------------------------------------------- /
3022  ! Post
3023  !
3024 #if defined(TEST_WMESMFMD) || defined(TEST_WMESMFMD_SETEXPORT)
3025  call nuopc_modelget(gcomp, exportstate=dumpstate, rc=rc)
3026  if (esmf_logfounderror(rc, passthru)) return
3027  call nuopc_write(dumpstate, overwrite=.true., &
3028  filenameprefix="field_"//trim(cname)//"_export_", &
3029  timeslice=timeslice, relaxedflag=.true., rc=rc)
3030  if (esmf_logfounderror(rc, passthru)) return
3031  timeslice = timeslice + 1
3032 #endif
3033  rc = esmf_success
3034  call esmf_vmwtime(wftime)
3035  wtime(iwt) = wtime(iwt) + wftime - wstime
3036  wtcnt(iwt) = wtcnt(iwt) + 1
3037  if (verbosity.gt.0) call esmf_logwrite(trim(cname)// &
3038  ': leaving SetExport', esmf_logmsg_info)
3039  !/
3040  !/ End of SetExport -------------------------------------------------- /
3041  !/

References file(), constants::lpdlib, w3odatmd::naproc, wmmdatmd::nmproc, w3gdatmd::ntri, w3gdatmd::nx, w3gdatmd::ny, and w3gdatmd::trigp.

◆ setservices()

subroutine, public wmesmfmd::setservices ( type(esmf_gridcomp)  gcomp,
integer, intent(out)  rc 
)

Wave model ESMF set services.

Parameters
gcompGridded component.
[out]rcReturn code.
Author
T. J. Campbell
Date
20-Jan-2017

Definition at line 352 of file wmesmfmd.F90.

352  !/
353  !/ +-----------------------------------+
354  !/ | WAVEWATCH III NOAA/NCEP |
355  !/ | T. J. Campbell, NRL |
356  !/ | FORTRAN 90 |
357  !/ | Last update : 20-Jan-2017 |
358  !/ +-----------------------------------+
359  !/
360  !/ 20-Jan-2017 : Origination. ( version 6.02 )
361  !/
362  ! 1. Purpose :
363  !
364  ! Wave model ESMF set services.
365  !
366  ! 2. Method :
367  !
368  ! 3. Parameters :
369  !
370  ! Parameter list
371  ! ----------------------------------------------------------------
372  ! gcomp Type I/O Gridded component
373  ! rc Int. O Return code
374  ! ----------------------------------------------------------------
375  !
376  ! 4. Subroutines used :
377  !
378  ! Name Type Module Description
379  ! ----------------------------------------------------------------
380  ! InitializeP0 Subr. WMESMFMD Wave model NUOPC/ESMF Initialize phase 0
381  ! InitializeP1 Subr. WMESMFMD Wave model NUOPC/ESMF Initialize phase 1
382  ! InitializeP3 Subr. WMESMFMD Wave model NUOPC/ESMF Initialize phase 3
383  ! Finalize Subr. WMESMFMD Wave model NUOPC/ESMF Finalize
384  ! DataInitialize Subr. WMESMFMD Wave model NUOPC/ESMF Data Initialize
385  ! ModelAdvance Subr. WMESMFMD Wave model NUOPC/ESMF Model Advance
386  ! ----------------------------------------------------------------
387  !
388  ! 5. Called by :
389  !
390  ! 6. Error messages :
391  !
392  ! 7. Remarks :
393  !
394  ! 8. Structure :
395  !
396  ! 9. Switches :
397  !
398  ! 10. Source code :
399  !
400  !/ ------------------------------------------------------------------- /
401  !/
402  !/ ------------------------------------------------------------------- /
403  !/ Parameter list
404  !/
405  implicit none
406  type(ESMF_GridComp) :: gcomp
407  integer,intent(out) :: rc
408  !/
409  !/ ------------------------------------------------------------------- /
410  !/ Local parameters
411  !/
412  !NONE
413  !
414  ! -------------------------------------------------------------------- /
415  ! Prep
416  !
417  rc = esmf_success
418 
419  ! --- Initialize wallclock timers
420 
421  wtnam( 1) = 'InitializeP0'
422  wtnam( 2) = 'InitializeP1'
423  wtnam( 3) = 'InitializeP3'
424  wtnam( 4) = 'DataInitialize'
425  wtnam( 5) = 'ModelAdvance'
426  wtnam( 6) = 'Finalize'
427  wtnam( 7) = 'GetImport'
428  wtnam( 8) = 'SetExport'
429  wtnam( 9) = 'FieldGather'
430  wtnam(10) = 'FieldFill'
431  wtcnt( :) = 0
432  wtime( :) = 0d0
433  !
434  ! -------------------------------------------------------------------- /
435  ! 1. NUOPC model component will register the generic methods
436  !
437  call nuopc_compderive(gcomp, parent_setservices, rc=rc)
438  if (esmf_logfounderror(rc, passthru)) return
439  !
440  ! -------------------------------------------------------------------- /
441  ! 2. Set model entry points
442  !
443  ! --- Initialize - phase 0 (requires use of ESMF method)
444 
445  call esmf_gridcompsetentrypoint(gcomp, esmf_method_initialize, &
446  userroutine=initializep0, phase=0, rc=rc)
447  if (esmf_logfounderror(rc, passthru)) return
448 
449  ! --- Set entry points for initialize methods
450 
451  ! >= IPDv03 supports satisfying inter-model data dependencies and
452  ! the transfer of ESMF Grid & Mesh objects between Model and/or
453  ! Mediator components during initialization
454  ! IPDv03p1: advertise import & export fields
455  call nuopc_compsetentrypoint(gcomp, esmf_method_initialize, &
456  phaselabellist=(/"IPDv03p1"/), userroutine=initializep1, rc=rc)
457  if (esmf_logfounderror(rc, passthru)) return
458  ! IPDv03p2: unspecified by NUOPC -- not required
459  ! IPDv03p3: realize import & export fields
460  call nuopc_compsetentrypoint(gcomp, esmf_method_initialize, &
461  phaselabellist=(/"IPDv03p3"/), userroutine=initializep3, rc=rc)
462  if (esmf_logfounderror(rc, passthru)) return
463  ! IPDv03p4: relevant for TransferActionGeomObject=="accept"
464  ! IPDv03p5: relevant for TransferActionGeomObject=="accept"
465  ! IPDv03p6: check compatibility of fields connected status
466  ! IPDv03p7: handle field data initialization
467 
468  !
469  ! -------------------------------------------------------------------- /
470  ! 3. Register specializing methods
471  !
472  ! --- Model initialize export data method
473 
474  call nuopc_compspecialize(gcomp, speclabel=label_datainitialize, &
475  specroutine=datainitialize, rc=rc)
476  if (esmf_logfounderror(rc, passthru)) return
477 
478  ! --- Model checkImport method (overriding default)
479 
480  call esmf_methodremove(gcomp, label_checkimport, rc=rc)
481  if (esmf_logfounderror(rc, passthru)) return
482  call nuopc_compspecialize(gcomp, speclabel=label_checkimport, &
483  specroutine=nuopc_noop, rc=rc)
484  if (esmf_logfounderror(rc, passthru)) return
485 
486  ! --- Model advance method
487 
488  call nuopc_compspecialize(gcomp, speclabel=label_advance, &
489  specroutine=modeladvance, rc=rc)
490  if (esmf_logfounderror(rc, passthru)) return
491 
492  ! --- Model finalize method
493 
494  call nuopc_compspecialize(gcomp, speclabel=label_finalize, &
495  specroutine=finalize, rc=rc)
496  if (esmf_logfounderror(rc, passthru)) return
497  !
498  ! -------------------------------------------------------------------- /
499  ! Post
500  !
501  rc = esmf_success
502  !/
503  !/ End of SetServices ------------------------------------------------ /
504  !/

References wmmdatmd::etime, w3gdatmd::grids, wmmdatmd::improc, constants::is_esmf_component, wmmdatmd::nmproc, wmmdatmd::nmpscr, wmmdatmd::stime, w3wdatmd::time, w3adatmd::us3d, and w3iogomd::w3outg().

◆ setupimpbmsk()

subroutine wmesmfmd::setupimpbmsk ( type(esmf_field)  bmskField,
type(esmf_field)  impField,
real(esmf_kind_rx)  missingVal,
integer, optional  rc 
)

Setup background blending mask field for an import field.

Parameters
bmskFieldBlending mask field
impFieldImport field
missingValMissing value
rcReturn code
Author
T. J. Campbell
Date
09-Aug-2017

Definition at line 5074 of file wmesmfmd.F90.

5074  !/
5075  !/ +-----------------------------------+
5076  !/ | WAVEWATCH III NOAA/NCEP |
5077  !/ | T. J. Campbell, NRL |
5078  !/ | FORTRAN 90 |
5079  !/ | Last update : 09-Aug-2017 |
5080  !/ +-----------------------------------+
5081  !/
5082  !/ 09-Aug-2017 : Origination. ( version 6.03 )
5083  !/
5084  ! 1. Purpose :
5085  !
5086  ! Setup background blending mask field for an import field
5087  !
5088  ! 2. Method :
5089  !
5090  ! 3. Parameters :
5091  !
5092  ! Parameter list
5093  ! ----------------------------------------------------------------
5094  ! bmskField Type I/O blending mask field
5095  ! impField Type I import field
5096  ! missingVal Real I missing value
5097  ! rc Int O Return code
5098  ! ----------------------------------------------------------------
5099  !
5100  ! 4. Subroutines used :
5101  !
5102  ! Name Type Module Description
5103  ! ----------------------------------------------------------------
5104  ! NONE
5105  ! ----------------------------------------------------------------
5106  !
5107  ! 5. Called by :
5108  !
5109  ! 6. Error messages :
5110  !
5111  ! 7. Remarks :
5112  !
5113  ! 8. Structure :
5114  !
5115  ! 9. Switches :
5116  !
5117  ! 10. Source code :
5118  !
5119  !/ ------------------------------------------------------------------- /
5120  !/
5121  !/ ------------------------------------------------------------------- /
5122  !/ Parameter list
5123  !/
5124  implicit none
5125  type(ESMF_Field) :: bmskField
5126  type(ESMF_Field) :: impField
5127  real(ESMF_KIND_RX) :: missingVal
5128  integer, optional :: rc
5129  !/
5130  !/ ------------------------------------------------------------------- /
5131  !/ Local parameters
5132  !/
5133  real(ESMF_KIND_RX), parameter :: zero = 0.0
5134  real(ESMF_KIND_RX), parameter :: half = 0.5
5135  real(ESMF_KIND_RX), parameter :: one = 1.0
5136  integer, parameter :: nsig = imphalowidth-1
5137  integer, parameter :: niter = 10
5138  integer, parameter :: iter0 = 1-niter
5139  integer, parameter :: lde = 0
5140  integer :: iter, i, j, ii, jj, k, l
5141  integer :: elb(2), eub(2)
5142  integer :: tlb(2), tub(2)
5143  real(ESMF_KIND_RX), pointer :: mptr(:,:), dptr(:,:)
5144  type(ESMF_Field) :: cmskField
5145  real(ESMF_KIND_RX), pointer :: bmsk(:,:), cmsk(:,:)
5146  real(ESMF_KIND_RX) :: bsum, wsum
5147  real(ESMF_KIND_RX) :: wflt(-nsig:nsig,-nsig:nsig)
5148  character(ESMF_MAXSTR) :: fnm
5149 #if defined(TEST_WMESMFMD) || defined(TEST_WMESMFMD_SETUPIMPBMSK)
5150  integer :: timeSlice
5151 #endif
5152  !
5153  ! -------------------------------------------------------------------- /
5154  ! Initialize filter
5155  !
5156  if (present(rc)) rc = esmf_success
5157 
5158  do l = -nsig,nsig
5159  do k = -nsig,nsig
5160  wflt(k,l) = exp( -half*( real(k,esmf_kind_rx)**2 &
5161  + real(l,esmf_kind_rx)**2 ) )
5162  enddo
5163  enddo
5164  !
5165  ! -------------------------------------------------------------------- /
5166  ! Set up fields and pointers
5167  !
5168  call esmf_fieldget( impfield, name=fnm, rc=rc )
5169  if (esmf_logfounderror(rc, passthru)) return
5170 
5171  cmskfield = esmf_fieldcreate( impgrid, imparrayspec2d, &
5172  totallwidth=imphalolwidth, totaluwidth=imphalouwidth, &
5173  staggerloc=impstaggerloc, indexflag=impindexflag, rc=rc )
5174  if (esmf_logfounderror(rc, passthru)) return
5175 
5176  if ( impgridislocal ) then
5177  call esmf_fieldget( impmask, localde=lde, farrayptr=mptr, &
5178  exclusivelbound=elb, exclusiveubound=eub, &
5179  totallbound=tlb, totalubound=tub, rc=rc )
5180  if (esmf_logfounderror(rc, passthru)) return
5181  call esmf_fieldget( impfield, localde=lde, farrayptr=dptr, rc=rc )
5182  if (esmf_logfounderror(rc, passthru)) return
5183  call esmf_fieldget( bmskfield, localde=lde, farrayptr=bmsk, rc=rc )
5184  if (esmf_logfounderror(rc, passthru)) return
5185  call esmf_fieldget( cmskfield, localde=lde, farrayptr=cmsk, rc=rc )
5186  if (esmf_logfounderror(rc, passthru)) return
5187  endif
5188  !
5189  ! -------------------------------------------------------------------- /
5190  ! Create blending mask
5191  !
5192  if ( impgridislocal ) then
5193  do j = elb(2),eub(2)
5194  do i = elb(1),eub(1)
5195  if ( nint(dptr(i,j)).eq.nint(missingval) ) then
5196  bmsk(i,j) = zero
5197  else
5198  bmsk(i,j) = one
5199  endif
5200  cmsk(i,j) = bmsk(i,j)
5201  enddo
5202  enddo
5203  endif
5204 #if defined(TEST_WMESMFMD) || defined(TEST_WMESMFMD_SETUPIMPBMSK)
5205  timeslice = 1
5206  call esmf_fieldwrite( bmskfield, &
5207  "wmesmfmd_setupimpbmsk_"//trim(fnm)//".nc", &
5208  overwrite=.true., timeslice=timeslice, rc=rc )
5209  if (esmf_logfounderror(rc, passthru)) return
5210 #endif
5211 
5212  iter_loop: do iter = iter0,niter
5213 
5214  call esmf_fieldhalo( bmskfield, imphalorh, rc=rc )
5215  if (esmf_logfounderror(rc, passthru)) return
5216 
5217  if ( impgridislocal ) then
5218 
5219  j_loop: do j = elb(2),eub(2)
5220  i_loop: do i = elb(1),eub(1)
5221  if ( nint(mptr(i,j)).eq.maskvalueland ) cycle i_loop
5222  if ( nint(dptr(i,j)).eq.nint(missingval) ) cycle i_loop
5223 
5224  if (iter.le.0) then
5225 
5226  ! initialize blending zone to zero
5227  l_loop0: do l = -1,1
5228  jj = j + l
5229  if ( jj.lt.tlb(2).or.jj.gt.tub(2) ) cycle l_loop0
5230  k_loop0: do k = -1,1
5231  ii = i + k
5232  if ( ii.lt.tlb(1).or.ii.gt.tub(1) ) cycle k_loop0
5233  if ( nint(mptr(ii,jj)).eq.maskvalueland ) cycle k_loop0
5234  if ( bmsk(ii,jj).eq.zero ) cmsk(i,j) = zero
5235  enddo k_loop0
5236  enddo l_loop0
5237 
5238  else
5239 
5240  ! iterate filter over blending zone
5241  bsum = zero
5242  wsum = zero
5243  l_loop: do l = -nsig,nsig
5244  jj = j + l
5245  if ( jj.lt.tlb(2).or.jj.gt.tub(2) ) cycle l_loop
5246  k_loop: do k = -nsig,nsig
5247  ii = i + k
5248  if ( ii.lt.tlb(1).or.ii.gt.tub(1) ) cycle k_loop
5249  if ( nint(mptr(ii,jj)).eq.maskvalueland ) cycle k_loop
5250  bsum = bsum + wflt(k,l)*bmsk(ii,jj)
5251  wsum = wsum + wflt(k,l)
5252  enddo k_loop
5253  enddo l_loop
5254  if ( wsum.gt.zero ) cmsk(i,j) = bsum/wsum
5255 
5256  endif
5257 
5258  enddo i_loop
5259  enddo j_loop
5260 
5261  do j = elb(2),eub(2)
5262  do i = elb(1),eub(1)
5263  if ( nint(mptr(i,j)).eq.maskvalueland ) cycle
5264  bmsk(i,j) = cmsk(i,j)
5265  enddo
5266  enddo
5267 
5268  endif
5269 #if defined(TEST_WMESMFMD) || defined(TEST_WMESMFMD_SETUPIMPBMSK)
5270  timeslice = timeslice + 1
5271  call esmf_fieldwrite( bmskfield, &
5272  "wmesmfmd_setupimpbmsk_"//trim(fnm)//".nc", &
5273  overwrite=.true., timeslice=timeslice, rc=rc )
5274  if (esmf_logfounderror(rc, passthru)) return
5275 #endif
5276 
5277  enddo iter_loop
5278  !
5279  ! -------------------------------------------------------------------- /
5280  ! Clean up
5281  !
5282  call esmf_fielddestroy( cmskfield, rc=rc )
5283  if (esmf_logfounderror(rc, passthru)) return
5284  !/
5285  !/ End of SetupImpBmsk ----------------------------------------------- /
5286  !/

References w3gdatmd::clgtype, w3gdatmd::gtype, yownodepool::iplg, constants::lpdlib, yownodepool::np, yowrankmodule::rank, w3gdatmd::rlgtype, and w3gdatmd::ungtype.

scrip_constants::half
real(kind=scrip_r8), parameter half
Definition: scrip_constants.f:50
yowelementpool::ielg
integer, dimension(:), allocatable, target, public ielg
global element array.
Definition: yowelementpool.F90:65
yowelementpool
Definition: yowelementpool.F90:38
yowrankmodule::rank
type(t_rank), dimension(:), allocatable, public rank
Provides access to some information of all threads e.g.
Definition: yowrankModule.F90:68
w3adatmd::us3d
real, dimension(:,:), pointer us3d
Definition: w3adatmd.F90:612
w3adatmd::cg
real, dimension(:,:), pointer cg
Definition: w3adatmd.F90:575
w3adatmd::dw
real, dimension(:), pointer dw
Definition: w3adatmd.F90:584
w3adatmd::ipass
integer, pointer ipass
Definition: w3adatmd.F90:686
w3odatmd::iaproc
integer, pointer iaproc
Definition: w3odatmd.F90:457
yownodepool::iplg
integer, dimension(:), allocatable, public iplg
Node local to global mapping.
Definition: yownodepool.F90:116
yownodepool::npa
integer, public npa
number of ghost + resident nodes this partition holds
Definition: yownodepool.F90:99
wmunitmd::wmuset
subroutine wmuset(NDSE, NDST, NDS, FLAG, TYPE, NAME, DESC)
Directly set information for a unit number in the data structure.
Definition: wmunitmd.F90:497
wmunitmd::wmuget
subroutine wmuget(NDSE, NDST, NDS, TYPE, NR)
Find a free unit number for a given file type.
Definition: wmunitmd.F90:667
scrip_constants::one
real(kind=scrip_r8), parameter one
Definition: scrip_constants.f:50
w3wdatmd::va
real, dimension(:,:), pointer va
Definition: w3wdatmd.F90:183
yowelementpool::ne
integer, public ne
number of local elements
Definition: yowelementpool.F90:48
w3adatmd::w3seta
subroutine w3seta(IMOD, NDSE, NDST)
Select one of the WAVEWATCH III grids / models.
Definition: w3adatmd.F90:2645
scrip_constants::zero
real(kind=scrip_r8), parameter zero
Definition: scrip_constants.f:50
yownodepool
Has data that belong to nodes.
Definition: yownodepool.F90:39
constants::lpdlib
logical lpdlib
LPDLIB Logical for using the PDLIB library.
Definition: constants.F90:101
w3wdatmd::w3setw
subroutine w3setw(IMOD, NDSE, NDST)
Select one of the WAVEWATCH III grids / models.
Definition: w3wdatmd.F90:660
w3idatmd::flagsc
logical, dimension(:), pointer flagsc
Definition: w3idatmd.F90:260
scrip_constants::two
real(kind=scrip_r8), parameter two
Definition: scrip_constants.f:50
w3odatmd::naproc
integer, pointer naproc
Definition: w3odatmd.F90:457
w3idatmd::w3seti
subroutine w3seti(IMOD, NDSE, NDST)
Select one of the WAVEWATCH III grids / models.
Definition: w3idatmd.F90:819
w3adatmd::wn
real, dimension(:,:), pointer wn
Definition: w3adatmd.F90:575
yowelementpool::ine
integer, dimension(:,:), allocatable, target, public ine
number of elements of the augmented domain
Definition: yowelementpool.F90:56
yownodepool::nodes_global
type(t_node), dimension(:), allocatable, target, public nodes_global
all nodes with their data.
Definition: yownodepool.F90:103
w3adatmd::iter
integer, dimension(:,:), pointer iter
Definition: w3adatmd.F90:654
w3fldsmd::w3fldg
subroutine w3fldg(INXOUT, IDFLD, NDS, NDST, NDSE, MX, MY, NX, NY, T0, TN, TF0, FX0, FY0, FA0, TFN, FXN, FYN, FAN, IERR, FLAGSC ifdef W3_OASIS
Definition: w3fldsmd.F90:958
w3fldsmd::w3fldo
subroutine w3fldo(INXOUT, IDFLD, NDS, NDST, NDSE, NX, NY, GTYPE, IERR, FEXT, FPRE, FHDR, TIDEFLAGIN)
Definition: w3fldsmd.F90:90
w3iogomd::w3outg
subroutine w3outg(A, FLPART, FLOUTG, FLOUTG2)
Fill necessary arrays with gridded data for output.
Definition: w3iogomd.F90:1198
check
subroutine check(status)
N/A.
Definition: ww3_systrk.F90:1299
wmunitmd
Dynamic assignement of unit numbers for the multi-grid wave model.
Definition: wmunitmd.F90:18
w3fldsmd
Definition: w3fldsmd.F90:3