Neko 0.9.99
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
tree_amg_utils.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!
35 use num_types
36 use utils
37 use math
38 use tree_amg, only : tamg_hierarchy_t
39 use gather_scatter, only : gs_op_add
40 implicit none
41 private
42
44
45contains
46
53 subroutine tamg_sample_matrix_val(val, amg, lvl, i, j)
54 real(kind=rp), intent(inout) :: val
55 type(tamg_hierarchy_t), intent(inout) :: amg
56 integer, intent(in) :: lvl, i, j
57 real(kind=rp), allocatable :: u(:), v(:)
58 integer :: n
59 n = amg%lvl(lvl+1)%fine_lvl_dofs
60 allocate( u(n) )
61 allocate( v(n) )
62 u = 0.0_rp
63 v = 0.0_rp
64
65 u(i) = 1.0_rp
66
68 if (lvl.eq.0) call amg%gs_h%op(u, n, gs_op_add)
69
70 call amg%matvec(v, u, lvl)
71
72 val = v(j)
73 deallocate( u )
74 deallocate( v )
75 end subroutine tamg_sample_matrix_val
76
79 subroutine tamg_print_matrix(amg, lvl)
80 type(tamg_hierarchy_t), intent(inout) :: amg
81 integer, intent(in) :: lvl
82 integer :: n, i, j
83 real(kind=rp) :: val
84 real(kind=rp), allocatable :: u(:), v(:)
85
86 n = amg%lvl(lvl+1)%fine_lvl_dofs
87 allocate( u(n) )
88 allocate( v(n) )
89
90 do i = 1, n
91 do j = 1, n
92 u = 0.0_rp
93 v = 0.0_rp
94
95 u(i) = 1.0_rp
97 if (lvl.eq.0) call amg%gs_h%op(u, n, gs_op_add)
98
99 call amg%matvec(v, u, lvl)
100
101 val = v(j)
102
103 print *, j, i, val
104 end do
105 end do
106
107 deallocate( u )
108 deallocate( v )
109 end subroutine tamg_print_matrix
110
113 subroutine tamg_print_restriction_matrix(amg, lvl)
114 type(tamg_hierarchy_t), intent(inout) :: amg
115 integer, intent(in) :: lvl
116 integer :: n, m, i, j
117 real(kind=rp) :: val
118 real(kind=rp), allocatable :: u(:), v(:)
119
120 n = amg%lvl(lvl+1)%fine_lvl_dofs
121 m = amg%lvl(lvl+2)%fine_lvl_dofs
122 allocate( u(n) )
123 allocate( v(m) )
124
125 do i = 1, n
126 do j = 1, m
127 u = 0.0_rp
128 v = 0.0_rp
129
130 u(i) = 1.0_rp
132 if (lvl.eq.0) call amg%gs_h%op(u, n, gs_op_add)
133
134 call amg%interp_f2c(v, u, lvl+1)
135
136 val = v(j)
137
138 print *, j, i, val
139 end do
140 end do
141
142 deallocate( u )
143 deallocate( v )
144 end subroutine tamg_print_restriction_matrix
145
149 type(tamg_hierarchy_t), intent(inout) :: amg
150 integer, intent(in) :: lvl
151 integer :: n, m, i, j
152 real(kind=rp) :: val
153 real(kind=rp), allocatable :: u(:), v(:)
154
155 n = amg%lvl(lvl+1)%fine_lvl_dofs
156 m = amg%lvl(lvl+2)%fine_lvl_dofs
157 allocate( u(m) )
158 allocate( v(n) )
159
160 do i = 1, m
161 do j = 1, n
162 u = 0.0_rp
163 v = 0.0_rp
164
165 u(i) = 1.0_rp
166 !---!> Get the ghosts to agree on the value of 1.0
167 !---if (lvl.eq.0) call amg%gs_h%op(u, n, GS_OP_ADD)
168
169 call amg%interp_c2f(v, u, lvl+1)
170
171 val = v(j)
172
173 print *, j, i, val
174 end do
175 end do
176
177 deallocate( u )
178 deallocate( v )
179 end subroutine tamg_print_prolongation_matrix
180
181
182end module tree_amg_utils
Gather-scatter.
Definition math.f90:60
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Implements utilities for the TreeAMG hierarchy structure.
subroutine tamg_print_restriction_matrix(amg, lvl)
never call this function also this does not account for duplicates/ghosts in dofs
subroutine tamg_print_matrix(amg, lvl)
never call this function also this does not account for duplicates/ghosts in dofs
subroutine tamg_print_prolongation_matrix(amg, lvl)
never call this function also this does not account for duplicates/ghosts in dofs
subroutine, public tamg_sample_matrix_val(val, amg, lvl, i, j)
Sample the values in a matix (expensive, use with caution)
Implements the base type for TreeAMG hierarchy structure.
Definition tree_amg.f90:34
Utilities.
Definition utils.f90:35
Type for a TreeAMG hierarchy.
Definition tree_amg.f90:83