public inbox for fortran@gcc.gnu.org
 help / color / mirror / Atom feed
From: Mikael Morin <mikael@gcc.gnu.org>
To: gcc-patches@gcc.gnu.org, fortran@gcc.gnu.org
Cc: Richard Biener <rguenther@suse.de>
Subject: [PATCH 4/4] fortran: Use pointer arithmetic to index arrays [PR102043]
Date: Sat, 16 Apr 2022 18:56:18 +0200	[thread overview]
Message-ID: <20220416165618.236666-5-mikael@gcc.gnu.org> (raw)
In-Reply-To: <20220416165618.236666-1-mikael@gcc.gnu.org>

The code generated for array references used to be ARRAY_REF trees as
could be expected.  However, the middle-end may conclude from those
trees that the indexes used are non-negative (more precisely not below
the lower bound), which is a wrong assumption in the case of "reversed-
order" arrays.

The problematic arrays are those with a descriptor and having a negative
stride for at least one dimension.  The descriptor data points to the
first element in array order (which is not the first in memory order in
that case), and the negative stride(s) makes walking the array backwards
(towards lower memory addresses), and we can access elements with
negative index wrt data pointer.

With this change, pointer arithmetic is generated by default for array
references, unless we are in a case where negative indexes can’t happen
(array descriptor’s dim element, substrings, explicit shape,
allocatable, or assumed shape contiguous).  A new flag is added to
choose between array indexing and pointer arithmetic, and it’s set
if the context can tell array indexing is safe (descriptor dim
element, substring, temporary array), or a new method is called
to decide on whether the flag should be set for one given array
expression.

	PR fortran/102043

gcc/fortran/ChangeLog:

	* trans.h (gfc_build_array_ref): Add non_negative_offset
	argument.
	* trans.cc (gfc_build_array_ref): Ditto. Use pointer arithmetic
	if non_negative_offset is false.
	* trans-expr.cc (gfc_conv_substring): Set flag in the call to
	gfc_build_array_ref.
	* trans-array.cc (gfc_get_cfi_dim_item,
	gfc_conv_descriptor_dimension): Same.
	(build_array_ref): Decide on whether to set the flag and update
	the call.
	(gfc_conv_scalarized_array_ref): Same.  New argument tmp_array.
	(gfc_conv_tmp_array_ref): Update call to
	gfc_conv_scalarized_ref.
	(non_negative_strides_array_p): New function.

gcc/testsuite/ChangeLog:

	* gfortran.dg/array_reference_3.f90: New.
	* gfortran.dg/negative_stride_1.f90: New.
	* gfortran.dg/vector_subscript_8.f90: New.
	* gfortran.dg/vector_subscript_9.f90: New.
	* gfortran.dg/c_loc_test_22.f90: Update dump patterns.
	* gfortran.dg/finalize_10.f90: Same.

Co-Authored-By: Richard Biener <rguenther@suse.de>
---
 gcc/fortran/trans-array.cc                    |  58 +++++-
 gcc/fortran/trans-expr.cc                     |   2 +-
 gcc/fortran/trans.cc                          |  42 +++-
 gcc/fortran/trans.h                           |   4 +-
 .../gfortran.dg/array_reference_3.f90         | 195 ++++++++++++++++++
 gcc/testsuite/gfortran.dg/c_loc_test_22.f90   |   4 +-
 gcc/testsuite/gfortran.dg/finalize_10.f90     |   2 +-
 .../gfortran.dg/negative_stride_1.f90         |  25 +++
 .../gfortran.dg/vector_subscript_8.f90        |  16 ++
 .../gfortran.dg/vector_subscript_9.f90        |  21 ++
 10 files changed, 354 insertions(+), 15 deletions(-)
 create mode 100644 gcc/testsuite/gfortran.dg/array_reference_3.f90
 create mode 100644 gcc/testsuite/gfortran.dg/negative_stride_1.f90
 create mode 100644 gcc/testsuite/gfortran.dg/vector_subscript_8.f90
 create mode 100644 gcc/testsuite/gfortran.dg/vector_subscript_9.f90

