-
Notifications
You must be signed in to change notification settings - Fork 62
/
analytic_analysis_routines.f90
1885 lines (1707 loc) · 109 KB
/
analytic_analysis_routines.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
!> \file
!> \author Ting Yu
!> \brief This module handles all analytic analysis routines.
!>
!> \section LICENSE
!>
!> Version: MPL 1.1/GPL 2.0/LGPL 2.1
!>
!> The contents of this file are subject to the Mozilla Public License
!> Version 1.1 (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.mozilla.org/MPL/
!>
!> Software distributed under the License is distributed on an "AS IS"
!> basis, WITHOUT WARRANTY OF ANY KIND, either express or implied. See the
!> License for the specific language governing rights and limitations
!> under the License.
!>
!> The Original Code is OpenCMISS
!>
!> The Initial Developer of the Original Code is University of Auckland,
!> Auckland, New Zealand, the University of Oxford, Oxford, United
!> Kingdom and King's College, London, United Kingdom. Portions created
!> by the University of Auckland, the University of Oxford and King's
!> College, London are Copyright (C) 2007-2010 by the University of
!> Auckland, the University of Oxford and King's College, London.
!> All Rights Reserved.
!>
!> Contributor(s):
!>
!> Alternatively, the contents of this file may be used under the terms of
!> either the GNU General Public License Version 2 or later (the "GPL"), or
!> the GNU Lesser General Public License Version 2.1 or later (the "LGPL"),
!> in which case the provisions of the GPL or the LGPL are applicable instead
!> of those above. If you wish to allow use of your version of this file only
!> under the terms of either the GPL or the LGPL, and not to allow others to
!> use your version of this file under the terms of the MPL, indicate your
!> decision by deleting the provisions above and replace them with the notice
!> and other provisions required by the GPL or the LGPL. If you do not delete
!> the provisions above, a recipient may use your version of this file under
!> the terms of any one of the MPL, the GPL or the LGPL.
!>
!>This module handles all analytic analysis routines.
MODULE ANALYTIC_ANALYSIS_ROUTINES
USE BasisRoutines
USE BasisAccessRoutines
USE CmissMPI
USE ComputationEnvironment
USE Constants
USE FIELD_ROUTINES
USE FieldAccessRoutines
USE INPUT_OUTPUT
USE ISO_VARYING_STRING
USE Kinds
USE MatrixVector
#ifndef NOMPIMOD
USE MPI
#endif
USE Strings
USE Timer
USE Types
#include "macros.h"
IMPLICIT NONE
PRIVATE
#ifdef NOMPIMOD
#include "mpif.h"
#endif
!Module parameters
!> \addtogroup ANALYTIC_ANALYSIS_ROUTINES_ErrorTypes ANALYTIC_ANALYSIS_ROUTINES::ErrorTypes
!> \brief errors definition type parameters
!> \see ANALYTIC_ANALYSIS_ROUTINES,OPENCMISS_ErrorTypes
!>@{
INTEGER(INTG), PARAMETER :: ABSOLUTE_ERROR_TYPE=1 !<The absolute type \see ANALYTIC_ANALYSIS_ROUTINES_ErrorTypes,ANALYTIC_ANALYSIS_ROUTINES
INTEGER(INTG), PARAMETER :: PERCENTAGE_ERROR_TYPE=2 !<The percentage type \see ANALYTIC_ANALYSIS_ROUTINES_ErrorTypes,ANALYTIC_ANALYSIS_ROUTINES
INTEGER(INTG), PARAMETER :: RELATIVE_ERROR_TYPE=3 !<The relative type \see ANALYTIC_ANALYSIS_ROUTINES_ErrorTypes,ANALYTIC_ANALYSIS_ROUTINES
!>@}
!Module types
!Module variables
!Interfaces
PUBLIC AnalyticAnalysis_Output
PUBLIC AnalyticAnalysis_AbsoluteErrorGetNode,AnalyticAnalysis_PercentageErrorGetNode, &
& AnalyticAnalysis_RelativeErrorGetNode,AnalyticAnalysis_RMSErrorGetNode
PUBLIC AnalyticAnalysis_AbsoluteErrorGetElement,AnalyticAnalysis_PercentageErrorGetElement, &
& AnalyticAnalysis_RelativeErrorGetElement,AnalyticAnalysis_RMSErrorGetElement
PUBLIC AnalyticAnalysis_AbsoluteErrorGetConstant,AnalyticAnalysis_PercentageErrorGetConstant, &
& AnalyticAnalysis_RelativeErrorGetConstant
PUBLIC AnalyticAnalysis_IntegralNumericalValueGet,AnalyticAnalysis_IntegralAnalyticValueGet, &
& AnalyticAnalysis_IntegralPercentageErrorGet,AnalyticAnalysis_IntegralAbsoluteErrorGet, &
& AnalyticAnalysis_IntegralRelativeErrorGet,AnalyticAnalysis_IntegralNIDNumericalValueGet, &
& AnalyticAnalysis_IntegralNIDErrorGet
CONTAINS
!
!================================================================================================================================
!
!>Output the analytic error analysis for a dependent field compared to the analytic values parameter set. \see OPENCMISS::CMISSAnalyticAnalytisOutput
SUBROUTINE AnalyticAnalysis_Output(FIELD,FILENAME,ERR,ERROR,*)
!Argument variables
TYPE(FIELD_TYPE), INTENT(IN), POINTER :: FIELD !<A pointer to the dependent field to calculate the analytic error analysis for
CHARACTER(LEN=*) :: FILENAME !<If not empty, the filename to output the analytic analysis to. If empty, the analysis will be output to the standard output
INTEGER(INTG), INTENT(OUT) :: ERR !<The error code
TYPE(VARYING_STRING), INTENT(OUT) :: ERROR !<The error string
!Local Variables
INTEGER(INTG) :: component_idx,deriv_idx,element_idx,GHOST_NUMBER(8),local_ny,MESH_COMPONENT,MPI_IERROR,node_idx, &
& NUMBER(8),OUTPUT_ID,var_idx,variable_type
REAL(DP) :: GHOST_RMS_ERROR_PER(8),GHOST_RMS_ERROR_ABS(8),GHOST_RMS_ERROR_REL(8),RMS_ERROR_PER(8),RMS_ERROR_ABS(8), &
& RMS_ERROR_REL(8),VALUES(5)
REAL(DP), POINTER :: ANALYTIC_VALUES(:),NUMERICAL_VALUES(:)
REAL(DP), ALLOCATABLE :: INTEGRAL_ERRORS(:,:),GHOST_INTEGRAL_ERRORS(:,:)
CHARACTER(LEN=40) :: FIRST_FORMAT
CHARACTER(LEN=MAXSTRLEN) :: FILE_NAME
TYPE(DECOMPOSITION_TYPE), POINTER :: DECOMPOSITION
TYPE(DECOMPOSITION_ELEMENTS_TYPE), POINTER :: ELEMENTS_DECOMPOSITION
TYPE(DECOMPOSITION_TOPOLOGY_TYPE), POINTER :: DECOMPOSITION_TOPOLOGY
TYPE(DOMAIN_TYPE), POINTER :: DOMAIN
TYPE(DOMAIN_ELEMENTS_TYPE), POINTER :: ELEMENTS_DOMAIN
TYPE(DOMAIN_NODES_TYPE), POINTER :: NODES_DOMAIN
TYPE(DOMAIN_TOPOLOGY_TYPE), POINTER :: DOMAIN_TOPOLOGY
TYPE(FIELD_VARIABLE_TYPE), POINTER :: FIELD_VARIABLE
TYPE(VARYING_STRING) :: LOCAL_ERROR,LOCAL_STRING
NULLIFY(ANALYTIC_VALUES)
NULLIFY(NUMERICAL_VALUES)
ENTERS("AnalyticAnalysis_Output",ERR,ERROR,*999)
IF(ASSOCIATED(FIELD)) THEN
IF(FIELD%FIELD_FINISHED) THEN
IF(FIELD%DEPENDENT_TYPE==FIELD_DEPENDENT_TYPE) THEN
IF(LEN_TRIM(FILENAME)>=1) THEN
!!TODO \todo have more general ascii file mechanism
IF(computationalEnvironment%numberOfComputationalNodes>1) THEN
WRITE(FILE_NAME,'(A,".opanal.",I0)') FILENAME(1:LEN_TRIM(FILENAME)),computationalEnvironment% &
& myComputationalNodeNumber
ELSE
FILE_NAME=FILENAME(1:LEN_TRIM(FILENAME))//".opanal"
ENDIF
OUTPUT_ID=IO1_FILE_UNIT
OPEN(UNIT=OUTPUT_ID,FILE=FILE_NAME(1:LEN_TRIM(FILE_NAME)),STATUS="REPLACE",FORM="FORMATTED",IOSTAT=ERR)
IF(ERR/=0) CALL FlagError("Error opening analysis output file.",ERR,ERROR,*999)
ELSE
OUTPUT_ID=GENERAL_OUTPUT_TYPE
ENDIF
DECOMPOSITION=>FIELD%DECOMPOSITION
IF(ASSOCIATED(DECOMPOSITION)) THEN
DECOMPOSITION_TOPOLOGY=>DECOMPOSITION%TOPOLOGY
IF(ASSOCIATED(DECOMPOSITION_TOPOLOGY)) THEN
CALL WRITE_STRING(OUTPUT_ID,"Analytic error analysis:",ERR,ERROR,*999)
CALL WRITE_STRING(OUTPUT_ID,"",ERR,ERROR,*999)
LOCAL_STRING="Field "//TRIM(NUMBER_TO_VSTRING(FIELD%USER_NUMBER,"*",ERR,ERROR))//" : "//FIELD%LABEL
IF(ERR/=0) GOTO 999
CALL WRITE_STRING(OUTPUT_ID,LOCAL_STRING,ERR,ERROR,*999)
NULLIFY(NUMERICAL_VALUES)
NULLIFY(ANALYTIC_VALUES)
!Loop over the variables
DO var_idx=1,FIELD%NUMBER_OF_VARIABLES
variable_type=FIELD%VARIABLES(var_idx)%VARIABLE_TYPE
FIELD_VARIABLE=>FIELD%VARIABLE_TYPE_MAP(variable_type)%PTR
IF(ASSOCIATED(FIELD_VARIABLE)) THEN
CALL WRITE_STRING(OUTPUT_ID,"",ERR,ERROR,*999)
LOCAL_STRING="Variable "//TRIM(NUMBER_TO_VSTRING(variable_type,"*",ERR,ERROR))//" : "// &
& FIELD_VARIABLE%VARIABLE_LABEL
IF(ERR/=0) GOTO 999
CALL WRITE_STRING(OUTPUT_ID,LOCAL_STRING,ERR,ERROR,*999)
CALL WRITE_STRING(OUTPUT_ID,"",ERR,ERROR,*999)
!Get the dependent and analytic parameter sets
CALL FIELD_PARAMETER_SET_DATA_GET(FIELD,variable_type,FIELD_VALUES_SET_TYPE,NUMERICAL_VALUES,ERR,ERROR,*999)
CALL FIELD_PARAMETER_SET_DATA_GET(FIELD,variable_type,FIELD_ANALYTIC_VALUES_SET_TYPE,ANALYTIC_VALUES, &
& ERR,ERROR,*999)
!Loop over the components
DO component_idx=1,FIELD%VARIABLES(var_idx)%NUMBER_OF_COMPONENTS
MESH_COMPONENT=FIELD_VARIABLE%COMPONENTS(component_idx)%MESH_COMPONENT_NUMBER
DOMAIN=>FIELD_VARIABLE%COMPONENTS(component_idx)%DOMAIN
IF(ASSOCIATED(DOMAIN)) THEN
DOMAIN_TOPOLOGY=>DOMAIN%TOPOLOGY
IF(ASSOCIATED(DOMAIN_TOPOLOGY)) THEN
LOCAL_STRING="Component "//TRIM(NUMBER_TO_VSTRING(component_idx,"*",ERR,ERROR))//" : "// &
& FIELD_VARIABLE%COMPONENTS(component_idx)%COMPONENT_LABEL
IF(ERR/=0) GOTO 999
CALL WRITE_STRING(OUTPUT_ID,LOCAL_STRING,ERR,ERROR,*999)
CALL WRITE_STRING(OUTPUT_ID,"",ERR,ERROR,*999)
SELECT CASE(FIELD%VARIABLES(var_idx)%COMPONENTS(component_idx)%INTERPOLATION_TYPE)
CASE(FIELD_CONSTANT_INTERPOLATION)
CALL WRITE_STRING(OUTPUT_ID,"Constant errors:",ERR,ERROR,*999)
LOCAL_STRING=" Numerical Analytic % error Absolute err Relative err"
CALL WRITE_STRING(OUTPUT_ID,LOCAL_STRING,ERR,ERROR,*999)
local_ny=FIELD_VARIABLE%COMPONENTS(component_idx)%PARAM_TO_DOF_MAP%CONSTANT_PARAM2DOF_MAP
VALUES(1)=NUMERICAL_VALUES(local_ny)
VALUES(2)=ANALYTIC_VALUES(local_ny)
VALUES(3)=ANALYTIC_ANALYSIS_PERCENTAGE_ERROR(VALUES(1),VALUES(2))
VALUES(4)=ANALYTIC_ANALYSIS_ABSOLUTE_ERROR(VALUES(1),VALUES(2))
VALUES(5)=ANALYTIC_ANALYSIS_RELATIVE_ERROR(VALUES(1),VALUES(2))
CALL WRITE_STRING_VECTOR(OUTPUT_ID,1,1,5,5,5,VALUES,"(20X,5(2X,E12.5))","(20X,5(2X,E12.5))", &
& ERR,ERROR,*999)
CASE(FIELD_ELEMENT_BASED_INTERPOLATION)
ELEMENTS_DOMAIN=>DOMAIN_TOPOLOGY%ELEMENTS
IF(ASSOCIATED(ELEMENTS_DOMAIN)) THEN
DECOMPOSITION=>DOMAIN%DECOMPOSITION
IF(ASSOCIATED(DECOMPOSITION)) THEN
DECOMPOSITION_TOPOLOGY=>DECOMPOSITION%TOPOLOGY
IF(ASSOCIATED(DECOMPOSITION_TOPOLOGY)) THEN
ELEMENTS_DECOMPOSITION=>DECOMPOSITION_TOPOLOGY%ELEMENTS
IF(ASSOCIATED(ELEMENTS_DECOMPOSITION)) THEN
NUMBER=0
RMS_ERROR_PER=0.0_DP
RMS_ERROR_ABS=0.0_DP
RMS_ERROR_REL=0.0_DP
GHOST_NUMBER=0
GHOST_RMS_ERROR_PER=0.0_DP
GHOST_RMS_ERROR_ABS=0.0_DP
GHOST_RMS_ERROR_REL=0.0_DP
CALL WRITE_STRING(OUTPUT_ID,"Element errors:",ERR,ERROR,*999)
LOCAL_STRING= &
& " Element# Numerical Analytic % error Absolute err Relative err"
CALL WRITE_STRING(OUTPUT_ID,LOCAL_STRING,ERR,ERROR,*999)
DO element_idx=1,ELEMENTS_DOMAIN%NUMBER_OF_ELEMENTS
local_ny=FIELD_VARIABLE%COMPONENTS(component_idx)%PARAM_TO_DOF_MAP% &
& ELEMENT_PARAM2DOF_MAP%ELEMENTS(element_idx)
VALUES(1)=NUMERICAL_VALUES(local_ny)
VALUES(2)=ANALYTIC_VALUES(local_ny)
VALUES(3)=ANALYTIC_ANALYSIS_PERCENTAGE_ERROR(VALUES(1),VALUES(2))
VALUES(4)=ANALYTIC_ANALYSIS_ABSOLUTE_ERROR(VALUES(1),VALUES(2))
VALUES(5)=ANALYTIC_ANALYSIS_RELATIVE_ERROR(VALUES(1),VALUES(2))
!Accumlate the RMS errors
NUMBER(1)=NUMBER(1)+1
RMS_ERROR_PER(1)=RMS_ERROR_PER(1)+VALUES(3)*VALUES(3)
RMS_ERROR_ABS(1)=RMS_ERROR_ABS(1)+VALUES(4)*VALUES(4)
RMS_ERROR_REL(1)=RMS_ERROR_REL(1)+VALUES(5)*VALUES(5)
WRITE(FIRST_FORMAT,"(A,I10,A)") "('",ELEMENTS_DECOMPOSITION%ELEMENTS(element_idx)%USER_NUMBER, &
& "',20X,3(2X,E12.5))"
CALL WRITE_STRING_VECTOR(OUTPUT_ID,1,1,5,5,5,VALUES,FIRST_FORMAT,"(20X,5(2X,E12.5))", &
& ERR,ERROR,*999)
ENDDO !element_idx
DO element_idx=ELEMENTS_DOMAIN%NUMBER_OF_ELEMENTS+1,ELEMENTS_DOMAIN%TOTAL_NUMBER_OF_ELEMENTS
local_ny=FIELD_VARIABLE%COMPONENTS(component_idx)%PARAM_TO_DOF_MAP% &
& ELEMENT_PARAM2DOF_MAP%ELEMENTS(element_idx)
VALUES(1)=NUMERICAL_VALUES(local_ny)
VALUES(2)=ANALYTIC_VALUES(local_ny)
VALUES(3)=ANALYTIC_ANALYSIS_PERCENTAGE_ERROR(VALUES(1),VALUES(2))
VALUES(4)=ANALYTIC_ANALYSIS_ABSOLUTE_ERROR(VALUES(1),VALUES(2))
VALUES(5)=ANALYTIC_ANALYSIS_RELATIVE_ERROR(VALUES(1),VALUES(2))
!Accumlate the RMS errors
GHOST_NUMBER(1)=GHOST_NUMBER(1)+1
GHOST_RMS_ERROR_PER(1)=GHOST_RMS_ERROR_PER(1)+VALUES(3)*VALUES(3)
GHOST_RMS_ERROR_ABS(1)=GHOST_RMS_ERROR_ABS(1)+VALUES(4)*VALUES(4)
GHOST_RMS_ERROR_REL(1)=GHOST_RMS_ERROR_REL(1)+VALUES(5)*VALUES(5)
WRITE(FIRST_FORMAT,"(A,I10,A)") "('",ELEMENTS_DECOMPOSITION%ELEMENTS(element_idx)%USER_NUMBER, &
& "',20X,3(2X,E12.5))"
CALL WRITE_STRING_VECTOR(OUTPUT_ID,1,1,5,5,5,VALUES,FIRST_FORMAT,"(20X,5(2X,E12.5))", &
& ERR,ERROR,*999)
ENDDO !node_idx
!Output RMS errors
CALL WRITE_STRING(OUTPUT_ID,"",ERR,ERROR,*999)
IF(NUMBER(1)>0) THEN
IF(computationalEnvironment%numberOfComputationalNodes>1) THEN
!Local elements only
CALL WRITE_STRING(OUTPUT_ID,"Local RMS errors:",ERR,ERROR,*999)
LOCAL_STRING= &
& " % error Absolute err Relative err"
CALL WRITE_STRING(OUTPUT_ID,LOCAL_STRING,ERR,ERROR,*999)
VALUES(1)=SQRT(RMS_ERROR_PER(deriv_idx)/NUMBER(deriv_idx))
VALUES(2)=SQRT(RMS_ERROR_ABS(deriv_idx)/NUMBER(deriv_idx))
VALUES(3)=SQRT(RMS_ERROR_REL(deriv_idx)/NUMBER(deriv_idx))
CALL WRITE_STRING_VECTOR(OUTPUT_ID,1,1,3,3,3,VALUES,"(46X,3(2X,E12.5))","(46X,3(2X,E12.5))", &
& ERR,ERROR,*999)
!Local and ghost nodes
CALL WRITE_STRING(OUTPUT_ID,"Local + Ghost RMS errors:",ERR,ERROR,*999)
LOCAL_STRING= &
& " % error Absolute err Relative err"
CALL WRITE_STRING(OUTPUT_ID,LOCAL_STRING,ERR,ERROR,*999)
VALUES(1)=SQRT((RMS_ERROR_PER(1)+GHOST_RMS_ERROR_PER(1))/(NUMBER(1)+GHOST_NUMBER(1)))
VALUES(2)=SQRT((RMS_ERROR_ABS(1)+GHOST_RMS_ERROR_ABS(1))/(NUMBER(1)+GHOST_NUMBER(1)))
VALUES(3)=SQRT((RMS_ERROR_REL(1)+GHOST_RMS_ERROR_REL(1))/(NUMBER(1)+GHOST_NUMBER(1)))
CALL WRITE_STRING_VECTOR(OUTPUT_ID,1,1,3,3,3,VALUES,"(46X,3(2X,E12.5))","(46X,3(2X,E12.5))", &
& ERR,ERROR,*999)
!Global RMS values
!Collect the values across the ranks
CALL MPI_ALLREDUCE(MPI_IN_PLACE,NUMBER,1,MPI_INTEGER,MPI_SUM, &
& computationalEnvironment%mpiCommunicator,MPI_IERROR)
CALL MPI_ERROR_CHECK("MPI_ALLREDUCE",MPI_IERROR,ERR,ERROR,*999)
CALL MPI_ALLREDUCE(MPI_IN_PLACE,RMS_ERROR_PER,1,MPI_DOUBLE_PRECISION,MPI_SUM, &
& computationalEnvironment%mpiCommunicator,MPI_IERROR)
CALL MPI_ERROR_CHECK("MPI_ALLREDUCE",MPI_IERROR,ERR,ERROR,*999)
CALL MPI_ALLREDUCE(MPI_IN_PLACE,RMS_ERROR_ABS,1,MPI_DOUBLE_PRECISION,MPI_SUM, &
& computationalEnvironment%mpiCommunicator,MPI_IERROR)
CALL MPI_ERROR_CHECK("MPI_ALLREDUCE",MPI_IERROR,ERR,ERROR,*999)
CALL MPI_ALLREDUCE(MPI_IN_PLACE,RMS_ERROR_REL,1,MPI_DOUBLE_PRECISION,MPI_SUM, &
& computationalEnvironment%mpiCommunicator,MPI_IERROR)
CALL MPI_ERROR_CHECK("MPI_ALLREDUCE",MPI_IERROR,ERR,ERROR,*999)
CALL WRITE_STRING(OUTPUT_ID,"Global RMS errors:",ERR,ERROR,*999)
LOCAL_STRING= &
& " % error Absolute err Relative err"
CALL WRITE_STRING(OUTPUT_ID,LOCAL_STRING,ERR,ERROR,*999)
VALUES(1)=SQRT(RMS_ERROR_PER(1)/NUMBER(1))
VALUES(2)=SQRT(RMS_ERROR_ABS(1)/NUMBER(1))
VALUES(3)=SQRT(RMS_ERROR_REL(1)/NUMBER(1))
CALL WRITE_STRING_VECTOR(OUTPUT_ID,1,1,3,3,3,VALUES,"(46X,3(2X,E12.5))","(46X,3(2X,E12.5))", &
& ERR,ERROR,*999)
ELSE
CALL WRITE_STRING(OUTPUT_ID,"RMS errors:",ERR,ERROR,*999)
LOCAL_STRING= &
& " % error Absolute err Relative err"
CALL WRITE_STRING(OUTPUT_ID,LOCAL_STRING,ERR,ERROR,*999)
VALUES(1)=SQRT(RMS_ERROR_PER(deriv_idx)/NUMBER(deriv_idx))
VALUES(2)=SQRT(RMS_ERROR_ABS(deriv_idx)/NUMBER(deriv_idx))
VALUES(3)=SQRT(RMS_ERROR_REL(deriv_idx)/NUMBER(deriv_idx))
CALL WRITE_STRING_VECTOR(OUTPUT_ID,1,1,3,3,3,VALUES,"(46X,3(2X,E12.5))","(46X,3(2X,E12.5))", &
& ERR,ERROR,*999)
ENDIF
ENDIF
ELSE
CALL FlagError("Decomposition topology elements is not associated.",ERR,ERROR,*999)
ENDIF
ELSE
CALL FlagError("Decomposition topology is not associated.",ERR,ERROR,*999)
ENDIF
ELSE
CALL FlagError("Domain decomposition is not associated.",ERR,ERROR,*999)
ENDIF
ELSE
CALL FlagError("Elements domain topology is not associated.",ERR,ERROR,*999)
ENDIF
CASE(FIELD_NODE_BASED_INTERPOLATION)
NODES_DOMAIN=>DOMAIN_TOPOLOGY%NODES
IF(ASSOCIATED(NODES_DOMAIN)) THEN
NUMBER=0
RMS_ERROR_PER=0.0_DP
RMS_ERROR_ABS=0.0_DP
RMS_ERROR_REL=0.0_DP
GHOST_NUMBER=0
GHOST_RMS_ERROR_PER=0.0_DP
GHOST_RMS_ERROR_ABS=0.0_DP
GHOST_RMS_ERROR_REL=0.0_DP
CALL WRITE_STRING(OUTPUT_ID,"Nodal errors:",ERR,ERROR,*999)
LOCAL_STRING=" Node# Deriv# Numerical Analytic % error Absolute err Relative err"
CALL WRITE_STRING(OUTPUT_ID,LOCAL_STRING,ERR,ERROR,*999)
DO node_idx=1,NODES_DOMAIN%NUMBER_OF_NODES
DO deriv_idx=1,NODES_DOMAIN%NODES(node_idx)%NUMBER_OF_DERIVATIVES
local_ny=FIELD_VARIABLE%COMPONENTS(component_idx)%PARAM_TO_DOF_MAP% &
& NODE_PARAM2DOF_MAP%NODES(node_idx)%DERIVATIVES(deriv_idx)%VERSIONS(1)
VALUES(1)=NUMERICAL_VALUES(local_ny)
VALUES(2)=ANALYTIC_VALUES(local_ny)
VALUES(3)=ANALYTIC_ANALYSIS_PERCENTAGE_ERROR(VALUES(1),VALUES(2))
VALUES(4)=ANALYTIC_ANALYSIS_ABSOLUTE_ERROR(VALUES(1),VALUES(2))
VALUES(5)=ANALYTIC_ANALYSIS_RELATIVE_ERROR(VALUES(1),VALUES(2))
!Accumlate the RMS errors
NUMBER(deriv_idx)=NUMBER(deriv_idx)+1
RMS_ERROR_PER(deriv_idx)=RMS_ERROR_PER(deriv_idx)+VALUES(3)*VALUES(3)
RMS_ERROR_ABS(deriv_idx)=RMS_ERROR_ABS(deriv_idx)+VALUES(4)*VALUES(4)
RMS_ERROR_REL(deriv_idx)=RMS_ERROR_REL(deriv_idx)+VALUES(5)*VALUES(5)
IF(deriv_idx==1) THEN
WRITE(FIRST_FORMAT,"(A,I10,A,I6,A)") "('",NODES_DOMAIN%NODES(node_idx)%USER_NUMBER,"',2X,'", &
& deriv_idx,"',5(2X,E12.5))"
ELSE
WRITE(FIRST_FORMAT,"(A,I6,A)") "(12X,'",deriv_idx,"',5(2X,E12.5))"
ENDIF
CALL WRITE_STRING_VECTOR(OUTPUT_ID,1,1,5,5,5,VALUES,FIRST_FORMAT,"(20X,5(2X,E12.5))",ERR,ERROR,*999)
ENDDO !deriv_idx
ENDDO !node_idx
DO node_idx=NODES_DOMAIN%NUMBER_OF_NODES+1,NODES_DOMAIN%TOTAL_NUMBER_OF_NODES
DO deriv_idx=1,NODES_DOMAIN%NODES(node_idx)%NUMBER_OF_DERIVATIVES
local_ny=FIELD_VARIABLE%COMPONENTS(component_idx)%PARAM_TO_DOF_MAP% &
& NODE_PARAM2DOF_MAP%NODES(node_idx)%DERIVATIVES(deriv_idx)%VERSIONS(1)
VALUES(1)=NUMERICAL_VALUES(local_ny)
VALUES(2)=ANALYTIC_VALUES(local_ny)
VALUES(3)=ANALYTIC_ANALYSIS_PERCENTAGE_ERROR(VALUES(1),VALUES(2))
VALUES(4)=ANALYTIC_ANALYSIS_ABSOLUTE_ERROR(VALUES(1),VALUES(2))
VALUES(5)=ANALYTIC_ANALYSIS_RELATIVE_ERROR(VALUES(1),VALUES(2))
!Accumlate the RMS errors
GHOST_NUMBER(deriv_idx)=GHOST_NUMBER(deriv_idx)+1
GHOST_RMS_ERROR_PER(deriv_idx)=GHOST_RMS_ERROR_PER(deriv_idx)+VALUES(3)*VALUES(3)
GHOST_RMS_ERROR_ABS(deriv_idx)=GHOST_RMS_ERROR_ABS(deriv_idx)+VALUES(4)*VALUES(4)
GHOST_RMS_ERROR_REL(deriv_idx)=GHOST_RMS_ERROR_REL(deriv_idx)+VALUES(5)*VALUES(5)
IF(deriv_idx==1) THEN
WRITE(FIRST_FORMAT,"(A,I10,A,I6,A)") "('",NODES_DOMAIN%NODES(node_idx)%USER_NUMBER, &
& "',2X,'",deriv_idx,"',5(2X,E12.5))"
ELSE
WRITE(FIRST_FORMAT,"(A,I6,A)") "(12X,'",deriv_idx,"',5(2X,E12.5))"
ENDIF
CALL WRITE_STRING_VECTOR(OUTPUT_ID,1,1,5,5,5,VALUES,FIRST_FORMAT,"(20X,5(2X,E12.5))",ERR,ERROR,*999)
ENDDO !deriv_idx
ENDDO !node_idx
!Output RMS errors
CALL WRITE_STRING(OUTPUT_ID,"",ERR,ERROR,*999)
IF(computationalEnvironment%numberOfComputationalNodes>1) THEN
IF(ANY(NUMBER>0)) THEN
!Local nodes only
CALL WRITE_STRING(OUTPUT_ID,"Local RMS errors:",ERR,ERROR,*999)
LOCAL_STRING= &
& " Deriv# % error Absolute err Relative err"
CALL WRITE_STRING(OUTPUT_ID,LOCAL_STRING,ERR,ERROR,*999)
DO deriv_idx=1,8
IF(NUMBER(deriv_idx)>0) THEN
VALUES(1)=SQRT(RMS_ERROR_PER(deriv_idx)/NUMBER(deriv_idx))
VALUES(2)=SQRT(RMS_ERROR_ABS(deriv_idx)/NUMBER(deriv_idx))
VALUES(3)=SQRT(RMS_ERROR_REL(deriv_idx)/NUMBER(deriv_idx))
WRITE(FIRST_FORMAT,"(A,I6,A)") "(12X,'",deriv_idx,"',28X,3(2X,E12.5))"
CALL WRITE_STRING_VECTOR(OUTPUT_ID,1,1,3,3,3,VALUES,FIRST_FORMAT,"(46X,3(2X,E12.5))", &
& ERR,ERROR,*999)
ENDIF
ENDDO !deriv_idx
!Local and ghost nodes
CALL WRITE_STRING(OUTPUT_ID,"Local + Ghost RMS errors:",ERR,ERROR,*999)
LOCAL_STRING= &
& " Deriv# % error Absolute err Relative err"
CALL WRITE_STRING(OUTPUT_ID,LOCAL_STRING,ERR,ERROR,*999)
DO deriv_idx=1,8
IF(NUMBER(deriv_idx)>0) THEN
VALUES(1)=SQRT((RMS_ERROR_PER(deriv_idx)+GHOST_RMS_ERROR_PER(deriv_idx))/ &
& (NUMBER(deriv_idx)+GHOST_NUMBER(deriv_idx)))
VALUES(2)=SQRT((RMS_ERROR_ABS(deriv_idx)+GHOST_RMS_ERROR_ABS(deriv_idx))/ &
& (NUMBER(deriv_idx)+GHOST_NUMBER(deriv_idx)))
VALUES(3)=SQRT((RMS_ERROR_REL(deriv_idx)+GHOST_RMS_ERROR_REL(deriv_idx))/ &
& (NUMBER(deriv_idx)+GHOST_NUMBER(deriv_idx)))
WRITE(FIRST_FORMAT,"(A,I6,A)") "(12X,'",deriv_idx,"',28X,3(2X,E12.5))"
CALL WRITE_STRING_VECTOR(OUTPUT_ID,1,1,3,3,3,VALUES,FIRST_FORMAT,"(46X,3(2X,E12.5))", &
& ERR,ERROR,*999)
ENDIF
ENDDO !deriv_idx
!Global RMS values
!Collect the values across the ranks
CALL MPI_ALLREDUCE(MPI_IN_PLACE,NUMBER,8,MPI_INTEGER,MPI_SUM, &
& computationalEnvironment%mpiCommunicator,MPI_IERROR)
CALL MPI_ERROR_CHECK("MPI_ALLREDUCE",MPI_IERROR,ERR,ERROR,*999)
CALL MPI_ALLREDUCE(MPI_IN_PLACE,RMS_ERROR_PER,8,MPI_DOUBLE_PRECISION,MPI_SUM, &
& computationalEnvironment%mpiCommunicator,MPI_IERROR)
CALL MPI_ERROR_CHECK("MPI_ALLREDUCE",MPI_IERROR,ERR,ERROR,*999)
CALL MPI_ALLREDUCE(MPI_IN_PLACE,RMS_ERROR_ABS,8,MPI_DOUBLE_PRECISION,MPI_SUM, &
& computationalEnvironment%mpiCommunicator,MPI_IERROR)
CALL MPI_ERROR_CHECK("MPI_ALLREDUCE",MPI_IERROR,ERR,ERROR,*999)
CALL MPI_ALLREDUCE(MPI_IN_PLACE,RMS_ERROR_REL,8,MPI_DOUBLE_PRECISION,MPI_SUM, &
& computationalEnvironment%mpiCommunicator,MPI_IERROR)
CALL MPI_ERROR_CHECK("MPI_ALLREDUCE",MPI_IERROR,ERR,ERROR,*999)
CALL WRITE_STRING(OUTPUT_ID,"Global RMS errors:",ERR,ERROR,*999)
LOCAL_STRING= &
& " Deriv# % error Absolute err Relative err"
CALL WRITE_STRING(OUTPUT_ID,LOCAL_STRING,ERR,ERROR,*999)
DO deriv_idx=1,8
IF(NUMBER(deriv_idx)>0) THEN
VALUES(1)=SQRT(RMS_ERROR_PER(deriv_idx)/NUMBER(deriv_idx))
VALUES(2)=SQRT(RMS_ERROR_ABS(deriv_idx)/NUMBER(deriv_idx))
VALUES(3)=SQRT(RMS_ERROR_REL(deriv_idx)/NUMBER(deriv_idx))
WRITE(FIRST_FORMAT,"(A,I6,A)") "(12X,'",deriv_idx,"',28X,3(2X,E12.5))"
CALL WRITE_STRING_VECTOR(OUTPUT_ID,1,1,3,3,3,VALUES,FIRST_FORMAT,"(46X,3(2X,E12.5))", &
& ERR,ERROR,*999)
ENDIF
ENDDO !deriv_idx
ENDIF
ELSE
IF(ANY(NUMBER>0)) THEN
CALL WRITE_STRING(OUTPUT_ID,"RMS errors:",ERR,ERROR,*999)
LOCAL_STRING= &
& " Deriv# % error Absolute err Relative err"
CALL WRITE_STRING(OUTPUT_ID,LOCAL_STRING,ERR,ERROR,*999)
DO deriv_idx=1,8
IF(NUMBER(deriv_idx)>0) THEN
VALUES(1)=SQRT(RMS_ERROR_PER(deriv_idx)/NUMBER(deriv_idx))
VALUES(2)=SQRT(RMS_ERROR_ABS(deriv_idx)/NUMBER(deriv_idx))
VALUES(3)=SQRT(RMS_ERROR_REL(deriv_idx)/NUMBER(deriv_idx))
WRITE(FIRST_FORMAT,"(A,I6,A)") "(12X,'",deriv_idx,"',28X,3(2X,E12.5))"
CALL WRITE_STRING_VECTOR(OUTPUT_ID,1,1,3,3,3,VALUES,FIRST_FORMAT,"(46X,3(2X,E12.5))", &
& ERR,ERROR,*999)
ENDIF
ENDDO !deriv_idx
ENDIF
ENDIF
ELSE
CALL FlagError("Nodes domain topology is not associated.",ERR,ERROR,*999)
ENDIF
CASE(FIELD_GRID_POINT_BASED_INTERPOLATION)
CALL FlagError("Not implemented.",ERR,ERROR,*999)
CASE(FIELD_GAUSS_POINT_BASED_INTERPOLATION)
CALL FlagError("Not implemented.",ERR,ERROR,*999)
CASE DEFAULT
LOCAL_ERROR="The interpolation type of "// &
& TRIM(NUMBER_TO_VSTRING(FIELD%VARIABLES(var_idx)%COMPONENTS(component_idx)% &
& INTERPOLATION_TYPE,"*",ERR,ERROR))//" for component number "// &
& TRIM(NUMBER_TO_VSTRING(component_idx,"*",ERR,ERROR))//" of variable type "// &
& TRIM(NUMBER_TO_VSTRING(variable_type,"*",ERR,ERROR))//" of field number "// &
& TRIM(NUMBER_TO_VSTRING(FIELD%USER_NUMBER,"*",ERR,ERROR))//" is invalid."
CALL FlagError(LOCAL_ERROR,ERR,ERROR,*999)
END SELECT
ELSE
CALL FlagError("Domain topology is not associated.",ERR,ERROR,*999)
ENDIF
ELSE
CALL FlagError("Domain is not associated.",ERR,ERROR,*999)
ENDIF
CALL WRITE_STRING(OUTPUT_ID,"",ERR,ERROR,*999)
ENDDO !component_idx
!Restore the dependent and analytic parameter sets
CALL FIELD_PARAMETER_SET_DATA_RESTORE(FIELD,variable_type,FIELD_VALUES_SET_TYPE,NUMERICAL_VALUES, &
& ERR,ERROR,*999)
CALL FIELD_PARAMETER_SET_DATA_RESTORE(FIELD,variable_type,FIELD_ANALYTIC_VALUES_SET_TYPE,ANALYTIC_VALUES, &
& ERR,ERROR,*999)
!Allocated the integral errors
ALLOCATE(INTEGRAL_ERRORS(6,FIELD_VARIABLE%NUMBER_OF_COMPONENTS),STAT=ERR)
IF(ERR/=0) CALL FlagError("Could not allocate integral errors.",ERR,ERROR,*999)
ALLOCATE(GHOST_INTEGRAL_ERRORS(6,FIELD_VARIABLE%NUMBER_OF_COMPONENTS),STAT=ERR)
IF(ERR/=0) CALL FlagError("Could not allocate ghost integral errors.",ERR,ERROR,*999)
CALL ANALYTIC_ANALYSIS_INTEGRAL_ERRORS(FIELD_VARIABLE,INTEGRAL_ERRORS,GHOST_INTEGRAL_ERRORS,ERR,ERROR,*999)
IF(computationalEnvironment%numberOfComputationalNodes>1) THEN
CALL WRITE_STRING(OUTPUT_ID,"Local Integral errors:",ERR,ERROR,*999)
LOCAL_STRING="Component# Numerical Analytic % error Absolute err Relative err"
CALL WRITE_STRING(OUTPUT_ID,LOCAL_STRING,ERR,ERROR,*999)
DO component_idx=1,FIELD_VARIABLE%NUMBER_OF_COMPONENTS
VALUES(1)=INTEGRAL_ERRORS(1,component_idx)
VALUES(2)=INTEGRAL_ERRORS(3,component_idx)
VALUES(3)=ANALYTIC_ANALYSIS_PERCENTAGE_ERROR(VALUES(1),VALUES(2))
VALUES(4)=ANALYTIC_ANALYSIS_ABSOLUTE_ERROR(VALUES(1),VALUES(2))
VALUES(5)=ANALYTIC_ANALYSIS_RELATIVE_ERROR(VALUES(1),VALUES(2))
WRITE(FIRST_FORMAT,"(A,I10,A)") "('",component_idx,"',2X,'Intg ',5(2X,E12.5))"
CALL WRITE_STRING_VECTOR(OUTPUT_ID,1,1,5,5,5,VALUES,FIRST_FORMAT,"(20X,5(2X,E12.5))",ERR,ERROR,*999)
VALUES(1)=INTEGRAL_ERRORS(2,component_idx)
VALUES(2)=INTEGRAL_ERRORS(4,component_idx)
VALUES(3)=ANALYTIC_ANALYSIS_PERCENTAGE_ERROR(VALUES(1),VALUES(2))
VALUES(4)=ANALYTIC_ANALYSIS_ABSOLUTE_ERROR(VALUES(1),VALUES(2))
VALUES(5)=ANALYTIC_ANALYSIS_RELATIVE_ERROR(VALUES(1),VALUES(2))
WRITE(FIRST_FORMAT,"(A)") "(12X,'Intg^2',5(2X,E12.5))"
CALL WRITE_STRING_VECTOR(OUTPUT_ID,1,1,5,5,5,VALUES,FIRST_FORMAT,"(20X,5(2X,E12.5))",ERR,ERROR,*999)
ENDDO !component_idx
LOCAL_STRING="Component# Numerical Analytic NID NID(%)"
CALL WRITE_STRING(OUTPUT_ID,LOCAL_STRING,ERR,ERROR,*999)
DO component_idx=1,FIELD_VARIABLE%NUMBER_OF_COMPONENTS
VALUES(1)=INTEGRAL_ERRORS(5,component_idx)
VALUES(2)=INTEGRAL_ERRORS(3,component_idx)
VALUES(3)=ANALYTIC_ANALYSIS_NID_ERROR(VALUES(1),VALUES(2))
VALUES(4)=VALUES(3)*100.0_DP
WRITE(FIRST_FORMAT,"(A,I10,A)") "('",component_idx,"',2X,'Diff ',4(2X,E12.5))"
CALL WRITE_STRING_VECTOR(OUTPUT_ID,1,1,4,4,4,VALUES,FIRST_FORMAT,"(20X,4(2X,E12.5))",ERR,ERROR,*999)
VALUES(1)=INTEGRAL_ERRORS(6,component_idx)
VALUES(2)=INTEGRAL_ERRORS(4,component_idx)
VALUES(3)=ANALYTIC_ANALYSIS_NID_ERROR(VALUES(1),VALUES(2))
VALUES(4)=VALUES(3)*100.0_DP
WRITE(FIRST_FORMAT,"(A)") "(12X,'Diff^2',4(2X,E12.5))"
CALL WRITE_STRING_VECTOR(OUTPUT_ID,1,1,4,4,4,VALUES,FIRST_FORMAT,"(20X,4(2X,E12.5))",ERR,ERROR,*999)
ENDDO !component_idx
CALL WRITE_STRING(OUTPUT_ID,"Local + Ghost Integral errors:",ERR,ERROR,*999)
LOCAL_STRING="Component# Numerical Analytic % error Absolute err Relative err"
CALL WRITE_STRING(OUTPUT_ID,LOCAL_STRING,ERR,ERROR,*999)
DO component_idx=1,FIELD_VARIABLE%NUMBER_OF_COMPONENTS
VALUES(1)=INTEGRAL_ERRORS(1,component_idx)+GHOST_INTEGRAL_ERRORS(1,component_idx)
VALUES(2)=INTEGRAL_ERRORS(3,component_idx)+GHOST_INTEGRAL_ERRORS(3,component_idx)
VALUES(3)=ANALYTIC_ANALYSIS_PERCENTAGE_ERROR(VALUES(1),VALUES(2))
VALUES(4)=ANALYTIC_ANALYSIS_ABSOLUTE_ERROR(VALUES(1),VALUES(2))
VALUES(5)=ANALYTIC_ANALYSIS_RELATIVE_ERROR(VALUES(1),VALUES(2))
WRITE(FIRST_FORMAT,"(A,I10,A)") "('",component_idx,"',2X,'Intg ',5(2X,E12.5))"
CALL WRITE_STRING_VECTOR(OUTPUT_ID,1,1,5,5,5,VALUES,FIRST_FORMAT,"(20X,5(2X,E12.5))",ERR,ERROR,*999)
VALUES(1)=INTEGRAL_ERRORS(2,component_idx)
VALUES(2)=INTEGRAL_ERRORS(4,component_idx)
VALUES(3)=ANALYTIC_ANALYSIS_PERCENTAGE_ERROR(VALUES(1),VALUES(2))
VALUES(4)=ANALYTIC_ANALYSIS_ABSOLUTE_ERROR(VALUES(1),VALUES(2))
VALUES(5)=ANALYTIC_ANALYSIS_RELATIVE_ERROR(VALUES(1),VALUES(2))
WRITE(FIRST_FORMAT,"(A)") "(12X,'Intg^2',5(2X,E12.5))"
CALL WRITE_STRING_VECTOR(OUTPUT_ID,1,1,5,5,5,VALUES,FIRST_FORMAT,"(20X,5(2X,E12.5))",ERR,ERROR,*999)
ENDDO !component_idx
LOCAL_STRING="Component# Numerical Analytic NID NID(%)"
CALL WRITE_STRING(OUTPUT_ID,LOCAL_STRING,ERR,ERROR,*999)
DO component_idx=1,FIELD_VARIABLE%NUMBER_OF_COMPONENTS
VALUES(1)=INTEGRAL_ERRORS(5,component_idx)+GHOST_INTEGRAL_ERRORS(5,component_idx)
VALUES(2)=INTEGRAL_ERRORS(3,component_idx)+GHOST_INTEGRAL_ERRORS(3,component_idx)
VALUES(3)=ANALYTIC_ANALYSIS_NID_ERROR(VALUES(1),VALUES(2))
VALUES(4)=VALUES(3)*100.0_DP
WRITE(FIRST_FORMAT,"(A,I10,A)") "('",component_idx,"',2X,'Diff ',4(2X,E12.5))"
CALL WRITE_STRING_VECTOR(OUTPUT_ID,1,1,4,4,4,VALUES,FIRST_FORMAT,"(20X,4(2X,E12.5))",ERR,ERROR,*999)
VALUES(1)=INTEGRAL_ERRORS(6,component_idx)+GHOST_INTEGRAL_ERRORS(6,component_idx)
VALUES(2)=INTEGRAL_ERRORS(4,component_idx)+GHOST_INTEGRAL_ERRORS(4,component_idx)
VALUES(3)=ANALYTIC_ANALYSIS_NID_ERROR(VALUES(1),VALUES(2))
VALUES(4)=VALUES(3)*100.0_DP
WRITE(FIRST_FORMAT,"(A)") "(12X,'Diff^2',4(2X,E12.5))"
CALL WRITE_STRING_VECTOR(OUTPUT_ID,1,1,4,4,4,VALUES,FIRST_FORMAT,"(20X,4(2X,E12.5))",ERR,ERROR,*999)
ENDDO !component_idx
!Collect the values across the ranks
CALL MPI_ALLREDUCE(MPI_IN_PLACE,INTEGRAL_ERRORS,6*FIELD_VARIABLE%NUMBER_OF_COMPONENTS,MPI_DOUBLE_PRECISION, &
& MPI_SUM,computationalEnvironment%mpiCommunicator,MPI_IERROR)
CALL WRITE_STRING(OUTPUT_ID,"Global Integral errors:",ERR,ERROR,*999)
LOCAL_STRING="Component# Numerical Analytic % error Absolute err Relative err"
CALL WRITE_STRING(OUTPUT_ID,LOCAL_STRING,ERR,ERROR,*999)
DO component_idx=1,FIELD_VARIABLE%NUMBER_OF_COMPONENTS
CALL MPI_ERROR_CHECK("MPI_ALLREDUCE",MPI_IERROR,ERR,ERROR,*999)
VALUES(1)=INTEGRAL_ERRORS(1,component_idx)
VALUES(2)=INTEGRAL_ERRORS(3,component_idx)
VALUES(4)=ANALYTIC_ANALYSIS_ABSOLUTE_ERROR(VALUES(1),VALUES(2))
VALUES(5)=ANALYTIC_ANALYSIS_RELATIVE_ERROR(VALUES(1),VALUES(2))
WRITE(FIRST_FORMAT,"(A,I10,A)") "('",component_idx,"',2X,'Intg ',5(2X,E12.5))"
CALL WRITE_STRING_VECTOR(OUTPUT_ID,1,1,5,5,5,VALUES,FIRST_FORMAT,"(20X,5(2X,E12.5))",ERR,ERROR,*999)
VALUES(1)=INTEGRAL_ERRORS(2,component_idx)
VALUES(2)=INTEGRAL_ERRORS(4,component_idx)
VALUES(3)=ANALYTIC_ANALYSIS_PERCENTAGE_ERROR(VALUES(1),VALUES(2))
VALUES(4)=ANALYTIC_ANALYSIS_ABSOLUTE_ERROR(VALUES(1),VALUES(2))
VALUES(5)=ANALYTIC_ANALYSIS_RELATIVE_ERROR(VALUES(1),VALUES(2))
WRITE(FIRST_FORMAT,"(A)") "(12X,'Intg^2',5(2X,E12.5))"
CALL WRITE_STRING_VECTOR(OUTPUT_ID,1,1,5,5,5,VALUES,FIRST_FORMAT,"(20X,5(2X,E12.5))",ERR,ERROR,*999)
ENDDO !component_idx
LOCAL_STRING="Component# Numerical Analytic NID NID(%)"
CALL WRITE_STRING(OUTPUT_ID,LOCAL_STRING,ERR,ERROR,*999)
DO component_idx=1,FIELD_VARIABLE%NUMBER_OF_COMPONENTS
VALUES(1)=INTEGRAL_ERRORS(5,component_idx)
VALUES(2)=INTEGRAL_ERRORS(3,component_idx)
VALUES(3)=ANALYTIC_ANALYSIS_NID_ERROR(VALUES(1),VALUES(2))
VALUES(4)=VALUES(3)*100.0_DP
WRITE(FIRST_FORMAT,"(A,I10,A)") "('",component_idx,"',2X,'Diff ',4(2X,E12.5))"
CALL WRITE_STRING_VECTOR(OUTPUT_ID,1,1,4,4,4,VALUES,FIRST_FORMAT,"(20X,4(2X,E12.5))",ERR,ERROR,*999)
VALUES(1)=INTEGRAL_ERRORS(6,component_idx)
VALUES(2)=INTEGRAL_ERRORS(4,component_idx)
VALUES(3)=ANALYTIC_ANALYSIS_NID_ERROR(VALUES(1),VALUES(2))
VALUES(4)=VALUES(3)*100.0_DP
WRITE(FIRST_FORMAT,"(A)") "(12X,'Diff^2',4(2X,E12.5))"
CALL WRITE_STRING_VECTOR(OUTPUT_ID,1,1,4,4,4,VALUES,FIRST_FORMAT,"(20X,4(2X,E12.5))",ERR,ERROR,*999)
ENDDO !component_idx
ELSE
CALL WRITE_STRING(OUTPUT_ID,"Integral errors:",ERR,ERROR,*999)
LOCAL_STRING="Component# Numerical Analytic % error Absolute err Relative err"
CALL WRITE_STRING(OUTPUT_ID,LOCAL_STRING,ERR,ERROR,*999)
DO component_idx=1,FIELD_VARIABLE%NUMBER_OF_COMPONENTS
VALUES(1)=INTEGRAL_ERRORS(1,component_idx)
VALUES(2)=INTEGRAL_ERRORS(3,component_idx)
VALUES(3)=ANALYTIC_ANALYSIS_PERCENTAGE_ERROR(VALUES(1),VALUES(2))
VALUES(4)=ANALYTIC_ANALYSIS_ABSOLUTE_ERROR(VALUES(1),VALUES(2))
VALUES(5)=ANALYTIC_ANALYSIS_RELATIVE_ERROR(VALUES(1),VALUES(2))
WRITE(FIRST_FORMAT,"(A,I10,A)") "('",component_idx,"',2X,'Intg ',5(2X,E12.5))"
CALL WRITE_STRING_VECTOR(OUTPUT_ID,1,1,5,5,5,VALUES,FIRST_FORMAT,"(20X,5(2X,E12.5))",ERR,ERROR,*999)
VALUES(1)=INTEGRAL_ERRORS(2,component_idx)
VALUES(2)=INTEGRAL_ERRORS(4,component_idx)
VALUES(3)=ANALYTIC_ANALYSIS_PERCENTAGE_ERROR(VALUES(1),VALUES(2))
VALUES(4)=ANALYTIC_ANALYSIS_ABSOLUTE_ERROR(VALUES(1),VALUES(2))
VALUES(5)=ANALYTIC_ANALYSIS_RELATIVE_ERROR(VALUES(1),VALUES(2))
WRITE(FIRST_FORMAT,"(A)") "(12X,'Intg^2',5(2X,E12.5))"
CALL WRITE_STRING_VECTOR(OUTPUT_ID,1,1,5,5,5,VALUES,FIRST_FORMAT,"(20X,5(2X,E12.5))",ERR,ERROR,*999)
ENDDO !component_idx
LOCAL_STRING="Component# Numerical Analytic NID NID(%)"
CALL WRITE_STRING(OUTPUT_ID,LOCAL_STRING,ERR,ERROR,*999)
DO component_idx=1,FIELD_VARIABLE%NUMBER_OF_COMPONENTS
VALUES(1)=INTEGRAL_ERRORS(5,component_idx)
VALUES(2)=INTEGRAL_ERRORS(3,component_idx)
VALUES(3)=ANALYTIC_ANALYSIS_NID_ERROR(VALUES(1),VALUES(2))
VALUES(4)=VALUES(3)*100.0_DP
WRITE(FIRST_FORMAT,"(A,I10,A)") "('",component_idx,"',2X,'Diff ',4(2X,E12.5))"
CALL WRITE_STRING_VECTOR(OUTPUT_ID,1,1,4,4,4,VALUES,FIRST_FORMAT,"(20X,4(2X,E12.5))",ERR,ERROR,*999)
VALUES(1)=INTEGRAL_ERRORS(6,component_idx)
VALUES(2)=INTEGRAL_ERRORS(4,component_idx)
VALUES(3)=ANALYTIC_ANALYSIS_NID_ERROR(VALUES(1),VALUES(2))
VALUES(4)=VALUES(3)*100.0_DP
WRITE(FIRST_FORMAT,"(A)") "(12X,'Diff^2',4(2X,E12.5))"
CALL WRITE_STRING_VECTOR(OUTPUT_ID,1,1,4,4,4,VALUES,FIRST_FORMAT,"(20X,4(2X,E12.5))",ERR,ERROR,*999)
ENDDO !component_idx
CALL WRITE_STRING(OUTPUT_ID,"",ERR,ERROR,*999)
ENDIF
IF(ALLOCATED(INTEGRAL_ERRORS)) DEALLOCATE(INTEGRAL_ERRORS)
IF(ALLOCATED(GHOST_INTEGRAL_ERRORS)) DEALLOCATE(GHOST_INTEGRAL_ERRORS)
ELSE
CALL FlagError("Field variable is not associated.",ERR,ERROR,*999)
ENDIF
ENDDO !var_idx
ELSE
CALL FlagError("Decomposition topology is not associated.",ERR,ERROR,*999)
ENDIF
ELSE
CALL FlagError("Field decomposition is not associated.",ERR,ERROR,*999)
ENDIF
IF(LEN_TRIM(FILENAME)>=1) THEN
CLOSE(UNIT=OUTPUT_ID,IOSTAT=ERR)
IF(ERR/=0) CALL FlagError("Error closing analysis output file.",ERR,ERROR,*999)
ENDIF
ELSE
LOCAL_ERROR="Field number "//TRIM(NUMBER_TO_VSTRING(FIELD%USER_NUMBER,"*",ERR,ERROR))// &
& " is not a dependent field."
CALL FlagError(LOCAL_ERROR,ERR,ERROR,*999)
ENDIF
ELSE
LOCAL_ERROR="Field number "//TRIM(NUMBER_TO_VSTRING(FIELD%USER_NUMBER,"*",ERR,ERROR))// &
& " has not been finished."
CALL FlagError(LOCAL_ERROR,ERR,ERROR,*999)
ENDIF
ELSE
CALL FlagError("Field is not associated.",ERR,ERROR,*999)
ENDIF
EXITS("AnalyticAnalysis_Output")
RETURN
999 IF(ALLOCATED(INTEGRAL_ERRORS)) DEALLOCATE(INTEGRAL_ERRORS)
IF(ALLOCATED(GHOST_INTEGRAL_ERRORS)) DEALLOCATE(GHOST_INTEGRAL_ERRORS)
ERRORSEXITS("AnalyticAnalysis_Output",ERR,ERROR)
RETURN 1
END SUBROUTINE AnalyticAnalysis_Output
!
!================================================================================================================================
!
!>Calculates the absolute error between a numeric value and an analytic value.
PURE FUNCTION ANALYTIC_ANALYSIS_ABSOLUTE_ERROR(NUMERIC_VALUE,ANALYTIC_VALUE)
!Argument variables
REAL(DP), INTENT(IN) :: NUMERIC_VALUE !<The numerical value for the error calculation
REAL(DP), INTENT(IN) :: ANALYTIC_VALUE !<The analytic value for the error calculation
!Function result
REAL(DP) :: ANALYTIC_ANALYSIS_ABSOLUTE_ERROR
ANALYTIC_ANALYSIS_ABSOLUTE_ERROR=ABS(ANALYTIC_VALUE-NUMERIC_VALUE)
END FUNCTION ANALYTIC_ANALYSIS_ABSOLUTE_ERROR
!
!================================================================================================================================
!
!>Calculates the percentage error between a numeric value and an analytic value.
PURE FUNCTION ANALYTIC_ANALYSIS_PERCENTAGE_ERROR(NUMERIC_VALUE,ANALYTIC_VALUE)
!Argument variables
REAL(DP), INTENT(IN) :: NUMERIC_VALUE !<The numerical value for the error calculation
REAL(DP), INTENT(IN) :: ANALYTIC_VALUE !<The analytic value for the error calculation
!Function result
REAL(DP) :: ANALYTIC_ANALYSIS_PERCENTAGE_ERROR
IF(ABS(ANALYTIC_VALUE)>ZERO_TOLERANCE) THEN
ANALYTIC_ANALYSIS_PERCENTAGE_ERROR=(ANALYTIC_VALUE-NUMERIC_VALUE)/ANALYTIC_VALUE*100.0_DP
ELSE
ANALYTIC_ANALYSIS_PERCENTAGE_ERROR=0.0_DP
ENDIF
END FUNCTION ANALYTIC_ANALYSIS_PERCENTAGE_ERROR
!
!================================================================================================================================
!
!>Calculates the relative error between a numeric value and an analytic value.
PURE FUNCTION ANALYTIC_ANALYSIS_RELATIVE_ERROR(NUMERIC_VALUE,ANALYTIC_VALUE)
!Argument variables
REAL(DP), INTENT(IN) :: NUMERIC_VALUE !<The numerical value for the error calculation
REAL(DP), INTENT(IN) :: ANALYTIC_VALUE !<The analytic value for the error calculation
!Function result
REAL(DP) :: ANALYTIC_ANALYSIS_RELATIVE_ERROR
IF(ABS(1.0_DP+ANALYTIC_VALUE)>ZERO_TOLERANCE) THEN
ANALYTIC_ANALYSIS_RELATIVE_ERROR=(ANALYTIC_VALUE-NUMERIC_VALUE)/(1.0_DP+ANALYTIC_VALUE)
ELSE
ANALYTIC_ANALYSIS_RELATIVE_ERROR=0.0_DP
ENDIF
END FUNCTION ANALYTIC_ANALYSIS_RELATIVE_ERROR
!
!================================================================================================================================
!
!>Calculates the Normalised Integral Difference (NID) error with a numeric value and an analytic value.
PURE FUNCTION ANALYTIC_ANALYSIS_NID_ERROR(NUMERIC_VALUE,ANALYTIC_VALUE)
!Argument variables
REAL(DP), INTENT(IN) :: NUMERIC_VALUE !<The numerical value for the error calculation
REAL(DP), INTENT(IN) :: ANALYTIC_VALUE !<The analytic value for the error calculation
!Function result
REAL(DP) :: ANALYTIC_ANALYSIS_NID_ERROR
IF(ABS(ANALYTIC_VALUE)>ZERO_TOLERANCE) THEN
ANALYTIC_ANALYSIS_NID_ERROR=(ANALYTIC_VALUE-NUMERIC_VALUE)/ANALYTIC_VALUE
ELSE
ANALYTIC_ANALYSIS_NID_ERROR=(ANALYTIC_VALUE-NUMERIC_VALUE)/(1.0_DP+ANALYTIC_VALUE)
ENDIF
END FUNCTION ANALYTIC_ANALYSIS_NID_ERROR
!
!================================================================================================================================
!
!>Calculates the relative error between a numeric value and an analytic value.
SUBROUTINE ANALYTIC_ANALYSIS_INTEGRAL_ERRORS(FIELD_VARIABLE,INTEGRAL_ERRORS,GHOST_INTEGRAL_ERRORS,ERR,ERROR,*)
!Argument variables
TYPE(FIELD_VARIABLE_TYPE), POINTER :: FIELD_VARIABLE !<A pointer to the field variable to calculate the integral errors for
REAL(DP), INTENT(OUT) :: INTEGRAL_ERRORS(:,:) !<On exit, the integral errors for the local elements
REAL(DP), INTENT(OUT) :: GHOST_INTEGRAL_ERRORS(:,:) !<On exit, the integral errors for the ghost elements
INTEGER(INTG), INTENT(OUT) :: ERR !<The error code
TYPE(VARYING_STRING), INTENT(OUT) :: ERROR !<The error string
!Local Variables
INTEGER(INTG) :: element_idx,component_idx,gauss_idx,parameter_idx,variable_type
REAL(DP) :: ANALYTIC_INT,NUMERICAL_INT,RWG
TYPE(BASIS_TYPE), POINTER :: BASIS,DEPENDENT_BASIS,GEOMETRIC_BASIS
TYPE(DECOMPOSITION_TYPE), POINTER :: DECOMPOSITION
TYPE(DOMAIN_ELEMENTS_TYPE), POINTER :: DOMAIN_ELEMENTS1,DOMAIN_ELEMENTS2,DOMAIN_ELEMENTS3
TYPE(FIELD_TYPE), POINTER :: DEPENDENT_FIELD,GEOMETRIC_FIELD
TYPE(FIELD_INTERPOLATED_POINT_PTR_TYPE), POINTER :: GEOMETRIC_INTERP_POINT(:)
TYPE(FIELD_INTERPOLATED_POINT_METRICS_PTR_TYPE), POINTER :: GEOMETRIC_INTERP_POINT_METRICS(:)
TYPE(FIELD_INTERPOLATION_PARAMETERS_PTR_TYPE), POINTER :: ANALYTIC_INTERP_PARAMETERS(:),GEOMETRIC_INTERP_PARAMETERS(:), &
& NUMERICAL_INTERP_PARAMETERS(:)
TYPE(FIELD_VARIABLE_TYPE), POINTER :: GEOMETRIC_VARIABLE
TYPE(QUADRATURE_SCHEME_TYPE), POINTER :: QUADRATURE_SCHEME
TYPE(VARYING_STRING) :: LOCAL_ERROR
NULLIFY(GEOMETRIC_INTERP_POINT)
NULLIFY(GEOMETRIC_INTERP_POINT_METRICS)
NULLIFY(ANALYTIC_INTERP_PARAMETERS)
NULLIFY(GEOMETRIC_INTERP_PARAMETERS)
NULLIFY(NUMERICAL_INTERP_PARAMETERS)
ENTERS("ANALYTIC_ANALYSIS_INTEGRAL_ERRORS",ERR,ERROR,*999)
INTEGRAL_ERRORS=0.0_DP
GHOST_INTEGRAL_ERRORS=0.0_DP
IF(ASSOCIATED(FIELD_VARIABLE)) THEN
IF(SIZE(INTEGRAL_ERRORS,1)>=6.AND.SIZE(INTEGRAL_ERRORS,2)>=FIELD_VARIABLE%NUMBER_OF_COMPONENTS) THEN
IF(SIZE(GHOST_INTEGRAL_ERRORS,1)>=6.AND.SIZE(GHOST_INTEGRAL_ERRORS,2)>=FIELD_VARIABLE%NUMBER_OF_COMPONENTS) THEN
variable_type=FIELD_VARIABLE%VARIABLE_TYPE
DEPENDENT_FIELD=>FIELD_VARIABLE%FIELD
IF(ASSOCIATED(DEPENDENT_FIELD)) THEN
DECOMPOSITION=>DEPENDENT_FIELD%DECOMPOSITION
IF(ASSOCIATED(DECOMPOSITION)) THEN
GEOMETRIC_FIELD=>DEPENDENT_FIELD%GEOMETRIC_FIELD
IF(ASSOCIATED(GEOMETRIC_FIELD)) THEN
GEOMETRIC_VARIABLE=>GEOMETRIC_FIELD%VARIABLE_TYPE_MAP(FIELD_U_VARIABLE_TYPE)%PTR
IF(ASSOCIATED(GEOMETRIC_VARIABLE)) THEN
CALL FIELD_INTERPOLATION_PARAMETERS_INITIALISE(GEOMETRIC_FIELD,GEOMETRIC_INTERP_PARAMETERS,ERR,ERROR,*999)
CALL FIELD_INTERPOLATION_PARAMETERS_INITIALISE(DEPENDENT_FIELD,NUMERICAL_INTERP_PARAMETERS,ERR,ERROR,*999)
CALL FIELD_INTERPOLATION_PARAMETERS_INITIALISE(DEPENDENT_FIELD,ANALYTIC_INTERP_PARAMETERS,ERR,ERROR,*999)
CALL FIELD_INTERPOLATED_POINTS_INITIALISE(GEOMETRIC_INTERP_PARAMETERS,GEOMETRIC_INTERP_POINT,ERR,ERROR,*999)
CALL Field_InterpolatedPointsMetricsInitialise(GEOMETRIC_INTERP_POINT,GEOMETRIC_INTERP_POINT_METRICS, &
& ERR,ERROR,*999)
DOMAIN_ELEMENTS1=>FIELD_VARIABLE%COMPONENTS(DECOMPOSITION%MESH_COMPONENT_NUMBER)%DOMAIN%TOPOLOGY%ELEMENTS
DOMAIN_ELEMENTS2=>GEOMETRIC_VARIABLE%COMPONENTS(DECOMPOSITION%MESH_COMPONENT_NUMBER)%DOMAIN%TOPOLOGY%ELEMENTS
DO element_idx=1,DOMAIN_ELEMENTS1%NUMBER_OF_ELEMENTS
DEPENDENT_BASIS=>DOMAIN_ELEMENTS1%ELEMENTS(element_idx)%BASIS
GEOMETRIC_BASIS=>DOMAIN_ELEMENTS2%ELEMENTS(element_idx)%BASIS
QUADRATURE_SCHEME=>DEPENDENT_BASIS%QUADRATURE%QUADRATURE_SCHEME_MAP(BASIS_DEFAULT_QUADRATURE_SCHEME)%PTR
CALL FIELD_INTERPOLATION_PARAMETERS_ELEMENT_GET(FIELD_VALUES_SET_TYPE,element_idx, &
& GEOMETRIC_INTERP_PARAMETERS(FIELD_U_VARIABLE_TYPE)%PTR,ERR,ERROR,*999)
CALL FIELD_INTERPOLATION_PARAMETERS_ELEMENT_GET(FIELD_VALUES_SET_TYPE,element_idx, &
& NUMERICAL_INTERP_PARAMETERS(variable_type)%PTR,ERR,ERROR,*999)
CALL FIELD_INTERPOLATION_PARAMETERS_ELEMENT_GET(FIELD_ANALYTIC_VALUES_SET_TYPE,element_idx, &
& ANALYTIC_INTERP_PARAMETERS(variable_type)%PTR,ERR,ERROR,*999)
DO gauss_idx=1,QUADRATURE_SCHEME%NUMBER_OF_GAUSS
CALL FIELD_INTERPOLATE_GAUSS(FIRST_PART_DERIV,BASIS_DEFAULT_QUADRATURE_SCHEME,gauss_idx, &
& GEOMETRIC_INTERP_POINT(FIELD_U_VARIABLE_TYPE)%PTR,ERR,ERROR,*999)
CALL FIELD_INTERPOLATED_POINT_METRICS_CALCULATE(GEOMETRIC_BASIS%NUMBER_OF_XI, &
& GEOMETRIC_INTERP_POINT_METRICS(FIELD_U_VARIABLE_TYPE)%PTR,ERR,ERROR,*999)
RWG=GEOMETRIC_INTERP_POINT_METRICS(FIELD_U_VARIABLE_TYPE)%PTR%JACOBIAN* &
& QUADRATURE_SCHEME%GAUSS_WEIGHTS(gauss_idx)
DO component_idx=1,FIELD_VARIABLE%NUMBER_OF_COMPONENTS
DOMAIN_ELEMENTS3=>FIELD_VARIABLE%COMPONENTS(component_idx)%DOMAIN%TOPOLOGY%ELEMENTS
BASIS=>DOMAIN_ELEMENTS3%ELEMENTS(element_idx)%BASIS
NUMERICAL_INT=0.0_DP
ANALYTIC_INT=0.0_DP
SELECT CASE(DEPENDENT_FIELD%SCALINGS%SCALING_TYPE)
CASE(FIELD_NO_SCALING)
DO parameter_idx=1,BASIS%NUMBER_OF_ELEMENT_PARAMETERS
NUMERICAL_INT=NUMERICAL_INT+QUADRATURE_SCHEME%GAUSS_BASIS_FNS(parameter_idx,NO_PART_DERIV,gauss_idx)* &
& NUMERICAL_INTERP_PARAMETERS(variable_type)%PTR%PARAMETERS(parameter_idx,component_idx)
ANALYTIC_INT=ANALYTIC_INT+QUADRATURE_SCHEME%GAUSS_BASIS_FNS(parameter_idx,NO_PART_DERIV,gauss_idx)* &
& ANALYTIC_INTERP_PARAMETERS(variable_type)%PTR%PARAMETERS(parameter_idx,component_idx)
ENDDO !parameter_idx
CASE(FIELD_UNIT_SCALING,FIELD_ARC_LENGTH_SCALING,FIELD_ARITHMETIC_MEAN_SCALING,FIELD_HARMONIC_MEAN_SCALING)
DO parameter_idx=1,BASIS%NUMBER_OF_ELEMENT_PARAMETERS
NUMERICAL_INT=NUMERICAL_INT+QUADRATURE_SCHEME%GAUSS_BASIS_FNS(parameter_idx,NO_PART_DERIV,gauss_idx)* &
& NUMERICAL_INTERP_PARAMETERS(variable_type)%PTR%PARAMETERS(parameter_idx,component_idx)* &
& NUMERICAL_INTERP_PARAMETERS(variable_type)%PTR%SCALE_FACTORS(parameter_idx,component_idx)
ANALYTIC_INT=ANALYTIC_INT+QUADRATURE_SCHEME%GAUSS_BASIS_FNS(parameter_idx,NO_PART_DERIV,gauss_idx)* &
& ANALYTIC_INTERP_PARAMETERS(variable_type)%PTR%PARAMETERS(parameter_idx,component_idx)* &
& ANALYTIC_INTERP_PARAMETERS(variable_type)%PTR%SCALE_FACTORS(parameter_idx,component_idx)
ENDDO !parameter_idx
CASE DEFAULT
LOCAL_ERROR="The dependent field scaling type of "// &
& TRIM(NUMBER_TO_VSTRING(DEPENDENT_FIELD%SCALINGS%SCALING_TYPE,"*",ERR,ERROR))//" is invalid."
CALL FlagError(LOCAL_ERROR,ERR,ERROR,*999)
END SELECT
INTEGRAL_ERRORS(1,component_idx)=INTEGRAL_ERRORS(1,component_idx)+NUMERICAL_INT*RWG
INTEGRAL_ERRORS(2,component_idx)=INTEGRAL_ERRORS(2,component_idx)+NUMERICAL_INT**2*RWG
INTEGRAL_ERRORS(3,component_idx)=INTEGRAL_ERRORS(3,component_idx)+ANALYTIC_INT*RWG
INTEGRAL_ERRORS(4,component_idx)=INTEGRAL_ERRORS(4,component_idx)+ANALYTIC_INT**2*RWG
INTEGRAL_ERRORS(5,component_idx)=INTEGRAL_ERRORS(5,component_idx)+(ANALYTIC_INT-NUMERICAL_INT)*RWG
INTEGRAL_ERRORS(6,component_idx)=INTEGRAL_ERRORS(6,component_idx)+(ANALYTIC_INT-NUMERICAL_INT)**2*RWG
ENDDO !component_idx
ENDDO !gauss_idx
ENDDO !element_idx
DO element_idx=DOMAIN_ELEMENTS1%NUMBER_OF_ELEMENTS+1,DOMAIN_ELEMENTS1%TOTAL_NUMBER_OF_ELEMENTS
DEPENDENT_BASIS=>DOMAIN_ELEMENTS1%ELEMENTS(element_idx)%BASIS
GEOMETRIC_BASIS=>DOMAIN_ELEMENTS2%ELEMENTS(element_idx)%BASIS
QUADRATURE_SCHEME=>DEPENDENT_BASIS%QUADRATURE%QUADRATURE_SCHEME_MAP(BASIS_DEFAULT_QUADRATURE_SCHEME)%PTR
CALL FIELD_INTERPOLATION_PARAMETERS_ELEMENT_GET(FIELD_VALUES_SET_TYPE,element_idx, &
& GEOMETRIC_INTERP_PARAMETERS(FIELD_U_VARIABLE_TYPE)%PTR,ERR,ERROR,*999)
CALL FIELD_INTERPOLATION_PARAMETERS_ELEMENT_GET(FIELD_VALUES_SET_TYPE,element_idx, &
& NUMERICAL_INTERP_PARAMETERS(variable_type)%PTR,ERR,ERROR,*999)
CALL FIELD_INTERPOLATION_PARAMETERS_ELEMENT_GET(FIELD_ANALYTIC_VALUES_SET_TYPE,element_idx, &
& ANALYTIC_INTERP_PARAMETERS(variable_type)%PTR,ERR,ERROR,*999)
DO gauss_idx=1,QUADRATURE_SCHEME%NUMBER_OF_GAUSS
CALL FIELD_INTERPOLATE_GAUSS(FIRST_PART_DERIV,BASIS_DEFAULT_QUADRATURE_SCHEME,gauss_idx, &
& GEOMETRIC_INTERP_POINT(FIELD_U_VARIABLE_TYPE)%PTR,ERR,ERROR,*999)
CALL FIELD_INTERPOLATED_POINT_METRICS_CALCULATE(GEOMETRIC_BASIS%NUMBER_OF_XI, &
& GEOMETRIC_INTERP_POINT_METRICS(FIELD_U_VARIABLE_TYPE)%PTR,ERR,ERROR,*999)
RWG=GEOMETRIC_INTERP_POINT_METRICS(FIELD_U_VARIABLE_TYPE)%PTR%JACOBIAN* &
& QUADRATURE_SCHEME%GAUSS_WEIGHTS(gauss_idx)
DO component_idx=1,FIELD_VARIABLE%NUMBER_OF_COMPONENTS
DOMAIN_ELEMENTS3=>FIELD_VARIABLE%COMPONENTS(component_idx)%DOMAIN%TOPOLOGY%ELEMENTS
BASIS=>DOMAIN_ELEMENTS3%ELEMENTS(element_idx)%BASIS
NUMERICAL_INT=0.0_DP
ANALYTIC_INT=0.0_DP
SELECT CASE(DEPENDENT_FIELD%SCALINGS%SCALING_TYPE)
CASE(FIELD_NO_SCALING)
DO parameter_idx=1,BASIS%NUMBER_OF_ELEMENT_PARAMETERS
NUMERICAL_INT=NUMERICAL_INT+QUADRATURE_SCHEME%GAUSS_BASIS_FNS(parameter_idx,NO_PART_DERIV,gauss_idx)* &
& NUMERICAL_INTERP_PARAMETERS(variable_type)%PTR%PARAMETERS(parameter_idx,component_idx)
ANALYTIC_INT=ANALYTIC_INT+QUADRATURE_SCHEME%GAUSS_BASIS_FNS(parameter_idx,NO_PART_DERIV,gauss_idx)* &
& ANALYTIC_INTERP_PARAMETERS(variable_type)%PTR%PARAMETERS(parameter_idx,component_idx)
ENDDO !parameter_idx
CASE(FIELD_UNIT_SCALING,FIELD_ARC_LENGTH_SCALING,FIELD_ARITHMETIC_MEAN_SCALING,FIELD_HARMONIC_MEAN_SCALING)
DO parameter_idx=1,BASIS%NUMBER_OF_ELEMENT_PARAMETERS
NUMERICAL_INT=NUMERICAL_INT+QUADRATURE_SCHEME%GAUSS_BASIS_FNS(parameter_idx,NO_PART_DERIV,gauss_idx)* &
& NUMERICAL_INTERP_PARAMETERS(variable_type)%PTR%PARAMETERS(parameter_idx,component_idx)* &
& NUMERICAL_INTERP_PARAMETERS(variable_type)%PTR%SCALE_FACTORS(parameter_idx,component_idx)
ANALYTIC_INT=ANALYTIC_INT+QUADRATURE_SCHEME%GAUSS_BASIS_FNS(parameter_idx,NO_PART_DERIV,gauss_idx)* &
& ANALYTIC_INTERP_PARAMETERS(variable_type)%PTR%PARAMETERS(parameter_idx,component_idx)* &
& ANALYTIC_INTERP_PARAMETERS(variable_type)%PTR%SCALE_FACTORS(parameter_idx,component_idx)
ENDDO !parameter_idx
CASE DEFAULT
LOCAL_ERROR="The dependent field scaling type of "// &
& TRIM(NUMBER_TO_VSTRING(DEPENDENT_FIELD%SCALINGS%SCALING_TYPE,"*",ERR,ERROR))//" is invalid."
CALL FlagError(LOCAL_ERROR,ERR,ERROR,*999)
END SELECT
GHOST_INTEGRAL_ERRORS(1,component_idx)=GHOST_INTEGRAL_ERRORS(1,component_idx)+NUMERICAL_INT*RWG
GHOST_INTEGRAL_ERRORS(2,component_idx)=GHOST_INTEGRAL_ERRORS(2,component_idx)+NUMERICAL_INT**2*RWG
GHOST_INTEGRAL_ERRORS(3,component_idx)=GHOST_INTEGRAL_ERRORS(3,component_idx)+ANALYTIC_INT*RWG
GHOST_INTEGRAL_ERRORS(4,component_idx)=GHOST_INTEGRAL_ERRORS(4,component_idx)+ANALYTIC_INT**2*RWG
GHOST_INTEGRAL_ERRORS(5,component_idx)=GHOST_INTEGRAL_ERRORS(5,component_idx)+ &
& (ANALYTIC_INT-NUMERICAL_INT)*RWG
GHOST_INTEGRAL_ERRORS(6,component_idx)=GHOST_INTEGRAL_ERRORS(6,component_idx)+ &
& (ANALYTIC_INT-NUMERICAL_INT)**2*RWG
ENDDO !component_idx
ENDDO !gauss_idx
ENDDO !element_idx
CALL Field_InterpolatedPointsMetricsFinalise(GEOMETRIC_INTERP_POINT_METRICS,ERR,ERROR,*999)
CALL FIELD_INTERPOLATED_POINTS_FINALISE(GEOMETRIC_INTERP_POINT,ERR,ERROR,*999)
CALL FIELD_INTERPOLATION_PARAMETERS_FINALISE(ANALYTIC_INTERP_PARAMETERS,ERR,ERROR,*999)
CALL FIELD_INTERPOLATION_PARAMETERS_FINALISE(NUMERICAL_INTERP_PARAMETERS,ERR,ERROR,*999)
CALL FIELD_INTERPOLATION_PARAMETERS_FINALISE(GEOMETRIC_INTERP_PARAMETERS,ERR,ERROR,*999)
ELSE
CALL FlagError("Geometric field variable is not associated.",ERR,ERROR,*999)
ENDIF
ELSE
CALL FlagError("Field geometric field is not associated.",ERR,ERROR,*999)
ENDIF
ELSE
CALL FlagError("Field decomposition is not associated.",ERR,ERROR,*999)
ENDIF
ELSE
CALL FlagError("Field variable field is not associated.",ERR,ERROR,*999)
ENDIF
ELSE
LOCAL_ERROR="Invalid size for GHOST_INTEGRAL_ERRORS. The size is ("// &
& TRIM(NUMBER_TO_VSTRING(SIZE(GHOST_INTEGRAL_ERRORS,1),"*",ERR,ERROR))//","// &
& TRIM(NUMBER_TO_VSTRING(SIZE(GHOST_INTEGRAL_ERRORS,2),"*",ERR,ERROR))//") and it needs to be at least (6,"// &
& TRIM(NUMBER_TO_VSTRING(FIELD_VARIABLE%NUMBER_OF_COMPONENTS,"*",ERR,ERROR))//")."
CALL FlagError(LOCAL_ERROR,ERR,ERROR,*999)
ENDIF
ELSE
LOCAL_ERROR="Invalid size for INTEGRAL_ERRORS. The size is ("// &
& TRIM(NUMBER_TO_VSTRING(SIZE(INTEGRAL_ERRORS,1),"*",ERR,ERROR))//","// &
& TRIM(NUMBER_TO_VSTRING(SIZE(INTEGRAL_ERRORS,2),"*",ERR,ERROR))//") and it needs to be at least (6,"// &
& TRIM(NUMBER_TO_VSTRING(FIELD_VARIABLE%NUMBER_OF_COMPONENTS,"*",ERR,ERROR))//")."
CALL FlagError(LOCAL_ERROR,ERR,ERROR,*999)
ENDIF
ELSE
CALL FlagError("Field variable is not associated.",ERR,ERROR,*999)
ENDIF
EXITS("ANALYTIC_ANALYSIS_INTEGRAL_ERRORS")
RETURN
999 ERRORSEXITS("ANALYTIC_ANALYSIS_INTEGRAL_ERRORS",ERR,ERROR)
RETURN 1
END SUBROUTINE ANALYTIC_ANALYSIS_INTEGRAL_ERRORS
!
!================================================================================================================================
!
!>Get integral absolute error value for the field
SUBROUTINE AnalyticAnalysis_IntegralAbsoluteErrorGet(FIELD,VARIABLE_TYPE,COMPONENT_NUMBER,INTEGRAL_ERROR, &
& GHOST_INTEGRAL_ERROR,ERR,ERROR,*)