HelloSodiumChlorideInFortran.f90 15.9 KB
Newer Older
Michael Sherman's avatar
Michael Sherman committed
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
! -----------------------------------------------------------------------------
!     OpenMM(tm) HelloSodiumChloride example in Fortran 95 (June 2009)
! ------------------------------------------------------------------------------
! This is a complete, self-contained "hello world" example demonstrating 
! GPU-accelerated constant temperature simulation of a very simple system with
! just nonbonded forces, consisting of several sodium (Na+) and chloride (Cl-) 
! ions in implicit solvent. A multi-frame PDB file is written to stdout which 
! can be read by VMD or other visualization tool to produce an animation of the 
! resulting trajectory.
!
! Pay particular attention to the handling of units in this example. Incorrect
! handling of units is a very common error; this example shows how you can
! continue to work with Amber-style units like Angstroms, kCals, and van der
! Waals radii while correctly communicating with OpenMM in nm, kJ, and sigma.
!
! This example is written entirely in Fortran 95, using a Fortran interface 
! module which is NOT official parts of the OpenMM distribution.
! ------------------------------------------------------------------------------
19
20

!-------------------------------------------------------------------------------
Michael Sherman's avatar
Michael Sherman committed
21
!                ATOM, FORCE FIELD, AND SIMULATION PARAMETERS
22
!-------------------------------------------------------------------------------
Michael Sherman's avatar
Michael Sherman committed
23
24
! We'll define this module as a simplified example of the kinds of data
! structures that may already be in an MD program that is to be converted
25
26
! to use OpenMM. Note that we're using data in Angstrom and kcal units; we'll
! show how to safely convert to and from OpenMM's internal units as we go.
Michael Sherman's avatar
Michael Sherman committed
27
MODULE MyAtomInfo
28
29
    ! Simulation parameters
    ! ---------------------
Michael Sherman's avatar
Michael Sherman committed
30
31
32
33
34
35
36
37
38
39
40
41
42
    real*8 Temperature, FrictionInPerPs, SolventDielectric, SoluteDielectric
    parameter(Temperature        = 300) !Kelvins
    parameter(FrictionInPerPs    = 91)  !collisions per picosecond
    parameter(SolventDielectric  = 80)  !typical for water
    parameter(SoluteDielectric   = 2)   !typical for protein
    
    real*8 StepSizeInFs, ReportIntervalInFs, SimulationTimeInPs
    parameter(StepSizeInFs       = 2)   !integration step size (fs)
    parameter(ReportIntervalInFs = 50)  !how often for PDB frame (fs)
    parameter(SimulationTimeInPs = 100) !total simulation time (ps)
    
    ! Currently energy calculation is not available in the GPU kernels so
    ! asking for it requires slow Reference Platform computation at 
43
    ! reporting intervals. If you have a big system you'll want this off.
Michael Sherman's avatar
Michael Sherman committed
44
    logical, parameter :: WantEnergy = .true.
45

46
47
    ! Atom and force field information
    ! --------------------------------
Michael Sherman's avatar
Michael Sherman committed
48
49
50
51
52
53
54
55
56
    type Atom
        character*4       pdb
        real*8            mass, charge, vdwRadiusInAng, vdwEnergyInKcal
        real*8            gbsaRadiusInAng, gbsaScaleFactor
        real*8            initPosInAng(3)
        real*8            posInAng(3) ! leave room for runtime state info
    end type
    integer, parameter :: NumAtoms = 6
    type (Atom) :: atoms(NumAtoms) = (/ &
57
    ! pdb   mass charge vdwRad vdwEnergy  gbRad gbScale   initPos      runtime
Michael Sherman's avatar
Michael Sherman committed
58
59
60
61
62
63
64
65
Atom(' NA ',22.99,  1,  1.8680, 0.00277,  1.992,  0.8, (/ 8, 0,  0/), (/0,0,0/)),&
Atom(' CL ',35.45, -1,  2.4700, 0.1000,   1.735,  0.8, (/-8, 0,  0/), (/0,0,0/)),&
Atom(' NA ',22.99,  1,  1.8680, 0.00277,  1.992,  0.8, (/ 0, 9,  0/), (/0,0,0/)),&
Atom(' CL ',35.45, -1,  2.4700, 0.1000,   1.735,  0.8, (/ 0,-9,  0/), (/0,0,0/)),&
Atom(' NA ',22.99,  1,  1.8680, 0.00277,  1.992,  0.8, (/ 0, 0,-10/), (/0,0,0/)),&
Atom(' CL ',35.45, -1,  2.4700, 0.1000,   1.735,  0.8, (/ 0, 0, 10/), (/0,0,0/)) &
    /)
