ROOT на GPU (OpenCL)

С самого начала своего развития ROOT поддерживал многопоточность и параллелизм. Однако обеспечение безопасности потоков требовало определенных усилий от пользователей. ROOT 6 и новый интерпретатор CLING упростили разработку многопоточных и параллельных приложений.

Современные вычислительные архитектуры часто имеют герерогенную структуру, то есть состоят из разнородных групп вычислительных узлов - центрального процессора (CPU) с несколькими ядрами и графических процессоров (GPU). Последние модели графических процессоров содержат сотни небольших специализированных процессоров, которые могут одновременно выполнять простые математические операции над входящими потоками данных и поэтому имеют возможность полноценно применяться для общих вычислений. GPGPU - (также GPGP, GP²U, англ. General-purpose computing on graphics processing units, неспециализированные вычисления на графических процессорах) - техника использования графического процессора видеокарты, предназначенного для компьютерной графики, для выполнения расчётов в приложениях для общих вычислений, которые обычно проводит центральный процессор. Использовать вычислительную мощь графических процессоров для неграфических задач возможно, но только в том случае, если решаемая задача допускает возможность распараллеливания алгоритмов на сотни исполнительных блоков, имеющихся в GPU. В частности, выполнение расчетов на GPU показывает отличные результаты в случае, когда одна и та же последовательность математических операций применяется к большому объему данных. При этом лучшие результаты достигаются, если отношение числа арифметических инструкций к числу обращений к памяти достаточно велико.

Существует несколько реализаций GPGPU, наиболее распространенные - OpenCL и CUDA.

В настоящее время CUDA является более быстрой технологией, но она применима только для GPU. Преимущество OpenCL состоит в том, что он является единым стандартом для написания приложений, которые должны исполняться в системе, где установлены различные по архитектуре процессоры, ускорители и платы расширения. Рассмотрим, как запускать макросы ROOT с использованием OpenCL. Подробно о этой технологии можно прочитать здесь, здесь и здесь.

Работа OpenCL-приложения.


Приложение, использующее OpenCL, состоит из Host - управляющей части и Kernel - загружаемой реализации OpenCL, в которой непосредственно происходят вычисления.

Первоначально OpenCL-приложение опрашивает имеющиеся OpenCL - платформы. Из списка найденных платформ приложение выбирает какую-то одну нужного типа (CPU, GPU и т.п.) и создаёт контекст - некое окружение, в котором будут запускаться на исполнение ядра. Контекст включает в себя информацию о наборе устройств, существующих в рамках платформы, памяти, доступной устройствам,а также очереди команд(command queues), используемые для организации исполнения ядер или операций над объектами памяти.

Взаимодействие хоста и устройств осуществляется с помощью команд, а для доставки этих команд устройствам и используются очереди команд. Одновременно с командой можно создать объект события (event). Такие объекты позволяют приложению проверять завершение исполнения команд, а потому могут использоваться для синхронизации. Для выделения памяти на устройствах создаются объекты памяти; свойства их (например, возможность чтения/записи) устанавливаются приложением. Программные объекты (program objects) создаются загрузкой исходного или бинарного представления одного или нескольких ядер и последующей процедурой построения (build) исполняемого кода. В результате возникают объекты ядер (kernel objects), которые обеспечивает параллелизм на уровне инструкций и на уровне данных.

На некоторых видеокартах по умолчанию отключен режим работы с числами типа double, что приводит к возникновению ошибки компиляции 5105. Для включения режима поддержки чисел double в текст OpenCL-программы нужно добавить директиву #pragma OPENCL EXTENSION cl_khr_fp64 : enable. Однако, если видеокарта не поддерживает double, то включение данной директивы не поможет.

Функции для выполнения программ на OpenCL


Функция Действие
CLHandleType - Возвращает тип OpenCL хендла в виде значения из перечисления ENUM_OPENCL_HANDLE_TYPE
CLGetInfoInteger - Возвращает значение целочисленного свойства для OpenCL-объекта или устройства
CLGetInfoString - Возвращает строковое значение свойства для OpenCL-объекта или устройства
CLContextCreate - Cоздает контекст OpenCL
CLContextFree - Удаляет контекст OpenCL
CLGetDeviceInfo - Получает свойство устройства из OpenCL драйвера
CLProgramCreate - Создает OpenCL программу из исходного кода
CLProgramFree - Удаляет OpenCL программу
CLKernelCreate - Создает функцию запуска OpenCL
CLKernelFree - Удаляет функцию запуска OpenCL
CLSetKernelArg - Выставляет параметр для функции OpenCL
CLSetKernelArgMem - Выставляет буфер OpenCL в качестве параметра функции OpenCL
CLSetKernelArgMemLocal - Задаёт локальный буфер в качестве аргумента kernel-функции
CLBufferCreate - Создает буфер OpenCL
CLBufferFree - Удаляет буфер OpenCL
CLBufferWrite - Записывает массив в буфер OpenCL и возвращает количество записанных элементов
CLBufferRead - Читает буфер OpenCL в массив и возвращает количество прочитанных элементов
CLExecute - Выполняет OpenCL программу
CLExecutionStatus - Возвращает состояние выполнения OpenCL программы

