-
Notifications
You must be signed in to change notification settings - Fork 24
/
Copy pathpot3d.F90
7202 lines (7200 loc) · 203 KB
/
pot3d.F90
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
!#######################################################################
! _____ ____ _______ ____ _____
! | __ \ / __ \__ __|___ \| __ \
! | |__) | | | | | | __) | | | |
! | ___/| | | | | | |__ <| | | |
! | | | |__| | | | ___) | |__| |
! |_| \____/ |_| |____/|_____/
!
! ****** POT3D: 3D potential magnetic field outside a sphere.
!
! ****** This program can find the classical potential field, the
! ****** fully open field, the source-surface field, and the
! ****** source-surface plus current-sheet field.
!
! Authors: Zoran Mikic
! Ronald M. Caplan
! Jon A. Linker
! Roberto Lionello
! Miko Stulajter
!
! Predictive Science Inc.
! www.predsci.com
! San Diego, California, USA 92121
!
!#######################################################################
! Copyright 2021 Predictive Science Inc.
!
! Licensed under the Apache License, Version 2.0 (the "License");
! you may not use this file except in compliance with the License.
! You may obtain a copy of the License at
!
! http://www.apache.org/licenses/LICENSE-2.0
!
! Unless required by applicable law or agreed to in writing, software
! distributed under the License is distributed on an "AS IS" BASIS,
! WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or
! implied.
! See the License for the specific language governing permissions and
! limitations under the License.
!#######################################################################
!
!#######################################################################
!
module ident
!
!-----------------------------------------------------------------------
!
implicit none
!
!-----------------------------------------------------------------------
! ****** Code name.
!-----------------------------------------------------------------------
!
character(*), parameter :: idcode='POT3D'
character(*), parameter :: vers ='4.3.1'
character(*), parameter :: update='12/05/2024'
!
end module
!#######################################################################
module number_types
!
!-----------------------------------------------------------------------
! ****** Basic number types.
! ****** This module is used to set the default precision for REALs.
!-----------------------------------------------------------------------
!
use iso_fortran_env
!
!-----------------------------------------------------------------------
!
implicit none
!
integer, parameter :: KIND_REAL_4=REAL32
integer, parameter :: KIND_REAL_8=REAL64
integer, parameter :: KIND_REAL_16=max(REAL128,REAL64)
!
integer, parameter :: r_typ=KIND_REAL_8
!
end module
!#######################################################################
module number_types_pc
!
!-----------------------------------------------------------------------
!
use number_types
use iso_fortran_env
!
!-----------------------------------------------------------------------
!
implicit none
!
!-----------------------------------------------------------------------
!
integer, parameter :: r_typ_pc=REAL32
!
end module
!#######################################################################
module constants
!
!-----------------------------------------------------------------------
!
use number_types
!
!-----------------------------------------------------------------------
!
implicit none
!
real(r_typ), parameter :: pi=3.1415926535897932_r_typ
!
end module
!#######################################################################
module global_dims
!
!-----------------------------------------------------------------------
! ****** Global number of mesh points.
!-----------------------------------------------------------------------
!
implicit none
!
! ****** Global mesh size.
!
integer :: nr_g,nrm1_g
integer :: nt_g,ntm1_g
integer :: np_g,npm1_g
!
! ****** Flag to indicate an axisymmetric run.
!
logical :: axisymmetric=.false.
!
end module
!#######################################################################
module local_dims
!
!-----------------------------------------------------------------------
! ****** Local number of mesh points.
!-----------------------------------------------------------------------
!
implicit none
!
! ****** Flags to indicate whether this processor has points
! ****** on the physical boundaries.
!
logical :: rb0,rb1,tb0,tb1
!
! ****** Local mesh size.
!
integer :: nr,nrm1
integer :: nt,ntm1
integer :: np,npm1
!
! ****** Dimensions of arrays on the "main" mesh.
!
integer :: nrm
integer :: ntm
integer :: npm
!
! ****** Indices of start and end points in the global mesh
! ****** belonging to this processor.
!
integer :: i0_g,i1_g
integer :: j0_g,j1_g
integer :: k0_g,k1_g
!
end module
!#######################################################################
module global_mesh
!
!-----------------------------------------------------------------------
! ****** Global mesh.
!-----------------------------------------------------------------------
!
use number_types
use constants
!
implicit none
!
real(r_typ), dimension(:), allocatable :: r_g,rh_g,dr_g,drh_g
real(r_typ), dimension(:), allocatable :: t_g,th_g,dt_g,dth_g
real(r_typ), dimension(:), allocatable :: p_g,ph_g,dp_g,dph_g
!
real(r_typ), dimension(:), allocatable :: st_g,ct_g,sth_g,cth_g
real(r_typ), dimension(:), allocatable :: sp_g,cp_g,sph_g,cph_g
!
! ****** Physical mesh size.
!
real(r_typ), parameter, private :: one=1._r_typ
real(r_typ), parameter, private :: two=2._r_typ
!
real(r_typ) :: r0=1._r_typ
real(r_typ) :: r1=30._r_typ
real(r_typ), parameter :: t0=0.
real(r_typ), parameter :: t1=pi
real(r_typ), parameter :: p0=0.
real(r_typ), parameter :: p1=two*pi
!
real(r_typ), parameter :: pl=p1-p0
real(r_typ), parameter :: pl_i=one/pl
!
end module
!#######################################################################
module local_mesh
!
!-----------------------------------------------------------------------
! ****** Local mesh.
!-----------------------------------------------------------------------
!
use number_types
!
implicit none
!
real(r_typ), dimension(:), allocatable :: r,r2,rh,dr,drh
real(r_typ) :: dr1,drn
!
real(r_typ), dimension(:), allocatable :: t,th,dt,dth
real(r_typ), dimension(:), allocatable :: p,ph,dp,dph
!
real(r_typ), dimension(:), allocatable :: st,ct,sth,cth
real(r_typ), dimension(:), allocatable :: sp,cp,sph,cph
!
! ****** Inverse quantities (for efficiency).
!
real(r_typ), dimension(:), allocatable :: r_i,rh_i
real(r_typ), dimension(:), allocatable :: dr_i,drh_i
real(r_typ), dimension(:), allocatable :: dt_i,dth_i
real(r_typ), dimension(:), allocatable :: st_i,sth_i
real(r_typ), dimension(:), allocatable :: dp_i,dph_i
!
end module
!#######################################################################
module mpidefs
!
!-----------------------------------------------------------------------
! ****** MPI variables, processor topology, and processor information.
!-----------------------------------------------------------------------
!
use mpi
!
implicit none
!
! ****** Total number of processors.
!
integer :: nproc
!
! ****** Total number of processors per node.
!
integer :: nprocsh
!
! ****** Processor rank of this process in communicator
! ****** MPI_COMM_WORLD.
!
integer :: iprocw
!
! ****** Processor rank of this process in communicator
! ****** comm_shared.
!
integer :: iprocsh
!
! ****** Flag to designate that this is the processor with
! ****** rank 0 in communicator MPI_COMM_WORLD.
!
logical :: iamp0
!
! ****** Communicator over all processors in the Cartesian topology.
!
integer :: comm_all
!
! ****** Processor rank of this process in communicator
! ****** COMM_ALL.
!
integer :: iproc
!
! ****** Processor rank in communicator COMM_ALL for the
! ****** processor that has rank 0 in MPI_COMM_WORLD.
!
integer :: iproc0
!
! ****** Communicators over all processors in the phi dimension.
!
integer :: comm_phi
!
! ****** Communicator over all shared processors on the node.
!
integer :: comm_shared
!
! ****** Communicators over all processors in the r dimension.
!
integer :: comm_r
!
! ****** Processor coordinate indices of this process
! ****** in the Cartesian topology.
!
integer :: iproc_r,iproc_t,iproc_p
!
! ****** Processor coordinate indices of the neighboring
! ****** processors in the Cartesian topology.
!
integer :: iproc_rm,iproc_rp
integer :: iproc_tm,iproc_tp
integer :: iproc_pm,iproc_pp
!
! ****** Number of processors along r, theta, and phi.
!
integer :: nproc_r,nproc_t,nproc_p
!
! ****** Number type for REALs to be used in MPI calls.
!
integer :: ntype_real
!
end module
!#######################################################################
module decomposition_params
!
!-----------------------------------------------------------------------
! ****** Input parameters that define the domain decomposition
! ****** among processors.
!-----------------------------------------------------------------------
!
implicit none
!
! ****** Number of processors per dimension.
!
integer, dimension(3) :: nprocs=(/-1,-1,-1/)
!
end module
!#######################################################################
module decomposition
!
!-----------------------------------------------------------------------
! ****** Information defining the domain decomposition.
!-----------------------------------------------------------------------
!
implicit none
!
! ****** Define the structure type for mapping local arrays
! ****** into global arrays.
!
type :: map_struct
integer :: n
integer :: i0
integer :: i1
integer :: offset
end type
!
! ****** Mapping structures for the different mesh types.
!
type(map_struct), dimension(:), pointer :: map_rh
type(map_struct), dimension(:), pointer :: map_rm
type(map_struct), dimension(:), pointer :: map_th
type(map_struct), dimension(:), pointer :: map_tm
type(map_struct), dimension(:), pointer :: map_ph
type(map_struct), dimension(:), pointer :: map_pm
!
end module
!#######################################################################
module meshdef
!
!-----------------------------------------------------------------------
! ****** Variables that define the mesh-point distributions.
!-----------------------------------------------------------------------
!
use number_types
!
implicit none
!
integer, parameter :: nmseg=30
!
real(r_typ), dimension(nmseg) :: drratio=0.
real(r_typ), dimension(nmseg) :: dtratio=0.
real(r_typ), dimension(nmseg) :: dpratio=0.
real(r_typ), dimension(nmseg) :: rfrac=0.
real(r_typ), dimension(nmseg) :: tfrac=0.
real(r_typ), dimension(nmseg) :: pfrac=0.
!
integer :: nfrmesh=0
integer :: nftmesh=0
integer :: nfpmesh=0
!
real(r_typ) :: phishift=0.
!
end module
!#######################################################################
module fields
!
!-----------------------------------------------------------------------
! ****** Local field arrays.
!-----------------------------------------------------------------------
!
use number_types
!
implicit none
!
! ****** Potential solution array and cg temp array.
!
real(r_typ), dimension(:,:,:), allocatable :: phi
real(r_typ), dimension(:,:,:), allocatable :: x_ax
!
! ****** Boundary radial magnetic field arrays.
!
real(r_typ), dimension(:,:), allocatable :: br0
real(r_typ), dimension(:,:), allocatable :: br1
!
! ****** Arrays used in polar boundary conditions.
!
real(r_typ), dimension(:), allocatable :: sum0,sum1
!
! ****** Arrays used for final magnetic field.
!
real(r_typ), dimension(:,:,:), allocatable :: br,bt,bp
!
end module
!#######################################################################
module cgcom
!
use number_types
!
implicit none
!
!-----------------------------------------------------------------------
! ****** Number of equations to solve in the CG solver.
!-----------------------------------------------------------------------
!
integer :: ncgeq
!
!-----------------------------------------------------------------------
! ****** CG field solver parameters.
!-----------------------------------------------------------------------
!
integer :: ifprec=1
integer :: ncgmax=1000000
integer :: ncghist=100
real(r_typ) :: epscg=1.e-9
!
!-----------------------------------------------------------------------
! ****** CG field solver variables.
!-----------------------------------------------------------------------
!
integer :: ncg
real(r_typ) :: epsn
!
! ****** Seam buffers.
!
real(r_typ), allocatable, dimension(:,:) :: sbuf_rt1,sbuf_rt2
real(r_typ), allocatable, dimension(:,:) :: sbuf_tp1,sbuf_tp2
real(r_typ), allocatable, dimension(:,:) :: sbuf_rp1,sbuf_rp2
!
real(r_typ), allocatable, dimension(:,:) :: rbuf_rt1,rbuf_rt2
real(r_typ), allocatable, dimension(:,:) :: rbuf_tp1,rbuf_tp2
real(r_typ), allocatable, dimension(:,:) :: rbuf_rp1,rbuf_rp2
!
end module
!#######################################################################
module vars
!
!-----------------------------------------------------------------------
! ****** Miscellaneous input variables.
!-----------------------------------------------------------------------
!
use number_types
!
implicit none
!
character(256) :: outfile='pot3d.out'
character(256) :: phifile='default'
character(256) :: br0file='default'
character(256) :: brfile='default'
character(256) :: btfile='default'
character(256) :: bpfile='default'
character(256) :: br_photo_file='default'
character(256) :: br_photo_original_file='default'
!
! ****** Type of field solution.
! ****** Select between 'potential', 'open', and 'source-surface'.
!
character(16) :: option='potential'
!
! ****** Interval at which to dump diagonstics during the
! ****** iteration for the source-surface plus current-sheet
! ****** solution.
!
integer :: ndump=0
!
! ****** Flag to skip the balancing of the flux (for PFSS and
! ****** OPEN field options only).
logical :: do_not_balance_flux=.false.
!
! ****** Set format for output files.
!
character(3) :: fmt='h5'
!
logical :: hdf32=.true.
!
! ***** Validation run (tilted dipole).
!
logical :: validation_run=.false.
!
real(r_typ) :: dipole_angle=0.7853981633974483_r_typ
!
end module
!#######################################################################
module solve_params
!
!-----------------------------------------------------------------------
! ****** Parameters used in the solver.
!-----------------------------------------------------------------------
!
use number_types
!
implicit none
!
! ****** Boundary condition switch at r=R1.
!
real(r_typ) :: pm_r1
!
end module
!#######################################################################
module timer
!
!-----------------------------------------------------------------------
! ****** Timer stack.
!-----------------------------------------------------------------------
!
use number_types
!
implicit none
!
integer, parameter :: nstack=10
integer :: istack=0
real(r_typ), dimension(nstack) :: tstart=0.
!
end module
!#######################################################################
module timing
!
!-----------------------------------------------------------------------
! ****** Timing variables.
!-----------------------------------------------------------------------
!
use number_types
!
implicit none
!
real(r_typ) :: t_startup=0.
real(r_typ) :: t_solve=0.
real(r_typ) :: t_pc_load=0.
real(r_typ) :: t_pc=0.
real(r_typ) :: t_ax=0.
real(r_typ) :: t_io=0.
real(r_typ) :: c_seam=0.
real(r_typ) :: c_cgdot=0.
real(r_typ) :: c_sumphi=0.
real(r_typ) :: t_wall=0.
!
end module
!#######################################################################
module debug
!
!-----------------------------------------------------------------------
! ****** Debugging level.
!-----------------------------------------------------------------------
!
implicit none
!
integer :: idebug=0
!
end module
!#######################################################################
module assemble_array_interface
interface
subroutine assemble_array (map_r,map_t,map_p,a,a_g)
use number_types
use decomposition
use mpidefs
implicit none
type(map_struct), dimension(0:nproc-1) :: map_r,map_t,map_p
real(r_typ), dimension(:,:,:) :: a,a_g
end subroutine
end interface
end module
!#######################################################################
module cusparse_interface
!
use, intrinsic :: iso_c_binding
!
#ifdef CUSPARSE
interface
subroutine load_lusol_cusparse(CSR_A,CSR_AI,CSR_AJ,N,M) &
BIND(C, name="load_lusol_cusparse")
use, intrinsic :: iso_c_binding
implicit none
integer(C_INT), value :: N,M
type(C_PTR), value :: CSR_A,CSR_AI,CSR_AJ
end subroutine load_lusol_cusparse
!
subroutine lusol_cusparse(x) &
BIND(C, name="lusol_cusparse")
use, intrinsic :: iso_c_binding
implicit none
type(C_PTR), value :: x
end subroutine lusol_cusparse
!
subroutine unload_lusol_cusparse() &
BIND(C, name="unload_lusol_cusparse")
end subroutine unload_lusol_cusparse
end interface
#endif
!
integer(c_int) :: cN,cM
!
end module
!#######################################################################
module matrix_storage_pot3d_solve
!
!-----------------------------------------------------------------------
! ****** Storage for the matrix/preconditioners of the solve.
!-----------------------------------------------------------------------
!
use number_types
use number_types_pc
!
implicit none
!
real(r_typ), dimension(:,:,:,:), allocatable :: a
real(r_typ_pc), dimension(:), allocatable :: a_i
!
integer, dimension(7) :: a_offsets
integer :: N,M
real(r_typ_pc), dimension(:), allocatable, target :: a_csr
real(r_typ_pc), dimension(:), allocatable :: lu_csr
real(r_typ_pc), dimension(:), allocatable :: a_csr_d
integer, dimension(:), allocatable :: lu_csr_ja
integer, dimension(:), allocatable, target :: a_csr_ia
integer, dimension(:), allocatable, target :: a_csr_ja
integer, dimension(:), allocatable :: a_N1
integer, dimension(:), allocatable :: a_N2
integer, dimension(:), allocatable :: a_csr_dptr
!
end module
!#######################################################################
program POT3D
!
!-----------------------------------------------------------------------
!
use ident
use mpidefs
use vars
use solve_params
use timing
!
!-----------------------------------------------------------------------
!
implicit none
!
!-----------------------------------------------------------------------
!
integer :: ierr
!
!-----------------------------------------------------------------------
!
! ****** Initialize MPI.
!
call init_mpi
!
! ****** Start the wall-clock timer.
!
call timer_on
!
! ****** Write the code name and version.
!
if (iamp0) then
write (*,*)
write (*,*) idcode,' Version ',vers,', updated on ',update
end if
!
call timer_on
!
! ****** Read the input file.
!
call read_input_file
!
! ****** Create the output file.
!
if (iamp0) then
call ffopen (9,outfile,'rw',ierr)
if (ierr.ne.0) then
write (*,*)
write (*,*) '### ERROR in POT3D:'
write (*,*) '### Could not create the output file.'
write (*,*) 'File name: ',trim(outfile)
end if
end if
call check_error_on_p0 (ierr)
!
! ****** Check the input parameters.
!
call check_input
!
! ****** Check the processor topology.
!
call check_proc_topology
!
! ****** Decompose the domain.
!
call decompose_domain
!
! ****** Allocate global arrays.
!
call allocate_global_arrays
!
! ****** Set the global meshes.
!
call set_global_mesh
!
! ****** Decompose the mesh.
!
call decompose_mesh_r
call decompose_mesh_tp
!
! ****** Allocate local arrays.
!
call allocate_local_arrays_r
call allocate_local_arrays_tp
!
! ****** Set the local meshes.
!
call set_local_mesh_r
call set_local_mesh_tp
!
! ****** Print decomposition diagnostics.
!
call decomp_diags
!
! ****** Set up seam for solver.
!
call seam_setup
!
! ****** Initialize the flux and if validating, write analytic solution.
!
if (validation_run) then
if (iamp0) then
write (*,*)
write (*,*) '### COMMENT from POT3D:'
write (*,*) '### Validation run requested.'
write (*,*) '### Ignoring br input file,'
write (*,*) '### setting HDF32=.FALSE.,'
write (*,*) '### and overriding output filenames.'
write (9,*)
write (9,*) '### COMMENT from POT3D:'
write (9,*) '### Validation run requested.'
write (9,*) '### Ignoring br input file'
write (9,*) '### setting HDF32=.FALSE.,'
write (9,*) '### and overriding output filenames.'
end if
hdf32=.false.
brfile='br_solved.'//trim(fmt)
btfile='bt_solved.'//trim(fmt)
bpfile='bp_solved.'//trim(fmt)
phifile='phi_solved.'//trim(fmt)
call set_validation_flux
call write_validation_solution
else
call set_flux
end if
!
call timer_off (t_startup)
!
! ****** Find the solution.
!
if (iamp0) then
write (*,*)
write (*,*) '### COMMENT from POT3D:'
write (*,*) '### Starting PCG solve.'
call FLUSH(OUTPUT_UNIT)
write (9,*)
write (9,*) '### COMMENT from POT3D:'
write (9,*) '### Starting PCG solve.'
call FLUSH(9)
end if
!
call timer_on
!
call potfld
!
call timer_off (t_solve)
!
! ****** Compute B.
!
call getb
!
! ****** Write solution to file.
!
call write_solution
!
! ****** Magnetic energy diagnostics.
!
call magnetic_energy
!
call MPI_Barrier(MPI_COMM_WORLD,ierr)
call timer_off (t_wall)
!
call write_timing
!
call endrun (.false.)
!
end program pot3d
!#######################################################################
subroutine read_input_file
!
!-----------------------------------------------------------------------
!
! ****** Read the input file.
!
!-----------------------------------------------------------------------
!
use global_dims
use global_mesh
use mpidefs
use meshdef
use cgcom
use debug
use vars
use decomposition_params
!
!-----------------------------------------------------------------------
!
implicit none
!
!-----------------------------------------------------------------------
!
! ****** Values for the global mesh size.
! ****** Since these names conflict with those in LOCAL_DIMS*, it is
! ****** important not to use these modules here.
!
integer :: nr=0
integer :: nt=0
integer :: np=0
!
!-----------------------------------------------------------------------
!
namelist /topology/ &
nprocs, &! MPI topology triplet. Default -1,-1,-1
! automatically sets "best" topology.
nr, &! Grid resolution in the `r` direction.
nt, &! Grid resolution in the `theta` direction.
np ! Grid resolution in the `phi` direction.
!
!-----------------------------------------------------------------------
!
namelist /inputvars/ &
r0, &! Lower radial boundary.
r1, &! Upper radial boundary.
drratio, &! Ratio of grid spacing at the end
! of each segment to that at the
! beginning for the radial grid
! [ length(drratio) = length(rfrac)-1 ].
dtratio, &! Ratio of grid spacing (theta)
dpratio, &! Ratio of grid spacing (phi)
rfrac, &! Normalized positions of grid segment
! boundaries (frac of domain size)
! for radial grid.
tfrac, &! Normalized positions of grid (theta)
pfrac, &! Normalized positions of grid (phi)
nfrmesh, &! Number of times to filter/smooth
! the radial grid spacing.
nftmesh, &! Number of times to filter (theta)
nfpmesh, &! Number of times to filter (phi)
phishift, &! Apply an optional phi shift
! (radians) to the input Br at r0.
ifprec, &! Preconditioner method:
! 1: Diagonal (use for GPU runs)
! 2: ILU0 (use for CPU runs or
! GPU runs when built with cusparse)
ncgmax, &! Maximum alowed solver iterations.
ncghist, &! Iteration information.
! 0: Only write # total iterations.
! >0: Write every NCGHIST iteration.
epscg, &! Solver convergence tolerance,
!|residual|/|right-hand-side|.
idebug, &! Output debugging info during run.
br0file, &! Filename of input 2D (tp) Br
phifile, &! Filename to write 3D PHI potential.
brfile, &! Filename to write 3D Br field.
btfile, &! Filename to write 3D Bt field.
bpfile, &! Filename to write 3D Bp field.
br_photo_file, &! Filename to write 2D Br@r=r0
!(after interp/flux-balance).
br_photo_original_file, &! Filename to write 2D Br@r=r0
! (option='open' only, writes Br
! before sign change).
option, &! 'ss' PFSS
! 'potential' PF with closed-wall
! (requires flux balance).
! 'open' Open field. Used for
! current sheet and fully
! open field runs.
! B field will be unsigned.
do_not_balance_flux, &! Do not balance flux of input Br@r0.
hdf32, &! Output precision:
! .true. Single (32-bit) output.
! .false. Double (64-bit) output.
validation_run, &! Set this to run a validation test
! with an analytic tilted dipole
! solution. Overrides other inputs.
dipole_angle ! Tilt angle for validation run dipole.
!
!-----------------------------------------------------------------------
!
integer :: ierr
character(80) :: infile='pot3d.dat'
!
!-----------------------------------------------------------------------
!
! ****** Read the input file.
!
call ffopen (8,infile,'r',ierr)
!
if (ierr.ne.0) then
if (iamp0) then
write (*,*)
write (*,*) '### ERROR in READ_INPUT_FILE:'
write (*,*) '### Could not open the input file.'
write (*,*) 'File name: ',trim(infile)
end if
call endrun (.true.)
end if
!
read (8,topology)
!
read (8,inputvars)
!
close (8)
!
if (trim(fmt).ne.'h5') then
if (iamp0) then
write (*,*)
write (*,*) '### ERROR in READ_INPUT_FILE:'
write (*,*) '### "fmt" must be "h5".'
write (*,*) 'fmt: ',trim(fmt)
end if
call endrun (.true.)
end if
!
! ****** Check if output names were overwritten.
! ****** If not, set default names with format fmt.
!
if (trim(phifile).eq.'default') then
phifile='phi.'//trim(fmt)
end if
if (trim(br0file).eq.'default') then
br0file='br0.'//trim(fmt)
end if
if (trim(brfile).eq.'default') then
brfile='br.'//trim(fmt)
end if
if (trim(btfile).eq.'default') then
btfile='bt.'//trim(fmt)
end if
if (trim(bpfile).eq.'default') then
bpfile='bp.'//trim(fmt)
end if
if (trim(br_photo_file).eq.'default') then
br_photo_file='br_photo.'//trim(fmt)
end if
if (trim(br_photo_original_file).eq.'default') then
br_photo_original_file='br_photo_original.'//trim(fmt)
end if
!
nr_g=nr
nt_g=nt
np_g=np
!
! ****** Check if the specified mesh dimensions are valid.
!
call check_mesh_dimensions (nr_g,nt_g,np_g)
!
nrm1_g=nr_g-1
ntm1_g=nt_g-1
npm1_g=np_g-1
!
end subroutine
!#######################################################################
subroutine check_error_on_p0 (ierr0)
!
!-----------------------------------------------------------------------
!
! ****** Check if the error flag IERR0 on processor 0 in
! ****** MPI_COMM_WORLD (i.e., processor IPROC0 in COMM_ALL)
! ****** indicates that the code should exit.
!
! ****** If IERR0 is non-zero, all the processors are directed
! ****** to call ENDRUN to terminate the code.
!
!-----------------------------------------------------------------------
!
use mpidefs
!
!-----------------------------------------------------------------------
!
implicit none