~reducedmodelling/fluidity/ROM_Non-intrusive-ann

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
!    Copyright (C) 2006 Imperial College London and others.
!    
!    Please see the AUTHORS file in the main source directory for a full list
!    of copyright holders.
!
!    Prof. C Pain
!    Applied Modelling and Computation Group
!    Department of Earth Science and Engineering
!    Imperial College London
!
!    amcgsoftware@imperial.ac.uk
!    
!    This library is free software; you can redistribute it and/or
!    modify it under the terms of the GNU Lesser General Public
!    License as published by the Free Software Foundation,
!    version 2.1 of the License.
!
!    This library is distributed in the hope that it will be useful,
!    but WITHOUT ANY WARRANTY; without even the implied warranty of
!    MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
!    Lesser General Public License for more details.
!
!    You should have received a copy of the GNU Lesser General Public
!    License along with this library; if not, write to the Free Software
!    Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307
!    USA

#include "fdebug.h"

module diagnostic_fields_wrapper
  !!< A module to link to diagnostic variable calculations.

  use global_parameters, only:FIELD_NAME_LEN 
  use fields
  use sparse_matrices_fields
  use field_derivatives
  use state_module
  use futils
  use fetools
  use spud
  use parallel_tools
  use diagnostic_fields, only: calculate_diagnostic_variable
  use multimaterial_module, only: calculate_material_mass, &
                                  calculate_bulk_material_pressure, &
                                  calculate_sum_material_volume_fractions, &
                                  calculate_material_volume
  use free_surface_module, only: calculate_diagnostic_free_surface, &
                                 calculate_diagnostic_wettingdrying_alpha
  use tidal_module, only: calculate_diagnostic_equilibrium_pressure
  use field_options, only: do_not_recalculate
  use vorticity_diagnostics
  use diagnostic_fields_matrices
  use equation_of_state
  use momentum_diagnostic_fields
  use spontaneous_potentials, only: calculate_formation_conductivity
  use sediment_diagnostics
  use geostrophic_pressure
  use multiphase_module
  
  implicit none

  private
  
  public :: calculate_diagnostic_variables