Рассмотрим пример, в котором данные считывается из дерева с 12 ветками. Затем производятся вычисления, заполнение гистограмм и результат записывается в новый root файл.

#include "TFile.h"
#include "TChain.h"
#include "TTree.h"
#include "TH1F.h"
#include "TCanvas.h"
#include "TStopwatch.h"
#include <iostream>
using namespace std;

int TreeNoPar(){
 TStopwatch timer;
 timer.Start();
 TFile *ff =new TFile("treeMy.root");
 TTree *treeMy=(TTree*)ff->Get("treeMy");
 const Int_t NumberBranch=12;
 Float_t* x=new Float_t[NumberBranch]; 
 Float_t* y=new Float_t[NumberBranch];
 TH1F *h[NumberBranch];
 for (Int_t i=0;i < NumberBranch; i++){
     treeMy->SetBranchAddress(Form("x%d",i),&x[i]);
     h[i]=new TH1F(Form("h%d",i),Form("title%d",i),150,-10,40);
 }
 Int_t i;
 for ( Int_t irow=0; irow < treeMy->GetEntries(); ++irow ) {
    treeMy->GetEntry(irow); 
    for (i=0;i <  NumberBranch; i++){
        y[i]=sqrt(10+x[i]*x[i]*x[i]*x[i])+sin(x[i])*cos(x[i])-tan(x[i]);
        h[i]->Fill(y[i]);
    }
 }
 TFile* f=new TFile("NoParallel.root","recreate");
 for(i=0;i < NumberBranch; i++)
     h[i]->Write();
 timer.Stop();
 cout<<"time: ";
 timer.Print("m");
 return 0;
 }

А теперь приведем программу, решающую ту же самую задачу, но уже с использованием OpenCL. Как отмечалось выше, программа в данном случае состоит из двух частей - хоста и ядра. В данном случаи имя первой, управляющей, части - TreeOpenCL.C. Текст хоста:

#include "stdio.h"
#include "math.h"
#include "stdlib.h"
#include "time.h"      
#include "CL/cl.h"
#include "TFile.h"
#include "TChain.h"
#include "TTree.h"
#include "TH1D.h"
#include "TCanvas.h"
#include "TStopwatch.h"
#include <iostream>
 
#define MAX_SOURCE_SIZE (0x100000)

