1 |
jgs |
149 |
# $Id$ |
2 |
|
|
|
3 |
|
|
""" |
4 |
|
|
Test suite for linearPDEs class |
5 |
|
|
|
6 |
|
|
The tests must be linked with a Domain class object in the setUp method: |
7 |
|
|
|
8 |
|
|
from esys.finley import Rectangle |
9 |
|
|
class Test_LinearPDEOnFinley(Test_LinearPDE): |
10 |
|
|
def setUp(self): |
11 |
|
|
self.domain = Rectangle(10,10,2) |
12 |
|
|
suite = unittest.TestSuite() |
13 |
|
|
suite.addTest(unittest.makeSuite(Test_LinearPDEOnFinley)) |
14 |
|
|
unittest.TextTestRunner(verbosity=2).run(suite) |
15 |
|
|
|
16 |
|
|
""" |
17 |
|
|
|
18 |
|
|
__author__="Lutz Gross, l.gross@uq.edu.au" |
19 |
|
|
__licence__="contact: esys@access.uq.edu.au" |
20 |
|
|
__url__="http://www.iservo.edu.au/esys/escript" |
21 |
|
|
__version__="$Revision$" |
22 |
|
|
__date__="$Date$" |
23 |
|
|
|
24 |
|
|
|
25 |
|
|
|
26 |
|
|
from esys.escript.util import Lsup,kronecker,interpolate |
27 |
|
|
from esys.escript import Function,FunctionOnBoundary,FunctionOnContactZero,Solution,ReducedSolution,Vector,ContinuousFunction,Scalar |
28 |
|
|
from esys.escript.linearPDEs import LinearPDE,IllegalCoefficientValue,Poisson |
29 |
|
|
import numarray |
30 |
|
|
import unittest |
31 |
|
|
|
32 |
|
|
class Test_linearPDEs(unittest.TestCase): |
33 |
|
|
TOL=1.e-6 |
34 |
|
|
DEBUG=False |
35 |
|
|
VERBOSE=False |
36 |
|
|
def check(self,arg,ref_arg,tol=None): |
37 |
|
|
""" |
38 |
|
|
checks if arg and ref_arg are nearly identical using the L{Lsup<esys.escript.util.Lsup>} |
39 |
|
|
""" |
40 |
|
|
if tol==None: tol=self.TOL |
41 |
|
|
return Lsup(arg-ref_arg)<=tol*Lsup(ref_arg) |
42 |
|
|
|
43 |
|
|
class Test_Poisson(Test_linearPDEs): |
44 |
|
|
|
45 |
|
|
def test_config(self): |
46 |
|
|
mypde=Poisson(self.domain,debug=self.DEBUG) |
47 |
|
|
self.failUnlessEqual((mypde.getNumEquations(),mypde.getNumSolutions(),mypde.isSymmetric()),(1,1,True),"set up incorrect") |
48 |
|
|
def test_setCoefficient_q(self): |
49 |
|
|
mypde=Poisson(self.domain,debug=self.DEBUG) |
50 |
|
|
x=self.domain.getX() |
51 |
|
|
q_ref=interpolate(x[0].whereZero(),Solution(self.domain)) |
52 |
|
|
A_ref=kronecker(self.domain) |
53 |
|
|
mypde.setValue(q=x[0].whereZero()) |
54 |
|
|
self.failUnless(self.check(mypde.getCoefficientOfGeneralPDE("A"),A_ref),"A is not kronecker") |
55 |
|
|
self.failUnless(mypde.getCoefficientOfGeneralPDE("B").isEmpty(),"B is not empty") |
56 |
|
|
self.failUnless(mypde.getCoefficientOfGeneralPDE("C").isEmpty(),"C is not empty") |
57 |
|
|
self.failUnless(mypde.getCoefficientOfGeneralPDE("D").isEmpty(),"D is not empty") |
58 |
|
|
self.failUnless(mypde.getCoefficientOfGeneralPDE("X").isEmpty(),"X is not empty") |
59 |
|
|
self.failUnless(mypde.getCoefficientOfGeneralPDE("Y").isEmpty(),"Y is not empty") |
60 |
|
|
self.failUnless(mypde.getCoefficientOfGeneralPDE("y").isEmpty(),"y is not empty") |
61 |
|
|
self.failUnless(mypde.getCoefficientOfGeneralPDE("d").isEmpty(),"d is not empty") |
62 |
|
|
self.failUnless(mypde.getCoefficientOfGeneralPDE("d_contact").isEmpty(),"d_contact is not empty") |
63 |
|
|
self.failUnless(mypde.getCoefficientOfGeneralPDE("y_contact").isEmpty(),"y_contact is not empty") |
64 |
|
|
self.failUnless(self.check(mypde.getCoefficientOfGeneralPDE("q"),q_ref),"q is not empty") |
65 |
|
|
self.failUnless(mypde.getCoefficientOfGeneralPDE("r").isEmpty(),"r is not empty") |
66 |
|
|
def test_setCoefficient_f(self): |
67 |
|
|
mypde=Poisson(self.domain,debug=self.DEBUG) |
68 |
|
|
x=self.domain.getX() |
69 |
|
|
Y_ref=interpolate(x[0],Function(self.domain)) |
70 |
|
|
A_ref=kronecker(self.domain) |
71 |
|
|
mypde.setValue(f=x[0]) |
72 |
|
|
self.failUnless(self.check(mypde.getCoefficientOfGeneralPDE("A"),A_ref),"A is not kronecker") |
73 |
|
|
self.failUnless(mypde.getCoefficientOfGeneralPDE("B").isEmpty(),"B is not empty") |
74 |
|
|
self.failUnless(mypde.getCoefficientOfGeneralPDE("C").isEmpty(),"C is not empty") |
75 |
|
|
self.failUnless(mypde.getCoefficientOfGeneralPDE("D").isEmpty(),"D is not empty") |
76 |
|
|
self.failUnless(mypde.getCoefficientOfGeneralPDE("X").isEmpty(),"X is not empty") |
77 |
|
|
self.failUnless(self.check(mypde.getCoefficientOfGeneralPDE("Y"),Y_ref),"Y is not x[0]") |
78 |
|
|
self.failUnless(mypde.getCoefficientOfGeneralPDE("y").isEmpty(),"y is not empty") |
79 |
|
|
self.failUnless(mypde.getCoefficientOfGeneralPDE("d").isEmpty(),"d is not empty") |
80 |
|
|
self.failUnless(mypde.getCoefficientOfGeneralPDE("d_contact").isEmpty(),"d_contact is not empty") |
81 |
|
|
self.failUnless(mypde.getCoefficientOfGeneralPDE("y_contact").isEmpty(),"y_contact is not empty") |
82 |
|
|
self.failUnless(mypde.getCoefficientOfGeneralPDE("q").isEmpty(),"q is not empty") |
83 |
|
|
self.failUnless(mypde.getCoefficientOfGeneralPDE("r").isEmpty(),"r is not empty") |
84 |
|
|
def test_solve(self): |
85 |
|
|
d=self.domain.getDim() |
86 |
|
|
cf=ContinuousFunction(self.domain) |
87 |
|
|
x=cf.getX() |
88 |
|
|
#construct exact solution: |
89 |
|
|
u_ex=Scalar(1.,cf) |
90 |
|
|
for i in range(d): |
91 |
|
|
u_ex*=x[i]*(2.-x[i]) |
92 |
|
|
#construct mask: |
93 |
|
|
msk=Scalar(0.,cf) |
94 |
|
|
for i in range(d): |
95 |
|
|
msk+=x[i].whereZero() |
96 |
|
|
#construct right hand side |
97 |
|
|
f=Scalar(0,cf) |
98 |
|
|
for i in range(d): |
99 |
|
|
f_p=Scalar(1,cf) |
100 |
|
|
for j in range(d): |
101 |
|
|
if i==j: |
102 |
|
|
f_p*=2. |
103 |
|
|
else: |
104 |
|
|
f_p*=x[j]*(2-x[j]) |
105 |
|
|
f+=f_p |
106 |
|
|
mypde=Poisson(self.domain) |
107 |
|
|
mypde.setValue(f=f,q=msk) |
108 |
|
|
u=mypde.getSolution() |
109 |
|
|
self.failUnless(self.check(u,u_ex,10*self.TOL),"incorrect solution") |
110 |
|
|
|
111 |
|
|
class Test_LinearPDE(Test_linearPDEs): |
112 |
|
|
N=4 |
113 |
|
|
def test_setCoefficient_WithIllegalFunctionSpace(self): |
114 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
115 |
|
|
try: |
116 |
|
|
success=True |
117 |
|
|
mypde.setValue(C=Vector(0.,FunctionOnBoundary(self.domain))) |
118 |
|
|
except IllegalCoefficientValue: |
119 |
|
|
success=False |
120 |
|
|
self.failUnless(not success,'inapropraite function space accepted') |
121 |
|
|
|
122 |
|
|
def test_resetCoefficient_WithWrongShape(self): |
123 |
|
|
mypde=LinearPDE(self.domain,numEquations=2,debug=self.DEBUG) |
124 |
|
|
try: |
125 |
|
|
success=True |
126 |
|
|
mypde.setValue(C=0.) |
127 |
|
|
except IllegalCoefficientValue: |
128 |
|
|
success=False |
129 |
|
|
self.failUnless(not success,'illegal shape accepted') |
130 |
|
|
def test_reducedOn(self): |
131 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
132 |
|
|
x=self.domain.getX() |
133 |
|
|
mypde.setReducedOrderOn() |
134 |
|
|
mypde.setValue(A=kronecker(self.domain),D=x[0],Y=x[0]) |
135 |
|
|
u=mypde.getSolution() |
136 |
|
|
self.failUnless(self.check(u,1.),'solution is wrong.') |
137 |
|
|
|
138 |
|
|
def test_attemptToChangeOrderAfterDefinedCoefficient(self): |
139 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
140 |
|
|
mypde.setValue(D=1.) |
141 |
|
|
try: |
142 |
|
|
success=True |
143 |
|
|
mypde.setReducedOrderOn() |
144 |
|
|
except RuntimeError: |
145 |
|
|
success=False |
146 |
|
|
self.failUnless(not success,'alterion of order after coefficient is changed not detected.') |
147 |
|
|
|
148 |
|
|
def test_reducedOnConfig(self): |
149 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
150 |
|
|
mypde.setReducedOrderOn() |
151 |
|
|
self.failUnlessEqual((mypde.getFunctionSpaceForSolution(),mypde.getFunctionSpaceForEquation()),(ReducedSolution(self.domain),ReducedSolution(self.domain)),"reduced function spaces expected.") |
152 |
|
|
# |
153 |
|
|
# set coefficients for scalars: |
154 |
|
|
# |
155 |
|
|
def test_setCoefficient_A_Scalar(self): |
156 |
|
|
d=self.domain.getDim() |
157 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
158 |
|
|
mypde.setValue(A=numarray.ones((d,d))) |
159 |
|
|
coeff=mypde.getCoefficientOfGeneralPDE("A") |
160 |
|
|
self.failUnlessEqual((coeff.getShape(),coeff.getFunctionSpace(),mypde.getNumSolutions(),mypde.getNumEquations()),((d,d),Function(self.domain),1,1)) |
161 |
|
|
def test_setCoefficient_B_Scalar(self): |
162 |
|
|
d=self.domain.getDim() |
163 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
164 |
|
|
mypde.setValue(B=numarray.ones((d,))) |
165 |
|
|
coeff=mypde.getCoefficientOfGeneralPDE("B") |
166 |
|
|
self.failUnlessEqual((coeff.getShape(),coeff.getFunctionSpace(),mypde.getNumSolutions(),mypde.getNumEquations()),((d,),Function(self.domain),1,1)) |
167 |
|
|
def test_setCoefficient_C_Scalar(self): |
168 |
|
|
d=self.domain.getDim() |
169 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
170 |
|
|
mypde.setValue(C=numarray.ones((d,))) |
171 |
|
|
coeff=mypde.getCoefficientOfGeneralPDE("C") |
172 |
|
|
self.failUnlessEqual((coeff.getShape(),coeff.getFunctionSpace(),mypde.getNumSolutions(),mypde.getNumEquations()),((d,),Function(self.domain),1,1)) |
173 |
|
|
def test_setCoefficient_D_Scalar(self): |
174 |
|
|
d=self.domain.getDim() |
175 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
176 |
|
|
mypde.setValue(D=1.) |
177 |
|
|
coeff=mypde.getCoefficientOfGeneralPDE("D") |
178 |
|
|
self.failUnlessEqual((coeff.getShape(),coeff.getFunctionSpace(),mypde.getNumSolutions(),mypde.getNumEquations()),((),Function(self.domain),1,1)) |
179 |
|
|
def test_setCoefficient_X_Scalar(self): |
180 |
|
|
d=self.domain.getDim() |
181 |
|
|
mypde=LinearPDE(self.domain,numSolutions=3,debug=self.DEBUG) |
182 |
|
|
mypde.setValue(X=numarray.ones((d,))) |
183 |
|
|
coeff=mypde.getCoefficientOfGeneralPDE("X") |
184 |
|
|
self.failUnlessEqual((coeff.getShape(),coeff.getFunctionSpace(),mypde.getNumEquations()),((d,),Function(self.domain),1)) |
185 |
|
|
def test_setCoefficient_Y_Scalar(self): |
186 |
|
|
d=self.domain.getDim() |
187 |
|
|
mypde=LinearPDE(self.domain,numSolutions=3,debug=self.DEBUG) |
188 |
|
|
mypde.setValue(Y=1.) |
189 |
|
|
coeff=mypde.getCoefficientOfGeneralPDE("Y") |
190 |
|
|
self.failUnlessEqual((coeff.getShape(),coeff.getFunctionSpace(),mypde.getNumEquations()),((),Function(self.domain),1)) |
191 |
|
|
def test_setCoefficient_y_Scalar(self): |
192 |
|
|
d=self.domain.getDim() |
193 |
|
|
mypde=LinearPDE(self.domain,numSolutions=3,debug=self.DEBUG) |
194 |
|
|
mypde.setValue(y=1.) |
195 |
|
|
coeff=mypde.getCoefficientOfGeneralPDE("y") |
196 |
|
|
self.failUnlessEqual((coeff.getShape(),coeff.getFunctionSpace(),mypde.getNumEquations()),((),FunctionOnBoundary(self.domain),1)) |
197 |
|
|
def test_setCoefficient_d_Scalar(self): |
198 |
|
|
d=self.domain.getDim() |
199 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
200 |
|
|
mypde.setValue(d=1.) |
201 |
|
|
coeff=mypde.getCoefficientOfGeneralPDE("d") |
202 |
|
|
self.failUnlessEqual((coeff.getShape(),coeff.getFunctionSpace(),mypde.getNumSolutions(),mypde.getNumEquations()),((),FunctionOnBoundary(self.domain),1,1)) |
203 |
|
|
def test_setCoefficient_d_contact_Scalar(self): |
204 |
|
|
d=self.domain.getDim() |
205 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
206 |
|
|
mypde.setValue(d_contact=1.) |
207 |
|
|
coeff=mypde.getCoefficientOfGeneralPDE("d_contact") |
208 |
|
|
self.failUnlessEqual((coeff.getShape(),coeff.getFunctionSpace(),mypde.getNumSolutions(),mypde.getNumEquations()),((),FunctionOnContactZero(self.domain),1,1)) |
209 |
|
|
def test_setCoefficient_y_contact_Scalar(self): |
210 |
|
|
d=self.domain.getDim() |
211 |
|
|
mypde=LinearPDE(self.domain,numSolutions=3,debug=self.DEBUG) |
212 |
|
|
mypde.setValue(y_contact=1.) |
213 |
|
|
coeff=mypde.getCoefficientOfGeneralPDE("y_contact") |
214 |
|
|
self.failUnlessEqual((coeff.getShape(),coeff.getFunctionSpace(),mypde.getNumEquations()),((),FunctionOnContactZero(self.domain),1)) |
215 |
|
|
def test_setCoefficient_r_Scalar(self): |
216 |
|
|
d=self.domain.getDim() |
217 |
|
|
mypde=LinearPDE(self.domain,numEquations=3,debug=self.DEBUG) |
218 |
|
|
mypde.setValue(r=1.) |
219 |
|
|
coeff=mypde.getCoefficientOfGeneralPDE("r") |
220 |
|
|
self.failUnlessEqual((coeff.getShape(),coeff.getFunctionSpace(),mypde.getNumSolutions()),((),Solution(self.domain),1)) |
221 |
|
|
def test_setCoefficient_q_Scalar(self): |
222 |
|
|
d=self.domain.getDim() |
223 |
|
|
mypde=LinearPDE(self.domain,numEquations=3,debug=self.DEBUG) |
224 |
|
|
mypde.setValue(q=1.) |
225 |
|
|
coeff=mypde.getCoefficientOfGeneralPDE("q") |
226 |
|
|
self.failUnlessEqual((coeff.getShape(),coeff.getFunctionSpace(),mypde.getNumSolutions()),((),Solution(self.domain),1)) |
227 |
|
|
def test_setCoefficient_r_Scalar_reducedOn(self): |
228 |
|
|
d=self.domain.getDim() |
229 |
|
|
mypde=LinearPDE(self.domain,numEquations=3,debug=self.DEBUG) |
230 |
|
|
mypde.setReducedOrderOn() |
231 |
|
|
mypde.setValue(r=1.) |
232 |
|
|
coeff=mypde.getCoefficientOfGeneralPDE("r") |
233 |
|
|
self.failUnlessEqual((coeff.getShape(),coeff.getFunctionSpace(),mypde.getNumSolutions()),((),ReducedSolution(self.domain),1)) |
234 |
|
|
def test_setCoefficient_q_Scalar_reducedOn(self): |
235 |
|
|
d=self.domain.getDim() |
236 |
|
|
mypde=LinearPDE(self.domain,numEquations=3,debug=self.DEBUG) |
237 |
|
|
mypde.setReducedOrderOn() |
238 |
|
|
mypde.setValue(q=1.) |
239 |
|
|
coeff=mypde.getCoefficientOfGeneralPDE("q") |
240 |
|
|
self.failUnlessEqual((coeff.getShape(),coeff.getFunctionSpace(),mypde.getNumSolutions()),((),ReducedSolution(self.domain),1)) |
241 |
|
|
|
242 |
|
|
# |
243 |
|
|
# set coefficients for systems: |
244 |
|
|
# |
245 |
|
|
def test_setCoefficient_A_System(self): |
246 |
|
|
d=self.domain.getDim() |
247 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
248 |
|
|
mypde.setValue(A=numarray.ones((self.N,d,self.N,d))) |
249 |
|
|
coeff=mypde.getCoefficientOfGeneralPDE("A") |
250 |
|
|
self.failUnlessEqual((coeff.getShape(),coeff.getFunctionSpace(),mypde.getNumSolutions(),mypde.getNumEquations()),((self.N,d,self.N,d),Function(self.domain),self.N,self.N)) |
251 |
|
|
def test_setCoefficient_B_System(self): |
252 |
|
|
d=self.domain.getDim() |
253 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
254 |
|
|
mypde.setValue(B=numarray.ones((self.N,d,self.N))) |
255 |
|
|
coeff=mypde.getCoefficientOfGeneralPDE("B") |
256 |
|
|
self.failUnlessEqual((coeff.getShape(),coeff.getFunctionSpace(),mypde.getNumSolutions(),mypde.getNumEquations()),((self.N,d,self.N),Function(self.domain),self.N,self.N)) |
257 |
|
|
def test_setCoefficient_C_System(self): |
258 |
|
|
d=self.domain.getDim() |
259 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
260 |
|
|
mypde.setValue(C=numarray.ones((self.N,self.N,d))) |
261 |
|
|
coeff=mypde.getCoefficientOfGeneralPDE("C") |
262 |
|
|
self.failUnlessEqual((coeff.getShape(),coeff.getFunctionSpace(),mypde.getNumSolutions(),mypde.getNumEquations()),((self.N,self.N,d),Function(self.domain),self.N,self.N)) |
263 |
|
|
def test_setCoefficient_D_System(self): |
264 |
|
|
d=self.domain.getDim() |
265 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
266 |
|
|
mypde.setValue(D=numarray.ones((self.N,self.N))) |
267 |
|
|
coeff=mypde.getCoefficientOfGeneralPDE("D") |
268 |
|
|
self.failUnlessEqual((coeff.getShape(),coeff.getFunctionSpace(),mypde.getNumSolutions(),mypde.getNumEquations()),((self.N,self.N),Function(self.domain),self.N,self.N)) |
269 |
|
|
def test_setCoefficient_X_System(self): |
270 |
|
|
d=self.domain.getDim() |
271 |
|
|
mypde=LinearPDE(self.domain,numSolutions=3,debug=self.DEBUG) |
272 |
|
|
mypde.setValue(X=numarray.ones((self.N,d))) |
273 |
|
|
coeff=mypde.getCoefficientOfGeneralPDE("X") |
274 |
|
|
self.failUnlessEqual((coeff.getShape(),coeff.getFunctionSpace(),mypde.getNumEquations()),((self.N,d),Function(self.domain),self.N)) |
275 |
|
|
def test_setCoefficient_Y_System(self): |
276 |
|
|
d=self.domain.getDim() |
277 |
|
|
mypde=LinearPDE(self.domain,numSolutions=3,debug=self.DEBUG) |
278 |
|
|
mypde.setValue(Y=numarray.ones((self.N,))) |
279 |
|
|
coeff=mypde.getCoefficientOfGeneralPDE("Y") |
280 |
|
|
self.failUnlessEqual((coeff.getShape(),coeff.getFunctionSpace(),mypde.getNumEquations()),((self.N,),Function(self.domain),self.N)) |
281 |
|
|
def test_setCoefficient_y_System(self): |
282 |
|
|
d=self.domain.getDim() |
283 |
|
|
mypde=LinearPDE(self.domain,numSolutions=3,debug=self.DEBUG) |
284 |
|
|
mypde.setValue(y=numarray.ones((self.N,))) |
285 |
|
|
coeff=mypde.getCoefficientOfGeneralPDE("y") |
286 |
|
|
self.failUnlessEqual((coeff.getShape(),coeff.getFunctionSpace(),mypde.getNumEquations()),((self.N,),FunctionOnBoundary(self.domain),self.N)) |
287 |
|
|
def test_setCoefficient_d_System(self): |
288 |
|
|
d=self.domain.getDim() |
289 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
290 |
|
|
mypde.setValue(d=numarray.ones((self.N,self.N))) |
291 |
|
|
coeff=mypde.getCoefficientOfGeneralPDE("d") |
292 |
|
|
self.failUnlessEqual((coeff.getShape(),coeff.getFunctionSpace(),mypde.getNumSolutions(),mypde.getNumEquations()),((self.N,self.N),FunctionOnBoundary(self.domain),self.N,self.N)) |
293 |
|
|
def test_setCoefficient_d_contact_System(self): |
294 |
|
|
d=self.domain.getDim() |
295 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
296 |
|
|
mypde.setValue(d_contact=numarray.ones((self.N,self.N))) |
297 |
|
|
coeff=mypde.getCoefficientOfGeneralPDE("d_contact") |
298 |
|
|
self.failUnlessEqual((coeff.getShape(),coeff.getFunctionSpace(),mypde.getNumSolutions(),mypde.getNumEquations()),((self.N,self.N),FunctionOnContactZero(self.domain),self.N,self.N)) |
299 |
|
|
def test_setCoefficient_y_contact_System(self): |
300 |
|
|
d=self.domain.getDim() |
301 |
|
|
mypde=LinearPDE(self.domain,numSolutions=3,debug=self.DEBUG) |
302 |
|
|
mypde.setValue(y_contact=numarray.ones((self.N,))) |
303 |
|
|
coeff=mypde.getCoefficientOfGeneralPDE("y_contact") |
304 |
|
|
self.failUnlessEqual((coeff.getShape(),coeff.getFunctionSpace(),mypde.getNumEquations()),((self.N,),FunctionOnContactZero(self.domain),self.N)) |
305 |
|
|
def test_setCoefficient_r_System(self): |
306 |
|
|
d=self.domain.getDim() |
307 |
|
|
mypde=LinearPDE(self.domain,numEquations=3,debug=self.DEBUG) |
308 |
|
|
mypde.setValue(r=numarray.ones((self.N,))) |
309 |
|
|
coeff=mypde.getCoefficientOfGeneralPDE("r") |
310 |
|
|
self.failUnlessEqual((coeff.getShape(),coeff.getFunctionSpace(),mypde.getNumSolutions()),((self.N,),Solution(self.domain),self.N)) |
311 |
|
|
def test_setCoefficient_q_System(self): |
312 |
|
|
d=self.domain.getDim() |
313 |
|
|
mypde=LinearPDE(self.domain,numEquations=3,debug=self.DEBUG) |
314 |
|
|
mypde.setValue(q=numarray.ones((self.N,))) |
315 |
|
|
coeff=mypde.getCoefficientOfGeneralPDE("q") |
316 |
|
|
self.failUnlessEqual((coeff.getShape(),coeff.getFunctionSpace(),mypde.getNumSolutions()),((self.N,),Solution(self.domain),self.N)) |
317 |
|
|
def test_setCoefficient_r_System_reducedOn(self): |
318 |
|
|
d=self.domain.getDim() |
319 |
|
|
mypde=LinearPDE(self.domain,numEquations=3,debug=self.DEBUG) |
320 |
|
|
mypde.setReducedOrderOn() |
321 |
|
|
mypde.setValue(r=numarray.ones((self.N,))) |
322 |
|
|
coeff=mypde.getCoefficientOfGeneralPDE("r") |
323 |
|
|
self.failUnlessEqual((coeff.getShape(),coeff.getFunctionSpace(),mypde.getNumSolutions()),((self.N,),ReducedSolution(self.domain),self.N)) |
324 |
|
|
def test_setCoefficient_q_System_reducedOn(self): |
325 |
|
|
d=self.domain.getDim() |
326 |
|
|
mypde=LinearPDE(self.domain,numEquations=3,debug=self.DEBUG) |
327 |
|
|
mypde.setReducedOrderOn() |
328 |
|
|
mypde.setValue(q=numarray.ones((self.N,))) |
329 |
|
|
coeff=mypde.getCoefficientOfGeneralPDE("q") |
330 |
|
|
self.failUnlessEqual((coeff.getShape(),coeff.getFunctionSpace(),mypde.getNumSolutions()),((self.N,),ReducedSolution(self.domain),self.N)) |
331 |
|
|
|
332 |
|
|
def test_resetCoefficient_HomogeneousConstraint(self): |
333 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
334 |
|
|
x=self.domain.getX() |
335 |
|
|
mypde.setValue(A=kronecker(self.domain),Y=1.,q=x[0].whereZero()) |
336 |
|
|
u1=mypde.getSolution() |
337 |
|
|
mypde.setValue(Y=2.) |
338 |
|
|
u2=mypde.getSolution() |
339 |
|
|
self.failUnless(self.check(u2,2*u1),'solution is wrong.') |
340 |
|
|
|
341 |
|
|
def test_resetCoefficient_InHomogeneousConstraint(self): |
342 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
343 |
|
|
x=self.domain.getX() |
344 |
|
|
mypde.setValue(A=kronecker(self.domain),D=1.,Y=1.,r=1,q=x[0].whereZero()) |
345 |
|
|
u1=mypde.getSolution() |
346 |
|
|
mypde.setValue(Y=2.,D=2) |
347 |
|
|
u2=mypde.getSolution() |
348 |
|
|
self.failUnless(self.check(u2,u1),'solution is wrong.') |
349 |
|
|
mypde.setValue(r=2,Y=4.) |
350 |
|
|
u2=mypde.getSolution() |
351 |
|
|
self.failUnless(self.check(u2,2*u1),'solution is wrong.') |
352 |
|
|
|
353 |
|
|
def test_symmetryCheckTrue_System(self): |
354 |
|
|
d=self.domain.getDim() |
355 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
356 |
|
|
A=numarray.ones((self.N,d,self.N,d)) |
357 |
|
|
C=2*numarray.ones((self.N,self.N,d)) |
358 |
|
|
B=2*numarray.ones((self.N,d,self.N)) |
359 |
|
|
D=3*numarray.ones((self.N,self.N)) |
360 |
|
|
d=4*numarray.ones((self.N,self.N)) |
361 |
|
|
d_contact=5*numarray.ones((self.N,self.N)) |
362 |
|
|
mypde.setValue(A=A,B=B,C=C,D=D,d=d,d_contact=d_contact) |
363 |
|
|
self.failUnless(mypde.checkSymmetry(verbose=False),"symmetry detected") |
364 |
|
|
|
365 |
|
|
def test_symmetryCheckFalse_A_System(self): |
366 |
|
|
d=self.domain.getDim() |
367 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
368 |
|
|
A=numarray.ones((self.N,d,self.N,d)) |
369 |
|
|
A[1,1,1,0]=0. |
370 |
|
|
mypde.setValue(A=A) |
371 |
|
|
self.failUnless(not mypde.checkSymmetry(verbose=False),"symmetry detected") |
372 |
|
|
def test_symmetryCheckFalse_BC_System(self): |
373 |
|
|
d=self.domain.getDim() |
374 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
375 |
|
|
C=2*numarray.ones((self.N,self.N,d)) |
376 |
|
|
B=2*numarray.ones((self.N,d,self.N)) |
377 |
|
|
B[0,0,1]=1. |
378 |
|
|
mypde.setValue(B=B,C=C) |
379 |
|
|
self.failUnless(not mypde.checkSymmetry(verbose=False),"symmetry detected") |
380 |
|
|
|
381 |
|
|
def test_symmetryCheckFalse_D_System(self): |
382 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
383 |
|
|
D=3*numarray.ones((self.N,self.N)) |
384 |
|
|
D[0,1]=0. |
385 |
|
|
mypde.setValue(D=D) |
386 |
|
|
self.failUnless(not mypde.checkSymmetry(verbose=False),"symmetry detected") |
387 |
|
|
|
388 |
|
|
def test_symmetryCheckFalse_d_System(self): |
389 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
390 |
|
|
d=4*numarray.ones((self.N,self.N)) |
391 |
|
|
d[0,1]=0. |
392 |
|
|
mypde.setValue(d=d) |
393 |
|
|
self.failUnless(not mypde.checkSymmetry(verbose=False),"symmetry detected") |
394 |
|
|
|
395 |
|
|
def test_symmetryCheckFalse_d_contact_System(self): |
396 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
397 |
|
|
d_contact=5*numarray.ones((self.N,self.N)) |
398 |
|
|
d_contact[0,1]=0. |
399 |
|
|
mypde.setValue(d_contact=d_contact) |
400 |
|
|
self.failUnless(not mypde.checkSymmetry(verbose=False),"symmetry detected") |
401 |
|
|
|
402 |
|
|
def test_symmetryCheckTrue_Scalar(self): |
403 |
|
|
d=self.domain.getDim() |
404 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
405 |
|
|
A=numarray.ones((d,d)) |
406 |
|
|
C=2*numarray.ones((d,)) |
407 |
|
|
B=2*numarray.ones((d,)) |
408 |
|
|
D=3 |
409 |
|
|
d=4 |
410 |
|
|
d_contact=5 |
411 |
|
|
mypde.setValue(A=A,B=B,C=C,D=D,d=d,d_contact=d_contact) |
412 |
|
|
self.failUnless(mypde.checkSymmetry(verbose=False),"symmetry detected") |
413 |
|
|
|
414 |
|
|
def test_symmetryCheckFalse_A_Scalar(self): |
415 |
|
|
d=self.domain.getDim() |
416 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
417 |
|
|
A=numarray.ones((d,d)) |
418 |
|
|
A[1,0]=0. |
419 |
|
|
mypde.setValue(A=A) |
420 |
|
|
self.failUnless(not mypde.checkSymmetry(verbose=False),"symmetry detected") |
421 |
|
|
def test_symmetryCheckFalse_BC_Scalar(self): |
422 |
|
|
d=self.domain.getDim() |
423 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
424 |
|
|
C=2*numarray.ones((d,)) |
425 |
|
|
B=2*numarray.ones((d,)) |
426 |
|
|
B[0]=1. |
427 |
|
|
mypde.setValue(B=B,C=C) |
428 |
|
|
self.failUnless(not mypde.checkSymmetry(verbose=False),"symmetry detected") |
429 |
|
|
# |
430 |
|
|
# solver checks: |
431 |
|
|
# |
432 |
|
|
def test_symmetryOnIterative(self): |
433 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
434 |
|
|
mypde.setValue(A=kronecker(self.domain),D=1.,Y=1.) |
435 |
|
|
u=mypde.getSolution(verbose=self.VERBOSE) |
436 |
|
|
self.failUnless(self.check(u,1.),'solution is wrong.') |
437 |
|
|
def test_symmetryOnDirect(self): |
438 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
439 |
|
|
mypde.setValue(A=kronecker(self.domain),D=1.,Y=1.) |
440 |
|
|
mypde.setSolverMethod(mypde.DIRECT) |
441 |
|
|
u=mypde.getSolution(verbose=self.VERBOSE) |
442 |
|
|
self.failUnless(self.check(u,1.),'solution is wrong.') |
443 |
|
|
def test_PCG_JACOBI(self): |
444 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
445 |
|
|
mypde.setValue(A=kronecker(self.domain),D=1.,Y=1.) |
446 |
|
|
mypde.setSolverMethod(mypde.PCG) |
447 |
|
|
u=mypde.getSolution(verbose=self.VERBOSE,preconditioner=mypde.JACOBI) |
448 |
|
|
self.failUnless(self.check(u,1.),'solution is wrong.') |
449 |
|
|
def test_PCG_ILU0(self): |
450 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
451 |
|
|
mypde.setValue(A=kronecker(self.domain),D=1.,Y=1.) |
452 |
|
|
mypde.setSolverMethod(mypde.PCG) |
453 |
|
|
u=mypde.getSolution(verbose=self.VERBOSE,preconditioner=mypde.ILU0) |
454 |
|
|
self.failUnless(self.check(u,1.),'solution is wrong.') |
455 |
|
|
def test_DIRECT(self): |
456 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
457 |
|
|
mypde.setValue(A=kronecker(self.domain),D=1.,Y=1.) |
458 |
|
|
mypde.setSolverMethod(mypde.DIRECT) |
459 |
|
|
u=mypde.getSolution(verbose=self.VERBOSE,preconditioner=mypde.ILU0) |
460 |
|
|
self.failUnless(self.check(u,1.),'solution is wrong.') |
461 |
|
|
def test_BICGSTAB_JACOBI(self): |
462 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
463 |
|
|
mypde.setSolverMethod(mypde.BICGSTAB) |
464 |
|
|
mypde.setValue(A=kronecker(self.domain),D=1.,Y=1.) |
465 |
|
|
u=mypde.getSolution(verbose=self.VERBOSE,preconditioner=mypde.JACOBI) |
466 |
|
|
self.failUnless(self.check(u,1.),'solution is wrong.') |
467 |
|
|
def test_BICGSTAB_ILU0(self): |
468 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
469 |
|
|
mypde.setValue(A=kronecker(self.domain),D=1.,Y=1.) |
470 |
|
|
mypde.setSolverMethod(mypde.BICGSTAB) |
471 |
|
|
u=mypde.getSolution(verbose=self.VERBOSE,preconditioner=mypde.ILU0) |
472 |
|
|
self.failUnless(self.check(u,1.),'solution is wrong.') |
473 |
|
|
def test_PRES20_JACOBI(self): |
474 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
475 |
|
|
mypde.setValue(A=kronecker(self.domain),D=1.,Y=1.) |
476 |
|
|
mypde.setSolverMethod(mypde.PRES20) |
477 |
|
|
u=mypde.getSolution(verbose=self.VERBOSE,preconditioner=mypde.JACOBI) |
478 |
|
|
self.failUnless(self.check(u,1.),'solution is wrong.') |
479 |
|
|
def test_PRES20_ILU0(self): |
480 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
481 |
|
|
mypde.setValue(A=kronecker(self.domain),D=1.,Y=1.) |
482 |
|
|
mypde.setSolverMethod(mypde.PRES20) |
483 |
|
|
u=mypde.getSolution(verbose=self.VERBOSE,preconditioner=mypde.ILU0) |
484 |
|
|
self.failUnless(self.check(u,1.),'solution is wrong.') |
485 |
|
|
def test_GMRESnoRestart_JACOBI(self): |
486 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
487 |
|
|
mypde.setValue(A=kronecker(self.domain),D=1.,Y=1.) |
488 |
|
|
mypde.setSolverMethod(mypde.GMRES) |
489 |
|
|
u=mypde.getSolution(verbose=self.VERBOSE,preconditioner=mypde.JACOBI) |
490 |
|
|
self.failUnless(self.check(u,1.),'solution is wrong.') |
491 |
|
|
def test_GMRESnoRestart_ILU0(self): |
492 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
493 |
|
|
mypde.setValue(A=kronecker(self.domain),D=1.,Y=1.) |
494 |
|
|
mypde.setSolverMethod(mypde.GMRES) |
495 |
|
|
u=mypde.getSolution(verbose=self.VERBOSE,preconditioner=mypde.ILU0) |
496 |
|
|
self.failUnless(self.check(u,1.),'solution is wrong.') |
497 |
|
|
def test_GMRES_JACOBI(self): |
498 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
499 |
|
|
mypde.setValue(A=kronecker(self.domain),D=1.,Y=1.) |
500 |
|
|
mypde.setSolverMethod(mypde.GMRES) |
501 |
|
|
u=mypde.getSolution(verbose=self.VERBOSE,preconditioner=mypde.JACOBI) |
502 |
|
|
self.failUnless(self.check(u,1.),'solution is wrong.') |
503 |
|
|
def test_GMRES_ILU0(self): |
504 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
505 |
|
|
mypde.setValue(A=kronecker(self.domain),D=1.,Y=1.) |
506 |
|
|
mypde.setSolverMethod(mypde.GMRES) |
507 |
|
|
u=mypde.getSolution(verbose=self.VERBOSE,preconditioner=mypde.ILU0) |
508 |
|
|
self.failUnless(self.check(u,1.),'solution is wrong.') |
509 |
|
|
def test_GMRES_truncation_restart_JACOBI(self): |
510 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
511 |
|
|
mypde.setValue(A=kronecker(self.domain),D=1.,Y=1.) |
512 |
|
|
mypde.setSolverMethod(mypde.GMRES) |
513 |
|
|
u=mypde.getSolution(verbose=self.VERBOSE,preconditioner=mypde.JACOBI,truncation=10,restart=20) |
514 |
|
|
self.failUnless(self.check(u,1.),'solution is wrong.') |
515 |
|
|
def test_GMRES_truncation_restart_ILU0(self): |
516 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
517 |
|
|
mypde.setValue(A=kronecker(self.domain),D=1.,Y=1.) |
518 |
|
|
mypde.setSolverMethod(mypde.GMRES) |
519 |
|
|
u=mypde.getSolution(verbose=self.VERBOSE,preconditioner=mypde.ILU0,truncation=10,restart=20) |
520 |
|
|
self.failUnless(self.check(u,1.),'solution is wrong.') |
521 |
|
|
def test_Lumping_attemptToSetA(self): |
522 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
523 |
|
|
try: |
524 |
|
|
success=True |
525 |
|
|
mypde.setSolverMethod(mypde.LUMPING) |
526 |
|
|
mypde.setValue(A=kronecker(self.domain)) |
527 |
|
|
u=mypde.getSolution(verbose=self.VERBOSE) |
528 |
|
|
except Warning: |
529 |
|
|
success=False |
530 |
|
|
self.failUnless(not success,'warning should be issued') |
531 |
|
|
def test_Lumping_attemptToSetB(self): |
532 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
533 |
|
|
try: |
534 |
|
|
success=True |
535 |
|
|
mypde.setSolverMethod(mypde.LUMPING) |
536 |
|
|
mypde.setValue(B=kronecker(self.domain)[0]) |
537 |
|
|
u=mypde.getSolution(verbose=self.VERBOSE) |
538 |
|
|
except Warning: |
539 |
|
|
success=False |
540 |
|
|
self.failUnless(not success,'warning should be issued') |
541 |
|
|
def test_Lumping_attemptToSetC(self): |
542 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
543 |
|
|
try: |
544 |
|
|
success=True |
545 |
|
|
mypde.setSolverMethod(mypde.LUMPING) |
546 |
|
|
mypde.setValue(C=kronecker(self.domain)[0]) |
547 |
|
|
u=mypde.getSolution(verbose=self.VERBOSE) |
548 |
|
|
except Warning: |
549 |
|
|
success=False |
550 |
|
|
self.failUnless(not success,'warning should be issued') |
551 |
|
|
|
552 |
|
|
def test_Lumping(self): |
553 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
554 |
|
|
mypde.setSolverMethod(mypde.LUMPING) |
555 |
|
|
mypde.setValue(D=1.,Y=1.) |
556 |
|
|
u=mypde.getSolution(verbose=self.VERBOSE,preconditioner=mypde.ILU0) |
557 |
|
|
self.failUnless(self.check(u,1.),'solution is wrong.') |
558 |
|
|
def test_Constrained_Lumping(self): |
559 |
|
|
x=self.domain.getX() |
560 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
561 |
|
|
mypde.setSolverMethod(mypde.LUMPING) |
562 |
|
|
mypde.setValue(D=1.,Y=1.,q=x[0].whereZero(),r=1.) |
563 |
|
|
u=mypde.getSolution(verbose=self.VERBOSE,preconditioner=mypde.ILU0) |
564 |
|
|
self.failUnless(self.check(u,1.),'solution is wrong.') |
565 |
|
|
def test_Lumping_updateRHS(self): |
566 |
|
|
x=self.domain.getX() |
567 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
568 |
|
|
mypde.setSolverMethod(mypde.LUMPING) |
569 |
|
|
mypde.setValue(D=1.,Y=1.) |
570 |
|
|
u=mypde.getSolution(verbose=self.VERBOSE,preconditioner=mypde.ILU0) |
571 |
|
|
self.failUnless(self.check(u,1.),'first solution is wrong.') |
572 |
|
|
mypde.setValue(Y=2.,q=x[0].whereZero(),r=2.) |
573 |
|
|
u=mypde.getSolution(verbose=self.VERBOSE,preconditioner=mypde.ILU0) |
574 |
|
|
self.failUnless(self.check(u,2.),'second solution is wrong.') |
575 |
|
|
def test_Lumping_updateOperator(self): |
576 |
|
|
x=self.domain.getX() |
577 |
|
|
mypde=LinearPDE(self.domain,debug=self.DEBUG) |
578 |
|
|
mypde.setSolverMethod(mypde.LUMPING) |
579 |
|
|
mypde.setValue(D=1.,Y=1.) |
580 |
|
|
u=mypde.getSolution(verbose=self.VERBOSE,preconditioner=mypde.ILU0) |
581 |
|
|
mypde.setValue(D=2.) |
582 |
|
|
u=mypde.getSolution(verbose=self.VERBOSE,preconditioner=mypde.ILU0) |
583 |
|
|
self.failUnless(self.check(u,0.5),'second solution is wrong.') |
584 |
|
|
|