38  use json_module, 
only : json_file
 
   70     generic :: init_from_components => &
 
   71          init_from_controllers, init_from_controllers_properties
 
   73     procedure, pass(this) :: init_from_controllers => &
 
   77     procedure, pass(this) :: init_from_controllers_properties => &
 
 
   92    type(json_file), 
intent(inout) :: json
 
   93    class(
case_t), 
intent(inout), 
target :: case
 
   94    character(len=:), 
allocatable :: field_name
 
   95    character(len=20) :: fields(3)
 
   96    character(len=:), 
allocatable :: computed_field
 
  100    call json_get(json, 
"field", field_name)
 
  102         "weak_gradient" // trim(field_name))
 
  104    fields(1) = computed_field // 
"_x" 
  105    fields(2) = computed_field // 
"_y" 
  106    fields(3) = computed_field // 
"_z" 
  108    call json%add(
"fields", fields)
 
  110    call this%init_base(json, 
case)
 
  111    call this%writer%init(json, 
case)
 
  113    call this%init_common(field_name, computed_field)
 
 
  119    character(len=*) :: field_name
 
  120    character(len=*) :: computed_field
 
  125         computed_field // 
"_x")
 
  127         computed_field // 
"_y")
 
  129         computed_field // 
"_z")
 
 
  146       preprocess_controller, compute_controller, output_controller, &
 
  147       field_name, computed_field, filename, precision)
 
  149    class(
case_t), 
intent(inout), 
target :: case
 
  154    character(len=*) :: field_name
 
  155    character(len=*) :: computed_field
 
  156    character(len=*), 
intent(in), 
optional :: filename
 
  157    integer, 
intent(in), 
optional :: precision
 
  159    character(len=20) :: fields(3)
 
  161    fields(1) = trim(computed_field) // 
"_x" 
  162    fields(2) = trim(computed_field) // 
"_y" 
  163    fields(3) = trim(computed_field) // 
"_z" 
  165    call this%init_base_from_components(
case, order, preprocess_controller, &
 
  167    call this%writer%init_from_components(
case, order, preprocess_controller, &
 
  169    call this%init_common(field_name, computed_field)
 
 
  190       case, order, preprocess_control, preprocess_value, compute_control, &
 
  191       compute_value, output_control, output_value, field_name, &
 
  192       computed_field, filename, precision)
 
  194    class(
case_t), 
intent(inout), 
target :: case
 
  196    character(len=*), 
intent(in) :: preprocess_control
 
  197    real(kind=
rp), 
intent(in) :: preprocess_value
 
  198    character(len=*), 
intent(in) :: compute_control
 
  199    real(kind=
rp), 
intent(in) :: compute_value
 
  200    character(len=*), 
intent(in) :: output_control
 
  201    real(kind=
rp), 
intent(in) :: output_value
 
  202    character(len=*) :: field_name
 
  203    character(len=*) :: computed_field
 
  204    character(len=*), 
intent(in), 
optional :: filename
 
  205    integer, 
intent(in), 
optional :: precision
 
  207    character(len=20) :: fields(3)
 
  209    fields(1) = trim(computed_field) // 
"_x" 
  210    fields(2) = trim(computed_field) // 
"_y" 
  211    fields(3) = trim(computed_field) // 
"_z" 
  213    call this%init_base_from_components(
case, order, preprocess_control, &
 
  214         preprocess_value, compute_control, compute_value, output_control, &
 
  216    call this%writer%init_from_components(
case, order, preprocess_control, &
 
  217         preprocess_value, compute_control, compute_value, output_control, &
 
  218         output_value, fields, filename, precision)
 
  219    call this%init_common(field_name, computed_field)
 
 
  226    call this%free_base()
 
  227    call this%writer%free()
 
  228    nullify(this%gradient_x)
 
  229    nullify(this%gradient_y)
 
  230    nullify(this%gradient_z)
 
 
  240    call opgrad(this%gradient_x%x, this%gradient_y%x, this%gradient_z%x, this%u%x,&
 
  241         this%case%fluid%c_Xh)
 
 
Retrieves a parameter by name or assigns a provided default value. In the latter case also adds the m...
 
Retrieves a parameter by name or throws an error.
 
Defines a simulation case.
 
Defines a registry for storing solution fields.
 
type(field_registry_t), target, public neko_field_registry
Global field registry.
 
Implements the field_writer_t type.
 
Implements fld_file_output_t.
 
Utilities for retrieving parameters from the case files.
 
integer, parameter, public dp
 
integer, parameter, public sp
 
integer, parameter, public rp
Global precision used in computations.
 
subroutine, public opgrad(ux, uy, uz, u, coef, es, ee)
Compute the weak gradient of a scalar field, i.e. the gradient multiplied by the mass matrix.
 
Implements output_controller_t
 
Simulation components are objects that encapsulate functionality that can be fit to a particular comp...
 
subroutine compute_(this, time)
Dummy compute function.
 
Contains the time_based_controller_t type.
 
Module with things related to the simulation time.
 
Implements the weak_gradient_t type.
 
subroutine weak_gradient_init_common(this, field_name, computed_field)
Common part of the constructors.
 
subroutine weak_gradient_free(this)
Destructor.
 
subroutine weak_gradient_init_from_controllers_properties(this, case, order, preprocess_control, preprocess_value, compute_control, compute_value, output_control, output_value, field_name, computed_field, filename, precision)
Constructor from components, passing properties to the time_based_controller` components in the base ...
 
subroutine weak_gradient_init_from_json(this, json, case)
Constructor from json.
 
subroutine weak_gradient_compute(this, time)
Compute the weak_gradient field.
 
subroutine weak_gradient_init_from_controllers(this, case, order, preprocess_controller, compute_controller, output_controller, field_name, computed_field, filename, precision)
Constructor from components, passing controllers.
 
A simulation component that writes a 3d field to a file.
 
A simple output saving a list of fields to a .fld file.
 
Base abstract class for simulation components.
 
A utility type for determining whether an action should be executed based on the current time value....
 
A struct that contains all info about the time, expand as needed.
 
A simulation component that computes the weak gradient of a field. Wraps the opgradient operator.