New LFRic BuiltIns - field_min_max and copy_field routines based on intrinsic Fortran precisions
Introduction: The release of LFRic version 1.1 has seen psykal_builtin_light module replaced by sci_psykal_builtin_light_mod, which, at present, provides support for field_min_max and two other routines (pending Psyclone issue #489), but does not include the copy_field invokes, which use similar parallelised methods, based on intrinsic Fortran precisions. Field_min_max and copy_field are covered by LFRic tickets #4109 and #4024.
1/ LFRic Core version 1.0
Ref: [log:LFRic/branches/dev/richardmalone/apps1.0_b6]
1.0/ The relevant modules are located at (apps1.0_b6):
./components/science/source/algorithm/field_minmax_alg_mod.f90
./components/science/source/algorithm/copy_field_alg_mod.f90
./components/science/source/psy/psykal_builtin_light_mod.f90
1.1/ Psykal_builtin_light_mod contained two new classes of invoke subroutines for field min_max and field copy operations, optimised for parallel processing, utlilising only intrinsic iso Fortran numeric precisions. These are distinct from similar functions based on science numeric classes (r_single, r_double, r_def, r_tran, r_solver).
1.2/ These new invoke functions, used in field_min_max_alg_mod and copy_field_alg_mod, are as follows:
invoke_r32_field_min_max
invoke_r64_field_min_max
invoke_copy_field_32_64
invoke_copy_field_64_32
invoke_copy_field_32_32
invoke_copy_field_64_64
where: [32, 64] are intrinsic iso Fortran real precisions
2/ LFRic Core version 1.1
Ref: [log:LFRic/branches/dev/richardmalone/core1.1_b6N]
2.0/ The relevant modules are located at (core1.1_b6N):
./components/science/source/algorithm/field_minmax_alg_mod.f90
./components/science/source/algorithm/copy_field_alg_mod.f90
./components/science/source/psy/sci_psykal_builtin_light_mod.f90
2.1/ In Core version 1.1, field_min_max invoke routines currently reside in sci_psykal_builtin_light_mod, with note that:
! This is a PSyKAl-lite implementation of a built-in that will be implemented
! under PSYclone issue #489. See that issue for further details.
2.2/ At present time there are no equivalent entries for copy_field invoke routines in sci_psykal_builtin_light_mod.f90
3.0/ PSYclone and structural similarities between field_minmax and field_copy
The field_minmax and field_copy code share much in common even though one returns two scalars and the other a whole field: They both have very similar assumptions and stucture, use the same omp decorators, span the same loop architecture, norms and numerical intrinsics. So I would suspect they would transform in Psyclone in similar fashion to PSYclone issue #489.
PSyclone could, in principle, optimise out the copy operations when it knows that the precision doesn't change.
Just wanted to add that since the copy routines appear in the linear, it may also be convenient to have adjoint copy routines as well when this work gets done.
An example algorithm that makes use of a copy might be:
module some_alg_mod
contains
subroutine some_alg(x, y)
type(r32_field_type) :: x
type(r64_field_type) :: y
call invoke(Copy(y, x))
end subroutine some_alg
end module some_alg_mod
and the metadata for the new built-in would be:
!> field2 = [cast] field1
type, public, extends(kernel_type) :: Copy
private
type(arg_type) :: meta_args(2) = (/ &
arg_type(GH_FIELD, GH_REAL, GH_WRITE, ANY_SPACE_1), &
arg_type(GH_FIELD, GH_REAL, GH_READ, ANY_SPACE_1) &
/)
integer :: operates_on = DOF
contains
procedure, nopass :: copy_code
end type Copy
Upon investigation, I realise we already have real_to_real_X which I think can be used for all of the copy_field cases listed in this issue. Please could @TeranIvy and/or @mo-alistairp confirm? If I'm right, then lfric_core/components/science/source/algorithm/copy_field_alg_mod.f90 can be converted to an x90 and the various call invoke copy_field_xx_yy(fsrc, fdest) can be replaced with call invoke(real_to_real_x(fdest, fsrc)) (note that the source and destination are re-ordered to match the LFRic built-in convention of having the modified field first). If that works then the places where these copy_field... subroutines are called can be updated.
An update, I've patched:
--- components/science/source/algorithm/solver/sci_r_solver_field_vector_mod.x90(revision 52054)
+++ components/science/source/algorithm/solver/sci_r_solver_field_vector_mod.x90(working copy)
@@ -12,7 +12,7 @@
use constants_mod, only : i_def, r_def, r_solver
use r_solver_field_mod, only : r_solver_field_type
use field_mod, only : field_type
- use copy_field_alg_mod, only : copy_field
+ !ARPDBG use copy_field_alg_mod, only : copy_field
use function_space_mod, only : function_space_type
use function_space_collection_mod, only : function_space_collection
use log_mod, only : log_event, LOG_LEVEL_ERROR, &
@@ -317,7 +317,7 @@
! Copy a r_def field into the field vector
subroutine import_rdef_field(self, rdef_field, position)
- use copy_field_alg_mod, only : copy_field
+ !ARPDBG use copy_field_alg_mod, only : copy_field
implicit none
class(r_solver_field_vector_type ), intent(inout) :: self
@@ -332,7 +332,8 @@
end if
self%vector_set=.true.
- call copy_field(rdef_field, self%vector(position))
+ !ARPDBG call copy_field(rdef_field, self%vector(position))
+ call invoke(real_to_real_x(self%vector(position), rdef_field))
end subroutine import_rdef_field
and gungho_model compiles and runs successfully.