Source code for openpnm.contrib._transient_multiphysics
importloggingimportnumpyasnpfromopenpnm.algorithmsimportAlgorithmfromopenpnm.algorithms._solutionimportSolutionContainer,TransientSolutionfromopenpnm.integratorsimportScipyRK45fromopenpnm.utilsimportDocoratorlogger=logging.getLogger(__name__)docstr=Docorator()__all__=['TransientMultiPhysics',]@docstr.dedentclassTransientMultiPhysicsSettings:r""" Parameters ---------- %(AlgorithmSettings.parameters)s algorithms: list List of transient algorithm objects to be solved in a coupled manner """algorithms=[]
[docs]@docstr.dedentclassTransientMultiPhysics(Algorithm):r""" A subclass for transient multiphysics simulations. """def__init__(self,algorithms,**kwargs):super().__init__(**kwargs)self.settings._update(TransientMultiPhysicsSettings())self.settings.algorithms=[alg.nameforalginalgorithms]self._algs=algorithms
[docs]defrun(self,x0,tspan,saveat=None,integrator=None):""" Runs all of the transient algorithms simultaneoulsy and returns the solution. Parameters ---------- x0 : ndarray or float Array (or scalar) containing initial condition values. tspan : array_like Tuple (or array) containing the integration time span. saveat : array_like or float, optional If an array is passed, it signifies the time points at which the solution is to be stored, and if a scalar is passed, it refers to the interval at which the solution is to be stored. integrator : Integrator, optional Integrator object which will be used to to the time stepping. Can be instantiated using openpnm.integrators module. Returns ------- TransientSolution The solution object, which is basically a numpy array with the added functionality that it can be called to return the solution at intermediate times (i.e., those not stored in the solution object). In the case of multiphysics, the solution object is a combined array of solutions for each physics. The solution for each physics is available on each algorithm object independently. """logger.info('Running TransientMultiphysics')ifnp.isscalar(saveat):saveat=np.arange(*tspan,saveat)if(saveatisnotNone)and(tspan[1]notinsaveat):saveat=np.hstack((saveat,[tspan[1]]))integrator=ScipyRK45()ifintegratorisNoneelseintegratorfori,alginenumerate(self._algs):# Perform pre-solve validationsalg._validate_settings()alg._validate_topology_health()alg._validate_linear_system()# Write x0 to algorithm the obj (needed by _update_iterative_props)x0_i=self._get_x0(x0,i)alg['pore.ic']=x0_i=np.ones(alg.Np,dtype=float)*x0_ialg._merge_inital_and_boundary_values()# Build RHS (dx/dt = RHS), then integrate the system of ODEsrhs=self._build_rhs()# Integrate RHS using the given solversoln=integrator.solve(rhs,x0,tspan,saveat)# Return dictionary containing solutionself.soln=SolutionContainer()fori,alginenumerate(self._algs):# Slice soln and attach as TransientSolution object to each algt=soln.tx=soln[i*alg.Np:(i+1)*alg.Np,:]alg.soln=TransientSolution(t,x)# Add solution of each alg to solution dictionaryself.soln[alg.settings['quantity']]=alg.soln
def_run_special(self,x0):...def_build_rhs(self):""" Returns a function handle, which calculates dy/dt = rhs(y, t). Notes ----- ``y`` is a composite array that contains ALL the variables that the multiphysics algorithm solves for, e.g., if the constituent algorithms are ``TransientFickianDiffusion``, and ``TransientFourierConduction``, ``y[0:Np-1]`` refers to the concentration, and ``[Np:2*Np-1]`` refers to the temperature values. """defode_func(t,y):# Initialize RHSrhs=[]fori,alginenumerate(self._algs):# Get x from y, assume alg.Np is same for all algsx=self._get_x0(y,i)# again use helper function# Store x onto algorithm,alg.x=x# Build A and balg._update_A_and_b()A=alg.A.tocsc()b=alg.b# Retrieve volumeV=alg.network[alg.settings["pore_volume"]]# Calcualte RHSrhs_alg=np.hstack(-A.dot(x)+b)/Vrhs=np.hstack((rhs,rhs_alg))returnrhsreturnode_funcdef_get_x0(self,x0,i):tmp=[alg.Npforalginself._algs]idx_end=np.cumsum(tmp)idx_start=np.hstack((0,idx_end[:-1]))x0=x0[idx_start[i]:idx_end[i]]returnx0