/* ************************************************************************
 * Copyright 2018-2021 Advanced Micro Devices, Inc.
 * ************************************************************************ */

#pragma once

#include "bytes.hpp"
#include "cblas_interface.hpp"
#include "flops.hpp"
#include "norm.hpp"
#include "rocblas.hpp"
#include "rocblas_init.hpp"
#include "rocblas_math.hpp"
#include "rocblas_random.hpp"
#include "rocblas_test.hpp"
#include "rocblas_vector.hpp"
#include "unit.hpp"
#include "utility.hpp"

template <typename T, bool CONJ>
void testing_ger_strided_batched_bad_arg(const Arguments& arg)
{
    auto rocblas_ger_strided_batched_fn
        = arg.fortran ? (CONJ ? rocblas_ger_strided_batched<T, true, true>
                              : rocblas_ger_strided_batched<T, false, true>)
                      : (CONJ ? rocblas_ger_strided_batched<T, true, false>
                              : rocblas_ger_strided_batched<T, false, false>);

    rocblas_int M           = 100;
    rocblas_int N           = 100;
    rocblas_int incx        = 1;
    rocblas_int incy        = 1;
    rocblas_int lda         = 100;
    T           alpha       = 0.6;
    rocblas_int abs_incx    = incx >= 0 ? incx : -incx;
    rocblas_int abs_incy    = incy >= 0 ? incy : -incy;
    rocblas_int stride_a    = lda * N;
    rocblas_int stride_x    = abs_incx * M;
    rocblas_int stride_y    = abs_incy * N;
    rocblas_int batch_count = 5;

    rocblas_local_handle handle{arg};

    size_t size_A = stride_a * batch_count;
    size_t size_x = stride_x * batch_count;
    size_t size_y = stride_y * batch_count;

    // allocate memory on device
    device_vector<T> dA_1(size_A);
    device_vector<T> dx(size_x);
    device_vector<T> dy(size_y);
    CHECK_DEVICE_ALLOCATION(dA_1.memcheck());
    CHECK_DEVICE_ALLOCATION(dx.memcheck());
    CHECK_DEVICE_ALLOCATION(dy.memcheck());

    EXPECT_ROCBLAS_STATUS((rocblas_ger_strided_batched_fn(handle,
                                                          M,
                                                          N,
                                                          &alpha,
                                                          nullptr,
                                                          incx,
                                                          stride_x,
                                                          dy,
                                                          incy,
                                                          stride_y,
                                                          dA_1,
                                                          lda,
                                                          stride_a,
                                                          batch_count)),
                          rocblas_status_invalid_pointer);
    EXPECT_ROCBLAS_STATUS((rocblas_ger_strided_batched_fn(handle,
                                                          M,
                                                          N,
                                                          &alpha,
                                                          dx,
                                                          incx,
                                                          stride_x,
                                                          nullptr,
                                                          incy,
                                                          stride_y,
                                                          dA_1,
                                                          lda,
                                                          stride_a,
                                                          batch_count)),
                          rocblas_status_invalid_pointer);
    EXPECT_ROCBLAS_STATUS((rocblas_ger_strided_batched_fn(handle,
                                                          M,
                                                          N,
                                                          &alpha,
                                                          dx,
                                                          incx,
                                                          stride_x,
                                                          dy,
                                                          incy,
                                                          stride_y,
                                                          nullptr,
                                                          lda,
                                                          stride_a,
                                                          batch_count)),
                          rocblas_status_invalid_pointer);
    EXPECT_ROCBLAS_STATUS((rocblas_ger_strided_batched_fn(nullptr,
                                                          M,
                                                          N,
                                                          &alpha,
                                                          dx,
                                                          incx,
                                                          stride_x,
                                                          dy,
                                                          incy,
                                                          stride_y,
                                                          dA_1,
                                                          lda,
                                                          stride_a,
                                                          batch_count)),
                          rocblas_status_invalid_handle);
}

