58 integer,
parameter :: rho_ndims = 3
63 real(
real64),
dimension(:,:,:),
allocatable,
intent(in) :: rho
64 real(
real64),
dimension(:),
allocatable,
intent(in) :: r
67 integer,
dimension(rho_ndims) :: rho_sizes, rho_dim_ids
68 integer :: r_size, r_dim_id
71 character(len=1),
dimension(rho_ndims) :: rho_dims=(/
"i",
"j",
"r"/)
72 character(len=3) :: r_dims =
"r_i"
75 character(len=*),
intent(in) :: filename
81 integer :: rho_id, r_id
84 rho_sizes = shape(rho)
88 call check(nf90_create(filename, nf90_clobber, file_id))
91 call check(nf90_put_att(file_id, nf90_global, &
93 call check(nf90_put_att(file_id, nf90_global, &
95 call check(nf90_put_att(file_id, nf90_global, &
97 call check(nf90_put_att(file_id, nf90_global, &
98 'Number of Species', setup%n_species))
99 call check(nf90_put_att(file_id, nf90_global, &
100 'Lattice Type', setup%lattice))
101 call check(nf90_put_att(file_id, nf90_global, &
102 'Interaction file', setup%interaction_file))
103 call check(nf90_put_att(file_id, nf90_global, &
104 'Concentrations', setup%species_concentrations))
105 call check(nf90_put_att(file_id, nf90_global, &
106 'Warren-Cowley Range', &
111 call check(nf90_def_dim(file_id, rho_dims(i), &
112 rho_sizes(i), rho_dim_ids(i)))
115 call check(nf90_def_var(file_id,
"rho data", nf90_double, &
116 rho_dim_ids, rho_id))
118 call check(nf90_def_dim(file_id, r_dims, r_size, r_dim_id))
120 call check(nf90_def_var(file_id,
"r data", nf90_double, &
124 call check(nf90_enddef(file_id))
127 call check(nf90_put_var(file_id, rho_id, rho))
128 call check(nf90_put_var(file_id, r_id, r))
131 call check(nf90_close(file_id))
152 integer,
parameter :: rho_ndims = 4
157 real(
real64),
dimension(:,:,:,:),
allocatable,
intent(in) :: rho
158 real(
real64),
dimension(:),
allocatable,
intent(in) :: r, t
159 real(
real64),
dimension(:),
allocatable,
intent(in) :: u_of_t
162 integer,
dimension(rho_ndims) :: rho_sizes, rho_dim_ids
163 integer :: r_size, r_dim_id, t_size, t_dim_id, u_size, u_dim_id
166 character(len=1),
dimension(rho_ndims) :: rho_dims=(/
"i",
"j",
"r",
"T"/)
167 character(len=3) :: r_dims =
"r_i"
168 character(len=3) :: t_dims =
"T_i"
169 character(len=3) :: u_dims =
"U_i"
172 character(len=*),
intent(in) :: filename
175 integer :: file_id, i
178 integer :: rho_id, r_id, t_id, u_id
181 rho_sizes = shape(rho)
184 u_size =
size(u_of_t)
187 call check(nf90_create(filename, nf90_clobber, file_id))
190 call check(nf90_put_att(file_id, nf90_global, &
192 call check(nf90_put_att(file_id, nf90_global, &
194 call check(nf90_put_att(file_id, nf90_global, &
196 call check(nf90_put_att(file_id, nf90_global, &
197 'Number of Species', setup%n_species))
198 call check(nf90_put_att(file_id, nf90_global, &
199 'Lattice Type', setup%lattice))
200 call check(nf90_put_att(file_id, nf90_global, &
201 'Interaction file', setup%interaction_file))
202 call check(nf90_put_att(file_id, nf90_global, &
203 'Concentrations', setup%species_concentrations))
204 call check(nf90_put_att(file_id, nf90_global, &
205 'Warren-Cowley Range', &
210 call check(nf90_def_dim(file_id, rho_dims(i), &
211 rho_sizes(i), rho_dim_ids(i)))
214 call check(nf90_def_var(file_id,
"rho data", nf90_double, &
215 rho_dim_ids, rho_id))
217 call check(nf90_def_dim(file_id, r_dims, r_size, r_dim_id))
219 call check(nf90_def_var(file_id,
"r data", nf90_double, &
222 call check(nf90_def_dim(file_id, t_dims, t_size, t_dim_id))
224 call check(nf90_def_var(file_id,
"T data", nf90_double, &
227 call check(nf90_def_dim(file_id, u_dims, u_size, u_dim_id))
229 call check(nf90_def_var(file_id,
"U data", nf90_double, &
233 call check(nf90_enddef(file_id))
236 call check(nf90_put_var(file_id, rho_id, rho))
237 call check(nf90_put_var(file_id, r_id, r))
238 call check(nf90_put_var(file_id, t_id, t))
239 call check(nf90_put_var(file_id, u_id, u_of_t))
242 call check(nf90_close(file_id))
262 integer,
parameter :: rho_ndims = 4
267 real(
real64),
dimension(:,:,:,:),
allocatable,
intent(in) :: rho
268 real(
real64),
dimension(:),
allocatable,
intent(in) :: r
269 real(
real64),
dimension(:),
allocatable,
intent(in) :: u
272 integer,
dimension(rho_ndims) :: rho_sizes, rho_dim_ids
273 integer :: r_size, r_dim_id, u_size, u_dim_id
276 character(len=1),
dimension(rho_ndims) :: rho_dims=(/
"i",
"j",
"r",
"U"/)
277 character(len=3) :: r_dims =
"r_i"
278 character(len=3) :: u_dims =
"U_i"
281 character(len=*),
intent(in) :: filename
284 integer :: file_id, i
287 integer :: rho_id, r_id, u_id
290 rho_sizes = shape(rho)
295 call check(nf90_create(filename, nf90_clobber, file_id))
298 call check(nf90_put_att(file_id, nf90_global, &
300 call check(nf90_put_att(file_id, nf90_global, &
302 call check(nf90_put_att(file_id, nf90_global, &
304 call check(nf90_put_att(file_id, nf90_global, &
305 'Number of Species', setup%n_species))
306 call check(nf90_put_att(file_id, nf90_global, &
307 'Lattice Type', setup%lattice))
308 call check(nf90_put_att(file_id, nf90_global, &
309 'Interaction file', setup%interaction_file))
310 call check(nf90_put_att(file_id, nf90_global, &
311 'Concentrations', setup%species_concentrations))
312 call check(nf90_put_att(file_id, nf90_global, &
313 'Warren-Cowley Range', &
318 call check(nf90_def_dim(file_id, rho_dims(i), &
319 rho_sizes(i), rho_dim_ids(i)))
322 call check(nf90_def_var(file_id,
"rho data", nf90_double, &
323 rho_dim_ids, rho_id))
325 call check(nf90_def_dim(file_id, r_dims, r_size, r_dim_id))
327 call check(nf90_def_var(file_id,
"r data", nf90_double, &
330 call check(nf90_def_dim(file_id, u_dims, u_size, u_dim_id))
332 call check(nf90_def_var(file_id,
"U data", nf90_double, &
336 call check(nf90_enddef(file_id))
339 call check(nf90_put_var(file_id, rho_id, rho))
340 call check(nf90_put_var(file_id, r_id, r))
341 call check(nf90_put_var(file_id, u_id, u))
344 call check(nf90_close(file_id))
364 integer,
parameter :: grid_ndims = 6
369 real(
real64),
dimension(:,:,:,:,:,:),
allocatable,
intent(in) :: order
370 real(
real64),
dimension(:),
allocatable,
intent(in) :: temperature
373 integer,
dimension(grid_ndims) :: grid_sizes, grid_dim_ids
374 integer :: temp_size, temp_dim_id
377 character(len=1),
dimension(grid_ndims) :: grid_dims=(/
"b",
"x",
"y" ,
"z",
"s",
"t"/)
378 character(len=4) :: temp_dims =
"temp"
381 character(len=*),
intent(in) :: filename
384 integer :: ierr, file_id, i
387 integer :: grid_id, temp_id
390 grid_sizes = shape(order)
391 temp_size =
size(temperature)
394 ierr = nf90_create(filename, nf90_clobber, file_id)
395 if (ierr /= nf90_noerr)
then
396 print*, trim(nf90_strerror(ierr))
401 call check(nf90_put_att(file_id, nf90_global, &
402 'N_Basis', setup%n_basis))
403 call check(nf90_put_att(file_id, nf90_global, &
405 call check(nf90_put_att(file_id, nf90_global, &
407 call check(nf90_put_att(file_id, nf90_global, &
409 call check(nf90_put_att(file_id, nf90_global, &
410 'Number of Species', setup%n_species))
411 call check(nf90_put_att(file_id, nf90_global, &
412 'Lattice Type', setup%lattice))
413 call check(nf90_put_att(file_id, nf90_global, &
414 'Interaction file', setup%interaction_file))
415 call check(nf90_put_att(file_id, nf90_global, &
416 'Concentrations', setup%species_concentrations))
417 call check(nf90_put_att(file_id, nf90_global, &
418 'Warren-Cowley Range', &
423 ierr = nf90_def_dim(file_id, grid_dims(i), grid_sizes(i), grid_dim_ids(i))
424 if (ierr /= nf90_noerr)
then
425 print*, trim(nf90_strerror(ierr))
430 ierr = nf90_def_var(file_id,
"grid data", nf90_double, grid_dim_ids, grid_id)
431 if (ierr /= nf90_noerr)
then
432 print*, trim(nf90_strerror(ierr))
436 if (
allocated(temperature))
then
438 ierr = nf90_def_dim(file_id, temp_dims, temp_size, temp_dim_id)
439 if (ierr /= nf90_noerr)
then
440 print*, trim(nf90_strerror(ierr))
444 ierr = nf90_def_var(file_id,
"temperature data", nf90_double, temp_dim_id, temp_id)
445 if (ierr /= nf90_noerr)
then
446 print*, trim(nf90_strerror(ierr))
452 ierr = nf90_enddef(file_id)
453 if (ierr /= nf90_noerr)
then
454 print*, trim(nf90_strerror(ierr))
459 ierr = nf90_put_var(file_id, grid_id, order)
460 if (ierr /= nf90_noerr)
then
461 print*, trim(nf90_strerror(ierr))
465 if (
allocated(temperature))
then
466 ierr = nf90_put_var(file_id, temp_id, temperature)
467 if (ierr /= nf90_noerr)
then
468 print*, trim(nf90_strerror(ierr))
474 ierr = nf90_close(file_id)
475 if (ierr /= nf90_noerr)
then
476 print*, trim(nf90_strerror(ierr))
497 integer,
parameter :: grid_ndims = 4
502 integer(array_int),
dimension(:,:,:,:),
allocatable,
intent(in) :: state
505 integer,
dimension(grid_ndims) :: grid_sizes, grid_dim_ids
508 character(len=1),
dimension(grid_ndims) :: &
509 grid_dims=(/
"b",
"x",
"y" ,
"z"/)
512 character(len=*),
intent(in) :: filename
515 integer :: ierr, file_id, i
521 grid_sizes = shape(state)
524 ierr = nf90_create(filename, nf90_clobber, file_id)
525 if (ierr /= nf90_noerr)
then
526 print*, trim(nf90_strerror(ierr))
531 call check(nf90_put_att(file_id, nf90_global, &
532 'N_basis', setup%n_basis))
534 call check(nf90_put_att(file_id, nf90_global, &
536 call check(nf90_put_att(file_id, nf90_global, &
538 call check(nf90_put_att(file_id, nf90_global, &
540 call check(nf90_put_att(file_id, nf90_global, &
541 'Number of Species', setup%n_species))
542 call check(nf90_put_att(file_id, nf90_global, &
543 'Lattice Type', setup%lattice))
544 call check(nf90_put_att(file_id, nf90_global, &
545 'Concentrations', setup%species_concentrations))
549 ierr = nf90_def_dim(file_id, grid_dims(i), grid_sizes(i), grid_dim_ids(i))
550 if (ierr /= nf90_noerr)
then
551 print*, trim(nf90_strerror(ierr))
556 ierr = nf90_def_var(file_id,
"configuration", nf90_short, grid_dim_ids, grid_id)
557 if (ierr /= nf90_noerr)
then
558 print*, trim(nf90_strerror(ierr))
563 ierr = nf90_enddef(file_id)
564 if (ierr /= nf90_noerr)
then
565 print*, trim(nf90_strerror(ierr))
570 ierr = nf90_put_var(file_id, grid_id, state)
571 if (ierr /= nf90_noerr)
then
572 print*, trim(nf90_strerror(ierr))
577 ierr = nf90_close(file_id)
578 if (ierr /= nf90_noerr)
then
579 print*, trim(nf90_strerror(ierr))
602 integer,
parameter :: grid_ndims = 4
607 integer(array_int),
dimension(:,:,:,:),
allocatable,
intent(in) :: states
608 real(
real64),
dimension(:),
allocatable,
intent(in) :: temperature
611 integer,
dimension(grid_ndims) :: grid_sizes, grid_dim_ids
612 integer :: temp_size, temp_dim_id
615 character(len=1),
dimension(grid_ndims) :: grid_dims=(/
"x",
"y" ,
"z",
"t"/)
616 character(len=4) :: temp_dims =
"temp"
619 character(len=*),
intent(in) :: filename
622 integer :: ierr, file_id, i
625 integer :: grid_id, temp_id
628 grid_sizes = shape(states)
629 temp_size =
size(temperature)
632 ierr = nf90_create(filename, nf90_clobber, file_id)
633 if (ierr /= nf90_noerr)
then
634 print*, trim(nf90_strerror(ierr))
639 call check(nf90_put_att(file_id, nf90_global, &
641 call check(nf90_put_att(file_id, nf90_global, &
643 call check(nf90_put_att(file_id, nf90_global, &
645 call check(nf90_put_att(file_id, nf90_global, &
646 'Number of Species', setup%n_species))
647 call check(nf90_put_att(file_id, nf90_global, &
648 'Lattice Type', setup%lattice))
649 call check(nf90_put_att(file_id, nf90_global, &
650 'Interaction file', setup%interaction_file))
651 call check(nf90_put_att(file_id, nf90_global, &
652 'Concentrations', setup%species_concentrations))
653 call check(nf90_put_att(file_id, nf90_global, &
654 'Warren-Cowley Range', &
659 ierr = nf90_def_dim(file_id, grid_dims(i), grid_sizes(i), grid_dim_ids(i))
660 if (ierr /= nf90_noerr)
then
661 print*, trim(nf90_strerror(ierr))
666 ierr = nf90_def_var(file_id,
"grid data", nf90_short, grid_dim_ids, grid_id)
667 if (ierr /= nf90_noerr)
then
668 print*, trim(nf90_strerror(ierr))
672 if (
allocated(temperature))
then
674 ierr = nf90_def_dim(file_id, temp_dims, temp_size, temp_dim_id)
675 if (ierr /= nf90_noerr)
then
676 print*, trim(nf90_strerror(ierr))
680 ierr = nf90_def_var(file_id,
"temperature data", nf90_double, temp_dim_id, temp_id)
681 if (ierr /= nf90_noerr)
then
682 print*, trim(nf90_strerror(ierr))
688 ierr = nf90_enddef(file_id)
689 if (ierr /= nf90_noerr)
then
690 print*, trim(nf90_strerror(ierr))
695 ierr = nf90_put_var(file_id, grid_id, states)
696 if (ierr /= nf90_noerr)
then
697 print*, trim(nf90_strerror(ierr))
701 if (
allocated(temperature))
then
702 ierr = nf90_put_var(file_id, temp_id, temperature)
703 if (ierr /= nf90_noerr)
then
704 print*, trim(nf90_strerror(ierr))
710 ierr = nf90_close(file_id)
711 if (ierr /= nf90_noerr)
then
712 print*, trim(nf90_strerror(ierr))
733 integer,
parameter :: grid_ndims = 1
736 real(
real64),
dimension(:),
allocatable,
intent(in) :: grid_data
738 integer,
dimension(grid_ndims) :: grid_sizes, grid_dim_ids
741 character(len=1),
dimension(grid_ndims) :: grid_dims=(/
"x"/)
744 character(len=*),
intent(in) :: filename
747 integer :: ierr, file_id, i
753 grid_sizes = shape(grid_data)
758 ierr = nf90_create(filename, nf90_clobber, file_id)
759 if (ierr /= nf90_noerr)
then
760 print*, trim(nf90_strerror(ierr))
768 ierr = nf90_def_dim(file_id, grid_dims(i), grid_sizes(i), grid_dim_ids(i))
769 if (ierr /= nf90_noerr)
then
770 print*, trim(nf90_strerror(ierr))
776 ierr = nf90_def_var(file_id,
"grid data", nf90_double, grid_dim_ids, grid_id)
777 if (ierr /= nf90_noerr)
then
778 print*, trim(nf90_strerror(ierr))
785 ierr = nf90_enddef(file_id)
786 if (ierr /= nf90_noerr)
then
787 print*, trim(nf90_strerror(ierr))
791 ierr = nf90_put_var(file_id, grid_id, grid_data)
792 if (ierr /= nf90_noerr)
then
793 print*, trim(nf90_strerror(ierr))
800 ierr = nf90_close(file_id)
801 if (ierr /= nf90_noerr)
then
802 print*, trim(nf90_strerror(ierr))
823 integer,
parameter :: grid_ndims = 2
826 real(
real64),
dimension(:,:),
allocatable,
intent(in) :: grid_data
828 integer,
dimension(grid_ndims) :: grid_sizes, grid_dim_ids
831 character(len=1),
dimension(grid_ndims) :: grid_dims=(/
"x",
"y" /)
834 character(len=*),
intent(in) :: filename
837 integer :: ierr, file_id, i
843 grid_sizes = shape(grid_data)
848 ierr = nf90_create(filename, nf90_clobber, file_id)
849 if (ierr /= nf90_noerr)
then
850 print*, trim(nf90_strerror(ierr))
858 ierr = nf90_def_dim(file_id, grid_dims(i), grid_sizes(i), grid_dim_ids(i))
859 if (ierr /= nf90_noerr)
then
860 print*, trim(nf90_strerror(ierr))
866 ierr = nf90_def_var(file_id,
"grid data", nf90_double, grid_dim_ids, grid_id)
867 if (ierr /= nf90_noerr)
then
868 print*, trim(nf90_strerror(ierr))
875 ierr = nf90_enddef(file_id)
876 if (ierr /= nf90_noerr)
then
877 print*, trim(nf90_strerror(ierr))
881 ierr = nf90_put_var(file_id, grid_id, grid_data)
882 if (ierr /= nf90_noerr)
then
883 print*, trim(nf90_strerror(ierr))
890 ierr = nf90_close(file_id)
891 if (ierr /= nf90_noerr)
then
892 print*, trim(nf90_strerror(ierr))
913 integer,
parameter :: grid_ndims = 5
916 real(
real64),
dimension(:,:,:,:,:),
allocatable,
intent(in) :: &
919 integer,
dimension(grid_ndims) :: grid_sizes, grid_dim_ids
922 character(len=2),
dimension(grid_ndims) :: &
923 grid_dims=(/
"s1",
"s2",
"xx",
"yy" ,
"zz"/)
926 character(len=*),
intent(in) :: filename
929 integer :: ierr, file_id, i
935 grid_sizes = shape(grid_data)
940 ierr = nf90_create(filename, nf90_clobber, file_id)
941 if (ierr /= nf90_noerr)
then
942 print*, trim(nf90_strerror(ierr))
950 ierr = nf90_def_dim(file_id, grid_dims(i), grid_sizes(i), grid_dim_ids(i))
951 if (ierr /= nf90_noerr)
then
952 print*, trim(nf90_strerror(ierr))
958 ierr = nf90_def_var(file_id,
"grid data", nf90_double, grid_dim_ids, grid_id)
959 if (ierr /= nf90_noerr)
then
960 print*, trim(nf90_strerror(ierr))
967 ierr = nf90_enddef(file_id)
968 if (ierr /= nf90_noerr)
then
969 print*, trim(nf90_strerror(ierr))
973 ierr = nf90_put_var(file_id, grid_id, grid_data)
974 if (ierr /= nf90_noerr)
then
975 print*, trim(nf90_strerror(ierr))
982 ierr = nf90_close(file_id)
983 if (ierr /= nf90_noerr)
then
984 print*, trim(nf90_strerror(ierr))
1005 integer,
parameter :: grid_ndims = 3
1008 real(
real64),
dimension(:,:,:),
allocatable,
intent(in) :: grid_data
1010 integer,
dimension(grid_ndims) :: grid_sizes, grid_dim_ids
1013 character(len=1),
dimension(grid_ndims) :: grid_dims=(/
"x",
"y" ,
"z"/)
1016 character(len=*),
intent(in) :: filename
1019 integer :: ierr, file_id, i
1025 grid_sizes = shape(grid_data)
1030 ierr = nf90_create(filename, nf90_clobber, file_id)
1031 if (ierr /= nf90_noerr)
then
1032 print*, trim(nf90_strerror(ierr))
1039 do i = 1, grid_ndims
1040 ierr = nf90_def_dim(file_id, grid_dims(i), grid_sizes(i), grid_dim_ids(i))
1041 if (ierr /= nf90_noerr)
then
1042 print*, trim(nf90_strerror(ierr))
1048 ierr = nf90_def_var(file_id,
"grid data", nf90_double, grid_dim_ids, grid_id)
1049 if (ierr /= nf90_noerr)
then
1050 print*, trim(nf90_strerror(ierr))
1057 ierr = nf90_enddef(file_id)
1058 if (ierr /= nf90_noerr)
then
1059 print*, trim(nf90_strerror(ierr))
1063 ierr = nf90_put_var(file_id, grid_id, grid_data)
1064 if (ierr /= nf90_noerr)
then
1065 print*, trim(nf90_strerror(ierr))
1072 ierr = nf90_close(file_id)
1073 if (ierr /= nf90_noerr)
then
1074 print*, trim(nf90_strerror(ierr))
1095 integer,
parameter :: grid_ndims = 4
1098 real(
real64),
dimension(:,:,:,:),
allocatable,
intent(in) :: grid_data
1100 integer,
dimension(grid_ndims) :: grid_sizes, grid_dim_ids
1103 character(len=1),
dimension(grid_ndims) :: grid_dims=(/
"x",
"y" ,
"z",
"t"/)
1106 character(len=*),
intent(in) :: filename
1109 integer :: ierr, file_id, i
1115 grid_sizes = shape(grid_data)
1120 ierr = nf90_create(filename, nf90_clobber, file_id)
1121 if (ierr /= nf90_noerr)
then
1122 print*, trim(nf90_strerror(ierr))
1129 do i = 1, grid_ndims
1130 ierr = nf90_def_dim(file_id, grid_dims(i), grid_sizes(i), grid_dim_ids(i))
1131 if (ierr /= nf90_noerr)
then
1132 print*, trim(nf90_strerror(ierr))
1138 ierr = nf90_def_var(file_id,
"grid data", nf90_double, grid_dim_ids, grid_id)
1139 if (ierr /= nf90_noerr)
then
1140 print*, trim(nf90_strerror(ierr))
1147 ierr = nf90_enddef(file_id)
1148 if (ierr /= nf90_noerr)
then
1149 print*, trim(nf90_strerror(ierr))
1153 ierr = nf90_put_var(file_id, grid_id, grid_data)
1154 if (ierr /= nf90_noerr)
then
1155 print*, trim(nf90_strerror(ierr))
1162 ierr = nf90_close(file_id)
1163 if (ierr /= nf90_noerr)
then
1164 print*, trim(nf90_strerror(ierr))
1187 real(
real64),
dimension(:),
allocatable,
intent(in) :: energies
1190 integer :: energy_size, energy_dim_id
1193 character(len=6) :: energy_dims =
"energy"
1196 integer :: energy_id
1198 integer,
parameter :: grid_ndims = 3
1201 integer(array_int),
dimension(:,:,:),
allocatable,
intent(in) :: grid_data
1203 integer,
dimension(grid_ndims) :: grid_sizes, grid_dim_ids
1206 character(len=1),
dimension(grid_ndims) :: grid_dims=(/
"x",
"y" ,
"z"/)
1209 character(len=*),
intent(in) :: filename
1212 integer :: ierr, file_id, i
1218 grid_sizes = shape(grid_data)
1219 energy_size =
size(energies)
1224 ierr = nf90_create(filename, nf90_clobber, file_id)
1225 if (ierr /= nf90_noerr)
then
1226 print*, trim(nf90_strerror(ierr))
1233 do i = 1, grid_ndims
1234 ierr = nf90_def_dim(file_id, grid_dims(i), grid_sizes(i), grid_dim_ids(i))
1235 if (ierr /= nf90_noerr)
then
1236 print*, trim(nf90_strerror(ierr))
1242 ierr = nf90_def_var(file_id,
"grid data", nf90_short, grid_dim_ids, grid_id)
1243 if (ierr /= nf90_noerr)
then
1244 print*, trim(nf90_strerror(ierr))
1249 ierr = nf90_def_dim(file_id, energy_dims, energy_size, energy_dim_id)
1250 if (ierr /= nf90_noerr)
then
1251 print*, trim(nf90_strerror(ierr))
1255 ierr = nf90_def_var(file_id,
"energy", nf90_double, energy_dim_id, energy_id)
1256 if (ierr /= nf90_noerr)
then
1257 print*, trim(nf90_strerror(ierr))
1264 ierr = nf90_enddef(file_id)
1265 if (ierr /= nf90_noerr)
then
1266 print*, trim(nf90_strerror(ierr))
1270 ierr = nf90_put_var(file_id, grid_id, grid_data)
1271 if (ierr /= nf90_noerr)
then
1272 print*, trim(nf90_strerror(ierr))
1276 ierr = nf90_put_var(file_id, energy_id, energies)
1277 if (ierr /= nf90_noerr)
then
1278 print*, trim(nf90_strerror(ierr))
1285 ierr = nf90_close(file_id)
1286 if (ierr /= nf90_noerr)
then
1287 print*, trim(nf90_strerror(ierr))
1304 integer,
intent ( in) :: stat
1307 if(stat /= nf90_noerr)
then
1308 print *, trim(nf90_strerror(stat))
1309 stop
"Stopped due to NetCDF error"
1311 end subroutine check
1325 character(len=*),
intent(in) :: filename
1326 character(len=*),
intent(in) :: varname
1327 real(
real64),
intent(out) :: array(:)
1328 integer :: ncid, varid, ierr
1331 ierr = nf90_open(filename, nf90_nowrite, ncid)
1332 if (ierr /= nf90_noerr)
then
1333 print *,
'Error opening file:', nf90_strerror(ierr)
1338 ierr = nf90_inq_varid(ncid, varname, varid)
1339 if (ierr /= nf90_noerr)
then
1340 print *,
'Error getting variable ID:', nf90_strerror(ierr)
1344 ierr = nf90_get_var(ncid, varid, array)
1345 if (ierr /= nf90_noerr)
then
1346 print *,
'Error reading variable:', nf90_strerror(ierr)
1350 ierr = nf90_close(ncid)
1351 if (ierr /= nf90_noerr)
then
1352 print *,
'Error closing file:', nf90_strerror(ierr)
1373 integer(array_int),
dimension(:,:,:,:),
allocatable,
intent(out) :: config
1376 character(len=*),
intent(in) :: filename
1379 character(len=20) :: lattice
1382 integer :: file_id, n_basis, n_1, n_2, n_3
1388 call check(nf90_open(filename, nf90_nowrite, file_id))
1391 call check(nf90_inquire_attribute(file_id, nf90_global,
"N_basis"))
1392 call check(nf90_inquire_attribute(file_id, nf90_global,
"N_1"))
1393 call check(nf90_inquire_attribute(file_id, nf90_global,
"N_2"))
1394 call check(nf90_inquire_attribute(file_id, nf90_global,
"N_3"))
1395 call check(nf90_inquire_attribute(file_id, nf90_global,
"Lattice Type"))
1398 call check(nf90_get_att(file_id, nf90_global,
"N_basis", n_basis))
1399 call check(nf90_get_att(file_id, nf90_global,
"N_1", n_1))
1400 call check(nf90_get_att(file_id, nf90_global,
"N_2", n_2))
1401 call check(nf90_get_att(file_id, nf90_global,
"N_3", n_3))
1402 call check(nf90_get_att(file_id, nf90_global,
"Lattice Type", lattice))
1405 if (trim(lattice) .ne. trim(setup%lattice))
then
1406 stop
"Lattice types do not match in ncdf_config_reader()"
1407 else if (n_basis .ne. setup%n_basis)
then
1408 stop
"n_basis values do not match in ncdf_config_reader()"
1409 else if (n_1 .ne. setup%n_1)
then
1410 stop
"n_1 values do not match in ncdf_config_reader()"
1411 else if (n_2 .ne. setup%n_2)
then
1412 stop
"n_2 values do not match in ncdf_config_reader()"
1413 else if (n_3 .ne. setup%n_3)
then
1414 stop
"n_3 values do not match in ncdf_config_reader()"
1418 call check(nf90_inq_varid(file_id,
"configuration", grid_id))
1421 allocate(config(n_basis, 2*n_1, 2*n_2, 2*n_3))
1424 call check(nf90_get_var(file_id, grid_id, config))
1427 call check(nf90_close(file_id))
1449 integer,
parameter :: rho_ndims = 4
1450 character(len=1),
dimension(rho_ndims) :: rho_dims=(/
"i",
"j",
"r",
"T"/)
1451 integer,
dimension(rho_ndims) :: rho_sizes, rho_dim_ids
1454 real(
real64),
dimension(:,:,:,:),
allocatable,
intent(out) :: asro
1455 real(
real64),
dimension(:),
allocatable,
intent(out) :: energy
1457 integer :: u_size, u_dim_id, u_id, file_id
1458 character(len=3) :: u_dims =
"U_i"
1460 integer :: n_steps, i
1462 character(len=20) :: lattice
1465 character(len=*),
intent(in) :: filename
1471 call check(nf90_open(filename, nf90_nowrite, file_id))
1474 call check(nf90_inquire_attribute(file_id, nf90_global,
"Lattice Type"))
1475 call check(nf90_get_att(file_id, nf90_global,
"Lattice Type", lattice))
1479 call check(nf90_inq_dimid(file_id, rho_dims(i), rho_dim_ids(i)))
1480 call check(nf90_inquire_dimension(file_id, rho_dim_ids(i), rho_dims(i), rho_sizes(i)))
1484 if (trim(lattice) .ne. trim(setup%lattice))
then
1485 stop
"Lattice types do not match in ncdf_radial_density_reader()"
1486 else if (rho_sizes(1) .ne. setup%n_species)
then
1487 stop
"n_species values do not match in ncdf_radial_density_reader()"
1488 else if (rho_sizes(2) .ne. setup%n_species)
then
1489 stop
"n_species values do not match in ncdf_radial_density_reader()"
1490 else if (rho_sizes(3) .ne. setup%wc_range)
then
1491 stop
"wc_range values do not match in ncdf_radial_density_reader()"
1492 else if (rho_sizes(4) .ne. n_steps)
then
1493 stop
"n_steps values do not match in ncdf_radial_density_reader()"
1497 allocate(asro(rho_sizes(1), rho_sizes(2), rho_sizes(3), rho_sizes(4)))
1500 call check(nf90_inq_varid(file_id,
"rho data", rho_id))
1503 call check(nf90_get_var(file_id, rho_id, asro))
1505 call check(nf90_inq_dimid(file_id, u_dims, u_dim_id))
1506 call check(nf90_inquire_dimension(file_id, u_dim_id, u_dims, u_size))
1508 allocate(energy(u_size))
1510 call check(nf90_inq_varid(file_id,
"U data", u_id))
1511 call check(nf90_get_var(file_id, u_id, energy))
1514 call check(nf90_close(file_id))
integer, parameter real64
Longer "double" (64 bit, approx -1.8e308 to 1.8e308 and covering values down to about 2e-308 magnitud...
subroutine, public read_1d_array(filename, varname, array)
Routine to read and parse bin edge NetCDF file.
subroutine, public ncdf_radial_density_writer_once(filename, rho, r, setup)
Routine to write radial density (calculated once) to file.
subroutine, public ncdf_grid_state_writer(filename, ierr, state, setup)
Routine to write current state of the grid to file.
subroutine, public ncdf_writer_2d(filename, ierr, grid_data)
Routine to write a 2D array to NetCDF file.
subroutine, public ncdf_config_reader(filename, config, setup)
Routine to read current state of the grid from file.
subroutine, public ncdf_radial_density_reader(filename, asro, energy, setup, n_steps)
Routine to read an array of calculated ASRO parameters and energies. (Mainly for testing purposes....
subroutine, public ncdf_writer_1d(filename, ierr, grid_data)
Routine to write a 1D array to NetCDF file.
subroutine, public ncdf_writer_3d(filename, ierr, grid_data)
Routine to write a 3D array to NetCDF file.
subroutine, public check(stat)
Routine to check NetCDF error codes.
subroutine, public ncdf_writer_3d_short(filename, ierr, grid_data, energies)
Routine to write a 3D array of shorts to NetCDF file.
subroutine, public ncdf_writer_4d(filename, ierr, grid_data)
Routine to write a 4D array to NetCDF file.
subroutine, public ncdf_order_writer(filename, ierr, order, temperature, setup)
Routine to write average atomic long-range order (ALRO) parameters as a function of temperature to fi...
subroutine, public ncdf_radial_density_writer(filename, rho, r, t, u_of_t, setup)
Routine to write radial densities as a function of temperature to file. Also writes internal energies...
subroutine, public ncdf_radial_density_writer_across_energy(filename, rho, r, u, setup)
Routine to write radial densities as a function of average energy to file (used in Wang-Landau).
subroutine, public ncdf_writer_5d(filename, ierr, grid_data)
Routine to write a 5D array to NetCDF file.
subroutine, public ncdf_grid_states_writer(filename, ierr, states, temperature, setup)
Routine to write list of states of the grid to file.
Derived type for parameters specifying general simulation parameters which are common to all sampling...