This module is based on the fields_local module (which is in term a re-implementation of the fields implicit functionality). fields_gf_local implements the implicit field calculation using a gf_lo based data decomposition to improve performance on large process counts. It utilises a smaller set of processes to reduce communication costs (at the expense of significant load imbalance in the fields calculations). It assumes gf_lo velocity space integration and stores field data in gf_lo format (i.e. ik,it are decomposed over processes, everthing else is local).
Implemented by Adrian Jackson, EPCC, The University of Edinburgh; 2014-2016. Funded by an EPSRC ARCHER eCSE project. Module level routines User set parameters
Made public for testing
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
logical, | private, | parameter | :: | debug | = | .false. |
Do we do debug stuff? |
logical, | private | :: | initialised | = | .false. |
Have we initialised yet? |
|
logical, | private | :: | reinit | = | .false. |
Are we reinitialising? |
|
logical, | public | :: | dump_response | = | .false. |
Do we dump the response matrix? |
|
logical, | public | :: | read_response | = | .false. |
Do we read the response matrix from dump? |
|
logical, | public | :: | field_local_allreduce_sub | = | .false. |
If true and field_local_allreduce true then do two sub comm all reduces rather than 1 |
|
type(pc_gf_type), | private, | save | :: | pc |
This is the parallel control object |
||
type(fieldmat_gf_type), | public, | save | :: | fieldmat |
This is the top level field matrix object |
This is the top level object, consisting of a collection of ky blocks.
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
type(ky_type), | public, | dimension(:), allocatable | :: | kyb |
The ky blocks |
||
integer, | public | :: | naky |
Number of ky blocks |
|||
integer, | public | :: | ntheta0 | ||||
integer, | public | :: | npts |
Total number of theta grid points on all extended domains |
|||
integer, | public | :: | nbound |
Number of ignored boundary points |
|||
logical, | public | :: | is_local |
Does this supercell have any data on this proc? |
|||
logical, | public | :: | is_empty |
Have we got any data for the fieldmat on this proc? |
|||
logical, | public | :: | is_all_local |
Is all of this supercells data on this proc? |
|||
type(comm_type), | public | :: | fm_sub_all |
Sub communicator involving all processors with fieldmat |
|||
type(comm_type), | public | :: | fm_sub_headsc_p0 |
Sub communicator involving the supercell heads and proc0 |
|||
integer, | public | :: | prepare_type | = | 0 |
What sort of preparation do we do to the field matrix |
|
integer, | public, | dimension(:,:), allocatable | :: | heads |
An array that holds the iproc (from comm world) of the supercell head for a given ik,it point |
||
logical, | public | :: | no_populate | = | .false. |
Advanced usage only, if true then don't populate response |
|
logical, | public | :: | no_prepare | = | .false. |
Advanced usage only, if true then don't prepare (invert) response |
|
integer, | public | :: | nfield |
Number of fields we're working with |
procedure , public , :: deallocate => fm_deallocate Subroutine | |
procedure , public , :: allocate => fm_allocate Subroutine | |
procedure , public , :: debug_print => fm_debug_print Subroutine | |
procedure , public , :: get_field_update => fm_get_field_update Subroutine | |
procedure , public , :: init => fm_init Subroutine | |
procedure , public , :: reset => fm_reset Subroutine | |
procedure , public , :: populate => fm_populate Subroutine | |
procedure , public , :: init_next_field_points => fm_init_next_field_points Subroutine | |
procedure , public , :: prepare => fm_prepare Subroutine | |
procedure , public , :: make_subcom_1 => fm_make_subcom_1 Subroutine | |
procedure , public , :: make_subcom_2 => fm_make_subcom_2 Subroutine | |
procedure , public , :: calc_sc_heads => fm_calc_sc_heads Subroutine | |
procedure , public , :: gather_fields => fm_gather_fields Subroutine | |
procedure , public , :: unpack_to_field => fm_unpack_to_field Subroutine | |
procedure , public , :: write_debug_data => fm_write_debug_data Subroutine | |
procedure , public , :: set_is_local => fm_set_is_local Subroutine | |
procedure , public , :: count_subcom => fm_count_subcom Subroutine | |
procedure , public , :: reduce_an => fm_reduce_an_head_send_broadcast Subroutine | |
procedure , public , :: check_an => fm_check_an Subroutine | |
procedure , public , :: getfieldeq1_nogath => fm_getfieldeq1_nogath Subroutine | |
procedure , public , :: update_fields => fm_update_fields Subroutine | |
procedure , public , :: get_condition_numbers => fm_get_condition_numbers Subroutine | |
procedure , public , :: gather_fields_gf => fm_gather_fields_gf Subroutine | |
procedure , public , :: gather_init_fields => fm_gather_init_fields Subroutine | |
procedure , public , :: getfieldeq_nogath => fm_getfieldeq_nogath Subroutine | |
procedure , public , nopass :: update_fields_gf => fm_update_fields_gf Subroutine | |
procedure , public , :: update_fields_newstep => fm_update_fields_newstep Subroutine |
This is the type which controls/organises the parallel data decomposition etc.
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
integer, | public, | dimension(:,:), allocatable | :: | is_local |
Is the given it,ik on this proc? |
||
type(comm_type), | public, | dimension(:,:), allocatable | :: | itik_subcom |
Cell level sub-communicators |
||
integer, | public, | dimension(:,:), allocatable | :: | nproc_per_cell |
How many procs have this cell |
||
integer, | public, | dimension(:,:), allocatable | :: | nresp_per_cell |
How many rows are we personally responsible for in each cell |
||
integer, | public, | dimension(:,:), allocatable | :: | navail_per_cell |
How many rows could we be responsible for in each cell |
||
integer, | public | :: | nresp_tot |
Total number of rows available/potentially available |
|||
integer, | public | :: | navail_tot |
Total number of rows available/potentially available |
|||
logical, | public | :: | force_local_invert | = | .false. |
If true then we force local inversion |
|
logical, | public | :: | has_to_gather | = | .true. |
If true we have to gather when calling getfieldeq, determined by decomp routine |
|
integer, | public | :: | decomp_type | = | 3 |
This sets what type of decomposition is done |
procedure , public , :: deallocate => pc_deallocate Subroutine | |
procedure , public , :: allocate => pc_allocate Subroutine | |
procedure , public , :: debug_print => pc_debug_print Subroutine | |
procedure , public , :: current_nresp => pc_current_nresp Function | |
procedure , public , :: count_avail => pc_count_avail Subroutine | |
procedure , public , :: has_ik => pc_has_ik Function | |
procedure , public , :: has_it => pc_has_it Function | |
procedure , public , :: init => pc_init Subroutine | |
procedure , public , :: reset => pc_reset Subroutine | |
procedure , public , :: decomp_all_serial_local => pc_decomp_all_serial_local Subroutine | |
procedure , public , :: decomp_own_serial_local => pc_decomp_own_serial_local Subroutine | |
procedure , public , :: decomp_owncells_serial_local => pc_decomp_owncells_serial_local Subroutine | |
procedure , public , :: decomp_owncells_simplempi => pc_decomp_owncells_simplempi Subroutine | |
procedure , public , :: find_locality => pc_find_locality Subroutine | |
procedure , public , :: make_decomp => pc_make_decomp Subroutine | |
procedure , public , :: decomp => pc_decomp Subroutine |
Returns true if GS2 was built in such a way as to allow this module to work. Currently does not work with PGI compilers. See online discussions.
Redistribute initial fields from gf_lo to g_lo. These are fields that have not been reduced to the supercell heads, just calculated on each process that owns gf_lo points
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(fieldmat_gf_type), | intent(inout) | :: | self | |||
complex, | intent(inout), | dimension(-ntgrid:,:,:) | :: | gf_phnew | ||
complex, | intent(inout), | dimension(-ntgrid:,:,:) | :: | gf_apnew | ||
complex, | intent(inout), | dimension(-ntgrid:,:,:) | :: | gf_bpnew | ||
complex, | intent(inout), | dimension(-ntgrid:,:,:) | :: | gf_ph | ||
complex, | intent(inout), | dimension(-ntgrid:,:,:) | :: | gf_ap | ||
complex, | intent(inout), | dimension(-ntgrid:,:,:) | :: | gf_bp | ||
complex, | intent(inout), | dimension(-ntgrid:,:,:) | :: | phnew | ||
complex, | intent(inout), | dimension(-ntgrid:,:,:) | :: | apnew | ||
complex, | intent(inout), | dimension(-ntgrid:,:,:) | :: | bpnew | ||
complex, | intent(inout), | dimension(-ntgrid:,:,:) | :: | ph | ||
complex, | intent(inout), | dimension(-ntgrid:,:,:) | :: | ap | ||
complex, | intent(inout), | dimension(-ntgrid:,:,:) | :: | bp |
Redistribute fields from gf_lo to g_lo. This first involves unpacking the calculated fields from tmpsum on the supercell heads to the field arrays and then sending them back to the processes that calculated them in gf_lo and the processes that need them in g_lo. DD>TAGGED:Worth looking at improving this bad memory pattern
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(fieldmat_gf_type), | intent(inout) | :: | self | |||
complex, | intent(inout), | dimension(-ntgrid:,:,:) | :: | gf_ph | ||
complex, | intent(inout), | dimension(-ntgrid:,:,:) | :: | gf_ap | ||
complex, | intent(inout), | dimension(-ntgrid:,:,:) | :: | gf_bp | ||
complex, | intent(inout), | dimension(-ntgrid:,:,:) | :: | ph | ||
complex, | intent(inout), | dimension(-ntgrid:,:,:) | :: | ap | ||
complex, | intent(inout), | dimension(-ntgrid:,:,:) | :: | bp |
Update the fields using calculated update We update both the g_lo (phi,apar,bpar) and gf_lo (gf_phi,gf_apar,gf_bpar) fields as these are needed to make sure we can progress the fields in g_lo and gf_lo both times and thereby not have to transfer the full fields for all processes every iteration only the required updates. With some refactoring it should be possible to do all the work on a single set of field arrays as they are the same size regardless of the layout being used.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
complex, | intent(inout), | dimension(-ntgrid:,:,:) | :: | phi | ||
complex, | intent(inout), | dimension(-ntgrid:,:,:) | :: | apar | ||
complex, | intent(inout), | dimension(-ntgrid:,:,:) | :: | bpar | ||
complex, | intent(inout), | dimension(-ntgrid:,:,:) | :: | gf_phi | ||
complex, | intent(inout), | dimension(-ntgrid:,:,:) | :: | gf_apar | ||
complex, | intent(inout), | dimension(-ntgrid:,:,:) | :: | gf_bpar |
Update the fields using the new fields
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(fieldmat_gf_type), | intent(in) | :: | self |
FIXME : Add documentation
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(fieldmat_gf_type), | intent(in) | :: | self | |||
complex, | intent(in), | dimension (-ntgrid:,:,:) | :: | phi | ||
complex, | intent(in), | dimension (-ntgrid:,:,:) | :: | apar | ||
complex, | intent(in), | dimension (-ntgrid:,:,:) | :: | bpar | ||
complex, | intent(inout), | dimension (-ntgrid:,:,:) | :: | fieldeq | ||
complex, | intent(inout), | dimension (-ntgrid:,:,:) | :: | fieldeqa | ||
complex, | intent(inout), | dimension (-ntgrid:,:,:) | :: | fieldeqp |
Work out which cells are local
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(pc_gf_type), | intent(inout) | :: | self |
A wrapper routine to do the decomposition
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(pc_gf_type), | intent(inout) | :: | self | |||
class(fieldmat_type), | intent(inout) | :: | fieldmat |
The process decomposition works on the gf_lo layout, so if you have a point in gf_lo that is in a supercell you are in that supercell and you own the whole cell associated with that point, otherwise you don't have anything in that supercell, i.e. each cell has a single owner.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(pc_gf_type), | intent(inout) | :: | self | |||
class(fieldmat_type), | intent(inout) | :: | fieldmat |
Routine to initialise
Reset the fields_gf_local module
Finish the fields_gf_local module
Calculate the update to the fields
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
complex, | intent(out), | dimension(-ntgrid:,:,:) | :: | phi | ||
complex, | intent(out), | dimension(-ntgrid:,:,:) | :: | apar | ||
complex, | intent(out), | dimension(-ntgrid:,:,:) | :: | bpar | ||
complex, | intent(out), | dimension(-ntgrid:,:,:) | :: | gf_phi | ||
complex, | intent(out), | dimension(-ntgrid:,:,:) | :: | gf_apar | ||
complex, | intent(out), | dimension(-ntgrid:,:,:) | :: | gf_bpar | ||
logical, | intent(in), | optional | :: | do_update_in |
Initialise the fields from the initial g, just uses the fields_implicit routine
This routine advances the solution by one full time step DD>TAGGED: Should we only this is fapar>0 as well?
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | istep | |||
logical, | intent(in) | :: | remove_zonal_flows_switch |
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
character(len=*), | intent(in), | optional | :: | suffix |
If passed then use as part of file suffix |