test_imexbdf.py 10 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
from dgpy import *
from math import *
import time, os, sys
import gmshPartition

exportData=0
#### Physical constants ###
gamma=1.4
Rd=287.0
p0=1.0e5
g=9.80616
Cv=Rd/(gamma-1.0)
Cp=Rd*(1.0+1.0/(gamma-1.0))
###########################
#### Begin and end time ###
Ti=0.0
Tf=200.0
###########################

#Period of the oscillations
period = 100


order=4
dimension=2

model = GModel()
name='line'
if Msg.GetCommSize()>1:
  partStr='_part_%i' % Msg.GetCommSize()
else:
  partStr=''
model.load(name + partStr + '_2d_XY.msh')

#############  for convergence plot################
###################################################
errp=[]
erru=[]
errv=[]
errthetap=[]
errall=[]
stepsizes=[]

groups = dgGroupCollection(model, dimension, order)
groups.splitGroupsByPhysicalTag();

claw = dgEulerAtmLaw(dimension)
####Load reference solution from file###
reloaded =dgDofContainer(groups,claw.getNbFields())
50
reloaded.importIdx("refSolERK"+str(period)+"/subwindow4/subwindow4.idx")
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
solution = dgDofContainer(groups, claw.getNbFields())

XYZ = groups.getFunctionCoordinates();
def hydrostaticState(z) :
    thetaHs=300.0
    exnerHs=1.0-g/(Cp*thetaHs)*z
    rhoHs=p0/(Rd*thetaHs)*exnerHs**(Cv/Rd)
    return rhoHs,exnerHs,thetaHs

def initialCondition(FCT,  XYZ) :
  FCT.setAll(0)

def boundaryForcing(FCT, XYZ, time) :
    thetac=5
    xc=500.0
    zc=0.0
    rc=250.0
    currentTime=time(0,0)
    # print currentTime 
    for i in range (0,XYZ.size1()) :
        x=XYZ(i,0)
        z=XYZ(i,1)
        rhoHs,exnerHs,thetaHs = hydrostaticState(z)
        r=sqrt((x-xc)**2+(z-zc)**2)
        if (r>rc):
            thetaPert=thetaHs
        else:
            thetaPert=thetaHs+thetac/2.0*(1.0+cos(pi*r/rc))*sin((currentTime)/period*2.0*pi)**2
        rhoPert=rhoHs#p0/(Rd*thetaPert)*exnerHs**(Cv/Rd)
        FCT.set(i,0,rhoPert-rhoHs)
        FCT.set(i,1,0.0)
        FCT.set(i,2,0.0)
        FCT.set(i,3,rhoPert*thetaPert-rhoHs*thetaHs)

def rhoHydrostatic(FCT,  XYZ) :
    for i in range (0,XYZ.size1()) :
        z=XYZ(i,1)
        rhoHs,exnerHs,thetaHs = hydrostaticState(z)
        FCT.set(i,0,rhoHs)

def rhoThetaHydrostatic(FCT,  XYZ) :
    for i in range (0,XYZ.size1()) :
        z=XYZ(i,1)
        rhoHs,exnerHs,thetaHs = hydrostaticState(z)
        FCT.set(i,0,rhoHs*thetaHs)
        
def getVelocity(FCT, sol, XYZ) :
    for i in range (0,sol.size1()) :
        z=XYZ(i,1)
        rhoHs,exnerHs,thetaHs = hydrostaticState(z)
        rho=rhoHs+sol(i,0)
        FCT.set(i,0,sol(i,1)/rho)
        FCT.set(i,1,sol(i,2)/rho)
        FCT.set(i,2,0)

def getRhop(FCT, sol) :
    for i in range (0,sol.size1()) :
        FCT.set(i,0,sol(i,0))

def getpp(FCT, sol, XYZ) :
    for i in range (0,sol.size1()) :
        z=XYZ(i,1)
        rhoHs,exnerHs,thetaHs = hydrostaticState(z)
        rho=rhoHs+sol(i,0)
        rhoTheta = rhoHs*thetaHs+sol(i,3)
        pHs=p0*(rhoHs*thetaHs*Rd/p0)**gamma
        p=p0*(rhoTheta*Rd/p0)**gamma
        FCT.set(i,0,p-pHs)

def getpHs(FCT, sol, XYZ) :
    for i in range (0,sol.size1()) :
        z=XYZ(i,1)
        rhoHs,exnerHs,thetaHs = hydrostaticState(z)
        pHs=p0*(rhoHs*thetaHs*Rd/p0)**gamma
        FCT.set(i,0,pHs)

