MODFLOW 6
version 6.6.0.dev0
USGS Modular Hydrologic Model
Main Page
Modules
Data Types List
Files
File List
File Members
SubcellRect.f90
Go to the documentation of this file.
1
module
subcellrectmodule
2
3
use
subcellmodule
,
only
:
subcelltype
4
implicit none
5
6
private
7
public
::
subcellrecttype
8
public
::
create_subcell_rect
9
10
type
,
extends
(
subcelltype
) ::
subcellrecttype
11
private
12
double precision
,
public
:: sinrot
!< sine of rotation angle for local (x, y)
13
double precision
,
public
:: cosrot
!< cosine of rotation angle for local (x, y)
14
double precision
,
public
:: xorigin
!< cell x origin for local (x, y)
15
double precision
,
public
:: yorigin
!< cell y origin for local (x, y)
16
double precision
,
public
:: zorigin
!< cell z origin for local z
17
double precision
,
public
:: dx, dy, dz
!< subcell dimensions
18
double precision
,
public
:: vx1, vx2, vy1, vy2, vz1, vz2
!< subcell face velocities
19
contains
20
procedure
,
public
::
destroy
=>
destroy_subcell_rect
!< destructor for the subcell
21
procedure
,
public
::
init
=>
init_subcell_rect
!< initializes the rectangular subcell
22
end type
subcellrecttype
23
24
contains
25
26
!> @brief Create a new rectangular subcell
27
subroutine
create_subcell_rect
(subcell)
28
type
(
subcellrecttype
),
pointer
:: subcell
29
allocate
(subcell)
30
allocate
(subcell%type)
31
subcell%type =
'subcellrect'
32
end subroutine
create_subcell_rect
33
34
!> @brief Destructor for a rectangular subcell
35
subroutine
destroy_subcell_rect
(this)
36
class
(
subcellrecttype
),
intent(inout)
:: this
37
deallocate
(this%type)
38
end subroutine
destroy_subcell_rect
39
40
!> @brief Initialize a rectangular subcell
41
subroutine
init_subcell_rect
(this)
42
class
(
subcellrecttype
),
intent(inout)
:: this
43
end subroutine
init_subcell_rect
44
45
end module
subcellrectmodule
subcellmodule::destroy
Definition:
Subcell.f90:20
subcellmodule::init
Definition:
Subcell.f90:25
subcellmodule
Definition:
Subcell.f90:1
subcellrectmodule
Definition:
SubcellRect.f90:1
subcellrectmodule::destroy_subcell_rect
subroutine destroy_subcell_rect(this)
Destructor for a rectangular subcell.
Definition:
SubcellRect.f90:36
subcellrectmodule::init_subcell_rect
subroutine init_subcell_rect(this)
Initialize a rectangular subcell.
Definition:
SubcellRect.f90:42
subcellrectmodule::create_subcell_rect
subroutine, public create_subcell_rect(subcell)
Create a new rectangular subcell.
Definition:
SubcellRect.f90:28
subcellmodule::subcelltype
A subcell of a cell.
Definition:
Subcell.f90:9
subcellrectmodule::subcellrecttype
Definition:
SubcellRect.f90:10
src
Solution
ParticleTracker
SubcellRect.f90
Generated on Wed Nov 13 2024 11:15:22 for MODFLOW 6 by
1.9.1