[flang] IEEE_REAL (#113948)

IEEE_REAL converts an integer or real argument to a real of a given
kind.
This commit is contained in:
vdonaldson 2024-10-30 09:56:42 -04:00 committed by GitHub
parent 7f498a865f
commit 8d406d882d
No known key found for this signature in database
GPG Key ID: B5690EEEBB952194
3 changed files with 451 additions and 2 deletions

View File

@ -289,6 +289,7 @@ struct IntrinsicLibrary {
template <mlir::arith::CmpFPredicate pred> template <mlir::arith::CmpFPredicate pred>
mlir::Value genIeeeQuietCompare(mlir::Type resultType, mlir::Value genIeeeQuietCompare(mlir::Type resultType,
llvm::ArrayRef<mlir::Value>); llvm::ArrayRef<mlir::Value>);
mlir::Value genIeeeReal(mlir::Type, llvm::ArrayRef<mlir::Value>);
mlir::Value genIeeeRint(mlir::Type, llvm::ArrayRef<mlir::Value>); mlir::Value genIeeeRint(mlir::Type, llvm::ArrayRef<mlir::Value>);
template <bool isFlag> template <bool isFlag>
void genIeeeSetFlagOrHaltingMode(llvm::ArrayRef<fir::ExtendedValue>); void genIeeeSetFlagOrHaltingMode(llvm::ArrayRef<fir::ExtendedValue>);

View File

@ -97,7 +97,6 @@ static bool isStaticallyPresent(const fir::ExtendedValue &exv) {
/// IEEE module procedure names not yet implemented for genModuleProcTODO. /// IEEE module procedure names not yet implemented for genModuleProcTODO.
static constexpr char ieee_get_underflow_mode[] = "ieee_get_underflow_mode"; static constexpr char ieee_get_underflow_mode[] = "ieee_get_underflow_mode";
static constexpr char ieee_real[] = "ieee_real";
static constexpr char ieee_rem[] = "ieee_rem"; static constexpr char ieee_rem[] = "ieee_rem";
static constexpr char ieee_set_underflow_mode[] = "ieee_set_underflow_mode"; static constexpr char ieee_set_underflow_mode[] = "ieee_set_underflow_mode";
@ -362,7 +361,7 @@ static constexpr IntrinsicHandler handlers[]{
{"ieee_quiet_le", &I::genIeeeQuietCompare<mlir::arith::CmpFPredicate::OLE>}, {"ieee_quiet_le", &I::genIeeeQuietCompare<mlir::arith::CmpFPredicate::OLE>},
{"ieee_quiet_lt", &I::genIeeeQuietCompare<mlir::arith::CmpFPredicate::OLT>}, {"ieee_quiet_lt", &I::genIeeeQuietCompare<mlir::arith::CmpFPredicate::OLT>},
{"ieee_quiet_ne", &I::genIeeeQuietCompare<mlir::arith::CmpFPredicate::UNE>}, {"ieee_quiet_ne", &I::genIeeeQuietCompare<mlir::arith::CmpFPredicate::UNE>},
{"ieee_real", &I::genModuleProcTODO<ieee_real>}, {"ieee_real", &I::genIeeeReal},
{"ieee_rem", &I::genModuleProcTODO<ieee_rem>}, {"ieee_rem", &I::genModuleProcTODO<ieee_rem>},
{"ieee_rint", &I::genIeeeRint}, {"ieee_rint", &I::genIeeeRint},
{"ieee_round_eq", &I::genIeeeTypeCompare<mlir::arith::CmpIPredicate::eq>}, {"ieee_round_eq", &I::genIeeeTypeCompare<mlir::arith::CmpIPredicate::eq>},
@ -4799,6 +4798,238 @@ IntrinsicLibrary::genIeeeQuietCompare(mlir::Type resultType,
return builder.create<fir::ConvertOp>(loc, resultType, res); return builder.create<fir::ConvertOp>(loc, resultType, res);
} }
// IEEE_REAL
mlir::Value IntrinsicLibrary::genIeeeReal(mlir::Type resultType,
llvm::ArrayRef<mlir::Value> args) {
// Convert integer or real argument A to a real of a specified kind.
// Round according to the current rounding mode.
// Signal IEEE_INVALID if A is an sNaN, and return a qNaN.
// Signal IEEE_UNDERFLOW for an inexact subnormal or zero result.
// Signal IEEE_OVERFLOW if A is finite and the result is infinite.
// Signal IEEE_INEXACT for an inexact result.
//
// if (type(a) == resultType) {
// // Conversion to the same type is a nop except for sNaN processing.
// result = a
// } else {
// result = r = real(a, kind(result))
// // Conversion to a larger type is exact.
// if (c_sizeof(a) >= c_sizeof(r)) {
// b = (a is integer) ? int(r, kind(a)) : real(r, kind(a))
// if (a == b || isNaN(a)) {
// // a is {-0, +0, -inf, +inf, NaN} or exact; result is r
// } else {
// // odd(r) is true if the low bit of significand(r) is 1
// // rounding mode ieee_other is an alias for mode ieee_nearest
// if (a < b) {
// if (mode == ieee_nearest && odd(r)) result = ieee_next_down(r)
// if (mode == ieee_other && odd(r)) result = ieee_next_down(r)
// if (mode == ieee_to_zero && a > 0) result = ieee_next_down(r)
// if (mode == ieee_away && a < 0) result = ieee_next_down(r)
// if (mode == ieee_down) result = ieee_next_down(r)
// } else { // a > b
// if (mode == ieee_nearest && odd(r)) result = ieee_next_up(r)
// if (mode == ieee_other && odd(r)) result = ieee_next_up(r)
// if (mode == ieee_to_zero && a < 0) result = ieee_next_up(r)
// if (mode == ieee_away && a > 0) result = ieee_next_up(r)
// if (mode == ieee_up) result = ieee_next_up(r)
// }
// }
// }
// }
assert(args.size() == 2);
mlir::Type i1Ty = builder.getI1Type();
mlir::Type f32Ty = mlir::FloatType::getF32(builder.getContext());
mlir::Value a = args[0];
mlir::Type aType = a.getType();
// If the argument is an sNaN, raise an invalid exception and return a qNaN.
// Otherwise return the argument.
auto processSnan = [&](mlir::Value x) {
fir::IfOp ifOp = builder.create<fir::IfOp>(loc, resultType,
genIsFPClass(i1Ty, x, snanTest),
/*withElseRegion=*/true);
builder.setInsertionPointToStart(&ifOp.getThenRegion().front());
genRaiseExcept(_FORTRAN_RUNTIME_IEEE_INVALID);
builder.create<fir::ResultOp>(loc, genQNan(resultType));
builder.setInsertionPointToStart(&ifOp.getElseRegion().front());
builder.create<fir::ResultOp>(loc, x);
builder.setInsertionPointAfter(ifOp);
return ifOp.getResult(0);
};
// Conversion is a nop, except that A may be an sNaN.
if (resultType == aType)
return processSnan(a);
// Can't directly convert between kind=2 and kind=3.
mlir::Value r, r1;
if ((aType.isBF16() && resultType.isF16()) ||
(aType.isF16() && resultType.isBF16())) {
a = builder.createConvert(loc, f32Ty, a);
aType = f32Ty;
}
r = builder.create<fir::ConvertOp>(loc, resultType, a);
mlir::IntegerType aIntType = mlir::dyn_cast<mlir::IntegerType>(aType);
mlir::FloatType aFloatType = mlir::dyn_cast<mlir::FloatType>(aType);
mlir::FloatType resultFloatType = mlir::dyn_cast<mlir::FloatType>(resultType);
// Conversion from a smaller type to a larger type is exact.
if ((aIntType ? aIntType.getWidth() : aFloatType.getWidth()) <
resultFloatType.getWidth())
return aIntType ? r : processSnan(r);
// A possibly inexact conversion result may need to be rounded up or down.
mlir::Value b = builder.create<fir::ConvertOp>(loc, aType, r);
mlir::Value aEqB;
if (aIntType)
aEqB = builder.create<mlir::arith::CmpIOp>(
loc, mlir::arith::CmpIPredicate::eq, a, b);
else
aEqB = builder.create<mlir::arith::CmpFOp>(
loc, mlir::arith::CmpFPredicate::UEQ, a, b);
// [a == b] a is a NaN or r is exact (a may be -0, +0, -inf, +inf) -- return r
fir::IfOp ifOp1 = builder.create<fir::IfOp>(loc, resultType, aEqB,
/*withElseRegion=*/true);
builder.setInsertionPointToStart(&ifOp1.getThenRegion().front());
builder.create<fir::ResultOp>(loc, aIntType ? r : processSnan(r));
// Code common to (a < b) and (a > b) branches.
builder.setInsertionPointToStart(&ifOp1.getElseRegion().front());
mlir::func::FuncOp getRound = fir::factory::getLlvmGetRounding(builder);
mlir::Value mode = builder.create<fir::CallOp>(loc, getRound).getResult(0);
mlir::Value aIsNegative, aIsPositive;
if (aIntType) {
mlir::Value zero = builder.createIntegerConstant(loc, aIntType, 0);
aIsNegative = builder.create<mlir::arith::CmpIOp>(
loc, mlir::arith::CmpIPredicate::slt, a, zero);
aIsPositive = builder.create<mlir::arith::CmpIOp>(
loc, mlir::arith::CmpIPredicate::sgt, a, zero);
} else {
mlir::Value zero = builder.createRealZeroConstant(loc, aFloatType);
aIsNegative = builder.create<mlir::arith::CmpFOp>(
loc, mlir::arith::CmpFPredicate::OLT, a, zero);
aIsPositive = builder.create<mlir::arith::CmpFOp>(
loc, mlir::arith::CmpFPredicate::OGT, a, zero);
}
mlir::Type resultIntType = builder.getIntegerType(resultFloatType.getWidth());
mlir::Value resultCast =
builder.create<mlir::arith::BitcastOp>(loc, resultIntType, r);
mlir::Value one = builder.createIntegerConstant(loc, resultIntType, 1);
mlir::Value rIsOdd = builder.create<fir::ConvertOp>(
loc, i1Ty, builder.create<mlir::arith::AndIOp>(loc, resultCast, one));
// Check for a rounding mode match.
auto match = [&](int m) {
return builder.create<mlir::arith::CmpIOp>(
loc, mlir::arith::CmpIPredicate::eq, mode,
builder.createIntegerConstant(loc, mode.getType(), m));
};
mlir::Value roundToNearestBit = builder.create<mlir::arith::OrIOp>(
loc,
// IEEE_OTHER is an alias for IEEE_NEAREST.
match(_FORTRAN_RUNTIME_IEEE_NEAREST), match(_FORTRAN_RUNTIME_IEEE_OTHER));
mlir::Value roundToNearest =
builder.create<mlir::arith::AndIOp>(loc, roundToNearestBit, rIsOdd);
mlir::Value roundToZeroBit = match(_FORTRAN_RUNTIME_IEEE_TO_ZERO);
mlir::Value roundAwayBit = match(_FORTRAN_RUNTIME_IEEE_AWAY);
mlir::Value roundToZero, roundAway, mustAdjust;
fir::IfOp adjustIfOp;
mlir::Value aLtB;
if (aIntType)
aLtB = builder.create<mlir::arith::CmpIOp>(
loc, mlir::arith::CmpIPredicate::slt, a, b);
else
aLtB = builder.create<mlir::arith::CmpFOp>(
loc, mlir::arith::CmpFPredicate::OLT, a, b);
mlir::Value upResult =
builder.create<mlir::arith::AddIOp>(loc, resultCast, one);
mlir::Value downResult =
builder.create<mlir::arith::SubIOp>(loc, resultCast, one);
// (a < b): r is inexact -- return r or ieee_next_down(r)
fir::IfOp ifOp2 = builder.create<fir::IfOp>(loc, resultType, aLtB,
/*withElseRegion=*/true);
builder.setInsertionPointToStart(&ifOp2.getThenRegion().front());
roundToZero =
builder.create<mlir::arith::AndIOp>(loc, roundToZeroBit, aIsPositive);
roundAway =
builder.create<mlir::arith::AndIOp>(loc, roundAwayBit, aIsNegative);
mlir::Value roundDown = match(_FORTRAN_RUNTIME_IEEE_DOWN);
mustAdjust =
builder.create<mlir::arith::OrIOp>(loc, roundToNearest, roundToZero);
mustAdjust = builder.create<mlir::arith::OrIOp>(loc, mustAdjust, roundAway);
mustAdjust = builder.create<mlir::arith::OrIOp>(loc, mustAdjust, roundDown);
adjustIfOp = builder.create<fir::IfOp>(loc, resultType, mustAdjust,
/*withElseRegion=*/true);
builder.setInsertionPointToStart(&adjustIfOp.getThenRegion().front());
if (resultType.isF80())
r1 = fir::runtime::genNearest(builder, loc, r,
builder.createBool(loc, false));
else
r1 = builder.create<mlir::arith::BitcastOp>(
loc, resultType,
builder.create<mlir::arith::SelectOp>(loc, aIsNegative, upResult,
downResult));
builder.create<fir::ResultOp>(loc, r1);
builder.setInsertionPointToStart(&adjustIfOp.getElseRegion().front());
builder.create<fir::ResultOp>(loc, r);
builder.setInsertionPointAfter(adjustIfOp);
builder.create<fir::ResultOp>(loc, adjustIfOp.getResult(0));
// (a > b): r is inexact -- return r or ieee_next_up(r)
builder.setInsertionPointToStart(&ifOp2.getElseRegion().front());
roundToZero =
builder.create<mlir::arith::AndIOp>(loc, roundToZeroBit, aIsNegative);
roundAway =
builder.create<mlir::arith::AndIOp>(loc, roundAwayBit, aIsPositive);
mlir::Value roundUp = match(_FORTRAN_RUNTIME_IEEE_UP);
mustAdjust =
builder.create<mlir::arith::OrIOp>(loc, roundToNearest, roundToZero);
mustAdjust = builder.create<mlir::arith::OrIOp>(loc, mustAdjust, roundAway);
mustAdjust = builder.create<mlir::arith::OrIOp>(loc, mustAdjust, roundUp);
adjustIfOp = builder.create<fir::IfOp>(loc, resultType, mustAdjust,
/*withElseRegion=*/true);
builder.setInsertionPointToStart(&adjustIfOp.getThenRegion().front());
if (resultType.isF80())
r1 = fir::runtime::genNearest(builder, loc, r,
builder.createBool(loc, true));
else
r1 = builder.create<mlir::arith::BitcastOp>(
loc, resultType,
builder.create<mlir::arith::SelectOp>(loc, aIsPositive, upResult,
downResult));
builder.create<fir::ResultOp>(loc, r1);
builder.setInsertionPointToStart(&adjustIfOp.getElseRegion().front());
builder.create<fir::ResultOp>(loc, r);
builder.setInsertionPointAfter(adjustIfOp);
builder.create<fir::ResultOp>(loc, adjustIfOp.getResult(0));
// Generate exceptions for (a < b) and (a > b) branches.
builder.setInsertionPointAfter(ifOp2);
r = ifOp2.getResult(0);
fir::IfOp exceptIfOp1 = builder.create<fir::IfOp>(
loc, genIsFPClass(i1Ty, r, infiniteTest), /*withElseRegion=*/true);
builder.setInsertionPointToStart(&exceptIfOp1.getThenRegion().front());
genRaiseExcept(_FORTRAN_RUNTIME_IEEE_OVERFLOW |
_FORTRAN_RUNTIME_IEEE_INEXACT);
builder.setInsertionPointToStart(&exceptIfOp1.getElseRegion().front());
fir::IfOp exceptIfOp2 = builder.create<fir::IfOp>(
loc, genIsFPClass(i1Ty, r, subnormalTest | zeroTest),
/*withElseRegion=*/true);
builder.setInsertionPointToStart(&exceptIfOp2.getThenRegion().front());
genRaiseExcept(_FORTRAN_RUNTIME_IEEE_UNDERFLOW |
_FORTRAN_RUNTIME_IEEE_INEXACT);
builder.setInsertionPointToStart(&exceptIfOp2.getElseRegion().front());
genRaiseExcept(_FORTRAN_RUNTIME_IEEE_INEXACT);
builder.setInsertionPointAfter(exceptIfOp1);
builder.create<fir::ResultOp>(loc, ifOp2.getResult(0));
builder.setInsertionPointAfter(ifOp1);
return ifOp1.getResult(0);
}
// IEEE_RINT // IEEE_RINT
mlir::Value IntrinsicLibrary::genIeeeRint(mlir::Type resultType, mlir::Value IntrinsicLibrary::genIeeeRint(mlir::Type resultType,
llvm::ArrayRef<mlir::Value> args) { llvm::ArrayRef<mlir::Value> args) {

View File

@ -0,0 +1,217 @@
! RUN: bbc -emit-hlfir -o - %s | FileCheck %s
! CHECK-LABEL: c.func @_QQmain
program p
use ieee_arithmetic, only: ieee_real
! CHECK: %[[V_0:[0-9]+]] = fir.alloca i16 {bindc_name = "j2", uniq_name = "_QFEj2"}
! CHECK: %[[V_1:[0-9]+]]:2 = hlfir.declare %[[V_0]] {uniq_name = "_QFEj2"} : (!fir.ref<i16>) -> (!fir.ref<i16>, !fir.ref<i16>)
! CHECK: %[[V_2:[0-9]+]] = fir.alloca i64 {bindc_name = "j8", uniq_name = "_QFEj8"}
! CHECK: %[[V_3:[0-9]+]]:2 = hlfir.declare %[[V_2]] {uniq_name = "_QFEj8"} : (!fir.ref<i64>) -> (!fir.ref<i64>, !fir.ref<i64>)
! CHECK: %[[V_4:[0-9]+]] = fir.alloca f16 {bindc_name = "x2", uniq_name = "_QFEx2"}
! CHECK: %[[V_5:[0-9]+]]:2 = hlfir.declare %[[V_4]] {uniq_name = "_QFEx2"} : (!fir.ref<f16>) -> (!fir.ref<f16>, !fir.ref<f16>)
! CHECK: %[[V_6:[0-9]+]] = fir.alloca f32 {bindc_name = "x4", uniq_name = "_QFEx4"}
! CHECK: %[[V_7:[0-9]+]]:2 = hlfir.declare %[[V_6]] {uniq_name = "_QFEx4"} : (!fir.ref<f32>) -> (!fir.ref<f32>, !fir.ref<f32>)
! CHECK: %[[V_8:[0-9]+]] = fir.alloca f64 {bindc_name = "x8", uniq_name = "_QFEx8"}
! CHECK: %[[V_9:[0-9]+]]:2 = hlfir.declare %[[V_8]] {uniq_name = "_QFEx8"} : (!fir.ref<f64>) -> (!fir.ref<f64>, !fir.ref<f64>)
integer(2) :: j2
integer(8) :: j8
real(2) :: x2
real(4) :: x4
real(8) :: x8
! CHECK: hlfir.assign %c-32768{{.*}} to %[[V_1]]#0 : i16, !fir.ref<i16>
j2 = -huge(j2) - 1
! CHECK: %[[V_10:[0-9]+]] = fir.load %[[V_1]]#0 : !fir.ref<i16>
! CHECK: %[[V_11:[0-9]+]] = fir.convert %[[V_10]] : (i16) -> f32
! CHECK: hlfir.assign %[[V_11]] to %[[V_7]]#0 : f32, !fir.ref<f32>
x4 = ieee_real(j2,4) ! exact
! print*, j2, ' -> ', x4
! CHECK: hlfir.assign %c33{{.*}} to %[[V_3]]#0 : i64, !fir.ref<i64>
j8 = 33
! CHECK: %[[V_12:[0-9]+]] = fir.load %[[V_3]]#0 : !fir.ref<i64>
! CHECK: %[[V_13:[0-9]+]] = fir.convert %[[V_12]] : (i64) -> f32
! CHECK: %[[V_14:[0-9]+]] = fir.convert %[[V_13]] : (f32) -> i64
! CHECK: %[[V_15:[0-9]+]] = arith.cmpi eq, %[[V_12]], %[[V_14]] : i64
! CHECK: %[[V_16:[0-9]+]] = fir.if %[[V_15]] -> (f32) {
! CHECK: fir.result %[[V_13]] : f32
! CHECK: } else {
! CHECK: %[[V_27:[0-9]+]] = fir.call @llvm.get.rounding() fastmath<contract> : () -> i32
! CHECK-DAG: %[[V_28:[0-9]+]] = arith.cmpi slt, %[[V_12]], %c0{{.*}} : i64
! CHECK-DAG: %[[V_29:[0-9]+]] = arith.cmpi sgt, %[[V_12]], %c0{{.*}} : i64
! CHECK-DAG: %[[V_30:[0-9]+]] = arith.bitcast %[[V_13]] : f32 to i32
! CHECK-DAG: %[[V_31:[0-9]+]] = arith.andi %[[V_30]], %c1{{.*}} : i32
! CHECK-DAG: %[[V_32:[0-9]+]] = fir.convert %[[V_31]] : (i32) -> i1
! CHECK-DAG: %[[V_33:[0-9]+]] = arith.cmpi eq, %[[V_27]], %c5{{.*}} : i32
! CHECK-DAG: %[[V_34:[0-9]+]] = arith.cmpi eq, %[[V_27]], %c1{{.*}} : i32
! CHECK-DAG: %[[V_35:[0-9]+]] = arith.ori %[[V_34]], %[[V_33]] : i1
! CHECK-DAG: %[[V_36:[0-9]+]] = arith.andi %[[V_35]], %[[V_32]] : i1
! CHECK-DAG: %[[V_37:[0-9]+]] = arith.cmpi eq, %[[V_27]], %c0{{.*}} : i32
! CHECK-DAG: %[[V_38:[0-9]+]] = arith.cmpi eq, %[[V_27]], %c4{{.*}} : i32
! CHECK-DAG: %[[V_39:[0-9]+]] = arith.cmpi slt, %[[V_12]], %[[V_14]] : i64
! CHECK-DAG: %[[V_40:[0-9]+]] = arith.addi %[[V_30]], %c1{{.*}} : i32
! CHECK-DAG: %[[V_41:[0-9]+]] = arith.subi %[[V_30]], %c1{{.*}} : i32
! CHECK: %[[V_42:[0-9]+]] = fir.if %[[V_39]] -> (f32) {
! CHECK-DAG: %[[V_44:[0-9]+]] = arith.andi %[[V_37]], %[[V_29]] : i1
! CHECK-DAG: %[[V_45:[0-9]+]] = arith.andi %[[V_38]], %[[V_28]] : i1
! CHECK-DAG: %[[V_46:[0-9]+]] = arith.cmpi eq, %[[V_27]], %c3{{.*}} : i32
! CHECK-DAG: %[[V_47:[0-9]+]] = arith.ori %[[V_36]], %[[V_44]] : i1
! CHECK-DAG: %[[V_48:[0-9]+]] = arith.ori %[[V_47]], %[[V_45]] : i1
! CHECK-DAG: %[[V_49:[0-9]+]] = arith.ori %[[V_48]], %[[V_46]] : i1
! CHECK: %[[V_50:[0-9]+]] = fir.if %[[V_49]] -> (f32) {
! CHECK: %[[V_51:[0-9]+]] = arith.select %[[V_28]], %[[V_40]], %[[V_41]] : i32
! CHECK: %[[V_52:[0-9]+]] = arith.bitcast %[[V_51]] : i32 to f32
! CHECK: fir.result %[[V_52]] : f32
! CHECK: } else {
! CHECK: fir.result %[[V_13]] : f32
! CHECK: }
! CHECK: fir.result %[[V_50]] : f32
! CHECK: } else {
! CHECK-DAG: %[[V_44:[0-9]+]] = arith.andi %[[V_37]], %[[V_28]] : i1
! CHECK-DAG: %[[V_45:[0-9]+]] = arith.andi %[[V_38]], %[[V_29]] : i1
! CHECK-DAG: %[[V_46:[0-9]+]] = arith.cmpi eq, %[[V_27]], %c2{{.*}} : i32
! CHECK-DAG: %[[V_47:[0-9]+]] = arith.ori %[[V_36]], %[[V_44]] : i1
! CHECK-DAG: %[[V_48:[0-9]+]] = arith.ori %[[V_47]], %[[V_45]] : i1
! CHECK-DAG: %[[V_49:[0-9]+]] = arith.ori %[[V_48]], %[[V_46]] : i1
! CHECK: %[[V_50:[0-9]+]] = fir.if %[[V_49]] -> (f32) {
! CHECK: %[[V_51:[0-9]+]] = arith.select %[[V_29]], %[[V_40]], %[[V_41]] : i32
! CHECK: %[[V_52:[0-9]+]] = arith.bitcast %[[V_51]] : i32 to f32
! CHECK: fir.result %[[V_52]] : f32
! CHECK: } else {
! CHECK: fir.result %[[V_13]] : f32
! CHECK: }
! CHECK: fir.result %[[V_50]] : f32
! CHECK: }
! CHECK: %[[V_43:[0-9]+]] = "llvm.intr.is.fpclass"(%[[V_42]]) <{bit = 516 : i32}> : (f32) -> i1
! CHECK: fir.if %[[V_43]] {
! CHECK: %[[V_44:[0-9]+]] = fir.call @_FortranAMapException(%c40{{.*}}) fastmath<contract> : (i32) -> i32
! CHECK: %[[V_45:[0-9]+]] = fir.call @feraiseexcept(%[[V_44]]) fastmath<contract> : (i32) -> i32
! CHECK: } else {
! CHECK: %[[V_44:[0-9]+]] = "llvm.intr.is.fpclass"(%[[V_42]]) <{bit = 240 : i32}> : (f32) -> i1
! CHECK: fir.if %[[V_44]] {
! CHECK: %[[V_45:[0-9]+]] = fir.call @_FortranAMapException(%c48{{.*}}) fastmath<contract> : (i32) -> i32
! CHECK: %[[V_46:[0-9]+]] = fir.call @feraiseexcept(%[[V_45]]) fastmath<contract> : (i32) -> i32
! CHECK: } else {
! CHECK: %[[V_45:[0-9]+]] = fir.call @_FortranAMapException(%c32{{.*}}) fastmath<contract> : (i32) -> i32
! CHECK: %[[V_46:[0-9]+]] = fir.call @feraiseexcept(%[[V_45]]) fastmath<contract> : (i32) -> i32
! CHECK: }
! CHECK: }
! CHECK: fir.result %[[V_42]] : f32
! CHECK: }
! CHECK: hlfir.assign %[[V_16]] to %[[V_7]]#0 : f32, !fir.ref<f32>
x4 = ieee_real(j8,4)
! print*, j8, ' -> ', x4
! CHECK: hlfir.assign %cst{{[_0-9]*}} to %[[V_5]]#0 : f16, !fir.ref<f16>
x2 = 3.33
! CHECK: %[[V_17:[0-9]+]] = fir.load %[[V_5]]#0 : !fir.ref<f16>
! CHECK: %[[V_18:[0-9]+]] = fir.convert %[[V_17]] : (f16) -> f32
! CHECK: %[[V_19:[0-9]+]] = "llvm.intr.is.fpclass"(%[[V_18]]) <{bit = 1 : i32}> : (f32) -> i1
! CHECK: %[[V_20:[0-9]+]] = fir.if %[[V_19]] -> (f32) {
! CHECK: %[[V_27:[0-9]+]] = fir.call @_FortranAMapException(%c1{{.*}}) fastmath<contract> : (i32) -> i32
! CHECK: %[[V_28:[0-9]+]] = fir.call @feraiseexcept(%[[V_27]]) fastmath<contract> : (i32) -> i32
! CHECK: %[[V_29:[0-9]+]] = fir.address_of(@_FortranAIeeeValueTable_4) : !fir.ref<!fir.array<12xi32>>
! CHECK: %[[V_30:[0-9]+]] = fir.coordinate_of %[[V_29]], %c2{{.*}} : (!fir.ref<!fir.array<12xi32>>, i8) -> !fir.ref<i32>
! CHECK: %[[V_31:[0-9]+]] = fir.load %[[V_30]] : !fir.ref<i32>
! CHECK: %[[V_32:[0-9]+]] = arith.bitcast %[[V_31]] : i32 to f32
! CHECK: fir.result %[[V_32]] : f32
! CHECK: } else {
! CHECK: fir.result %[[V_18]] : f32
! CHECK: }
! CHECK: %[[V_21:[0-9]+]] = fir.convert %[[V_20]] : (f32) -> f16
! CHECK: hlfir.assign %[[V_21]] to %[[V_5]]#0 : f16, !fir.ref<f16>
x2 = ieee_real(x2,4) ! exact
! print*, x2, ' -> ', x2
! CHECK: hlfir.assign %cst{{[_0-9]*}} to %[[V_9]]#0 : f64, !fir.ref<f64>
x8 = -0.
! CHECK: %[[V_22:[0-9]+]] = fir.load %[[V_9]]#0 : !fir.ref<f64>
! CHECK: %[[V_23:[0-9]+]] = fir.convert %[[V_22]] : (f64) -> f32
! CHECK: %[[V_24:[0-9]+]] = fir.convert %[[V_23]] : (f32) -> f64
! CHECK: %[[V_25:[0-9]+]] = arith.cmpf ueq, %[[V_22]], %[[V_24]] fastmath<contract> : f64
! CHECK: %[[V_26:[0-9]+]] = fir.if %[[V_25]] -> (f32) {
! CHECK: %[[V_27:[0-9]+]] = "llvm.intr.is.fpclass"(%[[V_23]]) <{bit = 1 : i32}> : (f32) -> i1
! CHECK: %[[V_28:[0-9]+]] = fir.if %[[V_27]] -> (f32) {
! CHECK: %[[V_29:[0-9]+]] = fir.call @_FortranAMapException(%c1{{.*}}) fastmath<contract> : (i32) -> i32
! CHECK: %[[V_30:[0-9]+]] = fir.call @feraiseexcept(%[[V_29]]) fastmath<contract> : (i32) -> i32
! CHECK: %[[V_31:[0-9]+]] = fir.address_of(@_FortranAIeeeValueTable_4) : !fir.ref<!fir.array<12xi32>>
! CHECK: %[[V_32:[0-9]+]] = fir.coordinate_of %[[V_31]], %c2{{.*}} : (!fir.ref<!fir.array<12xi32>>, i8) -> !fir.ref<i32>
! CHECK: %[[V_33:[0-9]+]] = fir.load %[[V_32]] : !fir.ref<i32>
! CHECK: %[[V_34:[0-9]+]] = arith.bitcast %[[V_33]] : i32 to f32
! CHECK: fir.result %[[V_34]] : f32
! CHECK: } else {
! CHECK: fir.result %[[V_23]] : f32
! CHECK: }
! CHECK: fir.result %[[V_28]] : f32
! CHECK: } else {
! CHECK-DAG: %[[V_27:[0-9]+]] = fir.call @llvm.get.rounding() fastmath<contract> : () -> i32
! CHECK-DAG: %[[V_28:[0-9]+]] = arith.cmpf olt, %[[V_22]], %cst{{[_0-9]*}} fastmath<contract> : f64
! CHECK-DAG: %[[V_29:[0-9]+]] = arith.cmpf ogt, %[[V_22]], %cst{{[_0-9]*}} fastmath<contract> : f64
! CHECK-DAG: %[[V_30:[0-9]+]] = arith.bitcast %[[V_23]] : f32 to i32
! CHECK-DAG: %[[V_31:[0-9]+]] = arith.andi %[[V_30]], %c1{{.*}} : i32
! CHECK-DAG: %[[V_32:[0-9]+]] = fir.convert %[[V_31]] : (i32) -> i1
! CHECK-DAG: %[[V_33:[0-9]+]] = arith.cmpi eq, %[[V_27]], %c5{{.*}} : i32
! CHECK-DAG: %[[V_34:[0-9]+]] = arith.cmpi eq, %[[V_27]], %c1{{.*}} : i32
! CHECK-DAG: %[[V_35:[0-9]+]] = arith.ori %[[V_34]], %[[V_33]] : i1
! CHECK-DAG: %[[V_36:[0-9]+]] = arith.andi %[[V_35]], %[[V_32]] : i1
! CHECK-DAG: %[[V_37:[0-9]+]] = arith.cmpi eq, %[[V_27]], %c0{{.*}} : i32
! CHECK-DAG: %[[V_38:[0-9]+]] = arith.cmpi eq, %[[V_27]], %c4{{.*}} : i32
! CHECK-DAG: %[[V_39:[0-9]+]] = arith.cmpf olt, %[[V_22]], %[[V_24]] fastmath<contract> : f64
! CHECK-DAG: %[[V_40:[0-9]+]] = arith.addi %[[V_30]], %c1{{.*}} : i32
! CHECK-DAG: %[[V_41:[0-9]+]] = arith.subi %[[V_30]], %c1{{.*}} : i32
! CHECK: %[[V_42:[0-9]+]] = fir.if %[[V_39]] -> (f32) {
! CHECK-DAG: %[[V_44:[0-9]+]] = arith.andi %[[V_37]], %[[V_29]] : i1
! CHECK-DAG: %[[V_45:[0-9]+]] = arith.andi %[[V_38]], %[[V_28]] : i1
! CHECK-DAG: %[[V_46:[0-9]+]] = arith.cmpi eq, %[[V_27]], %c3{{.*}} : i32
! CHECK-DAG: %[[V_47:[0-9]+]] = arith.ori %[[V_36]], %[[V_44]] : i1
! CHECK-DAG: %[[V_48:[0-9]+]] = arith.ori %[[V_47]], %[[V_45]] : i1
! CHECK-DAG: %[[V_49:[0-9]+]] = arith.ori %[[V_48]], %[[V_46]] : i1
! CHECK: %[[V_50:[0-9]+]] = fir.if %[[V_49]] -> (f32) {
! CHECK: %[[V_51:[0-9]+]] = arith.select %[[V_28]], %[[V_40]], %[[V_41]] : i32
! CHECK: %[[V_52:[0-9]+]] = arith.bitcast %[[V_51]] : i32 to f32
! CHECK: fir.result %[[V_52]] : f32
! CHECK: } else {
! CHECK: fir.result %[[V_23]] : f32
! CHECK: }
! CHECK: fir.result %[[V_50]] : f32
! CHECK: } else {
! CHECK-DAG: %[[V_44:[0-9]+]] = arith.andi %[[V_37]], %[[V_28]] : i1
! CHECK-DAG: %[[V_45:[0-9]+]] = arith.andi %[[V_38]], %[[V_29]] : i1
! CHECK-DAG: %[[V_46:[0-9]+]] = arith.cmpi eq, %[[V_27]], %c2{{.*}} : i32
! CHECK-DAG: %[[V_47:[0-9]+]] = arith.ori %[[V_36]], %[[V_44]] : i1
! CHECK-DAG: %[[V_48:[0-9]+]] = arith.ori %[[V_47]], %[[V_45]] : i1
! CHECK-DAG: %[[V_49:[0-9]+]] = arith.ori %[[V_48]], %[[V_46]] : i1
! CHECK: %[[V_50:[0-9]+]] = fir.if %[[V_49]] -> (f32) {
! CHECK: %[[V_51:[0-9]+]] = arith.select %[[V_29]], %[[V_40]], %[[V_41]] : i32
! CHECK: %[[V_52:[0-9]+]] = arith.bitcast %[[V_51]] : i32 to f32
! CHECK: fir.result %[[V_52]] : f32
! CHECK: } else {
! CHECK: fir.result %[[V_23]] : f32
! CHECK: }
! CHECK: fir.result %[[V_50]] : f32
! CHECK: }
! CHECK: %[[V_43:[0-9]+]] = "llvm.intr.is.fpclass"(%[[V_42]]) <{bit = 516 : i32}> : (f32) -> i1
! CHECK: fir.if %[[V_43]] {
! CHECK: %[[V_44:[0-9]+]] = fir.call @_FortranAMapException(%c40{{.*}}) fastmath<contract> : (i32) -> i32
! CHECK: %[[V_45:[0-9]+]] = fir.call @feraiseexcept(%[[V_44]]) fastmath<contract> : (i32) -> i32
! CHECK: } else {
! CHECK: %[[V_44:[0-9]+]] = "llvm.intr.is.fpclass"(%[[V_42]]) <{bit = 240 : i32}> : (f32) -> i1
! CHECK: fir.if %[[V_44]] {
! CHECK: %[[V_45:[0-9]+]] = fir.call @_FortranAMapException(%c48{{.*}}) fastmath<contract> : (i32) -> i32
! CHECK: %[[V_46:[0-9]+]] = fir.call @feraiseexcept(%[[V_45]]) fastmath<contract> : (i32) -> i32
! CHECK: } else {
! CHECK: %[[V_45:[0-9]+]] = fir.call @_FortranAMapException(%c32{{.*}}) fastmath<contract> : (i32) -> i32
! CHECK: %[[V_46:[0-9]+]] = fir.call @feraiseexcept(%[[V_45]]) fastmath<contract> : (i32) -> i32
! CHECK: }
! CHECK: }
! CHECK: fir.result %[[V_42]] : f32
! CHECK: }
! CHECK: hlfir.assign %[[V_26]] to %[[V_7]]#0 : f32, !fir.ref<f32>
x4 = ieee_real(x8,4)
! print*, x8, ' -> ', x4
end