Современные
информационные технологии/Вычислительная техника и программирование
Мясищев А.А.
Хмельницкий национальный
университет, Украина
СОПОСТАВЛЕНИЕ
ПРОИЗВОДИТЕЛЬНОСТЕЙ GPU TESLA C2075, GEFORCE 480 GTX С
6-И ЯДЕРНЫМ CPU AMD ПРИ РЕШЕНИИ СИСТЕМ ЛИНЕЙНЫХ УРАВНЕНИЙ
Задачи теории упругости, пластичности сводятся к системе дифференциальных уравнений в напряжениях и скоростях. Обычно
они решаются численно с использованием метода конечных элементов. В свою очередь этот метод сводится к построению так
называемых матриц
жесткости
и решению систем линейных уравнений, содержащей тысячи и даже десятки тысяч неизвестных. Для
упругопластических и пластических задач решение систем уравнений приходится выполнять сотни и тысячи раз для выявления конечной конфигурации деформируемого материала. Для
этих целей требуются вычислительные системы высокой производительности.
Для проведения вычислительных экспериментов рассмотрим две системы.
В каждой из них установлен процессор (CPU) AMD Phenom II X6 1090T с
частотой 3.2 ГГц. В первой системе установлен графический процессор (GPU) Tesla C2075, во второй GeForce GTX 480 фирмы
NVIDIA. GPU будут
выполнять роль вычислительных ускорителей для решения систем линейных
уравнений. Причем в первом случае расчет будет выполняться c учетом распараллеливания по 6-ядрам процессора с
использованием библиотек MPI, ScaLAPACK и ATLAS [1,2,3]. Во втором и третьем случаях – распараллеливанием
по ядрам GPU Tesla C2075 и GeForce GTX 480 с использованием
технологии CUDA[4]. Вычислительные системы работают
под управлением операционной системы Linux Ubuntu
ver. 10.10 desktop (64bit). На
них установлены компилятор фортран F77, gfortran с перечисленными выше
библиотеки для 6-и ядерного процессора. Для программирования на GPU Tesla C2075 и GeForce GTX 480 инсталлированы
видеодрайвер nvidia и программное обеспечение CUDA Toolkit
4.0 с сайта http://developer.nvidia.com/cuda-toolkit-archive.
Последовательность установки программного обеспечения для GPU представлена в
источнике [5]. Установка библиотек MPI, ScaLAPACK, ATLAS под ОС Linux
Ubuntu ver. 9.04 desktop
представлена в работе [6] для системы с
4-х ядерным процессором CORE
2 QUAD PENTIUM Q6600
2.4GHZ, поэтому здесь не рассматривается.
Рассмотрим лишь установку библиотеки BLACS, которая
в работе [6] используется готовой для 32-битной системы.
С сайта ftp.netlib.org копируем файл mpiblacs.tar.gz в каталог $(HOME)/fortran
Разархивируем его
gzip –d mpiblacs.tar.gz
tar xf mpiblacs.tar
Заходим в каталог BLACS
Из подкаталога BMAKES копируем файл Bmake.MPI-LINUX в каталог
BLACS под именем Bmake.inc
cp BMAKES/Bmake.MPI-LINUX ./Bmake.inc
Выборочно редактируем Bmake.inc (привязываем его к компилятору и каталогам библиотек MPI)
BTOPdir = $(HOME)/fortran/BLACS
MPIdir = /usr/local/mpi/ch_shmem/
F77 = f77
Редактируем Makefile. Вместо
строки
all : mpi cmmd mpl nx pvm tester
вставляем
all : mpi
Выполняем команду
make all
После компиляции в подкаталоге LIB должны появиться файлы библиотеки BLACS. При компиляции программ на CPU для решения систем уравнений использовались библиотеки scalapack-1.8.0
и atlas3.8.4.
Ниже
представлены командные строки для компиляции и запуска программы ur.f на 6-и ядерном процессоре для решения системы уравнений
методом LUP разложения.
f77 -O3 -I/usr/local/mpi/ch_shmem/include -f -o ur ur.f
libscalapack.a blacsF77.a blacs.a blacsF77.a -lgfortran libf77blas.a libatlas.a
/usr/local/mpi/ch_shmem/lib/libmpich.a
/usr/local/mpi/ch_shmem/bin/mpirun ur -np 6
machinefile machines
Текст
программы ur.f:
program solve
system
include
'mpif.h'
parameter
(nm = 11008, nxn = nm*nm)
double
precision a(nxn), b(nm), mem(nm)
double
precision aa(nm,nm),bb(nm,1)
double
precision ops1,ops2,total1,total2,t1,flop1,flop2
INTEGER
DESCA(9), DESCB(9)
integer
ipvt(nm)
CALL BLACS_PINFO( IAM, NPROCS )
NPROW =
INT(SQRT(REAL(NPROCS)))
NPCOL =
NPROCS/NPROW
IF(
IAM.EQ.0 ) THEN
print *,'
Input N and NB: '
read *,
N, NB
END IF
call
MPI_BCAST(N, 1, MPI_INTEGER, 0, MPI_COMM_WORLD, ierr)
call
MPI_BCAST(NB,1, MPI_INTEGER, 0, MPI_COMM_WORLD, ierr)
ops1 =
(2.0d0*dfloat(n)**3)/3.0d0
ops2 =
(2.0d0*(dfloat(n))**3)/3.0d0 + 2.0d0*(dfloat(n))**2
CALL
BLACS_GET( -1, 0, ICTXT )
CALL
BLACS_GRIDINIT( ICTXT, 'Row-major', NPROW, NPCOL )
CALL
BLACS_GRIDINFO( ICTXT, NPROW, NPCOL, MYROW, MYCOL )
IF(
MYROW.GE.NPROW .OR. MYCOL.GE.NPCOL ) GOTO 500
NP =
NUMROC( N, NB, MYROW, 0, NPROW )
NQ =
NUMROC( N, NB, MYCOL, 0, NPCOL )
CALL
DESCINIT( DESCA, N, N, NB,NB,0,0,ICTXT, NP, INFO )
CALL
DESCINIT( DESCB, N, 1, NB,1,0,0,ICTXT, NP, INFO )
lda =
DESCA(9)
c initial data
do 10
i=1,n
do 10
j=1,n
aa(i,j)=(sqrt(i/1.0d0))/11.0d0+(sqrt(j/1.0d0))/12.0d0
if(i.eq.j) aa(i,j)=10.0d0
10 continue
do 11
i=1,n
bb(i,1)=(sqrt(i/1.0d0))/20.0d0 - 1.0d0
11 continue
c initial matrix in processors
call
pmatgen(a,aa,DESCA,np,nq,b,bb,DESCB,nprow,npcol
&,myrow,mycol,nm)
t1 =
MPI_Wtime()
CALL
PDGETRF(N, N, A, 1, 1, DESCA,ipvt, info)
total1 =
MPI_Wtime() - t1
t1 =
MPI_Wtime()
CALL PDGETRS('n', n,1,a,1,1,DESCA,ipvt,b,1,1,DESCB,info)
total2 =
MPI_Wtime() - t1
if
(IAM.EQ.0) then
write(6,110) total1, total2
110
format(2x,'Time LU:',f12.4, ' sec.',' Time solv: = ',f12.4)
flop1=ops1/(1.0d6*total1)
flop2=ops2/(1.0d6*(total1+total2))
write(6,111)
flop1, flop2
111
format(2x,'Flops LU:',f12.4, ' Mflops', ' Flops solv: = ',f12.4)
end if
CALL
PDLAPRNT( 1, 1, A, n,n, DESCA, 0, 0, 'matrC', 6, MEM )
CALL
PDLAPRNT( 1, 1, B, 1,1, DESCB, 0, 0, 'matrX', 6, MEM )
CALL
PDLAPRNT( 1, 1, B, n,1, DESCB, 0, 0, 'matrX', 6, MEM )
CALL
BLACS_GRIDEXIT( ICTXT )
CALL
BLACS_EXIT(0)
500 continue
stop
end
subroutine pmatgen(a,aa,DESCA,np,nq,
&b,bb,DESCB,nprow, npcol,myrow,mycol,nm)
integer
i, j, DESCA(9), DESCB(9), nprow, npcol,myrow,mycol
double
precision a(*), b(*),aa(11008,11008),
&bb(11008,1)
nb = DESCA(5)
k = 0
do 250 j
= 1,nq
jc =
(j-1)/nb
jb =
mod(j-1,nb)
jg =
mycol*nb + jc*npcol*nb + jb + 1
do 240 i
= 1,np
ic =
(i-1)/nb
ib =
mod(i-1,nb)
ig =
myrow*nb + ic*nprow*nb + ib + 1
k = k + 1
a(k) = aa(ig,jg)
if(jg.eq.1) b(k) = bb(ig,1)
240 continue
250 continue
return
end
В
этой программе расчет проводиться для чисел с двойной точностью. Решение
системы уравнений разбито на два этапа:
-факторизация матрицы системы (A) [7];
-само решение системы с использованием результатов,
полученных на первом этапе.
Каждый этап выполняется
посредством обращения к соответствующим подпрограммам библиотек ScaLAPACK, BLACS, ATLAS. Первый этап использует подпрограмму PDGETRF, второй – PDGETRS. В
программе для каждой подпрограммы считается время ее выполнения. Рассчитываются количество операций
для первого (ops1) и второго этапа (ops2), и вычисляется производительность системы для LU-факторизации (flop1) и
решения системы в целом (flop2).
В подпрограмме PDGETRF
для квадратных заполненных матриц A применяется LU –факторизация
с выбором ведущего элемента по столбцам:
P*A = L*U,
где P – матрица перестановок; L – нижняя
треугольная матрица с единичными диагональными элементами; U – верхняя треугольная матрица. Например, для матрицы размером
3x3 перестановки будут выглядеть так:
2 1 3
11 7 5 11 7 5
11
7 5
=> 2 1
3 => 9
8 4
9 8 4 9 8 4 2 1 3
На первом этапе в первом столбце находим строку с наибольшим
по модулю значением (11) и меняем первую строку на найденную. На втором этапе во
втором столбце находим строку с наибольшим по модулю числом (8), однако рассматриваем только 2-ю и 3-ю строки и
меняем 2-ю строку на 3-ю. Если числа в столбцах одинаковые, то перестановка
выполняется с минимальным номером строки. Такие перестановки позволяют
исключить вероятность деления на ноль при LU - факторизации.
Подпрограмма pmatgen
распределяет матрицы A (размером n x n) и B(вектор
правой части размером n) по ядрам CPU.
Рассмотрим
процедуру решения этой же задачи на GPU. Для
обеспечения максимальной производительности будем также использовать готовые
библиотеки для решения систем уравнений методом LU – факторизации. С этой целью установим на вычислительных
системах дополнительно библиотеки lapack-3.4.0 и MAGMA (Matrix Algebra on GPU and Multicore Architectures)
версии 1.1.0.
С
сайта копируем lapack-3.4.0.tgz и разархивируем его:
gzip
–d lapack-3.4.0.tgz
tar xf lapack-3.4.0.tar
Заходим в каталог lapack-3.4.0
cd lapack-3.4.0/
и переименовываем файл компиляции make.inc.example на
make.inc
mv make.inc.example make.inc
Выполняем команды
make lapacklib
make blaslib
В текущем каталоге появятся файлы liblapack.a и
librefblas.a
Далее выполняем копирование исходника magma_1.1.0.tar.gz с сайта http://www.cs.utk.edu/~tomov/magma_1.1.0.tar.gz
и разархивируем
его
$gzip –d magma_1.1.0.tag.gz
$tar xf magma_1.1.0.tar
Переходим в каталог magma_1.1.0:
$cd magma_1.1.0
С помощью редактора nano создаем
файл с именем make.inc для
компиляции библиотеки magma совместно
с ранее установленными библиотеками:
GPU_TARGET = 1
CC = gcc
NVCC = nvcc
FORT = gfortran
ARCH = ar
ARCHFLAGS = cr
RANLIB = ranlib
OPTS = -O3 -DADD_
FOPTS = -O3 -DADD_ -x f95-cpp-input
NVOPTS = --compiler-options -fno-strict-aliasing -DUNIX -O3 -DADD_
LDOPTS = -fPIC -Xlinker
-zmuldefs
LIB =
/home/alex/cuda/ATLAS/my/lib/liblapack.a \
/home/alex/cuda/ATLAS/my/lib/libf77blas.a \
/home/alex/cuda/ATLAS/my/lib/libcblas.a \
/home/alex/cuda/ATLAS/my/lib/libatlas.a \
/home/alex/cuda/lapack-3.4.0/liblapack.a \
/home/alex/cuda/lapack-3.4.0/librefblas.a \
-lgfortran -lpthread -lcublas -lcudart -lm
CUDADIR = /usr/local/cuda
LIBDIR = -L/usr/local/cuda/lib64 -L/usr/lib
INC = -I$(CUDADIR)/include
LIBMAGMA = ../lib/libmagma.a
LIBMAGMABLAS = ../lib/libmagmablas.a
Запустим команду компиляции make all.
В этом случае будут создаваться библиотеки libmagma, libmagmablas и выполнится
компиляция тестовых примеров. Однако здесь возникают следующие проблемы.
../lib/libmagma.a(zlatrd.o):
In function `magma_zlatrd':
zlatrd.cpp:(.text+0x332): undefined reference to
`zdotc'
zlatrd.cpp:(.text+0xbe2): undefined reference to
`zdotc'
Для решения этой проблемы необходимо зайти в программу src/zlatrd.cpp и
сделать изменения
cblas_zdotc_sub(i_n,
W(i +1, i), ione, A(i +1, i), ione, &value);
// blasf77_zdotc(&value, &i_n,
W(i+1,i), &ione, A(i+1, i), &ione);
и
cblas_zdotc_sub(i, W(0, iw), ione, A(0, i), ione, &value);
// blasf77_zdotc(&value, &i, W(0,
iw), &ione, A(0, i), &ione);
Вторая проблема:
../lib/libmagma.a(clatrd.o): In function
`magma_clatrd':
clatrd.cpp:(.text+0x321): undefined reference to
`cdotc'
clatrd.cpp:(.text+0xb30): undefined reference to
`cdotc'
Аналогично
необходимо зайти в программу src/clatrd.cpp и сделать
изменения
//
blasf77_cdotc(&value, &i, W(0, iw), &ione, A(0, i),
&ione);
cblas_cdotc_sub(i, W(0, iw), ione, A(0, i), ione, &value);
и
cblas_cdotc_sub(i_n, W(i +1, i), ione, A(i +1, i), ione, &value);
//
blasf77_cdotc(&value, &i_n, W(i+1,i), &ione, A(i+1, i),
&ione);
После этого повторить компиляцию по команде
make all
Ниже представлена программа ur_f1.f90 для решения системы
уравнений методом LU-факторизации:
program testing_dgetrf_gpu_f
use magma
external cublas_init, cublas_set_matrix, cublas_get_matrix
external cublas_shutdown, cublas_alloc
integer cublas_alloc
double precision, allocatable :: h_A(:), h_B(:), h_X(:)
magma_devptr_t
:: devptrA, devptrB
integer, allocatable :: ipiv(:)
integer
:: i, n, info, stat, lda
integer
:: size_of_elt, nrhs
real(kind=8)
:: ops,ops1, t, t1
integer
:: tstart(2), tend(2)
integer
:: tstart1(2), tend1(2)
call cublas_init()
write (*,*) "Input n"
read (*,*) n
nrhs = 1
lda = n
ldda = ((n+31)/32)*32
size_of_elt = sizeof_double
!------ Allocate CPU memory
allocate(h_A(lda*n))
allocate(h_B(lda*nrhs))
allocate(h_X(lda*nrhs))
allocate(ipiv(n))
!------ Allocate GPU memory
stat = cublas_alloc(ldda*n, size_of_elt, devPtrA)
stat = cublas_alloc(ldda*nrhs, size_of_elt, devPtrB)
!---- Initializa the matrix
do 10 i=1,n
do 10 j=1,n
h_A(i+(j-1)*n)=(sqrt(i/1.0d0))/11.0d0+(sqrt(j/1.0d0))/12.0d0
if(i.eq.j) h_A(i+(j-1)*n)=10.0d0
10 continue
do 11 i=1,n
h_B(i)=(sqrt(i/1.0d0))/20.0d0 - 1.0d0
11 continue
h_X(:) = h_B(:)
!---- devPtrA = h_A
call cublas_set_matrix(n, n, size_of_elt, h_A, lda, devptrA,
ldda)
!---- devPtrB = h_B
call cublas_set_matrix(n, nrhs, size_of_elt, h_B, lda,
devptrB, ldda)
!---- Call magma LU
call magma_gettime_f(tstart)
call
magmaf_dgetrf_gpu(n, n, devptrA, ldda, ipiv, info)
call magma_gettime_f(tend)
call magma_gettimervalue_f(tstart, tend, t)
write(*,*) "time LU =", t
!---- Call magma solve -------------
call magma_gettime_f(tstart1)
call magmaf_dgetrs_gpu('n', n, nrhs, devptrA, ldda, ipiv,
devptrB, ldda, info)
call magma_gettime_f(tend1)
call magma_gettimervalue_f(tstart1, tend1, t1)
write(*,*) "time solve =", t1
!---- h_X = devptrB
call cublas_get_matrix (n, nrhs, size_of_elt, devptrB, ldda,
h_X, lda)
write(*,*)
"X(1)=",h_X(1),"X(",n,")=",h_X(n)
ops = 2.0d0 * (dfloat(n))**3 / 3.0d0
ops1 = 2.0d0 * (dfloat(n))**3 / 3.0d0 + 2.0d0*(dfloat(n))**2
write(*,*) ' Gflops LU
= ', ops /(t*1e6)
write(*,*) '
Gflops summa solve = ', ops1
/(1e6*(t+t1))
!---- Free CPU memory
deallocate(h_A, h_X, h_B, ipiv)
!---- Free GPU memory
call cublas_free(devPtrA)
call cublas_free(devPtrB)
call cublas_shutdown()
end
Здесь как и в предыдущей программе
используются подпрограммы для LU-факторизации
magmaf_dgetrf_gpu и решения системы magmaf_dgetrs_gpu. Для них также считается
время выполнения и производительности
вычислений отдельно для LU-факторизации
и решении системы в целом.
Компиляция программы выполняется в два этапа:
gfortran -O3 -DADD_ -x f95-cpp-input
-Dmagma_devptr_t="integer(kind=8)" -I/usr/local/cuda/include
-I../include -I../quark/include -c ur_f1.f90 -o ur_f1.o
и
gfortran -O3 -DADD_ -DGPUSHMEM=200 -fPIC -Xlinker -zmuldefs -DGPUSHMEM=200
ur_f1.o fortran.o -o ur_f1 lin/liblapacktest.a -L../lib -lcuda -lmagma -lmagmablas -lmagma
-L/usr/local/cuda/lib64 -L/usr/lib
/home/alex/cuda/ATLAS/my/lib/liblapack.a
/home/alex/cuda/ATLAS/my/lib/libf77blas.a /home/alex/cuda/ATLAS/my/lib/libcblas.a
/home/alex/cuda/ATLAS/my/lib/libatlas.a
/home/alex/cuda/lapack-3.4.0/liblapack.a
/home/alex/cuda/lapack-3.4.0/librefblas.a -lgfortran -lpthread -lcublas
-lcudart -lm
Рассмотрим
результаты вычислительных экспериментов. В
таблице 1 сведены результаты расчетов по представленным выше программам для CPU AMD и GPU Tesla C2075 при двойной точности вычислений. Для процессора даются два варианта: расчет
выполняется одним ядром и 6-ю ядрами. Представлена величина ускорения, которая
нигде не достигает шестикратной величины. Результаты приведены в виде дроби:
числитель – время счета в секундах, знаменатель – производительность в
Гигафлопсах в секунду. Отдельными столбцами представлены расчеты для LU-факторизации и
решения системы. В столбце “решение” в числителе дано время работы
подпрограммы PDGETRS,
т.е. обратного хода, а в знаменателе – суммарная производительность при решении
системы. Видно, что основная трудоемкость вычислений связана с LU-факторизацией. Для CPU она в 26.491/0.121=219
раз больше времени обратного хода, а для GPU в 3.718/0.071=52
раза для размера матрицы 11008x11008. Соотношение увеличивается при увеличении размера
уравнений. Таким образом, оценку производительности можно выполнять только по LU-факторизации.
Таблица 1.
Мат-рица NxN |
CPU AMD Phenom
II X6 1090T |
GPU Tesla C2075 |
|||||
6 ядер |
1 ядро |
Уско-рение |
LU-фактор. |
Решен-ние |
|||
LU-факт. |
Решен. |
LU-факт. |
Решен. |
||||
1920 |
0.306/15.4 |
0.011/14.9 |
0.598/7.89 |
0.010/7.77 |
1.92 |
0.076/61.7 |
0.010/54.4 |
3072 |
0.914/21.1 |
0.021/20.7 |
2.384/8.11 |
0.026/8.03 |
2.58 |
0.230/84.0 |
0.017/78.4 |
4992 |
2.620/31.7 |
0.028/31.3 |
8.508/9.75 |
0.064/9.68 |
3.24 |
0.518/160.0 |
0.029/151.7 |
7104 |
7.605/31.4 |
0.054/31.2 |
23.524/10.2 |
0.122/10.1 |
3.09 |
1.180/202.6 |
0.043/195.6 |
8064 |
10.880/32.1 |
0.067/31.9 |
33.915/10.3 |
0.152/10.3 |
3.10 |
1.684/207.6 |
0.050/201.7 |
9984 |
20.226/32.8 |
0.099/32.6 |
61.895/10.7 |
0.235/10.7 |
3.05 |
2.882/230.2 |
0.063/225.3 |
11008 |
26.491/33.6 |
0.121/33.4 |
82.275/10.8 |
0.281/10.8 |
3.09 |
3.718/239.2 |
0.071/234.7 |
25344 |
- |
- |
- |
- |
- |
38.665/280.7 |
0.204/279.2 |
В таблице 2
приведены результаты расчетов производительностей выполнения LU-факторизации в Гигафлопсах в секунду для GPU Tesla C2075, GeForce GTX
480 и CPU
AMD Phenom II X6
1090T для одинарной (real) и двойной (double) точностей.
Видно, что для размера уравнений до
3072, более дешевый GPU GeForce
480 обгоняет GPU Tesla C2075. При больших уравнениях
преимущество за Tesla C2075. В целом при двойной точности вычислений для размера
в 11008 уравнений GPU Tesla C2075 обгоняет процессор в 239.2/34.1=7.01
раз.
Таблица 2
Матрица NxN |
GPU Tesla C2075 |
GPU GeForce 480 |
CPU AMD |
|||
real |
double |
real |
double |
real |
double |
|
1920 |
77.9 |
61.7 |
90.9 |
66.4 |
18.3 |
15.4 |
3072 |
132.3 |
84.2 |
170.1 |
84.2 |
30.8 |
21.1 |
4992 |
150.9 |
160.3 |
161.2 |
128.5 |
50.0 |
31.7 |
7104 |
221.4 |
203.3 |
234.3 |
143.1 |
60.9 |
32.1 |
8064 |
251.4 |
208.1 |
263.4 |
145.0 |
58.5 |
32.7 |
9984 |
317.3 |
230.2 |
329.7 |
150.9 |
71.6 |
33.3 |
11008 |
351.4 |
239.2 |
365.3 |
153.2 |
73.9 |
34.1 |
В таблице 3 приведены результаты вычислений первого
аргумента системы уравнений. Здесь Xr(1)
соответствует одинарной точности, а Xd(1)
- двойной точности. Видно, что чем меньше размер уравнений, тем меньше
отличаются Xr(1) и Xd(1).
При размере системы в 3000 уравнений можно “верить” только первым трем значащим
цифрам первого аргумента, а при размере 10000 уравнений погрешность становиться
такой большой, что и первые значащие цифры становятся под вопросом. Таким образом,
для систем свыше 3000 уравнений целесообразно использовать при расчетах только
числа с двойной точностью.
Таблица 3
Матрица NxN |
Аргумент Xr(1), одинарная точн. |
Аргумент Xd(1), двойная точность |
Разность |Xr(1)-Xd(1)| |
100x100 |
-9.904994E-02 |
-9.904992E-02 |
2.0E-06 |
1000x1000 |
7.36474E-03 |
7.36434E-03 |
4.0E-04 |
2000x2000 |
2.50738E-03 |
2.50710E-03 |
2.80E-04 |
3000x3000 |
1.16325E-03 |
1.16273E-03 |
4.80E-04 |
4000x4000 |
4.93012E-04 |
5.01715E-04 |
8.70E-02 |
5000x5000 |
4.92534E-04 |
5.03876E-04 |
1.13E-01 |
10000X10000 |
7.34304E-04 |
7.90892E-04 |
5.66E-01 |
Решение системы в 11008 уравнений с числами двойной точности для Xd(1) показали следующие результаты:
GPU Tesla C2075 Xd(1)=9.9392662447128E-004
CPU AMD Phenom II X6=9.9392662452686E-004
Здесь можно предположить, что первые девять значащих цифр
являются правильными.
Выводы.
1. При решении систем уравнений оценку производительности
можно выполнять только по LU-факторизации.
2. Для систем свыше 3000 уравнений целесообразно при
расчетах использовать только числа с двойной точностью.
3. При решении систем до 3000 уравнений целесообразно
использовать более дешевый GPU GeForce 480 для вычислений с любой
точностью. Он работает в этом диапазоне уравнений быстрее CPU AMD Phenom II X6 1090T в 4 раза для двойной точности и в 5.5 раза для одинарной.
4.Для системы в 11008 уравнений при двойной точности вычислений GPU обгоняет CPU в 7 раз.
5. При решении систем уравнений на 6-и ядерном CPU удается достичь максимального ускорения в 3.24 раза по
сравнению с расчетами на 1-м ядре для используемых библиотек ScaLAPACK и ATLAS.
6. Преимущество GPU по сравнению с CPU
возрастает с увеличением числа уравнений. Например, GPU Tesla C2075
для 1920 уравнений в системе работает в 4 раза быстрее CPU AMD, а для 11008 уравнений в 7 раз.
Литература
1.Антонов А.С.
Параллельное программирование с использованием технологии MPI: Учебное пособие.
– М.: Изд-во МГУ, 2004. – 71 с.
2. The ScaLAPACK Project. [Electronic resource]. - Mode of access: http://www.netlib.org
3. Automatically Tuned Linear
Algebra Software (ATLAS). [Electronic resource]. - Mode of access: http://math-atlas.sourceforge.net/
4.Боресков А.В., Харламов
А.А. Основы работы с технологией CUDA.
- М.: ДМК Пресс, 2010. -232 с.: ил.
5.Установка CUDA Toolkit и драйвера NVIDIA для разработчиков.
[Electronic resource]. - Mode of access: http://forum.ubuntu.ru/index.php?topic=114802.0
6.Мясищев А.А. Достижение наибольшей
производительности перемножения матриц на системах с многоядерными
процессорами. - Кривий Ріг: Видавничний відділ НметАУ, 2010. – Т. 3: Теорія та
методика навчання інформатики. – 303 с.