diff --git a/gcc/fortran/trans-array.cc b/gcc/fortran/trans-array.cc
index 11e47c0c1ca..e4b6270ccf8 100644
--- a/gcc/fortran/trans-array.cc
+++ b/gcc/fortran/trans-array.cc
@@ -172,7 +172,7 @@ static tree
 gfc_get_cfi_dim_item (tree desc, tree idx, unsigned field_idx)
 {
   tree tmp = gfc_get_cfi_descriptor_field (desc, CFI_FIELD_DIM);
-  tmp = gfc_build_array_ref (tmp, idx, NULL);
+  tmp = gfc_build_array_ref (tmp, idx, NULL_TREE, true);
   tree field = gfc_advance_chain (TYPE_FIELDS (TREE_TYPE (tmp)), field_idx);
   gcc_assert (field != NULL_TREE);
   return fold_build3_loc (input_location, COMPONENT_REF, TREE_TYPE (field),
@@ -424,7 +424,7 @@ gfc_conv_descriptor_dimension (tree desc, tree dim)
 
   tmp = gfc_get_descriptor_dimension (desc);
 
-  return gfc_build_array_ref (tmp, dim, NULL);
+  return gfc_build_array_ref (tmp, dim, NULL_TREE, true);
 }
 
 
@@ -3664,10 +3664,51 @@ build_class_array_ref (gfc_se *se, tree base, tree index)
 }
 
 
+/* Indicates that the tree EXPR is a reference to an array that can’t
+   have any negative stride.  */
+
+static bool
+non_negative_strides_array_p (tree expr)
+{
+  if (expr == NULL_TREE)
+    return false;
+
+  tree type = TREE_TYPE (expr);
+  if (POINTER_TYPE_P (type))
+    type = TREE_TYPE (type);
+
+  if (TYPE_LANG_SPECIFIC (type))
+    {
+      gfc_array_kind array_kind = GFC_TYPE_ARRAY_AKIND (type);
+
+      if (array_kind == GFC_ARRAY_ALLOCATABLE
+	  || array_kind == GFC_ARRAY_ASSUMED_SHAPE_CONT)
+	return true;
+    }
+
+  /* An array with descriptor can have negative strides.
+     We try to be conservative and return false by default here
+     if we don’t recognize a contiguous array instead of
+     returning false if we can identify a non-contiguous one.  */
+  if (!GFC_ARRAY_TYPE_P (type))
+    return false;
+
+  /* If the array was originally a dummy with a descriptor, strides can be
+     negative.  */
+  if (DECL_P (expr)
+      && DECL_LANG_SPECIFIC (expr))
+    if (tree orig_decl = GFC_DECL_SAVED_DESCRIPTOR (expr))
+      return non_negative_strides_array_p (orig_decl);
+
+  return true;
+}
+
+
 /* Build a scalarized reference to an array.  */
 
 static void
-gfc_conv_scalarized_array_ref (gfc_se * se, gfc_array_ref * ar)
+gfc_conv_scalarized_array_ref (gfc_se * se, gfc_array_ref * ar,
+			       bool tmp_array = false)
 {
   gfc_array_info *info;
   tree decl = NULL_TREE;
@@ -3717,7 +3758,10 @@ gfc_conv_scalarized_array_ref (gfc_se * se, gfc_array_ref * ar)
 	decl = info->descriptor;
     }
 
-  se->expr = gfc_build_array_ref (base, index, decl);
+  bool non_negative_stride = tmp_array
+			     || non_negative_strides_array_p (info->descriptor);
+  se->expr = gfc_build_array_ref (base, index, decl,
+				  non_negative_stride);
 }
 
 
@@ -3727,7 +3771,7 @@ void
 gfc_conv_tmp_array_ref (gfc_se * se)
 {
   se->string_length = se->ss->info->string_length;
-  gfc_conv_scalarized_array_ref (se, NULL);
+  gfc_conv_scalarized_array_ref (se, NULL, true);
   gfc_advance_se_ss_chain (se);
 }
 
@@ -3779,7 +3823,9 @@ build_array_ref (tree desc, tree offset, tree decl, tree vptr)
 
   tmp = gfc_conv_array_data (desc);
   tmp = build_fold_indirect_ref_loc (input_location, tmp);
-  tmp = gfc_build_array_ref (tmp, offset, decl, vptr);
+  tmp = gfc_build_array_ref (tmp, offset, decl,
+			     non_negative_strides_array_p (desc),
+			     vptr);
   return tmp;
 }
 
