77
77
! SURFACE_ALBEDO, QTY_SURFACE_ALBEDO
78
78
! OCO2_SIF, QTY_SOLAR_INDUCED_FLUORESCENCE, COMMON_CODE
79
79
! ECOSTRESS_ET, QTY_LATENT_HEAT_FLUX, COMMON_CODE
80
+ ! HARMONIZED_SIF, QTY_SOLAR_INDUCED_FLUORESCENCE
80
81
! END DART PREPROCESS TYPE DEFINITIONS
81
82
82
83
!- ----------------------------------------------------------------------------
83
84
! BEGIN DART PREPROCESS USE OF SPECIAL OBS_DEF MODULE
84
85
! use obs_def_land_mod, only : calculate_albedo, &
85
86
! calculate_biomass, &
86
- ! calculate_fpar
87
+ ! calculate_fpar, &
88
+ ! calculate_sif, &
89
+ ! read_sif_wavelength, &
90
+ ! write_sif_wavelength, &
91
+ ! interactive_sif_wavelength
87
92
! END DART PREPROCESS USE OF SPECIAL OBS_DEF MODULE
88
93
!- ----------------------------------------------------------------------------
89
94
95
100
! call calculate_biomass(state_handle, ens_size, location, expected_obs, istatus)
96
101
! case(MODIS_FPAR)
97
102
! call calculate_fpar(state_handle, ens_size, location, expected_obs, istatus)
103
+ ! case(HARMONIZED_SIF)
104
+ ! call calculate_sif(state_handle, ens_size, location, expected_obs, istatus)
98
105
! END DART PREPROCESS GET_EXPECTED_OBS_FROM_DEF
99
106
!- ----------------------------------------------------------------------------
100
107
104
111
! BIOMASS, &
105
112
! MODIS_FPAR)
106
113
! continue
114
+ ! case(HARMONIZED_SIF)
115
+ ! call read_SIF_wavelength(obs_def%key, key, ifile, fform)
107
116
! END DART PREPROCESS READ_OBS_DEF
108
117
!- ----------------------------------------------------------------------------
109
118
113
122
! BIOMASS, &
114
123
! MODIS_FPAR)
115
124
! continue
125
+ ! case(HARMONIZED_SIF)
126
+ ! call write_SIF_wavelength(key, ifile, fform)
116
127
! END DART PREPROCESS WRITE_OBS_DEF
117
128
!- ----------------------------------------------------------------------------
118
129
122
133
! BIOMASS, &
123
134
! MODIS_FPAR)
124
135
! continue
136
+ ! case(HARMONIZED_SIF)
137
+ ! call interactive_SIF_wavelength(obs_def%key)
125
138
! END DART PREPROCESS INTERACTIVE_OBS_DEF
126
139
!- ----------------------------------------------------------------------------
127
140
@@ -138,10 +151,9 @@ module obs_def_land_mod
138
151
use location_mod, only : location_type, &
139
152
write_location
140
153
141
- use utilities_mod, only : register_module, &
142
- error_handler, &
154
+ use utilities_mod, only : error_handler, &
143
155
E_ERR, E_MSG, &
144
- do_output
156
+ do_output, ascii_file_format
145
157
146
158
use assim_model_mod, only : interpolate
147
159
@@ -157,18 +169,22 @@ module obs_def_land_mod
157
169
QTY_FRACTION_ABSORBED_PAR, &
158
170
QTY_PAR_DIRECT, &
159
171
QTY_PAR_DIFFUSE, &
160
- QTY_ABSORBED_PAR
172
+ QTY_ABSORBED_PAR, &
173
+ QTY_SOLAR_INDUCED_FLUORESCENCE
161
174
162
175
implicit none
163
176
private
164
177
165
178
public :: calculate_albedo, &
166
179
calculate_biomass, &
167
- calculate_fpar
180
+ calculate_fpar, &
181
+ calculate_sif, &
182
+ set_SIF_wavelength, &
183
+ read_SIF_wavelength, &
184
+ write_SIF_wavelength, &
185
+ interactive_SIF_wavelength
168
186
169
- character (len=* ), parameter :: source = ' obs_def_land_mod.f90'
170
- character (len=* ), parameter :: revision = ' '
171
- character (len=* ), parameter :: revdate = ' '
187
+ character (len=* ), parameter :: source = ' obs_def_land_mod.f90'
172
188
173
189
logical :: module_initialized = .false.
174
190
@@ -177,6 +193,11 @@ module obs_def_land_mod
177
193
! This might be useful, but not enough to warrant a namelist ... yet
178
194
logical :: debug = .false.
179
195
196
+ ! Bits and bobs for the solar-induced fluorescence metadata
197
+ integer :: max_num_sif_obs = 200000
198
+ integer :: sifkey = 0
199
+ integer , allocatable :: sif_wavelength(:)
200
+ character (len=* ), parameter :: SIF_STRING = ' lambda'
180
201
181
202
! ===============================================================================
182
203
contains
@@ -189,8 +210,7 @@ subroutine initialize_module()
189
210
190
211
module_initialized = .true.
191
212
192
- ! Log the version of this source file.
193
- call register_module(source, revision, revdate)
213
+ allocate (sif_wavelength(max_num_sif_obs))
194
214
195
215
end subroutine initialize_module
196
216
@@ -219,8 +239,7 @@ subroutine calculate_albedo(state_handle, ens_size, location, obs_val, istatus)
219
239
istatus = 1 ! 0 == success, anything else is a failure
220
240
obs_val = MISSING_R8
221
241
222
- call error_handler(E_ERR,' calculate_albedo' ,' routine untested - stopping.' , &
223
- source, revision, revdate)
242
+ call error_handler(E_ERR,' calculate_albedo' ,' routine untested - stopping.' , source )
224
243
225
244
if ( .not. module_initialized ) call initialize_module()
226
245
@@ -416,6 +435,210 @@ subroutine calculate_fpar(state_handle, ens_size, location, obs_val, istatus)
416
435
417
436
end subroutine calculate_fpar
418
437
438
+
439
+ !- ------------------------------------------------------------------------------
440
+
441
+
442
+ subroutine calculate_sif (state_handle , ens_size , location , obs_val , istatus )
443
+
444
+ type (ensemble_type), intent (in ) :: state_handle
445
+ integer , intent (in ) :: ens_size
446
+ type (location_type), intent (in ) :: location
447
+ real (r8 ), intent (out ) :: obs_val(ens_size)
448
+ integer , intent (out ) :: istatus(ens_size)
449
+
450
+ if ( .not. module_initialized ) call initialize_module()
451
+
452
+ ! If the model state has it directly, this is simple.
453
+ ! If it does not ... nothing else to try at the moment
454
+
455
+ call interpolate(state_handle, ens_size, location, &
456
+ QTY_SOLAR_INDUCED_FLUORESCENCE, obs_val, istatus)
457
+
458
+ end subroutine calculate_sif
459
+
460
+
461
+ !- ------------------------------------------------------------------------------
462
+ ! > stuff the value into the local metadata array
463
+
464
+ function set_SIF_wavelength (lambda ) result(key)
465
+
466
+ integer , intent (in ) :: lambda
467
+ integer :: key
468
+
469
+ if ( .not. module_initialized ) call initialize_module
470
+
471
+ ! update the index into the module array
472
+ sifkey = sifkey + 1
473
+
474
+ ! check that it fits
475
+ if (sifkey > max_num_sif_obs) call double_metadata()
476
+
477
+ sif_wavelength(sifkey) = lambda
478
+ key = sifkey
479
+
480
+ end function set_SIF_wavelength
481
+
482
+
483
+ !- ------------------------------------------------------------------------------
484
+ ! > writes the metadata for SIF observations.
485
+
486
+ subroutine read_SIF_wavelength (key , obsID , ifile , fform )
487
+
488
+ integer , intent (out ) :: key
489
+ integer , intent (in ) :: obsID
490
+ integer , intent (in ) :: ifile
491
+ character (len=* ), intent (in ), optional :: fform
492
+
493
+ character (len=* ), parameter :: routine = ' read_SIF_wavelength'
494
+
495
+ logical :: is_asciifile
496
+ integer :: lambda
497
+ character (len= 6 ) :: header
498
+ integer :: ierr
499
+ character (len= 512 ) :: msgstring
500
+
501
+ if ( .not. module_initialized ) call initialize_module
502
+
503
+ ! create string for error reporting
504
+ write (msgstring,* )' observation # ' ,obsID
505
+
506
+ ! given the index into the local metadata arrays - retrieve
507
+ ! the metadata for this particular observation.
508
+
509
+ is_asciifile = ascii_file_format(fform)
510
+
511
+ if (is_asciifile) then
512
+ read (ifile, * , iostat= ierr) header
513
+ call check_iostat(ierr,routine,' header' ,msgstring)
514
+ read (ifile, * , iostat= ierr) lambda
515
+ call check_iostat(ierr,routine,' lambda' ,msgstring)
516
+ read (ifile, * , iostat= ierr) key
517
+ call check_iostat(ierr,routine,' key' ,msgstring)
518
+ else
519
+ read (ifile , iostat= ierr) header
520
+ call check_iostat(ierr,routine,' header' ,msgstring)
521
+ read (ifile , iostat= ierr) lambda
522
+ call check_iostat(ierr,routine,' lambda' ,msgstring)
523
+ read (ifile , iostat= ierr) key
524
+ call check_iostat(ierr,routine,' key' ,msgstring)
525
+ endif
526
+
527
+ sifkey = sifkey + 1
528
+
529
+ ! check that it fits
530
+ if (sifkey > max_num_sif_obs) call double_metadata()
531
+
532
+ sif_wavelength(sifkey) = lambda
533
+ key = sifkey
534
+
535
+ end subroutine read_SIF_wavelength
536
+
537
+ !- ------------------------------------------------------------------------------
538
+ ! > writes the metadata for SIF observations.
539
+
540
+ subroutine write_SIF_wavelength (key , ifile , fform )
541
+
542
+ integer , intent (in ) :: key
543
+ integer , intent (in ) :: ifile
544
+ character (len=* ), intent (in ), optional :: fform
545
+
546
+ logical :: is_asciifile
547
+ integer :: lambda
548
+
549
+ if ( .not. module_initialized ) call initialize_module
550
+
551
+ ! given the index into the local metadata arrays - retrieve
552
+ ! the metadata for this particular observation.
553
+
554
+ lambda = sif_wavelength(key)
555
+
556
+ is_asciifile = ascii_file_format(fform)
557
+
558
+ if (is_asciifile) then
559
+ write (ifile, * ) trim (SIF_STRING)
560
+ write (ifile, * ) lambda
561
+ write (ifile, * ) key
562
+ else
563
+ write (ifile ) trim (SIF_STRING)
564
+ write (ifile ) lambda
565
+ write (ifile ) key
566
+ endif
567
+
568
+ end subroutine write_SIF_wavelength
569
+
570
+
571
+ !- ------------------------------------------------------------------------------
572
+ ! > interactively queries for metadata required for SIF observations.
573
+
574
+ subroutine interactive_SIF_wavelength (key )
575
+
576
+ integer , intent (out ) :: key
577
+
578
+ character (len=* ), parameter :: routine = ' interactive_SIF_wavelength'
579
+ integer :: lambda
580
+
581
+ if ( .not. module_initialized ) call initialize_module
582
+
583
+ write (* ,* ) ' Input wavelength of SIF (nm)'
584
+ read (* ,* )lambda
585
+
586
+ key = set_SIF_wavelength(lambda)
587
+
588
+ end subroutine interactive_SIF_wavelength
589
+
590
+
591
+ !- ------------------------------------------------------------------------------
592
+ ! > creates enough space for more SIF metadata
593
+
594
+ subroutine double_metadata ()
595
+
596
+ integer , allocatable :: temp_array(:)
597
+ integer :: existing_length
598
+ integer :: new_length
599
+
600
+ existing_length = size (sif_wavelength)
601
+ new_length = 2 * existing_length
602
+
603
+ write (string1,* )' increasing metadata length from ' ,existing_length, &
604
+ ' to ' ,new_length
605
+ call error_handler(E_MSG,' double_metadata' ,string1,source)
606
+
607
+ allocate (temp_array(existing_length))
608
+
609
+ temp_array = sif_wavelength
610
+
611
+ deallocate (sif_wavelength)
612
+ allocate ( sif_wavelength(new_length))
613
+
614
+ sif_wavelength(1 :existing_length) = temp_array
615
+
616
+ deallocate (temp_array)
617
+
618
+ max_num_sif_obs = new_length
619
+
620
+ end subroutine double_metadata
621
+
622
+
623
+ !- ------------------------------------------------------------------------------
624
+ ! > simple error handling routine
625
+
626
+ subroutine check_iostat (istat , routine , context , msgstring )
627
+
628
+ integer , intent (in ) :: istat
629
+ character (len=* ), intent (in ) :: routine
630
+ character (len=* ), intent (in ) :: context
631
+ character (len=* ), intent (in ) :: msgstring
632
+
633
+ if ( istat /= 0 ) then
634
+ write (string1,* )' istat should be 0 but is ' ,istat,' for ' // context
635
+ call error_handler(E_ERR, routine, string1, source, text2= msgstring)
636
+ end if
637
+
638
+ end subroutine check_iostat
639
+
640
+
641
+
419
642
end module obs_def_land_mod
420
643
421
644
! END DART PREPROCESS MODULE CODE
0 commit comments