source: issm/trunk/src/m/classes/diagnostic.py@ 15396

Last change on this file since 15396 was 15396, checked in by Mathieu Morlighem, 12 years ago

merged trunk-jpl and trunk for revision 15394

File size: 10.1 KB
Line 
1import numpy
2import sys
3import copy
4from fielddisplay import fielddisplay
5from EnumDefinitions import *
6from checkfield import *
7from WriteData import *
8
9class diagnostic(object):
10 """
11 DIAGNOSTIC class definition
12
13 Usage:
14 diagnostic=diagnostic();
15 """
16
17 def __init__(self): # {{{
18 self.spcvx = float('NaN')
19 self.spcvy = float('NaN')
20 self.spcvz = float('NaN')
21 self.restol = 0
22 self.reltol = 0
23 self.abstol = 0
24 self.isnewton = 0
25 self.stokesreconditioning = 0
26 self.viscosity_overshoot = 0
27 self.icefront = float('NaN')
28 self.maxiter = 0
29 self.shelf_dampening = 0
30 self.vertex_pairing = float('NaN')
31 self.penalty_factor = float('NaN')
32 self.rift_penalty_lock = float('NaN')
33 self.rift_penalty_threshold = 0
34 self.referential = float('NaN')
35 self.loadingforce = float('NaN')
36 self.requested_outputs = float('NaN')
37
38 #set defaults
39 self.setdefaultparameters()
40
41 #}}}
42 def __repr__(self): # {{{
43
44 string=' Diagnostic solution parameters:'
45 string="%s\n%s"%(string,' Convergence criteria:')
46 string="%s\n%s"%(string,fielddisplay(self,'restol','mechanical equilibrium residual convergence criterion'))
47 string="%s\n%s"%(string,fielddisplay(self,'reltol','velocity relative convergence criterion, NaN: not applied'))
48 string="%s\n%s"%(string,fielddisplay(self,'abstol','velocity absolute convergence criterion, NaN: not applied'))
49 string="%s\n%s"%(string,fielddisplay(self,'isnewton',"0: Picard's fixed point, 1: Newton's method, 2: hybrid"))
50 string="%s\n%s"%(string,fielddisplay(self,'maxiter','maximum number of nonlinear iterations'))
51 string="%s\n%s"%(string,fielddisplay(self,'viscosity_overshoot','over-shooting constant new=new+C*(new-old)'))
52
53 string="%s\n%s"%(string,'\n boundary conditions:')
54 string="%s\n%s"%(string,fielddisplay(self,'spcvx','x-axis velocity constraint (NaN means no constraint) [m/yr]'))
55 string="%s\n%s"%(string,fielddisplay(self,'spcvy','y-axis velocity constraint (NaN means no constraint) [m/yr]'))
56 string="%s\n%s"%(string,fielddisplay(self,'spcvz','z-axis velocity constraint (NaN means no constraint) [m/yr]'))
57 string="%s\n%s"%(string,fielddisplay(self,'icefront','segments on ice front list (last column 0: Air, 1: Water, 2: Ice'))
58
59 string="%s\n%s"%(string,'\n Rift options:')
60 string="%s\n%s"%(string,fielddisplay(self,'rift_penalty_threshold','threshold for instability of mechanical constraints'))
61 string="%s\n%s"%(string,fielddisplay(self,'rift_penalty_lock','number of iterations before rift penalties are locked'))
62
63 string="%s\n%s"%(string,'\n Penalty options:')
64 string="%s\n%s"%(string,fielddisplay(self,'penalty_factor','offset used by penalties: penalty = Kmax*10^offset'))
65 string="%s\n%s"%(string,fielddisplay(self,'vertex_pairing','pairs of vertices that are penalized'))
66
67 string="%s\n%s"%(string,'\n Other:')
68 string="%s\n%s"%(string,fielddisplay(self,'shelf_dampening','use dampening for floating ice ? Only for Stokes model'))
69 string="%s\n%s"%(string,fielddisplay(self,'stokesreconditioning','multiplier for incompressibility equation. Only for Stokes model'))
70 string="%s\n%s"%(string,fielddisplay(self,'referential','local referential'))
71 string="%s\n%s"%(string,fielddisplay(self,'loadingforce','loading force applied on each point [N/m^3]'))
72 string="%s\n%s"%(string,fielddisplay(self,'requested_outputs','additional outputs requested'))
73
74 return string
75 #}}}
76 def setdefaultparameters(self): # {{{
77 #maximum of non-linear iterations.
78 self.maxiter=100
79
80 #Convergence criterion: absolute, relative and residual
81 self.restol=10**-4
82 self.reltol=0.01
83 self.abstol=10
84
85 self.stokesreconditioning=10**13
86 self.shelf_dampening=0
87
88 #Penalty factor applied kappa=max(stiffness matrix)*10^penalty_factor
89 self.penalty_factor=3
90
91 #coefficient to update the viscosity between each iteration of
92 #a diagnostic according to the following formula
93 #viscosity(n)=viscosity(n)+viscosity_overshoot(viscosity(n)-viscosity(n-1))
94 self.viscosity_overshoot=0
95
96 #Stop the iterations of rift if below a threshold
97 self.rift_penalty_threshold=0
98
99 #in some solutions, it might be needed to stop a run when only
100 #a few constraints remain unstable. For thermal computation, this
101 #parameter is often used.
102 self.rift_penalty_lock=10
103
104 return self
105 #}}}
106 def checkconsistency(self,md,solution,analyses): # {{{
107
108 #Early return
109 if DiagnosticHorizAnalysisEnum() not in analyses:
110 return md
111 #if (DiagnosticHorizAnalysisEnum() not in analyses) | (solution==TransientSolutionEnum() and not md.transient.isdiagnostic):
112 # return md
113
114 md = checkfield(md,'diagnostic.spcvx','forcing',1)
115 md = checkfield(md,'diagnostic.spcvy','forcing',1)
116 if md.mesh.dimension==3:
117 md = checkfield(md,'diagnostic.spcvz','forcing',1)
118 md = checkfield(md,'diagnostic.restol','size',[1],'>',0)
119 md = checkfield(md,'diagnostic.reltol','size',[1])
120 md = checkfield(md,'diagnostic.abstol','size',[1])
121 md = checkfield(md,'diagnostic.isnewton','numel',[1],'values',[0,1,2])
122 md = checkfield(md,'diagnostic.stokesreconditioning','size',[1],'NaN',1)
123 md = checkfield(md,'diagnostic.viscosity_overshoot','size',[1],'NaN',1)
124 if md.mesh.dimension==2:
125 md = checkfield(md,'diagnostic.icefront','size',[float('NaN'),4],'NaN',1)
126 else:
127 md = checkfield(md,'diagnostic.icefront','size',[float('NaN'),6],'NaN',1)
128 md = checkfield(md,'diagnostic.icefront[:,-1]','values',[0,1,2])
129 md = checkfield(md,'diagnostic.maxiter','size',[1],'>=',1)
130 md = checkfield(md,'diagnostic.referential','size',[md.mesh.numberofvertices,6])
131 md = checkfield(md,'diagnostic.loadingforce','size',[md.mesh.numberofvertices,3])
132 if not md.diagnostic.requested_outputs:
133 md = checkfield(md,'diagnostic.requested_outputs','size',[float('NaN'),1])
134
135 #singular solution
136# if ~any((~isnan(md.diagnostic.spcvx)+~isnan(md.diagnostic.spcvy))==2),
137 if not numpy.any(numpy.logical_and(numpy.logical_not(numpy.isnan(md.diagnostic.spcvx)),numpy.logical_not(numpy.isnan(md.diagnostic.spcvy)))):
138 md.checkmessage("model is not well posed (singular). You need at least one node with fixed velocity!")
139 #CHECK THAT EACH LINES CONTAINS ONLY NAN VALUES OR NO NAN VALUES
140# if any(sum(isnan(md.diagnostic.referential),2)~=0 & sum(isnan(md.diagnostic.referential),2)~=6),
141 if numpy.any(numpy.logical_and(numpy.sum(numpy.isnan(md.diagnostic.referential),axis=1)!=0,numpy.sum(numpy.isnan(md.diagnostic.referential),axis=1)!=6)):
142 md.checkmessage("Each line of diagnostic.referential should contain either only NaN values or no NaN values")
143 #CHECK THAT THE TWO VECTORS PROVIDED ARE ORTHOGONAL
144# if any(sum(isnan(md.diagnostic.referential),2)==0),
145 if numpy.any(numpy.sum(numpy.isnan(md.diagnostic.referential),axis=1)==0):
146 pos=[i for i,item in enumerate(numpy.sum(numpy.isnan(md.diagnostic.referential),axis=1)) if item==0]
147# numpy.inner (and numpy.dot) calculate all the dot product permutations, resulting in a full matrix multiply
148# if numpy.any(numpy.abs(numpy.inner(md.diagnostic.referential[pos,0:2],md.diagnostic.referential[pos,3:5]).diagonal())>sys.float_info.epsilon):
149# md.checkmessage("Vectors in diagnostic.referential (columns 1 to 3 and 4 to 6) must be orthogonal")
150 for item in md.diagnostic.referential[pos,:]:
151 if numpy.abs(numpy.inner(item[0:2],item[3:5]))>sys.float_info.epsilon:
152 md.checkmessage("Vectors in diagnostic.referential (columns 1 to 3 and 4 to 6) must be orthogonal")
153 #CHECK THAT NO rotation specified for FS Grounded ice at base
154# if md.mesh.dimension==3 & md.flowequation.isstokes,
155 if md.mesh.dimension==3 and md.flowequation.isstokes:
156 pos=numpy.nonzero(numpy.logical_and(md.mask.vertexongroundedice,md.mesh.vertexonbed))
157 if numpy.any(numpy.logical_not(numpy.isnan(md.diagnostic.referential[pos,:]))):
158 md.checkmessage("no referential should be specified for basal vertices of grounded ice")
159
160 return md
161 # }}}
162 def marshall(self,md,fid): # {{{
163 WriteData(fid,'object',self,'fieldname','spcvx','format','DoubleMat','mattype',1,'forcinglength',md.mesh.numberofvertices+1)
164 WriteData(fid,'object',self,'fieldname','spcvy','format','DoubleMat','mattype',1,'forcinglength',md.mesh.numberofvertices+1)
165 WriteData(fid,'object',self,'fieldname','spcvz','format','DoubleMat','mattype',1,'forcinglength',md.mesh.numberofvertices+1)
166 WriteData(fid,'object',self,'fieldname','restol','format','Double')
167 WriteData(fid,'object',self,'fieldname','reltol','format','Double')
168 WriteData(fid,'object',self,'fieldname','abstol','format','Double')
169 WriteData(fid,'object',self,'fieldname','isnewton','format','Integer')
170 WriteData(fid,'object',self,'fieldname','stokesreconditioning','format','Double')
171 WriteData(fid,'object',self,'fieldname','viscosity_overshoot','format','Double')
172 WriteData(fid,'object',self,'fieldname','maxiter','format','Integer')
173 WriteData(fid,'object',self,'fieldname','shelf_dampening','format','Integer')
174 WriteData(fid,'object',self,'fieldname','vertex_pairing','format','DoubleMat','mattype',3)
175 WriteData(fid,'object',self,'fieldname','penalty_factor','format','Double')
176 WriteData(fid,'object',self,'fieldname','rift_penalty_lock','format','Integer')
177 WriteData(fid,'object',self,'fieldname','rift_penalty_threshold','format','Integer')
178 WriteData(fid,'object',self,'fieldname','referential','format','DoubleMat','mattype',1)
179 WriteData(fid,'data',self.loadingforce[:,0],'format','DoubleMat','mattype',1,'enum',LoadingforceXEnum())
180 WriteData(fid,'data',self.loadingforce[:,1],'format','DoubleMat','mattype',1,'enum',LoadingforceYEnum())
181 WriteData(fid,'data',self.loadingforce[:,2],'format','DoubleMat','mattype',1,'enum',LoadingforceZEnum())
182 WriteData(fid,'object',self,'fieldname','requested_outputs','format','DoubleMat','mattype',3)
183
184 #marshall ice front
185 data=copy.deepcopy(self.icefront)
186 data[numpy.nonzero(data[:,-1]==0),-1]=AirEnum()
187 data[numpy.nonzero(data[:,-1]==1),-1]=WaterEnum()
188 data[numpy.nonzero(data[:,-1]==2),-1]=IceEnum()
189 WriteData(fid,'data',data,'enum',DiagnosticIcefrontEnum(),'format','DoubleMat','mattype',3)
190 # }}}
Note: See TracBrowser for help on using the repository browser.