Neko  0.9.99
A portable framework for high-order spectral element flow simulations
weak_grad.f90
Go to the documentation of this file.
1 ! Copyright (c) 2024, The Neko Authors
2 ! All rights reserved.
3 !
4 ! Redistribution and use in source and binary forms, with or without
5 ! modification, are permitted provided that the following conditions
6 ! are met:
7 !
8 ! * Redistributions of source code must retain the above copyright
9 ! notice, this list of conditions and the following disclaimer.
10 !
11 ! * Redistributions in binary form must reproduce the above
12 ! copyright notice, this list of conditions and the following
13 ! disclaimer in the documentation and/or other materials provided
14 ! with the distribution.
15 !
16 ! * Neither the name of the authors nor the names of its
17 ! contributors may be used to endorse or promote products derived
18 ! from this software without specific prior written permission.
19 !
20 ! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21 ! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22 ! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23 ! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24 ! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25 ! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26 ! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27 ! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28 ! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29 ! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30 ! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31 ! POSSIBILITY OF SUCH DAMAGE.
32 !
33 !
35 
36 module weak_grad
37  use num_types, only : rp, dp, sp
38  use json_module, only : json_file
41  use field, only : field_t
42  use operators, only : opgrad
43  use case, only : case_t
46  use field_writer, only : field_writer_t
47  implicit none
48  private
49 
52  type, public, extends(simulation_component_t) :: weak_grad_t
54  type(field_t), pointer :: u
56  type(field_t), pointer :: grad_x
58  type(field_t), pointer :: grad_y
60  type(field_t), pointer :: grad_z
62  type(field_writer_t) :: writer
63 
64  contains
66  procedure, pass(this) :: init => weak_grad_init_from_json
68  procedure, pass(this) :: init_from_attributes => &
71  procedure, pass(this) :: free => weak_grad_free
73  procedure, pass(this) :: compute_ => weak_grad_compute
74  end type weak_grad_t
75 
76 contains
77 
79  subroutine weak_grad_init_from_json(this, json, case)
80  class(weak_grad_t), intent(inout) :: this
81  type(json_file), intent(inout) :: json
82  class(case_t), intent(inout), target :: case
83  character(len=:), allocatable :: fieldname
84  character(len=20) :: fields(3)
85 
86  ! Add fields keyword to the json so that the field_writer picks it up.
87  ! Will also add fields to the registry.
88  call json_get(json, "field", fieldname)
89 
90  fields(1) = "weak_grad_" // trim(fieldname) // "_x"
91  fields(2) = "weak_grad_" // trim(fieldname) // "_y"
92  fields(3) = "weak_grad_" // trim(fieldname) // "_z"
93 
94  call json%add("fields", fields)
95 
96  call this%init_base(json, case)
97  call this%writer%init(json, case)
98 
99  call weak_grad_init_from_attributes(this, fieldname)
100  end subroutine weak_grad_init_from_json
101 
103  subroutine weak_grad_init_from_attributes(this, fieldname)
104  class(weak_grad_t), intent(inout) :: this
105  character(len=*) :: fieldname
106 
107  this%u => neko_field_registry%get_field_by_name(trim(fieldname))
108 
109  this%grad_x => neko_field_registry%get_field_by_name(&
110  "weak_grad_" // fieldname // "_x")
111  this%grad_y => neko_field_registry%get_field_by_name(&
112  "weak_grad_" // fieldname // "_y")
113  this%grad_z => neko_field_registry%get_field_by_name(&
114  "weak_grad_" // fieldname // "_z")
115 
116 
117  end subroutine weak_grad_init_from_attributes
118 
120  subroutine weak_grad_free(this)
121  class(weak_grad_t), intent(inout) :: this
122  call this%free_base()
123  call this%writer%free()
124  nullify(this%grad_x)
125  nullify(this%grad_y)
126  nullify(this%grad_z)
127  nullify(this%u)
128  end subroutine weak_grad_free
129 
133  subroutine weak_grad_compute(this, t, tstep)
134  class(weak_grad_t), intent(inout) :: this
135  real(kind=rp), intent(in) :: t
136  integer, intent(in) :: tstep
137 
138  call opgrad(this%grad_x%x, this%grad_y%x, this%grad_z%x, this%u%x,&
139  this%case%fluid%c_Xh)
140  end subroutine weak_grad_compute
141 
142 end module weak_grad
Retrieves a parameter by name or assigns a provided default value. In the latter case also adds the m...
Definition: json_utils.f90:54
Retrieves a parameter by name or throws an error.
Definition: json_utils.f90:45
Defines a simulation case.
Definition: case.f90:34
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.
Defines a field.
Definition: field.f90:34
Implements fld_file_output_t.
Utilities for retrieving parameters from the case files.
Definition: json_utils.f90:34
integer, parameter, public dp
Definition: num_types.f90:9
integer, parameter, public sp
Definition: num_types.f90:8
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
Operators.
Definition: operators.f90:34
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.
Definition: operators.f90:171
Simulation components are objects that encapsulate functionality that can be fit to a particular comp...
Implements the weak_grad_t type.
Definition: weak_grad.f90:36
subroutine weak_grad_init_from_json(this, json, case)
Constructor from json.
Definition: weak_grad.f90:80
subroutine weak_grad_init_from_attributes(this, fieldname)
Actual constructor.
Definition: weak_grad.f90:104
subroutine weak_grad_compute(this, t, tstep)
Compute the weak_grad field.
Definition: weak_grad.f90:134
subroutine weak_grad_free(this)
Destructor.
Definition: weak_grad.f90:121
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 simulation component that computes the weak gradient of a field. Wraps the opgrad operator.
Definition: weak_grad.f90:52