END MODULE
66
67

!-------------------------------------------------------------------------------
Michael Sherman's avatar
Michael Sherman committed
68
!                                MAIN PROGRAM
69
!-------------------------------------------------------------------------------
70
71
72
73
74
75
76
77
78
79
! This makes use of four subroutines that encapsulate all the OpenMM calls:
!     myInitializeOpenMM
!     myStepWithOpenMM
!     myGetOpenMMState
!     myTerminateOpenMM
! and one minimalist PDB file writer that has nothing to do with OpenMM:
!     myWritePDBFrame
! All of these subroutines can be found later in this file. For use in a real
! MD code you would need to write your own interface routines along these lines.
! Note that the main program does NOT include the OpenMM module.
Michael Sherman's avatar
Michael Sherman committed
80
PROGRAM HelloSodiumChloride
81
    use MyAtomInfo
82

83
84
85
86
87
88
89
    ! Calculate the number of PDB frames we want to write out and how
    ! many steps to take on the GPU in between.
    integer NumReports, NumSilentSteps
    parameter(NumReports = (SimulationTimeInPs*1000 / ReportIntervalInFs + 0.5))
    parameter(NumSilentSteps = (ReportIntervalInFs / StepSizeInFs + 0.5))
    
    character*10 platformName; real*8 timeInPs, energyInKcal; integer frame
90

91
92
93
94
95
    ! This is an opaque handle to a container that holds the OpenMM runtime
    ! objects. You can use any type for this purpose as long as it is 
    ! big enough to hold a pointer. (If you use a pointer type you'll have
    ! to declare the subroutine interfaces before calling them.)
    integer*8 ommHandle
96

97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
    ! Set up OpenMM data structures; returns platform name and handle.
    call myInitializeOpenMM(ommHandle, platformName)
 
    ! Run the simulation:
    !  (1) Write the first line of the PDB file and the initial configuration.
    !  (2) Run silently entirely within OpenMM between reporting intervals.
    !  (3) Write a PDB frame when the time comes.
    print "('REMARK  Using OpenMM platform ', A)", platformName
    call myGetOpenMMState(ommHandle, timeInPs, energyInKcal)
    call myWritePDBFrame(0, timeInPs, energyInKcal)
    
    do frame = 1, NumReports
        call myStepWithOpenMM(ommHandle, NumSilentSteps)
        call myGetOpenMMState(ommHandle, timeInPs, energyInKcal)
        call myWritePDBFrame(frame, timeInPs, energyInKcal)
    end do
    
    ! Clean up OpenMM data structures.
    call myTerminateOpenMM(ommHandle)
Michael Sherman's avatar
Michael Sherman committed
116
END PROGRAM
117
118
119
120

!-------------------------------------------------------------------------------
!                                PDB FILE WRITER
!-------------------------------------------------------------------------------
121
122
123
! Given state data which was written into the atoms array of the MyAtomInfo
! module, output a single frame (pdb "model") of the trajectory. This has
! nothing to do with OpenMM.
Michael Sherman's avatar
Michael Sherman committed
124
SUBROUTINE myWritePDBFrame(frameNum, timeInPs, energyInKcal)
125
126
127
    use MyAtomInfo; implicit none
    integer frameNum; real*8 timeInPs, energyInKcal