def getp(FCT, sol, XYZ) :
    for i in range (0,sol.size1()) :
        z=XYZ(i,1)
        rhoHs,exnerHs,thetaHs = hydrostaticState(z)
        rhoTheta = rhoHs*thetaHs+sol(i,3)
        p=p0*(rhoTheta*Rd/p0)**gamma
        FCT.set(i,0,p)

def getThetap(FCT, sol, XYZ) :
    for i in range (0,sol.size1()) :
        z=XYZ(i,1)
        rhoHs,exnerHs,thetaHs = hydrostaticState(z)
        rho=rhoHs+sol(i,0)
        FCT.set(i,0,(sol(i,3)+(rhoHs*thetaHs))/rho-thetaHs)

uv=functionPython(3, getVelocity, [solution.getFunction(), XYZ])
rhop=functionPython(1, getRhop, [solution.getFunction()])
pp=functionPython(1, getpp, [solution.getFunction(), XYZ])
thetap=functionPython(1, getThetap, [solution.getFunction(), XYZ])

initF=functionPython(4, initialCondition, [XYZ])
solution.interpolate(initF)

rhoHs = dgDofContainer(groups, 1)
rhoHs.interpolate(functionPython(1, rhoHydrostatic, [XYZ]))
rhoThetaHs = dgDofContainer(groups, 1)
rhoThetaHs.interpolate(functionPython(1, rhoThetaHydrostatic, [XYZ]))

claw.setHydrostaticState(rhoHs,rhoThetaHs)
claw.setPhysicalConstants(gamma,Rd,p0,g)

boundaryWall = claw.newBoundaryWall()
timeFunction = function.getTime()
boundaryForcingF=functionPython(4, boundaryForcing, [XYZ, timeFunction])
toReplace = VectorFunctorConst([function.getSolution()])
replaceBy = VectorFunctorConst([boundaryForcingF])
boundaryHeated = claw.newOutsideValueBoundaryGeneric("",toReplace,replaceBy)
claw.addBoundaryCondition('bottom_domain', boundaryWall)
claw.addBoundaryCondition('top_domain', boundaryHeated) #top is actually bottom due to the inverse extrusion
claw.addBoundaryCondition('left', boundaryWall)
claw.addBoundaryCondition('right', boundaryWall)

claw.setFilterMode(FILTER_LINEAR)

solution0 = dgDofContainer(groups, claw.getNbFields())
solution0.copy(solution)

nbStepsList = [400, 600, 900, 1350, 2025, 3037, 4555, 6832, 10248]
ampf = 1.5

