-
Notifications
You must be signed in to change notification settings - Fork 2
/
Copy pathPMGsolver.F90
10223 lines (9840 loc) · 490 KB
/
PMGsolver.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
!********************************************************************
! PMGsolver
! A parallel multigrid solver for the 2D linear problems
! arising in the plasma fluid model GDB (see B. Zhu, et. al.,
! GDB: a global 3D two-fluid code for plasma turbulence and
! transport in the tokamak edge region, J. of Comp. Phys 2017).
!
! These problems are:
! 1. A Poisson solver to obtain the electrostatic potential phi
! from the vorticity w
! div(lambda*grad(u))=rho
! where lambda=plasma density, u=phi and rho=w.
! 2. An inhomogeneous modified Helmholtz equation to obtain the
! magnetic flux function psi from the variable rpsi which
! contains an electron inertia contribution
! u-ksq*L(u)=rho
! where u=psi, rho=rpsi and ksq=(electron skin depth)^2
! In the routines below we use the suffix 'ps' to indicate
! routines to obtain psi
! 3. Higher order variants of the (modified) Helmholtz equation
! arising from implied hyperdiffusion of order 2m
! u+((-1)^(m+1))*Diff*(L^2m)(u)=rho
! where Diff is the diffusion coefficient and L is the 2D
! Laplacian operator, although a non-Laplacian-based diffusion
! operator is also available respectively. For more details
! see M. Francisquez, et. al., Multigrid treatment of implicit
! continuum diffusion, J. of Comp. Phys 2017)
! 4. A simple Poisson solver
! div(grad(u))=rho
! In the routines below we use the suffix 'ph' to indicate
! routines to obtain phi, 'ps' for psi and 'hd' for hyperdiffusion.
!
! Though these are 2D solvers, we allow for a 3rd dimension so that
! PMGsolver effetively supports solving many xy planes of
! independent 2D problems.
!
! Manaure Francisquez
! 2017
!
!********************************************************************
! IMPORTANT NOTES
! i) Not all solvers support the same boundary conditions
! ii) Modified Helmholtz relaxation has damped relaxation commented out
! iii) Some restriction and prolongation operators have been removed
! for simplicity. I believe IR=0, ISR=1, IPRO=2 are the only options
! available in this file.
! ASSUMPTIONS
! a) Periodic along y and not periodic along x
! b) The length of the simulation domain is defined as 2*pi times
! a factor specified in the input file (fLx along x)
! c) The problem, the grids and the MPI decomposition are such that
! only two subdomains along x and/or y can be unified in a single
! coarsening.
!
! Poisson
! 1. Each process has at least nxmin and nymin points
! along x and y, respectively.
! 2. nxmin,nymin >= 4
! 3. smallest grid has >=4 points along x
! 4. Boundary condition on the right (along x) is Dirichlet(? check)
! 5. The boundary conditions of lambda are assumed even
! 6. Use boundary linear extrapolation when restricting lambda
!
! Hyperdiffusion
! 1. Hyperdiffusion may be applied to more than one quantity
! (e.g. density, temperature, etc.), and a different multigrid
! parameter set is given for each. nqhd = number of quantities
! hyperdiffusion is applied to. This variable is inferred from
! the size of iDiff, and the code looks for nqhd parameter sets
! in the input file.
! 2. HDop: type of diffusion, some options are (L=Laplacian)
! =2 u + Diff*(L^2)u=rho
! =3 u - Diff*(L^3)u=rho
! =4 u + Diff*(L^2)u=rho split into 2 eqns
! =6 u - Diff*(L^3)u=rho split into 3 eqns
! =8 u + Diffx*(d^4 u/dx^4) + Diffy*(d^4 u/dy^4)=rho
! =12 u - Diffx*(d^6 u/dx^6) - Diffy*(d^6 u/dy^6)=rho
!********************************************************************
module PMGsolver
implicit none
#include "gdbFLAGS.h"
!.Symbolic names for kind types of 4-byte integers and double precision reals:
integer, parameter :: I4B = SELECTED_INT_KIND(9)
integer, parameter :: DP = KIND(1.0D0)
!.Maximum number of grids to visit (=0 if no max)
integer(I4B) :: LGMAXph,LGMAXps
integer(I4B), allocatable :: LGMAXhd(:)
!.Semi-coarsening strategy: =2 doubling =4 quadrupling
integer(I4B) :: NISph,NISps
integer(I4B), allocatable :: NIShd(:)
!.Type of cycle: GAMMA=1 V-cycle, GAMMA=2 W-cycle
integer(I4B) :: GAMph,GAMps
integer(I4B), allocatable :: GAMhd(:)
!.Number of (GAMMA-)cycles (at every grid if using FMG)
integer(I4B) :: BETph,BETps
integer(I4B), allocatable :: BEThd(:)
!.Number of relaxation sweeps before and after the course grid correction
integer(I4B) :: NU1ph,NU2ph,NU1ps,NU2ps
integer(I4B), allocatable :: NU1hd(:),NU2hd(:)
!.Damping coefficients in relaxation
real(DP) :: omeph,Oomeph,omeps,Oomeps
real(DP), allocatable :: omehd(:),Oomehd(:)
integer(I4B), parameter :: HDopd2=HDop/2,HDopd4=HDop/4
!.Number of elements needed from diagonal neighbors
!.(only needed when HDop=2 or 3)
integer(I4B), parameter :: ndiag=1+(HDop-2)*2
!.Number of quantities hyperdiffusion is applied to
integer(I4B) :: nqhd
!.Define a type for an array 3D real arrays
type allo3D
real(DP), allocatable :: a(:,:,:)
end type allo3D
!.Define a type for an array 2D real arrays
type allo2D
real(DP), allocatable :: a(:,:)
end type allo2D
!.Define a type for an array 1D integer arrays
type Iallo1D
integer(I4B), allocatable :: a(:)
end type Iallo1D
!.Define a type for an array 2D integer arrays
type Iallo2D
integer(I4B), allocatable :: a(:,:)
end type Iallo2D
!.Define a type for an array 3D integer arrays
type Iallo3D
integer(I4B), allocatable :: a(:,:,:)
end type Iallo3D
!.Define a type for an array 2D boolean arrays
type Lallo2D
logical, allocatable :: a(:,:)
end type Lallo2D
!.Define a type for an array of arrays of 1D integer arrays
type IIallo2D
type(Iallo2D), allocatable :: a(:)
end type IIallo2D
!.For Poisson equation, one needs lambda coarsened to each
!.grid level. That is a 3D array at each grid level
type(allo3D), allocatable :: Diff(:)
!.Local (to this rank) fine and coarse dimensions of grids visited
integer(I4B), allocatable, dimension(:) :: nxFph,nyFph,nxFps,nyFps
type(Iallo1D), allocatable, dimension(:) :: nxFhd,nyFhd
integer(I4B), allocatable, dimension(:) :: nxCph,nyCph,nxCps,nyCps
type(Iallo1D), allocatable, dimension(:) :: nxChd,nyChd
!.Indicate which dimension was coarsened, =0 for both, =1 for x, =2 for y
integer(I4B), allocatable, dimension(:) :: wdcph,wdcps
type(Iallo1D), allocatable, dimension(:) :: wdchd
!.Number of grids to be visited
integer(I4B) :: lgph,lgps
integer(I4B), allocatable :: lghd(:)
!.Spatial grid details
real(DP) :: bLx,bLy !.Length of the domain box along x and y
integer(I4B) :: nxG,nyG !.Number of points in global grid
integer(I4B) :: nxL,nyL,nzL !.Number of points in local grid
!.Constant hyper-diffusion coefficient
real(DP), allocatable :: hDiff(:)
!.k^2 factor in Helmholtz equation = de^2 (electron skin depth squared) in GDB
real(DP) :: ksq
!.Resolution of current grids, fine and coarse
integer(I4B) :: nxc,nyc,nxf,nyf
!.Local grid dimensions after unification of subdomains takes place
integer(I4B) :: nxu,nyu
!........................................................!
!................ MPI related variables ................!
integer(I4B) :: wrnk0 !.initial rank number in MPI_WORLD_COMM
integer(I4B) :: iprocs0,jprocs0 !.Number of (finest grid) processes in x and y
integer(I4B) :: iID0,jID0 !.This process' initial ID in x and y COMMs
integer(I4B) :: myrank0 !.This process' initial ID in xy COMM
integer(I4B) :: Wcomm0,XYcomm0 !.Initial World and XY COMMs
logical :: Lhalf0,Bhalf0 !.Is this rank in Left 1/2 (along x) and/or Bottom 1/2 (along y)
integer(I4B) :: nxmin,nymin !.min number of cells per process
!.In the following arrays, entries correspond to differend grid levels
!.Note: x (East-West) and y (North-South) rank
!.rank ID in X and Y COMMs
integer(I4B), allocatable, dimension(:) :: iIDph,jIDph,iIDps,jIDps
type(Iallo1D), allocatable, dimension(:) :: iIDhd,jIDhd
!.Number of processors along each direction
integer(I4B), allocatable, dimension(:) :: iprocsph,jprocsph,iprocsps, &
jprocsps
type(Iallo1D), allocatable, dimension(:) :: iprocshd,jprocshd
!.Number of processors along each direction (minus 1)
integer(I4B) :: iprocsm1,jprocsm1
!.Rank of neighbors, starting with the one above (along y/North) and
!.proceeding clockwise: 1-North, 2-NE, 3-E, 4-SE, ..., 8-NW
integer(I4B), allocatable, dimension(:,:) :: nborph,nborps
type(Iallo2D), allocatable, dimension(:) :: nborhd
!.Array of communicators of procs still active at grid level ng.
!.COMM(ng,1) is the y (North-South) communicator, COMM(ng,2) is the x
!.(East-West) communicator and COMM(ng,3) the 2D communicator for
!.diagonal communication
integer(I4B), allocatable, dimension(:,:) :: COMMph,COMMps
type(Iallo2D), allocatable, dimension(:) :: COMMhd
!.This array contains a boolean for indicating whether this rank is active
!.and whether unification of subdomains takes place or not
!.acu(ng,1) TRUE if this process is still active
!.acu(ng,2) TRUE if it is time to unify subdomains
logical, allocatable, dimension(:,:) :: acuph,acups
type(Lallo2D), allocatable, dimension(:) :: acuhd
!.Ranks data is sent to/received from when subdomains are unified
type(Iallo2D), allocatable :: urnkph(:),urnkps(:)
type(IIallo2D), allocatable :: urnkhd(:)
!.Number of ranks to whose subdomains will be unified along x and y
integer(I4B), allocatable :: nurnkph(:,:),nurnkps(:,:)
type(Iallo2D), allocatable :: nurnkhd(:)
#IF (MGTIMER==1)
!.These are user-defined timers for timing different parts of the code
!.used in profiling
real(DP) :: initPMG,init2,init1,mg1,mg2,rstdiff1,rstdiff2, &
prolt1,prolt2,count1,count2,rest1,rest2
real(DP) :: mgtph,resDtph,postXtph,preXtph,proltph,rNrtph,restph
real(DP) :: mgtps,resDtps,postXtps,preXtps,proltps,rNrtps,restps
real(DP) :: mgthd,resDthd,postXthd,preXthd,prolthd,rNrthd,resthd
#ENDIF
contains
subroutine initialize_PMG(iWcomm,iXYcomm,iXcomm,iYcomm,iDiff)
!.Initialize the PMG solver arrays and communicators needed
!.INPUTS
!. integer iWcomm: WORLD communicator
!. integer iXYcomm: 2D XY communicator
!. integer iXcomm: 1D X communicator
!. integer iYcomm: 1D Y communicator
!. real(dp) dimension(:,2) iDiff: Diffusion coefficients along x
!. and y of nqhd quantities
implicit none
include 'mpif.h'
integer(I4B), intent(in) :: iWcomm,iXYcomm,iXcomm,iYcomm
real(DP), intent(in) :: iDiff(:,:)
integer(I4B) :: ng,i,j,k,ierr,jnk1(1),jnk2(1),IDs1D(1)
logical :: jnkl(1)
#IF (MGTIMER==1)
initPMG=0.0_dp
call CPU_TIME(init1)
#ENDIF
!.Use a duplicate of the input world COMM for safety/simplicity
call MPI_COMM_DUP(iWcomm,Wcomm0,ierr) ! World communicator
call MPI_COMM_RANK(Wcomm0,wrnk0,ierr) ! rank in world comm
nqhd = size(iDiff,1) !.Number of quantities hyperdiffused
call read_mginputs !.Read input multigrid parameter set(s)
!.Use a duplicate of the input 2D COMM for safety/simplicity
!.The 2D COMM is used in unifying subdomains, prolongation
!.and diagonal communication
call MPI_COMM_DUP(iXYcomm,XYcomm0,ierr)
!.Initial grid rank ID's in each COMM
call MPI_Comm_rank(iXcomm,iID0,ierr)
call MPI_Comm_rank(iYcomm,jID0,ierr)
call MPI_CART_RANK(iXYcomm,(/jID0,iID0/),myrank0,ierr)
!.Is this rank in the Left Half of the domain (along x)?
Lhalf0=.TRUE.
if (iID0>=iprocs0/2) Lhalf0=.FALSE.
!.Is this rank in the Bottom Half of the domain (along y)?
Bhalf0=.TRUE.
if (jID0>=jprocs0/2) Bhalf0=.FALSE.
!.Local number of grid points in each direction
nxL=nxG/iprocs0
nyL=nyG/jprocs0
!.Determine the grids to be visited by the Poisson and
!.Helmholtz solvers, when to unify subdomains, and initialize necessary arrays
call pgrids(1,LGMAXph,NISph,lgph,nxFph,nyFph,nxCph,nyCph,wdcph, &
acuph,iprocsph,jprocsph,urnkph,nurnkph)
call pgrids(2,LGMAXps,NISps,lgps,nxFps,nyFps,nxCps,nyCps,wdcps, &
acups,iprocsps,jprocsps,urnkps,nurnkps)
!.Need to allocate space for coarsened lambda at each grid (Diff)
allocate(Diff(lgph))
do j=lgph,1,-1
allocate(Diff(j)%a(nxFph(j),nyFph(j),nzL))
enddo
!.Communicators at each grid level. The 2D communicator (COMM(ng,3))
!.is used to deal with corner grid points in stencils.
allocate(COMMph(lgph,3),COMMps(lgps,3))
!.Copy initial (finest) 2D communicator into array of communicators
call MPI_COMM_DUP(XYcomm0,COMMph(lgph,3),ierr)
call MPI_COMM_DUP(XYcomm0,COMMps(lgps,3),ierr)
!.Create 2D xy communicators for each coarse grid
call COMMv2D(COMMph(:,3),XYcomm0,acuph(:,1),iprocsph,jprocsph,lgph)
call COMMv2D(COMMps(:,3),XYcomm0,acups(:,1),iprocsps,jprocsps,lgps)
!.Save the rank ID of the diagonal neighbors (corners) NE,SE,SW,NW
allocate(nborph(lgph,8),nborps(lgps,8))
call diagNeigh(nborph,COMMph(:,3),iprocsph,jprocsph,acuph(:,1),lgph)
call diagNeigh(nborps,COMMps(:,3),iprocsps,jprocsps,acups(:,1),lgps)
!.Copy initial x and y communicators into array of communicators
call MPI_COMM_DUP(iYcomm,COMMph(lgph,1),ierr)
call MPI_COMM_DUP(iXcomm,COMMph(lgph,2),ierr)
call MPI_COMM_DUP(iYcomm,COMMps(lgps,1),ierr)
call MPI_COMM_DUP(iXcomm,COMMps(lgps,2),ierr)
!.Create x and y communicators for coarse grids
call COMMv1D(COMMph(:,1),iYcomm,acuph(:,1),jprocsph,lgph,(/.TRUE./),jID0)
call COMMv1D(COMMph(:,2),iXcomm,acuph(:,1),iprocsph,lgph,(/.FALSE./),iID0)
call COMMv1D(COMMps(:,1),iYcomm,acups(:,1),jprocsps,lgps,(/.TRUE./),jID0)
call COMMv1D(COMMps(:,2),iXcomm,acups(:,1),iprocsps,lgps,(/.FALSE./),iID0)
allocate(iIDph(lgph),iIDps(lgps))
allocate(jIDph(lgph),jIDps(lgps))
do ng=lgph,1,-1 !.Loop over grids visited
if (acuph(ng,1)) then !.Check if this rank is active (not idle)
!.....This rank's ID in x (East-West) direction
call MPI_Comm_rank(COMMph(ng,2),iIDph(ng),ierr)
!.....For communication along x, save East and West neighbor rank
call MPI_CART_SHIFT(COMMph(ng,2),0,1,nborph(ng,7),nborph(ng,3),ierr)
!.....MF: These are needed on NERSC. Not sure why.
if (iIDph(ng)==0) nborph(ng,7)=MPI_PROC_NULL
if (iIDph(ng)==iprocsph(lgph)) nborph(ng,3)=MPI_PROC_NULL
!.....This rank's ID in y (North-South) direction
call MPI_Comm_rank(COMMph(ng,1),jIDph(ng), ierr)
!.....For communication along y, save North and South neighbor rank
call MPI_CART_SHIFT(COMMph(ng,1),0,1,nborph(ng,5),nborph(ng,1),ierr)
endif
enddo
do ng=lgps,1,-1 !.Loop over grids visited
if (acups(ng,1)) then !.Check if this rank is active (not idle)
!.....This rank's ID in x (East-West) direction
call MPI_Comm_rank(COMMps(ng,2),iIDps(ng), ierr)
!.....For communication along x, save East and West neighbor rank
call MPI_CART_SHIFT(COMMps(ng,2),0,1,nborps(ng,7),nborps(ng,3),ierr)
!.....MF: These are needed on NERSC. Not sure why.
if (iIDps(ng)==0) nborps(ng,7)=MPI_PROC_NULL
if (iIDps(ng)==iprocsps(lgps)) nborps(ng,3)=MPI_PROC_NULL
!.....This rank's ID in y (North-South) direction
call MPI_Comm_rank(COMMps(ng,1),jIDps(ng), ierr)
!.....For communication along y, save North and South neighbor rank
call MPI_CART_SHIFT(COMMps(ng,1),0,1,nborps(ng,5),nborps(ng,1),ierr)
endif
enddo
!.Now do the same for hyperdiffusion quantities
allocate(hDiff(2)) ! x-y diffusion coefficients used when solver is called
!.Check nxmin and nymin are greater than 8
if (nxmin<8) then
nxmin=8
if (wrnk0==0) print*,'Changing to nxmin=8 for PMG hyperdiff'
endif
if (nymin<8) then
nymin=8
if (wrnk0==0) print*,'Changing to nymin=8 for PMG hyperdiff'
endif
!.Allocate arrays for hyperdiffusion parameters of various quantities
allocate(lghd(nqhd),nxFhd(nqhd),nyFhd(nqhd),nxChd(nqhd),nyChd(nqhd), &
wdchd(nqhd),acuhd(nqhd),iprocshd(nqhd),jprocshd(nqhd), &
urnkhd(nqhd),nurnkhd(nqhd))
allocate(COMMhd(nqhd),nborhd(nqhd),iIDhd(nqhd),jIDhd(nqhd))
do j=1,nqhd !.Loop over quantities hyperdiffused
!...Determine the grids to be visited, when to unify subdomains,
!...and initialize necessary arrays
call pgrids(3,LGMAXhd(j),NIShd(j),lghd(j),nxFhd(j)%a,nyFhd(j)%a, &
nxChd(j)%a,nyChd(j)%a,wdchd(j)%a,acuhd(j)%a,iprocshd(j)%a, &
jprocshd(j)%a,urnkhd(j)%a,nurnkhd(j)%a,iDiff(j,:))
!...Communicators at each grid level. The 2D communicator (COMM(ng,3))
!...is used to deal with corner grid points in stencils.
allocate(COMMhd(j)%a(lghd(j),3))
!...Copy 2D communicator into array of communicators
call MPI_COMM_DUP(XYcomm0,COMMhd(j)%a(lghd(j),3),ierr)
!...Create 2D xy communicators for each coarse grid
call COMMv2D(COMMhd(j)%a(:,3),XYcomm0,acuhd(j)%a(:,1), &
iprocshd(j)%a,jprocshd(j)%a,lghd(j))
!.Save the rank ID of the diagonal neighbors (corners) NE,SE,SW,NW
allocate(nborhd(j)%a(lghd(j),8))
call diagNeigh(nborhd(j)%a,COMMhd(j)%a(:,3),iprocshd(j)%a, &
jprocshd(j)%a,acuhd(j)%a(:,1),lghd(j))
!...Copy initial x and y communicators into array of communicators
call MPI_COMM_DUP(iYcomm,COMMhd(j)%a(lghd(j),1),ierr)
call MPI_COMM_DUP(iXcomm,COMMhd(j)%a(lghd(j),2),ierr)
!...Create x and y communicators for coarse grids
call COMMv1D(COMMhd(j)%a(:,1),iYcomm,acuhd(j)%a(:,1), &
jprocshd(j)%a,lghd(j),(/.TRUE./),jID0)
call COMMv1D(COMMhd(j)%a(:,2),iXcomm,acuhd(j)%a(:,1), &
iprocshd(j)%a,lghd(j),(/.FALSE./),iID0)
allocate(iIDhd(j)%a(lghd(j)),jIDhd(j)%a(lghd(j)))
do ng=lghd(j),1,-1 !.Loop over grids visited
if (acuhd(j)%a(ng,1)) then !.Check if this rank is active (not idle)
!.......This rank's ID in x (East-West) direction
call MPI_Comm_rank(COMMhd(j)%a(ng,2),iIDhd(j)%a(ng),ierr)
!.......For communication along x, save East and West neighbor rank
call MPI_CART_SHIFT(COMMhd(j)%a(ng,2),0,1, &
nborhd(j)%a(ng,7),nborhd(j)%a(ng,3),ierr)
!.......MF: These are needed on NERSC. Not sure why.
if (iIDhd(j)%a(ng)==0) nborhd(j)%a(ng,7)=MPI_PROC_NULL
if (iIDhd(j)%a(ng)==iprocshd(j)%a(lghd(j))) nborhd(j)%a(ng,3)=MPI_PROC_NULL
!.......This rank's ID in y (North-South) direction
call MPI_Comm_rank(COMMhd(j)%a(ng,1),jIDhd(j)%a(ng), ierr)
!.......For communication along y, save North and South neighbor rank
call MPI_CART_SHIFT(COMMhd(j)%a(ng,1),0,1, &
nborhd(j)%a(ng,5),nborhd(j)%a(ng,1),ierr)
endif
enddo
enddo
#IF (MGTIMER==1)
call CPU_TIME(init2)
initPMG=init2-init1
#ENDIF
end subroutine initialize_PMG
!*********************************************************
subroutine read_mginputs
!.Read the input resolution, MPI decomposition, factors defining
!.the domain lengh along x and y, and the multigrid parameter set.
!.The multigrid parameters control which grids are visited, the
!.type and number of iterations to be performed, and when to
!.unify subdomains.
implicit none
include 'mpif.h'
integer(I4B) :: nx0G,ny0G,nz0G,xprocs,yprocs,zprocs
real(DP) :: fLx,fLy,fLz
integer(I4B) :: ierr
allocate(LGMAXhd(nqhd),NIShd(nqhd),GAMhd(nqhd),BEThd(nqhd), &
NU1hd(nqhd),NU2hd(nqhd),omehd(nqhd),Oomehd(nqhd))
!.Add here the namelists for parameters (first declared above) from input file:
namelist/dom/nx0G,ny0G,nz0G,xprocs,yprocs,zprocs,fLx,fLy,fLz
namelist/mgin/LGMAXph,NISph,GAMph,BETph,NU1ph,NU2ph,omeph, &
LGMAXps,NISps,GAMps,BETps,NU1ps,NU2ps,omeps, &
LGMAXhd,NIShd,GAMhd,BEThd,NU1hd,NU2hd,omehd, &
nxmin,nymin
if (wrnk0 .eq. 0) then
open(unit=10,file="gdb.in",status="old")
read (10,dom)
read (10,mgin)
close(unit=10)
nxG=nx0G
nyG=ny0G
iprocs0=xprocs
jprocs0=yprocs
nzL=nz0G/zprocs
bLx=fLx*8.0_dp*atan(1.0_8) !.fLx*2*pi
bLy=fLy*8.0_dp*atan(1.0_8) !.fLy*2*pi
endif
!.Since only wrnk0 read the inputs, broadcast them to World COMM
call MPI_BCAST(nxG,1,MPI_INTEGER,0,Wcomm0,ierr)
call MPI_BCAST(nyG,1,MPI_INTEGER,0,Wcomm0,ierr)
call MPI_BCAST(nzL,1,MPI_INTEGER,0,Wcomm0,ierr)
call MPI_BCAST(iprocs0,1,MPI_INTEGER,0,Wcomm0,ierr)
call MPI_BCAST(jprocs0,1,MPI_INTEGER,0,Wcomm0,ierr)
call MPI_BCAST(bLx,1,MPI_DOUBLE_PRECISION,0,Wcomm0,ierr)
call MPI_BCAST(bLy,1,MPI_DOUBLE_PRECISION,0,Wcomm0,ierr)
call MPI_BCAST(LGMAXph,1,MPI_INTEGER,0,Wcomm0,ierr)
call MPI_BCAST(NISph,1,MPI_INTEGER,0,Wcomm0,ierr)
call MPI_BCAST(GAMph,1,MPI_INTEGER,0,Wcomm0,ierr)
call MPI_BCAST(BETph,1,MPI_INTEGER,0,Wcomm0,ierr)
call MPI_BCAST(NU1ph,1,MPI_INTEGER,0,Wcomm0,ierr)
call MPI_BCAST(NU2ph,1,MPI_INTEGER,0,Wcomm0,ierr)
call MPI_BCAST(omeph,1,MPI_DOUBLE_PRECISION,0,Wcomm0,ierr)
call MPI_BCAST(LGMAXps,1,MPI_INTEGER,0,Wcomm0,ierr)
call MPI_BCAST(NISps,1,MPI_INTEGER,0,Wcomm0,ierr)
call MPI_BCAST(GAMps,1,MPI_INTEGER,0,Wcomm0,ierr)
call MPI_BCAST(BETps,1,MPI_INTEGER,0,Wcomm0,ierr)
call MPI_BCAST(NU1ps,1,MPI_INTEGER,0,Wcomm0,ierr)
call MPI_BCAST(NU2ps,1,MPI_INTEGER,0,Wcomm0,ierr)
call MPI_BCAST(omeps,1,MPI_DOUBLE_PRECISION,0,Wcomm0,ierr)
call MPI_BCAST(LGMAXhd,nqhd,MPI_INTEGER,0,Wcomm0,ierr)
call MPI_BCAST(NIShd,nqhd,MPI_INTEGER,0,Wcomm0,ierr)
call MPI_BCAST(GAMhd,nqhd,MPI_INTEGER,0,Wcomm0,ierr)
call MPI_BCAST(BEThd,nqhd,MPI_INTEGER,0,Wcomm0,ierr)
call MPI_BCAST(NU1hd,nqhd,MPI_INTEGER,0,Wcomm0,ierr)
call MPI_BCAST(NU2hd,nqhd,MPI_INTEGER,0,Wcomm0,ierr)
call MPI_BCAST(omehd,nqhd,MPI_DOUBLE_PRECISION,0,Wcomm0,ierr)
call MPI_BCAST(nxmin,1,MPI_INTEGER,0,Wcomm0,ierr)
call MPI_BCAST(nymin,1,MPI_INTEGER,0,Wcomm0,ierr)
!.One minus the damping parameter in relaxations. Useful in relax routines
Oomeph=1.0_dp-omeph
Oomeps=1.0_dp-omeps
Oomehd=1.0_dp-omehd
end subroutine read_mginputs
!*********************************************************
subroutine COMMv1D(COMMv,iCOMM,act,procs,lg,BCs,ID0)
!.Create communicators (along x or y) for each (coarse) level
!.visited and where this process is active. These communicators
!.are used to exchange information in the relaxation, residue,
!.prolongation and restriction routines.
!.INPUTS
!. integer iCOMM: initial (finest) 1D communicator
!. logical dimension(:) act: active status of this rank at each grid
!. integer dimension(:) procs: number of active processes at each grid
!. integer lg: number of grids to visit
!. integer BCs: boundary conditions, TRUE for periodic
!. integer ID0: rank ID in iCOMM
!.OUTPUT
!. integer dimension(:) COMMv: array of COMMs of active processes in each grid
implicit none
include 'mpif.h'
integer(I4B), intent(in) :: iCOMM,lg,ID0
integer(I4B), intent(in), allocatable :: procs(:)
integer(I4B), intent(inout) :: COMMv(:)
logical, intent(in) :: BCs(:)
logical, intent(in) :: act(:)
integer(I4B) :: ng,i,j,k,ierr,IDs1D(1),gap,procs0, &
iCOMMgroup,tmpgroup,tmpCOMM
logical :: reorder=.TRUE.
integer(I4B), allocatable :: ranks(:)
procs0=procs(lg) !.initial number of active processes
!.Obtain MPI group of original 1D communicator
call MPI_Comm_group(iCOMM,iCOMMgroup,ierr)
do ng=lg-1,1,-1 !.Loop over coarse grids
allocate(ranks(procs(ng))) !.Ranks in this grid
k=0 !.Number of ranks in this grid
if (procs(ng)>1) then
gap=procs0/procs(ng) !.gap (+1) between active processes
!.....This does the left or bottom half of the active processes
do i=gap,procs0/2,gap
k=k+1
! IDs1D=i-1
! call MPI_CART_RANK(iCOMM,IDs1D,ranks(k),ierr)
ranks(k)=i-1 !.Rank IDs are zero-indexed
enddo
!.....This does the right or top half of the active processes
do i=procs0/2+1,procs0,gap
k=k+1
! IDs1D=i-1
! call MPI_CART_RANK(iCOMM,IDs1D,ranks(k),ierr)
ranks(k)=i-1 !.Rank IDs are zero-indexed
enddo
else
! IDs1D=iID0
! call MPI_CART_RANK(iXcomm,IDs1D,xranks(1),ierr)
ranks=ID0 !.Single active rank
endif
!...Create new group with only active processes
call MPI_GROUP_incl(iCOMMgroup,procs(ng),ranks,tmpgroup,ierr)
!...Create new communicator for active processes
!...MPI_COMM_create group available in newer compilers only and is
!...collective only for processes in the group
call MPI_COMM_create_group(iCOMM,tmpgroup,0,tmpCOMM,ierr)
! call MPI_COMM_create(iCOMM,tmpgroup,tmpCOMM,ierr)
!...Make a cartesian communicator and save in array of communicators
if (act(ng)) call MPI_CART_create(tmpCOMM,1,(/procs(ng)/),BCs,reorder,COMMv(ng),ierr)
!...Destroy temporary group and communicator to be used in next iteration
call MPI_GROUP_free(tmpgroup,ierr)
if (act(ng)) call MPI_COMM_free(tmpCOMM,ierr)
deallocate(ranks)
enddo
end subroutine COMMv1D
!*********************************************************
subroutine COMMv2D(COMMv,iCOMM,act,iprocs,jprocs,lg)
!.Create 2D communicators for each (coarse) level visited
!.INPUTS
!. integer iCOMM: initial (finest) 2D communicator
!. logical dimension(:) act: active status of this rank at each grid
!. integer dimension(:) iprocs: number of active processes along x at each grid
!. integer dimension(:) jprocs: number of active processes along y at each grid
!. integer lg: number of grids to visit
!.OUTPUT
!. integer dimension(:) COMMv: array of 2D COMMs of active processes in each grid
implicit none
include 'mpif.h'
integer(I4B), intent(in) :: iCOMM,lg
integer(I4B), intent(in), allocatable :: iprocs(:),jprocs(:)
integer(I4B), intent(inout) :: COMMv(:)
logical :: BCs(2)
logical, intent(in) :: act(:)
integer(I4B) :: ng,i,j,k,ierr,IDs2D(2),igap,jgap,iprocs0, &
jprocs0,iCOMMgroup,tmpgroup,tmpCOMM
logical :: reorder=.TRUE.
integer(I4B), allocatable :: iranks(:),jranks(:),xyranks(:)
iprocs0=iprocs(lg) !.initial number of active processes
jprocs0=jprocs(lg) !.initial number of active processes
BCs=(/.TRUE.,.FALSE./) !.Boundary conditions along y and x
!.Obtain MPI group of original 1D communicator
call MPI_Comm_group(iCOMM,iCOMMgroup,ierr)
do ng=lg-1,1,-1 !.Loop over coarse grids
allocate(iranks(iprocs(ng))) !.Ranks along x in this grid
k=0 !.Number of ranks in this grid
if (iprocs(ng)>1) then
igap=iprocs0/iprocs(ng) !.gap (+1) between active processors
!.....This does the left half (along x) of the active processes
do i=igap,iprocs0/2,igap
k=k+1
iranks(k)=i-1 !.Rank IDs are zero-indexed
enddo
!.....This does the right half (along x) of the active processes
do i=iprocs0/2+1,iprocs0,igap
k=k+1
iranks(k)=i-1 !.Rank IDs are zero-indexed
enddo
else
iranks=iID0 !.Single active rank
endif
allocate(jranks(jprocs(ng))) !.Ranks along y in this grid
k=0
if (jprocs(ng)>1) then
jgap=jprocs0/jprocs(ng) !.gap (+1) between active processors
!.....This does the bottom half (along y) of the active processes
do j=jgap,jprocs0/2,jgap
k=k+1
jranks(k)=j-1 !.Rank IDs are zero-indexed
enddo
!.....This does the top half (along y) of the active processes
do j=jprocs0/2+1,jprocs0,jgap
k=k+1
jranks(k)=j-1 !.Rank IDs are zero-indexed
enddo
else
jranks=jID0 !.Single active rank
endif
!...Determine rank IDs of all active processes in 2D communicator
allocate(xyranks(iprocs(ng)*jprocs(ng)))
k=0
do j=1,jprocs(ng)
do i=1,iprocs(ng)
k=k+1
IDs2D=(/jranks(j),iranks(i)/)
call MPI_CART_RANK(iCOMM,IDs2D,xyranks(k),ierr)
enddo
enddo
!...Create new group with only active processes
call MPI_GROUP_incl(iCOMMgroup,iprocs(ng)*jprocs(ng),xyranks,tmpgroup,ierr)
!...Create new communicator for active processes
!....MPI_COMM_create group available in newer compilers only and is
!....collective only for processes in the group
call MPI_COMM_create_group(iCOMM,tmpgroup,0,tmpCOMM,ierr)
! call MPI_COMM_create(iCOMM,tmpgroup,tmpCOMM,ierr)
!...Make a cartesian communicator and save in array of communicators
if (act(ng)) call MPI_CART_create(tmpCOMM,2,(/jprocs(ng),iprocs(ng)/), &
BCs,reorder,COMMv(ng),ierr)
!...Destroy temporary group and communicator to be used in next iteration
call MPI_GROUP_free(tmpgroup,ierr)
if (act(ng)) call MPI_COMM_free(tmpCOMM,ierr)
deallocate(iranks,jranks,xyranks)
enddo
end subroutine COMMv2D
!*********************************************************
subroutine diagNeigh(nbors,XYcomm,iprocs,jprocs,act,lg)
!.Determine the diagonal (NE,SE,SW,NW) neighbors at every grid level
!.INPUTS
!. integer dimension(:) XYcomm: 2D COMM at each grid
!. integer dimension(:) iprocs: number of active processes along x at each grid
!. integer dimension(:) jprocs: number of active processes along y at each grid
!. logical dimension(:) act: active status of this rank at each grid
!. integer lg: number of grids to visit
!.OUTPUT
!. integer dimension(:,:) nbors: array of neibor rank IDs at each grid
!. only elements 2,4,6,8 are initialized here
implicit none
include 'mpif.h'
integer(I4B), intent(in) :: XYcomm(:),iprocs(:),jprocs(:)
logical, intent(in) :: act(:)
integer(I4B), intent(in) :: lg
integer(I4B), intent(inout) :: nbors(:,:)
integer(I4B) :: ng,i,ierr,coord(2),NE(2),SE(2),SW(2),NW(2)
logical :: jnkl(2)
do ng=lg,1,-1 !.Loop over grids visited
!...Initiate them all to NULL, so non-periodic boundaries are NULL
forall(i=2:8:2) nbors(ng,i)=MPI_PROC_NULL
if (act(ng)) then
call MPI_CART_GET(XYcomm(ng),2,(/jprocs(ng),iprocs(ng)/), &
jnkl,coord,ierr) !.This ranks coordinates in 2D COMM
NE=(/coord(1)+1,coord(2)+1/) !.NorthEast neighbor coordinates
SE=(/coord(1)-1,coord(2)+1/) !.SouthEast neighbor coordinates
SW=(/coord(1)-1,coord(2)-1/) !.SouthWest neighbor coordinates
NW=(/coord(1)+1,coord(2)-1/) !.NorthWest neighbor coordinates
if (coord(1)==0) then !.Enforce periodicity in bottom boundary
SE(1)=jprocs(ng)-1
SW(1)=jprocs(ng)-1
endif
if (coord(1)==jprocs(ng)-1) then !.Enforce periodicity in top boundary
NE(1)=0
NW(1)=0
endif
!.....For non-null diagonal neigbors conver coordinates to a rank ID
if (coord(2)/=0) then
call MPI_CART_rank(XYcomm(ng),SW,nbors(ng,6),ierr) !.SW
call MPI_CART_rank(XYcomm(ng),NW,nbors(ng,8),ierr) !.NW
endif
if (coord(2)/=iprocs(ng)-1) then
call MPI_CART_rank(XYcomm(ng),NE,nbors(ng,2),ierr) !.NE
call MPI_CART_rank(XYcomm(ng),SE,nbors(ng,4),ierr) !.SE
endif
endif
enddo
end subroutine diagNeigh
!*********************************************************
subroutine pgrids(eq,lgMAX,NIS,lg,nxv,nyv,Cnx,Cny,wdc,acu,iprocs, &
jprocs,urnk,nurnk,iDiff)
!.Determine the grids visited based on coupling coefficients, when processes
!.are active or not, when to unify subdomains and how many processes remain
!.active at each grid level
!.INPUTS
!. integer eq: equation, =1 Poisson, =2 modified Helmholtz, =3 hyperdiffusion
!. integer lgMAX: max number of grids
!. integer NIS: semi-coarsening factor, =2 double spacing, =4 quadruple spacing
!. real(dp) dimension(:) iDiff: diffusion coefficients along x and y
!.OUTPUTS
!. integer dimension(:) nxv: Number of points along x of grids visited
!. integer dimension(:) nyv: Number of points along y of grids visited
!. integer dimension(:) Cnx: Number of points along x of next coarser grid at each level
!. integer dimension(:) Cny: Number of points along y of next coarser grid at each level
!. integer dimension(:) wdc: "which dimension coarsened" last at each level
!. logical dimension(:,:) acu: is process active and is subdomain unification
!. happening at this level?
!. integer dimension(:) iprocs: number of processes along x at each level
!. integer dimension(:) jprocs: number of processes along y at each level
!. integer dimension(:) urnk: ranks data is sent to/received from when unifying
!. subdomains
!. integer dimension(:,:) nurnk: number of ranks data is sent to/received from
!. when unifying subdomains
implicit none
include 'mpif.h'
integer(I4B), intent(in) :: eq, lgMAX, NIS
integer(I4B), intent(out) :: lg
integer(I4B), intent(out), allocatable :: iprocs(:),jprocs(:)
integer(I4B), allocatable, dimension(:), intent(out) :: nxv,nyv,Cnx, &
Cny,wdc
integer(I4B), allocatable, intent(out) :: nurnk(:,:)
type(Iallo2D), allocatable, intent(out) :: urnk(:)
logical, allocatable, intent(out) :: acu(:,:)
logical, allocatable :: active(:),unify(:)
real(DP), intent(in), optional :: iDiff(:)
real(DP), allocatable :: hx(:),hy(:) !.cell lengths
real(DP) :: hxs,hys !.cell lengths squared
real(DP), allocatable :: cf(:,:) !.coupling coefficients
integer(I4B) :: ng !.current grid
integer(I4B) :: igap,jgap,igap0,jgap0 !.gap between processes along x and y
logical :: oneiproc,onejproc !.true when # of active procs=1
logical :: iactive,jactive !.is process active along x/y?
logical :: iunify,junify !.does subdomain unification happen along x/y?
logical :: iunisend,junisend !.does this process send its data during unification along x/y?
integer(I4B) :: minCnx,minCny !.minimum number of points in subdomain along x/y
integer(I4B) :: inurnk,jnurnk !.number of ranks to unify subdomains with along x/y
integer(I4B) :: i,j,k,ierror
!.temporary arrays
real(DP), allocatable :: tmp(:,:),tmphx(:),tmphy(:)
real(DP) :: tmpcfx,tmpcfy
integer(I4B), allocatable, dimension(:) :: tmpx,tmpy,tmpw,tmpproc,iugap,jugap
logical, allocatable :: tmpact(:),tmpuni(:)
type(Iallo2D), allocatable :: tmpurnk(:)
integer(I4B), allocatable :: tmpnurnk(:,:)
nxf=nxL !.Fine local resolution along x
nyf=nyL !.Fine local resolution along y
nxc=nxL !.Coarse local resolution along x
nyc=nyL !.Coarse local resolution along y
lg=1 !.Initial number of grids=1
iactive=.TRUE. !.All processes active along x
jactive=.TRUE. !.All processes active along y
allocate(hx(lg),hy(lg)) !.Current fine grid spacing
allocate(active(lg),unify(lg))
active(lg)=.TRUE. !.All processes active initially
!.No unification happens at first grid. Requires initial MPI decomposition
!.to be sufficiently coarse (see nxmin and nymin in input file)
unify(lg)=.FALSE.
!.Gap between active processors. igap=1 every processor along x
!.is active, igap=2 every other processor along x is active, etc.
igap=1
jgap=1
!.Number of active processors along each direction
allocate(iprocs(lg),jprocs(lg))
iprocs(lg)=iprocs0
jprocs(lg)=jprocs0
oneiproc=.FALSE. !.true when # of active procs=1
onejproc=.FALSE.
!.Where data is sent to/received from when subdomains are unified
allocate(urnk(lg),nurnk(lg,2))
nurnk(lg,:)=1 !.number of processes participating in unification
allocate(urnk(lg)%a(nurnk(lg,1),nurnk(lg,2)))
urnk(lg)%a=myrank0 !MPI_PROC_NULL !.lg entry is not used
!.These define the coarsest grid, which depends on the equation being solved
if (eq==1 .OR. eq==2) then
minCnx=2
minCny=2
elseif (eq==3) then
#IF (HDop==2 .OR. HDop==8)
minCnx=3
minCny=3
#ELIF (HDop==3 .OR. HDop==12)
minCnx=4
minCny=4
#ELIF (HDop==4 .OR. HDop==6)
minCnx=2
minCny=2
#ENDIF
endif
do
!...Begin assuming no unification happens then check later
iunify=.FALSE.
junify=.FALSE.
iunisend=.FALSE.
junisend=.FALSE.
!...Current cell size along each direction
hx(1)=bLx/dble(iprocs(1)*nxf)
hy(1)=bLy/dble(jprocs(1)*nyf)
if ((nxc*iprocs(1) == minCnx) .OR. (nyc*jprocs(1) == minCny) &
.OR. (lg==lgMAX)) exit !.reached last grid
lg=lg+1 !.one more coarser grid
!...Current fine grid spacing (squared)
hxs=hx(1)**2
hys=hy(1)**2
!...This allocating and deallocating is just so that quantities
!...are stored in the right order in memory, which is finest
!...last and coarsest first.
if (lg>2) then
!.....Temp coupling coefficients from previous grid
allocate(tmp(lg-1,2))
tmp=cf
deallocate(cf)
!.....Temp previous grid dimensions
allocate(tmpx(lg-1),tmpy(lg-1))
tmpx=nxv
tmpy=nyv
deallocate(nxv,nyv)
!.....Temp previous indicator of which dimension was coarsened
allocate(tmpw(lg-1))
tmpw=wdc
deallocate(wdc)
endif
!...Temp previous cell sizes
allocate(tmphx(lg-1),tmphy(lg-1))
tmphx=hx
tmphy=hy
deallocate(hx,hy)
!...Allocate arrays with new number of grids
allocate(cf(lg,2))
allocate(nxv(lg),nyv(lg))
allocate(wdc(lg))
allocate(hx(lg),hy(lg))
wdc=0 !.only affects the final wdc(1), which isn't used
!...Populate arrays with previous values stored in tmp arrays
if (lg>2) then
cf(3:lg,:)=tmp(2:lg-1,:)
deallocate(tmp)
nxv(3:lg)=tmpx(2:lg-1)
nyv(3:lg)=tmpy(2:lg-1)
deallocate(tmpx,tmpy)
wdc(2:lg)=tmpw
wdc(2)=wdc(3) !.temporarily
deallocate(tmpw)
endif
hx(2:lg)=tmphx
hy(2:lg)=tmphy
deallocate(tmphx,tmphy)
!...Store current fine grid dimensions
nxv(2)=nxf
nyv(2)=nyf
if (lg>2 .AND. wdc(2)==0) then ! .OR. (wdc() .ne. 0)) then
!.....If already reached a grid with balanced coupling coefficients
!.....coarsen using standard coarsening from now on
cf(2,1)=1.0_dp
cf(2,2)=1.0_dp
else
!.....Otherwise compute the coupling factors + determine coarsening
!.....This identifies along which dimension is the coupling in the
!.....elliptic equation strongest. It is a combination of the difference
!.....in dx vs dy and the coefficients in the equation.
!.....Two options here
!!____________________ Option A __________________________________________!
!!_____ For the Poisson eqn this seems like a better informed approach ____!
! allocate(tmp(nxf,nyf))
!........tmpD here is the coarsened Diff at every level (currently computed
!........after this subroutine is called)
! do i=2,nxf-1
!!.........Original, based on exact entries in matrix
!! tmp(i,:)=(tmpD(ng-lg+2)%a(i+1,:,1)+tmpD(ng-lg+2)%a(i-1,:,1) &
!! +2.0_dp*tmpD(ng-lg+2)%a(i,:,1))/hxs
!!.........This one requires less communication
! tmp(i,:)=(tmpD(ng-lg+2)%a(i+1,:,1)+tmpD(ng-lg+2)%a(i-1,:,1) &
! +2.0_dp*tmpD(ng-lg+2)%a(i,:,1))/hxs
! enddo
! cf(2,1)=maxval(tmp(2:nxf-1,:))
!! cf(2,1)=sum(tmp(2:nxf-1,:))/dble((nxf-2)*nyf)
! do j=2,nyf-1
! tmp(:,j)=(tmpD(ng-lg+2)%a(:,j+1,1)+tmpD(ng-lg+2)%a(:,j-1,1) &
! +2.0_dp*tmpD(ng-lg+2)%a(:,j,1))/hys
! enddo
! cf(2,2)=maxval(tmp(:,2:nyf-1))
!! cf(2,2)=sum(tmp(:,2:nyf-1))/dble(nxf*(nyf-2))
! deallocate(tmp)
!!____________________ Option B __________________________________________!
!.....It appears that using the same coupling coefficient for both the
!.....Poisson, modified Helmholtz and Laplacian-based hyperdiffusion
!.....works fine (see below). Then each processor can establish the
!.....communication patterns independently because the matrix coefficients
!.....are scalar constants.
#IF (HDop==8 .OR. HDop==12)
if (eq==3) then
cf(2,1)=iDiff(1)/(hxs**(HDop/4))
cf(2,2)=iDiff(2)/(hys**(HDop/4))
else
cf(2,1)=1.0_dp/hxs
cf(2,2)=1.0_dp/hys
endif
#ELSE
cf(2,1)=1.0_dp/hxs
cf(2,2)=1.0_dp/hys
#ENDIF
!!__________________________________________________________________________!
endif
if (NIS==2) then
!.....Doubling grid spacing at every level
if (cf(2,1)/cf(2,2)>1.30_dp) then !. 1.3 from Zubair 2007, requires testing
!.......Semi-coarsen along x
nxc=nxf/2
nyc=nyf
wdc(2)=1
elseif (cf(2,2)/cf(2,1)>1.30_dp) then !. 1.3 from Zubair 2007, requires testing
!.......Semi-coarsen along y
nxc=nxf
nyc=nyf/2
wdc(2)=2
else
!.......Standard coarsening
nxc=nxf/2
nyc=nyf/2
wdc(2)=0
endif
else !if (NIS==4) then
!.....Quadrupling grid spacing until problem is isotropic, then doubling
if (cf(2,1)/cf(2,2)>1.01_dp) then !. 1.3 from Zubair 2007, requires testing
!.......Semi-coarsen along x
nxc=nxf/4
nyc=nyf
wdc(2)=1
!.......If quadrupling the grid spacing overcoarsens the grid along that
!.......dimension so that the next step would require quadrupling along
!.......the other, then just double the grid spacing
#IF (HDop==8 .OR. HDop==12)
if (eq==3) then
tmpcfx=iDiff(1)/((bLx/(dble(nxc)*iprocs(1)))**HDopd4)
tmpcfy=iDiff(2)/((bLy/(dble(nyc)*jprocs(1)))**HDopd4)
else
tmpcfx=1.0_dp/(bLx/(dble(nxc)*iprocs(1)))**2
tmpcfy=1.0_dp/(bLy/(dble(nyc)*jprocs(1)))**2
endif
#ELSE
tmpcfx=1.0_dp/(bLx/(dble(nxc)*iprocs(1)))**2
tmpcfy=1.0_dp/(bLy/(dble(nyc)*jprocs(1)))**2
#ENDIF
if (tmpcfy/tmpcfx>1.01_dp) nxc=nxc*2
elseif (cf(2,2)/cf(2,1)>1.01_dp) then !. 1.3 from Zubair 2007, requires testing
!.......Semi-coarsen along y
nxc=nxf
nyc=nyf/4
wdc(2)=2
!.......If quadrupling the grid spacing overcoarsens the grid along that
!.......dimension so that the next step would require quadrupling along
!.......the other, then just double the grid spacing
#IF (HDop==8 .OR. HDop==12)
if (eq==3) then
tmpcfx=iDiff(1)/((bLx/(dble(nxc)*iprocs(1)))**HDopd4)
tmpcfy=iDiff(2)/((bLy/(dble(nyc)*jprocs(1)))**HDopd4)
else
tmpcfx=1.0_dp/(bLx/(dble(nxc)*iprocs(1)))**2
tmpcfy=1.0_dp/(bLy/(dble(nyc)*jprocs(1)))**2
endif
#ELSE
tmpcfx=1.0_dp/(bLx/(dble(nxc)*iprocs(1)))**2
tmpcfy=1.0_dp/(bLy/(dble(nyc)*jprocs(1)))**2
#ENDIF
if (tmpcfx/tmpcfy>1.01_dp) nyc=nyc*2
else
!.......Standard coarsening
nxc=nxf/2
nyc=nyf/2
wdc(2)=0
endif
endif
if (lg>2) then
!.....Temporarily repurpose tmphx/y to save coarse grid number of points
allocate(tmphx(lg-1),tmphy(lg-1))
tmphx=Cnx
tmphy=Cny
deallocate(Cnx,Cny)
endif
!...Save local coarse grid dimensions
allocate(Cnx(lg),Cny(lg))
Cnx(2)=nxc
Cny(2)=nyc