-
Notifications
You must be signed in to change notification settings - Fork 115
Adding new diagnostic fields
The name of each field for each material_phase in the options tree must be unique. Hence there can only be one field named "Gradient" in a single material_phase. The concept of a diagnostic algorithm is designed to solve this issue - multiple fields, each with their own name, can share the same algorithm. This, for example, allows the gradient of multiple fields to be calculated for a single material_phase.
Diagnostic algorithms are defined in two locations:
- The options tree, defining options for the algorithm
- The Fortran code implementing the algorithm
Each diagnostic field contains a named algorithm element, defining the type of algorithm and any other algorithm specific options. The possible algorithms are defined in tools/diagnostic_algorithms.rnc.
A new algorithm must, as a minimum, define:
some_algorithm =
(
## Documentation
element algorithm {
attribute name { "[algorithm name]" },
attribute material_phase_support { "[single or multiple]" }
}
)
The algorithm can then be made available for use by either adding it to the lists of generic diagnostic algorithms:
[scalar, vector of tensor]_diagnostic_algorithms |= some_algorithm
or by defining a new diagnostic field in the options tree using that algorithm. The latter case can be used for single instance algorithms (e.g. Coriolis).
Any fields on which the diagnostic algorithm depends (not necessarily themselves diagnostic fields) can be added via a "depends" attribute:
some_algorithm =
(
## Documentation
element algorithm {
attribute name { "[algorithm name]" },
attribute material_phase_support { "[single or multiple]" },
attribute depends { "[list of fields]" }
}
)
The list of fields must be comma delimited. For multiple material_phase dependencies these have the form: "[material_phase name]::[field name]".
In addition, for diagnostics that depend upon a source field, the source_field definitions can be used. For single source fields these take the form:
some_algorithm =
(
## Documentation
element algorithm {
attribute name { "[algorithm name]" },
attribute material_phase_support { "[single or multiple]" },
[scalar, vector, tensor or component]_source_field,
}
)
"component" source field indicates that the diagnostic algorithm expects a scalar field source, but can use a component of a vector field for that scalar field. For binary operator diagnostics these take the form:
some_algorithm =
(
## Documentation
element algorithm {
attribute name { "[algorithm name]" },
attribute material_phase_support { "[single or multiple]" },
[scalar, vector, tensor or component]_source_field_1,
[scalar, vector, tensor or component]_source_field_2,
}
)
Source fields added in this way are automatically treated as dependencies. The complete list of paths treated as dependency lists is: "depends", "source_field_name", "source_field_1_name" and "source_field_2_name", either as attributes for the "diagnostic" element or "algorithm" element of diagnostic fields.
Diagnostic field dependency resolution is implemented. Circular dependencies are checked for at simulation start.
Each diagnostic algorithm "some_algorithm" must be implemented as a (public moduled) subroutine "calculate_some_algorithm", contained in the diagnostics/ directory. This one subroutine is all that is required - the algorithm is automatically wired into the algorithm calculation and dependency checking routines.
These subroutines can only accept as arguments:
states | All system states (material_phase s) |
---|---|
state | The state containing the diagnostic field being calculated |
state_index | The index of state in states |
s_field / v_field / t_field | The diagnostic field being calculated (scalar, vector or tensor respectively) |
current_time | The current simulation time |
dt | The current simulation time step |
Obviously the entire options tree is also available to diagnostic algorithm implementations via libspud.
A module exists to facilitate the extraction of single source fields from the system state. This can be used via:
use diagnostic_source_fields
...
type([scalar, vector or tensor]_field), pointer :: source
...
source => [scalar, vector or tensor]_source_field(state, [diagnostic field to be calculated])
or:
source => [scalar, vector or tensor]_source_field(states, state_index, [diagnostic field to be calculated])
The scalar_field_field interface takes an optional "allocated" output argument for component_source_field support (see Dependencies above). If this argument is returned as .true., it indicates that the scalar source field is a vector field component, and the returned pointer has been allocated. In this case, the pointer "source" must be explicitly deallocated after use:
logical :: allocated
type([scalar, vector or tensor]_field), pointer :: source
...
source => scalar_source_field(state, state_index, [diagnostic field to be calculated], allocated = allocated)
...
if(allocated) deallocate(source)
scalar_source_field, vector_source_field and tensor_source_field all take an optional integer "index" argument, taking the value "1" or "2", for source fields of binary operator diagnostics.
The options in tools/diagnostic_algorithms.rnc:
scalar_copy_algorithm =
(
## Copies the source field. This is intended for testing purposes
## <b>only</b>.
element algorithm {
attribute name { "scalar_copy" },
attribute material_phase_support { "single" },
scalar_source_field
}
)
This is a generic diagnostic, in that we want to allow multiple fields using this algorithm in a single state. Hence we add it to the list of generic algorithms towards the bottom of tools/diagnostic_algorithms.rnc:
scalar_diagnostic_algorithms |= scalar_copy_algorithm
The implementation code in diagnostics/Field_Copies.F90:
#include "fdebug.h"
module field_copies_diagnostics
use diagnostic_source_fields
use fields
use fldebug
use state_module
implicit none
! ...
contains
subroutine calculate_scalar_copy(state, s_field)
type(state_type), intent(in) :: state
type(scalar_field), intent(inout) :: s_field
type(scalar_field), pointer :: source_field
source_field => scalar_source_field(state, s_field)
call remap_field(source_field, s_field)
end subroutine calculate_scalar_copy
! ...
end module field_copies_diagnostics
The options in tools/diagnostic_algorithms.rnc:
vector_difference_algorithm =
(
## Difference between two vector fields, field 1 - field 2
element algorithm {
attribute name { "vector_difference" },
attribute material_phase_support { "single" },
vector_source_field_1,
vector_source_field_2
}
)
This is a generic diagnostic, in that we want to allow multiple fields using this algorithm in a single state. Hence we add it to the list of generic algorithms towards the bottom of tools/diagnostic_algorithms.rnc:
vector_diagnostic_algorithms |= vector_difference_algorithm
The implementation code in diagnostics/Binary_Operators.F90:
#include "fdebug.h"
module binary_operators
use diagnostic_source_fields
use fields
use fldebug
use state_module
implicit none
! ...
contains
subroutine calculate_vector_difference(state, v_field)
type(state_type), intent(in) :: state
type(vector_field), intent(inout) :: v_field
type(vector_field), pointer :: source_field_1, source_field_2
source_field_1 => vector_source_field(state, v_field, index = 1)
source_field_2 => vector_source_field(state, v_field, index = 2)
call set(v_field, source_field_1)
call addto(v_field, source_field_2, scale = -1.0)
end subroutine calculate_vector_difference
! ...
end module binary_operators
All diagnostics calculated in the Fluidity internals must be marked as internal diagnostics, e.g.:
element tensor_field {
attribute rank { "2" },
attribute name { "GLSEddyViscosityKM" },
(
element diagnostic {
internal_algorithm,
...
To add new diagnostic fields using the new Fluidity options system:
- If the diagnostic field depends on assemble code, the routine should be added into assemble/Diagnostic_fields_wrapper.F90 and a call added directly into calculate_diagnostics (and steps 2 and 3 should be skipped). Otherwise, write a routine to calculate the new diagnostic field, preferably in femtools/Diagnostic_fields.F90. This should be of the form:
subroutine calculate_''NAME''(state, d_field, stat)
!!< Calculates the diagnostic field ''NAME''
type(state_type), [dimension(:),] intent(in) :: state
type(''{scalar | vector | tensor}''_field), intent(inout) :: d_field
integer, intent(out), optional :: stat
''!Routine''
end subroutine calculate_''NAME''
The use of an optional stat variable is not essential, but is recommended.
- In one of:
- calculate_diagnostic_scalar_variable_single_state
- calculate_diagnostic_scalar_variable_multiple_states
- calculate_diagnostic_vector_variable_single_state
- calculate_diagnostic_vector_variable_multiple_states
- calculate_diagnostic_tensor_variable_single_state
- calculate_diagnostic_tensor_variable_multiple_states (whichever is appropriate) add a new case corresponding to the new diagnostic field:
case("''NAME''")
call calculate_''NAME''(state, d_field, stat)
- Add a call to calculate_diagnostic_variable (in assemble/Diagnostic_fields_wrapper.F90), passing the state, diagnostic field name and output diagnostic field:
''{s | v | t}''_field => extract_''{scalar | vector | tensor}''_field(state[(i)], "''NAME''", stat)
if(stat == 0) then
call calculate_diagnostic_variable(state[(i)], "''NAME''", &
& ''{s | v | t}''_field)
end if
- Add the new diagnostic field to the options schema
In order for a diagnostic to be calculable off-line by fldiagnostics it must be entered into the wrappers in step 2. This means that diagnostics that depend on asemble code cannot currently be calculated off-line. Also, some future options for diagnostic fields may only be available for diagnostic fields calculated via calculate_diagnostics.