ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE_old/src/mdtools/libmdCode/force_globals.F90
Revision: 360
Committed: Tue Mar 18 16:46:47 2003 UTC (21 years, 5 months ago) by gezelter
File size: 6555 byte(s)
Log Message:
brought force_globals back from the dead to fix circular references

File Contents

# User Rev Content
1 gezelter 360 !! Fortran interface to C entry plug.
2    
3     module force_globals
4     use definitions
5     #ifdef IS_MPI
6     use mpiSimulation
7     #endif
8    
9     implicit none
10     PRIVATE
11    
12     logical, save :: force_globals_initialized = .false.
13     integer, save :: natoms
14    
15     #ifdef IS_MPI
16     real( kind = dp ), allocatable, dimension(:,:), public :: q_Row
17     real( kind = dp ), allocatable, dimension(:,:), public :: q_Col
18     real( kind = dp ), allocatable, dimension(:,:), public :: u_l_Row
19     real( kind = dp ), allocatable, dimension(:,:), public :: u_l_Col
20     real( kind = dp ), allocatable, dimension(:,:), public :: A_Row
21     real( kind = dp ), allocatable, dimension(:,:), public :: A_Col
22    
23     real( kind = dp ), allocatable, dimension(:), public :: pot_Row
24     real( kind = dp ), allocatable, dimension(:), public :: pot_Col
25     real( kind = dp ), allocatable, dimension(:), public :: pot_Temp
26     real( kind = dp ), allocatable, dimension(:,:), public :: f_Row
27     real( kind = dp ), allocatable, dimension(:,:), public :: f_Col
28     real( kind = dp ), allocatable, dimension(:,:), public :: f_Temp
29     real( kind = dp ), allocatable, dimension(:,:), public :: t_Row
30     real( kind = dp ), allocatable, dimension(:,:), public :: t_Col
31     real( kind = dp ), allocatable, dimension(:,:), public :: t_Temp
32     real( kind = dp ), allocatable, dimension(:,:), public :: rf_Row
33     real( kind = dp ), allocatable, dimension(:,:), public :: rf_Col
34     real( kind = dp ), allocatable, dimension(:,:), public :: rf_Temp
35    
36     integer, allocatable, dimension(:), public :: atid_Row
37     integer, allocatable, dimension(:), public :: atid_Col
38     #else
39     integer, allocatable, dimension(:), public :: atid
40     #endif
41     real( kind = dp ), allocatable, dimension(:,:), public :: rf
42     real(kind = dp), dimension(9), public :: tau_Temp = 0.0_dp
43     real(kind = dp), public :: virial_Temp = 0.0_dp
44    
45     public :: InitializeForceGlobals
46     public :: getNlocal
47    
48     contains
49    
50     subroutine InitializeForceGlobals(nlocal, thisStat)
51     integer, intent(out) :: thisStat
52     integer :: nrow, ncol
53     integer :: nlocal
54     integer :: ndim = 3
55     integer :: alloc_stat
56    
57     thisStat = 0
58     natoms = nlocal
59    
60     #ifdef IS_MPI
61     nrow = getNrow(plan_row)
62     ncol = getNcol(plan_col)
63     #endif
64    
65     call FreeForceGlobals()
66    
67     #ifdef IS_MPI
68    
69     allocate(q_Row(ndim,nrow),stat=alloc_stat)
70     if (alloc_stat /= 0 ) then
71     thisStat = -1
72     return
73     endif
74    
75     allocate(q_Col(ndim,ncol),stat=alloc_stat)
76     if (alloc_stat /= 0 ) then
77     thisStat = -1
78     return
79     endif
80    
81     allocate(u_l_Row(ndim,nrow),stat=alloc_stat)
82     if (alloc_stat /= 0 ) then
83     thisStat = -1
84     return
85     endif
86    
87     allocate(u_l_Col(ndim,ncol),stat=alloc_stat)
88     if (alloc_stat /= 0 ) then
89     thisStat = -1
90     return
91     endif
92    
93     allocate(A_row(9,nrow),stat=alloc_stat)
94     if (alloc_stat /= 0 ) then
95     thisStat = -1
96     return
97     endif
98    
99     allocate(A_Col(9,ncol),stat=alloc_stat)
100     if (alloc_stat /= 0 ) then
101     thisStat = -1
102     return
103     endif
104    
105     allocate(pot_row(nrow),stat=alloc_stat)
106     if (alloc_stat /= 0 ) then
107     thisStat = -1
108     return
109     endif
110    
111     allocate(pot_Col(ncol),stat=alloc_stat)
112     if (alloc_stat /= 0 ) then
113     thisStat = -1
114     return
115     endif
116    
117     allocate(pot_Temp(nlocal),stat=alloc_stat)
118     if (alloc_stat /= 0 ) then
119     thisStat = -1
120     return
121     endif
122    
123     allocate(f_Row(ndim,nrow),stat=alloc_stat)
124     if (alloc_stat /= 0 ) then
125     thisStat = -1
126     return
127     endif
128    
129     allocate(f_Col(ndim,ncol),stat=alloc_stat)
130     if (alloc_stat /= 0 ) then
131     thisStat = -1
132     return
133     endif
134    
135     allocate(f_Temp(ndim,nlocal),stat=alloc_stat)
136     if (alloc_stat /= 0 ) then
137     thisStat = -1
138     return
139     endif
140    
141     allocate(t_Row(ndim,nrow),stat=alloc_stat)
142     if (alloc_stat /= 0 ) then
143     thisStat = -1
144     return
145     endif
146    
147     allocate(t_Col(ndim,ncol),stat=alloc_stat)
148     if (alloc_stat /= 0 ) then
149     thisStat = -1
150     return
151     endif
152    
153     allocate(t_temp(ndim,nlocal),stat=alloc_stat)
154     if (alloc_stat /= 0 ) then
155     thisStat = -1
156     return
157     endif
158    
159     allocate(atid_Row(nrow),stat=alloc_stat)
160     if (alloc_stat /= 0 ) then
161     thisStat = -1
162     return
163     endif
164    
165     allocate(atid_Col(ncol),stat=alloc_stat)
166     if (alloc_stat /= 0 ) then
167     thisStat = -1
168     return
169     endif
170    
171     allocate(rf_Row(ndim,nrow),stat=alloc_stat)
172     if (alloc_stat /= 0 ) then
173     thisStat = -1
174     return
175     endif
176    
177     allocate(rf_Col(ndim,ncol),stat=alloc_stat)
178     if (alloc_stat /= 0 ) then
179     thisStat = -1
180     return
181     endif
182    
183     allocate(rf_Temp(ndim,nlocal),stat=alloc_stat)
184     if (alloc_stat /= 0 ) then
185     thisStat = -1
186     return
187     endif
188    
189     #else
190    
191     allocate(atid(nlocal),stat=alloc_stat)
192     if (alloc_stat /= 0 ) then
193     thisStat = -1
194     return
195     end if
196    
197     #endif
198    
199     allocate(rf(ndim,nlocal),stat=alloc_stat)
200     if (alloc_stat /= 0 ) then
201     thisStat = -1
202     return
203     endif
204    
205     force_globals_initialized = .true.
206    
207     end subroutine InitializeForceGlobals
208    
209     subroutine FreeForceGlobals()
210    
211     !We free in the opposite order in which we allocate in.
212    
213     if (allocated(rf)) deallocate(rf)
214     #ifdef IS_MPI
215     if (allocated(rf_Temp)) deallocate(rf_Temp)
216     if (allocated(rf_Col)) deallocate(rf_Col)
217     if (allocated(rf_Row)) deallocate(rf_Row)
218     if (allocated(atid_Col)) deallocate(atid_Col)
219     if (allocated(atid_Row)) deallocate(atid_Row)
220     if (allocated(t_Temp)) deallocate(t_Temp)
221     if (allocated(t_Col)) deallocate(t_Col)
222     if (allocated(t_Row)) deallocate(t_Row)
223     if (allocated(f_Temp)) deallocate(f_Temp)
224     if (allocated(f_Col)) deallocate(f_Col)
225     if (allocated(f_Row)) deallocate(f_Row)
226     if (allocated(pot_Temp)) deallocate(pot_Temp)
227     if (allocated(pot_Col)) deallocate(pot_Col)
228     if (allocated(pot_Row)) deallocate(pot_Row)
229     if (allocated(A_Col)) deallocate(A_Col)
230     if (allocated(A_Row)) deallocate(A_Row)
231     if (allocated(u_l_Col)) deallocate(u_l_Col)
232     if (allocated(u_l_Row)) deallocate(u_l_Row)
233     if (allocated(q_Col)) deallocate(q_Col)
234     if (allocated(q_Row)) deallocate(q_Row)
235     #else
236     if (allocated(atid)) deallocate(atid)
237     #endif
238    
239     end subroutine FreeForceGlobals
240    
241     pure function getNlocal() result(nlocal)
242     integer :: nlocal
243     nlocal = natoms
244     end function getNlocal
245    
246     end module force_globals