int main(void)
{
        TFile *ff =new TFile("treeMy.root");
        TTree *treeMy=(TTree*)ff->Get("treeMy");
	clock_t begin = clock();
        int Number=treeMy->GetEntries();
        const int NumberBranch=12;
	const int N = Number*NumberBranch;
	float *X = (float*)malloc(sizeof(float)*N);
        float *x = (float*)malloc(sizeof(float)*NumberBranch);
        TH1D *hh[NumberBranch];
	for (int i=0;i < NumberBranch;i++){
	       treeMy->SetBranchAddress(Form("x%d",i),&x[i]);
	       hh[i]=new TH1D(Form("hh%d",i),Form("title%d",i),150,-10,140);
        }

         float **p=new float * [N];
         for(int i=0;i < Number;i++)
	      p[i]=new float[NumberBranch];

//       Считывание переменных дерева из root файла и запись их в массив        
	 for (int irow=0; irow < Number; ++irow){
	       treeMy->GetEntry(irow);
	 for(int j=0;j < NumberBranch;j++){
	   p[irow][j]=x[j];
	 }
	 }
	  int k=0;
	  for (int i=0;i < Number;i++){
	        for(int j=0;j < NumberBranch;j++){
	              X[k]=p[i][j];
                      k++;
	  }
	}

//	Загрузка исходного кода ядра в массив source_str
	FILE *fp;
	char *source_str;
	size_t source_size;
 
	fp = fopen("TreeOpenCL_kernel.cl", "r");
	if (!fp)
	{
		fprintf(stderr, "Failed to load kernel.\n");
		exit(1);
	}
	source_str = (char*)malloc(MAX_SOURCE_SIZE);
	source_size = fread( source_str, 1, MAX_SOURCE_SIZE, fp);
	fclose( fp );
 
//	Получение информации о платформе и устройствах
	cl_platform_id platform_id = NULL;
	cl_device_id device_id = NULL;   
	cl_uint ret_num_devices;
	cl_uint ret_num_platforms;
	clGetPlatformIDs(1, &platform_id, &ret_num_platforms);
	clGetDeviceIDs( platform_id, CL_DEVICE_TYPE_DEFAULT, 1, &device_id, &ret_num_devices);
 
//      Создание OpenCL контекста
	cl_context context = clCreateContext( NULL, 1, &device_id, NULL, NULL, NULL);
 
//	Создание очереди  команд
	cl_command_queue command_queue = clCreateCommandQueue(context, device_id, 0, NULL);
 
//	Создание буферов памяти на устройстве для каждого вектора
	cl_mem X_dev = clCreateBuffer(context, CL_MEM_READ_ONLY,  N * sizeof(float), NULL, NULL);
	cl_mem S_dev = clCreateBuffer(context, CL_MEM_WRITE_ONLY, N * sizeof(float), NULL, NULL);
 
//	Копирование значений переменной в   буфер памяти
	clEnqueueWriteBuffer(command_queue, X_dev, CL_TRUE, 0, N * sizeof(float), X, 0, NULL, NULL);
 
//	Создание программы из исходного кода ядра
	cl_program program = clCreateProgramWithSource(context, 1, (const char **)&source_str, (const size_t *)&source_size, NULL);
 
//	Построение программы     
	clBuildProgram(program, 1, &device_id, NULL, NULL, NULL);
 
//	Создание OpenCL ядра
	cl_kernel kernel = clCreateKernel(program, "Tree_OpenCL", NULL);
 
//	Установка аргументов ядра
	clSetKernelArg(kernel, 0, sizeof(cl_mem), (void *)&X_dev);
	clSetKernelArg(kernel, 1, sizeof(cl_mem), (void *)&S_dev);
	
//	Выполнение ядра OpenCL 
	size_t global_item_size = N; // Process the entire lists
	size_t local_item_size = 64; // Divide work items into groups of 64
	clEnqueueNDRangeKernel(command_queue, kernel, 1, NULL, &global_item_size, &local_item_size, 0, NULL, NULL);

//	Считывание из буфера памяти S_dev на устройстве локальной  переменной S
	float *S = (float*)malloc(sizeof(float)*N);
	clEnqueueReadBuffer(command_queue, S_dev, CL_TRUE, 0, N * sizeof(float), S, 0, NULL, NULL);
      
        TFile *f=new TFile("TreeOpenCL.root","recreate");
	int m=0;
	for (int irow=0; irow < Number; ++irow){
	   for(int j=0;j < NumberBranch;j++){
              hh[j]->Fill(S[m]);
	      m++;
	     }
       	}
        for (int i=0; i<12; ++i){
	hh[i]->Write();
	}
	f->Close();

	clock_t end = clock();
	double elapsed_secs = (double)(end - begin) / CLOCKS_PER_SEC;	
        printf("Time = %e\n", elapsed_secs);

//      Очищение памяти	
	clFlush(command_queue);
	clFinish(command_queue);
	clReleaseKernel(kernel);
	clReleaseProgram(program);
	clReleaseMemObject(X_dev);
	clReleaseMemObject(S_dev);
	clReleaseCommandQueue(command_queue);
	clReleaseContext(context);
	free(X); free(S); free(x);     
	return 0;
}      

Текст ядра - код, который будет выполняться непосредственно на GPU (имя файла - TreeOpenCL_kernel.cl):

float func(float x) { return sqrt(10+x*x*x*x)+sin(x)*cos(x)-tan(x);
       		      	}
__kernel void Tree_OpenCL( __global float *X, __global float *Y)
{
//	Get the index of the current element
	int gid = get_global_id(0);

//	Do the operation
	Y[gid] = func(X[gid]);
}

Приложения, написанные на языке OpenCL, могут быть скомпилированы и запущены на различных платформах. Для запуска на hydra.jinr.ru необходимо подключить одну из доступных на кластере версий CUDA, используя пакетный модуль MODULES: Также необходимо подключить одну из версий ROOT. Например так:

module add ROOT/v6-13-02-1
module add cuda/v9.1-1

Для компиляции макроса необходимо выполнить команду

g++ TreeOpenCL.C -lOpenCL `root-config --cflags --glibs`

При успешной компиляции образуется исполняемый бинарный файл. По умолчанию имя бинарного файла - a.out. Подробно о том, как запускать OpenCL - приложения на кластере hybrilit можно посмотреть здесь. Приведем рекомендуемый вид script-файла для запуска на GPU. Имя скрипта - script_gpu.

#!/bin/sh
# Set the partition where the job will run:
#SBATCH -p gpu
# Set the number of GPUs per node
#SBATCH --gres=gpu:1
# Set time of work: 
#SBATCH -t 60
# Submit a job for execution:
srun  ./a.out
# End script

Для запуска на CPU можно использовать скрипт script_cpu:

#!/bin/sh
#SBATCH -p cpu
#SBATCH -c 12
#SBATCH t 60
./a.out

Для запуска приложения используется следующая команда:

sbatch script_gpu

или

sbatch script_cpu

Программирование GPGPU (OpenCL и CUDA) подходит для такой обработки данных, в которой происходят интенсивные АРИФМЕТИЧЕСКИЕ вычисления. Например это актуально для (почти) всех задач линейной алгебры. Эти проблемы довольно легко распараллеливаются и подходят для решения на графических процессоров.

По этой ссылке можно найти результаты тестов, в которых сравнивается эффективность распараллеливания программ с помощью OpenCL и PROOF.