BraWl
Loading...
Searching...
No Matches
netcdf_io.f90
Go to the documentation of this file.
1
15
16 use kinds
17 use constants
18 use netcdf
20
21 implicit none
22
23 private
24
40 check
41
42 contains
43
56 subroutine ncdf_radial_density_writer_once(filename, rho, r, setup)
57
58 integer, parameter :: rho_ndims = 3
59
60 type(run_params), intent(in) :: setup
61
62 ! Data to write to file
63 real(real64), dimension(:,:,:), allocatable, intent(in) :: rho
64 real(real64), dimension(:), allocatable, intent(in) :: r
65
66 ! Number of dimensions of my grid data
67 integer, dimension(rho_ndims) :: rho_sizes, rho_dim_ids
68 integer :: r_size, r_dim_id
69
70 ! Names of my dimensions
71 character(len=1), dimension(rho_ndims) :: rho_dims=(/"i", "j", "r"/)
72 character(len=3) :: r_dims = "r_i"
73
74 ! Filename to which to write
75 character(len=*), intent(in) :: filename
76
77 ! Variables used in writing process
78 integer :: file_id, i
79
80 ! Ids for variables
81 integer :: rho_id, r_id
82
83 ! Get the sizes of my incoming arrays
84 rho_sizes = shape(rho)
85 r_size = size(r)
86
87 ! Create the file
88 call check(nf90_create(filename, nf90_clobber, file_id))
89
90 ! Add information about global runtime data
91 call check(nf90_put_att(file_id, nf90_global, &
92 'N_1', setup%n_1))
93 call check(nf90_put_att(file_id, nf90_global, &
94 'N_2', setup%n_2))
95 call check(nf90_put_att(file_id, nf90_global, &
96 'N_3', setup%n_3))
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', &
107 setup%wc_range))
108
109 ! Define the 3D variables and dimensions
110 do i = 1, rho_ndims
111 call check(nf90_def_dim(file_id, rho_dims(i), &
112 rho_sizes(i), rho_dim_ids(i)))
113 end do
114
115 call check(nf90_def_var(file_id, "rho data", nf90_double, &
116 rho_dim_ids, rho_id))
117
118 call check(nf90_def_dim(file_id, r_dims, r_size, r_dim_id))
119
120 call check(nf90_def_var(file_id, "r data", nf90_double, &
121 r_dim_id, r_id))
122
123 ! Finish defining metadata
124 call check(nf90_enddef(file_id))
125
126 ! Dump the variables to file
127 call check(nf90_put_var(file_id, rho_id, rho))
128 call check(nf90_put_var(file_id, r_id, r))
129
130 ! Close the file
131 call check(nf90_close(file_id))
132
134
150 subroutine ncdf_radial_density_writer(filename, rho, r, T, U_of_T, setup)
151
152 integer, parameter :: rho_ndims = 4
153
154 type(run_params), intent(in) :: setup
155
156 ! Data to write to file
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
160
161 ! Number of dimensions of my grid data
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
164
165 ! Names of my dimensions
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"
170
171 ! Filename to which to write
172 character(len=*), intent(in) :: filename
173
174 ! Variables used in writing process
175 integer :: file_id, i
176
177 ! Ids for variables
178 integer :: rho_id, r_id, t_id, u_id
179
180 ! Get the sizes of my incoming arrays
181 rho_sizes = shape(rho)
182 r_size = size(r)
183 t_size = size(t)
184 u_size = size(u_of_t)
185
186 ! Create the file
187 call check(nf90_create(filename, nf90_clobber, file_id))
188
189 ! Add information about global runtime data
190 call check(nf90_put_att(file_id, nf90_global, &
191 'N_1', setup%n_1))
192 call check(nf90_put_att(file_id, nf90_global, &
193 'N_2', setup%n_2))
194 call check(nf90_put_att(file_id, nf90_global, &
195 'N_3', setup%n_3))
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', &
206 setup%wc_range))
207
208 ! Define the 3D variables and dimensions
209 do i = 1, rho_ndims
210 call check(nf90_def_dim(file_id, rho_dims(i), &
211 rho_sizes(i), rho_dim_ids(i)))
212 end do
213
214 call check(nf90_def_var(file_id, "rho data", nf90_double, &
215 rho_dim_ids, rho_id))
216
217 call check(nf90_def_dim(file_id, r_dims, r_size, r_dim_id))
218
219 call check(nf90_def_var(file_id, "r data", nf90_double, &
220 r_dim_id, r_id))
221
222 call check(nf90_def_dim(file_id, t_dims, t_size, t_dim_id))
223
224 call check(nf90_def_var(file_id, "T data", nf90_double, &
225 t_dim_id, t_id))
226
227 call check(nf90_def_dim(file_id, u_dims, u_size, u_dim_id))
228
229 call check(nf90_def_var(file_id, "U data", nf90_double, &
230 u_dim_id, u_id))
231
232 ! Finish defining metadata
233 call check(nf90_enddef(file_id))
234
235 ! Dump the variables to file
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))
240
241 ! Close the file
242 call check(nf90_close(file_id))
243
244 end subroutine ncdf_radial_density_writer
245
260 subroutine ncdf_radial_density_writer_across_energy(filename, rho, r, U, setup)
261
262 integer, parameter :: rho_ndims = 4
263
264 type(run_params), intent(in) :: setup
265
266 ! Data to write to file
267 real(real64), dimension(:,:,:,:), allocatable, intent(in) :: rho
268 real(real64), dimension(:), allocatable, intent(in) :: r
269 real(real64), dimension(:), allocatable, intent(in) :: u
270
271 ! Number of dimensions of my grid data
272 integer, dimension(rho_ndims) :: rho_sizes, rho_dim_ids
273 integer :: r_size, r_dim_id, u_size, u_dim_id
274
275 ! Names of my dimensions
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"
279
280 ! Filename to which to write
281 character(len=*), intent(in) :: filename
282
283 ! Variables used in writing process
284 integer :: file_id, i
285
286 ! Ids for variables
287 integer :: rho_id, r_id, u_id
288
289 ! Get the sizes of my incoming arrays
290 rho_sizes = shape(rho)
291 r_size = size(r)
292 u_size = size(u)
293
294 ! Create the file
295 call check(nf90_create(filename, nf90_clobber, file_id))
296
297 ! Add information about global runtime data
298 call check(nf90_put_att(file_id, nf90_global, &
299 'N_1', setup%n_1))
300 call check(nf90_put_att(file_id, nf90_global, &
301 'N_2', setup%n_2))
302 call check(nf90_put_att(file_id, nf90_global, &
303 'N_3', setup%n_3))
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', &
314 setup%wc_range))
315
316 ! Define the 3D variables and dimensions
317 do i = 1, rho_ndims
318 call check(nf90_def_dim(file_id, rho_dims(i), &
319 rho_sizes(i), rho_dim_ids(i)))
320 end do
321
322 call check(nf90_def_var(file_id, "rho data", nf90_double, &
323 rho_dim_ids, rho_id))
324
325 call check(nf90_def_dim(file_id, r_dims, r_size, r_dim_id))
326
327 call check(nf90_def_var(file_id, "r data", nf90_double, &
328 r_dim_id, r_id))
329
330 call check(nf90_def_dim(file_id, u_dims, u_size, u_dim_id))
331
332 call check(nf90_def_var(file_id, "U data", nf90_double, &
333 u_dim_id, u_id))
334
335 ! Finish defining metadata
336 call check(nf90_enddef(file_id))
337
338 ! Dump the variables to file
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))
342
343 ! Close the file
344 call check(nf90_close(file_id))
345
347
362 subroutine ncdf_order_writer(filename, ierr, order, temperature, setup)
363
364 integer, parameter :: grid_ndims = 6
365
366 type(run_params), intent(in) :: setup
367
368 ! Data to write to file
369 real(real64), dimension(:,:,:,:,:,:), allocatable, intent(in) :: order
370 real(real64), dimension(:), allocatable, intent(in) :: temperature
371
372 ! Number of dimensions of my grid data
373 integer, dimension(grid_ndims) :: grid_sizes, grid_dim_ids
374 integer :: temp_size, temp_dim_id
375
376 ! Names of my dimensions
377 character(len=1), dimension(grid_ndims) :: grid_dims=(/"b", "x", "y" , "z", "s", "t"/)
378 character(len=4) :: temp_dims = "temp"
379
380 ! Filename to which to write
381 character(len=*), intent(in) :: filename
382
383 ! Variables used in writing process
384 integer :: ierr, file_id, i
385
386 ! Ids for variables
387 integer :: grid_id, temp_id
388
389 ! Get the sizes of my incoming arrays
390 grid_sizes = shape(order)
391 temp_size = size(temperature)
392
393 ! Create the file, overwriting if it exists !
394 ierr = nf90_create(filename, nf90_clobber, file_id)
395 if (ierr /= nf90_noerr) then
396 print*, trim(nf90_strerror(ierr))
397 return
398 end if
399
400 ! Add information about global runtime data
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, &
404 'N_1', setup%n_1))
405 call check(nf90_put_att(file_id, nf90_global, &
406 'N_2', setup%n_2))
407 call check(nf90_put_att(file_id, nf90_global, &
408 'N_3', setup%n_3))
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', &
419 setup%wc_range))
420
421 ! Define the variables and dimensions
422 do i = 1, grid_ndims
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))
426 return
427 end if
428 end do
429
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))
433 return
434 end if
435
436 if (allocated(temperature)) then
437 ! Define the temperature variable and dimensions !
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))
441 return
442 end if
443
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))
447 return
448 end if
449 end if
450
451 ! Finish defining metadata !
452 ierr = nf90_enddef(file_id)
453 if (ierr /= nf90_noerr) then
454 print*, trim(nf90_strerror(ierr))
455 return
456 end if
457
458 ! Dump the variables to file !
459 ierr = nf90_put_var(file_id, grid_id, order)
460 if (ierr /= nf90_noerr) then
461 print*, trim(nf90_strerror(ierr))
462 return
463 end if
464
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))
469 return
470 end if
471 end if
472
473 ! Close the file !
474 ierr = nf90_close(file_id)
475 if (ierr /= nf90_noerr) then
476 print*, trim(nf90_strerror(ierr))
477 return
478 end if
479
480 end subroutine ncdf_order_writer
481
495 subroutine ncdf_grid_state_writer(filename, ierr, state, setup)
496
497 integer, parameter :: grid_ndims = 4
498
499 type(run_params), intent(in) :: setup
500
501 ! Data to write to file
502 integer(array_int), dimension(:,:,:,:), allocatable, intent(in) :: state
503
504 ! Number of dimensions of my grid data
505 integer, dimension(grid_ndims) :: grid_sizes, grid_dim_ids
506
507 ! Names of my dimensions
508 character(len=1), dimension(grid_ndims) :: &
509 grid_dims=(/"b", "x", "y" ,"z"/)
510
511 ! Filename to which to write
512 character(len=*), intent(in) :: filename
513
514 ! Variables used in writing process
515 integer :: ierr, file_id, i
516
517 ! Ids for variables
518 integer :: grid_id
519
520 ! Get the sizes of my incoming arrays
521 grid_sizes = shape(state)
522
523 ! Create the file, overwriting if it exists !
524 ierr = nf90_create(filename, nf90_clobber, file_id)
525 if (ierr /= nf90_noerr) then
526 print*, trim(nf90_strerror(ierr))
527 return
528 end if
529
531 call check(nf90_put_att(file_id, nf90_global, &
532 'N_basis', setup%n_basis))
533 ! Add information about global runtime data
534 call check(nf90_put_att(file_id, nf90_global, &
535 'N_1', setup%n_1))
536 call check(nf90_put_att(file_id, nf90_global, &
537 'N_2', setup%n_2))
538 call check(nf90_put_att(file_id, nf90_global, &
539 'N_3', setup%n_3))
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))
546
547 ! Define the 4D variables and dimensions !
548 do i = 1, grid_ndims
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))
552 return
553 end if
554 end do
555
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))
559 return
560 end if
561
562 ! Finish defining metadata !
563 ierr = nf90_enddef(file_id)
564 if (ierr /= nf90_noerr) then
565 print*, trim(nf90_strerror(ierr))
566 return
567 end if
568
569 ! Dump the variables to file !
570 ierr = nf90_put_var(file_id, grid_id, state)
571 if (ierr /= nf90_noerr) then
572 print*, trim(nf90_strerror(ierr))
573 return
574 end if
575
576 ! Close the file !
577 ierr = nf90_close(file_id)
578 if (ierr /= nf90_noerr) then
579 print*, trim(nf90_strerror(ierr))
580 return
581 end if
582
583 end subroutine ncdf_grid_state_writer
584
600 subroutine ncdf_grid_states_writer(filename, ierr, states, temperature, setup)
601
602 integer, parameter :: grid_ndims = 4
603
604 type(run_params), intent(in) :: setup
605
606 ! Data to write to file
607 integer(array_int), dimension(:,:,:,:), allocatable, intent(in) :: states
608 real(real64), dimension(:), allocatable, intent(in) :: temperature
609
610 ! Number of dimensions of my grid data
611 integer, dimension(grid_ndims) :: grid_sizes, grid_dim_ids
612 integer :: temp_size, temp_dim_id
613
614 ! Names of my dimensions
615 character(len=1), dimension(grid_ndims) :: grid_dims=(/"x", "y" , "z", "t"/)
616 character(len=4) :: temp_dims = "temp"
617
618 ! Filename to which to write
619 character(len=*), intent(in) :: filename
620
621 ! Variables used in writing process
622 integer :: ierr, file_id, i
623
624 ! Ids for variables
625 integer :: grid_id, temp_id
626
627 ! Get the sizes of my incoming arrays
628 grid_sizes = shape(states)
629 temp_size = size(temperature)
630
631 ! Create the file, overwriting if it exists !
632 ierr = nf90_create(filename, nf90_clobber, file_id)
633 if (ierr /= nf90_noerr) then
634 print*, trim(nf90_strerror(ierr))
635 return
636 end if
637
638 ! Add information about global runtime data
639 call check(nf90_put_att(file_id, nf90_global, &
640 'N_1', setup%n_1))
641 call check(nf90_put_att(file_id, nf90_global, &
642 'N_2', setup%n_2))
643 call check(nf90_put_att(file_id, nf90_global, &
644 'N_3', setup%n_3))
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', &
655 setup%wc_range))
656
657 ! Define the 4D variables and dimensions !
658 do i = 1, grid_ndims
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))
662 return
663 end if
664 end do
665
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))
669 return
670 end if
671
672 if (allocated(temperature)) then
673 ! Define the temperature variable and dimensions !
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))
677 return
678 end if
679
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))
683 return
684 end if
685 end if
686
687 ! Finish defining metadata !
688 ierr = nf90_enddef(file_id)
689 if (ierr /= nf90_noerr) then
690 print*, trim(nf90_strerror(ierr))
691 return
692 end if
693
694 ! Dump the variables to file !
695 ierr = nf90_put_var(file_id, grid_id, states)
696 if (ierr /= nf90_noerr) then
697 print*, trim(nf90_strerror(ierr))
698 return
699 end if
700
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))
705 return
706 end if
707 end if
708
709 ! Close the file !
710 ierr = nf90_close(file_id)
711 if (ierr /= nf90_noerr) then
712 print*, trim(nf90_strerror(ierr))
713 return
714 end if
715
716 end subroutine ncdf_grid_states_writer
717
731 subroutine ncdf_writer_1d(filename, ierr, grid_data)
732
733 integer, parameter :: grid_ndims = 1
734
735 ! Data to write to file
736 real(real64), dimension(:), allocatable, intent(in) :: grid_data
737 ! Number of dimensions of my rho data
738 integer, dimension(grid_ndims) :: grid_sizes, grid_dim_ids
739
740 ! Names of my dimensions
741 character(len=1), dimension(grid_ndims) :: grid_dims=(/"x"/)
742
743 ! Filename to which to write
744 character(len=*), intent(in) :: filename
745
746 ! Variables used in writing process
747 integer :: ierr, file_id, i
748
749 ! Ids for variables
750 integer :: grid_id
751
752 ! Get the sizes of my incoming arrays
753 grid_sizes = shape(grid_data)
754
755 !-------------------------------------------!
756 ! Create the file, overwriting if it exists !
757 !-------------------------------------------!
758 ierr = nf90_create(filename, nf90_clobber, file_id)
759 if (ierr /= nf90_noerr) then
760 print*, trim(nf90_strerror(ierr))
761 return
762 end if
763
764 !----------------------------------------!
765 ! Define the 2D variables and dimensions !
766 !----------------------------------------!
767 do i = 1, grid_ndims
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))
771 return
772 end if
773 end do
774
775 ! grid_data
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))
779 return
780 end if
781
782 !--------------------------!
783 ! Finish defining metadata !
784 !--------------------------!
785 ierr = nf90_enddef(file_id)
786 if (ierr /= nf90_noerr) then
787 print*, trim(nf90_strerror(ierr))
788 return
789 end if
790
791 ierr = nf90_put_var(file_id, grid_id, grid_data)
792 if (ierr /= nf90_noerr) then
793 print*, trim(nf90_strerror(ierr))
794 return
795 end if
796
797 !----------------!
798 ! Close the file !
799 !----------------!
800 ierr = nf90_close(file_id)
801 if (ierr /= nf90_noerr) then
802 print*, trim(nf90_strerror(ierr))
803 return
804 end if
805
806 end subroutine ncdf_writer_1d
807
821 subroutine ncdf_writer_2d(filename, ierr, grid_data)
822
823 integer, parameter :: grid_ndims = 2
824
825 ! Data to write to file
826 real(real64), dimension(:,:), allocatable, intent(in) :: grid_data
827 ! Number of dimensions of my rho data
828 integer, dimension(grid_ndims) :: grid_sizes, grid_dim_ids
829
830 ! Names of my dimensions
831 character(len=1), dimension(grid_ndims) :: grid_dims=(/"x", "y" /)
832
833 ! Filename to which to write
834 character(len=*), intent(in) :: filename
835
836 ! Variables used in writing process
837 integer :: ierr, file_id, i
838
839 ! Ids for variables
840 integer :: grid_id
841
842 ! Get the sizes of my incoming arrays
843 grid_sizes = shape(grid_data)
844
845 !-------------------------------------------!
846 ! Create the file, overwriting if it exists !
847 !-------------------------------------------!
848 ierr = nf90_create(filename, nf90_clobber, file_id)
849 if (ierr /= nf90_noerr) then
850 print*, trim(nf90_strerror(ierr))
851 return
852 end if
853
854 !----------------------------------------!
855 ! Define the 2D variables and dimensions !
856 !----------------------------------------!
857 do i = 1, grid_ndims
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))
861 return
862 end if
863 end do
864
865 ! grid_data
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))
869 return
870 end if
871
872 !--------------------------!
873 ! Finish defining metadata !
874 !--------------------------!
875 ierr = nf90_enddef(file_id)
876 if (ierr /= nf90_noerr) then
877 print*, trim(nf90_strerror(ierr))
878 return
879 end if
880
881 ierr = nf90_put_var(file_id, grid_id, grid_data)
882 if (ierr /= nf90_noerr) then
883 print*, trim(nf90_strerror(ierr))
884 return
885 end if
886
887 !----------------!
888 ! Close the file !
889 !----------------!
890 ierr = nf90_close(file_id)
891 if (ierr /= nf90_noerr) then
892 print*, trim(nf90_strerror(ierr))
893 return
894 end if
895
896 end subroutine ncdf_writer_2d
897
911 subroutine ncdf_writer_5d(filename, ierr, grid_data)
912
913 integer, parameter :: grid_ndims = 5
914
915 ! Data to write to file
916 real(real64), dimension(:,:,:,:,:), allocatable, intent(in) :: &
917 grid_data
918 ! Number of dimensions of my rho data
919 integer, dimension(grid_ndims) :: grid_sizes, grid_dim_ids
920
921 ! Names of my dimensions
922 character(len=2), dimension(grid_ndims) :: &
923 grid_dims=(/"s1", "s2", "xx", "yy" , "zz"/)
924
925 ! Filename to which to write
926 character(len=*), intent(in) :: filename
927
928 ! Variables used in writing process
929 integer :: ierr, file_id, i
930
931 ! Ids for variables
932 integer :: grid_id
933
934 ! Get the sizes of my incoming arrays
935 grid_sizes = shape(grid_data)
936
937 !-------------------------------------------!
938 ! Create the file, overwriting if it exists !
939 !-------------------------------------------!
940 ierr = nf90_create(filename, nf90_clobber, file_id)
941 if (ierr /= nf90_noerr) then
942 print*, trim(nf90_strerror(ierr))
943 return
944 end if
945
946 !----------------------------------------!
947 ! Define the 3D variables and dimensions !
948 !----------------------------------------!
949 do i = 1, grid_ndims
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))
953 return
954 end if
955 end do
956
957 ! grid_data
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))
961 return
962 end if
963
964 !--------------------------!
965 ! Finish defining metadata !
966 !--------------------------!
967 ierr = nf90_enddef(file_id)
968 if (ierr /= nf90_noerr) then
969 print*, trim(nf90_strerror(ierr))
970 return
971 end if
972
973 ierr = nf90_put_var(file_id, grid_id, grid_data)
974 if (ierr /= nf90_noerr) then
975 print*, trim(nf90_strerror(ierr))
976 return
977 end if
978
979 !----------------!
980 ! Close the file !
981 !----------------!
982 ierr = nf90_close(file_id)
983 if (ierr /= nf90_noerr) then
984 print*, trim(nf90_strerror(ierr))
985 return
986 end if
987
988 end subroutine ncdf_writer_5d
989
1003 subroutine ncdf_writer_3d(filename, ierr, grid_data)
1004
1005 integer, parameter :: grid_ndims = 3
1006
1007 ! Data to write to file
1008 real(real64), dimension(:,:,:), allocatable, intent(in) :: grid_data
1009 ! Number of dimensions of my rho data
1010 integer, dimension(grid_ndims) :: grid_sizes, grid_dim_ids
1011
1012 ! Names of my dimensions
1013 character(len=1), dimension(grid_ndims) :: grid_dims=(/"x", "y" , "z"/)
1014
1015 ! Filename to which to write
1016 character(len=*), intent(in) :: filename
1017
1018 ! Variables used in writing process
1019 integer :: ierr, file_id, i
1020
1021 ! Ids for variables
1022 integer :: grid_id
1023
1024 ! Get the sizes of my incoming arrays
1025 grid_sizes = shape(grid_data)
1026
1027 !-------------------------------------------!
1028 ! Create the file, overwriting if it exists !
1029 !-------------------------------------------!
1030 ierr = nf90_create(filename, nf90_clobber, file_id)
1031 if (ierr /= nf90_noerr) then
1032 print*, trim(nf90_strerror(ierr))
1033 return
1034 end if
1035
1036 !----------------------------------------!
1037 ! Define the 3D variables and dimensions !
1038 !----------------------------------------!
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))
1043 return
1044 end if
1045 end do
1046
1047 ! grid_data
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))
1051 return
1052 end if
1053
1054 !--------------------------!
1055 ! Finish defining metadata !
1056 !--------------------------!
1057 ierr = nf90_enddef(file_id)
1058 if (ierr /= nf90_noerr) then
1059 print*, trim(nf90_strerror(ierr))
1060 return
1061 end if
1062
1063 ierr = nf90_put_var(file_id, grid_id, grid_data)
1064 if (ierr /= nf90_noerr) then
1065 print*, trim(nf90_strerror(ierr))
1066 return
1067 end if
1068
1069 !----------------!
1070 ! Close the file !
1071 !----------------!
1072 ierr = nf90_close(file_id)
1073 if (ierr /= nf90_noerr) then
1074 print*, trim(nf90_strerror(ierr))
1075 return
1076 end if
1077
1078 end subroutine ncdf_writer_3d
1079
1093 subroutine ncdf_writer_4d(filename, ierr, grid_data)
1094
1095 integer, parameter :: grid_ndims = 4
1096
1097 ! Data to write to file
1098 real(real64), dimension(:,:,:,:), allocatable, intent(in) :: grid_data
1099 ! Number of dimensions of my rho data
1100 integer, dimension(grid_ndims) :: grid_sizes, grid_dim_ids
1101
1102 ! Names of my dimensions
1103 character(len=1), dimension(grid_ndims) :: grid_dims=(/"x", "y" , "z", "t"/)
1104
1105 ! Filename to which to write
1106 character(len=*), intent(in) :: filename
1107
1108 ! Variables used in writing process
1109 integer :: ierr, file_id, i
1110
1111 ! Ids for variables
1112 integer :: grid_id
1113
1114 ! Get the sizes of my incoming arrays
1115 grid_sizes = shape(grid_data)
1116
1117 !-------------------------------------------!
1118 ! Create the file, overwriting if it exists !
1119 !-------------------------------------------!
1120 ierr = nf90_create(filename, nf90_clobber, file_id)
1121 if (ierr /= nf90_noerr) then
1122 print*, trim(nf90_strerror(ierr))
1123 return
1124 end if
1125
1126 !----------------------------------------!
1127 ! Define the 3D variables and dimensions !
1128 !----------------------------------------!
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))
1133 return
1134 end if
1135 end do
1136
1137 ! grid_data
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))
1141 return
1142 end if
1143
1144 !--------------------------!
1145 ! Finish defining metadata !
1146 !--------------------------!
1147 ierr = nf90_enddef(file_id)
1148 if (ierr /= nf90_noerr) then
1149 print*, trim(nf90_strerror(ierr))
1150 return
1151 end if
1152
1153 ierr = nf90_put_var(file_id, grid_id, grid_data)
1154 if (ierr /= nf90_noerr) then
1155 print*, trim(nf90_strerror(ierr))
1156 return
1157 end if
1158
1159 !----------------!
1160 ! Close the file !
1161 !----------------!
1162 ierr = nf90_close(file_id)
1163 if (ierr /= nf90_noerr) then
1164 print*, trim(nf90_strerror(ierr))
1165 return
1166 end if
1167
1168 end subroutine ncdf_writer_4d
1169
1185 subroutine ncdf_writer_3d_short(filename, ierr, grid_data, energies)
1186
1187 real(real64), dimension(:), allocatable, intent(in) :: energies
1188
1189 ! Number of dimensions of my grid data
1190 integer :: energy_size, energy_dim_id
1191
1192 ! Names of my dimensions
1193 character(len=6) :: energy_dims = "energy"
1194
1195 ! Ids for variables
1196 integer :: energy_id
1197
1198 integer, parameter :: grid_ndims = 3
1199
1200 ! Data to write to file
1201 integer(array_int), dimension(:,:,:), allocatable, intent(in) :: grid_data
1202 ! Number of dimensions of my rho data
1203 integer, dimension(grid_ndims) :: grid_sizes, grid_dim_ids
1204
1205 ! Names of my dimensions
1206 character(len=1), dimension(grid_ndims) :: grid_dims=(/"x", "y" , "z"/)
1207
1208 ! Filename to which to write
1209 character(len=*), intent(in) :: filename
1210
1211 ! Variables used in writing process
1212 integer :: ierr, file_id, i
1213
1214 ! Ids for variables
1215 integer :: grid_id
1216
1217 ! Get the sizes of my incoming arrays
1218 grid_sizes = shape(grid_data)
1219 energy_size = size(energies)
1220
1221 !-------------------------------------------!
1222 ! Create the file, overwriting if it exists !
1223 !-------------------------------------------!
1224 ierr = nf90_create(filename, nf90_clobber, file_id)
1225 if (ierr /= nf90_noerr) then
1226 print*, trim(nf90_strerror(ierr))
1227 return
1228 end if
1229
1230 !----------------------------------------!
1231 ! Define the 3D variables and dimensions !
1232 !----------------------------------------!
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))
1237 return
1238 end if
1239 end do
1240
1241 ! grid_data
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))
1245 return
1246 end if
1247
1248 ! Define the energy variable and dimensions !
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))
1252 return
1253 end if
1254
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))
1258 return
1259 end if
1260
1261 !--------------------------!
1262 ! Finish defining metadata !
1263 !--------------------------!
1264 ierr = nf90_enddef(file_id)
1265 if (ierr /= nf90_noerr) then
1266 print*, trim(nf90_strerror(ierr))
1267 return
1268 end if
1269
1270 ierr = nf90_put_var(file_id, grid_id, grid_data)
1271 if (ierr /= nf90_noerr) then
1272 print*, trim(nf90_strerror(ierr))
1273 return
1274 end if
1275
1276 ierr = nf90_put_var(file_id, energy_id, energies)
1277 if (ierr /= nf90_noerr) then
1278 print*, trim(nf90_strerror(ierr))
1279 return
1280 end if
1281
1282 !----------------!
1283 ! Close the file !
1284 !----------------!
1285 ierr = nf90_close(file_id)
1286 if (ierr /= nf90_noerr) then
1287 print*, trim(nf90_strerror(ierr))
1288 return
1289 end if
1290
1291 end subroutine ncdf_writer_3d_short
1292
1302 subroutine check(stat)
1303
1304 integer, intent ( in) :: stat
1305
1307 if(stat /= nf90_noerr) then
1308 print *, trim(nf90_strerror(stat))
1309 stop "Stopped due to NetCDF error"
1310 end if
1311 end subroutine check
1312
1324 subroutine read_1d_array(filename, varname, array)
1325 character(len=*), intent(in) :: filename
1326 character(len=*), intent(in) :: varname
1327 real(real64), intent(out) :: array(:)
1328 integer :: ncid, varid, ierr
1329
1330 ! Open the NetCDF file
1331 ierr = nf90_open(filename, nf90_nowrite, ncid)
1332 if (ierr /= nf90_noerr) then
1333 print *, 'Error opening file:', nf90_strerror(ierr)
1334 return
1335 end if
1336
1337 ! Get the variable ID
1338 ierr = nf90_inq_varid(ncid, varname, varid)
1339 if (ierr /= nf90_noerr) then
1340 print *, 'Error getting variable ID:', nf90_strerror(ierr)
1341 end if
1342
1343 ! Read the data into the array
1344 ierr = nf90_get_var(ncid, varid, array)
1345 if (ierr /= nf90_noerr) then
1346 print *, 'Error reading variable:', nf90_strerror(ierr)
1347 end if
1348
1349 ! Close the NetCDF file
1350 ierr = nf90_close(ncid)
1351 if (ierr /= nf90_noerr) then
1352 print *, 'Error closing file:', nf90_strerror(ierr)
1353 end if
1354
1355 end subroutine read_1d_array
1356
1368 subroutine ncdf_config_reader(filename, config, setup)
1369
1370 type(run_params), intent(in) :: setup
1371
1372 ! Data to read from file
1373 integer(array_int), dimension(:,:,:,:), allocatable, intent(out) :: config
1374
1375 ! Filename from to which to read
1376 character(len=*), intent(in) :: filename
1377
1378 ! Filename from to which to read
1379 character(len=20) :: lattice
1380
1381 ! Variables used in writing process
1382 integer :: file_id, n_basis, n_1, n_2, n_3
1383
1384 ! Ids for variables
1385 integer :: grid_id
1386
1387 ! Open the file for read-in
1388 call check(nf90_open(filename, nf90_nowrite, file_id))
1389
1390 ! Get the IDs of the relevant attributes
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"))
1396
1397 ! Get the values of the relevant attributes
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))
1403
1404 ! Check that these match the current simulation
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()"
1415 end if
1416
1417 ! Check that the variable ID exists
1418 call check(nf90_inq_varid(file_id, "configuration", grid_id))
1419
1420 ! Allocate my grid for reading in
1421 allocate(config(n_basis, 2*n_1, 2*n_2, 2*n_3))
1422
1423 ! Read in
1424 call check(nf90_get_var(file_id, grid_id, config))
1425
1426 ! Close the file
1427 call check(nf90_close(file_id))
1428
1429 end subroutine ncdf_config_reader
1430
1443 subroutine ncdf_radial_density_reader(filename, asro, energy, setup, &
1444 n_steps)
1445
1446 type(run_params), intent(in) :: setup
1447
1448 ! Names of my dimensions
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
1452
1453 ! Data to read from file
1454 real(real64), dimension(:,:,:,:), allocatable, intent(out) :: asro
1455 real(real64), dimension(:), allocatable, intent(out) :: energy
1456
1457 integer :: u_size, u_dim_id, u_id, file_id
1458 character(len=3) :: u_dims = "U_i"
1459
1460 integer :: n_steps, i
1461
1462 character(len=20) :: lattice
1463
1464 ! Filename from to which to read
1465 character(len=*), intent(in) :: filename
1466
1467 ! Ids for variables
1468 integer :: rho_id
1469
1470 ! Open the file for read-in
1471 call check(nf90_open(filename, nf90_nowrite, file_id))
1472
1473 ! Double-check the lattice type
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))
1476
1477 ! Check the dimensions
1478 do i = 1, rho_ndims
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)))
1481 end do
1482
1483 ! Check that these match the current simulation
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()"
1494 end if
1495
1496 ! Allocate array for read-in
1497 allocate(asro(rho_sizes(1), rho_sizes(2), rho_sizes(3), rho_sizes(4)))
1498
1499 ! Get the ID of the variable
1500 call check(nf90_inq_varid(file_id, "rho data", rho_id))
1501
1502 ! Read it in
1503 call check(nf90_get_var(file_id, rho_id, asro))
1504
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))
1507
1508 allocate(energy(u_size))
1509
1510 call check(nf90_inq_varid(file_id, "U data", u_id))
1511 call check(nf90_get_var(file_id, u_id, energy))
1512
1513 ! Close the file
1514 call check(nf90_close(file_id))
1515
1516 end subroutine ncdf_radial_density_reader
1517
1518end module netcdf_io
integer, parameter real64
Longer "double" (64 bit, approx -1.8e308 to 1.8e308 and covering values down to about 2e-308 magnitud...
Definition kinds.f90:59
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.
Definition netcdf_io.f90:57
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...