Changeset 20266
- Timestamp:
- 02/27/16 19:14:19 (9 years ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
issm/trunk-jpl/src/m/classes/clusters/pfe.m
r20228 r20266 42 42 disp(sprintf(' name: %s',cluster.name)); 43 43 disp(sprintf(' login: %s',cluster.login)); 44 disp(sprintf(' modules: %s',strjoin(cluster.modules,', '))); 44 modules=''; for i=1:length(cluster.modules), modules=[modules cluster.modules{i} ',']; end; modules=modules(1:end-1); 45 disp(sprintf(' modules: %s',modules)); 45 46 disp(sprintf(' port: %i',cluster.port)); 46 47 disp(sprintf(' numnodes: %i',cluster.numnodes)); … … 190 191 end 191 192 end %}}} 193 function BuildQueueScriptMultipleModels(cluster,dirname,modelname,solution,dirnames,modelnames,nps) % {{{ 194 195 %some checks: 196 if isempty(modelname), error('BuildQueueScriptMultipleModels error message: need a non empty model name!');end 197 198 %what is the executable being called? 199 executable='issm_slr.exe'; 200 201 if ispc(), error('BuildQueueScriptMultipleModels not support yet on windows machines');end; 202 203 %write queuing script 204 fid=fopen([modelname '.queue'],'w'); 205 206 fprintf(fid,'#PBS -S /bin/bash\n'); 207 fprintf(fid,'#PBS -N %s\n',modelname); 208 fprintf(fid,'#PBS -l select=%i:ncpus=%i:model=%s\n',cluster.numnodes,cluster.cpuspernode,cluster.processor); 209 fprintf(fid,'#PBS -l walltime=%i\n',cluster.time*60); %walltime is in seconds. 210 fprintf(fid,'#PBS -q %s \n',cluster.queue); 211 fprintf(fid,'#PBS -W group_list=%s\n',cluster.grouplist); 212 fprintf(fid,'#PBS -m e\n'); 213 fprintf(fid,'#PBS -o %s.outlog \n',[cluster.executionpath '/' dirname '/' modelname]); 214 fprintf(fid,'#PBS -e %s.errlog \n\n',[cluster.executionpath '/' dirname '/' modelname]); 215 fprintf(fid,'. /usr/share/modules/init/bash\n\n'); 216 fprintf(fid,'module load comp-intel/2015.0.090\n'); 217 fprintf(fid,'module load mpi-sgi/mpt.2.12r16\n'); 218 fprintf(fid,'export PATH="$PATH:."\n\n'); 219 fprintf(fid,'export MPI_GROUP_MAX=64\n\n'); 220 fprintf(fid,'export ISSM_DIR="%s/../"\n',cluster.codepath); %FIXME 221 fprintf(fid,'source $ISSM_DIR/etc/environment.sh\n'); %FIXME 222 fprintf(fid,'cd %s/%s/\n\n',cluster.executionpath,dirname); 223 224 %number of cpus: 225 mpistring=sprintf('mpiexec -np %i ',cluster.numnodes*cluster.cpuspernode); 226 227 %executable: 228 mpistring=[mpistring sprintf('%s/%s ',cluster.codepath,executable)]; 229 230 %solution name: 231 mpistring=[mpistring sprintf('%s ',EnumToString(solution))]; 232 233 %execution directory and model name: 234 mpistring=[mpistring sprintf('%s/%s %s',cluster.executionpath,dirname,modelname)]; 235 236 %inform main executable of how many icecaps, glaciers and earth models are being run: 237 mpistring=[mpistring sprintf(' %i ',length(dirnames))]; 238 239 %icecaps, glaciers and earth location, names and number of processors associated: 240 for i=1:length(dirnames), 241 mpistring=[mpistring sprintf(' %s/%s %s %i ',cluster.executionpath,dirnames{i},modelnames{i},nps{i})]; 242 end 243 244 %write this long string to disk: 245 fprintf(fid,mpistring); 246 fclose(fid); 247 248 if cluster.interactive, 249 fid=fopen([modelname '.run'],'w'); 250 251 %number of cpus: 252 mpistring=sprintf('mpiexec -np %i ',cluster.numnodes*cluster.cpuspernode); 253 254 %executable: 255 mpistring=[mpistring sprintf('%s/%s ',cluster.codepath,executable)]; 256 257 %solution name: 258 mpistring=[mpistring sprintf('%s ',EnumToString(solution))]; 259 260 %execution directory and model name: 261 mpistring=[mpistring sprintf('%s/%s %s',cluster.executionpath,dirname,modelname)]; 262 263 %inform main executable of how many icecaps, glaciers and earth models are being run: 264 mpistring=[mpistring sprintf(' %i ',length(dirnames))]; 265 266 %icecaps, glaciers and earth location, names and number of processors associated: 267 for i=1:length(dirnames), 268 mpistring=[mpistring sprintf(' %s/Interactive%i %s %i ',cluster.executionpath,cluster.interactive,modelnames{i},nps{i})]; 269 end 270 271 %write this long string to disk: 272 fprintf(fid,mpistring); 273 fclose(fid); 274 275 fid=fopen([modelname '.errlog'],'w'); 276 fclose(fid); 277 fid=fopen([modelname '.outlog'],'w'); 278 fclose(fid); 279 end 280 end 281 %}}} 192 282 function BuildKrigingQueueScript(cluster,modelname,solution,io_gather,isvalgrind,isgprof) % {{{ 193 283
Note:
See TracChangeset
for help on using the changeset viewer.