Octopus
modelmb_1part.F90
Go to the documentation of this file.
1!! Copyright (C) 2009 N. Helbig and M. Verstraete
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 debug_oct_m
25 use global_oct_m
28 use mesh_oct_m
30
31 implicit none
32
33 private
34
35 public :: &
39
42 ! All components are public by default
43 integer :: ndim1part
44 integer, private :: npt_part
45 integer :: npt
46 real(real64) :: vol_elem_1part
47 real(real64), allocatable :: origin(:)
48 integer, allocatable :: enlarge_1part(:)
49 integer, allocatable :: nr_1part(:,:)
50 integer, allocatable, private :: ll(:)
51 real(real64), allocatable :: h_1part(:)
52 type(hypercube_t) :: hypercube_1part
53 end type modelmb_1part_t
54
55contains
56
57 subroutine modelmb_1part_init(this, mesh, ikeeppart, ndim1part)
58 type(modelmb_1part_t), intent(out) :: this
59 type(mesh_t), intent(in) :: mesh
60 integer, intent(in) :: ikeeppart
61 integer, intent(in) :: ndim1part
62
63 integer :: idir, irealdir
64
65 push_sub(modelmb_1part_init)
66
67 this%ndim1part = ndim1part
68
69! get full size of arrays for 1 particle only in ndim_modelmb dimensions
70 this%npt_part = 1
71 do idir = 1, ndim1part
72 this%npt_part = this%npt_part*(mesh%idx%nr(2,(ikeeppart - 1)*ndim1part + idir) &
73 - mesh%idx%nr(1,(ikeeppart - 1)*ndim1part + idir) + 1)
74 end do
75
76!real bounds for indices in
77 this%npt = 1
78 safe_allocate(this%ll(1:ndim1part))
79 do idir = 1, ndim1part
80 this%ll(idir) = mesh%idx%ll((ikeeppart-1)*ndim1part+idir)
81 this%npt = this%npt*this%ll(idir)
82 end do
83
84! volume element for the chosen particle
85 safe_allocate(this%h_1part(1:ndim1part))
86 this%vol_elem_1part = m_one
87 do idir = 1,ndim1part
88 irealdir = (ikeeppart-1)*ndim1part + idir
89 this%vol_elem_1part = this%vol_elem_1part*mesh%spacing(irealdir)
90 this%h_1part(idir) = mesh%spacing(irealdir)
91 end do
92
93! store start and end positions for the relevant dimensions for this particle
94 safe_allocate(this%nr_1part(1:2, 1:ndim1part))
95 this%nr_1part(:,:) = mesh%idx%nr(:,(ikeeppart-1)*ndim1part+1:ikeeppart*ndim1part)
96
97! initialize a hypercube for just this particle
98! NB: hypercube_* presume that enlarge is the same for all dimensions!
99 safe_allocate(this%enlarge_1part(1:ndim1part))
100 this%enlarge_1part = mesh%idx%enlarge((ikeeppart-1)*ndim1part+1:ikeeppart*ndim1part)
101 call hypercube_init(this%hypercube_1part, ndim1part, this%nr_1part, this%enlarge_1part(1))
102
103 ! not always the real origin if the box is shifted, no?
104 ! which happens to be my case...
105 ! only important for printout, so it is ok
106 safe_allocate(this%origin(1:ndim1part))
107 do idir = 1,ndim1part
108 irealdir = (ikeeppart-1)*ndim1part + idir
109 !origin(idir) = (npoints(irealdir)/2)*gr%mesh%spacing(irealdir)
110 this%origin(idir) = m_zero
111 end do
112
113 pop_sub(modelmb_1part_init)
114 end subroutine modelmb_1part_init
115
117 subroutine modelmb_1part_end(this)
118 type(modelmb_1part_t), intent(inout) :: this
119
120 push_sub(modelmb_1part_end)
121
122 safe_deallocate_a(this%origin)
123 safe_deallocate_a(this%enlarge_1part)
124 safe_deallocate_a(this%nr_1part)
125 safe_deallocate_a(this%ll)
126 safe_deallocate_a(this%h_1part)
127
128 call hypercube_end(this%hypercube_1part)
129
130 pop_sub(modelmb_1part_end)
131 end subroutine modelmb_1part_end
132
133
135
137!! Local Variables:
138!! mode: f90
139!! coding: utf-8
140!! End:
real(real64), parameter, public m_zero
Definition: global.F90:187
real(real64), parameter, public m_one
Definition: global.F90:188
subroutine, public hypercube_init(this, ndim, nr, enlarge)
Definition: hypercube.F90:141
subroutine, public hypercube_end(this)
Definition: hypercube.F90:186
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
general module for modelmb particles (eg 4 electrons in 1D equiv to 1 in 4D). Also calculate differen...
subroutine, public modelmb_1part_end(this)
subroutine, public modelmb_1part_init(this, mesh, ikeeppart, ndim1part)
container type for the position and dimensions for 1 particle