diff --git a/gcc/fortran/trans-expr.cc b/gcc/fortran/trans-expr.cc
index 3962b6848ce..7d8e7a7ffa4 100644
--- a/gcc/fortran/trans-expr.cc
+++ b/gcc/fortran/trans-expr.cc
@@ -2612,7 +2612,7 @@ gfc_conv_substring (gfc_se * se, gfc_ref * ref, int kind,
       /* For BIND(C), a BT_CHARACTER is not an ARRAY_TYPE.  */
       if (TREE_CODE (TREE_TYPE (tmp)) == ARRAY_TYPE)
 	{
-	  tmp = gfc_build_array_ref (tmp, start.expr, NULL);
+	  tmp = gfc_build_array_ref (tmp, start.expr, NULL_TREE, true);
 	  se->expr = gfc_build_addr_expr (type, tmp);
 	}
     }
diff --git a/gcc/fortran/trans.cc b/gcc/fortran/trans.cc
index 333dfa69642..f0a5dfb50fc 100644
--- a/gcc/fortran/trans.cc
+++ b/gcc/fortran/trans.cc
@@ -446,10 +446,14 @@ gfc_build_spanned_array_ref (tree base, tree offset, tree span)
 }
 
 
-/* Build an ARRAY_REF with its natural type.  */
+/* Build an ARRAY_REF with its natural type.
+   NON_NEGATIVE_OFFSET indicates if it’s true that OFFSET can’t be negative,
+   and thus that an ARRAY_REF can safely be generated.  If it’s false, we
+   have to play it safe and use pointer arithmetic.  */
 
 tree
-gfc_build_array_ref (tree base, tree offset, tree decl, tree vptr)
+gfc_build_array_ref (tree base, tree offset, tree decl,
+		     bool non_negative_offset, tree vptr)
 {
   tree type = TREE_TYPE (base);
   tree span = NULL_TREE;
@@ -495,10 +499,40 @@ gfc_build_array_ref (tree base, tree offset, tree decl, tree vptr)
      pointer arithmetic.  */
   if (span != NULL_TREE)
     return gfc_build_spanned_array_ref (base, offset, span);
-  /* Otherwise use a straightforward array reference.  */
-  else
+  /* Else use a straightforward array reference if possible.  */
+  else if (non_negative_offset)
     return build4_loc (input_location, ARRAY_REF, type, base, offset,
 		       NULL_TREE, NULL_TREE);
+  /* Otherwise use pointer arithmetic.  */
+  else
+    {
+      gcc_assert (TREE_CODE (TREE_TYPE (base)) == ARRAY_TYPE);
+      tree min = NULL_TREE;
+      if (TYPE_DOMAIN (TREE_TYPE (base))
+	  && !integer_zerop (TYPE_MIN_VALUE (TYPE_DOMAIN (TREE_TYPE (base)))))
+	min = TYPE_MIN_VALUE (TYPE_DOMAIN (TREE_TYPE (base)));
+
+      tree zero_based_index
+	   = min ? fold_build2_loc (input_location, MINUS_EXPR,
+				    gfc_array_index_type,
+				    fold_convert (gfc_array_index_type, offset),
+				    fold_convert (gfc_array_index_type, min))
+		 : fold_convert (gfc_array_index_type, offset);
+
+      tree elt_size = fold_convert (gfc_array_index_type,
+				    TYPE_SIZE_UNIT (type));
+
+      tree offset_bytes = fold_build2_loc (input_location, MULT_EXPR,
+					   gfc_array_index_type,
+					   zero_based_index, elt_size);
+
+      tree base_addr = gfc_build_addr_expr (pvoid_type_node, base);
+
+      tree ptr = fold_build_pointer_plus_loc (input_location, base_addr,
+					      offset_bytes);
+      return build1_loc (input_location, INDIRECT_REF, type,
+			 fold_convert (build_pointer_type (type), ptr));
+    }
 }
 
 
