This repository has been archived by the owner on Jul 20, 2021. It is now read-only.
-
Notifications
You must be signed in to change notification settings - Fork 5
/
Copy pathstella_save.fpp
1119 lines (900 loc) · 36.2 KB
/
stella_save.fpp
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
# include "define.inc"
module stella_save
use mp, only: mp_comm, mp_info
# ifdef NETCDF
! use netcdf, only: NF90_FLOAT, NF90_DOUBLE
# ifdef NETCDF_PARALLEL
! If using netcdf version 4.1.2 or older delete NF90_MPIIO
use netcdf, only: NF90_HDF5,NF90_MPIIO
use netcdf, only: nf90_var_par_access, NF90_COLLECTIVE
use netcdf, only: nf90_put_att, NF90_GLOBAL, nf90_get_att
# endif
use netcdf, only: NF90_NOWRITE, NF90_CLOBBER, NF90_NOERR
use netcdf, only: nf90_create, nf90_open, nf90_sync, nf90_close
use netcdf, only: nf90_def_dim, nf90_def_var, nf90_enddef
use netcdf, only: nf90_put_var, nf90_get_var, nf90_strerror
use netcdf, only: nf90_inq_dimid, nf90_inquire_dimension
use netcdf, only: nf90_inq_varid, nf90_inquire_variable
use netcdf, only: nf90_int
use netcdf_utils, only: get_netcdf_code_precision
use netcdf_utils, only: check_netcdf_file_precision
use netcdf_utils, only: netcdf_error
use netcdf_utils, only: netcdf_real, kind_nf
# endif
implicit none
public :: stella_restore, stella_save_for_restart
public :: read_many, save_many
public :: init_save, init_dt, init_tstart, finish_save
!# ifdef NETCDF
! public :: netcdf_real, kind_nf, get_netcdf_code_precision, netcdf_error
!# endif
interface stella_restore
module procedure stella_restore_many
end interface
logical :: read_many = .true., save_many = .true. ! Read and write single or multiple restart files
private
character (300), save :: restart_file
# ifdef NETCDF
real, allocatable, dimension (:,:,:) :: tmpr, tmpi
real, allocatable, dimension (:,:,:,:) :: ktmpr, ktmpi
real, allocatable, dimension (:,:,:,:) :: ptmpr, ptmpi
integer (kind_nf) :: ncid, zedid, vpaid, gloid, gvmuloid, kyid, kxid, muid, tubeid
integer (kind_nf) :: krookr_id, krooki_id, projr_id, proji_id
! integer (kind_nf) :: bparr_id, bpari_id
integer (kind_nf) :: t0id, gr_id, gi_id, delt0id, istep0id
integer (kind_nf) :: intkrook_id, intproj_id;
integer (kind_nf) :: shift_id
logical :: initialized = .false.
# endif
contains
!!----------------------------------------------------------------------!!
!!----------------------------------------------------------------------!!
!!--Save----------------------------------------------------------------!!
!!----------------------------------------------------------------------!!
!!----------------------------------------------------------------------!!
subroutine stella_save_for_restart &
(g, istep0, t0, delt0, istatus, exit_in, fileopt)
# ifdef NETCDF
use fields_arrays, only: shift_state
use dist_fn_arrays, only: g_krook, g_proj
use kt_grids, only: naky, nakx
# else
use mp, only: proc0
# endif
use mp, only: iproc, barrier
# ifdef NETCDF_PARALLEL
use zgrid, only: nztot
# endif
use zgrid, only: nzgrid, ntubes
! Must include kxkyz_layout_type here to avoid obscure bomb while compiling
! stella_diagnostics.f90 (which uses this module) with the Compaq F90 compiler:
use stella_layouts, only: kxkyz_lo, vmu_lo
use common_types, only: kxkyz_layout_type
use file_utils, only: error_unit
use vpamu_grids, only: nvpa, nmu
use dissipation, only: include_krook_operator, int_krook
use dissipation, only: remove_zero_projection, int_proj
use physics_flags, only: prp_shear_enabled
implicit none
complex, dimension (:,:,kxkyz_lo%llim_proc:), intent (in) :: g
real, intent (in) :: t0, delt0
integer, intent (in) :: istep0
integer, intent (out) :: istatus
logical, intent (in), optional :: exit_in
character (20), INTENT (in), optional :: fileopt
# ifdef NETCDF
character (306) :: file_proc
character (10) :: suffix
integer :: i, n_elements, nvmulo_elements, ierr
integer :: total_elements, total_vmulo_elements
logical :: has_vmulo
# ifdef NETCDF_PARALLEL
integer, dimension(3) :: start_pos, counts
# endif
logical :: exit
!*********-----------------------_**********************
istatus = 0
if (present(exit_in)) then
exit = exit_in
else
exit = .false.
end if
! if (proc0) then
! write (*,*) "Starting save_for_restart in ", restart_file
! write (*,*) "List restart files"
! call system("echo 'start' >> filelist.txt; ls nc/* >> filelist.txt; ")
! end if
n_elements = kxkyz_lo%ulim_proc-kxkyz_lo%llim_proc+1
total_elements = kxkyz_lo%ulim_world+1
nvmulo_elements = vmu_lo%ulim_proc-vmu_lo%llim_proc+1
total_vmulo_elements = vmu_lo%ulim_world+1
if (n_elements <= 0) return
has_vmulo = nvmulo_elements.gt.0.or..not.save_many
if (.not.initialized) then
initialized = .true.
file_proc = trim(restart_file)
!CMR, 5/4/2011: Add optional piece of filename
IF (PRESENT(fileopt)) THEN
file_proc=trim(file_proc)//trim(fileopt)
END IF
!CMRend
!</HL> The NETCDF_PARALLEL directives include code for parallel
! netcdf using HDF5 to write the output to a single restart file
! The read_many flag allows the old style multiple file output
# ifdef NETCDF_PARALLEL
if(save_many) then
# endif
WRITE (suffix,'(a1,i0)') '.', iproc
# ifdef NETCDF_PARALLEL
else
WRITE (suffix,*) ''
endif
# endif
file_proc = trim(trim(file_proc)//adjustl(suffix))
# ifdef NETCDF_PARALLEL
if(save_many) then
# endif
istatus = nf90_create (file_proc, NF90_CLOBBER, ncid)
# ifdef NETCDF_PARALLEL
else
call barrier
if(iproc .eq. 0) then
open(unit=tmpunit, file=file_proc)
close(unit=tmpunit, status='delete')
end if
call barrier
! If using netcdf version 4.1.2 or older replace NF90_MPIIO with NF90_CLOBBER
istatus = nf90_create (file_proc, IOR(NF90_HDF5,NF90_MPIIO), ncid, comm=mp_comm, info=mp_info)
end if
# endif
if (istatus /= NF90_NOERR) then
ierr = error_unit()
write(ierr,*) "nf90_create error: ", nf90_strerror(istatus)
goto 1
end if
# ifdef NETCDF_PARALLEL
if(.not.save_many) then
istatus = nf90_put_att(ncid, NF90_GLOBAL, 'xyzs_layout', xyzs_layout)
if (istatus /= NF90_NOERR) then
ierr = error_unit()
write(ierr,*) "nf90_put_attr error: ", nf90_strerror(istatus)
goto 1
end if
istatus = nf90_put_att(ncid, NF90_GLOBAL, 'vms_layout', vms_layout)
if (istatus /= NF90_NOERR) then
ierr = error_unit()
write(ierr,*) "nf90_put_attr error: ", nf90_strerror(istatus)
goto 1
end if
endif
# endif
if (n_elements > 0) then
istatus = nf90_def_dim (ncid, "tube", ntubes, tubeid)
if (istatus /= NF90_NOERR) then
ierr = error_unit()
write(ierr,*) "nf90_def_dim zed error: ", nf90_strerror(istatus)
goto 1
end if
istatus = nf90_def_dim (ncid, "zed", 2*nzgrid+1, zedid)
if (istatus /= NF90_NOERR) then
ierr = error_unit()
write(ierr,*) "nf90_def_dim zed error: ", nf90_strerror(istatus)
goto 1
end if
istatus = nf90_def_dim (ncid, "vpa", nvpa, vpaid)
if (istatus /= NF90_NOERR) then
ierr = error_unit()
write(ierr,*) "nf90_def_dim vpa error: ", nf90_strerror(istatus)
goto 1
end if
istatus = nf90_def_dim (ncid, "mu", nmu, muid)
if (istatus /= NF90_NOERR) then
ierr = error_unit()
write(ierr,*) "nf90_def_dim mu error: ", nf90_strerror(istatus)
goto 1
end if
# ifdef NETCDF_PARALLEL
if(save_many) then
# endif
istatus = nf90_def_dim (ncid, "glo", n_elements, gloid)
# ifdef NETCDF_PARALLEL
else
istatus = nf90_def_dim (ncid, "glo", total_elements, gloid)
endif
# endif
if (istatus /= NF90_NOERR) then
ierr = error_unit()
write(ierr,*) "nf90_def_dim glo error: ", nf90_strerror(istatus)
goto 1
end if
# ifdef NETCDF_PARALLEL
if(save_many) then
# endif
if(nvmulo_elements.gt.0) then
istatus = nf90_def_dim (ncid, "gvmulo", nvmulo_elements, gvmuloid)
if (istatus /= NF90_NOERR) then
ierr = error_unit()
write(ierr,*) "nf90_def_dim gvmulo error: ", nf90_strerror(istatus)
goto 1
endif
endif
# ifdef NETCDF_PARALLEL
else
istatus = nf90_def_dim (ncid, "gvmulo", total_vmulo_elements, gvmuloid)
if (istatus /= NF90_NOERR) then
ierr = error_unit()
write(ierr,*) "nf90_def_dim gvmulo error: ", nf90_strerror(istatus)
goto 1
endif
endif
# endif
istatus = nf90_def_dim (ncid, "aky", naky, kyid)
if (istatus /= NF90_NOERR) then
ierr = error_unit()
write(ierr,*) "nf90_def_dim aky error: ", nf90_strerror(istatus)
goto 1
end if
istatus = nf90_def_dim (ncid, "akx", nakx, kxid)
if (istatus /= NF90_NOERR) then
ierr = error_unit()
write(ierr,*) "nf90_def_dim akx error: ", nf90_strerror(istatus)
goto 1
end if
end if
if (netcdf_real == 0) netcdf_real = get_netcdf_code_precision()
istatus = nf90_def_var (ncid, "t0", netcdf_real, t0id)
if (istatus /= NF90_NOERR) then
ierr = error_unit()
write(ierr,*) "nf90_def_var t0 error: ", nf90_strerror(istatus)
goto 1
end if
istatus = nf90_def_var (ncid, "istep0", nf90_int, istep0id)
if (istatus /= NF90_NOERR) then
ierr = error_unit()
write(ierr,*) "nf90_def_var istep0 error: ", nf90_strerror(istatus)
goto 1
end if
istatus = nf90_def_var (ncid, "delt0", netcdf_real, delt0id)
if (istatus /= NF90_NOERR) then
ierr = error_unit()
write(ierr,*) "nf90_def_var delt0 error: ", nf90_strerror(istatus)
goto 1
end if
if (n_elements > 0) then
istatus = nf90_def_var (ncid, "gr", netcdf_real, &
(/ vpaid, muid, gloid /), gr_id)
if (istatus /= NF90_NOERR) then
ierr = error_unit()
write(ierr,*) "nf90_def_var g error: ", nf90_strerror(istatus)
goto 1
end if
istatus = nf90_def_var (ncid, "gi", netcdf_real, &
(/ vpaid, muid, gloid /), gi_id)
if (istatus /= NF90_NOERR) then
ierr = error_unit()
write(ierr,*) "nf90_def_var g error: ", nf90_strerror(istatus)
goto 1
end if
if (include_krook_operator.and.has_vmulo) then
istatus = nf90_def_var (ncid, "intkrook", netcdf_real, intkrook_id)
if (istatus /= NF90_NOERR) then
ierr = error_unit()
write(ierr,*) "nf90_def_var intkrook error: ", nf90_strerror(istatus)
goto 1
end if
istatus = nf90_def_var (ncid, "krookr", netcdf_real, &
(/ kxid, zedid, tubeid, gvmuloid /), krookr_id)
if (istatus /= NF90_NOERR) then
ierr = error_unit()
write(ierr,*) "nf90_def_var apar error: ", nf90_strerror(istatus)
goto 1
end if
istatus = nf90_def_var (ncid, "krooki", netcdf_real, &
(/ kxid, zedid, tubeid, gvmuloid /), krooki_id)
if (istatus /= NF90_NOERR) then
ierr = error_unit()
write(ierr,*) "nf90_def_var krooki error: ", nf90_strerror(istatus)
goto 1
end if
end if
if (remove_zero_projection.and.has_vmulo) then
istatus = nf90_def_var (ncid, "intproj", netcdf_real, intproj_id)
if (istatus /= NF90_NOERR) then
ierr = error_unit()
write(ierr,*) "nf90_def_var intproj error: ", nf90_strerror(istatus)
goto 1
end if
istatus = nf90_def_var (ncid, "projr", netcdf_real, &
(/ kxid, zedid, tubeid, gvmuloid /), projr_id)
if (istatus /= NF90_NOERR) then
ierr = error_unit()
write(ierr,*) "nf90_def_var projr error: ", nf90_strerror(istatus)
goto 1
end if
istatus = nf90_def_var (ncid, "proji", netcdf_real, &
(/ kxid, zedid, tubeid, gvmuloid /), proji_id)
if (istatus /= NF90_NOERR) then
ierr = error_unit()
write(ierr,*) "nf90_def_var proji error: ", nf90_strerror(istatus)
goto 1
end if
end if
if (prp_shear_enabled) then
istatus = nf90_def_var (ncid, "shiftstate", netcdf_real,&
(/ kyid /), shift_id)
if (istatus /= NF90_NOERR) then
ierr = error_unit()
write(ierr,*) "nf90_def_var shiftstate error: ", nf90_strerror(istatus)
goto 1
end if
endif
! if (fbpar > epsilon(0.)) then
! istatus = nf90_def_var (ncid, "bpar_r", netcdf_real, &
! (/ zedid, kxid, kyid /), bparr_id)
! if (istatus /= NF90_NOERR) then
! ierr = error_unit()
! write(ierr,*) "nf90_def_var bparr error: ", nf90_strerror(istatus)
! goto 1
! end if
! istatus = nf90_def_var (ncid, "bpar_i", netcdf_real, &
! (/ zedid, kxid, kyid /), bpari_id)
! if (istatus /= NF90_NOERR) then
! ierr = error_unit()
! write(ierr,*) "nf90_def_var bpari error: ", nf90_strerror(istatus)
! goto 1
! end if
! end if
end if
! remove allocated conditional because we want to be able to restart
! using exb shear from a case which does not have exb shear (i.e.
! we need kx_shift variable defined in netcdf file even if no exb
! shear present in simulation) -- MAB + CMR
! if (allocated(kx_shift)) then ! MR begin
! istatus = nf90_def_var (ncid, "kx_shift", netcdf_real, &
! (/ kyid /), kx_shift_id)
! if (istatus /= NF90_NOERR) then
! ierr = error_unit()
! write(ierr,*) "nf90_def_var kx_shift error: ", nf90_strerror(istatus)
! goto 1
! endif
! endif ! MR end
! if (proc0) then
! write (*,*) "Finished definitions"
! write (*,*) "List restart files"
! call system("echo 'defs' >> filelist.txt; ls nc/* >> filelist.txt; ")
! end if
istatus = nf90_enddef (ncid)
if (istatus /= NF90_NOERR) then
ierr = error_unit()
write (ierr,*) "nf90_enddef error: ", nf90_strerror(istatus)
goto 1
end if
end if
!!!-----------------------!!!
!!!-----------------------!!!
!!!-----------------------!!!
# ifdef NETCDF_PARALLEL
if(save_many .or. iproc == 0) then
# endif
istatus = nf90_put_var (ncid, delt0id, delt0)
if (istatus /= NF90_NOERR) then
ierr = error_unit()
write (ierr,*) "nf90_put_var delt0 error: ", nf90_strerror(istatus)
goto 1
end if
istatus = nf90_put_var (ncid, t0id, t0)
if (istatus /= NF90_NOERR) then
ierr = error_unit()
write (ierr,*) "nf90_put_var t0 error: ", nf90_strerror(istatus)
goto 1
end if
istatus = nf90_put_var (ncid, istep0id, istep0)
if (istatus /= NF90_NOERR) then
ierr = error_unit()
write (ierr,*) "nf90_put_var istep0 error: ", nf90_strerror(istatus)
goto 1
end if
# ifdef NETCDF_PARALLEL
endif
# endif
1 continue
if (istatus /= NF90_NOERR) then
i = nf90_close (ncid)
return
end if
if (n_elements > 0) then
if (.not. allocated(tmpr)) &
allocate (tmpr(nvpa,nmu,kxkyz_lo%llim_proc:kxkyz_lo%ulim_alloc))
tmpr = real(g)
# ifdef NETCDF_PARALLEL
if(save_many) then
# endif
istatus = nf90_put_var (ncid, gr_id, tmpr)
#ifdef NETCDF_PARALLEL
else
istatus = nf90_var_par_access(ncid, gr_id, NF90_COLLECTIVE)
istatus = nf90_var_par_access(ncid, gi_id, NF90_COLLECTIVE)
start_pos = (/1,1,kxkyz_lo%llim_proc+1/)
counts = (/nvpa, nmu, n_elements/)
istatus = nf90_put_var (ncid, gr_id, tmpr, start=start_pos, count=counts)
endif
# endif
if (istatus /= NF90_NOERR) call netcdf_error (istatus, ncid, gr_id)
tmpr = aimag(g)
# ifdef NETCDF_PARALLEL
if(save_many) then
# endif
istatus = nf90_put_var (ncid, gi_id, tmpr)
#ifdef NETCDF_PARALLEL
else
istatus = nf90_put_var (ncid, gi_id, tmpr, start=start_pos, count=counts)
endif
# endif
if (istatus /= NF90_NOERR) call netcdf_error (istatus, ncid, gi_id)
if (include_krook_operator) then
if (.not. allocated(ktmpr)) &
allocate (ktmpr(nakx,-nzgrid:nzgrid,ntubes,vmu_lo%llim_proc:vmu_lo%ulim_alloc))
if (.not. allocated(ktmpi)) &
allocate (ktmpi(nakx,-nzgrid:nzgrid,ntubes,vmu_lo%llim_proc:vmu_lo%ulim_alloc))
# ifdef NETCDF_PARALLEL
if(save_many .or. iproc == 0) then
# endif
istatus = nf90_put_var (ncid, intkrook_id, int_krook)
if (istatus /= NF90_NOERR) then
ierr = error_unit()
write (ierr,*) "nf90_put_var int_krook error: ", nf90_strerror(istatus)
goto 1
end if
# ifdef NETCDF_PARALLEL
endif
# endif
ktmpr = real(g_krook)
ktmpi = aimag(g_krook)
# ifdef NETCDF_PARALLEL
if(save_many) then
# endif
istatus = nf90_put_var (ncid, krookr_id, ktmpr)
#ifdef NETCDF_PARALLEL
else
istatus = nf90_var_par_access(ncid, krookr_id, NF90_COLLECTIVE)
istatus = nf90_var_par_access(ncid, krooki_id, NF90_COLLECTIVE)
start_pos = (/1,1,1,vmu_lo%llim_proc+1/)
counts = (/nakx, nztot, ntubes, nvmulo_elements/)
istatus = nf90_put_var (ncid, krookr_id, ktmpr, start=start_pos, count=counts)
endif
# endif
if (istatus /= NF90_NOERR) call netcdf_error (istatus, ncid, krookr_id)
# ifdef NETCDF_PARALLEL
if(save_many) then
# endif
istatus = nf90_put_var (ncid, krooki_id, ktmpi)
#ifdef NETCDF_PARALLEL
else
istatus = nf90_put_var (ncid, krooki_id, ktmpi, start=start_pos, count=counts)
endif
# endif
if (istatus /= NF90_NOERR) call netcdf_error (istatus, ncid, krooki_id)
end if
if (remove_zero_projection) then
if (.not. allocated(ptmpr)) &
allocate (ptmpr(nakx,-nzgrid:nzgrid,ntubes,vmu_lo%llim_proc:vmu_lo%ulim_alloc))
if (.not. allocated(ptmpi)) &
allocate (ptmpi(nakx,-nzgrid:nzgrid,ntubes,vmu_lo%llim_proc:vmu_lo%ulim_alloc))
# ifdef NETCDF_PARALLEL
if(save_many .or. iproc == 0) then
# endif
istatus = nf90_put_var (ncid, intproj_id, int_proj)
if (istatus /= NF90_NOERR) then
ierr = error_unit()
write (ierr,*) "nf90_put_var int_proj error: ", nf90_strerror(istatus)
goto 1
end if
# ifdef NETCDF_PARALLEL
endif
# endif
ptmpr = real(g_proj)
ptmpi = aimag(g_proj)
# ifdef NETCDF_PARALLEL
if(save_many) then
# endif
istatus = nf90_put_var (ncid, projr_id, ptmpr)
#ifdef NETCDF_PARALLEL
else
istatus = nf90_var_par_access(ncid, projr_id, NF90_COLLECTIVE)
istatus = nf90_var_par_access(ncid, proji_id, NF90_COLLECTIVE)
start_pos = (/1,1,1,vmu_lo%llim_proc+1/)
counts = (/nakx, nztot, ntubes, nvmulo_elements/)
istatus = nf90_put_var (ncid, projr_id, ptmpr, start=start_pos, count=counts)
endif
# endif
if (istatus /= NF90_NOERR) call netcdf_error (istatus, ncid, projr_id)
# ifdef NETCDF_PARALLEL
if(save_many) then
# endif
istatus = nf90_put_var (ncid, proji_id, ptmpi)
#ifdef NETCDF_PARALLEL
else
istatus = nf90_put_var (ncid, proji_id, ptmpi, start=start_pos, count=counts)
endif
# endif
if (istatus /= NF90_NOERR) call netcdf_error (istatus, ncid, proji_id)
end if
if (prp_shear_enabled) then
istatus = nf90_put_var (ncid, shift_id, shift_state)
if (istatus /= NF90_NOERR) call netcdf_error (istatus, ncid, shift_id)
end if
end if
if (exit) then
i = nf90_close (ncid)
if (i /= NF90_NOERR) &
call netcdf_error (istatus, message='nf90_close error')
else
i = nf90_sync (ncid)
if (i /= NF90_NOERR) &
call netcdf_error (istatus, message='nf90_sync error')
end if
# else
if (proc0) write (error_unit(),*) &
'WARNING: stella_save_for_restart is called without netcdf library'
# endif
if (allocated(tmpr)) deallocate (tmpr)
if (allocated(tmpi)) deallocate (tmpi)
if (allocated(ptmpr)) deallocate (ptmpr)
if (allocated(ptmpi)) deallocate (ptmpi)
if (allocated(ktmpr)) deallocate (ktmpr)
if (allocated(ktmpi)) deallocate (ktmpi)
end subroutine stella_save_for_restart
!!----------------------------------------------------------------------!!
!!----------------------------------------------------------------------!!
!!---Restart------------------------------------------------------------!!
!!----------------------------------------------------------------------!!
!!----------------------------------------------------------------------!!
subroutine stella_restore_many (g, scale, istatus)
# ifdef NETCDF
use mp, only: iproc
use fields_arrays, only: shift_state
use dist_fn_arrays, only: g_krook, g_proj
use kt_grids, only: naky, nakx
# endif
# ifdef NETCDF_PARALLEL
use zgrid, only: nztot
# endif
use zgrid, only: nzgrid, ntubes
use vpamu_grids, only: nvpa, nmu
use stella_layouts, only: kxkyz_lo, vmu_lo
use file_utils, only: error_unit
use dissipation, only: include_krook_operator, int_krook
use dissipation, only: remove_zero_projection, int_proj
use physics_flags, only: prp_shear_enabled
implicit none
complex, dimension (:,:,kxkyz_lo%llim_proc:), intent (out) :: g
real, intent (in) :: scale
integer, intent (out) :: istatus
# ifdef NETCDF
# ifdef NETCDF_PARALLEL
integer, dimension(3) :: counts, start_pos
# endif
character (306) :: file_proc
character (10) :: suffix
integer :: i, n_elements, nvmulo_elements, ierr
logical :: has_vmulo
n_elements = kxkyz_lo%ulim_proc-kxkyz_lo%llim_proc+1
nvmulo_elements = vmu_lo%ulim_proc-vmu_lo%llim_proc+1
if (n_elements <= 0) return
has_vmulo = nvmulo_elements.gt.0.or..not.read_many
if (.not.initialized) then
! initialized = .true.
file_proc = trim(restart_file)
# ifdef NETCDF_PARALLEL
if(read_many) then
# endif
write (suffix,'(a1,i0)') '.', iproc
file_proc = trim(trim(file_proc)//adjustl(suffix))
istatus = nf90_open (file_proc, NF90_NOWRITE, ncid)
# ifdef NETCDF_PARALLEL
else
! If using netcdf version 4.1.2 deleted NF90_MPIIO and the associated IOR
istatus = nf90_open (file_proc, IOR(NF90_NOWRITE, NF90_MPIIO), ncid, comm=mp_comm, info=mp_info)
endif
# endif
if (istatus /= NF90_NOERR) then
call netcdf_error (istatus, file=file_proc, abort=.true.)
endif
! check precision
if (netcdf_real == 0) netcdf_real = get_netcdf_code_precision()
call check_netcdf_file_precision (ncid)
istatus = nf90_inq_dimid (ncid, "tube", tubeid)
if (istatus /= NF90_NOERR) call netcdf_error (istatus, dim='tube')
istatus = nf90_inq_dimid (ncid, "zed", zedid)
if (istatus /= NF90_NOERR) call netcdf_error (istatus, dim='zed')
istatus = nf90_inq_dimid (ncid, "aky", kyid)
if (istatus /= NF90_NOERR) call netcdf_error (istatus, dim='aky')
istatus = nf90_inq_dimid (ncid, "akx", kxid)
if (istatus /= NF90_NOERR) call netcdf_error (istatus, dim='akx')
istatus = nf90_inq_dimid (ncid, "glo", gloid)
if (istatus /= NF90_NOERR) call netcdf_error (istatus, dim='glo')
if(has_vmulo) then
istatus = nf90_inq_dimid (ncid, "gvmulo", gvmuloid)
if (istatus /= NF90_NOERR) call netcdf_error (istatus, dim='gvmulo')
endif
istatus = nf90_inquire_dimension (ncid, tubeid, len=i)
if (istatus /= NF90_NOERR) call netcdf_error (istatus, ncid, dimid=tubeid)
if (i /= ntubes) write(*,*) 'Restart error: ntubes=? ',i,' : ',ntubes,' : ',iproc
istatus = nf90_inquire_dimension (ncid, zedid, len=i)
if (istatus /= NF90_NOERR) call netcdf_error (istatus, ncid, dimid=zedid)
if (i /= 2*nzgrid + 1) write(*,*) 'Restart error: nzgrid=? ',i,' : ',nzgrid,' : ',iproc
istatus = nf90_inquire_dimension (ncid, kyid, len=i)
if (istatus /= NF90_NOERR) call netcdf_error (istatus, ncid, dimid=kyid)
if (i /= naky) write(*,*) 'Restart error: naky=? ',i,' : ',naky,' : ',iproc
istatus = nf90_inquire_dimension (ncid, kxid, len=i)
if (istatus /= NF90_NOERR) call netcdf_error (istatus, ncid, dimid=kxid)
if (i /= nakx) write(*,*) 'Restart error: nakx=? ',i,' : ',nakx,' : ',iproc
istatus = nf90_inquire_dimension (ncid, gloid, len=i)
if (istatus /= NF90_NOERR) call netcdf_error (istatus, ncid, dimid=gloid)
#ifdef NETCDF_PARALLEL
if(read_many) then
#endif
if (i /= n_elements) write(*,*) 'Restart error: glo=? ',i,' : ',iproc
#ifdef NETCDF_PARALLEL
else
if (i /= kxkyz_lo%ulim_world+1) write(*,*) 'Restart error: glo=? ',i,' : ',iproc
endif
#endif
if(has_vmulo) then
istatus = nf90_inquire_dimension (ncid, gvmuloid, len=i)
if (istatus /= NF90_NOERR) call netcdf_error (istatus, ncid, dimid=gvmuloid)
#ifdef NETCDF_PARALLEL
if(read_many) then
#endif
if (i /= nvmulo_elements) write(*,*) 'Restart error: gvmulo=? ',i,' : ',iproc
#ifdef NETCDF_PARALLEL
else
if (i /= vmu_lo%ulim_world+1) write(*,*) 'Restart error: gvmulo=? ',i,' : ',iproc
endif
#endif
endif
if(include_krook_operator.and.has_vmulo) then
istatus = nf90_inq_varid (ncid, "intkrook", intkrook_id)
if (istatus /= NF90_NOERR) call netcdf_error (istatus, var='intkrook')
istatus = nf90_inq_varid (ncid, "krookr", krookr_id)
if (istatus /= NF90_NOERR) call netcdf_error (istatus, var='krookr')
istatus = nf90_inq_varid (ncid, "krooki", krooki_id)
if (istatus /= NF90_NOERR) call netcdf_error (istatus, var='krooki')
endif
if(remove_zero_projection.and.has_vmulo) then
istatus = nf90_inq_varid (ncid, "intproj", intproj_id)
if (istatus /= NF90_NOERR) call netcdf_error (istatus, var='intproj')
istatus = nf90_inq_varid (ncid, "projr", projr_id)
if (istatus /= NF90_NOERR) call netcdf_error (istatus, var='projr')
istatus = nf90_inq_varid (ncid, "proji", proji_id)
if (istatus /= NF90_NOERR) call netcdf_error (istatus, var='proji')
endif
if(prp_shear_enabled) then
istatus = nf90_inq_varid (ncid, "shiftstate", shift_id)
if (istatus /= NF90_NOERR) call netcdf_error (istatus, var='shiftstate')
endif
! if (fbpar > epsilon(0.)) then
! istatus = nf90_inq_varid (ncid, "bpar_r", bparr_id)
! if (istatus /= NF90_NOERR) call netcdf_error (istatus, var='bpar_r')
! istatus = nf90_inq_varid (ncid, "bpar_i", bpari_id)
! if (istatus /= NF90_NOERR) call netcdf_error (istatus, var='bpar_i')
! end if
! if (allocated(kx_shift)) then ! MR begin
! istatus = nf90_inq_varid (ncid, "kx_shift", kx_shift_id)
! if (istatus /= NF90_NOERR) call netcdf_error (istatus, var='kx_shift')
! endif ! MR end
istatus = nf90_inq_varid (ncid, "gr", gr_id)
if (istatus /= NF90_NOERR) call netcdf_error (istatus, var='gr')
istatus = nf90_inq_varid (ncid, "gi", gi_id)
if (istatus /= NF90_NOERR) call netcdf_error (istatus, var='gi')
end if
if (.not. allocated(tmpr)) &
allocate (tmpr(nvpa,nmu,kxkyz_lo%llim_proc:kxkyz_lo%ulim_alloc))
if (.not. allocated(tmpi)) &
allocate (tmpi(nvpa,nmu,kxkyz_lo%llim_proc:kxkyz_lo%ulim_alloc))
tmpr = 0.; tmpi = 0.
# ifdef NETCDF_PARALLEL
if(read_many) then
# endif
istatus = nf90_get_var (ncid, gr_id, tmpr)
#ifdef NETCDF_PARALLEL
else
start_pos = (/1,1,kxkyz_lo%llim_proc+1/)
counts = (/nvpa, nmu, n_elements/)
istatus = nf90_get_var (ncid, gr_id, tmpr, start=start_pos, count=counts)
end if
# endif
if (istatus /= NF90_NOERR) call netcdf_error (istatus, ncid, gr_id)
# ifdef NETCDF_PARALLEL
if(read_many) then
# endif
istatus = nf90_get_var (ncid, gi_id, tmpi)
#ifdef NETCDF_PARALLEL
else
istatus = nf90_get_var (ncid, gi_id, tmpi, start=start_pos, count=counts)
end if
# endif
if (istatus /= NF90_NOERR) call netcdf_error (istatus, ncid, gi_id)
g = cmplx(tmpr, tmpi)
if(include_krook_operator.and.has_vmulo) then
if (.not. allocated(ktmpr)) &
allocate (ktmpr(nakx,-nzgrid:nzgrid,ntubes,vmu_lo%llim_proc:vmu_lo%ulim_alloc))
if (.not. allocated(ktmpi)) &
allocate (ktmpi(nakx,-nzgrid:nzgrid,ntubes,vmu_lo%llim_proc:vmu_lo%ulim_alloc))
istatus = nf90_get_var (ncid, intkrook_id, int_krook)
if (istatus /= NF90_NOERR) call netcdf_error (istatus, ncid, intkrook_id)
ktmpr = 0.; ktmpi = 0.
# ifdef NETCDF_PARALLEL
if(read_many) then
# endif
istatus = nf90_get_var (ncid, krookr_id, ktmpr)
#ifdef NETCDF_PARALLEL
else
start_pos = (/1,1,1,vmu_lo%llim_proc+1/)
counts = (/nakx, nztot, ntubes, nvmulo_elements/)
istatus = nf90_get_var (ncid, krookr_id, ktmpr, start=start_pos, count=counts)
end if
# endif
if (istatus /= NF90_NOERR) call netcdf_error (istatus, ncid, krookr_id)
# ifdef NETCDF_PARALLEL
if(read_many) then
# endif
istatus = nf90_get_var (ncid, krooki_id, ktmpi)
#ifdef NETCDF_PARALLEL
else
istatus = nf90_get_var (ncid, krooki_id, ktmpi, start=start_pos, count=counts)
end if
# endif
if (istatus /= NF90_NOERR) call netcdf_error (istatus, ncid, krooki_id)
g_krook = cmplx(ktmpr, ktmpi)
endif
if(remove_zero_projection.and.has_vmulo) then
if (.not. allocated(ptmpr)) &
allocate (ptmpr(nakx,-nzgrid:nzgrid,ntubes,vmu_lo%llim_proc:vmu_lo%ulim_alloc))
if (.not. allocated(ptmpi)) &
allocate (ptmpi(nakx,-nzgrid:nzgrid,ntubes,vmu_lo%llim_proc:vmu_lo%ulim_alloc))
istatus = nf90_get_var (ncid, intproj_id, int_proj)
if (istatus /= NF90_NOERR) call netcdf_error (istatus, ncid, intproj_id)
ptmpr = 0.; ptmpi = 0.
# ifdef NETCDF_PARALLEL
if(read_many) then
# endif
istatus = nf90_get_var (ncid, projr_id, ptmpr)
#ifdef NETCDF_PARALLEL
else
start_pos = (/1,1,1,vmu_lo%llim_proc+1/)
counts = (/nakx, nztot, ntubes, nvmulo_elements/)
istatus = nf90_get_var (ncid, projr_id, ptmpr, start=start_pos, count=counts)
end if
# endif
if (istatus /= NF90_NOERR) call netcdf_error (istatus, ncid, projr_id)
# ifdef NETCDF_PARALLEL
if(read_many) then
# endif
istatus = nf90_get_var (ncid, proji_id, ptmpi)
#ifdef NETCDF_PARALLEL
else
istatus = nf90_get_var (ncid, proji_id, ptmpi, start=start_pos, count=counts)
end if
# endif
if (istatus /= NF90_NOERR) call netcdf_error (istatus, ncid, proji_id)
g_proj = cmplx(ptmpr, ptmpi)
endif
if(prp_shear_enabled) then
istatus = nf90_get_var (ncid, shift_id, shift_state)
if (istatus /= NF90_NOERR) call netcdf_error (istatus, ncid, shift_id)
endif
if (scale > 0.) then
g = g*scale
if(include_krook_operator) g_krook = g_krook*scale
if(remove_zero_projection) g_proj = g_proj*scale
endif
! RN 2008/05/23: this was commented out. why? HJL 2013/05/15 Because it stops future writing to the file
! istatus = nf90_close (ncid)
if (istatus /= NF90_NOERR) then
ierr = error_unit()
write(ierr,*) "nf90_close error: ", nf90_strerror(istatus),' ',iproc
end if
# else
write (error_unit(),*) &
'ERROR: stella_restore_many is called without netcdf'
# endif
if (allocated(tmpr)) deallocate (tmpr)
if (allocated(tmpi)) deallocate (tmpi)
if (allocated(ptmpr)) deallocate (ptmpr)
if (allocated(ptmpi)) deallocate (ptmpi)
if (allocated(ktmpr)) deallocate (ktmpr)
if (allocated(ktmpi)) deallocate (ktmpi)
end subroutine stella_restore_many
subroutine init_save (file)
character(300), intent (in) :: file
restart_file = file
end subroutine init_save