Octopus
mesh_batch.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
24 use accel_oct_m
26 use batch_oct_m
28 use blas_oct_m
29 use iso_c_binding
30 use comm_oct_m
31 use debug_oct_m
32 use global_oct_m
35 use math_oct_m
36 use mesh_oct_m
39 use mpi_oct_m
42#if defined(HAVE_OPENMP)
43 use omp_lib
44#endif
48 use types_oct_m
49
50 implicit none
51
52 private
53 public :: &
54 batch_p_t, &
72
73 type batch_p_t
74 class(batch_t), allocatable :: p
75 end type batch_p_t
76
77
78contains
79
80 ! -----------------------------------------------------
82 !
83 subroutine mesh_batch_nrm2(mesh, aa, nrm2, reduce)
84 class(mesh_t), intent(in) :: mesh
85 class(batch_t), intent(in) :: aa
86 real(real64), contiguous, intent(out) :: nrm2(:)
87 logical, optional, intent(in) :: reduce
88
89 push_sub(mesh_batch_nrm2)
90
91 if (aa%type() == type_float) then
92 call dpriv_mesh_batch_nrm2(mesh, aa, nrm2)
93 else
94 call zpriv_mesh_batch_nrm2(mesh, aa, nrm2)
95 end if
96
97 if (mesh%parallel_in_domains .and. optional_default(reduce, .true.)) then
98 nrm2(1:aa%nst) = nrm2(1:aa%nst)**2
99 call mesh%allreduce(nrm2, dim = aa%nst)
100 nrm2(1:aa%nst) = sqrt(nrm2(1:aa%nst))
101 end if
102
103 pop_sub(mesh_batch_nrm2)
104 end subroutine mesh_batch_nrm2
105
106#include "undef.F90"
107#include "real.F90"
108#include "mesh_batch_inc.F90"
109
110#include "undef.F90"
111#include "complex.F90"
112#include "mesh_batch_inc.F90"
113
114end module mesh_batch_oct_m
115
116!! Local Variables:
117!! mode: f90
118!! coding: utf-8
119!! End:
double sqrt(double __x) __attribute__((__nothrow__
This module implements batches of mesh functions.
Definition: batch.F90:133
This module implements common operations on batches of mesh functions.
Definition: batch_ops.F90:116
This module contains interfaces for BLAS routines You should not use these routines directly....
Definition: blas.F90:118
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
This module defines functions over batches of mesh functions.
Definition: mesh_batch.F90:116
subroutine, public dmesh_batch_dotp_matrix(mesh, aa, bb, dot, reduce)
Calculate the overlap matrix of two batches.
Definition: mesh_batch.F90:271
subroutine, public dmesh_batch_codensity(mesh, aa, psi, rho)
calculate the co-densities
Definition: mesh_batch.F90:965
subroutine, public dmesh_batch_dotp_self(mesh, aa, dot, reduce)
calculate the overlap matrix of a batch with itself
Definition: mesh_batch.F90:517
subroutine, public mesh_batch_nrm2(mesh, aa, nrm2, reduce)
Calculate the norms (norm2, not the square!) of a batch of mesh functions.
Definition: mesh_batch.F90:177
subroutine, public zmesh_batch_exchange_points(mesh, aa, forward_map, backward_map)
This functions exchanges points of a mesh according to a certain map. Two possible maps can be given....
subroutine, public zmesh_batch_dotp_vector(mesh, aa, bb, dot, reduce, cproduct)
calculate the vector of dot-products of mesh functions between two batches
subroutine zpriv_mesh_batch_nrm2(mesh, aa, nrm2)
This function should not be called directly, but through mesh_batch_nrm2.
subroutine, public zmesh_batch_dotp_matrix(mesh, aa, bb, dot, reduce)
Calculate the overlap matrix of two batches.
subroutine, public zmesh_batch_normalize(mesh, psib, norm)
Normalize a batch.
subroutine, public dmesh_batch_orthogonalization(mesh, nst, psib, phib, normalize, overlap, norm, gs_scheme, full_batch)
Orthonormalizes states of phib to the orbitals of nst batches of psi.
subroutine, public zmesh_batch_orthogonalization(mesh, nst, psib, phib, normalize, overlap, norm, gs_scheme, full_batch)
Orthonormalizes states of phib to the orbitals of nst batches of psi.
subroutine, public dmesh_batch_exchange_points(mesh, aa, forward_map, backward_map)
This functions exchanges points of a mesh according to a certain map. Two possible maps can be given....
subroutine, public zmesh_batch_codensity(mesh, aa, psi, rho)
calculate the co-densities
subroutine, public zmesh_batch_mf_dotp(mesh, aa, psi, dot, reduce, nst)
calculate the dot products between a batch and a vector of mesh functions
subroutine, public dmesh_batch_mf_dotp(mesh, aa, psi, dot, reduce, nst)
calculate the dot products between a batch and a vector of mesh functions
Definition: mesh_batch.F90:805
subroutine dpriv_mesh_batch_nrm2(mesh, aa, nrm2)
This function should not be called directly, but through mesh_batch_nrm2.
subroutine, public zmesh_batch_dotp_self(mesh, aa, dot, reduce)
calculate the overlap matrix of a batch with itself
subroutine, public dmesh_batch_normalize(mesh, psib, norm)
Normalize a batch.
subroutine, public dmesh_batch_dotp_vector(mesh, aa, bb, dot, reduce, cproduct)
calculate the vector of dot-products of mesh functions between two batches
Definition: mesh_batch.F90:632
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:145
Some general things and nomenclature:
Definition: par_vec.F90:171
type(type_t), public type_float
Definition: types.F90:133
int true(void)