diff --git a/gcc/fortran/trans.h b/gcc/fortran/trans.h
index 738c7487a56..623aceed520 100644
--- a/gcc/fortran/trans.h
+++ b/gcc/fortran/trans.h
@@ -619,7 +619,9 @@ tree gfc_get_extern_function_decl (gfc_symbol *,
 tree gfc_build_addr_expr (tree, tree);
 
 /* Build an ARRAY_REF.  */
-tree gfc_build_array_ref (tree, tree, tree, tree vptr = NULL_TREE);
+tree gfc_build_array_ref (tree, tree, tree,
+			  bool non_negative_offset = false,
+			  tree vptr = NULL_TREE);
 
 /* Build an array ref using pointer arithmetic.  */
 tree gfc_build_spanned_array_ref (tree base, tree offset, tree span);
diff --git a/gcc/testsuite/gfortran.dg/array_reference_3.f90 b/gcc/testsuite/gfortran.dg/array_reference_3.f90
new file mode 100644
index 00000000000..85fa3317d98
--- /dev/null
+++ b/gcc/testsuite/gfortran.dg/array_reference_3.f90
@@ -0,0 +1,195 @@
+! { dg-do run }
+! { dg-additional-options "-fdump-tree-original" }
+!
+! PR fortran/102043
+! Array indexing was causing the middle-end to conclude the index
+! to be non-negative, which can be wrong for arrays with a "reversed-order"
+! descriptor.  This was fixed by using pointer arithmetic when
+! the index can be negative.
+! 
+! This test checks the code generated for array references of various kinds
+! of arrays, using either array indexing or pointer arithmetic.
+
+program p
+  implicit none
+  call check_assumed_shape_elem
+  call check_assumed_shape_scalarized
+  call check_descriptor_dim
+  call check_cfi_dim
+  call check_substring
+  call check_ptr_elem
+  call check_ptr_scalarized
+  call check_explicit_shape_elem
+  call check_explicit_shape_scalarized
+  call check_tmp_array
+  call check_allocatable_array_elem
+  call check_allocatable_array_scalarized
+contains
+  subroutine cases(assumed_shape_x)
+    integer :: assumed_shape_x(:)
+    assumed_shape_x(2) = 10
+  end subroutine cases 
+  subroutine check_assumed_shape_elem
+    integer :: x(3)
+    x = 0
+    call cases(x)
+    if (any(x /= (/ 0, 10, 0 /))) stop 10
+    ! Assumed shape array are referenced with pointer arithmetic.
+    ! { dg-final { scan-tree-dump-times "\\*\\(\\(integer\\(kind=4\\) \\*\\) assumed_shape_x.\\d+ \\+ \\(sizetype\\) \\(\\(stride.\\d+ \\* 2 \\+ offset.\\d+\\) \\* 4\\)\\) = 10;" 1 "original" } }
+  end subroutine check_assumed_shape_elem
+  subroutine casss(assumed_shape_y)
+    integer :: assumed_shape_y(:)
+    assumed_shape_y = 11
+  end subroutine casss 
+  subroutine check_assumed_shape_scalarized
+    integer :: y(3)
+    call casss(y)
+    if (any(y /= 11)) stop 11
+    ! Assumed shape array are referenced with pointer arithmetic.
+    ! { dg-final { scan-tree-dump-times "\\*\\(\\(integer\\(kind=4\\) \\*\\) assumed_shape_y.\\d+ \\+ \\(sizetype\\) \\(\\(S.\\d+ \\* D.\\d+ \\+ D.\\d+\\) \\* 4\\)\\) = 11;" 1 "original" } }
+  end subroutine check_assumed_shape_scalarized
+  subroutine check_descriptor_dim
+    integer, allocatable :: descriptor(:)
+    allocate(descriptor(4))
+    descriptor(:) = 12
+    if (any(descriptor /= 12)) stop 12
+    ! The descriptor’s dim array is referenced with array indexing.
+    ! { dg-final { scan-tree-dump-times "descriptor\\.dim\\\[0\\\]\\.ubound = 4;" 1 "original" } }
+  end subroutine check_descriptor_dim
+  subroutine ccfis(cfi_descriptor) bind(c)
+    integer :: cfi_descriptor(:)
+    cfi_descriptor = 13
+  end subroutine ccfis 
+  subroutine check_cfi_dim 
+    integer :: x(5)
+    call ccfis(x)
+    if (any(x /= 13)) stop 13
+    ! The cfi descriptor’s dim array is referenced with array indexing.
+    ! { dg-final { scan-tree-dump-times "cfi_descriptor->dim\\\[idx.\\d+\\\]\\.ubound = _cfi_descriptor->dim\\\[idx.\\d+\\\]\\.extent \\+ \\(cfi_descriptor->dim\\\[idx.\\d+\\\]\\.lbound \\+ -1\\);" 1 "original" } }
+  end subroutine check_cfi_dim
+  subroutine css(c) bind(c)
+    character :: c
+    c = 'k'
+  end subroutine css
+  subroutine check_substring
+    character(5) :: x
+    x = 'abcde'
+    call css(x(3:3))
+    if (x /= 'abkde') stop 14
+    ! Substrings use array indexing
+    ! { dg-final { scan-tree-dump-times "css \\(\\(character\\(kind=1\\)\\\[\\d+:\\d+\\\] \\*\\) &x\\\[3\\\].lb: \\d+ sz: \\d+.\\);" 1 "original" } }
+  end subroutine check_substring
+  subroutine check_ptr_elem
+    integer, target :: x(7)
+    integer, pointer :: ptr_x(:)
+    x = 0
+    ptr_x => x
+    ptr_x(4) = 16
+    if (any(ptr_x /= (/ 0, 0, 0, 16, 0, 0, 0 /))) stop 16
+    ! pointers are referenced with pointer arithmetic.
+    ! { dg-final { scan-tree-dump-times "\\*\\(integer\\(kind=4\\) \\*\\) \\(ptr_x\\.data \\+ \\(sizetype\\) \\(\\(ptr_x\\.offset \\+ ptr_x\\.dim\\\[0\\\]\\.stride \\* 4\\) \\* ptr_x\\.span\\)\\) = 16;" 1 "original" } }
+  end subroutine check_ptr_elem
+  subroutine check_ptr_scalarized
+    integer, target :: y(8)
+    integer, pointer :: ptr_y(:)
+    y = 0
+    ptr_y => y
+    ptr_y = 17
+    if (any(ptr_y /= 17)) stop 17
+    ! pointers are referenced with pointer arithmetic.
+    ! { dg-final { scan-tree-dump-times "\\*\\(\\(integer\\(kind=4\\) \\*\\) D.\\d+ \\+ \\(sizetype\\) \\(\\(S.\\d+ \\* D.\\d+ \\+ D.\\d+\\) \\* ptr_y\\.span\\)\\) = 17;" 1 "original" } }
+  end subroutine check_ptr_scalarized
+  subroutine check_explicit_shape_elem
+    integer :: explicit_shape_x(9)
+    explicit_shape_x = 0
+    explicit_shape_x(5) = 18
+    if (any(explicit_shape_x /= (/ 0, 0, 0, 0, 18, 0, 0, 0, 0 /))) stop 18
+    ! Explicit shape arrays are referenced with array indexing.
+    ! { dg-final { scan-tree-dump-times "explicit_shape_x\\\[4\\\] = 18;" 1 "original" } }
+  end subroutine check_explicit_shape_elem
+  subroutine check_explicit_shape_scalarized
+    integer :: explicit_shape_y(3)
+    explicit_shape_y = 19
+    if (any(explicit_shape_y /= 19)) stop 19
+    ! Explicit shape arrays are referenced with array indexing.
+    ! { dg-final { scan-tree-dump-times "explicit_shape_y\\\[S.\\d+ \\+ -1\\\] = 19;" 1 "original" } }
+  end subroutine check_explicit_shape_scalarized
+  subroutine check_tmp_array
+    integer :: non_tmp(6)
+    non_tmp = 15
+    non_tmp(2:5) = non_tmp(1:4) + non_tmp(3:6)
+    if (any(non_tmp /= (/ 15, 30, 30, 30, 30, 15 /))) stop 15
+    ! temporary arrays use array indexing
+    ! { dg-final { scan-tree-dump-times "\\(*\\(integer\\(kind=4\\)\\\[4\\\] \\* restrict\\) atmp.\\d+\\.data\\)\\\[S.\\d+\\\] = non_tmp\\\[S.\\d+\\\] \\+ non_tmp\\\[S.\\d+ \\+ 2\\\];" 1 "original" } }
+    ! { dg-final { scan-tree-dump-times "non_tmp\\\[S.\\d+ \\+ 1\\\] = \\(\\*\\(integer\\(kind=4\\)\\\[4\\\] \\* restrict\\) atmp.\\d+\\.data\\)\\\[S.\\d+\\\];" 1 "original" } }
+  end subroutine check_tmp_array
+  subroutine check_allocatable_array_elem
+    integer, allocatable :: allocatable_x(:)
+    allocate(allocatable_x(4),source=0)
+    allocatable_x(2) = 20
+    if (any(allocatable_x /= (/ 0, 20, 0, 0 /))) stop 20
+    ! Allocatable arrays are referenced with array indexing.
+    ! { dg-final { scan-tree-dump-times "\\(\\*\\(integer\\(kind=4\\)\\\[0:\\\] \\* restrict\\) allocatable_x\\.data\\)\\\[allocatable_x\\.offset \\+ 2\\\] = 20;" 1 "original" } }
+  end subroutine check_allocatable_array_elem
+  subroutine check_allocatable_array_scalarized
+    integer, allocatable :: allocatable_y(:)
+    allocate(allocatable_y(5),source=0)
+    allocatable_y = 21
+    if (any(allocatable_y /= 21)) stop 21
+    ! Allocatable arrays are referenced with array indexing.
+    ! { dg-final { scan-tree-dump-times "\\(\\*D.\\d+\\)\\\[S.\\d+ \\+ \\D.\\d+\\\] = 21;" 1 "original" } }
+  end subroutine check_allocatable_array_scalarized
+  subroutine cares(assumed_rank_x)
+    integer :: assumed_rank_x(..)
+    select rank(rank_1_var_x => assumed_rank_x)
+      rank(1)
+        rank_1_var_x(3) = 22
+    end select
+  end subroutine cares 
+  subroutine check_assumed_rank_elem
+    integer :: x(6)
+    x = 0
+    call cares(x)
+    if (any(x /= (/ 0, 0, 22, 0, 0, 0 /))) stop 22
+    ! Assumed rank arrays are referenced with pointer arithmetic.
+    ! { dg-final { scan-tree-dump-times "\\*\\(\\(integer\\(kind=4\\) \\*\\) __tmp_INTEGER_4_rank_1\\.data \\+ \\(sizetype\\) \\(\\(__tmp_INTEGER_4_rank_1\\.offset \\+ __tmp_INTEGER_4_rank_1\\.dim\\\[0\\\]\\.stride \\* 3\\) \\* 4\\)\\) = 22;" 1 "original" } }
+  end subroutine check_assumed_rank_elem
+  subroutine carss(assumed_rank_y)
+    integer :: assumed_rank_y(..)
+    select rank(rank_1_var_y => assumed_rank_y)
+      rank(1)
+        rank_1_var_y = 23
+    end select
+  end subroutine carss 
+  subroutine check_assumed_rank_scalarized
+    integer :: y(7)
+    call carss(y)
+    if (any(y /= 23)) stop 23
+    ! Assumed rank arrays are referenced with pointer arithmetic.
+    ! { dg-final { scan-tree-dump-times "\\*\\(\\(integer\\(kind=4\\) \\*\\) D.\\d+ \\+ \\(sizetype\\) \\(\\(S.\\d+ \\* D.\\d+ \\+ D.\\d+\\) \\* 4\\)\\) = 23;" 1 "original" } }
+  end subroutine check_assumed_rank_scalarized
+  subroutine casces(assumed_shape_cont_x)
+    integer, dimension(:), contiguous :: assumed_shape_cont_x
+    assumed_shape_cont_x(4) = 24
+  end subroutine casces 
+  subroutine check_assumed_shape_cont_elem
+    integer :: x(8)
+    x = 0
+    call casces(x)
+    if (any(x /= (/ 0, 0, 0, 24, 0, 0, 0, 0 /))) stop 24
+    ! Contiguous assumed shape arrays are referenced with array indexing.
+    ! { dg-final { scan-tree-dump-times "\\(\\*assumed_shape_cont_x.\\d+\\)\\\[stride.\\d+ \\* 4 \\+ offset.\\d+\\\] = 24;" 1 "original" } }
+  end subroutine check_assumed_shape_cont_elem
+  subroutine cascss(assumed_shape_cont_y)
+    integer, dimension(:), contiguous :: assumed_shape_cont_y
+    assumed_shape_cont_y = 25
+  end subroutine cascss 
+  subroutine check_assumed_shape_cont_scalarized
+    integer :: y(9)
+    call cascss(y)
+    if (any(y /= 25)) stop 25
+    ! Contiguous assumed shape arrays are referenced with array indexing.
+    ! { dg-final { scan-tree-dump-times "\\(\\*assumed_shape_cont_y.\\d+\\)\\\[S.\\d+ \\* D.\\d+ \\+ D.\\d+\\\] = 25;" 1 "original" } }
+  end subroutine check_assumed_shape_cont_scalarized
+end program p
+
diff --git a/gcc/testsuite/gfortran.dg/c_loc_test_22.f90 b/gcc/testsuite/gfortran.dg/c_loc_test_22.f90
index 9c40b26d830..7b1149aaa45 100644
--- a/gcc/testsuite/gfortran.dg/c_loc_test_22.f90
+++ b/gcc/testsuite/gfortran.dg/c_loc_test_22.f90
@@ -17,7 +17,7 @@ end
 ! { dg-final { scan-tree-dump-not " _gfortran_internal_pack" "original" } }
 ! { dg-final { scan-tree-dump-times "parm.\[0-9\]+.data = \\(void .\\) &\\(.xxx.\[0-9\]+\\)\\\[0\\\];" 1 "original" } }
 ! { dg-final { scan-tree-dump-times "parm.\[0-9\]+.data = \\(void .\\) &\\(.xxx.\[0-9\]+\\)\\\[D.\[0-9\]+ \\* 4\\\];" 1 "original" } }
-! { dg-final { scan-tree-dump-times "parm.\[0-9\]+.data = \\(void .\\) &\\(.yyy.\[0-9\]+\\)\\\[0\\\];" 1 "original" } }
-! { dg-final { scan-tree-dump-times "parm.\[0-9\]+.data = \\(void .\\) &\\(.yyy.\[0-9\]+\\)\\\[D.\[0-9\]+ \\* 4\\\];" 1 "original" } }
+! { dg-final { scan-tree-dump-times "parm.\[0-9\]+.data = \\(void .\\) yyy.\[0-9\]+;" 1 "original" } }
+! { dg-final { scan-tree-dump-times "parm.\[0-9\]+.data = \\(void .\\) yyy.\[0-9\]+ \\+ \\(sizetype\\) \\(D.\[0-9\]+ \\* 16\\);" 1 "original" } }
 
 ! { dg-final { scan-tree-dump-times "D.\[0-9\]+ = parm.\[0-9\]+.data;\[^;]+ptr\[1-4\] = D.\[0-9\]+;" 4 "original" } }
diff --git a/gcc/testsuite/gfortran.dg/finalize_10.f90 b/gcc/testsuite/gfortran.dg/finalize_10.f90
index 937dff5a167..c0e4170fd66 100644
--- a/gcc/testsuite/gfortran.dg/finalize_10.f90
+++ b/gcc/testsuite/gfortran.dg/finalize_10.f90
@@ -31,7 +31,7 @@ end subroutine foo
 ! { dg-final { scan-tree-dump-times "x->_vptr->_copy \\(" 1 "original" } }
 
 ! FINALIZE TYPE:
-! { dg-final { scan-tree-dump-times "parm.\[0-9\]+.data = \\(void \\*\\) &\\(\\*aa.\[0-9\]+\\)\\\[0\\\];" 1 "original" } }
+! { dg-final { scan-tree-dump-times "parm.\[0-9\]+.data = \\(void \\*\\) aa.\[0-9\]+;" 1 "original" } }
 ! { dg-final { scan-tree-dump-times "__final_m_T2 \\(&parm.\[0-9\]+, 0, 0\\);" 1 "original" } }
 ! { dg-final { scan-tree-dump-times "desc.\[0-9\]+.data = \\(void \\* restrict\\) bb;" 1 "original" } }
 ! { dg-final { scan-tree-dump-times "__final_m_T2 \\(&desc.\[0-9\]+, 0, 0\\);" 1 "original" } }
diff --git a/gcc/testsuite/gfortran.dg/negative_stride_1.f90 b/gcc/testsuite/gfortran.dg/negative_stride_1.f90
new file mode 100644
index 00000000000..45da92a223f
--- /dev/null
+++ b/gcc/testsuite/gfortran.dg/negative_stride_1.f90
@@ -0,0 +1,25 @@
+! { dg-do run }
+!
+! PR fortran/102043
+! The middle-end used to conclude from array indexing that the index
+! should be non-negative and thus that array accesses to reversed arrays
+! (i.e. with negative strides) only access the last element of the array,
+! as the access involves a pointer to array that is initialized to point
+! to the last element in the case of a reversed array.
+
+program main
+   implicit none
+   integer :: a(3, 3)
+   integer :: i
+   a = 0
+   call s(a(3:1:-1,:))
+   if (any(a(:,1) /= (/  7,  5,  3 /))) stop 1
+   if (any(a(:,2) /= (/ 17, 13, 11 /))) stop 2
+   if (any(a(:,3) /= (/ 29, 23, 19 /))) stop 3
+contains
+  subroutine s(b)
+    implicit none
+    integer, dimension(:,:) :: b
+    b = reshape((/ 3, 5, 7, 11, 13, 17, 19, 23, 29 /), (/ 3, 3 /)) 
+  end subroutine s
+end program main
diff --git a/gcc/testsuite/gfortran.dg/vector_subscript_8.f90 b/gcc/testsuite/gfortran.dg/vector_subscript_8.f90
new file mode 100644
index 00000000000..e90450b2f1b
--- /dev/null
+++ b/gcc/testsuite/gfortran.dg/vector_subscript_8.f90
@@ -0,0 +1,16 @@
+! { dg-do run }
+!
+! PR fortran/102043
+! The middle-end used to conclude from array indexing that the index
+! should be non-negative and thus that array accesses to reversed arrays
+! (i.e. with negative strides) only access the last element of the array,
+! as the access involves a pointer to array that is initialized to point
+! to the last element in the case of a reversed array.
+
+program main
+   integer, dimension (4) :: idx, a, b
+   a = (/ 11, 13, 17, 19 /)
+   idx = (/ 1, 2, 3, 4 /)
+   a(idx(4:1:-1)) = idx
+   if (a(1).ne.4) STOP 1
+end program main
diff --git a/gcc/testsuite/gfortran.dg/vector_subscript_9.f90 b/gcc/testsuite/gfortran.dg/vector_subscript_9.f90
new file mode 100644
index 00000000000..d70a773287d
--- /dev/null
+++ b/gcc/testsuite/gfortran.dg/vector_subscript_9.f90
@@ -0,0 +1,21 @@
+! { dg-do run }
+!
+! PR fortran/102043
+! The middle-end used to conclude from array indexing that the index
+! should be non-negative and thus that array accesses to reversed arrays
+! (i.e. with negative strides) only access the last element of the array,
+! as the access involves a pointer to array that is initialized to point
+! to the last element in the case of a reversed array.
+
+program main
+   integer, dimension (2) :: idx, a, b
+   a = (/ 3, 4 /)
+   idx = (/ 1, 2 /)
+   call check_values(a(idx(2:1:-1)), (/ 4, 3 /))
+contains
+   subroutine check_values(values, expected)
+      integer, dimension(:) :: values, expected
+      if (size(values) /= size(expected)) stop 1
+      if (any(values /= expected)) stop 2
+   end subroutine check_values
+end program main
-- 
2.35.1


  parent reply	other threads:[~2022-04-16 16:56 UTC|newest]

Thread overview: 12+ messages / expand[flat|nested]  mbox.gz  Atom feed  top
2022-04-16 16:56 [PATCH 0/4] Use pointer arithmetic for array references [PR102043] Mikael Morin
2022-04-16 16:56 ` [PATCH 1/4] fortran: Pre-evaluate string pointers. [PR102043] Mikael Morin
2022-04-16 16:56 ` [PATCH 2/4] fortran: Update index extraction code. [PR102043] Mikael Morin
2022-04-16 16:56 ` [PATCH 3/4] fortran: Generate an array temporary reference [PR102043] Mikael Morin
2022-04-16 16:56 ` Mikael Morin [this message]
2022-04-19 15:05 ` [PATCH 0/4] Use pointer arithmetic for array references [PR102043] Richard Biener
2022-04-22 11:09 ` *PING* " Mikael Morin
2022-04-22 13:59   ` Thomas Koenig
2022-04-24  1:57     ` Jerry D
2022-04-26 14:40     ` Hans-Peter Nilsson
2022-04-27  5:48       ` Thomas Koenig
2022-05-02  7:16       ` Stefan Schulze Frielinghaus

Reply instructions:

You may reply publicly to this message via plain-text email
using any one of the following methods:

* Save the following mbox file, import it into your mail client,
  and reply-to-all from there: mbox

  Avoid top-posting and favor interleaved quoting:
  https://en.wikipedia.org/wiki/Posting_style#Interleaved_style

* Reply using the --to, --cc, and --in-reply-to
  switches of git-send-email(1):

  git send-email \
    --in-reply-to=20220416165618.236666-5-mikael@gcc.gnu.org \
    --to=mikael@gcc.gnu.org \
    --cc=fortran@gcc.gnu.org \
    --cc=gcc-patches@gcc.gnu.org \
    --cc=rguenther@suse.de \
    /path/to/YOUR_REPLY

  https://kernel.org/pub/software/scm/git/docs/git-send-email.html

* If your mail client supports setting the In-Reply-To header
  via mailto: links, try the mailto: link
Be sure your reply has a Subject: header at the top and a blank line before the message body.
This is a public inbox, see mirroring instructions
for how to clone and mirror all data and code used for this inbox;
as well as URLs for read-only IMAP folder(s) and NNTP newsgroup(s).