contains

  subroutine calculate_diagnostic_variables(state, exclude_nonrecalculated)
    !!< Updates diagnostic fields in the supplied states.

    type(state_type), dimension(:) :: state
    logical, intent(in), optional :: exclude_nonrecalculated

    integer :: i,stat
    type(scalar_field), pointer :: s_field
    type(vector_field), pointer :: v_field
    logical :: diagnostic

    ! An array of submaterials of the current phase in state(istate).
    type(state_type), dimension(:), pointer :: submaterials
    
    ewrite(1, *) "In calculate_diagnostic_variables"
 
    do i = 1, size(state)

       ! start of fields that can be called through the generic calculate_diagnostic_variable
       ! interface, i.e. - those that only require things available in f90modules

       s_field => extract_scalar_field(state(i), "CFLNumber", stat)
       if(stat == 0) then
         if(recalculate(trim(s_field%option_path))) then
           call calculate_diagnostic_variable(state(i), "CFLNumber", &
             & s_field)
         end if
       end if

       s_field => extract_scalar_field(state(i), "GridReynoldsNumber", stat)
       if(stat == 0) then
         if(recalculate(trim(s_field%option_path))) then
           call calculate_diagnostic_variable(state(i), "GridReynoldsNumber", &
             & s_field)
         end if
       end if

       s_field => extract_scalar_field(state(i), "GridPecletNumber", stat)
       if(stat == 0) then
         if(recalculate(trim(s_field%option_path))) then
           call calculate_diagnostic_variable(state(i), "GridPecletNumber", &
             & s_field)
         end if
       end if

       s_field => extract_scalar_field(state(i), "ControlVolumeCFLNumber", stat)
       if(stat == 0) then
         if(recalculate(trim(s_field%option_path))) then
           call calculate_diagnostic_variable(state(i), "ControlVolumeCFLNumber", &
             & s_field)
         end if
       end if

       s_field => extract_scalar_field(state(i), "DG_CourantNumber", stat)
       if(stat == 0) then
         if(recalculate(trim(s_field%option_path))) then
           call calculate_diagnostic_variable(state(i), "DG_CourantNumber", &
             & s_field)
         end if
       end if

       s_field => extract_scalar_field(state(i), "CVMaterialDensityCFLNumber", stat)
       if(stat == 0) then
         if(recalculate(trim(s_field%option_path))) then
           call calculate_diagnostic_variable(state(i), "CVMaterialDensityCFLNumber", &
             & s_field)
         end if
       end if

       s_field => extract_scalar_field(state(i), "InterstitialVelocityCGCourantNumber", stat)
       if(stat == 0) then
         if(recalculate(trim(s_field%option_path))) then
           call calculate_diagnostic_variable(state(i), "InterstitialVelocityCGCourantNumber", &
             & s_field)
         end if
       end if

       s_field => extract_scalar_field(state(i), "InterstitialVelocityDGCourantNumber", stat)
       if(stat == 0) then
         if(recalculate(trim(s_field%option_path))) then
           call calculate_diagnostic_variable(state(i), "InterstitialVelocityDGCourantNumber", &
             & s_field)
         end if
       end if

       s_field => extract_scalar_field(state(i), "InterstitialVelocityCVCourantNumber", stat)
       if(stat == 0) then
         if(recalculate(trim(s_field%option_path))) then
           call calculate_diagnostic_variable(state(i), "InterstitialVelocityCVCourantNumber", &
             & s_field)
         end if
       end if

       s_field => extract_scalar_field(state(i), "KineticEnergyDensity", stat)
       if(stat == 0) then
         if(recalculate(trim(s_field%option_path))) then
           call calculate_diagnostic_variable(state(i), "KineticEnergyDensity", &
             & s_field)
         end if
       end if

       s_field => extract_scalar_field(state(i),  "HorizontalVelocityDivergence", stat)
       if(stat == 0) then
         if(recalculate(trim(s_field%option_path))) then
           call calculate_diagnostic_variable(state(i), "HorizontalVelocityDivergence", s_field)
         end if
       end if

       s_field => extract_scalar_field(state(i), &
         & "VelocityDivergence", stat)
       if(stat == 0) then
         if(recalculate(trim(s_field%option_path))) then
           call calculate_diagnostic_variable(state(i), "VelocityDivergence", s_field)
         end if
       end if

       s_field => extract_scalar_field(state(i), &
         & "PerturbationDensity", stat)
       if(stat == 0) then
         if(recalculate(trim(s_field%option_path))) then
           ! this routine returns the density used in the buoyancy term, which we're not really interested in
           ! but it computes the PerturbationDensity as a side effect. Note that this means it will happen twice
           ! as it will be recalculated at the beginning of Momemtum_Equation after the Temperature and Salinity
           ! fields have been solved for.
           call calculate_perturbation_density(state(i), s_field)
         end if
       end if

       ! this diagnostic field depends on PerturbationDensity
       s_field => extract_scalar_field(state(i), "GravitationalPotentialEnergyDensity", stat)
       if(stat == 0) then
         if(recalculate(trim(s_field%option_path))) then
           call calculate_diagnostic_variable(state(i), "GravitationalPotentialEnergyDensity", s_field)
         end if
       end if
       
       ! this diagnostic field depends on PerturbationDensity
       s_field => extract_scalar_field(state(i), "IsopycnalCoordinate", stat)
       if(stat == 0) then
         if(recalculate(trim(s_field%option_path))) then
           call calculate_diagnostic_variable(state(i), "IsopycnalCoordinate", &
             & s_field)
         end if
       end if
       
       ! Must be calculated after IsopycnalCoordinate
       s_field => extract_scalar_field(state(i), "BackgroundPotentialEnergyDensity", stat)
       if(stat == 0) then
         if(recalculate(trim(s_field%option_path))) then
           call calculate_diagnostic_variable(state(i), "BackgroundPotentialEnergyDensity", s_field)
         end if
       end if
       
       v_field => extract_vector_field(state(i), "InnerElementFullVelocity", stat)
       if(stat == 0) then
         if(recalculate(trim(v_field%option_path))) then
           call calculate_diagnostic_variable(state(i), "InnerElementFullVelocity", &
             & v_field)
         end if
       end if

       ! Must be calculated after InnerElementFullVelocity
       v_field => extract_vector_field(state(i), "InnerElementFullVorticity", stat)
       if(stat == 0) then
         if(recalculate(trim(v_field%option_path))) then
           call calculate_diagnostic_variable(state(i), "InnerElementFullVorticity", &
             & v_field)
         end if
       end if

       v_field => extract_vector_field(state(i), "InnerElementVorticity", stat)
       if(stat == 0) then
         if(recalculate(trim(v_field%option_path))) then
           call calculate_diagnostic_variable(state(i), "InnerElementVorticity", &
             & v_field)
         end if
       end if

       v_field => extract_vector_field(state(i), "DgMappedVelocity", stat)
       if(stat == 0) then
         if(recalculate(trim(v_field%option_path))) then
           call calculate_diagnostic_variable(state(i), "DgMappedVelocity", &
             & v_field)
         end if
       end if

       v_field => extract_vector_field(state(i), "DgMappedVorticity", stat)
       if(stat == 0) then
         if(recalculate(trim(v_field%option_path))) then
           call calculate_diagnostic_variable(state(i), "DgMappedVorticity", &
             & v_field)
         end if
       end if
       
       s_field => extract_scalar_field(state(i), "HorizontalStreamFunction", stat)
       if(stat == 0) then
         if(recalculate(trim(s_field%option_path))) then
           call calculate_diagnostic_variable(state(i), "HorizontalStreamFunction", s_field)
         end if
       end if
       
       s_field => extract_scalar_field(state(i), "Speed", stat)
       if(stat == 0) then
         if(recalculate(trim(s_field%option_path))) then
           call calculate_diagnostic_variable(state(i), "Speed", &
             & s_field)
         end if
       end if
       
       s_field => extract_scalar_field(state(i), "DiffusiveDissipation", stat)
       if(stat == 0) then
         if(recalculate(trim(s_field%option_path))) then
           call calculate_diagnostic_variable(state(i), "DiffusiveDissipation", &
             & s_field)
         end if
       end if
       
       s_field => extract_scalar_field(state(i), "ViscousDissipation", stat)
       if(stat == 0) then
         if(recalculate(trim(s_field%option_path))) then
           call calculate_diagnostic_variable(state(i), "ViscousDissipation", &
             & s_field)
         end if
       end if

       s_field => extract_scalar_field(state(i), "RichardsonNumber", stat)
       if(stat == 0) then
         if(recalculate(trim(s_field%option_path))) then
           call calculate_diagnostic_variable(state(i), "RichardsonNumber", &
             & s_field)
         end if
       end if

       s_field => extract_scalar_field(state(i), "StreamFunction", stat)
       if(stat == 0) then
         if(recalculate(trim(s_field%option_path))) then
           call calculate_diagnostic_variable(state(i), "StreamFunction", s_field)
         end if
       end if

       s_field => extract_scalar_field(state(i), "MultiplyConnectedStreamFunction", stat)
       if(stat == 0) then
         if(recalculate(trim(s_field%option_path))) then
           call calculate_diagnostic_variable(state(i), "MultiplyConnectedStreamFunction", s_field)
         end if
       end if

       s_field => extract_scalar_field(state(i), "Time", stat)
       if(stat == 0) then
         if(recalculate(trim(s_field%option_path))) then
           call calculate_diagnostic_variable(state(i), "Time", s_field)
         end if
       end if

       v_field => extract_vector_field(state(i), "LinearMomentum", stat)
       if(stat == 0) then
         if(recalculate(trim(v_field%option_path))) then
           call calculate_diagnostic_variable(state(i), "LinearMomentum", v_field)
         end if
       end if

       v_field => extract_vector_field(state(i), "DiagnosticCoordinate", stat)
       if(stat == 0) then
         if(recalculate(trim(v_field%option_path))) then
           call calculate_diagnostic_variable(state(i), "DiagnosticCoordinate", v_field)
         end if
       end if

       v_field => extract_vector_field(state(i), "BedShearStress", stat)  
       if(stat == 0) then  
         if(recalculate(trim(v_field%option_path))) then
           call calculate_diagnostic_variable(state(i), "BedShearStress", v_field)  
         end if
       end if

       v_field => extract_vector_field(state(i), "MaxBedShearStress", stat)  
       if(stat == 0) then  
         if(recalculate(trim(v_field%option_path))) then
           call calculate_diagnostic_variable(state(i), "MaxBedShearStress", v_field)  
         end if
       end if

       s_field => extract_scalar_field(state(i), "GalerkinProjection", stat)
       if(stat == 0) then
         if(recalculate(trim(s_field%option_path))) then
           call calculate_diagnostic_variable(state(i), "GalerkinProjection", s_field)
         end if
       end if

       v_field => extract_vector_field(state(i), "GalerkinProjection", stat)
       if(stat == 0) then
         if(recalculate(trim(v_field%option_path))) then
           call calculate_diagnostic_variable(state(i), "GalerkinProjection", v_field)
         end if
       end if

       s_field => extract_scalar_field(state(i), "UniversalNumber", stat)
       if(stat == 0) then
          call calculate_diagnostic_variable(state(i), "UniversalNumber", s_field)
       end if

       s_field => extract_scalar_field(state(i), "NodeOwner", stat)
       if(stat == 0) then
          call calculate_diagnostic_variable(state(i), "NodeOwner", s_field)
       end if

       ! end of fields that can be called through the generic calculate_diagnostic_variable
       ! interface, i.e. - those that only require things available in f90modules

       ! start of fields that cannot be called through the generic calculate_diagnostic_variable
       ! interface, i.e. - those that need things from assemble
       s_field => extract_scalar_field(state(i), "ControlVolumeDivergence", stat)
       if(stat == 0) then
         if(recalculate(trim(s_field%option_path))) then
           call calculate_divergence_cv(state(i), s_field)
         end if
       end if

       s_field => extract_scalar_field(state(i), "FiniteElementDivergence", stat)
       if(stat == 0) then
         if(recalculate(trim(s_field%option_path))) then
           call calculate_divergence_fe(state(i), s_field)
         end if
       end if

       s_field => extract_scalar_field(state(i), "MaterialVolumeFraction", stat)
       if(stat == 0) then
         diagnostic = have_option(trim(s_field%option_path)//"/diagnostic")
         if(diagnostic .and. .not.(aliased(s_field))) then
           if(recalculate(trim(s_field%option_path))) then
             call calculate_sum_material_volume_fractions(state, s_field)
             call scale(s_field, -1.0)
             call addto(s_field, 1.0)
           end if
         end if
       end if
       
       s_field => extract_scalar_field(state(i), "MaterialMass", stat)
       if(stat == 0) then
         if(recalculate(trim(s_field%option_path))) then
           call calculate_material_mass(state(i), s_field)
         end if
       end if

       s_field => extract_scalar_field(state(i), "MaterialVolume", stat)
       if(stat == 0) then
         if(recalculate(trim(s_field%option_path))) then
           call calculate_material_volume(state(i), s_field)
         end if
       end if
       
       s_field => extract_scalar_field(state(i), "MaterialDensity", stat)
       if(stat == 0) then
         diagnostic = have_option(trim(s_field%option_path)//"/diagnostic")
         if(diagnostic .and. .not.(aliased(s_field))) then
           if(recalculate(trim(s_field%option_path))) then
            call calculate_densities(state(i), bulk_density=s_field)
           end if
         end if
       end if

       s_field => extract_scalar_field(state(i), "Density", stat)
       if(stat == 0) then
         diagnostic = have_option(trim(s_field%option_path)//"/diagnostic")
         if(diagnostic .and. .not.(aliased(s_field))) then
           if(recalculate(trim(s_field%option_path))) then
             if(option_count("/material_phase/vector_field::Velocity/prognostic") > 1) then 
               call get_phase_submaterials(state, i, submaterials)
               call calculate_densities(submaterials, bulk_density=s_field)
               deallocate(submaterials)
             else
               call calculate_densities(state, bulk_density=s_field)
             end if
           end if
         end if
       end if

       s_field => extract_scalar_field(state(i), "MaterialEOSDensity", stat)
       if(stat == 0) then
         if(recalculate(trim(s_field%option_path))) then
           call compressible_material_eos(state(i), materialdensity=s_field)
         end if
       end if

       s_field => extract_scalar_field(state(i), "MaterialPressure", stat)
       if(stat == 0) then
         if(recalculate(trim(s_field%option_path))) then
           call compressible_material_eos(state(i), materialpressure=s_field)
         end if
       end if

       s_field => extract_scalar_field(state(i), "BulkMaterialPressure", stat)
       if(stat == 0) then
         if(recalculate(trim(s_field%option_path))) then
           call calculate_bulk_material_pressure(state, s_field)
         end if
       end if

       s_field => extract_scalar_field(state(i), "SumMaterialVolumeFractions", stat)
       if(stat == 0) then
         if(recalculate(trim(s_field%option_path))) then
           call calculate_sum_material_volume_fractions(state, s_field)
         end if
       end if

       s_field => extract_scalar_field(state(i), "FreeSurface", stat)
       if(stat == 0) then
         if(recalculate(trim(s_field%option_path))) then
           call calculate_diagnostic_free_surface(state(i), s_field)
         end if
       end if
       
       s_field => extract_scalar_field(state(i), "WettingDryingAlpha", stat)
       if(stat == 0) then
         if(recalculate(trim(s_field%option_path))) then
           call calculate_diagnostic_wettingdrying_alpha(state(i), s_field)
         end if
       end if

       s_field => extract_scalar_field(state(i), "EquilibriumPressure", stat)
       if(stat == 0) then
         if(recalculate(trim(s_field%option_path))) then
           call calculate_diagnostic_equilibrium_pressure(state(i), s_field)
         end if
       end if

       v_field => extract_vector_field(state(i), "ControlVolumeDivergenceTransposed", stat)
       if(stat == 0) then
         if(recalculate(trim(v_field%option_path))) then
           call calculate_div_t_cv(state(i), v_field)
         end if
       end if

       v_field => extract_vector_field(state(i), "FiniteElementGradient", stat)
       if(stat == 0) then
         if(recalculate(trim(v_field%option_path))) then
           call calculate_grad_fe(state(i), v_field)
         end if
       end if

       v_field => extract_vector_field(state(i), "FiniteElementDivergenceTransposed", stat)
       if(stat == 0) then
         if(recalculate(trim(v_field%option_path))) then
           call calculate_div_t_fe(state(i), v_field)
         end if
       end if
              
       v_field => extract_vector_field(state(i), "PlanetaryVorticity", stat)
       if(stat == 0) then
         if(recalculate(trim(v_field%option_path))) then
           call calculate_planetary_vorticity(state(i), v_field)
         end if
       end if
       
       v_field => extract_vector_field(state(i), "AbsoluteVorticity", stat)
       if(stat == 0) then
         if(recalculate(trim(v_field%option_path))) then
           call calculate_absolute_vorticity(state(i), v_field)
         end if
       end if
       
       s_field => extract_scalar_field(state(i), "PotentialVorticity", stat)
       if(stat == 0) then
         diagnostic = have_option(trim(s_field%option_path)//"/diagnostic")
         if(diagnostic .and. recalculate(trim(s_field%option_path))) then
           call calculate_potential_vorticity(state(i), s_field)
         end if
       end if
       
       s_field => extract_scalar_field(state(i), "RelativePotentialVorticity", stat)
       if(stat == 0) then
         if(recalculate(trim(s_field%option_path))) then
           call calculate_relative_potential_vorticity(state(i), s_field)
         end if
       end if
       ! End of vorticity diagnostics

       ! Start of spontaneous potentials diagnostics
       if(i == 1) then
         s_field => extract_scalar_field(state(i), "ElectricalConductivity", stat)
         if(stat == 0) then
           diagnostic = have_option(trim(s_field%option_path)//"/diagnostic/algorithm::Internal")
           if(diagnostic .and. recalculate(trim(s_field%option_path))) then
             call calculate_formation_conductivity(state(i), i, s_field, stat)
           end if
         end if
       end if
       ! End of spontaneous potentials diagnostics

       ! Start of sediment diagnostics.
       if (have_option("/material_phase[0]/sediment")) then
          call calculate_sediment_flux(state(i))
       end if
       ! End of sediment diagnostics.

       ! Multiphase-related diagnostic fields
       s_field => extract_scalar_field(state(i), "PhaseVolumeFraction", stat)
       if(stat == 0) then
         diagnostic = have_option(trim(s_field%option_path)//"/diagnostic")
         if(diagnostic .and. .not.(aliased(s_field))) then
           if(recalculate(trim(s_field%option_path))) then
            call calculate_diagnostic_phase_volume_fraction(state)
           end if
         end if
       end if

       s_field => extract_scalar_field(state(i), "SumVelocityDivergence", stat)
       if(stat == 0) then
         ! Check that we are running a multiphase simulation
         if(option_count("/material_phase/vector_field::Velocity/prognostic") > 1) then 
            diagnostic = have_option(trim(s_field%option_path)//"/diagnostic")
            if(diagnostic .and. .not.(aliased(s_field))) then
               if(recalculate(trim(s_field%option_path))) then
                  call calculate_sum_velocity_divergence(state, s_field)
               end if
            end if
         else
            FLExit("The SumVelocityDivergence field is only used in multiphase simulations.")
         end if
       end if

       ! end of fields that cannot be called through the generic
       ! calculate_diagnostic_variable interface, i.e. - those that need things
       ! higher than femtools in the build
       
       ! the following fields need to be here in case they are taking the difference with
       ! other diagnostic fields
       s_field => extract_scalar_field(state(i), "ScalarAbsoluteDifference", stat)  
       if(stat == 0) then  
         if(recalculate(trim(s_field%option_path))) then
           call calculate_diagnostic_variable(state(i), "AbsoluteDifference", s_field)  
         end if
       end if
       
       v_field => extract_vector_field(state(i), "VectorAbsoluteDifference", stat)  
       if(stat == 0) then  
         if(recalculate(trim(v_field%option_path))) then
           call calculate_diagnostic_variable(state(i), "AbsoluteDifference", v_field)  
         end if
       end if

       s_field => extract_scalar_field(state(i), "AbsoluteDifference", stat)  
       if(stat == 0) then  
         if(recalculate(trim(s_field%option_path))) then
           call calculate_diagnostic_variable(state(i), "AbsoluteDifference", s_field)  
         end if
       end if
       
       v_field => extract_vector_field(state(i), "AbsoluteDifference", stat)  
       if(stat == 0) then  
         if(recalculate(trim(v_field%option_path))) then
           call calculate_diagnostic_variable(state(i), "AbsoluteDifference", v_field)  
         end if
       end if

    end do
    
    ewrite(1, *) "Exiting calculate_diagnostic_variables"
    
  contains
  
    logical function recalculate(option_path)
      character(len=*) :: option_path
      
      recalculate = ((.not.present_and_true(exclude_nonrecalculated)).or. &
           (.not.do_not_recalculate(option_path)))
    
    end function recalculate

  end subroutine calculate_diagnostic_variables

end module diagnostic_fields_wrapper