Michael Sherman's avatar
Michael Sherman committed
128
129
130
131
132
133
134
135
136
137
138
139
140
141
    integer n
    print "('MODEL',5X,I0)", frameNum
    print "('REMARK 250 time=', F0.3, ' picoseconds; Energy=', F0.3, ' kcal/mole')", &
        timeInPs, energyInKcal
    do n = 1,NumAtoms
        print "('ATOM  ', I5, ' ', A4, ' SLT     1    ', 3F8.3, '  1.00  0.00')",    &
            n, atoms(n)%pdb, atoms(n)%posInAng
    end do
    print "('ENDMDL')"
END SUBROUTINE 

!-------------------------------------------------------------------------------
!                            OpenMM-USING CODE
!-------------------------------------------------------------------------------
142
143
144
145
146
147
148
! The OpenMM Fortran interface module is included only at this point and below. 
! Normally these subroutines would be in a separate compilation module; we're 
! including them here for simplicity. We suggest that you write them in C++ 
! if possible, using extern "C" functions to make then callable from your
! Fortran main program. (See the C++ version of this example program for 
! an implementation of very similar routines.) However, these routines are 
! reimplemented entirely in Fortran 95 below in case you prefer.
149

150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167

! ------------------------------------------------------------------------------
!                      INITIALIZE OpenMM DATA STRUCTURES
! ------------------------------------------------------------------------------
! We take these actions here:
! (1) Load any available OpenMM plugins, e.g. Cuda and Brook.
! (2) Fill the OpenMM::System with the force field parameters we want to
!     use and the particular set of atoms to be simulated.
! (3) Create an Integrator and a Context associating the Integrator with
!     the System.
! (4) Select the OpenMM platform to be used.
! (5) Allocate a RuntimeObjects container to hang on to the System, 
!     Integrator, and Context.
! (6) Return an opaque handle to the RuntimeObjects and the name of the 
!     Platform in use.
!
! Note that this routine must understand the calling MD code's molecule and
! force field data structures so will need to be customized for each MD code.
Michael Sherman's avatar
Michael Sherman committed
168

169
170
171
172
173
174
175
SUBROUTINE myInitializeOpenMM(ommHandle, platformName)
    use OpenMM; use MyAtomInfo; implicit none
    integer*8,    intent(out) :: ommHandle
    character*10, intent(out) :: platformName

    ! This is the actual type of the opaque handle.
    type (OpenMM_RuntimeObjects) omm
176
177
178
179
180
    ! These are the objects we'll create here and store in the
    ! RuntimeObjects container for later access.
    type (OpenMM_System)             system
    type (OpenMM_LangevinIntegrator) langevin
    type (OpenMM_Context)            context
Michael Sherman's avatar
Michael Sherman committed
181
182

    ! These are temporary OpenMM objects used and discarded here.
183
184
185
    type(OpenMM_Vec3Array)      initialPosInNm
    type(OpenMM_NonbondedForce) nonbond
    type(OpenMM_GBSAOBCForce)   gbsa
186
    integer*4 n, ix
Michael Sherman's avatar
Michael Sherman committed
187
188

    type(OpenMM_String) dir
189
190
191

    ! Create a string, get the name of the default plugins directory, 
    ! and then load all the plugins found there.
Michael Sherman's avatar
Michael Sherman committed
192
193
194
    call OpenMM_String_create(dir, '')
    call OpenMM_Platform_getDefaultPluginsDirectory(dir)
    call OpenMM_Platform_loadPluginsFromDirectory(dir)
195
    call OpenMM_String_destroy(dir)
Michael Sherman's avatar
Michael Sherman committed
196
    
197
198
199
    ! Create a System and Force objects for it. The System will take
    ! over ownership of the Forces; don't destroy them yourself.
    call OpenMM_System_create(system)
Michael Sherman's avatar
Michael Sherman committed
200
201
202
203
204
    call OpenMM_NonbondedForce_create(nonbond)
    call OpenMM_GBSAOBCForce_create(gbsa)
    
    ! Convert specific force types to generic OpenMM_Force so that we can
    ! add them to the OpenMM_System.
