Octopus
derivatives.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
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
27
29 use accel_oct_m
31 use batch_oct_m
34 use debug_oct_m
40 use global_oct_m
41 use iso_c_binding
42 use, intrinsic :: iso_fortran_env
45 use math_oct_m
46 use mesh_oct_m
51 use parser_oct_m
53 use space_oct_m
61 use types_oct_m
62 use utils_oct_m
64
65 implicit none
66
67 private
68 public :: &
105
106
107 integer, parameter :: &
108 DER_BC_ZERO_F = 0, & !< function is zero at the boundaries
109 der_bc_zero_df = 1, &
110 der_bc_period = 2
111
112 integer, parameter, public :: &
113 DER_STAR = 1, &
114 der_variational = 2, &
115 der_cube = 3, &
116 der_starplus = 4, &
118
119 integer, parameter :: &
120 BLOCKING = 1, &
121 non_blocking = 2
122
125 type derivatives_t
126 ! Components are public by default
127 type(boundaries_t) :: boundaries
128 type(mesh_t), pointer :: mesh => null()
129 integer :: dim = 0
130 integer, private :: periodic_dim = 0
131 integer :: order = 0
132 integer :: stencil_type = 0
133
134
135 logical, private :: stencil_primitive_coordinates
136
137 class(coordinate_system_t), private, pointer :: coord_system
138 logical :: remove_zero_weight_points
139
140 real(real64), allocatable :: masses(:)
141
144 real(real64), private :: lapl_cutoff = m_zero
145
146 type(nl_operator_t), allocatable, private :: op(:)
148 type(nl_operator_t), pointer :: lapl => null()
149 type(nl_operator_t), pointer :: grad(:) => null()
150
151 integer, allocatable :: n_ghost(:)
152!#if defined(HAVE_MPI)
153 integer, private :: comm_method = 0
154!#endif
155 type(derivatives_t), pointer :: finer => null()
156 type(derivatives_t), pointer :: coarser => null()
157 type(transfer_table_t), pointer :: to_finer => null()
158 type(transfer_table_t), pointer :: to_coarser => null()
159 end type derivatives_t
160
163 private
164!#ifdef HAVE_MPI
165 type(par_vec_handle_batch_t) :: pv_h
166!#endif
167 type(derivatives_t), pointer :: der
168 type(nl_operator_t), pointer :: op
169 type(batch_t), pointer :: ff
170 type(batch_t), pointer :: opff
171 logical :: ghost_update
172 logical :: factor_present
173 real(real64) :: factor
175
176 type(accel_kernel_t) :: kernel_uvw_xyz, kernel_dcurl, kernel_zcurl
177
178contains
179
180 ! ---------------------------------------------------------
181 subroutine derivatives_init(der, namespace, space, coord_system, order)
182 type(derivatives_t), target, intent(inout) :: der
183 type(namespace_t), intent(in) :: namespace
184 class(space_t), intent(in) :: space
185 class(coordinate_system_t), target, intent(in) :: coord_system
186 integer, optional, intent(in) :: order
187
188 integer :: idir
189 integer :: default_stencil
190 character(len=40) :: flags
191 logical :: stencil_primitive_coordinates
192
193 push_sub(derivatives_init)
194
195 ! Non-orthogonal curvilinear coordinates are currently not implemented
196 assert(.not. coord_system%local_basis .or. coord_system%orthogonal)
197
198 ! copy this value to my structure
199 der%dim = space%dim
200 der%periodic_dim = space%periodic_dim
201 der%coord_system => coord_system
203 !%Variable DerivativesStencil
204 !%Type integer
205 !%Default stencil_star
206 !%Section Mesh::Derivatives
207 !%Description
208 !% Decides what kind of stencil is used, <i>i.e.</i> which points, around
209 !% each point in the mesh, are the neighboring points used in the
210 !% expression of the differential operator.
211 !%
212 !% If curvilinear coordinates are to be used, then only the <tt>stencil_starplus</tt>
213 !% or the <tt>stencil_cube</tt> may be used. We only recommend the <tt>stencil_starplus</tt>,
214 !% since the cube typically needs far too much memory.
215 !%Option stencil_star 1
216 !% A star around each point (<i>i.e.</i>, only points on the axis).
217 !%Option stencil_variational 2
218 !% Same as the star, but with coefficients built in a different way.
219 !%Option stencil_cube 3
220 !% A cube of points around each point.
221 !%Option stencil_starplus 4
222 !% The star, plus a number of off-axis points.
223 !%Option stencil_stargeneral 5
224 !% The general star. Default for non-orthogonal grids.
225 !%End
226 default_stencil = der_star
227 if (coord_system%local_basis) default_stencil = der_starplus
228 if (.not. coord_system%orthogonal) default_stencil = der_stargeneral
229
230 call parse_variable(namespace, 'DerivativesStencil', default_stencil, der%stencil_type)
231
232 if (.not. varinfo_valid_option('DerivativesStencil', der%stencil_type)) then
233 call messages_input_error(namespace, 'DerivativesStencil')
234 end if
235 call messages_print_var_option("DerivativesStencil", der%stencil_type, namespace=namespace)
236
237 if (coord_system%local_basis .and. der%stencil_type < der_cube) call messages_input_error(namespace, 'DerivativesStencil')
238 if (der%stencil_type == der_variational) then
239 !%Variable DerivativesLaplacianFilter
240 !%Type float
241 !%Default 1.0
242 !%Section Mesh::Derivatives
243 !%Description
244 !% Undocumented
245 !%End
246 call parse_variable(namespace, 'DerivativesLaplacianFilter', m_one, der%lapl_cutoff)
247 end if
249 !%Variable DerivativesOrder
250 !%Type integer
251 !%Default 4
252 !%Section Mesh::Derivatives
253 !%Description
254 !% This variable gives the discretization order for the approximation of
255 !% the differential operators. This means, basically, that
256 !% <tt>DerivativesOrder</tt> points are used in each positive/negative
257 !% spatial direction, <i>e.g.</i> <tt>DerivativesOrder = 1</tt> would give
258 !% the well-known three-point formula in 1D.
259 !% The number of points actually used for the Laplacian
260 !% depends on the stencil used. Let <math>O</math> = <tt>DerivativesOrder</tt>, and <math>d</math> = <tt>Dimensions</tt>.
261 !% <ul>
262 !% <li> <tt>stencil_star</tt>: <math>2 O d + 1</math>
263 !% <li> <tt>stencil_cube</tt>: <math>(2 O + 1)^d</math>
264 !% <li> <tt>stencil_starplus</tt>: <math>2 O d + 1 + n</math> with <i>n</i> being 8
265 !% in 2D and 24 in 3D.
266 !% </ul>
267 !%End
268 call parse_variable(namespace, 'DerivativesOrder', 4, der%order)
269 ! overwrite order if given as argument
270 if (present(order)) then
271 der%order = order
272 end if
273
274 !%Variable DerivativesRemoveZeroWeightPoints
275 !%Type logical
276 !%Default yes
277 !%Section Mesh::Derivatives
278 !%Description
279 !% By default, Octopus removes the zero-weight points in the stencils.
280 !% This debug variable gives the possibility to deactive this feature.
281 !%End
282 call parse_variable(namespace, 'DerivativesRemoveZeroWeightPoints', .true., der%remove_zero_weight_points)
283
284
285#ifdef HAVE_MPI
286 !%Variable ParallelizationOfDerivatives
287 !%Type integer
288 !%Default non_blocking
289 !%Section Execution::Parallelization
290 !%Description
291 !% This option selects how the communication of mesh boundaries is performed.
292 !%Option blocking 1
293 !% Blocking communication.
294 !%Option non_blocking 2
295 !% Communication is based on non-blocking point-to-point communication.
296 !%End
297
298 call parse_variable(namespace, 'ParallelizationOfDerivatives', non_blocking, der%comm_method)
299
300 if (.not. varinfo_valid_option('ParallelizationOfDerivatives', der%comm_method)) then
301 call messages_input_error(namespace, 'ParallelizationOfDerivatives')
302 end if
303
304 call messages_obsolete_variable(namespace, 'OverlapDerivatives', 'ParallelizationOfDerivatives')
305#endif
306
307 !%Variable StencilPrimitiveCoordinates
308 !%Type logical
309 !%Section Mesh::Derivatives
310 !%Description
311 !% This variable controls the method for generating the stencil weights for the Laplacian.
312 !% For the gradient, the weights are always computed for primitive coordinates.
313 !% If set to yes, primitive coordinates are used for the polynomials and the right-hand side
314 !% of the Laplacian is computed using the metric tensor and the trace of the Hessian.
315 !% If set to no, Cartesian coordinates are used for the polynomials and the right-hand side
316 !% of the Laplacian reduces to the Cartesian case.
317 !% For some non-orthogonal grids, using primitve coordinates is necessary becasuse the polynomials
318 !% become linearly dependent, thus the corresponding matrix cannot be inverted. For some curvilinear
319 !% coordinate systems, using Cartesian coordinates is more accurate.
320 !%
321 !% By default, use primitve coordinates except for curvilinear meshes.
322 !%End
323 stencil_primitive_coordinates = .not. coord_system%local_basis
324 call parse_variable(namespace, 'StencilPrimitiveCoordinates', stencil_primitive_coordinates, der%stencil_primitive_coordinates)
325
326 ! if needed, der%masses should be initialized in modelmb_particles_init
327 safe_allocate(der%masses(1:space%dim))
328 der%masses = m_one
329
330 ! construct lapl and grad structures
331 safe_allocate(der%op(1:der%dim + 1))
332 der%grad => der%op
333 der%lapl => der%op(der%dim + 1)
334
335 call derivatives_get_stencil_lapl(der, der%lapl, space, coord_system)
337
338 ! find out how many ghost points we need in each dimension
339 safe_allocate(der%n_ghost(1:der%dim))
340 der%n_ghost(:) = 0
341 do idir = 1, der%dim
342 der%n_ghost(idir) = maxval(abs(der%lapl%stencil%points(idir, :)))
343 end do
344
345 nullify(der%coarser)
346 nullify(der%finer)
347 nullify(der%to_coarser)
348 nullify(der%to_finer)
349
350 if (accel_is_enabled()) then
351 ! Check if we can build the uvw_to_xyz kernel
352 select type (coord_system)
353 type is (cartesian_t)
354 ! In this case one does not need to call the kernel, so all is fine
355 class default
356 if (der%dim > 3) then
357 message(1) = "Calculation of derivatives on the GPU with dimension > 3 are only implemented for Cartesian coordinates."
358 call messages_fatal(1, namespace=namespace)
359 else
360 write(flags, '(A,I1.1)') ' -DDIMENSION=', der%dim
361 call accel_kernel_build(kernel_uvw_xyz, 'uvw_to_xyz.cu', 'uvw_to_xyz', flags)
362 end if
363 end select
364 call accel_kernel_build(kernel_dcurl, 'curl.cu', 'dcurl', flags = '-DRTYPE_DOUBLE')
365 call accel_kernel_build(kernel_zcurl, 'curl.cu', 'zcurl', flags = '-DRTYPE_COMPLEX')
366 end if
367
368 pop_sub(derivatives_init)
369 end subroutine derivatives_init
370
371
372 ! ---------------------------------------------------------
373 subroutine derivatives_end(der)
374 type(derivatives_t), intent(inout) :: der
375
376 integer :: idim
377
378 push_sub(derivatives_end)
379
380 assert(allocated(der%op))
381
382 do idim = 1, der%dim+1
383 call nl_operator_end(der%op(idim))
384 end do
385
386 safe_deallocate_a(der%masses)
387
388 safe_deallocate_a(der%n_ghost)
389
390 safe_deallocate_a(der%op)
391 nullify(der%lapl, der%grad)
392
393 nullify(der%coarser)
394 nullify(der%finer)
395 nullify(der%to_coarser)
396 nullify(der%to_finer)
397
398 call boundaries_end(der%boundaries)
399
400 pop_sub(derivatives_end)
401 end subroutine derivatives_end
402
403
404 ! ---------------------------------------------------------
405 subroutine derivatives_get_stencil_lapl(der, lapl, space, coord_system, name, order)
406 type(derivatives_t), intent(in) :: der
407 type(nl_operator_t), intent(inout) :: lapl
408 class(space_t), intent(in) :: space
409 class(coordinate_system_t), intent(in) :: coord_system
410 character(len=80), optional, intent(in) :: name
411 integer, optional, intent(in) :: order
412
413 character(len=80) :: name_
414 integer :: order_
415
417
418 name_ = optional_default(name, "Laplacian")
419 order_ = optional_default(order, der%order)
420
421 ! initialize nl operator
422 call nl_operator_init(lapl, name_, symm=op_symmetric)
423
424 ! create stencil
425 select case (der%stencil_type)
426 case (der_star, der_variational)
427 call stencil_star_get_lapl(lapl%stencil, der%dim, order_)
428 case (der_cube)
429 call stencil_cube_get_lapl(lapl%stencil, der%dim, order_)
430 case (der_starplus)
431 call stencil_starplus_get_lapl(lapl%stencil, der%dim, order_)
432 case (der_stargeneral)
433 call stencil_stargeneral_get_arms(lapl%stencil, space%dim, coord_system)
434 call stencil_stargeneral_get_lapl(lapl%stencil, der%dim, order_)
435 end select
436
438 end subroutine derivatives_get_stencil_lapl
439
440
441 ! ---------------------------------------------------------
443 subroutine derivatives_lapl_diag(der, lapl)
444 type(derivatives_t), intent(in) :: der
445 real(real64), intent(out) :: lapl(:)
446
447 push_sub(derivatives_lapl_diag)
448
449 assert(ubound(lapl, dim=1) >= der%mesh%np)
450
451 ! the Laplacian is a real operator
452 call dnl_operator_operate_diag(der%lapl, lapl)
453
454 pop_sub(derivatives_lapl_diag)
455
456 end subroutine derivatives_lapl_diag
457
458
459 ! ---------------------------------------------------------
460 subroutine derivatives_get_stencil_grad(der)
461 type(derivatives_t), intent(inout) :: der
462
463
464 integer :: ii
465 character :: dir_label
466
468
469 assert(associated(der%grad))
470
471 ! initialize nl operator
472 do ii = 1, der%dim
473 dir_label = ' '
474 if (ii < 5) dir_label = index2axis(ii)
475
476 call nl_operator_init(der%grad(ii), "Gradient "//dir_label, symm=op_antisymmetric)
478 ! create stencil
479 select case (der%stencil_type)
480 case (der_star, der_variational)
481 call stencil_star_get_grad(der%grad(ii)%stencil, der%dim, ii, der%order)
482 case (der_cube)
483 call stencil_cube_get_grad(der%grad(ii)%stencil, der%dim, der%order)
484 case (der_starplus)
485 call stencil_starplus_get_grad(der%grad(ii)%stencil, der%dim, ii, der%order)
486 case (der_stargeneral)
487 ! use the simple star stencil
488 call stencil_star_get_grad(der%grad(ii)%stencil, der%dim, ii, der%order)
489 end select
490 end do
491
493
494 end subroutine derivatives_get_stencil_grad
495
496 ! ---------------------------------------------------------
542 subroutine derivatives_build(der, namespace, space, mesh, qvector, regenerate, verbose)
543 type(derivatives_t), intent(inout) :: der
544 type(namespace_t), intent(in) :: namespace
545 class(space_t), intent(in) :: space
546 class(mesh_t), target, intent(in) :: mesh
547 real(real64), optional, intent(in) :: qvector(:)
549 logical, optional, intent(in) :: regenerate
550 logical, optional, intent(in) :: verbose
551
552 integer :: i
553 logical :: const_w_
554 integer :: np_zero_bc
555
556 push_sub(derivatives_build)
557
558 if (.not. optional_default(regenerate, .false.)) then
559 call boundaries_init(der%boundaries, namespace, space, mesh, qvector)
560 end if
561
562 assert(allocated(der%op))
563 assert(der%stencil_type >= der_star .and. der%stencil_type <= der_stargeneral)
564 assert(.not.(der%stencil_type == der_variational .and. mesh%use_curvilinear))
565
566 der%mesh => mesh ! make a pointer to the underlying mesh
567
568 const_w_ = .true.
569
570 ! need non-constant weights for curvilinear and scattering meshes
571 if (mesh%use_curvilinear) const_w_ = .false.
572
573 np_zero_bc = 0
574
575 if (optional_default(regenerate, .false.)) then
576 do i = 1, der%dim+1
577 call nl_operator_end(der%op(i))
578 safe_deallocate_a(der%op(i)%w)
579 end do
580 call derivatives_get_stencil_lapl(der, der%lapl, space, der%coord_system)
582 end if
583
584 ! build operators
585 do i = 1, der%dim+1
586 call nl_operator_build(space, mesh, der%op(i), der%mesh%np, const_w = const_w_, &
587 regenerate=regenerate)
588 np_zero_bc = max(np_zero_bc, nl_operator_np_zero_bc(der%op(i)))
589 end do
590
591 assert(np_zero_bc > mesh%np .and. np_zero_bc <= mesh%np_part)
592
593 select case (der%stencil_type)
594
595 case (der_star, der_starplus, der_stargeneral) ! Laplacian and gradient have different stencils
596 do i = 1, der%dim + 1
597
598 call derivatives_make_discretization(namespace, der, 1, i, der%op(i:i), space, verbose=verbose)
599 end do
600
601 case (der_cube)
602 ! gradients have the same stencils, so use one call to derivatives_make_discretization
603 ! to solve the linear equation once for several right-hand sides
604 call derivatives_make_discretization(namespace, der, der%dim, der%dim, der%op(:), space, verbose=verbose)
605 ! the polynomials are evaluated differently for the Laplacian
606 call derivatives_make_discretization(namespace, der, 1, der%dim+1, der%op(der%dim+1:der%dim+1), space, &
607 verbose=verbose)
608
609 case (der_variational)
610 ! we have the explicit coefficients
611 call stencil_variational_coeff_lapl(der%dim, der%order, mesh%spacing, der%lapl, alpha = der%lapl_cutoff)
612 end select
613
615
616 end subroutine derivatives_build
617
618 ! ---------------------------------------------------------
619 subroutine stencil_pol_grad(stencil_type, dim, direction, order, polynomials)
620 integer, intent(in) :: stencil_type
621 integer, intent(in) :: dim
622 integer, intent(in) :: direction
623 integer, intent(in) :: order
624 integer, intent(inout) :: polynomials(:, :)
625
626 select case (stencil_type)
627 case (der_star, der_stargeneral)
628 call stencil_star_polynomials_grad(direction, order, polynomials)
629 case (der_starplus)
630 call stencil_starplus_pol_grad(dim, direction, order, polynomials)
631 case (der_cube)
632 ! same stencil for gradient and Laplacian in this case
633 call stencil_cube_polynomials_lapl(dim, order, polynomials)
634 end select
635 end subroutine stencil_pol_grad
636
637 ! ---------------------------------------------------------
638 subroutine stencil_pol_lapl(stencil_type, stencil, dim, order, polynomials)
639 integer, intent(in) :: stencil_type
640 type(stencil_t), intent(in) :: stencil
641 integer, intent(in) :: dim
642 integer, intent(in) :: order
643 integer, intent(inout) :: polynomials(:, :)
644
645 select case (stencil_type)
646 case (der_star)
647 call stencil_star_polynomials_lapl(dim, order, polynomials)
648 case (der_starplus)
649 call stencil_starplus_pol_lapl(dim, order, polynomials)
650 case (der_stargeneral)
651 call stencil_stargeneral_pol_lapl(stencil, dim, order, polynomials)
652 case (der_cube)
653 call stencil_cube_polynomials_lapl(dim, order, polynomials)
654 end select
655 end subroutine stencil_pol_lapl
656
657 ! ---------------------------------------------------------
658 subroutine get_rhs_lapl(coord_system, polynomials, chi, rhs, stencil_primitive_coordinates)
659 class(coordinate_system_t), intent(in) :: coord_system
660 integer, intent(in) :: polynomials(:,:)
661 real(real64), intent(in) :: chi(:)
662 real(real64), intent(out) :: rhs(:)
663 logical, intent(in) :: stencil_primitive_coordinates
664
665 integer :: i, j, k, dim
666 real(real64) :: metric_inverse(1:size(polynomials, dim=1), 1:size(polynomials, dim=1))
667 real(real64) :: hessian_trace(1:size(polynomials, dim=1))
668 integer :: powers(0:2)
669
670 push_sub(get_rhs_lapl)
671
672 dim = size(polynomials, dim=1)
673
674 if (stencil_primitive_coordinates) then
675 ! inverse metric and trace of Hessian in primitive coordinates
676 metric_inverse = coord_system%metric_inverse(chi)
677 hessian_trace = coord_system%trace_hessian(chi)
678 else
679 ! inverse metric and trace of Hessian in Cartesian coordinates
680 metric_inverse = m_zero
681 do i = 1, dim
682 metric_inverse(i, i) = m_one
683 end do
684 hessian_trace = m_zero
685 end if
686
687 ! find right-hand side for operator
688 rhs(:) = m_zero
689 do j = 1, size(polynomials, dim=2)
690 ! count the powers of the polynomials
691 powers = 0
692 do i = 1, dim
693 if (polynomials(i, j) <= 2) then
694 powers(polynomials(i, j)) = powers(polynomials(i, j)) + 1
695 end if
696 end do
697
698 ! find all polynomials for which exactly one term is quadratic
699 ! for these, the Laplacian on the polynomial is 2*metric_inverse(i, i)
700 if (powers(2) == 1 .and. powers(0) == dim - 1) then
701 do i = 1, dim
702 if (polynomials(i, j) == 2) then
703 rhs(j) = m_two*metric_inverse(i, i)
704 end if
705 end do
706 end if
707 ! find all polynomials for which exactly two terms are linear
708 ! for these, the Laplacian on the polynomial is metric_inverse(i, k) + metric_inverse(k, i)
709 if (powers(1) == 2 .and. powers(0) == dim - 2) then
710 do i = 1, dim
711 if (polynomials(i, j) == 1) then
712 do k = i+1, dim
713 if (polynomials(k, j) == 1) then
714 rhs(j) = metric_inverse(i, k) + metric_inverse(k, i)
715 end if
716 end do
717 end if
718 end do
719 end if
720 ! find all polynomials for which exactly one term is linear
721 ! for these, the Laplacian on the polynomial is hessian_trace(i)
722 ! this term is only non-zero for curvilinear coordinates with a varying metric tensor
723 if (powers(1) == 1 .and. powers(0) == dim - 1) then
724 do i = 1, dim
725 if (polynomials(i, j) == 1) then
726 rhs(j) = hessian_trace(i)
727 end if
728 end do
729 end if
730 end do
731
732 pop_sub(get_rhs_lapl)
733 end subroutine get_rhs_lapl
734
735 ! ---------------------------------------------------------
736 subroutine get_rhs_grad(polynomials, dir, rhs)
737 integer, intent(in) :: polynomials(:,:)
738 integer, intent(in) :: dir
739 real(real64), intent(out) :: rhs(:)
740
741 integer :: j, k, dim
742 logical :: this_one
743
744 push_sub(get_rhs_grad)
745
746 dim = size(polynomials, dim=1)
747
748 ! find right-hand side for operator
749 rhs(:) = m_zero
750 do j = 1, size(polynomials, dim=2)
751 this_one = .true.
752 do k = 1, dim
753 if (k == dir .and. polynomials(k, j) /= 1) this_one = .false.
754 if (k /= dir .and. polynomials(k, j) /= 0) this_one = .false.
755 end do
756 if (this_one) rhs(j) = m_one
757 end do
758
759 pop_sub(get_rhs_grad)
760 end subroutine get_rhs_grad
761
762
763 ! ---------------------------------------------------------
764 subroutine derivatives_make_discretization(namespace, der, nderiv, ideriv, op, &
765 space, name, verbose, order)
766 type(namespace_t), intent(in) :: namespace
767 type(derivatives_t), intent(in) :: der
768 integer, intent(in) :: nderiv
769 integer, intent(in) :: ideriv
770 type(nl_operator_t), intent(inout) :: op(:)
771 type(space_t), intent(in) :: space
772 character(len=80), optional, intent(in) :: name
773 logical, optional, intent(in) :: verbose
774 integer, optional, intent(in) :: order
775
776 integer :: p, p_max, i, j, k, pow_max, order_
777 real(real64) :: x(der%dim)
778 real(real64), allocatable :: mat(:,:), sol(:,:), powers(:,:)
779 integer, allocatable :: pol(:,:)
780 real(real64), allocatable :: rhs(:,:)
781 character(len=80) :: name_
782
784
785 order_ = optional_default(order, der%order)
786
787 safe_allocate(rhs(1:op(1)%stencil%size, nderiv))
788 safe_allocate(pol(1:der%dim, 1:op(1)%stencil%size))
789 ! get polynomials
790 if (nderiv == 1) then
791 if (ideriv <= der%dim) then ! gradient
792 call stencil_pol_grad(der%stencil_type, der%dim, ideriv, order_, pol)
793 name_ = index2axis(ideriv) // "-gradient"
794 else ! Laplacian
795 call stencil_pol_lapl(der%stencil_type, op(1)%stencil, der%dim, order_, pol)
796 name_ = "Laplacian"
797 end if
798 else if (nderiv == der%dim) then
799 ! in this case, we compute the weights for the gradients at once
800 ! because they have the same stencils
801 ! this is only used for the cube stencil
802 name_ = "gradients"
803 call stencil_pol_grad(der%stencil_type, der%dim, 1, order_, pol)
804 else
805 message(1) = "Error: derivatives_make_discretization can only be called with nderiv = 1 or der%dim."
806 call messages_fatal(1)
807 end if
809 safe_allocate(mat(1:op(1)%stencil%size, 1:op(1)%stencil%size))
810 safe_allocate(sol(1:op(1)%stencil%size, 1:nderiv))
811
812 if (optional_default(verbose, .true.)) then
813 name_ = optional_default(name, name_)
814 message(1) = 'Info: Generating weights for finite-difference discretization of ' // trim(name_)
815 call messages_info(1, namespace=namespace)
816 end if
817
818 ! use to generate power lookup table
819 pow_max = maxval(pol)
820 safe_allocate(powers(1:der%dim, 0:pow_max))
821 powers(:,:) = m_zero
822 powers(:,0) = m_one
823
824 p_max = op(1)%np
825 if (op(1)%const_w) p_max = 1
826
827 do p = 1, p_max
828 ! get polynomials and right-hand side
829 ! we need the current position for the right-hand side for curvilinear coordinates
830 if (nderiv == 1) then
831 if (ideriv <= der%dim) then ! gradient
832 call get_rhs_grad(pol, ideriv, rhs(:,1))
833 name_ = index2axis(ideriv) // "-gradient"
834 else ! Laplacian
835 call get_rhs_lapl(der%mesh%coord_system, pol, der%mesh%chi(:, p), rhs(:,1), der%stencil_primitive_coordinates)
836 end if
837 else if (nderiv == der%dim) then
838 ! in this case, we compute the weights for the gradients at once
839 ! because they have the same stencils
840 do i = 1, der%dim
841 call get_rhs_grad(pol, i, rhs(:,i))
842 end do
843 end if
844
845 ! first polynomial is just a constant
846 mat(1,:) = m_one
847 ! i indexes the point in the stencil
848 do i = 1, op(1)%stencil%size
849 if (ideriv <= der%dim) then
850 ! for gradients, we always need to evaluate the polynomials in primitive coordinates
851 x = real(op(1)%stencil%points(:, i), real64)*der%mesh%spacing
852 else
853 ! for the Laplacian, we can evaluate the polynomials in primitive or Cartesian coordinates
854 if (der%stencil_primitive_coordinates) then
855 x = real(op(1)%stencil%points(:, i), real64)*der%mesh%spacing
856 else
857 x = der%mesh%x(:, p + op(1)%ri(i, op(1)%rimap(p))) - der%mesh%x(:, p)
858 end if
859 end if
860
861 ! NB: these masses are applied on the cartesian directions. Should add a check for non-orthogonal axes
862 x = x*sqrt(der%masses)
863
864 ! calculate powers
865 powers(:, 1) = x
866 do k = 2, pow_max
867 powers(:, k) = x*powers(:, k-1)
868 end do
869
870 ! generate the matrix
871 ! j indexes the polynomial being used
872 do j = 2, op(1)%stencil%size
873 mat(j, i) = powers(1, pol(1, j))
874 do k = 2, der%dim
875 mat(j, i) = mat(j, i)*powers(k, pol(k, j))
876 end do
877 end do
878 end do ! loop over i = point in stencil
879
880 ! linear problem to solve for derivative weights:
881 ! mat * sol = rhs
882 call lalg_linsyssolve(op(1)%stencil%size, nderiv, mat, rhs, sol)
883
884 ! for the cube stencil, all derivatives are calculated at once, so assign
885 ! the correct solution to each operator
886 do i = 1, nderiv
887 op(i)%w(:, p) = sol(:, i)
888 end do
889
890 end do ! loop over points p
891
892 do i = 1, nderiv
893 if (der%remove_zero_weight_points) then
894 call nl_operator_remove_zero_weight_points(op(i), space, der%mesh)
895 end if
897 if (op(i)%const_w) then
899 end if
900 end do
901
902 ! In case of constant weights, we store the weights of the Laplacian on the GPU, as this
903 ! saves many unecessary transfers
904 do i = 1, nderiv
907 end do
908
909
910 safe_deallocate_a(mat)
911 safe_deallocate_a(sol)
912 safe_deallocate_a(powers)
913 safe_deallocate_a(pol)
914 safe_deallocate_a(rhs)
915
918
919#ifdef HAVE_MPI
920 ! ---------------------------------------------------------
921 logical function derivatives_overlap(this) result(overlap)
922 type(derivatives_t), intent(in) :: this
923
924 push_sub(derivatives_overlap)
925
926 overlap = this%comm_method /= blocking
927
928 pop_sub(derivatives_overlap)
929 end function derivatives_overlap
930#endif
931
932 ! ---------------------------------------------------------
933 subroutine derivatives_get_lapl(this, namespace, op, space, name, order)
934 type(derivatives_t), intent(in) :: this
935 type(namespace_t), intent(in) :: namespace
936 type(nl_operator_t), intent(inout) :: op(:)
937 class(space_t), intent(in) :: space
938 character(len=80), intent(in) :: name
939 integer, intent(in) :: order
940
941 push_sub(derivatives_get_lapl)
942
943 call derivatives_get_stencil_lapl(this, op(1), space, this%mesh%coord_system, name, order)
944 call nl_operator_build(space, this%mesh, op(1), this%mesh%np, const_w = .not. this%mesh%use_curvilinear)
945 call derivatives_make_discretization(namespace, this, 1, this%dim+1, op(1:1), space, name, order=order)
946
947 pop_sub(derivatives_get_lapl)
948 end subroutine derivatives_get_lapl
949
950 ! ---------------------------------------------------------
963 function derivatives_get_inner_boundary_mask(this) result(mask)
964 type(derivatives_t), intent(in) :: this
965
966 logical :: mask(1:this%mesh%np)
967 integer :: ip, is, index
968
969 mask = .false.
970 ! Loop through all points in the grid
971 do ip = 1, this%mesh%np
972 ! For each of them, loop through all points in the stencil
973 do is = 1, this%lapl%stencil%size
974 ! Get the index of the point obtained as: grid_point + displament_due_to_stencil
975 index = nl_operator_get_index(this%lapl, is, ip)
976 ! Check whether the displaced point if outsude the grid. Is so, it belongs to the mask
977 if (index > this%mesh%np + this%mesh%pv%np_ghost) then
978 mask(ip) = .true.
979 exit
980 end if
981 end do
982 end do
983
985
986 ! ---------------------------------------------------------
991 real(real64) function derivatives_lapl_get_max_eigenvalue(this)
992 type(derivatives_t), intent(in) :: this
993
994 integer :: i
995
996 push_sub(derivatives_lapl_get_max_eigenvalue)
997
998 derivatives_lapl_get_max_eigenvalue = m_zero
999 if ((this%stencil_type == der_star .or. this%stencil_type == der_stargeneral) &
1000 .and. .not. this%mesh%use_curvilinear) then
1001 ! use Fourier transform of stencil evaluated at the maximum phase
1002 do i = 1, this%lapl%stencil%size
1003 derivatives_lapl_get_max_eigenvalue = derivatives_lapl_get_max_eigenvalue + &
1004 (-1)**maxval(abs(this%lapl%stencil%points(:, i)))*this%lapl%w(i, 1)
1005 end do
1006 derivatives_lapl_get_max_eigenvalue = abs(derivatives_lapl_get_max_eigenvalue)
1007 else
1008 ! use upper bound from continuum for other stencils
1009 do i = 1, this%dim
1010 derivatives_lapl_get_max_eigenvalue = derivatives_lapl_get_max_eigenvalue + &
1011 m_pi**2/this%mesh%spacing(i)**2
1012 end do
1013 end if
1014
1015 pop_sub(derivatives_lapl_get_max_eigenvalue)
1017
1018 subroutine derivates_set_coordinates_system(der, coord_system)
1019 type(derivatives_t), target, intent(inout) :: der
1020 class(coordinate_system_t), target, intent(in) :: coord_system
1021
1022 der%coord_system => coord_system
1024
1025#include "undef.F90"
1026#include "real.F90"
1027#include "derivatives_inc.F90"
1028
1029#include "undef.F90"
1030#include "complex.F90"
1031#include "derivatives_inc.F90"
1032
1033end module derivatives_oct_m
1034
1035!! Local Variables:
1036!! mode: f90
1037!! coding: utf-8
1038!! End:
double sqrt(double __x) __attribute__((__nothrow__
subroutine, public accel_kernel_build(this, file_name, kernel_name, flags)
Definition: accel.F90:1227
pure logical function, public accel_is_enabled()
Definition: accel.F90:392
This module implements batches of mesh functions.
Definition: batch.F90:135
This module implements common operations on batches of mesh functions.
Definition: batch_ops.F90:118
Module implementing boundary conditions in Octopus.
Definition: boundaries.F90:124
subroutine, public boundaries_end(this)
Definition: boundaries.F90:430
subroutine, public boundaries_init(this, namespace, space, mesh, qvector)
initialize the boundary contitions
Definition: boundaries.F90:229
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 zderivatives_curl(der, ff, op_ff, ghost_update, set_bc)
apply the curl operator to a vector of mesh functions
subroutine, public dderivatives_perform(op, der, ff, op_ff, ghost_update, set_bc, factor)
apply a nl_operator to a mesh function
subroutine, public derivatives_lapl_diag(der, lapl)
Returns the diagonal elements of the Laplacian, needed for preconditioning.
integer, parameter, public der_cube
subroutine, public dderivatives_batch_curl_from_gradient(der, ffb, gradb)
calculate the curl from a batch and its gradient
subroutine derivatives_make_discretization(namespace, der, nderiv, ideriv, op, space, name, verbose, order)
subroutine, public derivatives_build(der, namespace, space, mesh, qvector, regenerate, verbose)
build the derivatives object:
subroutine, public dderivatives_batch_grad(der, ffb, opffb, ghost_update, set_bc, to_cartesian, factor)
apply the gradient to a batch of mesh functions
subroutine stencil_pol_grad(stencil_type, dim, direction, order, polynomials)
subroutine stencil_pol_lapl(stencil_type, stencil, dim, order, polynomials)
subroutine, public zderivatives_batch_curl_from_gradient(der, ffb, gradb)
calculate the curl from a batch and its gradient
subroutine, public zderivatives_batch_div(der, ffb, opffb, ghost_update, set_bc, to_cartesian)
calculate the divergence of a vector of batches
subroutine, public dderivatives_partial(der, ff, op_ff, dir, ghost_update, set_bc)
apply the partial derivative along dir to a mesh function
subroutine get_rhs_lapl(coord_system, polynomials, chi, rhs, stencil_primitive_coordinates)
subroutine, public dderivatives_batch_perform(op, der, ff, opff, ghost_update, set_bc, factor, async)
apply an operator to a bach of mesh functions
subroutine, public zderivatives_div(der, ff, op_ff, ghost_update, set_bc, to_cartesian)
apply the divergence operator to a vector of mesh functions
integer, parameter der_bc_period
boundary is periodic
subroutine, public dderivatives_batch_finish(handle, async)
apply a non-local operator to a batch (2nd part)
subroutine, public dderivatives_curl(der, ff, op_ff, ghost_update, set_bc)
apply the curl operator to a vector of mesh functions
subroutine derivatives_get_stencil_grad(der)
subroutine, public dderivatives_grad(der, ff, op_ff, ghost_update, set_bc, to_cartesian)
apply the gradient to a mesh function
subroutine derivatives_get_stencil_lapl(der, lapl, space, coord_system, name, order)
subroutine, public derivatives_init(der, namespace, space, coord_system, order)
subroutine, public derivatives_end(der)
integer, parameter der_bc_zero_df
first derivative of the function is zero
subroutine, public zderivatives_batch_grad(der, ffb, opffb, ghost_update, set_bc, to_cartesian, factor)
apply the gradient to a batch of mesh functions
subroutine, public dderivatives_batch_div(der, ffb, opffb, ghost_update, set_bc, to_cartesian)
calculate the divergence of a vector of batches
subroutine, public derivatives_get_lapl(this, namespace, op, space, name, order)
subroutine, public zderivatives_batch_perform(op, der, ff, opff, ghost_update, set_bc, factor, async)
apply an operator to a bach of mesh functions
subroutine, public dderivatives_batch_start(op, der, ff, opff, handle, ghost_update, set_bc, factor)
apply a non-local operator to a batch (1st part)
subroutine get_rhs_grad(polynomials, dir, rhs)
subroutine, public dderivatives_lapl(der, ff, op_ff, ghost_update, set_bc, factor)
apply the Laplacian to a mesh function
real(real64) function, public derivatives_lapl_get_max_eigenvalue(this)
Get maximum eigenvalue of discrete Laplacian. For the star and star_general stencils,...
subroutine, public dderivatives_batch_curl(der, ffb, ghost_update, set_bc)
apply the curl to a batch of mesh functions
subroutine, public zderivatives_grad(der, ff, op_ff, ghost_update, set_bc, to_cartesian)
apply the gradient to a mesh function
subroutine, public zderivatives_partial(der, ff, op_ff, dir, ghost_update, set_bc)
apply the partial derivative along dir to a mesh function
subroutine, public zderivatives_batch_start(op, der, ff, opff, handle, ghost_update, set_bc, factor)
apply a non-local operator to a batch (1st part)
logical function, dimension(1:this%mesh%np), public derivatives_get_inner_boundary_mask(this)
This function tells whether a point in the grid is contained in a layer of the width of the stencil b...
integer, parameter, public der_starplus
subroutine, public zderivatives_perform(op, der, ff, op_ff, ghost_update, set_bc, factor)
apply a nl_operator to a mesh function
subroutine, public zderivatives_batch_curl(der, ffb, ghost_update, set_bc)
apply the curl to a batch of mesh functions
integer, parameter, public der_variational
subroutine, public zderivatives_lapl(der, ff, op_ff, ghost_update, set_bc, factor)
apply the Laplacian to a mesh function
subroutine, public derivates_set_coordinates_system(der, coord_system)
subroutine, public zderivatives_batch_finish(handle, async)
apply a non-local operator to a batch (2nd part)
subroutine, public dderivatives_div(der, ff, op_ff, ghost_update, set_bc, to_cartesian)
apply the divergence operator to a vector of mesh functions
integer, parameter non_blocking
integer, parameter, public der_stargeneral
real(real64), parameter, public m_two
Definition: global.F90:193
real(real64), parameter, public m_zero
Definition: global.F90:191
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:189
real(real64), parameter, public m_one
Definition: global.F90:192
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:117
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1023
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_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:594
This module defines non-local operators.
subroutine, public dnl_operator_operate_diag(op, fo)
subroutine, public nl_operator_init(op, label, symm)
initialize an instance of a non-local operator by setting the label
subroutine, public nl_operator_build_symmetric_weights(op, max_size)
Builds (or rebuild) the necessary arrays for symmetric and antisymmetric stencils.
subroutine, public nl_operator_update_gpu_buffers(op)
subroutine, public nl_operator_output_weights(this)
subroutine, public nl_operator_end(op)
subroutine, public nl_operator_remove_zero_weight_points(op, space, mesh)
Removes the zero-weight points for constant weight stencils.
integer, parameter, public op_symmetric
subroutine, public nl_operator_build(space, mesh, op, np, const_w, regenerate)
Creates the nonlocal operators for the stencils used for finite differences.
integer pure function, public nl_operator_np_zero_bc(op)
integer, parameter, public op_antisymmetric
integer pure function, public nl_operator_get_index(op, is, ip)
subroutine, public nl_operator_allocate_gpu_buffers(op)
Some general things and nomenclature:
Definition: par_vec.F90:173
This module defines routines, generating operators for a cubic stencil.
subroutine, public stencil_cube_get_lapl(this, dim, order)
subroutine, public stencil_cube_polynomials_lapl(dim, order, pol)
subroutine, public stencil_cube_get_grad(this, dim, order)
This module defines stencils used in Octopus.
Definition: stencil.F90:137
This module defines routines, generating operators for a star stencil.
subroutine, public stencil_star_get_grad(this, dim, dir, order)
subroutine, public stencil_star_get_lapl(this, dim, order)
subroutine, public stencil_star_polynomials_grad(dir, order, pol)
subroutine, public stencil_star_polynomials_lapl(dim, order, pol)
This module defines routines, generating operators for a generalized star stencil.
subroutine, public stencil_stargeneral_get_arms(this, dim, coord_system)
Finds the direction of the arms of the star-general stencil as described in Natan et al....
subroutine, public stencil_stargeneral_pol_lapl(this, dim, order, pol)
subroutine, public stencil_stargeneral_get_lapl(this, dim, order)
This module defines routines, generating operators for a stencil consisting of a star and a cross....
subroutine, public stencil_starplus_pol_grad(dim, dir, order, pol)
subroutine, public stencil_starplus_pol_lapl(dim, order, pol)
subroutine, public stencil_starplus_get_lapl(this, dim, order)
subroutine, public stencil_starplus_get_grad(this, dim, dir, order)
Implements the variational discretization of the Laplacian as proposed by P. Maragakis,...
subroutine, public stencil_variational_coeff_lapl(dim, order, h, lapl, alpha)
This module is intended to contain simple general-purpose utility functions and procedures.
Definition: utils.F90:120
character pure function, public index2axis(idir)
Definition: utils.F90:205
abstract class to describe coordinate systems
handle to transfer data from the start() to finish() calls.
class representing derivatives
Describes mesh distribution to nodes.
Definition: mesh.F90:187
data type for non local operators
The class representing the stencil, which is used for non-local mesh operations.
Definition: stencil.F90:165
int true(void)