Octopus
stencil.F90
Go to the documentation of this file.
1!! Copyright (C) 2008 X. Andrade
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
41
42module stencil_oct_m
43 use debug_oct_m
44 use global_oct_m
47
48 implicit none
49
50 private
51 public :: &
52 stencil_t, &
58
60 ! Components are public by default
61 integer :: arms(1:3,1:3) = reshape([0,0,0,0,0,0,0,0,0], (/3, 3/))
62 integer :: narms = 0
63 end type stargeneral_arms_t
64
69
70 type stencil_t
71 ! Components are public by default
72 integer, private :: dim
73 integer :: center
74 integer :: size
75 integer, allocatable :: points(:, :)
76
77 ! The stargeneral arms
78 type(stargeneral_arms_t) :: stargeneral
79 end type stencil_t
80
81contains
82
83 !-------------------------------------------------------
84 subroutine stencil_allocate(this, dim, size)
85 type(stencil_t), intent(inout) :: this
86 integer, intent(in) :: dim
87 integer, intent(in) :: size
88
89 push_sub(stencil_allocate)
90
91 this%dim = dim
92 this%size = size
93
94 safe_allocate(this%points(1:this%dim, 1:size))
95
96 this%points = 0
97
98 pop_sub(stencil_allocate)
99 end subroutine stencil_allocate
100
101 !-------------------------------------------------------
102 subroutine stencil_copy(input, output)
103 type(stencil_t), intent(in) :: input
104 type(stencil_t), intent(out) :: output
105
106 push_sub(stencil_copy)
107
108 call stencil_allocate(output, input%dim, input%size)
109 output%points = input%points
110 output%center = input%center
111
112 output%stargeneral%narms = input%stargeneral%narms
113 output%stargeneral%arms = input%stargeneral%arms
114
115 pop_sub(stencil_copy)
116 end subroutine stencil_copy
117
118
119 !-------------------------------------------------------
120 subroutine stencil_end(this)
121 type(stencil_t), intent(inout) :: this
122
123 push_sub(stencil_end)
124
125 safe_deallocate_a(this%points)
126
127 pop_sub(stencil_end)
128 end subroutine stencil_end
129
130
131 !-------------------------------------------------------
132 subroutine stencil_union(st1, st2, stu)
133 type(stencil_t), intent(inout) :: st1
134 type(stencil_t), intent(inout) :: st2
135 type(stencil_t), intent(inout) :: stu
136
137 integer :: ii, jj, nstu
138 logical :: not_in_st1
139
140 push_sub(stencil_union)
141
142 assert(st1%dim == st2%dim)
143
144 call stencil_allocate(stu, st1%dim, st1%size + st2%size)
145
146 ! copy the first stencil
147 do ii = 1, st1%size
148 stu%points(:, ii) = st1%points(:, ii)
149 end do
150
151 nstu = st1%size
153 do ii = 1, st2%size
155 not_in_st1 = .true.
156
157 ! check whether that point was already in the stencil
158 do jj = 1, st1%size
159 if (all(st1%points(:, jj) == st2%points(:, ii))) then
160 not_in_st1 = .false.
161 exit
162 end if
163 end do
164
165 if (not_in_st1) then !add it
166 nstu = nstu + 1
167 stu%points(:, nstu) = st2%points(:, ii)
168 end if
169
170 end do
172 stu%size = nstu
173
174 call stencil_init_center(stu)
175
176 pop_sub(stencil_union)
177 end subroutine stencil_union
178
179
180 !-------------------------------------------------------
181 subroutine stencil_init_center(this)
182 type(stencil_t), intent(inout) :: this
183
184 integer :: ii
185
186 push_sub(stencil_init_center)
187
188 this%center = -1
189
190 do ii = 1, this%size
191 if (all(this%points(:, ii) == 0)) this%center = ii
192 end do
193
194 pop_sub(stencil_init_center)
195 end subroutine stencil_init_center
196
197end module stencil_oct_m
198
199!! Local Variables:
200!! mode: f90
201!! coding: utf-8
202!! End:
This module defines stencils used in Octopus.
Definition: stencil.F90:135
subroutine, public stencil_end(this)
Definition: stencil.F90:214
subroutine, public stencil_union(st1, st2, stu)
Definition: stencil.F90:226
subroutine, public stencil_copy(input, output)
Definition: stencil.F90:196
subroutine, public stencil_allocate(this, dim, size)
Definition: stencil.F90:178
subroutine, public stencil_init_center(this)
Definition: stencil.F90:275
The class representing the stencil, which is used for non-local mesh operations.
Definition: stencil.F90:163
int true(void)