205
206
    ix = OpenMM_System_addForce(system, transfer(nonbond, OpenMM_Force()))
    ix = OpenMM_System_addForce(system, transfer(gbsa,    OpenMM_Force()))
Michael Sherman's avatar
Michael Sherman committed
207
208
209
210
211

    ! Specify dielectrics for GBSA implicit solvation.
    call OpenMM_GBSAOBCForce_setSolventDielectric(gbsa, SolventDielectric)
    call OpenMM_GBSAOBCForce_setSoluteDielectric(gbsa, SoluteDielectric)
    
212
213
214
215
216
217
    ! Specify the atoms and their properties:
    !  (1) System needs to know the masses.
    !  (2) NonbondedForce needs charges,van der Waals properties (in MD units!).
    !  (3) GBSA needs charge, radius, and scale factor.
    !  (4) Collect default positions for initializing the simulation later.
    call OpenMM_Vec3Array_create(initialPosInNm, NumAtoms)
Michael Sherman's avatar
Michael Sherman committed
218
    do n=1,NumAtoms
Michael Sherman's avatar
Michael Sherman committed
219
        ix = OpenMM_System_addParticle(system, atoms(n)%mass)
Michael Sherman's avatar
Michael Sherman committed
220
    
221
        ix = OpenMM_NonbondedForce_addParticle(nonbond,                 &
Michael Sherman's avatar
Michael Sherman committed
222
223
224
225
226
                atoms(n)%charge,                                        &
                atoms(n)%vdwRadiusInAng  * OpenMM_NmPerAngstrom         &
                                         * OpenMM_SigmaPerVdwRadius,    &
                atoms(n)%vdwEnergyInKcal * OpenMM_KJPerKcal)
    
227
        ix = OpenMM_GBSAOBCForce_addParticle(gbsa,                      &
Michael Sherman's avatar
Michael Sherman committed
228
229
230
231
                atoms(n)%charge,                                        &
                atoms(n)%gbsaRadiusInAng * OpenMM_NmPerAngstrom,        &
                atoms(n)%gbsaScaleFactor)
    
232
233
234
235
        ! Sets initPos(n) = atoms(n)%initPos * nm/Angstrom.
        call OpenMM_Vec3Array_setScaled(initialPosInNm, n,              &
                                        atoms(n)%initPosInAng,          &
                                        OpenMM_NmPerAngstrom)
Michael Sherman's avatar
Michael Sherman committed
236
237
    end do
    
238
239
240
241
242
    ! Choose an Integrator for advancing time, and a Context connecting the
    ! System with the Integrator for simulation. Let the Context choose the
    ! best available Platform. Initialize the configuration from the default
    ! positions we collected above. Initial velocities will be zero but could
    ! have been set here.
Michael Sherman's avatar
Michael Sherman committed
243
244
245
    call OpenMM_LangevinIntegrator_create(langevin,                     &
                                          Temperature, FrictionInPerPs, &
                                          StepSizeInFs * OpenMM_PsPerFs)
246
247
248

    ! Convert LangevinIntegrator to generic Integrator type for this call.
    call OpenMM_Context_create(context, system,                         &
249
                               transfer(langevin, OpenMM_Integrator()))
250
251
252
253
    call OpenMM_Context_setPositions(context, initialPosInNm)

    ! Get the platform name to return.
    call OpenMM_Context_getPlatformName(context, platformName)
Michael Sherman's avatar
Michael Sherman committed
254
    
255
256
257
258
    ! Put the System, Integrator, and Context in the RuntimeObjects
    ! container and return an opaque reference to the container.
    call OpenMM_RuntimeObjects_create(omm)
    call OpenMM_RuntimeObjects_setSystem(omm, system)
259
    call OpenMM_RuntimeObjects_setIntegrator(omm, &
260
                               transfer(langevin, OpenMM_Integrator()))
