Run OpenSees Script at the CLI#

Most OpenSees analyses are run in non-interactive mode, where you execute a full input script all at once. Even if you’re in an interactive environment like Jupyter, running a script this way means the entire sequence of commands executes without waiting for further input. Jupyter Hub has an integrated editor to make things easy.

You can run therefore, run either OpenSees(Tcl) and OpenSeesPy (Python) in non-interactive batch mode at the terminal’s CLI by specifying the input script to execute. In both cases, the first argument provided after the executable is interpreted as the main script file to run.

Most OpenSees Analyses require that you work with an input script that you can edit as needed to fix errors and add features. Being able to work with this script in an interactive environment allows to have immediate feedback on your changes.

This is your starting place. Develop your script for a few test cases.

The following images show how to run OpenSees Tcl and OpenSeesPy at the terminal for both sequential and parallel analyses, for both Tcl and Python interpreters.

Sequential Execution#

Language

Command

Tcl

OpenSees inputfile.tcl

OpenSeesPy

python inputfile.py (import OpenSeesPy inside inputfile.py)

OpenSees-Tcl
1. At the Command-line

Type the following command at the command-line prompt:

 OpenSees Ex1a.Canti2D.Push.tcl

You will see the following response:

         OpenSees -- Open System For Earthquake Engineering Simulation
                 Pacific Earthquake Engineering Research Center
                        Version 3.7.1 64-Bit

      (c) Copyright 1999-2016 The Regents of the University of California
                              All Rights Reserved
  (Copyright and Disclaimer @ http://www.berkeley.edu/OpenSees/copyright.html)


Analysis-0 execution done
Analysis-1 execution done
Analysis-2 execution done
Analysis-3 execution done
Analysis-4 execution done
Analysis-5 execution done
Analysis-6 execution done
Analysis-7 execution done
ALL DONE!!!

The program exits and you are returned to the prompt:

(base) jovyan@3cd0f33abec1:~/BasicExamples$
2. From a Jupyter Notebook – shell command

Type the following command in a notebook cell:

 import os
 os.system('OpenSees Ex1a.Canti2D.Push.tcl')

You will see the following response in the output:

     OpenSees -- Open System For Earthquake Engineering Simulation
             Pacific Earthquake Engineering Research Center
                    Version 3.7.1 64-Bit

  (c) Copyright 1999-2016 The Regents of the University of California
                          All Rights Reserved
(Copyright and Disclaimer @ http://www.berkeley.edu/OpenSees/copyright.html)


Analysis-0 execution done
Analysis-1 execution done
Analysis-2 execution done
Analysis-3 execution done
Analysis-4 execution done
Analysis-5 execution done
Analysis-6 execution done
Analysis-7 execution done
ALL DONE!!!
OpenSeesPy
1. At the Command-line

Type the following command at the command-line prompt:

 python Ex1a.Canti2D.Push.py

You will see the following response:

Analysis-0 execution done
Analysis-1 execution done
Analysis-2 execution done
Analysis-3 execution done
Analysis-4 execution done
Analysis-5 execution done
Analysis-6 execution done
Analysis-7 execution done
ALL DONE!!!
Process 0 Terminating

The program exits and you are returned to the prompt:

(base) jovyan@3cd0f33abec1:~/BasicExamples$
2. From a Jupyter Notebook – shell command

Type the following command in a notebook cell:

 import os
 os.system('python Ex1a.Canti2D.Push.py')

You will see the following response:

Analysis-0 execution done
Analysis-1 execution done
Analysis-2 execution done
Analysis-3 execution done
Analysis-4 execution done
Analysis-5 execution done
Analysis-6 execution done
Analysis-7 execution done
ALL DONE!!!
Process 0 Terminating
3. Within a Jupyter Notebook

Type the following command in a notebook cell:

 %run 'Ex1a.Canti2D.Push.py'

You will see the following response:

Analysis-0 execution done
Analysis-1 execution done
Analysis-2 execution done
Analysis-3 execution done
Analysis-4 execution done
Analysis-5 execution done
Analysis-6 execution done
Analysis-7 execution done
ALL DONE!!!

The script was run in the kernel, so now you can continue with OpenSees commands or any python command. In the next cell type:

print(LColList)

You will see the following response:

[100, 120, 200, 240, 300, 360, 400, 480]

Type in the next cell:

ops.getNodeTags()

You will see the following response:

[1, 2]

Parallel Execution (N cores)#

Use mpiexec, mpirun, or ibrun to launch distributed memory jobs:

Language

Command

Tcl (MP)

mpiexec -np N OpenSeesMP inputfileMP.tcl

Tcl (SP)

mpiexec -np N OpenSeesSP inputfileSP.tcl

OpenSeesPy (MPI)

mpiexec -np N python inputfile.py *

  • (import OpenSeesPy inside. And script must use mpi4py or similar)

  • N is the number of concurrent processors

  • Use mpiexec in JupyterHub

  • Use ibrun in the HPC system, such as Stampede3

OpenSeesMP
1. At the Command-line prompt

Type the following command at the command-line prompt:

 mpiexec -np 3 OpenSeesMP Ex1a.Canti2D.Push.mp.tcl

You will see the following response:

         OpenSees -- Open System For Earthquake Engineering Simulation
                 Pacific Earthquake Engineering Research Center
                        Version 3.7.1 64-Bit

      (c) Copyright 1999-2016 The Regents of the University of California
                              All Rights Reserved
  (Copyright and Disclaimer @ http://www.berkeley.edu/OpenSees/copyright.html)


pid 1 of np=3  started
pid 0 of np=3  started
pid 2 of np=3  started
pid 2 of 3 Analysis-2 execution done
pid 1 of 3 Analysis-1 execution done
pid 0 of 3 Analysis-0 execution done
pid 2 of 3 Analysis-5 execution done
pid 2 ALL DONE!!!
Process Terminating 2
pid 1 of 3 Analysis-4 execution done
pid 0 of 3 Analysis-3 execution done
pid 1 of 3 Analysis-7 execution done
pid 1 ALL DONE!!!
Process Terminating 1
pid 0 of 3 Analysis-6 execution done
pid 0 ALL DONE!!!
Process Terminating 0

The program exits and you are returned to the prompt:

(base) jovyan@3cd0f33abec1:~/BasicExamples$
2. In a Jupyter Notebook, as a shell command

Type the following command in a notebook cell:

 import os
 os.system('mpiexec -np 3 OpenSeesMP Ex1a.Canti2D.Push.mp.tcl')

You will see the following response:

         OpenSees -- Open System For Earthquake Engineering Simulation
                 Pacific Earthquake Engineering Research Center
                        Version 3.7.1 64-Bit

      (c) Copyright 1999-2016 The Regents of the University of California
                              All Rights Reserved
  (Copyright and Disclaimer @ http://www.berkeley.edu/OpenSees/copyright.html)


pid 2 of np=3  started
pid 0 of np=3  started
pid 1 of np=3  started
pid 2 of 3 Analysis-2 execution done
pid 1 of 3 Analysis-1 execution done
pid 0 of 3 Analysis-0 execution done
pid 2 of 3 Analysis-5 execution done
pid 2 ALL DONE!!!
Process Terminating 2
pid 1 of 3 Analysis-4 execution done
pid 0 of 3 Analysis-3 execution done
pid 1 of 3 Analysis-7 execution done
pid 1 ALL DONE!!!
Process Terminating 1
pid 0 of 3 Analysis-6 execution done
pid 0 ALL DONE!!!
Process Terminating 0
OpenSeesPy Parallel
1a. At the Command-line prompt, OpenSeesPy mpi.

Type the following command at the command-line prompt:

 mpiexec -np 3 python Ex1a.Canti2D.Push.mpi.py

You will see the following response:

pid 0 of 1 started
pid 0 of 1 started
pid 0 of 1 started
pid 0 of np=1 Analysis-0 execution done
pid 0 of np=1 Analysis-0 execution done
pid 0 of np=1 Analysis-0 execution done
pid 0 of np=1 Analysis-1 execution done
pid 0 of np=1 Analysis-2 execution done
pid 0 of np=1 Analysis-1 execution done
pid 0 of np=1 Analysis-1 execution done
pid 0 of np=1 Analysis-3 execution done
pid 0 of np=1 Analysis-2 execution done
pid 0 of np=1 Analysis-2 execution done
pid 0 of np=1 Analysis-4 execution done
pid 0 of np=1 Analysis-3 execution done
pid 0 of np=1 Analysis-3 execution done
pid 0 of np=1 Analysis-5 execution done
pid 0 of np=1 Analysis-4 execution done
pid 0 of np=1 Analysis-4 execution done
pid 0 of np=1 Analysis-6 execution done
pid 0 of np=1 Analysis-5 execution done
pid 0 of np=1 Analysis-5 execution done
pid 0 of np=1 Analysis-7 execution done
pid 0 of np=1 ALL DONE!!!
pid 0 of np=1 Analysis-6 execution done
pid 0 of np=1 Analysis-6 execution done
pid 0 of np=1 Analysis-7 execution done
pid 0 of np=1 ALL DONE!!!
pid 0 of np=1 Analysis-7 execution done
pid 0 of np=1 ALL DONE!!!
Process 0 Terminating
Process 0 Terminating
Process 0 Terminating

The program exits and you are returned to the prompt:

(base) jovyan@3cd0f33abec1:~/BasicExamples$

You may run into the problem that the OpenSeesPy MPI does not work properly – it does run 3 instances of OpenSeesPy, but each instance is not aware of the others and identifies as pid=0 out of one instead of pid=0,1,or 2 out of 3. –> All three processors ran all the analyses!

1b. At the Command-line prompt, use mpi4py.

Type the following command in a notebook cell:

 mpiexec -np 3 python Ex1a.Canti2D.Push.mpi4py.py

You will see the following response:

mpi4py -- python pid 0 of 3 started
mpi4py -- python pid 1 of 3 started
mpi4py -- python pid 2 of 3 started
pid 2 of np=3 Analysis-2 execution done
pid 0 of np=3 Analysis-0 execution done
pid 1 of np=3 Analysis-1 execution done
pid 2 of np=3 Analysis-5 execution done
pid 2 of np=3 ALL DONE!!!
pid 1 of np=3 Analysis-4 execution done
pid 0 of np=3 Analysis-3 execution done
pid 1 of np=3 Analysis-7 execution done
pid 1 of np=3 ALL DONE!!!
pid 0 of np=3 Analysis-6 execution done
pid 0 of np=3 ALL DONE!!!
Process 0 Terminating
Process 0 Terminating
Process 0 Terminating

The program exits and you are returned to the prompt:

(base) jovyan@3cd0f33abec1:~/BasicExamples$

You now see the three pids out of 3, and each processor ran only the analysis it was assigned.

2a. In a Jupyter Notebook – shell command, OpenSeesPy-MPI

Type the following command in a notebook cell:

 import os
 os.system('mpiexec -np 3 python Ex1a.Canti2D.Push.mpi.py')

You will see the following response. Note that the program exits once the analysis has been run:

pid 0 of 1 started
pid 0 of 1 started
pid 0 of np=1 Analysis-0 execution done
pid 0 of 1 started
pid 0 of np=1 Analysis-0 execution done
pid 0 of np=1 Analysis-1 execution done
pid 0 of np=1 Analysis-0 execution done
pid 0 of np=1 Analysis-1 execution done
pid 0 of np=1 Analysis-2 execution done
pid 0 of np=1 Analysis-1 execution done
pid 0 of np=1 Analysis-2 execution done
pid 0 of np=1 Analysis-3 execution done
pid 0 of np=1 Analysis-2 execution done
pid 0 of np=1 Analysis-3 execution done
pid 0 of np=1 Analysis-4 execution done
pid 0 of np=1 Analysis-3 execution done
pid 0 of np=1 Analysis-4 execution done
pid 0 of np=1 Analysis-5 execution done
pid 0 of np=1 Analysis-4 execution done
pid 0 of np=1 Analysis-6 execution done
pid 0 of np=1 Analysis-5 execution done
pid 0 of np=1 Analysis-5 execution done
pid 0 of np=1 Analysis-7 execution done
pid 0 of np=1 ALL DONE!!!
pid 0 of np=1 Analysis-6 execution done
pid 0 of np=1 Analysis-6 execution done
pid 0 of np=1 Analysis-7 execution done
pid 0 of np=1 ALL DONE!!!
pid 0 of np=1 Analysis-7 execution done
pid 0 of np=1 ALL DONE!!!
Process 0 Terminating
Process 0 Terminating
Process 0 Terminating

As expected, you have the same problem as you had at the command line.

2. In a Jupyter Notebook – shell command, mpi4py

Type the following command in a notebook cell:

 import os
 os.system('mpiexec -np 3 python Ex1a.Canti2D.Push.mpi4py.py')

You will see the following response. Note that the program exits once the analysis has been run:

mpi4py -- python pid 0 of 3 started
mpi4py -- python pid 2 of 3 started
mpi4py -- python pid 1 of 3 started
pid 1 of np=3 Analysis-1 execution done
pid 0 of np=3 Analysis-0 execution done
pid 2 of np=3 Analysis-2 execution done
pid 1 of np=3 Analysis-4 execution done
pid 0 of np=3 Analysis-3 execution done
pid 2 of np=3 Analysis-5 execution done
pid 2 of np=3 ALL DONE!!!
pid 1 of np=3 Analysis-7 execution done
pid 1 of np=3 ALL DONE!!!
pid 0 of np=3 Analysis-6 execution done
pid 0 of np=3 ALL DONE!!!
Process Process 0 Terminating
0 Terminating
Process 0 Terminating

As expected, using mpi4pi has resolved the issue.

Example Files Used in this demo#

You can find these files in Community Data:

OpenSees-Tcl
Ex1a.Canti2D.Push.tcl
# OpenSees Ex1a.Canti2D.Push.tcl

############################################################
#  EXAMPLE: 
#       Ex1a.Canti2D.Push.tcl
#          for OpenSees.exe (tcl)
#  --------------------------------------------------------#
#  by: Silvia Mazzoni, 2020
#       silviamazzoni@yahoo.com
############################################################
# This file was obtained by updating the Tcl script in the original examples manual
# You can find the original Examples:
# https://opensees.berkeley.edu/wiki/index.php/Examples_Manual
# Original Examples by By Silvia Mazzoni & Frank McKenna, 2006, in Tcl
############################################################
# --------------------------------------------------------------------------------------------------
# Example 1. cantilever 2D
# static pushover analysis with gravity.
# all units are in kip, inch, second
# elasticBeamColumn ELEMENT
#			Silvia Mazzoni & Frank McKenna, 2006
#
#    ^Y
#    |
#    2       __ 
#    |         | 
#    |         | 
#    |         | 
#  (1)      36'
#    |         | 
#    |         | 
#    |         | 
#  =1=    ----  -------->X
#
#

if {[llength $argv]>0} {
    puts "Command-Line Arguments (argv): $argv"
}

set LcolList "100 120 200 240 300 360 400 480"

# ----------------------------------------------
set dataDir DataTCL;                # set up name of data directory
file mkdir $dataDir; # create data directory

set count 0;
foreach Lcol $LcolList {
    
    # SET UP ----------------------------------------------------------------------------
    wipe;						# clear opensees model
    model basic -ndm 2 -ndf 3;	# 2 dimensions, 3 dof per node

    # define GEOMETRY -------------------------------------------------------------
    # nodal coordinates:
    node 1 0 0;				# node#, X Y
    node 2 0 $Lcol
    
    # Single point constraints -- Boundary Conditions
    fix 1 1 1 1; 			# node DX DY RZ
    
    # nodal masses:
    mass 2 5.18 0. 0.;		# node#, Mx My Mz, Mass=Weight/g.
    
    # Define ELEMENTS -------------------------------------------------------------
    # define geometric transformation: performs a linear geometric transformation of beam stiffness 
    #	and resisting force from the basic system to the global-coordinate system
    geomTransf Linear 1;  		# associate a tag to transformation
    
    # connectivity: (make A very large, 10e6 times its actual value)
    element elasticBeamColumn 1 1 2 3600000000 4227 1080000 1;	
    
    # Define RECORDERS -------------------------------------------------------------
    recorder Node -file ${dataDir}/DFree_Lcol${Lcol}.out -time -node 2 -dof 1 2 3 disp;			# displacements of free nodes
    recorder Node -file ${dataDir}/DBase_Lcol${Lcol}.out -time -node 1 -dof 1 2 3 disp;			# displacements of support nodes
    recorder Node -file ${dataDir}/RBase_Lcol${Lcol}.out -time -node 1 -dof 1 2 3 reaction;		# support reaction
    recorder Element -file ${dataDir}/FCol_Lcol${Lcol}.out -time -ele 1 globalForce;			# element forces -- column
    recorder Element -file ${dataDir}/DCol_Lcol${Lcol}.out -time -ele 1 deformation;			# element deformations -- column
    
    
    # define GRAVITY -------------------------------------------------------------
    timeSeries Linear 1;
    pattern Plain 1 1 {
       load 2 0. -2000. 0.;			# node#, FX FY MZ --  superstructure-weight
    }
    constraints Plain;     			# how it handles boundary conditions
    numberer Plain;					# renumber dof's to minimize band-width (optimization), if you want to
    system BandGeneral;				# how to store and solve the system of equations in the analysis
    test NormDispIncr 1.0e-8 6 ; 	# determine if convergence has been achieved at the end of an iteration step
    algorithm Newton;				# use Newton's solution algorithm: updates tangent stiffness at every iteration
    integrator LoadControl 0.1;		# determine the next time step for an analysis, # apply gravity in 10 steps
    analysis Static					# define type of analysis static or transient
    analyze 10;						# perform gravity analysis
    loadConst -time 0.0;			# hold gravity constant and restart time
    
    # define LATERAL load -------------------------------------------------------------
    timeSeries Linear 2;
    pattern Plain 2 2 {
        load 2 2000. 0.0 0.0;		# node#, FX FY MZ -- representative lateral load at top node
    }
    
    # pushover: diplacement controlled static analysis
    integrator DisplacementControl 2 1 0.1;		# switch to displacement control, for node 11, dof 1, 0.1 increment
    analyze 1000;								# apply 100 steps of pushover analysis to a displacement of 10
    
    puts "Analysis-${count} execution done"

    incr count 1;
}

puts "ALL DONE!!!"
Ex1a.Canti2D.Push.mp.tcl
# mpiexec -np 3 OpenSeesMP Ex1a.Canti2D.Push.mp.tcl

############################################################
#  EXAMPLE: 
#       Ex1a.Canti2D.Push.tcl
#          for OpenSees.exe (tcl)
#  --------------------------------------------------------#
#  by: Silvia Mazzoni, 2020
#       silviamazzoni@yahoo.com
############################################################
# This file was obtained by updating the Tcl script in the original examples manual
# You can find the original Examples:
# https://opensees.berkeley.edu/wiki/index.php/Examples_Manual
# Original Examples by By Silvia Mazzoni & Frank McKenna, 2006, in Tcl
############################################################
# --------------------------------------------------------------------------------------------------
# Example 1. cantilever 2D
# static pushover analysis with gravity.
# all units are in kip, inch, second
# elasticBeamColumn ELEMENT
#			Silvia Mazzoni & Frank McKenna, 2006
#
#    ^Y
#    |
#    2       __ 
#    |         | 
#    |         | 
#    |         | 
#  (1)      36'
#    |         | 
#    |         | 
#    |         | 
#  =1=    ----  -------->X
#
#



set pid [getPID]
set np [getNP]
puts "pid $pid of np=$np  started"


if {[llength $argv]>0} {
    puts "pid $pid of np=$np  Command-Line Arguments (argv): $argv"
}

set LcolList "100 120 200 240 300 360 400 480"

# ----------------------------------------------
set dataDir DataTCLmp;                # set up name of data directory
file mkdir $dataDir; # create data directory

set count 0;
foreach Lcol $LcolList {
    # check if count is a multiple of pid : "is it its turn?"
    if {[expr $count % $np] == $pid} {
       
        # SET UP ----------------------------------------------------------------------------
        wipe;						# clear opensees model
        model basic -ndm 2 -ndf 3;	# 2 dimensions, 3 dof per node
        
        # define GEOMETRY -------------------------------------------------------------
        # nodal coordinates:
        node 1 0 0;				# node#, X Y
        node 2 0 $Lcol
        
        # Single point constraints -- Boundary Conditions
        fix 1 1 1 1; 			# node DX DY RZ
        
        # nodal masses:
        mass 2 5.18 0. 0.;		# node#, Mx My Mz, Mass=Weight/g.
        
        # Define ELEMENTS -------------------------------------------------------------
        # define geometric transformation: performs a linear geometric transformation of beam stiffness 
        #	and resisting force from the basic system to the global-coordinate system
        geomTransf Linear 1;  		# associate a tag to transformation
        
        # connectivity: (make A very large, 10e6 times its actual value)
        # element elasticBeamColumn $eleTag $iNode $jNode $A $E $Iz $transfTag
        element elasticBeamColumn 1 1 2 3600000000 4227 1080000 1;	
        
        # Define RECORDERS -------------------------------------------------------------
        recorder Node -file ${dataDir}/DFree_Lcol${Lcol}.out -time -node 2 -dof 1 2 3 disp;			# displacements of free nodes
        recorder Node -file ${dataDir}/DBase_Lcol${Lcol}.out -time -node 1 -dof 1 2 3 disp;			# displacements of support nodes
        recorder Node -file ${dataDir}/RBase_Lcol${Lcol}.out -time -node 1 -dof 1 2 3 reaction;		# support reaction
        recorder Element -file ${dataDir}/FCol_Lcol${Lcol}.out -time -ele 1 globalForce;			# element forces -- column
        recorder Element -file ${dataDir}/DCol_Lcol${Lcol}.out -time -ele 1 deformation;			# element deformations -- column
        
        
        # define GRAVITY -------------------------------------------------------------
        timeSeries Linear 1;
        pattern Plain 1 1 {
           load 2 0. -2000. 0.;			# node#, FX FY MZ --  superstructure-weight
        }
        constraints Plain;     			# how it handles boundary conditions
        numberer Plain;					# renumber dof's to minimize band-width (optimization), if you want to
        system BandGeneral;				# how to store and solve the system of equations in the analysis
        test NormDispIncr 1.0e-8 6 ; 	# determine if convergence has been achieved at the end of an iteration step
        algorithm Newton;				# use Newton's solution algorithm: updates tangent stiffness at every iteration
        integrator LoadControl 0.1;		# determine the next time step for an analysis, # apply gravity in 10 steps
        analysis Static					# define type of analysis static or transient
        analyze 10;						# perform gravity analysis
        loadConst -time 0.0;			# hold gravity constant and restart time
        
        # define LATERAL load -------------------------------------------------------------
        timeSeries Linear 2;
        pattern Plain 2 2 {
           load 2 2000. 0.0 0.0;		# node#, FX FY MZ -- representative lateral load at top node
        }
        
        # pushover: diplacement controlled static analysis
        integrator DisplacementControl 2 1 0.1;		# switch to displacement control, for node 11, dof 1, 0.1 increment
        analyze 1000;								# apply 100 steps of pushover analysis to a displacement of 10
        
        puts "pid $pid of $np Analysis-${count} execution done"
        
        
    }
    incr count 1;
}

puts "pid $pid ALL DONE!!!"
OpenSeesPy
Ex1a.Canti2D.Push.py
# python Ex1a.Canti2D.Push.py

############################################################
#  EXAMPLE: 
#       Ex1a.Canti2D.Push.py
#          for OpenSeesPy
#  --------------------------------------------------------#
#  by: Silvia Mazzoni, 2020
#       silviamazzoni@yahoo.com
############################################################
# This file was obtained from a conversion of the updated Tcl script
# You can find the original Examples:
# https://opensees.berkeley.edu/wiki/index.php/Examples_Manual
# Original Examples by By Silvia Mazzoni & Frank McKenna, 2006, in Tcl
# Converted to OpenSeesPy by SilviaMazzoni, 2020
############################################################
# --------------------------------------------------------------------------------------------------
# Example 1. cantilever 2D
# static pushover analysis with gravity.
# all units are in kip, inch, second
# elasticBeamColumn ELEMENT
#			Silvia Mazzoni & Frank McKenna, 2006
#
#    ^Y
#    |
#    2       __ 
#    |         | 
#    |         | 
#    |         | 
#  (1)      36'
#    |         | 
#    |         | 
#    |         | 
#  =1=    ----  -------->X
#
#

import openseespy.opensees as ops
import numpy as numpy
import matplotlib.pyplot as plt
import os
import sys

print(sys.argv)
if len(sys.argv)>1:
    print(f'Command-Line Arguments (argv): {sys.argv}')


LColList = [100,120,200,240,300,360,400,480]
#-----------------------------------------
dataDir=f'DataPY';                # set up name of data directory
os.makedirs(dataDir, exist_ok=True);    # create data directory

count = 0;
for Lcol in LColList:
    ops.wipe()

    # SET UP ----------------------------------------------------------------------------
    ops.wipe()     #  clear opensees model
    ops.model('basic','-ndm',2,'-ndf',3)     #  2 dimensions, 3 dof per node

    # define GEOMETRY -------------------------------------------------------------
    # nodal coordinates:
    ops.node(1,0,0)     #  node , X Y
    ops.node(2,0,Lcol)

    # Single point constraints -- Boundary Conditions
    ops.fix(1,1,1,1)     #  node DX DY RZ

    # nodal masses:
    ops.mass(2,5.18,0.,0.)     #  node , Mx My Mz, Mass=Weight/g.

    # Define ELEMENTS -------------------------------------------------------------
    # define geometric transformation: performs a linear geometric transformation of beam stiffness
    # and resisting force from the basic system to the global-coordinate system
    ops.geomTransf('Linear',1)     #  associate a tag to transformation

    # element elasticBeamColumn eleTag iNode jNode A E Iz transfTag
    ops.element('elasticBeamColumn',1,1,2,3600000000,4227,1080000,1)

    # Define RECORDERS -------------------------------------------------------------
    ops.recorder('Node','-file',f'{dataDir}/DFree_Lcol{Lcol}.out','-time','-node',2,'-dof',1,2,3,'disp')     #  displacements of free nodes
    ops.recorder('Node','-file',f'{dataDir}/DBase_Lcol{Lcol}.out','-time','-node',1,'-dof',1,2,3,'disp')     #  displacements of support nodes
    ops.recorder('Node','-file',f'{dataDir}/RBase_Lcol{Lcol}.out','-time','-node',1,'-dof',1,2,3,'reaction')     #  support reaction
    ops.recorder('Element','-file',f'{dataDir}/FCol_Lcol{Lcol}.out','-time','-ele',1,'globalForce')     #  element forces -- column
    ops.recorder('Element','-file',f'{dataDir}/DCol_Lcol{Lcol}.out','-time','-ele',1,'deformation')     #  element deformations -- column


    # define GRAVITY -------------------------------------------------------------
    ops.timeSeries('Linear',1)     # timeSeries Linear 1;
    ops.pattern('Plain',1,1) # 
    ops.load(2,0.,-2000.,0.)     #  node , FX FY MZ -- superstructure-weight
    ops.wipeAnalysis()     # adding this to clear Analysis module 
    ops.constraints('Plain')     #  how it handles boundary conditions
    ops.numberer('Plain')     #  renumber dofs to minimize band-width (optimization), if you want to
    ops.system('BandGeneral')     #  how to store and solve the system of equations in the analysis
    ops.test('NormDispIncr',1.0e-8,6)     #  determine if convergence has been achieved at the end of an iteration step
    ops.algorithm('Newton')     #  use Newtons solution algorithm: updates tangent stiffness at every iteration
    ops.integrator('LoadControl',0.1)     #  determine the next time step for an analysis,   apply gravity in 10 steps
    ops.analysis('Static')     #  define type of analysis static or transient
    ops.analyze(10)     #  perform gravity analysis
    ops.loadConst('-time',0.0)     #  hold gravity constant and restart time

    # define LATERAL load -------------------------------------------------------------
    ops.timeSeries('Linear',2)     # timeSeries Linear 2;
    ops.pattern('Plain',2,2) # 
    ops.load(2,2000.,0.0,0.0)     #  node , FX FY MZ -- representative lateral load at top node

    # pushover: diplacement controlled static analysis
    ops.integrator('DisplacementControl',2,1,0.1)     #  switch to displacement control, for node 11, dof 1, 0.1 increment
    ops.analyze(1000)     #  apply 100 steps of pushover analysis to a displacement of 10

    print(f'Analysis-{count} execution done')

    count +=1

print(f"ALL DONE!!!")
Ex1a.Canti2D.Push.mpi.py
# mpiexec -np 3 python Ex1a.Canti2D.Push.mpi.py

import openseespy.opensees as ops
import numpy as numpy
import matplotlib.pyplot as plt
import sys
import os

ops.start()
pid = ops.getPID()
np = ops.getNP()
print(f'pid {pid} of {np} started')

if len(sys.argv)>1:
    print(f'mpi -- python pid {pid} of {np} Command-Line Arguments (argv): {sys.argv}')

LColList = [100,120,200,240,300,360,400,480]
#-----------------------------------------
dataDir=f'DataPYmpi';                # set up name of data directory
os.makedirs(dataDir, exist_ok=True); # create data directory

count = 0;
for Lcol in LColList:
    # check if count is a multiple of pid : "is it its turn?"
    if count % np == pid:
        ops.wipe()
        
        # SET UP ----------------------------------------------------------------------------
        ops.wipe()     #  clear opensees model
        ops.model('basic','-ndm',2,'-ndf',3)     #  2 dimensions, 3 dof per node
        
        
        # define GEOMETRY -------------------------------------------------------------
        # nodal coordinates:
        ops.node(1,0,0)     #  node , X Y
        ops.node(2,0,Lcol)
        
        # Single point constraints -- Boundary Conditions
        ops.fix(1,1,1,1)     #  node DX DY RZ
        
        # nodal masses:
        ops.mass(2,5.18,0.,0.)     #  node , Mx My Mz, Mass=Weight/g.
        
        # Define ELEMENTS -------------------------------------------------------------
        # define geometric transformation: performs a linear geometric transformation of beam stiffness
        # and resisting force from the basic system to the global-coordinate system
        ops.geomTransf('Linear',1)     #  associate a tag to transformation
        
        # element elasticBeamColumn eleTag iNode jNode A E Iz transfTag
        ops.element('elasticBeamColumn',1,1,2,3600000000,4227,1080000,1)
        
        # Define RECORDERS -------------------------------------------------------------
        ops.recorder('Node','-file',f'{dataDir}/DFree_Lcol{Lcol}.out','-time','-node',2,'-dof',1,2,3,'disp')     #  displacements of free nodes
        ops.recorder('Node','-file',f'{dataDir}/DBase_Lcol{Lcol}.out','-time','-node',1,'-dof',1,2,3,'disp')     #  displacements of support nodes
        ops.recorder('Node','-file',f'{dataDir}/RBase_Lcol{Lcol}.out','-time','-node',1,'-dof',1,2,3,'reaction')     #  support reaction
        ops.recorder('Element','-file',f'{dataDir}/FCol_Lcol{Lcol}.out','-time','-ele',1,'globalForce')     #  element forces -- column
        ops.recorder('Element','-file',f'{dataDir}/DCol_Lcol{Lcol}.out','-time','-ele',1,'deformation')     #  element deformations -- column
        
        # define GRAVITY -------------------------------------------------------------
        ops.timeSeries('Linear',1)     # timeSeries Linear 1;
        ops.pattern('Plain',1,1) # 
        ops.load(2,0.,-2000.,0.)     #  node , FX FY MZ -- superstructure-weight
        ops.wipeAnalysis()     # adding this to clear Analysis module 
        ops.constraints('Plain')     #  how it handles boundary conditions
        ops.numberer('Plain')     #  renumber dofs to minimize band-width (optimization), if you want to
        ops.system('BandGeneral')     #  how to store and solve the system of equations in the analysis
        ops.test('NormDispIncr',1.0e-8,6)     #  determine if convergence has been achieved at the end of an iteration step
        ops.algorithm('Newton')     #  use Newtons solution algorithm: updates tangent stiffness at every iteration
        ops.integrator('LoadControl',0.1)     #  determine the next time step for an analysis,   apply gravity in 10 steps
        ops.analysis('Static')     #  define type of analysis static or transient
        ops.analyze(10)     #  perform gravity analysis
        ops.loadConst('-time',0.0)     #  hold gravity constant and restart time
        
        # define LATERAL load -------------------------------------------------------------
        ops.timeSeries('Linear',2)     # timeSeries Linear 2;
        ops.pattern('Plain',2,2) # 
        ops.load(2,2000.,0.0,0.0)     #  node , FX FY MZ -- representative lateral load at top node
        
        # pushover: diplacement controlled static analysis
        ops.integrator('DisplacementControl',2,1,0.1)     #  switch to displacement control, for node 11, dof 1, 0.1 increment
        ops.analyze(1000)     #  apply 100 steps of pushover analysis to a displacement of 10
        
        print(f'pid {pid} of np={np} Analysis-{count} execution done')

    count +=1

print(f"pid {pid} of np={np} ALL DONE!!!")
Ex1a.Canti2D.Push.mpi4py.py
# mpiexec -np 3 python Ex1a.Canti2D.Push.mpi4py.py

import openseespy.opensees as ops
import numpy as numpy
import matplotlib.pyplot as plt
import sys
import os

from mpi4py import MPI
comm = MPI.COMM_WORLD
pid = comm.Get_rank()
np=comm.Get_size()
print(f'mpi4py -- python pid {pid} of {np} started')

if len(sys.argv)>1:
    print(f'mpi4py -- python pid {pid} of {np} Command-Line Arguments (argv): {sys.argv}')


LColList = [100,120,200,240,300,360,400,480]
#-----------------------------------------
dataDir=f'DataPYmpi4py';                # set up name of data directory
os.makedirs(dataDir, exist_ok=True);    # create data directory

count = 0;
for Lcol in LColList:
    # check if count is a multiple of pid : "is it its turn?"
    if count % np == pid:
        ops.wipe()
        
        # SET UP ----------------------------------------------------------------------------
        ops.wipe()     #  clear opensees model
        ops.model('basic','-ndm',2,'-ndf',3)     #  2 dimensions, 3 dof per node
        
        
        # define GEOMETRY -------------------------------------------------------------
        # nodal coordinates:
        ops.node(1,0,0)     #  node , X Y
        ops.node(2,0,Lcol)
        
        # Single point constraints -- Boundary Conditions
        ops.fix(1,1,1,1)     #  node DX DY RZ
        
        # nodal masses:
        ops.mass(2,5.18,0.,0.)     #  node , Mx My Mz, Mass=Weight/g.
        
        # Define ELEMENTS -------------------------------------------------------------
        # define geometric transformation: performs a linear geometric transformation of beam stiffness
        # and resisting force from the basic system to the global-coordinate system
        ops.geomTransf('Linear',1)     #  associate a tag to transformation
        
        # connectivity: (make A very large, 10e6 times its actual value)
        ops.element('elasticBeamColumn',1,1,2,3600000000,4227,1080000,1)
        
        # Define RECORDERS -------------------------------------------------------------
        ops.recorder('Node','-file',f'{dataDir}/DFree_Lcol{Lcol}.out','-time','-node',2,'-dof',1,2,3,'disp')     #  displacements of free nodes
        ops.recorder('Node','-file',f'{dataDir}/DBase_Lcol{Lcol}.out','-time','-node',1,'-dof',1,2,3,'disp')     #  displacements of support nodes
        ops.recorder('Node','-file',f'{dataDir}/RBase_Lcol{Lcol}.out','-time','-node',1,'-dof',1,2,3,'reaction')     #  support reaction
        ops.recorder('Element','-file',f'{dataDir}/FCol_Lcol{Lcol}.out','-time','-ele',1,'globalForce')     #  element forces -- column
        ops.recorder('Element','-file',f'{dataDir}/DCol_Lcol{Lcol}.out','-time','-ele',1,'deformation')     #  element deformations -- column

        
        # define GRAVITY -------------------------------------------------------------
        ops.timeSeries('Linear',1)     # timeSeries Linear 1;
        ops.pattern('Plain',1,1) # 
        ops.load(2,0.,-2000.,0.)     #  node , FX FY MZ -- superstructure-weight
        ops.wipeAnalysis()     # adding this to clear Analysis module 
        ops.constraints('Plain')     #  how it handles boundary conditions
        ops.numberer('Plain')     #  renumber dofs to minimize band-width (optimization), if you want to
        ops.system('BandGeneral')     #  how to store and solve the system of equations in the analysis
        ops.test('NormDispIncr',1.0e-8,6)     #  determine if convergence has been achieved at the end of an iteration step
        ops.algorithm('Newton')     #  use Newtons solution algorithm: updates tangent stiffness at every iteration
        ops.integrator('LoadControl',0.1)     #  determine the next time step for an analysis,   apply gravity in 10 steps
        ops.analysis('Static')     #  define type of analysis static or transient
        ops.analyze(10)     #  perform gravity analysis
        ops.loadConst('-time',0.0)     #  hold gravity constant and restart time
        
        # define LATERAL load -------------------------------------------------------------
        ops.timeSeries('Linear',2)     # timeSeries Linear 2;
        ops.pattern('Plain',2,2) # 
        ops.load(2,2000.,0.0,0.0)     #  node , FX FY MZ -- representative lateral load at top node
        
        # pushover: diplacement controlled static analysis
        ops.integrator('DisplacementControl',2,1,0.1)     #  switch to displacement control, for node 11, dof 1, 0.1 increment
        ops.analyze(1000)     #  apply 100 steps of pushover analysis to a displacement of 10
        
        print(f'pid {pid} of np={np} Analysis-{count} execution done')

    count +=1

print(f"pid {pid} of np={np} ALL DONE!!!")