Class for a generic grid defined in bands of latitude. More...
Data Types | |
type | band_grid |
The grid is defined as bands. We store the positions of the boundaries Latitude is defined in degree in [-90,90] and the longitude is in [0,360] Also, the ghost cells are already contained in the grid (and not in a separate array) More... | |
Public Member Functions | |
subroutine | new_band_grid (this, n_lat, dlat, nb_lon, dlon, lat_min, lon_min, ghost_cells) |
Default constructor. Does not set number of ghosts. More... | |
subroutine | new_band_grid_copy (this, other) |
Copy constructor. More... | |
subroutine | allocate_winds (size_merid, this) |
Must compte size of meridional winds before (need partition info) More... | |
subroutine | allocate_advection_data_band_grid (size_merid, this) |
Allocate data used in advection : slope, fluxes, mass. More... | |
integer function | size_winds_zonal (this) |
Size of zonal winds. More... | |
integer function | sum_nb_interfaces_merid (i, this) |
Sums the number of meridional interfaces up to line i (south latitude) More... | |
integer function | sum_nb_interfaces_zonal (i, this) |
Sums the number of zonal interfaces up to line i (south latitude) More... | |
integer function | nb_interfaces_lat (i, this) |
Number of cell interfaces between latitude line i and i+1. Only inside the current grid. Assume i and i+1 exists. More... | |
logical function | is_latitude_just_before_equator_band_grid (i, this) |
Returns true if the latitude line i is just before the equator Concurrent of is_latitude_just_before_equator Uses float. TODO : optimize. More... | |
logical function | is_first_on_band_band_grid (this) |
Concurrent version of is_first_on_band for a grid. Uses float. TODO : optimise. More... | |
logical function | is_last_on_band_band_grid (this) |
Concurrent version of is_last_on_band for a grid. Uses float. TODO : optimise. More... | |
subroutine | init_advection_data (this) |
Initializes arrays for advection. More... | |
integer function | nb_cells_band_grid (this) |
Number of cells except ghost cells. More... | |
subroutine | set_first_cells_band_grid (this) |
Set the position of first cells on each latitude line. More... | |
subroutine | allocate_most_data (this, nb_lat) |
Allocate most of the data. More... | |
subroutine | new_band_grid_from_pos (this, i_pos, j_pos, first_lat, height, width, nb_resized, zone, ghost_cells) |
Create a grid based on its index in a (lat-lon) matrix. Grid is in most case square but the last band with a square number of partition can sometimes be rectangular hence the height and width. However, only the height varies, while the width is constant. We then need the first latitude. More... | |
subroutine | set_data_for_neighbours (this, list_firsts, first_lat, zone) |
Store the position of the first north and south neighbour line, and the number of east and west ghost cells. More... | |
subroutine | set_north_south_offsets (i_lat, k, this, zone, list_firsts) |
Set the north and south offset. More... | |
logical function | offsets_sector_extrem (not_first, not_last, k, this, list_firsts, zone) |
subroutine | set_nb_east_west_ghosts (i_lat, k, this, zone, list_firsts) |
Set the number of east and west ghost cells at line k in the partition. The number of ghost cells must be initialized to 0. More... | |
integer function | nb_ghosts_side (is_first, k, list_firsts, i_lat, zone, this) |
Compute the number of ghosts on a side (either east or west) for north/south. More... | |
integer function | cur_cell_from_first (k, is_first, list_firsts, this) |
Returns either the first or last cell. More... | |
integer function | north_extrema_cell (cur, is_first, list_firsts, i_lat, zone) |
Return either the leftmost or rightmost north neighbour. More... | |
integer function | south_extrema_cell (cur, is_first, list_firsts, i_lat, zone) |
Return either the leftmost or rightmost north neighbour. More... | |
subroutine | add_east_west_ghosts (this) |
Correct nb_lon and lon_min for east/west ghost cells. More... | |
subroutine | new_band_grid_leftover (this, i_pos, j_pos, height, width, n_left, zone, ghost_cells) |
Create a grid at position i on the leftover band. Also finishes to set some. More... | |
integer function | global_to_local_cell (j_cur, i, sector, this) |
Convert global cell position to local Float operations and comparison. More... | |
subroutine | get_tracer_pointer (this, ptr, datatype, tracer) |
subroutine | get_temperature_pointer (this, ptr, array) |
subroutine | get_winds_pointer (this, ptr, datatype, array) |
integer function | get_nb_lon (i, this) |
Returns number of cells at latitude i (in the grid), without the ghost cells. More... | |
logical function | has_north_ghost_cells (this) |
Tell us if there are ghost cells on the north boundary line. More... | |
logical function | has_east_ghost_cells (this) |
Tell us if there are ghost cells on the east boundary line. More... | |
logical function | has_south_ghost_cells (this) |
Tell us if there are ghost cells on the south boundary line. More... | |
logical function | has_west_ghost_cells (this) |
Tell us if there are ghost cells on the west boundary line. More... | |
subroutine | set_longitude_data (this, k, i_pos, j_pos, cur_dlon, nb_resized, width, nb_cells, zone, list_firsts) |
Set longitude data (nb_lon, d_lon, lon_min). But does not contains east/west ghost cells Does not work for the leftover band. Another function does that : set_north_south_ghost in partition_class Also, weast/east ghost cells are added later. More... | |
subroutine | set_longitude_data_ghost (k, this, nb_resized, width, i_pos, j_corr, nb_cells, zone, list_firsts) |
Compute the lon_min and nb_lon for ghost cells. The strategy is simple, we compute the neighbours for the extremities Only works on the first sector. More... | |
subroutine | set_longitude_data_ghost_leftover (k, this, width, i_pos, j_corr, nb_cells, n_left, zone, list_firsts) |
Compute the lon_min and nb_lon for ghost cells on leftover band. Same strategy as "normal" bands. Only works on the first sector. More... | |
subroutine | find_first_length_cell_leftover (first, length, this, width, j_pos, nb_cells, n_left) |
Does not work for ghost cells ! (but used by its computation) Only works on the first sector. More... | |
subroutine | compute_longitude_data_ghost (k, first, last, this, width, i_pos, j_corr, nb_cells, zone, k_next, is_north, nb_cells_next, list_firsts) |
Does the true computing for nb_lon and lon_min (ghost cells). Needs the neighbour first and last cells. More... | |
subroutine | set_longitude_data_leftover (this, k, i_pos, j_pos, cur_dlon, width, nb_cells, n_left, zone, list_firsts) |
Set longitude data for leftover band. More... | |
subroutine | find_first_length_cell (first, length, k, this, nb_resized, width, i_pos, j_corr, nb_cells, zone) |
Compute the first and last cell on latitude k. It does not take west/ghost cells into account. Does not work for leftover band Does not work for ghost cells ! (but used by its computation) Only works on the first sector. More... | |
subroutine | correct_with_east_west_ghost_cells (this) |
Correct the partition extremities by adding east and west ghost cells Does not work for north/south ghost cells ! (but used by its computation) Only works on the first sector. More... | |
integer function | nb_ghost_cells_lat (k, this) |
Computes the total number of ghost cells laterally (east-west) More... | |
logical function | is_latitude_just_after_equator_float (i, this, zone) |
Returns true if the latitude line i is just after the equator Use float ! More... | |
subroutine | find_east_weast_corr (corr_east, corr_west, k, this, zone, i, j) |
Find the correction for minimal longitude. For non-north/south ghost cells, the test is simple. For north-south ghost cells, there are 2 cases. More... | |
subroutine | find_lon_min_south_ghost_last (lon_min, grid, i, width, dlon, i_lat, zone) |
Compute the minimal longitude for south ghost cells on the leftover band (south hemisphere). For that we compute the north neighbour of the west extremity. More... | |
subroutine | allocate_mpi_borders (this, nb) |
Allocate array for mpi types for borders. More... | |
subroutine | new_mpi_global_type (val, offsets, blocklens, this, nb) |
Set the MPI type for ratio or slope field as array have the same dimensions. More... | |
subroutine | new_mpi_ghost_cells (val, offsets, blocklens, this, nb, direction, start_prev, end_prev) |
Create an MPI type for ghost cells and assing it to val. More... | |
subroutine | create_mpi_border_zwinds (this, direction) |
subroutine | start_indice_ghost (start, this, side) |
Update/search starting indice of the data array for the ghost cells, with the following disposition. North/south are converted it from the local position (on a band) to a global (for the grid), while east/west are created
N | |W | | E | || |S | |. More... | |
subroutine | get_start_border_in (this, start, side) |
Get starting indice of the data array for the border inside the partition. More... | |
subroutine | new_indexed_mpi (nb, blocklens, offsets, mpi_type) |
Creates and commit an new indexed mpi type for an array of double precision. More... | |
subroutine | get_center_band_grid (this, mid_lat, mid_lon) |
Find the center of the grid, that is the middle latitude and the middle longitude at this latitude. More... | |
subroutine | new_mpi_band_grid (grid) |
MPI user-defined datatype : create the datatype from a grid The grid must have all of its array allocated ! More... | |
subroutine | print_band_grid (this) |
Print grid features (for debug) More... | |
subroutine | print_nb_ghosts (i, this) |
Print conditionnaly the current number of east/west ghost cells at line i. More... | |
subroutine | print_offsets (i, this) |
Print conditionnaly the current north/south offset at line i. More... | |
subroutine | write_ratio_band_grid (this, tracer, id_file) |
Write the cell ratio and coordinates of the cell center dlat is supposed constant Writing is done from the north. More... | |
subroutine | write_id_band_grid (this, id_part, id_file) |
Write the partition id and the cell corners. More... | |
subroutine | interior_lat_indices (i_start, i_end, this) |
Give the starting indices and minimal latitude for the interior cells (i.e without the ghost cells) More... | |
integer function | first_interior_lat (this) |
First latitude indice for interior cells. More... | |
integer function | last_interior_lat (this) |
Last latitude indice for interior cells. More... | |
double precision function | cell_interface_length (i, j, i_neighb, j_neighb, this) |
Compute interface length between cell (i,j) and (i_neighb, j_neighb) in degree It is a length on the sphere, not a difference of longitudes. More... | |
double precision function | cell_interface_middle (i, j, i_neighb, j_neighb, this) |
Compute middle of the interface between cell (i,j) and (i_neighb, j_neighb) More... | |
double precision function | get_lat_min (this) |
Minimal north latitude (without ghost cells) More... | |
double precision function | cell_center_lat (i, this) |
Latitude (in degrees) of the cell center (with ghost cells) More... | |
double precision function | cell_south_lat (i, this) |
Latitude (in degrees) of the cell south border (with ghost cells) More... | |
double precision function | cell_north_lat (i, this) |
Latitude (in degrees) of the cell north border (with ghost cells) More... | |
double precision function | cell_center_lon (i, j, this) |
Longitude (in degrees) of the cell center (with ghost cells) More... | |
subroutine | interior_lon_indices (j_start, j_end, this, k) |
Give the starting indices for the interior cells (i.e without the ghost cells). If we are on a ghost line, return all cells. More... | |
integer function | first_interior_lon (k, this) |
integer function | last_interior_lon (k, this) |
integer function | nb_lon_band_grid (i, this) |
Give the number of interior cells at line i inside the grid (without the ghost cells) More... | |
integer function | true_nb_lon_band_grid (this) |
Give the number of cells at line i inside the grid (with the ghost cells) More... | |
integer function | get_nb_ghosts_west (i, this) |
integer function | get_nb_ghosts_east (i, this) |
integer function | get_nb_ghosts_north (this) |
integer function | get_nb_ghosts_south (this) |
integer function | get_nb_lat_band_grid (this) |
Give the number of interior latitude inside the grid (without the ghost cells) More... | |
integer function | first_lat_band_grid (this) |
Give the first interior latitude inside the grid (without the ghost cells) More... | |
integer function | last_lat_band_grid (this) |
Give the last interior latitude inside the grid (without the ghost cells) More... | |
integer function | first_lon_band_grid (i, this) |
Give the first interior longitude at latitude i inside the grid (without the ghost cells) More... | |
integer function | last_lon_band_grid (i, this) |
Give the last interior longitude at latitude i inside the grid (without the ghost cells) More... | |
integer function | pos_lon_line_band_grid (this, i) |
subroutine | north_neighbour_cell_band_grid (j_start, j_end, i, j, this, i_neighb, offset, i_pos, j_pos, zone) |
Returns the north neighbours for cells (i,j) inside the partition. Returns the position inside the grid. Used by north_neighbour_cell_Partition. More... | |
subroutine | south_neighbour_cell_band_grid (j_start, j_end, i, j, this, i_neighb, offset, i_pos, j_pos, zone) |
Returns the south neighbours for cells (i,j) inside the partition. Returns also the ghost cells. Used by south_neighbour_cell_Partition. More... | |
subroutine | correct_indices (i_neighb, j, this, i, i_pos, j_pos) |
If j is outside the grid, takes the closest indice. More... | |
subroutine | correct_indices_other (i_neighb, j, this, i_pos, j_pos) |
Neighbour indices must be in 1 and the number of cells. More... | |
subroutine | correct_indices_ghost (j, i, this, i_neighb, i_pos, j_pos) |
If we are on north/south ghost cells and the partition is on a boundary, we must adapt the indices (as there are no east/west ghost cells on north/south ghost cells) More... | |
logical function | is_westmost_on_sector (i, j, i_pos, j_pos, this) |
Returns true if the cell is the westmost on the current sector. More... | |
logical function | is_north_south_ghost_for_boundary (i, this, i_pos, j_pos) |
Returns false if we can adjust north or south ghost cells. We cannot for a boundary partition (i.e the first on the sector). More... | |
logical function | is_first_part_sector (i_pos, j_pos) |
Returns true if the current partition is the first on the sector. More... | |
logical function | is_second_cell_sector (i, j, i_pos, j_pos, this) |
Returns true if the current cell is the second on the sector. Returns false for ghost cells. More... | |
logical function | is_beforelast_cell_sector (i, j, i_pos, j_pos, this) |
Returns true if the current cell is before the last on the sector. Returns _false for ghost cells. More... | |
logical function | is_last_part_sector (i_pos, j_pos) |
Returns true if the current partition is the last one on the sector. More... | |
logical function | is_ghost_cell (i, j, this) |
logical function | to_or_from_interior_cell (cur, i, j, i_neighb, this) |
Returns true if the cell (i,j) is a ghost cell interfacing with an interior cell, or if (i,j) is an interior cell Assume there is a neighbouring latitude. More... | |
logical function | is_interior_cell (i, j, this) |
logical function | is_north_south_ghost (i, this) |
Returns true if we are on north or south ghost cells. More... | |
logical function | is_north_ghost (i, this) |
Returns true if we are on north ghost cells. More... | |
logical function | is_south_ghost (i, this) |
Returns true if we are on south ghost cells. More... | |
logical function | is_eastmost_on_sector (i, j, i_pos, j_pos, this) |
Returns true if the cell is the westmost on the current sector. More... | |
integer function | global_cell_position_from_local (this, i, j) |
Compute the global position of the cell on the line i with float. More... | |
integer function | middle_cell (i, grid, zone) |
Returns the indice of the middle cell on the line i for the current grid. If the middle is not in the current partition, returns -1 if the middle is to the east, -2 otherwise (to the west). The middle cell is the one for the current sector. More... | |
logical function | is_indice_inside (k, this, i) |
subroutine | north_partition_reduced (j, nb_cur, nb_prev, k_start, k_end) |
subroutine | south_partition_reduced (j, nb_cur, nb_next, k_start, k_end) |
integer function | find_global_lat (lat, this) |
Find the cell whose south latitude is lat Assume global latitude is constant. More... | |
subroutine | numerical_concentration (this) |
Copies the data from the global grid (which contains numerical data). Global is an integer checking we parse all the cells of the grid. We need the previous latitude line (cannotget it analytically if dlat is not constant). Data must be allocated. global_ratio is a pointer to the global grid ratio data. More... | |
double precision function | get_lon_min (k, this) |
Minimal longitude without the ghost cells. More... | |
double precision function | cell_west_lon (i, j, this) |
Returns longitude of the west boundary for the cell at position (i,j), We suppose the caller know how to avoid the ghost cells. More... | |
double precision function | cell_west_lon2 (i, j, this) |
double precision function | cell_east_lon (i, j, this) |
Returns longitude of the east boundary for the cell at position (i,j) We suppose the caller know how to avoid the ghost cells. More... | |
double precision function | get_dlon (i, this) |
Get the variation of longitude in cells at latitude i. Constant at a given latitude. More... | |
double precision function | get_dlat (this) |
Get the variation of latitude in cells at latitude i. Constant for all cells for now. More... | |
subroutine | set_cell_ratio (val, i, j, tracer, this) |
Set the cell tracer ratio at the position (i,j), the indices being strictly positive. More... | |
subroutine | set_cell_slope (val, i, j, tracer, this) |
Set the cell slope at the position (i,j), the indices being strictly positive. More... | |
subroutine | set_cell_gradient (val, i, j, tracer, this) |
Set the cell gradient at the position (i,j), the indices being strictly positive. More... | |
subroutine | set_merid_winds (val, k, array, this) |
Set a value for meridional winds at position k in the wind array. More... | |
subroutine | set_zonal_winds (val, k, array, this) |
Set a value for zonal winds at position k in the wind array. More... | |
subroutine | set_cell_air_mass (val, i, j, this) |
Set air mass for cell (i,j), i, j > 0. More... | |
double precision function | get_cell_ratio (i, j, tracer, this) |
Get a tracer ratio for cell (i,j), i, j > 0. More... | |
double precision function | get_cell_slope (i, j, tracer, this) |
Get a tracer slope for cell (i,j), i, j > 0. More... | |
double precision function | get_cell_gradient (i, j, tracer, this) |
Get a tracer gradient for cell (i,j), i, j > 0. More... | |
double precision function | get_cell_air_mass (i, j, this) |
Get air mass for cell (i,j), i, j > 0. More... | |
double precision function | get_temperature (i, j, this) |
Get air temperature for cell (i,j), i, j > 0. More... | |
double precision function | cell_area (i, j, this) |
Area of the cell (i,j) : suppose dlat/2 << 1 for approximation Returns an area in degrees. More... | |
logical function | at_the_poles (i, this) |
logical function | is_middle_float (i, j, this, zone) |
Check if the current cell is at the middle of the sector (float computation) More... | |
integer function | local_to_global (i, j, this) |
Assuming the cell in inside a sector, return the global position. More... | |
integer function | first_global_lat (lat_min, grid) |
Assume constant dlat and round to nearest integer (the cast to int does not work in some cases) More... | |
double precision function | mass_tracer (tracer, this) |
Compute a tracer mass inside the grid (interior cells only) More... | |
subroutine | set_merid_flux (this, val, k_f, tracer) |
subroutine | set_zonal_flux (this, val, k_f, tracer) |
double precision function | get_merid_flux (this, k_f, tracer) |
double precision function | get_zonal_flux (this, k_f, tracer) |
subroutine | unset_gradients (this) |
subroutine | init_tracer_ratio (this) |
subroutine | free_band_grid (this) |
Public Attributes | |
character(*), parameter | fname_band = "band_grid_class.F90" |
Class for a generic grid defined in bands of latitude.
subroutine band_grid_class::add_east_west_ghosts | ( | type (band_grid) | this | ) |
Correct nb_lon and lon_min for east/west ghost cells.
subroutine band_grid_class::allocate_advection_data_band_grid | ( | integer, intent(in) | size_merid, |
type (band_grid) | this | ||
) |
Allocate data used in advection : slope, fluxes, mass.
subroutine band_grid_class::allocate_most_data | ( | type(band_grid) | this, |
integer, intent(in) | nb_lat | ||
) |
Allocate most of the data.
subroutine band_grid_class::allocate_mpi_borders | ( | type (band_grid) | this, |
integer | nb | ||
) |
Allocate array for mpi types for borders.
subroutine band_grid_class::allocate_winds | ( | integer, intent(in) | size_merid, |
type (band_grid) | this | ||
) |
Must compte size of meridional winds before (need partition info)
logical function band_grid_class::at_the_poles | ( | integer, intent(in) | i, |
type (band_grid) | this | ||
) |
double precision function band_grid_class::cell_area | ( | integer, intent(in) | i, |
integer, intent(in) | j, | ||
type (band_grid) | this | ||
) |
Area of the cell (i,j) : suppose dlat/2 << 1 for approximation Returns an area in degrees.
double precision function band_grid_class::cell_center_lat | ( | integer, intent(in) | i, |
type (band_grid) | this | ||
) |
Latitude (in degrees) of the cell center (with ghost cells)
double precision function band_grid_class::cell_center_lon | ( | integer, intent(in) | i, |
integer, intent(in) | j, | ||
type (band_grid) | this | ||
) |
Longitude (in degrees) of the cell center (with ghost cells)
double precision function band_grid_class::cell_east_lon | ( | integer | i, |
integer | j, | ||
type (band_grid) | this | ||
) |
Returns longitude of the east boundary for the cell at position (i,j) We suppose the caller know how to avoid the ghost cells.
double precision function band_grid_class::cell_interface_length | ( | integer, intent(in) | i, |
integer, intent(in) | j, | ||
integer, intent(in) | i_neighb, | ||
integer, intent(in) | j_neighb, | ||
type (band_grid) | this | ||
) |
Compute interface length between cell (i,j) and (i_neighb, j_neighb) in degree It is a length on the sphere, not a difference of longitudes.
double precision function band_grid_class::cell_interface_middle | ( | integer, intent(in) | i, |
integer, intent(in) | j, | ||
integer, intent(in) | i_neighb, | ||
integer, intent(in) | j_neighb, | ||
type (band_grid) | this | ||
) |
Compute middle of the interface between cell (i,j) and (i_neighb, j_neighb)
double precision function band_grid_class::cell_north_lat | ( | integer, intent(in) | i, |
type (band_grid) | this | ||
) |
Latitude (in degrees) of the cell north border (with ghost cells)
double precision function band_grid_class::cell_south_lat | ( | integer, intent(in) | i, |
type (band_grid) | this | ||
) |
Latitude (in degrees) of the cell south border (with ghost cells)
double precision function band_grid_class::cell_west_lon | ( | integer | i, |
integer | j, | ||
type (band_grid) | this | ||
) |
Returns longitude of the west boundary for the cell at position (i,j), We suppose the caller know how to avoid the ghost cells.
double precision function band_grid_class::cell_west_lon2 | ( | integer | i, |
integer | j, | ||
type (band_grid) | this | ||
) |
subroutine band_grid_class::compute_longitude_data_ghost | ( | integer | k, |
integer | first, | ||
integer | last, | ||
type (band_grid) | this, | ||
integer | width, | ||
integer | i_pos, | ||
integer | j_corr, | ||
integer, intent(in) | nb_cells, | ||
integer, intent(in) | zone, | ||
integer | k_next, | ||
logical | is_north, | ||
integer | nb_cells_next, | ||
integer, dimension(:) | list_firsts | ||
) |
Does the true computing for nb_lon and lon_min (ghost cells). Needs the neighbour first and last cells.
first,last | : neighbours inside the grid |
subroutine band_grid_class::correct_indices | ( | integer, intent(in) | i_neighb, |
integer | j, | ||
type (band_grid) | this, | ||
integer, intent(in) | i, | ||
integer, intent(in) | i_pos, | ||
integer, intent(in) | j_pos | ||
) |
If j is outside the grid, takes the closest indice.
subroutine band_grid_class::correct_indices_ghost | ( | integer | j, |
integer, intent(in) | i, | ||
type (band_grid) | this, | ||
integer, intent(in) | i_neighb, | ||
integer, intent(in) | i_pos, | ||
integer, intent(in) | j_pos | ||
) |
If we are on north/south ghost cells and the partition is on a boundary, we must adapt the indices (as there are no east/west ghost cells on north/south ghost cells)
subroutine band_grid_class::correct_indices_other | ( | integer, intent(in) | i_neighb, |
integer | j, | ||
type (band_grid) | this, | ||
integer, intent(in) | i_pos, | ||
integer, intent(in) | j_pos | ||
) |
Neighbour indices must be in 1 and the number of cells.
subroutine band_grid_class::correct_with_east_west_ghost_cells | ( | type (band_grid) | this | ) |
Correct the partition extremities by adding east and west ghost cells Does not work for north/south ghost cells ! (but used by its computation) Only works on the first sector.
subroutine band_grid_class::create_mpi_border_zwinds | ( | type (band_grid) | this, |
integer, intent(in) | direction | ||
) |
integer function band_grid_class::cur_cell_from_first | ( | integer, intent(in) | k, |
logical, intent(in) | is_first, | ||
integer, dimension(:), intent(in) | list_firsts, | ||
type (band_grid) | this | ||
) |
Returns either the first or last cell.
subroutine band_grid_class::find_east_weast_corr | ( | integer | corr_east, |
integer | corr_west, | ||
integer, intent(in) | k, | ||
type(band_grid) | this, | ||
integer, intent(in) | zone, | ||
integer, intent(in) | i, | ||
integer, intent(in) | j | ||
) |
Find the correction for minimal longitude. For non-north/south ghost cells, the test is simple. For north-south ghost cells, there are 2 cases.
k | : position in the array |
corr_east,corr_west | : longitude correction |
zone | : partition zone |
i,j | : partition position |
subroutine band_grid_class::find_first_length_cell | ( | integer | first, |
integer | length, | ||
integer, intent(in) | k, | ||
type(band_grid) | this, | ||
integer, intent(in) | nb_resized, | ||
integer, intent(in) | width, | ||
integer, intent(in) | i_pos, | ||
integer, intent(in) | j_corr, | ||
integer, intent(in) | nb_cells, | ||
integer, intent(in) | zone | ||
) |
Compute the first and last cell on latitude k. It does not take west/ghost cells into account. Does not work for leftover band Does not work for ghost cells ! (but used by its computation) Only works on the first sector.
k | : latitude on the grid |
first,length | : first cell on the latitude line, and the partition length (first is NOT an offset) |
subroutine band_grid_class::find_first_length_cell_leftover | ( | integer | first, |
integer | length, | ||
type (band_grid) | this, | ||
integer, intent(in) | width, | ||
integer, intent(in) | j_pos, | ||
integer, intent(in) | nb_cells, | ||
integer, intent(in) | n_left | ||
) |
Does not work for ghost cells ! (but used by its computation) Only works on the first sector.
j_pos | : partition position on the band (first sector) |
integer function band_grid_class::find_global_lat | ( | double precision, intent(in) | lat, |
type(band_grid) | this | ||
) |
Find the cell whose south latitude is lat Assume global latitude is constant.
subroutine band_grid_class::find_lon_min_south_ghost_last | ( | double precision, dimension(:) | lon_min, |
type (band_grid) | grid, | ||
integer, intent(in) | i, | ||
integer, intent(in) | width, | ||
double precision, intent(in) | dlon, | ||
integer, intent(in) | i_lat, | ||
integer, intent(in) | zone | ||
) |
Compute the minimal longitude for south ghost cells on the leftover band (south hemisphere). For that we compute the north neighbour of the west extremity.
i_lat : latitude from the equator |
integer function band_grid_class::first_global_lat | ( | double precision, intent(in) | lat_min, |
type(band_grid) | grid | ||
) |
Assume constant dlat and round to nearest integer (the cast to int does not work in some cases)
integer function band_grid_class::first_interior_lat | ( | type (band_grid) | this | ) |
First latitude indice for interior cells.
integer function band_grid_class::first_interior_lon | ( | integer | k, |
type (band_grid) | this | ||
) |
integer function band_grid_class::first_lat_band_grid | ( | type (band_grid) | this | ) |
Give the first interior latitude inside the grid (without the ghost cells)
integer function band_grid_class::first_lon_band_grid | ( | integer, intent(in) | i, |
type (band_grid) | this | ||
) |
Give the first interior longitude at latitude i inside the grid (without the ghost cells)
subroutine band_grid_class::free_band_grid | ( | type(band_grid) | this | ) |
double precision function band_grid_class::get_cell_air_mass | ( | integer, intent(in) | i, |
integer, intent(in) | j, | ||
type (band_grid) | this | ||
) |
Get air mass for cell (i,j), i, j > 0.
double precision function band_grid_class::get_cell_gradient | ( | integer, intent(in) | i, |
integer, intent(in) | j, | ||
integer, intent(in) | tracer, | ||
type (band_grid) | this | ||
) |
Get a tracer gradient for cell (i,j), i, j > 0.
double precision function band_grid_class::get_cell_ratio | ( | integer, intent(in) | i, |
integer, intent(in) | j, | ||
integer, intent(in) | tracer, | ||
type (band_grid) | this | ||
) |
Get a tracer ratio for cell (i,j), i, j > 0.
double precision function band_grid_class::get_cell_slope | ( | integer, intent(in) | i, |
integer, intent(in) | j, | ||
integer, intent(in) | tracer, | ||
type (band_grid) | this | ||
) |
Get a tracer slope for cell (i,j), i, j > 0.
subroutine band_grid_class::get_center_band_grid | ( | type(band_grid) | this, |
double precision | mid_lat, | ||
double precision | mid_lon | ||
) |
Find the center of the grid, that is the middle latitude and the middle longitude at this latitude.
double precision function band_grid_class::get_dlat | ( | type (band_grid) | this | ) |
Get the variation of latitude in cells at latitude i. Constant for all cells for now.
double precision function band_grid_class::get_dlon | ( | integer | i, |
type (band_grid) | this | ||
) |
Get the variation of longitude in cells at latitude i. Constant at a given latitude.
double precision function band_grid_class::get_lat_min | ( | type (band_grid) | this | ) |
Minimal north latitude (without ghost cells)
double precision function band_grid_class::get_lon_min | ( | integer | k, |
type (band_grid) | this | ||
) |
Minimal longitude without the ghost cells.
double precision function band_grid_class::get_merid_flux | ( | type(band_grid) | this, |
integer | k_f, | ||
integer | tracer | ||
) |
integer function band_grid_class::get_nb_ghosts_east | ( | integer | i, |
type (band_grid) | this | ||
) |
integer function band_grid_class::get_nb_ghosts_north | ( | type (band_grid) | this | ) |
integer function band_grid_class::get_nb_ghosts_south | ( | type (band_grid) | this | ) |
integer function band_grid_class::get_nb_ghosts_west | ( | integer | i, |
type (band_grid) | this | ||
) |
integer function band_grid_class::get_nb_lat_band_grid | ( | type (band_grid) | this | ) |
Give the number of interior latitude inside the grid (without the ghost cells)
integer function band_grid_class::get_nb_lon | ( | integer | i, |
type (band_grid) | this | ||
) |
Returns number of cells at latitude i (in the grid), without the ghost cells.
subroutine band_grid_class::get_start_border_in | ( | type (band_grid) | this, |
integer | start, | ||
character(*) | side | ||
) |
Get starting indice of the data array for the border inside the partition.
double precision function band_grid_class::get_temperature | ( | integer, intent(in) | i, |
integer, intent(in) | j, | ||
type (band_grid) | this | ||
) |
Get air temperature for cell (i,j), i, j > 0.
subroutine band_grid_class::get_temperature_pointer | ( | type(band_grid), target | this, |
double precision, dimension(:), pointer | ptr, | ||
integer, intent(in) | array | ||
) |
array | : IS_PREV, IS_CUR, IS_NEXT for correct array |
subroutine band_grid_class::get_tracer_pointer | ( | type(band_grid), target | this, |
double precision, dimension(:), pointer | ptr, | ||
integer, intent(in) | datatype, | ||
integer, intent(in) | tracer | ||
) |
subroutine band_grid_class::get_winds_pointer | ( | type(band_grid), target | this, |
double precision, dimension(:), pointer | ptr, | ||
integer, intent(in) | datatype, | ||
integer, intent(in) | array | ||
) |
datatype | : IS_ZWINDS, IS_MWINDS |
array | : IS_PREV, IS_CUR, IS_NEXT for correct array |
double precision function band_grid_class::get_zonal_flux | ( | type(band_grid) | this, |
integer | k_f, | ||
integer | tracer | ||
) |
integer function band_grid_class::global_cell_position_from_local | ( | type (band_grid) | this, |
integer, intent(in) | i, | ||
integer, intent(in) | j | ||
) |
Compute the global position of the cell on the line i with float.
integer function band_grid_class::global_to_local_cell | ( | integer, intent(in) | j_cur, |
integer, intent(in) | i, | ||
integer, intent(in) | sector, | ||
type (band_grid) | this | ||
) |
Convert global cell position to local Float operations and comparison.
j_cur | : cell in global grid (first sector) |
sector | : current sector (starts from 1) |
logical function band_grid_class::has_east_ghost_cells | ( | type (band_grid) | this | ) |
Tell us if there are ghost cells on the east boundary line.
logical function band_grid_class::has_north_ghost_cells | ( | type (band_grid) | this | ) |
Tell us if there are ghost cells on the north boundary line.
logical function band_grid_class::has_south_ghost_cells | ( | type (band_grid) | this | ) |
Tell us if there are ghost cells on the south boundary line.
logical function band_grid_class::has_west_ghost_cells | ( | type (band_grid) | this | ) |
Tell us if there are ghost cells on the west boundary line.
subroutine band_grid_class::init_advection_data | ( | type(band_grid) | this | ) |
Initializes arrays for advection.
subroutine band_grid_class::init_tracer_ratio | ( | type(band_grid) | this | ) |
subroutine band_grid_class::interior_lat_indices | ( | integer | i_start, |
integer | i_end, | ||
type (band_grid) | this | ||
) |
Give the starting indices and minimal latitude for the interior cells (i.e without the ghost cells)
subroutine band_grid_class::interior_lon_indices | ( | integer | j_start, |
integer | j_end, | ||
type (band_grid) | this, | ||
integer | k | ||
) |
Give the starting indices for the interior cells (i.e without the ghost cells). If we are on a ghost line, return all cells.
k | : latitude line |
logical function band_grid_class::is_beforelast_cell_sector | ( | integer, intent(in) | i, |
integer, intent(in) | j, | ||
integer, intent(in) | i_pos, | ||
integer, intent(in) | j_pos, | ||
type (band_grid) | this | ||
) |
Returns true if the current cell is before the last on the sector. Returns _false for ghost cells.
i,j | : cell position in the grid |
logical function band_grid_class::is_eastmost_on_sector | ( | integer, intent(in) | i, |
integer, intent(in) | j, | ||
integer, intent(in) | i_pos, | ||
integer, intent(in) | j_pos, | ||
type (band_grid) | this | ||
) |
Returns true if the cell is the westmost on the current sector.
i,j | : cell position in the partition |
i_pos,j_pos | : partition position |
logical function band_grid_class::is_first_on_band_band_grid | ( | type (band_grid) | this | ) |
Concurrent version of is_first_on_band for a grid. Uses float. TODO : optimise.
logical function band_grid_class::is_first_part_sector | ( | integer, intent(in) | i_pos, |
integer, intent(in) | j_pos | ||
) |
Returns true if the current partition is the first on the sector.
i_pos,j_pos | : partition position in the matrix array |
logical function band_grid_class::is_ghost_cell | ( | integer, intent(in) | i, |
integer, intent(in) | j, | ||
type (band_grid) | this | ||
) |
logical function band_grid_class::is_indice_inside | ( | integer, intent(in) | k, |
type(band_grid) | this, | ||
integer, intent(in) | i | ||
) |
logical function band_grid_class::is_interior_cell | ( | integer, intent(in) | i, |
integer, intent(in) | j, | ||
type (band_grid) | this | ||
) |
logical function band_grid_class::is_last_on_band_band_grid | ( | type (band_grid) | this | ) |
Concurrent version of is_last_on_band for a grid. Uses float. TODO : optimise.
logical function band_grid_class::is_last_part_sector | ( | integer, intent(in) | i_pos, |
integer, intent(in) | j_pos | ||
) |
Returns true if the current partition is the last one on the sector.
i_pos,j_pos | : partition position in the matrix array |
logical function band_grid_class::is_latitude_just_after_equator_float | ( | integer, intent(in) | i, |
type (band_grid) | this, | ||
integer, intent(in) | zone | ||
) |
Returns true if the latitude line i is just after the equator Use float !
logical function band_grid_class::is_latitude_just_before_equator_band_grid | ( | integer, intent(in) | i, |
type (band_grid) | this | ||
) |
Returns true if the latitude line i is just before the equator Concurrent of is_latitude_just_before_equator Uses float. TODO : optimize.
logical function band_grid_class::is_middle_float | ( | integer, intent(in) | i, |
integer, intent(in) | j, | ||
type (band_grid) | this, | ||
integer, intent(in) | zone | ||
) |
Check if the current cell is at the middle of the sector (float computation)
i, j : cell position in grid |
logical function band_grid_class::is_north_ghost | ( | integer, intent(in) | i, |
type (band_grid) | this | ||
) |
Returns true if we are on north ghost cells.
logical function band_grid_class::is_north_south_ghost | ( | integer, intent(in) | i, |
type (band_grid) | this | ||
) |
Returns true if we are on north or south ghost cells.
logical function band_grid_class::is_north_south_ghost_for_boundary | ( | integer, intent(in) | i, |
type (band_grid) | this, | ||
integer, intent(in) | i_pos, | ||
integer, intent(in) | j_pos | ||
) |
Returns false if we can adjust north or south ghost cells. We cannot for a boundary partition (i.e the first on the sector).
i_pos,j_pos | : partition position in the matrix array |
logical function band_grid_class::is_second_cell_sector | ( | integer, intent(in) | i, |
integer, intent(in) | j, | ||
integer, intent(in) | i_pos, | ||
integer, intent(in) | j_pos, | ||
type (band_grid) | this | ||
) |
Returns true if the current cell is the second on the sector. Returns false for ghost cells.
i,j | : cell position in the grid |
logical function band_grid_class::is_south_ghost | ( | integer, intent(in) | i, |
type (band_grid) | this | ||
) |
Returns true if we are on south ghost cells.
logical function band_grid_class::is_westmost_on_sector | ( | integer, intent(in) | i, |
integer, intent(in) | j, | ||
integer, intent(in) | i_pos, | ||
integer, intent(in) | j_pos, | ||
type (band_grid) | this | ||
) |
Returns true if the cell is the westmost on the current sector.
i,j | : cell position in the partition |
i_pos,j_pos | : partition position |
integer function band_grid_class::last_interior_lat | ( | type (band_grid) | this | ) |
Last latitude indice for interior cells.
integer function band_grid_class::last_interior_lon | ( | integer | k, |
type (band_grid) | this | ||
) |
integer function band_grid_class::last_lat_band_grid | ( | type (band_grid) | this | ) |
Give the last interior latitude inside the grid (without the ghost cells)
integer function band_grid_class::last_lon_band_grid | ( | integer, intent(in) | i, |
type (band_grid) | this | ||
) |
Give the last interior longitude at latitude i inside the grid (without the ghost cells)
integer function band_grid_class::local_to_global | ( | integer, intent(in) | i, |
integer, intent(in) | j, | ||
type (band_grid) | this | ||
) |
Assuming the cell in inside a sector, return the global position.
double precision function band_grid_class::mass_tracer | ( | integer | tracer, |
type (band_grid) | this | ||
) |
Compute a tracer mass inside the grid (interior cells only)
integer function band_grid_class::middle_cell | ( | integer, intent(in) | i, |
type (band_grid) | grid, | ||
integer, intent(in) | zone | ||
) |
Returns the indice of the middle cell on the line i for the current grid. If the middle is not in the current partition, returns -1 if the middle is to the east, -2 otherwise (to the west). The middle cell is the one for the current sector.
i | : latitude inside the grid |
zone | : partition zone |
integer function band_grid_class::nb_cells_band_grid | ( | type(band_grid) | this | ) |
Number of cells except ghost cells.
integer function band_grid_class::nb_ghost_cells_lat | ( | integer, intent(in) | k, |
type (band_grid) | this | ||
) |
Computes the total number of ghost cells laterally (east-west)
integer function band_grid_class::nb_ghosts_side | ( | logical, intent(in) | is_first, |
integer, intent(in) | k, | ||
integer, dimension(:), intent(in) | list_firsts, | ||
integer, intent(in) | i_lat, | ||
integer, intent(in) | zone, | ||
type (band_grid) | this | ||
) |
Compute the number of ghosts on a side (either east or west) for north/south.
integer function band_grid_class::nb_interfaces_lat | ( | integer, intent(in) | i, |
type (band_grid) | this | ||
) |
Number of cell interfaces between latitude line i and i+1. Only inside the current grid. Assume i and i+1 exists.
integer function band_grid_class::nb_lon_band_grid | ( | integer | i, |
type (band_grid) | this | ||
) |
Give the number of interior cells at line i inside the grid (without the ghost cells)
subroutine band_grid_class::new_band_grid | ( | type(band_grid) | this, |
integer | n_lat, | ||
double precision, dimension(:) | dlat, | ||
integer, dimension(:) | nb_lon, | ||
double precision, dimension(:) | dlon, | ||
double precision | lat_min, | ||
double precision, dimension(:) | lon_min, | ||
integer, optional | ghost_cells | ||
) |
Default constructor. Does not set number of ghosts.
Copy constructor.
subroutine band_grid_class::new_band_grid_from_pos | ( | type(band_grid) | this, |
integer, intent(in) | i_pos, | ||
integer, intent(in) | j_pos, | ||
integer, intent(in) | first_lat, | ||
integer, intent(in) | height, | ||
integer, intent(in) | width, | ||
integer | nb_resized, | ||
integer, intent(in) | zone, | ||
integer, intent(in) | ghost_cells | ||
) |
Create a grid based on its index in a (lat-lon) matrix. Grid is in most case square but the last band with a square number of partition can sometimes be rectangular hence the height and width. However, only the height varies, while the width is constant. We then need the first latitude.
zon | : partition zone |
first_lat | : starting latitude from the closest pole (beware of the minimal latitude) |
nb_resized | : number of resized partition on the current band (only one one half) |
width | : width of a square partition (before resizing) |
i_pos | : position in the partition matrix (latitude) |
j_pos | : position in the partition matrix (longitude) |
subroutine band_grid_class::new_band_grid_leftover | ( | type (band_grid) | this, |
integer, intent(in) | i_pos, | ||
integer, intent(in) | j_pos, | ||
integer | height, | ||
integer, intent(in) | width, | ||
integer | n_left, | ||
integer, intent(in) | zone, | ||
integer, intent(in) | ghost_cells | ||
) |
Create a grid at position i on the leftover band. Also finishes to set some.
height | : partition height |
width | : partition width (except for the last one) |
j | : partition number on the band (no need to correct it) |
n_left | : number of partitions on the band |
zone | : zone indice |
subroutine band_grid_class::new_indexed_mpi | ( | integer | nb, |
integer, dimension(:) | blocklens, | ||
integer, dimension(:) | offsets, | ||
integer | mpi_type | ||
) |
Creates and commit an new indexed mpi type for an array of double precision.
mpi_type | : storage in the class for the new MPI type |
subroutine band_grid_class::new_mpi_band_grid | ( | type(band_grid) | grid | ) |
MPI user-defined datatype : create the datatype from a grid The grid must have all of its array allocated !
subroutine band_grid_class::new_mpi_ghost_cells | ( | integer | val, |
integer, dimension(:), allocatable | offsets, | ||
integer, dimension(:), allocatable | blocklens, | ||
type (band_grid) | this, | ||
integer, intent(in) | nb, | ||
integer, intent(in) | direction, | ||
integer | start_prev, | ||
integer | end_prev | ||
) |
Create an MPI type for ghost cells and assing it to val.
subroutine band_grid_class::new_mpi_global_type | ( | integer | val, |
integer, dimension(:), allocatable | offsets, | ||
integer, dimension(:), allocatable | blocklens, | ||
type (band_grid) | this, | ||
integer, intent(in) | nb | ||
) |
Set the MPI type for ratio or slope field as array have the same dimensions.
val | : value to set |
integer function band_grid_class::north_extrema_cell | ( | integer, intent(in) | cur, |
logical, intent(in) | is_first, | ||
integer, dimension(:), intent(in) | list_firsts, | ||
integer, intent(in) | i_lat, | ||
integer, intent(in) | zone | ||
) |
Return either the leftmost or rightmost north neighbour.
i_lat | : current latitude |
subroutine band_grid_class::north_neighbour_cell_band_grid | ( | integer | j_start, |
integer | j_end, | ||
integer, intent(in) | i, | ||
integer, intent(in) | j, | ||
type (band_grid) | this, | ||
integer, intent(in) | i_neighb, | ||
integer, intent(in) | offset, | ||
integer, intent(in) | i_pos, | ||
integer, intent(in) | j_pos, | ||
integer, intent(in) | zone | ||
) |
Returns the north neighbours for cells (i,j) inside the partition. Returns the position inside the grid. Used by north_neighbour_cell_Partition.
i | : latitude line inside the partition |
j | : longitude on the latitude line inside the partition |
subroutine band_grid_class::north_partition_reduced | ( | integer, intent(in) | j, |
integer, intent(in) | nb_cur, | ||
integer, intent(in) | nb_prev, | ||
integer | k_start, | ||
integer | k_end | ||
) |
subroutine band_grid_class::numerical_concentration | ( | type(band_grid) | this | ) |
Copies the data from the global grid (which contains numerical data). Global is an integer checking we parse all the cells of the grid. We need the previous latitude line (cannotget it analytically if dlat is not constant). Data must be allocated. global_ratio is a pointer to the global grid ratio data.
logical function band_grid_class::offsets_sector_extrem | ( | logical, intent(in) | not_first, |
logical, intent(in) | not_last, | ||
integer, intent(in) | k, | ||
type (band_grid) | this, | ||
integer, dimension(:), intent(in) | list_firsts, | ||
integer, intent(in) | zone | ||
) |
integer function band_grid_class::pos_lon_line_band_grid | ( | type(band_grid) | this, |
integer | i | ||
) |
subroutine band_grid_class::print_band_grid | ( | type(band_grid) | this | ) |
Print grid features (for debug)
subroutine band_grid_class::print_nb_ghosts | ( | integer | i, |
type(band_grid) | this | ||
) |
Print conditionnaly the current number of east/west ghost cells at line i.
subroutine band_grid_class::print_offsets | ( | integer | i, |
type(band_grid) | this | ||
) |
Print conditionnaly the current north/south offset at line i.
subroutine band_grid_class::set_cell_air_mass | ( | double precision, intent(in) | val, |
integer, intent(in) | i, | ||
integer, intent(in) | j, | ||
type (band_grid) | this | ||
) |
Set air mass for cell (i,j), i, j > 0.
subroutine band_grid_class::set_cell_gradient | ( | double precision, intent(in) | val, |
integer, intent(in) | i, | ||
integer, intent(in) | j, | ||
integer, intent(in) | tracer, | ||
type (band_grid) | this | ||
) |
Set the cell gradient at the position (i,j), the indices being strictly positive.
subroutine band_grid_class::set_cell_ratio | ( | double precision, intent(in) | val, |
integer, intent(in) | i, | ||
integer, intent(in) | j, | ||
integer, intent(in) | tracer, | ||
type (band_grid) | this | ||
) |
Set the cell tracer ratio at the position (i,j), the indices being strictly positive.
subroutine band_grid_class::set_cell_slope | ( | double precision, intent(in) | val, |
integer, intent(in) | i, | ||
integer, intent(in) | j, | ||
integer, intent(in) | tracer, | ||
type (band_grid) | this | ||
) |
Set the cell slope at the position (i,j), the indices being strictly positive.
subroutine band_grid_class::set_data_for_neighbours | ( | type(band_grid) | this, |
integer, dimension(:), intent(in) | list_firsts, | ||
integer, intent(in) | first_lat, | ||
integer, intent(in) | zone | ||
) |
Store the position of the first north and south neighbour line, and the number of east and west ghost cells.
list_firsts | : list of the first cells on each latitude line |
first_lat | : first latitude (frome the north pole) of the partition |
subroutine band_grid_class::set_first_cells_band_grid | ( | type(band_grid) | this | ) |
Set the position of first cells on each latitude line.
subroutine band_grid_class::set_longitude_data | ( | type (band_grid) | this, |
integer, intent(in) | k, | ||
integer, intent(in) | i_pos, | ||
integer, intent(in) | j_pos, | ||
double precision, intent(in) | cur_dlon, | ||
integer, intent(in) | nb_resized, | ||
integer, intent(in) | width, | ||
integer, intent(in) | nb_cells, | ||
integer, intent(in) | zone, | ||
integer, dimension(:) | list_firsts | ||
) |
Set longitude data (nb_lon, d_lon, lon_min). But does not contains east/west ghost cells Does not work for the leftover band. Another function does that : set_north_south_ghost in partition_class Also, weast/east ghost cells are added later.
k | : position in the array |
cur_dlon | : current dlon |
zone | : partition zone |
nb_resized | : number of resized partition on the current band (only one one half) |
width | : width of a square partition (before resizing) |
i_pos,j_pos | : partition position in the array. j can be on any sector |
nb_cells | : number of cells on a sector for this latitude |
list_firsts | : store the first cell position (used later) |
subroutine band_grid_class::set_longitude_data_ghost | ( | integer, intent(in) | k, |
type (band_grid) | this, | ||
integer, intent(in) | nb_resized, | ||
integer, intent(in) | width, | ||
integer, intent(in) | i_pos, | ||
integer, intent(in) | j_corr, | ||
integer, intent(in) | nb_cells, | ||
integer, intent(in) | zone, | ||
integer, dimension(:) | list_firsts | ||
) |
Compute the lon_min and nb_lon for ghost cells. The strategy is simple, we compute the neighbours for the extremities Only works on the first sector.
k | : latitude line in the grid |
nb_cells | : number of cells on a sector for this latitude |
list_firsts | : store the first cell position (used later) |
subroutine band_grid_class::set_longitude_data_ghost_leftover | ( | integer, intent(in) | k, |
type (band_grid) | this, | ||
integer, intent(in) | width, | ||
integer, intent(in) | i_pos, | ||
integer, intent(in) | j_corr, | ||
integer, intent(in) | nb_cells, | ||
integer, intent(in) | n_left, | ||
integer, intent(in) | zone, | ||
integer, dimension(:) | list_firsts | ||
) |
Compute the lon_min and nb_lon for ghost cells on leftover band. Same strategy as "normal" bands. Only works on the first sector.
k | : latitude line in the grid |
nb_cells | : number of cells on a sector for this latitude |
subroutine band_grid_class::set_longitude_data_leftover | ( | type (band_grid) | this, |
integer, intent(in) | k, | ||
integer, intent(in) | i_pos, | ||
integer, intent(in) | j_pos, | ||
double precision, intent(in) | cur_dlon, | ||
integer, intent(in) | width, | ||
integer, intent(in) | nb_cells, | ||
integer, intent(in) | n_left, | ||
integer, intent(in) | zone, | ||
integer, dimension(:) | list_firsts | ||
) |
Set longitude data for leftover band.
subroutine band_grid_class::set_merid_flux | ( | type(band_grid) | this, |
double precision | val, | ||
integer | k_f, | ||
integer | tracer | ||
) |
subroutine band_grid_class::set_merid_winds | ( | double precision | val, |
integer, intent(in) | k, | ||
integer, intent(in) | array, | ||
type(band_grid) | this | ||
) |
Set a value for meridional winds at position k in the wind array.
array | : IS_PREV, IS_CUR, IS_NEXT for correct array |
subroutine band_grid_class::set_nb_east_west_ghosts | ( | integer, intent(in) | i_lat, |
integer, intent(in) | k, | ||
type (band_grid) | this, | ||
integer, intent(in) | zone, | ||
integer, dimension(:), intent(in) | list_firsts | ||
) |
Set the number of east and west ghost cells at line k in the partition. The number of ghost cells must be initialized to 0.
k | : line inside the partition |
subroutine band_grid_class::set_north_south_offsets | ( | integer, intent(in) | i_lat, |
integer, intent(in) | k, | ||
type (band_grid) | this, | ||
integer, intent(in) | zone, | ||
integer, dimension(:), intent(in) | list_firsts | ||
) |
Set the north and south offset.
i_lat | : latitude line from the north pole |
k | : line inside the grid |
subroutine band_grid_class::set_zonal_flux | ( | type(band_grid) | this, |
double precision | val, | ||
integer | k_f, | ||
integer | tracer | ||
) |
subroutine band_grid_class::set_zonal_winds | ( | double precision | val, |
integer, intent(in) | k, | ||
integer, intent(in) | array, | ||
type(band_grid) | this | ||
) |
Set a value for zonal winds at position k in the wind array.
array | : IS_PREV, IS_CUR, IS_NEXT for correct array |
integer function band_grid_class::size_winds_zonal | ( | type (band_grid) | this | ) |
Size of zonal winds.
integer function band_grid_class::south_extrema_cell | ( | integer, intent(in) | cur, |
logical, intent(in) | is_first, | ||
integer, dimension(:), intent(in) | list_firsts, | ||
integer, intent(in) | i_lat, | ||
integer, intent(in) | zone | ||
) |
Return either the leftmost or rightmost north neighbour.
subroutine band_grid_class::south_neighbour_cell_band_grid | ( | integer | j_start, |
integer | j_end, | ||
integer, intent(in) | i, | ||
integer, intent(in) | j, | ||
type (band_grid) | this, | ||
integer, intent(in) | i_neighb, | ||
integer, intent(in) | offset, | ||
integer, intent(in) | i_pos, | ||
integer, intent(in) | j_pos, | ||
integer, intent(in) | zone | ||
) |
Returns the south neighbours for cells (i,j) inside the partition. Returns also the ghost cells. Used by south_neighbour_cell_Partition.
i_pos,j_pos | : partition position in the matrix |
i_neighb | : neighbour latitude line (not always i+1, as we can be on the south hemisphere) |
subroutine band_grid_class::south_partition_reduced | ( | integer, intent(in) | j, |
integer, intent(in) | nb_cur, | ||
integer, intent(in) | nb_next, | ||
integer | k_start, | ||
integer | k_end | ||
) |
subroutine band_grid_class::start_indice_ghost | ( | integer | start, |
type (band_grid) | this, | ||
character(*) | side | ||
) |
Update/search starting indice of the data array for the ghost cells, with the following disposition. North/south are converted it from the local position (on a band) to a global (for the grid), while east/west are created
integer function band_grid_class::sum_nb_interfaces_merid | ( | integer, intent(in) | i, |
type (band_grid) | this | ||
) |
Sums the number of meridional interfaces up to line i (south latitude)
integer function band_grid_class::sum_nb_interfaces_zonal | ( | integer, intent(in) | i, |
type (band_grid) | this | ||
) |
Sums the number of zonal interfaces up to line i (south latitude)
logical function band_grid_class::to_or_from_interior_cell | ( | integer, intent(in) | cur, |
integer, intent(in) | i, | ||
integer, intent(in) | j, | ||
integer, intent(in) | i_neighb, | ||
type (band_grid) | this | ||
) |
Returns true if the cell (i,j) is a ghost cell interfacing with an interior cell, or if (i,j) is an interior cell Assume there is a neighbouring latitude.
cur | : cell position at i_neighb |
integer function band_grid_class::true_nb_lon_band_grid | ( | type (band_grid) | this | ) |
Give the number of cells at line i inside the grid (with the ghost cells)
subroutine band_grid_class::unset_gradients | ( | type(band_grid) | this | ) |
subroutine band_grid_class::write_id_band_grid | ( | type(band_grid) | this, |
integer, intent(in) | id_part, | ||
integer, intent(in) | id_file | ||
) |
Write the partition id and the cell corners.
id_part | : partition id |
subroutine band_grid_class::write_ratio_band_grid | ( | type(band_grid) | this, |
integer, intent(in) | tracer, | ||
integer, intent(in) | id_file | ||
) |
Write the cell ratio and coordinates of the cell center dlat is supposed constant Writing is done from the north.
character(*), parameter band_grid_class::fname_band = "band_grid_class.F90" |