261
262
    call OpenMM_RuntimeObjects_setContext(omm, context)
    ommHandle = transfer(omm, ommHandle)
Michael Sherman's avatar
Michael Sherman committed
263
264
265
266
267
268
END SUBROUTINE


!-------------------------------------------------------------------------------
!                     COPY STATE BACK TO CPU FROM OpenMM
!-------------------------------------------------------------------------------
269
270
271
272
273
274
SUBROUTINE myGetOpenMMState(ommHandle, timeInPs, energyInKcal)
    use OpenMM; use MyAtomInfo; implicit none
    integer*8, intent(in) :: ommHandle
    real*8, intent(out)   :: timeInPs, energyInKcal
    
    type (OpenMM_State)     state
275
    type (OpenMM_Vec3Array) posArrayInNm
276
    integer                 infoMask, n
277

278
279
    type (OpenMM_RuntimeObjects)  omm
    type (OpenMM_Context)         context
280

281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
    omm = transfer(ommHandle, omm)
    call OpenMM_RuntimeObjects_getContext(omm, context)
    
    infoMask = OpenMM_State_Positions
    if (WantEnergy) then
        infoMask = infoMask + OpenMM_State_Velocities ! for KE (cheap)
        infoMask = infoMask + OpenMM_State_Energy     ! for PE (very expensive)
    end if
    ! Forces are also available (and cheap).
    
    ! Don't forget to destroy this State when you're done with it.
    call OpenMM_Context_createState(context, infoMask, state) 
    timeInPs = OpenMM_State_getTime(state) ! OpenMM time is in ps already.
    
    ! Positions are maintained as a Vec3Array inside the State. This will give
    ! us access, but don't destroy it yourself -- it will go away with the State.
297
    call OpenMM_State_getPositions(state, posArrayInNm)
298
    do n = 1, NumAtoms
299
300
301
        ! Sets atoms(n)%pos = posArray(n) * Angstroms/nm.
        call OpenMM_Vec3Array_getScaled(posArrayInNm, n, OpenMM_AngstromsPerNm, &
                                        atoms(n)%posInAng)
302
303
304
305
306
307
308
309
310
311
312
313
    end do
    
    energyInKcal = 0
    if (WantEnergy) then
        energyInKcal = (  OpenMM_State_getPotentialEnergy(state)   &
                        + OpenMM_State_getKineticEnergy(state))    &
                       * OpenMM_KcalPerKJ
    end if
    
    ! Clean up the State memory
    call OpenMM_State_destroy(state)
END SUBROUTINE
314
315


316
317
318
319
320
321
322
!-------------------------------------------------------------------------------
!                     TAKE MULTIPLE STEPS USING OpenMM
!-------------------------------------------------------------------------------
SUBROUTINE myStepWithOpenMM(ommHandle, numSteps)
    use OpenMM; implicit none
    integer*8, intent(in) :: ommHandle
    integer,   intent(in) :: numSteps
Michael Sherman's avatar
Michael Sherman committed
323

324
325
326
327
    type (OpenMM_RuntimeObjects) omm
    type (OpenMM_Integrator)     integrator
    omm = transfer(ommHandle, omm)
    call OpenMM_RuntimeObjects_getIntegrator(omm, integrator)
Michael Sherman's avatar
Michael Sherman committed
328
    
329
    call OpenMM_Integrator_step(integrator, numSteps)
Michael Sherman's avatar
Michael Sherman committed
330
331
END SUBROUTINE

332

333
334
335
336
337
338
339
340
341
!-------------------------------------------------------------------------------
!                      DEALLOCATE ALL OpenMM OBJECTS
!-------------------------------------------------------------------------------
SUBROUTINE myTerminateOpenMM(ommHandle)
    use OpenMM; implicit none
    integer*8, intent(in) :: ommHandle

    type (OpenMM_RuntimeObjects) omm
    omm = transfer(ommHandle, omm)
342

343
    call OpenMM_RuntimeObjects_destroy(omm)
Michael Sherman's avatar
Michael Sherman committed
344
END SUBROUTINE