Сайт Георгия Мурого

Создание файла NetCDF на Fortran90

Первую программу для работы с NetCDF-файлами мне пришлось создавать методом проб и ошибок — готового примера для начинающих найти не удалось. Надеюсь, сей труд послужит таким примером.

Описывается создание файла со значениями в узлах регулярной сетки географических координат. Компилятор — GNU Fortran, операционная система — Linux Debian версии wheeze или новее.

Для упрощения создаваемый файл имеет версию 3 формата NetCDF, а программа не содержит дополнительных проверок (существование и формат входного файла, правильность указания даты и т.п.).

Исходный код форматирован для экрана шириной 80 символов.

Перед началом работы

Убедитесь, что на вашем компьютере установлены пакеты gfortran, libnetcdf-dev, libnetcdff-dev, netcdf-bin, cdo.

Обзаведитесь руководством по API (NetCDF Fortran90 Interface Guide): скачайте с веб-страницы https://www.unidata.ucar.edu/software/netcdf/netcdf-4/newdocs/index.html либо установите пакет libnetcdff-doc.

В качестве образца на будущее выберите готовый файл формата NetCDF. Например, скачайте с FTP-сайта ESRL (ftp://ftp.cdc.noaa.gpv) подходящий файл реанализа со среднемесячными данными и в программе CDO выберите данные за один месяц одного года. Желательно использовать гидрометеорологический параметр с данными на единственном горизонте (приземное атмосферное давление, температура поверхности океана, атмосферное давление на уровне моря и т. п.).

На основе получившегося файла (input.nc) сформируйте две шпаргалки:

CDL-файл с заголовками и значениями координат (header.cdl) формируется командой:

$ ncdump -c input.nc > header.cdl

Исходный текст программы на Fortran77 формируется командами:

$ ncdump input.nc > temp.cdl
$ ncgen -lf77 temp.cdl > create.for

При чтении дальнейшего текста, а также при разработке собственной программы, сверяйтесь с шпаргалками и руководством по API.

Постановка задачи

Вы работаете в институте «СуперНИИ».

Ваше подразделение получает от научного учреждения «DataSource» данные объективного анализа гидрометеорологического параметра «kq на уровне моря» и добавляет их в институтский массив научных данных, официально называющийся «Database of SuperNII».

Данные поставляются в текстовых файлах (ASCII text) в узлах регулярной сетки с шагом 5° по широте и долготе, по одному файлу на каждый срок наблюдения (осреднения). Файл — таблица из 19 строк и 72 колонок. Каждая строка таблицы содержит значения на одной широте, каждый столбец — на одной долготе.

Широты уменьшаются с шагом 5° от Северного полюса (первая строка) до экватора (последняя строка). Долготы увеличиваются в восточном направлении с шагом 5° от Гринвичского меридиана (первая колонка). Последняя колонка (355° в. д.) соответствует 5° з. д.

Единица измерения kq — безразмерный коэффициент. Пропуску в данных соответствует значение, равное -999.9.

Институт переводит свой массив данных в формат NetCDF. Необходимо разработать программу для перевода данных за один срок в файл формата NetCDF. Полученный файл будет объединяться с файлами за другие сроки.

Методика формирования массива «Database of SuperNII» не публиковалась, так что библиографические ссылки на массив отсутствуют.

Формирование тестовых данных

Файл тестовых данных (test_in.dat) сформируем в следующей программе:

program testdata
implicit none
integer :: ilat, ilon
open (1, file = 'test_in.dat')
write (1, 800) ( 90., ilon = 0, 71 )
do ilat =  17, 0, -1
  write (1, fmt = 800 ) (float (ilat * 5) + float (ilon * 5) / 1000., &
   ilon = 0, 71 )
end do
close ( 1 )
800 format (72f8.3)
stop 'Файл сформирован'
end program testdata

Проектирование файла NetCDF

Строго говоря, требования к структуре NetCDF-файлов должны содержаться в документации на «Database of SuperNII». Однако для примера разработаем структуру с начала.

Файл создадим в соответствии с конвенцией «Climate and Forecast Metadata Conventions» (CF) версии 1.6 (http://cfconventions.org).

Глобальные атрибуты:

Атрибут Назначение Значение
Conventions Конвенция (типовой или местный стандарт структуры файла) CF-1.6
title Название массива данных Database of SuperNII
institution Организация, подготовившая данные DataSource
source Источник данных objective analyses of hydrometeorological fields
references Библиографические ссылки на массив данных None
history История изменений NetCDF-файла Формируется автоматически

В файле будут 4 переменные: три координаты — широта (lat), долгота (lon), время (time) — и гидрометеорологический параметр (kq).

Широта и долгота — одномерные массивы типа REAL*4. Размерность массива: для широты — 19, для долготы — 72. Значения широты уменьшаются с 90 до 0 с шагом 5, долготы — возрастают от 0 до 355 с шагом 5.

Атрибуты широты:

Атрибут Назначение Значение
standard_name Краткое название latitude
long_name Полное название latitude
units Единица измерения degrees_north

Атрибуты долготы:

Атрибут Назначение Значение
standard_name Краткое название longitude
long_name Полное название longitude
units Единица измерения degrees_east

Время — одномерный массив типа REAL*8 единичной длины.

Чтобы объединение получившегося NetCDF-файла с файлами за другие сроки было простым, размерность переменной времени — неограниченной длины (nf90_unlimited).

Для простоты использована единица времени «day as %Y%m%d.%f», при которой значение для полуночи 2 января 2018 года равно 20180102,0. В документации на CDO это называется абсолютной шкалой времени. Однако для работы с реальными данными такая единица непригодна — из-за особенности обработки чисел в памяти ЭВМ при операциях с большими числами неминуемо происходит потеря точности и, в результате, искажение значений. Необходимо использовать относительную шкалу времени (время, прошедшее с определённой даты).

Атрибуты времени:

Атрибут Назначение Значение
standard_name Краткое название time
long_name Полное название time
units Единица измерения day as %Y%m%d.%f

Переменная kq — трехмерный массив типа REAL*4 ( time, lat, lon ). Размерность по широте — 19, по долготе — 72, по времени — единичный. Единица измерения — безразмерный коэффициент (counts). Признак отсутствующего (забракованного) значения — значение, равное -999.9.

Атрибут Назначение Значение
name Краткое название kq
long_name Полное название Sea Level KQ
units Единица измерения counts
_FillValue Значение — признак отбраковки или отсутствия данных -999.9

Вспомогательные процедуры

Программы для формата NetCDF изобилуют элементарными операциями, исходный текст безобразно распухает, его трудно понимать.

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

Подпрограмму check_err возьмём из сформированного нами файла create.for и адаптируем под Fortran90:

subroutine check_err ( stat )
use netcdf
implicit none
integer :: stat
if ( stat /= NF90_NOERR ) then
  print *, nf90_strerror ( stat )
  stop
endif
return
end subroutine check_err 

Подпрограмма load_input_data считывает файл с исходными данными и возвращает массив c значениями параметра в узлах регулярной сетки:

subroutine load_input_data ( file, array )
implicit none
character ( len = * ), intent ( in ) :: file
real, dimension ( 19, 72 ), intent ( out ) :: array
integer :: ilat
open ( 1, file = file, action='read', status = 'old' )
do ilat =  19, 1, -1
  read (1, '(72f8.3)', err=99, end=99) array ( ilat, : )
end do
close ( 1 )
return
99 stop 'Ошибка! Не удалось считать файл' // trim( file )
end subroutine load_input_data

Подпрограмма fill_coord заполняет одномерный массив значениями арифметической прогрессии.

subroutine fill_coord ( array, n, start, step )
implicit none
real, dimension ( n ), intent ( out ) :: array
integer, intent ( in ) :: n
real, intent ( in ) :: start, step
integer :: i
do i = 1, n
 array ( i ) = start + step * float ( i - 1 )
end do
return
end subroutine fill_coord

Функция define_time вычисляет значение переменной времени для даты, указанной в командной строке программы:

real ( kind = 8 ) function define_time ( year, month, day )
implicit none
integer, intent ( in ) :: year, month, day
define_time = real ( year, 8 ) * 1d4 + real ( month, 8 ) * 1d2 + real &
  ( day, 8 )
return
end function define_time

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

subroutine define_history ( history )
implicit none
character ( len = * ), intent ( out ) :: history
character ( len = 128 ) :: param
integer :: numparam, iparam
call ctime ( time8 (), history )
history = trim ( history ) // ' : '
numparam = command_argument_count () 
do iparam = 0, numparam
 call get_command_argument ( iparam, param )
 history = trim ( history ) // ' ' // trim ( param )
end do
return
end subroutine define_history

Алгоритм работы программы

Опишем алгоритм нашей программы. Для сокращения текста однотипные места в исходном коде опускаем.

Вначале из командной строки считываем имена входного и выходного файлов и дату:

call get_command_argument ( 3, chyear )
read ( chyear, '(i4)' ) year

Считываем исходный файл:

call load_input_data ( infile, kq )

Формируем значение атрибута истории:

call define_history ( history )

Вычисляем значение времени:

outtime = define_time ( year, month, day )

Входной файл содержит двумерный массив kq ( 19, 72 ), а в NetCDF-файл выводится трёхмерный kq_out ( 72, 19, 1 ). Формируем выходной массив:

do ilon = 1, 72
 do ilat = 19, 1, -1
  kq_out ( ilon, 20 - ilat , 1 ) = kq ( ilat, ilon )
 end do
end do

Формируем массивы с значениями координат:

call fill_coord ( lat, n = 19, start = 90., step = -5. )

Все выводимые данные готовы, создаём NetCDF-файл:

status = nf90_create ( path = trim ( outfile ), cmode = nf90_noclobber, &
  ncid =  ncid )
call check_err ( status )

Объявляем размерности переменных координат:

status = nf90_def_dim ( ncid = ncid, name = 'lon', len = 72, &
  dimid = lon_dimid )
call check_err ( status )

Объявляем переменные координат:

status = nf90_def_var ( ncid = ncid, name = 'lon', xtype = nf90_float,  &
 dimids = lon_dimid, varid =  lon_varid )
call check_err ( status )

Формируем массив с размерностями переменной гидрометпараметра:

hmparam_dimids ( 1 )  = lon_dimid
hmparam_dimids ( 2 )  = lat_dimid
hmparam_dimids ( 3 )  = time_dimid

Объявляем переменную гидрометпараметра:

status = nf90_def_var ( ncid = ncid, name = 'kq', xtype = nf90_float,  &
 dimids = hmparam_dimids, varid = hmparam_varid )
call check_err ( status )

Записываем атрибуты координат:

status = nf90_put_att ( ncid = ncid, varid = lat_varid, name = &
  'standard_name', values = 'latitude' )
call check_err ( status )

Теперь — атрибуты параметра:

status = nf90_put_att ( ncid = ncid, varid = hmparam_varid, name =&
  'long_name', values = 'Sea Level KQ' )
call check_err ( status )

И, наконец, глобальные атрибуты:

status = nf90_put_att ( ncid = ncid, varid = nf90_global, name = &
  'history', values = history )
call check_err ( status )

Записывать атрибуты можно в любое время после объявления их переменных. В нашей программы атрибуты записываются, когда структура NetCDF-файла уже готова, но можно было записать глобальные атрибуты сразу после создания файла, а атрибуты координат и гидрометпараметра — сразу после создания соответствующих переменных.

Атрибуты записываются в файл в порядке их объявления. Этот порядок, теоретически, определяется конвенцией, но главный критерий — ваше удобство. Например, я привык, что в данных реанализа ESRL история — второй глобальный атрибут, и в своих NetCDF-файлах также ставлю историю на второе место.

Кстати, ничего не мешает «застолбить» место для нужного атрибута, объявив его хоть с пустым значением, а позже присвоить нужное значение повторным вызовом функции nf90_put_att.

Когда объявление структуры файла и всех атрибутов закончено, переходим в режим записи значений:

status = nf90_enddef ( ncid )
call check_err ( status )

Записываем в файл значения координат и параметра:

status = nf90_put_var ( ncid = ncid, varid = time_varid, values = outtime )
call check_err ( status )
status = nf90_put_var ( ncid = ncid, varid = hmparam_varid, values = &
  kq_out )
call check_err ( status )

Осталось только закрыть файл:

status = nf90_close ( ncid )
call check_err ( status )

Завершаем программу.

Исходный текст программы

program netcdf_create
use netcdf
implicit none
!
! Рабочие переменные программы
character ( len = 80 ) :: infile, outfile
integer :: status
!
! Переменные, относящиеся к широте
real, dimension ( 19 ) :: lat
integer :: lat_dimid, lat_varid
integer:: ilat
!
! Переменные, относящиеся к долготе
real, dimension ( 72 ) :: lon
integer :: lon_dimid, lon_varid
integer :: ilon
!
! Переменные, относящиеся к времени
character ( len = 4 ) :: chyear
character ( len = 2 ) :: chmonth, chday
real ( kind = 8 ) :: outtime 
integer :: time_dimid, time_varid
integer :: year, month, day
!
! Переменные, относящиеся к гидрометпараметру
real, dimension ( 19, 72 ) :: kq
real, dimension ( 72, 19, 1 ) :: kq_out
integer, dimension ( 3 ) :: hmparam_dimids
integer :: hmparam_varid
!
! Переменные, относящиеся к глобальным атрибутам
! и файлу NetCDF в целом
character ( len = 128 ) :: history
integer :: ncid
!
! Обрабатываем командную строку
call get_command_argument ( 1, infile )
call get_command_argument ( 2, outfile )
call get_command_argument ( 3, chyear )
call get_command_argument ( 4, chmonth )
call get_command_argument ( 5, chday )
read ( chyear, '(i4)' ) year
read ( chmonth, '(i2)' ) month
read ( chday, '(i2)' ) day
!
! Считываем входной файл
call  load_input_data ( infile, kq )
!
! Вычисляем значение времени
outtime = define_time ( year, month, day )
!
! Формируем значение истории
call define_history ( history )
!
! Формируем выходной массив:
do ilon = 1, 72
 do ilat = 19, 1, -1
  kq_out ( ilon, 20 - ilat , 1 ) = kq ( ilat, ilon )
 end do
end do
!
! Формируем массивы координат:
call fill_coord ( lat, n = 19, start = 90., step = -5. )
call fill_coord ( lon, n = 72, start = 0., step = 5. )
!
! Создаём файл NetCDF
status = nf90_create ( path = trim ( outfile ), cmode = nf90_noclobber, &
  ncid =  ncid )
call check_err ( status )
!
! Объявляем размерности координат
status = nf90_def_dim ( ncid = ncid, name = 'lon', len = 72,  dimid = &
  lon_dimid )
call check_err ( status )
status = nf90_def_dim ( ncid = ncid, name = 'lat', len = 19,  dimid = &
  lat_dimid )
call check_err ( status )
status = nf90_def_dim ( ncid = ncid, name = 'time', len = nf90_unlimited, &
  dimid = time_dimid )
call check_err ( status )
!
! Объявляем переменные координат
status = nf90_def_var ( ncid = ncid, name = 'lon', xtype = nf90_float,  &
 dimids = lon_dimid, varid =  lon_varid )
call check_err ( status )
status = nf90_def_var ( ncid = ncid, name = 'lat', xtype = nf90_float,  &
 dimids = lat_dimid, varid =  lat_varid )
call check_err ( status )
status = nf90_def_var ( ncid = ncid, name = 'time', xtype = nf90_double,  &
 dimids = time_dimid, varid = time_varid )
call check_err ( status )
!
! Объявляем переменную параметра
hmparam_dimids ( 1 )  = lon_dimid
hmparam_dimids ( 2 )  = lat_dimid
hmparam_dimids ( 3 )  = time_dimid
status = nf90_def_var ( ncid = ncid, name = 'kq', xtype = nf90_float,  &
 dimids = hmparam_dimids, varid = hmparam_varid )
call check_err ( status )
!
! Записываем атрибуты широты
status = nf90_put_att ( ncid = ncid, varid = lat_varid, name = &
  'standard_name', values = 'latitude' )
call check_err ( status )
status = nf90_put_att ( ncid = ncid, varid = lat_varid, name = 'long_name', &
  values = 'latitude' )
call check_err ( status )
status = nf90_put_att ( ncid = ncid, varid = lat_varid, name = 'units', &
  values = 'degrees_north' )
call check_err ( status )
!
! Записываем атрибуты долготы
status = nf90_put_att ( ncid = ncid, varid = lon_varid, name = &
  'standard_name', values = 'longitude' )
call check_err ( status )
status = nf90_put_att ( ncid = ncid, varid = lon_varid, name = 'long_name',&
  values = 'longitude' )
call check_err ( status )
status = nf90_put_att ( ncid = ncid, varid = lon_varid, name = 'units', &
  values = 'degrees_east' )
call check_err ( status )
!
! Записываем атрибуты времени
status = nf90_put_att ( ncid = ncid, varid = time_varid, name =&
  'standard_name', values = 'time' )
call check_err ( status )
status = nf90_put_att ( ncid = ncid, varid = time_varid, name = 'long_name',&
  values = 'time' )
call check_err ( status )
status = nf90_put_att ( ncid = ncid, varid = time_varid, name = 'units', &
  values = 'day as %Y%m%d.%f' )
call check_err ( status )
!
! Записываем атрибуты гидрометпараметра
status = nf90_put_att ( ncid = ncid, varid = hmparam_varid, name = &
  'standard_name', values = 'kq' )
call check_err ( status )
status = nf90_put_att ( ncid = ncid, varid = hmparam_varid, name =&
  'long_name', values = 'Sea Level KQ' )
call check_err ( status )
status = nf90_put_att ( ncid = ncid, varid = hmparam_varid, name = 'units', &
  values = 'counts' )
call check_err ( status )
status = nf90_put_att ( ncid = ncid, varid = hmparam_varid, name = &
 '_FillValue', values = -999.9 )
call check_err ( status )
!
! Записываем глобальные атрибуты.
status = nf90_put_att ( ncid = ncid, varid = nf90_global, name =&
  'Conventions', values = 'CF-1.6' )
call check_err ( status )
status = nf90_put_att ( ncid = ncid, varid = nf90_global, name = 'history', &
  values = history )
call check_err ( status )
status = nf90_put_att ( ncid = ncid, varid = nf90_global, name = 'title', &
  values = 'Database of SuperNII' )
call check_err ( status )
status = nf90_put_att ( ncid = ncid, varid = nf90_global, name =&
 'institution', values = 'DataSource' )
call check_err ( status )
status = nf90_put_att ( ncid = ncid, varid = nf90_global, name = 'source', &
  values = 'objective analyses of hydrometeorological fields' )
call check_err ( status )
status = nf90_put_att ( ncid = ncid, varid = nf90_global, name = & 
  'references', values = 'None' )
call check_err ( status )
!
! Переводим файл NetCDF в режим записи значений:
status = nf90_enddef ( ncid )
call check_err ( status )
!
! Записываем в файл значения координат и параметра:
status = nf90_put_var ( ncid = ncid, varid = lat_varid, values = lat )
call check_err ( status )
status = nf90_put_var ( ncid = ncid, varid = lon_varid, values = lon )
call check_err ( status )
status = nf90_put_var ( ncid = ncid, varid = time_varid, values = outtime )
call check_err ( status )
status = nf90_put_var ( ncid = ncid, varid = hmparam_varid, values = kq_out )
call check_err ( status )
!
! Закрываем выходной файл:
status = nf90_close ( ncid )
call check_err ( status )
!
stop 'Программа доработала до конца'
!
contains 
!
subroutine check_err ( stat )
use netcdf
implicit none
integer :: stat
if (stat /= NF90_NOERR) then
  print *, nf90_strerror ( stat )
  stop
endif
return
end subroutine check_err 
!
subroutine load_input_data ( file, array )
implicit none
character ( len = * ), intent ( in ) :: file
real, dimension (19, 72), intent ( out ) :: array
integer :: ilat
open (1, file =  file, action='read', status = 'old' )
do ilat =  19, 1, -1
  read (1, '(72f8.3)', err = 99, end = 99 ) array ( ilat, : )
end do
close (1)
return
99 stop 'Ошибка! Не удалось считать файл' // trim( file )
end subroutine load_input_data
!
subroutine fill_coord ( array, n, start, step )
implicit none
real, dimension ( n ), intent ( out ) :: array
integer, intent ( in ) :: n
real, intent ( in ) :: start, step
integer :: i
do i = 1, n
 array ( i ) = start + step * float ( i - 1 )
end do
return
end subroutine fill_coord
!
real ( kind = 8 ) function define_time ( year, month, day )
implicit none
integer, intent ( in ) :: year, month, day
define_time = real ( year, 8 ) * 1d4 + real ( month, 8 ) * 1d2 + &
  real ( day, 8 )
return
end function define_time
!
subroutine define_history ( history )
implicit none
character ( len = * ), intent ( out ) :: history
character ( len = 128 ) :: param
integer :: numparam, iparam
call ctime ( time8 (), history )
history = trim ( history ) // ' : '
numparam = command_argument_count () 
do iparam = 0, numparam
 call get_command_argument ( iparam, param )
 history = trim ( history ) // ' ' // trim ( param )
end do
return
end subroutine define_history
!
end program netcdf_create
! 

Компиляция и отладка

Программу компилируем командой:

$ gfortran netcdf_create.f90 -o netcdf_create -l netcdff -I/usr/include

Обращаю внимание, что программа линкуется с библиотекой libnetcdff, а не с libnetcdf.

Выполняем программу:

$ ./netcdf_create test_in.dat testout.nc 2018 2 15

Проверяем заголовки и значения координат в получившемся файле:

$ ncdump -c testout.nc > testout.cdl

Получившийся CDL-файл приведён ниже. Ваш файл будет отличаться от него только датой и временем выполнения программы в атрибуте истории.

netcdf test_out {
dimensions:
	lon = 72 ;
	lat = 19 ;
	time = UNLIMITED ; // (1 currently)
variables:
	float lon(lon) ;
		lon:standard_name = "longitude" ;
		lon:long_name = "longitude" ;
		lon:units = "degrees_east" ;
	float lat(lat) ;
		lat:standard_name = "latitude" ;
		lat:long_name = "latitude" ;
		lat:units = "degrees_north" ;
	double time(time) ;
		time:standard_name = "time" ;
		time:long_name = "time" ;
		time:units = "day as %Y%m%d.%f" ;
	float kq(time, lat, lon) ;
		kq:standard_name = "kq" ;
		kq:long_name = "Sea Level KQ" ;
		kq:units = "counts" ;
		kq:_FillValue = -999.9f ;

// global attributes:
		:Conventions = "CF-1.6" ;
		:history = "Fri Oct 12 00:55:47 2018 : ./netcdf_create 
test_in.dat test_out.nc 2018 2 15" ;
		:title = "Database of SuperNII" ;
		:institution = "DataSource" ;
		:source = "objective analyses of hydrometeorological fields" ;
		:references = "None" ;
data:

 lon = 0, 5, 10, 15, 20, 25, 30, 35, 40, 45, 50, 55, 60, 65, 70, 75, 80, 85, 
    90, 95, 100, 105, 110, 115, 120, 125, 130, 135, 140, 145, 150, 155, 160, 
    165, 170, 175, 180, 185, 190, 195, 200, 205, 210, 215, 220, 225, 230, 
    235, 240, 245, 250, 255, 260, 265, 270, 275, 280, 285, 290, 295, 300, 
    305, 310, 315, 320, 325, 330, 335, 340, 345, 350, 355 ;

 lat = 90, 85, 80, 75, 70, 65, 60, 55, 50, 45, 40, 35, 30, 25, 20, 15, 10, 5, 0 ;

 time = 20180215 ;
}

Выберем точку (допустим, 35° с. ш. 125° в. д.). Значение kq в этой точке должно быть 35 + 125/1000 = 35,125. Проверяем:

$ cdo sellonlatbox,125,125,35,35 test_out.nc test_35_125.nc
$ cdo outputtab,date,value test_35_125.nc 

Вывод программы должен быть таким:

#      date    value 
 2018-02-15   35.125 

Дополнение про многомерные массивы

Легко заметить, что размерности массива гидрометеорологического параметра в программе и в готовом NetCDF-файле различаются.

В программе в массиве размерностей переменной kq (hmparam_dimids) координаты перечислены в следующем порядке: долгота, широта, время:

hmparam_dimids ( 1 )  = lon_dimid
hmparam_dimids ( 2 )  = lat_dimid
hmparam_dimids ( 3 )  = time_dimid

Однако в NetCDF-файле у этой переменной порядок следования координат — обратный:

	float kq(time, lat, lon) ;

Это — одно из проявлений специфики многомерных переменных в NetCDF-файлах.

При разработке программ для чтения и записи NetCDF-файлов более сложной структуры (несколько многомерных переменных, связанных между собой координат и т. п.) необходимо учитывать, что на структуру многомерных переменных может влиять порядок, в котором объявлялись предыдущие переменные и их размерности.

Таким образом, при освоении формата NetCDF целесообразно рассчитывать программы только на структуру файлов, с которыми работаете постоянно. К программам для более сложных файлов, а, тем более, универсальным, переходите после получения необходимого опыта.

Заключение

При подготовке данного текста я воспроизвёл разработку программы для решения реальной задачи, максимально эту программу упростив.

Я постарался подробно описать вопросы, возникшие передо мной в начале работы с NetCDF: проектирование структуры файла, дополнительные процедуры, алгоритм работы программы, процесс отладки программы и некоторые, не очевидные для новичков, вещи.

Надеюсь, этот текст окажется вам полезен.

Программа для считывания получившегося NetCDF-файла описана здесь.

Автор: Георгий Мурый, site@moury.ru, страница создана 15 октября 2018 года, изменена 2 марта 2023 года.