Octopus
grid.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
2!! Copyright (C) 2021 S. Ohlmann
3!!
4!! This program is free software; you can redistribute it and/or modify
5!! it under the terms of the GNU General Public License as published by
6!! the Free Software Foundation; either version 2, or (at your option)
7!! any later version.
8!!
9!! This program is distributed in the hope that it will be useful,
10!! but WITHOUT ANY WARRANTY; without even the implied warranty of
11!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12!! GNU General Public License for more details.
13!!
14!! You should have received a copy of the GNU General Public License
15!! along with this program; if not, write to the Free Software
16!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
17!! 02110-1301, USA.
18!!
19
20#include "global.h"
21
23
24module grid_oct_m
25 use accel_oct_m
27 use box_oct_m
32 use cube_oct_m
37 use global_oct_m
39 use mesh_oct_m
42 use mpi_oct_m
46 use parser_oct_m
48 use space_oct_m
53 use unit_oct_m
56
57 implicit none
58
59 private
60 public :: &
61 grid_t, &
65 grid_end, &
74
76 type, extends(mesh_t) :: grid_t
77 ! Components are public by default
78 type(derivatives_t) :: der
79 type(stencil_t) :: stencil
80 type(symmetries_t) :: symm
81 type(symmetrizer_t) :: symmetrizer
82 end type grid_t
83
84 integer, parameter :: &
85 CURV_AFFINE = 1, &
86 curv_gygi = 2, &
87 curv_briggs = 3, &
88 curv_modine = 4
89
90contains
91
92 !-------------------------------------------------------------------
101 subroutine grid_init_stage_1(gr, namespace, space, grp, symm, latt, n_sites, site_position)
102 type(grid_t), intent(inout) :: gr
103 type(namespace_t), intent(in) :: namespace
104 class(space_t), intent(in) :: space
105 type(mpi_grp_t), intent(in) :: grp
106 type(symmetries_t), optional, intent(in) :: symm
107 type(lattice_vectors_t), optional, intent(in) :: latt
108 integer, optional, intent(in) :: n_sites
109 real(real64), optional, intent(in) :: site_position(:,:)
110
111 type(stencil_t) :: cube
112 integer :: enlarge(space%dim)
113 type(block_t) :: blk
114 integer :: idir
115 real(real64) :: grid_spacing(space%dim)
116
117 push_sub(grid_init_stage_1)
118
119 gr%box => box_factory_create(namespace, space, latt=latt, n_sites=n_sites, site_position=site_position)
120
121 if (present(symm)) then
122 gr%symm = symm
123 end if
124
125 !%Variable Spacing
126 !%Type float
127 !%Section Mesh
128 !%Description
129 !% The spacing between the points in the mesh. This controls the
130 !% quality of the discretization: smaller spacing gives more
131 !% precise results but increased computational cost.
132 !%
133 !% When using curvilinear coordinates, this is a canonical spacing
134 !% that will be changed locally by the transformation. In periodic
135 !% directions, your spacing may be slightly different than what
136 !% you request here, since the box size must be an integer
137 !% multiple of the spacing.
138 !%
139 !% The default value is defined by the image resolution if <tt>BoxShape =
140 !% box_image</tt>. Othewise here is no default otherwise.
141 !%
142 !% It is possible to have a different spacing in each one of the Cartesian directions
143 !% if we define <tt>Spacing</tt> as block of the form
144 !%
145 !% <tt>%Spacing
146 !% <br>&nbsp;&nbsp;spacing_x | spacing_y | spacing_z
147 !% <br>%</tt>
148 !%End
149
150 if(parse_block(namespace, 'Spacing', blk) == 0) then
151 if(parse_block_cols(blk,0) < space%dim) call messages_input_error(namespace, 'Spacing')
152 do idir = 1, space%dim
153 call parse_block_float(blk, 0, idir - 1, grid_spacing(idir), units_inp%length)
154 end do
155 call parse_block_end(blk)
156 else
157 call parse_variable(namespace, 'Spacing', -m_one, grid_spacing(1), units_inp%length)
158 grid_spacing = grid_spacing(1)
159 end if
160
161 if (associated(gr%box)) then
162 select type (box => gr%box)
163 type is (box_image_t)
164 do idir = 1, space%dim
165 ! default grid_spacing is determined from the pixel size such that one grid point = one pixel.
166 if (grid_spacing(idir) < m_zero) then
167 grid_spacing(idir) = box%pixel_size(idir)
168 end if
169 end do
170 end select
171 end if
172
173 if (any(grid_spacing(1:space%dim) < m_epsilon)) then
174 message(1) = " Your input for 'Spacing' is negative or zero."
175 call messages_fatal(1, namespace=namespace)
176 end if
177
178 !%Variable PeriodicBoundaryMask
179 !%Type block
180 !%Section Mesh
181 !%Description
182 !% (Experimental) Defines a mask for which periodic boundaries are replaced by zero boundary conditions.
183 !%End
184 if (parse_block(namespace, 'PeriodicBoundaryMask', blk) < 0) then
185 gr%masked_periodic_boundaries = .false.
186 else
187 gr%masked_periodic_boundaries = .true.
188 call parse_block_string(blk, 0, 0, gr%periodic_boundary_mask)
189 call messages_experimental('PeriodicBoundaryMask')
190 call parse_block_end(blk)
191 end if
192
193 ! Initialize coordinate system
194 call initialize_coordinate_system(gr, namespace, space, grid_spacing, latt, n_sites, site_position)
195
196 ! initialize derivatives
197 call nl_operator_global_init(namespace)
198 call derivatives_init(gr%der, namespace, space, gr%coord_system)
199 ! the stencil used to generate the grid is a union of a cube (for
200 ! multigrid) and the Laplacian
201 call stencil_cube_get_lapl(cube, space%dim, order = 2)
202 call stencil_union(cube, gr%der%lapl%stencil, gr%stencil)
203 call stencil_end(cube)
204
205 enlarge = 2
206 enlarge = max(enlarge, gr%der%n_ghost)
207
208 call mesh_init_stage_1(gr, namespace, space, gr%box, gr%coord_system, grid_spacing, enlarge)
209 call mesh_init_stage_2(gr, namespace, space, gr%box, gr%stencil, grp)
210
211 pop_sub(grid_init_stage_1)
212 end subroutine grid_init_stage_1
213
214 !--------------------------------------------------------------------
216 subroutine grid_init_from_grid_stage_1(gr, original_gr, namespace, space, grp, &
217 box, spacing_prefactor, latt, n_sites, site_position)
218 type(grid_t), intent(inout) :: gr
219 type(grid_t), intent(in) :: original_gr
220 type(namespace_t), intent(in) :: namespace
221 class(space_t), intent(in) :: space
222 type(mpi_grp_t), intent(in) :: grp
223 class(box_t), target, optional, intent(in) :: box
224 real(real64), optional, intent(in) :: spacing_prefactor(:)
225 type(lattice_vectors_t), optional, intent(in) :: latt
226 integer, optional, intent(in) :: n_sites
227 real(real64), optional, intent(in) :: site_position(:,:)
228
229 type(stencil_t) :: cube
230 integer :: enlarge(space%dim)
231
233
234 ! First of all, create and associate the box that contains the grid
235 if (.not. present(box)) then
236 gr%box => box_factory_create(namespace, space)
237 else
238 gr%box => box
239 end if
240
241 !!! Copy all data !!!
242 ! Symmetries
243 gr%symm = original_gr%symm
244 ! Spacing
245 if (present(spacing_prefactor)) then
246 gr%spacing = spacing_prefactor * original_gr%spacing
247 if (any(abs(spacing_prefactor - m_one) > m_epsilon)) then
248 write(message(1), '(a)') "The two grids have different spacings, it will not be possible to establish a map between them"
249 call messages_warning(1)
250 end if
251 else
252 gr%spacing = original_gr%spacing
253 end if
254 ! Periodic boundaries
255 gr%masked_periodic_boundaries = original_gr%masked_periodic_boundaries
256 if (gr%masked_periodic_boundaries) gr%periodic_boundary_mask = original_gr%periodic_boundary_mask
257 ! Coordinate System
258 call initialize_coordinate_system(gr, namespace, space, gr%spacing, latt, n_sites, site_position)
259
260 ! initialize derivatives
261 call nl_operator_global_init(namespace)
262 call derivatives_init(gr%der, namespace, space, gr%coord_system)
263 ! the stencil used to generate the grid is a union of a cube (for
264 ! multigrid) and the Laplacian
265 call stencil_cube_get_lapl(cube, space%dim, order = 2)
266 call stencil_union(cube, gr%der%lapl%stencil, gr%stencil)
267 call stencil_end(cube)
268
269 enlarge = 2
270 enlarge = max(enlarge, gr%der%n_ghost)
271
272 call mesh_init_stage_1(gr, namespace, space, gr%box, gr%coord_system, gr%spacing, enlarge)
273 call mesh_init_stage_2(gr, namespace, space, gr%box, gr%stencil, grp)
274
276 end subroutine grid_init_from_grid_stage_1
277
278
279 !--------------------------------------------------------------------
283 subroutine initialize_coordinate_system(gr, namespace, space, grid_spacing, latt, n_sites, site_position)
284 type(grid_t), intent(inout) :: gr
285 type(namespace_t), intent(in) :: namespace
286 class(space_t), intent(in) :: space
287 real(real64), intent(in) :: grid_spacing(:)
288 type(lattice_vectors_t), optional, intent(in) :: latt
289 integer, optional, intent(in) :: n_sites
290 real(real64), optional, intent(in) :: site_position(:,:)
291
292 integer :: cv_method
293
295
296 !%Variable CurvMethod
297 !%Type integer
298 !%Default curv_uniform
299 !%Section Mesh::Curvilinear
300 !%Description
301 !% The relevant functions in octopus are represented on a mesh in real space.
302 !% This mesh may be an evenly spaced regular rectangular grid (standard mode),
303 !% or else an adaptive or curvilinear grid. We have implemented
304 !% three kinds of adaptive meshes, although only one is currently working,
305 !% the one invented by F. Gygi (<tt>curv_gygi</tt>). The code will stop if any of
306 !% the other two is invoked. All are experimental with domain parallelization.
307 !%Option curv_affine 1
308 !% Regular, uniform rectangular grid.
309 !%Option curv_gygi 2
310 !% The deformation of the grid is done according to the scheme described by
311 !% F. Gygi [F. Gygi and G. Galli, <i>Phys. Rev. B</i> <b>52</b>, R2229 (1995)].
312 !%Option curv_briggs 3
313 !% The deformation of the grid is done according to the scheme described by
314 !% Briggs [E.L. Briggs, D.J. Sullivan, and J. Bernholc, <i>Phys. Rev. B</i> <b>54</b> 14362 (1996)]
315 !%Option curv_modine 4
316 !% The deformation of the grid is done according to the scheme described by
317 !% Modine [N.A. Modine, G. Zumbach and E. Kaxiras, <i>Phys. Rev. B</i> <b>55</b>, 10289 (1997)]
318 !%End
319 call parse_variable(namespace, 'CurvMethod', curv_affine, cv_method)
320 if (.not. varinfo_valid_option('CurvMethod', cv_method)) call messages_input_error(namespace, 'CurvMethod')
321 if (cv_method /= curv_affine .and. accel_is_enabled()) then
322 call messages_write('Curvilinear coordinates on GPUs is not implemented')
323 call messages_fatal(namespace=namespace)
324 end if
325 if (present(latt)) then
326 if (cv_method /= curv_affine .and. latt%nonorthogonal) then
327 call messages_input_error(namespace, 'CurvMethod', 'Curvilinear coordinates with non-orthogonal cells are not implemented')
328 end if
329 end if
330 call messages_print_var_option("CurvMethod", cv_method, namespace=namespace)
331
332 ! FIXME: The other two methods are apparently not working
333 if (cv_method > curv_gygi) then
334 call messages_experimental('Selected curvilinear coordinates method')
335 end if
336
337 select case (cv_method)
338 case (curv_briggs)
339 gr%coord_system => curv_briggs_t(namespace, space%dim, &
340 gr%box%bounding_box_l(1:space%dim), grid_spacing(1:space%dim))
341 case (curv_gygi)
342 if (present(site_position) .and. present(n_sites)) then
343 gr%coord_system => curv_gygi_t(namespace, space%dim, n_sites, site_position)
344 else
345 message(1) = "Option CurvMethod = curv_gygi is not currently implemented without ions present."
346 call messages_fatal(1, namespace=namespace)
347 end if
348 case (curv_modine)
349 if (present(site_position) .and. present(n_sites)) then
350 gr%coord_system => curv_modine_t(namespace, space%dim, n_sites, site_position, &
351 gr%box%bounding_box_l(1:space%dim), grid_spacing(1:space%dim))
352 else
353 message(1) = "Option CurvMethod = curv_modine is not currently implemented without ions present."
354 call messages_fatal(1, namespace=namespace)
355 end if
356 case (curv_affine)
357 if (present(latt)) then
358 if (latt%nonorthogonal) then
359 gr%coord_system => affine_coordinates_t(namespace, space%dim, latt%rlattice_primitive)
360 else
361 gr%coord_system => cartesian_t(namespace, space%dim)
362 end if
363 else
364 gr%coord_system => cartesian_t(namespace, space%dim)
365 end if
366 end select
367
369 end subroutine initialize_coordinate_system
370
371
372 !-------------------------------------------------------------------
380 subroutine grid_init_stage_2(gr, namespace, space, mc, qvector)
381 type(grid_t), target, intent(inout) :: gr
382 type(namespace_t), intent(in) :: namespace
383 class(space_t), intent(in) :: space
384 type(multicomm_t), intent(in) :: mc
385 real(real64), optional, intent(in) :: qvector(:)
386
387 push_sub(grid_init_stage_2)
388
389 call mesh_init_stage_3(gr, namespace, space, gr%stencil, mc)
390
391 call derivatives_build(gr%der, namespace, space, gr, qvector)
392
393 ! print info concerning the grid
394 call grid_write_info(gr, namespace=namespace)
395
396 if(space%dim == 3) then
397 call symmetrizer_init(gr%symmetrizer, gr, gr%symm)
398 end if
399
400 pop_sub(grid_init_stage_2)
401 end subroutine grid_init_stage_2
402
403
404 !-------------------------------------------------------------------
407 subroutine grid_end(gr)
408 type(grid_t), intent(inout) :: gr
409
410 class(coordinate_system_t), pointer :: coord_system
411 class(box_t), pointer :: box
412
413 push_sub(grid_end)
414
416
417 call derivatives_end(gr%der)
418
419 ! We need to take a pointer here, otherwise we run into a gfortran bug.
420 coord_system => gr%coord_system
421 safe_deallocate_p(coord_system)
422 box => gr%box
423 safe_deallocate_p(box)
424
425 call mesh_end(gr)
426
427 call stencil_end(gr%stencil)
428
429 call symmetrizer_end(gr%symmetrizer)
430
431 pop_sub(grid_end)
432 end subroutine grid_end
433
434
435 !-------------------------------------------------------------------
436 subroutine grid_write_info(gr, iunit, namespace)
437 type(grid_t), intent(in) :: gr
438 integer, optional, intent(in) :: iunit
439 type(namespace_t), optional, intent(in) :: namespace
440
441 push_sub(grid_write_info)
442
443 call messages_print_with_emphasis(msg="Grid", iunit=iunit, namespace=namespace)
444
445 message(1) = 'Simulation Box:'
446 call messages_info(1, iunit=iunit, namespace=namespace)
447 call gr%box%write_info(iunit, namespace)
448
449 message(1) = "Main mesh:"
450 call messages_info(1, iunit=iunit, namespace=namespace)
451 call mesh_write_info(gr, iunit, namespace)
452
453 if (gr%use_curvilinear) then
454 call gr%coord_system%write_info(iunit, namespace)
455 end if
456
457 call messages_print_with_emphasis(iunit=iunit, namespace=namespace)
458
459 pop_sub(grid_write_info)
460 end subroutine grid_write_info
461
462 !-------------------------------------------------------------------
466 subroutine grid_lattice_vectors_update(gr, space, namespace, mc, new_latt)
467 type(grid_t), intent(inout) :: gr
468 type(space_t), intent(in) :: space
469 type(namespace_t), intent(in) :: namespace
470 type(multicomm_t), intent(in) :: mc
471 type(lattice_vectors_t), intent(in) :: new_latt
472
473 real(real64) :: new_box_bounds(2, space%dim)
474
476
477 call new_latt%write_info(namespace)
478
479 ! Regenerate the coordinate system
480 select type (coord_system=>gr%coord_system)
481 type is (affine_coordinates_t)
482 if (new_latt%nonorthogonal) then
483 deallocate(gr%coord_system)
484 gr%coord_system => affine_coordinates_t(namespace, space%dim, new_latt%rlattice_primitive)
485 end if
486 end select
487
488 ! Regenerate the grid keeping the same number of points
489 select type (coord_system=>gr%coord_system)
490 class is (affine_coordinates_t)
491 new_box_bounds = gr%box%bounds(coord_system%basis)
492 gr%spacing = (new_box_bounds(2,:)-new_box_bounds(1,:))/gr%idx%ll(:)
493 write(message(1), '(a)') trim(gr%box%short_info(unit_angstrom))
494 call messages_info(1, namespace=namespace)
495 call messages_new_line()
496 end select
497
498 safe_deallocate_a(gr%x)
499 safe_deallocate_a(gr%x_t)
500 safe_deallocate_a(gr%chi)
501 safe_deallocate_a(gr%vol_pp)
502 safe_deallocate_a(gr%jacobian_inverse)
503 call mesh_init_stage_3(gr, namespace, space, gr%stencil, mc, regenerate=.true.)
504
505 call derivates_set_coordinates_system(gr%der, gr%coord_system)
506 call derivatives_build(gr%der, namespace, space, gr, regenerate=.true., verbose=.false.)
507
509 end subroutine grid_lattice_vectors_update
510
511#include "undef.F90"
512#include "real.F90"
513#include "grid_inc.F90"
514
515#include "undef.F90"
516#include "complex.F90"
517#include "grid_inc.F90"
518
519end module grid_oct_m
520
521!! Local Variables:
522!! mode: f90
523!! coding: utf-8
524!! End:
pure logical function, public accel_is_enabled()
Definition: accel.F90:403
Module, implementing a factory for boxes.
class(box_t) function, pointer, public box_factory_create(namespace, space, latt, n_sites, site_position)
initialize a box of any type
This module implements the curvilinear coordinates given in E.L. Briggs, D.J. Sullivan,...
This module implements the curvilinear coordinates given in F. Gygi and G. Galli, PRB 52 R2229 (1996)...
Definition: curv_gygi.F90:120
This module implements the curvilinear coordinates given in N. A. Modine, G. Zumbach,...
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
subroutine, public derivatives_build(der, namespace, space, mesh, qvector, regenerate, verbose)
build the derivatives object:
subroutine, public derivatives_init(der, namespace, space, coord_system, order)
subroutine, public derivatives_end(der)
subroutine, public derivates_set_coordinates_system(der, coord_system)
real(real64), parameter, public m_zero
Definition: global.F90:200
real(real64), parameter, public m_epsilon
Definition: global.F90:216
real(real64), parameter, public m_one
Definition: global.F90:201
This module implements the underlying real-space grid.
Definition: grid.F90:119
subroutine initialize_coordinate_system(gr, namespace, space, grid_spacing, latt, n_sites, site_position)
this subroutine initializes the coordinate system
Definition: grid.F90:379
integer, parameter curv_gygi
Definition: grid.F90:179
integer, parameter curv_briggs
Definition: grid.F90:179
subroutine, public zgrid_symmetrize_scalar_field(gr, field, suppress_warning)
Definition: grid.F90:814
subroutine, public grid_init_from_grid_stage_1(gr, original_gr, namespace, space, grp, box, spacing_prefactor, latt, n_sites, site_position)
this subroutine allows to create a grid from an existing grid
Definition: grid.F90:313
subroutine, public grid_init_stage_1(gr, namespace, space, grp, symm, latt, n_sites, site_position)
First stage of the grid initialization.
Definition: grid.F90:197
subroutine, public zgrid_symmetrize_single(gr, iop, field, symm_field)
Definition: grid.F90:868
subroutine, public zgrid_symmetrize_vector_field(gr, field, suppress_warning)
Definition: grid.F90:840
subroutine, public grid_init_stage_2(gr, namespace, space, mc, qvector)
Second stage of the grid initialization.
Definition: grid.F90:476
subroutine, public dgrid_symmetrize_vector_field(gr, field, suppress_warning)
Definition: grid.F90:701
subroutine, public grid_write_info(gr, iunit, namespace)
Definition: grid.F90:532
subroutine, public grid_lattice_vectors_update(gr, space, namespace, mc, new_latt)
Regenerate the grid information after update of the lattice vectors.
Definition: grid.F90:562
subroutine, public dgrid_symmetrize_single(gr, iop, field, symm_field)
Definition: grid.F90:729
subroutine, public grid_end(gr)
finalize a grid object
Definition: grid.F90:503
subroutine, public dgrid_symmetrize_scalar_field(gr, field, suppress_warning)
Definition: grid.F90:675
integer, parameter curv_modine
Definition: grid.F90:179
This module contains subroutines, related to the initialization of the mesh.
Definition: mesh_init.F90:119
subroutine, public mesh_init_stage_3(mesh, namespace, space, stencil, mc, parent, regenerate)
When running parallel in domains, stencil and np_stencil are needed to compute the ghost points....
Definition: mesh_init.F90:531
subroutine, public mesh_init_stage_1(mesh, namespace, space, box, coord_system, spacing, enlarge)
First stage mesh initialization.
Definition: mesh_init.F90:172
subroutine, public mesh_init_stage_2(mesh, namespace, space, box, stencil, grp, regenerate)
This subroutine creates the global array of spatial indices and the inverse mapping.
Definition: mesh_init.F90:292
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
subroutine, public mesh_write_info(this, iunit, namespace)
Definition: mesh.F90:310
recursive subroutine, public mesh_end(this)
Definition: mesh.F90:686
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
Definition: messages.F90:898
character(len=512), private msg
Definition: messages.F90:167
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:525
subroutine, public messages_new_line()
Definition: messages.F90:1089
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:162
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:410
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:691
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1040
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:594
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:147
This module defines non-local operators.
subroutine, public nl_operator_global_init(namespace)
initialize global settings for non-local operators
subroutine, public nl_operator_global_end()
subroutine, public parse_block_string(blk, l, c, res, convert_to_c)
Definition: parser.F90:818
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:623
This module defines routines, generating operators for a cubic stencil.
subroutine, public stencil_cube_get_lapl(this, dim, order)
This module defines stencils used in Octopus.
Definition: stencil.F90:137
subroutine, public stencil_end(this)
Definition: stencil.F90:217
subroutine, public stencil_union(st1, st2, stu)
Definition: stencil.F90:229
subroutine, public symmetrizer_end(this)
subroutine, public symmetrizer_init(this, mesh, symm)
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:134
This module defines the unit system, used for input and output.
type(unit_t), public unit_angstrom
For XYZ files.
type(unit_system_t), public units_inp
the units systems for reading and writing
Class implementing a box generated from a 2D image.
Definition: box_image.F90:134
class to tell whether a point is inside or outside
Definition: box.F90:143
abstract class to describe coordinate systems
Description of the grid, containing information on derivatives, stencil, and symmetries.
Definition: grid.F90:171
Describes mesh distribution to nodes.
Definition: mesh.F90:187
This is defined even when running serial.
Definition: mpi.F90:144
Stores all communicators and groups.
Definition: multicomm.F90:208
The class representing the stencil, which is used for non-local mesh operations.
Definition: stencil.F90:165
int true(void)