-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathPlaneRemesh.f90
1103 lines (1011 loc) · 43.1 KB
/
PlaneRemesh.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
module PlaneRemeshModule
!------------------------------------------------------------------------------
! Lagrangian Particle / Panel Method
!------------------------------------------------------------------------------
!
!> @author
!> Peter Bosler, Department of Mathematics, University of Michigan
!
!> @defgroup PlaneRemesh Plane Remesh
!> Data types, subroutines, and functions for adaptive remeshing planar particle/panel meshes.
!
!
! DESCRIPTION:
!> @file
!> Data types, subroutines, and functions for adaptive remeshing planar particle/panel meshes.
!
!------------------------------------------------------------------------------
use NumberKindsModule
use LoggerModule
use PlaneGeomModule
use ParticlesModule
use EdgesModule
use PanelsModule
use PlaneMeshModule
use PlaneVorticityModule
use PlaneTracerModule
use BIVARInterfaceModule
use BIVARModule
implicit none
private
public RemeshSetup
public New, Delete
public InitialRefinement
public LagrangianRemeshToInitialTime, LagrangianRemeshToReferenceTime
!
!----------------
! Types and module constants
!----------------
!
type RemeshSetup
logical(klog) :: uniformMesh ! if TRUE, no AMR
logical(klog) :: vorticityRefine ! if TRUE, remeshing will apply vorticity AMR criteria
real(kreal) :: maxCircTol
real(kreal) :: vortVarTol
logical(klog) :: flowMapRefine ! if TRUE, remeshing will apply flow map AMR criteria
real(kreal) :: lagVarTol
logical(klog) :: tracerRefine ! if TRUE, remeshing will apply tracer AMR criteria
integer(kint) :: tracerID
real(kreal) :: maxMassTol
real(kreal) :: massVarTol
integer(kint) :: refinementLimit
end type
!
!----------------
! Interfaces
!----------------
!
interface New
module procedure NewPrivateAll
module procedure NewPrivateVorticity
module procedure NewPrivateVorticityAndFlowMap
module procedure NewPrivateTracer
module procedure NewPrivateNoAMR
end interface
interface Delete
module procedure DeletePrivate
end interface
interface InitialRefinement
module procedure InitialRefinementPrivate
end interface
interface LagrangianRemeshToInitialTime
module procedure LagrangianRemeshToInitialTimePrivate
end interface
interface LagrangianRemeshToReferenceTime
module procedure LagrangianRemeshToReferenceTimePrivate
end interface
interface
subroutine SetVorticityOnMesh( genMesh, genVort)
use PlaneMeshModule
use PlaneVorticityModule
implicit none
type(PlaneMesh), intent(inout) :: genMesh
type(VorticitySetup), intent(in) :: genVort
end subroutine
end interface
interface
subroutine SetTracerOnMesh( genMesh, genTracer )
use PlaneMeshModule
use PlaneTracerModule
implicit none
type(PlaneMesh), intent(inout) :: genMesh
type(TracerSetup), intent(in) :: genTracer
end subroutine
end interface
!
!----------------
! Logging
!----------------
!
logical(klog), save :: logInit = .FALSE.
type(Logger) :: log
character(len=28), save :: logKey = 'PlaneRemesh'
integer(kint), parameter :: logLevel = TRACE_LOGGING_LEVEL
character(len=128) :: logString
character(len=24) :: formatString
contains
!
!----------------
! Standard methods : Constructor / Destructor
!----------------
!
subroutine NewPrivateAll(self, maxCircTol, vortVarTol, lagVarTol, tracerID, maxMassTol, massVarTol, limit)
type(RemeshSetup), intent(out) :: self
real(kreal), intent(in) :: maxCircTol, vortVarTol, lagVarTol
integer(kint), intent(in) :: tracerID
real(kreal), intent(in) :: maxMassTol, massVarTol
integer(kint), intent(in) :: limit
if ( .NOT. logInit ) call InitLogger(log, procRank)
self%uniformMesh = .FALSE.
self%vorticityRefine = .TRUE.
self%maxCircTol = maxCircTol
self%vortVarTol = vortVarTol
self%flowMapRefine = .TRUE.
self%lagVarTol = lagVarTol
self%tracerRefine = .TRUE.
self%tracerID = tracerID
self%maxMassTol = maxMassTol
self%massVarTol = massVarTol
self%refinementLimit = limit
end subroutine
subroutine NewPrivateVorticity(self, maxCircTol, vortVarTol, limit)
type(RemeshSetup), intent(out) :: self
real(kreal), intent(in) :: maxCircTol, vortVarTol
integer(kint), intent(in) :: limit
if ( .NOT. logInit ) call InitLogger(log, procRank)
self%uniformMesh = .FALSE.
self%vorticityRefine = .TRUE.
self%maxCircTol = maxCircTol
self%vortVarTol = vortVarTol
self%flowMapRefine = .FALSE.
self%lagVarTol = 0.0_kreal
self%tracerRefine = .FALSE.
self%tracerID = 0
self%maxMassTol = 0.0_kreal
self%massVarTol = 0.0_kreal
self%refinementLimit = limit
end subroutine
subroutine NewPrivateVorticityAndFlowMap(self, maxCircTol, vortVarTol, lagVarTol, limit)
type(RemeshSetup), intent(out) :: self
real(kreal), intent(in) :: maxCircTol, vortVarTol, lagVarTol
integer(kint), intent(in) :: limit
if ( .NOT. logInit ) call InitLogger(log, procRank)
self%uniformMesh = .FALSE.
self%vorticityRefine = .TRUE.
self%maxCircTol = maxCircTol
self%vortVarTol = vortVarTol
self%flowMapRefine = .TRUE.
self%lagVarTol = lagVarTol
self%tracerRefine = .FALSE.
self%tracerID = 0
self%maxMassTol = 0.0_kreal
self%massVarTol = 0.0_kreal
self%refinementLimit = limit
end subroutine
subroutine NewPrivateTracer(self, tracerID, maxMassTol, massVarTol, limit)
type(RemeshSetup), intent(inout) :: self
integer(kint), intent(in) :: tracerID
real(kreal), intent(in) :: maxMassTol, massVarTol
integer(kint), intent(in) :: limit
if ( .NOT. logInit ) call InitLogger(log, procRank)
self%uniformMesh = .FALSE.
self%vorticityRefine = .FALSE.
self%maxCircTol = 0.0_kreal
self%vortVarTol = 0.0_kreal
self%flowMapRefine = .FALSE.
self%lagVarTol = 0.0_kreal
self%tracerRefine = .TRUE.
self%tracerID = tracerID
self%maxMassTol = maxMassTol
self%massVarTol = massVarTol
self%refinementLimit = limit
end subroutine
subroutine NewPrivateNoAMR(self)
type(RemeshSetup), intent(inout) :: self
if ( .NOT. logInit ) call InitLogger(log, procRank)
self%uniformMesh = .TRUE.
self%vorticityRefine = .FALSE.
self%maxCircTol = 0.0_kreal
self%vortVarTol = 0.0_kreal
self%flowmapRefine = .FALSE.
self%lagVarTol = 0.0_kreal
self%tracerRefine = .FALSE.
self%tracerID = 0
self%maxMassTol = 0.0_kreal
self%massVarTol = 0.0_kreal
self%refinementLimit = 0
end subroutine
subroutine DeletePrivate(self)
type(RemeshSetup), intent(inout) :: self
self%uniformMesh = .TRUE.
self%vorticityRefine = .FALSE.
self%flowMapRefine = .FALSE.
self%tracerRefine = .FALSE.
end subroutine
!
!----------------
! Public functions
!----------------
!
subroutine InitialRefinementPrivate( aMesh, remeshData, updateTracerOnMesh, tracerDefinition, &
updateVorticityOnMesh, vorticityDefinition)
type(PlaneMesh), intent(inout) :: aMesh
type(RemeshSetup), intent(in) :: remeshData
procedure(SetTracerOnMesh) :: updateTracerOnMesh
type(TracerSetup), intent(in) :: tracerDefinition
procedure(SetVorticityOnMesh) :: updateVorticityOnMesh
type(VorticitySetup), intent(in) :: vorticityDefinition
!
integer(kint) :: refineCount, spaceLeft, counters(5), j
integer(kint) :: startIndex, nOldPanels, amrLoopCounter, limit
logical(klog) :: keepGoing
logical(klog), allocatable :: refineFlag(:)
type(Panels), pointer :: aPanels
if ( remeshData%uniformMesh ) then
call LogMessage(log, WARNING_LOGGING_LEVEL, logKey, ' InitialRefinement WARNING : uniform mesh = no AMR.')
return
endif
if ( remeshData%refinementLimit == 0 ) then
call LogMessage(log, WARNING_LOGGING_LEVEL, logKey, ' InitialRefinement WARNING : amr limit = 0.')
return
endif
call StartSection(log, 'InitialRefinement')
aPanels => aMesh%panels
allocate(refineFlag(aPanels%N_Max))
refineFlag = .FALSE.
keepGoing = .FALSE.
!
! apply refinement criteria
!
limit = remeshData%refinementLimit
startIndex = 1
counters = 0
if ( remeshData%vorticityRefine ) then
!limit = max( remeshData%refinementLimit, limit)
call FlagPanelsForMaxCirculationRefinement(refineFlag, aMesh, remeshData, startIndex, counters(1))
call FlagPanelsForVorticityVariationRefinement(refineFlag, aMesh, remeshData, startIndex, counters(2))
endif
if ( remeshData%flowMapRefine ) then
!limit = max( remeshData%refinementLimit, limit)
call FlagPanelsForFlowMapRefinement(refineFlag, aMesh, remeshData, startIndex, counters(3))
endif
if ( remeshData%tracerRefine) then
!limit = max( remeshData%refinementLimit, limit)
call FlagPanelsForTracerMassRefinement(refineFlag, aMesh, remeshData, startIndex, counters(4))
call FlagPanelsForTracerVariationRefinement(refineFlag, aMesh, remeshData, startIndex, counters(5))
endif
refineCount = count(refineFlag)
spaceLeft = aPanels%N_Max - aPanels%N
if ( refineCount > 0 ) then
if ( spaceleft/4 > refineCount ) then
keepGoing = .TRUE.
amrLoopCounter = 0
do while (keepGoing)
amrLoopCounter = amrLoopCounter + 1
write(logstring,'(A,I2,A,I8,A)') 'AMR loop ', amrLoopCounter, ' : refining ', refineCount,' panels.'
call LogMessage(log,TRACE_LOGGING_LEVEL,'InitRefine : ',logString)
!
! divide flagged panels
!
nOldPanels = aPanels%N
do j = startIndex, nOldPanels
if ( refineFlag(j) ) then
call DividePanel(aMesh, j)
refineFlag(j) = .FALSE.
endif
enddo
!
! TO DO : ensure adjacent panels do not differ by more than one level of refinement
!
! set data on refined panels
if ( aMesh%nTracer > 0 ) call UpdateTracerOnMesh(aMesh, tracerDefinition)
call UpdateVorticityOnMesh(amesh, vorticityDefinition)
!
! prevent too much refinement
!
if ( amrLoopCounter >= remeshData%refinementLimit ) then
keepGoing = .FALSE.
call LogMessage(log, WARNING_LOGGING_LEVEL, 'InitRefine WARNING : ',' refinement limit reached.')
call LogMessage(log, TRACE_LOGGING_LEVEL,'InitRefine : nActive = ', aPanels%N_Active)
write(logstring,'(A, I8, A)') ' max circulation criterion triggered ', counters(1), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'InitRefine : ',trim(logstring))
write(logstring,'(A, I8, A)') ' vort. variation criterion triggered ', counters(2), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'InitRefine : ',trim(logstring))
write(logstring,'(A, I8, A)') 'flwmap variation criterion triggered ', counters(3), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'InitRefine : ',trim(logstring))
write(logstring,'(A, I8, A)') ' tracermax criterion triggered ', counters(4), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'InitRefine : ',trim(logstring))
write(logstring,'(A, I8, A)') 'tracer variation criterion triggered ', counters(5), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'InitRefine : ',trim(logstring))
else
!
! apply refinement criteria to new panels
!
startIndex = nOldPanels + 1
nOldPanels = aPanels%N
if ( remeshData%vorticityRefine ) then
call FlagPanelsForMaxCirculationRefinement(refineFlag, &
aMesh, remeshData, startIndex, counters(1))
call FlagPanelsForVorticityVariationRefinement(refineFlag, &
aMesh, remeshData, startIndex, counters(2))
endif
if ( remeshData%flowMapRefine ) then
call FlagPanelsForFlowMapRefinement(refineFlag, aMesh, remeshData, startIndex, counters(3))
endif
if ( remeshData%tracerRefine ) then
call FlagPanelsForTracerMassRefinement(refineFlag, &
aMesh, remeshData, startIndex, counters(4))
call FlagPanelsForTracerVariationRefinement(refineFlag, &
aMesh, remeshData, startIndex, counters(5))
endif
refineCount = count(refineFlag)
spaceleft = aPanels%N_Max - aPanels%N
if ( refineCount == 0 ) then
!
! stop if refinement is not needed
!
call LogMessage(log, TRACE_LOGGING_LEVEL, 'InitRefine : ', 'refinement converged.')
call LogMessage(log, TRACE_LOGGING_LEVEL,'InitRefine : nActive = ', aPanels%N_Active)
write(logstring,'(A, I8, A)') ' max circulation criterion triggered ', counters(1), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'InitRefine : ',trim(logstring))
write(logstring,'(A, I8, A)') ' vort. variation criterion triggered ', counters(2), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'InitRefine : ',trim(logstring))
write(logstring,'(A, I8, A)') 'flwmap variation criterion triggered ', counters(3), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'InitRefine : ',trim(logstring))
write(logstring,'(A, I8, A)') ' tracermax criterion triggered ', counters(4), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'InitRefine : ',trim(logstring))
write(logstring,'(A, I8, A)') 'tracer variation criterion triggered ', counters(5), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'InitRefine : ',trim(logstring))
keepGoing = .FALSE.
elseif ( spaceLeft / 4 < refineCount ) then
!
! stop if not enough memory
!
call LogMessage(log, WARNING_LOGGING_LEVEL, 'InitRefine WARNING : ',&
'not enough memory to continue AMR.')
call LogMessage(log, TRACE_LOGGING_LEVEL,'InitRefine : nActive = ', aPanels%N_Active)
write(logstring,'(A, I8, A)') ' max circulation criterion triggered ', counters(1), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'InitRefine : ',trim(logstring))
write(logstring,'(A, I8, A)') ' vort. variation criterion triggered ', counters(2), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'InitRefine : ',trim(logstring))
write(logstring,'(A, I8, A)') 'flwmap variation criterion triggered ', counters(3), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'InitRefine : ',trim(logstring))
write(logstring,'(A, I8, A)') ' tracermax criterion triggered ', counters(4), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'InitRefine : ',trim(logstring))
write(logstring,'(A, I8, A)') 'tracer variation criterion triggered ', counters(5), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'InitRefine : ',trim(logstring))
keepGoing = .FALSE.
endif ! stopping criteria met
endif! below refinement limit
enddo! while keepgoing
else ! not enough memory
call LogMessage(log, WARNING_LOGGING_LEVEL, 'InitRefine WARNING : ', 'not enough memory.')
endif
else ! no refinement needed
call LogMessage(log, TRACE_LOGGING_LEVEL, 'InitRefine : ',' no refinement necessary.')
endif
call EndSection(log)
deallocate(refineFlag)
end subroutine
subroutine LagrangianRemeshToInitialTimePrivate(aMesh, remesh, setVorticity, vorticityDef, setTracer, tracerDef)
type(PlaneMesh), intent(inout) :: aMesh
type(RemeshSetup), intent(in) :: remesh
procedure(SetVorticityOnMesh) :: setVorticity
type(VorticitySetup), intent(in) :: vorticityDef
procedure(SetTracerOnMesh) :: setTracer
type(TracerSetup), intent(in) :: tracerDef
!
type(BIVARSetup) :: bivarData
type(PlaneMesh) :: newMesh
type(Particles), pointer :: newParticles
type(Edges), pointer :: newEdges
type(Panels), pointer :: newPanels
integer(kint) :: j, amrLoopCounter, counters(5), md
logical(klog), allocatable :: refineFlag(:)
logical(klog) :: keepGoing
integer(kint) :: startIndex, nOldPanels, nOldParticles, refineCount, spaceLeft
integer(kint), allocatable :: integerWorkParticles(:), integerWorkPanels(:)
nullify(newParticles)
nullify(newEdges)
nullify(newPanels)
keepGoing = .FALSE.
amrLoopCounter = 0
counters = 0
call LogMessage(log, DEBUG_LOGGING_LEVEL, logkey, ' entering LagrangianRemeshToInitialTime.')
!
! build a new uniform mesh
!
call New(newMesh, aMesh%initNest, aMesh%amr, aMesh%nTracer)
call InitializeRectangle(newMesh, MinX(aMesh), MaxX(aMesh), MinY(aMesh), MaxY(aMesh), aMesh%boundaryType)
newParticles => newMesh%particles
newEdges => newMesh%edges
newPanels => newMesh%panels
!
! set old mesh as interpolation data source
!
call New(bivarData, aMesh)
allocate(integerWorkParticles( 31*bivarData%n + newParticles%N_Max))
allocate(integerWorkPanels( 31*bivarData%n + newPanels%N_Max))
integerWorkParticles = 0
integerWorkPanels = 0
!
! interpolate Lagrangian parameter from old mesh to new mesh
!
md = 1
call IDBVIP( md, bivarData%n, bivarData%x, bivarData%y, bivarData%x0, &
newParticles%N, newParticles%x(1,1:newParticles%N), newParticles%x(2,1:newParticles%N), &
newParticles%x0(1,1:newParticles%N), integerWorkParticles, bivarData%realWork)
md = 3
call IDBVIP( md, bivarData%n, bivarData%x, bivarData%y, bivarData%y0, &
newParticles%N, newParticles%x(1,1:newParticles%N), newParticles%x(2,1:newParticles%N), &
newParticles%x0(2,1:newParticles%N), integerWorkParticles, bivarData%realWork)
md = 1
call IDBVIP( md, bivarData%n, bivarData%x, bivarData%y, bivarData%x0, &
newPanels%N, newPanels%x(1,1:newPanels%N), newPanels%x(2,1:newPanels%N), &
newPanels%x0(1,1:newPanels%N), integerWorkPanels, bivarData%realWork)
md = 3
call IDBVIP( md, bivarData%n, bivarData%x, bivarData%y, bivarData%y0, &
newPanels%N, newPanels%x(1,1:newPanels%N), newPanels%x(2,1:newPanels%N), &
newPanels%x0(2,1:newPanels%N), integerWorkPanels, bivarData%realWork)
!
! set tracer values on new mesh
!
if ( aMesh%nTracer > 0 ) call SetTracer(newMesh, tracerDef)
!
! set vorticity values on new mesh
!
call SetVorticity(newMesh, vorticityDef)
call LogMessage(log, DEBUG_LOGGING_LEVEL, logkey, ' LagRemeshInitTime : new uniform mesh ready.')
!
! AMR
!
if ( aMesh%AMR > 0 ) then
call StartSection(log,'LagrangianRemeshToInitialTime AMR :')
allocate(refineFlag(newPanels%N_Max))
refineFlag = .FALSE.
startIndex = 1
if ( remesh%vorticityRefine ) then
call FlagPanelsForMaxCirculationRefinement(refineFlag, newMesh, remesh, startIndex, counters(1) )
call FlagPanelsForVorticityVariationRefinement(refineFlag, newMesh, remesh, startIndex, counters(2) )
endif
if ( remesh%flowMapRefine ) then
call FlagPanelsForFlowMapRefinement(refineFlag, newMesh, remesh, startIndex, counters(3))
endif
if ( remesh%tracerRefine ) then
call FlagPanelsForTracerMassRefinement(refineFlag, newMesh, remesh, startIndex, counters(4))
call FlagPanelsForTracerVariationRefinement(refineFlag, newMesh, remesh, startIndex, counters(5))
endif
refineCount = count(refineFlag)
spaceLeft = newPanels%N_Max - newPanels%N
if ( refineCount > 0 ) then
if ( spaceLeft/4 > refineCount ) then
keepGoing = .TRUE.
do while (keepGoing)
amrLoopCounter = amrLoopCounter + 1
write(logString,'(A,I2,A,I8,A)') 'AMR Loop ', amrLoopCounter,&
' : refining ', refineCount, ' panels.'
call LogMessage(log, TRACE_LOGGING_LEVEL, 'LagRemeshInitTime : ',trim(logstring))
!
! divide flagged panels
!
nOldPanels = newPanels%N
nOldParticles = newParticles%N
do j = startIndex, nOldPanels
if ( refineFlag(j) ) then
call DividePanel(newMesh, j)
refineFlag(j) = .FALSE.
endif
enddo
!
! TO DO : ensure adjacent panels differ by no more than 1 level of refinement
!
!
! interpolate Lagrangian parameter to new panels
!
md = 1
call IDBVIP( md, bivarData%n, bivarData%x, bivarData%y, bivarData%x0, &
newParticles%N - nOldParticles, newParticles%x(1,nOldParticles+1:newParticles%N), &
newParticles%x(2,nOldParticles+1:newParticles%N), &
newParticles%x0(1,nOldParticles+1:newParticles%N), integerWorkParticles, bivarData%realWork)
md = 3
call IDBVIP( md, bivarData%n, bivarData%x, bivarData%y, bivarData%y0, &
newParticles%N - nOldParticles, newParticles%x(1,nOldParticles+1:newParticles%N), &
newParticles%x(2,nOldParticles+1:newParticles%N), &
newParticles%x0(2,nOldParticles+1:newParticles%N), integerWorkParticles, bivarData%realWork)
md = 1
call IDBVIP( md, bivarData%n, bivarData%x, bivarData%y, bivarData%x0, &
newPanels%N - nOldPanels, newPanels%x(1,nOldPanels+1:newPanels%N), &
newPanels%x(2,nOldPanels+1:newPanels%N), &
newPanels%x0(1,nOldPanels+1:newPanels%N), integerWorkPanels, bivarData%realWork)
md = 3
call IDBVIP( md, bivarData%n, bivarData%x, bivarData%y, bivarData%y0, &
newPanels%N - nOldPanels, newPanels%x(1,nOldPanels+1:newPanels%N), &
newPanels%x(2,nOldPanels+1:newPanels%N), &
newPanels%x0(2,nOldPanels+1:newPanels%N), integerWorkPanels, bivarData%realWork)
!
! set data on new mesh
!
if ( amesh%nTracer > 0 ) call SetTracer(newMesh, tracerDef)
call SetVorticity(newMesh, vorticityDef)
if ( amrLoopCounter >= remesh%refinementLimit ) then
!
! prevent too much refinement
!
keepGoing = .FALSE.
call LogMessage(log, WARNING_LOGGING_LEVEL, 'LagRemeshInitTime WARNING : ',&
'refinement limit reached.')
call LogMessage(log, TRACE_LOGGING_LEVEL, 'LagRemeshInitTime : nActive = ', newPanels%N_Active)
write(logstring,'(A, I8, A)') ' max circulation criterion triggered ', counters(1), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'LagRemeshInitTime : ',trim(logstring))
write(logstring,'(A, I8, A)') ' vort. variation criterion triggered ', counters(2), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'LagRemeshInitTime : ',trim(logstring))
write(logstring,'(A, I8, A)') 'flwmap variation criterion triggered ', counters(3), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'LagRemeshInitTime : ',trim(logstring))
write(logstring,'(A, I8, A)') ' tracermax criterion triggered ', counters(4), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'LagRemeshInitTime : ',trim(logstring))
write(logstring,'(A, I8, A)') 'tracer variation criterion triggered ', counters(5), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'LagRemeshInitTime : ',trim(logstring))
else
!
! apply refinement criteria to new panels
!
startIndex = nOldPanels + 1
nOldPanels = newPanels%N
if ( remesh%vorticityRefine ) then
call FlagPanelsForMaxCirculationRefinement(refineFlag, &
newMesh, remesh, startIndex, counters(1) )
call FlagPanelsForVorticityVariationRefinement(refineFlag, &
newMesh, remesh, startIndex, counters(2) )
endif
if ( remesh%flowMapRefine ) then
call FlagPanelsForFlowMapRefinement(refineFlag, newMesh, remesh, startIndex, counters(3))
endif
if ( remesh%tracerRefine ) then
call FlagPanelsForTracerMassRefinement(refineFlag, &
newMesh, remesh, startIndex, counters(4))
call FlagPanelsForTracerVariationRefinement(refineFlag, &
newMesh, remesh, startIndex, counters(5))
endif
refineCount = count(refineFlag)
spaceLeft = newPanels%N_Max - newPanels%N
if ( refineCount == 0 ) then
keepGoing = .FALSE.
call LogMessage(log, TRACE_LOGGING_LEVEL, 'LagRemeshInitTime : ', &
' refinement converged.')
call LogMessage(log, TRACE_LOGGING_LEVEL, 'LagRemeshInitTime : nActive = ',&
newPanels%N_Active)
write(logstring,'(A, I8, A)') ' max circulation criterion triggered ', &
counters(1), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'LagRemeshInitTime : ',trim(logstring))
write(logstring,'(A, I8, A)') ' vort. variation criterion triggered ',&
counters(2), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'LagRemeshInitTime : ',trim(logstring))
write(logstring,'(A, I8, A)') 'flwmap variation criterion triggered ', &
counters(3), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'LagRemeshInitTime : ',trim(logstring))
write(logstring,'(A, I8, A)') ' tracermax criterion triggered ',&
counters(4), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'LagRemeshInitTime : ',trim(logstring))
write(logstring,'(A, I8, A)') 'tracer variation criterion triggered ', &
counters(5), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'LagRemeshInitTime : ',trim(logstring))
elseif ( spaceLeft / 4 < refineCount ) then
keepGoing = .FALSE.
call LogMessage(log, WARNING_LOGGING_LEVEL, 'LagRemeshInitTime WARNING : ', &
'not enough memory to continue AMR.')
call LogMessage(log, TRACE_LOGGING_LEVEL, 'LagRemeshInitTime : nActive = ',&
newPanels%N_Active)
write(logstring,'(A, I8, A)') ' max circulation criterion triggered ',&
counters(1), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'LagRemeshInitTime : ',trim(logstring))
write(logstring,'(A, I8, A)') ' vort. variation criterion triggered ', &
counters(2), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'LagRemeshInitTime : ',trim(logstring))
write(logstring,'(A, I8, A)') 'flwmap variation criterion triggered ',&
counters(3), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'LagRemeshInitTime : ',trim(logstring))
write(logstring,'(A, I8, A)') ' tracermax criterion triggered ', &
counters(4), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'LagRemeshInitTime : ',trim(logstring))
write(logstring,'(A, I8, A)') 'tracer variation criterion triggered ',&
counters(5), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'LagRemeshInitTime : ',trim(logstring))
endif ! stopping criteria met
endif ! limit reached
enddo! while keepGoing
else ! not enough space left in memory
call LogMessage(log, WARNING_LOGGING_LEVEL, 'LagRemeshInitTime WARNING : ',&
'not enough memory for AMR.')
endif
else ! refineCount == 0
call LogMessage(log, TRACE_LOGGING_LEVEL, 'LagRemeshInitTime : ', ' no refinement necessary.')
endif
call EndSection(log)
deallocate(refineFlag)
endif! AMR
!
! replace old mesh with new mesh
!
call Copy(aMesh, newMesh)
!
! clean up
!
deallocate(integerWorkParticles)
deallocate(integerWorkPanels)
call Delete(bivarData)
call Delete(newMesh)
end subroutine
subroutine LagrangianRemeshToReferenceTimePrivate(aMesh, reference, remesh)
type(PlaneMesh), intent(inout) :: aMesh
type(BIVARSetup), intent(in) :: reference
type(RemeshSetup), intent(in) :: remesh
!
type(PlaneMesh) :: newMesh
type(Particles), pointer :: newParticles
type(Edges), pointer :: newEdges
type(Panels), pointer :: newPanels
type(BIVARSetup) :: bivarData
integer(kint) :: j, k, amrLoopCounter, counters(5), md
logical(klog), allocatable :: refineFlag(:)
logical(klog) :: keepGoing
integer(kint) :: startIndex, nOldPanels, nOldParticles, refineCount, spaceLeft
integer(kint), allocatable :: integerWorkParticles(:), integerWorkPanels(:)
nullify(newparticles)
nullify(newEdges)
nullify(newPanels)
amrLoopCounter = 0
counters = 0
md = 1
startIndex = 1
keepGoing = .FALSE.
call LogMessage(log, DEBUG_LOGGING_LEVEL, logkey, ' entering LagrangianRemeshToReferenceTime')
!
! build a new uniform mesh
!
call New(newMesh, aMesh%initNest, aMesh%amr, aMesh%ntracer)
call InitializeRectangle(newMesh, MinX(aMesh), MaxX(aMesh), MinY(aMesh), MaxY(aMesh), aMesh%boundaryType)
call LogMessage(log, DEBUG_LOGGING_LEVEL, logkey, ' base mesh returned.')
newParticles => newMesh%particles
newEdges => newMesh%edges
newPanels => newMesh%panels
!
! setup old mesh as interpolation source
!
call New(bivarData, amesh)
allocate(integerWorkParticles( 31*bivarData%n + newParticles%N_Max))
allocate(integerWorkPanels( 31*bivarData%n + newPanels%N_Max))
integerWorkParticles = 0
integerWorkPanels = 0
call LogMessage(log, DEBUG_LOGGING_LEVEL, logkey, ' alpha source data ready.')
!
! interpolate Lagrangian parameter from old mesh to new mesh
!
md = 1
call IDBVIP( md, bivarData%n, bivarData%x, bivarData%y, bivarData%x0, &
newParticles%N, newParticles%x(1,1:newParticles%N), newParticles%x(2,1:newParticles%N), &
newParticles%x0(1,1:newParticles%N), integerWorkParticles, bivarData%realWork)
md = 3
call IDBVIP( md, bivarData%n, bivarData%x, bivarData%y, bivarData%y0, &
newParticles%N, newParticles%x(1,1:newParticles%N), newParticles%x(2,1:newParticles%N), &
newParticles%x0(2,1:newParticles%N), integerWorkParticles, bivarData%realWork)
md = 1
call IDBVIP( md, bivarData%N, bivarData%x, bivarData%y, bivarData%x0, &
newPanels%N, newPanels%x(1,1:newPanels%n), newPanels%x(2,1:newPanels%n), &
newPanels%x0(1,1:newPanels%N), integerWorkPanels, bivarData%realWork)
md = 3
call IDBVIP( md, bivarData%N, bivarData%x, bivarData%y, bivarData%y0, &
newPanels%N, newPanels%x(1,1:newPanels%n), newPanels%x(2,1:newPanels%n), &
newPanels%x0(2,1:newPanels%N), integerWorkPanels, bivarData%realWork)
call LogMessage(log, DEBUG_LOGGING_LEVEL, logkey, ' base mesh interpolation complete.')
!
! set tracer values on new mesh
!
do k = 1, aMesh%nTracer
call AssignTracer(newParticles%tracer(1:newParticles%N,k), newParticles%x0(:,1:newParticles%N), reference, k)
call AssignTracer(newPanels%tracer(1:newPanels%N, k), newPanels%x0(:,1:newPanels%N), reference, k)
enddo
!
! set vorticity values on new mesh
!
call AssignVorticity(newParticles%relVort(1:newParticles%N), newParticles%x0(:,1:newParticles%N), reference)
call AssignVorticity(newPanels%relVort(1:newPanels%N), newPanels%x0(:,1:newPanels%N), reference)
call LogMessage(log, DEBUG_LOGGING_LEVEL, logkey, ' LagRemeshReference : new uniform mesh ready.')
if ( aMesh%AMR > 0 ) then
call StartSection(log,'LagrangianRemeshToReference AMR :')
allocate(refineFlag(newPanels%N_Max))
refineFlag = .FALSE.
if ( remesh%vorticityRefine ) then
call FlagPanelsForMaxCirculationRefinement(refineFlag, newMesh, remesh, startIndex, counters(1))
call FlagPanelsForVorticityVariationRefinement(refineFlag, newMesh, remesh, startIndex, counters(2))
endif
if ( remesh%flowMapRefine ) then
call FlagPanelsForFlowMapRefinement(refineFlag, newMesh, remesh, startIndex, counters(3))
endif
if ( remesh%tracerRefine ) then
call FlagPanelsForTracerMassRefinement(refineFlag, newMesh, remesh, startIndex, counters(4))
call FlagPanelsForTracerVariationRefinement(refineFlag, newMesh, remesh, startIndex, counters(5))
endif
refineCount = count(refineFlag)
spaceLeft = newPanels%N_Max - newPanels%N
if ( refineCount > 0 ) then
if ( spaceLeft / 4 > refineCount ) then
keepGoing = .TRUE.
do while (keepGoing)
amrLoopCounter = amrLoopCounter + 1
write(logString,'(A,I2,A,I8,A)') 'AMR loop ', amrLoopCounter,&
' : refining ', refineCount,' panels.'
call LogMessage(log, TRACE_LOGGING_LEVEL, 'LagRemeshReference : ', trim(logstring))
!
! divide flagged panels
!
nOldPanels = newPanels%N
nOldParticles = newParticles%N
do j = startIndex, nOldPanels
if ( refineFlag(j) ) then
call DividePanel(newMesh, j)
refineFlag(j) = .FALSE.
endif
enddo
!
! TO DO : ensure adjacent panels differ by no more than 1 level of refinement
!
!
! interpolate Lagrangian parameter to new panels
!
md = 1
call IDBVIP(md, bivarData%n, bivarData%x, bivarData%y, bivarData%x0, &
newParticles%N - nOldParticles, &
newParticles%x(1,nOldParticles+1:newParticles%N), newParticles%x(2,nOldParticles+1:newParticles%N), &
newParticles%x0(1, nOldParticles+1:newParticles%N), integerWorkParticles, bivarData%realwork)
md = 3
call IDBVIP(md, bivarData%n, bivarData%x, bivarData%y, bivarData%y0, &
newParticles%N - nOldParticles, &
newParticles%x(1,nOldParticles+1:newParticles%N), newParticles%x(2,nOldParticles+1:newParticles%N), &
newParticles%x0(2, nOldParticles+1:newParticles%N), integerWorkParticles, bivarData%realwork)
md = 1
call IDBVIP(md, bivarData%n, bivarData%x, bivarData%y, bivarData%x0, &
newPanels%N - nOldPanels, &
newPanels%x(1,nOldPanels+1:newPanels%N), newPanels%x(2,nOldPanels+1:newPanels%N), &
newPanels%x0(1,nOldPanels+1:newPanels%N), integerWorkPanels, bivarData%realWork)
md = 3
call IDBVIP(md, bivarData%n, bivarData%x, bivarData%y, bivarData%y0, &
newPanels%N - nOldPanels, &
newPanels%x(1,nOldPanels+1:newPanels%N), newPanels%x(2,nOldPanels+1:newPanels%N), &
newPanels%x0(2,nOldPanels+1:newPanels%N), integerWorkPanels, bivarData%realWork)
!
! set data on new particles, panels
!
do k = 1, aMesh%nTracer
call AssignTracer(newParticles%tracer(nOldParticles+1:newParticles%N,k), &
newParticles%x0(:,nOldParticles+1:newParticles%N), reference, k)
call AssignTracer(newPanels%tracer(nOldPanels+1:newPanels%N,k), &
newPanels%x0(:,nOldPanels+1:newPanels%N), reference, k)
enddo
call AssignVorticity(newParticles%relVort(nOldParticles+1:newParticles%N), &
newParticles%x0(:, nOldParticles+1:newParticles%N), reference)
call AssignVorticity(newPanels%relVort(nOldPanels+1:newPanels%N), &
newPanels%x0(:,nOldPanels+1:newPanels%N), reference)
!
! prevent too much refinement
!
if ( amrLoopCounter >= remesh%refinementLimit ) then
keepGoing = .FALSE.
call LogMessage(log, WARNING_LOGGING_LEVEL, 'LagRemeshReference WARNING : ', 'refinement limit reached.')
call LogMessage(log, TRACE_LOGGING_LEVEL, 'LagRemeshReference N_Active = ', newPanels%N_Active)
write(logstring,'(A, I8, A)') ' max circulation criterion triggered ', counters(1), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'LagRemeshReference : ',trim(logstring))
write(logstring,'(A, I8, A)') ' vort. variation criterion triggered ', counters(2), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'LagRemeshReference : ',trim(logstring))
write(logstring,'(A, I8, A)') 'flwmap variation criterion triggered ', counters(3), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'LagRemeshReference : ',trim(logstring))
write(logstring,'(A, I8, A)') ' tracermax criterion triggered ', counters(4), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'LagRemeshReference : ',trim(logstring))
write(logstring,'(A, I8, A)') 'tracer variation criterion triggered ', counters(5), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'LagRemeshReference : ',trim(logstring))
else
!
! apply refinement criteria to new panels
!
startIndex = nOldPanels+1
nOldPanels = newPanels%N
if ( remesh%vorticityRefine ) then
call FlagPanelsForMaxCirculationRefinement(refineFlag, newMesh, remesh, startIndex, counters(1) )
call FlagPanelsForVorticityVariationRefinement(refineFlag, newMesh, remesh, startIndex, counters(2) )
endif
if ( remesh%flowMapRefine ) then
call FlagPanelsForFlowMapRefinement(refineFlag, newMesh, remesh, startIndex, counters(3))
endif
if ( remesh%tracerRefine ) then
call FlagPanelsForTracerMassRefinement(refineFlag, newMesh, remesh, startIndex, counters(4))
call FlagPanelsForTracerVariationRefinement(refineFlag, newMesh, remesh, startIndex, counters(5))
endif
refineCount = count(refineFlag)
spaceLeft = newPanels%N_Max - newPanels%N
!
! check stopping criteria
!
if ( refineCount == 0 ) then
keepGoing = .FALSE.
call LogMessage(log, TRACE_LOGGING_LEVEL, 'LagRemeshReference : ', ' refinement converged.')
call LogMessage(log, TRACE_LOGGING_LEVEL, 'LagRemeshReference N_Active = ', newPanels%N_Active)
write(logstring,'(A, I8, A)') ' max circulation criterion triggered ', counters(1), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'LagRemeshReference : ',trim(logstring))
write(logstring,'(A, I8, A)') ' vort. variation criterion triggered ', counters(2), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'LagRemeshReference : ',trim(logstring))
write(logstring,'(A, I8, A)') 'flwmap variation criterion triggered ', counters(3), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'LagRemeshReference : ',trim(logstring))
write(logstring,'(A, I8, A)') ' tracermax criterion triggered ', counters(4), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'LagRemeshReference : ',trim(logstring))
write(logstring,'(A, I8, A)') 'tracer variation criterion triggered ', counters(5), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'LagRemeshReference : ',trim(logstring))
elseif ( spaceLeft / 4 < refineCount ) then
keepGoing = .FALSE.
call LogMessage(log, WARNING_LOGGING_LEVEL, 'LagRemeshReference WARNING: ', ' not enough memory to continue AMR.')
call LogMessage(log, TRACE_LOGGING_LEVEL, 'LagRemeshReference N_Active = ', newPanels%N_Active)
write(logstring,'(A, I8, A)') ' max circulation criterion triggered ', counters(1), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'LagRemeshReference : ',trim(logstring))
write(logstring,'(A, I8, A)') ' vort. variation criterion triggered ', counters(2), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'LagRemeshReference : ',trim(logstring))
write(logstring,'(A, I8, A)') 'flwmap variation criterion triggered ', counters(3), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'LagRemeshReference : ',trim(logstring))
write(logstring,'(A, I8, A)') ' tracermax criterion triggered ', counters(4), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'LagRemeshReference : ',trim(logstring))
write(logstring,'(A, I8, A)') 'tracer variation criterion triggered ', counters(5), ' times.'
call LogMessage(log, TRACE_LOGGING_LEVEL,'LagRemeshReference : ',trim(logstring))
endif ! stopping criteria met
endif! refinement limit reached
enddo!while keepgoing
else
call LogMessage(log, WARNING_LOGGING_LEVEL,'LagRemeshReference : ',' not enough memory for AMR.')
endif ! spaceleft
else ! refinecount == 0
call LogMessage(log, TRACE_LOGGING_LEVEL, 'LagRemeshReference : ',' no refinement necessary.')
endif
call EndSection(log)
deallocate(refineFlag)
endif!AMR
!
! replace old mesh with new mesh
!
call Copy(aMesh, newMesh)
!
! clean up
!
deallocate(integerWorkParticles)
deallocate(integerWorkPanels)
call Delete(bivarData)
call Delete(newMesh)
end subroutine
!
!----------------
! Module methods : module- or type-specific private functions
!----------------
!
subroutine FlagPanelsForMaxCirculationRefinement(refineFlag, aMesh, remeshData, startIndex, counter)
logical(klog), intent(inout) :: refineFlag(:)
type(PlaneMesh), intent(in) :: aMesh
type(RemeshSetup), intent(in) :: remeshData
integer(kint), intent(in) :: startIndex
integer(kint), intent(inout) :: counter
!
type(Particles), pointer :: aParticles
type(Panels), pointer :: aPanels
integer(kint) :: j
aPanels => aMesh%panels
aParticles => aMesh%particles
do j = startIndex, aPanels%N
if ( .NOT. aPanels%hasChildren(j) ) then
if ( abs(aPanels%relVort(j)) * aPanels%area(j) > remeshData%maxCircTol ) then
refineFlag(j) = .TRUE.
counter = counter + 1
endif
endif
enddo
end subroutine
subroutine FlagPanelsForVorticityVariationRefinement(refineFlag, aMesh, remeshData, startIndex, counter)
logical(klog), intent(inout) :: refineFlag(:)
type(PlaneMesh), intent(in) :: aMesh
type(RemeshSetup), intent(in) :: remeshData
integer(kint), intent(in) :: startIndex
integer(kint), intent(inout) :: counter
!
type(Particles), pointer :: aParticles
type(Panels), pointer :: aPanels
integer(kint) :: edgeList(8), vertList(8), nverts, j, k
real(kreal) :: maxVort, minVort, vortVar
aPanels => aMesh%panels
aParticles => aMesh%particles
do j = startIndex, aPanels%N
if ( .NOT. aPanels%hasChildren(j) ) then
maxVort = aPanels%relVort(j)
minVort = aPanels%relVort(j)
call CCWEdgesAndParticlesAroundPanel(edgeList, vertLIst, nVerts, aMesh, j)
do k = 1, nVerts
if ( aParticles%relVort(vertList(k)) > maxVort ) maxVort = aParticles%relVort(vertList(k))
if ( aParticles%relVort(vertList(k)) < minVort ) minVort = aParticles%relVort(vertList(k))
enddo
vortVar = maxVort - minVort
if ( vortVar > remeshData%vortVarTol ) then
refineFlag(j) = .TRUE.
counter = counter + 1
endif
endif
enddo
end subroutine
subroutine FlagPanelsForFlowMapRefinement(refineFlag, aMesh, remeshData, startIndex, counter)
logical(klog), intent(inout) :: refineFlag(:)
type(PlaneMesh), intent(in) :: aMesh
type(RemeshSetup), intent(in) :: remeshData
integer(kint), intent(in) :: startIndex
integer(kint), intent(inout) :: counter