-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathLipidSys.sh
executable file
·235 lines (211 loc) · 6.54 KB
/
LipidSys.sh
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
#!/bin/bash
set -e
#export AMBERHOME=/home/thomas/amber18
#source $AMBERHOME/amber.sh
#Unpacked structure file:
#pdb=unpacked.pdb
## packmol-memgen --pdb ${pdb} --lipids DOPE:DOPG --ratio 3:1
pdbfilename=GakA.pdb
runname=GakA_mem
# Generate the lipid bilayer system:
packmol-memgen --pdb ${pdbfilename} --lipids DOPE:DOPG --ratio 3:1 --overwrite
## Calculate box size:
dims=$(Rscript SysSize.R bilayer_${pdbfilename})
## Generate AMBER input topology (prmtop) and input coordinates (inpcrd)
cat << EOF > system.in
source leaprc.protein.ff19SB
loadamberparams frcmod.ions1lm_126_iod_opc
source leaprc.water.opc
source leaprc.lipid17
system = loadpdb bilayer_${pdbfilename}
set system box { ${dims} }
saveamberparm system ${runname}.prmtop ${runname}.inpcrd
quit
EOF
tleap -f system.in
#a=$(grep -n -m 1 "TER" bilayer_${pdbfilename} |cut -f1 -d:)
#atomcount=$(expr $a - 3)
#echo $atomcount
## minimization
cat <<EOF > min.in
Lipid minimization $(date)
&cntrl
imin=1, ! Minimize the initial structure
maxcyc=10000, ! Maximum number of cycles for minimization
ncyc=5000, ! Switch from steepest descent to conjugate gradient minimization after ncyc cycles
ntb=1, ! Constant volume
ntp=0, ! No pressure scaling
ntf=1, ! Complete force evaluation
ntc=1, ! No SHAKE
ntpr=50, ! Print to mdout every ntpr steps
ntwr=2000, ! Write a restart file every ntwr steps
cut=10.0, ! Nonbonded cutoff in Angstroms
/
EOF
date +"%T"
echo 'Minimization'
mpirun -np 4 pmemd.MPI -O -i min.in -o minmdout -p ${runname}.prmtop -c ${runname}.inpcrd -r min.rst
## heat.in
cat <<EOF > heat1.in
Lipid heating $(date)
&cntrl
imin=0, ! Molecular dynamics
ntx=1, ! Positions read formatted with no initial velocities
irest=0, ! No restart
ntc=2, ! SHAKE on for bonds with hydrogen
ntf=2, ! No force evaluation for bonds with hydrogen
tol=0.0000001, ! SHAKE tolerance
nstlim=2500, ! Number of MD steps
ntt=3, ! Langevin thermostat
gamma_ln=1.0, ! Collision frequency for Langevin thermostat
ntr=1, ! Restrain atoms using a harmonic potential
restraint_wt=10 ! degree of restraint
restraintmask=':OL,PC,PGR', ! Restrain all except oleoyl lipid tails / acyl chains
ig=-1, ! Random seed for Langevin thermostat
ntpr=100,
ntwr=10000,
ntwx=100, ! Write to trajectory file every ntwx steps
dt=0.002, ! Timestep (ps)
nmropt=1, ! NMR restraints will be read (See TEMP0 control below)
ntb=1,
ntp=0,
cut=10.0,
ioutfm=1, ! Write a binary (netcdf) trajectory
ntxo=2, ! Write binary restart files
/
&wt
type='TEMP0', ! Varies the target temperature TEMP0
istep1=0, ! Initial step
istep2=2500, ! Final step
value1=0.0, ! Initial temp0 (K)
value2=100.0 / ! final temp0 (K)
&wt type='END' / ! End of varying conditions
Hold lipid fixed
RES :OL,PC,PGR
END
END ! End GROUP input
EOF
date +"%T"
echo 'Heating 1'
pmemd.cuda -O -i heat1.in -o heatmdout -p ${runname}.prmtop -c min.rst -r heat.rst -ref min.rst -x heat.nc
## heat2.in
cat <<EOF > heat2.in
heat2 ${date}
Lipid 128 heating 303K
&cntrl
imin=0,
ntx=5, ! Positions and velocities read formatted
irest=1, ! Restart calculation
ntc=2,
ntf=2,
tol=0.0000001,
nstlim=50000, ! Number of MD steps
ntt=3,
gamma_ln=1.0,
ntr=1,
restraint_wt=10 ! degree of restraint
restraintmask=':OL,PC,PGR', ! Restrain all except oleoyl lipid tails / acyl chains
ig=-1,
ntpr=100,
ntwr=10000,
ntwx=100,
dt=0.002,
nmropt=1,
ntb=2, ! Constant pressure periodic boundary conditions
ntp=2, ! Anisotropic pressure coupling
taup=2.0, ! Pressure relaxation time (ps)
cut=10.0,
ioutfm=1,
ntxo=2,
/
&wt
type='TEMP0',
istep1=0,
istep2=50000,
value1=100.0,
value2=303.0 /
&wt type='END' /
Hold lipids fixed
RES :OL,PC,PGR
END
END
EOF
date +"%T"
echo 'Heating 2'
pmemd.cuda -O -i heat2.in -o heat2mdout -p ${runname}.prmtop -c heat.rst -r heat2.rst -ref heat.rst -x heat2.nc
cat <<EOF > equil.in
Boundary dimensions equilibriation
&cntrl
imin=0,
ntx=5,
irest=1,
ntc=2,
ntf=2,
tol=0.0000001,
nstlim=250000,
ntt=3,
gamma_ln=1.0,
temp0=303.0,
ntpr=5000,
ntwr=5000,
ntwx=5000,
dt=0.002,
ig=-1,
ntb=2,
ntp=2,
cut=10.0,
ioutfm=1,
ntxo=2,
/
/
&ewald
skinnb=5, ! Increase skinnb to avoid skinnb errors
/
EOF
pmemd.cuda -O -i equil.in -o equil1mdout -p ${runname}.prmtop -c heat2.rst -r equil1.rst -x equil1.nc
# Repeat equilibriation
skinnb=5
j=2
while [ $j -lt 11 ] && [ $skinnb -lt 20 ]
do
date +"%T"
echo 'Equilibriation run' $j
pmemd.cuda -O -i equil.in -o equil"$j"mdout -p ${runname}.prmtop -c equil"$(($j-1))".rst -r equil"$j".rst -x equil"$j".nc
if [ $? -eq 0 ]
then
j=$(expr ${j} + 1)
else
grep -rl skinnb=${skinnb} equil.in | xargs sed -i 's/skinnb=${skinnb}/skinnb="$(($j+1))"/g'
skinnb=$(expr ${skinnb} + 1)
fi
done
echo "Equilibriation finished with a skinnb value of $skinnb."
cat <<EOF > prod.in
Lipid production 303K 125ns
&cntrl
imin=0, ! Molecular dynamics
ntx=5, ! Positions and velocities read formatted
irest=1, ! Restart calculation
ntc=2, ! SHAKE on for bonds with hydrogen
ntf=2, ! No force evaluation for bonds with hydrogen
tol=0.0000001, ! SHAKE tolerance
nstlim=62500000, ! Number of MD steps
ntt=3, ! Langevin thermostat
gamma_ln=1.0, ! Collision frequency for thermostat
temp0=303.0, ! Simulation temperature (K)
ntpr=5000, ! Print to mdout every ntpr steps
ntwr=500000, ! Write a restart file every ntwr steps
ntwx=5000, ! Write to trajectory file every ntwx steps
dt=0.002, ! Timestep (ps)
ig=-1, ! Random seed for Langevin thermostat
ntb=2, ! Constant pressure periodic boundary conditions
ntp=2, ! Anisotropic pressure coupling
cut=10.0, ! Nonbonded cutoff (Angstroms)
ioutfm=1, ! Write binary NetCDF trajectory
ntxo=2, ! Write binary restart file
/
EOF
date +"%T"
echo 'Production run'
pmemd.cuda -O -i prod.in -o prod1mdout -p ${runname}.prmtop -c equil10.rst -r prod1.rst -x prod1.nc
date +"%T"