ddx 0.6.8
Libary for domain-decomposition methods for polarizable continuum models
ddx_cosmo.f90
Go to the documentation of this file.
1
11
14! Get ddx-operators
17implicit none
18
21
22contains
23
32subroutine cosmo_energy(constants, state, esolv, ddx_error)
33 implicit none
34 type(ddx_constants_type), intent(in) :: constants
35 type(ddx_state_type), intent(in) :: state
36 type(ddx_error_type), intent(inout) :: ddx_error
37 real(dp), intent(out) :: esolv
38 real(dp), external :: ddot
39
40 ! dummy operation on unused interface arguments
41 if (ddx_error % flag .eq. 0) continue
42
43 esolv = pt5*ddot(constants % n, state % xs, 1, state % psi, 1)
44end subroutine cosmo_energy
45
58subroutine cosmo_setup(params, constants, workspace, state, phi_cav, psi, ddx_error)
59 implicit none
60 type(ddx_params_type), intent(in) :: params
61 type(ddx_constants_type), intent(in) :: constants
62 type(ddx_workspace_type), intent(inout) :: workspace
63 type(ddx_state_type), intent(inout) :: state
64 type(ddx_error_type), intent(inout) :: ddx_error
65 real(dp), intent(in) :: phi_cav(constants % ncav)
66 real(dp), intent(in) :: psi(constants % nbasis, params % nsph)
67
68 ! dummy operation on unused interface arguments
69 if (ddx_error % flag .eq. 0) continue
70
71 call cav_to_spherical(params, constants, workspace, phi_cav, &
72 & state % phi)
73 state % phi = - state % phi
74 state % phi_cav = phi_cav
75 state % psi = psi
76 state % rhs_done = .true.
77 state % adjoint_rhs_done = .true.
78end subroutine cosmo_setup
79
89subroutine cosmo_guess(params, constants, workspace, state, ddx_error)
90 implicit none
91 type(ddx_params_type), intent(in) :: params
92 type(ddx_constants_type), intent(in) :: constants
93 type(ddx_workspace_type), intent(inout) :: workspace
94 type(ddx_state_type), intent(inout) :: state
95 type(ddx_error_type), intent(inout) :: ddx_error
96
97 ! apply the diagonal preconditioner as a guess
98 call ldm1x(params, constants, workspace, state % phi, state % xs, ddx_error)
99
100end subroutine cosmo_guess
101
111subroutine cosmo_guess_adjoint(params, constants, workspace, state, ddx_error)
112 implicit none
113 type(ddx_params_type), intent(in) :: params
114 type(ddx_constants_type), intent(in) :: constants
115 type(ddx_workspace_type), intent(inout) :: workspace
116 type(ddx_state_type), intent(inout) :: state
117 type(ddx_error_type), intent(inout) :: ddx_error
118
119 ! apply the diagonal preconditioner as a guess
120 call ldm1x(params, constants, workspace, state % psi, state % s, ddx_error)
121
122end subroutine cosmo_guess_adjoint
123
134subroutine cosmo_solve(params, constants, workspace, state, tol, ddx_error)
135 implicit none
136 type(ddx_params_type), intent(in) :: params
137 type(ddx_constants_type), intent(in) :: constants
138 type(ddx_workspace_type), intent(inout) :: workspace
139 type(ddx_state_type), intent(inout) :: state
140 real(dp), intent(in) :: tol
141 type(ddx_error_type), intent(inout) :: ddx_error
142 ! local variables
143 real(dp) :: start_time, finish_time
144
145 state % xs_niter = params % maxiter
146 start_time = omp_get_wtime()
147 call jacobi_diis(params, constants, workspace, tol, state % phi, &
148 & state % xs, state % xs_niter, state % xs_rel_diff, lx, ldm1x, &
149 & hnorm, ddx_error)
150 finish_time = omp_get_wtime()
151 state % xs_time = finish_time - start_time
152
153end subroutine cosmo_solve
154
165subroutine cosmo_solve_adjoint(params, constants, workspace, state, tol, &
166 & ddx_error)
167 implicit none
168 type(ddx_params_type), intent(in) :: params
169 type(ddx_constants_type), intent(in) :: constants
170 type(ddx_workspace_type), intent(inout) :: workspace
171 type(ddx_state_type), intent(inout) :: state
172 real(dp), intent(in) :: tol
173 type(ddx_error_type), intent(inout) :: ddx_error
174 ! local variables
175 real(dp) :: start_time, finish_time
176
177 state % s_niter = params % maxiter
178 start_time = omp_get_wtime()
179 call jacobi_diis(params, constants, workspace, tol, state % psi, &
180 & state % s, state % s_niter, state % s_rel_diff, lstarx, ldm1x, &
181 & hnorm, ddx_error)
182 finish_time = omp_get_wtime()
183 state % s_time = finish_time - start_time
184
185 state % q = state % s
186
187 call cosmo_derivative_setup(params, constants, workspace, state)
188
189end subroutine cosmo_solve_adjoint
190
203subroutine cosmo_solvation_force_terms(params, constants, workspace, &
204 & state, e_cav, force, ddx_error)
205 implicit none
206 type(ddx_params_type), intent(in) :: params
207 type(ddx_constants_type), intent(in) :: constants
208 type(ddx_workspace_type), intent(inout) :: workspace
209 type(ddx_state_type), intent(inout) :: state
210 type(ddx_error_type), intent(inout) :: ddx_error
211 real(dp), intent(in) :: e_cav(3, constants % ncav)
212 real(dp), intent(inout) :: force(3, params % nsph)
213 ! local variables
214 real(dp) :: start_time, finish_time
215 integer :: isph
216
217 ! dummy operation on unused interface arguments
218 if (ddx_error % flag .eq. 0) continue
219
220 start_time = omp_get_wtime()
221
222 force = zero
223 do isph = 1, params % nsph
224 call contract_grad_l(params, constants, isph, state % xs, &
225 & state % sgrid, workspace % tmp_vylm(:, 1), &
226 & workspace % tmp_vdylm(:, :, 1), workspace % tmp_vplm(:, 1), &
227 & workspace % tmp_vcos(:, 1), workspace % tmp_vsin(:, 1), &
228 & force(:, isph))
229 call contract_grad_u(params, constants, isph, state % sgrid, &
230 & state % phi_grid, force(:, isph))
231 end do
232
233 force = -pt5*force
234
235 call zeta_grad(params, constants, state, e_cav, force)
236
237 finish_time = omp_get_wtime()
238 state % force_time = finish_time - start_time
239
240end subroutine cosmo_solvation_force_terms
241
250subroutine cosmo_derivative_setup(params, constants, workspace, state)
251 type(ddx_params_type), intent(in) :: params
252 type(ddx_constants_type), intent(in) :: constants
253 type(ddx_workspace_type), intent(inout) :: workspace
254 type(ddx_state_type), intent(inout) :: state
255
256 real(dp), external :: ddot
257 integer :: icav, isph, igrid
258
259 ! dummy operation on unused interface arguments
260 if (allocated(workspace % tmp_pot)) continue
261
262 ! Get values of S on the grid
263 call ddeval_grid_work(constants % nbasis, params % ngrid, params % nsph, &
264 & constants % vgrid, constants % vgrid_nbasis, one, state % s, zero, &
265 & state % sgrid)
266 ! Get the values of phi on the grid
267 call ddcav_to_grid_work(params % ngrid, params % nsph, constants % ncav, &
268 & constants % icav_ia, constants % icav_ja, state % phi_cav, &
269 & state % phi_grid)
270
271 ! assemble the intermediate zeta: S weighted by U evaluated on the
272 ! exposed grid points.
273 icav = 0
274 do isph = 1, params % nsph
275 do igrid = 1, params % ngrid
276 if (constants % ui(igrid, isph) .ne. zero) then
277 icav = icav + 1
278 state % zeta(icav) = constants % wgrid(igrid) * &
279 & constants % ui(igrid, isph) * ddot(constants % nbasis, &
280 & constants % vgrid(1, igrid), 1, state % s(1, isph), 1)
281 end if
282 end do
283 end do
284
285end subroutine cosmo_derivative_setup
286
287end module ddx_cosmo
288
subroutine cosmo_guess(params, constants, workspace, state, ddx_error)
Do a guess for the primal ddCOSMO linear system.
Definition: ddx_cosmo.f90:90
subroutine cosmo_solve(params, constants, workspace, state, tol, ddx_error)
Solve the primal ddCOSMO linear system.
Definition: ddx_cosmo.f90:135
subroutine cosmo_solve_adjoint(params, constants, workspace, state, tol, ddx_error)
Solve the adjoint ddCOSMO linear system.
Definition: ddx_cosmo.f90:167
subroutine cosmo_solvation_force_terms(params, constants, workspace, state, e_cav, force, ddx_error)
Compute the solvation term of the forces (solute aspecific). This must be summed to the solute specif...
Definition: ddx_cosmo.f90:205
subroutine cosmo_setup(params, constants, workspace, state, phi_cav, psi, ddx_error)
Given the potential at the cavity points, assemble the RHS for ddCOSMO or for ddPCM.
Definition: ddx_cosmo.f90:59
subroutine cosmo_guess_adjoint(params, constants, workspace, state, ddx_error)
Do a guess for the adjoint ddCOSMO linear system.
Definition: ddx_cosmo.f90:112
subroutine cosmo_energy(constants, state, esolv, ddx_error)
Compute the ddCOSMO energy.
Definition: ddx_cosmo.f90:33
High-level subroutines for ddcosmo.
Definition: ddx_cosmo.f90:13
subroutine cosmo_derivative_setup(params, constants, workspace, state)
This routines precomputes the intermediates to be used in the evaluation of ddCOSMO analytical deriva...
Definition: ddx_cosmo.f90:251
Routines to build rhs (phi and psi)
Operators shared among ddX methods.
subroutine lx(params, constants, workspace, x, y, ddx_error)
Single layer operator matvec.
subroutine ldm1x(params, constants, workspace, x, y, ddx_error)
Diagonal preconditioning for Lx operator.
subroutine lstarx(params, constants, workspace, x, y, ddx_error)
Adjoint single layer operator matvec.