Page MenuHomec4science

laplace_disk_gaussian.py
No OneTemporary

File Metadata

Created
Wed, May 15, 20:41

laplace_disk_gaussian.py

# Copyright (C) 2018-2020 by the RROMPy authors
#
# This file is part of RROMPy.
#
# RROMPy is free software: you can redistribute it and/or modify
# it under the terms of the GNU Lesser General Public License as published by
# the Free Software Foundation, either version 3 of the License, or
# (at your option) any later version.
#
# RROMPy is distributed in the hope that it will be useful,
# but WITHOUT ANY WARRANTY; without even the implied warranty of
# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
# GNU Lesser General Public License for more details.
#
# You should have received a copy of the GNU Lesser General Public License
# along with RROMPy. If not, see <http://www.gnu.org/licenses/>.
#
import numpy as np
from numpy.polynomial.polynomial import polyfit as fit
import fenics as fen
from rrompy.utilities.base.types import paramVal
from rrompy.hfengines.fenics_engines import LaplaceBaseProblemEngine
from rrompy.solver.fenics import fenZERO
from rrompy.utilities.base import verbosityManager as vbMng
from rrompy.utilities.poly_fitting.polynomial import (
PolynomialInterpolator as PI)
from rrompy.solver.fenics import fenics2Vector
class LaplaceDiskGaussian(LaplaceBaseProblemEngine):
"""
Solver for disk Laplace problems with parametric forcing term center.
- \Delta u = C exp(-.5 * ||\cdot - (mu, 0)||^2) in \Omega = B(0, 5)
u = 0 on \partial\Omega.
"""
def __init__(self, n:int, mu0 : paramVal = [0.], *args, **kwargs):
super().__init__(mu0, *args, **kwargs)
self.nbs = 19
import mshr
mesh = mshr.generate_mesh(mshr.Circle(fen.Point(0., 0.), 5.), 3 * n)
self.V = fen.FunctionSpace(mesh, "P", 1)
def buildb(self):
"""Build terms of operator of linear system."""
if self.thbs[0] is None:
self.thbs = self.getMonomialWeights(self.nbs)
bDEIMCoeffs = None
for j in range(self.nbs):
if self.bs[j] is None:
vbMng(self, "INIT", "Assembling forcing term b{}.".format(j),
20)
if bDEIMCoeffs is None:
muDEIM = 3. * np.linspace(0., 1., self.nbs // 2 + 1) ** 2.
muDEIM = np.concatenate((-muDEIM[:0:-1], muDEIM))
for jj, muD in enumerate(muDEIM):
x, y = fen.SpatialCoordinate(self.V.mesh())[:]
C = np.exp(-.5 * muD ** 2.)
CR, CI = np.real(C), np.imag(C)
f0 = ((2 * np.pi) ** -.5
* fen.exp(-.5 * (x ** 2. + y ** 2.)))
muR, muI = np.real(muD), np.imag(muD)
f1R = fen.exp(muR * x) * fen.cos(muI * x)
f1I = fen.exp(muR * x) * fen.sin(muI * x)
fRe = f0 * (CR * f1R - CI * f1I) + fenZERO
fIm = f0 * (CR * f1I + CI * f1R) + fenZERO
parsRe = self.iterReduceQuadratureDegree(zip([fRe],
["forcingTerm{}Real".format(jj)]))
parsIm = self.iterReduceQuadratureDegree(zip([fIm],
["forcingTerm{}Imag".format(jj)]))
LR = fen.dot(fRe, self.v) * fen.dx
LI = fen.dot(fIm, self.v) * fen.dx
DBC0 = fen.DirichletBC(self.V, fenZERO,
self.DirichletBoundary)
bjj = (fenics2Vector(LR, parsRe, DBC0, 1)
+ 1.j * fenics2Vector(LI, parsIm, DBC0, 1))
if jj == 0:
bDEIM = np.empty((self.nbs, len(bjj)),
dtype = np.complex)
bDEIM[jj] = bjj
bDEIMCoeffs = (fit(muDEIM / 3., bDEIM, self.nbs - 1).T
* np.power(3., - np.arange(self.nbs))).T
self.bs[j] = bDEIMCoeffs[j]
vbMng(self, "DEL", "Done assembling forcing term.", 20)
class LaplaceDiskGaussian2(LaplaceDiskGaussian):
"""
Solver for disk Laplace problems with parametric forcing term center.
- \Delta u = C exp(-.5 * ||\cdot - (mu1, mu2)||^2) in \Omega = B(0, 5)
u = 0 on \partial\Omega.
"""
def __init__(self, n:int, mu0 : paramVal = [0., 0.], *args, **kwargs):
super().__init__(n = n, mu0 = mu0, *args, **kwargs)
self.nbs = 16
self.npar = 2
def buildb(self):
"""Build terms of operator of linear system."""
if self.thbs[0] is None:
self.thbs = [None] * self.nbs
bDEIMCoeffs = None
for j in range(self.nbs):
j1, j2 = j % int(self.nbs ** .5), j // int(self.nbs ** .5)
if self.bs[j] is None:
vbMng(self, "INIT", "Assembling forcing term b{}.".format(j),
20)
if bDEIMCoeffs is None:
muD1 = np.linspace(-2., 2., int(self.nbs ** .5))
muDEIM = np.empty((self.nbs, 2))
muDEIM[:, 0] = np.repeat(muD1, int(self.nbs ** .5))
muDEIM[:, 1] = np.tile(muD1, int(self.nbs ** .5))
for jj, muD in enumerate(muDEIM):
x, y = fen.SpatialCoordinate(self.V.mesh())[:]
C = np.exp(-.5 * (muD[0] ** 2. + muD[1] ** 2.))
CR, CI = np.real(C), np.imag(C)
f0 = ((2 * np.pi) ** -.5
* fen.exp(-.5 * (x ** 2. + y ** 2.)))
muxR, muxI = np.real(muD[0]), np.imag(muD[0])
muyR, muyI = np.real(muD[1]), np.imag(muD[1])
f1R = (fen.exp(muxR * x + muyR * y)
* fen.cos(muxI * x + muyI * y))
f1I = (fen.exp(muxR * x + muyR * y)
* fen.sin(muxI * x + muyI * y))
fRe = f0 * (CR * f1R - CI * f1I) + fenZERO
fIm = f0 * (CR * f1I + CI * f1R) + fenZERO
parsRe = self.iterReduceQuadratureDegree(zip([fRe],
["forcingTerm{}Real".format(jj)]))
parsIm = self.iterReduceQuadratureDegree(zip([fIm],
["forcingTerm{}Imag".format(jj)]))
LR = fen.dot(fRe, self.v) * fen.dx
LI = fen.dot(fIm, self.v) * fen.dx
DBC0 = fen.DirichletBC(self.V, fenZERO,
self.DirichletBoundary)
bjj = (fenics2Vector(LR, parsRe, DBC0, 1)
+ 1.j * fenics2Vector(LI, parsIm, DBC0, 1))
if jj == 0:
bDEIM = np.empty((self.nbs, len(bjj)),
dtype = np.complex)
bDEIM[jj] = bjj
p = PI()
wellCond, _ = p.setupByInterpolation(muDEIM, bDEIM,
int(self.nbs ** .5) - 1,
"MONOMIAL", False, False)
bDEIMCoeffs = p.coeffs
self.bs[j1 + int(self.nbs ** .5) * j2] = bDEIMCoeffs[j1, j2]
self.thbs[j1 + int(self.nbs ** .5) * j2] = (
self.getMonomialSingleWeight([j1, j2]))
vbMng(self, "DEL", "Done assembling forcing term.", 20)

Event Timeline