template <typename T, bool CONJ>
void testing_ger_strided_batched(const Arguments& arg)
{
    auto rocblas_ger_strided_batched_fn
        = arg.fortran ? (CONJ ? rocblas_ger_strided_batched<T, true, true>
                              : rocblas_ger_strided_batched<T, false, true>)
                      : (CONJ ? rocblas_ger_strided_batched<T, true, false>
                              : rocblas_ger_strided_batched<T, false, false>);

    rocblas_int M           = arg.M;
    rocblas_int N           = arg.N;
    rocblas_int incx        = arg.incx;
    rocblas_int incy        = arg.incy;
    rocblas_int lda         = arg.lda;
    T           h_alpha     = arg.get_alpha<T>();
    rocblas_int stride_x    = arg.stride_x;
    rocblas_int stride_y    = arg.stride_y;
    rocblas_int stride_a    = arg.stride_a;
    rocblas_int batch_count = arg.batch_count;

    rocblas_local_handle handle{arg};

    size_t abs_incx = incx >= 0 ? incx : -incx;
    size_t abs_incy = incy >= 0 ? incy : -incy;
    size_t size_A   = size_t(lda) * N;
    size_t size_x   = M * abs_incx;
    size_t size_y   = N * abs_incy;

    // argument check before allocating invalid memory
    bool invalid_size = M < 0 || N < 0 || lda < M || lda < 1 || !incx || !incy || batch_count < 0;
    if(invalid_size || !M || !N || !batch_count)
    {
        EXPECT_ROCBLAS_STATUS((rocblas_ger_strided_batched_fn(handle,
                                                              M,
                                                              N,
                                                              nullptr,
                                                              nullptr,
                                                              incx,
                                                              stride_x,
                                                              nullptr,
                                                              incy,
                                                              stride_y,
                                                              nullptr,
                                                              lda,
                                                              stride_a,
                                                              batch_count)),
                              invalid_size || batch_count < 0 ? rocblas_status_invalid_size
                                                              : rocblas_status_success);

        return;
    }

    size_A += size_t(stride_a) * size_t(batch_count - 1);
    size_x += size_t(stride_x) * size_t(batch_count - 1);
    size_y += size_t(stride_y) * size_t(batch_count - 1);

    // Naming: dK is in GPU (device) memory. hK is in CPU (host) memory
    host_vector<T> hA_1(size_A);
    host_vector<T> hA_2(size_A);
    host_vector<T> hA_gold(size_A);
    host_vector<T> hx(size_x);
    host_vector<T> hy(size_y);

    // allocate memory on device
    device_vector<T> dA_1(size_A);
    device_vector<T> dA_2(size_A);
    device_vector<T> dx(size_x);
    device_vector<T> dy(size_y);
    device_vector<T> d_alpha(1);
    CHECK_DEVICE_ALLOCATION(dA_1.memcheck());
    CHECK_DEVICE_ALLOCATION(dA_2.memcheck());
    CHECK_DEVICE_ALLOCATION(dx.memcheck());
    CHECK_DEVICE_ALLOCATION(dy.memcheck());
    CHECK_DEVICE_ALLOCATION(d_alpha.memcheck());

    double gpu_time_used, cpu_time_used;
    double rocblas_error_1;
    double rocblas_error_2;

    // Initialize data on host memory
    rocblas_init_matrix(
        hA_1, arg, M, N, lda, stride_a, batch_count, rocblas_client_never_set_nan, true);
    rocblas_init_vector(
        hx, arg, M, abs_incx, stride_x, batch_count, rocblas_client_alpha_sets_nan, false, true);
    rocblas_init_vector(hy, arg, N, abs_incy, stride_y, batch_count, rocblas_client_alpha_sets_nan);

    // copy matrix is easy in STL; hA_gold = hA_1: save a copy in hA_gold which will be output of
    // CPU BLAS
    hA_gold = hA_1;
    hA_2    = hA_1;

    // copy data from CPU to device
    CHECK_HIP_ERROR(hipMemcpy(dA_1, hA_1, sizeof(T) * size_A, hipMemcpyHostToDevice));
    CHECK_HIP_ERROR(hipMemcpy(dx, hx, sizeof(T) * size_x, hipMemcpyHostToDevice));
    CHECK_HIP_ERROR(hipMemcpy(dy, hy, sizeof(T) * size_y, hipMemcpyHostToDevice));

    if(arg.unit_check || arg.norm_check)
    {
        // copy data from CPU to device
        CHECK_HIP_ERROR(hipMemcpy(dA_2, hA_2, sizeof(T) * size_A, hipMemcpyHostToDevice));
        CHECK_HIP_ERROR(hipMemcpy(d_alpha, &h_alpha, sizeof(T), hipMemcpyHostToDevice));
        CHECK_ROCBLAS_ERROR(rocblas_set_pointer_mode(handle, rocblas_pointer_mode_host));
        CHECK_ROCBLAS_ERROR((rocblas_ger_strided_batched_fn(handle,
                                                            M,
                                                            N,
                                                            &h_alpha,
                                                            dx,
                                                            incx,
                                                            stride_x,
                                                            dy,
                                                            incy,
                                                            stride_y,
                                                            dA_1,
                                                            lda,
                                                            stride_a,
                                                            batch_count)));
        CHECK_ROCBLAS_ERROR(rocblas_set_pointer_mode(handle, rocblas_pointer_mode_device));
        CHECK_ROCBLAS_ERROR((rocblas_ger_strided_batched_fn(handle,
                                                            M,
                                                            N,
                                                            d_alpha,
                                                            dx,
                                                            incx,
                                                            stride_x,
                                                            dy,
                                                            incy,
                                                            stride_y,
                                                            dA_2,
                                                            lda,
                                                            stride_a,
                                                            batch_count)));

        // CPU BLAS
        cpu_time_used = get_time_us_no_sync();

        for(int b = 0; b < batch_count; ++b)
        {
            cblas_ger<T, CONJ>(M,
                               N,
                               h_alpha,
                               hx + b * stride_x,
                               incx,
                               hy + b * stride_y,
                               incy,
                               hA_gold + b * stride_a,
                               lda);
        }

        cpu_time_used = get_time_us_no_sync() - cpu_time_used;

        // copy output from device to CPU
        hipMemcpy(hA_1, dA_1, sizeof(T) * size_A, hipMemcpyDeviceToHost);
        hipMemcpy(hA_2, dA_2, sizeof(T) * size_A, hipMemcpyDeviceToHost);

        if(arg.unit_check)
        {
            if(std::is_same<T, float>{} || std::is_same<T, double>{})
            {
                unit_check_general<T>(M, N, lda, stride_a, hA_gold, hA_1, batch_count);
                unit_check_general<T>(M, N, lda, stride_a, hA_gold, hA_2, batch_count);
            }
            else
            {
                const double tol = N * sum_error_tolerance<T>;
                near_check_general<T>(M, N, lda, stride_a, hA_gold, hA_1, batch_count, tol);
                near_check_general<T>(M, N, lda, stride_a, hA_gold, hA_2, batch_count, tol);
            }
        }

        if(arg.norm_check)
        {
            rocblas_error_1
                = norm_check_general<T>('F', M, N, lda, stride_a, hA_gold, hA_1, batch_count);
            rocblas_error_2
                = norm_check_general<T>('F', M, N, lda, stride_a, hA_gold, hA_2, batch_count);
        }
    }

    if(arg.timing)
    {
        int number_cold_calls = arg.cold_iters;
        int number_hot_calls  = arg.iters;
        CHECK_ROCBLAS_ERROR(rocblas_set_pointer_mode(handle, rocblas_pointer_mode_host));

        for(int iter = 0; iter < number_cold_calls; iter++)
        {
            rocblas_ger_strided_batched_fn(handle,
                                           M,
                                           N,
                                           &h_alpha,
                                           dx,
                                           incx,
                                           stride_x,
                                           dy,
                                           incy,
                                           stride_y,
                                           dA_1,
                                           lda,
                                           stride_a,
                                           batch_count);
        }

        hipStream_t stream;
        CHECK_ROCBLAS_ERROR(rocblas_get_stream(handle, &stream));
        gpu_time_used = get_time_us_sync(stream); // in microseconds

        for(int iter = 0; iter < number_hot_calls; iter++)
        {
            rocblas_ger_strided_batched_fn(handle,
                                           M,
                                           N,
                                           &h_alpha,
                                           dx,
                                           incx,
                                           stride_x,
                                           dy,
                                           incy,
                                           stride_y,
                                           dA_1,
                                           lda,
                                           stride_a,
                                           batch_count);
        }

        gpu_time_used = get_time_us_sync(stream) - gpu_time_used;

        ArgumentModel<e_M,
                      e_N,
                      e_alpha,
                      e_lda,
                      e_stride_a,
                      e_incx,
                      e_stride_x,
                      e_incy,
                      e_stride_y,
                      e_batch_count>{}
            .log_args<T>(rocblas_cout,
                         arg,
                         gpu_time_used,
                         ger_gflop_count<T>(M, N),
                         ger_gbyte_count<T>(M, N),
                         cpu_time_used,
                         rocblas_error_1,
                         rocblas_error_2);
    }
}
