MODFLOW 6  version 6.6.0.dev0
USGS Modular Hydrologic Model
CellRect.f90
Go to the documentation of this file.
2 
3  use cellmodule, only: celltype
5  implicit none
6 
7  private
8  public :: cellrecttype
9  public :: create_cell_rect
10 
11  type, extends(celltype) :: cellrecttype
12  private
13  double precision, public :: dx ! dimension of cell in local x direction
14  double precision, public :: dy ! dimension of cell in local y direction
15  double precision, public :: dz ! dimension of cell in z direction
16 
17  double precision, public :: sinrot ! sine of rotation angle for local (x, y)
18  double precision, public :: cosrot ! cosine of rotation angle for local (x, y)
19 
20  double precision, public :: vx1 ! west-boundary local-x velocity
21  double precision, public :: vx2 ! east-boundary local-x velocity
22  double precision, public :: vy1 ! south-boundary local-y velocity
23  double precision, public :: vy2 ! north-boundary local-y velocity
24  double precision, public :: vz1 ! bottom-boundary z velocity
25  double precision, public :: vz2 ! top-boundary z velocity
26 
27  integer, public :: ipvorigin ! origin vertex
28  double precision, public :: xorigin ! model x origin for local (x, y)
29  double precision, public :: yorigin ! model y origin for local (x, y)
30  double precision, public :: zorigin ! model z origin for local z
31 
32  contains
33  procedure :: destroy => destroy_rect ! destructor for the cell
34  end type cellrecttype
35 
36 contains
37 
38  !> @brief Create a new rectangular cell
39  subroutine create_cell_rect(cell)
40  type(cellrecttype), pointer :: cell
41  allocate (cell)
42  call create_defn(cell%defn)
43  allocate (cell%type)
44  cell%type = 'rect'
45  end subroutine create_cell_rect
46 
47  !> @brief Destroy the rectangular cell
48  subroutine destroy_rect(this)
49  class(cellrecttype), intent(inout) :: this
50  deallocate (this%defn)
51  deallocate (this%type)
52  end subroutine destroy_rect
53 
54 end module cellrectmodule
subroutine, public create_defn(cellDefn)
Create a new cell definition object.
Definition: CellDefn.f90:41
subroutine, public create_cell_rect(cell)
Create a new rectangular cell.
Definition: CellRect.f90:40
subroutine destroy_rect(this)
Destroy the rectangular cell.
Definition: CellRect.f90:49
Base grid cell definition.
Definition: CellDefn.f90:10
Base type for grid cells of a concrete type. Contains a cell-definition which is information shared b...
Definition: Cell.f90:13