for iruns in range(2,9) :

    petscIm = linearSystemPETScDouble()
    petscIm.setParameter("petscOptions","-ksp_type preonly -pc_type lu")
    #petscIm.setParameter("petscOptions",  "-pc_type lu -ksp_type preonly -pc_factor_mat_solver_package mumps")
    dofIm = dgDofManager.newDG(groups, claw, petscIm)
    petscImStarter = linearSystemPETScDouble()
    petscImStarter.setParameter("petscOptions",  "")
    dofImStarter = dgDofManager.newDG(groups,claw,petscImStarter)
    #timeIter = dgIMEXLMM(claw, dofIm, IMEX_BDF2)
    timeIter = dgIMEXLMM(claw, dofIm, dofImStarter, IMEX_BDF2,None,True)
    timeIter.getNewton().setVerb(1)     #Verbosity

    solution.copy(solution0)
    #dt=claw.getMinOfTimeSteps(solution)*4.0/(2*timeOrder+1) #infinite because no initial velocity
    nbSteps = nbStepsList[iruns]
    dt=(Tf-Ti)/nbSteps

    if (Msg.GetCommRank()==0):
        print ("Time step:",dt,"No. of steps:",nbSteps)
   
    # Compute relative error
    def getRefSolSq(FCT,  refSol) :
        for i in range(FCT.size1()):
            FCT.set(i,0, refSol(i,0)**2)
            FCT.set(i,1, refSol(i,1)**2)
            FCT.set(i,2, refSol(i,2)**2)
            FCT.set(i,3, refSol(i,3)**2)
    solRefSq=functionPython(4, getRefSolSq, [reloaded.getFunction()])

    def getError(FCT, refSol, compSol) :
        for i in range(FCT.size1()):
            FCT.set(i,0,(refSol(i,0)-compSol(i,0))**2)
            FCT.set(i,1,(refSol(i,1)-compSol(i,1))**2)
            FCT.set(i,2,(refSol(i,2)-compSol(i,2))**2)
            FCT.set(i,3,(refSol(i,3)-compSol(i,3))**2)

    error = functionPython(4, getError, [reloaded.getFunction(), solution.getFunction()])
    integratorError = dgFunctionIntegrator(groups, error)
    intErr = fullMatrixDouble(4,1)
    integratorRefSq = dgFunctionIntegrator(groups, solRefSq)
    intRefSq = fullMatrixDouble(4,1)

    #Export data
    def getExp(FCT, uv, rhop, thetap, pp) :
        for i in range (0,uv.size1()) :
            FCT.set(i,0,uv(i,0))
            FCT.set(i,1,uv(i,1))
            FCT.set(i,2,rhop(i,0))
            FCT.set(i,3,thetap(i,0))
            FCT.set(i,4,pp(i,0))
    Exp=functionPython(5, getExp, [uv, rhop, thetap, pp])
    nCompExp=[2,1,1,1]
    namesExp=["uv","rhop","thetap","pp"]

    prefixFileName = 'output_imexdimsim'+str(iruns)+'/export'

    # sBV = 2.0/3*(order+1.0)
    # muBV = 0.2
    # filterBV = dgFilterBoydVandeven(groups, "", sBV, muBV, True, False)
    # solutionBV = dgDofContainer(solution);
    # minDiff = fullMatrixDouble(claw.getNbFields(),1)
    # maxDiff = fullMatrixDouble(claw.getNbFields(),1)

    t=Ti
    # starter for IMEX DIMSIM
    #timeIter.starter(solution,dt/4.0,dt,t,4) 
    # starter for IMEX LMM
    t=timeIter.starter(solution,dt,t,3)

    n_export=0
    timeStart=time.clock();

    start = time.clock()
    #IMEX LMM
    for i in range(1,nbSteps):
    #for i in range(0,nbSteps):
        if (i%(250*2**iruns) == 0 and exportData):
            solution.exportFunctionVtk(Exp,prefixFileName, t*1000, i,"solution",nCompExp,namesExp)
            if (Msg.GetCommRank()==0):
                print ('\nWriting output',n_export,'at time',t*1000,'ms and step',i,'over',nbSteps)
                print("CFL: %.1e"%(dt/(claw.getMinOfTimeSteps(solution)*4.0)))
            n_export=n_export+1
        #We filter the dof manager and evaluate the max difference with unfiltered solution
        # solutionBV.copy(solution)
        # filterBV.apply(solutionBV)
        # solutionBV.axpy(solution, -1)
        # solutionBV.minmax(minDiff, maxDiff);
        # if (Msg.GetCommRank()==0):
        #     #Print the max difference on density
        #    print("[%.1e]"%max(maxDiff(0,0),-minDiff(0,0)))
        #    print("CFL dt %.1e"%(claw.getMinOfTimeSteps(solution)*4.0/(2*timeOrder+1)))
    ##
        timeIter.iterate(solution, dt, t)
        t=t+dt
        
        if (Msg.GetCommRank()==0):
            sys.stdout.write('.')
            sys.stdout.flush()
    if (Msg.GetCommRank()==0):
        print ('')
        print 'Tf=',t,
    elapsed = (time.clock() - start)
    if (Msg.GetCommRank()==0):
        print ('Time elapsed: ',elapsed,' s')
    if(exportData) :
        solution.exportFunctionVtk(Exp,prefixFileName, Tf*1000, nbSteps,"solution",nCompExp,namesExp)

    integratorError.compute(intErr)
    integratorRefSq.compute(intRefSq)
    errp.append(sqrt(intErr(0,0)/intRefSq(0,0)))
    erru.append(sqrt(intErr(1,0)/intRefSq(1,0)))
    errv.append(sqrt(intErr(2,0)/intRefSq(2,0)))
    errthetap.append(sqrt(intErr(3,0)/intRefSq(3,0)))
    errall.append( sqrt((intErr(0,0)+intErr(1,0)+intErr(2,0)+intErr(3,0))/(intRefSq(0,0)+intRefSq(1,0)+intRefSq(2,0)+intRefSq(3,0))) )
    stepsizes.append(dt)

    if (Msg.GetCommRank()==0):
        print "\nError on rhop:",errp[-1],"rhou:",erru[-1],"rhov:", errv[-1],"rhothetap:",errthetap[-1],"all:",errall[-1]
        if(iruns>2):
            print "rhop order:",log(errp[-2]/errp[-1])/log(ampf)
            print "rhou order:",log(erru[-2]/erru[-1])/log(ampf)
            print "rhov order:",log(errv[-2]/errv[-1])/log(ampf)
            print "rhothetap order:",log(errthetap[-2]/errthetap[-1])/log(ampf)
            print "overall order:",log(errall[-2]/errall[-1])/log(ampf)
    
Msg.Exit(0)