Page MenuHomec4science

all_forcing_engine.py
No OneTemporary

File Metadata

Created
Sun, May 5, 05:29

all_forcing_engine.py

import numpy as np
import fenics as fen
from rrompy.hfengines.linear_problem import LaplaceBaseProblemEngine as LBPE
from rrompy.utilities.base import verbosityManager as vbMng
from rrompy.solver.fenics import fenics2Vector
class AllForcingEngine(LBPE):
def __init__(self, mu0:float, n : int = 30,
degree_threshold : int = np.inf,
verbosity : int = 10, timestamp : bool = True):
super().__init__(mu0 = mu0, degree_threshold = degree_threshold,
verbosity = verbosity, timestamp = timestamp)
mesh = fen.RectangleMesh(fen.Point(-5., -5.), fen.Point(5., 5.), n, n)
self.V = fen.FunctionSpace(mesh, "P", 1)
self.nAs, self.nbs = 1, 4
x, y = fen.SpatialCoordinate(mesh)[:]
scaling = (2. * np.pi) ** -1.
r2 = x ** 2. + y ** 2.
self.forcingCoeffs = [
scaling * fen.exp(- (r2 + 1. - 2. * x + 1. - 2. * y) / 2. / 4.) / 2.,
scaling * fen.exp(- (r2 + 1. + 2. * x + 1. + 2. * y) / 2. / 16.) / 4.,
- scaling * fen.exp(- (r2 + 1. + 2. * x + 1. - 2. * y) / 2. / 9.) / 30.,
scaling * fen.exp(- (r2 + 1. - 2. * x + 1. + 2. * y) / 2. / 25.) / 120.]
def b(self, mu = [], der = 0, homogeneized = False):
mu = self.checkParameter(mu)
if not hasattr(der, "__len__"): der = [der] * self.npar
derI = der[0]
nbsTot = self.nbsH if homogeneized else self.nbs
bs = self.bsH if homogeneized else self.bs
if homogeneized and self.mu0 != self.mu0BC:
self.liftDirichletData(self.mu0)
for j in range(derI, nbsTot):
if bs[j] is None:
vbMng(self, "INIT", "Assembling forcing term b{}.".format(j),
20)
parsRe = self.iterReduceQuadratureDegree([(
self.forcingCoeffs[j],
"forcingCoefficient")])
u0Re = self.DirichletDatum[0]
L0Re = fen.dot(self.forcingCoeffs[j], self.v) * fen.dx
DBCR = fen.DirichletBC(self.V, u0Re, self.DirichletBoundary)
b = fenics2Vector(L0Re, parsRe, DBCR, 1)
if homogeneized:
self.bsH[j] = b
else:
self.bs[j] = b
vbMng(self, "DEL", "Done assembling forcing term.", 20)
return self._assembleb(mu, der, derI, homogeneized)

Event Timeline