FE-Project
|
module FElib / Mesh / Cubed-sphere 2D domain More...
Data Types | |
type | meshcubedspheredom2d |
Functions/Subroutines | |
subroutine | meshcubedspheredom2d_init (this, negx, negy, rplanet, refelem, nlocalmeshperprc, nproc, myrank) |
subroutine, public | meshcubedspheredom2d_check_division_params (nprcx_lc, nprcy_lc, prc_num, local_mesh_num_global, call_prc_abort) |
subroutine, public | meshcubedspheredom2d_setuplocaldom (lcmesh, tileid, panelid, i, j, nprcx, nprcy, dom_xmin, dom_xmax, dom_ymin, dom_ymax, planet_radius, nex, ney) |
module FElib / Mesh / Cubed-sphere 2D domain
subroutine scale_mesh_cubedspheredom2d::meshcubedspheredom2d_init | ( | class(meshcubedspheredom2d), intent(inout) | this, |
integer, intent(in) | negx, | ||
integer, intent(in) | negy, | ||
real(rp), intent(in) | rplanet, | ||
type(quadrilateralelement), intent(in), target | refelem, | ||
integer, intent(in) | nlocalmeshperprc, | ||
integer, intent(in), optional | nproc, | ||
integer, intent(in), optional | myrank ) |
Definition at line 78 of file scale_mesh_cubedspheredom2d.F90.
References scale_mesh_base2d::meshbase2d_dimtypeid_x, scale_mesh_base2d::meshbase2d_dimtypeid_xy, scale_mesh_base2d::meshbase2d_dimtypeid_xyt, scale_mesh_base2d::meshbase2d_dimtypeid_y, scale_mesh_base2d::meshbase2d_final(), scale_mesh_base2d::meshbase2d_init(), meshcubedspheredom2d_check_division_params(), and meshcubedspheredom2d_setuplocaldom().
subroutine, public scale_mesh_cubedspheredom2d::meshcubedspheredom2d_check_division_params | ( | integer, intent(out) | nprcx_lc, |
integer, intent(out) | nprcy_lc, | ||
integer, intent(in) | prc_num, | ||
integer, intent(in) | local_mesh_num_global, | ||
logical, intent(in), optional | call_prc_abort ) |
Definition at line 193 of file scale_mesh_cubedspheredom2d.F90.
Referenced by scale_mesh_cubedspheredom2d::meshcubedspheredom2d::assigndomid(), meshcubedspheredom2d_init(), and scale_mesh_cubedspheredom3d::meshcubedspheredom3d_init().
subroutine, public scale_mesh_cubedspheredom2d::meshcubedspheredom2d_setuplocaldom | ( | type(localmesh2d), intent(inout) | lcmesh, |
integer, intent(in) | tileid, | ||
integer, intent(in) | panelid, | ||
integer, intent(in) | i, | ||
integer, intent(in) | j, | ||
integer, intent(in) | nprcx, | ||
integer, intent(in) | nprcy, | ||
real(rp), intent(in) | dom_xmin, | ||
real(rp), intent(in) | dom_xmax, | ||
real(rp), intent(in) | dom_ymin, | ||
real(rp), intent(in) | dom_ymax, | ||
real(rp), intent(in) | planet_radius, | ||
integer, intent(in) | nex, | ||
integer, intent(in) | ney ) |
Definition at line 249 of file scale_mesh_cubedspheredom2d.F90.
References scale_localmesh_base::bctype_interior, scale_cubedsphere_coord_cnv::cubedspherecoordcnv_cs2lonlatpos(), scale_cubedsphere_coord_cnv::cubedspherecoordcnv_getmetric(), scale_mesh_base2d::meshbase2d_setgeometricinfo(), and scale_meshutil_cubedsphere2d::meshutilcubedsphere2d_buildglobalmap().
Referenced by scale_mesh_cubedspheredom2d::meshcubedspheredom2d::assigndomid(), meshcubedspheredom2d_init(), and scale_mesh_cubedspheredom3d::meshcubedspheredom3d_init().